simple.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 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 
30 #include "simple.H"
31 #include "findRefCell.H"
32 #include "constrainHbyA.H"
33 #include "constrainPressure.H"
34 #include "adjustPhi.H"
35 #include "Time.H"
36 #include "fvOptions.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
44  addToRunTimeSelectionTable(incompressiblePrimalSolver, simple, dictionary);
45 }
46 
47 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
48 
50 {
52  return getIncoVars();
53 }
54 
55 
57 {
58  if (incoVars_.useSolverNameForFields())
59  {
61  << "useSolverNameForFields is set to true for primalSolver "
62  << solverName() << nl << tab
63  << "Appending variable names with the solver name" << nl << tab
64  << "Please adjust the necessary entries in fvSchemes and fvSolution"
65  << nl << endl;
66  }
67 }
68 
69 
71 {
72  surfaceScalarField& phi = incoVars_.phiInst();
74 
75  scalar sumLocalContErr = mesh_.time().deltaTValue()*
76  mag(contErr)().weightedAverage(mesh_.V()).value();
77 
78  scalar globalContErr = mesh_.time().deltaTValue()*
79  contErr.weightedAverage(mesh_.V()).value();
80  cumulativeContErr_ += globalContErr;
81 
82  Info<< "time step continuity errors : sum local = " << sumLocalContErr
83  << ", global = " << globalContErr
84  << ", cumulative = " << cumulativeContErr_
85  << endl;
86 }
87 
88 
89 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
90 
91 Foam::simple::simple
92 (
93  fvMesh& mesh,
94  const word& managerType,
95  const dictionary& dict
96 )
97 :
98  incompressiblePrimalSolver(mesh, managerType, dict),
99  solverControl_(SIMPLEControl::New(mesh, managerType, *this)),
100  incoVars_(allocateVars()),
101  MRF_(mesh, word(useSolverNameForFields() ? solverName() : word::null)),
102  cumulativeContErr_(Zero),
103  objectives_(0)
104 {
105  addExtraSchemes();
106  setRefCell
107  (
108  incoVars_.pInst(),
109  solverControl_().dict(),
110  solverControl_().pRefCell(),
111  solverControl_().pRefValue()
112  );
113 }
114 
115 
116 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117 
119 {
121  {
122  return true;
123  }
124 
125  return false;
126 }
127 
129 {
130  preIter();
131  mainIter();
132  postIter();
133 }
134 
135 
137 {
138  Info<< "Time = " << mesh_.time().timeName() << "\n" << endl;
139 }
140 
141 
143 {
144  // Grab references
145  volScalarField& p = incoVars_.pInst();
146  volVectorField& U = incoVars_.UInst();
147  surfaceScalarField& phi = incoVars_.phiInst();
149  incoVars_.turbulence();
150  label& pRefCell = solverControl_().pRefCell();
151  scalar& pRefValue = solverControl_().pRefValue();
152  fv::options& fvOptions(fv::options::New(this->mesh_));
153 
154  // Momentum predictor
155  //~~~~~~~~~~~~~~~~~~~
156 
157  MRF_.correctBoundaryVelocity(U);
158 
160  (
161  fvm::div(phi, U)
162  + MRF_.DDt(U)
163  + turbulence->divDevReff(U)
164  ==
165  fvOptions(U)
166  );
167  fvVectorMatrix& UEqn = tUEqn.ref();
168 
169  UEqn.relax();
170 
172 
173  if (solverControl_().momentumPredictor())
174  {
175  Foam::solve(UEqn == -fvc::grad(p));
176 
178  }
179 
180  // Pressure Eq
181  //~~~~~~~~~~~~
182  {
183  volScalarField rAU(1.0/UEqn.A());
186  MRF_.makeRelative(phiHbyA);
187  adjustPhi(phiHbyA, U, p);
188 
190 
191  if (solverControl_().consistent())
192  {
193  rAtU = 1.0/(1.0/rAU - UEqn.H1());
194  phiHbyA +=
195  fvc::interpolate(rAtU() - rAU)*fvc::snGrad(p)*mesh_.magSf();
196  HbyA -= (rAU - rAtU())*fvc::grad(p);
197  }
198 
199  tUEqn.clear();
200 
201  // Update the pressure BCs to ensure flux consistency
202  constrainPressure(p, U, phiHbyA, rAtU(), MRF_);
203 
204  // Non-orthogonal pressure corrector loop
205  while (solverControl_().correctNonOrthogonal())
206  {
207  fvScalarMatrix pEqn
208  (
210  );
211 
213 
214  pEqn.solve();
215 
216  if (solverControl_().finalNonOrthogonalIter())
217  {
218  phi = phiHbyA - pEqn.flux();
219  }
220  }
221 
222  continuityErrors();
223 
224  // Explicitly relax pressure for momentum corrector
225  p.relax();
226 
227  // Momentum corrector
228  U = HbyA - rAtU()*fvc::grad(p);
229  U.correctBoundaryConditions();
231  }
232 
233  incoVars_.laminarTransport().correct();
234  turbulence->correct();
235 }
236 
237 
239 {
240  solverControl_().write();
241 
242  // Print objective values to screen and compute mean value
243  Info<< endl;
244  for (objective* obj : objectives_)
245  {
246  Info<< obj->objectiveName() << " : " << obj->J() << endl;
247  obj->accumulateJMean(solverControl_());
248  obj->writeInstantaneousValue();
249  }
250 
251  // Average fields if necessary
252  incoVars_.computeMeanFields();
253 
254  // Print execution time
255  mesh_.time().printExecutionTime(Info);
256 }
257 
258 
260 {
261  // Iterate
262  if (active_)
263  {
264  preLoop();
265  while (solverControl_().loop())
266  {
267  solveIter();
268  }
269  postLoop();
270  }
271 }
272 
273 
275 {
276  return solverControl_().loop();
277 }
278 
279 
281 {
282  incoVars_.restoreInitValues();
283 }
284 
285 
287 {
288  // Get the objectives for this solver
289  if (objectives_.empty())
290  {
291  objectives_ = getObjectiveFunctions();
292  }
293 
294  // Reset initial and mean fields before solving
295  restoreInitValues();
296  incoVars_.resetMeanFields();
297 
298  // Validate turbulence model fields
299  incoVars_.turbulence()->validate();
300 }
301 
302 
304 {
305  for (objective* obj : objectives_)
306  {
307  obj->writeInstantaneousSeparator();
308  }
309 
310  // Safety
311  objectives_.clear();
312 }
313 
314 
316 {
317  os.writeEntry("averageIter", solverControl_().averageIter());
318 
319  return true;
320 }
321 
322 
323 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::constrainPressure
void constrainPressure(volScalarField &p, const RhoType &rho, const volVectorField &U, const surfaceScalarField &phiHbyA, const RAUType &rhorAU, const MRFType &MRF)
Definition: constrainPressure.C:39
Foam::fvMatrix::H
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition: fvMatrix.C:1423
Foam::constrainHbyA
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:35
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:47
sumLocalContErr
scalar sumLocalContErr
Definition: compressibleContinuityErrs.H:34
rAtU
volScalarField rAtU(1.0/(1.0/rAU - UEqn.H1()))
contErr
volScalarField::Internal contErr((fvc::ddt(rho)+fvc::div(rhoPhi) -(fvOptions(alpha1, mixture.thermo1().rho())&rho1) -(fvOptions(alpha2, mixture.thermo2().rho())&rho2))())
simple.H
Foam::adjustPhi
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:37
simple
const dictionary & simple
Definition: readFluidMultiRegionSIMPLEControls.H:1
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::fvc::flux
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
momentumPredictor
const bool momentumPredictor
Definition: readFluidMultiRegionSIMPLEControls.H:6
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
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
Foam::simple::solverControl_
autoPtr< SIMPLEControl > solverControl_
Solver control.
Definition: simple.H:72
Foam::fv::optionList::correct
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
Definition: fvOptionListTemplates.C:355
Foam::simple::mainIter
virtual void mainIter()
The main SIMPLE iter.
Definition: simple.C:142
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::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
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
findRefCell.H
Find the reference cell nearest (in index) to the given cell but which is not on a cyclic,...
Foam::simple::postIter
virtual void postIter()
Steps to be executed before each main SIMPLE iteration.
Definition: simple.C:238
Foam::solver::mesh_
fvMesh & mesh_
Reference to the mesh database.
Definition: solver.H:60
constrainPressure.H
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
HbyA
HbyA
Definition: pcEqn.H:74
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::fvMatrix::setReference
void setReference(const label celli, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition: fvMatrix.C:1092
Foam::simple::addExtraSchemes
void addExtraSchemes()
Definition: simple.C:56
adjustPhi.H
For cases which do no have a pressure boundary adjust the balance of fluxes to obey continuity....
Foam::simple::continuityErrors
void continuityErrors()
Compute continuity errors.
Definition: simple.C:70
Foam::simple::solveIter
virtual void solveIter()
Execute one iteration of the solution algorithm.
Definition: simple.C:128
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
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::simple::writeData
virtual bool writeData(Ostream &os) const
Write average iteration.
Definition: simple.C:315
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::H1
tmp< volScalarField > H1() const
Return H(1)
Definition: fvMatrix.C:1485
Foam::fvMatrix::relax
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1183
Foam::simple::preLoop
virtual void preLoop()
Functions to be called before loop.
Definition: simple.C:286
Foam::incompressiblePrimalSolver::readDict
virtual bool readDict(const dictionary &dict)
Read dict if updated.
Definition: incompressiblePrimalSolver.C:109
Foam::setRefCell
bool setRefCell(const volScalarField &field, const volScalarField &fieldRef, const dictionary &dict, label &refCelli, scalar &refValue, const bool forceReference=false)
If the field fieldRef needs referencing find the reference cell nearest.
Definition: findRefCell.C:34
Foam::solve
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::simple::restoreInitValues
virtual void restoreInitValues()
Restore initial field values if necessary.
Definition: simple.C:280
Foam::fv::options
Finite-volume options.
Definition: fvOptions.H:55
Foam::incompressiblePrimalSolver
Base class for primal incompressible solvers.
Definition: incompressiblePrimalSolver.H:53
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
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMatrix::A
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:1394
Foam::simple::loop
virtual bool loop()
Looper (advances iters, time step)
Definition: simple.C:274
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::fvc::interpolate
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Foam::SIMPLEControl::New
static autoPtr< SIMPLEControl > New(fvMesh &mesh, const word &managerType, const solver &solver)
Return a reference to the selected turbulence model.
Definition: SIMPLEControl.C:64
globalContErr
scalar globalContErr
Definition: compressibleContinuityErrs.H:37
Foam::simple::preIter
virtual void preIter()
Steps to be executed before each main SIMPLE iteration.
Definition: simple.C:136
Time.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
U
U
Definition: pEqn.H:72
Foam::tab
constexpr char tab
Definition: Ostream.H:403
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::solver::vars_
autoPtr< variablesSet > vars_
Base variableSet pointer.
Definition: solver.H:82
Foam::DimensionedField::weightedAverage
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &weightField) const
Calculate and return weighted average.
Definition: DimensionedField.C:461
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::simple::allocateVars
incompressibleVars & allocateVars()
Protected Member Functions.
Definition: simple.C:49
Foam::word::null
static const word null
An empty word.
Definition: word.H:80
Foam::GeometricField::relax
void relax(const scalar alpha)
Relax field (for steady-state solution).
Definition: GeometricField.C:972
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:319
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
Foam::simple::readDict
virtual bool readDict(const dictionary &dict)
Read dict if updated.
Definition: simple.C:118
UEqn
fvVectorMatrix & UEqn
Definition: UEqn.H:13
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::incompressiblePrimalSolver::getIncoVars
const incompressibleVars & getIncoVars() const
Access to the incompressible variables set.
Definition: incompressiblePrimalSolver.C:154
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >
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)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
constrainHbyA.H
Foam::fvMatrix::flux
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:1534
Foam::incompressibleVars
Base class for solution control classes.
Definition: incompressibleVars.H:54
pRefCell
const label pRefCell
Definition: setRegionFluidFields.H:36
Foam::simple::solve
virtual void solve()
Main control loop.
Definition: simple.C:259
Foam::simple::postLoop
virtual void postLoop()
Functions to be called after loop.
Definition: simple.C:303
pRefValue
const scalar pRefValue
Definition: setRegionFluidFields.H:37