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-2023 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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
46 {
47  maxCo =
48  runTime.controlDict().lookupOrDefault<scalar>("maxCo", 1.0);
49 
50  maxDeltaT_ =
51  runTime.controlDict().lookupOrDefault<scalar>("maxDeltaT", vGreat);
52 
54  (
55  "correctPhi",
56  mesh.dynamic()
57  );
58 
59  checkMeshCourantNo = pimple.dict().lookupOrDefault
60  (
61  "checkMeshCourantNo",
62  false
63  );
64 }
65 
66 
68 {
70  {
71  const scalarField sumPhi
72  (
73  fvc::surfaceSum(mag(mesh.phi()))().primitiveField()
74  );
75 
76  const scalar meshCoNum
77  (
78  0.5*gMax(sumPhi/mesh.V().field())*runTime.deltaTValue()
79  );
80 
81  const scalar meanMeshCoNum
82  (
83  0.5*(gSum(sumPhi)/gSum(mesh.V().field()))*runTime.deltaTValue()
84  );
85 
86  Info<< "Mesh Courant Number mean: " << meanMeshCoNum
87  << " max: " << meshCoNum << endl;
88  }
89 }
90 
91 
92 template<class RhoType>
93 void Foam::solvers::fluidSolver::correctCoNum
94 (
95  const RhoType& rho,
96  const surfaceScalarField& phi
97 )
98 {
99  const scalarField sumPhi
100  (
101  fvc::surfaceSum(mag(phi))().primitiveField()/rho.primitiveField()
102  );
103 
104  CoNum = 0.5*gMax(sumPhi/mesh.V().field())*runTime.deltaTValue();
105 
106  const scalar meanCoNum =
107  0.5*(gSum(sumPhi)/gSum(mesh.V().field()))*runTime.deltaTValue();
108 
109  Info<< "Courant Number mean: " << meanCoNum
110  << " max: " << CoNum << endl;
111 }
112 
113 
114 void Foam::solvers::fluidSolver::correctCoNum(const surfaceScalarField& phi)
115 {
116  correctCoNum(geometricOneField(), phi);
117 }
118 
119 
120 void Foam::solvers::fluidSolver::correctCoNum
121 (
122  const volScalarField& rho,
123  const surfaceScalarField& phi
124 )
125 {
126  correctCoNum<volScalarField>(rho, phi);
127 }
128 
129 
131 (
132  const surfaceScalarField& phi
133 )
134 {
135  const volScalarField contErr(fvc::div(phi));
136 
137  const scalar sumLocalContErr =
138  runTime.deltaTValue()
139  *mag(contErr)().weightedAverage(mesh.V()).value();
140 
141  const scalar globalContErr =
142  runTime.deltaTValue()
143  *contErr.weightedAverage(mesh.V()).value();
144 
145  Info<< "time step continuity errors : sum local = " << sumLocalContErr
146  << ", global = " << globalContErr;
147 
148  if (pimple.finalPisoIter() && pimple.finalIter())
149  {
151 
152  Info<< ", cumulative = " << cumulativeContErr;
153  }
154 
155  Info<< endl;
156 }
157 
158 
160 (
161  const volScalarField& rho,
162  const volScalarField& thermoRho,
163  const surfaceScalarField& phi
164 )
165 {
166  if (mesh.schemes().steady())
167  {
168  continuityErrors(phi);
169  }
170  else
171  {
172  const dimensionedScalar totalMass = fvc::domainIntegrate(rho);
173 
174  const scalar sumLocalContErr =
175  (fvc::domainIntegrate(mag(rho - thermoRho))/totalMass).value();
176 
177  const scalar globalContErr =
178  (fvc::domainIntegrate(rho - thermoRho)/totalMass).value();
179 
181 
182  Info<< "time step continuity errors : sum local = " << sumLocalContErr
183  << ", global = " << globalContErr
184  << ", cumulative = " << cumulativeContErr
185  << endl;
186  }
187 }
188 
189 
190 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
191 
193 :
194  solver(mesh),
196  CoNum(0)
197 {
198  // Read the controls
199  readControls();
200 }
201 
202 
203 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
204 
206 {}
207 
208 
209 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
210 
212 {
213  scalar deltaT = min(fvModels().maxDeltaT(), maxDeltaT_);
214 
215  if (CoNum > small)
216  {
217  deltaT = min(deltaT, maxCo/CoNum*runTime.deltaTValue());
218  }
219 
220  return deltaT;
221 }
222 
223 
224 // ************************************************************************* //
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &) const
Calculate and return weighted average.
Generic GeometricField class.
const dictionary & controlDict() const
Return the control dict.
Definition: Time.H:265
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
bool dynamic() const
Is this mesh dynamic?
Definition: fvMesh.C:641
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
virtual const dictionary & dict() const
Return the solution dictionary.
Abstract base class for run-time selectable region solvers.
Definition: solver.H:55
pimpleNoLoopControl pimple
PIMPLE inner-loop controls.
Definition: solver.H:100
const Time & runTime
Time.
Definition: solver.H:97
const fvMesh & mesh
Region mesh.
Definition: solver.H:94
Base solver module for fluid solvers.
Definition: fluidSolver.H:62
fluidSolver(fvMesh &mesh)
Construct from region mesh.
Definition: fluidSolver.C:192
bool correctPhi
Switch to correct the flux after mesh change.
Definition: fluidSolver.H:95
void continuityErrors(const surfaceScalarField &phi)
Calculate and print the continuity errors.
Definition: fluidSolver.C:131
virtual scalar maxDeltaT() const
Return the current maximum time-step for stable solution.
Definition: fluidSolver.C:211
void meshCourantNo() const
Check mesh Courant numbers for moving mesh cases.
Definition: fluidSolver.C:67
virtual ~fluidSolver()
Destructor.
Definition: fluidSolver.C:205
void readControls()
Read controls.
Definition: fluidSolver.C:45
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:251
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensioned< scalar > mag(const dimensioned< Type > &)
Type gMax(const FieldField< Field, Type > &f)
checkMeshCourantNo
Foam::surfaceFields.
scalar sumLocalContErr
Definition: volContinuity.H:9
scalar globalContErr
Definition: volContinuity.H:12