Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | 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)
 
virtual ~ProcessStreamFunction ()
 
virtual void Process (po::variables_map &vm)
 Write mesh to output file. More...
 
virtual std::string GetModuleName ()
 
virtual std::string GetModuleDescription ()
 
virtual ModulePriority GetModulePriority ()
 
- 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)
 
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 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
 

Additional Inherited Members

- Protected Member Functions inherited from Nektar::FieldUtils::Module
 Module ()
 
- Protected Attributes inherited from Nektar::FieldUtils::Module
FieldSharedPtr m_f
 Field object. More...
 
std::map< std::string, ConfigOptionm_config
 List of configuration values. 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 48 of file ProcessStreamFunction.h.

Constructor & Destructor Documentation

◆ ProcessStreamFunction()

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

Definition at line 52 of file ProcessStreamFunction.cpp.

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

53  : ProcessModule(f)
54 {
55  m_f->m_declareExpansionAsContField = true;
56 }
FieldSharedPtr m_f
Field object.

◆ ~ProcessStreamFunction()

Nektar::FieldUtils::ProcessStreamFunction::~ProcessStreamFunction ( )
virtual

Definition at line 58 of file ProcessStreamFunction.cpp.

59 {
60 }

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 52 of file ProcessStreamFunction.h.

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

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

◆ GetModuleDescription()

virtual std::string Nektar::FieldUtils::ProcessStreamFunction::GetModuleDescription ( )
inlinevirtual

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 69 of file ProcessStreamFunction.h.

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

◆ GetModuleName()

virtual std::string Nektar::FieldUtils::ProcessStreamFunction::GetModuleName ( )
inlinevirtual

Implements Nektar::FieldUtils::Module.

Definition at line 64 of file ProcessStreamFunction.h.

65  {
66  return "ProcessStreamFunction";
67  }

◆ GetModulePriority()

virtual ModulePriority Nektar::FieldUtils::ProcessStreamFunction::GetModulePriority ( )
inlinevirtual

Implements Nektar::FieldUtils::Module.

Definition at line 74 of file ProcessStreamFunction.h.

References Nektar::FieldUtils::eModifyExp.

◆ Process()

void Nektar::FieldUtils::ProcessStreamFunction::Process ( po::variables_map &  vm)
virtual

Write mesh to output file.

Implements Nektar::FieldUtils::Module.

Definition at line 62 of file ProcessStreamFunction.cpp.

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

63 {
64  boost::ignore_unused(vm);
65 
66  int nfields = m_f->m_variables.size();
67  // Append field name
68  m_f->m_variables.push_back("StreamFunc");
69 
70  // Skip in case of empty partition
71  if (m_f->m_exp[0]->GetNumElmts() == 0)
72  {
73  return;
74  }
75 
76  // Check if dimension is correct
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.");
80 
81  // Allocate arrays
82  int npoints = m_f->m_exp[0]->GetNpoints();
83  Array<OneD, NekDouble> vx(npoints);
84  Array<OneD, NekDouble> uy(npoints);
85  Array<OneD, NekDouble> vortNeg(npoints);
86 
87  // Resize expansion
88  m_f->m_exp.resize(nfields + 1);
89  m_f->m_exp[nfields] = m_f->AppendExpList(m_f->m_numHomogeneousDir);
90 
91  // Calculate vorticity: -W_z = Uy - Vx
92  m_f->m_exp[0]->PhysDeriv(MultiRegions::DirCartesianMap[1],
93  m_f->m_exp[0]->GetPhys(),uy);
94  m_f->m_exp[1]->PhysDeriv(MultiRegions::DirCartesianMap[0],
95  m_f->m_exp[1]->GetPhys(),vx);
96 
97  Vmath::Vsub(npoints, uy, 1, vx, 1, vortNeg, 1);
98 
99  // Calculate the Stream Function as the solution of the
100  // Poisson equation: \nabla^2 StreamFunction = -Vorticity
101  // Use same boundary conditions as v
103  factor[StdRegions::eFactorLambda] = 0.0;
104 
105  m_f->m_exp[1]->HelmSolve(vortNeg,
106  m_f->m_exp[nfields]->UpdateCoeffs(),
107  NullFlagList, factor);
108  m_f->m_exp[1]->BwdTrans(m_f->m_exp[nfields]->GetCoeffs(),
109  m_f->m_exp[nfields]->UpdatePhys());
110 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:294
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.cpp:346
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:88
static FlagList NullFlagList
An empty flag list.
FieldSharedPtr m_f
Field object.

Member Data Documentation

◆ className

ModuleKey Nektar::FieldUtils::ProcessStreamFunction::className
static
Initial value:
=
ModuleKey(eProcessModule, "streamfunction"),
"Computes stream-function for 2D field.")

Definition at line 56 of file ProcessStreamFunction.h.