Loading...
Searching...
No Matches
RBFMotionSolver.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 4.0
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7-------------------------------------------------------------------------------
8License
9 This file is part of foam-extend.
10
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.
15
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.
20
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/>.
23
24\*---------------------------------------------------------------------------*/
25
26#include "RBFMotionSolver.H"
27#include "addToRunTimeSelectionTable.H"
28
29// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30
31namespace Foam
32{
34
36(
37 motionSolver,
39 dictionary
40);
41}
42
43// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44
45void Foam::RBFMotionSolver::makeControlIDs()
46{
47 // Points that are neither on moving nor on static patches
48 // will be marked with 0
49 labelList markedPoints(mesh().nPoints(), 0);
50 // Mark all points on moving patches with 1
51 label nMovingPoints = 0;
52 forAll (movingPatches_, patchI)
53 {
54 // Find the patch in boundary
55 label patchIndex =
56 mesh().boundaryMesh().findPatchID(movingPatches_[patchI]);
57
58 if (patchIndex < 0)
59 {
60 FatalErrorIn("void RBFMotionSolver::makeControlIDs()")
61 << "Patch " << movingPatches_[patchI] << " not found. "
62 << "valid patch names: " << mesh().boundaryMesh().names()
63 << abort(FatalError);
64 }
65
66 const labelList& mp = mesh().boundaryMesh()[patchIndex].meshPoints();
67 forAll (mp, i)
68 {
69 markedPoints[mp[i]] = 1;
70 nMovingPoints++;
71 }
72 }
73
74 // Mark moving points and select control points from moving patches
75 movingIDs_.setSize(nMovingPoints);
76 Info << "Total points on moving boundaries: " << nMovingPoints << endl;
77 const pointField& points = mesh().points();
78 // Re-use counter to count moving points
79 // Note: the control points also hold static points in the second part
80 // of the list if static patches are included in the RBF
81 // HJ, 24/Mar/2011
82 nMovingPoints = 0;
83 // Count moving points first
84 forAll (markedPoints, i)
85 {
86 if (markedPoints[i] == 1)
87 {
88 // Grab internal point
89 movingIDs_[nMovingPoints] = i;
90 nMovingPoints++;
91 }
92 }
93
94 movingIDs_.setSize(nMovingPoints);
95 // Actual location of moving points will be set later on request
96 // HJ, 19/Dec/2008
97 movingPoints_.setSize(nMovingPoints, vector::zero);
98 // Mark all points on static patches with -1
99 label nStaticPoints = 0;
100 forAll (staticPatches_, patchI)
101 {
102 // Find the patch in boundary
103 label patchIndex =
104 mesh().boundaryMesh().findPatchID(staticPatches_[patchI]);
105
106 if (patchIndex < 0)
107 {
108 FatalErrorIn("void RBFMotionSolver::makeControlPoints()")
109 << "Patch " << staticPatches_[patchI] << " not found. "
110 << "valid patch names: " << mesh().boundaryMesh().names()
111 << abort(FatalError);
112 }
113
114 const labelList& mp = mesh().boundaryMesh()[patchIndex].meshPoints();
115 forAll (mp, i)
116 {
117 markedPoints[mp[i]] = -1;
118 nStaticPoints++;
119 }
120 }
121
122 Info << "Total points on static boundaries: " << nStaticPoints << endl;
123 staticIDs_.setSize(nStaticPoints);
124 // Re-use counter
125 nStaticPoints = 0;
126 // Count total number of control points
127 forAll (markedPoints, i)
128 {
129 if (markedPoints[i] == -1)
130 {
131 staticIDs_[nStaticPoints] = i;
132 nStaticPoints++;
133 }
134 }
135
136 staticIDs_.setSize(nStaticPoints);
137 // Control IDs also potentially include points on static patches
138 // HJ, 24/Mar/2011
139 controlIDs_.setSize(movingIDs_.size() + staticIDs_.size());
140 motion_.setSize(controlIDs_.size(), vector::zero);
141 label nControlPoints = 0;
142 forAll (movingPatches_, patchI)
143 {
144 // Find the patch in boundary
145 label patchIndex =
146 mesh().boundaryMesh().findPatchID(movingPatches_[patchI]);
147 const labelList& mp = mesh().boundaryMesh()[patchIndex].meshPoints();
148
149 for
150 (
151 label pickedPoint = 0;
152 pickedPoint < mp.size();
153 pickedPoint += coarseningRatio_
154 )
155 {
156 // Pick point as control point
157 controlIDs_[nControlPoints] = mp[pickedPoint];
158 // Mark the point as picked
159 markedPoints[mp[pickedPoint]] = 2;
160 nControlPoints++;
161 }
162 }
163
164 Info << "Selected " << nControlPoints
165 << " control points on moving boundaries" << endl;
166
167 if (includeStaticPatches_)
168 {
169 forAll (staticPatches_, patchI)
170 {
171 // Find the patch in boundary
172 label patchIndex =
173 mesh().boundaryMesh().findPatchID(staticPatches_[patchI]);
174 const labelList& mp =
175 mesh().boundaryMesh()[patchIndex].meshPoints();
176
177 for
178 (
179 label pickedPoint = 0;
180 pickedPoint < mp.size();
181 pickedPoint += coarseningRatio_
182 )
183 {
184 // Pick point as control point
185 controlIDs_[nControlPoints] = mp[pickedPoint];
186 // Mark the point as picked
187 markedPoints[mp[pickedPoint]] = 2;
188 nControlPoints++;
189 }
190 }
191
192 Info << "Selected " << nControlPoints
193 << " total control points" << endl;
194 }
195
196 // Resize control IDs
197 controlIDs_.setSize(nControlPoints);
198 // Pick up point locations
199 controlPoints_.setSize(nControlPoints);
200 // Set control points
201 forAll (controlIDs_, i)
202 {
203 controlPoints_[i] = points[controlIDs_[i]];
204 }
205
206 // Pick up all internal points
207 internalIDs_.setSize(points.size());
208 internalPoints_.setSize(points.size());
209 // Count internal points
210 label nInternalPoints = 0;
211 forAll (markedPoints, i)
212 {
213 if (markedPoints[i] == 0)
214 {
215 // Grab internal point
216 internalIDs_[nInternalPoints] = i;
217 internalPoints_[nInternalPoints] = points[i];
218 nInternalPoints++;
219 }
220 }
221
222 Info << "Number of internal points: " << nInternalPoints << endl;
223 // Resize the lists
224 internalIDs_.setSize(nInternalPoints);
225 internalPoints_.setSize(nInternalPoints);
226}
227
228
229void Foam::RBFMotionSolver::setMovingPoints() const
230{
231 const pointField& points = mesh().points();
232 // Set moving points
233 forAll (movingIDs_, i)
234 {
235 movingPoints_[i] = points[movingIDs_[i]];
236 }
237}
238
239
240// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
241
242Foam::RBFMotionSolver::RBFMotionSolver
243(
244 const polyMesh& mesh,
245 const IOdictionary& dict
246)
247 :
248 displacementMotionSolver(mesh, dict, typeName),
249 movingPatches_(lookup("movingPatches")),
250 staticPatches_(lookup("staticPatches")),
251 coarseningRatio_(readLabel(lookup("coarseningRatio"))),
252 includeStaticPatches_(lookup("includeStaticPatches")),
253 frozenInterpolation_(lookup("frozenInterpolation")),
254 movingIDs_(0),
255 movingPoints_(0),
256 staticIDs_(0),
257 controlIDs_(0),
258 controlPoints_(0),
259 internalIDs_(0),
260 internalPoints_(0),
261 motion_(0),
262 interpolation_
263 (
264 subDict("interpolation"),
265 controlPoints_,
266 internalPoints_
267 )
268{
269 makeControlIDs();
270}
271
272
273// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
274
277
278
279// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
280
281void Foam::RBFMotionSolver::setMotion(const vectorField& m)
282{
283 if (m.size() != movingIDs_.size())
284 {
285 FatalErrorIn
286 (
287 "void RBFMotionSolver::setMotion(const vectorField& m)"
288 ) << "Incorrect size of motion points: m = " << m.size()
289 << " movingIDs = " << movingIDs_.size()
290 << abort(FatalError);
291 }
292
293 // Motion of static points is zero and moving points are first
294 // in the list. HJ, 24/Mar/2011
295 motion_ = vector::zero;
296 forAll (m, i)
297 {
298 motion_[i] = m[i];
299 }
300
301 if (!frozenInterpolation_)
302 {
303 // Set control points
304 const pointField& points = mesh().points();
305 forAll (controlIDs_, i)
306 {
307 controlPoints_[i] = points[controlIDs_[i]];
308 }
309
310 // Re-calculate interpolation
311 interpolation_.movePoints();
312 }
313}
314
315
316const Foam::vectorField& Foam::RBFMotionSolver::movingPoints() const
317{
318 // Update moving points based on current mesh
319 setMovingPoints();
320 return movingPoints_;
321}
322
323void Foam::RBFMotionSolver::movePoints(const pointField&)
324{}
325
326
327Foam::tmp<Foam::pointField> Foam::RBFMotionSolver::curPoints() const
328{
329 // Prepare new points: same as old point
330 tmp<pointField> tcurPoints
331 (
332 new vectorField(mesh().nPoints(), vector::zero)
333 );
334 //fvMesh& mesh = const_cast<fvMesh&>(T.mesh());
335 pointField& curPoints = const_cast<pointField&>(tcurPoints());
336 //pointField& curPoints = tcurPoints();
337 // Add motion to existing points
338 // 1. Insert prescribed motion of moving points
339 forAll (movingIDs_, i)
340 {
341 curPoints[movingIDs_[i]] = motion_[i];
342 }
343
344 // 2. Insert zero motion of static points
345 forAll (staticIDs_, i)
346 {
347 curPoints[staticIDs_[i]] = vector::zero;
348 }
349
350 // Set motion of control
351 vectorField motionOfControl(controlIDs_.size());
352 // 2. Capture positions of control points
353 forAll (controlIDs_, i)
354 {
355 motionOfControl[i] = curPoints[controlIDs_[i]];
356 }
357
358 // Call interpolation
359 vectorField interpolatedMotion =
360 interpolation_.interpolate(motionOfControl).ref();
361 // 3. Insert RBF interpolated motion
362 forAll (internalIDs_, i)
363 {
364 curPoints[internalIDs_[i]] = interpolatedMotion[i];
365 }
366
367 // 4. Add old point positions
368 curPoints += mesh().points();
369 twoDCorrectPoints(const_cast<pointField&>(tcurPoints()));
370 return tcurPoints;
371}
372
373
376
377
378void Foam::RBFMotionSolver::updateMesh(const mapPolyMesh&)
379{
380 // Recalculate control point IDs
381 makeControlIDs();
382}
383
384
385// ************************************************************************* //
forAll(example_CG.gList, solutionI)
Definition CGtest.H:21
Foam::fvMesh & mesh
Definition createMesh.H:47
virtual void updateMesh(const mapPolyMesh &)
const vectorField & movingPoints() const
void setMotion(const vectorField &)
virtual tmp< pointField > curPoints() const
addToRunTimeSelectionTable(Filter, IntegralFilter, dictionary)
defineTypeNameAndDebug(Filter, 0)
label i
Definition pEqn.H:46