Loading...
Searching...
No Matches
incompressiblePrimalSolver.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2007-2023 PCOpt/NTUA
9 Copyright (C) 2013-2023 FOSS GP
10 Copyright (C) 2019-2021 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
31#include "adjustPhi.H"
32#include "adjointSolver.H"
33#include "fvOptions.H"
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
44 (
47 primalSolver
48 );
49}
50
51
52// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53
54Foam::incompressiblePrimalSolver::incompressiblePrimalSolver
55(
56 fvMesh& mesh,
57 const word& managerType,
58 const dictionary& dict,
59 const word& solverName
60)
61:
62 primalSolver(mesh, managerType, dict, solverName),
63 phiReconstructionTol_
64 (
65 dict.subOrEmptyDict("fieldReconstruction").
66 getOrDefault<scalar>("tolerance", 5.e-5)
67 ),
68 phiReconstructionIters_
69 (
70 dict.subOrEmptyDict("fieldReconstruction").
71 getOrDefault<label>("iters", 10)
72 )
73{}
74
75
76// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
77
80(
81 fvMesh& mesh,
82 const word& managerType,
83 const dictionary& dict,
84 const word& solverName
85)
86{
87 const word solverType(dict.get<word>("solver"));
88 auto* ctorPtr = dictionaryConstructorTable(solverType);
89
90 if (!ctorPtr)
91 {
93 (
94 dict,
95 "incompressiblePrimalSolver",
96 solverType,
97 *dictionaryConstructorTablePtr_
98 ) << exit(FatalIOError);
99 }
100
101 return
102 autoPtr<incompressiblePrimalSolver>
103 (
105 );
106}
107
108
109// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
110
112{
114 {
115 return true;
117
118 return false;
119}
120
121
124{
125 DynamicList<objective*> objectives;
126
127 for (auto& adjoint : mesh_.sorted<adjointSolver>())
128 {
129 if (adjoint.primalSolverName() == solverName_)
130 {
131 PtrList<objective>& managerObjectives =
132 adjoint.getObjectiveManager().getObjectiveFunctions();
133
134 for (objective& obj : managerObjectives)
135 {
136 objectives.push_back(&obj);
137 }
138 }
140
142}
143
144
147{
148 const incompressibleVars& incoVars =
150 return incoVars;
151}
152
153
157 incompressibleVars& incoVars =
159 return incoVars;
160}
161
162
164{
165 incompressibleVars& vars = getIncoVars();
166 // Update boundary conditions for all primal volFields,
167 // including averaged ones, if present
169
170 // phi cannot be updated through correctBoundaryConditions.
171 // Re-compute based on the Rhie-Chow interpolation scheme.
172 // This is a non-linear process
173 // (phi depends on UEqn().A() which depends on phi)
174 // so iterations are required
175
176 volScalarField& p = vars.p();
177 volVectorField& U = vars.U();
178 surfaceScalarField& phi = vars.phi();
181
182 scalar contError(GREAT), diff(GREAT);
183 for (label iter = 0; iter < phiReconstructionIters_; ++iter)
184 {
185 Info<< "phi correction iteration " << iter << endl;
186
187 // Form momentum equations to grab A
188 //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
190 (
191 fvm::div(phi, U)
192 + turbulence->divDevReff(U)
193 ==
194 fvOptions(U)
195 );
196 fvVectorMatrix& UEqn = tUEqn.ref();
197 UEqn.relax();
198 fvOptions.constrain(UEqn);
199
200 // Pressure equation will give the Rhie-Chow correction
201 //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
202 volScalarField rAU(1.0/UEqn.A());
203 volVectorField HbyA("HbyA", U);
204 HbyA = rAU*UEqn.H();
205 tUEqn.clear();
206
208 adjustPhi(phiHbyA, U, p);
209
210 //fvOptions.makeRelative(phiHbyA);
211
212 fvScalarMatrix pEqn
213 (
215 );
216 phi = phiHbyA - pEqn.flux();
217
218 // Check convergence
219 // ~~~~~~~~~~~~~~~~~
220 scalar contErrorNew =
221 mesh_.time().deltaTValue()*
222 mag(fvc::div(phi)())().weightedAverage(mesh_.V()).value();
223
224 Info<< "sum local = " << contErrorNew << endl;
225 diff = mag(contErrorNew - contError)/contError;
226 contError = contErrorNew;
227
228 if (diff < phiReconstructionTol_) break;
229
230 Info<< endl;
231 }
232}
233
234
235// ************************************************************************* //
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....
fv::options & fvOptions
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition UPtrList.H:101
Base class for adjoint solvers.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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.
Definition dictionary.C:521
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.
Definition fvMatrix.C:1415
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Finite-volume options, which is an IOdictionary of values and a fv::optionList.
Definition fvOptions.H:69
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
Definition fvOptions.C:116
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...
Definition objective.H:58
Base class for primal solvers.
virtual bool readDict(const dictionary &dict)
autoPtr< variablesSet > vars_
Base variableSet pointer.
Definition solver.H:95
const word & managerType() const
Return the manager type.
Definition solverI.H:72
const dictionary & dict() const
Return the solver dictionary.
Definition solverI.H:54
const word solverName_
Solver name.
Definition solver.H:71
const fvMesh & mesh() const
Return the solver mesh.
Definition solverI.H:24
const word & solverName() const
Return the solver name.
Definition solverI.H:30
fvMesh & mesh_
Reference to the mesh database.
Definition solver.H:56
A class for managing temporary objects.
Definition tmp.H:75
Base class for creating a set of variables.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
volScalarField & p
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
fvVectorMatrix & UEqn
Definition UEqn.H:13
phiHbyA
Definition pcEqn.H:73
HbyA
Definition pcEqn.H:74
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
tmp< volScalarField > rAU
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition fvcFlux.C:27
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)
Definition fvmDiv.C:41
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition adjustPhi.C:30
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.
Definition Ostream.H:519
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.
Definition triad.C:373
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)
Definition errorManip.H:125
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict