Nektar++
ProcessVorticity.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: ProcessVorticity.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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Computes vorticity field.
33 //
34 ////////////////////////////////////////////////////////////////////////////////
35 
36 #include <string>
37 #include <iostream>
38 using namespace std;
39 
40 #include "ProcessVorticity.h"
41 
44 
45 namespace Nektar
46 {
47 namespace Utilities
48 {
49 
50 ModuleKey ProcessVorticity::className =
52  ModuleKey(eProcessModule, "vorticity"),
53  ProcessVorticity::create, "Computes vorticity field.");
54 
55 ProcessVorticity::ProcessVorticity(FieldSharedPtr f) : ProcessModule(f)
56 {
57 }
58 
60 {
61 }
62 
63 void ProcessVorticity::Process(po::variables_map &vm)
64 {
65  if (m_f->m_verbose)
66  {
67  cout << "ProcessVorticity: Calculating vorticity..." << endl;
68  }
69 
70  int i, j, s;
71  int expdim = m_f->m_graph->GetMeshDimension();
72  int spacedim = expdim;
73  if ((m_f->m_fielddef[0]->m_numHomogeneousDir) == 1 ||
74  (m_f->m_fielddef[0]->m_numHomogeneousDir) == 2)
75  {
76  spacedim = 3;
77  }
78  int nfields = m_f->m_fielddef[0]->m_fields.size();
79  if (spacedim == 1)
80  {
81  ASSERTL0(false, "Error: Vorticity for a 1D problem cannot "
82  "be computed")
83  }
84  int addfields = (spacedim == 2)? 1:3;
85 
86  int npoints = m_f->m_exp[0]->GetNpoints();
87  Array<OneD, Array<OneD, NekDouble> > grad(nfields*nfields);
88  Array<OneD, Array<OneD, NekDouble> > outfield(addfields);
89 
90  int nstrips;
91 
92  m_f->m_session->LoadParameter("Strip_Z",nstrips,1);
93 
94  m_f->m_exp.resize(nfields*nstrips);
95 
96  for (i = 0; i < nfields*nfields; ++i)
97  {
98  grad[i] = Array<OneD, NekDouble>(npoints);
99  }
100 
101  for (i = 0; i < addfields; ++i)
102  {
103  outfield[i] = Array<OneD, NekDouble>(npoints);
104  }
105 
106  vector<MultiRegions::ExpListSharedPtr> Exp(nstrips*addfields);
107 
108  for(s = 0; s < nstrips; ++s) //homogeneous strip varient
109  {
110  // Calculate Gradient & Vorticity
111  if (spacedim == 2)
112  {
113  for (i = 0; i < nfields; ++i)
114  {
115  m_f->m_exp[s*nfields+i]->PhysDeriv(m_f->m_exp[s*nfields+i]->GetPhys(),
116  grad[i*nfields],
117  grad[i*nfields+1]);
118  }
119  // W_z = Vx - Uy
120  Vmath::Vsub(npoints, grad[1*nfields+0], 1,
121  grad[0*nfields+1], 1,
122  outfield[0], 1);
123  }
124  else
125  {
126  for (i = 0; i < nfields; ++i)
127  {
128 
129  m_f->m_exp[s*nfields+i]->PhysDeriv(m_f->m_exp[s*nfields+i]->GetPhys(),
130  grad[i*nfields],
131  grad[i*nfields+1],
132  grad[i*nfields+2]);
133  }
134 
135  // W_x = Wy - Vz
136  Vmath::Vsub(npoints, grad[2*nfields+1], 1, grad[1*nfields+2], 1,
137  outfield[0],1);
138  // W_y = Uz - Wx
139  Vmath::Vsub(npoints, grad[0*nfields+2], 1, grad[2*nfields+0], 1,
140  outfield[1], 1);
141  // W_z = Vx - Uy
142  Vmath::Vsub(npoints, grad[1*nfields+0], 1, grad[0*nfields+1], 1,
143  outfield[2], 1);
144  }
145 
146  for (i = 0; i < addfields; ++i)
147  {
148  int n = s*addfields + i;
149  Exp[n] = m_f->AppendExpList(m_f->m_fielddef[0]->m_numHomogeneousDir);
150  Exp[n]->UpdatePhys() = outfield[i];
151  Exp[n]->FwdTrans_IterPerExp(outfield[i],
152  Exp[n]->UpdateCoeffs());
153  }
154  }
155 
157  for(s = 0; s < nstrips; ++s)
158  {
159  for(i = 0; i < addfields; ++i)
160  {
161  it = m_f->m_exp.begin()+s*(nfields+addfields)+nfields+i;
162  m_f->m_exp.insert(it, Exp[s*addfields+i]);
163  }
164  }
165 
166  vector<string > outname;
167  if (addfields == 1)
168  {
169  outname.push_back("W_z");
170  }
171  else
172  {
173  outname.push_back("W_x");
174  outname.push_back("W_y");
175  outname.push_back("W_z");
176  }
177 
178  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef
179  = m_f->m_exp[0]->GetFieldDefinitions();
180  std::vector<std::vector<NekDouble> > FieldData(FieldDef.size());
181 
182  for(s = 0; s < nstrips; ++s) //homogeneous strip varient
183  {
184  for (j = 0; j < nfields + addfields; ++j)
185  {
186  for (i = 0; i < FieldDef.size()/nstrips; ++i)
187  {
188  int n = s * FieldDef.size()/nstrips + i;
189 
190  if (j >= nfields)
191  {
192  FieldDef[n]->m_fields.push_back(outname[j-nfields]);
193  }
194  else
195  {
196  FieldDef[n]->m_fields.push_back(m_f->m_fielddef[0]->m_fields[j]);
197  }
198  m_f->m_exp[s*(nfields + addfields)+j]->AppendFieldData(FieldDef[n], FieldData[n]);
199  }
200  }
201  }
202 
203  m_f->m_fielddef = FieldDef;
204  m_f->m_data = FieldData;
205 }
206 
207 }
208 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
pair< ModuleType, string > ModuleKey
virtual void Process()=0
STL namespace.
FieldSharedPtr m_f
Field object.
boost::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:677
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:329
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
ModuleFactory & GetModuleFactory()
Abstract base class for processing modules.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215