CorrectPhi.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "CorrectPhi.H"
27 #include "volFields.H"
28 #include "surfaceFields.H"
29 #include "fvScalarMatrix.H"
30 #include "fvmDdt.H"
31 #include "fvmLaplacian.H"
32 #include "fvcDiv.H"
35 #include "adjustPhi.H"
36 #include "fvcMeshPhi.H"
37 #include "pressureReference.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 template<class RAUfType, class DivUType>
44 (
45  surfaceScalarField& phi,
46  const volVectorField& U,
47  const volScalarField& p,
48  const RAUfType& rAUf,
49  const DivUType& divU,
50  const pressureReference& pressureReference,
51  nonOrthogonalSolutionControl& pcorrControl
52 )
53 {
54  const fvMesh& mesh = phi.mesh();
55  const Time& runTime = mesh.time();
56 
57  // Initialise BCs list for pcorr to zero-gradient
59  (
60  p.boundaryField().size(),
61  zeroGradientFvPatchScalarField::typeName
62  );
63 
64  // Set BCs of pcorr to fixed-value for patches at which p is fixed
66  {
67  if (p.boundaryField()[patchi].fixesValue())
68  {
69  pcorrTypes[patchi] = fixedValueFvPatchScalarField::typeName;
70  }
71  }
72 
73  volScalarField pcorr
74  (
75  IOobject
76  (
77  "pcorr",
78  runTime.timeName(),
79  mesh
80  ),
81  mesh,
84  );
85 
86  if (pcorr.needReference())
87  {
88  fvc::makeRelative(phi, U);
89  adjustPhi(phi, U, pcorr);
90  fvc::makeAbsolute(phi, U);
91  }
92 
93  mesh.schemes().setFluxRequired(pcorr.name());
94 
95  while (pcorrControl.correctNonOrthogonal())
96  {
97  // Solve for pcorr such that the divergence of the corrected flux
98  // matches the divU provided (from previous iteration, time-step...)
99  fvScalarMatrix pcorrEqn
100  (
101  fvm::laplacian(rAUf, pcorr) == fvc::div(phi) - divU
102  );
103 
104  pcorrEqn.setReference(pressureReference.refCell(), 0);
105 
106  pcorrEqn.solve();
107 
108  if (pcorrControl.finalNonOrthogonalIter())
109  {
110  phi -= pcorrEqn.flux();
111  }
112  }
113 }
114 
115 
116 template<class RAUfType, class DivRhoUType>
117 void Foam::CorrectPhi
118 (
119  surfaceScalarField& phi,
120  const volScalarField& p,
121  const volScalarField& rho,
122  const volScalarField& psi,
123  const RAUfType& rAUf,
124  const DivRhoUType& divRhoU,
125  nonOrthogonalSolutionControl& pcorrControl
126 )
127 {
128  const fvMesh& mesh = phi.mesh();
129  const Time& runTime = mesh.time();
130 
131  // Initialise BCs list for pcorr to zero-gradient
133  (
134  p.boundaryField().size(),
135  zeroGradientFvPatchScalarField::typeName
136  );
137 
138  // Set BCs of pcorr to fixed-value for patches at which p is fixed
140  {
141  if (p.boundaryField()[patchi].fixesValue())
142  {
143  pcorrTypes[patchi] = fixedValueFvPatchScalarField::typeName;
144  }
145  }
146 
147  volScalarField pcorr
148  (
149  IOobject
150  (
151  "pcorr",
152  runTime.timeName(),
153  mesh
154  ),
155  mesh,
157  pcorrTypes
158  );
159 
160  mesh.schemes().setFluxRequired(pcorr.name());
161 
162  while (pcorrControl.correctNonOrthogonal())
163  {
164  // Solve for pcorr such that the divergence of the corrected flux
165  // matches the divRhoU provided (from previous iteration, time-step...)
166  fvScalarMatrix pcorrEqn
167  (
168  fvm::ddt(psi, pcorr)
169  + fvc::div(phi)
170  - fvm::laplacian(rAUf, pcorr)
171  ==
172  divRhoU
173  );
174 
175  pcorrEqn.solve();
176 
177  if (pcorrControl.finalNonOrthogonalIter())
178  {
179  phi += pcorrEqn.flux();
180  }
181  }
182 }
183 
184 
185 // ************************************************************************* //
Foam::surfaceFields.
bool finalNonOrthogonalIter() const
Flag to indicate the last non-orthogonal iteration.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const word & name() const
Return name.
Definition: IOobject.H:315
bool correctNonOrthogonal()
Non-orthogonal corrector loop.
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1672
const Boundary & boundaryField() const
Return const-reference to the boundary field.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
void setReference(const label celli, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition: fvMatrix.C:563
Calculate the matrix for the laplacian of the field.
void setFluxRequired(const word &name) const
Definition: fvSchemes.C:487
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
fvMesh & mesh
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:35
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
bool needReference() const
Does the field need a reference level for solution.
const dimensionSet & dimensions() const
Return dimensions.
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
Calculate the matrix for the first temporal derivative.
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:72
Calculate the divergence of the given field.
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
const Mesh & mesh() const
Return mesh.
label refCell() const
Return the cell in which the reference pressure is set.
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
label patchi
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A scalar instance of fvMatrix.
wordList pcorrTypes(p.boundaryField().size(), zeroGradientFvPatchScalarField::typeName)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
void makeAbsolute(surfaceScalarField &phi, const volVectorField &U)
Make the given flux absolute.
Definition: fvcMeshPhi.C:128
For cases which do no have a pressure boundary adjust the balance of fluxes to obey continuity...
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:927
void CorrectPhi(surfaceScalarField &phi, const volVectorField &U, const volScalarField &p, const RAUfType &rAUf, const DivUType &divU, const pressureReference &pressureReference, nonOrthogonalSolutionControl &pcorrControl)
Definition: CorrectPhi.C:44
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:89
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
Non-orthogonal solution control class. Provides non-orthogonal-loop control methods.