54Foam::incompressiblePrimalSolver::incompressiblePrimalSolver
57 const word& managerType,
59 const word& solverName
65 dict.subOrEmptyDict(
"fieldReconstruction").
66 getOrDefault<scalar>(
"tolerance", 5.e-5)
68 phiReconstructionIters_
70 dict.subOrEmptyDict(
"fieldReconstruction").
71 getOrDefault<label>(
"iters", 10)
82 const word& managerType,
84 const word& solverName
88 auto* ctorPtr = dictionaryConstructorTable(solverType);
95 "incompressiblePrimalSolver",
97 *dictionaryConstructorTablePtr_
102 autoPtr<incompressiblePrimalSolver>
125 DynamicList<objective*> objectives;
127 for (
auto& adjoint : mesh_.sorted<adjointSolver>())
129 if (adjoint.primalSolverName() == solverName_)
131 PtrList<objective>& managerObjectives =
132 adjoint.getObjectiveManager().getObjectiveFunctions();
134 for (objective& obj : managerObjectives)
136 objectives.push_back(&obj);
182 scalar contError(GREAT),
diff(GREAT);
183 for (label iter = 0; iter < phiReconstructionIters_; ++iter)
185 Info<<
"phi correction iteration " << iter <<
endl;
220 scalar contErrorNew =
221 mesh_.time().deltaTValue()*
224 Info<<
"sum local = " << contErrorNew <<
endl;
225 diff =
mag(contErrorNew - contError)/contError;
226 contError = contErrorNew;
228 if (
diff < phiReconstructionTol_)
break;
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
For cases which do no have a pressure boundary adjust the balance of fluxes to obey continuity....
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Base class for adjoint solvers.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Find and return a sub-dictionary as a copy, otherwise return an empty dictionary.
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...
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Mesh data needed to do the Finite Volume discretisation.
Finite-volume options, which is an IOdictionary of values and a fv::optionList.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
Base class for primal incompressible solvers.
static autoPtr< incompressiblePrimalSolver > New(fvMesh &mesh, const word &managerType, const dictionary &dict, const word &solverName)
Return a reference to the selected incompressible primal solver.
const incompressibleVars & getIncoVars() const
Access to the incompressible variables set.
virtual bool readDict(const dictionary &dict)
Read dict if updated.
label phiReconstructionIters_
Max iterations for reconstructing phi from U and p.
scalar phiReconstructionTol_
Convergence criterion for reconstructing phi from U and p.
UPtrList< objective > getObjectiveFunctions() const
Return the list of objectives assodicated with this solver.
virtual void correctBoundaryConditions()
Update boundary conditions.
Base class for solution control classes.
const volVectorField & U() const
Return const reference to velocity.
const volScalarField & p() const
Return const reference to pressure.
const surfaceScalarField & phi() const
Return const reference to volume flux.
const autoPtr< incompressible::turbulenceModel > & turbulence() const
Return const reference to the turbulence model.
void correctBoundaryConditions()
correct boundaryconditions for all volFields
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Base class for primal solvers.
virtual bool readDict(const dictionary &dict)
autoPtr< variablesSet > vars_
Base variableSet pointer.
const word & managerType() const
Return the manager type.
const dictionary & dict() const
Return the solver dictionary.
const word solverName_
Solver name.
const fvMesh & mesh() const
Return the solver mesh.
const word & solverName() const
Return the solver name.
fvMesh & mesh_
Reference to the mesh database.
A class for managing temporary objects.
Base class for creating a set of variables.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
tmp< volScalarField > rAU
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
GeometricField< vector, fvPatchField, volMesh > volVectorField
fvMatrix< scalar > fvScalarMatrix
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
Type weightedAverage(const UList< scalar > &weights, const UList< Type > &fld)
The local weighted average of a field, using the mag() of the weights.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fvMatrix< vector > fvVectorMatrix
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.