20 integer,
parameter ::
dp = kind(1.0d0)
33 integer,
parameter :: PREC = 8
34 integer,
intent(in) :: N
35 real(PREC),
intent(in) :: A(N,N)
36 real(PREC),
intent(out) :: AI(N,N)
37 real(PREC) :: WORK(N+N*N)
38 integer :: LDA, LWORK, M, INFO, IPIV(N)
49 call dgetrf(m, n, ai, lda, ipiv, info)
50 call dgetri(n, ai, n, ipiv, work, lwork, info)
84 subroutine prg_kernel_multirank_latte(KRes,KK0_bml,Res,FelTol,L,LMAX,NUMRANK,HO_bml,mu,beta,RXYZ,Box,Hubbard_U,Element_Pointer, &
85 Nr_atoms,HDIM,Max_Nr_Neigh,Coulomb_acc,nebcoul,totnebcoul,Hinxlist, &
86 S_bml,Z_bml,Nocc,Inv_bml,H1_bml,DO_bml,D1_bml,m_rec,threshold,Nr_elem,mdstep,update)
93 integer,
intent(in) :: nr_atoms, hdim, max_nr_neigh,m_rec,nr_elem,numrank,mdstep,update
94 real(
dp),
parameter :: one = 1.d0, two = 2.d0, zero = 0.d0
95 real(
dp),
intent(in) :: coulomb_acc, feltol, beta, mu, nocc
96 real(
dp),
intent(in) :: rxyz(3,nr_atoms), box(3,3)
97 real(
dp),
intent(in) :: res(nr_atoms)
98 integer,
intent(in) :: hinxlist(hdim),element_pointer(nr_atoms)
99 real(
dp),
intent(in) :: hubbard_u(nr_elem)
100 type(bml_matrix_t),
intent(inout) :: ho_bml, s_bml, z_bml, inv_bml(m_rec), kk0_bml
101 real(
dp) :: k0res(nr_atoms)
102 type(bml_matrix_t),
intent(inout) :: h1_bml, do_bml,d1_bml
103 real(
dp),
intent(in) :: threshold
104 integer,
intent(in) :: lmax
105 character(20) :: bml_type
106 integer,
intent(in) :: totnebcoul(nr_atoms), nebcoul(4,max_nr_neigh,nr_atoms)
107 real(
dp) :: coulomb_pot_real(nr_atoms), coulomb_pot(nr_atoms)
108 real(
dp) :: coulomb_pot_real_i, coulomb_pot_k(nr_atoms),dq_dv(nr_atoms)
109 real(
dp) :: coulomb_pot_real_dq_v(nr_atoms), coulomb_pot_dq_v(nr_atoms)
110 real(
dp) :: dq_v(nr_atoms)
111 real(
dp),
intent(inout) :: kres(nr_atoms)
112 integer :: i,j,k,it,n,mn
113 integer,
intent(out) :: l
114 real(
dp) :: fel, proj_tmp, start, finish
115 real(
dp) :: vi(nr_atoms,lmax), fi(nr_atoms,lmax), v(nr_atoms)
116 real(
dp) :: dr(nr_atoms), identres(nr_atoms)
117 real(
dp),
allocatable :: o(:,:),m(:,:),row1(:),row2(:),row_na(:)
118 type(bml_matrix_t) :: kk0t_bml, k0res_bml, res_bml, zt_bml, t_bml, x_bml, y_bml
121 bml_type = bml_get_type(ho_bml)
122 n = bml_get_n(ho_bml)
123 mn = bml_get_m(ho_bml)
124 allocate(row1(hdim));
allocate(row2(hdim));
allocate(row_na(nr_atoms))
125 call bml_zero_matrix(bml_type,bml_element_real,
dp,nr_atoms,nr_atoms,kk0t_bml)
126 call bml_zero_matrix(bml_type,bml_element_real,
dp,nr_atoms,nr_atoms,k0res_bml)
127 call bml_zero_matrix(bml_type,bml_element_real,
dp,nr_atoms,nr_atoms,res_bml)
128 call bml_zero_matrix(bml_type,bml_element_real,
dp,n,mn,zt_bml)
129 call bml_zero_matrix(bml_type,bml_element_real,
dp,n,mn,x_bml)
130 call bml_zero_matrix(bml_type,bml_element_real,
dp,n,mn,y_bml)
132 call bml_transpose(z_bml,zt_bml)
135 call bml_set_row(res_bml,1,res,1.0_dp*1e-10)
136 call bml_transpose(kk0_bml,kk0t_bml)
137 call bml_multiply(res_bml,kk0t_bml,k0res_bml,1.0_dp,0.0_dp,1.0_dp*1e-10)
138 call bml_get_row(k0res_bml,1,row_na)
142 write(*,*)
'resnorm =', norm2(res),
'kresnorm =', norm2(dr)
145 do while ((fel > feltol).AND.(i < (lmax)).AND.(i < numrank))
148 vi(:,i) = dr/norm2(dr)
150 vi(:,i) = vi(:,i) - dot_product(vi(:,i),vi(:,j))*vi(:,j)
152 vi(:,i) = vi(:,i)/norm2(vi(:,i))
161 dq_v,hubbard_u,element_pointer,nr_atoms,coulomb_acc,nebcoul,totnebcoul,hdim,max_nr_neigh,nr_elem)
162 coulomb_pot_real(j) = coulomb_pot_real_i
167 coulomb_pot_dq_v = coulomb_pot_real+coulomb_pot_k
171 call bml_deallocate(h1_bml)
172 call bml_zero_matrix(bml_type,bml_element_real,
dp,hdim,mn,h1_bml)
176 do k = hinxlist(it)+1,hinxlist(it+1)
177 row1(k) = hubbard_u(element_pointer(it))*dq_v(it) + coulomb_pot_dq_v(it)
181 do k = hinxlist(nr_atoms)+1,hdim
182 row1(k) = hubbard_u(element_pointer(nr_atoms))*dq_v(nr_atoms) + coulomb_pot_dq_v(nr_atoms)
186 call bml_set_diagonal(h1_bml,row1,threshold)
187 call bml_multiply(s_bml,h1_bml,x_bml,1.0_dp,0.0_dp,threshold)
188 call bml_multiply(h1_bml,s_bml,x_bml,0.5_dp,0.5_dp,threshold)
189 call bml_copy(x_bml,h1_bml)
192 call bml_multiply(zt_bml,h1_bml,x_bml,1.0_dp,0.0_dp,threshold)
193 call bml_multiply(x_bml,z_bml,h1_bml,1.0_dp,0.0_dp,threshold)
197 m_rec, mu, beta, real(nocc,
dp), threshold)
200 call bml_multiply(z_bml,d1_bml,y_bml,2.0_dp,0.0_dp,threshold)
201 call bml_multiply(y_bml,zt_bml,d1_bml,1.0_dp,0.0_dp,threshold)
204 call bml_multiply(d1_bml,s_bml,x_bml, 1.0_dp,0.0_dp,threshold)
205 call bml_get_diagonal(x_bml,row1)
207 do it = 1, nr_atoms-1
209 do k = hinxlist(it)+1,hinxlist(it+1)
210 dq_dv(it) = dq_dv(it) + row1(k)
215 do k = hinxlist(nr_atoms)+1,hdim
216 dq_dv(nr_atoms) = dq_dv(nr_atoms) + row1(k)
221 call bml_set_row(res_bml,1,dr,1.0_dp*1e-10)
222 call bml_multiply(res_bml,kk0t_bml,k0res_bml,1.0_dp,0.0_dp,1.0_dp*1e-10)
223 call bml_get_row(k0res_bml,1,row_na)
228 allocate(o(l,l), m(l,l))
231 o(k,j) = dot_product(fi(:,k),fi(:,j))
235 identres = 0.d0*k0res
239 proj_tmp = m(k,j)*dot_product(fi(:,j),k0res)
240 identres = identres + proj_tmp*fi(:,k)
241 kres = kres + proj_tmp*vi(:,k)
244 fel = norm2(identres-k0res)/norm2(identres)
245 write(*,*)
'# I, L, Fel = ',i,l,fel
250 if (modulo(mdstep,update) .eq. 0)
then
254 deallocate(row1);
deallocate(row2);
deallocate(row_na)
255 call bml_deallocate(kk0t_bml)
256 call bml_deallocate(k0res_bml)
257 call bml_deallocate(res_bml)
258 call bml_deallocate(zt_bml)
259 call bml_deallocate(x_bml)
260 call bml_deallocate(y_bml)
268 type(bml_matrix_t),
intent(inout) :: K0
269 real(dp),
intent(in) :: v(Nr_atoms),fv(Nr_atoms),threshold
270 integer,
intent(in) :: Nr_atoms
271 type(bml_matrix_t) :: v_bml,fv_bml,vt_bml,tmp1_bml,tmp2_bml
272 type(bml_matrix_t) :: ones_bml,onest_bml,K0_update,K0_T
273 real(dp) :: const,kthresh
274 real(dp) :: ones(Nr_atoms)
275 real(dp),
allocatable :: row(:),norm(:)
277 character(20) :: bml_type
279 bml_type = bml_get_type(k0)
280 kthresh = 1.0_dp*1e-5
281 allocate(row(nr_atoms))
282 allocate(norm(nr_atoms))
284 call bml_zero_matrix(bml_type,bml_element_real,dp,nr_atoms,nr_atoms,v_bml)
285 call bml_zero_matrix(bml_type,bml_element_real,dp,nr_atoms,nr_atoms,vt_bml)
286 call bml_zero_matrix(bml_type,bml_element_real,dp,nr_atoms,nr_atoms,fv_bml)
287 call bml_zero_matrix(bml_type,bml_element_real,dp,nr_atoms,nr_atoms,tmp1_bml)
288 call bml_zero_matrix(bml_type,bml_element_real,dp,nr_atoms,nr_atoms,tmp2_bml)
289 call bml_zero_matrix(bml_type,bml_element_real,dp,nr_atoms,nr_atoms,ones_bml)
290 call bml_zero_matrix(bml_type,bml_element_real,dp,nr_atoms,nr_atoms,onest_bml)
291 call bml_zero_matrix(bml_type,bml_element_real,dp,nr_atoms,nr_atoms,k0_update)
292 call bml_zero_matrix(bml_type,bml_element_real,dp,nr_atoms,nr_atoms,k0_t)
295 call bml_set_row(onest_bml,1,ones,0.0_dp)
296 call bml_multiply(k0,ones_bml,tmp1_bml,1.0_dp,0.0_dp,threshold)
297 const = 1.0_dp/(bml_trace_mult(onest_bml,tmp1_bml) + 1)
299 call bml_set_row(vt_bml,1,v,1.0_dp*1e-10)
300 call bml_set_row(fv_bml,1,fv,1.0_dp*1e-10)
301 call bml_transpose(vt_bml,v_bml)
303 call bml_multiply(k0,v_bml,tmp1_bml,1.0_dp,0.0_dp,threshold)
304 call bml_copy(vt_bml,tmp2_bml)
305 call bml_multiply(fv_bml,k0,tmp2_bml,1.0_dp,-1.0_dp,threshold)
306 call bml_multiply(tmp1_bml,tmp2_bml,k0_update,const,0.0_dp,kthresh)
307 call bml_add(k0,k0_update,1.0_dp,-1.0_dp,kthresh)
309 call bml_multiply(onest_bml,k0,tmp1_bml,1.0_dp,0.0_dp,0.0_dp)
310 call bml_get_row(tmp1_bml,1,norm)
311 call bml_transpose(k0,k0_t)
313 call bml_get_row(k0_t,1,row)
314 call bml_set_row(k0_t,1,row/norm(i),1.0_dp*1e-10)
316 call bml_transpose(k0_t,k0)
320 call bml_deallocate(v_bml)
321 call bml_deallocate(vt_bml)
322 call bml_deallocate(fv_bml)
323 call bml_deallocate(tmp1_bml)
324 call bml_deallocate(tmp2_bml)
325 call bml_deallocate(onest_bml)
326 call bml_deallocate(ones_bml)
327 call bml_deallocate(k0_update)
328 call bml_deallocate(k0_t)
334 subroutine prg_kernel_multirank(KRes,KK0_bml,Res,FelTol,L,LMAX,HO_bml,mu,beta,RX,RY,RZ,LBox,Hubbard_U,Element_Type, &
335 Nr_atoms,HDIM,Max_Nr_Neigh,Coulomb_acc,TIMERATIO,nnRx,nnRy,nnRz,nrnnlist,nnType,H_INDEX_START,H_INDEX_END, &
336 S_bml,Z_bml,Nocc,Znuc,Inv_bml,H1_bml,X_bml,Y_bml,DO_bml,D1_bml,m_rec,threshold)
343 integer,
intent(in) :: nr_atoms, hdim, nocc,max_nr_neigh,m_rec
344 real(
dp),
parameter :: one = 1.d0, two = 2.d0, zero = 0.d0
345 real(
dp),
intent(in) :: coulomb_acc, timeratio, feltol, beta, mu
346 real(
dp),
intent(in) :: rx(nr_atoms), ry(nr_atoms), rz(nr_atoms), lbox(3)
347 real(
dp),
intent(in) :: res(nr_atoms)
348 integer,
intent(in) :: h_index_start(nr_atoms), h_index_end(nr_atoms)
349 real(
dp),
intent(in) :: znuc(nr_atoms), hubbard_u(nr_atoms)
350 type(bml_matrix_t),
intent(in) :: ho_bml, s_bml, z_bml, inv_bml(m_rec), kk0_bml
351 real(
dp) :: k0res(nr_atoms)
352 type(bml_matrix_t),
intent(inout) :: h1_bml, x_bml,y_bml,do_bml,d1_bml
353 real(
dp),
intent(in) :: threshold
354 integer,
intent(in) :: lmax
355 character(10),
intent(in) :: element_type(nr_atoms)
356 integer,
intent(in) :: nrnnlist(nr_atoms), nntype(nr_atoms,max_nr_neigh)
357 real(
dp),
intent(in) :: nnrx(nr_atoms,max_nr_neigh)
358 real(
dp),
intent(in) :: nnry(nr_atoms,max_nr_neigh), nnrz(nr_atoms,max_nr_neigh)
359 real(
dp) :: coulomb_pot_real(nr_atoms), coulomb_pot(nr_atoms)
360 real(
dp) :: coulomb_pot_real_i, coulomb_pot_k(nr_atoms),dq_dv(nr_atoms)
361 real(
dp) :: coulomb_pot_real_dq_v(nr_atoms), coulomb_pot_dq_v(nr_atoms)
362 real(
dp) :: coulomb_force_real_i(3),coulomb_force_k(3,nr_atoms)
363 real(
dp) :: dq_v(nr_atoms)
364 real(
dp),
intent(inout) :: kres(nr_atoms)
366 integer,
intent(out) :: l
367 real(
dp) :: fel, proj_tmp
368 real(
dp) :: vi(nr_atoms,lmax), fi(nr_atoms,lmax), v(nr_atoms)
369 real(
dp) :: dr(nr_atoms), identres(nr_atoms)
370 real(
dp),
allocatable :: o(:,:),m(:,:),row1(:),row2(:),row_na(:)
371 type(bml_matrix_t) :: kk0t_bml, k0res_bml, res_bml
373 allocate(row1(hdim));
allocate(row2(hdim));
allocate(row_na(nr_atoms))
374 call bml_zero_matrix(
"ellpack",bml_element_real,
dp,nr_atoms,nr_atoms,kk0t_bml)
375 call bml_zero_matrix(
"ellpack",bml_element_real,
dp,nr_atoms,nr_atoms,k0res_bml)
376 call bml_zero_matrix(
"ellpack",bml_element_real,
dp,nr_atoms,nr_atoms,res_bml)
379 call bml_transpose(kk0_bml,kk0t_bml)
380 call bml_set_row(res_bml,1,res,one*1e-14)
381 call bml_multiply(res_bml,kk0t_bml,k0res_bml,1.0_dp,0.0_dp,one*1e-14)
382 call bml_get_row(k0res_bml,1,row_na)
388 do while ((fel > feltol).AND.(i < (lmax)))
390 vi(:,i) = dr/norm2(dr)
392 vi(:,i) = vi(:,i) - dot_product(vi(:,i),vi(:,j))*vi(:,j)
394 vi(:,i) = vi(:,i)/norm2(vi(:,i))
401 call ewald_real_space(coulomb_pot_real_i,coulomb_force_real_i,j,rx,ry,rz,lbox, &
402 dq_v,hubbard_u,element_type,nr_atoms,coulomb_acc,timeratio,nnrx,nnry,nnrz,nrnnlist,nntype,hdim,max_nr_neigh)
403 coulomb_pot_real(j) = coulomb_pot_real_i
412 call bml_set_row(h1_bml,j,row1,threshold)
417 do k = h_index_start(j),h_index_end(j)
418 row1(k) = hubbard_u(j)*dq_v(j) + coulomb_pot_dq_v(j)
423 call bml_set_diagonal(h1_bml,row1,threshold)
424 call bml_multiply(s_bml,h1_bml,d1_bml,1.0_dp,0.0_dp,threshold)
425 call bml_multiply(h1_bml,s_bml,d1_bml,0.5_dp,0.5_dp,threshold)
426 call bml_copy(d1_bml,h1_bml)
428 call bml_transpose(z_bml,x_bml)
429 call bml_multiply(x_bml,h1_bml,d1_bml,1.0_dp,0.0_dp,threshold)
430 call bml_multiply(d1_bml,z_bml,h1_bml,1.0_dp,0.0_dp,threshold)
433 m_rec, mu, beta, real(nocc,
dp), threshold)
435 call bml_transpose(z_bml,x_bml)
436 call bml_multiply(z_bml,d1_bml,y_bml,2.0_dp,0.0_dp,threshold)
437 call bml_multiply(y_bml,x_bml,d1_bml,1.0_dp,0.0_dp,threshold)
442 do k = h_index_start(it), h_index_end(it)
443 call bml_get_row(s_bml,k,row1)
444 call bml_get_row(d1_bml,k,row2)
445 dq_dv(it) = dq_dv(it) + dot_product(row1,row2)
452 call bml_set_row(res_bml,1,dr,threshold)
453 call bml_multiply(res_bml,kk0t_bml,k0res_bml,1.0_dp,0.0_dp,threshold)
454 call bml_get_row(k0res_bml,1,row_na)
459 allocate(o(l,l), m(l,l))
462 o(k,j) = dot_product(fi(:,k),fi(:,j))
466 identres = 0.d0*k0res
470 proj_tmp = m(k,j)*dot_product(fi(:,j),k0res)
471 identres = identres + proj_tmp*fi(:,k)
472 kres = kres + proj_tmp*vi(:,k)
475 fel = norm2(identres-k0res)/norm2(identres)
476 write(*,*)
'# I, L, Fel = ',i,l,fel
481 deallocate(row1);
deallocate(row2);
deallocate(row_na)
482 call bml_deallocate(kk0t_bml)
483 call bml_deallocate(k0res_bml)
484 call bml_deallocate(res_bml)
512 Element_Pointer, Nr_atoms,HDIM, Max_Nr_Neigh,Coulomb_acc, &
513 Hinxlist,S_bml,Z_bml,Inv_bml,D1_bml,H1_bml,HO_bml, &
514 Nocc,m_rec,threshold,beta,Nr_elem,nebcoul,totnebcoul)
519 integer,
parameter :: prec = 8,
dp = kind(1.0d0)
520 integer,
intent(in) :: nr_atoms, nr_elem,hdim, max_nr_neigh,m_rec
521 real(prec),
parameter :: one = 1.d0, two = 2.d0, zero = 0.d0
522 real(prec),
parameter :: kb = 8.61739d-5
523 real(prec),
intent(in) :: coulomb_acc, threshold,beta,nocc
524 real(prec) :: v(nr_atoms)
525 real(prec),
intent(in) :: rxyz(3,nr_atoms),box(3,3)
526 integer,
intent(in) :: hinxlist(hdim)
527 real(prec),
intent(in) :: hubbard_u(nr_elem)
528 type(bml_matrix_t),
intent(in) :: inv_bml(m_rec)
529 type(bml_matrix_t),
intent(in) :: s_bml,z_bml,ho_bml
530 type(bml_matrix_t),
intent(inout) :: h1_bml,do_bml,d1_bml
531 real(prec),
intent(inout) :: mu0
532 integer,
intent(in) :: element_pointer(nr_atoms),nebcoul(4,max_nr_neigh,nr_atoms),totnebcoul(nr_atoms)
533 real(prec) :: coulomb_pot_real(nr_atoms), coulomb_pot(nr_atoms)
534 real(prec) :: coulomb_pot_real_i, coulomb_pot_k(nr_atoms)
535 real(prec) :: coulomb_pot_real_dq_v(nr_atoms), coulomb_pot_dq_v(nr_atoms)
536 real(prec) :: dq_v(nr_atoms)
537 real(prec) :: dq_dv(nr_atoms), err,tol,start,ewaldk,ewaldr,ewaldkacc,ewaldracc,response,respacc
538 real(prec),
intent(out) :: kk(nr_atoms,nr_atoms)
539 integer :: i,j,k, iter, mm,it,n,mn
540 real(prec),
allocatable :: row1(:),row2(:),jj(:,:),e(:,:),ereal(:,:),ekspace(:,:)
541 type(bml_matrix_t) :: zt_bml, r_bml, jji_bml, jj_bml, tmp_bml, x_bml, y_bml
542 character(20) :: bml_type
544 bml_type = bml_get_type(ho_bml)
545 n = bml_get_n(ho_bml)
546 mn = bml_get_m(ho_bml)
547 call bml_zero_matrix(bml_type,bml_element_real,
dp,n,mn,zt_bml)
548 call bml_zero_matrix(bml_type,bml_element_real,
dp,n,mn,x_bml)
549 call bml_zero_matrix(bml_type,bml_element_real,
dp,n,mn,y_bml)
550 call bml_transpose(z_bml,zt_bml)
551 allocate(row1(hdim));
allocate(row2(hdim));
allocate(jj(nr_atoms,nr_atoms))
552 allocate(ereal(nr_atoms,nr_atoms));
allocate(ekspace(nr_atoms,nr_atoms))
553 allocate(e(nr_atoms,nr_atoms))
554 coulomb_pot_dq_v = zero
564 write(*,*)
'Beginning response calculations'
565 ewaldracc = 0.0; ewaldkacc = 0.0; respacc = 0.0;
573 dq_v,j,hubbard_u,element_pointer,nr_atoms,coulomb_acc,hdim,max_nr_neigh)
574 coulomb_pot_real(i) = coulomb_pot_real_i
577 call cpu_time(ewaldr)
578 ewaldracc = ewaldracc + ewaldr-start
580 coulomb_pot_dq_v = coulomb_pot_real+coulomb_pot_k
581 call cpu_time(ewaldk)
582 ewaldkacc = ewaldkacc + ewaldk-ewaldr
585 call bml_deallocate(h1_bml)
586 call bml_zero_matrix(bml_type,bml_element_real,
dp,n,mn,h1_bml)
590 do k = hinxlist(i)+1,hinxlist(i+1)
591 row1(k) = hubbard_u(element_pointer(i))*dq_v(i) + coulomb_pot_dq_v(i)
595 do k = hinxlist(nr_atoms)+1,hdim
596 row1(k) = hubbard_u(element_pointer(nr_atoms))*dq_v(nr_atoms) + coulomb_pot_dq_v(nr_atoms)
599 call bml_set_diagonal(h1_bml,row1,threshold)
600 call bml_multiply(s_bml,h1_bml,x_bml,1.0_dp,0.0_dp,threshold)
601 call bml_multiply(h1_bml,s_bml,x_bml,0.5_dp,0.5_dp,threshold)
602 call bml_copy(x_bml,h1_bml)
604 call bml_multiply(zt_bml,h1_bml,x_bml,1.0_dp,0.0_dp,threshold)
605 call bml_multiply(x_bml,z_bml,h1_bml,1.0_dp,0.0_dp,threshold)
608 m_rec,mu0,beta,nocc,threshold)
609 call cpu_time(response)
610 respacc = respacc + response-ewaldk
612 call bml_multiply(z_bml,d1_bml,y_bml,2.0_dp,0.0_dp,threshold)
613 call bml_multiply(y_bml,zt_bml,d1_bml,1.0_dp,0.0_dp,threshold)
615 call bml_multiply(d1_bml,s_bml,x_bml,1.0_dp,0.0_dp,threshold)
616 call bml_get_diagonal(x_bml,row1)
618 do it = 1, nr_atoms-1
620 do k = hinxlist(it)+1,hinxlist(it+1)
621 dq_dv(it) = dq_dv(it) + row1(k)
627 do k = hinxlist(nr_atoms)+1,hdim
628 dq_dv(nr_atoms) = dq_dv(nr_atoms) + row1(k)
630 jj(nr_atoms,j) = dq_dv(nr_atoms)
635 jj(i,i) = jj(i,i) - one
638 write(*,*)
'ewaldracc =', ewaldracc
639 write(*,*)
'ewaldkacc =', ewaldkacc
640 write(*,*)
'implcit response time =', respacc
641 call invert(jj,kk,nr_atoms)
643 deallocate(row1);
deallocate(row2);
deallocate(jj)
644 deallocate(e);
deallocate(ereal);
deallocate(ekspace)
645 call bml_deallocate(zt_bml)
646 call bml_deallocate(x_bml)
647 call bml_deallocate(y_bml)
681 Nr_atoms,HDIM, Max_Nr_Neigh,Coulomb_acc,TIMERATIO,nnRx,nnRy,nnRz, &
682 nrnnlist,nnType,H_INDEX_START,H_INDEX_END,S_bml,Z_bml,Inv_bml,D1_bml,H1_bml,HO_bml,Y_bml,X_bml, &
683 Nocc,Znuc,m_rec,threshold,beta,diagonal)
688 integer,
parameter :: prec = 8,
dp = kind(1.0d0)
689 integer,
intent(in) :: nr_atoms, hdim, nocc, max_nr_neigh,m_rec
690 real(prec),
parameter :: one = 1.d0, two = 2.d0, zero = 0.d0
691 real(prec),
parameter :: kb = 8.61739d-5
692 real(prec),
intent(in) :: coulomb_acc, timeratio,threshold,beta
693 real(prec) :: v(nr_atoms)
694 real(prec),
intent(in) :: rx(nr_atoms), ry(nr_atoms), rz(nr_atoms), lbox(3)
695 integer,
intent(in) :: h_index_start(nr_atoms), h_index_end(nr_atoms)
696 real(prec),
intent(in) :: znuc(nr_atoms), hubbard_u(nr_atoms)
697 type(bml_matrix_t),
intent(in) :: inv_bml(m_rec)
698 type(bml_matrix_t),
intent(in) :: s_bml,z_bml,ho_bml
699 type(bml_matrix_t),
intent(inout) :: h1_bml,do_bml,d1_bml,y_bml,x_bml
700 real(prec),
intent(inout) :: mu0, diagonal(hdim)
701 character(10),
intent(in) :: element_type(nr_atoms)
702 integer,
intent(in) :: nrnnlist(nr_atoms), nntype(nr_atoms,max_nr_neigh)
703 real(prec),
intent(in) :: nnrx(nr_atoms,max_nr_neigh)
704 real(prec),
intent(in) :: nnry(nr_atoms,max_nr_neigh), nnrz(nr_atoms,max_nr_neigh)
705 real(prec) :: coulomb_pot_real(nr_atoms), coulomb_pot(nr_atoms)
706 real(prec) :: coulomb_pot_real_i, coulomb_pot_k(nr_atoms)
707 real(prec) :: coulomb_pot_real_dq_v(nr_atoms), coulomb_pot_dq_v(nr_atoms)
708 real(prec) :: coulomb_force_real_i(3),coulomb_force_k(3,nr_atoms)
709 real(prec) :: dq_v(nr_atoms)
710 real(prec) :: dq_dv(nr_atoms)
711 real(prec),
intent(out) :: kk(nr_atoms,nr_atoms)
712 integer :: i,j,k, iter, mm,it
713 real(prec),
allocatable :: row1(:),row2(:),jj(:,:)
715 allocate(row1(hdim));
allocate(row2(hdim));
allocate(jj(nr_atoms,nr_atoms))
716 coulomb_pot_dq_v = zero
727 dq_v,j,hubbard_u,element_type,nr_atoms,coulomb_acc,timeratio,hdim,max_nr_neigh)
728 coulomb_pot_real(i) = coulomb_pot_real_i
733 coulomb_pot_dq_v = coulomb_pot_real+coulomb_pot_k
737 call bml_set_row(h1_bml,i,diagonal,threshold)
742 do k = h_index_start(i),h_index_end(i)
743 diagonal(k) = hubbard_u(i)*dq_v(i) + coulomb_pot_dq_v(i)
748 call bml_set_diagonal(h1_bml,diagonal,threshold)
749 call bml_multiply(s_bml,h1_bml,d1_bml,1.0_dp,0.0_dp,threshold)
750 call bml_multiply(h1_bml,s_bml,d1_bml,0.5_dp,0.5_dp,threshold)
751 call bml_copy(d1_bml,h1_bml)
753 call bml_transpose(z_bml,x_bml)
754 call bml_multiply(x_bml,h1_bml,d1_bml,1.0_dp,0.0_dp,threshold)
755 call bml_multiply(d1_bml,z_bml,h1_bml,1.0_dp,0.0_dp,threshold)
758 m_rec,mu0,beta,real(nocc,prec),threshold)
760 call bml_transpose(z_bml,x_bml)
761 call bml_multiply(z_bml,d1_bml,y_bml,2.0_dp,0.0_dp,threshold)
762 call bml_multiply(y_bml,x_bml,d1_bml,1.0_dp,0.0_dp,threshold)
766 do k = h_index_start(it), h_index_end(it)
767 call bml_get_row(s_bml,k,row1)
768 call bml_get_row(d1_bml,k,row2)
769 dq_dv(it) = dq_dv(it) + dot_product(row1,row2)
777 jj(i,i) = jj(i,i) - one
779 call invert(jj,kk,nr_atoms)
780 deallocate(row1);
deallocate(row2);
deallocate(jj)
786 subroutine prg_kernel_matrix_multirank(KRes,KK0_bml,Res,FelTol,L,LMAX,HO_bml,mu,beta,RX,RY,RZ,LBox,Hubbard_U,Element_Type, &
787 Nr_atoms,HDIM,Max_Nr_Neigh,Coulomb_acc,TIMERATIO,nnRx,nnRy,nnRz,nrnnlist,nnType,H_INDEX_START,H_INDEX_END, &
788 S_bml,Z_bml,Nocc,Znuc,Inv_bml,H1_bml,X_bml,Y_bml,DO_bml,D1_bml,m_rec,threshold)
795 integer,
intent(in) :: nr_atoms, hdim, max_nr_neigh,m_rec
796 real(
dp),
parameter :: one = 1.d0, two = 2.d0, zero = 0.d0
797 real(
dp),
intent(in) :: coulomb_acc, timeratio, feltol, beta, mu,nocc
798 real(
dp),
intent(in) :: rx(nr_atoms), ry(nr_atoms), rz(nr_atoms), lbox(3)
799 real(
dp),
intent(in) :: res(nr_atoms)
800 integer,
intent(in) :: h_index_start(nr_atoms), h_index_end(nr_atoms)
801 real(
dp),
intent(in) :: znuc(nr_atoms), hubbard_u(nr_atoms)
802 type(bml_matrix_t),
intent(in) :: ho_bml, s_bml, z_bml, inv_bml(m_rec)
803 real(
dp) :: k0res(nr_atoms)
804 type(bml_matrix_t),
intent(inout) :: h1_bml, x_bml,y_bml,do_bml,d1_bml
805 type(bml_matrix_t),
intent(inout) :: kk0_bml
806 real(
dp),
intent(in) :: threshold
807 integer,
intent(in) :: lmax
808 character(10),
intent(in) :: element_type(nr_atoms)
809 integer,
intent(in) :: nrnnlist(nr_atoms), nntype(nr_atoms,max_nr_neigh)
810 real(
dp),
intent(in) :: nnrx(nr_atoms,max_nr_neigh)
811 real(
dp),
intent(in) :: nnry(nr_atoms,max_nr_neigh), nnrz(nr_atoms,max_nr_neigh)
812 real(
dp) :: coulomb_pot_real(nr_atoms), coulomb_pot(nr_atoms)
813 real(
dp) :: coulomb_pot_real_i, coulomb_pot_k(nr_atoms),dq_dv(nr_atoms)
814 real(
dp) :: coulomb_pot_real_dq_v(nr_atoms), coulomb_pot_dq_v(nr_atoms)
815 real(
dp) :: coulomb_force_real_i(3),coulomb_force_k(3,nr_atoms)
816 real(
dp) :: dq_v(nr_atoms)
817 real(
dp),
intent(out) :: kres(nr_atoms)
818 integer :: i,j,k,it,col,row
819 integer,
intent(out) :: l
820 real(
dp) :: fel, proj_tmp,elem
821 real(
dp) :: vi(nr_atoms,lmax), fi(nr_atoms,lmax), v(nr_atoms)
822 real(
dp) :: dr(nr_atoms), identres(nr_atoms)
823 real(
dp),
allocatable :: o(:,:),m(:,:),row1(:),row2(:)
825 allocate(row1(hdim));
allocate(row2(hdim));
830 do while ((fel > feltol).AND.(i < (lmax)))
832 vi(:,i) = dr/norm2(dr)
834 vi(:,i) = vi(:,i) - dot_product(vi(:,i),vi(:,j))*vi(:,j)
836 vi(:,i) = vi(:,i)/norm2(vi(:,i))
843 call ewald_real_space(coulomb_pot_real_i,coulomb_force_real_i,j,rx,ry,rz,lbox, &
844 dq_v,hubbard_u,element_type,nr_atoms,coulomb_acc,timeratio,nnrx,nnry,nnrz,nrnnlist,nntype,hdim,max_nr_neigh)
845 coulomb_pot_real(j) = coulomb_pot_real_i
854 call bml_set_row(h1_bml,j,row1,threshold)
859 do k = h_index_start(j),h_index_end(j)
860 row1(k) = hubbard_u(j)*dq_v(j) + coulomb_pot_dq_v(j)
865 call bml_set_diagonal(h1_bml,row1,threshold)
866 call bml_multiply(s_bml,h1_bml,d1_bml,1.0_dp,0.0_dp,threshold)
867 call bml_multiply(h1_bml,s_bml,d1_bml,0.5_dp,0.5_dp,threshold)
868 call bml_copy(d1_bml,h1_bml)
870 call bml_transpose(z_bml,x_bml)
871 call bml_multiply(x_bml,h1_bml,d1_bml,1.0_dp,0.0_dp,threshold)
872 call bml_multiply(d1_bml,z_bml,h1_bml,1.0_dp,0.0_dp,threshold)
875 m_rec, mu, beta, nocc, threshold)
877 call bml_transpose(z_bml,x_bml)
878 call bml_multiply(z_bml,d1_bml,y_bml,2.0_dp,0.0_dp,threshold)
879 call bml_multiply(y_bml,x_bml,d1_bml,1.0_dp,0.0_dp,threshold)
884 do k = h_index_start(it), h_index_end(it)
885 call bml_get_row(s_bml,k,row1)
886 call bml_get_row(d1_bml,k,row2)
887 dq_dv(it) = dq_dv(it) + dot_product(row1,row2)
896 allocate(o(l,l), m(l,l))
899 o(k,j) = dot_product(fi(:,k),fi(:,j))
907 proj_tmp = m(k,j)*dot_product(fi(:,j),res)
908 identres = identres + proj_tmp*fi(:,k)
909 kres = kres + proj_tmp*vi(:,k)
912 fel = norm2(identres-res)/norm2(identres)
913 write(*,*)
'# I, L, Fel = ',i,l,fel
916 if ((fel > feltol).AND.(i < (lmax)))
then
924 elem = elem + m(j,k)*vi(row,j)*fi(col,k)
927 if (abs(elem) > threshold)
then
928 call bml_set_element(kk0_bml,row,col,elem)
935 deallocate(row1);
deallocate(row2);
Module to obtain the density matrix by diagonalizing an orthogonalized Hamiltonian.
subroutine, public ewald_k_space_latte(COULOMBV, RXYZ, Box, DELTAQ, Nr_atoms, COULACC, Max_Nr_Neigh)
subroutine, public ewald_k_space_latte_single(COULOMBV, J, RXYZ, Box, DELTAQ, Nr_atoms, COULACC)
subroutine, public ewald_real_space_single_latte(COULOMBV, I, RXYZ, Box, Nr_elem, DELTAQ, J, U, Element_Pointer, Nr_atoms, COULACC, HDIM, Max_Nr_Neigh)
Find Coulomb potential on site I from single charge at site J.
subroutine, public ewald_real_space(COULOMBV, FCOUL, I, RX, RY, RZ, LBox, DELTAQ, U, Element_Type, Nr_atoms, COULACC, TIMERATIO, nnRx, nnRy, nnRz, nrnnlist, nnType, HDIM, Max_Nr_Neigh)
subroutine, public ewald_real_space_latte(COULOMBV, I, RXYZ, Box, DELTAQ, U, Element_Pointer, Nr_atoms, COULACC, nebcoul, totnebcoul, HDIM, Max_Nr_Neigh, Nr_Elem)
subroutine, public ewald_real_space_single(COULOMBV, FCOUL, I, RX, RY, RZ, LBox, DELTAQ, J, U, Element_Type, Nr_atoms, COULACC, TIMERATIO, HDIM, Max_Nr_Neigh)
subroutine, public prg_implicit_fermi_first_order_response(H0_bml, H1_bml, P0_bml, P1_bml, Inv_bml, nsteps, mu0, beta, nocc, threshold)
Calculate first order density matrix response to perturbations using Implicit Fermi Dirac.
The prg_normalize module.
subroutine, public prg_timer_shutdown()
Done with timers.
subroutine, public timer_prg_init()
Initialize timers.
subroutine, public prg_timer_start(itimer, tag)
Start Timing.
subroutine, public prg_timer_stop(itimer, verbose)
Stop timing.
Pre-conditioned O(N) calculation of the kernel for XL-BOMD.
subroutine, public prg_kernel_multirank(KRes, KK0_bml, Res, FelTol, L, LMAX, HO_bml, mu, beta, RX, RY, RZ, LBox, Hubbard_U, Element_Type, Nr_atoms, HDIM, Max_Nr_Neigh, Coulomb_acc, TIMERATIO, nnRx, nnRy, nnRz, nrnnlist, nnType, H_INDEX_START, H_INDEX_END, S_bml, Z_bml, Nocc, Znuc, Inv_bml, H1_bml, X_bml, Y_bml, DO_bml, D1_bml, m_rec, threshold)
subroutine, public prg_kernel_matrix_multirank(KRes, KK0_bml, Res, FelTol, L, LMAX, HO_bml, mu, beta, RX, RY, RZ, LBox, Hubbard_U, Element_Type, Nr_atoms, HDIM, Max_Nr_Neigh, Coulomb_acc, TIMERATIO, nnRx, nnRy, nnRz, nrnnlist, nnType, H_INDEX_START, H_INDEX_END, S_bml, Z_bml, Nocc, Znuc, Inv_bml, H1_bml, X_bml, Y_bml, DO_bml, D1_bml, m_rec, threshold)
subroutine, public prg_full_kernel(KK, DO_bml, mu0, RX, RY, RZ, LBox, Hubbard_U, Element_Type, Nr_atoms, HDIM, Max_Nr_Neigh, Coulomb_acc, TIMERATIO, nnRx, nnRy, nnRz, nrnnlist, nnType, H_INDEX_START, H_INDEX_END, S_bml, Z_bml, Inv_bml, D1_bml, H1_bml, HO_bml, Y_bml, X_bml, Nocc, Znuc, m_rec, threshold, beta, diagonal)
Compute full inverse Jacobian of q[n]-n (for development code)
subroutine, public prg_kernel_multirank_latte(KRes, KK0_bml, Res, FelTol, L, LMAX, NUMRANK, HO_bml, mu, beta, RXYZ, Box, Hubbard_U, Element_Pointer, Nr_atoms, HDIM, Max_Nr_Neigh, Coulomb_acc, nebcoul, totnebcoul, Hinxlist, S_bml, Z_bml, Nocc, Inv_bml, H1_bml, DO_bml, D1_bml, m_rec, threshold, Nr_elem, mdstep, update)
Compute low rank approximation of (K0*J)^(-1)*K0*(q[n]-n)(for LATTE)
subroutine invert(A, AI, N)
subroutine, public prg_full_kernel_latte(KK, DO_bml, mu0, RXYZ, Box, Hubbard_U, Element_Pointer, Nr_atoms, HDIM, Max_Nr_Neigh, Coulomb_acc, Hinxlist, S_bml, Z_bml, Inv_bml, D1_bml, H1_bml, HO_bml, Nocc, m_rec, threshold, beta, Nr_elem, nebcoul, totnebcoul)
Compute full inverse Jacobian of q[n]-n (for LATTE)
subroutine prg_update_preconditioner(K0, v, fv, Nr_atoms, threshold)