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-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 "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 
47 {
48  maxDi =
49  runTime.controlDict().lookupOrDefault<scalar>("maxDi", 1.0);
50 
51  maxDeltaT_ =
52  runTime.controlDict().lookupOrDefault<scalar>("maxDeltaT", vGreat);
53 }
54 
55 
56 void Foam::solvers::solid::correctDiNum()
57 {
58  const volScalarField kappa
59  (
60  thermo.isotropic()
61  ? thermo.kappa()
62  : mag(thermo.Kappa())()
63  );
64 
65  const volScalarField::Internal DiNumvf
66  (
68  (
69  mesh.magSf()
71  *mesh.surfaceInterpolation::deltaCoeffs()
72  )()()
73  /(mesh.V()*thermo.rho()()*thermo.Cp()())
74  *runTime.deltaT()
75  );
76 
77  const scalar meanDiNum = gAverage(DiNumvf);
78  const scalar maxDiNum = gMax(DiNumvf);
79 
80  Info<< "Diffusion Number mean: " << meanDiNum
81  << " max: " << maxDiNum << endl;
82 
83  DiNum = maxDiNum;
84 }
85 
86 
87 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
88 
90 (
91  fvMesh& mesh,
92  autoPtr<solidThermo> thermoPtr
93 )
94 :
95  solver(mesh),
96 
97  thermoPtr_(thermoPtr),
98  thermo_(thermoPtr_()),
99 
100  T_(thermo_.T()),
101 
103 
104  DiNum(0),
105 
106  thermo(thermo_),
107  T(T_)
108 {
109  thermo.validate("solid", "h", "e");
110 
111  if (transient())
112  {
113  correctDiNum();
114  }
115  else if (LTS)
116  {
117  FatalError
118  << type()
119  << " solver does not support LTS, use 'steadyState' ddtScheme"
120  << exit(FatalError);
121  }
122 }
123 
124 
125 
127 :
128  solid(mesh, solidThermo::New(mesh))
129 {
130  // Read the controls
131  readControls();
132 }
133 
134 
135 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
136 
138 {}
139 
140 
141 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
142 
144 {
145  scalar deltaT = min(fvModels().maxDeltaT(), maxDeltaT_);
146 
147  if (DiNum > small)
148  {
149  deltaT = min(deltaT, maxDi/DiNum*runTime.deltaTValue());
150  }
151 
152  return deltaT;
153 }
154 
155 
157 {
158  // Read the controls
159  readControls();
160 
162 
163  // Update the mesh for topology change, mesh to mesh mapping
164  mesh_.update();
165 
166  if (transient())
167  {
168  correctDiNum();
169  }
170 }
171 
172 
174 {
175  if (pimple.firstIter() || pimple.moveMeshOuterCorrectors())
176  {
177  if (!mesh_.mover().solidBody())
178  {
180  << "Region " << name() << " of type " << type()
181  << " does not support non-solid body mesh motion"
182  << exit(FatalError);
183  }
184 
185  mesh_.move();
186  }
187 }
188 
189 
191 {
192  if (pimple.predictTransport())
193  {
194  thermophysicalTransport->predict();
195  }
196 }
197 
198 
200 {}
201 
202 
204 {
205  volScalarField& e = thermo_.he();
206  const volScalarField& rho = thermo_.rho();
207 
208  while (pimple.correctNonOrthogonal())
209  {
210  fvScalarMatrix eEqn
211  (
212  fvm::ddt(rho, e)
213  + thermophysicalTransport->divq(e)
214  ==
215  fvModels().source(rho, e)
216  );
217 
218  eEqn.relax();
219 
220  fvConstraints().constrain(eEqn);
221 
222  eEqn.solve();
223 
225  }
226 
227  thermo_.correct();
228 }
229 
230 
232 {}
233 
234 
236 {
237  if (pimple.correctTransport())
238  {
239  thermophysicalTransport->correct();
240  }
241 }
242 
243 
245 {}
246 
247 
248 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
const dictionary & controlDict() const
Return the control dict.
Definition: Time.H:265
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
void validate(const string &app, const word &) const
Check that the thermodynamics package is consistent.
Definition: basicThermo.C:331
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
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:604
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:101
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:56
Abstract base class for solid thermophysical transport models.
Abstract base class for run-time selectable region solvers.
Definition: solver.H:55
bool LTS
Switch for local time step transient operation.
Definition: solver.H:69
const Time & runTime
Time.
Definition: solver.H:97
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:203
virtual void prePredictor()
Called at the beginning of the PIMPLE loop.
Definition: solid.C:190
virtual void postSolve()
Called after the PIMPLE loop at the end of the time-step.
Definition: solid.C:244
virtual void moveMesh()
Called at the start of the PIMPLE loop to move the mesh.
Definition: solid.C:173
virtual scalar maxDeltaT() const
Return the current maximum time-step for stable solution.
Definition: solid.C:143
const solidThermo & thermo
Reference to the solid thermophysical properties.
Definition: solid.H:106
virtual void pressureCorrector()
Construct and solve the pressure equation in the PISO loop.
Definition: solid.C:231
virtual void postCorrector()
Correct the thermophysical transport modelling.
Definition: solid.C:235
virtual void momentumPredictor()
Construct and optionally solve the momentum equation.
Definition: solid.C:199
virtual void readControls()
Read controls.
Definition: solid.C:46
solid(fvMesh &mesh, autoPtr< solidThermo >)
Construct from region mesh and thermophysical properties.
Definition: solid.C:90
virtual void preSolve()
Called at the start of the time-step, before the PIMPLE loop.
Definition: solid.C:156
scalar maxDeltaT_
Definition: solid.H:64
virtual ~solid()
Destructor.
Definition: solid.C:137
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:306
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:105
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 > &)
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
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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