LoAdSG
|
#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[]) |
This class provides a Matrix-Vector-Multiplication on Dirichlet-SparseGrids
|
inline |
Constructor for MatrixVectorHomogen.
grid | AdaptiveSparseGrid, make sure that the grid is completed, i.e. fullfills Definition of locally adaptive sparse grid |
mgrid | locally adaptive full subgrids |
|
inline |
Constructor for MatrixVectorHomogen.
grid | AdaptiveSparseGrid, make sure that the grid is completed, i.e. fullfills Definition of locally adaptive sparse grid |
mgrid | locally adaptive full subgrids |
symmetry | make use of symmetry of the bilinear form |
|
inline |
Constructor for MatrixVectorHomogen.
grid | AdaptiveSparseGrid, make sure that the grid is completed, i.e. fullfills Definition of locally adaptive sparse grid |
mgrid | locally adaptive full subgrids |
symmetry | make use of symmetry of the bilinear form |
mpi_on | distribute local stiffness matrices with MPI |
|
private |
Apply a stencil on all grid points with depth less or equal to Depth T
StencilClass |
input | Nodal values |
output | a(input, v_p) |
stencil | |
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.)
StencilClass |
input | Nodal values |
output | a(input, v_p) |
stencil | |
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.)
StencilClass |
input | Nodal values |
output | a(input, v_p) |
stencil | |
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.)
StencilClass |
input | Nodal values |
output | a(input, v_p) |
stencil | |
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.
StencilClass |
input | Nodal values |
output | a(input, v_p) |
stencil | |
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.
StencilClass |
Index | grid point |
input | nodal values |
stencil | |
T |
|
private |
Calculates nodal values of a function given of prewavelets (prew) on a fixed level (depth)
prew | |
depth |
|
private |
Calculates nodal values of a function given of prewavelets (prew) on a fixed level (depth)
prew | |
depth |
|
private |
Calculates one of the 2^d cases.
StencilClass | Stencil shall be given by a class |
prew | Function given in prewavelet basis |
restrictions | 2^d cases of prolongations (=0) and restrictions=1) |
stencil |
|
private |
Calculates one of the 2^d cases.
prew | Function given in prewavelet basis |
restrictions | 2^d cases of prolongations (=0) and restrictions=1) |
stencil |
void MatrixVectorHomogen::casefunction_mpi | ( | VectorSparseG & | prew, |
VectorSparseG & | Ax, | ||
Problem & | stencilClass | ||
) |
Distribute the cases with dynamic mpi_scheduling.
Problem |
prew | |
Ax | |
stencilClass |
|
inlineprivate |
Convert dualspace to prewavlet basis i.e. a(u,nodal_p) -> a(u,prewavelet_p)
Ax_nodal | a(u,nodal_p) |
Ax | a(u,prewavlet_p) |
T | T(p) \leq T \forall p in D |
|
inlineprivate |
Converts an integer to a binary and stores the resulting binary as a boolean.
num | |
binary |
|
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.
num_workers | |
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.
|
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.
|
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.
num_workers | |
total_tasks |
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.
StencilClass |
prew | |
Ax | |
stencilClass |
|
inlineprivate |
Prolongation operator --> prolongation on primal space
values | |
t | depth |
d | direction |
|
inlineprivate |
Restriction operator -> Restriction of dual space
values | |
t | |
d |
|
private |
worker function for dynamic mpi scheduling of the cases
StencilClass |
prew | |
Ax | |
stencilClass |
bool MatrixVectorHomogen::option_mpi = false |
Distribute local stiffness matrices with MPI. (No speedup yet!)
bool MatrixVectorHomogen::option_symmetry = false |
make use of symmetry of the bilinear form
int MatrixVectorHomogen::stencil_count =0 |
counter for stencil evaluations, just needed for testing