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