4#include <thrust/host_vector.h>
6#include "../backend/config.h"
53template<
class real_type,
size_t Nd>
93template<
class real_type,
size_t Nd>
106 constexpr static unsigned ndim() {
return Nd;}
118 <<u<<
" while Nd is "<<Nd);
119 return m_n[u]*m_N[u];
135 <<u<<
" while Nd is "<<Nd);
139 for(
unsigned i=0; i<m_N[u]; i++)
140 for(
unsigned j=0; j<m_n[u]; j++)
142 real_type xmiddle = DG_FMA( hu, (real_type)(i), m_x0[u]);
143 real_type h2 = hu/2.;
144 real_type absj = 1.+aa[j];
145 abs[i*m_n[u]+j] = DG_FMA( h2, absj, xmiddle);
158 throw Error(
Message(
_ping_)<<
"u>Nd not allowed! You typed: "<<u<<
" while Nd is "<<Nd);
162 for(
unsigned i=0; i<m_N[u]; i++)
163 for(
unsigned j=0; j<m_n[u]; j++)
164 v[i*m_n[u] + j] = hu/2.*ww[j];
174 std::array<unsigned,Nd> ss;
175 for(
unsigned u=0; u<Nd; u++)
184 std::array<host_vector,Nd> abs;
185 for(
unsigned u=0; u<Nd; u++)
194 std::array<host_vector,Nd> w;
195 for(
unsigned u=0; u<Nd; u++)
202 std::array<real_type,Nd>
get_p()
const{
207 std::array<real_type,Nd>
get_q()
const{
212 std::array<real_type,Nd>
get_l()
const{
213 std::array<real_type, Nd>
p;
214 for(
unsigned u=0; u<Nd; u++)
220 std::array<real_type,Nd>
get_h()
const{
221 std::array<real_type, Nd> hh;
222 for(
unsigned u=0; u<Nd; u++)
228 std::array<unsigned, Nd>
get_N()
const
234 std::array<unsigned, Nd>
get_n()
const
250 real_type
p(
unsigned u=0)
const {
return m_x0.at(u);}
253 real_type
q(
unsigned u=0)
const {
return m_x1.at(u);}
256 real_type
h(
unsigned u=0)
const {
return (m_x1.at(u) - m_x0.at(u))/(real_type)m_N.at(u);}
259 real_type
l(
unsigned u=0)
const {
return m_x1.at(u) - m_x0.at(u);}
262 unsigned n(
unsigned u=0)
const {
return m_n.at(u);}
265 unsigned N(
unsigned u=0)
const {
return m_N.at(u);}
268 dg::bc bc(
unsigned u=0)
const {
return m_bcs.at(u);}
278 throw Error(
Message(
_ping_)<<
"u>Nd not allowed! You typed: "<<u<<
" while Nd is "<<Nd);
284 template<
size_t Md = Nd>
285 real_type
x0()
const {
return std::get<0>(m_x0);}
287 template<
size_t Md = Nd>
288 real_type
x1()
const {
return std::get<0>(m_x1);}
290 template<
size_t Md = Nd>
291 real_type
y0()
const {
return std::get<1>(m_x0);}
293 template<
size_t Md = Nd>
294 real_type
y1()
const {
return std::get<1>(m_x1);}
296 template<
size_t Md = Nd>
297 real_type
z0()
const {
return std::get<2>(m_x0);}
299 template<
size_t Md = Nd>
300 real_type
z1()
const {
return std::get<2>(m_x1);}
303 template<
size_t Md = Nd>
304 real_type
lx()
const {
return std::get<0>(
get_l());}
306 template<
size_t Md = Nd>
307 real_type
ly()
const {
return std::get<1>(
get_l());}
309 template<
size_t Md = Nd>
310 real_type
lz()
const {
return std::get<2>(
get_l());}
313 template<
size_t Md = Nd>
314 real_type
hx()
const {
return std::get<0>(
get_h());}
316 template<
size_t Md = Nd>
317 real_type
hy()
const {
return std::get<1>(
get_h());}
319 template<
size_t Md = Nd>
320 real_type
hz()
const {
return std::get<2>(
get_h());}
323 template<
size_t Md = Nd>
324 unsigned nx()
const {
return std::get<0>(m_n);}
326 template<
size_t Md = Nd>
327 unsigned ny()
const {
return std::get<1>(m_n);}
329 template<
size_t Md = Nd>
330 unsigned nz()
const {
return std::get<2>(m_n);}
333 template<
size_t Md = Nd>
334 unsigned Nx()
const {
return std::get<0>(m_N);}
336 template<
size_t Md = Nd>
337 unsigned Ny()
const {
return std::get<1>(m_N);}
339 template<
size_t Md = Nd>
340 unsigned Nz()
const {
return std::get<2>(m_N);}
343 template<
size_t Md = Nd>
346 template<
size_t Md = Nd>
349 template<
size_t Md = Nd>
353 template<
size_t Md = Nd>
355 static_assert( Nd > 0);
359 template<
size_t Md = Nd>
361 static_assert( Nd > 1);
365 template<
size_t Md = Nd>
367 static_assert( Nd > 2);
376 std::array<unsigned, Nd>
start()
const {
return {0};}
385 std::array<unsigned, Nd>
count()
const
387 std::array<unsigned,Nd> ss;
388 for(
unsigned u=0; u<Nd; u++)
389 ss[Nd-1-u] = m_n[u]*m_N[u];
403 template<
size_t Md = Nd>
406 Ns[0] = round(fx*(real_type)m_N[0]);
407 Ns[1] = round(fy*(real_type)m_N[1]);
408 if( fx != 1 or fy != 1)
417 template<
size_t Md = Nd>
418 std::enable_if_t<(Md == 1),
void>
set(
unsigned new_n,
unsigned new_Nx)
420 set(std::array{new_n}, std::array{new_Nx});
430 template<
size_t Md = Nd>
431 std::enable_if_t<(Md == 2),
void>
set(
unsigned new_n,
unsigned new_Nx,
434 set({new_n,new_n}, {new_Nx,new_Ny});
446 template<
size_t Md = Nd>
447 std::enable_if_t<(Md == 3),
void>
set(
unsigned new_n,
unsigned new_Nx,
448 unsigned new_Ny,
unsigned new_Nz)
450 set({new_n,new_n,1}, {new_Nx,new_Ny,new_Nz});
453 void set(
unsigned new_n, std::array<unsigned,Nd> new_N)
455 std::array<unsigned , Nd> tmp;
456 for(
unsigned u=0; u<Nd; u++)
467 void set_axis(
unsigned coord,
unsigned new_n ,
unsigned new_N)
469 std::array<unsigned,Nd>
n = m_n,
N = m_N;
481 void set( std::array<unsigned,Nd> new_n, std::array<unsigned,Nd> new_N)
483 if( new_n==m_n && new_N == m_N)
494 void set_pq( std::array<real_type,Nd> new_p, std::array<real_type,Nd>
518 void set( std::array<real_type,Nd> new_p, std::array<real_type,Nd> new_q,
519 std::array<unsigned,Nd> new_n, std::array<unsigned,Nd> new_N,
520 std::array<dg::bc,Nd> new_bcs)
534 for(
unsigned u=0; u<Nd; u++)
535 size *= m_n[u]*m_N[u];
545 void display( std::ostream& os = std::cout)
const
547 for(
unsigned u=0; u<Nd; u++)
549 os <<
"Topology parameters for Grid "<<u<<
" are: \n"
550 <<
" n = "<<m_n[u]<<
"\n"
551 <<
" N = "<<m_N[u]<<
"\n"
552 <<
" x0 = "<<m_x0[u]<<
"\n"
553 <<
" x1 = "<<m_x1[u]<<
"\n"
554 <<
" h = "<<
h(u)<<
"\n"
555 <<
" l = "<<
l(u)<<
"\n"
556 <<
" bc = "<<
bc2str(m_bcs[u])<<
"\n";
561 template<
size_t Md = Nd>
562 std::enable_if_t<(Md == 1),
bool>
contains( real_type x)
const
564 return contains( std::array<real_type,1>{
x});
577 template<
class Vector>
580 for(
unsigned u=0; u<Nd; u++)
582 if( !std::isfinite(
x[u]) )
return false;
584 if(
x[u] < m_x0[u])
return false;
585 if(
x[u] > m_x1[u])
return false;
605 std::array<real_type,Nd>
p,
606 std::array<real_type,Nd>
q,
607 std::array<unsigned,Nd>
n,
608 std::array<unsigned,Nd>
N,
609 std::array<dg::bc, Nd> bcs) : m_x0(
p), m_x1(
q), m_n(
n), m_N(
N), m_bcs(bcs)
618 for(
unsigned u=0; u<Nd; u++)
620 m_x0[u] = axes[u].p();
621 m_x1[u] = axes[u].q();
622 m_n[u] = axes[u].n();
623 m_N[u] = axes[u].N();
624 m_bcs[u] = axes[u].bc();
636 virtual void do_set(std::array<unsigned,Nd> new_n, std::array<unsigned,Nd> new_N) = 0;
638 virtual void do_set_pq( std::array<real_type, Nd> new_p, std::array<real_type,Nd> new_q) = 0;
640 virtual void do_set( std::array<dg::bc, Nd> new_bcs) = 0;
674 std::array<real_type,Nd> m_x0;
675 std::array<real_type,Nd> m_x1;
676 std::array<unsigned,Nd> m_n;
677 std::array<unsigned,Nd> m_N;
678 std::array<dg::bc,Nd> m_bcs;
682template<
class real_type,
size_t Nd>
688template<
class real_type,
size_t Nd>
694template<
class real_type,
size_t Nd>
708template<
class real_type,
size_t Nd>
726 template<
size_t Md = Nd>
733 template<
size_t Md = Nd>
735 aRealTopology<real_type,2>({
x0,
y0}, {
x1,
y1}, {
n,
n}, {
Nx,
Ny}, {
bcx,
bcy})
740 template<
size_t Md = Nd>
741 RealGrid( real_type
x0, real_type
x1, real_type
y0, real_type
y1, real_type
z0, real_type
z1,
unsigned n,
unsigned Nx,
unsigned Ny,
unsigned Nz,
dg::bc bcx =
PER,
dg::bc bcy =
PER,
dg::bc bcz=
PER):
742 aRealTopology<real_type,3>({
x0,
y0,
z0}, {
x1,
y1,
z1}, {
n,
n,1}, {
Nx,
Ny,
Nz}, {
bcx,
bcy,
bcz})
755 template<
class ...Grid1ds>
760 RealGrid( std::array<real_type,Nd>
p, std::array<real_type,Nd>
q,
761 std::array<unsigned,Nd>
n, std::array<unsigned,Nd>
N,
769 virtual void do_set( std::array<unsigned,Nd> new_n, std::array<unsigned,Nd> new_N)
override final{
772 virtual void do_set_pq( std::array<real_type,Nd> new_x0, std::array<real_type,Nd> new_x1)
override final{
775 virtual void do_set( std::array<dg::bc,Nd> new_bcs)
override final{
class intended for the use in throw statements
Definition exceptions.h:83
small class holding a stringstream
Definition exceptions.h:29
The discrete legendre trafo class.
#define _ping_
Definition exceptions.h:12
std::string bc2str(bc bcx)
write a string describing boundary condition to an output stream
Definition enums.h:37
bc
Switch between boundary conditions.
Definition enums.h:15
@ PER
periodic boundaries
Definition enums.h:16
aMPITopology3d aTopology3d
Definition mpi_grid.h:909
MPIGrid3d Grid3d
Definition mpi_grid.h:905
MPIGrid1d Grid1d
Definition mpi_grid.h:903
MPIGrid0d Grid0d
Definition mpi_grid.h:902
MPIGrid2d Grid2d
Definition mpi_grid.h:904
aMPITopology2d aTopology2d
Definition mpi_grid.h:908
This is the namespace for all functions and classes defined and used by the discontinuous Galerkin li...
static std::vector< real_type > abscissas(unsigned n)
Return Gauss-Legendre nodes on the interval [-1,1].
Definition dlt.h:27
static std::vector< real_type > weights(unsigned n)
Return Gauss-Legendre weights.
Definition dlt.h:83
The simplest implementation of aRealTopology.
Definition grid.h:710
RealGrid(real_type x0, real_type x1, unsigned n, unsigned Nx, dg::bc bcx=dg::PER)
1D grid
Definition grid.h:727
RealGrid(const RealGrid< real_type, 1 > &g0, const Grid1ds &...gs)
Construct from given 1d grids Equivalent to RealGrid( std::array{g0,gs...})
Definition grid.h:756
RealGrid(const aRealTopology< real_type, Nd > &src)
allow explicit type conversion from any other topology
Definition grid.h:767
RealGrid(const std::array< RealGrid< real_type, 1 >, Nd > &axes)
Construct a topology as the product of 1d axes grids.
Definition grid.h:746
RealGrid(std::array< real_type, Nd > p, std::array< real_type, Nd > q, std::array< unsigned, Nd > n, std::array< unsigned, Nd > N, std::array< dg::bc, Nd > bcs)
Construct a topology directly from points and dimensions.
Definition grid.h:760
RealGrid()=default
construct an empty grid this leaves the access functions undefined
RealGrid(real_type x0, real_type x1, real_type y0, real_type y1, real_type z0, real_type z1, unsigned n, unsigned Nx, unsigned Ny, unsigned Nz, dg::bc bcx=PER, dg::bc bcy=PER, dg::bc bcz=PER)
Construct with equal polynomial coefficients.
Definition grid.h:741
RealGrid(real_type x0, real_type x1, real_type y0, real_type y1, unsigned n, unsigned Nx, unsigned Ny, dg::bc bcx=PER, dg::bc bcy=PER)
Construct with equal polynomial coefficients.
Definition grid.h:734
An abstract base class for Nd-dimensional dG grids.
Definition grid.h:95
dg::bc bcz() const
Equivalent to bc(2)
Definition grid.h:350
RealGrid< real_type, 1 > gy() const
Equivalent to grid(1)
Definition grid.h:360
real_type h(unsigned u=0) const
Get grid constant for axis u.
Definition grid.h:256
real_type hy() const
Equivalent to h(1)
Definition grid.h:317
std::array< unsigned, Nd > get_n() const
Get number of polynomial coefficients for all axes.
Definition grid.h:234
static constexpr unsigned ndim()
Dimensionality == Nd.
Definition grid.h:106
real_type q(unsigned u=0) const
Get right boundary point for axis u.
Definition grid.h:253
real_type hz() const
Equivalent to h(2)
Definition grid.h:320
std::array< unsigned, Nd > start() const
Start coordinate in C-order for dg::file::NcHyperslab.
Definition grid.h:376
unsigned nz() const
Equivalent to n(2)
Definition grid.h:330
real_type lx() const
Equivalent to l(0)
Definition grid.h:304
real_type l(unsigned u=0) const
Get grid length for axis u.
Definition grid.h:259
real_type hx() const
Equivalent to h(0)
Definition grid.h:314
std::array< real_type, Nd > get_q() const
Get right boundary point .
Definition grid.h:207
real_type value_type
value type of abscissas and weights
Definition grid.h:98
real_type x0() const
Equivalent to p(0)
Definition grid.h:285
host_vector weights(unsigned u=0) const
Get the weights of the u axis.
Definition grid.h:155
unsigned ny() const
Equivalent to n(1)
Definition grid.h:327
RealGrid< real_type, 1 > gx() const
Equivalent to grid(0)
Definition grid.h:354
void display(std::ostream &os=std::cout) const
Display.
Definition grid.h:545
real_type z1() const
Equivalent to q(2)
Definition grid.h:300
unsigned size() const
The total number of points.
Definition grid.h:532
std::array< dg::bc, Nd > get_bc() const
Get boundary condition for all axes.
Definition grid.h:240
host_vector abscissas(unsigned u=0) const
Construct grid abscissas of the u axis.
Definition grid.h:128
std::array< host_vector, Nd > get_weights() const
Construct weights for all axes.
Definition grid.h:193
dg::bc bcy() const
Equivalent to bc(1)
Definition grid.h:347
real_type z0() const
Equivalent to q(1)
Definition grid.h:297
unsigned n(unsigned u=0) const
Get number of polynomial coefficients for axis u.
Definition grid.h:262
void set(unsigned new_n, std::array< unsigned, Nd > new_N)
Same as set( {new_n, new_n,...}, new_N);
Definition grid.h:453
unsigned Nx() const
Equivalent to N(0)
Definition grid.h:334
~aRealTopology()=default
disallow deletion through base class pointer
unsigned nx() const
Equivalent to n(0)
Definition grid.h:324
aRealTopology()=default
default constructor
real_type ly() const
Equivalent to l(1)
Definition grid.h:307
std::array< real_type, Nd > get_p() const
Get left boundary point .
Definition grid.h:202
RealGrid< real_type, 1 > grid(unsigned u) const
Get axis u as a 1d grid.
Definition grid.h:274
virtual void do_set_pq(std::array< real_type, Nd > new_p, std::array< real_type, Nd > new_q)=0
Reset the boundaries of the grid.
std::array< unsigned, Nd > get_N() const
Get number of cells for all axes.
Definition grid.h:228
std::enable_if_t<(Md==1), void > set(unsigned new_n, unsigned new_Nx)
Set n and N in a 1-dimensional grid.
Definition grid.h:418
real_type y0() const
Equivalent to p(2)
Definition grid.h:291
virtual void do_set(std::array< unsigned, Nd > new_n, std::array< unsigned, Nd > new_N)=0
Set the number of polynomials and cells.
RealGrid< real_type, 1 > gz() const
Equivalent to grid(2)
Definition grid.h:366
std::array< unsigned, Nd > count() const
Count vector in C-order for dg::file::NcHyperslab.
Definition grid.h:385
void set(std::array< unsigned, Nd > new_n, std::array< unsigned, Nd > new_N)
Set the number of polynomials and cells.
Definition grid.h:481
std::array< unsigned, Nd > get_shape() const
the total number of points of an axis
Definition grid.h:173
thrust::host_vector< real_type > host_vector
Definition grid.h:102
aRealTopology(const aRealTopology &src)=default
aRealTopology & operator=(const aRealTopology &src)=default
virtual void do_set(std::array< dg::bc, Nd > new_bcs)=0
Reset the boundary conditions of the grid.
real_type lz() const
Equivalent to l(2)
Definition grid.h:310
real_type p(unsigned u=0) const
Get left boundary point for axis u.
Definition grid.h:250
std::enable_if_t<(Md >=2), void > multiplyCellNumbers(real_type fx, real_type fy)
Multiply the number of cells in the first two dimensions with a given factor.
Definition grid.h:404
void set_pq(std::array< real_type, Nd > new_p, std::array< real_type, Nd > new_q)
Reset the boundaries of the grid.
Definition grid.h:494
std::array< real_type, Nd > get_h() const
Get grid constant for all axes.
Definition grid.h:220
void set_axis(unsigned coord, unsigned new_n, unsigned new_N)
Set n and N for axis coord.
Definition grid.h:467
real_type y1() const
Equivalent to q(0)
Definition grid.h:294
std::enable_if_t<(Md==3), void > set(unsigned new_n, unsigned new_Nx, unsigned new_Ny, unsigned new_Nz)
Set n and N in a 3-dimensional grid.
Definition grid.h:447
dg::bc bc(unsigned u=0) const
Get boundary condition for axis u.
Definition grid.h:268
std::array< host_vector, Nd > get_abscissas() const
Construct abscissas for all axes.
Definition grid.h:183
RealGrid< real_type, 1 > axis(unsigned u) const
An alias for "grid".
Definition grid.h:282
real_type x1() const
Equivalent to p(1)
Definition grid.h:288
bool contains(const Vector &x) const
Check if the grid contains a point.
Definition grid.h:578
std::array< real_type, Nd > get_l() const
Get grid length for all axes.
Definition grid.h:212
void set_bcs(std::array< dg::bc, Nd > new_bcs)
Reset the boundary conditions of the grid.
Definition grid.h:504
dg::bc bcx() const
Equivalent to bc(0)
Definition grid.h:344
aRealTopology(std::array< real_type, Nd > p, std::array< real_type, Nd > q, std::array< unsigned, Nd > n, std::array< unsigned, Nd > N, std::array< dg::bc, Nd > bcs)
Construct a topology directly from points and dimensions.
Definition grid.h:604
void set(std::array< real_type, Nd > new_p, std::array< real_type, Nd > new_q, std::array< unsigned, Nd > new_n, std::array< unsigned, Nd > new_N, std::array< dg::bc, Nd > new_bcs)
Reset the entire grid.
Definition grid.h:518
unsigned Ny() const
Equivalent to N(1)
Definition grid.h:337
unsigned Nz() const
Equivalent to N(2)
Definition grid.h:340
std::enable_if_t<(Md==1), bool > contains(real_type x) const
Check if the grid contains a point.
Definition grid.h:562
unsigned N(unsigned u=0) const
Get number of cells for axis u.
Definition grid.h:265
unsigned shape(unsigned u=0) const
the total number of points of an axis
Definition grid.h:114
aRealTopology(const std::array< RealGrid< real_type, 1 >, Nd > &axes)
Construct a topology as the product of 1d axes grids.
Definition grid.h:616
std::enable_if_t<(Md==2), void > set(unsigned new_n, unsigned new_Nx, unsigned new_Ny)
Set n and N in a 2-dimensional grid.
Definition grid.h:431