(git:374b731)
Loading...
Searching...
No Matches
particle_methods.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 Define methods related to particle_type
10!> \par History
11!> 10.2014 Move routines out of particle_types.F [Ole Schuett]
12!> \author Ole Schuett
13! **************************************************************************************************
18 USE cell_methods, ONLY: cell_create,&
20 USE cell_types, ONLY: cell_clone,&
22 cell_type,&
23 get_cell,&
24 pbc,&
34 USE input_constants, ONLY: dump_atomic,&
35 dump_dcd,&
37 dump_pdb,&
42 USE kinds, ONLY: default_string_length,&
43 dp,&
44 sp
45 USE mathconstants, ONLY: degree
46 USE mathlib, ONLY: angle,&
51 USE physcon, ONLY: massunit
53 USE qs_kind_types, ONLY: get_qs_kind,&
58 USE util, ONLY: sort,&
60#include "./base/base_uses.f90"
61
62 IMPLICIT NONE
63
64 PRIVATE
65
66 ! Public subroutines
67
75
76 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'particle_methods'
77
78CONTAINS
79
80! **************************************************************************************************
81!> \brief Get the components of a particle set.
82!> \param particle_set ...
83!> \param qs_kind_set ...
84!> \param first_sgf ...
85!> \param last_sgf ...
86!> \param nsgf ...
87!> \param nmao ...
88!> \param basis ...
89!> \date 14.01.2002
90!> \par History
91!> - particle type cleaned (13.10.2003,MK)
92!> - refactoring and add basis set option (17.08.2010,jhu)
93!> \author MK
94!> \version 1.0
95! **************************************************************************************************
96 SUBROUTINE get_particle_set(particle_set, qs_kind_set, first_sgf, last_sgf, nsgf, &
97 nmao, basis)
98
99 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
100 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
101 INTEGER, DIMENSION(:), INTENT(INOUT), OPTIONAL :: first_sgf, last_sgf, nsgf, nmao
102 TYPE(gto_basis_set_p_type), DIMENSION(:), OPTIONAL :: basis
103
104 INTEGER :: ikind, iparticle, isgf, nparticle, ns
105
106 cpassert(ASSOCIATED(particle_set))
107
108 nparticle = SIZE(particle_set)
109 IF (PRESENT(first_sgf)) THEN
110 cpassert(SIZE(first_sgf) >= nparticle)
111 END IF
112 IF (PRESENT(last_sgf)) THEN
113 cpassert(SIZE(last_sgf) >= nparticle)
114 END IF
115 IF (PRESENT(nsgf)) THEN
116 cpassert(SIZE(nsgf) >= nparticle)
117 END IF
118 IF (PRESENT(nmao)) THEN
119 cpassert(SIZE(nmao) >= nparticle)
120 END IF
121
122 IF (PRESENT(first_sgf) .OR. PRESENT(last_sgf) .OR. PRESENT(nsgf)) THEN
123 isgf = 0
124 DO iparticle = 1, nparticle
125 CALL get_atomic_kind(particle_set(iparticle)%atomic_kind, kind_number=ikind)
126 IF (PRESENT(basis)) THEN
127 IF (ASSOCIATED(basis(ikind)%gto_basis_set)) THEN
128 CALL get_gto_basis_set(gto_basis_set=basis(ikind)%gto_basis_set, nsgf=ns)
129 ELSE
130 ns = 0
131 END IF
132 ELSE
133 CALL get_qs_kind(qs_kind_set(ikind), nsgf=ns)
134 END IF
135 IF (PRESENT(nsgf)) nsgf(iparticle) = ns
136 IF (PRESENT(first_sgf)) first_sgf(iparticle) = isgf + 1
137 isgf = isgf + ns
138 IF (PRESENT(last_sgf)) last_sgf(iparticle) = isgf
139 END DO
140 END IF
141
142 IF (PRESENT(first_sgf)) THEN
143 IF (SIZE(first_sgf) > nparticle) first_sgf(nparticle + 1) = isgf + 1
144 END IF
145
146 IF (PRESENT(nmao)) THEN
147 DO iparticle = 1, nparticle
148 CALL get_atomic_kind(particle_set(iparticle)%atomic_kind, kind_number=ikind)
149 CALL get_qs_kind(qs_kind_set(ikind), mao=ns)
150 nmao(iparticle) = ns
151 END DO
152 END IF
153
154 END SUBROUTINE get_particle_set
155
156! **************************************************************************************************
157!> \brief Should be able to write a few formats e.g. xmol, and some binary
158!> format (dcd) some format can be used for x,v,f
159!>
160!> FORMAT CONTENT UNITS x, v, f
161!> XMOL POS, VEL, FORCE, POS_VEL, POS_VEL_FORCE Angstrom, a.u., a.u.
162!>
163!> \param particle_set ...
164!> \param iunit ...
165!> \param output_format ...
166!> \param content ...
167!> \param title ...
168!> \param cell ...
169!> \param array ...
170!> \param unit_conv ...
171!> \param charge_occup ...
172!> \param charge_beta ...
173!> \param charge_extended ...
174!> \param print_kind ...
175!> \date 14.01.2002
176!> \author MK
177!> \version 1.0
178! **************************************************************************************************
179 SUBROUTINE write_particle_coordinates(particle_set, iunit, output_format, &
180 content, title, cell, array, unit_conv, &
181 charge_occup, charge_beta, &
182 charge_extended, print_kind)
183
184 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
185 INTEGER :: iunit, output_format
186 CHARACTER(LEN=*) :: content, title
187 TYPE(cell_type), OPTIONAL, POINTER :: cell
188 REAL(kind=dp), DIMENSION(:), INTENT(IN), OPTIONAL :: array
189 REAL(kind=dp), INTENT(IN), OPTIONAL :: unit_conv
190 LOGICAL, INTENT(IN), OPTIONAL :: charge_occup, charge_beta, &
191 charge_extended, print_kind
192
193 CHARACTER(len=*), PARAMETER :: routinen = 'write_particle_coordinates'
194
195 CHARACTER(LEN=120) :: line
196 CHARACTER(LEN=2) :: element_symbol
197 CHARACTER(LEN=4) :: name
198 CHARACTER(LEN=default_string_length) :: atm_name, my_format
199 INTEGER :: handle, iatom, natom
200 LOGICAL :: dummy, my_charge_beta, &
201 my_charge_extended, my_charge_occup, &
202 my_print_kind
203 REAL(kind=dp) :: angle_alpha, angle_beta, angle_gamma, &
204 factor, qeff
205 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: arr
206 REAL(kind=dp), DIMENSION(3) :: abc, angles, f, r, v
207 REAL(kind=dp), DIMENSION(3, 3) :: h
208 REAL(kind=sp), ALLOCATABLE, DIMENSION(:) :: x4, y4, z4
209 TYPE(cell_type), POINTER :: cell_dcd
210 TYPE(fist_potential_type), POINTER :: fist_potential
211 TYPE(shell_kind_type), POINTER :: shell
212
213 CALL timeset(routinen, handle)
214
215 natom = SIZE(particle_set)
216 IF (PRESENT(array)) THEN
217 SELECT CASE (trim(content))
218 CASE ("POS_VEL", "POS_VEL_FORCE")
219 cpabort("Illegal usage")
220 END SELECT
221 END IF
222 factor = 1.0_dp
223 IF (PRESENT(unit_conv)) THEN
224 factor = unit_conv
225 END IF
226 SELECT CASE (output_format)
227 CASE (dump_xmol)
228 my_print_kind = .false.
229 IF (PRESENT(print_kind)) my_print_kind = print_kind
230 WRITE (iunit, "(I8)") natom
231 WRITE (iunit, "(A)") trim(title)
232 DO iatom = 1, natom
233 CALL get_atomic_kind(atomic_kind=particle_set(iatom)%atomic_kind, &
234 element_symbol=element_symbol)
235 IF (len_trim(element_symbol) == 0 .OR. my_print_kind) THEN
236 CALL get_atomic_kind(atomic_kind=particle_set(iatom)%atomic_kind, &
237 name=atm_name)
238 dummy = qmmm_ff_precond_only_qm(id1=atm_name)
239 my_format = "(A,"
240 atm_name = trim(atm_name)
241 ELSE
242 my_format = "(T2,A2,"
243 atm_name = trim(element_symbol)
244 END IF
245 SELECT CASE (trim(content))
246 CASE ("POS")
247 IF (PRESENT(array)) THEN
248 r(1:3) = get_particle_pos_or_vel(iatom, particle_set, array)
249 ELSE
250 r(:) = particle_set(iatom)%r(:)
251 END IF
252 WRITE (iunit, trim(my_format)//"1X,3F20.10)") trim(atm_name), r(1:3)*factor
253 CASE ("VEL")
254 IF (PRESENT(array)) THEN
255 v(1:3) = get_particle_pos_or_vel(iatom, particle_set, array)
256 ELSE
257 v(:) = particle_set(iatom)%v(:)
258 END IF
259 WRITE (iunit, trim(my_format)//"1X,3F20.10)") trim(atm_name), v(1:3)*factor
260 CASE ("FORCE")
261 IF (PRESENT(array)) THEN
262 f(:) = array((iatom - 1)*3 + 1:(iatom - 1)*3 + 3)
263 ELSE
264 f(:) = particle_set(iatom)%f(:)
265 END IF
266 WRITE (iunit, trim(my_format)//"1X,3F20.10)") trim(atm_name), f(1:3)*factor
267 CASE ("FORCE_MIXING_LABELS")
268 IF (PRESENT(array)) THEN
269 f(:) = array((iatom - 1)*3 + 1:(iatom - 1)*3 + 3)
270 ELSE
271 f(:) = particle_set(iatom)%f(:)
272 END IF
273 WRITE (iunit, trim(my_format)//"1X,3F20.10)") trim(atm_name), f(1:3)*factor
274 END SELECT
275 END DO
276 CASE (dump_atomic)
277 DO iatom = 1, natom
278 SELECT CASE (trim(content))
279 CASE ("POS")
280 IF (PRESENT(array)) THEN
281 r(1:3) = get_particle_pos_or_vel(iatom, particle_set, array)
282 ELSE
283 r(:) = particle_set(iatom)%r(:)
284 END IF
285 WRITE (iunit, "(3F20.10)") r(1:3)*factor
286 CASE ("VEL")
287 IF (PRESENT(array)) THEN
288 v(1:3) = get_particle_pos_or_vel(iatom, particle_set, array)
289 ELSE
290 v(:) = particle_set(iatom)%v(:)
291 END IF
292 WRITE (iunit, "(3F20.10)") v(1:3)*factor
293 CASE ("FORCE")
294 IF (PRESENT(array)) THEN
295 f(:) = array((iatom - 1)*3 + 1:(iatom - 1)*3 + 3)
296 ELSE
297 f(:) = particle_set(iatom)%f(:)
298 END IF
299 WRITE (iunit, "(3F20.10)") f(1:3)*factor
300 CASE ("FORCE_MIXING_LABELS")
301 IF (PRESENT(array)) THEN
302 f(:) = array((iatom - 1)*3 + 1:(iatom - 1)*3 + 3)
303 ELSE
304 f(:) = particle_set(iatom)%f(:)
305 END IF
306 WRITE (iunit, "(3F20.10)") f(1:3)*factor
307 END SELECT
308 END DO
310 IF (.NOT. (PRESENT(cell))) THEN
311 cpabort("Cell is not present! Report this bug!")
312 END IF
313 CALL get_cell(cell, alpha=angle_alpha, beta=angle_beta, gamma=angle_gamma, &
314 abc=abc)
315 IF (.NOT. cell%orthorhombic .AND. (output_format == dump_dcd_aligned_cell)) THEN
316 ! In the case of a non-orthorhombic cell adopt a common convention
317 ! for the orientation of the cell with respect to the Cartesian axes:
318 ! Cell vector a is aligned with the x axis and the cell vector b lies
319 ! in the xy plane.
320 NULLIFY (cell_dcd)
321 CALL cell_create(cell_dcd)
322 CALL cell_clone(cell, cell_dcd, tag="CELL_DCD")
323 angles(1) = angle_alpha/degree
324 angles(2) = angle_beta/degree
325 angles(3) = angle_gamma/degree
326 CALL set_cell_param(cell_dcd, abc, angles, &
327 do_init_cell=.true.)
328 h(1:3, 1:3) = matmul(cell_dcd%hmat(1:3, 1:3), cell%h_inv(1:3, 1:3))
329 CALL cell_release(cell_dcd)
330 END IF
331 ALLOCATE (arr(3, natom))
332 IF (PRESENT(array)) THEN
333 arr(1:3, 1:natom) = reshape(array, (/3, natom/))
334 ELSE
335 SELECT CASE (trim(content))
336 CASE ("POS")
337 DO iatom = 1, natom
338 arr(1:3, iatom) = particle_set(iatom)%r(1:3)
339 END DO
340 CASE ("VEL")
341 DO iatom = 1, natom
342 arr(1:3, iatom) = particle_set(iatom)%v(1:3)
343 END DO
344 CASE ("FORCE")
345 DO iatom = 1, natom
346 arr(1:3, iatom) = particle_set(iatom)%f(1:3)
347 END DO
348 CASE DEFAULT
349 cpabort("Illegal DCD dump type")
350 END SELECT
351 END IF
352 ALLOCATE (x4(natom))
353 ALLOCATE (y4(natom))
354 ALLOCATE (z4(natom))
355 IF (.NOT. cell%orthorhombic .AND. (output_format == dump_dcd_aligned_cell)) THEN
356 x4(1:natom) = real(matmul(h(1, 1:3), arr(1:3, 1:natom)), kind=sp)
357 y4(1:natom) = real(matmul(h(2, 1:3), arr(1:3, 1:natom)), kind=sp)
358 z4(1:natom) = real(matmul(h(3, 1:3), arr(1:3, 1:natom)), kind=sp)
359 ELSE
360 x4(1:natom) = real(arr(1, 1:natom), kind=sp)
361 y4(1:natom) = real(arr(2, 1:natom), kind=sp)
362 z4(1:natom) = real(arr(3, 1:natom), kind=sp)
363 END IF
364 WRITE (iunit) abc(1)*factor, angle_gamma, abc(2)*factor, &
365 angle_beta, angle_alpha, abc(3)*factor
366 WRITE (iunit) x4*real(factor, kind=sp)
367 WRITE (iunit) y4*real(factor, kind=sp)
368 WRITE (iunit) z4*real(factor, kind=sp)
369 ! Release work storage
370 DEALLOCATE (arr)
371 DEALLOCATE (x4)
372 DEALLOCATE (y4)
373 DEALLOCATE (z4)
374 CASE (dump_pdb)
375 my_charge_occup = .false.
376 IF (PRESENT(charge_occup)) my_charge_occup = charge_occup
377 my_charge_beta = .false.
378 IF (PRESENT(charge_beta)) my_charge_beta = charge_beta
379 my_charge_extended = .false.
380 IF (PRESENT(charge_extended)) my_charge_extended = charge_extended
381 IF (len_trim(title) > 0) THEN
382 WRITE (unit=iunit, fmt="(A6,T11,A)") &
383 "REMARK", trim(title)
384 END IF
385 CALL get_cell(cell, alpha=angle_alpha, beta=angle_beta, gamma=angle_gamma, abc=abc)
386 ! COLUMNS DATA TYPE CONTENTS
387 ! --------------------------------------------------
388 ! 1 - 6 Record name "CRYST1"
389 ! 7 - 15 Real(9.3) a (Angstroms)
390 ! 16 - 24 Real(9.3) b (Angstroms)
391 ! 25 - 33 Real(9.3) c (Angstroms)
392 ! 34 - 40 Real(7.2) alpha (degrees)
393 ! 41 - 47 Real(7.2) beta (degrees)
394 ! 48 - 54 Real(7.2) gamma (degrees)
395 ! 56 - 66 LString Space group
396 ! 67 - 70 Integer Z value
397 WRITE (unit=iunit, fmt="(A6,3F9.3,3F7.2)") &
398 "CRYST1", abc(1:3)*factor, angle_alpha, angle_beta, angle_gamma
399 WRITE (unit=line(1:6), fmt="(A6)") "ATOM "
400 DO iatom = 1, natom
401 line = ""
402 ! COLUMNS DATA TYPE CONTENTS
403 ! 1 - 6 Record name "ATOM "
404 ! 7 - 11 Integer Atom serial number
405 ! 13 - 16 Atom Atom name
406 ! 17 Character Alternate location indicator
407 ! 18 - 20 Residue name Residue name
408 ! 22 Character Chain identifier
409 ! 23 - 26 Integer Residue sequence number
410 ! 27 AChar Code for insertion of residues
411 ! 31 - 38 Real(8.3) Orthogonal coordinates for X in Angstrom
412 ! 39 - 46 Real(8.3) Orthogonal coordinates for Y in Angstrom
413 ! 47 - 54 Real(8.3) Orthogonal coordinates for Z in Angstrom
414 ! 55 - 60 Real(6.2) Occupancy
415 ! 61 - 66 Real(6.2) Temperature factor (Default = 0.0)
416 ! 73 - 76 LString(4) Segment identifier, left-justified
417 ! 77 - 78 LString(2) Element symbol, right-justified
418 ! 79 - 80 LString(2) Charge on the atom
419 CALL get_atomic_kind(atomic_kind=particle_set(iatom)%atomic_kind, &
420 element_symbol=element_symbol, name=atm_name, &
421 fist_potential=fist_potential, shell=shell)
422 IF (len_trim(element_symbol) == 0) THEN
423 dummy = qmmm_ff_precond_only_qm(id1=atm_name)
424 END IF
425 name = trim(atm_name)
426 IF (ASSOCIATED(fist_potential)) THEN
427 CALL get_potential(potential=fist_potential, qeff=qeff)
428 ELSE
429 qeff = 0.0_dp
430 END IF
431 IF (ASSOCIATED(shell)) CALL get_shell(shell=shell, charge=qeff)
432 WRITE (unit=line(1:6), fmt="(A6)") "ATOM "
433 WRITE (unit=line(7:11), fmt="(I5)") modulo(iatom, 100000)
434 WRITE (unit=line(13:16), fmt="(A4)") adjustl(name)
435 ! WRITE (UNIT=line(18:20),FMT="(A3)") TRIM(resname)
436 ! WRITE (UNIT=line(23:26),FMT="(I4)") MODULO(idres,10000)
437 SELECT CASE (trim(content))
438 CASE ("POS")
439 IF (PRESENT(array)) THEN
440 r(1:3) = get_particle_pos_or_vel(iatom, particle_set, array)
441 ELSE
442 r(:) = particle_set(iatom)%r(:)
443 END IF
444 WRITE (unit=line(31:54), fmt="(3F8.3)") r(1:3)*factor
445 CASE DEFAULT
446 cpabort("PDB dump only for trajectory available")
447 END SELECT
448 IF (my_charge_occup) THEN
449 WRITE (unit=line(55:60), fmt="(F6.2)") qeff
450 ELSE
451 WRITE (unit=line(55:60), fmt="(F6.2)") 0.0_dp
452 END IF
453 IF (my_charge_beta) THEN
454 WRITE (unit=line(61:66), fmt="(F6.2)") qeff
455 ELSE
456 WRITE (unit=line(61:66), fmt="(F6.2)") 0.0_dp
457 END IF
458 ! WRITE (UNIT=line(73:76),FMT="(A4)") ADJUSTL(TRIM(molname))
459 WRITE (unit=line(77:78), fmt="(A2)") adjustr(trim(element_symbol))
460 IF (my_charge_extended) THEN
461 WRITE (unit=line(81:), fmt="(SP,F0.8)") qeff
462 END IF
463 WRITE (unit=iunit, fmt="(A)") trim(line)
464 END DO
465 WRITE (unit=iunit, fmt="(A)") "END"
466 CASE DEFAULT
467 cpabort("Illegal dump type")
468 END SELECT
469
470 CALL timestop(handle)
471
472 END SUBROUTINE write_particle_coordinates
473
474! **************************************************************************************************
475!> \brief Write the atomic coordinates to the output unit.
476!> \param particle_set ...
477!> \param subsys_section ...
478!> \param charges ...
479!> \date 05.06.2000
480!> \author MK
481!> \version 1.0
482! **************************************************************************************************
483 SUBROUTINE write_fist_particle_coordinates(particle_set, subsys_section, charges)
484 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
485 TYPE(section_vals_type), POINTER :: subsys_section
486 REAL(kind=dp), DIMENSION(:), OPTIONAL :: charges
487
488 CHARACTER(LEN=default_string_length) :: name, unit_str
489 INTEGER :: iatom, ikind, iw, natom
490 REAL(kind=dp) :: conv, mass, qcore, qeff, qshell
491 TYPE(cp_logger_type), POINTER :: logger
492 TYPE(shell_kind_type), POINTER :: shell_kind
493
494 NULLIFY (logger)
495 NULLIFY (shell_kind)
496
497 logger => cp_get_default_logger()
498 iw = cp_print_key_unit_nr(logger, subsys_section, &
499 "PRINT%ATOMIC_COORDINATES", extension=".coordLog")
500
501 CALL section_vals_val_get(subsys_section, "PRINT%ATOMIC_COORDINATES%UNIT", c_val=unit_str)
502 conv = cp_unit_from_cp2k(1.0_dp, trim(unit_str))
503 CALL uppercase(unit_str)
504 IF (iw > 0) THEN
505 WRITE (unit=iw, fmt="(/,/,T2,A)") &
506 "MODULE FIST: ATOMIC COORDINATES IN "//trim(unit_str)
507 WRITE (unit=iw, fmt="(/,T4,A,T30,A,T44,A,T58,A,T66,A,T77,A)") &
508 "Atom Kind Name", "X", "Y", "Z", "q(eff)", "Mass"
509 natom = SIZE(particle_set)
510 DO iatom = 1, natom
511 CALL get_atomic_kind(atomic_kind=particle_set(iatom)%atomic_kind, &
512 kind_number=ikind, &
513 name=name, &
514 mass=mass, &
515 qeff=qeff, &
516 shell=shell_kind)
517 IF (PRESENT(charges)) qeff = charges(iatom)
518 IF (ASSOCIATED(shell_kind)) THEN
519 CALL get_shell(shell=shell_kind, &
520 charge_core=qcore, &
521 charge_shell=qshell)
522 qeff = qcore + qshell
523 END IF
524 WRITE (unit=iw, fmt="(T2,I6,1X,I4,1X,A7,3(1X,F13.6),2(1X,F8.4))") &
525 iatom, ikind, name, particle_set(iatom)%r(1:3)*conv, qeff, mass/massunit
526 END DO
527 WRITE (iw, "(A)") ""
528 END IF
529
530 CALL cp_print_key_finished_output(iw, logger, subsys_section, &
531 "PRINT%ATOMIC_COORDINATES")
532
534
535! **************************************************************************************************
536!> \brief Write the atomic coordinates to the output unit.
537!> \param particle_set ...
538!> \param qs_kind_set ...
539!> \param subsys_section ...
540!> \param label ...
541!> \date 05.06.2000
542!> \author MK
543!> \version 1.0
544! **************************************************************************************************
545 SUBROUTINE write_qs_particle_coordinates(particle_set, qs_kind_set, subsys_section, label)
546
547 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
548 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
549 TYPE(section_vals_type), POINTER :: subsys_section
550 CHARACTER(LEN=*), INTENT(IN) :: label
551
552 CHARACTER(len=*), PARAMETER :: routinen = 'write_qs_particle_coordinates'
553
554 CHARACTER(LEN=2) :: element_symbol
555 CHARACTER(LEN=default_string_length) :: unit_str
556 INTEGER :: handle, iatom, ikind, iw, natom, z
557 REAL(kind=dp) :: conv, mass, zeff
558 TYPE(cp_logger_type), POINTER :: logger
559
560 CALL timeset(routinen, handle)
561
562 NULLIFY (logger)
563 logger => cp_get_default_logger()
564 iw = cp_print_key_unit_nr(logger, subsys_section, &
565 "PRINT%ATOMIC_COORDINATES", extension=".coordLog")
566
567 CALL section_vals_val_get(subsys_section, "PRINT%ATOMIC_COORDINATES%UNIT", c_val=unit_str)
568 conv = cp_unit_from_cp2k(1.0_dp, trim(unit_str))
569 CALL uppercase(unit_str)
570 IF (iw > 0) THEN
571 WRITE (unit=iw, fmt="(/,/,T2,A)") &
572 "MODULE "//trim(label)//": ATOMIC COORDINATES IN "//trim(unit_str)
573 WRITE (unit=iw, fmt="(/,T4,A,T30,A,T44,A,T58,A,T66,A,T77,A)") &
574 "Atom Kind Element", "X", "Y", "Z", "Z(eff)", "Mass"
575 natom = SIZE(particle_set)
576 DO iatom = 1, natom
577 CALL get_atomic_kind(atomic_kind=particle_set(iatom)%atomic_kind, &
578 kind_number=ikind, &
579 element_symbol=element_symbol, &
580 mass=mass, &
581 z=z)
582 CALL get_qs_kind(qs_kind_set(ikind), zeff=zeff)
583 WRITE (unit=iw, fmt="(T2,I6,1X,I4,1X,A2,1X,I4,3(1X,F13.6),2(1X,F8.4))") &
584 iatom, ikind, element_symbol, z, particle_set(iatom)%r(1:3)*conv, zeff, mass/massunit
585 END DO
586 WRITE (iw, "(A)") ""
587 END IF
588
589 CALL cp_print_key_finished_output(iw, logger, subsys_section, &
590 "PRINT%ATOMIC_COORDINATES")
591
592 CALL timestop(handle)
593
594 END SUBROUTINE write_qs_particle_coordinates
595
596! **************************************************************************************************
597!> \brief Write the matrix of the particle distances to the output unit.
598!> \param particle_set ...
599!> \param cell ...
600!> \param subsys_section ...
601!> \date 06.10.2000
602!> \author Matthias Krack
603!> \version 1.0
604! **************************************************************************************************
605 SUBROUTINE write_particle_distances(particle_set, cell, subsys_section)
606
607 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
608 TYPE(cell_type), POINTER :: cell
609 TYPE(section_vals_type), POINTER :: subsys_section
610
611 CHARACTER(len=*), PARAMETER :: routinen = 'write_particle_distances'
612
613 CHARACTER(LEN=default_string_length) :: unit_str
614 INTEGER :: handle, iatom, iw, jatom, natom
615 INTEGER, DIMENSION(3) :: periodic
616 LOGICAL :: explicit
617 REAL(kind=dp) :: conv, dab, dab_abort, dab_min, dab_warn
618 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: distance_matrix
619 REAL(kind=dp), DIMENSION(3) :: rab
620 TYPE(cp_logger_type), POINTER :: logger
621
622 CALL timeset(routinen, handle)
623
624 cpassert(ASSOCIATED(particle_set))
625 cpassert(ASSOCIATED(cell))
626 cpassert(ASSOCIATED(subsys_section))
627
628 NULLIFY (logger)
629 logger => cp_get_default_logger()
630 iw = cp_print_key_unit_nr(logger, subsys_section, &
631 "PRINT%INTERATOMIC_DISTANCES", extension=".distLog")
632
633 CALL section_vals_val_get(subsys_section, "PRINT%INTERATOMIC_DISTANCES%UNIT", c_val=unit_str)
634 conv = cp_unit_from_cp2k(1.0_dp, trim(unit_str))
635 CALL section_vals_val_get(subsys_section, "PRINT%INTERATOMIC_DISTANCES%CHECK_INTERATOMIC_DISTANCES", &
636 r_val=dab_min, explicit=explicit)
637
638 dab_abort = 0.0_dp
639 dab_warn = 0.0_dp
640 natom = SIZE(particle_set)
641
642 ! Compute interatomic distances only if their printout or check is explicitly requested
643 ! Disable the default check for systems with more than 3000 atoms
644 IF (explicit .OR. (iw > 0) .OR. (natom <= 2000)) THEN
645 IF (dab_min > 0.0_dp) THEN
646 dab_warn = dab_min*conv
647 ELSE IF (dab_min < 0.0_dp) THEN
648 dab_abort = abs(dab_min)*conv
649 END IF
650 END IF
651
652 IF ((iw > 0) .OR. (dab_abort > 0.0_dp) .OR. (dab_warn > 0.0_dp)) THEN
653 CALL get_cell(cell=cell, periodic=periodic)
654 IF (iw > 0) THEN
655 ALLOCATE (distance_matrix(natom, natom))
656 distance_matrix(:, :) = 0.0_dp
657 END IF
658 DO iatom = 1, natom
659 DO jatom = iatom + 1, natom
660 rab(:) = pbc(particle_set(iatom)%r(:), &
661 particle_set(jatom)%r(:), cell)
662 dab = sqrt(rab(1)*rab(1) + rab(2)*rab(2) + rab(3)*rab(3))*conv
663 IF (dab_abort > 0.0_dp) THEN
664 ! Stop the run for interatomic distances smaller than the requested threshold
665 IF (dab < dab_abort) THEN
666 CALL cp_abort(__location__, "The distance between the atoms "// &
667 trim(adjustl(cp_to_string(iatom, fmt="(I8)")))//" and "// &
668 trim(adjustl(cp_to_string(jatom, fmt="(I8)")))//" is only "// &
669 trim(adjustl(cp_to_string(dab, fmt="(F6.3)")))//" "// &
670 trim(adjustl(unit_str))//" and thus smaller than the requested threshold of "// &
671 trim(adjustl(cp_to_string(dab_abort, fmt="(F6.3)")))//" "// &
672 trim(adjustl(unit_str)))
673 END IF
674 END IF
675 IF (dab < dab_warn) THEN
676 ! Print warning for interatomic distances smaller than the requested threshold
677 CALL cp_warn(__location__, "The distance between the atoms "// &
678 trim(adjustl(cp_to_string(iatom, fmt="(I8)")))//" and "// &
679 trim(adjustl(cp_to_string(jatom, fmt="(I8)")))//" is only "// &
680 trim(adjustl(cp_to_string(dab, fmt="(F6.3)")))//" "// &
681 trim(adjustl(unit_str))//" and thus smaller than the threshold of "// &
682 trim(adjustl(cp_to_string(dab_warn, fmt="(F6.3)")))//" "// &
683 trim(adjustl(unit_str)))
684 END IF
685 IF (iw > 0) THEN
686 distance_matrix(iatom, jatom) = dab
687 distance_matrix(jatom, iatom) = distance_matrix(iatom, jatom)
688 END IF
689 END DO
690 END DO
691 IF (iw > 0) THEN
692 ! Print the distance matrix
693 WRITE (unit=iw, fmt="(/,/,T2,A)") &
694 "INTERATOMIC DISTANCES IN "//trim(unit_str)
695 CALL write_particle_matrix(distance_matrix, particle_set, iw)
696 IF (ALLOCATED(distance_matrix)) DEALLOCATE (distance_matrix)
697 END IF
698 CALL cp_print_key_finished_output(iw, logger, subsys_section, &
699 "PRINT%INTERATOMIC_DISTANCES")
700 END IF
701
702 CALL timestop(handle)
703
704 END SUBROUTINE write_particle_distances
705
706! **************************************************************************************************
707!> \brief ...
708!> \param matrix ...
709!> \param particle_set ...
710!> \param iw ...
711!> \param el_per_part ...
712!> \param Ilist ...
713!> \param parts_per_line : number of particle columns to be printed in one line
714! **************************************************************************************************
715 SUBROUTINE write_particle_matrix(matrix, particle_set, iw, el_per_part, Ilist, parts_per_line)
716 REAL(kind=dp), DIMENSION(:, :) :: matrix
717 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
718 INTEGER, INTENT(IN) :: iw
719 INTEGER, INTENT(IN), OPTIONAL :: el_per_part
720 INTEGER, DIMENSION(:), OPTIONAL, POINTER :: ilist
721 INTEGER, INTENT(IN), OPTIONAL :: parts_per_line
722
723 CHARACTER(LEN=2) :: element_symbol
724 CHARACTER(LEN=default_string_length) :: fmt_string1, fmt_string2
725 INTEGER :: from, i, iatom, icol, jatom, katom, &
726 my_el_per_part, my_parts_per_line, &
727 natom, to
728 INTEGER, DIMENSION(:), POINTER :: my_list
729
730 my_el_per_part = 1
731 IF (PRESENT(el_per_part)) my_el_per_part = el_per_part
732 my_parts_per_line = 5
733 IF (PRESENT(parts_per_line)) my_parts_per_line = max(parts_per_line, 1)
734 WRITE (fmt_string1, fmt='(A,I0,A)') &
735 "(/,T2,11X,", my_parts_per_line, "(4X,I5,4X))"
736 WRITE (fmt_string2, fmt='(A,I0,A)') &
737 "(T2,I5,2X,A2,2X,", my_parts_per_line, "(1X,F12.6))"
738 IF (PRESENT(ilist)) THEN
739 natom = SIZE(ilist)
740 ELSE
741 natom = SIZE(particle_set)
742 END IF
743 ALLOCATE (my_list(natom))
744 IF (PRESENT(ilist)) THEN
745 my_list = ilist
746 ELSE
747 DO i = 1, natom
748 my_list(i) = i
749 END DO
750 END IF
751 natom = natom*my_el_per_part
752 DO jatom = 1, natom, my_parts_per_line
753 from = jatom
754 to = min(from + my_parts_per_line - 1, natom)
755 WRITE (unit=iw, fmt=trim(fmt_string1)) (icol, icol=from, to)
756 DO iatom = 1, natom
757 katom = iatom/my_el_per_part
758 IF (mod(iatom, my_el_per_part) /= 0) katom = katom + 1
759 CALL get_atomic_kind(atomic_kind=particle_set(my_list(katom))%atomic_kind, &
760 element_symbol=element_symbol)
761 WRITE (unit=iw, fmt=trim(fmt_string2)) &
762 iatom, element_symbol, &
763 (matrix(iatom, icol), icol=from, to)
764 END DO
765 END DO
766
767 DEALLOCATE (my_list)
768
769 END SUBROUTINE write_particle_matrix
770
771! **************************************************************************************************
772!> \brief Write structure data requested by a separate structure data input
773!> section to the output unit.
774!> input_section can be either motion_section or subsys_section.
775!>
776!> \param particle_set ...
777!> \param cell ...
778!> \param input_section ...
779!> \date 11.03.04
780!> \par History
781!> Recovered (23.03.06,MK)
782!> \author MK
783!> \version 1.0
784! **************************************************************************************************
785 SUBROUTINE write_structure_data(particle_set, cell, input_section)
786 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
787 TYPE(cell_type), POINTER :: cell
788 TYPE(section_vals_type), POINTER :: input_section
789
790 CHARACTER(LEN=*), PARAMETER :: routinen = 'write_structure_data'
791
792 CHARACTER(LEN=default_string_length) :: string, unit_str
793 INTEGER :: handle, i, i_rep, iw, n, n_rep, n_vals, &
794 natom, new_size, old_size, wrk2(2), &
795 wrk3(3), wrk4(4)
796 INTEGER, ALLOCATABLE, DIMENSION(:) :: work
797 INTEGER, DIMENSION(:), POINTER :: atomic_indices, index_list
798 LOGICAL :: unique
799 REAL(kind=dp) :: conv, dab
800 REAL(kind=dp), DIMENSION(3) :: r, rab, rbc, rcd, s
801 TYPE(cp_logger_type), POINTER :: logger
802 TYPE(section_vals_type), POINTER :: section
803
804 CALL timeset(routinen, handle)
805 NULLIFY (atomic_indices)
806 NULLIFY (index_list)
807 NULLIFY (logger)
808 NULLIFY (section)
809 string = ""
810
811 logger => cp_get_default_logger()
812 iw = cp_print_key_unit_nr(logger=logger, &
813 basis_section=input_section, &
814 print_key_path="PRINT%STRUCTURE_DATA", &
815 extension=".coordLog")
816
817 CALL section_vals_val_get(input_section, "PRINT%STRUCTURE_DATA%UNIT", c_val=unit_str)
818 conv = cp_unit_from_cp2k(1.0_dp, trim(unit_str))
819 CALL uppercase(unit_str)
820 IF (iw > 0) THEN
821 natom = SIZE(particle_set)
822 section => section_vals_get_subs_vals(section_vals=input_section, &
823 subsection_name="PRINT%STRUCTURE_DATA")
824
825 WRITE (unit=iw, fmt="(/,T2,A)") "REQUESTED STRUCTURE DATA"
826 ! Print the requested atomic position vectors
827 CALL section_vals_val_get(section_vals=section, &
828 keyword_name="POSITION", &
829 n_rep_val=n_rep)
830 IF (n_rep > 0) THEN
831 WRITE (unit=iw, fmt="(/,T3,A,/)") &
832 "Position vectors r(i) of the atoms i in "//trim(unit_str)
833 old_size = 0
834 DO i_rep = 1, n_rep
835 CALL section_vals_val_get(section_vals=section, &
836 keyword_name="POSITION", &
837 i_rep_val=i_rep, &
838 i_vals=atomic_indices)
839 n_vals = SIZE(atomic_indices)
840 new_size = old_size + n_vals
841 CALL reallocate(index_list, 1, new_size)
842 index_list(old_size + 1:new_size) = atomic_indices(1:n_vals)
843 old_size = new_size
844 END DO
845 ALLOCATE (work(new_size))
846 CALL sort(index_list, new_size, work)
847 DEALLOCATE (work)
848 DO i = 1, new_size
849 WRITE (unit=string, fmt="(A,I0,A)") "(", index_list(i), ")"
850 IF ((index_list(i) < 1) .OR. (index_list(i) > natom)) THEN
851 WRITE (unit=iw, fmt="(T3,A)") &
852 "Invalid atomic index "//trim(string)//" specified. Print request is ignored."
853 cycle
854 END IF
855 IF (i > 1) THEN
856 ! Skip redundant indices
857 IF (index_list(i) == index_list(i - 1)) cycle
858 END IF
859 WRITE (unit=iw, fmt="(T3,A,T20,A,3F13.6)") &
860 "r"//trim(string), "=", pbc(particle_set(index_list(i))%r(1:3), cell)*conv
861 END DO
862 DEALLOCATE (index_list)
863 END IF
864
865 ! Print the requested atomic position vectors in scaled coordinates
866 CALL section_vals_val_get(section_vals=section, &
867 keyword_name="POSITION_SCALED", &
868 n_rep_val=n_rep)
869 IF (n_rep > 0) THEN
870 WRITE (unit=iw, fmt="(/,T3,A,/)") &
871 "Position vectors s(i) of the atoms i in scaled coordinates"
872 old_size = 0
873 DO i_rep = 1, n_rep
874 CALL section_vals_val_get(section_vals=section, &
875 keyword_name="POSITION_SCALED", &
876 i_rep_val=i_rep, &
877 i_vals=atomic_indices)
878 n_vals = SIZE(atomic_indices)
879 new_size = old_size + n_vals
880 CALL reallocate(index_list, 1, new_size)
881 index_list(old_size + 1:new_size) = atomic_indices(1:n_vals)
882 old_size = new_size
883 END DO
884 ALLOCATE (work(new_size))
885 CALL sort(index_list, new_size, work)
886 DEALLOCATE (work)
887 DO i = 1, new_size
888 WRITE (unit=string, fmt="(A,I0,A)") "(", index_list(i), ")"
889 IF ((index_list(i) < 1) .OR. (index_list(i) > natom)) THEN
890 WRITE (unit=iw, fmt="(T3,A)") &
891 "Invalid atomic index "//trim(string)//" specified. Print request is ignored."
892 cycle
893 END IF
894 IF (i > 1) THEN
895 ! Skip redundant indices
896 IF (index_list(i) == index_list(i - 1)) cycle
897 END IF
898 r(1:3) = pbc(particle_set(index_list(i))%r(1:3), cell)
899 CALL real_to_scaled(s, r, cell)
900 WRITE (unit=iw, fmt="(T3,A,T20,A,3F13.6)") &
901 "s"//trim(string), "=", s(1:3)
902 END DO
903 DEALLOCATE (index_list)
904 END IF
905
906 ! Print the requested distances
907 CALL section_vals_val_get(section_vals=section, &
908 keyword_name="DISTANCE", &
909 n_rep_val=n)
910 IF (n > 0) THEN
911 WRITE (unit=iw, fmt="(/,T3,A,/)") &
912 "Distance vector r(i,j) between the atom i and j in "// &
913 trim(unit_str)
914 DO i = 1, n
915 CALL section_vals_val_get(section_vals=section, &
916 keyword_name="DISTANCE", &
917 i_rep_val=i, &
918 i_vals=atomic_indices)
919 string = ""
920 WRITE (unit=string, fmt="(A,2(I0,A))") &
921 "(", atomic_indices(1), ",", atomic_indices(2), ")"
922 wrk2 = atomic_indices
923 CALL sort_unique(wrk2, unique)
924 IF (((wrk2(1) >= 1) .AND. (wrk2(SIZE(wrk2)) <= natom)) .AND. unique) THEN
925 rab(:) = pbc(particle_set(atomic_indices(1))%r(:), &
926 particle_set(atomic_indices(2))%r(:), cell)
927 dab = sqrt(rab(1)*rab(1) + rab(2)*rab(2) + rab(3)*rab(3))
928 WRITE (unit=iw, fmt="(T3,A,T20,A,3F13.6,3X,A,F13.6)") &
929 "r"//trim(string), "=", rab(:)*conv, &
930 "|r| =", dab*conv
931 ELSE
932 WRITE (unit=iw, fmt="(T3,A)") &
933 "Invalid atomic indices "//trim(string)//" specified. Print request is ignored."
934 END IF
935 END DO
936 END IF
937
938 ! Print the requested angles
939 CALL section_vals_val_get(section_vals=section, &
940 keyword_name="ANGLE", &
941 n_rep_val=n)
942 IF (n > 0) THEN
943 WRITE (unit=iw, fmt="(/,T3,A,/)") &
944 "Angle a(i,j,k) between the atomic distance vectors r(j,i) and "// &
945 "r(j,k) in DEGREE"
946 DO i = 1, n
947 CALL section_vals_val_get(section_vals=section, &
948 keyword_name="ANGLE", &
949 i_rep_val=i, &
950 i_vals=atomic_indices)
951 string = ""
952 WRITE (unit=string, fmt="(A,3(I0,A))") &
953 "(", atomic_indices(1), ",", atomic_indices(2), ",", atomic_indices(3), ")"
954 wrk3 = atomic_indices
955 CALL sort_unique(wrk3, unique)
956 IF (((wrk3(1) >= 1) .AND. (wrk3(SIZE(wrk3)) <= natom)) .AND. unique) THEN
957 rab(:) = pbc(particle_set(atomic_indices(1))%r(:), &
958 particle_set(atomic_indices(2))%r(:), cell)
959 rbc(:) = pbc(particle_set(atomic_indices(2))%r(:), &
960 particle_set(atomic_indices(3))%r(:), cell)
961 WRITE (unit=iw, fmt="(T3,A,T26,A,F9.3)") &
962 "a"//trim(string), "=", angle(-rab, rbc)*degree
963 ELSE
964 WRITE (unit=iw, fmt="(T3,A)") &
965 "Invalid atomic indices "//trim(string)//" specified. Print request is ignored."
966 END IF
967 END DO
968 END IF
969
970 ! Print the requested dihedral angles
971 CALL section_vals_val_get(section_vals=section, &
972 keyword_name="DIHEDRAL_ANGLE", &
973 n_rep_val=n)
974 IF (n > 0) THEN
975 WRITE (unit=iw, fmt="(/,T3,A,/)") &
976 "Dihedral angle d(i,j,k,l) between the planes (i,j,k) and (j,k,l) "// &
977 "in DEGREE"
978 DO i = 1, n
979 CALL section_vals_val_get(section_vals=section, &
980 keyword_name="DIHEDRAL_ANGLE", &
981 i_rep_val=i, &
982 i_vals=atomic_indices)
983 string = ""
984 WRITE (unit=string, fmt="(A,4(I0,A))") &
985 "(", atomic_indices(1), ",", atomic_indices(2), ",", &
986 atomic_indices(3), ",", atomic_indices(4), ")"
987 wrk4 = atomic_indices
988 CALL sort_unique(wrk4, unique)
989 IF (((wrk4(1) >= 1) .AND. (wrk4(SIZE(wrk4)) <= natom)) .AND. unique) THEN
990 rab(:) = pbc(particle_set(atomic_indices(1))%r(:), &
991 particle_set(atomic_indices(2))%r(:), cell)
992 rbc(:) = pbc(particle_set(atomic_indices(2))%r(:), &
993 particle_set(atomic_indices(3))%r(:), cell)
994 rcd(:) = pbc(particle_set(atomic_indices(3))%r(:), &
995 particle_set(atomic_indices(4))%r(:), cell)
996 WRITE (unit=iw, fmt="(T3,A,T26,A,F9.3)") &
997 "d"//trim(string), "=", dihedral_angle(rab, rbc, rcd)*degree
998 ELSE
999 WRITE (unit=iw, fmt="(T3,A)") &
1000 "Invalid atomic indices "//trim(string)//" specified. Print request is ignored."
1001 END IF
1002 END DO
1003 END IF
1004 END IF
1005 CALL cp_print_key_finished_output(iw, logger, input_section, &
1006 "PRINT%STRUCTURE_DATA")
1007
1008 CALL timestop(handle)
1009
1010 END SUBROUTINE write_structure_data
1011
1012END MODULE particle_methods
static GRID_HOST_DEVICE int modulo(int a, int m)
Equivalent of Fortran's MODULO, which always return a positive number. https://gcc....
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind(atomic_kind, fist_potential, element_symbol, name, mass, kind_number, natom, atom_list, rcov, rvdw, z, qeff, apol, cpol, mm_radius, shell, shell_active, damping)
Get attributes of an atomic kind.
subroutine, public get_gto_basis_set(gto_basis_set, name, aliases, norm_type, kind_radius, ncgf, nset, nsgf, cgf_symbol, sgf_symbol, norm_cgf, set_radius, lmax, lmin, lx, ly, lz, m, ncgf_set, npgf, nsgf_set, nshell, cphi, pgf_radius, sphi, scon, zet, first_cgf, first_sgf, l, last_cgf, last_sgf, n, gcc, maxco, maxl, maxpgf, maxsgf_set, maxshell, maxso, nco_sum, npgf_sum, nshell_sum, maxder, short_kind_radius)
...
Handles all functions related to the CELL.
subroutine, public set_cell_param(cell, cell_length, cell_angle, periodic, do_init_cell)
Sets the cell using the internal parameters (a,b,c) (alpha,beta,gamma) using the convention: a parall...
subroutine, public cell_create(cell, hmat, periodic, tag)
allocates and initializes a cell
Handles all functions related to the CELL.
Definition cell_types.F:15
subroutine, public real_to_scaled(s, r, cell)
Transform real to scaled cell coordinates. s=h_inv*r.
Definition cell_types.F:486
subroutine, public cell_release(cell)
releases the given cell (see doc/ReferenceCounting.html)
Definition cell_types.F:559
subroutine, public cell_clone(cell_in, cell_out, tag)
Clone cell variable.
Definition cell_types.F:107
subroutine, public get_cell(cell, alpha, beta, gamma, deth, orthorhombic, abc, periodic, h, h_inv, symmetry_id, tag)
Get informations about a simulation cell.
Definition cell_types.F:195
various routines to log and control the output. The idea is that decisions about where to log should ...
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer function, public cp_print_key_unit_nr(logger, basis_section, print_key_path, extension, middle_name, local, log_filename, ignore_should_output, file_form, file_position, file_action, file_status, do_backup, on_file, is_new_file, mpi_io, fout)
...
subroutine, public cp_print_key_finished_output(unit_nr, logger, basis_section, print_key_path, local, ignore_should_output, on_file, mpi_io)
should be called after you finish working with a unit obtained with cp_print_key_unit_nr,...
unit conversion facility
Definition cp_units.F:30
real(kind=dp) function, public cp_unit_from_cp2k(value, unit_str, defaults, power)
converts from the internal cp2k units to the given unit
Definition cp_units.F:1179
Definition of the atomic potential types.
Calculation of the incomplete Gamma function F_n(t) for multi-center integrals over Cartesian Gaussia...
Definition gamma.F:15
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public dump_xmol
integer, parameter, public dump_pdb
integer, parameter, public dump_atomic
integer, parameter, public dump_dcd_aligned_cell
integer, parameter, public dump_dcd
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_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
integer, parameter, public sp
Definition kinds.F:33
Definition of mathematical constants and functions.
real(kind=dp), parameter, public degree
Collection of simple mathematical functions and subroutines.
Definition mathlib.F:15
pure real(kind=dp) function, public angle(a, b)
Calculation of the angle between the vectors a and b. The angle is returned in radians.
Definition mathlib.F:175
pure real(kind=dp) function, public dihedral_angle(ab, bc, cd)
Returns the dihedral angle, i.e. the angle between the planes defined by the vectors (-ab,...
Definition mathlib.F:468
Utility routines for the memory handling.
Define methods related to particle_type.
subroutine, public write_qs_particle_coordinates(particle_set, qs_kind_set, subsys_section, label)
Write the atomic coordinates to the output unit.
subroutine, public write_fist_particle_coordinates(particle_set, subsys_section, charges)
Write the atomic coordinates to the output unit.
subroutine, public write_particle_matrix(matrix, particle_set, iw, el_per_part, ilist, parts_per_line)
...
subroutine, public write_structure_data(particle_set, cell, input_section)
Write structure data requested by a separate structure data input section to the output unit....
subroutine, public write_particle_distances(particle_set, cell, subsys_section)
Write the matrix of the particle distances to the output unit.
subroutine, public get_particle_set(particle_set, qs_kind_set, first_sgf, last_sgf, nsgf, nmao, basis)
Get the components of a particle set.
subroutine, public write_particle_coordinates(particle_set, iunit, output_format, content, title, cell, array, unit_conv, charge_occup, charge_beta, charge_extended, print_kind)
Should be able to write a few formats e.g. xmol, and some binary format (dcd) some format can be used...
Define the data structure for the particle information.
pure real(kind=dp) function, dimension(3), public get_particle_pos_or_vel(iatom, particle_set, vector)
Return the atomic position or velocity of atom iatom in x from a packed vector even if core-shell par...
Definition of physical constants:
Definition physcon.F:68
real(kind=dp), parameter, public massunit
Definition physcon.F:141
logical function, public qmmm_ff_precond_only_qm(id1, id2, id3, id4, is_link)
This function handles the atom names and modifies the "_QM_" prefix, in order to find the parameters ...
Define the quickstep kind type and their sub types.
subroutine, public get_qs_kind(qs_kind, basis_set, basis_type, ncgf, nsgf, all_potential, tnadd_potential, gth_potential, sgp_potential, upf_potential, se_parameter, dftb_parameter, xtb_parameter, dftb3_param, zeff, elec_conf, mao, lmax_dftb, alpha_core_charge, ccore_charge, core_charge, core_charge_radius, paw_proj_set, paw_atom, hard_radius, hard0_radius, max_rad_local, covalent_radius, vdw_radius, gpw_r3d_rs_type_forced, harmonics, max_iso_not0, max_s_harm, grid_atom, ngrid_ang, ngrid_rad, lmax_rho0, dft_plus_u_atom, l_of_dft_plus_u, n_of_dft_plus_u, u_minus_j, u_of_dft_plus_u, j_of_dft_plus_u, alpha_of_dft_plus_u, beta_of_dft_plus_u, j0_of_dft_plus_u, occupation_of_dft_plus_u, dispersion, bs_occupation, magnetization, no_optimize, addel, laddel, naddel, orbitals, max_scf, eps_scf, smear, u_ramping, u_minus_j_target, eps_u_ramping, init_u_ramping_each_scf, reltmat, ghost, floating, name, element_symbol, pao_basis_size, pao_potentials, pao_descriptors, nelec)
Get attributes of an atomic kind.
elemental subroutine, public get_shell(shell, charge, charge_core, charge_shell, mass_core, mass_shell, k2_spring, k4_spring, max_dist, shell_cutoff)
...
Utilities for string manipulations.
elemental subroutine, public uppercase(string)
Convert all lower case characters in a string to upper case.
All kind of helpful little routines.
Definition util.F:14
Type defining parameters related to the simulation cell.
Definition cell_types.F:55
type of a logger, at the moment it contains just a print level starting at which level it should be l...
Provides all information about a quickstep kind.