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-2023 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 
223 {
224  clearGeom();
226 
227  clearAddressing();
228 
229  // Clear mesh motion flux
230  deleteDemandDrivenData(phiPtr_);
231 }
232 
233 
235 {
236  clearOutLocal();
238 }
239 
240 
241 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
242 
243 Foam::fvMesh::fvMesh(const IOobject& io, const bool doInit)
244 :
245  polyMesh(io, doInit),
246  fvSchemes(static_cast<const objectRegistry&>(*this)),
247  surfaceInterpolation(*this),
248  fvSolution(static_cast<const objectRegistry&>(*this)),
249  data(static_cast<const objectRegistry&>(*this)),
250  boundary_(*this, boundaryMesh()),
251  lduPtr_(nullptr),
252  curTimeIndex_(time().timeIndex()),
253  VPtr_(nullptr),
254  V0Ptr_(nullptr),
255  V00Ptr_(nullptr),
256  SfPtr_(nullptr),
257  magSfPtr_(nullptr),
258  CPtr_(nullptr),
259  CfPtr_(nullptr),
260  phiPtr_(nullptr)
261 {
262  DebugInFunction << "Constructing fvMesh from IOobject" << endl;
263 
264  if (doInit)
265  {
266  fvMesh::init(false); // do not initialise lower levels
267  }
268 }
269 
270 
271 bool Foam::fvMesh::init(const bool doInit)
272 {
273  if (doInit)
274  {
275  // Construct basic geometry calculation engine. Note: do before
276  // doing anything with primitiveMesh::cellCentres etc.
277  (void)geometry();
278 
279  // Initialise my data
280  polyMesh::init(doInit);
281  }
282 
283  // Read some optional fields
284  // ~~~~~~~~~~~~~~~~~~~~~~~~~
285  // For redistributePar + fileHandler can have master processor
286  // find the file but not the subprocs.
287 
288  IOobject rio
289  (
290  "none",
291  this->time().timeName(),
292  *this,
296  );
297 
298 
299  // Read old cell volumes (if present) and set the storage of V00
300  rio.resetHeader("V0");
301  if (returnReduceOr(rio.typeHeaderOk<regIOobject>(false)))
302  {
304  << "Detected V0: " << rio.objectRelPath() << nl;
305 
306  V0Ptr_ = new DimensionedField<scalar, volMesh>
307  (
308  rio,
309  *this,
311  );
312 
313  // Set the moving flag early in case demand-driven geometry
314  // construction checks for it
315  moving(true);
316 
317  (void)V00();
318  }
319 
320 
321  // Read mesh fluxes (if present) and set the mesh to be moving
322  rio.resetHeader("meshPhi");
323  if (returnReduceOr(rio.typeHeaderOk<regIOobject>(false)))
324  {
326  << "Detected meshPhi: " << rio.objectRelPath() << nl;
327 
328  phiPtr_ = new surfaceScalarField
329  (
330  rio,
331  *this,
333  );
334 
335  // Set the moving flag early in case demand-driven geometry
336  // construction checks for it
337  moving(true);
338 
339  // The mesh is now considered moving so the old-time cell volumes
340  // will be required for the time derivatives so if they haven't been
341  // read initialise to the current cell volumes
342  if (!V0Ptr_)
343  {
344  V0Ptr_ = new DimensionedField<scalar, volMesh>
345  (
346  IOobject
347  (
348  "V0",
349  this->time().timeName(),
350  *this,
354  ),
355  V()
356  );
357  }
358  }
360  // Assume something changed
361  return true;
362 }
363 
364 
366 (
367  const IOobject& io,
368  pointField&& points,
369  faceList&& faces,
370  labelList&& allOwner,
371  labelList&& allNeighbour,
372  const bool syncPar
373 )
374 :
375  polyMesh
376  (
377  io,
378  std::move(points),
379  std::move(faces),
380  std::move(allOwner),
381  std::move(allNeighbour),
382  syncPar
383  ),
384  fvSchemes(static_cast<const objectRegistry&>(*this)),
385  surfaceInterpolation(*this),
386  fvSolution(static_cast<const objectRegistry&>(*this)),
387  data(static_cast<const objectRegistry&>(*this)),
388  boundary_(*this),
389  lduPtr_(nullptr),
390  curTimeIndex_(time().timeIndex()),
391  VPtr_(nullptr),
392  V0Ptr_(nullptr),
393  V00Ptr_(nullptr),
394  SfPtr_(nullptr),
395  magSfPtr_(nullptr),
396  CPtr_(nullptr),
397  CfPtr_(nullptr),
398  phiPtr_(nullptr)
399 {
400  DebugInFunction << "Constructing fvMesh from components" << endl;
401 }
402 
403 
405 (
406  const IOobject& io,
407  pointField&& points,
408  faceList&& faces,
409  cellList&& cells,
410  const bool syncPar
411 )
412 :
413  polyMesh
414  (
415  io,
416  std::move(points),
417  std::move(faces),
418  std::move(cells),
419  syncPar
420  ),
421  fvSchemes(static_cast<const objectRegistry&>(*this)),
422  surfaceInterpolation(*this),
423  fvSolution(static_cast<const objectRegistry&>(*this)),
424  data(static_cast<const objectRegistry&>(*this)),
425  boundary_(*this),
426  lduPtr_(nullptr),
427  curTimeIndex_(time().timeIndex()),
428  VPtr_(nullptr),
429  V0Ptr_(nullptr),
430  V00Ptr_(nullptr),
431  SfPtr_(nullptr),
432  magSfPtr_(nullptr),
433  CPtr_(nullptr),
434  CfPtr_(nullptr),
435  phiPtr_(nullptr)
436 {
437  DebugInFunction << "Constructing fvMesh from components" << endl;
438 }
439 
440 
441 Foam::fvMesh::fvMesh(const IOobject& io, const Foam::zero, const bool syncPar)
442 :
443  fvMesh(io, pointField(), faceList(), labelList(), labelList(), syncPar)
444 {}
445 
446 
448 (
449  const IOobject& io,
450  const fvMesh& baseMesh,
451  const Foam::zero,
452  const bool syncPar
453 )
454 :
455  fvMesh
456  (
457  io,
458  baseMesh,
459  pointField(),
460  faceList(),
461  labelList(), // owner
462  labelList(), // neighbour
463  syncPar
464  )
465 {}
466 
467 
469 (
470  const IOobject& io,
471  const fvMesh& baseMesh,
472  pointField&& points,
473  faceList&& faces,
474  labelList&& allOwner,
475  labelList&& allNeighbour,
476  const bool syncPar
477 )
478 :
479  polyMesh
480  (
481  io,
482  std::move(points),
483  std::move(faces),
484  std::move(allOwner),
485  std::move(allNeighbour),
486  syncPar
487  ),
488  fvSchemes
489  (
490  static_cast<const objectRegistry&>(*this),
491  static_cast<const fvSchemes&>(baseMesh)
492  ),
493  surfaceInterpolation(*this),
494  fvSolution
495  (
496  static_cast<const objectRegistry&>(*this),
497  static_cast<const fvSolution&>(baseMesh)
498  ),
499  data
500  (
501  static_cast<const objectRegistry&>(*this),
502  static_cast<const data&>(baseMesh)
503  ),
504  boundary_(*this),
505  lduPtr_(nullptr),
506  curTimeIndex_(time().timeIndex()),
507  VPtr_(nullptr),
508  V0Ptr_(nullptr),
509  V00Ptr_(nullptr),
510  SfPtr_(nullptr),
511  magSfPtr_(nullptr),
512  CPtr_(nullptr),
513  CfPtr_(nullptr),
514  phiPtr_(nullptr)
515 {
516  DebugInFunction << "Constructing fvMesh as copy and primitives" << endl;
517 }
518 
519 
521 (
522  const IOobject& io,
523  const fvMesh& baseMesh,
524  pointField&& points,
525  faceList&& faces,
526  cellList&& cells,
527  const bool syncPar
528 )
529 :
530  polyMesh
531  (
532  io,
533  std::move(points),
534  std::move(faces),
535  std::move(cells),
536  syncPar
537  ),
538  fvSchemes
539  (
540  static_cast<const objectRegistry&>(*this),
541  static_cast<const fvSchemes&>(baseMesh)
542  ),
543  surfaceInterpolation(*this),
544  fvSolution
545  (
546  static_cast<const objectRegistry&>(*this),
547  static_cast<const fvSolution&>(baseMesh)
548  ),
549  data
550  (
551  static_cast<const objectRegistry&>(*this),
552  static_cast<const data&>(baseMesh)
553  ),
554  boundary_(*this),
555  lduPtr_(nullptr),
556  curTimeIndex_(time().timeIndex()),
557  VPtr_(nullptr),
558  V0Ptr_(nullptr),
559  V00Ptr_(nullptr),
560  SfPtr_(nullptr),
561  magSfPtr_(nullptr),
562  CPtr_(nullptr),
563  CfPtr_(nullptr),
564  phiPtr_(nullptr)
565 {
566  DebugInFunction << "Constructing fvMesh as copy and primitives" << endl;
567 }
568 
569 
570 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
571 
573 {
574  clearOut();
575 }
576 
577 
578 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
579 
581 (
582  fvMatrix<scalar>& m,
583  const dictionary& dict
584 ) const
585 {
586  // Redirect to fvMatrix solver
587  return m.solveSegregatedOrCoupled(dict);
588 }
589 
590 
592 (
593  fvMatrix<vector>& m,
594  const dictionary& dict
595 ) const
596 {
597  // Redirect to fvMatrix solver
598  return m.solveSegregatedOrCoupled(dict);
599 }
600 
601 
603 (
605  const dictionary& dict
606 ) const
607 {
608  // Redirect to fvMatrix solver
609  return m.solveSegregatedOrCoupled(dict);
610 }
611 
612 
614 (
616  const dictionary& dict
617 ) const
618 {
619  // Redirect to fvMatrix solver
620  return m.solveSegregatedOrCoupled(dict);
621 }
622 
623 
625 (
626  fvMatrix<tensor>& m,
627  const dictionary& dict
628 ) const
629 {
630  // Redirect to fvMatrix solver
631  return m.solveSegregatedOrCoupled(dict);
632 }
633 
634 
636 (
637  polyPatchList& plist,
638  const bool validBoundary
639 )
640 {
641  if (boundary().size())
642  {
644  << " boundary already exists"
645  << abort(FatalError);
646  }
648  addPatches(plist, validBoundary);
649  boundary_.addPatches(boundaryMesh());
650 }
651 
652 
654 (
655  const List<polyPatch*>& p,
656  const bool validBoundary
657 )
658 {
659  // Acquire ownership of the pointers
660  polyPatchList plist(const_cast<List<polyPatch*>&>(p));
661 
662  addFvPatches(plist, validBoundary);
663 }
664 
665 
667 {
668  DebugInFunction << "Removing boundary patches." << endl;
669 
670  // Remove fvBoundaryMesh data first.
671  boundary_.clear();
673 
674  clearOut();
675 }
676 
677 
679 {
680  DebugInFunction << "Updating fvMesh";
681 
683 
684  if (state == polyMesh::TOPO_PATCH_CHANGE)
685  {
686  DebugInfo << "Boundary and topological update" << endl;
687 
688  boundary_.readUpdate(boundaryMesh());
689 
690  clearOut();
691 
692  }
693  else if (state == polyMesh::TOPO_CHANGE)
694  {
695  DebugInfo << "Topological update" << endl;
696 
697  // fvMesh::clearOut() but without the polyMesh::clearOut
698  clearOutLocal();
699  }
700  else if (state == polyMesh::POINTS_MOVED)
701  {
702  DebugInfo << "Point motion update" << endl;
703 
704  clearGeom();
705  }
706  else
707  {
708  DebugInfo << "No update" << endl;
709  }
710 
711  return state;
712 }
713 
714 
716 {
717  if (!lduPtr_)
718  {
720  << "Calculating fvMeshLduAddressing from nFaces:"
721  << nFaces() << endl;
722 
723  lduPtr_ = new fvMeshLduAddressing(*this);
724 
725  return *lduPtr_;
726  }
727 
728  return *lduPtr_;
729 }
730 
733 {
734  return boundary().interfaces();
735 }
736 
737 
738 void Foam::fvMesh::mapFields(const mapPolyMesh& meshMap)
739 {
741  << " nOldCells:" << meshMap.nOldCells()
742  << " nCells:" << nCells()
743  << " nOldFaces:" << meshMap.nOldFaces()
744  << " nFaces:" << nFaces()
745  << endl;
746 
747  // We require geometric properties valid for the old mesh
748  if
749  (
750  meshMap.cellMap().size() != nCells()
751  || meshMap.faceMap().size() != nFaces()
752  )
753  {
755  << "mapPolyMesh does not correspond to the old mesh."
756  << " nCells:" << nCells()
757  << " cellMap:" << meshMap.cellMap().size()
758  << " nOldCells:" << meshMap.nOldCells()
759  << " nFaces:" << nFaces()
760  << " faceMap:" << meshMap.faceMap().size()
761  << " nOldFaces:" << meshMap.nOldFaces()
762  << exit(FatalError);
763  }
764 
765  // Create a mapper
766  const fvMeshMapper mapper(*this, meshMap);
767 
768  // Map all the volFields in the objectRegistry
769  MapGeometricFields<scalar, fvPatchField, fvMeshMapper, volMesh>
770  (mapper);
771  MapGeometricFields<vector, fvPatchField, fvMeshMapper, volMesh>
772  (mapper);
773  MapGeometricFields<sphericalTensor, fvPatchField, fvMeshMapper, volMesh>
774  (mapper);
775  MapGeometricFields<symmTensor, fvPatchField, fvMeshMapper, volMesh>
776  (mapper);
777  MapGeometricFields<tensor, fvPatchField, fvMeshMapper, volMesh>
778  (mapper);
779 
780  // Map all the surfaceFields in the objectRegistry
781  MapGeometricFields<scalar, fvsPatchField, fvMeshMapper, surfaceMesh>
782  (mapper);
783  MapGeometricFields<vector, fvsPatchField, fvMeshMapper, surfaceMesh>
784  (mapper);
786  <
787  sphericalTensor, fvsPatchField, fvMeshMapper, surfaceMesh
788  >
789  (mapper);
790  MapGeometricFields<symmTensor, fvsPatchField, fvMeshMapper, surfaceMesh>
791  (mapper);
792  MapGeometricFields<tensor, fvsPatchField, fvMeshMapper, surfaceMesh>
793  (mapper);
794 
795  // Map all the dimensionedFields in the objectRegistry
796  MapDimensionedFields<scalar, fvMeshMapper, volMesh>(mapper);
797  MapDimensionedFields<vector, fvMeshMapper, volMesh>(mapper);
798  MapDimensionedFields<sphericalTensor, fvMeshMapper, volMesh>(mapper);
799  MapDimensionedFields<symmTensor, fvMeshMapper, volMesh>(mapper);
800  MapDimensionedFields<tensor, fvMeshMapper, volMesh>(mapper);
801 
802  // Map all the clouds in the objectRegistry
803  mapClouds(*this, meshMap);
804 
805 
806  const labelList& cellMap = meshMap.cellMap();
807 
808  // Map the old volume. Just map to new cell labels.
809  if (V0Ptr_)
810  {
811  scalarField& V0 = *V0Ptr_;
812 
813  scalarField savedV0(V0);
814  V0.resize_nocopy(nCells());
815 
816  forAll(V0, i)
817  {
818  if (cellMap[i] > -1)
819  {
820  V0[i] = savedV0[cellMap[i]];
821  }
822  else
823  {
824  V0[i] = 0.0;
825  }
826  }
827 
828  // Inject volume of merged cells
829  label nMerged = 0;
830  forAll(meshMap.reverseCellMap(), oldCelli)
831  {
832  label index = meshMap.reverseCellMap()[oldCelli];
833 
834  if (index < -1)
835  {
836  label celli = -index-2;
837 
838  V0[celli] += savedV0[oldCelli];
839 
840  nMerged++;
841  }
842  }
843 
844  DebugInfo
845  << "Mapping old time volume V0. Merged "
846  << nMerged << " out of " << nCells() << " cells" << endl;
847  }
848 
849 
850  // Map the old-old volume. Just map to new cell labels.
851  if (V00Ptr_)
852  {
853  scalarField& V00 = *V00Ptr_;
854 
855  scalarField savedV00(V00);
856  V00.resize_nocopy(nCells());
857 
858  forAll(V00, i)
859  {
860  if (cellMap[i] > -1)
861  {
862  V00[i] = savedV00[cellMap[i]];
863  }
864  else
865  {
866  V00[i] = 0.0;
867  }
868  }
869 
870  // Inject volume of merged cells
871  label nMerged = 0;
872  forAll(meshMap.reverseCellMap(), oldCelli)
873  {
874  label index = meshMap.reverseCellMap()[oldCelli];
875 
876  if (index < -1)
877  {
878  label celli = -index-2;
879 
880  V00[celli] += savedV00[oldCelli];
881  nMerged++;
882  }
883  }
884 
886  << "Mapping old time volume V00. Merged "
887  << nMerged << " out of " << nCells() << " cells" << endl;
888  }
889 }
890 
891 
893 {
895 
896  // Grab old time volumes if the time has been incremented
897  // This will update V0, V00
898  if (curTimeIndex_ < time().timeIndex())
899  {
900  storeOldVol(V());
901  }
902 
903 
904  // Move the polyMesh and initialise the mesh motion fluxes field
905  // Note: mesh flux updated by the fvGeometryScheme
906 
907  if (!phiPtr_)
908  {
909  DebugInFunction<< "Creating initial meshPhi field" << endl;
910 
911  // Create mesh motion flux
912  phiPtr_ = new surfaceScalarField
913  (
914  IOobject
915  (
916  "meshPhi",
917  this->time().timeName(),
918  *this,
922  ),
923  *this,
925  );
926  }
927  else
928  {
929  // Grab old time mesh motion fluxes if the time has been incremented
930  if (phiPtr_->timeIndex() != time().timeIndex())
931  {
932  DebugInFunction<< "Accessing old-time meshPhi field" << endl;
933  phiPtr_->oldTime();
934  }
935  }
936 
938 
939  // Update or delete the local geometric properties as early as possible so
940  // they can be used if necessary. These get recreated here instead of
941  // demand driven since they might do parallel transfers which can conflict
942  // with when they're actually being used.
943  // Note that between above "polyMesh::movePoints(p)" and here nothing
944  // should use the local geometric properties.
945  updateGeomNotOldVol();
946 
947  // Update other local data
948  boundary_.movePoints();
949 
950  // Clear weights, deltaCoeffs, nonOrthoDeltaCoeffs, nonOrthCorrectionVectors
952 
953  meshObject::movePoints<fvMesh>(*this);
954  meshObject::movePoints<lduMesh>(*this);
955 }
956 
957 
959 {
962  // Let surfaceInterpolation handle geometry calculation. Note: this does
963  // lower levels updateGeom
965 }
966 
967 
968 void Foam::fvMesh::updateMesh(const mapPolyMesh& mpm)
969 {
971 
972  // Update polyMesh. This needs to keep volume existent!
974 
975  // Our slice of the addressing is no longer valid
976  deleteDemandDrivenData(lduPtr_);
977 
978  if (VPtr_)
979  {
980  // Grab old time volumes if the time has been incremented
981  // This will update V0, V00
982  storeOldVol(mpm.oldCellVolumes());
983 
984  // Few checks
985  if (VPtr_ && (VPtr_->size() != mpm.nOldCells()))
986  {
988  << "V:" << VPtr_->size()
989  << " not equal to the number of old cells "
990  << mpm.nOldCells()
991  << exit(FatalError);
992  }
993  if (V0Ptr_ && (V0Ptr_->size() != mpm.nOldCells()))
994  {
996  << "V0:" << V0Ptr_->size()
997  << " not equal to the number of old cells "
998  << mpm.nOldCells()
999  << exit(FatalError);
1000  }
1001  if (V00Ptr_ && (V00Ptr_->size() != mpm.nOldCells()))
1002  {
1004  << "V0:" << V00Ptr_->size()
1005  << " not equal to the number of old cells "
1006  << mpm.nOldCells()
1007  << exit(FatalError);
1008  }
1009  }
1010 
1011 
1012  // Clear mesh motion flux (note: could instead save & map like volumes)
1013  if (phiPtr_)
1014  {
1015  // Mesh moving and topology change. Recreate meshPhi
1016  deleteDemandDrivenData(phiPtr_);
1017 
1018  // Create mesh motion flux
1019  phiPtr_ = new surfaceScalarField
1020  (
1021  IOobject
1022  (
1023  "meshPhi",
1024  this->time().timeName(),
1025  *this,
1029  ),
1030  *this,
1032  );
1033  }
1034 
1035  // Clear the sliced fields
1036  clearGeomNotOldVol();
1037 
1038  // Map all fields
1039  mapFields(mpm);
1040 
1041  // Clear the current volume and other geometry factors
1043 
1044  // Clear any non-updateable addressing
1045  clearAddressing(true);
1047  meshObject::updateMesh<fvMesh>(*this, mpm);
1048  meshObject::updateMesh<lduMesh>(*this, mpm);
1049 }
1050 
1051 
1053 (
1054  IOstreamOption streamOpt,
1055  const bool writeOnProc
1056 ) const
1057 {
1058  bool ok = true;
1059  if (phiPtr_)
1060  {
1061  ok = phiPtr_->write(writeOnProc);
1062  // NOTE: The old old time mesh phi might be necessary for certain
1063  // solver smooth restart using second order time schemes.
1064  //ok = phiPtr_->oldTime().write();
1065  }
1066  if (V0Ptr_ && V0Ptr_->writeOpt() == IOobject::AUTO_WRITE)
1067  {
1068  // For second order restarts we need to write V0
1069  ok = V0Ptr_->write(writeOnProc);
1070  }
1071 
1072  return ok && polyMesh::writeObject(streamOpt, writeOnProc);
1073 }
1074 
1075 
1076 bool Foam::fvMesh::write(const bool writeOnProc) const
1078  return polyMesh::write(writeOnProc);
1079 }
1080 
1081 
1082 template<>
1084 Foam::fvMesh::validComponents<Foam::sphericalTensor>() const
1087 }
1088 
1089 
1090 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1092 bool Foam::fvMesh::operator!=(const fvMesh& rhs) const
1093 {
1094  return &rhs != this;
1095 }
1096 
1097 
1098 bool Foam::fvMesh::operator==(const fvMesh& rhs) const
1099 {
1100  return &rhs == this;
1101 }
1102 
1103 
1104 // ************************************************************************* //
Foam::surfaceFields.
faceListList boundary
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:227
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
slicedVolVectorField * CPtr_
Cell centres.
Definition: fvMesh.H:139
void clearAddressing()
Clear topological data.
void removeFvBoundary()
Remove boundary patches. Warning: fvPatchFields hold ref to these fvPatches.
Definition: fvMesh.C:659
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
Definition: polyMesh.C:359
label nOldCells() const
Number of old cells.
Definition: mapPolyMesh.H:469
const labelList & cellMap() const
Old cell map.
Definition: mapPolyMesh.H:534
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
virtual void movePoints(const pointField &)
Move points.
Definition: polyMesh.C:1163
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc) const
Write the underlying polyMesh and other data.
Definition: fvMesh.C:1046
static void clearUpto(objectRegistry &obr)
Clear all meshObject derived from FromType up to (but not including) ToType.
Definition: MeshObject.C:257
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:120
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
void clearOut()
Clear all geometry and addressing unnecessary for CFD.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
void clearOutLocal()
Clear local-only storage (geometry, addressing etc)
Definition: fvMesh.C:215
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:487
A traits class, which is primarily used for primitives.
Definition: pTraits.H:51
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
const dimensionSet dimVol(dimVolume)
Older spelling for dimVolume.
Definition: dimensionSets.H:65
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:157
void MapGeometricFields(const MeshMapper &mapper)
Generic Geometric field mapper.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc) const
Write the objects using stream options.
void clearGeomNotOldVol()
Clear geometry but not the old-time cell volumes.
Definition: fvMesh.C:47
virtual ~fvMesh()
Destructor.
Definition: fvMesh.C:565
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:961
virtual void updateGeom()
Update all geometric data.
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:671
const cellShapeList & cells
const pointField & points
virtual void updateMesh(const mapPolyMesh &mpm)
Update the mesh corresponding to given map.
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:1085
#define DebugInFunction
Report an information message using Foam::Info.
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:708
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
void clearGeom()
Clear local geometry.
Definition: fvMesh.C:109
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:574
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:1069
#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:1091
const labelList & reverseCellMap() const
Reverse cell map.
Definition: mapPolyMesh.H:653
int debug
Static debugging option.
defineTypeNameAndDebug(combustionModel, 0)
Database for solution data, solver performance and other reduced data.
Definition: data.H:52
volScalarField & C
void addFvPatches(polyPatchList &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:629
const labelList & faceMap() const
Old face map.
Definition: mapPolyMesh.H:500
Generic Geometric field mapper. For "real" mapping, add template specialisations for mapping of inter...
surfaceScalarField * magSfPtr_
Mag face area vectors.
Definition: fvMesh.H:134
Template functions to aid in the implementation of demand driven data.
const scalarField & oldCellVolumes() const
Definition: mapPolyMesh.H:802
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: fvMesh.C:725
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:79
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
Definition: fvMesh.C:264
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:129
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].
virtual void updateGeom()
Update all geometric data. This gets redirected up from primitiveMesh level.
Definition: fvMesh.C:951
The class contains the addressing required by the lduMatrix: upper, lower and losort.
SlicedDimensionedField< scalar, volMesh > * VPtr_
Cell volumes.
Definition: fvMesh.H:114
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:73
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:89
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
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
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:171
label nOldFaces() const
Number of old faces.
Definition: mapPolyMesh.H:461
slicedSurfaceVectorField * CfPtr_
Face centres.
Definition: fvMesh.H:144
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:731
#define InfoInFunction
Report an information message using Foam::Info.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133