LoAdSG
MatrixVectorHomogen Class Reference

#include <MatrixVectorHomogen.h>

Public Member Functions

 MatrixVectorHomogen (AdaptiveSparseGrid &grid, MultiLevelAdaptiveSparseGrid &mgrid)
 Constructor for MatrixVectorHomogen.
 
 MatrixVectorHomogen (AdaptiveSparseGrid &grid, MultiLevelAdaptiveSparseGrid &mgrid, bool symmetry)
 Constructor for MatrixVectorHomogen.
 
 MatrixVectorHomogen (AdaptiveSparseGrid &grid, MultiLevelAdaptiveSparseGrid &mgrid, bool symmetry, bool mpi_on)
 Constructor for MatrixVectorHomogen.
 
template<class StencilClass >
void multiplication (VectorSparseG &prew, VectorSparseG &Ax, StencilClass &stencilClass)
 
template<class Problem >
void casefunction_mpi (VectorSparseG &prew, VectorSparseG &Ax, Problem &stencilClass)
 Distribute the cases with dynamic mpi_scheduling.
 

Public Attributes

bool option_symmetry = false
 
bool option_mpi = false
 
int stencil_count =0
 

Private Member Functions

void calcNodal (VectorSparseG &prew, Depth &depth)
 
void calcNodal (VectorSparseG &nodal_u, VectorSparseG &prew, Depth &depth)
 
template<class StencilClass >
void caseFunction (VectorSparseG &prew, bool *restrictions, StencilClass &stencil)
 
template<class Problem >
void caseFunction_LocalStiffnessMatrix (VectorSparseG &prew, bool *restrictions, Problem &localStiffnessMatrixAllDepths)
 
template<class StencilClass >
void applyStencilGlobal (const VectorSparseG &input, VectorSparseG &output, StencilClass &stencil, Depth &T)
 
template<class StencilClass >
void applyStencilGlobal_symmetry (VectorSparseG &input, VectorSparseG &output, StencilClass &stencil, Depth &T)
 
template<class StencilClass >
void applyStencilGlobal_MPI (VectorSparseG &input, VectorSparseG &output, StencilClass &stencil, Depth &T)
 
template<class StencilClass >
void applyStencilGlobal_MPI_2 (VectorSparseG &input, VectorSparseG &output, StencilClass &stencil, Depth &T)
 
template<class StencilClass >
void applyStencilGlobal_MPI_3 (VectorSparseG &input, VectorSparseG &output, StencilClass &stencil, Depth &T)
 
template<class StencilClass >
double applyStencilLocal (const IndexDimension &Index, const VectorSparseG &input, StencilClass &stencil, Depth &T)
 
void prolongation1D_inplace (VectorSparseG &values, int &t, int &d)
 
void restriction1D_inverted_inplace (const VectorSparseG &values, int t, int d)
 
void ConvertToPrewavelet2 (VectorSparseG &Ax_nodal, VectorSparseG &Ax, Depth &T)
 
void master (int num_workers, int total_tasks)
 
void master_start (int num_workers)
 
void master_distribute ()
 
void master_end_distribute ()
 
template<class StencilClass >
void worker (VectorSparseG &prew, VectorSparseG &Ax, StencilClass &stencilClass)
 
void intToBinary (int num, bool binary[])
 

Detailed Description

This class provides a Matrix-Vector-Multiplication on Dirichlet-SparseGrids

Constructor & Destructor Documentation

◆ MatrixVectorHomogen() [1/3]

MatrixVectorHomogen::MatrixVectorHomogen ( AdaptiveSparseGrid grid,
MultiLevelAdaptiveSparseGrid &  mgrid 
)
inline

Constructor for MatrixVectorHomogen.

Parameters
gridAdaptiveSparseGrid, make sure that the grid is completed, i.e. fullfills Definition of locally adaptive sparse grid
mgridlocally adaptive full subgrids

◆ MatrixVectorHomogen() [2/3]

MatrixVectorHomogen::MatrixVectorHomogen ( AdaptiveSparseGrid grid,
MultiLevelAdaptiveSparseGrid &  mgrid,
bool  symmetry 
)
inline

Constructor for MatrixVectorHomogen.

Parameters
gridAdaptiveSparseGrid, make sure that the grid is completed, i.e. fullfills Definition of locally adaptive sparse grid
mgridlocally adaptive full subgrids
symmetrymake use of symmetry of the bilinear form

◆ MatrixVectorHomogen() [3/3]

MatrixVectorHomogen::MatrixVectorHomogen ( AdaptiveSparseGrid grid,
MultiLevelAdaptiveSparseGrid &  mgrid,
bool  symmetry,
bool  mpi_on 
)
inline

