Loading...
Searching...
No Matches
singleLayerRegion.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2020-2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "singleLayerRegion.H"
30#include "fvMesh.H"
31#include "Time.H"
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace regionModels
39{
41}
42}
43
44// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
45
46void Foam::regionModels::singleLayerRegion::constructMeshObjects()
47{
48 // construct patch normal vectors
49 nHatPtr_.reset
50 (
52 (
54 (
55 "nHat",
57 regionMesh(),
61 ),
62 regionMesh(),
65 )
66 );
67
68 // construct patch areas
69 magSfPtr_.reset
70 (
72 (
74 (
75 "magSf",
77 regionMesh(),
81 ),
82 regionMesh(),
85 )
86 );
87}
88
89
90void Foam::regionModels::singleLayerRegion::initialise()
91{
92 if (debug)
93 {
94 Pout<< "singleLayerRegion::initialise()" << endl;
95 }
96
97 label nBoundaryFaces = 0;
98 const polyBoundaryMesh& rbm = regionMesh().boundaryMesh();
99 volVectorField& nHat = nHatPtr_();
100 volScalarField& magSf = magSfPtr_();
101 forAll(intCoupledPatchIDs_, i)
102 {
103 const label patchi = intCoupledPatchIDs_[i];
104 const polyPatch& pp = rbm[patchi];
105 const labelUList& fCells = pp.faceCells();
106
107 nBoundaryFaces += fCells.size();
108
109 UIndirectList<vector>(nHat, fCells) = pp.faceNormals();
110 UIndirectList<scalar>(magSf, fCells) = mag(pp.faceAreas());
111 }
112 nHat.correctBoundaryConditions();
113 magSf.correctBoundaryConditions();
114
115 if (nBoundaryFaces != regionMesh().nCells())
116 {
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
122 << abort(FatalError);
123 }
124
125 scalarField passiveMagSf(magSf.size(), Zero);
126 passivePatchIDs_.setSize(intCoupledPatchIDs_.size(), -1);
127 forAll(intCoupledPatchIDs_, i)
128 {
129 const label patchi = intCoupledPatchIDs_[i];
130 const polyPatch& ppIntCoupled = rbm[patchi];
131 if (ppIntCoupled.size() > 0)
132 {
133 label cellId = rbm[patchi].faceCells()[0];
134 const cell& cFaces = regionMesh().cells()[cellId];
135
136 label facei = ppIntCoupled.start();
137 label faceO = cFaces.opposingFaceLabel(facei, regionMesh().faces());
138
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());
144 }
145 }
146
147 Pstream::listReduce(passivePatchIDs_, maxOp<label>());
148
149 magSf.field() = 0.5*(magSf + passiveMagSf);
150 magSf.correctBoundaryConditions();
151}
152
153
154// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
155
158 return regionModel::read();
159}
160
161
162// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
163
164Foam::regionModels::singleLayerRegion::singleLayerRegion
165(
166 const fvMesh& mesh,
167 const word& regionType
168)
169:
170 regionModel(mesh, regionType),
171 nHatPtr_(nullptr),
172 magSfPtr_(nullptr),
174{}
175
176
177Foam::regionModels::singleLayerRegion::singleLayerRegion
178(
179 const fvMesh& mesh,
180 const word& regionType,
181 const word& modelName,
182 bool readFields
183)
184:
185 regionModel(mesh, regionType, modelName, false),
186 nHatPtr_(nullptr),
187 magSfPtr_(nullptr),
188 passivePatchIDs_()
189{
190 if (active_)
191 {
192 constructMeshObjects();
193 initialise();
194
195 if (readFields)
196 {
197 read();
199 }
200}
201
202
203// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
206{}
207
208
209// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
210
212{
213 if (!nHatPtr_)
214 {
216 << "Region patch normal vectors not available"
218 }
219
220 return *nHatPtr_;
221}
222
223
225{
226 if (!magSfPtr_)
227 {
229 << "Region patch areas not available"
230 << abort(FatalError);
232
233 return *magSfPtr_;
234}
235
236
237const Foam::labelList&
239{
240 return passivePatchIDs_;
241}
242
243
244// ************************************************************************* //
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,...
Definition IOobject.H:191
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.
Definition Time.C:714
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Base class for region models.
Definition regionModel.H:59
const Time & time_
Reference to the time database.
Definition regionModel.H:97
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 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.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
label cellId
Namespace for OpenFOAM.
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.
Definition List.H:62
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.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition errorManip.H:139
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
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.
Definition UList.H:75
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299