Loading...
Searching...
No Matches
oversetFvMeshBase.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-2022,2024 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::oversetFvMeshBase
28
29Description
30 Support for overset functionality.
31
32SourceFiles
33 oversetFvMeshBase.C
34
35\*---------------------------------------------------------------------------*/
36
37#ifndef oversetFvMeshBase_H
38#define oversetFvMeshBase_H
39
42#include "volFieldsFwd.H"
43
44// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45
46namespace Foam
47{
48
52/*---------------------------------------------------------------------------*\
53 Class oversetFvMeshBase Declaration
54\*---------------------------------------------------------------------------*/
55
56class oversetFvMeshBase
57{
58 // Private Member Functions
59
60 //- No copy construct
61 oversetFvMeshBase(const oversetFvMeshBase&) = delete;
62
63 //- No copy assignment
64 void operator=(const oversetFvMeshBase&) = delete;
65
66
67protected:
68
69 // Protected Data
70
71 //- Reference to mesh
72 const fvMesh& mesh_;
74 //- Select base addressing (false) or locally stored extended
75 // lduAddressing (true)
76 mutable bool active_;
77
78 //- Extended addressing (extended with local interpolation stencils)
81 //- Added (processor)lduInterfaces for remote bits of stencil.
84
85 //- Interfaces for above mesh. Contains both original and
86 //- above added processorLduInterfaces
88
89 //- Corresponding faces (in above lduPtr) to the stencil
92 //- Corresponding patches (in above lduPtr) to the stencil
94
95 //- From old to new face labels
98
99 // Protected Member Functions
100
101 //- Calculate the extended lduAddressing
102 virtual bool updateAddressing() const;
103
104 //- Debug: print matrix
105 template<class Type>
106 void write
108 Ostream&,
109 const fvMatrix<Type>&,
110 const lduAddressing&,
112 ) const;
113
114 //- Freeze values at holes
115 //template<class Type>
116 //void freezeHoles(fvMatrix<Type>&) const;
117
118 //- Scale coefficient depending on cell type
119 template<class Type>
120 void scaleConnection
121 (
122 Field<Type>& coeffs,
123 const labelUList& types,
124 const scalarList& factor,
125 const bool setHoleCellValue,
126 const label celli,
127 const label facei
128 ) const;
129
130 //- Solve given dictionary with settings
131 template<class Type>
133 (
135 const dictionary&
136 ) const;
137
138 //- Debug: correct coupled bc
139 template<class GeoField>
140 static void correctCoupledBoundaryConditions(GeoField& fld);
141
142 //- Average norm of valid neighbours
143 scalar cellAverage
144 (
145 const labelUList& types,
146 const labelUList& nbrTypes,
147 const scalarField& norm,
148 const scalarField& nbrNorm,
149 const label celli,
150 bitSet& isFront
151 ) const;
152
153 //- Debug: dump agglomeration
155 (
156 const GAMGAgglomeration& agglom
157 ) const;
158
160public:
161
162 //- Runtime type information
163 TypeName("oversetFvMeshBase");
164
165
166 // Constructors
167
168 //- Construct from IOobject
169 oversetFvMeshBase(const fvMesh& mesh, const bool doInit=true);
170
171
172 //- Destructor
173 virtual ~oversetFvMeshBase();
174
175
176 // Member Functions
177
178 // Extended addressing
179
180 //- Return extended ldu addressing
182
183 //- Return ldu addressing. If active: is (extended)
184 // primitiveLduAddr
185 virtual const lduAddressing& lduAddr() const;
186
187 //- Return a list of pointers for each patch
188 // with only those pointing to interfaces being set. If active:
189 // return additional remoteStencilInterfaces_
190 virtual lduInterfacePtrsList interfaces() const;
191
192 //- Return old to new face addressing
193 const labelList& reverseFaceMap() const
194 {
195 return reverseFaceMap_;
196 }
197
198 //- Return true if using extended addressing
199 bool active() const
200 {
201 return active_;
202 }
203
204 //- Enable/disable extended addressing
205 void active(const bool f) const
206 {
207 active_ = f;
208
209 if (active_)
210 {
211 DebugInfo<< "Switching to extended addressing with nFaces:"
213 << " nInterfaces:" << allInterfaces_.size()
214 << endl;
215 }
216 else
217 {
218 DebugInfo<< "Switching to base addressing with nFaces:"
219 << mesh_.fvMesh::lduAddr().lowerAddr().size()
220 << " nInterfaces:" << mesh_.fvMesh::interfaces().size()
221 << endl;
222 }
223 }
224
225
226 // Overset
227
228 //- Manipulate the matrix to add the interpolation/set hole
229 // values
230 template<class Type>
232 (
233 fvMatrix<Type>& m,
235 const bool setHoleCellValue,
236 const Type& holeCellValue
237 ) const;
238
239
240 //- Clear out local storage
241 void clearOut();
243 //- Update the mesh for both mesh motion and topology change
244 virtual bool update();
245
246 //- Update fields when mesh is updated
247 virtual bool interpolateFields();
248
249 //- Write using stream options
250 virtual bool writeObject
251 (
252 IOstreamOption streamOpt,
253 const bool writeOnProc
254 ) const;
255
256 //- Debug: check halo swap is ok
257 template<class GeoField>
258 static void checkCoupledBC(const GeoField& fld);
259
260 //- Correct boundary conditions of certain type (TypeOnly = true)
261 //- or explicitly not of the type (TypeOnly = false)
262 template<class GeoField, class PatchType, bool TypeOnly>
263 static void correctBoundaryConditions
264 (
265 typename GeoField::Boundary& bfld
266 );
267
268 //- Determine normalisation for interpolation. This equals the
269 //- original diagonal except stabilised for zero diagonals (possible
270 //- in hole cells)
271 template<class Type>
273};
274
275
276// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
277
278} // End namespace Foam
279
280// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281
282#ifdef NoRepository
284#endif
285
286// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
287
288#endif
289
290// ************************************************************************* //
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Geometric agglomerated algebraic multigrid agglomeration class.
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
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:59
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
Variant of fvMeshLduAddressing that contains addressing instead of slices.
virtual const labelUList & lowerAddr() const noexcept
Return lower addressing (i.e. lower label = upper triangle).
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
The class contains the addressing required by the lduMatrix: upper, lower and losort.
Concrete implementation of processor interface. Used to temporarily store settings.
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc) const
Write using stream options.
bool active_
Select base addressing (false) or locally stored extended.
void writeAgglomeration(const GAMGAgglomeration &agglom) const
Debug: dump agglomeration.
void scaleConnection(Field< Type > &coeffs, const labelUList &types, const scalarList &factor, const bool setHoleCellValue, const label celli, const label facei) const
Freeze values at holes.
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
virtual bool interpolateFields()
Update fields when mesh is updated.
const fvMesh & mesh_
Reference to mesh.
SolverPerformance< Type > solveOverset(fvMatrix< Type > &m, const dictionary &) const
Solve given dictionary with settings.
virtual const lduAddressing & lduAddr() const
Return ldu addressing. If active: is (extended).
PtrList< const lduPrimitiveProcessorInterface > remoteStencilInterfaces_
Added (processor)lduInterfaces for remote bits of stencil.
autoPtr< fvMeshPrimitiveLduAddressing > lduPtr_
Extended addressing (extended with local interpolation stencils).
static void checkCoupledBC(const GeoField &fld)
Debug: check halo swap is ok.
labelListList stencilFaces_
Corresponding faces (in above lduPtr) to the stencil.
void active(const bool f) const
Enable/disable extended addressing.
virtual bool updateAddressing() const
Calculate the extended lduAddressing.
labelListList stencilPatches_
Corresponding patches (in above lduPtr) to the stencil.
const fvMeshPrimitiveLduAddressing & primitiveLduAddr() const
Return extended ldu addressing.
virtual bool update()
Update the mesh for both mesh motion and topology change.
void addInterpolation(fvMatrix< Type > &m, const scalarField &normalisation, const bool setHoleCellValue, const Type &holeCellValue) const
Manipulate the matrix to add the interpolation/set hole.
TypeName("oversetFvMeshBase")
Runtime type information.
lduInterfacePtrsList allInterfaces_
Interfaces for above mesh. Contains both original and above added processorLduInterfaces.
virtual ~oversetFvMeshBase()
Destructor.
tmp< scalarField > normalisation(const fvMatrix< Type > &m) const
Determine normalisation for interpolation. This equals the original diagonal except stabilised for ze...
static void correctCoupledBoundaryConditions(GeoField &fld)
Debug: correct coupled bc.
scalar cellAverage(const labelUList &types, const labelUList &nbrTypes, const scalarField &norm, const scalarField &nbrNorm, const label celli, bitSet &isFront) const
Average norm of valid neighbours.
bool active() const
Return true if using extended addressing.
const labelList & reverseFaceMap() const
Return old to new face addressing.
void clearOut()
Clear out local storage.
labelList reverseFaceMap_
From old to new face labels.
A class for managing temporary objects.
Definition tmp.H:75
dynamicFvMesh & mesh
#define DebugInfo
Report an information message using Foam::Info.
Namespace for OpenFOAM.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
UPtrList< const lduInterface > lduInterfacePtrsList
Store lists of lduInterface as a UPtrList.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
UList< label > labelUList
A UList of labels.
Definition UList.H:75
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
runTime write()
labelList f(nPoints)
cellMask correctBoundaryConditions()
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68
Forwards and collection of common volume field types.