(git:b1f098b)
Loading...
Searching...
No Matches
min_basis_set.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief generate or use from input minimal basis set
10!> \par History
11!> 03.2023 created [JGH]
12!> \author JGH
13! **************************************************************************************************
16 USE basis_set_types, ONLY: &
21 USE kinds, ONLY: default_string_length,&
22 dp
24 ptable
28 USE qs_kind_types, ONLY: get_qs_kind,&
30#include "./base/base_uses.f90"
31
32 IMPLICIT NONE
33 PRIVATE
34
35 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'min_basis_set'
36
37 PUBLIC :: create_minbas_set
38
39! **************************************************************************************************
40
41CONTAINS
42
43! **************************************************************************************************
44!> \brief ...
45!> \param qs_env ...
46!> \param unit_nr ...
47!> \param basis_type ...
48!> \param primitive ...
49! **************************************************************************************************
50 SUBROUTINE create_minbas_set(qs_env, unit_nr, basis_type, primitive)
51 TYPE(qs_environment_type), POINTER :: qs_env
52 INTEGER, INTENT(IN) :: unit_nr
53 CHARACTER(LEN=*), INTENT(IN), OPTIONAL :: basis_type
54 INTEGER, INTENT(IN), OPTIONAL :: primitive
55
56 CHARACTER(LEN=2) :: element_symbol
57 CHARACTER(LEN=default_string_length) :: bname, btype
58 INTEGER :: ikind, lb, mao, ne, ngau, nkind, nprim, &
59 nsgf, ub, z
60 INTEGER, DIMENSION(0:3) :: elcon
61 INTEGER, DIMENSION(:), POINTER :: econf
62 REAL(kind=dp) :: zval
63 TYPE(dft_control_type), POINTER :: dft_control
64 TYPE(gto_basis_set_type), POINTER :: pbasis, refbasis
65 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
66 TYPE(qs_kind_type), POINTER :: qs_kind
67
68 IF (PRESENT(basis_type)) THEN
69 btype = trim(basis_type)
70 ELSE
71 btype = "MIN"
72 END IF
73 IF (PRESENT(primitive)) THEN
74 nprim = primitive
75 ELSE
76 nprim = -1
77 END IF
78
79 ! check for or generate reference basis
80 CALL get_qs_env(qs_env=qs_env, qs_kind_set=qs_kind_set)
81 CALL get_qs_env(qs_env, dft_control=dft_control)
82 nkind = SIZE(qs_kind_set)
83 DO ikind = 1, nkind
84 qs_kind => qs_kind_set(ikind)
85 CALL get_qs_kind(qs_kind, zeff=zval, elec_conf=econf, element_symbol=element_symbol)
86 CALL get_ptable_info(element_symbol, ielement=z)
87 NULLIFY (refbasis, pbasis)
88 CALL get_qs_kind(qs_kind=qs_kind, basis_set=refbasis, basis_type=btype)
89 IF (.NOT. ASSOCIATED(refbasis)) THEN
90 CALL get_qs_kind(qs_kind=qs_kind, mao=mao)
91 ! generate a minimal basis set
92 elcon = 0
93 lb = lbound(econf, 1)
94 ub = ubound(econf, 1)
95 ne = ub - lb
96 elcon(0:ne) = econf(lb:ub)
97 IF (nprim > 0) THEN
98 ngau = nprim
99 CALL create_min_basis(refbasis, z, elcon, mao, ngau)
100 CALL create_primitive_basis_set(refbasis, pbasis)
101 CALL init_interaction_radii_orb_basis(pbasis, dft_control%qs_control%eps_pgf_orb)
102 CALL add_basis_set_to_container(qs_kind%basis_sets, pbasis, btype)
103 CALL deallocate_gto_basis_set(refbasis)
104 ELSE
105 ! STO-3G
106 ngau = 3
107 CALL create_min_basis(refbasis, z, elcon, mao, ngau)
108 CALL init_interaction_radii_orb_basis(refbasis, dft_control%qs_control%eps_pgf_orb)
109 CALL add_basis_set_to_container(qs_kind%basis_sets, refbasis, btype)
110 END IF
111 CALL get_qs_kind(qs_kind=qs_kind, basis_set=refbasis, basis_type=btype)
112 END IF
113 IF (unit_nr > 0) THEN
114 CALL get_gto_basis_set(refbasis, name=bname, nsgf=nsgf)
115 WRITE (unit_nr, '(T2,A,A,T14,A,I4,T40,A,A24)') "Kind: ", element_symbol, "NBasFun: ", nsgf, &
116 "Reference Basis: ", adjustl(trim(bname))
117 END IF
118 END DO
119
120 END SUBROUTINE create_minbas_set
121
122! **************************************************************************************************
123!> \brief Creates a minimal basis set based on Slater rules
124!> \param min_basis_set ...
125!> \param zval ...
126!> \param econf ...
127!> \param mao ...
128!> \param ngau ...
129!> \par History
130!> 03.2023 created [JGH]
131! **************************************************************************************************
132 SUBROUTINE create_min_basis(min_basis_set, zval, econf, mao, ngau)
133 TYPE(gto_basis_set_type), POINTER :: min_basis_set
134 INTEGER, INTENT(IN) :: zval
135 INTEGER, DIMENSION(0:3) :: econf
136 INTEGER, INTENT(IN) :: mao, ngau
137
138 CHARACTER(len=1), DIMENSION(0:3), PARAMETER :: lnam = (/"S", "P", "D", "F"/)
139
140 CHARACTER(len=6) :: str
141 CHARACTER(len=6), DIMENSION(:), POINTER :: sym
142 CHARACTER(LEN=default_string_length) :: bname
143 INTEGER :: i, iss, l, lm, n, nmao, nn, nss
144 INTEGER, DIMENSION(0:3) :: nae, nco, npe
145 INTEGER, DIMENSION(4, 7) :: ne
146 INTEGER, DIMENSION(:), POINTER :: lq, nq
147 REAL(kind=dp), DIMENSION(:), POINTER :: zet
148 TYPE(sto_basis_set_type), POINTER :: sto_basis_set
149
150 cpassert(.NOT. ASSOCIATED(min_basis_set))
151 NULLIFY (sto_basis_set)
152
153 ! electronic configuration
154 ne = 0
155 DO l = 1, 4
156 nn = 2*(l - 1) + 1
157 DO i = l, 7
158 ne(l, i) = ptable(zval)%e_conv(l - 1) - 2*nn*(i - l)
159 ne(l, i) = max(ne(l, i), 0)
160 ne(l, i) = min(ne(l, i), 2*nn)
161 END DO
162 END DO
163 ! STO definition
164 nae = 0
165 npe = 0
166 DO l = 0, 3
167 nn = 2*(2*l + 1)
168 nae(l) = ptable(zval)%e_conv(l)/nn
169 IF (mod(ptable(zval)%e_conv(l), nn) /= 0) nae(l) = nae(l) + 1
170 npe(l) = econf(l)/nn
171 IF (mod(econf(l), nn) /= 0) npe(l) = npe(l) + 1
172 END DO
173 cpassert(all(nae - npe >= 0))
174 nco = nae - npe
175 ! MAO count
176 nmao = 0
177 DO l = 0, 3
178 nmao = nmao + npe(l)*(2*l + 1)
179 END DO
180
181 IF (mao > nmao) THEN
182 nmao = mao - nmao
183 SELECT CASE (nmao)
184 CASE (1)
185 npe(0) = npe(0) + 1
186 CASE (2)
187 npe(0) = npe(0) + 2
188 CASE (3)
189 npe(1) = npe(1) + 1
190 CASE (4)
191 npe(0) = npe(0) + 1
192 npe(1) = npe(1) + 1
193 CASE (5)
194 IF (npe(2) == 0) THEN
195 npe(2) = npe(2) + 1
196 ELSE
197 npe(0) = npe(0) + 2
198 npe(1) = npe(1) + 1
199 END IF
200 CASE (6)
201 npe(0) = npe(0) + 1
202 npe(2) = npe(2) + 1
203 CASE (7)
204 npe(0) = npe(0) + 2
205 npe(2) = npe(2) + 1
206 CASE (8)
207 npe(0) = npe(0) + 3
208 npe(2) = npe(2) + 1
209 CASE (9)
210 npe(0) = npe(0) + 1
211 npe(1) = npe(1) + 1
212 npe(2) = npe(2) + 1
213 CASE DEFAULT
214 cpabort("Default setting of minimal basis failed")
215 END SELECT
216 CALL cp_warn(__location__, "Reference basis has been adjusted according to MAO value.")
217 END IF
218
219 ! All shells should have at least 1 function
220 lm = 0
221 DO l = 0, 3
222 IF (npe(l) > 0) lm = l
223 END DO
224 DO l = 0, lm
225 IF (npe(l) == 0) npe(l) = 1
226 END DO
227
228 nss = sum(npe)
229 ALLOCATE (sym(nss), lq(nss), nq(nss), zet(nss))
230 iss = 0
231 DO l = 0, 3
232 DO i = 1, npe(l)
233 iss = iss + 1
234 lq(iss) = l
235 n = nco(l) + l
236 nq(iss) = n + i
237 str = " "
238 WRITE (str(5:5), fmt='(I1)') nq(iss)
239 str(6:6) = lnam(l)
240 sym(iss) = str
241 zet(iss) = srules(zval, ne, nq(iss), lq(iss))
242 END DO
243 END DO
244
245 bname = adjustr(ptable(zval)%symbol)//"_MBas"
246 CALL allocate_sto_basis_set(sto_basis_set)
247 CALL set_sto_basis_set(sto_basis_set, name=bname, nshell=nss, nq=nq, &
248 lq=lq, zet=zet, symbol=sym)
249 CALL create_gto_from_sto_basis(sto_basis_set, min_basis_set, ngau)
250 min_basis_set%norm_type = 2
252 CALL deallocate_sto_basis_set(sto_basis_set)
253
254 DEALLOCATE (sym, lq, nq, zet)
255
256 END SUBROUTINE create_min_basis
257
258! **************************************************************************************************
259
260END MODULE min_basis_set
subroutine, public add_basis_set_to_container(container, basis_set, basis_set_type)
...
subroutine, public create_primitive_basis_set(basis_set, pbasis)
...
subroutine, public deallocate_gto_basis_set(gto_basis_set)
...
pure real(dp) function, public srules(z, ne, n, l)
...
subroutine, public deallocate_sto_basis_set(sto_basis_set)
...
subroutine, public allocate_sto_basis_set(sto_basis_set)
...
subroutine, public create_gto_from_sto_basis(sto_basis_set, gto_basis_set, ngauss, ortho)
...
subroutine, public set_sto_basis_set(sto_basis_set, name, nshell, symbol, nq, lq, zet)
...
subroutine, public init_orb_basis_set(gto_basis_set)
Initialise a Gaussian-type orbital (GTO) basis set data set.
subroutine, public get_gto_basis_set(gto_basis_set, name, aliases, norm_type, kind_radius, ncgf, nset, nsgf, cgf_symbol, sgf_symbol, norm_cgf, set_radius, lmax, lmin, lx, ly, lz, m, ncgf_set, npgf, nsgf_set, nshell, cphi, pgf_radius, sphi, scon, zet, first_cgf, first_sgf, l, last_cgf, last_sgf, n, gcc, maxco, maxl, maxpgf, maxsgf_set, maxshell, maxso, nco_sum, npgf_sum, nshell_sum, maxder, short_kind_radius)
...
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
integer, parameter, public default_string_length
Definition kinds.F:57
generate or use from input minimal basis set
subroutine, public create_minbas_set(qs_env, unit_nr, basis_type, primitive)
...
Periodic Table related data definitions.
type(atom), dimension(0:nelem), public ptable
subroutine, public get_ptable_info(symbol, number, amass, ielement, covalent_radius, metallic_radius, vdw_radius, found)
Pass information about the kind given the element symbol.
subroutine, public get_qs_env(qs_env, atomic_kind_set, qs_kind_set, cell, super_cell, cell_ref, use_ref_cell, kpoints, dft_control, mos, sab_orb, sab_all, qmmm, qmmm_periodic, sac_ae, sac_ppl, sac_lri, sap_ppnl, sab_vdw, sab_scp, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, ecoul_1c, rho0_s_rs, rho0_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, wanniercentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, rhs)
Get the QUICKSTEP environment.
Calculate the interaction radii for the operator matrix calculation.
subroutine, public init_interaction_radii_orb_basis(orb_basis_set, eps_pgf_orb, eps_pgf_short)
...
Define the quickstep kind type and their sub types.
subroutine, public get_qs_kind(qs_kind, basis_set, basis_type, ncgf, nsgf, all_potential, tnadd_potential, gth_potential, sgp_potential, upf_potential, se_parameter, dftb_parameter, xtb_parameter, dftb3_param, zeff, elec_conf, mao, lmax_dftb, alpha_core_charge, ccore_charge, core_charge, core_charge_radius, paw_proj_set, paw_atom, hard_radius, hard0_radius, max_rad_local, covalent_radius, vdw_radius, gpw_r3d_rs_type_forced, harmonics, max_iso_not0, max_s_harm, grid_atom, ngrid_ang, ngrid_rad, lmax_rho0, dft_plus_u_atom, l_of_dft_plus_u, n_of_dft_plus_u, u_minus_j, u_of_dft_plus_u, j_of_dft_plus_u, alpha_of_dft_plus_u, beta_of_dft_plus_u, j0_of_dft_plus_u, occupation_of_dft_plus_u, dispersion, bs_occupation, magnetization, no_optimize, addel, laddel, naddel, orbitals, max_scf, eps_scf, smear, u_ramping, u_minus_j_target, eps_u_ramping, init_u_ramping_each_scf, reltmat, ghost, floating, name, element_symbol, pao_basis_size, pao_potentials, pao_descriptors, nelec)
Get attributes of an atomic kind.
Provides all information about a quickstep kind.