(git:a1d4336)
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
25 cp_p_file,&
40 USE input_constants, ONLY: &
49 USE kinds, ONLY: default_string_length,&
50 dp
51 USE machine, ONLY: m_walltime
52 USE md_ener_types, ONLY: create_md_ener,&
56 md_energy,&
64 USE md_util, ONLY: md_output,&
66 USE md_vel_utils, ONLY: angvel_control,&
71 USE mdctrl_types, ONLY: mdctrl_type
80 USE reftraj_types, ONLY: create_reftraj,&
96 USE virial_types, ONLY: virial_type
99#include "../base/base_uses.f90"
100
101 IMPLICIT NONE
102
103 PRIVATE
104
105 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'md_run'
106
107 PUBLIC :: qs_mol_dyn
108
109CONTAINS
110
111! **************************************************************************************************
112!> \brief Main driver module for Molecular Dynamics
113!> \param force_env ...
114!> \param globenv ...
115!> \param averages ...
116!> \param rm_restart_info ...
117!> \param hmc_e_initial ...
118!> \param hmc_e_final ...
119!> \param mdctrl ...
120! **************************************************************************************************
121 SUBROUTINE qs_mol_dyn(force_env, globenv, averages, rm_restart_info, hmc_e_initial, hmc_e_final, mdctrl)
122
123 TYPE(force_env_type), POINTER :: force_env
124 TYPE(global_environment_type), POINTER :: globenv
125 TYPE(average_quantities_type), OPTIONAL, POINTER :: averages
126 LOGICAL, INTENT(IN), OPTIONAL :: rm_restart_info
127 REAL(kind=dp), OPTIONAL :: hmc_e_initial, hmc_e_final
128 TYPE(mdctrl_type), OPTIONAL, POINTER :: mdctrl
129
130 LOGICAL :: my_rm_restart_info
131 TYPE(md_environment_type), POINTER :: md_env
132 TYPE(mp_para_env_type), POINTER :: para_env
133 TYPE(section_vals_type), POINTER :: md_section, motion_section
134
135 my_rm_restart_info = .true.
136 IF (PRESENT(rm_restart_info)) my_rm_restart_info = rm_restart_info
137 NULLIFY (md_env, para_env)
138 para_env => force_env%para_env
139 motion_section => section_vals_get_subs_vals(force_env%root_section, "MOTION")
140 md_section => section_vals_get_subs_vals(motion_section, "MD")
141
142 ! Real call to MD driver - Low Level
143 ALLOCATE (md_env)
144 CALL md_env_create(md_env, md_section, para_env, force_env=force_env)
145 CALL set_md_env(md_env, averages=averages)
146 IF (PRESENT(hmc_e_initial) .AND. PRESENT(hmc_e_final)) THEN
147 CALL qs_mol_dyn_low(md_env, md_section, motion_section, force_env, globenv, &
148 hmc_e_initial=hmc_e_initial, hmc_e_final=hmc_e_final)
149 ELSE
150 CALL qs_mol_dyn_low(md_env, md_section, motion_section, force_env, globenv, mdctrl=mdctrl)
151 END IF
152 CALL md_env_release(md_env)
153 DEALLOCATE (md_env)
154
155 ! Clean restartable sections..
156 IF (my_rm_restart_info) CALL remove_restart_info(force_env%root_section)
157 END SUBROUTINE qs_mol_dyn
158
159! **************************************************************************************************
160!> \brief Purpose: Driver routine for MD run using QUICKSTEP.
161!> \param md_env ...
162!> \param md_section ...
163!> \param motion_section ...
164!> \param force_env ...
165!> \param globenv ...
166!> \param hmc_e_initial ...
167!> \param hmc_e_final ...
168!> \param mdctrl ...
169!> \par History
170!> - Cleaning (09.2007) Teodoro Laino [tlaino] - University of Zurich
171!> - Added lines to print out langevin regions (2014/02/04, LT)
172!> \author Creation (07.11.2002,MK)
173! **************************************************************************************************
174 SUBROUTINE qs_mol_dyn_low(md_env, md_section, motion_section, force_env, globenv, hmc_e_initial, hmc_e_final, mdctrl)
175
176 TYPE(md_environment_type), POINTER :: md_env
177 TYPE(section_vals_type), POINTER :: md_section, motion_section
178 TYPE(force_env_type), POINTER :: force_env
179 TYPE(global_environment_type), POINTER :: globenv
180 REAL(kind=dp), OPTIONAL :: hmc_e_initial, hmc_e_final
181 TYPE(mdctrl_type), OPTIONAL, POINTER :: mdctrl
182
183 CHARACTER(LEN=*), PARAMETER :: routinen = 'qs_mol_dyn_low'
184
185 CHARACTER(LEN=default_string_length) :: my_act, my_pos
186 INTEGER :: handle, i, istep, md_stride, run_type_id
187 INTEGER, POINTER :: itimes
188 LOGICAL :: check, ehrenfest_md, save_mem, &
189 should_stop, write_binary_restart_file
190 REAL(kind=dp) :: dummy, time_iter_start, time_iter_stop
191 REAL(kind=dp), POINTER :: constant, time, used_time
192 TYPE(atomic_kind_list_type), POINTER :: atomic_kinds
193 TYPE(barostat_type), POINTER :: barostat
194 TYPE(cell_type), POINTER :: cell
195 TYPE(cp_logger_type), POINTER :: logger
196 TYPE(cp_subsys_type), POINTER :: subsys, subsys_i
197 TYPE(distribution_1d_type), POINTER :: local_particles
198 TYPE(free_energy_type), POINTER :: fe_env
199 TYPE(md_ener_type), POINTER :: md_ener
200 TYPE(mp_para_env_type), POINTER :: para_env
201 TYPE(particle_list_type), POINTER :: particles
202 TYPE(reftraj_type), POINTER :: reftraj
203 TYPE(section_vals_type), POINTER :: constraint_section, force_env_section, &
204 free_energy_section, global_section, reftraj_section, subsys_section, work_section
205 TYPE(simpar_type), POINTER :: simpar
206 TYPE(thermal_regions_type), POINTER :: thermal_regions
207 TYPE(thermostats_type), POINTER :: thermostats
208 TYPE(virial_type), POINTER :: virial
209
210 CALL timeset(routinen, handle)
211 cpassert(ASSOCIATED(globenv))
212 cpassert(ASSOCIATED(force_env))
213
214 NULLIFY (particles, cell, simpar, itimes, used_time, subsys, &
215 md_ener, thermostats, barostat, reftraj, force_env_section, &
216 reftraj_section, work_section, atomic_kinds, &
217 local_particles, time, fe_env, free_energy_section, &
218 constraint_section, thermal_regions, virial, subsys_i)
219 logger => cp_get_default_logger()
220 para_env => force_env%para_env
221
222 global_section => section_vals_get_subs_vals(force_env%root_section, "GLOBAL")
223 free_energy_section => section_vals_get_subs_vals(motion_section, "FREE_ENERGY")
224 constraint_section => section_vals_get_subs_vals(motion_section, "CONSTRAINT")
225 CALL section_vals_val_get(global_section, "SAVE_MEM", l_val=save_mem)
226
227 CALL section_vals_val_get(global_section, "RUN_TYPE", i_val=run_type_id)
228 IF (run_type_id == ehrenfest) CALL set_md_env(md_env, ehrenfest_md=.true.)
229
230 CALL create_simpar_type(simpar)
231 force_env_section => force_env%force_env_section
232 subsys_section => section_vals_get_subs_vals(force_env_section, "SUBSYS")
233 CALL cp_add_iter_level(logger%iter_info, "MD")
234 CALL cp_iterate(logger%iter_info, iter_nr=0)
235 ! Read MD section
236 CALL read_md_section(simpar, motion_section, md_section)
237 ! Setup print_keys
238 simpar%info_constraint = cp_print_key_unit_nr(logger, constraint_section, &
239 "CONSTRAINT_INFO", extension=".shakeLog", log_filename=.false.)
240 simpar%lagrange_multipliers = cp_print_key_unit_nr(logger, constraint_section, &
241 "LAGRANGE_MULTIPLIERS", extension=".LagrangeMultLog", log_filename=.false.)
242 simpar%dump_lm = btest(cp_print_key_should_output(logger%iter_info, constraint_section, &
243 "LAGRANGE_MULTIPLIERS"), cp_p_file)
244
245 ! Create the structure for the md energies
246 ALLOCATE (md_ener)
247 CALL create_md_ener(md_ener)
248 CALL set_md_env(md_env, md_ener=md_ener)
249 NULLIFY (md_ener)
250
251 ! If requested setup Thermostats
252 CALL create_thermostats(thermostats, md_section, force_env, simpar, para_env, &
253 globenv, global_section)
254
255 ! If requested setup Barostat
256 CALL create_barostat_type(barostat, md_section, force_env, simpar, globenv)
257
258 ! If requested setup different thermal regions
259 CALL create_thermal_regions(thermal_regions, md_section, simpar, force_env)
260
261 ! If doing langevin_ensemble, then print out langevin_regions information upon request
262 IF (simpar%ensemble == langevin_ensemble) THEN
263 my_pos = "REWIND"
264 my_act = "WRITE"
265 CALL print_thermal_regions_langevin(thermal_regions, simpar, &
266 pos=my_pos, act=my_act)
267 END IF
268
269 CALL set_md_env(md_env, thermostats=thermostats, barostat=barostat, thermal_regions=thermal_regions)
270
271 CALL get_md_env(md_env, ehrenfest_md=ehrenfest_md)
272
273 !If requested set up the REFTRAJ run
274 IF (simpar%ensemble == reftraj_ensemble .AND. ehrenfest_md) &
275 cpabort("Ehrenfest MD does not support reftraj ensemble ")
276 IF (simpar%ensemble == reftraj_ensemble) THEN
277 reftraj_section => section_vals_get_subs_vals(md_section, "REFTRAJ")
278 ALLOCATE (reftraj)
279 CALL create_reftraj(reftraj, reftraj_section, para_env)
280 CALL set_md_env(md_env, reftraj=reftraj)
281 END IF
282
283 CALL force_env_get(force_env, subsys=subsys, cell=cell, &
284 force_env_section=force_env_section)
285 CALL cp_subsys_get(subsys, virial=virial)
286
287 ! Set V0 if needed
288 IF (simpar%ensemble == nph_uniaxial_ensemble .OR. simpar%ensemble == nph_uniaxial_damped_ensemble) THEN
289 IF (simpar%v0 == 0._dp) simpar%v0 = cell%deth
290 END IF
291
292 ! Initialize velocities possibly applying constraints at the zeroth MD step
293 CALL section_vals_val_get(motion_section, "PRINT%RESTART%SPLIT_RESTART_FILE", &
294 l_val=write_binary_restart_file)
295 CALL setup_velocities(force_env, simpar, globenv, md_env, md_section, constraint_section, &
296 write_binary_restart_file)
297
298 ! Setup Free Energy Calculation (if required)
299 CALL fe_env_create(fe_env, free_energy_section)
300
301 CALL set_md_env(md_env=md_env, simpar=simpar, fe_env=fe_env, cell=cell, &
302 force_env=force_env)
303
304 ! Possibly initialize Wiener processes
305 ![NB] Tested again within create_wiener_process. Why??
307
308 time_iter_start = m_walltime()
309
310 CALL get_md_env(md_env, force_env=force_env, itimes=itimes, constant=constant, &
311 md_ener=md_ener, t=time, used_time=used_time)
312
313 ! Attach the time counter of the meta_env to the one of the MD
314 CALL set_meta_env(force_env%meta_env, time=time)
315
316 ! Initialize the md_ener structure
317 CALL initialize_md_ener(md_ener, force_env, simpar)
318
319 ! Check for ensembles requiring the stress tensor - takes into account the possibility for
320 ! multiple force_evals
321 IF ((simpar%ensemble == npt_i_ensemble) .OR. &
322 (simpar%ensemble == npt_ia_ensemble) .OR. &
323 (simpar%ensemble == npt_f_ensemble) .OR. &
324 (simpar%ensemble == npe_f_ensemble) .OR. &
325 (simpar%ensemble == npe_i_ensemble) .OR. &
326 (simpar%ensemble == nph_uniaxial_ensemble) .OR. &
327 (simpar%ensemble == nph_uniaxial_damped_ensemble)) THEN
328 check = virial%pv_availability
329 IF (.NOT. check) &
330 CALL cp_abort(__location__, &
331 "Virial evaluation not requested for this run in the input file!"// &
332 " You may consider to switch on the virial evaluation with the keyword: STRESS_TENSOR."// &
333 " Be sure the method you are using can compute the virial!")
334 IF (ASSOCIATED(force_env%sub_force_env)) THEN
335 DO i = 1, SIZE(force_env%sub_force_env)
336 IF (ASSOCIATED(force_env%sub_force_env(i)%force_env)) THEN
337 CALL force_env_get(force_env%sub_force_env(i)%force_env, subsys=subsys_i)
338 CALL cp_subsys_get(subsys_i, virial=virial)
339 check = check .AND. virial%pv_availability
340 END IF
341 END DO
342 END IF
343 IF (.NOT. check) &
344 CALL cp_abort(__location__, &
345 "Virial evaluation not requested for all the force_eval sections present in"// &
346 " the input file! You have to switch on the virial evaluation with the keyword: STRESS_TENSOR"// &
347 " in each force_eval section. Be sure the method you are using can compute the virial!")
348 END IF
349
350 ! Computing Forces at zero MD step
351 IF (simpar%ensemble /= reftraj_ensemble) THEN
352 CALL section_vals_val_get(md_section, "STEP_START_VAL", i_val=itimes)
353 CALL section_vals_val_get(md_section, "TIME_START_VAL", r_val=time)
354 CALL section_vals_val_get(md_section, "ECONS_START_VAL", r_val=constant)
355 CALL cp_iterate(logger%iter_info, iter_nr=itimes)
356 IF (save_mem) THEN
357 work_section => section_vals_get_subs_vals(subsys_section, "VELOCITY")
358 CALL section_vals_remove_values(work_section)
359 work_section => section_vals_get_subs_vals(subsys_section, "SHELL_VELOCITY")
360 CALL section_vals_remove_values(work_section)
361 work_section => section_vals_get_subs_vals(subsys_section, "CORE_VELOCITY")
362 CALL section_vals_remove_values(work_section)
363 END IF
364
365 IF (ehrenfest_md) THEN
366 CALL rt_prop_setup(force_env)
367 force_env%qs_env%rtp%dt = simpar%dt
368 ELSE
369 ![NB] Lets let all methods, even ones without consistent energies, succeed here.
370 ! They'll fail in actual integrator if needed
371 ! consistent_energies=.FALSE. by default
372 CALL force_env_calc_energy_force(force_env, calc_force=.true.)
373 END IF
374
375 IF (ASSOCIATED(force_env%qs_env)) THEN
376 CALL qs_env_time_update(force_env%qs_env, time, itimes)
377 END IF
378 ! Warm-up engines for metadynamics
379 IF (ASSOCIATED(force_env%meta_env)) THEN
380 ! Setup stuff for plumed if needed
381 IF (force_env%meta_env%use_plumed .EQV. .true.) THEN
382 CALL metadyn_initialise_plumed(force_env, simpar, itimes)
383 ELSE
384 IF (force_env%meta_env%langevin) THEN
385 CALL create_wiener_process_cv(force_env%meta_env)
386 END IF
387 IF (force_env%meta_env%well_tempered) THEN
388 force_env%meta_env%wttemperature = simpar%temp_ext
389 IF (force_env%meta_env%wtgamma > epsilon(1._dp)) THEN
390 dummy = force_env%meta_env%wttemperature*(force_env%meta_env%wtgamma - 1._dp)
391 IF (force_env%meta_env%delta_t > epsilon(1._dp)) THEN
392 check = abs(force_env%meta_env%delta_t - dummy) < 1.e+3_dp*epsilon(1._dp)
393 IF (.NOT. check) &
394 CALL cp_abort(__location__, &
395 "Inconsistency between DELTA_T and WTGAMMA (both specified):"// &
396 " please, verify that DELTA_T=(WTGAMMA-1)*TEMPERATURE")
397 ELSE
398 force_env%meta_env%delta_t = dummy
399 END IF
400 ELSE
401 force_env%meta_env%wtgamma = 1._dp &
402 + force_env%meta_env%delta_t/force_env%meta_env%wttemperature
403 END IF
404 force_env%meta_env%invdt = 1._dp/force_env%meta_env%delta_t
405 END IF
406 CALL metadyn_forces(force_env)
407 CALL metadyn_write_colvar(force_env)
408 END IF
409 END IF
410
411 IF (simpar%do_respa) THEN
412 CALL force_env_calc_energy_force(force_env%sub_force_env(1)%force_env, &
413 calc_force=.true.)
414 END IF
415
416 CALL force_env_get(force_env, subsys=subsys)
417
418 CALL cp_subsys_get(subsys, atomic_kinds=atomic_kinds, local_particles=local_particles, &
419 particles=particles, virial=virial)
420
421 CALL virial_evaluate(atomic_kinds%els, particles%els, local_particles, &
422 virial, force_env%para_env)
423
424 CALL md_energy(md_env, md_ener)
425 CALL md_write_output(md_env) !inits the print env at itimes == 0 also writes trajectories
426 md_stride = 1
427 ELSE
428 CALL get_md_env(md_env, reftraj=reftraj)
429 CALL initialize_reftraj(reftraj, reftraj_section, md_env)
430 itimes = reftraj%info%first_snapshot - 1
431 md_stride = reftraj%info%stride
432 IF (ASSOCIATED(force_env%meta_env)) THEN
433 IF (force_env%meta_env%use_plumed .EQV. .true.) THEN
434 CALL metadyn_initialise_plumed(force_env, simpar, itimes)
435 END IF
436 END IF
437 END IF
438
439 CALL cp_print_key_finished_output(simpar%info_constraint, logger, &
440 constraint_section, "CONSTRAINT_INFO")
441 CALL cp_print_key_finished_output(simpar%lagrange_multipliers, logger, &
442 constraint_section, "LAGRANGE_MULTIPLIERS")
443
444! if we need the initial kinetic energy for Hybrid Monte Carlo
445 IF (PRESENT(hmc_e_initial)) hmc_e_initial = md_ener%ekin
446
447 IF (itimes >= simpar%max_steps) CALL cp_abort(__location__, &
448 "maximum step number smaller than initial step value")
449
450 ! Real MD Loop
451 DO istep = 1, simpar%nsteps, md_stride
452 ! Increase counters
453 itimes = itimes + 1
454 time = time + simpar%dt
455 !needed when electric field fields are applied
456 IF (ASSOCIATED(force_env%qs_env)) THEN
457 CALL qs_env_time_update(force_env%qs_env, time, itimes)
458 END IF
459 IF (ehrenfest_md) force_env%qs_env%rtp%istep = istep
460
461 IF (.NOT. logger%iter_info%last_iter(logger%iter_info%n_rlevel)) THEN
462 CALL cp_iterate(logger%iter_info, last=(istep == simpar%nsteps), iter_nr=itimes)
463 ELSE
464 CALL cp_iterate(logger%iter_info, last=.true., iter_nr=itimes)
465 END IF
466
467 ! Open possible Shake output units
468 simpar%info_constraint = cp_print_key_unit_nr(logger, constraint_section, "CONSTRAINT_INFO", &
469 extension=".shakeLog", log_filename=.false.)
470 simpar%lagrange_multipliers = cp_print_key_unit_nr( &
471 logger, constraint_section, &
472 "LAGRANGE_MULTIPLIERS", extension=".LagrangeMultLog", log_filename=.false.)
473 simpar%dump_lm = btest(cp_print_key_should_output(logger%iter_info, constraint_section, &
474 "LAGRANGE_MULTIPLIERS"), cp_p_file)
475
476 ! Update temperature for thermal regions and thermostat regions
477 CALL update_expected_temperature(md_env)
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
562END 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, 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....
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:153
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:122
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.