Loading...
Searching...
No Matches
CloudIO.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017, 2020 OpenFOAM Foundation
9 Copyright (C) 2017-2024 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "Cloud.H"
30#include "Time.H"
31#include "IOPosition.H"
32#include "IOdictionary.H"
33#include "IOobjectList.H"
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37template<class ParticleType>
39
40
41// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
42
43template<class ParticleType>
44void Foam::Cloud<ParticleType>::readCloudUniformProperties()
45{
46 IOobject dictObj
47 (
48 cloudPropertiesName,
49 time().timeName(),
50 "uniform"/cloud::prefix/name(),
51 db(),
52 IOobject::MUST_READ,
53 IOobject::NO_WRITE,
54 IOobject::NO_REGISTER
55 );
56
57 if (dictObj.typeHeaderOk<IOdictionary>(true))
58 {
59 const IOdictionary uniformPropsDict(dictObj);
60
61 // Fall back to positions mode if the entry is not present for
62 // backwards compatibility
63 geometryType_ =
64 cloud::geometryTypeNames.getOrDefault
65 (
66 "geometry",
67 uniformPropsDict,
68 cloud::geometryType::POSITIONS
69 );
70
71 const word procName("processor" + Foam::name(Pstream::myProcNo()));
72
73 const dictionary* dictptr = uniformPropsDict.findDict(procName);
74
75 if (dictptr)
76 {
77 dictptr->readEntry("particleCount", ParticleType::particleCount_);
78 }
79 }
80 else
81 {
82 ParticleType::particleCount_ = 0;
83 }
84}
85
86
87template<class ParticleType>
88void Foam::Cloud<ParticleType>::writeCloudUniformProperties() const
89{
90 IOdictionary uniformPropsDict
91 (
92 IOobject
93 (
94 cloudPropertiesName,
95 time().timeName(),
96 "uniform"/cloud::prefix/name(),
97 db(),
98 IOobject::NO_READ,
99 IOobject::NO_WRITE,
100 IOobject::NO_REGISTER
101 )
102 );
103
104 labelList np(UPstream::nProcs(), Foam::zero{});
105 np[UPstream::myProcNo()] = ParticleType::particleCount_;
106 Pstream::allGatherList(np);
107
108 uniformPropsDict.add
109 (
110 "geometry",
111 cloud::geometryTypeNames[geometryType_]
112 );
113
114 forAll(np, i)
115 {
116 const word procName("processor" + Foam::name(i));
117 uniformPropsDict.subDictOrAdd(procName).add("particleCount", np[i]);
118 }
119
120 uniformPropsDict.writeObject
121 (
122 IOstreamOption(IOstreamOption::ASCII, time().writeCompression()),
123 true
124 );
125}
126
127
128template<class ParticleType>
129void Foam::Cloud<ParticleType>::initCloud(const bool checkClass)
130{
131 readCloudUniformProperties();
132
133 IOPosition<Cloud<ParticleType>> ioP(*this, geometryType_);
134
135 const bool haveFile = ioP.headerOk();
136 Istream& is = ioP.readStream(checkClass ? typeName : word::null, haveFile);
137 if (haveFile)
138 {
139 ioP.readData(is, *this);
140 ioP.close();
141 }
142
143 if (!haveFile && debug)
144 {
145 Pout<< "Not reading particle positions file: "
146 << ioP.objectRelPath() << nl
147 << "Assuming the initial cloud contains 0 particles." << endl;
148 }
149
150 // Always operate in coordinates mode after reading
151 geometryType_ = cloud::geometryType::COORDINATES;
152
153 // Ask for the tetBasePtIs to trigger all processors to build
154 // them, otherwise, if some processors have no particles then
155 // there is a comms mismatch.
156 (void)polyMesh_.tetBasePtIs();
157}
158
159
160// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
161
162template<class ParticleType>
164(
165 const polyMesh& pMesh,
166 const word& cloudName,
167 const bool checkClass
168)
169:
170 Cloud<ParticleType>(pMesh, Foam::zero{}, cloudName)
171{
172 initCloud(checkClass);
174
175
176// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
177
178template<class ParticleType>
179template<class DataType>
181(
182 const Cloud<ParticleType>& c,
183 const IOField<DataType>& data
184) const
185{
186 if (data.size() != c.size())
187 {
189 << "Size of " << data.name()
190 << " field " << data.size()
191 << " does not match the number of particles " << c.size()
193 }
194}
195
196
197template<class ParticleType>
198template<class DataType>
200(
201 const Cloud<ParticleType>& c,
202 const CompactIOField<Field<DataType>>& data
203) const
204{
205 if (data.size() != c.size())
206 {
208 << "Size of " << data.name()
209 << " field " << data.size()
210 << " does not match the number of particles " << c.size()
212 }
213}
214
215
216template<class ParticleType>
217template<class Type>
219(
220 const IOobject& io,
221 const IOobject& ioNew
222) const
223{
224 if (io.isHeaderClass<IOField<Type>>())
225 {
226 IOField<Type> fld(io);
227 auto* fldNewPtr = new IOField<Type>(ioNew, std::move(fld));
228 return fldNewPtr->store();
230
231 return false;
232}
233
234
235template<class ParticleType>
237(
238 objectRegistry& obr,
239 const wordRes& selectFields,
240 const wordRes& excludeFields
241) const
242{
243 IOobjectList cloudObjects
244 (
245 *this,
246 time().timeName(),
248 );
249
250 const wordRes::filter pred(selectFields, excludeFields);
251
252 forAllConstIters(cloudObjects, iter)
253 {
254 const IOobject& io = *(iter.val());
255 const word& fldName = io.name();
256
257 if (!pred(fldName))
258 {
259 continue; // reject
260 }
261
262 IOobject ioNew
263 (
264 fldName,
265 time().timeName(),
266 obr,
269 );
270
271 const bool stored
272 (
273 readStoreFile<label>(io, ioNew)
274 || readStoreFile<scalar>(io, ioNew)
275 || readStoreFile<vector>(io, ioNew)
276 || readStoreFile<sphericalTensor>(io, ioNew)
277 || readStoreFile<symmTensor>(io, ioNew)
278 || readStoreFile<tensor>(io, ioNew)
279 );
280
281 if (!stored)
282 {
284 << "Unhandled field:" << fldName
285 << " type:" << io.headerClassName() << endl;
286 }
287 }
288}
289
290
291template<class ParticleType>
294 ParticleType::writeFields(*this);
295}
296
297
298template<class ParticleType>
300(
301 IOstreamOption streamOpt,
302 const bool /* writeOnProc */
303) const
304{
305 writeCloudUniformProperties();
306
307 writeFields();
308 return cloud::writeObject(streamOpt, (this->size() > 0));
309}
310
311
312// * * * * * * * * * * * * * * * Ostream Operators * * * * * * * * * * * * * //
313
314template<class ParticleType>
316{
317 c.writeData(os);
318
319 os.check(FUNCTION_NAME);
320 return os;
321}
322
323
324// ************************************************************************* //
wordRes excludeFields
const word cloudName(propsDict.get< word >("cloud"))
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Base cloud calls templated on particle type.
Definition Cloud.H:64
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc) const
Write using stream options.
Definition CloudIO.C:293
bool readStoreFile(const IOobject &io, const IOobject &ioNew) const
Helper function to store a cloud field on its registry.
Definition CloudIO.C:212
static word cloudPropertiesName
Name of cloud properties dictionary.
Definition Cloud.H:149
void readFromFiles(objectRegistry &obr, const wordRes &selectFields, const wordRes &excludeFields=wordRes::null()) const
Read from files into objectRegistry.
Definition CloudIO.C:230
void checkFieldIOobject(const Cloud< ParticleType > &c, const IOField< DataType > &data) const
Check lagrangian data field.
Definition CloudIO.C:174
Cloud(const polyMesh &mesh, const Foam::zero, const word &cloudName)
Construct without particles.
Definition Cloud.C:57
virtual void writeFields() const
Write the field data for the cloud of particles Dummy at.
Definition CloudIO.C:285
const polyMesh & pMesh() const noexcept
Return the polyMesh reference.
Definition Cloud.H:193
void checkFieldFieldIOobject(const Cloud< ParticleType > &c, const CompactIOField< Field< DataType > > &data) const
Check lagrangian data fieldfield.
Definition CloudIO.C:193
A Field of objects of type <T> with automated input and output using a compact storage....
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
A primitive field of type <T> with automated input and output.
Definition IOField.H:53
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable,...
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
A simple container for options an IOstream can normally have.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Registry of regIOobjects.
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc) const
Write the objects using stream options.
const Time & time() const noexcept
Return time registry.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const auto & io
auto & name
word timeName
Definition getTimeIndex.H:3
#define DebugInfo
Report an information message using Foam::Info.
#define FUNCTION_NAME
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition List.H:62
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
errorManip< error > abort(error &err)
Definition errorManip.H:139
void writeFields(const fvMesh &mesh, const wordHashSet &selectedFields, const bool writeFaceFields)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235
Functor wrapper of allow/deny lists of wordRe for filtering.
Definition wordRes.H:275