56#include "./base/base_uses.f90"
61 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'metadynamics_utils'
83 SUBROUTINE metadyn_read(meta_env, force_env, root_section, para_env, fe_section)
90 CHARACTER(len=*),
PARAMETER :: routinen =
'metadyn_read'
92 CHARACTER(LEN=default_path_length) :: walkers_file_name
93 INTEGER :: handle, i, id_method, n_colvar, n_rep, &
94 number_allocated_colvars
95 INTEGER,
DIMENSION(:),
POINTER :: walkers_status
96 LOGICAL :: check, explicit
100 metavar_section, walkers_section
103 CALL timeset(routinen, handle)
107 number_allocated_colvars = 0
109 IF (
ASSOCIATED(subsys%colvar_p))
THEN
110 number_allocated_colvars =
SIZE(subsys%colvar_p)
114 CALL timestop(handle)
118 cpassert(.NOT.
ASSOCIATED(meta_env))
127 dt=dt, para_env=para_env, metadyn_section=metadyn_section)
131 IF (meta_env%use_plumed .EQV. .true.)
THEN
132 CALL section_vals_val_get(metadyn_section,
"PLUMED_INPUT_FILE", c_val=meta_env%plumed_input_file)
133 meta_env%plumed_input_file = trim(meta_env%plumed_input_file)//char(0)
134 meta_env%langevin = .false.
135 CALL timestop(handle)
142 IF (meta_env%TAMCSteps < 0)
THEN
143 cpabort(
"TAMCSteps must be positive!")
146 IF (meta_env%zdt <= 0.0_dp)
THEN
147 cpabort(
"Timestep must be positive!")
151 CALL section_vals_val_get(metadyn_section,
"MIN_NT_HILLS", i_val=meta_env%hills_env%min_nt_hills)
152 IF (meta_env%hills_env%nt_hills <= 0)
THEN
153 meta_env%hills_env%min_nt_hills = meta_env%hills_env%nt_hills
154 CALL cp_warn(__location__, &
155 "NT_HILLS has a value <= 0; "// &
156 "Setting MIN_NT_HILLS to the same value! "// &
157 "Overriding input specification!")
159 check = meta_env%hills_env%nt_hills >= meta_env%hills_env%min_nt_hills
161 CALL cp_abort(__location__,
"MIN_NT_HILLS must have a value smaller or equal to NT_HILLS! "// &
162 "Cross check with the input reference!")
165 CALL section_vals_val_get(metadyn_section,
"OLD_HILL_NUMBER", i_val=meta_env%hills_env%old_hill_number)
166 CALL section_vals_val_get(metadyn_section,
"OLD_HILL_STEP", i_val=meta_env%hills_env%old_hill_step)
169 CALL section_vals_val_get(metadyn_section,
"HILL_TAIL_CUTOFF", r_val=meta_env%hills_env%tail_cutoff)
177 cpassert(meta_env%n_steps >= 0)
179 i_val=meta_env%hills_env%n_hills)
183 r_val=meta_env%toll_temp)
186 r_val=meta_env%delta_t)
188 r_val=meta_env%wtgamma)
189 IF (meta_env%well_tempered)
THEN
190 meta_env%hills_env%wtcontrol = meta_env%hills_env%wtcontrol .OR. check
191 check = meta_env%hills_env%wtcontrol
193 CALL cp_abort(__location__,
"When using Well-Tempered metadynamics, "// &
194 "DELTA_T (or WTGAMMA) should be explicitly specified.")
195 IF (meta_env%extended_lagrange) &
196 CALL cp_abort(__location__, &
197 "Well-Tempered metadynamics not possible with extended-lagrangian formulation.")
198 IF (meta_env%hills_env%min_disp > 0.0_dp) &
199 CALL cp_abort(__location__, &
200 "Well-Tempered metadynamics not possible with Adaptive hills.")
204 r_val=meta_env%avg_temp)
207 CALL metavar_read(meta_env%metavar(i), meta_env%extended_lagrange, &
208 meta_env%langevin, i, metavar_section)
209 check = (meta_env%metavar(i)%icolvar <= number_allocated_colvars)
211 CALL cp_abort(__location__, &
212 "An error occurred in the specification of COLVAR for METAVAR. "// &
213 "Specified COLVAR #("//trim(adjustl(
cp_to_string(meta_env%metavar(i)%icolvar)))//
") "// &
214 "is larger than the maximum number of COLVARS defined in the SUBSYS ("// &
215 trim(adjustl(
cp_to_string(number_allocated_colvars)))//
") !")
219 IF (meta_env%do_multiple_walkers)
THEN
220 NULLIFY (walkers_status)
225 i_val=meta_env%multiple_walkers%walker_id)
227 i_val=meta_env%multiple_walkers%walkers_tot_nr)
229 i_val=meta_env%multiple_walkers%walkers_freq_comm)
232 ALLOCATE (meta_env%multiple_walkers%walkers_status(meta_env%multiple_walkers%walkers_tot_nr))
233 ALLOCATE (meta_env%multiple_walkers%walkers_file_name(meta_env%multiple_walkers%walkers_tot_nr))
237 check = (
SIZE(walkers_status) == meta_env%multiple_walkers%walkers_tot_nr)
239 CALL cp_abort(__location__, &
240 "Number of Walkers specified in the input does not match with the "// &
241 "size of the WALKERS_STATUS. Please check your input and in case "// &
242 "this is a restart run consider the possibility to switch off the "// &
243 "RESTART_WALKERS in the EXT_RESTART section! ")
244 meta_env%multiple_walkers%walkers_status = walkers_status
246 meta_env%multiple_walkers%walkers_status = 0
248 meta_env%multiple_walkers%n_hills_local = &
249 meta_env%multiple_walkers%walkers_status(meta_env%multiple_walkers%walker_id)
253 check = (n_rep == meta_env%multiple_walkers%walkers_tot_nr)
255 CALL cp_abort(__location__, &
256 "Number of Walkers specified in the input does not match with the "// &
257 "number of Walkers File names provided. Please check your input and in case "// &
258 "this is a restart run consider the possibility to switch off the "// &
259 "RESTART_WALKERS in the EXT_RESTART section! ")
262 i_rep_val=i, c_val=walkers_file_name)
263 meta_env%multiple_walkers%walkers_file_name(i) = walkers_file_name
268 CALL print_metadyn_info(meta_env, n_colvar, metadyn_section)
271 CALL timestop(handle)
282 SUBROUTINE print_metadyn_info(meta_env, n_colvar, metadyn_section)
284 INTEGER,
INTENT(IN) :: n_colvar
287 CHARACTER(len=*),
PARAMETER :: routinen =
'print_metadyn_info'
289 CHARACTER(LEN=10) :: my_id, my_tag
290 INTEGER :: handle, i, iw, j
294 TYPE(
section_type),
POINTER :: section, wall_section, work_section
296 CALL timeset(routinen, handle)
300 "PRINT%PROGRAM_RUN_INFO", extension=
".metadynLog")
301 NULLIFY (section, enum, keyword)
305 WRITE (iw,
'( /A )')
' METADYN| Meta Dynamics Protocol '
306 WRITE (iw,
'( A,T71,I10)')
' METADYN| Number of interval time steps to spawn hills', &
307 meta_env%hills_env%nt_hills
308 WRITE (iw,
'( A,T71,I10)')
' METADYN| Number of previously spawned hills', &
309 meta_env%hills_env%n_hills
310 IF (meta_env%extended_lagrange)
THEN
311 WRITE (iw,
'( A )')
' METADYN| Extended Lagrangian Scheme '
312 IF (meta_env%tempcontrol)
WRITE (iw,
'( A,T71,F10.2)') &
313 ' METADYN| Collective Variables Temperature control', meta_env%toll_temp
314 IF (meta_env%langevin)
THEN
315 WRITE (iw,
'(A,T71)')
' METADYN| Langevin Thermostat in use for COLVAR '
316 WRITE (iw,
'(A,T71,F10.4)')
' METADYN| Langevin Thermostat. Target Temperature = ', &
317 meta_env%temp_wanted*
kelvin
319 WRITE (iw,
'(A,T71,F10.4)')
' METADYN| COLVARS restarted average temperature ', &
322 IF (meta_env%do_hills)
THEN
323 WRITE (iw,
'( A )')
' METADYN| Spawning the Hills '
324 WRITE (iw,
'( A,T71,F10.3)')
' METADYN| Height of the Spawned Gaussian', meta_env%hills_env%ww
326 IF (meta_env%hills_env%min_disp .GT. 0.0_dp)
THEN
327 WRITE (iw,
'(A)')
' METADYN| Adapative meta time step is activated'
328 WRITE (iw,
'(A,T71,F10.4)')
' METADYN| Minimum displacement for next hill', &
329 meta_env%hills_env%min_disp
334 IF (meta_env%well_tempered)
THEN
335 WRITE (iw,
'( A )')
' METADYN| Well-Tempered metadynamics '
336 IF (meta_env%delta_t > epsilon(1._dp))
THEN
337 WRITE (iw,
'( A,T71,F10.3)')
' METADYN| Temperature parameter (Delta T) [K]', meta_env%delta_t*
kelvin
339 WRITE (iw,
'( A,T71,F10.3)')
' METADYN| Temperature parameter (gamma)', meta_env%wtgamma
343 IF (meta_env%do_multiple_walkers)
THEN
344 WRITE (iw,
'( A,T71,A10)')
' METADYN| Multiple Walkers',
' ENABLED'
345 WRITE (iw,
'( A,T71,I10)')
' METADYN| Number of Multiple Walkers', &
346 meta_env%multiple_walkers%walkers_tot_nr
347 WRITE (iw,
'( A,T71,I10)')
' METADYN| Local Walker ID', &
348 meta_env%multiple_walkers%walker_id
349 WRITE (iw,
'( A,T71,I10)')
' METADYN| Walker Communication Frequency', &
350 meta_env%multiple_walkers%walkers_freq_comm
351 DO i = 1, meta_env%multiple_walkers%walkers_tot_nr
353 IF (i == meta_env%multiple_walkers%walker_id) my_tag =
" ( Local )"
355 WRITE (iw,
'(/,A,T71,A10)')
' WALKERS| Walker ID'//trim(my_tag), adjustr(my_id)
356 WRITE (iw,
'( A,T71,I10)')
' WALKERS| Number of Hills communicated', &
357 meta_env%multiple_walkers%walkers_status(i)
358 WRITE (iw,
'( A,T24,A57)')
' WALKERS| Base Filename', &
359 adjustr(meta_env%multiple_walkers%walkers_file_name(i) (1:57))
364 WRITE (iw,
'( A,T71,I10)')
' METADYN| Number of collective variables', meta_env%n_colvar
366 WRITE (iw,
'( A )')
' '//
'----------------------------------------------------------------------'
367 WRITE (iw,
'( A,T71,I10)')
' METAVARS| Collective Variable Number', meta_env%metavar(i)%icolvar
368 IF (meta_env%extended_lagrange)
THEN
369 WRITE (iw,
'( A,T71,F10.6)')
' METAVARS| Lambda Parameter', meta_env%metavar(i)%lambda
370 WRITE (iw,
'( A,T66,F15.6)')
' METAVARS| Collective Variable Mass', meta_env%metavar(i)%mass
372 WRITE (iw,
'( A,T71,F10.6)')
' METAVARS| Scaling factor', meta_env%metavar(i)%delta_s
373 IF (meta_env%langevin)
WRITE (iw,
'( A,T71,F10.6)')
' METAVARS| Friction for Langevin Thermostat', &
374 meta_env%metavar(i)%gamma
375 IF (meta_env%metavar(i)%do_wall)
THEN
376 WRITE (iw,
'( A,T71,I10)')
' METAVARS| Number of Walls present',
SIZE(meta_env%metavar(i)%walls)
377 DO j = 1,
SIZE(meta_env%metavar(i)%walls)
380 WRITE (iw,
'(/,A,5X,I10,T50,A,T70,A11)')
' METAVARS| Wall Number:', j,
'Type of Wall:', &
381 adjustr(trim(
enum_i2c(enum, meta_env%metavar(i)%walls(j)%id_type)))
383 SELECT CASE (meta_env%metavar(i)%walls(j)%id_type)
391 WRITE (iw,
'(A,T70,A11)')
' METAVARS| Wall direction', &
392 adjustr(trim(
enum_i2c(enum, meta_env%metavar(i)%walls(j)%id_direction)))
397 WRITE (iw,
'(A,T70,A11)')
' METAVARS| Wall direction', &
398 adjustr(trim(
enum_i2c(enum, meta_env%metavar(i)%walls(j)%id_direction)))
399 WRITE (iw,
'(A,T70,F11.6)')
' METAVARS| Constant K of the quadratic potential', &
400 meta_env%metavar(i)%walls(j)%k_quadratic
402 WRITE (iw,
'(A,T70,F11.6)')
' METAVARS| Height of the Wall Gaussian', &
403 meta_env%metavar(i)%walls(j)%ww_gauss
404 WRITE (iw,
'(A,T70,F11.6)')
' METAVARS| Scale of the Wall Gaussian', &
405 meta_env%metavar(i)%walls(j)%sigma_gauss
407 WRITE (iw,
'(A,T70,F11.6)')
' METAVARS| Wall location', &
408 meta_env%metavar(i)%walls(j)%pos
411 WRITE (iw,
'( A )')
' '//
'----------------------------------------------------------------------'
417 CALL timestop(handle)
419 END SUBROUTINE print_metadyn_info
433 SUBROUTINE metavar_read(metavar, extended_lagrange, langevin, icol, metavar_section)
435 LOGICAL,
INTENT(IN) :: extended_lagrange, langevin
436 INTEGER,
INTENT(IN) :: icol
439 CHARACTER(len=*),
PARAMETER :: routinen =
'metavar_read'
441 INTEGER :: handle, i, n_walls
444 CALL timeset(routinen, handle)
451 IF (n_walls /= 0)
THEN
452 metavar%do_wall = .true.
453 ALLOCATE (metavar%walls(n_walls))
455 CALL section_vals_val_get(wall_section,
"TYPE", i_rep_section=i, i_val=metavar%walls(i)%id_type)
456 CALL section_vals_val_get(wall_section,
"POSITION", i_rep_section=i, r_val=metavar%walls(i)%pos)
457 SELECT CASE (metavar%walls(i)%id_type)
472 SELECT CASE (metavar%walls(i)%id_direction)
474 metavar%walls(i)%pos0 = metavar%walls(i)%pos + (0.05_dp/metavar%walls(i)%k_quartic**(0.25_dp))
476 metavar%walls(i)%pos0 = metavar%walls(i)%pos - (0.05_dp/metavar%walls(i)%k_quartic**(0.25_dp))
486 IF (extended_lagrange)
THEN
494 CALL timestop(handle)
496 END SUBROUTINE metavar_read
508 n_colvar, metadyn_section)
512 INTEGER,
INTENT(IN) :: n_colvar
515 CHARACTER(len=*),
PARAMETER :: routinen =
'synchronize_multiple_walkers'
517 CHARACTER(LEN=default_path_length) :: filename, tmpname
518 INTEGER :: delta_hills, handle, i, i_hills, ih, iw, &
521 REAL(kind=
dp) :: invdt, ww
522 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: delta_s_save, ss0_save
523 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: delta_s_ss0_buf
527 CALL timeset(routinen, handle)
530 para_env => logger%para_env
533 IF (para_env%is_source())
THEN
535 i = multiple_walkers%walker_id
536 filename = trim(multiple_walkers%walkers_file_name(i))//
"_"// &
537 trim(adjustl(
cp_to_string(multiple_walkers%n_hills_local)))
538 tmpname = trim(filename)//
".tmp"
539 CALL open_file(file_name=tmpname, file_status=
"UNKNOWN", &
540 file_form=
"FORMATTED", file_action=
"WRITE", &
541 file_position=
"APPEND", unit_number=unit_nr)
542 WRITE (unit_nr, *) hills_env%ww_history(hills_env%n_hills)
544 WRITE (unit_nr, *) hills_env%ss_history(ih, hills_env%n_hills)
545 WRITE (unit_nr, *) hills_env%delta_s_history(ih, hills_env%n_hills)
547 IF (hills_env%wtcontrol)
WRITE (unit_nr, *) hills_env%invdt_history(hills_env%n_hills)
549 CALL m_mov(tmpname, filename)
552 IF (
modulo(multiple_walkers%n_hills_local, multiple_walkers%walkers_freq_comm) == 0)
THEN
554 ALLOCATE (ss0_save(n_colvar))
555 ALLOCATE (delta_s_save(n_colvar))
556 ALLOCATE (delta_s_ss0_buf(2, 0:n_colvar))
559 ss0_save(i) = colvars(i)%ss0
560 delta_s_save(i) = colvars(i)%delta_s
564 DO i = 1, multiple_walkers%walkers_tot_nr
565 IF (i == multiple_walkers%walker_id)
THEN
567 multiple_walkers%walkers_status(i) = multiple_walkers%n_hills_local
571 i_hills = multiple_walkers%walkers_status(i) + 1
572 filename = trim(multiple_walkers%walkers_file_name(i))//
"_"// &
575 IF (para_env%is_source())
THEN
576 INQUIRE (file=trim(filename), exist=exist)
578 CALL para_env%bcast(exist)
583 IF (para_env%is_source())
THEN
584 CALL open_file(file_name=filename, file_status=
"OLD", &
585 file_form=
"FORMATTED", file_action=
"READ", &
586 file_position=
"REWIND", unit_number=unit_nr)
587 READ (unit_nr, *) delta_s_ss0_buf(1, 0)
589 READ (unit_nr, *) delta_s_ss0_buf(1, ih)
590 READ (unit_nr, *) delta_s_ss0_buf(2, ih)
592 IF (hills_env%wtcontrol)
READ (unit_nr, *) delta_s_ss0_buf(2, 0)
595 CALL para_env%bcast(delta_s_ss0_buf)
596 ww = delta_s_ss0_buf(1, 0)
597 IF (hills_env%wtcontrol) invdt = delta_s_ss0_buf(2, 0)
599 colvars(ih)%ss0 = delta_s_ss0_buf(1, ih)
600 colvars(ih)%delta_s = delta_s_ss0_buf(2, ih)
604 IF (hills_env%wtcontrol)
THEN
605 CALL add_hill_single(hills_env, colvars, ww, hills_env%n_hills, n_colvar, invdt=invdt)
607 CALL add_hill_single(hills_env, colvars, ww, hills_env%n_hills, n_colvar)
610 i_hills = i_hills + 1
611 filename = trim(multiple_walkers%walkers_file_name(i))//
"_"// &
613 IF (para_env%is_source())
THEN
614 INQUIRE (file=trim(filename), exist=exist)
616 CALL para_env%bcast(exist)
619 delta_hills = i_hills - 1 - multiple_walkers%walkers_status(i)
620 multiple_walkers%walkers_status(i) = i_hills - 1
622 extension=
".metadynLog")
624 WRITE (iw,
'(T2,A,I0,A,I0,A,I0,A)')
'WALKERS| Walker #', i,
'. Reading [', delta_hills, &
625 '] Hills. Total number of Hills acquired [', multiple_walkers%walkers_status(i),
']'
628 "PRINT%PROGRAM_RUN_INFO")
633 colvars(i)%ss0 = ss0_save(i)
634 colvars(i)%delta_s = delta_s_save(i)
636 DEALLOCATE (ss0_save)
637 DEALLOCATE (delta_s_save)
640 CALL timestop(handle)
657 REAL(kind=
dp),
INTENT(IN) :: ww
658 INTEGER,
INTENT(INOUT) :: n_hills
659 INTEGER,
INTENT(IN) :: n_colvar
660 REAL(kind=
dp),
INTENT(IN),
OPTIONAL :: invdt
662 CHARACTER(len=*),
PARAMETER :: routinen =
'add_hill_single'
666 REAL(kind=
dp),
DIMENSION(:),
POINTER :: tnp
667 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: tmp
669 CALL timeset(routinen, handle)
671 wtcontrol =
PRESENT(invdt)
673 IF (
SIZE(hills_env%ss_history, 2) < n_hills + 1)
THEN
674 ALLOCATE (tmp(n_colvar, n_hills + 100))
675 tmp(:, :n_hills) = hills_env%ss_history
676 tmp(:, n_hills + 1:) = 0.0_dp
677 DEALLOCATE (hills_env%ss_history)
678 hills_env%ss_history => tmp
681 IF (
SIZE(hills_env%delta_s_history, 2) < n_hills + 1)
THEN
682 ALLOCATE (tmp(n_colvar, n_hills + 100))
683 tmp(:, :n_hills) = hills_env%delta_s_history
684 tmp(:, n_hills + 1:) = 0.0_dp
685 DEALLOCATE (hills_env%delta_s_history)
686 hills_env%delta_s_history => tmp
689 IF (
SIZE(hills_env%ww_history) < n_hills + 1)
THEN
690 ALLOCATE (tnp(n_hills + 100))
691 tnp(1:n_hills) = hills_env%ww_history
692 tnp(n_hills + 1:) = 0.0_dp
693 DEALLOCATE (hills_env%ww_history)
694 hills_env%ww_history => tnp
698 IF (
SIZE(hills_env%invdt_history) < n_hills + 1)
THEN
699 ALLOCATE (tnp(n_hills + 100))
700 tnp(1:n_hills) = hills_env%invdt_history
701 tnp(n_hills + 1:) = 0.0_dp
702 DEALLOCATE (hills_env%invdt_history)
703 hills_env%invdt_history => tnp
707 n_hills = n_hills + 1
710 hills_env%ss_history(i, n_hills) = colvars(i)%ss0
711 hills_env%delta_s_history(i, n_hills) = colvars(i)%delta_s
713 hills_env%ww_history(n_hills) = ww
714 IF (wtcontrol) hills_env%invdt_history(n_hills) = invdt
716 CALL timestop(handle)
734 n_hills, n_colvar, colvars, metadyn_section, invdt_history)
735 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: ss_history, delta_s_history
736 REAL(kind=
dp),
DIMENSION(:),
POINTER :: ww_history
738 INTEGER,
INTENT(IN) :: n_hills, n_colvar
741 REAL(kind=
dp),
DIMENSION(:),
OPTIONAL,
POINTER :: invdt_history
743 CHARACTER(len=*),
PARAMETER :: routinen =
'restart_hills'
745 INTEGER :: handle, i, j, ndum
746 LOGICAL :: explicit, wtcontrol
747 REAL(kind=
dp) :: rval
748 REAL(kind=
dp),
DIMENSION(:),
POINTER :: rvals
751 CALL timeset(routinen, handle)
753 wtcontrol =
PRESENT(invdt_history)
760 DEALLOCATE (ss_history)
761 DEALLOCATE (delta_s_history)
762 DEALLOCATE (ww_history)
764 DEALLOCATE (invdt_history)
767 cpassert(n_hills == ndum)
768 ALLOCATE (ss_history(n_colvar, n_hills))
769 ALLOCATE (delta_s_history(n_colvar, n_hills))
770 ALLOCATE (ww_history(n_hills))
772 ALLOCATE (invdt_history(n_hills))
777 i_rep_val=i, r_vals=rvals)
778 cpassert(
SIZE(rvals) == n_colvar)
779 ss_history(1:n_colvar, i) = rvals
787 cpassert(n_hills == ndum)
790 i_rep_val=i, r_vals=rvals)
791 cpassert(
SIZE(rvals) == n_colvar)
792 delta_s_history(1:n_colvar, i) = rvals
795 CALL cp_warn(__location__, &
796 "Section SPAWNED_HILLS_SCALE is not present! Setting the scales of the "// &
797 "restarted hills according the parameters specified in the input file.")
800 delta_s_history(j, i) = colvars(i)%delta_s
811 cpassert(n_hills == ndum)
814 i_rep_val=i, r_val=rval)
815 cpassert(
SIZE(rvals) == n_colvar)
819 CALL cp_warn(__location__, &
820 "Section SPAWNED_HILLS_HEIGHT is not present! Setting the height of the"// &
821 " restarted hills according the parameters specified in the input file. ")
832 cpassert(n_hills == ndum)
835 i_rep_val=i, r_val=rval)
836 cpassert(
SIZE(rvals) == n_colvar)
837 invdt_history(i) = rval
840 CALL cp_warn(__location__, &
841 "Section SPAWNED_HILLS_INVDT is not present! Restarting from standard"// &
842 " metadynamics run i.e. setting 1/(Delta T) equal to zero. ")
843 invdt_history = 0._dp
847 CALL cp_abort(__location__, &
848 "Found section SPAWNED_HILLS_INVDT while restarting a standard metadynamics run..."// &
849 " Cannot restart metadynamics from well-tempered MetaD runs. ")
854 CALL timestop(handle)
866 INTEGER,
INTENT(OUT) :: iter_nr
868 IF (meta_env%do_multiple_walkers)
THEN
869 iter_nr = meta_env%multiple_walkers%n_hills_local
871 iter_nr = meta_env%hills_env%n_hills
887 REAL(
dp) :: ddp, delta_s, dfunc, diff_ss, dp2, &
891 colvars => meta_env%metavar
893 DO ih = 1,
SIZE(colvars)
894 IF (colvars(ih)%do_wall)
THEN
895 colvars(ih)%epot_walls = 0.0_dp
896 colvars(ih)%ff_walls = 0.0_dp
897 DO iwall = 1,
SIZE(colvars(ih)%walls)
898 SELECT CASE (colvars(ih)%walls(iwall)%id_type)
903 diff_ss = colvars(ih)%ss0 - colvars(ih)%walls(iwall)%pos
904 IF (colvars(ih)%periodic)
THEN
906 diff_ss = sign(1.0_dp, asin(sin(diff_ss)))*acos(cos(diff_ss))
908 efunc = colvars(ih)%walls(iwall)%k_quadratic*diff_ss**2
909 dfunc = 2.0_dp*colvars(ih)%walls(iwall)%k_quadratic*diff_ss
910 SELECT CASE (colvars(ih)%walls(iwall)%id_direction)
912 IF (diff_ss > 0.0_dp)
THEN
913 colvars(ih)%ff_walls = colvars(ih)%ff_walls - dfunc
914 colvars(ih)%epot_walls = colvars(ih)%epot_walls + efunc
917 IF (diff_ss < 0.0_dp)
THEN
918 colvars(ih)%ff_walls = colvars(ih)%ff_walls - dfunc
919 colvars(ih)%epot_walls = colvars(ih)%epot_walls + efunc
923 diff_ss = colvars(ih)%ss0 - colvars(ih)%walls(iwall)%pos0
924 IF (colvars(ih)%periodic)
THEN
926 diff_ss = sign(1.0_dp, asin(sin(diff_ss)))*acos(cos(diff_ss))
928 efunc = colvars(ih)%walls(iwall)%k_quartic*diff_ss*diff_ss**4
929 dfunc = 4.0_dp*colvars(ih)%walls(iwall)%k_quartic*diff_ss**3
930 SELECT CASE (colvars(ih)%walls(iwall)%id_direction)
932 IF (diff_ss > 0.0_dp)
THEN
933 colvars(ih)%ff_walls = colvars(ih)%ff_walls - dfunc
934 colvars(ih)%epot_walls = colvars(ih)%epot_walls + efunc
937 IF (diff_ss < 0.0_dp)
THEN
938 colvars(ih)%ff_walls = colvars(ih)%ff_walls - dfunc
939 colvars(ih)%epot_walls = colvars(ih)%epot_walls + efunc
943 diff_ss = colvars(ih)%ss0 - colvars(ih)%walls(iwall)%pos
944 IF (colvars(ih)%periodic)
THEN
946 diff_ss = sign(1.0_dp, asin(sin(diff_ss)))*acos(cos(diff_ss))
948 ww = colvars(ih)%walls(iwall)%ww_gauss
949 delta_s = colvars(ih)%walls(iwall)%sigma_gauss
950 ddp = (diff_ss)/delta_s
952 efunc = ww*exp(-0.5_dp*dp2)
953 dfunc = -efunc*ddp/delta_s
954 colvars(ih)%ff_walls = colvars(ih)%ff_walls - dfunc
955 colvars(ih)%epot_walls = colvars(ih)%epot_walls + efunc
static GRID_HOST_DEVICE int modulo(int a, int m)
Equivalent of Fortran's MODULO, which always return a positive number. https://gcc....
Utility routines to open and close files. Tracking of preconnections.
subroutine, public open_file(file_name, file_status, file_form, file_action, file_position, file_pad, unit_number, debug, skip_get_unit_number, file_access)
Opens the requested file using a free unit number.
subroutine, public close_file(unit_number, file_status, keep_preconnection)
Close an open file given by its logical unit number. Optionally, keep the file and unit preconnected.
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,...
types that represent a subsys, i.e. a part of the system
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
Defines the basic variable types.
integer, parameter, public dp
integer, parameter, public default_path_length
Machine interface based on Fortran 2003 and POSIX.
subroutine, public m_mov(source, target)
...
Interface to the message passing library MPI.
Definition of physical constants:
real(kind=dp), parameter, public kelvin
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,...
wrapper to abstract the force evaluation of the various methods
stores all the informations relevant to an mpi environment