Skip to content

GradientStencil

The LaplacianStencil is the storage container for storing the Laplacian stencil components. It stores only non-zero stencil components for 1D, 2D and 3D stencils and offers iterators allowing convenient evaluation of the Laplacian operators which is beneficial for compact stencils.

The LaplacianStencil can be initialized to use the following stencils:

  • Simple - simple stencils in 1D, 2D and 3D.
  • Isotropic - isotropic stencils in 1D, 2D and 3D.
  • LB - isotropic stencils based on lattice Boltzmann stencils.

The available stencils for 1D, 2D and 3D cases are given below. Note that 1D and 2D stencils for convenience defined in 3D arrays but should only be used with inactive dimensions suppressed.

Simple Laplacian stencils

cpp
const double LaplacianStencil1D_3[3][3][3] = 
          {{{0.0,   0.0,   0.0},
            {0.0,   1.0,   0.0},
            {0.0,   0.0,   0.0}},

           {{0.0,   1.0,   0.0},
            {1.0,  -2.0,   1.0},
            {0.0,   1.0,   0.0}},

           {{0.0,   0.0,   0.0},
            {0.0,   1.0,   0.0},
            {0.0,   0.0,   0.0}}};

3-point 1D Laplacian stencil.

cpp
const double LaplacianStencil2D_5[3][3][3] = 
          {{{0.0,   0.0,   0.0},
            {0.0,   1.0,   0.0},
            {0.0,   0.0,   0.0}},

           {{0.0,   1.0,   0.0},
            {1.0,  -4.0,   1.0},
            {0.0,   1.0,   0.0}},

           {{0.0,   0.0,   0.0},
            {0.0,   1.0,   0.0},
            {0.0,   0.0,   0.0}}};

5-point simple 2D Laplacian stencil.

cpp
const double LaplacianStencil3D_7[3][3][3] = 
          {{{0.0,   0.0,  0.0},
            {0.0,   1.0,  0.0},
            {0.0,   0.0,  0.0}},

           {{0.0,   1.0,  0.0},
            {1.0,  -6.0,  1.0},
            {0.0,   1.0,  0.0}},

           {{0.0,   0.0,  0.0},
            {0.0,   1.0,  0.0},
            {0.0,   0.0,  0.0}}};

7-point simple 3D Laplacian stencil

Isotropic Laplacian stencils

cpp
const double LaplacianStencil2D_9[3][3][3] = 
          {{{    0.0,   1.0/6.0,     0.0},
            {1.0/6.0,   2.0/3.0, 1.0/6.0},
            {    0.0,   1.0/6.0,     0.0}},

           {{1.0/6.0,   2.0/3.0, 1.0/6.0},
            {2.0/3.0, -10.0/3.0, 2.0/3.0},
            {1.0/6.0,   2.0/3.0, 1.0/6.0}},

           {{    0.0,   1.0/6.0,     0.0},
            {1.0/6.0,   2.0/3.0, 1.0/6.0},
            {    0.0,   1.0/6.0,     0.0}}};

9-point 2D Laplacian stencil by Dave Hale.

cpp
const double LaplacianStencil3D_19[3][3][3] = 
          {{{     0.0,   1.0/6.0,      0.0},
            { 1.0/6.0,   1.0/3.0,  1.0/6.0},
            {     0.0,   1.0/6.0,      0.0}},

           {{ 1.0/6.0,   1.0/3.0,  1.0/6.0},
            { 1.0/3.0,      -4.0,  1.0/3.0},
            { 1.0/6.0,   1.0/3.0,  1.0/6.0}},

           {{     0.0,   1.0/6.0,      0.0},
            { 1.0/6.0,   1.0/3.0,  1.0/6.0},
            {     0.0,   1.0/6.0,      0.0}}};

19-point 3D Laplacian stencil by Patra and Karttunen.

cpp
const double LaplacianStencil3D_27a[3][3][3] = 
          {{{1.0/30.0,   1.0/10.0, 1.0/30.0},
            {1.0/10.0,   7.0/15.0, 1.0/10.0},
            {1.0/30.0,   1.0/10.0, 1.0/30.0}},

           {{1.0/10.0,   7.0/15.0, 1.0/10.0},
            {7.0/15.0, -64.0/15.0, 7.0/15.0},
            {1.0/10.0,   7.0/15.0, 1.0/10.0}},

           {{1.0/30.0,   1.0/10.0, 1.0/30.0},
            {1.0/10.0,   7.0/15.0, 1.0/10.0},
            {1.0/30.0,   1.0/10.0, 1.0/30.0}}};

