Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Protected Member Functions | Protected Attributes | Static Private Attributes | List of all members
Nektar::Extrapolate Class Reference

#include <Extrapolate.h>

Inheritance diagram for Nektar::Extrapolate:
Inheritance graph
[legend]
Collaboration diagram for Nektar::Extrapolate:
Collaboration graph
[legend]

Public Member Functions

 Extrapolate (const LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields, MultiRegions::ExpListSharedPtr pPressure, const Array< OneD, int > pVel, const AdvectionTermSharedPtr advObject)
virtual ~Extrapolate ()
void GenerateHOPBCMap ()
void SubSteppingTimeIntegration (const int intMethod, const LibUtilities::TimeIntegrationWrapperSharedPtr &IntegrationScheme)
void SubStepSaveFields (const int nstep)
void SubStepSetPressureBCs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, const NekDouble Aii_DT, NekDouble kinvis)
void SubStepAdvance (const LibUtilities::TimeIntegrationSolutionSharedPtr &integrationSoln, const int nstep, NekDouble time)
void MountHOPBCs (int HBCdata, NekDouble kinvis, Array< OneD, NekDouble > &Q, Array< OneD, const NekDouble > &Advection)
void EvaluatePressureBCs (const Array< OneD, const Array< OneD, NekDouble > > &fields, const Array< OneD, const Array< OneD, NekDouble > > &N, NekDouble kinvis)
Array< OneD, NekDoubleGetMaxStdVelocity (const Array< OneD, Array< OneD, NekDouble > > inarray)

Protected Member Functions

virtual void v_SubSteppingTimeIntegration (int intMethod, const LibUtilities::TimeIntegrationWrapperSharedPtr &IntegrationScheme)=0
virtual void v_SubStepSaveFields (int nstep)=0
virtual void v_SubStepSetPressureBCs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, NekDouble Aii_DT, NekDouble kinvis)=0
virtual void v_SubStepAdvance (const LibUtilities::TimeIntegrationSolutionSharedPtr &integrationSoln, int nstep, NekDouble time)=0
virtual void v_MountHOPBCs (int HBCdata, NekDouble kinvis, Array< OneD, NekDouble > &Q, Array< OneD, const NekDouble > &Advection)=0
void CalcNeumannPressureBCs (const Array< OneD, const Array< OneD, NekDouble > > &fields, const Array< OneD, const Array< OneD, NekDouble > > &N, NekDouble kinvis)
void CalcOutflowBCs (const Array< OneD, const Array< OneD, NekDouble > > &fields, const Array< OneD, const Array< OneD, NekDouble > > &N, NekDouble kinvis)
void RollOver (Array< OneD, Array< OneD, NekDouble > > &input)
void CurlCurl (Array< OneD, Array< OneD, const NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q, const int j)

Protected Attributes

LibUtilities::SessionReaderSharedPtr m_session
LibUtilities::CommSharedPtr m_comm
Array< OneD,
MultiRegions::ExpListSharedPtr
m_fields
MultiRegions::ExpListSharedPtr m_pressure
 Pointer to field holding pressure field.
Array< OneD, int > m_velocity
 int which identifies which components of m_fields contains the velocity (u,v,w);
AdvectionTermSharedPtr m_advObject
Array< OneD, Array< OneD,
NekDouble > > 
m_previousVelFields
int m_curl_dim
 Curl-curl dimensionality.
int m_bnd_dim
 bounday dimensionality
Array< OneD, const
SpatialDomains::BoundaryConditionShPtr
m_PBndConds
 pressure boundary conditions container
Array< OneD,
MultiRegions::ExpListSharedPtr
m_PBndExp
 pressure boundary conditions expansion container
int m_pressureCalls
 number of times the high-order pressure BCs have been called
int m_pressureBCsMaxPts
 Maximum points used in pressure BC evaluation.
int m_pressureBCsElmtMaxPts
 Maximum points used in Element adjacent to pressure BC evaluation.
int m_intSteps
 Maximum points used in pressure BC evaluation.
NekDouble m_timestep
bool m_SingleMode
 Flag to determine if single homogeneous mode is used.
bool m_HalfMode
 Flag to determine if half homogeneous mode is used.
bool m_MultipleModes
 Flag to determine if use multiple homogenenous modes are used.
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous)
int m_npointsX
 number of points in X direction (if homogeneous)
int m_npointsY
 number of points in Y direction (if homogeneous)
int m_npointsZ
 number of points in Z direction (if homogeneous)
Array< OneD, int > m_pressureBCtoElmtID
 Id of element to which pressure boundary condition belongs.
Array< OneD, int > m_pressureBCtoTraceID
 Id of edge (2D) or face (3D) to which pressure boundary condition belongs.
Array< OneD, Array< OneD,
NekDouble > > 
m_pressureHBCs
 Storage for current and previous levels of high order pressure boundary conditions.
Array< OneD, Array< OneD,
NekDouble > > 
m_acceleration
 Storage for current and previous levels of the acceleration term.
Array< OneD, HBCInfom_HBCdata
 data structure to old all the information regarding High order pressure BCs
Array< OneD, NekDoublem_wavenumber
 wave number 2 pi k /Lz
Array< OneD, NekDoublem_negWavenumberSq
 minus Square of wavenumber
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_outflowVel
 Storage for current and previous velocity fields at the otuflow for high order outflow BCs.

Static Private Attributes

static std::string def
static NekDouble StifflyStable_Betaq_Coeffs [3][3]
static NekDouble StifflyStable_Alpha_Coeffs [3][3]
static NekDouble StifflyStable_Gamma0_Coeffs [3]

Detailed Description

Definition at line 83 of file Extrapolate.h.

Constructor & Destructor Documentation

Nektar::Extrapolate::Extrapolate ( const LibUtilities::SessionReaderSharedPtr  pSession,
Array< OneD, MultiRegions::ExpListSharedPtr pFields,
MultiRegions::ExpListSharedPtr  pPressure,
const Array< OneD, int >  pVel,
const AdvectionTermSharedPtr  advObject 
)

Definition at line 56 of file Extrapolate.cpp.

References m_comm, m_session, and m_timestep.

: m_session(pSession),
m_fields(pFields),
m_pressure(pPressure),
m_velocity(pVel),
m_advObject(advObject)
{
m_session->LoadParameter("TimeStep", m_timestep, 0.01);
m_comm = m_session->GetComm();
}
Nektar::Extrapolate::~Extrapolate ( )
virtual

Definition at line 72 of file Extrapolate.cpp.

{
}

Member Function Documentation

void Nektar::Extrapolate::CalcNeumannPressureBCs ( const Array< OneD, const Array< OneD, NekDouble > > &  fields,
const Array< OneD, const Array< OneD, NekDouble > > &  N,
NekDouble  kinvis 
)
protected

Unified routine for calculation high-oder terms

Casting the boundary expansion to the specific case

Picking up the element where the HOPBc is located

Assigning

Calculating the curl-curl and storing it in Q

Definition at line 189 of file Extrapolate.cpp.

References ASSERTL0, CurlCurl(), Nektar::MultiRegions::e2D, Nektar::MultiRegions::e3D, Nektar::MultiRegions::e3DH1D, Nektar::MultiRegions::e3DH2D, m_acceleration, m_bnd_dim, m_curl_dim, m_HBCdata, m_PBndExp, m_pressure, m_pressureBCsElmtMaxPts, m_pressureBCsMaxPts, and MountHOPBCs().

Referenced by EvaluatePressureBCs().

