40 #define checkField(fld1, fld2, op) \ 41 if (&(fld1).mesh() != &(fld2).mesh()) \ 43 FatalErrorInFunction \ 44 << "Different mesh for fields " \ 45 << (fld1).name() << " and " << (fld2).name() \ 46 << " during operation " << op \ 47 << abort(FatalError); \ 53 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
56 const dictionary&
dict 59 Internal::readField(
dict,
"internalField");
61 boundaryField_.readField(*
this,
dict.subDict(
"boundaryField"));
65 if (
dict.readIfPresent(
"referenceLevel", refLevel))
67 Field<Type>::operator+=(refLevel);
69 forAll(boundaryField_, patchi)
71 boundaryField_[patchi] == boundaryField_[patchi] + refLevel;
77 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
80 const localIOdictionary
dict 101 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
104 if (this->isReadRequired())
107 <<
"read option IOobject::MUST_READ or MUST_READ_IF_MODIFIED" 108 <<
" suggests that a read constructor for field " << this->
name()
109 <<
" would be more appropriate." <<
endl;
113 this->isReadOptional()
114 && this->
template typeHeaderOk<GeometricField<Type, PatchField, GeoMesh>>
123 if (this->size() != GeoMesh::size(this->
mesh()))
126 <<
" number of field elements = " << this->size()
127 <<
" number of mesh elements = " 128 << GeoMesh::size(this->
mesh())
132 readOldTimeIfPresent();
141 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
150 IOobject::READ_IF_PRESENT,
151 IOobject::AUTO_WRITE,
152 this->registerObject()
157 field0.template typeHeaderOk<GeometricField<Type, PatchField, GeoMesh>>
164 <<
"Reading old time level for field" <<
nl << this->info() <<
endl;
166 field0Ptr_ =
new GeometricField<Type, PatchField, GeoMesh>
175 field0Ptr_->oriented() = this->oriented();
177 field0Ptr_->timeIndex_ = timeIndex_ - 1;
179 if (!field0Ptr_->readOldTimeIfPresent())
181 field0Ptr_->oldTime();
193 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
199 const word& patchFieldType
202 Internal(
io,
mesh, dims, false),
205 fieldPrevIterPtr_(nullptr),
209 <<
"Creating" <<
nl << this->info() <<
endl;
215 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
225 Internal(
io,
mesh, dims, false),
228 fieldPrevIterPtr_(nullptr),
229 boundaryField_(
mesh.
boundary(), *this, patchFieldTypes, actualPatchTypes)
232 <<
"Creating" <<
nl << this->info() <<
endl;
239 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
246 const word& patchFieldType
249 Internal(
io,
mesh, value, dims, false),
252 fieldPrevIterPtr_(nullptr),
256 <<
"Creating" <<
nl << this->info() <<
endl;
258 boundaryField_ == value;
264 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
275 Internal(
io,
mesh, value, dims, false),
278 fieldPrevIterPtr_(nullptr),
279 boundaryField_(
mesh.
boundary(), *this, patchFieldTypes, actualPatchTypes)
282 <<
"Creating" <<
nl << this->info() <<
endl;
284 boundaryField_ == value;
290 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
296 const word& patchFieldType
310 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
332 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
336 const Internal& diField,
337 const PtrList<PatchField<Type>>& ptfl
340 Internal(
io, diField),
343 fieldPrevIterPtr_(nullptr),
347 <<
"Copy construct from components" <<
nl << this->info() <<
endl;
353 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
358 const PtrList<PatchField<Type>>& ptfl
361 Internal(
io,
std::move(diField)),
364 fieldPrevIterPtr_(nullptr),
368 <<
"Move construct from components" <<
nl << this->info() <<
endl;
374 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
379 const PtrList<PatchField<Type>>& ptfl
382 Internal(
io, tfield),
385 fieldPrevIterPtr_(nullptr),
389 <<
"Construct from tmp internalField" <<
nl << this->info() <<
endl;
395 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
398 const Internal& diField,
399 const PtrList<PatchField<Type>>& ptfl
405 fieldPrevIterPtr_(nullptr),
409 <<
"Copy construct from components" <<
nl << this->info() <<
endl;
415 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
419 const PtrList<PatchField<Type>>& ptfl
422 Internal(
std::move(diField)),
425 fieldPrevIterPtr_(nullptr),
429 <<
"Move construct from components" <<
nl << this->info() <<
endl;
435 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
442 const word& patchFieldType
445 Internal(
io,
mesh, dims, iField),
448 fieldPrevIterPtr_(nullptr),
452 <<
"Copy construct from internal field" <<
nl << this->info() <<
endl;
458 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
465 const word& patchFieldType
468 Internal(
io,
mesh, dims,
std::move(iField)),
471 fieldPrevIterPtr_(nullptr),
475 <<
"Move construct from internal field" <<
nl << this->info() <<
endl;
481 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
488 const PtrList<PatchField<Type>>& ptfl
491 Internal(
io,
mesh, dims, iField),
494 fieldPrevIterPtr_(nullptr),
498 <<
"Copy construct from components" <<
nl << this->info() <<
endl;
504 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
511 const PtrList<PatchField<Type>>& ptfl
514 Internal(
io,
mesh, dims,
std::move(iField)),
517 fieldPrevIterPtr_(nullptr),
521 <<
"Move construct from components" <<
nl << this->info() <<
endl;
527 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
534 const PtrList<PatchField<Type>>& ptfl
537 Internal(
io,
mesh, dims, tfield),
540 fieldPrevIterPtr_(nullptr),
544 <<
"Construct from tmp internalField" <<
nl << this->info() <<
endl;
550 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
555 const bool readOldTime
561 fieldPrevIterPtr_(nullptr),
568 if (this->size() != GeoMesh::size(this->
mesh()))
571 <<
" number of field elements = " << this->size()
572 <<
" number of mesh elements = " << GeoMesh::size(this->
mesh())
578 readOldTimeIfPresent();
582 <<
"Finishing read-construction" <<
nl << this->info() <<
endl;
586 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
597 fieldPrevIterPtr_(nullptr),
604 if (this->size() != GeoMesh::size(this->
mesh()))
607 <<
" number of field elements = " << this->size()
608 <<
" number of mesh elements = " << GeoMesh::size(this->
mesh())
613 <<
"Finishing dictionary-construct" <<
nl << this->info() <<
endl;
617 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
626 fieldPrevIterPtr_(nullptr),
627 boundaryField_(*this, gf.boundaryField_)
630 <<
"Copy construct" <<
nl << this->info() <<
endl;
640 this->writeOpt(IOobject::NO_WRITE);
644 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
650 Internal(tgf.constCast(), tgf.movable()),
653 fieldPrevIterPtr_(nullptr),
654 boundaryField_(*this, tgf().boundaryField_)
657 <<
"Constructing from tmp" <<
nl << this->info() <<
endl;
659 this->writeOpt(IOobject::NO_WRITE);
665 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
675 fieldPrevIterPtr_(nullptr),
676 boundaryField_(*this, gf.boundaryField_)
679 <<
"Copy construct, resetting IO params" <<
nl 680 << this->info() <<
endl;
693 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
697 const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf
700 Internal(
io, tgf.constCast(), tgf.movable()),
703 fieldPrevIterPtr_(nullptr),
704 boundaryField_(*this, tgf().boundaryField_)
707 <<
"Constructing from tmp resetting IO params" <<
nl 708 << this->info() <<
endl;
716 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
723 Internal(newName, gf),
726 fieldPrevIterPtr_(nullptr),
727 boundaryField_(*this, gf.boundaryField_)
730 <<
"Copy construct, resetting name" <<
nl 731 << this->info() <<
endl;
744 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
748 const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf
751 Internal(newName, tgf.constCast(), tgf.movable()),
754 fieldPrevIterPtr_(nullptr),
755 boundaryField_(*this, tgf().boundaryField_)
758 <<
"Constructing from tmp resetting name" <<
nl 759 << this->info() <<
endl;
765 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
770 const word& patchFieldType
776 fieldPrevIterPtr_(nullptr),
777 boundaryField_(this->
mesh().
boundary(), *this, patchFieldType)
780 <<
"Copy construct, resetting IO params" <<
nl 781 << this->info() <<
endl;
783 boundaryField_ == gf.boundaryField_;
796 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
800 const GeometricField<Type, PatchField, GeoMesh>& gf,
808 fieldPrevIterPtr_(nullptr),
818 <<
"Copy construct, resetting IO params and patch types" <<
nl 819 << this->info() <<
endl;
821 boundaryField_ == gf.boundaryField_;
825 field0Ptr_ =
new GeometricField<Type, PatchField, GeoMesh>
834 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
838 const GeometricField<Type, PatchField, GeoMesh>& gf,
840 const word& patchFieldType
846 fieldPrevIterPtr_(nullptr),
847 boundaryField_(*this, gf.boundaryField_,
patchIDs, patchFieldType)
850 <<
"Copy construct, resetting IO params and setting patchFieldType " 851 <<
"for patchIDs" <<
nl 852 << this->info() <<
endl;
856 field0Ptr_ =
new GeometricField<Type, PatchField, GeoMesh>
865 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
869 const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf,
874 Internal(
io, tgf.constCast(), tgf.movable()),
877 fieldPrevIterPtr_(nullptr),
887 <<
"Constructing from tmp resetting IO params and patch types" <<
nl 888 << this->info() <<
endl;
890 boundaryField_ == tgf().boundaryField_;
896 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
906 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
916 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
921 const bool updateAccessTime
924 if (updateAccessTime)
933 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
938 const bool updateAccessTime
941 if (updateAccessTime)
950 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
955 const bool updateAccessTime
958 if (updateAccessTime)
963 return boundaryField_;
967 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
973 && timeIndex_ != this->time().
timeIndex()
974 && !this->
name().ends_with(
"_0")
978 timeIndex_ = this->time().timeIndex();
986 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
991 field0Ptr_->storeOldTime();
994 <<
"Storing old time field for field" <<
nl << this->info() <<
endl;
996 *field0Ptr_ == *
this;
997 field0Ptr_->timeIndex_ = timeIndex_;
999 if (field0Ptr_->field0Ptr_)
1001 field0Ptr_->writeOpt(this->writeOpt());
1007 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1012 return field0Ptr_->nOldTimes() + 1;
1019 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1029 this->
name() +
"_0",
1034 this->registerObject()
1042 <<
"created old time field " << field0Ptr_->info() <<
endl;
1046 error::printStack(
Info);
1059 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1070 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1073 if (!fieldPrevIterPtr_)
1076 <<
"Allocating previous iteration field" <<
nl 1077 << this->info() <<
endl;
1081 this->
name() +
"PrevIter",
1087 *fieldPrevIterPtr_ == *
this;
1092 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1096 if (!fieldPrevIterPtr_)
1099 <<
"previous iteration field" <<
endl << this->info() <<
endl 1101 <<
" Use field.storePrevIter() at start of iteration." 1105 return *fieldPrevIterPtr_;
1109 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1113 this->setUpToDate();
1115 boundaryField_.evaluate();
1119 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1125 bool needRef =
true;
1127 for (
const auto& pf : boundaryField_)
1129 if (pf.fixesValue())
1140 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1150 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1157 this->
mesh().data::template getOrDefault<bool>
1174 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1182 return this->
name() +
"Final";
1185 return this->
name();
1189 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1202 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1213 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1219 this->
name() +
".T()",
1224 Foam::T(tresult.ref().primitiveFieldRef(), primitiveField());
1225 Foam::T(tresult.ref().boundaryFieldRef(), boundaryField());
1231 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1253 Foam::component(tresult.ref().primitiveFieldRef(), primitiveField(), d);
1254 Foam::component(tresult.ref().boundaryFieldRef(), boundaryField(), d);
1260 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1264 const GeometricField
1266 typename GeometricField<Type, PatchField, GeoMesh>::cmptType,
1277 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1289 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1300 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1311 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1321 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1327 this->clamp_max(
upper.value());
1331 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1343 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1354 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1361 this->clamp_range(
lower.value(),
upper.value());
1365 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1371 this->clamp_range(
range.value());
1375 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1383 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1393 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1413 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1419 const auto& gf = tgf();
1449 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1455 internalFieldRef() = dt;
1460 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1466 const auto& gf = tgf();
1479 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1485 internalFieldRef() = dt;
1490 #define COMPUTED_ASSIGNMENT(TYPE, op) \ 1492 template<class Type, template<class> class PatchField, class GeoMesh> \ 1493 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \ 1495 const GeometricField<TYPE, PatchField, GeoMesh>& gf \ 1498 checkField(*this, gf, #op); \ 1500 internalFieldRef() op gf.internalField(); \ 1501 boundaryFieldRef() op gf.boundaryField(); \ 1504 template<class Type, template<class> class PatchField, class GeoMesh> \ 1505 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \ 1507 const tmp<GeometricField<TYPE, PatchField, GeoMesh>>& tgf \ 1510 operator op(tgf()); \ 1514 template<class Type, template<class> class PatchField, class GeoMesh> \ 1515 void Foam::GeometricField<Type, PatchField, GeoMesh>::operator op \ 1517 const dimensioned<TYPE>& dt \ 1520 internalFieldRef() op dt; \ 1521 boundaryFieldRef() op dt.value(); \ 1529 #undef COMPUTED_ASSIGNMENT 1534 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1550 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
void clamp_min(const Type &lower)
Impose lower (floor) clamp on the field values (in-place)
void clamp_range(const dimensioned< MinMax< Type >> &range)
Clamp field values (in-place) to the specified range.
#define COMPUTED_ASSIGNMENT(TYPE, op)
const Type & value() const noexcept
Return const reference to value.
const labelList patchIDs(pbm.patchSet(polyPatchNames, false, true).sortedToc())
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
errorManipArg< error, int > exit(error &err, const int errNo=1)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Y [inertIndex] clamp_min(0)
bool writeData(Ostream &) const
WriteData member function required by regIOobject.
Field< Type >::cmptType cmptType
Component type of the field elements.
const Internal & internalField() const noexcept
Return a const-reference to the dimensioned internal field.
constexpr char nl
The newline '\n' character (0x0a)
A min/max value pair with additional methods. In addition to conveniently storing values...
GeometricField(const IOobject &io, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=PatchField< Type >::calculatedType())
Construct given IOobject, mesh, dimensions and patch type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
orientedType oriented() const noexcept
Return oriented type.
string upper(const std::string &s)
Return string copy transformed with std::toupper on each character.
Generic GeometricField class.
Generic dimensioned Type class.
const dimensionSet dimless
Dimensionless.
label nOldTimes() const
Return the number of old time fields stored.
word select(bool final) const
Select the final iteration parameters if final is true by returning the field name + "Final" otherwis...
conserve primitiveFieldRef()+
MinMax< label > minMax(const labelHashSet &set)
Find the min/max values of labelHashSet.
#define forAll(list, i)
Loop across all elements in list.
bool needReference() const
Does the field need a reference level for solution.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
const T & min() const noexcept
The min value (first)
bool returnReduceAnd(const bool value, const label comm=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
void storeOldTimes() const
Store the old-time fields.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Generic templated field type.
A class for handling words, derived from Foam::string.
#define DebugInFunction
Report an information message using Foam::Info.
void storeOldTime() const
Store the old-time field.
void writeMinMax(Ostream &os) const
Helper function to write the min and max to an Ostream.
errorManip< error > abort(error &err)
Us boundaryFieldRef().evaluateCoupled< coupledFaPatch >()
void negate()
Negate the field inplace. See notes in Field.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Internal & internalFieldRef(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
int debug
Static debugging option.
void normalise()
Normalise the field inplace. See notes in Field.
const GeometricField< Type, PatchField, GeoMesh > & prevIter() const
Return previous iteration field.
OBJstream os(runTime.globalPath()/outputName)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
virtual ~GeometricField()
Destructor.
Generic GeometricBoundaryField class.
void clamp_max(const Type &upper)
Impose upper (ceiling) clamp on the field values (in-place)
List< word > wordList
List of word.
Template functions to aid in the implementation of demand driven data.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
tmp< GeometricField< Type, PatchField, GeoMesh > > clone() const
Clone.
string lower(const std::string &s)
Return string copy transformed with std::tolower on each character.
#define WarningInFunction
Report a warning using Foam::Warning.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
meshDefDict readIfPresent("polyMeshPatches", polyPatchNames)
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
messageStream Info
Information stream (stdout output on master, null elsewhere)
void storePrevIter() const
Store the field as the previous iteration value.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the templated type.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
List< label > labelList
A List of labels.
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
A class for managing temporary objects.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
void deleteDemandDrivenData(DataPtr &dataPtr)
Defines the attributes of an object for which implicit objectRegistry management is supported...
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
void relax()
Relax field (for steady-state solution).
const dimensionSet & dimensions() const noexcept
Return dimensions.
#define checkField(fld1, fld2, op)
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
#define InfoInFunction
Report an information message using Foam::Info.
void replace(const direction d, const GeometricField< cmptType, PatchField, GeoMesh > &gcf)
Replace specified field component with content from another field.