diff --git a/src/casscf_cipsi/gradient.irp.f b/src/casscf_cipsi/gradient.irp.f index 961d260d..cd608fb3 100644 --- a/src/casscf_cipsi/gradient.irp.f +++ b/src/casscf_cipsi/gradient.irp.f @@ -14,8 +14,8 @@ END_PROVIDER implicit none n_c_a_prov = n_core_inact_orb * n_act_orb n_c_v_prov = n_core_inact_orb * n_virt_orb - n_a_v_prov = n_act_orb * n_virt_orb - END_PROVIDER + n_a_v_prov = n_act_orb * n_virt_orb + END_PROVIDER BEGIN_PROVIDER [integer, excit, (2,nMonoEx)] &BEGIN_PROVIDER [character*3, excit_class, (nMonoEx)] @@ -28,7 +28,7 @@ END_PROVIDER BEGIN_DOC ! a list of the orbitals involved in the excitation END_DOC - + implicit none integer :: i,t,a,ii,tt,aa,indx,indx_tmp indx=0 @@ -48,7 +48,7 @@ END_PROVIDER mat_idx_c_a(ii,tt) = indx end do end do - + indx_tmp = 0 do ii=1,n_core_inact_orb i=list_core_inact(ii) @@ -61,11 +61,11 @@ END_PROVIDER indx_tmp += 1 list_idx_c_v(1,indx_tmp) = indx list_idx_c_v(2,indx_tmp) = ii - list_idx_c_v(3,indx_tmp) = aa + list_idx_c_v(3,indx_tmp) = aa mat_idx_c_v(ii,aa) = indx end do end do - + indx_tmp = 0 do tt=1,n_act_orb t=list_act(tt) @@ -82,7 +82,7 @@ END_PROVIDER mat_idx_a_v(tt,aa) = indx end do end do - + if (bavard) then write(6,*) ' Filled the table of the Monoexcitations ' do indx=1,nMonoEx @@ -90,7 +90,7 @@ END_PROVIDER ,excit(2,indx),' ',excit_class(indx) end do end if - + END_PROVIDER BEGIN_PROVIDER [real*8, gradvec2, (nMonoEx)] @@ -104,7 +104,7 @@ END_PROVIDER implicit none integer :: i,t,a,indx real*8 :: gradvec_it,gradvec_ia,gradvec_ta - + indx=0 norm_grad_vec2_tab = 0.d0 do i=1,n_core_inact_orb @@ -114,7 +114,7 @@ END_PROVIDER norm_grad_vec2_tab(1) += gradvec2(indx)*gradvec2(indx) end do end do - + do i=1,n_core_inact_orb do a=1,n_virt_orb indx+=1 @@ -122,7 +122,7 @@ END_PROVIDER norm_grad_vec2_tab(2) += gradvec2(indx)*gradvec2(indx) end do end do - + do t=1,n_act_orb do a=1,n_virt_orb indx+=1 @@ -130,7 +130,7 @@ END_PROVIDER norm_grad_vec2_tab(3) += gradvec2(indx)*gradvec2(indx) end do end do - + norm_grad_vec2=0.d0 do indx=1,nMonoEx norm_grad_vec2+=gradvec2(indx)*gradvec2(indx) @@ -144,7 +144,7 @@ END_PROVIDER write(6,*) ' Norm of the orbital gradient (via D, P and integrals): ', norm_grad_vec2 write(6,*) endif - + END_PROVIDER real*8 function gradvec_it(i,t) @@ -154,23 +154,30 @@ real*8 function gradvec_it(i,t) END_DOC implicit none integer :: i,t - + integer :: ii,tt,v,vv,x,y integer :: x3,y3 double precision :: bielec_PQxx_no - + ii=list_core_inact(i) tt=list_act(t) gradvec_it=2.D0*(Fipq(tt,ii)+Fapq(tt,ii)) gradvec_it-=occnum(tt)*Fipq(ii,tt) - do v=1,n_act_orb ! active - vv=list_act(v) - do x=1,n_act_orb ! active - x3=x+n_core_inact_orb ! list_act(x) - do y=1,n_act_orb ! active - y3=y+n_core_inact_orb ! list_act(y) + do y=1,n_act_orb ! active +! y3=y+n_core_inact_orb ! list_act(y) + do x=1,n_act_orb ! active +! x3=x+n_core_inact_orb ! list_act(x) + do v=1,n_act_orb ! active + vv=list_act(v) ! Gamma(2) a a a a 1/r12 i a a a - gradvec_it-=2.D0*P0tuvx_no(t,v,x,y)*bielec_PQxx_no(ii,vv,x3,y3) +! gradvec_it-=2.D0*P0tuvx_no(t,v,x,y)*bielec_PQxx_no(ii,vv,x3,y3) + integer :: ichol + double precision :: tmp + tmp = 0.d0 + do ichol=1,cholesky_mo_num + tmp = tmp + cholesky_no_total_transp(ichol,vv,ii) * cholesky_no_total_transp(ichol,list_act(x),list_act(y)) + enddo + gradvec_it = gradvec_it - 2.D0*P0tuvx_no(t,v,x,y)*tmp end do end do end do @@ -183,12 +190,12 @@ real*8 function gradvec_ia(i,a) END_DOC implicit none integer :: i,a,ii,aa - + ii=list_core_inact(i) aa=list_virt(a) gradvec_ia=2.D0*(Fipq(aa,ii)+Fapq(aa,ii)) gradvec_ia*=2.D0 - + end function gradvec_ia real*8 function gradvec_ta(t,a) @@ -198,7 +205,7 @@ real*8 function gradvec_ta(t,a) END_DOC implicit none integer :: t,a,tt,aa,v,vv,x,y - + tt=list_act(t) aa=list_virt(a) gradvec_ta=0.D0 @@ -211,6 +218,6 @@ real*8 function gradvec_ta(t,a) end do end do gradvec_ta*=2.D0 - + end function gradvec_ta diff --git a/src/casscf_cipsi/hessian.irp.f b/src/casscf_cipsi/hessian.irp.f index 1ee073d2..431a6979 100644 --- a/src/casscf_cipsi/hessian.irp.f +++ b/src/casscf_cipsi/hessian.irp.f @@ -11,13 +11,14 @@ real*8 function hessmat_itju(i,t,j,u) integer :: i,t,j,u,ii,tt,uu,v,vv,x,xx,y,jj real*8 :: term,t2 double precision :: bielec_pqxx_no,bielec_pxxq_no - + ii=list_core_inact(i) tt=list_act(t) if (i.eq.j) then if (t.eq.u) then ! diagonal element - term=occnum(tt)*Fipq(ii,ii)+2.D0*(Fipq(tt,tt)+Fapq(tt,tt)) & + term = occnum(tt)*Fipq(ii,ii) + & + 2.D0*(Fipq(tt,tt)+Fapq(tt,tt)) & -2.D0*(Fipq(ii,ii)+Fapq(ii,ii)) term+=2.D0*(3.D0*bielec_pxxq_no(tt,i,i,tt)-bielec_pqxx_no(tt,tt,i,i)) term-=2.D0*occnum(tt)*(3.D0*bielec_pxxq_no(tt,i,i,tt) & @@ -83,10 +84,10 @@ real*8 function hessmat_itju(i,t,j,u) end do end do end if - + term*=2.D0 hessmat_itju=term - + end function hessmat_itju real*8 function hessmat_itja(i,t,j,a) @@ -97,7 +98,7 @@ real*8 function hessmat_itja(i,t,j,a) integer :: i,t,j,a,ii,tt,jj,aa,v,vv,x,y real*8 :: term double precision :: bielec_pqxx_no,bielec_pxxq_no - + ! it/ja ii=list_core_inact(i) tt=list_act(t) @@ -120,7 +121,7 @@ real*8 function hessmat_itja(i,t,j,a) end if term*=2.D0 hessmat_itja=term - + end function hessmat_itja real*8 function hessmat_itua(i,t,u,a) @@ -131,7 +132,7 @@ real*8 function hessmat_itua(i,t,u,a) integer :: i,t,u,a,ii,tt,uu,aa,v,vv,x,xx,u3,t3,v3 real*8 :: term double precision :: bielec_pqxx_no,bielec_pxxq_no - + ii=list_core_inact(i) tt=list_act(t) t3=t+n_core_inact_orb @@ -162,7 +163,7 @@ real*8 function hessmat_itua(i,t,u,a) end if term*=2.D0 hessmat_itua=term - + end function hessmat_itua real*8 function hessmat_iajb(i,a,j,b) @@ -173,7 +174,7 @@ real*8 function hessmat_iajb(i,a,j,b) integer :: i,a,j,b,ii,aa,jj,bb real*8 :: term double precision :: bielec_pqxx_no,bielec_pxxq_no - + ii=list_core_inact(i) aa=list_virt(a) if (i.eq.j) then @@ -199,7 +200,7 @@ real*8 function hessmat_iajb(i,a,j,b) end if term*=2.D0 hessmat_iajb=term - + end function hessmat_iajb real*8 function hessmat_iatb(i,a,t,b) @@ -210,7 +211,7 @@ real*8 function hessmat_iatb(i,a,t,b) integer :: i,a,t,b,ii,aa,tt,bb,v,vv,x,y,v3,t3 real*8 :: term double precision :: bielec_pqxx_no,bielec_pxxq_no - + ii=list_core_inact(i) aa=list_virt(a) tt=list_act(t) @@ -231,7 +232,7 @@ real*8 function hessmat_iatb(i,a,t,b) end if term*=2.D0 hessmat_iatb=term - + end function hessmat_iatb real*8 function hessmat_taub(t,a,u,b) @@ -243,7 +244,7 @@ real*8 function hessmat_taub(t,a,u,b) integer :: v3,x3 real*8 :: term,t1,t2,t3 double precision :: bielec_pqxx_no,bielec_pxxq_no - + tt=list_act(t) aa=list_virt(a) if (t == u) then @@ -311,12 +312,12 @@ real*8 function hessmat_taub(t,a,u,b) end do end do end if - + end if - + term*=2.D0 hessmat_taub=term - + end function hessmat_taub BEGIN_PROVIDER [real*8, hessdiag, (nMonoEx)] @@ -326,7 +327,7 @@ BEGIN_PROVIDER [real*8, hessdiag, (nMonoEx)] implicit none integer :: i,t,a,indx,indx_shift real*8 :: hessmat_itju,hessmat_iajb,hessmat_taub - + !$OMP PARALLEL DEFAULT(NONE) & !$OMP SHARED(hessdiag,n_core_inact_orb,n_act_orb,n_virt_orb,nMonoEx) & !$OMP PRIVATE(i,indx,t,a,indx_shift) @@ -339,9 +340,9 @@ BEGIN_PROVIDER [real*8, hessdiag, (nMonoEx)] end do end do !$OMP END DO NOWAIT - + indx_shift = n_core_inact_orb*n_act_orb - !$OMP DO + !$OMP DO do a=1,n_virt_orb do i=1,n_core_inact_orb indx = a + (i-1)*n_virt_orb + indx_shift @@ -349,9 +350,9 @@ BEGIN_PROVIDER [real*8, hessdiag, (nMonoEx)] end do end do !$OMP END DO NOWAIT - + indx_shift += n_core_inact_orb*n_virt_orb - !$OMP DO + !$OMP DO do a=1,n_virt_orb do t=1,n_act_orb indx = a + (t-1)*n_virt_orb + indx_shift @@ -360,7 +361,7 @@ BEGIN_PROVIDER [real*8, hessdiag, (nMonoEx)] end do !$OMP END DO !$OMP END PARALLEL - + END_PROVIDER @@ -377,7 +378,7 @@ BEGIN_PROVIDER [double precision, hessmat, (nMonoEx,nMonoEx)] real*8 :: hessmat_taub ! c-a c-v a-v ! c-a | X X X - ! c-v | X X + ! c-v | X X ! a-v | X provide all_mo_integrals @@ -390,12 +391,12 @@ BEGIN_PROVIDER [double precision, hessmat, (nMonoEx,nMonoEx)] !$OMP DO !!!! < Core-active| H |Core-active > - ! Core-active excitations + ! Core-active excitations do indx_tmp = 1, n_c_a_prov indx = list_idx_c_a(1,indx_tmp) i = list_idx_c_a(2,indx_tmp) t = list_idx_c_a(3,indx_tmp) - ! Core-active excitations + ! Core-active excitations do j = 1, n_core_inact_orb if (i.eq.j) then ustart=t @@ -418,12 +419,12 @@ BEGIN_PROVIDER [double precision, hessmat, (nMonoEx,nMonoEx)] !$OMP DO !!!! < Core-active| H |Core-VIRTUAL > - ! Core-active excitations + ! Core-active excitations do indx_tmp = 1, n_c_a_prov indx = list_idx_c_a(1,indx_tmp) i = list_idx_c_a(2,indx_tmp) t = list_idx_c_a(3,indx_tmp) - ! Core-VIRTUAL excitations + ! Core-VIRTUAL excitations do jndx_tmp = 1, n_c_v_prov jndx = list_idx_c_v(1,jndx_tmp) j = list_idx_c_v(2,jndx_tmp) @@ -441,12 +442,12 @@ BEGIN_PROVIDER [double precision, hessmat, (nMonoEx,nMonoEx)] !$OMP DO !!!! < Core-active| H |ACTIVE-VIRTUAL > - ! Core-active excitations + ! Core-active excitations do indx_tmp = 1, n_c_a_prov indx = list_idx_c_a(1,indx_tmp) i = list_idx_c_a(2,indx_tmp) t = list_idx_c_a(3,indx_tmp) - ! ACTIVE-VIRTUAL excitations + ! ACTIVE-VIRTUAL excitations do jndx_tmp = 1, n_a_v_prov jndx = list_idx_a_v(1,jndx_tmp) u = list_idx_a_v(2,jndx_tmp) @@ -466,12 +467,12 @@ BEGIN_PROVIDER [double precision, hessmat, (nMonoEx,nMonoEx)] !$OMP PRIVATE(indx_tmp,indx,i,a,j,b,bstart,jndx) !$OMP DO !!!!! < Core-VIRTUAL | H |Core-VIRTUAL > - ! Core-VIRTUAL excitations + ! Core-VIRTUAL excitations do indx_tmp = 1, n_c_v_prov indx = list_idx_c_v(1,indx_tmp) i = list_idx_c_v(2,indx_tmp) a = list_idx_c_v(3,indx_tmp) - ! Core-VIRTUAL excitations + ! Core-VIRTUAL excitations do j = 1, n_core_inact_orb if (i.eq.j) then bstart=a @@ -485,7 +486,7 @@ BEGIN_PROVIDER [double precision, hessmat, (nMonoEx,nMonoEx)] enddo enddo enddo - + !$OMP END DO NOWAIT !$OMP END PARALLEL endif @@ -496,12 +497,12 @@ BEGIN_PROVIDER [double precision, hessmat, (nMonoEx,nMonoEx)] !$OMP DO !!!! < Core-VIRTUAL | H |Active-VIRTUAL > - ! Core-VIRTUAL excitations + ! Core-VIRTUAL excitations do indx_tmp = 1, n_c_v_prov indx = list_idx_c_v(1,indx_tmp) i = list_idx_c_v(2,indx_tmp) a = list_idx_c_v(3,indx_tmp) - ! Active-VIRTUAL excitations + ! Active-VIRTUAL excitations do jndx_tmp = 1, n_a_v_prov jndx = list_idx_a_v(1,jndx_tmp) t = list_idx_a_v(2,jndx_tmp) @@ -520,12 +521,12 @@ BEGIN_PROVIDER [double precision, hessmat, (nMonoEx,nMonoEx)] !$OMP DO !!!! < Active-VIRTUAL | H |Active-VIRTUAL > - ! Active-VIRTUAL excitations + ! Active-VIRTUAL excitations do indx_tmp = 1, n_a_v_prov indx = list_idx_a_v(1,indx_tmp) t = list_idx_a_v(2,indx_tmp) a = list_idx_a_v(3,indx_tmp) - ! Active-VIRTUAL excitations + ! Active-VIRTUAL excitations do u=t,n_act_orb if (t.eq.u) then bstart=a @@ -542,4 +543,4 @@ BEGIN_PROVIDER [double precision, hessmat, (nMonoEx,nMonoEx)] !$OMP END DO NOWAIT !$OMP END PARALLEL -END_PROVIDER +END_PROVIDER diff --git a/src/casscf_cipsi/natorb.irp.f b/src/casscf_cipsi/natorb.irp.f index 9ce90304..6376308d 100644 --- a/src/casscf_cipsi/natorb.irp.f +++ b/src/casscf_cipsi/natorb.irp.f @@ -72,84 +72,27 @@ BEGIN_PROVIDER [real*8, P0tuvx_no, (n_act_orb,n_act_orb,n_act_orb,n_act_orb)] BEGIN_DOC ! 4-index transformation of 2part matrices END_DOC - integer :: i,j,k,l,p,q - real*8 :: d(n_act_orb) - ! index per index - ! first quarter - P0tuvx_no(:,:,:,:) = P0tuvx(:,:,:,:) + double precision, allocatable :: tmp(:,:,:,:) + allocate(tmp(n_act_orb,n_act_orb,n_act_orb,n_act_orb)) - do j=1,n_act_orb - do k=1,n_act_orb - do l=1,n_act_orb - do p=1,n_act_orb - d(p)=0.D0 - end do - do p=1,n_act_orb - do q=1,n_act_orb - d(p)+=P0tuvx_no(q,j,k,l)*natorbsCI(q,p) - end do - end do - do p=1,n_act_orb - P0tuvx_no(p,j,k,l)=d(p) - end do - end do - end do - end do - ! 2nd quarter - do j=1,n_act_orb - do k=1,n_act_orb - do l=1,n_act_orb - do p=1,n_act_orb - d(p)=0.D0 - end do - do p=1,n_act_orb - do q=1,n_act_orb - d(p)+=P0tuvx_no(j,q,k,l)*natorbsCI(q,p) - end do - end do - do p=1,n_act_orb - P0tuvx_no(j,p,k,l)=d(p) - end do - end do - end do - end do - ! 3rd quarter - do j=1,n_act_orb - do k=1,n_act_orb - do l=1,n_act_orb - do p=1,n_act_orb - d(p)=0.D0 - end do - do p=1,n_act_orb - do q=1,n_act_orb - d(p)+=P0tuvx_no(j,k,q,l)*natorbsCI(q,p) - end do - end do - do p=1,n_act_orb - P0tuvx_no(j,k,p,l)=d(p) - end do - end do - end do - end do - ! 4th quarter - do j=1,n_act_orb - do k=1,n_act_orb - do l=1,n_act_orb - do p=1,n_act_orb - d(p)=0.D0 - end do - do p=1,n_act_orb - do q=1,n_act_orb - d(p)+=P0tuvx_no(j,k,l,q)*natorbsCI(q,p) - end do - end do - do p=1,n_act_orb - P0tuvx_no(j,k,l,p)=d(p) - end do - end do - end do - end do + call dgemm('T','N',(n_act_orb*n_act_orb*n_act_orb), n_act_orb, n_act_orb, 1.d0, & + P0tuvx, n_act_orb, natorbsCI, n_act_orb, 0.d0, & + tmp, (n_act_orb*n_act_orb*n_act_orb)) + + call dgemm('T','N',(n_act_orb*n_act_orb*n_act_orb), n_act_orb, n_act_orb, 1.d0, & + tmp, n_act_orb, natorbsCI, n_act_orb, 0.d0, & + P0tuvx_no, (n_act_orb*n_act_orb*n_act_orb)) + + call dgemm('T','N',(n_act_orb*n_act_orb*n_act_orb), n_act_orb, n_act_orb, 1.d0, & + P0tuvx_no, n_act_orb, natorbsCI, n_act_orb, 0.d0, & + tmp, (n_act_orb*n_act_orb*n_act_orb)) + + call dgemm('T','N',(n_act_orb*n_act_orb*n_act_orb), n_act_orb, n_act_orb, 1.d0, & + tmp, n_act_orb, natorbsCI, n_act_orb, 0.d0, & + P0tuvx_no, (n_act_orb*n_act_orb*n_act_orb)) + + deallocate(tmp) END_PROVIDER @@ -160,6 +103,7 @@ BEGIN_PROVIDER [real*8, one_ints_no, (mo_num,mo_num)] BEGIN_DOC ! Transformed one-e integrals END_DOC + integer :: i,j, p, q real*8 :: d(n_act_orb) one_ints_no(:,:)=mo_one_e_integrals(:,:) @@ -168,10 +112,8 @@ BEGIN_PROVIDER [real*8, one_ints_no, (mo_num,mo_num)] do j=1,mo_num do p=1,n_act_orb d(p)=0.D0 - end do - do p=1,n_act_orb do q=1,n_act_orb - d(p)+=one_ints_no(list_act(q),j)*natorbsCI(q,p) + d(p) = d(p) + one_ints_no(list_act(q),j)*natorbsCI(q,p) end do end do do p=1,n_act_orb @@ -183,8 +125,6 @@ BEGIN_PROVIDER [real*8, one_ints_no, (mo_num,mo_num)] do j=1,mo_num do p=1,n_act_orb d(p)=0.D0 - end do - do p=1,n_act_orb do q=1,n_act_orb d(p)+=one_ints_no(j,list_act(q))*natorbsCI(q,p) end do diff --git a/src/mo_two_e_ints/map_integrals.irp.f b/src/mo_two_e_ints/map_integrals.irp.f index b5f78b7b..5a005d7b 100644 --- a/src/mo_two_e_ints/map_integrals.irp.f +++ b/src/mo_two_e_ints/map_integrals.irp.f @@ -86,10 +86,8 @@ BEGIN_PROVIDER [ double precision, mo_integrals_cache, (0_8:mo_integrals_cache_s call set_multiple_levels_omp(.False.) - !$OMP PARALLEL DO PRIVATE(k,l,ii) SCHEDULE(dynamic) do l=mo_integrals_cache_min,mo_integrals_cache_max - print *, l do k=mo_integrals_cache_min,mo_integrals_cache_max ii = int(l-mo_integrals_cache_min,8) ii = ior( shiftl(ii,mo_integrals_cache_shift), int(k-mo_integrals_cache_min,8))