get_gcp_param Subroutine

public subroutine get_gcp_param(param, mol, method, basis, eta)

Arguments

Type IntentOptional Attributes Name
type(gcp_param), intent(out) :: param
type(structure_type), intent(in) :: mol
character(len=*), intent(in), optional :: method
character(len=*), intent(in), optional :: basis
real(kind=wp), intent(in), optional :: eta

Source Code

subroutine get_gcp_param(param, mol, method, basis, eta)
   type(gcp_param), intent(out) :: param
   type(structure_type), intent(in) :: mol
   character(len=*), intent(in), optional :: method
   character(len=*), intent(in), optional :: basis
   real(wp), intent(in), optional :: eta

   real(wp) :: eta_, eta_spec
   real(wp), allocatable :: nbas(:)
   integer :: isp, izp, jsp, basis_id, method_id
   integer, allocatable :: nel(:)
   logical :: valence_minimal_basis

   method_id = p_unknown_method
   if (present(method)) then
      method_id = get_method_id(method)
      basis_id = get_basis_id(method)
   else
      basis_id = p_unknown_bas
   end if
   if (present(basis)) basis_id = get_basis_id(basis)
   if (basis_id == p_unknown_bas .and. method_id == p_unknown_method) return

   eta_ = 0.0_wp
   if (present(eta)) eta_ = eta
   eta_spec = 0.0_wp

   allocate(param%zeff(mol%nid))
   do isp = 1, mol%nid
      izp = mol%num(isp)
      param%zeff(isp) = effective_atomic_number(izp)
   end do

   allocate(param%rvdw(mol%nid, mol%nid))
   do isp = 1, mol%nid
      do jsp = 1, mol%nid
         param%rvdw(isp, jsp) = get_vdw_rad(param%zeff(isp), param%zeff(jsp))
      end do
   end do

   allocate(param%rvdw_srb(mol%nid, mol%nid))
   do isp = 1, mol%nid
      do jsp = 1, mol%nid
         param%rvdw_srb(isp, jsp) = get_vdw_rad(mol%num(isp), mol%num(jsp))
      end do
   end do

   valence_minimal_basis = basis_id == p_vmb_bas

   allocate(nel(mol%nid))
   do isp = 1, mol%nid
      izp = param%zeff(isp)
      nel(isp) = number_of_electrons(izp, valence_minimal_basis)
   end do

   select case(basis_id)
   case(p_vmb_bas)
      param%emiss = emiss_hf_vmb(param%zeff)
      nbas = nbas_vmb(param%zeff)

   case(p_sv_bas)
      param%emiss = emiss_hf_sv(param%zeff)
      nbas = nbas_sv(param%zeff)

      select case(method_id)
      case (p_hf_method) ! RMS=0.3218975
         param%sigma = 0.1724_wp
         eta_ = 1.2804_wp
         param%alpha = 0.8568_wp
         param%beta = 1.2342_wp
      case (p_hyb_method, p_b3lyp_method, p_pw6b95_method) ! RMS= 0.557
         param%sigma = 0.4048_wp
         eta_ = 1.1626_wp
         param%alpha = 0.8652_wp
         param%beta = 1.2375_wp
      case (p_gga_method, p_tpss_method, p_blyp_method) ! RMS = 0.6652
         param%sigma = 0.2727_wp
         eta_ = 1.4022_wp
         param%alpha = 0.8055_wp
         param%beta = 1.3000_wp
      end select

   case(p_sv_p_bas)
      param%emiss = emiss_hf_sv_p(param%zeff)
      nbas = nbas_sv_p(param%zeff)

      if (method_id == p_hf_method) then !RMS=0.3502
         param%sigma = 0.1373_wp
         eta_ = 1.4271_wp
         param%alpha = 0.8141_wp
         param%beta = 1.2760_wp
      elseif (is_dft_method(method_id)) then ! RMS= 0.57 ! def2-SV(P)
         param%sigma = 0.2424_wp
         eta_ = 1.2371_wp
         param%alpha = 0.6076_wp
         param%beta = 1.4078_wp
      end if

   case(p_svx_bas) ! RMS=  ! def2-SV(P/h,c)  = SV at h,c
      param%emiss = emiss_hf_svx(param%zeff)
      nbas = nbas_svx(param%zeff)

      if (is_dft_method(method_id)) then ! RMS=  0.56 ! def2-SV(P/h,c)  = SV at h,c
         param%sigma = 0.1861_wp
         eta_ = 1.3200_wp
         param%alpha = 0.6171_wp
         param%beta = 1.4019_wp
      end if

   case(p_svp_bas)
      param%emiss = emiss_hf_svp(param%zeff)
      nbas = nbas_svp(param%zeff)

      if (method_id == p_hf_method) then ! RMS=0.4065
         param%sigma = 0.2054_wp
         eta_ = 1.3157_wp
         param%alpha = 0.8136_wp
         param%beta = 1.2572_wp
      elseif (method_id == p_tpss_method) then ! RMS=  0.618
         param%sigma = 0.6647_wp
         eta_ = 1.3306_wp
         param%alpha = 1.0792_wp
         param%beta = 1.1651_wp
      elseif (method_id == p_pw6b95_method) then  ! RMS = 0.58312
         param%sigma = 0.3098_wp
         eta_ = 1.2373_wp
         param%alpha = 0.6896_wp
         param%beta = 1.3347_wp
      elseif (is_hyb_method(method_id)) then ! RMS=0.6498
         param%sigma = 0.2990_wp
         eta_ = 1.2605_wp
         param%alpha = 0.6438_wp
         param%beta = 1.3694_wp
      elseif (is_gga_method(method_id)) then ! RMS=
         param%sigma = 0.6823_wp
         eta_ = 1.2491_wp
         param%alpha = 0.8225_wp
         param%beta = 1.2811_wp
      end if

   case(p_svp_old_bas)
      param%emiss = emiss_hf_svp_old(param%zeff)
      nbas = nbas_svp_old(param%zeff)

      if (method_id == p_hf_method) then ! RMS=0.4065
         param%sigma = 0.2054_wp
         eta_ = 1.3157_wp
         param%alpha = 0.8136_wp
         param%beta = 1.2572_wp
      elseif (is_dft_method(method_id)) then ! RMS=0.6498
         param%sigma = 0.2990_wp
         eta_ = 1.2605_wp
         param%alpha = 0.6438_wp
         param%beta = 1.3694_wp
      end if

   case(p_minis_bas)
      param%emiss = emiss_hf_minis(param%zeff)
      nbas = nbas_minis(param%zeff)

      if (method_id == p_hf_method) then ! RMS= 0.3040
         param%sigma = 0.1290_wp
         eta_ = 1.1526_wp
         param%alpha = 1.1549_wp
         param%beta = 1.1763_wp
      elseif (method_id == p_tpss_method) then ! RMS=
         param%sigma = 0.22982_wp
         eta_ = 1.35401_wp
         param%alpha = 1.47633_wp
         param%beta = 1.11300_wp
      elseif (method_id == p_pw6b95_method) then  ! RMS = 0.3279929
         param%sigma = 0.21054_wp
         eta_ = 1.25458_wp
         param%alpha = 1.35003_wp
         param%beta = 1.14061_wp
      elseif (is_gga_method(method_id)) then ! RMS= 0.3462
         param%sigma = 0.1566_wp
         eta_ = 1.0271_wp
         param%alpha = 1.0732_wp
         param%beta = 1.1968_wp
      elseif (is_hyb_method(method_id)) then ! RMS= 0.3400
         param%sigma = 0.2059_wp
         eta_ = 0.9722_wp
         param%alpha = 1.1961_wp
         param%beta = 1.1456_wp
      end if

   case(p_631gd_bas)
      param%emiss = emiss_hf_631gd(param%zeff)
      nbas = nbas_631gd(param%zeff)

      if (method_id == p_hf_method) then ! RMS= 0.40476
         param%sigma = 0.2048_wp
         eta_ = 1.5652_wp
         param%alpha = 0.9447_wp
         param%beta = 1.2100_wp
      elseif (is_dft_method(method_id)) then ! RMS=  0.47856
         param%sigma = 0.3405_wp
         eta_ = 1.6127_wp
         param%alpha = 0.8589_wp
         param%beta = 1.2830_wp
      end if

   case(p_tz_bas)
      param%emiss = emiss_hf_tz(param%zeff)
      nbas = nbas_tz(param%zeff)

      if (method_id == p_hf_method) then !  RMS= 0.1150
         param%sigma = 0.3127_wp
         eta_ = 1.9914_wp
         param%alpha = 1.0216_wp
         param%beta = 1.2833_wp
      elseif (is_hyb_method(method_id)) then ! RMS=0.19648
         param%sigma = 0.2905_wp
         eta_ = 2.2495_wp
         param%alpha = 0.8120_wp
         param%beta = 1.4412_wp
      elseif (is_gga_method(method_id)) then !RMS = 0.21408
         param%sigma = 0.1182_wp
         eta_ = 1.0631_wp
         param%alpha = 1.0510_wp
         param%beta = 1.1287_wp
      end if

   case(p_deftzvp_bas)
      param%emiss = emiss_hf_def1tzvp(param%zeff)
      nbas = nbas_def1tzvp(param%zeff)

      if (method_id == p_hf_method) then ! RMS=0.209
         param%sigma = 0.2600_wp
         eta_ = 2.2448_wp
         param%alpha = 0.7998_wp
         param%beta = 1.4381_wp
      elseif (is_dft_method(method_id)) then ! RMS=0.1817
         param%sigma = 0.2393_wp
         eta_ = 2.2247_wp
         param%alpha = 0.8185_wp
         param%beta = 1.4298_wp
      end if

   case(p_ccdz_bas)
      param%emiss = emiss_hf_ccdz(param%zeff)
      nbas = nbas_ccdz(param%zeff)

      if (method_id == p_hf_method) then ! RMS=0.4968
         param%sigma = 0.4416_wp
         eta_ = 1.5185_wp
         param%alpha = 0.6902_wp
         param%beta = 1.3713_wp
      elseif (is_dft_method(method_id)) then ! RMS=0.7610
         param%sigma = 0.5383_wp
         eta_ = 1.6482_wp
         param%alpha = 0.6230_wp
         param%beta = 1.4523_wp
      end if

   case(p_accdz_bas)
      param%emiss = emiss_hf_accdz(param%zeff)
      nbas = nbas_accdz(param%zeff)

      if (method_id == p_hf_method) then !RMS=0.2222
         param%sigma = 0.0748_wp
         eta_ = 0.0663_wp
         param%alpha = 0.3811_wp
         param%beta = 1.0155_wp
      elseif (is_dft_method(method_id)) then ! RMS=0.1840
         param%sigma = 0.1465_wp
         eta_ = 0.0500_wp
         param%alpha = 0.6003_wp
         param%beta = 0.8761_wp
      end if

   case(p_pobtz_bas)
      param%emiss = emiss_hf_pobtz(param%zeff)
      nbas = nbas_pobtz(param%zeff)

      if (is_dft_method(method_id)) then
         param%sigma = 0.1300_wp
         eta_ = 1.3743_wp
         param%alpha = 0.4792_wp
         param%beta = 1.3962_wp
      end if

   ! case(p_pobdzvp_bas)
   !    param%emiss = emiss_hf_pobdzvp(param%zeff)
   !    nbas = nbas_pobdzvp(param%zeff)

   case(p_minix_bas)
      param%emiss = emiss_hf_minix(param%zeff)
      nbas = nbas_minix(param%zeff)

      if (method_id == p_hf_method .or. method_id == p_hf3c_method) then
         param%sigma = 0.1290_wp
         eta_ = 1.1526_wp
         param%alpha = 1.1549_wp
         param%beta = 1.1763_wp
         param%base = method_id == p_hf3c_method
         if (param%base) then
            param%rscal = 0.7_wp
            param%qscal = 0.03_wp
         end if
      elseif (is_dft_method(method_id)) then
         param%sigma = 0.2059_wp
         eta_ = 0.9722_wp
         param%alpha = 1.1961_wp
         param%beta = 1.1456_wp
      end if

   case(p_gcore_bas)
      param%emiss = emiss_hf_2gcore(param%zeff)
      nbas = nbas_2gcore(param%zeff)

   case(p_2g_bas)
      param%emiss = emiss_hf_2g(param%zeff)
      nbas = nbas_2g(param%zeff)

      if (method_id == p_hf_method) then
         param%sigma = 0.2461_wp
         eta_ = 1.1616_wp
         param%alpha = 0.7335_wp
         param%beta = 1.4709_wp
      end if

   case(p_dzp_bas)
      param%emiss = emiss_hf_dzp(param%zeff)
      nbas = nbas_dzp(param%zeff)

      if (method_id == p_hf_method) then !RMS=0.4571
         param%sigma = 0.1443_wp
         eta_ = 1.4547_wp
         param%alpha = 0.3711_wp
         param%beta = 1.6300_wp
      elseif (is_dft_method(method_id)) then !RMS=0.7184
         param%sigma = 0.2687_wp
         eta_ = 1.4634_wp
         param%alpha = 0.3513_wp
         param%beta = 1.6880_wp
      end if

   ! case(p_hsv_bas)
   !    param%emiss = emiss_hf_hsv(param%zeff)
   !    nbas = nbas_hsv(param%zeff)

   case(p_dz_bas)
      param%emiss = emiss_hf_dz(param%zeff)
      nbas = nbas_dz(param%zeff)

      if (method_id == p_hf_method) then  !RMS=0.3754
         param%sigma = 0.1059_wp
         eta_ = 1.4554_wp
         param%alpha = 0.3711_wp
         param%beta = 1.6342_wp
      elseif (is_dft_method(method_id)) then
         param%sigma = 0.2687_wp
         eta_ = 1.4634_wp
         param%alpha = 0.3513_wp
         param%beta = 1.6880_wp
      end if

   case(p_msvp_bas)
      param%emiss = emiss_hf_msvp(param%zeff)
      nbas = nbas_msvp(param%zeff)

   case(p_def2mtzvp_bas)
      param%emiss = emiss_hf_def2mtzvp(param%zeff)
      nbas = nbas_def2mtzvp(param%zeff)

      if (method_id == p_b3pbe3c_method) then
         param%sigma = 1.0000_wp
         eta_ = 2.98561_wp
         param%alpha = 0.3011_wp
         param%beta = 2.4405_wp
      end if

   case(p_pbeh3c_bas)
      param%emiss = emiss_hf_pbeh3c(param%zeff)
      nbas = nbas_msvp(param%zeff)

      if (method_id == p_pbeh3c_method) then
         param%sigma = 1.00000_wp
         eta_ = 1.32492_wp
         param%alpha = 0.27649_wp
         param%beta = 1.95600_wp
         param%damp = .true.
      elseif (method_id == p_hse3c_method) then
         param%sigma = 1.00000_wp
         eta_ = 1.32378_wp
         param%alpha = 0.28314_wp
         param%beta = 1.94527_wp
         param%damp = .true.
      end if

   case(p_lanl2_bas)
      param%emiss = emiss_hf_lanl2(param%zeff)
      nbas = nbas_lanl2(param%zeff)

      if (is_dft_method(method_id)) then
         param%sigma = 0.3405_wp
         eta_ = 1.6127_wp
         param%alpha = 0.8589_wp
         param%beta = 1.2830_wp
      end if

   case(p_def2mtzvpp_bas)
      param%emiss = emiss_hf_def2mtzvpp(param%zeff)
      nbas = nbas_def2mtzvpp(param%zeff)

      if (method_id == p_r2scan3c_method) then
         param%sigma = 1.0000_wp
         eta_ = 1.3150_wp
         eta_spec = 1.15_wp
         param%alpha = 0.9410_wp
         param%beta = 1.4636_wp
         param%damp = .true.
      end if

   end select
   if (allocated(nbas)) then
      param%xv = nbas - 0.5_wp * nel
      if (basis_id == p_def2mtzvpp_bas) then
         param%xv(:) = 1.0_wp
         where(param%zeff == 6) param%xv = 3.0_wp
         where(param%zeff == 7) param%xv = 0.5_wp
         where(param%zeff == 8) param%xv = 0.5_wp
      end if
   end if

   if (eta_ > 0.0_wp) then
      param%slater = eta_ * slater_exp(param%zeff)
      if (eta_spec > 0.0_wp) then
         where(param%zeff > 10) param%slater = eta_spec * param%slater
      end if
   end if

   if (method_id == p_b973c_method) then
      param%srb=.true.
      param%rscal=10.00_wp
      param%qscal=0.08_wp
   end if
end subroutine get_gcp_param