LagrangianSpScheme.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 Class
25  Foam::Lagrangian::SpScheme
26 
27 Description
28  Abstract base class for Lagrangian Sp schemes
29 
30 SourceFiles
31  LagrangianSpScheme.C
32  LagrangianSpSchemes.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef LagrangianSpScheme_H
37 #define LagrangianSpScheme_H
38 
39 #include "tmp.H"
40 #include "dimensionedType.H"
41 #include "LagrangianFieldsFwd.H"
42 #include "LagrangianSubFieldsFwd.H"
43 #include "runTimeSelectionTables.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 class LagrangianMesh;
51 
52 template<class Type>
53 class LagrangianEqn;
54 
55 namespace Lagrangian
56 {
57 
58 /*---------------------------------------------------------------------------*\
59  Class SpScheme Declaration
60 \*---------------------------------------------------------------------------*/
61 
62 template<class Type, class SpType>
63 class SpScheme
64 :
65  public tmp<SpScheme<Type, SpType>>::refCount
66 {
67 protected:
68 
69  // Protected Data
70 
71  //- Reference to the mesh
72  const LagrangianMesh& mesh_;
73 
74 
75  // Protected Member Functions
76 
77  //- Generate the inner product of a scalar coefficient with a field
79  (
82  );
83 
84  //- Generate the inner product of a tensor coefficient with a field
85  static tmp
86  <
88  <
90  >
91  > inner
92  (
95  );
96 
97 
98 public:
99 
100  //- Runtime type information
101  TypeName("SpScheme");
102 
103 
104  //- Declare run-time constructor selection tables
106  (
107  tmp,
108  SpScheme,
109  Istream,
110  (const LagrangianMesh& mesh, Istream& is),
111  (mesh, is)
112  );
113 
114 
115  // Constructors
116 
117  //- Construct from a mesh
119  :
120  mesh_(mesh)
121  {}
122 
123  //- Construct from a mesh and a stream
125  :
126  mesh_(mesh)
127  {}
128 
129  //- Disallow default bitwise copy construction
130  SpScheme(const SpScheme&) = delete;
131 
132 
133  // Selectors
134 
135  //- Return a pointer to a new SpScheme
137  (
138  const LagrangianMesh& mesh,
139  Istream& is
140  );
141 
142 
143  //- Destructor
144  virtual ~SpScheme();
145 
146 
147  // Member Functions
148 
149  //- Return mesh reference
150  const LagrangianMesh& mesh() const
151  {
152  return mesh_;
153  }
154 
155  //- Return the source matrix
157  (
160  ) = 0;
161 
162 
163  // Member Operators
164 
165  //- Disallow default bitwise assignment
166  void operator=(const SpScheme&) = delete;
167 };
168 
169 
170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171 
172 } // End namespace Lagrangian
173 } // End namespace Foam
174 
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 
177 #define defineLagrangianSpScheme(Type, SpType) \
178  \
179  namespace Foam \
180  { \
181  namespace Lagrangian \
182  { \
183  typedef SpScheme<Type, SpType> Type##SpType##SpScheme; \
184  \
185  defineNamedTemplateTypeNameAndDebug \
186  ( \
187  Type##SpType##SpScheme, \
188  0 \
189  ); \
190  \
191  defineTemplateRunTimeSelectionTable \
192  ( \
193  Type##SpType##SpScheme, \
194  Istream \
195  ); \
196  } \
197  }
198 
199 
200 #define makeLagrangianSpScheme(Type, SpType, SpSchemeType) \
201  \
202  namespace Foam \
203  { \
204  namespace Lagrangian \
205  { \
206  typedef SpScheme<Type, SpType> Type##SpType##SpScheme; \
207  \
208  namespace SpSchemes \
209  { \
210  typedef SpSchemeType<Type, SpType> Type##SpType##SpSchemeType; \
211  \
212  defineNamedTemplateTypeNameAndDebug \
213  ( \
214  Type##SpType##SpSchemeType, \
215  0 \
216  ); \
217  \
218  addToRunTimeSelectionTable \
219  ( \
220  Type##SpType##SpScheme, \
221  Type##SpType##SpSchemeType, \
222  Istream \
223  ); \
224  } \
225  } \
226  }
227 
228 
229 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230 
231 #ifdef NoRepository
232  #include "LagrangianSpScheme.C"
233 #endif
234 
235 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236 
237 #endif
238 
239 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
Class containing Lagrangian geometry and topology.
Abstract base class for Lagrangian Sp schemes.
declareRunTimeSelectionTable(tmp, SpScheme, Istream,(const LagrangianMesh &mesh, Istream &is),(mesh, is))
Declare run-time constructor selection tables.
void operator=(const SpScheme &)=delete
Disallow default bitwise assignment.
TypeName("SpScheme")
Runtime type information.
static tmp< SpScheme< Type, SpType > > New(const LagrangianMesh &mesh, Istream &is)
Return a pointer to a new SpScheme.
static tmp< LagrangianSubField< Type > > inner(const LagrangianSubScalarField &Sp, const LagrangianSubSubField< Type > &psi)
Generate the inner product of a scalar coefficient with a field.
SpScheme(const LagrangianMesh &mesh)
Construct from a mesh.
virtual ~SpScheme()
Destructor.
const LagrangianMesh & mesh_
Reference to the mesh.
virtual tmp< LagrangianEqn< Type > > LagrangianmSp(const LagrangianSubField< SpType > &Sp, const LagrangianSubSubField< Type > &psi)=0
Return the source matrix.
const LagrangianMesh & mesh() const
Return mesh reference.
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) - 2 >::type type
Definition: products.H:115
Reference counter for various OpenFOAM components.
Definition: refCount.H:50
A class for managing temporary objects.
Definition: tmp.H:55
const volScalarField & psi
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
Definition: fvcSup.C:67
Namespace for OpenFOAM.
Macros to ease declaration of run-time selection tables.