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 -------------------------------------------------------------------------------
11 License
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 
36 namespace Foam
37 {
38 namespace regionModels
39 {
40  defineTypeNameAndDebug(singleLayerRegion, 0);
41 }
42 }
43 
44 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
45 
46 void Foam::regionModels::singleLayerRegion::constructMeshObjects()
47 {
48  // construct patch normal vectors
49  nHatPtr_.reset
50  (
51  new volVectorField
52  (
53  IOobject
54  (
55  "nHat",
56  time_.timeName(),
57  regionMesh(),
61  ),
62  regionMesh(),
65  )
66  );
67 
68  // construct patch areas
69  magSfPtr_.reset
70  (
71  new volScalarField
72  (
73  IOobject
74  (
75  "magSf",
76  time_.timeName(),
77  regionMesh(),
81  ),
82  regionMesh(),
85  )
86  );
87 }
88 
89 
90 void 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 labelList& 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::listCombineReduce(passivePatchIDs_, maxEqOp<label>());
148 
149  magSf.field() = 0.5*(magSf + passiveMagSf);
150  magSf.correctBoundaryConditions();
151 }
152 
153 
154 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
155 
157 {
159 }
160 
161 
162 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
163 
164 Foam::regionModels::singleLayerRegion::singleLayerRegion
165 (
166  const fvMesh& mesh,
167  const word& regionType
168 )
169 :
170  regionModel(mesh, regionType),
171  nHatPtr_(nullptr),
172  magSfPtr_(nullptr),
173  passivePatchIDs_()
174 {}
175 
176 
177 Foam::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();
198  }
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);
231  }
232 
233  return *magSfPtr_;
234 }
235 
236 
237 const Foam::labelList&
239 {
240  return passivePatchIDs_;
241 }
242 
243 
244 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
virtual const labelList & passivePatchIDs() const
Return the list of patch IDs opposite to internally.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
void setSize(const label n)
Alias for resize()
Definition: List.H:320
dynamicFvMesh & mesh
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.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual bool read()
Read control parameters from dictionary.
Definition: regionModel.C:145
Reading is optional [identical to LAZY_READ].
virtual const volVectorField & nHat() const
Return the patch normal vectors.
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:26
errorManip< error > abort(error &err)
Definition: errorManip.H:139
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
Switch active_
Active flag.
Definition: regionModel.H:102
int debug
Static debugging option.
autoPtr< volScalarField > magSfPtr_
Face area magnitudes / [m2].
virtual const volScalarField & magSf() const
Return the face area magnitudes / [m2].
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Time & time_
Reference to the time database.
Definition: regionModel.H:97
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
label cellId
autoPtr< volVectorField > nHatPtr_
Patch normal vectors.
Base class for region models.
Definition: regionModel.H:56
List< label > labelList
A List of labels.
Definition: List.H:62
defineTypeNameAndDebug(KirchhoffShell, 0)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Do not request registration (bool: false)
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
static void listCombineReduce(List< T > &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Combines List elements. After completion all processors have the same data.
virtual bool read()
Read control parameters from dictionary.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127