Forward compatibility: flex
[foam-extend-3.2.git] / src / dynamicMesh / topoChangerFvMesh / movingBodyTopoFvMesh / movingBodyTopoFvMesh.C
blob18bbaa85cedce11398e86a53bbd207cf72de401f
1 /*---------------------------------------------------------------------------*\
2   =========                 |
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 -------------------------------------------------------------------------------
8 License
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"
27 #include "foamTime.H"
28 #include "mapPolyMesh.H"
29 #include "layerAdditionRemoval.H"
30 #include "volMesh.H"
31 #include "transformField.H"
32 #include "addToRunTimeSelectionTable.H"
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 namespace Foam
38     defineTypeNameAndDebug(movingBodyTopoFvMesh, 0);
40     addToRunTimeSelectionTable
41     (
42         topoChangerFvMesh,
43         movingBodyTopoFvMesh,
44         IOobject
45     );
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())
70     {
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)
78         {
79             const labelList& curCp = cp[movingCells[cellI]];
81             forAll (curCp, pointI)
82             {
83                 vertexMarkup[curCp[pointI]] = 1;
84             }
85         }
86     }
88     faceZoneID frontFacesID(frontFacesName_, faceZones());
90     if (frontFacesID.active())
91     {
92         const faceZone& frontFaces = faceZones()[frontFacesID.index()];
94         const labelList& mp = frontFaces().meshPoints();
96         forAll (mp, mpI)
97         {
98             vertexMarkup[mp[mpI]] = 1;
99         }
100     }
102     faceZoneID backFacesID(backFacesName_, faceZones());
104     if (backFacesID.active())
105     {
106         const faceZone& backFaces = faceZones()[backFacesID.index()];
108         const labelList& mp = backFaces().meshPoints();
110         forAll (mp, mpI)
111         {
112             vertexMarkup[mp[mpI]] = 1;
113         }
114     }
116     return tvertexMarkup;
120 void Foam::movingBodyTopoFvMesh::addZonesAndModifiers()
122     // Add zones and modifiers for motion action
124     if (topoChanger_.size() > 0)
125     {
126         Info<< "void movingBodyTopoFvMesh::addZonesAndModifiers() : "
127             << "Zones and modifiers already present.  Skipping."
128             << endl;
130         return;
131     }
133     // Add layer addition/removal interfaces
134     topoChanger_.setSize(2);
135     label nMods = 0;
138     faceZoneID frontFacesID(frontFacesName_, faceZones());
139     faceZoneID backFacesID(backFacesName_, faceZones());
141     if (frontFacesID.active())
142     {
143         const faceZone& frontFaces = faceZones()[frontFacesID.index()];
145         if (!frontFaces.empty())
146         {
147             topoChanger_.set
148             (
149                 nMods,
150                 new layerAdditionRemoval
151                 (
152                     frontFacesName_ + "Layer",
153                     nMods,
154                     topoChanger_,
155                     frontFacesName_,
156                     readScalar
157                     (
158                         dict_.subDict("front").lookup("minThickness")
159                     ),
160                     readScalar
161                     (
162                         dict_.subDict("front").lookup("maxThickness")
163                     )
164                 )
165             );
167             nMods++;
168         }
169     }
171     if (backFacesID.active())
172     {
173         const faceZone& backFaces = faceZones()[backFacesID.index()];
175         if (!backFaces.empty())
176         {
177             topoChanger_.set
178             (
179                 nMods,
180                 new layerAdditionRemoval
181                 (
182                     backFacesName_ + "Layer",
183                     nMods,
184                     topoChanger_,
185                     backFacesName_,
186                     readScalar
187                     (
188                         dict_.subDict("back").lookup("minThickness")
189                     ),
190                     readScalar
191                     (
192                         dict_.subDict("back").lookup("maxThickness")
193                     )
194                 )
195             );
197             nMods++;
198         }
199     }
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.
211     // HJ, 18/Feb/2011
212 //     write();
216 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
218 // Construct from components
219 Foam::movingBodyTopoFvMesh::movingBodyTopoFvMesh(const IOobject& io)
221     topoChangerFvMesh(io),
222     dict_
223     (
224         IOdictionary
225         (
226             IOobject
227             (
228                 "dynamicMeshDict",
229                 time().constant(),
230                 *this,
231                 IOobject::MUST_READ,
232                 IOobject::NO_WRITE
233             )
234         ).subDict(typeName + "Coeffs")
235     ),
236     movingCellsName_(dict_.lookup("movingCells")),
237     frontFacesName_(dict_.lookup("frontFaces")),
238     backFacesName_(dict_.lookup("backFaces")),
239     SBMFPtr_(solidBodyMotionFunction::New(dict_, time())),
240     motionMask_()
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)
268     {
269         Pout<< "Topology change. Calculating motion point mask" << endl;
270         motionMask_ = calcMotionMask();
271     }
273     if (localMeshChanged)
274     {
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);
282 //         resetMotion();
283 //         setV0();
285         // Get new points from preMotion
286         newPoints = topoChangeMap().preMotionPoints();
287     }
288 //     else
289 //     {
290 //         // No change, use old points
291 //         // Note: using setOldPoints instead of movePoints.
292 //         // HJ, 23/Aug/2015
293 //         setOldPoints(oldPointsNew);
294 //         resetMotion();
295 //         setV0();
296 //     }
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 // ************************************************************************* //