Octopus
stencil_cube.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
26 use debug_oct_m
27 use global_oct_m
29
30 implicit none
31
32 private
33 public :: &
40
41
42contains
43
44 ! ---------------------------------------------------------
45 integer function stencil_cube_size_lapl(dim, order)
46 integer, intent(in) :: dim
47 integer, intent(in) :: order
48
50
51 stencil_cube_size_lapl = (2*order+1)**dim
52
54 end function stencil_cube_size_lapl
55
56
57 ! ---------------------------------------------------------
58 subroutine stencil_cube_get_lapl(this, dim, order)
59 type(stencil_t), intent(inout) :: this
60 integer, intent(in) :: dim
61 integer, intent(in) :: order
62
63 integer :: i, j, k, n
64
65 push_sub(stencil_cube_get_lapl)
66
67 call stencil_allocate(this, dim, stencil_cube_size_lapl(dim, order))
68
69 n = 1
70 select case (dim)
71 case (1)
72 do i = -order, order
73 this%points(1, n) = i
74 n = n + 1
75 end do
76 case (2)
77 do i = -order, order
78 do j = -order, order
79 this%points(1, n) = i
80 this%points(2, n) = j
81 n = n + 1
82 end do
83 end do
84 case (3)
85 do i = -order, order
86 do j = -order, order
87 do k = -order, order
88 this%points(1, n) = i
89 this%points(2, n) = j
90 this%points(3, n) = k
91 n = n + 1
92 end do
93 end do
94 end do
95 end select
96
97 call stencil_init_center(this)
98
100 end subroutine stencil_cube_get_lapl
101
102
103 ! ---------------------------------------------------------
104 subroutine stencil_cube_polynomials_lapl(dim, order, pol)
105 integer, intent(in) :: dim
106 integer, intent(in) :: order
107 integer, intent(out) :: pol(:,:)
108
109 integer :: i, j, k, n
110
112
113 n = 1
114 select case (dim)
115 case (1)
116 do i = 0, 2*order
117 pol(1, n) = i
118 n = n + 1
119 end do
120 case (2)
121 do i = 0, 2*order
122 do j = 0, 2*order
123 pol(1, n) = i
124 pol(2, n) = j
125 n = n + 1
126 end do
127 end do
128 case (3)
129 do i = 0, 2*order
130 do j = 0, 2*order
131 do k = 0, 2*order
132 pol(1, n) = i
133 pol(2, n) = j
134 pol(3, n) = k
135 n = n + 1
136 end do
137 end do
138 end do
139 end select
140
142 end subroutine stencil_cube_polynomials_lapl
143
144
148
149 ! ---------------------------------------------------------
150 integer function stencil_cube_size_grad(dim, order)
151 integer, intent(in) :: dim
152 integer, intent(in) :: order
153
154 push_sub(stencil_cube_size_grad)
155
157
159 end function stencil_cube_size_grad
160
161
162 ! ---------------------------------------------------------
163 subroutine stencil_cube_get_grad(this, dim, order)
164 type(stencil_t), intent(inout) :: this
165 integer, intent(in) :: dim
166 integer, intent(in) :: order
167
168 push_sub(stencil_cube_get_grad)
169
170 call stencil_cube_get_lapl(this, dim, order)
171
172 pop_sub(stencil_cube_get_grad)
173 end subroutine stencil_cube_get_grad
174
175
176 ! ---------------------------------------------------------
177 subroutine stencil_cube_polynomials_grad(dim, order, pol)
178 integer, intent(in) :: dim
179 integer, intent(in) :: order
180 integer, intent(out) :: pol(:,:)
181
183
184 call stencil_cube_polynomials_lapl(dim, order, pol)
185
187 end subroutine stencil_cube_polynomials_grad
188
189end module stencil_cube_oct_m
190
191!! Local Variables:
192!! mode: f90
193!! coding: utf-8
194!! End:
This module defines routines, generating operators for a cubic stencil.
subroutine, public stencil_cube_get_lapl(this, dim, order)
subroutine, public stencil_cube_polynomials_lapl(dim, order, pol)
integer function, public stencil_cube_size_grad(dim, order)
Now come the gradient routines. As this stencil is the same for the laplacian and the gradient,...
integer function, public stencil_cube_size_lapl(dim, order)
subroutine, public stencil_cube_polynomials_grad(dim, order, pol)
subroutine, public stencil_cube_get_grad(this, dim, order)
This module defines stencils used in Octopus.
Definition: stencil.F90:135
subroutine, public stencil_allocate(this, dim, size)
Definition: stencil.F90:178
subroutine, public stencil_init_center(this)
Definition: stencil.F90:275
The class representing the stencil, which is used for non-local mesh operations.
Definition: stencil.F90:163