LagrangianmDdt.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 
26 #include "LagrangianmDdt.H"
27 #include "LagrangianDdtScheme.H"
28 #include "LagrangianMesh.H"
29 #include "LagrangianSubFields.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 template<class Type>
35 (
36  const dimensionSet& mDims,
38  const bool instantaneousDdt
39 )
40 {
41  const LagrangianMesh& mesh = psi.mesh().mesh();
42 
43  return
45  (
46  mesh,
47  mesh.schemes().ddt("ddt(" + psi.name() + ')')
48  ).ref().LagrangianmInitDdt(mDims, psi, instantaneousDdt);
49 }
50 
51 
52 template<class Type, template<class> class PrimitiveField>
54 (
55  const dimensionSet& mDims,
57  const bool instantaneousDdt
58 )
59 {
60  return Lagrangianm::initDdt(mDims, toSubField(psi)(), instantaneousDdt);
61 }
62 
63 
64 template<class Type>
66 (
67  const LagrangianSubScalarField& deltaT,
68  const dimensionSet& mDims,
70 )
71 {
72  const LagrangianMesh& mesh = psi.mesh().mesh();
73 
74  return
76  (
77  mesh,
78  mesh.schemes().ddt("ddt(" + psi.name() + ')')
79  ).ref().LagrangianmNoDdt(deltaT, mDims, psi);
80 }
81 
82 
83 template<class Type, template<class> class PrimitiveField>
85 (
86  const LagrangianSubScalarField& deltaT,
87  const dimensionSet& mDims,
88  const LagrangianSubField<Type, PrimitiveField>& psi
89 )
90 {
91  return Lagrangianm::noDdt(deltaT, mDims, toSubField(psi)());
92 }
93 
94 
95 template<class Type>
97 (
98  const LagrangianSubScalarField& deltaT,
99  LagrangianSubSubField<Type>& psi
100 )
101 {
102  const LagrangianMesh& mesh = psi.mesh().mesh();
103 
104  return
106  (
107  mesh,
108  mesh.schemes().ddt("ddt(" + psi.name() + ')')
109  ).ref().LagrangianmDdt(deltaT, psi);
110 }
111 
112 
113 template<class Type>
115 (
116  const LagrangianSubScalarField& deltaT,
118  LagrangianSubSubField<Type>& psi
119 )
120 {
121  const LagrangianMesh& mesh = psi.mesh().mesh();
122 
123  return
125  (
126  mesh,
127  mesh.schemes().ddt("ddt(" + m.name() + ',' + psi.name() + ')')
128  ).ref().LagrangianmDdt(deltaT, m, psi);
129 }
130 
131 
132 template<class Type>
134 (
135  const LagrangianSubScalarField& deltaT,
136  const LagrangianSubSubField<Type>& psi
137 )
138 {
139  return Lagrangian::ddtScheme<Type>::Lagrangianmddt(deltaT, psi);
140 }
141 
142 
143 template<class Type>
145 (
146  const LagrangianSubScalarField& deltaT,
148  const LagrangianSubSubField<Type>& psi
149 )
150 {
151  return Lagrangian::ddtScheme<Type>::Lagrangianmddt(deltaT, m, psi);
152 }
153 
154 
155 template<class Type>
157 (
158  const LagrangianSubScalarField& deltaT,
159  const LagrangianSubSubField<Type>& psi
160 )
161 {
162  return Lagrangian::ddtScheme<Type>::Lagrangianmddt0(deltaT, psi);
163 }
164 
165 
166 template<class Type>
168 (
169  const LagrangianSubScalarField& deltaT,
171  const LagrangianSubSubField<Type>& psi
172 )
173 {
174  return Lagrangian::ddtScheme<Type>::Lagrangianmddt0(deltaT, m, psi);
175 }
176 
177 
178 // ************************************************************************* //
Functions for calculating the time derivative for a Lagrangian equation.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Class containing Lagrangian geometry and topology.
Dimension set for the base types.
Definition: dimensionSet.H:125
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1799
ITstream & ddt(const word &name) const
Definition: fvSchemes.C:400
A class for managing temporary objects.
Definition: tmp.H:55
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
const volScalarField & psi
tmp< LagrangianEqn< Type > > Ddt(const LagrangianSubScalarField &deltaT, LagrangianSubSubField< Type > &psi)
tmp< LagrangianEqn< Type > > ddt0(const LagrangianSubScalarField &deltaT, const LagrangianSubSubField< Type > &psi)
tmp< LagrangianEqn< Type > > ddt(const LagrangianSubScalarField &deltaT, const LagrangianSubSubField< Type > &psi)
bool initDdt(const dimensionSet &mDims, const LagrangianSubSubField< Type > &psi, const bool instantaneousDdt=false)
tmp< LagrangianEqn< Type > > noDdt(const LagrangianSubScalarField &deltaT, const dimensionSet &mDims, const LagrangianSubSubField< Type > &psi)
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
LagrangianSubSubField< scalar > LagrangianSubScalarSubField
LagrangianSubField< scalar > LagrangianSubScalarField
tmp< DimensionedField< Type, GeoMesh, SubField > > toSubField(const DimensionedField< Type, GeoMesh, Field > &)
Return a temporary sub-field from a reference to a field.