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-2025 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 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
45 
47 {
48  return runTime.controlDict().modified();
49 }
50 
51 
53 {
54  solver::read();
55 
56  maxDeltaT_ =
57  runTime.controlDict().found("maxDeltaT")
58  ? runTime.controlDict().lookup<scalar>("maxDeltaT", runTime.userUnits())
59  : vGreat;
60 
61  return true;
62 }
63 
64 
65 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
66 
68 (
69  fvMesh& mesh,
70  autoPtr<solidThermo> thermoPtr
71 )
72 :
73  solver(mesh),
74 
75  thermoPtr_(thermoPtr),
76  thermo_(thermoPtr_()),
77 
78  T_(thermo_.T()),
79 
81 
82  thermo(thermo_),
83  T(T_)
84 {
85  thermo.validate("solid", "h", "e");
86 
87  if (LTS)
88  {
90  << type()
91  << " solver does not support LTS, use 'steadyState' ddtScheme"
92  << exit(FatalError);
93  }
94 }
95 
96 
97 
99 :
101 {
102  // Read the controls
103  read();
104 }
105 
106 
107 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
108 
110 {}
111 
112 
113 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
114 
116 {
117  return min(fvModels().maxDeltaT(), maxDeltaT_);
118 }
119 
120 
122 {
124 
125  // Update the mesh for topology change, mesh to mesh mapping
126  mesh_.update();
127 }
128 
129 
131 {
132  if (pimple.firstIter() || pimple.moveMeshOuterCorrectors())
133  {
134  if (!mesh_.mover().solidBody())
135  {
137  << "Region " << name() << " of type " << type()
138  << " does not support non-solid body mesh motion"
139  << exit(FatalError);
140  }
141 
142  mesh_.move();
143  }
144 }
145 
146 
148 {}
149 
150 
152 {}
153 
154 
156 {}
157 
158 
160 {
161  thermophysicalTransport->predict();
162 }
163 
164 
166 {}
167 
168 
170 {
171  volScalarField& e = thermo_.he();
172  const volScalarField& rho = thermo_.rho();
173 
174  while (pimple.correctNonOrthogonal())
175  {
176  fvScalarMatrix eEqn
177  (
178  fvm::ddt(rho, e)
179  + thermophysicalTransport->divq(e)
180  ==
181  fvModels().source(rho, e)
182  );
183 
184  eEqn.relax();
185 
186  fvConstraints().constrain(eEqn);
187 
188  eEqn.solve();
189 
191 
192  thermo_.correct();
193  }
194 }
195 
196 
198 {}
199 
200 
202 {}
203 
204 
206 {
207  thermophysicalTransport->correct();
208 }
209 
210 
212 {}
213 
214 
215 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
const IOdictionary & controlDict() const
Return the control dict.
Definition: Time.H:269
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:308
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:96
Finite volume models.
Definition: fvModels.H:65
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: fvModels.C:260
virtual bool modified() const
Return true if the object's file (or files for objectRegistry)
Base-class for solid thermodynamic properties.
Definition: solidThermo.H:59
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
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:169
virtual void momentumTransportCorrector()
Correct the momentum transport.
Definition: solid.C:201
virtual void prePredictor()
Called at the beginning of the PIMPLE loop.
Definition: solid.C:151
virtual void postSolve()
Called after the PIMPLE loop at the end of the time-step.
Definition: solid.C:211
virtual void momentumTransportPredictor()
Predict the momentum transport.
Definition: solid.C:155
virtual bool dependenciesModified() const
Return true if the solver's dependencies have been modified.
Definition: solid.C:46
virtual void moveMesh()
Called at the start of the PIMPLE loop to move the mesh.
Definition: solid.C:130
virtual scalar maxDeltaT() const
Return the current maximum time-step for stable solution.
Definition: solid.C:115
virtual void motionCorrector()
Corrections that follow mesh motion.
Definition: solid.C:147
const solidThermo & thermo
Reference to the solid thermophysical properties.
Definition: solid.H:94
virtual void pressureCorrector()
Construct and solve the pressure equation in the PISO loop.
Definition: solid.C:197
virtual void momentumPredictor()
Construct and optionally solve the momentum equation.
Definition: solid.C:165
virtual void thermophysicalTransportCorrector()
Correct the thermophysical transport.
Definition: solid.C:205
solid(fvMesh &mesh, autoPtr< solidThermo >)
Construct from region mesh and thermophysical properties.
Definition: solid.C:68
virtual void preSolve()
Called at the start of the time-step, before the PIMPLE loop.
Definition: solid.C:121
virtual void thermophysicalTransportPredictor()
Predict thermophysical transport.
Definition: solid.C:159
virtual bool read()
Read controls.
Definition: solid.C:52
virtual ~solid()
Destructor.
Definition: solid.C:109
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
pimpleControl pimple(mesh)
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
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)
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
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
tmp< DimensionedField< TypeR, GeoMesh, Field > > New(const tmp< DimensionedField< TypeR, GeoMesh, Field >> &tdf1, const word &name, const dimensionSet &dimensions)
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