44     namespace MultiRegions
 
   47         ExpList2DHomogeneous1D::ExpList2DHomogeneous1D():
 
   58             const bool                                  dealiasing,
 
   65                      "Size of basis number of points and number" 
   66                      "of planes are not the same");
 
   71                     planes.num_elements() * planes[0]->GetExpSize());
 
   73             for(cnt = n = 0; n < planes.num_elements(); ++n)
 
   76                 for (i = 0; i < planes[n]->GetExpSize(); ++i)
 
   78                     (*m_exp)[cnt++] = planes[n]->GetExp(i);
 
   96             const bool dealiasing,
 
  111             for (j = 0; j < nel; ++j)
 
  116             for (n = 1; n < 
m_planes.num_elements(); ++n)
 
  120                 for(j = 0; j < nel; ++j)
 
  122                     (*m_exp).push_back((*
m_exp)[j]);
 
  145             for (
int n = 0; n < 
m_planes.num_elements(); ++n)
 
  164             int ncoeffs_per_plane = 
m_planes[0]->GetNcoeffs();
 
  165             int npoints_per_plane = 
m_planes[0]->GetTotPoints();
 
  167             int nzplanes = 
m_planes.num_elements();
 
  176             int nel = 
m_planes[0]->GetExpSize();
 
  182             for (cnt  = n = 0; n < nzplanes; ++n)
 
  185                     tmparray= 
m_coeffs + ncoeffs_per_plane*n);
 
  187                     tmparray = 
m_phys + npoints_per_plane*n);
 
  189                 for(i = 0; i < nel; ++i)
 
  192                         + n*ncoeffs_per_plane;
 
  193                     m_phys_offset[cnt] =  
m_planes[n]->GetPhys_Offset(i)
 
  194                         + n*npoints_per_plane;
 
  195                     m_offset_elmt_id[cnt++] = 
m_planes[n]->GetOffset_Elmt_Id(i)
 
  209             int nyplanes = 
m_planes.num_elements();
 
  216                     (*m_exp)[eid]->GetCoords(xc0);
 
  222                          "output coord_1 is not defined");
 
  224                     (*m_exp)[eid]->GetCoords(xc0,xc1);
 
  229                 ASSERTL0(
false, 
"Cannot have coordim being three dimensions" 
  230                          "in a homogeneous field");
 
  240             for(n = 0; n < 
m_planes.num_elements(); n++)
 
  248             for (n = 0; n < nyplanes; ++n)
 
  250                 Vmath::Fill(npoints,z[n],tmp_xc = xhom + npoints*n,1);
 
  287             int nyplanes = 
m_planes.num_elements();
 
  288             int npoints = 
m_planes[0]->GetTotPoints();
 
  306             for(n = 0; n < 
m_planes.num_elements(); n++)
 
  314             for(n = 0; n < nyplanes; ++n)
 
  316                 Vmath::Fill(npoints,z[n],tmp_xc = xhom + npoints*n,1);
 
  335             std::ostream &outfile,
 
  340             int nquad0 = (*m_exp)[expansion]->GetNumPoints(0);
 
  341             int nquad1 = 
m_planes.num_elements();
 
  346             coords[1] = coords[0] + nquad0*nquad1;
 
  347             coords[2] = coords[1] + nquad0*nquad1;
 
  349             GetCoords(expansion,coords[0],coords[1],coords[2]);
 
  351             outfile << 
"Zone, I=" << nquad0 << 
", J=" << nquad1
 
  352                     << 
", F=Block" << std::endl;
 
  356                 for(i = 0; i < nquad0*nquad1; ++i)
 
  358                     outfile << coords[j][i] << 
" ";
 
  360                 outfile << std::endl;
 
  366             std::ostream &outfile,
 
  371             int nquad0 = (*m_exp)[expansion]->GetNumPoints(0);
 
  372             int nquad1 = 
m_planes.num_elements();
 
  373             int ntot = nquad0*nquad1;
 
  374             int ntotminus = (nquad0-1)*(nquad1-1);
 
  380             GetCoords(expansion,coords[0],coords[1],coords[2]);
 
  382             outfile << 
"    <Piece NumberOfPoints=\"" 
  383                     << ntot << 
"\" NumberOfCells=\"" 
  384                     << ntotminus << 
"\">" << endl;
 
  385             outfile << 
"      <Points>" << endl;
 
  386             outfile << 
"        <DataArray type=\"Float64\" " 
  387                     << 
