1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 3.2
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
9 This file is part of foam-extend.
11 foam-extend is free software: you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by the
13 Free Software Foundation, either version 3 of the License, or (at your
14 option) any later version.
16 foam-extend is distributed in the hope that it will be useful, but
17 WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 General Public License for more details.
21 You should have received a copy of the GNU General Public License
22 along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "movingBodyTopoFvMesh.H"
28 #include "mapPolyMesh.H"
29 #include "layerAdditionRemoval.H"
31 #include "transformField.H"
32 #include "addToRunTimeSelectionTable.H"
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 defineTypeNameAndDebug(movingBodyTopoFvMesh, 0);
40 addToRunTimeSelectionTable
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 Foam::tmp<Foam::scalarField>
52 Foam::movingBodyTopoFvMesh::calcMotionMask() const
54 Info<< "Updating vertex markup" << endl;
56 tmp<scalarField> tvertexMarkup(new scalarField(allPoints().size(), 0));
57 scalarField& vertexMarkup = tvertexMarkup();
59 cellZoneID movingCellsID(movingCellsName_, cellZones());
61 // In order to do a correct update on a mask on processor boundaries,
62 // Detection of moving cells should use patchNeighbourField for
63 // processor (not coupled!) boundaries. This is done by expanding
64 // a moving cell set into a field and making sure that processor patch
65 // points move in sync. Not done at the moment, probably best to do
66 // using parallel update of pointFields. HJ, 19/Feb/2011
68 // If moving cells are found, perform mark-up
69 if (movingCellsID.active())
71 // Get cell-point addressing
72 const labelListList& cp = cellPoints();
74 // Get labels of all moving cells
75 const labelList& movingCells = cellZones()[movingCellsID.index()];
77 forAll (movingCells, cellI)
79 const labelList& curCp = cp[movingCells[cellI]];
81 forAll (curCp, pointI)
83 vertexMarkup[curCp[pointI]] = 1;
88 faceZoneID frontFacesID(frontFacesName_, faceZones());
90 if (frontFacesID.active())
92 const faceZone& frontFaces = faceZones()[frontFacesID.index()];
94 const labelList& mp = frontFaces().meshPoints();
98 vertexMarkup[mp[mpI]] = 1;
102 faceZoneID backFacesID(backFacesName_, faceZones());
104 if (backFacesID.active())
106 const faceZone& backFaces = faceZones()[backFacesID.index()];
108 const labelList& mp = backFaces().meshPoints();
112 vertexMarkup[mp[mpI]] = 1;
116 return tvertexMarkup;
120 void Foam::movingBodyTopoFvMesh::addZonesAndModifiers()
122 // Add zones and modifiers for motion action
124 if (topoChanger_.size() > 0)
126 Info<< "void movingBodyTopoFvMesh::addZonesAndModifiers() : "
127 << "Zones and modifiers already present. Skipping."
133 // Add layer addition/removal interfaces
134 topoChanger_.setSize(2);
138 faceZoneID frontFacesID(frontFacesName_, faceZones());
139 faceZoneID backFacesID(backFacesName_, faceZones());
141 if (frontFacesID.active())
143 const faceZone& frontFaces = faceZones()[frontFacesID.index()];
145 if (!frontFaces.empty())
150 new layerAdditionRemoval
152 frontFacesName_ + "Layer",
158 dict_.subDict("front").lookup("minThickness")
162 dict_.subDict("front").lookup("maxThickness")
171 if (backFacesID.active())
173 const faceZone& backFaces = faceZones()[backFacesID.index()];
175 if (!backFaces.empty())
180 new layerAdditionRemoval
182 backFacesName_ + "Layer",
188 dict_.subDict("back").lookup("minThickness")
192 dict_.subDict("back").lookup("maxThickness")
201 topoChanger_.setSize(nMods);
203 reduce(nMods, sumOp<label>());
205 Info << "Adding " << nMods << " mesh modifiers" << endl;
207 // Write mesh and modifiers
208 topoChanger_.write();
210 // No need to write the mesh - only modifiers are added.
216 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
218 // Construct from components
219 Foam::movingBodyTopoFvMesh::movingBodyTopoFvMesh(const IOobject& io)
221 topoChangerFvMesh(io),
234 ).subDict(typeName + "Coeffs")
236 movingCellsName_(dict_.lookup("movingCells")),
237 frontFacesName_(dict_.lookup("frontFaces")),
238 backFacesName_(dict_.lookup("backFaces")),
239 SBMFPtr_(solidBodyMotionFunction::New(dict_, time())),
242 addZonesAndModifiers();
243 motionMask_ = calcMotionMask();
247 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
249 Foam::movingBodyTopoFvMesh::~movingBodyTopoFvMesh()
253 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
255 bool Foam::movingBodyTopoFvMesh::update()
257 // Store points to recreate mesh motion
258 pointField oldPointsNew = allPoints();
259 pointField newPoints = allPoints();
261 autoPtr<mapPolyMesh> topoChangeMap = topoChanger_.changeMesh();
263 bool localMeshChanged = topoChangeMap->morphing();
264 bool globalMeshChanged = localMeshChanged;
265 reduce(globalMeshChanged, orOp<bool>());
267 if (globalMeshChanged)
269 Pout<< "Topology change. Calculating motion point mask" << endl;
270 motionMask_ = calcMotionMask();
273 if (localMeshChanged)
275 // // Map old points onto the new mesh
276 // pointField mappedOldPointsNew(allPoints().size());
277 // mappedOldPointsNew.map(oldPointsNew, topoChangeMap->pointMap());
279 // // Note: using setOldPoints instead of movePoints.
280 // // HJ, 23/Aug/2015
281 // setOldPoints(mappedOldPointsNew);
285 // Get new points from preMotion
286 newPoints = topoChangeMap().preMotionPoints();
290 // // No change, use old points
291 // // Note: using setOldPoints instead of movePoints.
292 // // HJ, 23/Aug/2015
293 // setOldPoints(oldPointsNew);
298 // Calculate new points using a velocity transformation
299 newPoints += motionMask_*
300 transform(SBMFPtr_().velocity(), newPoints)*time().deltaT().value();
302 Info << "Executing mesh motion" << endl;
303 movePoints(newPoints);
305 return localMeshChanged;
309 // ************************************************************************* //