singleCellFvMesh.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "singleCellFvMesh.H"
30 #include "syncTools.H"
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 void Foam::singleCellFvMesh::agglomerateMesh
36 (
37  const fvMesh& mesh,
38  const labelListList& agglom
39 )
40 {
41  // Conversion is a two step process:
42  // - from original (fine) patch faces to agglomerations (aggloms might not
43  // be in correct patch order)
44  // - from agglomerations to coarse patch faces
45 
46  const polyBoundaryMesh& oldPatches = mesh.boundaryMesh();
47 
48  // Check agglomeration within patch face range and continuous
49  labelList nAgglom(oldPatches.size(), Zero);
50 
51  forAll(oldPatches, patchi)
52  {
53  const polyPatch& pp = oldPatches[patchi];
54  if (pp.size() > 0)
55  {
56  nAgglom[patchi] = max(agglom[patchi])+1;
57 
58  forAll(pp, i)
59  {
60  if (agglom[patchi][i] < 0 || agglom[patchi][i] >= pp.size())
61  {
63  << "agglomeration on patch " << patchi
64  << " is out of range 0.." << pp.size()-1
65  << exit(FatalError);
66  }
67  }
68  }
69  }
70 
71  // Check agglomeration is sync
72  {
73  // Get neighbouring agglomeration
74  labelList nbrAgglom(mesh.nBoundaryFaces());
75  forAll(oldPatches, patchi)
76  {
77  const polyPatch& pp = oldPatches[patchi];
78 
79  if (pp.coupled())
80  {
81  label offset = pp.start()-mesh.nInternalFaces();
82  forAll(pp, i)
83  {
84  nbrAgglom[offset+i] = agglom[patchi][i];
85  }
86  }
87  }
89 
90 
91  // Get correspondence between this agglomeration and remote one
92  Map<label> localToNbr(nbrAgglom.size()/10);
93 
94  forAll(oldPatches, patchi)
95  {
96  const polyPatch& pp = oldPatches[patchi];
97 
98  if (pp.coupled())
99  {
100  label offset = pp.start()-mesh.nInternalFaces();
101 
102  forAll(pp, i)
103  {
104  label bFacei = offset+i;
105  label myZone = agglom[patchi][i];
106  label nbrZone = nbrAgglom[bFacei];
107 
108  const auto iter = localToNbr.cfind(myZone);
109 
110  if (iter.found())
111  {
112  // Check that zone numbers are still the same.
113  if (iter.val() != nbrZone)
114  {
116  << "agglomeration is not synchronised across"
117  << " coupled patch " << pp.name()
118  << endl
119  << "Local agglomeration " << myZone
120  << ". Remote agglomeration " << nbrZone
121  << exit(FatalError);
122  }
123  }
124  else
125  {
126  // First occurrence of this zone. Store correspondence
127  // to remote zone number.
128  localToNbr.insert(myZone, nbrZone);
129  }
130  }
131  }
132  }
133  }
134 
135 
136  label coarseI = 0;
137  forAll(nAgglom, patchi)
138  {
139  coarseI += nAgglom[patchi];
140  }
141  // New faces
142  faceList patchFaces(coarseI);
143  // New patch start and size
144  labelList patchStarts(oldPatches.size());
145  labelList patchSizes(oldPatches.size());
146 
147  // From new patch face back to agglomeration
148  patchFaceMap_.setSize(oldPatches.size());
149 
150  // From fine face to coarse face (or -1)
151  reverseFaceMap_.setSize(mesh.nFaces());
152  reverseFaceMap_.labelList::operator=(-1);
153 
154  // Face counter
155  coarseI = 0;
156 
157 
158  forAll(oldPatches, patchi)
159  {
160  patchStarts[patchi] = coarseI;
161 
162  const polyPatch& pp = oldPatches[patchi];
163 
164  if (pp.size() > 0)
165  {
166  patchFaceMap_[patchi].setSize(nAgglom[patchi]);
167 
168  // Patchfaces per agglomeration
169  labelListList agglomToPatch
170  (
171  invertOneToMany(nAgglom[patchi], agglom[patchi])
172  );
173 
174  // From agglomeration to compact patch face
175  labelList agglomToFace(nAgglom[patchi], -1);
176 
177  forAll(pp, i)
178  {
179  label myAgglom = agglom[patchi][i];
180 
181  if (agglomToFace[myAgglom] == -1)
182  {
183  // Agglomeration not yet done. We now have:
184  // - coarseI : current coarse mesh face
185  // - patchStarts[patchi] : coarse mesh patch start
186  // - myAgglom : agglomeration
187  // - agglomToPatch[myAgglom] : fine mesh faces for zone
188  label coarsePatchFacei = coarseI - patchStarts[patchi];
189  patchFaceMap_[patchi][coarsePatchFacei] = myAgglom;
190  agglomToFace[myAgglom] = coarsePatchFacei;
191 
192  const labelList& fineFaces = agglomToPatch[myAgglom];
193 
194  // Create overall map from fine mesh faces to coarseI.
195  forAll(fineFaces, fineI)
196  {
197  reverseFaceMap_[pp.start()+fineFaces[fineI]] = coarseI;
198  }
199 
200  // Construct single face
202  (
203  UIndirectList<face>(pp, fineFaces),
204  pp.points()
205  );
206 
207  if (upp.edgeLoops().size() != 1)
208  {
210  << "agglomeration does not create a"
211  << " single, non-manifold"
212  << " face for agglomeration " << myAgglom
213  << " on patch " << patchi
214  << exit(FatalError);
215  }
216 
217  patchFaces[coarseI++] = face
218  (
219  renumber
220  (
221  upp.meshPoints(),
222  upp.edgeLoops()[0]
223  )
224  );
225  }
226  }
227  }
228 
229  patchSizes[patchi] = coarseI-patchStarts[patchi];
230  }
231 
232  //Pout<< "patchStarts:" << patchStarts << endl;
233  //Pout<< "patchSizes:" << patchSizes << endl;
234 
235  // Compact numbering for points
236  reversePointMap_.setSize(mesh.nPoints());
237  reversePointMap_.labelList::operator=(-1);
238  label newI = 0;
239 
240  forAll(patchFaces, coarseI)
241  {
242  face& f = patchFaces[coarseI];
243 
244  forAll(f, fp)
245  {
246  if (reversePointMap_[f[fp]] == -1)
247  {
248  reversePointMap_[f[fp]] = newI++;
249  }
250 
251  f[fp] = reversePointMap_[f[fp]];
252  }
253  }
254 
255  pointMap_ = invert(newI, reversePointMap_);
256 
257  // Subset used points
258  pointField boundaryPoints(mesh.points(), pointMap_);
259 
260  // Add patches (on still zero sized mesh)
261  List<polyPatch*> newPatches(oldPatches.size());
262  forAll(oldPatches, patchi)
263  {
264  newPatches[patchi] = oldPatches[patchi].clone
265  (
266  boundaryMesh(),
267  patchi,
268  0,
269  0
270  ).ptr();
271  }
272  addFvPatches(newPatches);
273 
274  const label nFace = patchFaces.size();
275 
276  // Actually change the mesh. // Owner, neighbour is trivial
278  (
279  autoPtr<pointField>::New(std::move(boundaryPoints)),
280  autoPtr<faceList>::New(std::move(patchFaces)),
281  autoPtr<labelList>::New(nFace, Zero), // owner
282  autoPtr<labelList>::New(), // neighbour
283  patchSizes,
284  patchStarts,
285  true // syncPar
286  );
287 
288  // Adapt the zones
289  cellZones().clear();
290  cellZones().setSize(mesh.cellZones().size());
291  {
292  forAll(mesh.cellZones(), zoneI)
293  {
294  const cellZone& oldFz = mesh.cellZones()[zoneI];
295 
296  DynamicList<label> newAddressing;
297 
298  //Note: uncomment if you think it makes sense. Note that value
299  // of cell0 is the average.
301  //if (oldFz.localID(0) != -1)
302  //{
303  // newAddressing.append(0);
304  //}
305 
306  cellZones().set
307  (
308  zoneI,
309  oldFz.clone
310  (
311  newAddressing,
312  zoneI,
313  cellZones()
314  )
315  );
316  }
317  }
318 
319  faceZones().clear();
320  faceZones().setSize(mesh.faceZones().size());
321  {
322  forAll(mesh.faceZones(), zoneI)
323  {
324  const faceZone& oldFz = mesh.faceZones()[zoneI];
325 
326  DynamicList<label> newAddressing(oldFz.size());
327  DynamicList<bool> newFlipMap(oldFz.size());
328 
329  forAll(oldFz, i)
330  {
331  label newFacei = reverseFaceMap_[oldFz[i]];
332 
333  if (newFacei != -1)
334  {
335  newAddressing.append(newFacei);
336  newFlipMap.append(oldFz.flipMap()[i]);
337  }
338  }
339 
340  faceZones().set
341  (
342  zoneI,
343  oldFz.clone
344  (
345  newAddressing,
346  newFlipMap,
347  zoneI,
348  faceZones()
349  )
350  );
351  }
352  }
353 
354 
355  pointZones().clear();
356  pointZones().setSize(mesh.pointZones().size());
357  {
358  forAll(mesh.pointZones(), zoneI)
359  {
360  const pointZone& oldFz = mesh.pointZones()[zoneI];
361 
362  DynamicList<label> newAddressing(oldFz.size());
363 
364  forAll(oldFz, i)
365  {
366  label newPointi = reversePointMap_[oldFz[i]];
367  if (newPointi != -1)
368  {
369  newAddressing.append(newPointi);
370  }
371  }
372 
373  pointZones().set
374  (
375  zoneI,
376  oldFz.clone
377  (
378  pointZones(),
379  zoneI,
380  newAddressing
381  )
382  );
383  }
384  }
385 
386  // Make sure we don't start dumping mesh every timestep (since
387  // resetPrimitives sets AUTO_WRITE)
389 }
390 
391 
392 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
393 
394 Foam::singleCellFvMesh::singleCellFvMesh
395 (
396  const IOobject& io,
397  const fvMesh& mesh
398 )
399 :
400  fvMesh(io, Zero, false),
401  patchFaceAgglomeration_
402  (
403  IOobject
404  (
405  "patchFaceAgglomeration",
406  io.instance(),
408  *this,
409  io.readOpt(),
410  io.writeOpt()
411  ),
412  0
413  ),
414  patchFaceMap_
415  (
416  IOobject
417  (
418  "patchFaceMap",
419  io.instance(),
421  *this,
422  io.readOpt(),
423  io.writeOpt()
424  ),
425  mesh.boundaryMesh().size()
426  ),
427  reverseFaceMap_
428  (
429  IOobject
430  (
431  "reverseFaceMap",
432  io.instance(),
434  *this,
435  io.readOpt(),
436  io.writeOpt()
437  ),
438  mesh.nFaces()
439  ),
440  pointMap_
441  (
442  IOobject
443  (
444  "pointMap",
445  io.instance(),
447  *this,
448  io.readOpt(),
449  io.writeOpt()
450  ),
451  mesh.nPoints()
452  ),
453  reversePointMap_
454  (
455  IOobject
456  (
457  "reversePointMap",
458  io.instance(),
460  *this,
461  io.readOpt(),
462  io.writeOpt()
463  ),
464  mesh.nPoints()
465  )
466 {
467  const polyBoundaryMesh& oldPatches = mesh.boundaryMesh();
468 
469  labelListList agglom(oldPatches.size());
470 
471  forAll(oldPatches, patchi)
472  {
473  agglom[patchi] = identity(oldPatches[patchi].size());
474  }
475 
476  agglomerateMesh(mesh, agglom);
477 }
478 
479 
480 Foam::singleCellFvMesh::singleCellFvMesh
481 (
482  const IOobject& io,
483  const fvMesh& mesh,
484  const labelListList& patchFaceAgglomeration
485 )
486 :
487  fvMesh(io, Zero, false),
488  patchFaceAgglomeration_
489  (
490  IOobject
491  (
492  "patchFaceAgglomeration",
493  io.instance(),
495  *this,
496  io.readOpt(),
497  io.writeOpt()
498  ),
499  patchFaceAgglomeration
500  ),
501  patchFaceMap_
502  (
503  IOobject
504  (
505  "patchFaceMap",
506  io.instance(),
508  *this,
509  io.readOpt(),
510  io.writeOpt()
511  ),
512  mesh.boundaryMesh().size()
513  ),
514  reverseFaceMap_
515  (
516  IOobject
517  (
518  "reverseFaceMap",
519  io.instance(),
521  *this,
522  io.readOpt(),
523  io.writeOpt()
524  ),
525  mesh.nFaces()
526  ),
527  pointMap_
528  (
529  IOobject
530  (
531  "pointMap",
532  io.instance(),
534  *this,
535  io.readOpt(),
536  io.writeOpt()
537  ),
538  mesh.nPoints()
539  ),
540  reversePointMap_
541  (
542  IOobject
543  (
544  "reversePointMap",
545  io.instance(),
547  *this,
548  io.readOpt(),
549  io.writeOpt()
550  ),
551  mesh.nPoints()
552  )
553 {
554  agglomerateMesh(mesh, patchFaceAgglomeration);
555 }
556 
557 
558 Foam::singleCellFvMesh::singleCellFvMesh(const IOobject& io)
559 :
560  fvMesh(io),
561  patchFaceAgglomeration_
562  (
563  IOobject
564  (
565  "patchFaceAgglomeration",
566  io.instance(),
567  fvMesh::meshSubDir,
568  *this,
569  io.readOpt(),
570  io.writeOpt()
571  )
572  ),
573  patchFaceMap_
574  (
575  IOobject
576  (
577  "patchFaceMap",
578  io.instance(),
579  fvMesh::meshSubDir,
580  *this,
581  io.readOpt(),
582  io.writeOpt()
583  )
584  ),
585  reverseFaceMap_
586  (
587  IOobject
588  (
589  "reverseFaceMap",
590  io.instance(),
591  fvMesh::meshSubDir,
592  *this,
593  io.readOpt(),
594  io.writeOpt()
595  )
596  ),
597  pointMap_
598  (
599  IOobject
600  (
601  "pointMap",
602  io.instance(),
603  fvMesh::meshSubDir,
604  *this,
605  io.readOpt(),
606  io.writeOpt()
607  )
608  ),
609  reversePointMap_
610  (
611  IOobject
612  (
613  "reversePointMap",
614  io.instance(),
615  fvMesh::meshSubDir,
616  *this,
617  io.readOpt(),
618  io.writeOpt()
619  )
620  )
621 {}
622 
623 
624 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::autoPtr::New
static autoPtr< T > New(Args &&... args)
Construct autoPtr of T with forwarding arguments.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::ZoneMesh::clear
void clear()
Clear the zones.
Definition: ZoneMesh.C:724
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
uindirectPrimitivePatch.H
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:63
Foam::IOobject::instance
const fileName & instance() const noexcept
Definition: IOobjectI.H:196
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::polyMesh::meshSubDir
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:321
Foam::syncTools::swapBoundaryFaceList
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues)
Swap coupled boundary face values. Uses eqOp.
Definition: syncTools.H:445
Foam::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::IOobject::writeOpt
writeOption writeOpt() const noexcept
The write option.
Definition: IOobjectI.H:179
syncTools.H
Foam::invert
labelList invert(const label len, const labelUList &map)
Create an inverse one-to-one mapping.
Definition: ListOps.C:36
Foam::primitiveMesh::nPoints
label nPoints() const noexcept
Number of mesh points.
Definition: primitiveMeshI.H:37
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:492
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::polyMesh::resetPrimitives
void resetPrimitives(autoPtr< pointField > &&points, autoPtr< faceList > &&faces, autoPtr< labelList > &&owner, autoPtr< labelList > &&neighbour, const labelUList &patchSizes, const labelUList &patchStarts, const bool validBoundary=true)
Reset mesh primitive data. Assumes all patch info correct.
Definition: polyMesh.C:718
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:486
singleCellFvMesh.H
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
newPointi
label newPointi
Definition: readKivaGrid.H:496
Foam::FatalError
error FatalError
Foam::fvMesh::addFvPatches
void addFvPatches(PtrList< polyPatch > &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:605
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::IOobject::readOpt
readOption readOpt() const noexcept
The read option.
Definition: IOobjectI.H:164
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Foam::renumber
IntListType renumber(const labelUList &oldToNew, const IntListType &input)
Renumber the values (not the indices) of a list.
Definition: ListOpsTemplates.C:37
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
f
labelList f(nPoints)
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::List< labelList >
Foam::identity
labelList identity(const label len, label start=0)
Create identity map of the given length with (map[i] == i)
Definition: labelList.C:38
Foam::polyMesh::pointZones
const pointZoneMesh & pointZones() const noexcept
Return point zone mesh.
Definition: polyMesh.H:480
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::primitiveMesh::nFaces
label nFaces() const noexcept
Number of mesh faces.
Definition: primitiveMeshI.H:90
Foam::uindirectPrimitivePatch
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
Definition: uindirectPrimitivePatch.H:49
Foam::invertOneToMany
labelListList invertOneToMany(const label len, const labelUList &map)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:114
constant
constant condensation/saturation model.
Foam::polyMesh::setInstance
void setInstance(const fileName &instance, const IOobject::writeOption wOpt=IOobject::AUTO_WRITE)
Set the instance for mesh files.
Definition: polyMeshIO.C:36