meshToMesh0.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) 2011-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 "meshToMesh0.H"
27 #include "processorFvPatch.H"
28 #include "demandDrivenData.H"
29 #include "treeDataCell.H"
30 #include "treeDataFace.H"
31 #include "tetOverlapVolume.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 defineTypeNameAndDebug(meshToMesh0, 0);
38 }
39 
40 const Foam::scalar Foam::meshToMesh0::directHitTol = 1e-5;
41 
42 
43 
44 
45 
46 
47 void Foam::meshToMesh0::calcAddressing()
48 {
49  if (debug)
50  {
52  << "Calculating mesh-to-mesh cell addressing" << endl;
53  }
54 
55  // Set reference to cells
56  const cellList& fromCells = fromMesh_.cells();
57  const pointField& fromPoints = fromMesh_.points();
58 
59  // In an attempt to preserve the efficiency of linear search and prevent
60  // failure, a RESCUE mechanism will be set up. Here, we shall mark all
61  // cells next to the solid boundaries. If such a cell is found as the
62  // closest, the relationship between the origin and cell will be examined.
63  // If the origin is outside the cell, a global n-squared search is
64  // triggered.
65 
66  // SETTING UP RESCUE
67 
68  // Visit all boundaries and mark the cell next to the boundary.
69 
70  if (debug)
71  {
72  InfoInFunction << "Setting up rescue" << endl;
73  }
74 
75  List<bool> boundaryCell(fromCells.size(), false);
76 
77  // Set reference to boundary
78  const polyPatchList& patchesFrom = fromMesh_.boundaryMesh();
79 
80  forAll(patchesFrom, patchi)
81  {
82  // Get reference to cells next to the boundary
83  const labelUList& bCells = patchesFrom[patchi].faceCells();
84 
85  forAll(bCells, facei)
86  {
87  boundaryCell[bCells[facei]] = true;
88  }
89  }
90 
91  treeBoundBox meshBb(fromPoints);
92 
93  scalar typDim = meshBb.avgDim()/(2.0*cbrt(scalar(fromCells.size())));
94 
95  treeBoundBox shiftedBb
96  (
97  meshBb.min(),
98  meshBb.max() + vector(typDim, typDim, typDim)
99  );
100 
101  if (debug)
102  {
103  Info<< "\nMesh\n"
104  << " bounding box : " << meshBb << nl
105  << " bounding box (shifted) : " << shiftedBb << nl
106  << " typical dimension :" << shiftedBb.typDim() << endl;
107  }
108 
109  indexedOctree<treeDataCell> oc
110  (
111  treeDataCell(false, fromMesh_, polyMesh::CELL_TETS),
112  shiftedBb, // overall bounding box
113  8, // maxLevel
114  10, // leafsize
115  6.0 // duplicity
116  );
117 
118  if (debug)
119  {
120  oc.print(Pout, false, 0);
121  }
122 
123  cellAddresses
124  (
125  cellAddressing_,
126  toMesh_.cellCentres(),
127  fromMesh_,
128  boundaryCell,
129  oc
130  );
131 
132  forAll(toMesh_.boundaryMesh(), patchi)
133  {
134  const polyPatch& toPatch = toMesh_.boundaryMesh()[patchi];
135 
136  if (cuttingPatches_.found(toPatch.name()))
137  {
138  boundaryAddressing_[patchi].setSize(toPatch.size());
139 
140  cellAddresses
141  (
142  boundaryAddressing_[patchi],
143  toPatch.faceCentres(),
144  fromMesh_,
145  boundaryCell,
146  oc
147  );
148  }
149  else if
150  (
151  patchMap_.found(toPatch.name())
152  && fromMeshPatches_.found(patchMap_.find(toPatch.name())())
153  )
154  {
155  const polyPatch& fromPatch = fromMesh_.boundaryMesh()
156  [
157  fromMeshPatches_.find(patchMap_.find(toPatch.name())())()
158  ];
159 
160  if (fromPatch.empty())
161  {
163  << "Source patch " << fromPatch.name()
164  << " has no faces. Not performing mapping for it."
165  << endl;
166  boundaryAddressing_[patchi].setSize(toPatch.size());
167  boundaryAddressing_[patchi] = -1;
168  }
169  else
170  {
171  treeBoundBox wallBb(fromPatch.localPoints());
172  scalar typDim =
173  wallBb.avgDim()/(2.0*sqrt(scalar(fromPatch.size())));
174 
175  treeBoundBox shiftedBb
176  (
177  wallBb.min(),
178  wallBb.max() + vector(typDim, typDim, typDim)
179  );
180 
181  // Note: allow more levels than in meshSearch. Assume patch
182  // is not as big as all boundary faces
183  indexedOctree<treeDataFace> oc
184  (
185  treeDataFace(false, fromPatch),
186  shiftedBb, // overall search domain
187  12, // maxLevel
188  10, // leafsize
189  6.0 // duplicity
190  );
191 
192  const vectorField::subField centresToBoundary =
193  toPatch.faceCentres();
194 
195  boundaryAddressing_[patchi].setSize(toPatch.size());
196 
197  const scalar distSqr = sqr(wallBb.mag());
198 
199  forAll(toPatch, toi)
200  {
201  boundaryAddressing_[patchi][toi] = oc.findNearest
202  (
203  centresToBoundary[toi],
204  distSqr
205  ).index();
206  }
207  }
208  }
209  }
210 
211  if (debug)
212  {
214  << "Finished calculating mesh-to-mesh cell addressing" << endl;
215  }
216 }
217 
218 
219 void Foam::meshToMesh0::cellAddresses
220 (
221  labelList& cellAddressing_,
222  const pointField& points,
223  const fvMesh& fromMesh,
224  const List<bool>& boundaryCell,
225  const indexedOctree<treeDataCell>& oc
226 ) const
227 {
228  // The implemented search method is a simple neighbour array search.
229  // It starts from a cell zero, searches its neighbours and finds one
230  // which is nearer to the target point than the current position.
231  // The location of the "current position" is reset to that cell and
232  // search through the neighbours continues. The search is finished
233  // when all the neighbours of the cell are farther from the target
234  // point than the current cell
235 
236  // Set curCell label to zero (start)
237  label curCell = 0;
238 
239  // Set reference to cell to cell addressing
240  const vectorField& centresFrom = fromMesh.cellCentres();
241  const labelListList& cc = fromMesh.cellCells();
242 
243  forAll(points, toi)
244  {
245  // Pick up target position
246  const vector& p = points[toi];
247 
248  // Set the sqr-distance
249  scalar distSqr = magSqr(p - centresFrom[curCell]);
250 
251  bool closer;
252 
253  do
254  {
255  closer = false;
256 
257  // Set the current list of neighbouring cells
258  const labelList& neighbours = cc[curCell];
259 
260  forAll(neighbours, ni)
261  {
262  const scalar curDistSqr =
263  magSqr(p - centresFrom[neighbours[ni]]);
264 
265  // Search through all the neighbours.
266  // If the cell is closer, reset current cell and distance
267  if (curDistSqr < (1 - small)*distSqr)
268  {
269  curCell = neighbours[ni];
270  distSqr = curDistSqr;
271  closer = true; // A closer neighbour has been found
272  }
273  }
274  } while (closer);
275 
276  cellAddressing_[toi] = -1;
277 
278  // Check point is actually in the nearest cell
279  if (fromMesh.pointInCell(p, curCell))
280  {
281  cellAddressing_[toi] = curCell;
282  }
283  else
284  {
285  // If curCell is a boundary cell then the point maybe either outside
286  // the domain or in an other region of the domain, either way use
287  // the octree search to find it.
288  if (boundaryCell[curCell])
289  {
290  cellAddressing_[toi] = oc.findInside(p);
291  if (cellAddressing_[toi] != -1)
292  {
293  curCell = cellAddressing_[toi];
294  }
295  }
296  else
297  {
298  // If not on the boundary search the neighbours
299  bool found = false;
300 
301  // Set the current list of neighbouring cells
302  const labelList& neighbours = cc[curCell];
303 
304  forAll(neighbours, ni)
305  {
306  // Search through all the neighbours.
307  // If point is in neighbour reset current cell
308  if (fromMesh.pointInCell(p, neighbours[ni]))
309  {
310  cellAddressing_[toi] = neighbours[ni];
311  found = true;
312  break;
313  }
314  }
315 
316  if (!found)
317  {
318  // If still not found search the neighbour-neighbours
319 
320  // Set the current list of neighbouring cells
321  const labelList& neighbours = cc[curCell];
322 
323  forAll(neighbours, ni)
324  {
325  // Set the current list of neighbour-neighbouring cells
326  const labelList& nn = cc[neighbours[ni]];
327 
328  forAll(nn, ni)
329  {
330  // Search through all the neighbours.
331  // If point is in neighbour reset current cell
332  if (fromMesh.pointInCell(p, nn[ni]))
333  {
334  cellAddressing_[toi] = nn[ni];
335  found = true;
336  break;
337  }
338  }
339  if (found) break;
340  }
341  }
342 
343  if (!found)
344  {
345  // Still not found so use the octree
346  cellAddressing_[toi] = oc.findInside(p);
347 
348  if (cellAddressing_[toi] != -1)
349  {
350  curCell = cellAddressing_[toi];
351  }
352  }
353  }
354  }
355  }
356 }
357 
358 
359 void Foam::meshToMesh0::calculateInverseDistanceWeights() const
360 {
361  if (debug)
362  {
364  << "Calculating inverse distance weighting factors" << endl;
365  }
366 
367  if (inverseDistanceWeightsPtr_)
368  {
370  << "weighting factors already calculated"
371  << exit(FatalError);
372  }
373 
374  //- Initialise overlap volume to zero
375  V_ = 0.0;
376 
377  inverseDistanceWeightsPtr_ = new scalarListList(toMesh_.nCells());
378  scalarListList& invDistCoeffs = *inverseDistanceWeightsPtr_;
379 
380  // get reference to source mesh data
381  const labelListList& cc = fromMesh_.cellCells();
382  const vectorField& centreFrom = fromMesh_.C();
383  const vectorField& centreTo = toMesh_.C();
384 
385  forAll(cellAddressing_, celli)
386  {
387  if (cellAddressing_[celli] != -1)
388  {
389  const vector& target = centreTo[celli];
390  scalar m = mag(target - centreFrom[cellAddressing_[celli]]);
391 
392  const labelList& neighbours = cc[cellAddressing_[celli]];
393 
394  // if the nearest cell is a boundary cell or there is a direct hit,
395  // pick up the value
396  label directCelli = -1;
397  if (m < directHitTol || neighbours.empty())
398  {
399  directCelli = celli;
400  }
401  else
402  {
403  forAll(neighbours, ni)
404  {
405  scalar nm = mag(target - centreFrom[neighbours[ni]]);
406  if (nm < directHitTol)
407  {
408  directCelli = neighbours[ni];
409  break;
410  }
411  }
412  }
413 
414 
415  if (directCelli != -1)
416  {
417  // Direct hit
418  invDistCoeffs[directCelli].setSize(1);
419  invDistCoeffs[directCelli][0] = 1.0;
420  V_ += fromMesh_.V()[cellAddressing_[directCelli]];
421  }
422  else
423  {
424  invDistCoeffs[celli].setSize(neighbours.size() + 1);
425 
426  // The first coefficient corresponds to the centre cell.
427  // The rest is ordered in the same way as the cellCells list.
428  scalar invDist = 1.0/m;
429  invDistCoeffs[celli][0] = invDist;
430  scalar sumInvDist = invDist;
431 
432  // now add the neighbours
433  forAll(neighbours, ni)
434  {
435  invDist = 1.0/mag(target - centreFrom[neighbours[ni]]);
436  invDistCoeffs[celli][ni + 1] = invDist;
437  sumInvDist += invDist;
438  }
439 
440  // divide by the total inverse-distance
441  forAll(invDistCoeffs[celli], i)
442  {
443  invDistCoeffs[celli][i] /= sumInvDist;
444  }
445 
446 
447  V_ +=
448  invDistCoeffs[celli][0]
449  *fromMesh_.V()[cellAddressing_[celli]];
450  for (label i = 1; i < invDistCoeffs[celli].size(); i++)
451  {
452  V_ +=
453  invDistCoeffs[celli][i]*fromMesh_.V()[neighbours[i-1]];
454  }
455  }
456  }
457  }
458 }
459 
460 
461 void Foam::meshToMesh0::calculateInverseVolumeWeights() const
462 {
463  if (debug)
464  {
466  << "Calculating inverse volume weighting factors" << endl;
467  }
468 
469  if (inverseVolumeWeightsPtr_)
470  {
472  << "weighting factors already calculated"
473  << exit(FatalError);
474  }
475 
476  //- Initialise overlap volume to zero
477  V_ = 0.0;
478 
479  inverseVolumeWeightsPtr_ = new scalarListList(toMesh_.nCells());
480  scalarListList& invVolCoeffs = *inverseVolumeWeightsPtr_;
481 
482  const labelListList& cellToCell = cellToCellAddressing();
483 
484  tetOverlapVolume overlapEngine;
485 
486  forAll(cellToCell, celli)
487  {
488  const labelList& overlapCells = cellToCell[celli];
489 
490  if (overlapCells.size() > 0)
491  {
492  invVolCoeffs[celli].setSize(overlapCells.size());
493 
494  forAll(overlapCells, j)
495  {
496  label cellFrom = overlapCells[j];
497  treeBoundBox bbFromMesh
498  (
499  pointField
500  (
501  fromMesh_.points(),
502  fromMesh_.cellPoints()[cellFrom]
503  )
504  );
505 
506  scalar v = overlapEngine.cellCellOverlapVolumeMinDecomp
507  (
508  toMesh_,
509  celli,
510 
511  fromMesh_,
512  cellFrom,
513  bbFromMesh
514  );
515  invVolCoeffs[celli][j] = v/toMesh_.V()[celli];
516 
517  V_ += v;
518  }
519  }
520  }
521 }
522 
523 
524 void Foam::meshToMesh0::calculateCellToCellAddressing() const
525 {
526  if (debug)
527  {
529  << "Calculating cell to cell addressing" << endl;
530  }
531 
532  if (cellToCellAddressingPtr_)
533  {
535  << "addressing already calculated"
536  << exit(FatalError);
537  }
538 
539  //- Initialise overlap volume to zero
540  V_ = 0.0;
541 
542  tetOverlapVolume overlapEngine;
543 
544  cellToCellAddressingPtr_ = new labelListList(toMesh_.nCells());
545  labelListList& cellToCell = *cellToCellAddressingPtr_;
546 
547 
548  forAll(cellToCell, iTo)
549  {
550  const labelList overLapCells =
551  overlapEngine.overlappingCells(fromMesh_, toMesh_, iTo);
552  if (overLapCells.size() > 0)
553  {
554  cellToCell[iTo].setSize(overLapCells.size());
555  forAll(overLapCells, j)
556  {
557  cellToCell[iTo][j] = overLapCells[j];
558  V_ += fromMesh_.V()[overLapCells[j]];
559  }
560  }
561  }
562 }
563 
564 
565 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
566 
568 (
569  const fvMesh& meshFrom,
570  const fvMesh& meshTo,
571  const HashTable<word>& patchMap,
572  const wordList& cuttingPatchNames
573 )
574 :
575  fromMesh_(meshFrom),
576  toMesh_(meshTo),
577  patchMap_(patchMap),
578  cellAddressing_(toMesh_.nCells()),
579  boundaryAddressing_(toMesh_.boundaryMesh().size()),
580  inverseDistanceWeightsPtr_(nullptr),
581  inverseVolumeWeightsPtr_(nullptr),
582  cellToCellAddressingPtr_(nullptr),
583  V_(0.0)
584 {
585  forAll(fromMesh_.boundaryMesh(), patchi)
586  {
587  fromMeshPatches_.insert
588  (
589  fromMesh_.boundaryMesh()[patchi].name(),
590  patchi
591  );
592  }
593 
594  forAll(toMesh_.boundaryMesh(), patchi)
595  {
596  toMeshPatches_.insert
597  (
598  toMesh_.boundaryMesh()[patchi].name(),
599  patchi
600  );
601  }
602 
603  forAll(cuttingPatchNames, i)
604  {
605  if (toMeshPatches_.found(cuttingPatchNames[i]))
606  {
607  cuttingPatches_.insert
608  (
609  cuttingPatchNames[i],
610  toMeshPatches_.find(cuttingPatchNames[i])()
611  );
612  }
613  else
614  {
616  << "Cannot find cutting-patch " << cuttingPatchNames[i]
617  << " in destination mesh" << endl;
618  }
619  }
620 
621  forAll(toMesh_.boundaryMesh(), patchi)
622  {
623  // Add the processor patches in the toMesh to the cuttingPatches list
624  if (isA<processorPolyPatch>(toMesh_.boundaryMesh()[patchi]))
625  {
626  cuttingPatches_.insert
627  (
628  toMesh_.boundaryMesh()[patchi].name(),
629  patchi
630  );
631  }
632  }
633 
634  calcAddressing();
635 }
636 
637 
639 (
640  const fvMesh& meshFrom,
641  const fvMesh& meshTo
642 )
643 :
644  fromMesh_(meshFrom),
645  toMesh_(meshTo),
646  cellAddressing_(toMesh_.nCells()),
647  boundaryAddressing_(toMesh_.boundaryMesh().size()),
648  inverseDistanceWeightsPtr_(nullptr),
649  inverseVolumeWeightsPtr_(nullptr),
650  cellToCellAddressingPtr_(nullptr),
651  V_(0.0)
652 {
653  // check whether both meshes have got the same number
654  // of boundary patches
655  if (fromMesh_.boundary().size() != toMesh_.boundary().size())
656  {
658  << "Incompatible meshes: different number of patches, "
659  << "fromMesh = " << fromMesh_.boundary().size()
660  << ", toMesh = " << toMesh_.boundary().size()
661  << exit(FatalError);
662  }
663 
664  forAll(fromMesh_.boundaryMesh(), patchi)
665  {
666  if
667  (
668  fromMesh_.boundaryMesh()[patchi].name()
669  != toMesh_.boundaryMesh()[patchi].name()
670  )
671  {
673  << "Incompatible meshes: different patch names for patch "
674  << patchi
675  << ", fromMesh = " << fromMesh_.boundary()[patchi].name()
676  << ", toMesh = " << toMesh_.boundary()[patchi].name()
677  << exit(FatalError);
678  }
679 
680  if
681  (
682  fromMesh_.boundaryMesh()[patchi].type()
683  != toMesh_.boundaryMesh()[patchi].type()
684  )
685  {
687  << "Incompatible meshes: different patch types for patch "
688  << patchi
689  << ", fromMesh = " << fromMesh_.boundary()[patchi].type()
690  << ", toMesh = " << toMesh_.boundary()[patchi].type()
691  << exit(FatalError);
692  }
693 
694  fromMeshPatches_.insert
695  (
696  fromMesh_.boundaryMesh()[patchi].name(),
697  patchi
698  );
699 
700  toMeshPatches_.insert
701  (
702  toMesh_.boundaryMesh()[patchi].name(),
703  patchi
704  );
705 
706  patchMap_.insert
707  (
708  toMesh_.boundaryMesh()[patchi].name(),
709  fromMesh_.boundaryMesh()[patchi].name()
710  );
711  }
712 
713  calcAddressing();
714 }
715 
716 
717 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
718 
720 {
721  deleteDemandDrivenData(inverseDistanceWeightsPtr_);
722  deleteDemandDrivenData(inverseVolumeWeightsPtr_);
723  deleteDemandDrivenData(cellToCellAddressingPtr_);
724 }
725 
726 
727 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
728 
729 const Foam::scalarListList& Foam::meshToMesh0::inverseDistanceWeights() const
730 {
731  if (!inverseDistanceWeightsPtr_)
732  {
733  calculateInverseDistanceWeights();
734  }
735 
736  return *inverseDistanceWeightsPtr_;
737 }
738 
739 
740 const Foam::scalarListList& Foam::meshToMesh0::inverseVolumeWeights() const
741 {
742  if (!inverseVolumeWeightsPtr_)
743  {
744  calculateInverseVolumeWeights();
745  }
746 
747  return *inverseVolumeWeightsPtr_;
748 }
749 
750 
751 const Foam::labelListList& Foam::meshToMesh0::cellToCellAddressing() const
752 {
753  if (!cellToCellAddressingPtr_)
754  {
755  calculateCellToCellAddressing();
756  }
757 
758  return *cellToCellAddressingPtr_;
759 }
760 
761 
762 // ************************************************************************* //
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
SubField< vector > subField
Declare type of subField.
Definition: Field.H:100
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
void setSize(const label)
Reset size of List.
Definition: List.C:281
meshToMesh0(const fvMesh &fromMesh, const fvMesh &toMesh, const HashTable< word > &patchMap, const wordList &cuttingPatchNames)
Construct from the two meshes, the patch name map for the patches.
~meshToMesh0()
Destructor.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:405
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1360
const vectorField & cellCentres() const
const cellList & cells() const
Template functions to aid in the implementation of demand driven data.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
label patchi
const pointField & points
#define WarningInFunction
Report a warning using Foam::Warning.
#define InfoInFunction
Report an information message using Foam::Info.
const dimensionedScalar e
Elementary charge.
Namespace for OpenFOAM.
List< scalarList > scalarListList
Definition: scalarList.H:51
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
List< word > wordList
A List of words.
Definition: fileName.H:54
List< label > labelList
A List of labels.
Definition: labelList.H:56
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)
List< cell > cellList
list of cells
Definition: cellList.H:42
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
void deleteDemandDrivenData(DataPtr &dataPtr)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
dimensionedScalar sqrt(const dimensionedScalar &ds)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
PtrList< polyPatch > polyPatchList
container classes for polyPatch
Definition: polyPatchList.H:45
Field< vector > vectorField
Specialisation of Field<T> for vector.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
error FatalError
dimensionedScalar cbrt(const dimensionedScalar &ds)
UList< label > labelUList
Definition: UList.H:65
static const char nl
Definition: Ostream.H:260
dimensioned< scalar > magSqr(const dimensioned< Type > &)
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(1e-3))
volScalarField & p