fvMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017,2022 OpenFOAM Foundation
9  Copyright (C) 2016-2024 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "fvMesh.H"
30 #include "volFields.H"
31 #include "surfaceFields.H"
32 #include "slicedVolFields.H"
33 #include "slicedSurfaceFields.H"
34 #include "SubField.H"
35 #include "demandDrivenData.H"
36 #include "fvMeshLduAddressing.H"
37 #include "mapPolyMesh.H"
38 #include "MapFvFields.H"
39 #include "fvMeshMapper.H"
40 #include "mapClouds.H"
41 #include "MeshObject.H"
42 #include "fvMatrix.H"
43 
44 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
49 }
50 
51 
52 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53 
55 {
57  <
58  fvMesh,
61  >(*this);
62 
64  <
65  lduMesh,
68  >(*this);
69 
75 }
76 
77 
79 {
80  const bool haveV = (VPtr_ != nullptr);
81  const bool haveSf = (SfPtr_ != nullptr);
82  const bool haveMagSf = (magSfPtr_ != nullptr);
83  const bool haveCP = (CPtr_ != nullptr);
84  const bool haveCf = (CfPtr_ != nullptr);
85 
86  clearGeomNotOldVol();
87 
88  // Now recreate the fields
89  if (haveV)
90  {
91  (void)V();
92  }
93 
94  if (haveSf)
95  {
96  (void)Sf();
97  }
98 
99  if (haveMagSf)
100  {
101  (void)magSf();
102  }
103 
104  if (haveCP)
105  {
106  (void)C();
107  }
108 
109  if (haveCf)
110  {
111  (void)Cf();
112  }
113 }
114 
115 
117 {
118  clearGeomNotOldVol();
119 
120  deleteDemandDrivenData(V0Ptr_);
122 
123  // Mesh motion flux cannot be deleted here because the old-time flux
124  // needs to be saved.
125 }
126 
127 
128 void Foam::fvMesh::clearAddressing(const bool isMeshUpdate)
129 {
130  DebugInFunction << "isMeshUpdate: " << isMeshUpdate << endl;
131 
132  if (isMeshUpdate)
133  {
134  // Part of a mesh update. Keep meshObjects that have an updateMesh
135  // callback
137  <
138  fvMesh,
141  >
142  (
143  *this
144  );
146  <
147  lduMesh,
150  >
151  (
152  *this
153  );
154  }
155  else
156  {
157  meshObject::clear<fvMesh, TopologicalMeshObject>(*this);
158  meshObject::clear<lduMesh, TopologicalMeshObject>(*this);
159  }
160  deleteDemandDrivenData(lduPtr_);
161 }
162 
163 
165 {
166  if (curTimeIndex_ < time().timeIndex())
167  {
169  << " Storing old time volumes since from time " << curTimeIndex_
170  << " and time now " << time().timeIndex()
171  << " V:" << V.size() << endl;
172 
173  if (V00Ptr_ && V0Ptr_)
174  {
175  // Copy V0 into V00 storage
176  *V00Ptr_ = *V0Ptr_;
177  }
178 
179  if (!V0Ptr_)
180  {
181  V0Ptr_ = new DimensionedField<scalar, volMesh>
182  (
183  IOobject
184  (
185  "V0",
186  this->time().timeName(),
187  *this,
191  ),
192  *this,
193  dimVolume
194  );
195  // Note: V0 now sized with current mesh, not with (potentially
196  // different size) V.
197  V0Ptr_->scalarField::resize_nocopy(V.size());
198  }
199 
200  // Copy V into V0 storage
201  V0Ptr_->scalarField::operator=(V);
202 
203  curTimeIndex_ = time().timeIndex();
204 
205  if (debug)
206  {
208  << " Stored old time volumes V0:" << V0Ptr_->size()
209  << endl;
210 
211  if (V00Ptr_)
212  {
214  << " Stored oldold time volumes V00:" << V00Ptr_->size()
215  << endl;
216  }
217  }
218  }
219 }
220 
221 
222 void Foam::fvMesh::clearOutLocal(const bool isMeshUpdate)
223 {
224  clearGeom();
226 
227  clearAddressing(isMeshUpdate);
228 
229  // Clear mesh motion flux
230  phiPtr_.reset(nullptr);
231 }
232 
233 
234 void Foam::fvMesh::clearOut(const bool isMeshUpdate)
235 {
236  clearOutLocal(isMeshUpdate);
237 
238  polyMesh::clearOut(isMeshUpdate);
239 }
240 
241 
243 {
244  // Clear mesh motion flux
245  phiPtr_.reset(nullptr);
246 }
247 
248 
249 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
250 
251 Foam::fvMesh::fvMesh(const IOobject& io, const bool doInit)
252 :
253  polyMesh(io, doInit),
254  fvSchemes(static_cast<const objectRegistry&>(*this)),
255  surfaceInterpolation(*this),
256  fvSolution(static_cast<const objectRegistry&>(*this)),
257  boundary_(*this, boundaryMesh()),
258  lduPtr_(nullptr),
259  curTimeIndex_(time().timeIndex()),
260  VPtr_(nullptr),
261  V0Ptr_(nullptr),
262  V00Ptr_(nullptr),
263  SfPtr_(nullptr),
264  magSfPtr_(nullptr),
265  CPtr_(nullptr),
266  CfPtr_(nullptr),
267  phiPtr_(nullptr)
268 {
269  DebugInFunction << "Constructing fvMesh from IOobject" << endl;
270 
271  if (doInit)
272  {
273  fvMesh::init(false); // do not initialise lower levels
274  }
275 }
276 
277 
278 bool Foam::fvMesh::init(const bool doInit)
279 {
280  if (doInit)
281  {
282  // Construct basic geometry calculation engine. Note: do before
283  // doing anything with primitiveMesh::cellCentres etc.
284  (void)geometry();
285 
286  // Initialise my data
287  polyMesh::init(doInit);
288  }
289 
290  // Read some optional fields
291  // ~~~~~~~~~~~~~~~~~~~~~~~~~
292  // For redistributePar + fileHandler can have master processor
293  // find the file but not the subprocs.
294 
295  IOobject rio
296  (
297  "none",
298  this->time().timeName(),
299  *this,
303  );
304 
305 
306  // Read old cell volumes (if present) and set the storage of V00
307  rio.resetHeader("V0");
308  if (returnReduceOr(rio.typeHeaderOk<regIOobject>(false)))
309  {
311  << "Detected V0: " << rio.objectRelPath() << nl;
312 
313  V0Ptr_ = new DimensionedField<scalar, volMesh>
314  (
315  rio,
316  *this,
318  );
319 
320  // Set the moving flag early in case demand-driven geometry
321  // construction checks for it
322  moving(true);
323 
324  (void)V00();
325  }
326 
327 
328  // Read mesh fluxes (if present) and set the mesh to be moving
329  rio.resetHeader("meshPhi");
330  if (returnReduceOr(rio.typeHeaderOk<regIOobject>(false)))
331  {
333  << "Detected meshPhi: " << rio.objectRelPath() << nl;
334 
335  // Clear mesh motion flux
336  phiPtr_.reset(nullptr);
337 
338  phiPtr_ = std::make_unique<surfaceScalarField>
339  (
340  rio,
341  *this,
343  );
344 
345  // Set the moving flag early in case demand-driven geometry
346  // construction checks for it
347  moving(true);
348 
349  // The mesh is now considered moving so the old-time cell volumes
350  // will be required for the time derivatives so if they haven't been
351  // read initialise to the current cell volumes
352  if (!V0Ptr_)
353  {
354  V0Ptr_ = new DimensionedField<scalar, volMesh>
355  (
356  IOobject
357  (
358  "V0",
359  this->time().timeName(),
360  *this,
364  ),
365  V()
366  );
367  }
368  }
370  // Assume something changed
371  return true;
372 }
373 
374 
376 (
377  const IOobject& io,
378  pointField&& points,
379  faceList&& faces,
380  labelList&& allOwner,
381  labelList&& allNeighbour,
382  const bool syncPar
383 )
384 :
385  polyMesh
386  (
387  io,
388  std::move(points),
389  std::move(faces),
390  std::move(allOwner),
391  std::move(allNeighbour),
392  syncPar
393  ),
394  fvSchemes(static_cast<const objectRegistry&>(*this)),
395  surfaceInterpolation(*this),
396  fvSolution(static_cast<const objectRegistry&>(*this)),
397  boundary_(*this),
398  lduPtr_(nullptr),
399  curTimeIndex_(time().timeIndex()),
400  VPtr_(nullptr),
401  V0Ptr_(nullptr),
402  V00Ptr_(nullptr),
403  SfPtr_(nullptr),
404  magSfPtr_(nullptr),
405  CPtr_(nullptr),
406  CfPtr_(nullptr),
407  phiPtr_(nullptr)
408 {
409  DebugInFunction << "Constructing fvMesh from components" << endl;
410 }
411 
412 
414 (
415  const IOobject& io,
416  pointField&& points,
417  faceList&& faces,
418  cellList&& cells,
419  const bool syncPar
420 )
421 :
422  polyMesh
423  (
424  io,
425  std::move(points),
426  std::move(faces),
427  std::move(cells),
428  syncPar
429  ),
430  fvSchemes(static_cast<const objectRegistry&>(*this)),
431  surfaceInterpolation(*this),
432  fvSolution(static_cast<const objectRegistry&>(*this)),
433  boundary_(*this),
434  lduPtr_(nullptr),
435  curTimeIndex_(time().timeIndex()),
436  VPtr_(nullptr),
437  V0Ptr_(nullptr),
438  V00Ptr_(nullptr),
439  SfPtr_(nullptr),
440  magSfPtr_(nullptr),
441  CPtr_(nullptr),
442  CfPtr_(nullptr),
443  phiPtr_(nullptr)
444 {
445  DebugInFunction << "Constructing fvMesh from components" << endl;
446 }
447 
448 
449 Foam::fvMesh::fvMesh(const IOobject& io, const Foam::zero, const bool syncPar)
450 :
451  fvMesh(io, pointField(), faceList(), labelList(), labelList(), syncPar)
452 {}
453 
454 
456 (
457  const IOobject& io,
458  const fvMesh& baseMesh,
459  const Foam::zero,
460  const bool syncPar
461 )
462 :
463  fvMesh
464  (
465  io,
466  baseMesh,
467  pointField(),
468  faceList(),
469  labelList(), // owner
470  labelList(), // neighbour
471  syncPar
472  )
473 {}
474 
475 
477 (
478  const IOobject& io,
479  const fvMesh& baseMesh,
480  pointField&& points,
481  faceList&& faces,
482  labelList&& allOwner,
483  labelList&& allNeighbour,
484  const bool syncPar
485 )
486 :
487  polyMesh
488  (
489  io,
490  std::move(points),
491  std::move(faces),
492  std::move(allOwner),
493  std::move(allNeighbour),
494  syncPar
495  ),
496  fvSchemes
497  (
498  static_cast<const objectRegistry&>(*this),
499  io.readOpt(),
500  static_cast<const dictionary*>(baseMesh.hasSchemes())
501  ),
502  surfaceInterpolation(*this),
503  fvSolution
504  (
505  static_cast<const objectRegistry&>(*this),
506  io.readOpt(),
507  static_cast<const dictionary*>(baseMesh.hasSolution())
508  ),
509  boundary_(*this),
510  lduPtr_(nullptr),
511  curTimeIndex_(time().timeIndex()),
512  VPtr_(nullptr),
513  V0Ptr_(nullptr),
514  V00Ptr_(nullptr),
515  SfPtr_(nullptr),
516  magSfPtr_(nullptr),
517  CPtr_(nullptr),
518  CfPtr_(nullptr),
519  phiPtr_(nullptr)
520 {
521  DebugInFunction << "Constructing fvMesh as copy and primitives" << endl;
523  // Reset mesh data
524  data().reset(baseMesh.data());
525 }
526 
527 
529 (
530  const IOobject& io,
531  const fvMesh& baseMesh,
532  pointField&& points,
533  faceList&& faces,
534  cellList&& cells,
535  const bool syncPar
536 )
537 :
538  polyMesh
539  (
540  io,
541  std::move(points),
542  std::move(faces),
543  std::move(cells),
544  syncPar
545  ),
546  fvSchemes
547  (
548  static_cast<const objectRegistry&>(*this),
549  io.readOpt(),
550  static_cast<const dictionary*>(baseMesh.hasSchemes())
551  ),
552  surfaceInterpolation(*this),
553  fvSolution
554  (
555  static_cast<const objectRegistry&>(*this),
556  io.readOpt(),
557  static_cast<const dictionary*>(baseMesh.hasSolution())
558  ),
559  boundary_(*this),
560  lduPtr_(nullptr),
561  curTimeIndex_(time().timeIndex()),
562  VPtr_(nullptr),
563  V0Ptr_(nullptr),
564  V00Ptr_(nullptr),
565  SfPtr_(nullptr),
566  magSfPtr_(nullptr),
567  CPtr_(nullptr),
568  CfPtr_(nullptr),
569  phiPtr_(nullptr)
570 {
571  DebugInFunction << "Constructing fvMesh as copy and primitives" << endl;
572 
573  // Reset mesh data
574  data().reset(baseMesh.data());
575 }
576 
577 
578 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
579 
581 {
582  clearOut();
583 }
584 
585 
586 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
589 {
590  return static_cast<const fvSchemes*>(this);
591 }
592 
595 {
596  return static_cast<const fvSolution*>(this);
597 }
598 
601 {
602  return static_cast<const fvSchemes&>(*this);
603 }
604 
607 {
608  return static_cast<fvSchemes&>(*this);
609 }
610 
613 {
614  return static_cast<const fvSolution&>(*this);
615 }
616 
617 
619 {
620  return static_cast<fvSolution&>(*this);
621 }
622 
623 
625 (
626  fvMatrix<scalar>& m,
627  const dictionary& dict
628 ) const
629 {
630  // Redirect to fvMatrix solver
631  return m.solveSegregatedOrCoupled(dict);
632 }
633 
634 
636 (
637  fvMatrix<vector>& m,
638  const dictionary& dict
639 ) const
640 {
641  // Redirect to fvMatrix solver
642  return m.solveSegregatedOrCoupled(dict);
643 }
644 
645 
647 (
649  const dictionary& dict
650 ) const
651 {
652  // Redirect to fvMatrix solver
653  return m.solveSegregatedOrCoupled(dict);
654 }
655 
656 
658 (
660  const dictionary& dict
661 ) const
662 {
663  // Redirect to fvMatrix solver
664  return m.solveSegregatedOrCoupled(dict);
665 }
666 
667 
669 (
670  fvMatrix<tensor>& m,
671  const dictionary& dict
672 ) const
673 {
674  // Redirect to fvMatrix solver
675  return m.solveSegregatedOrCoupled(dict);
676 }
677 
678 
680 (
681  polyPatchList& plist,
682  const bool validBoundary
683 )
684 {
685  if (boundary().size())
686  {
688  << " boundary already exists"
689  << abort(FatalError);
690  }
692  addPatches(plist, validBoundary);
693  boundary_.addPatches(boundaryMesh());
694 }
695 
696 
698 (
699  const List<polyPatch*>& p,
700  const bool validBoundary
701 )
702 {
703  // Acquire ownership of the pointers
704  polyPatchList plist(const_cast<List<polyPatch*>&>(p));
705 
706  addFvPatches(plist, validBoundary);
707 }
708 
709 
711 {
712  DebugInFunction << "Removing boundary patches." << endl;
713 
714  // Remove fvBoundaryMesh data first.
715  boundary_.clear();
717 
718  clearOut();
719 }
720 
721 
723 {
724  DebugInFunction << "Updating fvMesh";
725 
727 
728  if (state == polyMesh::TOPO_PATCH_CHANGE)
729  {
730  DebugInfo << "Boundary and topological update" << endl;
731 
732  boundary_.readUpdate(boundaryMesh());
733 
734  clearOut();
735 
736  }
737  else if (state == polyMesh::TOPO_CHANGE)
738  {
739  DebugInfo << "Topological update" << endl;
740 
741  // fvMesh::clearOut() but without the polyMesh::clearOut
742  clearOutLocal();
743  }
744  else if (state == polyMesh::POINTS_MOVED)
745  {
746  DebugInfo << "Point motion update" << endl;
747 
748  clearGeom();
749  }
750  else
751  {
752  DebugInfo << "No update" << endl;
753  }
754 
755  return state;
756 }
757 
758 
760 {
761  if (!lduPtr_)
762  {
764  << "Calculating fvMeshLduAddressing from nFaces:"
765  << nFaces() << endl;
766 
767  lduPtr_ = new fvMeshLduAddressing(*this);
768 
769  return *lduPtr_;
770  }
771 
772  return *lduPtr_;
773 }
774 
777 {
778  return boundary().interfaces();
779 }
780 
781 
782 void Foam::fvMesh::mapFields(const mapPolyMesh& meshMap)
783 {
785  << " nOldCells:" << meshMap.nOldCells()
786  << " nCells:" << nCells()
787  << " nOldFaces:" << meshMap.nOldFaces()
788  << " nFaces:" << nFaces()
789  << endl;
790 
791  // We require geometric properties valid for the old mesh
792  if
793  (
794  meshMap.cellMap().size() != nCells()
795  || meshMap.faceMap().size() != nFaces()
796  )
797  {
799  << "mapPolyMesh does not correspond to the old mesh."
800  << " nCells:" << nCells()
801  << " cellMap:" << meshMap.cellMap().size()
802  << " nOldCells:" << meshMap.nOldCells()
803  << " nFaces:" << nFaces()
804  << " faceMap:" << meshMap.faceMap().size()
805  << " nOldFaces:" << meshMap.nOldFaces()
806  << exit(FatalError);
807  }
808 
809  // Create a mapper
810  const fvMeshMapper mapper(*this, meshMap);
811 
812  // Map all the volFields in the objectRegistry
813  MapGeometricFields<scalar, fvPatchField, fvMeshMapper, volMesh>
814  (mapper);
815  MapGeometricFields<vector, fvPatchField, fvMeshMapper, volMesh>
816  (mapper);
817  MapGeometricFields<sphericalTensor, fvPatchField, fvMeshMapper, volMesh>
818  (mapper);
819  MapGeometricFields<symmTensor, fvPatchField, fvMeshMapper, volMesh>
820  (mapper);
821  MapGeometricFields<tensor, fvPatchField, fvMeshMapper, volMesh>
822  (mapper);
823 
824  // Map all the surfaceFields in the objectRegistry
825  MapGeometricFields<scalar, fvsPatchField, fvMeshMapper, surfaceMesh>
826  (mapper);
827  MapGeometricFields<vector, fvsPatchField, fvMeshMapper, surfaceMesh>
828  (mapper);
830  <
831  sphericalTensor, fvsPatchField, fvMeshMapper, surfaceMesh
832  >
833  (mapper);
834  MapGeometricFields<symmTensor, fvsPatchField, fvMeshMapper, surfaceMesh>
835  (mapper);
836  MapGeometricFields<tensor, fvsPatchField, fvMeshMapper, surfaceMesh>
837  (mapper);
838 
839  // Map all the dimensionedFields in the objectRegistry
840  MapDimensionedFields<scalar, fvMeshMapper, volMesh>(mapper);
841  MapDimensionedFields<vector, fvMeshMapper, volMesh>(mapper);
842  MapDimensionedFields<sphericalTensor, fvMeshMapper, volMesh>(mapper);
843  MapDimensionedFields<symmTensor, fvMeshMapper, volMesh>(mapper);
844  MapDimensionedFields<tensor, fvMeshMapper, volMesh>(mapper);
845 
846  // Map all the clouds in the objectRegistry
847  mapClouds(*this, meshMap);
848 
849 
850  const labelList& cellMap = meshMap.cellMap();
851 
852  // Map the old volume. Just map to new cell labels.
853  if (V0Ptr_)
854  {
855  scalarField& V0 = *V0Ptr_;
856 
857  scalarField savedV0(V0);
858  V0.resize_nocopy(nCells());
859 
860  forAll(V0, i)
861  {
862  if (cellMap[i] > -1)
863  {
864  V0[i] = savedV0[cellMap[i]];
865  }
866  else
867  {
868  V0[i] = 0.0;
869  }
870  }
871 
872  // Inject volume of merged cells
873  label nMerged = 0;
874  forAll(meshMap.reverseCellMap(), oldCelli)
875  {
876  label index = meshMap.reverseCellMap()[oldCelli];
877 
878  if (index < -1)
879  {
880  label celli = -index-2;
881 
882  V0[celli] += savedV0[oldCelli];
883 
884  nMerged++;
885  }
886  }
887 
888  DebugInfo
889  << "Mapping old time volume V0. Merged "
890  << nMerged << " out of " << nCells() << " cells" << endl;
891  }
892 
893 
894  // Map the old-old volume. Just map to new cell labels.
895  if (V00Ptr_)
896  {
897  scalarField& V00 = *V00Ptr_;
898 
899  scalarField savedV00(V00);
900  V00.resize_nocopy(nCells());
901 
902  forAll(V00, i)
903  {
904  if (cellMap[i] > -1)
905  {
906  V00[i] = savedV00[cellMap[i]];
907  }
908  else
909  {
910  V00[i] = 0.0;
911  }
912  }
913 
914  // Inject volume of merged cells
915  label nMerged = 0;
916  forAll(meshMap.reverseCellMap(), oldCelli)
917  {
918  label index = meshMap.reverseCellMap()[oldCelli];
919 
920  if (index < -1)
921  {
922  label celli = -index-2;
923 
924  V00[celli] += savedV00[oldCelli];
925  nMerged++;
926  }
927  }
928 
930  << "Mapping old time volume V00. Merged "
931  << nMerged << " out of " << nCells() << " cells" << endl;
932  }
933 }
934 
935 
937 {
939 
940  // Grab old time volumes if the time has been incremented
941  // This will update V0, V00
942  if (curTimeIndex_ < time().timeIndex())
943  {
944  storeOldVol(V());
945  }
946 
947 
948  // Move the polyMesh and initialise the mesh motion fluxes field
949  // Note: mesh flux updated by the fvGeometryScheme
950 
951  if (!phiPtr_)
952  {
953  DebugInFunction<< "Creating initial meshPhi field" << endl;
954 
955  // Create mesh motion flux
956  phiPtr_ = std::make_unique<surfaceScalarField>
957  (
958  IOobject
959  (
960  "meshPhi",
961  this->time().timeName(),
962  *this,
966  ),
967  *this,
969  );
970  }
971  else
972  {
973  // Grab old time mesh motion fluxes if the time has been incremented
974  if (phiPtr_->timeIndex() != time().timeIndex())
975  {
976  DebugInFunction<< "Accessing old-time meshPhi field" << endl;
977  phiPtr_->oldTime();
978  }
979  }
980 
982 
983  // Update or delete the local geometric properties as early as possible so
984  // they can be used if necessary. These get recreated here instead of
985  // demand driven since they might do parallel transfers which can conflict
986  // with when they're actually being used.
987  // Note that between above "polyMesh::movePoints(p)" and here nothing
988  // should use the local geometric properties.
989  updateGeomNotOldVol();
990 
991  // Update other local data
992  boundary_.movePoints();
993 
994  // Clear weights, deltaCoeffs, nonOrthoDeltaCoeffs, nonOrthCorrectionVectors
996 
997  meshObject::movePoints<fvMesh>(*this);
998  meshObject::movePoints<lduMesh>(*this);
999 }
1000 
1001 
1003 {
1004  DebugInFunction << endl;
1006  // Let surfaceInterpolation handle geometry calculation. Note: this does
1007  // lower levels updateGeom
1009 }
1010 
1011 
1012 void Foam::fvMesh::updateMesh(const mapPolyMesh& mpm)
1013 {
1014  DebugInFunction << endl;
1015 
1016  // Update polyMesh. This needs to keep volume existent!
1017  polyMesh::updateMesh(mpm);
1018 
1019  // Our slice of the addressing is no longer valid
1020  deleteDemandDrivenData(lduPtr_);
1021 
1022  if (VPtr_)
1023  {
1024  // Grab old time volumes if the time has been incremented
1025  // This will update V0, V00
1026  storeOldVol(mpm.oldCellVolumes());
1027 
1028  // Few checks
1029  if (VPtr_ && (VPtr_->size() != mpm.nOldCells()))
1030  {
1032  << "V:" << VPtr_->size()
1033  << " not equal to the number of old cells "
1034  << mpm.nOldCells()
1035  << exit(FatalError);
1036  }
1037  if (V0Ptr_ && (V0Ptr_->size() != mpm.nOldCells()))
1038  {
1040  << "V0:" << V0Ptr_->size()
1041  << " not equal to the number of old cells "
1042  << mpm.nOldCells()
1043  << exit(FatalError);
1044  }
1045  if (V00Ptr_ && (V00Ptr_->size() != mpm.nOldCells()))
1046  {
1048  << "V0:" << V00Ptr_->size()
1049  << " not equal to the number of old cells "
1050  << mpm.nOldCells()
1051  << exit(FatalError);
1052  }
1053  }
1054 
1055 
1056  // Clear mesh motion flux (note: could instead save & map like volumes)
1057  if (phiPtr_)
1058  {
1059  // Mesh moving and topology change. Recreate meshPhi
1060  phiPtr_.reset(nullptr);
1061 
1062  // Create mesh motion flux
1063  phiPtr_ = std::make_unique<surfaceScalarField>
1064  (
1065  IOobject
1066  (
1067  "meshPhi",
1068  this->time().timeName(),
1069  *this,
1073  ),
1074  *this,
1076  );
1077  }
1078 
1079  // Clear the sliced fields
1080  clearGeomNotOldVol();
1081 
1082  // Map all fields
1083  mapFields(mpm);
1084 
1085  // Clear the current volume and other geometry factors
1087 
1088  // Clear any non-updateable addressing
1089  clearAddressing(true);
1091  meshObject::updateMesh<fvMesh>(*this, mpm);
1092  meshObject::updateMesh<lduMesh>(*this, mpm);
1093 }
1094 
1095 
1097 (
1098  IOstreamOption streamOpt,
1099  const bool writeOnProc
1100 ) const
1101 {
1102  bool ok = true;
1103  if (phiPtr_)
1104  {
1105  ok = phiPtr_->write(writeOnProc);
1106  // NOTE: The old old time mesh phi might be necessary for certain
1107  // solver smooth restart using second order time schemes.
1108  //ok = phiPtr_->oldTime().write();
1109  }
1110  if (V0Ptr_ && V0Ptr_->writeOpt() == IOobject::AUTO_WRITE)
1111  {
1112  // For second order restarts we need to write V0
1113  ok = V0Ptr_->write(writeOnProc);
1114  }
1115 
1116  return ok && polyMesh::writeObject(streamOpt, writeOnProc);
1117 }
1118 
1119 
1120 bool Foam::fvMesh::write(const bool writeOnProc) const
1122  return polyMesh::write(writeOnProc);
1123 }
1124 
1125 
1126 template<>
1128 Foam::fvMesh::validComponents<Foam::sphericalTensor>() const
1131 }
1132 
1133 
1134 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1136 bool Foam::fvMesh::operator!=(const fvMesh& rhs) const
1137 {
1138  return &rhs != this;
1139 }
1140 
1141 
1142 bool Foam::fvMesh::operator==(const fvMesh& rhs) const
1143 {
1144  return &rhs == this;
1145 }
1146 
1147 
1148 // ************************************************************************* //
Foam::surfaceFields.
faceListList boundary
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
slicedVolVectorField * CPtr_
Cell centres.
Definition: fvMesh.H:137
void clearAddressing()
Clear topological data.
void removeFvBoundary()
Remove boundary patches. Warning: fvPatchFields hold ref to these fvPatches.
Definition: fvMesh.C:703
void clearMeshPhi()
Clear cell face motion fluxes.
Definition: fvMesh.C:235
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
Definition: polyMesh.C:346
const fvSchemes & schemes() const
Read-access to the fvSchemes controls.
Definition: fvMesh.C:593
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
virtual void movePoints(const pointField &)
Move points.
Definition: polyMesh.C:1168
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc) const
Write the underlying polyMesh and other data.
Definition: fvMesh.C:1090
static void clearUpto(objectRegistry &obr)
Clear all meshObject derived from FromType up to (but not including) ToType.
Definition: MeshObject.C:339
void mapClouds(const objectRegistry &db, const mapPolyMesh &mapper)
Generic Geometric field mapper.
Definition: mapClouds.H:48
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
void clearOutLocal(const bool isMeshUpdate=false)
Clear local-only storage (geometry, addressing etc)
Definition: fvMesh.C:215
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:58
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A traits class, which is primarily used for primitives and vector-space.
Definition: pTraits.H:75
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:53
virtual bool movePoints()
Do what is necessary if the mesh has moved.
void storeOldVol(const scalarField &)
Preserve old volume(s)
Definition: fvMesh.C:157
virtual const meshState & data() const noexcept
Const reference to the mesh and solver state data.
Definition: polyMesh.H:559
const dimensionSet dimVol(dimVolume)
Older spelling for dimVolume.
Definition: dimensionSets.H:65
const fvSolution & solution() const
Read-access to the fvSolution controls.
Definition: fvMesh.C:605
A simple container for options an IOstream can normally have.
Cell to surface interpolation scheme. Included in fvMesh.
Ignore writing from objectRegistry::writeObject()
SolverPerformance< Type > solveSegregatedOrCoupled(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:55
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
void MapGeometricFields(const MeshMapper &mapper)
Generic Geometric field mapper.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
void clearGeomNotOldVol()
Clear geometry but not the old-time cell volumes.
Definition: fvMesh.C:47
virtual ~fvMesh()
Destructor.
Definition: fvMesh.C:573
word timeName
Definition: getTimeIndex.H:3
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:1005
virtual void updateGeom()
Update all geometric data.
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:715
const cellShapeList & cells
const pointField & points
const labelList & cellMap() const noexcept
Old cell map.
Definition: mapPolyMesh.H:537
virtual void updateMesh(const mapPolyMesh &mpm)
Update the mesh corresponding to given map.
void clearOut()
Clear all geometry and addressing unnecessary for CFD.
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh for topology changes.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
bool operator!=(const fvMesh &rhs) const
Compares addresses.
Definition: fvMesh.C:1129
#define DebugInFunction
Report an information message using Foam::Info.
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:752
const labelList & reverseCellMap() const noexcept
Reverse cell map.
Definition: mapPolyMesh.H:656
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
void clearGeom()
Clear local geometry.
Definition: fvMesh.C:109
void reset(const meshState &ms)
Reset the dictionary.
Definition: meshState.C:95
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
virtual SolverPerformance< scalar > solve(fvMatrix< scalar > &, const dictionary &) const
Solve returning the solution statistics given convergence tolerance. Use the given solver controls...
Definition: fvMesh.C:618
errorManip< error > abort(error &err)
Definition: errorManip.H:139
fvMesh(const fvMesh &)=delete
No copy construct.
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1113
#define DebugInfo
Report an information message using Foam::Info.
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:32
bool operator==(const fvMesh &rhs) const
Compares addresses.
Definition: fvMesh.C:1135
int debug
Static debugging option.
defineTypeNameAndDebug(combustionModel, 0)
volScalarField & C
void addFvPatches(polyPatchList &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:673
label nOldCells() const noexcept
Number of old cells.
Definition: mapPolyMesh.H:472
Generic Geometric field mapper. For "real" mapping, add template specialisations for mapping of inter...
surfaceScalarField * magSfPtr_
Mag face area vectors.
Definition: fvMesh.H:132
Template functions to aid in the implementation of demand driven data.
const fvSchemes * hasSchemes() const
Non-null if fvSchemes exists (can test as bool).
Definition: fvMesh.C:581
const scalarField & oldCellVolumes() const
Definition: mapPolyMesh.H:805
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: fvMesh.C:769
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Selector class for finite volume solution solution. fvMesh is derived from fvSolution so that all fie...
Definition: fvSolution.H:51
void clearOut()
Clear all geometry and addressing.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Selector class for finite volume differencing schemes. fvMesh is derived from fvSchemes so that all f...
Definition: fvSchemes.H:51
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const fvSolution * hasSolution() const
Non-null if fvSolution exists (can test as bool).
Definition: fvMesh.C:587
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
Definition: fvMesh.C:271
Nothing to be read.
Automatically write from objectRegistry::writeObject()
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
slicedSurfaceVectorField * SfPtr_
Face area vectors.
Definition: fvMesh.H:127
void updateGeomNotOldVol()
Clear geometry like clearGeomNotOldVol but recreate any.
Definition: fvMesh.C:71
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
Reading is optional [identical to READ_IF_PRESENT].
label nOldFaces() const noexcept
Number of old faces.
Definition: mapPolyMesh.H:464
virtual void updateGeom()
Update all geometric data. This gets redirected up from primitiveMesh level.
Definition: fvMesh.C:995
The class contains the addressing required by the lduMatrix: upper, lower and losort.
SlicedDimensionedField< scalar, volMesh > * VPtr_
Cell volumes.
Definition: fvMesh.H:112
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:75
const labelList & faceMap() const noexcept
Old face map.
Definition: mapPolyMesh.H:503
SphericalTensor< scalar > sphericalTensor
SphericalTensor of scalars, i.e. SphericalTensor<scalar>.
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:91
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in.
Definition: polyMeshIO.C:68
Registry of regIOobjects.
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
void deleteDemandDrivenData(DataPtr &dataPtr)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
slicedSurfaceVectorField * CfPtr_
Face centres.
Definition: fvMesh.H:142
Do not request registration (bool: false)
Namespace for OpenFOAM.
label timeIndex
Definition: getTimeIndex.H:24
virtual void mapFields(const mapPolyMesh &mpm)
Map all fields in time using given map.
Definition: fvMesh.C:775
#define InfoInFunction
Report an information message using Foam::Info.
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc=true) const
Write items held in the objectRegistry. Normally includes mesh components (points, faces, etc) and any registered fields.
Definition: polyMesh.C:1592