{
Array<OneD, NekDouble> Pvals;
Array<OneD, NekDouble> Uvals;
Array<OneD, Array<OneD, const NekDouble> > Velocity(m_curl_dim);
Array<OneD, Array<OneD, const NekDouble> > Advection(m_bnd_dim);
Array<OneD, Array<OneD, NekDouble> > BndValues(m_bnd_dim);
Array<OneD, Array<OneD, NekDouble> > Q(m_bnd_dim);
for(int i = 0; i < m_bnd_dim; i++)
{
BndValues[i] = Array<OneD, NekDouble> (m_pressureBCsMaxPts,0.0);
Q[i] = Array<OneD, NekDouble> (m_pressureBCsElmtMaxPts,0.0);
}
for(int j = 0 ; j < m_HBCdata.num_elements() ; j++)
{
/// Casting the boundary expansion to the specific case
Pbc = boost::dynamic_pointer_cast<StdRegions::StdExpansion>
(m_PBndExp[m_HBCdata[j].m_bndryElmtID]
->GetExp(m_HBCdata[j].m_bndElmtOffset));
/// Picking up the element where the HOPBc is located
elmt = m_pressure->GetExp(m_HBCdata[j].m_globalElmtID);
/// Assigning
for(int i = 0; i < m_bnd_dim; i++)
{
Velocity[i] = fields[i] + m_HBCdata[j].m_physOffset;
Advection[i] = N[i] + m_HBCdata[j].m_physOffset;
}
// for the 3DH1D case we need to grab the conjugate mode
if(m_pressure->GetExpType() == MultiRegions::e3DH1D)
{
Velocity[2] = fields[2] + m_HBCdata[j].m_assPhysOffset;
}
/// Calculating the curl-curl and storing it in Q
CurlCurl(Velocity,Q,j);
// Mounting advection component into the high-order condition
for(int i = 0; i < m_bnd_dim; i++)
{
MountHOPBCs(m_HBCdata[j].m_ptsInElmt,kinvis,Q[i],Advection[i]);
}
Pvals = m_PBndExp[m_HBCdata[j].m_bndryElmtID]->UpdateCoeffs()
+ m_PBndExp[m_HBCdata[j].m_bndryElmtID]
->GetCoeff_Offset(m_HBCdata[j].m_bndElmtOffset);
Uvals = (m_acceleration[0]) + m_HBCdata[j].m_coeffOffset;
// Getting values on the edge and filling the pressure boundary
// expansion and the acceleration term. Multiplication by the
// normal is required
switch(m_pressure->GetExpType())
{
{
elmt->GetEdgePhysVals(m_HBCdata[j].m_elmtTraceID, Pbc,
Q[0], BndValues[0]);
elmt->GetEdgePhysVals(m_HBCdata[j].m_elmtTraceID, Pbc,
Q[1], BndValues[1]);
Pbc->NormVectorIProductWRTBase(BndValues[0], BndValues[1],
Pvals);
elmt->GetEdgePhysVals(m_HBCdata[j].m_elmtTraceID, Pbc,
Velocity[0], BndValues[0]);
elmt->GetEdgePhysVals(m_HBCdata[j].m_elmtTraceID, Pbc,
Velocity[1], BndValues[1]);
Pbc->NormVectorIProductWRTBase(BndValues[0], BndValues[1],
Uvals);
}
break;
{
if(m_HBCdata[j].m_elmtTraceID == 0)
{
(m_PBndExp[m_HBCdata[j].m_bndryElmtID]->UpdateCoeffs()
+ m_PBndExp[m_HBCdata[j].m_bndryElmtID]
->GetCoeff_Offset(
m_HBCdata[j].m_bndElmtOffset))[0]
= -1.0*Q[0][0];
}
else if (m_HBCdata[j].m_elmtTraceID == 1)
{
(m_PBndExp[m_HBCdata[j].m_bndryElmtID]->UpdateCoeffs()
+ m_PBndExp[m_HBCdata[j].m_bndryElmtID]
->GetCoeff_Offset(
m_HBCdata[j].m_bndElmtOffset))[0]
= Q[0][m_HBCdata[j].m_ptsInElmt-1];
}
else
{
ASSERTL0(false,
"In the 3D homogeneous 2D approach BCs edge "
"ID can be just 0 or 1 ");
}
}
break;
{
elmt->GetFacePhysVals(m_HBCdata[j].m_elmtTraceID, Pbc,
Q[0], BndValues[0]);
elmt->GetFacePhysVals(m_HBCdata[j].m_elmtTraceID, Pbc,
Q[1], BndValues[1]);
elmt->GetFacePhysVals(m_HBCdata[j].m_elmtTraceID, Pbc,
Q[2], BndValues[2]);
Pbc->NormVectorIProductWRTBase(BndValues[0], BndValues[1],
BndValues[2], Pvals);
elmt->GetFacePhysVals(m_HBCdata[j].m_elmtTraceID, Pbc,
Velocity[0], BndValues[0]);
elmt->GetFacePhysVals(m_HBCdata[j].m_elmtTraceID, Pbc,
Velocity[1], BndValues[1]);
elmt->GetFacePhysVals(m_HBCdata[j].m_elmtTraceID, Pbc,
Velocity[2], BndValues[2]);
Pbc->NormVectorIProductWRTBase(BndValues[0], BndValues[1],
BndValues[2], Uvals);
}
break;
default:
ASSERTL0(0,"Dimension not supported");
break;
}
}
}
void Nektar::Extrapolate::CalcOutflowBCs ( const Array< OneD, const Array< OneD, NekDouble > > &  fields,
const Array< OneD, const Array< OneD, NekDouble > > &  N,
NekDouble  kinvis 
)
protected

Definition at line 327 of file Extrapolate.cpp.

References Nektar::MultiRegions::DirCartesianMap, Nektar::SpatialDomains::eHighOutflow, m_bnd_dim, m_curl_dim, m_fields, m_intSteps, m_outflowVel, m_PBndConds, m_PBndExp, m_pressureBCsElmtMaxPts, m_pressureBCsMaxPts, m_pressureBCtoElmtID, m_pressureBCtoTraceID, m_pressureCalls, m_session, m_velocity, RollOver(), Vmath::Smul(), Vmath::Svtvp(), Vmath::Vadd(), Vmath::Vcopy(), Vmath::Vvtvp(), and Vmath::Zero().

Referenced by EvaluatePressureBCs().

