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-2022 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.reset
46 (
47 new GeoField
48 (
50 (
51 fieldName,
53 mesh_,
56 false
57 ),
58 mesh_
59 )
60 );
61 }
62 else
63 {
64 // Slightly paranoid here
65 tfield.cref(mesh_.cfindObject<GeoField>(fieldName));
66 }
67
68 return tfield;
69}
70
71
72template<class Type>
73void Foam::sampledSets::writeCoordSet
74(
76 const Field<Type>& values,
77 const word& fieldName
78)
79{
81 if (Pstream::master())
82 {
83 outputName = writer.write(fieldName, values);
84 }
86
87 if (outputName.size())
88 {
89 // Case-local file name with "<case>" to make relocatable
90
91 dictionary propsDict;
93 (
94 "file",
95 time_.relativePath(outputName, true)
96 );
97 setProperty(fieldName, propsDict);
98 }
99}
100
101
102template<class Type>
103void Foam::sampledSets::performAction
104(
105 const VolumeField<Type>& fld,
106 unsigned request
107)
108{
109 const word& fieldName = fld.name();
110 const scalar timeValue = fld.time().timeOutputValue();
111
112 // The interpolator for this field
113 autoPtr<interpolation<Type>> interpPtr;
114
115 if (!samplePointScheme_.empty() && samplePointScheme_ != "cell")
116 {
117 interpPtr.reset(interpolation<Type>::New(samplePointScheme_, fld));
118 }
119
120 const unsigned int width(IOstream::defaultPrecision() + 7);
121 OFstream* osptr = nullptr;
122
123 if (writeAsProbes_ && (request & ACTION_WRITE))
124 {
125 osptr = createProbeFile(fieldName);
126
127 if (Pstream::master() && osptr)
128 {
129 (*osptr) << setw(width) << timeValue;
130 }
131 }
132
133 // Ensemble min/max/avg values
134 Type avgEnsemble = Zero;
135 label sizeEnsemble = 0;
136 MinMax<Type> limitsEnsemble;
137
138 forAll(*this, seti)
139 {
140 const sampledSet& s = (*this)[seti];
141 const globalIndex& globIdx = globalIndices_[seti];
142 const labelList& globOrder = gatheredSorting_[seti];
143
144 const word& setName = s.name();
145 Field<Type> values(s.size());
146
147 if (interpPtr)
148 {
149 forAll(s, samplei)
150 {
151 const point& p = s[samplei];
152 const label celli = s.cells()[samplei];
153 const label facei = s.faces()[samplei];
154
155 if (celli == -1 && facei == -1)
156 {
157 // Special condition for illegal sampling points
158 values[samplei] = pTraits<Type>::max;
159 }
160 else
161 {
162 values[samplei] = interpPtr().interpolate(p, celli, facei);
163 }
164 }
165 }
166 else
167 {
168 forAll(s, samplei)
169 {
170 const label celli = s.cells()[samplei];
171
172 if (celli == -1)
173 {
174 values[samplei] = pTraits<Type>::max;
175 }
176 else
177 {
178 values[samplei] = fld[celli];
179 }
180 }
181 }
182
183 // Collect data from all processors
184 globIdx.gatherInplace(values);
185
186 // Local min/max/avg values - calculate on master
187 Type avgValue = Zero;
188 label sizeValue = 0;
189 MinMax<Type> limits;
190
191 if (Pstream::master())
192 {
193 avgValue = sum(values);
194 sizeValue = values.size();
195 limits = MinMax<Type>(values);
196
197 // Ensemble values
198 avgEnsemble += avgValue;
199 sizeEnsemble += sizeValue;
200 limitsEnsemble += limits;
201
202 if (sizeValue)
203 {
204 avgValue /= sizeValue;
205 }
206
207 // Use sorted order
208 values = UIndirectList<Type>(values, globOrder)();
209 }
210 Pstream::broadcasts(UPstream::worldComm, avgValue, sizeValue, limits);
211
212 // Store results: min/max/average/size with the name of the set
213 // for scoping.
214 // Eg, average(lines,T) ...
215 const word resultArg('(' + setName + ',' + fieldName + ')');
216
217 this->setResult("average" + resultArg, avgValue);
218 this->setResult("min" + resultArg, limits.min());
219 this->setResult("max" + resultArg, limits.max());
220 this->setResult("size" + resultArg, sizeValue);
221
222 if (verbose_)
223 {
224 Info<< name() << ' ' << setName << " : " << fieldName << nl
225 << " avg: " << avgValue << nl
226 << " min: " << limits.min() << nl
227 << " max: " << limits.max() << nl << nl;
228 }
229
230 if ((request & ACTION_WRITE) != 0)
231 {
232 if (writeAsProbes_)
233 {
234 if (osptr)
235 {
236 for (const Type& val : values)
237 {
238 (*osptr) << ' ' << setw(width) << val;
239 }
240 }
241 }
242 else
243 {
244 writeCoordSet<Type>(writers_[seti], values, fieldName);
245 }
246 }
247 }
248
249 // Finish probes write
250 if (Pstream::master() && osptr)
251 {
252 (*osptr) << endl;
253 }
254
255 if (sizeEnsemble)
256 {
257 avgEnsemble /= sizeEnsemble;
258 }
259
260 if (size())
261 {
263 (
265 avgEnsemble,
266 sizeEnsemble,
267 limitsEnsemble
268 );
269
270 // Store results: min/max/average/size for the ensemble
271 // Eg, average(T) ...
272 const word resultArg('(' + fieldName + ')');
273
274 this->setResult("average" + resultArg, avgEnsemble);
275 this->setResult("min" + resultArg, limitsEnsemble.min());
276 this->setResult("max" + resultArg, limitsEnsemble.max());
277 this->setResult("size" + resultArg, sizeEnsemble);
278 }
279}
280
281
282template<class GeoField>
283void Foam::sampledSets::performAction
284(
285 const IOobjectList& objects,
286 unsigned request
287)
288{
289 wordList fieldNames;
290 if (loadFromFiles_)
291 {
292 fieldNames = objects.sortedNames<GeoField>(fieldSelection_);
293 }
294 else
295 {
296 fieldNames = mesh_.thisDb().sortedNames<GeoField>(fieldSelection_);
297 }
298
299 for (const word& fieldName : fieldNames)
300 {
301 tmp<GeoField> tfield = getOrLoadField<GeoField>(fieldName);
302
303 if (tfield)
304 {
305 performAction<typename GeoField::value_type>(tfield(), request);
306 }
307 }
308}
309
310
311// ************************************************************************* //
vtk::internalMeshWriter writer(topoMesh, topoCells, vtk::formatType::INLINE_ASCII, runTime.path()/"blockTopology")
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Generic templated field type.
Definition: Field.H:82
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition: IOstream.H:342
static void broadcast(Type &value, const label comm=UPstream::worldComm)
static void broadcasts(const label comm, Type &arg1, Args &&... args)
Broadcast multiple items to all processes in communicator.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
static label worldComm
Default communicator (all processors)
Definition: UPstream.H:293
Base class for writing coordSet(s) and tracks with fields.
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:640
A class for handling file names.
Definition: fileName.H:76
const fvMesh & mesh_
Reference to the fvMesh.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
const Type * cfindObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
splitCell * master() const
Definition: splitCell.H:113
A class for managing temporary objects.
Definition: tmp.H:65
const T & cref() const
Definition: tmpI.H:213
void reset(tmp< T > &&other) noexcept
Clear existing and transfer ownership.
Definition: tmpI.H:314
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
word outputName("finiteArea-edges.obj")
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:149
List< word > wordList
A List of words.
Definition: fileName.H:63
List< label > labelList
A List of labels.
Definition: List.H:66
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
messageStream Info
Information stream (stdout output on master, null elsewhere)
vector point
Point is a vector.
Definition: point.H:43
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
IOdictionary propsDict(dictIO)