Loading...
Searching...
No Matches
sampledSetsImpl.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 OpenFOAM Foundation
9 Copyright (C) 2015-2023 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 "sampledSets.H"
30#include "globalIndex.H"
31#include "interpolation.H"
32#include "volFields.H"
34
35// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36
37template<class GeoField>
39Foam::sampledSets::getOrLoadField(const word& fieldName) const
40{
41 tmp<GeoField> tfield;
42
43 if (loadFromFiles_)
44 {
45 tfield.emplace
46 (
48 (
49 fieldName,
50 mesh_.time().timeName(),
51 mesh_.thisDb(),
55 ),
56 mesh_
57 );
58 }
59 else
60 {
61 tfield.cref(mesh_.cfindObject<GeoField>(fieldName));
62 }
63
64 return tfield;
65}
66
67
68template<class Type>
69void Foam::sampledSets::writeCoordSet
70(
72 const Field<Type>& values,
73 const word& fieldName
74)
75{
77 if (Pstream::master())
78 {
79 outputName = writer.write(fieldName, values);
80 }
82
83 if (outputName.size())
84 {
85 // Case-local file name with "<case>" to make relocatable
86
88 propsDict.add
89 (
90 "file",
91 time_.relativePath(outputName, true)
92 );
93 setProperty(fieldName, propsDict);
94 }
95}
96
97
98template<class Type>
99void Foam::sampledSets::performAction
100(
101 const VolumeField<Type>& fld,
102 unsigned request
103)
104{
105 const word& fieldName = fld.name();
106 const scalar timeValue = fld.time().timeOutputValue();
107
108 // The interpolator for this field
110
111 if (!samplePointScheme_.empty() && samplePointScheme_ != "cell")
112 {
113 interpPtr.reset(interpolation<Type>::New(samplePointScheme_, fld));
114 }
115
116 const unsigned int width(IOstream::defaultPrecision() + 7);
117 OFstream* osptr = nullptr;
118
119 if (writeAsProbes_ && (request & ACTION_WRITE))
120 {
121 osptr = createProbeFile(fieldName);
122
123 if (Pstream::master() && osptr)
124 {
125 (*osptr) << setw(width) << timeValue;
126 }
127 }
128
129 // Ensemble min/max/avg values
130 Type avgEnsemble = Zero;
131 label sizeEnsemble = 0;
132 MinMax<Type> limitsEnsemble;
133
134 forAll(*this, seti)
135 {
136 const sampledSet& s = (*this)[seti];
137 const globalIndex& globIdx = globalIndices_[seti];
138 const labelList& globOrder = gatheredSorting_[seti];
139
140 const word& setName = s.name();
141 Field<Type> values(s.size());
142
143 if (interpPtr)
144 {
145 forAll(s, samplei)
146 {
147 const point& p = s[samplei];
148 const label celli = s.cells()[samplei];
149 const label facei = s.faces()[samplei];
150
151 if (celli == -1 && facei == -1)
152 {
153 // Special condition for illegal sampling points
154 values[samplei] = pTraits<Type>::max;
155 }
156 else
157 {
158 values[samplei] = interpPtr().interpolate(p, celli, facei);
159 }
160 }
161 }
162 else
163 {
164 forAll(s, samplei)
165 {
166 const label celli = s.cells()[samplei];
167
168 if (celli == -1)
169 {
170 values[samplei] = pTraits<Type>::max;
171 }
172 else
173 {
174 values[samplei] = fld[celli];
175 }
176 }
177 }
178
179 // Collect data from all processors
180 globIdx.gatherInplace(values);
181
182 // Local min/max/avg values - calculate on master
183 Type avgValue = Zero;
184 label sizeValue = 0;
186
187 if (Pstream::master())
188 {
189 avgValue = sum(values);
190 sizeValue = values.size();
191 limits = MinMax<Type>(values);
192
193 // Ensemble values
194 avgEnsemble += avgValue;
195 sizeEnsemble += sizeValue;
196 limitsEnsemble += limits;
197
198 if (sizeValue)
199 {
200 avgValue /= sizeValue;
201 }
202
203 // Use sorted order
204 values = UIndirectList<Type>(values, globOrder)();
205 }
206 Pstream::broadcasts(UPstream::worldComm, avgValue, sizeValue, limits);
207
208 // Store results: min/max/average/size with the name of the set
209 // for scoping.
210 // Eg, average(lines,T) ...
211 const word resultArg('(' + setName + ',' + fieldName + ')');
212
213 this->setResult("average" + resultArg, avgValue);
214 this->setResult("min" + resultArg, limits.min());
215 this->setResult("max" + resultArg, limits.max());
216 this->setResult("size" + resultArg, sizeValue);
217
218 if (verbose_)
219 {
220 Info<< name() << ' ' << setName << " : " << fieldName << nl
221 << " avg: " << avgValue << nl
222 << " min: " << limits.min() << nl
223 << " max: " << limits.max() << nl << nl;
224 }
225
226 if ((request & ACTION_WRITE) != 0)
227 {
228 if (writeAsProbes_)
229 {
230 if (osptr)
231 {
232 for (const Type& val : values)
233 {
234 (*osptr) << ' ' << setw(width) << val;
235 }
236 }
237 }
238 else
239 {
240 writeCoordSet<Type>(writers_[seti], values, fieldName);
241 }
242 }
243 }
244
245 // Finish probes write
246 if (Pstream::master() && osptr)
247 {
248 (*osptr) << endl;
249 }
250
251 if (sizeEnsemble)
252 {
253 avgEnsemble /= sizeEnsemble;
254 }
255
256 if (size())
257 {
259 (
261 avgEnsemble,
262 sizeEnsemble,
263 limitsEnsemble
264 );
265
266 // Store results: min/max/average/size for the ensemble
267 // Eg, average(T) ...
268 const word resultArg('(' + fieldName + ')');
269
270 this->setResult("average" + resultArg, avgEnsemble);
271 this->setResult("min" + resultArg, limitsEnsemble.min());
272 this->setResult("max" + resultArg, limitsEnsemble.max());
273 this->setResult("size" + resultArg, sizeEnsemble);
274 }
275}
276
277
278template<class GeoField>
279void Foam::sampledSets::performAction
280(
281 const IOobjectList& objects,
282 unsigned request
283)
284{
285 wordList fieldNames;
286 if (loadFromFiles_)
287 {
288 fieldNames = objects.sortedNames<GeoField>(fieldSelection_);
289 }
290 else
291 {
292 fieldNames = mesh_.thisDb().sortedNames<GeoField>(fieldSelection_);
293 }
294
295 for (const word& fieldName : fieldNames)
296 {
297 tmp<GeoField> tfield = getOrLoadField<GeoField>(fieldName);
298
299 if (tfield)
300 {
301 performAction<typename GeoField::value_type>(tfield(), request);
302 }
303 }
304}
305
306
307// ************************************************************************* //
IOdictionary propsDict(dictIO)
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))
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edgesCentres")))
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable,...
@ NO_REGISTER
Do not request registration (bool: false).
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition IOstream.H:437
A min/max value pair with additional methods. In addition to conveniently storing values,...
Definition MinMax.H:106
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
static void broadcasts(const int communicator, Type &value, Args &&... values)
Broadcast multiple items to all communicator ranks. Does nothing in non-parallel.
A List with indirect addressing. Like IndirectList but does not store addressing.
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static label worldComm
Communicator for all ranks. May differ from commGlobal() if local worlds are in use.
Definition UPstream.H:1069
@ broadcast
broadcast [MPI]
Definition UPstream.H:189
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Base class for writing coordSet(s) and tracks with fields.
A class for handling file names.
Definition fileName.H:75
const fvMesh & mesh_
Reference to the fvMesh.
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
static autoPtr< interpolation< Type > > New(const word &interpolationType, const GeometricField< Type, fvPatchField, volMesh > &psi)
Return a reference to the specified interpolation scheme.
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition sampledSet.H:82
A class for managing temporary objects.
Definition tmp.H:75
T & emplace(Args &&... args)
Reset with emplace construction. Return reference to the new content.
Definition tmpI.H:357
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition tmpI.H:221
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
auto limits
Definition setRDeltaT.H:186
word outputName("finiteArea-edges.obj")
auto & name
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition HashOps.H:164
List< word > wordList
List of word.
Definition fileName.H:60
GeometricField< Type, fvPatchField, volMesh > VolumeField
A volume field for a given type.
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
Omanip< int > setw(const int i)
Definition IOmanip.H:199
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
vector point
Point is a vector.
Definition point.H:37
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
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