externalCoupledMixedFvPatchField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "fvPatchFieldMapper.H"
28 #include "volFields.H"
29 #include "IFstream.H"
30 #include "globalIndex.H"
31 #include "OSspecific.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 template<class Type>
37 
38 template<class Type>
41 
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
45 template<class Type>
47 {
48  word regionName(this->internalField().mesh().name());
49  if (regionName == polyMesh::defaultRegion)
50  {
51  regionName = ".";
52  }
53 
54  fileName result(commsDir_/regionName);
55  result.clean();
56 
57  return result;
58 }
59 
60 
61 template<class Type>
63 (
64  const labelList& patchIDs
65 )
66 {
67  const VolField<Type>& cvf =
68  static_cast<const VolField<Type>&>(this->internalField());
69 
70  VolField<Type>& vf = const_cast<VolField<Type>&>(cvf);
71 
72  typename VolField<Type>::Boundary& bf = vf.boundaryFieldRef();
73 
74  // number of patches can be different in parallel...
75  label nPatch = bf.size();
76  reduce(nPatch, maxOp<label>());
77 
78  offsets_.setSize(nPatch);
79  forAll(offsets_, i)
80  {
81  offsets_[i].setSize(Pstream::nProcs());
82  offsets_[i] = 0;
83  }
84 
85  // set the master patch
86  forAll(patchIDs, i)
87  {
88  label patchi = patchIDs[i];
89 
90  externalCoupledMixedFvPatchField<Type>& pf =
91  refCast<externalCoupledMixedFvPatchField<Type>>(bf[patchi]);
92 
93  offsets_[patchi][Pstream::myProcNo()] = pf.size();
94 
95  if (i == 0)
96  {
97  pf.master() = true;
98  }
99  else
100  {
101  pf.master() = false;
102  }
103  }
104 
105  // set the patch offsets
106  int tag = Pstream::msgType() + 1;
107  forAll(offsets_, i)
108  {
109  Pstream::gatherList(offsets_[i], tag);
110  Pstream::scatterList(offsets_[i], tag);
111  }
112 
113  label patchOffset = 0;
114  forAll(offsets_, patchi)
115  {
116  label sumOffset = 0;
117  List<label>& procOffsets = offsets_[patchi];
118 
119  forAll(procOffsets, proci)
120  {
121  label o = procOffsets[proci];
122  if (o > 0)
123  {
124  procOffsets[proci] = patchOffset + sumOffset;
125  sumOffset += o;
126  }
127  }
128  patchOffset += sumOffset;
129  }
130 }
131 
132 
133 template<class Type>
135 (
136  OFstream& osPoints,
137  OFstream& osFaces
138 ) const
139 {
140  int tag = Pstream::msgType() + 1;
141 
142  const label proci = Pstream::myProcNo();
143  const polyPatch& p = this->patch().patch();
144  const polyMesh& mesh = p.boundaryMesh().mesh();
145 
146  labelList pointToGlobal;
147  labelList uniquePointIDs;
148  (void)mesh.globalData().mergePoints
149  (
150  p.meshPoints(),
151  p.meshPointMap(),
152  pointToGlobal,
153  uniquePointIDs
154  );
155 
156  List<pointField> allPoints(Pstream::nProcs());
157  allPoints[proci] = pointField(mesh.points(), uniquePointIDs);
158  Pstream::gatherList(allPoints, tag);
159 
160  List<faceList> allFaces(Pstream::nProcs());
161  faceList& patchFaces = allFaces[proci];
162  patchFaces = p.localFaces();
163  forAll(patchFaces, facei)
164  {
165  inplaceRenumber(pointToGlobal, patchFaces[facei]);
166  }
167 
168  Pstream::gatherList(allFaces, tag);
169 
170  if (Pstream::master())
171  {
172  pointField pts
173  (
174  ListListOps::combine<pointField>(allPoints, accessOp<pointField>())
175  );
176 
177  // write points
178  osPoints << patchKey.c_str() << this->patch().name() << pts << endl;
179 
180  faceList fcs
181  (
182  ListListOps::combine<faceList>(allFaces, accessOp<faceList>())
183  );
184 
185  // write faces
186  osFaces<< patchKey.c_str() << this->patch().name() << fcs << endl;
187  }
188 }
189 
190 
191 template<class Type>
193 {
194  return fileName(baseDir()/(lockName + ".lock"));
195 }
196 
197 
198 template<class Type>
200 {
201  if (!master_ || !Pstream::master())
202  {
203  return;
204  }
205 
206  const fileName fName(lockFile());
207  IFstream is(fName);
208 
209  // only create lock file if it doesn't already exist
210  if (!is.good())
211  {
212  if (log_)
213  {
214  Info<< type() << ": creating lock file" << endl;
215  }
216 
217  OFstream os(fName);
218  os << "lock file";
219  os.flush();
220  }
221 }
222 
223 
224 template<class Type>
226 {
227  if (!master_ || !Pstream::master())
228  {
229  return;
230  }
231 
232  if (log_)
233  {
234  Info<< type() << ": removing lock file" << endl;
235  }
236 
237  rm(lockFile());
238 }
239 
240 
241 template<class Type>
243 {
244  // only wait on master patch
245 
246  const VolField<Type>& cvf =
247  static_cast<const VolField<Type>&>(this->internalField());
248 
249  const typename VolField<Type>::Boundary& bf =
250  cvf.boundaryField();
251 
252  forAll(coupledPatchIDs_, i)
253  {
254  label patchi = coupledPatchIDs_[i];
255 
257  refCast<const externalCoupledMixedFvPatchField<Type>>(bf[patchi]);
258 
259  if (pf.master())
260  {
261  pf.wait();
262  break;
263  }
264  }
265 }
266 
267 
268 template<class Type>
270 {
271  const fileName fName(lockFile());
272  label found = 0;
273  label totalTime = 0;
274 
275  if (log_)
276  {
277  Info<< type() << ": beginning wait for lock file " << fName << endl;
278  }
279 
280  while (found == 0)
281  {
282  if (Pstream::master())
283  {
284  if (totalTime > timeOut_)
285  {
287  << "Wait time exceeded time out time of " << timeOut_
288  << " s" << abort(FatalError);
289  }
290 
291  IFstream is(fName);
292 
293  if (is.good())
294  {
295  found++;
296 
297  if (log_)
298  {
299  Info<< type() << ": found lock file " << fName << endl;
300  }
301  }
302  else
303  {
304  sleep(waitInterval_);
305  totalTime += waitInterval_;
306 
307  if (log_)
308  {
309  Info<< type() << ": wait time = " << totalTime << endl;
310  }
311  }
312  }
313 
314  // prevent other procs from racing ahead
315  reduce(found, sumOp<label>());
316  }
317 }
318 
319 
320 template<class Type>
322 (
323  IFstream& is
324 ) const
325 {
326  if (!is.good())
327  {
329  << "Unable to open data transfer file " << is.name()
330  << " for patch " << this->patch().name()
331  << exit(FatalError);
332  }
333 
334  label offset = offsets_[this->patch().index()][Pstream::myProcNo()];
335 
336  // scan forward to start reading at relevant line/position
337  string line;
338  for (label i = 0; i < offset; i++)
339  {
340  if (is.good())
341  {
342  is.getLine(line);
343  }
344  else
345  {
347  << "Unable to scan forward to appropriate read position for "
348  << "data transfer file " << is.name()
349  << " for patch " << this->patch().name()
350  << exit(FatalError);
351  }
352  }
353 }
354 
355 
356 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
357 
358 template<class Type>
360 (
361  const fileName& transferFile
362 )
363 {
364  // read data passed back from external source
365  IFstream is(transferFile + ".in");
366 
367  // pre-process the input transfer file
368  initialiseRead(is);
369 
370  // read data from file
371  forAll(this->patch(), facei)
372  {
373  if (is.good())
374  {
375  is >> this->refValue()[facei]
376  >> this->refGrad()[facei]
377  >> this->valueFraction()[facei];
378  }
379  else
380  {
382  << "Insufficient data for patch "
383  << this->patch().name()
384  << " in file " << is.name()
385  << exit(FatalError);
386  }
387  }
388 
389  initialised_ = true;
390 
391  // update the value from the mixed condition
393 }
394 
395 
396 template<class Type>
398 (
399  const fileName& transferFile
400 ) const
401 {
402  if (!master_)
403  {
404  return;
405  }
406 
407  OFstream os(transferFile);
408 
409  writeHeader(os);
410 
411  const VolField<Type>& cvf =
412  static_cast<const VolField<Type>&>(this->internalField());
413 
414  const typename VolField<Type>::Boundary& bf =
415  cvf.boundaryField();
416 
417  forAll(coupledPatchIDs_, i)
418  {
419  label patchi = coupledPatchIDs_[i];
420 
422  refCast<const externalCoupledMixedFvPatchField<Type>>(bf[patchi]);
423 
424  pf.transferData(os);
425  }
426 }
427 
428 
429 template<class Type>
431 (
432  OFstream& os
433 ) const
434 {
435  os << "# Values: magSf value snGrad" << endl;
436 }
437 
438 
439 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
440 
441 template<class Type>
444 (
445  const fvPatch& p,
447  const dictionary& dict
448 )
449 :
450  mixedFvPatchField<Type>(p, iF, dict, false),
451  commsDir_(dict.lookup("commsDir")),
452  fName_(dict.lookup("file")),
453  waitInterval_(dict.lookupOrDefault("waitInterval", 1)),
454  timeOut_(dict.lookupOrDefault("timeOut", 100*waitInterval_)),
455  calcFrequency_(dict.lookupOrDefault("calcFrequency", 1)),
456  initByExternal_(readBool(dict.lookup("initByExternal"))),
457  log_(dict.lookupOrDefault("log", false)),
458  master_(true),
459  offsets_(),
460  initialised_(false),
461  coupledPatchIDs_()
462 {
463  if (dict.found("value"))
464  {
466  (
467  Field<Type>("value", dict, p.size())
468  );
469  }
470  else
471  {
473  }
474 
475  commsDir_.expand();
476 
477  if (Pstream::master())
478  {
479  mkDir(baseDir());
480  }
481 
482  if (!initByExternal_)
483  {
484  createLockFile();
485  }
486 
487  // initialise as a fixed value
488  this->refValue() = *this;
489  this->refGrad() = Zero;
490  this->valueFraction() = 1.0;
491 }
492 
493 
494 template<class Type>
497 (
499  const fvPatch& p,
501  const fvPatchFieldMapper& mapper
502 )
503 :
504  mixedFvPatchField<Type>(ptf, p, iF, mapper),
505  commsDir_(ptf.commsDir_),
506  fName_(ptf.fName_),
507  waitInterval_(ptf.waitInterval_),
508  timeOut_(ptf.timeOut_),
509  calcFrequency_(ptf.calcFrequency_),
510  initByExternal_(ptf.initByExternal_),
511  log_(ptf.log_),
512  master_(ptf.master_),
513  offsets_(ptf.offsets_),
514  initialised_(ptf.initialised_),
515  coupledPatchIDs_(ptf.coupledPatchIDs_)
516 {}
517 
518 
519 template<class Type>
522 (
525 )
526 :
527  mixedFvPatchField<Type>(ecmpf, iF),
528  commsDir_(ecmpf.commsDir_),
529  fName_(ecmpf.fName_),
530  waitInterval_(ecmpf.waitInterval_),
531  timeOut_(ecmpf.timeOut_),
532  calcFrequency_(ecmpf.calcFrequency_),
533  initByExternal_(ecmpf.initByExternal_),
534  log_(ecmpf.log_),
535  master_(ecmpf.master_),
536  offsets_(ecmpf.offsets_),
537  initialised_(ecmpf.initialised_),
538  coupledPatchIDs_(ecmpf.coupledPatchIDs_)
539 {}
540 
541 
542 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
543 
544 template<class Type>
547 {
548  removeLockFile();
549 }
550 
551 
552 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
553 
554 template<class Type>
556 (
557  const fileName& transferFile
558 )
559 {
560  if (initialised_)
561  {
562  return;
563  }
564 
565  const VolField<Type>& cvf =
566  static_cast<const VolField<Type>&>(this->internalField());
567 
568  VolField<Type>& vf = const_cast<VolField<Type>&>(cvf);
569 
570  typename VolField<Type>::Boundary& bf = vf.boundaryFieldRef();
571 
572  // identify all coupled patches
573  DynamicList<label> coupledPatchIDs(bf.size());
574 
575  forAll(bf, patchi)
576  {
577  if (isA<externalCoupledMixedFvPatchField<Type>>(bf[patchi]))
578  {
579  coupledPatchIDs.append(patchi);
580  }
581  }
582 
583  coupledPatchIDs_.transfer(coupledPatchIDs);
584 
585 
586  // initialise by external solver, or just set the master patch
587  if (initByExternal_)
588  {
589  // remove lock file, signalling external source to execute
590 // removeLockFile();
591 
592  forAll(coupledPatchIDs_, i)
593  {
594  label patchi = coupledPatchIDs_[i];
595 
596  externalCoupledMixedFvPatchField<Type>& pf =
597  refCast<externalCoupledMixedFvPatchField<Type>>(bf[patchi]);
598 
599  pf.setMaster(coupledPatchIDs_);
600  }
601 
602 
603  // wait for initial data to be made available
604  startWait();
605 
606  // read the initial data
607  if (master_)
608  {
609  forAll(coupledPatchIDs_, i)
610  {
611  label patchi = coupledPatchIDs_[i];
612 
613  externalCoupledMixedFvPatchField<Type>& pf =
614  refCast<externalCoupledMixedFvPatchField<Type>>(bf[patchi]);
615 
616  pf.readData(transferFile);
617  }
618  }
619  }
620  else
621  {
622  setMaster(coupledPatchIDs_);
623  }
624 
625  initialised_ = true;
626 }
627 
628 
629 template<class Type>
631 (
632  const Pstream::commsTypes comms
633 )
634 {
635  if (!initialised_ || this->db().time().timeIndex() % calcFrequency_ == 0)
636  {
637  const fileName transferFile(baseDir()/fName_);
638 
639  // initialise the coupling
640  initialise(transferFile);
641 
642  // write data for external source
643  writeData(transferFile + ".out");
644 
645  // remove lock file, signalling external source to execute
646  removeLockFile();
647 
648  // wait for response
649  startWait();
650 
651  if (master_ && Pstream::master())
652  {
653  // remove old data file from OpenFOAM
654  rm(transferFile + ".out");
655  }
656 
657  // read data passed back from external source
658  readData(transferFile);
659 
660  // create lock file for external source
661  createLockFile();
662  }
663 }
664 
665 
666 template<class Type>
668 (
669  OFstream& os
670 ) const
671 {
672  if (log_)
673  {
674  Info<< type() << ": writing data to " << os.name() << endl;
675  }
676 
677  if (Pstream::parRun())
678  {
679  int tag = Pstream::msgType() + 1;
680 
682  magSfs[Pstream::myProcNo()].setSize(this->patch().size());
683  magSfs[Pstream::myProcNo()] = this->patch().magSf();
684  Pstream::gatherList(magSfs, tag);
685 
687  values[Pstream::myProcNo()].setSize(this->patch().size());
688  values[Pstream::myProcNo()] = this->refValue();
689  Pstream::gatherList(values, tag);
690 
691  List<Field<Type>> snGrads(Pstream::nProcs());
692  snGrads[Pstream::myProcNo()].setSize(this->patch().size());
693  snGrads[Pstream::myProcNo()] = this->snGrad();
694  Pstream::gatherList(snGrads, tag);
695 
696  if (Pstream::master())
697  {
698  forAll(values, proci)
699  {
700  const Field<scalar>& magSf = magSfs[proci];
701  const Field<Type>& value = values[proci];
702  const Field<Type>& snGrad = snGrads[proci];
703 
704  forAll(magSf, facei)
705  {
706  os << magSf[facei] << token::SPACE
707  << value[facei] << token::SPACE
708  << snGrad[facei] << nl;
709  }
710  }
711 
712  os.flush();
713  }
714  }
715  else
716  {
717  const Field<scalar>& magSf(this->patch().magSf());
718  const Field<Type>& value(this->refValue());
719  const Field<Type> snGrad(this->snGrad());
720 
721  forAll(magSf, facei)
722  {
723  os << magSf[facei] << token::SPACE
724  << value[facei] << token::SPACE
725  << snGrad[facei] << nl;
726  }
727 
728  os.flush();
729  }
730 }
731 
732 
733 template<class Type>
735 {
736  const VolField<Type>& cvf =
737  static_cast<const VolField<Type>&>(this->internalField());
738 
739  const typename VolField<Type>::Boundary& bf =
740  cvf.boundaryField();
741 
742  OFstream osPoints(baseDir()/"patchPoints");
743  OFstream osFaces(baseDir()/"patchFaces");
744 
745  if (log_)
746  {
747  Info<< "writing collated patch points to: " << osPoints.name() << nl
748  << "writing collated patch faces to: " << osFaces.name() << endl;
749  }
750 
751  forAll(bf, patchi)
752  {
754  {
756  refCast<const externalCoupledMixedFvPatchField<Type>>
757  (
758  bf[patchi]
759  );
760 
761  pf.writeGeometry(osPoints, osFaces);
762  }
763  }
764 }
765 
766 
767 template<class Type>
769 {
771 
772  writeEntry(os, "commsDir", commsDir_);
773  writeEntry(os, "file", fName_);
774  writeEntry(os, "waitInterval", waitInterval_);
775  writeEntry(os, "timeOut", timeOut_);
776  writeEntry(os, "calcFrequency", calcFrequency_);
777  writeEntry(os, "initByExternal", initByExternal_);
778  writeEntry(os, "log", log_);
779 
780  writeEntry(os, "value", *this);
781 }
782 
783 
784 // ************************************************************************* //
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Pre-declare SubField and related Field type.
Definition: Field.H:82
Generic GeometricBoundaryField class.
Generic GeometricField class.
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Input from file stream.
Definition: IFstream.H:85
const fileName & name() const
Return the name of the stream.
Definition: IFstream.H:116
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:330
ISstream & getLine(string &, const bool continuation=true)
Read line into a string.
Definition: ISstream.C:692
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void setSize(const label)
Reset size of List.
Definition: List.C:281
Output to file stream.
Definition: OFstream.H:86
const fileName & name() const
Return the name of the stream.
Definition: OFstream.H:120
virtual void flush()
Flush stream.
Definition: OSstream.C:207
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
commsTypes
Types of communications.
Definition: UPstream.H:65
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:476
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
This boundary condition provides an interface to an external application. Values are transferred as p...
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
virtual void writeHeader(OFstream &os) const
Write header to transfer file.
externalCoupledMixedFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void writeData(const fileName &transferFile) const
Write data for external source - calls transferData.
void writeGeometry() const
Write the geometry to the comms dir.
virtual void transferData(OFstream &os) const
Transfer data for external source.
virtual void readData(const fileName &transferFile)
Read data from external source.
A class for handling file names.
Definition: fileName.H:82
word name() const
Return file name (part beyond last /)
Definition: fileName.C:195
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:172
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
A line primitive.
Definition: line.H:71
This boundary condition provides a base class for 'mixed' type boundary conditions,...
virtual void write(Ostream &) const
Write.
virtual Field< Type > & refValue()
virtual scalarField & valueFraction()
virtual Field< Type > & refGrad()
A class for handling character strings derived from std::string.
Definition: string.H:79
string & expand(const bool allowEmpty=false)
Expand initial tildes and all occurrences of environment variables.
Definition: string.C:125
A class for handling words, derived from string.
Definition: word.H:62
Foam::word regionName
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
label patchi
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
void writeHeader(std::ostream &, const bool isBinary, const std::string &title)
Write header.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
bool readBool(Istream &)
Definition: boolIO.C:60
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
bool rm(const fileName &)
Remove a file, returning true if successful otherwise false.
Definition: POSIX.C:1017
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:290
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
errorManip< error > abort(error &err)
Definition: errorManip.H:131
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:139
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
error FatalError
void offset(label &lst, const label o)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
List< face > faceList
Definition: faceListFwd.H:43
unsigned int sleep(const unsigned int)
Sleep for the specified number of seconds.
Definition: POSIX.C:1130
static const char nl
Definition: Ostream.H:260
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
label timeIndex
Definition: getTimeIndex.H:4
dictionary dict
volScalarField & p
const bool writeData(readBool(pdfDictionary.lookup("writeData")))