39 template<
class Type,
template<
class>
class GeoField>
40 void Foam::fvMeshStitcher::resizePatchFields()
42 UPtrList<GeoField<Type>>
fields(mesh_.
fields<GeoField<Type>>());
48 typename GeoField<Type>::Patch& pf =
51 if (isA<nonConformalFvPatch>(pf.patch()))
53 pf.map(pf, setSizeFieldMapper(pf.patch().size()));
60 template<
template<
class>
class GeoField>
61 void Foam::fvMeshStitcher::resizePatchFields()
63 #define ResizePatchFields(Type, nullArg) \
64 resizePatchFields<Type, GeoField>();
66 #undef ResizePatchFields
71 void Foam::fvMeshStitcher::preConformVolFields()
73 UPtrList<VolField<Type>>
fields(mesh_.curFields<VolField<Type>>());
77 VolField<Type>& field =
fields[i];
79 for (
label ti=0; ti<=field.nOldTimes(
false); ti++)
83 boundaryFieldRefNoUpdate(field.oldTime(ti))
91 void Foam::fvMeshStitcher::preConformSurfaceFields()
93 UPtrList<SurfaceField<Type>>
fields(mesh_.curFields<SurfaceField<Type>>());
97 SurfaceField<Type>& field =
fields[i];
99 for (
label ti=0; ti<=field.nOldTimes(
false); ti++)
103 boundaryFieldRefNoUpdate(field.oldTime(ti))
111 void Foam::fvMeshStitcher::postUnconformVolFields()
113 UPtrList<VolField<Type>>
fields(mesh_.curFields<VolField<Type>>());
117 VolField<Type>& field =
fields[i];
119 for (
label ti=0; ti<=field.nOldTimes(
false); ti++)
123 boundaryFieldRefNoUpdate(field.oldTime(ti))
131 void Foam::fvMeshStitcher::postUnconformEvaluateVolFields()
137 isA<nonConformalFvPatch>(pf.patch())
138 && pf.type() == pf.patch().patch().type()
141 || isA<nonConformalErrorFvPatch>(pf.patch());
144 UPtrList<VolField<Type>>
fields(mesh_.fields<VolField<Type>>());
185 void Foam::fvMeshStitcher::postUnconformSurfaceFields()
187 UPtrList<SurfaceField<Type>>
fields(mesh_.curFields<SurfaceField<Type>>());
191 SurfaceField<Type>& field =
fields[i];
193 for (
label ti=0; ti<=field.nOldTimes(
false); ti++)
197 boundaryFieldRefNoUpdate(field.oldTime(ti))
206 template<
class GeoField>
212 return const_cast<typename GeoField::Boundary&
>(
fld.boundaryField());
#define forAll(list, i)
Loop across all elements in list.
PatchField< Type > Patch
Type of the patch field of which the Boundary is composed.
static label nRequests()
Get number of outstanding requests.
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
static bool & parRun()
Is this a parallel run?
static commsTypes defaultCommsType
Default commsType.
static GeoField::Boundary & boundaryFieldRefNoUpdate(GeoField &fld)
Access the boundary field reference of a field, without updating.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
UPtrList< GeoField > fields(const bool strict=false) const
Return the list of fields of type GeoField.
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
#define ResizePatchFields(Type, nullArg)
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)