CrankNicolson_LagrangianDdtScheme.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) 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 
27 #include "LagrangianFields.H"
28 #include "LagrangianEqn.H"
29 #include "LagrangianModels.H"
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 template<class Type>
37 (
38  const LagrangianSubSubField<Type>& psi
39 ) const
40 {
41  return typedName("S0(" + LagrangianModel::fieldName(psi) + ")");
42 }
43 
44 
45 template<class Type>
47 (
48  const LagrangianSubSubField<Type>& psi
49 ) const
50 {
51  return typedName("deltaTSp0(" + LagrangianModel::fieldName(psi) + ")");
52 }
53 
54 
55 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
56 
57 template<class Type>
59 (
60  const dimensionSet& mDims,
62  const bool instantaneousDdt
63 )
64 {
65  const LagrangianSubMesh& subMesh = psi.mesh();
66  const LagrangianMesh& mesh = subMesh.mesh();
67 
68  const word S0Name = this->S0Name(psi);
69 
71  {
72  regIOobject::store
73  (
75  (
76  IOobject
77  (
78  S0Name,
79  mesh.time().name(),
80  mesh,
81  IOobject::READ_IF_PRESENT,
82  IOobject::AUTO_WRITE,
83  false
84  ),
85  mesh,
86  dimensioned<Type>(mDims*psi.dimensions()/dimTime, Zero),
87  wordList
88  (
89  mesh.boundary().size(),
91  ),
92  wordList::null(),
93  LagrangianModels::New(mesh).modelTypeFieldSourceTypes
94  <
99  >()
100  )
101  );
102  }
103 
104  if (!instantaneousDdt) return true;
105 
106  const word deltaTSp0Name = this->deltaTSp0Name(psi);
107 
108  if (!mesh.foundObject<LagrangianDynamicField<scalar>>(deltaTSp0Name))
109  {
110  regIOobject::store
111  (
113  (
114  IOobject
115  (
116  deltaTSp0Name,
117  mesh.time().name(),
118  mesh,
119  IOobject::READ_IF_PRESENT,
120  IOobject::AUTO_WRITE,
121  false
122  ),
123  mesh,
124  dimensioned<scalar>(mDims, Zero),
125  wordList
126  (
127  mesh.boundary().size(),
129  ),
130  wordList::null(),
131  LagrangianModels::New(mesh).modelTypeFieldSourceTypes
132  <
137  >()
138  )
139  );
140  }
141 
142  return true;
143 }
144 
145 
146 template<class Type>
149 (
150  const LagrangianSubScalarField& deltaT,
151  const dimensionSet& mDims,
153 )
154 {
155  const LagrangianSubMesh& subMesh = deltaT.mesh();
156  const LagrangianMesh& mesh = subMesh.mesh();
157 
158  const word S0Name = this->S0Name(psi);
159  const word deltaTSp0Name = this->deltaTSp0Name(psi);
160 
163  LagrangianDynamicField<scalar>& deltaTSp0 =
166  : NullObjectNonConstRef<LagrangianDynamicField<scalar>>();
167 
168  const bool isNone = subMesh.group() == LagrangianGroup::none;
169 
171  (
173  (
174  !isNone ? deltaT/2 : tmp<LagrangianSubScalarField>(deltaT),
175  psi,
176  deltaTSp0,
177  S0
178  )
179  );
180 
181  if (!isNone) tEqn.ref().Su += subMesh.sub(S0);
182 
183  return tEqn;
184 }
185 
186 
187 template<class Type>
190 (
191  const LagrangianSubScalarField& deltaT,
193 )
194 {
197 
198  tEqn.ref().deltaTSp += 1;
199  tEqn.ref().deltaTSu -= psi.oldTime();
200 
201  // Set the equation name and the non-const field reference
202  tEqn->op(LagrangianEqn<Type>(psi));
203 
204  return tEqn;
205 }
206 
207 
208 template<class Type>
211 (
212  const LagrangianSubScalarField& deltaT,
215 )
216 {
219 
220  tEqn.ref().deltaTSp += m;
221  tEqn.ref().deltaTSu -= m.oldTime()*psi.oldTime();
222 
223  // Set the equation name and the non-const field reference
224  tEqn->op(LagrangianEqn<Type>(psi));
225 
226  return tEqn;
227 }
228 
229 
230 template<class Type>
233 (
235 )
236 {
237  const LagrangianSubMesh& subMesh = psi.mesh();
238  const LagrangianMesh& mesh = subMesh.mesh();
239 
240  return
241  - subMesh.sub
242  (
244  (
245  S0Name(psi)
246  )
247  )
248  /subMesh.sub
249  (
251  (
252  deltaTSp0Name(psi)
253  )
254  );
255 }
256 
257 
258 template<class Type>
261 (
264 )
265 {
266  // This is an approximation if the time-coefficient (deltaTSp0) is not
267  // "just" m. Some models may add to the time-coefficient; a virtual mass
268  // model, for example, adds the added mass. The time-derivative of these
269  // additions is not known and cannot therefore be factored out. This
270  // approximation assumes that those additional terms are proportional to m
271  // and neglects the time-derivative of the coefficient. For a
272  // constant-coefficient virtual mass and constant relative density this
273  // approach is exact.
274 
275  return m*LagrangiancDdt(psi);
276 }
277 
278 
279 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
const Mesh & mesh() const
Return mesh.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
This class stores the coefficients of a Lagrangian equation, and facilitates solving that equation an...
Definition: LagrangianEqn.H:56
Base class for Lagrangian injections. Minimal wrapper over LagrangianSource. Implements some utility ...
Class containing Lagrangian geometry and topology.
Base class for Lagrangian sources. Minimal wrapper over LagrangianModel that provides an interface to...
Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate...
SubList< Type > sub(const List< Type > &list) const
Return a sub-list corresponding to this sub-mesh.
LagrangianGroup group() const
Return the group.
const LagrangianMesh & mesh() const
Return the mesh.
Basic first-order implicit CrankNicolson Lagrangian ddt scheme.
virtual tmp< LagrangianEqn< Type > > LagrangianmNoDdt(const LagrangianSubScalarField &deltaT, const dimensionSet &mDims, const LagrangianSubSubField< Type > &psi)
Return the no-time-derivative matrix.
virtual bool LagrangianmInitDdt(const dimensionSet &mDims, const LagrangianSubSubField< Type > &psi, const bool instantaneousDdt)
Initialise time-derivative information.
virtual tmp< LagrangianEqn< Type > > LagrangianmDdt(const LagrangianSubScalarField &deltaT, LagrangianSubSubField< Type > &psi)
Return the time-derivative matrix.
virtual tmp< LagrangianSubField< Type > > LagrangiancDdt(const LagrangianSubSubField< Type > &psi)
Return the instantaneous time-derivative.
This source condition provides a NaN value.
const Field0Type & oldTime() const
Return the old-time field.
Definition: OldTimeField.C:322
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
A calculated boundary condition for Lagrangian. This condition does not provide any evaluation functi...
Dimension set for the base types.
Definition: dimensionSet.H:125
Generic dimensioned Type class.
const word & name() const
Return const reference to name.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:962
This source condition retains the internal value.
Type & lookupObjectRef(const word &name) const
Lookup and return the object reference of the given Type.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
bool foundObject(const word &name) const
Is the named Type in registry.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
const volScalarField & psi
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
static const zero Zero
Definition: zero.H:97
const dimensionSet dimless
const dimensionSet dimTime
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:181