46void Foam::regionModels::singleLayerRegion::constructMeshObjects()
90void Foam::regionModels::singleLayerRegion::initialise()
94 Pout<<
"singleLayerRegion::initialise()" <<
endl;
97 label nBoundaryFaces = 0;
98 const polyBoundaryMesh& rbm = regionMesh().boundaryMesh();
101 forAll(intCoupledPatchIDs_, i)
103 const label patchi = intCoupledPatchIDs_[i];
104 const polyPatch&
pp = rbm[patchi];
107 nBoundaryFaces += fCells.
size();
109 UIndirectList<vector>(nHat, fCells) =
pp.faceNormals();
110 UIndirectList<scalar>(magSf, fCells) =
mag(
pp.faceAreas());
112 nHat.correctBoundaryConditions();
113 magSf.correctBoundaryConditions();
115 if (nBoundaryFaces != regionMesh().nCells())
118 <<
"Number of primary region coupled boundary faces not equal to "
119 <<
"the number of cells in the local region" <<
nl <<
nl
120 <<
"Number of cells = " << regionMesh().nCells() <<
nl
121 <<
"Boundary faces = " << nBoundaryFaces <<
nl
126 passivePatchIDs_.setSize(intCoupledPatchIDs_.size(), -1);
127 forAll(intCoupledPatchIDs_, i)
129 const label patchi = intCoupledPatchIDs_[i];
130 const polyPatch& ppIntCoupled = rbm[patchi];
131 if (ppIntCoupled.size() > 0)
133 label
cellId = rbm[patchi].faceCells()[0];
134 const cell& cFaces = regionMesh().cells()[
cellId];
136 label facei = ppIntCoupled.start();
137 label faceO = cFaces.opposingFaceLabel(facei, regionMesh().faces());
139 label passivePatchi = rbm.whichPatch(faceO);
140 passivePatchIDs_[i] = passivePatchi;
141 const polyPatch& ppPassive = rbm[passivePatchi];
142 UIndirectList<scalar>(passiveMagSf, ppPassive.faceCells()) =
143 mag(ppPassive.faceAreas());
150 magSf.correctBoundaryConditions();
164Foam::regionModels::singleLayerRegion::singleLayerRegion
167 const word& regionType
177Foam::regionModels::singleLayerRegion::singleLayerRegion
180 const word& regionType,
181 const word& modelName,
192 constructMeshObjects();
216 <<
"Region patch normal vectors not available"
229 <<
"Region patch areas not available"
240 return passivePatchIDs_;
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
@ NO_REGISTER
Do not request registration (bool: false).
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
static void listReduce(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce list elements (list must be equal size on all ranks), applying bop to each list element.
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
void size(const label n)
Older name for setAddressableSize.
Mesh data needed to do the Finite Volume discretisation.
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Base class for region models.
Switch active_
Active flag.
const Time & time_
Reference to the time database.
const fvMesh & regionMesh() const
Return the region mesh database.
const word & modelName() const noexcept
Return the model name.
virtual bool read()
Read control parameters from dictionary.
Base class for single layer region models.
autoPtr< volVectorField > nHatPtr_
Patch normal vectors.
virtual const labelList & passivePatchIDs() const
Return the list of patch IDs opposite to internally.
virtual ~singleLayerRegion()
Destructor.
virtual const volScalarField & magSf() const
Return the face area magnitudes / [m2].
labelList passivePatchIDs_
List of patch IDs opposite to internally coupled patches.
virtual const volVectorField & nHat() const
Return the patch normal vectors.
autoPtr< volScalarField > magSfPtr_
Face area magnitudes / [m2].
virtual bool read()
Read control parameters from dictionary.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
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.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
UList< label > labelUList
A UList of labels.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.