32 template<
class GeoField>
33 bool Foam::LagrangianFieldDecomposer::decomposes(
const IOobjectList&
objects)
35 return !
objects.lookupClass(GeoField::typeName).empty();
39 template<
class Type,
template<
class>
class PrimitiveField>
44 Foam::LagrangianFieldDecomposer::decomposeLagrangianField
46 const DimensionedField<Type, LagrangianMesh, PrimitiveField>& completeField
49 PtrList<DimensionedField<Type, LagrangianMesh, PrimitiveField>>
50 procFields(procMeshes_.size());
57 new DimensionedField<Type, LagrangianMesh, PrimitiveField>
62 procMeshes_[proci].time().name(),
68 completeField.dimensions(),
72 particleProcAddressing_[proci]
82 template<
class Type,
template<
class>
class PrimitiveField>
87 Foam::LagrangianFieldDecomposer::decomposeLagrangianField
89 const GeometricField<Type, LagrangianMesh, PrimitiveField>& completeField
92 PtrList<GeometricField<Type, LagrangianMesh, PrimitiveField>>
93 procFields(procMeshes_.size());
102 PtrList<LagrangianPatchField<Type>> procPatchFields
104 procMeshes_[proci].
boundary().size()
106 forAll(completeField.boundaryField(), completePatchi)
111 completeField.boundaryField()[completePatchi].clone
119 label procPatchi = completeMesh_.boundary().size();
120 procPatchi < procMeshes_[proci].boundary().size();
130 procMeshes_[proci].
boundary()[procPatchi],
140 new GeometricField<Type, LagrangianMesh, PrimitiveField>
144 completeField.name(),
145 procMeshes_[proci].time().name(),
151 completeField.dimensions(),
155 particleProcAddressing_[proci]
158 completeField.sources().table()
169 template<
class GeoField>
177 decomposeLagrangianField
184 completeMesh_.time().name(),
196 template<
class GeoField>
206 Info<<
nl <<
" Decomposing " << GeoField::typeName <<
"s"
211 Info<<
" " << fieldIter()->name() <<
endl;
214 decomposeField<GeoField>(*fieldIter());
218 procFields[proci].write();
#define forAll(list, i)
Loop across all elements in list.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
List of IOobjects with searching and retrieving facilities.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
const word & name() const
Return name.
void decomposeFields(const IOobjectList &objects) const
Read, decompose and write all Lagrangian fields.
PtrList< GeoField > decomposeField(const IOobject &) const
Read and decompose a field.
static autoPtr< LagrangianPatchField< Type > > New(const word &patchFieldType, const word &actualPatchType, const LagrangianPatch &p, const regIOobject &)
Return a pointer to a new LagrangianPatchField with a given type.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
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.
Ostream & endl(Ostream &os)
Add newline and flush stream.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
faceListList boundary(nPatches)