40 template<
class FieldType>
41 bool Foam::fvFieldReconstructor::reconstructs
44 const HashSet<word>& selectedFields
47 IOobjectList
fields =
objects.lookupClass(FieldType::typeName);
49 if (
fields.size() && selectedFields.empty())
56 if (selectedFields.found(fieldIter()->
name()))
67 void Foam::fvFieldReconstructor::rmapFaceToFace
70 const Field<Type>& fromField,
77 toField[
mag(addressing[i]) - 1] =
78 (
isFlux && addressing[i] < 0 ? -1 : +1)*fromField[i];
85 Foam::fvFieldReconstructor::reconstructVolInternalField
87 const IOobject& fieldIoObject
91 PtrList<DimensionedField<Type, volMesh>> procFields(procMeshes_.size());
97 new DimensionedField<Type, volMesh>
101 fieldIoObject.name(),
102 procMeshes_[proci].time().name(),
114 Field<Type> internalField(completeMesh_.nCells());
116 forAll(procMeshes_, proci)
118 const DimensionedField<Type, volMesh>& procField = procFields[proci];
123 procField.primitiveField(),
124 cellProcAddressing_[proci]
128 return tmp<DimensionedField<Type, volMesh>>
130 new DimensionedField<Type, volMesh>
134 fieldIoObject.name(),
135 completeMesh_.time().name(),
142 procFields[0].dimensions(),
151 Foam::fvFieldReconstructor::reconstructVolField
153 const IOobject& fieldIoObject
157 PtrList<VolField<Type>> procFields(procMeshes_.size());
158 forAll(procMeshes_, proci)
167 fieldIoObject.name(),
168 procMeshes_[proci].time().name(),
180 Field<Type> internalField(completeMesh_.nCells());
183 PtrList<fvPatchField<Type>> patchFields(completeMesh_.boundary().size());
187 const VolField<Type>& procField =
193 procField.primitiveField(),
194 cellProcAddressing_[proci]
198 forAll(procField.boundaryField(), procPatchi)
200 const fvPatch& procPatch =
201 procMeshes_[proci].boundary()[procPatchi];
203 const label completePatchi = completePatchID(proci, procPatchi);
205 if (completePatchi == procPatchi)
207 if (!patchFields(completePatchi))
214 procField.boundaryField()[procPatchi],
215 completeMesh_.boundary()[completePatchi],
219 completeMesh_.boundary()[completePatchi].size()
225 patchFields[completePatchi].map
227 procField.boundaryField()[procPatchi],
230 faceProcAddressingBf_[proci][procPatchi] - 1
234 else if (isA<processorCyclicFvPatch>(procPatch))
236 if (!patchFields(completePatchi))
243 completeMesh_.boundary()[completePatchi].type(),
244 completeMesh_.boundary()[completePatchi],
250 if (patchFields[completePatchi].overridesConstraint())
253 str <<
"\nThe field \"" << procFields[0].name()
254 <<
"\" on cyclic patch \""
255 << patchFields[completePatchi].patch().name()
256 <<
"\" cannot be reconstructed as it is not a cyclic "
257 <<
"patch field. A \"patchType cyclic;\" setting has "
258 <<
"been used to override the cyclic patch type.\n\n"
259 <<
"Cyclic patches like this with non-cyclic boundary "
260 <<
"conditions should be confined to a single "
261 <<
"processor using decomposition constraints.";
267 patchFields[completePatchi].map
269 procField.boundaryField()[procPatchi],
272 faceProcAddressingBf_[proci][procPatchi] - 1
280 return tmp<VolField<Type>>
286 fieldIoObject.name(),
287 completeMesh_.time().name(),
294 procFields[0].dimensions(),
297 procFields[0].sources().table()
305 Foam::fvFieldReconstructor::reconstructFvSurfaceField
307 const IOobject& fieldIoObject
311 PtrList<SurfaceField<Type>> procFields(procMeshes_.size());
312 forAll(procMeshes_, proci)
317 new SurfaceField<Type>
321 fieldIoObject.name(),
322 procMeshes_[proci].time().name(),
334 Field<Type> internalField(completeMesh_.nInternalFaces());
337 PtrList<fvsPatchField<Type>> patchFields(completeMesh_.boundary().size());
339 forAll(procMeshes_, proci)
341 const SurfaceField<Type>& procField =
348 procField.primitiveField(),
351 faceProcAddressing_[proci],
352 procMeshes_[proci].nInternalFaces()
358 forAll(procField.boundaryField(), procPatchi)
360 const fvPatch& procPatch =
361 procMeshes_[proci].boundary()[procPatchi];
363 const label completePatchi = completePatchID(proci, procPatchi);
365 if (completePatchi == procPatchi)
367 if (!patchFields(completePatchi))
374 procField.boundaryField()[procPatchi],
375 completeMesh_.boundary()[completePatchi],
379 completeMesh_.boundary()[completePatchi].size()
385 patchFields[completePatchi].map
387 procField.boundaryField()[procPatchi],
390 faceProcAddressingBf_[proci][procPatchi] - 1
394 else if (isA<processorCyclicFvPatch>(procPatch))
396 if (!patchFields(completePatchi))
403 completeMesh_.boundary()[completePatchi].type(),
404 completeMesh_.boundary()[completePatchi],
410 patchFields[completePatchi].map
412 procField.boundaryField()[procPatchi],
415 faceProcAddressingBf_[proci][procPatchi] - 1
419 else if (isA<processorFvPatch>(procPatch))
424 procField.boundaryField()[procPatchi],
425 faceProcAddressingBf_[proci][procPatchi],
433 return tmp<SurfaceField<Type>>
435 new SurfaceField<Type>
439 fieldIoObject.name(),
440 completeMesh_.time().name(),
447 procFields[0].dimensions(),
461 const HashSet<word>& selectedFields
470 Info<<
nl <<
" Reconstructing " << fieldClassName <<
"s"
477 selectedFields.empty()
478 || selectedFields.found(fieldIter()->
name())
481 Info<<
" " << fieldIter()->name() <<
endl;
483 reconstructVolInternalField<Type>(*fieldIter())().
write();
494 const HashSet<word>& selectedFields
497 const word& fieldClassName =
504 Info<<
nl <<
" Reconstructing " << fieldClassName <<
"s"
511 selectedFields.empty()
512 || selectedFields.found(fieldIter()->
name())
515 Info<<
" " << fieldIter()->name() <<
endl;
517 reconstructVolField<Type>(*fieldIter())().
write();
528 const HashSet<word>& selectedFields
531 const word& fieldClassName =
538 Info<<
nl <<
" Reconstructing " << fieldClassName <<
"s"
545 selectedFields.empty()
546 || selectedFields.found(fieldIter()->
name())
549 Info<<
" " << fieldIter()->name() <<
endl;
551 reconstructFvSurfaceField<Type>(*fieldIter())().
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.
static const DimensionedField< Type, GeoMesh > & null()
Return a null DimensionedField.
static const char *const typeName
void reconstructVolInternalFields(const IOobjectList &objects, const HashSet< word > &selectedFields)
Read, reconstruct and write all/selected volume internal fields.
void reconstructVolFields(const IOobjectList &objects, const HashSet< word > &selectedFields)
Read, reconstruct and write all/selected volume fields.
void reconstructFvSurfaceFields(const IOobjectList &objects, const HashSet< word > &selectedFields)
Read, reconstruct and write all/selected surface fields.
static tmp< fvPatchField< Type > > New(const word &, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
static tmp< fvsPatchField< Type > > New(const word &, const fvPatch &, const DimensionedField< Type, surfaceMesh > &)
Return a pointer to a new patchField created on freestore given.
A class for managing temporary objects.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
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
string breakIntoIndentedLines(const string &str, const string::size_type nLength=80, const string::size_type nIndent=0)
Break a string up into indented lines.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
errorManipArg< error, int > exit(error &err, const int errNo=1)
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.
word name(const bool)
Return a word representation of a bool.
bool isFlux(const DimensionedField< Type, surfaceMesh > &df)
Check if surfaceField is a flux.
dimensioned< scalar > mag(const dimensioned< Type > &)
UList< label > labelUList