primitiveMeshCellCentresAndVols.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 Description
27  Efficient cell-centre calculation using face-addressing, face-centres and
28  face-areas.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "primitiveMesh.H"
33 #include "PrecisionAdaptor.H"
34 
35 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 
38 {
39  if (debug)
40  {
41  Pout<< "primitiveMesh::calcCellCentresAndVols() : "
42  << "Calculating cell centres and cell volumes"
43  << endl;
44  }
45 
46  // It is an error to attempt to recalculate cellCentres
47  // if the pointer is already set
48  if (cellCentresPtr_ || cellVolumesPtr_)
49  {
51  << "Cell centres or cell volumes already calculated"
52  << abort(FatalError);
53  }
54 
55  // set the accumulated cell centre to zero vector
56  cellCentresPtr_ = new vectorField(nCells());
57  vectorField& cellCtrs = *cellCentresPtr_;
58 
59  // Initialise cell volumes to 0
60  cellVolumesPtr_ = new scalarField(nCells());
61  scalarField& cellVols = *cellVolumesPtr_;
62 
63  // Make centres and volumes
65 
66  if (debug)
67  {
68  Pout<< "primitiveMesh::calcCellCentresAndVols() : "
69  << "Finished calculating cell centres and cell volumes"
70  << endl;
71  }
72 }
73 
74 
76 (
77  const vectorField& fCtrs,
78  const vectorField& fAreas,
79  vectorField& cellCtrs_s,
80  scalarField& cellVols_s
81 ) const
82 {
83  typedef Vector<solveScalar> solveVector;
84 
85  PrecisionAdaptor<solveVector, vector> tcellCtrs(cellCtrs_s);
86  Field<solveVector>& cellCtrs = tcellCtrs.ref();
87  PrecisionAdaptor<solveScalar, scalar> tcellVols(cellVols_s);
88  Field<solveScalar>& cellVols = tcellVols.ref();
89 
90  // Clear the fields for accumulation
91  cellCtrs = Zero;
92  cellVols = 0.0;
93 
94  const labelList& own = faceOwner();
95  const labelList& nei = faceNeighbour();
96 
97  // first estimate the approximate cell centre as the average of
98  // face centres
99 
100  Field<solveVector> cEst(nCells(), Zero);
101  labelField nCellFaces(nCells(), Zero);
102 
103  forAll(own, facei)
104  {
105  cEst[own[facei]] += solveVector(fCtrs[facei]);
106  ++nCellFaces[own[facei]];
107  }
108 
109  forAll(nei, facei)
110  {
111  cEst[nei[facei]] += solveVector(fCtrs[facei]);
112  ++nCellFaces[nei[facei]];
113  }
114 
115  forAll(cEst, celli)
116  {
117  cEst[celli] /= nCellFaces[celli];
118  }
119 
120  forAll(own, facei)
121  {
122  const solveVector fc(fCtrs[facei]);
123  const solveVector fA(fAreas[facei]);
124 
125  // Calculate 3*face-pyramid volume
126  solveScalar pyr3Vol =
127  fA & (fc - cEst[own[facei]]);
128 
129  // Calculate face-pyramid centre
130  solveVector pc = (3.0/4.0)*fc + (1.0/4.0)*cEst[own[facei]];
131 
132  // Accumulate volume-weighted face-pyramid centre
133  cellCtrs[own[facei]] += pyr3Vol*pc;
134 
135  // Accumulate face-pyramid volume
136  cellVols[own[facei]] += pyr3Vol;
137  }
138 
139  forAll(nei, facei)
140  {
141  const solveVector fc(fCtrs[facei]);
142  const solveVector fA(fAreas[facei]);
143 
144  // Calculate 3*face-pyramid volume
145  solveScalar pyr3Vol =
146  fA & (cEst[nei[facei]] - fc);
147 
148  // Calculate face-pyramid centre
149  solveVector pc = (3.0/4.0)*fc + (1.0/4.0)*cEst[nei[facei]];
150 
151  // Accumulate volume-weighted face-pyramid centre
152  cellCtrs[nei[facei]] += pyr3Vol*pc;
153 
154  // Accumulate face-pyramid volume
155  cellVols[nei[facei]] += pyr3Vol;
156  }
157 
158  forAll(cellCtrs, celli)
159  {
160  if (mag(cellVols[celli]) > VSMALL)
161  {
162  cellCtrs[celli] /= cellVols[celli];
163  }
164  else
165  {
166  cellCtrs[celli] = cEst[celli];
167  }
168  }
169 
170  cellVols *= (1.0/3.0);
171 }
172 
173 
174 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
175 
177 {
178  if (!cellCentresPtr_)
179  {
180  calcCellCentresAndVols();
181  }
182 
183  return *cellCentresPtr_;
184 }
185 
186 
188 {
189  if (!cellVolumesPtr_)
190  {
191  calcCellCentresAndVols();
192  }
193 
194  return *cellVolumesPtr_;
195 }
196 
197 
198 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
PrecisionAdaptor.H
primitiveMesh.H
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::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::primitiveMesh::nCells
label nCells() const
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::Field< vector >
Foam::primitiveMesh::calcCellCentresAndVols
void calcCellCentresAndVols() const
Calculate cell centres and volumes.
Definition: primitiveMeshCellCentresAndVols.C:37
cellVols
const scalarField & cellVols
Definition: temperatureAndPressureVariables.H:51
Foam::FatalError
error FatalError
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::primitiveMesh::cellVolumes
const scalarField & cellVolumes() const
Definition: primitiveMeshCellCentresAndVols.C:187
Foam::PrecisionAdaptor
A Field wrapper with possible data conversion.
Definition: PrecisionAdaptor.H:158
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:176
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::primitiveMesh::makeCellCentresAndVols
void makeCellCentresAndVols(const vectorField &fCtrs, const vectorField &fAreas, vectorField &cellCtrs, scalarField &cellVols) const
Definition: primitiveMeshCellCentresAndVols.C:76
Foam::PrecisionAdaptor::ref
FieldType & ref()
Allow modification without const-ref check.
Definition: PrecisionAdaptor.H:219
Foam::Vector
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition: Vector.H:62
Foam::List< label >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::primitiveMesh::faceCentres
const vectorField & faceCentres() const
Definition: primitiveMeshFaceCentresAndAreas.C:145
Foam::primitiveMesh::faceAreas
const vectorField & faceAreas() const
Definition: primitiveMeshFaceCentresAndAreas.C:156