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 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::oversetFvMeshBase
28 
29 Description
30  Support for overset functionality.
31 
32 SourceFiles
33  oversetFvMeshBase.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef oversetFvMeshBase_H
38 #define oversetFvMeshBase_H
39 
43 #include "volFieldsFwd.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 class lduPrimitiveProcessorInterface;
51 class GAMGAgglomeration;
52 
53 /*---------------------------------------------------------------------------*\
54  Class oversetFvMeshBase Declaration
55 \*---------------------------------------------------------------------------*/
56 
58 {
59  // Private Member Functions
60 
61  //- No copy construct
62  oversetFvMeshBase(const oversetFvMeshBase&) = delete;
63 
64  //- No copy assignment
65  void operator=(const oversetFvMeshBase&) = delete;
66 
67 
68 protected:
69 
70  // Protected Data
71 
72  //- Reference to mesh
73  const fvMesh& mesh_;
74 
75  //- Select base addressing (false) or locally stored extended
76  // lduAddressing (true)
77  mutable bool active_;
78 
79  //- Extended addressing (extended with local interpolation stencils)
81 
82  //- Added (processor)lduInterfaces for remote bits of stencil.
85 
86  //- Interfaces for above mesh. Contains both original and
87  //- above added processorLduInterfaces
89 
90  //- Corresponding faces (in above lduPtr) to the stencil
92 
93  //- Corresponding patches (in above lduPtr) to the stencil
95 
96  //- From old to new face labels
97  mutable labelList reverseFaceMap_;
98 
99 
100  // Protected Member Functions
101 
102  //- Calculate the extended lduAddressing
103  virtual bool updateAddressing() const;
104 
105  //- Debug: print matrix
106  template<class Type>
107  void write
108  (
109  Ostream&,
110  const fvMatrix<Type>&,
111  const lduAddressing&,
112  const lduInterfacePtrsList&
113  ) const;
114 
115  //- Freeze values at holes
116  //template<class Type>
117  //void freezeHoles(fvMatrix<Type>&) const;
118 
119  //- Scale coefficient depending on cell type
120  template<class Type>
121  void scaleConnection
122  (
123  Field<Type>& coeffs,
124  const labelUList& types,
125  const scalarList& factor,
126  const bool setHoleCellValue,
127  const label celli,
128  const label facei
129  ) const;
130 
131  //- Solve given dictionary with settings
132  template<class Type>
134  (
135  fvMatrix<Type>& m,
136  const dictionary&
137  ) const;
138 
139  //- Debug: correct coupled bc
140  template<class GeoField>
141  static void correctCoupledBoundaryConditions(GeoField& fld);
142 
143  //- Average norm of valid neighbours
144  scalar cellAverage
145  (
146  const labelList& types,
147  const labelList& nbrTypes,
148  const scalarField& norm,
149  const scalarField& nbrNorm,
150  const label celli,
151  bitSet& isFront
152  ) const;
153 
154  //- Debug: dump agglomeration
155  void writeAgglomeration
156  (
157  const GAMGAgglomeration& agglom
158  ) const;
159 
160 
161 public:
162 
163  //- Runtime type information
164  TypeName("oversetFvMeshBase");
165 
166 
167  // Constructors
168 
169  //- Construct from IOobject
170  oversetFvMeshBase(const fvMesh& mesh, const bool doInit=true);
171 
172 
173  //- Destructor
174  virtual ~oversetFvMeshBase();
175 
176 
177  // Member Functions
178 
179  // Extended addressing
180 
181  //- Return extended ldu addressing
183 
184  //- Return ldu addressing. If active: is (extended)
185  // primitiveLduAddr
186  virtual const lduAddressing& lduAddr() const;
187 
188  //- Return a list of pointers for each patch
189  // with only those pointing to interfaces being set. If active:
190  // return additional remoteStencilInterfaces_
191  virtual lduInterfacePtrsList interfaces() const;
192 
193  //- Return old to new face addressing
194  const labelList& reverseFaceMap() const
195  {
196  return reverseFaceMap_;
197  }
198 
199  //- Return true if using extended addressing
200  bool active() const
201  {
202  return active_;
203  }
204 
205  //- Enable/disable extended addressing
206  void active(const bool f) const
207  {
208  active_ = f;
209 
210  if (active_)
211  {
212  DebugInfo<< "Switching to extended addressing with nFaces:"
214  << " nInterfaces:" << allInterfaces_.size()
215  << endl;
216  }
217  else
218  {
219  DebugInfo<< "Switching to base addressing with nFaces:"
220  << mesh_.fvMesh::lduAddr().lowerAddr().size()
221  << " nInterfaces:" << mesh_.fvMesh::interfaces().size()
222  << endl;
223  }
224  }
225 
226 
227  // Overset
228 
229  //- Manipulate the matrix to add the interpolation/set hole
230  // values
231  template<class Type>
232  void addInterpolation
233  (
234  fvMatrix<Type>& m,
235  const scalarField& normalisation,
236  const bool setHoleCellValue,
237  const Type& holeCellValue
238  ) const;
239 
240 
241  //- Update the mesh for both mesh motion and topology change
242  virtual bool update();
244  //- Update fields when mesh is updated
245  virtual bool interpolateFields();
246 
247  //- Write using stream options
248  virtual bool writeObject
249  (
250  IOstreamOption streamOpt,
251  const bool writeOnProc
252  ) const;
253 
254  //- Debug: check halo swap is ok
255  template<class GeoField>
256  static void checkCoupledBC(const GeoField& fld);
257 
258  //- Correct boundary conditions of certain type (typeOnly = true)
259  // or explicitly not of the type (typeOnly = false)
260  template<class GeoField, class PatchType>
261  static void correctBoundaryConditions
262  (
263  typename GeoField::Boundary& bfld,
264  const bool typeOnly
265  );
266 
267  //- Determine normalisation for interpolation. This equals the
268  //- original diagonal except stabilised for zero diagonals (possible
269  //- in hole cells)
270  template<class Type>
272 };
273 
274 
275 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
276 
277 } // End namespace Foam
278 
279 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
280 
281 #ifdef NoRepository
283 #endif
284 
285 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
286 
287 #endif
288 
289 // ************************************************************************* //
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
TypeName("oversetFvMeshBase")
Runtime type information.
Forwards and collection of common volume field types.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
static void checkCoupledBC(const GeoField &fld)
Debug: check halo swap is ok.
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc) const
Write using stream options.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const fvMeshPrimitiveLduAddressing & primitiveLduAddr() const
Return extended ldu addressing.
lduInterfacePtrsList allInterfaces_
Interfaces for above mesh. Contains both original and above added processorLduInterfaces.
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.
A simple container for options an IOstream can normally have.
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.
SolverPerformance< Type > solveOverset(fvMatrix< Type > &m, const dictionary &) const
Solve given dictionary with settings.
labelList reverseFaceMap_
From old to new face labels.
labelListList stencilFaces_
Corresponding faces (in above lduPtr) to the stencil.
virtual const labelUList & lowerAddr() const noexcept
Return lower addressing (i.e. lower label = upper triangle)
label size() const noexcept
The number of elements in table.
Definition: HashTable.H:342
dynamicFvMesh & mesh
Generic templated field type.
Definition: Field.H:62
bool active() const
Return true if using extended addressing.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
void addInterpolation(fvMatrix< Type > &m, const scalarField &normalisation, const bool setHoleCellValue, const Type &holeCellValue) const
Manipulate the matrix to add the interpolation/set hole.
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
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.
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
#define DebugInfo
Report an information message using Foam::Info.
static void correctBoundaryConditions(typename GeoField::Boundary &bfld, const bool typeOnly)
Correct boundary conditions of certain type (typeOnly = true)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
labelList f(nPoints)
virtual bool interpolateFields()
Update fields when mesh is updated.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Support for overset functionality.
virtual bool update()
Update the mesh for both mesh motion and topology change.
autoPtr< fvMeshPrimitiveLduAddressing > lduPtr_
Extended addressing (extended with local interpolation stencils)
labelListList stencilPatches_
Corresponding patches (in above lduPtr) to the stencil.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
virtual const lduAddressing & lduAddr() const
Return ldu addressing. If active: is (extended)
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
PtrList< const lduPrimitiveProcessorInterface > remoteStencilInterfaces_
Added (processor)lduInterfaces for remote bits of stencil.
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
void writeAgglomeration(const GAMGAgglomeration &agglom) const
Debug: dump agglomeration.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
The class contains the addressing required by the lduMatrix: upper, lower and losort.
const labelList & reverseFaceMap() const
Return old to new face addressing.
virtual ~oversetFvMeshBase()
Destructor.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
void write(Ostream &, const fvMatrix< Type > &, const lduAddressing &, const lduInterfacePtrsList &) const
Debug: print matrix.
virtual bool updateAddressing() const
Calculate the extended lduAddressing.
Geometric agglomerated algebraic multigrid agglomeration class.
bool active_
Select base addressing (false) or locally stored extended.
const fvMesh & mesh_
Reference to mesh.
Variant of fvMeshLduAddressing that contains addressing instead of slices.
Namespace for OpenFOAM.