Skip to content

Commit

Permalink
further cleanup for clarity, comments, and V.T. removal of use allmod
Browse files Browse the repository at this point in the history
  • Loading branch information
etiennepalos committed May 14, 2024
1 parent 665fb8a commit c6955ab
Showing 1 changed file with 32 additions and 19 deletions.
51 changes: 32 additions & 19 deletions src/modules/quick_oeproperties_module.f90
Original file line number Diff line number Diff line change
Expand Up @@ -5,7 +5,8 @@
#endif

!---------------------------------------------------------------------!
! Created by Etienne Palos on xx/xx/2024 !
! Created by Etienne Palos on 01/20/2024 !
! Contributor: Vikrant Tripathy !
! !
! Copyright (C) 2024-2025 Götz lab !
! !
Expand All @@ -22,12 +23,22 @@ module quick_oeproperties_module
contains

!----------------------------------------------------------------------------!
! This is the main subroutine that computes the Electrostatic Potential (ESP)!
! at a given point , V(r) = V_nuc(r) + V_elec(r), and prints it to file.prop !
! This is the subroutine that "computes" the Electrostatic Potential (ESP) !
! at a given point , V(r) = V_nuc(r) + V_elec(r), and prints it to file.prop !
! !
! This subroutine is called from the main program. !
! It calls the following subroutines: !
! 1. esp_nuc: Computes the nuclear contribution to the ESP !
! 2. esp_shell_pair: Computes the electronic contribution to the ESP !
! 3. print_esp: Prints the ESP to the output file.prop !
!----------------------------------------------------------------------------!
subroutine compute_esp(ierr)
use allmod

use quick_mpi_module, only: master
use quick_timer_module, only : timer_begin, timer_end, timer_cumer
use quick_molspec_module, only : quick_molspec
use quick_files_module, only : iPropFile, propFileName
use quick_basis_module, only: jshell

#ifdef MPIV
use mpi
#endif
Expand Down Expand Up @@ -76,15 +87,15 @@ subroutine compute_esp(ierr)
deallocate(esp_nuclear)

close(iPropFile)

end subroutine compute_esp

!---------------------------------------------------------------------------------------------!
! This subroutine formats and prints the Electrostatic Potential (ESP) to the output file.prop!
!---------------------------------------------------------------------------------------------!
subroutine print_esp(esp_nuclear, esp_electronic, ierr)
use quick_molspec_module
use quick_method_module
use quick_mpi_module, only: master
use quick_molspec_module, only: quick_molspec
use quick_method_module, only: quick_method
use quick_files_module, only: ioutfile, iPropFile
use quick_mpi_module, only: master
use quick_constants_module, only: BOHRS_TO_A
Expand Down Expand Up @@ -167,16 +178,15 @@ end subroutine print_esp
! This subroutine calculates V_nuc(r) = sum Z_k/|r-Rk| !
!-----------------------------------------------------------------------!
subroutine esp_nuc(ierr, igridpoint, esp_nuclear_term)
use quick_molspec_module
!, only: extxyz, chg
!use allmod
use quick_molspec_module, only: natom, quick_molspec, xyz

implicit none

integer, intent(inout) :: ierr

double precision :: distance
double precision, external :: rootSquare
integer i,j
integer inucleus

double precision, intent(inout) :: esp_nuclear_term
integer ,intent(in) :: igridpoint
Expand All @@ -185,9 +195,9 @@ subroutine esp_nuc(ierr, igridpoint, esp_nuclear_term)

esp_nuclear_term = 0.d0

do i=1,natom
distance = rootSquare(xyz(1:3,i), quick_molspec%extxyz(1:3,igridpoint), 3)
esp_nuclear_term = esp_nuclear_term + quick_molspec%chg(i) / distance
do inucleus=1,natom
distance = rootSquare(xyz(1:3,inucleus), quick_molspec%extxyz(1:3,igridpoint), 3)
esp_nuclear_term = esp_nuclear_term + quick_molspec%chg(inucleus) / distance
enddo
end subroutine esp_nuc

Expand All @@ -199,7 +209,10 @@ end subroutine esp_nuc
! Then, P_{mu nu} * 〈 phi_mu | 1/|r-C| | phi_nu 〉
!-----------------------------------------------------------------------------------------
subroutine esp_1pdm(Ips,Jps,IIsh,JJsh,NIJ1,Ax,Ay,Az,Bx,By,Bz,Cx,Cy,Cz,Px,Py,Pz,esp)
use allmod
! use allmod
use quick_params_module, only: trans
use quick_calculated_module, only : quick_qm_struct
use quick_basis_module, only: attraxiao,quick_basis
use quick_files_module, only: ioutfile

implicit double precision(a-h,o-z)
Expand Down Expand Up @@ -438,12 +451,12 @@ End subroutine esp_1pdm
! Then, P_{mu nu} * 〈 phi_mu | 1/|r-C| | phi_nu 〉
!-----------------------------------------------------------------------------------------
subroutine esp_shell_pair(IIsh, JJsh, esp_electronic)
use allmod ! revisit, avoid allmod
use quick_molspec_module
use quick_method_module, only: quick_method
use quick_basis_module, only: quick_basis, attraxiao
use quick_molspec_module, only: quick_molspec, xyz
use quick_overlap_module, only: gpt, opf, overlap_core
use quick_files_module, only: ioutfile, iPropFile
use quick_mpi_module, only: master
implicit double precision(a-h,o-z)

dimension aux(0:20)
double precision AA(3),BB(3),CC(3),PP(3)
Expand Down

0 comments on commit c6955ab

Please sign in to comment.