Octopus
stencil_stargeneral.F90
Go to the documentation of this file.
1!! Copyright (C) 2016 U. De Giovannini, H Huebener; 2019 N. Tancogne-Dejean
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
28 use debug_oct_m
29 use global_oct_m
30 use, intrinsic :: iso_fortran_env
33
34 implicit none
35
36 private
37 public :: &
43
44
45
46contains
47
48 ! ---------------------------------------------------------
51 subroutine stencil_stargeneral_get_arms(this, dim, coord_system)
52 type(stencil_t), intent(inout) :: this
53 integer, intent(in) :: dim
54 class(coordinate_system_t), intent(in) :: coord_system
55
56 real(real64) :: theta(dim), basis_vectors(dim, dim)
57 integer :: arm(dim)
58 integer :: ii, jj, kk
59 real(real64) :: norm, min_norm(dim)
60
61 real(real64), parameter :: tol = 1e-14_real64
62
64
65 assert(dim <= 3)
66
67 this%stargeneral%narms = 0
68 this%stargeneral%arms = 0
69
70 if (dim == 1) then
71 !we are done
73 return
74 end if
75
76 select type (coord_system)
77 class is (affine_coordinates_t)
78 basis_vectors = coord_system%basis%vectors
79 class default
80 message(1) = "Stencil stargeneral can currently only be used with affine coordinate systems."
81 call messages_fatal(1)
82 end select
83
84 if (dim == 2) then
85 !get the angle between the primitive vectors
86 theta(1) = acos(dot_product(basis_vectors(:, 1), basis_vectors(:, 2)))
87
88 if (theta(1) < m_pi*m_half) then
89 this%stargeneral%narms = this%stargeneral%narms + 1
90 arm = (/1,-1/)
91 this%stargeneral%arms(this%stargeneral%narms, 1:dim) = arm
92 else if (theta(1) > m_pi*m_half) then
93 this%stargeneral%narms = this%stargeneral%narms + 1
94 arm = (/1,+1/)
95 this%stargeneral%arms(this%stargeneral%narms, 1:dim) = arm
96 end if
97 !if theta == pi/2 we do not need additional arms
98
99 ! NB: you have supposed the axis of the 2D system is along z.
100 !we are done
102 return
103 end if
104
105 ! dim>2
106
107 !We first count how many arms we need
108 theta(1) = acos(dot_product(basis_vectors(:, 1), basis_vectors(:, 2)))
109 if (abs(theta(1) - m_pi*m_half) > 1e-8_real64) this%stargeneral%narms = this%stargeneral%narms + 1
110
111 theta(2) = acos(dot_product(basis_vectors(:, 2), basis_vectors(:, 3)))
112 if (abs(theta(2)-m_pi*m_half) > 1e-8_real64) this%stargeneral%narms = this%stargeneral%narms + 1
113
114 theta(3) = acos(dot_product(basis_vectors(:, 3), basis_vectors(:, 1)))
115 if (abs(theta(3) - m_pi*m_half) > 1e-8_real64) this%stargeneral%narms = this%stargeneral%narms + 1
116
117
118 !Cubic cell, no extra arms
119 if (this%stargeneral%narms == 0) then
121 return
122 end if
123
124
125 !We find the three closest neighbours, to define the direction of the arms
126 min_norm = m_huge
127 do ii= -1,1
128 do jj = -1,1
129 do kk = -1,1
130 !The nearest points could not be along the primitive axis
131 if (ii == 0 .and. jj == 0) cycle
132 if (ii == 0 .and. kk == 0) cycle
133 if (kk == 0 .and. jj == 0) cycle
134
135 !Nor in a plane formed by orthogonal vectors
136 if (abs(theta(1) - m_pi*m_half) <= 1e-8_real64 .and. kk == 0) cycle
137 if (abs(theta(2) - m_pi*m_half) <= 1e-8_real64 .and. ii == 0) cycle
138 if (abs(theta(3) - m_pi*m_half) <= 1e-8_real64 .and. jj == 0) cycle
139
140 norm = sum((ii*basis_vectors(:, 1) + jj*basis_vectors(:, 2) + kk*basis_vectors(:, 3))**2)
141
142 if (min_norm(1) - norm > tol) then
143 if (min_norm(2) - min_norm(1) > tol .and. this%stargeneral%narms >= 2) then
144 if (min_norm(3) - min_norm(2) > tol .and. this%stargeneral%narms == 3) then
145 this%stargeneral%arms(3, 1:3) = this%stargeneral%arms(2, 1:3)
146 min_norm(3) = min_norm(2)
147 end if
148 this%stargeneral%arms(2, 1:3) = this%stargeneral%arms(1, 1:3)
149 min_norm(2) = min_norm(1)
150 end if
151 min_norm(1) = norm
152 this%stargeneral%arms(1, 1:3) = (/ii, jj, kk/)
153 cycle
154 end if
155
156 ! Note: this should not be needed since we now have tolerances
157 if (this%stargeneral%arms(1, 1) == -ii &
158 .and. this%stargeneral%arms(1, 2) == -jj &
159 .and. this%stargeneral%arms(1, 3) == -kk) cycle
160
161 if (this%stargeneral%narms == 1) cycle
162
163 if (min_norm(2) - norm > tol) then
164 if (min_norm(3) - min_norm(2) > tol .and. this%stargeneral%narms == 3) then
165 this%stargeneral%arms(3, 1:3) = this%stargeneral%arms(2, 1:3)
166 min_norm(3) = min_norm(2)
167 end if
168 min_norm(2) = norm
169 this%stargeneral%arms(2, 1:3) = (/ii, jj, kk/)
170 cycle
171 end if
172
173 if (this%stargeneral%narms == 2) cycle
174
175 ! Note: this should not be needed since we now have tolerances
176 if (this%stargeneral%arms(2, 1) == -ii &
177 .and. this%stargeneral%arms(2, 2) == -jj &
178 .and. this%stargeneral%arms(2, 3) == -kk) cycle
179
180 if (min_norm(3) - norm > tol) then
181 min_norm(3) = norm
182 this%stargeneral%arms(3, 1:3) = (/ii, jj, kk/)
183 end if
184 end do
185 end do
186 end do
187
189 end subroutine stencil_stargeneral_get_arms
190
191
192 ! ---------------------------------------------------------
193 integer function stencil_stargeneral_size_lapl(this, dim, order) result(n)
194 type(stencil_t), intent(inout) :: this
195 integer, intent(in) :: dim
196 integer, intent(in) :: order
197
199
200 !normal star
201 n = 2*dim*order + 1
202
203 ! star general
204 n = n + 2 * order * this%stargeneral%narms
205
206
209
210
211 ! ---------------------------------------------------------
214 integer function stencil_stargeneral_extent(dir, order)
215 integer, intent(in) :: dir
216 integer, intent(in) :: order
217
218 integer :: extent
219
221
222 extent = 0
223 if (dir >= 1 .or. dir <= 3) then
224 if (order <= 2) then
225 extent = 2
226 else
227 extent = order
228 end if
229 end if
231
233 end function stencil_stargeneral_extent
234
235
236
237 ! ---------------------------------------------------------
238 subroutine stencil_stargeneral_get_lapl(this, dim, order)
239 type(stencil_t), intent(inout) :: this
240 integer, intent(in) :: dim
241 integer, intent(in) :: order
242
243 integer :: i, j, n
244 logical :: got_center
245
247
248 call stencil_allocate(this, dim, stencil_stargeneral_size_lapl(this, dim, order))
249
250 n = 1
251 select case (dim)
252 case (1)
253 n = 1
254 do i = 1, dim
255 do j = -order, order
256 if (j == 0) cycle
257 n = n + 1
258 this%points(i, n) = j
259 end do
260 end do
261 case (2)
262 n = 1
263 do i = 1, dim
264 do j = -order, order
265 if (j == 0) cycle
266 n = n + 1
267 this%points(i, n) = j
268 end do
269 end do
270
271 do j = -order, order
272 if (j == 0) cycle
273 do i = 1, this%stargeneral%narms
274 n = n + 1
275 this%points(1:2, n) = this%stargeneral%arms(i, 1:2)*j
276 end do
277 end do
278
279 case (3)
280 got_center = .false.
281
282 n = 0
283 do i = 1, dim
284 do j = -order, order
285
286 ! count center only once
287 if (j == 0) then
288 if (got_center) then
289 cycle
290 else
291 got_center = .true.
292 end if
293
294 end if
295 n = n + 1
296 this%points(i, n) = j
297 end do
298 end do
299
300 do j = -order, order
301 if (j == 0) cycle
302 do i = 1, this%stargeneral%narms
303 n = n + 1
304 this%points(1:3, n) = this%stargeneral%arms(i, 1:3)*j
305 end do
306 end do
308 end select
309
310 call stencil_init_center(this)
311
313 end subroutine stencil_stargeneral_get_lapl
314
315
316
317
318 ! ---------------------------------------------------------
319 subroutine stencil_stargeneral_pol_lapl(this, dim, order, pol)
320 type(stencil_t), intent(in) :: this
321 integer, intent(in) :: dim
322 integer, intent(in) :: order
323 integer, intent(out) :: pol(:,:)
324
325 integer :: i, j, n, j1
326 logical :: zeros_treated(dim)
327
329
330 n = 1
331 select case (dim)
332 case (1)
333 n = 1
334 pol(:,:) = 0
335 do i = 1, dim
336 do j = 1, 2*order
337 n = n + 1
338 pol(i, n) = j
339 end do
340 end do
341 case (2)
342 n = 1
343 pol(:,:) = 0
344 do i = 1, dim
345 do j = 1, 2*order
346 n = n + 1
347 pol(i, n) = j
348 end do
349 end do
350
351 do j = 1, 2*order
352 do i = 1, this%stargeneral%narms
353 n = n + 1
354 if (sum(this%stargeneral%arms(i,1:dim)) == 0) then
355 pol(1:2, n) = (/j,1/)
356 else
357 pol(1:2, n) = (/1,j/)
358 end if
359 end do
360 end do
361
362 case (3)
363 n = 1
364 pol(:,:) = 0
365 do i = 1, dim
366 do j = 1, 2*order
367 n = n + 1
368 pol(i, n) = j
369 end do
370 end do
371
372 ! first, treat arms with zeros: they should have polynomials with zeros in the
373 ! corresponding dimension
374 zeros_treated = .false.
375 do i = 1, this%stargeneral%narms
376 if (this%stargeneral%arms(i, 1) == 0) then
377 do j1 = 1, 2*order
378 n = n + 1
379 pol(1:3, n) = (/0, 1, j1/)
380 end do
381 zeros_treated(1) = .true.
382 else if (this%stargeneral%arms(i, 2) == 0) then
383 do j1 = 1, 2*order
384 n = n + 1
385 pol(1:3, n) = (/j1, 0, 1/)
386 end do
387 zeros_treated(2) = .true.
388 else if (this%stargeneral%arms(i, 3) == 0) then
389 do j1 = 1, 2*order
390 n = n + 1
391 pol(1:3, n) = (/1, j1, 0/)
392 end do
393 zeros_treated(3) = .true.
394 end if
395 end do ! end loop on number of arms
396 ! then treat the remaining arms, use the polynomials that have not been used yet
397 do i = 1, this%stargeneral%narms - count(zeros_treated)
398 if (.not. zeros_treated(1)) then
399 do j1 = 1, 2*order
400 n = n + 1
401 pol(1:3, n) = (/0, 1, j1/)
402 end do
403 zeros_treated(1) = .true.
404 else if (.not. zeros_treated(2)) then
405 do j1 = 1, 2*order
406 n = n + 1
407 pol(1:3, n) = (/j1, 0, 1/)
408 end do
409 zeros_treated(2) = .true.
410 else if (.not. zeros_treated(3)) then
411 do j1 = 1, 2*order
412 n = n + 1
413 pol(1:3, n) = (/1, j1, 0/)
414 end do
415 zeros_treated(3) = .true.
416 end if
417 end do ! end loop on number of arms
418 end select
419
421 end subroutine stencil_stargeneral_pol_lapl
422
423
425
426!! Local Variables:
427!! mode: f90
428!! coding: utf-8
429!! End:
double acos(double __x) __attribute__((__nothrow__
real(real64), parameter, public m_huge
Definition: global.F90:205
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:185
real(real64), parameter, public m_half
Definition: global.F90:193
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
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
This module defines routines, generating operators for a generalized star stencil.
subroutine, public stencil_stargeneral_get_arms(this, dim, coord_system)
Finds the direction of the arms of the star-general stencil as described in Natan et al....
subroutine, public stencil_stargeneral_pol_lapl(this, dim, order, pol)
integer function, public stencil_stargeneral_extent(dir, order)
Returns maximum extension of the stencil in spatial direction dir = 1, 2, 3 for a given discretizatio...
subroutine, public stencil_stargeneral_get_lapl(this, dim, order)
integer function, public stencil_stargeneral_size_lapl(this, dim, order)
The class representing the stencil, which is used for non-local mesh operations.
Definition: stencil.F90:163
int true(void)