{
static bool init = true;
static bool noHOBC = false;
if(noHOBC == true)
{
return;
}
if(init) // set up storage for boundary velocity at outflow
{
init = false;
int totbndpts = 0;
for(int n = 0; n < m_PBndConds.num_elements(); ++n)
{
if(m_PBndConds[n]->GetUserDefined()
{
totbndpts += m_PBndExp[n]->GetTotPoints();
}
}
if(totbndpts == 0)
{
noHOBC = true;
return;
}
m_outflowVel = Array<OneD, Array<OneD, Array<OneD, NekDouble> > > (m_bnd_dim);
for(int i = 0; i < m_bnd_dim; ++i)
{
m_outflowVel[i] = Array<OneD, Array<OneD, NekDouble> >(m_curl_dim);
for(int j = 0; j < m_curl_dim; ++j)
{
// currently just set up for 2nd order extrapolation
m_outflowVel[i][j] = Array<OneD, NekDouble>(totbndpts,0.0);
}
}
}
Array<OneD, Array<OneD, const SpatialDomains::BoundaryConditionShPtr> >
UBndConds(m_curl_dim);
Array<OneD, Array<OneD, MultiRegions::ExpListSharedPtr> >
UBndExp(m_curl_dim);
for (int i = 0; i < m_curl_dim; ++i)
{
UBndConds[i] = m_fields[m_velocity[i]]->GetBndConditions();
UBndExp[i] = m_fields[m_velocity[i]]->GetBndCondExpansions();
}
Array<OneD, Array<OneD, NekDouble> > BndValues(m_bnd_dim);
Array<OneD, Array<OneD, NekDouble> > BndElmt (m_bnd_dim);
Array<OneD, Array<OneD, NekDouble> > nGradu(m_bnd_dim);
Array<OneD, NekDouble > gradtmp (m_pressureBCsElmtMaxPts),
nGradu[0] = Array<OneD, NekDouble>(m_bnd_dim*m_pressureBCsMaxPts);
for(int i = 0; i < m_bnd_dim; ++i)
{
BndElmt[i] = Array<OneD, NekDouble> (m_pressureBCsElmtMaxPts,
0.0);
BndValues[i] = Array<OneD, NekDouble> (m_pressureBCsMaxPts,0.0);
nGradu[i] = nGradu[0] + i*m_pressureBCsMaxPts;
}
int nbc,cnt,cnt_start;
int veloffset = 0;
int nint = min(m_pressureCalls,m_intSteps);
Array<OneD, NekDouble> PBCvals, UBCvals;
Array<OneD, Array<OneD, NekDouble> > ubc(m_curl_dim);
Array<OneD, Array<OneD, NekDouble> > normals;
cnt = 0;
for(int n = 0; n < m_PBndConds.num_elements(); ++n)
{
// Do outflow boundary conditions if they exist
if(m_PBndConds[n]->GetUserDefined() == SpatialDomains::eHighOutflow)
{
for(int i = 0; i < m_PBndExp[n]->GetExpSize(); ++i,cnt++)
{
cnt = max(cnt,m_PBndExp[n]->GetTotPoints());
}
}
}
for(int i =0; i < m_curl_dim; ++i)
{
ubc[i] = Array<OneD, NekDouble>(cnt);
}
NekDouble U0,delta;
m_session->LoadParameter("U0_HighOrderBC",U0,1.0);
m_session->LoadParameter("Delta_HighOrderBC",delta,1/20.0);
cnt = 0;
for(int n = 0; n < m_PBndConds.num_elements(); ++n)
{
cnt_start = cnt;
// Do outflow boundary conditions if they exist
if(m_PBndConds[n]->GetUserDefined() == SpatialDomains::eHighOutflow)
{
for(int i = 0; i < m_PBndExp[n]->GetExpSize(); ++i,cnt++)
{
// find element and edge of this expansion.
Bc = boost::dynamic_pointer_cast<StdRegions::StdExpansion>
(m_PBndExp[n]->GetExp(i));
int elmtid = m_pressureBCtoElmtID[cnt];
elmt = m_fields[0]->GetExp(elmtid);
int offset = m_fields[0]->GetPhys_Offset(elmtid);
// Determine extrapolated U,V values
int nq = elmt->GetTotPoints();
// currently just using first order approximation here.
// previously have obtained value from m_integrationSoln
for(int j = 0; j < m_bnd_dim; ++j)
{
Vmath::Vcopy(nq, &fields[m_velocity[j]][offset], 1,
&BndElmt[j][0], 1);
}
int nbc = m_PBndExp[n]->GetExp(i)->GetTotPoints();
int boundary = m_pressureBCtoTraceID[cnt];
Array<OneD, NekDouble> veltmp, ptmp(nbc,0.0),
normDotu(nbc,0.0), utot(nbc,0.0),
divU(nbc,0.0);
normals=elmt->GetSurfaceNormal(boundary);
Vmath::Zero(m_bnd_dim*m_pressureBCsMaxPts,nGradu[0],1);
for (int j = 0; j < m_bnd_dim; j++)
{
// Calculate Grad u = du/dx, du/dy, du/dz, etc.
for (int k = 0; k< m_bnd_dim; k++)
{
elmt->PhysDeriv(MultiRegions::DirCartesianMap[k],
BndElmt[j], gradtmp);
elmt->GetTracePhysVals(boundary, Bc, gradtmp,
fgradtmp);
Vmath::Vvtvp(nbc,normals[k], 1, fgradtmp, 1,
nGradu[j], 1, nGradu[j],1);
if(j == k)
{
Vmath::Vadd(nbc,fgradtmp, 1, divU, 1, divU, 1);
}
}
}
// extract velocity and store
for(int j = 0; j < m_bnd_dim; ++j)
{
elmt->GetTracePhysVals(boundary,Bc,BndElmt[j],
veltmp = m_outflowVel[j][0] + veloffset);
}
// extrapolate velocity
if(nint <= 1)
{
for(int j = 0; j < m_bnd_dim; ++j)
{
veltmp = m_outflowVel[j][0] +veloffset, 1,
BndValues[j], 1);
}
}
else // only set up for 2nd order extrapolation
{
for(int j = 0; j < m_bnd_dim; ++j)
{
Vmath::Smul(nbc, 2.0,
veltmp = m_outflowVel[j][0] + veloffset, 1,
BndValues[j], 1);
Vmath::Svtvp(nbc, -1.0,
veltmp = m_outflowVel[j][1] + veloffset, 1,
BndValues[j], 1,
BndValues[j], 1);
}
}
// Set up |u|^2, n.u, div(u), and (n.grad(u) . n) for
// pressure condition
for(int j = 0; j < m_bnd_dim; ++j)
{
Vmath::Vvtvp(nbc, BndValues[j], 1, BndValues[j], 1,
utot, 1, utot, 1);
Vmath::Vvtvp(nbc, normals[j], 1, BndValues[j], 1,
normDotu, 1, normDotu, 1);
Vmath::Vvtvp(nbc, normals[j], 1, nGradu[j], 1,
ptmp, 1, ptmp, 1);
}
PBCvals = m_PBndExp[n]->GetPhys() +
m_PBndExp[n]->GetPhys_Offset(i);
for(int k = 0; k < nbc; ++k)
{
NekDouble fac = 0.5*(1.0-tanh(normDotu[k]/(U0*delta)));
// Set up Dirichlet pressure condition and
// store in ptmp (PBCvals contains a
// function from the input file )
ptmp[k] = kinvis * ptmp[k] - 0.5 * utot[k] * fac
- PBCvals[k];
}
int u_offset = UBndExp[0][n]->GetPhys_Offset(i);
for(int j = 0; j < m_bnd_dim; ++j)
{
UBCvals = UBndExp[j][n]->GetPhys()
+ UBndExp[j][n]->GetPhys_Offset(i);
for(int k = 0; k < nbc; ++k)
{
NekDouble fac = 0.5 * (1.0 - tanh(normDotu[k]
/ (U0 * delta)));
ubc[j][k + u_offset] = (1.0 / kinvis)
* (UBCvals[k] + 0.5 * utot[k] * fac
* normals[j][k]);
}
}
// set up pressure boundary condition
PBCvals = m_PBndExp[n]->UpdateCoeffs()
+ m_PBndExp[n]->GetCoeff_Offset(i);
m_PBndExp[n]->GetExp(i)->FwdTrans(ptmp,PBCvals);
veloffset += nbc;
}
// Now set up Velocity conditions.
for(int j = 0; j < m_bnd_dim; j++)
{
if(UBndConds[j][n]->GetUserDefined()
{
cnt = cnt_start;
for(int i = 0; i < UBndExp[0][n]->GetExpSize();
++i, cnt++)
{
(m_PBndExp[n]->GetExp(i));
(UBndExp[0][n]->GetExp(i));
nbc = UBndExp[0][n]->GetExp(i)->GetTotPoints();
int boundary = m_pressureBCtoTraceID[cnt];
Array<OneD, NekDouble> pb(nbc), ub(nbc);
int elmtid = m_pressureBCtoElmtID[cnt];
elmt = m_fields[0]->GetExp(elmtid);
normals = elmt->GetSurfaceNormal(boundary);
// Get p from projected boundary condition
PBCvals = m_PBndExp[n]->UpdateCoeffs()
+ m_PBndExp[n]->GetCoeff_Offset(i);
Pbc->BwdTrans(PBCvals,pb);
int u_offset = UBndExp[j][n]->GetPhys_Offset(i);
for(int k = 0; k < nbc; ++k)
{
ub[k] = ubc[j][k + u_offset]
+ pb[k] * normals[j][k] / kinvis;
}
UBCvals = UBndExp[j][n]->UpdateCoeffs()
+ UBndExp[j][n]->GetCoeff_Offset(i);
Bc->IProductWRTBase(ub,UBCvals);
}
}
}
}
else
{
cnt += m_PBndExp[n]->GetExpSize();
}
}
}
void Nektar::Extrapolate::CurlCurl ( Array< OneD, Array< OneD, const NekDouble > > &  Vel,
Array< OneD, Array< OneD, NekDouble > > &  Q,
const int  j 
)
protected

Curl Curl routine - dimension dependent

Definition at line 627 of file Extrapolate.cpp.

References ASSERTL0, Nektar::MultiRegions::DirCartesianMap, Nektar::MultiRegions::e2D, Nektar::MultiRegions::e3D, Nektar::MultiRegions::e3DH1D, Nektar::MultiRegions::e3DH2D, m_fields, m_HBCdata, m_negWavenumberSq, m_pressureBCsElmtMaxPts, m_wavenumber, Vmath::Smul(), Vmath::Vadd(), and Vmath::Vsub().

Referenced by CalcNeumannPressureBCs().

{
= m_fields[0]->GetExp(m_HBCdata[j].m_globalElmtID);
Array<OneD,NekDouble> Vx(m_pressureBCsElmtMaxPts);
Array<OneD,NekDouble> Uy(m_pressureBCsElmtMaxPts);
switch(m_fields[0]->GetExpType())
{
{
Array<OneD,NekDouble> Dummy(m_pressureBCsElmtMaxPts);
elmt->PhysDeriv(MultiRegions::DirCartesianMap[0], Vel[1], Vx);
elmt->PhysDeriv(MultiRegions::DirCartesianMap[1], Vel[0], Uy);
Vmath::Vsub(m_HBCdata[j].m_ptsInElmt, Vx, 1, Uy, 1, Dummy, 1);
elmt->PhysDeriv(Dummy,Q[1],Q[0]);
Vmath::Smul(m_HBCdata[j].m_ptsInElmt, -1.0, Q[1], 1, Q[1], 1);
}
break;
{
Array<OneD,NekDouble> Wz(m_pressureBCsElmtMaxPts);
Array<OneD,NekDouble> Dummy1(m_pressureBCsElmtMaxPts);
Array<OneD,NekDouble> Dummy2(m_pressureBCsElmtMaxPts);
elmt->PhysDeriv(MultiRegions::DirCartesianMap[0], Vel[1], Vx);
elmt->PhysDeriv(MultiRegions::DirCartesianMap[1], Vel[0], Uy);
Vmath::Smul(m_HBCdata[j].m_ptsInElmt, m_wavenumber[j],
Vel[2], 1, Wz, 1);
elmt->PhysDeriv(MultiRegions::DirCartesianMap[1], Vx, Dummy1);
elmt->PhysDeriv(MultiRegions::DirCartesianMap[1], Uy, Dummy2);
Vmath::Vsub(m_HBCdata[j].m_ptsInElmt, Dummy1, 1, Dummy2, 1,
Q[0], 1);
Vel[0], 1, Dummy1, 1);
elmt->PhysDeriv(MultiRegions::DirCartesianMap[0], Wz, Dummy2);
Vmath::Vsub(m_HBCdata[j].m_ptsInElmt, Q[0], 1, Dummy1, 1,
Q[0], 1);
Vmath::Vadd(m_HBCdata[j].m_ptsInElmt, Q[0], 1, Dummy2, 1,
Q[0], 1);
elmt->PhysDeriv(MultiRegions::DirCartesianMap[1], Wz, Dummy1);
Vel[1], 1, Dummy2, 1);
Vmath::Vsub(m_HBCdata[j].m_ptsInElmt, Dummy1, 1, Dummy2, 1,
Q[1], 1);
elmt->PhysDeriv(MultiRegions::DirCartesianMap[0], Vx, Dummy1);
elmt->PhysDeriv(MultiRegions::DirCartesianMap[0], Uy, Dummy2);
Vmath::Vsub(m_HBCdata[j].m_ptsInElmt, Q[1], 1, Dummy1, 1,
Q[1], 1);
Vmath::Vadd(m_HBCdata[j].m_ptsInElmt, Q[1], 1, Dummy2, 1,
Q[1], 1);
}
break;
{
Array<OneD,NekDouble> Wx(m_pressureBCsElmtMaxPts);
Array<OneD,NekDouble> Wz(m_pressureBCsElmtMaxPts);
Array<OneD,NekDouble> Uz(m_pressureBCsElmtMaxPts);
Array<OneD,NekDouble> qz(m_pressureBCsElmtMaxPts);
Array<OneD,NekDouble> qy(m_pressureBCsElmtMaxPts);
elmt->PhysDeriv(MultiRegions::DirCartesianMap[0], Vel[2], Wx);
elmt->PhysDeriv(MultiRegions::DirCartesianMap[0], Vel[1], Vx);
Vel[0], 1, Uy, 1);
Vmath::Smul(m_HBCdata[j].m_ptsInElmt, m_wavenumber[j],
Vel[0], 1, Uz, 1);
Vmath::Vsub(m_HBCdata[j].m_ptsInElmt, Wz, 1, Wx, 1,
qy, 1);
Vmath::Vsub(m_HBCdata[j].m_ptsInElmt, Vx, 1, Uy, 1,
qz, 1);
qz, 1, Uy, 1);
Vmath::Smul(m_HBCdata[j].m_ptsInElmt, m_wavenumber[j],
qy, 1, Uz, 1);
Vmath::Vsub(m_HBCdata[j].m_ptsInElmt, Uy, 1, Uz, 1,
Q[0], 1);
}
break;
{
Array<OneD,NekDouble> Dummy(m_pressureBCsElmtMaxPts);
Array<OneD,NekDouble> Vz(m_pressureBCsElmtMaxPts);
Array<OneD,NekDouble> Uz(m_pressureBCsElmtMaxPts);
Array<OneD,NekDouble> Wx(m_pressureBCsElmtMaxPts);
Array<OneD,NekDouble> Wy(m_pressureBCsElmtMaxPts);
elmt->PhysDeriv(Vel[0], Dummy, Uy, Uz);
elmt->PhysDeriv(Vel[1], Vx, Dummy, Vz);
elmt->PhysDeriv(Vel[2], Wx, Wy, Dummy);
Vmath::Vsub(m_HBCdata[j].m_ptsInElmt, Wy, 1, Vz, 1, Q[0], 1);
Vmath::Vsub(m_HBCdata[j].m_ptsInElmt, Uz, 1, Wx, 1, Q[1], 1);
Vmath::Vsub(m_HBCdata[j].m_ptsInElmt, Vx, 1, Uy, 1, Q[2], 1);
elmt->PhysDeriv(Q[0], Dummy, Wy, Vx);
elmt->PhysDeriv(Q[1], Wx, Dummy, Uz);
elmt->PhysDeriv(Q[2], Vz, Uy, Dummy);
Vmath::Vsub(m_HBCdata[j].m_ptsInElmt, Uy, 1, Uz, 1, Q[0], 1);
Vmath::Vsub(m_HBCdata[j].m_ptsInElmt, Vx, 1, Vz, 1, Q[1], 1);
Vmath::Vsub(m_HBCdata[j].m_ptsInElmt, Wx, 1, Wy, 1, Q[2], 1);
}
break;
default:
ASSERTL0(0,"Dimension not supported");
break;
}
}
void Nektar::Extrapolate::EvaluatePressureBCs ( const Array< OneD, const Array< OneD, NekDouble > > &  fields,
const Array< OneD, const Array< OneD, NekDouble > > &  N,
NekDouble  kinvis 
)