Constructor for MatrixVectorHomogen.

Parameters
gridAdaptiveSparseGrid, make sure that the grid is completed, i.e. fullfills Definition of locally adaptive sparse grid
mgridlocally adaptive full subgrids
symmetrymake use of symmetry of the bilinear form
mpi_ondistribute local stiffness matrices with MPI

Member Function Documentation

◆ applyStencilGlobal()

template<class StencilClass >
void MatrixVectorHomogen::applyStencilGlobal ( const VectorSparseG &  input,
VectorSparseG &  output,
StencilClass &  stencil,
Depth &  T 
)
private

Apply a stencil on all grid points with depth less or equal to Depth T

Template Parameters
StencilClass
Parameters
inputNodal values
outputa(input, v_p)
stencil
T
Here is the caller graph for this function:

◆ applyStencilGlobal_MPI()

template<class StencilClass >
void MatrixVectorHomogen::applyStencilGlobal_MPI ( VectorSparseG &  input,
VectorSparseG &  output,
StencilClass &  stencil,
Depth &  T 
)
private

Apply a stencil on all grid points with depth less or equal to Depth T. Make use of Symmetry of the bilinearform by calculating local stiffness matrices on cells. Distribute the local stiffness matrices with MPI. (No SpeedUp yet.)

Template Parameters
StencilClass
Parameters
inputNodal values
outputa(input, v_p)
stencil
T

◆ applyStencilGlobal_MPI_2()

template<class StencilClass >
void MatrixVectorHomogen::applyStencilGlobal_MPI_2 ( VectorSparseG &  input,
VectorSparseG &  output,
StencilClass &  stencil,
Depth &  T 
)
private

Apply a stencil on all grid points with depth less or equal to Depth T. Make use of Symmetry of the bilinearform by calculating local stiffness matrices on cells. Distribute the local stiffness matrices with MPI. (No SpeedUp yet.)

Template Parameters
StencilClass
Parameters
inputNodal values
outputa(input, v_p)
stencil
T

◆ applyStencilGlobal_MPI_3()

template<class StencilClass >
void MatrixVectorHomogen::applyStencilGlobal_MPI_3 ( VectorSparseG &  input,
VectorSparseG &  output,
StencilClass &  stencil,
Depth &  T 
)
private

Apply a stencil on all grid points with depth less or equal to Depth T. Make use of Symmetry of the bilinearform by calculating local stiffness matrices on cells. Distribute the local stiffness matrices with MPI. (No SpeedUp yet.)

Template Parameters
StencilClass
Parameters
inputNodal values
outputa(input, v_p)
stencil
T
Here is the caller graph for this function:

◆ applyStencilGlobal_symmetry()

template<class StencilClass >
void MatrixVectorHomogen::applyStencilGlobal_symmetry ( VectorSparseG &  input,
VectorSparseG &  output,
StencilClass &  stencil,
Depth &  T 
)
private

Apply a stencil on all grid points with depth less or equal to Depth T. Make use of Symmetry of the bilinearform by calculating local stiffness matrices on cells.

Template Parameters
StencilClass
Parameters
inputNodal values
outputa(input, v_p)
stencil
T
Here is the caller graph for this function:

◆ applyStencilLocal()

template<class StencilClass >
double MatrixVectorHomogen::applyStencilLocal ( const IndexDimension &  Index,
const VectorSparseG &  input,
StencilClass &  stencil,
Depth &  T 
)
private

Apply a stencil on a given grid point (Index) w.r.t. Depth T. I.e. a(u,v_p) with u = nodal, p = Index.

Template Parameters
StencilClass
Parameters
Indexgrid point
inputnodal values
stencil
T
Returns

◆ calcNodal() [1/2]

void MatrixVectorHomogen::calcNodal ( VectorSparseG &  nodal_u,
VectorSparseG &  prew,
Depth &  depth 
)
private

Calculates nodal values of a function given of prewavelets (prew) on a fixed level (depth)

Parameters
prew
depth

◆ calcNodal() [2/2]

void MatrixVectorHomogen::calcNodal ( VectorSparseG &  prew,
Depth &  depth 
)
private

Calculates nodal values of a function given of prewavelets (prew) on a fixed level (depth)

Parameters
prew
depth

◆ caseFunction()

template<class StencilClass >
void MatrixVectorHomogen::caseFunction ( VectorSparseG &  prew,
bool *  restrictions,
StencilClass &  stencil 
)
private

Calculates one of the 2^d cases.

Template Parameters
StencilClassStencil shall be given by a class
Parameters
prewFunction given in prewavelet basis
restrictions2^d cases of prolongations (=0) and restrictions=1)
stencil
Here is the caller graph for this function:

