Peano
Loading...
Searching...
No Matches
ElementMatrix.h
Go to the documentation of this file.
1// This file is part of the Peano project. For conditions of distribution and
2// use, please see the copyright notice at www.peano-framework.org
3#pragma once
4
5
6#include "Stencil.h"
7#include <bitset>
8
9
10namespace toolbox {
71 namespace finiteelements {
96
104
105
123
133
169
174
217 ElementWiseAssemblyMatrix getElementWiseAssemblyMatrix( const VectorOfStencils& vectorOfStencils, const std::bitset<ThreePowerD>& cellIsInside );
218
222 ComplexElementWiseAssemblyMatrix getElementWiseAssemblyMatrix( const VectorOfComplexStencils& vectorOfStencils, const std::bitset<ThreePowerD>& cellIsInside );
223
238 int mapElementMatrixEntryOntoStencilEntry(int row, int col);
239
241
242 double getDiagonalElement( const ElementWiseAssemblyMatrix& matrix );
243 double getDiagonalElement( const Stencil& stencil );
244
245 template<int StencilSize>
249 );
250
251 template<int StencilSize>
256 );
257
262 const ElementWiseAssemblyMatrix& matrix,
264 double weight = 1.0
265 );
267 const ElementWiseAssemblyMatrix& matrix,
269 double weight = 1.0
270 );
271 }
272}
273
274#include "ElementMatrix.cpph"
275
const float const float const float struct part *restrict struct part *restrict const float a
Definition hydro_iact.h:55
tarch::la::Vector< ThreePowerD, std::complex< double > > ComplexStencil
Definition Stencil.h:30
tarch::la::Matrix< TwoPowerD, TwoPowerD, double > ElementWiseAssemblyMatrix
Element-wise Matrix.
Definition Stencil.h:18
tarch::la::Vector< ThreePowerD *TwoPowerD, std::complex< double > > VectorOfComplexStencils
Definition Stencil.h:36
double getDiagonalElement(const ElementWiseAssemblyMatrix &matrix)
ElementWiseAssemblyMatrix inverseHierarchicalTransform(const ElementWiseAssemblyMatrix &matrix, const tarch::la::Vector< Dimensions, double > &h, double weight=1.0)
tarch::la::Vector< ThreePowerD, double > Stencil
Stencil.
Definition Stencil.h:29
int mapElementMatrixEntryOntoStencilEntry(int row, int col)
An element-wise assembly matrix is distilled from a sequence of stencils which are in turn double a...
ElementWiseAssemblyMatrix hierarchicalTransform(const ElementWiseAssemblyMatrix &matrix, const tarch::la::Vector< Dimensions, double > &h, double weight=1.0)
Stencil reconstructUniformStencilFragments(const ElementWiseAssemblyMatrix &matrix)
Reconstruct stencils from element matrices.
tarch::la::Vector< ThreePowerD *TwoPowerD, double > VectorOfStencils
Vectors.
Definition Stencil.h:35
tarch::la::Matrix< TwoPowerD, TwoPowerD, std::complex< double > > ComplexElementWiseAssemblyMatrix
Definition Stencil.h:20
tarch::la::Vector< StencilSize *StencilSize, double > stencilProduct(const tarch::la::Vector< StencilSize, double > &a, const tarch::la::Vector< StencilSize, double > &b)
ElementWiseAssemblyMatrix getElementWiseAssemblyMatrix(const Stencil &stencil)
Derive element-wise matrix from stencil.
tarch::la::Vector< TwoPowerD *ThreePowerD, double > reconstructStencilFragments(const ElementWiseAssemblyMatrix &matrix)
Reconstruct stencils from element matrices.
Simple vector class.
Definition Vector.h:150