27-point 3D Laplacian stencil by Spotz and Carey (1995).

cpp
const double LaplacianStencil3D_27b[3][3][3] = 
          {{{1.0/48.0,   1.0/8.0, 1.0/48.0},
            {1.0/8.0,   5.0/12.0, 1.0/8.0},
            {1.0/48.0,   1.0/8.0, 1.0/48.0}},

           {{1.0/8.0,   5.0/12.0, 1.0/8.0},
            {5.0/12.0, -25.0/6.0, 5.0/12.0},
            {1.0/8.0,   5.0/12.0, 1.0/8.0}},

           {{1.0/48.0,   1.0/8.0, 1.0/48.0},
            { 1.0/8.0,  5.0/12.0, 1.0/8.0},
            {1.0/48.0,   1.0/8.0, 1.0/48.0}}};

27-point 3D Laplacian stencil by Dave Hale.

Laplacian stencils based on lattice Boltzmann stencils

cpp
const double LaplacianStencil2D_LB[3][3][3] = 
          {{{    0.0, 1.0/6.0,      0.0},
            {1.0/6.0, 2.0/3.0,  1.0/6.0},
            {    0.0, 1.0/6.0,      0.0}},

           {{1.0/6.0,   2.0/3.0, 1.0/6.0},
            {2.0/3.0, -10.0/3.0, 2.0/3.0},
            {1.0/6.0,   2.0/3.0, 1.0/6.0}},

           {{    0.0, 1.0/6.0,     0.0},
            {1.0/6.0, 2.0/3.0, 1.0/6.0},
            {    0.0, 1.0/6.0,     0.0}}};

9-point 2D Laplacian stencil based on D2Q9 lattice Boltzmann stencil. It is the same as 2D stencil by Dave Hale.

cpp
const double LaplacianStencil3D_LB[3][3][3] = 
          {{{1.0/36.0, 1.0/9.0, 1.0/36.0},
            {1.0/9.0,  4.0/9.0, 1.0/9.0 },
            {1.0/36.0, 1.0/9.0, 1.0/36.0}},

           {{1.0/9.0,   4.0/9.0, 1.0/9.0},
            {4.0/9.0, -38.0/9.0, 4.0/9.0},
            {1.0/9.0,   4.0/9.0, 1.0/9.0}},

           {{1.0/36.0, 1.0/9.0, 1.0/36.0},
            {1.0/9.0,  4.0/9.0, 1.0/9.0 },
            {1.0/36.0, 1.0/9.0, 1.0/36.0}}};

27-point Laplacian stencil based on D3Q27 lattice Boltzmann stencil.

LaplacianStencilEntry

The stencil components inside of the LaplacianStencil are stored using the corresponding LaplacianStencilEntry structure which contains three stencil direction coordinates and the weight of the stencil element:

cpp
int di;

X coordinate of the stencil element.

cpp
int dj;

Y coordinate of the stencil element.

cpp
int dk;

Z coordinate of the stencil element.

cpp
double weight;

Weight associated with the stencil element.

For convenience the weight of the zero's element of the stencil is stored additionally in the separate parameter:

cpp
double weight_zero;

Copy of the weight of the center element of the Laplacian stencil.

Laplacian stencil construction methods

cpp
void SetNoCenter(const double UserStencil[3][3][3], double dx,
                 int dNx = 1, int dNy = 1, int dNz = 1);

Sets the Laplacian stencil for diffusion operators without center element using user specified Laplacian stencil.

cpp
void SetNoCenter(const double UserStencil[3][3][3], GridParameters& Dimensions);

Sets the Laplacian stencil for diffusion operators without center element using user specified stencil.

cpp
void Set(const double UserStencil[3][3][3], double dx,
         int dNx = 1, int dNy = 1, int dNz = 1);

Sets the Laplacian stencil using user specified stencil.

cpp
void Set(const double UserStencil[3][3][3], GridParameters& Dimensions);

Sets the Laplacian stencil using user specified stencil.

Container properties

cpp
size_t size() const;

Returns the size of the stencil storage (considers only non-zero entries).

Iterators

cpp
iterator  begin();

Iterator to the begin of the stencil.

cpp
iterator  end();

Iterator to the end of the stencil.

cpp
citerator cbegin() const;

Constant iterator to the begin of the stencil.

cpp
citerator cend() const;

Constant iterator to the end of the stencil.

Data manipulation methods

cpp
iterator erase(iterator it);

Erases the entry pointed to by the given iterator it.

Released under the GNU GPLv3 License.