2024-07-25 10:27:17 +02:00

440 lines
19 KiB
Fortran
Executable File

!***********************************************************************
subroutine pcm_core(pcm,epol,verbosity,iout,dcore,icore,imem,nbasis,focka,fockb,pa,pb,scftype,dero,grad)
!***********************************************************************
! Calculate solvent effects using an interface to PCMSolver
!***********************************************************************
implicit none
integer natoms,verbosity,iout,imem,ncent,icoord,dblalloc,iatchg,i
integer nbasis,scftype,nangmax,ncontrmax,nprimmax,nbset,dero
integer ncartmax,nsphermax,icore(*),inang,incontr,inprim,igexp
integer igcoef,ictostr,inangmin,igcn,inshrange,iboysval,iccf,iamat
integer intalloc,nmboys,idpre,iindarr
real*8 epol,dcore(*),focka(nbasis,nbasis),fockb(nbasis,nbasis)
real*8 pa(nbasis,nbasis),pb(nbasis,nbasis),itol,grad(*)
character(len=4) citol
character(len=32) pcm
logical cartg
!
if(verbosity.ge.2) then
if(dero.eq.0) then
write(iout,*) 'Calculating solvent polarization...'
else
write(iout,*)
write(iout,*) &
'Calculating solvent contribution to Cartesian gradient...'
endif
endif
! Read variables and allocate memory for grid-independent arrays
call getkey('itol',4,citol,4)
read(citol,*) i
itol=10.d0**(-i)
call getvar('nbset ',nbset)
call getvar('natoms ',natoms)
call getvar('ncent ',ncent)
call getvar('nangmax ',nangmax)
call getvar('ncontrmax ',ncontrmax)
call getvar('nprimmax ',nprimmax)
call getvar('ncartmax ',ncartmax)
call getvar('nsphermax ',nsphermax)
call getvar('cartg ',cartg)
call getvar('nmboys ',nmboys)
inang=intalloc(natoms*nbset)
call getvar('nang ',icore(inang))
incontr=intalloc((nangmax+1)*natoms*nbset)
call getvar('ncontr ',icore(incontr))
inprim=intalloc((nangmax+1)*natoms*nbset)
call getvar('nprim ',icore(inprim))
igexp=dblalloc((nangmax+1)*nprimmax*natoms*nbset)
call getvar('gexp ',dcore(igexp))
igcoef=dblalloc(nprimmax*ncontrmax*(nangmax+1)*natoms*nbset)
call getvar('gcoef ',dcore(igcoef))
icoord=dblalloc(3*ncent)
call getvar('coord ',dcore(icoord))
ictostr=dblalloc((nangmax+1)*ncartmax**2)
call getvar('ctostr ',dcore(ictostr))
inangmin=intalloc(natoms*nbset)
call getvar('nangmin ',icore(inangmin))
igcn=intalloc(2*ncontrmax*(nangmax+1)*natoms*nbset)
call getvar('gcn ',icore(igcn))
iatchg=dblalloc(ncent)
call getvar('atchg ',dcore(iatchg))
inshrange=intalloc(2*(nangmax+1)*natoms*nbset)
call getvar('nshrange ',icore(inshrange))
iboysval=dblalloc((nmboys+1)*1481)
call getvar('boysval ',dcore(iboysval))
iccf=dblalloc(nmboys+1)
call getvar('cf ',dcore(iccf))
iamat=dblalloc((2*dero+1)*nbasis**2)
if(dero.gt.0) then
idpre=dblalloc((nangmax+1)**2*natoms**2)
iindarr=intalloc(nsphermax*ncontrmax*(nangmax+1)*natoms*nbset)
call getvar('indarr ',icore(iindarr))
else
idpre=imem
iindarr=imem
endif
!
if(scftype.eq.2) call daxpy(nbasis**2,1.d0,pb,1,pa,1)
call pcmifc(pcm,epol,verbosity,natoms,dcore(icoord),dcore(iatchg),dcore,imem,iout,nbasis,focka,fockb,pa,scftype,nangmax,ncontrmax,nprimmax,ncartmax,nsphermax,icore(inang),icore(inangmin),icore(incontr),icore(inprim),dcore(igexp),dcore(igcoef),dcore(ictostr),dcore(iamat),icore(igcn),icore(inshrange),cartg,dcore(iboysval),nmboys,dcore(iccf),itol,dero,grad,dcore(idpre),icore(iindarr),ncent)
if(scftype.eq.2) call daxpy(nbasis**2,-1.d0,pb,1,pa,1)
!
call dbldealloc(inang)
if(verbosity.ge.2) call timer
write(iout,*)
!
end subroutine pcm_core
!***********************************************************************
subroutine pcmifc(pcm,epol,verbosity,natoms,coord,atchg,dcore,imem,iout,nbasis,focka,fockb,pa,scftype,nangmax,ncontrmax,nprimmax,ncartmax,nsphermax,nang,nangmin,ncontr,nprim,gexp,gcoef,ctostr,amat,gcn,nshrange,cartg,boysval,nmboys,ccf,itol,dero,grad,dpre,indarr,ncent)
!***********************************************************************
! Calculate solvent effects using an interface to PCMSolver
!***********************************************************************
use, intrinsic :: iso_c_binding
use pcmsolver
implicit none
integer natoms,ngrid,verbosity,igrid,imep,iasc,iout,nbasis,dero,ncent
integer scftype,imem,dblalloc,nangmax,ncontrmax,nprimmax,ncartmax
integer nang(natoms,*),ncontr(0:nangmax,natoms,*),nsphermax,nmboys
integer nangmin(natoms,*),nprim(0:nangmax,natoms,*),indarr
integer gcn(2,ncontrmax,0:nangmax,natoms,*)
integer nshrange(2,0:nangmax,natoms,*)
real*8 gexp(nprimmax,0:nangmax,natoms,*),ctostr,itol,ddot,grad(*)
real*8 gcoef(nprimmax,ncontrmax,0:nangmax,natoms,*),dpre(*)
real*8 amat(nbasis,nbasis),boysval,ccf,pa(nbasis,nbasis)
real*8 focka(nbasis,nbasis),fockb(nbasis,nbasis)
real(c_double) epol,atchg(ncent),coord(3*ncent),dcore(*)
character(kind=c_char,len=*), parameter :: mep_lbl='NucMEP'
character(kind=c_char,len=*), parameter :: asc_lbl='NucASC'
character(len=32) pcm
logical cartg
integer(c_int) ngrid_c
type(c_ptr) :: pcm_context
common/pcm/ pcm_context,ngrid_c
interface
subroutine host_writer(message) bind(c)
end subroutine host_writer
function fs2cs(fs) result(cs)
use, intrinsic :: iso_c_binding
implicit none
character(len=*), intent(in) :: fs
character(kind=c_char, len=1) :: cs(len(fs)+1)
end function
end interface
ngrid=ngrid_c
igrid=dblalloc(3*ngrid)
imep=dblalloc(ngrid)
iasc=dblalloc(ngrid)
dcore(igrid:igrid+3*ngrid-1)=0.d0
dcore(iasc :iasc + ngrid-1)=0.d0
call pcmsolver_get_centers(pcm_context,dcore(igrid:igrid+3*ngrid-1))
call mepcalc(natoms,atchg,coord,ngrid,dcore(igrid),dcore(imep),nangmax,ncontrmax,nprimmax,ncartmax,nsphermax,nang,nangmin,ncontr,nprim,gexp,gcoef,ctostr,nbasis,amat,dcore(imem),gcn,nshrange,cartg,boysval,nmboys,ccf,itol,pa,'mep',focka,fockb,scftype,dcore(iasc),dero,grad,dpre,indarr,ncent)
! call pcmsolver_set_surface_function(pcm_context,ngrid_c,dcore(imep:imep+ngrid-1),fs2cs(mep_lbl))
! call pcmsolver_compute_asc(pcm_context,fs2cs(mep_lbl),fs2cs(asc_lbl),irrep=0_c_int)
! call pcmsolver_get_surface_function(pcm_context,ngrid_c,dcore(iasc:iasc+ngrid-1),fs2cs(asc_lbl))
call pcmsolver_set_surface_function(pcm_context,ngrid_c,dcore(imep:imep+ngrid-1),mep_lbl)
call pcmsolver_compute_asc(pcm_context,mep_lbl,asc_lbl,irrep=0_c_int)
call pcmsolver_get_surface_function(pcm_context,ngrid_c,dcore(iasc:iasc+ngrid-1),asc_lbl)
call mepcalc(natoms,atchg,coord,ngrid,dcore(igrid),dcore(imep),nangmax,ncontrmax,nprimmax,ncartmax,nsphermax,nang,nangmin,ncontr,nprim,gexp,gcoef,ctostr,nbasis,amat,dcore(imem),gcn,nshrange,cartg,boysval,nmboys,ccf,itol,pa,'asc',focka,fockb,scftype,dcore(iasc),dero,grad,dpre,indarr,ncent)
epol=0.5d0*ddot(ngrid,dcore(imep),1,dcore(iasc),1)
if(verbosity.ge.3) write(iout,'(a, f20.12)') ' Polarization energy [au]: ',epol
! call pcmsolver_delete(pcm_context)
end subroutine pcmifc
!***********************************************************************
function pcmsolver_input(pcm,pcm1,minpfile) result(host_input)
!***********************************************************************
! Parse PCMSolver input and fill the data structure holding input data
!***********************************************************************
use, intrinsic :: iso_c_binding
use pcmsolver, only:PCMInput
!
implicit none
type(PCMInput) :: host_input
integer minpfile,pcm_der_order,i
character(len=1) pcm1(32)
character(len=32) pcm,c32
character(kind=c_char,len=32) :: PCM_Cavity_Type
character(kind=c_char,len=32) :: PCM_Cavity_RadiiSet
character(kind=c_char,len=32) :: PCM_Cavity_NpzFile
character(kind=c_char,len=32) :: PCM_Medium_SolverType
character(kind=c_char,len=32) :: PCM_Medium_Solvent
character(kind=c_char,len=32) :: pcm_equation_type
character(kind=c_char,len=32) :: pcm_inside_type
character(kind=c_char,len=32) :: PCM_Green_Type
logical(c_bool) :: PCM_Cavity_Scaling
real(c_double) :: PCM_Cavity_Area
real(c_double) :: pcm_min_distance
real(c_double) :: PCM_Cavity_MinRadius
real(c_double) :: PCM_Medium_Correction
real(c_double) :: PCM_Medium_ProbeRadius
real(c_double) :: PCM_Green_Eps
interface
function fs2cs(fs) result(cs)
use, intrinsic :: iso_c_binding
implicit none
character(len=*), intent(in) :: fs
character(kind=c_char, len=1) :: cs(len(fs)+1)
end function
end interface
!
do i=1,32
if(pcm1(i).eq.'_') pcm1(i)=' '
enddo
open(minpfile,file='MINP')
! Cavity/Type, dfault: GePol
PCM_Cavity_Type='gepol'//c_null_char
call getkeym('pcm_cavity_type',15,c32,32)
if(trim(c32).ne.'') PCM_Cavity_Type=trim(c32)//c_null_char
host_input%cavity_type=fs2cs(PCM_Cavity_Type)
! Cavity/Area, dfault: 0.3 A^2
PCM_Cavity_Area=0.3
call getkeym('pcm_cavity_area',15,c32,32)
if(trim(c32).ne.'') read(c32,*) PCM_Cavity_Area
host_input%area=PCM_Cavity_Area
! ???
pcm_min_distance=0.1
host_input%min_distance=pcm_min_distance
! ???
pcm_der_order=4
host_input%der_order=int(pcm_der_order, kind=c_int)
! Cavity/Scaling, dfault: False
PCM_Cavity_Scaling=.false.
call getkeym('pcm_cavity_scaling',18,c32,32)
if(trim(c32).eq.'false'.or.trim(c32).eq.'.false.') PCM_Cavity_Scaling=.false.
host_input%scaling=PCM_Cavity_Scaling
! Cavity/RadiiSet, dfault: UFF
PCM_Cavity_RadiiSet='uff'//c_null_char
call getkeym('pcm_cavity_radiiset',19,c32,32)
if(trim(c32).ne.'') PCM_Cavity_RadiiSet=trim(c32)//c_null_char
host_input%radii_set=fs2cs(PCM_Cavity_RadiiSet)
! Cavity/NpzFile, dfault: <empty string>
PCM_Cavity_NpzFile=''//c_null_char
call getkeym('pcm_cavity_npzfile',18,c32,32)
if(trim(c32).ne.'') PCM_Cavity_NpzFile=trim(c32)//c_null_char
host_input%restart_name=fs2cs(PCM_Cavity_NpzFile)
! Cavity/MinRadius, dfault: 100.0 A
PCM_Cavity_MinRadius=100.0
call getkeym('pcm_cavity_minradius',20,c32,32)
if(trim(c32).ne.'') read(c32,*) PCM_Cavity_MinRadius
host_input%min_radius=PCM_Cavity_MinRadius
! Medium/SolverType, dfault: IEFPCM
PCM_Medium_SolverType='iefpcm'//c_null_char
call getkeym('pcm_medium_solvertype',21,c32,32)
if(trim(c32).ne.'') PCM_Medium_SolverType=trim(c32)//c_null_char
host_input%solver_type=fs2cs(PCM_Medium_SolverType)
! Medium/Solvent, dfault: none
PCM_Medium_Solvent=trim(pcm)//c_null_char
call getkeym('pcm_medium_solvent',18,c32,32)
if(trim(c32).ne.'') PCM_Medium_Solvent=trim(c32)//c_null_char
host_input%solvent=fs2cs(PCM_Medium_Solvent)
! ???
pcm_equation_type='secondkind'//c_null_char
host_input%equation_type=fs2cs(pcm_equation_type)
! Medium/Correction, dfault: 0.0
PCM_Medium_Correction=0.0
call getkeym('pcm_medium_correction',21,c32,32)
if(trim(c32).ne.'') read(c32,*) PCM_Medium_Correction
host_input%correction=PCM_Medium_Correction
! Medium/ProbeRadius, dfault: 1.0 A
PCM_Medium_ProbeRadius=1.0
call getkeym('pcm_medium_proberadius',22,c32,32)
if(trim(c32).ne.'') read(c32,*) PCM_Medium_ProbeRadius
host_input%probe_radius=PCM_Medium_ProbeRadius
! ???
pcm_inside_type='vacuum'//c_null_char
host_input%inside_type=fs2cs(pcm_inside_type)
! Green/Eps, default: 1.0
PCM_Green_Eps=1.0
call getkeym('pcm_green_eps',13,c32,32)
if(trim(c32).ne.'') read(c32,*) PCM_Green_Eps
host_input%outside_epsilon=PCM_Green_Eps
! Green/Type, default: Vacuum
PCM_Green_Type='UniformDielectric'//c_null_char
call getkeym('pcm_green_type',14,c32,32)
if(trim(c32).ne.'') PCM_Green_Type=trim(c32)//c_null_char
host_input%outside_type=fs2cs(PCM_Green_Type)
!
close(minpfile)
end function pcmsolver_input
!***********************************************************************
subroutine mepcalc(natoms,atchg,coord,ngrid,grid,mep,nangmax,ncontrmax,nprimmax,ncartmax,nsphermax,nang,nangmin,ncontr,nprim,gexp,gcoef,ctostr,nbasis,amat,work,gcn,nshrange,cartg,boysval,nmboys,ccf,itol,pa,route,focka,fockb,scftype,asc,dero,grad,dpre,indarr,ncent)
!***********************************************************************
! Calculate molecular electrostatic potential
!***********************************************************************
implicit none
integer natoms,iatoms,ngrid,igrid,nbasis,scftype,dero,i,ncent
integer nangmax,ncontrmax,nprimmax,ncartmax,nsphermax,nmboys
integer nang(natoms,*),ncontr(0:nangmax,natoms,*)
integer nangmin(natoms,*),nprim(0:nangmax,natoms,*)
integer gcn(2,ncontrmax,0:nangmax,natoms,*)
integer nshrange(2,0:nangmax,natoms,*)
integer indarr(natoms,0:nangmax,ncontrmax,nsphermax)
real*8 dist,atchg(ncent),coord(3,ncent),grid(3,ngrid),mep(ngrid)
real*8 gexp(nprimmax,0:nangmax,natoms,*),ctostr,tmp,itol,expval
real*8 gcoef(nprimmax,ncontrmax,0:nangmax,natoms,*),asc(ngrid)
real*8 amat(nbasis,nbasis,3),work(*),boysval,ccf,pa(nbasis,nbasis)
real*8 focka(nbasis,nbasis),fockb(nbasis,nbasis),grad(3,natoms)
real*8 dpre(natoms,0:nangmax,natoms,0:nangmax),g1,g2,g3
logical cartg
character(len=3) route
!
if(route.eq.'mep') then
! Calculate nuclear contribution to MEP
!$OMP PARALLEL DO &
!$OMP SHARED(ngrid,ncent,coord,grid,atchg,mep) &
!$OMP PRIVATE(igrid,tmp,iatoms,dist)
do igrid=1,ngrid
tmp=0.d0
do iatoms=1,ncent
dist=dsqrt((coord(1,iatoms)-grid(1,igrid))**2 + &
(coord(2,iatoms)-grid(2,igrid))**2 + &
(coord(3,iatoms)-grid(3,igrid))**2)
tmp=tmp+atchg(iatoms)/dist
enddo
mep(igrid)=tmp
enddo
!$OMP END PARALLEL DO
! Calculate electronic contribution to MEP
do igrid=1,ngrid
call nucatt(natoms,nangmax,ncontrmax,nprimmax,ncartmax,&
nsphermax,nang,nangmin,ncontr,nprim,gexp,gcoef,coord,ctostr,&
nbasis,amat,work,gcn,nshrange,cartg,boysval,nmboys,ccf,itol,&
1.d0,1,grid(1,igrid))
mep(igrid)=mep(igrid)+expval(amat,pa,nbasis)
enddo
else
if(dero.eq.0) then
! Calculate contribution to Fock matrix
call nucatt(natoms,nangmax,ncontrmax,nprimmax,ncartmax,&
nsphermax,nang,nangmin,ncontr,nprim,gexp,gcoef,coord,ctostr,&
nbasis,amat,work,gcn,nshrange,cartg,boysval,nmboys,ccf,itol,&
asc,ngrid,grid)
call fillup(amat,nbasis)
call daxpy(nbasis**2,1.d0,amat,1,focka,1)
if(scftype.eq.2) call daxpy(nbasis**2,1.d0,amat,1,fockb,1)
else
! Calculate nuclear contribution to Cartesian gradient
!$OMP PARALLEL DO &
!$OMP SHARED(ngrid,ncent,coord,grid,atchg,asc,grad) &
!$OMP PRIVATE(igrid,tmp,iatoms,g1,g2,g3)
do iatoms=1,ncent
g1=0.d0
g2=0.d0
g3=0.d0
do igrid=1,ngrid
tmp=asc(igrid)/dsqrt((coord(1,iatoms)-grid(1,igrid))**2 + &
(coord(2,iatoms)-grid(2,igrid))**2 + &
(coord(3,iatoms)-grid(3,igrid))**2)**3
g1=g1+tmp*(coord(1,iatoms)-grid(1,igrid))
g2=g2+tmp*(coord(2,iatoms)-grid(2,igrid))
g3=g3+tmp*(coord(3,iatoms)-grid(3,igrid))
enddo
tmp=-0.5d0*atchg(iatoms)
grad(1,iatoms)=tmp*g1
grad(2,iatoms)=tmp*g2
grad(3,iatoms)=tmp*g3
enddo
!$OMP END PARALLEL DO
! Calculate electronic contribution to Cartesian gradient
call denspre(pa,pa,dpre,natoms,nangmax,ncontrmax,nsphermax, &
nbasis,cartg,indarr,nang,ncontr,'rhf ',.false.,.false.,.true.)
do iatoms=1,natoms
amat=0.d0
call nucatx(natoms,nangmax,ncontrmax,nprimmax,ncartmax, &
nsphermax,nang,nangmin,ncontr,nprim,gexp,gcoef,coord,ctostr, &
nbasis,amat(1,1,1),amat(1,1,2),amat(1,1,3),work,gcn,nshrange, &
cartg,boysval,nmboys,ccf,itol,asc,ngrid,iatoms,dpre,grid)
do i=1,3
grad(i,iatoms)= &
grad(i,iatoms)+0.5d0*expval(pa,amat(1,1,i),nbasis)
enddo
! write(6,"(3f12.5)") grad(1:3,iatoms)
enddo
endif
endif
!
end subroutine mepcalc
!***********************************************************************
subroutine host_writer(message) bind(c)
!***********************************************************************
! Writer for host code
!***********************************************************************
use, intrinsic :: iso_c_binding, only:c_char,c_null_char
implicit none
character(kind=c_char), intent(in) :: message(*)
integer n,i
n=0
do
if(message(n+1).eq.c_null_char) exit
n=n+1
enddo
write(6,'(1x,10000A)') (message(i),i=1,n)
end subroutine host_writer
!***********************************************************************
subroutine pcmini(pcm,verbosity,natoms,iout,minpfile)
!***********************************************************************
! Initialize PCM calculation
!***********************************************************************
use, intrinsic :: iso_c_binding
use pcmsolver
implicit none
integer natoms,verbosity,iout,minpfile,dblalloc
integer(c_int) natoms_c,ngrid_c
integer(c_int) :: symmetry_info(4)=(/0,0,0,0/)
real(c_double), allocatable :: atchg(:),coord(:)
type(PCMInput) :: host_input,pcmsolver_input
character(len=32) pcm
type(c_ptr) :: pcm_context
common/pcm/ pcm_context,ngrid_c
interface
subroutine host_writer(message) bind(c)
end subroutine host_writer
function fs2cs(fs) result(cs)
use, intrinsic :: iso_c_binding
implicit none
character(len=*), intent(in) :: fs
character(kind=c_char, len=1) :: cs(len(fs)+1)
end function
end interface
!
write(iout,*) 'Initializing PCM calculation...'
if(.not.pcmsolver_is_compatible_library()) then
write(iout,*) 'PCMSolver library is not compatible!'
call mrccend(1)
end if
!
allocate(atchg(natoms),coord(3*natoms))
call getvar('coord ',coord)
call getvar('atchg ',atchg)
natoms_c=natoms
host_input=pcmsolver_input(pcm,pcm,minpfile)
pcm_context=pcmsolver_new(PCMSOLVER_READER_HOST,natoms_c,atchg,coord,&
symmetry_info,host_input,c_funloc(host_writer)) ! Uncomment this for newer versions
ngrid_c=pcmsolver_get_cavity_size(pcm_context)
if(verbosity.ge.3) call pcmsolver_print(pcm_context)
deallocate(atchg,coord)
!
call timer
write(iout,*)
end subroutine pcmini
!***********************************************************************
function fs2cs(fs) result(cs)
!***********************************************************************
use, intrinsic :: iso_c_binding
implicit none
character(len=*), intent(in) :: fs
character(kind=c_char, len=1) :: cs(len(fs)+1)
integer :: i
do i = 1, len(fs)
cs(i) = fs(i:i)
end do
cs(i) = c_null_char
end function