Go to the documentation of this file.
54 regionNames[i] = meshes[i].dbDir();
60 if (Pstream::master())
64 groupDir(commsDir_, compositeName(regionNames), groupName)
70 Info<<
type() <<
": reading data from " << transferFile <<
endl;
74 if (!masterFilePtr().good())
77 <<
"Cannot open file for region " << compositeName(regionNames)
78 <<
", field " << fieldName
91 if (!
mesh.foundObject<volFieldType>(fieldName))
98 const volFieldType& cvf =
mesh.lookupObject<volFieldType>(fieldName);
99 const typename volFieldType::GeometricBoundaryField& bf =
106 mesh.boundaryMesh().patchSet
115 label patchI = patchIDs[i];
117 if (isA<patchFieldType>(bf[patchI]))
122 patchFieldType& pf =
const_cast<patchFieldType&
>
124 refCast<const patchFieldType>
144 pf.patchFieldType::evaluate();
162 refCast<const mixedFvPatchField<Type> >
176 cmpt < pTraits<Type>::nComponents;
186 cmpt < pTraits<Type>::nComponents;
196 pf.mixedFvPatchField<Type>::evaluate();
213 refCast<const fixedGradientFvPatchField<Type> >
224 cmpt < pTraits<Type>::nComponents;
237 pf.fixedGradientFvPatchField<Type>::evaluate();
256 cmpt < pTraits<Type>::nComponents;
266 refCast<const fixedValueFvPatchField<Type> >
276 pf.fixedValueFvPatchField<Type>::evaluate();
281 <<
"Unsupported boundary condition " << bf[patchI].type()
282 <<
" for patch " << bf[patchI].patch().name()
283 <<
" in region " <<
mesh.name()
304 gatheredValues[Pstream::myProcNo()] =
fld;
305 Pstream::gatherList(gatheredValues);
311 if (Pstream::master())
314 label globalElemI = 0;
316 forAll(gatheredValues, lstI)
318 globalElemI += gatheredValues[lstI].
size();
321 result.setSize(globalElemI);
325 forAll(gatheredValues, lstI)
331 result[globalElemI++] = sub[elemI];
345 const word& fieldName
354 regionNames[i] = meshes[i].dbDir();
360 if (Pstream::master())
364 groupDir(commsDir_, compositeName(regionNames), groupName)
370 Info<<
type() <<
": writing data to " << transferFile <<
endl;
374 if (!masterFilePtr().good())
377 <<
"Cannot open file for region " << compositeName(regionNames)
378 <<
", field " << fieldName
384 bool headerDone =
false;
392 if (!
mesh.foundObject<volFieldType>(fieldName))
399 const volFieldType& cvf =
mesh.lookupObject<volFieldType>(fieldName);
400 const typename volFieldType::GeometricBoundaryField& bf =
407 mesh.boundaryMesh().patchSet
416 label patchI = patchIDs[i];
420 if (isA<patchFieldType>(bf[patchI]))
425 const patchFieldType& pf = refCast<const patchFieldType>
436 if (Pstream::master())
441 pf.writeHeader(masterFilePtr());
444 masterFilePtr() << os.
str().c_str();
446 for (
label procI = 1; procI < Pstream::nProcs(); procI++)
448 IPstream fromSlave(Pstream::scheduled, procI);
449 string str(fromSlave);
450 masterFilePtr() << str.c_str();
455 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
456 toMaster << os.
str();
462 refCast<const mixedFvPatchField<Type> >(bf[patchI]);
470 if (Pstream::master())
475 << value[faceI] << token::SPACE
476 <<
snGrad[faceI] << token::SPACE
477 << refValue[faceI] << token::SPACE
478 << refGrad[faceI] << token::SPACE
479 << valueFraction[faceI] <<
nl;
488 if (Pstream::master())
493 << value[faceI] << token::SPACE
void writeData() const
Write data for all regions, all fields.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
virtual scalarField & valueFraction()
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
A class for handling words, derived from string.
A class for handling file names.
virtual Field< Type > & refValue()
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
Output inter-processor communications stream.
static tmp< Field< Type > > gatherAndCombine(const Field< Type > &fld)
Helper: append data from all processors onto master.
virtual Field< Type > & refGrad()
Ostream & endl(Ostream &os)
Add newline and flush stream.
virtual Field< Type > & gradient()
Return gradient at boundary.
A wordRe is a word, but can also have a regular expression for matching words.
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
string str() const
Return the string.
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
This boundary condition extends the mixed boundary condition with serialisation through.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Mesh data needed to do the Finite Volume discretisation.
Input from memory buffer stream.
This boundary condition provides a base class for 'mixed' type boundary conditions,...
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
errorManipArg< error, int > exit(error &err, const int errNo=1)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Traits class for primitives.
Output to memory buffer stream.
void reset(T *=0)
If object pointer already set, delete object and set to given.
Input inter-processor communications stream.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
void size(const label)
Override size to be inconsistent with allocated storage.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
This boundary condition supplies a fixed gradient condition, such that the patch values are calculate...
Generic GeometricField class.
Database for solution data, solver performance and other reduced data.
void readData()
Read data for all regions, all fields.
label size() const
Return the number of elements in the UPtrList.