Go to the documentation of this file.
31 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
46 Info<<
"GeometricField<Type, PatchField, GeoMesh>::"
47 "GeometricBoundaryField::readField"
49 "const DimensionedField<Type, GeoMesh>&, "
56 label nUnset = this->size();
62 if (iter().isDict() && !iter().keyword().isPattern())
64 label patchi = bmesh_.findPatchID(iter().keyword());
105 if (
e.isDict() && !
e.keyword().isPattern())
107 const labelList patchIDs = bmesh_.findIndices
141 if (bmesh_[
patchi].
type() == emptyPolyPatch::typeName)
148 emptyPolyPatch::typeName,
181 if (bmesh_[
patchi].
type() == cyclicPolyPatch::typeName)
186 ) <<
"Cannot find patchField entry for cyclic "
188 <<
"Is your field uptodate with split cyclics?" <<
endl
189 <<
"Run foamUpgradeCyclics to convert mesh and fields"
197 ) <<
"Cannot find patchField entry for "
207 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
219 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
225 const word& patchFieldType
233 Info<<
"GeometricField<Type, PatchField, GeoMesh>::"
234 "GeometricBoundaryField::"
235 "GeometricBoundaryField(const BoundaryMesh&, "
236 "const DimensionedField<Type>&, const word&)"
256 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
261 const DimensionedField<Type, GeoMesh>& field,
266 FieldField<PatchField, Type>(bmesh.size()),
271 Info<<
"GeometricField<Type, PatchField, GeoMesh>::"
272 "GeometricBoundaryField::"
273 "GeometricBoundaryField"
275 "const BoundaryMesh&, "
276 "const DimensionedField<Type>&, "
285 patchFieldTypes.size() != this->size()
286 || (constraintTypes.size() && (constraintTypes.size() != this->size()))
290 <<
"Incorrect number of patch type specifications given" <<
nl
291 <<
" Number of patches in mesh = " << bmesh.size()
292 <<
" number of patch type specifications = "
293 << patchFieldTypes.size()
297 if (constraintTypes.size())
333 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
338 const DimensionedField<Type, GeoMesh>& field,
339 const PtrList<PatchField<Type> >& ptfl
342 FieldField<PatchField, Type>(bmesh.size()),
347 Info<<
"GeometricField<Type, PatchField, GeoMesh>::"
348 "GeometricBoundaryField::"
349 "GeometricBoundaryField"
351 "const BoundaryMesh&, "
352 "const DimensionedField<Type, GeoMesh>&, "
353 "const PtrLIst<PatchField<Type> >&"
365 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
369 const DimensionedField<Type, GeoMesh>& field,
370 const typename GeometricField<Type, PatchField, GeoMesh>::
371 GeometricBoundaryField& btf
374 FieldField<PatchField, Type>(btf.size()),
379 Info<<
"GeometricField<Type, PatchField, GeoMesh>::"
380 "GeometricBoundaryField::"
381 "GeometricBoundaryField"
383 "const DimensionedField<Type, GeoMesh>&, "
384 "const typename GeometricField<Type, PatchField, GeoMesh>::"
385 "GeometricBoundaryField&"
402 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
406 const typename GeometricField<Type, PatchField, GeoMesh>::
407 GeometricBoundaryField& btf
410 FieldField<PatchField, Type>(btf),
415 Info<<
"GeometricField<Type, PatchField, GeoMesh>::"
416 "GeometricBoundaryField::"
417 "GeometricBoundaryField"
419 "const GeometricField<Type, PatchField, GeoMesh>::"
420 "GeometricBoundaryField&"
427 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
432 const DimensionedField<Type, GeoMesh>& field,
433 const dictionary&
dict
436 FieldField<PatchField, Type>(bmesh.size()),
445 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
451 Info<<
"GeometricField<Type, PatchField, GeoMesh>::"
452 "GeometricBoundaryField::"
453 "updateCoeffs()" <<
endl;
458 this->operator[](
patchi).updateCoeffs();
463 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
469 Info<<
"GeometricField<Type, PatchField, GeoMesh>::"
470 "GeometricBoundaryField::"
471 "evaluate()" <<
endl;
505 bmesh_.mesh().globalData().patchSchedule();
507 forAll(patchSchedule, patchEvali)
509 if (patchSchedule[patchEvali].init)
511 this->operator[](patchSchedule[patchEvali].patch)
516 this->operator[](patchSchedule[patchEvali].patch)
524 <<
"Unsuported communications type "
531 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
549 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
555 BoundaryInternalField(*
this);
559 BoundaryInternalField[
patchi] ==
560 this->operator[](
patchi).patchInternalField();
563 return BoundaryInternalField;
567 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
593 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
602 if (isA<lduInterfaceField>(this->
operator[](
patchi)))
607 &refCast<const lduInterfaceField>
619 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
638 "GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField::"
639 "writeEntry(const word& keyword, Ostream& os) const"
646 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
651 GeometricBoundaryField& bf
658 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
669 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
681 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
685 const typename GeometricField<Type, PatchField, GeoMesh>::
686 GeometricBoundaryField& bf
691 this->operator[](patchI) == bf[patchI];
696 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
700 const FieldField<PatchField, Type>& ptff
705 this->operator[](patchI) == ptff[patchI];
710 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
719 this->operator[](patchI) == t;
726 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
730 const typename GeometricField<Type, PatchField, GeoMesh>::
731 GeometricBoundaryField& bf
734 os << static_cast<const FieldField<PatchField, Type>&>(bf);
A keyword and a list of tokens is an 'entry'.
LduInterfaceFieldPtrsList< Type > interfaces() const
Return a list of pointers for each patch field with only those.
void evaluate()
Evaluate boundary conditions.
points setSize(newPointi)
void updateCoeffs()
Update the boundary condition coefficients.
A class for handling words, derived from string.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
#define forAll(list, i)
Loop across all elements in list.
static bool & parRun()
Is this a parallel run?
static const NamedEnum< commsTypes, 3 > commsTypeNames
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
void readField(const DimensionedField< Type, GeoMesh > &field, const dictionary &dict)
Read the boundary field.
Ostream & endl(Ostream &os)
Add newline and flush stream.
GeoMesh::BoundaryMesh BoundaryMesh
tmp< FieldField< Field, Type > > clone() const
Clone.
static commsTypes defaultCommsType
Default commsType.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
const BoundaryMesh & bmesh_
Reference to BoundaryMesh for which this field is defined.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
lduInterfaceFieldPtrsList scalarInterfaces() const
Return a list of pointers for each patch field with only those.
List< word > wordList
A List of words.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Intrusive doubly-linked list.
GeometricBoundaryField(const BoundaryMesh &)
Construct from a BoundaryMesh.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
An abstract base class for implicitly-coupled interface fields e.g. processor and cyclic patch fields...
void operator=(const FieldField< Field, Type > &)
virtual bool check(const char *operation) const
Check IOstream status for given operation.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
errorManip< error > abort(error &err)
Ostream & indent(Ostream &os)
Indent stream.
To & refCast(From &r)
Reference type cast template function.
errorManipArg< error, int > exit(error &err, const int errNo=1)
wordList types() const
Return a list of the patch types.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
static label nRequests()
Get number of outstanding requests.
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
bool set(const label) const
Is element set.
const dimensionedScalar e
Elementary charge.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
GeometricBoundaryField boundaryInternalField() const
Return BoundaryField of the cell values neighbouring.
friend Ostream & operator(Ostream &, const GeometricField< Type, PatchField, GeoMesh > &)
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Generic GeometricField class.
virtual const fileName & name() const
Return the name of the stream.
void writeEntry(const word &keyword, Ostream &os) const
Write boundary field as dictionary entry.
word name(const complex &)
Return a string representation of a complex.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...