Loading...
Searching...
No Matches
externalCoupled.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) 2015-2024 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 "externalCoupled.H"
29#include "stringListOps.H"
31#include "OSspecific.H"
32#include "Fstream.H"
33#include "SpanStream.H"
34#include "volFields.H"
35#include "globalIndex.H"
36#include "fvMesh.H"
37#include "DynamicField.H"
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41namespace Foam
42{
43namespace functionObjects
44{
47}
48}
49
51
52
53namespace Foam
54{
56//- Write list content with size, bracket, content, bracket one-per-line.
57// This makes for consistent for parsing, regardless of the list length.
58template <class T>
59static void writeList(Ostream& os, const string& header, const UList<T>& L)
60{
61 // Header string
62 os << header.c_str() << nl;
63
64 // Write size and start delimiter
65 os << L.size() << nl
67
68 // Write contents
69 forAll(L, i)
70 {
71 os << L[i] << nl;
72 }
73
74 // Write end delimiter
75 os << token::END_LIST << nl << endl;
76}
78
79} // End namespace Foam
80
81
82// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
83
84Foam::fileName Foam::functionObjects::externalCoupled::groupDir
85(
86 const fileName& commsDir,
87 const word& regionGroupName,
88 const wordRe& groupName
89)
90{
91 fileName result
92 (
93 commsDir
94 / regionGroupName
95 / word::validate(groupName)
96 );
97 result.clean(); // Remove unneeded ".."
98
99 return result;
100}
101
102
103void Foam::functionObjects::externalCoupled::readColumns
104(
105 const label nRows,
106 const label nColumns,
107 autoPtr<IFstream>& masterFilePtr,
108 List<scalarField>& data
109) const
110{
111 // Get sizes for all processors
112 const globalIndex globalFaces(globalIndex::gatherOnly{}, nRows);
113
114 PstreamBuffers pBufs;
115
116 if (UPstream::master())
117 {
118 // Read data from file and send to destination processor
119 auto& ifile = masterFilePtr();
120 string line;
121
122 for (const int proci : UPstream::allProcs())
123 {
124 // Temporary storage
125 List<scalarField> values(nColumns);
126
127 // Number of rows to read for processor proci
128 const label procNRows = globalFaces.localSize(proci);
129
130 forAll(values, columni)
131 {
132 values[columni].setSize(procNRows);
133 }
134
135 for (label rowi = 0; rowi < procNRows; ++rowi)
136 {
137 // Get a line
138 do
139 {
140 if (!ifile.good())
141 {
143 << "Trying to read data for processor " << proci
144 << " row " << rowi
145 << ". Does your file have as many rows as there are"
146 << " patch faces (" << globalFaces.totalSize()
147 << ") ?" << exit(FatalIOError);
148 }
149
150 ifile.getLine(line);
151 }
152 while (line.empty() || line[0] == '#');
153
154 ISpanStream isstr(line);
155
156 for (label columni = 0; columni < nColumns; ++columni)
157 {
158 isstr >> values[columni][rowi];
159 }
160 }
161
162 // Send to proci
163 UOPstream toProc(proci, pBufs);
164 toProc << values;
165 }
166 }
167 pBufs.finishedScatters();
168
169 // Get scattered data
170 UIPstream fromMaster(UPstream::masterNo(), pBufs);
171 fromMaster >> data;
172}
173
174
175void Foam::functionObjects::externalCoupled::readLines
176(
177 const label nRows,
178 autoPtr<IFstream>& masterFilePtr,
179 std::string& lines
180) const
181{
182 // Get sizes for all processors
183 const globalIndex globalFaces(globalIndex::gatherOnly{}, nRows);
184
185 PstreamBuffers pBufs;
186
187 if (UPstream::master())
188 {
189 // Read line from file and send to destination processor
190 auto& ifile = masterFilePtr();
191 string line;
192
193 for (const int proci : UPstream::allProcs())
194 {
195 // Number of rows to read for processor proci
196 const label procNRows = globalFaces.localSize(proci);
197
198 UOPstream toProc(proci, pBufs);
199
200 for (label rowi = 0; rowi < procNRows; ++rowi)
201 {
202 // Get a line
203 do
204 {
205 if (!ifile.good())
206 {
208 << "Trying to read data for processor " << proci
209 << " row " << rowi
210 << ". Does your file have as many rows as there are"
211 << " patch faces (" << globalFaces.totalSize()
212 << ") ?" << exit(FatalIOError);
213 }
214
215 ifile.getLine(line);
216 }
217 while (line.empty() || line[0] == '#');
218
219 // Send line (without newline) to the destination processor
220 toProc << line;
221 }
222 }
223 }
224
225 pBufs.finishedScatters();
226
227 // Sizing is approximate (slightly too large)
228 lines.reserve(pBufs.recvDataCount(UPstream::masterNo()));
229
230 // Get scattered data
231 UIPstream fromMaster(UPstream::masterNo(), pBufs);
232 for (label rowi = 0; rowi < nRows; ++rowi)
233 {
234 string line(fromMaster);
235 lines += line;
236 lines += '\n';
238}
239
240
242(
244 const fileName& commsDir,
245 const wordRe& groupName
246)
247{
249 forAll(meshes, i)
250 {
251 regionNames[i] = meshes[i].dbDir();
252 }
253
254 // Make sure meshes are provided in sorted order
255 checkOrder(regionNames);
256
257 fileName dir(groupDir(commsDir, compositeName(regionNames), groupName));
258
259 autoPtr<OFstream> osPointsPtr;
260 autoPtr<OFstream> osFacesPtr;
261 if (UPstream::master())
262 {
263 Foam::mkDir(dir);
264 osPointsPtr.reset(new OFstream(dir/"patchPoints"));
265 osFacesPtr.reset(new OFstream(dir/"patchFaces"));
266
267 osPointsPtr() << "// Group: " << groupName << endl;
268 osFacesPtr() << "// Group: " << groupName << endl;
269
270 Info<< typeName << ": writing geometry to " << dir << endl;
271 }
272
273 // Individual region/patch entries
274
275 DynamicList<face> allFaces;
276 DynamicField<point> allPoints;
277
278 labelList pointToGlobal;
279 labelList uniquePointIDs;
280 for (const fvMesh& mesh : meshes)
281 {
282 const labelList patchIDs
283 (
284 mesh.boundaryMesh().patchSet
285 (
286 wordRes(Foam::one{}, groupName)
287 ).sortedToc()
288 );
289
290 for (const label patchi : patchIDs)
291 {
292 const polyPatch& p = mesh.boundaryMesh()[patchi];
293
294 mesh.globalData().mergePoints
295 (
296 p.meshPoints(),
297 p.meshPointMap(),
298 pointToGlobal,
299 uniquePointIDs
300 );
301
302 label proci = Pstream::myProcNo();
303
304 List<pointField> collectedPoints(Pstream::nProcs());
305 collectedPoints[proci] = pointField(mesh.points(), uniquePointIDs);
306 Pstream::gatherList(collectedPoints);
307
308 List<faceList> collectedFaces(Pstream::nProcs());
309 faceList& patchFaces = collectedFaces[proci];
310 patchFaces = p.localFaces();
311 for (auto& f : patchFaces)
312 {
313 inplaceRenumber(pointToGlobal, f);
314 }
315 Pstream::gatherList(collectedFaces);
316
317 if (UPstream::master())
318 {
319 allPoints.clear();
320 allFaces.clear();
321
322 for (const int proci : Pstream::allProcs())
323 {
324 allPoints.append(collectedPoints[proci]);
325 allFaces.append(collectedFaces[proci]);
326 }
327
328 Info<< typeName << ": mesh " << mesh.name()
329 << ", patch " << p.name()
330 << ": writing " << allPoints.size() << " points to "
331 << osPointsPtr().name() << nl
332 << typeName << ": mesh " << mesh.name()
333 << ", patch " << p.name()
334 << ": writing " << allFaces.size() << " faces to "
335 << osFacesPtr().name() << endl;
336
337 // The entry name (region / patch)
338 const string entryHeader =
339 patchKey + ' ' + mesh.name() + ' ' + p.name();
340
341 writeList(osPointsPtr(), entryHeader, allPoints);
342 writeList(osFacesPtr(), entryHeader, allFaces);
343 }
344 }
346}
347
348
350(
351 const wordList& regionNames
352)
353{
354 if (regionNames.size() == 0)
355 {
357 << "Empty regionNames" << abort(FatalError);
358 return word::null;
359 }
360 else if (regionNames.size() == 1)
361 {
362 // For compatibility with single region cases
363 // - suppress single region name
365 }
366
367 // Enforce lexical ordering
368 checkOrder(regionNames);
369
370 word composite(regionNames[0]);
371 for (label i = 1; i < regionNames.size(); ++i)
372 {
373 composite += "_" + regionNames[i];
374 }
375
376 return composite;
377}
378
379
380void Foam::functionObjects::externalCoupled::checkOrder
381(
382 const wordList& regionNames
383)
384{
386 if (order != identity(regionNames.size()))
387 {
389 << "regionNames " << regionNames << " not in alphabetical order :"
390 << order << exit(FatalError);
391 }
392}
393
394
395void Foam::functionObjects::externalCoupled::initCoupling()
396{
397 if (initialisedCoupling_)
398 {
399 return;
400 }
401
402 // Write the geometry if not already there
403 forAll(regionGroupNames_, regioni)
404 {
405 const word& compName = regionGroupNames_[regioni];
406 const wordList& regionNames = regionGroupRegions_[regioni];
407
408 // Get the meshes for the region-group
409 UPtrList<const fvMesh> meshes(regionNames.size());
410 forAll(regionNames, regi)
411 {
412 meshes.set(regi, time_.findObject<fvMesh>(regionNames[regi]));
413 }
414
415 const labelList& groups = regionToGroups_[compName];
416
417 for (const label groupi : groups)
418 {
419 const wordRe& groupName = groupNames_[groupi];
420
421 bool geomExists = false;
422 if (Pstream::master())
423 {
424 fileName dir(groupDir(commDirectory(), compName, groupName));
425
426 geomExists =
427 isFile(dir/"patchPoints")
428 || isFile(dir/"patchFaces");
429 }
430
431 Pstream::broadcast(geomExists);
432
433 if (!geomExists)
434 {
435 writeGeometry(meshes, commDirectory(), groupName);
436 }
437 }
438 }
439
440 if (slaveFirst())
441 {
442 // Wait for initial data to be made available
443 waitForSlave();
444
445 // Read data passed back from external source
446 readDataMaster();
447 }
448
449 initialisedCoupling_ = true;
450}
451
452
453void Foam::functionObjects::externalCoupled::performCoupling()
454{
455 // Ensure coupling has been initialised
456 initCoupling();
457
458 // Write data for external source
459 writeDataMaster();
460
461 // Signal external source to execute (by removing lock file)
462 // - Wait for slave to provide data
463 useSlave();
464
465 // Wait for response - and catch any abort information sent from slave
466 const auto action = waitForSlave();
467
468 // Remove old data files from OpenFOAM
469 removeDataMaster();
470
471 // Read data passed back from external source
472 readDataMaster();
473
474 // Signal external source to wait (by creating the lock file)
475 useMaster();
476
477 // Update information about last triggering
478 lastTrigger_ = time_.timeIndex();
479
480 // Process any abort information sent from slave
481 if
482 (
483 action != time_.stopAt()
485 )
486 {
487 Info<< type() << ": slave requested action "
488 << Time::stopAtControlNames[action] << endl;
489
490 time_.stopAt(action);
491 }
492}
494
495// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
496
498(
499 const word& name,
500 const Time& runTime,
501 const dictionary& dict
502)
503:
505 externalFileCoupler(),
506 calcFrequency_(-1),
507 lastTrigger_(-1),
508 initialisedCoupling_(false)
509{
510 read(dict);
511
512 if (!slaveFirst())
513 {
514 useMaster();
515 }
516}
518
519// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
520
522{
523 // Not initialized or overdue
524 if
525 (
526 !initialisedCoupling_
527 || (time_.timeIndex() >= lastTrigger_ + calcFrequency_)
528 )
529 {
530 performCoupling();
531 }
532
533 return false;
534}
535
536
538{
539 performCoupling();
540
541 return true;
542}
543
544
546{
548
549 // Remove old data files
550 removeDataMaster();
551 removeDataSlave();
552 shutdown();
553
554 return true;
555}
556
557
559{
562
563 calcFrequency_ =
564 dict.getCheckOrDefault("calcFrequency", 1, labelMinMax::ge(1));
565
566 // Leave trigger intact
567
568 // Get names of all fvMeshes (and derived types)
569 wordList allRegionNames(time_.sortedNames<fvMesh>());
570
571 const dictionary& allRegionsDict = dict.subDict("regions");
572 for (const entry& dEntry : allRegionsDict)
573 {
574 if (!dEntry.isDict())
575 {
576 FatalIOErrorInFunction(allRegionsDict)
577 << "Regions must be specified in dictionary format"
578 << exit(FatalIOError);
579 }
580
581 const wordRe regionGroupName(dEntry.keyword());
582 const dictionary& regionDict = dEntry.dict();
583
584 labelList regionIDs
585 (
586 wordRes::matching(regionGroupName, allRegionNames)
587 );
588
589 const wordList regionNames(allRegionNames, regionIDs);
590
591 regionGroupNames_.append(compositeName(regionNames));
592 regionGroupRegions_.append(regionNames);
593
594 for (const entry& dEntry : regionDict)
595 {
596 if (!dEntry.isDict())
597 {
598 FatalIOErrorInFunction(regionDict)
599 << "Regions must be specified in dictionary format"
600 << exit(FatalIOError);
601 }
602
603 const wordRe groupName(dEntry.keyword());
604 const dictionary& groupDict = dEntry.dict();
605
606 const label nGroups = groupNames_.size();
607 const wordList readFields(groupDict.get<wordList>("readFields"));
608 const wordList writeFields(groupDict.get<wordList>("writeFields"));
609
610 auto fnd = regionToGroups_.find(regionGroupNames_.last());
611 if (fnd.good())
612 {
613 fnd().append(nGroups);
614 }
615 else
616 {
617 regionToGroups_.insert
618 (
619 regionGroupNames_.last(),
620 labelList(one{}, nGroups)
621 );
622 }
623 groupNames_.append(groupName);
624 groupReadFields_.append(readFields);
625 groupWriteFields_.append(writeFields);
626 }
627 }
628
629
630 Info<< type() << ": Communicating with regions:" << endl;
631 for (const word& compName : regionGroupNames_)
632 {
633 Info<< "Region: " << compName << nl << incrIndent;
634 const labelList& groups = regionToGroups_[compName];
635 for (const label groupi : groups)
636 {
637 const wordRe& groupName = groupNames_[groupi];
638
639 Info<< indent << "patchGroup: " << groupName << "\t"
640 << nl
641 << incrIndent
642 << indent << "Reading fields: "
643 << groupReadFields_[groupi]
644 << nl
645 << indent << "Writing fields: "
646 << groupWriteFields_[groupi]
647 << nl
648 << decrIndent;
649 }
651 }
652 Info<< endl;
653
654
655 // Note: we should not have to make directories since the geometry
656 // should already be written - but just make sure
657 if (Pstream::master())
658 {
659 for (const word& compName : regionGroupNames_)
660 {
661 const labelList& groups = regionToGroups_[compName];
662 for (const label groupi : groups)
663 {
664 const wordRe& groupName = groupNames_[groupi];
665
666 fileName dir(groupDir(commDirectory(), compName, groupName));
667
668 if (!isDir(dir))
669 {
670 Log << type() << ": creating communications directory "
671 << dir << endl;
672 mkDir(dir);
673 }
674 }
675 }
676 }
677
678 return true;
679}
680
681
683{
684 forAll(regionGroupNames_, regioni)
685 {
686 const word& compName = regionGroupNames_[regioni];
687 const wordList& regionNames = regionGroupRegions_[regioni];
688
689 // Get the meshes for the region-group
691 forAll(regionNames, regi)
692 {
693 meshes.set(regi, time_.findObject<fvMesh>(regionNames[regi]));
694 }
695
696 const labelList& groups = regionToGroups_[compName];
697
698 for (const label groupi : groups)
699 {
700 const wordRe& groupName = groupNames_[groupi];
701 const wordList& fieldNames = groupReadFields_[groupi];
702
703 for (const word& fieldName : fieldNames)
704 {
705 const bool ok =
706 (
707 readData<scalar>(meshes, groupName, fieldName)
708 || readData<vector>(meshes, groupName, fieldName)
709 || readData<sphericalTensor>(meshes, groupName, fieldName)
710 || readData<symmTensor>(meshes, groupName, fieldName)
711 || readData<tensor>(meshes, groupName, fieldName)
712 );
713
714 if (!ok)
715 {
717 << "Field " << fieldName << " in regions " << compName
718 << " was not found." << endl;
719 }
720 }
721 }
723}
724
725
727{
728 forAll(regionGroupNames_, regioni)
729 {
730 const word& compName = regionGroupNames_[regioni];
731 const wordList& regionNames = regionGroupRegions_[regioni];
732
733 // Get the meshes for the region-group
735 forAll(regionNames, regi)
736 {
737 meshes.set(regi, time_.findObject<fvMesh>(regionNames[regi]));
738 }
739
740 const labelList& groups = regionToGroups_[compName];
741
742 for (const label groupi : groups)
743 {
744 const wordRe& groupName = groupNames_[groupi];
745 const wordList& fieldNames = groupWriteFields_[groupi];
746
747 for (const word& fieldName : fieldNames)
748 {
749 const bool ok =
750 (
751 writeData<scalar>(meshes, groupName, fieldName)
752 || writeData<vector>(meshes, groupName, fieldName)
753 || writeData<sphericalTensor>(meshes, groupName, fieldName)
754 || writeData<symmTensor>(meshes, groupName, fieldName)
755 || writeData<tensor>(meshes, groupName, fieldName)
756 );
757
758 if (!ok)
759 {
761 << "Field " << fieldName << " in regions " << compName
762 << " was not found." << endl;
763 }
764 }
765 }
767}
768
769
771{
772 if (!Pstream::master())
773 {
774 return;
775 }
776
777 Log << type() << ": removing data files written by master" << nl;
778
779 for (const word& compName : regionGroupNames_)
780 {
781 const labelList& groups = regionToGroups_[compName];
782 for (const label groupi : groups)
783 {
784 const wordRe& groupName = groupNames_[groupi];
785 const wordList& fieldNames = groupReadFields_[groupi];
786
787 for (const word& fieldName : fieldNames)
788 {
790 (
791 groupDir(commDirectory(), compName, groupName)
792 / fieldName + ".out"
793 );
794 }
795 }
797}
798
799
801{
802 if (!Pstream::master())
803 {
804 return;
805 }
806
807 Log << type() << ": removing data files written by slave" << nl;
808
809 for (const word& compName : regionGroupNames_)
810 {
811 const labelList& groups = regionToGroups_[compName];
812 for (const label groupi : groups)
813 {
814 const wordRe& groupName = groupNames_[groupi];
815 const wordList& fieldNames = groupReadFields_[groupi];
816
817 for (const word& fieldName : fieldNames)
818 {
820 (
821 groupDir(commDirectory(), compName, groupName)
822 / fieldName + ".in"
823 );
824 }
825 }
827}
828
829
831{
832 return true;
833}
834
835
836// ************************************************************************* //
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
#define Log
Definition PDRblock.C:28
Input/output streams with (internal or external) character storage.
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
writer writeGeometry()
labelList patchIDs
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
Dynamically sized Field. Similar to DynamicList, but inheriting from a Field instead of a List.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
void append(const T &val)
Copy append an element to the end of this list.
virtual const fileName & name() const
The name of the stream.
Definition IOstream.C:33
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
static MinMax< label > ge(const label &minVal)
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
static const Enum< stopAtControls > stopAtControlNames
Names for stopAtControls.
Definition Time.H:113
@ saUnknown
Dummy no-op. Do not change current value.
Definition Time.H:102
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
static constexpr int masterNo() noexcept
Relative rank for the master process - is always 0.
Definition UPstream.H:1691
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
Definition UPstream.H:1697
static rangeType allProcs(const label communicator=worldComm)
Range of process indices for all processes.
Definition UPstream.H:1857
@ gatherList
gatherList [manual algorithm]
Definition UPstream.H:194
@ broadcast
broadcast [MPI]
Definition UPstream.H:189
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition UPtrList.H:101
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition autoPtrI.H:37
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
A keyword and a list of tokens is an 'entry'.
Definition entry.H:66
const fileName & commDirectory() const
Return the file path to the base communications directory.
enum Time::stopAtControls useMaster(const bool wait=false) const
Create lock file to indicate that OpenFOAM is in charge.
bool slaveFirst() const
External application provides initial values.
bool readDict(const dictionary &dict)
Read communication settings from dictionary.
void shutdown() const
Generate status=done in lock (only when run-state = master).
A class for handling file names.
Definition fileName.H:75
static std::string name(const std::string &str)
Return basename (part beyond last /), including its extension.
Definition fileNameI.H:192
Abstract base-class for Time/database function objects.
const word & name() const noexcept
Return the name of this functionObject.
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
virtual bool end()
Called when Time::run() determines that the time-loop exits.
Provides a simple file-based communication interface for explicit coupling with an external applicati...
static void writeGeometry(const UPtrList< const fvMesh > &meshes, const fileName &commsDir, const wordRe &groupName)
Write geometry for the group as region/patch.
virtual void writeDataMaster() const
Write data files (all regions, all fields) from master (OpenFOAM).
virtual void removeDataSlave() const
Remove data files written by slave (external code).
static string patchKey
Name of patch key, e.g. '// Patch:' when looking for start of patch data.
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
virtual void readDataMaster()
Read data files (all regions, all fields) on master (OpenFOAM).
virtual void removeDataMaster() const
Remove data files written by master (OpenFOAM).
externalCoupled(const word &name, const Time &runTime, const dictionary &dict)
Construct given time and dictionary.
virtual bool execute()
Called at each ++ or += of the time-loop.
virtual bool write()
Write, currently a no-op.
virtual bool end()
Called when Time::run() determines that the time-loop exits.
static word compositeName(const wordList &)
Create single name by appending words (in sorted order), separated by '_'.
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition readFields.H:146
const Time & time_
Reference to the time database.
timeFunctionObject(const timeFunctionObject &)=delete
No copy construct.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition one.H:57
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition polyMesh.C:832
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
A class for handling character strings derived from std::string.
Definition string.H:76
@ BEGIN_LIST
Begin list [isseparator].
Definition token.H:174
@ END_LIST
End list [isseparator].
Definition token.H:175
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
Definition wordRe.H:81
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
static labelList matching(const wordRe &select, const UList< StringType > &input, const bool invert=false)
Determine the list indices for all matches.
A class for handling words, derived from Foam::string.
Definition word.H:66
static word validate(const std::string &s, const bool prefix=false)
Construct validated word (no invalid characters).
Definition word.C:39
static const word null
An empty word.
Definition word.H:84
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
dynamicFvMesh & mesh
engineTime & runTime
Foam::PtrList< Foam::fvMesh > meshes(regionNames.size())
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
wordList regionNames
#define WarningInFunction
Report a warning using Foam::Warning.
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition HashOps.H:164
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
void writeList(vtk::formatter &fmt, const UList< uint8_t > &values)
Write a list of uint8_t values.
Namespace for OpenFOAM.
bool rm(const fileName &file)
Remove a file (or its gz equivalent), returning true if successful.
Definition POSIX.C:1406
List< word > wordList
List of word.
Definition fileName.H:60
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values within a list.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
List< label > labelList
A List of labels.
Definition List.H:62
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
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
Definition faceListFwd.H:41
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition Ostream.H:490
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Ostream & indent(Ostream &os)
Indent stream.
Definition Ostream.H:481
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
errorManip< error > abort(error &err)
Definition errorManip.H:139
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
void writeFields(const fvMesh &mesh, const wordHashSet &selectedFields, const bool writeFaceFields)
bool isFile(const fileName &name, const bool checkGzip=true, const bool followLink=true)
Does the name exist as a FILE in the file system?
Definition POSIX.C:879
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
vectorField pointField
pointField is a vectorField.
static void writeData(Ostream &os, const Type &val)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition Ostream.H:499
bool isDir(const fileName &name, const bool followLink=true)
Does the name exist as a DIRECTORY in the file system?
Definition POSIX.C:862
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Operations on lists of strings.
const vector L(dict.get< vector >("L"))