CorrectPhi.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 | Copyright (C) 2015-2016 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 "pimpleControl.H"
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 template<class RAUfType, class DivUType>
43 (
44  volVectorField& U,
45  surfaceScalarField& phi,
46  const volScalarField& p,
47  const RAUfType& rAUf,
48  const DivUType& divU,
49  pimpleControl& pimple
50 )
51 {
52  const fvMesh& mesh = U.mesh();
53  const Time& runTime = mesh.time();
54 
55  correctUphiBCs(U, phi);
56 
57  // Initialize 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,
82  dimensionedScalar("pcorr", p.dimensions(), 0.0),
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.setFluxRequired(pcorr.name());
94 
95  while (pimple.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(0, 0);
105  pcorrEqn.solve();
106 
107  if (pimple.finalNonOrthogonalIter())
108  {
109  phi -= pcorrEqn.flux();
110  }
111  }
112 }
113 
114 
115 template<class RAUfType, class DivRhoUType>
116 void Foam::CorrectPhi
117 (
118  volVectorField& U,
119  surfaceScalarField& phi,
120  const volScalarField& p,
121  const volScalarField& rho,
122  const volScalarField& psi,
123  const RAUfType& rAUf,
124  const DivRhoUType& divRhoU,
125  pimpleControl& pimple
126 )
127 {
128  const fvMesh& mesh = U.mesh();
129  const Time& runTime = mesh.time();
130 
131  correctUphiBCs(rho, U, phi);
132 
133  // Initialize BCs list for pcorr to zero-gradient
135  (
136  p.boundaryField().size(),
137  zeroGradientFvPatchScalarField::typeName
138  );
139 
140  // Set BCs of pcorr to fixed-value for patches at which p is fixed
142  {
143  if (p.boundaryField()[patchi].fixesValue())
144  {
145  pcorrTypes[patchi] = fixedValueFvPatchScalarField::typeName;
146  }
147  }
148 
149  volScalarField pcorr
150  (
151  IOobject
152  (
153  "pcorr",
154  runTime.timeName(),
155  mesh
156  ),
157  mesh,
158  dimensionedScalar("pcorr", p.dimensions(), 0.0),
159  pcorrTypes
160  );
161 
162  mesh.setFluxRequired(pcorr.name());
163 
164  while (pimple.correctNonOrthogonal())
165  {
166  // Solve for pcorr such that the divergence of the corrected flux
167  // matches the divRhoU provided (from previous iteration, time-step...)
168  fvScalarMatrix pcorrEqn
169  (
170  fvm::ddt(psi, pcorr)
171  + fvc::div(phi)
172  - fvm::laplacian(rAUf, pcorr)
173  ==
174  divRhoU
175  );
176 
177  pcorrEqn.solve();
178 
179  if (pimple.finalNonOrthogonalIter())
180  {
181  phi += pcorrEqn.flux();
182  }
183  }
184 }
185 
186 
187 // ************************************************************************* //
void CorrectPhi(volVectorField &U, surfaceScalarField &phi, const volScalarField &p, const RAUfType &rAUf, const DivUType &divU, pimpleControl &pimple)
Definition: CorrectPhi.C:43
Foam::surfaceFields.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
void setFluxRequired(const word &name) const
Definition: fvSchemes.C:497
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:496
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:864
Calculate the matrix for the laplacian of the field.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
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:68
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:34
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
dynamicFvMesh & mesh
const Boundary & boundaryField() const
Return const-reference to the boundary field.
bool correctNonOrthogonal()
Non-orthogonal corrector loop.
Calulate the matrix for the first temporal derivative.
bool finalNonOrthogonalIter() const
Helper function to identify final non-orthogonal iteration.
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:71
const dimensionSet & dimensions() const
Return dimensions.
Calculate the divergence of the given field.
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:57
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
correctUphiBCs(U, phi)
label patchi
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Mesh & mesh() const
Return mesh.
A scalar instance of fvMatrix.
PIMPLE control class to supply convergence information/checks for the PIMPLE loop.
Definition: pimpleControl.H:52
wordList pcorrTypes(p.boundaryField().size(), zeroGradientFvPatchScalarField::typeName)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void makeAbsolute(surfaceScalarField &phi, const volVectorField &U)
Make the given flux absolute.
Definition: fvcMeshPhi.C:114
For cases which do no have a pressure boundary adjust the balance of fluxes to obey continuity...
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:75
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243