44     namespace MultiRegions
 
   47         ExpList3DHomogeneous1D::ExpList3DHomogeneous1D():
 
   58             const bool dealiasing):
 
   70             const bool dealiasing,
 
   72             const std::string &var):
 
   86             const bool dealiasing,
 
  107             for(j = 0; j < nel; ++j)
 
  112             for(n = 1; n < 
m_planes.num_elements(); ++n)
 
  115                 for(j = 0; j < nel; ++j)
 
  117                     (*m_exp).push_back((*
m_exp)[j]);
 
  138             if(DeclarePlanesSetCoeffPhys)
 
  143                 for(
int n = 0; n < 
m_planes.num_elements(); ++n)
 
  157                 const std::vector<unsigned int> &eIDs,
 
  158                 bool DeclarePlanesSetCoeffPhys):
 
  163             if(DeclarePlanesSetCoeffPhys)
 
  166                 std::vector<unsigned int> eIDsPlane;
 
  167                 int nel = eIDs.size()/
m_planes.num_elements();
 
  169                 for (
int i = 0; i < nel; ++i)
 
  171                     eIDsPlane.push_back(eIDs[i]);
 
  180                 for(
int n = 0; n < 
m_planes.num_elements(); ++n)
 
  199             int ncoeffs_per_plane = 
m_planes[0]->GetNcoeffs();
 
  200             int npoints_per_plane = 
m_planes[0]->GetTotPoints();
 
  202             int nzplanes = 
m_planes.num_elements();
 
  211             int nel = 
m_planes[0]->GetExpSize();
 
  217             for(cnt  = n = 0; n < nzplanes; ++n)
 
  220                 m_planes[n]->SetPhysArray(tmparray = 
m_phys + npoints_per_plane*n);
 
  222                 for(i = 0; i < nel; ++i)
 
  225                     m_phys_offset[cnt] =  
m_planes[n]->GetPhys_Offset(i) + n*npoints_per_plane;
 
  226                     m_offset_elmt_id[cnt++] = 
m_planes[n]->GetOffset_Elmt_Id(i) + n*nel;
 
  238             int nzplanes = 
m_planes.num_elements();
 
  241             (*m_exp)[eid]->GetCoords(xc0,xc1);
 
  247             for(n = 0; n < 
m_planes.num_elements(); n++)
 
  257             for(n = 0; n < nzplanes; ++n)
 
  259                 Vmath::Fill(npoints,z[n],tmp_xc = xc2 + npoints*n,1);
 
  290             int nzplanes = 
m_planes.num_elements();
 
  291             int npoints = 
m_planes[0]->GetTotPoints();
 
  300             for(n = 0; n < 
m_planes.num_elements(); n++)
 
  310             for(n = 0; n < nzplanes; ++n)
 
  312                 Vmath::Fill(npoints,z[n],tmp_xc = xc2 + npoints*n,1);
 
  323             ASSERTL0(expansion == -1, 
"Multi-zone output not supported for homogeneous expansions.");
 
  325             const int nPtsPlane = 
m_planes[0]->GetNpoints();
 
  326             const int nElmt     = 
m_planes[0]->GetExpSize();
 
  327             const int nPlanes   = 
m_planes.num_elements();
 
  331             for (
int i = 0; i < nElmt; ++i)
 
  333                 const int np0 = (*m_exp)[i]->GetNumPoints(0);
 
  334                 const int np1 = (*m_exp)[i]->GetNumPoints(1);
 
  336                 for (
int n = 1; n < nPlanes; ++n)
 
  338                     const int o1 = (n-1) * nPtsPlane;
 
  339                     const int o2 =  n    * nPtsPlane;
 
  340                     for (
int j = 1; j < np1; ++j)
 
  342                         for(
int k = 1; k < np0; ++k)
 
  344                             outfile << cnt + (j-1)*np0 + (k-1) + o1 + 1 << 
" ";
 
  345                             outfile << cnt + (j-1)*np0 + (k-1) + o2 + 1 << 
