dynamicOversetFvMesh.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2022 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::dynamicOversetFvMesh
28 
29 Description
30  dynamicFvMesh with support for overset meshes.
31 
32 SourceFiles
33  dynamicOversetFvMesh.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef dynamicOversetFvMesh_H
38 #define dynamicOversetFvMesh_H
39 
41 #include "oversetFvMeshBase.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 /*---------------------------------------------------------------------------*\
49  Class dynamicOversetFvMesh Declaration
50 \*---------------------------------------------------------------------------*/
51 
53 :
55  public oversetFvMeshBase
56 {
57  // Private Member Functions
58 
59  //- No copy construct
61 
62  //- No copy assignment
63  void operator=(const dynamicOversetFvMesh&) = delete;
64 
65 
66 public:
67 
68  //- Runtime type information
69  TypeName("dynamicOversetFvMesh");
70 
71 
72  // Constructors
73 
74  //- Construct from IOobject
75  dynamicOversetFvMesh(const IOobject& io, const bool doInit=true);
76 
77 
78  //- Destructor
79  virtual ~dynamicOversetFvMesh();
80 
81 
82  // Member Functions
83 
84  //- Override ldu addressing
85  virtual const lduAddressing& lduAddr() const
86  {
88  }
89 
90  //- Override ldu addressing
91  virtual lduInterfacePtrsList interfaces() const
92  {
94  }
95 
96 
97  // Overset
98 
99  // Explicit interpolation
101  virtual void interpolate(scalarField& psi) const
102  {
103  Stencil::New(*this).interpolate<scalar>(*this, psi);
104  }
105 
106  virtual void interpolate(vectorField& psi) const
107  {
108  Stencil::New(*this).interpolate<vector>(*this, psi);
109  }
111  virtual void interpolate(sphericalTensorField& psi) const
112  {
114  (
115  *this
116  ).interpolate<sphericalTensor>(*this, psi);
117  }
118 
119  virtual void interpolate(symmTensorField& psi) const
120  {
121  Stencil::New(*this).interpolate<symmTensor>(*this, psi);
122  }
123 
124  virtual void interpolate(tensorField& psi) const
125  {
126  Stencil::New(*this).interpolate<tensor>(*this, psi);
127  }
129  virtual void interpolate(volScalarField& psi) const
130  {
131  Stencil::New(*this).interpolate<volScalarField>(psi);
132  }
134  virtual void interpolate(volVectorField& psi) const
135  {
136  Stencil::New(*this).interpolate<volVectorField>(psi);
137  }
139  virtual void interpolate(volSphericalTensorField& psi) const
140  {
141  Stencil::New(*this).interpolate
142  <
144  >(psi);
145  }
146 
147  virtual void interpolate(volSymmTensorField& psi) const
148  {
149  Stencil::New(*this).interpolate<volSymmTensorField>(psi);
150  }
151 
152  virtual void interpolate(volTensorField& psi) const
153  {
154  Stencil::New(*this).interpolate<volTensorField>(psi);
155  }
157 
158  // Implicit interpolation (matrix manipulation)
159 
160  //- Solve returning the solution statistics given convergence
161  //- tolerance. Use the given solver controls
163  (
164  fvMatrix<scalar>& m,
165  const dictionary& dict
166  ) const
167  {
168  return oversetFvMeshBase::solveOverset<scalar>(m, dict);
169  }
170 
171  //- Solve returning the solution statistics given convergence
172  //- tolerance. Use the given solver controls
174  (
175  fvMatrix<vector>& m,
176  const dictionary& dict
177  ) const
178  {
179  return oversetFvMeshBase::solveOverset<vector>(m, dict);
180  }
181 
182  //- Solve returning the solution statistics given convergence
183  //- tolerance. Use the given solver controls
185  (
188  ) const
189  {
190  return oversetFvMeshBase::solveOverset<symmTensor>
191  (m, dict);
192  }
193 
194  //- Solve returning the solution statistics given convergence
195  //- tolerance. Use the given solver controls
197  (
198  fvMatrix<tensor>& m,
199  const dictionary& dict
200  ) const
201  {
202  return oversetFvMeshBase::solveOverset<tensor>(m, dict);
203  }
204 
205 
206  //- Update the mesh for both mesh motion and topology change
207  virtual bool update();
208 
209  //- Write using given format, version and compression
210  virtual bool writeObject
211  (
212  IOstreamOption streamOpt,
213  const bool writeOnProc
214  ) const;
215 };
216 
217 
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 
220 } // End namespace Foam
221 
222 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223 
224 #endif
225 
226 // ************************************************************************* //
dictionary dict
virtual lduInterfacePtrsList interfaces() const
Override ldu addressing.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
static const Type & New(const Mesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
A simple container for options an IOstream can normally have.
virtual const lduAddressing & lduAddr() const
Override ldu addressing.
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
virtual bool update()
Update the mesh for both mesh motion and topology change.
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
dynamicFvMesh with support for overset meshes.
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc) const
Write using given format, version and compression.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Dynamic mesh able to handle multiple motion solvers. NOTE: If the word entry "solvers" is not found i...
TypeName("dynamicOversetFvMesh")
Runtime type information.
Support for overset functionality.
virtual const lduAddressing & lduAddr() const
Return ldu addressing. If active: is (extended)
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
const volScalarField & psi
The class contains the addressing required by the lduMatrix: upper, lower and losort.
virtual ~dynamicOversetFvMesh()
Destructor.
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
Tensor of scalars, i.e. Tensor<scalar>.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
virtual SolverPerformance< scalar > solve(fvMatrix< scalar > &m, const dictionary &dict) const
Solve returning the solution statistics given convergence tolerance. Use the given solver controls...
Namespace for OpenFOAM.
virtual void interpolate(scalarField &psi) const
Interpolate interpolationCells only. No bcs.