Line data Source code
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
10 : !> \author
11 : ! **************************************************************************************************
12 : MODULE shell_opt
13 :
14 : USE atomic_kind_list_types, ONLY: atomic_kind_list_type
15 : USE atomic_kind_types, ONLY: atomic_kind_type,&
16 : get_atomic_kind
17 : USE cg_optimizer, ONLY: geoopt_cg
18 : USE cp_log_handling, ONLY: cp_get_default_logger,&
19 : cp_logger_type
20 : USE cp_output_handling, ONLY: cp_add_iter_level,&
21 : cp_rm_iter_level
22 : USE cp_subsys_types, ONLY: cp_subsys_get,&
23 : cp_subsys_type
24 : USE distribution_1d_types, ONLY: distribution_1d_type
25 : USE force_env_types, ONLY: force_env_get,&
26 : force_env_type
27 : USE global_types, ONLY: global_environment_type
28 : USE gopt_f_types, ONLY: gopt_f_create,&
29 : gopt_f_release,&
30 : gopt_f_type
31 : USE gopt_param_types, ONLY: gopt_param_read,&
32 : gopt_param_type
33 : USE input_constants, ONLY: default_shellcore_method_id
34 : USE input_section_types, ONLY: section_vals_get,&
35 : section_vals_get_subs_vals,&
36 : section_vals_type
37 : USE integrator_utils, ONLY: tmp_variables_type
38 : USE kinds, ONLY: dp
39 : USE message_passing, ONLY: mp_para_env_type
40 : USE particle_types, ONLY: particle_type
41 : USE shell_potential_types, ONLY: shell_kind_type
42 : #include "../base/base_uses.f90"
43 :
44 : IMPLICIT NONE
45 : PRIVATE
46 : PUBLIC :: optimize_shell_core
47 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'shell_opt'
48 :
49 : CONTAINS
50 :
51 : ! **************************************************************************************************
52 : !> \brief Optimize shell-core positions along an MD run
53 : !> \param force_env ...
54 : !> \param particle_set ...
55 : !> \param shell_particle_set ...
56 : !> \param core_particle_set ...
57 : !> \param globenv ...
58 : !> \param tmp ...
59 : !> \param check ...
60 : !> \author
61 : ! **************************************************************************************************
62 :
63 714 : SUBROUTINE optimize_shell_core(force_env, particle_set, shell_particle_set, core_particle_set, globenv, tmp, check)
64 : TYPE(force_env_type), POINTER :: force_env
65 : TYPE(particle_type), DIMENSION(:), POINTER :: particle_set, shell_particle_set, &
66 : core_particle_set
67 : TYPE(global_environment_type), POINTER :: globenv
68 : TYPE(tmp_variables_type), OPTIONAL, POINTER :: tmp
69 : LOGICAL, INTENT(IN), OPTIONAL :: check
70 :
71 : CHARACTER(len=*), PARAMETER :: routineN = 'optimize_shell_core'
72 :
73 : INTEGER :: handle, i, iat, nshell
74 : LOGICAL :: do_update, explicit, my_check, optimize
75 714 : REAL(dp), DIMENSION(:), POINTER :: dvec_sc, dvec_sc_0
76 : TYPE(atomic_kind_list_type), POINTER :: atomic_kinds
77 : TYPE(cp_logger_type), POINTER :: logger
78 : TYPE(cp_subsys_type), POINTER :: subsys
79 : TYPE(distribution_1d_type), POINTER :: local_particles
80 : TYPE(gopt_f_type), POINTER :: gopt_env
81 : TYPE(gopt_param_type), POINTER :: gopt_param
82 : TYPE(mp_para_env_type), POINTER :: para_env
83 : TYPE(section_vals_type), POINTER :: force_env_section, geo_section, &
84 : root_section
85 :
86 714 : NULLIFY (logger)
87 1428 : logger => cp_get_default_logger()
88 :
89 714 : CPASSERT(ASSOCIATED(force_env))
90 714 : CPASSERT(ASSOCIATED(globenv))
91 :
92 714 : NULLIFY (gopt_param, force_env_section, gopt_env, dvec_sc, dvec_sc_0, root_section, geo_section)
93 714 : root_section => force_env%root_section
94 714 : force_env_section => force_env%force_env_section
95 714 : geo_section => section_vals_get_subs_vals(root_section, "MOTION%SHELL_OPT")
96 :
97 714 : CALL section_vals_get(geo_section, explicit=explicit)
98 824 : IF (.NOT. explicit) RETURN
99 :
100 130 : CALL timeset(routineN, handle)
101 :
102 130 : optimize = .FALSE.
103 130 : my_check = .FALSE.
104 130 : IF (PRESENT(check)) my_check = check
105 120 : IF (my_check) THEN
106 120 : NULLIFY (subsys, para_env, atomic_kinds, local_particles)
107 120 : CALL force_env_get(force_env=force_env, subsys=subsys, para_env=para_env)
108 120 : CALL cp_subsys_get(subsys=subsys, atomic_kinds=atomic_kinds, local_particles=local_particles)
109 : CALL check_shell_core_distance(atomic_kinds, local_particles, particle_set, shell_particle_set, &
110 120 : core_particle_set, para_env, optimize)
111 :
112 120 : IF (.NOT. optimize) THEN
113 110 : CALL timestop(handle)
114 110 : RETURN
115 : END IF
116 : END IF
117 :
118 20 : nshell = SIZE(shell_particle_set)
119 60 : ALLOCATE (dvec_sc(3*nshell))
120 40 : ALLOCATE (dvec_sc_0(3*nshell))
121 1940 : DO i = 1, nshell
122 1920 : dvec_sc(1 + 3*(i - 1)) = core_particle_set(i)%r(1) - shell_particle_set(i)%r(1)
123 1920 : dvec_sc(2 + 3*(i - 1)) = core_particle_set(i)%r(2) - shell_particle_set(i)%r(2)
124 1940 : dvec_sc(3 + 3*(i - 1)) = core_particle_set(i)%r(3) - shell_particle_set(i)%r(3)
125 : END DO
126 11540 : dvec_sc_0 = dvec_sc
127 :
128 20 : ALLOCATE (gopt_param)
129 20 : CALL gopt_param_read(gopt_param, geo_section, type_id=default_shellcore_method_id)
130 : CALL gopt_f_create(gopt_env, gopt_param, force_env=force_env, globenv=globenv, &
131 20 : geo_opt_section=geo_section)
132 :
133 20 : CALL cp_add_iter_level(logger%iter_info, "SHELL_OPT")
134 20 : gopt_env%eval_opt_geo = .FALSE.
135 : CALL geoopt_cg(force_env, gopt_param, globenv, &
136 20 : geo_section, gopt_env, dvec_sc, do_update=do_update)
137 20 : IF (.NOT. do_update) THEN
138 0 : DO i = 1, nshell
139 0 : shell_particle_set(i)%r(1) = -dvec_sc_0(1 + 3*(i - 1)) + core_particle_set(i)%r(1)
140 0 : shell_particle_set(i)%r(2) = -dvec_sc_0(2 + 3*(i - 1)) + core_particle_set(i)%r(2)
141 0 : shell_particle_set(i)%r(3) = -dvec_sc_0(3 + 3*(i - 1)) + core_particle_set(i)%r(3)
142 : END DO
143 : END IF
144 20 : CALL cp_rm_iter_level(logger%iter_info, "SHELL_OPT")
145 :
146 20 : CALL gopt_f_release(gopt_env)
147 20 : DEALLOCATE (gopt_param)
148 20 : DEALLOCATE (dvec_sc)
149 20 : DEALLOCATE (dvec_sc_0)
150 :
151 20 : IF (PRESENT(tmp)) THEN
152 970 : DO i = 1, nshell
153 960 : iat = shell_particle_set(i)%atom_index
154 3840 : tmp%shell_vel(1:3, i) = tmp%vel(1:3, iat)
155 3850 : tmp%core_vel(1:3, i) = tmp%vel(1:3, iat)
156 : END DO
157 : ELSE
158 970 : DO i = 1, nshell
159 960 : iat = shell_particle_set(i)%atom_index
160 6720 : shell_particle_set(i)%v(1:3) = particle_set(iat)%v(1:3)
161 6730 : core_particle_set(i)%v(1:3) = particle_set(iat)%v(1:3)
162 : END DO
163 : END IF
164 :
165 20 : CALL timestop(handle)
166 :
167 1578 : END SUBROUTINE optimize_shell_core
168 :
169 : ! **************************************************************************************************
170 : !> \brief Check shell_core_distance
171 : !> \param atomic_kinds ...
172 : !> \param local_particles ...
173 : !> \param particle_set ...
174 : !> \param shell_particle_set ...
175 : !> \param core_particle_set ...
176 : !> \param para_env ...
177 : !> \param optimize ...
178 : !> \par History
179 : !> none
180 : !> \author MI (October 2008)
181 : !> I soliti ignoti
182 : ! **************************************************************************************************
183 120 : SUBROUTINE check_shell_core_distance(atomic_kinds, local_particles, particle_set, &
184 : shell_particle_set, core_particle_set, para_env, optimize)
185 :
186 : TYPE(atomic_kind_list_type), POINTER :: atomic_kinds
187 : TYPE(distribution_1d_type), POINTER :: local_particles
188 : TYPE(particle_type), DIMENSION(:), POINTER :: particle_set, shell_particle_set, &
189 : core_particle_set
190 : TYPE(mp_para_env_type), POINTER :: para_env
191 : LOGICAL, INTENT(INOUT) :: optimize
192 :
193 : INTEGER :: ikind, iparticle, iparticle_local, &
194 : itest, nkind, nparticle_local, &
195 : shell_index
196 : LOGICAL :: is_shell
197 : REAL(dp) :: dsc, rc(3), rs(3)
198 : TYPE(atomic_kind_type), POINTER :: atomic_kind
199 : TYPE(shell_kind_type), POINTER :: shell
200 :
201 120 : nkind = atomic_kinds%n_els
202 120 : itest = 0
203 360 : DO ikind = 1, nkind
204 240 : NULLIFY (atomic_kind)
205 240 : atomic_kind => atomic_kinds%els(ikind)
206 240 : CALL get_atomic_kind(atomic_kind=atomic_kind, shell_active=is_shell, shell=shell)
207 360 : IF (is_shell) THEN
208 240 : IF (shell%max_dist > 0.0_dp) THEN
209 240 : nparticle_local = local_particles%n_el(ikind)
210 6000 : DO iparticle_local = 1, nparticle_local
211 5760 : iparticle = local_particles%list(ikind)%array(iparticle_local)
212 5760 : shell_index = particle_set(iparticle)%shell_index
213 :
214 23040 : rc(:) = core_particle_set(shell_index)%r(:)
215 23040 : rs(:) = shell_particle_set(shell_index)%r(:)
216 5760 : dsc = SQRT((rc(1) - rs(1))**2 + (rc(2) - rs(2))**2 + (rc(3) - rs(3))**2)
217 6000 : IF (dsc > shell%max_dist) THEN
218 15 : itest = 1
219 : END IF
220 : END DO
221 : END IF
222 : END IF
223 : END DO
224 :
225 120 : CALL para_env%sum(itest)
226 120 : IF (itest > 0) optimize = .TRUE.
227 :
228 120 : END SUBROUTINE check_shell_core_distance
229 : END MODULE shell_opt
|