Peano
Loading...
Searching...
No Matches
StencilFactory.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
8#include <bitset>
9#include <functional>
10
11
12namespace toolbox {
23 namespace finiteelements {
57 void preprocessBoundaryStencil( Stencil& stencil, const std::bitset<Dimensions*2>& boundaryFaceNormals );
58
65 Stencil exchangeCoordinates( const Stencil& stencil, int coord0, int coord1 );
66
71
72
79
86
91
96
103
115 );
119 );
120
128 );
133 );
134
140 );
146 );
147
154 );
161 );
162
173 const tarch::la::Vector<Dimensions,double>& scaling = 1.0
174 );
176 const tarch::la::Vector<Dimensions,std::complex<double> >& h,
177 const tarch::la::Vector<Dimensions,std::complex<double> >& scaling = std::complex<double>(1.0)
178 );
181 const tarch::la::Vector<Dimensions,std::complex<double> >& scaling
182 );
183
197 );
198
245 void applyNeumannBC(
247 double& rhs,
248 const std::bitset<Dimensions*2>& boundaryFaceNormals,
250 double derivative
251 );
252
261 tarch::la::Vector<ThreePowerD,std::complex<double> > getMassMatrix(const tarch::la::Vector<Dimensions,std::complex<double> >& h);
262
263
271 const std::complex<double>& phi
272 );
273
275 const tarch::la::Vector<Dimensions,std::complex<double> >& h,
276 const std::complex<double>& phi
277 );
278
279
289 tarch::la::Vector<ThreePowerD,std::complex<double> > getIdentity(const tarch::la::Vector<Dimensions,std::complex<double> >& h);
290
295
309 Stencil extractElementStencil( const Stencil& stencil, const tarch::la::Vector<Dimensions,int>& adjacentCell );
310
311
324 const tarch::la::Vector<Dimensions,double>& cellCentre,
326 int integrationPointsPerAxis,
327 std::function<double(const tarch::la::Vector<Dimensions,double>&)>
328 );
329 }
330}
331
332
333
const float const float const float struct part *restrict struct part *restrict const float a
Definition hydro_iact.h:55
ElementWiseAssemblyMatrix getPoissonMatrixWithJumpingCoefficient(const tarch::la::Vector< Dimensions, double > &cellCentre, const tarch::la::Vector< Dimensions, double > &h, int integrationPointsPerAxis, std::function< double(const tarch::la::Vector< Dimensions, double > &)>)
Integrate over element with a subvoxel scheme (cmp.
tarch::la::Matrix< TwoPowerD, TwoPowerD, double > ElementWiseAssemblyMatrix
Element-wise Matrix.
Definition Stencil.h:18
tarch::la::Vector< 3, double > get1DIdentity()
Stencil is not scaled at all with any mesh width.
void applyNeumannBC(tarch::la::Vector< ThreePowerD, double > &stencil, double &rhs, const std::bitset< Dimensions *2 > &boundaryFaceNormals, const tarch::la::Vector< Dimensions, double > &h, double derivative)
tarch::la::Vector< FivePowerD, double > getDLinearInterpolation()
Computes the stencil for a d-linear interpolation.
tarch::la::Vector< 5, double > get1DLinearInterpolationStencil()
Stencil is not scaled at all with any mesh width.
tarch::la::Vector< ThreePowerD, double > getUpwindDiscretisedConvection(const tarch::la::Vector< Dimensions, double > &h, const tarch::la::Vector< Dimensions, double > &convCoeff)
Computes a convection operator, defined by the convection coefficients, discretised with the first or...
Stencil exchangeCoordinates(const Stencil &stencil, int coord0, int coord1)
Exchanges the coordinates of a stencil.
tarch::la::Vector< ThreePowerD, double > getMassMatrix(const tarch::la::Vector< Dimensions, double > &h)
Computes the mass matrix.
tarch::la::Vector< ThreePowerD, double > Stencil
Stencil.
Definition Stencil.h:29
tarch::la::Vector< ThreePowerD, std::complex< double > > getHelmholtzShiftMatrix(const tarch::la::Vector< Dimensions, double > &h, const std::complex< double > &phi)
Computes the Helmholtz shift term matrix.
tarch::la::Vector< 3, double > get1DCentralDifferencesStencil()
tarch::la::Vector< ThreePowerD, double > getLaplacian(const tarch::la::Vector< Dimensions, double > &h, const tarch::la::Vector< Dimensions, double > &scaling=1.0)
Computes the Laplacian.
tarch::la::Vector< 3, double > get1DRightIdentity()
void preprocessBoundaryStencil(Stencil &stencil, const std::bitset< Dimensions *2 > &boundaryFaceNormals)
Preprocess boundary stencils held by boundary vertices.
tarch::la::Vector< 3, double > get1DMeanValueStencil()
Stencil is not scaled at all with any mesh width.
tarch::la::Vector< StencilSize *StencilSize, double > stencilProduct(const tarch::la::Vector< StencilSize, double > &a, const tarch::la::Vector< StencilSize, double > &b)
tarch::la::Vector< 3, double > get1DDownwindStencil()
Stencil is not scaled at all with any mesh width.
tarch::la::Vector< 3, double > get1DLaplaceStencil()
Stencil is not scaled at all with any mesh width.
tarch::la::Vector< 3, double > get1DUpwindStencil()
tarch::la::Vector< 3, double > get1DLeftIdentity()
tarch::la::Vector< 3, double > get1DMassStencil()
tarch::la::Vector< ThreePowerD, double > getIdentity(const tarch::la::Vector< Dimensions, double > &h)
Computes the mass matrix.
Stencil extractElementStencil(const Stencil &stencil, const tarch::la::Vector< Dimensions, int > &adjacentCell)
The standard stencil spans 2^d cells adjacent to a vertex.
Simple vector class.
Definition Vector.h:150