46using namespace MultiRegions;
95 SetUpExpansions<ContField>(nvel);
99 SetUpExpansions<DisContField>(nvel);
106 SetUpExpansions<ContField>(nvel);
110 SetUpExpansions<DisContField>(nvel);
119 SetUpExpansions<ContField>(nvel);
123 SetUpExpansions<ContField3DHomogeneous1D>(nvel);
128 SetUpExpansions<ContField3DHomogeneous2D>(nvel);
135 SetUpExpansions<DisContField>(nvel);
139 SetUpExpansions<DisContField3DHomogeneous1D>(nvel);
144 SetUpExpansions<DisContField3DHomogeneous2D>(nvel);
152 "ShapeFunction must be defined in the session file.")
156 size_t physTot =
m_phi->GetTotPoints();
168 for (
size_t i = 0; i < nvel; ++i)
193 boost::iequals(timeInt.
method,
"IMEX") && 1 <= timeInt.
order &&
195 "The TimeIntegrationMethod scheme must be IMEX with order '1' to '4'.")
197 switch (timeInt.
order)
218 for (
size_t i = 0; i <
m_session->GetFilters().size(); ++i)
220 if (boost::iequals(
m_session->GetFilters()[i].first,
"AeroForcesSPM"))
237 "Smoothed Profile Method (SPM)");
272 static_pointer_cast<FilterAeroForcesSPM>(
305 size_t physTot =
m_fs[0]->GetNpoints();
315 for (
size_t i = 1; i < nvel; ++i)
361 size_t physTot =
m_phi->GetNpoints();
376 for (
size_t i = 0; i < nvel; ++i)
380 if (
m_session->DefinesSolverInfo(
"ForceBoundary") &&
381 boost::iequals(
m_session->GetSolverInfo(
"ForceBoundary"),
"True"))
415 for (
size_t b = 0; b < BndExp.size(); ++b)
418 if (BndCond[b]->GetBoundaryConditionType() ==
424 for (
size_t i = 0; i < nvel; ++i)
426 f_s[i] =
m_fs[0]->GetBndCondExpansions()[b]->GetPhys();
430 BndExp[b]->NormVectorIProductWRTBase(f_s, BndExp[b]->UpdatePhys());
488 size_t nq =
m_phi->GetNpoints();
490 for (
size_t i = 0; i < nvel; ++i)
499 m_fields[ind]->HomogeneousBwdTrans(nq, fields[i], tmpField);
500 m_fs[i]->HomogeneousBwdTrans(nq,
m_fs[i]->GetPhys(),
501 m_fs[i]->UpdatePhys());
505 tmpField = fields[i];
510 m_fs[i]->UpdatePhys(), 1);
512 m_fs[i]->UpdatePhys(), 1);
518 m_fs[i]->HomogeneousFwdTrans(nq,
m_fs[i]->GetPhys(),
519 m_fs[i]->UpdatePhys());
539 TiXmlElement *functionDef = function->FirstChildElement();
540 ASSERTL0(functionDef,
"At least one element must be defined in " + funcName)
543 string varName = functionDef->Attribute(
"VAR");
544 while (functionDef && !boost::iequals(varName, elemName))
546 functionDef = functionDef->NextSiblingElement();
547 varName = functionDef->Attribute(
"VAR");
551 "Variable " + elemName +
" must be defined in " + funcName +
".");
555 int err = functionDef->QueryStringAttribute(
"USERDEFINEDTYPE", &attr);
556 bool output = boost::iequals(attr,
"TimeDependent");
558 ASSERTL0((err == TIXML_NO_ATTRIBUTE) || (err == TIXML_SUCCESS && output),
559 "USERDEFINEDTYPE in " + elemName +
560 " must be TimeDependent if defined");
575 TiXmlElement *conds =
m_session->GetElement(
"Nektar/Conditions");
576 TiXmlElement *function = conds->FirstChildElement(
"FUNCTION");
579 string functionType = function->Attribute(
"NAME");
580 while (function && !boost::iequals(functionType, functionName))
582 function = function->NextSiblingElement(
"FUNCTION");
583 functionType = function->Attribute(
"NAME");
595 TiXmlElement *child = function->FirstChildElement();
599 if (boost::iequals(child->ValueStr(),
"F"))
603 int status = child->QueryStringAttribute(
"FILE", &fileName);
605 "An FLD file with the values "
606 "of the phi function has to be supplied.")
607 ASSERTL0(boost::iequals(fileName.substr(fileName.length() - 4),
".fld"),
608 "A valid FLD file must be supplied in the "
609 "'ShapeFunction' field.")
613 std::vector<LibUtilities::FieldDefinitionsSharedPtr> fieldDef;
614 std::vector<std::vector<NekDouble>> fieldData;
618 phiFile->Import(fileName, fieldDef, fieldData, fieldMetaData);
621 ASSERTL0(fieldData.size() == 1,
"Only one field (phi) must be "
622 "defined in the FLD file.")
626 m_phi->ExtractDataToCoeffs(fieldDef[0], fieldData[0], tmp,
627 m_phi->UpdateCoeffs());
#define ASSERTL0(condition, msg)
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field.
NekDouble m_kinvis
Kinematic viscosity.
ExtrapolateSharedPtr m_extrapolation
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
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::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
void UpdateForcing(const Array< OneD, const Array< OneD, NekDouble > > &fields, NekDouble dt)
For a body with a velocity , the force applied to the fluid ensures that the IBC are met:
bool m_filePhi
Flag indicating that phi was defined in a file.
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Generates the summary of the current simulation.
void SetCorrectionPressureBCs()
Updates the BCs for boundaries with Neumann or periodic BCs in the pressure:
void v_SolveUnsteadyStokesSystem(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, NekDouble a_iixDt) override
Linear terms due to pressure and visosity are calculated here. After solving the velocity filed witho...
Array< OneD, Array< OneD, NekDouble > > m_up
Velocity of the immersed body(ies)
std::vector< std::string > m_velName
Vector storing the names of the components of \u_p.
NekDouble m_gamma0
Stiffly-stable scheme coefficient.
int m_forcesFilter
Position of "AeroForcesSPM" filter in 'm_session->GetFilters()'.
void v_InitObject(bool DeclareField=true) override
Initialisation object for EquationSystem.
bool m_timeDependentPhi
Flag that is true when phi depends on time.
~SmoothedProfileMethod() override
Destroy the Smoothed Profile Method object.
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
SmoothedProfileMethod(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Construct a new Smoothed Profile Method object.
void UpdatePhiUp(NekDouble time)
Calculates the values of the shape function.
SolverUtils::SessionFunctionSharedPtr m_phiEvaluator
Function that evaluates the values of \Phi.
void SolveCorrectedVelocity(Array< OneD, Array< OneD, NekDouble > > &Forcing, Array< OneD, Array< OneD, NekDouble > > &fields, NekDouble dt)
Corrects the velocity field so that the IBs are taken into account. Solves the explicit equation:
TiXmlElement * GetFunctionHdl(std::string functionName)
Returns a handle to the requested function. Returns NULL if it does not exist.
static std::string className
Name of class.
bool m_timeDependentUp
Flag signaling if depends on time.
MultiRegions::ExpListSharedPtr m_pressureP
Correction pressure field for SPM.
bool GetVarTimeDependence(std::string funcName, std::string attrName)
True if the function is timedependent, false otherwise.
Array< OneD, Array< OneD, NekDouble > > m_upPrev
Array< OneD, MultiRegions::ExpListSharedPtr > m_fs
Forcing function 'f_s'.
MultiRegions::ExpListSharedPtr m_phi
Shape function 'phi' as expansion list.
static std::string solverTypeLookupId
void SolveCorrectionPressure(const Array< OneD, NekDouble > &Forcing)
Solves the Poisson equation for the correction pressure :
void SetUpCorrectionPressure(const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing)
Sets the forcing term of the equation for the correction pressure :
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
enum HomogeneousType m_HomogeneousType
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
Defines a forcing term to be explicitly applied.
Base class for unsteady solvers.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
Array< OneD, Array< OneD, NekDouble > > m_F
virtual void v_SolveUnsteadyStokesSystem(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time, const NekDouble a_iixDt)
void v_InitObject(bool DeclareField=true) override
Initialisation object for EquationSystem.
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
std::shared_ptr< FieldIO > FieldIOSharedPtr
std::map< std::string, std::string > FieldMetaDataMap
std::shared_ptr< SessionReader > SessionReaderSharedPtr
MultiRegions::Direction const DirCartesianMap[]
std::vector< std::pair< std::string, std::string > > SummaryList
EquationSystemFactory & GetEquationSystemFactory()
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
std::map< ConstFactorType, NekDouble > ConstFactorMap
StdRegions::ConstFactorMap factors
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = 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 Vvtvm(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvm (vector times vector minus vector): z = w*x - y
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
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.