" ";
 
  346                             outfile << cnt + (j-1)*np0 +  k    + o2 + 1 << 
" ";
 
  347                             outfile << cnt + (j-1)*np0 +  k    + o1 + 1 << 
" ";
 
  348                             outfile << cnt +  j   *np0 + (k-1) + o1 + 1 << 
" ";
 
  349                             outfile << cnt +  j   *np0 + (k-1) + o2 + 1 << 
" ";
 
  350                             outfile << cnt +  j   *np0 +  k    + o2 + 1 << 
" ";
 
  351                             outfile << cnt +  j   *np0 +  k    + o1 + 1 << endl;
 
  363                 std::ostream    &outfile,
 
  370                 m_planes[0]->WriteVtkPieceHeader(outfile, expansion);
 
  375             int outputExtraPlane = 0;
 
  380                 outputExtraPlane = 1;
 
  384             int nq0 = (*m_exp)[expansion]->GetNumPoints(0);
 
  385             int nq1 = (*m_exp)[expansion]->GetNumPoints(1);
 
  386             int nq2 = 
m_planes.num_elements() + outputExtraPlane;
 
  387             int ntot = nq0*nq1*nq2;
 
  388             int ntotminus = (nq0-1)*(nq1-1)*(nq2-1);
 
  394             GetCoords(expansion,coords[0],coords[1],coords[2]);
 
  396             if (outputExtraPlane)
 
  401                                       tmp = coords[0] + (nq2-1)*nq0*nq1, 1);
 
  403                                       tmp = coords[1] + (nq2-1)*nq0*nq1, 1);
 
  406                               (coords[2][nq0*nq1] - coords[2][0]);
 
  407                 Vmath::Fill(nq0*nq1, z, tmp = coords[2] + (nq2-1)*nq0*nq1, 1);
 
  411             m_session->LoadParameter(
"DistStrip", DistStrip, 0);
 
  413             for(
int i = 0; i < ntot; i++)
 
  415                 coords[2][i] += istrip*DistStrip;
 
  418             outfile << 
"    <Piece NumberOfPoints=\"" 
  419                     << ntot << 
"\" NumberOfCells=\"" 
  420                     << ntotminus << 
"\">" << endl;
 
  421             outfile << 
"      <Points>" << endl;
 
  422             outfile << 
"        <DataArray type=\"Float64\" " 
  423                     << 
"NumberOfComponents=\"3\" format=\"ascii\">" << endl;
 
  425             for (i = 0; i < ntot; ++i)
 
  427                 for (j = 0; j < 3; ++j)
 
  429                     outfile << coords[j][i] << 
" ";
 
  434             outfile << 
"        </DataArray>" << endl;
 
  435             outfile << 
"      </Points>" << endl;
 
  436             outfile << 
"      <Cells>" << endl;
 
  437             outfile << 
"        <DataArray type=\"Int32\" " 
  438                     << 
"Name=\"connectivity\" format=\"ascii\">" << endl;
 
  439             for (i = 0; i < nq0-1; ++i)
 
  441                 for (j = 0; j < nq1-1; ++j)
 
  443                     for (k = 0; k < nq2-1; ++k)
 
  445                         outfile << k*nq0*nq1     + j*nq0     + i     << 
" " 
  446                                 << k*nq0*nq1     + j*nq0     + i + 1 << 
" " 
  447                                 << k*nq0*nq1     + (j+1)*nq0 + i + 1 << 
" " 
  448                                 << k*nq0*nq1     + (j+1)*nq0 + i     << 
" " 
  449                                 << (k+1)*nq0*nq1 + j*nq0     + i     << 
" " 
  450                                 << (k+1)*nq0*nq1 + j*nq0     + i + 1 << 
" " 
  451                                 << (k+1)*nq0*nq1 + (j+1)*nq0 + i + 1 << 
" " 
  452                                 << (k+1)*nq0*nq1 + (j+1)*nq0 + i     << endl;
 
  457             outfile << 
