MappedFile.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) 2018-2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "polyMesh.H"
29#include "rawIOField.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class Type>
35(
36 const polyPatch& pp,
37 const word& redirectType,
38 const word& entryName,
39 const dictionary& dict,
40 const bool faceValues
41)
42:
43 PatchFunction1<Type>(pp, entryName, dict, faceValues),
44 dictConstructed_(true),
45 setAverage_(dict.getOrDefault("setAverage", false)),
46 fieldTableName_(dict.getOrDefault<word>("fieldTable", entryName)),
47 perturb_(dict.getOrDefault<scalar>("perturb", 1e-5)),
48 pointsName_(dict.getOrDefault<word>("points", "points")),
49 mapMethod_
50 (
51 dict.getOrDefault<word>
52 (
53 "mapMethod",
54 "planarInterpolation"
55 )
56 ),
57 mapperPtr_(nullptr),
58 sampleTimes_(0),
59 startSampleTime_(-1),
60 startSampledValues_(0),
61 startAverage_(Zero),
62 endSampleTime_(-1),
63 endSampledValues_(0),
64 endAverage_(Zero),
65 offset_(Function1<Type>::NewIfPresent("offset", dict))
66{
67 if
68 (
69 mapMethod_ != "planarInterpolation"
70 && mapMethod_ != "nearest"
71 )
72 {
74 << "mapMethod should be one of 'planarInterpolation'"
75 << ", 'nearest'" << exit(FatalIOError);
76 }
77}
78
79
80template<class Type>
82(
83 const polyPatch& pp,
84 const word& entryName,
85 const dictionary& dict,
86 const word& fieldTableName,
87 const bool faceValues
88)
89:
90 PatchFunction1<Type>(pp, entryName, dict, faceValues),
91 dictConstructed_(false),
92 setAverage_(dict.getOrDefault("setAverage", false)),
93 fieldTableName_(fieldTableName),
94 perturb_(dict.getOrDefault<scalar>("perturb", 1e-5)),
95 pointsName_(dict.getOrDefault<word>("points", "points")),
96 mapMethod_
97 (
98 dict.getOrDefault<word>
99 (
100 "mapMethod",
101 "planarInterpolation"
102 )
103 ),
104 mapperPtr_(nullptr),
105 sampleTimes_(0),
106 startSampleTime_(-1),
107 startSampledValues_(0),
108 startAverage_(Zero),
109 endSampleTime_(-1),
110 endSampledValues_(0),
111 endAverage_(Zero),
112 offset_(Function1<Type>::NewIfPresent("offset", dict))
113{
114 if
115 (
116 mapMethod_ != "planarInterpolation"
117 && mapMethod_ != "nearest"
118 )
119 {
121 << "mapMethod should be one of 'planarInterpolation'"
122 << ", 'nearest'" << exit(FatalIOError);
123 }
124}
125
126
127template<class Type>
129(
130 const MappedFile<Type>& rhs
131)
132:
133 MappedFile<Type>(rhs, rhs.patch())
134{}
135
136
137template<class Type>
139(
140 const MappedFile<Type>& rhs,
141 const polyPatch& pp
142)
143:
144 PatchFunction1<Type>(rhs, pp),
145 dictConstructed_(rhs.dictConstructed_),
146 setAverage_(rhs.setAverage_),
147 fieldTableName_(rhs.fieldTableName_),
148 perturb_(rhs.perturb_),
149 pointsName_(rhs.pointsName_),
150 mapMethod_(rhs.mapMethod_),
151 mapperPtr_(rhs.mapperPtr_.clone()),
152 sampleTimes_(rhs.sampleTimes_),
153 startSampleTime_(rhs.startSampleTime_),
154 startSampledValues_(rhs.startSampledValues_),
155 startAverage_(rhs.startAverage_),
156 endSampleTime_(rhs.endSampleTime_),
157 endSampledValues_(rhs.endSampledValues_),
158 endAverage_(rhs.endAverage_),
159 offset_(rhs.offset_.clone())
160{}
161
162
163// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
164
165template<class Type>
167(
168 const FieldMapper& mapper
169)
170{
172
173 if (startSampledValues_.size())
174 {
175 startSampledValues_.autoMap(mapper);
176 }
177
178 if (endSampledValues_.size())
179 {
180 endSampledValues_.autoMap(mapper);
181 }
182
183 // Clear interpolator
184 mapperPtr_.clear();
185 startSampleTime_ = -1;
186 endSampleTime_ = -1;
187}
188
189
190template<class Type>
192(
193 const PatchFunction1<Type>& pf1,
194 const labelList& addr
195)
196{
198
200 refCast<const PatchFunction1Types::MappedFile<Type>>(pf1);
201
202 if (tiptf.startSampledValues_.size())
203 {
204 startSampledValues_.setSize(this->size());
205 startSampledValues_.rmap(tiptf.startSampledValues_, addr);
206 }
207
208 if (tiptf.endSampledValues_.size())
209 {
210 endSampledValues_.setSize(this->size());
211 endSampledValues_.rmap(tiptf.endSampledValues_, addr);
212 }
213
214 // Clear interpolator
215 mapperPtr_.clear();
216 startSampleTime_ = -1;
217 endSampleTime_ = -1;
218}
219
221template<class Type>
224 const scalar t
225) const
226{
227 const polyMesh& mesh = this->patch_.boundaryMesh().mesh();
228 const Time& time = mesh.time();
229
230 // Initialise
231 if (!mapperPtr_)
232 {
233 // Reread values and interpolate
234 const fileName samplePointsFile
235 (
236 time.globalPath()
237 /time.constant() // instance
238 /mesh.dbDir() // region
239 /"boundaryData"
240 /this->patch_.name()
241 /pointsName_
242 );
243
245 (
246 samplePointsFile, // absolute path
247 time,
250 false, // no need to register
251 true // is global object (currently not used)
252 );
253
254 // Read data
255 const rawIOField<point> samplePoints(io, false);
256
258 << "Read " << samplePoints.size() << " sample points from "
259 << samplePointsFile << endl;
260
261
262 // tbd: run-time selection
263 bool nearestOnly =
264 (
265 !mapMethod_.empty()
266 && mapMethod_ != "planarInterpolation"
267 );
268
269 // Allocate the interpolator
270 if (this->faceValues())
271 {
272 mapperPtr_.reset
275 (
276 samplePoints,
277 this->localPosition(this->patch_.faceCentres()),
278 perturb_,
279 nearestOnly
280 )
281 );
282 }
283 else
284 {
285 mapperPtr_.reset
288 (
289 samplePoints,
290 this->localPosition(this->patch_.localPoints()),
291 perturb_,
292 nearestOnly
293 )
294 );
295 }
297
298 // Read the times for which data is available
299 const fileName samplePointsDir = samplePointsFile.path();
300 sampleTimes_ = Time::findTimes(samplePointsDir);
301
303 << "In directory "
304 << samplePointsDir << " found times "
306 << endl;
307 }
308
309
310 // Find current time in sampleTimes
311 label lo = -1;
312 label hi = -1;
313
314 bool foundTime = mapperPtr_().findTime
315 (
316 sampleTimes_,
317 startSampleTime_,
318 t, //mesh.time().value(),
319 lo,
320 hi
321 );
322
323 if (!foundTime)
324 {
326 << "Cannot find starting sampling values for index "
327 << t << nl
328 << "Have sampling values for "
330 << "In directory "
331 << time.constant()/mesh.dbDir()/"boundaryData"/this->patch_.name()
332 << "\n on patch " << this->patch_.name()
333 << " of field " << fieldTableName_
334 << exit(FatalError);
335 }
336
337
338 // Update sampled data fields.
339
340 if (lo != startSampleTime_)
341 {
342 startSampleTime_ = lo;
343
344 if (startSampleTime_ == endSampleTime_)
345 {
346 // No need to reread since are end values
347 if (debug)
348 {
349 Pout<< "checkTable : Setting startValues to (already read) "
350 << "boundaryData"
351 /this->patch_.name()
352 /sampleTimes_[startSampleTime_].name()
353 << endl;
354 }
355 startSampledValues_ = endSampledValues_;
356 startAverage_ = endAverage_;
357 }
358 else
359 {
360 if (debug)
361 {
362 Pout<< "checkTable : Reading startValues from "
363 << "boundaryData"
364 /this->patch_.name()
365 /sampleTimes_[lo].name()
366 /fieldTableName_
367 << endl;
368 }
369
370 // Reread values and interpolate
371 const fileName valsFile
372 (
373 time.globalPath()
374 /time.constant()
375 /mesh.dbDir() // region
376 /"boundaryData"
377 /this->patch_.name()
378 /sampleTimes_[startSampleTime_].name()
379 /fieldTableName_
380 );
381
382 IOobject io
383 (
384 valsFile, // absolute path
385 time,
388 false, // no need to register
389 true // is global object (currently not used)
390 );
391
392 const rawIOField<Type> vals(io, setAverage_);
393 if (setAverage_)
394 {
395 startAverage_ = vals.average();
396 }
397
398 if (vals.size() != mapperPtr_().sourceSize())
399 {
401 << "Number of values (" << vals.size()
402 << ") differs from the number of points ("
403 << mapperPtr_().sourceSize()
404 << ") in file " << valsFile << exit(FatalError);
405 }
406
407 startSampledValues_ = mapperPtr_().interpolate(vals);
408 }
409 }
410
411 if (hi != endSampleTime_)
412 {
413 endSampleTime_ = hi;
414
415 if (endSampleTime_ == -1)
416 {
417 // endTime no longer valid. Might as well clear endValues.
418 if (debug)
419 {
420 Pout<< "checkTable : Clearing endValues" << endl;
421 }
422 endSampledValues_.clear();
423 }
424 else
425 {
426 if (debug)
427 {
428 Pout<< "checkTable : Reading endValues from "
429 << "boundaryData"
430 /this->patch_.name()
431 /sampleTimes_[endSampleTime_].name()
432 << endl;
433 }
434
435 // Reread values and interpolate
436 fileName valsFile
437 (
438 time.globalPath()
439 /time.constant()
440 /mesh.dbDir() // region
441 /"boundaryData"
442 /this->patch_.name()
443 /sampleTimes_[endSampleTime_].name()
444 /fieldTableName_
445 );
446
447 IOobject io
448 (
449 valsFile, // absolute path
450 time,
453 false, // no need to register
454 true // is global object (currently not used)
455 );
456
457 const rawIOField<Type> vals(io, setAverage_);
458 if (setAverage_)
459 {
460 endAverage_ = vals.average();
461 }
462
463 if (vals.size() != mapperPtr_().sourceSize())
464 {
466 << "Number of values (" << vals.size()
467 << ") differs from the number of points ("
468 << mapperPtr_().sourceSize()
469 << ") in file " << valsFile << exit(FatalError);
470 }
471
472 endSampledValues_ = mapperPtr_().interpolate(vals);
473 }
474 }
475}
476
477
478template<class Type>
481(
482 const scalar x
483) const
484{
485 checkTable(x);
486
487 auto tfld = tmp<Field<Type>>::New(startSampledValues_.size());
488 auto& fld = tfld.ref();
489 Type wantedAverage;
490
491 if (endSampleTime_ == -1)
492 {
493 // Only start value
494 if (debug)
495 {
496 Pout<< "MappedFile<Type>::value : Sampled, non-interpolated values"
497 << " from start time:"
498 << sampleTimes_[startSampleTime_].name() << nl;
499 }
500
501 fld = startSampledValues_;
502 wantedAverage = startAverage_;
503 }
504 else
505 {
506 scalar start = sampleTimes_[startSampleTime_].value();
507 scalar end = sampleTimes_[endSampleTime_].value();
508
509 scalar s = (x - start)/(end - start);
510
511 if (debug)
512 {
513 Pout<< "MappedFile<Type>::value : Sampled, interpolated values"
514 << " between start time:"
515 << sampleTimes_[startSampleTime_].name()
516 << " and end time:" << sampleTimes_[endSampleTime_].name()
517 << " with weight:" << s << endl;
518 }
519
520 fld = ((1 - s)*startSampledValues_ + s*endSampledValues_);
521 wantedAverage = (1 - s)*startAverage_ + s*endAverage_;
522 }
523
524 // Enforce average. Either by scaling (if scaling factor > 0.5) or by
525 // offsetting.
526 if (setAverage_)
527 {
528 Type averagePsi;
529 if (this->faceValues())
530 {
531 const scalarField magSf(mag(this->patch_.faceAreas()));
532 averagePsi = gSum(magSf*fld)/gSum(magSf);
533 }
534 else
535 {
536 averagePsi = gAverage(fld);
537 }
538
539 if (debug)
540 {
541 Pout<< "MappedFile<Type>::value :"
542 << " actual average:" << averagePsi
543 << " wanted average:" << wantedAverage
544 << endl;
545 }
546
547 if (mag(averagePsi) < VSMALL)
548 {
549 // Field too small to scale. Offset instead.
550 const Type offset = wantedAverage - averagePsi;
551 if (debug)
552 {
553 Pout<< "MappedFile<Type>::value :"
554 << " offsetting with:" << offset << endl;
555 }
556 fld += offset;
557 }
558 else
559 {
560 const scalar scale = mag(wantedAverage)/mag(averagePsi);
561
562 if (debug)
563 {
564 Pout<< "MappedFile<Type>::value :"
565 << " scaling with:" << scale << endl;
566 }
567 fld *= scale;
568 }
569 }
570
571 // Apply offset to mapped values
572 if (offset_)
573 {
574 fld += offset_->value(x);
575 }
576
577 if (debug)
578 {
579 Pout<< "MappedFile<Type>::value : set fixedValue to min:" << gMin(fld)
580 << " max:" << gMax(fld)
581 << " avg:" << gAverage(fld) << endl;
582 }
583
584 return this->transform(tfld);
585}
586
587
588template<class Type>
591(
592 const scalar x1,
593 const scalar x2
594) const
595{
597 return nullptr;
598}
599
600
601template<class Type>
603(
604 Ostream& os
605) const
606{
607 os.writeEntryIfDifferent
608 (
609 "fieldTable",
610 this->name(),
611 fieldTableName_
612 );
613
614 if (setAverage_)
615 {
616 os.writeEntry("setAverage", setAverage_);
617 }
618
619 os.writeEntryIfDifferent<scalar>("perturb", 1e-5, perturb_);
620
621 os.writeEntryIfDifferent<word>("points", "points", pointsName_);
622
623 os.writeEntryIfDifferent<word>
624 (
625 "mapMethod",
626 "planarInterpolation",
627 mapMethod_
628 );
629
630 if (offset_)
631 {
632 offset_->writeData(os);
633 }
634}
635
636
637template<class Type>
639(
640 Ostream& os
641) const
642{
644
645 // Check if field name explicitly provided
646 // (e.g. through timeVaryingMapped bc)
647 if (dictConstructed_)
648 {
649 os.writeEntry(this->name(), type());
650
651 os.beginBlock(word(this->name() + "Coeffs"));
652 writeEntries(os);
653 os.endBlock();
654 }
655 else
656 {
657 // Note that usually dictConstructed = true. The
658 // construct-from-dictionary (dictConstructed_ = false)
659 // should only be used if there is only
660 // a single potential MappedFile in the local scope.
661 writeEntries(os);
662 }
663}
664
665
666// ************************************************************************* //
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))
Abstract base class to hold the Field mapping addressing and weights.
Definition: FieldMapper.H:50
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:96
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
virtual const fileName & name() const
Get the name of the stream.
Definition: OSstream.H:107
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
virtual Ostream & beginBlock(const keyType &kw)
Write begin block group with the given name.
Definition: Ostream.C:87
Patch value mapping from a set of values stored in a file and a set of unstructured points using the ...
Definition: MappedFile.H:132
virtual tmp< Field< Type > > integrate(const scalar x1, const scalar x2) const
Integrate between two values.
Definition: MappedFile.C:591
virtual void rmap(const PatchFunction1< Type > &pf1, const labelList &addr)
Reverse map the given PatchFunction1 onto this PatchFunction1.
Definition: MappedFile.C:192
virtual void autoMap(const FieldMapper &mapper)
Map (and resize as needed) from self given a mapping object.
Definition: MappedFile.C:167
virtual void writeEntries(Ostream &os) const
Write coefficient entries in dictionary format.
Definition: MappedFile.C:603
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
virtual void rmap(const PatchFunction1< Type > &rhs, const labelList &addr)
Reverse map the given PatchFunction1 onto this PatchFunction1.
virtual void autoMap(const FieldMapper &mapper)
Map (and resize as needed) from self given a mapping object.
const polyPatch const word const word const dictionary & dict
static instantList findTimes(const fileName &directory, const word &constantName="constant")
Search a given directory for valid time directories.
Definition: TimePaths.C:140
const word & constant() const
Return constant name.
Definition: TimePathsI.H:96
fileName globalPath() const
Return global path for the case.
Definition: TimePathsI.H:80
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A class for handling file names.
Definition: fileName.H:76
static std::string path(const std::string &str)
Return directory path name (part before last /)
Definition: fileNameI.H:176
static std::string name(const std::string &str)
Return basename (part beyond last /), including its extension.
Definition: fileNameI.H:199
Interpolates between two sets of unstructured points using 2D Delaunay triangulation....
static wordList timeNames(const instantList &)
Helper: extract words of times.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
Like IOField but falls back to raw IFstream if no header found. Optionally reads average value....
Definition: rawIOField.H:56
void checkTable()
Find boundary data inbetween current time and interpolate.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
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))
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
#define DebugInfo
Report an information message using Foam::Info.
Type gSum(const FieldField< Field, Type > &f)
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
IOerror FatalIOError
Type gAverage(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
Type gMin(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Type gMax(const FieldField< Field, Type > &f)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
volScalarField & e
Definition: createFields.H:11