1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | cfMesh: A library for mesh generation
5 \\ / A nd | Copyright held by the original author
7 -------------------------------------------------------------------------------
9 This file is part of cfMesh.
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.
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
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/>.
28 Efficient cell-centre calculation using face-addressing, face-centres and
31 \*---------------------------------------------------------------------------*/
33 #include "polyMeshGenAddressing.H"
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 void polyMeshGenAddressing::calcCellCentresAndVols() const
48 if( cellCentresPtr_ || cellVolumesPtr_ )
50 FatalErrorIn("polyMeshGenAddressing::calcCellCentresAndVols() const")
51 << "Cell centres or cell volumes already calculated"
55 const cellListPMG& cells = mesh_.cells();
57 // set the accumulated cell centre to zero vector
58 cellCentresPtr_ = new vectorField(cells.size());
59 vectorField& cellCtrs = *cellCentresPtr_;
61 // Initialise cell volumes to 0
62 cellVolumesPtr_ = new scalarField(cells.size());
63 scalarField& cellVols = *cellVolumesPtr_;
65 // Make centres and volumes
66 makeCellCentresAndVols(faceCentres(), faceAreas(), cellCtrs, cellVols);
69 void polyMeshGenAddressing::makeCellCentresAndVols
71 const vectorField& fCtrs,
72 const vectorField& fAreas,
73 vectorField& cellCtrs,
77 const labelList& own = mesh_.owner();
78 const cellListPMG& cells = mesh_.cells();
79 const label nCells = cells.size();
82 # pragma omp parallel for if( nCells > 1000 )
84 for(label cellI=0;cellI<nCells;++cellI)
86 const cell& c = cells[cellI];
88 //- approximate the centre first
89 vector cEst(vector::zero);
95 //- start evaluating the volume and the cell centre
96 vector cellCentre(vector::zero);
101 // Calculate 3*face-pyramid volume
102 scalar pyr3Vol = (fAreas[c[fI]] & (fCtrs[c[fI]] - cEst));
104 if( own[c[fI]] != cellI )
107 pyr3Vol = Foam::max(pyr3Vol, VSMALL);
109 // Calculate face-pyramid centre
110 const vector pc = (3.0/4.0)*fCtrs[c[fI]] + (1.0/4.0)*cEst;
112 // Accumulate volume-weighted face-pyramid centre
113 cellCentre += pyr3Vol*pc;
115 // Accumulate face-pyramid volume
119 cellCtrs[cellI] = cellCentre / cellVol;
120 cellVols[cellI] = cellVol / 3.0;
124 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
126 const vectorField& polyMeshGenAddressing::cellCentres() const
128 if( !cellCentresPtr_ )
131 if( omp_in_parallel() )
134 "const vectorField& polyMeshGenAddressing::cellCentres() const"
135 ) << "Calculating addressing inside a parallel region."
136 << " This is not thread safe" << exit(FatalError);
139 calcCellCentresAndVols();
142 return *cellCentresPtr_;
145 const scalarField& polyMeshGenAddressing::cellVolumes() const
147 if( !cellVolumesPtr_ )
150 if( omp_in_parallel() )
153 "const scalarField& polyMeshGenAddressing::cellVolumes() const"
154 ) << "Calculating addressing inside a parallel region."
155 << " This is not thread safe" << exit(FatalError);
158 calcCellCentresAndVols();
161 return *cellVolumesPtr_;
164 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166 } // End namespace Foam
168 // ************************************************************************* //