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-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 FOSS GP
10  Copyright (C) 2019 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"
34 
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(incompressiblePrimalSolver, 0);
41  defineRunTimeSelectionTable(incompressiblePrimalSolver, dictionary);
43  (
44  primalSolver,
45  incompressiblePrimalSolver,
46  primalSolver
47  );
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
53 Foam::incompressiblePrimalSolver::incompressiblePrimalSolver
54 (
55  fvMesh& mesh,
56  const word& managerType,
57  const dictionary& dict
58 )
59 :
60  primalSolver(mesh, managerType, dict),
61  phiReconstructionTol_
62  (
63  dict.subOrEmptyDict("fieldReconstruction").
64  lookupOrDefault<scalar>("tolerance", scalar(5.e-5))
65  ),
66  phiReconstructionIters_
67  (
68  dict.subOrEmptyDict("fieldReconstruction").
69  lookupOrDefault<label>("iters", label(10))
70  ),
71  fvOptions_(nullptr)
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 cstrIter = dictionaryConstructorTablePtr_->cfind(solverType);
87 
88  if (!cstrIter.found())
89  {
91  (
92  dict,
93  "incompressiblePrimalSolver",
94  solverType,
95  *dictionaryConstructorTablePtr_
96  ) << exit(FatalIOError);
97  }
98 
99  return
101  (
102  cstrIter()(mesh, managerType, dict)
103  );
104 }
105 
106 
107 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
108 
110 {
112  {
113  fvOptions_().read(dict.subOrEmptyDict("fvOptions"));
114 
115  return true;
116  }
117 
118  return false;
119 }
120 
121 
124 {
125  DynamicList<objective*> objectives(10);
126 
127  auto adjointSolvers = mesh_.lookupClass<adjointSolver>();
128 
129  for (adjointSolver* adjointPtr : adjointSolvers)
130  {
131  adjointSolver& adjoint = *adjointPtr;
132 
133  if (adjoint.active() && adjoint.primalSolverName() == solverName_)
134  {
135  PtrList<objective>& managerObjectives =
137 
138  for (objective& obj : managerObjectives)
139  {
140  objectives.append(&obj);
141  }
142  }
143  }
144 
145  return objectives;
146 }
147 
148 
150 {
151  return vars_().useSolverNameForFields();
152 }
153 
154 
157 {
158  const incompressibleVars& incoVars =
159  refCast<incompressibleVars>(const_cast<variablesSet&>(vars_()));
160  return incoVars;
161 }
162 
163 
166 {
167  incompressibleVars& incoVars =
168  refCast<incompressibleVars>(const_cast<variablesSet&>(vars_()));
169  return incoVars;
170 }
171 
172 
174 {
175  incompressibleVars& vars = getIncoVars();
176  // Update boundary conditions for all primal volFields,
177  // including averaged ones, if present
179 
180  // phi cannot be updated through correctBoundayrConditions.
181  // Re-compute based on the Rhie-Chow interpolation scheme.
182  // This is a non-linear process
183  // (phi depends on UEqn().A() which depends on phi)
184  // so iterations are required
185 
186  volScalarField& p = vars.p();
187  volVectorField& U = vars.U();
188  surfaceScalarField& phi = vars.phi();
190 
191  scalar contError(GREAT), diff(GREAT);
192  for (label iter = 0; iter < phiReconstructionIters_; ++iter)
193  {
194  Info<< "phi correction iteration " << iter << endl;
195 
196  // Form momentum equations to grab A
197  //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
199  (
200  fvm::div(phi, U)
201  + turbulence->divDevReff(U)
202  ==
203  fvOptions_()(U)
204  );
205  fvVectorMatrix& UEqn = tUEqn.ref();
206  UEqn.relax();
207  fvOptions_().constrain(UEqn);
208 
209  // Pressure equation will give the Rhie-Chow correction
210  //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
211  volScalarField rAU(1.0/UEqn.A());
212  volVectorField HbyA("HbyA", U);
213  HbyA = rAU*UEqn.H();
214  tUEqn.clear();
215 
217  adjustPhi(phiHbyA, U, p);
218 
219  //fvOptions_().makeRelative(phiHbyA);
220 
221  fvScalarMatrix pEqn
222  (
224  );
225  phi = phiHbyA - pEqn.flux();
226 
227  // Check convergence
228  // ~~~~~~~~~~~~~~~~~
229  scalar contErrorNew =
230  mesh_.time().deltaTValue()*
231  mag(fvc::div(phi)())().weightedAverage(mesh_.V()).value();
232 
233  Info<< "sum local = " << contErrorNew << endl;
234  diff = mag(contErrorNew - contError)/contError;
235  contError = contErrorNew;
236 
237  if (diff < phiReconstructionTol_) break;
238 
239  Info<< endl;
240  }
241 }
242 
243 
244 // ************************************************************************* //
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:817
Foam::incompressiblePrimalSolver::getObjectiveFunctions
List< objective * > getObjectiveFunctions() const
Return the list of objectives assodicated with this solver.
Definition: incompressiblePrimalSolver.C:123
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::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
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:59
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:57
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Foam::incompressiblePrimalSolver::fvOptions_
autoPtr< fv::optionList > fvOptions_
optionList for source terms addition
Definition: incompressiblePrimalSolver.H:80
Foam::solver::active
virtual bool active()
Return state of solver.
Definition: solver.C:107
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:337
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:380
tUEqn
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
phiHbyA
phiHbyA
Definition: pcEqn.H:73
adjointSolver.H
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
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:574
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:472
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:65
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
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:121
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:149
Foam::incompressiblePrimalSolver::correctBoundaryConditions
virtual void correctBoundaryConditions()
Update boundary conditions.
Definition: incompressiblePrimalSolver.C:173
Foam::fvMatrix::A
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:788
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
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::dictionary::subOrEmptyDict
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Definition: dictionary.C:603
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:76
rAU
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
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:156
Foam::fvc::flux
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:34
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:58
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:927
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:243
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