"NumberOfComponents=\"3\" format=\"ascii\">" << endl;
 
  389             for (i = 0; i < ntot; ++i)
 
  391                 for (j = 0; j < 3; ++j)
 
  393                     outfile << coords[j][i] << 
" ";
 
  398             outfile << 
"        </DataArray>" << endl;
 
  399             outfile << 
"      </Points>" << endl;
 
  400             outfile << 
"      <Cells>" << endl;
 
  401             outfile << 
"        <DataArray type=\"Int32\" " 
  402                     << 
"Name=\"connectivity\" format=\"ascii\">" << endl;
 
  403             for (i = 0; i < nquad0-1; ++i)
 
  405                 for (j = 0; j < nquad1-1; ++j)
 
  407                     outfile << j*nquad0 + i << 
" " 
  408                             << j*nquad0 + i + 1 << 
" " 
  409                             << (j+1)*nquad0 + i + 1 << 
" " 
  410                             << (j+1)*nquad0 + i << endl;
 
  414             outfile << 
"        </DataArray>" << endl;
 
  415             outfile << 
"        <DataArray type=\"Int32\" " 
  416                     << 
"Name=\"offsets\" format=\"ascii\">" << endl;
 
  417             for (i = 0; i < ntotminus; ++i)
 
  419                 outfile << i*4+4 << 
" ";
 
  422             outfile << 
"        </DataArray>" << endl;
 
  423             outfile << 
"        <DataArray type=\"UInt8\" " 
  424                     << 
"Name=\"types\" format=\"ascii\">" << endl;
 
  425             for (i = 0; i < ntotminus; ++i)
 
  430             outfile << 
"        </DataArray>" << endl;
 
  431             outfile << 
"      </Cells>" << endl;
 
  432             outfile << 
"      <PointData>" << endl;
 
  438             int nPlanes   = 
m_planes.num_elements();
 
  439             int nPtsPlane = 
m_planes[0]->GetNpoints();
 
  442             ASSERTL1(normals.num_elements() >= nDim,
 
  443                      "Output vector does not have sufficient dimensions to" 
  445             ASSERTL1(normals[0].num_elements() >= nPtsPlane,
 
  446                      "Output vector does not have sufficient dimensions to" 
  453             for (
int i = 0; i < nDim; ++i)
 
  455                 for (
int n = 1; n < nPlanes; ++n)
 
  458                                  &normals[i][n*nPtsPlane], 1);
 
Abstraction of a two-dimensional multi-elemental expansion which is merely a collection of local expa...
#define ASSERTL0(condition, msg)
virtual ~ExpList2DHomogeneous1D()
Destructor. 
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool. 
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
LibUtilities::TranspositionSharedPtr m_transposition
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value. 
const boost::shared_ptr< LocalRegions::ExpansionVector > GetExp() const 
This function returns the vector of elements in the expansion. 
virtual void v_WriteTecplotZone(std::ostream &outfile, int expansion)
NekDouble m_lhom
Width of homogeneous direction. 
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points. 
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients. 
int GetExpSize(void)
This function returns the number of elements in the expansion. 
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Abstraction of a two-dimensional multi-elemental expansion which is merely a collection of local expa...
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs. 
LibUtilities::BasisSharedPtr m_homogeneousBasis
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y. 
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions. 
int GetTotPoints(void) const 
Returns the total number of quadrature points m_npoints . 
boost::shared_ptr< ExpList1D > ExpList1DSharedPtr
Shared pointer to an ExpList1D object. 
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys. 
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs . 
Array< OneD, int > m_offset_elmt_id
Array containing the element id m_offset_elmt_id[n] that the n^th consecutive block of data in m_coef...
LibUtilities::SessionReaderSharedPtr m_session
Session. 
virtual void v_GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1, Array< OneD, NekDouble > &coord_2)
Array< OneD, ExpListSharedPtr > m_planes
virtual void v_GetNormals(Array< OneD, Array< OneD, NekDouble > > &normals)
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x. 
void SetCoeffPhys(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
ExpList2DHomogeneous1D()
Default constructor. 
This class is the abstraction of a one-dimensional multi-elemental expansions which is merely a colle...
int GetCoordim(int eid)
This function returns the dimension of the coordinates of the element eid. 
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 . ...
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
boost::shared_ptr< MeshGraph > MeshGraphSharedPtr
virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip)
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Describes the specification for a Basis.