Peano
Loading...
Searching...
No Matches
toolbox::finiteelements::BSplinesStencilFactory Class Reference

Stencil Factory. More...

#include <BSplinesStencilFactory.h>

Static Public Member Functions

static tarch::la::Vector< 3, doubleget1DMassStencilP1 ()
 
static tarch::la::Vector< 5, doubleget1DMassStencilP2 ()
 
static tarch::la::Vector< 7, doubleget1DMassStencilP3 ()
 
static tarch::la::Vector< 9, doubleget1DMassStencilP4 ()
 
static tarch::la::Vector< 11, doubleget1DMassStencilP5 ()
 
static tarch::la::Vector< 13, doubleget1DMassStencilP6 ()
 
static tarch::la::Vector< 15, doubleget1DMassStencilP7 ()
 
static tarch::la::Vector< 17, doubleget1DMassStencilP8 ()
 
static tarch::la::Vector< 3, doubleget1DMassStencilP1 (int elementOfSupport)
 A stencil of order p decomposes into 2p stencils among its support.
 
static tarch::la::Vector< 5, doubleget1DMassStencilP2 (int elementOfSupport)
 Spans 1.5 elements in each direction, i.e.
 
static tarch::la::Vector< 7, doubleget1DMassStencilP3 (int elementOfSupport)
 Spans 2 elements in each direction.
 
static tarch::la::Vector< 9, doubleget1DMassStencilP4 (int elementOfSupport)
 Spans 2.5 elements in each direction, i.e.
 
static tarch::la::Vector< 3, doubleget1DLaplaceStencilP1 ()
 Stencil is not scaled at all with any mesh width.
 
static tarch::la::Vector< 5, doubleget1DLaplaceStencilP2 ()
 
static tarch::la::Vector< 7, doubleget1DLaplaceStencilP3 ()
 
static tarch::la::Vector< 9, doubleget1DLaplaceStencilP4 ()
 
static tarch::la::Vector< 11, doubleget1DLaplaceStencilP5 ()
 
static tarch::la::Vector< 13, doubleget1DLaplaceStencilP6 ()
 
static tarch::la::Vector< 15, doubleget1DLaplaceStencilP7 ()
 
static tarch::la::Vector< 17, doubleget1DLaplaceStencilP8 ()
 
static tarch::la::Vector< 3, doubleget1DLaplaceStencilP1 (int elementOfSupport)
 
static tarch::la::Vector< 5, doubleget1DLaplaceStencilP2 (int elementOfSupport)
 
static tarch::la::Vector< 7, doubleget1DLaplaceStencilP3 (int elementOfSupport)
 
static tarch::la::Vector< 9, doubleget1DLaplaceStencilP4 (int elementOfSupport)
 
static tarch::la::Vector< ThreePowerD, doublegetLaplacianStencilP1 (const tarch::la::Vector< Dimensions, double > &scaling=1.0)
 
static tarch::la::Vector< FivePowerD, doublegetLaplacianStencilP2 (const tarch::la::Vector< Dimensions, double > &scaling=1.0)
 
static tarch::la::Vector< SevenPowerD, doublegetLaplacianStencilP3 (const tarch::la::Vector< Dimensions, double > &scaling=1.0)
 
static tarch::la::Vector< NinePowerD, doublegetLaplacianStencilP4 (const tarch::la::Vector< Dimensions, double > &scaling=1.0)
 
static tarch::la::Vector< ThreePowerD, doublegetMassStencilP1 ()
 
static tarch::la::Vector< FivePowerD, doublegetMassStencilP2 ()
 
static tarch::la::Vector< SevenPowerD, doublegetMassStencilP3 ()
 
static tarch::la::Vector< NinePowerD, doublegetMassStencilP4 ()
 
template<int p, int N>
static tarch::la::Matrix< N *N, N *N, doublegetElementWiseAssemblyMatrix (tarch::la::Vector< 2 *p+1, double > stencilX[N], tarch::la::Vector< 2 *p+1, double > stencilY[N])
 N means NumberOfInfluencingVerticesAroundCurrentElementAlongEachAxis.
 
