Nektar++
ForcingMovingReferenceFrame.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: ForcingMovingReferenceFrame.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: Solving the absolute flow in a moving body frame,
32// by adding (U0 + Omega X (x - x0)) . grad u - Omega X u
33// as the body force.
34// U0 is the translational velocity of the body frame.
35// Omega is the angular velocity.
36// x0 is the rotation pivot in the body frame.
37// All vectors use the basis of the body frame.
38// Translational motion is allowed for all dimensions.
39// Rotation is not allowed for 1D, 2DH1D, 3DH2D.
40// Rotation in z direction is allowed for 2D and 3DH1D.
41// Rotation in 3 directions are allowed for 3D.
42// TODO: add suport for 3D rotation using Quaternion
43///////////////////////////////////////////////////////////////////////////////
44
50#include <boost/format.hpp>
51
52namespace Nektar::SolverUtils
53{
54
57 "MovingReferenceFrame", ForcingMovingReferenceFrame::create,
58 "Moving Frame");
59
60/**
61 * @brief
62 * @param pSession
63 * @param pEquation
64 */
67 const std::weak_ptr<EquationSystem> &pEquation)
68 : Forcing(pSession, pEquation)
69{
70}
71
73{
74 if (m_isRoot)
75 {
76 m_outputStream.close();
77 }
78}
79
80/**
81 * @brief Initialise the forcing module
82 * @param pFields
83 * @param pNumForcingFields
84 * @param pForce
85 */
88 [[maybe_unused]] const unsigned int &pNumForcingFields,
89 const TiXmlElement *pForce)
90{
91 m_session->MatchSolverInfo("Homogeneous", "1D", m_isH1d, false);
92 m_session->MatchSolverInfo("Homogeneous", "2D", m_isH2d, false);
93 bool singleMode, halfMode;
94 m_session->MatchSolverInfo("ModeType", "SingleMode", singleMode, false);
95 m_session->MatchSolverInfo("ModeType", "HalfMode", halfMode, false);
96 if (singleMode || halfMode)
97 {
98 m_isH1d = false;
99 }
100 m_expdim = m_isH2d ? 1 : pFields[0]->GetGraph()->GetMeshDimension();
101 m_spacedim = m_expdim + (m_isH1d ? 1 : 0) + (m_isH2d ? 2 : 0);
102 m_isRoot = pFields[0]->GetComm()->TreatAsRankZero();
103 m_hasPlane0 = true;
104 if (m_isH1d)
105 {
106 m_hasPlane0 = pFields[0]->GetZIDs()[0] == 0;
107 }
108 // initialize variables
111 m_body.extForceXYZ = Array<OneD, NekDouble>(6, 0.0);
112 m_hasVel = Array<OneD, bool>(3, false);
113 m_hasOmega = Array<OneD, bool>(3, false);
116 m_index = 0;
117 m_currentTime = -1.;
118 LoadParameters(pForce);
119 InitBodySolver(pForce);
120 if (m_body.isCircular || m_expdim == 1)
121 {
122 for (int i = 0; i < 3; ++i)
123 {
124 m_hasOmega[i] = false;
125 }
126 m_hasRotation = false;
127 }
128 // account for the effect of rotation
129 // Omega_X results in having v and w even if not defined by user
130 // Omega_Y results in having u and w even if not defined by user
131 // Omega_Z results in having u and v even if not defined by user
132 for (int i = 0; i < 3; ++i)
133 {
134 int j = (i + 1) % 3;
135 int k = (i + 2) % 3;
136 if (m_hasOmega[i])
137 {
138 m_hasVel[j] = true;
139 m_hasVel[k] = true;
140 }
141 }
142 if (m_hasRotation)
143 {
144 int npoints = pFields[0]->GetNpoints();
146 for (int j = 0; j < m_spacedim; ++j)
147 {
148 m_coords[j] = Array<OneD, NekDouble>(npoints);
149 }
150 pFields[0]->GetCoords(m_coords[0], m_coords[1], m_coords[2]);
151 // move the origin to the pivot point
152 for (int i = 0; i < m_spacedim; ++i)
153 {
154 Vmath::Sadd(npoints, -m_pivotPoint[i], m_coords[i], 1, m_coords[i],
155 1);
156 }
157 }
158 // initialise pivot point for fluid interface
159 auto equ = m_equ.lock();
160 ASSERTL0(equ, "Weak pointer to the equation system is expired");
161 auto FluidEq = std::dynamic_pointer_cast<FluidInterface>(equ);
162 FluidEq->SetMovingFramePivot(m_pivotPoint);
163 // initialize aeroforce filter
164 if (m_body.hasFreeMotion)
165 {
166 InitialiseFilter(m_session, pFields, pForce);
167 }
168}
169
171 std::string expression)
172{
173 NekDouble value = 0.;
174 try
175 {
176 LibUtilities::Equation expession(m_session->GetInterpreter(),
177 expression);
178 value = expession.Evaluate();
179 }
180 catch (const std::runtime_error &)
181 {
182 NEKERROR(ErrorUtil::efatal, "Error evaluating expression" + expression);
183 }
184 return value;
185}
186
187void ForcingMovingReferenceFrame::LoadParameters(const TiXmlElement *pForce)
188{
189 const TiXmlElement *funcNameElmt;
190 // body is a circular cylinder
191 const TiXmlElement *mssgTagSpecial =
192 pForce->FirstChildElement("CIRCULARCYLINDER");
193 if (mssgTagSpecial)
194 {
195 m_body.isCircular = true;
196 }
197 else
198 {
199 m_body.isCircular = false;
200 }
201 // load frame velocity
202 funcNameElmt = pForce->FirstChildElement("FRAMEVELOCITY");
203 if (funcNameElmt)
204 {
205 std::string FuncName = funcNameElmt->GetText();
206 ASSERTL0(m_session->DefinesFunction(FuncName),
207 "Function '" + FuncName + "' is not defined in the session.");
208 // linear velocity
209 for (int i = 0; i < m_spacedim; ++i)
210 {
211 std::string var = m_session->GetVariable(i);
212 if (m_session->DefinesFunction(FuncName, var))
213 {
214 m_frameVelFunction[i] = m_session->GetFunction(FuncName, var);
215 m_hasVel[i] = true;
216 }
217 }
218 // linear displacement and acceleration
219 std::vector<std::string> linearDispVar = {"X", "Y", "Z"};
220 std::vector<std::string> linearAcceVar = {"A_x", "A_y", "A_z"};
221 for (int i = 0; i < m_spacedim; ++i)
222 {
223 if (m_session->DefinesFunction(FuncName, linearDispVar[i]))
224 {
225 m_frameVelFunction[i + 6] =
226 m_session->GetFunction(FuncName, linearDispVar[i]);
227 }
228 if (m_session->DefinesFunction(FuncName, linearAcceVar[i]))
229 {
230 m_frameVelFunction[i + 12] =
231 m_session->GetFunction(FuncName, linearAcceVar[i]);
232 }
233 }
234 // angular velocities
235 m_hasRotation = false;
236 std::vector<std::string> angularVar = {"Omega_x", "Omega_y", "Omega_z"};
237 for (int i = 0; i < 3; ++i)
238 {
239 std::string var = angularVar[i];
240 if (m_session->DefinesFunction(FuncName, var))
241 {
242 m_frameVelFunction[i + 3] =
243 m_session->GetFunction(FuncName, var);
244 m_hasOmega[i] = true;
245 m_hasRotation = true;
246 }
247 }
248 // angular displacement and acceleration
249 std::vector<std::string> angularDispVar = {"Theta_x", "Theta_y",
250 "Theta_z"};
251 std::vector<std::string> angularAccepVar = {"DOmega_x", "DOmega_y",
252 "DOmega_z"};
253 for (int i = 0; i < 3; ++i)
254 {
255 if (m_session->DefinesFunction(FuncName, angularDispVar[i]))
256 {
257 m_frameVelFunction[i + 3 + 6] =
258 m_session->GetFunction(FuncName, angularDispVar[i]);
259 }
260 if (m_session->DefinesFunction(FuncName, angularAccepVar[i]))
261 {
262 m_frameVelFunction[i + 3 + 12] =
263 m_session->GetFunction(FuncName, angularAccepVar[i]);
264 }
265 }
266 // TODO: add the support for general rotation
267 for (int i = 0; i < 2; ++i)
268 {
269 ASSERTL0(!m_hasOmega[i], "Currently only Omega_z is supported");
270 }
271 }
272
273 // load external force
274 funcNameElmt = pForce->FirstChildElement("EXTERNALFORCE");
275 if (funcNameElmt)
276 {
277 std::string FuncName = funcNameElmt->GetText();
278 if (m_session->DefinesFunction(FuncName))
279 {
280 std::vector<std::string> forceVar = {"Fx", "Fy", "Fz"};
281 for (int i = 0; i < m_spacedim; ++i)
282 {
283 std::string var = forceVar[i];
284 if (m_session->DefinesFunction(FuncName, var))
285 {
286 m_body.extForceFunction[i] =
287 m_session->GetFunction(FuncName, var);
288 }
289 }
290 std::vector<std::string> momentVar = {"Mx", "My", "Mz"};
291 for (int i = 0; i < 3; ++i)
292 {
293 std::string var = momentVar[i];
294 if (m_session->DefinesFunction(FuncName, var))
295 {
296 m_body.extForceFunction[i + 3] =
297 m_session->GetFunction(FuncName, var);
298 }
299 }
300 }
301 }
302
303 // load pitching pivot
304 const TiXmlElement *pivotElmt = pForce->FirstChildElement("PIVOTPOINT");
305 if (pivotElmt)
306 {
307 std::vector<std::string> values;
308 std::string mssgStr = pivotElmt->GetText();
309 ParseUtils::GenerateVector(mssgStr, values);
310 for (int j = 0; j < m_spacedim; ++j)
311 {
312 m_pivotPoint[j] = EvaluateExpression(values[j]);
313 }
314 }
315
316 // load travelling wave speed
317 const TiXmlElement *TWSElmt =
318 pForce->FirstChildElement("TRAVELINGWAVESPEED");
319 if (TWSElmt)
320 {
321 std::vector<std::string> values;
322 std::string mssgStr = TWSElmt->GetText();
323 ParseUtils::GenerateVector(mssgStr, values);
324 for (int j = 0; j < m_spacedim; ++j)
325 {
326 NekDouble tmp = EvaluateExpression(values[j]);
327 if (fabs(tmp) > NekConstants::kNekMachineEpsilon)
328 {
329 m_travelWave[j] = tmp;
330 m_hasVel[j] = true;
331 }
332 }
333 }
334
335 // OutputFile
336 const TiXmlElement *mssgTag = pForce->FirstChildElement("OutputFile");
337 std::string filename;
338 if (mssgTag)
339 {
340 filename = mssgTag->GetText();
341 }
342 else
343 {
344 filename = m_session->GetSessionName();
345 }
346 if (!(filename.length() >= 4 &&
347 filename.substr(filename.length() - 4) == ".mrf"))
348 {
349 filename += ".mrf";
350 }
351 if (m_isRoot)
352 {
353 m_outputStream.open(filename.c_str());
354 if (m_spacedim == 2)
355 {
357 << "Variables = t, x, ux, ax, y, uy, ay, theta, omega, domega"
358 << std::endl;
359 }
360 else if (m_spacedim == 3)
361 {
362 m_outputStream << "Variables = t, x, ux, ax, y, uy, ay, z, uz, az, "
363 "theta, omega, domega"
364 << std::endl;
365 }
366 }
367
368 // output frequency
370 mssgTag = pForce->FirstChildElement("OutputFrequency");
371 if (mssgTag)
372 {
373 std::vector<std::string> values;
374 std::string mssgStr = mssgTag->GetText();
375 m_outputFrequency = round(EvaluateExpression(mssgStr));
376 }
378 "OutputFrequency should be greater than zero.");
379}
380
381void ForcingMovingReferenceFrame::InitBodySolver(const TiXmlElement *pForce)
382{
383 int NumDof = m_spacedim + 1;
384 const TiXmlElement *mssgTag;
385 std::string mssgStr;
386 // allocate memory and initialise
388 for (size_t i = 0; i < 3; ++i)
389 {
390 m_body.vel[i] = Array<OneD, NekDouble>(NumDof, 0.);
391 }
394 // read free motion DoFs
395 m_body.dirDoFs.clear();
396 for (int i = 0; i < NumDof; ++i)
397 {
398 m_body.dirDoFs.insert(i);
399 }
400 mssgTag = pForce->FirstChildElement("MOTIONPRESCRIBED");
401 if (mssgTag)
402 {
403 std::vector<std::string> values;
404 mssgStr = mssgTag->GetText();
405 ParseUtils::GenerateVector(mssgStr, values);
406 ASSERTL0(values.size() == NumDof,
407 "MOTIONPRESCRIBED vector should be of size " +
408 std::to_string(NumDof));
409 for (int i = 0; i < NumDof; ++i)
410 {
411 if (EvaluateExpression(values[i]) == 0)
412 {
413 m_body.dirDoFs.erase(i);
414 if (i < m_spacedim)
415 {
416 m_hasVel[i] = true;
417 }
418 else if (i == m_spacedim)
419 {
420 m_hasOmega[2] = true;
421 m_hasRotation = true;
422 }
423 }
424 }
425 }
426 m_body.hasFreeMotion = m_body.dirDoFs.size() < NumDof;
427 // read mass matrix
428 m_body.M = Array<OneD, NekDouble>(NumDof * NumDof, 0.);
429 mssgTag = pForce->FirstChildElement("MASS");
430 ASSERTL0(m_body.dirDoFs.size() == NumDof || mssgTag,
431 "Mass matrix is required.");
432 if (mssgTag)
433 {
434 std::vector<std::string> values;
435 mssgStr = mssgTag->GetText();
436 ParseUtils::GenerateVector(mssgStr, values);
437 ASSERTL0(values.size() == NumDof * NumDof,
438 "Mass matrix should be of size " + std::to_string(NumDof) +
439 "X" + std::to_string(NumDof));
440 int count = 0;
441 for (int i = 0; i < NumDof; ++i)
442 {
443 for (int j = 0; j < NumDof; ++j)
444 {
445 m_body.M[count] = EvaluateExpression(values[count]);
446 ++count;
447 }
448 }
449 }
450 // read damping matrix
451 m_body.C = Array<OneD, NekDouble>(NumDof * NumDof, 0.);
452 mssgTag = pForce->FirstChildElement("DAMPING");
453 if (mssgTag)
454 {
455 std::vector<std::string> values;
456 mssgStr = mssgTag->GetText();
457 ParseUtils::GenerateVector(mssgStr, values);
458 ASSERTL0(values.size() == NumDof * NumDof,
459 "Damping matrix should be of size " + std::to_string(NumDof) +
460 "X" + std::to_string(NumDof));
461 int count = 0;
462 for (int i = 0; i < NumDof; ++i)
463 {
464 for (int j = 0; j < NumDof; ++j)
465 {
466 m_body.C[count] = EvaluateExpression(values[count]);
467 ++count;
468 }
469 }
470 }
471 // read rigidity matrix
472 m_body.K = Array<OneD, NekDouble>(NumDof * NumDof, 0.);
473 mssgTag = pForce->FirstChildElement("RIGIDITY");
474 if (mssgTag)
475 {
476 std::vector<std::string> values;
477 mssgStr = mssgTag->GetText();
478 ParseUtils::GenerateVector(mssgStr, values);
479 ASSERTL0(values.size() == NumDof * NumDof,
480 "Rigidity matrix should be of size " + std::to_string(NumDof) +
481 "X" + std::to_string(NumDof));
482 int count = 0;
483 for (int i = 0; i < NumDof; ++i)
484 {
485 for (int j = 0; j < NumDof; ++j)
486 {
487 m_body.K[count] = EvaluateExpression(values[count]);
488 ++count;
489 }
490 }
491 }
492 // read Newmark Beta paramters
493 m_timestep = m_session->GetParameter("TimeStep");
494 NekDouble beta = 0.25;
495 NekDouble gamma = 0.51;
496 if (m_session->DefinesParameter("NewmarkBeta"))
497 {
498 beta = m_session->GetParameter("NewmarkBeta");
499 }
500 if (m_session->DefinesParameter("NewmarkGamma"))
501 {
502 gamma = m_session->GetParameter("NewmarkGamma");
503 }
505 m_body.K, m_body.dirDoFs);
506}
507
509 const NekDouble &time, std::map<int, NekDouble> &Dirs)
510{
511 int NumDof = m_spacedim + 1;
512 for (auto it : m_frameVelFunction)
513 {
514 if (it.first < 3)
515 {
516 Dirs[it.first] = it.second->Evaluate(0., 0., 0., time);
517 }
518 else if (it.first == 5)
519 {
520 Dirs[m_spacedim] = it.second->Evaluate(0., 0., 0., time);
521 }
522 else if (it.first < 9)
523 {
524 Dirs[NumDof + it.first - 6] = it.second->Evaluate(0., 0., 0., time);
525 }
526 else if (it.first == 11)
527 {
528 Dirs[NumDof + m_spacedim] = it.second->Evaluate(0., 0., 0., time);
529 }
530 else if (it.first < 15)
531 {
532 Dirs[(NumDof << 1) + it.first - 12] =
533 it.second->Evaluate(0., 0., 0., time);
534 }
535 else if (it.first == 17)
536 {
537 Dirs[(NumDof << 1) + m_spacedim] =
538 it.second->Evaluate(0., 0., 0., time);
539 }
540 }
541 for (auto i : m_body.dirDoFs)
542 {
543 if (Dirs.find(i) == Dirs.end())
544 {
545 Dirs[i] = 0.;
546 Dirs[i + NumDof] = 0.;
547 Dirs[i + (NumDof << 1)] = 0.;
548 }
549 }
550}
551/**
552 * @brief Updates the forcing array with the current required forcing.
553 * @param pFields
554 * @param time
555 */
558 const NekDouble &time)
559{
560 if (m_currentTime >= time)
561 {
562 return;
563 }
564 m_currentTime = time;
565 std::map<int, NekDouble> Dirs;
566 UpdatePrescribed(time, Dirs);
567 if (m_index == 0)
568 {
570 }
571 else
572 {
573 // compute the velocites whoes functions are provided in inertial frame
574 Array<OneD, NekDouble> forcebody(6, 0.); // fluid force
575 if (m_body.hasFreeMotion)
576 {
577 for (auto it : m_body.extForceFunction)
578 {
579 m_body.extForceXYZ[it.first] =
580 it.second->Evaluate(0., 0., 0., time);
581 }
582 m_body.aeroforceFilter->GetForces(pFields, NullNekDouble1DArray,
583 time);
584 auto equ = m_equ.lock();
585 ASSERTL0(equ, "Weak pointer to the equation system is expired");
586 auto FluidEq = std::dynamic_pointer_cast<FluidInterface>(equ);
587 FluidEq->GetAeroForce(forcebody);
588 }
589 SolveBodyMotion(m_body.vel, forcebody, Dirs);
590 }
591 if (m_isRoot && m_index % m_outputFrequency == 0)
592 {
593 m_outputStream << boost::format("%25.19e") % time << " ";
594 for (size_t i = 0; i < m_body.vel[0].size(); ++i)
595 {
596 m_outputStream << boost::format("%25.19e") % m_body.vel[0][i] << " "
597 << boost::format("%25.19e") % m_body.vel[1][i] << " "
598 << boost::format("%25.19e") % m_body.vel[2][i]
599 << " ";
600 }
601 m_outputStream << std::endl;
602 }
603 // extract values and transform to body frame
604 for (int i = 0; i < m_spacedim; ++i)
605 {
606 m_velxyz[i] = m_body.vel[1][i];
607 }
608 if (m_hasRotation)
609 {
610 m_omegaxyz[2] = m_body.vel[1][m_spacedim];
611 Array<OneD, NekDouble> angle(3, 0.);
612 angle[2] = m_body.vel[0][m_spacedim];
613 m_frame.SetAngle(angle);
615 }
616 // set displacements at the current time step
618 ++m_index;
619}
620
622 Array<OneD, Array<OneD, NekDouble>> &bodyVel, const int step)
623{
624 // set displacements at the current or next time step
625 Array<OneD, NekDouble> disp(6, 0.);
627 for (int i = 0; i < m_spacedim; ++i)
628 {
629 disp[i] = bodyVel[0][i];
630 }
631 if (!m_body.isCircular)
632 {
633 disp[5] = bodyVel[0][m_spacedim];
634 }
635 // to set the boundary condition of the next time step
636 // update the frame velocities and accelerations
637 for (int i = 0; i < m_spacedim; ++i)
638 {
639 vel[i] = bodyVel[1][i];
640 vel[i + 6] = bodyVel[2][i];
641 }
642 vel[5] = bodyVel[1][m_spacedim];
643 vel[11] = bodyVel[2][m_spacedim];
645 if (step && m_hasRotation)
646 {
647 m_frame.SetAngle(disp + 3);
649 m_frame.IneritalToBody(3, vel + 6, tmp = vel + 6);
650 }
651 auto equ = m_equ.lock();
652 ASSERTL0(equ, "Weak pointer to the equation system is expired");
653 auto FluidEq = std::dynamic_pointer_cast<FluidInterface>(equ);
654 FluidEq->SetMovingFrameVelocities(vel, step);
655 FluidEq->SetMovingFrameDisp(disp, step);
656}
657
660 const Array<OneD, NekDouble> &forcebody, std::map<int, NekDouble> &Dirs)
661{
662 if (!m_body.hasFreeMotion)
663 {
664 m_bodySolver.SolvePrescribed(bodyVel, Dirs);
665 }
666 else if (!m_hasRotation || Dirs.find(m_spacedim) != Dirs.end())
667 {
668 m_bodySolver.SolvePrescribed(bodyVel, Dirs);
669 Array<OneD, NekDouble> force(6, 0.), tmp;
670 if (m_hasRotation && Dirs.find(m_spacedim) != Dirs.end())
671 {
672 Array<OneD, NekDouble> angle(3, 0.);
673 angle[2] = bodyVel[0][m_spacedim];
674 m_frame.SetAngle(angle);
675 m_frame.BodyToInerital(3, forcebody, force);
676 m_frame.BodyToInerital(3, forcebody + 3, tmp = force + 3);
677 }
678 else
679 {
680 Vmath::Vcopy(6, forcebody, 1, force, 1);
681 }
682 for (int i = 0; i < m_spacedim; ++i)
683 {
684 force[i] += m_body.extForceXYZ[i];
685 }
686 force[m_spacedim] = force[5] + m_body.extForceXYZ[5];
687 m_bodySolver.SolveFree(bodyVel, force);
688 }
689 else
690 {
691 Array<OneD, Array<OneD, NekDouble>> tmpbodyVel(bodyVel.size());
692 for (size_t i = 0; i < bodyVel.size(); ++i)
693 {
694 tmpbodyVel[i] = Array<OneD, NekDouble>(bodyVel[i].size());
695 }
696 Array<OneD, NekDouble> angle(3, 0.);
697 angle[2] = bodyVel[0][m_spacedim];
698 for (int iter = 0; iter < 2; ++iter)
699 {
700 // copy initial condition
701 for (size_t i = 0; i < bodyVel.size(); ++i)
702 {
703 Vmath::Vcopy(bodyVel[i].size(), bodyVel[i], 1, tmpbodyVel[i],
704 1);
705 }
706 Array<OneD, NekDouble> force(6, 0.), tmp;
707 m_frame.SetAngle(angle);
708 m_frame.BodyToInerital(3, forcebody, force);
709 m_frame.BodyToInerital(3, forcebody + 3, tmp = force + 3);
710 for (int i = 0; i < m_spacedim; ++i)
711 {
712 force[i] += m_body.extForceXYZ[i];
713 }
714 force[m_spacedim] = force[5] + m_body.extForceXYZ[5];
715 m_bodySolver.Solve(tmpbodyVel, force, Dirs);
716 angle[2] = tmpbodyVel[0][m_spacedim];
717 }
718 // copy final results
719 for (size_t i = 0; i < bodyVel.size(); ++i)
720 {
721 Vmath::Vcopy(bodyVel[i].size(), tmpbodyVel[i], 1, bodyVel[i], 1);
722 }
723 }
724}
725
726/**
727 * @brief Adds the body force, -Omega X u.
728 * @param fields
729 * @param inarray
730 * @param outarray
731 * @param time
732 */
735 const Array<OneD, Array<OneD, NekDouble>> &inarray,
736 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble &time)
737{
738 // If there is no rotation, body force is zero,
739 // nothing needs to be done here.
740 if (m_hasRotation)
741 {
742 int npoints = fields[0]->GetNpoints();
743 addRotation(npoints, outarray, -1., inarray, outarray);
744 }
746}
747
748/**
749 * @brief Set velocity boundary condition at the next time step
750 */
752{
753 time += m_timestep;
754 // compute the velocities whose functions are provided in inertial frame
755 std::map<int, NekDouble> Dirs;
756 UpdatePrescribed(time, Dirs);
757 Array<OneD, Array<OneD, NekDouble>> bodyVel(m_body.vel.size());
758 for (size_t i = 0; i < m_body.vel.size(); ++i)
759 {
760 bodyVel[i] = Array<OneD, NekDouble>(m_body.vel[i].size());
761 Vmath::Vcopy(m_body.vel[i].size(), m_body.vel[i], 1, bodyVel[i], 1);
762 }
763 m_bodySolver.SolvePrescribed(bodyVel, Dirs);
764 // set displacements at the next time step
765 UpdateFluidInterface(bodyVel, 1);
766}
767
768/**
769 * @brief outarray = inarray0 + angVelScale Omega x inarray1
770 */
772 int nPnts, // number of points
773 const Array<OneD, Array<OneD, NekDouble>> &inarray0, NekDouble angVelScale,
774 const Array<OneD, Array<OneD, NekDouble>> &inarray1,
776{
777 ASSERTL0(&inarray1 != &outarray, "inarray1 and outarray "
778 "should not be the same.");
779
780 // TODO: In case of having support for all three components of Omega,
781 // they should be transformed into the rotating frame first!
782
783 // In case that the inarray0 and outarry are different, to avoid using
784 // un-initialized array, copy the array first
785 if (&inarray0 != &outarray)
786 {
787 ASSERTL0(inarray0.size() == outarray.size(),
788 "inarray0 and outarray must have same dimentions");
789 for (int i = 0; i < inarray0.size(); ++i)
790 {
791 Vmath::Vcopy(nPnts, inarray0[i], 1, outarray[i], 1);
792 }
793 }
794
795 if (m_spacedim >= 2 && m_hasOmega[2])
796 {
797 NekDouble cp = m_omegaxyz[2] * angVelScale;
798 NekDouble cm = -1. * cp;
799
800 Vmath::Svtvp(nPnts, cm, inarray1[1], 1, outarray[0], 1, outarray[0], 1);
801 Vmath::Svtvp(nPnts, cp, inarray1[0], 1, outarray[1], 1, outarray[1], 1);
802 }
803
804 if (m_spacedim == 3 && m_hasOmega[0])
805 {
806 NekDouble cp = m_omegaxyz[0] * angVelScale;
807 NekDouble cm = -1. * cp;
808
809 Vmath::Svtvp(nPnts, cp, inarray1[1], 1, outarray[2], 1, outarray[2], 1);
810 Vmath::Svtvp(nPnts, cm, inarray1[2], 1, outarray[1], 1, outarray[1], 1);
811 }
812
813 if (m_spacedim == 3 && m_hasOmega[1])
814 {
815 NekDouble cp = m_omegaxyz[1] * angVelScale;
816 NekDouble cm = -1. * cp;
817
818 Vmath::Svtvp(nPnts, cp, inarray1[2], 1, outarray[0], 1, outarray[0], 1);
819 Vmath::Svtvp(nPnts, cm, inarray1[0], 1, outarray[2], 1, outarray[2], 1);
820 }
821}
822
823/**
824 * @brief Compute the moving frame velocity at given time
825 */
828 const Array<OneD, Array<OneD, NekDouble>> &inarray,
829 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble &time)
830{
831 UpdateFrameVelocity(fields, time);
832 int npoints = fields[0]->GetNpoints();
833 if (m_isH2d && fields[0]->GetWaveSpace())
834 {
835 for (int i = 0; i < m_spacedim; ++i)
836 {
837 if (m_hasVel[i])
838 {
839 Array<OneD, NekDouble> tmpphys(npoints,
840 -m_velxyz[i] - m_travelWave[i]);
841 Array<OneD, NekDouble> tmpcoef(npoints);
842 fields[0]->HomogeneousFwdTrans(npoints, tmpphys, tmpcoef);
843 Vmath::Vadd(npoints, tmpcoef, 1, inarray[i], 1, outarray[i], 1);
844 }
845 else if (&inarray != &outarray)
846 {
847 Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
848 }
849 }
850 }
851 else
852 {
853 int npoints0 = npoints;
854 if (m_isH1d && fields[0]->GetWaveSpace())
855 {
856 npoints0 = m_hasPlane0 ? fields[0]->GetPlane(0)->GetNpoints() : 0;
857 }
858 for (int i = 0; i < m_spacedim; ++i)
859 {
860 if (m_hasVel[i])
861 {
862 Vmath::Sadd(npoints0, -m_velxyz[i] - m_travelWave[i],
863 inarray[i], 1, outarray[i], 1);
864 if (&inarray != &outarray && npoints != npoints0)
865 {
866 Array<OneD, NekDouble> tmp = outarray[i] + npoints0;
867 Vmath::Vcopy(npoints - npoints0, inarray[i] + npoints0, 1,
868 tmp, 1);
869 }
870 }
871 else if (&inarray != &outarray)
872 {
873 Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
874 }
875 }
876 if (m_hasRotation)
877 {
878 addRotation(npoints0, outarray, -1., m_coords, outarray);
879 }
880 }
881}
882
884{
885 NekDouble time = 0.;
886 std::map<std::string, std::string> fieldMetaDataMap;
887 std::vector<std::string> strFrameData = {
888 "X", "Y", "Z", "Theta_x", "Theta_y", "Theta_z",
889 "U", "V", "W", "Omega_x", "Omega_y", "Omega_z",
890 "A_x", "A_y", "A_z", "DOmega_x", "DOmega_y", "DOmega_z"};
891 std::map<std::string, NekDouble> fileData;
892 if (m_session->DefinesFunction("InitialConditions"))
893 {
894 for (int i = 0; i < m_session->GetVariables().size(); ++i)
895 {
896 if (m_session->GetFunctionType("InitialConditions",
897 m_session->GetVariable(i)) ==
899 {
900 std::string filename = m_session->GetFunctionFilename(
901 "InitialConditions", m_session->GetVariable(i));
902 fs::path pfilename(filename);
903 // redefine path for parallel file which is in directory
904 if (fs::is_directory(pfilename))
905 {
906 fs::path metafile("Info.xml");
907 fs::path fullpath = pfilename / metafile;
908 filename = LibUtilities::PortablePath(fullpath);
909 }
912 fld->ImportFieldMetaData(filename, fieldMetaDataMap);
913
914 // check to see if time is defined
915 if (fieldMetaDataMap != LibUtilities::NullFieldMetaDataMap)
916 {
917 if (fieldMetaDataMap.find("Time") != fieldMetaDataMap.end())
918 {
919 time = std::stod(fieldMetaDataMap["Time"]);
920 }
921 fileData.clear();
922 for (auto &var : strFrameData)
923 {
924 if (fieldMetaDataMap.find(var) !=
925 fieldMetaDataMap.end())
926 {
927 fileData[var] = std::stod(fieldMetaDataMap[var]);
928 }
929 }
930 if (fileData.size() == strFrameData.size())
931 {
932 break;
933 }
934 }
935 }
936 }
937 }
938 if (m_session->DefinesCmdLineArgument("set-start-time"))
939 {
940 time = std::stod(
941 m_session->GetCmdLineArgument<std::string>("set-start-time"));
942 }
943 if (fileData.size() == strFrameData.size())
944 {
945 int NumDofm1 = m_body.vel[0].size() - 1;
946 for (int i = 0; i < m_spacedim; ++i)
947 {
948 m_body.vel[0][i] = fileData[strFrameData[i]];
949 m_body.vel[1][i] = fileData[strFrameData[i + 6]];
950 m_body.vel[2][i] = fileData[strFrameData[i + 12]];
951 }
952 m_body.vel[0][NumDofm1] = fileData[strFrameData[5]];
953 m_body.vel[1][NumDofm1] = fileData[strFrameData[11]];
954 m_body.vel[2][NumDofm1] = fileData[strFrameData[17]];
955 }
956 std::map<int, NekDouble> Dirs;
957 UpdatePrescribed(time, Dirs);
959}
960
962 std::map<int, NekDouble> &Dirs)
963{
964 for (auto it : Dirs)
965 {
966 int NumDof = m_body.vel[0].size();
967 int NumDof2 = NumDof << 1;
968 if (it.first < m_body.vel[0].size())
969 {
970 m_body.vel[1][it.first] = it.second;
971 }
972 else if (it.first < NumDof2)
973 {
974 m_body.vel[0][it.first - NumDof] = it.second;
975 }
976 else
977 {
978 m_body.vel[2][it.first - NumDof2] = it.second;
979 }
980 }
981}
982
986 const TiXmlElement *pForce)
987{
988 std::map<std::string, std::string> vParams;
989 vParams["OutputFile"] = ".dummyMRFForceFile";
990 const TiXmlElement *param = pForce->FirstChildElement("BOUNDARY");
991 ASSERTL0(param, "Body surface should be assigned");
992
993 vParams["Boundary"] = param->GetText();
994 const TiXmlElement *pivotElmt = pForce->FirstChildElement("PIVOTPOINT");
995 if (pivotElmt)
996 {
997 std::string pstr = pivotElmt->GetText();
998 std::replace(pstr.begin(), pstr.end(), ',', ' ');
999 vParams["MomentPoint"] = pstr;
1000 }
1002 pSession, m_equ.lock(), vParams);
1003
1004 m_body.aeroforceFilter->Initialise(pFields, 0.0);
1005}
1006
1011 std::set<int> DirDoFs)
1012{
1014 m_coeffs[0] = 1. / (gamma * dt);
1015 m_coeffs[1] = 1. / gamma - 1.;
1016 m_coeffs[2] = beta * dt / gamma;
1017 m_coeffs[3] = dt * (1. - beta / gamma);
1018 m_coeffs[4] = (0.5 - beta / gamma) * dt * dt;
1019
1020 m_rows = sqrt(M.size());
1021 m_index.resize(m_rows, -1);
1022 m_motionDofs = 0;
1023 for (int i = 0; i < m_rows; ++i)
1024 {
1025 if (DirDoFs.find(i) == DirDoFs.end())
1026 {
1027 m_index[m_motionDofs++] = i;
1028 }
1029 }
1030 for (int i = 0, count = m_motionDofs; i < m_rows; ++i)
1031 {
1032 if (DirDoFs.find(i) != DirDoFs.end())
1033 {
1034 m_index[count++] = i;
1035 }
1036 }
1037 if (m_motionDofs)
1038 {
1044 DNekMatSharedPtr inverseMatrix =
1046 m_motionDofs, 0.0, eFULL);
1047 for (int i = 0; i < m_motionDofs; ++i)
1048 {
1053 int offset = m_index[i] * m_rows;
1054 for (int j = 0; j < m_rows; ++j)
1055 {
1056 int ind = offset + m_index[j];
1057 m_M[i][j] = M[ind];
1058 m_C[i][j] = C[ind];
1059 m_K[i][j] = K[ind];
1060 NekDouble value =
1061 m_coeffs[0] * M[ind] + C[ind] + m_coeffs[2] * K[ind];
1062 if (j < m_motionDofs)
1063 {
1064 inverseMatrix->SetValue(i, j, value);
1065 }
1066 }
1067 }
1068 inverseMatrix->Invert();
1069 for (int i = 0; i < m_motionDofs; ++i)
1070 {
1071 for (int j = 0; j < m_rows; ++j)
1072 {
1073 if (j < m_motionDofs)
1074 {
1075 m_Matrix[i][j] = inverseMatrix->GetValue(i, j);
1076 }
1077 }
1078 }
1079 }
1080}
1081
1084 std::map<int, NekDouble> motionPrescribed)
1085{
1086 for (int i = 0; i < m_rows; ++i)
1087 {
1088 if (motionPrescribed.find(i) != motionPrescribed.end())
1089 {
1090 int i0 = i + m_rows;
1091 int i2 = i0 + m_rows;
1092 NekDouble bm = 0., bk = 0.;
1093 if (motionPrescribed.find(i2) == motionPrescribed.end())
1094 {
1095 bm = m_coeffs[0] * u[1][i] + m_coeffs[1] * u[2][i];
1096 }
1097 if (motionPrescribed.find(i0) == motionPrescribed.end())
1098 {
1099 bk = u[0][i] + m_coeffs[3] * u[1][i] + m_coeffs[4] * u[2][i];
1100 }
1101
1102 u[1][i] = motionPrescribed[i];
1103 if (motionPrescribed.find(i2) == motionPrescribed.end())
1104 {
1105 u[2][i] = m_coeffs[0] * u[1][i] - bm;
1106 }
1107 else
1108 {
1109 u[2][i] = motionPrescribed[i2];
1110 }
1111 if (motionPrescribed.find(i0) == motionPrescribed.end())
1112 {
1113 u[0][i] = m_coeffs[2] * u[1][i] + bk;
1114 }
1115 else
1116 {
1117 u[0][i] = motionPrescribed[i0];
1118 }
1119 }
1120 }
1121}
1122
1125 std::map<int, NekDouble> motionPrescribed)
1126{
1127 SolvePrescribed(u, motionPrescribed);
1128 SolveFree(u, force);
1129}
1130
1133{
1134 if (m_motionDofs)
1135 {
1138 for (int j = 0; j < m_motionDofs; ++j)
1139 {
1140 int j1 = m_index[j];
1141 bm[j] = m_coeffs[0] * u[1][j1] + m_coeffs[1] * u[2][j1];
1142 bk[j] = u[0][j1] + m_coeffs[3] * u[1][j1] + m_coeffs[4] * u[2][j1];
1143 }
1145 for (int i = 0; i < m_motionDofs; ++i)
1146 {
1147 rhs[i] = force[m_index[i]];
1148 for (int j = 0; j < m_motionDofs; ++j)
1149 {
1150 rhs[i] += m_M[i][j] * bm[j] - m_K[i][j] * bk[j];
1151 }
1152 for (int j = m_motionDofs; j < m_rows; ++j)
1153 {
1154 int j1 = m_index[j];
1155 rhs[i] -= m_M[i][j] * u[2][j1] + m_C[i][j] * u[1][j1] +
1156 m_K[i][j] * u[0][j1];
1157 }
1158 }
1159 for (int j = 0; j < m_motionDofs; ++j)
1160 {
1161 int j1 = m_index[j];
1162 u[1][j1] = Vmath::Dot(m_motionDofs, m_Matrix[j], 1, rhs, 1);
1163 u[0][j1] = m_coeffs[2] * u[1][j1] + bk[j];
1164 u[2][j1] = m_coeffs[0] * u[1][j1] - bm[j];
1165 }
1166 }
1167}
1168
1170{
1172}
1173
1175{
1176 m_matrix[0] = cos(theta[2]);
1177 m_matrix[1] = sin(theta[2]);
1178}
1179
1181 const Array<OneD, NekDouble> &body,
1182 Array<OneD, NekDouble> &inertial)
1183{
1184 NekDouble xi = body[0], eta = body[1];
1185 inertial[0] = xi * m_matrix[0] - eta * m_matrix[1];
1186 inertial[1] = eta * m_matrix[0] + xi * m_matrix[1];
1187 if (body != inertial && dim > 2)
1188 {
1189 for (int i = 2; i < dim; ++i)
1190 {
1191 inertial[i] = body[i];
1192 }
1193 }
1194}
1195
1197 const Array<OneD, NekDouble> &inertial,
1199{
1200 NekDouble x = inertial[0], y = inertial[1];
1201 body[0] = x * m_matrix[0] + y * m_matrix[1];
1202 body[1] = y * m_matrix[0] - x * m_matrix[1];
1203 if (body != inertial && dim > 2)
1204 {
1205 for (int i = 2; i < dim; ++i)
1206 {
1207 body[i] = inertial[i];
1208 }
1209 }
1210}
1211
1212} // namespace Nektar::SolverUtils
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
static std::shared_ptr< FieldIO > CreateForFile(const LibUtilities::SessionReaderSharedPtr session, const std::string &filename)
Construct a FieldIO object for a given input filename.
Definition: FieldIO.cpp:223
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
static bool GenerateVector(const std::string &str, std::vector< T > &out)
Takes a comma-separated string and converts it to entries in a vector.
Definition: ParseUtils.cpp:130
Defines a forcing term to be explicitly applied.
Definition: Forcing.h:71
const std::weak_ptr< EquationSystem > m_equ
Weak pointer to equation system using this forcing.
Definition: Forcing.h:125
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
Definition: Forcing.h:123
SOLVER_UTILS_EXPORT void v_PreApply(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time) override
Compute the moving frame velocity at given time.
void addRotation(int npoints, const Array< OneD, Array< OneD, NekDouble > > &inarray0, NekDouble angVelScale, const Array< OneD, Array< OneD, NekDouble > > &inarray1, Array< OneD, Array< OneD, NekDouble > > &outarray)
outarray = inarray0 + angVelScale Omega x inarray1
std::map< int, LibUtilities::EquationSharedPtr > m_frameVelFunction
void SolveBodyMotion(Array< OneD, Array< OneD, NekDouble > > &bodyVel, const Array< OneD, NekDouble > &forcebody, std::map< int, NekDouble > &Dirs)
void UpdateFluidInterface(Array< OneD, Array< OneD, NekDouble > > &bodyVel, const int step)
void UpdatePrescribed(const NekDouble &time, std::map< int, NekDouble > &Dirs)
ForcingMovingReferenceFrame(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation)
SOLVER_UTILS_EXPORT 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
Adds the body force, -Omega X u.
SOLVER_UTILS_EXPORT void v_InitObject(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields, const TiXmlElement *pForce) override
Initialise the forcing module.
void UpdateFrameVelocity(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
Updates the forcing array with the current required forcing.
void UpdateBoundaryConditions(NekDouble time)
Set velocity boundary condition at the next time step.
static SOLVER_UTILS_EXPORT ForcingSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields, const TiXmlElement *pForce)
Creates an instance of this class.
struct Nektar::SolverUtils::ForcingMovingReferenceFrame::@2 m_body
void InitialiseFilter(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pForce)
void SetAngle(const Array< OneD, NekDouble > theta)
void BodyToInerital(const int dim, const Array< OneD, NekDouble > &body, Array< OneD, NekDouble > &inertial)
void IneritalToBody(const int dim, const Array< OneD, NekDouble > &inertial, Array< OneD, NekDouble > &body)
void SetNewmarkBeta(NekDouble beta, NekDouble gamma, NekDouble dt, Array< OneD, NekDouble > M, Array< OneD, NekDouble > C, Array< OneD, NekDouble > K, std::set< int > DirDoFs)
Array< OneD, Array< OneD, NekDouble > > m_K
void SolveFree(Array< OneD, Array< OneD, NekDouble > > u, Array< OneD, NekDouble > force)
void SolvePrescribed(Array< OneD, Array< OneD, NekDouble > > u, std::map< int, NekDouble > motionPrescribed)
Array< OneD, Array< OneD, NekDouble > > m_Matrix
Array< OneD, Array< OneD, NekDouble > > m_C
void Solve(Array< OneD, Array< OneD, NekDouble > > u, Array< OneD, NekDouble > force, std::map< int, NekDouble > motionPrescribed)
Array< OneD, Array< OneD, NekDouble > > m_M
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:322
static std::string PortablePath(const fs::path &path)
create portable path on different platforms for std::filesystem path.
Definition: Filesystem.hpp:66
std::shared_ptr< SessionReader > SessionReaderSharedPtr
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:51
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59
static const NekDouble kNekMachineEpsilon
ForcingFactory & GetForcingFactory()
Declaration of the forcing factory singleton.
Definition: Forcing.cpp:41
static Array< OneD, NekDouble > NullNekDouble1DArray
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
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.hpp:396
T Dot(int n, const T *w, const T *x)
dot product
Definition: Vmath.hpp:761
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.hpp:180
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.hpp:194
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:285