ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / regionModels / regionModel / singleLayerRegion / singleLayerRegion.C
blob87118d9d2e10925cac1bbda2cf1ab06ee0003770
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "singleLayerRegion.H"
27 #include "fvMesh.H"
28 #include "Time.H"
29 #include "directMappedWallPolyPatch.H"
30 #include "zeroGradientFvPatchFields.H"
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 namespace Foam
36 namespace regionModels
38     defineTypeNameAndDebug(singleLayerRegion, 0);
42 // * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * * //
44 void Foam::regionModels::singleLayerRegion::constructMeshObjects()
46     // construct patch normal vectors
47     nHatPtr_.reset
48     (
49         new volVectorField
50         (
51             IOobject
52             (
53                 "nHat",
54                 time_.timeName(),
55                 regionMesh(),
56                 IOobject::READ_IF_PRESENT,
57                 NO_WRITE
58             ),
59             regionMesh(),
60             dimensionedVector("zero", dimless, vector::zero),
61             zeroGradientFvPatchField<vector>::typeName
62         )
63     );
65     // construct patch areas
66     magSfPtr_.reset
67     (
68         new volScalarField
69         (
70             IOobject
71             (
72                 "magSf",
73                 time_.timeName(),
74                 regionMesh(),
75                 IOobject::READ_IF_PRESENT,
76                 NO_WRITE
77             ),
78             regionMesh(),
79             dimensionedScalar("zero", dimArea, 0.0),
80             zeroGradientFvPatchField<scalar>::typeName
81         )
82     );
86 void Foam::regionModels::singleLayerRegion::initialise()
88     if (debug)
89     {
90         Pout<< "singleLayerRegion::initialise()" << endl;
91     }
93     label nBoundaryFaces = 0;
94     const polyBoundaryMesh& rbm = regionMesh().boundaryMesh();
95     volVectorField& nHat = nHatPtr_();
96     volScalarField& magSf = magSfPtr_();
97     forAll(intCoupledPatchIDs_, i)
98     {
99         const label patchI = intCoupledPatchIDs_[i];
100         const polyPatch& pp = rbm[patchI];
101         const labelList& fCells = pp.faceCells();
103         nBoundaryFaces += fCells.size();
105         UIndirectList<vector>(nHat, fCells) = pp.faceNormals();
106         UIndirectList<scalar>(magSf, fCells) = mag(pp.faceAreas());
107     }
108     nHat.correctBoundaryConditions();
109     magSf.correctBoundaryConditions();
111     if (nBoundaryFaces != regionMesh().nCells())
112     {
113         FatalErrorIn("singleLayerRegion::initialise()")
114             << "Number of primary region coupled boundary faces not equal to "
115             << "the number of cells in the local region" << nl << nl
116             << "Number of cells = " << regionMesh().nCells() << nl
117             << "Boundary faces  = " << nBoundaryFaces << nl
118             << abort(FatalError);
119     }
121     scalarField passiveMagSf(magSf.size(), 0.0);
122     passivePatchIDs_.setSize(intCoupledPatchIDs_.size(), -1);
123     forAll(intCoupledPatchIDs_, i)
124     {
125         const label patchI = intCoupledPatchIDs_[i];
126         const polyPatch& ppIntCoupled = rbm[patchI];
127         if (ppIntCoupled.size() > 0)
128         {
129             label cellId = rbm[patchI].faceCells()[0];
130             const cell& cFaces = regionMesh().cells()[cellId];
132             label faceI = ppIntCoupled.start();
133             label faceO = cFaces.opposingFaceLabel(faceI, regionMesh().faces());
135             label passivePatchI = rbm.whichPatch(faceO);
136             passivePatchIDs_[i] = passivePatchI;
137             const polyPatch& ppPassive = rbm[passivePatchI];
138             UIndirectList<scalar>(passiveMagSf, ppPassive.faceCells()) =
139                 mag(ppPassive.faceAreas());
140         }
141     }
143     Pstream::listCombineGather(passivePatchIDs_, maxEqOp<label>());
144     Pstream::listCombineScatter(passivePatchIDs_);
146     magSf.field() = 0.5*(magSf + passiveMagSf);
147     magSf.correctBoundaryConditions();
151 // * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //
153 bool Foam::regionModels::singleLayerRegion::read()
155     return regionModel::read();
159 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
161 Foam::regionModels::singleLayerRegion::singleLayerRegion(const fvMesh& mesh)
163     regionModel(mesh),
164     nHatPtr_(NULL),
165     magSfPtr_(NULL),
166     passivePatchIDs_()
170 Foam::regionModels::singleLayerRegion::singleLayerRegion
172     const fvMesh& mesh,
173     const word& regionType,
174     const word& modelName,
175     bool readFields
178     regionModel(mesh, regionType, modelName, false),
179     nHatPtr_(NULL),
180     magSfPtr_(NULL),
181     passivePatchIDs_()
183     if (active_)
184     {
185         constructMeshObjects();
186         initialise();
188         if (readFields)
189         {
190             read();
191         }
192     }
196 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
198 Foam::regionModels::singleLayerRegion::~singleLayerRegion()
202 // * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * * //
204 const Foam::volVectorField& Foam::regionModels::singleLayerRegion::nHat() const
206     if (!nHatPtr_.valid())
207     {
208         FatalErrorIn("const fvMesh& singleLayerRegion::nHat() const")
209             << "Region patch normal vectors not available"
210             << abort(FatalError);
211     }
213     return nHatPtr_();
217 const Foam::volScalarField& Foam::regionModels::singleLayerRegion::magSf() const
219     if (!magSfPtr_.valid())
220     {
221         FatalErrorIn("const fvMesh& singleLayerRegion::magSf() const")
222             << "Region patch areas not available"
223             << abort(FatalError);
224     }
226     return magSfPtr_();
230 const Foam::labelList&
231 Foam::regionModels::singleLayerRegion::passivePatchIDs() const
233     return passivePatchIDs_;
237 // ************************************************************************* //