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 
35 #include <boost/core/ignore_unused.hpp>
36 
37 #include "ProcessStreamFunction.h"
38 
39 using namespace std;
40 
41 namespace Nektar
42 {
43 namespace FieldUtils
44 {
45 
46 ModuleKey ProcessStreamFunction::className =
48  ModuleKey(eProcessModule, "streamfunction"),
49  ProcessStreamFunction::create,
50  "Computes stream-function for 2D field.");
51 
52 ProcessStreamFunction::ProcessStreamFunction(FieldSharedPtr f)
53  : ProcessModule(f)
54 {
55  m_f->m_declareExpansionAsContField = true;
56 }
57 
59 {
60 }
61 
62 void ProcessStreamFunction::Process(po::variables_map &vm)
63 {
64  m_f->SetUpExp(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(), factor);
107  m_f->m_exp[1]->BwdTrans(m_f->m_exp[nfields]->GetCoeffs(),
108  m_f->m_exp[nfields]->UpdatePhys());
109 }
110 
111 }
112 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
FieldSharedPtr m_f
Field object.
Definition: Module.h:230
Abstract base class for processing modules.
Definition: Module.h:265
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.
Definition: NekFactory.hpp:200
std::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:989
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:290
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:49
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:90
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:314
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
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:372