Discontinuous Galerkin Library
#include "dg/algorithm.h"
|
Modules | |
Timer class | |
t.tic(); t.toc(); t.diff(); | |
Functions and functors | |
Lowlevel helper functions and classes | |
Classes | |
class | dg::Message |
small class holding a stringstream More... | |
class | dg::Error |
class intended for the use in throw statements More... | |
struct | dg::Fail |
Indicate failure to converge. More... | |
struct | dg::MultiMatrix< MatrixType, ContainerType > |
Struct that applies given matrices one after the other. More... | |
class | dg::NoRoot1d |
Exception class, that stores boundaries for 1D root finding. More... | |
Functions | |
template<class ContainerType > | |
void | dg::transpose (unsigned nx, unsigned ny, const ContainerType &in, ContainerType &out) |
Transpose vector. More... | |
template<class ContainerType > | |
void | dg::extend_line (unsigned nx, unsigned ny, const ContainerType &in, ContainerType &out) |
Copy a line into rows of output vector. More... | |
template<class ContainerType > | |
void | dg::extend_column (unsigned nx, unsigned ny, const ContainerType &in, ContainerType &out) |
Copy a line into columns of output vector. More... | |
static void | dg::abort_program (int code=-1) |
Abort program (both MPI and non-MPI) More... | |
static void | dg::mpi_init (int argc, char *argv[]) |
Convencience shortcut: Calls MPI_Init or MPI_Init_thread. More... | |
static void | dg::mpi_init2d (dg::bc bcx, dg::bc bcy, MPI_Comm &comm, std::istream &is=std::cin, bool verbose=true) |
Read in number of processses and create Cartesian MPI communicator. More... | |
static void | dg::mpi_init2d (unsigned &n, unsigned &Nx, unsigned &Ny, MPI_Comm comm, std::istream &is=std::cin, bool verbose=true) |
Read in number of processes and broadcast to process group. More... | |
static void | dg::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) |
Read in number of processses and grid size and create Cartesian MPI communicator. More... | |
static void | dg::mpi_init3d (dg::bc bcx, dg::bc bcy, dg::bc bcz, MPI_Comm &comm, std::istream &is=std::cin, bool verbose=true) |
Read in number of processses and create Cartesian MPI communicator. More... | |
static void | dg::mpi_init3d (unsigned &n, unsigned &Nx, unsigned &Ny, unsigned &Nz, MPI_Comm comm, std::istream &is=std::cin, bool verbose=true) |
Read in number of processes and broadcast to process group. More... | |
static void | dg::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) |
Read in number of processses and grid size and create Cartesian MPI communicator. More... | |
template<class real_type > | |
thrust::host_vector< real_type > | dg::forward_transform (const thrust::host_vector< real_type > &in, const aRealTopology2d< real_type > &g) |
Transform a vector from dg::xspace (nodal values) to dg::lspace (modal values) More... | |
template<class UnaryOp , class real_type > | |
dg::Operator< real_type > | dg::create::modal_filter (UnaryOp op, const DLT< real_type > &dlt) |
Create a modal filter block \( V D V^{-1}\). More... | |
template<class T > | |
T | dg::gcd (T a, T b) |
Greatest common divisor. More... | |
template<class T > | |
T | dg::lcm (T a, T b) |
Least common multiple. More... | |
bool | dg::is_same (double x, double y, double eps=1e-15) |
bool | dg::is_same (float x, float y, float eps=1e-6) |
bool | dg::is_divisable (double a, double b, double eps=1e-15) |
bool | dg::is_divisable (float a, float b, float eps=1e-6) |
|
inlinestatic |
Abort program (both MPI and non-MPI)
code | The abortion code that will be signalled to the caller (of the program) |
void dg::extend_column | ( | unsigned | nx, |
unsigned | ny, | ||
const ContainerType & | in, | ||
ContainerType & | out | ||
) |
Copy a line into columns of output vector.
out
[i*nx+j] = in[i]
ContainerType | Any class for which a specialization of TensorTraits exists and which fulfills the requirements of the there defined data and execution policies derived from AnyVectorTag and AnyPolicyTag . Among others
ContainerTypes in the argument list, then TensorTraits must exist for all of them |
nx | number of columns (size of contiguous chunks) in output vector |
ny | size of input vector, number of rows in output vector |
in | input (size ny ) |
out | output (size nx*ny may not alias in) |
void dg::extend_line | ( | unsigned | nx, |
unsigned | ny, | ||
const ContainerType & | in, | ||
ContainerType & | out | ||
) |
Copy a line into rows of output vector.
Just append the line of size nx
until the output vector size of nx*ny
is reached: out
[i*nx+j] = in[j]
ContainerType | Any class for which a specialization of TensorTraits exists and which fulfills the requirements of the there defined data and execution policies derived from AnyVectorTag and AnyPolicyTag . Among others
ContainerTypes in the argument list, then TensorTraits must exist for all of them |
nx | size of the input vector/ number of columns (size of contiguous chunks) in output vector |
ny | number of rows in output vector |
in | input (size nx ) |
out | output (size nx*ny ; may not alias in) |
thrust::host_vector< real_type > dg::forward_transform | ( | const thrust::host_vector< real_type > & | in, |
const aRealTopology2d< real_type > & | g | ||
) |
Transform a vector from dg::xspace (nodal values) to dg::lspace (modal values)
in | input |
g | grid |
T dg::gcd | ( | T | a, |
T | b | ||
) |
Greatest common divisor.
a | First number |
b | Second number |
|
inline |
Checks if two number are integer divisable \(a/b \in \mathbb{Z}\) within accuracy
|
inline |
Checks if two number are integer divisable \(a/b \in \mathbb{Z}\) within accuracy
|
inline |
Checks if two number are same within accuracy
|
inline |
Checks if two number are same within accuracy
T dg::lcm | ( | T | a, |
T | b | ||
) |
Least common multiple.
a | Fist number |
b | Second number |
dg::Operator< real_type > dg::create::modal_filter | ( | UnaryOp | op, |
const DLT< real_type > & | dlt | ||
) |
Create a modal filter block \( V D V^{-1}\).
where \( V\) is the Vandermonde matrix (the backward transformation matrix) and \( D \) is a diagonal matrix with \( D_{ii} = \sigma(i)\)
UnaryOp | Model of Unary Function real_type sigma(unsigned) The input will be the modal number i where \( i=0,...,n-1\) and n is the number of polynomial coefficients in use. The output is the filter strength for the given mode number |
op | the unary function |
dlt | provide the forward, backward transformation and number of coefficients n |
dg::create::fast_transform()
to create a matrix that applies the filter to vectors. For example to create a modal filter that acts in two dimensions:
|
inlinestatic |
Convencience shortcut: Calls MPI_Init or MPI_Init_thread.
Shortcut for
argc | command line argument number |
argv | command line arguments |
|
inlinestatic |
Read in number of processses and create Cartesian MPI communicator.
Also sets the GPU a process should use via rank%
num_devices_per_node if THRUST_DEVICE_SYSTEM==THRUST_DEVICE_SYSTEM_CUDA
bcx | if bcx==dg::PER then the communicator is periodic in x |
bcy | if bcy==dg::PER then the communicator is periodic in y |
comm | (write only) MPI_COMM_WORLD as a 2d Cartesian MPI communicator |
is | Input stream rank 0 reads parameters (npx , npy ) |
verbose | If true, rank 0 prints queries and information on std::cout |
MPI_Comm_free
. (Using mpi_init2d
and mpi_init3d
in the same program has sometimes led to Segmentation faults in the past)
|
inlinestatic |
Read in number of processses and grid size and create Cartesian MPI communicator.
Also sets the GPU a process should use via rank%
num_devices_per_node if THRUST_DEVICE_SYSTEM==THRUST_DEVICE_SYSTEM_CUDA
bcx | if bcx==dg::PER then the communicator is periodic in x |
bcy | if bcy==dg::PER then the communicator is periodic in y |
n | rank 0 reads in from is and broadcasts to all processes in MPI_COMM_WORLD |
Nx | rank 0 reads in from is and broadcasts to all processes in MPI_COMM_WORLD |
Ny | rank 0 reads in from is and broadcasts to all processes in MPI_COMM_WORLD |
comm | (write only) MPI_COMM_WORLD as a 2d Cartesian MPI communicator |
is | Input stream rank 0 reads parameters (npx , npy , n , Nx , Ny ) |
verbose | If true, rank 0 prints queries and information on std::cout |
MPI_Comm_free
. (Using mpi_init2d
and mpi_init3d
in the same program has sometimes led to Segmentation faults in the past)
|
inlinestatic |
Read in number of processes and broadcast to process group.
n | rank 0 reads in from is and broadcasts to all processes in comm |
Nx | rank 0 reads in from is and broadcasts to all processes in comm |
Ny | rank 0 reads in from is and broadcasts to all processes in comm |
comm | (read only) a 2d Cartesian MPI communicator |
is | Input stream rank 0 reads parameters (n , Nx , Ny ) |
verbose | If true, rank 0 prints queries and information on std::cout |
|
inlinestatic |
Read in number of processses and create Cartesian MPI communicator.
Also sets the GPU a process should use via rank%
num_devices_per_node if THRUST_DEVICE_SYSTEM==THRUST_DEVICE_SYSTEM_CUDA
bcx | if bcx==dg::PER then the communicator is periodic in x |
bcy | if bcy==dg::PER then the communicator is periodic in y |
bcz | if bcz==dg::PER then the communicator is periodic in z |
comm | (write only) MPI_COMM_WORLD as a 3d Cartesian MPI communicator |
is | Input stream rank 0 reads parameters (npx , npy , npz ) |
verbose | If true, rank 0 prints queries and information on std::cout |
MPI_Comm_free
. (Using mpi_init2d
and mpi_init3d
in the same program has sometimes led to Segmentation faults in the past)
|
inlinestatic |
Read in number of processses and grid size and create Cartesian MPI communicator.
Also sets the GPU via rank%
num_devices_per_node a process should use if THRUST_DEVICE_SYSTEM==THRUST_DEVICE_SYSTEM_CUDA
bcx | if bcx==dg::PER then the communicator is periodic in x |
bcy | if bcy==dg::PER then the communicator is periodic in y |
bcz | if bcz==dg::PER then the communicator is periodic in z |
n | rank 0 reads in from is and broadcasts to all processes in MPI_COMM_WORLD |
Nx | rank 0 reads in from is and broadcasts to all processes in MPI_COMM_WORLD |
Ny | rank 0 reads in from is and broadcasts to all processes in MPI_COMM_WORLD |
Nz | rank 0 reads in from is and broadcasts to all processes in MPI_COMM_WORLD |
comm | (write only) MPI_COMM_WORLD as a 3d Cartesian MPI communicator |
is | Input stream rank 0 reads parameters (npx , npy , npz , n , Nx , Ny , Nz ) |
verbose | If true, rank 0 prints queries and information on std::cout |
MPI_Comm_free
. (Using mpi_init2d
and mpi_init3d
in the same program has sometimes led to Segmentation faults in the past)
|
inlinestatic |
Read in number of processes and broadcast to process group.
n | rank 0 reads in from is and broadcasts to all processes in comm |
Nx | rank 0 reads in from is and broadcasts to all processes in comm |
Ny | rank 0 reads in from is and broadcasts to all processes in comm |
Nz | rank 0 reads in from is and broadcasts to all processes in comm |
comm | (read only) a 3d Cartesian MPI communicator |
is | Input stream rank 0 reads parameters (n , Nx , Ny , Nz ) |
verbose | If true, rank 0 prints queries and information on std::cout |
void dg::transpose | ( | unsigned | nx, |
unsigned | ny, | ||
const ContainerType & | in, | ||
ContainerType & | out | ||
) |
Transpose vector.
The equivalent of out[j*ny+i] = in[i*nx+j];
ContainerType | Any class for which a specialization of TensorTraits exists and which fulfills the requirements of the there defined data and execution policies derived from AnyVectorTag and AnyPolicyTag . Among others
ContainerTypes in the argument list, then TensorTraits must exist for all of them |
nx | number of columns in input vector (size of contiguous chunks) /rows in output vector |
ny | number of rows in input vector /columns in output vector |
in | input |
out | output (may not alias in) |