◆ caseFunction_LocalStiffnessMatrix()

template<class Problem >
void MatrixVectorHomogen::caseFunction_LocalStiffnessMatrix ( VectorSparseG &  prew,
bool *  restrictions,
Problem &  localStiffnessMatrixAllDepths 
)
private

Calculates one of the 2^d cases.

Parameters
prewFunction given in prewavelet basis
restrictions2^d cases of prolongations (=0) and restrictions=1)
stencil

◆ casefunction_mpi()

template<class Problem >
void MatrixVectorHomogen::casefunction_mpi ( VectorSparseG &  prew,
VectorSparseG &  Ax,
Problem &  stencilClass 
)

Distribute the cases with dynamic mpi_scheduling.

Template Parameters
Problem
Parameters
prew
Ax
stencilClass

◆ ConvertToPrewavelet2()

void MatrixVectorHomogen::ConvertToPrewavelet2 ( VectorSparseG &  Ax_nodal,
VectorSparseG &  Ax,
Depth &  T 
)
inlineprivate

Convert dualspace to prewavlet basis i.e. a(u,nodal_p) -> a(u,prewavelet_p)

Parameters
Ax_nodala(u,nodal_p)
Axa(u,prewavlet_p)
TT(p) \leq T \forall p in D
Here is the caller graph for this function:

◆ intToBinary()

void MatrixVectorHomogen::intToBinary ( int  num,
bool  binary[] 
)
inlineprivate

Converts an integer to a binary and stores the resulting binary as a boolean.

Parameters
num
binary
Here is the caller graph for this function:

◆ master()

void MatrixVectorHomogen::master ( int  num_workers,
int  total_tasks 
)
private

This function is used for dynamic scheduling of the 2^d different cases of prolongations and restrictions. The master sends a case (bool restrictions) to the worker. The worker then starts calculation the required partial sum of the bilinear form.

Parameters
num_workers
total_tasks

◆ master_distribute()

void MatrixVectorHomogen::master_distribute ( )
private

This function is used for dynamic scheduling of the 2^d different cases of prolongations and restrictions. The master sends a case (bool restrictions) to the worker. The worker then starts calculation the required partial sum of the bilinear form.

◆ master_end_distribute()

void MatrixVectorHomogen::master_end_distribute ( )
private

This function is used for dynamic scheduling of the 2^d different cases of prolongations and restrictions. The master sends a case (bool restrictions) to the worker. The worker then starts calculation the required partial sum of the bilinear form.

Here is the caller graph for this function:

◆ master_start()

void MatrixVectorHomogen::master_start ( int  num_workers)
private

This function is used for dynamic scheduling of the 2^d different cases of prolongations and restrictions. The master sends a case (bool restrictions) to the worker. The worker then starts calculation the required partial sum of the bilinear form.

Parameters
num_workers
total_tasks
Here is the caller graph for this function:

◆ multiplication()

template<class StencilClass >
void MatrixVectorHomogen::multiplication ( VectorSparseG &  prew,
VectorSparseG &  Ax,
StencilClass &  stencilClass 
)

Calculates \( a(prew, \varphi_p) \quad \forall p \in D \) and stores it in Ax, with prew given in prewavelet basis.

Template Parameters
StencilClass
Parameters
prew
Ax
stencilClass

◆ prolongation1D_inplace()

void MatrixVectorHomogen::prolongation1D_inplace ( VectorSparseG &  values,
int &  t,
int &  d 
)
inlineprivate

Prolongation operator --> prolongation on primal space

Parameters
values
tdepth
ddirection
Here is the caller graph for this function:

◆ restriction1D_inverted_inplace()

void MatrixVectorHomogen::restriction1D_inverted_inplace ( const VectorSparseG &  values,
int  t,
int  d 
)
inlineprivate

Restriction operator -> Restriction of dual space

Parameters
values
t
d
Here is the caller graph for this function:

◆ worker()

template<class StencilClass >
void MatrixVectorHomogen::worker ( VectorSparseG &  prew,
VectorSparseG &  Ax,
StencilClass &  stencilClass 
)
private

worker function for dynamic mpi scheduling of the cases

Template Parameters
StencilClass
Parameters
prew
Ax
stencilClass
Here is the caller graph for this function:

Member Data Documentation

◆ option_mpi

bool MatrixVectorHomogen::option_mpi = false

Distribute local stiffness matrices with MPI. (No speedup yet!)

◆ option_symmetry

bool MatrixVectorHomogen::option_symmetry = false

make use of symmetry of the bilinear form

◆ stencil_count

int MatrixVectorHomogen::stencil_count =0

counter for stencil evaluations, just needed for testing


The documentation for this class was generated from the following files: