42#include "../base/base_uses.f90"
47 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'shell_opt'
63 SUBROUTINE optimize_shell_core(force_env, particle_set, shell_particle_set, core_particle_set, globenv, tmp, check)
65 TYPE(
particle_type),
DIMENSION(:),
POINTER :: particle_set, shell_particle_set, &
69 LOGICAL,
INTENT(IN),
OPTIONAL :: check
71 CHARACTER(len=*),
PARAMETER :: routinen =
'optimize_shell_core'
73 INTEGER :: handle, i, iat, nshell
74 LOGICAL :: do_update, explicit, my_check, optimize
75 REAL(
dp),
DIMENSION(:),
POINTER :: dvec_sc, dvec_sc_0
89 cpassert(
ASSOCIATED(force_env))
90 cpassert(
ASSOCIATED(globenv))
92 NULLIFY (gopt_param, force_env_section, gopt_env, dvec_sc, dvec_sc_0, root_section, geo_section)
93 root_section => force_env%root_section
94 force_env_section => force_env%force_env_section
98 IF (.NOT. explicit)
RETURN
100 CALL timeset(routinen, handle)
104 IF (
PRESENT(check)) my_check = check
106 NULLIFY (subsys, para_env, atomic_kinds, local_particles)
107 CALL force_env_get(force_env=force_env, subsys=subsys, para_env=para_env)
108 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 core_particle_set, para_env, optimize)
112 IF (.NOT. optimize)
THEN
113 CALL timestop(handle)
118 nshell =
SIZE(shell_particle_set)
119 ALLOCATE (dvec_sc(3*nshell))
120 ALLOCATE (dvec_sc_0(3*nshell))
122 dvec_sc(1 + 3*(i - 1)) = core_particle_set(i)%r(1) - shell_particle_set(i)%r(1)
123 dvec_sc(2 + 3*(i - 1)) = core_particle_set(i)%r(2) - shell_particle_set(i)%r(2)
124 dvec_sc(3 + 3*(i - 1)) = core_particle_set(i)%r(3) - shell_particle_set(i)%r(3)
128 ALLOCATE (gopt_param)
130 CALL gopt_f_create(gopt_env, gopt_param, force_env=force_env, globenv=globenv, &
131 geo_opt_section=geo_section)
134 gopt_env%eval_opt_geo = .false.
135 CALL geoopt_cg(force_env, gopt_param, globenv, &
136 geo_section, gopt_env, dvec_sc, do_update=do_update)
137 IF (.NOT. do_update)
THEN
139 shell_particle_set(i)%r(1) = -dvec_sc_0(1 + 3*(i - 1)) + core_particle_set(i)%r(1)
140 shell_particle_set(i)%r(2) = -dvec_sc_0(2 + 3*(i - 1)) + core_particle_set(i)%r(2)
141 shell_particle_set(i)%r(3) = -dvec_sc_0(3 + 3*(i - 1)) + core_particle_set(i)%r(3)
147 DEALLOCATE (gopt_param)
149 DEALLOCATE (dvec_sc_0)
151 IF (
PRESENT(tmp))
THEN
153 iat = shell_particle_set(i)%atom_index
154 tmp%shell_vel(1:3, i) = tmp%vel(1:3, iat)
155 tmp%core_vel(1:3, i) = tmp%vel(1:3, iat)
159 iat = shell_particle_set(i)%atom_index
160 shell_particle_set(i)%v(1:3) = particle_set(iat)%v(1:3)
161 core_particle_set(i)%v(1:3) = particle_set(iat)%v(1:3)
165 CALL timestop(handle)
183 SUBROUTINE check_shell_core_distance(atomic_kinds, local_particles, particle_set, &
184 shell_particle_set, core_particle_set, para_env, optimize)
188 TYPE(
particle_type),
DIMENSION(:),
POINTER :: particle_set, shell_particle_set, &
191 LOGICAL,
INTENT(INOUT) :: optimize
193 INTEGER :: ikind, iparticle, iparticle_local, &
194 itest, nkind, nparticle_local, &
197 REAL(
dp) :: dsc, rc(3), rs(3)
201 nkind = atomic_kinds%n_els
204 NULLIFY (atomic_kind)
205 atomic_kind => atomic_kinds%els(ikind)
206 CALL get_atomic_kind(atomic_kind=atomic_kind, shell_active=is_shell, shell=shell)
208 IF (shell%max_dist > 0.0_dp)
THEN
209 nparticle_local = local_particles%n_el(ikind)
210 DO iparticle_local = 1, nparticle_local
211 iparticle = local_particles%list(ikind)%array(iparticle_local)
212 shell_index = particle_set(iparticle)%shell_index
214 rc(:) = core_particle_set(shell_index)%r(:)
215 rs(:) = shell_particle_set(shell_index)%r(:)
216 dsc = sqrt((rc(1) - rs(1))**2 + (rc(2) - rs(2))**2 + (rc(3) - rs(3))**2)
217 IF (dsc > shell%max_dist)
THEN
225 CALL para_env%sum(itest)
226 IF (itest > 0) optimize = .true.
228 END SUBROUTINE check_shell_core_distance
represent a simple array based list of the given type
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.
Routines for Geometry optimization using Conjugate Gradients.
recursive subroutine, public geoopt_cg(force_env, gopt_param, globenv, geo_section, gopt_env, x0, do_update)
Driver for conjugate gradient optimization technique.
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...
subroutine, public cp_rm_iter_level(iteration_info, level_name, n_rlevel_att)
Removes an iteration level.
subroutine, public cp_add_iter_level(iteration_info, level_name, n_rlevel_new)
Adds an iteration level.
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, ipi_env)
returns various attributes about the force environment
Define type storing the global information of a run. Keep the amount of stored data small....
contains a functional that calculates the energy and its derivatives for the geometry optimizer
recursive subroutine, public gopt_f_create(gopt_env, gopt_param, force_env, globenv, geo_opt_section, eval_opt_geo)
...
recursive subroutine, public gopt_f_release(gopt_env)
...
contains typo and related routines to handle parameters controlling the GEO_OPT module
subroutine, public gopt_param_read(gopt_param, gopt_section, type_id)
reads the parameters of the geopmetry optimizer
Provides integrator utility routines for the integrators.
Defines the basic variable types.
integer, parameter, public dp
Interface to the message passing library MPI.
Define the data structure for the particle information.
subroutine, public optimize_shell_core(force_env, particle_set, shell_particle_set, core_particle_set, globenv, tmp, check)
Optimize shell-core positions along an MD run.
represent a list of objects
Provides all information about an atomic kind.
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,...
structure to store local (to a processor) ordered lists of integers.
wrapper to abstract the force evaluation of the various methods
contains the initially parsed file and the initial parallel environment
calculates the potential energy of a system, and its derivatives
stores all the informations relevant to an mpi environment