50#include <boost/format.hpp>
67 const std::weak_ptr<EquationSystem> &pEquation)
88 [[maybe_unused]]
const unsigned int &pNumForcingFields,
89 const TiXmlElement *pForce)
93 bool singleMode, halfMode;
94 m_session->MatchSolverInfo(
"ModeType",
"SingleMode", singleMode,
false);
95 m_session->MatchSolverInfo(
"ModeType",
"HalfMode", halfMode,
false);
96 if (singleMode || halfMode)
102 m_isRoot = pFields[0]->GetComm()->TreatAsRankZero();
122 for (
int i = 0; i < 3; ++i)
132 for (
int i = 0; i < 3; ++i)
144 int npoints = pFields[0]->GetNpoints();
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);
171 std::string expression)
180 catch (
const std::runtime_error &)
189 const TiXmlElement *funcNameElmt;
191 const TiXmlElement *mssgTagSpecial =
192 pForce->FirstChildElement(
"CIRCULARCYLINDER");
199 m_body.isCircular =
false;
202 funcNameElmt = pForce->FirstChildElement(
"FRAMEVELOCITY");
205 std::string FuncName = funcNameElmt->GetText();
207 "Function '" + FuncName +
"' is not defined in the session.");
211 std::string var =
m_session->GetVariable(i);
212 if (
m_session->DefinesFunction(FuncName, var))
219 std::vector<std::string> linearDispVar = {
"X",
"Y",
"Z"};
220 std::vector<std::string> linearAcceVar = {
"A_x",
"A_y",
"A_z"};
223 if (
m_session->DefinesFunction(FuncName, linearDispVar[i]))
226 m_session->GetFunction(FuncName, linearDispVar[i]);
228 if (
m_session->DefinesFunction(FuncName, linearAcceVar[i]))
231 m_session->GetFunction(FuncName, linearAcceVar[i]);
236 std::vector<std::string> angularVar = {
"Omega_x",
"Omega_y",
"Omega_z"};
237 for (
int i = 0; i < 3; ++i)
239 std::string var = angularVar[i];
240 if (
m_session->DefinesFunction(FuncName, var))
249 std::vector<std::string> angularDispVar = {
"Theta_x",
"Theta_y",
251 std::vector<std::string> angularAccepVar = {
"DOmega_x",
"DOmega_y",
253 for (
int i = 0; i < 3; ++i)
255 if (
m_session->DefinesFunction(FuncName, angularDispVar[i]))
258 m_session->GetFunction(FuncName, angularDispVar[i]);
260 if (
m_session->DefinesFunction(FuncName, angularAccepVar[i]))
263 m_session->GetFunction(FuncName, angularAccepVar[i]);
267 for (
int i = 0; i < 2; ++i)
274 funcNameElmt = pForce->FirstChildElement(
"EXTERNALFORCE");
277 std::string FuncName = funcNameElmt->GetText();
278 if (
m_session->DefinesFunction(FuncName))
280 std::vector<std::string> forceVar = {
"Fx",
"Fy",
"Fz"};
283 std::string var = forceVar[i];
284 if (
m_session->DefinesFunction(FuncName, var))
286 m_body.extForceFunction[i] =
290 std::vector<std::string> momentVar = {
"Mx",
"My",
"Mz"};
291 for (
int i = 0; i < 3; ++i)
293 std::string var = momentVar[i];
294 if (
m_session->DefinesFunction(FuncName, var))
296 m_body.extForceFunction[i + 3] =
304 const TiXmlElement *pivotElmt = pForce->FirstChildElement(
"PIVOTPOINT");
307 std::vector<std::string> values;
308 std::string mssgStr = pivotElmt->GetText();
317 const TiXmlElement *TWSElmt =
318 pForce->FirstChildElement(
"TRAVELINGWAVESPEED");
321 std::vector<std::string> values;
322 std::string mssgStr = TWSElmt->GetText();
336 const TiXmlElement *mssgTag = pForce->FirstChildElement(
"OutputFile");
337 std::string filename;
340 filename = mssgTag->GetText();
346 if (!(filename.length() >= 4 &&
347 filename.substr(filename.length() - 4) ==
".mrf"))
357 <<
"Variables = t, x, ux, ax, y, uy, ay, theta, omega, domega"
362 m_outputStream <<
"Variables = t, x, ux, ax, y, uy, ay, z, uz, az, "
363 "theta, omega, domega"
370 mssgTag = pForce->FirstChildElement(
"OutputFrequency");
373 std::vector<std::string> values;
374 std::string mssgStr = mssgTag->GetText();
378 "OutputFrequency should be greater than zero.");
384 const TiXmlElement *mssgTag;
388 for (
size_t i = 0; i < 3; ++i)
396 for (
int i = 0; i < NumDof; ++i)
400 mssgTag = pForce->FirstChildElement(
"MOTIONPRESCRIBED");
403 std::vector<std::string> values;
404 mssgStr = mssgTag->GetText();
407 "MOTIONPRESCRIBED vector should be of size " +
408 std::to_string(NumDof));
409 for (
int i = 0; i < NumDof; ++i)
429 mssgTag = pForce->FirstChildElement(
"MASS");
431 "Mass matrix is required.");
434 std::vector<std::string> values;
435 mssgStr = mssgTag->GetText();
437 ASSERTL0(values.size() == NumDof * NumDof,
438 "Mass matrix should be of size " + std::to_string(NumDof) +
439 "X" + std::to_string(NumDof));
441 for (
int i = 0; i < NumDof; ++i)
443 for (
int j = 0; j < NumDof; ++j)
452 mssgTag = pForce->FirstChildElement(
"DAMPING");
455 std::vector<std::string> values;
456 mssgStr = mssgTag->GetText();
458 ASSERTL0(values.size() == NumDof * NumDof,
459 "Damping matrix should be of size " + std::to_string(NumDof) +
460 "X" + std::to_string(NumDof));
462 for (
int i = 0; i < NumDof; ++i)
464 for (
int j = 0; j < NumDof; ++j)
473 mssgTag = pForce->FirstChildElement(
"RIGIDITY");
476 std::vector<std::string> values;
477 mssgStr = mssgTag->GetText();
479 ASSERTL0(values.size() == NumDof * NumDof,
480 "Rigidity matrix should be of size " + std::to_string(NumDof) +
481 "X" + std::to_string(NumDof));
483 for (
int i = 0; i < NumDof; ++i)
485 for (
int j = 0; j < NumDof; ++j)
496 if (
m_session->DefinesParameter(
"NewmarkBeta"))
500 if (
m_session->DefinesParameter(
"NewmarkGamma"))
502 gamma =
m_session->GetParameter(
"NewmarkGamma");
509 const NekDouble &time, std::map<int, NekDouble> &Dirs)
516 Dirs[it.first] = it.second->Evaluate(0., 0., 0., time);
518 else if (it.first == 5)
520 Dirs[
m_spacedim] = it.second->Evaluate(0., 0., 0., time);
522 else if (it.first < 9)
524 Dirs[NumDof + it.first - 6] = it.second->Evaluate(0., 0., 0., time);
526 else if (it.first == 11)
528 Dirs[NumDof +
m_spacedim] = it.second->Evaluate(0., 0., 0., time);
530 else if (it.first < 15)
532 Dirs[(NumDof << 1) + it.first - 12] =
533 it.second->Evaluate(0., 0., 0., time);
535 else if (it.first == 17)
538 it.second->Evaluate(0., 0., 0., time);
541 for (
auto i :
m_body.dirDoFs)
543 if (Dirs.find(i) == Dirs.end())
546 Dirs[i + NumDof] = 0.;
547 Dirs[i + (NumDof << 1)] = 0.;
565 std::map<int, NekDouble> Dirs;
577 for (
auto it :
m_body.extForceFunction)
579 m_body.extForceXYZ[it.first] =
580 it.second->Evaluate(0., 0., 0., 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);
594 for (
size_t i = 0; i <
m_body.vel[0].size(); ++i)
629 disp[i] = bodyVel[0][i];
639 vel[i] = bodyVel[1][i];
640 vel[i + 6] = bodyVel[2][i];
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);
662 if (!
m_body.hasFreeMotion)
684 force[i] +=
m_body.extForceXYZ[i];
692 for (
size_t i = 0; i < bodyVel.size(); ++i)
698 for (
int iter = 0; iter < 2; ++iter)
701 for (
size_t i = 0; i < bodyVel.size(); ++i)
703 Vmath::Vcopy(bodyVel[i].size(), bodyVel[i], 1, tmpbodyVel[i],
712 force[i] +=
m_body.extForceXYZ[i];
719 for (
size_t i = 0; i < bodyVel.size(); ++i)
721 Vmath::Vcopy(bodyVel[i].size(), tmpbodyVel[i], 1, bodyVel[i], 1);
742 int npoints = fields[0]->GetNpoints();
743 addRotation(npoints, outarray, -1., inarray, outarray);
755 std::map<int, NekDouble> Dirs;
758 for (
size_t i = 0; i <
m_body.vel.size(); ++i)
777 ASSERTL0(&inarray1 != &outarray,
"inarray1 and outarray "
778 "should not be the same.");
785 if (&inarray0 != &outarray)
787 ASSERTL0(inarray0.size() == outarray.size(),
788 "inarray0 and outarray must have same dimentions");
789 for (
int i = 0; i < inarray0.size(); ++i)
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);
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);
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);
832 int npoints = fields[0]->GetNpoints();
833 if (
m_isH2d && fields[0]->GetWaveSpace())
842 fields[0]->HomogeneousFwdTrans(npoints, tmpphys, tmpcoef);
843 Vmath::Vadd(npoints, tmpcoef, 1, inarray[i], 1, outarray[i], 1);
845 else if (&inarray != &outarray)
853 int npoints0 = npoints;
854 if (
m_isH1d && fields[0]->GetWaveSpace())
856 npoints0 =
m_hasPlane0 ? fields[0]->GetPlane(0)->GetNpoints() : 0;
863 inarray[i], 1, outarray[i], 1);
864 if (&inarray != &outarray && npoints != npoints0)
867 Vmath::Vcopy(npoints - npoints0, inarray[i] + npoints0, 1,
871 else if (&inarray != &outarray)
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"))
894 for (
int i = 0; i <
m_session->GetVariables().size(); ++i)
896 if (
m_session->GetFunctionType(
"InitialConditions",
900 std::string filename =
m_session->GetFunctionFilename(
901 "InitialConditions",
m_session->GetVariable(i));
902 fs::path pfilename(filename);
904 if (fs::is_directory(pfilename))
906 fs::path metafile(
"Info.xml");
907 fs::path fullpath = pfilename / metafile;
912 fld->ImportFieldMetaData(filename, fieldMetaDataMap);
917 if (fieldMetaDataMap.find(
"Time") != fieldMetaDataMap.end())
919 time = std::stod(fieldMetaDataMap[
"Time"]);
922 for (
auto &var : strFrameData)
924 if (fieldMetaDataMap.find(var) !=
925 fieldMetaDataMap.end())
927 fileData[var] = std::stod(fieldMetaDataMap[var]);
930 if (fileData.size() == strFrameData.size())
938 if (
m_session->DefinesCmdLineArgument(
"set-start-time"))
941 m_session->GetCmdLineArgument<std::string>(
"set-start-time"));
943 if (fileData.size() == strFrameData.size())
945 int NumDofm1 =
m_body.vel[0].size() - 1;
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]];
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]];
956 std::map<int, NekDouble> Dirs;
962 std::map<int, NekDouble> &Dirs)
966 int NumDof =
m_body.vel[0].size();
967 int NumDof2 = NumDof << 1;
968 if (it.first <
m_body.vel[0].size())
970 m_body.vel[1][it.first] = it.second;
972 else if (it.first < NumDof2)
974 m_body.vel[0][it.first - NumDof] = it.second;
978 m_body.vel[2][it.first - NumDof2] = it.second;
986 const TiXmlElement *pForce)
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");
993 vParams[
"Boundary"] = param->GetText();
994 const TiXmlElement *pivotElmt = pForce->FirstChildElement(
"PIVOTPOINT");
997 std::string pstr = pivotElmt->GetText();
998 std::replace(pstr.begin(), pstr.end(),
',',
' ');
999 vParams[
"MomentPoint"] = pstr;
1002 pSession,
m_equ.lock(), vParams);
1004 m_body.aeroforceFilter->Initialise(pFields, 0.0);
1011 std::set<int> DirDoFs)
1023 for (
int i = 0; i <
m_rows; ++i)
1025 if (DirDoFs.find(i) == DirDoFs.end())
1032 if (DirDoFs.find(i) != DirDoFs.end())
1054 for (
int j = 0; j <
m_rows; ++j)
1056 int ind = offset +
m_index[j];
1064 inverseMatrix->SetValue(i, j, value);
1068 inverseMatrix->Invert();
1071 for (
int j = 0; j <
m_rows; ++j)
1075 m_Matrix[i][j] = inverseMatrix->GetValue(i, j);
1084 std::map<int, NekDouble> motionPrescribed)
1086 for (
int i = 0; i <
m_rows; ++i)
1088 if (motionPrescribed.find(i) != motionPrescribed.end())
1093 if (motionPrescribed.find(i2) == motionPrescribed.end())
1097 if (motionPrescribed.find(i0) == motionPrescribed.end())
1102 u[1][i] = motionPrescribed[i];
1103 if (motionPrescribed.find(i2) == motionPrescribed.end())
1105 u[2][i] =
m_coeffs[0] * u[1][i] - bm;
1109 u[2][i] = motionPrescribed[i2];
1111 if (motionPrescribed.find(i0) == motionPrescribed.end())
1113 u[0][i] =
m_coeffs[2] * u[1][i] + bk;
1117 u[0][i] = motionPrescribed[i0];
1125 std::map<int, NekDouble> motionPrescribed)
1150 rhs[i] +=
m_M[i][j] * bm[j] -
m_K[i][j] * bk[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];
1163 u[0][j1] =
m_coeffs[2] * u[1][j1] + bk[j];
1164 u[2][j1] =
m_coeffs[0] * u[1][j1] - bm[j];
1187 if (body != inertial && dim > 2)
1189 for (
int i = 2; i < dim; ++i)
1191 inertial[i] = body[i];
1200 NekDouble x = inertial[0], y = inertial[1];
1203 if (body != inertial && dim > 2)
1205 for (
int i = 2; i < dim; ++i)
1207 body[i] = inertial[i];
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
NekDouble Evaluate() const
static std::shared_ptr< FieldIO > CreateForFile(const LibUtilities::SessionReaderSharedPtr session, const std::string &filename)
Construct a FieldIO object for a given input filename.
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.
Defines a forcing term to be explicitly applied.
const std::weak_ptr< EquationSystem > m_equ
Weak pointer to equation system using this forcing.
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
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 InitBodySolver(const TiXmlElement *pForce)
Array< OneD, NekDouble > m_velxyz
static std::string classNameBody
Name of the class.
Array< OneD, bool > m_hasVel
~ForcingMovingReferenceFrame(void) override
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
NekDouble EvaluateExpression(std::string expression)
Newmark_BetaSolver m_bodySolver
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)
unsigned int m_outputFrequency
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)
Array< OneD, NekDouble > m_omegaxyz
void SetInitialConditions()
std::ofstream m_outputStream
Array< OneD, NekDouble > m_pivotPoint
Array< OneD, Array< OneD, NekDouble > > vel
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.
void LoadParameters(const TiXmlElement *pForce)
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
Array< OneD, bool > m_hasOmega
Array< OneD, Array< OneD, NekDouble > > m_coords
Array< OneD, NekDouble > m_travelWave
void InitialiseFilter(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pForce)
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)
std::vector< int > m_index
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, NekDouble > m_coeffs
Array< OneD, Array< OneD, NekDouble > > m_M
std::shared_ptr< FieldIO > FieldIOSharedPtr
static std::string PortablePath(const fs::path &path)
create portable path on different platforms for std::filesystem path.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
static FieldMetaDataMap NullFieldMetaDataMap
@ beta
Gauss Radau pinned at x=-1,.
static const NekDouble kNekMachineEpsilon
ForcingFactory & GetForcingFactory()
Declaration of the forcing factory singleton.
static Array< OneD, NekDouble > NullNekDouble1DArray
std::shared_ptr< DNekMat > DNekMatSharedPtr
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.
T Dot(int n, const T *w, const T *x)
dot product
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.
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
scalarT< T > sqrt(scalarT< T > in)