Loading...
Searching...
No Matches
polyMeshUpdate.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-2016, 2020 OpenFOAM Foundation
9 Copyright (C) 2020 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 "polyMesh.H"
30#include "mapPolyMesh.H"
31#include "Time.H"
33#include "pointMesh.H"
34#include "indexedOctree.H"
35#include "treeDataCell.H"
36
37// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
38
40{
42 << "Updating addressing and (optional) pointMesh/pointFields" << endl;
43
44 // Update boundaryMesh (note that patches themselves already ok)
45// boundary_.updateMesh(mpm);
46 boundary_.updateMesh();
47
48 // Update zones
49 pointZones_.clearAddressing();
50 faceZones_.clearAddressing();
51 cellZones_.clearAddressing();
52
53 // Remove the stored tet base points
54 tetBasePtIsPtr_.reset(nullptr);
55
56 // Remove the cell tree
57 cellTreePtr_.reset(nullptr);
58
59 // Update parallel data
60 if (globalMeshDataPtr_)
61 {
62 globalMeshDataPtr_->updateMesh();
63 }
64
66
67 // Map the old motion points if present
68 if (oldPointsPtr_)
69 {
70 // Make a copy of the original points
71 pointField oldMotionPoints(*oldPointsPtr_);
72
73 pointField& newMotionPoints = *oldPointsPtr_;
74
75 // Resize the list to new size
76 newMotionPoints.resize(points_.size());
77
78 // Map the list
79 if (mpm.hasMotionPoints())
80 {
81 newMotionPoints.map(oldMotionPoints, mpm.pointMap());
82
83 // Any points created out-of-nothing get set to the current
84 // coordinate for lack of anything better.
86 {
87 if (mpm.pointMap()[newPointi] == -1)
88 {
89 newMotionPoints[newPointi] = points_[newPointi];
90 }
91 }
92 }
93 else
94 {
95 const labelList& pointMap = mpm.pointMap();
96 const labelList& revPointMap = mpm.reversePointMap();
97
98 forAll(pointMap, newPointi)
99 {
100 label oldPointi = pointMap[newPointi];
101 if (oldPointi >= 0)
102 {
103 if (revPointMap[oldPointi] == newPointi) // master point
104 {
105 newMotionPoints[newPointi] = oldMotionPoints[oldPointi];
106 }
107 else
108 {
109 newMotionPoints[newPointi] = points_[newPointi];
110 }
111 }
112 else
113 {
114 newMotionPoints[newPointi] = points_[newPointi];
115 }
116 }
117 }
118 }
119
120 // Map the old motion cell-centres if present
121 if (oldCellCentresPtr_)
122 {
123 // Make a copy of the original cell-centres
124 pointField oldMotionCellCentres(*oldCellCentresPtr_);
125
126 pointField& newMotionCellCentres = *oldCellCentresPtr_;
127
128 // Resize the list to new size
129 newMotionCellCentres.resize(cellCentres().size());
130
131 // Map the list
132 newMotionCellCentres.map(oldMotionCellCentres, mpm.cellMap());
133
134 // Any points created out-of-nothing get set to the current coordinate
135 // for lack of anything better.
136 forAll(mpm.cellMap(), newCelli)
137 {
138 if (mpm.cellMap()[newCelli] == -1)
139 {
140 newMotionCellCentres[newCelli] = cellCentres()[newCelli];
141 }
142 }
143 }
144
147
148 // Reset valid directions (could change by faces put into empty patches)
149 geometricD_ = Zero;
150 solutionD_ = Zero;
151
152 const_cast<Time&>(time()).functionObjects().updateMesh(mpm);
153}
154
155
156// ************************************************************************* //
void map(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 map from the given field
Definition Field.C:302
label size() const noexcept
Definition HashTable.H:358
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
bool hasMotionPoints() const noexcept
Has valid preMotionPoints?
const labelList & cellMap() const noexcept
Old cell map.
const labelList & reversePointMap() const noexcept
Reverse point map.
const labelList & pointMap() const noexcept
Old point map.
static void updateMesh(objectRegistry &obr, const mapPolyMesh &mpm)
Update topology using the given map.
const Time & time() const noexcept
Return time registry.
void setInstance(const fileName &instance, const IOobjectOption::writeOption wOpt=IOobject::AUTO_WRITE)
Set the instance for mesh files.
Definition polyMeshIO.C:29
virtual void updateMesh(const mapPolyMesh &mpm)
Update the mesh corresponding to given map.
const vectorField & cellCentres() const
word timeName
Definition getTimeIndex.H:3
#define DebugInFunction
Report an information message using Foam::Info.
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
List< label > labelList
A List of labels.
Definition List.H:62
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
vectorField pointField
pointField is a vectorField.
label newPointi
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299