50#include <boost/format.hpp>
68 const std::weak_ptr<EquationSystem> &pEquation)
89 [[maybe_unused]]
const unsigned int &pNumForcingFields,
90 const TiXmlElement *pForce)
94 bool singleMode, halfMode;
95 m_session->MatchSolverInfo(
"ModeType",
"SingleMode", singleMode,
false);
96 m_session->MatchSolverInfo(
"ModeType",
"HalfMode", halfMode,
false);
97 if (singleMode || halfMode)
102 int npoints = pFields[0]->GetNpoints();
148 m_rank = pFields[0]->GetComm()->GetRank();
154 auto equ =
m_equ.lock();
155 ASSERTL0(equ,
"Weak pointer to the equation system is expired");
156 auto FluidEq = std::dynamic_pointer_cast<FluidInterface>(equ);
157 FluidEq->SetMovingFrameDisp(
m_disp);
168 m_ProjMatZ = bn::ublas::zero_matrix<NekDouble>(3, 3);
186 for (
int i = 0; i < 3; ++i)
205 std::string expression)
214 catch (
const std::runtime_error &)
224 const TiXmlElement *funcNameElmt;
226 funcNameElmt = pForce->FirstChildElement(
"FRAMEVELOCITY");
229 std::string FuncName = funcNameElmt->GetText();
231 "Function '" + FuncName +
"' is not defined in the session.");
235 std::string var =
m_session->GetVariable(i);
236 if (
m_session->DefinesFunction(FuncName, var))
244 std::vector<std::string> angularVar = {
"Omega_x",
"Omega_y",
"Omega_z"};
245 for (
int i = 0; i < 3; ++i)
247 std::string var = angularVar[i];
248 if (
m_session->DefinesFunction(FuncName, var))
257 for (
int i = 0; i < 2; ++i)
268 funcNameElmt = pForce->FirstChildElement(
"EXTERNALFORCE");
271 std::string FuncName = funcNameElmt->GetText();
272 if (
m_session->DefinesFunction(FuncName))
274 std::vector<std::string> forceVar = {
"fx",
"fy",
"fz"};
277 std::string var = forceVar[i];
278 if (
m_session->DefinesFunction(FuncName, var))
284 std::vector<std::string> momentVar = {
"Mx",
"My",
"Mz"};
285 for (
int i = 0; i < 3; ++i)
287 std::string var = momentVar[i];
288 if (
m_session->DefinesFunction(FuncName, var))
298 const TiXmlElement *pivotElmt = pForce->FirstChildElement(
"PIVOTPOINT");
301 std::vector<std::string> values;
302 std::string mssgStr = pivotElmt->GetText();
312 const TiXmlElement *TWSElmt =
313 pForce->FirstChildElement(
"TRAVELINGWAVESPEED");
316 std::vector<std::string> values;
317 std::string mssgStr = TWSElmt->GetText();
331 std::map<int, LibUtilities::EquationSharedPtr> initDisplacement;
332 funcNameElmt = pForce->FirstChildElement(
"INITIALDISPLACEMENT");
335 std::string FuncName = funcNameElmt->GetText();
337 "Function '" + FuncName +
"' is not defined in the session.");
339 std::vector<std::string> dispVar = {
"x",
"y",
"z"};
342 std::string var = dispVar[i];
343 if (
m_session->DefinesFunction(FuncName, var))
346 ->Evaluate(0., 0., 0., time);
350 std::vector<std::string> angleVar = {
"theta_x",
"theta_y",
"theta_z"};
351 for (
int i = 0; i < 3; ++i)
353 std::string var = angleVar[i];
354 if (
m_session->DefinesFunction(FuncName, var))
357 m_session->GetFunction(FuncName, var)->Evaluate();
364 const int dim,
const int rank,
367 int NumDof = dim + 1;
368 std::set<int> DirBCs;
369 const TiXmlElement *mssgTag;
373 for (
int i = 0; i < NumDof; ++i)
377 mssgTag = pForce->FirstChildElement(
"MOTIONPRESCRIBED");
380 std::vector<std::string> values;
381 mssgStr = mssgTag->GetText();
384 "MOTIONPRESCRIBED vector should be of size " +
385 std::to_string(NumDof));
386 for (
int i = 0; i < NumDof; ++i)
406 mssgTag = pForce->FirstChildElement(
"MASS");
410 std::vector<std::string> values;
411 mssgStr = mssgTag->GetText();
413 ASSERTL0(values.size() >= NumDof * NumDof,
414 "Mass matrix should be of size " + std::to_string(NumDof) +
415 "X" + std::to_string(NumDof));
416 for (
int i = 0; i < NumDof; ++i)
418 for (
int j = 0; j < NumDof; ++j)
427 mssgTag = pForce->FirstChildElement(
"DAMPING");
430 std::vector<std::string> values;
431 mssgStr = mssgTag->GetText();
433 ASSERTL0(values.size() >= NumDof * NumDof,
434 "Damping matrix should be of size " + std::to_string(NumDof) +
435 "X" + std::to_string(NumDof));
436 for (
int i = 0; i < NumDof; ++i)
438 for (
int j = 0; j < NumDof; ++j)
447 mssgTag = pForce->FirstChildElement(
"RIGIDITY");
450 std::vector<std::string> values;
451 mssgStr = mssgTag->GetText();
453 ASSERTL0(values.size() >= NumDof * NumDof,
454 "Rigidity matrix should be of size " + std::to_string(NumDof) +
455 "X" + std::to_string(NumDof));
456 for (
int i = 0; i < NumDof; ++i)
458 for (
int j = 0; j < NumDof; ++j)
468 if (
m_session->DefinesParameter(
"NewmarkBeta"))
472 if (
m_session->DefinesParameter(
"NewmarkGamma"))
474 gamma =
m_session->GetParameter(
"NewmarkGamma");
479 mssgTag = pForce->FirstChildElement(
"OutputFile");
482 filename = mssgTag->GetText();
488 if (!(filename.length() >= 4 &&
489 filename.substr(filename.length() - 4) ==
".mrf"))
499 <<
"Variables = t, x, ux, ax, y, uy, ay, theta, omega, domega"
504 m_outputStream <<
"Variables = t, x, ux, ax, y, uy, ay, z, uz, az, "
505 "theta, omega, domega"
512 mssgTag = pForce->FirstChildElement(
"OutputFrequency");
515 std::vector<std::string> values;
516 mssgStr = mssgTag->GetText();
520 "OutputFrequency should be greater than zero.");
523 for (
size_t i = 0; i < 3; ++i)
536 m_bodyVel[1][it.first] = it.second->Evaluate(0., 0., 0., time);
538 else if (it.first == 5)
540 m_bodyVel[1][NumDof - 1] = it.second->Evaluate(0., 0., 0., time);
559 m_velXYZ[it.first] = it.second->Evaluate(0., 0., 0., time);
563 m_omegaXYZ[it.first - 3] = it.second->Evaluate(0., 0., 0., time);
568 m_extForceXYZ[it.first] = it.second->Evaluate(0., 0., 0., time);
570 auto equ =
m_equ.lock();
571 ASSERTL0(equ,
"Weak pointer to the equation system is expired");
572 auto FluidEq = std::dynamic_pointer_cast<FluidInterface>(equ);
576 std::map<int, NekDouble> Dirs;
592 FluidEq->GetAeroForce(forcebody);
598 for (
size_t i = 0; i <
m_bodyVel[0].size(); ++i)
668 for (
size_t i = 0; i < bodyVel.size(); ++i)
674 for (
int iter = 0; iter < 2; ++iter)
677 for (
size_t i = 0; i < bodyVel.size(); ++i)
679 Vmath::Vcopy(bodyVel[i].size(), bodyVel[i], 1, tmpbodyVel[i],
695 for (
size_t i = 0; i < bodyVel.size(); ++i)
697 Vmath::Vcopy(bodyVel[i].size(), tmpbodyVel[i], 1, bodyVel[i], 1);
732 int npoints = fields[0]->GetNpoints();
733 addRotation(npoints, outarray, -1., inarray, outarray);
749 m_velXYZ[it.first] = it.second->Evaluate(0., 0., 0., time);
753 m_omegaXYZ[it.first - 3] = it.second->Evaluate(0., 0., 0., time);
758 m_extForceXYZ[it.first] = it.second->Evaluate(0., 0., 0., time);
760 std::map<int, NekDouble> Dirs;
773 auto equ =
m_equ.lock();
774 ASSERTL0(equ,
"Weak pointer to the equation system is expired");
775 auto FluidEq = std::dynamic_pointer_cast<FluidInterface>(equ);
777 FluidEq->GetAeroForce(forcebody);
780 for (
size_t i = 0; i <
m_bodyVel.size(); ++i)
785 for (
size_t i = 0; i <
m_bodyVel.size(); ++i)
792 for (
size_t i = 0; i < bodyVel[0].size(); ++i)
794 if (Dirs.find(i) == Dirs.end())
796 for (
size_t k = 0; k < bodyVel.size(); ++k)
813 vel[i] = bodyVel[1][i];
814 vel[i + 6] = bodyVel[2][i];
828 FluidEq->SetMovingFrameVelocities(vel);
830 FluidEq->SetMovingFrameProjectionMat(
m_ProjMatZ);
834 FluidEq->SetMovingFrameDisp(
m_disp);
846 ASSERTL0(&inarray1 != &outarray,
"inarray1 and outarray "
847 "should not be the same.");
854 if (&inarray0 != &outarray)
856 ASSERTL0(inarray0.size() == outarray.size(),
857 "inarray0 and outarray must have same dimentions");
858 for (
int i = 0; i < inarray0.size(); ++i)
869 Vmath::Svtvp(nPnts, cm, inarray1[1], 1, outarray[0], 1, outarray[0], 1);
870 Vmath::Svtvp(nPnts, cp, inarray1[0], 1, outarray[1], 1, outarray[1], 1);
878 Vmath::Svtvp(nPnts, cp, inarray1[1], 1, outarray[2], 1, outarray[2], 1);
879 Vmath::Svtvp(nPnts, cm, inarray1[2], 1, outarray[1], 1, outarray[1], 1);
887 Vmath::Svtvp(nPnts, cp, inarray1[2], 1, outarray[0], 1, outarray[0], 1);
888 Vmath::Svtvp(nPnts, cm, inarray1[0], 1, outarray[2], 1, outarray[2], 1);
898 int npoints = fields[0]->GetNpoints();
899 if (
m_isH2d && fields[0]->GetWaveSpace())
908 fields[0]->HomogeneousFwdTrans(npoints, tmpphys, tmpcoef);
909 Vmath::Vadd(npoints, tmpcoef, 1, inarray[i], 1, outarray[i], 1);
911 else if (&inarray != &outarray)
919 int npoints0 = npoints;
920 if (
m_isH1d && fields[0]->GetWaveSpace())
922 npoints0 =
m_hasPlane0 ? fields[0]->GetPlane(0)->GetNpoints() : 0;
929 inarray[i], 1, outarray[i], 1);
930 if (&inarray != &outarray && npoints != npoints0)
933 Vmath::Vcopy(npoints - npoints0, inarray[i] + npoints0, 1,
937 else if (&inarray != &outarray)
952 std::map<std::string, std::string> fieldMetaDataMap;
953 if (
m_session->DefinesFunction(
"InitialConditions"))
955 for (
int i = 0; i < pFields.size(); ++i)
959 vType =
m_session->GetFunctionType(
"InitialConditions",
964 std::string filename =
m_session->GetFunctionFilename(
965 "InitialConditions",
m_session->GetVariable(i));
967 fs::path pfilename(filename);
970 if (fs::is_directory(pfilename))
972 fs::path metafile(
"Info.xml");
973 fs::path fullpath = pfilename / metafile;
978 fld->ImportFieldMetaData(filename, fieldMetaDataMap);
983 std::string strt =
"Time";
984 if (fieldMetaDataMap.find(strt) != fieldMetaDataMap.end())
986 time = boost::lexical_cast<NekDouble>(
987 fieldMetaDataMap[strt]);
999 const TiXmlElement *pForce)
1001 std::map<std::string, std::string> vParams;
1002 vParams[
"OutputFile"] =
".dummyfilename";
1003 const TiXmlElement *param = pForce->FirstChildElement(
"BOUNDARY");
1004 ASSERTL0(param,
"Body surface should be assigned");
1005 vParams[
"Boundary"] = param->GetText();
1007 pSession,
m_equ.lock(), vParams);
1014 std::set<int> DirDoFs)
1028 for (
int i = 0; i <
m_rows; ++i)
1030 for (
int j = 0; j <
m_rows; ++j)
1036 if (DirDoFs.find(i) != DirDoFs.end() ||
1037 DirDoFs.find(j) != DirDoFs.end())
1039 value = (i == j) ? 1. : 0.;
1052 std::map<int, NekDouble> motionPrescribed)
1055 for (
int i = 0; i <
m_rows; ++i)
1062 for (
int i = 0; i <
m_rows; ++i)
1069 for (
int i = 0; i <
m_rows; ++i)
1071 rhs[i] = force[i] - fbk[i] + fbm[i];
1074 for (
int i = 0; i < 3; ++i)
1076 if (motionPrescribed.find(i) != motionPrescribed.end())
1078 rhs[i] = motionPrescribed[i];
1082 for (
auto it : motionPrescribed)
1091 for (
int i = 0; i <
m_rows; ++i)
1094 u[0][i] =
m_coeffs[2] * u[1][i] + bk[i];
1095 u[2][i] =
m_coeffs[0] * u[1][i] - bm[i];
1116 if (body != inertial && dim > 2)
1118 for (
int i = 2; i < dim; ++i)
1120 inertial[i] = body[i];
1127 NekDouble x = inertial[0], y = inertial[1];
1130 if (body != inertial && dim > 2)
1132 for (
int i = 2; i < dim; ++i)
1134 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.
int m_NumVariable
Number of variables.
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
Array< OneD, NekDouble > m_velxyz
static std::string classNameBody
Name of the class.
FilterAeroForcesSharedPtr m_aeroforceFilter
Array< OneD, bool > m_hasVel
Array< OneD, NekDouble > m_velXYZ
bn::ublas::matrix< NekDouble > m_ProjMatZ
~ForcingMovingReferenceFrame(void) override
Array< OneD, Array< OneD, NekDouble > > m_bodyVel
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)
void CheckForRestartTime(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, NekDouble &time)
void InitBodySolver(const TiXmlElement *pForce, const int dim, const int rank, const NekDouble time)
Array< OneD, NekDouble > m_disp
unsigned int m_outputFrequency
ForcingMovingReferenceFrame(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation)
Array< OneD, NekDouble > m_omegaXYZ
Array< OneD, NekDouble > m_omegaxyz
std::map< int, LibUtilities::EquationSharedPtr > m_extForceFunction
std::ofstream m_outputStream
void LoadParameters(const TiXmlElement *pForce, const NekDouble time)
Array< OneD, NekDouble > m_pivotPoint
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 Update(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.
Array< OneD, NekDouble > m_extForceXYZ
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.
std::set< int > m_DirDoFs
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, DNekMatSharedPtr M, DNekMatSharedPtr C, DNekMatSharedPtr K, std::set< int > DirDoFs)
DNekMatSharedPtr m_inverseMatrix
void Solve(Array< OneD, Array< OneD, NekDouble > > u, Array< OneD, NekDouble > force, std::map< int, NekDouble > motionPrescribed)
DNekMatSharedPtr m_coeffMatrix
Array< OneD, NekDouble > m_coeffs
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.
void Multiply(NekMatrix< ResultDataType, StandardMatrixTag > &result, const NekMatrix< LhsDataType, LhsMatrixType > &lhs, const ResultDataType &rhs)
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.
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)