33 template<
class ParticleType>
42 template<
class ParticleType>
45 typeIOobject<timeIOdictionary> dictObj
49 "uniform"/cloud::prefix/
name(),
51 IOobject::MUST_READ_IF_MODIFIED,
56 if (dictObj.headerOk())
58 const timeIOdictionary uniformPropsDict(dictObj);
60 const word procName(
"processor" +
Foam::name(Pstream::myProcNo()));
61 if (uniformPropsDict.found(procName))
63 uniformPropsDict.subDict(procName).lookup(
"particleCount")
64 >> ParticleType::particleCount;
69 ParticleType::particleCount = 0;
74 template<
class ParticleType>
77 timeIOdictionary uniformPropsDict
83 "uniform"/cloud::prefix/
name(),
92 np[Pstream::myProcNo()] = ParticleType::particleCount;
94 Pstream::listCombineGather(np, maxEqOp<label>());
95 Pstream::listCombineScatter(np);
100 uniformPropsDict.add(procName, dictionary());
101 uniformPropsDict.subDict(procName).add(
"particleCount", np[i]);
104 uniformPropsDict.writeObject
107 IOstream::currentVersion,
108 time().writeCompression(),
114 template<
class ParticleType>
117 readCloudUniformProperties();
119 IOPosition<Cloud<ParticleType>> ioP(*
this);
121 bool valid = ioP.headerOk();
122 Istream& is = ioP.readStream(checkClass ? typeName :
"",
valid);
125 ioP.readData(is, *
this);
131 Pout<<
"Cannot read particle positions file:" <<
nl
132 <<
" " << ioP.objectPath() <<
nl
133 <<
"Assuming the initial cloud contains 0 particles." <<
endl;
139 pMesh_.tetBasePtIs();
145 template<
class ParticleType>
150 const bool checkClass
155 patchNbrProc_(patchNbrProc(pMesh)),
156 patchNbrProcPatch_(patchNbrProcPatch(pMesh)),
157 patchNonConformalCyclicPatches_(patchNonConformalCyclicPatches(pMesh)),
158 globalPositionsPtr_()
161 pMesh_.tetBasePtIs();
162 if (!ParticleType::instantaneous) pMesh_.oldCellCentres();
164 initCloud(checkClass);
170 template<
class ParticleType>
173 const word& fieldName,
189 template<
class ParticleType>
190 template<
class DataType>
197 if (data.
size() !=
c.size())
200 <<
"Size of " << data.
name()
201 <<
" field " << data.
size()
202 <<
" does not match the number of particles " <<
c.size()
208 template<
class ParticleType>
209 template<
class DataType>
216 if (data.size() !=
c.size())
219 <<
"Size of " << data.name()
220 <<
" field " << data.size()
221 <<
" does not match the number of particles " <<
c.size()
227 template<
class ParticleType>
230 ParticleType::writeFields(*
this);
234 template<
class ParticleType>
243 writeCloudUniformProperties();
246 return cloud::writeObject(fmt, ver, cmp, this->size());
252 template<
class ParticleType>
262 os.
check(
"Ostream& operator<<(Ostream&, const Cloud<ParticleType>&)");
#define forAll(list, i)
Loop across all elements in list.
A Field of objects of type <Type> with automated input and output using a compact storage....
Pre-declare SubField and related Field type.
A primitive field of type <Type> with automated input and output.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
readOption
Enumeration defining the read options.
const word & name() const
Return name.
virtual bool check(const char *operation) const
Check IOstream status for given operation.
streamFormat
Enumeration for the format of data in the stream.
compressionType
Enumeration for the format of data in the stream.
void size(const label)
Override size to be inconsistent with allocated storage.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
static word cloudPropertiesName
Name of cloud properties dictionary.
void checkFieldIOobject(const Cloud< ParticleType > &c, const IOField< DataType > &data) const
Check lagrangian data field.
virtual void writeFields() const
Write the field data for the cloud of particles Dummy at.
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp, const bool write=true) const
Write using given format, version and compression.
IOobject fieldIOobject(const word &fieldName, const IOobject::readOption r) const
Helper to construct IOobject for field and current time.
Cloud(const polyMesh &mesh, const word &cloudName, const IDLList< ParticleType > &particles)
Construct from mesh and a list of particles.
void checkFieldFieldIOobject(const Cloud< ParticleType > &c, const CompactIOField< Field< DataType >> &data) const
Check lagrangian data fieldfield.
Mesh consisting of general polyhedral cells.
A class for handling words, derived from string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
bool valid(const PtrList< ModelType > &l)
const dimensionedScalar c
Speed of light in a vacuum.
List< label > labelList
A List of labels.
Ostream & endl(Ostream &os)
Add newline and flush stream.
errorManip< error > abort(error &err)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
prefixOSstream Pout(cout, "Pout")
const word cloudName(propsDict.lookup("cloudName"))