54#include "../base/base_uses.f90"
59 LOGICAL,
PRIVATE,
PARAMETER :: debug_this_module = .true.
60 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'input_cp2k_md'
77 cpassert(.NOT.
ASSOCIATED(section))
79 description=
"This section defines the whole set of parameters needed perform an MD run.", &
80 n_keywords=13, n_subsections=6, repeats=.false.)
82 NULLIFY (keyword, subsection)
84 description=
"The ensemble/integrator that you want to use for MD propagation", &
85 usage=
"ensemble nve", &
87 enum_c_vals=
s2a(
"NVE",
"NVT",
"NPT_I",
"NPT_F",
"MSST",
"MSST_DAMPED", &
88 "HYDROSTATICSHOCK",
"ISOKIN",
"REFTRAJ",
"LANGEVIN",
"NPE_F", &
89 "NPE_I",
"NVT_ADIABATIC",
"NPT_IA"), &
90 enum_desc=
s2a(
"constant energy (microcanonical)", &
91 "constant temperature and volume (canonical)", &
92 "constant temperature and pressure using an isotropic cell", &
93 "constant temperature and pressure using a flexible cell", &
94 "simulate steady shock (uniaxial)", &
95 "simulate steady shock (uniaxial) with extra viscosity", &
96 "simulate steady shock with hydrostatic pressure", &
97 "constant kinetic energy", &
98 "reading frames from a file called reftraj.xyz (e.g. for property calculation)", &
99 "langevin dynamics (constant temperature)", &
100 "constant pressure ensemble (no thermostat)", &
101 "constant pressure ensemble using an isotropic cell (no thermostat)", &
102 "adiabatic dynamics in constant temperature and volume ensemble (CAFES)", &
103 "NPT_I ensemble with frozen atoms in absolute coordinate"), &
113 description=
"The number of MD steps to perform, counting from step_start_val. ", &
114 usage=
"STEPS 100", default_i_val=3)
119 description=
"The number of MD steps to perform, counting from step 1", &
120 usage=
"MAX_STEPS 100", default_i_val=1000000000)
125 description=
"The length of an integration step (in case RESPA the large TIMESTEP)", &
126 usage=
"TIMESTEP 1.0", default_r_val=
cp_unit_to_cp2k(
value=0.5_dp, unit_str=
"fs"), &
131 CALL keyword_create(keyword, __location__, name=
"STEP_START_VAL", &
132 description=
"The starting step value for the MD", usage=
"STEP_START_VAL <integer>", &
137 CALL keyword_create(keyword, __location__, name=
"TIME_START_VAL", &
138 description=
"The starting timer value for the MD", &
139 usage=
"TIME_START_VAL <real>", default_r_val=
cp_unit_to_cp2k(
value=0.0_dp, unit_str=
"fs"), &
144 CALL keyword_create(keyword, __location__, name=
"ECONS_START_VAL", &
145 description=
"The starting value of the conserved quantity", &
146 usage=
"ECONS_START_VAL <real>", default_r_val=0.0_dp, &
152 description=
"The temperature in K used to initialize "// &
153 "the velocities with init and pos restart, and in the NPT/NVT simulations", &
154 usage=
"TEMPERATURE 325.0", default_r_val=
cp_unit_to_cp2k(
value=300.0_dp, unit_str=
"K"), &
160 variants=
s2a(
"temp_to",
"temperature_tolerance"), &
161 description=
"The maximum accepted deviation of the (global) temperature "// &
162 "from the desired target temperature before a rescaling of the velocities "// &
163 "is performed. If it is 0 no rescaling is performed. NOTE: This keyword is "// &
164 "obsolescent; Using a CSVR thermostat with a short timeconstant is "// &
165 "recommended as a better alternative.", &
166 usage=
"TEMP_TOL 0.0", default_r_val=0.0_dp, unit_str=
'K')
171 description=
"Compute the temperature per each kind separately", &
172 usage=
"TEMP_KIND LOGICAL", &
173 default_l_val=.false., lone_keyword_l_val=.true.)
177 CALL keyword_create(keyword, __location__, name=
"SCALE_TEMP_KIND", &
178 description=
"When necessary rescale the temperature per each kind separately", &
179 usage=
"SCALE_TEMP_KIND LOGICAL", &
180 default_l_val=.false., lone_keyword_l_val=.true.)
185 description=
"The maximum accepted velocity of the center of mass. "// &
186 "With Shell-Model, comvel may drift if MD%THERMOSTAT%REGION /= GLOBAL ", &
187 usage=
"COMVEL_TOL 0.1", type_of_var=
real_t, n_var=1, unit_str=
"bohr*au_t^-1")
192 description=
"The maximum accepted angular velocity. This option is ignored "// &
193 "when the system is periodic. Removes the components of the velocities that "// &
194 "project on the external rotational degrees of freedom.", &
195 usage=
"ANGVEL_TOL 0.1", type_of_var=
real_t, n_var=1, unit_str=
"bohr*au_t^-1")
200 description=
"Set the initial angular velocity to zero. This option is ignored "// &
201 "when the system is periodic or when initial velocities are defined. Technically, "// &
202 "the part of the random initial velocities that projects on the external "// &
203 "rotational degrees of freedom is subtracted.", &
204 usage=
"ANGVEL_ZERO LOGICAL", &
205 default_l_val=.false., lone_keyword_l_val=.true.)
210 description=
"Specifies the rescaling factor for annealing velocities. "// &
211 "Automatically enables the annealing procedure. This scheme works only for ensembles "// &
212 "that do not have thermostats on particles.", &
213 usage=
"annealing <REAL>", default_r_val=1.0_dp)
217 CALL keyword_create(keyword, __location__, name=
"ANNEALING_CELL", &
218 description=
"Specifies the rescaling factor for annealing velocities of the CELL "// &
219 "Automatically enables the annealing procedure for the CELL. This scheme works only "// &
220 "for ensambles that do not have thermostat on CELLS velocities.", &
221 usage=
"ANNEALING_CELL <REAL>", default_r_val=1.0_dp)
225 CALL keyword_create(keyword, __location__, name=
"TEMPERATURE_ANNEALING", &
226 description=
"Specifies the rescaling factor for the external temperature. "// &
227 "This scheme works only for the Langevin ensemble.", &
228 usage=
"TEMPERATURE_ANNEALING <REAL>", default_r_val=1.0_dp)
232 CALL keyword_create(keyword, __location__, name=
"DISPLACEMENT_TOL", &
233 description=
"This keyword sets a maximum atomic displacement "// &
234 "in each Cartesian direction. "// &
235 "The maximum velocity is evaluated and if it is too large to remain "// &
236 "within the assigned limit, the time step is rescaled accordingly, "// &
237 "and the first half step of the velocity verlet is repeated.", &
238 usage=
"DISPLACEMENT_TOL <REAL>", default_r_val=100.0_dp, &
243 CALL keyword_create(keyword, __location__, name=
"INITIALIZATION_METHOD", &
244 description=
"This keyword selects which method to use to initialize MD. "// &
245 "If velecities are not set explicitly, DEFAULT optioin will assign "// &
246 "random velocities and then scale according to TEMPERATURE; VIBRATIONAL "// &
247 "option will then use previously calculated vibrational modes to "// &
248 "initialise both the atomic positions and velocities so that the "// &
249 "starting point for MD is as close to canonical ensemble as possible, "// &
250 "without the need for lengthy equilibration steps. See PRL 96, 115504 "// &
251 "(2006). The user input atomic positions in this case are expected to "// &
252 "be already geometry optimised. Further options for VIBRATIONAL mode "// &
253 "is can be set in INITIAL_VIBRATION subsection. If unspecified, then "// &
254 "the DEFAULT mode will be used.", &
255 usage=
"INITIALIZATION_METHOD DEFAULT", &
257 enum_c_vals=
s2a(
"DEFAULT",
"VIBRATIONAL"), &
258 enum_desc=
s2a(
"Assign random velocities and then scale according to "// &
260 "Initialise positions and velocities to give canonical ensemble "// &
261 "with TEMPERATURE, using the method described in PRL 96, 115504 (2006)"), &
266 CALL create_langevin_section(subsection)
270 CALL create_msst_section(subsection)
282 CALL create_respa_section(subsection)
286 CALL create_shell_section(subsection)
290 CALL create_adiabatic_section(subsection)
294 CALL create_softening_section(subsection)
298 CALL create_reftraj_section(subsection)
302 CALL create_avgs_section(subsection)
306 CALL create_thermal_region_section(subsection)
310 CALL create_md_print_section(subsection)
314 CALL create_cascade_section(subsection)
318 CALL create_vib_init_section(subsection)
329 SUBROUTINE create_langevin_section(section)
334 cpassert(.NOT.
ASSOCIATED(section))
336 description=
"Controls the set of parameters to run a Langevin MD. "// &
337 "The integrator used follows that given in the article by Ricci et al. "// &
338 "The user can define regions in the system where the atoms inside "// &
339 "undergoes Langevin MD, while those outside the regions undergoes "// &
340 "NVE Born Oppenheimer MD. To define the regions, the user should "// &
341 "use THERMAL_REGION subsection of MOTION%MD. ", &
343 n_keywords=0, n_subsections=1, repeats=.false.)
347 description=
"Gamma parameter for the Langevin dynamics (LD)", &
348 usage=
"gamma 0.001", &
349 default_r_val=0.0_dp, unit_str=
'fs^-1')
354 variants=[
"NoisyGamma"], &
355 description=
"Imaginary Langevin Friction term for LD with noisy forces.", &
357 usage=
"Noisy_Gamma 4.0E-5", default_r_val=0.0_dp, unit_str=
'fs^-1')
362 variants=[
"ShadowGamma"], &
363 description=
"Shadow Langevin Friction term for LD with noisy forces in order to adjust Noisy_Gamma.", &
365 usage=
"Shadow_Gamma 0.001", default_r_val=0.0_dp, unit_str=
'fs^-1')
368 END SUBROUTINE create_langevin_section
375 SUBROUTINE create_md_print_section(section)
381 cpassert(.NOT.
ASSOCIATED(section))
383 description=
"Controls the printing properties during an MD run", &
384 n_keywords=0, n_subsections=1, repeats=.false.)
385 NULLIFY (print_key, keyword)
388 description=
"Print the output and restart file if walltime is reached or "// &
389 "if an external EXIT command is given. It still requires the keyword LAST "// &
390 "to be present for the specific print key (in case the last step should not "// &
391 "match with the print_key iteration number).", &
392 usage=
"FORCE_LAST LOGICAL", &
393 default_l_val=.false., lone_keyword_l_val=.true.)
398 description=
"Controls the output the ener file", &
405 description=
"Controls the output of the shell-energy file (only if shell-model)", &
412 description=
"Controls the output of the temperature"// &
413 " computed separately for each kind", &
420 description=
"Controls the output of the temperature of the"// &
421 " shell-core motion computed separately for each kind", &
428 description=
"Controls the printing of COM velocity during an MD", &
430 filename=
"__STD_OUT__")
435 description=
"Controls the printing of coefficients during an MD run.", &
442 description=
"Controls the printing basic info during the calculation of the "// &
443 "translational/rotational degrees of freedom.", print_level=
low_print_level, &
446 description=
"Prints atomic coordinates in the standard orientation. "// &
447 "Coordinates are not affected during the calculation.", &
448 default_l_val=.false., lone_keyword_l_val=.true.)
455 description=
"Controls the printing of basic and summary information during the"// &
456 " Molecular Dynamics", &
460 END SUBROUTINE create_md_print_section
467 SUBROUTINE create_respa_section(section)
472 cpassert(.NOT.
ASSOCIATED(section))
475 description=
"Multiple timestep integration based on RESPA (implemented for NVE only)."// &
476 " RESPA exploits multiple force_eval."// &
477 " In this case the order of the force_eval maps"// &
478 " the order of the respa shells from the slowest to the fastest force evaluation."// &
479 " If force_evals share the same subsys, it's enough then to specify the"// &
480 " subsys in the force_eval corresponding at the first index in the multiple_force_eval list."// &
481 " Can be used to speedup classical and ab initio MD simulations.", &
482 n_keywords=1, n_subsections=0, repeats=.false., &
487 description=
"The number of reference MD steps between two RESPA corrections.", &
488 usage=
"FREQUENCY <INTEGER>", default_i_val=5)
492 END SUBROUTINE create_respa_section
499 SUBROUTINE create_reftraj_section(section)
505 cpassert(.NOT.
ASSOCIATED(section))
508 description=
"Loads an external trajectory file and performs analysis on the"// &
509 " loaded snapshots.", &
510 n_keywords=1, n_subsections=1, repeats=.false.)
512 NULLIFY (keyword, print_key, subsection)
514 CALL keyword_create(keyword, __location__, name=
"TRAJ_FILE_NAME", &
515 description=
"Specify the filename where the trajectory is stored. "// &
516 "If you built your own trajectory file make sure it has the trajectory format. "// &
517 'In particular, each structure has to be enumerated using " i = ..."', &
519 usage=
"TRAJ_FILE_NAME <CHARACTER>", default_lc_val=
"reftraj.xyz")
523 CALL keyword_create(keyword, __location__, name=
"CELL_FILE_NAME", &
524 description=
"Specify the filename where the cell is stored "// &
525 "(for trajectories generated within variable cell ensembles).", repeats=.false., &
526 usage=
"CELL_FILE_NAME <CHARACTER>", default_lc_val=
"reftraj.cell")
531 keyword, __location__, name=
"VARIABLE_VOLUME", &
532 description=
"Enables the possibility to read a CELL file with information on the CELL size during the MD.", &
533 repeats=.false., default_l_val=.false., lone_keyword_l_val=.true.)
537 CALL keyword_create(keyword, __location__, name=
"FIRST_SNAPSHOT", &
538 description=
"Index of the snapshot stored in the trajectory file "// &
539 "from which to start a REFTRAJ run", &
540 repeats=.false., usage=
"FIRST_SNAPSHOT <INTEGER>", default_i_val=1)
544 CALL keyword_create(keyword, __location__, name=
"LAST_SNAPSHOT", &
545 description=
"Index of the last snapshot stored in the trajectory file "// &
546 "that is read along a REFTRAJ run. Must be specified as default is 0. "// &
547 "Must be specified exactly as LAST_SNAPSHOT = FIRST_SNAPSHOT + STRIDE*(number of strides) "// &
548 "to avoid an error 'Unexpected EOF'. Note that STRIDE*(number of strides) "// &
549 "is simply the number of steps between the first to last snapshot.", &
550 repeats=.false., usage=
"LAST_SNAPSHOT", default_i_val=0)
555 description=
" Stride in number of snapshot for the reftraj analysis", &
556 repeats=.false., usage=
"STRIDE", default_i_val=1)
561 description=
"Selects the properties to evaluate for each retrieved snapshot during a REFTRAJ run", &
564 enum_c_vals=
s2a(
"NONE",
"ENERGY",
"ENERGY_FORCES"), &
565 enum_desc=
s2a(
"Evaluate nothing",
"Evaluate only the energy",
"Evaluate energy and forces"))
569 CALL keyword_create(keyword, __location__, name=
"eval_energy_forces", &
570 description=
"Evaluate energy and forces for each retrieved snapshot during a REFTRAJ run", &
571 repeats=.false., default_l_val=.false., lone_keyword_l_val=.true., &
572 deprecation_notice=
"Please use MOTION / MD / REFTRAJ / EVAL instead.")
577 description=
"Evaluate the forces for each retrieved snapshot during a REFTRAJ run", &
578 repeats=.false., default_l_val=.false., lone_keyword_l_val=.true., &
579 deprecation_notice=
"Please use MOTION / MD / REFTRAJ / EVAL instead.")
583 CALL create_msd_section(subsection)
588 description=
"The section that controls the output of a reftraj run", &
589 n_keywords=1, n_subsections=0, repeats=.false.)
593 description=
"Controls the output of msd per kind", &
600 description=
"Controls the output of msd per molecule kind", &
607 description=
"Controls the output of index and dislacement of "// &
608 "atoms that moved away from the initial position of more than a "// &
609 "given distance (see msd%disp_tol)", &
618 END SUBROUTINE create_reftraj_section
625 SUBROUTINE create_msd_section(section)
631 cpassert(.NOT.
ASSOCIATED(section))
634 description=
"Loads an external trajectory file and performs analysis on the"// &
635 " loaded snapshots.", &
636 n_keywords=3, n_subsections=0, repeats=.false.)
638 NULLIFY (keyword, subsection)
640 CALL keyword_create(keyword, __location__, name=
"_SECTION_PARAMETERS_", &
641 description=
"controls the activation of core-level spectroscopy simulations", &
643 default_l_val=.false., &
644 lone_keyword_l_val=.true.)
648 CALL keyword_create(keyword, __location__, name=
"REF0_FILENAME", &
649 description=
"Specify the filename where the initial reference configuration is stored.", &
650 repeats=.false., usage=
"REF0_FILENAME <CHARACTER>", default_lc_val=
"")
655 description=
"Set up the calculation of the MSD for each atomic kind", &
656 usage=
"MSD_PER_KIND <LOGICAL>", repeats=.false., &
657 default_l_val=.false., lone_keyword_l_val=.true.)
661 CALL keyword_create(keyword, __location__, name=
"MSD_PER_MOLKIND", &
662 description=
"Set up the calculation of the MSD for each molecule kind. "// &
663 "The position of the center of mass of the molecule is considered.", &
664 usage=
"MSD_PER_MOLKIND <LOGICAL>", repeats=.false., &
665 default_l_val=.false., lone_keyword_l_val=.true.)
669 CALL keyword_create(keyword, __location__, name=
"MSD_PER_REGION", &
670 description=
"Set up the calculation of the MSD for each defined region.", &
671 usage=
"MSD_PER_REGION <LOGICAL>", repeats=.false., &
672 default_l_val=.false., lone_keyword_l_val=.true.)
680 CALL keyword_create(keyword, __location__, name=
"DISPLACED_ATOM", &
681 description=
"Identify the atoms that moved from their initial "// &
682 "position of a distance larger than a given tolerance (see msd%displacement_tol).", &
683 usage=
"DISPLACED_ATOM <LOGICAL>", repeats=.false., &
684 default_l_val=.false., lone_keyword_l_val=.true.)
688 CALL keyword_create(keyword, __location__, name=
"displacement_tol", &
689 description=
"Lower limit to define displaced atoms", &
690 usage=
"DISPLACEMENT_TOL real", &
691 default_r_val=0._dp, n_var=1, unit_str=
'bohr')
695 END SUBROUTINE create_msd_section
702 SUBROUTINE create_msst_section(section)
707 cpassert(.NOT.
ASSOCIATED(section))
710 description=
"Parameters for Multi-Scale Shock Technique (MSST) "// &
711 "which simulate the effect of a steady planar shock on a unit cell. "// &
712 "Reed et. al. Physical Review Letters 90, 235503 (2003).", &
713 n_keywords=1, n_subsections=0, repeats=.false.)
717 description=
"Initial pressure", &
718 usage=
"PRESSURE real", &
719 default_r_val=0._dp, n_var=1, unit_str=
'bar')
724 description=
"Initial energy", &
725 usage=
"ENERGY real", &
726 default_r_val=0._dp, n_var=1, unit_str=
'hartree')
731 description=
"Initial volume", &
732 usage=
"VOLUME real", &
733 default_r_val=0._dp, n_var=1, unit_str=
'angstrom^3')
738 description=
"Effective cell mass", &
739 usage=
"CMASS real", &
740 default_r_val=0._dp, n_var=1, unit_str=
'au_m')
744 CALL keyword_create(keyword, __location__, name=
"VSHOCK", variants=[
"V_SHOCK"], &
745 description=
"Velocity shock", &
746 usage=
"VSHOCK real", &
747 default_r_val=0._dp, n_var=1, unit_str=
'm/s')
752 description=
"Damping coefficient for cell volume", &
753 usage=
"GAMMA real", &
755 default_r_val=0.0_dp)
759 END SUBROUTINE create_msst_section
765 SUBROUTINE create_shell_section(section)
771 cpassert(.NOT.
ASSOCIATED(section))
774 description=
"Parameters of shell model in adiabatic dynamics.", &
775 n_keywords=4, n_subsections=1, repeats=.false.)
777 NULLIFY (keyword, thermo_section)
780 description=
"Temperature in K used to control "// &
781 "the internal velocities of the core-shell motion ", &
782 usage=
"temperature 5.0", &
789 description=
"Maximum accepted temperature deviation"// &
790 " from the expected value, for the internal core-shell motion."// &
791 " If 0, no rescaling is performed", &
792 usage=
"temp_tol 0.0", default_r_val=0.0_dp, unit_str=
'K')
796 CALL keyword_create(keyword, __location__, name=
"nose_particle", &
797 description=
"If nvt or npt, the core and shell velocities are controlled "// &
798 "by the same thermostat used for the particle. This might favour heat exchange "// &
799 "and additional rescaling of the internal core-shell velocity is needed (TEMP_TOL)", &
800 default_l_val=.false., lone_keyword_l_val=.true.)
804 CALL keyword_create(keyword, __location__, name=
"DISPLACEMENT_SHELL_TOL", &
805 description=
"This keyword sets a maximum variation of the shell"// &
806 " core distance in each Cartesian direction."// &
807 " The maximum internal core-shell velocity is evaluated and"// &
808 " if it is too large to remain"// &
809 " within the assigned limit, the time step is rescaled accordingly,"// &
810 " and the first half step of the velocity verlet is repeated.", &
811 usage=
"DISPLACEMENT_SHELL_TOL <REAL>", default_r_val=100.0_dp, &
820 END SUBROUTINE create_shell_section
826 SUBROUTINE create_adiabatic_section(section)
830 TYPE(
section_type),
POINTER :: thermo_fast_section, thermo_slow_section
832 cpassert(.NOT.
ASSOCIATED(section))
834 CALL section_create(section, __location__, name=
"ADIABATIC_DYNAMICS", &
835 description=
"Parameters used in canonical adiabatic free energy sampling (CAFES).", &
836 n_keywords=5, n_subsections=2, repeats=.false., &
839 NULLIFY (keyword, thermo_fast_section, thermo_slow_section)
842 description=
"Temperature in K used to control "// &
843 "the fast degrees of freedom ", &
844 usage=
"temp_fast 5.0", &
851 description=
"Temperature in K used to control "// &
852 "the slow degrees of freedom ", &
853 usage=
"temp_slow 5.0", &
859 CALL keyword_create(keyword, __location__, name=
"temp_tol_fast", &
860 description=
"Maximum accepted temperature deviation"// &
861 " from the expected value, for the fast motion."// &
862 " If 0, no rescaling is performed", &
863 usage=
"temp_tol_fast 0.0", default_r_val=0.0_dp, unit_str=
'K')
867 CALL keyword_create(keyword, __location__, name=
"temp_tol_slow", &
868 description=
"Maximum accepted temperature deviation"// &
869 " from the expected value, for the slow motion."// &
870 " If 0, no rescaling is performed", &
871 usage=
"temp_tol_slow 0.0", default_r_val=0.0_dp, unit_str=
'K')
876 description=
"number of respa steps for fast degrees of freedom", &
877 repeats=.false., default_i_val=1)
889 END SUBROUTINE create_adiabatic_section
896 SUBROUTINE create_softening_section(section)
901 CALL section_create(section, __location__, name=
"VELOCITY_SOFTENING", &
902 description=
"A method to initialize the velocities along low-curvature " &
903 //
"directions in order to favors MD trajectories to cross rapidly over " &
904 //
"small energy barriers into neighboring basins. " &
905 //
"In each iteration the forces are calculated at a point y, which " &
906 //
"is slightly displaced from the current positions x in the direction " &
907 //
"of the original velocities v. The velocities are then updated with " &
908 //
"the force component F_t, which is perpendicular to N. " &
909 //
"N = v / |v|; y = x + delta * N; F_t = F(y) - ⟨ F(y) | N ⟩ * N; " &
910 //
"v' = v + alpha * F_t")
914 description=
"Number of softening iterations performed. " &
915 //
"Typical values are around 40 steps.", &
921 description=
"Displacement used to obtain y.", &
922 default_r_val=0.1_dp)
927 description=
"Mixing factor used for updating velocities.", &
928 default_r_val=0.15_dp)
932 END SUBROUTINE create_softening_section
942 SUBROUTINE create_thermal_region_section(section)
946 TYPE(
section_type),
POINTER :: print_key, region_section, subsection, &
949 cpassert(.NOT.
ASSOCIATED(section))
951 CALL section_create(section, __location__, name=
"THERMAL_REGION", &
952 description=
"Define regions where different initialization and control "// &
953 "of the temperature is used. When MOTION%MD%ENSEMBLE is set to LANGEVIN, "// &
954 "this section controls if the atoms defined inside and outside the "// &
955 "thermal regions should undergo Langevin MD or NVE Born-Oppenheimer MD. "// &
956 "The theory behind Langevin MD using different regions can be found in "// &
957 "articles by Kantorovitch et al. listed below.", &
959 n_keywords=0, n_subsections=1, repeats=.false.)
961 NULLIFY (region_section)
962 NULLIFY (keyword, subsection, tfunc_section)
964 CALL keyword_create(keyword, __location__, name=
"FORCE_RESCALING", &
965 description=
"Control the rescaling ot the velocities in all the regions, "// &
966 "according to the temperature assigned to each reagion, when "// &
967 "RESTART_VELOCITY in EXT_RESTART is active.", &
968 default_l_val=.false., lone_keyword_l_val=.true.)
972 CALL keyword_create(keyword, __location__, name=
"DO_LANGEVIN_DEFAULT", &
973 description=
"If ENSEMBLE is set to LANGEVIN, controls whether the "// &
974 "atoms NOT defined in the thermal regions to undergo langevin MD "// &
975 "or not. If not, then the atoms will undergo NVE Born-Oppenheimer MD.", &
976 usage=
"DO_LANGEVIN_DEFAULT .FALSE.", &
977 default_l_val=.false., lone_keyword_l_val=.true.)
981 CALL section_create(region_section, __location__, name=
"DEFINE_REGION", &
982 description=
"Define arbitrary thermal region with distinct "// &
983 "temperature controls. For thermostat regions, see section "// &
984 "MOTION%MD%THERMOSTAT%DEFINE_REGION; when the keyword "// &
985 "MOTION%MD%THERMOSTAT%REGION is set to THERMAL, the "// &
986 "definition of thermostat regions will be the same as "// &
987 "that of thermal regions.", &
988 n_keywords=3, n_subsections=1, repeats=.true.)
992 description=
"Specifies a list of atoms belonging to the region.", &
993 usage=
"LIST {integer} {integer} .. {integer}", &
994 repeats=.true., n_var=-1, type_of_var=
integer_t)
999 variants=[
"SEGNAME"], &
1000 description=
"Specifies the name of the molecules belonging to the region.", &
1001 usage=
"MOLNAME WAT MEOH", repeats=.true., &
1002 n_var=-1, type_of_var=
char_t)
1007 variants=[
"PROTEIN"], &
1008 description=
"Not yet implemented: In a QM/MM run all "// &
1009 "MM atoms are specified as a whole region", &
1010 usage=
"MM_SUBSYS (NONE|ATOMIC|MOLECULAR)", &
1011 enum_c_vals=
s2a(
"NONE",
"ATOMIC",
"MOLECULAR"), &
1013 enum_desc=
s2a(
"Nothing in the region", &
1014 "Only the MM atoms itself", &
1015 "The full molecule/residue that contains a MM atom"), &
1021 description=
"Not yet implemented: In a QM/MM run all "// &
1022 "QM atoms are specified as a whole region", &
1023 usage=
"QM_SUBSYS (NONE|ATOMIC|MOLECULAR)", &
1024 enum_c_vals=
s2a(
"NONE",
"ATOMIC",
"MOLECULAR"), &
1025 enum_desc=
s2a(
"Nothing in the region", &
1026 "Only the QM atoms itself", &
1027 "The full molecule/residue that contains a QM atom"), &
1034 description=
"The temperature in K used to initialize the velocities "// &
1035 "of the atoms in this region. It overrides the value evaluated by the "// &
1036 "function defined in TFUNC section.", &
1037 usage=
"TEMPERATURE 5.0", &
1044 description=
"Maximum accepted temperature deviation from the target "// &
1045 "value for this region, which is defined by TEMPERATURE or evaluated "// &
1046 "with the function defined in TFUNC section. Default value is 0 that "// &
1047 "turns off rescaling.", &
1048 usage=
"TEMP_TOL 0.0", &
1049 default_r_val=0.0_dp, unit_str=
'K')
1054 description=
"When ENSEMBLE is set to LANGEVIN, Controls whether "// &
1055 "the atoms in the thermal region should undergo Langevin MD. If "// &
1056 "not, then they will undergo NVE Born-Oppenheimer MD.", &
1057 usage=
"DO_LANGEVIN .TRUE.", &
1058 default_l_val=.true., lone_keyword_l_val=.true.)
1062 CALL keyword_create(keyword, __location__, name=
"NOISY_GAMMA_REGION", &
1063 description=
"Special imaginary Langevin Friction term"// &
1064 " for Langevin Dynamics with noisy forces for the atoms in this region."// &
1065 " When set, overrides the general value set by NOISY_GAMMA in the MOTION%MD%LANGEVIN section."// &
1066 " When unset for a defined region, the general NOISY_GAMMA value applies.", &
1067 citations=[
kuhne2007], usage=
"NOISY_GAMMA_REGION 4.0E-5", &
1075 description=
"This section specifies the target temperature by using "// &
1076 "the time step and the temperature of the last step, allowing for "// &
1077 "simple region-wise simulated annealing (heating/cooling) controls. "// &
1078 "When the absolute difference between target and actual temperature "// &
1079 "is larger than TEMP_TOL value, the atomic velocities in this region "// &
1080 "are rescaled to match the target temperature. However, note that "// &
1081 "the actual temperature is subject to thermostat coupling as well "// &
1082 "and a small restrictive TEMP_TOL together with a rapidly varying "// &
1083 "target temperature function leaves little time for the system to "// &
1084 "adjust itself, possibly showing unphysical relaxation behaviors. "// &
1086 n_keywords=4, n_subsections=0, repeats=.false.)
1090 description=
"Specifies the functional form of the temperature "// &
1091 "in Kelvin in mathematical notation. Two variables are available "// &
1092 "with one denoted S for the integer time step and another denoted "// &
1093 "T for the temperature of the last step in Kelvin respectively. When "// &
1094 "DEFINE_REGION%TEMPERATURE is not explicitly set, the function "// &
1095 "will be evaluated at S set in MD%STEP_START_VAL and T = 0.0 "// &
1096 "to initialize the velocities of atoms in the region.", &
1097 usage=
"FUNCTION A*T+(1-A)*B*(2/(1+EXP(-0.001*S))-1)", &
1103 description=
"Defines the parameters of the functional form", &
1104 usage=
"PARAMETERS A B", type_of_var=
char_t, &
1105 n_var=-1, repeats=.true.)
1110 description=
"Optionally, allows to define valid CP2K unit strings for each parameter value. "// &
1111 "It is assumed that the corresponding parameter value is specified in this unit.", &
1112 usage=
"UNITS ", type_of_var=
char_t, &
1113 n_var=-1, repeats=.true.)
1118 description=
"Defines the values of the parameter of the functional form", &
1119 usage=
"VALUES 0.75 3.00E+2", type_of_var=
real_t, &
1120 n_var=-1, repeats=.true., unit_str=
"internal_cp2k")
1133 description=
"Collects all print_keys for thermal_regions", &
1134 n_keywords=1, n_subsections=0, repeats=.false.)
1137 description=
"Controls the output of the temperature per "// &
1138 "region. The first three columns are: step number, time "// &
1139 "(in fs), and the temperature of atoms not in any temperature "// &
1140 "region (or zero if all atoms are in some temperature region). "// &
1141 "Subsequent columns are, in pairs, the target value "// &
1142 "(calculated by TFUNC or specified by TEMPERATURE) "// &
1143 "and the actual value of the temperature for each region.", &
1150 description=
"Controls output of information on which atoms "// &
1151 "underwent Langevin MD and which atoms did not.", &
1160 END SUBROUTINE create_thermal_region_section
1169 SUBROUTINE create_cascade_section(section)
1177 NULLIFY (subsection)
1178 cpassert(.NOT.
ASSOCIATED(section))
1181 description=
"Defines the parameters for the setup of a cascade simulation.", &
1187 name=
"_SECTION_PARAMETERS_", &
1188 description=
"Controls the activation of the CASCADE section.", &
1189 usage=
"&CASCADE on", &
1190 default_l_val=.false., &
1191 lone_keyword_l_val=.true.)
1196 description=
"Total energy transferred to the system during the cascade event.", &
1197 usage=
"ENERGY 20.0", &
1198 default_r_val=0.0_dp, &
1203 CALL section_create(subsection, __location__, name=
"ATOM_LIST", &
1204 description=
"Defines a list of atoms for which the initial velocities are modified", &
1209 CALL keyword_create(keyword, __location__, name=
"_DEFAULT_KEYWORD_", &
1210 description=
"Defines the list of atoms for which the velocities are modified. "// &
1211 "Each record consists of the atomic index, the velocity vector, and "// &
1212 "a weight to define which fraction of the total energy is assigned "// &
1213 "to the current atom: `Atomic_index v_x v_y v_z Weight`", &
1214 usage=
"{{Integer} {Real} {Real} {Real} {Real}}", &
1223 END SUBROUTINE create_cascade_section
1230 SUBROUTINE create_avgs_section(section)
1236 cpassert(.NOT.
ASSOCIATED(section))
1238 description=
"Controls the calculation of the averages during an MD run.", &
1239 n_keywords=1, n_subsections=1, repeats=.false.)
1240 NULLIFY (keyword, print_key, subsection)
1242 CALL keyword_create(keyword, __location__, name=
"_SECTION_PARAMETERS_", &
1243 description=
"Controls the calculations of the averages.", &
1244 usage=
"&AVERAGES T", default_l_val=.true., lone_keyword_l_val=.true.)
1248 CALL keyword_create(keyword, __location__, name=
"ACQUISITION_START_TIME", &
1249 description=
"Setup up the simulation time when the acquisition process to compute"// &
1250 " averages is started.", &
1251 usage=
"ACQUISITION_START_TIME <REAL>", &
1252 default_r_val=0.0_dp, unit_str=
'fs')
1256 CALL keyword_create(keyword, __location__, name=
"AVERAGE_COLVAR", &
1257 description=
"Switch for computing the averages of COLVARs.", &
1258 usage=
"AVERAGE_COLVAR <LOGICAL>", default_l_val=.false., &
1259 lone_keyword_l_val=.true.)
1264 description=
"Controls the output the averaged quantities", &
1270 CALL create_avgs_restart_section(subsection)
1273 END SUBROUTINE create_avgs_section
1280 SUBROUTINE create_avgs_restart_section(section)
1285 cpassert(.NOT.
ASSOCIATED(section))
1286 CALL section_create(section, __location__, name=
"RESTART_AVERAGES", &
1287 description=
"Stores information for restarting averages.", &
1288 n_keywords=1, n_subsections=1, repeats=.false.)
1291 CALL keyword_create(keyword, __location__, name=
"ITIMES_START", &
1292 description=
"TIME STEP starting the evaluation of averages", &
1293 usage=
"ITIMES_START <INTEGER>", type_of_var=
integer_t, n_var=1)
1298 description=
"CPU average", usage=
"AVECPU <REAL>", &
1299 type_of_var=
real_t, n_var=1)
1304 description=
"HUGONIOT average", usage=
"AVEHUGONIOT <REAL>", &
1305 type_of_var=
real_t, n_var=1)
1309 CALL keyword_create(keyword, __location__, name=
"AVETEMP_BARO", &
1310 description=
"BAROSTAT TEMPERATURE average", usage=
"AVETEMP_BARO <REAL>", &
1311 type_of_var=
real_t, n_var=1)
1316 description=
"POTENTIAL ENERGY average", usage=
"AVEPOT <REAL>", &
1317 type_of_var=
real_t, n_var=1)
1322 description=
"KINETIC ENERGY average", usage=
"AVEKIN <REAL>", &
1323 type_of_var=
real_t, n_var=1)
1328 description=
"TEMPERATURE average", usage=
"AVETEMP <REAL>", &
1329 type_of_var=
real_t, n_var=1)
1334 description=
"QM KINETIC ENERGY average in QMMM runs", usage=
"AVEKIN_QM <REAL>", &
1335 type_of_var=
real_t, n_var=1)
1340 description=
"QM TEMPERATURE average in QMMM runs", usage=
"AVETEMP_QM <REAL>", &
1341 type_of_var=
real_t, n_var=1)
1346 description=
"VOLUME average", usage=
"AVEVOL <REAL>", &
1347 type_of_var=
real_t, n_var=1)
1352 description=
"CELL VECTOR A average", usage=
"AVECELL_A <REAL>", &
1353 type_of_var=
real_t, n_var=1)
1358 description=
"CELL VECTOR B average", usage=
"AVECELL_B <REAL>", &
1359 type_of_var=
real_t, n_var=1)
1364 description=
"CELL VECTOR C average", usage=
"AVECELL_C <REAL>", &
1365 type_of_var=
real_t, n_var=1)
1370 description=
"ALPHA cell angle average", usage=
"AVEALPHA <REAL>", &
1371 type_of_var=
real_t, n_var=1)
1376 description=
"BETA cell angle average", usage=
"AVEBETA <REAL>", &
1377 type_of_var=
real_t, n_var=1)
1382 description=
"GAMMA cell angle average", usage=
"AVEGAMMA <REAL>", &
1383 type_of_var=
real_t, n_var=1)
1388 description=
"CONSTANT ENERGY average", usage=
"AVE_ECONS <REAL>", &
1389 type_of_var=
real_t, n_var=1)
1394 description=
"PRESSURE average", usage=
"AVE_PRESS <REAL>", &
1395 type_of_var=
real_t, n_var=1)
1400 description=
"P_{XX} average", usage=
"AVE_PXX <REAL>", &
1401 type_of_var=
real_t, n_var=1)
1406 description=
"PV VIRIAL average", usage=
"AVE_PV_VIR <REAL> .. <REAL>", &
1407 type_of_var=
real_t, n_var=9)
1412 description=
"PV TOTAL average", usage=
"AVE_PV_TOT <REAL> .. <REAL>", &
1413 type_of_var=
real_t, n_var=9)
1418 description=
"PV KINETIC average", usage=
"AVE_PV_KIN <REAL> .. <REAL>", &
1419 type_of_var=
real_t, n_var=9)
1423 CALL keyword_create(keyword, __location__, name=
"AVE_PV_CNSTR", &
1424 description=
"PV CONSTRAINTS average", usage=
"AVE_PV_CNSTR <REAL> .. <REAL>", &
1425 type_of_var=
real_t, n_var=9)
1430 description=
"PV XC average", usage=
"AVE_PV_XC <REAL> .. <REAL>", &
1431 type_of_var=
real_t, n_var=9)
1435 CALL keyword_create(keyword, __location__, name=
"AVE_PV_FOCK_4C", &
1436 description=
"PV XC average", usage=
"AVE_PV_FOCK_4C <REAL> .. <REAL>", &
1437 type_of_var=
real_t, n_var=9)
1442 description=
"COLVARS averages", usage=
"AVE_COLVARS <REAL> .. <REAL>", &
1443 type_of_var=
real_t, n_var=-1)
1448 description=
"METRIC TENSOR averages", usage=
"AVE_MMATRIX <REAL> .. <REAL>", &
1449 type_of_var=
real_t, n_var=-1)
1452 END SUBROUTINE create_avgs_restart_section
1459 SUBROUTINE create_vib_init_section(section)
1464 cpassert(.NOT.
ASSOCIATED(section))
1465 CALL section_create(section, __location__, name=
"INITIAL_VIBRATION", &
1466 description=
"Controls the set of parameters for MD initialisation "// &
1467 "based on vibration analysis data. The starting atomic coordinates "// &
1468 "should be based on the relaxed positions obtained from a previous "// &
1469 "geometry/cell optimisation calculation, and the vibrational "// &
1470 "frequencies and displacements data should be obtained from a "// &
1471 "vibrational analysis calculation done based on the relaxed "// &
1472 "coordinates. The MD initialisation process expects the user has "// &
1473 "performed both geometry optimisation and vibrational analysis "// &
1474 "before hand, and won't perform those calculations automatically ", &
1476 n_keywords=0, n_subsections=1, repeats=.false.)
1478 CALL keyword_create(keyword, __location__, name=
"VIB_EIGS_FILE_NAME", &
1479 description=
"The file name of vibrational frequency (eigenvalue) "// &
1480 "and displacement (eigenvector) data calculated from the a "// &
1481 "vibrational analysis calculation done previously. This file must "// &
1482 "be the same unformatted binary file as referred to by "// &
1483 "VIBRATIONAL_ANALYSIS%PRINT%CARTESIAN_EIGS keyword. If this keyword "// &
1484 "is not explicitly defined by the user, then the default file "// &
1485 "name of: <project_name>-<CARTESIAN_EIGS_FILENAME>.eig will be used", &
1486 usage=
"VIB_EIGS_FILE_NAME <FILENAME>", &
1492 description=
"Controls the initial ratio of potential and kinetic "// &
1493 "contribution to the total energy. The contribution is determined by "// &
1494 "COS**2(2*pi*PHASE) for potential energy, and SIN**2(2*pi*PHASE) "// &
1495 "for kinetic energy. If PHASE is negative, then for each vibration "// &
1496 "mode the phase is determined randomly. Otherwise, PHASE must be "// &
1497 "between 0.0 and 1.0 and will be the same for all vibration modes. "// &
1498 "If value > 1.0 it will just be treated as 1.0. "// &
1499 "For example, setting PHASE = 0.25 would set all modes to "// &
1500 "begin with entirely kinetic energy --- in other words, the initial "// &
1501 "atomic positions will remain at their optimised location", &
1502 default_r_val=-1.0_dp, &
1503 usage=
"PHASE <REAL>")
1506 END SUBROUTINE create_vib_init_section
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public rengaraj2020
integer, save, public vandevondele2002
integer, save, public west2006
integer, save, public tuckerman1992
integer, save, public kuhne2007
integer, save, public ricci2003
integer, save, public guidon2008
integer, save, public kantorovich2008
integer, save, public minary2003
integer, save, public kantorovich2008a
integer, save, public evans1983
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer, parameter, public debug_print_level
integer, parameter, public low_print_level
integer, parameter, public medium_print_level
integer, parameter, public high_print_level
integer, parameter, public add_last_numeric
subroutine, public cp_print_key_section_create(print_key_section, location, name, description, print_level, each_iter_names, each_iter_values, add_last, filename, common_iter_levels, citations, unit_str)
creates a print_key section
real(kind=dp) function, public cp_unit_to_cp2k(value, unit_str, defaults, power)
converts to the internal cp2k units to the given unit
This public domain function parser module is intended for applications where a set of mathematical ex...
character(len=:) function, allocatable, public docf()
...
Defines the basic variable types.
integer, parameter, public dp
initialization of the reftraj structure used to analyse previously generated trajectories
integer, parameter, public reftraj_eval_energy
integer, parameter, public reftraj_eval_energy_forces
integer, parameter, public reftraj_eval_none
Utilities for string manipulations.
character(len=1), parameter, public newline