Octopus
scaling_function.F90
Go to the documentation of this file.
1!! Copyright (C) 2011 X. Andrade
2!! Copyright (C) Luigi Genovese, Thierry Deutsch, CEA Grenoble, 2006
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include <global.h>
21
23 use, intrinsic :: iso_fortran_env
24 use blas_oct_m
26
27 implicit none
28
29 private
30
31 public :: &
34
35contains
36
37 !!****h* BigDFT/scaling_function
38 !! NAME
39 !! scaling_function
40 !!
41 !! FUNCTION
42 !! Calculate the values of a scaling function in real uniform grid
43 !!
44 !! SOURCE
45 !!
46 subroutine scaling_function(itype,nd,nrange,a,x)
47 integer, intent(in) :: itype
48 integer, intent(in) :: nd
49 integer, intent(out) :: nrange
50 real(real64), dimension(0:nd), intent(out) :: a,x
51
52 real(real64), dimension(:), allocatable :: y
53 integer :: i,nt,ni
54
55 !Only itype=8,14,16,20,24,30,40,50,60,100
56 select case (itype)
57 case (8)
58 !O.K.
59 case default
60 message(1) = "Only interpolating functions 8, 14, 16, 20, 24, 30, 40, 50, 60, 100."
61 call messages_fatal(1)
62 end select
63!!$ write(unit=*,fmt="(1x,a,i0,a)") &
64!!$ "Use interpolating scaling functions of ",itype," order"
65
66 !Give the range of the scaling function
67 !from -itype to itype
68 ni=2*itype
69 nrange = ni
70 allocate(y(0:nd))
71
72 ! plot scaling function
73 x = 0.0_8
74 y = 0.0_8
75
76 nt=ni
77 x(nt/2-1)=1._real64
78 loop1: do
79 nt=2*nt
80 ! write(6,*) 'nd,nt',nd,nt
81 select case (itype)
82 case (8)
83 call back_trans_8(nd,nt,x,y)
84 end select
85 call blas_copy(nt, y(0), 1, x(0) ,1)
86 if (nt == nd) then
87 exit loop1
88 end if
89 end do loop1
90
91 !open (unit=1,file='scfunction',status='unknown')
92 do i=0,nd
93 a(i) = 1._real64 * i * ni / nd - (.5_real64 * ni - 1._real64)
94 !write(1,*) 1._real64*i*ni/nd-(.5_real64*ni-1._real64),x(i)
95 end do
96 !close(1)
97
98 deallocate(y)
99 end subroutine scaling_function
100 !!***
101
102 !!****h* BigDFT/scf_recursion
103 !! NAME
104 !! scf_recursion
105 !!
106 !! FUNCTION
107 !! Do iterations to go from p0gauss to pgauss
108 !! order interpolating scaling function
109 !!
110 !! SOURCE
111 !!
112 subroutine scf_recursion(itype,n_iter,n_range,kernel_scf,kern_1_scf)
113 integer, intent(in) :: itype,n_iter,n_range
114 real(kind=8), intent(inout) :: kernel_scf(-n_range:n_range)
115 real(kind=8), intent(out) :: kern_1_scf(-n_range:n_range)
116
117 !Only itype=8,14,16,20,24,30,40,50,60,100
118 select case (itype)
119 case (8)
120 !O.K.
121 case default
122 message(1) = "Only interpolating functions 8, 14, 16, 20, 24, 30, 40, 50, 60, 100."
123 call messages_fatal(1)
124 end select
125
126 select case (itype)
127 case (8)
128 call scf_recursion_8(n_iter,n_range,kernel_scf,kern_1_scf)
129 end select
130
131 end subroutine scf_recursion
132 !!***
133
134 !!****h* BigDFT/back_trans_8
135 !! NAME
136 !! back_trans_8
137 !!
138 !! FUNCTION
139 !!
140 !! SOURCE
141 !!
142 ! backward wavelet transform
143 ! nd: length of data set
144 ! nt length of data in data set to be transformed
145 ! m filter length (m has to be even!)
146 ! x input data, y output data
147 subroutine back_trans_8(nd,nt,x,y)
148 integer, intent(in) :: nd,nt
149 real(kind=8), intent(in) :: x(0:nd-1)
150 real(kind=8), intent(out) :: y(0:nd-1)
151
152 integer :: i,j,ind
153
154 ! TODO(Alex) Issue #872. Move lazy_8_inc into a function with scaling_function.F90
155#include "lazy_8_inc.F90"
156
157 do i=0,nt/2-1
158 y(2*i+0)=0._real64
159 y(2*i+1)=0._real64
160
161 do j=-m/2,m/2-1
162
163 ! periodically wrap index if necessary
164 ind=i-j
165 loop99: do
166 if (ind < 0) then
167 ind=ind+nt/2
168 cycle loop99
169 end if
170 if (ind >= nt/2) then
171 ind=ind-nt/2
172 cycle loop99
173 end if
174 exit loop99
175 end do loop99
176
177 y(2*i+0)=y(2*i+0) + ch(2*j-0)*x(ind)+cg(2*j-0)*x(ind+nt/2)
178 y(2*i+1)=y(2*i+1) + ch(2*j+1)*x(ind)+cg(2*j+1)*x(ind+nt/2)
179 end do
180
181 end do
182
183 end subroutine back_trans_8
184 !!***
185
186 !!****h* BigDFT/scf_recursion_8
187 !! NAME
188 !! scf_recursion_8
189 !!
190 !! FUNCTION
191 !! Do iterations to go from p0gauss to pgauss
192 !! 8th-order interpolating scaling function
193 !!
194 !! SOURCE
195 !!
196 subroutine scf_recursion_8(n_iter,n_range,kernel_scf,kern_1_scf)
197 integer, intent(in) :: n_iter,n_range
198 real(kind=8), intent(inout) :: kernel_scf(-n_range:n_range)
199 real(kind=8), intent(out) :: kern_1_scf(-n_range:n_range)
200
201 real(kind=8) :: kern,kern_tot
202 integer :: i_iter,i,j,ind
203
204 ! TODO(Alex) Issue #872. Move lazy_8_inc into a function with scaling_function.F90
205#include "lazy_8_inc.F90"
206
207 !Start the iteration to go from p0gauss to pgauss
208 loop_iter_scf: do i_iter=1,n_iter
209 kern_1_scf(:) = kernel_scf(:)
210 kernel_scf(:) = 0._real64
211 loop_iter_i: do i=0,n_range
212 kern_tot = 0._real64
213 do j=-m,m
214 ind = 2*i-j
215 if (abs(ind) > n_range) then
216 kern = 0._real64
217 else
218 kern = kern_1_scf(ind)
219 end if
220 kern_tot = kern_tot + ch(j)*kern
221 end do
222 if (abs(kern_tot) <= 1d-150) then
223 !zero after (be sure because strictly == 0._real64)
224 exit loop_iter_i
225 else
226 kernel_scf( i) = 0.5_real64*kern_tot
227 kernel_scf(-i) = kernel_scf(i)
228 end if
229 end do loop_iter_i
230 end do loop_iter_scf
231 end subroutine scf_recursion_8
232 !!***
233
234end module scaling_function_oct_m
235
236!! Local Variables:
237!! mode: f90
238!! coding: utf-8
239!! End:
--------------— copy ---------------— Copies a vector, x, to a vector, y.
Definition: blas.F90:205
This module contains interfaces for BLAS routines You should not use these routines directly....
Definition: blas.F90:118
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
subroutine scf_recursion_8(n_iter, n_range, kernel_scf, kern_1_scf)
subroutine back_trans_8(nd, nt, x, y)
subroutine, public scaling_function(itype, nd, nrange, a, x)
subroutine, public scf_recursion(itype, n_iter, n_range, kernel_scf, kern_1_scf)