Nektar++
ForcingMovingBody.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: ForcingMovingBody.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Moving Body m_forcing (movement of a body in a domain is
32// achieved via a m_forcing term which is the results of a coordinate system
33// motion)
34//
35///////////////////////////////////////////////////////////////////////////////
36
39
40using namespace std;
41
42namespace Nektar
43{
46 "MovingBody", ForcingMovingBody::create, "Moving Body Forcing");
47
50 const std::weak_ptr<SolverUtils::EquationSystem> &pEquation)
51 : Forcing(pSession, pEquation)
52{
53}
54
57 const unsigned int &pNumForcingFields, const TiXmlElement *pForce)
58{
59 boost::ignore_unused(pNumForcingFields);
60
61 // Just 3D homogenous 1D problems can use this techinque
62 ASSERTL0(pFields[0]->GetExpType() == MultiRegions::e3DH1D,
63 "Moving body implemented just for 3D Homogenous 1D expansions.");
64
65 // At this point we know in the xml file where those quantities
66 // are declared (equation or file) - via a function name which is now
67 // stored in funcNameD etc. We need now to fill in with this info the
68 // m_zta and m_eta vectors (actuallythey are matrices) Array to control if
69 // the motion is determined by an equation or is from a file.(not Nektar++)
70 // check if we need to load a file or we have an equation
71 CheckIsFromFile(pForce);
72
73 // Initialise movingbody filter
74 InitialiseFilter(m_session, pFields, pForce);
75
76 // Initialise the cable model
78
79 // Load mapping
81 m_mapping->SetTimeDependent(true);
82
83 if (m_vdim > 0)
84 {
85 m_mapping->SetFromFunction(false);
86 }
87
90 // What are this bi-dimensional vectors
91 // ------------------------------------------ m_zta[0] = m_zta | m_eta[0] =
92 // m_eta | m_zta[1] = d(m_zta)/dt |
93 // m_eta[1] = d(m_eta)/dt | m_zta[2] = dd(m_zta)/ddtt |
94 // m_eta[2] = dd(m_eta)/ddtt |
95 //--------------------------------------------------------------------------------
96 size_t phystot = pFields[0]->GetTotPoints();
97 for (size_t i = 0; i < m_zta.size(); i++)
98 {
99 m_zta[i] = Array<OneD, NekDouble>(phystot, 0.0);
100 m_eta[i] = Array<OneD, NekDouble>(phystot, 0.0);
101 }
102}
103
106 const Array<OneD, Array<OneD, NekDouble>> &inarray,
107 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble &time)
108{
109 boost::ignore_unused(inarray, outarray);
110
111 // Update the forces from the calculation of fluid field, which is
112 // implemented in the movingbody filter
113 Array<OneD, NekDouble> Hydroforces(2 * m_np, 0.0);
114 m_MovBodyfilter->UpdateForce(m_session, pFields, Hydroforces, time);
115
116 // for "free" type (m_vdim = 2), the cable vibrates both in streamwise and
117 // crossflow dimections, for "constrained" type (m_vdim = 1), the cable only
118 // vibrates in crossflow dimection, and for "forced" type (m_vdim = 0), the
119 // calbe vibrates specifically along a given function or file
120 if (m_vdim == 1 || m_vdim == 2)
121 {
122 // For free vibration case, displacements, velocities and acceleartions
123 // are obtained through solving structure dynamic model
124 EvaluateStructDynModel(pFields, Hydroforces, time);
125
126 // Convert result to format required by mapping
127 size_t physTot = pFields[0]->GetTotPoints();
130 for (size_t i = 0; i < 3; i++)
131 {
132 coords[i] = Array<OneD, NekDouble>(physTot, 0.0);
133 coordsVel[i] = Array<OneD, NekDouble>(physTot, 0.0);
134 }
135 // Get original coordinates
136 pFields[0]->GetCoords(coords[0], coords[1], coords[2]);
137
138 // Add displacement to coordinates
139 Vmath::Vadd(physTot, coords[0], 1, m_zta[0], 1, coords[0], 1);
140 Vmath::Vadd(physTot, coords[1], 1, m_eta[0], 1, coords[1], 1);
141 // Copy velocities
142 Vmath::Vcopy(physTot, m_zta[1], 1, coordsVel[0], 1);
143 Vmath::Vcopy(physTot, m_eta[1], 1, coordsVel[1], 1);
144
145 // Update mapping
146 m_mapping->UpdateMapping(time, coords, coordsVel);
147 }
148 else if (m_vdim == 0)
149 {
150 // For forced vibration case, load from specified file or function
151 size_t cnt = 0;
152 for (size_t j = 0; j < m_funcName.size(); j++)
153 {
154 if (m_IsFromFile[cnt] && m_IsFromFile[cnt + 1])
155 {
156 ASSERTL0(false, "Motion loading from file needs specific "
157 "implementation: Work in Progress!");
158 }
159 else
160 {
161 GetFunction(pFields, m_session, m_funcName[j], true)
162 ->Evaluate(m_motion[0], m_zta[j], time);
163 GetFunction(pFields, m_session, m_funcName[j], true)
164 ->Evaluate(m_motion[1], m_eta[j], time);
165 cnt = cnt + 2;
166 }
167 }
168
169 // Update mapping
170 m_mapping->UpdateMapping(time);
171
172 // Convert result from mapping
173 size_t physTot = pFields[0]->GetTotPoints();
176 for (size_t i = 0; i < 3; i++)
177 {
178 coords[i] = Array<OneD, NekDouble>(physTot, 0.0);
179 coordsVel[i] = Array<OneD, NekDouble>(physTot, 0.0);
180 }
181 // Get original coordinates
182 pFields[0]->GetCoords(coords[0], coords[1], coords[2]);
183
184 // Get Coordinates and coord velocity from mapping
185 m_mapping->GetCartesianCoordinates(m_zta[0], m_eta[0], coords[2]);
186 m_mapping->GetCoordVelocity(coordsVel);
187
188 // Calculate displacement
189 Vmath::Vsub(physTot, m_zta[0], 1, coords[0], 1, m_zta[0], 1);
190 Vmath::Vsub(physTot, m_eta[0], 1, coords[1], 1, m_eta[0], 1);
191
192 Vmath::Vcopy(physTot, coordsVel[0], 1, m_zta[1], 1);
193 Vmath::Vcopy(physTot, coordsVel[1], 1, m_eta[1], 1);
194
195 for (size_t var = 0; var < 3; var++)
196 {
197 for (size_t plane = 0; plane < m_np; plane++)
198 {
199 size_t n = pFields[0]->GetPlane(plane)->GetTotPoints();
200 size_t offset = plane * n;
201 size_t Offset = var * m_np + plane;
202
203 m_MotionVars[0][Offset] = m_zta[var][offset];
204 m_MotionVars[1][Offset] = m_eta[var][offset];
205 }
206 }
207 }
208 else
209 {
210 ASSERTL0(false, "Unrecogized vibration type for cable's dynamic model");
211 }
212
213 LibUtilities::CommSharedPtr vcomm = pFields[0]->GetComm();
214 size_t colrank = vcomm->GetColumnComm()->GetRank();
215 // Pass the variables of the cable's motion to the movingbody filter
216 if (colrank == 0)
217 {
218 size_t n = m_MotionVars[0].size();
219 Array<OneD, NekDouble> tmpArray(2 * n), tmp(n);
220 Vmath::Vcopy(n, m_MotionVars[0], 1, tmpArray, 1);
221 Vmath::Vcopy(n, m_MotionVars[1], 1, tmp = tmpArray + n, 1);
222 m_MovBodyfilter->UpdateMotion(m_session, pFields, tmpArray, time);
223 }
224}
225
226/**
227 *
228 */
231 Array<OneD, NekDouble> &Hydroforces, NekDouble time)
232{
233 boost::ignore_unused(time);
234
235 LibUtilities::CommSharedPtr vcomm = pFields[0]->GetComm();
236 size_t colrank = vcomm->GetColumnComm()->GetRank();
237 size_t nproc = vcomm->GetColumnComm()->GetSize();
238
239 bool homostrip;
240 m_session->MatchSolverInfo("HomoStrip", "True", homostrip, false);
241
242 // number of structral modes and number of strips
243 size_t npts, nstrips;
244 if (!homostrip) // full resolutions
245 {
246 npts = m_session->GetParameter("HomModesZ");
247 }
248 else
249 {
250 m_session->LoadParameter("HomStructModesZ", npts);
251 m_session->LoadParameter("Strip_Z", nstrips);
252
253 // ASSERTL0(nstrips < = npts,
254 // "the number of struct. modes must be larger than that of the
255 // strips.");
256 }
257
258 // the hydrodynamic forces
260 // forces in x-direction
261 fces[0] = Array<OneD, NekDouble>(npts, 0.0);
262 // forces in y-direction
263 fces[1] = Array<OneD, NekDouble>(npts, 0.0);
264
265 // fill the force vectors
266 if (colrank == 0)
267 {
268 if (!homostrip) // full resolutions
269 {
270 Vmath::Vcopy(m_np, Hydroforces, 1, fces[0], 1);
271 Vmath::Vcopy(m_np, Hydroforces + m_np, 1, fces[1], 1);
272 }
273 else // strip modelling
274 {
275 fces[0][0] = Hydroforces[0];
276 fces[1][0] = Hydroforces[m_np];
277 }
278
279 if (!homostrip) // full resolutions
280 {
282 for (size_t i = 1; i < nproc; ++i)
283 {
284 vcomm->GetColumnComm()->Recv(i, tmp);
285 for (size_t n = 0; n < m_np; n++)
286 {
287 for (size_t j = 0; j < 2; ++j)
288 {
289 fces[j][i * m_np + n] = tmp[j * m_np + n];
290 }
291 }
292 }
293 }
294 else // strip modelling
295 // if the body is submerged partly, the fces are filled partly
296 // by the flow induced forces
297 {
299 for (size_t i = 1; i < nstrips; ++i)
300 {
301 vcomm->GetColumnComm()->Recv(i, tmp);
302
303 for (size_t j = 0; j < 2; ++j)
304 {
305 fces[j][i] = tmp[j];
306 }
307 }
308 }
309 }
310 else
311 {
312 if (!homostrip) // full resolutions
313 {
314 vcomm->GetColumnComm()->Send(0, Hydroforces);
315 }
316 else // strip modelling
317 {
318 for (size_t i = 1; i < nstrips; ++i)
319 {
320 if (colrank == i)
321 {
323 tmp[0] = Hydroforces[0];
324 tmp[1] = Hydroforces[m_np];
325 vcomm->GetColumnComm()->Send(0, tmp);
326 }
327 }
328 }
329 }
330
331 if (colrank == 0)
332 {
333 // Fictitious mass method used to stablize the explicit coupling at
334 // relatively lower mass ratio
335 bool fictmass;
336 m_session->MatchSolverInfo("FictitiousMassMethod", "True", fictmass,
337 false);
338 if (fictmass)
339 {
340 NekDouble fictrho, fictdamp;
341 m_session->LoadParameter("FictMass", fictrho);
342 m_session->LoadParameter("FictDamp", fictdamp);
343
344 static NekDouble Betaq_Coeffs[2][2] = {{1.0, 0.0}, {2.0, -1.0}};
345
346 // only consider second order approximation for fictitious variables
347 size_t intOrder = 2;
348 size_t nint = min(m_movingBodyCalls + 1, intOrder);
349 size_t nlevels = m_fV[0].size();
350
351 for (size_t i = 0; i < m_motion.size(); ++i)
352 {
353 RollOver(m_fV[i]);
354 RollOver(m_fA[i]);
355
356 size_t Voffset = npts;
357 size_t Aoffset = 2 * npts;
358
359 Vmath::Vcopy(npts, m_MotionVars[i] + Voffset, 1, m_fV[i][0], 1);
360 Vmath::Vcopy(npts, m_MotionVars[i] + Aoffset, 1, m_fA[i][0], 1);
361
362 // Extrapolate to n+1
363 Vmath::Smul(npts, Betaq_Coeffs[nint - 1][nint - 1],
364 m_fV[i][nint - 1], 1, m_fV[i][nlevels - 1], 1);
365 Vmath::Smul(npts, Betaq_Coeffs[nint - 1][nint - 1],
366 m_fA[i][nint - 1], 1, m_fA[i][nlevels - 1], 1);
367
368 for (size_t n = 0; n < nint - 1; ++n)
369 {
370 Vmath::Svtvp(npts, Betaq_Coeffs[nint - 1][n], m_fV[i][n], 1,
371 m_fV[i][nlevels - 1], 1, m_fV[i][nlevels - 1],
372 1);
373 Vmath::Svtvp(npts, Betaq_Coeffs[nint - 1][n], m_fA[i][n], 1,
374 m_fA[i][nlevels - 1], 1, m_fA[i][nlevels - 1],
375 1);
376 }
377
378 // Add the fictitious forces on the RHS of the equation
379 Vmath::Svtvp(npts, fictdamp, m_fV[i][nlevels - 1], 1, fces[i],
380 1, fces[i], 1);
381 Vmath::Svtvp(npts, fictrho, m_fA[i][nlevels - 1], 1, fces[i], 1,
382 fces[i], 1);
383 }
384 }
385 }
386 // structural solver is implemented on the root process
387 if (colrank == 0)
388 {
389 // Tensioned cable model is evaluated in wave space
390 for (size_t n = 0, cn = 1; n < m_vdim; n++, cn--)
391 {
392 Newmark_betaSolver(pFields, fces[cn], m_MotionVars[cn]);
393 }
394 }
395
396 Array<OneD, NekDouble> Motvars(2 * 2 * m_np);
397 // send physical coeffients to all planes of each processor
398 if (!homostrip) // full resolutions
399 {
400 Array<OneD, NekDouble> tmp(2 * 2 * m_np);
401
402 if (colrank != 0)
403 {
404 vcomm->GetColumnComm()->Recv(0, tmp);
405 Vmath::Vcopy(2 * 2 * m_np, tmp, 1, Motvars, 1);
406 }
407 else
408 {
409 for (size_t i = 1; i < nproc; ++i)
410 {
411 for (size_t j = 0; j < 2; j++) // moving dimensions
412 {
413 for (size_t k = 0; k < 2; k++) // disp. and vel.
414 {
415 for (size_t n = 0; n < m_np; n++)
416 {
417 tmp[j * 2 * m_np + k * m_np + n] =
418 m_MotionVars[j][k * npts + i * m_np + n];
419 }
420 }
421 }
422 vcomm->GetColumnComm()->Send(i, tmp);
423 }
424
425 for (size_t j = 0; j < 2; j++)
426 {
427 for (size_t k = 0; k < 2; k++)
428 {
429 for (size_t n = 0; n < m_np; n++)
430 {
431 tmp[j * 2 * m_np + k * m_np + n] =
432 m_MotionVars[j][k * npts + n];
433 }
434 Vmath::Vcopy(2 * 2 * m_np, tmp, 1, Motvars, 1);
435 }
436 }
437 }
438 }
439 else // strip modelling
440 {
441 Array<OneD, NekDouble> tmp(2 * 2);
442
443 if (colrank != 0)
444 {
445 for (size_t j = 1; j < nproc / nstrips; j++)
446 {
447 if (colrank == j * nstrips)
448 {
449 vcomm->GetColumnComm()->Recv(0, tmp);
450
451 for (size_t plane = 0; plane < m_np; plane++)
452 {
453 for (size_t var = 0; var < 2; var++)
454 {
455 for (size_t k = 0; k < 2; k++)
456 {
457 Motvars[var * 2 * m_np + k * m_np + plane] =
458 tmp[var * 2 + k];
459 }
460 }
461 }
462 }
463 }
464
465 for (size_t i = 1; i < nstrips; i++)
466 {
467 for (size_t j = 0; j < nproc / nstrips; j++)
468 {
469 if (colrank == i + j * nstrips)
470 {
471 vcomm->GetColumnComm()->Recv(0, tmp);
472
473 for (size_t plane = 0; plane < m_np; plane++)
474 {
475 for (size_t var = 0; var < 2; var++)
476 {
477 for (size_t k = 0; k < 2; k++)
478 {
479 Motvars[var * 2 * m_np + k * m_np + plane] =
480 tmp[var * 2 + k];
481 }
482 }
483 }
484 }
485 }
486 }
487 }
488 else
489 {
490 for (size_t j = 0; j < 2; ++j)
491 {
492 for (size_t k = 0; k < 2; ++k)
493 {
494 tmp[j * 2 + k] = m_MotionVars[j][k * npts];
495 }
496 }
497
498 for (size_t j = 1; j < nproc / nstrips; j++)
499 {
500 vcomm->GetColumnComm()->Send(j * nstrips, tmp);
501 }
502
503 for (size_t plane = 0; plane < m_np; plane++)
504 {
505 for (size_t j = 0; j < 2; j++)
506 {
507 for (size_t k = 0; k < 2; ++k)
508 {
509 Motvars[j * 2 * m_np + k * m_np + plane] =
510 m_MotionVars[j][k * npts];
511 }
512 }
513 }
514
515 for (size_t i = 1; i < nstrips; ++i)
516 {
517 for (size_t j = 0; j < 2; ++j)
518 {
519 for (size_t k = 0; k < 2; ++k)
520 {
521 tmp[j * 2 + k] = m_MotionVars[j][i + k * npts];
522 }
523 }
524
525 for (size_t j = 0; j < nproc / nstrips; j++)
526 {
527 vcomm->GetColumnComm()->Send(i + j * nstrips, tmp);
528 }
529 }
530 }
531 }
532
533 // Set the m_forcing term based on the motion of the cable
534 for (size_t var = 0; var < 2; var++)
535 {
536 for (size_t plane = 0; plane < m_np; plane++)
537 {
538 size_t n = pFields[0]->GetPlane(plane)->GetTotPoints();
539
541
542 size_t offset = plane * n;
543 size_t xoffset = var * m_np + plane;
544 size_t yoffset = 2 * m_np + xoffset;
545
546 Vmath::Fill(n, Motvars[xoffset], tmp = m_zta[var] + offset, 1);
547 Vmath::Fill(n, Motvars[yoffset], tmp = m_eta[var] + offset, 1);
548 }
549 }
550}
551
552/**
553 *
554 */
557 Array<OneD, NekDouble> &HydroForces, Array<OneD, NekDouble> &BodyMotions)
558{
559 boost::ignore_unused(pFields);
560
561 std::string supptype = m_session->GetSolverInfo("SupportType");
562
563 size_t npts = HydroForces.size();
564
567
568 for (size_t i = 0; i < 4; i++)
569 {
570 fft_i[i] = Array<OneD, NekDouble>(npts, 0.0);
571 fft_o[i] = Array<OneD, NekDouble>(npts, 0.0);
572 }
573
574 Vmath::Vcopy(npts, HydroForces, 1, fft_i[0], 1);
575 for (size_t i = 0; i < 3; i++)
576 {
577 Vmath::Vcopy(npts, BodyMotions + i * npts, 1, fft_i[i + 1], 1);
578 }
579
580 // Implement Fourier transformation of the motion variables
581 if (boost::iequals(supptype, "Free-Free"))
582 {
583 for (size_t j = 0; j < 4; ++j)
584 {
585 m_FFT->FFTFwdTrans(fft_i[j], fft_o[j]);
586 }
587 }
588 else if (boost::iequals(supptype, "Pinned-Pinned"))
589 {
590 // TODO:
591 size_t N = fft_i[0].size();
592
593 for (size_t var = 0; var < 4; var++)
594 {
595 for (size_t i = 0; i < N; i++)
596 {
597 fft_o[var][i] = 0;
598
599 for (size_t k = 0; k < N; k++)
600 {
601 fft_o[var][i] +=
602 fft_i[var][k] * sin(M_PI / (N) * (k + 1 / 2) * (i + 1));
603 }
604 }
605 }
606 }
607 else
608 {
609 ASSERTL0(false, "Unrecognized support type for cable's motion");
610 }
611
612 // solve the ODE in the wave space
613 for (size_t i = 0; i < npts; ++i)
614 {
615 size_t nrows = 3;
616
617 Array<OneD, NekDouble> tmp0, tmp1, tmp2;
618 tmp0 = Array<OneD, NekDouble>(3, 0.0);
619 tmp1 = Array<OneD, NekDouble>(3, 0.0);
620 tmp2 = Array<OneD, NekDouble>(3, 0.0);
621
622 for (size_t var = 0; var < 3; var++)
623 {
624 tmp0[var] = fft_o[var + 1][i];
625 }
626
627 tmp2[0] = fft_o[0][i];
628
629 Blas::Dgemv('N', nrows, nrows, 1.0, &(m_CoeffMat_B[i]->GetPtr())[0],
630 nrows, &tmp0[0], 1, 0.0, &tmp1[0], 1);
631 Blas::Dgemv('N', nrows, nrows, 1.0 / m_structrho,
632 &(m_CoeffMat_A[i]->GetPtr())[0], nrows, &tmp2[0], 1, 1.0,
633 &tmp1[0], 1);
634
635 for (size_t var = 0; var < 3; var++)
636 {
637 fft_i[var][i] = tmp1[var];
638 }
639 }
640
641 // get physical coeffients via Backward fourier transformation of wave
642 // coefficients
643 if (boost::iequals(supptype, "Free-Free"))
644 {
645 for (size_t var = 0; var < 3; var++)
646 {
647 m_FFT->FFTBwdTrans(fft_i[var], fft_o[var]);
648 }
649 }
650 else if (boost::iequals(supptype, "Pinned-Pinned"))
651 {
652 // TODO:
653 size_t N = fft_i[0].size();
654
655 for (size_t var = 0; var < 3; var++)
656 {
657 for (size_t i = 0; i < N; i++)
658 {
659 fft_o[var][i] = 0;
660
661 for (size_t k = 0; k < N; k++)
662 {
663 fft_o[var][i] += fft_i[var][k] *
664 sin(M_PI / (N) * (k + 1) * (i + 1 / 2)) *
665 2 / N;
666 }
667 }
668 }
669 }
670 else
671 {
672 ASSERTL0(false, "Unrecognized support type for cable's motion");
673 }
674
675 for (size_t i = 0; i < 3; i++)
676 {
677 Array<OneD, NekDouble> tmp(npts, 0.0);
678 Vmath::Vcopy(npts, fft_o[i], 1, tmp = BodyMotions + i * npts, 1);
679 }
680}
681
682/**
683 *
684 */
688{
689 boost::ignore_unused(pSession, pFields);
690
692 m_session->LoadParameter("Kinvis", m_kinvis);
693 m_session->LoadParameter("TimeStep", m_timestep, 0.01);
694
695 LibUtilities::CommSharedPtr vcomm = pFields[0]->GetComm();
696 size_t colrank = vcomm->GetColumnComm()->GetRank();
697 size_t nproc = vcomm->GetColumnComm()->GetSize();
698
699 // number of structral modes
700 size_t npts;
701 bool homostrip;
702 m_session->MatchSolverInfo("HomoStrip", "True", homostrip, false);
703
704 if (!homostrip) // full resolutions
705 {
706 npts = m_session->GetParameter("HomModesZ");
707 }
708 else
709 {
710 m_session->LoadParameter("HomStructModesZ", npts);
711 }
712
714 m_MotionVars[0] = Array<OneD, NekDouble>(3 * npts, 0.0);
715 m_MotionVars[1] = Array<OneD, NekDouble>(3 * npts, 0.0);
716
718 ZIDs = pFields[0]->GetZIDs();
719 m_np = ZIDs.size();
720
721 std::string vibtype = m_session->GetSolverInfo("VibrationType");
722
723 if (boost::iequals(vibtype, "Constrained"))
724 {
725 m_vdim = 1;
726 }
727 else if (boost::iequals(vibtype, "Free"))
728 {
729 m_vdim = 2;
730 }
731 else if (boost::iequals(vibtype, "Forced"))
732 {
733 m_vdim = 0;
734 return;
735 }
736
737 if (!homostrip)
738 {
739 m_session->LoadParameter("LZ", m_lhom);
740 int nplanes = m_session->GetParameter("HomModesZ");
742 nplanes);
743 }
744 else
745 {
746 int nstrips;
747 NekDouble DistStrip;
748
749 m_session->LoadParameter("DistStrip", DistStrip);
750 m_session->LoadParameter("Strip_Z", nstrips);
751 m_lhom = nstrips * DistStrip;
753 nstrips);
754 }
755
756 // load the structural dynamic parameters from xml file
757 m_session->LoadParameter("StructRho", m_structrho);
758 m_session->LoadParameter("StructDamp", m_structdamp, 0.0);
759
760 // Identify whether the fictitious mass method is active for explicit
761 // coupling of fluid solver and structural dynamics solver
762 bool fictmass;
763 m_session->MatchSolverInfo("FictitiousMassMethod", "True", fictmass, false);
764 if (fictmass)
765 {
766 NekDouble fictrho, fictdamp;
767 m_session->LoadParameter("FictMass", fictrho);
768 m_session->LoadParameter("FictDamp", fictdamp);
769 m_structrho += fictrho;
770 m_structdamp += fictdamp;
771
772 // Storage array of Struct Velocity and Acceleration used for
773 // extrapolation of fictitious force
776 for (size_t i = 0; i < m_motion.size(); ++i)
777 {
780
781 for (size_t n = 0; n < 2; ++n)
782 {
783 m_fV[i][n] = Array<OneD, NekDouble>(npts, 0.0);
784 m_fA[i][n] = Array<OneD, NekDouble>(npts, 0.0);
785 }
786 }
787 }
788
789 // Setting the coefficient matrices for solving structural dynamic ODEs
790 SetDynEqCoeffMatrix(pFields);
791
792 // Set initial condition for cable's motion
793 size_t cnt = 0;
794
795 for (size_t j = 0; j < m_funcName.size(); j++)
796 {
797 // loading from the specified files through inputstream
798 if (m_IsFromFile[cnt] && m_IsFromFile[cnt + 1])
799 {
800 std::ifstream inputStream;
801 size_t nzpoints;
802
803 if (homostrip)
804 {
805 m_session->LoadParameter("HomStructModesZ", nzpoints);
806 }
807 else
808 {
809 nzpoints = pFields[0]->GetHomogeneousBasis()->GetNumModes();
810 }
811
812 if (vcomm->GetRank() == 0)
813 {
814 std::string filename =
815 m_session->GetFunctionFilename(m_funcName[j], m_motion[0]);
816
817 // Open intputstream for cable motions
818 inputStream.open(filename.c_str());
819
820 // Import the head string from the file
822 for (size_t n = 0; n < tmp.size(); n++)
823 {
824 inputStream >> tmp[n];
825 }
826
827 NekDouble time, z_cds;
828 // Import the motion variables from the file
829 for (size_t n = 0; n < nzpoints; n++)
830 {
831 inputStream >> setprecision(6) >> time;
832 inputStream >> setprecision(6) >> z_cds;
833 inputStream >> setprecision(8) >> m_MotionVars[0][n];
834 inputStream >> setprecision(8) >>
835 m_MotionVars[0][n + nzpoints];
836 inputStream >> setprecision(8) >>
837 m_MotionVars[0][n + 2 * nzpoints];
838 inputStream >> setprecision(8) >> m_MotionVars[1][n];
839 inputStream >> setprecision(8) >>
840 m_MotionVars[1][n + nzpoints];
841 inputStream >> setprecision(8) >>
842 m_MotionVars[1][n + 2 * nzpoints];
843 }
844 // Close inputstream for cable motions
845 inputStream.close();
846 }
847 cnt = cnt + 2;
848 }
849 else // Evaluate from the functions specified in xml file
850 {
851 if (!homostrip)
852 {
853 size_t nzpoints =
854 pFields[0]->GetHomogeneousBasis()->GetNumModes();
855 Array<OneD, NekDouble> z_coords(nzpoints, 0.0);
857 pFields[0]->GetHomogeneousBasis()->GetZ();
858
859 Vmath::Smul(nzpoints, m_lhom / 2.0, pts, 1, z_coords, 1);
860 Vmath::Sadd(nzpoints, m_lhom / 2.0, z_coords, 1, z_coords, 1);
861
865 for (size_t plane = 0; plane < m_np; plane++)
866 {
867 x2[plane] = z_coords[ZIDs[plane]];
868 }
869
871 Array<OneD, NekDouble> tmp0(m_np, 0.0);
872 Array<OneD, NekDouble> tmp1(m_np, 0.0);
873 LibUtilities::EquationSharedPtr ffunc0, ffunc1;
874
875 ffunc0 = m_session->GetFunction(m_funcName[j], m_motion[0]);
876 ffunc1 = m_session->GetFunction(m_funcName[j], m_motion[1]);
877 ffunc0->Evaluate(x0, x1, x2, 0.0, tmp0);
878 ffunc1->Evaluate(x0, x1, x2, 0.0, tmp1);
879
880 size_t offset = j * npts;
881 Vmath::Vcopy(m_np, tmp0, 1, tmp = m_MotionVars[0] + offset, 1);
882 Vmath::Vcopy(m_np, tmp1, 1, tmp = m_MotionVars[1] + offset, 1);
883
884 if (colrank == 0)
885 {
886 for (size_t i = 1; i < nproc; ++i)
887 {
888 vcomm->GetColumnComm()->Recv(i, tmp0);
889 vcomm->GetColumnComm()->Recv(i, tmp1);
890 Vmath::Vcopy(m_np, tmp0, 1,
891 tmp = m_MotionVars[0] + offset + i * m_np,
892 1);
893 Vmath::Vcopy(m_np, tmp1, 1,
894 tmp = m_MotionVars[1] + offset + i * m_np,
895 1);
896 }
897 }
898 else
899 {
900 vcomm->GetColumnComm()->Send(0, tmp0);
901 vcomm->GetColumnComm()->Send(0, tmp1);
902 }
903 }
904 else
905 {
906 if (colrank == 0)
907 {
908 size_t nstrips;
909 m_session->LoadParameter("Strip_Z", nstrips);
910
911 ASSERTL0(
912 m_session->DefinesSolverInfo("USEFFT"),
913 "Fourier transformation of cable motion is currently "
914 "implemented only for FFTW module.");
915
916 NekDouble DistStrip;
917 m_session->LoadParameter("DistStrip", DistStrip);
918
919 Array<OneD, NekDouble> x0(npts, 0.0);
920 Array<OneD, NekDouble> x1(npts, 0.0);
921 Array<OneD, NekDouble> x2(npts, 0.0);
922 Array<OneD, NekDouble> tmp(npts, 0.0);
923 Array<OneD, NekDouble> tmp0(npts, 0.0);
924 Array<OneD, NekDouble> tmp1(npts, 0.0);
925 for (size_t plane = 0; plane < npts; plane++)
926 {
927 x2[plane] = plane * DistStrip;
928 }
929 LibUtilities::EquationSharedPtr ffunc0, ffunc1;
930 ffunc0 = m_session->GetFunction(m_funcName[j], m_motion[0]);
931 ffunc1 = m_session->GetFunction(m_funcName[j], m_motion[1]);
932 ffunc0->Evaluate(x0, x1, x2, 0.0, tmp0);
933 ffunc1->Evaluate(x0, x1, x2, 0.0, tmp1);
934
935 size_t offset = j * npts;
936 Vmath::Vcopy(npts, tmp0, 1, tmp = m_MotionVars[0] + offset,
937 1);
938 Vmath::Vcopy(npts, tmp1, 1, tmp = m_MotionVars[1] + offset,
939 1);
940 }
941 }
942
943 cnt = cnt + 2;
944 }
945 }
946}
947
948/**
949 *
950 */
953{
954 boost::ignore_unused(pFields);
955
956 size_t nplanes;
957
958 bool homostrip;
959 m_session->MatchSolverInfo("HomoStrip", "True", homostrip, false);
960
961 if (!homostrip)
962 {
963 nplanes = m_session->GetParameter("HomModesZ");
964 }
965 else
966 {
967 m_session->LoadParameter("Strip_Z", nplanes);
968 }
969
972
973 NekDouble tmp1, tmp2, tmp3;
974 NekDouble tmp4, tmp5, tmp6, tmp7;
975
976 // load the structural dynamic parameters from xml file
977 NekDouble cabletension;
978 NekDouble bendingstiff;
979 NekDouble structstiff;
980 m_session->LoadParameter("StructStiff", structstiff, 0.0);
981 m_session->LoadParameter("CableTension", cabletension, 0.0);
982 m_session->LoadParameter("BendingStiff", bendingstiff, 0.0);
983
984 tmp1 = m_timestep * m_timestep;
985 tmp2 = structstiff / m_structrho;
986 tmp3 = m_structdamp / m_structrho;
987 tmp4 = cabletension / m_structrho;
988 tmp5 = bendingstiff / m_structrho;
989
990 // solve the ODE in the wave space for cable motion to obtain disp, vel and
991 // accel
992
993 std::string supptype = m_session->GetSolverInfo("SupportType");
994
995 for (size_t plane = 0; plane < nplanes; plane++)
996 {
997 int nel = 3;
998 m_CoeffMat_A[plane] =
1000 m_CoeffMat_B[plane] =
1002
1003 // Initialised to avoid compiler warnings.
1004 unsigned int K = 0;
1005 NekDouble beta = 0.0;
1006
1007 if (boost::iequals(supptype, "Free-Free"))
1008 {
1009 K = plane / 2;
1010 beta = 2.0 * M_PI / m_lhom;
1011 }
1012 else if (boost::iequals(supptype, "Pinned-Pinned"))
1013 {
1014 K = plane + 1;
1015 beta = M_PI / m_lhom;
1016 }
1017 else
1018 {
1019 ASSERTL0(false, "Unrecognized support type for cable's motion");
1020 }
1021
1022 tmp6 = beta * K;
1023 tmp6 = tmp6 * tmp6;
1024 tmp7 = tmp6 * tmp6;
1025 tmp7 = tmp2 + tmp4 * tmp6 + tmp5 * tmp7;
1026
1027 (*m_CoeffMat_A[plane])(0, 0) = tmp7;
1028 (*m_CoeffMat_A[plane])(0, 1) = tmp3;
1029 (*m_CoeffMat_A[plane])(0, 2) = 1.0;
1030 (*m_CoeffMat_A[plane])(1, 0) = 0.0;
1031 (*m_CoeffMat_A[plane])(1, 1) = 1.0;
1032 (*m_CoeffMat_A[plane])(1, 2) = -m_timestep / 2.0;
1033 (*m_CoeffMat_A[plane])(2, 0) = 1.0;
1034 (*m_CoeffMat_A[plane])(2, 1) = 0.0;
1035 (*m_CoeffMat_A[plane])(2, 2) = -tmp1 / 4.0;
1036 (*m_CoeffMat_B[plane])(0, 0) = 0.0;
1037 (*m_CoeffMat_B[plane])(0, 1) = 0.0;
1038 (*m_CoeffMat_B[plane])(0, 2) = 0.0;
1039 (*m_CoeffMat_B[plane])(1, 0) = 0.0;
1040 (*m_CoeffMat_B[plane])(1, 1) = 1.0;
1041 (*m_CoeffMat_B[plane])(1, 2) = m_timestep / 2.0;
1042 (*m_CoeffMat_B[plane])(2, 0) = 1.0;
1043 (*m_CoeffMat_B[plane])(2, 1) = m_timestep;
1044 (*m_CoeffMat_B[plane])(2, 2) = tmp1 / 4.0;
1045
1046 m_CoeffMat_A[plane]->Invert();
1047 (*m_CoeffMat_B[plane]) =
1048 (*m_CoeffMat_A[plane]) * (*m_CoeffMat_B[plane]);
1049 }
1050}
1051
1052/**
1053 * Function to roll time-level storages to the next step layout.
1054 * The stored data associated with the oldest time-level
1055 * (not required anymore) are moved to the top, where they will
1056 * be overwritten as the solution process progresses.
1057 */
1059{
1060 size_t nlevels = input.size();
1061
1063 tmp = input[nlevels - 1];
1064
1065 for (size_t n = nlevels - 1; n > 0; --n)
1066 {
1067 input[n] = input[n - 1];
1068 }
1069
1070 input[0] = tmp;
1071}
1072
1073/**
1074 *
1075 */
1076void ForcingMovingBody::CheckIsFromFile(const TiXmlElement *pForce)
1077{
1078
1081 m_motion[0] = "x";
1082 m_motion[1] = "y";
1083
1085 // Loading the x-dispalcement (m_zta) and the y-displacement (m_eta)
1086 // Those two variables are bith functions of z and t and the may come
1087 // from an equation (forced vibration) or from another solver which, given
1088 // the aerodynamic forces at the previous step, calculates the
1089 // displacements.
1090
1091 // Get the body displacement: m_zta and m_eta
1092 const TiXmlElement *funcNameElmt_D =
1093 pForce->FirstChildElement("DISPLACEMENTS");
1094 ASSERTL0(funcNameElmt_D,
1095 "MOVINGBODYFORCE tag has to specify a function name which "
1096 "prescribes the body displacement as d(z,t).");
1097
1098 m_funcName[0] = funcNameElmt_D->GetText();
1099 ASSERTL0(m_session->DefinesFunction(m_funcName[0]),
1100 "Function '" + m_funcName[0] + "' not defined.");
1101
1102 // Get the body velocity of movement: d(m_zta)/dt and d(m_eta)/dt
1103 const TiXmlElement *funcNameElmt_V =
1104 pForce->FirstChildElement("VELOCITIES");
1105 ASSERTL0(funcNameElmt_D,
1106 "MOVINGBODYFORCE tag has to specify a function name which "
1107 "prescribes the body velocity of movement as v(z,t).");
1108
1109 m_funcName[1] = funcNameElmt_V->GetText();
1110 ASSERTL0(m_session->DefinesFunction(m_funcName[1]),
1111 "Function '" + m_funcName[1] + "' not defined.");
1112
1113 // Get the body acceleration: dd(m_zta)/ddt and dd(m_eta)/ddt
1114 const TiXmlElement *funcNameElmt_A =
1115 pForce->FirstChildElement("ACCELERATIONS");
1116 ASSERTL0(funcNameElmt_A,
1117 "MOVINGBODYFORCE tag has to specify a function name which "
1118 "prescribes the body acceleration as a(z,t).");
1119
1120 m_funcName[2] = funcNameElmt_A->GetText();
1121 ASSERTL0(m_session->DefinesFunction(m_funcName[2]),
1122 "Function '" + m_funcName[2] + "' not defined.");
1123
1125
1126 // Check Displacement x
1127 vType = m_session->GetFunctionType(m_funcName[0], m_motion[0]);
1129 {
1130 m_IsFromFile[0] = false;
1131 }
1132 else if (vType == LibUtilities::eFunctionTypeFile)
1133 {
1134 m_IsFromFile[0] = true;
1135 }
1136 else
1137 {
1138 ASSERTL0(false, "The displacements in x must be specified via an "
1139 "equation <E> or a file <F>");
1140 }
1141
1142 // Check Displacement y
1143 vType = m_session->GetFunctionType(m_funcName[0], m_motion[1]);
1145 {
1146 m_IsFromFile[1] = false;
1147 }
1148 else if (vType == LibUtilities::eFunctionTypeFile)
1149 {
1150 m_IsFromFile[1] = true;
1151 }
1152 else
1153 {
1154 ASSERTL0(false, "The displacements in y must be specified via an "
1155 "equation <E> or a file <F>");
1156 }
1157
1158 // Check Velocity x
1159 vType = m_session->GetFunctionType(m_funcName[1], m_motion[0]);
1161 {
1162 m_IsFromFile[2] = false;
1163 }
1164 else if (vType == LibUtilities::eFunctionTypeFile)
1165 {
1166 m_IsFromFile[2] = true;
1167 }
1168 else
1169 {
1170 ASSERTL0(false, "The velocities in x must be specified via an equation "
1171 "<E> or a file <F>");
1172 }
1173
1174 // Check Velocity y
1175 vType = m_session->GetFunctionType(m_funcName[1], m_motion[1]);
1177 {
1178 m_IsFromFile[3] = false;
1179 }
1180 else if (vType == LibUtilities::eFunctionTypeFile)
1181 {
1182 m_IsFromFile[3] = true;
1183 }
1184 else
1185 {
1186 ASSERTL0(false, "The velocities in y must be specified via an equation "
1187 "<E> or a file <F>");
1188 }
1189
1190 // Check Acceleration x
1191 vType = m_session->GetFunctionType(m_funcName[2], m_motion[0]);
1193 {
1194 m_IsFromFile[4] = false;
1195 }
1196 else if (vType == LibUtilities::eFunctionTypeFile)
1197 {
1198 m_IsFromFile[4] = true;
1199 }
1200 else
1201 {
1202 ASSERTL0(false, "The accelerations in x must be specified via an "
1203 "equation <E> or a file <F>");
1204 }
1205
1206 // Check Acceleration y
1207 vType = m_session->GetFunctionType(m_funcName[2], m_motion[1]);
1209 {
1210 m_IsFromFile[5] = false;
1211 }
1212 else if (vType == LibUtilities::eFunctionTypeFile)
1213 {
1214 m_IsFromFile[5] = true;
1215 }
1216 else
1217 {
1218 ASSERTL0(false, "The accelerations in y must be specified via an "
1219 "equation <E> or a file <F>");
1220 }
1221}
1222
1223/**
1224 *
1225 */
1229 const TiXmlElement *pForce)
1230{
1231 // Get the outputfile name, output frequency and
1232 // the boundary's ID for the cable's wall
1233 std::string typeStr = pForce->Attribute("TYPE");
1234 std::map<std::string, std::string> vParams;
1235
1236 const TiXmlElement *param = pForce->FirstChildElement("PARAM");
1237 while (param)
1238 {
1239 ASSERTL0(param->Attribute("NAME"),
1240 "Missing attribute 'NAME' for parameter in filter " + typeStr +
1241 "'.");
1242 std::string nameStr = param->Attribute("NAME");
1243
1244 ASSERTL0(param->GetText(), "Empty value string for param.");
1245 std::string valueStr = param->GetText();
1246
1247 vParams[nameStr] = valueStr;
1248
1249 param = param->NextSiblingElement("PARAM");
1250 }
1251
1252 // Creat the filter for MovingBody, where we performed the calculation of
1253 // fluid forces and write both the aerodynamic forces and motion variables
1254 // into the output files
1256 pSession, m_equ, vParams);
1257
1258 // Initialise the object of MovingBody filter
1259 m_MovBodyfilter->Initialise(pFields, 0.0);
1260}
1261
1262} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
GlobalMapping::MappingSharedPtr m_mapping
Array< OneD, std::string > m_funcName
[0] is displacements, [1] is velocities, [2] is accelerations
ForcingMovingBody(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< SolverUtils::EquationSystem > &pEquation)
Array< OneD, Array< OneD, NekDouble > > m_zta
Store the derivatives of motion variables in x-direction.
NekDouble m_lhom
length ratio of the cable
NekDouble m_kinvis
fluid viscous
size_t m_np
number of planes per processors
Array< OneD, DNekMatSharedPtr > m_CoeffMat_A
matrices in Newmart-beta method
Array< OneD, Array< OneD, NekDouble > > m_eta
Store the derivatives of motion variables in y-direction.
void CheckIsFromFile(const TiXmlElement *pForce)
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_fV
fictitious velocity storage
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_fA
fictitious acceleration storage
void SetDynEqCoeffMatrix(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
void RollOver(Array< OneD, Array< OneD, NekDouble > > &input)
void InitialiseCableModel(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
NekDouble m_structrho
mass of the cable per unit length
size_t m_vdim
vibration dimension
NekDouble m_timestep
time step
virtual void v_Apply(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time) override
Array< OneD, bool > m_IsFromFile
do determine if the the body motion come from an extern file
FilterMovingBodySharedPtr m_MovBodyfilter
Array< OneD, Array< OneD, NekDouble > > m_MotionVars
storage for the cable's motion(x,y) variables
static SolverUtils::ForcingSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< SolverUtils::EquationSystem > &pEquation, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields, const TiXmlElement *pForce)
Creates an instance of this class.
Array< OneD, DNekMatSharedPtr > m_CoeffMat_B
matrices in Newmart-beta method
void EvaluateStructDynModel(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, Array< OneD, NekDouble > &Hydroforces, NekDouble time)
size_t m_movingBodyCalls
number of times the movbody have been called
void Newmark_betaSolver(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, Array< OneD, NekDouble > &FcePhysinArray, Array< OneD, NekDouble > &MotPhysinArray)
void InitialiseFilter(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pForce)
Array< OneD, std::string > m_motion
motion direction: [0] is 'x' and [1] is 'y'
LibUtilities::NektarFFTSharedPtr m_FFT
static std::string className
Name of the class.
virtual void v_InitObject(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields, const TiXmlElement *pForce) override
NekDouble m_structdamp
damping ratio of the cable
static GLOBAL_MAPPING_EXPORT MappingSharedPtr Load(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
Return a pointer to the mapping, creating it on first call.
Definition: Mapping.cpp:272
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
Defines a forcing term to be explicitly applied.
Definition: Forcing.h:73
const std::weak_ptr< EquationSystem > m_equ
Weak pointer to equation system using this forcing.
Definition: Forcing.h:119
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const LibUtilities::SessionReaderSharedPtr &pSession, std::string pName, bool pCache=false)
Get a SessionFunction by name.
Definition: Forcing.cpp:194
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
Definition: Forcing.h:117
static void Dgemv(const char &trans, const int &m, const int &n, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = alpha A x plus beta y where A[m x n].
Definition: Blas.hpp:213
NektarFFTFactory & GetNektarFFTFactory()
Definition: NektarFFT.cpp:69
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:129
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:61
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:57
ForcingFactory & GetForcingFactory()
Declaration of the forcing factory singleton.
Definition: Forcing.cpp:44
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:617
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:354
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:245
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:43
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add scalar y = alpha + x.
Definition: Vmath.cpp:379
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:414