37#ifndef dynamicOversetFvMesh_H
38#define dynamicOversetFvMesh_H
52class lduPrimitiveProcessorInterface;
53class GAMGAgglomeration;
114 template<
class GeoField>
122 template<
class GeoField>
148 template<
class GeoField>
226 void active(
const bool f)
const
232 DebugInfo<<
"Switching to extended addressing with nFaces:"
238 DebugInfo<<
"Switching to base addressing with nFaces:"
251 interpolate<scalar>(
psi);
256 interpolate<vector>(
psi);
261 interpolate<sphericalTensor>(
psi);
266 interpolate<symmTensor>(
psi);
271 interpolate<tensor>(
psi);
276 interpolate<volScalarField>(
psi);
281 interpolate<volVectorField>(
psi);
286 interpolate<volSphericalTensorField>(
psi);
291 interpolate<volSymmTensorField>(
psi);
296 interpolate<volTensorField>(
psi);
310 return solve<scalar>(m,
dict);
321 return solve<vector>(m,
dict);
332 return solve<symmTensor>(m,
dict);
343 return solve<tensor>(m,
dict);
348 virtual bool init(
const bool doInit);
364 template<
class GeoField>
369 template<
class GeoField,
class PatchType>
372 typename GeoField::Boundary& bfld,
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.
Geometric agglomerated algebraic multigrid agglomeration class.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
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,...
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
void size(const label n)
Older name for setAddressableSize.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
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....
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.
const word & name() const
Return reference to name.
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.
A class for handling words, derived from Foam::string.
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.
Ostream & endl(Ostream &os)
Add newline and flush stream.
cellMask correctBoundaryConditions()
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.