Nektar++
ProcessStreamFunction.cpp
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: ProcessStreamFunction.cpp
4//
5// For more information, please see: http://www.nektar.info/
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Computes stream-function for 2D field.
32//
33////////////////////////////////////////////////////////////////////////////////
34
36
37using namespace std;
38
39namespace Nektar::FieldUtils
40{
41
44 ModuleKey(eProcessModule, "streamfunction"),
46 "Computes stream-function for 2D field.");
47
49 : ProcessModule(f)
50{
51 m_f->m_declareExpansionAsContField = true;
52}
53
55{
56}
57
58void ProcessStreamFunction::v_Process(po::variables_map &vm)
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}
106
107} // namespace Nektar::FieldUtils
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
FieldSharedPtr m_f
Field object.
Definition: Module.h:239
Abstract base class for processing modules.
Definition: Module.h:301
static std::shared_ptr< Module > create(FieldSharedPtr f)
Creates an instance of this class.
void v_Process(po::variables_map &vm) override
Write mesh to output file.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
std::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:990
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:180
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:47
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:402
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