(git:374b731)
Loading...
Searching...
No Matches
free_energy_types.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 defines types for metadynamics calculation
10!> \par History
11!> 01.2007 created [tlaino] Teodoro Laino
12! **************************************************************************************************
14
15 USE input_constants, ONLY: do_fe_ac,&
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, &
35
36! **************************************************************************************************
38 REAL(kind=dp), DIMENSION(:), POINTER :: ss
39 INTEGER :: icolvar
40 END TYPE ui_var_type
41
42! **************************************************************************************************
43 TYPE ui_conv_type
44 ! Specifying convergence parameters
45 INTEGER :: cg_width, max_cg_width
46 INTEGER :: cg_points
47 REAL(kind=dp) :: eps_conv
48 REAL(kind=dp) :: k_conf_lm
49 REAL(kind=dp) :: sw_conf_lm
50 REAL(kind=dp) :: vn_conf_lm
51 LOGICAL :: test_k, &
52 test_sw, &
53 test_vn
54 END TYPE ui_conv_type
55
56! **************************************************************************************************
57 TYPE statistical_type
58 ! Collecting coarse grained data
59 REAL(kind=dp), DIMENSION(:), POINTER :: avg
60 REAL(kind=dp), DIMENSION(:, :), POINTER :: var
61 END TYPE statistical_type
62
63! **************************************************************************************************
65 INTEGER :: ncolvar
66 INTEGER :: type
67 INTEGER :: nr_points, &
68 nr_rejected
69 TYPE(ui_conv_type), POINTER :: conv_par
70 TYPE(ui_var_type), POINTER, DIMENSION(:) :: uivar
71 TYPE(statistical_type), DIMENSION(:), POINTER :: cg_data
72 ! Old data
73 REAL(kind=dp) :: eps_conv
74 REAL(kind=dp), DIMENSION(:, :), POINTER :: covmx
75
76 CHARACTER(len=default_string_length) :: plumed_input_file
77 END TYPE free_energy_type
78
79CONTAINS
80
81! **************************************************************************************************
82!> \brief creates the fe_env
83!> \param fe_env ...
84!> \param fe_section ...
85!> \author Teodoro Laino 01.2007
86! **************************************************************************************************
87 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 cpassert(.NOT. ASSOCIATED(fe_env))
96
97 CALL section_vals_get(fe_section, explicit=explicit)
98 IF (explicit) THEN
99 CALL section_vals_val_get(fe_section, "METHOD", i_val=id_method)
100 SELECT CASE (id_method)
101 CASE (do_fe_ui)
102 ALLOCATE (fe_env)
103 NULLIFY (fe_env%covmx, fe_env%uivar, fe_env%conv_par, fe_env%cg_data)
104 fe_env%type = id_method
105 fe_env%nr_points = 0
106 fe_env%nr_rejected = 0
107 NULLIFY (fe_env%cg_data)
108 ui_section => section_vals_get_subs_vals(fe_section, "UMBRELLA_INTEGRATION")
109 ui_var_section => section_vals_get_subs_vals(ui_section, "UVAR")
110 CALL section_vals_get(ui_var_section, n_repetition=fe_env%ncolvar)
111 ! Convergence controlling parameters
112 ALLOCATE (fe_env%conv_par)
113 fe_env%conv_par%test_k = .false.
114 fe_env%conv_par%test_sw = .false.
115 fe_env%conv_par%test_vn = .false.
116 CALL section_vals_val_get(ui_section, "CONVERGENCE_CONTROL%COARSE_GRAINED_WIDTH", &
117 i_val=fe_env%conv_par%cg_width)
118 CALL section_vals_val_get(ui_section, "CONVERGENCE_CONTROL%MAX_COARSE_GRAINED_WIDTH", &
119 i_val=fe_env%conv_par%max_cg_width)
120 CALL section_vals_val_get(ui_section, "CONVERGENCE_CONTROL%COARSE_GRAINED_POINTS", &
121 i_val=fe_env%conv_par%cg_points)
122 CALL section_vals_val_get(ui_section, "CONVERGENCE_CONTROL%EPS_CONV", &
123 r_val=fe_env%conv_par%eps_conv)
124 CALL section_vals_val_get(ui_section, "CONVERGENCE_CONTROL%K_CONFIDENCE_LIMIT", &
125 r_val=fe_env%conv_par%k_conf_lm)
126 CALL section_vals_val_get(ui_section, "CONVERGENCE_CONTROL%SW_CONFIDENCE_LIMIT", &
127 r_val=fe_env%conv_par%sw_conf_lm)
128 CALL section_vals_val_get(ui_section, "CONVERGENCE_CONTROL%VN_CONFIDENCE_LIMIT", &
129 r_val=fe_env%conv_par%vn_conf_lm)
130 ! Umbrella Integration variables
131 ALLOCATE (fe_env%uivar(fe_env%ncolvar))
132 DO i = 1, fe_env%ncolvar
133 ! Read Umbrella Integration Variable definition
134 CALL section_vals_val_get(ui_var_section, "COLVAR", &
135 i_val=fe_env%uivar(i)%icolvar, i_rep_section=i)
136 NULLIFY (fe_env%uivar(i)%ss)
137 END DO
138 CASE (do_fe_ac)
139 ALLOCATE (fe_env)
140 NULLIFY (fe_env%covmx, fe_env%uivar, fe_env%conv_par, fe_env%cg_data)
141 ALLOCATE (fe_env%covmx(3, 0))
142 fe_env%type = id_method
143 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 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 SUBROUTINE fe_env_release(fe_env)
156 TYPE(free_energy_type), POINTER :: fe_env
157
158 INTEGER :: i
159
160 IF (ASSOCIATED(fe_env)) THEN
161 IF (ASSOCIATED(fe_env%covmx)) THEN
162 DEALLOCATE (fe_env%covmx)
163 END IF
164 IF (ASSOCIATED(fe_env%cg_data)) THEN
165 DO i = 1, SIZE(fe_env%cg_data)
166 IF (ASSOCIATED(fe_env%cg_data(i)%avg)) THEN
167 DEALLOCATE (fe_env%cg_data(i)%avg)
168 END IF
169 IF (ASSOCIATED(fe_env%cg_data(i)%var)) THEN
170 DEALLOCATE (fe_env%cg_data(i)%var)
171 END IF
172 END DO
173 DEALLOCATE (fe_env%cg_data)
174 END IF
175 IF (ASSOCIATED(fe_env%conv_par)) THEN
176 DEALLOCATE (fe_env%conv_par)
177 END IF
178 IF (ASSOCIATED(fe_env%uivar)) THEN
179 DO i = 1, SIZE(fe_env%uivar)
180 IF (ASSOCIATED(fe_env%uivar(i)%ss)) THEN
181 DEALLOCATE (fe_env%uivar(i)%ss)
182 END IF
183 END DO
184 DEALLOCATE (fe_env%uivar)
185 END IF
186 DEALLOCATE (fe_env)
187 END IF
188 END SUBROUTINE fe_env_release
189
190END MODULE free_energy_types
defines types for metadynamics calculation
subroutine, public fe_env_create(fe_env, fe_section)
creates the fe_env
subroutine, public fe_env_release(fe_env)
releases the fe_env
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public do_fe_ac
integer, parameter, public do_fe_ui
objects that represent the structure of input sections and the data contained in an input 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
integer, parameter, public default_string_length
Definition kinds.F:57