-
Notifications
You must be signed in to change notification settings - Fork 0
/
module_operator_radial.f90
333 lines (279 loc) · 12.2 KB
/
module_operator_radial.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
!======================================================================================
! calculate the matrix elment of kinetic and nuclei attractive potential operator
! in radial part
!======================================================================================
module operator_radial
use global
use fedvr3d_basis_set
use not_full_operator_index
implicit none
! the matrix element of kinetic operator in radial part
real(kind=k1),allocatable,save :: tmat_radial(:)
! the matrix element of v potential in radial part
real(kind=k1),allocatable,save :: vmat_radial(:)
contains
!==============================================================================
! calc. the kinetic energy operator in fedvr_radial part space
!calc. the kinetic operator in fe-dvr. e f
! d x (x) d x (x)
! e 1 d f 1 i j
! < x (x) | - ~~~ ~~~ | x (x) > = ~~~~~ (delta(e,f) + delta(e,f-1) + delta(e,f+1))* Int { dx * ~~~~~~~ * ~~~~~~~}
! i 2 dx j 2 dx dx
!
!==============================================================================
subroutine calc_tmat_radial
implicit none
integer :: idicp,e,i,f,j,ee,ii
integer :: ibasis,jbasis
real(kind=k1) :: sumtemp
allocate(tmat_radial(n_total_kinetic))
tmat_radial = 0.0_k1
do idicp =1, n_total_kinetic
ibasis = index_kinetic_basis(idicp,1); jbasis = index_kinetic_basis(idicp,2)
! from global to local e,i,f,j
e = which_element(ibasis)
i = which_basis(ibasis)
f = which_element(jbasis)
j = which_basis(jbasis)
!===================================================================================================
if( e==f ) then
if((i/= fedvr3d%fedvr_nb(e)).and.(j/=fedvr3d%fedvr_nb(f))) then
sumtemp = 0.0_k1
do ii =1,fedvr3d%fedvr_nb(e)
sumtemp = sumtemp + fedvr_w (e,ii) * derivate_x(e,i,e,ii) * derivate_x(f,j,e,ii)
enddo
elseif(((i==fedvr3d%fedvr_nb(e)).and.(j/=fedvr3d%fedvr_nb(f))) &
.or.((i/=fedvr3d%fedvr_nb(e)).and.(j==fedvr3d%fedvr_nb(f))) ) then
sumtemp = 0.0_k1
do ii =1,fedvr3d%fedvr_nb(e)
sumtemp = sumtemp + fedvr_w (e,ii) * derivate_x(e,i,e,ii) * derivate_x(f,j,e,ii)
enddo
elseif((i== fedvr3d%fedvr_nb(e)).and.(j==fedvr3d%fedvr_nb(f))) then
sumtemp = 0.0_k1
do ee =e, e+1
do ii =1, fedvr3d%fedvr_nb(ee)
sumtemp = sumtemp + fedvr_w (ee,ii) * derivate_x(e,i,ee,ii) * derivate_x(f,j,ee,ii)
enddo
enddo
endif
!==================================================================================================
elseif(e==(f-1)) then
sumtemp = 0.0_k1
do ii =1,fedvr3d%fedvr_nb(f)
sumtemp = sumtemp + fedvr_w (f,ii) * derivate_x(e,i,f,ii) * derivate_x(f,j,f,ii)
enddo
!===================================================================================================
elseif(e==(f+1)) then
sumtemp = 0.0_k1
do ii =1,fedvr3d%fedvr_nb(e)
sumtemp = sumtemp + fedvr_w (e,ii) * derivate_x(e,i,e,ii) * derivate_x(f,j,e,ii)
enddo
!====================================================================================================
endif
tmat_radial(idicp) = sumtemp*0.5d0
enddo
return
end subroutine calc_tmat_radial
!==============================================================================
! calc. the kinetic energy operator in fedvr_radial part space
!==============================================================================
!====================================================================================
! the fedvr1d , fedvr3d is the same one
! e f
! < x | v | x > == v ( x(e,i) ) * delta(e,f)*delta(i,j)
! i j
!
! e 1 f e 1 f e 1 f
! for fedvr-1d,< x | ~~~~~ | x >, for,fedvr-3d <x |~~~~~~~~~~~| x > = <x |~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~|x >
! i r j i --> j i sqrt((rsin[theta]cos[phi])^2+(rsin[theta]sin[phi])^2) + (r*cos[theta])^2 j
! r
! e 1 f
! = <x |~~~~~~ | x >
! i r j
!====================================================================================
subroutine calc_vmat_radial
implicit none
integer :: idicp, e,i
allocate(vmat_radial(fedvr3d%nb_r))
do idicp =1,fedvr3d%nb_r
e = which_element(idicp)
i = which_basis(idicp)
vmat_radial(idicp) = - system%zz/( fedvr_x(e,i ) )
enddo
return
end subroutine calc_vmat_radial
!==============================================================================================
! assume : num_of_element, the number of the element
! n(e), the number of the basis function in e-th element
!
! e
! find derivate of x
! i
!
!
! for element 1,2,3....., num_of_element -1, there are i =2,3,4,...,ne
!
!
! e e e+1
! f (x) f (x) + f (x)
! e i e ne 1
! x (x) =~~~~~~~~~~~~~ , x = ~~~~~~~~~~~~~~~~~~~ , {i=2, n(e) -1 , ne }
! i e ne e e+1
! sqrt( w ) sqrt( w + w )
! i ne 1
!
! e e
! d x (x) d f (x)
! i i 1
! ~~~~~~~~~~~ == ~~~~~~~~~~~~ * ~~~~~~~~~~~ , for i = 2, ne -1
! d x dx e
! sqrt( w )
! i
!
! e e e+1
! d x (x) d f (x) d f (x)
! ne ne 1 1
! ~~~~~~~~~~~ == [ ~~~~~~~~~~~~~ + ~~~~~~~~~~ ] * ~~~~~~~~~~~~~ , for i = ne
! d x dx dx e e+1
! sqrt(w + w )
! ne 1
!
! For num_of_element-th, the last element have no brige function
!
!
!
! e
! f (x)
! e i
! x (x) =~~~~~~~~~~~~~ {i =2, n(num_of_element) - 1 }
! i e
! sqrt( w )
! i
!
!
! e e
! d x (x) d f (x)
! i i 1
! ~~~~~~~~~~~ == ~~~~~~~~~~~~ * ~~~~~~~~~~~ , for i = 2, ne -1
! d x dx e
! sqrt( w )
! i
! there are no ne !!!!! Attention here
!=================================================================================================
!
!
! e
! d x (x) |
! i | === derivate_x(e,i,f,j)
!~~~~~~~~~~~ | f
! dx | x = x
! j
!
!=================================================================================================
double precision function derivate_x(e,i,f,j)
implicit none
integer :: e,i,f,j,ne_here
if( e<1 .or. e>fedvr3d%number_of_element .or. i<2 .or. i > fedvr3d%fedvr_nb(e)) then
write(6,*) 'error happen in function deriate_x'
stop
endif
if((e==fedvr3d%number_of_element) .and. (i > fedvr3d%fedvr_nb(e)-1)) then
write(6,*) 'error happen in function deriate_x'
stop
endif
if( f<1 .or. f>fedvr3d%number_of_element .or. i<1 .or. i > fedvr3d%fedvr_nb(f)) then
write(6,*) 'error happen in function deriate_x'
stop
endif
!
! 1< e < num_of_element - 1
!
!
ne_here = fedvr3d%fedvr_nb(e)
if((e>=1) .and. (e<= fedvr3d%number_of_element -1)) then
if((i>=2) .and. (i<= fedvr3d%fedvr_nb(e) -1 )) then
derivate_x = derivate_lobatto(e,i,f,j)/(sqrt(fedvr_w(e,i)))
elseif(i ==fedvr3d%fedvr_nb(e) ) then
derivate_x = (derivate_lobatto(e,ne_here,f,j) + derivate_lobatto(e+1,1,f,j) ) &
/(sqrt( fedvr_w(e,ne_here) + fedvr_w(e+1,1) ))
else
write(6,*) 'error happen in function derivate_x'
endif
elseif(e == fedvr3d%number_of_element) then
if((i>=2) .and. (i<= fedvr3d%fedvr_nb(e) -1 )) then
derivate_x = derivate_lobatto(e,i,f,j)/(sqrt(fedvr_w(e,i)))
else
write(6,*) 'error happen in function derivate_x'
endif
else
write(6,*) 'error element in function derviate_x'
endif
return
end function derivate_x
! derivative of Lobatto shape function
! a speed up is possible by circumventing somehow the calculation of indices
! df_{e,i}(x_{f,j})
!=============================================================================================
!calc. the derivative of Lobatto shape function
!
! e
! d f (x) |
! i | = 0, e/=f
!~~~~~~~~~ | f
! d x | x = x
! j
! 1
! = ~~~~~~~~~~*( delta(i, n_e) - delta(i, 1) ) , i==j , e==f
! e
! 2* w
! i
! e e
! 1 ____ x - x
! = ~~~~~~~~~~~ * || j k
! e e || ~~~~~~~~~~~~ , i /= j , e==f
! x x k/=i,j e e
! i j x - x
! i k
!===============================================================================================
double precision function derivate_lobatto(e,i,f,j)
use global
implicit none
integer :: e,i,f,j,k
real(kind=k1) :: sum_temp,product_temp
if (e/=f) then
derivate_lobatto=0.0_k1
else
if (i==j) then
sum_temp = 0.0_k1
if(i == fedvr3d%fedvr_nb(e)) then
sum_temp = sum_temp + 1.0_k1
endif
if(i == 1) then
sum_temp = sum_temp - 1.0_k1
endif
derivate_lobatto=1.0d0/(2.0d0*fedvr_w(e,i))*(sum_temp)
else
product_temp=1.0_k1
do k=1,fedvr3d%fedvr_nb(e)
if ((k/=i).and.(k/=j)) then
product_temp=product_temp*(fedvr_x(e,j)-fedvr_x(e,k))/(fedvr_x(e,i)-fedvr_x(e,k))
endif
enddo
derivate_lobatto=product_temp/(fedvr_x(e,i)-fedvr_x(e,j))
endif
endif
end function derivate_lobatto
end module operator_radial
!
! drive subroutine, find the information about the radial part and angular part
!
subroutine drive_operator_radial
use operator_radial
implicit none
!
! radial part, kinetic energy operator and nuclei attractive potential
!
call calc_tmat_radial
call calc_vmat_radial
return
end subroutine drive_operator_radial