(git:ccc2433)
tamc_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 temperature accelarated hybrid monte carlo (TAHMC) run using QUICKSTEP
10 !> \par History
11 !> none
12 !> \author Alin M Elena
13 ! **************************************************************************************************
14 MODULE tamc_run
15 
16  USE atomic_kind_list_types, ONLY: atomic_kind_list_type
17  USE atomic_kind_types, ONLY: atomic_kind_type
18  USE averages_types, ONLY: average_quantities_type
19  USE barostat_types, ONLY: barostat_type,&
21  USE bibliography, ONLY: vandencic2006
22  USE cell_types, ONLY: cell_type
24  USE colvar_types, ONLY: hbp_colvar_id,&
25  wc_colvar_id,&
26  colvar_p_type
31  cp_logger_type
33  cp_iterate,&
34  cp_p_file,&
39  USE cp_subsys_types, ONLY: cp_subsys_get,&
40  cp_subsys_type
41  USE cp_units, ONLY: cp_unit_from_cp2k
42  USE distribution_1d_types, ONLY: distribution_1d_type
44  USE force_env_types, ONLY: force_env_get,&
45  force_env_type
46  USE free_energy_types, ONLY: fe_env_create,&
47  free_energy_type
48  USE global_types, ONLY: global_environment_type
49  USE input_constants, ONLY: &
57  section_vals_type,&
60  USE kinds, ONLY: dp
61  USE machine, ONLY: m_walltime
62  USE mc_environment_types, ONLY: get_mc_env,&
65  mc_environment_type,&
67  USE mc_misc, ONLY: mc_averages_create,&
69  USE mc_move_control, ONLY: init_mc_moves,&
71  USE mc_types, ONLY: get_mc_par,&
72  mc_averages_type,&
73  mc_ekin_type,&
74  mc_moves_type,&
75  mc_simpar_type,&
77  USE md_ener_types, ONLY: create_md_ener,&
78  md_ener_type
79  USE md_energies, ONLY: initialize_md_ener,&
80  md_energy
81  USE md_environment_types, ONLY: get_md_env,&
84  md_environment_type,&
86  USE md_run, ONLY: qs_mol_dyn
87  USE message_passing, ONLY: mp_comm_type,&
88  mp_para_env_type
89  USE metadynamics_types, ONLY: meta_env_type,&
90  metavar_type,&
92  USE molecule_kind_list_types, ONLY: molecule_kind_list_type
93  USE molecule_kind_types, ONLY: molecule_kind_type
94  USE molecule_list_types, ONLY: molecule_list_type
95  USE molecule_types, ONLY: global_constraint_type,&
96  molecule_type
97  USE parallel_rng_types, ONLY: uniform,&
98  rng_stream_type
99  USE particle_list_types, ONLY: particle_list_type
100  USE particle_types, ONLY: particle_type
101  USE physcon, ONLY: boltzmann,&
102  femtoseconds,&
103  joule,&
104  kelvin
109  USE reftraj_types, ONLY: create_reftraj,&
110  reftraj_type
113  USE simpar_types, ONLY: create_simpar_type,&
115  simpar_type
116  USE string_utilities, ONLY: str_comp
117  USE thermal_region_types, ONLY: thermal_regions_type
120  USE thermostat_types, ONLY: thermostats_type
122  USE virial_types, ONLY: virial_type
125 !!!!! monte carlo part
126 #include "../../base/base_uses.f90"
127 
128  IMPLICIT NONE
129 
130  PRIVATE
131 
132  CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'tamc_run'
133 
134  PUBLIC :: qs_tamc
135 
136 CONTAINS
137 
138 ! **************************************************************************************************
139 !> \brief Driver routine for TAHMC
140 !> \param force_env ...
141 !> \param globenv ...
142 !> \param averages ...
143 !> \author Alin M Elena
144 !> \note it computes the forces using QuickStep.
145 ! **************************************************************************************************
146  SUBROUTINE qs_tamc(force_env, globenv, averages)
147 
148  TYPE(force_env_type), POINTER :: force_env
149  TYPE(global_environment_type), POINTER :: globenv
150  TYPE(average_quantities_type), OPTIONAL, POINTER :: averages
151 
152  CHARACTER(LEN=*), PARAMETER :: routinen = 'qs_tamc'
153 
154  CHARACTER(LEN=20) :: ensemble
155  INTEGER :: handle, i, initialstep, iprint, isos, &
156  istep, j, md_stride, nmccycles, &
157  output_unit, rand2skip, run_type_id
158  INTEGER, POINTER :: itimes
159  LOGICAL :: check, explicit, my_rm_restart_info, &
160  save_mem, should_stop
161  REAL(kind=dp) :: auxrandom, inittime, rval, temp, &
162  time_iter_start, time_iter_stop
163  REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: an, fz, xieta, zbuff
164  REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: r
165  REAL(kind=dp), POINTER :: constant, time, used_time
166  TYPE(atomic_kind_list_type), POINTER :: atomic_kinds
167  TYPE(barostat_type), POINTER :: barostat
168  TYPE(cell_type), POINTER :: cell
169  TYPE(cp_logger_type), POINTER :: logger
170  TYPE(cp_subsys_type), POINTER :: subsys, subsys_i
171  TYPE(distribution_1d_type), POINTER :: local_particles
172  TYPE(free_energy_type), POINTER :: fe_env
173  TYPE(mc_averages_type), POINTER :: mcaverages
174  TYPE(mc_environment_type), POINTER :: mc_env
175  TYPE(mc_moves_type), POINTER :: gmoves, moves
176  TYPE(mc_simpar_type), POINTER :: mc_par
177  TYPE(md_ener_type), POINTER :: md_ener
178  TYPE(md_environment_type), POINTER :: md_env
179  TYPE(meta_env_type), POINTER :: meta_env_saved
180  TYPE(mp_para_env_type), POINTER :: para_env
181  TYPE(particle_list_type), POINTER :: particles
182  TYPE(reftraj_type), POINTER :: reftraj
183  TYPE(rng_stream_type) :: rng_stream_mc
184  TYPE(section_vals_type), POINTER :: constraint_section, force_env_section, &
185  free_energy_section, fs_section, global_section, mc_section, md_section, motion_section, &
186  reftraj_section, subsys_section, work_section
187  TYPE(simpar_type), POINTER :: simpar
188  TYPE(thermal_regions_type), POINTER :: thermal_regions
189  TYPE(thermostats_type), POINTER :: thermostats
190  TYPE(virial_type), POINTER :: virial
191 
192  initialstep = 0
193  inittime = 0.0_dp
194 
195  CALL timeset(routinen, handle)
196  my_rm_restart_info = .true.
197  NULLIFY (para_env, fs_section, virial)
198  para_env => force_env%para_env
199  motion_section => section_vals_get_subs_vals(force_env%root_section, "MOTION")
200  md_section => section_vals_get_subs_vals(motion_section, "MD")
201 
202  ! Real call to MD driver - Low Level
203  ALLOCATE (md_env)
204  CALL md_env_create(md_env, md_section, para_env, force_env=force_env)
205  IF (PRESENT(averages)) CALL set_md_env(md_env, averages=averages)
206 
207  cpassert(ASSOCIATED(globenv))
208  cpassert(ASSOCIATED(force_env))
209 
210  NULLIFY (particles, cell, simpar, itimes, used_time, subsys, &
211  md_ener, thermostats, barostat, reftraj, force_env_section, &
212  reftraj_section, work_section, atomic_kinds, &
213  local_particles, time, fe_env, free_energy_section, &
214  constraint_section, thermal_regions, subsys_i)
215  logger => cp_get_default_logger()
216  para_env => force_env%para_env
217 
218  global_section => section_vals_get_subs_vals(force_env%root_section, "GLOBAL")
219  free_energy_section => section_vals_get_subs_vals(motion_section, "FREE_ENERGY")
220  constraint_section => section_vals_get_subs_vals(motion_section, "CONSTRAINT")
221  CALL section_vals_val_get(global_section, "SAVE_MEM", l_val=save_mem)
222 
223  CALL section_vals_val_get(global_section, "RUN_TYPE", i_val=run_type_id)
224 
225  CALL create_simpar_type(simpar)
226  force_env_section => force_env%force_env_section
227  subsys_section => section_vals_get_subs_vals(force_env_section, "SUBSYS")
228  CALL cp_add_iter_level(logger%iter_info, "MD")
229  CALL cp_iterate(logger%iter_info, iter_nr=initialstep)
230  ! Read MD section
231  CALL read_md_section(simpar, motion_section, md_section)
232  ! Setup print_keys
233  simpar%info_constraint = cp_print_key_unit_nr(logger, constraint_section, &
234  "CONSTRAINT_INFO", extension=".shakeLog", log_filename=.false.)
235  simpar%lagrange_multipliers = cp_print_key_unit_nr(logger, constraint_section, &
236  "LAGRANGE_MULTIPLIERS", extension=".LagrangeMultLog", log_filename=.false.)
237  simpar%dump_lm = btest(cp_print_key_should_output(logger%iter_info, constraint_section, &
238  "LAGRANGE_MULTIPLIERS"), cp_p_file)
239 
240  ! Create the structure for the md energies
241  ALLOCATE (md_ener)
242  CALL create_md_ener(md_ener)
243  CALL set_md_env(md_env, md_ener=md_ener)
244 
245  ! If requested setup Thermostats
246  CALL create_thermostats(thermostats, md_section, force_env, simpar, para_env, &
247  globenv, global_section)
248 
249  ! If requested setup Barostat
250  CALL create_barostat_type(barostat, md_section, force_env, simpar, globenv)
251 
252  ! If requested setup different thermal regions
253  CALL create_thermal_regions(thermal_regions, md_section, simpar, force_env)
254 
255  CALL set_md_env(md_env, thermostats=thermostats, barostat=barostat, thermal_regions=thermal_regions)
256 
257  IF (simpar%ensemble == reftraj_ensemble) THEN
258  reftraj_section => section_vals_get_subs_vals(md_section, "REFTRAJ")
259  ALLOCATE (reftraj)
260  CALL create_reftraj(reftraj, reftraj_section, para_env)
261  CALL set_md_env(md_env, reftraj=reftraj)
262  END IF
263 
264  CALL force_env_get(force_env, subsys=subsys, cell=cell, &
265  force_env_section=force_env_section)
266 
267  ! Set V0 if needed
268  IF (simpar%ensemble == nph_uniaxial_ensemble .OR. simpar%ensemble == nph_uniaxial_damped_ensemble) THEN
269  IF (simpar%v0 == 0._dp) simpar%v0 = cell%deth
270  END IF
271 
272  ! Initialize velocities possibly applying constraints at the zeroth MD step
273 ! ! ! CALL section_vals_val_get(motion_section,"PRINT%RESTART%SPLIT_RESTART_FILE",&
274 ! ! ! l_val=write_binary_restart_file)
275 !! let us see if this created all the trouble
276 ! CALL setup_velocities(force_env,simpar,globenv,md_env,md_section,constraint_section, &
277 ! write_binary_restart_file)
278 
279  ! Setup Free Energy Calculation (if required)
280  CALL fe_env_create(fe_env, free_energy_section)
281  CALL set_md_env(md_env=md_env, simpar=simpar, fe_env=fe_env, cell=cell, &
282  force_env=force_env)
283 
284  ! Possibly initialize Wiener processes
285  IF (simpar%ensemble == langevin_ensemble) CALL create_wiener_process(md_env)
286  time_iter_start = m_walltime()
287 
288  CALL get_md_env(md_env, force_env=force_env, itimes=itimes, constant=constant, &
289  md_ener=md_ener, t=time, used_time=used_time)
290 
291  ! Attach the time counter of the meta_env to the one of the MD
292  CALL set_meta_env(force_env%meta_env, time=time)
293  ! Initialize the md_ener structure
294 
295  force_env%meta_env%dt = force_env%meta_env%zdt
296  CALL initialize_md_ener(md_ener, force_env, simpar)
297 ! force_env%meta_env%dt=force_env%meta_env%zdt
298 
299 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! MC setup up
300 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
301 
302  NULLIFY (mc_env, mc_par, mcaverages)
303 
304  CALL section_vals_get(force_env_section, n_repetition=isos)
305  cpassert(isos == 1)
306 ! set some values...will use get_globenv if that ever comes around
307 
308 ! initialize the random numbers
309 ! IF (para_env%is_source()) THEN
310  rng_stream_mc = rng_stream_type(name="Random numbers for monte carlo acc/rej", &
311  distribution_type=uniform)
312 ! ENDIF
313 !!!!! this should go in a routine hmc_read
314 
315  NULLIFY (mc_section)
316  ALLOCATE (mc_par)
317 
318  mc_section => section_vals_get_subs_vals(force_env%root_section, &
319  "MOTION%MC")
320  CALL section_vals_val_get(mc_section, "ENSEMBLE", &
321  c_val=ensemble)
322  cpassert(str_comp(ensemble, "TRADITIONAL"))
323  CALL section_vals_val_get(mc_section, "NSTEP", &
324  i_val=nmccycles)
325  cpassert(nmccycles > 0)
326  CALL section_vals_val_get(mc_section, "IPRINT", &
327  i_val=iprint)
328  CALL section_vals_val_get(mc_section, "RANDOMTOSKIP", i_val=rand2skip)
329  cpassert(rand2skip >= 0)
330  temp = cp_unit_from_cp2k(simpar%temp_ext, "K")
331 !
332 
333  CALL set_mc_par(mc_par, ensemble=ensemble, nstep=nmccycles, iprint=iprint, temperature=temp, &
334  beta=1.0_dp/temp/boltzmann*joule, exp_max_val=0.9_dp*log(huge(0.0_dp)), &
335  exp_min_val=0.9_dp*log(tiny(0.0_dp)), max_val=huge(0.0_dp), min_val=0.0_dp, &
336  source=para_env%source, group=para_env, ionode=para_env%is_source(), rand2skip=rand2skip)
337 
338  output_unit = cp_logger_get_default_io_unit(logger)
339  IF (output_unit > 0) THEN
340  WRITE (output_unit, '(a,a)') "HMC| Hybrid Monte Carlo Scheme "
341  WRITE (output_unit, '(a,a)') "HMC| Ensemble ", adjustl(ensemble)
342  WRITE (output_unit, '(a,i0)') "HMC| MC Cycles ", nmccycles
343  WRITE (output_unit, '(a,i0,a)') "HMC| Print every ", iprint, " cycles"
344  WRITE (output_unit, '(a,i0)') "HMC| Number of random numbers to skip ", rand2skip
345  WRITE (output_unit, '(a,f16.8,a)') "HMC| Temperature ", temp, "K"
346  END IF
347 
348  CALL force_env_get(force_env, subsys=subsys)
349 
350  CALL cp_subsys_get(subsys, atomic_kinds=atomic_kinds, local_particles=local_particles, &
351  particles=particles, virial=virial)
352 
353  DO i = 1, rand2skip
354  auxrandom = rng_stream_mc%next()
355  DO j = 1, 3*SIZE(particles%els)
356  auxrandom = globenv%gaussian_rng_stream%next()
357  END DO
358  END DO
359 
360  ALLOCATE (mc_env)
361  CALL mc_env_create(mc_env)
362  CALL set_mc_env(mc_env, mc_par=mc_par, force_env=force_env)
363 !!!!!!!end mc setup
364 
365  ! Check for ensembles requiring the stress tensor - takes into account the possibility for
366  ! multiple force_evals
367  IF ((simpar%ensemble == npt_i_ensemble) .OR. &
368  (simpar%ensemble == npt_ia_ensemble) .OR. &
369  (simpar%ensemble == npt_f_ensemble) .OR. &
370  (simpar%ensemble == npe_f_ensemble) .OR. &
371  (simpar%ensemble == npe_i_ensemble) .OR. &
372  (simpar%ensemble == nph_uniaxial_ensemble) .OR. &
373  (simpar%ensemble == nph_uniaxial_damped_ensemble)) THEN
374  check = virial%pv_availability
375  IF (.NOT. check) &
376  CALL cp_abort(__location__, &
377  "Virial evaluation not requested for this run in the input file! "// &
378  "You may consider to switch on the virial evaluation with the keyword: STRESS_TENSOR. "// &
379  "Be sure the method you are using can compute the virial!")
380  IF (ASSOCIATED(force_env%sub_force_env)) THEN
381  DO i = 1, SIZE(force_env%sub_force_env)
382  IF (ASSOCIATED(force_env%sub_force_env(i)%force_env)) THEN
383  CALL force_env_get(force_env%sub_force_env(i)%force_env, subsys=subsys_i)
384  CALL cp_subsys_get(subsys_i, virial=virial)
385  check = check .AND. virial%pv_availability
386  END IF
387  END DO
388  END IF
389  IF (.NOT. check) &
390  CALL cp_abort(__location__, &
391  "Virial evaluation not requested for all the force_eval sections present in"// &
392  " the input file! You have to switch on the virial evaluation with the keyword: STRESS_TENSOR"// &
393  " in each force_eval section. Be sure the method you are using can compute the virial!")
394  END IF
395 
396  ! Computing Forces at zero MD step
397  IF (simpar%ensemble /= reftraj_ensemble) THEN
398  CALL section_vals_val_get(md_section, "STEP_START_VAL", i_val=itimes)
399  CALL section_vals_val_get(md_section, "TIME_START_VAL", r_val=time)
400  CALL section_vals_val_get(md_section, "ECONS_START_VAL", r_val=constant)
401  CALL section_vals_val_set(md_section, "STEP_START_VAL", i_val=initialstep)
402  CALL section_vals_val_set(md_section, "TIME_START_VAL", r_val=inittime)
403  initialstep = itimes
404  CALL cp_iterate(logger%iter_info, iter_nr=itimes)
405  IF (save_mem) THEN
406  work_section => section_vals_get_subs_vals(subsys_section, "VELOCITY")
407  CALL section_vals_remove_values(work_section)
408  work_section => section_vals_get_subs_vals(subsys_section, "SHELL_VELOCITY")
409  CALL section_vals_remove_values(work_section)
410  work_section => section_vals_get_subs_vals(subsys_section, "CORE_VELOCITY")
411  CALL section_vals_remove_values(work_section)
412  END IF
413 
414 ! CALL force_env_calc_energy_force (force_env, calc_force=.TRUE.)
415  meta_env_saved => force_env%meta_env
416  NULLIFY (force_env%meta_env)
417  CALL force_env_calc_energy_force(force_env, calc_force=.false.)
418  force_env%meta_env => meta_env_saved
419 
420  IF (ASSOCIATED(force_env%qs_env)) THEN
421 ! force_env%qs_env%sim_time=time
422 ! force_env%qs_env%sim_step=itimes
423  force_env%qs_env%sim_time = 0.0_dp
424  force_env%qs_env%sim_step = 0
425  END IF
426  ! Warm-up engines for metadynamics
427  IF (ASSOCIATED(force_env%meta_env)) THEN
428  IF (force_env%meta_env%langevin) THEN
429  CALL create_wiener_process_cv(force_env%meta_env)
430  DO j = 1, (rand2skip - 1)/nmccycles
431  DO i = 1, force_env%meta_env%n_colvar
432  auxrandom = force_env%meta_env%rng(i)%next()
433  auxrandom = force_env%meta_env%rng(i)%next()
434  END DO
435  END DO
436  END IF
437 ! IF (force_env%meta_env%well_tempered) THEN
438 ! force_env%meta_env%wttemperature = simpar%temp_ext
439 ! IF (force_env%meta_env%wtgamma>EPSILON(1._dp)) THEN
440 ! dummy=force_env%meta_env%wttemperature*(force_env%meta_env%wtgamma-1._dp)
441 ! IF (force_env%meta_env%delta_t>EPSILON(1._dp)) THEN
442 ! check=ABS(force_env%meta_env%delta_t-dummy)<1.E+3_dp*EPSILON(1._dp)
443 ! IF(.NOT.check) CALL cp_abort(__LOCATION__,&
444 ! "Inconsistency between DELTA_T and WTGAMMA (both specified):"//&
445 ! " please, verify that DELTA_T=(WTGAMMA-1)*TEMPERATURE")
446 ! ELSE
447 ! force_env%meta_env%delta_t = dummy
448 ! ENDIF
449 ! ELSE
450 ! force_env%meta_env%wtgamma = 1._dp &
451 ! + force_env%meta_env%delta_t/force_env%meta_env%wttemperature
452 ! ENDIF
453 ! force_env%meta_env%invdt = 1._dp/force_env%meta_env%delta_t
454 ! ENDIF
455  CALL tamc_force(force_env)
456 ! CALL metadyn_write_colvar(force_env)
457  END IF
458 
459  IF (simpar%do_respa) THEN
460  CALL force_env_calc_energy_force(force_env%sub_force_env(1)%force_env, &
461  calc_force=.true.)
462  END IF
463 
464 ! CALL force_env_get( force_env, subsys=subsys)
465 !
466 ! CALL cp_subsys_get(subsys,atomic_kinds=atomic_kinds,local_particles=local_particles,&
467 ! particles=particles)
468 
469  CALL virial_evaluate(atomic_kinds%els, particles%els, local_particles, &
470  virial, force_env%para_env)
471 
472  CALL md_energy(md_env, md_ener)
473 ! CALL md_write_output(md_env) !inits the print env at itimes == 0 also writes trajectories
474  md_stride = 1
475  ELSE
476  CALL get_md_env(md_env, reftraj=reftraj)
477  CALL initialize_reftraj(reftraj, reftraj_section, md_env)
478  itimes = reftraj%info%first_snapshot - 1
479  md_stride = reftraj%info%stride
480  END IF
481 
482  CALL cp_print_key_finished_output(simpar%info_constraint, logger, &
483  constraint_section, "CONSTRAINT_INFO")
484  CALL cp_print_key_finished_output(simpar%lagrange_multipliers, logger, &
485  constraint_section, "LAGRANGE_MULTIPLIERS")
486  CALL init_mc_moves(moves)
487  CALL init_mc_moves(gmoves)
488  ALLOCATE (r(1:3, SIZE(particles%els)))
489 ! ALLOCATE (r_old(1:3,size(particles%els)))
490  CALL mc_averages_create(mcaverages)
491  !!!!! some more buffers
492  ! Allocate random number for Langevin Thermostat acting on COLVARS
493  ALLOCATE (xieta(2*force_env%meta_env%n_colvar))
494  xieta(:) = 0.0_dp
495  ALLOCATE (an(force_env%meta_env%n_colvar))
496  an(:) = 0.0_dp
497  ALLOCATE (fz(force_env%meta_env%n_colvar))
498  fz(:) = 0.0_dp
499  ALLOCATE (zbuff(2*force_env%meta_env%n_colvar))
500  zbuff(:) = 0.0_dp
501 
502  IF (output_unit > 0) THEN
503  WRITE (output_unit, '(a)') "HMC|==== Initial average forces"
504  END IF
505  CALL metadyn_write_colvar_header(force_env)
506  moves%hmc%attempts = 0
507  moves%hmc%successes = 0
508  gmoves%hmc%attempts = 0
509  gmoves%hmc%successes = 0
510  IF (initialstep == 0) THEN
511 !!! if we come from a restart we shall properly compute the average force
512 !!! read the average force up to now
513  DO i = 1, force_env%meta_env%n_colvar
514  fs_section => section_vals_get_subs_vals(force_env%meta_env%metadyn_section, "EXT_LAGRANGE_FS")
515  CALL section_vals_get(fs_section, explicit=explicit)
516  IF (explicit) THEN
517  CALL section_vals_val_get(fs_section, "_DEFAULT_KEYWORD_", &
518  i_rep_val=i, r_val=rval)
519  fz(i) = rval*rand2skip
520  END IF
521  END DO
522 
523  CALL hmcsampler(globenv, force_env, mcaverages, r, mc_par, moves, gmoves, rng_stream_mc, output_unit, &
524  fz, zbuff, nskip=rand2skip)
525  CALL cp_iterate(logger%iter_info, last=.false., iter_nr=0)
526  CALL section_vals_val_set(mc_section, "RANDOMTOSKIP", i_val=rand2skip + nmccycles)
527  CALL write_restart(md_env=md_env, root_section=force_env%root_section)
528  END IF
529  IF (output_unit > 0) THEN
530  WRITE (output_unit, '(a)') "HMC|==== end initial average forces"
531  END IF
532 ! call set_md_env(md_env, init=.FALSE.)
533 
534  CALL metadyn_write_colvar(force_env)
535 
536  DO istep = 1, force_env%meta_env%TAMCSteps
537  ! Increase counters
538  itimes = itimes + 1
539  time = time + force_env%meta_env%dt
540  IF (output_unit > 0) THEN
541  WRITE (output_unit, '(a)') "HMC|==================================="
542  WRITE (output_unit, '(a,1x,i0)') "HMC| on z step ", istep
543  END IF
544  !needed when electric field fields are applied
545  IF (ASSOCIATED(force_env%qs_env)) THEN
546  force_env%qs_env%sim_time = time
547  force_env%qs_env%sim_step = itimes
548  force_env%meta_env%time = force_env%qs_env%sim_time
549  END IF
550 
551  CALL cp_iterate(logger%iter_info, last=(istep == force_env%meta_env%TAMCSteps), iter_nr=itimes)
552  ! Open possible Shake output units
553  simpar%info_constraint = cp_print_key_unit_nr(logger, constraint_section, "CONSTRAINT_INFO", &
554  extension=".shakeLog", log_filename=.false.)
555  simpar%lagrange_multipliers = cp_print_key_unit_nr( &
556  logger, constraint_section, &
557  "LAGRANGE_MULTIPLIERS", extension=".LagrangeMultLog", log_filename=.false.)
558  simpar%dump_lm = btest(cp_print_key_should_output(logger%iter_info, constraint_section, &
559  "LAGRANGE_MULTIPLIERS"), cp_p_file)
560 
561  ! Velocity Verlet Integrator
562 
563  moves%hmc%attempts = 0
564  moves%hmc%successes = 0
565  CALL langevinvec(md_env, globenv, mc_env, moves, gmoves, r, &
566  rng_stream_mc, xieta, an, fz, mcaverages, zbuff)
567 
568  ! Close Shake output if requested...
569  CALL cp_print_key_finished_output(simpar%info_constraint, logger, &
570  constraint_section, "CONSTRAINT_INFO")
571  CALL cp_print_key_finished_output(simpar%lagrange_multipliers, logger, &
572  constraint_section, "LAGRANGE_MULTIPLIERS")
573  CALL cp_iterate(logger%iter_info, iter_nr=initialstep)
574  CALL metadyn_write_colvar(force_env)
575  ! Free Energy calculation
576 ! CALL free_energy_evaluate(md_env,should_stop,free_energy_section)
577 
578  ![AME:UB] IF (should_stop) EXIT
579 
580  ! Test for <PROJECT_NAME>.EXIT_MD or for WALL_TIME to exit
581  ! Default:
582  ! IF so we don't overwrite the restart or append to the trajectory
583  ! because the execution could in principle stop inside the SCF where energy
584  ! and forces are not converged.
585  ! But:
586  ! You can force to print the last step (for example if the method used
587  ! to compute energy and forces is not SCF based) activating the print_key
588  ! MOTION%MD%PRINT%FORCE_LAST.
589  CALL external_control(should_stop, "MD", globenv=globenv)
590  IF (should_stop) THEN
591  CALL cp_iterate(logger%iter_info, last=.true., iter_nr=itimes)
592 ! CALL md_output(md_env,md_section,force_env%root_section,should_stop)
593  EXIT
594  END IF
595 
596 ! IF(simpar%ensemble /= reftraj_ensemble) THEN
597 ! CALL md_energy(md_env, md_ener)
598 ! CALL temperature_control(simpar, md_env, md_ener, force_env, logger)
599 ! CALL comvel_control(md_ener, force_env, md_section, logger)
600 ! CALL angvel_control(md_ener, force_env, md_section, logger)
601 ! ELSE
602 ! CALL md_ener_reftraj(md_env, md_ener)
603 ! END IF
604 
605  time_iter_stop = m_walltime()
606  used_time = time_iter_stop - time_iter_start
607  time_iter_start = time_iter_stop
608 
609 !!!!! this writes the restart...
610 ! CALL md_output(md_env,md_section,force_env%root_section,should_stop)
611 
612 ! IF(simpar%ensemble == reftraj_ensemble ) THEN
613 ! CALL write_output_reftraj(md_env)
614 ! END IF
615 
616  IF (output_unit > 0) THEN
617  WRITE (output_unit, '(a,1x,i0)') "HMC| end z step ", istep
618  WRITE (output_unit, '(a)') "HMC|==================================="
619  END IF
620  END DO
621  CALL cp_iterate(logger%iter_info, last=.true., iter_nr=itimes)
622  force_env%qs_env%sim_time = 0.0_dp
623  force_env%qs_env%sim_step = 0
624  rand2skip = rand2skip + nmccycles*force_env%meta_env%TAMCSteps
625  IF (initialstep == 0) rand2skip = rand2skip + nmccycles
626  CALL section_vals_val_set(mc_section, "RANDOMTOSKIP", i_val=rand2skip)
627 
628  CALL write_restart(md_env=md_env, root_section=force_env%root_section)
629 ! if we need the final kinetic energy for Hybrid Monte Carlo
630 ! hmc_ekin%final_ekin=md_ener%ekin
631 
632  ! Remove the iteration level
633  CALL cp_rm_iter_level(logger%iter_info, "MD")
634 
635  ! Deallocate Thermostats and Barostats
636  CALL release_simpar_type(simpar)
637 
638  CALL md_env_release(md_env)
639  DEALLOCATE (md_env)
640  ! Clean restartable sections..
641  IF (my_rm_restart_info) CALL remove_restart_info(force_env%root_section)
642 ! IF (para_env%is_source()) THEN
643 ! ENDIF
644  CALL mc_env_release(mc_env)
645  DEALLOCATE (mc_env)
646  DEALLOCATE (mc_par)
647  CALL mc_moves_release(moves)
648  CALL mc_moves_release(gmoves)
649  DEALLOCATE (r)
650 ! DEALLOCATE(r_old)
651  DEALLOCATE (xieta)
652  DEALLOCATE (an)
653  DEALLOCATE (fz)
654  DEALLOCATE (zbuff)
655  CALL mc_averages_release(mcaverages)
656  CALL timestop(handle)
657 
658  END SUBROUTINE qs_tamc
659 
660 ! **************************************************************************************************
661 !> \brief Propagates velocities for z half a step
662 !> \param force_env ...
663 !> \param An ...
664 !> \author Alin M Elena
665 !> \note Vanden-Eijnden Ciccotti C.Phys.Letter 429 (2006) 310-316
666 ! **************************************************************************************************
667  SUBROUTINE tamc_velocities_colvar(force_env, An)
668  TYPE(force_env_type), POINTER :: force_env
669  REAL(kind=dp), DIMENSION(:), INTENT(INOUT) :: an
670 
671  CHARACTER(len=*), PARAMETER :: routinen = 'tamc_velocities_colvar'
672 
673  INTEGER :: handle, i_c
674  REAL(kind=dp) :: dt, fft, sigma
675  TYPE(cp_logger_type), POINTER :: logger
676  TYPE(meta_env_type), POINTER :: meta_env
677  TYPE(metavar_type), POINTER :: cv
678 
679  NULLIFY (logger, meta_env, cv)
680  meta_env => force_env%meta_env
681  CALL timeset(routinen, handle)
682  logger => cp_get_default_logger()
683  ! Add citation
684  IF (meta_env%langevin) CALL cite_reference(vandencic2006)
685  dt = meta_env%dt
686 
687  ! Evolve Velocities
688  meta_env%epot_walls = 0.0_dp
689  DO i_c = 1, meta_env%n_colvar
690  cv => meta_env%metavar(i_c)
691  fft = cv%ff_s + cv%ff_hills
692  sigma = sqrt((meta_env%temp_wanted*kelvin)*2.0_dp*(boltzmann/joule)*cv%gamma/cv%mass)
693  cv%vvp = cv%vvp + 0.5_dp*dt*(fft/cv%mass - cv%gamma*cv%vvp)*(1.0_dp - 0.25_dp*dt*cv%gamma) + an(i_c)
694  meta_env%epot_walls = meta_env%epot_walls + cv%epot_walls
695  END DO
696  CALL timestop(handle)
697  END SUBROUTINE tamc_velocities_colvar
698 
699 ! **************************************************************************************************
700 !> \brief propagates z one step
701 !> \param force_env ...
702 !> \param xieta ...
703 !> \author Alin M Elena
704 !> \note Vanden-Eijnden Ciccotti C.Phys.Letter 429 (2006) 310-316
705 ! **************************************************************************************************
706  SUBROUTINE tamc_position_colvar(force_env, xieta)
707  TYPE(force_env_type), POINTER :: force_env
708  REAL(kind=dp), DIMENSION(:), INTENT(INOUT) :: xieta
709 
710  CHARACTER(len=*), PARAMETER :: routinen = 'tamc_position_colvar'
711 
712  INTEGER :: handle, i_c
713  REAL(kind=dp) :: dt, sigma
714  TYPE(cp_logger_type), POINTER :: logger
715  TYPE(meta_env_type), POINTER :: meta_env
716  TYPE(metavar_type), POINTER :: cv
717 
718  NULLIFY (logger, meta_env, cv)
719  meta_env => force_env%meta_env
720 ! IF (.NOT.ASSOCIATED(meta_env)) RETURN
721 
722  CALL timeset(routinen, handle)
723  logger => cp_get_default_logger()
724 
725  ! Add citation
726  IF (meta_env%langevin) CALL cite_reference(vandencic2006)
727  dt = meta_env%dt
728 
729  ! Update of ss0
730  DO i_c = 1, meta_env%n_colvar
731  cv => meta_env%metavar(i_c)
732  sigma = sqrt((meta_env%temp_wanted*kelvin)*2.0_dp*(boltzmann/joule)*cv%gamma/cv%mass)
733 ! cv%ss0 =cv%ss0 +dt*cv%vvp
734  cv%ss0 = cv%ss0 + dt*cv%vvp + dt*sqrt(dt/12.0_dp)*sigma*xieta(i_c + meta_env%n_colvar)
735  IF (cv%periodic) THEN
736  ! A periodic COLVAR is always within [-pi,pi]
737  cv%ss0 = sign(1.0_dp, asin(sin(cv%ss0)))*acos(cos(cv%ss0))
738  END IF
739  END DO
740  CALL timestop(handle)
741 
742  END SUBROUTINE tamc_position_colvar
743 
744 ! **************************************************************************************************
745 !> \brief Computes forces on z
746 !> #details also can be used to get the potenzial evergy of z
747 !> \param force_env ...
748 !> \param zpot ...
749 !> \author Alin M Elena
750 ! **************************************************************************************************
751  SUBROUTINE tamc_force(force_env, zpot)
752  TYPE(force_env_type), POINTER :: force_env
753  REAL(kind=dp), INTENT(inout), OPTIONAL :: zpot
754 
755  CHARACTER(len=*), PARAMETER :: routinen = 'tamc_force'
756 
757  INTEGER :: handle, i, i_c, icolvar, ii
758  LOGICAL :: explicit
759  REAL(kind=dp) :: diff_ss, dt, rval
760  TYPE(colvar_p_type), DIMENSION(:), POINTER :: colvar_p
761  TYPE(cp_logger_type), POINTER :: logger
762  TYPE(cp_subsys_type), POINTER :: subsys
763  TYPE(meta_env_type), POINTER :: meta_env
764  TYPE(metavar_type), POINTER :: cv
765  TYPE(particle_list_type), POINTER :: particles
766  TYPE(section_vals_type), POINTER :: ss0_section, ss_section, vvp_section
767 
768  NULLIFY (logger, meta_env)
769  meta_env => force_env%meta_env
770 ! IF (.NOT.ASSOCIATED(meta_env)) RETURN
771 
772  CALL timeset(routinen, handle)
773  logger => cp_get_default_logger()
774  NULLIFY (colvar_p, subsys, cv, ss0_section, vvp_section, ss_section)
775  CALL force_env_get(force_env, subsys=subsys)
776 
777  dt = meta_env%dt
778  IF (.NOT. meta_env%restart) meta_env%n_steps = meta_env%n_steps + 1
779  ! compute ss and the derivative of ss with respect to the atomic positions
780  DO i_c = 1, meta_env%n_colvar
781  cv => meta_env%metavar(i_c)
782  icolvar = cv%icolvar
783  CALL colvar_eval_glob_f(icolvar, force_env)
784  cv%ss = subsys%colvar_p(icolvar)%colvar%ss
785  ! Restart for Extended Lagrangian Metadynamics
786  IF (meta_env%restart) THEN
787  ! Initialize the position of the collective variable in the extended lagrange
788  ss0_section => section_vals_get_subs_vals(meta_env%metadyn_section, "EXT_LAGRANGE_SS0")
789  CALL section_vals_get(ss0_section, explicit=explicit)
790  IF (explicit) THEN
791  CALL section_vals_val_get(ss0_section, "_DEFAULT_KEYWORD_", &
792  i_rep_val=i_c, r_val=rval)
793  cv%ss0 = rval
794  ELSE
795  cv%ss0 = cv%ss
796  END IF
797  vvp_section => section_vals_get_subs_vals(meta_env%metadyn_section, "EXT_LAGRANGE_VVP")
798  CALL section_vals_get(vvp_section, explicit=explicit)
799  IF (explicit) THEN
800  CALL setup_velocities_z(force_env)
801  CALL section_vals_val_get(vvp_section, "_DEFAULT_KEYWORD_", &
802  i_rep_val=i_c, r_val=rval)
803  cv%vvp = rval
804  ELSE
805  CALL setup_velocities_z(force_env)
806  END IF
807  ss_section => section_vals_get_subs_vals(meta_env%metadyn_section, "EXT_LAGRANGE_SS")
808  CALL section_vals_get(ss_section, explicit=explicit)
809  IF (explicit) THEN
810  CALL section_vals_val_get(ss_section, "_DEFAULT_KEYWORD_", &
811  i_rep_val=i_c, r_val=rval)
812  cv%ss = rval
813  END IF
814  END IF
815  !
816  END DO
817  ! forces on the atoms
818  NULLIFY (particles)
819  CALL cp_subsys_get(subsys, colvar_p=colvar_p, &
820  particles=particles)
821 
822  meta_env%restart = .false.
823  meta_env%epot_s = 0.0_dp
824  meta_env%epot_walls = 0.0_dp
825  DO i_c = 1, meta_env%n_colvar
826  cv => meta_env%metavar(i_c)
827  diff_ss = cv%ss - cv%ss0
828  IF (cv%periodic) THEN
829  ! The difference of a periodic COLVAR is always within [-pi,pi]
830  diff_ss = sign(1.0_dp, asin(sin(diff_ss)))*acos(cos(diff_ss))
831  END IF
832  cv%epot_s = 0.5_dp*cv%lambda*diff_ss*diff_ss
833  cv%ff_s = cv%lambda*(diff_ss)
834  meta_env%epot_s = meta_env%epot_s + cv%epot_s
835  icolvar = cv%icolvar
836 
837  DO ii = 1, colvar_p(icolvar)%colvar%n_atom_s
838  i = colvar_p(icolvar)%colvar%i_atom(ii)
839  particles%els(i)%f = particles%els(i)%f - cv%ff_s*colvar_p(icolvar)%colvar%dsdr(:, ii)
840  END DO
841 
842  END DO
843  IF (PRESENT(zpot)) zpot = meta_env%epot_s
844  CALL fix_atom_control(force_env)
845 
846  CALL timestop(handle)
847  END SUBROUTINE tamc_force
848 
849 ! **************************************************************************************************
850 !> \brief propagates one time step both z systems and samples the x system
851 !> \param md_env ...
852 !> \param globenv ...
853 !> \param mc_env ...
854 !> \param moves ...
855 !> \param gmoves ...
856 !> \param r ...
857 !> \param rng_stream_mc ...
858 !> \param xieta ...
859 !> \param An ...
860 !> \param fz ...
861 !> \param averages ...
862 !> \param zbuff ...
863 !> \author Alin M Elena
864 ! **************************************************************************************************
865  SUBROUTINE langevinvec(md_env, globenv, mc_env, moves, gmoves, r, &
866  rng_stream_mc, xieta, An, fz, averages, zbuff)
867 
868  TYPE(md_environment_type), POINTER :: md_env
869  TYPE(global_environment_type), POINTER :: globenv
870  TYPE(mc_environment_type), POINTER :: mc_env
871  TYPE(mc_moves_type), POINTER :: moves, gmoves
872  REAL(kind=dp), DIMENSION(:, :), INTENT(INOUT) :: r
873  TYPE(rng_stream_type), INTENT(INOUT) :: rng_stream_mc
874  REAL(kind=dp), DIMENSION(:), INTENT(INOUT) :: xieta, an, fz
875  TYPE(mc_averages_type), INTENT(INOUT), POINTER :: averages
876  REAL(kind=dp), DIMENSION(:), INTENT(INOUT) :: zbuff
877 
878  INTEGER :: iprint, ivar, nparticle, nparticle_kind, &
879  nstep, output_unit
880  REAL(kind=dp) :: dt, gamma, mass, sigma
881  TYPE(atomic_kind_list_type), POINTER :: atomic_kinds
882  TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
883  TYPE(cell_type), POINTER :: cell
884  TYPE(cp_logger_type), POINTER :: logger
885  TYPE(cp_subsys_type), POINTER :: subsys
886  TYPE(distribution_1d_type), POINTER :: local_molecules, local_particles
887  TYPE(force_env_type), POINTER :: force_env
888  TYPE(global_constraint_type), POINTER :: gci
889  TYPE(mc_simpar_type), POINTER :: mc_par
890  TYPE(molecule_kind_list_type), POINTER :: molecule_kinds
891  TYPE(molecule_kind_type), DIMENSION(:), POINTER :: molecule_kind_set
892  TYPE(molecule_list_type), POINTER :: molecules
893  TYPE(molecule_type), DIMENSION(:), POINTER :: molecule_set
894  TYPE(mp_para_env_type), POINTER :: para_env
895  TYPE(particle_list_type), POINTER :: particles
896  TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
897  TYPE(simpar_type), POINTER :: simpar
898  TYPE(virial_type), POINTER :: virial
899 
900  NULLIFY (logger, mc_par)
901  logger => cp_get_default_logger()
902  output_unit = cp_logger_get_default_io_unit(logger)
903 
904 ! quantitites to be nullified for the get_md_env
905  NULLIFY (simpar, force_env, para_env)
906 ! quantities to be nullified for the force_env_get environment
907  NULLIFY (subsys, cell)
908 ! quantitites to be nullified for the cp_subsys_get
909  NULLIFY (atomic_kinds, local_particles, particles, local_molecules, molecules, molecule_kinds, gci)
910 
911  CALL get_md_env(md_env=md_env, simpar=simpar, force_env=force_env, &
912  para_env=para_env)
913  CALL get_mc_env(mc_env, mc_par=mc_par)
914  CALL get_mc_par(mc_par, nstep=nstep, iprint=iprint)
915 
916  dt = simpar%dt
917  CALL force_env_get(force_env=force_env, subsys=subsys, cell=cell)
918 
919 !!!! this bit should vanish once I understand what the hell is with it
920 
921 ! ! Do some checks on coordinates and box
922  CALL apply_qmmm_walls_reflective(force_env)
923 
924  CALL cp_subsys_get(subsys=subsys, atomic_kinds=atomic_kinds, local_particles=local_particles, &
925  particles=particles, local_molecules=local_molecules, molecules=molecules, &
926  molecule_kinds=molecule_kinds, gci=gci, virial=virial)
927 
928  nparticle_kind = atomic_kinds%n_els
929  atomic_kind_set => atomic_kinds%els
930  molecule_kind_set => molecule_kinds%els
931 
932  nparticle = particles%n_els
933  particle_set => particles%els
934  molecule_set => molecules%els
935  cpassert(ASSOCIATED(force_env%meta_env))
936  cpassert(force_env%meta_env%langevin)
937  ! *** Velocity Verlet for Langevin *** v(t)--> v(t+1/2)
938  !!!!!! noise xi is in the first half, eta in the second half
939  DO ivar = 1, force_env%meta_env%n_colvar
940  xieta(ivar) = force_env%meta_env%rng(ivar)%next()
941  xieta(ivar + force_env%meta_env%n_colvar) = force_env%meta_env%rng(ivar)%next()
942  gamma = force_env%meta_env%metavar(ivar)%gamma
943  mass = force_env%meta_env%metavar(ivar)%mass
944  sigma = sqrt((force_env%meta_env%temp_wanted*kelvin)*2.0_dp*(boltzmann/joule)*gamma/mass)
945  an(ivar) = 0.5_dp*sqrt(dt)*sigma*(xieta(ivar)*(1.0_dp - 0.5_dp*dt*gamma) - &
946  dt*gamma*xieta(ivar + force_env%meta_env%n_colvar)/sqrt(12.0_dp))
947  END DO
948 ! *** Velocity Verlet for Langeving *** v(t)--> v(t+1/2)
949  CALL tamc_velocities_colvar(force_env, an)
950 ! *** Velocity Verlet for Langevin S(t)->S(t+1)
951  CALL tamc_position_colvar(force_env, xieta)
952 !!!!! start zHMC sampler
953  CALL hmcsampler(globenv, force_env, averages, r, mc_par, moves, gmoves, rng_stream_mc, output_unit, fz, zbuff)
954 
955 ! CALL final_mc_write(mc_par,tmp_moves,&
956 ! output_unit,energy_check,&
957 ! initial_energy,final_energy,&
958 ! averages)
959 
960 !!!!!!!!!!!!!!!!!!!! end zHMC sampler
961  ! *** Velocity Verlet for Langeving *** v(t+1/2)--> v(t+1)
962  CALL tamc_velocities_colvar(force_env, an)
963 ! CALL virial_evaluate ( atomic_kind_set, particle_set, &
964 ! local_particles, virial, para_env)
965 
966  END SUBROUTINE langevinvec
967 
968 ! **************************************************************************************************
969 !> \brief Driver routin for the canonical sampler using modified HMC
970 !> \param globenv ...
971 !> \param force_env ...
972 !> \param averages ...
973 !> \param r ...
974 !> \param mc_par ...
975 !> \param moves ...
976 !> \param gmoves ...
977 !> \param rng_stream_mc ...
978 !> \param output_unit ...
979 !> \param fz ...
980 !> \param zbuff ...
981 !> \param nskip ...
982 !> \author Alin M Elena
983 !> \note at the end of this routine %ff_s will contain mean force
984 ! **************************************************************************************************
985 
986  SUBROUTINE hmcsampler(globenv, force_env, averages, r, mc_par, moves, gmoves, rng_stream_mc, output_unit, &
987  fz, zbuff, nskip)
988  TYPE(global_environment_type), POINTER :: globenv
989  TYPE(force_env_type), POINTER :: force_env
990  TYPE(mc_averages_type), POINTER :: averages
991  REAL(kind=dp), DIMENSION(:, :), INTENT(INOUT) :: r
992  TYPE(mc_simpar_type), POINTER :: mc_par
993  TYPE(mc_moves_type), POINTER :: moves, gmoves
994  TYPE(rng_stream_type), INTENT(INOUT) :: rng_stream_mc
995  INTEGER, INTENT(IN) :: output_unit
996  REAL(kind=dp), DIMENSION(:), INTENT(INOUT) :: fz, zbuff
997  INTEGER, INTENT(IN), OPTIONAL :: nskip
998 
999  INTEGER :: i, iprint, ishift, it1, j, nsamples, &
1000  nstep
1001  REAL(kind=dp) :: energy_check, old_epx, old_epz, t1
1002  TYPE(meta_env_type), POINTER :: meta_env_saved
1003 
1004  IF (PRESENT(nskip)) THEN
1005  nsamples = nskip
1006  ishift = nskip
1007  ELSE
1008  nsamples = 0
1009  fz = 0.0_dp
1010  ishift = 0
1011  END IF
1012 ! lrestart = .false.
1013 ! if (present(logger) .and. present(iter)) THEN
1014 ! lrestart=.true.
1015 ! ENDIF
1016  CALL get_mc_par(mc_par, nstep=nstep, iprint=iprint)
1017  meta_env_saved => force_env%meta_env
1018  NULLIFY (force_env%meta_env)
1019  CALL force_env_get(force_env, potential_energy=old_epx)
1020  force_env%meta_env => meta_env_saved
1021 
1022  old_epz = force_env%meta_env%epot_s
1023 !!! average energy will be wrong on restarts
1024  averages%ave_energy = 0.0_dp
1025  t1 = force_env%qs_env%sim_time
1026  it1 = force_env%qs_env%sim_step
1027  IF (output_unit > 0) THEN
1028  WRITE (output_unit, '(a,l4)') "HMC|restart? ", force_env%meta_env%restart
1029  WRITE (output_unit, '(a,3(f16.8,1x))') &
1030  "HMC|Ep, Epx, Epz ", old_epx + force_env%meta_env%epot_s, old_epx, force_env%meta_env%epot_s
1031  WRITE (output_unit, '(a)') "#HMC| No | z.. | theta.. | ff_z... | ff_z/n |"
1032  END IF
1033  DO i = 1, nstep
1034  IF (mod(i, iprint) == 0 .AND. (output_unit > 0)) THEN
1035  WRITE (output_unit, '(a,1x,i0)') "HMC|========== On Monte Carlo cycle ", i + ishift
1036  WRITE (output_unit, '(a)') "HMC| Attempting a minitrajectory move"
1037  WRITE (output_unit, '(a,1x,i0)') "HMC| start mini-trajectory", i + ishift
1038  WRITE (output_unit, '(a,1x,i0,1x)', advance="no") "#HMC|0 ", i + ishift
1039  DO j = 1, force_env%meta_env%n_colvar
1040  WRITE (output_unit, '(f16.8,1x,f16.8,1x,f16.8)', advance="no") force_env%meta_env%metavar(j)%ss0, &
1041  force_env%meta_env%metavar(j)%ss, &
1042  force_env%meta_env%metavar(j)%ff_s !,fz(j)/real(i+ishift,dp)
1043  END DO
1044  WRITE (output_unit, *)
1045  END IF
1046 
1047  CALL mc_hmc_move(mc_par, force_env, globenv, moves, gmoves, old_epx, old_epz, energy_check, &
1048  r, output_unit, rng_stream_mc, zbuff)
1049  ! check averages...
1050  ! force average for z needed too...
1051  averages%ave_energy = averages%ave_energy*real(i - 1, dp)/real(i, dp) + &
1052  old_epx/real(i, dp)
1053  DO j = 1, force_env%meta_env%n_colvar
1054  fz(j) = fz(j) + force_env%meta_env%metavar(j)%ff_s
1055  END DO
1056  IF (output_unit > 0) THEN
1057  WRITE (output_unit, '(a,1x,i0)') "HMC|end mini-trajectory", i + ishift
1058 !!!!!!!! this prints z and theta(x) --ss0,ss-- needed to determine an acceptable k then
1059  ! the instanteneous force and some instanteneous average for force
1060  WRITE (output_unit, '(a,1x,i0,1x)', advance="no") "#HMC|1 ", i + ishift
1061  DO j = 1, force_env%meta_env%n_colvar
1062  WRITE (output_unit, '(f16.8,1x,f16.8,1x,f16.8,1x,f16.8)', advance="no") force_env%meta_env%metavar(j)%ss0, &
1063  force_env%meta_env%metavar(j)%ss, &
1064  force_env%meta_env%metavar(j)%ff_s, fz(j)/real(i + ishift, dp)
1065  END DO
1066  WRITE (output_unit, *)
1067  END IF
1068  nsamples = nsamples + 1
1069  IF (mod(i, iprint) == 0 .AND. (output_unit > 0)) THEN
1070  WRITE (output_unit, '(a,f16.8)') "HMC| Running average for potential energy ", averages%ave_energy
1071  WRITE (output_unit, '(a,1x,i0)') "HMC|======== End Monte Carlo cycle ", i + ishift
1072  END IF
1073 ! IF (lrestart) THEN
1074 ! k=nstep/5
1075 ! IF(MOD(i,k) == 0) THEN
1076 ! force_env%qs_env%sim_time=t1
1077 ! force_env%qs_env%sim_step=it1
1078 ! DO j=1,force_env%meta_env%n_colvar
1079 ! force_env%meta_env%metavar(j)%ff_s=fz(j)/real(i+ishift,dp)
1080 ! ENDDO
1081 ! ! CALL cp_iterate(logger%iter_info,last=.FALSE.,iter_nr=-1)
1082 ! CALL section_vals_val_set(mcsec,"RANDOMTOSKIP",i_val=i+ishift)
1083 ! CALL write_restart(md_env=mdenv,root_section=force_env%root_section)
1084 ! ! CALL cp_iterate(logger%iter_info,last=.FALSE.,iter_nr=iter)
1085 ! ENDIF
1086 ! ENDIF
1087  END DO
1088  force_env%qs_env%sim_time = t1
1089  force_env%qs_env%sim_step = it1
1090  IF (output_unit > 0) THEN
1091  WRITE (output_unit, '(a,i0,a,i0,a,f16.8)') "HMC| local acceptance ratio: ", moves%hmc%successes, "/", &
1092  moves%hmc%attempts, "=", real(moves%hmc%successes, dp)/real(moves%hmc%attempts, dp)
1093  WRITE (output_unit, '(a,i0,a,i0,a,f16.8)') "HMC| global acceptance ratio: ", gmoves%hmc%successes, "/", &
1094  gmoves%hmc%attempts, "=", real(gmoves%hmc%successes, dp)/real(gmoves%hmc%attempts, dp)
1095  END IF
1096  !average force
1097  DO j = 1, force_env%meta_env%n_colvar
1098  force_env%meta_env%metavar(j)%ff_s = fz(j)/nsamples
1099  END DO
1100  END SUBROUTINE hmcsampler
1101 
1102 ! **************************************************************************************************
1103 !> \brief performs a hybrid Monte Carlo move
1104 !> \param mc_par ...
1105 !> \param force_env ...
1106 !> \param globenv ...
1107 !> \param moves ...
1108 !> \param gmoves ...
1109 !> \param old_epx ...
1110 !> \param old_epz ...
1111 !> \param energy_check ...
1112 !> \param r ...
1113 !> \param output_unit ...
1114 !> \param rng_stream ...
1115 !> \param zbuff ...
1116 !> \author Alin M Elena
1117 !> \note It runs a NVE trajectory, followed by localisation and accepts rejects
1118 !> using the biased Hamiltonian, rather than the traditional guiding Hamiltonian
1119 ! **************************************************************************************************
1120  SUBROUTINE mc_hmc_move(mc_par, force_env, globenv, moves, gmoves, old_epx, old_epz, &
1121  energy_check, r, output_unit, rng_stream, zbuff)
1122 
1123  TYPE(mc_simpar_type), POINTER :: mc_par
1124  TYPE(force_env_type), POINTER :: force_env
1125  TYPE(global_environment_type), POINTER :: globenv
1126  TYPE(mc_moves_type), POINTER :: moves, gmoves
1127  REAL(kind=dp), INTENT(INOUT) :: old_epx, old_epz, energy_check
1128  REAL(kind=dp), DIMENSION(:, :), INTENT(INOUT) :: r
1129  INTEGER, INTENT(IN) :: output_unit
1130  TYPE(rng_stream_type), INTENT(INOUT) :: rng_stream
1131  REAL(kind=dp), DIMENSION(:), INTENT(INOUT) :: zbuff
1132 
1133  CHARACTER(LEN=*), PARAMETER :: routinen = 'mc_hmc_move'
1134 
1135  INTEGER :: handle, iatom, j, natoms, source
1136  LOGICAL :: ionode, localise
1137  REAL(kind=dp) :: beta, energy_term, exp_max_val, &
1138  exp_min_val, new_energy, new_epx, &
1139  new_epz, rand, value, w
1140  TYPE(cp_subsys_type), POINTER :: oldsys
1141  TYPE(mc_ekin_type), POINTER :: hmc_ekin
1142  TYPE(meta_env_type), POINTER :: meta_env_saved
1143  TYPE(mp_comm_type) :: group
1144  TYPE(particle_list_type), POINTER :: particles_set
1145  TYPE(section_vals_type), POINTER :: dft_section, input
1146 
1147 ! begin the timing of the subroutine
1148 
1149  CALL timeset(routinen, handle)
1150 
1151  CALL get_qs_env(force_env%qs_env, input=input)
1152  dft_section => section_vals_get_subs_vals(input, "DFT")
1153 
1154 ! get a bunch of stuff from mc_par
1155  CALL get_mc_par(mc_par, ionode=ionode, &
1156  beta=beta, exp_max_val=exp_max_val, &
1157  exp_min_val=exp_min_val, source=source, group=group)
1158 
1159 ! nullify some pointers
1160 ! NULLIFY(particles_set,oldsys,hmc_ekin)
1161  NULLIFY (particles_set, oldsys, meta_env_saved, hmc_ekin)
1162  ! now let's grab the particle positions
1163  CALL force_env_get(force_env, subsys=oldsys)
1164  CALL cp_subsys_get(oldsys, particles=particles_set)
1165  natoms = SIZE(particles_set%els)
1166 ! do some allocation
1167 
1168  ALLOCATE (hmc_ekin)
1169 
1170 ! record the attempt
1171  moves%hmc%attempts = moves%hmc%attempts + 1
1172  gmoves%hmc%attempts = gmoves%hmc%attempts + 1
1173 
1174 ! save the old coordinates just in case we need to go back
1175  DO iatom = 1, natoms
1176  r(1:3, iatom) = particles_set%els(iatom)%r(1:3)
1177  END DO
1178  localise = .true.
1179 ! the same for collective variables data should be made,ss first half and ff_s the last half
1180  DO j = 1, force_env%meta_env%n_colvar
1181  zbuff(j) = force_env%meta_env%metavar(j)%ss
1182  zbuff(j + force_env%meta_env%n_colvar) = force_env%meta_env%metavar(j)%ff_s
1183  IF ((oldsys%colvar_p(force_env%meta_env%metavar(j)%icolvar)%colvar%type_id == hbp_colvar_id) .OR. &
1184  (oldsys%colvar_p(force_env%meta_env%metavar(j)%icolvar)%colvar%type_id == wc_colvar_id)) THEN
1185  localise = .false.
1186  END IF
1187  END DO
1188 
1189 ! now run the MD simulation
1190  meta_env_saved => force_env%meta_env
1191  NULLIFY (force_env%meta_env)
1192  force_env%qs_env%sim_time = 0.0_dp
1193  force_env%qs_env%sim_step = 0
1194  IF (.NOT. localise) THEN
1195  CALL section_vals_val_set(dft_section, "LOCALIZE%_SECTION_PARAMETERS_", l_val=.false.)
1196  END IF
1197  CALL qs_mol_dyn(force_env, globenv, hmc_e_initial=hmc_ekin%initial_ekin, hmc_e_final=hmc_ekin%final_ekin)
1198  IF (.NOT. localise) THEN
1199  CALL section_vals_val_set(dft_section, "LOCALIZE%_SECTION_PARAMETERS_", l_val=.true.)
1200  CALL scf_post_calculation_gpw(force_env%qs_env)
1201  END IF
1202 
1203  CALL force_env_get(force_env, potential_energy=new_epx)
1204 
1205  force_env%meta_env => meta_env_saved
1206  CALL tamc_force(force_env, zpot=new_epz)
1207  new_energy = new_epx + new_epz
1208  IF (output_unit > 0) THEN
1209  WRITE (output_unit, '(a,4(f16.8,1x))') &
1210  "HMC|old Ep, Ekx, Epz, Epx ", old_epx + old_epz, hmc_ekin%initial_ekin, old_epz, old_epx
1211  WRITE (output_unit, '(a,4(f16.8,1x))') "HMC|new Ep, Ekx, Epz, Epx ", new_energy, hmc_ekin%final_ekin, new_epz, new_epx
1212  END IF
1213  energy_term = new_energy - old_epx - old_epz + hmc_ekin%final_ekin - hmc_ekin%initial_ekin
1214 
1215  value = -beta*(energy_term)
1216 ! to prevent overflows
1217  IF (value > exp_max_val) THEN
1218  w = 10.0_dp
1219  ELSEIF (value < exp_min_val) THEN
1220  w = 0.0_dp
1221  ELSE
1222  w = exp(value)
1223  END IF
1224 
1225  rand = rng_stream%next()
1226  IF (rand < w) THEN
1227 ! accept the move
1228  moves%hmc%successes = moves%hmc%successes + 1
1229  gmoves%hmc%successes = gmoves%hmc%successes + 1
1230 ! update energies
1231  energy_check = energy_check + (new_energy - old_epx - old_epz)
1232  old_epx = new_epx
1233  old_epz = new_epz
1234  ELSE
1235 ! reset the cell and particle positions
1236  DO iatom = 1, natoms
1237  particles_set%els(iatom)%r(1:3) = r(1:3, iatom)
1238  END DO
1239  DO j = 1, force_env%meta_env%n_colvar
1240  force_env%meta_env%metavar(j)%ss = zbuff(j)
1241  force_env%meta_env%metavar(j)%ff_s = zbuff(j + force_env%meta_env%n_colvar)
1242  END DO
1243 
1244  END IF
1245 
1246  DEALLOCATE (hmc_ekin)
1247 
1248 ! end the timing
1249  CALL timestop(handle)
1250 
1251  END SUBROUTINE mc_hmc_move
1252 
1253 ! **************************************************************************************************
1254 !> \brief ...
1255 !> \param force_env ...
1256 ! **************************************************************************************************
1257  SUBROUTINE metadyn_write_colvar_header(force_env)
1258  TYPE(force_env_type), POINTER :: force_env
1259 
1260  CHARACTER(len=*), PARAMETER :: routinen = 'metadyn_write_colvar_header'
1261 
1262  CHARACTER(len=100) :: aux, fmt
1263  CHARACTER(len=255) :: label1, label2, label3, label4, label5, &
1264  label6
1265  INTEGER :: handle, i, iw, m
1266  TYPE(cp_logger_type), POINTER :: logger
1267  TYPE(meta_env_type), POINTER :: meta_env
1268 
1269  NULLIFY (logger, meta_env)
1270  meta_env => force_env%meta_env
1271  IF (.NOT. ASSOCIATED(meta_env)) RETURN
1272 
1273  CALL timeset(routinen, handle)
1274  logger => cp_get_default_logger()
1275 
1276  iw = cp_print_key_unit_nr(logger, meta_env%metadyn_section, &
1277  "PRINT%COLVAR", extension=".metadynLog")
1278  IF (iw > 0) THEN
1279  label1 = ""
1280  label2 = ""
1281  label3 = ""
1282  label4 = ""
1283  label5 = ""
1284  label6 = ""
1285  DO i = 1, meta_env%n_colvar
1286  WRITE (aux, '(a,i0)') "z_", i
1287  label1 = trim(label1)//trim(aux)
1288  m = 15*i - len_trim(label1) - 1
1289  label1 = trim(label1)//repeat(" ", m)//"|"
1290  WRITE (aux, '(a,i0)') "Theta_", i
1291  label2 = trim(label2)//trim(aux)
1292  m = 15*i - len_trim(label2) - 1
1293  label2 = trim(label2)//repeat(" ", m)//"|"
1294  WRITE (aux, '(a,i0)') "F_z", i
1295  label3 = trim(label3)//trim(aux)
1296  m = 15*i - len_trim(label3) - 1
1297  label3 = trim(label3)//repeat(" ", m)//"|"
1298  WRITE (aux, '(a,i0)') "F_h", i
1299  label4 = trim(label4)//trim(aux)
1300  m = 15*i - len_trim(label4) - 1
1301  label4 = trim(label4)//repeat(" ", m)//"|"
1302  WRITE (aux, '(a,i0)') "F_w", i
1303  label5 = trim(label5)//trim(aux)
1304  m = 15*i - len_trim(label5) - 1
1305  label5 = trim(label5)//repeat(" ", m)//"|"
1306  WRITE (aux, '(a,i0)') "v_z", i
1307  label6 = trim(label6)//trim(aux)
1308  m = 15*i - len_trim(label6) - 1
1309  label6 = trim(label6)//repeat(" ", m)//"|"
1310  END DO
1311  WRITE (fmt, '("(a17,6a",i0 ,",4a15)")') meta_env%n_colvar*15
1312  WRITE (iw, trim(fmt)) "#Time[fs] |", &
1313  trim(label1), &
1314  trim(label2), &
1315  trim(label3), &
1316  trim(label4), &
1317  trim(label5), &
1318  trim(label6), &
1319  "Epot_z |", &
1320  "Ene hills |", &
1321  "Epot walls |", &
1322  "Temperature |"
1323 
1324  END IF
1325  CALL cp_print_key_finished_output(iw, logger, meta_env%metadyn_section, &
1326  "PRINT%COLVAR")
1327 
1328  CALL timestop(handle)
1329 
1330  END SUBROUTINE metadyn_write_colvar_header
1331 
1332 ! **************************************************************************************************
1333 !> \brief ...
1334 !> \param force_env ...
1335 ! **************************************************************************************************
1336  SUBROUTINE metadyn_write_colvar(force_env)
1337  TYPE(force_env_type), POINTER :: force_env
1338 
1339  CHARACTER(len=*), PARAMETER :: routinen = 'metadyn_write_colvar'
1340 
1341  INTEGER :: handle, i, i_c, iw
1342  REAL(kind=dp) :: temp
1343  TYPE(cp_logger_type), POINTER :: logger
1344  TYPE(meta_env_type), POINTER :: meta_env
1345  TYPE(metavar_type), POINTER :: cv
1346 
1347  NULLIFY (logger, meta_env, cv)
1348  meta_env => force_env%meta_env
1349  IF (.NOT. ASSOCIATED(meta_env)) RETURN
1350 
1351  CALL timeset(routinen, handle)
1352  logger => cp_get_default_logger()
1353 
1354  IF (meta_env%langevin) THEN
1355  meta_env%ekin_s = 0.0_dp
1356 ! meta_env%epot_s = 0.0_dp
1357  DO i_c = 1, meta_env%n_colvar
1358  cv => meta_env%metavar(i_c)
1359  meta_env%ekin_s = meta_env%ekin_s + 0.5_dp*cv%mass*cv%vvp**2
1360  END DO
1361  END IF
1362 
1363  ! write COLVAR file
1364  iw = cp_print_key_unit_nr(logger, meta_env%metadyn_section, &
1365  "PRINT%COLVAR", extension=".metadynLog")
1366  IF (iw > 0) THEN
1367  IF (meta_env%extended_lagrange) THEN
1368  WRITE (iw, '(f16.8,70f15.8)') meta_env%time*femtoseconds, &
1369  (meta_env%metavar(i)%ss0, i=1, meta_env%n_colvar), &
1370  (meta_env%metavar(i)%ss, i=1, meta_env%n_colvar), &
1371  (meta_env%metavar(i)%ff_s, i=1, meta_env%n_colvar), &
1372  (meta_env%metavar(i)%ff_hills, i=1, meta_env%n_colvar), &
1373  (meta_env%metavar(i)%ff_walls, i=1, meta_env%n_colvar), &
1374  (meta_env%metavar(i)%vvp, i=1, meta_env%n_colvar), &
1375  meta_env%epot_s, &
1376  meta_env%hills_env%energy, &
1377  meta_env%epot_walls, &
1378  (meta_env%ekin_s)*2.0_dp/(real(meta_env%n_colvar, kind=dp))*kelvin
1379  ELSE
1380  WRITE (iw, '(f16.8,40f13.5)') meta_env%time*femtoseconds, &
1381  (meta_env%metavar(i)%ss0, i=1, meta_env%n_colvar), &
1382  (meta_env%metavar(i)%ff_hills, i=1, meta_env%n_colvar), &
1383  (meta_env%metavar(i)%ff_walls, i=1, meta_env%n_colvar), &
1384  meta_env%hills_env%energy, &
1385  meta_env%epot_walls
1386  END IF
1387  END IF
1388  CALL cp_print_key_finished_output(iw, logger, meta_env%metadyn_section, &
1389  "PRINT%COLVAR")
1390  ! Temperature for COLVAR
1391  IF (meta_env%extended_lagrange) THEN
1392  temp = meta_env%ekin_s*2.0_dp/(real(meta_env%n_colvar, kind=dp))*kelvin
1393  meta_env%avg_temp = (meta_env%avg_temp*real(meta_env%n_steps, kind=dp) + &
1394  temp)/real(meta_env%n_steps + 1, kind=dp)
1395  iw = cp_print_key_unit_nr(logger, meta_env%metadyn_section, &
1396  "PRINT%TEMPERATURE_COLVAR", extension=".metadynLog")
1397  IF (iw > 0) THEN
1398  WRITE (iw, '(T2,79("-"))')
1399  WRITE (iw, '( A,T51,f10.2,T71,f10.2)') ' COLVARS INSTANTANEOUS/AVERAGE TEMPERATURE ', &
1400  temp, meta_env%avg_temp
1401  WRITE (iw, '(T2,79("-"))')
1402  END IF
1403  CALL cp_print_key_finished_output(iw, logger, meta_env%metadyn_section, &
1404  "PRINT%TEMPERATURE_COLVAR")
1405  END IF
1406  CALL timestop(handle)
1407 
1408  END SUBROUTINE metadyn_write_colvar
1409 
1410 ! **************************************************************************************************
1411 !> \brief ...
1412 !> \param force_env ...
1413 ! **************************************************************************************************
1414  SUBROUTINE setup_velocities_z(force_env)
1415  TYPE(force_env_type), POINTER :: force_env
1416 
1417  INTEGER :: i_c
1418  REAL(kind=dp) :: ekin_w, fac_t
1419  TYPE(meta_env_type), POINTER :: meta_env
1420  TYPE(metavar_type), POINTER :: cv
1421 
1422  NULLIFY (meta_env)
1423  meta_env => force_env%meta_env
1424  meta_env%ekin_s = 0.0_dp
1425  DO i_c = 1, meta_env%n_colvar
1426  cv => meta_env%metavar(i_c)
1427  cv%vvp = force_env%globenv%gaussian_rng_stream%next()
1428  meta_env%ekin_s = meta_env%ekin_s + 0.5_dp*cv%mass*cv%vvp**2
1429  END DO
1430  ekin_w = 0.5_dp*meta_env%temp_wanted*real(meta_env%n_colvar, kind=dp)
1431  fac_t = sqrt(ekin_w/max(meta_env%ekin_s, 1.0e-8_dp))
1432  DO i_c = 1, meta_env%n_colvar
1433  cv => meta_env%metavar(i_c)
1434  cv%vvp = cv%vvp*fac_t
1435  END DO
1436  END SUBROUTINE setup_velocities_z
1437 END MODULE tamc_run
represent a simple array based list of the given type
Define the atomic kind types and their sub types.
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)
...
collects all references to literature in CP2K as new algorithms / method are included from literature...
Definition: bibliography.F:28
integer, save, public vandencic2006
Definition: bibliography.F:43
Handles all functions related to the CELL.
Definition: cell_types.F:15
defines collective variables s({R}) and the derivative of this variable wrt R these can then be used ...
subroutine, public colvar_eval_glob_f(icolvar, force_env)
evaluates the derivatives (dsdr) given and due to the given colvar
Initialize the collective variables types.
Definition: colvar_types.F:15
integer, parameter, public wc_colvar_id
Definition: colvar_types.F:56
integer, parameter, public hbp_colvar_id
Definition: colvar_types.F:56
subroutine, public fix_atom_control(force_env, w)
allows for fix atom constraints
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)
returns information about various attributes of the given subsys
unit conversion facility
Definition: cp_units.F:30
real(kind=dp) function, public cp_unit_from_cp2k(value, unit_str, defaults, power)
converts from the internal cp2k units to the given unit
Definition: cp_units.F:1179
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
defines types for metadynamics calculation
subroutine, public fe_env_create(fe_env, fe_section)
creates the fe_env
Calculation of the incomplete Gamma function F_n(t) for multi-center integrals over Cartesian Gaussia...
Definition: gamma.F:15
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 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.
Set of routines to dump the restart file of CP2K.
subroutine, public write_restart(md_env, force_env, root_section, coords, vels, pint_env, helium_env)
checks if a restart needs to be written and does so, updating all necessary fields in the input file....
objects that represent the structure of input sections and the data contained in an input section
subroutine, public section_vals_val_set(section_vals, keyword_name, i_rep_section, i_rep_val, val, l_val, i_val, r_val, c_val, l_vals_ptr, i_vals_ptr, r_vals_ptr, c_vals_ptr)
sets the requested value
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_get(section_vals, ref_count, n_repetition, n_subs_vals_rep, section, explicit)
returns various attributes about the section_vals
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
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
contains the subroutines for dealing with the mc_env
subroutine, public get_mc_env(mc_env, mc_par, force_env)
provides a method for getting the various structures attached to an mc_env
subroutine, public mc_env_create(mc_env)
creates and initializes an mc_env
subroutine, public set_mc_env(mc_env, mc_par, force_env)
provides a method for attaching various structures to an mc_env
subroutine, public mc_env_release(mc_env)
releases the given mc env
contains miscellaneous subroutines used in the Monte Carlo runs, mostly I/O stuff
Definition: mc_misc.F:13
subroutine, public mc_averages_release(averages)
deallocates the structure that holds running averages of MC variables
Definition: mc_misc.F:77
subroutine, public mc_averages_create(averages)
initializes the structure that holds running averages of MC variables
Definition: mc_misc.F:46
control the handling of the move data in Monte Carlo (MC) simulations
subroutine, public mc_moves_release(moves)
deallocates all the structures and nullifies the pointer
subroutine, public init_mc_moves(moves)
allocates and initializes the structure to record all move attempts/successes
holds all the structure types needed for Monte Carlo, except the mc_environment_type
Definition: mc_types.F:15
subroutine, public set_mc_par(mc_par, rm, cl, diff, nstart, rmvolume, rmcltrans, rmbond, rmangle, rmdihedral, rmrot, rmtrans, PROGRAM, nmoves, nswapmoves, lstop, temperature, pressure, rclus, iuptrans, iupcltrans, iupvolume, pmswap, pmvolume, pmtraion, pmtrans, pmcltrans, BETA, rcut, iprint, lbias, nstep, lrestart, ldiscrete, discrete_step, pmavbmc, mc_molecule_info, pmavbmc_mol, pmtrans_mol, pmrot_mol, pmtraion_mol, pmswap_mol, avbmc_rmin, avbmc_rmax, avbmc_atom, pbias, ensemble, pmvol_box, pmclus_box, eta, mc_input_file, mc_bias_file, exp_max_val, exp_min_val, min_val, max_val, pmhmc, pmhmc_box, lhmc, ionode, source, group, rand2skip)
changes the private elements of the mc_parameters_type
Definition: mc_types.F:667
subroutine, public get_mc_par(mc_par, nstep, nvirial, iuptrans, iupcltrans, iupvolume, nmoves, nswapmoves, rm, cl, diff, nstart, source, group, lbias, ionode, lrestart, lstop, rmvolume, rmcltrans, rmbond, rmangle, rmrot, rmtrans, temperature, pressure, rclus, BETA, pmswap, pmvolume, pmtraion, pmtrans, pmcltrans, ensemble, PROGRAM, restart_file_name, molecules_file, moves_file, coords_file, energy_file, displacement_file, cell_file, dat_file, data_file, box2_file, fft_lib, iprint, rcut, ldiscrete, discrete_step, pmavbmc, pbias, avbmc_atom, avbmc_rmin, avbmc_rmax, rmdihedral, input_file, mc_molecule_info, pmswap_mol, pmavbmc_mol, pmtrans_mol, pmrot_mol, pmtraion_mol, mc_input_file, mc_bias_file, pmvol_box, pmclus_box, virial_temps, exp_min_val, exp_max_val, min_val, max_val, eta, pmhmc, pmhmc_box, lhmc, rand2skip)
...
Definition: mc_types.F:405
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 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 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
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
Interface to the message passing library MPI.
defines types for metadynamics calculation
subroutine, public set_meta_env(meta_env, time)
sets the meta_env
represent a simple array based list of the given type
Define the molecule kind structure types and the corresponding functionality.
represent a simple array based list of the given type
Define the data structure for the molecule information.
Parallel (pseudo)random number generator (RNG) for multiple streams and substreams of random numbers.
integer, parameter, public uniform
represent a simple array based list of the given type
Define the data structure for the particle information.
Definition of physical constants:
Definition: physcon.F:68
real(kind=dp), parameter, public boltzmann
Definition: physcon.F:129
real(kind=dp), parameter, public femtoseconds
Definition: physcon.F:153
real(kind=dp), parameter, public joule
Definition: physcon.F:159
real(kind=dp), parameter, public kelvin
Definition: physcon.F:165
subroutine, public apply_qmmm_walls_reflective(force_env)
Apply reflective QM walls in order to avoid QM atoms escaping from the QM Box.
Definition: qmmm_util.F:99
subroutine, public get_qs_env(qs_env, atomic_kind_set, qs_kind_set, cell, super_cell, cell_ref, use_ref_cell, kpoints, dft_control, mos, sab_orb, sab_all, qmmm, qmmm_periodic, sac_ae, sac_ppl, sac_lri, sap_ppnl, sab_vdw, sab_scp, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_RI_aux_kp, matrix_s, matrix_s_RI_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, ecoul_1c, rho0_s_rs, rho0_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, WannierCentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, rhs)
Get the QUICKSTEP environment.
Does all kind of post scf calculations for GPW/GAPW.
subroutine, public scf_post_calculation_gpw(qs_env, wf_type, do_mp2)
collects possible post - scf calculations and prints info / computes properties.
provides a uniform framework to add references to CP2K cite and output these
subroutine, public cite_reference(key)
marks a given reference as cited.
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 initialize_reftraj(reftraj, reftraj_section, md_env)
...
Definition: reftraj_util.F:82
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
Utilities for string manipulations.
elemental logical function, public str_comp(str1, str2)
...
Perform a temperature accelarated hybrid monte carlo (TAHMC) run using QUICKSTEP.
Definition: tamc_run.F:14
subroutine, public qs_tamc(force_env, globenv, averages)
Driver routine for TAHMC.
Definition: tamc_run.F:147
Thermal regions type: to initialize and control the temperature of different regions.
Setup of regions with different temperature.
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.
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...