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