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) 2015-2020 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
26Class
27 Foam::dynamicOversetFvMesh
28
29Description
30 dynamicFvMesh with support for overset meshes.
31
32SourceFiles
33 dynamicOversetFvMesh.C
34
35\*---------------------------------------------------------------------------*/
36
37#ifndef dynamicOversetFvMesh_H
38#define dynamicOversetFvMesh_H
39
41#include "labelIOList.H"
44#include "volFields.H"
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48namespace Foam
49{
50
51class mapDistribute;
52class lduPrimitiveProcessorInterface;
53class GAMGAgglomeration;
54
55/*---------------------------------------------------------------------------*\
56 Class dynamicOversetFvMesh Declaration
57\*---------------------------------------------------------------------------*/
60:
62{
63protected:
64
65 // Protected data
66
67 //- Select base addressing (false) or locally stored extended
68 // lduAddressing (true)
69 mutable bool active_;
70
71 //- Extended addressing (extended with local interpolation stencils)
73
74 //- Added (processor)lduInterfaces for remote bits of stencil.
75 //PtrList<const lduInterface> remoteStencilInterfaces_;
78
79 //- Interfaces for above mesh. Contains both original and
80 // above added processorLduInterfaces
82
83 //- Corresponding faces (in above lduPtr) to the stencil
85
86 //- Corresponding patches (in above lduPtr) to the stencil
88
89 //- From old to new face labels
91
92
93 // Protected Member Functions
94
95 //- Calculate the extended lduAddressing
96 virtual bool updateAddressing() const;
97
98 //- Debug: print matrix
99 template<class Type>
100 void write
101 (
102 Ostream&,
103 const fvMatrix<Type>&,
104 const lduAddressing&,
106 ) const;
107
108 //- Explicit interpolation of acceptor cells from donor cells
109 template<class T>
110 void interpolate(Field<T>& psi) const;
111
112 //- Explicit interpolation of acceptor cells from donor cells with
113 // boundary condition handling
114 template<class GeoField>
115 void interpolate(GeoField& psi) const;
116
117 //- Helper: strip off trailing _0
118 static word baseName(const word& name);
119
120 //- Explicit interpolation of all registered fields. Excludes
121 // selected fields (and their old-time fields)
122 template<class GeoField>
123 void interpolate(const wordHashSet& suppressed);
124
125 //- Freeze values at holes
126 //template<class Type>
127 //void freezeHoles(fvMatrix<Type>&) const;
128
129 //- Get scalar interfaces of certain type
130 //template<class GeoField, class PatchType>
131 //lduInterfaceFieldPtrsList scalarInterfaces(const GeoField& psi) const;
132
133 //- Determine normalisation for interpolation. This equals the
134 // original diagonal except stabilised for zero diagonals (possible
135 // in hole cells)
136 template<class Type>
138
139 //- Add interpolation to matrix (coefficients)
140 template<class Type>
141 void addInterpolation(fvMatrix<Type>&, const scalarField& norm) const;
142
143 //- Solve given dictionary with settings
144 template<class Type>
146
147 //- Debug: correct coupled bc
148 template<class GeoField>
149 static void correctCoupledBoundaryConditions(GeoField& fld);
150
151 //- Average norm of valid neighbours
152 scalar cellAverage
153 (
154 const labelList& types,
155 const labelList& nbrTypes,
156 const scalarField& norm,
157 const scalarField& nbrNorm,
158 const label celli,
159 bitSet& isFront
160 ) const;
161
162 //- Debug: dump agglomeration
164 (
165 const GAMGAgglomeration& agglom
166 ) const;
167
168
169private:
170
171 // Private Member Functions
172
173 //- No copy construct
175
176 //- No copy assignment
177 void operator=(const dynamicOversetFvMesh&) = delete;
178
179
180public:
181
182 //- Runtime type information
183 TypeName("dynamicOversetFvMesh");
184
185
186 // Constructors
187
188 //- Construct from IOobject
189 dynamicOversetFvMesh(const IOobject& io, const bool doInit=true);
190
191
192 //- Destructor
193 virtual ~dynamicOversetFvMesh();
194
195
196 // Member Functions
197
198
199 // Extended addressing
200
201 //- Return extended ldu addressing
203
204 //- Return ldu addressing. If active: is (extended)
205 // primitiveLduAddr
206 virtual const lduAddressing& lduAddr() const;
207
208 //- Return a list of pointers for each patch
209 // with only those pointing to interfaces being set. If active:
210 // return additional remoteStencilInterfaces_
211 virtual lduInterfacePtrsList interfaces() const;
212
213 //- Return old to new face addressing
214 const labelList& reverseFaceMap() const
215 {
216 return reverseFaceMap_;
217 }
218
219 //- Return true if using extended addressing
220 bool active() const
221 {
222 return active_;
223 }
224
225 //- Enable/disable extended addressing
226 void active(const bool f) const
227 {
228 active_ = f;
229
230 if (active_)
231 {
232 DebugInfo<< "Switching to extended addressing with nFaces:"
234 << endl;
235 }
236 else
237 {
238 DebugInfo<< "Switching to base addressing with nFaces:"
240 << endl;
241 }
242 }
243
244
245 // Overset
246
247 // Explicit interpolation
249 virtual void interpolate(scalarField& psi) const
250 {
251 interpolate<scalar>(psi);
252 }
254 virtual void interpolate(vectorField& psi) const
255 {
256 interpolate<vector>(psi);
257 }
259 virtual void interpolate(sphericalTensorField& psi) const
260 {
261 interpolate<sphericalTensor>(psi);
262 }
264 virtual void interpolate(symmTensorField& psi) const
265 {
266 interpolate<symmTensor>(psi);
267 }
269 virtual void interpolate(tensorField& psi) const
270 {
271 interpolate<tensor>(psi);
272 }
274 virtual void interpolate(volScalarField& psi) const
275 {
276 interpolate<volScalarField>(psi);
277 }
279 virtual void interpolate(volVectorField& psi) const
280 {
281 interpolate<volVectorField>(psi);
282 }
284 virtual void interpolate(volSphericalTensorField& psi) const
285 {
286 interpolate<volSphericalTensorField>(psi);
287 }
289 virtual void interpolate(volSymmTensorField& psi) const
290 {
291 interpolate<volSymmTensorField>(psi);
292 }
294 virtual void interpolate(volTensorField& psi) const
295 {
296 interpolate<volTensorField>(psi);
297 }
298
299
300 // Implicit interpolation (matrix manipulation)
301
302 //- Solve returning the solution statistics given convergence
303 // tolerance. Use the given solver controls
305 (
307 const dictionary& dict
308 ) const
309 {
310 return solve<scalar>(m, dict);
311 }
312
313 //- Solve returning the solution statistics given convergence
314 // tolerance. Use the given solver controls
316 (
318 const dictionary& dict
319 ) const
320 {
321 return solve<vector>(m, dict);
322 }
323
324 //- Solve returning the solution statistics given convergence
325 // tolerance. Use the given solver controls
327 (
329 const dictionary& dict
330 ) const
331 {
332 return solve<symmTensor>(m, dict);
333 }
334
335 //- Solve returning the solution statistics given convergence
336 // tolerance. Use the given solver controls
338 (
340 const dictionary& dict
341 ) const
342 {
343 return solve<tensor>(m, dict);
344 }
345
346
347 //- Initialise all non-demand-driven data
348 virtual bool init(const bool doInit);
349
350 //- Update the mesh for both mesh motion and topology change
351 virtual bool update();
352
353 //- Update fields when mesh is updated
354 virtual bool interpolateFields();
355
356 //- Write using stream options
357 virtual bool writeObject
358 (
359 IOstreamOption streamOpt,
360 const bool valid
361 ) const;
362
363 //- Debug: check halo swap is ok
364 template<class GeoField>
365 static void checkCoupledBC(const GeoField& fld);
366
367 //- Correct boundary conditions of certain type (typeOnly = true)
368 // or explicitly not of the type (typeOnly = false)
369 template<class GeoField, class PatchType>
370 static void correctBoundaryConditions
371 (
372 typename GeoField::Boundary& bfld,
373 const bool typeOnly
374 );
375};
376
377
378// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
379
380} // End namespace Foam
381
382// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
383
384#ifdef NoRepository
386#endif
387
388// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
389
390#endif
391
392// ************************************************************************* //
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Generic templated field type.
Definition: Field.H:82
Geometric agglomerated algebraic multigrid agglomeration class.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
The IOstreamOption is a simple container for options an IOstream can normally have.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Dynamic mesh able to handle multiple motion solvers. NOTE: If the word entry "solvers" is not found i...
dynamicFvMesh with support for overset meshes.
static word baseName(const word &name)
Helper: strip off trailing _0.
bool active_
Select base addressing (false) or locally stored extended.
void writeAgglomeration(const GAMGAgglomeration &agglom) const
Debug: dump agglomeration.
virtual void interpolate(volTensorField &psi) const
Interpolate interpolationCells only.
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
virtual bool interpolateFields()
Update fields when mesh is updated.
virtual ~dynamicOversetFvMesh()
Destructor.
virtual void interpolate(volVectorField &psi) const
Interpolate interpolationCells only.
virtual const lduAddressing & lduAddr() const
Return ldu addressing. If active: is (extended)
virtual void interpolate(tensorField &psi) const
Interpolate interpolationCells only. No bcs.
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
PtrList< const lduPrimitiveProcessorInterface > remoteStencilInterfaces_
Added (processor)lduInterfaces for remote bits of stencil.
autoPtr< fvMeshPrimitiveLduAddressing > lduPtr_
Extended addressing (extended with local interpolation stencils)
virtual void interpolate(sphericalTensorField &psi) const
Interpolate interpolationCells only. No bcs.
static void checkCoupledBC(const GeoField &fld)
Debug: check halo swap is ok.
TypeName("dynamicOversetFvMesh")
Runtime type information.
labelListList stencilFaces_
Corresponding faces (in above lduPtr) to the stencil.
void active(const bool f) const
Enable/disable extended addressing.
virtual void interpolate(scalarField &psi) const
Interpolate interpolationCells only. No bcs.
virtual void interpolate(symmTensorField &psi) const
Interpolate interpolationCells only. No bcs.
void addInterpolation(fvMatrix< Type > &, const scalarField &norm) const
Add interpolation to matrix (coefficients)
virtual bool updateAddressing() const
Calculate the extended lduAddressing.
labelListList stencilPatches_
Corresponding patches (in above lduPtr) to the stencil.
virtual void interpolate(vectorField &psi) const
Interpolate interpolationCells only. No bcs.
const fvMeshPrimitiveLduAddressing & primitiveLduAddr() const
Return extended ldu addressing.
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write using stream options.
virtual bool update()
Update the mesh for both mesh motion and topology change.
virtual void interpolate(volScalarField &psi) const
Interpolate interpolationCells only.
virtual void interpolate(volSymmTensorField &psi) const
Interpolate interpolationCells only.
virtual void interpolate(volSphericalTensorField &psi) const
Interpolate interpolationCells only.
lduInterfacePtrsList allInterfaces_
Interfaces for above mesh. Contains both original and.
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &) const
Solve given dictionary with settings.
void interpolate(Field< T > &psi) const
Explicit interpolation of acceptor cells from donor cells.
tmp< scalarField > normalisation(const fvMatrix< Type > &m) const
Freeze values at holes.
static void correctCoupledBoundaryConditions(GeoField &fld)
Debug: correct coupled bc.
bool active() const
Return true if using extended addressing.
const labelList & reverseFaceMap() const
Return old to new face addressing.
scalar cellAverage(const labelList &types, const labelList &nbrTypes, const scalarField &norm, const scalarField &nbrNorm, const label celli, bitSet &isFront) const
Average norm of valid neighbours.
labelList reverseFaceMap_
From old to new face labels.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
Variant of fvMeshLduAddressing that contains addressing instead of slices.
virtual const labelUList & lowerAddr() const noexcept
Return lower addressing (i.e. lower label = upper triangle)
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:718
const word & name() const
Return reference to name.
Definition: fvMesh.H:310
The class contains the addressing required by the lduMatrix: upper, lower and losort.
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
const volScalarField & psi
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
#define DebugInfo
Report an information message using Foam::Info.
Namespace for OpenFOAM.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
runTime write()
labelList f(nPoints)
dictionary dict
cellMask correctBoundaryConditions()
CEqn solve()
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73