fvMeshStitcher.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2021-2023 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 "fvMeshStitcher.H"
27 #include "globalIndex.H"
28 #include "fvcSurfaceIntegrate.H"
29 #include "meshObjects.H"
30 #include "polyTopoChangeMap.H"
31 #include "syncTools.H"
32 #include "surfaceToVolVelocity.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 template<class FaceList, class PointField>
41 (
43  const label edgei
44 )
45 {
46  return edge
47  (
48  p.meshPoints()[p.edges()[edgei][0]],
49  p.meshPoints()[p.edges()[edgei][1]]
50  );
51 }
52 
53 }
54 
55 
56 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
57 
58 namespace Foam
59 {
62 }
63 
64 
65 const Foam::word Foam::fvMeshStitcher::nccFieldPrefix_ =
66  fvMeshStitcher::typeName + ":";
67 
68 
69 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
70 
71 void Foam::fvMeshStitcher::intersectNonConformalCyclic
72 (
73  const nonConformalCyclicFvPatch& nccFvp,
74  surfaceLabelField::Boundary& polyFacesBf,
75  surfaceVectorField::Boundary& SfBf,
76  surfaceVectorField::Boundary& CfBf,
77  const tmp<surfaceLabelField::Boundary>& tOrigFacesNbrBf,
78  const tmp<surfaceVectorField::Boundary>& tOrigSfNbrBf,
79  const tmp<surfaceVectorField::Boundary>& tOrigCfNbrBf,
80  List<part>& origEdgeParts
81 ) const
82 {
83  const nonConformalCyclicFvPatch& nbrNccFvp = nccFvp.nbrPatch();
84 
85  // Alias the original poly patches
86  const polyPatch& origPp = nccFvp.origPatch().patch();
87  const polyPatch& nbrOrigPp = nbrNccFvp.origPatch().patch();
88 
89  // Get the indices of the non-conformal patches
90  labelList patchis(Pstream::nProcs(), -1);
91  labelList nbrPatchis(Pstream::nProcs(), -1);
92  patchis[Pstream::myProcNo()] = nccFvp.index();
93  nbrPatchis[Pstream::myProcNo()] = nbrNccFvp.index();
94  forAll(mesh_.boundary(), patchj)
95  {
96  const fvPatch& fvp = mesh_.boundary()[patchj];
97 
98  if (isA<nonConformalProcessorCyclicFvPatch>(fvp))
99  {
100  const nonConformalProcessorCyclicFvPatch& ncpcFvp =
101  refCast<const nonConformalProcessorCyclicFvPatch>(fvp);
102 
103  if (ncpcFvp.referPatchID() == nccFvp.index())
104  {
105  patchis[ncpcFvp.neighbProcNo()] = patchj;
106  }
107  if (ncpcFvp.referPatchID() == nbrNccFvp.index())
108  {
109  nbrPatchis[ncpcFvp.neighbProcNo()] = patchj;
110  }
111  }
112  }
113 
114  // Get the intersection geometry
115  const patchToPatches::intersection& intersection =
116  nccFvp.nonConformalCyclicPatch().intersection();
117 
118  // Unpack the patchToPatch addressing into lists of indices (fixed lists of
119  // 3 labels; owner face, neighbour face, couple index). These will be used
120  // to create the non-conformal faces, so sort them to make sure the
121  // non-conformal interfaces are ordered.
122  auto procFacesToIndices = []
123  (
124  const List<List<remote>>& faceOtherProcFaces,
125  const bool owner
126  )
127  {
128  // Compute interface sizes
129  labelList is(Pstream::nProcs(), 0);
130  forAll(faceOtherProcFaces, facei)
131  {
132  forAll(faceOtherProcFaces[facei], i)
133  {
134  const remote& otherProcFacei = faceOtherProcFaces[facei][i];
135 
136  is[otherProcFacei.proci] ++;
137  }
138  }
139 
140  // Allocate the indices
141  List<List<FixedList<label, 3>>> indices(Pstream::nProcs());
142  forAll(indices, proci)
143  {
144  indices[proci].resize(is[proci]);
145  }
146 
147  // Set the indices
148  is = 0;
149  forAll(faceOtherProcFaces, facei)
150  {
151  forAll(faceOtherProcFaces[facei], i)
152  {
153  const remote& otherProcFacei = faceOtherProcFaces[facei][i];
154 
155  FixedList<label, 3>& index =
156  indices[otherProcFacei.proci][is[otherProcFacei.proci] ++];
157 
158  index = {facei, otherProcFacei.elementi, i};
159 
160  if (!owner) Swap(index[0], index[1]);
161  }
162  }
163 
164  // Sort to ensure ordering
165  forAll(indices, proci)
166  {
167  sort(indices[proci]);
168  }
169 
170  return indices;
171  };
172 
173  List<List<FixedList<label, 3>>> indices =
174  procFacesToIndices(intersection.srcTgtProcFaces(), true);
175 
176  List<List<FixedList<label, 3>>> nbrIndices =
177  procFacesToIndices(intersection.tgtSrcProcFaces(), false);
178 
179  // If addressing has been provided, then modify the indices to match. When
180  // a coupling has to be added, the couple index is set to -1. This
181  // indicates that there is no geometry in the patchToPatch engine for this
182  // coupling, and for a small stabilisation value to be used instead.
183  auto matchIndices = [&polyFacesBf, &tOrigFacesNbrBf]
184  (
185  List<List<FixedList<label, 3>>>& indices,
186  const nonConformalCyclicFvPatch& nccFvp,
187  const polyPatch& origPp,
188  const labelList& patchis,
189  const bool owner
190  )
191  {
192  const surfaceLabelField::Boundary& origFacesNbrBf = tOrigFacesNbrBf();
193 
194  // Create what the indices should be
195  List<List<FixedList<label, 3>>> indicesRef(indices.size());
196  forAll(indices, proci)
197  {
198  const label patchi = patchis[proci];
199 
200  if (patchi != -1)
201  {
202  indicesRef[proci].resize(polyFacesBf[patchi].size());
203 
204  forAll(polyFacesBf[patchi], patchFacei)
205  {
206  FixedList<label, 3>& indexRef =
207  indicesRef[proci][patchFacei];
208 
209  indexRef =
210  {
211  polyFacesBf[patchi][patchFacei] - origPp.start(),
212  origFacesNbrBf[patchi][patchFacei],
213  -1
214  };
215 
216  if (!owner) Swap(indexRef[0], indexRef[1]);
217  }
218  }
219  }
220 
221  // Populate the indices with the index of the coupling
222  label nCouplesRemoved = 0, nCouplesAdded = 0;
223  forAll(indices, proci)
224  {
225  const label patchi = patchis[proci];
226 
227  if (patchi != -1)
228  {
229  label patchFacei = 0, i = 0;
230 
231  while
232  (
233  patchFacei < polyFacesBf[patchi].size()
234  && i < indices[proci].size()
235  )
236  {
237  const FixedList<label, 3> index
238  ({
239  indices[proci][i][0],
240  indices[proci][i][1],
241  -1
242  });
243 
244  FixedList<label, 3>& indexRef =
245  indicesRef[proci][patchFacei];
246 
247  if (index < indexRef)
248  {
249  nCouplesRemoved ++;
250  i ++;
251  }
252  else if (index == indexRef)
253  {
254  indexRef[2] = indices[proci][i][2];
255  patchFacei ++;
256  i ++;
257  }
258  else // (index > indexRef)
259  {
260  nCouplesAdded ++;
261  patchFacei ++;
262  }
263  }
264 
265  nCouplesRemoved +=
266  min(indices[proci].size() - i, 0);
267  nCouplesAdded +=
268  min(polyFacesBf[patchi].size() - patchFacei, 0);
269  }
270  }
271 
272  // Report if changes have been made
273  reduce(nCouplesRemoved, sumOp<label>());
274  reduce(nCouplesAdded, sumOp<label>());
275  if ((nCouplesRemoved || nCouplesAdded) && nccFvp.owner())
276  {
277  Info<< indent << nCouplesRemoved << '/' << nCouplesAdded
278  << " small couplings removed/added to " << nccFvp.name()
279  << endl;
280  }
281 
282  // Set the indices to the correct values
283  Swap(indices, indicesRef);
284  };
285 
286  if (tOrigFacesNbrBf.valid())
287  {
288  matchIndices(indices, nccFvp, origPp, patchis, true);
289  matchIndices(nbrIndices, nbrNccFvp, nbrOrigPp, nbrPatchis, false);
290  }
291 
292  // Create couplings by transferring geometry from the original to the
293  // non-conformal patches
294  auto createCouplings =
295  [
296  &polyFacesBf,
297  &tOrigSfNbrBf,
298  &tOrigCfNbrBf,
299  &SfBf,
300  &CfBf
301  ]
302  (
303  const List<List<FixedList<label, 3>>>& indices,
304  const List<DynamicList<couple>>& couples,
305  const nonConformalCyclicFvPatch& nccFvp,
306  const polyPatch& origPp,
307  const labelList& patchis,
308  const bool owner
309  )
310  {
311  forAll(patchis, proci)
312  {
313  const label patchi = patchis[proci];
314  const label patchSize = indices[proci].size();
315 
316  if (patchi == -1 && patchSize)
317  {
319  << "Intersection generated coupling through "
320  << nccFvp.name() << " between processes "
321  << Pstream::myProcNo() << " and " << proci
322  << " but a corresponding processor interface was not found"
323  << exit(FatalError);
324  }
325 
326  if (patchi != -1)
327  {
328  polyFacesBf[patchi].resize(patchSize);
329  SfBf[patchi].resize(patchSize);
330  CfBf[patchi].resize(patchSize);
331 
332  forAll(indices[proci], patchFacei)
333  {
334  const label origFacei = indices[proci][patchFacei][!owner];
335  const label i = indices[proci][patchFacei][2];
336 
337  polyFacesBf[patchi][patchFacei] =
338  origFacei + origPp.start();
339 
340  couple c;
341  if (i != -1)
342  {
343  c = couples[origFacei][i];
344  }
345  else
346  {
347  c =
348  couple
349  (
350  part
351  (
352  small*origPp.faceAreas()[origFacei],
353  origPp.faceCentres()[origFacei]
354  ),
355  part
356  (
357  small*tOrigSfNbrBf()[patchi][patchFacei],
358  tOrigCfNbrBf()[patchi][patchFacei]
359  )
360  );
361  }
362 
363  if (!owner)
364  {
365  c.nbr = c;
366  }
367 
368  SfBf[patchi][patchFacei] = c.nbr.area;
369  CfBf[patchi][patchFacei] = c.nbr.centre;
370 
371  if (i != -1)
372  {
373  part origP
374  (
375  SfBf[origPp.index()][origFacei],
376  CfBf[origPp.index()][origFacei]
377  );
378  origP -= c;
379 
380  SfBf[origPp.index()][origFacei] = origP.area;
381  CfBf[origPp.index()][origFacei] = origP.centre;
382  }
383  }
384  }
385  }
386  };
387 
388  createCouplings
389  (
390  indices,
391  intersection.srcCouples(),
392  nccFvp,
393  origPp,
394  patchis,
395  true
396  );
397 
398  createCouplings
399  (
400  nbrIndices,
401  intersection.tgtCouples(),
402  nbrNccFvp,
403  nbrOrigPp,
404  nbrPatchis,
405  false
406  );
407 
408  // Add error geometry to the original patches
409  forAll(intersection.srcErrorParts(), origFacei)
410  {
411  part origP
412  (
413  SfBf[origPp.index()][origFacei],
414  CfBf[origPp.index()][origFacei]
415  );
416  origP += intersection.srcErrorParts()[origFacei];
417 
418  SfBf[origPp.index()][origFacei] = origP.area;
419  CfBf[origPp.index()][origFacei] = origP.centre;
420  }
421 
422  // Store the orig edge parts
423  forAll(intersection.srcEdgeParts(), origEdgei)
424  {
425  origEdgeParts[origEdgei] += intersection.srcEdgeParts()[origEdgei];
426  }
427 }
428 
429 
431 Foam::fvMeshStitcher::calculateOwnerOrigBoundaryEdgeParts
432 (
433  const List<List<part>>& patchEdgeParts
434 ) const
435 {
436  const polyBoundaryMesh& pbMesh = mesh_.boundaryMesh();
437 
438  const nonConformalBoundary& ncb = nonConformalBoundary::New(mesh_);
439  const labelList& ownerOrigBoundaryPointMeshPoint =
440  ncb.ownerOrigBoundaryPointMeshPoint();
441  const labelList& ownerOrigBoundaryEdgeMeshEdge =
442  ncb.ownerOrigBoundaryEdgeMeshEdge();
443  const edgeList& ownerOrigBoundaryEdges = ncb.ownerOrigBoundaryEdges();
444  const edgeList& ownerOrigBoundaryMeshEdges =
445  ncb.ownerOrigBoundaryMeshEdges();
446 
447  // Sum up boundary edge parts, being careful with signs
448  labelList ownerOrigBoundaryEdgeNParts(ownerOrigBoundaryEdges.size(), 0);
449  List<part> ownerOrigBoundaryEdgeParts(ownerOrigBoundaryEdges.size());
450  forAll(patchEdgeParts, patchi)
451  {
452  if (patchEdgeParts[patchi].empty()) continue;
453 
454  const polyPatch& patch = pbMesh[patchi];
455 
456  const labelList& patchEdgeOwnerOrigBoundaryEdges =
457  ncb.patchEdgeOwnerOrigBoundaryEdges(patchi);
458 
459  forAll(patchEdgeParts[patchi], patchEdgei)
460  {
461  const label ownerOrigBoundaryEdgei =
462  patchEdgeOwnerOrigBoundaryEdges[patchEdgei];
463 
464  const label sign =
466  (
467  meshEdge(patch, patchEdgei),
468  ownerOrigBoundaryMeshEdges[ownerOrigBoundaryEdgei]
469  );
470 
471  const part& pU = patchEdgeParts[patchi][patchEdgei];
472  const part p = sign > 0 ? pU : -pU;
473 
474  ownerOrigBoundaryEdgeNParts[ownerOrigBoundaryEdgei] ++;
475  ownerOrigBoundaryEdgeParts[ownerOrigBoundaryEdgei] += p;
476  }
477  }
478 
479  // Give every point in the all boundary a global unique index
480  const globalIndex globalOwnerOrigBoundaryPoints
481  (
482  ownerOrigBoundaryPointMeshPoint.size()
483  );
484  labelList ownerOrigBoundaryPointIndices
485  (
486  ownerOrigBoundaryPointMeshPoint.size()
487  );
488  forAll(ownerOrigBoundaryPointIndices, ownerOrigBoundaryPointi)
489  {
490  ownerOrigBoundaryPointIndices[ownerOrigBoundaryPointi] =
491  globalOwnerOrigBoundaryPoints.toGlobal(ownerOrigBoundaryPointi);
492  }
494  (
495  mesh_,
496  ownerOrigBoundaryPointMeshPoint,
497  ownerOrigBoundaryPointIndices,
498  minEqOp<label>(),
499  labelMax
500  );
501 
502  // Synchronise the sign of the edge parts so that they can be
503  // meaningfully summed. This is done using the sign of the global point
504  // indices; if they are not in order, then reverse the part area.
505  forAll(ownerOrigBoundaryEdgeParts, ownerOrigBoundaryEdgei)
506  {
507  const edge& e = ownerOrigBoundaryEdges[ownerOrigBoundaryEdgei];
508 
509  if
510  (
511  ownerOrigBoundaryPointIndices[e.start()]
512  > ownerOrigBoundaryPointIndices[e.end()])
513  {
514  ownerOrigBoundaryEdgeParts[ownerOrigBoundaryEdgei] =
515  - ownerOrigBoundaryEdgeParts[ownerOrigBoundaryEdgei];
516  }
517  }
518 
519  // Average the boundary edge parts across couplings
521  (
522  mesh_,
523  ownerOrigBoundaryEdgeMeshEdge,
524  ownerOrigBoundaryEdgeNParts,
525  plusEqOp<label>(),
526  label(0)
527  );
529  (
530  mesh_,
531  ownerOrigBoundaryEdgeMeshEdge,
532  ownerOrigBoundaryEdgeParts,
533  plusEqOp<part>(),
534  part(),
535  []
536  (
537  const transformer& vt,
538  const bool forward,
539  List<part>& fld
540  )
541  {
542  if (forward)
543  {
544  forAll(fld, i)
545  {
546  if (magSqr(fld[i].area) != 0)
547  {
548  fld[i].area = vt.transform(fld[i].area);
549  fld[i].centre = vt.transformPosition(fld[i].centre);
550  }
551  }
552  }
553  else
554  {
555  forAll(fld, i)
556  {
557  if (magSqr(fld[i].area) != 0)
558  {
559  fld[i].area = vt.invTransform(fld[i].area);
560  fld[i].centre = vt.invTransformPosition(fld[i].centre);
561  }
562  }
563  }
564  }
565  );
566  forAll(ownerOrigBoundaryEdgeParts, ownerOrigBoundaryEdgei)
567  {
568  if (ownerOrigBoundaryEdgeNParts[ownerOrigBoundaryEdgei] != 0)
569  {
570  ownerOrigBoundaryEdgeParts[ownerOrigBoundaryEdgei].area /=
571  ownerOrigBoundaryEdgeNParts[ownerOrigBoundaryEdgei];
572  }
573  }
574 
575  // Put the edge parts back to their original orientations
576  forAll(ownerOrigBoundaryEdgeParts, ownerOrigBoundaryEdgei)
577  {
578  const edge& e = ownerOrigBoundaryEdges[ownerOrigBoundaryEdgei];
579 
580  if
581  (
582  ownerOrigBoundaryPointIndices[e.start()]
583  > ownerOrigBoundaryPointIndices[e.end()]
584  )
585  {
586  ownerOrigBoundaryEdgeParts[ownerOrigBoundaryEdgei] =
587  - ownerOrigBoundaryEdgeParts[ownerOrigBoundaryEdgei];
588  }
589  }
590 
591  return ownerOrigBoundaryEdgeParts;
592 }
593 
594 
595 
596 void Foam::fvMeshStitcher::applyOwnerOrigBoundaryEdgeParts
597 (
598  surfaceVectorField& SfSf,
599  surfaceVectorField& CfSf,
600  const List<part>& ownerOrigBoundaryEdgeParts
601 ) const
602 {
603  const polyBoundaryMesh& pbMesh = mesh_.boundaryMesh();
604 
605  const nonConformalBoundary& ncb = nonConformalBoundary::New(mesh_);
606  const labelList ownerOrigPatchIDs = ncb.ownerOrigPatchIDs();
607  const labelList& ownerOrigBoundaryEdgeMeshEdge =
608  ncb.ownerOrigBoundaryEdgeMeshEdge();
609  const edgeList& ownerOrigBoundaryMeshEdges =
610  ncb.ownerOrigBoundaryMeshEdges();
611 
612  boolList patchIsOwnerOrig(pbMesh.size(), false);
613  UIndirectList<bool>(patchIsOwnerOrig, ownerOrigPatchIDs) = true;
614 
615  // Count the number of original faces attached to each boundary edge
616  labelList ownerOrigBoundaryEdgeNOrigFaces
617  (
618  ownerOrigBoundaryEdgeMeshEdge.size(),
619  0
620  );
621  forAll(ownerOrigBoundaryEdgeMeshEdge, ownerOrigBoundaryEdgei)
622  {
623  const label meshEdgei =
624  ownerOrigBoundaryEdgeMeshEdge[ownerOrigBoundaryEdgei];
625 
626  forAll(mesh_.edgeFaces()[meshEdgei], edgeFacei)
627  {
628  const label facei = mesh_.edgeFaces()[meshEdgei][edgeFacei];
629 
630  const label patchi =
631  mesh_.isInternalFace(facei)
632  ? -1 : pbMesh.patchID()[facei - mesh_.nInternalFaces()];
633 
634  if (patchi != -1 && patchIsOwnerOrig[patchi])
635  {
636  ownerOrigBoundaryEdgeNOrigFaces[ownerOrigBoundaryEdgei] ++;
637  }
638  }
639  }
640 
641  // Synchronise the boundary edge original face count
643  (
644  mesh_,
645  ownerOrigBoundaryEdgeMeshEdge,
646  ownerOrigBoundaryEdgeNOrigFaces,
647  plusEqOp<label>(),
648  label(0)
649  );
650 
651  // If debugging, check that face changes are in sync
652  tmp<surfaceLabelField> tChanged =
653  debug
654  ? surfaceLabelField::New("changed", mesh_, dimensioned<label>(dimless, 0))
655  : tmp<surfaceLabelField>(nullptr);
656 
657  // Modify faces connected to the boundary edges
658  forAll(ownerOrigBoundaryEdgeMeshEdge, ownerOrigBoundaryEdgei)
659  {
660  const label meshEdgei =
661  ownerOrigBoundaryEdgeMeshEdge[ownerOrigBoundaryEdgei];
662 
663  const part& ownerOrigBoundaryEdgeP =
664  ownerOrigBoundaryEdgeParts[ownerOrigBoundaryEdgei];
665 
666  switch (ownerOrigBoundaryEdgeNOrigFaces[ownerOrigBoundaryEdgei])
667  {
668  case 0:
669  {
670  continue;
671  }
672 
673  // If this edge is connected to just one original face then add any
674  // edge area into that face
675  case 1:
676  {
677  label origFacei = -1, origPatchi = -1, origPatchFacei = -1;
678  forAll(mesh_.edgeFaces()[meshEdgei], edgeFacei)
679  {
680  const label facei = mesh_.edgeFaces()[meshEdgei][edgeFacei];
681 
682  const label patchi =
683  mesh_.isInternalFace(facei)
684  ? -1 : pbMesh.patchID()[facei - mesh_.nInternalFaces()];
685 
686  if (patchi != -1 && patchIsOwnerOrig[patchi])
687  {
688  origFacei = facei;
689  origPatchi = patchi;
690  origPatchFacei = facei - pbMesh[patchi].start();
691  break;
692  }
693  }
694 
695  if (origFacei != -1)
696  {
697  vector& area =
698  SfSf.boundaryFieldRef()[origPatchi][origPatchFacei];
699  point& centre =
700  CfSf.boundaryFieldRef()[origPatchi][origPatchFacei];
701 
702  const label sign =
703  mesh_.faces()[origFacei].edgeDirection
704  (
705  ownerOrigBoundaryMeshEdges[ownerOrigBoundaryEdgei]
706  );
707 
708  part p(area, centre);
709  p +=
710  sign > 0
711  ? ownerOrigBoundaryEdgeP
712  : -ownerOrigBoundaryEdgeP;
713 
714  area = p.area;
715  centre = p.centre;
716 
717  if (debug)
718  {
719  tChanged->boundaryFieldRef()
720  [origPatchi][origPatchFacei] = label(1);
721  }
722  }
723 
724  break;
725  }
726 
727  // If this edge is connected to two original faces then add any
728  // edge area into non-original connected faces
729  case 2:
730  {
731  forAll(mesh_.edgeFaces()[meshEdgei], edgeFacei)
732  {
733  const label facei = mesh_.edgeFaces()[meshEdgei][edgeFacei];
734 
735  const label patchi =
736  mesh_.isInternalFace(facei)
737  ? -1 : pbMesh.patchID()[facei - mesh_.nInternalFaces()];
738 
739  if (patchi != -1 && patchIsOwnerOrig[patchi])
740  {
741  continue;
742  }
743 
744  const label patchFacei =
745  patchi == -1 ? -1 : facei - pbMesh[patchi].start();
746 
747  vector& area =
748  patchi == -1
749  ? SfSf[facei]
750  : SfSf.boundaryFieldRef()[patchi][patchFacei];
751  point& centre =
752  patchi == -1
753  ? CfSf[facei]
754  : CfSf.boundaryFieldRef()[patchi][patchFacei];
755 
756  const label sign =
757  mesh_.faces()[facei].edgeDirection
758  (
759  ownerOrigBoundaryMeshEdges[ownerOrigBoundaryEdgei]
760  );
761 
762  part p(area, centre);
763  p +=
764  sign < 0
765  ? ownerOrigBoundaryEdgeP
766  : -ownerOrigBoundaryEdgeP;
767 
768  area = p.area;
769  centre = p.centre;
770 
771  if (debug && patchi != -1)
772  {
773  tChanged->boundaryFieldRef()
774  [patchi][patchFacei] = label(1);
775  }
776  }
777 
778  break;
779  }
780 
781  default:
782  {
784  << "Boundary is not manifold"
785  << exit(FatalError);
786  }
787  }
788  }
789 
790  // If debugging, check that face changes are in sync
791  if (debug)
792  {
793  const surfaceLabelField::Boundary& changedBf =
794  tChanged->boundaryField();
795 
796  const surfaceLabelField::Boundary changedNbrBf
797  (
799  changedBf.boundaryNeighbourField()
800  );
801 
802  bool synchronised = true;
803 
804  forAll(changedBf, patchi)
805  {
806  const fvsPatchLabelField& changedPf = changedBf[patchi];
807  const fvsPatchLabelField& changedNbrPf = changedNbrBf[patchi];
808 
809  if (!isA<processorFvPatch>(changedPf.patch())) continue;
810 
811  forAll(changedPf, patchFacei)
812  {
813  if (changedPf[patchFacei] != changedNbrPf[patchFacei])
814  {
815  const label facei =
816  changedPf.patch().start() + patchFacei;
817 
819  << "Face not synchronised with centre at "
820  << mesh_.faceCentres()[facei] << endl;
821 
822  synchronised = false;
823  }
824  }
825  }
826 
827  if (!synchronised)
828  {
830  << exit(FatalError);
831  }
832  }
833 }
834 
835 
836 void Foam::fvMeshStitcher::stabiliseOrigPatchFaces
837 (
840 ) const
841 {
842  const nonConformalBoundary& ncb = nonConformalBoundary::New(mesh_);
843  const labelList allOrigPatchIDs = ncb.allOrigPatchIDs();
844 
845  forAll(allOrigPatchIDs, i)
846  {
847  const label origPatchi = allOrigPatchIDs[i];
848  const polyPatch& origPp = mesh_.boundaryMesh()[origPatchi];
849 
850  forAll(origPp, origPatchFacei)
851  {
852  part p
853  (
854  SfBf[origPatchi][origPatchFacei],
855  CfBf[origPatchi][origPatchFacei]
856  );
857 
858  const part smallP
859  (
860  small*origPp.faceAreas()[origPatchFacei],
861  origPp.faceCentres()[origPatchFacei]
862  );
863 
864  p += smallP;
865 
866  SfBf[origPatchi][origPatchFacei] = p.area;
867  CfBf[origPatchi][origPatchFacei] = p.centre;
868  }
869  }
870 }
871 
872 
873 void Foam::fvMeshStitcher::intersectNonConformalCyclics
874 (
875  surfaceLabelField::Boundary& polyFacesBf,
876  surfaceVectorField& SfSf,
877  surfaceVectorField& CfSf,
878  const bool haveTopology
879 ) const
880 {
881  const polyBoundaryMesh& pbMesh = mesh_.boundaryMesh();
882 
883  const nonConformalBoundary& ncb = nonConformalBoundary::New(mesh_);
884  const labelList ownerOrigPatchIDs = ncb.ownerOrigPatchIDs();
885 
886  // Alias the boundary geometry fields
887  surfaceVectorField::Boundary& SfBf = SfSf.boundaryFieldRef();
888  surfaceVectorField::Boundary& CfBf = CfSf.boundaryFieldRef();
889 
890  // Create storage for and initialise the edge parts of source patches
891  List<List<part>> patchEdgeParts(mesh_.boundary().size());
892  forAll(ownerOrigPatchIDs, i)
893  {
894  const label origPatchi = ownerOrigPatchIDs[i];
895 
896  patchEdgeParts[origPatchi].resize
897  (
898  mesh_.boundaryMesh()[origPatchi].nEdges(),
899  part(Zero)
900  );
901  }
902 
903  // If topology is specified, also create a boundary field of the original
904  // patch indices on the neighbouring interface, as well as the
905  // corresponding areas and centres for stabilisation purposes
906  tmp<surfaceLabelField::Boundary> tOrigFacesNbrBf;
907  tmp<surfaceVectorField::Boundary> tOrigSfNbrBf;
908  tmp<surfaceVectorField::Boundary> tOrigCfNbrBf;
909  if (haveTopology)
910  {
911  surfaceLabelField origFaces
912  (
914  (
915  "origFaces",
916  mesh_,
917  dimensioned<label>(dimless, -1)
918  )
919  );
920  surfaceVectorField origSf("origSf", mesh_.Sf());
921  surfaceVectorField origCf("origCf", mesh_.Cf());
922 
923  forAll(mesh_.boundary(), patchi)
924  {
925  const fvPatch& fvp = mesh_.boundary()[patchi];
926 
927  if (!isA<nonConformalFvPatch>(fvp)) continue;
928 
929  const nonConformalFvPatch& ncFvp =
930  refCast<const nonConformalFvPatch>(fvp);
931 
932  origFaces.boundaryFieldRef()[patchi] =
933  polyFacesBf[patchi] - ncFvp.origPatch().start();
934  origSf.boundaryFieldRef()[patchi] =
935  vectorField(mesh_.faceAreas(), polyFacesBf[patchi]);
936  origCf.boundaryFieldRef()[patchi] =
937  pointField(mesh_.faceCentres(), polyFacesBf[patchi]);
938  }
939 
940  tOrigFacesNbrBf =
942  (
944  origFaces.boundaryField().boundaryNeighbourField()
945  );
946  tOrigSfNbrBf =
948  (
950  origSf.boundaryField().boundaryNeighbourField()
951  );
952  tOrigCfNbrBf =
954  (
956  origCf.boundaryField().boundaryNeighbourField()
957  );
958 
959  // See note in fvMesh::unconform regarding the transformation
960  // of cell centres. The same applies here.
961  forAll(tOrigCfNbrBf(), patchi)
962  {
963  fvsPatchVectorField& Cfp = tOrigCfNbrBf.ref()[patchi];
964  if (Cfp.patch().coupled())
965  {
966  const transformer& t =
967  refCast<const coupledFvPatch>(Cfp.patch())
968  .transform();
969  t.invTransform(Cfp, Cfp);
970  t.transformPosition(Cfp, Cfp);
971  }
972  }
973  }
974 
975  // Do the intersections
976  forAll(mesh_.boundary(), patchi)
977  {
978  const fvPatch& fvp = mesh_.boundary()[patchi];
979 
980  if (!isA<nonConformalCyclicFvPatch>(fvp)) continue;
981 
982  const nonConformalCyclicFvPatch& nccFvp =
983  refCast<const nonConformalCyclicFvPatch>(fvp);
984 
985  if (!nccFvp.owner()) continue;
986 
987  intersectNonConformalCyclic
988  (
989  nccFvp,
990  polyFacesBf,
991  SfBf,
992  CfBf,
993  tOrigFacesNbrBf,
994  tOrigSfNbrBf,
995  tOrigCfNbrBf,
996  patchEdgeParts[nccFvp.origPatchID()]
997  );
998  }
999 
1000  // Construct the boundary edge geometry
1001  List<part> ownerOrigBoundaryEdgeParts =
1002  calculateOwnerOrigBoundaryEdgeParts(patchEdgeParts);
1003 
1004  // Add the difference between patch edge parts and all boundary
1005  // edge parts to the adjacent patch faces. This is an error part.
1006  forAll(ownerOrigPatchIDs, i)
1007  {
1008  const label origPatchi = ownerOrigPatchIDs[i];
1009  const polyPatch& origPatch = pbMesh[origPatchi];
1010 
1011  const labelList& origPatchEdgeOwnerOrigBoundaryEdges =
1012  ncb.patchEdgeOwnerOrigBoundaryEdges(origPatchi);
1013 
1014  forAll(patchEdgeParts[origPatchi], origPatchEdgei)
1015  {
1016  const label ownerOrigBoundaryEdgei =
1017  origPatchEdgeOwnerOrigBoundaryEdges[origPatchEdgei];
1018 
1019  part errorP =
1020  patchEdgeParts[origPatchi][origPatchEdgei];
1021  errorP -= ownerOrigBoundaryEdgeParts[ownerOrigBoundaryEdgei];
1022 
1023  forAll(origPatch.edgeFaces()[origPatchEdgei], patchEdgeFacei)
1024  {
1025  const label patchFacei =
1026  origPatch.edgeFaces()[origPatchEdgei][patchEdgeFacei];
1027 
1028  part p
1029  (
1030  SfBf[origPatchi][patchFacei],
1031  CfBf[origPatchi][patchFacei]
1032  );
1033  p += errorP;
1034 
1035  SfBf[origPatchi][patchFacei] = p.area;
1036  CfBf[origPatchi][patchFacei] = p.centre;
1037  }
1038  }
1039  }
1040 
1041  // Use the boundary edge geometry to correct all edge-connected faces
1042  applyOwnerOrigBoundaryEdgeParts(SfSf, CfSf, ownerOrigBoundaryEdgeParts);
1043 
1044  // Stabilise
1045  stabiliseOrigPatchFaces(SfBf, CfBf);
1046 }
1047 
1048 
1049 template<class NonConformalFvPatch>
1050 inline void Foam::fvMeshStitcher::createNonConformalStabilisationGeometry
1051 (
1052  const surfaceLabelField::Boundary& polyFacesBf,
1053  surfaceVectorField& SfSf,
1054  surfaceVectorField& CfSf
1055 ) const
1056 {
1057  // Alias the boundary geometry fields
1060 
1061  // Create small stabilisation geometry
1062  forAll(mesh_.boundary(), patchi)
1063  {
1064  const fvPatch& fvp = mesh_.boundary()[patchi];
1065 
1066  if (!isA<NonConformalFvPatch>(fvp)) continue;
1067 
1068  const polyPatch& origPp =
1069  refCast<const nonConformalFvPatch>(fvp).origPatch().patch();
1070 
1071  SfBf[patchi] ==
1072  vectorField
1073  (
1074  small*origPp.faceAreas(),
1075  polyFacesBf[patchi] - origPp.start()
1076  );
1077  CfBf[patchi] ==
1078  vectorField
1079  (
1080  origPp.faceCentres(),
1081  polyFacesBf[patchi] - origPp.start()
1082  );
1083  }
1084 }
1085 
1086 
1087 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
1088 
1090 {
1091  bool result = false;
1092 
1093  forAll(mesh_.boundary(), patchi)
1094  {
1095  const fvPatch& fvp = mesh_.boundary()[patchi];
1096 
1097  if (!isA<nonConformalFvPatch>(fvp)) continue;
1098 
1099  const fvsPatchScalarField& magSfp =
1100  mesh_.magSf().boundaryField()[patchi];
1101 
1102  const polyPatch& origPp =
1103  refCast<const nonConformalFvPatch>(fvp).origPatch().patch();
1104 
1105  const scalarField origMagSfp
1106  (
1107  origPp.magFaceAreas(),
1108  mesh_.polyFacesBf()[patchi] - origPp.start()
1109  );
1110 
1111  if (max(magSfp/origMagSfp) > rootSmall)
1112  {
1113  result = true;
1114  }
1115  }
1116 
1117  reduce(result, orOp<bool>());
1118 
1119  return returnReduce(result, orOp<bool>());
1120 }
1121 
1122 
1125 {
1126  return
1127  mag(fvc::surfaceIntegrate(mesh_.Sf()))()*mesh_.V()
1128  /max
1129  (
1130  mag(fvc::surfaceSum(cmptMag(mesh_.Sf())))(),
1131  (small*sqr(cbrt(mesh_.V())))()
1132  )();
1133 }
1134 
1135 
1138 {
1139  if (0 > n || n > 1)
1140  {
1142  << "Can only compute volume conservation error for this time, or "
1143  << "the previous time" << exit(FatalError);
1144  }
1145 
1146  const surfaceScalarField& phi = mesh_.phi().oldTime(n);
1147 
1148  const dimensionedScalar deltaT =
1149  n == 0 ? mesh_.time().deltaT() : mesh_.time().deltaT0();
1150 
1151  const volScalarField::Internal& V = n == 0 ? mesh_.V() : mesh_.V0();
1152  const volScalarField::Internal& V0 = n == 0 ? mesh_.V0() : mesh_.V00();
1153 
1154  return fvc::surfaceIntegrate(phi*deltaT)() - (V - V0)/mesh_.V();
1155 }
1156 
1157 
1158 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1159 
1161 :
1162  mesh_(mesh)
1163 {}
1164 
1165 
1166 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1167 
1169 {}
1170 
1171 
1172 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1173 
1175 {
1176  forAll(mesh_.boundary(), patchi)
1177  {
1178  if (isA<nonConformalFvPatch>(mesh_.boundary()[patchi]))
1179  {
1180  return true;
1181  }
1182  }
1183 
1184  return false;
1185 }
1186 
1187 
1189 {}
1190 
1191 
1193 {}
1194 
1195 
1197 (
1198  const bool changing,
1199  const bool geometric
1200 )
1201 {
1202  if (!stitches() || (changing && !mesh_.dynamic()))
1203  {
1204  return false;
1205  }
1206 
1207  // Is this mesh coupled?
1208  const bool coupled = Pstream::parRun() || !mesh_.time().processorCase();
1209 
1210  if (coupled && geometric)
1211  {
1212  Info<< indent << typeName << ": Disconnecting" << incrIndent << endl;
1213  }
1214 
1215  if (changing)
1216  {
1217  // Pre-conform surface fields. This splits the original and cyclic
1218  // parts of the interface fields into separate boundary fields, with
1219  // both sets of values store on the original faces. The original field
1220  // overwrites the existing boundary values, whilst the cyclic field is
1221  // stored as a separate field for use later.
1222  preConformSurfaceFields();
1223  }
1224 
1225  // Conform the mesh
1226  if (mesh_.moving())
1227  {
1228  surfaceScalarField phi(mesh_.phi());
1229  conformCorrectMeshPhi(phi);
1230  mesh_.conform(phi);
1231  }
1232  else
1233  {
1234  mesh_.conform();
1235  }
1236 
1237  // Resize all the affected patch fields
1238  resizePatchFields<VolField>();
1239  resizePatchFields<SurfaceField>();
1240 
1241  // Prevent hangs caused by processor cyclic patches using mesh geometry
1242  mesh_.deltaCoeffs();
1243 
1244  if (coupled && geometric)
1245  {
1246  const volScalarField::Internal o(openness());
1247  Info<< indent << "Cell min/avg/max openness = "
1248  << gMin(o) << '/' << gAverage(o) << '/' << gMax(o) << endl;
1249 
1250  for (label i = 0; mesh_.moving() && i <= mesh_.phi().nOldTimes(); ++ i)
1251  {
1252  const volScalarField::Internal vce(volumeConservationError(i));
1253  Info<< indent << "Cell min/avg/max ";
1254  for (label j = 0; j < i; ++ j) Info<< "old-";
1255  Info<< (i ? "time " : "") << "volume conservation error = "
1256  << gMin(vce) << '/' << gAverage(vce) << '/' << gMax(vce)
1257  << endl;
1258  }
1259 
1260  }
1261 
1262  if (coupled && geometric)
1263  {
1264  Info<< decrIndent;
1265  }
1266 
1267  // Create null polyTopoChangeMap
1268  const polyTopoChangeMap map(mesh_);
1269 
1270  meshObjects::topoChange<fvMesh>(mesh_, map);
1271  meshObjects::topoChange<lduMesh>(mesh_, map);
1272 
1273  const_cast<Time&>(mesh_.time()).functionObjects().topoChange(map);
1274 
1275  return true;
1276 }
1277 
1278 
1280 (
1281  const bool changing,
1282  const bool geometric,
1283  const bool load
1284 )
1285 {
1286  if (!stitches() || (changing && !mesh_.dynamic()))
1287  {
1288  return false;
1289  }
1290 
1291  // Is this mesh coupled?
1292  const bool coupled = Pstream::parRun() || !mesh_.time().processorCase();
1293 
1294  // Create a copy of the conformal poly face addressing
1295  surfaceLabelField::Boundary polyFacesBf
1296  (
1298  mesh_.polyFacesBf()
1299  );
1300 
1301  // If starting up then load topology from disk, if it is available
1302  bool haveTopology = false;
1303  if (load)
1304  {
1305  const IOobject polyFacesBfIO(mesh_.polyFacesBfIO(IOobject::MUST_READ));
1306 
1307  if (fileHandler().isFile(polyFacesBfIO.objectPath(false)))
1308  {
1309  polyFacesBf.reset
1310  (
1311  surfaceLabelField(polyFacesBfIO, mesh_).boundaryField()
1312  );
1313 
1314  haveTopology = true;
1315  }
1316  }
1317 
1318  // Access all the intersections in advance. Makes the log nicer.
1319  if (coupled && (geometric || !haveTopology))
1320  {
1321  forAll(mesh_.boundary(), patchi)
1322  {
1323  const fvPatch& fvp = mesh_.boundary()[patchi];
1324 
1325  if (!isA<nonConformalCyclicFvPatch>(fvp)) continue;
1326 
1327  const nonConformalCyclicFvPatch& nccFvp =
1328  refCast<const nonConformalCyclicFvPatch>(fvp);
1329 
1330  if (!nccFvp.owner()) continue;
1331 
1333  }
1334  }
1335 
1336  if (coupled && (geometric || !haveTopology))
1337  {
1338  Info<< indent << typeName << ": Connecting" << incrIndent << endl;
1339  }
1340 
1341  // Create copies of geometry fields to be modified
1342  surfaceVectorField Sf(mesh_.Sf().cloneUnSliced()());
1343  surfaceVectorField Cf(mesh_.Cf().cloneUnSliced()());
1344 
1345  // Construct non-conformal geometry
1346  if (coupled && (geometric || !haveTopology))
1347  {
1348  // Do the intersection and create the non-conformal cyclic faces
1349  intersectNonConformalCyclics(polyFacesBf, Sf, Cf, haveTopology);
1350 
1351  // Create stabilisation geometry on any error faces specified in the
1352  // addressing
1353  createNonConformalStabilisationGeometry<nonConformalErrorFvPatch>
1354  (
1355  polyFacesBf,
1356  Sf,
1357  Cf
1358  );
1359  }
1360  else
1361  {
1362  // If not coupled, then create stabilisation geometry on all
1363  // non-conformal faces specified in the addressing
1364  createNonConformalStabilisationGeometry<nonConformalFvPatch>
1365  (
1366  polyFacesBf,
1367  Sf,
1368  Cf
1369  );
1370  }
1371 
1372  // If the mesh is moving then create any additional non-conformal geometry
1373  // necessary to correct the mesh fluxes
1374  if (mesh_.moving())
1375  {
1376  createNonConformalCorrectMeshPhiGeometry(polyFacesBf, Sf, Cf);
1377  }
1378 
1379  // Make the mesh non-conformal
1380  if (mesh_.moving())
1381  {
1382  surfaceScalarField phi(mesh_.phi());
1383  unconformCorrectMeshPhi(polyFacesBf, Sf, Cf, phi);
1384  mesh_.unconform(polyFacesBf, Sf, Cf, phi);
1385  }
1386  else
1387  {
1388  mesh_.unconform(polyFacesBf, Sf, Cf);
1389  }
1390 
1391  // Resize all the affected patch fields
1392  resizePatchFields<VolField>();
1393  resizePatchFields<SurfaceField>();
1394 
1395  if (changing)
1396  {
1397  // Post-non-conform surface fields. This reconstructs the original and
1398  // cyclic parts of the interface fields from separate original and
1399  // cyclic parts. The original part was store in the same field, whilst
1400  // the cyclic part was separately registered.
1401  postNonConformSurfaceFields();
1402 
1403  // Volume fields are assumed to be intensive. So, the value on a face
1404  // which has changed in size can be retained without modification. New
1405  // faces need values to be set. This is done by evaluating all the
1406  // nonConformalCoupled patch fields.
1407  evaluateVolFields();
1408 
1409  // Do special post-non-conformation for surface velocities.
1410  postNonConformSurfaceVelocities();
1411  }
1412 
1413  // Prevent hangs caused by processor cyclic patches using mesh geometry
1414  mesh_.deltaCoeffs();
1415 
1416  if (coupled && geometric)
1417  {
1418  const volScalarField::Internal o(openness());
1419  Info<< indent << "Cell min/avg/max openness = "
1420  << gMin(o) << '/' << gAverage(o) << '/' << gMax(o) << endl;
1421 
1422  for (label i = 0; mesh_.moving() && i <= mesh_.phi().nOldTimes(); ++ i)
1423  {
1424  const volScalarField::Internal vce(volumeConservationError(i));
1425  Info<< indent << "Cell min/avg/max ";
1426  for (label j = 0; j < i; ++ j) Info<< "old-";
1427  Info<< (i ? "time " : "") << "volume conservation error = "
1428  << gMin(vce) << '/' << gAverage(vce) << '/' << gMax(vce)
1429  << endl;
1430  }
1431 
1432  if (mesh_.moving())
1433  {
1434  // Create a boundary field of the imbalance between the mesh fluxes
1435  // on either side of interfaces (like phiErrorb above)
1437  (
1439  mesh_.phi().boundaryField()
1440  );
1441  mfe += mesh_.phi().boundaryField().boundaryNeighbourField();
1442 
1443  // Determine the number of non-processor patches
1444  label nNonProcPatches = 0;
1445  forAll(mesh_.boundary(), patchi)
1446  {
1447  const fvPatch& fvp = mesh_.boundary()[patchi];
1448 
1449  if (!isA<processorFvPatch>(fvp))
1450  {
1451  if (nNonProcPatches != patchi)
1452  {
1454  << "Processor patches do not follow non-processor "
1455  << "patches" << exit(FatalError);
1456  }
1457 
1458  nNonProcPatches ++;
1459  }
1460  }
1461 
1462  // Work out the min, avg and max error for all non-conformal
1463  // cyclics, including any errors on referred processor cyclics
1464  scalarList minMfe(nNonProcPatches, vGreat);
1465  scalarField sumMfe(nNonProcPatches, 0), nSumMfe(nNonProcPatches, 0);
1466  scalarList maxMfe(nNonProcPatches, -vGreat);
1467  forAll(mfe, patchi)
1468  {
1469  const fvPatch& fvp = mesh_.boundary()[patchi];
1470 
1471  const label nccPatchi =
1472  isA<nonConformalCyclicFvPatch>(fvp)
1473  ? refCast<const nonConformalCyclicFvPatch>(fvp)
1474  .index()
1475  : isA<nonConformalProcessorCyclicFvPatch>(fvp)
1476  ? refCast<const nonConformalProcessorCyclicFvPatch>(fvp)
1477  .referPatchID()
1478  : -1;
1479 
1480  if (nccPatchi != -1)
1481  {
1482  minMfe[nccPatchi] =
1483  min(minMfe[nccPatchi], min(mfe[patchi]));
1484  sumMfe[nccPatchi] += sum(mfe[patchi]);
1485  nSumMfe[nccPatchi] += mfe[patchi].size();
1486  maxMfe[nccPatchi] =
1487  max(maxMfe[nccPatchi], max(mfe[patchi]));
1488  }
1489  }
1490  reduce(minMfe, ListOp<minOp<scalar>>());
1491  reduce(sumMfe, ListOp<minOp<scalar>>());
1492  reduce(nSumMfe, ListOp<minOp<scalar>>());
1493  reduce(maxMfe, ListOp<maxOp<scalar>>());
1494 
1495  // Report
1496  forAll(minMfe, patchi)
1497  {
1498  const fvPatch& fvp = mesh_.boundary()[patchi];
1499 
1500  if
1501  (
1502  isA<nonConformalCyclicFvPatch>(fvp)
1503  && refCast<const nonConformalCyclicFvPatch>(fvp).owner()
1504  && nSumMfe[patchi] > 0
1505  )
1506  {
1507  Info<< indent << fvp.name()
1508  << " min/avg/max mesh flux error = " << minMfe[patchi]
1509  << '/' << sumMfe[patchi]/(nSumMfe[patchi] + vSmall)
1510  << '/' << maxMfe[patchi] << endl;
1511  }
1512  }
1513  }
1514  }
1515 
1516  if (coupled && (geometric || !haveTopology))
1517  {
1518  Info<< decrIndent;
1519  }
1520 
1521  // Create null polyTopoChangeMap
1522  const polyTopoChangeMap map(mesh_);
1523 
1524  meshObjects::topoChange<fvMesh>(mesh_, map);
1525  meshObjects::topoChange<lduMesh>(mesh_, map);
1526 
1527  const_cast<Time&>(mesh_.time()).functionObjects().topoChange(map);
1528 
1529  return true;
1530 }
1531 
1532 
1533 void Foam::fvMeshStitcher::reconnect(const bool geometric) const
1534 {
1535  if (mesh_.conformal() || geometric == this->geometric())
1536  {
1537  return;
1538  }
1539 
1540  // Create a copy of the conformal poly face addressing
1541  surfaceLabelField::Boundary polyFacesBf
1542  (
1544  mesh_.polyFacesBf()
1545  );
1546 
1547  // Undo all non-conformal changes and clear all geometry and topology
1548  mesh_.conform();
1549 
1550  // Create copies of geometry fields to be modified
1551  surfaceVectorField Sf(mesh_.Sf().cloneUnSliced()());
1552  surfaceVectorField Cf(mesh_.Cf().cloneUnSliced()());
1553 
1554  // Construct non-conformal geometry
1555  if (geometric)
1556  {
1557  // Do the intersection and create the non-conformal cyclic faces
1558  intersectNonConformalCyclics(polyFacesBf, Sf, Cf, true);
1559 
1560  // Create stabilisation geometry on any error faces specified in the
1561  // addressing
1562  createNonConformalStabilisationGeometry<nonConformalErrorFvPatch>
1563  (
1564  polyFacesBf,
1565  Sf,
1566  Cf
1567  );
1568  }
1569  else
1570  {
1571  // If not coupled, then create stabilisation geometry on all
1572  // non-conformal faces specified in the addressing
1573  createNonConformalStabilisationGeometry<nonConformalFvPatch>
1574  (
1575  polyFacesBf,
1576  Sf,
1577  Cf
1578  );
1579  }
1580 
1581  // Apply changes to the mesh
1582  mesh_.unconform(polyFacesBf, Sf, Cf);
1583 
1584  // Prevent hangs caused by processor cyclic patches using mesh geometry
1585  mesh_.deltaCoeffs();
1586 
1587  // Create null polyTopoChangeMap
1588  const polyTopoChangeMap map(mesh_);
1589 
1590  meshObjects::topoChange<fvMesh>(mesh_, map);
1591  meshObjects::topoChange<lduMesh>(mesh_, map);
1592 
1593  const_cast<Time&>(mesh_.time()).functionObjects().topoChange(map);
1594 }
1595 
1596 
1597 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static const DimensionedField< Type, GeoMesh > & null()
Return a null DimensionedField.
Generic GeometricBoundaryField class.
tmp< GeometricBoundaryField > boundaryNeighbourField() const
Return BoundaryField of the values on the other side of couples.
void reset(const GeometricBoundaryField< Type, PatchField, GeoMesh > &)
Reset the boundary field contents to the given field.
Generic GeometricField class.
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of the boundary field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
fileName objectPath(const bool global) const
Return complete path + object name including the processor.
Definition: IOobject.H:411
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.H:294
A list of faces which address into the list of points.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
static int compare(const edge &, const edge &)
Compare edges.
Definition: edgeI.H:36
Mesh manipulator that uses the intersection provided by the cyclic non-conformal poly patches to crea...
virtual bool disconnect(const bool changing, const bool geometric)
Disconnect the mesh by removing faces from the nonConformalCyclics.
virtual void updateMesh(const polyTopoChangeMap &)
Update local data for topology changes.
virtual bool connect(const bool changing, const bool geometric, const bool load)
Connect the mesh by adding faces into the nonConformalCyclics.
void reconnect(const bool geometric) const
Re-compute the connection. Topology is preserved. Permits a change.
bool geometric() const
Is the connection "geometric", or has the topology just been.
fvMeshStitcher(fvMesh &)
Construct from fvMesh.
virtual void movePoints()
Update local data for mesh motion.
bool stitches() const
Does this stitcher do anything?
tmp< DimensionedField< scalar, volMesh > > openness() const
Return the non-dimensional cell openness for debugging/checking.
patchToPatches::intersection::part part
Alias the patchToPatch intersection part struct.
tmp< DimensionedField< scalar, volMesh > > volumeConservationError(const label n) const
Return the non-dimensional old-time volume conservation error.
patchToPatches::intersection::couple couple
Alias the patchToPatch intersection couple struct.
virtual ~fvMeshStitcher()
Destructor.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:893
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:175
virtual const word & name() const
Return name.
Definition: fvPatch.H:145
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:82
Non-conformal cyclic FV patch. As nonConformalCoupledFvPatch, but the neighbouring patch is local and...
virtual bool owner() const
...
const nonConformalCyclicPolyPatch & nonConformalCyclicPatch() const
Poly patch.
const patchToPatches::intersection & intersection() const
Access the intersection engine.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
const vectorField::subField faceAreas() const
Return face areas.
Definition: polyPatch.C:282
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:276
const scalarField::subField magFaceAreas() const
Return face area magnitudes.
Definition: polyPatch.C:288
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronise values on all mesh points.
static void syncEdgeList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronise values on all mesh edges.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
label patchi
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const dimensionedScalar c
Speed of light in a vacuum.
tmp< VolField< Type > > surfaceSum(const SurfaceField< Type > &ssf)
void surfaceIntegrate(Field< Type > &ivf, const SurfaceField< Type > &ssf)
Namespace for OpenFOAM.
bool isFile(const fileName &, const bool checkVariants=true, const bool followLink=true)
Does the name exist as a file in the file system?
Definition: POSIX.C:555
const fileOperation & fileHandler()
Get current file handler.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
const doubleScalar e
Definition: doubleScalar.H:105
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
List< label > labelList
A List of labels.
Definition: labelList.H:56
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
dimensionedScalar sign(const dimensionedScalar &ds)
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
const dimensionSet dimless
messageStream Info
SurfaceField< label > surfaceLabelField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
vector point
Point is a vector.
Definition: point.H:41
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
fvsPatchField< label > fvsPatchLabelField
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
void sort(UList< T > &)
Definition: UList.C:115
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
edge meshEdge(const PrimitivePatch< FaceList, PointField > &p, const label edgei)
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensionSet cmptMag(const dimensionSet &)
Definition: dimensionSet.C:275
Type gAverage(const FieldField< Field, Type > &f)
error FatalError
void Swap(T &a, T &b)
Definition: Swap.H:43
List< edge > edgeList
Definition: edgeList.H:38
dimensionedScalar cbrt(const dimensionedScalar &ds)
Type gMin(const FieldField< Field, Type > &f)
static const label labelMax
Definition: label.H:62
SurfaceField< vector > surfaceVectorField
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
fvsPatchField< vector > fvsPatchVectorField
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Type gMax(const FieldField< Field, Type > &f)
volScalarField & p
Operator to apply a binary operation to a pair of lists.
Definition: ListOps.H:284
Return the vol-field velocity corresponding to a given surface-field velocity.