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-2020 PCOpt/NTUA
9 Copyright (C) 2013-2020 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"
35
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
44 (
48 );
49}
50
51
52// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53
55(
56 fvMesh& mesh,
57 const word& managerType,
58 const dictionary& dict
59)
60:
61 primalSolver(mesh, managerType, dict),
62 phiReconstructionTol_
63 (
64 dict.subOrEmptyDict("fieldReconstruction").
65 getOrDefault<scalar>("tolerance", 5.e-5)
66 ),
67 phiReconstructionIters_
68 (
69 dict.subOrEmptyDict("fieldReconstruction").
70 getOrDefault<label>("iters", 10)
71 )
72{}
73
74
75// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
76
79(
80 fvMesh& mesh,
81 const word& managerType,
82 const dictionary& dict
83)
84{
85 const word solverType(dict.get<word>("solver"));
86 auto* ctorPtr = dictionaryConstructorTable(solverType);
87
88 if (!ctorPtr)
89 {
91 (
92 dict,
93 "incompressiblePrimalSolver",
94 solverType,
95 *dictionaryConstructorTablePtr_
96 ) << exit(FatalIOError);
97 }
98
99 return
101 (
102 ctorPtr(mesh, managerType, dict)
103 );
104}
105
106
107// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
108
110{
111 if (primalSolver::readDict(dict))
112 {
113 return true;
114 }
115
116 return false;
117}
118
119
122{
123 DynamicList<objective*> objectives(10);
124
125 auto adjointSolvers = mesh_.lookupClass<adjointSolver>();
126
127 for (adjointSolver* adjointPtr : adjointSolvers)
128 {
129 adjointSolver& adjoint = *adjointPtr;
130
131 if (adjoint.primalSolverName() == solverName_)
132 {
133 PtrList<objective>& managerObjectives =
135
136 for (objective& obj : managerObjectives)
137 {
138 objectives.append(&obj);
139 }
140 }
141 }
142
143 return objectives;
144}
145
146
148{
149 return vars_().useSolverNameForFields();
150}
151
152
155{
156 const incompressibleVars& incoVars =
157 refCast<incompressibleVars>(const_cast<variablesSet&>(vars_()));
158 return incoVars;
159}
160
161
164{
165 incompressibleVars& incoVars =
166 refCast<incompressibleVars>(const_cast<variablesSet&>(vars_()));
167 return incoVars;
168}
169
170
172{
173 incompressibleVars& vars = getIncoVars();
174 // Update boundary conditions for all primal volFields,
175 // including averaged ones, if present
177
178 // phi cannot be updated through correctBoundayrConditions.
179 // Re-compute based on the Rhie-Chow interpolation scheme.
180 // This is a non-linear process
181 // (phi depends on UEqn().A() which depends on phi)
182 // so iterations are required
183
184 volScalarField& p = vars.p();
185 volVectorField& U = vars.U();
186 surfaceScalarField& phi = vars.phi();
189
190 scalar contError(GREAT), diff(GREAT);
191 for (label iter = 0; iter < phiReconstructionIters_; ++iter)
192 {
193 Info<< "phi correction iteration " << iter << endl;
194
195 // Form momentum equations to grab A
196 //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
198 (
199 fvm::div(phi, U)
200 + turbulence->divDevReff(U)
201 ==
202 fvOptions(U)
203 );
204 fvVectorMatrix& UEqn = tUEqn.ref();
205 UEqn.relax();
207
208 // Pressure equation will give the Rhie-Chow correction
209 //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
210 volScalarField rAU(1.0/UEqn.A());
211 volVectorField HbyA("HbyA", U);
212 HbyA = rAU*UEqn.H();
213 tUEqn.clear();
214
216 adjustPhi(phiHbyA, U, p);
217
218 //fvOptions.makeRelative(phiHbyA);
219
220 fvScalarMatrix pEqn
221 (
223 );
224 phi = phiHbyA - pEqn.flux();
225
226 // Check convergence
227 // ~~~~~~~~~~~~~~~~~
228 scalar contErrorNew =
229 mesh_.time().deltaTValue()*
230 mag(fvc::div(phi)())().weightedAverage(mesh_.V()).value();
231
232 Info<< "sum local = " << contErrorNew << endl;
233 diff = mag(contErrorNew - contError)/contError;
234 contError = contErrorNew;
235
236 if (diff < phiReconstructionTol_) break;
237
238 Info<< endl;
239 }
240}
241
242
243// ************************************************************************* //
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
surfaceScalarField & phi
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Base class for adjoint solvers.
Definition: adjointSolver.H:60
const objectiveManager & getObjectiveManager() const
Return a const reference to the objective manager.
const word & primalSolverName() const
Return the primal solver name.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1092
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:1443
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:1303
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition: fvMatrix.C:1332
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
Finite-volume options.
Definition: fvOptions.H:59
Base class for primal incompressible solvers.
List< objective * > getObjectiveFunctions() const
Return the list of objectives assodicated with this solver.
const incompressibleVars & getIncoVars() const
Access to the incompressible variables set.
virtual bool readDict(const dictionary &dict)
Read dict if updated.
bool useSolverNameForFields() const
Should solver name be appended to fields.
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
PtrList< objective > & getObjectiveFunctions()
Return reference to objective functions.
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition: objective.H:62
Base class for primal solvers.
Definition: primalSolver.H:55
A class for managing temporary objects.
Definition: tmp.H:65
Base class for creating a set of variables.
Definition: variablesSet.H:50
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
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:478
compressible::turbulenceModel & turbulence
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:34
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
Namespace for OpenFOAM.
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:37
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:378
IOerror FatalIOError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict