Loading...
Searching...
No Matches
sampledSets.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 "dictionary.H"
31#include "Time.H"
32#include "globalIndex.H"
33#include "volFields.H"
34#include "mapPolyMesh.H"
35#include "IOmanip.H"
36#include "IOobjectList.H"
37#include "IndirectList.H"
38#include "ListOps.H"
41// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42
43namespace Foam
44{
46
48 (
52 );
53}
54
55#include "sampledSetsImpl.C"
56
57// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
58
59Foam::autoPtr<Foam::coordSetWriter> Foam::sampledSets::newWriter
60(
61 word writerType,
62 const dictionary& topDict,
63 const dictionary& setDict
64)
65{
66 // Per-set adjustment
67 setDict.readIfPresent<word>("setFormat", writerType);
68
70 (
71 writerType,
72 // Top-level/set-specific "formatOptions"
73 coordSetWriter::formatOptions(topDict, setDict, writerType)
74 );
75}
76
77
78Foam::OFstream* Foam::sampledSets::createProbeFile(const word& fieldName)
79{
80 // Open new output stream
81
82 OFstream* osptr = probeFilePtrs_.lookup(fieldName, nullptr);
83
84 if (!osptr && Pstream::master())
85 {
86 // Put in undecomposed case
87 // (Note: gives problems for distributed data running)
88
89 fileName probeDir
90 (
91 mesh_.time().globalPath()
93 / name()/mesh_.regionName()
94 // Use startTime as the instance for output files
95 / mesh_.time().timeName(mesh_.time().startTime().value())
96 );
97 probeDir.clean(); // Remove unneeded ".."
98
99 // Create directory if needed
100 Foam::mkDir(probeDir);
101
102 probeFilePtrs_.insert
103 (
104 fieldName,
105 autoPtr<OFstream>::New(probeDir/fieldName)
106 );
107 osptr = probeFilePtrs_.lookup(fieldName, nullptr);
108
109 if (osptr)
110 {
111 auto& os = *osptr;
112
113 DebugInfo<< "open probe stream: " << os.name() << endl;
114
115 const unsigned int width(IOstream::defaultPrecision() + 7);
116
117 label nPoints = 0;
118 forAll(*this, seti)
119 {
120 const coordSet& s = gatheredSets_[seti];
121
122 const pointField& pts = static_cast<const pointField&>(s);
123
124 for (const point& p : pts)
125 {
126 os << "# Probe " << nPoints++ << ' ' << p << nl;
127 }
128 }
129
130 os << '#' << setw(IOstream::defaultPrecision() + 6)
131 << "Probe";
132
133 for (label probei = 0; probei < nPoints; ++probei)
134 {
135 os << ' ' << setw(width) << probei;
136 }
137 os << nl;
138
139 os << '#' << setw(IOstream::defaultPrecision() + 6)
140 << "Time" << endl;
141 }
142 }
143
144 return osptr;
145}
146
147
148void Foam::sampledSets::gatherAllSets()
149{
150 // Any writer references will become invalid
151 for (auto& writer : writers_)
152 {
153 writer.expire();
154 }
155
156 const PtrList<sampledSet>& localSets = *this;
157
158 gatheredSets_.resize_null(localSets.size());
159 gatheredSorting_.resize_nocopy(localSets.size());
160 globalIndices_.resize_nocopy(localSets.size());
161
162 forAll(localSets, seti)
163 {
164 const coordSet& coords = localSets[seti];
165
166 globalIndices_[seti].reset(globalIndex::gatherOnly{}, coords.size());
167 gatheredSets_.set(seti, coords.gatherSort(gatheredSorting_[seti]));
168 }
169}
170
171
172Foam::IOobjectList Foam::sampledSets::preCheckFields(unsigned request)
173{
174 wordList allFields; // Just needed for warnings
175 HashTable<wordHashSet> selected;
176
177 IOobjectList objects;
178
179 if (loadFromFiles_)
180 {
181 // Check files for a particular time
182 objects = IOobjectList(mesh_, mesh_.time().timeName());
183
184 allFields = objects.names();
185 selected = objects.classes(fieldSelection_);
186 }
187 else
188 {
189 // Check currently available fields
190 allFields = mesh_.names();
191 selected = mesh_.classes(fieldSelection_);
192 }
193
194 // Parallel consistency (no-op in serial)
195 // Probably not needed...
197
198
199 DynamicList<label> missed(fieldSelection_.size());
200
201 // Detect missing fields
202 forAll(fieldSelection_, i)
203 {
204 if
205 (
206 fieldSelection_[i].isLiteral()
207 && !ListOps::found(allFields, fieldSelection_[i])
208 )
209 {
210 missed.append(i);
211 }
212 }
213
214 if (missed.size() && (request != ACTION_NONE))
215 {
217 << nl
218 << "Cannot find "
219 << (loadFromFiles_ ? "field file" : "registered field")
220 << " matching "
221 << UIndirectList<wordRe>(fieldSelection_, missed) << endl;
222 }
223
224
225 // The selected field names, ordered by (scalar, vector, ...)
226 // with internal sorting
227
228 selectedFieldNames_.clear();
229
230 do
231 {
232 #undef doLocalCode
233 #define doLocalCode(InputType) \
234 { \
235 const auto iter = selected.find(InputType::typeName); \
236 if (iter.good()) \
237 { \
238 selectedFieldNames_.append(iter.val().sortedToc()); \
239 } \
240 }
241
247 #undef doLocalCode
248 }
249 while (false);
250
251
252 // Now propagate field counts (per surface)
253 // - can update writer even when not writing without problem
254
255 const label nFields = selectedFieldNames_.size();
256
257 if (writeAsProbes_)
258 {
259 // Close streams for fields that no longer exist
260 forAllIters(probeFilePtrs_, iter)
261 {
262 if (!selectedFieldNames_.found(iter.key()))
263 {
265 << "close probe stream: "
266 << iter()->name() << endl;
267
268 probeFilePtrs_.remove(iter);
269 }
270 }
271 }
272 else if ((request & ACTION_WRITE) != 0)
273 {
274 forAll(writers_, seti)
275 {
276 coordSetWriter& writer = writers_[seti];
277
278 writer.nFields(nFields);
279 }
280 }
281
282 return objects;
283}
284
285
286void Foam::sampledSets::initDict(const dictionary& dict, const bool initial)
287{
289 if (initial)
290 {
291 writers_.clear();
292 }
293
294 const entry* eptr = dict.findEntry("sets", keyType::LITERAL);
295
296 if (eptr && eptr->isDict())
297 {
298 PtrList<sampledSet> sampSets(eptr->dict().size());
299 if (initial && !writeAsProbes_)
300 {
301 writers_.resize(sampSets.size());
302 }
303
304 label seti = 0;
305
306 for (const entry& dEntry : eptr->dict())
307 {
308 if (!dEntry.isDict())
309 {
310 continue;
311 }
312
313 const dictionary& subDict = dEntry.dict();
314
315 autoPtr<sampledSet> sampSet =
317 (
318 dEntry.keyword(),
319 mesh_,
320 searchEngine_,
321 subDict
322 );
323
324 // if (!sampSet || !sampSet->enabled())
325 // {
326 // continue;
327 // }
328
329 // Define the set
330 sampSets.set(seti, sampSet);
331
332 // Define writer, but do not attached
333 if (initial && !writeAsProbes_)
334 {
335 writers_.set
336 (
337 seti,
338 newWriter(writeFormat_, dict_, subDict)
339 );
340
341 // Use outputDir/TIME/set-name
342 writers_[seti].useTimeDir(true);
343 writers_[seti].verbose(verbose_);
344 }
345 ++seti;
346 }
347
348 sampSets.resize(seti);
349 if (initial && !writeAsProbes_)
350 {
351 writers_.resize(seti);
352 }
353 static_cast<PtrList<sampledSet>&>(*this).transfer(sampSets);
354 }
355 else if (eptr)
356 {
357 // This is slightly trickier.
358 // We want access to the individual dictionaries used for construction
359
361
363 (
364 eptr->stream(),
365 sampledSet::iNewCapture(mesh_, searchEngine_, capture)
366 );
367
368 PtrList<sampledSet> sampSets(input.size());
369 if (initial && !writeAsProbes_)
370 {
371 writers_.resize(sampSets.size());
372 }
373
374 label seti = 0;
375
376 forAll(input, inputi)
377 {
378 const dictionary& subDict = capture[inputi];
379
380 autoPtr<sampledSet> sampSet = input.release(inputi);
381
382 // if (!sampSet || !sampSet->enabled())
383 // {
384 // continue;
385 // }
386
387 // Define the set
388 sampSets.set(seti, sampSet);
389
390 // Define writer, but do not attached
391 if (initial && !writeAsProbes_)
392 {
393 writers_.set
394 (
395 seti,
396 newWriter(writeFormat_, dict_, subDict)
397 );
398
399 // Use outputDir/TIME/set-name
400 writers_[seti].useTimeDir(true);
401 writers_[seti].verbose(verbose_);
402 }
403 ++seti;
404 }
405
406 sampSets.resize(seti);
407 if (initial && !writeAsProbes_)
408 {
409 writers_.resize(seti);
410 }
411
412 static_cast<PtrList<sampledSet>&>(*this).transfer(sampSets);
413 }
414
415 gatherAllSets();
417 needsCorrect_ = false;
418}
419
420
421// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
422
423Foam::sampledSets::sampledSets
424(
425 const word& name,
426 const Time& runTime,
427 const dictionary& dict
428)
429:
430 functionObjects::fvMeshFunctionObject(name, runTime, dict),
432 dict_(dict),
433 loadFromFiles_(false),
434 verbose_(false),
435 onExecute_(false),
436 needsCorrect_(false),
437 writeAsProbes_(false),
438 outputPath_
439 (
440 time_.globalPath()/functionObject::outputPrefix
441 / name/mesh_.regionName()
442 ),
443 searchEngine_(mesh_),
444 samplePointScheme_(),
445 writeFormat_(),
446 selectedFieldNames_(),
447 writers_(),
448 probeFilePtrs_(),
449 gatheredSets_(),
450 gatheredSorting_(),
451 globalIndices_()
453 outputPath_.clean(); // Remove unneeded ".."
454
455 read(dict);
456}
457
458
459Foam::sampledSets::sampledSets
460(
461 const word& name,
462 const objectRegistry& obr,
463 const dictionary& dict,
464 const bool loadFromFiles
465)
466:
467 functionObjects::fvMeshFunctionObject(name, obr, dict),
469 dict_(dict),
470 loadFromFiles_(loadFromFiles),
471 verbose_(false),
472 onExecute_(false),
473 needsCorrect_(false),
474 writeAsProbes_(false),
475 outputPath_
476 (
477 time_.globalPath()/functionObject::outputPrefix
478 / name/mesh_.regionName()
479 ),
480 searchEngine_(mesh_),
481 samplePointScheme_(),
482 writeFormat_(),
483 selectedFieldNames_(),
484 writers_(),
485 probeFilePtrs_(),
486 gatheredSets_(),
487 gatheredSorting_(),
488 globalIndices_()
489{
490 outputPath_.clean(); // Remove unneeded ".."
492 read(dict);
493}
494
495
496// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
497
498bool Foam::sampledSets::verbose(const bool on) noexcept
500 bool old(verbose_);
501 verbose_ = on;
502 return old;
503}
504
505
507{
508 if (&dict_ != &dict)
509 {
510 // Update local copy of dictionary
511 dict_ = dict;
512 }
513
514 fvMeshFunctionObject::read(dict);
515
517 writers_.clear();
518 fieldSelection_.clear();
519 selectedFieldNames_.clear();
520
521 gatheredSets_.clear();
522 gatheredSorting_.clear();
523 globalIndices_.clear();
524
525 verbose_ = dict.getOrDefault("verbose", false);
526 onExecute_ = dict.getOrDefault("sampleOnExecute", false);
527
528 samplePointScheme_ =
529 dict.getOrDefault<word>("interpolationScheme", "cellPoint");
530
531 const entry* eptr = dict.findEntry("sets", keyType::LITERAL);
532
533 if (eptr)
534 {
535 dict.readEntry("setFormat", writeFormat_);
536 }
537
538 // Hard-coded handling of ensemble 'probes' writer
539 writeAsProbes_ = ("probes" == writeFormat_);
540 if (!writeAsProbes_)
541 {
542 // Close all streams
543 probeFilePtrs_.clear();
544 }
545
546 initDict(dict, true);
547
548 // Have some sets, so sort out which fields are needed and report
549
550 if (this->size())
551 {
552 dict_.readEntry("fields", fieldSelection_);
553 fieldSelection_.uniq();
554
555 // Report
556 if (writeAsProbes_)
557 {
558 Info<< "Sampled set as probes ensemble:" << nl;
559
560 forAll(*this, seti)
561 {
562 const sampledSet& s = (*this)[seti];
563 Info<< " " << s.name();
564 }
565 Info<< nl;
566 }
567 else
568 {
569 Info<< "Sampled set:" << nl;
570
571 forAll(*this, seti)
572 {
573 const sampledSet& s = (*this)[seti];
574
575 Info<< " " << s.name() << " -> "
576 << writers_[seti].type() << nl;
577 }
578 }
579
580 Info<< endl;
581 }
582
583 if (debug && Pstream::master())
584 {
585 Pout<< "sample fields:" << flatOutput(fieldSelection_) << nl
586 << "sample sets:" << nl << '(' << nl;
587
588 for
589 (
590 const sampledSet& s
591 : static_cast<const PtrList<sampledSet>&>(*this)
592 )
593 {
594 Pout<< " " << s << endl;
595 }
596 Pout<< ')' << endl;
597 }
598
599 if (writeAsProbes_)
600 {
601 (void) preCheckFields(ACTION_NONE);
602 }
603
604 // FUTURE:
605 // Ensure all sets and merge information are expired
606 // expire(true);
607
608 return true;
609}
610
611
612bool Foam::sampledSets::performAction(unsigned request)
613{
614 if (empty())
615 {
616 // Nothing to do
617 return true;
618 }
619 else if (needsCorrect_)
620 {
621 searchEngine_.correct();
622 initDict(dict_, false);
623 }
624
625 // FUTURE:
626 // Update sets and store
627 // ...
628
629 // Determine availability of fields.
630 // Count number of fields (only seems to be needed for VTK legacy)
631
632 IOobjectList objects = preCheckFields(request);
633
634 const label nFields = selectedFieldNames_.size();
635
636 if (!nFields)
637 {
638 // Nothing to do
639 return true;
640 }
641
642 // Update writers
643 if (!writeAsProbes_)
644 {
645 forAll(*this, seti)
646 {
647 const coordSet& s = gatheredSets_[seti];
648
649 if ((request & ACTION_WRITE) != 0)
650 {
651 coordSetWriter& writer = writers_[seti];
652
653 if (writer.needsUpdate())
654 {
655 writer.setCoordinates(s);
656 }
657
658 if (writer.buffering())
659 {
661 (
662 outputPath_
663 / word
664 (
665 s.name()
666 + coordSetWriter::suffix(selectedFieldNames_)
667 )
668 );
669 }
670 else
671 {
672 writer.open(outputPath_/s.name());
673 }
674
675 writer.beginTime(mesh_.time());
676 }
677 }
678 }
679
680 // Sample fields
681
682 performAction<VolumeField<scalar>>(objects, request);
683 performAction<VolumeField<vector>>(objects, request);
684 performAction<VolumeField<sphericalTensor>>(objects, request);
685 performAction<VolumeField<symmTensor>>(objects, request);
686 performAction<VolumeField<tensor>>(objects, request);
687
688
689 // Finish this time step
690 if (!writeAsProbes_)
691 {
692 forAll(writers_, seti)
693 {
694 // Write geometry if no fields were written so that we still
695 // can have something to look at
696
697 if ((request & ACTION_WRITE) != 0)
698 {
703
704 writers_[seti].endTime();
705 }
707 }
708
709 return true;
710}
711
712
714{
715 if (onExecute_)
716 {
717 return performAction(ACTION_ALL & ~ACTION_WRITE);
718 }
719
720 return true;
721}
722
725{
726 return performAction(ACTION_ALL);
727}
728
731{
732 needsCorrect_ = true;
733}
734
735
737{
738 if (&mpm.mesh() == &mesh_)
739 {
740 correct();
741 }
742}
743
744
745void Foam::sampledSets::movePoints(const polyMesh& mesh)
746{
747 if (&mesh == &mesh_)
748 {
749 correct();
750 }
751}
752
753
755{
756 if (state != polyMesh::UNCHANGED)
757 {
758 correct();
759 }
760}
761
762
763// ************************************************************************* //
Istream and Ostream manipulators taking arguments.
Various functions to operate on Lists.
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edgesCentres")))
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
A HashTable similar to std::unordered_map.
Definition HashTable.H:124
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable,...
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition IOstream.H:437
bool set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition List.H:469
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
virtual const fileName & name() const override
Read/write access to the name of the stream.
Definition OSstream.H:134
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
void transfer(PtrList< T > &list)
Transfer into this list and annul the argument list.
Definition PtrListI.H:289
void resize_null(const label newLen)
Set the addressed list to the given size, deleting all existing entries. Afterwards the list contains...
Definition PtrListI.H:113
void clear()
Clear the PtrList. Delete allocated entries and set size to zero.
Definition PtrListI.H:98
constexpr PtrList() noexcept
Definition PtrListI.H:29
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
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
label size() const noexcept
Definition UPtrListI.H:106
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
T * release() noexcept
Release ownership and return the pointer.
Definition autoPtrI.H:28
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition autoPtr.H:178
Base class for writing coordSet(s) and tracks with fields.
static autoPtr< coordSetWriter > New(const word &writeFormat)
Return a reference to the selected writer.
static word suffix(const word &fldName, const word &fileExt=word::null)
Name suffix based on fieldName (underscore separator).
static dictionary formatOptions(const dictionary &dict, const word &formatName, const word &entryName="formatOptions")
Same as fileFormats::getFormatOptions.
Holds list of sampling positions.
Definition coordSet.H:52
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const entry * findEntry(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Definition dictionaryI.H:84
A keyword and a list of tokens is an 'entry'.
Definition entry.H:66
virtual const dictionary & dict() const =0
Return dictionary, if entry is a dictionary, otherwise Fatal.
A class for handling file names.
Definition fileName.H:75
static bool clean(std::string &str)
Cleanup filename string, possibly applies other transformations such as changing the path separator e...
Definition fileName.C:192
Abstract base-class for Time/database function objects.
static word outputPrefix
Directory prefix.
const fvMesh & mesh_
Reference to the fvMesh.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
virtual const objectRegistry & obr() const
The region or sub-region registry being used.
const Time & time_
Reference to the time database.
@ LITERAL
String literal.
Definition keyType.H:82
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const polyMesh & mesh() const noexcept
Return polyMesh.
Registry of regIOobjects.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition polyMesh.H:92
PtrList read-construction helper that captures dictionaries used during creation.
Definition sampledSet.H:270
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition sampledSet.H:82
static autoPtr< sampledSet > New(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Return a reference to the selected sampledSet.
Definition sampledSet.C:477
Set of sets to sample.
bool verbose(const bool on) noexcept
Enable/disable verbose output.
void correct()
Correct for mesh changes.
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
virtual void readUpdate(const polyMesh::readUpdateState state)
Update for changes of mesh due to readUpdate.
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
virtual bool execute()
Execute, currently does nothing.
virtual bool write()
Sample and write.
virtual bool read(const dictionary &)
Read the sampledSets.
virtual bool open(const fileName &file, bool parallel=UPstream::parRun())
Open file for writing (creates parent directory).
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
thermo correct()
dynamicFvMesh & mesh
engineTime & runTime
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
OBJstream os(runTime.globalPath()/outputName)
#define doLocalCode(FieldType, Variable)
label nPoints
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))
#define DebugInfo
Report an information message using Foam::Info.
#define WarningInFunction
Report a warning using Foam::Warning.
bool found(const ListType &input, const UnaryPredicate &pred, const label start=0)
Same as found_if.
Definition ListOps.H:886
Namespace for handling debugging switches.
Definition debug.C:45
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
GeometricField< vector, fvPatchField, volMesh > volVectorField
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition POSIX.C:616
GeometricField< scalar, fvPatchField, volMesh > volScalarField
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
GeometricField< tensor, fvPatchField, volMesh > volTensorField
GeometricField< sphericalTensor, fvPatchField, volMesh > volSphericalTensorField
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition FlatOutput.H:217
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
vector point
Point is a vector.
Definition point.H:37
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.
vectorField pointField
pointField is a vectorField.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
const pointField & pts
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition stdFoam.H:214
Dispatch tag: Construct 'one-sided' from local sizes, using gather but no broadcast.