faceAgglomerate.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 | Copyright (C) 2011-2013 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 Application
25  faceAgglomerate
26 
27 Description
28  Agglomerate boundary faces using the pairPatchAgglomeration algorithm.
29 
30  It writes a map from the fine to coarse grid.
31 
32 SeeAlso
33  pairPatchAgglomeration.H
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #include "argList.H"
38 #include "fvMesh.H"
39 #include "Time.H"
40 #include "volFields.H"
41 #include "CompactListList.H"
42 #include "unitConversion.H"
43 #include "pairPatchAgglomeration.H"
44 #include "labelListIOList.H"
45 #include "syncTools.H"
46 #include "globalIndex.H"
47 
48 using namespace Foam;
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 int main(int argc, char *argv[])
53 {
54  #include "addRegionOption.H"
55  #include "addDictOption.H"
56  #include "setRootCase.H"
57  #include "createTime.H"
58  #include "createNamedMesh.H"
59 
60  const word dictName("viewFactorsDict");
61 
63 
64  // Read control dictionary
65  const IOdictionary agglomDict(dictIO);
66 
67  bool writeAgglom = readBool(agglomDict.lookup("writeFacesAgglomeration"));
68 
69 
71 
72  labelListIOList finalAgglom
73  (
74  IOobject
75  (
76  "finalAgglom",
78  mesh,
81  false
82  ),
83  boundary.size()
84  );
85 
86  label nCoarseFaces = 0;
87 
88  forAllConstIter(dictionary, agglomDict, iter)
89  {
90  labelList patchIds = boundary.findIndices(iter().keyword());
91  forAll(patchIds, i)
92  {
93  label patchI = patchIds[i];
94  const polyPatch& pp = boundary[patchI];
95  if (!pp.coupled())
96  {
97  Info << "\nAgglomerating patch : " << pp.name() << endl;
98  pairPatchAgglomeration agglomObject
99  (
100  pp,
101  agglomDict.subDict(pp.name())
102  );
103  agglomObject.agglomerate();
104  finalAgglom[patchI] =
105  agglomObject.restrictTopBottomAddressing();
106 
107  if (finalAgglom[patchI].size())
108  {
109  nCoarseFaces += max(finalAgglom[patchI] + 1);
110  }
111  }
112  }
113  }
114 
115 
116  // - All patches which are not agglomarated are identity for finalAgglom
118  {
119  if (finalAgglom[patchId].size() == 0)
120  {
121  finalAgglom[patchId] = identity(boundary[patchId].size());
122  }
123  }
124 
125  // Sync agglomeration across coupled patches
126  labelList nbrAgglom(mesh.nFaces() - mesh.nInternalFaces(), -1);
127 
129  {
130  const polyPatch& pp = boundary[patchId];
131  if (pp.coupled())
132  {
133  finalAgglom[patchId] = identity(pp.size());
134  forAll(pp, i)
135  {
136  nbrAgglom[pp.start() - mesh.nInternalFaces() + i] =
137  finalAgglom[patchId][i];
138  }
139  }
140  }
141 
144  {
145  const polyPatch& pp = boundary[patchId];
146  if (pp.coupled() && !refCast<const coupledPolyPatch>(pp).owner())
147  {
148  forAll(pp, i)
149  {
150  finalAgglom[patchId][i] =
151  nbrAgglom[pp.start() - mesh.nInternalFaces() + i];
152  }
153  }
154  }
155 
156  finalAgglom.write();
157 
158  if (writeAgglom)
159  {
160  globalIndex index(nCoarseFaces);
161  volScalarField facesAgglomeration
162  (
163  IOobject
164  (
165  "facesAgglomeration",
166  mesh.time().timeName(),
167  mesh,
170  ),
171  mesh,
172  dimensionedScalar("facesAgglomeration", dimless, 0)
173  );
174 
175  label coarsePatchIndex = 0;
177  {
178  const polyPatch& pp = boundary[patchId];
179  if (pp.size() > 0)
180  {
181  fvPatchScalarField& bFacesAgglomeration =
182  facesAgglomeration.boundaryField()[patchId];
183 
184  forAll(bFacesAgglomeration, j)
185  {
186  bFacesAgglomeration[j] =
187  index.toGlobal
188  (
190  finalAgglom[patchId][j] + coarsePatchIndex
191  );
192  }
193 
194  coarsePatchIndex += max(finalAgglom[patchId]) + 1;
195  }
196  }
197 
198  Info<< "\nWriting facesAgglomeration" << endl;
199  facesAgglomeration.write();
200  }
201 
202  Info<< "End\n" << endl;
203  return 0;
204 }
205 
206 
207 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
volFields.H
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
globalIndex.H
Foam::polyPatch::coupled
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:322
boundary
faceListList boundary(nPatches)
Foam::pairPatchAgglomeration
Primitive patch pair agglomerate method.
Definition: pairPatchAgglomeration.H:57
Foam::polyMesh::facesInstance
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:778
unitConversion.H
Unit conversion functions.
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
syncTools.H
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
addDictOption.H
labelListIOList.H
dictName
const word dictName("particleTrackDict")
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::Ostream::write
virtual Ostream & write(const token &)=0
Write next token to stream.
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
argList.H
addRegionOption.H
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Foam::identity
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Definition: primitiveMeshI.H:52
createNamedMesh.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:312
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:405
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
setRootCase.H
Foam::syncTools::swapBoundaryFaceList
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled boundary face values.
Definition: syncTools.H:430
Foam::primitiveMesh::nFaces
label nFaces() const
Definition: primitiveMeshI.H:58
pairPatchAgglomeration.H
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
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::IOList< labelList >
createTime.H
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
patchId
label patchId(-1)
CompactListList.H
Foam::readBool
bool readBool(Istream &)
Definition: boolIO.C:60
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::patchIdentifier::name
const word & name() const
Return name.
Definition: patchIdentifier.H:109
setConstantMeshDictionaryIO.H
dictIO
IOobject dictIO(dictName, runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE)