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 // Mark moving points and select control points from moving patches
74 movingIDs_.setSize(nMovingPoints);
75 Info << "Total points on moving boundaries: " << nMovingPoints << endl;
76 const pointField& points = mesh().points();
77 // Re-use counter to count moving points
78 // Note: the control points also hold static points in the second part
79 // of the list if static patches are included in the RBF
80 // HJ, 24/Mar/2011
81 nMovingPoints = 0;
82 // Count moving points first
83 forAll (markedPoints, i)
84 {
85 if (markedPoints[i] == 1)
86 {
87 // Grab internal point
88 movingIDs_[nMovingPoints] = i;
89 nMovingPoints++;
90 }
91 }
92 movingIDs_.setSize(nMovingPoints);
93 // Actual location of moving points will be set later on request
94 // HJ, 19/Dec/2008
95 movingPoints_.setSize(nMovingPoints, vector::zero);
96 // Mark all points on static patches with -1
97 label nStaticPoints = 0;
98 forAll (staticPatches_, patchI)
99 {
100 // Find the patch in boundary
101 label patchIndex =
102 mesh().boundaryMesh().findPatchID(staticPatches_[patchI]);
103
104 if (patchIndex < 0)
105 {
106 FatalErrorIn("void RBFMotionSolver::makeControlPoints()")
107 << "Patch " << staticPatches_[patchI] << " not found. "
108 << "valid patch names: " << mesh().boundaryMesh().names()
109 << abort(FatalError);
110 }
111
112 const labelList& mp = mesh().boundaryMesh()[patchIndex].meshPoints();
113 forAll (mp, i)
114 {
115 markedPoints[mp[i]] = -1;
116 nStaticPoints++;
117 }
118 }
119 Info << "Total points on static boundaries: " << nStaticPoints << endl;
120 staticIDs_.setSize(nStaticPoints);
121 // Re-use counter
122 nStaticPoints = 0;
123 // Count total number of control points
124 forAll (markedPoints, i)
125 {
126 if (markedPoints[i] == -1)
127 {
128 staticIDs_[nStaticPoints] = i;
129 nStaticPoints++;
130 }
131 }
132 staticIDs_.setSize(nStaticPoints);
133 // Control IDs also potentially include points on static patches
134 // HJ, 24/Mar/2011
135 controlIDs_.setSize(movingIDs_.size() + staticIDs_.size());
136 motion_.setSize(controlIDs_.size(), vector::zero);
137 label nControlPoints = 0;
138 forAll (movingPatches_, patchI)
139 {
140 // Find the patch in boundary
141 label patchIndex =
142 mesh().boundaryMesh().findPatchID(movingPatches_[patchI]);
143 const labelList& mp = mesh().boundaryMesh()[patchIndex].meshPoints();
144
145 for
146 (
147 label pickedPoint = 0;
148 pickedPoint < mp.size();
149 pickedPoint += coarseningRatio_
150 )
151 {
152 // Pick point as control point
153 controlIDs_[nControlPoints] = mp[pickedPoint];
154 // Mark the point as picked
155 markedPoints[mp[pickedPoint]] = 2;
156 nControlPoints++;
157 }
158 }
159 Info << "Selected " << nControlPoints
160 << " control points on moving boundaries" << endl;
161
162 if (includeStaticPatches_)
163 {
164 forAll (staticPatches_, patchI)
165 {
166 // Find the patch in boundary
167 label patchIndex =
168 mesh().boundaryMesh().findPatchID(staticPatches_[patchI]);
169 const labelList& mp =
170 mesh().boundaryMesh()[patchIndex].meshPoints();
171
172 for
173 (
174 label pickedPoint = 0;
175 pickedPoint < mp.size();
176 pickedPoint += coarseningRatio_
177 )
178 {
179 // Pick point as control point
180 controlIDs_[nControlPoints] = mp[pickedPoint];
181 // Mark the point as picked
182 markedPoints[mp[pickedPoint]] = 2;
183 nControlPoints++;
184 }
185 }
186 Info << "Selected " << nControlPoints
187 << " total control points" << endl;
188 }
189
190 // Resize control IDs
191 controlIDs_.setSize(nControlPoints);
192 // Pick up point locations
193 controlPoints_.setSize(nControlPoints);
194 // Set control points
195 forAll (controlIDs_, i)
196 {
197 controlPoints_[i] = points[controlIDs_[i]];
198 }
199 // Pick up all internal points
200 internalIDs_.setSize(points.size());
201 internalPoints_.setSize(points.size());
202 // Count internal points
203 label nInternalPoints = 0;
204 forAll (markedPoints, i)
205 {
206 if (markedPoints[i] == 0)
207 {
208 // Grab internal point
209 internalIDs_[nInternalPoints] = i;
210 internalPoints_[nInternalPoints] = points[i];
211 nInternalPoints++;
212 }
213 }
214 Info << "Number of internal points: " << nInternalPoints << endl;
215 // Resize the lists
216 internalIDs_.setSize(nInternalPoints);
217 internalPoints_.setSize(nInternalPoints);
218}
219
220
221void Foam::RBFMotionSolver::setMovingPoints() const
222{
223 const pointField& points = mesh().points();
224 // Set moving points
225 forAll (movingIDs_, i)
226 {
227 movingPoints_[i] = points[movingIDs_[i]];
228 }
229}
230
231
232// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
233
234Foam::RBFMotionSolver::RBFMotionSolver
235(
236 const polyMesh& mesh,
237 const IOdictionary& dict
238)
239 :
240 displacementMotionSolver(mesh, dict, typeName),
241 movingPatches_(lookup("movingPatches")),
242 staticPatches_(lookup("staticPatches")),
243 coarseningRatio_(readLabel(lookup("coarseningRatio"))),
244 includeStaticPatches_(lookup("includeStaticPatches")),
245 frozenInterpolation_(lookup("frozenInterpolation")),
246 movingIDs_(0),
247 movingPoints_(0),
248 staticIDs_(0),
249 controlIDs_(0),
250 controlPoints_(0),
251 internalIDs_(0),
252 internalPoints_(0),
253 motion_(0),
254 interpolation_
255 (
256 subDict("interpolation"),
257 controlPoints_,
258 internalPoints_
259 )
260{
261 makeControlIDs();
262}
263
264
265// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
266
269
270
271// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
272
273void Foam::RBFMotionSolver::setMotion(const vectorField& m)
274{
275 if (m.size() != movingIDs_.size())
276 {
277 FatalErrorIn
278 (
279 "void RBFMotionSolver::setMotion(const vectorField& m)"
280 ) << "Incorrect size of motion points: m = " << m.size()
281 << " movingIDs = " << movingIDs_.size()
282 << abort(FatalError);
283 }
284
285 // Motion of static points is zero and moving points are first
286 // in the list. HJ, 24/Mar/2011
287 motion_ = vector::zero;
288 forAll (m, i)
289 {
290 motion_[i] = m[i];
291 }
292
293 if (!frozenInterpolation_)
294 {
295 // Set control points
296 const pointField& points = mesh().points();
297 forAll (controlIDs_, i)
298 {
299 controlPoints_[i] = points[controlIDs_[i]];
300 }
301 // Re-calculate interpolation
302 interpolation_.movePoints();
303 }
304}
305
306
307const Foam::vectorField& Foam::RBFMotionSolver::movingPoints() const
308{
309 // Update moving points based on current mesh
310 setMovingPoints();
311 return movingPoints_;
312}
313
314void Foam::RBFMotionSolver::movePoints(const pointField&)
315{}
316
317
318Foam::tmp<Foam::pointField> Foam::RBFMotionSolver::curPoints() const
319{
320 // Prepare new points: same as old point
321 tmp<pointField> tcurPoints
322 (
323 new vectorField(mesh().nPoints(), vector::zero)
324 );
325 //fvMesh& mesh = const_cast<fvMesh&>(T.mesh());
326 pointField& curPoints = const_cast<pointField&>(tcurPoints());
327 //pointField& curPoints = tcurPoints();
328 // Add motion to existing points
329 // 1. Insert prescribed motion of moving points
330 forAll (movingIDs_, i)
331 {
332 curPoints[movingIDs_[i]] = motion_[i];
333 }
334 // 2. Insert zero motion of static points
335 forAll (staticIDs_, i)
336 {
337 curPoints[staticIDs_[i]] = vector::zero;
338 }
339 // Set motion of control
340 vectorField motionOfControl(controlIDs_.size());
341 // 2. Capture positions of control points
342 forAll (controlIDs_, i)
343 {
344 motionOfControl[i] = curPoints[controlIDs_[i]];
345 }
346 // Call interpolation
347 vectorField interpolatedMotion =
348 interpolation_.interpolate(motionOfControl).ref();
349 // 3. Insert RBF interpolated motion
350 forAll (internalIDs_, i)
351 {
352 curPoints[internalIDs_[i]] = interpolatedMotion[i];
353 }
354 // 4. Add old point positions
355 curPoints += mesh().points();
356 twoDCorrectPoints(const_cast<pointField&>(tcurPoints()));
357 return tcurPoints;
358}
359
360
363
364
365void Foam::RBFMotionSolver::updateMesh(const mapPolyMesh&)
366{
367 // Recalculate control point IDs
368 makeControlIDs();
369}
370
371
372// ************************************************************************* //
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