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 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"
34 #include "adjustPhi.H"
35 #include "fvcMeshPhi.H"
36 #include "pimpleControl.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 template<class RAUfType, class DivUType>
42 (
43  volVectorField& U,
44  surfaceScalarField& phi,
45  const volScalarField& p,
46  const RAUfType& rAUf,
47  const DivUType& divU,
48  pimpleControl& pimple
49 )
50 {
51  const fvMesh& mesh = U.mesh();
52  const Time& runTime = mesh.time();
53 
54  correctUphiBCs(U, phi);
55 
56  // Initialize BCs list for pcorr to zero-gradient
58  (
59  p.boundaryField().size(),
60  zeroGradientFvPatchScalarField::typeName
61  );
62 
63  // Set BCs of pcorr to fixed-value for patches at which p is fixed
65  {
66  if (p.boundaryField()[patchi].fixesValue())
67  {
68  pcorrTypes[patchi] = fixedValueFvPatchScalarField::typeName;
69  }
70  }
71 
72  volScalarField pcorr
73  (
74  IOobject
75  (
76  "pcorr",
77  runTime.timeName(),
78  mesh
79  ),
80  mesh,
81  dimensionedScalar("pcorr", p.dimensions(), 0.0),
83  );
84 
85  if (pcorr.needReference())
86  {
87  fvc::makeRelative(phi, U);
88  adjustPhi(phi, U, pcorr);
89  fvc::makeAbsolute(phi, U);
90  }
91 
92  mesh.setFluxRequired(pcorr.name());
93 
94  while (pimple.correctNonOrthogonal())
95  {
96  // Solve for pcorr such that the divergence of the corrected flux
97  // matches the divU provided (from previous iteration, time-step...)
98  fvScalarMatrix pcorrEqn
99  (
100  fvm::laplacian(rAUf, pcorr) == fvc::div(phi) - divU
101  );
102 
103  pcorrEqn.setReference(0, 0);
104  pcorrEqn.solve();
105 
106  if (pimple.finalNonOrthogonalIter())
107  {
108  phi -= pcorrEqn.flux();
109  }
110  }
111 }
112 
113 
114 template<class RAUfType, class DivRhoUType>
115 void Foam::CorrectPhi
116 (
117  volVectorField& U,
118  surfaceScalarField& phi,
119  const volScalarField& p,
120  const volScalarField& rho,
121  const volScalarField& psi,
122  const RAUfType& rAUf,
123  const DivRhoUType& divRhoU,
124  pimpleControl& pimple
125 )
126 {
127  const fvMesh& mesh = U.mesh();
128  const Time& runTime = mesh.time();
129 
130  correctUphiBCs(rho, U, phi);
131 
132  // Initialize BCs list for pcorr to zero-gradient
134  (
135  p.boundaryField().size(),
136  zeroGradientFvPatchScalarField::typeName
137  );
138 
139  // Set BCs of pcorr to fixed-value for patches at which p is fixed
141  {
142  if (p.boundaryField()[patchi].fixesValue())
143  {
144  pcorrTypes[patchi] = fixedValueFvPatchScalarField::typeName;
145  }
146  }
147 
148  volScalarField pcorr
149  (
150  IOobject
151  (
152  "pcorr",
153  runTime.timeName(),
154  mesh
155  ),
156  mesh,
157  dimensionedScalar("pcorr", p.dimensions(), 0.0),
158  pcorrTypes
159  );
160 
161  mesh.setFluxRequired(pcorr.name());
162 
163  while (pimple.correctNonOrthogonal())
164  {
165  // Solve for pcorr such that the divergence of the corrected flux
166  // matches the divRhoU provided (from previous iteration, time-step...)
167  fvScalarMatrix pcorrEqn
168  (
169  fvm::ddt(psi, pcorr)
170  + fvc::div(phi)
171  - fvm::laplacian(rAUf, pcorr)
172  ==
173  divRhoU
174  );
175 
176  pcorrEqn.solve();
177 
178  if (pimple.finalNonOrthogonalIter())
179  {
180  phi += pcorrEqn.flux();
181  }
182  }
183 }
184 
185 
186 // ************************************************************************* //
bool correctNonOrthogonal()
Non-orthogonal corrector loop.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::surfaceFields.
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
void makeAbsolute(surfaceScalarField &phi, const volVectorField &U)
Make the given flux absolute.
Definition: fvcMeshPhi.C:113
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:883
void setFluxRequired(const word &name) const
Definition: fvSchemes.C:497
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:34
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
Calulate the matrix for the first temporal derivative.
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
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:68
const Mesh & mesh() const
Return mesh.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
dynamicFvMesh & mesh
correctUphiBCs(U, phi)
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
PIMPLE control class to supply convergence information/checks for the PIMPLE loop.
Definition: pimpleControl.H:49
bool finalNonOrthogonalIter() const
Helper function to identify final non-orthogonal iteration.
#define forAll(list, i)
Definition: UList.H:421
label patchi
const dimensionSet & dimensions() const
Return dimensions.
solverPerformance solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:57
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:74
wordList pcorrTypes(p.boundaryField().size(), zeroGradientFvPatchScalarField::typeName)
Calculate the matrix for the laplacian of the field.
Calculate the divergence of the given field.
For cases which do no have a pressure boundary adjust the balance of fluxes to obey continuity...
void setReference(const label celli, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition: fvMatrix.C:510
A scalar instance of fvMatrix.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void CorrectPhi(volVectorField &U, surfaceScalarField &phi, const volScalarField &p, const RAUfType &rAUf, const DivUType &divU, pimpleControl &pimple)
Definition: CorrectPhi.C:42