Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 
43 #include <boost/math/special_functions/fpclassify.hpp>
44 
45 namespace Nektar
46 {
47 namespace Utilities
48 {
49 
50 ModuleKey ProcessEquiSpacedOutput::className =
52  ModuleKey(eProcessModule, "equispacedoutput"),
53  ProcessEquiSpacedOutput::create,
54  "Write data as equi-spaced output using simplices to represent the data for connecting points");
55 
56 
57 ProcessEquiSpacedOutput::ProcessEquiSpacedOutput(FieldSharedPtr f)
58  : ProcessModule(f)
59 {
61  f->m_setUpEquiSpacedFields = true;
62 
63  m_config["tetonly"] = ConfigOption(true, "NotSet",
64  "Only process tetrahedral elements");
65 
66 }
67 
69 {
70 }
71 
72 void ProcessEquiSpacedOutput::Process(po::variables_map &vm)
73 {
75 }
76 
77 
79 {
80 
81  if(m_f->m_verbose)
82  {
83  cout << "Interpolating fields to equispaced" << endl;
84  }
85 
86  int coordim = m_f->m_exp[0]->GetCoordim(0);
87  int shapedim = m_f->m_exp[0]->GetExp(0)->GetShapeDimension();
88  int npts = m_f->m_exp[0]->GetTotPoints();
89  Array<OneD, Array<OneD, NekDouble> > coords(3);
90 
91  int nel = m_f->m_exp[0]->GetExpSize();
92 
93  // set up the number of points in each element
94  int newpoints;
95  int newtotpoints = 0;
96 
97  Array<OneD,int> conn;
98  int prevNcoeffs = 0;
99  int prevNpoints = 0;
100  int cnt = 0;
101 
102  // identify face 1 connectivity for prisms
103  map<int,StdRegions::Orientation > face0orient;
104  set<int> prismorient;
106 
107  for(int i = 0; i < nel; ++i)
108  {
109  e = m_f->m_exp[0]->GetExp(i);
110  if(e->DetShapeType() == LibUtilities::ePrism)
111  {
112  StdRegions::Orientation forient = e->GetForient(0);
113  int fid = e->GetGeom()->GetFid(0);
114  if(face0orient.count(fid))
115  { // face 1 meeting face 1 so reverse this id
116  prismorient.insert(i);
117  }
118  else
119  {
120  // just store if Dir 1 is fwd or bwd
121  if((forient == StdRegions::eDir1BwdDir1_Dir2FwdDir2) ||
125  {
126  face0orient[fid] = StdRegions::eBwd;
127  }
128  else
129  {
130  face0orient[fid] = StdRegions::eFwd;
131  }
132  }
133  }
134  }
135 
136  for(int i = 0; i < nel; ++i)
137  {
138  e = m_f->m_exp[0]->GetExp(i);
139  if(e->DetShapeType() == LibUtilities::ePrism)
140  {
141  int fid = e->GetGeom()->GetFid(2);
142  // check to see if face 2 meets face 1
143  if(face0orient.count(fid))
144  {
145  // check to see how face 2 is orientated
146  StdRegions::Orientation forient2 = e->GetForient(2);
147  StdRegions::Orientation forient0 = face0orient[fid];
148 
149  // If dir 1 or forient2 is bwd then check agains
150  // face 1 value
151  if((forient2 == StdRegions::eDir1BwdDir1_Dir2FwdDir2) ||
152  (forient2 == StdRegions::eDir1BwdDir1_Dir2BwdDir2) ||
153  (forient2 == StdRegions::eDir1BwdDir2_Dir2FwdDir1) ||
155  {
156  if(forient0 == StdRegions::eFwd)
157  {
158  prismorient.insert(i);
159  }
160  }
161  else
162  {
163  if(forient0 == StdRegions::eBwd)
164  {
165  prismorient.insert(i);
166  }
167  }
168  }
169  }
170  }
171 
172  for(int i = 0; i < nel; ++i)
173  {
174  e = m_f->m_exp[0]->GetExp(i);
175  if(m_config["tetonly"].m_beenSet)
176  {
177  if(m_f->m_exp[0]->GetExp(i)->DetShapeType() !=
179  {
180  if(i == nel-1)
181  {
182  break;
183  }
184  else
185  {
186  continue;
187  }
188  }
189  }
190 
191  switch(e->DetShapeType())
192  {
194  {
195  int npoints0 = e->GetBasis(0)->GetNumPoints();
196 
197  newpoints = LibUtilities::StdSegData::
198  getNumberOfCoefficients(npoints0);
199  m_f->m_fieldPts->m_npts.push_back(newpoints);
200  newtotpoints += newpoints;
201  }
202  break;
204  {
205  int np0 = e->GetBasis(0)->GetNumPoints();
206  int np1 = e->GetBasis(1)->GetNumPoints();
207  int np = max(np0,np1);
208  newpoints = LibUtilities::StdTriData::
210  m_f->m_fieldPts->m_npts.push_back(newpoints);
211  newtotpoints += newpoints;
212  }
213  break;
215  {
216  int np0 = e->GetBasis(0)->GetNumPoints();
217  int np1 = e->GetBasis(1)->GetNumPoints();
218  int np = max(np0,np1);
219 
220  newpoints = LibUtilities::StdQuadData::
222  m_f->m_fieldPts->m_npts.push_back(newpoints);
223  newtotpoints += newpoints;
224  }
225  break;
227  {
228  int np0 = e->GetBasis(0)->GetNumPoints();
229  int np1 = e->GetBasis(1)->GetNumPoints();
230  int np2 = e->GetBasis(2)->GetNumPoints();
231  int np = max(np0,max(np1,np2));
232 
233  newpoints = LibUtilities::StdTetData::
234  getNumberOfCoefficients(np,np,np);
235  m_f->m_fieldPts->m_npts.push_back(newpoints);
236  newtotpoints += newpoints;
237  }
238  break;
240  {
241  int np0 = e->GetBasis(0)->GetNumPoints();
242  int np1 = e->GetBasis(1)->GetNumPoints();
243  int np2 = e->GetBasis(2)->GetNumPoints();
244  int np = max(np0,max(np1,np2));
245 
246  newpoints = LibUtilities::StdPrismData::
247  getNumberOfCoefficients(np,np,np);
248  m_f->m_fieldPts->m_npts.push_back(newpoints);
249  newtotpoints += newpoints;
250  }
251  break;
253  {
254  int np0 = e->GetBasis(0)->GetNumPoints();
255  int np1 = e->GetBasis(1)->GetNumPoints();
256  int np2 = e->GetBasis(2)->GetNumPoints();
257  int np = max(np0,max(np1,np2));
258 
259  newpoints = LibUtilities::StdPyrData::
260  getNumberOfCoefficients(np,np,np);
261  m_f->m_fieldPts->m_npts.push_back(newpoints);
262  newtotpoints += newpoints;
263  }
264  break;
266  {
267  int np0 = e->GetBasis(0)->GetNumPoints();
268  int np1 = e->GetBasis(1)->GetNumPoints();
269  int np2 = e->GetBasis(2)->GetNumPoints();
270  int np = max(np0,max(np1,np2));
271 
272  newpoints = LibUtilities::StdPyrData::
273  getNumberOfCoefficients(np,np,np);
274  m_f->m_fieldPts->m_npts.push_back(newpoints);
275  newtotpoints += newpoints;
276  }
277  break;
278  default:
279  {
280  ASSERTL0(false,"Points not known");
281  }
282  }
283 
284  if(e->DetShapeType() == LibUtilities::ePrism)
285  {
286  bool standard = true;
287 
288  if(prismorient.count(i))
289  {
290  standard = false; // reverse direction
291  }
292 
293  e->GetSimplexEquiSpacedConnectivity(conn,standard);
294  }
295  else
296  {
297 
298  if((prevNcoeffs != e->GetNcoeffs()) ||
299  (prevNpoints != e->GetTotPoints()))
300  {
301  prevNcoeffs = e->GetNcoeffs();
302  prevNpoints = e->GetTotPoints();
303 
304  e->GetSimplexEquiSpacedConnectivity(conn);
305  }
306  }
307  Array<OneD, int> newconn(conn.num_elements());
308  for(int j = 0; j < conn.num_elements(); ++j)
309  {
310  newconn[j] = conn[j] + cnt;
311  }
312 
313  m_f->m_fieldPts->m_ptsConn.push_back(newconn);
314  cnt += newpoints;
315  }
316 
317  m_f->m_fieldPts->m_ptsDim = coordim;
318  if(m_f->m_fielddef.size())
319  {
320  m_f->m_fieldPts->m_nFields = m_f->m_exp.size();
321  }
322  else // just the mesh points
323  {
324  m_f->m_fieldPts->m_nFields = 0;
325  }
326 
327  m_f->m_fieldPts->m_pts = Array<OneD, Array<OneD, NekDouble> >(
328  m_f->m_fieldPts->m_nFields + coordim);
329 
330  for(int i = 0; i < m_f->m_fieldPts->m_nFields + coordim; ++i)
331  {
332  m_f->m_fieldPts->m_pts[i] = Array<OneD, NekDouble>(newtotpoints);
333  }
334 
335  // Interpolate coordinates
336  for(int i = 0; i < coordim; ++i)
337  {
338  coords[i] = Array<OneD, NekDouble>(npts);
339  }
340 
341  for(int i = coordim; i < 3; ++i)
342  {
343  coords[i] = NullNekDouble1DArray;
344  }
345 
346  m_f->m_exp[0]->GetCoords(coords[0],coords[1],coords[2]);
347 
348  int nq1 = m_f->m_exp[0]->GetTotPoints();
349 
350  Array<OneD, NekDouble> x1(nq1);
351  Array<OneD, NekDouble> y1(nq1);
352  Array<OneD, NekDouble> z1(nq1);
353 
354  m_f->m_exp[0]->GetCoords(x1, y1, z1);
355 
356  if (shapedim == 2)
357  {
358  m_f->m_fieldPts->m_ptype = ePtsTriBlock;
359  }
360  else if (shapedim == 3)
361  {
362  m_f->m_fieldPts->m_ptype = ePtsTetBlock;
363  }
364 
365  Array<OneD, NekDouble> tmp;
366 
367  for(int n = 0; n < coordim; ++n)
368  {
369  cnt = 0;
370  int cnt1 = 0;
371  for(int i = 0; i < nel; ++i)
372  {
373  m_f->m_exp[0]->GetExp(i)->PhysInterpToSimplexEquiSpaced(
374  coords[n] + cnt,
375  tmp = m_f->m_fieldPts->m_pts[n] + cnt1);
376  cnt1 += m_f->m_fieldPts->m_npts[i];
377  cnt += m_f->m_exp[0]->GetExp(i)->GetTotPoints();
378  }
379  }
380 
381  if(m_f->m_fielddef.size())
382  {
383  ASSERTL0(m_f->m_fielddef[0]->m_fields.size() == m_f->m_exp.size(),
384  "More expansion defined than fields");
385 
386  for(int n = 0; n < m_f->m_exp.size(); ++n)
387  {
388  cnt = 0;
389  int cnt1 = 0;
390  Array<OneD, const NekDouble> phys = m_f->m_exp[n]->GetPhys();
391  for(int i = 0; i < nel; ++i)
392  {
393  m_f->m_exp[0]->GetExp(i)->PhysInterpToSimplexEquiSpaced(
394  phys + cnt,
395  tmp = m_f->m_fieldPts->m_pts[coordim + n] + cnt1);
396  cnt1 += m_f->m_fieldPts->m_npts[i];
397  cnt += m_f->m_exp[0]->GetExp(i)->GetTotPoints();
398  }
399 
400  // Set up Variable string.
401  m_f->m_fieldPts->m_fields.push_back(
402  m_f->m_fielddef[0]->m_fields[n]);
403  }
404  }
405 }
406 }
407 }
408 
409