Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | List of all members
Nektar::FieldUtils::ProcessStreamFunction Class Reference

This processing module calculates the stream function of a 2D field and adds it as an extra-field to the output file. More...

#include <ProcessStreamFunction.h>

Inheritance diagram for Nektar::FieldUtils::ProcessStreamFunction:
[legend]

Public Member Functions

 ProcessStreamFunction (FieldSharedPtr f)
 
 ~ProcessStreamFunction () override
 
- Public Member Functions inherited from Nektar::FieldUtils::ProcessModule
 ProcessModule ()
 
 ProcessModule (FieldSharedPtr p_f)
 
- Public Member Functions inherited from Nektar::FieldUtils::Module
FIELD_UTILS_EXPORT Module (FieldSharedPtr p_f)
 
virtual ~Module ()=default
 
void Process (po::variables_map &vm)
 
std::string GetModuleName ()
 
std::string GetModuleDescription ()
 
const ConfigOptionGetConfigOption (const std::string &key) const
 
ModulePriority GetModulePriority ()
 
std::vector< ModuleKeyGetModulePrerequisites ()
 
FIELD_UTILS_EXPORT void RegisterConfig (std::string key, std::string value="")
 Register a configuration option with a module. More...
 
FIELD_UTILS_EXPORT void PrintConfig ()
 Print out all configuration options for a module. More...
 
FIELD_UTILS_EXPORT void SetDefaults ()
 Sets default configuration options for those which have not been set. More...
 
FIELD_UTILS_EXPORT void AddFile (std::string fileType, std::string fileName)
 
FIELD_UTILS_EXPORT void EvaluateTriFieldAtEquiSpacedPts (LocalRegions::ExpansionSharedPtr &exp, const Array< OneD, const NekDouble > &infield, Array< OneD, NekDouble > &outfield)
 

Static Public Member Functions

static std::shared_ptr< Modulecreate (FieldSharedPtr f)
 Creates an instance of this class. More...
 

Static Public Attributes

static ModuleKey className
 

Protected Member Functions

void v_Process (po::variables_map &vm) override
 Write mesh to output file. More...
 
std::string v_GetModuleName () override
 
std::string v_GetModuleDescription () override
 
ModulePriority v_GetModulePriority () override
 
- Protected Member Functions inherited from Nektar::FieldUtils::Module
 Module ()
 
virtual void v_Process (po::variables_map &vm)
 
virtual std::string v_GetModuleName ()
 
virtual std::string v_GetModuleDescription ()
 
virtual ModulePriority v_GetModulePriority ()
 
virtual std::vector< ModuleKeyv_GetModulePrerequisites ()
 

Additional Inherited Members

- Public Attributes inherited from Nektar::FieldUtils::Module
FieldSharedPtr m_f
 Field object. More...
 
- Protected Attributes inherited from Nektar::FieldUtils::Module
std::map< std::string, ConfigOptionm_config
 List of configuration values. More...
 
std::set< std::string > m_allowedFiles
 List of allowed file formats. More...
 

Detailed Description

This processing module calculates the stream function of a 2D field and adds it as an extra-field to the output file.

Definition at line 46 of file ProcessStreamFunction.h.

Constructor & Destructor Documentation

◆ ProcessStreamFunction()

Nektar::FieldUtils::ProcessStreamFunction::ProcessStreamFunction ( FieldSharedPtr  f)

Definition at line 48 of file ProcessStreamFunction.cpp.

49 : ProcessModule(f)
50{
51 m_f->m_declareExpansionAsContField = true;
52}
FieldSharedPtr m_f
Field object.
Definition: Module.h:239

References Nektar::FieldUtils::Module::m_f.

◆ ~ProcessStreamFunction()

Nektar::FieldUtils::ProcessStreamFunction::~ProcessStreamFunction ( )
override

Definition at line 54 of file ProcessStreamFunction.cpp.

55{
56}

Member Function Documentation

◆ create()

static std::shared_ptr< Module > Nektar::FieldUtils::ProcessStreamFunction::create ( FieldSharedPtr  f)
inlinestatic

