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 defines types for metadynamics calculation
10 : !> \par History
11 : !> 01.2007 created [tlaino] Teodoro Laino
12 : ! **************************************************************************************************
13 : MODULE free_energy_types
14 :
15 : USE input_constants, ONLY: do_fe_ac,&
16 : do_fe_ui
17 : USE input_section_types, ONLY: section_vals_get,&
18 : section_vals_get_subs_vals,&
19 : section_vals_type,&
20 : section_vals_val_get
21 : USE kinds, ONLY: default_string_length,&
22 : dp
23 : #include "./base/base_uses.f90"
24 :
25 : IMPLICIT NONE
26 :
27 : PRIVATE
28 :
29 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'free_energy_types'
30 :
31 : PUBLIC :: free_energy_type, &
32 : fe_env_release, &
33 : fe_env_create, &
34 : ui_var_type
35 :
36 : ! **************************************************************************************************
37 : TYPE ui_var_type
38 : REAL(KIND=dp), DIMENSION(:), POINTER :: ss => NULL()
39 : INTEGER :: icolvar = -1
40 : END TYPE ui_var_type
41 :
42 : ! **************************************************************************************************
43 : TYPE ui_conv_type
44 : ! Specifying convergence parameters
45 : INTEGER :: cg_width = -1, max_cg_width = -1
46 : INTEGER :: cg_points = -1
47 : REAL(KIND=dp) :: eps_conv = 0.0_dp
48 : REAL(KIND=dp) :: k_conf_lm = 0.0_dp
49 : REAL(KIND=dp) :: sw_conf_lm = 0.0_dp
50 : REAL(KIND=dp) :: vn_conf_lm = 0.0_dp
51 : LOGICAL :: test_k = .FALSE., &
52 : test_sw = .FALSE., &
53 : test_vn = .FALSE.
54 : END TYPE ui_conv_type
55 :
56 : ! **************************************************************************************************
57 : TYPE statistical_type
58 : ! Collecting coarse grained data
59 : REAL(KIND=dp), DIMENSION(:), POINTER :: avg => NULL()
60 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: var => NULL()
61 : END TYPE statistical_type
62 :
63 : ! **************************************************************************************************
64 : TYPE free_energy_type
65 : INTEGER :: ncolvar = -1
66 : INTEGER :: TYPE = -1
67 : INTEGER :: nr_points = -1, &
68 : nr_rejected = -1
69 : TYPE(ui_conv_type), POINTER :: conv_par => NULL()
70 : TYPE(ui_var_type), POINTER, DIMENSION(:) :: uivar => NULL()
71 : TYPE(statistical_type), DIMENSION(:), POINTER :: cg_data => NULL()
72 : ! Old data
73 : REAL(KIND=dp) :: eps_conv = 0.0_dp
74 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: covmx => NULL()
75 :
76 : CHARACTER(len=default_string_length) :: plumed_input_file = ""
77 : END TYPE free_energy_type
78 :
79 : CONTAINS
80 :
81 : ! **************************************************************************************************
82 : !> \brief creates the fe_env
83 : !> \param fe_env ...
84 : !> \param fe_section ...
85 : !> \author Teodoro Laino 01.2007
86 : ! **************************************************************************************************
87 3572 : SUBROUTINE fe_env_create(fe_env, fe_section)
88 : TYPE(free_energy_type), POINTER :: fe_env
89 : TYPE(section_vals_type), POINTER :: fe_section
90 :
91 : INTEGER :: i, id_method
92 : LOGICAL :: explicit
93 : TYPE(section_vals_type), POINTER :: ui_section, ui_var_section
94 :
95 1786 : CPASSERT(.NOT. ASSOCIATED(fe_env))
96 :
97 1786 : CALL section_vals_get(fe_section, explicit=explicit)
98 1786 : IF (explicit) THEN
99 172 : CALL section_vals_val_get(fe_section, "METHOD", i_val=id_method)
100 2 : SELECT CASE (id_method)
101 : CASE (do_fe_ui)
102 2 : ALLOCATE (fe_env)
103 : NULLIFY (fe_env%covmx, fe_env%uivar, fe_env%conv_par, fe_env%cg_data)
104 2 : fe_env%type = id_method
105 2 : fe_env%nr_points = 0
106 2 : fe_env%nr_rejected = 0
107 : NULLIFY (fe_env%cg_data)
108 2 : ui_section => section_vals_get_subs_vals(fe_section, "UMBRELLA_INTEGRATION")
109 2 : ui_var_section => section_vals_get_subs_vals(ui_section, "UVAR")
110 2 : CALL section_vals_get(ui_var_section, n_repetition=fe_env%ncolvar)
111 : ! Convergence controlling parameters
112 2 : ALLOCATE (fe_env%conv_par)
113 2 : fe_env%conv_par%test_k = .FALSE.
114 2 : fe_env%conv_par%test_sw = .FALSE.
115 2 : fe_env%conv_par%test_vn = .FALSE.
116 : CALL section_vals_val_get(ui_section, "CONVERGENCE_CONTROL%COARSE_GRAINED_WIDTH", &
117 2 : i_val=fe_env%conv_par%cg_width)
118 : CALL section_vals_val_get(ui_section, "CONVERGENCE_CONTROL%MAX_COARSE_GRAINED_WIDTH", &
119 2 : i_val=fe_env%conv_par%max_cg_width)
120 : CALL section_vals_val_get(ui_section, "CONVERGENCE_CONTROL%COARSE_GRAINED_POINTS", &
121 2 : i_val=fe_env%conv_par%cg_points)
122 : CALL section_vals_val_get(ui_section, "CONVERGENCE_CONTROL%EPS_CONV", &
123 2 : r_val=fe_env%conv_par%eps_conv)
124 : CALL section_vals_val_get(ui_section, "CONVERGENCE_CONTROL%K_CONFIDENCE_LIMIT", &
125 2 : r_val=fe_env%conv_par%k_conf_lm)
126 : CALL section_vals_val_get(ui_section, "CONVERGENCE_CONTROL%SW_CONFIDENCE_LIMIT", &
127 2 : r_val=fe_env%conv_par%sw_conf_lm)
128 : CALL section_vals_val_get(ui_section, "CONVERGENCE_CONTROL%VN_CONFIDENCE_LIMIT", &
129 2 : r_val=fe_env%conv_par%vn_conf_lm)
130 : ! Umbrella Integration variables
131 8 : ALLOCATE (fe_env%uivar(fe_env%ncolvar))
132 4 : DO i = 1, fe_env%ncolvar
133 : ! Read Umbrella Integration Variable definition
134 : CALL section_vals_val_get(ui_var_section, "COLVAR", &
135 2 : i_val=fe_env%uivar(i)%icolvar, i_rep_section=i)
136 4 : NULLIFY (fe_env%uivar(i)%ss)
137 : END DO
138 : CASE (do_fe_ac)
139 12 : ALLOCATE (fe_env)
140 : NULLIFY (fe_env%covmx, fe_env%uivar, fe_env%conv_par, fe_env%cg_data)
141 12 : ALLOCATE (fe_env%covmx(3, 0))
142 12 : fe_env%type = id_method
143 184 : CALL section_vals_val_get(fe_section, "ALCHEMICAL_CHANGE%EPS_CONV", r_val=fe_env%eps_conv)
144 : CASE DEFAULT
145 : ! Do Nothing
146 : END SELECT
147 : END IF
148 1786 : END SUBROUTINE fe_env_create
149 :
150 : ! **************************************************************************************************
151 : !> \brief releases the fe_env
152 : !> \param fe_env ...
153 : !> \author Laino Teodoro 01.2007
154 : ! **************************************************************************************************
155 1786 : SUBROUTINE fe_env_release(fe_env)
156 : TYPE(free_energy_type), POINTER :: fe_env
157 :
158 : INTEGER :: i
159 :
160 1786 : IF (ASSOCIATED(fe_env)) THEN
161 14 : IF (ASSOCIATED(fe_env%covmx)) THEN
162 12 : DEALLOCATE (fe_env%covmx)
163 : END IF
164 14 : IF (ASSOCIATED(fe_env%cg_data)) THEN
165 0 : DO i = 1, SIZE(fe_env%cg_data)
166 0 : IF (ASSOCIATED(fe_env%cg_data(i)%avg)) THEN
167 0 : DEALLOCATE (fe_env%cg_data(i)%avg)
168 : END IF
169 0 : IF (ASSOCIATED(fe_env%cg_data(i)%var)) THEN
170 0 : DEALLOCATE (fe_env%cg_data(i)%var)
171 : END IF
172 : END DO
173 0 : DEALLOCATE (fe_env%cg_data)
174 : END IF
175 14 : IF (ASSOCIATED(fe_env%conv_par)) THEN
176 2 : DEALLOCATE (fe_env%conv_par)
177 : END IF
178 14 : IF (ASSOCIATED(fe_env%uivar)) THEN
179 4 : DO i = 1, SIZE(fe_env%uivar)
180 4 : IF (ASSOCIATED(fe_env%uivar(i)%ss)) THEN
181 2 : DEALLOCATE (fe_env%uivar(i)%ss)
182 : END IF
183 : END DO
184 2 : DEALLOCATE (fe_env%uivar)
185 : END IF
186 14 : DEALLOCATE (fe_env)
187 : END IF
188 1786 : END SUBROUTINE fe_env_release
189 :
190 0 : END MODULE free_energy_types
|