Function to extrapolate the new pressure boundary condition. Based on the velocity field and on the advection term. Acceleration term is also computed. This routine is a general one for 2d and 3D application and it can be called directly from velocity correction scheme. Specialisation on dimensionality is redirected to the CalcPressureBCs method.

Definition at line 89 of file Extrapolate.cpp.

References CalcNeumannPressureBCs(), CalcOutflowBCs(), Nektar::SpatialDomains::eHigh, m_acceleration, m_HBCdata, m_intSteps, m_PBndConds, m_PBndExp, m_pressureCalls, m_pressureHBCs, m_timestep, RollOver(), Vmath::Smul(), StifflyStable_Alpha_Coeffs, StifflyStable_Betaq_Coeffs, StifflyStable_Gamma0_Coeffs, Vmath::Svtvp(), and Vmath::Vcopy().

Referenced by Nektar::SubSteppingExtrapolate::v_SubStepSetPressureBCs().

{
if(m_HBCdata.num_elements()>0)
{
Array<OneD, NekDouble> tmp;
Array<OneD, NekDouble> accelerationTerm;
int n,cnt;
int nint = min(m_pressureCalls,m_intSteps);
int nlevels = m_pressureHBCs.num_elements();
int acc_order = 0;
accelerationTerm =
Array<OneD, NekDouble>(m_acceleration[0].num_elements(), 0.0);
// Rotate HOPBCs storage
// Rotate acceleration term
// Calculate Neumann BCs at current level
CalcNeumannPressureBCs(fields, N, kinvis);
// Copy High order values into storage array
for(cnt = n = 0; n < m_PBndConds.num_elements(); ++n)
{
// High order boundary condition;
if(m_PBndConds[n]->GetUserDefined() == SpatialDomains::eHigh)
{
int nq = m_PBndExp[n]->GetNcoeffs();
Vmath::Vcopy(nq, &(m_PBndExp[n]->GetCoeffs()[0]), 1,
&(m_pressureHBCs[0])[cnt], 1);
cnt += nq;
}
}
//Calculate acceleration term at level n based on previous steps
if (m_pressureCalls > 2)
{
acc_order = min(m_pressureCalls-2,m_intSteps);
accelerationTerm, 1);
for(int i = 0; i < acc_order; i++)
{
-1*StifflyStable_Alpha_Coeffs[acc_order-1][i],
m_acceleration[i+1], 1,
accelerationTerm, 1,
accelerationTerm, 1);
}
}
// Adding acceleration term to HOPBCs
accelerationTerm, 1,
// Extrapolate to n+1
m_pressureHBCs[nint-1], 1,
m_pressureHBCs[nlevels-1], 1);
for(n = 0; n < nint-1; ++n)
{
m_pressureHBCs[n],1,m_pressureHBCs[nlevels-1],1,
m_pressureHBCs[nlevels-1],1);
}
// Copy values of [dP/dn]^{n+1} in the pressure bcs storage.
// m_pressureHBCS[nlevels-1] will be cancelled at next time step
for(cnt = n = 0; n < m_PBndConds.num_elements(); ++n)
{
if(m_PBndConds[n]->GetUserDefined() == SpatialDomains::eHigh)
{
int nq = m_PBndExp[n]->GetNcoeffs();
Vmath::Vcopy(nq, &(m_pressureHBCs[nlevels-1])[cnt], 1,
&(m_PBndExp[n]->UpdateCoeffs()[0]), 1);
cnt += nq;
}
}
}
CalcOutflowBCs(fields, N, kinvis);
}
void Nektar::Extrapolate::GenerateHOPBCMap ( )