Creates an instance of this class.

Definition at line 50 of file ProcessStreamFunction.h.

51 {
53 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

◆ v_GetModuleDescription()

std::string Nektar::FieldUtils::ProcessStreamFunction::v_GetModuleDescription ( )
inlineoverrideprotectedvirtual

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 68 of file ProcessStreamFunction.h.

69 {
70 return "Calculating stream-function";
71 }

◆ v_GetModuleName()

std::string Nektar::FieldUtils::ProcessStreamFunction::v_GetModuleName ( )
inlineoverrideprotectedvirtual

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 63 of file ProcessStreamFunction.h.

64 {
65 return "ProcessStreamFunction";
66 }

◆ v_GetModulePriority()

ModulePriority Nektar::FieldUtils::ProcessStreamFunction::v_GetModulePriority ( )
inlineoverrideprotectedvirtual

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 73 of file ProcessStreamFunction.h.

74 {
75 return eModifyExp;
76 }

References Nektar::FieldUtils::eModifyExp.

◆ v_Process()

void Nektar::FieldUtils::ProcessStreamFunction::v_Process ( po::variables_map &  vm)
overrideprotectedvirtual

Write mesh to output file.

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 58 of file ProcessStreamFunction.cpp.

59{
60 m_f->SetUpExp(vm);
61
62 int nfields = m_f->m_variables.size();
63 // Append field name
64 m_f->m_variables.push_back("StreamFunc");
65
66 // Skip in case of empty partition
67 if (m_f->m_exp[0]->GetNumElmts() == 0)
68 {
69 return;
70 }
71
72 // Check if dimension is correct
73 int expdim = m_f->m_graph->GetMeshDimension();
74 int spacedim = expdim + m_f->m_numHomogeneousDir;
75 ASSERTL0(spacedim == 2, "Stream function can only be obtained in 2D.");
76
77 // Allocate arrays
78 int npoints = m_f->m_exp[0]->GetNpoints();
79 Array<OneD, NekDouble> vx(npoints);
80 Array<OneD, NekDouble> uy(npoints);
81 Array<OneD, NekDouble> vortNeg(npoints);
82
83 // Resize expansion
84 m_f->m_exp.resize(nfields + 1);
85 m_f->m_exp[nfields] = m_f->AppendExpList(m_f->m_numHomogeneousDir);
86
87 // Calculate vorticity: -W_z = Uy - Vx
88 m_f->m_exp[0]->PhysDeriv(MultiRegions::DirCartesianMap[1],
89 m_f->m_exp[0]->GetPhys(), uy);
90 m_f->m_exp[1]->PhysDeriv(MultiRegions::DirCartesianMap[0],
91 m_f->m_exp[1]->GetPhys(), vx);
92
93 Vmath::Vsub(npoints, uy, 1, vx, 1, vortNeg, 1);
94
95 // Calculate the Stream Function as the solution of the
96 // Poisson equation: \nabla^2 StreamFunction = -Vorticity
97 // Use same boundary conditions as v
99 factor[StdRegions::eFactorLambda] = 0.0;
100
101 m_f->m_exp[1]->HelmSolve(vortNeg, m_f->m_exp[nfields]->UpdateCoeffs(),
102 factor);
103 m_f->m_exp[1]->BwdTrans(m_f->m_exp[nfields]->GetCoeffs(),
104 m_f->m_exp[nfields]->UpdatePhys());
105}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:87
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:430
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.
Definition: Vmath.hpp:220

References ASSERTL0, Nektar::MultiRegions::DirCartesianMap, Nektar::StdRegions::eFactorLambda, Nektar::FieldUtils::Module::m_f, and Vmath::Vsub().

Member Data Documentation

◆ className

ModuleKey Nektar::FieldUtils::ProcessStreamFunction::className
static
Initial value:
=
ModuleKey(eProcessModule, "streamfunction"),
"Computes stream-function for 2D field.")
static std::shared_ptr< Module > create(FieldSharedPtr f)
Creates an instance of this class.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:180
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:47

Definition at line 54 of file ProcessStreamFunction.h.