template<int p, int N>
static tarch::la::Matrix< N *N *N, N *N *N, doublegetElementWiseAssemblyMatrix (tarch::la::Vector< 2 *p+1, double > stencilX[N], tarch::la::Vector< 2 *p+1, double > stencilY[N], tarch::la::Vector< 2 *p+1, double > stencilZ[N])
 

Detailed Description

Stencil Factory.

Many components in the Peano framework switch from a stencil notation to element-wise assembly matrices. This class provides some operations to construct the stencils. To convert these stencil into element-wise assembly matrices, use the class ElementMatrix.

Author
Tobias Weinzierl

Definition at line 28 of file BSplinesStencilFactory.h.

Member Function Documentation

◆ get1DLaplaceStencilP1() [1/2]

tarch::la::Vector< 3, double > toolbox::finiteelements::BSplinesStencilFactory::get1DLaplaceStencilP1 ( )
static

Stencil is not scaled at all with any mesh width.

Returns
\( [-1, 2, -1] \)

Definition at line 61 of file BSplinesStencilFactory.cpp.

◆ get1DLaplaceStencilP1() [2/2]

tarch::la::Vector< 3, double > toolbox::finiteelements::BSplinesStencilFactory::get1DLaplaceStencilP1 ( int elementOfSupport)
static

Definition at line 203 of file BSplinesStencilFactory.cpp.

References assertion, and assertion1.

◆ get1DLaplaceStencilP2() [1/2]

tarch::la::Vector< 5, double > toolbox::finiteelements::BSplinesStencilFactory::get1DLaplaceStencilP2 ( )
static

Definition at line 68 of file BSplinesStencilFactory.cpp.

◆ get1DLaplaceStencilP2() [2/2]

tarch::la::Vector< 5, double > toolbox::finiteelements::BSplinesStencilFactory::get1DLaplaceStencilP2 ( int elementOfSupport)
static

Definition at line 218 of file BSplinesStencilFactory.cpp.

References assertion, and assertion1.

◆ get1DLaplaceStencilP3() [1/2]

tarch::la::Vector< 7, double > toolbox::finiteelements::BSplinesStencilFactory::get1DLaplaceStencilP3 ( )
static

Definition at line 75 of file BSplinesStencilFactory.cpp.

◆ get1DLaplaceStencilP3() [2/2]

tarch::la::Vector< 7, double > toolbox::finiteelements::BSplinesStencilFactory::get1DLaplaceStencilP3 ( int elementOfSupport)
static

Definition at line 237 of file BSplinesStencilFactory.cpp.

References assertion, and assertion1.

◆ get1DLaplaceStencilP4() [1/2]

tarch::la::Vector< 9, double > toolbox::finiteelements::BSplinesStencilFactory::get1DLaplaceStencilP4 ( )
static

Definition at line 82 of file BSplinesStencilFactory.cpp.

◆ get1DLaplaceStencilP4() [2/2]

tarch::la::Vector< 9, double > toolbox::finiteelements::BSplinesStencilFactory::get1DLaplaceStencilP4 ( int elementOfSupport)
static

Definition at line 260 of file BSplinesStencilFactory.cpp.

References assertion, and assertion1.

◆ get1DLaplaceStencilP5()

tarch::la::Vector< 11, double > toolbox::finiteelements::BSplinesStencilFactory::get1DLaplaceStencilP5 ( )
static

Definition at line 89 of file BSplinesStencilFactory.cpp.

◆ get1DLaplaceStencilP6()

tarch::la::Vector< 13, double > toolbox::finiteelements::BSplinesStencilFactory::get1DLaplaceStencilP6 ( )
static

Definition at line 96 of file BSplinesStencilFactory.cpp.

◆ get1DLaplaceStencilP7()

tarch::la::Vector< 15, double > toolbox::finiteelements::BSplinesStencilFactory::get1DLaplaceStencilP7 ( )
static

Definition at line 103 of file BSplinesStencilFactory.cpp.

◆ get1DLaplaceStencilP8()