Map to directly locate HOPBCs position and offsets in all scenarios

Definition at line 784 of file Extrapolate.cpp.

References ASSERTL0, Nektar::MultiRegions::e2D, Nektar::MultiRegions::e3D, Nektar::MultiRegions::e3DH1D, Nektar::MultiRegions::e3DH2D, Nektar::SpatialDomains::eHigh, m_acceleration, m_bnd_dim, m_comm, m_curl_dim, m_HalfMode, m_HBCdata, m_intSteps, m_LhomZ, m_MultipleModes, m_negWavenumberSq, m_npointsY, m_npointsZ, m_PBndConds, m_PBndExp, m_pressure, m_pressureBCsElmtMaxPts, m_pressureBCsMaxPts, m_pressureBCtoElmtID, m_pressureBCtoTraceID, m_pressureCalls, m_pressureHBCs, m_session, m_SingleMode, m_wavenumber, Nektar::LibUtilities::ReduceSum, and sign.

{
m_PBndConds = m_pressure->GetBndConditions();
m_PBndExp = m_pressure->GetBndCondExpansions();
// Set up mapping from pressure boundary condition to pressure element
// details.
m_pressure->GetBoundaryToElmtMap(m_pressureBCtoElmtID,
// find the maximum values of points for pressure BC evaluation
int cnt, n;
for(cnt = n = 0; n < m_PBndConds.num_elements(); ++n)
{
for(int i = 0; i < m_PBndExp[n]->GetExpSize(); ++i)
{
m_PBndExp[n]->GetExp(i)->GetTotPoints());
->GetTotPoints());
}
}
// Storage array for high order pressure BCs
m_pressureHBCs = Array<OneD, Array<OneD, NekDouble> > (m_intSteps);
m_acceleration = Array<OneD, Array<OneD, NekDouble> > (m_intSteps + 1);
int HBCnumber = 0;
for(cnt = n = 0; n < m_PBndConds.num_elements(); ++n)
{
// High order boundary condition;
if(m_PBndConds[n]->GetUserDefined() == SpatialDomains::eHigh)
{
cnt += m_PBndExp[n]->GetNcoeffs();
HBCnumber += m_PBndExp[n]->GetExpSize();
}
}
int checkHBC = HBCnumber;
m_comm->AllReduce(checkHBC,LibUtilities::ReduceSum);
//ASSERTL0(checkHBC > 0 ,"At least one high-order pressure boundary "
// "condition is required for scheme "
// "consistency");
m_acceleration[0] = Array<OneD, NekDouble>(cnt, 0.0);
for(n = 0; n < m_intSteps; ++n)
{
m_pressureHBCs[n] = Array<OneD, NekDouble>(cnt, 0.0);
m_acceleration[n+1] = Array<OneD, NekDouble>(cnt, 0.0);
}
switch(m_pressure->GetExpType())
{
{
m_bnd_dim = 2;
}
break;
{
m_bnd_dim = 2;
}
break;
{
m_bnd_dim = 1;
}
break;
{
m_bnd_dim = 3;
}
break;
default:
ASSERTL0(0,"Dimension not supported");
break;
}
m_HBCdata = Array<OneD, HBCInfo>(HBCnumber);
switch(m_pressure->GetExpType())
{
{
int coeff_count = 0;
int exp_size;
int j=0;
int cnt = 0;
for(int n = 0 ; n < m_PBndConds.num_elements(); ++n)
{
exp_size = m_PBndExp[n]->GetExpSize();
if(m_PBndConds[n]->GetUserDefined()
{
for(int i = 0; i < exp_size; ++i,cnt++)
{
m_HBCdata[j].m_globalElmtID =
elmt = m_pressure->GetExp(
m_HBCdata[j].m_globalElmtID);
m_HBCdata[j].m_ptsInElmt =
elmt->GetTotPoints();
m_HBCdata[j].m_physOffset =
m_pressure->GetPhys_Offset(
m_HBCdata[j].m_globalElmtID);
m_HBCdata[j].m_bndElmtOffset = i;
m_HBCdata[j].m_elmtTraceID =
m_HBCdata[j].m_bndryElmtID = n;
m_HBCdata[j].m_coeffOffset = coeff_count;
coeff_count += elmt->GetEdgeNcoeffs(
m_HBCdata[j].m_elmtTraceID);
j = j+1;
}
}
else // setting if just standard BC no High order
{
cnt += exp_size;
}
}
}
break;
{
Array<OneD, unsigned int> planes;
planes = m_pressure->GetZIDs();
int num_planes = planes.num_elements();
int num_elm_per_plane = (m_pressure->GetExpSize())/num_planes;
m_wavenumber = Array<OneD, NekDouble>(HBCnumber);
m_negWavenumberSq = Array<OneD, NekDouble>(HBCnumber);
int coeff_count = 0;
int exp_size, exp_size_per_plane;
int j=0;
int K;
NekDouble sign = -1.0;
int cnt = 0;
m_session->MatchSolverInfo("ModeType", "SingleMode",
m_SingleMode, false);
m_session->MatchSolverInfo("ModeType", "HalfMode",
m_HalfMode, false);
m_session->MatchSolverInfo("ModeType", "MultipleModes",
m_MultipleModes, false);
m_session->LoadParameter("LZ", m_LhomZ);
// Stability Analysis flags
if(m_session->DefinesSolverInfo("ModeType"))
{
{
}
else if(m_HalfMode)
{
}
else if(m_MultipleModes)
{
m_npointsZ = m_session->GetParameter("HomModesZ");
}
else
{
ASSERTL0(false, "SolverInfo ModeType not valid");
}
}
else
{
m_npointsZ = m_session->GetParameter("HomModesZ");
}
for(int k = 0; k < num_planes; k++)
{
K = planes[k]/2;
for(int n = 0 ; n < m_PBndConds.num_elements(); ++n)
{
exp_size = m_PBndExp[n]->GetExpSize();
exp_size_per_plane = exp_size/num_planes;
if(m_PBndConds[n]->GetUserDefined() == SpatialDomains::eHigh)
{
for(int i = 0; i < exp_size_per_plane; ++i,cnt++)
{
m_HBCdata[j].m_globalElmtID = m_pressureBCtoElmtID[cnt];
elmt = m_pressure->GetExp(m_HBCdata[j].m_globalElmtID);
m_HBCdata[j].m_ptsInElmt = elmt->GetTotPoints();
m_HBCdata[j].m_physOffset = m_pressure->GetPhys_Offset(m_HBCdata[j].m_globalElmtID);
m_HBCdata[j].m_bndElmtOffset = i+k*exp_size_per_plane;
m_HBCdata[j].m_elmtTraceID = m_pressureBCtoTraceID[cnt];
m_HBCdata[j].m_bndryElmtID = n;
m_HBCdata[j].m_coeffOffset = coeff_count;
coeff_count += elmt->GetEdgeNcoeffs(m_HBCdata[j].m_elmtTraceID);
{
m_wavenumber[j] = -2*M_PI/m_LhomZ;
m_negWavenumberSq[j] = -1.0*m_wavenumber[j]*m_wavenumber[j];
}
{
m_wavenumber[j] = 2*M_PI/m_LhomZ;
m_negWavenumberSq[j] = -1.0*m_wavenumber[j]*m_wavenumber[j];
}
else
{
m_wavenumber[j] = 2*M_PI*sign*(NekDouble(K))/m_LhomZ;
m_negWavenumberSq[j] = -1.0*m_wavenumber[j]*m_wavenumber[j];
}
int assElmtID;
if(k%2==0)
{
{
assElmtID = m_HBCdata[j].m_globalElmtID;
}
else
{
assElmtID = m_HBCdata[j].m_globalElmtID + num_elm_per_plane;
}
}
else
{
assElmtID = m_HBCdata[j].m_globalElmtID - num_elm_per_plane;
}
m_HBCdata[j].m_assPhysOffset = m_pressure->GetPhys_Offset(assElmtID);
j = j+1;
}
}
else // setting if just standard BC no High order
{
cnt += exp_size_per_plane;
}
}
sign = -1.0*sign;
}
}
break;
{
int cnt = 0;
int exp_size, exp_size_per_line;
int j=0;
for(int k1 = 0; k1 < m_npointsZ; k1++)
{
for(int k2 = 0; k2 < m_npointsY; k2++)
{
for(int n = 0 ; n < m_PBndConds.num_elements(); ++n)
{
exp_size = m_PBndExp[n]->GetExpSize();
exp_size_per_line = exp_size/(m_npointsZ*m_npointsY);
if(m_PBndConds[n]->GetUserDefined() == SpatialDomains::eHigh)
{
for(int i = 0; i < exp_size_per_line; ++i,cnt++)
{
// find element and edge of this expansion.
m_HBCdata[j].m_globalElmtID = m_pressureBCtoElmtID[cnt];
elmt = m_pressure->GetExp(m_HBCdata[j].m_globalElmtID);
m_HBCdata[j].m_ptsInElmt = elmt->GetTotPoints();
m_HBCdata[j].m_physOffset = m_pressure->GetPhys_Offset(m_HBCdata[j].m_globalElmtID);
m_HBCdata[j].m_bndElmtOffset = i+(k1*m_npointsY+k2)*exp_size_per_line;
m_HBCdata[j].m_elmtTraceID = m_pressureBCtoTraceID[cnt];
m_HBCdata[j].m_bndryElmtID = n;
}
}
else
{
cnt += exp_size_per_line;
}
}
}
}
}
break;
default:
ASSERTL0(0,"Dimension not supported");
break;
}
}
Array< OneD, NekDouble > Nektar::Extrapolate::GetMaxStdVelocity ( const Array< OneD, Array< OneD, NekDouble > >  inarray)

