Line data Source code
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 Provides an interface to the velocity-verlet based integrator
10 : !> routines for all ensembles
11 : !> \author CJM (11-SEPT-2002)
12 : ! **************************************************************************************************
13 : MODULE velocity_verlet_control
14 :
15 : USE force_env_types, ONLY: force_env_type
16 : USE global_types, ONLY: global_environment_type
17 : USE input_constants, ONLY: &
18 : isokin_ensemble, langevin_ensemble, npe_f_ensemble, npe_i_ensemble, &
19 : nph_uniaxial_damped_ensemble, nph_uniaxial_ensemble, npt_f_ensemble, npt_i_ensemble, &
20 : npt_ia_ensemble, nve_ensemble, nvt_adiabatic_ensemble, nvt_ensemble, reftraj_ensemble
21 : USE integrator, ONLY: &
22 : isokin, langevin, nph_uniaxial, nph_uniaxial_damped, npt_f, npt_i, nve, nve_respa, nvt, &
23 : nvt_adiabatic, reftraj
24 : USE md_environment_types, ONLY: get_md_env,&
25 : md_environment_type
26 : USE simpar_types, ONLY: simpar_type
27 : #include "../base/base_uses.f90"
28 :
29 : IMPLICIT NONE
30 :
31 : PRIVATE
32 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'velocity_verlet_control'
33 : PUBLIC :: velocity_verlet
34 :
35 : CONTAINS
36 :
37 : ! **************************************************************************************************
38 : !> \brief ...
39 : !> \param md_env ...
40 : !> \param globenv ...
41 : !> \par History
42 : !> none
43 : !> \author CJM
44 : ! **************************************************************************************************
45 41595 : SUBROUTINE velocity_verlet(md_env, globenv)
46 :
47 : TYPE(md_environment_type), POINTER :: md_env
48 : TYPE(global_environment_type), POINTER :: globenv
49 :
50 : CHARACTER(LEN=*), PARAMETER :: routineN = 'velocity_verlet'
51 :
52 : INTEGER :: handle
53 : TYPE(force_env_type), POINTER :: force_env
54 : TYPE(simpar_type), POINTER :: simpar
55 :
56 41595 : CALL timeset(routineN, handle)
57 :
58 : ! Get force environment
59 41595 : CALL get_md_env(md_env, force_env=force_env, simpar=simpar)
60 :
61 : ! RESPA implemented only for NVE
62 41595 : IF (simpar%do_respa .AND. nve_ensemble .NE. simpar%ensemble) THEN
63 0 : CPABORT("RESPA integrator not implemented for this ensemble")
64 : END IF
65 :
66 : ! Choice of the ensemble
67 41595 : SELECT CASE (simpar%ensemble)
68 : CASE DEFAULT
69 0 : CPABORT("Integrator not implemented")
70 : CASE (nve_ensemble)
71 31199 : IF (simpar%do_respa) THEN
72 14 : CALL nve_respa(md_env)
73 : ELSE
74 31185 : CALL nve(md_env, globenv)
75 : END IF
76 : CASE (nvt_ensemble)
77 7360 : CALL nvt(md_env, globenv)
78 : CASE (nvt_adiabatic_ensemble)
79 0 : CALL nvt_adiabatic(md_env, globenv)
80 : CASE (isokin_ensemble)
81 6 : CALL isokin(md_env)
82 : CASE (npt_i_ensemble, npt_ia_ensemble, npe_i_ensemble)
83 1564 : CALL npt_i(md_env, globenv)
84 : CASE (npt_f_ensemble)
85 656 : CALL npt_f(md_env, globenv)
86 : CASE (nph_uniaxial_ensemble)
87 40 : CALL nph_uniaxial(md_env)
88 : CASE (nph_uniaxial_damped_ensemble)
89 20 : CALL nph_uniaxial_damped(md_env)
90 : CASE (reftraj_ensemble)
91 288 : CALL reftraj(md_env)
92 : CASE (langevin_ensemble)
93 222 : CALL langevin(md_env)
94 : CASE (npe_f_ensemble)
95 41595 : CALL npt_f(md_env, globenv)
96 : END SELECT
97 :
98 41595 : CALL timestop(handle)
99 :
100 41595 : END SUBROUTINE velocity_verlet
101 :
102 : END MODULE velocity_verlet_control
|