easyconfigs-it4i/m/MRCC/mrcc_files/2/qmmod_handler.f
2024-07-25 10:27:17 +02:00

533 lines
24 KiB
Fortran

***********************************************************************
module qmmod_memory_handler
***********************************************************************
* Upgrades the memory model of dballoc for qmmod
* 2023-05-12 Bence Hégely
***********************************************************************
#include "MRCCCOMMON"
save
integer, external :: dblalloc,intalloc
integer, parameter :: stack_max_size = 64
type :: Stack_memory
integer :: capacity
integer :: isize
integer, pointer :: address( : )
integer, pointer :: memsize( : )
integer, pointer :: active ( : )
contains
procedure, pass :: init => memory_init
procedure, pass :: push => memory_push
procedure, pass :: a_top => memory_top_address
procedure, pass :: s_top => memory_top_memsize
procedure, pass :: top_r => memory_top_r
procedure, pass :: top_i => memory_top_i
procedure, pass :: pop => memory_pop
procedure, pass :: erase => memory_erase
procedure, pass :: erase_from => memory_erase_from
procedure, pass :: get_size
end type Stack_memory
!
contains
!
subroutine memory_init( this , capacity )
class(Stack_memory) :: this
integer, intent(in) :: capacity
call memalloc
this % capacity = capacity
this % isize = 0
call get_reference_i
$ ( icore( intalloc( capacity ) ) , this % address, capacity)
call get_reference_i
$ ( icore( intalloc( capacity ) ) , this % memsize, capacity)
call get_reference_i
$ ( icore( intalloc( capacity ) ) , this % active , capacity)
this % address(1:capacity) = 0
this % memsize(1:capacity) = 0
this % active (1:capacity) = 0
end subroutine memory_init
!
subroutine memory_push( this , n )
class(Stack_memory) :: this
integer, intent(in) :: n
integer :: iarr
if( stack_max_size .eq. this % isize ) then
write(iout,'(a)') ' Cannot manage more memory.'
call mrccend(1)
endif
iarr = dblalloc( n )
this % isize = this % isize + 1
this % address ( this % isize ) = iarr ! address
this % memsize ( this % isize ) = n ! size
this % active ( this % isize ) = 1 ! used
end subroutine memory_push
!
integer function memory_top_address( this ) result( res )
class(Stack_memory) :: this
res = this % address( this % isize )
end function memory_top_address
!
integer function memory_top_memsize( this ) result( res )
class(Stack_memory) :: this
res = this % memsize( this % isize )
end function memory_top_memsize
!
function memory_top_r( this ) result( array )
class(Stack_memory) :: this
double precision, pointer :: array(:)
integer :: iarr,n
iarr = this % a_top()
n = this % s_top()
call get_reference_r( dcore( iarr ) , array , n )
end function memory_top_r
!
function memory_top_i( this ) result( array )
class(Stack_memory) :: this
integer, pointer :: array(:)
integer :: iarr,n
iarr = this % a_top()
n = this % s_top()
call get_reference_i( icore( iarr ) , array , n )
end function memory_top_i
!
subroutine get_reference_r
$ ( data_target , data_pointer , nsize )
integer, intent(in) :: nsize
double precision, target :: data_target( nsize )
double precision, pointer :: data_pointer(:)
data_pointer => data_target
end subroutine get_reference_r
!
subroutine get_reference_i
$ ( data_target , data_pointer , nsize )
integer, intent(in) :: nsize
integer, target :: data_target( nsize )
integer, pointer :: data_pointer(:)
data_pointer => data_target
end subroutine get_reference_i
!
subroutine memory_pop( this )
class(Stack_memory) :: this
integer :: address
address = this % a_top()
this % address( this % isize ) = 0
this % memsize( this % isize ) = 0
this % active ( this % isize ) = 0
this % isize = this % isize - 1
call dbldealloc( address )
end subroutine memory_pop
!
subroutine memory_erase( this )
class(Stack_memory) :: this
do while( this % isize .gt. 0 )
call this % pop()
enddo
end subroutine memory_erase
!
subroutine memory_erase_from( this , address )
class(Stack_memory) :: this
integer, intent(in) :: address
do while( this % address( this % isize ) .ne. address )
call this % pop()
enddo
end subroutine memory_erase_from
!
integer function get_size( this ) result(res)
class(Stack_memory) :: this
res = this % isize
end function get_size
!
end module qmmod_memory_handler
!
***********************************************************************
module qmmod_handler
***********************************************************************
* Manages common features of qmmod-related procedures, which mostly
* builds around the execution of the Boughton-Pulay algorithm
* and the analysis of the BP-related results.
* 2023-05-12 Bence Hégely
***********************************************************************
use error_handler
use qmmod_memory_handler, only: Stack_memory, stack_max_size
use, intrinsic :: ISO_C_BINDING
#include "MRCCCOMMON"
save
type :: Keywords
character(len=6) :: core
character(len=8) :: corembed
character(len=8) :: embed
character(len=16) :: bpcompo
character(len=16) :: bpcompv
character(len=16) :: orbloce
character(len=16) :: orbloce_special(3)
character(len=8) :: qmreg
character(len=5) :: scftype
character(len=4) :: verbosity
end type Keywords
type :: BP_algorithm
integer :: nMO
integer :: nMO_nocc
double precision :: bpcompo
double precision :: bpcompv
integer :: max_natmo
logical :: return_mocomp
character(len=64) :: determinant
double precision, allocatable :: s_aux(:)
double precision, pointer :: s(:)
double precision, allocatable :: c_aux(:)
double precision, pointer :: c(:,:)
double precision, allocatable :: sa_aux(:)
double precision, pointer :: sa(:)
double precision, allocatable :: ai_aux(:)
double precision, pointer :: ai(:)
double precision, allocatable :: vi_aux(:)
double precision, pointer :: vi(:)
double precision, allocatable :: v_aux(:)
double precision, pointer :: v(:,:)
double precision, allocatable :: mocomp_aux(:)
double precision, pointer :: mocomp(:)
integer , allocatable :: natmo_aux(:)
integer , pointer :: natmo(:)
integer , allocatable :: atmo_aux(:)
integer , pointer :: atmo(:,:)
integer , allocatable :: nmoat_aux(:)
integer , pointer :: nmoat(:)
integer , allocatable :: moat_aux(:)
integer , pointer :: moat(:)
integer , allocatable :: ind_aux(:)
integer , pointer :: ind(:)
integer , allocatable :: natrange_aux(:)
integer , pointer :: natrange(:,:,:)
double precision, allocatable :: work_aux(:)
double precision, pointer :: work(:)
double precision, pointer :: work2d(:,:,:)
contains
procedure, pass :: get_resources => get_memory_for_bp_alg
procedure, pass :: init_data => init_data_for_bp_alg
procedure, pass :: execute => execute_bp_algorithm
end type BP_algorithm
class( Stack_memory ), allocatable :: smem
type( Keywords ) :: keywd
integer :: natoms
integer :: nbas
integer :: index_offset
logical :: do_multilevel_localcorrelation
logical :: do_projection_based_embedding
logical :: do_oniom_assist
logical :: do_frozen_MO_extraction_for_QMMM
logical :: do_projector_build_for_QMMM
integer :: iverbosity
integer :: bp_type_size
logical :: use_beta
character(len=32) :: route
integer, parameter :: buffer_size = 1024
! aux
logical :: lfound
integer :: istat
integer, allocatable :: ncorb(:)
contains
subroutine qmmod_handler_init
allocate( smem )
call smem % init( stack_max_size )
call get_route
call get_keywords
call set_task
call set_parameters
call read_var_arrays
end subroutine qmmod_handler_init
!
subroutine qmmod_handler_clear
call smem % erase()
deallocate( smem )
end subroutine qmmod_handler_clear
!
subroutine get_route
open(scrfile1,file='QMMOD_ROUTE',position='REWIND',
$ status='OLD',action='READ',iostat=istat)
if( istat .ne. 0 ) call io_error
$('Cannot found QMMOD_ROUTE','get_route @ qmmod_handler.f')
read(scrfile1,'(a)') route
close(scrfile1,status='DELETE')
end subroutine get_route
!
subroutine get_keywords
call getkey('scftype' , 7, keywd % scftype , 5)
call getkey('qmreg' , 5, keywd % qmreg , 8)
call getkey('embed' , 5, keywd % embed , 8)
call getkey('corembed' , 8, keywd % corembed, 8)
call getkey('core' , 4 ,keywd % core , 6)
call getkey('bpcompo' , 7 ,keywd % bpcompo , 16)
call getkey('bpcompv' , 7 ,keywd % bpcompv , 16)
call getkey('verbosity', 9 , keywd % verbosity, 4)
read( keywd % verbosity , * ) iverbosity
keywd % orbloce_special = ''
call getkey('orbloce' , 7 , keywd % orbloce ,16 )
if(keywd % orbloce.eq.'special') then
open(minpfile,file='MINP',form='FORMATTED',iostat=istat)
if(istat.ne.0) call io_error
$('Cannot open MINP file
$ in embedorb % MO_selection_init % read_keywords','bopu.f')
call getkeym('orbloce',7,keywd % orbloce,16)
read(minpfile,'(a)') keywd % orbloce_special(1)
read(minpfile,'(a)') keywd % orbloce_special(2)
read(minpfile,'(a)') keywd % orbloce_special(3)
close(minpfile, iostat=istat)
if(istat.ne.0) call io_error
$('Cannot close MINP file in
$ embedorb % MO_selection_init % read_keywords','bopu.f')
endif
end subroutine get_keywords
!
subroutine set_task
integer :: oniomcalc, qmreg
oniomcalc = -1
qmreg = -1
do_oniom_assist = .false.
inquire(file='ONIOMROUTE',exist=lfound)
if(lfound) then
open(scrfile1,file='ONIOMROUTE',form='UNFORMATTED')
read(scrfile1) oniomcalc
close(scrfile1)
do_oniom_assist = oniomcalc .eq. 1
endif
if( len_trim(keywd % qmreg) .ne. 0 )
$ read( keywd % qmreg , * ) qmreg
do_frozen_MO_extraction_for_QMMM
$ = len_trim(keywd % qmreg) .ne. 0
$ .and.
$ qmreg .ge. 1
do_projection_based_embedding = keywd % embed .ne. 'off'
do_multilevel_localcorrelation = keywd % corembed .ne. 'off'
if( route .eq. 'build_proj' ) then
do_projection_based_embedding = .false.
do_multilevel_localcorrelation = .false.
endif
do_projector_build_for_QMMM = keywd % qmreg .eq. '0'
$ .and.
$ .not.do_projection_based_embedding
end subroutine set_task
!
subroutine set_parameters
type( BP_algorithm ), allocatable :: bp
call getvar('nbasis ',nbasis)
call getvar('nal ',nal)
call getvar('nbe ',nbe)
call getvar('ncore ',ncore)
call getvar('nbf ',nbf)
call getvar('natoms ',natoms)
call getvar('nbset ',nbset)
nbas = nbf(1)
index_offset = 0
if( keywd % corembed .ne. 'off' .and.
$ keywd % core .ne. 'corr' ) then
index_offset = ncore
nbasis = nbasis - ncore
nal = nal - ncore
nbe = nbe - ncore
endif
nvirtal=nbasis-nal
nvirtbe=nbasis-nbe
nocc=max(nal,nbe)
nvirt=max(nvirtal,nvirtbe)
bp_type_size = sizeof( bp )
use_beta = .false.
end subroutine set_parameters
!
subroutine read_var_arrays
call smem % push( natoms )
ncorb = smem % top_i()
call getvar('ncorb ',ncorb)
end subroutine read_var_arrays
!
C
subroutine get_memory_for_bp_alg( this , bp_alpha )
class( BP_algorithm ), target :: this
class( BP_algorithm ), optional, target :: bp_alpha
read(keywd % bpcompo,*) this % bpcompo
read(keywd % bpcompv,*) this % bpcompv
if ( route .eq. 'select_only_occ' ) then
this % nMO = nal
this % nMO_nocc = nal
if( present( bp_alpha ) ) then
this % nMO = nbe
this % nMO_nocc = nbe
endif
else if( route .eq. 'select_only_virt') then
this % nMO = nvirtal
if( present( bp_alpha ) ) this % nMO = nvirtbe
else if( route .eq. 'select_occ+virt') then
this % nMO = nbas
this % nMO_nocc = nal
if( present( bp_alpha ) ) then
this % nMO_nocc = nbe
endif
endif
this % max_natmo = natoms
this % return_mocomp = .false.
call smem % push( nbas**2 )
this % c_aux = smem % top_r()
call c_f_pointer
$( c_loc( this % c_aux ), this % c , (/nbas,nbas/) )
call smem % push( nbas**2 )
this % v_aux = smem % top_r()
call c_f_pointer
$( c_loc( this % v_aux ) , this % v , (/ nbas, nbas /) )
call smem % push( natoms * nbas )
this % atmo_aux = smem % top_i()
call c_f_pointer
$( c_loc( this % atmo_aux ) , this % atmo , (/ natoms,nbas /) )
call smem % push( nbas )
this % natmo_aux = smem % top_i()
call c_f_pointer
$( c_loc( this % natmo_aux ) , this % natmo , (/ nbas /) )
! common
if( .not. present( bp_alpha ) ) then
call smem % push( nbas**2 )
this % s_aux = smem % top_r()
call c_f_pointer
$( c_loc( this % s_aux ) , this % s , (/ nbas**2 /) )
call smem % push( nbas**2 )
this % sa_aux = smem % top_r()
call c_f_pointer
$( c_loc( this % sa_aux ) , this % sa , (/ nbas**2 /) )
call smem % push( natoms )
this % ind_aux = smem % top_i()
call c_f_pointer
$( c_loc( this % ind_aux ) , this % ind , (/ natoms /) )
call smem % push( nbas )
this % vi_aux = smem % top_r()
call c_f_pointer
$( c_loc( this % vi_aux ) , this % vi , (/ nbas /) )
call smem % push( nbas )
this % ai_aux = smem % top_r()
call c_f_pointer
$( c_loc( this % ai_aux ) , this % ai , (/ nbas /) )
call smem % push( natoms )
this % nmoat_aux = smem % top_i()
call c_f_pointer
$( c_loc( this % nmoat_aux ) , this % nmoat , (/ natoms /) )
call smem % push( natoms * nbas )
this % moat_aux = smem % top_i()
call c_f_pointer
$( c_loc( this % moat_aux ) , this % moat , (/ natoms * nbas /) )
call smem % push( nbas )
this % mocomp_aux = smem % top_r()
call c_f_pointer
$( c_loc( this % mocomp_aux ) , this % mocomp , (/ nbas /) )
call smem % push( 2*natoms*nbset )
this % natrange_aux = smem % top_i()
call c_f_pointer
$( c_loc( this % natrange_aux ),this % natrange,(/2,natoms,nbset/))
call smem % push( 2 * nbas**2 )
this % work_aux = smem % top_r()
call c_f_pointer
$( c_loc( this % work_aux ),this % work,(/nbas*nbas*2/))
call c_f_pointer
$( c_loc( this % work_aux ),this % work2d,(/nbas,nbas,2/))
else
this % s => bp_alpha % s
this % sa => bp_alpha % sa
this % ind => bp_alpha % ind
this % vi => bp_alpha % vi
this % ai => bp_alpha % ai
this % nmoat => bp_alpha % nmoat
this % moat => bp_alpha % moat
this % mocomp => bp_alpha % mocomp
this % natrange => bp_alpha % natrange
this % work => bp_alpha % work
this % work2d => bp_alpha % work2d
endif
end subroutine get_memory_for_bp_alg
!
subroutine init_data_for_bp_alg( this , bp_alpha )
class( BP_algorithm ) :: this
class( BP_algorithm ), optional :: bp_alpha
if( .not. present( bp_alpha ) ) then
if( keywd % scftype .eq. 'rhf' ) this % determinant =
$ ' Restricted HF/KS'
if( keywd % scftype .eq. 'uhf' ) this % determinant =
$ ' Unrestricted HF/KS - ALPHA PARTICLES'
open(mocoeffile,file='MOCOEF',form='UNFORMATTED',
$ position='REWIND',status='OLD')
call readmo
$( this % work , this % work , this % c , mocoeffile , nbas , nbas)
close(mocoeffile)
open(oeintfile ,file='OEINT',form='UNFORMATTED',
$ position='REWIND',status='OLD')
call roeint
$( this % work , this % work, this % s , oeintfile, nbas )
close(oeintfile)
call dsymm
$('l','l',nbas,nbas,1.d0, this%s,nbas,this%c,nbas,0.d0,this%v,nbas)
call getvar('natrange ',this % natrange)
else
this % determinant =
$ ' Unrestricted HF/KS - BETA PARTICLES'
open(mocoeffile,file='MOCOEF',form='UNFORMATTED',
$ position='REWIND',status='OLD')
call readmo
$( this % work , this % work , this % c , mocoeffile , nbas , nbas)
call readmo
$( this % work , this % work , this % c , mocoeffile , nbas , nbas)
close(mocoeffile)
call dsymm
$('l','l',nbas,nbas,1.d0, this%s,nbas,this%c,nbas,0.d0,this%v,nbas)
endif
end subroutine init_data_for_bp_alg
!
subroutine execute_bp_algorithm( this )
class( BP_algorithm ) :: this
write(iout,'(a)')
write(iout,'(a)')
$' %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%'
write(iout,*)
write(iout,'(a)')
$' Executing the Boughton-Pulay algorithm '
write(iout,*)
write(iout,'(a)')
$' %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%'
write(iout,*)
write(iout,'(a,f6.4)')
$' with occupied MO completeness: ', this % bpcompo
write(iout,'(a,f6.4)')
$' with virtual MO completeness: ', this % bpcompv
write(iout,'(a,a)')
$' Determinant type : ', this % determinant
call bopu( ! INTERFACE
$nbas, ! integer :: nbasis
$this%nMO, ! integer :: nmatrix
$natoms, ! integer :: natoms
$this%natrange, ! integer :: natrange(2,natoms)
$this%c(1,1+index_offset), ! real*8 :: c(nbasis,nmatrix)
$this%v(1,1+index_offset), ! real*8 :: v(nbasis,nmatrix)
$this%ind, ! integer :: ind(natoms)
$this%work, ! real*8 :: r8heap(*)
$this%s, ! real*8 :: s(nbasis,nbasis))
$this%sa, ! real*8 :: sa(nbasis,nbasis)
$this%vi, ! real*8 :: vi(nbasis)
$this%ai, ! real*8 :: ai(nbasis)
$this%atmo, ! integer :: atmo(natoms,nmatrix)
$this%nmoat, ! integer :: nmoat(natoms)
$this%moat, ! integer :: moat(natoms,nmatrix)
$this%natmo, ! integer :: natmo(nmatrix)
$this%bpcompo, ! real*8 :: bpcompo
$this%bpcompv, ! real*8 :: bpcompv
$this%nMO_nocc, ! integer :: nocc
$this%max_natmo, ! integer :: natmaxo
$'bopu', ! character(len=4) :: meth
$this%work, ! integer :: i8heap(*)
$this%natrange, ! integer :: mnatrange(2,natoms)
$nbas, ! integer :: mnbasis
$this%c(1,1+index_offset), ! real*8 :: cc(nmatrix,nmatrix)
$this%work, ! real*8 :: distc(natoms,natoms)
$0, ! integer :: ndao
$0.d0, ! real*8 :: ovltol
$this%mocomp, ! real*8 :: mocomp(nmatrix)
$this%return_mocomp, ! logical :: pmul
$iout, ! integer :: iout
$0.05d0) ! real*8 :: minlc
write(iout,*)
write(iout,'(a)')
$' %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%'
end subroutine execute_bp_algorithm
end module qmmod_handler