tarch::la::Vector< 17, double > toolbox::finiteelements::BSplinesStencilFactory::get1DLaplaceStencilP8 ( )
static

Definition at line 110 of file BSplinesStencilFactory.cpp.

◆ get1DMassStencilP1() [1/2]

tarch::la::Vector< 3, double > toolbox::finiteelements::BSplinesStencilFactory::get1DMassStencilP1 ( )
static
Returns
\( [\frac{1}{6}, \frac{2}{3}, \frac{1}{6}] \)

Definition at line 4 of file BSplinesStencilFactory.cpp.

◆ get1DMassStencilP1() [2/2]

tarch::la::Vector< 3, double > toolbox::finiteelements::BSplinesStencilFactory::get1DMassStencilP1 ( int elementOfSupport)
static

A stencil of order p decomposes into 2p stencils among its support.

With the present operation, you can access individual substencils.

Definition at line 117 of file BSplinesStencilFactory.cpp.

References assertion, and assertion1.

◆ get1DMassStencilP2() [1/2]

tarch::la::Vector< 5, double > toolbox::finiteelements::BSplinesStencilFactory::get1DMassStencilP2 ( )
static

Definition at line 11 of file BSplinesStencilFactory.cpp.

◆ get1DMassStencilP2() [2/2]

tarch::la::Vector< 5, double > toolbox::finiteelements::BSplinesStencilFactory::get1DMassStencilP2 ( int elementOfSupport)
static

Spans 1.5 elements in each direction, i.e.

4 in total.

Definition at line 134 of file BSplinesStencilFactory.cpp.

References assertion, and assertion1.

◆ get1DMassStencilP3() [1/2]

tarch::la::Vector< 7, double > toolbox::finiteelements::BSplinesStencilFactory::get1DMassStencilP3 ( )
static

Definition at line 18 of file BSplinesStencilFactory.cpp.

◆ get1DMassStencilP3() [2/2]

tarch::la::Vector< 7, double > toolbox::finiteelements::BSplinesStencilFactory::get1DMassStencilP3 ( int elementOfSupport)
static

Spans 2 elements in each direction.

Definition at line 157 of file BSplinesStencilFactory.cpp.

References assertion, and assertion1.

◆ get1DMassStencilP4() [1/2]

tarch::la::Vector< 9, double > toolbox::finiteelements::BSplinesStencilFactory::get1DMassStencilP4 ( )
static

Definition at line 26 of file BSplinesStencilFactory.cpp.

◆ get1DMassStencilP4() [2/2]

tarch::la::Vector< 9, double > toolbox::finiteelements::BSplinesStencilFactory::get1DMassStencilP4 ( int elementOfSupport)
static

Spans 2.5 elements in each direction, i.e.

6 in total.

Definition at line 180 of file BSplinesStencilFactory.cpp.

References assertion, and assertion1.

◆ get1DMassStencilP5()

tarch::la::Vector< 11, double > toolbox::finiteelements::BSplinesStencilFactory::get1DMassStencilP5 ( )
static

Definition at line 33 of file BSplinesStencilFactory.cpp.

◆ get1DMassStencilP6()

tarch::la::Vector< 13, double > toolbox::finiteelements::BSplinesStencilFactory::get1DMassStencilP6 ( )
static

Definition at line 40 of file BSplinesStencilFactory.cpp.

◆ get1DMassStencilP7()

tarch::la::Vector< 15, double > toolbox::finiteelements::BSplinesStencilFactory::get1DMassStencilP7 ( )
static

Definition at line 47 of file BSplinesStencilFactory.cpp.

◆ get1DMassStencilP8()

tarch::la::Vector< 17, double > toolbox::finiteelements::BSplinesStencilFactory::get1DMassStencilP8 ( )
static

Definition at line 54 of file BSplinesStencilFactory.cpp.

◆ getElementWiseAssemblyMatrix() [1/2]

template<int p, int N>
tarch::la::Matrix< N *N, N *N, double > toolbox::finiteelements::BSplinesStencilFactory::getElementWiseAssemblyMatrix ( tarch::la::Vector< 2 *p+1, double > stencilX[N],
tarch::la::Vector< 2 *p+1, double > stencilY[N] )
static

