fvMesh.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017,2022 OpenFOAM Foundation
9 Copyright (C) 2016-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
46namespace Foam
47{
49}
50
51
52// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53
55{
57 <
58 fvMesh,
61 >(*this);
62
64 <
65 lduMesh,
68 >(*this);
69
75}
76
77
79{
80 const bool haveV = (VPtr_ != nullptr);
81 const bool haveSf = (SfPtr_ != nullptr);
82 const bool haveMagSf = (magSfPtr_ != nullptr);
83 const bool haveCP = (CPtr_ != nullptr);
84 const bool haveCf = (CfPtr_ != nullptr);
85
86 clearGeomNotOldVol();
87
88 // Now recreate the fields
89 if (haveV)
90 {
91 (void)V();
92 }
93
94 if (haveSf)
95 {
96 (void)Sf();
97 }
98
99 if (haveMagSf)
100 {
101 (void)magSf();
102 }
103
104 if (haveCP)
105 {
106 (void)C();
107 }
108
109 if (haveCf)
110 {
111 (void)Cf();
112 }
113}
114
115
117{
118 clearGeomNotOldVol();
119
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
128void 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 (
190 (
191 "V0",
192 time().timeName(),
193 *this,
196 false
197 ),
198 *this,
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}
237
238
240{
241 clearOutLocal();
243}
244
245
246// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
247
248Foam::fvMesh::fvMesh(const IOobject& io, const bool doInit)
249:
250 polyMesh(io, doInit),
251 fvSchemes(static_cast<const objectRegistry&>(*this)),
253 fvSolution(static_cast<const objectRegistry&>(*this)),
254 data(static_cast<const objectRegistry&>(*this)),
255 boundary_(*this, boundaryMesh()),
256 lduPtr_(nullptr),
257 curTimeIndex_(time().timeIndex()),
258 VPtr_(nullptr),
259 V0Ptr_(nullptr),
260 V00Ptr_(nullptr),
261 SfPtr_(nullptr),
262 magSfPtr_(nullptr),
263 CPtr_(nullptr),
264 CfPtr_(nullptr),
265 phiPtr_(nullptr)
266{
267 DebugInFunction << "Constructing fvMesh from IOobject" << endl;
268
269 if (doInit)
270 {
271 fvMesh::init(false); // do not initialise lower levels
272 }
273}
274
275
276bool Foam::fvMesh::init(const bool doInit)
277{
278 if (doInit)
279 {
280 // Construct basic geometry calculation engine. Note: do before
281 // doing anything with primitiveMesh::cellCentres etc.
282 (void)geometry();
283
284 // Initialise my data
285 polyMesh::init(doInit);
286 }
287
288 // Check the existence of the cell volumes and read if present
289 // and set the storage of V00
290 if (fileHandler().isFile(time().timePath()/dbDir()/"V0"))
291 {
292 // Set the moving flag early in case the demand-driven geometry
293 // construction checks for it
294 moving(true);
295
297 (
299 (
300 "V0",
301 time().timeName(),
302 *this,
305 false
306 ),
307 *this
308 );
309
310 V00();
311 }
312
313 // Check the existence of the mesh fluxes, read if present and set the
314 // mesh to be moving
315 if (fileHandler().isFile(time().timePath()/dbDir()/"meshPhi"))
316 {
317 // Set the moving flag early in case the demand-driven geometry
318 // construction checks for it
319 moving(true);
320
321 phiPtr_ = new surfaceScalarField
322 (
324 (
325 "meshPhi",
326 time().timeName(),
327 *this,
330 false
331 ),
332 *this
333 );
334
335 // The mesh is now considered moving so the old-time cell volumes
336 // will be required for the time derivatives so if they haven't been
337 // read initialise to the current cell volumes
338 if (!V0Ptr_)
339 {
341 (
343 (
344 "V0",
345 time().timeName(),
346 *this,
349 false
350 ),
351 V()
352 );
353 }
354 }
355
356 // Assume something changed
357 return true;
358}
359
360
362(
363 const IOobject& io,
365 faceList&& faces,
366 labelList&& allOwner,
367 labelList&& allNeighbour,
368 const bool syncPar
369)
370:
372 (
373 io,
374 std::move(points),
375 std::move(faces),
376 std::move(allOwner),
377 std::move(allNeighbour),
378 syncPar
379 ),
380 fvSchemes(static_cast<const objectRegistry&>(*this)),
382 fvSolution(static_cast<const objectRegistry&>(*this)),
383 data(static_cast<const objectRegistry&>(*this)),
384 boundary_(*this),
385 lduPtr_(nullptr),
386 curTimeIndex_(time().timeIndex()),
387 VPtr_(nullptr),
388 V0Ptr_(nullptr),
389 V00Ptr_(nullptr),
390 SfPtr_(nullptr),
391 magSfPtr_(nullptr),
392 CPtr_(nullptr),
393 CfPtr_(nullptr),
394 phiPtr_(nullptr)
395{
396 DebugInFunction << "Constructing fvMesh from components" << endl;
397}
398
399
401(
402 const IOobject& io,
404 faceList&& faces,
405 cellList&& cells,
406 const bool syncPar
407)
408:
410 (
411 io,
412 std::move(points),
413 std::move(faces),
414 std::move(cells),
415 syncPar
416 ),
417 fvSchemes(static_cast<const objectRegistry&>(*this)),
419 fvSolution(static_cast<const objectRegistry&>(*this)),
420 data(static_cast<const objectRegistry&>(*this)),
421 boundary_(*this),
422 lduPtr_(nullptr),
423 curTimeIndex_(time().timeIndex()),
424 VPtr_(nullptr),
425 V0Ptr_(nullptr),
426 V00Ptr_(nullptr),
427 SfPtr_(nullptr),
428 magSfPtr_(nullptr),
429 CPtr_(nullptr),
430 CfPtr_(nullptr),
431 phiPtr_(nullptr)
432{
433 DebugInFunction << "Constructing fvMesh from components" << endl;
434}
435
436
437Foam::fvMesh::fvMesh(const IOobject& io, const Foam::zero, const bool syncPar)
438:
439 fvMesh(io, pointField(), faceList(), labelList(), labelList(), syncPar)
440{}
441
442
444(
445 const IOobject& io,
446 const fvMesh& baseMesh,
447 const Foam::zero,
448 const bool syncPar
449)
450:
451 fvMesh
452 (
453 io,
454 baseMesh,
455 pointField(),
456 faceList(),
457 labelList(), // owner
458 labelList(), // neighbour
459 syncPar
460 )
461{}
462
463
465(
466 const IOobject& io,
467 const fvMesh& baseMesh,
469 faceList&& faces,
470 labelList&& allOwner,
471 labelList&& allNeighbour,
472 const bool syncPar
473)
474:
476 (
477 io,
478 std::move(points),
479 std::move(faces),
480 std::move(allOwner),
481 std::move(allNeighbour),
482 syncPar
483 ),
485 (
486 static_cast<const objectRegistry&>(*this),
487 static_cast<const fvSchemes&>(baseMesh)
488 ),
491 (
492 static_cast<const objectRegistry&>(*this),
493 static_cast<const fvSolution&>(baseMesh)
494 ),
495 data
496 (
497 static_cast<const objectRegistry&>(*this),
498 static_cast<const data&>(baseMesh)
499 ),
500 boundary_(*this),
501 lduPtr_(nullptr),
502 curTimeIndex_(time().timeIndex()),
503 VPtr_(nullptr),
504 V0Ptr_(nullptr),
505 V00Ptr_(nullptr),
506 SfPtr_(nullptr),
507 magSfPtr_(nullptr),
508 CPtr_(nullptr),
509 CfPtr_(nullptr),
510 phiPtr_(nullptr)
511{
512 DebugInFunction << "Constructing fvMesh as copy and primitives" << endl;
513}
514
515
517(
518 const IOobject& io,
519 const fvMesh& baseMesh,
521 faceList&& faces,
522 cellList&& cells,
523 const bool syncPar
524)
525:
527 (
528 io,
529 std::move(points),
530 std::move(faces),
531 std::move(cells),
532 syncPar
533 ),
535 (
536 static_cast<const objectRegistry&>(*this),
537 static_cast<const fvSchemes&>(baseMesh)
538 ),
541 (
542 static_cast<const objectRegistry&>(*this),
543 static_cast<const fvSolution&>(baseMesh)
544 ),
545 data
546 (
547 static_cast<const objectRegistry&>(*this),
548 static_cast<const data&>(baseMesh)
549 ),
550 boundary_(*this),
551 lduPtr_(nullptr),
552 curTimeIndex_(time().timeIndex()),
553 VPtr_(nullptr),
554 V0Ptr_(nullptr),
555 V00Ptr_(nullptr),
556 SfPtr_(nullptr),
557 magSfPtr_(nullptr),
558 CPtr_(nullptr),
559 CfPtr_(nullptr),
560 phiPtr_(nullptr)
561{
562 DebugInFunction << "Constructing fvMesh as copy and primitives" << endl;
563}
564
565
566// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
567
569{
570 clearOut();
571}
572
573
574// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
575
577(
579 const dictionary& dict
580) const
581{
582 // Redirect to fvMatrix solver
584}
585
586
588(
590 const dictionary& dict
591) const
592{
593 // Redirect to fvMatrix solver
595}
596
597
599(
601 const dictionary& dict
602) const
603{
604 // Redirect to fvMatrix solver
606}
607
608
610(
612 const dictionary& dict
613) const
614{
615 // Redirect to fvMatrix solver
617}
618
619
621(
623 const dictionary& dict
624) const
625{
626 // Redirect to fvMatrix solver
628}
629
630
632(
633 polyPatchList& plist,
634 const bool validBoundary
635)
636{
637 if (boundary().size())
638 {
640 << " boundary already exists"
641 << abort(FatalError);
642 }
643
644 addPatches(plist, validBoundary);
645 boundary_.addPatches(boundaryMesh());
646}
647
648
650(
651 const List<polyPatch*>& p,
652 const bool validBoundary
653)
654{
655 // Acquire ownership of the pointers
656 polyPatchList plist(const_cast<List<polyPatch*>&>(p));
657
658 addFvPatches(plist, validBoundary);
659}
660
661
663{
664 DebugInFunction << "Removing boundary patches." << endl;
665
666 // Remove fvBoundaryMesh data first.
667 boundary_.clear();
668 boundary_.setSize(0);
670
671 clearOut();
672}
673
674
676{
677 DebugInFunction << "Updating fvMesh";
678
680
681 if (state == polyMesh::TOPO_PATCH_CHANGE)
682 {
683 DebugInfo << "Boundary and topological update" << endl;
684
685 boundary_.readUpdate(boundaryMesh());
686
687 clearOut();
688
689 }
690 else if (state == polyMesh::TOPO_CHANGE)
691 {
692 DebugInfo << "Topological update" << endl;
693
694 // fvMesh::clearOut() but without the polyMesh::clearOut
695 clearOutLocal();
696 }
697 else if (state == polyMesh::POINTS_MOVED)
698 {
699 DebugInfo << "Point motion update" << endl;
700
701 clearGeom();
702 }
703 else
704 {
705 DebugInfo << "No update" << endl;
706 }
707
708 return state;
709}
710
711
713{
714 return boundary_;
715}
716
717
719{
720 if (!lduPtr_)
721 {
723 << "Calculating fvMeshLduAddressing from nFaces:"
724 << nFaces() << endl;
725
726 lduPtr_ = new fvMeshLduAddressing(*this);
727
728 return *lduPtr_;
729 }
730
731 return *lduPtr_;
732}
733
734
736{
737 return boundary().interfaces();
738}
739
740
742{
744 << " nOldCells:" << meshMap.nOldCells()
745 << " nCells:" << nCells()
746 << " nOldFaces:" << meshMap.nOldFaces()
747 << " nFaces:" << nFaces()
748 << endl;
749
750 // We require geometric properties valid for the old mesh
751 if
752 (
753 meshMap.cellMap().size() != nCells()
754 || meshMap.faceMap().size() != nFaces()
755 )
756 {
758 << "mapPolyMesh does not correspond to the old mesh."
759 << " nCells:" << nCells()
760 << " cellMap:" << meshMap.cellMap().size()
761 << " nOldCells:" << meshMap.nOldCells()
762 << " nFaces:" << nFaces()
763 << " faceMap:" << meshMap.faceMap().size()
764 << " nOldFaces:" << meshMap.nOldFaces()
765 << exit(FatalError);
766 }
767
768 // Create a mapper
769 const fvMeshMapper mapper(*this, meshMap);
770
771 // Map all the volFields in the objectRegistry
772 MapGeometricFields<scalar, fvPatchField, fvMeshMapper, volMesh>
773 (mapper);
774 MapGeometricFields<vector, fvPatchField, fvMeshMapper, volMesh>
775 (mapper);
776 MapGeometricFields<sphericalTensor, fvPatchField, fvMeshMapper, volMesh>
777 (mapper);
778 MapGeometricFields<symmTensor, fvPatchField, fvMeshMapper, volMesh>
779 (mapper);
780 MapGeometricFields<tensor, fvPatchField, fvMeshMapper, volMesh>
781 (mapper);
782
783 // Map all the surfaceFields in the objectRegistry
784 MapGeometricFields<scalar, fvsPatchField, fvMeshMapper, surfaceMesh>
785 (mapper);
786 MapGeometricFields<vector, fvsPatchField, fvMeshMapper, surfaceMesh>
787 (mapper);
789 <
791 >
792 (mapper);
793 MapGeometricFields<symmTensor, fvsPatchField, fvMeshMapper, surfaceMesh>
794 (mapper);
795 MapGeometricFields<tensor, fvsPatchField, fvMeshMapper, surfaceMesh>
796 (mapper);
797
798 // Map all the dimensionedFields in the objectRegistry
799 MapDimensionedFields<scalar, fvMeshMapper, volMesh>(mapper);
800 MapDimensionedFields<vector, fvMeshMapper, volMesh>(mapper);
801 MapDimensionedFields<sphericalTensor, fvMeshMapper, volMesh>(mapper);
802 MapDimensionedFields<symmTensor, fvMeshMapper, volMesh>(mapper);
803 MapDimensionedFields<tensor, fvMeshMapper, volMesh>(mapper);
804
805 // Map all the clouds in the objectRegistry
806 mapClouds(*this, meshMap);
807
808
809 const labelList& cellMap = meshMap.cellMap();
810
811 // Map the old volume. Just map to new cell labels.
812 if (V0Ptr_)
813 {
814 scalarField& V0 = *V0Ptr_;
815
816 scalarField savedV0(V0);
817 V0.setSize(nCells());
818
819 forAll(V0, i)
820 {
821 if (cellMap[i] > -1)
822 {
823 V0[i] = savedV0[cellMap[i]];
824 }
825 else
826 {
827 V0[i] = 0.0;
828 }
829 }
830
831 // Inject volume of merged cells
832 label nMerged = 0;
833 forAll(meshMap.reverseCellMap(), oldCelli)
834 {
835 label index = meshMap.reverseCellMap()[oldCelli];
836
837 if (index < -1)
838 {
839 label celli = -index-2;
840
841 V0[celli] += savedV0[oldCelli];
842
843 nMerged++;
844 }
845 }
846
848 << "Mapping old time volume V0. Merged "
849 << nMerged << " out of " << nCells() << " cells" << endl;
850 }
851
852
853 // Map the old-old volume. Just map to new cell labels.
854 if (V00Ptr_)
855 {
856 scalarField& V00 = *V00Ptr_;
857
858 scalarField savedV00(V00);
859 V00.setSize(nCells());
860
861 forAll(V00, i)
862 {
863 if (cellMap[i] > -1)
864 {
865 V00[i] = savedV00[cellMap[i]];
866 }
867 else
868 {
869 V00[i] = 0.0;
870 }
871 }
872
873 // Inject volume of merged cells
874 label nMerged = 0;
875 forAll(meshMap.reverseCellMap(), oldCelli)
876 {
877 label index = meshMap.reverseCellMap()[oldCelli];
878
879 if (index < -1)
880 {
881 label celli = -index-2;
882
883 V00[celli] += savedV00[oldCelli];
884 nMerged++;
885 }
886 }
887
889 << "Mapping old time volume V00. Merged "
890 << nMerged << " out of " << nCells() << " cells" << endl;
891 }
892}
893
894
896{
898
899 // Grab old time volumes if the time has been incremented
900 // This will update V0, V00
901 if (curTimeIndex_ < time().timeIndex())
902 {
903 storeOldVol(V());
904 }
905
906
907 // Move the polyMesh and initialise the mesh motion fluxes field
908 // Note: mesh flux updated by the fvGeometryScheme
909
910 if (!phiPtr_)
911 {
912 DebugInFunction<< "Creating initial meshPhi field" << endl;
913
914 // Create mesh motion flux
915 phiPtr_ = new surfaceScalarField
916 (
918 (
919 "meshPhi",
920 this->time().timeName(),
921 *this,
924 false
925 ),
926 *this,
928 );
929 }
930 else
931 {
932 // Grab old time mesh motion fluxes if the time has been incremented
933 if (phiPtr_->timeIndex() != time().timeIndex())
934 {
935 DebugInFunction<< "Accessing old-time meshPhi field" << endl;
936 phiPtr_->oldTime();
937 }
938 }
939
941
942 // Update or delete the local geometric properties as early as possible so
943 // they can be used if necessary. These get recreated here instead of
944 // demand driven since they might do parallel transfers which can conflict
945 // with when they're actually being used.
946 // Note that between above "polyMesh::movePoints(p)" and here nothing
947 // should use the local geometric properties.
948 updateGeomNotOldVol();
949
950 // Update other local data
951 boundary_.movePoints();
952
953 // Clear weights, deltaCoeffs, nonOrthoDeltaCoeffs, nonOrthCorrectionVectors
955
956 meshObject::movePoints<fvMesh>(*this);
957 meshObject::movePoints<lduMesh>(*this);
958}
959
960
962{
964
965 // Let surfaceInterpolation handle geometry calculation. Note: this does
966 // lower levels updateGeom
968}
969
970
972{
974
975 // Update polyMesh. This needs to keep volume existent!
977
978 // Our slice of the addressing is no longer valid
979 deleteDemandDrivenData(lduPtr_);
980
981 if (VPtr_)
982 {
983 // Grab old time volumes if the time has been incremented
984 // This will update V0, V00
985 storeOldVol(mpm.oldCellVolumes());
986
987 // Few checks
988 if (VPtr_ && (VPtr_->size() != mpm.nOldCells()))
989 {
991 << "V:" << VPtr_->size()
992 << " not equal to the number of old cells "
993 << mpm.nOldCells()
994 << exit(FatalError);
995 }
996 if (V0Ptr_ && (V0Ptr_->size() != mpm.nOldCells()))
997 {
999 << "V0:" << V0Ptr_->size()
1000 << " not equal to the number of old cells "
1001 << mpm.nOldCells()
1002 << exit(FatalError);
1003 }
1004 if (V00Ptr_ && (V00Ptr_->size() != mpm.nOldCells()))
1005 {
1007 << "V0:" << V00Ptr_->size()
1008 << " not equal to the number of old cells "
1009 << mpm.nOldCells()
1010 << exit(FatalError);
1011 }
1012 }
1013
1014
1015 // Clear mesh motion flux (note: could instead save & map like volumes)
1016 if (phiPtr_)
1017 {
1018 // Mesh moving and topology change. Recreate meshPhi
1019 deleteDemandDrivenData(phiPtr_);
1020
1021 // Create mesh motion flux
1022 phiPtr_ = new surfaceScalarField
1023 (
1024 IOobject
1025 (
1026 "meshPhi",
1027 this->time().timeName(),
1028 *this,
1031 false
1032 ),
1033 *this,
1035 );
1036 }
1037
1038 // Clear the sliced fields
1039 clearGeomNotOldVol();
1040
1041 // Map all fields
1042 mapFields(mpm);
1043
1044 // Clear the current volume and other geometry factors
1046
1047 // Clear any non-updateable addressing
1048 clearAddressing(true);
1049
1050 meshObject::updateMesh<fvMesh>(*this, mpm);
1051 meshObject::updateMesh<lduMesh>(*this, mpm);
1052}
1053
1054
1056(
1057 IOstreamOption streamOpt,
1058 const bool valid
1059) const
1060{
1061 bool ok = true;
1062 if (phiPtr_)
1063 {
1064 ok = phiPtr_->write(valid);
1065 // NOTE: The old old time mesh phi might be necessary for certain
1066 // solver smooth restart using second order time schemes.
1067 //ok = phiPtr_->oldTime().write();
1068 }
1069 if (V0Ptr_ && V0Ptr_->writeOpt() == IOobject::AUTO_WRITE)
1070 {
1071 // For second order restarts we need to write V0
1072 ok = V0Ptr_->write(valid);
1073 }
1074
1075 return ok && polyMesh::writeObject(streamOpt, valid);
1076}
1077
1078
1079bool Foam::fvMesh::write(const bool valid) const
1080{
1081 return polyMesh::write(valid);
1082}
1083
1084
1085template<>
1086typename Foam::pTraits<Foam::sphericalTensor>::labelType
1088{
1089 return Foam::pTraits<Foam::sphericalTensor>::labelType(1);
1090}
1091
1092
1093// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1094
1095bool Foam::fvMesh::operator!=(const fvMesh& rhs) const
1096{
1097 return &rhs != this;
1098}
1099
1100
1101bool Foam::fvMesh::operator==(const fvMesh& rhs) const
1102{
1103 return &rhs == this;
1104}
1105
1106
1107// ************************************************************************* //
Graphite solid properties.
Definition: C.H:53
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
The IOstreamOption is a simple container for options an IOstream can normally have.
void setSize(const label n)
Alias for resize()
Definition: List.H:218
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:63
Database for solution data, solver performance and other reduced data.
Definition: data.H:58
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write using stream options.
virtual bool write()
Write the output fields.
Maps input fields from local mesh to secondary mesh at runtime.
Definition: mapFields.H:225
Foam::fvBoundaryMesh.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
SolverPerformance< Type > solveSegregatedOrCoupled(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:62
Foam::fvMeshLduAddressing.
Class holds all the necessary information for mapping fields associated with fvMesh.
Definition: fvMeshMapper.H:58
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
virtual bool movePoints()
Avoid masking surfaceInterpolation method.
slicedSurfaceVectorField * SfPtr_
Face area vectors.
Definition: fvMesh.H:119
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:712
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: fvMesh.C:735
void addFvPatches(polyPatchList &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:632
slicedSurfaceVectorField * CfPtr_
Face centres.
Definition: fvMesh.H:128
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:718
void clearGeom()
Clear local geometry.
Definition: fvMesh.C:116
SlicedDimensionedField< scalar, volMesh > * VPtr_
Cell volumes.
Definition: fvMesh.H:110
void clearOutLocal()
Clear local-only storage (geometry, addressing etc)
Definition: fvMesh.C:227
virtual ~fvMesh()
Destructor.
Definition: fvMesh.C:568
void updateGeomNotOldVol()
Clear geometry like clearGeomNotOldVol but recreate any.
Definition: fvMesh.C:78
void clearGeomNotOldVol()
Clear geometry but not the old-time cell volumes.
Definition: fvMesh.C:54
void storeOldVol(const scalarField &)
Preserve old volume(s)
Definition: fvMesh.C:164
void clearAddressing(const bool isMeshUpdate=false)
Clear local addressing.
Definition: fvMesh.C:128
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:675
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write the underlying polyMesh and other data.
Definition: fvMesh.C:1056
void removeFvBoundary()
Definition: fvMesh.C:662
slicedVolVectorField * CPtr_
Cell centres.
Definition: fvMesh.H:125
surfaceScalarField * magSfPtr_
Mag face area vectors.
Definition: fvMesh.H:122
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:239
virtual void updateGeom()
Definition: fvMesh.C:961
Selector class for finite volume differencing schemes. fvMesh is derived from fvSchemes so that all f...
Definition: fvSchemes.H:60
Selector class for finite volume solution solution. fvMesh is derived from fvSolution so that all fie...
Definition: fvSolution.H:60
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:79
The class contains the addressing required by the lduMatrix: upper, lower and losort.
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:63
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
label nOldCells() const
Number of old cells.
Definition: mapPolyMesh.H:387
const labelList & faceMap() const
Old face map.
Definition: mapPolyMesh.H:410
const labelList & reverseCellMap() const
Reverse cell map.
Definition: mapPolyMesh.H:532
const labelList & cellMap() const
Old cell map.
Definition: mapPolyMesh.H:435
label nOldFaces() const
Number of old faces.
Definition: mapPolyMesh.H:381
const scalarField & oldCellVolumes() const
Definition: mapPolyMesh.H:651
static void clearUpto(objectRegistry &obr)
Definition: MeshObject.C:217
void movePoints()
Update for new mesh geometry.
void updateMesh()
Update for new mesh topology.
Registry of regIOobjects.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:91
@ TOPO_PATCH_CHANGE
Definition: polyMesh.H:95
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in.
Definition: polyMeshIO.C:75
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:39
void clearOut()
Clear all geometry and addressing unnecessary for CFD.
friend bool operator!=(const refineCell &rc1, const refineCell &rc2)
Definition: refineCell.H:106
friend bool operator==(const refineCell &rc1, const refineCell &rc2)
Definition: refineCell.H:97
Cell to surface interpolation scheme. Included in fvMesh.
void clearOut()
Clear all geometry and addressing.
virtual void updateGeom()
Update all geometric data.
Mesh data needed to do the Finite Volume discretisation.
Definition: surfaceMesh.H:52
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:63
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
faceListList boundary
Template functions to aid in the implementation of demand driven data.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
Foam::pTraits< Foam::sphericalTensor >::labelType Foam::fvMesh::validComponents< Foam::sphericalTensor >() const
Definition: fvMesh.C:1087
const pointField & points
const cellShapeList & cells
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
Generic Geometric field mapper. For "real" mapping, add template specialisations for mapping of inter...
word timeName
Definition: getTimeIndex.H:3
#define DebugInfo
Report an information message using Foam::Info.
#define DebugInFunction
Report an information message using Foam::Info.
#define InfoInFunction
Report an information message using Foam::Info.
Namespace for OpenFOAM.
const fileOperation & fileHandler()
Get current file handler.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void mapClouds(const objectRegistry &db, const mapPolyMesh &mapper)
Generic Geometric field mapper.
Definition: mapClouds.H:51
SphericalTensor< scalar > sphericalTensor
SphericalTensor of scalars, i.e. SphericalTensor<scalar>.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
errorManip< error > abort(error &err)
Definition: errorManip.H:144
void MapGeometricFields(const MeshMapper &mapper)
Generic Geometric field mapper.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
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:666
error FatalError
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
void deleteDemandDrivenData(DataPtr &dataPtr)
label timeIndex
Definition: getTimeIndex.H:30
dictionary dict
CEqn solve()
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Foam::surfaceFields.