(git:ccc2433)
qmmmx_update.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 Update a QM/MM calculations with force mixing
10 !> \par History
11 !> 5.2004 created [fawzi]
12 !> \author Fawzi Mohamed
13 ! **************************************************************************************************
15  USE atomic_kind_list_types, ONLY: atomic_kind_list_type
16  USE cp_subsys_types, ONLY: cp_subsys_get,&
17  cp_subsys_type
18  USE distribution_1d_types, ONLY: distribution_1d_type
19  USE force_env_types, ONLY: force_env_get,&
20  force_env_type
25  section_vals_type
26  USE qmmm_create, ONLY: qmmm_env_create
27  USE qmmm_types, ONLY: qmmm_env_get
28  USE qmmmx_types, ONLY: qmmmx_env_release,&
29  qmmmx_env_type
32 #include "./base/base_uses.f90"
33 
34  IMPLICIT NONE
35  PRIVATE
36 
37  LOGICAL, PRIVATE, PARAMETER :: debug_this_module = .true.
38  CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qmmmx_update'
39 
40  PUBLIC :: qmmmx_update_force_env
41 
42 CONTAINS
43 
44 ! **************************************************************************************************
45 !> \brief ...
46 !> \param force_env ...
47 !> \param root_section ...
48 ! **************************************************************************************************
49  SUBROUTINE qmmmx_update_force_env(force_env, root_section)
50  TYPE(force_env_type), POINTER :: force_env
51  TYPE(section_vals_type), POINTER :: root_section
52 
53  LOGICAL :: force_mixing_active, labels_changed
54  TYPE(atomic_kind_list_type), POINTER :: atomic_kinds, new_atomic_kinds
55  TYPE(cp_subsys_type), POINTER :: subsys, subsys_new
56  TYPE(distribution_1d_type), POINTER :: local_particles, new_local_particles
57  TYPE(qmmmx_env_type) :: new_qmmmx_env
58  TYPE(section_vals_type), POINTER :: qmmm_core_section, &
59  qmmm_extended_section, &
60  qmmm_force_mixing, qmmm_section, &
61  subsys_section
62 
63 ! check everything for not null, because sometimes (e.g. metadynamics in parallel) it happens
64 
65  IF (.NOT. ASSOCIATED(force_env)) RETURN
66  IF (.NOT. ASSOCIATED(force_env%force_env_section)) RETURN
67  ! these two should never happen, because the sections exist, but just in case...
68  qmmm_section => section_vals_get_subs_vals(force_env%force_env_section, "QMMM", can_return_null=.true.)
69  IF (.NOT. ASSOCIATED(qmmm_section)) RETURN
70  qmmm_force_mixing => section_vals_get_subs_vals(qmmm_section, "FORCE_MIXING", can_return_null=.true.)
71  IF (.NOT. ASSOCIATED(qmmm_force_mixing)) RETURN
72  CALL section_vals_get(qmmm_force_mixing, explicit=force_mixing_active)
73  IF (.NOT. force_mixing_active) RETURN
74  IF (.NOT. ASSOCIATED(force_env%qmmmx_env)) cpabort("force_env%qmmmx not associated")
75 
76  CALL force_env_get(force_env, subsys=subsys)
77  CALL update_force_mixing_labels(subsys, qmmm_section, labels_changed=labels_changed)
78  IF (.NOT. labels_changed) RETURN
79  cpwarn("Adaptive force-mixing labels changed, rebuilding QM/MM calculations! ")
80 
81  CALL update_force_eval(force_env, root_section, .false.)
82 
83  ! using CUR_INDICES and CUR_LABELS, create appropriate QM_KIND sections for two QM/MM calculations
84  CALL setup_force_mixing_qmmm_sections(subsys, qmmm_section, qmmm_core_section, qmmm_extended_section)
85 
86  subsys_section => section_vals_get_subs_vals(force_env%force_env_section, "SUBSYS")
87  ![ADAPT] no sure about use_motion_section
88  ALLOCATE (new_qmmmx_env%core)
89  CALL qmmm_env_create(new_qmmmx_env%core, &
90  force_env%root_section, force_env%para_env, force_env%globenv, &
91  force_env%force_env_section, qmmm_core_section, subsys_section, use_motion_section=.true., &
92  prev_subsys=subsys, ignore_outside_box=.true.)
93  ALLOCATE (new_qmmmx_env%ext)
94  CALL qmmm_env_create(new_qmmmx_env%ext, &
95  force_env%root_section, force_env%para_env, force_env%globenv, &
96  force_env%force_env_section, qmmm_extended_section, subsys_section, use_motion_section=.true., &
97  prev_subsys=subsys, ignore_outside_box=.true.)
98 
99  ! [NB] need to copy wiener process data, since it's not recreated when
100  ! fist subsys is recreated by qmmm_env_create
101  CALL qmmm_env_get(force_env%qmmmx_env%core, subsys=subsys)
102  CALL cp_subsys_get(subsys, atomic_kinds=atomic_kinds, local_particles=local_particles)
103  CALL qmmm_env_get(new_qmmmx_env%core, subsys=subsys_new)
104  CALL cp_subsys_get(subsys_new, atomic_kinds=new_atomic_kinds, local_particles=new_local_particles)
105  IF (ASSOCIATED(local_particles%local_particle_set)) THEN
106  CALL copy_wiener_process(atomic_kinds, local_particles, new_atomic_kinds, new_local_particles)
107  END IF
108 
109  CALL qmmm_env_get(force_env%qmmmx_env%ext, subsys=subsys)
110  CALL cp_subsys_get(subsys, atomic_kinds=atomic_kinds, local_particles=local_particles)
111  CALL qmmm_env_get(new_qmmmx_env%ext, subsys=subsys_new)
112  CALL cp_subsys_get(subsys_new, atomic_kinds=new_atomic_kinds, local_particles=new_local_particles)
113  IF (ASSOCIATED(local_particles%local_particle_set)) THEN
114  CALL copy_wiener_process(atomic_kinds, local_particles, new_atomic_kinds, new_local_particles)
115  END IF
116 
117  CALL section_vals_release(qmmm_core_section)
118  CALL section_vals_release(qmmm_extended_section)
119 
120  ! release old qmmmx_env and point to new one
121  CALL qmmmx_env_release(force_env%qmmmx_env)
122  force_env%qmmmx_env = new_qmmmx_env
123 
124  END SUBROUTINE qmmmx_update_force_env
125 
126 ! **************************************************************************************************
127 !> \brief ...
128 !> \param from_local_particle_kinds ...
129 !> \param from_local_particles ...
130 !> \param to_local_particle_kinds ...
131 !> \param to_local_particles ...
132 ! **************************************************************************************************
133  SUBROUTINE copy_wiener_process(from_local_particle_kinds, from_local_particles, &
134  to_local_particle_kinds, to_local_particles)
135  TYPE(atomic_kind_list_type), POINTER :: from_local_particle_kinds
136  TYPE(distribution_1d_type), POINTER :: from_local_particles
137  TYPE(atomic_kind_list_type), POINTER :: to_local_particle_kinds
138  TYPE(distribution_1d_type), POINTER :: to_local_particles
139 
140  CHARACTER(LEN=*), PARAMETER :: routinen = 'copy_wiener_process'
141 
142  INTEGER :: from_iparticle_kind, from_iparticle_local(1), from_nparticle_kind, &
143  from_nparticle_local, handle, to_iparticle_global, to_iparticle_kind, to_iparticle_local, &
144  to_nparticle_kind, to_nparticle_local, tot_from_nparticle_local, tot_to_nparticle_local
145  LOGICAL :: found_it
146 
147  CALL timeset(routinen, handle)
148  cpassert(ASSOCIATED(from_local_particles))
149  cpassert(ASSOCIATED(to_local_particles))
150 
151  IF (.NOT. ASSOCIATED(from_local_particles%local_particle_set)) RETURN
152  cpassert(.NOT. ASSOCIATED(to_local_particles%local_particle_set))
153 
154  from_nparticle_kind = from_local_particle_kinds%n_els
155  to_nparticle_kind = to_local_particle_kinds%n_els
156 
157  ! make sure total number of particles hasn't changed, even if particle kinds have
158  tot_from_nparticle_local = 0
159  DO from_iparticle_kind = 1, from_nparticle_kind
160  tot_from_nparticle_local = tot_from_nparticle_local + from_local_particles%n_el(from_iparticle_kind)
161  END DO
162  tot_to_nparticle_local = 0
163  DO to_iparticle_kind = 1, to_nparticle_kind
164  tot_to_nparticle_local = tot_to_nparticle_local + to_local_particles%n_el(to_iparticle_kind)
165  END DO
166  cpassert(tot_from_nparticle_local == tot_to_nparticle_local)
167 
168  ALLOCATE (to_local_particles%local_particle_set(to_nparticle_kind))
169  DO to_iparticle_kind = 1, to_nparticle_kind
170 
171  to_nparticle_local = to_local_particles%n_el(to_iparticle_kind)
172  ALLOCATE (to_local_particles%local_particle_set(to_iparticle_kind)%rng(to_nparticle_local))
173 
174  DO to_iparticle_local = 1, to_nparticle_local
175  to_iparticle_global = to_local_particles%list(to_iparticle_kind)%array(to_iparticle_local)
176  ALLOCATE (to_local_particles%local_particle_set(to_iparticle_kind)%rng(to_iparticle_local)%stream)
177 
178  found_it = .false.
179  ! find the matching kind/index where this particle was before
180  DO from_iparticle_kind = 1, from_nparticle_kind
181  from_nparticle_local = from_local_particles%n_el(from_iparticle_kind)
182  IF (minval(abs(from_local_particles%list(from_iparticle_kind)%array(1:from_nparticle_local) - &
183  to_iparticle_global)) == 0) THEN
184  from_iparticle_local = &
185  minloc(abs(from_local_particles%list(from_iparticle_kind)%array(1:from_nparticle_local) - &
186  to_iparticle_global))
187  to_local_particles%local_particle_set(to_iparticle_kind)%rng(to_iparticle_local)%stream = &
188  from_local_particles%local_particle_set(from_iparticle_kind)%rng(from_iparticle_local(1))%stream
189  found_it = .true.
190  EXIT
191  END IF
192  END DO
193  cpassert(found_it)
194 
195  END DO ! to_iparticle_local
196 
197  END DO ! to_iparticle_kind
198  CALL timestop(handle)
199 
200  END SUBROUTINE copy_wiener_process
201 
202 END MODULE qmmmx_update
represent a simple array based list of the given type
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)
returns information about various attributes of the given subsys
stores a lists of integer that are local to a processor. The idea is that these integers represent ob...
Interface for the force calculations.
recursive subroutine, public force_env_get(force_env, in_use, fist_env, qs_env, meta_env, fp_env, subsys, para_env, potential_energy, additional_potential, kinetic_energy, harmonic_shell, kinetic_shell, cell, sub_force_env, qmmm_env, qmmmx_env, eip_env, pwdft_env, globenv, input, force_env_section, method_name_id, root_section, mixed_env, nnp_env, embed_env)
returns various attributes about the force environment
subroutine, public update_force_eval(force_env, root_section, write_binary_restart_file, respa)
Updates the force_eval section of the input file.
objects that represent the structure of input sections and the data contained in an input section
recursive type(section_vals_type) function, pointer, public section_vals_get_subs_vals(section_vals, subsection_name, i_rep_section, can_return_null)
returns the values of the requested subsection
subroutine, public section_vals_get(section_vals, ref_count, n_repetition, n_subs_vals_rep, section, explicit)
returns various attributes about the section_vals
recursive subroutine, public section_vals_release(section_vals)
releases the given object
Initialize a QM/MM calculation.
Definition: qmmm_create.F:14
subroutine, public qmmm_env_create(qmmm_env, root_section, para_env, globenv, force_env_section, qmmm_section, subsys_section, use_motion_section, prev_subsys, ignore_outside_box)
...
Definition: qmmm_create.F:104
Basic container type for QM/MM.
Definition: qmmm_types.F:12
subroutine, public qmmm_env_get(qmmm_env, subsys, potential_energy, kinetic_energy)
...
Definition: qmmm_types.F:50
Basic container type for QM/MM with force mixing.
Definition: qmmmx_types.F:12
subroutine, public qmmmx_env_release(qmmmx_env)
releases the given qmmmx_env (see doc/ReferenceCounting.html)
Definition: qmmmx_types.F:61
Update a QM/MM calculations with force mixing.
Definition: qmmmx_update.F:14
subroutine, public qmmmx_update_force_env(force_env, root_section)
...
Definition: qmmmx_update.F:50
Routines used for force-mixing QM/MM calculations.
Definition: qmmmx_util.F:14
subroutine, public setup_force_mixing_qmmm_sections(subsys, qmmm_section, qmmm_core_section, qmmm_extended_section)
...
Definition: qmmmx_util.F:587
subroutine, public update_force_mixing_labels(subsys, qmmm_section, labels_changed)
...
Definition: qmmmx_util.F:116