Nektar++
ProcessEquiSpacedOutput.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: ProcessEquiSpacedOutput.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: Set up fields as interpolation to equispaced output
33 //
34 ////////////////////////////////////////////////////////////////////////////////
35 #include <string>
36 #include <iostream>
37 using namespace std;
38 
40 
44 #include <boost/math/special_functions/fpclassify.hpp>
45 #include <StdRegions/StdQuadExp.h>
46 #include <StdRegions/StdTriExp.h>
47 
48 namespace Nektar
49 {
50 namespace Utilities
51 {
52 
53 ModuleKey ProcessEquiSpacedOutput::className =
55  ModuleKey(eProcessModule, "equispacedoutput"),
56  ProcessEquiSpacedOutput::create,
57  "Write data as equi-spaced output using simplices to represent the data for connecting points");
58 
59 
60 ProcessEquiSpacedOutput::ProcessEquiSpacedOutput(FieldSharedPtr f)
61  : ProcessModule(f)
62 {
63  f->m_setUpEquiSpacedFields = true;
64 
65  m_config["tetonly"] = ConfigOption(true, "NotSet",
66  "Only process tetrahedral elements");
67 
68  m_config["modalenergy"] = ConfigOption(true,"NotSet","Write output as modal energy");
69 
70 }
71 
73 {
74 }
75 
76 void ProcessEquiSpacedOutput::Process(po::variables_map &vm)
77 {
79 }
80 
81 
83 {
84 
85  if(m_f->m_verbose)
86  {
87  cout << "Interpolating fields to equispaced" << endl;
88  }
89 
90  int coordim = m_f->m_exp[0]->GetCoordim(0);
91  int shapedim = m_f->m_exp[0]->GetExp(0)->GetShapeDimension();
92  int npts = m_f->m_exp[0]->GetTotPoints();
94 
95  int nel = m_f->m_exp[0]->GetExpSize();
96 
97  // set up the number of points in each element
98  int newpoints = 0;
99  int newtotpoints = 0;
100 
101  Array<OneD,int> conn;
102  int prevNcoeffs = 0;
103  int prevNpoints = 0;
104  int cnt = 0;
105 
106  // identify face 1 connectivity for prisms
107  map<int,StdRegions::Orientation > face0orient;
108  set<int> prismorient;
110 
111  // prepare PtsField
112  vector<std::string> fieldNames;
113  vector<int> ppe;
114  vector<Array<OneD, int> > ptsConn;
115  int nfields;
116 
117  for(int i = 0; i < nel; ++i)
118  {
119  e = m_f->m_exp[0]->GetExp(i);
120  if(e->DetShapeType() == LibUtilities::ePrism)
121  {
122  StdRegions::Orientation forient = e->GetForient(0);
123  int fid = e->GetGeom()->GetFid(0);
124  if(face0orient.count(fid))
125  { // face 1 meeting face 1 so reverse this id
126  prismorient.insert(i);
127  }
128  else
129  {
130  // just store if Dir 1 is fwd or bwd
131  if((forient == StdRegions::eDir1BwdDir1_Dir2FwdDir2) ||
135  {
136  face0orient[fid] = StdRegions::eBwd;
137  }
138  else
139  {
140  face0orient[fid] = StdRegions::eFwd;
141  }
142  }
143  }
144  }
145 
146  for(int i = 0; i < nel; ++i)
147  {
148  e = m_f->m_exp[0]->GetExp(i);
149  if(e->DetShapeType() == LibUtilities::ePrism)
150  {
151  int fid = e->GetGeom()->GetFid(2);
152  // check to see if face 2 meets face 1
153  if(face0orient.count(fid))
154  {
155  // check to see how face 2 is orientated
156  StdRegions::Orientation forient2 = e->GetForient(2);
157  StdRegions::Orientation forient0 = face0orient[fid];
158 
159  // If dir 1 or forient2 is bwd then check agains
160  // face 1 value
161  if((forient2 == StdRegions::eDir1BwdDir1_Dir2FwdDir2) ||
162  (forient2 == StdRegions::eDir1BwdDir1_Dir2BwdDir2) ||
163  (forient2 == StdRegions::eDir1BwdDir2_Dir2FwdDir1) ||
165  {
166  if(forient0 == StdRegions::eFwd)
167  {
168  prismorient.insert(i);
169  }
170  }
171  else
172  {
173  if(forient0 == StdRegions::eBwd)
174  {
175  prismorient.insert(i);
176  }
177  }
178  }
179  }
180  }
181 
182  for(int i = 0; i < nel; ++i)
183  {
184  e = m_f->m_exp[0]->GetExp(i);
185  if(m_config["tetonly"].m_beenSet)
186  {
187  if(m_f->m_exp[0]->GetExp(i)->DetShapeType() !=
189  {
190  continue;
191  }
192  }
193 
194  switch(e->DetShapeType())
195  {
197  {
198  int npoints0 = e->GetBasis(0)->GetNumPoints();
199 
200  newpoints = LibUtilities::StdSegData::
201  getNumberOfCoefficients(npoints0);
202  }
203  break;
205  {
206  int np0 = e->GetBasis(0)->GetNumPoints();
207  int np1 = e->GetBasis(1)->GetNumPoints();
208  int np = max(np0,np1);
209  newpoints = LibUtilities::StdTriData::
211  }
212  break;
214  {
215  int np0 = e->GetBasis(0)->GetNumPoints();
216  int np1 = e->GetBasis(1)->GetNumPoints();
217  int np = max(np0,np1);
218 
219  newpoints = LibUtilities::StdQuadData::
221  }
222  break;
224  {
225  int np0 = e->GetBasis(0)->GetNumPoints();
226  int np1 = e->GetBasis(1)->GetNumPoints();
227  int np2 = e->GetBasis(2)->GetNumPoints();
228  int np = max(np0,max(np1,np2));
229 
230  newpoints = LibUtilities::StdTetData::
231  getNumberOfCoefficients(np,np,np);
232  }
233  break;
235  {
236  int np0 = e->GetBasis(0)->GetNumPoints();
237  int np1 = e->GetBasis(1)->GetNumPoints();
238  int np2 = e->GetBasis(2)->GetNumPoints();
239  int np = max(np0,max(np1,np2));
240 
241  newpoints = LibUtilities::StdPrismData::
242  getNumberOfCoefficients(np,np,np);
243  }
244  break;
246  {
247  int np0 = e->GetBasis(0)->GetNumPoints();
248  int np1 = e->GetBasis(1)->GetNumPoints();
249  int np2 = e->GetBasis(2)->GetNumPoints();
250  int np = max(np0,max(np1,np2));
251 
252  newpoints = LibUtilities::StdPyrData::
253  getNumberOfCoefficients(np,np,np);
254  }
255  break;
257  {
258  int np0 = e->GetBasis(0)->GetNumPoints();
259  int np1 = e->GetBasis(1)->GetNumPoints();
260  int np2 = e->GetBasis(2)->GetNumPoints();
261  int np = max(np0,max(np1,np2));
262 
263  newpoints = LibUtilities::StdPyrData::
264  getNumberOfCoefficients(np,np,np);
265  }
266  break;
267  default:
268  {
269  ASSERTL0(false,"Points not known");
270  }
271  }
272 
273  ppe.push_back(newpoints);
274  newtotpoints += newpoints;
275 
276 
277  if(e->DetShapeType() == LibUtilities::ePrism)
278  {
279  bool standard = true;
280 
281  if(prismorient.count(i))
282  {
283  standard = false; // reverse direction
284  }
285 
286  e->GetSimplexEquiSpacedConnectivity(conn,standard);
287  }
288  else
289  {
290 
291  if((prevNcoeffs != e->GetNcoeffs()) ||
292  (prevNpoints != e->GetTotPoints()))
293  {
294  prevNcoeffs = e->GetNcoeffs();
295  prevNpoints = e->GetTotPoints();
296 
297  e->GetSimplexEquiSpacedConnectivity(conn);
298  }
299  }
300  Array<OneD, int> newconn(conn.num_elements());
301  for(int j = 0; j < conn.num_elements(); ++j)
302  {
303  newconn[j] = conn[j] + cnt;
304  }
305 
306  ptsConn.push_back(newconn);
307  cnt += newpoints;
308  }
309 
310  if(m_f->m_fielddef.size())
311  {
312  nfields = m_f->m_exp.size();
313  }
314  else // just the mesh points
315  {
316  nfields = 0;
317  }
318 
319  Array<OneD, Array<OneD, NekDouble> > pts(nfields + coordim);
320 
321  for(int i = 0; i < nfields + coordim; ++i)
322  {
323  pts[i] = Array<OneD, NekDouble>(newtotpoints);
324  }
325 
326  // Interpolate coordinates
327  for(int i = 0; i < coordim; ++i)
328  {
329  coords[i] = Array<OneD, NekDouble>(npts);
330  }
331 
332  for(int i = coordim; i < 3; ++i)
333  {
334  coords[i] = NullNekDouble1DArray;
335  }
336 
337  m_f->m_exp[0]->GetCoords(coords[0],coords[1],coords[2]);
338 
339  int nq1 = m_f->m_exp[0]->GetTotPoints();
340 
341  Array<OneD, NekDouble> x1(nq1);
342  Array<OneD, NekDouble> y1(nq1);
343  Array<OneD, NekDouble> z1(nq1);
344 
345  m_f->m_exp[0]->GetCoords(x1, y1, z1);
346 
347 
349 
350  for(int n = 0; n < coordim; ++n)
351  {
352  cnt = 0;
353  int cnt1 = 0;
354  for(int i = 0; i < nel; ++i)
355  {
356  m_f->m_exp[0]->GetExp(i)->PhysInterpToSimplexEquiSpaced(
357  coords[n] + cnt,
358  tmp = pts[n] + cnt1);
359  cnt1 += ppe[i];
360  cnt += m_f->m_exp[0]->GetExp(i)->GetTotPoints();
361  }
362  }
363 
364  if(m_f->m_fielddef.size())
365  {
366  ASSERTL0(m_f->m_fielddef[0]->m_fields.size() == m_f->m_exp.size(),
367  "More expansion defined than fields");
368 
369  for(int n = 0; n < m_f->m_exp.size(); ++n)
370  {
371  cnt = 0;
372  int cnt1 = 0;
373 
374  if(m_config["modalenergy"].m_beenSet)
375  {
376  Array<OneD, const NekDouble> phys = m_f->m_exp[n]->GetPhys();
377  for(int i = 0; i < nel; ++i)
378  {
379  GenOrthoModes(i,phys+cnt,tmp = pts[coordim + n] + cnt1);
380  cnt1 += ppe[i];
381  cnt += m_f->m_exp[0]->GetExp(i)->GetTotPoints();
382  }
383  }
384  else
385  {
386  Array<OneD, const NekDouble> phys = m_f->m_exp[n]->GetPhys();
387  for(int i = 0; i < nel; ++i)
388  {
389  m_f->m_exp[0]->GetExp(i)->PhysInterpToSimplexEquiSpaced(
390  phys + cnt,
391  tmp = pts[coordim + n] + cnt1);
392  cnt1 += ppe[i];
393  cnt += m_f->m_exp[0]->GetExp(i)->GetTotPoints();
394  }
395  }
396 
397  // Set up Variable string.
398  fieldNames.push_back(m_f->m_fielddef[0]->m_fields[n]);
399  }
400  }
401 
402  m_f->m_fieldPts = MemoryManager<LibUtilities::PtsField>::AllocateSharedPtr(coordim, fieldNames, pts);
403  if (shapedim == 2)
404  {
405  m_f->m_fieldPts->SetPtsType(LibUtilities::ePtsTriBlock);
406  }
407  else if (shapedim == 3)
408  {
409  m_f->m_fieldPts->SetPtsType(LibUtilities::ePtsTetBlock);
410  }
411  m_f->m_fieldPts->SetConnectivity(ptsConn);
412 }
413 
414 
416  int n,
417  const Array<OneD,const NekDouble> &phys,
418  Array<OneD, NekDouble> &coeffs)
419  {
421  e = m_f->m_exp[0]->GetExp(n);
422 
423  switch(e->DetShapeType())
424  {
426  {
427  int np0 = e->GetBasis(0)->GetNumPoints();
428  int np1 = e->GetBasis(1)->GetNumPoints();
429  int np = max(np0,np1);
430 
431  // to ensure points are correctly projected to np need
432  // to increase the order slightly of coordinates
433  LibUtilities::PointsKey pa(np+1,e->GetPointsType(0));
434  LibUtilities::PointsKey pb(np,e->GetPointsType(1));
435  Array<OneD, NekDouble> tophys(np*(np+1));
436 
439  StdRegions::StdTriExp OrthoExp(Ba,Bb);
440 
441  // interpolate points to new phys points!
442  LibUtilities::Interp2D(e->GetBasis(0)->GetBasisKey(),
443  e->GetBasis(1)->GetBasisKey(),
444  phys,Ba,Bb,tophys);
445 
446  OrthoExp.FwdTrans(tophys,coeffs);
447  break;
448  }
450  {
451  int np0 = e->GetBasis(0)->GetNumPoints();
452  int np1 = e->GetBasis(1)->GetNumPoints();
453  int np = max(np0,np1);
454 
455  LibUtilities::PointsKey pa(np+1,e->GetPointsType(0));
456  LibUtilities::PointsKey pb(np+1,e->GetPointsType(1));
457  Array<OneD, NekDouble> tophys((np+1)*(np+1));
458 
461  StdRegions::StdQuadExp OrthoExp(Ba,Bb);
462 
463  // interpolate points to new phys points!
464  LibUtilities::Interp2D(e->GetBasis(0)->GetBasisKey(),
465  e->GetBasis(1)->GetBasisKey(),
466  phys,Ba,Bb,tophys);
467 
468  OrthoExp.FwdTrans(phys,coeffs);
469  break;
470  }
471  default:
472  ASSERTL0(false,"Shape needs setting up");
473  break;
474  }
475  }
476 
477 }
478 }
479 
480 
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
int getNumberOfCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:286
pair< ModuleType, string > ModuleKey
static Array< OneD, NekDouble > NullNekDouble1DArray
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
virtual void Process()=0
map< string, ConfigOption > m_config
List of configuration values.
STL namespace.
void GenOrthoModes(int n, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &coeffs)
int getNumberOfCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:186
FieldSharedPtr m_f
Field object.
Principle Orthogonal Functions .
Definition: BasisType.h:47
void Interp2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
this function interpolates a 2D function evaluated at the quadrature points of the 2D basis...
Definition: Interp.cpp:116
static std::string npts
Definition: InputFld.cpp:43
int getNumberOfCoefficients(int Na, int Nb)
Definition: ShapeType.hpp:111
Principle Orthogonal Functions .
Definition: BasisType.h:46
Defines a specification for a set of points.
Definition: Points.h:58
boost::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:677
Represents a command-line configuration option.
int getNumberOfCoefficients(int Na, int Nb, int Nc)
Definition: ShapeType.hpp:232
int getNumberOfCoefficients(int Na, int Nb)
Definition: ShapeType.hpp:132
boost::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
Describes the specification for a Basis.
Definition: Basis.h:50
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
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Forward transformation from physical space to coefficient space...