Definition at line 1090 of file Extrapolate.cpp.

References ASSERTL0, Nektar::SpatialDomains::eDeformed, m_curl_dim, and m_fields.

Referenced by Nektar::SubSteppingExtrapolate::GetSubstepTimeStep().

{
// Checking if the problem is 2D
ASSERTL0(m_curl_dim >= 2, "Method not implemented for 1D");
int n_points_0 = m_fields[0]->GetExp(0)->GetTotPoints();
int n_element = m_fields[0]->GetExpSize();
int nvel = inarray.num_elements();
int cnt;
NekDouble pntVelocity;
// Getting the standard velocity vector on the 2D normal space
Array<OneD, Array<OneD, NekDouble> > stdVelocity(nvel);
Array<OneD, NekDouble> maxV(n_element, 0.0);
for (int i = 0; i < nvel; ++i)
{
stdVelocity[i] = Array<OneD, NekDouble>(n_points_0);
}
if (nvel == 2)
{
cnt = 0.0;
for (int el = 0; el < n_element; ++el)
{
int n_points = m_fields[0]->GetExp(el)->GetTotPoints();
ptsKeys = m_fields[0]->GetExp(el)->GetPointsKeys();
// reset local space if necessary
if(n_points != n_points_0)
{
for (int j = 0; j < nvel; ++j)
{
stdVelocity[j] = Array<OneD, NekDouble>(n_points);
}
n_points_0 = n_points;
}
Array<TwoD, const NekDouble> gmat =
m_fields[0]->GetExp(el)->GetGeom()->GetMetricInfo()->GetDerivFactors(ptsKeys);
if (m_fields[0]->GetExp(el)->GetGeom()->GetMetricInfo()->GetGtype()
{
for (int i = 0; i < n_points; i++)
{
stdVelocity[0][i] = gmat[0][i]*inarray[0][i+cnt]
+ gmat[2][i]*inarray[1][i+cnt];
stdVelocity[1][i] = gmat[1][i]*inarray[0][i+cnt]
+ gmat[3][i]*inarray[1][i+cnt];
}
}
else
{
for (int i = 0; i < n_points; i++)
{
stdVelocity[0][i] = gmat[0][0]*inarray[0][i+cnt]
+ gmat[2][0]*inarray[1][i+cnt];
stdVelocity[1][i] = gmat[1][0]*inarray[0][i+cnt]
+ gmat[3][0]*inarray[1][i+cnt];
}
}
cnt += n_points;
for (int i = 0; i < n_points; i++)
{
pntVelocity = stdVelocity[0][i]*stdVelocity[0][i]
+ stdVelocity[1][i]*stdVelocity[1][i];
if (pntVelocity>maxV[el])
{
maxV[el] = pntVelocity;
}
}
maxV[el] = sqrt(maxV[el]);
}
}
else
{
cnt = 0;
for (int el = 0; el < n_element; ++el)
{
int n_points = m_fields[0]->GetExp(el)->GetTotPoints();
ptsKeys = m_fields[0]->GetExp(el)->GetPointsKeys();
// reset local space if necessary
if(n_points != n_points_0)
{
for (int j = 0; j < nvel; ++j)
{
stdVelocity[j] = Array<OneD, NekDouble>(n_points);
}
n_points_0 = n_points;
}
Array<TwoD, const NekDouble> gmat =
m_fields[0]->GetExp(el)->GetGeom()->GetMetricInfo()->GetDerivFactors(ptsKeys);
if (m_fields[0]->GetExp(el)->GetGeom()->GetMetricInfo()->GetGtype()
{
for (int i = 0; i < n_points; i++)
{
stdVelocity[0][i] = gmat[0][i]*inarray[0][i+cnt]
+ gmat[3][i]*inarray[1][i+cnt]
+ gmat[6][i]*inarray[2][i+cnt];
stdVelocity[1][i] = gmat[1][i]*inarray[0][i+cnt]
+ gmat[4][i]*inarray[1][i+cnt]
+ gmat[7][i]*inarray[2][i+cnt];
stdVelocity[2][i] = gmat[2][i]*inarray[0][i+cnt]
+ gmat[5][i]*inarray[1][i+cnt]
+ gmat[8][i]*inarray[2][i+cnt];
}
}
else
{
for (int i = 0; i < n_points; i++)
{
stdVelocity[0][i] = gmat[0][0]*inarray[0][i+cnt]
+ gmat[3][0]*inarray[1][i+cnt]
+ gmat[6][0]*inarray[2][i+cnt];
stdVelocity[1][i] = gmat[1][0]*inarray[0][i+cnt]
+ gmat[4][0]*inarray[1][i+cnt]
+ gmat[7][0]*inarray[2][i+cnt];
stdVelocity[2][i] = gmat[2][0]*inarray[0][i+cnt]
+ gmat[5][0]*inarray[1][i+cnt]
+ gmat[8][0]*inarray[2][i+cnt];
}
}
cnt += n_points;
for (int i = 0; i < n_points; i++)
{
pntVelocity = stdVelocity[0][i]*stdVelocity[0][i]
+ stdVelocity[1][i]*stdVelocity[1][i]
+ stdVelocity[2][i]*stdVelocity[2][i];
if (pntVelocity > maxV[el])
{
maxV[el] = pntVelocity;
}
}
maxV[el] = sqrt(maxV[el]);
//cout << maxV[el]*maxV[el] << endl;
}
}
return maxV;
}
void Nektar::Extrapolate::MountHOPBCs ( int  HBCdata,
NekDouble  kinvis,
Array< OneD, NekDouble > &  Q,
Array< OneD, const NekDouble > &  Advection 
)
inline

Definition at line 306 of file Extrapolate.h.

References v_MountHOPBCs().

Referenced by CalcNeumannPressureBCs().

{
v_MountHOPBCs(HBCdata,kinvis,Q,Advection);
}
void Nektar::Extrapolate::RollOver ( Array< OneD, Array< OneD, NekDouble > > &  input)
protected

Function to roll time-level storages to the next step layout. The stored data associated with the oldest time-level (not required anymore) are moved to the top, where they will be overwritten as the solution process progresses.

Definition at line 764 of file Extrapolate.cpp.

Referenced by CalcOutflowBCs(), and EvaluatePressureBCs().

{
int nlevels = input.num_elements();
Array<OneD, NekDouble> tmp;
tmp = input[nlevels-1];
for(int n = nlevels-1; n > 0; --n)
{
input[n] = input[n-1];
}
input[0] = tmp;
}
void Nektar::Extrapolate::SubStepAdvance ( const LibUtilities::TimeIntegrationSolutionSharedPtr integrationSoln,
const int  nstep,
NekDouble  time 
)
inline

Definition at line 295 of file Extrapolate.h.

References v_SubStepAdvance().

{
v_SubStepAdvance(integrationSoln,nstep, time);
}
void Nektar::Extrapolate::SubSteppingTimeIntegration ( const int  intMethod,
const LibUtilities::TimeIntegrationWrapperSharedPtr IntegrationScheme 
)
inline

Definition at line 265 of file Extrapolate.h.

References v_SubSteppingTimeIntegration().

{
v_SubSteppingTimeIntegration(intMethod, IntegrationScheme);
}
void Nektar::Extrapolate::SubStepSaveFields ( const int  nstep)
inline

Definition at line 275 of file Extrapolate.h.

References v_SubStepSaveFields().

{
}
void Nektar::Extrapolate::SubStepSetPressureBCs ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
const NekDouble  Aii_DT,
NekDouble  kinvis 
)
inline

Definition at line 284 of file Extrapolate.h.

References v_SubStepSetPressureBCs().

{
v_SubStepSetPressureBCs(inarray,Aii_DT,kinvis);
}
virtual void Nektar::Extrapolate::v_MountHOPBCs ( int  HBCdata,
NekDouble  kinvis,
Array< OneD, NekDouble > &  Q,
Array< OneD, const NekDouble > &  Advection 
)
protectedpure virtual
virtual void Nektar::Extrapolate::v_SubStepAdvance ( const LibUtilities::TimeIntegrationSolutionSharedPtr integrationSoln,
int  nstep,
NekDouble  time 
)
protectedpure virtual
virtual void Nektar::Extrapolate::v_SubSteppingTimeIntegration ( int  intMethod,
const LibUtilities::TimeIntegrationWrapperSharedPtr IntegrationScheme 
)
protectedpure virtual
virtual void Nektar::Extrapolate::v_SubStepSaveFields ( int  nstep)
protectedpure virtual
virtual void Nektar::Extrapolate::v_SubStepSetPressureBCs ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
NekDouble  Aii_DT,
NekDouble  kinvis 
)
protectedpure virtual

