36template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
51 this->
resize(bmesh_.size());
53 label nUnset = this->size();
60 if (dEntry.isDict() && dEntry.keyword().isLiteral())
62 const label patchi = bmesh_.findPatchID(dEntry.keyword());
92 for (
auto iter =
dict.crbegin(); iter !=
dict.crend(); ++iter)
94 const entry& dEntry = *iter;
99 bmesh_.indices(dEntry.
keyword(),
true);
103 if (!this->set(patchi))
108 PatchField<Type>::New
124 if (!this->set(patchi))
126 if (bmesh_[patchi].type() == emptyPolyPatch::typeName)
131 PatchField<Type>::New
133 emptyPolyPatch::typeName,
148 PatchField<Type>::New
152 dict.subDict(bmesh_[patchi].name())
164 if (!this->set(patchi))
166 if (bmesh_[patchi].type() == cyclicPolyPatch::typeName)
169 <<
"Cannot find patchField entry for cyclic "
170 << bmesh_[patchi].name() <<
endl
171 <<
"Is your field uptodate with split cyclics?" <<
endl
172 <<
"Run foamUpgradeCyclics to convert mesh and fields"
173 <<
" to split cyclics." <<
exit(FatalIOError);
178 <<
"Cannot find patchField entry for "
179 << bmesh_[patchi].name() <<
exit(FatalIOError);
188template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
199template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
204 const word& patchFieldType
220 PatchField<Type>::New
231template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
250 patchFieldTypes.
size() != this->size()
251 || (constraintTypes.
size() && (constraintTypes.
size() != this->size()))
255 <<
"Incorrect number of patch type specifications given" <<
nl
256 <<
" Number of patches in mesh = " << bmesh.size()
257 <<
" number of patch type specifications = "
258 << patchFieldTypes.
size()
262 if (constraintTypes.
size())
269 PatchField<Type>::New
271 patchFieldTypes[patchi],
272 constraintTypes[patchi],
286 PatchField<Type>::New
288 patchFieldTypes[patchi],
298template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
303 const PtrList<PatchField<Type>>& ptfl
321template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
343template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
349 const word& patchFieldType
360 for (
const label patchi : patchIDs)
365 PatchField<Type>::New
376 if (!this->
set(patchi))
384template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
400template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
417template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
425 for (
auto& pfld : *
this)
432template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
450 for (
auto& pfld : *
this)
452 pfld.initEvaluate(commsType);
465 for (
auto& pfld : *
this)
467 pfld.evaluate(commsType);
473 bmesh_.mesh().globalData().patchSchedule();
475 for (
const auto& schedEval : patchSchedule)
477 const label patchi = schedEval.patch;
478 auto& pfld = (*this)[patchi];
482 pfld.initEvaluate(commsType);
486 pfld.evaluate(commsType);
493 <<
"Unsupported communications type "
500template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
501template<
class CoupledPatchType>
519 for (
auto& pfld : *
this)
521 const auto* cpp = isA<CoupledPatchType>(pfld.patch());
523 if (cpp && cpp->coupled())
525 pfld.initEvaluate(commsType);
539 for (
auto& pfld : *
this)
541 const auto* cpp = isA<CoupledPatchType>(pfld.patch());
543 if (cpp && cpp->coupled())
545 pfld.evaluate(commsType);
552 bmesh_.mesh().globalData().patchSchedule();
554 for (
const auto& schedEval : patchSchedule)
556 const label patchi = schedEval.patch;
557 auto& pfld = (*this)[patchi];
559 const auto* cpp = isA<CoupledPatchType>(pfld.patch());
561 if (cpp && cpp->coupled())
565 pfld.initEvaluate(commsType);
569 pfld.evaluate(commsType);
577 <<
"Unsupported communications type "
584template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
594 list[patchi] = pff[patchi].type();
601template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
610 result[patchi] == this->operator[](patchi).patchInternalField();
617template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
626 isA<LduInterfaceField<Type>>(this->operator[](patchi));
630 list.
set(patchi, lduPtr);
638template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
648 isA<lduInterfaceField>(this->
operator[](patchi));
652 list.
set(patchi, lduPtr);
660template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
668 this->writeEntries(
os);
675template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
681 for (
const auto& pfld : *
this)
692template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
702template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
712template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
722template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
730 this->operator[](patchi) == bf[patchi];
735template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
743 this->operator[](patchi) == bf[patchi];
748template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
756 this->operator[](patchi) == val;
763template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
770 os << static_cast<const FieldField<PatchField, Type>&>(bf);
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A field of fields is a PtrList of fields with reference counting.
tmp< FieldField< PatchField, Type > > clone() const
Clone.
void operator=(const FieldField< Field, Type > &)
Copy assignment.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Generic GeometricBoundaryField class.
lduInterfaceFieldPtrsList scalarInterfaces() const
void readField(const DimensionedField< Type, GeoMesh > &field, const dictionary &dict)
Read the boundary field.
wordList types() const
Return a list of the patch types.
void evaluateCoupled()
Evaluate boundary conditions on a subset of coupled patches.
void evaluate()
Evaluate boundary conditions.
LduInterfaceFieldPtrsList< Type > interfaces() const
void writeEntries(Ostream &os) const
Write dictionary entries of the individual boundary fields.
void writeEntry(const word &keyword, Ostream &os) const
Write boundary field as dictionary entry.
void updateCoeffs()
Update the boundary condition coefficients.
GeoMesh::BoundaryMesh BoundaryMesh
The boundary mesh type for the boundary fields.
GeometricBoundaryField boundaryInternalField() const
Return boundary field of values neighbouring the boundary.
virtual bool check(const char *operation) const
Check IOstream status for given operation.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
virtual Ostream & endBlock()
Write end block group.
virtual Ostream & beginBlock(const keyType &kw)
Write begin block group with the given name.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
const T * set(const label i) const
void size(const label n)
Older name for setAddressableSize.
commsTypes
Types of communications.
@ nonBlocking
"nonBlocking"
static const Enum< commsTypes > commsTypeNames
Names of the communication types.
static label nRequests()
Get number of outstanding requests.
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
static commsTypes defaultCommsType
Default commsType.
static bool & parRun() noexcept
Test if this a parallel run.
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
const T * set(const label i) const
label size() const noexcept
The number of elements in the list.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
A keyword and a list of tokens is an 'entry'.
virtual bool isDict() const noexcept
Return true if this entry is a dictionary.
const keyType & keyword() const noexcept
Return keyword.
virtual const dictionary & dict() const =0
Return dictionary, if entry is a dictionary.
bool isLiteral() const noexcept
The keyType is treated as literal, not as pattern.
A class for handling words, derived from Foam::string.
patchWriters resize(patchIds.size())
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
OBJstream os(runTime.globalPath()/outputName)
Ostream & endl(Ostream &os)
Add newline and flush stream.
errorManip< error > abort(error &err)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
errorManipArg< error, int > exit(error &err, const int errNo=1)
constexpr char nl
The newline '\n' character (0x0a)
#define forAll(list, i)
Loop across all elements in list.