Octopus
splines.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
21!/*----------------------------------------------------------------------------
22! This module contains a data type (spline_t) to contain 1D functions,
23! along with a series of procedures to manage them (to define them, to obtain
24! its values, to operate on them, etc). The internal representation of the
25! functions is done through cubic splines, handled by the GSL library. For
26! the user of the module, this internal representation is hidden; one just
27! works with what are called hereafter "spline functions".
28!
29! To define a function, one must supply a set {x(i),y(i)} of pairs of values
30! -- the abscissa and the value of the function.
31!
32! [*] DATA TYPE:
33! To define a spline function:
34! type(spline_t) :: f
35!
36! [1] INITIALIZATION:
37! Before using any function, one should initialize it:
38!
39! Interface:
40! subroutine spline_init(spl)
41! type(spline_t), intent(out) :: spl [or spl(:) or spl(:, :)]
42! end subroutine spline_init
43!
44! Usage:
45! call spline_init(f)
46!
47! [2] FINALIZATION:
48! To empty any allocated space, one should finalize the function:
49!
50! Interface:
51! subroutine spline_end(spl)
52! type(spline_t), intent(inout) :: spl [or spl(:) or spl(:, :)]
53! end subroutine spline_end
54!
55! Usage
56! type(spline_t) :: f
57! call spline_end(f)
58!
59! [3] TO DEFINE A FUNCTION:
60! To "fill" an initialized function f,use spline_fit
61!
62! Interface:
63! subroutine spline_fit(n, x, y, spl, threshold)
64! integer, intent(in) :: nrc
65! real(X), intent(in) :: x(n), y(n)
66! type(spline_t), intent(out) :: spl
67! real(real64), optional, intent(in) :: threshold
68! end subroutine spline_fit
69!
70! (X may be 4 or eight, for single or double precision)
71!
72! Usage:
73! call spline_fit(n, x, y, f)
74! n is the number of values that are supplied, x the abscissas, and y
75! the value of the function to represent at each point.
76! The optional argument threshold is used to find the radius associated with the spline.
77!
78! [4] FUNCTION VALUES:
79! To retrieve the value of a function at a given point:
80!
81! Interface:
82!
83! real(real64) function spline_eval(spl, x)
84! type(spline_t), intent(in) :: spl
85! real(real64), intent(in) :: x
86! end function spline_eval
87!
88! Usage:
89! To retrieve the value of function f at point x, and place it into
90! real value val:
91! val = spline_eval(f, x)
92!
93! [5] MULTIPLICATION BY A SCALAR
94! You may multiply a given spline-represented spline by a real number:
95!
96! Interface:
97! subroutine spline_times(a, spl)
98! type(spline_t), intent(inout) :: spl
99! real(real64), intent(in) :: a
100! real(real64), optional, intent(in) :: threshold
101! end subroutine spline_times
102!
103! Usage:
104! call spline_init(f)
105! call spline_fit(npoints, x, y, f) ! Fill f with y values at x points
106! call spline_times(a, f) ! Now f contains a*y values at x points.
107!
108! [6] INTEGRAL:
109! Given a defined function, the function spline_integral returns its
110! integral. The interval of integration may or may not be supplied.
111!
112! Interface:
113! real(real64) function spline_integral(spl [,a,b])
114! type(spline_integral), intent(in) :: spl
115! real(real64), intent(in), optional :: a, b
116! end function spline_integral
117!
118! Note: The following routines, spline_3dft, spline_cut and spline_filter,
119! assume that the spline functions are the radial part of a 3 dimensional function with
120! spherical symmetry. This is why the Fourier transform of F(\vec{r}) = f(r), is:
121! F(\vec{g}) = f(g) = \frac{4\pi}{g} \int_{0}^{\infty} { r*sin(g*r)*f(r) }
122! which coincides with the inverse Fourier transform, except that the inverse Fourier
123! transform should be multiplied by a (2*\pi)^{-3} factor.
124!
125! [7] FOURIER TRANSFORM:
126! If a spline function f is representing the radial part of a spherically
127! symmetric function F(\vec{r}), its Fourier transform is:
128! F(\vec{g}) = f(g) = \frac{4\pi}{g} \int_{0}^{\infty} { r*sin(g*r)*f(r) }
129! It is assumed that f is defined in some interval [0,rmax], and that it is
130! null at rmax and beyond. One may obtain f(g) by using spline_3dft.
131! The result is placed on the spline data type splw. This has to be initialized,
132! and may or may not be filled. If it is filled, the abscissas that it has
133! are used to define the function. If it is not filled, an equally spaced
134! grid is constructed to define the function, in the interval [0, gmax], where
135! gmax has to be supplied by the caller.
136!
137! Interface:
138! subroutine spline_3dft(spl, splw, gmax)
139! type(spline_t), intent(in) :: spl
140! type(spline_t), intent(inout) :: splw
141! real(real64), intent(in), optional :: gmax
142! end subroutine spline_3dft
143!
144! [8] BESSEL TRANSFORM:
145!
146!
147! [9] CUTTING A FUNCTION:
148! spline_cut multiplies a given function by a cutoff-function, which
149! is defined to be one in [0, cutoff], and \exp\{-beta*(x/cutoff-1)^2\}
150!
151! Interface:
152! subroutine spline_cut(spl, cutoff, beta)
153! type(spline_t), intent(in) :: spl
154! real(real64), intent(in) :: cutoff, beta
155! end subroutine spline_cut
157! [10] PRINTING A FUNCTION:
158! It prints to a file the (x,y) values that were used to define a function.
159! The file is pointed to by its Fortran unit given by argument iunit.
161! Interface:
162! subroutine spline_print(spl, iunit)
163! type(spline_t), intent(in) :: spl [ or spl(:) or spl(:, :)]
164! integer, intent(in) :: iunit
165! end subroutine spline_print
167!----------------------------------------------------------------------------*/!
168module splines_oct_m
169 use debug_oct_m
170 use global_oct_m
171 use iso_c_binding
172 use, intrinsic :: iso_fortran_env
176
177 implicit none
178
179 ! Define which routines can be seen from the outside.
180 private
181 public :: &
182 spline_t, & ! [*]
183 spline_init, & ! [1]
184 spline_end, & ! [2]
185 spline_fit, & ! [3]
186 spline_eval, & ! [4]
187 spline_eval_vec, & ! [4]
188 spline_times, & ! [5]
189 spline_integral, & ! [6]
190 spline_3dft, & ! [7]
192 spline_cut, & ! [9]
193 spline_print, & ! [10]
194 spline_der, &
195 spline_der2, &
196 spline_div, &
197 spline_mult, &
203
205 type spline_t
206 private
207 real(real64) :: x_limit(2)
208 type(c_ptr) :: spl, acc
209 real(real64), public :: x_threshold
211 end type spline_t
212
215 interface spline_eval_vec
216 module procedure spline_eval8_array
217 module procedure spline_evalz_array
218 end interface spline_eval_vec
219
222 interface spline_integral
223 module procedure spline_integral_full
224 module procedure spline_integral_limits
225 end interface spline_integral
228 interface spline_init
229 module procedure spline_init_0
230 module procedure spline_init_1
231 module procedure spline_init_2
232 end interface spline_init
233
234 interface spline_end
235 module procedure spline_end_0
236 module procedure spline_end_1
237 module procedure spline_end_2
238 end interface spline_end
239
240 interface spline_print
241 module procedure spline_print_0
242 module procedure spline_print_1
243 module procedure spline_print_2
244 end interface spline_print
245
246 ! These are interfaces to functions defined in oct_gsl_f.c, which in turn
247 ! take care of calling the GSL library.
248 interface
249 subroutine oct_spline_end(spl, acc)
250 use iso_c_binding
251 implicit none
252 type(c_ptr), intent(inout) :: spl, acc
253 end subroutine oct_spline_end
254
255 subroutine oct_spline_fit(nrc, x, y, spl, acc)
256 use iso_c_binding
257 use, intrinsic :: iso_fortran_env
258 implicit none
259 integer, intent(in) :: nrc
260 real(real64), intent(in) :: x
261 real(real64), intent(in) :: y
262 type(c_ptr), intent(inout) :: spl
263 type(c_ptr), intent(inout) :: acc
264 end subroutine oct_spline_fit
265
266 real(real64) function oct_spline_eval(x, spl, acc)
267 use iso_c_binding
268 use, intrinsic :: iso_fortran_env
269 implicit none
270 real(real64), intent(in) :: x
271 type(c_ptr), intent(in) :: spl
272 type(c_ptr), intent(in) :: acc
273 end function oct_spline_eval
274
275 subroutine oct_spline_eval_array(nn, xf, spl, acc)
276 use iso_c_binding
277 use, intrinsic :: iso_fortran_env
278 implicit none
279 integer, intent(in) :: nn
280 real(real64), intent(inout) :: xf
281 type(c_ptr), intent(in) :: spl
282 type(c_ptr), intent(in) :: acc
283 end subroutine oct_spline_eval_array
284
285 subroutine oct_spline_eval_arrayz(nn, xf, spl, acc)
286 use iso_c_binding
287 use, intrinsic :: iso_fortran_env
288 implicit none
289 integer, intent(in) :: nn
290 complex(real64), intent(inout) :: xf
291 type(c_ptr), intent(in) :: spl
292 type(c_ptr), intent(in) :: acc
293 end subroutine oct_spline_eval_arrayz
294
295 real(real64) function oct_spline_eval_der(x, spl, acc)
296 use iso_c_binding
297 use, intrinsic :: iso_fortran_env
298 implicit none
299 real(real64), intent(in) :: x
300 type(c_ptr), intent(in) :: spl
301 type(c_ptr), intent(in) :: acc
302 end function oct_spline_eval_der
303
304 real(real64) function oct_spline_eval_der2(x, spl, acc)
305 use iso_c_binding
306 use, intrinsic :: iso_fortran_env
307 implicit none
308 real(real64), intent(in) :: x
309 type(c_ptr), intent(in) :: spl
310 type(c_ptr), intent(in) :: acc
312
313 integer pure function oct_spline_npoints(spl, acc)
314 use iso_c_binding
315 implicit none
316 type(c_ptr), intent(in) :: spl
317 type(c_ptr), intent(in) :: acc
318 end function oct_spline_npoints
319
320 pure subroutine oct_spline_x(spl, x)
321 use iso_c_binding
322 use, intrinsic :: iso_fortran_env
323 implicit none
324 type(c_ptr), intent(in) :: spl
325 real(real64), intent(out) :: x
326 end subroutine oct_spline_x
327
328 subroutine oct_spline_y(spl, y)
329 use iso_c_binding
330 use, intrinsic :: iso_fortran_env
331 implicit none
332 type(c_ptr), intent(in) :: spl
333 real(real64), intent(out) :: y
334 end subroutine oct_spline_y
335
336 real(real64) pure function oct_spline_eval_integ(spl, a, b, acc)
337 use iso_c_binding
338 use, intrinsic :: iso_fortran_env
339 implicit none
340 type(c_ptr), intent(in) :: spl
341 real(real64), intent(in) :: a
342 real(real64), intent(in) :: b
343 type(c_ptr), intent(in) :: acc
344 end function oct_spline_eval_integ
346 real(real64) pure function oct_spline_eval_integ_full(spl, acc)
347 use iso_c_binding
348 use, intrinsic :: iso_fortran_env
349 implicit none
350 type(c_ptr), intent(in) :: spl
351 type(c_ptr), intent(in) :: acc
352 end function oct_spline_eval_integ_full
353 end interface
354
355 real(real64), parameter :: tol_x = 1e-14_real64
356
357contains
358
359 !------------------------------------------------------------
360 subroutine spline_init_0(spl)
361 type(spline_t), intent(out) :: spl
362
363 ! No PUSH SUB, called too often
364
365 spl%spl = c_null_ptr
366 spl%acc = c_null_ptr
367
368 ! deliberately illegal values, for checking
369 spl%x_limit(1) = -1_real64
370 spl%x_limit(2) = -2_real64
371
372 spl%x_threshold = m_zero
373
374 end subroutine spline_init_0
375
376
377 !------------------------------------------------------------
378 subroutine spline_init_1(spl)
379 type(spline_t), intent(out) :: spl(:)
380
381 integer :: i
382
383 push_sub(spline_init_1)
384
385 do i = 1, size(spl)
386 call spline_init_0(spl(i))
387 end do
388
389 pop_sub(spline_init_1)
390 end subroutine spline_init_1
391
392
393 !------------------------------------------------------------
394 subroutine spline_init_2(spl)
395 type(spline_t), intent(out) :: spl(:, :)
396
397 integer :: i, j
398
399 push_sub(spline_init_2)
400
401 do i = 1, size(spl, 1)
402 do j = 1, size(spl, 2)
403 call spline_init_0(spl(i, j))
404 end do
405 end do
406
407 pop_sub(spline_init_2)
408 end subroutine spline_init_2
409
410
411 !------------------------------------------------------------
412 subroutine spline_end_0(spl)
413 type(spline_t), intent(inout) :: spl
414
415 ! No PUSH SUB, called too often.
416
417 if (c_associated(spl%spl) .and. c_associated(spl%acc)) then
418 call oct_spline_end(spl%spl, spl%acc)
419 end if
420 spl%spl = c_null_ptr
421 spl%acc = c_null_ptr
422
423 end subroutine spline_end_0
424
425
426 !------------------------------------------------------------
427 subroutine spline_end_1(spl)
428 type(spline_t), intent(inout) :: spl(:)
429
430 integer :: i
431
432 push_sub(spline_end_1)
433
434 do i = 1, size(spl)
435 call spline_end_0(spl(i))
436 end do
437
438 pop_sub(spline_end_1)
439 end subroutine spline_end_1
441
442 !------------------------------------------------------------
443 subroutine spline_end_2(spl)
444 type(spline_t), intent(inout) :: spl(:, :)
445
446 integer :: i, j
447
448 push_sub(spline_end_2)
450 do i = 1, size(spl, 1)
451 do j = 1, size(spl, 2)
452 call spline_end_0(spl(i, j))
453 end do
454 end do
455
456 pop_sub(spline_end_2)
457 end subroutine spline_end_2
458
459
460 !------------------------------------------------------------
461 subroutine spline_fit(nrc, rofi, ffit, spl, threshold)
462 integer, intent(in) :: nrc
463 real(real64), intent(in) :: rofi(:)
464 real(real64), intent(in) :: ffit(:)
465 type(spline_t), intent(inout) :: spl
466 real(real64), optional, intent(in) :: threshold
467
468 !No PUSH SUB, called too often
469
470 assert(nrc > 0)
471
472 spl%x_limit(1) = rofi(1)
473 spl%x_limit(2) = rofi(nrc)
474 call oct_spline_fit(nrc, rofi(1), ffit(1), spl%spl, spl%acc)
475
476 if (present(threshold)) then
477 if (threshold > m_zero) then
478 spl%x_threshold = spline_x_threshold(spl, threshold)
479 else
480 spl%x_threshold = spl%x_limit(2)
481 end if
482 else
483 spl%x_threshold = m_zero
484 end if
485
486 end subroutine spline_fit
487
488 !------------------------------------------------------------
489 real(real64) function spline_eval(spl, x)
490 type(spline_t), intent(in) :: spl
491 real(real64), intent(in) :: x
492
493 spline_eval = oct_spline_eval(x, spl%spl, spl%acc)
494 end function spline_eval
495
496
497 !------------------------------------------------------------
498 subroutine spline_eval8_array(spl, nn, xf)
499 type(spline_t), intent(in) :: spl
500 integer, intent(in) :: nn
501 real(real64), intent(inout) :: xf(:)
502
503 assert(nn > 0)
504
505 call oct_spline_eval_array(nn, xf(1), spl%spl, spl%acc)
506 end subroutine spline_eval8_array
507
508
509 !------------------------------------------------------------
510 subroutine spline_evalz_array(spl, nn, xf)
511 type(spline_t), intent(in) :: spl
512 integer, intent(in) :: nn
513 complex(real64), intent(inout) :: xf(:)
514
515 assert(nn > 0)
516
517 call oct_spline_eval_arrayz(nn, xf(1), spl%spl, spl%acc)
518 end subroutine spline_evalz_array
519
520 !------------------------------------------------------------
521 subroutine spline_times(a, spl, threshold)
522 real(real64), intent(in) :: a
523 type(spline_t), intent(inout) :: spl
524 real(real64), intent(in) :: threshold
525
526 integer :: npoints, i
527 real(real64), allocatable :: x(:), y(:)
528
529 push_sub(spline_times)
530
531 npoints = oct_spline_npoints(spl%spl, spl%acc)
532 assert(npoints > 0)
533
534 safe_allocate(x(1:npoints))
535 safe_allocate(y(1:npoints))
536
537 call oct_spline_x(spl%spl, x(1))
538 call oct_spline_y(spl%spl, y(1))
539 call oct_spline_end(spl%spl, spl%acc)
540 !$omp parallel do simd
541 do i = 1, npoints
542 y(i) = a*y(i)
543 end do
544 call spline_fit(npoints, x, y, spl, threshold)
545
546 safe_deallocate_a(x)
547 safe_deallocate_a(y)
548
549 pop_sub(spline_times)
550 end subroutine spline_times
551
552
553 !------------------------------------------------------------
554 real(real64) function spline_integral_full(spl) result(res)
555 type(spline_t), intent(in) :: spl
556
557 push_sub(spline_integral_full)
558
559 res = oct_spline_eval_integ_full(spl%spl, spl%acc)
560
561 pop_sub(spline_integral_full)
562 end function spline_integral_full
563
564
565 !------------------------------------------------------------
566 real(real64) pure function spline_integral_limits(spl, a, b) result(res)
567 type(spline_t), intent(in) :: spl
568 real(real64), intent(in) :: a
569 real(real64), intent(in) :: b
570
571 res = oct_spline_eval_integ(spl%spl, a, b, spl%acc)
572 end function spline_integral_limits
573
574 !------------------------------------------------------------
575 subroutine spline_3dft(spl, splw, threshold, gmax)
576 type(spline_t), intent(in) :: spl
577 type(spline_t), intent(inout) :: splw
578 real(real64), intent(in) :: threshold
579 real(real64), optional, intent(in) :: gmax
580
581 type(spline_t) :: aux
582 real(real64) :: g, dg
583 integer :: np
584 integer :: npoints, i, j
585 real(real64), allocatable :: x(:), y(:), y2(:), xw(:), yw(:)
586
587 push_sub(spline_3dft)
588
589 npoints = oct_spline_npoints(spl%spl, spl%acc)
590 assert(npoints > 0)
591
592 safe_allocate( x(1:npoints))
593 safe_allocate( y(1:npoints))
594 safe_allocate(y2(1:npoints))
595
596 call oct_spline_x(spl%spl, x(1))
597 call oct_spline_y(spl%spl, y(1))
598
599 ! Check whether splw comes with a defined grid, or else build it.
600 if (c_associated(splw%spl)) then
601 np = oct_spline_npoints(splw%spl, splw%acc)
602 safe_allocate(xw(1:np))
603 safe_allocate(yw(1:np))
604 call oct_spline_x(splw%spl, xw(1))
605 ! But now we need to kill the input:
606 call spline_end(splw)
607 else
608 np = 200 ! hard coded value
609 dg = gmax/(np-1)
610 safe_allocate(xw(1:np))
611 safe_allocate(yw(1:np))
612 do i = 1, np
613 g = (i-1)*dg
614 xw(i) = g
615 end do
616 end if
617
618 ! The first point, xw(1) = 0.0 and it has to be treated separately.
619 !$omp parallel do simd
620 do j = 1, npoints
621 y2(j) = m_four*m_pi*y(j)*x(j)**2
622 end do
623 call spline_init(aux)
624 call spline_fit(npoints, x, y2, aux, m_zero)
625 yw(1) = oct_spline_eval_integ_full(aux%spl, aux%acc)
626 call spline_end(aux)
627
628 do i = 2, np
629 !$omp parallel do simd
630 do j = 1, npoints
631 y2(j) = (m_four*m_pi/xw(i))*y(j)*x(j)*sin(xw(i)*x(j))
632 end do
633 call spline_init(aux)
634 call spline_fit(npoints, x, y2, aux, m_zero)
635 yw(i) = oct_spline_eval_integ_full(aux%spl, aux%acc)
636 call spline_end(aux)
637 end do
638
639 call spline_init(splw)
640 call spline_fit(np, xw, yw, splw, threshold)
641
642 safe_deallocate_a(x)
643 safe_deallocate_a(y)
644 safe_deallocate_a(y2)
645 safe_deallocate_a(xw)
646 safe_deallocate_a(yw)
647
648 pop_sub(spline_3dft)
649 end subroutine spline_3dft
650
651
652 !------------------------------------------------------------
653 subroutine spline_besselft(spl, splw, l, threshold, gmax)
654 type(spline_t), intent(in) :: spl
655 type(spline_t), intent(inout) :: splw
656 integer, intent(in) :: l
657 real(real64), optional, intent(in) :: threshold
658 real(real64), optional, intent(in) :: gmax
659
660 type(spline_t) :: aux
661 real(real64) :: dg
662 integer :: np
663 integer :: npoints, i, j
664 real(real64), allocatable :: x(:), y(:), y2(:), xw(:), yw(:)
665
666 push_sub(spline_besselft)
667
668 npoints = oct_spline_npoints(spl%spl, spl%acc)
669 assert(npoints > 0)
670
671 safe_allocate( x(1:npoints))
672 safe_allocate( y(1:npoints))
673 safe_allocate(y2(1:npoints))
674
675 call oct_spline_x(spl%spl, x(1))
676 call oct_spline_y(spl%spl, y(1))
677
678 ! Check whether splw comes with a defined grid, or else build it.
679 if (c_associated(splw%spl)) then
680 np = oct_spline_npoints(splw%spl, splw%acc)
681 safe_allocate(xw(1:np))
682 safe_allocate(yw(1:np))
683 call oct_spline_x(splw%spl, xw(1))
684 ! But now we need to kill the input:
685 call spline_end(splw)
686 else
687 assert(present(gmax))
688 np = 1000 ! hard coded value
689 dg = gmax/(np-1)
690 safe_allocate(xw(1:np))
691 safe_allocate(yw(1:np))
692 !$omp parallel do
693 do i = 1, np
694 xw(i) = (i-1)*dg
695 end do
696 end if
697
698 do i = 1, np
699 !$omp parallel do
700 do j = 1, npoints
701 y2(j) = y(j) * x(j)**2 * loct_sph_bessel(l, x(j)*xw(i))
702 end do
703 call spline_init(aux)
704 call spline_fit(npoints, x, y2, aux)
705 yw(i) = sqrt(m_two/m_pi)*oct_spline_eval_integ_full(aux%spl, aux%acc)
706
707 call spline_end(aux)
708 end do
709
710 call spline_init(splw)
711 call spline_fit(np, xw, yw, splw, threshold)
712
713 safe_deallocate_a(x)
714 safe_deallocate_a(y)
715 safe_deallocate_a(y2)
716 safe_deallocate_a(xw)
717 safe_deallocate_a(yw)
718
719 pop_sub(spline_besselft)
720 end subroutine spline_besselft
721
722
723 !------------------------------------------------------------
724 subroutine spline_cut(spl, cutoff, beta, threshold)
725 type(spline_t), intent(inout) :: spl
726 real(real64), intent(in) :: cutoff
727 real(real64), intent(in) :: beta
728 real(real64), optional,intent(in) :: threshold
729
730 integer :: npoints, i
731 real(real64), allocatable :: x(:), y(:)
732 real(real64) :: exp_arg
733
734 push_sub(spline_cut)
735
736 npoints = oct_spline_npoints(spl%spl, spl%acc)
737 assert(npoints > 0)
738
739 safe_allocate(x(1:npoints))
740 safe_allocate(y(1:npoints))
741
742 call oct_spline_x(spl%spl, x(1))
743 call oct_spline_y(spl%spl, y(1))
744 call oct_spline_end(spl%spl, spl%acc)
745
746 do i = npoints, 1, -1
747 if (x(i)<cutoff) then
748 exit
749 end if
750
751 exp_arg = -beta*(x(i)/cutoff - m_one)**2
752 !To avoid underflows
753 if (exp_arg > m_min_exp_arg) then
754 y(i) = y(i) * exp(exp_arg)
755 else
756 y(i) = m_zero
757 end if
758 end do
759 call spline_fit(npoints, x, y, spl, threshold)
760
761 safe_deallocate_a(x)
762 safe_deallocate_a(y)
763
764 pop_sub(spline_cut)
765 end subroutine spline_cut
766
767
768 !------------------------------------------------------------
773 subroutine spline_div(spla, splb, threshold)
774 type(spline_t), intent(inout) :: spla
775 type(spline_t), intent(in) :: splb
776 real(real64), optional, intent(in) :: threshold
777
778 integer :: npoints, i, new_npoints
779 real(real64), allocatable :: x(:), y(:)
780 real(real64) :: aa
781
782 push_sub(spline_div)
783
784 npoints = oct_spline_npoints(spla%spl, spla%acc)
785 assert(npoints > 0)
786
787 safe_allocate(x(1:npoints))
788 safe_allocate(y(1:npoints))
789
790 call oct_spline_x(spla%spl, x(1))
791 call oct_spline_y(spla%spl, y(1))
792 call oct_spline_end(spla%spl, spla%acc)
793
794 assert(splb%x_limit(2) >= splb%x_limit(1))
795
796 ! As an output, the new spline only contains the new number of points
797 ! that goes up to the cutoff of the splb spline
798 new_npoints = 0
799 do i = npoints, 1, -1
800 if (x(i) > splb%x_limit(2)-tol_x) cycle
801 aa = spline_eval(splb, x(i))
802 y(i) = y(i)/aa
803 new_npoints = new_npoints + 1
804 end do
805
806 call spline_fit(new_npoints, x, y, spla, threshold)
807
808 safe_deallocate_a(x)
809 safe_deallocate_a(y)
810
811 pop_sub(spline_div)
812 end subroutine spline_div
813
814 !------------------------------------------------------------
815 !Force all the values of the spline to be positive
816 !------------------------------------------------------------
817 subroutine spline_force_pos(spl, threshold)
818 type(spline_t), intent(inout) :: spl
819 real(real64), intent(in) :: threshold
820
821 integer :: npoints, i
822 real(real64), allocatable :: x(:), y(:)
823
824 push_sub(spline_force_pos)
825
826 npoints = oct_spline_npoints(spl%spl, spl%acc)
827 assert(npoints > 0)
828
829 safe_allocate(x(1:npoints))
830 safe_allocate(y(1:npoints))
831
832 call oct_spline_x(spl%spl, x(1))
833 call oct_spline_y(spl%spl, y(1))
834 call oct_spline_end(spl%spl, spl%acc)
835
836 !$omp parallel do simd
837 do i = npoints, 1, -1
838 y(i) = max(y(i),m_zero)
839 end do
840
841 call spline_fit(npoints, x, y, spl, threshold)
842
843 safe_deallocate_a(x)
844 safe_deallocate_a(y)
845
846 pop_sub(spline_force_pos)
847 end subroutine spline_force_pos
849
850 !------------------------------------------------------------
851 subroutine spline_mult(spla, splb, threshold)
852 type(spline_t), intent(inout) :: spla
853 type(spline_t), intent(in) :: splb
854 real(real64), intent(in) :: threshold
855
856 integer :: npoints, i
857 real(real64), allocatable :: x(:), y(:)
858 real(real64) :: aa
859
860 push_sub(spline_mult)
861
862 npoints = oct_spline_npoints(spla%spl, spla%acc)
863 if(npoints <= 0) then
864 pop_sub(spline_mult)
865 return
866 end if
867
868 safe_allocate(x(1:npoints))
869 safe_allocate(y(1:npoints))
870
871 call oct_spline_x(spla%spl, x(1))
872 call oct_spline_y(spla%spl, y(1))
873 call oct_spline_end(spla%spl, spla%acc)
874
875 assert(splb%x_limit(2) >= splb%x_limit(1))
876
877 !$omp parallel do simd private(aa)
878 do i = npoints, 1, -1
879 if (x(i) > splb%x_limit(2) - tol_x) then
880 aa = m_zero
881 else
882 aa = spline_eval(splb, x(i))
883 end if
884 y(i) = y(i)*aa
885 end do
887 call spline_fit(npoints, x, y, spla, threshold)
888
889 safe_deallocate_a(x)
890 safe_deallocate_a(y)
891
892 pop_sub(spline_mult)
893 end subroutine spline_mult
894
895
896 !------------------------------------------------------------
897 subroutine spline_der(spl, dspl, threshold)
898 type(spline_t), intent(in) :: spl
899 type(spline_t), intent(inout) :: dspl
900 real(real64), intent(in) :: threshold
901
902 integer :: npoints, i
903 real(real64), allocatable :: x(:), y(:)
904
905 push_sub(spline_der)
906
907 ! Use the grid of dspl if it is present, otherwise use the same one of spl.
908 if (.not. c_associated(dspl%spl)) then ! use the grid of spl
909 npoints = oct_spline_npoints(spl%spl, spl%acc)
910 assert(npoints > 0)
911 safe_allocate(x(1:npoints))
912 safe_allocate(y(1:npoints))
913 call oct_spline_x(spl%spl, x(1))
914 else ! use the grid of dspl
915 npoints = oct_spline_npoints(dspl%spl, dspl%acc)
916 assert(npoints > 0)
917 safe_allocate(x(1:npoints))
918 safe_allocate(y(1:npoints))
919 call oct_spline_x(dspl%spl, x(1))
920 end if
921 do i = 1, npoints
922 y(i) = oct_spline_eval_der(x(i), spl%spl, spl%acc)
923 end do
924 call spline_fit(npoints, x, y, dspl, threshold)
925
926 safe_deallocate_a(x)
927 safe_deallocate_a(y)
929 pop_sub(spline_der)
930 end subroutine spline_der
931
932
933 !------------------------------------------------------------
935 subroutine spline_der2(spl, dspl, threshold, grid)
936 type(spline_t), intent(in) :: spl
937 type(spline_t), intent(inout) :: dspl
938 real(real64), intent(in) :: threshold
939 real(real64), optional, intent(in) :: grid(:)
940
941 integer :: npoints, i
942 real(real64), allocatable :: x(:), y(:)
943
944 push_sub(spline_der2)
945
946 ! Use the grid of dspl if it is present, otherwise use the same one of spl.
947 if (.not. c_associated(dspl%spl)) then ! use the grid of spl
948 npoints = oct_spline_npoints(spl%spl, spl%acc)
949 assert(npoints > 0)
950 safe_allocate(x(1:npoints))
951 safe_allocate(y(1:npoints))
952 if (present(grid)) then
953 x(:) = grid(:)
954 else
955 call oct_spline_x(spl%spl, x(1))
956 end if
957 else ! use the grid of dspl
958 npoints = oct_spline_npoints(dspl%spl, dspl%acc)
959 assert(npoints > 0)
960 safe_allocate(x(1:npoints))
961 safe_allocate(y(1:npoints))
962 call oct_spline_x(dspl%spl, x(1))
963 end if
964 do i = 1, npoints
965 y(i) = oct_spline_eval_der2(x(i), spl%spl, spl%acc)
966 end do
967 call spline_fit(npoints, x, y, dspl, threshold)
968
969 safe_deallocate_a(x)
970 safe_deallocate_a(y)
971
972 pop_sub(spline_der2)
973 end subroutine spline_der2
974
975
976 !------------------------------------------------------------
977 subroutine spline_print_0(spl, iunit)
978 type(spline_t), intent(in) :: spl
979 integer, intent(in) :: iunit
980
981 integer :: np, i
982 real(real64), allocatable :: x(:), y(:)
983
984 push_sub(spline_print_0)
985
986 np = oct_spline_npoints(spl%spl, spl%acc)
987 assert(np > 0)
988 safe_allocate(x(1:np))
989 safe_allocate(y(1:np))
990
991 call oct_spline_x(spl%spl, x(1))
992 call oct_spline_y(spl%spl, y(1))
993 do i = 1, np
994 write(iunit, '(2es18.8E3)') x(i), y(i)
995 end do
996
997 safe_deallocate_a(x)
998 safe_deallocate_a(y)
999
1000 pop_sub(spline_print_0)
1001 end subroutine spline_print_0
1002
1003
1004 !------------------------------------------------------------
1006 subroutine spline_print_1(spl, iunit)
1007 type(spline_t), intent(in) :: spl(:)
1008 integer, intent(in) :: iunit
1009
1010 character(len=4) :: fm
1011 integer :: np, i, n, j
1012 real(real64), allocatable :: x(:)
1013 real(real64) :: x_limit
1014
1015 push_sub(spline_print_1)
1016
1017 n = size(spl)
1018 if (n <= 0) then
1019 pop_sub(spline_print_1)
1020 return
1021 end if
1022
1023 write(fm,'(i4)') n + 1
1024 fm = adjustl(fm)
1025 np = oct_spline_npoints(spl(1)%spl, spl(1)%acc)
1026 x_limit = spl(1)%x_limit(2)
1027 safe_allocate(x(1:np))
1028
1029 call oct_spline_x(spl(1)%spl, x(1))
1030
1031 ! Not all splines have the same number of points
1032 do i = 2, n
1033 np = min(np, oct_spline_npoints(spl(i)%spl, spl(i)%acc))
1034 x_limit = min(x_limit, spl(i)%x_limit(2))
1035 end do
1036
1037 do i = 1, np
1038 if (x(i) > x_limit - tol_x) cycle
1039 write(iunit, '('//trim(fm)//'es18.8E3)') x(i), (spline_eval(spl(j), x(i)), j = 1, size(spl))
1040 end do
1041
1042 safe_deallocate_a(x)
1043
1044 pop_sub(spline_print_1)
1045 end subroutine spline_print_1
1046
1047
1048 !------------------------------------------------------------
1050 subroutine spline_print_2(spl, iunit)
1051 type(spline_t), intent(in) :: spl(:, :)
1052 integer, intent(in) :: iunit
1053
1054 character(len=4) :: fm
1055 integer :: np, i, n1, n2, j, k
1056 real(real64), allocatable :: x(:)
1057 real(real64) :: x_limit
1058
1059 push_sub(spline_print_2)
1060
1061 n1 = size(spl, 1)
1062 n2 = size(spl, 2)
1063 if (n1 * n2 <= 0) then
1064 pop_sub(spline_print_2)
1065 return
1066 end if
1067
1068 write(fm,'(i4)') n1*n2 + 1
1069 fm = adjustl(fm)
1070
1071 np = oct_spline_npoints(spl(1, 1)%spl, spl(1, 1)%acc)
1072 safe_allocate(x(1:np))
1073
1074 call oct_spline_x(spl(1, 1)%spl, x(1))
1075 x_limit = spl(1,1)%x_limit(2)
1076
1077 ! Not all splines have the same number of points
1078 do i = 2, n1
1079 do j = 1, n2
1080 np = min(np, oct_spline_npoints(spl(i, j)%spl, spl(i, j)%acc))
1081 x_limit = min(x_limit, spl(i,j)%x_limit(2))
1082 end do
1083 end do
1084
1085 do i = 1, np
1086 if (x(i) > x_limit - tol_x) cycle
1087 write(iunit, '('//trim(fm)//'es18.8E3)') x(i), &
1088 ((spline_eval(spl(j, k), x(i)), j = 1, n1), k = 1, n2)
1089 end do
1091 safe_deallocate_a(x)
1092
1093 pop_sub(spline_print_2)
1094 end subroutine spline_print_2
1095
1096
1097 !------------------------------------------------------------
1100 real(real64) function spline_x_threshold(spl, threshold) result(r)
1101 type(spline_t), intent(in) :: spl
1102 real(real64), intent(in) :: threshold
1103
1104 integer :: ii, jj
1105 real(real64), parameter :: dx = 0.01_real64
1106
1107 ! No PUSH SUB, called too often.
1108 assert(spl%x_limit(2) >= spl%x_limit(1))
1109
1110 call profiling_in('SPLINE_CUTOFF_RADIUS')
1111
1112 jj = floor(spl%x_limit(2)/dx) + 1
1113 ! Numerical errors can lead to jj to be overestimated by 1
1114 ! This test prevents this to occur
1115 do while(dx*(jj-1) > spl%x_limit(2) - m_half * dx)
1116 jj = jj-1
1117 end do
1118
1119 do ii = jj, 1, -1
1120 r = dx*(ii-1)
1121 if (abs(spline_eval(spl, r)) > threshold) exit
1122 end do
1123
1124 call profiling_out('SPLINE_CUTOFF_RADIUS')
1125
1126 end function spline_x_threshold
1127
1128 ! -------------------------------------------------------
1129
1130 real(real64) pure function spline_range_min(this) result(range)
1131 type(spline_t), intent(in) :: this
1132
1133 range = this%x_limit(1)
1134
1135 end function spline_range_min
1136
1137 ! -------------------------------------------------------
1138
1139 real(real64) pure function spline_range_max(this) result(range)
1140 type(spline_t), intent(in) :: this
1141
1142 range = this%x_limit(2)
1143
1144 end function spline_range_max
1145
1146 ! -------------------------------------------------------
1147 subroutine spline_generate_shifted_grid(this, x_shifted)
1148 type(spline_t), intent(in) :: this
1149 real(real64), allocatable, intent(inout) :: x_shifted(:)
1150
1151 real(real64), allocatable :: x(:)
1152 integer :: npoints, i
1153
1155
1156 ! We shift the grid by half of a spacing to force interpolation of the function
1157 ! This is useful for getting back the derivative on the same grid
1158 npoints = oct_spline_npoints(this%spl, this%acc)
1159 safe_allocate(x(1:npoints))
1160 safe_allocate(x_shifted(1:npoints))
1161 call oct_spline_x(this%spl, x(1))
1162 do i = 1, npoints -1
1163 x_shifted(i) = m_half*(x(i) + x(i+1))
1164 end do
1165 x_shifted(npoints) = x(npoints)
1166 safe_deallocate_a(x)
1167
1169 end subroutine spline_generate_shifted_grid
1170
1171end module splines_oct_m
1172
1173!! Local Variables:
1174!! mode: f90
1175!! coding: utf-8
1176!! End:
Both the filling of the function, and the retrieval of the values may be done using single- or double...
Definition: splines.F90:166
Some operations may be done for one spline-function, or for an array of them.
Definition: splines.F90:179
The integral may be done with or without integration limits, but we want the interface to be common.
Definition: splines.F90:173
double exp(double __x) __attribute__((__nothrow__
double sin(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
double floor(double __x) __attribute__((__nothrow__
subroutine, public spline_der(spl, dspl, threshold)
Definition: splines.F90:849
subroutine spline_end_0(spl)
Definition: splines.F90:364
subroutine, public spline_force_pos(spl, threshold)
Definition: splines.F90:769
real(real64) pure function spline_integral_limits(spl, a, b)
Definition: splines.F90:518
subroutine, public spline_fit(nrc, rofi, ffit, spl, threshold)
Definition: splines.F90:413
real(real64) function, public spline_x_threshold(spl, threshold)
Determines the largest value of x for which the spline values are above the threshold.
Definition: splines.F90:1052
subroutine, public spline_3dft(spl, splw, threshold, gmax)
Definition: splines.F90:527
subroutine, public spline_der2(spl, dspl, threshold, grid)
Returns a spline that contains the second derivative of the original spline.
Definition: splines.F90:887
subroutine spline_end_2(spl)
Definition: splines.F90:395
real(real64), parameter tol_x
Tolerance for checking if x is inside or outside the limits of the splines.
Definition: splines.F90:306
subroutine, public spline_besselft(spl, splw, l, threshold, gmax)
Definition: splines.F90:605
real(real64) pure function, public spline_range_min(this)
Definition: splines.F90:1082
subroutine spline_end_1(spl)
Definition: splines.F90:379
subroutine spline_print_1(spl, iunit)
Print debug information for a 1D array of splines.
Definition: splines.F90:958
real(real64) function spline_integral_full(spl)
Definition: splines.F90:506
real(real64) function, public spline_eval(spl, x)
Definition: splines.F90:441
subroutine spline_init_0(spl)
Definition: splines.F90:312
subroutine spline_print_0(spl, iunit)
Definition: splines.F90:929
subroutine spline_print_2(spl, iunit)
Print debug information for a 2D array of splines.
Definition: splines.F90:1002
subroutine spline_eval8_array(spl, nn, xf)
Definition: splines.F90:450
subroutine spline_init_1(spl)
Definition: splines.F90:330
subroutine, public spline_cut(spl, cutoff, beta, threshold)
Definition: splines.F90:676
subroutine, public spline_div(spla, splb, threshold)
Returns the values of spla divided by the values of splb.
Definition: splines.F90:725
real(real64) pure function, public spline_range_max(this)
Definition: splines.F90:1091
subroutine, public spline_times(a, spl, threshold)
Definition: splines.F90:473
subroutine spline_init_2(spl)
Definition: splines.F90:346
subroutine, public spline_mult(spla, splb, threshold)
Definition: splines.F90:803
subroutine spline_evalz_array(spl, nn, xf)
Definition: splines.F90:462
subroutine, public spline_generate_shifted_grid(this, x_shifted)
Definition: splines.F90:1099
the basic spline datatype
Definition: splines.F90:156