43void Foam::regionModels::regionModel1D::constructMeshObjects()
65void Foam::regionModels::regionModel1D::initialise()
69 Pout<<
"regionModel1D::initialise()" <<
endl;
74 DynamicList<label> faceIDs;
77 label localPyrolysisFacei = 0;
79 const polyBoundaryMesh& rbm = regionMesh().boundaryMesh();
81 forAll(intCoupledPatchIDs_, i)
83 const label patchi = intCoupledPatchIDs_[i];
84 const polyPatch& ppCoupled = rbm[patchi];
85 localPyrolysisFacei += ppCoupled.size();
88 boundaryFaceOppositeFace_.setSize(localPyrolysisFacei);
89 boundaryFaceFaces_.setSize(localPyrolysisFacei);
90 boundaryFaceCells_.setSize(localPyrolysisFacei);
92 localPyrolysisFacei = 0;
94 forAll(intCoupledPatchIDs_, i)
96 const label patchi = intCoupledPatchIDs_[i];
97 const polyPatch& ppCoupled = rbm[patchi];
98 forAll(ppCoupled, localFacei)
100 label facei = ppCoupled.start() + localFacei;
105 label ownCelli = regionMesh().faceOwner()[facei];
106 if (ownCelli != celli)
112 celli = regionMesh().faceNeighbour()[facei];
116 const cell& cFaces = regionMesh().cells()[celli];
117 faceIDs.append(facei);
119 cFaces.opposingFaceLabel(facei, regionMesh().faces());
121 }
while (regionMesh().isInternalFace(facei));
123 boundaryFaceOppositeFace_[localPyrolysisFacei] = facei;
126 boundaryFaceFaces_[localPyrolysisFacei].transfer(faceIDs);
127 boundaryFaceCells_[localPyrolysisFacei].transfer(
cellIDs);
129 localPyrolysisFacei++;
140 localPyrolysisFacei = 0;
142 forAll(intCoupledPatchIDs_, i)
144 const label patchi = intCoupledPatchIDs_[i];
145 const polyPatch& ppCoupled = rbm[patchi];
146 const vectorField& pNormals = ppCoupled.faceNormals();
148 nMagSfBf[patchi] = regionMesh().Sf().boundaryField()[patchi] & pNormals;
150 forAll(pNormals, localFacei)
152 const vector n = pNormals[localFacei];
153 const labelList& faces = boundaryFaceFaces_[localPyrolysisFacei++];
159 const label faceID = faces[facei];
160 nMagSf[faceID] = regionMesh().Sf()[faceID] &
n;
185 moveMesh_.readIfPresent(
"moveMesh", coeffs_);
197 const scalar minDelta
213 label totalFacei = 0;
214 forAll(intCoupledPatchIDs_, localPatchi)
216 label patchi = intCoupledPatchIDs_[localPatchi];
221 const labelList& faces = boundaryFaceFaces_[totalFacei];
223 const label oFace = boundaryFaceOppositeFace_[totalFacei];
233 oldCf[i] = regionMesh().faceCentres()[faces[i]];
236 oldCf[faces.
size()] = regionMesh().faceCentres()[oFace];
240 const label celli =
cells[i];
242 if (
mag(oldCf[i + 1] - oldCf[i]) < minDelta)
245 cellMoveMap[celli] = 1;
254 const label celli =
cells[i];
255 newDelta[j+1] = (deltaV[celli]/
mag(sf))*
n + newDelta[j];
260 const face of = regionMesh().faces()[oFace];
262 scalar omagV =
mag(newDelta[newDelta.
size()-1]);
264 if (!frozen[
cells.size()-1] && (omagV > ROOTVSMALL))
268 const label pointi = of[pti];
270 oldPoints[pointi] - newDelta[newDelta.
size()-1];
275 for (label i=0; i < faces.
size(); i++)
277 const label facei = faces[i];
278 const face f = regionMesh().faces()[facei];
280 scalar magV =
mag(newDelta[i]);
281 if (!frozen[i] && magV > 0)
285 const label pointi =
f[pti];
286 newPoints[pointi] = oldPoints[pointi] - newDelta[i];
296 regionMesh().movePoints(newPoints);
307 const word& regionType
311 boundaryFaceFaces_(),
312 boundaryFaceCells_(),
313 boundaryFaceOppositeFace_(),
323 const word& regionType,
324 const word& modelName,
329 boundaryFaceFaces_(regionMesh().nCells()),
330 boundaryFaceCells_(regionMesh().nCells()),
331 boundaryFaceOppositeFace_(regionMesh().nCells()),
338 constructMeshObjects();
351 const word& regionType,
352 const word& modelName,
358 boundaryFaceFaces_(regionMesh().nCells()),
359 boundaryFaceCells_(regionMesh().nCells()),
360 boundaryFaceOppositeFace_(regionMesh().nCells()),
367 constructMeshObjects();
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
Type of boundary fields.
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
bool readIfPresent(const word &key, const dictionary &dict)
Update the value of the Switch if it is found in the dictionary.
void size(const label n)
Older name for setAddressableSize.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
A face is a list of labels corresponding to mesh vertices.
Mesh data needed to do the Finite Volume discretisation.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
A patch is a list of labels that address the faces in the global face list.
const vectorField::subField faceAreas() const
Return face normals.
Base class for 1-D region models.
autoPtr< surfaceScalarField > nMagSfPtr_
Face area magnitude normal to patch.
virtual ~regionModel1D()
Destructor.
Switch moveMesh_
Flag to allow mesh movement.
virtual bool read()
Read control parameters from dictionary.
Base class for region models.
Switch active_
Active flag.
const Time & time() const
Return the reference to the time database.
dictionary coeffs_
Model coefficients dictionary.
const fvMesh & regionMesh() const
Return the region mesh database.
virtual bool read()
Read control parameters from dictionary.
A class for managing temporary objects.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
List< label > labelList
A List of labels.
const dimensionSet dimArea(sqr(dimLength))
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
Field< label > labelField
Specialisation of Field<T> for label.
static constexpr const zero Zero
Global zero (0)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
#define forAll(list, i)
Loop across all elements in list.