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 -------------------------------------------------------------------------------
12 License
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 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(incompressiblePrimalSolver, 0);
42  defineRunTimeSelectionTable(incompressiblePrimalSolver, dictionary);
44  (
45  primalSolver,
46  incompressiblePrimalSolver,
47  primalSolver
48  );
49 }
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
54 Foam::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;
116  }
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  }
139  }
140 
141  return UPtrList<objective>(objectives);
142 }
143 
144 
147 {
148  const incompressibleVars& incoVars =
149  refCast<incompressibleVars>(const_cast<variablesSet&>(vars_()));
150  return incoVars;
151 }
152 
153 
156 {
157  incompressibleVars& incoVars =
158  refCast<incompressibleVars>(const_cast<variablesSet&>(vars_()));
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();
180  fv::options& fvOptions(fv::options::New(this->mesh_));
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();
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 // ************************************************************************* //
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:373
dictionary dict
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
static autoPtr< incompressiblePrimalSolver > New(fvMesh &mesh, const word &managerType, const dictionary &dict, const word &solverName)
Return a reference to the selected incompressible primal solver.
Base class for primal solvers.
Definition: primalSolver.H:46
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
const fvMesh & mesh() const
Return the solver mesh.
Definition: solverI.H:24
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition: fvMatrix.C:1332
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
Finite-volume options.
Definition: fvOptions.H:51
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
HbyA
Definition: pcEqn.H:74
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:30
Macros for easy insertion into run-time selection tables.
fv::options & fvOptions
phiHbyA
Definition: pcEqn.H:73
const surfaceScalarField & phi() const
Return const reference to volume flux.
Base class for solution control classes.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
const volScalarField & p() const
Return const reference to pressure.
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
virtual void correctBoundaryConditions()
Update boundary conditions.
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 volVectorField & U() const
Return const reference to velocity.
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
defineTypeNameAndDebug(combustionModel, 0)
UPtrList< objective > getObjectiveFunctions() const
Return the list of objectives assodicated with this solver.
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1101
Base class for creating a set of variables.
Definition: variablesSet.H:43
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:1313
virtual bool readDict(const dictionary &dict)
Read dict if updated.
U
Definition: pEqn.H:72
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
const word & solverName() const
Return the solver name.
Definition: solverI.H:30
fvVectorMatrix & UEqn
Definition: UEqn.H:13
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void correctBoundaryConditions()
correct boundaryconditions for all volFields
For cases which do no have a pressure boundary adjust the balance of fluxes to obey continuity...
messageStream Info
Information stream (stdout output on master, null elsewhere)
const incompressibleVars & getIncoVars() const
Access to the incompressible variables set.
const autoPtr< incompressible::turbulenceModel > & turbulence() const
Return const reference to the turbulence model.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:96
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:635
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
virtual bool readDict(const dictionary &dict)
Definition: primalSolver.C:82