72#include "./base/base_uses.f90"
77 LOGICAL,
PRIVATE,
PARAMETER :: debug_this_module = .true.
78 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'qmmm_init'
112 mm_el_pot_radius_corr, mm_atom_index, mm_link_atoms, &
113 mm_link_scale_factor, added_shells, shell_model)
115 REAL(kind=
dp),
DIMENSION(:),
POINTER :: charges
116 REAL(
dp),
DIMENSION(:),
POINTER :: mm_atom_chrg, mm_el_pot_radius, &
117 mm_el_pot_radius_corr
118 INTEGER,
DIMENSION(:),
POINTER :: mm_atom_index, mm_link_atoms
119 REAL(
dp),
DIMENSION(:),
POINTER :: mm_link_scale_factor
121 LOGICAL :: shell_model
123 INTEGER :: i, ilink, indmm, indshell, ishell
125 REAL(
dp) :: qcore, qi, qshell, rc, ri
130 TYPE(
particle_type),
DIMENSION(:),
POINTER :: core_set, particle_set, shell_set
133 NULLIFY (particle_set, my_kind, added_shells)
134 CALL cp_subsys_get(subsys=subsys, particles=particles, core_particles=core_particles, &
135 shell_particles=shell_particles)
136 particle_set => particles%els
138 IF (all(particle_set(:)%shell_index == 0))
THEN
139 shell_model = .false.
145 IF (shell_model)
THEN
146 shell_set => shell_particles%els
147 core_set => core_particles%els
148 ishell =
SIZE(shell_set)
150 added_shells%added_particles => shell_set
151 added_shells%added_cores => core_set
154 DO i = 1,
SIZE(mm_atom_index)
155 indmm = mm_atom_index(i)
156 my_kind => particle_set(indmm)%atomic_kind
157 CALL get_atomic_kind(atomic_kind=my_kind, fist_potential=my_potential, &
158 shell_active=is_shell, shell=shell_kind)
163 IF (
ASSOCIATED(charges)) qi = charges(indmm)
165 mm_el_pot_radius(i) = ri
166 mm_el_pot_radius_corr(i) = rc
168 indshell = particle_set(indmm)%shell_index
169 IF (
ASSOCIATED(shell_kind))
THEN
173 mm_atom_chrg(i) = qcore
175 added_shells%mm_core_index(indshell) = indmm
176 added_shells%mm_core_chrg(indshell) = qshell
177 added_shells%mm_el_pot_radius(indshell) = ri*1.0_dp
178 added_shells%mm_el_pot_radius_corr(indshell) = rc*1.0_dp
182 IF (
ASSOCIATED(mm_link_atoms))
THEN
183 DO ilink = 1,
SIZE(mm_link_atoms)
184 DO i = 1,
SIZE(mm_atom_index)
185 IF (mm_atom_index(i) == mm_link_atoms(ilink))
EXIT
187 indmm = mm_atom_index(i)
188 mm_atom_chrg(i) = mm_atom_chrg(i)*mm_link_scale_factor(ilink)
209 SUBROUTINE print_qmmm_charges(mm_atom_index, mm_atom_chrg, mm_el_pot_radius, mm_el_pot_radius_corr, &
210 added_charges, added_shells, qmmm_section, nocompatibility, shell_model)
211 INTEGER,
DIMENSION(:),
POINTER :: mm_atom_index
212 REAL(
dp),
DIMENSION(:),
POINTER :: mm_atom_chrg, mm_el_pot_radius, &
213 mm_el_pot_radius_corr
217 LOGICAL,
INTENT(IN) :: nocompatibility, shell_model
219 INTEGER :: i, ind1, ind2, indmm, iw
220 REAL(kind=
dp) :: qi, qtot, rc, ri
228 WRITE (iw, fmt=
"(/,T2,A)") repeat(
"-", 79)
229 WRITE (iw, fmt=
'(/5X,A)')
"MM POINT CHARGES GENERATING THE QM/MM ELECTROSTATIC POTENTIAL"
230 WRITE (iw, fmt=
"(/,T2,A)") repeat(
"-", 79)
231 DO i = 1,
SIZE(mm_atom_index)
232 indmm = mm_atom_index(i)
235 ri = mm_el_pot_radius(i)
236 rc = mm_el_pot_radius_corr(i)
237 IF (nocompatibility)
THEN
238 WRITE (iw,
'(5X,A9,T15,I5,T28,A8,T38,F12.6,T60,A8,T69,F12.6)')
' MM ATOM:', indmm,
' RADIUS:', ri, &
241 WRITE (iw,
'(5X,A9,T15,I5,T28,A8,T38,F12.6,T60,A8,T69,F12.6,/,T56,A12,T69,F12.6)') &
242 ' MM ATOM:', indmm,
' RADIUS:', ri,
' CHARGE:', qi,
'CORR. RADIUS', rc
245 IF (added_charges%num_mm_atoms /= 0)
THEN
246 WRITE (iw, fmt=
"(/,T2,A)") repeat(
"-", 79)
247 WRITE (iw,
'(/5X,A)')
"ADDED POINT CHARGES GENERATING THE QM/MM ELECTROSTATIC POTENTIAL"
248 WRITE (iw, fmt=
"(/,T2,A)") repeat(
"-", 79)
249 DO i = 1,
SIZE(added_charges%mm_atom_index)
250 indmm = added_charges%mm_atom_index(i)
251 qi = added_charges%mm_atom_chrg(i)
253 ri = added_charges%mm_el_pot_radius(i)
254 ind1 = added_charges%add_env(i)%Index1
255 ind2 = added_charges%add_env(i)%Index2
256 IF (nocompatibility)
THEN
257 WRITE (iw,
'(5X,A9,I5,T25,A8,T35,F12.6,T50,A8,T59,F12.6,I5,I5)')
'MM POINT:', indmm,
' RADIUS:', ri, &
258 ' CHARGE:', qi, ind1, ind2
260 WRITE (iw,
'(5X,A9,I5,T25,A8,T35,F12.6,T50,A8,T59,F12.6,I5,I5,/,T56,A12,T69,F12.6)') &
261 'MM POINT:', indmm,
' RADIUS:', ri,
' CHARGE:', qi, ind1, ind2,
'CORR. RADIUS', rc
267 IF (shell_model)
THEN
268 WRITE (iw, fmt=
"(/,T2,A)") repeat(
"-", 73)
269 WRITE (iw,
'(/5X,A)')
"ADDED SHELL CHARGES GENERATING THE QM/MM ELECTROSTATIC POTENTIAL"
270 WRITE (iw, fmt=
"(/,T2,A)") repeat(
"-", 73)
272 DO i = 1,
SIZE(added_shells%mm_core_index)
273 indmm = added_shells%mm_core_index(i)
274 qi = added_shells%mm_core_chrg(i)
276 ri = added_shells%mm_el_pot_radius(i)
277 IF (nocompatibility)
THEN
278 WRITE (iw,
'(7X,A,I5,A8,F12.6,A8,F12.6,3F12.6)')
'SHELL:', indmm,
' RADIUS:', ri, &
279 ' CHARGE:', qi, added_shells%added_particles(i)%r
281 WRITE (iw,
'(7X,A,I5,A8,F12.6,A8,F12.6,A,F12.6)')
'SHELL:', indmm,
' RADIUS:', ri, &
282 ' CHARGE:', qi,
' CORR. RADIUS', rc
289 WRITE (iw, fmt=
"(/,T2,A)") repeat(
"-", 79)
290 WRITE (iw,
'(/,T50,A,T69,F12.6)')
' TOTAL CHARGE:', qtot
291 WRITE (iw, fmt=
"(/,T2,A,/)") repeat(
"-", 79)
294 "PRINT%QMMM_CHARGES")
309 INTEGER :: i, iw, mm_index, qm_index
310 REAL(kind=
dp) :: alpha
316 IF (
ASSOCIATED(qmmm_links))
THEN
317 WRITE (iw, fmt=
"(/,T2, A)") repeat(
"-", 73)
318 WRITE (iw, fmt=
"(/,T31,A)")
" QM/MM LINKS "
319 WRITE (iw, fmt=
"(/,T2, A)") repeat(
"-", 73)
320 IF (
ASSOCIATED(qmmm_links%imomm))
THEN
321 WRITE (iw, fmt=
"(/,T31,A)")
" IMOMM TYPE LINK "
322 DO i = 1,
SIZE(qmmm_links%imomm)
323 qm_index = qmmm_links%imomm(i)%link%qm_index
324 mm_index = qmmm_links%imomm(i)%link%mm_index
325 alpha = qmmm_links%imomm(i)%link%alpha
326 WRITE (iw, fmt=
"(T2,A,T20,A9,I8,1X,A9,I8,T62,A6,F12.6)")
"TYPE: IMOMM", &
327 "QM INDEX:", qm_index,
"MM INDEX:", mm_index,
"ALPHA:", alpha
330 IF (
ASSOCIATED(qmmm_links%pseudo))
THEN
331 WRITE (iw, fmt=
"(/,T31,A)")
" PSEUDO TYPE LINK "
332 DO i = 1,
SIZE(qmmm_links%pseudo)
333 qm_index = qmmm_links%pseudo(i)%link%qm_index
334 mm_index = qmmm_links%pseudo(i)%link%mm_index
335 WRITE (iw, fmt=
"(T2,A,T20,A9,I8,1X,A9,I8)")
"TYPE: PSEUDO", &
336 "QM INDEX:", qm_index,
"MM INDEX:", mm_index
339 WRITE (iw, fmt=
"(/,T2,A,/)") repeat(
"-", 73)
341 WRITE (iw, fmt=
"(/,T2, A)") repeat(
"-", 73)
342 WRITE (iw, fmt=
"(/,T26,A)")
" NO QM/MM LINKS DETECTED"
343 WRITE (iw, fmt=
"(/,T2, A)") repeat(
"-", 73)
347 "PRINT%qmmm_link_info")
365 mm_atom_chrg, qs_env, added_charges, added_shells, &
366 print_section, qmmm_section)
369 REAL(kind=
dp),
DIMENSION(:),
POINTER :: mm_atom_chrg
376 REAL(kind=
dp) :: maxchrg
377 REAL(kind=
dp),
DIMENSION(:),
POINTER :: maxradius, maxradius2
380 NULLIFY (maxradius, maxradius2, pw_env)
382 maxchrg = maxval(abs(mm_atom_chrg(:)))
384 IF (qmmm_env_qm%add_mm_charges) maxchrg = max(maxchrg, maxval(abs(added_charges%mm_atom_chrg(:))))
388 mm_el_pot_radius=qmmm_env_qm%mm_el_pot_radius, &
389 mm_el_pot_radius_corr=qmmm_env_qm%mm_el_pot_radius_corr, &
390 qmmm_coupl_type=qmmm_env_qm%qmmm_coupl_type, &
391 eps_mm_rspace=qmmm_env_qm%eps_mm_rspace, &
392 maxradius=maxradius, &
394 compatibility=qmmm_env_qm%compatibility, &
395 print_section=print_section, &
396 qmmm_section=qmmm_section)
398 IF (qmmm_env_qm%move_mm_charges .OR. qmmm_env_qm%add_mm_charges)
THEN
402 mm_el_pot_radius=added_charges%mm_el_pot_radius, &
403 mm_el_pot_radius_corr=added_charges%mm_el_pot_radius_corr, &
404 qmmm_coupl_type=qmmm_env_qm%qmmm_coupl_type, &
405 eps_mm_rspace=qmmm_env_qm%eps_mm_rspace, &
406 maxradius=maxradius2, &
408 compatibility=qmmm_env_qm%compatibility, &
409 print_section=print_section, &
410 qmmm_section=qmmm_section)
412 SELECT CASE (qmmm_env_qm%qmmm_coupl_type)
414 DO i = 1,
SIZE(maxradius)
415 maxradius(i) = max(maxradius(i), maxradius2(i))
419 IF (
ASSOCIATED(maxradius2))
DEALLOCATE (maxradius2)
422 IF (qmmm_env_qm%added_shells%num_mm_atoms > 0)
THEN
424 maxchrg = maxval(abs(added_shells%mm_core_chrg(:)))
428 mm_el_pot_radius=added_shells%mm_el_pot_radius, &
429 mm_el_pot_radius_corr=added_shells%mm_el_pot_radius_corr, &
430 qmmm_coupl_type=qmmm_env_qm%qmmm_coupl_type, &
431 eps_mm_rspace=qmmm_env_qm%eps_mm_rspace, &
432 maxradius=maxradius2, &
434 compatibility=qmmm_env_qm%compatibility, &
435 print_section=print_section, &
436 qmmm_section=qmmm_section)
438 SELECT CASE (qmmm_env_qm%qmmm_coupl_type)
440 DO i = 1,
SIZE(maxradius)
441 maxradius(i) = max(maxradius(i), maxradius2(i))
445 IF (
ASSOCIATED(maxradius2))
DEALLOCATE (maxradius2)
449 qmmm_env_qm%maxradius => maxradius
465 added_charges, added_shells, print_section)
473 mm_el_pot_radius=qmmm_env_qm%mm_el_pot_radius, &
474 potentials=qmmm_env_qm%potentials, &
475 pgfs=qmmm_env_qm%pgfs, &
477 compatibility=qmmm_env_qm%compatibility, &
478 print_section=print_section)
480 IF (qmmm_env_qm%move_mm_charges .OR. qmmm_env_qm%add_mm_charges)
THEN
483 mm_el_pot_radius=added_charges%mm_el_pot_radius, &
484 potentials=added_charges%potentials, &
485 pgfs=added_charges%pgfs, &
487 compatibility=qmmm_env_qm%compatibility, &
488 print_section=print_section)
491 IF (qmmm_env_qm%added_shells%num_mm_atoms > 0)
THEN
494 mm_el_pot_radius=added_shells%mm_el_pot_radius, &
495 potentials=added_shells%potentials, &
496 pgfs=added_shells%pgfs, &
498 compatibility=qmmm_env_qm%compatibility, &
499 print_section=print_section)
521 added_charges, added_shells, qmmm_periodic, print_section, mm_atom_chrg)
523 TYPE(
cell_type),
POINTER :: qm_cell_small, mm_cell
529 REAL(kind=
dp),
DIMENSION(:),
POINTER :: mm_atom_chrg
531 REAL(kind=
dp) :: maxchrg
534 IF (qmmm_env_qm%periodic)
THEN
536 NULLIFY (dft_control)
537 CALL get_qs_env(qs_env, dft_control=dft_control)
539 IF (dft_control%qs_control%semi_empirical)
THEN
540 cpabort(
"QM/MM periodic calculations not implemented for semi empirical methods")
541 ELSE IF (dft_control%qs_control%dftb)
THEN
543 qmmm_coupl_type=qmmm_env_qm%qmmm_coupl_type, mm_cell=mm_cell, &
544 para_env=para_env, qmmm_periodic=qmmm_periodic, print_section=print_section)
545 ELSE IF (dft_control%qs_control%xtb)
THEN
547 qmmm_coupl_type=qmmm_env_qm%qmmm_coupl_type, mm_cell=mm_cell, &
548 para_env=para_env, qmmm_periodic=qmmm_periodic, print_section=print_section)
552 maxchrg = maxval(abs(mm_atom_chrg(:)))
553 IF (qmmm_env_qm%add_mm_charges) maxchrg = max(maxchrg, maxval(abs(added_charges%mm_atom_chrg(:))))
556 per_potentials=qmmm_env_qm%per_potentials, &
557 potentials=qmmm_env_qm%potentials, &
558 pgfs=qmmm_env_qm%pgfs, &
559 qm_cell_small=qm_cell_small, &
561 compatibility=qmmm_env_qm%compatibility, &
562 qmmm_periodic=qmmm_periodic, &
563 print_section=print_section, &
564 eps_mm_rspace=qmmm_env_qm%eps_mm_rspace, &
566 ncp=qmmm_env_qm%aug_pools(
SIZE(qmmm_env_qm%aug_pools))%pool%pw_grid%npts, &
567 ncpl=qmmm_env_qm%aug_pools(
SIZE(qmmm_env_qm%aug_pools))%pool%pw_grid%npts_local)
569 IF (qmmm_env_qm%move_mm_charges .OR. qmmm_env_qm%add_mm_charges)
THEN
572 per_potentials=added_charges%per_potentials, &
573 potentials=added_charges%potentials, &
574 pgfs=added_charges%pgfs, &
575 qm_cell_small=qm_cell_small, &
577 compatibility=qmmm_env_qm%compatibility, &
578 qmmm_periodic=qmmm_periodic, &
579 print_section=print_section, &
580 eps_mm_rspace=qmmm_env_qm%eps_mm_rspace, &
582 ncp=qmmm_env_qm%aug_pools(
SIZE(qmmm_env_qm%aug_pools))%pool%pw_grid%npts, &
583 ncpl=qmmm_env_qm%aug_pools(
SIZE(qmmm_env_qm%aug_pools))%pool%pw_grid%npts_local)
586 IF (qmmm_env_qm%added_shells%num_mm_atoms > 0)
THEN
589 per_potentials=added_shells%per_potentials, &
590 potentials=added_shells%potentials, &
591 pgfs=added_shells%pgfs, &
592 qm_cell_small=qm_cell_small, &
594 compatibility=qmmm_env_qm%compatibility, &
595 qmmm_periodic=qmmm_periodic, &
596 print_section=print_section, &
597 eps_mm_rspace=qmmm_env_qm%eps_mm_rspace, &
599 ncp=qmmm_env_qm%aug_pools(
SIZE(qmmm_env_qm%aug_pools))%pool%pw_grid%npts, &
600 ncpl=qmmm_env_qm%aug_pools(
SIZE(qmmm_env_qm%aug_pools))%pool%pw_grid%npts_local)
627 qm_atom_index, mm_atom_index, qm_cell_small, qmmm_coupl_type, eps_mm_rspace, &
632 CHARACTER(len=default_string_length), &
633 DIMENSION(:),
POINTER :: qm_atom_type
634 INTEGER,
DIMENSION(:),
POINTER :: qm_atom_index, mm_atom_index
635 TYPE(
cell_type),
POINTER :: qm_cell_small
636 INTEGER,
INTENT(OUT) :: qmmm_coupl_type
637 REAL(kind=
dp),
INTENT(OUT) :: eps_mm_rspace
638 LOGICAL,
INTENT(OUT) :: qmmm_link
641 CHARACTER(len=default_string_length) :: atmname, mm_atom_kind
642 INTEGER :: i, icount, ikind, ikindr, my_type, &
643 n_rep_val, nkind, size_mm_system
644 INTEGER,
DIMENSION(:),
POINTER :: mm_link_atoms
645 LOGICAL :: explicit, is_mm, is_qm
646 REAL(kind=
dp) :: tmp_radius, tmp_radius_c
647 REAL(kind=
dp),
DIMENSION(:),
POINTER :: tmp_sph_cut
654 NULLIFY (mm_link_atoms, tmp_sph_cut)
655 NULLIFY (image_charge_section)
663 cpassert(
SIZE(tmp_sph_cut) == 2)
664 qmmm_env%spherical_cutoff = tmp_sph_cut
665 IF (qmmm_env%spherical_cutoff(1) <= 0.0_dp)
THEN
666 qmmm_env%spherical_cutoff(2) = 0.0_dp
668 IF (qmmm_env%spherical_cutoff(2) <= 0.0_dp) qmmm_env%spherical_cutoff(2) = epsilon(0.0_dp)
669 tmp_radius = qmmm_env%spherical_cutoff(1) - 20.0_dp*qmmm_env%spherical_cutoff(2)
670 IF (tmp_radius <= 0.0_dp) &
671 CALL cp_abort(__location__, &
672 "SPHERICAL_CUTOFF(1) > 20*SPHERICAL_CUTOFF(1)! Please correct parameters for "// &
673 "the Spherical Cutoff in order to satisfy the previous condition!")
679 CALL cp_subsys_get(subsys=subsys_mm, atomic_kinds=atomic_kinds)
680 DO ikind = 1,
SIZE(atomic_kinds%els)
681 atomic_kind => atomic_kinds%els(ikind)
683 fist_potential=fist_potential)
685 qmmm_radius=tmp_radius, &
686 qmmm_corr_radius=tmp_radius)
688 fist_potential=fist_potential)
690 CALL setup_qm_atom_list(qmmm_section=qmmm_section, &
691 qm_atom_index=qm_atom_index, &
692 qm_atom_type=qm_atom_type, &
693 mm_link_atoms=mm_link_atoms, &
704 set_radius_pot_0:
DO ikindr = 1,
SIZE(atomic_kinds%els)
705 atomic_kind => atomic_kinds%els(ikindr)
708 fist_potential=fist_potential)
709 CALL set_potential(potential=fist_potential, qmmm_radius=tmp_radius, &
710 qmmm_corr_radius=tmp_radius)
712 fist_potential=fist_potential)
713 END DO set_radius_pot_0
722 tmp_radius_c = tmp_radius
726 set_radius_pot_1:
DO ikindr = 1,
SIZE(atomic_kinds%els)
727 atomic_kind => atomic_kinds%els(ikindr)
730 IF (trim(mm_atom_kind) == atmname)
THEN
732 fist_potential=fist_potential)
734 qmmm_radius=tmp_radius, &
735 qmmm_corr_radius=tmp_radius_c)
737 fist_potential=fist_potential)
739 END DO set_radius_pot_1
749 cpabort(
"QMMM section not present in input file!")
754 size_mm_system =
SIZE(subsys_mm%particles%els) -
SIZE(qm_atom_index)
755 IF (qmmm_link .AND.
ASSOCIATED(mm_link_atoms)) size_mm_system = size_mm_system +
SIZE(mm_link_atoms)
756 ALLOCATE (mm_atom_index(size_mm_system))
759 DO i = 1,
SIZE(subsys_mm%particles%els)
761 IF (any(qm_atom_index == i))
THEN
764 IF (
ASSOCIATED(mm_link_atoms))
THEN
765 IF (any(mm_link_atoms == i) .AND. qmmm_link) is_mm = .true.
769 IF (icount <= size_mm_system) mm_atom_index(icount) = i
772 cpassert(icount == size_mm_system)
773 IF (
ASSOCIATED(mm_link_atoms))
THEN
774 DEALLOCATE (mm_link_atoms)
779 IF (qmmm_env%image_charge)
THEN
782 IF (explicit) qmmm_env%image_charge_pot%all_mm = .false.
784 IF (qmmm_env%image_charge_pot%all_mm)
THEN
785 qmmm_env%image_charge_pot%image_mm_list => mm_atom_index
787 CALL setup_image_atom_list(image_charge_section, qmmm_env, &
788 qm_atom_index, subsys_mm)
791 qmmm_env%image_charge_pot%particles_all => subsys_mm%particles%els
794 r_val=qmmm_env%image_charge_pot%V0)
796 r_val=qmmm_env%image_charge_pot%eta)
799 SELECT CASE (my_type)
801 qmmm_env%image_charge_pot%coeff_iterative = .false.
803 qmmm_env%image_charge_pot%coeff_iterative = .true.
807 l_val=qmmm_env%image_charge_pot%image_restart)
810 i_val=qmmm_env%image_charge_pot%image_matrix_method)
812 IF (qmmm_env%image_charge_pot%image_matrix_method ==
do_eri_mme)
THEN
816 hmat=qm_cell_small%hmat, is_ortho=qm_cell_small%orthorhombic, &
817 zet_min=qmmm_env%image_charge_pot%eta, &
818 zet_max=qmmm_env%image_charge_pot%eta, &
843 mm_link_atoms, mm_link_scale_factor, &
844 fist_scale_charge_link, qmmm_coupl_type, &
848 INTEGER,
DIMENSION(:),
POINTER :: qm_atom_index, mm_link_atoms
849 REAL(kind=
dp),
DIMENSION(:),
POINTER :: mm_link_scale_factor, &
850 fist_scale_charge_link
851 INTEGER,
INTENT(OUT) :: qmmm_coupl_type
852 LOGICAL,
INTENT(OUT) :: qmmm_link
857 NULLIFY (qmmm_ff_section)
862 CALL setup_qm_atom_list(qmmm_section, qm_atom_index=qm_atom_index, qmmm_link=qmmm_link, &
863 mm_link_atoms=mm_link_atoms, mm_link_scale_factor=mm_link_scale_factor, &
864 fist_scale_charge_link=fist_scale_charge_link)
870 IF (qmmm_env%use_qmmm_ff)
THEN
872 l_val=qmmm_env%multiple_potential)
873 CALL read_qmmm_ff_section(qmmm_ff_section, qmmm_env%inp_info)
887 SUBROUTINE read_qmmm_ff_section(qmmm_ff_section, inp_info)
891 INTEGER :: n_gd, n_gp, n_lj, n_wl, np
916 np = n_lj + n_wl + n_gd
919 CALL read_gd_section(inp_info%nonbonded, gd_section, start=n_lj + n_wl)
922 np = n_lj + n_wl + n_gd + n_gp
925 CALL read_gp_section(inp_info%nonbonded, gp_section, start=n_lj + n_wl + n_gd)
947 np = n_lj + n_wl + n_gd
950 CALL read_gd_section(inp_info%nonbonded14, gd_section, start=n_lj + n_wl)
953 np = n_lj + n_wl + n_gd + n_gp
956 CALL read_gp_section(inp_info%nonbonded14, gp_section, start=n_lj + n_wl + n_gd)
959 END SUBROUTINE read_qmmm_ff_section
974 SUBROUTINE setup_qm_atom_list(qmmm_section, qm_atom_index, qm_atom_type, &
975 mm_link_atoms, mm_link_scale_factor, qmmm_link, fist_scale_charge_link)
977 INTEGER,
DIMENSION(:),
OPTIONAL,
POINTER :: qm_atom_index
978 CHARACTER(len=default_string_length), &
979 DIMENSION(:),
OPTIONAL,
POINTER :: qm_atom_type
980 INTEGER,
DIMENSION(:),
OPTIONAL,
POINTER :: mm_link_atoms
981 REAL(kind=
dp),
DIMENSION(:),
OPTIONAL,
POINTER :: mm_link_scale_factor
982 LOGICAL,
INTENT(OUT),
OPTIONAL :: qmmm_link
983 REAL(kind=
dp),
DIMENSION(:),
OPTIONAL,
POINTER :: fist_scale_charge_link
985 CHARACTER(len=default_string_length) :: qm_atom_kind, qm_link_element
986 INTEGER :: ikind, k, link_involv_mm, link_type, &
987 mm_index, n_var, nkind, nlinks, &
989 INTEGER,
DIMENSION(:),
POINTER :: mm_indexes
991 REAL(kind=
dp) :: scale_f
1007 num_qm_atom_tot = num_qm_atom_tot +
SIZE(mm_indexes)
1019 DO ikind = 1, nlinks
1022 SELECT CASE (link_type)
1024 num_qm_atom_tot = num_qm_atom_tot + 1
1025 link_involv_mm = link_involv_mm + 1
1027 num_qm_atom_tot = num_qm_atom_tot + 1
1035 IF (
PRESENT(mm_link_scale_factor) .AND. (link_involv_mm /= 0)) &
1036 ALLOCATE (mm_link_scale_factor(link_involv_mm))
1037 IF (
PRESENT(fist_scale_charge_link) .AND. (link_involv_mm /= 0)) &
1038 ALLOCATE (fist_scale_charge_link(link_involv_mm))
1039 IF (
PRESENT(mm_link_atoms) .AND. (link_involv_mm /= 0)) &
1040 ALLOCATE (mm_link_atoms(link_involv_mm))
1041 IF (
PRESENT(qm_atom_index))
ALLOCATE (qm_atom_index(num_qm_atom_tot))
1042 IF (
PRESENT(qm_atom_type))
ALLOCATE (qm_atom_type(num_qm_atom_tot))
1043 IF (
PRESENT(qm_atom_index)) qm_atom_index = 0
1044 IF (
PRESENT(qm_atom_type)) qm_atom_type =
" "
1051 IF (
PRESENT(qm_atom_index))
THEN
1052 qm_atom_index(num_qm_atom_tot:num_qm_atom_tot +
SIZE(mm_indexes) - 1) = mm_indexes(:)
1054 IF (
PRESENT(qm_atom_type))
THEN
1057 qm_atom_type(num_qm_atom_tot:num_qm_atom_tot +
SIZE(mm_indexes) - 1) = qm_atom_kind
1059 num_qm_atom_tot = num_qm_atom_tot +
SIZE(mm_indexes)
1062 IF (
PRESENT(mm_link_scale_factor) .AND. (link_involv_mm /= 0)) mm_link_scale_factor = 0.0_dp
1063 IF (
PRESENT(fist_scale_charge_link) .AND. (link_involv_mm /= 0)) fist_scale_charge_link = 0.0_dp
1064 IF (
PRESENT(mm_link_atoms) .AND. (link_involv_mm /= 0)) mm_link_atoms = 0
1066 DO ikind = 1, nlinks
1067 IF (
PRESENT(qm_atom_type))
THEN
1069 qm_atom_type(num_qm_atom_tot:num_qm_atom_tot) = trim(qm_link_element)//
"_LINK"
1071 IF (
PRESENT(qm_atom_index))
THEN
1073 cpassert(all(qm_atom_index /= mm_index))
1074 qm_atom_index(num_qm_atom_tot:num_qm_atom_tot) = mm_index
1075 num_qm_atom_tot = num_qm_atom_tot + 1
1077 IF (
PRESENT(mm_link_atoms) .AND. (link_involv_mm /= 0))
THEN
1079 mm_link_atoms(ikind) = mm_index
1081 IF (
PRESENT(mm_link_scale_factor) .AND. (link_involv_mm /= 0))
THEN
1083 mm_link_scale_factor(ikind) = scale_f
1085 IF (
PRESENT(fist_scale_charge_link) .AND. (link_involv_mm /= 0))
THEN
1087 fist_scale_charge_link(ikind) = scale_f
1091 cpassert(num_qm_atom_tot - 1 ==
SIZE(qm_atom_index))
1093 END SUBROUTINE setup_qm_atom_list
1107 SUBROUTINE setup_qmmm_links(qmmm_section, qmmm_links, mm_el_pot_radius, mm_el_pot_radius_corr, &
1111 REAL(kind=
dp),
DIMENSION(:),
POINTER :: mm_el_pot_radius, mm_el_pot_radius_corr
1112 INTEGER,
DIMENSION(:),
POINTER :: mm_atom_index
1113 INTEGER,
INTENT(IN) :: iw
1115 INTEGER :: ikind, link_type, mm_index, n_gho, &
1116 n_imomm, n_pseudo, n_rep_val, n_tot, &
1118 INTEGER,
DIMENSION(:),
POINTER :: wrk_tmp
1119 REAL(kind=
dp) :: alpha, my_radius
1128 cpassert(nlinks /= 0)
1129 DO ikind = 1, nlinks
1130 CALL section_vals_val_get(qmmm_link_section,
"LINK_TYPE", i_rep_section=ikind, i_val=link_type)
1135 n_tot = n_imomm + n_gho + n_pseudo
1136 cpassert(n_tot /= 0)
1137 ALLOCATE (qmmm_links)
1138 NULLIFY (qmmm_links%imomm, &
1141 IF (n_imomm /= 0)
THEN
1142 ALLOCATE (qmmm_links%imomm(n_imomm))
1143 ALLOCATE (wrk_tmp(n_imomm))
1144 DO ikind = 1, n_imomm
1145 NULLIFY (qmmm_links%imomm(ikind)%link)
1146 ALLOCATE (qmmm_links%imomm(ikind)%link)
1149 DO ikind = 1, nlinks
1150 CALL section_vals_val_get(qmmm_link_section,
"LINK_TYPE", i_rep_section=ikind, i_val=link_type)
1152 n_imomm = n_imomm + 1
1156 CALL section_vals_val_get(qmmm_link_section,
"RADIUS", i_rep_section=ikind, n_rep_val=n_rep_val)
1157 qmmm_links%imomm(n_imomm)%link%qm_index = qm_index
1158 qmmm_links%imomm(n_imomm)%link%mm_index = mm_index
1159 qmmm_links%imomm(n_imomm)%link%alpha = alpha
1160 wrk_tmp(n_imomm) = mm_index
1161 IF (n_rep_val == 1)
THEN
1163 WHERE (mm_atom_index == mm_index) mm_el_pot_radius = my_radius
1164 WHERE (mm_atom_index == mm_index) mm_el_pot_radius_corr = my_radius
1166 CALL section_vals_val_get(qmmm_link_section,
"CORR_RADIUS", i_rep_section=ikind, n_rep_val=n_rep_val)
1167 IF (n_rep_val == 1)
THEN
1168 CALL section_vals_val_get(qmmm_link_section,
"CORR_RADIUS", i_rep_section=ikind, r_val=my_radius)
1169 WHERE (mm_atom_index == mm_index) mm_el_pot_radius_corr = my_radius
1176 DO ikind = 1,
SIZE(wrk_tmp)
1177 IF (count(wrk_tmp == wrk_tmp(ikind)) > 1)
THEN
1178 WRITE (iw,
'(/A)')
"In the IMOMM scheme no more than one QM atom can be bounded to the same MM atom."
1179 WRITE (iw,
'(A)')
"Multiple link MM atom not allowed. Check your link sections."
1183 DEALLOCATE (wrk_tmp)
1186 IF (n_pseudo /= 0)
THEN
1187 ALLOCATE (qmmm_links%pseudo(n_pseudo))
1188 DO ikind = 1, n_pseudo
1189 NULLIFY (qmmm_links%pseudo(ikind)%link)
1190 ALLOCATE (qmmm_links%pseudo(ikind)%link)
1193 DO ikind = 1, nlinks
1194 CALL section_vals_val_get(qmmm_link_section,
"LINK_TYPE", i_rep_section=ikind, i_val=link_type)
1196 n_pseudo = n_pseudo + 1
1199 qmmm_links%pseudo(n_pseudo)%link%qm_index = qm_index
1200 qmmm_links%pseudo(n_pseudo)%link%mm_index = mm_index
1205 IF (n_gho /= 0)
THEN
1227 mm_atom_chrg, mm_el_pot_radius, mm_el_pot_radius_corr, &
1228 added_charges, mm_atom_index)
1230 LOGICAL,
INTENT(OUT) :: move_mm_charges, add_mm_charges
1231 REAL(kind=
dp),
DIMENSION(:),
POINTER :: mm_atom_chrg, mm_el_pot_radius, &
1232 mm_el_pot_radius_corr
1234 INTEGER,
DIMENSION(:),
POINTER :: mm_atom_index
1236 INTEGER :: i_add, icount, ikind, ind1, index1, &
1237 index2, n_add_tot, n_adds, n_move_tot, &
1238 n_moves, n_rep_val, nlinks
1240 REAL(kind=
dp) :: alpha, c_radius, charge, radius
1245 move_mm_charges = .false.
1246 add_mm_charges = .false.
1247 NULLIFY (qmmm_link_section, move_section, add_section)
1250 cpassert(nlinks /= 0)
1254 DO ikind = 1, nlinks
1256 i_rep_section=ikind)
1259 i_rep_section=ikind)
1261 n_move_tot = n_move_tot + n_moves
1262 n_add_tot = n_add_tot + n_adds
1264 icount = n_move_tot + n_add_tot
1265 IF (n_add_tot /= 0) add_mm_charges = .true.
1266 IF (n_move_tot /= 0) move_mm_charges = .true.
1275 DO ikind = 1, nlinks
1277 i_rep_section=ikind)
1278 CALL section_vals_get(move_section, explicit=explicit, n_repetition=n_moves)
1283 DO i_add = 1, n_moves
1289 CALL section_vals_val_get(move_section,
"CORR_RADIUS", n_rep_val=n_rep_val, i_rep_section=i_add)
1291 IF (n_rep_val == 1) &
1294 CALL set_add_set_type(added_charges, icount, index1, index2, alpha, radius, c_radius, &
1295 mm_atom_chrg=mm_atom_chrg, mm_el_pot_radius=mm_el_pot_radius, &
1296 mm_el_pot_radius_corr=mm_el_pot_radius_corr, &
1297 mm_atom_index=mm_atom_index, move=n_moves, ind1=ind1)
1299 mm_atom_chrg(ind1) = 0.0_dp
1303 i_rep_section=ikind)
1309 DO i_add = 1, n_adds
1316 CALL section_vals_val_get(add_section,
"CORR_RADIUS", n_rep_val=n_rep_val, i_rep_section=i_add)
1318 IF (n_rep_val == 1) &
1321 CALL set_add_set_type(added_charges, icount, index1, index2, alpha, radius, c_radius, charge, &
1322 mm_atom_chrg=mm_atom_chrg, mm_el_pot_radius=mm_el_pot_radius, &
1323 mm_el_pot_radius_corr=mm_el_pot_radius_corr, &
1324 mm_atom_index=mm_atom_index)
1351 SUBROUTINE set_add_set_type(added_charges, icount, Index1, Index2, alpha, radius, c_radius, charge, &
1352 mm_atom_chrg, mm_el_pot_radius, mm_el_pot_radius_corr, mm_atom_index, move, ind1)
1354 INTEGER,
INTENT(IN) :: icount, index1, index2
1355 REAL(kind=
dp),
INTENT(IN) :: alpha, radius, c_radius
1356 REAL(kind=
dp),
INTENT(IN),
OPTIONAL :: charge
1357 REAL(kind=
dp),
DIMENSION(:),
POINTER :: mm_atom_chrg, mm_el_pot_radius, &
1358 mm_el_pot_radius_corr
1359 INTEGER,
DIMENSION(:),
POINTER :: mm_atom_index
1360 INTEGER,
INTENT(in),
OPTIONAL :: move
1361 INTEGER,
INTENT(OUT),
OPTIONAL :: ind1
1363 INTEGER :: i, my_move
1364 REAL(kind=
dp) :: my_c_radius, my_charge, my_radius
1368 my_c_radius = c_radius
1369 IF (
PRESENT(charge)) my_charge = charge
1370 IF (
PRESENT(move)) my_move = move
1372 getid:
DO WHILE (i <=
SIZE(mm_atom_index))
1373 IF (index1 == mm_atom_index(i))
EXIT getid
1376 IF (
PRESENT(ind1)) ind1 = i
1377 cpassert(i <=
SIZE(mm_atom_index))
1378 IF (.NOT.
PRESENT(charge)) my_charge = mm_atom_chrg(i)/real(my_move, kind=
dp)
1379 IF (my_radius == 0.0_dp) my_radius = mm_el_pot_radius(i)
1380 IF (my_c_radius == 0.0_dp) my_c_radius = mm_el_pot_radius_corr(i)
1382 added_charges%add_env(icount)%Index1 = index1
1383 added_charges%add_env(icount)%Index2 = index2
1384 added_charges%add_env(icount)%alpha = alpha
1385 added_charges%mm_atom_index(icount) = icount
1386 added_charges%mm_atom_chrg(icount) = my_charge
1387 added_charges%mm_el_pot_radius(icount) = my_radius
1388 added_charges%mm_el_pot_radius_corr(icount) = my_c_radius
1389 END SUBROUTINE set_add_set_type
1407 TYPE(
cell_type),
POINTER :: qm_cell_small
1408 REAL(kind=
dp),
DIMENSION(3),
INTENT(in) :: dr
1410 LOGICAL :: center_grid
1411 REAL(kind=
dp),
DIMENSION(3) :: tmp
1412 REAL(kind=
dp),
DIMENSION(:),
POINTER :: vec
1416 tmp(1) = qm_cell_small%hmat(1, 1)
1417 tmp(2) = qm_cell_small%hmat(2, 2)
1418 tmp(3) = qm_cell_small%hmat(3, 3)
1419 cpassert(all(tmp > 0))
1420 qmmm_env%dOmmOqm = tmp/2.0_dp
1424 IF (center_grid)
THEN
1425 qmmm_env%utrasl = dr
1427 qmmm_env%utrasl = 1.0_dp
1430 qmmm_env%transl_v = vec
1443 SUBROUTINE setup_image_atom_list(image_charge_section, qmmm_env, &
1444 qm_atom_index, subsys_mm)
1448 INTEGER,
DIMENSION(:),
POINTER :: qm_atom_index
1451 INTEGER :: atom_a, atom_b, i, j, k, max_index, &
1452 n_var, num_const_atom, &
1454 INTEGER,
DIMENSION(:),
POINTER :: mm_indexes
1455 LOGICAL :: fix_xyz, imageind_in_range
1458 NULLIFY (mm_indexes, molecule_kind)
1459 imageind_in_range = .false.
1460 num_image_mm_atom = 0
1467 i_rep_val=i, i_vals=mm_indexes)
1468 num_image_mm_atom = num_image_mm_atom +
SIZE(mm_indexes)
1471 ALLOCATE (qmmm_env%image_charge_pot%image_mm_list(num_image_mm_atom))
1473 qmmm_env%image_charge_pot%image_mm_list = 0
1474 num_image_mm_atom = 1
1478 i_rep_val=i, i_vals=mm_indexes)
1479 qmmm_env%image_charge_pot%image_mm_list(num_image_mm_atom:num_image_mm_atom &
1480 +
SIZE(mm_indexes) - 1) = mm_indexes(:)
1481 num_image_mm_atom = num_image_mm_atom +
SIZE(mm_indexes)
1485 num_image_mm_atom = num_image_mm_atom - 1
1487 max_index =
SIZE(subsys_mm%particles%els)
1489 cpassert(
SIZE(qmmm_env%image_charge_pot%image_mm_list) /= 0)
1490 imageind_in_range = (maxval(qmmm_env%image_charge_pot%image_mm_list) <= max_index) &
1491 .AND. (minval(qmmm_env%image_charge_pot%image_mm_list) > 0)
1492 cpassert(imageind_in_range)
1494 DO i = 1, num_image_mm_atom
1495 atom_a = qmmm_env%image_charge_pot%image_mm_list(i)
1496 IF (any(qm_atom_index == atom_a))
THEN
1497 cpabort(
"Image atom list must only contain MM atoms")
1499 DO j = i + 1, num_image_mm_atom
1500 atom_b = qmmm_env%image_charge_pot%image_mm_list(j)
1501 IF (atom_a == atom_b) &
1502 cpabort(
"There are atoms doubled in image list.")
1509 IF (
ASSOCIATED(subsys_mm%molecule_kinds))
THEN
1510 IF (
ASSOCIATED(subsys_mm%molecule_kinds%els))
THEN
1511 molecule_kind => subsys_mm%molecule_kinds%els
1512 DO i = 1,
SIZE(molecule_kind)
1513 IF (.NOT.
ASSOCIATED(molecule_kind(i)%fixd_list))
EXIT
1514 IF (.NOT. fix_xyz)
EXIT
1515 DO j = 1,
SIZE(molecule_kind(i)%fixd_list)
1516 IF (.NOT. fix_xyz)
EXIT
1517 DO k = 1, num_image_mm_atom
1518 atom_a = qmmm_env%image_charge_pot%image_mm_list(k)
1519 IF (atom_a == molecule_kind(i)%fixd_list(j)%fixd)
THEN
1520 num_const_atom = num_const_atom + 1
1521 IF (molecule_kind(i)%fixd_list(j)%itype /=
use_perd_xyz)
THEN
1534 IF (num_const_atom == num_image_mm_atom .AND. fix_xyz)
THEN
1535 qmmm_env%image_charge_pot%state_image_matrix =
calc_once
1537 qmmm_env%image_charge_pot%state_image_matrix =
calc_always
1540 END SUBROUTINE setup_image_atom_list
1556 REAL(kind=
dp) :: eta, eta_conv, v0, v0_conv
1562 eta = qmmm_env%image_charge_pot%eta
1564 v0 = qmmm_env%image_charge_pot%V0
1568 WRITE (iw, fmt=
"(T25,A)")
"IMAGE CHARGE PARAMETERS"
1569 WRITE (iw, fmt=
"(T25,A)") repeat(
"-", 23)
1570 WRITE (iw, fmt=
"(/)")
1571 WRITE (iw, fmt=
"(T2,A)")
"INDEX OF MM ATOMS CARRYING AN IMAGE CHARGE:"
1572 WRITE (iw, fmt=
"(/)")
1574 WRITE (iw,
"(7X,10I6)") qmmm_env%image_charge_pot%image_mm_list
1575 WRITE (iw, fmt=
"(/)")
1576 WRITE (iw,
"(T2,A52,T69,F12.8)") &
1577 "WIDTH OF GAUSSIAN CHARGE DISTRIBUTION [angstrom^-2]:", eta_conv
1578 WRITE (iw,
"(T2,A26,T69,F12.8)")
"EXTERNAL POTENTIAL [volt]:", v0_conv
1579 WRITE (iw, fmt=
"(/,T2,A,/)") repeat(
"-", 79)
1582 "PRINT%PROGRAM_RUN_INFO")
represent a simple array based list of the given type
Define the atomic kind types and their sub types.
subroutine, public set_atomic_kind(atomic_kind, element_symbol, name, mass, kind_number, natom, atom_list, fist_potential, shell, shell_active, damping)
Set the components of an atomic kind data set.
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.
integer, parameter, public use_perd_xyz
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
Interface to Minimax-Ewald method for periodic ERI's to be used in CP2K.
subroutine, public cp_eri_mme_init_read_input(mme_section, param)
Read input and initialize parameter type.
various routines to log and control the output. The idea is that decisions about where to log should ...
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer function, public cp_print_key_unit_nr(logger, basis_section, print_key_path, extension, middle_name, local, log_filename, ignore_should_output, file_form, file_position, file_action, file_status, do_backup, on_file, is_new_file, mpi_io, fout)
...
subroutine, public cp_print_key_finished_output(unit_nr, logger, basis_section, print_key_path, local, ignore_should_output, on_file, mpi_io)
should be called after you finish working with a unit obtained with cp_print_key_unit_nr,...
types that represent a subsys, i.e. a part of the system
subroutine, public cp_subsys_get(subsys, ref_count, atomic_kinds, atomic_kind_set, particles, particle_set, local_particles, molecules, molecule_set, molecule_kinds, molecule_kind_set, local_molecules, para_env, colvar_p, shell_particles, core_particles, gci, multipoles, natom, nparticle, ncore, nshell, nkind, atprop, virial, results, cell, cell_ref, use_ref_cell)
returns information about various attributes of the given subsys
real(kind=dp) function, public cp_unit_from_cp2k(value, unit_str, defaults, power)
converts from the internal cp2k units to the given unit
real(kind=dp) function, public cp_unit_to_cp2k(value, unit_str, defaults, power)
converts to the internal cp2k units to the given unit
Definition of the atomic potential types.
Define all structures types related to force_fields.
Defines the basic variable types.
integer, parameter, public dp
integer, parameter, public default_string_length
Interface to the message passing library MPI.
Define the molecule kind structure types and the corresponding functionality.
subroutine, public pair_potential_reallocate(p, lb1_new, ub1_new, lj, lj_charmm, williams, goodwin, eam, nequip, bmhft, bmhftd, ipbv, buck4r, buckmo, gp, tersoff, siepmann, gal, gal21, tab, deepmd, ace)
Cleans the potential parameter type.
represent a simple array based list of the given type
Define the data structure for the particle information.
container for various plainwaves related things
subroutine, public qmmm_potential_init(qmmm_coupl_type, mm_el_pot_radius, potentials, pgfs, mm_cell, compatibility, print_section)
Initialize the QMMM potential stored on vector, according the qmmm_coupl_type.
logical function, public qmmm_ff_precond_only_qm(id1, id2, id3, id4, is_link)
This function handles the atom names and modifies the "_QM_" prefix, in order to find the parameters ...
Initialize the use of the gaussians to treat the QMMM coupling potential.
subroutine, public qmmm_gaussian_initialize(qmmm_gaussian_fns, para_env, pw_env, mm_el_pot_radius, mm_el_pot_radius_corr, qmmm_coupl_type, eps_mm_rspace, maxradius, maxchrg, compatibility, print_section, qmmm_section)
Initialize the Gaussian QMMM Environment.
Initialize a QM/MM calculation.
subroutine, public setup_qmmm_vars_qm(qmmm_section, qmmm_env, subsys_mm, qm_atom_type, qm_atom_index, mm_atom_index, qm_cell_small, qmmm_coupl_type, eps_mm_rspace, qmmm_link, para_env)
...
subroutine, public setup_qmmm_links(qmmm_section, qmmm_links, mm_el_pot_radius, mm_el_pot_radius_corr, mm_atom_index, iw)
this routine sets up all variables to treat qmmm links
subroutine, public qmmm_init_gaussian_type(qmmm_env_qm, para_env, mm_atom_chrg, qs_env, added_charges, added_shells, print_section, qmmm_section)
...
subroutine, public move_or_add_atoms(qmmm_section, move_mm_charges, add_mm_charges, mm_atom_chrg, mm_el_pot_radius, mm_el_pot_radius_corr, added_charges, mm_atom_index)
this routine sets up all variables to treat qmmm links
subroutine, public assign_mm_charges_and_radius(subsys, charges, mm_atom_chrg, mm_el_pot_radius, mm_el_pot_radius_corr, mm_atom_index, mm_link_atoms, mm_link_scale_factor, added_shells, shell_model)
Assigns charges and radius to evaluate the MM electrostatic potential.
subroutine, public setup_origin_mm_cell(qmmm_section, qmmm_env, qm_cell_small, dr)
this routine sets up the origin of the MM cell respect to the origin of the QM cell....
subroutine, public qmmm_init_periodic_potential(qmmm_env_qm, qm_cell_small, mm_cell, para_env, qs_env, added_charges, added_shells, qmmm_periodic, print_section, mm_atom_chrg)
...
subroutine, public print_qmmm_charges(mm_atom_index, mm_atom_chrg, mm_el_pot_radius, mm_el_pot_radius_corr, added_charges, added_shells, qmmm_section, nocompatibility, shell_model)
Print info on charges generating the qmmm potential..
subroutine, public print_qmmm_links(qmmm_section, qmmm_links)
Print info on qm/mm links.
subroutine, public setup_qmmm_vars_mm(qmmm_section, qmmm_env, qm_atom_index, mm_link_atoms, mm_link_scale_factor, fist_scale_charge_link, qmmm_coupl_type, qmmm_link)
...
subroutine, public print_image_charge_info(qmmm_env, qmmm_section)
Print info on image charges.
subroutine, public qmmm_init_potential(qmmm_env_qm, mm_cell, added_charges, added_shells, print_section)
...
Setting up the potential for QM/MM periodic boundary conditions calculations.
subroutine, public qmmm_per_potential_init(qmmm_coupl_type, per_potentials, potentials, pgfs, qm_cell_small, mm_cell, compatibility, qmmm_periodic, print_section, eps_mm_rspace, maxchrg, ncp, ncpl)
Initialize the QMMM potential stored on vector, according the qmmm_coupl_type.
subroutine, public qmmm_ewald_potential_init(ewald_env, ewald_pw, qmmm_coupl_type, mm_cell, para_env, qmmm_periodic, print_section)
Initialize the QMMM Ewald potential needed for QM-QM Coupling using point charges.
subroutine, public create_add_shell_type(added_shells, ndim)
creates the add_shell_type structure
subroutine, public create_add_set_type(added_charges, ndim)
creates the add_set_type structure
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, mimic, 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, sab_cneo, 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, xcint_weights, 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, rhoz_cneo_set, ecoul_1c, rho0_s_rs, rho0_s_gs, rhoz_cneo_s_rs, rhoz_cneo_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, do_rixs, tb_tblite)
Get the QUICKSTEP environment.
elemental subroutine, public get_shell(shell, charge, charge_core, charge_shell, mass_core, mass_shell, k2_spring, k4_spring, max_dist, shell_cutoff)
...
represent a list of objects
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
type of a logger, at the moment it contains just a print level starting at which level it should be l...
represents a system: atoms, molecules, their pos,vel,...
stores all the informations relevant to an mpi environment
represent a list of objects
contained for different pw related things
parameters for core-shell model potentials