35 #include <boost/core/ignore_unused.hpp> 46 ModuleKey ProcessStreamFunction::className =
49 ProcessStreamFunction::create,
50 "Computes stream-function for 2D field.");
55 m_f->m_declareExpansionAsContField =
true;
64 boost::ignore_unused(vm);
66 int nfields =
m_f->m_variables.size();
68 m_f->m_variables.push_back(
"StreamFunc");
71 if (
m_f->m_exp[0]->GetNumElmts() == 0)
77 int expdim =
m_f->m_graph->GetMeshDimension();
78 int spacedim = expdim +
m_f->m_numHomogeneousDir;
79 ASSERTL0(spacedim == 2,
"Stream function can only be obtained in 2D.");
82 int npoints =
m_f->m_exp[0]->GetNpoints();
88 m_f->m_exp.resize(nfields + 1);
89 m_f->m_exp[nfields] =
m_f->AppendExpList(
m_f->m_numHomogeneousDir);
93 m_f->m_exp[0]->GetPhys(),uy);
95 m_f->m_exp[1]->GetPhys(),vx);
105 m_f->m_exp[1]->HelmSolve(vortNeg,
106 m_f->m_exp[nfields]->UpdateCoeffs(),
108 m_f->m_exp[1]->BwdTrans(
m_f->m_exp[nfields]->GetCoeffs(),
109 m_f->m_exp[nfields]->UpdatePhys());
#define ASSERTL0(condition, msg)
std::map< ConstFactorType, NekDouble > ConstFactorMap
std::shared_ptr< Field > FieldSharedPtr
std::pair< ModuleType, std::string > ModuleKey
virtual ~ProcessStreamFunction()
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.
MultiRegions::Direction const DirCartesianMap[]
virtual void Process(po::variables_map &vm)
Write mesh to output file.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Abstract base class for processing modules.
static FlagList NullFlagList
An empty flag list.
ModuleFactory & GetModuleFactory()
FieldSharedPtr m_f
Field object.