LagrangianmDdt.H
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 InNamespace
25  Foam::Lagrangianm
26 
27 Description
28  Functions for calculating the time derivative for a Lagrangian equation
29 
30 SourceFiles
31  LagrangianmDdt.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef LagrangianmDdt_H
36 #define LagrangianmDdt_H
37 
38 #include "LagrangianFieldsFwd.H"
39 #include "LagrangianSubFieldsFwd.H"
40 #include "dimensionSet.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 template<class Type>
48 class LagrangianEqn;
49 
50 /*---------------------------------------------------------------------------*\
51  Namespace Lagrangianm functions Declaration
52 \*---------------------------------------------------------------------------*/
53 
54 namespace Lagrangianm
55 {
56  template<class Type>
57  bool initDdt
58  (
59  const dimensionSet& mDims,
61  const bool instantaneousDdt = false
62  );
63 
64  template<class Type, template<class> class PrimitiveField>
65  bool initDdt
66  (
67  const dimensionSet& mDims,
69  const bool instantaneousDdt = false
70  );
71 
72  template<class Type>
74  (
75  const LagrangianSubScalarField& deltaT,
76  const dimensionSet& mDims,
78  );
79 
80  template<class Type, template<class> class PrimitiveField>
82  (
83  const LagrangianSubScalarField& deltaT,
84  const dimensionSet& mDims,
86  );
87 
88  template<class Type>
90  (
91  const LagrangianSubScalarField& deltaT,
93  );
94 
95  template<class Type>
97  (
98  const LagrangianSubScalarField& deltaT,
101  );
102 
103  template<class Type>
105  (
106  const LagrangianSubScalarField& deltaT,
108  );
109 
110  template<class Type>
112  (
113  const LagrangianSubScalarField& deltaT,
116  );
117 
118  template<class Type>
120  (
121  const LagrangianSubScalarField& deltaT,
123  );
124 
125  template<class Type>
127  (
128  const LagrangianSubScalarField& deltaT,
131  );
132 }
133 
134 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
135 
136 } // End namespace Foam
137 
138 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139 
140 #ifdef NoRepository
141  #include "LagrangianmDdt.C"
142 #endif
143 
144 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
145 
146 #endif
147 
148 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Dimension set for the base types.
Definition: dimensionSet.H:125
A class for managing temporary objects.
Definition: tmp.H:55
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)
Namespace for OpenFOAM.