fvMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 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 "pointFields.H"
30 #include "slicedVolFields.H"
31 #include "slicedSurfaceFields.H"
32 #include "SubField.H"
33 #include "demandDrivenData.H"
34 #include "fvMeshLduAddressing.H"
35 #include "fvMeshTopoChanger.H"
36 #include "fvMeshDistributor.H"
37 #include "fvMeshMover.H"
38 #include "fvMeshStitcher.H"
39 #include "nonConformalFvPatch.H"
41 #include "polyTopoChangeMap.H"
42 #include "MapFvFields.H"
43 #include "fvMeshMapper.H"
44 #include "pointMesh.H"
45 #include "pointMeshMapper.H"
46 #include "MapPointField.H"
47 #include "meshObjects.H"
48 #include "HashPtrTable.H"
49 #include "CompactListList.H"
50 
51 #include "fvcSurfaceIntegrate.H"
52 #include "fvcReconstruct.H"
53 #include "surfaceInterpolate.H"
54 
55 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
56 
57 namespace Foam
58 {
60 }
61 
63 {
64  "Vc",
65  "Vc0",
66  "Vc00",
67  "Sf",
68  "magSf",
69  "Cc",
70  "Cf",
71  "meshPhi"
72 };
73 
74 
75 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
76 
77 void Foam::fvMesh::clearGeomNotOldVol()
78 {
79  if (debug)
80  {
81  Pout<< FUNCTION_NAME << "clearGeomNotOldVol" << endl;
82  }
83 
85  <
86  fvMesh,
87  DeletableMeshObject,
88  MoveableMeshObject
89  >(*this);
90 
92  <
93  lduMesh,
94  DeletableMeshObject,
95  MoveableMeshObject
96  >(*this);
97 
99  deleteDemandDrivenData(SfSlicePtr_);
100  deleteDemandDrivenData(SfPtr_);
101  deleteDemandDrivenData(magSfSlicePtr_);
102  deleteDemandDrivenData(magSfPtr_);
103  deleteDemandDrivenData(CSlicePtr_);
104  deleteDemandDrivenData(CPtr_);
105  deleteDemandDrivenData(CfSlicePtr_);
106  deleteDemandDrivenData(CfPtr_);
107 }
108 
109 
110 void Foam::fvMesh::updateGeomNotOldVol()
111 {
112  bool haveV = (VPtr_ != nullptr);
113  bool haveSf = (SfSlicePtr_ != nullptr || SfPtr_ != nullptr);
114  bool haveMagSf = (magSfSlicePtr_ != nullptr || magSfPtr_ != nullptr);
115  bool haveCP = (CSlicePtr_ != nullptr || CPtr_ != nullptr);
116  bool haveCf = (CfSlicePtr_ != nullptr || CfPtr_ != nullptr);
117 
118  clearGeomNotOldVol();
119 
120  // Now recreate the fields
121  if (haveV)
122  {
123  (void)V();
124  }
125  if (haveSf)
126  {
127  (void)Sf();
128  }
129  if (haveMagSf)
130  {
131  (void)magSf();
132  }
133  if (haveCP)
134  {
135  (void)C();
136  }
137  if (haveCf)
138  {
139  (void)Cf();
140  }
141 }
142 
143 
144 void Foam::fvMesh::clearGeom()
145 {
146  if (debug)
147  {
148  Pout<< FUNCTION_NAME << "Clearing geometric data" << endl;
149  }
150 
151  clearGeomNotOldVol();
152 
153  deleteDemandDrivenData(phiPtr_);
154  deleteDemandDrivenData(V0Ptr_);
155  deleteDemandDrivenData(V00Ptr_);
156 }
157 
158 
159 void Foam::fvMesh::clearAddressing(const bool isMeshUpdate)
160 {
161  if (debug)
162  {
163  Pout<< FUNCTION_NAME << "isMeshUpdate: " << isMeshUpdate << endl;
164  }
165 
166  if (isMeshUpdate)
167  {
168  // Part of a mesh update. Keep meshObjects that have an topoChange
169  // callback
171  <
172  fvMesh,
173  DeletableMeshObject,
174  TopoChangeableMeshObject
175  >
176  (
177  *this
178  );
180  <
181  lduMesh,
182  DeletableMeshObject,
183  TopoChangeableMeshObject
184  >
185  (
186  *this
187  );
188  }
189  else
190  {
191  meshObjects::clear<fvMesh, DeletableMeshObject>(*this);
192  meshObjects::clear<lduMesh, DeletableMeshObject>(*this);
193  }
194 
195  deleteDemandDrivenData(lduPtr_);
196  deleteDemandDrivenData(polyFacesBfIOPtr_);
197  deleteDemandDrivenData(polyFacesBfPtr_);
198  deleteDemandDrivenData(polyBFaceOffsetsPtr_);
199  deleteDemandDrivenData(polyBFaceOffsetPatchesPtr_);
200  deleteDemandDrivenData(polyBFaceOffsetPatchFacesPtr_);
201  deleteDemandDrivenData(polyBFacePatchesPtr_);
202  deleteDemandDrivenData(polyBFacePatchFacesPtr_);
203  deleteDemandDrivenData(ownerBfPtr_);
204 }
205 
206 
207 void Foam::fvMesh::storeOldTimeFields()
208 {
209  storeOldTimeFields<PointField>();
210  storeOldTimeFields<VolField>();
211  storeOldTimeFields<SurfaceField>();
212 }
213 
214 
215 void Foam::fvMesh::nullOldestTimeFields()
216 {
217  nullOldestTimeFields<PointField>();
218  nullOldestTimeFields<VolField>();
219  nullOldestTimeFields<SurfaceField>();
220 }
221 
222 
224 {
225  clearGeom();
226 
228 
229  clearAddressing();
230 
232 }
233 
234 
235 Foam::wordList Foam::fvMesh::polyFacesPatchTypes() const
236 {
237  wordList wantedPatchTypes
238  (
239  boundary().size(),
241  );
242 
244  {
245  const fvPatch& fvp = boundary()[patchi];
246 
247  if (isA<nonConformalFvPatch>(fvp))
248  {
249  wantedPatchTypes[patchi] =
250  refCast<const nonConformalFvPatch>(fvp).polyFacesType();
251  }
252  }
253 
254  return wantedPatchTypes;
255 }
256 
257 
258 Foam::surfaceLabelField::Boundary& Foam::fvMesh::polyFacesBfRef()
259 {
260  if (!polyFacesBfPtr_)
261  {
262  polyFacesBf();
263  }
264 
265  setPolyFacesBfInstance(time().name());
266 
267  return *polyFacesBfPtr_;
268 }
269 
270 
271 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
272 
273 Foam::fvMesh::fvMesh(const IOobject& io, const bool doPost)
274 :
275  polyMesh(io),
276  surfaceInterpolation(*this),
277  boundary_(*this, boundaryMesh()),
278  stitcher_(nullptr),
279  topoChanger_(nullptr),
280  distributor_(nullptr),
281  mover_(nullptr),
282  lduPtr_(nullptr),
283  polyFacesBfIOPtr_(nullptr),
284  polyFacesBfPtr_(nullptr),
285  polyBFaceOffsetsPtr_(nullptr),
286  polyBFaceOffsetPatchesPtr_(nullptr),
287  polyBFaceOffsetPatchFacesPtr_(nullptr),
288  polyBFacePatchesPtr_(nullptr),
289  polyBFacePatchFacesPtr_(nullptr),
290  ownerBfPtr_(nullptr),
291  curTimeIndex_(time().timeIndex()),
292  VPtr_(nullptr),
293  V0Ptr_(nullptr),
294  V00Ptr_(nullptr),
295  SfSlicePtr_(nullptr),
296  SfPtr_(nullptr),
297  magSfSlicePtr_(nullptr),
298  magSfPtr_(nullptr),
299  CSlicePtr_(nullptr),
300  CPtr_(nullptr),
301  CfSlicePtr_(nullptr),
302  CfPtr_(nullptr),
303  phiPtr_(nullptr)
304 {
305  if (debug)
306  {
307  Pout<< FUNCTION_NAME << "Constructing fvMesh from IOobject" << endl;
308  }
309 
310  if (doPost)
311  {
313  }
314 }
315 
316 
318 (
319  const IOobject& io,
320  pointField&& points,
321  const cellShapeList& shapes,
322  const faceListList& boundaryFaces,
323  const wordList& boundaryPatchNames,
324  const PtrList<dictionary>& boundaryDicts,
325  const word& defaultBoundaryPatchName,
326  const word& defaultBoundaryPatchType,
327  const bool syncPar
328 )
329 :
330  polyMesh
331  (
332  io,
333  std::move(points),
334  shapes,
335  boundaryFaces,
336  boundaryPatchNames,
337  boundaryDicts,
338  defaultBoundaryPatchName,
339  defaultBoundaryPatchType,
340  syncPar
341  ),
342  surfaceInterpolation(*this),
343  boundary_(*this, boundaryMesh()),
344  stitcher_(nullptr),
345  topoChanger_(nullptr),
346  distributor_(nullptr),
347  mover_(nullptr),
348  lduPtr_(nullptr),
349  polyFacesBfIOPtr_(nullptr),
350  polyFacesBfPtr_(nullptr),
351  polyBFaceOffsetsPtr_(nullptr),
352  polyBFaceOffsetPatchesPtr_(nullptr),
353  polyBFaceOffsetPatchFacesPtr_(nullptr),
354  polyBFacePatchesPtr_(nullptr),
355  polyBFacePatchFacesPtr_(nullptr),
356  ownerBfPtr_(nullptr),
357  curTimeIndex_(time().timeIndex()),
358  VPtr_(nullptr),
359  V0Ptr_(nullptr),
360  V00Ptr_(nullptr),
361  SfSlicePtr_(nullptr),
362  SfPtr_(nullptr),
363  magSfSlicePtr_(nullptr),
364  magSfPtr_(nullptr),
365  CSlicePtr_(nullptr),
366  CPtr_(nullptr),
367  CfSlicePtr_(nullptr),
368  CfPtr_(nullptr),
369  phiPtr_(nullptr)
370 {
371  if (debug)
372  {
373  Pout<< FUNCTION_NAME << "Constructing fvMesh from cellShapes" << endl;
374  }
375 }
376 
377 
379 (
380  const IOobject& io,
381  pointField&& points,
382  faceList&& faces,
383  labelList&& allOwner,
384  labelList&& allNeighbour,
385  const bool syncPar
386 )
387 :
388  polyMesh
389  (
390  io,
391  std::move(points),
392  std::move(faces),
393  std::move(allOwner),
394  std::move(allNeighbour),
395  syncPar
396  ),
397  surfaceInterpolation(*this),
398  boundary_(*this, boundaryMesh()),
399  stitcher_(nullptr),
400  topoChanger_(nullptr),
401  distributor_(nullptr),
402  mover_(nullptr),
403  lduPtr_(nullptr),
404  polyFacesBfIOPtr_(nullptr),
405  polyFacesBfPtr_(nullptr),
406  polyBFaceOffsetsPtr_(nullptr),
407  polyBFaceOffsetPatchesPtr_(nullptr),
408  polyBFaceOffsetPatchFacesPtr_(nullptr),
409  polyBFacePatchesPtr_(nullptr),
410  polyBFacePatchFacesPtr_(nullptr),
411  ownerBfPtr_(nullptr),
412  curTimeIndex_(time().timeIndex()),
413  VPtr_(nullptr),
414  V0Ptr_(nullptr),
415  V00Ptr_(nullptr),
416  SfSlicePtr_(nullptr),
417  SfPtr_(nullptr),
418  magSfSlicePtr_(nullptr),
419  magSfPtr_(nullptr),
420  CSlicePtr_(nullptr),
421  CPtr_(nullptr),
422  CfSlicePtr_(nullptr),
423  CfPtr_(nullptr),
424  phiPtr_(nullptr)
425 {
426  if (debug)
427  {
428  Pout<< FUNCTION_NAME << "Constructing fvMesh from components" << endl;
429  }
430 }
431 
432 
434 (
435  const IOobject& io,
436  pointField&& points,
437  faceList&& faces,
438  cellList&& cells,
439  const bool syncPar
440 )
441 :
442  polyMesh
443  (
444  io,
445  std::move(points),
446  std::move(faces),
447  std::move(cells),
448  syncPar
449  ),
450  surfaceInterpolation(*this),
451  boundary_(*this),
452  stitcher_(nullptr),
453  topoChanger_(nullptr),
454  distributor_(nullptr),
455  mover_(nullptr),
456  lduPtr_(nullptr),
457  polyFacesBfIOPtr_(nullptr),
458  polyFacesBfPtr_(nullptr),
459  polyBFaceOffsetsPtr_(nullptr),
460  polyBFaceOffsetPatchesPtr_(nullptr),
461  polyBFaceOffsetPatchFacesPtr_(nullptr),
462  polyBFacePatchesPtr_(nullptr),
463  polyBFacePatchFacesPtr_(nullptr),
464  ownerBfPtr_(nullptr),
465  curTimeIndex_(time().timeIndex()),
466  VPtr_(nullptr),
467  V0Ptr_(nullptr),
468  V00Ptr_(nullptr),
469  SfSlicePtr_(nullptr),
470  SfPtr_(nullptr),
471  magSfSlicePtr_(nullptr),
472  magSfPtr_(nullptr),
473  CSlicePtr_(nullptr),
474  CPtr_(nullptr),
475  CfSlicePtr_(nullptr),
476  CfPtr_(nullptr),
477  phiPtr_(nullptr)
478 {
479  if (debug)
480  {
481  Pout<< FUNCTION_NAME << "Constructing fvMesh from components" << endl;
482  }
483 }
484 
485 
487 :
488  polyMesh(Foam::move(mesh)),
489  surfaceInterpolation(Foam::move(mesh)),
490  boundary_(Foam::move(mesh.boundary_)),
491  stitcher_(Foam::move(mesh.stitcher_)),
492  topoChanger_(Foam::move(mesh.topoChanger_)),
493  distributor_(Foam::move(mesh.distributor_)),
494  mover_(Foam::move(mesh.mover_)),
495  lduPtr_(Foam::move(mesh.lduPtr_)),
496  polyFacesBfIOPtr_(Foam::move(mesh.polyFacesBfIOPtr_)),
497  polyFacesBfPtr_(Foam::move(mesh.polyFacesBfPtr_)),
498  polyBFaceOffsetsPtr_(Foam::move(mesh.polyBFaceOffsetsPtr_)),
499  polyBFaceOffsetPatchesPtr_(Foam::move(mesh.polyBFaceOffsetPatchesPtr_)),
500  polyBFaceOffsetPatchFacesPtr_
501  (
502  Foam::move(mesh.polyBFaceOffsetPatchFacesPtr_)
503  ),
504  polyBFacePatchesPtr_(Foam::move(mesh.polyBFacePatchesPtr_)),
505  polyBFacePatchFacesPtr_(Foam::move(mesh.polyBFacePatchFacesPtr_)),
506  ownerBfPtr_(Foam::move(mesh.ownerBfPtr_)),
507  curTimeIndex_(mesh.curTimeIndex_),
508  VPtr_(Foam::move(mesh.VPtr_)),
509  V0Ptr_(Foam::move(mesh.V0Ptr_)),
510  V00Ptr_(Foam::move(mesh.V00Ptr_)),
511  SfSlicePtr_(Foam::move(mesh.SfSlicePtr_)),
512  SfPtr_(Foam::move(mesh.SfPtr_)),
513  magSfSlicePtr_(Foam::move(mesh.magSfSlicePtr_)),
514  magSfPtr_(Foam::move(mesh.magSfPtr_)),
515  CSlicePtr_(Foam::move(mesh.CSlicePtr_)),
516  CPtr_(Foam::move(mesh.CPtr_)),
517  CfSlicePtr_(Foam::move(mesh.CfSlicePtr_)),
518  CfPtr_(Foam::move(mesh.CfPtr_)),
519  phiPtr_(Foam::move(mesh.phiPtr_))
520 {
521  if (debug)
522  {
523  Pout<< FUNCTION_NAME << "Moving fvMesh" << endl;
524  }
525 }
526 
527 
528 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
529 
531 {
532  clearOut();
533 }
534 
535 
536 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
537 
538 void Foam::fvMesh::postConstruct(const bool changers, const stitchType stitch)
539 {
540  // Construct the stitcher
541  stitcher_.set(fvMeshStitcher::New(*this, changers).ptr());
542 
543  // Stitch or Re-stitch if necessary
544  if (stitch != stitchType::none)
545  {
546  stitcher_->connect(false, stitch == stitchType::geometric, true);
547  }
548 
549  // Construct changers
550  if (changers)
551  {
552  topoChanger_.set(fvMeshTopoChanger::New(*this).ptr());
553  distributor_.set(fvMeshDistributor::New(*this).ptr());
554  mover_.set(fvMeshMover::New(*this).ptr());
555 
556  // Check the existence of the cell volumes and read if present
557  // and set the storage of V00
558  if (fileHandler().isFile(time().timePath()/"Vc0"))
559  {
561  (
562  IOobject
563  (
564  "Vc0",
565  time().name(),
566  *this,
569  true
570  ),
571  *this
572  );
573 
574  V00();
575  }
576 
577  // Check the existence of the mesh fluxes and read if present
578  if (fileHandler().isFile(time().timePath()/"meshPhi"))
579  {
580  phiPtr_ = new surfaceScalarField
581  (
582  IOobject
583  (
584  "meshPhi",
585  time().name(),
586  *this,
589  true
590  ),
591  *this
592  );
593  }
594  }
595 }
596 
597 
599 {
600  return topoChanger_.valid() && topoChanger_->dynamic();
601 }
602 
603 
605 {
606  return distributor_.valid() && distributor_->dynamic();
607 }
608 
609 
611 {
612  return
613  (topoChanger_.valid() && topoChanger_->dynamic())
614  || (mover_.valid() && mover_->dynamic());
615 }
616 
617 
619 {
620  if
621  (
622  stitcher_->stitches()
623  || topoChanger_->dynamic()
624  || distributor_->dynamic()
625  )
626  {
628  }
629 
630  if (!conformal()) stitcher_->disconnect(true, true);
631 
632  // Remove the oldest cell volume field
633  if (V00Ptr_)
634  {
635  nullDemandDrivenData(V00Ptr_);
636  }
637  else
638  {
639  nullDemandDrivenData(V0Ptr_);
640  }
641 
642  // Set topoChanged_ false before any mesh change
643  topoChanged_ = false;
644  bool updated = topoChanger_->update();
645  topoChanged_ = updated;
646 
647  updated = distributor_->update() || updated;
648 
649  return updated;
650 }
651 
652 
654 {
655  if (!conformal()) stitcher_->disconnect(true, true);
656 
657  if (curTimeIndex_ < time().timeIndex() && stitcher_->stitches())
658  {
659  // Store all old-time fields. If we don't do this then we risk
660  // triggering a store in the middle of mapping and potentially
661  // overwriting a mapped old-time field with a not-yet-mapped new-time
662  // field.
663  storeOldTimeFields();
664  }
665 
666  // Do not set moving false
667  // Once the mesh starts moving it is considered to be moving
668  // for the rest of the run
669  const bool moved = mover_->update();
670 
671  curTimeIndex_ = time().timeIndex();
672 
673  stitcher_->connect(true, true, false);
674 
675  return moved;
676 }
677 
678 
680 (
681  const List<polyPatch*>& p,
682  const bool validBoundary
683 )
684 {
685  if (boundary().size())
686  {
688  << " boundary already exists"
689  << abort(FatalError);
690  }
691 
692  // first add polyPatches
693  addPatches(p, validBoundary);
694  boundary_.addPatches(boundaryMesh());
695 }
696 
697 
699 {
700  if (debug)
701  {
702  Pout<< FUNCTION_NAME << "Removing boundary patches." << endl;
703  }
704 
705  // Remove fvBoundaryMesh data first.
706  boundary_.clear();
707  boundary_.setSize(0);
709 
710  clearOut();
711 }
712 
713 
714 void Foam::fvMesh::swap(fvMesh& otherMesh)
715 {
716  // Clear the sliced fields
717  clearGeom();
718 
719  // Clear the current volume and other geometry factors
721 
722  // Clear any non-updateable addressing
723  clearAddressing(true);
724 
725  polyMesh::swap(otherMesh);
726 
727  auto updatePatches = []
728  (
729  const polyPatchList& patches,
730  fvBoundaryMesh& boundaryMesh
731  )
732  {
733  boundaryMesh.setSize(patches.size());
734 
736  {
737  // Construct new processor patches, as the decomposition may have
738  // changed. Leave other patches as is.
739 
740  if (isA<processorPolyPatch>(patches[patchi]))
741  {
742  boundaryMesh.set
743  (
744  patchi,
746  (
747  patches[patchi],
748  boundaryMesh
749  )
750  );
751  }
752  }
753  };
754 
755  updatePatches(boundaryMesh(), boundary_);
756  updatePatches(otherMesh.boundaryMesh(), otherMesh.boundary_);
757 }
758 
759 
761 (
762  const stitchType stitch
763 )
764 {
765  if (debug)
766  {
767  Pout<< FUNCTION_NAME << "Updating fvMesh. ";
768  }
769 
771 
772  const fileName polyFacesInst =
773  time().findInstance
774  (
775  dbDir()/typeName,
776  "polyFaces",
778  );
779 
780  const bool reStitch =
781  stitcher_.valid()
782  && stitcher_->stitches()
783  && stitch != stitchType::none
784  && (!conformal() || polyFacesInst != polyFacesBfIOPtr_->instance());
785 
786  if (reStitch)
787  {
788  stitcher_->disconnect(false, false);
789  }
790  else if (state != polyMesh::UNCHANGED)
791  {
792  conform();
793  }
794 
795  if (state == polyMesh::TOPO_PATCH_CHANGE)
796  {
797  boundary_.readUpdate(boundaryMesh());
798  }
799 
800  if (state == polyMesh::TOPO_PATCH_CHANGE)
801  {
802  if (debug)
803  {
804  Info<< "Boundary and topological update" << endl;
805  }
806 
807  clearOut();
808  }
809  else if (state == polyMesh::TOPO_CHANGE)
810  {
811  if (debug)
812  {
813  Info<< "Topological update" << endl;
814  }
815 
816  clearOut();
817  }
818  else if (state == polyMesh::POINTS_MOVED)
819  {
820  if (debug)
821  {
822  Info<< "Point motion update" << endl;
823  }
824 
825  clearGeom();
826  }
827  else
828  {
829  if (debug)
830  {
831  Info<< "No update" << endl;
832  }
833  }
834 
835  if (reStitch && stitch != stitchType::none)
836  {
837  stitcher_->connect(false, stitch == stitchType::geometric, true);
838  }
839 
840  // Return the corresponding fvMesh read update state
841  switch (state)
842  {
843  case polyMesh::UNCHANGED:
844  return reStitch ? STITCHED : UNCHANGED;
846  return reStitch ? STITCHED : POINTS_MOVED;
848  return TOPO_CHANGE;
850  return TOPO_PATCH_CHANGE;
851  }
852 
853  return UNCHANGED;
854 }
855 
856 
858 {
859  return boundary_;
860 }
861 
862 
864 {
865  if (!lduPtr_)
866  {
867  lduPtr_ = new fvMeshLduAddressing(*this);
868  }
869 
870  return *lduPtr_;
871 }
872 
873 
875 {
876  return !(polyFacesBfPtr_ && SfPtr_);
877 }
878 
879 
881 {
882  if (!polyFacesBfPtr_)
883  {
884  polyFacesBfIOPtr_ =
885  new IOobject
886  (
887  "polyFaces",
888  facesInstance(),
889  typeName,
890  *this,
893  false
894  );
895 
896  polyFacesBfPtr_ =
898  (
899  boundary(),
901  polyFacesPatchTypes(),
902  boundaryMesh().types()
903  );
904  }
905 
906  return *polyFacesBfPtr_;
907 }
908 
909 
912 {
913  if (!polyBFacePatchesPtr_)
914  {
915  const label nPolyBFaces = nFaces() - nInternalFaces();
916 
917  // Count face-poly-bFaces to get the offsets
918  polyBFaceOffsetsPtr_ = new labelList(nPolyBFaces + 1, 0);
919  labelList& offsets = *polyBFaceOffsetsPtr_;
921  {
922  forAll(boundary()[patchi], patchFacei)
923  {
924  const label polyBFacei =
925  (
926  polyFacesBfPtr_
927  ? (*polyFacesBfPtr_)[patchi][patchFacei]
928  : boundary()[patchi].start() + patchFacei
929  )
930  - nInternalFaces();
931 
932  offsets[polyBFacei + 1] ++;
933  }
934  }
935  for (label polyBFacei = 0; polyBFacei < nPolyBFaces; ++ polyBFacei)
936  {
937  offsets[polyBFacei + 1] += offsets[polyBFacei];
938  }
939 
940  // Set the poly-bFace patches and patch-faces, using the offsets as
941  // counters
942  polyBFaceOffsetPatchesPtr_ = new labelList(offsets.last());
943  polyBFaceOffsetPatchFacesPtr_ = new labelList(offsets.last());
944  labelUList& patches = *polyBFaceOffsetPatchesPtr_;
945  labelUList& patchFaces = *polyBFaceOffsetPatchFacesPtr_;
947  {
948  forAll(boundary()[patchi], patchFacei)
949  {
950  const label polyBFacei =
951  (
952  polyFacesBfPtr_
953  ? (*polyFacesBfPtr_)[patchi][patchFacei]
954  : boundary()[patchi].start() + patchFacei
955  )
956  - nInternalFaces();
957  patches[offsets[polyBFacei]] = patchi;
958  patchFaces[offsets[polyBFacei]] = patchFacei;
959  offsets[polyBFacei] ++;
960  }
961  }
962 
963  // Restore the offsets by removing the count
964  for
965  (
966  label polyBFacei = nPolyBFaces - 1;
967  polyBFacei >= 0;
968  -- polyBFacei
969  )
970  {
971  offsets[polyBFacei + 1] = offsets[polyBFacei];
972  }
973  offsets[0] = 0;
974 
975  // List-lists
976  polyBFacePatchesPtr_ =
977  new UCompactListList<label>(offsets, patches);
978  polyBFacePatchFacesPtr_ =
979  new UCompactListList<label>(offsets, patchFaces);
980  }
981 
982  return *polyBFacePatchesPtr_;
983 }
984 
985 
988 {
989  if (!polyBFacePatchFacesPtr_)
990  {
991  polyBFacePatches();
992  }
993 
994  return *polyBFacePatchFacesPtr_;
995 }
996 
997 
999 {
1000  if (!ownerBfPtr_)
1001  {
1002  ownerBfPtr_ =
1004  (
1005  boundary(),
1007  wordList
1008  (
1009  boundary().size(),
1010  calculatedFvsPatchLabelField::typeName
1011  ),
1012  boundaryMesh().types()
1013  );
1014 
1015  forAll(boundary(), patchi)
1016  {
1017  (*ownerBfPtr_)[patchi] =
1018  labelField(faceOwner(), polyFacesBf()[patchi]);
1019  }
1020  }
1021 
1022  return *ownerBfPtr_;
1023 }
1024 
1025 
1027 {
1028  return stitcher_();
1029 }
1030 
1031 
1033 {
1034  return stitcher_();
1035 }
1036 
1037 
1039 {
1040  return topoChanger_();
1041 }
1042 
1043 
1045 {
1046  return distributor_();
1047 }
1048 
1049 
1051 {
1052  return mover_();
1053 }
1054 
1055 
1057 {
1058  if (debug)
1059  {
1061  << " nOldCells:" << map.nOldCells()
1062  << " nCells:" << nCells()
1063  << " nOldFaces:" << map.nOldFaces()
1064  << " nFaces:" << nFaces()
1065  << endl;
1066  }
1067 
1068 
1069  // We require geometric properties valid for the old mesh
1070  if
1071  (
1072  map.cellMap().size() != nCells()
1073  || map.faceMap().size() != nFaces()
1074  )
1075  {
1077  << "polyTopoChangeMap does not correspond to the old mesh."
1078  << " nCells:" << nCells()
1079  << " cellMap:" << map.cellMap().size()
1080  << " nOldCells:" << map.nOldCells()
1081  << " nFaces:" << nFaces()
1082  << " faceMap:" << map.faceMap().size()
1083  << " nOldFaces:" << map.nOldFaces()
1084  << exit(FatalError);
1085  }
1086 
1087  // Create a fv mapper
1088  const fvMeshMapper fvMap(*this, map);
1089 
1090  // Map all the volFields in the objectRegistry
1091  #define mapVolFieldType(Type, nullArg) \
1092  MapGeometricFields<Type, fvPatchField, fvMeshMapper, volMesh>(fvMap);
1094 
1095  // Map all the surfaceFields in the objectRegistry
1096  #define mapSurfaceFieldType(Type, nullArg) \
1097  MapGeometricFields<Type, fvsPatchField, fvMeshMapper, surfaceMesh> \
1098  (fvMap);
1100 
1101  // Map all the dimensionedFields in the objectRegistry
1102  #define mapVolInternalFieldType(Type, nullArg) \
1103  MapDimensionedFields<Type, fvMeshMapper, volMesh>(fvMap);
1105 
1106  if (pointMesh::found(*this))
1107  {
1108  // Create the pointMesh mapper
1109  const pointMeshMapper mapper(pointMesh::New(*this), map);
1110 
1111  #define mapPointFieldType(Type, nullArg) \
1112  MapGeometricFields \
1113  < \
1114  Type, \
1115  pointPatchField, \
1116  pointMeshMapper, \
1117  pointMesh \
1118  > \
1119  (mapper);
1121  }
1122 }
1123 
1124 
1126 {
1128 
1129  clearGeom();
1130 
1131  // Update other local data
1133 
1134  meshObjects::movePoints<fvMesh>(*this);
1135  meshObjects::movePoints<lduMesh>(*this);
1136 
1137  const_cast<Time&>(time()).functionObjects().movePoints(*this);
1138 }
1139 
1140 
1142 {
1143  // Set the mesh to be moving. This remains true for the rest of the run.
1144  moving_ = true;
1145 
1146  // Create old-time volumes, if necessary, at the start of a new timestep
1147  if (curTimeIndex_ < time().timeIndex())
1148  {
1149  if (V00Ptr_ && notNull(V00Ptr_))
1150  {
1152  << "Old-old volumes should not be maintained across mesh "
1153  << "changes" << exit(FatalError);
1154  }
1155 
1156  // If old-old-volumes are necessary then copy them from the old-volumes
1157  if (Foam::isNull(V00Ptr_))
1158  {
1159  V00Ptr_ = new DimensionedField<scalar, volMesh>
1160  (
1161  IOobject
1162  (
1163  "Vc00",
1164  time().name(),
1165  *this,
1168  true
1169  ),
1170  V0()
1171  );
1172  }
1173 
1174  // Copy old-volumes from the volumes
1175  if (!V0Ptr_ || Foam::isNull(V0Ptr_))
1176  {
1178  (
1179  IOobject
1180  (
1181  "Vc0",
1182  time().name(),
1183  *this,
1186  true
1187  ),
1188  V()
1189  );
1190  }
1191  else
1192  {
1193  V0Ptr_->scalarField::operator=(V());
1194  }
1195  }
1196 
1197  // Create mesh motion flux, if necessary
1198  if (!phiPtr_)
1199  {
1200  phiPtr_ = new surfaceScalarField
1201  (
1202  IOobject
1203  (
1204  "meshPhi",
1205  this->time().name(),
1206  *this,
1209  true
1210  ),
1211  *this,
1213  );
1214  }
1215  else
1216  {
1217  phiPtr_->storeOldTimes();
1218  }
1219 
1220  surfaceScalarField& phi = *phiPtr_;
1221 
1222  // Move the polyMesh and set the mesh motion fluxes to the swept-volumes
1223 
1224  scalar rDeltaT = 1.0/time().deltaTValue();
1225 
1226  tmp<scalarField> tsweptVols = polyMesh::movePoints(p);
1227 
1228  scalarField& sweptVols = tsweptVols.ref();
1229 
1230  phi.primitiveFieldRef() =
1231  scalarField::subField(sweptVols, nInternalFaces());
1232  phi.primitiveFieldRef() *= rDeltaT;
1233 
1234  const fvPatchList& patches = boundary();
1235 
1237 
1239  {
1240  phibf[patchi] = patches[patchi].patchSlice(sweptVols);
1241  phibf[patchi] *= rDeltaT;
1242  }
1243 
1244  // Update or delete the local geometric properties as early as possible so
1245  // they can be used if necessary. These get recreated here instead of
1246  // demand driven since they might do parallel transfers which can conflict
1247  // with when they're actually being used.
1248  // Note that between above "polyMesh::movePoints(p)" and here nothing
1249  // should use the local geometric properties.
1250  updateGeomNotOldVol();
1251 
1252  // Update other local data
1254 
1255  meshObjects::movePoints<fvMesh>(*this);
1256  meshObjects::movePoints<lduMesh>(*this);
1257 
1258  const_cast<Time&>(time()).functionObjects().movePoints(*this);
1259 
1260  return tsweptVols;
1261 }
1262 
1263 
1265 {
1266  // Update polyMesh. This needs to keep volume existent!
1267  polyMesh::topoChange(map);
1268 
1269  // Clear the sliced fields
1270  clearGeomNotOldVol();
1271 
1272  // Check that we're not trying to maintain old-time mesh geometry
1273  if (V0Ptr_ && Foam::notNull(V0Ptr_))
1274  {
1276  << "It is not possible to use mesh motion, topology change, and "
1277  << "second order time schemes simultaneously"
1278  << exit(FatalError);
1279  }
1280 
1281  // Map all fields
1282  mapFields(map);
1283 
1284  // Clear the current volume and other geometry factors
1286 
1287  // Clear any non-updateable addressing
1288  clearAddressing(true);
1289 
1290  meshObjects::topoChange<fvMesh>(*this, map);
1291  meshObjects::topoChange<lduMesh>(*this, map);
1292 
1293  const_cast<Time&>(time()).functionObjects().topoChange(map);
1294 
1295  if (stitcher_.valid())
1296  {
1297  stitcher_->topoChange(map);
1298  }
1299 
1300  if (topoChanger_.valid())
1301  {
1302  topoChanger_->topoChange(map);
1303  }
1304 
1305  if (distributor_.valid())
1306  {
1307  distributor_->topoChange(map);
1308  }
1309 
1310  if (mover_.valid())
1311  {
1312  mover_->topoChange(map);
1313  }
1314 }
1315 
1316 
1318 {
1319  // Distribute polyMesh data
1320  polyMesh::mapMesh(map);
1321 
1322  // Clear the sliced fields
1323  clearGeomNotOldVol();
1324 
1325  // Clear the current volume and other geometry factors
1327 
1328  // Clear any non-updateable addressing
1329  clearAddressing(true);
1330 
1331  meshObjects::mapMesh<fvMesh>(*this, map);
1332  meshObjects::mapMesh<lduMesh>(*this, map);
1333 
1334  const_cast<Time&>(time()).functionObjects().mapMesh(map);
1335 
1336  stitcher_->mapMesh(map);
1337  topoChanger_->mapMesh(map);
1338  distributor_->mapMesh(map);
1339  mover_->mapMesh(map);
1340 }
1341 
1342 
1344 {
1345  // Distribute polyMesh data
1346  polyMesh::distribute(map);
1347 
1348  // Clear the sliced fields
1349  clearGeomNotOldVol();
1350 
1351  // Clear the current volume and other geometry factors
1353 
1354  // Clear any non-updateable addressing
1355  clearAddressing(true);
1356 
1357  meshObjects::distribute<fvMesh>(*this, map);
1358  meshObjects::distribute<lduMesh>(*this, map);
1359 
1360  const_cast<Time&>(time()).functionObjects().distribute(map);
1361 
1362  stitcher_->distribute(map);
1363  topoChanger_->distribute(map);
1364  distributor_->distribute(map);
1365  mover_->distribute(map);
1366 }
1367 
1368 
1370 {
1371  if (!polyFacesBfPtr_)
1372  {
1373  polyFacesBf();
1374  }
1375 
1376  polyFacesBfIOPtr_->instance() = inst;
1377  polyFacesBfIOPtr_->writeOpt() = IOobject::AUTO_WRITE;
1378 }
1379 
1380 
1382 {
1383  if (!polyFacesBfPtr_)
1384  {
1385  return facesInstance();
1386  }
1387 
1388  return polyFacesBfIOPtr_->instance();
1389 }
1390 
1391 
1393 {
1394  // Clear the geometry fields
1395  clearGeomNotOldVol();
1396 
1397  // Clear the current volume and other geometry factors
1399 
1400  // Clear any non-updateable addressing
1401  clearAddressing(true);
1402 
1403  // Modify the mesh fluxes, if necessary
1404  if (notNull(phi) && phiPtr_)
1405  {
1406  for (label i = 0; i <= phi.nOldTimes(); ++ i)
1407  {
1408  phiRef().oldTimeRef(i) == phi.oldTime(i);
1409  }
1410  }
1411 }
1412 
1413 
1415 (
1416  const surfaceLabelField::Boundary& polyFacesBf,
1417  const surfaceVectorField& Sf,
1418  const surfaceVectorField& Cf,
1419  const surfaceScalarField& phi,
1420  const bool sync
1421 )
1422 {
1423  // Clear the geometry fields
1424  clearGeomNotOldVol();
1425 
1426  // Clear the current volume and other geometry factors
1428 
1429  // Clear any non-updateable addressing
1430  clearAddressing(true);
1431 
1432  // Create non-sliced copies of geometry fields
1433  SfRef();
1434  magSfRef();
1435  CRef();
1436  CfRef();
1437 
1438  // Set the topology
1439  forAll(polyFacesBf, patchi)
1440  {
1441  polyFacesBfRef()[patchi].reset(polyFacesBf[patchi]);
1442  }
1443 
1444  // Set the face geometry
1445  SfRef() == Sf;
1446  magSfRef() == max(mag(Sf), dimensionedScalar(dimArea, rootVSmall));
1447  CRef().boundaryFieldRef() == Cf.boundaryField();
1448  CfRef() == Cf;
1449 
1450  // Communicate processor-coupled cell geometry. Cell-centre processor patch
1451  // fields must contain the (transformed) cell-centre locations on the other
1452  // side of the coupling. This is so that non-conformal patches can
1453  // construct weights and deltas without reference to the poly mesh
1454  // geometry.
1455  //
1456  // Note that the initEvaluate/evaluate communication does a transformation,
1457  // but it is wrong in this instance. A vector field gets transformed as if
1458  // it were a displacement, but the cell-centres need a positional
1459  // transformation. That's why there's the un-transform and re-transform bit
1460  // below just after the evaluate call.
1461  //
1462  // This transform handling is a bit of a hack. It would be nicer to have a
1463  // field attribute which identifies a field as needing a positional
1464  // transformation, and for it to apply automatically within the coupled
1465  // patch field. However, at the moment, the cell centres field is the only
1466  // vol-field containing an absolute position, so the hack is functionally
1467  // sufficient for now.
1468  if (sync && (Pstream::parRun() || !time().processorCase()))
1469  {
1470  volVectorField::Boundary& CBf = CRef().boundaryFieldRef();
1471 
1472  const label nReq = Pstream::nRequests();
1473 
1474  forAll(CBf, patchi)
1475  {
1476  if (isA<processorFvPatch>(CBf[patchi].patch()))
1477  {
1478  CBf[patchi].initEvaluate(Pstream::defaultCommsType);
1479  }
1480  }
1481 
1482  if
1483  (
1484  Pstream::parRun()
1486  )
1487  {
1488  Pstream::waitRequests(nReq);
1489  }
1490 
1491  forAll(CBf, patchi)
1492  {
1493  if (isA<processorFvPatch>(CBf[patchi].patch()))
1494  {
1496 
1497  const transformer& t =
1498  refCast<const processorFvPatch>(CBf[patchi].patch())
1499  .transform();
1500 
1501  t.invTransform(CBf[patchi], CBf[patchi]);
1502  t.transformPosition(CBf[patchi], CBf[patchi]);
1503  }
1504  }
1505  }
1506 
1507  // Modify the mesh fluxes, if necessary
1508  if (notNull(phi) && phiPtr_)
1509  {
1510  for (label i = 0; i <= phi.nOldTimes(); ++ i)
1511  {
1512  phiRef().oldTimeRef(i) == phi.oldTime(i);
1513  }
1514  }
1515 }
1516 
1517 
1519 (
1520  const label insertPatchi,
1521  const polyPatch& patch,
1522  const dictionary& patchFieldDict,
1523  const word& defaultPatchFieldType,
1524  const bool validBoundary
1525 )
1526 {
1527  // Remove my local data (see topoChange)
1528  // Clear mesh motion flux
1529  deleteDemandDrivenData(phiPtr_);
1530 
1531  // Clear the sliced fields
1532  clearGeomNotOldVol();
1533 
1534  // Clear the current volume and other geometry factors
1536 
1537  // Clear any non-updateable addressing
1538  clearAddressing(true);
1539 
1540 
1541  const label sz = boundary_.size();
1542 
1544  (
1545  insertPatchi,
1546  patch,
1547  patchFieldDict,
1548  defaultPatchFieldType,
1549  validBoundary
1550  );
1551 
1552  boundary_.setSize(sz+1);
1553  boundary_.set
1554  (
1555  insertPatchi,
1556  fvPatch::New
1557  (
1558  boundaryMesh()[insertPatchi],
1559  boundary_
1560  )
1561  );
1562 
1563  objectRegistry& db = const_cast<objectRegistry&>(thisDb());
1564  AddPatchFields<volScalarField>
1565  (
1566  db,
1567  insertPatchi,
1568  patchFieldDict,
1569  defaultPatchFieldType,
1570  Zero
1571  );
1572  AddPatchFields<volVectorField>
1573  (
1574  db,
1575  insertPatchi,
1576  patchFieldDict,
1577  defaultPatchFieldType,
1578  Zero
1579  );
1580  AddPatchFields<volSphericalTensorField>
1581  (
1582  db,
1583  insertPatchi,
1584  patchFieldDict,
1585  defaultPatchFieldType,
1586  Zero
1587  );
1588  AddPatchFields<volSymmTensorField>
1589  (
1590  db,
1591  insertPatchi,
1592  patchFieldDict,
1593  defaultPatchFieldType,
1594  Zero
1595  );
1596  AddPatchFields<volTensorField>
1597  (
1598  db,
1599  insertPatchi,
1600  patchFieldDict,
1601  defaultPatchFieldType,
1602  Zero
1603  );
1604 
1605  // Surface fields
1606 
1607  AddPatchFields<surfaceScalarField>
1608  (
1609  db,
1610  insertPatchi,
1611  patchFieldDict,
1612  defaultPatchFieldType,
1613  Zero
1614  );
1615  AddPatchFields<surfaceVectorField>
1616  (
1617  db,
1618  insertPatchi,
1619  patchFieldDict,
1620  defaultPatchFieldType,
1621  Zero
1622  );
1623  AddPatchFields<surfaceSphericalTensorField>
1624  (
1625  db,
1626  insertPatchi,
1627  patchFieldDict,
1628  defaultPatchFieldType,
1629  Zero
1630  );
1631  AddPatchFields<surfaceSymmTensorField>
1632  (
1633  db,
1634  insertPatchi,
1635  patchFieldDict,
1636  defaultPatchFieldType,
1637  Zero
1638  );
1639  AddPatchFields<surfaceTensorField>
1640  (
1641  db,
1642  insertPatchi,
1643  patchFieldDict,
1644  defaultPatchFieldType,
1645  Zero
1646  );
1647 }
1648 
1649 
1651 (
1652  const labelUList& newToOld,
1653  const bool validBoundary
1654 )
1655 {
1656  polyMesh::reorderPatches(newToOld, validBoundary);
1657 
1658  boundary_.shuffle(newToOld, validBoundary);
1659 
1660  objectRegistry& db = const_cast<objectRegistry&>(thisDb());
1661  ReorderPatchFields<volScalarField>(db, newToOld);
1662  ReorderPatchFields<volVectorField>(db, newToOld);
1663  ReorderPatchFields<volSphericalTensorField>(db, newToOld);
1664  ReorderPatchFields<volSymmTensorField>(db, newToOld);
1665  ReorderPatchFields<volTensorField>(db, newToOld);
1666 
1667  ReorderPatchFields<surfaceScalarField>(db, newToOld);
1668  ReorderPatchFields<surfaceVectorField>(db, newToOld);
1669  ReorderPatchFields<surfaceSphericalTensorField>(db, newToOld);
1670  ReorderPatchFields<surfaceSymmTensorField>(db, newToOld);
1671  ReorderPatchFields<surfaceTensorField>(db, newToOld);
1672 }
1673 
1674 
1676 (
1680  const bool write
1681 ) const
1682 {
1683  bool ok = true;
1684 
1685  if (!conformal() && polyFacesBfIOPtr_->writeOpt() == IOobject::AUTO_WRITE)
1686  {
1687  // Create a full surface field with the polyFacesBf boundary field to
1688  // write to disk. Make the internal field uniform to save disk space.
1689 
1690  surfaceLabelField polyFaces
1691  (
1692  *polyFacesBfIOPtr_,
1693  *this,
1694  dimless,
1695  labelField(nInternalFaces(), -1),
1696  *polyFacesBfPtr_
1697  );
1698 
1699  ok = ok & polyFaces.write(write);
1700  }
1701 
1702  if (phiPtr_)
1703  {
1704  ok = ok && phiPtr_->write(write);
1705  }
1706 
1707  // Write V0 only if V00 exists
1708  if (V00Ptr_)
1709  {
1710  ok = ok && V0Ptr_->write(write);
1711  }
1712 
1713  if (stitcher_.valid())
1714  {
1715  stitcher_->write(write);
1716  }
1717 
1718  if (topoChanger_.valid())
1719  {
1720  topoChanger_->write(write);
1721  }
1722 
1723  if (distributor_.valid())
1724  {
1725  distributor_->write(write);
1726  }
1727 
1728  if (mover_.valid())
1729  {
1730  mover_->write(write);
1731  }
1732 
1733  return ok && polyMesh::writeObject(fmt, ver, cmp, write);
1734 }
1735 
1736 
1737 bool Foam::fvMesh::write(const bool write) const
1738 {
1739  return polyMesh::write(write);
1740 }
1741 
1742 
1743 template<>
1745 Foam::fvMesh::validComponents<Foam::sphericalTensor>() const
1746 {
1748 }
1749 
1750 
1752 {
1753  if (!fvSchemes_.valid())
1754  {
1755  fvSchemes_ = new fvSchemes(*this);
1756  }
1757 
1758  return fvSchemes_;
1759 }
1760 
1761 
1763 {
1764  if (!fvSolution_.valid())
1765  {
1766  fvSolution_ = new fvSolution(*this);
1767  }
1768 
1769  return fvSolution_;
1770 }
1771 
1772 
1773 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1774 
1775 bool Foam::fvMesh::operator!=(const fvMesh& bm) const
1776 {
1777  return &bm != this;
1778 }
1779 
1780 
1781 bool Foam::fvMesh::operator==(const fvMesh& bm) const
1782 {
1783  return &bm == this;
1784 }
1785 
1786 
1787 // ************************************************************************* //
static const Foam::dimensionedScalar C("C", Foam::dimTemperature, 234.5)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
static bool found(const word &name, const polyMesh &mesh)
Return true if the DemandDrivenMeshObject with the given name.
static pointMesh & New(const word &name, const polyMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
SubField< scalar > subField
Declare type of subField.
Definition: Field.H:101
static const char *const typeName
Definition: Field.H:106
Generic GeometricBoundaryField class.
void evaluate()
Evaluate boundary conditions.
Generic GeometricField class.
Internal::FieldType & primitiveFieldRef()
Return a reference to the primitive field.
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Version number type.
Definition: IOstream.H:97
streamFormat
Enumeration for the format of data in the stream.
Definition: IOstream.H:87
compressionType
Enumeration for the format of data in the stream.
Definition: IOstream.H:194
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
const FieldType & oldTime() const
Return the old time field.
Definition: OldTimeField.C:315
label nOldTimes(const bool includeNull=true) const
Return the number of old time fields stored.
Definition: OldTimeField.C:292
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
Unallocated base class of CompactListList.
T & last()
Return the last element of the list.
Definition: UListI.H:128
static label nRequests()
Get number of outstanding requests.
Definition: UPstream.C:137
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
Definition: UPstream.C:147
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:272
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
A class for handling file names.
Definition: fileName.H:82
static bool valid(char)
Is this character valid for a fileName?
Definition: fileNameI.H:102
Foam::fvBoundaryMesh.
Abstract base class for fvMesh movers.
static autoPtr< fvMeshDistributor > New(fvMesh &)
Select, construct and return the fvMeshDistributor.
Foam::fvMeshLduAddressing.
Class holds all the necessary information for mapping fields associated with fvMesh.
Definition: fvMeshMapper.H:56
Abstract base class for fvMesh movers.
Definition: fvMeshMover.H:53
static autoPtr< fvMeshMover > New(fvMesh &)
Select, construct and return the fvMeshMover.
Mesh manipulator that uses the intersection provided by the cyclic non-conformal poly patches to crea...
static autoPtr< fvMeshStitcher > New(fvMesh &, bool changing)
Select, construct and return the fvMeshStitcher.
Abstract base class for fvMesh topology changers.
static autoPtr< fvMeshTopoChanger > New(fvMesh &, const dictionary &dict)
Select, construct and return the fvMeshTopoChanger.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const fileName & polyFacesBfInstance() const
Get the current instance for the poly faces boundary field.
Definition: fvMesh.C:1381
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:857
stitchType
Extent to which to stitch on read and readUpdate. By default, a.
Definition: fvMesh.H:123
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:863
virtual void addPatch(const label insertPatchi, const polyPatch &patch, const dictionary &patchFieldDict, const word &defaultPatchFieldType, const bool validBoundary)
Add/insert single patch. If validBoundary the new situation.
Definition: fvMesh.C:1519
const UCompactListList< label > & polyBFacePatches() const
Return poly-bFace-patch addressing.
Definition: fvMesh.C:911
const fvMeshDistributor & distributor() const
Return the distributor function class.
Definition: fvMesh.C:1044
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1751
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: fvMesh.H:108
virtual ~fvMesh()
Destructor.
Definition: fvMesh.C:530
const GeometricBoundaryField< label, fvsPatchField, surfaceMesh > & ownerBf() const
Return face-owner addressing.
Definition: fvMesh.C:998
static const HashSet< word > geometryFields
Set of names of registered geometric fields.
Definition: fvMesh.H:331
void unconform(const GeometricBoundaryField< label, fvsPatchField, surfaceMesh > &polyFacesBf, const surfaceVectorField &Sf, const surfaceVectorField &Cf, const surfaceScalarField &phi=NullObjectRef< surfaceScalarField >(), const bool sync=true)
Unconform the fvMesh from the polyMesh.
Definition: fvMesh.C:1415
const fvMeshStitcher & stitcher() const
Return the stitcher function class.
Definition: fvMesh.C:1026
const fvSolution & solution() const
Return the fvSolution.
Definition: fvMesh.C:1762
void addFvPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:680
void swap(fvMesh &)
Swap mesh.
Definition: fvMesh.C:714
const GeometricBoundaryField< label, fvsPatchField, surfaceMesh > & polyFacesBf() const
Return face-poly-face addressing.
Definition: fvMesh.C:880
bool operator==(const fvMesh &) const
Definition: fvMesh.C:1781
const UCompactListList< label > & polyBFacePatchFaces() const
Return poly-bFace-patch-face addressing.
Definition: fvMesh.C:987
bool operator!=(const fvMesh &) const
Definition: fvMesh.C:1775
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp, const bool write=true) const
Write the underlying polyMesh and other data.
Definition: fvMesh.C:1676
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: fvMesh.C:1317
void mapFields(const polyTopoChangeMap &map)
Map all fields in time using given map.
Definition: fvMesh.C:1056
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1737
bool update()
Update the mesh for topology change, mesh to mesh mapping.
Definition: fvMesh.C:618
void setPolyFacesBfInstance(const fileName &)
Set the instance for the poly faces boundary field.
Definition: fvMesh.C:1369
virtual void distribute(const polyDistributionMap &map)
Redistribute or update using the given distribution map.
Definition: fvMesh.C:1343
void removeFvBoundary()
Remove boundary patches. Warning: fvPatchFields hold ref to.
Definition: fvMesh.C:698
bool distributing() const
Does the mesh get redistributed?
Definition: fvMesh.C:604
virtual void reorderPatches(const labelUList &newToOld, const bool validBoundary)
Reorder and trim existing patches. If validBoundary the new.
Definition: fvMesh.C:1651
bool dynamic() const
Is this mesh dynamic?
Definition: fvMesh.C:610
fvMesh(const IOobject &io, const bool doPost=true)
Construct from IOobject. Optionally prevent post-construction, so.
Definition: fvMesh.C:273
const fvMeshTopoChanger & topoChanger() const
Return the topo-changer function class.
Definition: fvMesh.C:1038
virtual void topoChange(const polyTopoChangeMap &map)
Update mesh corresponding to the given map.
Definition: fvMesh.C:1264
bool topoChanging() const
Does the mesh topology change?
Definition: fvMesh.C:598
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:223
bool move()
Move the mesh.
Definition: fvMesh.C:653
virtual void setPoints(const pointField &)
Reset the points.
Definition: fvMesh.C:1125
void postConstruct(const bool changers, const stitchType stitch)
Complete construction of the mesh.
Definition: fvMesh.C:538
const fvMeshMover & mover() const
Return the mover function class.
Definition: fvMesh.C:1050
void conform(const surfaceScalarField &phi=NullObjectRef< surfaceScalarField >())
Conform the fvMesh to the polyMesh.
Definition: fvMesh.C:1392
bool conformal() const
Return whether the fvMesh is conformal with the polyMesh.
Definition: fvMesh.C:874
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
static autoPtr< fvPatch > New(const polyPatch &, const fvBoundaryMesh &)
Return a pointer to a new patch created on freestore from polyPatch.
Definition: fvPatchNew.C:32
Selector class for finite volume differencing schemes. fvMesh is derived from fvSchemes so that all f...
Definition: fvSchemes.H:53
Selector class for finite volume solution solution. fvMesh is derived from fvSolution so that all fie...
Definition: fvSolution.H:50
The class contains the addressing required by the lduMatrix: upper, lower and losort.
static void clearUpto(objectRegistry &)
Clear all meshObjects derived from FromType up to (but not including)
Registry of regIOobjects.
Traits class for primitives.
Definition: pTraits.H:53
Class holds all the necessary information for mapping fields associated with pointMesh.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: polyMesh.C:1444
virtual void addPatch(const label insertPatchi, const polyPatch &patch, const dictionary &patchFieldDict, const word &defaultPatchFieldType, const bool validBoundary)
Add/insert single patch. If validBoundary the new situation.
Definition: polyMesh.C:1231
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:90
@ TOPO_PATCH_CHANGE
Definition: polyMesh.H:94
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
void swap(polyMesh &)
Swap mesh.
Definition: polyMesh.C:804
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp, const bool write=true) const
Write the underlying polyMesh.
Definition: polyMeshIO.C:460
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
readUpdateState readUpdate()
Update the mesh based on the mesh files saved in.
Definition: polyMeshIO.C:99
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
virtual void distribute(const polyDistributionMap &map)
Redistribute or update using the given distribution map.
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:36
virtual void reorderPatches(const labelUList &newToOld, const bool validBoundary)
Reorder and trim existing patches. If validBoundary the new.
Definition: polyMesh.C:1186
void clearOut()
Clear all geometry and addressing unnecessary for CFD.
virtual void setPoints(const pointField &)
Reset the points.
Definition: polyMesh.C:1401
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
label nOldCells() const
Number of old cells.
const labelList & cellMap() const
Old cell map.
label nOldFaces() const
Number of old faces.
const labelList & faceMap() const
Old face map.
void clearAddressing()
Clear topological data.
virtual bool write(const bool write=true) const
Write using setting from DB.
Cell to surface interpolation scheme. Included in fvMesh.
bool movePoints()
Do what is necessary if the mesh has moved.
void clearOut()
Clear all geometry and addressing.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space.
Definition: transformer.H:84
Type invTransform(const Type &) const
Inverse transform the given type.
vector transformPosition(const vector &v) const
Transform the given position.
Definition: transformerI.H:153
Type transform(const Type &) const
Transform the given type.
A class for handling words, derived from string.
Definition: word.H:62
Template functions to aid in the implementation of demand driven data.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
#define nullOldestTimeFields(Type, nullArg)
#define mapVolInternalFieldType(Type, nullArg)
#define mapPointFieldType(Type, nullArg)
#define mapVolFieldType(Type, nullArg)
#define mapSurfaceFieldType(Type, nullArg)
Reconstruct volField from a face flux field.
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
label patchi
const pointField & points
const cellShapeList & cells
const fvPatchList & patches
#define FUNCTION_NAME
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
bool isFile(const fileName &, const bool checkVariants=true, const bool followLink=true)
Does the name exist as a file in the file system?
Definition: POSIX.C:555
const fileOperation & fileHandler()
Get current file handler.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
void deleteDemandDrivenData(DataType *&dataPtr)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
SurfaceField< scalar > surfaceScalarField
messageStream Info
const dimensionSet dimTime
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:64
const dimensionSet dimVolume
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:49
bool isNull(const T &t)
Return true if t is a reference to the nullObject of type T.
Definition: nullObjectI.H:58
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
error FatalError
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
const dimensionSet dimArea
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void nullDemandDrivenData(DataType *&dataPtr)
label timeIndex
Definition: getTimeIndex.H:4
faceListList boundary(nPatches)
volScalarField & p
Foam::surfaceFields.