51#include "./base/base_uses.f90"
57 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'qmmm_se_forces'
58 LOGICAL,
PARAMETER,
PRIVATE :: debug_this_module = .false.
77 calc_force, Forces, Forces_added_charges)
84 LOGICAL,
INTENT(in),
OPTIONAL :: calc_force
85 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: forces, forces_added_charges
87 CHARACTER(len=*),
PARAMETER :: routinen =
'deriv_se_qmmm_matrix'
89 INTEGER :: handle, i, iatom, ikind, iqm, ispin, &
90 itype, natom, natorb_a, nkind, &
92 INTEGER,
DIMENSION(:),
POINTER ::
list
93 LOGICAL :: anag, defined, found
94 REAL(kind=
dp) :: delta, enuclear
95 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: forces_qm, p_block_a
100 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
107 CALL timeset(routinen, handle)
109 NULLIFY (rho, atomic_kind_set, qs_kind_set, se_taper)
110 NULLIFY (se_kind_a, se_kind_mm, particles_qm)
114 atomic_kind_set=atomic_kind_set, &
115 qs_kind_set=qs_kind_set, &
116 ks_qmmm_env=ks_qmmm_env_loc, &
117 dft_control=dft_control, &
118 particle_set=particles_qm, &
119 natom=number_qm_atoms)
120 SELECT CASE (dft_control%qs_control%method_id)
126 cpabort(
"Method not available")
128 anag = dft_control%qs_control%se_control%analytical_gradients
129 delta = dft_control%qs_control%se_control%delta
132 se_int_control, shortrange=.false., do_ewald_r3=.false., &
133 do_ewald_gks=.false., integral_screening=dft_control%qs_control%se_control%integral_screening, &
137 ALLOCATE (forces_qm(3, number_qm_atoms))
141 nkind =
SIZE(atomic_kind_set)
146 DO ispin = 1, dft_control%nspins
148 kinds:
DO ikind = 1, nkind
150 CALL get_qs_kind(qs_kind_set(ikind), se_parameter=se_kind_a)
154 IF (.NOT. defined .OR. natorb_a < 1) cycle
155 atoms:
DO i = 1,
SIZE(
list)
161 row=iatom, col=iatom, block=p_block_a, found=found)
163 IF (
ASSOCIATED(p_block_a))
THEN
165 CALL deriv_se_qmmm_matrix_low(p_block_a, &
168 qmmm_env%Potentials, &
170 qmmm_env%mm_atom_chrg, &
171 qmmm_env%mm_atom_index, &
181 qmmm_env%spherical_cutoff, &
184 IF (qmmm_env%move_mm_charges .OR. qmmm_env%add_mm_charges)
THEN
185 CALL deriv_se_qmmm_matrix_low(p_block_a, &
188 qmmm_env%added_charges%potentials, &
189 qmmm_env%added_charges%added_particles, &
190 qmmm_env%added_charges%mm_atom_chrg, &
191 qmmm_env%added_charges%mm_atom_index, &
195 forces_added_charges, &
201 qmmm_env%spherical_cutoff, &
208 cpassert(iqm == number_qm_atoms)
210 CALL para_env%sum(forces_qm)
214 CALL get_qs_kind(qs_kind_set(ikind), se_parameter=se_kind_a)
218 IF (.NOT. defined .OR. natorb_a < 1) cycle
221 iatom = qmmm_env%qm_atom_index(
list(i))
222 particles_mm(iatom)%f(:) = particles_mm(iatom)%f(:) + forces_qm(:, iqm)
227 DEALLOCATE (forces_qm)
231 CALL timestop(handle)
257 SUBROUTINE deriv_se_qmmm_matrix_low(p_block_a, se_kind_a, se_kind_mm, &
258 potentials, particles_mm, mm_charges, mm_atom_index, &
259 mm_cell, IndQM, itype, forces, forces_qm, se_taper, &
260 se_int_control, anag, delta, qmmm_spherical_cutoff, particles_qm)
262 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: p_block_a
266 REAL(kind=
dp),
DIMENSION(:),
POINTER :: mm_charges
267 INTEGER,
DIMENSION(:),
POINTER :: mm_atom_index
269 INTEGER,
INTENT(IN) :: indqm, itype
270 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: forces
271 REAL(kind=
dp),
DIMENSION(:),
INTENT(INOUT) :: forces_qm
274 LOGICAL,
INTENT(IN) :: anag
275 REAL(kind=
dp),
INTENT(IN) :: delta, qmmm_spherical_cutoff(2)
278 CHARACTER(len=*),
PARAMETER :: routinen =
'deriv_se_qmmm_matrix_low'
280 INTEGER :: handle, i1, i1l, i2, imm, imp, indmm, &
282 REAL(kind=
dp) :: rt1, rt2, rt3, sph_chrg_factor
283 REAL(kind=
dp),
DIMENSION(3) :: denuc, force_ab, r_pbc, rij
284 REAL(kind=
dp),
DIMENSION(3, 45) :: de1b
287 CALL timeset(routinen, handle)
290 mainlooppot:
DO ipot = 1,
SIZE(potentials)
291 pot => potentials(ipot)%Pot
293 loopmm:
DO imp = 1,
SIZE(pot%mm_atom_index)
294 imm = pot%mm_atom_index(imp)
295 indmm = mm_atom_index(imm)
296 r_pbc =
pbc(particles_mm(indmm)%r - particles_qm(indqm)%r, mm_cell)
300 rij = (/rt1, rt2, rt3/)
301 se_kind_mm%zeff = mm_charges(imm)
303 IF (qmmm_spherical_cutoff(1) > 0.0_dp)
THEN
305 se_kind_mm%zeff = se_kind_mm%zeff*sph_chrg_factor
307 IF (abs(se_kind_mm%zeff) <= epsilon(0.0_dp)) cycle
309 CALL drotnuc(se_kind_a, se_kind_mm, rij, itype=itype, de1b=de1b, &
310 se_int_control=se_int_control, anag=anag, delta=delta, &
312 CALL dcorecore(se_kind_a, se_kind_mm, rij, itype=itype, denuc=denuc, &
313 se_int_control=se_int_control, anag=anag, delta=delta, &
316 force_ab(1:3) = -denuc(1:3)
319 DO i1l = 1, se_kind_a%natorb
324 force_ab = force_ab - 2.0_dp*de1b(:, i2)*p_block_a(i1, j1)
328 force_ab = force_ab - de1b(:, i2)*p_block_a(i1, j1)
331 forces_qm(:) = forces_qm(:) - force_ab
333 forces(:, imm) = forces(:, imm) - force_ab
336 CALL timestop(handle)
337 END SUBROUTINE deriv_se_qmmm_matrix_low
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind(atomic_kind, fist_potential, element_symbol, name, mass, kind_number, natom, atom_list, rcov, rvdw, z, qeff, apol, cpol, mm_radius, shell, shell_active, damping)
Get attributes of an atomic kind.
Handles all functions related to the CELL.
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
subroutine, public dbcsr_get_block_p(matrix, row, col, block, found, row_size, col_size)
...
Defines the basic variable types.
integer, parameter, public dp
An array-based list which grows on demand. When the internal array is full, a new array of twice the ...
Interface to the message passing library MPI.
Multipole structure: for multipole (fixed and induced) in FF based MD.
integer, parameter, public do_multipole_none
Define the data structure for the particle information.
Calculation of the derivative of the QMMM Hamiltonian integral matrix <a|\sum_i q_i|b> for semi-empir...
subroutine, public deriv_se_qmmm_matrix(qs_env, qmmm_env, particles_mm, mm_cell, para_env, calc_force, forces, forces_added_charges)
Constructs the derivative w.r.t. 1-el semi-empirical hamiltonian QMMM terms.
subroutine, public spherical_cutoff_factor(spherical_cutoff, rij, factor)
Computes a spherical cutoff factor for the QMMM interactions.
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_pp, 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, harris_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, eeq, rhs)
Get the QUICKSTEP environment.
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, zatom, 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_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_model_file, pao_potentials, pao_descriptors, nelec)
Get attributes of an atomic kind.
superstucture that hold various representations of the density and keeps track of which ones are vali...
subroutine, public qs_rho_get(rho_struct, rho_ao, rho_ao_im, rho_ao_kp, rho_ao_im_kp, rho_r, drho_r, rho_g, drho_g, tau_r, tau_g, rho_r_valid, drho_r_valid, rho_g_valid, drho_g_valid, tau_r_valid, tau_g_valid, tot_rho_r, tot_rho_g, rho_r_sccs, soft_valid, complex_rho_ao)
returns info about the density described by this object. If some representation is not available an e...
Arrays of parameters used in the semi-empirical calculations \References Everywhere in this module TC...
integer, dimension(9), public se_orbital_pointer
Set of wrappers for semi-empirical analytical/numerical Integrals routines.
subroutine, public dcorecore(sepi, sepj, rij, denuc, itype, delta, anag, se_int_control, se_taper)
wrapper for numerical/analytical routines
subroutine, public drotnuc(sepi, sepj, rij, de1b, de2a, itype, delta, anag, se_int_control, se_taper)
wrapper for numerical/analytical routines
Definition of the semi empirical parameter types.
subroutine, public semi_empirical_create(sep)
Allocate semi-empirical type.
subroutine, public setup_se_int_control_type(se_int_control, shortrange, do_ewald_r3, do_ewald_gks, integral_screening, max_multipole, pc_coulomb_int)
Setup the Semiempirical integral control type.
subroutine, public get_se_param(sep, name, typ, defined, z, zeff, natorb, eheat, beta, sto_exponents, uss, upp, udd, uff, alp, eisol, gss, gsp, gpp, gp2, acoul, nr, de, ass, asp, app, hsp, gsd, gpd, gdd, ppddg, dpddg, ngauss)
Get info from the semi-empirical type.
subroutine, public semi_empirical_release(sep)
Deallocate the semi-empirical type.
Working with the semi empirical parameter types.
integer function, public get_se_type(se_method)
Gives back the unique semi_empirical METHOD type.
subroutine, public se_param_set_default(sep, z, method)
Initialize parameter for a semi_empirival type.
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
stores all the informations relevant to an mpi environment
Provides all information about a quickstep kind.
calculation environment to calculate the ks_qmmm matrix, holds the QM/MM potential and all the needed...
keeps the density in various representations, keeping track of which ones are valid.
Taper type use in semi-empirical calculations.