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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "fvMesh.H"
27 #include "volFields.H"
28 #include "surfaceFields.H"
29 #include "slicedVolFields.H"
30 #include "slicedSurfaceFields.H"
31 #include "SubField.H"
32 #include "demandDrivenData.H"
33 #include "fvMeshLduAddressing.H"
34 #include "mapPolyMesh.H"
35 #include "MapFvFields.H"
36 #include "fvMeshMapper.H"
37 #include "mapClouds.H"
38 #include "MeshObject.H"
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44  defineTypeNameAndDebug(fvMesh, 0);
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 void Foam::fvMesh::clearGeomNotOldVol()
51 {
53  <
54  fvMesh,
55  GeometricMeshObject,
56  MoveableMeshObject
57  >(*this);
58 
60  <
61  lduMesh,
62  GeometricMeshObject,
63  MoveableMeshObject
64  >(*this);
65 
69  VPtr_ = NULL;
70 
71  deleteDemandDrivenData(SfPtr_);
72  deleteDemandDrivenData(magSfPtr_);
74  deleteDemandDrivenData(CfPtr_);
75 }
76 
77 
78 void Foam::fvMesh::updateGeomNotOldVol()
79 {
80  bool haveV = (VPtr_ != NULL);
81  bool haveSf = (SfPtr_ != NULL);
82  bool haveMagSf = (magSfPtr_ != NULL);
83  bool haveCP = (CPtr_ != NULL);
84  bool haveCf = (CfPtr_ != NULL);
85 
86  clearGeomNotOldVol();
87 
88  // Now recreate the fields
89  if (haveV)
90  {
91  (void)V();
92  }
93  if (haveSf)
94  {
95  (void)Sf();
96  }
97  if (haveMagSf)
98  {
99  (void)magSf();
100  }
101  if (haveCP)
102  {
103  (void)C();
104  }
105  if (haveCf)
106  {
107  (void)Cf();
108  }
109 }
110 
111 
112 void Foam::fvMesh::clearGeom()
113 {
114  clearGeomNotOldVol();
115 
116  deleteDemandDrivenData(V0Ptr_);
117  deleteDemandDrivenData(V00Ptr_);
118 
119  // Mesh motion flux cannot be deleted here because the old-time flux
120  // needs to be saved.
121 }
122 
123 
124 void Foam::fvMesh::clearAddressing(const bool isMeshUpdate)
125 {
126  if (debug)
127  {
128  Info<< "fvMesh::clearAddressing(const bool) :"
129  << " isMeshUpdate:" << isMeshUpdate << endl;
130  }
131 
132  if (isMeshUpdate)
133  {
134  // Part of a mesh update. Keep meshObjects that have an updateMesh
135  // callback
137  <
138  fvMesh,
139  TopologicalMeshObject,
140  UpdateableMeshObject
141  >
142  (
143  *this
144  );
146  <
147  lduMesh,
148  TopologicalMeshObject,
149  UpdateableMeshObject
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 
164 void Foam::fvMesh::storeOldVol(const scalarField& V)
165 {
166  if (curTimeIndex_ < time().timeIndex())
167  {
168  if (debug)
169  {
170  Info<< "fvMesh::storeOldVol(const scalarField&) :"
171  << " Storing old time volumes since from time " << curTimeIndex_
172  << " and time now " << time().timeIndex()
173  << " V:" << V.size()
174  << endl;
175  }
176 
177 
178  if (V00Ptr_ && V0Ptr_)
179  {
180  // Copy V0 into V00 storage
181  *V00Ptr_ = *V0Ptr_;
182  }
183 
184  if (V0Ptr_)
185  {
186  // Copy V into V0 storage
187  V0Ptr_->scalarField::operator=(V);
188  }
189  else
190  {
191  // Allocate V0 storage, fill with V
192  V0Ptr_ = new DimensionedField<scalar, volMesh>
193  (
194  IOobject
195  (
196  "V0",
197  time().timeName(),
198  *this,
201  false
202  ),
203  *this,
204  dimVolume
205  );
206  scalarField& V0 = *V0Ptr_;
207  // Note: V0 now sized with current mesh, not with (potentially
208  // different size) V.
209  V0.setSize(V.size());
210  V0 = V;
211  }
212 
213  curTimeIndex_ = time().timeIndex();
214 
215  if (debug)
216  {
217  Info<< "fvMesh::storeOldVol() :"
218  << " Stored old time volumes V0:" << V0Ptr_->size()
219  << endl;
220  if (V00Ptr_)
221  {
222  Info<< "fvMesh::storeOldVol() :"
223  << " Stored oldold time volumes V00:" << V00Ptr_->size()
224  << endl;
225  }
226  }
227  }
228 }
229 
230 
232 {
233  clearGeom();
235 
236  clearAddressing();
237 
238  // Clear mesh motion flux
239  deleteDemandDrivenData(phiPtr_);
240 
242 }
243 
244 
245 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
246 
247 Foam::fvMesh::fvMesh(const IOobject& io)
248 :
249  polyMesh(io),
250  surfaceInterpolation(*this),
251  fvSchemes(static_cast<const objectRegistry&>(*this)),
252  fvSolution(static_cast<const objectRegistry&>(*this)),
253  data(static_cast<const objectRegistry&>(*this)),
254  boundary_(*this, boundaryMesh()),
255  lduPtr_(NULL),
256  curTimeIndex_(time().timeIndex()),
257  VPtr_(NULL),
258  V0Ptr_(NULL),
259  V00Ptr_(NULL),
260  SfPtr_(NULL),
261  magSfPtr_(NULL),
262  CPtr_(NULL),
263  CfPtr_(NULL),
264  phiPtr_(NULL)
265 {
266  if (debug)
267  {
268  Info<< "Constructing fvMesh from IOobject"
269  << endl;
270  }
271 
272  // Check the existance of the cell volumes and read if present
273  // and set the storage of V00
274  if (isFile(time().timePath()/"V0"))
275  {
277  (
278  IOobject
279  (
280  "V0",
281  time().timeName(),
282  *this,
285  false
286  ),
287  *this
288  );
289 
290  V00();
291  }
292 
293  // Check the existance of the mesh fluxes, read if present and set the
294  // mesh to be moving
295  if (isFile(time().timePath()/"meshPhi"))
296  {
297  phiPtr_ = new surfaceScalarField
298  (
299  IOobject
300  (
301  "meshPhi",
302  time().timeName(),
303  *this,
306  false
307  ),
308  *this
309  );
310 
311  // The mesh is now considered moving so the old-time cell volumes
312  // will be required for the time derivatives so if they haven't been
313  // read initialise to the current cell volumes
314  if (!V0Ptr_)
315  {
317  (
318  IOobject
319  (
320  "V0",
321  time().timeName(),
322  *this,
325  false
326  ),
327  V()
328  );
329  }
330 
331  moving(true);
332  }
333 }
334 
335 
336 Foam::fvMesh::fvMesh
337 (
338  const IOobject& io,
339  const Xfer<pointField>& points,
340  const cellShapeList& shapes,
341  const faceListList& boundaryFaces,
342  const wordList& boundaryPatchNames,
343  const PtrList<dictionary>& boundaryDicts,
344  const word& defaultBoundaryPatchName,
345  const word& defaultBoundaryPatchType,
346  const bool syncPar
347 )
348 :
349  polyMesh
350  (
351  io,
352  points,
353  shapes,
354  boundaryFaces,
355  boundaryPatchNames,
356  boundaryDicts,
357  defaultBoundaryPatchName,
358  defaultBoundaryPatchType,
359  syncPar
360  ),
361  surfaceInterpolation(*this),
362  fvSchemes(static_cast<const objectRegistry&>(*this)),
363  fvSolution(static_cast<const objectRegistry&>(*this)),
364  data(static_cast<const objectRegistry&>(*this)),
365  boundary_(*this, boundaryMesh()),
366  lduPtr_(NULL),
367  curTimeIndex_(time().timeIndex()),
368  VPtr_(NULL),
369  V0Ptr_(NULL),
370  V00Ptr_(NULL),
371  SfPtr_(NULL),
372  magSfPtr_(NULL),
373  CPtr_(NULL),
374  CfPtr_(NULL),
375  phiPtr_(NULL)
376 {
377  if (debug)
378  {
379  Info<< "Constructing fvMesh from cellShapes" << endl;
380  }
381 }
382 
383 
384 Foam::fvMesh::fvMesh
385 (
386  const IOobject& io,
387  const Xfer<pointField>& points,
388  const Xfer<faceList>& faces,
389  const Xfer<labelList>& allOwner,
390  const Xfer<labelList>& allNeighbour,
391  const bool syncPar
392 )
393 :
394  polyMesh(io, points, faces, allOwner, allNeighbour, syncPar),
395  surfaceInterpolation(*this),
396  fvSchemes(static_cast<const objectRegistry&>(*this)),
397  fvSolution(static_cast<const objectRegistry&>(*this)),
398  data(static_cast<const objectRegistry&>(*this)),
399  boundary_(*this, boundaryMesh()),
400  lduPtr_(NULL),
401  curTimeIndex_(time().timeIndex()),
402  VPtr_(NULL),
403  V0Ptr_(NULL),
404  V00Ptr_(NULL),
405  SfPtr_(NULL),
406  magSfPtr_(NULL),
407  CPtr_(NULL),
408  CfPtr_(NULL),
409  phiPtr_(NULL)
410 {
411  if (debug)
412  {
413  Info<< "Constructing fvMesh from components" << endl;
414  }
415 }
416 
417 
418 Foam::fvMesh::fvMesh
419 (
420  const IOobject& io,
421  const Xfer<pointField>& points,
422  const Xfer<faceList>& faces,
423  const Xfer<cellList>& cells,
424  const bool syncPar
425 )
426 :
427  polyMesh(io, points, faces, cells, syncPar),
428  surfaceInterpolation(*this),
429  fvSchemes(static_cast<const objectRegistry&>(*this)),
430  fvSolution(static_cast<const objectRegistry&>(*this)),
431  data(static_cast<const objectRegistry&>(*this)),
432  boundary_(*this),
433  lduPtr_(NULL),
434  curTimeIndex_(time().timeIndex()),
435  VPtr_(NULL),
436  V0Ptr_(NULL),
437  V00Ptr_(NULL),
438  SfPtr_(NULL),
439  magSfPtr_(NULL),
440  CPtr_(NULL),
441  CfPtr_(NULL),
442  phiPtr_(NULL)
443 {
444  if (debug)
445  {
446  Info<< "Constructing fvMesh from components" << endl;
447  }
448 }
449 
450 
451 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
452 
454 {
455  clearOut();
456 }
457 
458 
459 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
460 
462 (
463  const List<polyPatch*> & p,
464  const bool validBoundary
465 )
466 {
467  if (boundary().size())
468  {
470  (
471  "fvMesh::addFvPatches(const List<polyPatch*>&, const bool)"
472  ) << " boundary already exists"
473  << abort(FatalError);
474  }
475 
476  // first add polyPatches
477  addPatches(p, validBoundary);
478  boundary_.addPatches(boundaryMesh());
479 }
480 
481 
483 {
484  if (debug)
485  {
486  Info<< "void fvMesh::removeFvBoundary(): "
487  << "Removing boundary patches."
488  << endl;
489  }
490 
491  // Remove fvBoundaryMesh data first.
492  boundary_.clear();
493  boundary_.setSize(0);
495 
496  clearOut();
497 }
498 
499 
501 {
502  if (debug)
503  {
504  Info<< "polyMesh::readUpdateState fvMesh::readUpdate() : "
505  << "Updating fvMesh. ";
506  }
507 
509 
510  if (state == polyMesh::TOPO_PATCH_CHANGE)
511  {
512  if (debug)
513  {
514  Info<< "Boundary and topological update" << endl;
515  }
516 
517  boundary_.readUpdate(boundaryMesh());
518 
519  clearOut();
520 
521  }
522  else if (state == polyMesh::TOPO_CHANGE)
523  {
524  if (debug)
525  {
526  Info<< "Topological update" << endl;
527  }
528 
529  clearOut();
530  }
531  else if (state == polyMesh::POINTS_MOVED)
532  {
533  if (debug)
534  {
535  Info<< "Point motion update" << endl;
536  }
537 
538  clearGeom();
539  }
540  else
541  {
542  if (debug)
543  {
544  Info<< "No update" << endl;
545  }
546  }
547 
548  return state;
549 }
550 
551 
553 {
554  return boundary_;
555 }
556 
557 
559 {
560  if (!lduPtr_)
561  {
562  lduPtr_ = new fvMeshLduAddressing(*this);
563  }
564 
565  return *lduPtr_;
566 }
567 
568 
570 {
571  if (debug)
572  {
573  Info<< "fvMesh::mapFields :"
574  << " nOldCells:" << meshMap.nOldCells()
575  << " nCells:" << nCells()
576  << " nOldFaces:" << meshMap.nOldFaces()
577  << " nFaces:" << nFaces()
578  << endl;
579  }
580 
581 
582  // We require geometric properties valid for the old mesh
583  if
584  (
585  meshMap.cellMap().size() != nCells()
586  || meshMap.faceMap().size() != nFaces()
587  )
588  {
589  FatalErrorIn("fvMesh::mapFields(const mapPolyMesh&)")
590  << "mapPolyMesh does not correspond to the old mesh."
591  << " nCells:" << nCells()
592  << " cellMap:" << meshMap.cellMap().size()
593  << " nOldCells:" << meshMap.nOldCells()
594  << " nFaces:" << nFaces()
595  << " faceMap:" << meshMap.faceMap().size()
596  << " nOldFaces:" << meshMap.nOldFaces()
597  << exit(FatalError);
598  }
599 
600  // Create a mapper
601  const fvMeshMapper mapper(*this, meshMap);
602 
603  // Map all the volFields in the objectRegistry
604  MapGeometricFields<scalar, fvPatchField, fvMeshMapper, volMesh>
605  (mapper);
606  MapGeometricFields<vector, fvPatchField, fvMeshMapper, volMesh>
607  (mapper);
608  MapGeometricFields<sphericalTensor, fvPatchField, fvMeshMapper, volMesh>
609  (mapper);
610  MapGeometricFields<symmTensor, fvPatchField, fvMeshMapper, volMesh>
611  (mapper);
612  MapGeometricFields<tensor, fvPatchField, fvMeshMapper, volMesh>
613  (mapper);
614 
615  // Map all the surfaceFields in the objectRegistry
616  MapGeometricFields<scalar, fvsPatchField, fvMeshMapper, surfaceMesh>
617  (mapper);
618  MapGeometricFields<vector, fvsPatchField, fvMeshMapper, surfaceMesh>
619  (mapper);
620  MapGeometricFields<symmTensor, fvsPatchField, fvMeshMapper, surfaceMesh>
621  (mapper);
622  MapGeometricFields<symmTensor, fvsPatchField, fvMeshMapper, surfaceMesh>
623  (mapper);
624  MapGeometricFields<tensor, fvsPatchField, fvMeshMapper, surfaceMesh>
625  (mapper);
626 
627  // Map all the dimensionedFields in the objectRegistry
628  MapDimensionedFields<scalar, fvMeshMapper, volMesh>(mapper);
629  MapDimensionedFields<vector, fvMeshMapper, volMesh>(mapper);
630  MapDimensionedFields<sphericalTensor, fvMeshMapper, volMesh>(mapper);
631  MapDimensionedFields<symmTensor, fvMeshMapper, volMesh>(mapper);
632  MapDimensionedFields<tensor, fvMeshMapper, volMesh>(mapper);
633 
634  // Map all the clouds in the objectRegistry
635  mapClouds(*this, meshMap);
636 
637 
638  const labelList& cellMap = meshMap.cellMap();
639 
640  // Map the old volume. Just map to new cell labels.
641  if (V0Ptr_)
642  {
643  scalarField& V0 = *V0Ptr_;
644 
645  scalarField savedV0(V0);
646  V0.setSize(nCells());
647 
648  forAll(V0, i)
649  {
650  if (cellMap[i] > -1)
651  {
652  V0[i] = savedV0[cellMap[i]];
653  }
654  else
655  {
656  V0[i] = 0.0;
657  }
658  }
659 
660  // Inject volume of merged cells
661  label nMerged = 0;
662  forAll(meshMap.reverseCellMap(), oldCellI)
663  {
664  label index = meshMap.reverseCellMap()[oldCellI];
665 
666  if (index < -1)
667  {
668  label cellI = -index-2;
669 
670  V0[cellI] += savedV0[oldCellI];
671 
672  nMerged++;
673  }
674  }
675 
676  if (debug)
677  {
678  Info<< "Mapping old time volume V0. Merged "
679  << nMerged << " out of " << nCells() << " cells" << endl;
680  }
681  }
682 
683 
684  // Map the old-old volume. Just map to new cell labels.
685  if (V00Ptr_)
686  {
687  scalarField& V00 = *V00Ptr_;
688 
689  scalarField savedV00(V00);
690  V00.setSize(nCells());
691 
692  forAll(V00, i)
693  {
694  if (cellMap[i] > -1)
695  {
696  V00[i] = savedV00[cellMap[i]];
697  }
698  else
699  {
700  V00[i] = 0.0;
701  }
702  }
703 
704  // Inject volume of merged cells
705  label nMerged = 0;
706  forAll(meshMap.reverseCellMap(), oldCellI)
707  {
708  label index = meshMap.reverseCellMap()[oldCellI];
709 
710  if (index < -1)
711  {
712  label cellI = -index-2;
713 
714  V00[cellI] += savedV00[oldCellI];
715  nMerged++;
716  }
717  }
718 
719  if (debug)
720  {
721  Info<< "Mapping old time volume V00. Merged "
722  << nMerged << " out of " << nCells() << " cells" << endl;
723  }
724  }
725 }
726 
727 
729 {
730  // Grab old time volumes if the time has been incremented
731  // This will update V0, V00
732  if (curTimeIndex_ < time().timeIndex())
733  {
734  storeOldVol(V());
735  }
736 
737  if (!phiPtr_)
738  {
739  // Create mesh motion flux
740  phiPtr_ = new surfaceScalarField
741  (
742  IOobject
743  (
744  "meshPhi",
745  this->time().timeName(),
746  *this,
749  false
750  ),
751  *this,
753  );
754  }
755  else
756  {
757  // Grab old time mesh motion fluxes if the time has been incremented
758  if (phiPtr_->timeIndex() != time().timeIndex())
759  {
760  phiPtr_->oldTime();
761  }
762  }
763 
764  surfaceScalarField& phi = *phiPtr_;
765 
766  // Move the polyMesh and set the mesh motion fluxes to the swept-volumes
767 
768  scalar rDeltaT = 1.0/time().deltaTValue();
769 
770  tmp<scalarField> tsweptVols = polyMesh::movePoints(p);
771  scalarField& sweptVols = tsweptVols();
772 
773  phi.internalField() = scalarField::subField(sweptVols, nInternalFaces());
774  phi.internalField() *= rDeltaT;
775 
776  const fvPatchList& patches = boundary();
777 
778  forAll(patches, patchI)
779  {
780  phi.boundaryField()[patchI] = patches[patchI].patchSlice(sweptVols);
781  phi.boundaryField()[patchI] *= rDeltaT;
782  }
783 
784  // Update or delete the local geometric properties as early as possible so
785  // they can be used if necessary. These get recreated here instead of
786  // demand driven since they might do parallel transfers which can conflict
787  // with when they're actually being used.
788  // Note that between above "polyMesh::movePoints(p)" and here nothing
789  // should use the local geometric properties.
790  updateGeomNotOldVol();
791 
792 
793  // Update other local data
794  boundary_.movePoints();
796 
797  meshObject::movePoints<fvMesh>(*this);
798  meshObject::movePoints<lduMesh>(*this);
799 
800  return tsweptVols;
801 }
802 
803 
805 {
806  // Update polyMesh. This needs to keep volume existent!
808 
809  if (VPtr_)
810  {
811  // Grab old time volumes if the time has been incremented
812  // This will update V0, V00
813  storeOldVol(mpm.oldCellVolumes());
814 
815  // Few checks
816  if (VPtr_ && (V().size() != mpm.nOldCells()))
817  {
818  FatalErrorIn("fvMesh::updateMesh(const mapPolyMesh&)")
819  << "V:" << V().size()
820  << " not equal to the number of old cells "
821  << mpm.nOldCells()
822  << exit(FatalError);
823  }
824  if (V0Ptr_ && (V0Ptr_->size() != mpm.nOldCells()))
825  {
826  FatalErrorIn("fvMesh::updateMesh(const mapPolyMesh&)")
827  << "V0:" << V0Ptr_->size()
828  << " not equal to the number of old cells "
829  << mpm.nOldCells()
830  << exit(FatalError);
831  }
832  if (V00Ptr_ && (V00Ptr_->size() != mpm.nOldCells()))
833  {
834  FatalErrorIn("fvMesh::updateMesh(const mapPolyMesh&)")
835  << "V0:" << V00Ptr_->size()
836  << " not equal to the number of old cells "
837  << mpm.nOldCells()
838  << exit(FatalError);
839  }
840  }
841 
842 
843  // Clear mesh motion flux (note: could instead save & map like volumes)
844  deleteDemandDrivenData(phiPtr_);
845 
846  // Clear the sliced fields
847  clearGeomNotOldVol();
848 
849  // Map all fields
850  mapFields(mpm);
851 
852  // Clear the current volume and other geometry factors
854 
855  // Clear any non-updateable addressing
856  clearAddressing(true);
857 
858  meshObject::updateMesh<fvMesh>(*this, mpm);
859  meshObject::updateMesh<lduMesh>(*this, mpm);
860 }
861 
862 
864 (
868 ) const
869 {
870  return polyMesh::writeObject(fmt, ver, cmp);
871 }
872 
873 
874 //- Write mesh using IO settings from the time
876 {
877  bool ok = true;
878  if (phiPtr_)
879  {
880  ok = phiPtr_->write();
881  }
882 
883  return ok && polyMesh::write();
884 }
885 
886 
887 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
888 
889 bool Foam::fvMesh::operator!=(const fvMesh& bm) const
890 {
891  return &bm != this;
892 }
893 
894 
895 bool Foam::fvMesh::operator==(const fvMesh& bm) const
896 {
897  return &bm == this;
898 }
899 
900 
901 // ************************************************************************* //
virtual bool write() const
Write mesh using IO settings from time.
Definition: fvMesh.C:875
virtual void updateMesh(const mapPolyMesh &mpm)
Update the mesh corresponding to given map.
void addFvPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:462
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: polyMesh.C:1105
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:142
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::surfaceFields.
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp) const
Write the objects.
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:36
label nFaces() const
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual bool writeObjects(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp) const
Write the underlying polyMesh and other data.
Definition: fvMesh.C:864
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in.
Definition: polyMeshIO.C:66
void clearAddressing()
Clear topological data.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
void clear()
Clear the PtrList, i.e. set size to zero deleting all the.
Definition: PtrList.C:185
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const surfaceScalarField & phi() const
Return cell face motion fluxes.
A simple container for copying or transferring objects of type <T>.
Definition: Xfer.H:85
Foam::fvBoundaryMesh.
void deleteDemandDrivenData(DataPtr &dataPtr)
label timeIndex
Definition: getTimeIndex.H:4
A class for handling words, derived from string.
Definition: word.H:59
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:59
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
InternalField & internalField()
Return internal field.
label nCells() const
Database for solution data, solver performance and other reduced data.
Definition: data.H:52
void removeFvBoundary()
Remove boundary patches. Warning: fvPatchFields hold ref to.
Definition: fvMesh.C:482
messageStream Info
Selector class for finite volume differencing schemes. fvMesh is derived from fvShemes so that all fi...
Definition: fvSchemes.H:50
void addPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:925
static void clearUpto(objectRegistry &)
Clear all meshObject derived from FromType up to (but not including)
Definition: MeshObject.C:401
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:804
Namespace for OpenFOAM.
label nOldFaces() const
Number of old faces.
Definition: mapPolyMesh.H:377
const scalarField & oldCellVolumes() const
Definition: mapPolyMesh.H:647
patches[0]
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
virtual bool write() const
Write using setting from DB.
void setSize(const label)
Reset size of List.
Definition: List.C:318
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
scalar deltaTValue() const
Return time step value.
Definition: TimeState.H:100
void movePoints()
Correct patches after moving points.
label timeIndex() const
Return current time index.
Definition: TimeState.C:73
streamFormat
Enumeration for the format of data in the stream.
Definition: IOstream.H:86
bool operator!=(const fvMesh &) const
Definition: fvMesh.C:889
Selector class for finite volume solution solution. fvMesh is derived from fvSolution so that all fie...
Definition: fvSolution.H:47
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:88
virtual void mapFields(const mapPolyMesh &mpm)
Map all fields in time using given map.
Definition: fvMesh.C:569
#define forAll(list, i)
Definition: UList.H:421
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const labelList & faceMap() const
Old face map.
Definition: mapPolyMesh.H:406
Template functions to aid in the implementation of demand driven data.
label size() const
Return number of elements in table.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:500
void clearOut()
Clear all geometry and addressing unnecessary for CFD.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
compressionType
Enumeration for the format of data in the stream.
Definition: IOstream.H:193
label nInternalFaces() const
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:558
void clearOut()
Clear all geometry and addressing.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
error FatalError
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:131
SubField< scalar > subField
Declare type of subField.
Definition: Field.H:89
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Registry of regIOobjects.
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:231
DimensionedField< Type, GeoMesh > DimensionedInternalField
Class holds all the necessary information for mapping fields associated with fvMesh.
Definition: fvMeshMapper.H:55
label timeIndex() const
Return the time index of the field.
const labelList & reverseCellMap() const
Reverse cell map.
Definition: mapPolyMesh.H:528
surfaceInterpolation(const fvMesh &)
Construct given an fvMesh.
The class contains the addressing required by the lduMatrix: upper, lower and losort.
Cell to surface interpolation scheme. Included in fvMesh.
bool isFile(const fileName &, const bool checkGzip=true)
Does the name exist as a FILE in the file system?
Definition: POSIX.C:623
Version number type.
Definition: IOstream.H:96
virtual ~fvMesh()
Destructor.
Definition: fvMesh.C:453
bool operator==(const fvMesh &) const
Definition: fvMesh.C:895
bool moving() const
Is mesh moving.
Definition: polyMesh.H:493
volScalarField & C
A class for managing temporary objects.
Definition: PtrList.H:118
void readUpdate(const polyBoundaryMesh &)
Update boundary based on new polyBoundaryMesh.
bool movePoints()
Do what is neccessary if the mesh has moved.
defineTypeNameAndDebug(combustionModel, 0)
word timeName
Definition: getTimeIndex.H:3
label nOldCells() const
Number of old cells.
Definition: mapPolyMesh.H:383
Generic Geometric field mapper. For "real" mapping, add template specialisations for mapping of inter...
void mapClouds(const objectRegistry &db, const mapPolyMesh &mapper)
Generic Geometric field mapper.
Definition: mapClouds.H:48
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:552
Foam::fvMeshLduAddressing.
const labelList & cellMap() const
Old cell map.
Definition: mapPolyMesh.H:431
static int debug
Debug switch.
Definition: fvSchemes.H:103