Octopus
spline_filter.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
22 use debug_oct_m
23 use global_oct_m
24 use io_oct_m
25 use, intrinsic :: iso_fortran_env
29
30 implicit none
31
32 private
33
34 public :: &
38
39 integer, parameter :: mask_n = 201
40 real(real64) :: mask_x(mask_n), mask_y(mask_n)
41
42contains
43
44 !----------------------------------------------------------------------------
45 subroutine spline_filter_bessel(spl, l, qmax, alpha, beta_fs, rcut, beta_rs, threshold)
46 type(spline_t), intent(inout) :: spl
47 integer, intent(in) :: l
48 real(real64), intent(in) :: qmax, alpha, beta_fs, rcut, beta_rs
49 real(real64), intent(in) :: threshold
50
51 type(spline_t) :: splw
52
53 push_sub(spline_filter_bessel)
54
55 call spline_init(splw)
56 call spline_besselft(spl, splw, l, gmax=m_four*qmax)
57 call spline_cut(splw, alpha*qmax, beta_fs)
58 call spline_besselft(splw, spl, l)
59 call spline_end(splw)
60 call spline_cut(spl, rcut, beta_rs, threshold)
61
63 end subroutine spline_filter_bessel
64
65
66 !----------------------------------------------------------------------------
67 subroutine spline_filter_mask_init()
68
69 integer :: iunit, i
70
72
73 iunit = io_open(trim(conf%share)//"/filter_mask.data", action='read', status='old')
74
75 do i = 1, mask_n
76 read(iunit, *) mask_x(i), mask_y(i)
77 end do
78
79 call io_close(iunit)
80
82 end subroutine spline_filter_mask_init
83
84
85 !----------------------------------------------------------------------------
86 subroutine spline_filter_mask(spl, l, rmax, qmax, alpha, gamma, threshold)
87 type(spline_t), intent(inout) :: spl
88 integer, intent(in) :: l
89 real(real64), intent(in) :: rmax
90 real(real64), intent(in) :: qmax
91 real(real64), intent(in) :: alpha
92 real(real64), intent(in) :: gamma
93 real(real64), intent(in) :: threshold
94
95 type(spline_t) :: mask, splw
96 real(real64) :: local_mask_x(mask_n), rcut, beta
97
98 push_sub(spline_filter_mask)
99
100 rcut = gamma*rmax
101
102 ! If rut is zero, the code below fails
103 assert(rcut > m_epsilon)
104
105 ! we define the mask function as f(r/rcut)
106 local_mask_x = mask_x*rcut
107 call spline_fit(mask_n, local_mask_x, mask_y, mask)
108
109 ! apply the mask function
110 call spline_div(spl, mask)
111
112 ! transform to bessel space
113 call spline_init(splw)
114 call spline_besselft(spl, splw, l, gmax=m_four*qmax)
115
116 ! apply a cutoff
117 beta = log(1.e5_real64)/(alpha - m_one)**2
118 call spline_cut(splw, qmax/alpha, beta)
119
120 ! transform back to real space
121 call spline_besselft(splw, spl, l)
122 call spline_end(splw)
123
124 ! remove the mask function
125 call spline_mult(spl, mask, threshold)
126
127 call spline_end(mask)
128
129 pop_sub(spline_filter_mask)
130
131 end subroutine spline_filter_mask
132
133end module spline_filter_oct_m
135!! Local Variables:
136!! mode: f90
137!! coding: utf-8
138!! End:
Some operations may be done for one spline-function, or for an array of them.
Definition: splines.F90:179
double log(double __x) __attribute__((__nothrow__
real(real64), parameter, public m_four
Definition: global.F90:195
real(real64), parameter, public m_epsilon
Definition: global.F90:207
type(conf_t), public conf
Global instance of Octopus configuration.
Definition: global.F90:181
real(real64), parameter, public m_one
Definition: global.F90:192
Definition: io.F90:116
subroutine, public io_close(iunit, grp)
Definition: io.F90:467
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:402
subroutine, public spline_filter_mask(spl, l, rmax, qmax, alpha, gamma, threshold)
subroutine, public spline_filter_mask_init()
subroutine, public spline_filter_bessel(spl, l, qmax, alpha, beta_fs, rcut, beta_rs, threshold)
subroutine, public spline_fit(nrc, rofi, ffit, spl, threshold)
Definition: splines.F90:413
subroutine, public spline_besselft(spl, splw, l, threshold, gmax)
Definition: splines.F90:605
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
subroutine, public spline_mult(spla, splb, threshold)
Definition: splines.F90:803