This is the namespace for all functions and classes defined and used by the discontinuous Galerkin library.
More...
|
struct | ABS |
| \( f(x) = |x|\) More...
|
|
struct | AbsMax |
| \( f(x,y) = \max(|x|,|y|)\) More...
|
|
struct | AbsMin |
| \( f(x,y) = \min(|x|,|y|)\) More...
|
|
struct | Adaptive |
| Driver class for adaptive timestep ODE integration. More...
|
|
struct | AdaptiveTimeloop |
| Integrate using a while loop. More...
|
|
struct | Advection |
| Upwind discretization of advection operator \( \vec v\cdot\nabla f\) More...
|
|
struct | AndersonAcceleration |
| Anderson Acceleration of Fixed Point/Richardson Iteration for the nonlinear equation \( f(x) = b\). More...
|
|
struct | AnyMatrixTag |
| tensor_category base class More...
|
|
struct | AnyPolicyTag |
| Execution Policy base class. More...
|
|
struct | AnyScalarTag |
| Scalar Tag base class. More...
|
|
struct | AnyVectorTag |
| Vector Tag base class, indicates the basic Vector/container concept. More...
|
|
struct | ArakawaX |
| Arakawa's scheme for Poisson bracket \( \{ f,g\} \). More...
|
|
struct | aRealGeometry2d |
| This is the abstract interface class for a two-dimensional Geometry. More...
|
|
struct | aRealGeometry3d |
| This is the abstract interface class for a three-dimensional Geometry. More...
|
|
struct | aRealGeometryX2d |
| This is the abstract interface class for a two-dimensional RealGeometryX. More...
|
|
struct | aRealGeometryX3d |
| This is the abstract interface class for a three-dimensional RealGeometryX. More...
|
|
struct | aRealMPIGeometry2d |
| This is the abstract interface class for a two-dimensional Geometry. More...
|
|
struct | aRealMPIGeometry3d |
| This is the abstract interface class for a three-dimensional MPIGeometry. More...
|
|
struct | aRealMPITopology |
| An abstract base class for MPI distributed Nd-dimensional dG grids. More...
|
|
struct | aRealProductGeometry3d |
| A 3d product space Geometry \( g_{2d} \otimes g_{1d}\). More...
|
|
struct | aRealProductMPIGeometry3d |
| A 3d product space Geometry \( g_{2d} \otimes g_{1d}\). More...
|
|
struct | aRealRefinement1d |
| Abstract base class for 1d grid refinement that increases the number of grid cells of a fixed basis grid. More...
|
|
struct | aRealRefinementX2d |
| Abstract base class for 2d grid refinement that increases the number of grid cells of a fixed basis grid. More...
|
|
struct | aRealTopology |
| An abstract base class for Nd-dimensional dG grids. More...
|
|
struct | aRealTopologyX2d |
| A 2D grid class with X-point topology. More...
|
|
struct | aRealTopologyX3d |
| A 3D grid class with X-point topology. More...
|
|
struct | ArithmeticTag |
| Types where std::is_arithmetic is true. More...
|
|
struct | ARKStep |
| Additive Runge Kutta (semi-implicit) time-step with error estimate following The ARKode library More...
|
|
struct | ArrayScalarTag |
|
struct | ArrayVectorTag |
|
struct | aTimeloop |
| Abstract timeloop independent of stepper and ODE. More...
|
|
struct | Average |
| Topological average computations in a Cartesian topology. More...
|
|
struct | BathRZ |
| \(f(R,Z) = A B \sum_\vec{k} \sqrt{E_k} \alpha_k \cos{\left(k \kappa_k + \theta_k \right)}
\) More...
|
|
class | BICGSTABl |
| Preconditioned BICGSTAB(l) method to solve \( Ax=b\). More...
|
|
struct | ButcherTableau |
| Manage coefficients of a (extended) Butcher tableau. More...
|
|
struct | Cauchy |
| \(
f(x,y) = \begin{cases}
Ae^{1 + \left(\frac{(x-x_0)^2}{\sigma_x^2} + \frac{(y-y_0)^2}{\sigma_y^2} - 1\right)^{-1}} \text{ if } \frac{(x-x_0)^2}{\sigma_x^2} + \frac{(y-y_0)^2}{\sigma_y^2} < 1\\
0 \text{ else}
\end{cases}
\) More...
|
|
struct | CauchyX |
| \(
f(x,y) = \begin{cases}
Ae^{1 + \left(\frac{(x-x_0)^2}{\sigma_x^2} - 1\right)^{-1}} \text{ if } \frac{(x-x_0)^2}{\sigma_x^2} < 1\\
0 \text{ else}
\end{cases}
\) More...
|
|
class | ChebyshevIteration |
| Preconditioned Chebyshev iteration for solving \( PAx=Pb\). More...
|
|
struct | ChebyshevPreconditioner |
| Chebyshev Polynomial Preconditioner \( C( A)\). More...
|
|
struct | ClonePtr |
| Manager class that invokes the clone() method on the managed ptr when copied. More...
|
|
struct | ComplexTag |
| complex number type More...
|
|
struct | Composite |
|
struct | CONSTANT |
| \( f(x,...) = c\) More...
|
|
struct | ConvertsToButcherTableau |
| Convert identifiers to their corresponding dg::ButcherTableau . More...
|
|
struct | ConvertsToMultistepTableau |
| Convert identifiers to their corresponding dg::MultistepTableau . More...
|
|
struct | ConvertsToShuOsherTableau |
| Convert identifiers to their corresponding dg::ShuOsherTableau . More...
|
|
struct | CooSparseBlockMat |
| Coo Sparse Block Matrix format. More...
|
|
struct | CosXCosY |
| \( f(x,y) =B+ A \cos(k_x x) \cos(k_y y) \) More...
|
|
struct | CosY |
| \( f(x,y) =B+ A \cos(k_y y) \) More...
|
|
struct | CSRAverageFilter |
| Average filter that computes the average of all points in the stencil More...
|
|
struct | CSRMedianFilter |
| Compute (lower) Median of input numbers. More...
|
|
struct | CSRSlopeLimiter |
| Generalized slope limiter for dG methods. More...
|
|
struct | CSRSWMFilter |
| Switching median filter. More...
|
|
struct | CSRSymvFilter |
| Test filter that computes the symv csr matrix-vector product if used. More...
|
|
struct | CudaTag |
| CUDA implementation. More...
|
|
struct | CuspMatrixTag |
| One of cusp's matrices, for these only the blas2 transfer and the symv( m,x,y) are implemented. More...
|
|
struct | CuspVectorTag |
| special tag for cusp arrays More...
|
|
struct | DefaultSolver |
| PCG Solver class for solving \( (y-\alpha\hat I(t,y)) = \rho\). More...
|
|
struct | DenseMatrixTag |
| indicate our dense matrix format More...
|
|
struct | DIRKStep |
| Embedded diagonally implicit Runge Kutta time-step with error estimate \(
\begin{align}
k_i = f\left( t^n + c_i \Delta t, u^n + \Delta t \sum_{j=1}^{i} a_{ij} k_j\right) \\
u^{n+1} = u^{n} + \Delta t\sum_{j=1}^s b_j k_j \\
\tilde u^{n+1} = u^{n} + \Delta t\sum_{j=1}^s \tilde b_j k_j \\
\delta^{n+1} = u^{n+1} - \tilde u^{n+1} = \Delta t\sum_{j=1}^s (b_j-\tilde b_j) k_j
\end{align}
\). More...
|
|
struct | Distance |
| \( f(x,y) = \sqrt{ (x-x_0)^2 + (y-y_0)^2} \) More...
|
|
struct | divides_equals |
| \( y = y/x\) More...
|
|
struct | DLT |
| Struct providing coefficients for Discrete Legendre Transformation (DLT) related operations. More...
|
|
struct | DPolynomialHeaviside |
| \( f(x) = \begin{cases}
0 \text{ if } x < x_b-a || x > x_b+a \\
(35 (a + x - x_b)^3 (a - x + x_b)^3)/(32 a^7)
\text{ if } |x-x_b| < a
\end{cases}\) The derivative of PolynomialHeaviside approximates delta(x) More...
|
|
class | Elliptic1d |
| A 1d negative elliptic differential operator \( -\partial_x ( \chi \partial_x ) \). More...
|
|
class | Elliptic2d |
| A 2d negative elliptic differential operator \( -\nabla \cdot ( \mathbf{\chi}\cdot \nabla ) \). More...
|
|
class | Elliptic3d |
| A 3d negative elliptic differential operator \( -\nabla \cdot ( \mathbf{\chi}\cdot \nabla ) \). More...
|
|
struct | EllSparseBlockMat |
| Ell Sparse Block Matrix format. More...
|
|
struct | EmbeddedPairSum |
| \( y = \sum_i b_i x_i + b_0 y,\quad \tilde y = \sum_i \tilde b_i x_i + \tilde b_0 y \) More...
|
|
struct | EntireDomain |
| The domain that contains all points. More...
|
|
struct | equals |
| \( y=x\) More...
|
|
struct | ERKStep |
| Embedded Runge Kutta explicit time-step with error estimate \(
\begin{align}
k_i = f\left( t^n + c_i \Delta t, u^n + \Delta t \sum_{j=1}^{i-1} a_{ij} k_j\right) \\
u^{n+1} = u^{n} + \Delta t\sum_{j=1}^s b_j k_j \\
\tilde u^{n+1} = u^{n} + \Delta t\sum_{j=1}^s \tilde b_j k_j \\
\delta^{n+1} = u^{n+1} - \tilde u^{n+1} = \Delta t\sum_{j=1}^s (b_j - \tilde b_j) k_j
\end{align}
\). More...
|
|
class | Error |
| class intended for the use in throw statements More...
|
|
class | EVE |
| The Eigen-Value-Estimator (EVE) finds largest Eigenvalue of \( M^{-1}Ax = \lambda_\max x\). More...
|
|
struct | EXP |
| \( f(x) = \exp( x)\) More...
|
|
struct | ExplicitMultistep |
| General explicit linear multistep ODE integrator \(
\begin{align}
v^{n+1} = \sum_{j=0}^{s-1} a_j v^{n-j} + \Delta t\left(\sum_{j=0}^{s-1}b_j \hat f\left(t^{n}-j\Delta t, v^{n-j}\right)\right)
\end{align}
\). More...
|
|
struct | ExponentialFilter |
| \( f(i) = \begin{cases}
1 \text{ if } \eta < \eta_c \\
\exp\left( -\alpha \left(\frac{\eta-\eta_c}{1-\eta_c} \right)^{2s}\right) \text { if } \eta \geq \eta_c \\
0 \text{ else} \\
\eta=\frac{i}{1-n}
\end{cases}\) More...
|
|
struct | ExpProfX |
| \( f(x) = f(x,y) = f(x,y,z) = A\exp(-x/L_n) + B \) More...
|
|
struct | Extrapolation |
| Extrapolate a polynomial passing through up to three points. More...
|
|
struct | Fail |
| Indicate failure to converge. More...
|
|
struct | FilteredERKStep |
| EXPERIMENTAL: Filtered Embedded Runge Kutta explicit time-step with error estimate \(
\begin{align}
k_i = f\left( t^n + c_i \Delta t, \Lambda\Pi \left[u^n + \Delta t \sum_{j=1}^{i-1} a_{ij} k_j\right]\right) \\
u^{n+1} = \Lambda\Pi\left[u^{n} + \Delta t\sum_{j=1}^s b_j k_j\right] \\
\delta^{n+1} = \Delta t\sum_{j=1}^s (\tilde b_j - b_j) k_j
\end{align}
\). More...
|
|
struct | FilteredExplicitMultistep |
| EXPERIMENTAL: General explicit linear multistep ODE integrator with Limiter / Filter \(
\begin{align}
\tilde v &= \sum_{j=0}^{s-1} a_j v^{n-j} + \Delta t\left(\sum_{j=0}^{s-1}b_j \hat f\left(t^{n}-j\Delta t, v^{n-j}\right)\right) \\
v^{n+1} &= \Lambda\Pi \left( \tilde v\right)
\end{align}
\). More...
|
|
struct | FloatingPointTag |
| Types where std::is_floating_point is true. More...
|
|
struct | Gaussian |
| \(
f(x,y) = Ae^{-\left(\frac{(x-x_0)^2}{2\sigma_x^2} + \frac{(y-y_0)^2}{2\sigma_y^2}\right)}
\) More...
|
|
struct | Gaussian3d |
| \(
f(x,y,z) = Ae^{-\left(\frac{(x-x_0)^2}{2\sigma_x^2} + \frac{(y-y_0)^2}{2\sigma_y^2} + \frac{(z-z_0)^2}{2\sigma_z^2}\right)}
\) More...
|
|
struct | GaussianDamping |
| \( f(\psi) = \begin{cases}
1 \text{ if } \psi < \psi_{\max}\\
0 \text{ if } \psi > (\psi_{\max} + 4\alpha) \\
\exp\left( - \frac{(\psi - \psi_{\max})^2}{2\alpha^2}\right), \text{ else}
\end{cases}
\) More...
|
|
struct | GaussianX |
| \(
f(x,y) = Ae^{-\frac{(x-x_0)^2}{2\sigma_x^2} }
\) More...
|
|
struct | GaussianY |
| \(
f(x,y) = Ae^{-\frac{(y-y_0)^2}{2\sigma_y^2}}
\) More...
|
|
struct | GaussianZ |
| \(
f(x,y,z) = Ae^{-\frac{(z-z_0)^2}{2\sigma_z^2}}
\) More...
|
|
struct | GeneralHelmholtz |
| A general Helmholtz-type operator \( (\chi-\alpha F) \). More...
|
|
struct | Heaviside |
| \( f(x) = \begin{cases}
0 \text{ if } x < x_b \\
1 \text{ else}
\end{cases}\) More...
|
|
struct | Helmholtz2 |
| DEPRECATED, Matrix class that represents a more general Helmholtz-type operator. More...
|
|
struct | Histogram |
| Compute a histogram on a 1D grid. More...
|
|
struct | Histogram2D |
| Compute a histogram on a 2D grid. More...
|
|
struct | Horner2d |
| \( f(x,y) = \sum_{i=0}^{M-1} \sum_{j=0}^{N-1} c_{iN+j} x^i y^j \) More...
|
|
struct | IDENTITY |
| \( f(x) = x\) More...
|
|
struct | IdentityFilter |
| A filter that does nothing. More...
|
|
struct | ImExMultistep |
| Semi-implicit multistep ODE integrator \(
\begin{align}
v^{n+1} = \sum_{q=0}^{s-1} a_q v^{n-q} + \Delta t\left[\left(\sum_{q=0}^{s-1}b_q \hat E(t^{n}-q\Delta t, v^{n-q}) + \sum_{q=1}^{s} c_q \hat I( t^n - q\Delta t, v^{n-q})\right) + c_0\hat I(t^{n}+\Delta t, v^{n+1})\right]
\end{align}
\). More...
|
|
struct | ImplicitMultistep |
| Implicit multistep ODE integrator \(
\begin{align}
v^{n+1} &= \sum_{i=0}^{s-1} a_i v^{n-i} + \Delta t \sum_{i=1}^{s} c_i\hat I(t^{n+1-i}, v^{n+1-i}) + \Delta t c_{0} \hat I (t + \Delta t, v^{n+1}) \\
\end{align}
\). More...
|
|
struct | IntegralTag |
| Types where std::is_integral is true. More...
|
|
struct | InvCoshXsq |
| \( f(x,y) =A/\cosh^2(k_x x) \) More...
|
|
struct | InverseKroneckerTriDiagonal2d |
| Fast inverse tridiagonal sparse matrix in 2d \( T_y^{-1}\otimes T_x^{-1}\). More...
|
|
struct | InverseTriDiagonal |
| DEPRECATED/UNTESTED Fast inverse tridiagonal sparse matrix. More...
|
|
struct | INVERT |
| \( f(x) = 1/x \) More...
|
|
struct | InvSqrt |
| \( f(x) = \frac{1}{\sqrt{x}}\) More...
|
|
struct | IPolynomialHeaviside |
| \( f(x) = \begin{cases}
x_b \text{ if } x < x_b-a \\
x_b + ((35 a^3 - 47 a^2 (x - x_b) + 25 a (x - x_b)^2 - 5 (x - x_b)^3) (a + x - x_b)^5)/(256 a^7)
\text{ if } |x-x_b| < a \\
x \text{ if } x > x_b + a
\end{cases}\) The integral of PolynomialHeaviside approximates xH(x) More...
|
|
struct | Iris |
| \( f(\psi) = \begin{cases}
1 \text{ if } \psi_{\min} < \psi < \psi_{\max}\\
0 \text{ else}
\end{cases}\) More...
|
|
struct | IslandXY |
| \( f(x,y) = \lambda \ln{(\cosh{(x/\lambda) } +\epsilon \cos(y/\lambda)) } \) More...
|
|
struct | ISNFINITE |
| \( f(x) = \mathrm{!std::isfinite(x)}\) More...
|
|
struct | ISNSANE |
| \( f(x) =\begin{cases} \mathrm{true\ if}\ |x| > 10^{100}\\
\mathrm{false\ else}
\end{cases}\) More...
|
|
struct | KroneckerTriDiagonal2d |
| Fast tridiagonal sparse matrix in 2d \( T_y\otimes T_x\). More...
|
|
struct | Lamb |
| \( f(x,y) = \begin{cases} 2\lambda U J_1(\lambda r) / J_0(\gamma)\cos(\theta) \text{ for } r<R \\
0 \text{ else}
\end{cases}
\) More...
|
|
struct | LeastSquaresExtrapolation |
| Evaluate a least squares fit More...
|
|
struct | LeastSquaresPreconditioner |
| Least Squares Polynomial Preconditioner \( M^{-1} s( AM^{-1})\). More...
|
|
class | LGMRES |
| Functor class for the right preconditioned LGMRES method to solve \( Ax=b\). More...
|
|
struct | Line |
| \( f(x) = y_1\frac{x-x_0}{x_1-x_0} + y_0\frac{x-x_1}{x_0-x_1}\) More...
|
|
struct | LinearX |
| \( f(x) = f(x,y) = f(x,y,z) = ax+b \) More...
|
|
struct | LinearY |
| \( f(x,y) = f(x,y,z) = ay+b \) More...
|
|
struct | LinearZ |
| \( f(x,y,z) = az+b \) More...
|
|
struct | LN |
| \( f(x) = \ln(x)\) More...
|
|
class | Message |
| small class holding a stringstream More...
|
|
struct | MinMod |
| \( f(x_1, x_2, ...) = \begin{cases}
\min(x_1, x_2, ...) &\text{ for } x_1, x_2, ... >0 \\
\max(x_1, x_2, ...) &\text{ for } x_1, x_2, ... <0 \\
0 &\text{ else}
\end{cases}
\) More...
|
|
struct | minus_equals |
| \( y=y-x\) More...
|
|
struct | MOD |
| \( f(x) = \) x mod m > 0 ? x mod m : x mod m + m More...
|
|
struct | ModifiedChebyshevPreconditioner |
| Approximate inverse Chebyshev Polynomial Preconditioner \( A^{-1} = \frac{c_0}{2} I + \sum_{k=1}^{r}c_kT_k( Z)\). More...
|
|
struct | MPI_Vector |
| A simple wrapper around a container object and an MPI_Comm. More...
|
|
struct | MPIDistMat |
| Distributed memory matrix class, asynchronous communication. More...
|
|
struct | MPIGather |
| Optimized MPI Gather operation. More...
|
|
struct | MPIKroneckerGather |
| Communicator for asynchronous communication of MPISparseBlockMat . More...
|
|
struct | MPIMatrixTag |
| indicate one of our mpi matrices More...
|
|
struct | MPISparseBlockMat |
| Distributed memory Sparse block matrix class, asynchronous communication. More...
|
|
struct | MPIVectorTag |
| A distributed vector contains a data container and a MPI communicator. More...
|
|
struct | MultigridCG2d |
| Solve. More...
|
|
struct | MultiMatrix |
| Struct that applies given matrices one after the other. More...
|
|
struct | MultistepTableau |
| Manage coefficients of Multistep methods. More...
|
|
struct | MultistepTimeloop |
| Integrate using a for loop and a fixed non-changeable time-step. More...
|
|
struct | NestedGrids |
| Hold nested grids and provide dg fast interpolation and projection matrices. More...
|
|
struct | NoPolicyTag |
| Indicate that a type does not have an execution policy. More...
|
|
class | NoRoot1d |
| Exception class, that stores boundaries for 1D root finding. More...
|
|
struct | NotATensorTag |
| Indicate that a type is not a tensor. More...
|
|
struct | OmpTag |
| OpenMP parallel execution. More...
|
|
struct | ONE |
| \( f(x,...) = 1\) More...
|
|
struct | PairSum |
| \( y = \sum_i a_i x_i \) More...
|
|
class | PCG |
| Preconditioned conjugate gradient method to solve \( Ax=b\). More...
|
|
struct | PLUS |
| \( f(x) = x + c\) More...
|
|
struct | plus_equals |
| \( y=y+x\) More...
|
|
struct | Poisson |
| Direct discretization of Poisson bracket \( \{ f,g\} \). More...
|
|
struct | PolynomialHeaviside |
| \( f(x) = \begin{cases}
0 \text{ if } x < x_b-a \\
((16 a^3 - 29 a^2 (x - x_b) + 20 a (x - x_b)^2 - 5 (x - x_b)^3) (a + x -
x_b)^4)/(32 a^7) \text{ if } |x-x_b| < a \\
1 \text{ if } x > x_b + a
\end{cases}\) More...
|
|
struct | PolynomialRectangle |
| \( f(x) = \begin{cases}
0 \text{ if } x < x_l-a_l \\
((16 a_l^3 - 29 a_l^2 (x - x_l) + 20 a_l (x - x_l)^2 - 5 (x - x_l)^3) (a_l + x -
x_l)^4)/(32 a_l^7) \text{ if } |x-x_l| < a_l \\
1 \text{ if } x_l + a_l < x < x_r-a_r \\
((16 a_r^3 - 29 a_r^2 (x - x_r) + 20 a_r (x - x_r)^2 - 5 (x - x_r)^3) (a_r + x -
x_l)^4)/(32 a_r^7) \text{ if } |x-x_r| < a_r \\
0 \text{ if } x > x_r + a_r
\end{cases}\) More...
|
|
struct | POSVALUE |
| \( f(x) = \begin{cases}
x \text{ for } x>0 \\
0 \text{ else}
\end{cases}
\) More...
|
|
struct | Product |
| \( y = \prod_i x_i \) More...
|
|
struct | PsiPupil |
| \( f(\psi) = \begin{cases}
\psi_{\max} \text{ if } \psi > \psi_{\max} \\
\psi \text{ else}
\end{cases}\) More...
|
|
struct | Pupil |
| \( f(\psi) = \begin{cases}
0 \text{ if } \psi > \psi_{\max} \\
1 \text{ else}
\end{cases}\) More...
|
|
struct | RealCartesianGrid2d |
| Two-dimensional Grid with Cartesian metric. More...
|
|
struct | RealCartesianGrid3d |
| Three-dimensional Grid with Cartesian metric. More...
|
|
struct | RealCartesianGridX2d |
| two-dimensional GridX with RealCartesian metric More...
|
|
struct | RealCartesianGridX3d |
| three-dimensional GridX with RealCartesian metric More...
|
|
struct | RealCartesianMPIGrid2d |
| The mpi version of RealCartesianGrid2d. More...
|
|
struct | RealCartesianMPIGrid3d |
| The mpi version of RealCartesianGrid3d. More...
|
|
struct | RealCartesianRefinedGrid2d |
| Refined RealCartesian grid. More...
|
|
struct | RealCartesianRefinedGrid3d |
| Refined RealCartesian grid. More...
|
|
struct | RealCartesianRefinedGridX2d |
| Refined X-point grid. More...
|
|
struct | RealCartesianRefinedGridX3d |
| Refined X-point grid. More...
|
|
struct | RealCylindricalGrid3d |
| three-dimensional Grid with Cylindrical metric More...
|
|
struct | RealCylindricalMPIGrid3d |
| the mpi version of RealCylindricalGrid3d More...
|
|
struct | RealEquidistRefinement |
| Cell refinement around a given node. More...
|
|
struct | RealEquidistXRefinement |
| RealEquidistant cell refinement around the X-point. More...
|
|
struct | RealExponentialRefinement |
| The exponential refinement around a node. More...
|
|
struct | RealExponentialXRefinement |
| The exponential refinement around the X-point. More...
|
|
struct | RealFemRefinement |
| Insert equidistant points in between dG nodes. More...
|
|
struct | RealGrid |
| The simplest implementation of aRealTopology. More...
|
|
struct | RealGridX1d |
| 1D grid for X-point topology More...
|
|
struct | RealGridX2d |
| The simplest implementation of aRealTopologyX2d. More...
|
|
struct | RealGridX3d |
| The simplest implementation of aRealTopologyX3d. More...
|
|
struct | RealIdentityRefinement |
| No refinement. More...
|
|
struct | RealIdentityXRefinement |
| No refinement. More...
|
|
struct | RealLinearRefinement |
| Multiply every cell in the grid by a factor. More...
|
|
struct | RealMPIGrid |
| The simplest implementation of aRealMPITopology3d. More...
|
|
struct | RecursiveVectorTag |
| This tag indicates composition/recursion. More...
|
|
class | RefinedElliptic |
| The refined version of Elliptic . More...
|
|
struct | SelfMadeMatrixTag |
| Indicates that the type has a member function with the same name and interface (up to the matrix itself of course) as the corresponding blas2 member function, for example void symv( const ContainerType1&, ContainerType2& ); More...
|
|
struct | SerialTag |
| Indicate sequential execution. More...
|
|
struct | SharedVectorTag |
| Indicate a contiguous chunk of shared memory. More...
|
|
struct | ShuOsher |
| Shu-Osher fixed-step explicit ODE integrator with Slope Limiter / Filter \(
\begin{align}
u_0 &= u_n \\
u_i &= \Lambda\Pi \left(\sum_{j=0}^{i-1}\left[ \alpha_{ij} u_j + \Delta t \beta_{ij} f( t_j, u_j)\right]\right)\\
u^{n+1} &= u_s
\end{align}
\). More...
|
|
struct | ShuOsherTableau |
| Manage coefficients in Shu-Osher form. More...
|
|
struct | Sign |
| \( f(x) = \text{sgn}(x) = \begin{cases}
-1 \text{ for } x < 0 \\
0 \text{ for } x = 0 \\
+1 \text{ for } x > 0
\end{cases}\) More...
|
|
struct | Simpsons |
| Time integration based on Simpson's rule. More...
|
|
struct | SinglestepTimeloop |
| Integrate using a for loop and a fixed time-step. More...
|
|
struct | SinProfX |
| \( f(x) = f(x,y) = f(x,y,z) = B + A(1 - \sin(k_xx )) \) More...
|
|
struct | SinX |
| \( f(x) = f(x,y) = f(x,y,z) =B+ A \sin(k_x x) \) More...
|
|
struct | SinXCosY |
| \( f(x,y) =B+ A \sin(k_x x) \cos(k_y y) \) More...
|
|
struct | SinXSinY |
| \( f(x,y) =B+ A \sin(k_x x) \sin(k_y y) \) More...
|
|
struct | SinY |
| \( f(x,y) =B+ A \sin(k_y y) \) More...
|
|
struct | SlopeLimiter |
| \( \text{up}(v, g_m, g_0, g_p, h_m, h_p ) = \begin{cases} +h_m \Lambda( g_0, g_m) &\text{ if } v \geq 0 \\
-h_p \Lambda( g_p, g_0) &\text{ else}
\end{cases}
\) More...
|
|
struct | SlopeLimiterProduct |
| \( \text{up}(v, g_m, g_0, g_p, h_m, h_p ) = v \begin{cases} +h_m \Lambda( g_0, g_m) &\text{ if } v \geq 0 \\
-h_p \Lambda( g_p, g_0) &\text{ else}
\end{cases}
\) More...
|
|
struct | SparseBlockMatrixTag |
| indicate our sparse block matrix format More...
|
|
struct | SparseMatrix |
| A CSR formatted sparse matrix. More...
|
|
struct | SparseMatrixTag |
| indicate our sparse matrix format More...
|
|
struct | SparseTensor |
| Class for 2x2 and 3x3 matrices sharing elements. More...
|
|
struct | SQRT |
| \( f(x) = \sqrt{x}\) More...
|
|
struct | Square |
| \( f(x) = x^2\) More...
|
|
class | SquareMatrix |
| A square nxn matrix. More...
|
|
struct | StdMapTag |
|
struct | Sum |
| \( y = \sum_i x_i \) More...
|
|
struct | TanhProfX |
| \( f(x) = B + 0.5 A(1+ \text{sign} \tanh((x-x_b)/\alpha ) ) \) More...
|
|
struct | TensorTraits |
| The vector traits. More...
|
|
struct | TensorTraits< CooSparseBlockMat< T, V > > |
|
struct | TensorTraits< cusp::coo_matrix< I, V, M > > |
|
struct | TensorTraits< cusp::csr_matrix< I, V, M > > |
|
struct | TensorTraits< cusp::dia_matrix< I, V, M > > |
|
struct | TensorTraits< cusp::ell_matrix< I, V, M > > |
|
struct | TensorTraits< cusp::hyb_matrix< I, V, M > > |
|
struct | TensorTraits< EllSparseBlockMat< T, V > > |
|
struct | TensorTraits< mat::PolChargeN< G, M, V > > |
|
struct | TensorTraits< MPI_Vector< container > > |
| prototypical MPI vector More...
|
|
struct | TensorTraits< MPIDistMat< V, LI, LO > > |
|
struct | TensorTraits< MPISparseBlockMat< V, LI, LO > > |
|
struct | TensorTraits< SparseMatrix< I, T, V > > |
|
struct | TensorTraits< SquareMatrix< T > > |
|
struct | TensorTraits< std::array< T, N >, std::enable_if_t< !dg::is_scalar< T >::value > > |
|
struct | TensorTraits< std::array< T, N >, std::enable_if_t< dg::is_scalar< T >::value > > |
|
struct | TensorTraits< std::complex< T > > |
| Enable std::complex as a Scalar. More...
|
|
struct | TensorTraits< std::map< Key, T > > |
| Behaves like a RecursiveVector. More...
|
|
struct | TensorTraits< std::vector< T >, std::enable_if_t< !dg::is_scalar< T >::value > > |
| Prototypical Recursive Vector (unless vector of scalars) More...
|
|
struct | TensorTraits< std::vector< T >, std::enable_if_t< dg::is_scalar< T >::value > > |
|
struct | TensorTraits< T, std::enable_if_t< std::is_floating_point_v< T > > > |
| Enable double and float as a floating point. More...
|
|
struct | TensorTraits< T, std::enable_if_t< std::is_integral_v< T > > > |
| Enable integers and anything promotable to integer (such as bool and char) as integral. More...
|
|
struct | TensorTraits< thrust::complex< T > > |
|
struct | TensorTraits< thrust::device_vector< T > > |
| prototypical Shared Vector with Cuda or Omp Tag More...
|
|
struct | TensorTraits< thrust::host_vector< T > > |
| prototypical Shared Vector with Serial Tag More...
|
|
struct | TensorTraits< View< ThrustVector > > |
| A View has identical value_type and execution_policy as the underlying container. More...
|
|
struct | ThrustVectorTag |
| Indicate thrust/std - like behaviour. More...
|
|
struct | Timer |
| Simple tool for performance measuring. More...
|
|
struct | times_equals |
| \( y=xy\) More...
|
|
struct | TriDiagonal |
| Fast (shared memory) tridiagonal sparse matrix. More...
|
|
struct | TripletSum |
| \( y = \sum_i a_i x_i y_i \) More...
|
|
struct | Upwind |
| \( \text{up}(v, b, f ) = \begin{cases} b &\text{ if } v \geq 0 \\
f &\text{ else}
\end{cases}
\) More...
|
|
struct | UpwindProduct |
| \( \text{up}(v, b, f ) = v \begin{cases} b &\text{ if } v \geq 0 \\
f &\text{ else}
\end{cases}
\) More...
|
|
struct | VanLeer |
| \( f(x_1,x_2) = 2\begin{cases}
\frac{x_1x_2}{x_1+x_2} &\text{ if } x_1x_2 > 0 \\
0 & \text { else }
\end{cases}
\) More...
|
|
struct | View |
| A vector view class, usable in dg functions. More...
|
|
struct | Vortex |
| \(f(x,y) =\begin{cases}
\frac{u_d}{1.2965125} \left(
r\left(1+\frac{\beta_i^2}{g_i^2}\right)
- R \frac{\beta_i^2}{g_i^2} \frac{J_1(g_ir/R)}{J_1(g_i)}\right)\cos(\theta) \text{ if } r < R \\
\frac{u_d}{1.2965125} R \frac{K_1(\beta_i {r}/{R})}{K_1(\beta)} \cos(\theta) \text{ else }
\end{cases}
\) More...
|
|
struct | WallDistance |
| Shortest Distance to a collection of vertical and horizontal lines. More...
|
|
struct | WhichType |
|
struct | ZERO |
| \( f(x, ...) = 0\) More...
|
|
|
template<class T , class Tag = AnyScalarTag> |
using | is_scalar = typename std::is_base_of<Tag, get_tensor_category<T>>::type |
| Does a type have a tensor_category derived from Tag .
|
|
template<class T , class Tag = AnyVectorTag> |
using | is_vector = typename std::is_base_of<Tag, get_tensor_category<T>>::type |
| Does a type have a tensor_category derived from Tag .
|
|
template<class T , class Tag = AnyMatrixTag> |
using | is_matrix = typename std::is_base_of<Tag, get_tensor_category<T>>::type |
| Does a type have a tensor_category derived from Tag .
|
|
template<class T , class Tag = AnyPolicyTag> |
using | has_policy = std::is_same<Tag, get_execution_policy<T>> |
| Does a type have an execution_policy equal to Tag .
|
|
template<class Vector > |
using | get_value_type = typename TensorTraits<std::decay_t<Vector>>::value_type |
|
template<class Vector > |
using | get_tensor_category = typename TensorTraits< std::decay_t<Vector>>::tensor_category |
|
template<class Vector > |
using | get_execution_policy = typename TensorTraits<std::decay_t<Vector>>::execution_policy |
|
template<class T > |
using | HVec_t = thrust::host_vector<T> |
| Host Vector.
|
|
using | HVec = thrust::host_vector<double> |
| Host Vector.
|
|
using | cHVec = thrust::host_vector<thrust::complex<double>> |
| complex Host Vector
|
|
using | iHVec = thrust::host_vector<int> |
| integer Host Vector
|
|
using | fHVec = thrust::host_vector<float> |
| Host Vector.
|
|
using | DVec = thrust::device_vector<double> |
| Device Vector. The device can be an OpenMP parallelized cpu or a gpu. This depends on the value of the macro THRUST_DEVICE_SYSTEM, which can be either THRUST_DEVICE_SYSTEM_OMP for openMP or THRUST_DEVICE_SYSTEM_CUDA for a gpu or THRUST_DEVICE_SYSTEM_CPP for a cpu.
|
|
using | iDVec = thrust::device_vector<int> |
| integer Device Vector
|
|
using | cDVec = thrust::device_vector<thrust::complex<double>> |
| complex Device Vector
|
|
using | fDVec = thrust::device_vector<float> |
| Device Vector. The device can be an OpenMP parallelized cpu or a gpu. This depends on the value of the macro THRUST_DEVICE_SYSTEM, which can be either THRUST_DEVICE_SYSTEM_OMP for openMP or THRUST_DEVICE_SYSTEM_CUDA for a gpu or THRUST_DEVICE_SYSTEM_CPP for a cpu.
|
|
template<class T > |
using | HMatrix_t = EllSparseBlockMat<T, thrust::host_vector> |
|
using | HMatrix = EllSparseBlockMat<double, thrust::host_vector> |
| Host Matrix for derivatives.
|
|
using | fHMatrix = EllSparseBlockMat<float, thrust::host_vector> |
| Host Matrix for derivatives.
|
|
using | DMatrix = EllSparseBlockMat<double, thrust::device_vector> |
| Device Matrix for derivatives.
|
|
using | fDMatrix = EllSparseBlockMat<float, thrust::device_vector> |
| Device Matrix for derivatives.
|
|
template<class real_type > |
using | IHMatrix_t = dg::SparseMatrix<int, real_type, thrust::host_vector> |
|
template<class real_type > |
using | IDMatrix_t = dg::SparseMatrix<int, real_type, thrust::device_vector> |
|
using | IHMatrix = IHMatrix_t<double> |
|
using | IDMatrix = IDMatrix_t<double> |
|
template<class T > |
using | MHVec_t = dg::MPI_Vector<dg::HVec_t<T> > |
| MPI Host Vector s.a. dg::HVec_t.
|
|
using | MHVec = dg::MPI_Vector<dg::HVec > |
| MPI Host Vector s.a. dg::HVec.
|
|
using | cMHVec = dg::MPI_Vector<dg::cHVec > |
| MPI Host Vector s.a. dg::cHVec.
|
|
using | fMHVec = dg::MPI_Vector<dg::fHVec > |
| MPI Host Vector s.a. dg::fHVec.
|
|
using | MDVec = dg::MPI_Vector<dg::DVec > |
| MPI Device Vector s.a. dg::DVec.
|
|
using | cMDVec = dg::MPI_Vector<dg::cDVec > |
| MPI Device Vector s.a. dg::cDVec.
|
|
using | fMDVec = dg::MPI_Vector<dg::fDVec > |
| MPI Device Vector s.a. dg::fDVec.
|
|
template<class T > |
using | CooMat_t = dg::CooSparseBlockMat<T, thrust::host_vector> |
|
using | CooMat = dg::CooSparseBlockMat<double, thrust::host_vector> |
|
using | fCooMat = dg::CooSparseBlockMat<float, thrust::host_vector> |
|
using | DCooMat = dg::CooSparseBlockMat<double, thrust::device_vector> |
|
using | fDCooMat = dg::CooSparseBlockMat<float, thrust::device_vector> |
|
template<class T > |
using | MHMatrix_t = dg::MPISparseBlockMat<thrust::host_vector, dg::HMatrix_t<T>, dg::CooMat_t<T>> |
| MPI Host Matrix for derivatives.
|
|
using | MHMatrix = dg::MPISparseBlockMat<thrust::host_vector, dg::HMatrix, dg::CooMat> |
| MPI Host Matrix for derivatives.
|
|
using | fMHMatrix = dg::MPISparseBlockMat<thrust::host_vector, dg::fHMatrix, dg::fCooMat> |
| MPI Host Matrix for derivatives.
|
|
using | MDMatrix = dg::MPISparseBlockMat<thrust::device_vector, dg::DMatrix, dg::DCooMat> |
| MPI Device Matrix for derivatives.
|
|
using | fMDMatrix = dg::MPISparseBlockMat<thrust::device_vector, dg::fDMatrix, dg::fDCooMat> |
| MPI Device Matrix for derivatives.
|
|
template<class real_type > |
using | MIHMatrix_t = MPIDistMat< thrust::host_vector, IHMatrix_t<real_type> > |
|
template<class real_type > |
using | MIDMatrix_t = MPIDistMat< thrust::device_vector, IDMatrix_t<real_type> > |
|
using | MIHMatrix = MIHMatrix_t<double> |
|
using | MIDMatrix = MIDMatrix_t<double> |
|
using | aGeometry2d = dg::aRealGeometry2d<double> |
|
using | aGeometry3d = dg::aRealGeometry3d<double> |
|
using | aProductGeometry3d = dg::aRealProductGeometry3d<double> |
|
using | CartesianGrid2d = dg::RealCartesianGrid2d<double> |
|
using | CartesianGrid3d = dg::RealCartesianGrid3d<double> |
|
using | CylindricalGrid3d = dg::RealCylindricalGrid3d<double> |
|
using | CartesianGridX2d = dg::RealCartesianGridX2d<double> |
|
using | CartesianGridX3d = dg::RealCartesianGridX3d<double> |
|
using | aGeometryX2d = dg::aRealGeometryX2d<double> |
|
using | aGeometryX3d = dg::aRealGeometryX3d<double> |
|
using | Grid0d = dg::RealGrid<double,0> |
|
using | Grid1d = dg::RealGrid<double,1> |
|
using | Grid2d = dg::RealGrid<double,2> |
|
using | Grid3d = dg::RealGrid<double,3> |
|
template<size_t Nd> |
using | Grid = dg::RealGrid<double,Nd> |
|
using | aTopology2d = dg::aRealTopology<double,2> |
|
using | aTopology3d = dg::aRealTopology<double,3> |
|
template<class T > |
using | aRealTopology2d = dg::aRealTopology<T,2> |
|
template<class T > |
using | aRealTopology3d = dg::aRealTopology<T,3> |
|
template<class T > |
using | RealGrid0d = dg::RealGrid<T,0> |
|
template<class T > |
using | RealGrid1d = dg::RealGrid<T,1> |
|
template<class T > |
using | RealGrid2d = dg::RealGrid<T,2> |
|
template<class T > |
using | RealGrid3d = dg::RealGrid<T,3> |
|
using | GridX1d = dg::RealGridX1d<double> |
|
using | GridX2d = dg::RealGridX2d<double> |
|
using | GridX3d = dg::RealGridX3d<double> |
|
using | aTopologyX2d = dg::aRealTopologyX2d<double> |
|
using | aTopologyX3d = dg::aRealTopologyX3d<double> |
|
using | aMPIGeometry2d = dg::aRealMPIGeometry2d<double> |
|
using | aMPIGeometry3d = dg::aRealMPIGeometry3d<double> |
|
using | aProductMPIGeometry3d = dg::aRealProductMPIGeometry3d<double> |
|
using | CartesianMPIGrid2d = dg::RealCartesianMPIGrid2d<double> |
|
using | CartesianMPIGrid3d = dg::RealCartesianMPIGrid3d<double> |
|
using | CylindricalMPIGrid3d = dg::RealCylindricalMPIGrid3d<double> |
|
using | MPIGrid0d = dg::RealMPIGrid<double,0> |
|
using | MPIGrid1d = dg::RealMPIGrid<double,1> |
|
using | MPIGrid2d = dg::RealMPIGrid<double,2> |
|
using | MPIGrid3d = dg::RealMPIGrid<double,3> |
|
template<size_t Nd> |
using | MPIGrid = dg::RealMPIGrid<double,Nd> |
|
using | aMPITopology2d = dg::aRealMPITopology<double,2> |
|
using | aMPITopology3d = dg::aRealMPITopology<double,3> |
|
template<class T > |
using | aRealMPITopology2d = dg::aRealMPITopology<T,2> |
|
template<class T > |
using | aRealMPITopology3d = dg::aRealMPITopology<T,3> |
|
template<class T > |
using | RealMPIGrid0d = dg::RealMPIGrid<T,0> |
|
template<class T > |
using | RealMPIGrid1d = dg::RealMPIGrid<T,1> |
|
template<class T > |
using | RealMPIGrid2d = dg::RealMPIGrid<T,2> |
|
template<class T > |
using | RealMPIGrid3d = dg::RealMPIGrid<T,3> |
|
template<class T > |
using | Operator = SquareMatrix<T> |
| The old name for SquareMatrix was Operator.
|
|
using | aRefinement1d = dg::aRealRefinement1d<double> |
|
using | IdentityRefinement = dg::RealIdentityRefinement<double> |
|
using | FemRefinement = dg::RealFemRefinement<double> |
|
using | LinearRefinement = dg::RealLinearRefinement<double> |
|
using | EquidistRefinement = dg::RealEquidistRefinement<double> |
|
using | ExponentialRefinement = dg::RealExponentialRefinement<double> |
|
using | CartesianRefinedGrid2d = dg::RealCartesianRefinedGrid2d<double> |
|
using | CartesianRefinedGrid3d = dg::RealCartesianRefinedGrid3d<double> |
|
using | aRefinementX2d = dg::aRealRefinementX2d<double> |
|
using | IdentityXRefinement = dg::RealIdentityXRefinement<double> |
|
using | EquidistXRefinement = dg::RealEquidistXRefinement<double> |
|
using | ExponentialXRefinement = dg::RealExponentialXRefinement<double> |
|
using | CartesianRefinedGridX2d = dg::RealCartesianRefinedGridX2d<double> |
|
using | CartesianRefinedGridX3d = dg::RealCartesianRefinedGridX3d<double> |
|
template<class MPIContainer > |
using | get_mpi_view_type |
|
template<class ContainerType > |
using | FixedPointIteration = AndersonAcceleration<ContainerType> |
| If you are looking for fixed point iteration: it is a special case of Anderson Acceleration.
|
|
template<class Geometry , class Matrix , class Container > |
using | Elliptic = Elliptic2d<Geometry, Matrix, Container> |
| A 2d negative elliptic differential operator \( -\nabla \cdot ( \mathbf{\chi}\cdot \nabla ) \).
|
|
template<class Geometry , class Matrix , class Container > |
using | Helmholtz = GeneralHelmholtz<dg::Elliptic2d<Geometry,Matrix,Container>, Container> |
| a 2d Helmholtz opereator \( (\chi - \alpha F)\) with \( F = -\Delta\)
|
|
template<class Geometry , class Matrix , class Container > |
using | Helmholtz1d = GeneralHelmholtz<dg::Elliptic1d<Geometry,Matrix,Container>, Container> |
| a 1d Helmholtz opereator \( (\chi - \alpha F)\) with \( F = -\partial_x^2\)
|
|
template<class Geometry , class Matrix , class Container > |
using | Helmholtz2d = GeneralHelmholtz<dg::Elliptic2d<Geometry,Matrix,Container>, Container> |
| a 2d Helmholtz opereator \( (\chi - \alpha F)\) with \( F = -\Delta\)
|
|
template<class Geometry , class Matrix , class Container > |
using | Helmholtz3d = GeneralHelmholtz<dg::Elliptic3d<Geometry,Matrix,Container>, Container> |
| a 3d Helmholtz opereator \( (\chi - \alpha F)\) with \( F = -\Delta\)
|
|
template<class ContainerType > |
using | RungeKutta = ERKStep<ContainerType> |
| Runge-Kutta fixed-step explicit ODE integrator \(
\begin{align}
k_i = f\left( t^n + c_i \Delta t, u^n + \Delta t \sum_{j=1}^{i-1} a_{ij} k_j\right) \\
u^{n+1} = u^{n} + \Delta t\sum_{j=1}^s b_j k_j
\end{align}
\).
|
|
template<class ContainerType > |
using | FilteredRungeKutta = FilteredERKStep<ContainerType> |
| Filtered Runge-Kutta fixed-step explicit ODE integrator \(
\begin{align}
k_i = f\left( t^n + c_i \Delta t, \Lambda\Pi \left[u^n + \Delta t \sum_{j=1}^{i-1} a_{ij} k_j\right]\right) \\
u^{n+1} = \Lambda\Pi\left[u^{n} + \Delta t\sum_{j=1}^s b_j k_j\right]
\end{align}
\).
|
|
template<class ContainerType > |
using | ImplicitRungeKutta = DIRKStep<ContainerType> |
| Runge-Kutta fixed-step implicit ODE integrator \(
\begin{align}
k_i = f\left( t^n + c_i \Delta t, u^n + \Delta t \sum_{j=1}^{s} a_{ij} k_j\right) \\
u^{n+1} = u^{n} + \Delta t\sum_{j=1}^s b_j k_j
\end{align}
\).
|
|
|
template<class ContainerType > |
auto | asDenseMatrix (const std::vector< const ContainerType * > &in) |
| Lightweight DenseMatrix for dg::blas2::gemv .
|
|
template<class ContainerType > |
auto | asDenseMatrix (const std::vector< const ContainerType * > &in, unsigned size) |
| Lightweight DenseMatrix for dg::blas2::gemv .
|
|
template<class ContainerType > |
std::vector< const ContainerType * > | asPointers (const std::vector< ContainerType > &in) |
| Convert a vector of vectors to a vector of pointers.
|
|
void | abort_program (int code=-1) |
| Abort program (both MPI and non-MPI)
|
|
void | mpi_init (int argc, char *argv[]) |
| Convencience shortcut: Calls MPI_Init or MPI_Init_thread and inits CUDA devices.
|
|
template<class T > |
std::vector< T > | mpi_read_as (unsigned num, MPI_Comm comm, std::istream &is=std::cin) |
| Read num values from is and broadcast to all processes as type T .
|
|
void | mpi_read_grid (unsigned &n, std::vector< unsigned > &N, MPI_Comm comm, std::istream &is=std::cin, bool verbose=true, std::ostream &os=std::cout) |
| Read in grid sizes from is .
|
|
void | mpi_read_grid (unsigned &n, std::vector< unsigned * > N, MPI_Comm comm, std::istream &is=std::cin, bool verbose=true, std::ostream &os=std::cout) |
| Convenience shortcut allowing a call like.
|
|
MPI_Comm | mpi_cart_create (MPI_Comm comm_old, std::vector< int > dims, std::vector< int > periods, bool reorder=true) |
| Convenience call to MPI_Cart_create preceded by MPI_Dims_create .
|
|
MPI_Comm | mpi_cart_create (std::vector< dg::bc > bcs, std::istream &is=std::cin, MPI_Comm comm_old=MPI_COMM_WORLD, bool reorder=true, bool verbose=true, std::ostream &os=std::cout) |
| Convenience call: read in number of processses from istream and create Cartesian MPI communicator.
|
|
void | mpi_init1d (dg::bc bcx, MPI_Comm &comm, std::istream &is=std::cin, bool verbose=true) |
| DEPRECATED: Short for.
|
|
void | mpi_init2d (dg::bc bcx, dg::bc bcy, MPI_Comm &comm, std::istream &is=std::cin, bool verbose=true) |
| DEPRECATED: Short for.
|
|
void | mpi_init3d (dg::bc bcx, dg::bc bcy, dg::bc bcz, MPI_Comm &comm, std::istream &is=std::cin, bool verbose=true) |
| DEPRECATED: Short for.
|
|
void | mpi_init1d (dg::bc bcx, unsigned &n, unsigned &N, MPI_Comm &comm, std::istream &is=std::cin, bool verbose=true) |
| DEPRECATED.
|
|
void | mpi_init2d (dg::bc bcx, dg::bc bcy, unsigned &n, unsigned &Nx, unsigned &Ny, MPI_Comm &comm, std::istream &is=std::cin, bool verbose=true) |
| DEPRECATED.
|
|
void | mpi_init3d (dg::bc bcx, dg::bc bcy, dg::bc bcz, unsigned &n, unsigned &Nx, unsigned &Ny, unsigned &Nz, MPI_Comm &comm, std::istream &is=std::cin, bool verbose=true) |
| DEPRECATED.
|
|
MPI_Comm | mpi_cart_sub (MPI_Comm comm, std::vector< int > remain_dims, bool duplicate=false) |
| Call and register a call to MPI_Cart_sub with the dg library.
|
|
MPI_Comm | mpi_cart_kron (std::vector< MPI_Comm > comms) |
| Form a Kronecker product among Cartesian communicators.
|
|
template<class Vector > |
MPI_Comm | mpi_cart_kron (Vector comms) |
| Convenience shortcut for return mpi_cart_kron( std::vector<MPI_Comm>(comms.begin(), comms.end()));
|
|
std::vector< MPI_Comm > | mpi_cart_split (MPI_Comm comm) |
| Split a Cartesian communicator along each dimensions.
|
|
template<size_t Nd> |
std::array< MPI_Comm, Nd > | mpi_cart_split_as (MPI_Comm comm) |
| Same as mpi_cart_split but differen return type.
|
|
template<class real_type , class ConversionPolicyRows , class ConversionPolicyCols > |
auto | make_mpi_sparseblockmat (const EllSparseBlockMat< real_type, thrust::host_vector > &src, const ConversionPolicyRows &g_rows, const ConversionPolicyCols &g_cols) |
| Split given EllSparseBlockMat into computation and communication part.
|
|
template<class MessageType > |
std::map< int, MessageType > | mpi_permute (const std::map< int, MessageType > &messages, MPI_Comm comm) |
| Exchange messages between processes in a communicator.
|
|
template<class ContainerType > |
void | mpi_gather (const thrust::host_vector< std::array< int, 2 > > &gather_map, const ContainerType &gatherFrom, ContainerType &result, MPI_Comm comm) |
| Un-optimized distributed gather operation.
|
|
template<class ContainerType > |
void | mpi_scatter (const thrust::host_vector< std::array< int, 2 > > &scatter_map, const ContainerType &toScatter, ContainerType &result, MPI_Comm comm, bool resize_result=false) |
| Un-optimized distributed scatter operation.
|
|
template<class Integer > |
thrust::host_vector< std::array< Integer, 2 > > | mpi_invert_permutation (const thrust::host_vector< std::array< Integer, 2 > > &p, MPI_Comm comm) |
| Invert a globally bijective index map.
|
|
template<class real_type , class value_type > |
void | ell_cpu_multiply_kernel (value_type alpha, value_type beta, const real_type *RESTRICT data, const int *RESTRICT cols_idx, const int *RESTRICT data_idx, const int num_rows, const int num_cols, const int blocks_per_line, const int n, const int left_size, const int right_size, const int *RESTRICT right_range, const value_type *RESTRICT x, value_type *RESTRICT y) |
|
template<class real_type , class value_type , int n, int blocks_per_line> |
void | ell_cpu_multiply_kernel (value_type alpha, value_type beta, const real_type *RESTRICT data, const int *RESTRICT cols_idx, const int *RESTRICT data_idx, const int num_rows, const int num_cols, const int left_size, const int right_size, const int *RESTRICT right_range, const value_type *RESTRICT x, value_type *RESTRICT y) |
|
template<class real_type , class value_type , int n> |
void | call_ell_cpu_multiply_kernel (value_type alpha, value_type beta, const real_type *RESTRICT data_ptr, const int *RESTRICT cols_ptr, const int *RESTRICT block_ptr, const int num_rows, const int num_cols, const int blocks_per_line, const int left_size, const int right_size, const int *RESTRICT right_range_ptr, const value_type *RESTRICT x_ptr, value_type *RESTRICT y_ptr) |
|
template<class real_type , class value_type , template< class > class Vector> |
void | coo_cpu_multiply_kernel (value_type alpha, const value_type **x, value_type beta, value_type *RESTRICT y, const CooSparseBlockMat< real_type, Vector > &m) |
|
template<class real_type , class value_type , int n, template< class > class Vector> |
void | coo_cpu_multiply_kernel (value_type alpha, const value_type **x, value_type beta, value_type *RESTRICT y, const CooSparseBlockMat< real_type, Vector > &m) |
|
template<class Functor , class Topology > |
auto | evaluate (Functor &&f, const Topology &g) |
| Evaluate a function on grid coordinates
|
|
template<class real_type > |
thrust::host_vector< real_type > | integrate (const thrust::host_vector< real_type > &in, const RealGrid< real_type, 1 > &g, dg::direction dir=dg::forward) |
| Indefinite integral of a function on a grid.
|
|
template<class UnaryOp , class real_type > |
thrust::host_vector< real_type > | integrate (UnaryOp f, const RealGrid< real_type, 1 > &g, dg::direction dir=dg::forward) |
| Untility shortcut.
|
|
template<class Topology > |
Topology::host_vector | forward_transform (const typename Topology::host_vector &in, const Topology &g) |
| Transform a vector from dg::xspace (nodal values) to dg::lspace (modal values)
|
|
template<class T , class ... Ts> |
DG_DEVICE T | zero (T x, Ts ...xs) |
| This enum can be used in dg::evaluate .
|
|
template<class T , class ... Ts> |
DG_DEVICE T | one (T x, Ts ...xs) |
| \( f(x, ...) = 1\)
|
|
DG_DEVICE double | cooX1d (double x) |
| \( f(x) = x\)
|
|
DG_DEVICE double | cooX2d (double x, double y) |
| \( f(x,y) = x\)
|
|
DG_DEVICE double | cooX3d (double x, double y, double z) |
| \( f(x,y,z) = x\)
|
|
DG_DEVICE double | cooY2d (double x, double y) |
| \( f(x,y) = y\)
|
|
DG_DEVICE double | cooY3d (double x, double y, double z) |
| \( f(x,y,z) = y\)
|
|
DG_DEVICE double | cooZ3d (double x, double y, double z) |
| \( f(x,y,z) = z\)
|
|
DG_DEVICE double | cooRZP2X (double R, double Z, double P) |
| \( x = R\sin(\varphi)\)
|
|
DG_DEVICE double | cooRZP2Y (double R, double Z, double P) |
| \( y = R\cos(\varphi)\)
|
|
DG_DEVICE double | cooRZP2Z (double R, double Z, double P) |
| \( z = Z\)
|
|
template<class host_vector , class real_type > |
real_type | interpolate (dg::space sp, const host_vector &v, real_type x, const RealGrid1d< real_type > &g, dg::bc bcx=dg::NEU) |
| Interpolate a vector on a single point on a 1d Grid.
|
|
template<class host_vector , class real_type > |
real_type | interpolate (dg::space sp, const host_vector &v, real_type x, real_type y, const aRealTopology< real_type, 2 > &g, dg::bc bcx=dg::NEU, dg::bc bcy=dg::NEU) |
| Interpolate a vector on a single point on a 2d Grid.
|
|
template<class real_type > |
real_type | interpolate (dg::space sp, const thrust::host_vector< real_type > &v, real_type x, real_type y, const aRealTopologyX2d< real_type > &g, dg::bc bcx=dg::NEU, dg::bc bcy=dg::NEU) |
| Interpolate a vector on a single point on a 2d Grid.
|
|
template<class real_type , class MPITopology > |
MPI_Vector< thrust::host_vector< real_type > > | global2local (const thrust::host_vector< real_type > &global, const MPITopology &g) |
| Take the relevant local part of a global vector.
|
|
template<class ConversionPolicy , class real_type > |
dg::MIHMatrix_t< real_type > | make_mpi_matrix (const dg::IHMatrix_t< real_type > &global_cols, const ConversionPolicy &col_policy) |
| Convert a (row-distributed) matrix with local row and global column indices to a row distributed MPI matrix.
|
|
template<class ConversionPolicy , class real_type > |
dg::IHMatrix_t< real_type > | convertGlobal2LocalRows (const dg::IHMatrix_t< real_type > &global, const ConversionPolicy &row_policy) |
| Convert a (column-distributed) matrix with global row and column indices to a row distributed matrix.
|
|
template<class ConversionPolicy , class real_type > |
void | convertLocal2GlobalCols (dg::IHMatrix_t< real_type > &local, const ConversionPolicy &policy) |
| Convert a matrix with local column indices to a matrix with global column indices.
|
|
template<class T > |
void | lu_solve (const dg::SquareMatrix< T > &lu, const std::vector< unsigned > &p, std::vector< T > &b) |
| Solve the linear system with the LU decomposition.
|
|
template<class T > |
dg::SquareMatrix< T > | invert (const dg::SquareMatrix< T > &in) |
| Compute inverse of square matrix (alias for dg::create::inverse )
|
|
template<class T > |
SquareMatrix< T > | tensorproduct (const SquareMatrix< T > &op1, const SquareMatrix< T > &op2) |
| Form the tensor product between two operators.
|
|
template<class T > |
dg::SparseMatrix< int, T, thrust::host_vector > | tensorproduct (unsigned N, const SquareMatrix< T > &op) |
| Tensor product between Identity matrix and an operator.
|
|
template<class T > |
dg::SparseMatrix< int, T, thrust::host_vector > | sandwich (const SquareMatrix< T > &left, const dg::SparseMatrix< int, T, thrust::host_vector > &m, const SquareMatrix< T > &right) |
| Multiply 1d matrices by diagonal block matrices from left and right.
|
|
template<class T > |
T | gcd (T a, T b) |
| Greatest common divisor.
|
|
template<class T > |
T | lcm (T a, T b) |
| Least common multiple.
|
|
template<class SharedContainer , class real_type > |
void | split (SharedContainer &in, std::vector< View< SharedContainer > > &out, const aRealTopology3d< real_type > &grid) |
| Split a vector into planes along the last dimension (fast version)
|
|
template<class SharedContainer , class real_type > |
std::vector< View< SharedContainer > > | split (SharedContainer &in, const aRealTopology3d< real_type > &grid) |
| Split a vector into planes along the last dimension (construct version)
|
|
template<class Container , class host_vector , class Topology > |
void | assign3dfrom2d (const host_vector &in2d, Container &out, const Topology &grid) |
| Construct a 3d vector given a 2d host vector.
|
|
template<class MPIContainer , class real_type > |
void | split (MPIContainer &in, std::vector< get_mpi_view_type< MPIContainer > > &out, const aRealMPITopology3d< real_type > &grid) |
| MPI Version of split (fast version)
|
|
template<class MPIContainer , class real_type > |
std::vector< get_mpi_view_type< MPIContainer > > | split (MPIContainer &in, const aRealMPITopology3d< real_type > &grid) |
| MPI Version of split (construct version)
|
|
template<class Functor , class Geometry > |
Geometry::host_vector | pullback (const Functor &f, const Geometry &g) |
| \( f_i = f( x(\zeta_i, \eta_i), y(\zeta_i, \eta_i))\)
|
|
template<class Functor1 , class Functor2 , class container , class Geometry > |
void | pushForwardPerp (const Functor1 &vR, const Functor2 &vZ, container &vx, container &vy, const Geometry &g) |
| \( \bar v = J v\)
|
|
template<class Functor1 , class Functor2 , class Functor3 , class container , class Geometry > |
void | pushForward (const Functor1 &vR, const Functor2 &vZ, const Functor3 &vPhi, container &vx, container &vy, container &vz, const Geometry &g) |
| \( {\bar v} = J v\)
|
|
template<class FunctorRR , class FunctorRZ , class FunctorZZ , class container , class Geometry > |
void | pushForwardPerp (const FunctorRR &chiRR, const FunctorRZ &chiRZ, const FunctorZZ &chiZZ, SparseTensor< container > &chi, const Geometry &g) |
| \( \bar \chi = J \chi J^T\)
|
|
template<class T > |
dg::SparseMatrix< int, T, thrust::host_vector > | tensorproduct (const dg::SparseMatrix< int, T, thrust::host_vector > &lhs, const dg::SparseMatrix< int, T, thrust::host_vector > &rhs) |
| \( L\otimes R\) Form the tensor (Kronecker) product between two matrices
|
|
template<class T > |
dg::SparseMatrix< int, T, thrust::host_vector > | tensorproduct_cols (const dg::SparseMatrix< int, T, thrust::host_vector > &lhs, const dg::SparseMatrix< int, T, thrust::host_vector > &rhs) |
| \( L\otimes R\) Form the tensor (Kronecker) product between two matrices in the column index
|
|
template<class from_ContainerType , class ContainerType , class ... Params> |
void | assign (const from_ContainerType &from, ContainerType &to, Params &&... ps) |
| Generic way to assign the contents of a from_ContainerType object to a ContainerType object optionally given additional parameters.
|
|
template<class ContainerType , class from_ContainerType , class ... Params> |
ContainerType | construct (const from_ContainerType &from, Params &&... ps) |
| Generic way to construct an object of ContainerType given a from_ContainerType object and optional additional parameters.
|
|
template<class ContainerType , class Functor , class ... ContainerTypes> |
auto | kronecker (Functor &&f, const ContainerType &x0, const ContainerTypes &... xs) |
| \( y_I = f(x_{0i_0}, x_{1i_1}, ...) \) Memory allocating version of dg::blas1::kronecker
|
|
template<class MatrixType , class ContainerType1 , class ContainerType2 > |
void | apply (get_value_type< ContainerType1 > alpha, MatrixType &&M, const ContainerType1 &x, get_value_type< ContainerType1 > beta, ContainerType2 &y) |
| Alias for dg::blas2::symv \( y = \alpha M(x) + \beta y \);.
|
|
template<class MatrixType , class ContainerType1 , class ContainerType2 > |
void | apply (MatrixType &&M, const ContainerType1 &x, ContainerType2 &y) |
| Alias for dg::blas2::symv \( y = M( x)\);.
|
|
std::string | bc2str (bc bcx) |
| write a string describing boundary condition to an output stream
|
|
bc | str2bc (std::string s) |
| convert a string to a bc
|
|
bc | inverse (bc bound) |
| invert boundary condition
|
|
direction | str2direction (std::string s) |
| convert a string to a direction
|
|
std::string | direction2str (enum direction dir) |
| convert a direciton to string
|
|
direction | inverse (direction dir) |
| invert direction
|
|
template<class ContainerType0 , class ContainerType1 , class ContainerType2 > |
std::vector< double > | least_squares (const std::vector< ContainerType0 > &bs, const ContainerType1 &b, const ContainerType2 &weights) |
| Compute \( a = (B^T W B)^{-1} WB^T b\) for given \( B\), weights \( W\) and right hand side \( b\).
|
|
template<class ContainerType0 , class ContainerType1 > |
std::vector< double > | least_squares (const std::vector< ContainerType0 > &bs, const ContainerType1 &b) |
| An alias for least_squares( bs, b, 1.)
|
|
template<class MatrixType0 , class ContainerType0 , class ContainerType1 , class MatrixType1 , class NestedGrids > |
void | nested_iterations (std::vector< MatrixType0 > &ops, ContainerType0 &x, const ContainerType1 &b, std::vector< MatrixType1 > &inverse_ops, NestedGrids &nested_grids) |
| Full approximation nested iterations.
|
|
template<class NestedGrids , class MatrixType0 , class MatrixType1 , class MatrixType2 > |
void | multigrid_cycle (std::vector< MatrixType0 > &ops, std::vector< MatrixType1 > &inverse_ops_down, std::vector< MatrixType2 > &inverse_ops_up, NestedGrids &nested_grids, unsigned gamma, unsigned p) |
| EXPERIMENTAL Full approximation multigrid cycle.
|
|
template<class MatrixType0 , class MatrixType1 , class MatrixType2 , class NestedGrids , class ContainerType0 , class ContainerType1 > |
void | full_multigrid (std::vector< MatrixType0 > &ops, ContainerType0 &x, const ContainerType1 &b, std::vector< MatrixType1 > &inverse_ops_down, std::vector< MatrixType2 > &inverse_ops_up, NestedGrids &nested_grids, unsigned gamma, unsigned mu) |
| EXPERIMENTAL One Full multigrid cycle.
|
|
template<class NestedGrids , class MatrixType0 , class MatrixType1 , class MatrixType2 , class ContainerType0 , class ContainerType1 , class ContainerType2 > |
void | fmg_solve (std::vector< MatrixType0 > &ops, ContainerType0 &x, const ContainerType1 &b, std::vector< MatrixType1 > &inverse_ops_down, std::vector< MatrixType2 > &inverse_ops_up, NestedGrids &nested_grids, const ContainerType2 &weights, double eps, unsigned gamma) |
| EXPERIMENTAL Full multigrid cycles.
|
|
template<typename UnaryOp > |
int | bisection1d (UnaryOp &&op, double &x_min, double &x_max, const double eps) |
| Find a root of a 1d function \( f(x) = 0\).
|
|
std::string | to2str (enum to mode) |
| Convert integration mode to string.
|
|
bool | is_same (double x, double y, double eps=1e-15) |
|
bool | is_same (float x, float y, float eps=1e-6) |
|
template<class T > |
bool | is_same (T x, T y) |
|
bool | is_divisable (double a, double b, double eps=1e-15) |
|
bool | is_divisable (float a, float b, float eps=1e-6) |
|
template<class UnaryOp , class Functor > |
auto | compose (UnaryOp f, Functor g) |
| Create Composition functor \( f(g(x_0,x_1,...)) \).
|
|
template<class UnaryOp , typename... Functors> |
auto | compose (UnaryOp f0, Functors... fs) |
| Create Composition funtor of an arbitrary number of functions \( f_0(f_1(f_2( ... f_s(x_0, x_1, ...)))\).
|
|