singleLayerRegion.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "singleLayerRegion.H"
27 #include "fvMesh.H"
28 #include "surfaceFields.H"
29 #include "Time.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace regionModels
37 {
38  defineTypeNameAndDebug(singleLayerRegion, 0);
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
46 (
47  const fvMesh& mesh,
48  const word& regionType,
49  const word& modelName,
50  bool readFields
51 )
52 :
53  regionModel(mesh, regionType, modelName, false),
54  nHat_
55  (
56  IOobject
57  (
58  "nHat",
59  time_.timeName(),
60  regionMesh()
61  ),
62  regionMesh(),
65  ),
66  magSf_
67  (
68  IOobject
69  (
70  "magSf",
71  time_.timeName(),
72  regionMesh()
73  ),
74  regionMesh(),
76  ),
77  VbyA_
78  (
79  IOobject
80  (
81  "VbyA",
82  time_.timeName(),
83  regionMesh()
84  ),
85  regionMesh(),
88  ),
89  passivePatchIDs_()
90 {
91  label nBoundaryFaces = 0;
92  const polyBoundaryMesh& rbm = regionMesh().boundaryMesh();
93  forAll(intCoupledPatchIDs_, i)
94  {
95  const label patchi = intCoupledPatchIDs_[i];
96  const polyPatch& pp = rbm[patchi];
97  const labelList& fCells = pp.faceCells();
98 
99  nBoundaryFaces += fCells.size();
100 
101  UIndirectList<vector>(nHat_, fCells) = pp.faceNormals();
102  UIndirectList<scalar>(magSf_, fCells) = pp.magFaceAreas();
103  }
104  nHat_.correctBoundaryConditions();
105 
106  if (nBoundaryFaces != regionMesh().nCells())
107  {
109  << "Number of primary region coupled boundary faces not equal to "
110  << "the number of cells in the local region" << nl << nl
111  << "Number of cells = " << regionMesh().nCells() << nl
112  << "Boundary faces = " << nBoundaryFaces << nl
113  << abort(FatalError);
114  }
115 
116  passivePatchIDs_.setSize(intCoupledPatchIDs_.size(), -1);
117  forAll(intCoupledPatchIDs_, i)
118  {
119  const label patchi = intCoupledPatchIDs_[i];
120  const polyPatch& ppIntCoupled = rbm[patchi];
121  if (ppIntCoupled.size() > 0)
122  {
123  const label cellId = rbm[patchi].faceCells()[0];
124  const cell& cFaces = regionMesh().cells()[cellId];
125  const label faceO
126  (
127  cFaces.opposingFaceLabel
128  (
129  ppIntCoupled.start(), regionMesh().faces()
130  )
131  );
132  passivePatchIDs_[i] = rbm.whichPatch(faceO);
133  }
134  }
135 
136  Pstream::listCombineGather(passivePatchIDs_, maxEqOp<label>());
137  Pstream::listCombineScatter(passivePatchIDs_);
138 
139  VbyA_.primitiveFieldRef() = regionMesh().V()/magSf_;
140  VbyA_.correctBoundaryConditions();
141 }
142 
143 
144 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
145 
147 {}
148 
149 
150 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
151 
153 {
154  return nHat_;
155 }
156 
157 
160 {
161  return magSf_;
162 }
163 
164 
166 {
167  return VbyA_;
168 }
169 
170 
171 const Foam::labelList&
173 {
174  return passivePatchIDs_;
175 }
176 
177 
178 // ************************************************************************* //
Foam::surfaceFields.
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const volScalarField & VbyA() const
Return the cell layer volume/area [m].
const dimensionSet dimArea
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
const labelList & passivePatchIDs() const
Return the list of patch IDs opposite to internally.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
const dimensionSet dimless
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
const dimensionSet dimLength
defineTypeNameAndDebug(regionModel, 0)
label opposingFaceLabel(const label masterFaceLabel, const faceUList &meshFaces) const
Return index of opposite face.
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:340
A class for handling words, derived from string.
Definition: word.H:59
const volVectorField & nHat() const
Return the patch normal vectors.
static const zero Zero
Definition: zero.H:97
const Field< PointType > & faceNormals() const
Return face normals for patch.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::polyBoundaryMesh.
static const char nl
Definition: Ostream.H:260
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
label patchi
const scalarField::subField magFaceAreas() const
Return face area magnitudes.
Definition: polyPatch.C:315
const volScalarField::Internal & magSf() const
Return the face area magnitudes [m^2].
singleLayerRegion(const fvMesh &mesh, const word &regionType, const word &modelName, bool readFields=true)
Construct from mesh, region type and name.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
A List with indirect addressing.
Definition: fvMatrix.H:106
label cellId
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:306
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Base class for region models.
Definition: regionModel.H:55
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Namespace for OpenFOAM.