Go to the documentation of this file.
36 template<
class Container>
43 const word fieldClassName(Container::typeName);
50 if (!selectedFields.
empty())
84 if (objectNames.
size())
88 Info<<
" Redistributing lagrangian "
89 << fieldClassName <<
"s\n" <<
endl;
93 Info<<
" " << objectNames[i] <<
endl;
101 srcMesh_.time().timeName(),
104 IOobject::READ_IF_PRESENT,
121 tgtMesh_.time().timeName(),
166 objectNames.
append(ioFieldNames);
170 if (objectNames.
size())
174 Info<<
" Redistributing lagrangian "
175 << fieldClassName <<
"s\n" <<
endl;
179 Info<<
" " << objectNames[i] <<
endl;
187 srcMesh_.time().timeName(),
190 IOobject::READ_IF_PRESENT,
208 tgtMesh_.time().timeName(),
215 xferMove<Field<Field<Type> > >(field)
223 template<
class Container>
233 filterObjects<Container>
240 if (objectNames.
size())
242 const word fieldClassName(Container::typeName);
244 Info<<
" Reading lagrangian "
245 << fieldClassName <<
"s\n" <<
endl;
249 Info<<
" " << objectNames[i] <<
endl;
252 Container* fieldPtr =
new Container
259 IOobject::READ_IF_PRESENT,
271 template<
class Container>
285 const word fieldClassName(Container::typeName);
287 Info<<
" Redistributing lagrangian "
288 << fieldClassName <<
"s\n" <<
endl;
292 Container& field = *iter();
305 tgtMesh_.time().timeName(),
void redistributeStoredLagrangianFields(const mapDistributeBase &map, passiveParticleCloud &cloud) const
Redistribute and write stored lagrangian fields.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
static void readLagrangianFields(const passiveParticleCloud &cloud, const IOobjectList &objects, const HashSet< word > &selectedFields)
Read and store all fields of a cloud.
A class for handling words, derived from string.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
A primitive field of type <T> with automated input and output.
#define forAll(list, i)
Loop across all elements in list.
static void distribute(const Pstream::commsTypes commsType, const List< labelPair > &schedule, const label constructSize, const labelListList &subMap, const bool subHasFlip, const labelListList &constructMap, const bool constructHasFlip, List< T > &, const negateOp &negOp, const int tag=UPstream::msgType())
Distribute data. Note:schedule only used for Pstream::scheduled.
Info<< "Creating field dpdt\n"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar("dpdt", p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);p_rgh=p - rho *gh;mesh.setFluxRequired(p_rgh.name());multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
const Time & time() const
Return time.
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
virtual bool write() const
Write using setting from DB.
Ostream & endl(Ostream &os)
Add newline and flush stream.
A Field of objects of type <T> with automated input and output using a compact storage....
wordList names() const
Return the list of names of the IOobjects.
A HashTable with keys but without contents.
void redistributeLagrangianFields(const mapDistributeBase &map, const word &cloudName, const IOobjectList &objects, const HashSet< word > &selectedFields) const
Read, redistribute and write all/selected lagrangian fields.
HashTable< const Type * > lookupClass(const bool strict=false) const
Lookup and return all objects of the given Type.
void append(const T &)
Append an element at the end of the list.
bool empty() const
Return true if the hash table is empty.
const word & name() const
Return name.
A Cloud of passive particles.
static wordList filterObjects(const IOobjectList &objects, const HashSet< word > &selectedFields)
Pick up any fields of a given type.
static List< word > fieldNames
bool found(const Key &) const
Return true if hashedEntry is found in table.
List of IOobjects with searching and retrieving facilities.
Helper class for list to append unique elelements of y onto the end of x.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
An STL-conforming hash table.
A cloud is a collection of lagrangian particles.
Xfer< T > xferMove(T &)
Construct by transferring the contents of the arg.
void redistributeLagrangianFieldFields(const mapDistributeBase &map, const word &cloudName, const IOobjectList &objects, const HashSet< word > &selectedFields) const
Read, redistribute and write all/selected lagrangian fieldFields.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
const word cloudName(propsDict.lookup("cloudName"))
Class containing processor-to-processor mapping information.
void size(const label)
Override size to be inconsistent with allocated storage.