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 OpenFOAM Foundation
9  Copyright (C) 2016-2020 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 {
48  defineTypeNameAndDebug(fvMesh, 0);
49 }
50 
51 
52 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53 
55 {
57  <
58  fvMesh,
61  >(*this);
62 
64  <
65  lduMesh,
68  >(*this);
69 
73  VPtr_ = nullptr;
74 
79 }
80 
81 
83 {
84  bool haveV = (VPtr_ != nullptr);
85  bool haveSf = (SfPtr_ != nullptr);
86  bool haveMagSf = (magSfPtr_ != nullptr);
87  bool haveCP = (CPtr_ != nullptr);
88  bool haveCf = (CfPtr_ != nullptr);
89 
90  clearGeomNotOldVol();
91 
92  // Now recreate the fields
93  if (haveV)
94  {
95  (void)V();
96  }
97  if (haveSf)
98  {
99  (void)Sf();
100  }
101  if (haveMagSf)
102  {
103  (void)magSf();
104  }
105  if (haveCP)
106  {
107  (void)C();
108  }
109  if (haveCf)
110  {
111  (void)Cf();
112  }
113 }
114 
115 
117 {
118  clearGeomNotOldVol();
119 
120  deleteDemandDrivenData(V0Ptr_);
121  deleteDemandDrivenData(V00Ptr_);
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  // Copy V into V0 storage
182  V0Ptr_->scalarField::operator=(V);
183  }
184  else
185  {
186  // Allocate V0 storage, fill with V
188  (
189  IOobject
190  (
191  "V0",
192  time().timeName(),
193  *this,
196  false
197  ),
198  *this,
199  dimVolume
200  );
201  scalarField& V0 = *V0Ptr_;
202  // Note: V0 now sized with current mesh, not with (potentially
203  // different size) V.
204  V0.setSize(V.size());
205  V0 = V;
206  }
207 
208  curTimeIndex_ = time().timeIndex();
209 
210  if (debug)
211  {
213  << " Stored old time volumes V0:" << V0Ptr_->size()
214  << endl;
215 
216  if (V00Ptr_)
217  {
219  << " Stored oldold time volumes V00:" << V00Ptr_->size()
220  << endl;
221  }
222  }
223  }
224 }
225 
226 
228 {
229  clearGeom();
231 
232  clearAddressing();
233 
234  // Clear mesh motion flux
235  deleteDemandDrivenData(phiPtr_);
236 
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  // Intialise my data
280  polyMesh::init(doInit);
281  }
282 
283  // Check the existence of the cell volumes and read if present
284  // and set the storage of V00
285  if (fileHandler().isFile(time().timePath()/dbDir()/"V0"))
286  {
288  (
289  IOobject
290  (
291  "V0",
292  time().timeName(),
293  *this,
296  false
297  ),
298  *this
299  );
300 
301  V00();
302  }
303 
304  // Check the existence of the mesh fluxes, read if present and set the
305  // mesh to be moving
306  if (fileHandler().isFile(time().timePath()/dbDir()/"meshPhi"))
307  {
308  phiPtr_ = new surfaceScalarField
309  (
310  IOobject
311  (
312  "meshPhi",
313  time().timeName(),
314  *this,
317  false
318  ),
319  *this
320  );
321 
322  // The mesh is now considered moving so the old-time cell volumes
323  // will be required for the time derivatives so if they haven't been
324  // read initialise to the current cell volumes
325  if (!V0Ptr_)
326  {
328  (
329  IOobject
330  (
331  "V0",
332  time().timeName(),
333  *this,
336  false
337  ),
338  V()
339  );
340  }
341 
342  moving(true);
343  }
344 
345  // Assume something changed
346  return true;
347 }
348 
349 
351 (
352  const IOobject& io,
353  pointField&& points,
354  faceList&& faces,
355  labelList&& allOwner,
356  labelList&& allNeighbour,
357  const bool syncPar
358 )
359 :
360  polyMesh
361  (
362  io,
363  std::move(points),
364  std::move(faces),
365  std::move(allOwner),
366  std::move(allNeighbour),
367  syncPar
368  ),
369  fvSchemes(static_cast<const objectRegistry&>(*this)),
370  surfaceInterpolation(*this),
371  fvSolution(static_cast<const objectRegistry&>(*this)),
372  data(static_cast<const objectRegistry&>(*this)),
373  boundary_(*this),
374  lduPtr_(nullptr),
375  curTimeIndex_(time().timeIndex()),
376  VPtr_(nullptr),
377  V0Ptr_(nullptr),
378  V00Ptr_(nullptr),
379  SfPtr_(nullptr),
380  magSfPtr_(nullptr),
381  CPtr_(nullptr),
382  CfPtr_(nullptr),
383  phiPtr_(nullptr)
384 {
385  DebugInFunction << "Constructing fvMesh from components" << endl;
386 }
387 
388 
390 (
391  const IOobject& io,
392  pointField&& points,
393  faceList&& faces,
394  cellList&& cells,
395  const bool syncPar
396 )
397 :
398  polyMesh
399  (
400  io,
401  std::move(points),
402  std::move(faces),
403  std::move(cells),
404  syncPar
405  ),
406  fvSchemes(static_cast<const objectRegistry&>(*this)),
407  surfaceInterpolation(*this),
408  fvSolution(static_cast<const objectRegistry&>(*this)),
409  data(static_cast<const objectRegistry&>(*this)),
410  boundary_(*this),
411  lduPtr_(nullptr),
412  curTimeIndex_(time().timeIndex()),
413  VPtr_(nullptr),
414  V0Ptr_(nullptr),
415  V00Ptr_(nullptr),
416  SfPtr_(nullptr),
417  magSfPtr_(nullptr),
418  CPtr_(nullptr),
419  CfPtr_(nullptr),
420  phiPtr_(nullptr)
421 {
422  DebugInFunction << "Constructing fvMesh from components" << endl;
423 }
424 
425 
426 Foam::fvMesh::fvMesh(const IOobject& io, const zero, const bool syncPar)
427 :
428  fvMesh(io, pointField(), faceList(), labelList(), labelList(), syncPar)
429 {}
430 
431 
433 (
434  const IOobject& io,
435  const fvMesh& baseMesh,
436  pointField&& points,
437  faceList&& faces,
438  labelList&& allOwner,
439  labelList&& allNeighbour,
440  const bool syncPar
441 )
442 :
443  polyMesh
444  (
445  io,
446  std::move(points),
447  std::move(faces),
448  std::move(allOwner),
449  std::move(allNeighbour),
450  syncPar
451  ),
452  fvSchemes
453  (
454  static_cast<const objectRegistry&>(*this),
455  static_cast<const fvSchemes&>(baseMesh)
456  ),
457  surfaceInterpolation(*this),
458  fvSolution
459  (
460  static_cast<const objectRegistry&>(*this),
461  static_cast<const fvSolution&>(baseMesh)
462  ),
463  data
464  (
465  static_cast<const objectRegistry&>(*this),
466  static_cast<const data&>(baseMesh)
467  ),
468  boundary_(*this),
469  lduPtr_(nullptr),
470  curTimeIndex_(time().timeIndex()),
471  VPtr_(nullptr),
472  V0Ptr_(nullptr),
473  V00Ptr_(nullptr),
474  SfPtr_(nullptr),
475  magSfPtr_(nullptr),
476  CPtr_(nullptr),
477  CfPtr_(nullptr),
478  phiPtr_(nullptr)
479 {
480  DebugInFunction << "Constructing fvMesh as copy and primitives" << endl;
481 }
482 
483 
485 (
486  const IOobject& io,
487  const fvMesh& baseMesh,
488  pointField&& points,
489  faceList&& faces,
490  cellList&& cells,
491  const bool syncPar
492 )
493 :
494  polyMesh
495  (
496  io,
497  std::move(points),
498  std::move(faces),
499  std::move(cells),
500  syncPar
501  ),
502  fvSchemes
503  (
504  static_cast<const objectRegistry&>(*this),
505  static_cast<const fvSchemes&>(baseMesh)
506  ),
507  surfaceInterpolation(*this),
508  fvSolution
509  (
510  static_cast<const objectRegistry&>(*this),
511  static_cast<const fvSolution&>(baseMesh)
512  ),
513  data
514  (
515  static_cast<const objectRegistry&>(*this),
516  static_cast<const data&>(baseMesh)
517  ),
518  boundary_(*this),
519  lduPtr_(nullptr),
520  curTimeIndex_(time().timeIndex()),
521  VPtr_(nullptr),
522  V0Ptr_(nullptr),
523  V00Ptr_(nullptr),
524  SfPtr_(nullptr),
525  magSfPtr_(nullptr),
526  CPtr_(nullptr),
527  CfPtr_(nullptr),
528  phiPtr_(nullptr)
529 {
530  DebugInFunction << "Constructing fvMesh as copy and primitives" << endl;
531 }
532 
533 
534 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
535 
537 {
538  clearOut();
539 }
540 
541 
542 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
543 
545 (
546  fvMatrix<scalar>& m,
547  const dictionary& dict
548 ) const
549 {
550  // Redirect to fvMatrix solver
551  return m.solveSegregatedOrCoupled(dict);
552 }
553 
554 
556 (
557  fvMatrix<vector>& m,
558  const dictionary& dict
559 ) const
560 {
561  // Redirect to fvMatrix solver
562  return m.solveSegregatedOrCoupled(dict);
563 }
564 
565 
567 (
569  const dictionary& dict
570 ) const
571 {
572  // Redirect to fvMatrix solver
573  return m.solveSegregatedOrCoupled(dict);
574 }
575 
576 
578 (
580  const dictionary& dict
581 ) const
582 {
583  // Redirect to fvMatrix solver
584  return m.solveSegregatedOrCoupled(dict);
585 }
586 
587 
589 (
590  fvMatrix<tensor>& m,
591  const dictionary& dict
592 ) const
593 {
594  // Redirect to fvMatrix solver
595  return m.solveSegregatedOrCoupled(dict);
596 }
597 
598 
600 (
601  PtrList<polyPatch>& plist,
602  const bool validBoundary
603 )
604 {
605  if (boundary().size())
606  {
608  << " boundary already exists"
609  << abort(FatalError);
610  }
611 
612  addPatches(plist, validBoundary);
613  boundary_.addPatches(boundaryMesh());
614 }
615 
616 
618 (
619  const List<polyPatch*>& p,
620  const bool validBoundary
621 )
622 {
623  // Acquire ownership of the pointers
624  PtrList<polyPatch> plist(const_cast<List<polyPatch*>&>(p));
625 
626  addFvPatches(plist, validBoundary);
627 }
628 
629 
631 {
632  DebugInFunction << "Removing boundary patches." << endl;
633 
634  // Remove fvBoundaryMesh data first.
635  boundary_.clear();
636  boundary_.setSize(0);
638 
639  clearOut();
640 }
641 
642 
644 {
645  DebugInFunction << "Updating fvMesh. ";
646 
648 
649  if (state == polyMesh::TOPO_PATCH_CHANGE)
650  {
651  DebugInfo << "Boundary and topological update" << endl;
652 
653  boundary_.readUpdate(boundaryMesh());
654 
655  clearOut();
656 
657  }
658  else if (state == polyMesh::TOPO_CHANGE)
659  {
660  DebugInfo << "Topological update" << endl;
661 
662  clearOut();
663  }
664  else if (state == polyMesh::POINTS_MOVED)
665  {
666  DebugInfo << "Point motion update" << endl;
667 
668  clearGeom();
669  }
670  else
671  {
672  DebugInfo << "No update" << endl;
673  }
674 
675  return state;
676 }
677 
678 
680 {
681  return boundary_;
682 }
683 
684 
686 {
687  if (!lduPtr_)
688  {
690  << "Calculating fvMeshLduAddressing from nFaces:"
691  << nFaces() << endl;
692 
693  lduPtr_ = new fvMeshLduAddressing(*this);
694  }
695 
696  return *lduPtr_;
697 }
698 
699 
701 {
703  << " nOldCells:" << meshMap.nOldCells()
704  << " nCells:" << nCells()
705  << " nOldFaces:" << meshMap.nOldFaces()
706  << " nFaces:" << nFaces()
707  << endl;
708 
709  // We require geometric properties valid for the old mesh
710  if
711  (
712  meshMap.cellMap().size() != nCells()
713  || meshMap.faceMap().size() != nFaces()
714  )
715  {
717  << "mapPolyMesh does not correspond to the old mesh."
718  << " nCells:" << nCells()
719  << " cellMap:" << meshMap.cellMap().size()
720  << " nOldCells:" << meshMap.nOldCells()
721  << " nFaces:" << nFaces()
722  << " faceMap:" << meshMap.faceMap().size()
723  << " nOldFaces:" << meshMap.nOldFaces()
724  << exit(FatalError);
725  }
726 
727  // Create a mapper
728  const fvMeshMapper mapper(*this, meshMap);
729 
730  // Map all the volFields in the objectRegistry
731  MapGeometricFields<scalar, fvPatchField, fvMeshMapper, volMesh>
732  (mapper);
733  MapGeometricFields<vector, fvPatchField, fvMeshMapper, volMesh>
734  (mapper);
735  MapGeometricFields<sphericalTensor, fvPatchField, fvMeshMapper, volMesh>
736  (mapper);
737  MapGeometricFields<symmTensor, fvPatchField, fvMeshMapper, volMesh>
738  (mapper);
739  MapGeometricFields<tensor, fvPatchField, fvMeshMapper, volMesh>
740  (mapper);
741 
742  // Map all the surfaceFields in the objectRegistry
743  MapGeometricFields<scalar, fvsPatchField, fvMeshMapper, surfaceMesh>
744  (mapper);
745  MapGeometricFields<vector, fvsPatchField, fvMeshMapper, surfaceMesh>
746  (mapper);
747  MapGeometricFields<symmTensor, fvsPatchField, fvMeshMapper, surfaceMesh>
748  (mapper);
749  MapGeometricFields<symmTensor, fvsPatchField, fvMeshMapper, surfaceMesh>
750  (mapper);
751  MapGeometricFields<tensor, fvsPatchField, fvMeshMapper, surfaceMesh>
752  (mapper);
753 
754  // Map all the dimensionedFields in the objectRegistry
755  MapDimensionedFields<scalar, fvMeshMapper, volMesh>(mapper);
756  MapDimensionedFields<vector, fvMeshMapper, volMesh>(mapper);
757  MapDimensionedFields<sphericalTensor, fvMeshMapper, volMesh>(mapper);
758  MapDimensionedFields<symmTensor, fvMeshMapper, volMesh>(mapper);
759  MapDimensionedFields<tensor, fvMeshMapper, volMesh>(mapper);
760 
761  // Map all the clouds in the objectRegistry
762  mapClouds(*this, meshMap);
763 
764 
765  const labelList& cellMap = meshMap.cellMap();
766 
767  // Map the old volume. Just map to new cell labels.
768  if (V0Ptr_)
769  {
770  scalarField& V0 = *V0Ptr_;
771 
772  scalarField savedV0(V0);
773  V0.setSize(nCells());
774 
775  forAll(V0, i)
776  {
777  if (cellMap[i] > -1)
778  {
779  V0[i] = savedV0[cellMap[i]];
780  }
781  else
782  {
783  V0[i] = 0.0;
784  }
785  }
786 
787  // Inject volume of merged cells
788  label nMerged = 0;
789  forAll(meshMap.reverseCellMap(), oldCelli)
790  {
791  label index = meshMap.reverseCellMap()[oldCelli];
792 
793  if (index < -1)
794  {
795  label celli = -index-2;
796 
797  V0[celli] += savedV0[oldCelli];
798 
799  nMerged++;
800  }
801  }
802 
803  DebugInfo
804  << "Mapping old time volume V0. Merged "
805  << nMerged << " out of " << nCells() << " cells" << endl;
806  }
807 
808 
809  // Map the old-old volume. Just map to new cell labels.
810  if (V00Ptr_)
811  {
812  scalarField& V00 = *V00Ptr_;
813 
814  scalarField savedV00(V00);
815  V00.setSize(nCells());
816 
817  forAll(V00, i)
818  {
819  if (cellMap[i] > -1)
820  {
821  V00[i] = savedV00[cellMap[i]];
822  }
823  else
824  {
825  V00[i] = 0.0;
826  }
827  }
828 
829  // Inject volume of merged cells
830  label nMerged = 0;
831  forAll(meshMap.reverseCellMap(), oldCelli)
832  {
833  label index = meshMap.reverseCellMap()[oldCelli];
834 
835  if (index < -1)
836  {
837  label celli = -index-2;
838 
839  V00[celli] += savedV00[oldCelli];
840  nMerged++;
841  }
842  }
843 
844  DebugInfo
845  << "Mapping old time volume V00. Merged "
846  << nMerged << " out of " << nCells() << " cells" << endl;
847  }
848 }
849 
850 
852 {
854 
855  // Grab old time volumes if the time has been incremented
856  // This will update V0, V00
857  if (curTimeIndex_ < time().timeIndex())
858  {
859  storeOldVol(V());
860  }
861 
862 
863  // Move the polyMesh and set the mesh motion fluxes to the swept-volumes
864 
865  scalar rDeltaT = 1.0/time().deltaTValue();
866 
868  scalarField& sweptVols = tsweptVols.ref();
869 
870  if (!phiPtr_)
871  {
872  // Create mesh motion flux
873  phiPtr_ = new surfaceScalarField
874  (
875  IOobject
876  (
877  "meshPhi",
878  this->time().timeName(),
879  *this,
882  false
883  ),
884  *this,
886  );
887  }
888  else
889  {
890  // Grab old time mesh motion fluxes if the time has been incremented
891  if (phiPtr_->timeIndex() != time().timeIndex())
892  {
893  phiPtr_->oldTime();
894  }
895  }
896 
897  surfaceScalarField& phi = *phiPtr_;
899  scalarField::subField(sweptVols, nInternalFaces());
900  phi.primitiveFieldRef() *= rDeltaT;
901 
902  const fvPatchList& patches = boundary();
903 
904  surfaceScalarField::Boundary& phibf = phi.boundaryFieldRef();
905  forAll(patches, patchi)
906  {
907  phibf[patchi] = patches[patchi].patchSlice(sweptVols);
908  phibf[patchi] *= rDeltaT;
909  }
910 
911  // Update or delete the local geometric properties as early as possible so
912  // they can be used if necessary. These get recreated here instead of
913  // demand driven since they might do parallel transfers which can conflict
914  // with when they're actually being used.
915  // Note that between above "polyMesh::movePoints(p)" and here nothing
916  // should use the local geometric properties.
917  updateGeomNotOldVol();
918 
919 
920  // Update other local data
921  boundary_.movePoints();
923 
924  meshObject::movePoints<fvMesh>(*this);
925  meshObject::movePoints<lduMesh>(*this);
926 
927  return tsweptVols;
928 }
929 
930 
932 {
933  // Let surfaceInterpolation handle geometry calculation. Note: this does
934  // lower levels updateGeom
936 }
937 
938 
940 {
942 
943  // Update polyMesh. This needs to keep volume existent!
945 
946  // Our slice of the addressing is no longer valid
947  deleteDemandDrivenData(lduPtr_);
948 
949  if (VPtr_)
950  {
951  // Grab old time volumes if the time has been incremented
952  // This will update V0, V00
953  storeOldVol(mpm.oldCellVolumes());
954 
955  // Few checks
956  if (VPtr_ && (V().size() != mpm.nOldCells()))
957  {
959  << "V:" << V().size()
960  << " not equal to the number of old cells "
961  << mpm.nOldCells()
962  << exit(FatalError);
963  }
964  if (V0Ptr_ && (V0Ptr_->size() != mpm.nOldCells()))
965  {
967  << "V0:" << V0Ptr_->size()
968  << " not equal to the number of old cells "
969  << mpm.nOldCells()
970  << exit(FatalError);
971  }
972  if (V00Ptr_ && (V00Ptr_->size() != mpm.nOldCells()))
973  {
975  << "V0:" << V00Ptr_->size()
976  << " not equal to the number of old cells "
977  << mpm.nOldCells()
978  << exit(FatalError);
979  }
980  }
981 
982 
983  // Clear mesh motion flux (note: could instead save & map like volumes)
984  deleteDemandDrivenData(phiPtr_);
985 
986  // Clear the sliced fields
987  clearGeomNotOldVol();
988 
989  // Map all fields
990  mapFields(mpm);
991 
992  // Clear the current volume and other geometry factors
994 
995  // Clear any non-updateable addressing
996  clearAddressing(true);
997 
998  meshObject::updateMesh<fvMesh>(*this, mpm);
999  meshObject::updateMesh<lduMesh>(*this, mpm);
1000 }
1001 
1002 
1005  IOstreamOption streamOpt,
1006  const bool valid
1007 ) const
1008 {
1009  bool ok = true;
1010  if (phiPtr_)
1011  {
1012  ok = phiPtr_->write(valid);
1013  // NOTE: The old old time mesh phi might be necessary for certain
1014  // solver smooth restart using second order time schemes.
1015  //ok = phiPtr_->oldTime().write();
1016  }
1017  if (V0Ptr_ && V0Ptr_->writeOpt() == IOobject::AUTO_WRITE)
1018  {
1019  // For second order restarts we need to write V0
1020  ok = V0Ptr_->write(valid);
1021  }
1022 
1023  return ok && polyMesh::writeObject(streamOpt, valid);
1024 }
1025 
1026 
1027 bool Foam::fvMesh::write(const bool valid) const
1028 {
1029  return polyMesh::write(valid);
1030 }
1031 
1032 
1033 template<>
1035 Foam::fvMesh::validComponents<Foam::sphericalTensor>() const
1036 {
1038 }
1039 
1040 
1041 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1042 
1043 bool Foam::fvMesh::operator!=(const fvMesh& rhs) const
1044 {
1045  return &rhs != this;
1046 }
1047 
1048 
1049 bool Foam::fvMesh::operator==(const fvMesh& rhs) const
1050 {
1051  return &rhs == this;
1052 }
1053 
1054 
1055 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::lduAddressing
The class contains the addressing required by the lduMatrix: upper, lower and losort.
Definition: lduAddressing.H:114
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::surfaceInterpolation::clearOut
void clearOut()
Clear all geometry and addressing.
Definition: surfaceInterpolation.C:49
Foam::fvMesh::storeOldVol
void storeOldVol(const scalarField &)
Preserve old volume(s)
Definition: fvMesh.C:164
volFields.H
Foam::fvMesh::~fvMesh
virtual ~fvMesh()
Destructor.
Definition: fvMesh.C:536
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:325
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
SubField.H
Foam::fvMesh::write
virtual bool write(const bool valid=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1027
Foam::polyMesh::POINTS_MOVED
Definition: polyMesh.H:93
Foam::primitiveMesh::clearAddressing
void clearAddressing()
Clear topological data.
Definition: primitiveMeshClear.C:157
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::fvMesh::updateGeomNotOldVol
void updateGeomNotOldVol()
Clear geometry like clearGeomNotOldVol but recreate any.
Definition: fvMesh.C:82
Foam::polyMesh::movePoints
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: polyMesh.C:1160
Foam::fvMatrix::solveSegregatedOrCoupled
SolverPerformance< Type > solveSegregatedOrCoupled(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:62
Foam::fvSchemes
Selector class for finite volume differencing schemes. fvMesh is derived from fvShemes so that all fi...
Definition: fvSchemes.H:53
slicedVolFields.H
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::GeometricMeshObject
Definition: MeshObject.H:201
mapPolyMesh.H
Foam::fvMesh::VPtr_
void * VPtr_
Cell volumes old time level.
Definition: fvMesh.H:110
Foam::fvMesh::mapFields
virtual void mapFields(const mapPolyMesh &mpm)
Map all fields in time using given map.
Definition: fvMesh.C:700
Foam::fvMesh::clearGeomNotOldVol
void clearGeomNotOldVol()
Clear geometry but not the old-time cell volumes.
Definition: fvMesh.C:54
Foam::fvMesh::clearGeom
void clearGeom()
Clear geometry.
Definition: fvMesh.C:116
Foam::mapPolyMesh::nOldFaces
label nOldFaces() const
Number of old faces.
Definition: mapPolyMesh.H:381
Foam::MoveableMeshObject
Definition: MeshObject.H:220
Foam::fvMesh::removeFvBoundary
void removeFvBoundary()
Definition: fvMesh.C:630
Foam::isFile
bool isFile(const fileName &name, const bool checkGzip=true, const bool followLink=true)
Does the name exist as a FILE in the file system?
Definition: MSwindows.C:658
Foam::fileHandler
const fileOperation & fileHandler()
Get current file handler.
Definition: fileOperation.C:1354
fvMatrix.H
Foam::polyMesh::clearOut
void clearOut()
Clear all geometry and addressing unnecessary for CFD.
Definition: polyMeshClear.C:222
Foam::polyMesh::init
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
Definition: polyMesh.C:367
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
surfaceFields.H
Foam::surfaceFields.
Foam::fvMesh::operator!=
bool operator!=(const fvMesh &rhs) const
Compares addresses.
Definition: fvMesh.C:1043
Foam::fvMesh::operator==
bool operator==(const fvMesh &rhs) const
Compares addresses.
Definition: fvMesh.C:1049
Foam::UpdateableMeshObject
Definition: MeshObject.H:241
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
C
volScalarField & C
Definition: readThermalProperties.H:102
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
mapClouds.H
Generic Geometric field mapper. For "real" mapping, add template specialisations for mapping of inter...
Foam::fvMeshLduAddressing
Foam::fvMeshLduAddressing.
Definition: fvMeshLduAddressing.H:51
Foam::fvMesh::magSfPtr_
surfaceScalarField * magSfPtr_
Mag face area vectors.
Definition: fvMesh.H:122
Foam::fvBoundaryMesh
Foam::fvBoundaryMesh.
Definition: fvBoundaryMesh.H:56
Foam::fvMesh::SfPtr_
slicedSurfaceVectorField * SfPtr_
Face area vectors.
Definition: fvMesh.H:119
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::fvMesh::readUpdate
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:643
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::polyMesh::TOPO_PATCH_CHANGE
Definition: polyMesh.H:95
Foam::fvMesh::solve
virtual SolverPerformance< scalar > solve(fvMatrix< scalar > &, const dictionary &) const
Definition: fvMesh.C:545
Foam::fvMesh::CfPtr_
slicedSurfaceVectorField * CfPtr_
Face centres.
Definition: fvMesh.H:128
Foam::fvMesh::updateGeom
virtual void updateGeom()
Definition: fvMesh.C:931
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
fvMeshMapper.H
Foam::Field< scalar >
Foam::regIOobject::write
virtual bool write(const bool valid=true) const
Write using setting from DB.
Definition: regIOobjectWrite.C:165
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:365
Foam::fvMesh::lduAddr
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:685
Foam::polyMesh::updateMesh
virtual void updateMesh(const mapPolyMesh &mpm)
Update the mesh corresponding to given map.
Definition: polyMeshUpdate.C:42
Foam::fvSolution
Selector class for finite volume solution solution. fvMesh is derived from fvSolution so that all fie...
Definition: fvSolution.H:50
Foam::IOstreamOption
The IOstreamOption is a simple container for options an IOstream can normally have.
Definition: IOstreamOption.H:63
Foam::fvMesh::updateMesh
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:939
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:62
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::polyMesh::TOPO_CHANGE
Definition: polyMesh.H:94
timeName
word timeName
Definition: getTimeIndex.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::objectRegistry::writeObject
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write the objects using stream options.
Definition: objectRegistry.C:475
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::fvMesh::addFvPatches
void addFvPatches(PtrList< polyPatch > &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:600
Foam::surfaceInterpolation::movePoints
virtual bool movePoints()
Do what is necessary if the mesh has moved.
Definition: surfaceInterpolation.C:151
MapFvFields.H
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::fvMesh::fvMesh
fvMesh(const fvMesh &)=delete
No copy construct.
Foam::mapPolyMesh::oldCellVolumes
const scalarField & oldCellVolumes() const
Definition: mapPolyMesh.H:651
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::polyMesh::removeBoundary
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:39
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:766
Foam::polyMesh::readUpdateState
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:90
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:679
Foam::fvMesh::CPtr_
slicedVolVectorField * CPtr_
Cell centres.
Definition: fvMesh.H:125
Foam::mapPolyMesh::reverseCellMap
const labelList & reverseCellMap() const
Reverse cell map.
Definition: mapPolyMesh.H:532
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::fvMesh::init
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
Definition: fvMesh.C:271
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:359
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::fvMeshMapper
Class holds all the necessary information for mapping fields associated with fvMesh.
Definition: fvMeshMapper.H:57
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::List< face >
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:54
Foam::mapClouds
void mapClouds(const objectRegistry &db, const mapPolyMesh &mapper)
Generic Geometric field mapper.
Definition: mapClouds.H:51
Foam::mapPolyMesh::nOldCells
label nOldCells() const
Number of old cells.
Definition: mapPolyMesh.H:387
Foam::meshObject::clearUpto
static void clearUpto(objectRegistry &obr)
Definition: MeshObject.C:217
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::surfaceInterpolation::updateGeom
virtual void updateGeom()
Update all geometric data.
Definition: surfaceInterpolation.C:172
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::mapPolyMesh::faceMap
const labelList & faceMap() const
Old face map.
Definition: mapPolyMesh.H:410
Foam::surfaceInterpolation
Cell to surface interpolation scheme. Included in fvMesh.
Definition: surfaceInterpolation.H:58
Foam::SlicedGeometricField::Internal
The internalField of a SlicedGeometricField.
Definition: SlicedGeometricField.H:196
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
Foam::TopologicalMeshObject
Definition: MeshObject.H:182
slicedSurfaceFields.H
MeshObject.H
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:62
Foam::Field< scalar >::subField
SubField< scalar > subField
Declare type of subField.
Definition: Field.H:89
timeIndex
label timeIndex
Definition: getTimeIndex.H:4
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::SolverPerformance
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
Definition: SolverPerformance.H:51
fvMeshLduAddressing.H
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::data
Database for solution data, solver performance and other reduced data.
Definition: data.H:55
boundary
faceListList boundary
Definition: createBlockMesh.H:4
Foam::mapPolyMesh::cellMap
const labelList & cellMap() const
Old cell map.
Definition: mapPolyMesh.H:435
Foam::lduMesh
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:62
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::fvMesh::clearOut
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:227
Foam::polyMesh::readUpdate
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in.
Definition: polyMeshIO.C:75
Foam::zero
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:62
Foam::IOobject::MUST_READ
Definition: IOobject.H:120
Foam::fvMesh::writeObject
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write the underlying polyMesh and other data.
Definition: fvMesh.C:1004