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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "singleLayerRegion.H"
29 #include "fvMesh.H"
30 #include "Time.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace regionModels
38 {
39  defineTypeNameAndDebug(singleLayerRegion, 0);
40 }
41 }
42 
43 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
44 
45 void Foam::regionModels::singleLayerRegion::constructMeshObjects()
46 {
47  // construct patch normal vectors
48  nHatPtr_.reset
49  (
50  new volVectorField
51  (
52  IOobject
53  (
54  "nHat",
55  time_.timeName(),
56  regionMesh(),
58  NO_WRITE
59  ),
60  regionMesh(),
62  zeroGradientFvPatchField<vector>::typeName
63  )
64  );
65 
66  // construct patch areas
67  magSfPtr_.reset
68  (
69  new volScalarField
70  (
71  IOobject
72  (
73  "magSf",
74  time_.timeName(),
75  regionMesh(),
77  NO_WRITE
78  ),
79  regionMesh(),
81  zeroGradientFvPatchField<scalar>::typeName
82  )
83  );
84 }
85 
86 
87 void Foam::regionModels::singleLayerRegion::initialise()
88 {
89  if (debug)
90  {
91  Pout<< "singleLayerRegion::initialise()" << endl;
92  }
93 
94  label nBoundaryFaces = 0;
95  const polyBoundaryMesh& rbm = regionMesh().boundaryMesh();
96  volVectorField& nHat = nHatPtr_();
97  volScalarField& magSf = magSfPtr_();
98  forAll(intCoupledPatchIDs_, i)
99  {
100  const label patchi = intCoupledPatchIDs_[i];
101  const polyPatch& pp = rbm[patchi];
102  const labelList& fCells = pp.faceCells();
103 
104  nBoundaryFaces += fCells.size();
105 
106  UIndirectList<vector>(nHat, fCells) = pp.faceNormals();
107  UIndirectList<scalar>(magSf, fCells) = mag(pp.faceAreas());
108  }
109  nHat.correctBoundaryConditions();
110  magSf.correctBoundaryConditions();
111 
112  if (nBoundaryFaces != regionMesh().nCells())
113  {
115  << "Number of primary region coupled boundary faces not equal to "
116  << "the number of cells in the local region" << nl << nl
117  << "Number of cells = " << regionMesh().nCells() << nl
118  << "Boundary faces = " << nBoundaryFaces << nl
119  << abort(FatalError);
120  }
121 
122  scalarField passiveMagSf(magSf.size(), Zero);
123  passivePatchIDs_.setSize(intCoupledPatchIDs_.size(), -1);
124  forAll(intCoupledPatchIDs_, i)
125  {
126  const label patchi = intCoupledPatchIDs_[i];
127  const polyPatch& ppIntCoupled = rbm[patchi];
128  if (ppIntCoupled.size() > 0)
129  {
130  label cellId = rbm[patchi].faceCells()[0];
131  const cell& cFaces = regionMesh().cells()[cellId];
132 
133  label facei = ppIntCoupled.start();
134  label faceO = cFaces.opposingFaceLabel(facei, regionMesh().faces());
135 
136  label passivePatchi = rbm.whichPatch(faceO);
137  passivePatchIDs_[i] = passivePatchi;
138  const polyPatch& ppPassive = rbm[passivePatchi];
139  UIndirectList<scalar>(passiveMagSf, ppPassive.faceCells()) =
140  mag(ppPassive.faceAreas());
141  }
142  }
143 
144  Pstream::listCombineGather(passivePatchIDs_, maxEqOp<label>());
145  Pstream::listCombineScatter(passivePatchIDs_);
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 
162 Foam::regionModels::singleLayerRegion::singleLayerRegion
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 
175 Foam::regionModels::singleLayerRegion::singleLayerRegion
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_.valid())
212  {
214  << "Region patch normal vectors not available"
215  << abort(FatalError);
216  }
217 
218  return *nHatPtr_;
219 }
220 
221 
223 {
224  if (!magSfPtr_.valid())
225  {
227  << "Region patch areas not available"
228  << abort(FatalError);
229  }
230 
231  return *magSfPtr_;
232 }
233 
234 
235 const Foam::labelList&
237 {
238  return passivePatchIDs_;
239 }
240 
241 
242 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::regionModels::singleLayerRegion::passivePatchIDs
virtual const labelList & passivePatchIDs() const
Return the list of patch IDs opposite to internally.
Definition: singleLayerRegion.C:236
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::regionModels::singleLayerRegion::~singleLayerRegion
virtual ~singleLayerRegion()
Destructor.
Definition: singleLayerRegion.C:203
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:764
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::regionModels::regionModel
Base class for region models.
Definition: regionModel.H:59
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:60
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::regionModels::singleLayerRegion::magSfPtr_
autoPtr< volScalarField > magSfPtr_
Face area magnitudes / [m2].
Definition: singleLayerRegion.H:81
singleLayerRegion.H
Foam::regionModels::singleLayerRegion::magSf
virtual const volScalarField & magSf() const
Return the face area magnitudes / [m2].
Definition: singleLayerRegion.C:222
Foam::regionModels::regionModel::regionMesh
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:63
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:122
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::regionModels::regionModel::read
virtual bool read()
Read control parameters from dictionary.
Definition: regionModel.C:149
Foam::regionModels::defineTypeNameAndDebug
defineTypeNameAndDebug(regionModel, 0)
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
cellId
label cellId
Definition: interrogateWallPatches.H:67
Foam::readFields
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Definition: ReadFieldsTemplates.C:312
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:60
Foam::Pstream::listCombineGather
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:290
Time.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::regionModels::regionModel::time_
const Time & time_
Reference to the time database.
Definition: regionModel.H:89
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::IOobject::IOobject
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:276
Foam::List< label >
Foam::regionModels::singleLayerRegion::nHatPtr_
autoPtr< volVectorField > nHatPtr_
Patch normal vectors.
Definition: singleLayerRegion.H:78
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::regionModels::singleLayerRegion::nHat
virtual const volVectorField & nHat() const
Return the patch normal vectors.
Definition: singleLayerRegion.C:209
Foam::Pstream::listCombineScatter
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:432
Foam::GeometricField< vector, fvPatchField, volMesh >
zeroGradientFvPatchFields.H
Foam::regionModels::singleLayerRegion::read
virtual bool read()
Read control parameters from dictionary.
Definition: singleLayerRegion.C:154