(git:34ef472)
md_run.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 Perform a molecular dynamics (MD) run using QUICKSTEP
10 !> \par History
11 !> - Added support for Langevin regions (2014/02/05, LT)
12 !> \author Matthias Krack (07.11.2002)
13 ! **************************************************************************************************
14 MODULE md_run
15  USE atomic_kind_list_types, ONLY: atomic_kind_list_type
16  USE averages_types, ONLY: average_quantities_type
17  USE barostat_types, ONLY: barostat_type,&
19  USE cell_types, ONLY: cell_type
22  cp_logger_type
24  cp_iterate,&
25  cp_p_file,&
30  USE cp_subsys_types, ONLY: cp_subsys_get,&
31  cp_subsys_type
32  USE distribution_1d_types, ONLY: distribution_1d_type
34  USE force_env_types, ONLY: force_env_get,&
35  force_env_type
37  USE free_energy_types, ONLY: fe_env_create,&
38  free_energy_type
39  USE global_types, ONLY: global_environment_type
40  USE input_constants, ONLY: &
47  section_vals_type,&
49  USE kinds, ONLY: default_string_length,&
50  dp
51  USE machine, ONLY: m_walltime
52  USE md_ener_types, ONLY: create_md_ener,&
53  md_ener_type
54  USE md_energies, ONLY: initialize_md_ener,&
56  md_energy,&
58  USE md_environment_types, ONLY: get_md_env,&
61  md_environment_type,&
64  USE md_util, ONLY: md_output
65  USE md_vel_utils, ONLY: angvel_control,&
70  USE mdctrl_types, ONLY: mdctrl_type
71  USE message_passing, ONLY: mp_para_env_type
77  USE particle_list_types, ONLY: particle_list_type
79  USE reftraj_types, ONLY: create_reftraj,&
80  reftraj_type
81  USE reftraj_util, ONLY: initialize_reftraj,&
83  USE rt_propagation, ONLY: rt_prop_setup
85  USE simpar_types, ONLY: create_simpar_type,&
87  simpar_type
88  USE thermal_region_types, ONLY: thermal_regions_type
92  USE thermostat_types, ONLY: thermostats_type
95  USE virial_types, ONLY: virial_type
98 #include "../base/base_uses.f90"
99 
100  IMPLICIT NONE
101 
102  PRIVATE
103 
104  CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'md_run'
105 
106  PUBLIC :: qs_mol_dyn
107 
108 CONTAINS
109 
110 ! **************************************************************************************************
111 !> \brief Main driver module for Molecular Dynamics
112 !> \param force_env ...
113 !> \param globenv ...
114 !> \param averages ...
115 !> \param rm_restart_info ...
116 !> \param hmc_e_initial ...
117 !> \param hmc_e_final ...
118 !> \param mdctrl ...
119 ! **************************************************************************************************
120  SUBROUTINE qs_mol_dyn(force_env, globenv, averages, rm_restart_info, hmc_e_initial, hmc_e_final, mdctrl)
121 
122  TYPE(force_env_type), POINTER :: force_env
123  TYPE(global_environment_type), POINTER :: globenv
124  TYPE(average_quantities_type), OPTIONAL, POINTER :: averages
125  LOGICAL, INTENT(IN), OPTIONAL :: rm_restart_info
126  REAL(kind=dp), OPTIONAL :: hmc_e_initial, hmc_e_final
127  TYPE(mdctrl_type), OPTIONAL, POINTER :: mdctrl
128 
129  LOGICAL :: my_rm_restart_info
130  TYPE(md_environment_type), POINTER :: md_env
131  TYPE(mp_para_env_type), POINTER :: para_env
132  TYPE(section_vals_type), POINTER :: md_section, motion_section
133 
134  my_rm_restart_info = .true.
135  IF (PRESENT(rm_restart_info)) my_rm_restart_info = rm_restart_info
136  NULLIFY (md_env, para_env)
137  para_env => force_env%para_env
138  motion_section => section_vals_get_subs_vals(force_env%root_section, "MOTION")
139  md_section => section_vals_get_subs_vals(motion_section, "MD")
140 
141  ! Real call to MD driver - Low Level
142  ALLOCATE (md_env)
143  CALL md_env_create(md_env, md_section, para_env, force_env=force_env)
144  CALL set_md_env(md_env, averages=averages)
145  IF (PRESENT(hmc_e_initial) .AND. PRESENT(hmc_e_final)) THEN
146  CALL qs_mol_dyn_low(md_env, md_section, motion_section, force_env, globenv, &
147  hmc_e_initial=hmc_e_initial, hmc_e_final=hmc_e_final)
148  ELSE
149  CALL qs_mol_dyn_low(md_env, md_section, motion_section, force_env, globenv, mdctrl=mdctrl)
150  END IF
151  CALL md_env_release(md_env)
152  DEALLOCATE (md_env)
153 
154  ! Clean restartable sections..
155  IF (my_rm_restart_info) CALL remove_restart_info(force_env%root_section)
156  END SUBROUTINE qs_mol_dyn
157 
158 ! **************************************************************************************************
159 !> \brief Purpose: Driver routine for MD run using QUICKSTEP.
160 !> \param md_env ...
161 !> \param md_section ...
162 !> \param motion_section ...
163 !> \param force_env ...
164 !> \param globenv ...
165 !> \param hmc_e_initial ...
166 !> \param hmc_e_final ...
167 !> \param mdctrl ...
168 !> \par History
169 !> - Cleaning (09.2007) Teodoro Laino [tlaino] - University of Zurich
170 !> - Added lines to print out langevin regions (2014/02/04, LT)
171 !> \author Creation (07.11.2002,MK)
172 ! **************************************************************************************************
173  SUBROUTINE qs_mol_dyn_low(md_env, md_section, motion_section, force_env, globenv, hmc_e_initial, hmc_e_final, mdctrl)
174 
175  TYPE(md_environment_type), POINTER :: md_env
176  TYPE(section_vals_type), POINTER :: md_section, motion_section
177  TYPE(force_env_type), POINTER :: force_env
178  TYPE(global_environment_type), POINTER :: globenv
179  REAL(kind=dp), OPTIONAL :: hmc_e_initial, hmc_e_final
180  TYPE(mdctrl_type), OPTIONAL, POINTER :: mdctrl
181 
182  CHARACTER(LEN=*), PARAMETER :: routinen = 'qs_mol_dyn_low'
183 
184  CHARACTER(LEN=default_string_length) :: my_act, my_pos
185  INTEGER :: handle, i, istep, md_stride, run_type_id
186  INTEGER, POINTER :: itimes
187  LOGICAL :: check, ehrenfest_md, save_mem, &
188  should_stop, write_binary_restart_file
189  REAL(kind=dp) :: dummy, time_iter_start, time_iter_stop
190  REAL(kind=dp), POINTER :: constant, time, used_time
191  TYPE(atomic_kind_list_type), POINTER :: atomic_kinds
192  TYPE(barostat_type), POINTER :: barostat
193  TYPE(cell_type), POINTER :: cell
194  TYPE(cp_logger_type), POINTER :: logger
195  TYPE(cp_subsys_type), POINTER :: subsys, subsys_i
196  TYPE(distribution_1d_type), POINTER :: local_particles
197  TYPE(free_energy_type), POINTER :: fe_env
198  TYPE(md_ener_type), POINTER :: md_ener
199  TYPE(mp_para_env_type), POINTER :: para_env
200  TYPE(particle_list_type), POINTER :: particles
201  TYPE(reftraj_type), POINTER :: reftraj
202  TYPE(section_vals_type), POINTER :: constraint_section, force_env_section, &
203  free_energy_section, global_section, reftraj_section, subsys_section, work_section
204  TYPE(simpar_type), POINTER :: simpar
205  TYPE(thermal_regions_type), POINTER :: thermal_regions
206  TYPE(thermostats_type), POINTER :: thermostats
207  TYPE(virial_type), POINTER :: virial
208 
209  CALL timeset(routinen, handle)
210  cpassert(ASSOCIATED(globenv))
211  cpassert(ASSOCIATED(force_env))
212 
213  NULLIFY (particles, cell, simpar, itimes, used_time, subsys, &
214  md_ener, thermostats, barostat, reftraj, force_env_section, &
215  reftraj_section, work_section, atomic_kinds, &
216  local_particles, time, fe_env, free_energy_section, &
217  constraint_section, thermal_regions, virial, subsys_i)
218  logger => cp_get_default_logger()
219  para_env => force_env%para_env
220 
221  global_section => section_vals_get_subs_vals(force_env%root_section, "GLOBAL")
222  free_energy_section => section_vals_get_subs_vals(motion_section, "FREE_ENERGY")
223  constraint_section => section_vals_get_subs_vals(motion_section, "CONSTRAINT")
224  CALL section_vals_val_get(global_section, "SAVE_MEM", l_val=save_mem)
225 
226  CALL section_vals_val_get(global_section, "RUN_TYPE", i_val=run_type_id)
227  IF (run_type_id == ehrenfest) CALL set_md_env(md_env, ehrenfest_md=.true.)
228 
229  CALL create_simpar_type(simpar)
230  force_env_section => force_env%force_env_section
231  subsys_section => section_vals_get_subs_vals(force_env_section, "SUBSYS")
232  CALL cp_add_iter_level(logger%iter_info, "MD")
233  CALL cp_iterate(logger%iter_info, iter_nr=0)
234  ! Read MD section
235  CALL read_md_section(simpar, motion_section, md_section)
236  ! Setup print_keys
237  simpar%info_constraint = cp_print_key_unit_nr(logger, constraint_section, &
238  "CONSTRAINT_INFO", extension=".shakeLog", log_filename=.false.)
239  simpar%lagrange_multipliers = cp_print_key_unit_nr(logger, constraint_section, &
240  "LAGRANGE_MULTIPLIERS", extension=".LagrangeMultLog", log_filename=.false.)
241  simpar%dump_lm = btest(cp_print_key_should_output(logger%iter_info, constraint_section, &
242  "LAGRANGE_MULTIPLIERS"), cp_p_file)
243 
244  ! Create the structure for the md energies
245  ALLOCATE (md_ener)
246  CALL create_md_ener(md_ener)
247  CALL set_md_env(md_env, md_ener=md_ener)
248  NULLIFY (md_ener)
249 
250  ! If requested setup Thermostats
251  CALL create_thermostats(thermostats, md_section, force_env, simpar, para_env, &
252  globenv, global_section)
253 
254  ! If requested setup Barostat
255  CALL create_barostat_type(barostat, md_section, force_env, simpar, globenv)
256 
257  ! If requested setup different thermal regions
258  CALL create_thermal_regions(thermal_regions, md_section, simpar, force_env)
259 
260  ! If doing langevin_ensemble, then print out langevin_regions information upon request
261  IF (simpar%ensemble == langevin_ensemble) THEN
262  my_pos = "REWIND"
263  my_act = "WRITE"
264  CALL print_thermal_regions_langevin(thermal_regions, simpar, &
265  pos=my_pos, act=my_act)
266  END IF
267 
268  CALL set_md_env(md_env, thermostats=thermostats, barostat=barostat, thermal_regions=thermal_regions)
269 
270  CALL get_md_env(md_env, ehrenfest_md=ehrenfest_md)
271 
272  !If requested set up the REFTRAJ run
273  IF (simpar%ensemble == reftraj_ensemble .AND. ehrenfest_md) &
274  cpabort("Ehrenfest MD does not support reftraj ensemble ")
275  IF (simpar%ensemble == reftraj_ensemble) THEN
276  reftraj_section => section_vals_get_subs_vals(md_section, "REFTRAJ")
277  ALLOCATE (reftraj)
278  CALL create_reftraj(reftraj, reftraj_section, para_env)
279  CALL set_md_env(md_env, reftraj=reftraj)
280  END IF
281 
282  CALL force_env_get(force_env, subsys=subsys, cell=cell, &
283  force_env_section=force_env_section)
284  CALL cp_subsys_get(subsys, virial=virial)
285 
286  ! Set V0 if needed
287  IF (simpar%ensemble == nph_uniaxial_ensemble .OR. simpar%ensemble == nph_uniaxial_damped_ensemble) THEN
288  IF (simpar%v0 == 0._dp) simpar%v0 = cell%deth
289  END IF
290 
291  ! Initialize velocities possibly applying constraints at the zeroth MD step
292  CALL section_vals_val_get(motion_section, "PRINT%RESTART%SPLIT_RESTART_FILE", &
293  l_val=write_binary_restart_file)
294  CALL setup_velocities(force_env, simpar, globenv, md_env, md_section, constraint_section, &
295  write_binary_restart_file)
296 
297  ! Setup Free Energy Calculation (if required)
298  CALL fe_env_create(fe_env, free_energy_section)
299 
300  CALL set_md_env(md_env=md_env, simpar=simpar, fe_env=fe_env, cell=cell, &
301  force_env=force_env)
302 
303  ! Possibly initialize Wiener processes
304  ![NB] Tested again within create_wiener_process. Why??
305  IF (need_per_atom_wiener_process(md_env)) CALL create_wiener_process(md_env)
306 
307  time_iter_start = m_walltime()
308 
309  CALL get_md_env(md_env, force_env=force_env, itimes=itimes, constant=constant, &
310  md_ener=md_ener, t=time, used_time=used_time)
311 
312  ! Attach the time counter of the meta_env to the one of the MD
313  CALL set_meta_env(force_env%meta_env, time=time)
314 
315  ! Initialize the md_ener structure
316  CALL initialize_md_ener(md_ener, force_env, simpar)
317 
318  ! Check for ensembles requiring the stress tensor - takes into account the possibility for
319  ! multiple force_evals
320  IF ((simpar%ensemble == npt_i_ensemble) .OR. &
321  (simpar%ensemble == npt_ia_ensemble) .OR. &
322  (simpar%ensemble == npt_f_ensemble) .OR. &
323  (simpar%ensemble == npe_f_ensemble) .OR. &
324  (simpar%ensemble == npe_i_ensemble) .OR. &
325  (simpar%ensemble == nph_uniaxial_ensemble) .OR. &
326  (simpar%ensemble == nph_uniaxial_damped_ensemble)) THEN
327  check = virial%pv_availability
328  IF (.NOT. check) &
329  CALL cp_abort(__location__, &
330  "Virial evaluation not requested for this run in the input file!"// &
331  " You may consider to switch on the virial evaluation with the keyword: STRESS_TENSOR."// &
332  " Be sure the method you are using can compute the virial!")
333  IF (ASSOCIATED(force_env%sub_force_env)) THEN
334  DO i = 1, SIZE(force_env%sub_force_env)
335  IF (ASSOCIATED(force_env%sub_force_env(i)%force_env)) THEN
336  CALL force_env_get(force_env%sub_force_env(i)%force_env, subsys=subsys_i)
337  CALL cp_subsys_get(subsys_i, virial=virial)
338  check = check .AND. virial%pv_availability
339  END IF
340  END DO
341  END IF
342  IF (.NOT. check) &
343  CALL cp_abort(__location__, &
344  "Virial evaluation not requested for all the force_eval sections present in"// &
345  " the input file! You have to switch on the virial evaluation with the keyword: STRESS_TENSOR"// &
346  " in each force_eval section. Be sure the method you are using can compute the virial!")
347  END IF
348 
349  ! Computing Forces at zero MD step
350  IF (simpar%ensemble /= reftraj_ensemble) THEN
351  CALL section_vals_val_get(md_section, "STEP_START_VAL", i_val=itimes)
352  CALL section_vals_val_get(md_section, "TIME_START_VAL", r_val=time)
353  CALL section_vals_val_get(md_section, "ECONS_START_VAL", r_val=constant)
354  CALL cp_iterate(logger%iter_info, iter_nr=itimes)
355  IF (save_mem) THEN
356  work_section => section_vals_get_subs_vals(subsys_section, "VELOCITY")
357  CALL section_vals_remove_values(work_section)
358  work_section => section_vals_get_subs_vals(subsys_section, "SHELL_VELOCITY")
359  CALL section_vals_remove_values(work_section)
360  work_section => section_vals_get_subs_vals(subsys_section, "CORE_VELOCITY")
361  CALL section_vals_remove_values(work_section)
362  END IF
363 
364  IF (ehrenfest_md) THEN
365  CALL rt_prop_setup(force_env)
366  force_env%qs_env%rtp%dt = simpar%dt
367  ELSE
368  ![NB] Lets let all methods, even ones without consistent energies, succeed here.
369  ! They'll fail in actual integrator if needed
370  ! consistent_energies=.FALSE. by default
371  CALL force_env_calc_energy_force(force_env, calc_force=.true.)
372  END IF
373 
374  IF (ASSOCIATED(force_env%qs_env)) THEN
375  CALL qs_env_time_update(force_env%qs_env, time, itimes)
376 !deb force_env%qs_env%sim_time = time
377 !deb force_env%qs_env%sim_step = itimes
378  END IF
379  ! Warm-up engines for metadynamics
380  IF (ASSOCIATED(force_env%meta_env)) THEN
381  ! Setup stuff for plumed if needed
382  IF (force_env%meta_env%use_plumed .EQV. .true.) THEN
383  CALL metadyn_initialise_plumed(force_env, simpar, itimes)
384  ELSE
385  IF (force_env%meta_env%langevin) THEN
386  CALL create_wiener_process_cv(force_env%meta_env)
387  END IF
388  IF (force_env%meta_env%well_tempered) THEN
389  force_env%meta_env%wttemperature = simpar%temp_ext
390  IF (force_env%meta_env%wtgamma > epsilon(1._dp)) THEN
391  dummy = force_env%meta_env%wttemperature*(force_env%meta_env%wtgamma - 1._dp)
392  IF (force_env%meta_env%delta_t > epsilon(1._dp)) THEN
393  check = abs(force_env%meta_env%delta_t - dummy) < 1.e+3_dp*epsilon(1._dp)
394  IF (.NOT. check) &
395  CALL cp_abort(__location__, &
396  "Inconsistency between DELTA_T and WTGAMMA (both specified):"// &
397  " please, verify that DELTA_T=(WTGAMMA-1)*TEMPERATURE")
398  ELSE
399  force_env%meta_env%delta_t = dummy
400  END IF
401  ELSE
402  force_env%meta_env%wtgamma = 1._dp &
403  + force_env%meta_env%delta_t/force_env%meta_env%wttemperature
404  END IF
405  force_env%meta_env%invdt = 1._dp/force_env%meta_env%delta_t
406  END IF
407  CALL metadyn_forces(force_env)
408  CALL metadyn_write_colvar(force_env)
409  END IF
410  END IF
411 
412  IF (simpar%do_respa) THEN
413  CALL force_env_calc_energy_force(force_env%sub_force_env(1)%force_env, &
414  calc_force=.true.)
415  END IF
416 
417  CALL force_env_get(force_env, subsys=subsys)
418 
419  CALL cp_subsys_get(subsys, atomic_kinds=atomic_kinds, local_particles=local_particles, &
420  particles=particles, virial=virial)
421 
422  CALL virial_evaluate(atomic_kinds%els, particles%els, local_particles, &
423  virial, force_env%para_env)
424 
425  CALL md_energy(md_env, md_ener)
426  CALL md_write_output(md_env) !inits the print env at itimes == 0 also writes trajectories
427  md_stride = 1
428  ELSE
429  CALL get_md_env(md_env, reftraj=reftraj)
430  CALL initialize_reftraj(reftraj, reftraj_section, md_env)
431  itimes = reftraj%info%first_snapshot - 1
432  md_stride = reftraj%info%stride
433  IF (ASSOCIATED(force_env%meta_env)) THEN
434  IF (force_env%meta_env%use_plumed .EQV. .true.) THEN
435  CALL metadyn_initialise_plumed(force_env, simpar, itimes)
436  END IF
437  END IF
438  END IF
439 
440  CALL cp_print_key_finished_output(simpar%info_constraint, logger, &
441  constraint_section, "CONSTRAINT_INFO")
442  CALL cp_print_key_finished_output(simpar%lagrange_multipliers, logger, &
443  constraint_section, "LAGRANGE_MULTIPLIERS")
444 
445 ! if we need the initial kinetic energy for Hybrid Monte Carlo
446  IF (PRESENT(hmc_e_initial)) hmc_e_initial = md_ener%ekin
447 
448  IF (itimes >= simpar%max_steps) CALL cp_abort(__location__, &
449  "maximum step number smaller than initial step value")
450 
451  ! Real MD Loop
452  DO istep = 1, simpar%nsteps, md_stride
453  ! Increase counters
454  itimes = itimes + 1
455  time = time + simpar%dt
456  !needed when electric field fields are applied
457  IF (ASSOCIATED(force_env%qs_env)) THEN
458  CALL qs_env_time_update(force_env%qs_env, time, itimes)
459 !deb force_env%qs_env%sim_time = time
460 !deb force_env%qs_env%sim_step = itimes
461  END IF
462  IF (ehrenfest_md) force_env%qs_env%rtp%istep = istep
463 
464  IF (.NOT. logger%iter_info%last_iter(logger%iter_info%n_rlevel)) THEN
465  CALL cp_iterate(logger%iter_info, last=(istep == simpar%nsteps), iter_nr=itimes)
466  ELSE
467  CALL cp_iterate(logger%iter_info, last=.true., iter_nr=itimes)
468  END IF
469 
470  ! Open possible Shake output units
471  simpar%info_constraint = cp_print_key_unit_nr(logger, constraint_section, "CONSTRAINT_INFO", &
472  extension=".shakeLog", log_filename=.false.)
473  simpar%lagrange_multipliers = cp_print_key_unit_nr( &
474  logger, constraint_section, &
475  "LAGRANGE_MULTIPLIERS", extension=".LagrangeMultLog", log_filename=.false.)
476  simpar%dump_lm = btest(cp_print_key_should_output(logger%iter_info, constraint_section, &
477  "LAGRANGE_MULTIPLIERS"), cp_p_file)
478 
479  ! Velocity Verlet Integrator
480  CALL velocity_verlet(md_env, globenv)
481 
482  ! Close Shake output if requested...
483  CALL cp_print_key_finished_output(simpar%info_constraint, logger, &
484  constraint_section, "CONSTRAINT_INFO")
485  CALL cp_print_key_finished_output(simpar%lagrange_multipliers, logger, &
486  constraint_section, "LAGRANGE_MULTIPLIERS")
487 
488  ! Free Energy calculation
489  CALL free_energy_evaluate(md_env, should_stop, free_energy_section)
490 
491  IF (should_stop) EXIT
492 
493  ! Test for <PROJECT_NAME>.EXIT_MD or for WALL_TIME to exit
494  ! Default:
495  ! IF so we don't overwrite the restart or append to the trajectory
496  ! because the execution could in principle stop inside the SCF where energy
497  ! and forces are not converged.
498  ! But:
499  ! You can force to print the last step (for example if the method used
500  ! to compute energy and forces is not SCF based) activating the print_key
501  ! MOTION%MD%PRINT%FORCE_LAST.
502  CALL external_control(should_stop, "MD", globenv=globenv)
503 
504  !check if upper bound of total steps has been reached
505  IF (.NOT. (istep == simpar%nsteps) .AND. logger%iter_info%last_iter(logger%iter_info%n_rlevel)) should_stop = .true.
506  IF (itimes >= simpar%max_steps) should_stop = .true.
507 
508  ! call external hook e.g. from global optimization
509  IF (PRESENT(mdctrl)) &
510  CALL mdctrl_callback(mdctrl, md_env, should_stop)
511 
512  IF (should_stop) THEN
513  CALL cp_iterate(logger%iter_info, last=.true., iter_nr=itimes)
514  !In Ehrenfest molecular dynamics the external control is only checked after a converged propagation
515  !The restart needs to be written in order to be consistent with the mos/density matrix for the restart
516  IF (run_type_id == ehrenfest) THEN
517  CALL md_output(md_env, md_section, force_env%root_section, .false.)
518  ELSE
519  CALL md_output(md_env, md_section, force_env%root_section, should_stop)
520  END IF
521  EXIT
522  END IF
523 
524  IF (simpar%ensemble /= reftraj_ensemble) THEN
525  CALL md_energy(md_env, md_ener)
526  CALL temperature_control(simpar, md_env, md_ener, force_env, logger)
527  CALL comvel_control(md_ener, force_env, md_section, logger)
528  CALL angvel_control(md_ener, force_env, md_section, logger)
529  ELSE
530  CALL md_ener_reftraj(md_env, md_ener)
531  END IF
532 
533  time_iter_stop = m_walltime()
534  used_time = time_iter_stop - time_iter_start
535  time_iter_start = time_iter_stop
536 
537  CALL md_output(md_env, md_section, force_env%root_section, should_stop)
538  IF (simpar%ensemble == reftraj_ensemble) THEN
539  CALL write_output_reftraj(md_env)
540  END IF
541  END DO
542 
543 ! if we need the final kinetic energy for Hybrid Monte Carlo
544  IF (PRESENT(hmc_e_final)) hmc_e_final = md_ener%ekin
545 
546  ! Remove the iteration level
547  CALL cp_rm_iter_level(logger%iter_info, "MD")
548 
549  ! Clean up PLUMED
550  IF (ASSOCIATED(force_env%meta_env)) THEN
551  IF (force_env%meta_env%use_plumed .EQV. .true.) THEN
553  END IF
554  END IF
555 
556  ! Deallocate Thermostats and Barostats
557  CALL release_simpar_type(simpar)
558  CALL timestop(handle)
559 
560  END SUBROUTINE qs_mol_dyn_low
561 
562 END MODULE md_run
represent a simple array based list of the given type
Handles the type to compute averages during an MD.
Barostat structure: module containing barostat available for MD.
subroutine, public create_barostat_type(barostat, md_section, force_env, simpar, globenv)
...
Handles all functions related to the CELL.
Definition: cell_types.F:15
Routines to handle the external control of CP2K.
subroutine, public external_control(should_stop, flag, globenv, target_time, start_time, force_check)
External manipulations during a run : when the <PROJECT_NAME>.EXIT_$runtype command is sent the progr...
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,...
integer, parameter, public cp_p_file
subroutine, public cp_iterate(iteration_info, last, iter_nr, increment, iter_nr_out)
adds one to the actual iteration
subroutine, public cp_rm_iter_level(iteration_info, level_name, n_rlevel_att)
Removes an iteration level.
integer function, public cp_print_key_should_output(iteration_info, basis_section, print_key_path, used_print_key, first_time)
returns what should be done with the given property if btest(res,cp_p_store) then the property should...
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_calc_energy_force(force_env, calc_force, consistent_energies, skip_external_control, eval_energy_forces, require_consistent_energy_force, linres, calc_stress_tensor)
Interface routine for force and energy calculations.
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
Methods to perform free energy and free energy derivatives calculations.
subroutine, public free_energy_evaluate(md_env, converged, fe_section)
Main driver for free energy calculations In this routine we handle specifically biased MD.
defines types for metadynamics calculation
subroutine, public fe_env_create(fe_env, fe_section)
creates the fe_env
Define type storing the global information of a run. Keep the amount of stored data small....
Definition: global_types.F:21
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public ehrenfest
integer, parameter, public nph_uniaxial_ensemble
integer, parameter, public npt_i_ensemble
integer, parameter, public nph_uniaxial_damped_ensemble
integer, parameter, public npe_f_ensemble
integer, parameter, public langevin_ensemble
integer, parameter, public npe_i_ensemble
integer, parameter, public npt_ia_ensemble
integer, parameter, public npt_f_ensemble
integer, parameter, public reftraj_ensemble
checks the input and perform some automatic "magic" on it
subroutine, public remove_restart_info(input_file)
Removes section used to restart a calculation from an input file in memory.
objects that represent the structure of input sections and the data contained in an input section
subroutine, public section_vals_remove_values(section_vals)
removes the values of a repetition of the 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_val_get(section_vals, keyword_name, i_rep_section, i_rep_val, n_rep_val, val, l_val, i_val, r_val, c_val, l_vals, i_vals, r_vals, c_vals, explicit)
returns the requested value
Defines the basic variable types.
Definition: kinds.F:23
integer, parameter, public dp
Definition: kinds.F:34
integer, parameter, public default_string_length
Definition: kinds.F:57
Machine interface based on Fortran 2003 and POSIX.
Definition: machine.F:17
real(kind=dp) function, public m_walltime()
returns time from a real-time clock, protected against rolling early/easily
Definition: machine.F:123
Split md_ener module from md_environment_type.
Definition: md_ener_types.F:12
subroutine, public create_md_ener(md_ener)
retains the given md_ener structure
Definition: md_ener_types.F:66
prints all energy info per timestep to the screen or to user defined output files
Definition: md_energies.F:16
subroutine, public md_ener_reftraj(md_env, md_ener)
...
Definition: md_energies.F:223
subroutine, public initialize_md_ener(md_ener, force_env, simpar)
...
Definition: md_energies.F:115
subroutine, public md_energy(md_env, md_ener)
...
Definition: md_energies.F:177
subroutine, public md_write_output(md_env)
This routine computes the conserved quantity, temperature and things like that and prints them out.
Definition: md_energies.F:252
subroutine, public set_md_env(md_env, itimes, constant, cell, simpar, fe_env, force_env, para_env, init, first_time, thermostats, barostat, reftraj, md_ener, averages, thermal_regions, ehrenfest_md)
Set the integrator environment to the correct program.
subroutine, public md_env_create(md_env, md_section, para_env, force_env)
Creates MD environment Purpose: Initialise the integrator environment. retain the para_env for this e...
subroutine, public md_env_release(md_env)
releases the given md env
pure logical function, public need_per_atom_wiener_process(md_env)
...
subroutine, public get_md_env(md_env, itimes, constant, used_time, cell, simpar, npt, force_env, para_env, reftraj, t, init, first_time, fe_env, thermostats, barostat, thermostat_coeff, thermostat_part, thermostat_shell, thermostat_baro, thermostat_fast, thermostat_slow, md_ener, averages, thermal_regions, ehrenfest_md)
get components of MD environment type
Perform a molecular dynamics (MD) run using QUICKSTEP.
Definition: md_run.F:14
subroutine, public qs_mol_dyn(force_env, globenv, averages, rm_restart_info, hmc_e_initial, hmc_e_final, mdctrl)
Main driver module for Molecular Dynamics.
Definition: md_run.F:121
Utilities for Molecular Dynamics.
Definition: md_util.F:12
subroutine, public md_output(md_env, md_section, root_section, forced_io)
collects the part of the MD that, basically, does the output
Definition: md_util.F:53
Collection of utilities for setting-up and handle velocities in MD runs.
Definition: md_vel_utils.F:17
subroutine, public setup_velocities(force_env, simpar, globenv, md_env, md_section, constraint_section, write_binary_restart_file)
Initialize Velocities for MD runs.
subroutine, public angvel_control(md_ener, force_env, md_section, logger)
Set to 0 the angular velocity along MD runs, if required.
subroutine, public temperature_control(simpar, md_env, md_ener, force_env, logger)
Perform all temperature manipulations during a QS MD run.
subroutine, public comvel_control(md_ener, force_env, md_section, logger)
Set to 0 the velocity of the COM along MD runs, if required.
A common interface (wrapper) for a callback into the md_run loop. Currently this is only used by the ...
subroutine, public mdctrl_callback(mdctrl, md_env, should_stop)
This is called by md_run for each step during during its main-loop.
A common interface for passing a callback into the md_run loop.
Definition: mdctrl_types.F:13
Interface to the message passing library MPI.
defines types for metadynamics calculation
subroutine, public set_meta_env(meta_env, time)
sets the meta_env
Performs the metadynamics calculation.
Definition: metadynamics.F:14
subroutine, public metadyn_finalise_plumed()
makes sure PLUMED is shut down cleanly
Definition: metadynamics.F:200
subroutine, public metadyn_initialise_plumed(force_env, simpar, itimes)
...
Definition: metadynamics.F:129
subroutine, public metadyn_forces(force_env, vel)
add forces to the subsys due to the metadynamics run possibly modifies the velocites (if reflective w...
Definition: metadynamics.F:369
subroutine, public metadyn_write_colvar(force_env)
Write down COLVAR information evolved according to Vanden-Eijnden Ciccotti C.Phys....
Definition: metadynamics.F:730
represent a simple array based list of the given type
qs_environment methods that use many other modules
subroutine, public qs_env_time_update(qs_env, time, itimes)
...
initialization of the reftraj structure used to analyse previously generated trajectories
Definition: reftraj_types.F:15
subroutine, public create_reftraj(reftraj, reftraj_section, para_env)
...
Definition: reftraj_types.F:89
Initialize the analysis of trajectories to be done by activating the REFTRAJ ensemble.
Definition: reftraj_util.F:15
subroutine, public write_output_reftraj(md_env)
...
Definition: reftraj_util.F:493
subroutine, public initialize_reftraj(reftraj, reftraj_section, md_env)
...
Definition: reftraj_util.F:82
Routines for the real time propagation.
subroutine, public rt_prop_setup(force_env)
creates rtp_type, gets the initial state, either by reading MO's from file or calling SCF run
Methods for storing MD parameters type.
subroutine, public read_md_section(simpar, motion_section, md_section)
Reads the MD section and setup the simulation parameters type.
Type for storing MD parameters.
Definition: simpar_types.F:14
subroutine, public release_simpar_type(simpar)
Releases the simulation parameters type.
Definition: simpar_types.F:118
subroutine, public create_simpar_type(simpar)
Creates the simulation parameters type.
Definition: simpar_types.F:106
Thermal regions type: to initialize and control the temperature of different regions.
Setup of regions with different temperature.
subroutine, public print_thermal_regions_langevin(thermal_regions, simpar, pos, act)
print out information regarding to langevin regions defined in thermal_regions section
subroutine, public create_thermal_regions(thermal_regions, md_section, simpar, force_env)
create thermal_regions
Methods for Thermostats.
subroutine, public create_thermostats(thermostats, md_section, force_env, simpar, para_env, globenv, global_section)
...
Thermostat structure: module containing thermostat available for MD.
Provides an interface to the velocity-verlet based integrator routines for all ensembles.
subroutine, public velocity_verlet(md_env, globenv)
...
subroutine, public virial_evaluate(atomic_kind_set, particle_set, local_particles, virial, igroup)
Computes the kinetic part of the pressure tensor and updates the full VIRIAL (PV)
Handling of the Wiener process currently employed in turn of the Langevin dynamics.
subroutine, public create_wiener_process(md_env)
Create a Wiener process for Langevin dynamics and initialize an independent random number generator f...
subroutine, public create_wiener_process_cv(meta_env)
Create a Wiener process for Langevin dynamics used for metadynamics and initialize an independent ran...