fvMeshStitchersMoving.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-2024 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 "fvMeshStitchersMoving.H"
27 #include "FvFaceCellWave.H"
28 #include "fvMeshSubset.H"
29 #include "fvmLaplacian.H"
30 #include "meshPhiCorrectInfo.H"
31 #include "meshPhiPreCorrectInfo.H"
32 #include "nonConformalBoundary.H"
36 #include "regionSplit.H"
37 #include "solutionControl.H"
38 #include "syncTools.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 template<class Type>
48 inline List<Type> repeat(const UList<Type>& a, const UList<Type>& b)
49 {
50  List<Type> result(a.size() + b.size());
51  forAll(a, i)
52  {
53  result[2*i] = a[i];
54  result[2*i + 1] = b[i];
55  }
56  return result;
57 }
58 
59 
60 template<class Type>
61 inline List<Type> repeat(const UList<Type>& l)
62 {
63  return repeat(l, l);
64 }
65 
66 
68 {
70  scalar weight;
71 
72  static const layerAndWeight min, max;
73 
74  typedef nil cmptType;
75 
76  friend bool operator==(const layerAndWeight& a, const layerAndWeight& b)
77  {
78  return a.layer == b.layer && a.weight == b.weight;
79  }
80 
81  friend bool operator!=(const layerAndWeight& a, const layerAndWeight& b)
82  {
83  return !(a == b);
84  }
85 
86  friend Ostream& operator<<(Ostream& os, const layerAndWeight& l)
87  {
88  return os << l.layer << token::SPACE << l.weight;
89  }
90 
92  {
93  return is >> l.layer >> l.weight;
94  }
95 };
96 
97 
98 const layerAndWeight layerAndWeight::min({-labelMax, NaN});
99 
100 
101 const layerAndWeight layerAndWeight::max({labelMax, NaN});
102 
103 
105 {
106  return a.layer > b.layer ? a : b;
107 }
108 
109 
111 {
112  return a.layer < b.layer ? a : b;
113 }
114 
115 }
116 
117 
118 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
119 
120 namespace Foam
121 {
122 namespace fvMeshStitchers
123 {
126 }
127 }
128 
129 
130 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
131 
132 void Foam::fvMeshStitchers::moving::conformCorrectMeshPhi
133 (
134  surfaceScalarField& phi
135 )
136 {
137  // Add the non-conformal parts of the mesh flux into the original faces
138  surfaceScalarField::Boundary& phiBf = phi.boundaryFieldRef();
139  forAll(phiBf, nccPatchi)
140  {
141  if (isA<nonConformalFvPatch>(phiBf[nccPatchi].patch()))
142  {
143  const nonConformalFvPatch& ncFvp =
144  refCast<const nonConformalFvPatch>(phiBf[nccPatchi].patch());
145 
146  const label origPatchi = ncFvp.origPatchIndex();
147  const fvPatch& origFvp = ncFvp.origPatch();
148 
149  for (label i = 0; i <= phi.nOldTimes(false); ++ i)
150  {
151  phi.oldTimeRef(i).boundaryFieldRef()[origPatchi] +=
153  (
154  phi.oldTime(i).boundaryField()[nccPatchi],
155  origFvp.size(),
156  ncFvp.polyFaces(),
157  origFvp.start()
158  );
159 
160  phi.oldTimeRef(i).boundaryFieldRef()[nccPatchi].clear();
161  }
162  }
163  }
164 }
165 
166 
167 void Foam::fvMeshStitchers::moving::createNonConformalCorrectMeshPhiGeometry
168 (
169  surfaceLabelField::Boundary& polyFacesBf,
170  surfaceVectorField& SfSf,
171  surfaceVectorField& CfSf
172 )
173 {
174  const nonConformalBoundary& ncb = nonConformalBoundary::New(mesh());
175  const labelList origPatchIndices = ncb.allOrigPatchIndices();
176  const labelList errorPatchIndices = ncb.allErrorPatchIndices();
177 
178  forAll(errorPatchIndices, i)
179  {
180  const label origPatchi = origPatchIndices[i];
181  const polyPatch& origPp = mesh().boundaryMesh()[origPatchi];
182 
183  const label errorPatchi = errorPatchIndices[i];
184 
185  polyFacesBf[errorPatchi] =
186  repeat((identityMap(origPp.size()) + origPp.start())());
187 
188  SfSf.boundaryFieldRef()[errorPatchi] =
189  repeat
190  (
191  (rootVSmall*origPp.faceNormals())(),
192  (-rootVSmall*origPp.faceNormals())()
193  );
194  CfSf.boundaryFieldRef()[errorPatchi] =
195  repeat(origPp.faceCentres());
196  }
197 }
198 
199 
200 Foam::labelHashSet Foam::fvMeshStitchers::moving::ownerCoupledCellSet()
201 {
202  const nonConformalBoundary& ncb = nonConformalBoundary::New(mesh());
203 
204  // For every boundary face, count how many owner-orig faces it is connected
205  // to (will be 1 or 0)
206  labelList bFaceNSet(mesh().nFaces() - mesh().nInternalFaces(), 0);
207  forAll(mesh().boundary(), nccPatchi)
208  {
209  const fvPatch& fvp = mesh().boundary()[nccPatchi];
210 
211  if (!isA<nonConformalCoupledFvPatch>(fvp)) continue;
212 
213  const nonConformalCoupledFvPatch& nccFvp =
214  refCast<const nonConformalCoupledFvPatch>(fvp);
215 
216  if (!nccFvp.owner()) continue;
217 
218  forAll(nccFvp, nccPatchFacei)
219  {
220  bFaceNSet
221  [
222  mesh().polyFacesBf()[nccPatchi][nccPatchFacei]
223  - mesh().nInternalFaces()
224  ] = 1;
225  }
226  }
227 
228  // For every boundary edge, count how many owner-orig faces it is connected
229  // to (should be 0, 1 or 2)
230  labelList ownerOrigBoundaryEdgeNSet
231  (
232  ncb.ownerOrigBoundaryEdgeMeshEdge().size(),
233  0
234  );
235  forAll(ncb.ownerOrigBoundaryEdgeMeshEdge(), ownerOrigBoundaryEdgei)
236  {
237  const label meshEdgei =
238  ncb.ownerOrigBoundaryEdgeMeshEdge()[ownerOrigBoundaryEdgei];
239 
240  forAll(mesh().edgeFaces()[meshEdgei], edgeFacei)
241  {
242  const label facei = mesh().edgeFaces()[meshEdgei][edgeFacei];
243 
244  if (mesh().isInternalFace(facei)) continue;
245 
246  const label bFacei = facei - mesh().nInternalFaces();
247 
248  if (bFaceNSet[bFacei])
249  {
250  ownerOrigBoundaryEdgeNSet[ownerOrigBoundaryEdgei] += 1;
251  }
252  }
253  }
255  (
256  mesh(),
257  ncb.ownerOrigBoundaryEdgeMeshEdge(),
258  ownerOrigBoundaryEdgeNSet,
259  plusEqOp<label>(),
260  label(0)
261  );
262 
263  // Isolate the cells that are edge connected to the owner patches. This
264  // will form the sub-mesh in which the correction is computed.
265  labelHashSet set;
266  forAll(ncb.ownerOrigBoundaryEdgeMeshEdge(), ownerOrigBoundaryEdgei)
267  {
268  if (ownerOrigBoundaryEdgeNSet[ownerOrigBoundaryEdgei] < 2) continue;
269 
270  const label meshEdgei =
271  ncb.ownerOrigBoundaryEdgeMeshEdge()[ownerOrigBoundaryEdgei];
272 
273  forAll(mesh().edgeFaces()[meshEdgei], edgeFacei)
274  {
275  const label facei = mesh().edgeFaces()[meshEdgei][edgeFacei];
276 
277  set.insert(mesh().faceOwner()[facei]);
278 
279  if (facei < mesh().nInternalFaces())
280  {
281  set.insert(mesh().faceNeighbour()[facei]);
282  }
283  }
284  }
285 
286  return set;
287 }
288 
289 
290 void Foam::fvMeshStitchers::moving::unconformInternalFaceCorrectMeshPhi
291 (
292  surfaceScalarField& phi
293 )
294 {
295  const surfaceScalarField::Boundary& magSfBf =
296  mesh().magSf().boundaryField();
297 
298  surfaceScalarField::Boundary& phiBf = phi.boundaryFieldRef();
299 
300  // Step 1: Construct some boundary information
301 
302  // For each boundary face mark if it is an owner orig face and sum the total
303  // non-conformal coupled area if so
304  boolList bFaceIsOwnerOrig(mesh().nFaces() - mesh().nInternalFaces(), false);
305  scalarList bFaceNccMagSf(mesh().nFaces() - mesh().nInternalFaces(), Zero);
306  forAll(mesh().boundary(), nccPatchi)
307  {
308  const fvPatch& fvp = mesh().boundary()[nccPatchi];
309 
310  if (!isA<nonConformalCoupledFvPatch>(fvp)) continue;
311 
312  const nonConformalCoupledFvPatch& nccFvp =
313  refCast<const nonConformalCoupledFvPatch>(fvp);
314 
315  if (!nccFvp.owner()) continue;
316 
317  forAll(nccFvp, nccPatchFacei)
318  {
319  const label bFacei =
320  nccFvp.polyFaces()[nccPatchFacei] - mesh().nInternalFaces();
321 
322  bFaceIsOwnerOrig[bFacei] = true;
323  bFaceNccMagSf[bFacei] += magSfBf[nccPatchi][nccPatchFacei];
324  }
325  }
326 
327  // Step 2: Construct a sub-mesh for all cells that are connected to the
328  // owner patches.
329 
330  // Isolate the cells that are edge connected to the owner patches. This
331  // will form the sub-mesh in which the correction is computed.
332  labelHashSet subCellSet = ownerCoupledCellSet();
333 
334  // Create a mesh for the cells that are edge connected to coupled faces
335  // of the owner patches. This is where the correction will be computed.
336  fvMeshSubset subsetter(mesh());
337  subsetter.setLargeCellSubset(ownerCoupledCellSet());
338  const fvMesh& subMesh = subsetter.subMesh();
339  subMesh.deltaCoeffs();
340 
341  // Determine the disconnected regions of the sub mesh
342  const regionSplit subMeshRegions(subMesh);
343  const label subNRegions = subMeshRegions.nRegions();
344 
345  // Map from mesh boundary face to sub-mesh region
346  labelList bFaceSubRegion(mesh().nFaces() - mesh().nInternalFaces(), -1);
347  forAll(subsetter.faceMap(), subFacei)
348  {
349  const label facei = subsetter.faceMap()[subFacei];
350 
351  if (mesh().isInternalFace(facei)) continue;
352 
353  const label bFacei = facei - mesh().nInternalFaces();
354 
355  bFaceSubRegion[bFacei] =
356  subMeshRegions[subMesh.faceOwner()[subFacei]];
357  }
358 
359  // Get a single reference cell for each region
360  labelList subMeshRegionRefCells(subNRegions, -1);
361  {
362  static const label proci = Pstream::myProcNo();
363 
364  labelList subMeshRegionRefProcs(subNRegions, labelMax);
365  forAll(subMeshRegions, subCelli)
366  {
367  subMeshRegionRefProcs[subMeshRegions[subCelli]] = proci;
368  }
369  reduce(subMeshRegionRefProcs, ListOp<minOp<label>>());
370 
371  forAll(subMeshRegions, subCelli)
372  {
373  if
374  (
375  subMeshRegionRefProcs[subMeshRegions[subCelli]] == proci
376  && subMeshRegionRefCells[subMeshRegions[subCelli]] == -1
377  )
378  {
379  subMeshRegionRefCells[subMeshRegions[subCelli]] = subCelli;
380  }
381  }
382  }
383 
384  // Step 3: Synchronise the flux across the interfaces
385 
386  // Create a synchronised mesh flux for the coupled patches. Set both sides
387  // to the neighbour value so that all the error is on the owner side.
389  (
392  );
393 
394  // Determine the total mesh flux error and area magnitude for each region
395  scalarList regionPhiError(subNRegions, scalar(0));
396  scalarList regionMagSf(subNRegions, vSmall);
397  forAll(phiBf, nccPatchi)
398  {
399  const fvPatch& fvp = mesh().boundary()[nccPatchi];
400 
401  if (!isA<nonConformalCoupledFvPatch>(fvp)) continue;
402 
403  const nonConformalCoupledFvPatch& nccFvp =
404  refCast<const nonConformalCoupledFvPatch>(fvp);
405 
406  if (!nccFvp.owner()) continue;
407 
408  forAll(nccFvp, nccPatchFacei)
409  {
410  const label subRegioni =
411  bFaceSubRegion
412  [
413  mesh().polyFacesBf()[nccPatchi][nccPatchFacei]
414  - mesh().nInternalFaces()
415  ];
416 
417  regionPhiError[subRegioni] +=
418  syncPhiBf[nccPatchi][nccPatchFacei]
419  - phiBf[nccPatchi][nccPatchFacei];
420 
421  regionMagSf[subRegioni] +=
422  magSfBf[nccPatchi][nccPatchFacei];
423  }
424  }
425  reduce(regionPhiError, ListOp<sumOp<scalar>>());
426  reduce(regionMagSf, ListOp<sumOp<scalar>>());
427 
428  // Synchronise the mesh fluxes, but offset so that the total flux for each
429  // region is the same as for the non-synchronised mesh fluxes
430  forAll(phiBf, nccPatchi)
431  {
432  const fvPatch& fvp = mesh().boundary()[nccPatchi];
433 
434  if (!isA<nonConformalCoupledFvPatch>(fvp)) continue;
435 
436  const nonConformalCoupledFvPatch& nccFvp =
437  refCast<const nonConformalCoupledFvPatch>(fvp);
438 
439  if (!nccFvp.owner()) continue;
440 
441  forAll(nccFvp, nccPatchFacei)
442  {
443  const label subRegioni =
444  bFaceSubRegion
445  [
446  mesh().polyFacesBf()[nccPatchi][nccPatchFacei]
447  - mesh().nInternalFaces()
448  ];
449 
450  phiBf[nccPatchi][nccPatchFacei] =
451  syncPhiBf[nccPatchi][nccPatchFacei]
452  - magSfBf[nccPatchi][nccPatchFacei]
453  /regionMagSf[subRegioni]
454  *regionPhiError[subRegioni];
455  }
456  }
457 
458  // Step 4: Set up a system on the sub-mesh with which to solve for a flux
459  // that corrects the volume conservation error in the cells connected to
460  // synchronised faces
461 
462  // Map volumes to the sub mesh
463  const volScalarField::Internal subV(subsetter.interpolate(mesh().V()));
464  const volScalarField::Internal subV0(subsetter.interpolate(mesh().V0()));
465 
466  // Map mesh flux to the sub mesh, accumulating when a face is split into
467  // multiple non-conformal parts
468  surfaceScalarField subPhi
469  (
471  (
472  "phi",
473  subMesh,
475  )
476  );
477  forAll(subPhi, subFacei)
478  {
479  const label facei = subsetter.faceMap()[subFacei];
480 
481  subPhi[subFacei] = phi[facei];
482  }
483  forAll(subPhi.boundaryField(), subPatchi)
484  {
485  const fvPatch& subFvp = subPhi.boundaryField()[subPatchi].patch();
486 
487  forAll(subPhi.boundaryField()[subPatchi], subPatchFacei)
488  {
489  const label subFacei = subFvp.start() + subPatchFacei;
490  const label facei = subsetter.faceMap()[subFacei];
491 
492  if (mesh().isInternalFace(facei))
493  {
494  const label s = sign(subsetter.faceFlipMap()[subFacei]);
495 
496  subPhi.boundaryFieldRef()[subPatchi][subPatchFacei] =
497  s*phi[facei];
498  }
499  else
500  {
501  const label bFacei = facei - mesh().nInternalFaces();
502 
503  const labelUList patches =
504  mesh().polyBFacePatches()[bFacei];
505  const labelUList patchFaces =
506  mesh().polyBFacePatchFaces()[bFacei];
507 
508  forAll(patches, i)
509  {
510  subPhi.boundaryFieldRef()[subPatchi][subPatchFacei] +=
511  phiBf[patches[i]][patchFaces[i]];
512  }
513  }
514  }
515  }
516 
517  // Calculate the volume conservation error for the sub mesh
518  const volScalarField::Internal subVce
519  (
520  fvc::surfaceIntegrate(subPhi*subMesh.time().deltaT())()
521  - (subV - subV0)/subV
522  );
523 
524  // Construct boundary conditions for the sub-mesh potential. Zero gradient
525  // is used for all non-constraint boundaries (i.e., this is a closed
526  // domain) and for exposed internal faces, which therefore require a patch
527  // type override
528  wordList MeshPhiPatchTypes(subMesh.boundary().size());
529  forAll(subMesh.boundary(), patchi)
530  {
531  const fvPatch& subFvp = subMesh.boundary()[patchi];
532  MeshPhiPatchTypes[patchi] =
533  polyPatch::constraintType(subFvp.type())
534  && !isA<internalFvPatch>(subFvp)
535  ? subFvp.type()
537  }
538 
539  // Solve
540  volScalarField MeshPhi
541  (
542  IOobject
543  (
544  "MeshPhi",
545  subMesh.time().name(),
546  subMesh,
549  ),
550  subMesh,
552  MeshPhiPatchTypes,
553  subMesh.boundaryMesh().types()
554  );
555 
556  subMesh.schemes().setFluxRequired(MeshPhi.name());
557 
558  fvScalarMatrix MeshPhiEqn
559  (
560  fvm::laplacian(MeshPhi) + subVce
561  );
562 
563  forAll(subMeshRegionRefCells, i)
564  {
565  MeshPhiEqn.setReference(subMeshRegionRefCells[i], 0);
566  }
567 
568  MeshPhiEqn.solve();
569 
570  // Create a correction for the sub-mesh face fluxes
571  surfaceScalarField subDeltaPhi(MeshPhiEqn.flux()/mesh().time().deltaT());
572 
573  // Step 5: Currently, the computed correction fixes the volume conservation
574  // error to the level of the solve tolerance. We would rather the volume
575  // conservation error be at round off and the error be in the mismatch
576  // across the couplings, as the mismatch will be fixed later when area and
577  // mesh flux is added to the error faces. So, here, we correct the
578  // correction (!) by propagating the error from the cells furthest from the
579  // interface back to the interface. This is done using waves.
580 
581  // Dynamic memory for wave initialisation
582  DynamicList<labelPair> subChangedPatchAndFaces;
583  DynamicList<meshPhiPreCorrectInfo> subChangedFacePci;
584  DynamicList<meshPhiCorrectInfo> subChangedFaceCi;
585 
586  // Allocate data for the pre-correction wave
587  List<meshPhiPreCorrectInfo> subInternalFacePci(subMesh.nInternalFaces());
588  List<List<meshPhiPreCorrectInfo>> subPatchFacePci
589  (
590  FvFaceCellWave<meshPhiPreCorrectInfo>::template
591  sizesListList<List<List<meshPhiPreCorrectInfo>>>
592  (
593  FvFaceCellWave<meshPhiPreCorrectInfo>::template
594  listListSizes(subMesh.boundary()),
595  meshPhiPreCorrectInfo()
596  )
597  );
598  List<meshPhiPreCorrectInfo> subCellPci(subMesh.nCells());
599 
600  // Initialisation for the pre-correction wave
601  subChangedPatchAndFaces.clear();
602  subChangedFacePci.clear();
603  forAll(subMesh.boundary(), subPatchi)
604  {
605  const fvPatch& subFvp = subMesh.boundary()[subPatchi];
606 
607  forAll(subFvp, subPatchFacei)
608  {
609  const label subFacei = subFvp.start() + subPatchFacei;
610  const label facei = subsetter.faceMap()[subFacei];
611  const label bFacei = facei - mesh().nInternalFaces();
612 
613  if (bFacei >= 0 && bFaceIsOwnerOrig[bFacei])
614  {
615  subChangedPatchAndFaces.append({subPatchi, subPatchFacei});
616  subChangedFacePci.append
617  (
618  meshPhiPreCorrectInfo(0, bFaceNccMagSf[bFacei])
619  );
620  }
621  }
622  }
623 
624  // Pre-correction wave
625  FvFaceCellWave<meshPhiPreCorrectInfo> preWave
626  (
627  subMesh,
628  subInternalFacePci,
629  subPatchFacePci,
630  subCellPci
631  );
632  preWave.setFaceInfo(subChangedPatchAndFaces, subChangedFacePci);
633  const label nWaveLayers =
634  preWave.iterate(subMesh.globalData().nTotalCells() + 1);
635 
636  // Allocate data for the correction wave
637  List<meshPhiCorrectInfo> subInternalFaceCi(subMesh.nInternalFaces());
638  List<List<meshPhiCorrectInfo>> subPatchFaceCi
639  (
640  FvFaceCellWave<meshPhiCorrectInfo>::template
641  sizesListList<List<List<meshPhiCorrectInfo>>>
642  (
643  FvFaceCellWave<meshPhiCorrectInfo>::template
644  listListSizes(subMesh.boundary()),
645  meshPhiCorrectInfo()
646  )
647  );
648  List<meshPhiCorrectInfo> subCellCi(subMesh.nCells());
649 
650  // Calculate the current error in the rate of change of volume
651  const volScalarField::Internal dVdtError
652  (
653  (subV - subV0)/subMesh.time().deltaT()
654  - fvc::surfaceIntegrate(subPhi + subDeltaPhi)()*subV
655  );
656 
657  // Construct track data for the correction wave
658  meshPhiCorrectInfo::trackData td
659  (
660  subInternalFacePci,
661  subPatchFacePci,
662  subCellPci,
663  dVdtError
664  );
665 
666  // Wave backwards through the layers to generate the corrections. Note that
667  // this has to be done in stages, so that later layers complete in their
668  // entirety before the earlier layers begin. Otherwise they interfere.
669  for (label waveLayeri = nWaveLayers - 1; waveLayeri >= 0; waveLayeri --)
670  {
671  // The layer indices on the faces that we want to wave from
672  const label faceLayeri = (waveLayeri + 1)*2;
673 
674  // Initialisation for the correction wave
675  subChangedPatchAndFaces.clear();
676  subChangedFaceCi.clear();
677  forAll(subInternalFacePci, subFacei)
678  {
679  if (subInternalFacePci[subFacei].layer() == faceLayeri)
680  {
681  subChangedPatchAndFaces.append({-1, subFacei});
682  subChangedFaceCi.append
683  (
684  subInternalFaceCi[subFacei].valid(td)
685  ? subInternalFaceCi[subFacei]
686  : meshPhiCorrectInfo(meshPhiCorrectInfo::shape::face)
687  );
688  }
689  }
690  forAll(subPatchFacePci, subPatchi)
691  {
692  forAll(subPatchFacePci[subPatchi], subPatchFacei)
693  {
694  if
695  (
696  subPatchFacePci[subPatchi][subPatchFacei].layer()
697  == faceLayeri
698  )
699  {
700  subChangedPatchAndFaces.append({subPatchi, subPatchFacei});
701  subChangedFaceCi.append
702  (
703  subPatchFaceCi[subPatchi][subPatchFacei].valid(td)
704  ? subPatchFaceCi[subPatchi][subPatchFacei]
705  : meshPhiCorrectInfo(meshPhiCorrectInfo::shape::face)
706  );
707  }
708  }
709  }
710 
711  // Correction wave
712  FvFaceCellWave<meshPhiCorrectInfo, meshPhiCorrectInfo::trackData> wave
713  (
714  subMesh,
715  subInternalFaceCi,
716  subPatchFaceCi,
717  subCellCi,
718  td
719  );
720  wave.setFaceInfo(subChangedPatchAndFaces, subChangedFaceCi);
721  wave.iterate(1);
722  }
723 
724  // Apply corrections
725  forAll(subInternalFaceCi, subFacei)
726  {
727  subDeltaPhi.primitiveFieldRef()[subFacei] +=
728  subInternalFaceCi[subFacei].deltaPhi();
729  }
730  forAll(subPatchFacePci, subPatchi)
731  {
732  forAll(subPatchFacePci[subPatchi], subPatchFacei)
733  {
734  subDeltaPhi.boundaryFieldRef()[subPatchi][subPatchFacei] +=
735  subPatchFaceCi[subPatchi][subPatchFacei].deltaPhi();
736  }
737  }
738 
739  // Step 6: Apply the corrections to the mesh flux
740 
741  // Correct the internal mesh face fluxes
742  forAll(subPhi, subFacei)
743  {
744  phi[subsetter.faceMap()[subFacei]] =
745  subPhi[subFacei] + subDeltaPhi[subFacei];
746  }
747 
748  // Map the sub-mesh flux changes to the conformal mesh boundary
750  (
751  mesh().boundary(),
754  );
755  deltaPhiBf = 0;
756  forAll(subMesh.boundary(), subPatchi)
757  {
758  const label patchi = subsetter.patchMap()[subPatchi];
759 
760  if (patchi == -1) continue;
761 
762  const fvPatch& subFvp = subMesh.boundary()[subPatchi];
763  const fvPatch& fvp = mesh().boundary()[patchi];
764 
765  const bool coupled = subFvp.coupled();
766 
767  forAll(subMesh.boundary()[subPatchi], subPatchFacei)
768  {
769  const label facei =
770  subsetter.faceMap()[subFvp.start() + subPatchFacei];
771 
772  const label patchFacei = facei - fvp.start();
773 
774  if (coupled)
775  {
776  deltaPhiBf[patchi][patchFacei] =
777  subPhi.boundaryField()[subPatchi][subPatchFacei]
778  + subDeltaPhi.boundaryField()[subPatchi][subPatchFacei]
779  - phiBf[patchi][patchFacei];
780  }
781  else
782  {
783  deltaPhiBf[patchi][patchFacei] =
784  subDeltaPhi.boundaryField()[subPatchi][subPatchFacei];
785  }
786  }
787  }
788 
789  // Move the changes from the owner-orig to the non-conformal coupled faces
790  forAll(mesh().boundary(), nccPatchi)
791  {
792  const fvPatch& fvp = mesh().boundary()[nccPatchi];
793 
794  if (!isA<nonConformalCoupledFvPatch>(fvp)) continue;
795 
796  const nonConformalCoupledFvPatch& nccFvp =
797  refCast<const nonConformalCoupledFvPatch>(fvp);
798 
799  if (!nccFvp.owner()) continue;
800 
801  const label origPatchi = nccFvp.origPatchIndex();
802  const fvPatch& origFvp = nccFvp.origPatch();
803 
804  forAll(nccFvp, nccPatchFacei)
805  {
806  const label bFacei =
807  nccFvp.polyFaces()[nccPatchFacei] - mesh().nInternalFaces();
808 
809  const label origPatchFacei =
810  nccFvp.polyFaces()[nccPatchFacei] - origFvp.start();
811 
812  const scalar deltaPhi =
813  magSfBf[nccPatchi][nccPatchFacei]
814  /bFaceNccMagSf[bFacei]
815  *deltaPhiBf[origPatchi][origPatchFacei];
816 
817  deltaPhiBf[nccPatchi][nccPatchFacei] = deltaPhi;
818  }
819  }
820  forAll(mesh().boundary(), nccPatchi)
821  {
822  const fvPatch& fvp = mesh().boundary()[nccPatchi];
823 
824  if (!isA<nonConformalCoupledFvPatch>(fvp)) continue;
825 
826  const nonConformalCoupledFvPatch& nccFvp =
827  refCast<const nonConformalCoupledFvPatch>(fvp);
828 
829  if (!nccFvp.owner()) continue;
830 
831  const label origPatchi = nccFvp.origPatchIndex();
832  const fvPatch& origFvp = nccFvp.origPatch();
833 
834  forAll(nccFvp, nccPatchFacei)
835  {
836  const label origPatchFacei =
837  nccFvp.polyFaces()[nccPatchFacei] - origFvp.start();
838 
839  deltaPhiBf[origPatchi][origPatchFacei] = 0;
840  }
841  }
842 
843  // Correct the boundary mesh face fluxes
844  phiBf += deltaPhiBf;
845 }
846 
847 
848 void Foam::fvMeshStitchers::moving::unconformErrorFaceCorrectMeshPhi
849 (
850  const surfaceLabelField::Boundary& polyFacesBf,
851  surfaceVectorField& SfSf,
852  surfaceVectorField& CfSf,
853  surfaceScalarField& phi
854 )
855 {
856  const nonConformalBoundary& ncb = nonConformalBoundary::New(mesh());
857  const labelList origPatchIndices = ncb.allOrigPatchIndices();
858  const labelList errorPatchIndices = ncb.allErrorPatchIndices();
859 
860  // Synchronise the mesh fluxes on both sides of coupled patches. Store
861  // the change made to the mesh flux as an error.
862  PtrList<surfaceScalarField::Boundary> phiErrorbs(phi.nOldTimes(false) + 1);
863  for (label i = 0; i <= phi.nOldTimes(false); ++ i)
864  {
865  tmp<surfaceScalarField::Boundary> tphib =
867  (
868  phi.oldTime(i).boundaryField()
869  );
870 
871  phiErrorbs.set
872  (
873  i,
875  (
877  phi.oldTime(i).boundaryField()
878  )
879  );
880  phiErrorbs[i] = phi.oldTime(i).boundaryField() - tphib();
881 
882  phi.oldTimeRef(i).boundaryFieldRef() = tphib;
883  }
884 
885  // Add the mesh flux error into the error patch so that the mesh fluxes
886  // and volume changes still match
887  forAll(mesh().boundary(), nccPatchi)
888  {
889  const fvPatch& fvp = mesh().boundary()[nccPatchi];
890 
891  if (isA<nonConformalCoupledFvPatch>(fvp))
892  {
893  const nonConformalCoupledFvPatch& nccFvp =
894  refCast<const nonConformalCoupledFvPatch>(fvp);
895 
896  const label origPatchi = nccFvp.origPatchIndex();
897  const polyPatch& origPp = mesh().boundaryMesh()[origPatchi];
898 
899  const label errorPatchi = nccFvp.errorPatchIndex();
900 
901  forAll(nccFvp, nccPatchFacei)
902  {
903  const label origPatchFacei =
904  nccFvp.polyFaces()[nccPatchFacei] - origPp.start();
905 
906  const label errorPatchFacei0 = 2*origPatchFacei;
907  const label errorPatchFacei1 = 2*origPatchFacei + 1;
908 
909  for (label i = 0; i <= phi.nOldTimes(false); ++ i)
910  {
911  fvsPatchField<scalar>& phip =
912  phi.oldTimeRef(i).boundaryFieldRef()[errorPatchi];
913  phip[errorPatchFacei0] +=
914  phiErrorbs[i][nccPatchi][nccPatchFacei]/2;
915  phip[errorPatchFacei1] +=
916  phiErrorbs[i][nccPatchi][nccPatchFacei]/2;
917  }
918  }
919  }
920  }
921 
922  // Create a boundary field with a mesh velocity magnitude. Take the
923  // maximum mesh velocity on either side of non-conformal-coupled faces,
924  // and then average that into the original faces. This means even
925  // stationary original faces have a velocity magnitude stored that is
926  // representative of the interface motion.
927  tmp<surfaceScalarField> tnccMeshMagUf =
929  (
930  "nccMeshMagUf",
931  mesh(),
933  );
934  surfaceScalarField::Boundary& tnccMeshMagUfb =
935  tnccMeshMagUf.ref().boundaryFieldRef();
936  forAll(mesh().boundary(), nccPatchi)
937  {
938  const fvPatch& fvp = mesh().boundary()[nccPatchi];
939 
940  if (isA<nonConformalCoupledFvPatch>(fvp))
941  {
942  const nonConformalCoupledFvPatch& nccFvp =
943  refCast<const nonConformalCoupledFvPatch>(fvp);
944 
945  const fvPatch& origFvp = nccFvp.origPatch();
946 
947  forAll(nccFvp, nccPatchFacei)
948  {
949  const label origPatchFacei =
950  nccFvp.polyFaces()[nccPatchFacei]
951  - origFvp.start();
952 
953  const point& origC =
954  origFvp.patch().faceCentres()[origPatchFacei];
955  const point origC0 =
956  origFvp.patch()[origPatchFacei]
957  .centre(mesh().oldPoints());
958 
959  tnccMeshMagUfb[nccPatchi][nccPatchFacei] =
960  mag(origC - origC0)/mesh().time().deltaTValue();
961  }
962  }
963  }
964  tnccMeshMagUfb =
965  max
966  (
967  tnccMeshMagUfb,
968  tnccMeshMagUfb.boundaryNeighbourField()()
969  );
971  (
974  );
975  tnccMeshMagUf.clear();
976 
977  // Resize the error patch faces so that mesh flux divided by area
978  // results in a velocity equal to the mesh velocity of the original
979  // patch face
980  forAll(errorPatchIndices, i)
981  {
982  const label origPatchi = origPatchIndices[i];
983  const polyPatch& origPp = mesh().boundaryMesh()[origPatchi];
984 
985  const label errorPatchi = errorPatchIndices[i];
986 
987  forAll(origPp, origPatchFacei)
988  {
989  const label errorPatchFacei0 = 2*origPatchFacei;
990  const label errorPatchFacei1 = 2*origPatchFacei + 1;
991 
992  const vector errorSf =
993  min
994  (
995  mag(phi.boundaryField()[errorPatchi][errorPatchFacei0])
996  /max(meshMagUfb[origPatchi][origPatchFacei], vSmall),
997  origPp.magFaceAreas()[origPatchFacei]
998  )
999  *origPp.faceNormals()[origPatchFacei];
1000 
1001  fvsPatchField<vector>& Sfp =
1002  SfSf.boundaryFieldRef()[errorPatchi];
1003  Sfp[errorPatchFacei0] += errorSf;
1004  Sfp[errorPatchFacei1] -= errorSf;
1005  }
1006  }
1007 }
1008 
1009 
1010 void Foam::fvMeshStitchers::moving::unconformCorrectMeshPhi
1011 (
1012  const SurfaceFieldBoundary<label>& polyFacesBf,
1013  surfaceVectorField& SfSf,
1014  surfaceVectorField& CfSf,
1015  surfaceScalarField& phi
1016 )
1017 {
1018  // !!! At present, the correction procedures that follow require the mesh
1019  // to be unconformed to its final topological state. This means we have to
1020  // call fvMesh::unconform twice; once to change the topology to allow for
1021  // calculation of the error areas and mesh flux corrections, and once again
1022  // to actually apply these error areas and mesh flux corrections. This is
1023  // OK for now (the operation is relatively quick), but ideally everything
1024  // would just be calculated in advance and applied to the mesh in one go,
1025  // and this function would only modify its arguments and leave calling
1026  // fvMesh::unconform to the base class.
1027  mesh().unconform(polyFacesBf, SfSf, CfSf);
1028 
1029  // Resize the patched in the flux field
1030  for (label i = 0; i <= phi.nOldTimes(false); ++ i)
1031  {
1033  boundaryFieldRefNoUpdate(phi.oldTime(i));
1034 
1035  forAll(polyFacesBf, ncPatchi)
1036  {
1037  if (!isA<nonConformalFvPatch>(polyFacesBf[ncPatchi].patch()))
1038  {
1039  phi0Bf[ncPatchi].map
1040  (
1041  phi0Bf[ncPatchi],
1042  setSizeFieldMapper(polyFacesBf[ncPatchi].size())
1043  );
1044  }
1045  }
1046  }
1047 
1048  // Set mesh fluxes on the original and cyclic faces as a proportion of
1049  // the area taken from the old original faces
1050  for (label i = 0; i <= phi.nOldTimes(false); ++ i)
1051  {
1052  phi.oldTimeRef(i).boundaryFieldRef() =
1054  (
1055  phi.oldTime(i).boundaryField(),
1056  phi.oldTime(i).boundaryField()
1057  );
1058  }
1059 
1060  // Correct the mesh flux error by modifying values on the internal
1061  // faces that are edge-connected to the owner orig patches.
1062  //
1063  // !!! This only corrects the new-time flux. For smooth operation with
1064  // second order time schemes on wonky meshes, this will probably need to
1065  // correct old time mesh fluxes, too.
1066  //
1067  if
1068  (
1069  mesh().foundObject<solutionControl>(solutionControl::typeName)
1070  && mesh().lookupObject<solutionControl>(solutionControl::typeName)
1071  .dict().lookup<Switch>("correctMeshPhi")
1072  )
1073  {
1074  unconformInternalFaceCorrectMeshPhi(phi);
1075  }
1076 
1077  // Correct the mesh flux error by adding flux and area to the error
1078  // patch faces
1079  unconformErrorFaceCorrectMeshPhi(polyFacesBf, SfSf, CfSf, phi);
1080 }
1081 
1082 
1083 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1084 
1086 :
1087  fvMeshStitcher(mesh)
1088 {}
1089 
1090 
1091 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1092 
1094 {}
1095 
1096 
1097 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
static nonConformalBoundary & New(const word &name, const polyMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
static const DimensionedField< Type, GeoMesh > & null()
Return a null DimensionedField.
static const char *const typeName
Definition: Field.H:106
DimensionedField< Type, GeoMesh > Internal
Type of the internal field from which this GeometricField is derived.
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 >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
const Vector< Cmpt > & centre(const Foam::List< Vector< Cmpt >> &) const
Return *this (used for point which is a typedef to Vector<scalar>.
Definition: VectorI.H:116
Mesh manipulator that uses the intersection provided by the cyclic non-conformal poly patches to crea...
Mesh stitcher for moving meshes.
moving(fvMesh &)
Construct from fvMesh.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
A zero-sized class without any storage. Used, for example, in HashSet.
Definition: nil.H:59
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: polyPatch.C:238
static void syncEdgeList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronise values on all mesh edges.
Calculate the matrix for the laplacian of the field.
label patchi
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){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const fvPatchList & patches
volScalarField & b
Definition: createFields.H:25
bool valid(const PtrList< ModelType > &l)
tmp< SurfaceFieldBoundary< Type > > conformedNcBoundaryField(const SurfaceFieldBoundary< Type > &fieldb)
Extract the non-conformal parts of the boundary field and store it on the.
tmp< Field< Type > > fieldRMapSum(const Field< Type > &f, const label size, const labelUList &addr, const label addr0=0)
Reverse map a field with an (optional) addressing offset, initialising the.
tmp< SurfaceFieldBoundary< Type > > synchronisedBoundaryField(const SurfaceFieldBoundary< Type > &fieldb, const bool flip, const scalar ownerWeight, const scalar neighbourWeight)
Synchronise the boundary field by combining corresponding.
tmp< SurfaceFieldBoundary< Type > > unconformedBoundaryField(const SurfaceFieldBoundary< Type > &ncFieldb, const SurfaceFieldBoundary< Type > &origFieldb)
Combine non-conformal and original parts of the boundary field from the.
addToRunTimeSelectionTable(fvMeshStitcher, stationary, fvMesh)
defineTypeNameAndDebug(stationary, 0)
label wave(const fvMesh &mesh, const List< labelPair > &changedPatchAndFaces, const label nCorrections, GeometricField< scalar, PatchField, GeoMesh > &distance, TrackingData &td, GeometricField< DataType, PatchField, GeoMesh > &... data)
Wave distance (and maybe additional) data from faces. If nCorrections is.
void surfaceIntegrate(Field< Type > &ivf, const SurfaceField< Type > &ssf)
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
List< word > wordList
A List of words.
Definition: fileName.H:54
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
SurfaceField< scalar > surfaceScalarField
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
vector point
Point is a vector.
Definition: point.H:41
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
const dimensionSet dimVolume
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:213
const dimensionSet dimVelocity
const dimensionSet dimArea
static const label labelMax
Definition: label.H:62
SurfaceField< vector > surfaceVectorField
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
List< Type > repeat(const UList< Type > &a, const UList< Type > &b)
UList< label > labelUList
Definition: UList.H:65
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
faceListList boundary(nPatches)
dictionary dict
friend Ostream & operator<<(Ostream &os, const layerAndWeight &l)
static const layerAndWeight max
friend bool operator!=(const layerAndWeight &a, const layerAndWeight &b)
friend Istream & operator>>(Istream &is, layerAndWeight &l)
friend bool operator==(const layerAndWeight &a, const layerAndWeight &b)
static const layerAndWeight min