78 rho_.setSize(cValues_.size() + 1, raa0_);
81 oldObjectiveValue_ = objectiveValue_;
82 oldCValues_ = cValues_;
86 const scalarField x(designVars_().getVars(), activeDesignVars_);
88 (designVars_().lowerBounds()(), activeDesignVars_);
90 (designVars_().upperBounds()(), activeDesignVars_);
93 (objectiveDerivatives_, activeDesignVars_);
98 /globalSum(
x.size())*globalSum(
mag(activeObjDerivs)*span);
101 forAll(constraintDerivatives_, i)
105 constraintDerivatives_[i], activeDesignVars_
109 /globalSum(
x.size())*globalSum(
mag(activeDerivs)*span);
116 if (dynamicRhoInitialisation_)
118 Info<<
"-----------------------------------------" <<
endl;
119 Info<<
"Solving sub problems for initializing rho" <<
endl;
120 Info<<
"-----------------------------------------" <<
endl;
131 designVars_().storeDesignVariables();
132 x0_.map(designVars_().getVars(), activeDesignVars_);
134 Info<<
"Initial pseudo update of the design variables" <<
endl;
139 designVars_().update(correction_);
145 while (!converged() && iter++ < 10)
147 Info<<
nl <<
"Dynamic rho initialisation iteration " << iter
149 designVars_().resetDesignVariables();
152 designVars_().update(correction_);
155 Info<<
"-----------------------------------------" <<
endl;
156 Info<<
"Dynamic rho initialisation converged in " << iter
157 <<
" iterations " <<
endl;
158 Info<<
"-----------------------------------------" <<
endl;
165 designVars_().resetDesignVariables();
170 rho_ *= dynamicRhoMult_;
175 rho_ =
max(rho_, scalar(1.e-6));
178 <<
"Computed r values " << rho_ <<
endl;
185 const scalarField x(designVars_().getVars(), activeDesignVars_);
186 const scalarField xMin(designVars_().lowerBoundsRef(), activeDesignVars_);
187 const scalarField xMax(designVars_().upperBoundsRef(), activeDesignVars_);
196 (
x[i] - x0_[i])*(x0_[i] - x00_[i]) < scalar(0)
197 ? asymDecr_ : asymIncr_
199 lower_[i] =
x[i] -
gamma*(x0_[i] - lower_[i]);
200 upper_[i] =
x[i] +
gamma*(upper_[i] - x0_[i]);
211 lower_ =
x - sInit_*span;
212 upper_ =
x + sInit_*span;
224 const label m(cValues_.size());
228 alpha_.setSize(m,
Zero);
231 Info<<
"Reading constraint penalty function type from dict" <<
endl;
232 c_.reset(Function1<scalar>::New(
"c", coeffsDict(
typeName)));
236 Info<<
"Setting constant penalty factor" <<
endl;
237 c_.reset(
new Function1Types::Constant<scalar>(
"c", 100));
239 d_.setSize(m, coeffsDict(
typeName).getOrDefault<scalar>(
"d", 1));
240 deltaLamda_.setSize(m,
Zero);
241 deltaY_.setSize(m,
Zero);
242 deltaS_.setSize(m,
Zero);
243 deltaMu_.setSize(m,
Zero);
245 const scalar t = mesh_.time().timeOutputValue();
246 Info<<
"Penalty multiplier (c):: " << c_->value(t) <<
endl;
254 xNew_ = 0.5*(a_ + b_);
255 ksi_ =
max(scalar(1), 1./(xNew_ - a_));
256 Eta_ =
max(scalar(1), 1./(b_ - xNew_));
259 y_.setSize(m, scalar(1));
261 s_.setSize(m, scalar(1));
263 mu_ =
max(scalar(1), 0.5*
c_->value(t));
274 x0_.map(designVars_().getVars(), activeDesignVars_);
287 designVars_().lowerBoundsRef(), activeDesignVars_
291 designVars_().upperBoundsRef(), activeDesignVars_
299 1.001*
max(derivs, scalar(0))
300 + 0.001*
max(-derivs, scalar(0))
301 + r/(upperBounds - lowerBounds)
327 designVars_().lowerBoundsRef(), activeDesignVars_
331 designVars_().upperBoundsRef(), activeDesignVars_
339 0.001*
max(derivs, scalar(0))
340 + 1.001*
max(-derivs, scalar(0))
341 + r/(upperBounds - lowerBounds)
369 constraintDerivatives_[i], activeDesignVars_
391 const scalarField activeObjDerivs(objectiveDerivatives_, activeDesignVars_);
394 forAll(constraintDerivatives_, cI)
398 constraintDerivatives_[cI], activeDesignVars_
409 const scalarField activeObjDerivs(objectiveDerivatives_, activeDesignVars_);
410 tmp<scalarField> tres(q(activeObjDerivs, rho_[0]));
412 forAll(constraintDerivatives_, cI)
416 constraintDerivatives_[cI], activeDesignVars_
491 field = defaultField;
502 const label
n(xNew_.size());
503 const label m(cValues_.size());
504 const scalarField x(designVars_().getVars(), activeDesignVars_);
507 tmp<scalarField> tpL(pLamda());
509 tmp<scalarField> tqL(qLamda());
519 scalar(2)*pL/
pow3(upper_ - xNew_)
520 + scalar(2)*qL/
pow3(xNew_ - lower_)
528 for (label i = 0; i < m; ++i)
532 constraintDerivatives_[i], activeDesignVars_
534 tmp<scalarField> tp =
p(activeDerivs, rho_[i + 1]);
535 tmp<scalarField> tq = q(activeDerivs, rho_[i + 1]);
546 pL/
sqr(upper_ - xNew_) - qL/
sqr(xNew_ - lower_)
550 const scalar t = mesh_.time().timeOutputValue();
553 c_->value(t) + d_*y_ - lamda_ - eps_/y_
555 scalar deltaZTilda = alpha0_ -
sum(lamda_*alpha_) - eps_/z_;
574 scalar zRatio(z_/zeta_);
578 for (label i = 0; i < m; ++i)
581 for (label j = 0; j < m; ++j)
583 for (label l = 0; l <
n; ++l)
586 lhs[i][j] += DxInv[l]*
G[i][l]*
G[j][l];
590 for (label l = 0; l <
n; ++l)
593 rhs[i] -= DxInv[l]*
G[i][l]*deltaXTilda[l];
604 rhs += deltaLamdaYTilda + alpha_*deltaZTilda*zRatio;
605 for (label i = 0; i < m; ++i)
608 for (label j = 0; j < m; ++j)
611 lhs[i][j] += alpha_[i]*alpha_[j]*zRatio;
614 lhs[i][i] += DLamdaY[i];
617 deltaZ_ = (
sum(alpha_*deltaLamda_) - deltaZTilda)*zRatio;
623 deltaX_ = - DxInv*(
G.Tmul(deltaLamda_) + deltaXTilda);
626 deltaY_ = DYInv*(deltaLamda_ - deltaYTilda);
629 deltaS_ = (-s_*deltaLamda_ + eps_)/lamda_ - s_;
632 deltaZeta_ = -zeta_/z_*deltaZ_ - zeta_ + eps_/z_;
635 deltaMu_ = (-mu_*deltaY_ + eps_)/y_ - mu_;
647 const label
n(xNew_.size());
648 const label m(cValues_.size());
652 for (label i = 0; i <
n; ++i)
656 xNew_[i] + step*deltaX_[i] - a_[i] - 0.01*(xNew_[i] - a_[i])
660 step = -0.99*(xNew_[i] - a_[i])/deltaX_[i];
665 -xNew_[i] - step*deltaX_[i] + b_[i] - 0.01*(b_[i] - xNew_[i])
669 step = 0.99*(b_[i] - xNew_[i])/deltaX_[i];
672 adjustStep(step, ksi_[i], deltaKsi_[i]);
673 adjustStep(step, Eta_[i], deltaEta_[i]);
676 for (label i = 0; i < m; ++i)
678 adjustStep(step, y_[i], deltaY_[i]);
679 adjustStep(step, lamda_[i], deltaLamda_[i]);
680 adjustStep(step, mu_[i], deltaMu_[i]);
681 adjustStep(step, s_[i], deltaS_[i]);
684 adjustStep(step, z_, deltaZ_);
685 adjustStep(step, zeta_, deltaZeta_);
696 Info<<
"Step before line search is " << step <<
endl;
700 scalar normResOld =
sqrt(globalSum(
magSqr(computeResiduals())));
701 scalar maxRes(GREAT);
703 for (label i = 0; i < maxLineSearchIters_ ; ++i)
706 updateSolution(step);
710 scalar normResNew =
sqrt(globalSum(
magSqr(resNew)));
713 if (normResNew < normResOld)
716 <<
"Initial residual = " << normResOld <<
", "
717 <<
"Final residual = " << normResNew <<
", "
718 <<
"No of LineSearch Iterations = " << i + 1
725 updateSolution(-step);
729 if (i == maxLineSearchIters_ - 1)
732 Info<<
"Line search could not find a step that reduced"
733 <<
" residuals while satisfying the constraints" <<
nl
734 <<
"Increasing eps to " << eps_ <<
endl;
755 if (0.99*value + step*
update < scalar(0))
757 step = -0.99*value/
update;
764 xNew_ += step*deltaX_;
767 lamda_ += step*deltaLamda_;
768 ksi_ += step*deltaKsi_;
778 const label
n(xNew_.size());
779 const label m(cValues_.size());
787 (pLamda()/
sqr(upper_ - xNew_) - qLamda()/
sqr(xNew_ - lower_));
788 for (label i = 0; i <
n; ++i)
790 res[iRes++] = dPsidx[i] - ksi_[i] + Eta_[i];
794 const scalar t = mesh_.time().timeOutputValue();
795 for (label i = 0; i < m; ++i)
797 res[iRes++] = c_->value(t) + d_[i]*y_[i] - lamda_[i] - mu_[i];
801 res[iRes++] = alpha0_ - zeta_ -
sum(lamda_*alpha_);
805 for (label i = 0; i < m; ++i)
807 res[iRes++] = gb[i] - alpha_[i]*z_ - y_[i] + s_[i];
811 for (label i = 0; i <
n; ++i)
813 res[iRes++] = ksi_[i]*(xNew_[i] - a_[i]) - eps_;
817 for (label i = 0; i <
n; ++i)
819 res[iRes++] = Eta_[i]*(b_[i] - xNew_[i]) - eps_;
823 for (label i = 0; i < m; ++i)
825 res[iRes++] = mu_[i]*y_[i] - eps_;
829 res[iRes++] = z_*zeta_ - eps_;
832 for (label i = 0; i < m; ++i)
849 || (continuousNormalisation_ && counter_ < lastNormalisationStep_)
852 scalarField activeSens(objectiveDerivatives_, activeDesignVars_);
853 Jnorm_.reset(
new scalar(
Foam::sqrt(globalSum(
sqr(activeSens)))));
856 forAll(constraintDerivatives_, cI)
859 (constraintDerivatives_[cI], activeDesignVars_);
862 Info<<
"Computed normalisation factors " <<
nl
863 <<
"J: " << Jnorm_() <<
nl
864 <<
"C: " << Cnorm_() <<
endl;
868 objectiveValue_ /= Jnorm_() + SMALL;
869 objectiveDerivatives_ /= Jnorm_() + SMALL;
872 cValues_ *= cw_/(Cnorm_() + SMALL);
873 forAll(constraintDerivatives_, cI)
875 constraintDerivatives_[cI] *= cw_/(Cnorm_()[cI] + SMALL);
888 const label nConstraints,
894 x0_(getOrDefaultScalarField(
"x0", activeDesignVars_.size())),
895 x00_(getOrDefaultScalarField(
"x00", activeDesignVars_.size())),
896 xNew_(activeDesignVars_.size()),
897 oldObjectiveValue_(getOrDefault<scalar>(
"J0",
Zero)),
898 oldCValues_(getOrDefaultScalarField(
"C0", nConstraints_)),
900 z_(coeffsDict(
type).getOrDefault<scalar>(
"z", 1.)),
901 alpha0_(coeffsDict(
type).getOrDefault<scalar>(
"alpha0", 1.)),
906 lower_(xNew_.size(), -GREAT),
907 upper_(xNew_.size(), GREAT),
911 boundRho_(coeffsDict(
type).getOrDefault<bool>(
"boundRho", false)),
912 correctDVs_(coeffsDict(
type).getOrDefault<bool>(
"correct", true)),
917 zeta_(coeffsDict(
type).getOrDefault<scalar>(
"zeta", 1.)),
921 deltaX_(xNew_.size(),
Zero),
926 deltaEta_(xNew_.size(),
Zero),
927 deltaKsi_(xNew_.size(),
Zero),
929 maxNewtonIters_(coeffsDict(
type).getOrDefault<label>(
"maxIters", 6000)),
932 coeffsDict(
type).getOrDefault<label>(
"maxLineSearchIters", 10)
934 initializeEverySubproblem_
935 (coeffsDict(
type).getOrDefault<bool>(
"initializeEverySubproblem", false)),
936 asymDecr_(coeffsDict(
type).getOrDefault<scalar>(
"asymptoteDecrease", 0.7)),
937 asymIncr_(coeffsDict(
type).getOrDefault<scalar>(
"asymptoteIncrease", 1.2)),
938 sInit_(coeffsDict(
type).getOrDefault<scalar>(
"sInit", 0.5)),
939 move_(coeffsDict(
type).getOrDefault<scalar>(
"move", 0.5)),
940 raa0_(coeffsDict(
type).getOrDefault<scalar>(
"raa0", 1.e-05)),
941 maxInitRhoMult_(coeffsDict(
type).getOrDefault<scalar>(
"maxInitRhoMult", 0.1)),
942 maxRhoMult_(coeffsDict(
type).getOrDefault<scalar>(
"maxRhoMult", 10)),
943 variableRho_(coeffsDict(
type).getOrDefault<bool>(
"variableRho", false)),
944 dynamicRhoInitialisation_
945 (coeffsDict(
type).getOrDefault<bool>(
"dynamicRhoInitialisation", false)),
947 (coeffsDict(
type).getOrDefault<scalar>(
"dynamicRhoMult", 0.1)),
948 normalise_(coeffsDict(
type).getOrDefault<bool>(
"normalise", false)),
949 continuousNormalisation_
950 (coeffsDict(
type).getOrDefault<bool>(
"continuousNormalisation", false)),
953 cw_(coeffsDict(
type).getOrDefault<scalar>(
"constraintWeight", 1)),
954 lastNormalisationStep_
956 coeffsDict(
type).getOrDefault<label>(
"lastNormalisationStep", 20)
960 if (!designVars().lowerBounds() || !designVars().upperBounds())
963 <<
"Lower and upper design variable bounds must be set for MMA"
1022 if (initializeEverySubproblem_)
1027 scalar resMax(
gMax(
mag(computeResiduals())));
1032 <<
"Newton iter " << iter <<
nl <<
endl;
1035 if (resMax < 0.9*eps_)
1041 computeNewtonDirection();
1048 <<
"max residual = " << resMax <<
", "
1049 <<
"eps = " << eps_ <<
nl <<
endl;
1051 mesh_.time().printExecutionTime(
Info);
1052 }
while (iter++ < maxNewtonIters_ && (eps_ > 1.e-07 || resMax > 0.9*eps_));
1054 Info<<
"Solved the MMA Newton problem in " << iter <<
" iterations "
1058 <<
"Constraint penalization variables (y) " << y_ <<
endl;
1061 const scalarField& oldVars = designVars_().getVars();
1062 forAll(activeDesignVars_, avI)
1074 label m(cValues_.size());
1081 const scalarField activeObjDerivs(objectiveDerivatives_, activeDesignVars_);
1082 vals[0] = objectiveValue_;
1083 approx[0] = fTilda(x0_, xNew_, oldObjectiveValue_, activeObjDerivs, rho_[0]);
1086 forAll(constraintDerivatives_, i)
1090 constraintDerivatives_[i], activeDesignVars_
1092 vals[i + 1] = cValues_[i];
1112 const scalarField xMin(designVars_().lowerBounds()(), activeDesignVars_);
1113 const scalarField xMax(designVars_().upperBounds()(), activeDesignVars_);
1115 scalar d = globalSum
1117 (upper_ - lower_)*
sqr(xNew_ - x0_)
1118 /(upper_ - xNew_)/(xNew_ - lower_)/(
xMax -
xMin)
1124 const scalar
delta = (vals[i] - approx[i])/d;
1127 rho_[i] =
min(1.1*(rho_[i] +
delta), maxRhoMult_*rho_[i]);
1131 <<
"Updated rho values to " << rho_ <<
endl;
1149 updateValuesAndApproximations();
1153 bool isConverged(
true);
1157 <<
nl <<
"MMA, objective/constraint " << i <<
nl
1158 <<
"Approximation " << approx[i]
1159 <<
" | old value " << vals[i] <<
nl <<
endl;
1160 isConverged = isConverged && (approx[i] - vals[i] > 0.);
1169 x0_.writeEntry(
"x0",
os);
1170 x00_.writeEntry(
"x00",
os);
1171 lower_.writeEntry(
"lower",
os);
1172 upper_.writeEntry(
"upper",
os);
1173 os.writeEntry(
"J0", oldObjectiveValue_);
1174 oldCValues_.writeEntry(
"C0",
os);
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
label size() const noexcept
The number of elements in list.
Templated function that returns a constant value.
void setSize(label n)
Alias for resize().
Update design variables using the Method of Moving Assymptotes. Can handle inequality constraints.
tmp< scalarField > gConstr(const scalarField &vars)
g of all constraint functions
scalarField x0_
The previous design variables. Used to compute new bounds.
bool normalise_
Perform the normalisation.
scalar lineSearch()
Perform line search and return max residual corresponding to the updated solution.
void updateSolution(const scalar step)
Update the current solution using the known direction and the given step length.
void solveSubproblem()
Solve subproblem using a Newton optimiser.
void computeCorrection()
Compute design variables correction.
scalar maxInitRhoMult_
Multiplier of the mean derivative used for the initialisation of rho.
bool converged()
Checks convergence of GCMMA.
void initialize()
Allocate fields related to constraints.
scalar zeta_
Lagrange multiplier for the z constraint.
PtrList< scalarField > valsAndApproxs_
Objective/Constraint values and approximations in the new point.
scalarField mu_
Lagrange multipliers for the y constraints.
autoPtr< Function1< scalar > > c_
y multipliers in the objective function
scalarField alpha_
Terms multyplying z in the constraint functions.
scalar oldObjectiveValue_
Old objective value.
void adjustStep(scalar &step, const scalar value, const scalar update)
Adjust step to satisfy cireteria.
scalar maxRhoMult_
Multiplier for the max. rho value during its update.
virtual bool writeData(Ostream &os) const
Write useful quantities to files.
scalarField xNew_
The set of design variables being updated during the subproblem.
scalar sInit_
Used to update the assymptotes in the first optimisation cycle.
tmp< scalarField > p(const scalarField &derivs, const scalar r, const scalarField &x)
p and q functions, used to approximate the objective and contraints
tmp< scalarField > q(const scalarField &derivs, const scalar r, const scalarField &x)
bool initializeEverySubproblem_
Initialize every subproblem with x = 0.5*(a + b) and reset Lagrange multipliers.
scalarField x00_
The twice previous design variables. Used to compute new bounds.
void updateValuesAndApproximations()
Compute objective/constraint values and their approximations.
void updateBounds()
Update the bounds associated with the design variables.
bool correctDVs_
Correct the design variables.
void computeNewtonDirection()
Compute direction for the Newton problem.
scalar move_
Movement of the unatainable upper and lower bounds.
scalarField ksi_
Lagrange multipliers for the lower limits constraints.
tmp< scalarField > b()
The rhs of the contraint approximations.
void initializeRho()
Initialize rho constants in the (p, q) functions.
tmp< scalarField > getOrDefaultScalarField(const word &name, const label size, const scalar value=Zero)
Read in scalarField from self (potentially in binary), or allocate field with the size of the active ...
scalar z_
Second term in the approximation function.
scalar fTilda(const scalarField &xInit, const scalarField &x, const scalar f, const scalarField &dfdx, const scalar rho)
Computation of the approximate function.
bool dynamicRhoInitialisation_
Change rho constants in each iteration?
scalarField b_
Upper design variables constraints.
bool continuousNormalisation_
Perform the normalisation in each optimisation cycle.
void zeroUpdates()
Zero the updates computed in the previous optimisation cycle.
tmp< scalarField > pLamda()
p and q functions, using the dual lamda
scalar eps_
Infinitesimal quantity.
scalarField d_
y^2 multipliers in the objective function
scalarField lamda_
Constraint Lagrange multipliers.
void updateRho()
Update rho value if needed.
label maxLineSearchIters_
Maxmimum number of line search iterations for each iteration of the subproblem.
scalarField upper_
Upper design variables asymptotes.
scalarField oldCValues_
Old constraint values.
scalar asymIncr_
Upper assymprote increase rate.
scalarField s_
Slack variables for the inequality constraints.
label lastNormalisationStep_
Constaint weight after the normalisation.
scalar cw_
Constaint weight after the normalisation.
scalarField a_
Lower design variables constraints.
scalar raa0_
Constant in p, q functions.
autoPtr< scalar > Jnorm_
Normalisation factor for the objective.
scalar asymDecr_
Lower assymprote reduction rate.
void setOrDefaultScalarField(scalarField &field, const word &name, const label size, const scalarField &defaultField)
Read in scalarField from self (potentially in binary), or allocate field with the size of the active ...
const scalarField & getRho() const
Get rho value if needed.
scalarField lower_
Lower design variables asymptotes.
bool variableRho_
Change rho constants in each iteration?
tmp< scalarField > qLamda()
void normalise()
Normalise the objective and constraints if necessary.
scalarField y_
y in the constraints functions
void storeOldValues()
Store old objcective and constraint values.
bool boundRho_
Bound the rho value with an upper value?
scalarField rho_
Constants in the (p,q) functions.
scalar dynamicRhoMult_
The rho values obtained by the dynamic rho initialisation might be too conservative....
tmp< scalarField > computeResiduals()
Compute and return residuals based on the current solution.
void updateSizes()
Update sizes of fields related to the active design variables.
tmp< scalarField > Cnorm_
Normalisation factors for the constraints.
scalar alpha0_
Term multiplying z in the objective function.
void setVariableRho(bool varRho=true)
Set variable rho.
label maxNewtonIters_
Maxmimum number of Newton iterations for the subproblem.
scalarField Eta_
Lagrange multipliers for the upper limits constraints.
const PtrList< scalarField > & getValuesAndApproximations() const
Get objective/constraint values and their approximations.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
static void listReduce(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce list elements (list must be equal size on all ranks), applying bop to each list element.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
void size(const label n)
Older name for setAddressableSize.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Abstract base class for optimisation methods supporting constraints. Does not add functionality to up...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
Mesh data needed to do the Finite Volume discretisation.
Abstract base class for line search methods.
A class for managing temporary objects.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Abstract base class for optimisation methods.
scalarField objectiveDerivatives_
Derivatives of the objective functions.
label nConstraints_
Number of constraints.
scalarField correction_
Design variables correction.
bool globalSum_
Whether to use gSum or sum in the inner products.
dictionary coeffsDict(const word &type) const
Return optional dictionary with parameters specific to each method.
scalar globalSum(const scalarField &field)
Compute either global or local sum, based on globalSum flag.
scalarField cValues_
Constraint values.
scalar objectiveValue_
Objective value.
PtrList< scalarField > constraintDerivatives_
Derivatives of the constraints.
virtual bool writeData(Ostream &os) const
Write continuation data under uniform.
const labelList & activeDesignVars_
Map to active design variables.
label counter_
Optimisation cycle count.
label nConstraints() const
Get the number of constraints.
autoPtr< designVariables > & designVars_
Design variables.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
OBJstream os(runTime.globalPath()/outputName)
#define DebugInfo
Report an information message using Foam::Info.
const dimensionedScalar G
Newtonian constant of gravitation.
Namespace for handling debugging switches.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
RectangularMatrix< scalar > scalarRectangularMatrix
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
errorManip< error > abort(error &err)
SquareMatrix< scalar > scalarSquareMatrix
void rhs(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Type gMax(const FieldField< Field, Type > &f)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.