Nektar++
ExpList2DHomogeneous1D.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File ExpList2DHomogeneous1D.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: An ExpList2D which is homogeneous in 1 direction and so
32 // uses much of the functionality from a ExpList and its daughters
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
37 #include <MultiRegions/ExpList.h>
38 
39 using namespace std;
40 
41 namespace Nektar
42 {
43  namespace MultiRegions
44  {
45  // Forward declaration for typedefs
46  ExpList2DHomogeneous1D::ExpList2DHomogeneous1D():
48  {
49  }
50 
51  // Constructor for ExpList2DHomogeneous1D to act as a Explist2D field
54  const LibUtilities::BasisKey &HomoBasis,
55  const NekDouble lhom,
56  const bool useFFT,
57  const bool dealiasing,
58  const Array<OneD, ExpListSharedPtr> &planes,
59  const LibUtilities::CommSharedPtr comm)
60  : ExpListHomogeneous1D(e2DH1D,pSession,HomoBasis,
61  lhom,useFFT,dealiasing)
62  {
63  int i, n, cnt;
64 
65  if (comm)
66  {
67  m_comm = comm;
68  }
69 
70  ASSERTL1(m_planes.size() == planes.size(),
71  "Size of basis number of points and number"
72  "of planes are not the same");
73 
74  // Set up expansion list with elements from all planes.
77  planes.size() * planes[0]->GetExpSize());
78 
79  for(cnt = n = 0; n < planes.size(); ++n)
80  {
81  m_planes[n] = planes[n];
82  for (i = 0; i < planes[n]->GetExpSize(); ++i)
83  {
84  (*m_exp)[cnt++] = planes[n]->GetExp(i);
85  }
86  }
87 
88  SetCoeffPhys();
89  }
90 
91  // Constructor for ExpList2DHomogeneous1D to act as a Explist2D field
94  const LibUtilities::BasisKey &HomoBasis,
95  const NekDouble lhom,
96  const bool useFFT,
97  const bool dealiasing,
99  const Collections::ImplementationType ImpType):
100  ExpListHomogeneous1D(e2DH1D, pSession,HomoBasis,
101  lhom,useFFT,dealiasing)
102  {
103  int n, j, nel;
104  ExpListSharedPtr plane_zero;
105 
106  // note that nzplanes can be larger than nzmodes
107  m_planes[0] = plane_zero = MemoryManager<ExpList>::
108  AllocateSharedPtr(m_session, graph1D, false, "DefaultVar", ImpType);
109 
111 
112  nel = m_planes[0]->GetExpSize();
113 
114  for (j = 0; j < nel; ++j)
115  {
116  (*m_exp).push_back(m_planes[0]->GetExp(j));
117  }
118 
119  for (n = 1; n < m_planes.size(); ++n)
120  {
122  AllocateSharedPtr(*plane_zero, false);
123  for(j = 0; j < nel; ++j)
124  {
125  (*m_exp).push_back((*m_exp)[j]);
126  }
127  }
128 
129  SetCoeffPhys();
130  }
131 
132 
133  /**
134  * @param In ExpList2DHomogeneous1D object to copy.
135  */
137  const ExpList2DHomogeneous1D &In):
139  {
140  ExpListSharedPtr zero_plane =
141  std::dynamic_pointer_cast<ExpList> (In.m_planes[0]);
142 
143  for (int n = 0; n < m_planes.size(); ++n)
144  {
146  AllocateSharedPtr(*zero_plane, false);
147  }
148 
149  SetCoeffPhys();
150  }
151 
152  /**
153  * Destructor
154  */
156  {
157  }
158 
160  {
161  int i, n, cnt;
162  int ncoeffs_per_plane = m_planes[0]->GetNcoeffs();
163  int npoints_per_plane = m_planes[0]->GetTotPoints();
164 
165  int nzplanes = m_planes.size();
166 
167  // Set total coefficients and points
168  m_ncoeffs = ncoeffs_per_plane*nzplanes;
169  m_npoints = npoints_per_plane*nzplanes;
170 
173 
174  int nel = m_planes[0]->GetExpSize();
175  m_coeff_offset = Array<OneD,int>(nel*nzplanes);
176  m_phys_offset = Array<OneD,int>(nel*nzplanes);
177  Array<OneD, NekDouble> tmparray;
178 
179  for (cnt = n = 0; n < nzplanes; ++n)
180  {
181  m_planes[n]->SetCoeffsArray(
182  tmparray= m_coeffs + ncoeffs_per_plane*n);
183  m_planes[n]->SetPhysArray(
184  tmparray = m_phys + npoints_per_plane*n);
185 
186  for(i = 0; i < nel; ++i)
187  {
188  m_coeff_offset[cnt] = m_planes[n]->GetCoeff_Offset(i)
189  + n*ncoeffs_per_plane;
190  m_phys_offset[cnt++] = m_planes[n]->GetPhys_Offset(i)
191  + n*npoints_per_plane;
192  }
193  }
194  }
195 
197  const int eid,
201  {
202  int n, coordim;
203  Array<OneD, NekDouble> tmp_xc,xhom;
204  int nyplanes = m_planes.size();
205  int npoints = GetTotPoints(eid);
206 
207  switch(coordim = GetCoordim(0))
208  {
209  case 1:
210  {
211  (*m_exp)[eid]->GetCoords(xc0);
212  xhom = xc1;
213  }
214  break;
215  case 2:
216  ASSERTL0(xc1.size() != 0,
217  "output coord_1 is not defined");
218  {
219  (*m_exp)[eid]->GetCoords(xc0,xc1);
220  xhom = xc2;
221  }
222  break;
223  default:
224  ASSERTL0(false, "Cannot have coordim being three dimensions"
225  "in a homogeneous field");
226  break;
227  }
228 
229  // Fill homogeneous-direction
231  Array<OneD, NekDouble> z(nyplanes);
232 
233  Array<OneD, NekDouble> local_pts(m_planes.size());
234 
235  for(n = 0; n < m_planes.size(); n++)
236  {
237  local_pts[n] = pts[m_transposition->GetPlaneID(n)];
238  }
239 
240  Vmath::Smul(nyplanes,m_lhom/2.0,local_pts,1,z,1);
241  Vmath::Sadd(nyplanes,m_lhom/2.0,z,1,z,1);
242 
243  for (n = 0; n < nyplanes; ++n)
244  {
245  Vmath::Fill(npoints,z[n],tmp_xc = xhom + npoints*n,1);
246  if(n)
247  {
248  Vmath::Vcopy(npoints,xc0,1,tmp_xc = xc0+npoints*n,1);
249  if(coordim == 2)
250  {
251  Vmath::Vcopy(npoints,xc1,1,tmp_xc = xc1+npoints*n,1);
252  }
253  }
254  }
255  }
256 
257  /**
258  * The operation calls the 2D plane coordinates through the
259  * function ExpList#GetCoords and then evaluated the third
260  * coordinate using the member \a m_lhom
261  *
262  * @param coord_0 After calculation, the \f$x_1\f$ coordinate
263  * will be stored in this array.
264  *
265  * @param coord_1 After calculation, the \f$x_2\f$ coordinate
266  * will be stored in this array. This
267  * coordinate might be evaluated using the
268  * predefined value \a m_lhom
269  *
270  * @param coord_2 After calculation, the \f$x_3\f$ coordinate
271  * will be stored in this array. This
272  * coordinate is evaluated using the
273  * predefined value \a m_lhom
274  */
279  {
280  int n,coordim;
281  Array<OneD, NekDouble> tmp_xc, xhom;
282  int nyplanes = m_planes.size();
283  int npoints = m_planes[0]->GetTotPoints();
284 
285  m_planes[0]->GetCoords(xc0,xc1);
286 
287  if((coordim = GetCoordim(0)) == 1)
288  {
289  xhom = xc1;
290  }
291  else
292  {
293  xhom = xc2;
294  }
295 
296  // Fill z-direction
298  Array<OneD, NekDouble> z(nyplanes);
299  Array<OneD, NekDouble> local_pts(m_planes.size());
300 
301  for(n = 0; n < m_planes.size(); n++)
302  {
303  local_pts[n] = pts[m_transposition->GetPlaneID(n)];
304  }
305 
306  Vmath::Smul(nyplanes,m_lhom/2.0,local_pts,1,z,1);
307  Vmath::Sadd(nyplanes,m_lhom/2.0,z,1,z,1);
308 
309  for(n = 0; n < nyplanes; ++n)
310  {
311  Vmath::Fill(npoints,z[n],tmp_xc = xhom + npoints*n,1);
312  if(n)
313  {
314  Vmath::Vcopy(npoints,xc0,1,tmp_xc = xc0+npoints*n,1);
315  if(coordim == 2)
316  {
317  Vmath::Vcopy(npoints,xc1,1,tmp_xc = xc1+npoints*n,1);
318  }
319  }
320  }
321  }
322 
323 
324  /**
325  * Write Tecplot Files Zone
326  * @param outfile Output file name.
327  * @param expansion Expansion that is considered
328  */
330  std::ostream &outfile,
331  int expansion)
332  {
333  int i, j;
334 
335  int nquad0 = (*m_exp)[expansion]->GetNumPoints(0);
336  int nquad1 = m_planes.size();
337 
338  Array<OneD,NekDouble> coords[3];
339 
340  coords[0] = Array<OneD,NekDouble>(3*nquad0*nquad1);
341  coords[1] = coords[0] + nquad0*nquad1;
342  coords[2] = coords[1] + nquad0*nquad1;
343 
344  GetCoords(expansion,coords[0],coords[1],coords[2]);
345 
346  outfile << "Zone, I=" << nquad0 << ", J=" << nquad1
347  << ", F=Block" << std::endl;
348 
349  for(j = 0; j < GetCoordim(0)+1; ++j)
350  {
351  for(i = 0; i < nquad0*nquad1; ++i)
352  {
353  outfile << coords[j][i] << " ";
354  }
355  outfile << std::endl;
356  }
357  }
358 
359 
361  std::ostream &outfile,
362  int expansion,
363  int istrip)
364  {
365  boost::ignore_unused(istrip);
366 
367  // If there is only one plane (e.g. HalfMode), we write a 2D plane.
368  if (m_planes.size() == 1)
369  {
370  m_planes[0]->WriteVtkPieceHeader(outfile, expansion);
371  return;
372  }
373 
374  // If we are using Fourier points, output extra plane to fill domain
375  int outputExtraPlane = 0;
376  if ( m_homogeneousBasis->GetBasisType() == LibUtilities::eFourier
377  && m_homogeneousBasis->GetPointsType() ==
379  {
380  outputExtraPlane = 1;
381  }
382  int i, j;
383  int nquad0 = (*m_exp)[expansion]->GetNumPoints(0);
384  int nquad1 = m_planes.size() + outputExtraPlane;
385  int ntot = nquad0*nquad1;
386  int ntotminus = (nquad0-1)*(nquad1-1);
387 
388  Array<OneD,NekDouble> coords[3];
389  coords[0] = Array<OneD,NekDouble>(ntot);
390  coords[1] = Array<OneD,NekDouble>(ntot);
391  coords[2] = Array<OneD,NekDouble>(ntot);
392  GetCoords(expansion,coords[0],coords[1],coords[2]);
393 
394  if (outputExtraPlane)
395  {
396  // Copy coords[0] and coords[1] to extra plane
398  Vmath::Vcopy (nquad0, coords[0], 1,
399  tmp = coords[0] + (nquad1-1)*nquad0, 1);
400  Vmath::Vcopy (nquad0, coords[1], 1,
401  tmp = coords[1] + (nquad1-1)*nquad0, 1);
402  // Fill coords[2] for extra plane
403  NekDouble z = coords[2][nquad0*m_planes.size()-1] +
404  (coords[2][nquad0] - coords[2][0]);
405  Vmath::Fill(nquad0, z, tmp = coords[2] + (nquad1-1)*nquad0, 1);
406  }
407 
408  outfile << " <Piece NumberOfPoints=\""
409  << ntot << "\" NumberOfCells=\""
410  << ntotminus << "\">" << endl;
411  outfile << " <Points>" << endl;
412  outfile << " <DataArray type=\"Float64\" "
413  << "NumberOfComponents=\"3\" format=\"ascii\">" << endl;
414  outfile << " ";
415  for (i = 0; i < ntot; ++i)
416  {
417  for (j = 0; j < 3; ++j)
418  {
419  outfile << coords[j][i] << " ";
420  }
421  outfile << endl;
422  }
423  outfile << endl;
424  outfile << " </DataArray>" << endl;
425  outfile << " </Points>" << endl;
426  outfile << " <Cells>" << endl;
427  outfile << " <DataArray type=\"Int32\" "
428  << "Name=\"connectivity\" format=\"ascii\">" << endl;
429  for (i = 0; i < nquad0-1; ++i)
430  {
431  for (j = 0; j < nquad1-1; ++j)
432  {
433  outfile << j*nquad0 + i << " "
434  << j*nquad0 + i + 1 << " "
435  << (j+1)*nquad0 + i + 1 << " "
436  << (j+1)*nquad0 + i << endl;
437  }
438  }
439  outfile << endl;
440  outfile << " </DataArray>" << endl;
441  outfile << " <DataArray type=\"Int32\" "
442  << "Name=\"offsets\" format=\"ascii\">" << endl;
443  for (i = 0; i < ntotminus; ++i)
444  {
445  outfile << i*4+4 << " ";
446  }
447  outfile << endl;
448  outfile << " </DataArray>" << endl;
449  outfile << " <DataArray type=\"UInt8\" "
450  << "Name=\"types\" format=\"ascii\">" << endl;
451  for (i = 0; i < ntotminus; ++i)
452  {
453  outfile << "9 ";
454  }
455  outfile << endl;
456  outfile << " </DataArray>" << endl;
457  outfile << " </Cells>" << endl;
458  outfile << " <PointData>" << endl;
459  }
460 
462  Array<OneD, Array<OneD, NekDouble> > &normals)
463  {
464  int nPlanes = m_planes.size();
465  int nPtsPlane = m_planes[0]->GetNpoints();
466  int nDim = GetCoordim(0) + 1;
467 
468  ASSERTL1(normals.size() >= nDim,
469  "Output vector does not have sufficient dimensions to"
470  "match coordim");
471  ASSERTL1(normals[0].size() >= nPtsPlane,
472  "Output vector does not have sufficient dimensions to"
473  "match coordim");
474 
475  // Calculate normals from plane 0.
476  m_planes[0]->GetNormals(normals);
477 
478  // Copy remaining planes.
479  for (int i = 0; i < nDim; ++i)
480  {
481  for (int n = 1; n < nPlanes; ++n)
482  {
483  Vmath::Vcopy(nPtsPlane, &normals[i][0], 1,
484  &normals[i][n*nPtsPlane], 1);
485  }
486  }
487  }
488  } //end of namespace
489 } //end of namespace
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
Describes the specification for a Basis.
Definition: Basis.h:50
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
Abstraction of a two-dimensional multi-elemental expansion which is merely a collection of local expa...
void SetCoeffPhys(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
virtual void v_GetNormals(Array< OneD, Array< OneD, NekDouble > > &normals)
Populate normals with the normals of all expansions.
virtual void v_WriteTecplotZone(std::ostream &outfile, int expansion)
virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip)
void GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1=NullNekDouble1DArray, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
This function calculates the coordinates of all the elemental quadrature points .
virtual void v_GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1, Array< OneD, NekDouble > &coord_2)
Abstraction of a two-dimensional multi-elemental expansion which is merely a collection of local expa...
NekDouble m_lhom
Width of homogeneous direction.
LibUtilities::TranspositionSharedPtr m_transposition
Array< OneD, ExpListSharedPtr > m_planes
LibUtilities::BasisSharedPtr m_homogeneousBasis
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:1252
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1304
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: ExpList.h:1220
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1290
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:1230
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:2422
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:1223
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:1307
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:1269
int GetTotPoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1850
int GetCoordim(int eid)
This function returns the dimension of the coordinates of the element eid.
Definition: ExpList.h:2244
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:65
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:54
@ eFourier
Fourier Expansion .
Definition: BasisType.h:53
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
double NekDouble
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:225
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha - x.
Definition: Vmath.cpp:341
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199