Octopus
math.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
22module math_oct_m
23 use debug_oct_m
24 use global_oct_m
25 use, intrinsic :: iso_c_binding, only: c_double
26 use, intrinsic :: iso_fortran_env
31
32 implicit none
33
34 private
35 public :: &
36 is_close, &
38 ylmr_cmplx, &
39 ylmr_real, &
40 weights, &
41 factorial, &
42 hermite, &
47 ddelta, &
50 even, &
51 odd, &
55 pad, &
56 pad_pow2, &
57 log2, &
58 expm1, &
60 phi1, &
61 phi2, &
63 is_prime, &
75
76 interface is_close
77 module procedure dis_close_scalar, zis_close_scalar
78 end interface
79
80 interface diagonal_matrix
82 end interface diagonal_matrix
83
87
91
94 interface interpolate
97 end interface interpolate
98
99 interface symmetrize_matrix
100 module procedure dsymmetrize_matrix, zsymmetrize_matrix
101 end interface symmetrize_matrix
102
103 interface log2
104 module procedure dlog2, ilog2, llog2
105 end interface log2
106
107 interface pad
108 module procedure pad4, pad8, pad48, pad88
109 end interface pad
110
111 interface
112 real(c_double) function expm1(x) bind(c, name='expm1')
113 import c_double
114 real(c_double), intent(in), value :: x
115 end function expm1
116 end interface
118contains
119
138 elemental logical function dis_close_scalar(x, y, rtol, atol)
139 real(real64), intent(in) :: x, y
140 real(real64), optional, intent(in) :: rtol
141 real(real64), optional, intent(in) :: atol
142 real(real64) :: atol_, rtol_
143
144 atol_ = optional_default(atol, 1.e-8_real64)
145 rtol_ = optional_default(rtol, 1.e-5_real64)
146
147 dis_close_scalar = abs(x - y) <= (atol_ + rtol_ * abs(y))
148 end function dis_close_scalar
149
151 elemental logical function zis_close_scalar(x, y, rtol, atol)
152 complex(real64), intent(in) :: x, y
153 real(real64), optional, intent(in) :: rtol
154 real(real64), optional, intent(in) :: atol
155 real(real64) :: atol_, rtol_
156
157 atol_ = optional_default(atol, 1.e-8_real64)
158 rtol_ = optional_default(rtol, 1.e-5_real64)
159
160 zis_close_scalar = abs(x - y) <= (atol_ + rtol_ * abs(y))
161 end function zis_close_scalar
162
163
164 ! ---------------------------------------------------------
167 pure function idiagonal_matrix(dim, diag) result(matrix)
168 integer, intent(in) :: dim
169 integer, intent(in) :: diag
170 integer :: matrix(dim, dim)
172 integer :: ii
173
174 matrix = 0
175 do ii = 1, dim
176 matrix(ii, ii) = diag
177 end do
178
179 end function idiagonal_matrix
180
181 ! ---------------------------------------------------------
182 recursive function hermite(n, x) result (h)
183 integer, intent(in) :: n
184 real(real64), intent(in) :: x
185
186 real(real64) :: h
187
188 ! no push_sub, called too frequently
190 if (n <= 0) then
191 h = m_one
192 elseif (n == 1) then
193 h = m_two*x
194 else
195 h = m_two*x*hermite(n-1,x) - m_two*(n-1)*hermite(n-2,x)
196 end if
197
198 end function hermite
199
200
201 ! ---------------------------------------------------------
202 recursive function factorial (n) RESULT (fac)
203 integer, intent(in) :: n
204
205 integer :: fac
206
207 ! no push_sub for recursive function
208
209 if (n <= 1) then
210 fac = 1
211 else
212 fac = n * factorial(n-1)
213 end if
214 end function factorial
215
216 ! ---------------------------------------------------------
218 subroutine ylmr_cmplx(xx, li, mi, ylm)
219 real(real64), intent(in) :: xx(3)
220 integer, intent(in) :: li, mi
221 complex(real64), intent(out) :: ylm
222
223 integer :: i
224 real(real64) :: dx, dy, dz, r, plm, cosm, sinm, cosmm1, sinmm1, cosphi, sinphi
225 real(real64), parameter :: tiny = 1.e-15_real64
226
227 ! no push_sub: this routine is called too many times
228
229 ! if l=0, no calculations are required
230 if (li == 0) then
231 ylm = cmplx(sqrt(m_one/(m_four*m_pi)), m_zero, real64)
232 return
233 end if
234
235 r = sqrt(xx(1)**2 + xx(2)**2 + xx(3)**2)
236
237 ! if r=0, direction is undefined => make ylm=0 except for l=0
238 if (r <= tiny) then
239 ylm = m_z0
240 return
241 end if
242 dx = xx(1)/r
243 dy = xx(2)/r
244 dz = xx(3)/r
245
246 ! get the associated Legendre polynomial (including the normalization factor)
247 ! We need to avoid getting outside [-1:1] due to rounding errors
248 if(abs(dz) > m_one) dz = sign(m_one, dz)
249 plm = loct_legendre_sphplm(li, abs(mi), dz)
250
251 ! compute sin(|m|*phi) and cos(|m|*phi)
252 r = sqrt(dx*dx + dy*dy)
253 if(r > tiny) then
254 cosphi = dx/r
255 sinphi = dy/r
256 else ! In this case the values are ill defined so we choose \phi=0
257 cosphi = m_zero
258 sinphi = m_zero
259 end if
260
261 cosm = m_one
262 sinm = m_zero
263 do i = 1, abs(mi)
264 cosmm1 = cosm
265 sinmm1 = sinm
266 cosm = cosmm1*cosphi - sinmm1*sinphi
267 sinm = cosmm1*sinphi + sinmm1*cosphi
268 end do
269
270 !And now ylm
271 ylm = plm*cmplx(cosm, sinm, real64)
272
273 if (mi < 0) then
274 ylm = conjg(ylm)
275 do i = 1, abs(mi)
276 ylm = -ylm
277 end do
278 end if
279
280 end subroutine ylmr_cmplx
281
282
283 ! ---------------------------------------------------------
289 subroutine ylmr_real(xx, li, mi, ylm)
290 real(real64), intent(in) :: xx(3)
291 integer, intent(in) :: li, mi
292 real(real64), intent(out) :: ylm
293
294 integer, parameter :: lmaxd = 20
295 real(real64), parameter :: tiny = 1.e-15_real64
296 integer :: i, ilm0, l, m, mabs
297 integer, save :: lmax = -1
298
299 real(real64) :: cmi, cosm, cosmm1, cosphi, dphase, dplg, fac, &
300 fourpi, plgndr, phase, pll, pmm, pmmp1, sinm, &
301 sinmm1, sinphi, rsize, rx, ry, rz, xysize
302 real(real64), save :: c(0:(lmaxd+1)*(lmaxd+1))
303
304 ! no push_sub: called too frequently
305
306 ! evaluate normalization constants once and for all
307 if (li > lmax) then
308 fourpi = m_four*m_pi
309 do l = 0, li
310 ilm0 = l*l + l
311 do m = 0, l
312 fac = (2*l+1)/fourpi
313 do i = l - m + 1, l + m
314 fac = fac/i
315 end do
316 c(ilm0 + m) = sqrt(fac)
317 ! next line because ylm are real combinations of m and -m
318 if (m /= 0) c(ilm0 + m) = c(ilm0 + m)*sqrt(m_two)
319 c(ilm0 - m) = c(ilm0 + m)
320 end do
321 end do
322 lmax = li
323 end if
324
325 ! if l=0, no calculations are required
326 if (li == 0) then
327 ylm = c(0)
328 return
329 end if
330
331 ! if r=0, direction is undefined => make ylm=0 except for l=0
332 rsize = sqrt(xx(1)**2 + xx(2)**2 + xx(3)**2)
333 if (rsize < tiny) then
334 ylm = m_zero
335 return
336 end if
337
338 rx = xx(1)/rsize
339 ry = xx(2)/rsize
340 rz = xx(3)/rsize
341
342 ! explicit formulas for l=1 and l=2
343 if (li == 1) then
344 select case (mi)
345 case (-1)
346 ylm = (-c(1))*ry
347 case (0)
348 ylm = c(2)*rz
349 case (1)
350 ylm = (-c(3))*rx
351 end select
352 return
353 end if
354
355 if (li == 2) then
356 select case (mi)
357 case (-2)
358 ylm = c(4)*6.0_real64*rx*ry
359 case (-1)
360 ylm = (-c(5))*m_three*ry*rz
361 case (0)
362 ylm = c(6)*m_half*(m_three*rz*rz - m_one)
363 case (1)
364 ylm = (-c(7))*m_three*rx*rz
365 case (2)
366 ylm = c(8)*m_three*(rx*rx - ry*ry)
367 end select
368 return
369 end if
370
371 ! general algorithm based on routine plgndr of numerical recipes
372 mabs = abs(mi)
373 xysize = sqrt(max(rx*rx + ry*ry, tiny))
374 cosphi = rx/xysize
375 sinphi = ry/xysize
376 cosm = m_one
377 sinm = m_zero
378 do m = 1, mabs
379 cosmm1 = cosm
380 sinmm1 = sinm
381 cosm = cosmm1*cosphi - sinmm1*sinphi
382 sinm = cosmm1*sinphi + sinmm1*cosphi
383 end do
385 if (mi < 0) then
386 phase = sinm
387 dphase = mabs*cosm
388 else
389 phase = cosm
390 dphase = (-mabs)*sinm
391 end if
392
393 pmm = m_one
394 fac = m_one
395
396 if (mabs > m_zero) then
397 do i = 1, mabs
398 pmm = (-pmm)*fac*xysize
399 fac = fac + m_two
400 end do
401 end if
402
403 if (li == mabs) then
404 plgndr = pmm
405 dplg = (-li)*rz*pmm/(xysize**2)
406 else
407 pmmp1 = rz*(2*mabs + 1)*pmm
408 if (li == mabs + 1) then
409 plgndr = pmmp1
410 dplg = -((li*rz*pmmp1 - (mabs + li)*pmm)/(xysize**2))
411 else
412 do l = mabs + 2, li
413 pll = (rz*(2*l - 1)*pmmp1 - (l + mabs - 1)*pmm)/(l - mabs)
414 pmm = pmmp1
415 pmmp1 = pll
416 end do
417 plgndr = pll
418 dplg = -((li*rz*pll - (l + mabs - 1)*pmm)/(xysize**2))
419 end if
420 end if
421
422 ilm0 = li*li + li
423 cmi = c(ilm0 + mi)
424 ylm = cmi*plgndr*phase
425
426 end subroutine ylmr_real
427
428
429 ! ---------------------------------------------------------
439 subroutine weights(N, M, cc, side)
440 integer, intent(in) :: n, m
441 real(real64), intent(out) :: cc(0:,0:,0:)
442 integer, optional, intent(in) :: side
443
444 integer :: i, j, k, mn, side_
445 real(real64) :: c1, c2, c3, c4, c5, xi
446 real(real64), allocatable :: x(:)
447
448 push_sub(weights)
449
450 safe_allocate(x(0:m))
451
452 if (present(side)) then
453 side_ = side
454 else
455 side_ = 0
456 end if
457
458 select case (side_)
459 case (-1)
460 ! grid-points for left-side finite-difference formulas on an equi-spaced grid
461 mn = m
462 x(:) = (/(-i,i=0,mn)/)
463 case (+1)
464 ! grid-points for right-side finite-difference formulas on an equi-spaced grid
465 mn = m
466 x(:) = (/(i,i=0,mn)/)
467 case default
468 ! grid-points for centered finite-difference formulas on an equi-spaced grid
469 mn = m/2
470 x(:) = (/0,(-i,i,i=1,mn)/)
471 end select
472
473
474
475 xi = m_zero ! point at which the approx. are to be accurate
476
477 cc = m_zero
478 cc(0,0,0) = m_one
479
480 c1 = m_one
481 c4 = x(0) - xi
482
483 do j = 1, m
484 mn = min(j,n)
485 c2 = m_one
486 c5 = c4
487 c4 = x(j) - xi
488
489 do k = 0, j - 1
490 c3 = x(j) - x(k)
491 c2 = c2*c3
492
493 if (j <= n) cc(k, j - 1, j) = m_zero
494 cc(k, j, 0) = c4*cc(k, j - 1, 0)/c3
495
496 do i = 1, mn
497 cc(k, j, i) = (c4*cc(k, j - 1, i) - i*cc(k, j - 1, i - 1))/c3
498 end do
499
500 cc(j, j, 0) = -c1*c5*cc(j - 1, j - 1, 0) / c2
501 end do
502
503 do i = 1, mn
504 cc(j, j, i) = c1*(i*cc(j - 1, j - 1, i - 1) - c5*cc(j - 1, j - 1, i))/c2
505 end do
506
507 c1 = c2
508 end do
509
510 safe_deallocate_a(x)
511 pop_sub(weights)
512 end subroutine weights
513
514 ! ---------------------------------------------------------
515 real(real64) pure function ddelta(i, j)
516 integer, intent(in) :: i
517 integer, intent(in) :: j
518
519 ! no push_sub in pure function
520
521 if (i == j) then
522 ddelta = m_one
523 else
524 ddelta = m_zero
525 end if
526
527 end function ddelta
528
529 ! ---------------------------------------------------------
530 subroutine interpolation_coefficients(nn, xa, xx, cc)
531 integer, intent(in) :: nn
532 real(real64), intent(in) :: xa(:)
533 real(real64), intent(in) :: xx
534 real(real64), intent(out) :: cc(:)
535
536 integer :: ii, kk
537
538 ! no push_sub, called too frequently
539
540 do ii = 1, nn
541 cc(ii) = m_one
542 do kk = 1, nn
543 if (kk == ii) cycle
544 cc(ii) = cc(ii)*(xx - xa(kk))/(xa(ii) - xa(kk))
545 end do
546 end do
547
548 end subroutine interpolation_coefficients
549
550
551 ! ---------------------------------------------------------
553 logical pure function even(n)
554 integer, intent(in) :: n
555
556 even = (mod(n, 2) == 0)
557
558 end function even
559
560
561 ! ---------------------------------------------------------
563 logical pure function odd(n)
564 integer, intent(in) :: n
565
566 odd = .not. even(n)
567
568 end function odd
569
570 ! ---------------------------------------------------------
573 subroutine cartesian2hyperspherical(x, u)
574 real(real64), intent(in) :: x(:)
575 real(real64), intent(out) :: u(:)
576
577 integer :: n, k, j
578 real(real64) :: sumx2
579 real(real64), allocatable :: xx(:)
580
582
583 n = size(x)
584 assert(n>1)
585 assert(size(u) == n-1)
586
587 ! These lines make the code less machine-dependent.
588 safe_allocate(xx(1:n))
589 do j = 1, n
590 if (abs(x(j))<1.0e-8_real64) then
591 xx(j) = m_zero
592 else
593 xx(j) = x(j)
594 end if
595 end do
596
597 u = m_zero
598 do k = 1, n-1
599 sumx2 = m_zero
600 do j = k+1, n
601 sumx2 = sumx2 + xx(j)**2
602 end do
603 if (abs(sumx2) <= m_epsilon .and. abs(xx(k)) <= m_epsilon) exit
604 if (k < n - 1) then
605 u(k) = atan2(sqrt(sumx2), xx(k))
606 else
607 u(n-1) = atan2(xx(n), xx(n-1))
608 end if
609 end do
611
613 end subroutine cartesian2hyperspherical
614 ! ---------------------------------------------------------
615
616
617 ! ---------------------------------------------------------
620 subroutine hyperspherical2cartesian(u, x)
621 real(real64), intent(in) :: u(:)
622 real(real64), intent(out) :: x(:)
623
624 integer :: n, j, k
627
628 n = size(x)
629 assert(n>1)
630 assert(size(u) == n-1)
631
632 if (n == 2) then
633 x(1) = cos(u(1))
634 x(2) = sin(u(1))
635 elseif (n == 3) then
636 x(1) = cos(u(1))
637 x(2) = sin(u(1))*cos(u(2))
638 x(3) = sin(u(1))*sin(u(2))
639 else
640 x(1) = cos(u(1))
641 x(2) = sin(u(1))*cos(u(2))
642 x(3) = sin(u(1))*sin(u(2))*cos(u(3))
643 do j = 4, n - 1
644 x(j) = m_one
645 do k = 1, j - 1
646 x(j) = x(j) * sin(u(k))
647 end do
648 x(j) = x(j) * cos(u(j))
649 end do
650 x(n) = m_one
651 do k = 1, n - 2
652 x(n) = x(n) * sin(u(k))
653 end do
654 x(n) = x(n) * sin(u(n-1))
655 end if
656
659 ! ---------------------------------------------------------
660
661
662 ! ---------------------------------------------------------
666 subroutine hypersphere_grad_matrix(grad_matrix, r, x)
667 real(real64), intent(out) :: grad_matrix(:,:)
668 real(real64), intent(in) :: r
669 real(real64), intent(in) :: x(:)
670
671 integer :: n, l, m
672
674
675 n = size(x)+1 ! the dimension of the matrix is (n-1)x(n)
676
677 ! --- l=1 ---
678 grad_matrix = m_one
679 grad_matrix(1,1) = -r*sin(x(1))
680 do m = 2, n-1
681 grad_matrix(m,1) = m_zero
682 end do
683
684 ! --- l=2..(n-1) ---
685 do l=2, n-1
686 do m=1, n-1
687 if (m == l) then
688 grad_matrix(m,l) = -r*grad_matrix(m,l)*product(sin(x(1:l)))
689 elseif (m < l) then
690 grad_matrix(m,l) = grad_matrix(m,l)*r*cos(x(m))*cos(x(l))*product(sin(x(1:m-1)))*product(sin(x(m+1:l-1)))
691 else
692 grad_matrix(m,l) = m_zero
693 end if
694 end do
695 end do
696
697 ! --- l=n ---
698 do m=1, n-1
699 grad_matrix(m,n) = r*cos(x(m))*grad_matrix(m,n)*product(sin(x(1:m-1)))*product(sin(x(m+1:n-1)))
700 end do
701
703 end subroutine hypersphere_grad_matrix
704
705 ! ---------------------------------------------------------
706 integer(int64) pure function pad88(size, blk)
707 integer(int64), intent(in) :: size
708 integer(int64), intent(in) :: blk
709
710 integer(int64) :: mm
711
712 mm = mod(size, blk)
713 if (mm == 0) then
714 pad88 = size
715 else
716 pad88 = size + blk - mm
717 end if
718 end function pad88
719
720 ! ---------------------------------------------------------
721 integer(int64) pure function pad48(size, blk)
722 integer, intent(in) :: size
723 integer(int64), intent(in) :: blk
724
725 pad48 = pad88(int(size, int64), blk)
726 end function pad48
727
728 ! ---------------------------------------------------------
729 integer(int64) pure function pad8(size, blk)
730 integer(int64), intent(in) :: size
731 integer, intent(in) :: blk
732
733 pad8 = pad88(size, int(blk, int64))
734 end function pad8
735
736 ! ---------------------------------------------------------
737 integer pure function pad4(size, blk)
738 integer, intent(in) :: size
739 integer, intent(in) :: blk
740
741 pad4 = int(pad88(int(size, int64), int(blk, int64)), int32)
742 end function pad4
743
744 ! ---------------------------------------------------------
745
750 integer pure function pad_pow2(size)
751 integer, intent(in) :: size
752
753 integer :: mm, mm2
754
755 mm = size
756 pad_pow2 = 1
757
758 ! loop below never terminates otherwise! just pick 1 as value.
759 if (size == 0) return
760
761 ! first we divide by two all the times we can, so we catch exactly
762 ! the case when size is already a power of two
763 do
764 mm2 = mm/2
765 if (mm - 2*mm2 /= 0) exit
767 mm = mm2
768 end do
769
770 ! the rest is handled by log2
771 if (mm /= 1) then
772 pad_pow2 = pad_pow2*2**(ceiling(log(real(mm, real64) )/log(2.0_8)))
773 end if
774
775 end function pad_pow2
776
777 ! -------------------------------------------------------
778
779 real(real64) pure function dlog2(xx)
780 real(real64), intent(in) :: xx
781
782 dlog2 = log(xx)/log(m_two)
783 end function dlog2
784
785 ! -------------------------------------------------------
786
787 integer pure function ilog2(xx)
788 integer, intent(in) :: xx
789
790 ilog2 = nint(log2(real(xx, real64) ))
791 end function ilog2
792
793 ! -------------------------------------------------------
794
795 integer(int64) pure function llog2(xx)
796 integer(int64), intent(in) :: xx
797
798 llog2 = nint(log2(real(xx, real64) ), kind=int64)
799 end function llog2
800
801 ! -------------------------------------------------------
803 complex(real64) pure function exponential(z)
804 complex(real64), intent(in) :: z
805
806 exponential = exp(z)
807 end function exponential
808
809 ! -------------------------------------------------------
815 complex(real64) pure function phi1(z)
816 complex(real64), intent(in) :: z
817
818 real(real64), parameter :: cut = 0.002_real64
819
820 if (abs(z) < cut) then
821 phi1 = m_one + m_half*z*(m_one + z/m_three*(m_one + m_fourth*z*(m_one + z/m_five)))
822 else
823 phi1 = (exp(z) - m_z1)/z
824 end if
825 end function phi1
826
827 ! -------------------------------------------------------
834 complex(real64) pure function phi2(z)
835 complex(real64), intent(in) :: z
836
837 real(real64), parameter :: cut = 0.002_real64
838
839 if (abs(z) < cut) then
840 phi2 = m_half + z/6.0_real64*(m_one + m_fourth*z*(m_one + z/m_five*(m_one + z/6.0_real64)))
841 else
842 phi2 = (exp(z) - z - m_z1)/z**2
843 end if
844 end function phi2
846 ! -------------------------------------------------------
848 real(real64) pure function square_root(x)
849 real(real64), intent(in) :: x
850
851 square_root = sqrt(x)
852 end function square_root
853
854 ! -------------------------------------------------------
855
856 logical function is_prime(n)
857 integer, intent(in) :: n
858
859 integer :: i, root
860
861 push_sub(is_prime)
862
863 if (n < 1) then
864 message(1) = "Internal error: is_prime does not take negative numbers."
865 call messages_fatal(1)
866 end if
867 if (n == 1) then
868 is_prime = .false.
869 pop_sub(is_prime)
870 return
871 end if
872
873 root = nint(sqrt(real(n, real64) ))
874 do i = 2, root
875 if (mod(n,i) == 0) then
876 is_prime = .false.
877 pop_sub(is_prime)
878 return
879 end if
880 end do
881
883 pop_sub(is_prime)
884 end function is_prime
885
886 ! ---------------------------------------------------------
891 subroutine generate_rotation_matrix(R, ff, tt)
892 real(real64), intent(out) :: r(:,:)
893 real(real64), intent(in) :: ff(:)
894 real(real64), intent(in) :: tt(:)
895
896 integer :: dim, i, j
897 real(real64) :: th, uv, uu, vv, ft
898 real(real64), allocatable :: axis(:), u(:), v(:), f(:), t(:), p(:)
899
901
902 dim = size(ff,1)
903
904 assert((dim < 3) .or. (dim > 2))
905 assert(size(tt,1) == dim)
906 assert((size(r,1) == dim) .and. (size(r,2) == dim))
907
908 safe_allocate(u(1:dim))
909 safe_allocate(v(1:dim))
910 safe_allocate(f(1:dim))
911 safe_allocate(t(1:dim))
912
913
914 !normalize
915 f = ff / norm2(ff)
916 t = tt / norm2(tt)
917
918 ft = dot_product(f,t)
919
920 if (abs(ft) < m_one) then
921 select case (dim)
922 case (2)
923 th = acos(ft)
924 r(1,1) = cos(th)
925 r(1,2) = -sin(th)
926
927 r(2,1) = sin(th)
928 r(2,2) = cos(th)
930 case (3)
931 if (.false.) then
932 !Old implementation
933 safe_allocate(axis(1:dim))
934 th = acos(ft)
935
936 u = f / dot_product(f,f)
937 v = t /dot_product(t,t)
938
939 axis = dcross_product(u,v)
940 axis = axis / norm2(axis)
941
942 r(1,1) = cos(th) + axis(1)**2 * (1 - cos(th))
943 r(1,2) = axis(1)*axis(2)*(1-cos(th)) + axis(3)*sin(th)
944 r(1,3) = axis(1)*axis(3)*(1-cos(th)) - axis(2)*sin(th)
945
946 r(2,1) = axis(2)*axis(1)*(1-cos(th)) - axis(3)*sin(th)
947 r(2,2) = cos(th) + axis(2)**2 * (1 - cos(th))
948 r(2,3) = axis(2)*axis(3)*(1-cos(th)) + axis(1)*sin(th)
949
950 r(3,1) = axis(3)*axis(1)*(1-cos(th)) + axis(2)*sin(th)
951 r(3,2) = axis(3)*axis(2)*(1-cos(th)) - axis(1)*sin(th)
952 r(3,3) = cos(th) + axis(3)**2 * (1 - cos(th))
953
954 safe_deallocate_a(axis)
955 end if
956
957 if (.true.) then
958 !Naive implementation
959 th = acos(ft)
960 u = dcross_product(f,t)
961 u = u / norm2(u)
962
963 r(1,1) = u(1)**2 + (1-u(1)**2)*cos(th)
964 r(1,2) = u(1)*u(2)*(1-cos(th)) - u(3)*sin(th)
965 r(1,3) = u(1)*u(3) + u(2)*sin(th)
966
967 r(2,1) = u(1)*u(2)*(1-cos(th)) + u(3)*sin(th)
968 r(2,2) = u(2)**2 + (1-u(2)**2)*cos(th)
969 r(2,3) = u(2)*u(3)*(1-cos(th)) - u(1)*sin(th)
970
971 r(3,1) = u(1)*u(3)*(1-cos(th)) - u(2)*sin(th)
972 r(3,2) = u(2)*u(3)*(1-cos(th)) + u(1)*sin(th)
973 r(3,3) = u(3)**2 + (1-u(3)**2)*cos(th)
974 end if
975
976 if (.false.) then
977 !Fast
978 safe_allocate(p(1:dim))
979
980 if (abs(f(1)) <= abs(f(2)) .and. abs(f(1)) < abs(f(3))) then
981 p = (/m_one, m_zero, m_zero/)
982 else if (abs(f(2)) < abs(f(1)) .and. abs(f(2)) <= abs(f(3))) then
983 p = (/m_zero, m_one, m_zero/)
984 else if (abs(f(3)) <= abs(f(1)) .and. abs(f(3)) < abs(f(2))) then
985 p = (/m_zero, m_zero, m_one/)
986 end if
987
988 u = p - f
989 v = p - t
990
991 uu = dot_product(u,u)
992 vv = dot_product(v,v)
993 uv = dot_product(u,v)
994
995 do i = 1,3
996 do j = 1,3
997
998 r(i,j) = ddelta(i,j) - m_two * u(i)*u(j)/uu - m_two * v(i)*v(j)/vv &
999 + 4_real64*uv * v(i)*u(j) /(uu*vv)
1000
1001 end do
1002 end do
1003
1004 safe_deallocate_a(p)
1005 end if
1006
1007
1008 end select
1009
1010 else
1011
1012 r = m_zero
1013 do i = 1,dim
1014 r(i,i) = m_one
1015 end do
1016
1017 end if
1018
1019
1020 safe_deallocate_a(u)
1021 safe_deallocate_a(v)
1022
1024 end subroutine generate_rotation_matrix
1025
1026 ! ---------------------------------------------------------
1035 subroutine numder_ridders(x, h, res, err, f)
1036 implicit none
1037 real(real64), intent(in) :: x, h
1038 real(real64), intent(out) :: res, err
1039 interface
1040 subroutine f(x, fx)
1041 use, intrinsic :: iso_fortran_env
1042 implicit none
1043 real(real64), intent(in) :: x
1044 real(real64), intent(inout) :: fx
1045 end subroutine f
1046 end interface
1047
1048 real(real64), parameter :: con = 1.4_real64, big = 1.0e30_real64, safe = m_two
1049 integer, parameter :: ntab = 20
1050 integer :: i, j
1051 real(real64) :: errt, fac, hh, fx1, fx2
1052 real(real64), allocatable :: a(:, :)
1053
1054 push_sub(numder_ridders)
1055
1056 if (abs(h) <= m_epsilon) then
1057 message(1) = "h must be nonzero in numder_ridders"
1058 call messages_fatal(1)
1059 end if
1060
1061 safe_allocate(a(1:ntab, 1:ntab))
1062
1063 hh = h
1064 call f(x+hh, fx1)
1065 call f(x-hh, fx2)
1066 a(1,1) = (fx1-fx2) / (m_two*hh)
1067 err = big
1068 do i = 2, ntab
1069 hh = hh / con
1070 call f(x+hh, fx1)
1071 call f(x-hh, fx2)
1072 a(1,i) = (fx1-fx2) / (m_two*hh)
1073 fac = con**2
1074 do j = 2, i
1075 a(j,i) = (a(j-1,i)*fac-a(j-1,i-1)) / (fac-m_one)
1076 fac = con**2*fac
1077 errt = max(abs(a(j,i)-a(j-1,i)),abs(a(j,i)-a(j-1,i-1)))
1078 if (errt .le. err) then
1079 err = errt
1080 res = a(j,i)
1081 end if
1082 end do
1083 if (abs(a(i,i)-a(i-1,i-1)) .ge. safe*err) return
1084 end do
1085
1086 safe_deallocate_a(a)
1087 pop_sub(numder_ridders)
1088 end subroutine numder_ridders
1089
1090 ! ---------------------------------------------------------
1091 pure function dzcross_product(a, b) result(c)
1092 real(real64), intent(in) :: a(:)
1093 complex(real64), intent(in) :: b(:)
1094
1095 complex(real64) :: c(1:3)
1096
1097 c(1) = a(2)*b(3) - a(3)*b(2)
1098 c(2) = a(3)*b(1) - a(1)*b(3)
1099 c(3) = a(1)*b(2) - a(2)*b(1)
1100
1101 end function dzcross_product
1102
1103 ! ---------------------------------------------------------
1104 pure function zdcross_product(a, b) result(c)
1105 complex(real64), intent(in) :: a(:)
1106 real(real64), intent(in) :: b(:)
1107
1108 complex(real64) :: c(1:3)
1109
1110 c(1) = a(2)*b(3) - a(3)*b(2)
1111 c(2) = a(3)*b(1) - a(1)*b(3)
1112 c(3) = a(1)*b(2) - a(2)*b(1)
1113
1114 end function zdcross_product
1115
1116
1117!*****************************************************************************80
1118!
1119!! LM_POLYNOMIAL evaluates Laguerre polynomials Lm(n,m,x).
1120!
1121! First terms:
1122!
1123! M = 0
1124!
1125! Lm(0,0,X) = 1
1126! Lm(1,0,X) = -X + 1
1127! Lm(2,0,X) = X^2 - 4 X + 2
1128! Lm(3,0,X) = -X^3 + 9 X^2 - 18 X + 6
1129! Lm(4,0,X) = X^4 - 16 X^3 + 72 X^2 - 96 X + 24
1130! Lm(5,0,X) = -X^5 + 25 X^4 - 200 X^3 + 600 X^2 - 600 x + 120
1131! Lm(6,0,X) = X^6 - 36 X^5 + 450 X^4 - 2400 X^3 + 5400 X^2 - 4320 X + 720
1132!
1133! M = 1
1134!
1135! Lm(0,1,X) = 0
1136! Lm(1,1,X) = -1,
1137! Lm(2,1,X) = 2 X - 4,
1138! Lm(3,1,X) = -3 X^2 + 18 X - 18,
1139! Lm(4,1,X) = 4 X^3 - 48 X^2 + 144 X - 96
1140!
1141! M = 2
1142!
1143! Lm(0,2,X) = 0
1144! Lm(1,2,X) = 0,
1145! Lm(2,2,X) = 2,
1146! Lm(3,2,X) = -6 X + 18,
1147! Lm(4,2,X) = 12 X^2 - 96 X + 144
1148!
1149! M = 3
1150!
1151! Lm(0,3,X) = 0
1152! Lm(1,3,X) = 0,
1153! Lm(2,3,X) = 0,
1154! Lm(3,3,X) = -6,
1155! Lm(4,3,X) = 24 X - 96
1156!
1157! M = 4
1158!
1159! Lm(0,4,X) = 0
1160! Lm(1,4,X) = 0
1161! Lm(2,4,X) = 0
1162! Lm(3,4,X) = 0
1163! Lm(4,4,X) = 24
1164!
1165! Recursion:
1166!
1167! Lm(0,M,X) = 1
1168! Lm(1,M,X) = (M+1-X)
1169!
1170! if 2 <= N:
1171!
1172! Lm(N,M,X) = ( (M+2*N-1-X) * Lm(N-1,M,X)
1173! + (1-M-N) * Lm(N-2,M,X) ) / N
1174!
1175! Special values:
1176!
1177! For M = 0, the associated Laguerre polynomials Lm(N,M,X) are equal
1178! to the Laguerre polynomials L(N,X).
1179!
1180! Licensing:
1181!
1182! This code is distributed under the GNU LGPL license.
1183!
1184! Modified:
1185!
1186! 08 February 2003
1187!
1188! Author:
1189!
1190! John Burkardt
1191!
1192! Reference:
1193!
1194! Milton Abramowitz, Irene Stegun,
1195! Handbook of Mathematical Functions,
1196! National Bureau of Standards, 1964,
1197! ISBN: 0-486-61272-4,
1198! LC: QA47.A34.
1200! Parameters:
1201!
1202! Input, integer ( kind = 4 ) MM, the number of evaluation points.
1203!
1204! Input, integer ( kind = 4 ) N, the highest order polynomial to compute.
1205! Note that polynomials 0 through N will be computed.
1206!
1207! Input, integer ( kind = 4 ) M, the parameter. M must be nonnegative.
1208!
1209! Input, real ( kind = rk ) X(MM), the evaluation points.
1210!
1211! Output, real ( kind = rk ) CX(MM,0:N), the associated Laguerre polynomials
1212! of degrees 0 through N evaluated at the evaluation points.
1213!
1214! Taken from
1215! https://people.sc.fsu.edu/~jburkardt/f_src/laguerre_polynomial/laguerre_polynomial.html
1216! and adapted to Octopus by N. Tancogne-Dejean
1217 subroutine generalized_laguerre_polynomial ( np, nn, mm, xx, cx )
1218 integer, intent(in) :: np
1219 integer, intent(in) :: nn, mm
1220 real(real64), intent(in) :: xx(np)
1221 real(real64), intent(out) :: cx(np)
1222
1223 integer ii
1224 real(real64), allocatable :: cx_tmp(:,:)
1225
1226 assert(mm >= 0)
1227
1228 safe_allocate(cx_tmp(1:np, 0:nn))
1229
1230 if (nn < 0) then
1231 cx = m_zero
1232 return
1233 end if
1234
1235 if (nn == 0) then
1236 cx(1:np) = m_one
1237 return
1238 end if
1239
1240 cx_tmp(1:np,0) = m_one
1241 cx_tmp(1:np,1) = real(mm + 1, real64) - xx(1:np)
1242
1243 do ii = 2, nn
1244 cx_tmp(1:np, ii) = &
1245 (( real(mm + 2*ii - 1, real64) - xx(1:np)) * cx_tmp(1:np,ii-1) &
1246 + real(-mm -ii + 1, real64) * cx_tmp(1:np,ii-2)) / real(ii, real64)
1247 end do
1248
1249 cx(1:np) = cx_tmp(1:np, nn)
1250
1251 safe_deallocate_a(cx_tmp)
1252
1253 end subroutine generalized_laguerre_polynomial
1254
1255 subroutine dupper_triangular_to_hermitian(nn, aa)
1256 integer, intent(in) :: nn
1257 real(real64), intent(inout) :: aa(:, :)
1258
1259 integer :: ii, jj
1260
1261 do jj = 1, nn
1262 do ii = 1, jj-1
1263 aa(jj, ii) = aa(ii, jj)
1264 end do
1265 end do
1266 end subroutine dupper_triangular_to_hermitian
1267
1268 subroutine zupper_triangular_to_hermitian(nn, aa)
1269 integer, intent(in) :: nn
1270 complex(real64), intent(inout) :: aa(:, :)
1271
1272 integer :: ii, jj
1273
1274 do jj = 1, nn
1275 do ii = 1, jj-1
1276 aa(jj, ii) = conjg(aa(ii, jj))
1277 end do
1278 aa(jj, jj) = real(aa(jj, jj), real64)
1279 end do
1280 end subroutine zupper_triangular_to_hermitian
1281
1282 subroutine dlower_triangular_to_hermitian(nn, aa)
1283 integer, intent(in) :: nn
1284 real(real64), intent(inout) :: aa(:, :)
1285
1286 integer :: ii, jj
1287
1288 do jj = 1, nn
1289 do ii = 1, jj-1
1290 aa(ii, jj) = aa(jj, ii)
1291 end do
1292 end do
1293 end subroutine dlower_triangular_to_hermitian
1294
1295 subroutine zlower_triangular_to_hermitian(nn, aa)
1296 integer, intent(in) :: nn
1297 complex(real64), intent(inout) :: aa(:, :)
1298
1299 integer :: ii, jj
1300
1301 do jj = 1, nn
1302 do ii = 1, jj-1
1303 aa(ii, jj) = conjg(aa(jj, ii))
1304 end do
1305 aa(jj, jj) = real(aa(jj, jj), real64)
1306 end do
1307 end subroutine zlower_triangular_to_hermitian
1308
1310 subroutine dsymmetrize_matrix(nn, aa)
1311 integer, intent(in) :: nn
1312 real(real64), intent(inout) :: aa(:, :)
1313
1314 integer :: ii, jj
1315
1316 do ii = 1, nn
1317 do jj = ii+1, nn
1318 aa(ii, jj) = m_half * (aa(ii, jj) + aa(jj, ii))
1319 aa(jj, ii) = aa(ii, jj)
1320 end do
1321 end do
1322 end subroutine dsymmetrize_matrix
1323
1325 subroutine zsymmetrize_matrix(n, a)
1326 integer, intent(in) :: n
1327 complex(real64), intent(inout) :: a(:, :)
1328
1329 integer :: i,j
1330
1331 do i = 1, n
1332 do j = i+1, n
1333 a(i, j) = m_half * (a(i, j) + conjg(a(j, i)))
1334 a(j, i) = conjg(a(i, j))
1335 end do
1336 end do
1337 end subroutine zsymmetrize_matrix
1338
1339
1340 subroutine dzero_small_elements_matrix(aa, tol)
1341 real(real64), intent(inout) :: aa(:, :)
1342 real(real64) :: tol
1343
1344 where(abs(aa)<tol)
1345 aa = m_zero
1346 end where
1347 end subroutine dzero_small_elements_matrix
1348
1350 logical function is_symmetric(a, tol)
1351 real(real64), intent(in) :: a(:, :)
1352 real(real64), optional, intent(in) :: tol
1353
1354 real(real64) :: tolerance
1355
1356 push_sub(is_symmetric)
1357
1358 tolerance = optional_default(tol, 1.0e-8_real64)
1359 ! a must be square
1360 assert(size(a, 1) == size(a, 2))
1361
1362 is_symmetric = all(abs(a - transpose(a)) < tolerance)
1364 pop_sub(is_symmetric)
1365
1366 end function is_symmetric
1367
1369 logical pure function is_diagonal(dim, matrix, tol)
1370 integer, intent(in) :: dim
1371 real(real64), intent(in) :: matrix(:, :)
1372 real(real64), intent(in) :: tol
1373
1374 integer :: ii, jj
1375 real(real64) :: max_diag
1376
1377 max_diag = m_zero
1378 do ii = 1, dim
1379 max_diag = max(abs(matrix(ii, ii)), max_diag)
1380 end do
1381
1382 is_diagonal = .true.
1383 do ii = 1, dim
1384 do jj = 1, dim
1385 if (ii == jj) cycle
1386 if (abs(matrix(ii, jj)) > tol*max_diag) is_diagonal = .false.
1387 end do
1388 end do
1389
1390 end function is_diagonal
1391
1393 subroutine convert_to_base(num, base, converted)
1394 integer, intent(in) :: num
1395 integer, intent(in) :: base
1396 integer, intent(out) :: converted(:)
1397
1398 integer :: i, num_, length
1399
1400 ! no PUSH_SUB, called too often
1401
1402 converted = 0
1403 if (num > 0) then
1404 ! length of converted number
1405 length = floor(log(real(num, real64)) / log(real(base, real64))) + 1
1406 assert(size(converted) >= length)
1407
1408 num_ = num
1409 ! do conversion by repeated division
1410 do i = 1, size(converted)
1411 converted(i) = mod(num_, base)
1412 num_ = num_ / base
1413 end do
1414 end if
1415
1416 end subroutine convert_to_base
1417
1419 subroutine convert_from_base(converted, base, num)
1420 integer, intent(in) :: converted(:)
1421 integer, intent(in) :: base
1422 integer, intent(out) :: num
1423
1424 integer :: i
1425
1426 ! no PUSH_SUB, called too often
1427
1428 num = m_zero
1429 do i = 1, size(converted)
1430 num = num + converted(i) * base**(i-1)
1431 end do
1432
1433 end subroutine convert_from_base
1434
1435#include "undef.F90"
1436#include "complex.F90"
1437#include "math_inc.F90"
1438
1439#include "undef.F90"
1440#include "real.F90"
1441#include "math_inc.F90"
1442
1443end module math_oct_m
1444
1445!! Local Variables:
1446!! mode: f90
1447!! coding: utf-8
1448!! End:
This is the common interface to a simple-minded polynomical interpolation procedure (simple use of th...
Definition: math.F90:189
double acos(double __x) __attribute__((__nothrow__
double log(double __x) __attribute__((__nothrow__
double exp(double __x) __attribute__((__nothrow__
double sin(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
double cos(double __x) __attribute__((__nothrow__
double atan2(double __y, double __x) __attribute__((__nothrow__
double floor(double __x) __attribute__((__nothrow__
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
subroutine zsymmetrize_matrix(n, a)
make a complex matrix Hermitian
Definition: math.F90:1421
complex(real64) pure function, public phi2(z)
Compute phi2(z) = (phi1(z)-1)/z = (exp(z) - z - 1)/z^2.
Definition: math.F90:930
integer pure function ilog2(xx)
Definition: math.F90:883
logical pure function, public even(n)
Returns if n is even.
Definition: math.F90:649
pure real(real64) function, dimension(dim, dim) ddiagonal_matrix(dim, diag)
Currently only returns a matrix whose diagonal elements are all the same. Note that the integer versi...
Definition: math.F90:1792
subroutine, public generalized_laguerre_polynomial(np, nn, mm, xx, cx)
Definition: math.F90:1313
complex(real64) pure function, public exponential(z)
Wrapper for exponential.
Definition: math.F90:899
subroutine dinterpolate_0(xa, ya, x, y)
Definition: math.F90:1870
subroutine, public ylmr_cmplx(xx, li, mi, ylm)
Computes spherical harmonics ylm at position (x, y, z)
Definition: math.F90:314
subroutine zinterpolate_0(xa, ya, x, y)
Definition: math.F90:1679
subroutine, public convert_to_base(num, base, converted)
convert an integer in base 10 to a different base
Definition: math.F90:1489
subroutine dinterpolate_1(xa, ya, x, y)
Definition: math.F90:1842
real(real64) pure function dlog2(xx)
Definition: math.F90:875
integer(int64) pure function pad88(size, blk)
Definition: math.F90:802
real(real64) pure function, public square_root(x)
Wrapper for sqrt.
Definition: math.F90:944
integer pure function pad4(size, blk)
Definition: math.F90:833
subroutine dlower_triangular_to_hermitian(nn, aa)
Definition: math.F90:1378
pure complex(real64) function, dimension(1:3), public dzcross_product(a, b)
Definition: math.F90:1187
subroutine, public cartesian2hyperspherical(x, u)
Performs a transformation of an n-dimensional vector from Cartesian coordinates to hyperspherical coo...
Definition: math.F90:669
integer pure function, public pad_pow2(size)
create array size, which is padded to powers of 2
Definition: math.F90:846
subroutine, public dzero_small_elements_matrix(aa, tol)
Definition: math.F90:1436
subroutine zinterpolate_1(xa, ya, x, y)
Definition: math.F90:1651
pure real(real64) function, dimension(1:3), public dcross_product(a, b)
Definition: math.F90:1895
integer(int64) pure function pad48(size, blk)
Definition: math.F90:817
subroutine dsymmetrize_matrix(nn, aa)
symmetrize a real matrix
Definition: math.F90:1406
subroutine, public ylmr_real(xx, li, mi, ylm)
This is a Numerical Recipes-based subroutine computes real spherical harmonics ylm at position (x,...
Definition: math.F90:385
subroutine dupper_triangular_to_hermitian(nn, aa)
Definition: math.F90:1351
subroutine, public interpolation_coefficients(nn, xa, xx, cc)
Definition: math.F90:626
logical function, public is_symmetric(a, tol)
Check if a 2D array is symmetric.
Definition: math.F90:1446
subroutine, public generate_rotation_matrix(R, ff, tt)
Generates a rotation matrix R to rotate a vector f to t.
Definition: math.F90:987
elemental logical function zis_close_scalar(x, y, rtol, atol)
Same as dis_close_scalar for complex numbers.
Definition: math.F90:247
logical pure function, public odd(n)
Returns if n is odd.
Definition: math.F90:659
real(real64) pure function, public ddelta(i, j)
Definition: math.F90:611
logical function, public is_prime(n)
Definition: math.F90:952
complex(real64) pure function, public phi1(z)
Compute phi1(z) = (exp(z)-1)/z.
Definition: math.F90:911
pure complex(real64) function, dimension(1:3), public zcross_product(a, b)
Definition: math.F90:1704
subroutine, public weights(N, M, cc, side)
Compute the weights for finite-difference calculations:
Definition: math.F90:535
subroutine zinterpolate_2(xa, ya, x, y)
Definition: math.F90:1616
pure complex(real64) function, dimension(dim, dim) zdiagonal_matrix(dim, diag)
Currently only returns a matrix whose diagonal elements are all the same. Note that the integer versi...
Definition: math.F90:1601
subroutine zlower_triangular_to_hermitian(nn, aa)
Definition: math.F90:1391
logical pure function, public is_diagonal(dim, matrix, tol)
Returns true is the matrix of size dim x dim is diagonal, given a relative tolerance.
Definition: math.F90:1465
subroutine zupper_triangular_to_hermitian(nn, aa)
Definition: math.F90:1364
pure integer function, dimension(dim, dim) idiagonal_matrix(dim, diag)
Currently only returns a matrix whose diagonal elements are all the same. Note that the real and comp...
Definition: math.F90:263
integer(int64) pure function pad8(size, blk)
Definition: math.F90:825
subroutine dinterpolate_2(xa, ya, x, y)
Definition: math.F90:1807
subroutine, public numder_ridders(x, h, res, err, f)
Numerical derivative (Ridder`s algorithm).
Definition: math.F90:1131
integer(int64) pure function llog2(xx)
Definition: math.F90:891
recursive real(real64) function, public hermite(n, x)
Definition: math.F90:278
elemental logical function dis_close_scalar(x, y, rtol, atol)
Are and equal within a tolerance.
Definition: math.F90:234
pure complex(real64) function, dimension(1:3), public zdcross_product(a, b)
Definition: math.F90:1200
subroutine, public convert_from_base(converted, base, num)
convert an integer to base 10 from a different base
Definition: math.F90:1515
subroutine, public hyperspherical2cartesian(u, x)
Performs the inverse transformation of cartesian2hyperspherical.
Definition: math.F90:716
subroutine, public hypersphere_grad_matrix(grad_matrix, r, x)
Gives the hyperspherical gradient matrix, which contains the derivatives of the Cartesian coordinates...
Definition: math.F90:762
recursive integer function, public factorial(n)
Definition: math.F90:298
static double f(double w, void *p)
int true(void)