Test-GAMGAgglomeration.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-2015 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  Test-GAMGAgglomeration
26 
27 Description
28  Test application for GAMG agglomeration. Hardcoded to expect GAMG on p.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "fvCFD.H"
33 #include "GAMGAgglomeration.H"
34 #include "OFstream.H"
35 #include "meshTools.H"
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 // Main program:
39 
40 int main(int argc, char *argv[])
41 {
42  argList::addBoolOption
43  (
44  "writeObj",
45  "write obj files of agglomeration"
46  );
47  argList::addBoolOption
48  (
49  "normalise",
50  "normalise agglomeration (0..1)"
51  );
52 
53  #include "setRootCase.H"
54  #include "createTime.H"
55 
56  bool writeObj = args.optionFound("writeObj");
57  bool normalise = args.optionFound("normalise");
58 
59  #include "createMesh.H"
60 
61  const fvSolution& sol = static_cast<const fvSolution&>(mesh);
62  const dictionary& pDict = sol.subDict("solvers").subDict("p");
63 
65  (
66  mesh,
67  pDict
68  );
69 
70  labelList cellToCoarse(identity(mesh.nCells()));
71  labelListList coarseToCell(invertOneToMany(mesh.nCells(), cellToCoarse));
72 
73  runTime++;
74 
75  // Write initial agglomeration
76  {
77  volScalarField scalarAgglomeration
78  (
79  IOobject
80  (
81  "agglomeration",
82  runTime.timeName(),
83  mesh,
84  IOobject::NO_READ,
85  IOobject::AUTO_WRITE
86  ),
87  mesh,
88  dimensionedScalar("aggomeration", dimless, 0.0)
89  );
90  scalarField& fld = scalarAgglomeration.internalField();
91  forAll(fld, cellI)
92  {
93  fld[cellI] = cellToCoarse[cellI];
94  }
95  fld /= max(fld);
96  scalarAgglomeration.correctBoundaryConditions();
97  scalarAgglomeration.write();
98 
99  Info<< "Writing initial cell distribution to "
100  << runTime.timeName() << endl;
101  }
102 
103 
104  for (label level = 0; level < agglom.size(); level++)
105  {
106  runTime++;
107 
108  Info<< "Time = " << runTime.timeName() << nl << endl;
109 
110  const labelList& addr = agglom.restrictAddressing(level);
111  label coarseSize = max(addr)+1;
112 
113  Info<< "Level : " << level << endl
114  << returnReduce(addr.size(), sumOp<label>()) << endl
115  << " current size : "
116  << returnReduce(addr.size(), sumOp<label>()) << endl
117  << " agglomerated size : "
118  << returnReduce(coarseSize, sumOp<label>()) << endl;
119 
120  labelList newAddr;
121  label newCoarseSize = 0;
122  bool ok = GAMGAgglomeration::checkRestriction
123  (
124  newAddr,
125  newCoarseSize,
126 
127  agglom.meshLevel(level).lduAddr(),
128  addr,
129  coarseSize
130  );
131  if (!ok)
132  {
134  << "At level " << level
135  << " there are " << coarseSize
136  << " agglomerated cells but " << newCoarseSize
137  << " disconnected regions" << endl
138  << " This means that some agglomerations (coarse cells)"
139  << " consist of multiple disconnected regions."
140  << endl;
141  }
142 
143 
144  forAll(addr, fineI)
145  {
146  const labelList& cellLabels = coarseToCell[fineI];
147  forAll(cellLabels, i)
148  {
149  cellToCoarse[cellLabels[i]] = addr[fineI];
150  }
151  }
152  coarseToCell = invertOneToMany(coarseSize, cellToCoarse);
153 
154  // Write agglomeration
155  {
156  volScalarField scalarAgglomeration
157  (
158  IOobject
159  (
160  "agglomeration",
161  runTime.timeName(),
162  mesh,
163  IOobject::NO_READ,
164  IOobject::AUTO_WRITE
165  ),
166  mesh,
167  dimensionedScalar("aggomeration", dimless, 0.0)
168  );
169  scalarField& fld = scalarAgglomeration.internalField();
170  forAll(fld, cellI)
171  {
172  fld[cellI] = cellToCoarse[cellI];
173  }
174  if (normalise)
175  {
176  fld /= max(fld);
177  }
178  scalarAgglomeration.correctBoundaryConditions();
179  scalarAgglomeration.write();
180  }
181 
182  if (writeObj)
183  {
184  OFstream str(runTime.path()/runTime.timeName()/"aggomeration.obj");
185  label vertI = 0;
186 
187  // Write all mesh cc
188  forAll(mesh.cellCentres(), cellI)
189  {
190  meshTools::writeOBJ(str, mesh.cellCentres()[cellI]);
191  vertI++;
192  }
193 
194  // Determine coarse cc
195  forAll(coarseToCell, coarseI)
196  {
197  const labelList& cellLabels = coarseToCell[coarseI];
198 
199  point coarseCc = average
200  (
201  pointField(mesh.cellCentres(), cellLabels)
202  );
203  meshTools::writeOBJ(str, coarseCc);
204  vertI++;
205 
206  forAll(cellLabels, i)
207  {
208  label cellI = cellLabels[i];
209 
210  str << "l " << cellI+1 << ' ' << vertI << nl;
211  }
212  }
213  }
214 
215  Info<< endl;
216  }
217 
218 
219  Info<< "End\n" << endl;
220 
221  return 0;
222 }
223 
224 
225 // ************************************************************************* //
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
meshTools.H
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::GAMGAgglomeration
Geometric agglomerated algebraic multigrid agglomeration class.
Definition: GAMGAgglomeration.H:61
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
writeOBJ
void writeOBJ(Ostream &os, label &vertI, const tetPoints &tet)
Definition: Test-tetTetOverlap.C:38
Foam::invertOneToMany
labelListList invertOneToMany(const label len, const labelUList &)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:67
Foam::GAMGAgglomeration::meshLevel
const lduMesh & meshLevel(const label leveli) const
Return LDU mesh of given level.
Definition: GAMGAgglomeration.C:436
GAMGAgglomeration.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::fvc::average
tmp< GeometricField< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:46
Foam::GAMGAgglomeration::size
label size() const
Definition: GAMGAgglomeration.H:320
OFstream.H
Foam::primitiveMesh::nCells
label nCells() const
Definition: primitiveMeshI.H:64
Foam::GAMGAgglomeration::restrictAddressing
const labelField & restrictAddressing(const label leveli) const
Return cell restrict addressing of given level.
Definition: GAMGAgglomeration.H:338
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::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::fvSolution
Selector class for finite volume solution solution. fvMesh is derived from fvSolution so that all fie...
Definition: fvSolution.H:47
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
main
int main(int argc, char *argv[])
Definition: Test-GAMGAgglomeration.C:37
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
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
fld
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
writeObj
void writeObj(Ostream &os, const pointField &points)
Definition: Test-PrimitivePatch.C:44
Foam::OFstream
Output to file stream.
Definition: OFstream.H:81
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:885
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
setRootCase.H
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:211
Foam::sumOp
Definition: ops.H:162
Foam::Vector< scalar >
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
createMesh.H
createTime.H
Foam::argList::optionFound
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
fvCFD.H
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::dictionary::subDict
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:631
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
args
Foam::argList args(argc, argv)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::lduMesh::lduAddr
virtual const lduAddressing & lduAddr() const =0
Return ldu addressing.