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