fluidSolver.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) 2022-2024 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 "fluidSolver.H"
27 #include "surfaceFields.H"
28 #include "fvcDiv.H"
29 #include "fvcSurfaceIntegrate.H"
30 #include "fvcVolumeIntegrate.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace solvers
37 {
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
44 
46 {
48 }
49 
50 
52 {
53  solver::read();
54 
55  maxCo =
56  runTime.controlDict().lookupOrDefault<scalar>("maxCo", vGreat);
57 
58  maxDeltaT_ =
59  runTime.controlDict().found("maxDeltaT")
60  ? runTime.controlDict().lookup<scalar>("maxDeltaT", runTime.userUnits())
61  : vGreat;
62 
63  correctPhi = pimple.dict().lookupOrDefault
64  (
65  "correctPhi",
66  mesh.dynamic()
67  );
68 
69  checkMeshCourantNo = pimple.dict().lookupOrDefault
70  (
71  "checkMeshCourantNo",
72  false
73  );
74 
75  return true;
76 }
77 
79 {
81  {
82  const scalarField sumPhi
83  (
84  fvc::surfaceSum(mag(mesh.phi()))().primitiveField()
85  );
86 
87  const scalar meshCoNum
88  (
89  0.5*gMax(sumPhi/mesh.V().primitiveField())*runTime.deltaTValue()
90  );
91 
92  const scalar meanMeshCoNum
93  (
94  0.5
95  *(gSum(sumPhi)/gSum(mesh.V().primitiveField()))
96  *runTime.deltaTValue()
97  );
98 
99  Info<< "Mesh Courant Number mean: " << meanMeshCoNum
100  << " max: " << meshCoNum << endl;
101  }
102 }
103 
104 
105 template<class RhoType>
106 void Foam::solvers::fluidSolver::correctCoNum
107 (
108  const RhoType& rho,
109  const surfaceScalarField& phi
110 )
111 {
112  const scalarField sumPhi
113  (
114  fvc::surfaceSum(mag(phi))().primitiveField()/rho.primitiveField()
115  );
116 
117  CoNum_ = 0.5*gMax(sumPhi/mesh.V().primitiveField())*runTime.deltaTValue();
118 
119  const scalar meanCoNum =
120  0.5
121  *(gSum(sumPhi)/gSum(mesh.V().primitiveField()))
122  *runTime.deltaTValue();
123 
124  Info<< "Courant Number mean: " << meanCoNum
125  << " max: " << CoNum << endl;
126 }
127 
128 
129 void Foam::solvers::fluidSolver::correctCoNum(const surfaceScalarField& phi)
130 {
131  correctCoNum(geometricOneField(), phi);
132 }
133 
134 
135 void Foam::solvers::fluidSolver::correctCoNum
136 (
137  const volScalarField& rho,
138  const surfaceScalarField& phi
139 )
140 {
141  correctCoNum<volScalarField>(rho, phi);
142 }
143 
144 
146 (
147  const surfaceScalarField& phi
148 )
149 {
150  const volScalarField contErr(fvc::div(phi));
151 
152  const scalar sumLocalContErr =
153  runTime.deltaTValue()
154  *mag(contErr)().weightedAverage(mesh.V()).value();
155 
156  const scalar globalContErr =
157  runTime.deltaTValue()
158  *contErr.weightedAverage(mesh.V()).value();
159 
160  Info<< "time step continuity errors : sum local = " << sumLocalContErr
161  << ", global = " << globalContErr;
162 
163  if (pimple.finalPisoIter() && pimple.finalIter())
164  {
166 
167  Info<< ", cumulative = " << cumulativeContErr;
168  }
169 
170  Info<< endl;
171 }
172 
173 
175 (
176  const volScalarField& rho,
177  const volScalarField& thermoRho,
178  const surfaceScalarField& phi
179 )
180 {
181  if (mesh.schemes().steady())
182  {
183  continuityErrors(phi);
184  }
185  else
186  {
187  const dimensionedScalar totalMass = fvc::domainIntegrate(rho);
188 
189  const scalar sumLocalContErr =
190  (fvc::domainIntegrate(mag(rho - thermoRho))/totalMass).value();
191 
192  const scalar globalContErr =
193  (fvc::domainIntegrate(rho - thermoRho)/totalMass).value();
194 
196 
197  Info<< "time step continuity errors : sum local = " << sumLocalContErr
198  << ", global = " << globalContErr
199  << ", cumulative = " << cumulativeContErr
200  << endl;
201  }
202 }
203 
204 
205 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
206 
208 :
209  solver(mesh),
210  maxCo(0),
211  maxDeltaT_(0),
213  CoNum_(0),
214  CoNum(CoNum_)
215 {
216  // Read the controls
217  read();
218 }
219 
220 
221 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
222 
224 {}
225 
226 
227 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
228 
230 {
231  scalar deltaT = min(fvModels().maxDeltaT(), maxDeltaT_);
232 
233  if (maxCo < vGreat && CoNum > small)
234  {
235  deltaT = min(deltaT, maxCo/CoNum*runTime.deltaTValue());
236  }
237 
238  return deltaT;
239 }
240 
241 
242 // ************************************************************************* //
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &) const
Calculate and return weighted average.
Generic GeometricField class.
const IOdictionary & controlDict() const
Return the control dict.
Definition: Time.H:269
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const fvSolution & solution() const
Return the fvSolution.
Definition: fvMesh.C:1762
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
virtual bool modified() const
Return true if the object's file (or files for objectRegistry)
Abstract base class for run-time selectable region solvers.
Definition: solver.H:56
const Time & runTime
Time.
Definition: solver.H:104
const fvMesh & mesh
Region mesh.
Definition: solver.H:101
virtual bool read()
Read controls.
Definition: solver.C:52
Base solver module for fluid solvers.
Definition: fluidSolver.H:62
fluidSolver(fvMesh &mesh)
Construct from region mesh.
Definition: fluidSolver.C:207
void continuityErrors(const surfaceScalarField &phi)
Calculate and print the continuity errors.
Definition: fluidSolver.C:146
virtual bool dependenciesModified() const
Return true if the solver's dependencies have been modified.
Definition: fluidSolver.C:45
virtual scalar maxDeltaT() const
Return the current maximum time-step for stable solution.
Definition: fluidSolver.C:229
void meshCourantNo() const
Check mesh Courant numbers for moving mesh cases.
Definition: fluidSolver.C:78
virtual ~fluidSolver()
Destructor.
Definition: fluidSolver.C:223
virtual bool read()
Read controls.
Definition: fluidSolver.C:51
scalar CoNum
scalar meanCoNum
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
pimpleControl pimple(mesh)
scalar maxDeltaT
scalar maxCo
Calculate the divergence of the given field.
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
Volume integrate volField creating a volField.
scalar cumulativeContErr
tmp< VolField< Type > > surfaceSum(const SurfaceField< Type > &ssf)
dimensioned< Type > domainIntegrate(const VolField< Type > &vf)
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
defineTypeNameAndDebug(compressibleMultiphaseVoF, 0)
Namespace for OpenFOAM.
Type gSum(const FieldField< Field, Type > &f)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensioned< scalar > mag(const dimensioned< Type > &)
Type gMax(const FieldField< Field, Type > &f)
checkMeshCourantNo
correctPhi
Foam::surfaceFields.
scalar sumLocalContErr
Definition: volContinuity.H:9
scalar globalContErr
Definition: volContinuity.H:12