(git:4cf809f)
Loading...
Searching...
No Matches
md_run.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2026 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! **************************************************************************************************
14MODULE md_run
17 USE barostat_types, ONLY: barostat_type,&
19 USE cell_types, ONLY: cell_type
26 cp_p_file,&
42 USE input_constants, ONLY: &
51 USE kinds, ONLY: default_string_length,&
52 dp
53 USE machine, ONLY: m_flush,&
55 USE md_ener_types, ONLY: create_md_ener,&
59 md_energy,&
67 USE md_util, ONLY: md_output,&
69 USE md_vel_utils, ONLY: angvel_control,&
74 USE mdctrl_types, ONLY: mdctrl_type
83 USE reftraj_types, ONLY: create_reftraj,&
99 USE virial_types, ONLY: virial_type
102#include "../base/base_uses.f90"
103
104 IMPLICIT NONE
105
106 PRIVATE
107
108 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'md_run'
109
110 PUBLIC :: qs_mol_dyn
111
112CONTAINS
113
114! **************************************************************************************************
115!> \brief Main driver module for Molecular Dynamics
116!> \param force_env ...
117!> \param globenv ...
118!> \param averages ...
119!> \param rm_restart_info ...
120!> \param hmc_e_initial ...
121!> \param hmc_e_final ...
122!> \param mdctrl ...
123! **************************************************************************************************
124 SUBROUTINE qs_mol_dyn(force_env, globenv, averages, rm_restart_info, hmc_e_initial, hmc_e_final, mdctrl)
125
126 TYPE(force_env_type), POINTER :: force_env
127 TYPE(global_environment_type), POINTER :: globenv
128 TYPE(average_quantities_type), OPTIONAL, POINTER :: averages
129 LOGICAL, INTENT(IN), OPTIONAL :: rm_restart_info
130 REAL(kind=dp), OPTIONAL :: hmc_e_initial, hmc_e_final
131 TYPE(mdctrl_type), OPTIONAL, POINTER :: mdctrl
132
133 LOGICAL :: my_rm_restart_info
134 TYPE(md_environment_type), POINTER :: md_env
135 TYPE(mp_para_env_type), POINTER :: para_env
136 TYPE(section_vals_type), POINTER :: md_section, motion_section
137
138 my_rm_restart_info = .true.
139 IF (PRESENT(rm_restart_info)) my_rm_restart_info = rm_restart_info
140 NULLIFY (md_env, para_env)
141
142 ! Tell ACE that this is a dynamic run: Bypass C will use full HFX
143 ! for the entire first MD step so wavefunction propagation delivers
144 ! a near-converged C_occ to step 1, making the ACE BUILD accurate.
145 CALL hfx_ace_set_dynamic_mode(.true.)
146
147 para_env => force_env%para_env
148 motion_section => section_vals_get_subs_vals(force_env%root_section, "MOTION")
149 md_section => section_vals_get_subs_vals(motion_section, "MD")
150
151 ! Real call to MD driver - Low Level
152 ALLOCATE (md_env)
153 CALL md_env_create(md_env, md_section, para_env, force_env=force_env)
154 CALL set_md_env(md_env, averages=averages)
155 IF (PRESENT(hmc_e_initial) .AND. PRESENT(hmc_e_final)) THEN
156 CALL qs_mol_dyn_low(md_env, md_section, motion_section, force_env, globenv, &
157 hmc_e_initial=hmc_e_initial, hmc_e_final=hmc_e_final)
158 ELSE
159 CALL qs_mol_dyn_low(md_env, md_section, motion_section, force_env, globenv, mdctrl=mdctrl)
160 END IF
161 CALL md_env_release(md_env)
162 DEALLOCATE (md_env)
163
164 ! Clean restartable sections..
165 IF (my_rm_restart_info) CALL remove_restart_info(force_env%root_section)
166 END SUBROUTINE qs_mol_dyn
167
168! **************************************************************************************************
169!> \brief Purpose: Driver routine for MD run using QUICKSTEP.
170!> \param md_env ...
171!> \param md_section ...
172!> \param motion_section ...
173!> \param force_env ...
174!> \param globenv ...
175!> \param hmc_e_initial ...
176!> \param hmc_e_final ...
177!> \param mdctrl ...
178!> \par History
179!> - Cleaning (09.2007) Teodoro Laino [tlaino] - University of Zurich
180!> - Added lines to print out langevin regions (2014/02/04, LT)
181!> \author Creation (07.11.2002,MK)
182! **************************************************************************************************
183 SUBROUTINE qs_mol_dyn_low(md_env, md_section, motion_section, force_env, globenv, hmc_e_initial, hmc_e_final, mdctrl)
184
185 TYPE(md_environment_type), POINTER :: md_env
186 TYPE(section_vals_type), POINTER :: md_section, motion_section
187 TYPE(force_env_type), POINTER :: force_env
188 TYPE(global_environment_type), POINTER :: globenv
189 REAL(kind=dp), OPTIONAL :: hmc_e_initial, hmc_e_final
190 TYPE(mdctrl_type), OPTIONAL, POINTER :: mdctrl
191
192 CHARACTER(LEN=*), PARAMETER :: routinen = 'qs_mol_dyn_low'
193
194 CHARACTER(LEN=80) :: md_step_line
195 CHARACTER(LEN=default_string_length) :: my_act, my_pos
196 INTEGER :: handle, i, istep, md_stride, &
197 output_unit, run_type_id
198 INTEGER, POINTER :: itimes
199 LOGICAL :: check, ehrenfest_md, save_mem, &
200 should_stop, write_binary_restart_file
201 REAL(kind=dp) :: dummy, time_iter_start, time_iter_stop
202 REAL(kind=dp), POINTER :: constant, time, used_time
203 TYPE(atomic_kind_list_type), POINTER :: atomic_kinds
204 TYPE(barostat_type), POINTER :: barostat
205 TYPE(cell_type), POINTER :: cell
206 TYPE(cp_logger_type), POINTER :: logger
207 TYPE(cp_subsys_type), POINTER :: subsys, subsys_i
208 TYPE(distribution_1d_type), POINTER :: local_particles
209 TYPE(free_energy_type), POINTER :: fe_env
210 TYPE(md_ener_type), POINTER :: md_ener
211 TYPE(mp_para_env_type), POINTER :: para_env
212 TYPE(particle_list_type), POINTER :: particles
213 TYPE(reftraj_type), POINTER :: reftraj
214 TYPE(section_vals_type), POINTER :: constraint_section, force_env_section, &
215 free_energy_section, global_section, reftraj_section, subsys_section, work_section
216 TYPE(simpar_type), POINTER :: simpar
217 TYPE(thermal_regions_type), POINTER :: thermal_regions
218 TYPE(thermostats_type), POINTER :: thermostats
219 TYPE(virial_type), POINTER :: virial
220
221 CALL timeset(routinen, handle)
222 cpassert(ASSOCIATED(globenv))
223 cpassert(ASSOCIATED(force_env))
224
225 NULLIFY (particles, cell, simpar, itimes, used_time, subsys, &
226 md_ener, thermostats, barostat, reftraj, force_env_section, &
227 reftraj_section, work_section, atomic_kinds, &
228 local_particles, time, fe_env, free_energy_section, &
229 constraint_section, thermal_regions, virial, subsys_i)
230 logger => cp_get_default_logger()
231 para_env => force_env%para_env
232 output_unit = cp_logger_get_default_io_unit(logger)
233
234 global_section => section_vals_get_subs_vals(force_env%root_section, "GLOBAL")
235 free_energy_section => section_vals_get_subs_vals(motion_section, "FREE_ENERGY")
236 constraint_section => section_vals_get_subs_vals(motion_section, "CONSTRAINT")
237 CALL section_vals_val_get(global_section, "SAVE_MEM", l_val=save_mem)
238
239 CALL section_vals_val_get(global_section, "RUN_TYPE", i_val=run_type_id)
240 IF (run_type_id == ehrenfest) CALL set_md_env(md_env, ehrenfest_md=.true.)
241
242 CALL create_simpar_type(simpar)
243 force_env_section => force_env%force_env_section
244 subsys_section => section_vals_get_subs_vals(force_env_section, "SUBSYS")
245 CALL cp_add_iter_level(logger%iter_info, "MD")
246 CALL cp_iterate(logger%iter_info, iter_nr=0)
247 ! Read MD section
248 CALL read_md_section(simpar, motion_section, md_section)
249 ! Setup print_keys
250 simpar%info_constraint = cp_print_key_unit_nr(logger, constraint_section, &
251 "CONSTRAINT_INFO", extension=".shakeLog", log_filename=.false.)
252 simpar%lagrange_multipliers = cp_print_key_unit_nr(logger, constraint_section, &
253 "LAGRANGE_MULTIPLIERS", extension=".LagrangeMultLog", log_filename=.false.)
254 simpar%dump_lm = btest(cp_print_key_should_output(logger%iter_info, constraint_section, &
255 "LAGRANGE_MULTIPLIERS"), cp_p_file)
256
257 ! Create the structure for the md energies
258 ALLOCATE (md_ener)
259 CALL create_md_ener(md_ener)
260 CALL set_md_env(md_env, md_ener=md_ener)
261 NULLIFY (md_ener)
262
263 ! If requested setup Thermostats
264 CALL create_thermostats(thermostats, md_section, force_env, simpar, para_env, &
265 globenv, global_section)
266
267 ! If requested setup Barostat
268 CALL create_barostat_type(barostat, md_section, force_env, simpar, globenv)
269
270 ! If requested setup different thermal regions
271 CALL create_thermal_regions(thermal_regions, md_section, simpar, force_env)
272
273 ! If doing langevin_ensemble, then print out langevin_regions information upon request
274 IF (simpar%ensemble == langevin_ensemble) THEN
275 my_pos = "REWIND"
276 my_act = "WRITE"
277 CALL print_thermal_regions_langevin(thermal_regions, simpar, &
278 pos=my_pos, act=my_act)
279 END IF
280
281 CALL set_md_env(md_env, thermostats=thermostats, barostat=barostat, thermal_regions=thermal_regions)
282
283 CALL get_md_env(md_env, ehrenfest_md=ehrenfest_md)
284
285 !If requested set up the REFTRAJ run
286 IF (simpar%ensemble == reftraj_ensemble .AND. ehrenfest_md) &
287 cpabort("Ehrenfest MD does not support reftraj ensemble ")
288 IF (simpar%ensemble == reftraj_ensemble) THEN
289 reftraj_section => section_vals_get_subs_vals(md_section, "REFTRAJ")
290 ALLOCATE (reftraj)
291 CALL create_reftraj(reftraj, reftraj_section, para_env)
292 CALL set_md_env(md_env, reftraj=reftraj)
293 END IF
294
295 CALL force_env_get(force_env, subsys=subsys, cell=cell, &
296 force_env_section=force_env_section)
297 CALL cp_subsys_get(subsys, virial=virial)
298
299 ! Set V0 if needed
300 IF (simpar%ensemble == nph_uniaxial_ensemble .OR. simpar%ensemble == nph_uniaxial_damped_ensemble) THEN
301 IF (simpar%v0 == 0._dp) simpar%v0 = cell%deth
302 END IF
303
304 ! Initialize velocities possibly applying constraints at the zeroth MD step
305 CALL section_vals_val_get(motion_section, "PRINT%RESTART%SPLIT_RESTART_FILE", &
306 l_val=write_binary_restart_file)
307 CALL setup_velocities(force_env, simpar, globenv, md_env, md_section, constraint_section, &
308 write_binary_restart_file)
309
310 ! Setup Free Energy Calculation (if required)
311 CALL fe_env_create(fe_env, free_energy_section)
312
313 CALL set_md_env(md_env=md_env, simpar=simpar, fe_env=fe_env, cell=cell, &
314 force_env=force_env)
315
316 ! Possibly initialize Wiener processes
317 ![NB] Tested again within create_wiener_process. Why??
319
320 time_iter_start = m_walltime()
321
322 CALL get_md_env(md_env, force_env=force_env, itimes=itimes, constant=constant, &
323 md_ener=md_ener, t=time, used_time=used_time)
324
325 ! Attach the time counter of the meta_env to the one of the MD
326 CALL set_meta_env(force_env%meta_env, time=time)
327
328 ! Initialize the md_ener structure
329 CALL initialize_md_ener(md_ener, force_env, simpar)
330
331 ! Check for ensembles requiring the stress tensor - takes into account the possibility for
332 ! multiple force_evals
333 IF ((simpar%ensemble == npt_i_ensemble) .OR. &
334 (simpar%ensemble == npt_ia_ensemble) .OR. &
335 (simpar%ensemble == npt_f_ensemble) .OR. &
336 (simpar%ensemble == npe_f_ensemble) .OR. &
337 (simpar%ensemble == npe_i_ensemble) .OR. &
338 (simpar%ensemble == nph_uniaxial_ensemble) .OR. &
339 (simpar%ensemble == nph_uniaxial_damped_ensemble)) THEN
340 check = virial%pv_availability
341 IF (.NOT. check) &
342 CALL cp_abort(__location__, &
343 "Virial evaluation not requested for this run in the input file!"// &
344 " You may consider to switch on the virial evaluation with the keyword: STRESS_TENSOR."// &
345 " Be sure the method you are using can compute the virial!")
346 IF (ASSOCIATED(force_env%sub_force_env)) THEN
347 DO i = 1, SIZE(force_env%sub_force_env)
348 IF (ASSOCIATED(force_env%sub_force_env(i)%force_env)) THEN
349 CALL force_env_get(force_env%sub_force_env(i)%force_env, subsys=subsys_i)
350 CALL cp_subsys_get(subsys_i, virial=virial)
351 check = check .AND. virial%pv_availability
352 END IF
353 END DO
354 END IF
355 IF (.NOT. check) &
356 CALL cp_abort(__location__, &
357 "Virial evaluation not requested for all the force_eval sections present in"// &
358 " the input file! You have to switch on the virial evaluation with the keyword: STRESS_TENSOR"// &
359 " in each force_eval section. Be sure the method you are using can compute the virial!")
360 END IF
361
362 ! Computing Forces at zero MD step
363 IF (simpar%ensemble /= reftraj_ensemble) THEN
364 CALL section_vals_val_get(md_section, "STEP_START_VAL", i_val=itimes)
365 CALL section_vals_val_get(md_section, "TIME_START_VAL", r_val=time)
366 CALL section_vals_val_get(md_section, "ECONS_START_VAL", r_val=constant)
367 CALL cp_iterate(logger%iter_info, iter_nr=itimes)
368 IF (save_mem) THEN
369 work_section => section_vals_get_subs_vals(subsys_section, "VELOCITY")
370 CALL section_vals_remove_values(work_section)
371 work_section => section_vals_get_subs_vals(subsys_section, "SHELL_VELOCITY")
372 CALL section_vals_remove_values(work_section)
373 work_section => section_vals_get_subs_vals(subsys_section, "CORE_VELOCITY")
374 CALL section_vals_remove_values(work_section)
375 END IF
376
377 IF (ehrenfest_md) THEN
378 CALL rt_prop_setup(force_env)
379 force_env%qs_env%rtp%dt = simpar%dt
380 ELSE
381 ![NB] Lets let all methods, even ones without consistent energies, succeed here.
382 ! They'll fail in actual integrator if needed
383 ! consistent_energies=.FALSE. by default
384 CALL force_env_calc_energy_force(force_env, calc_force=.true.)
385 END IF
386
387 IF (ASSOCIATED(force_env%qs_env)) THEN
388 CALL qs_env_time_update(force_env%qs_env, time, itimes)
389 END IF
390 ! Warm-up engines for metadynamics
391 IF (ASSOCIATED(force_env%meta_env)) THEN
392 ! Setup stuff for plumed if needed
393 IF (force_env%meta_env%use_plumed .EQV. .true.) THEN
394 CALL metadyn_initialise_plumed(force_env, simpar, itimes)
395 ELSE
396 IF (force_env%meta_env%langevin) THEN
397 CALL create_wiener_process_cv(force_env%meta_env)
398 END IF
399 IF (force_env%meta_env%well_tempered) THEN
400 force_env%meta_env%wttemperature = simpar%temp_ext
401 IF (force_env%meta_env%wtgamma > epsilon(1._dp)) THEN
402 dummy = force_env%meta_env%wttemperature*(force_env%meta_env%wtgamma - 1._dp)
403 IF (force_env%meta_env%delta_t > epsilon(1._dp)) THEN
404 check = abs(force_env%meta_env%delta_t - dummy) < 1.e+3_dp*epsilon(1._dp)
405 IF (.NOT. check) &
406 CALL cp_abort(__location__, &
407 "Inconsistency between DELTA_T and WTGAMMA (both specified):"// &
408 " please, verify that DELTA_T=(WTGAMMA-1)*TEMPERATURE")
409 ELSE
410 force_env%meta_env%delta_t = dummy
411 END IF
412 ELSE
413 force_env%meta_env%wtgamma = 1._dp &
414 + force_env%meta_env%delta_t/force_env%meta_env%wttemperature
415 END IF
416 force_env%meta_env%invdt = 1._dp/force_env%meta_env%delta_t
417 END IF
418 CALL metadyn_forces(force_env)
419 CALL metadyn_write_colvar(force_env)
420 END IF
421 END IF
422
423 IF (simpar%do_respa) THEN
424 CALL force_env_calc_energy_force(force_env%sub_force_env(1)%force_env, &
425 calc_force=.true.)
426 END IF
427
428 CALL force_env_get(force_env, subsys=subsys)
429
430 CALL cp_subsys_get(subsys, atomic_kinds=atomic_kinds, local_particles=local_particles, &
431 particles=particles, virial=virial)
432
433 CALL virial_evaluate(atomic_kinds%els, particles%els, local_particles, &
434 virial, force_env%para_env)
435
436 CALL md_energy(md_env, md_ener)
437 CALL md_write_output(md_env) !inits the print env at itimes == 0 also writes trajectories
438 md_stride = 1
439 ELSE
440 CALL get_md_env(md_env, reftraj=reftraj)
441 CALL initialize_reftraj(reftraj, reftraj_section, md_env)
442 itimes = reftraj%info%first_snapshot - 1
443 md_stride = reftraj%info%stride
444 IF (ASSOCIATED(force_env%meta_env)) THEN
445 IF (force_env%meta_env%use_plumed .EQV. .true.) THEN
446 CALL metadyn_initialise_plumed(force_env, simpar, itimes)
447 END IF
448 END IF
449 END IF
450
451 CALL cp_print_key_finished_output(simpar%info_constraint, logger, &
452 constraint_section, "CONSTRAINT_INFO")
453 CALL cp_print_key_finished_output(simpar%lagrange_multipliers, logger, &
454 constraint_section, "LAGRANGE_MULTIPLIERS")
455
456! if we need the initial kinetic energy for Hybrid Monte Carlo
457 IF (PRESENT(hmc_e_initial)) hmc_e_initial = md_ener%ekin
458
459 IF (itimes >= simpar%max_steps) CALL cp_abort(__location__, &
460 "maximum step number smaller than initial step value")
461
462 ! Real MD Loop
463 DO istep = 1, simpar%nsteps, md_stride
464 IF (output_unit > 0) THEN
465 WRITE (md_step_line, fmt="(A,I12,A,I0)") "MD STEP: ", istep, " / ", simpar%nsteps
466
467 WRITE (unit=output_unit, fmt="(/,T2,A)") repeat("-", len_trim(md_step_line))
468 WRITE (unit=output_unit, fmt="(T2,A)") trim(md_step_line)
469 WRITE (unit=output_unit, fmt="(T2,A)") repeat("-", len_trim(md_step_line))
470 CALL m_flush(output_unit)
471 END IF
472 ! Increase counters
473 itimes = itimes + 1
474 time = time + simpar%dt
475 !needed when electric field fields are applied
476 IF (ASSOCIATED(force_env%qs_env)) THEN
477 CALL qs_env_time_update(force_env%qs_env, time, itimes)
478 END IF
479 IF (ehrenfest_md) force_env%qs_env%rtp%istep = istep
480
481 IF (.NOT. logger%iter_info%last_iter(logger%iter_info%n_rlevel)) THEN
482 CALL cp_iterate(logger%iter_info, last=(istep == simpar%nsteps), iter_nr=itimes)
483 ELSE
484 CALL cp_iterate(logger%iter_info, last=.true., iter_nr=itimes)
485 END IF
486
487 ! Open possible Shake output units
488 simpar%info_constraint = cp_print_key_unit_nr(logger, constraint_section, "CONSTRAINT_INFO", &
489 extension=".shakeLog", log_filename=.false.)
490 simpar%lagrange_multipliers = cp_print_key_unit_nr( &
491 logger, constraint_section, &
492 "LAGRANGE_MULTIPLIERS", extension=".LagrangeMultLog", log_filename=.false.)
493 simpar%dump_lm = btest(cp_print_key_should_output(logger%iter_info, constraint_section, &
494 "LAGRANGE_MULTIPLIERS"), cp_p_file)
495
496 ! Update temperature for thermal regions and thermostat regions
497 CALL update_expected_temperature(md_env)
498
499 ! Velocity Verlet Integrator
500 CALL velocity_verlet(md_env, globenv)
501
502 ! Close Shake output if requested...
503 CALL cp_print_key_finished_output(simpar%info_constraint, logger, &
504 constraint_section, "CONSTRAINT_INFO")
505 CALL cp_print_key_finished_output(simpar%lagrange_multipliers, logger, &
506 constraint_section, "LAGRANGE_MULTIPLIERS")
507
508 ! Free Energy calculation
509 CALL free_energy_evaluate(md_env, should_stop, free_energy_section)
510
511 IF (should_stop) EXIT
512
513 ! Test for <PROJECT_NAME>.EXIT_MD or for WALL_TIME to exit
514 ! Default:
515 ! IF so we don't overwrite the restart or append to the trajectory
516 ! because the execution could in principle stop inside the SCF where energy
517 ! and forces are not converged.
518 ! But:
519 ! You can force to print the last step (for example if the method used
520 ! to compute energy and forces is not SCF based) activating the print_key
521 ! MOTION%MD%PRINT%FORCE_LAST.
522 CALL external_control(should_stop, "MD", globenv=globenv)
523
524 !check if upper bound of total steps has been reached
525 IF (.NOT. (istep == simpar%nsteps) .AND. logger%iter_info%last_iter(logger%iter_info%n_rlevel)) should_stop = .true.
526 IF (itimes >= simpar%max_steps) should_stop = .true.
527
528 ! call external hook e.g. from global optimization
529 IF (PRESENT(mdctrl)) &
530 CALL mdctrl_callback(mdctrl, md_env, should_stop)
531
532 IF (should_stop) THEN
533 CALL cp_iterate(logger%iter_info, last=.true., iter_nr=itimes)
534 !In Ehrenfest molecular dynamics the external control is only checked after a converged propagation
535 !The restart needs to be written in order to be consistent with the mos/density matrix for the restart
536 IF (run_type_id == ehrenfest) THEN
537 CALL md_output(md_env, md_section, force_env%root_section, .false.)
538 ELSE
539 CALL md_output(md_env, md_section, force_env%root_section, should_stop)
540 END IF
541 EXIT
542 END IF
543
544 IF (simpar%ensemble /= reftraj_ensemble) THEN
545 CALL md_energy(md_env, md_ener)
546 CALL temperature_control(simpar, md_env, md_ener, force_env, logger)
547 CALL comvel_control(md_ener, force_env, md_section, logger)
548 CALL angvel_control(md_ener, force_env, md_section, logger)
549 ELSE
550 CALL md_ener_reftraj(md_env, md_ener)
551 END IF
552
553 time_iter_stop = m_walltime()
554 used_time = time_iter_stop - time_iter_start
555 time_iter_start = time_iter_stop
556
557 CALL md_output(md_env, md_section, force_env%root_section, should_stop)
558 IF (simpar%ensemble == reftraj_ensemble) THEN
559 CALL write_output_reftraj(md_env)
560 END IF
561 END DO
562
563! if we need the final kinetic energy for Hybrid Monte Carlo
564 IF (PRESENT(hmc_e_final)) hmc_e_final = md_ener%ekin
565
566 ! Remove the iteration level
567 CALL cp_rm_iter_level(logger%iter_info, "MD")
568
569 ! Clean up PLUMED
570 IF (ASSOCIATED(force_env%meta_env)) THEN
571 IF (force_env%meta_env%use_plumed .EQV. .true.) THEN
573 END IF
574 END IF
575
576 ! Deallocate Thermostats and Barostats
577 CALL release_simpar_type(simpar)
578 CALL timestop(handle)
579
580 END SUBROUTINE qs_mol_dyn_low
581
582END 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 ...
integer function, public cp_logger_get_default_io_unit(logger)
returns the unit nr for the ionode (-1 on all other processors) skips as well checks if the procs cal...
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, cell_ref, use_ref_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, ipi_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....
Adaptively Compressed Exchange (ACE) operator for HFX. Reference: Lin, J. Chem. Theory Comput....
subroutine, public hfx_ace_set_dynamic_mode(is_dynamic)
Mark this run as dynamic (GEO_OPT/MD) so Bypass C fires for geo step 0. Call this once from the geo_o...
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
subroutine, public m_flush(lunit)
flushes units if the &GLOBAL flag is set accordingly
Definition machine.F:124
real(kind=dp) function, public m_walltime()
returns time from a real-time clock, protected against rolling early/easily
Definition machine.F:141
Split md_ener module from md_environment_type.
subroutine, public create_md_ener(md_ener)
retains the given md_ener structure
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)
...
subroutine, public initialize_md_ener(md_ener, force_env, simpar)
...
subroutine, public md_energy(md_env, md_ener)
...
subroutine, public md_write_output(md_env)
This routine computes the conserved quantity, temperature and things like that and prints them out.
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:125
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:125
subroutine, public update_expected_temperature(md_env)
update_expected_temperature
Definition md_util.F:71
Collection of utilities for setting-up and handle velocities in MD runs.
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.
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.
subroutine, public metadyn_finalise_plumed()
makes sure PLUMED is shut down cleanly
subroutine, public metadyn_initialise_plumed(force_env, simpar, itimes)
...
subroutine, public metadyn_forces(force_env, vel)
add forces to the subsys due to the metadynamics run possibly modifies the velocites (if reflective w...
subroutine, public metadyn_write_colvar(force_env)
Write down COLVAR information evolved according to Vanden-Eijnden Ciccotti C.Phys....
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
subroutine, public create_reftraj(reftraj, reftraj_section, para_env)
...
Initialize the analysis of trajectories to be done by activating the REFTRAJ ensemble.
subroutine, public write_output_reftraj(md_env)
...
subroutine, public initialize_reftraj(reftraj, reftraj_section, md_env)
...
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.
subroutine, public release_simpar_type(simpar)
Releases the simulation parameters type.
subroutine, public create_simpar_type(simpar)
Creates the simulation parameters type.
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...
Type defining parameters related to the simulation cell.
Definition cell_types.F:60
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
stores all the informations relevant to an mpi environment
Simulation parameter type for molecular dynamics.