(git:374b731)
Loading...
Searching...
No Matches
tmc_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 module handles definition of the tree nodes for the global and
10!> the subtrees binary tree
11!> parent element
12!> / \
13!> accepted (acc) / \ not accepted (nacc)
14!> / \
15!> child child
16!> / \ / \
17!>
18!> tree creation assuming acceptance (acc) AND rejectance (nacc)
19!> of configuration
20!> if configuration is accepted: new configuration (child on acc) on basis
21!> of last configuration (one level up)
22!> if configuration is rejected: child on nacc on basis of last accepted
23!> element (last element which is on acc brach of its parent element)
24!> The global tree handles all configurations of different subtrees.
25!> The structure element "conf" is an array related to the temperature
26!> (sorted) and points to the subtree elements.
27!> \par History
28!> 11.2012 created [Mandes Schoenherr]
29!> \author Mandes
30! **************************************************************************************************
31
33 USE cell_types, ONLY: cell_type
34 USE kinds, ONLY: default_path_length,&
36 dp
41 USE tmc_stati, ONLY: task_type_mc
42 USE tmc_tree_types, ONLY: clean_list,&
46#include "../base/base_uses.f90"
47
48 IMPLICIT NONE
49
50 PRIVATE
51
52 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'tmc_types'
53
58 PUBLIC :: tmc_atom_type
60
61 ! global environment
63 TYPE(tmc_comp_set_type), POINTER :: tmc_comp_set => null()
64 TYPE(tmc_param_type), POINTER :: params => null()
65 TYPE(rng_stream_type), ALLOCATABLE :: rng_stream
66 TYPE(master_env_type), POINTER :: m_env => null()
67 TYPE(worker_env_type), POINTER :: w_env => null()
68 END TYPE tmc_env_type
69
70 ! structure for remembering the main values used for reordering MPI communicators
71 ! \param group_nr: the first group_ener_nr groups are for energy calculation,
72 ! then the Configurational Change groups,
73 ! the global master has group 0 and unused cores have negative group numbering
75 INTEGER :: group_ener_size = 0
76 INTEGER :: group_ener_nr = 0
77 INTEGER :: group_cc_size = 0
78 INTEGER :: group_cc_nr = 0
79 INTEGER :: group_nr = 0
80 INTEGER :: ana_on_the_fly = -1
81 ! the communicators (para_env)
82 TYPE(mp_para_env_type), POINTER :: para_env_m_w => null()
83 TYPE(mp_para_env_type), POINTER :: para_env_sub_group => null()
84 TYPE(mp_para_env_type), POINTER :: para_env_m_first_w => null()
85 TYPE(mp_para_env_type), POINTER :: para_env_m_ana => null()
86 TYPE(mp_para_env_type), POINTER :: para_env_m_only => null()
87 END TYPE tmc_comp_set_type
88
89 ! struct for TMC global variables
91 INTEGER :: task_type = task_type_mc
92 INTEGER :: dim_per_elem = 3
93 INTEGER :: nr_temp = -1
94 REAL(kind=dp), DIMENSION(:), POINTER :: temp => null()
95 TYPE(cell_type), POINTER :: cell => null()
96 REAL(kind=dp), DIMENSION(:), POINTER :: sub_box_size => null()
97 TYPE(tmc_atom_type), DIMENSION(:), POINTER :: atoms => null()
98
99 INTEGER :: nr_elem_mv = -1
100 TYPE(tmc_move_type), POINTER :: move_types => null()
101 TYPE(tmc_move_type), POINTER :: nmc_move_types => null()
102 REAL(kind=dp) :: pressure = 0.0_dp
103 LOGICAL :: v_isotropic = .false.
104 LOGICAL :: mv_cen_of_mass = .false.
105 LOGICAL :: esimate_acc_prob = .false.
106 LOGICAL :: speculative_canceling = .false.
107 LOGICAL :: use_scf_energy_info = .false.
108 LOGICAL :: use_reduced_tree = .false.
109 CHARACTER(LEN=default_path_length) :: energy_inp_file = ""
110 CHARACTER(LEN=default_path_length) :: nmc_inp_file = ""
111 LOGICAL :: draw_tree = .false.
112 CHARACTER(LEN=default_path_length) :: dot_file_name = ""
113 CHARACTER(LEN=default_path_length) :: all_conf_file_name = ""
114 LOGICAL :: print_only_diff_conf = .false.
115 LOGICAL :: print_trajectory = .false.
116 LOGICAL :: print_dipole = .false.
117 LOGICAL :: print_forces = .false.
118 LOGICAL :: print_cell = .false.
119 LOGICAL :: print_energies = .false.
120 TYPE(prior_estimate_acceptance_type), POINTER :: prior_nmc_acc => null()
121 LOGICAL :: print_test_output = .false.
122 END TYPE tmc_param_type
123
125 CHARACTER(LEN=default_string_length) :: name = ""
126 REAL(kind=dp) :: mass = 0.0_dp
127 END TYPE
128
129 ! to estimate the prior acceptance
130 TYPE prior_estimate_acceptance_type
131 INTEGER :: counter = 0
132 REAL(kind=dp) :: aver = 0.0_dp, aver_2 = 0.0_dp
133 END TYPE prior_estimate_acceptance_type
134
135 ! environments for the master
136 TYPE master_env_type
137 INTEGER :: num_mc_elem = 0! the specified number of Markov Chain elements, to be reached
138 CHARACTER(LEN=default_path_length) :: restart_in_file_name = ""
139 CHARACTER(LEN=default_path_length) :: restart_out_file_name = ""
140 INTEGER :: restart_out_step = 0
141 INTEGER :: io_unit = -1
142 INTEGER :: info_out_step_size = 0
143 REAL(kind=dp) :: walltime = 0.0_dp
144 INTEGER :: rnd_init = 0
145 REAL(kind=dp) :: temp_decrease = 0.0_dp ! for simulated annealing
146 TYPE(elem_list_type), POINTER :: cancelation_list => null()
147 INTEGER :: count_cancel_ener = 0
148 INTEGER :: count_cancel_nmc = 0
149 ! masters tree stuff
150 TYPE(global_tree_type), POINTER :: gt_head => null(), gt_act => null()
151 INTEGER, DIMENSION(:), POINTER :: tree_node_count => null()
152 INTEGER, DIMENSION(:), POINTER :: result_count => null()
153 TYPE(elem_array_type), DIMENSION(:), &
154 POINTER :: result_list => null(), &
155 st_heads => null(), &
156 st_clean_ends => null()
157 TYPE(global_tree_type), POINTER :: gt_clean_end => null()
158 INTEGER, DIMENSION(4) :: estim_corr_wrong = 0
159 TYPE(elem_list_type), POINTER :: analysis_list => null()
160 END TYPE master_env_type
161
162 ! environment for the worker
163 TYPE worker_env_type
164 INTEGER :: env_id_ener = -1, env_id_approx = -1
165 INTEGER :: io_unit = -1
166 REAL(kind=dp) :: act_temp = 0.0_dp
167 END TYPE worker_env_type
168
169CONTAINS
170
171! **************************************************************************************************
172!> \brief creates a new structure environment for TMC
173!> \param tmc_env structure with parameters for TMC
174!> \author Mandes 11.2012
175! **************************************************************************************************
176 SUBROUTINE tmc_env_create(tmc_env)
177 TYPE(tmc_env_type), POINTER :: tmc_env
178
179 CHARACTER(LEN=*), PARAMETER :: routinen = 'tmc_env_create'
180
181 INTEGER :: handle
182
183 CALL timeset(routinen, handle)
184
185 cpassert(.NOT. ASSOCIATED(tmc_env))
186
187 ALLOCATE (tmc_env)
188
189 ALLOCATE (tmc_env%tmc_comp_set)
190
191 ! initialize the parameter section
192 ALLOCATE (tmc_env%params)
193
194 ALLOCATE (tmc_env%params%sub_box_size(tmc_env%params%dim_per_elem))
195 tmc_env%params%sub_box_size(:) = -1.0_dp
196
197 CALL timestop(handle)
198
199 END SUBROUTINE tmc_env_create
200
201! **************************************************************************************************
202!> \brief releases the structure environment for TMC
203!> \param tmc_env structure with parameters for TMC
204!> \author Mandes 11.2012
205! **************************************************************************************************
206 SUBROUTINE tmc_env_release(tmc_env)
207 TYPE(tmc_env_type), POINTER :: tmc_env
208
209 CHARACTER(LEN=*), PARAMETER :: routinen = 'tmc_env_release'
210
211 INTEGER :: handle
212
213 CALL timeset(routinen, handle)
214
215 cpassert(ASSOCIATED(tmc_env))
216 cpassert(ASSOCIATED(tmc_env%params))
217
218 DEALLOCATE (tmc_env%params%sub_box_size)
219 IF (ASSOCIATED(tmc_env%params%Temp)) &
220 DEALLOCATE (tmc_env%params%Temp)
221 IF (ASSOCIATED(tmc_env%params%cell)) &
222 DEALLOCATE (tmc_env%params%cell)
223 IF (ASSOCIATED(tmc_env%params%atoms)) &
224 CALL deallocate_tmc_atom_type(tmc_env%params%atoms)
225 DEALLOCATE (tmc_env%params)
226
227 CALL mp_para_env_release(tmc_env%tmc_comp_set%para_env_sub_group)
228 CALL mp_para_env_release(tmc_env%tmc_comp_set%para_env_m_w)
229 IF (ASSOCIATED(tmc_env%tmc_comp_set%para_env_m_first_w)) &
230 CALL mp_para_env_release(tmc_env%tmc_comp_set%para_env_m_first_w)
231 IF (ASSOCIATED(tmc_env%tmc_comp_set%para_env_m_ana)) &
232 CALL mp_para_env_release(tmc_env%tmc_comp_set%para_env_m_ana)
233 IF (ASSOCIATED(tmc_env%tmc_comp_set%para_env_m_only)) &
234 CALL mp_para_env_release(tmc_env%tmc_comp_set%para_env_m_only)
235
236 DEALLOCATE (tmc_env%tmc_comp_set)
237
238 DEALLOCATE (tmc_env)
239
240 CALL timestop(handle)
241
242 END SUBROUTINE tmc_env_release
243
244! **************************************************************************************************
245!> \brief creates a new structure environment for TMC master
246!> \param tmc_env structure with parameters for TMC
247!> \author Mandes 11.2012
248! **************************************************************************************************
249 SUBROUTINE tmc_master_env_create(tmc_env)
250 TYPE(tmc_env_type), POINTER :: tmc_env
251
252 CHARACTER(LEN=*), PARAMETER :: routinen = 'tmc_master_env_create'
253
254 INTEGER :: handle, i
255
256 CALL timeset(routinen, handle)
257
258 cpassert(ASSOCIATED(tmc_env))
259 cpassert(ASSOCIATED(tmc_env%params))
260 cpassert(tmc_env%params%nr_temp .GT. 0)
261
262 cpassert(.NOT. ASSOCIATED(tmc_env%m_env))
263
264 ALLOCATE (tmc_env%m_env)
265 NULLIFY (tmc_env%m_env%gt_head, tmc_env%m_env%gt_act, tmc_env%m_env%tree_node_count, &
266 tmc_env%m_env%result_count, tmc_env%m_env%result_list, &
267 tmc_env%m_env%st_heads, tmc_env%m_env%st_clean_ends, &
268 tmc_env%m_env%gt_clean_end, tmc_env%m_env%cancelation_list, tmc_env%m_env%analysis_list)
269
270 tmc_env%m_env%restart_in_file_name = ""
271 tmc_env%m_env%restart_out_file_name = ""
272 ALLOCATE (tmc_env%m_env%tree_node_count(0:tmc_env%params%nr_temp))
273 tmc_env%m_env%tree_node_count(:) = 0
274 ALLOCATE (tmc_env%m_env%result_count(0:tmc_env%params%nr_temp))
275 tmc_env%m_env%result_count(:) = 0
276 ALLOCATE (tmc_env%m_env%st_heads(tmc_env%params%nr_temp))
277 ALLOCATE (tmc_env%m_env%st_clean_ends(tmc_env%params%nr_temp))
278
279 IF (tmc_env%params%USE_REDUCED_TREE) ALLOCATE (tmc_env%m_env%result_list(tmc_env%params%nr_temp))
280
281 DO i = 1, tmc_env%params%nr_temp
282 tmc_env%m_env%st_heads(i)%elem => null()
283 tmc_env%m_env%st_clean_ends(i)%elem => null()
284 IF (tmc_env%params%USE_REDUCED_TREE) &
285 tmc_env%m_env%result_list(i)%elem => null()
286 END DO
287 tmc_env%m_env%gt_head => null()
288 tmc_env%m_env%gt_clean_end => null()
289 tmc_env%m_env%temp_decrease = 1.0_dp
290 tmc_env%m_env%count_cancel_ener = 0
291 tmc_env%m_env%count_cancel_NMC = 0
292 tmc_env%m_env%estim_corr_wrong(:) = 0
293
294 ALLOCATE (tmc_env%params%prior_NMC_acc)
295 tmc_env%params%prior_NMC_acc%counter = 0
296 tmc_env%params%prior_NMC_acc%aver = 0.0_dp
297 tmc_env%params%prior_NMC_acc%aver_2 = 0.0_dp
298
299 CALL timestop(handle)
300
301 END SUBROUTINE tmc_master_env_create
302
303! **************************************************************************************************
304!> \brief releases the structure environment for TMC master
305!> \param tmc_env structure with parameters for TMC
306!> \author Mandes 11.2012
307! **************************************************************************************************
308 SUBROUTINE tmc_master_env_release(tmc_env)
309 TYPE(tmc_env_type), POINTER :: tmc_env
310
311 CHARACTER(LEN=*), PARAMETER :: routinen = 'tmc_master_env_release'
312
313 INTEGER :: handle
314
315 CALL timeset(routinen, handle)
316
317 cpassert(ASSOCIATED(tmc_env))
318 cpassert(ASSOCIATED(tmc_env%m_env))
319
320 CALL clean_list(tmc_env%m_env%analysis_list)
321 CALL clean_list(tmc_env%m_env%cancelation_list)
322
323 DEALLOCATE (tmc_env%m_env%tree_node_count)
324 DEALLOCATE (tmc_env%m_env%result_count)
325 DEALLOCATE (tmc_env%m_env%st_heads)
326 DEALLOCATE (tmc_env%m_env%st_clean_ends)
327 IF (tmc_env%params%USE_REDUCED_TREE) DEALLOCATE (tmc_env%m_env%result_list)
328 DEALLOCATE (tmc_env%params%prior_NMC_acc)
329
330 DEALLOCATE (tmc_env%m_env)
331
332 CALL timestop(handle)
333
334 END SUBROUTINE tmc_master_env_release
335
336! **************************************************************************************************
337!> \brief creates a new structure environment for TMC master
338!> \param tmc_env structure with parameters for TMC
339!> \author Mandes 11.2012
340! **************************************************************************************************
341 SUBROUTINE tmc_worker_env_create(tmc_env)
342 TYPE(tmc_env_type), POINTER :: tmc_env
343
344 CHARACTER(LEN=*), PARAMETER :: routinen = 'tmc_worker_env_create'
345
346 INTEGER :: handle
347
348 CALL timeset(routinen, handle)
349
350 cpassert(ASSOCIATED(tmc_env))
351 cpassert(.NOT. ASSOCIATED(tmc_env%w_env))
352
353 ALLOCATE (tmc_env%w_env)
354
355 tmc_env%w_env%env_id_ener = -1
356 tmc_env%w_env%env_id_approx = -1
357 tmc_env%w_env%io_unit = -1
358 tmc_env%w_env%act_temp = -1.0_dp
359
360 CALL timestop(handle)
361
362 END SUBROUTINE tmc_worker_env_create
363
364! **************************************************************************************************
365!> \brief releases the structure environment for TMC master
366!> \param tmc_env structure with parameters for TMC
367!> \author Mandes 11.2012
368! **************************************************************************************************
369 SUBROUTINE tmc_worker_env_release(tmc_env)
370 TYPE(tmc_env_type), POINTER :: tmc_env
371
372 CHARACTER(LEN=*), PARAMETER :: routinen = 'tmc_worker_env_release'
373
374 INTEGER :: handle
375
376 CALL timeset(routinen, handle)
377
378 cpassert(ASSOCIATED(tmc_env))
379 cpassert(ASSOCIATED(tmc_env%w_env))
380
381 DEALLOCATE (tmc_env%w_env)
382
383 CALL timestop(handle)
384
385 END SUBROUTINE tmc_worker_env_release
386
387! **************************************************************************************************
388!> \brief creates a structure for storing the atom informations
389!> \param atoms pointer to a list of tmc_atoms_type
390!> \param nr_atoms the amount of atoms
391!> \author Mandes 01.2013
392! **************************************************************************************************
393 SUBROUTINE allocate_tmc_atom_type(atoms, nr_atoms)
394 TYPE(tmc_atom_type), DIMENSION(:), POINTER :: atoms
395 INTEGER, INTENT(IN) :: nr_atoms
396
397 cpassert(.NOT. ASSOCIATED(atoms))
398 cpassert(nr_atoms .GT. 0)
399
400 ALLOCATE (atoms(nr_atoms))
401
402 cpassert(ASSOCIATED(atoms))
403
404 END SUBROUTINE allocate_tmc_atom_type
405
406! **************************************************************************************************
407!> \brief releases the structure for storing the atom informations
408!> \param atoms pointer to a list of tmc_atoms_type
409!> \author Mandes 01.2013
410! **************************************************************************************************
411 SUBROUTINE deallocate_tmc_atom_type(atoms)
412 TYPE(tmc_atom_type), DIMENSION(:), POINTER :: atoms
413
414 cpassert(ASSOCIATED(atoms))
415
416 DEALLOCATE (atoms)
417
418 cpassert(.NOT. ASSOCIATED(atoms))
419 END SUBROUTINE deallocate_tmc_atom_type
420
421END MODULE tmc_types
Handles all functions related to the CELL.
Definition cell_types.F:15
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
integer, parameter, public default_path_length
Definition kinds.F:58
Interface to the message passing library MPI.
subroutine, public mp_para_env_release(para_env)
releases the para object (to be called when you don't want anymore the shared copy of this object)
Parallel (pseudo)random number generator (RNG) for multiple streams and substreams of random numbers.
tree nodes creation, searching, deallocation, references etc.
tree nodes creation, searching, deallocation, references etc.
Definition tmc_stati.F:15
integer, parameter, public task_type_mc
Definition tmc_stati.F:44
module handles definition of the tree nodes for the global and the subtrees binary tree parent elemen...
subroutine, public clean_list(list)
clean a certain element element list
module handles definition of the tree nodes for the global and the subtrees binary tree parent elemen...
Definition tmc_types.F:32
subroutine, public allocate_tmc_atom_type(atoms, nr_atoms)
creates a structure for storing the atom informations
Definition tmc_types.F:394
subroutine, public tmc_env_release(tmc_env)
releases the structure environment for TMC
Definition tmc_types.F:207
subroutine, public tmc_worker_env_create(tmc_env)
creates a new structure environment for TMC master
Definition tmc_types.F:342
subroutine, public tmc_master_env_create(tmc_env)
creates a new structure environment for TMC master
Definition tmc_types.F:250
subroutine, public tmc_master_env_release(tmc_env)
releases the structure environment for TMC master
Definition tmc_types.F:309
subroutine, public tmc_worker_env_release(tmc_env)
releases the structure environment for TMC master
Definition tmc_types.F:370
subroutine, public tmc_env_create(tmc_env)
creates a new structure environment for TMC
Definition tmc_types.F:177
Type defining parameters related to the simulation cell.
Definition cell_types.F:55
stores all the informations relevant to an mpi environment