Go to the documentation of this file.
26 #include "parLagrangianRedistributor.H"
36 template<
class Container>
39 const IOobjectList& objects,
40 const HashSet<word>& selectedFields
43 const word fieldClassName(Container::typeName);
50 if (!selectedFields.empty())
69 const mapDistributeBase& map,
71 const IOobjectList& objects,
72 const HashSet<word>& selectedFields
77 filterObjects<IOField<Type> >
84 if (objectNames.size())
88 Info<<
" Redistributing lagrangian "
89 << fieldClassName <<
"s\n" <<
endl;
93 Info<<
" " << objectNames[i] <<
endl;
101 srcMesh_.time().timeName(),
111 map.distribute(field);
121 tgtMesh_.time().timeName(),
141 const mapDistributeBase& map,
143 const IOobjectList& objects,
144 const HashSet<word>& selectedFields
149 filterObjects<CompactIOField<Field<Type>, Type> >
160 filterObjects<IOField<Field<Type> > >
166 objectNames.append(ioFieldNames);
170 if (objectNames.size())
172 const word fieldClassName(CompactIOField<Field<Type>, Type>::typeName);
174 Info<<
" Redistributing lagrangian "
175 << fieldClassName <<
"s\n" <<
endl;
179 Info<<
" " << objectNames[i] <<
endl;
182 CompactIOField<Field<Type>, Type > field
187 srcMesh_.time().timeName(),
198 map.distribute(field);
203 CompactIOField<Field<Type>, Type>
208 tgtMesh_.time().timeName(),
215 xferMove<Field<Field<Type> > >(field)
223 template<
class Container>
226 const passiveParticleCloud& cloud,
227 const IOobjectList& objects,
228 const HashSet<word>& selectedFields
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
257 cloud.time().timeName(),
271 template<
class Container>
274 const mapDistributeBase& map,
275 passiveParticleCloud& cloud
278 HashTable<Container*>
fields
280 cloud.lookupClass<Container >()
285 const word fieldClassName(Container::typeName);
287 Info<<
" Redistributing lagrangian "
288 << fieldClassName <<
"s\n" <<
endl;
292 Container& field = *iter();
296 map.distribute(field);
305 tgtMesh_.time().timeName(),
312 xferMove<Field<typename Container::value_type> >(field)
void redistributeStoredLagrangianFields(const mapDistributeBase &map, passiveParticleCloud &cloud) const
Redistribute and write stored lagrangian fields.
static void readLagrangianFields(const passiveParticleCloud &cloud, const IOobjectList &objects, const HashSet< word > &selectedFields)
Read and store all fields of a cloud.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
#define forAll(list, i)
Loop across all elements in list.
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
static const char *const typeName
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Ostream & endl(Ostream &os)
Add newline and flush stream.
void redistributeLagrangianFields(const mapDistributeBase &map, const word &cloudName, const IOobjectList &objects, const HashSet< word > &selectedFields) const
Read, redistribute and write all/selected lagrangian fields.
List< word > wordList
A List of words.
virtual Ostream & write(const token &)=0
Write next token to stream.
static wordList filterObjects(const IOobjectList &objects, const HashSet< word > &selectedFields)
Pick up any fields of a given type.
static void combineGather(const List< commsStruct > &comms, T &Value, const CombineOp &cop, const int tag, const label comm)
static const word prefix
The prefix to local: lagrangian.
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.
wordList fieldNames(const IOobjectList &objects, const bool syncPar)
Get sorted names of fields of type. If syncPar and running in parallel.
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"))
void size(const label)
Override size to be inconsistent with allocated storage.
static void combineScatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.