1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 option) any later version.
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
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "movingBodyTopoFvMesh.H"
29 #include "mapPolyMesh.H"
30 #include "layerAdditionRemoval.H"
32 #include "transformField.H"
33 #include "addToRunTimeSelectionTable.H"
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 defineTypeNameAndDebug(movingBodyTopoFvMesh, 0);
41 addToRunTimeSelectionTable
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 Foam::tmp<Foam::scalarField>
53 Foam::movingBodyTopoFvMesh::calcMotionMask() const
55 Info<< "Updating vertex markup" << endl;
57 tmp<scalarField> tvertexMarkup(new scalarField(allPoints().size(), 0));
58 scalarField& vertexMarkup = tvertexMarkup();
60 cellZoneID movingCellsID(movingCellsName_, cellZones());
62 // In order to do a correct update on a mask on processor boundaries,
63 // Detection of moving cells should use patchNeighbourField for
64 // processor (not coupled!) boundaries. This is done by expanding
65 // a moving cell set into a field and making sure that processor patch
66 // points move in sync. Not done at the moment, probably best to do
67 // using parallel update of pointFields. HJ, 19/Feb/2011
69 // If moving cells are found, perform mark-up
70 if (movingCellsID.active())
72 // Get cell-point addressing
73 const labelListList& cp = cellPoints();
75 // Get labels of all moving cells
76 const labelList& movingCells = cellZones()[movingCellsID.index()];
78 forAll (movingCells, cellI)
80 const labelList& curCp = cp[movingCells[cellI]];
82 forAll (curCp, pointI)
84 vertexMarkup[curCp[pointI]] = 1;
89 faceZoneID frontFacesID(frontFacesName_, faceZones());
91 if (frontFacesID.active())
93 const faceZone& frontFaces = faceZones()[frontFacesID.index()];
95 const labelList& mp = frontFaces().meshPoints();
99 vertexMarkup[mp[mpI]] = 1;
103 faceZoneID backFacesID(backFacesName_, faceZones());
105 if (backFacesID.active())
107 const faceZone& backFaces = faceZones()[backFacesID.index()];
109 const labelList& mp = backFaces().meshPoints();
113 vertexMarkup[mp[mpI]] = 1;
117 return tvertexMarkup;
121 void Foam::movingBodyTopoFvMesh::addZonesAndModifiers()
123 // Add zones and modifiers for motion action
125 if (topoChanger_.size() > 0)
127 Info<< "void movingBodyTopoFvMesh::addZonesAndModifiers() : "
128 << "Zones and modifiers already present. Skipping."
134 // Add layer addition/removal interfaces
135 topoChanger_.setSize(2);
139 faceZoneID frontFacesID(frontFacesName_, faceZones());
140 faceZoneID backFacesID(backFacesName_, faceZones());
142 if (frontFacesID.active())
144 const faceZone& frontFaces = faceZones()[frontFacesID.index()];
146 if (!frontFaces.empty())
151 new layerAdditionRemoval
153 frontFacesName_ + "Layer",
159 dict_.subDict("front").lookup("minThickness")
163 dict_.subDict("front").lookup("maxThickness")
172 if (backFacesID.active())
174 const faceZone& backFaces = faceZones()[backFacesID.index()];
176 if (!backFaces.empty())
181 new layerAdditionRemoval
183 backFacesName_ + "Layer",
189 dict_.subDict("back").lookup("minThickness")
193 dict_.subDict("back").lookup("maxThickness")
202 topoChanger_.setSize(nMods);
204 reduce(nMods, sumOp<label>());
206 Info << "Adding " << nMods << " mesh modifiers" << endl;
208 // Write mesh and modifiers
209 topoChanger_.write();
211 // No need to write the mesh - only modifiers are added.
217 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
219 // Construct from components
220 Foam::movingBodyTopoFvMesh::movingBodyTopoFvMesh(const IOobject& io)
222 topoChangerFvMesh(io),
235 ).subDict(typeName + "Coeffs")
237 movingCellsName_(dict_.lookup("movingCells")),
238 frontFacesName_(dict_.lookup("frontFaces")),
239 backFacesName_(dict_.lookup("backFaces")),
240 SBMFPtr_(solidBodyMotionFunction::New(dict_, time())),
243 addZonesAndModifiers();
244 motionMask_ = calcMotionMask();
248 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
250 Foam::movingBodyTopoFvMesh::~movingBodyTopoFvMesh()
254 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
256 bool Foam::movingBodyTopoFvMesh::update()
258 // Store points to recreate mesh motion
259 pointField oldPointsNew = allPoints();
260 pointField newPoints = allPoints();
262 autoPtr<mapPolyMesh> topoChangeMap = topoChanger_.changeMesh();
264 bool localMeshChanged = topoChangeMap->morphing();
265 bool globalMeshChanged = localMeshChanged;
266 reduce(globalMeshChanged, orOp<bool>());
268 if (globalMeshChanged)
270 Pout<< "Topology change. Calculating motion point mask" << endl;
271 motionMask_ = calcMotionMask();
274 if (localMeshChanged)
276 // // Map old points onto the new mesh
277 // pointField mappedOldPointsNew(allPoints().size());
278 // mappedOldPointsNew.map(oldPointsNew, topoChangeMap->pointMap());
280 // movePoints(mappedOldPointsNew);
284 // Get new points from preMotion
285 newPoints = topoChangeMap().preMotionPoints();
289 // // No change, use old points
290 // movePoints(oldPointsNew);
295 // Calculate new points using a velocity transformation
296 newPoints += motionMask_*
297 transform(SBMFPtr_().velocity(), newPoints)*time().deltaT().value();
299 Info << "Executing mesh motion" << endl;
300 movePoints(newPoints);
302 return localMeshChanged;
306 // ************************************************************************* //