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 | 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 \*---------------------------------------------------------------------------*/
25 
26 #include "singleCellFvMesh.H"
27 #include "syncTools.H"
28 #include "uindirectPrimitivePatch.H"
29 
30 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31 
32 // Conversion is a two step process:
33 // - from original (fine) patch faces to agglomerations (aggloms might not
34 // be in correct patch order)
35 // - from agglomerations to coarse patch faces
37 (
38  const fvMesh& mesh,
39  const labelListList& agglom
40 )
41 {
42  const polyBoundaryMesh& oldPatches = mesh.boundaryMesh();
43 
44  // Check agglomeration within patch face range and continuous
45  labelList nAgglom(oldPatches.size(), 0);
46 
47  forAll(oldPatches, patchI)
48  {
49 //Info << "now the patchI is = " << patchI <<endl;
50  const polyPatch& pp = oldPatches[patchI];
51 //Info << "the patch name is = " << pp.name();
52  if (pp.size() > 0)
53  {
54  nAgglom[patchI] = max(agglom[patchI])+1;
55 
56  forAll(pp, i)
57  {
58  if (agglom[patchI][i] < 0 || agglom[patchI][i] >= pp.size())
59  {
61  << "agglomeration on patch " << patchI
62  << " is out of range 0.." << pp.size()-1
63  << exit(FatalError);
64  }
65  }
66  }
67  }
68 
69  // Check agglomeration is sync
70  {
71  // Get neighbouring agglomeration
72  labelList nbrAgglom(mesh.nFaces()-mesh.nInternalFaces());
73  forAll(oldPatches, patchI)
74  {
75  const polyPatch& pp = oldPatches[patchI];
76 
77  if (pp.coupled())
78  {
79  label offset = pp.start()-mesh.nInternalFaces();
80  forAll(pp, i)
81  {
82  nbrAgglom[offset+i] = agglom[patchI][i];
83  }
84  }
85  }
86  syncTools::swapBoundaryFaceList(mesh, nbrAgglom);
87 
88 
89  // Get correspondence between this agglomeration and remote one
90  Map<label> localToNbr(nbrAgglom.size()/10);
91 
92  forAll(oldPatches, patchI)
93  {
94  const polyPatch& pp = oldPatches[patchI];
95 
96  if (pp.coupled())
97  {
98  label offset = pp.start()-mesh.nInternalFaces();
99 
100  forAll(pp, i)
101  {
102  label bFaceI = offset+i;
103  label myZone = agglom[patchI][i];
104  label nbrZone = nbrAgglom[bFaceI];
105 
106  Map<label>::const_iterator iter = localToNbr.find(myZone);
107 
108  if (iter == localToNbr.end())
109  {
110  // First occurence of this zone. Store correspondence
111  // to remote zone number.
112  localToNbr.insert(myZone, nbrZone);
113  }
114  else
115  {
116  // Check that zone numbers are still the same.
117  if (iter() != nbrZone)
118  {
120  << "agglomeration is not synchronised across"
121  << " coupled patch " << pp.name()
122  << endl
123  << "Local agglomeration " << myZone
124  << ". Remote agglomeration " << nbrZone
125  << exit(FatalError);
126  }
127  }
128  }
129  }
130  }
131  }
132 
133 
134  label coarseI = 0;
135  forAll(nAgglom, patchI)
136  {
137  coarseI += nAgglom[patchI];
138  }
139  // New faces
140  faceList patchFaces(coarseI);
141  // New patch start and size
142  labelList patchStarts(oldPatches.size());
143  labelList patchSizes(oldPatches.size());
144 
145  // From new patch face back to agglomeration
146  patchFaceMap_.setSize(oldPatches.size());
147 
148  // From fine face to coarse face (or -1)
149  reverseFaceMap_.setSize(mesh.nFaces());
150  reverseFaceMap_.labelList::operator=(-1);
151 
152  // Face counter
153  coarseI = 0;
154 
155 
156  forAll(oldPatches, patchI)
157  {
158  patchStarts[patchI] = coarseI;
159 
160  const polyPatch& pp = oldPatches[patchI];
161 
162  if (pp.size() > 0)
163  {
164  patchFaceMap_[patchI].setSize(nAgglom[patchI]);
165 
166  // Patchfaces per agglomeration
167  labelListList agglomToPatch
168  (
169  invertOneToMany(nAgglom[patchI], agglom[patchI])
170  );
171 
172  // From agglomeration to compact patch face
173  labelList agglomToFace(nAgglom[patchI], -1);
174 
175  forAll(pp, i)
176  {
177  label myAgglom = agglom[patchI][i];
178 
179  if (agglomToFace[myAgglom] == -1)
180  {
181  // Agglomeration not yet done. We now have:
182  // - coarseI : current coarse mesh face
183  // - patchStarts[patchI] : coarse mesh patch start
184  // - myAgglom : agglomeration
185  // - agglomToPatch[myAgglom] : fine mesh faces for zone
186  label coarsePatchFaceI = coarseI - patchStarts[patchI];
187  patchFaceMap_[patchI][coarsePatchFaceI] = myAgglom;
188  agglomToFace[myAgglom] = coarsePatchFaceI;
189 
190  const labelList& fineFaces = agglomToPatch[myAgglom];
191 
192  // Create overall map from fine mesh faces to coarseI.
193  forAll(fineFaces, fineI)
194  {
195  reverseFaceMap_[pp.start()+fineFaces[fineI]] = coarseI;
196  }
197 
198  // Construct single face
200  (
201  UIndirectList<face>(pp, fineFaces),
202  pp.points()
203  );
204 
205  if (upp.edgeLoops().size() != 1)
206  {
208  << "agglomeration does not create a"
209  << " single, non-manifold"
210  << " face for agglomeration " << myAgglom
211  << " on patch " << patchI
212  << exit(FatalError);
213  }
214 
215  patchFaces[coarseI++] = face
216  (
217  renumber
218  (
219  upp.meshPoints(),
220  upp.edgeLoops()[0]
221  )
222  );
223  }
224  }
225  }
226 
227  patchSizes[patchI] = coarseI-patchStarts[patchI];
228  }
229 
230  //Pout<< "patchStarts:" << patchStarts << endl;
231  //Pout<< "patchSizes:" << patchSizes << endl;
232 
233  // Compact numbering for points
234  reversePointMap_.setSize(mesh.nPoints());
235  reversePointMap_.labelList::operator=(-1);
236  label newI = 0;
237 
238  forAll(patchFaces, coarseI)
239  {
240  face& f = patchFaces[coarseI];
241 
242  forAll(f, fp)
243  {
244  if (reversePointMap_[f[fp]] == -1)
245  {
246  reversePointMap_[f[fp]] = newI++;
247  }
248 
249  f[fp] = reversePointMap_[f[fp]];
250  }
251  }
252 
253  pointMap_ = invert(newI, reversePointMap_);
254 
255  // Subset used points
256  pointField boundaryPoints(mesh.points(), pointMap_);
257 
258  // Add patches (on still zero sized mesh)
259  List<polyPatch*> newPatches(oldPatches.size());
260  forAll(oldPatches, patchI)
261  {
262  newPatches[patchI] = oldPatches[patchI].clone
263  (
264  boundaryMesh(),
265  patchI,
266  0,
267  0
268  ).ptr();
269  }
270  addFvPatches(newPatches);
271 
272  // Owner, neighbour is trivial
273  labelList owner(patchFaces.size(), 0);
274  labelList neighbour(0);
275 
276 
277  // actually change the mesh
278  resetPrimitives
279  (
280  xferMove(boundaryPoints),
282  xferMove(owner),
283  xferMove(neighbour),
284  patchSizes,
285  patchStarts,
286  true //syncPar
287  );
288 
289 
290  // Adapt the zones
291  cellZones().clear();
292  cellZones().setSize(mesh.cellZones().size());
293  {
294  forAll(mesh.cellZones(), zoneI)
295  {
296  const cellZone& oldFz = mesh.cellZones()[zoneI];
297 
298  DynamicList<label> newAddressing;
299 
300  //Note: uncomment if you think it makes sense. Note that value
301  // of cell0 is the average.
303  //if (oldFz.localID(0) != -1)
304  //{
305  // newAddressing.append(0);
306  //}
307 
308  cellZones().set
309  (
310  zoneI,
311  oldFz.clone
312  (
313  newAddressing,
314  zoneI,
315  cellZones()
316  )
317  );
318  }
319  }
320 
321  faceZones().clear();
322  faceZones().setSize(mesh.faceZones().size());
323  {
324  forAll(mesh.faceZones(), zoneI)
325  {
326  const faceZone& oldFz = mesh.faceZones()[zoneI];
327 
328  DynamicList<label> newAddressing(oldFz.size());
329  DynamicList<bool> newFlipMap(oldFz.size());
330 
331  forAll(oldFz, i)
332  {
333  label newFaceI = reverseFaceMap_[oldFz[i]];
334 
335  if (newFaceI != -1)
336  {
337  newAddressing.append(newFaceI);
338  newFlipMap.append(oldFz.flipMap()[i]);
339  }
340  }
341 
342  faceZones().set
343  (
344  zoneI,
345  oldFz.clone
346  (
347  newAddressing,
348  newFlipMap,
349  zoneI,
350  faceZones()
351  )
352  );
353  }
354  }
355 
356 
357  pointZones().clear();
358  pointZones().setSize(mesh.pointZones().size());
359  {
360  forAll(mesh.pointZones(), zoneI)
361  {
362  const pointZone& oldFz = mesh.pointZones()[zoneI];
363 
364  DynamicList<label> newAddressing(oldFz.size());
365 
366  forAll(oldFz, i)
367  {
368  label newPointI = reversePointMap_[oldFz[i]];
369  if (newPointI != -1)
370  {
371  newAddressing.append(newPointI);
372  }
373  }
374 
375  pointZones().set
376  (
377  zoneI,
378  oldFz.clone
379  (
380  pointZones(),
381  zoneI,
382  newAddressing
383  )
384  );
385  }
386  }
387 }
388 
389 
390 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
391 
393 (
394  const IOobject& io,
395  const fvMesh& mesh
396 )
397 :
398  fvMesh
399  (
400  io,
401  xferCopy(pointField()), //points
402  xferCopy(faceList()), //faces
403  xferCopy(labelList()), //allOwner
404  xferCopy(labelList()), //allNeighbour
405  false //syncPar
406  ),
407  patchFaceAgglomeration_
408  (
409  IOobject
410  (
411  "patchFaceAgglomeration",
412  io.instance(),
413  fvMesh::meshSubDir,
414  *this,
415  io.readOpt(),
416  io.writeOpt()
417  ),
418  0
419  ),
420  patchFaceMap_
421  (
422  IOobject
423  (
424  "patchFaceMap",
425  io.instance(),
426  fvMesh::meshSubDir,
427  *this,
428  io.readOpt(),
429  io.writeOpt()
430  ),
431  mesh.boundaryMesh().size()
432  ),
433  reverseFaceMap_
434  (
435  IOobject
436  (
437  "reverseFaceMap",
438  io.instance(),
439  fvMesh::meshSubDir,
440  *this,
441  io.readOpt(),
442  io.writeOpt()
443  ),
444  mesh.nFaces()
445  ),
446  pointMap_
447  (
448  IOobject
449  (
450  "pointMap",
451  io.instance(),
452  fvMesh::meshSubDir,
453  *this,
454  io.readOpt(),
455  io.writeOpt()
456  ),
457  mesh.nPoints()
458  ),
459  reversePointMap_
460  (
461  IOobject
462  (
463  "reversePointMap",
464  io.instance(),
465  fvMesh::meshSubDir,
466  *this,
467  io.readOpt(),
468  io.writeOpt()
469  ),
470  mesh.nPoints()
471  )
472 {
473  const polyBoundaryMesh& oldPatches = mesh.boundaryMesh();
474 
475  labelListList agglom(oldPatches.size());
476 
477  forAll(oldPatches, patchI)
478  {
479  agglom[patchI] = identity(oldPatches[patchI].size());
480  }
481 
482  agglomerateMesh(mesh, agglom);
483 }
484 
485 
487 (
488  const IOobject& io,
489  const fvMesh& mesh,
490  const labelListList& patchFaceAgglomeration
491 )
492 :
493  fvMesh
494  (
495  io,
496  xferCopy(pointField()), //points
497  xferCopy(faceList()), //faces
498  xferCopy(labelList()), //allOwner
499  xferCopy(labelList()), //allNeighbour
500  false //syncPar
501  ),
502  patchFaceAgglomeration_
503  (
504  IOobject
505  (
506  "patchFaceAgglomeration",
507  io.instance(),
508  fvMesh::meshSubDir,
509  *this,
510  io.readOpt(),
511  io.writeOpt()
512  ),
513  patchFaceAgglomeration
514  ),
515  patchFaceMap_
516  (
517  IOobject
518  (
519  "patchFaceMap",
520  io.instance(),
521  fvMesh::meshSubDir,
522  *this,
523  io.readOpt(),
524  io.writeOpt()
525  ),
526  mesh.boundaryMesh().size()
527  ),
528  reverseFaceMap_
529  (
530  IOobject
531  (
532  "reverseFaceMap",
533  io.instance(),
534  fvMesh::meshSubDir,
535  *this,
536  io.readOpt(),
537  io.writeOpt()
538  ),
539  mesh.nFaces()
540  ),
541  pointMap_
542  (
543  IOobject
544  (
545  "pointMap",
546  io.instance(),
547  fvMesh::meshSubDir,
548  *this,
549  io.readOpt(),
550  io.writeOpt()
551  ),
552  mesh.nPoints()
553  ),
554  reversePointMap_
555  (
556  IOobject
557  (
558  "reversePointMap",
559  io.instance(),
560  fvMesh::meshSubDir,
561  *this,
562  io.readOpt(),
563  io.writeOpt()
564  ),
565  mesh.nPoints()
566  )
567 {
568  agglomerateMesh(mesh, patchFaceAgglomeration);
569 }
570 
571 
573 :
574  fvMesh(io),
575  patchFaceAgglomeration_
576  (
577  IOobject
578  (
579  "patchFaceAgglomeration",
580  io.instance(),
581  fvMesh::meshSubDir,
582  *this,
583  io.readOpt(),
584  io.writeOpt()
585  )
586  ),
587  patchFaceMap_
588  (
589  IOobject
590  (
591  "patchFaceMap",
592  io.instance(),
593  fvMesh::meshSubDir,
594  *this,
595  io.readOpt(),
596  io.writeOpt()
597  )
598  ),
599  reverseFaceMap_
600  (
601  IOobject
602  (
603  "reverseFaceMap",
604  io.instance(),
605  fvMesh::meshSubDir,
606  *this,
607  io.readOpt(),
608  io.writeOpt()
609  )
610  ),
611  pointMap_
612  (
613  IOobject
614  (
615  "pointMap",
616  io.instance(),
617  fvMesh::meshSubDir,
618  *this,
619  io.readOpt(),
620  io.writeOpt()
621  )
622  ),
623  reversePointMap_
624  (
625  IOobject
626  (
627  "reversePointMap",
628  io.instance(),
629  fvMesh::meshSubDir,
630  *this,
631  io.readOpt(),
632  io.writeOpt()
633  )
634  )
635 {}
636 
637 
638 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
639 
640 
641 
642 // ************************************************************************* //
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
Foam::PrimitivePatch::points
const Field< PointType > & points() const
Return reference to global points.
Definition: PrimitivePatchTemplate.H:282
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::DynamicList< label >
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::cellZone::clone
virtual autoPtr< cellZone > clone(const cellZoneMesh &zm) const
Construct and return a clone, resetting the zone mesh.
Definition: cellZone.H:159
Foam::pointZone::clone
virtual autoPtr< pointZone > clone(const pointZoneMesh &zm) const
Construct and return a clone, resetting the zone mesh.
Definition: pointZone.H:158
Foam::pointZone
A subset of mesh points. The labels of points in the zone can be obtained from the addressing() list.
Definition: pointZone.H:62
Foam::polyPatch::coupled
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:322
Foam::singleCellFvMesh::singleCellFvMesh
singleCellFvMesh(const singleCellFvMesh &)
Disallow default bitwise copy construct.
Foam::Map< label >
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::IOobject::instance
const fileName & instance() const
Definition: IOobject.H:350
Foam::IOobject::writeOpt
writeOption writeOpt() const
Definition: IOobject.H:327
patchFaces
labelList patchFaces(const polyBoundaryMesh &patches, const wordList &names)
Definition: extrudeMesh.C:148
syncTools.H
Foam::cellZone
A subset of mesh cells.
Definition: cellZone.H:61
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::IOobject::readOpt
readOption readOpt() const
Definition: IOobject.H:317
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::DynamicList::clear
void clear()
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:242
singleCellFvMesh.H
Foam::identity
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:312
Foam::singleCellFvMesh::agglomerateMesh
void agglomerateMesh(const fvMesh &, const labelListList &)
Calculate agglomerated mesh.
Definition: singleCellFvMesh.C:37
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::List::setSize
void setSize(const label)
Reset size of List.
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::xferMove
Xfer< T > xferMove(T &)
Construct by transferring the contents of the arg.
f
labelList f(nPoints)
Foam::faceList
List< face > faceList
Definition: faceListFwd.H:43
Foam::xferCopy
Xfer< T > xferCopy(const T &)
Construct by copying the contents of the arg.
Foam::PrimitivePatch::edgeLoops
const labelListList & edgeLoops() const
Return list of closed loops of boundary vertices.
Definition: PrimitivePatchEdgeLoops.C:176
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::IOobject::clone
Foam::autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:252
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
Foam::PtrList::setSize
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:142
addFvPatches
mesh addFvPatches(patches)
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:106
Foam::invert
labelList invert(const label len, const labelUList &)
Invert one-to-one map. Unmapped elements will be -1.
Definition: ListOps.C:37
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:59
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatchTemplate.C:392
Foam::patchIdentifier::name
const word & name() const
Return name.
Definition: patchIdentifier.H:109
Foam::faceZone::clone
virtual autoPtr< faceZone > clone(const faceZoneMesh &zm) const
Construct and return a clone, resetting the zone mesh.
Definition: faceZone.H:204
renumber
static void renumber(const labelList &map, labelList &elems)
Definition: reconstructParMesh.C:61
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatchTemplate.H:88