Fixed URL for libccmio-2.6.1 (bug report #5 by Thomas Oliveira)
[foam-extend-3.2.git] / src / mesh / cfMesh / utilities / meshes / polyMeshGenAddressing / polyMeshGenAddressingCentresAndVols.C
blob7f08cd751dbdb7aad4ed026a83b6bb87a016af0b
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.
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
19     for more details.
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/>.
24 Class
25     polyMeshGenAddressing
27 Description
28     Efficient cell-centre calculation using face-addressing, face-centres and
29     face-areas.
31 \*---------------------------------------------------------------------------*/
33 #include "polyMeshGenAddressing.H"
35 # ifdef USE_OMP
36 #include <omp.h>
37 # endif
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 namespace Foam
44 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
46 void polyMeshGenAddressing::calcCellCentresAndVols() const
48     if( cellCentresPtr_ || cellVolumesPtr_ )
49     {
50         FatalErrorIn("polyMeshGenAddressing::calcCellCentresAndVols() const")
51             << "Cell centres or cell volumes already calculated"
52             << abort(FatalError);
53     }
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,
74     scalarField& cellVols
75 ) const
77     const labelList& own = mesh_.owner();
78     const cellListPMG& cells = mesh_.cells();
79     const label nCells = cells.size();
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];
88         //- approximate the centre first
89         vector cEst(vector::zero);
90         forAll(c, fI)
91             cEst += fCtrs[c[fI]];
93         cEst /= c.size();
95         //- start evaluating the volume and the cell centre
96         vector cellCentre(vector::zero);
97         scalar cellVol(0.0);
99         forAll(c, fI)
100         {
101             // Calculate 3*face-pyramid volume
102             scalar pyr3Vol = (fAreas[c[fI]] & (fCtrs[c[fI]] - cEst));
104             if( own[c[fI]] != cellI )
105                 pyr3Vol *= -1.0;
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
116             cellVol += pyr3Vol;
117         }
119         cellCtrs[cellI] = cellCentre / cellVol;
120         cellVols[cellI] = cellVol / 3.0;
121     }
124 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
126 const vectorField& polyMeshGenAddressing::cellCentres() const
128     if( !cellCentresPtr_ )
129     {
130         # ifdef USE_OMP
131         if( omp_in_parallel() )
132             FatalErrorIn
133             (
134                 "const vectorField& polyMeshGenAddressing::cellCentres() const"
135             ) << "Calculating addressing inside a parallel region."
136                 << " This is not thread safe" << exit(FatalError);
137         # endif
139         calcCellCentresAndVols();
140     }
142     return *cellCentresPtr_;
145 const scalarField& polyMeshGenAddressing::cellVolumes() const
147     if( !cellVolumesPtr_ )
148     {
149         # ifdef USE_OMP
150         if( omp_in_parallel() )
151             FatalErrorIn
152             (
153                 "const scalarField& polyMeshGenAddressing::cellVolumes() const"
154             ) << "Calculating addressing inside a parallel region."
155                 << " This is not thread safe" << exit(FatalError);
156         # endif
158         calcCellCentresAndVols();
159     }
161     return *cellVolumesPtr_;
164 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166 } // End namespace Foam
168 // ************************************************************************* //