"        </DataArray>" << endl;
 
  458             outfile << 
"        <DataArray type=\"Int32\" " 
  459                     << 
"Name=\"offsets\" format=\"ascii\">" << endl;
 
  460             for (i = 0; i < ntotminus; ++i)
 
  462                 outfile << i*8+8 << 
" ";
 
  465             outfile << 
"        </DataArray>" << endl;
 
  466             outfile << 
"        <DataArray type=\"UInt8\" " 
  467                     << 
"Name=\"types\" format=\"ascii\">" << endl;
 
  468             for (i = 0; i < ntotminus; ++i)
 
  473             outfile << 
"        </DataArray>" << endl;
 
  474             outfile << 
"      </Cells>" << endl;
 
  475             outfile << 
"      <PointData>" << endl;
 
  487             for(
int n = 0; n < 
m_planes.num_elements(); n++)
 
  494                 for(
int n = 0; n < 
m_planes.num_elements(); ++n)
 
  496                     errL2 = 
m_planes[n]->L2(inarray + cnt);
 
  498                     err  += errL2*errL2*local_w[n]*
m_lhom*0.5;
 
  503                 for(
int n = 0; n < 
m_planes.num_elements(); ++n)
 
  505                     errL2 = 
m_planes[n]->L2(inarray + cnt, soln + cnt);
 
  507                     err  += errL2*errL2*local_w[n]*
m_lhom*0.5;
 
  521             for (
int n = 0; n < 
m_planes[0]->GetExpSize(); ++n)
 
  524                 area += 
m_planes[0]->GetExp(n)->Integral(inarray);
 
  530             for (
int n = 0; n < 
m_planes.num_elements(); n += 2)
 
  536                 for(
int i = 0; i < 
m_planes[n]->GetExpSize(); ++i)
 
  543                     energy[n/2] += err*err;
 
  549                     energy[n/2] += err*err;
 
  553                 energy[n/2] /= 2.0*area;
 
Abstraction of a two-dimensional multi-elemental expansion which is merely a collection of local expa...
 
const Array< OneD, const NekDouble > & GetCoeffs() const 
This function returns (a reference to) the array  (implemented as m_coeffs) containing all local expa...
 
virtual void v_WriteTecplotConnectivity(std::ostream &outfile, int expansion)
 
#define ASSERTL0(condition, msg)
 
int GetCoeff_Offset(int n) const 
Get the start offset position for a global list of m_coeffs correspoinding to element n...
 
virtual Array< OneD, const NekDouble > v_HomogeneousEnergy(void)
 
static Array< OneD, NekDouble > NullNekDouble1DArray
 
virtual ~ExpList3DHomogeneous1D()
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. 
 
void GenExpList3DHomogeneous1D(const SpatialDomains::ExpansionMap &expansions)
 
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
 
virtual void v_GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1, Array< OneD, NekDouble > &coord_2)
 
void SetCoeffPhys(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
 
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...
 
1D Evenly-spaced points using Fourier Fit 
 
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y. 
 
ExpList3DHomogeneous1D()
Default constructor. 
 
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 . 
 
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. 
 
Array< OneD, ExpListSharedPtr > m_planes
 
boost::shared_ptr< ExpList2D > ExpList2DSharedPtr
Shared pointer to an ExpList2D object. 
 
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x. 
 
virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip)
 
LibUtilities::CommSharedPtr m_comm
Communicator. 
 
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 . ...
 
Abstraction of a two-dimensional multi-elemental expansion which is merely a collection of local expa...
 
virtual NekDouble v_L2(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
 
void SetExpType(ExpansionType Type)
Returns the type of the expansion. 
 
boost::shared_ptr< StdExpansion > StdExpansionSharedPtr
 
boost::shared_ptr< MeshGraph > MeshGraphSharedPtr
 
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
 
Abstraction of a two-dimensional multi-elemental expansion which is merely a collection of local expa...
 
Describes the specification for a Basis. 
 
std::map< int, ExpansionShPtr > ExpansionMap