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 -------------------------------------------------------------------------------
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 )
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 {
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();
188  fv::options& fvOptions(fv::options::New(this->mesh_));
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 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::adjointSolver
Base class for adjoint solvers.
Definition: adjointSolver.H:57
Foam::fvMatrix::H
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition: fvMatrix.C:1423
Foam::incompressiblePrimalSolver::getObjectiveFunctions
List< objective * > getObjectiveFunctions() const
Return the list of objectives assodicated with this solver.
Definition: incompressiblePrimalSolver.C:121
Foam::adjustPhi
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:37
p
volScalarField & p
Definition: createFieldRefs.H:8
incompressiblePrimalSolver.H
Foam::fvc::flux
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
rAU
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
fvOptions.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
turbulence
Info<< "Reading field U\n"<< 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\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::adjointSolver::primalSolverName
const word & primalSolverName() const
Return the primal solver name.
Definition: adjointSolver.H:152
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:55
Foam::fv::optionList::constrain
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
Definition: fvOptionListTemplates.C:314
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:103
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Foam::incompressibleVars::phi
const surfaceScalarField & phi() const
Return const reference to volume flux.
Definition: incompressibleVars.C:357
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
Foam::FatalIOError
IOerror FatalIOError
HbyA
HbyA
Definition: pcEqn.H:74
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::incompressibleVars::p
const volScalarField & p() const
Return const reference to pressure.
Definition: incompressibleVars.C:305
adjustPhi.H
For cases which do no have a pressure boundary adjust the balance of fluxes to obey continuity....
Foam::variablesSet
Base class for creating a set of variables.
Definition: variablesSet.H:49
Foam::diff
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:378
FatalIOErrorInLookup
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
tUEqn
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
phiHbyA
phiHbyA
Definition: pcEqn.H:73
adjointSolver.H
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
Foam::fvMatrix::relax
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1183
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:511
Foam::incompressiblePrimalSolver::readDict
virtual bool readDict(const dictionary &dict)
Read dict if updated.
Definition: incompressiblePrimalSolver.C:109
Foam::incompressibleVars::U
const volVectorField & U() const
Return const reference to velocity.
Definition: incompressibleVars.C:331
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::fv::options
Finite-volume options.
Definition: fvOptions.H:55
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::incompressiblePrimalSolver::useSolverNameForFields
bool useSolverNameForFields() const
Should solver name be appended to fields.
Definition: incompressiblePrimalSolver.C:147
Foam::incompressiblePrimalSolver::correctBoundaryConditions
virtual void correctBoundaryConditions()
Update boundary conditions.
Definition: incompressiblePrimalSolver.C:171
Foam::fvMatrix::A
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:1394
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::incompressibleVars::turbulence
const autoPtr< incompressible::turbulenceModel > & turbulence() const
Return const reference to the turbulence model.
Definition: incompressibleVars.C:431
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::adjointSolver::getObjectiveManager
const objectiveManager & getObjectiveManager() const
Return a const reference to the objective manager.
Definition: adjointSolver.C:142
U
U
Definition: pEqn.H:72
Foam::primalSolver::readDict
virtual bool readDict(const dictionary &dict)
Definition: primalSolver.C:86
Foam::List< Foam::objective * >
Foam::primalSolver
Base class for primal solvers.
Definition: primalSolver.H:52
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
UEqn
fvVectorMatrix & UEqn
Definition: UEqn.H:13
Foam::incompressibleVars::correctBoundaryConditions
void correctBoundaryConditions()
correct boundaryconditions for all volFields
Definition: incompressibleVars.C:506
Foam::incompressiblePrimalSolver::getIncoVars
const incompressibleVars & getIncoVars() const
Access to the incompressible variables set.
Definition: incompressiblePrimalSolver.C:154
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
Foam::objective
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition: objective.H:59
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::fvMatrix::flux
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:1534
Foam::solver::dict
virtual const dictionary & dict() const
Return the solver dictionary.
Definition: solver.C:113
Foam::incompressibleVars
Base class for solution control classes.
Definition: incompressibleVars.H:54
Foam::objectiveManager::getObjectiveFunctions
PtrList< objective > & getObjectiveFunctions()
Return reference to objective functions.
Definition: objectiveManager.C:296
Foam::incompressiblePrimalSolver::New
static autoPtr< incompressiblePrimalSolver > New(fvMesh &mesh, const word &managerType, const dictionary &dict)
Return a reference to the selected incompressible primal solver.
Definition: incompressiblePrimalSolver.C:79