51 return exp->GetExp(i);
56 std::ofstream out(filename.c_str());
57 exp->WriteVtkHeader(out);
58 for (
int i = 0; i < exp->GetExpSize(); ++i)
60 exp->WriteVtkPieceHeader(out, i);
61 exp->WriteVtkPieceFooter(out, i);
63 exp->WriteVtkFooter(out);
71 exp->FwdTrans(in, out);
80 exp->BwdTrans(in, out);
89 exp->IProductWRTBase(in, out);
98 exp->MultiplyByInvMassMatrix(in, out);
105 const py::object constFactorMap,
106 const py::object varCoeffMap)
113 if (!constFactorMap.is_none())
115 facMap = py::extract<StdRegions::ConstFactorMap>(constFactorMap);
117 if (!varCoeffMap.is_none())
119 coeffMap = py::extract<StdRegions::VarCoeffMap>(varCoeffMap);
122 exp->HelmSolve(in, out, facMap, coeffMap);
138 return exp->L2(in, err);
145 return exp->Linf(in);
153 return exp->Linf(in, err);
158 int nPhys = exp->GetNpoints();
159 int coordim = exp->GetCoordim(0);
161 std::vector<Array<OneD, NekDouble> > coords(coordim);
162 for (
int i = 0; i < coordim; ++i)
170 exp->GetCoords(coords[0]);
171 return py::make_tuple(coords[0]);
174 exp->GetCoords(coords[0], coords[1]);
175 return py::make_tuple(coords[0], coords[1]);
178 exp->GetCoords(coords[0], coords[1], coords[2]);
179 return py::make_tuple(coords[0], coords[1], coords[2]);
190 exp->SetPhysArray(inarray);
197 exp->SetPhys(inarray);
202 return exp->GetPhys();
207 return exp->Integral();
212 std::stringstream ss;
213 ss << static_cast<const void*>(&(exp->GetPhys()[0]));
219 exp->ClearGlobalLinSysManager();
230 int nExp = exp->GetExpSize();
234 for (
int i = 0; i < nExp; ++i)
236 elementGIDs[i] = exp->GetExp(i)->GetGeom()->GetGlobalID();
239 std::vector<LibUtilities::FieldDefinitionsSharedPtr> def;
240 std::vector<std::vector<NekDouble> > data;
248 Vmath::Zero(exp->GetNcoeffs(), exp->UpdateCoeffs(), 1);
251 for (
int i = 0; i < def.size(); ++i)
254 for (
int j = 0; j < def[i]->m_fields.size(); ++j)
256 if (def[i]->m_fields[j] == varName)
264 exp->ExtractDataToCoeffs(
265 def[i], data[i], def[i]->m_fields[idx], exp->UpdateCoeffs());
269 std::cout <<
"Field " + varName +
" not found." << std::endl;
273 exp->BwdTrans_IterPerExp(exp->GetCoeffs(), exp->UpdatePhys());
281 std::shared_ptr<ExpList>,
291 .def(
"GetNcoeffs", GetNcoeffs)
307 py::arg(
"constFactorMap") = py::object(),
308 py::arg(
"varCoeffMap") = py::object()
NekDouble ExpList_Linf_Error(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in, const Array< OneD, const NekDouble > &err)
void ExpList_ResetManagers(ExpListSharedPtr exp)
NekDouble ExpList_Linf(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
Array< OneD, NekDouble > ExpList_IProductWRTBase(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
NekDouble ExpList_Integral(ExpListSharedPtr exp)
ExpansionSharedPtr ExpList_GetExp(ExpListSharedPtr exp, int i)
void ExpList_SetPhysArray(ExpListSharedPtr exp, Array< OneD, NekDouble > inarray)
const Array< OneD, const NekDouble > ExpList_GetPhys(ExpListSharedPtr exp)
void ExpList_LoadField(ExpListSharedPtr exp, std::string filename, std::string varName)
Array< OneD, NekDouble > ExpList_MultiplyByInvMassMatrix(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
void ExpList_WriteVTK(ExpListSharedPtr exp, std::string filename)
Array< OneD, NekDouble > ExpList_FwdTrans(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
void ExpList_SetPhys(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &inarray)
py::tuple ExpList_GetCoords(ExpListSharedPtr exp)
NekDouble ExpList_L2(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
Array< OneD, NekDouble > ExpList_HelmSolve(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in, const py::object constFactorMap, const py::object varCoeffMap)
Array< OneD, NekDouble > ExpList_BwdTrans(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
std::string ExpList_GetPhysAddress(ExpListSharedPtr exp)
NekDouble ExpList_L2_Error(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in, const Array< OneD, const NekDouble > &err)
static std::shared_ptr< FieldIO > CreateForFile(const LibUtilities::SessionReaderSharedPtr session, const std::string &filename)
Construct a FieldIO object for a given input filename.
Base class for all multi-elemental spectral/hp expansions.
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
int GetCoeff_Offset(int n) const
Get the start offset position for a global list of m_coeffs correspoinding to element n.
int GetExpSize(void)
This function returns the number of elements in the expansion.
int GetNpoints(void) const
Returns the total number of quadrature points m_npoints .
NekDouble PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &phys)
bool GetPhysState(void) const
This function indicates whether the array of physical values (implemented as m_phys) is filled or no...
void SetPhysState(const bool physState)
This function manually sets whether the array of physical values (implemented as m_phys) is filled o...
int GetPhys_Offset(int n) const
Get the start offset position for a global list of m_phys correspoinding to element n.
std::shared_ptr< FieldIO > FieldIOSharedPtr
std::shared_ptr< SessionReader > SessionReaderSharedPtr
static FieldMetaDataMap NullFieldMetaDataMap
std::shared_ptr< Expansion > ExpansionSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
The namespace associated with the the StdRegions library (StdRegions introduction)
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
std::map< ConstFactorType, NekDouble > ConstFactorMap
static ConstFactorMap NullConstFactorMap
static VarCoeffMap NullVarCoeffMap
The above copyright notice and this permission notice shall be included.
NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag > DNekScalMat
NekMatrix< NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag >, BlockMatrixTag > DNekScalBlkMat
void Zero(int n, T *x, const int incx)
Zero vector.
Used to lookup the create function in NekManager.