N means NumberOfInfluencingVerticesAroundCurrentElementAlongEachAxis.

p | N

1 | 2 2 | 4 3 | 4 4 | 6

Definition at line 6 of file BSplinesStencilFactory.cpph.

References tarch::la::allGreaterEquals(), tarch::la::allSmaller(), d2for, peano4::utils::d2Linearised(), enddforx, and toolbox::finiteelements::stencilProduct().

Here is the call graph for this function:

◆ getElementWiseAssemblyMatrix() [2/2]

template<int p, int N>
tarch::la::Matrix< N *N *N, N *N *N, double > toolbox::finiteelements::BSplinesStencilFactory::getElementWiseAssemblyMatrix ( tarch::la::Vector< 2 *p+1, double > stencilX[N],
tarch::la::Vector< 2 *p+1, double > stencilY[N],
tarch::la::Vector< 2 *p+1, double > stencilZ[N] )
static

◆ getLaplacianStencilP1()

tarch::la::Vector< ThreePowerD, double > toolbox::finiteelements::BSplinesStencilFactory::getLaplacianStencilP1 ( const tarch::la::Vector< Dimensions, double > & scaling = 1.0)
static

Definition at line 411 of file BSplinesStencilFactory.cpp.

References assertionMsg, and toolbox::finiteelements::stencilProduct().

Here is the call graph for this function:

◆ getLaplacianStencilP2()

tarch::la::Vector< FivePowerD, double > toolbox::finiteelements::BSplinesStencilFactory::getLaplacianStencilP2 ( const tarch::la::Vector< Dimensions, double > & scaling = 1.0)
static

Definition at line 288 of file BSplinesStencilFactory.cpp.

References assertionMsg, and toolbox::finiteelements::stencilProduct().

Here is the call graph for this function:

◆ getLaplacianStencilP3()

tarch::la::Vector< SevenPowerD, double > toolbox::finiteelements::BSplinesStencilFactory::getLaplacianStencilP3 ( const tarch::la::Vector< Dimensions, double > & scaling = 1.0)
static

Definition at line 329 of file BSplinesStencilFactory.cpp.

References assertionMsg, and toolbox::finiteelements::stencilProduct().

Here is the call graph for this function:

◆ getLaplacianStencilP4()

tarch::la::Vector< NinePowerD, double > toolbox::finiteelements::BSplinesStencilFactory::getLaplacianStencilP4 ( const tarch::la::Vector< Dimensions, double > & scaling = 1.0)
static

Definition at line 370 of file BSplinesStencilFactory.cpp.

References assertionMsg, and toolbox::finiteelements::stencilProduct().

Here is the call graph for this function:

◆ getMassStencilP1()

tarch::la::Vector< ThreePowerD, double > toolbox::finiteelements::BSplinesStencilFactory::getMassStencilP1 ( )
static

Definition at line 452 of file BSplinesStencilFactory.cpp.

References assertionMsg, and toolbox::finiteelements::stencilProduct().

Here is the call graph for this function:

◆ getMassStencilP2()

tarch::la::Vector< FivePowerD, double > toolbox::finiteelements::BSplinesStencilFactory::getMassStencilP2 ( )
static

Definition at line 474 of file BSplinesStencilFactory.cpp.

References assertionMsg, and toolbox::finiteelements::stencilProduct().

Here is the call graph for this function:

◆ getMassStencilP3()

tarch::la::Vector< SevenPowerD, double > toolbox::finiteelements::BSplinesStencilFactory::getMassStencilP3 ( )
static

Definition at line 496 of file BSplinesStencilFactory.cpp.

References assertionMsg, and toolbox::finiteelements::stencilProduct().

Here is the call graph for this function:

◆ getMassStencilP4()

tarch::la::Vector< NinePowerD, double > toolbox::finiteelements::BSplinesStencilFactory::getMassStencilP4 ( )
static

Definition at line 518 of file BSplinesStencilFactory.cpp.

References assertionMsg, and toolbox::finiteelements::stencilProduct().

Here is the call graph for this function:

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