polyMeshGenAddressingCentresAndVols.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | cfMesh: A library for mesh generation
4  \\ / O peration |
5  \\ / A nd | Copyright held by the original author
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of cfMesh.
10 
11  cfMesh is free software; you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by the
13  Free Software Foundation; either version 3 of the License, or (at your
14  option) any later version.
15 
16  cfMesh 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 cfMesh. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  polyMeshGenAddressing
26 
27 Description
28  Efficient cell-centre calculation using face-addressing, face-centres and
29  face-areas.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "polyMeshGenAddressing.H"
34 
35 # ifdef USE_OMP
36 #include <omp.h>
37 # endif
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
47 {
49  {
50  FatalErrorIn("polyMeshGenAddressing::calcCellCentresAndVols() const")
51  << "Cell centres or cell volumes already calculated"
52  << abort(FatalError);
53  }
54 
55  const cellListPMG& cells = mesh_.cells();
56 
57  // set the accumulated cell centre to zero vector
58  cellCentresPtr_ = new vectorField(cells.size());
59  vectorField& cellCtrs = *cellCentresPtr_;
60 
61  // Initialise cell volumes to 0
62  cellVolumesPtr_ = new scalarField(cells.size());
64 
65  // Make centres and volumes
67 }
68 
70 (
71  const vectorField& fCtrs,
72  const vectorField& fAreas,
73  vectorField& cellCtrs,
75 ) const
76 {
77  const labelList& own = mesh_.owner();
78  const cellListPMG& cells = mesh_.cells();
79  const label nCells = cells.size();
80 
81  # ifdef USE_OMP
82  # pragma omp parallel for if( nCells > 1000 )
83  # endif
84  for(label cellI=0;cellI<nCells;++cellI)
85  {
86  const cell& c = cells[cellI];
87 
88  //- approximate the centre first
89  vector cEst(vector::zero);
90  forAll(c, fI)
91  cEst += fCtrs[c[fI]];
92 
93  cEst /= c.size();
94 
95  //- start evaluating the volume and the cell centre
96  vector cellCentre(vector::zero);
97  scalar cellVol(0.0);
98 
99  forAll(c, fI)
100  {
101  // Calculate 3*face-pyramid volume
102  scalar pyr3Vol = (fAreas[c[fI]] & (fCtrs[c[fI]] - cEst));
103 
104  if( own[c[fI]] != cellI )
105  pyr3Vol *= -1.0;
106 
107  pyr3Vol = Foam::max(pyr3Vol, VSMALL);
108 
109  // Calculate face-pyramid centre
110  const vector pc = (3.0/4.0)*fCtrs[c[fI]] + (1.0/4.0)*cEst;
111 
112  // Accumulate volume-weighted face-pyramid centre
113  cellCentre += pyr3Vol*pc;
114 
115  // Accumulate face-pyramid volume
116  cellVol += pyr3Vol;
117  }
118 
119  cellCtrs[cellI] = cellCentre / cellVol;
120  cellVols[cellI] = cellVol / 3.0;
121  }
122 }
123 
124 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
125 
127 {
128  if( !cellCentresPtr_ )
129  {
130  # ifdef USE_OMP
131  if( omp_in_parallel() )
133  (
134  "const vectorField& polyMeshGenAddressing::cellCentres() const"
135  ) << "Calculating addressing inside a parallel region."
136  << " This is not thread safe" << exit(FatalError);
137  # endif
138 
140  }
141 
142  return *cellCentresPtr_;
143 }
144 
146 {
147  if( !cellVolumesPtr_ )
148  {
149  # ifdef USE_OMP
150  if( omp_in_parallel() )
152  (
153  "const scalarField& polyMeshGenAddressing::cellVolumes() const"
154  ) << "Calculating addressing inside a parallel region."
155  << " This is not thread safe" << exit(FatalError);
156  # endif
157 
159  }
160 
161  return *cellVolumesPtr_;
162 }
163 
164 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
165 
166 } // End namespace Foam
167 
168 // ************************************************************************* //
Foam::polyMeshGenFaces::owner
const labelList & owner() const
owner and neighbour cells for faces
Definition: polyMeshGenFacesI.H:67
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
Foam::polyMeshGenAddressing::cellCentresPtr_
vectorField * cellCentresPtr_
Cell centres.
Definition: polyMeshGenAddressing.H:118
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
Foam::polyMeshGenAddressing::faceCentres
const vectorField & faceCentres() const
Definition: polyMeshGenAddressingCentresAndAreas.C:123
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::polyMeshGenAddressing::cellCentres
const vectorField & cellCentres() const
Definition: polyMeshGenAddressingCentresAndVols.C:125
Foam::polyMeshGenAddressing::mesh_
const polyMeshGenCells & mesh_
reference to the mesh
Definition: polyMeshGenAddressing.H:74
Foam::polyMeshGenAddressing::cellVolumesPtr_
scalarField * cellVolumesPtr_
Cell volumes.
Definition: polyMeshGenAddressing.H:124
Foam::cellListPMG
Definition: cellListPMG.H:49
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:49
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:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::polyMeshGenCells::cells
const cellListPMG & cells() const
access to cells
Definition: polyMeshGenCellsI.H:39
Foam::polyMeshGenAddressing::calcCellCentresAndVols
void calcCellCentresAndVols() const
Calculate cell centres and volumes.
Definition: polyMeshGenAddressingCentresAndVols.C:45
cellVols
const scalarField & cellVols
Definition: temperatureAndPressureVariables.H:49
Foam::FatalError
error FatalError
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::polyMeshGenAddressing::faceAreas
const vectorField & faceAreas() const
Definition: polyMeshGenAddressingCentresAndAreas.C:142
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::polyMeshGenAddressing::makeCellCentresAndVols
void makeCellCentresAndVols(const vectorField &fCtrs, const vectorField &fAreas, vectorField &cellCtrs, scalarField &cellVols) const
Definition: polyMeshGenAddressingCentresAndVols.C:69
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
FatalErrorIn
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:313
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::polyMeshGenAddressing::cellVolumes
const scalarField & cellVolumes() const
Definition: polyMeshGenAddressingCentresAndVols.C:144
polyMeshGenAddressing.H
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56