solid.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 "solid.H"
27 #include "fvcSurfaceIntegrate.H"
28 #include "fvMeshMover.H"
29 #include "localEulerDdtScheme.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace solvers
37 {
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 void Foam::solvers::solid::correctDiNum()
47 {
48  const volScalarField kappa
49  (
51  ? thermo.kappa()
52  : mag(thermo.Kappa())()
53  );
54 
55  const volScalarField::Internal DiNumvf
56  (
58  (
59  mesh.magSf()
61  *mesh.surfaceInterpolation::deltaCoeffs()
62  )()()
63  /(mesh.V()*thermo.rho()()*thermo.Cp()())
64  *runTime.deltaT()
65  );
66 
67  const scalar meanDiNum = gAverage(DiNumvf);
68  const scalar maxDiNum = gMax(DiNumvf);
69 
70  Info<< "Diffusion Number mean: " << meanDiNum
71  << " max: " << maxDiNum << endl;
72 
73  DiNum = maxDiNum;
74 }
75 
76 
77 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
78 
80 {
81  return runTime.controlDict().modified();
82 }
83 
84 
86 {
87  solver::read();
88 
89  maxDi =
90  runTime.controlDict().lookupOrDefault<scalar>("maxDi", 1.0);
91 
92  maxDeltaT_ =
93  runTime.controlDict().found("maxDeltaT")
94  ? runTime.controlDict().lookup<scalar>("maxDeltaT", runTime.userUnits())
95  : vGreat;
96 
97  return true;
98 }
99 
100 
101 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
102 
104 (
105  fvMesh& mesh,
106  autoPtr<solidThermo> thermoPtr
107 )
108 :
109  solver(mesh),
110 
111  thermoPtr_(thermoPtr),
112  thermo_(thermoPtr_()),
113 
114  T_(thermo_.T()),
115 
117 
118  DiNum(0),
119 
120  thermo(thermo_),
121  T(T_)
122 {
123  thermo.validate("solid", "h", "e");
124 
125  if (transient())
126  {
127  correctDiNum();
128  }
129  else if (LTS)
130  {
131  FatalError
132  << type()
133  << " solver does not support LTS, use 'steadyState' ddtScheme"
134  << exit(FatalError);
135  }
136 }
137 
138 
139 
141 :
142  solid(mesh, solidThermo::New(mesh))
143 {
144  // Read the controls
145  read();
146 }
147 
148 
149 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
150 
152 {}
153 
154 
155 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
156 
158 {
159  scalar deltaT = min(fvModels().maxDeltaT(), maxDeltaT_);
160 
161  if (DiNum > small)
162  {
163  deltaT = min(deltaT, maxDi/DiNum*runTime.deltaTValue());
164  }
165 
166  return deltaT;
167 }
168 
169 
171 {
173 
174  if (transient())
175  {
176  correctDiNum();
177  }
178 
179  // Update the mesh for topology change, mesh to mesh mapping
180  mesh_.update();
181 }
182 
183 
185 {
186  if (pimple.firstIter() || pimple.moveMeshOuterCorrectors())
187  {
188  if (!mesh_.mover().solidBody())
189  {
191  << "Region " << name() << " of type " << type()
192  << " does not support non-solid body mesh motion"
193  << exit(FatalError);
194  }
195 
196  mesh_.move();
197  }
198 }
199 
200 
202 {}
203 
204 
206 {
207  if (pimple.predictTransport())
208  {
209  thermophysicalTransport->predict();
210  }
211 }
212 
213 
215 {}
216 
217 
219 {
220  volScalarField& e = thermo_.he();
221  const volScalarField& rho = thermo_.rho();
222 
223  while (pimple.correctNonOrthogonal())
224  {
225  fvScalarMatrix eEqn
226  (
227  fvm::ddt(rho, e)
228  + thermophysicalTransport->divq(e)
229  ==
230  fvModels().source(rho, e)
231  );
232 
233  eEqn.relax();
234 
235  fvConstraints().constrain(eEqn);
236 
237  eEqn.solve();
238 
240 
241  thermo_.correct();
242  }
243 }
244 
245 
247 {}
248 
249 
251 {
252  if (pimple.correctTransport())
253  {
254  thermophysicalTransport->correct();
255  }
256 }
257 
258 
260 {}
261 
262 
263 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
DimensionedField< Type, GeoMesh > Internal
Type of the internal field from which this GeometricField is derived.
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:46
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
virtual const volScalarField & Cp() const =0
Heat capacity at constant pressure [J/kg/K].
void validate(const string &app, const word &) const
Check that the thermodynamics package is consistent.
Definition: basicThermo.C:308
virtual const volScalarField & kappa() const =0
Thermal conductivity of mixture [W/m/K].
virtual tmp< volScalarField > rho() const =0
Density [kg/m^3].
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:603
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Finite volume models.
Definition: fvModels.H:65
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: fvModels.C:260
Base-class for solid thermodynamic properties.
Definition: solidThermo.H:59
virtual const volVectorField & Kappa() const =0
Anisotropic thermal conductivity [W/m/K].
virtual bool isotropic() const =0
Return true if thermal conductivity is isotropic.
Abstract base class for solid thermophysical transport models.
Abstract base class for run-time selectable region solvers.
Definition: solver.H:56
bool LTS
Switch for local time step transient operation.
Definition: solver.H:70
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
Solver module for thermal transport in solid domains and regions for conjugate heat transfer,...
Definition: solid.H:56
virtual void thermophysicalPredictor()
Construct and solve the energy equation,.
Definition: solid.C:218
virtual void prePredictor()
Called at the beginning of the PIMPLE loop.
Definition: solid.C:205
virtual void postSolve()
Called after the PIMPLE loop at the end of the time-step.
Definition: solid.C:259
virtual bool dependenciesModified() const
Return true if the solver's dependencies have been modified.
Definition: solid.C:79
virtual void moveMesh()
Called at the start of the PIMPLE loop to move the mesh.
Definition: solid.C:184
virtual scalar maxDeltaT() const
Return the current maximum time-step for stable solution.
Definition: solid.C:157
virtual void motionCorrector()
Corrections that follow mesh motion.
Definition: solid.C:201
const solidThermo & thermo
Reference to the solid thermophysical properties.
Definition: solid.H:109
virtual void pressureCorrector()
Construct and solve the pressure equation in the PISO loop.
Definition: solid.C:246
virtual void postCorrector()
Correct the thermophysical transport modelling.
Definition: solid.C:250
virtual void momentumPredictor()
Construct and optionally solve the momentum equation.
Definition: solid.C:214
solid(fvMesh &mesh, autoPtr< solidThermo >)
Construct from region mesh and thermophysical properties.
Definition: solid.C:104
virtual void preSolve()
Called at the start of the time-step, before the PIMPLE loop.
Definition: solid.C:170
virtual bool read()
Read controls.
Definition: solid.C:85
virtual ~solid()
Destructor.
Definition: solid.C:151
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
pimpleControl pimple(mesh)
scalar maxDeltaT
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
Info<< "Creating thermophysical transport model\n"<< endl;turbulenceThermophysicalTransportModels::unityLewisEddyDiffusivity< RASThermophysicalTransportModel< ThermophysicalTransportModel< compressibleMomentumTransportModel, fluidThermo > >> thermophysicalTransport(turbulence(), thermo, true)
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
tmp< VolField< Type > > surfaceSum(const SurfaceField< Type > &ssf)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
addToRunTimeSelectionTable(solver, compressibleMultiphaseVoF, fvMesh)
defineTypeNameAndDebug(compressibleMultiphaseVoF, 0)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const doubleScalar e
Definition: doubleScalar.H:106
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Type gAverage(const FieldField< Field, Type > &f)
error FatalError
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Type gMax(const FieldField< Field, Type > &f)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
fluidMulticomponentThermo & thermo
Definition: createFields.H:31