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 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"
33
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 (
53 IOobject
54 (
55 "nHat",
57 regionMesh(),
60 ),
61 regionMesh(),
64 )
65 );
66
67 // construct patch areas
68 magSfPtr_.reset
69 (
71 (
72 IOobject
73 (
74 "magSf",
76 regionMesh(),
79 ),
80 regionMesh(),
82 zeroGradientFvPatchField<scalar>::typeName
83 )
84 );
85}
86
87
88void Foam::regionModels::singleLayerRegion::initialise()
89{
90 if (debug)
91 {
92 Pout<< "singleLayerRegion::initialise()" << endl;
93 }
94
95 label nBoundaryFaces = 0;
96 const polyBoundaryMesh& rbm = regionMesh().boundaryMesh();
97 volVectorField& nHat = nHatPtr_();
98 volScalarField& magSf = magSfPtr_();
99 forAll(intCoupledPatchIDs_, i)
100 {
101 const label patchi = intCoupledPatchIDs_[i];
102 const polyPatch& pp = rbm[patchi];
103 const labelList& fCells = pp.faceCells();
104
105 nBoundaryFaces += fCells.size();
106
107 UIndirectList<vector>(nHat, fCells) = pp.faceNormals();
108 UIndirectList<scalar>(magSf, fCells) = mag(pp.faceAreas());
109 }
110 nHat.correctBoundaryConditions();
111 magSf.correctBoundaryConditions();
112
113 if (nBoundaryFaces != regionMesh().nCells())
114 {
116 << "Number of primary region coupled boundary faces not equal to "
117 << "the number of cells in the local region" << nl << nl
118 << "Number of cells = " << regionMesh().nCells() << nl
119 << "Boundary faces = " << nBoundaryFaces << nl
120 << abort(FatalError);
121 }
122
123 scalarField passiveMagSf(magSf.size(), Zero);
124 passivePatchIDs_.setSize(intCoupledPatchIDs_.size(), -1);
125 forAll(intCoupledPatchIDs_, i)
126 {
127 const label patchi = intCoupledPatchIDs_[i];
128 const polyPatch& ppIntCoupled = rbm[patchi];
129 if (ppIntCoupled.size() > 0)
130 {
131 label cellId = rbm[patchi].faceCells()[0];
132 const cell& cFaces = regionMesh().cells()[cellId];
133
134 label facei = ppIntCoupled.start();
135 label faceO = cFaces.opposingFaceLabel(facei, regionMesh().faces());
136
137 label passivePatchi = rbm.whichPatch(faceO);
138 passivePatchIDs_[i] = passivePatchi;
139 const polyPatch& ppPassive = rbm[passivePatchi];
140 UIndirectList<scalar>(passiveMagSf, ppPassive.faceCells()) =
141 mag(ppPassive.faceAreas());
142 }
143 }
144
145 Pstream::listCombineAllGather(passivePatchIDs_, maxEqOp<label>());
146
147 magSf.field() = 0.5*(magSf + passiveMagSf);
148 magSf.correctBoundaryConditions();
149}
150
151
152// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
153
155{
156 return regionModel::read();
157}
158
159
160// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
161
163(
164 const fvMesh& mesh,
165 const word& regionType
166)
167:
168 regionModel(mesh, regionType),
169 nHatPtr_(nullptr),
170 magSfPtr_(nullptr),
171 passivePatchIDs_()
172{}
173
174
176(
177 const fvMesh& mesh,
178 const word& regionType,
179 const word& modelName,
180 bool readFields
181)
182:
183 regionModel(mesh, regionType, modelName, false),
184 nHatPtr_(nullptr),
185 magSfPtr_(nullptr),
186 passivePatchIDs_()
187{
188 if (active_)
189 {
190 constructMeshObjects();
191 initialise();
192
193 if (readFields)
194 {
195 read();
196 }
197 }
198}
199
200
201// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
202
204{}
205
206
207// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
208
210{
211 if (!nHatPtr_)
212 {
214 << "Region patch normal vectors not available"
215 << abort(FatalError);
216 }
217
218 return *nHatPtr_;
219}
220
221
223{
224 if (!magSfPtr_)
225 {
227 << "Region patch areas not available"
228 << abort(FatalError);
229 }
230
231 return *magSfPtr_;
232}
233
234
235const Foam::labelList&
237{
238 return passivePatchIDs_;
239}
240
241
242// ************************************************************************* //
static void listCombineAllGather(const List< commsStruct > &comms, List< T > &values, const CombineOp &cop, const int tag, const label comm)
After completion all processors have the same data.
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Base class for region models.
Definition: regionModel.H:63
Switch active_
Active flag.
Definition: regionModel.H:93
const Time & time_
Reference to the time database.
Definition: regionModel.H:90
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:64
virtual bool read()
Read control parameters from dictionary.
Definition: regionModel.C:149
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].
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:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
label cellId
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
const dimensionSet dimless
Dimensionless.
List< label > labelList
A List of labels.
Definition: List.H:66
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
static const char *const typeName
The type name used in ensight case files.