Member Data Documentation

std::string Nektar::Extrapolate::def
staticprivate
Initial value:
"StandardExtrapolate", "StandardExtrapolate")

Definition at line 254 of file Extrapolate.h.

Array<OneD, Array<OneD, NekDouble> > Nektar::Extrapolate::m_acceleration
protected

Storage for current and previous levels of the acceleration term.

Definition at line 238 of file Extrapolate.h.

Referenced by CalcNeumannPressureBCs(), EvaluatePressureBCs(), and GenerateHOPBCMap().

AdvectionTermSharedPtr Nektar::Extrapolate::m_advObject
protected

Definition at line 183 of file Extrapolate.h.

Referenced by Nektar::SubSteppingExtrapolate::SubStepAdvection().

int Nektar::Extrapolate::m_bnd_dim
protected

bounday dimensionality

Definition at line 191 of file Extrapolate.h.

Referenced by CalcNeumannPressureBCs(), CalcOutflowBCs(), and GenerateHOPBCMap().

LibUtilities::CommSharedPtr Nektar::Extrapolate::m_comm
protected
int Nektar::Extrapolate::m_curl_dim
protected
Array<OneD, MultiRegions::ExpListSharedPtr> Nektar::Extrapolate::m_fields
protected
bool Nektar::Extrapolate::m_HalfMode
protected

