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
34 use sort_oct_m
35
36 implicit none
37
38 private
39 public :: &
45
46
47
48contains
49
50 ! ---------------------------------------------------------
53 subroutine stencil_stargeneral_get_arms(this, dim, coord_system)
54 type(stencil_t), intent(inout) :: this
55 integer, intent(in) :: dim
56 class(coordinate_system_t), intent(in) :: coord_system
57
58 real(real64) :: theta(dim), basis_vectors(dim, dim)
59 integer :: arm(dim)
60 integer :: ii, jj, kk
61 real(real64) :: norm
62
63 integer, allocatable :: cand_vecs(:,:), idx(:)
64 real(real64), allocatable :: cand_norms(:)
65 integer :: n_cand, max_cand
66 logical :: active(3)
67 integer :: n_found
68 integer :: det
69 integer :: v1_xy, v1_xz, v1_yz
70 integer :: v2_xy, v2_xz, v2_yz
71 integer :: v3_xy, v3_xz, v3_yz
72 integer :: i1, i2, i3, i, j, k
73
75
76 assert(dim <= 3)
77
78 this%stargeneral%narms = 0
79 this%stargeneral%arms = 0
80
81 if (dim == 1) then
82 !we are done
84 return
85 end if
86
87 select type (coord_system)
88 class is (affine_coordinates_t)
89 basis_vectors = coord_system%basis%vectors
90 class default
91 message(1) = "Stencil stargeneral can currently only be used with affine coordinate systems."
92 call messages_fatal(1)
93 end select
94
95 if (dim == 2) then
96 !get the angle between the primitive vectors
97 theta(1) = acos(dot_product(basis_vectors(:, 1), basis_vectors(:, 2)))
98
99 if (theta(1) < m_pi*m_half) then
100 this%stargeneral%narms = this%stargeneral%narms + 1
101 arm = (/1,-1/)
102 this%stargeneral%arms(this%stargeneral%narms, 1:dim) = arm
103 else if (theta(1) > m_pi*m_half) then
104 this%stargeneral%narms = this%stargeneral%narms + 1
105 arm = (/1,+1/)
106 this%stargeneral%arms(this%stargeneral%narms, 1:dim) = arm
107 end if
108 !if theta == pi/2 we do not need additional arms
109
111 return
112 end if
113
114 ! dim>2
115
116 !We first count how many arms we need
117 theta(1) = acos(dot_product(basis_vectors(:, 1), basis_vectors(:, 2)))
118 if (abs(theta(1) - m_pi*m_half) > 1e-8_real64) this%stargeneral%narms = this%stargeneral%narms + 1
119
120 theta(2) = acos(dot_product(basis_vectors(:, 2), basis_vectors(:, 3)))
121 if (abs(theta(2)-m_pi*m_half) > 1e-8_real64) this%stargeneral%narms = this%stargeneral%narms + 1
122
123 theta(3) = acos(dot_product(basis_vectors(:, 3), basis_vectors(:, 1)))
124 if (abs(theta(3) - m_pi*m_half) > 1e-8_real64) this%stargeneral%narms = this%stargeneral%narms + 1
125
126
127 !Cubic cell, no extra arms
128 if (this%stargeneral%narms == 0) then
130 return
131 end if
132
133 ! Determine active off-diagonal components
134 ! Mapping:
135 ! active(1): xy (theta 12) -> theta(1)
136 ! active(2): xz (theta 13) -> theta(3)
137 ! active(3): yz (theta 23) -> theta(2)
138 active = .false.
139 if (abs(theta(1) - m_pi*m_half) > 1e-8_real64) active(1) = .true.
140 if (abs(theta(3) - m_pi*m_half) > 1e-8_real64) active(2) = .true.
141 if (abs(theta(2) - m_pi*m_half) > 1e-8_real64) active(3) = .true.
142
143 ! Generate candidates in range [-2, 2]
144 max_cand = 5**3
145 safe_allocate(cand_vecs(3, max_cand))
146 safe_allocate(cand_norms(max_cand))
147 safe_allocate(idx(max_cand))
148 n_cand = 0
149
150 do ii= -2,2
151 do jj = -2,2
152 do kk = -2,2
153 ! Skip origin
154 if (ii == 0 .and. jj == 0 .and. kk == 0) cycle
155
156 ! Skip axes
157 if (ii /= 0 .and. jj == 0 .and. kk == 0) cycle
158 if (ii == 0 .and. jj /= 0 .and. kk == 0) cycle
159 if (ii == 0 .and. jj == 0 .and. kk /= 0) cycle
160
161 ! Skip orthogonal planes if angle is 90
162 if (.not. active(1) .and. kk == 0) cycle ! xy plane, but x-y orthogonal
163 if (.not. active(3) .and. ii == 0) cycle ! yz plane, but y-z orthogonal
164 if (.not. active(2) .and. jj == 0) cycle ! xz plane, but x-z orthogonal
165
166 ! We only store one vector from the pair (v, -v) as they define the same line.
167 ! To ensure uniqueness, we enforce that the first non-zero component is negative.
168 if (ii > 0) cycle
169 if (ii == 0 .and. jj > 0) cycle
170 if (ii == 0 .and. jj == 0 .and. kk > 0) cycle
171
172 norm = sum((ii*basis_vectors(:, 1) + jj*basis_vectors(:, 2) + kk*basis_vectors(:, 3))**2)
173
174 n_cand = n_cand + 1
175 cand_vecs(1, n_cand) = ii
176 cand_vecs(2, n_cand) = jj
177 cand_vecs(3, n_cand) = kk
178 cand_norms(n_cand) = norm
179 end do
180 end do
181 end do
182
183 ! Sort candidates by norm
184 call robust_sort_by_abs(cand_norms(1:n_cand), cand_vecs(:, 1:n_cand), idx(1:n_cand))
185
186 ! Select narms vectors
187 n_found = 0
188
189 ! We iterate through sorted candidates
190 ! Case 1: narms = 1
191 ! This is a 3D system with two orthogonal directions.
192 ! We only need to have a non-zero component along the non-orthogonal direction.
193 if (this%stargeneral%narms == 1) then
194 do i = 1, n_cand
195 i1 = idx(i)
196 ! Check if this vector has non-zero component for the active direction
197 if (active(1) .and. cand_vecs(1,i1)*cand_vecs(2,i1) /= 0) n_found = 1
198 if (active(2) .and. cand_vecs(1,i1)*cand_vecs(3,i1) /= 0) n_found = 1
199 if (active(3) .and. cand_vecs(2,i1)*cand_vecs(3,i1) /= 0) n_found = 1
200
201 if (n_found == 1) then
202 this%stargeneral%arms(1, 1:3) = cand_vecs(:, i1)
203 exit
204 end if
205 end do
206 end if
207
208 ! Case 2: narms = 2
209 ! This corresponds to a 3D system, with one orthogonal direction.
210 ! We just need to pick two linearly independent vectors.
211 if (this%stargeneral%narms == 2) then
212 do i = 1, n_cand
213 i1 = idx(i)
214 do j = i + 1, n_cand
215 i2 = idx(j)
216
217 ! Check linear independence (cross product non-zero)
218 if ( (cand_vecs(2, i1) * cand_vecs(3, i2) - cand_vecs(3, i1) * cand_vecs(2, i2)) /= 0 .or. &
219 (cand_vecs(3, i1) * cand_vecs(1, i2) - cand_vecs(1, i1) * cand_vecs(3, i2)) /= 0 .or. &
220 (cand_vecs(1, i1) * cand_vecs(2, i2) - cand_vecs(2, i1) * cand_vecs(1, i2)) /= 0 ) then
221
222 this%stargeneral%arms(1, 1:3) = cand_vecs(:, i1)
223 this%stargeneral%arms(2, 1:3) = cand_vecs(:, i2)
224 n_found = 2
225 exit
226 end if
227 end do
228 if (n_found == 2) exit
229 end do
230 end if
231
232 ! Case 3: narms = 3 (Triclinic)
233 if (this%stargeneral%narms == 3) then
234 do i = 1, n_cand
235 i1 = idx(i)
236 do j = i + 1, n_cand
237 i2 = idx(j)
238 do k = j + 1, n_cand
239 i3 = idx(k)
240
241 ! We check if the off-diagonal terms of the outer product of the vectors
242 ! form a linearly independent set. This is required for the matrix of the
243 ! linear system to determine the weights to be non-singular.
244 ! See Natan et al., PRB 78, 075109 (2008), eq. 19.
245 ! The matrix is M_ij = (n_l)_i * (n_l)_j, where n_l are the arm vectors.
246 ! The off-diagonal terms are xy, xz, yz.
247 ! We form the 3x3 matrix where each row corresponds to one arm vector,
248 ! and the columns correspond to the xy, xz, yz components.
249
250 v1_xy = cand_vecs(1, i1) * cand_vecs(2, i1)
251 v1_xz = cand_vecs(1, i1) * cand_vecs(3, i1)
252 v1_yz = cand_vecs(2, i1) * cand_vecs(3, i1)
253
254 v2_xy = cand_vecs(1, i2) * cand_vecs(2, i2)
255 v2_xz = cand_vecs(1, i2) * cand_vecs(3, i2)
256 v2_yz = cand_vecs(2, i2) * cand_vecs(3, i2)
257
258 v3_xy = cand_vecs(1, i3) * cand_vecs(2, i3)
259 v3_xz = cand_vecs(1, i3) * cand_vecs(3, i3)
260 v3_yz = cand_vecs(2, i3) * cand_vecs(3, i3)
261
262 det = v1_xy * ( v2_xz * v3_yz - v3_xz * v2_yz ) &
263 - v1_xz * ( v2_xy * v3_yz - v3_xy * v2_yz ) &
264 + v1_yz * ( v2_xy * v3_xz - v3_xy * v2_xz )
265
266 if (det /= 0) then
267 this%stargeneral%arms(1, 1:3) = cand_vecs(:, i1)
268 this%stargeneral%arms(2, 1:3) = cand_vecs(:, i2)
269 this%stargeneral%arms(3, 1:3) = cand_vecs(:, i3)
270 n_found = 3
271 exit
272 end if
273 end do
274 if (n_found == 3) exit
275 end do
276 if (n_found == 3) exit
277 end do
278 end if
279
280 if (n_found < this%stargeneral%narms) then
281 write(message(1), '(a, i0, a, i0)') 'Stencil stargeneral could not find enough independent arms. Found ', &
282 n_found, ' needed ', this%stargeneral%narms
283 call messages_fatal(1)
284 end if
285
286 safe_deallocate_a(cand_vecs)
287 safe_deallocate_a(cand_norms)
288 safe_deallocate_a(idx)
289
291 end subroutine stencil_stargeneral_get_arms
292
293
294 ! ---------------------------------------------------------
295 integer function stencil_stargeneral_size_lapl(this, dim, order) result(n)
296 type(stencil_t), intent(inout) :: this
297 integer, intent(in) :: dim
298 integer, intent(in) :: order
299
301
302 !normal star
303 n = 2*dim*order + 1
304
305 ! star general
306 n = n + 2 * order * this%stargeneral%narms
307
308
311
312
313 ! ---------------------------------------------------------
316 integer function stencil_stargeneral_extent(dir, order)
317 integer, intent(in) :: dir
318 integer, intent(in) :: order
319
320 integer :: extent
321
323
324 extent = 0
325 if (dir >= 1 .or. dir <= 3) then
326 if (order <= 2) then
327 extent = 2
328 else
329 extent = order
330 end if
331 end if
333
335 end function stencil_stargeneral_extent
336
337
338
339 ! ---------------------------------------------------------
340 subroutine stencil_stargeneral_get_lapl(this, dim, order)
341 type(stencil_t), intent(inout) :: this
342 integer, intent(in) :: dim
343 integer, intent(in) :: order
344
345 integer :: i, j, n
346 logical :: got_center
347
349
350 call stencil_allocate(this, dim, stencil_stargeneral_size_lapl(this, dim, order))
351
352 n = 1
353 select case (dim)
354 case (1)
355 n = 1
356 do i = 1, dim
357 do j = -order, order
358 if (j == 0) cycle
359 n = n + 1
360 this%points(i, n) = j
361 end do
362 end do
363 case (2)
364 n = 1
365 do i = 1, dim
366 do j = -order, order
367 if (j == 0) cycle
368 n = n + 1
369 this%points(i, n) = j
370 end do
371 end do
372
373 do j = -order, order
374 if (j == 0) cycle
375 do i = 1, this%stargeneral%narms
376 n = n + 1
377 this%points(1:2, n) = this%stargeneral%arms(i, 1:2)*j
378 end do
379 end do
380
381 case (3)
382 got_center = .false.
383
384 n = 0
385 do i = 1, dim
386 do j = -order, order
387
388 ! count center only once
389 if (j == 0) then
390 if (got_center) then
391 cycle
392 else
393 got_center = .true.
394 end if
395
396 end if
397 n = n + 1
398 this%points(i, n) = j
399 end do
400 end do
401
402 do j = -order, order
403 if (j == 0) cycle
404 do i = 1, this%stargeneral%narms
405 n = n + 1
406 this%points(1:3, n) = this%stargeneral%arms(i, 1:3)*j
407 end do
408 end do
409
410 end select
412 call stencil_init_center(this)
413
415 end subroutine stencil_stargeneral_get_lapl
416
417
418
419
420 ! ---------------------------------------------------------
421 subroutine stencil_stargeneral_pol_lapl(this, dim, order, pol)
422 type(stencil_t), intent(in) :: this
423 integer, intent(in) :: dim
424 integer, intent(in) :: order
425 integer, intent(out) :: pol(:,:)
426
427 integer :: i, j, n, j1
428 logical :: zeros_treated(dim)
429
431
432 n = 1
433 select case (dim)
434 case (1)
435 n = 1
436 pol(:,:) = 0
437 do i = 1, dim
438 do j = 1, 2*order
439 n = n + 1
440 pol(i, n) = j
441 end do
442 end do
443 case (2)
444 n = 1
445 pol(:,:) = 0
446 do i = 1, dim
447 do j = 1, 2*order
448 n = n + 1
449 pol(i, n) = j
450 end do
451 end do
452
453 do j = 1, 2*order
454 do i = 1, this%stargeneral%narms
455 n = n + 1
456 if (sum(this%stargeneral%arms(i,1:dim)) == 0) then
457 pol(1:2, n) = (/j,1/)
458 else
459 pol(1:2, n) = (/1,j/)
460 end if
461 end do
462 end do
463
464 case (3)
465 n = 1
466 pol(:,:) = 0
467 do i = 1, dim
468 do j = 1, 2*order
469 n = n + 1
470 pol(i, n) = j
471 end do
472 end do
473
474 ! first, treat arms with zeros: they should have polynomials with zeros in the
475 ! corresponding dimension
476 zeros_treated = .false.
477 do i = 1, this%stargeneral%narms
478 if (this%stargeneral%arms(i, 1) == 0) then
479 do j1 = 1, 2*order
480 n = n + 1
481 pol(1:3, n) = (/0, 1, j1/)
482 end do
483 zeros_treated(1) = .true.
484 else if (this%stargeneral%arms(i, 2) == 0) then
485 do j1 = 1, 2*order
486 n = n + 1
487 pol(1:3, n) = (/j1, 0, 1/)
488 end do
489 zeros_treated(2) = .true.
490 else if (this%stargeneral%arms(i, 3) == 0) then
491 do j1 = 1, 2*order
492 n = n + 1
493 pol(1:3, n) = (/1, j1, 0/)
494 end do
495 zeros_treated(3) = .true.
496 end if
497 end do ! end loop on number of arms
498 ! then treat the remaining arms, use the polynomials that have not been used yet
499 do i = 1, this%stargeneral%narms - count(zeros_treated)
500 if (.not. zeros_treated(1)) then
501 do j1 = 1, 2*order
502 n = n + 1
503 pol(1:3, n) = (/0, 1, j1/)
504 end do
505 zeros_treated(1) = .true.
506 else if (.not. zeros_treated(2)) then
507 do j1 = 1, 2*order
508 n = n + 1
509 pol(1:3, n) = (/j1, 0, 1/)
510 end do
511 zeros_treated(2) = .true.
512 else if (.not. zeros_treated(3)) then
513 do j1 = 1, 2*order
514 n = n + 1
515 pol(1:3, n) = (/1, j1, 0/)
516 end do
517 zeros_treated(3) = .true.
518 end if
519 end do ! end loop on number of arms
520 end select
521
523 end subroutine stencil_stargeneral_pol_lapl
524
525
527
528!! Local Variables:
529!! mode: f90
530!! coding: utf-8
531!! End:
double acos(double __x) __attribute__((__nothrow__
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:189
real(real64), parameter, public m_half
Definition: global.F90:197
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:119
This module defines stencils used in Octopus.
Definition: stencil.F90:137
subroutine, public stencil_allocate(this, dim, size)
Definition: stencil.F90:180
subroutine, public stencil_init_center(this)
Definition: stencil.F90:277
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:165
int true(void)