(git:b31c023)
Loading...
Searching...
No Matches
fist_environment_types.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2025 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \par History
10!> gt sept-23-02 added atomic_kind_set to replica_environment_type
11!> to allow use of kind_based neighbor list
12!> CJM rewrite
13!> \author CJM SEPT-01-02
14! **************************************************************************************************
20 USE cell_types, ONLY: cell_release,&
61#include "./base/base_uses.f90"
62
63 IMPLICIT NONE
64 PRIVATE
65
66! **************************************************************************************************
67!> \par History
68!> 11/03
69!> \author CJM
70! **************************************************************************************************
72 PRIVATE
73 LOGICAL :: qmmm = .false.
74 LOGICAL :: shell_model = .false., shell_model_ad = .false.
75 TYPE(qmmm_env_mm_type), POINTER :: qmmm_env => null()
76 TYPE(cell_type), POINTER :: cell_ref => null()
77 TYPE(ewald_environment_type), POINTER :: ewald_env => null()
78 TYPE(ewald_pw_type), POINTER :: ewald_pw => null()
79 TYPE(fist_energy_type), POINTER :: thermo => null()
80 TYPE(mp_para_env_type), POINTER :: para_env => null()
81 TYPE(cp_subsys_type), POINTER :: subsys => null()
82 TYPE(fist_nonbond_env_type), POINTER :: fist_nonbond_env => null()
83 TYPE(section_vals_type), POINTER :: input => null()
84 TYPE(exclusion_type), DIMENSION(:), POINTER :: exclusions => null()
85 TYPE(fist_efield_type), POINTER :: efield => null()
87
88! *** Public data types ***
89 PUBLIC :: fist_environment_type
90
91! *** Public subroutines ***
92 PUBLIC :: fist_env_get, &
96
97 CHARACTER(len=*), PARAMETER, PRIVATE :: modulen = 'fist_environment_types'
98
99!***
100
101CONTAINS
102
103! **************************************************************************************************
104!> \brief Purpose: Get the FIST environment.
105!> \param fist_env the pointer to the fist_env
106!> \param atomic_kind_set ...
107!> \param particle_set ...
108!> \param ewald_pw ...
109!> \param local_particles ...
110!> \param local_molecules ...
111!> \param molecule_kind_set ...
112!> \param molecule_set ...
113!> \param cell ...
114!> \param cell_ref ...
115!> \param ewald_env ...
116!> \param fist_nonbond_env ...
117!> \param thermo ...
118!> \param para_env ...
119!> \param subsys ...
120!> \param qmmm ...
121!> \param qmmm_env ...
122!> \param input ...
123!> \param shell_model ...
124!> \param shell_model_ad ...
125!> \param shell_particle_set ...
126!> \param core_particle_set ...
127!> \param multipoles ...
128!> \param results ...
129!> \param exclusions ...
130!> \param efield ...
131!> \par History
132!> 11/03
133!> \author CJM
134! **************************************************************************************************
135 SUBROUTINE fist_env_get(fist_env, atomic_kind_set, particle_set, ewald_pw, &
136 local_particles, local_molecules, molecule_kind_set, molecule_set, cell, &
137 cell_ref, ewald_env, fist_nonbond_env, thermo, para_env, subsys, qmmm, &
138 qmmm_env, input, shell_model, shell_model_ad, shell_particle_set, &
139 core_particle_set, multipoles, results, exclusions, efield)
140
141 TYPE(fist_environment_type), INTENT(IN) :: fist_env
142 TYPE(atomic_kind_type), OPTIONAL, POINTER :: atomic_kind_set(:)
143 TYPE(particle_type), OPTIONAL, POINTER :: particle_set(:)
144 TYPE(ewald_pw_type), OPTIONAL, POINTER :: ewald_pw
145 TYPE(distribution_1d_type), OPTIONAL, POINTER :: local_particles, local_molecules
146 TYPE(molecule_kind_type), OPTIONAL, POINTER :: molecule_kind_set(:)
147 TYPE(molecule_type), OPTIONAL, POINTER :: molecule_set(:)
148 TYPE(cell_type), OPTIONAL, POINTER :: cell, cell_ref
149 TYPE(ewald_environment_type), OPTIONAL, POINTER :: ewald_env
150 TYPE(fist_nonbond_env_type), OPTIONAL, POINTER :: fist_nonbond_env
151 TYPE(fist_energy_type), OPTIONAL, POINTER :: thermo
152 TYPE(mp_para_env_type), OPTIONAL, POINTER :: para_env
153 TYPE(cp_subsys_type), OPTIONAL, POINTER :: subsys
154 LOGICAL, OPTIONAL :: qmmm
155 TYPE(qmmm_env_mm_type), OPTIONAL, POINTER :: qmmm_env
156 TYPE(section_vals_type), OPTIONAL, POINTER :: input
157 LOGICAL, OPTIONAL :: shell_model, shell_model_ad
158 TYPE(particle_type), OPTIONAL, POINTER :: shell_particle_set(:), &
159 core_particle_set(:)
160 TYPE(multipole_type), OPTIONAL, POINTER :: multipoles
161 TYPE(cp_result_type), OPTIONAL, POINTER :: results
162 TYPE(exclusion_type), DIMENSION(:), OPTIONAL, &
163 POINTER :: exclusions
164 TYPE(fist_efield_type), OPTIONAL, POINTER :: efield
165
166 TYPE(atomic_kind_list_type), POINTER :: atomic_kinds
167 TYPE(molecule_kind_list_type), POINTER :: molecule_kinds
168 TYPE(molecule_list_type), POINTER :: molecules
169 TYPE(multipole_type), POINTER :: fist_multipoles
170 TYPE(particle_list_type), POINTER :: core_particles, particles, &
171 shell_particles
172
173 NULLIFY (atomic_kinds, particles, molecules, molecule_kinds, fist_multipoles)
174
175 IF (PRESENT(input)) input => fist_env%input
176 IF (PRESENT(qmmm)) qmmm = fist_env%qmmm
177 IF (PRESENT(qmmm_env)) qmmm_env => fist_env%qmmm_env
178 IF (PRESENT(cell_ref)) cell_ref => fist_env%cell_ref
179 IF (PRESENT(ewald_env)) ewald_env => fist_env%ewald_env
180 IF (PRESENT(thermo)) thermo => fist_env%thermo
181 IF (PRESENT(exclusions)) exclusions => fist_env%exclusions
182 IF (PRESENT(para_env)) para_env => fist_env%para_env
183 IF (PRESENT(ewald_pw)) ewald_pw => fist_env%ewald_pw
184 IF (PRESENT(fist_nonbond_env)) fist_nonbond_env => fist_env%fist_nonbond_env
185 IF (PRESENT(shell_model)) shell_model = fist_env%shell_model
186 IF (PRESENT(shell_model_ad)) shell_model_ad = fist_env%shell_model_ad
187 IF (PRESENT(subsys)) subsys => fist_env%subsys
188 IF (PRESENT(efield)) efield => fist_env%efield
189
190 IF (ASSOCIATED(fist_env%subsys)) &
191 CALL cp_subsys_get(fist_env%subsys, &
192 atomic_kinds=atomic_kinds, &
193 local_molecules=local_molecules, &
194 local_particles=local_particles, &
195 particles=particles, &
196 molecule_kinds=molecule_kinds, &
197 molecules=molecules, &
198 shell_particles=shell_particles, &
199 core_particles=core_particles, &
200 multipoles=fist_multipoles, &
201 results=results, &
202 cell=cell)
203 IF (PRESENT(atomic_kind_set)) atomic_kind_set => atomic_kinds%els
204 IF (PRESENT(particle_set)) particle_set => particles%els
205 IF (PRESENT(molecule_kind_set)) molecule_kind_set => molecule_kinds%els
206 IF (PRESENT(molecule_set)) molecule_set => molecules%els
207 IF (PRESENT(shell_particle_set)) shell_particle_set => shell_particles%els
208 IF (PRESENT(core_particle_set)) core_particle_set => core_particles%els
209 IF (PRESENT(multipoles)) multipoles => fist_multipoles
210 END SUBROUTINE fist_env_get
211
212! **************************************************************************************************
213!> \brief Initialise the FIST environment.
214!> \param fist_env the pointer to the fist_env
215!> \param para_env ...
216!> \par History
217!> 11/03
218!> \author CJM
219! **************************************************************************************************
220 SUBROUTINE init_fist_env(fist_env, para_env)
221
222 TYPE(fist_environment_type), INTENT(OUT) :: fist_env
223 TYPE(mp_para_env_type), POINTER :: para_env
224
225 NULLIFY (fist_env%input)
226 NULLIFY (fist_env%qmmm_env)
227 NULLIFY (fist_env%cell_ref)
228 NULLIFY (fist_env%ewald_env)
229 NULLIFY (fist_env%ewald_pw)
230 NULLIFY (fist_env%thermo)
231 NULLIFY (fist_env%fist_nonbond_env)
232 NULLIFY (fist_env%subsys)
233 NULLIFY (fist_env%exclusions)
234 NULLIFY (fist_env%efield)
235 fist_env%qmmm = .false.
236 fist_env%shell_model = .false.
237 fist_env%shell_model_ad = .false.
238 ALLOCATE (fist_env%qmmm_env)
239 CALL qmmm_env_mm_create(fist_env%qmmm_env)
240 NULLIFY (fist_env%subsys)
241 CALL para_env%retain()
242 fist_env%para_env => para_env
243
244 END SUBROUTINE init_fist_env
245
246! **************************************************************************************************
247!> \brief Set the FIST environment.
248!> \param fist_env the pointer to the fist_env
249!> \param atomic_kind_set ...
250!> \param particle_set ...
251!> \param ewald_pw ...
252!> \param local_particles ...
253!> \param local_molecules ...
254!> \param molecule_kind_set ...
255!> \param molecule_set ...
256!> \param cell_ref ...
257!> \param ewald_env ...
258!> \param fist_nonbond_env ...
259!> \param thermo ...
260!> \param subsys ...
261!> \param qmmm ...
262!> \param qmmm_env ...
263!> \param input ...
264!> \param shell_model ...
265!> \param shell_model_ad ...
266!> \param exclusions ...
267!> \param efield ...
268!> \par History
269!> 11/03
270!> \author CJM
271! **************************************************************************************************
272 SUBROUTINE fist_env_set(fist_env, atomic_kind_set, particle_set, ewald_pw, &
273 local_particles, local_molecules, molecule_kind_set, &
274 molecule_set, cell_ref, ewald_env, &
275 fist_nonbond_env, thermo, subsys, qmmm, qmmm_env, &
276 input, shell_model, shell_model_ad, exclusions, efield)
277
278 TYPE(fist_environment_type), INTENT(INOUT) :: fist_env
279 TYPE(atomic_kind_type), OPTIONAL, POINTER :: atomic_kind_set(:)
280 TYPE(particle_type), OPTIONAL, POINTER :: particle_set(:)
281 TYPE(ewald_pw_type), OPTIONAL, POINTER :: ewald_pw
282 TYPE(distribution_1d_type), OPTIONAL, POINTER :: local_particles, local_molecules
283 TYPE(molecule_kind_type), OPTIONAL, POINTER :: molecule_kind_set(:)
284 TYPE(molecule_type), OPTIONAL, POINTER :: molecule_set(:)
285 TYPE(cell_type), OPTIONAL, POINTER :: cell_ref
286 TYPE(ewald_environment_type), OPTIONAL, POINTER :: ewald_env
287 TYPE(fist_nonbond_env_type), OPTIONAL, POINTER :: fist_nonbond_env
288 TYPE(fist_energy_type), OPTIONAL, POINTER :: thermo
289 TYPE(cp_subsys_type), OPTIONAL, POINTER :: subsys
290 LOGICAL, OPTIONAL :: qmmm
291 TYPE(qmmm_env_mm_type), OPTIONAL, POINTER :: qmmm_env
292 TYPE(section_vals_type), OPTIONAL, POINTER :: input
293 LOGICAL, OPTIONAL :: shell_model, shell_model_ad
294 TYPE(exclusion_type), DIMENSION(:), OPTIONAL, &
295 POINTER :: exclusions
296 TYPE(fist_efield_type), OPTIONAL, POINTER :: efield
297
298 TYPE(atomic_kind_list_type), POINTER :: atomic_kinds
299 TYPE(molecule_kind_list_type), POINTER :: molecule_kinds
300 TYPE(molecule_list_type), POINTER :: molecules
301 TYPE(particle_list_type), POINTER :: particles
302
303 IF (PRESENT(qmmm)) fist_env%qmmm = qmmm
304 IF (PRESENT(qmmm_env)) THEN
305 IF (ASSOCIATED(fist_env%qmmm_env)) THEN
306 CALL qmmm_env_mm_release(fist_env%qmmm_env)
307 DEALLOCATE (fist_env%qmmm_env)
308 END IF
309 fist_env%qmmm_env => qmmm_env
310 END IF
311 IF (PRESENT(ewald_env)) THEN
312 IF (ASSOCIATED(fist_env%ewald_env)) THEN
313 IF (.NOT. ASSOCIATED(fist_env%ewald_env, ewald_env)) THEN
314 CALL ewald_env_release(fist_env%ewald_env)
315 DEALLOCATE (fist_env%ewald_env)
316 END IF
317 END IF
318 fist_env%ewald_env => ewald_env
319 END IF
320 IF (PRESENT(ewald_pw)) THEN
321 IF (ASSOCIATED(fist_env%ewald_pw)) THEN
322 IF (.NOT. ASSOCIATED(fist_env%ewald_pw, ewald_pw)) THEN
323 CALL ewald_pw_release(fist_env%ewald_pw)
324 DEALLOCATE (fist_env%ewald_pw)
325 END IF
326 END IF
327 fist_env%ewald_pw => ewald_pw
328 END IF
329 IF (PRESENT(cell_ref)) THEN
330 CALL cell_retain(cell_ref)
331 CALL cell_release(fist_env%cell_ref)
332 fist_env%cell_ref => cell_ref
333 END IF
334 IF (PRESENT(fist_nonbond_env)) THEN
335 IF (ASSOCIATED(fist_env%fist_nonbond_env)) THEN
336 IF (.NOT. ASSOCIATED(fist_env%fist_nonbond_env, fist_nonbond_env)) THEN
337 CALL fist_nonbond_env_release(fist_env%fist_nonbond_env)
338 DEALLOCATE (fist_env%fist_nonbond_env)
339 END IF
340 END IF
341 fist_env%fist_nonbond_env => fist_nonbond_env
342 END IF
343 IF (PRESENT(input)) THEN
344 CALL section_vals_retain(input)
345 CALL section_vals_release(fist_env%input)
346 fist_env%input => input
347 END IF
348 IF (PRESENT(thermo)) fist_env%thermo => thermo
349 IF (PRESENT(subsys)) THEN
350 IF (ASSOCIATED(fist_env%subsys)) THEN
351 IF (.NOT. ASSOCIATED(fist_env%subsys, subsys)) THEN
352 CALL cp_subsys_release(fist_env%subsys)
353 END IF
354 END IF
355 fist_env%subsys => subsys
356 END IF
357 IF (PRESENT(atomic_kind_set)) THEN
358 CALL atomic_kind_list_create(atomic_kinds, &
359 els_ptr=atomic_kind_set)
360 CALL cp_subsys_set(fist_env%subsys, &
361 atomic_kinds=atomic_kinds)
362 CALL atomic_kind_list_release(atomic_kinds)
363 END IF
364 IF (PRESENT(particle_set)) THEN
365 CALL particle_list_create(particles, &
366 els_ptr=particle_set)
367 CALL cp_subsys_set(fist_env%subsys, &
368 particles=particles)
369 CALL particle_list_release(particles)
370 END IF
371 IF (PRESENT(local_particles)) THEN
372 CALL cp_subsys_set(fist_env%subsys, &
373 local_particles=local_particles)
374 END IF
375 IF (PRESENT(local_molecules)) THEN
376 CALL cp_subsys_set(fist_env%subsys, &
377 local_molecules=local_molecules)
378 END IF
379 IF (PRESENT(molecule_kind_set)) THEN
380 CALL molecule_kind_list_create(molecule_kinds, &
381 els_ptr=molecule_kind_set)
382 CALL cp_subsys_set(fist_env%subsys, &
383 molecule_kinds=molecule_kinds)
384 CALL molecule_kind_list_release(molecule_kinds)
385 END IF
386 IF (PRESENT(molecule_set)) THEN
387 CALL molecule_list_create(molecules, &
388 els_ptr=molecule_set)
389 CALL cp_subsys_set(fist_env%subsys, &
390 molecules=molecules)
391 CALL molecule_list_release(molecules)
392 END IF
393 IF (PRESENT(exclusions)) fist_env%exclusions => exclusions
394 IF (PRESENT(shell_model)) THEN
395 fist_env%shell_model = shell_model
396 END IF
397 IF (PRESENT(shell_model_ad)) THEN
398 fist_env%shell_model_ad = shell_model_ad
399 END IF
400 IF (PRESENT(efield)) fist_env%efield => efield
401
402 END SUBROUTINE fist_env_set
403
404! **************************************************************************************************
405!> \brief allocates and intitializes a fist_env
406!> \param fist_env the object to create
407!> \param para_env the parallel environment for the qs_env
408!> \par History
409!> 12.2002 created [fawzi]
410!> \author Fawzi Mohamed
411! **************************************************************************************************
412 SUBROUTINE fist_env_create(fist_env, para_env)
413 TYPE(fist_environment_type), INTENT(OUT) :: fist_env
414 TYPE(mp_para_env_type), INTENT(IN), POINTER :: para_env
415
416 CALL init_fist_env(fist_env, para_env=para_env)
417 END SUBROUTINE fist_env_create
418
419! **************************************************************************************************
420!> \brief releases the given fist_env (see doc/ReferenceCounting.html)
421!> \param fist_env the object to release
422!> \par History
423!> 12.2002 created [fawzi]
424!> \author Fawzi Mohamed
425! **************************************************************************************************
426 SUBROUTINE fist_env_release(fist_env)
427 TYPE(fist_environment_type), INTENT(INOUT) :: fist_env
428
429 IF (ASSOCIATED(fist_env%qmmm_env)) THEN
430 CALL qmmm_env_mm_release(fist_env%qmmm_env)
431 DEALLOCATE (fist_env%qmmm_env)
432 END IF
433 CALL cell_release(fist_env%cell_ref)
434 IF (ASSOCIATED(fist_env%ewald_pw)) THEN
435 CALL ewald_pw_release(fist_env%ewald_pw)
436 DEALLOCATE (fist_env%ewald_pw)
437 END IF
438 IF (ASSOCIATED(fist_env%ewald_env)) THEN
439 CALL ewald_env_release(fist_env%ewald_env)
440 DEALLOCATE (fist_env%ewald_env)
441 END IF
442 CALL mp_para_env_release(fist_env%para_env)
443 CALL deallocate_fist_energy(fist_env%thermo)
444
445 IF (ASSOCIATED(fist_env%fist_nonbond_env)) THEN
446 CALL fist_nonbond_env_release(fist_env%fist_nonbond_env)
447 DEALLOCATE (fist_env%fist_nonbond_env)
448 END IF
449 CALL cp_subsys_release(fist_env%subsys)
450 CALL section_vals_release(fist_env%input)
451 CALL exclusion_release(fist_env%exclusions)
452
453 IF (ASSOCIATED(fist_env%efield)) THEN
454 DEALLOCATE (fist_env%efield)
455 END IF
456
457 END SUBROUTINE fist_env_release
458
459END MODULE fist_environment_types
represent a simple array based list of the given type
subroutine, public atomic_kind_list_release(list)
releases a list (see doc/ReferenceCounting.html)
subroutine, public atomic_kind_list_create(list, els_ptr, owns_els, n_els)
creates a list
Define the atomic kind types and their sub types.
Handles all functions related to the CELL.
Definition cell_types.F:15
subroutine, public cell_release(cell)
releases the given cell (see doc/ReferenceCounting.html)
Definition cell_types.F:559
subroutine, public cell_retain(cell)
retains the given cell (see doc/ReferenceCounting.html)
Definition cell_types.F:542
set of type/routines to handle the storage of results in force_envs
types that represent a subsys, i.e. a part of the system
subroutine, public cp_subsys_release(subsys)
releases a subsys (see doc/ReferenceCounting.html)
subroutine, public cp_subsys_set(subsys, atomic_kinds, particles, local_particles, molecules, molecule_kinds, local_molecules, para_env, colvar_p, shell_particles, core_particles, gci, multipoles, results, cell)
sets various propreties of the subsys
subroutine, public cp_subsys_get(subsys, ref_count, atomic_kinds, atomic_kind_set, particles, particle_set, local_particles, molecules, molecule_set, molecule_kinds, molecule_kind_set, local_molecules, para_env, colvar_p, shell_particles, core_particles, gci, multipoles, natom, nparticle, ncore, nshell, nkind, atprop, virial, results, cell)
returns information about various attributes of the given subsys
stores a lists of integer that are local to a processor. The idea is that these integers represent ob...
subroutine, public ewald_env_release(ewald_env)
releases the given ewald_env (see doc/ReferenceCounting.html)
subroutine, public ewald_pw_release(ewald_pw)
releases the memory used by the ewald_pw
an exclusion type
subroutine, public exclusion_release(exclusions)
Release exclusion type.
subroutine, public deallocate_fist_energy(fist_energy)
Deallocate a Fist energy data structure.
subroutine, public fist_env_release(fist_env)
releases the given fist_env (see doc/ReferenceCounting.html)
subroutine, public fist_env_get(fist_env, atomic_kind_set, particle_set, ewald_pw, local_particles, local_molecules, molecule_kind_set, molecule_set, cell, cell_ref, ewald_env, fist_nonbond_env, thermo, para_env, subsys, qmmm, qmmm_env, input, shell_model, shell_model_ad, shell_particle_set, core_particle_set, multipoles, results, exclusions, efield)
Purpose: Get the FIST environment.
subroutine, public fist_env_create(fist_env, para_env)
allocates and intitializes a fist_env
subroutine, public fist_env_set(fist_env, atomic_kind_set, particle_set, ewald_pw, local_particles, local_molecules, molecule_kind_set, molecule_set, cell_ref, ewald_env, fist_nonbond_env, thermo, subsys, qmmm, qmmm_env, input, shell_model, shell_model_ad, exclusions, efield)
Set the FIST environment.
subroutine, public fist_nonbond_env_release(fist_nonbond_env)
releases the given fist_nonbond_env (see doc/ReferenceCounting.html)
objects that represent the structure of input sections and the data contained in an input section
subroutine, public section_vals_retain(section_vals)
retains the given section values (see doc/ReferenceCounting.html)
recursive subroutine, public section_vals_release(section_vals)
releases the given object
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)
represent a simple array based list of the given type
subroutine, public molecule_kind_list_create(list, els_ptr, owns_els, n_els)
creates a list
subroutine, public molecule_kind_list_release(list)
releases a list (see doc/ReferenceCounting.html)
Define the molecule kind structure types and the corresponding functionality.
represent a simple array based list of the given type
subroutine, public molecule_list_create(list, els_ptr, owns_els, n_els)
creates a list
subroutine, public molecule_list_release(list)
releases a list (see doc/ReferenceCounting.html)
Define the data structure for the molecule information.
Multipole structure: for multipole (fixed and induced) in FF based MD.
represent a simple array based list of the given type
subroutine, public particle_list_create(list, els_ptr, owns_els, n_els)
creates a list
subroutine, public particle_list_release(list)
releases a list (see doc/ReferenceCounting.html)
Define the data structure for the particle information.
subroutine, public qmmm_env_mm_create(qmmm_env)
...
subroutine, public qmmm_env_mm_release(qmmm_env)
releases the given qmmm_env (see doc/ReferenceCounting.html)
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
Definition cell_types.F:55
contains arbitrary information which need to be stored
represents a system: atoms, molecules, their pos,vel,...
structure to store local (to a processor) ordered lists of integers.
A type used to store lists of exclusions and onfos.
stores all the informations relevant to an mpi environment