Flag to determine if half homogeneous mode is used.

Definition at line 216 of file Extrapolate.h.

Referenced by GenerateHOPBCMap().

Array<OneD, HBCInfo > Nektar::Extrapolate::m_HBCdata
protected

data structure to old all the information regarding High order pressure BCs

Definition at line 241 of file Extrapolate.h.

Referenced by CalcNeumannPressureBCs(), CurlCurl(), EvaluatePressureBCs(), and GenerateHOPBCMap().

int Nektar::Extrapolate::m_intSteps
protected
NekDouble Nektar::Extrapolate::m_LhomZ
protected

physical length in Z direction (if homogeneous)

Definition at line 220 of file Extrapolate.h.

Referenced by GenerateHOPBCMap().

bool Nektar::Extrapolate::m_MultipleModes
protected

Flag to determine if use multiple homogenenous modes are used.

Definition at line 218 of file Extrapolate.h.

Referenced by GenerateHOPBCMap().

Array<OneD, NekDouble> Nektar::Extrapolate::m_negWavenumberSq
protected

minus Square of wavenumber

Definition at line 247 of file Extrapolate.h.

Referenced by CurlCurl(), and GenerateHOPBCMap().

int Nektar::Extrapolate::m_npointsX
protected

number of points in X direction (if homogeneous)

Definition at line 222 of file Extrapolate.h.

int Nektar::Extrapolate::m_npointsY
protected

number of points in Y direction (if homogeneous)

Definition at line 223 of file Extrapolate.h.

Referenced by GenerateHOPBCMap().

int Nektar::Extrapolate::m_npointsZ
protected

number of points in Z direction (if homogeneous)

Definition at line 224 of file Extrapolate.h.

Referenced by GenerateHOPBCMap().

Array<OneD, Array<OneD, Array<OneD, NekDouble > > > Nektar::Extrapolate::m_outflowVel
protected

Storage for current and previous velocity fields at the otuflow for high order outflow BCs.

Definition at line 250 of file Extrapolate.h.

Referenced by CalcOutflowBCs().

Array<OneD, const SpatialDomains::BoundaryConditionShPtr> Nektar::Extrapolate::m_PBndConds
protected

pressure boundary conditions container

Definition at line 194 of file Extrapolate.h.

Referenced by CalcOutflowBCs(), EvaluatePressureBCs(), and GenerateHOPBCMap().

Array<OneD, MultiRegions::ExpListSharedPtr> Nektar::Extrapolate::m_PBndExp
protected

pressure boundary conditions expansion container

Definition at line 197 of file Extrapolate.h.

Referenced by CalcNeumannPressureBCs(), CalcOutflowBCs(), EvaluatePressureBCs(), and GenerateHOPBCMap().

MultiRegions::ExpListSharedPtr Nektar::Extrapolate::m_pressure
protected

Pointer to field holding pressure field.

Definition at line 177 of file Extrapolate.h.

Referenced by CalcNeumannPressureBCs(), and GenerateHOPBCMap().

int Nektar::Extrapolate::m_pressureBCsElmtMaxPts
protected

Maximum points used in Element adjacent to pressure BC evaluation.

Definition at line 206 of file Extrapolate.h.

Referenced by CalcNeumannPressureBCs(), CalcOutflowBCs(), CurlCurl(), and GenerateHOPBCMap().

int Nektar::Extrapolate::m_pressureBCsMaxPts
protected
Array<OneD, int> Nektar::Extrapolate::m_pressureBCtoElmtID
protected

Id of element to which pressure boundary condition belongs.

Definition at line 229 of file Extrapolate.h.

Referenced by Nektar::SubSteppingExtrapolate::AddDuDt2D(), Nektar::SubSteppingExtrapolate::AddDuDt3D(), CalcOutflowBCs(), and GenerateHOPBCMap().

Array<OneD, int> Nektar::Extrapolate::m_pressureBCtoTraceID
protected

Id of edge (2D) or face (3D) to which pressure boundary condition belongs.

Definition at line 232 of file Extrapolate.h.

Referenced by Nektar::SubSteppingExtrapolate::AddDuDt2D(), Nektar::SubSteppingExtrapolate::AddDuDt3D(), CalcOutflowBCs(), and GenerateHOPBCMap().

int Nektar::Extrapolate::m_pressureCalls
protected

number of times the high-order pressure BCs have been called

Definition at line 200 of file Extrapolate.h.

Referenced by CalcOutflowBCs(), EvaluatePressureBCs(), and GenerateHOPBCMap().

Array<OneD, Array<OneD, NekDouble> > Nektar::Extrapolate::m_pressureHBCs
protected

Storage for current and previous levels of high order pressure boundary conditions.

Definition at line 235 of file Extrapolate.h.

Referenced by EvaluatePressureBCs(), and GenerateHOPBCMap().

Array<OneD, Array<OneD, NekDouble> > Nektar::Extrapolate::m_previousVelFields
protected

Definition at line 185 of file Extrapolate.h.

LibUtilities::SessionReaderSharedPtr Nektar::Extrapolate::m_session
protected
bool Nektar::Extrapolate::m_SingleMode
protected

Flag to determine if single homogeneous mode is used.

Definition at line 214 of file Extrapolate.h.

Referenced by GenerateHOPBCMap().

NekDouble Nektar::Extrapolate::m_timestep
protected
Array<OneD, int> Nektar::Extrapolate::m_velocity
protected
Array<OneD, NekDouble> Nektar::Extrapolate::m_wavenumber
protected

wave number 2 pi k /Lz

Definition at line 244 of file Extrapolate.h.

Referenced by CurlCurl(), and GenerateHOPBCMap().

NekDouble Nektar::Extrapolate::StifflyStable_Alpha_Coeffs
staticprivate
Initial value:
{
{ 1.0, 0.0, 0.0},{ 2.0, -0.5, 0.0},{ 3.0, -1.5, 1.0/3.0}}

Definition at line 258 of file Extrapolate.h.

Referenced by EvaluatePressureBCs().

NekDouble Nektar::Extrapolate::StifflyStable_Betaq_Coeffs
staticprivate
Initial value:
{
{ 1.0, 0.0, 0.0},{ 2.0, -1.0, 0.0},{ 3.0, -3.0, 1.0}}

Definition at line 257 of file Extrapolate.h.

Referenced by EvaluatePressureBCs().

NekDouble Nektar::Extrapolate::StifflyStable_Gamma0_Coeffs
staticprivate
Initial value:
{
1.0, 1.5, 11.0/6.0}

Definition at line 259 of file Extrapolate.h.

Referenced by EvaluatePressureBCs().