(git:8dea62c)
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-2026 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,&
26 USE cp2k_info, ONLY: compile_revision,&
28 r_cwd,&
41 USE input_constants, ONLY: dump_atomic,&
42 dump_dcd,&
45 dump_pdb,&
58 USE kinds, ONLY: default_path_length,&
60 dp,&
61 sp
62 USE machine, ONLY: m_timestamp,&
64 USE mathconstants, ONLY: degree
65 USE mathlib, ONLY: angle,&
67 gcd
71 USE periodic_table, ONLY: nelem
72 USE physcon, ONLY: massunit
74 USE qs_kind_types, ONLY: get_qs_kind,&
79 USE util, ONLY: sort,&
81#include "./base/base_uses.f90"
82
83 IMPLICIT NONE
84
85 PRIVATE
86
87 ! Public subroutines
88
97
98 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'particle_methods'
99
100CONTAINS
101
102! **************************************************************************************************
103!> \brief Get the components of a particle set.
104!> \param particle_set ...
105!> \param qs_kind_set ...
106!> \param first_sgf ...
107!> \param last_sgf ...
108!> \param nsgf ...
109!> \param nmao ...
110!> \param basis ...
111!> \param ncgf ...
112!> \date 14.01.2002
113!> \par History
114!> - particle type cleaned (13.10.2003,MK)
115!> - refactoring and add basis set option (17.08.2010,jhu)
116!> \author MK
117!> \version 1.0
118! **************************************************************************************************
119 SUBROUTINE get_particle_set(particle_set, qs_kind_set, first_sgf, last_sgf, nsgf, &
120 nmao, basis, ncgf)
121
122 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
123 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
124 INTEGER, DIMENSION(:), INTENT(INOUT), OPTIONAL :: first_sgf, last_sgf, nsgf, nmao
125 TYPE(gto_basis_set_p_type), DIMENSION(:), OPTIONAL :: basis
126 INTEGER, DIMENSION(:), INTENT(INOUT), OPTIONAL :: ncgf
127
128 INTEGER :: ikind, iparticle, isgf, nparticle, ns
129
130 cpassert(ASSOCIATED(particle_set))
131
132 nparticle = SIZE(particle_set)
133 IF (PRESENT(first_sgf)) THEN
134 cpassert(SIZE(first_sgf) >= nparticle)
135 END IF
136 IF (PRESENT(last_sgf)) THEN
137 cpassert(SIZE(last_sgf) >= nparticle)
138 END IF
139 IF (PRESENT(nsgf)) THEN
140 cpassert(SIZE(nsgf) >= nparticle)
141 END IF
142 IF (PRESENT(nmao)) THEN
143 cpassert(SIZE(nmao) >= nparticle)
144 END IF
145 IF (PRESENT(ncgf)) THEN
146 cpassert(SIZE(ncgf) >= nparticle)
147 END IF
148
149 IF (PRESENT(first_sgf) .OR. PRESENT(last_sgf) .OR. PRESENT(nsgf)) THEN
150 isgf = 0
151 DO iparticle = 1, nparticle
152 CALL get_atomic_kind(particle_set(iparticle)%atomic_kind, kind_number=ikind)
153 IF (PRESENT(basis)) THEN
154 IF (ASSOCIATED(basis(ikind)%gto_basis_set)) THEN
155 CALL get_gto_basis_set(gto_basis_set=basis(ikind)%gto_basis_set, nsgf=ns)
156 ELSE
157 ns = 0
158 END IF
159 ELSE
160 CALL get_qs_kind(qs_kind_set(ikind), nsgf=ns)
161 END IF
162 IF (PRESENT(nsgf)) nsgf(iparticle) = ns
163 IF (PRESENT(first_sgf)) first_sgf(iparticle) = isgf + 1
164 isgf = isgf + ns
165 IF (PRESENT(last_sgf)) last_sgf(iparticle) = isgf
166 END DO
167 END IF
168
169 IF (PRESENT(ncgf)) THEN
170 DO iparticle = 1, nparticle
171 CALL get_atomic_kind(particle_set(iparticle)%atomic_kind, kind_number=ikind)
172 IF (PRESENT(basis)) THEN
173 IF (ASSOCIATED(basis(ikind)%gto_basis_set)) THEN
174 CALL get_gto_basis_set(gto_basis_set=basis(ikind)%gto_basis_set, ncgf=ns)
175 ELSE
176 ns = 0
177 END IF
178 ELSE
179 CALL get_qs_kind(qs_kind_set(ikind), ncgf=ns)
180 END IF
181 ncgf(iparticle) = ns
182 END DO
183 END IF
184
185 IF (PRESENT(first_sgf)) THEN
186 IF (SIZE(first_sgf) > nparticle) first_sgf(nparticle + 1) = isgf + 1
187 END IF
188
189 IF (PRESENT(nmao)) THEN
190 DO iparticle = 1, nparticle
191 CALL get_atomic_kind(particle_set(iparticle)%atomic_kind, kind_number=ikind)
192 CALL get_qs_kind(qs_kind_set(ikind), mao=ns)
193 nmao(iparticle) = ns
194 END DO
195 END IF
196
197 END SUBROUTINE get_particle_set
198
199! **************************************************************************************************
200!> \brief Should be able to write a few formats e.g. xmol, and some binary
201!> format (dcd) some format can be used for x,v,f
202!>
203!> FORMAT CONTENT UNITS x, v, f
204!> XMOL POS, VEL, FORCE, POS_VEL, POS_VEL_FORCE Angstrom, a.u., a.u.
205!>
206!> \param particle_set ...
207!> \param iunit ...
208!> \param output_format ...
209!> \param content ...
210!> \param title ...
211!> \param cell ...
212!> \param array ...
213!> \param unit_conv ...
214!> \param charge_occup ...
215!> \param charge_beta ...
216!> \param charge_extended ...
217!> \param print_kind ...
218!> \date 14.01.2002
219!> \author MK
220!> \version 1.0
221! **************************************************************************************************
222 SUBROUTINE write_particle_coordinates(particle_set, iunit, output_format, &
223 content, title, cell, array, unit_conv, &
224 charge_occup, charge_beta, &
225 charge_extended, print_kind)
226
227 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
228 INTEGER :: iunit, output_format
229 CHARACTER(LEN=*) :: content, title
230 TYPE(cell_type), OPTIONAL, POINTER :: cell
231 REAL(kind=dp), DIMENSION(:), INTENT(IN), OPTIONAL :: array
232 REAL(kind=dp), INTENT(IN), OPTIONAL :: unit_conv
233 LOGICAL, INTENT(IN), OPTIONAL :: charge_occup, charge_beta, &
234 charge_extended, print_kind
235
236 CHARACTER(len=*), PARAMETER :: routinen = 'write_particle_coordinates'
237
238 CHARACTER(LEN=120) :: line
239 CHARACTER(LEN=2) :: element_symbol
240 CHARACTER(LEN=4) :: name
241 CHARACTER(LEN=default_string_length) :: atm_name, my_format
242 INTEGER :: handle, iatom, natom
243 LOGICAL :: dummy, my_charge_beta, &
244 my_charge_extended, my_charge_occup, &
245 my_print_kind
246 REAL(kind=dp) :: angle_alpha, angle_beta, angle_gamma, &
247 factor, qeff
248 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: arr
249 REAL(kind=dp), DIMENSION(3) :: abc, angles, f, r, v
250 REAL(kind=dp), DIMENSION(3, 3) :: h
251 REAL(kind=sp), ALLOCATABLE, DIMENSION(:) :: x4, y4, z4
252 TYPE(cell_type), POINTER :: cell_dcd
253 TYPE(fist_potential_type), POINTER :: fist_potential
254 TYPE(shell_kind_type), POINTER :: shell
255
256 CALL timeset(routinen, handle)
257
258 natom = SIZE(particle_set)
259 IF (PRESENT(array)) THEN
260 SELECT CASE (trim(content))
261 CASE ("POS_VEL", "POS_VEL_FORCE")
262 cpabort("Illegal usage")
263 END SELECT
264 END IF
265 factor = 1.0_dp
266 IF (PRESENT(unit_conv)) THEN
267 factor = unit_conv
268 END IF
269 SELECT CASE (output_format)
270 CASE (dump_xmol, dump_extxyz)
271 my_print_kind = .false.
272 IF (PRESENT(print_kind)) my_print_kind = print_kind
273 WRITE (iunit, "(I8)") natom
274 WRITE (iunit, "(A)") trim(title)
275 DO iatom = 1, natom
276 CALL get_atomic_kind(atomic_kind=particle_set(iatom)%atomic_kind, &
277 element_symbol=element_symbol)
278 IF (len_trim(element_symbol) == 0 .OR. my_print_kind) THEN
279 CALL get_atomic_kind(atomic_kind=particle_set(iatom)%atomic_kind, &
280 name=atm_name)
281 dummy = qmmm_ff_precond_only_qm(id1=atm_name)
282 my_format = "(A,"
283 atm_name = trim(atm_name)
284 ELSE
285 my_format = "(T2,A2,"
286 atm_name = trim(element_symbol)
287 END IF
288 SELECT CASE (trim(content))
289 CASE ("POS")
290 IF (PRESENT(array)) THEN
291 r(1:3) = get_particle_pos_or_vel(iatom, particle_set, array)
292 ELSE
293 r(:) = particle_set(iatom)%r(:)
294 END IF
295 WRITE (iunit, trim(my_format)//"1X,3F20.10)") trim(atm_name), r(1:3)*factor
296 CASE ("VEL")
297 IF (PRESENT(array)) THEN
298 v(1:3) = get_particle_pos_or_vel(iatom, particle_set, array)
299 ELSE
300 v(:) = particle_set(iatom)%v(:)
301 END IF
302 WRITE (iunit, trim(my_format)//"1X,3F20.10)") trim(atm_name), v(1:3)*factor
303 CASE ("FORCE")
304 IF (PRESENT(array)) THEN
305 f(:) = array((iatom - 1)*3 + 1:(iatom - 1)*3 + 3)
306 ELSE
307 f(:) = particle_set(iatom)%f(:)
308 END IF
309 WRITE (iunit, trim(my_format)//"1X,3F20.10)") trim(atm_name), f(1:3)*factor
310 CASE ("FORCE_MIXING_LABELS")
311 IF (PRESENT(array)) THEN
312 f(:) = array((iatom - 1)*3 + 1:(iatom - 1)*3 + 3)
313 ELSE
314 f(:) = particle_set(iatom)%f(:)
315 END IF
316 WRITE (iunit, trim(my_format)//"1X,3F20.10)") trim(atm_name), f(1:3)*factor
317 END SELECT
318 END DO
319 CASE (dump_atomic)
320 DO iatom = 1, natom
321 SELECT CASE (trim(content))
322 CASE ("POS")
323 IF (PRESENT(array)) THEN
324 r(1:3) = get_particle_pos_or_vel(iatom, particle_set, array)
325 ELSE
326 r(:) = particle_set(iatom)%r(:)
327 END IF
328 WRITE (iunit, "(3F20.10)") r(1:3)*factor
329 CASE ("VEL")
330 IF (PRESENT(array)) THEN
331 v(1:3) = get_particle_pos_or_vel(iatom, particle_set, array)
332 ELSE
333 v(:) = particle_set(iatom)%v(:)
334 END IF
335 WRITE (iunit, "(3F20.10)") v(1:3)*factor
336 CASE ("FORCE")
337 IF (PRESENT(array)) THEN
338 f(:) = array((iatom - 1)*3 + 1:(iatom - 1)*3 + 3)
339 ELSE
340 f(:) = particle_set(iatom)%f(:)
341 END IF
342 WRITE (iunit, "(3F20.10)") f(1:3)*factor
343 CASE ("FORCE_MIXING_LABELS")
344 IF (PRESENT(array)) THEN
345 f(:) = array((iatom - 1)*3 + 1:(iatom - 1)*3 + 3)
346 ELSE
347 f(:) = particle_set(iatom)%f(:)
348 END IF
349 WRITE (iunit, "(3F20.10)") f(1:3)*factor
350 END SELECT
351 END DO
353 IF (.NOT. (PRESENT(cell))) THEN
354 cpabort("Cell is not present! Report this bug!")
355 END IF
356 CALL get_cell(cell, alpha=angle_alpha, beta=angle_beta, gamma=angle_gamma, &
357 abc=abc)
358 IF (.NOT. cell%orthorhombic .AND. (output_format == dump_dcd_aligned_cell)) THEN
359 ! In the case of a non-orthorhombic cell adopt a common convention
360 ! for the orientation of the cell with respect to the Cartesian axes:
361 ! Cell vector a is aligned with the x axis and the cell vector b lies
362 ! in the xy plane.
363 NULLIFY (cell_dcd)
364 CALL cell_create(cell_dcd)
365 CALL cell_clone(cell, cell_dcd, tag="CELL_DCD")
366 angles(1) = angle_alpha/degree
367 angles(2) = angle_beta/degree
368 angles(3) = angle_gamma/degree
369 CALL set_cell_param(cell_dcd, abc, angles, &
370 do_init_cell=.true.)
371 h(1:3, 1:3) = matmul(cell_dcd%hmat(1:3, 1:3), cell%h_inv(1:3, 1:3))
372 CALL cell_release(cell_dcd)
373 END IF
374 ALLOCATE (arr(3, natom))
375 IF (PRESENT(array)) THEN
376 arr(1:3, 1:natom) = reshape(array, [3, natom])
377 ELSE
378 SELECT CASE (trim(content))
379 CASE ("POS")
380 DO iatom = 1, natom
381 arr(1:3, iatom) = particle_set(iatom)%r(1:3)
382 END DO
383 CASE ("VEL")
384 DO iatom = 1, natom
385 arr(1:3, iatom) = particle_set(iatom)%v(1:3)
386 END DO
387 CASE ("FORCE")
388 DO iatom = 1, natom
389 arr(1:3, iatom) = particle_set(iatom)%f(1:3)
390 END DO
391 CASE DEFAULT
392 cpabort("Illegal DCD dump type")
393 END SELECT
394 END IF
395 ALLOCATE (x4(natom))
396 ALLOCATE (y4(natom))
397 ALLOCATE (z4(natom))
398 IF (.NOT. cell%orthorhombic .AND. (output_format == dump_dcd_aligned_cell)) THEN
399 x4(1:natom) = real(matmul(h(1, 1:3), arr(1:3, 1:natom)), kind=sp)
400 y4(1:natom) = real(matmul(h(2, 1:3), arr(1:3, 1:natom)), kind=sp)
401 z4(1:natom) = real(matmul(h(3, 1:3), arr(1:3, 1:natom)), kind=sp)
402 ELSE
403 x4(1:natom) = real(arr(1, 1:natom), kind=sp)
404 y4(1:natom) = real(arr(2, 1:natom), kind=sp)
405 z4(1:natom) = real(arr(3, 1:natom), kind=sp)
406 END IF
407 WRITE (iunit) abc(1)*factor, angle_gamma, abc(2)*factor, &
408 angle_beta, angle_alpha, abc(3)*factor
409 WRITE (iunit) x4*real(factor, kind=sp)
410 WRITE (iunit) y4*real(factor, kind=sp)
411 WRITE (iunit) z4*real(factor, kind=sp)
412 ! Release work storage
413 DEALLOCATE (arr)
414 DEALLOCATE (x4)
415 DEALLOCATE (y4)
416 DEALLOCATE (z4)
417 CASE (dump_pdb)
418 my_charge_occup = .false.
419 IF (PRESENT(charge_occup)) my_charge_occup = charge_occup
420 my_charge_beta = .false.
421 IF (PRESENT(charge_beta)) my_charge_beta = charge_beta
422 my_charge_extended = .false.
423 IF (PRESENT(charge_extended)) my_charge_extended = charge_extended
424 IF (len_trim(title) > 0) THEN
425 WRITE (unit=iunit, fmt="(A6,T11,A)") &
426 "REMARK", trim(title)
427 END IF
428 CALL get_cell(cell, alpha=angle_alpha, beta=angle_beta, gamma=angle_gamma, abc=abc)
429 ! COLUMNS DATA TYPE CONTENTS
430 ! --------------------------------------------------
431 ! 1 - 6 Record name "CRYST1"
432 ! 7 - 15 Real(9.3) a (Angstroms)
433 ! 16 - 24 Real(9.3) b (Angstroms)
434 ! 25 - 33 Real(9.3) c (Angstroms)
435 ! 34 - 40 Real(7.2) alpha (degrees)
436 ! 41 - 47 Real(7.2) beta (degrees)
437 ! 48 - 54 Real(7.2) gamma (degrees)
438 ! 56 - 66 LString Space group
439 ! 67 - 70 Integer Z value
440 WRITE (unit=iunit, fmt="(A6,3F9.3,3F7.2)") &
441 "CRYST1", abc(1:3)*factor, angle_alpha, angle_beta, angle_gamma
442 WRITE (unit=line(1:6), fmt="(A6)") "ATOM "
443 DO iatom = 1, natom
444 line = ""
445 ! COLUMNS DATA TYPE CONTENTS
446 ! 1 - 6 Record name "ATOM "
447 ! 7 - 11 Integer Atom serial number
448 ! 13 - 16 Atom Atom name
449 ! 17 Character Alternate location indicator
450 ! 18 - 20 Residue name Residue name
451 ! 22 Character Chain identifier
452 ! 23 - 26 Integer Residue sequence number
453 ! 27 AChar Code for insertion of residues
454 ! 31 - 38 Real(8.3) Orthogonal coordinates for X in Angstrom
455 ! 39 - 46 Real(8.3) Orthogonal coordinates for Y in Angstrom
456 ! 47 - 54 Real(8.3) Orthogonal coordinates for Z in Angstrom
457 ! 55 - 60 Real(6.2) Occupancy
458 ! 61 - 66 Real(6.2) Temperature factor (Default = 0.0)
459 ! 73 - 76 LString(4) Segment identifier, left-justified
460 ! 77 - 78 LString(2) Element symbol, right-justified
461 ! 79 - 80 LString(2) Charge on the atom
462 CALL get_atomic_kind(atomic_kind=particle_set(iatom)%atomic_kind, &
463 element_symbol=element_symbol, name=atm_name, &
464 fist_potential=fist_potential, shell=shell)
465 IF (len_trim(element_symbol) == 0) THEN
466 dummy = qmmm_ff_precond_only_qm(id1=atm_name)
467 END IF
468 name = trim(atm_name)
469 IF (ASSOCIATED(fist_potential)) THEN
470 CALL get_potential(potential=fist_potential, qeff=qeff)
471 ELSE
472 qeff = 0.0_dp
473 END IF
474 IF (ASSOCIATED(shell)) CALL get_shell(shell=shell, charge=qeff)
475 WRITE (unit=line(1:6), fmt="(A6)") "ATOM "
476 WRITE (unit=line(7:11), fmt="(I5)") modulo(iatom, 100000)
477 WRITE (unit=line(13:16), fmt="(A4)") adjustl(name)
478 ! WRITE (UNIT=line(18:20),FMT="(A3)") TRIM(resname)
479 ! WRITE (UNIT=line(23:26),FMT="(I4)") MODULO(idres,10000)
480 SELECT CASE (trim(content))
481 CASE ("POS")
482 IF (PRESENT(array)) THEN
483 r(1:3) = get_particle_pos_or_vel(iatom, particle_set, array)
484 ELSE
485 r(:) = particle_set(iatom)%r(:)
486 END IF
487 WRITE (unit=line(31:54), fmt="(3F8.3)") r(1:3)*factor
488 CASE DEFAULT
489 cpabort("PDB dump only for trajectory available")
490 END SELECT
491 IF (my_charge_occup) THEN
492 WRITE (unit=line(55:60), fmt="(F6.2)") qeff
493 ELSE
494 WRITE (unit=line(55:60), fmt="(F6.2)") 0.0_dp
495 END IF
496 IF (my_charge_beta) THEN
497 WRITE (unit=line(61:66), fmt="(F6.2)") qeff
498 ELSE
499 WRITE (unit=line(61:66), fmt="(F6.2)") 0.0_dp
500 END IF
501 ! WRITE (UNIT=line(73:76),FMT="(A4)") ADJUSTL(TRIM(molname))
502 WRITE (unit=line(77:78), fmt="(A2)") adjustr(trim(element_symbol))
503 IF (my_charge_extended) THEN
504 WRITE (unit=line(81:), fmt="(SP,F0.8)") qeff
505 END IF
506 WRITE (unit=iunit, fmt="(A)") trim(line)
507 END DO
508 WRITE (unit=iunit, fmt="(A)") "END"
509 CASE DEFAULT
510 cpabort("Illegal dump type")
511 END SELECT
512
513 CALL timestop(handle)
514
515 END SUBROUTINE write_particle_coordinates
516
517! **************************************************************************************************
518!> \brief Write the atomic coordinates to the output unit.
519!> \param particle_set ...
520!> \param subsys_section ...
521!> \param charges ...
522!> \date 05.06.2000
523!> \author MK
524!> \version 1.0
525! **************************************************************************************************
526 SUBROUTINE write_fist_particle_coordinates(particle_set, subsys_section, charges)
527 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
528 TYPE(section_vals_type), POINTER :: subsys_section
529 REAL(kind=dp), DIMENSION(:), OPTIONAL :: charges
530
531 CHARACTER(LEN=default_string_length) :: name, unit_str
532 INTEGER :: iatom, ikind, iw, natom
533 REAL(kind=dp) :: conv, mass, qcore, qeff, qshell
534 TYPE(cp_logger_type), POINTER :: logger
535 TYPE(shell_kind_type), POINTER :: shell_kind
536
537 NULLIFY (logger)
538 NULLIFY (shell_kind)
539
540 logger => cp_get_default_logger()
541 iw = cp_print_key_unit_nr(logger, subsys_section, &
542 "PRINT%ATOMIC_COORDINATES", extension=".coordLog")
543
544 CALL section_vals_val_get(subsys_section, "PRINT%ATOMIC_COORDINATES%UNIT", c_val=unit_str)
545 conv = cp_unit_from_cp2k(1.0_dp, trim(unit_str))
546 CALL uppercase(unit_str)
547 IF (iw > 0) THEN
548 WRITE (unit=iw, fmt="(/,/,T2,A)") &
549 "MODULE FIST: ATOMIC COORDINATES IN "//trim(unit_str)
550 WRITE (unit=iw, fmt="(/,T4,A,T30,A,T44,A,T58,A,T66,A,T77,A)") &
551 "Atom Kind Name", "X", "Y", "Z", "q(eff)", "Mass"
552 natom = SIZE(particle_set)
553 DO iatom = 1, natom
554 CALL get_atomic_kind(atomic_kind=particle_set(iatom)%atomic_kind, &
555 kind_number=ikind, &
556 name=name, &
557 mass=mass, &
558 qeff=qeff, &
559 shell=shell_kind)
560 IF (PRESENT(charges)) qeff = charges(iatom)
561 IF (ASSOCIATED(shell_kind)) THEN
562 CALL get_shell(shell=shell_kind, &
563 charge_core=qcore, &
564 charge_shell=qshell)
565 qeff = qcore + qshell
566 END IF
567 WRITE (unit=iw, fmt="(T2,I6,1X,I4,1X,A7,3(1X,F13.6),2(1X,F8.4))") &
568 iatom, ikind, name, particle_set(iatom)%r(1:3)*conv, qeff, mass/massunit
569 END DO
570 WRITE (iw, "(A)") ""
571 END IF
572
573 CALL cp_print_key_finished_output(iw, logger, subsys_section, &
574 "PRINT%ATOMIC_COORDINATES")
575
577
578! **************************************************************************************************
579!> \brief Write the atomic coordinates to the output unit.
580!> \param particle_set ...
581!> \param qs_kind_set ...
582!> \param subsys_section ...
583!> \param label ...
584!> \date 05.06.2000
585!> \author MK
586!> \version 1.0
587! **************************************************************************************************
588 SUBROUTINE write_qs_particle_coordinates(particle_set, qs_kind_set, subsys_section, label)
589
590 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
591 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
592 TYPE(section_vals_type), POINTER :: subsys_section
593 CHARACTER(LEN=*), INTENT(IN) :: label
594
595 CHARACTER(len=*), PARAMETER :: routinen = 'write_qs_particle_coordinates'
596
597 CHARACTER(LEN=2) :: element_symbol
598 CHARACTER(LEN=default_string_length) :: unit_str
599 INTEGER :: handle, iatom, ikind, iw, natom, z
600 REAL(kind=dp) :: conv, mass, zeff
601 TYPE(cp_logger_type), POINTER :: logger
602
603 CALL timeset(routinen, handle)
604
605 NULLIFY (logger)
606 logger => cp_get_default_logger()
607 iw = cp_print_key_unit_nr(logger, subsys_section, &
608 "PRINT%ATOMIC_COORDINATES", extension=".coordLog")
609
610 CALL section_vals_val_get(subsys_section, "PRINT%ATOMIC_COORDINATES%UNIT", c_val=unit_str)
611 conv = cp_unit_from_cp2k(1.0_dp, trim(unit_str))
612 CALL uppercase(unit_str)
613 IF (iw > 0) THEN
614 WRITE (unit=iw, fmt="(/,/,T2,A)") &
615 "MODULE "//trim(label)//": ATOMIC COORDINATES IN "//trim(unit_str)
616 WRITE (unit=iw, fmt="(/,T4,A,T30,A,T44,A,T58,A,T66,A,T77,A)") &
617 "Atom Kind Element", "X", "Y", "Z", "Z(eff)", "Mass"
618 natom = SIZE(particle_set)
619 DO iatom = 1, natom
620 CALL get_atomic_kind(atomic_kind=particle_set(iatom)%atomic_kind, &
621 kind_number=ikind, &
622 element_symbol=element_symbol, &
623 mass=mass, &
624 z=z)
625 CALL get_qs_kind(qs_kind_set(ikind), zeff=zeff)
626 WRITE (unit=iw, fmt="(T2,I6,1X,I4,1X,A2,1X,I4,3(1X,F13.6),2(1X,F8.4))") &
627 iatom, ikind, element_symbol, z, particle_set(iatom)%r(1:3)*conv, zeff, mass/massunit
628 END DO
629 WRITE (iw, "(A)") ""
630 END IF
631
632 CALL cp_print_key_finished_output(iw, logger, subsys_section, &
633 "PRINT%ATOMIC_COORDINATES")
634
635 CALL timestop(handle)
636
637 END SUBROUTINE write_qs_particle_coordinates
638
639! **************************************************************************************************
640!> \brief Write the matrix of the particle distances to the output unit.
641!> \param particle_set ...
642!> \param cell ...
643!> \param subsys_section ...
644!> \date 06.10.2000
645!> \author Matthias Krack
646!> \version 1.0
647! **************************************************************************************************
648 SUBROUTINE write_particle_distances(particle_set, cell, subsys_section)
649
650 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
651 TYPE(cell_type), POINTER :: cell
652 TYPE(section_vals_type), POINTER :: subsys_section
653
654 CHARACTER(len=*), PARAMETER :: routinen = 'write_particle_distances'
655
656 CHARACTER(LEN=default_string_length) :: unit_str
657 INTEGER :: handle, iatom, iw, jatom, natom
658 INTEGER, DIMENSION(3) :: periodic
659 LOGICAL :: explicit
660 REAL(kind=dp) :: conv, dab, dab_abort, dab_min, dab_warn
661 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: distance_matrix
662 REAL(kind=dp), DIMENSION(3) :: rab
663 TYPE(cp_logger_type), POINTER :: logger
664
665 CALL timeset(routinen, handle)
666
667 cpassert(ASSOCIATED(particle_set))
668 cpassert(ASSOCIATED(cell))
669 cpassert(ASSOCIATED(subsys_section))
670
671 NULLIFY (logger)
672 logger => cp_get_default_logger()
673 iw = cp_print_key_unit_nr(logger, subsys_section, &
674 "PRINT%INTERATOMIC_DISTANCES", extension=".distLog")
675
676 CALL section_vals_val_get(subsys_section, "PRINT%INTERATOMIC_DISTANCES%UNIT", c_val=unit_str)
677 conv = cp_unit_from_cp2k(1.0_dp, trim(unit_str))
678 CALL section_vals_val_get(subsys_section, "PRINT%INTERATOMIC_DISTANCES%CHECK_INTERATOMIC_DISTANCES", &
679 r_val=dab_min, explicit=explicit)
680
681 dab_abort = 0.0_dp
682 dab_warn = 0.0_dp
683 natom = SIZE(particle_set)
684
685 ! Compute interatomic distances only if their printout or check is explicitly requested
686 ! Disable the default check for systems with more than 3000 atoms
687 IF (explicit .OR. (iw > 0) .OR. (natom <= 2000)) THEN
688 IF (dab_min > 0.0_dp) THEN
689 dab_warn = dab_min*conv
690 ELSE IF (dab_min < 0.0_dp) THEN
691 dab_abort = abs(dab_min)*conv
692 END IF
693 END IF
694
695 IF ((iw > 0) .OR. (dab_abort > 0.0_dp) .OR. (dab_warn > 0.0_dp)) THEN
696 CALL get_cell(cell=cell, periodic=periodic)
697 IF (iw > 0) THEN
698 ALLOCATE (distance_matrix(natom, natom))
699 distance_matrix(:, :) = 0.0_dp
700 END IF
701 DO iatom = 1, natom
702 DO jatom = iatom + 1, natom
703 rab(:) = pbc(particle_set(iatom)%r(:), &
704 particle_set(jatom)%r(:), cell)
705 dab = sqrt(rab(1)*rab(1) + rab(2)*rab(2) + rab(3)*rab(3))*conv
706 IF (dab_abort > 0.0_dp) THEN
707 ! Stop the run for interatomic distances smaller than the requested threshold
708 IF (dab < dab_abort) THEN
709 CALL cp_abort(__location__, "The distance between the atoms "// &
710 trim(adjustl(cp_to_string(iatom, fmt="(I8)")))//" and "// &
711 trim(adjustl(cp_to_string(jatom, fmt="(I8)")))//" is only "// &
712 trim(adjustl(cp_to_string(dab, fmt="(F6.3)")))//" "// &
713 trim(adjustl(unit_str))//" and thus smaller than the requested threshold of "// &
714 trim(adjustl(cp_to_string(dab_abort, fmt="(F6.3)")))//" "// &
715 trim(adjustl(unit_str)))
716 END IF
717 END IF
718 IF (dab < dab_warn) THEN
719 ! Print warning for interatomic distances smaller than the requested threshold
720 CALL cp_warn(__location__, "The distance between the atoms "// &
721 trim(adjustl(cp_to_string(iatom, fmt="(I8)")))//" and "// &
722 trim(adjustl(cp_to_string(jatom, fmt="(I8)")))//" is only "// &
723 trim(adjustl(cp_to_string(dab, fmt="(F6.3)")))//" "// &
724 trim(adjustl(unit_str))//" and thus smaller than the threshold of "// &
725 trim(adjustl(cp_to_string(dab_warn, fmt="(F6.3)")))//" "// &
726 trim(adjustl(unit_str)))
727 END IF
728 IF (iw > 0) THEN
729 distance_matrix(iatom, jatom) = dab
730 distance_matrix(jatom, iatom) = distance_matrix(iatom, jatom)
731 END IF
732 END DO
733 END DO
734 IF (iw > 0) THEN
735 ! Print the distance matrix
736 WRITE (unit=iw, fmt="(/,/,T2,A)") &
737 "INTERATOMIC DISTANCES IN "//trim(unit_str)
738 CALL write_particle_matrix(distance_matrix, particle_set, iw)
739 IF (ALLOCATED(distance_matrix)) DEALLOCATE (distance_matrix)
740 END IF
741 CALL cp_print_key_finished_output(iw, logger, subsys_section, &
742 "PRINT%INTERATOMIC_DISTANCES")
743 END IF
744
745 CALL timestop(handle)
746
747 END SUBROUTINE write_particle_distances
748
749! **************************************************************************************************
750!> \brief ...
751!> \param matrix ...
752!> \param particle_set ...
753!> \param iw ...
754!> \param el_per_part ...
755!> \param Ilist ...
756!> \param parts_per_line : number of particle columns to be printed in one line
757! **************************************************************************************************
758 SUBROUTINE write_particle_matrix(matrix, particle_set, iw, el_per_part, Ilist, parts_per_line)
759 REAL(kind=dp), DIMENSION(:, :) :: matrix
760 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
761 INTEGER, INTENT(IN) :: iw
762 INTEGER, INTENT(IN), OPTIONAL :: el_per_part
763 INTEGER, DIMENSION(:), OPTIONAL, POINTER :: ilist
764 INTEGER, INTENT(IN), OPTIONAL :: parts_per_line
765
766 CHARACTER(LEN=2) :: element_symbol
767 CHARACTER(LEN=default_string_length) :: fmt_string1, fmt_string2
768 INTEGER :: from, i, iatom, icol, jatom, katom, &
769 my_el_per_part, my_parts_per_line, &
770 natom, to
771 INTEGER, DIMENSION(:), POINTER :: my_list
772
773 my_el_per_part = 1
774 IF (PRESENT(el_per_part)) my_el_per_part = el_per_part
775 my_parts_per_line = 5
776 IF (PRESENT(parts_per_line)) my_parts_per_line = max(parts_per_line, 1)
777 WRITE (fmt_string1, fmt='(A,I0,A)') &
778 "(/,T2,11X,", my_parts_per_line, "(4X,I5,4X))"
779 WRITE (fmt_string2, fmt='(A,I0,A)') &
780 "(T2,I5,2X,A2,2X,", my_parts_per_line, "(1X,F12.6))"
781 IF (PRESENT(ilist)) THEN
782 natom = SIZE(ilist)
783 ELSE
784 natom = SIZE(particle_set)
785 END IF
786 ALLOCATE (my_list(natom))
787 IF (PRESENT(ilist)) THEN
788 my_list = ilist
789 ELSE
790 DO i = 1, natom
791 my_list(i) = i
792 END DO
793 END IF
794 natom = natom*my_el_per_part
795 DO jatom = 1, natom, my_parts_per_line
796 from = jatom
797 to = min(from + my_parts_per_line - 1, natom)
798 WRITE (unit=iw, fmt=trim(fmt_string1)) (icol, icol=from, to)
799 DO iatom = 1, natom
800 katom = iatom/my_el_per_part
801 IF (mod(iatom, my_el_per_part) /= 0) katom = katom + 1
802 CALL get_atomic_kind(atomic_kind=particle_set(my_list(katom))%atomic_kind, &
803 element_symbol=element_symbol)
804 WRITE (unit=iw, fmt=trim(fmt_string2)) &
805 iatom, element_symbol, &
806 (matrix(iatom, icol), icol=from, to)
807 END DO
808 END DO
809
810 DEALLOCATE (my_list)
811
812 END SUBROUTINE write_particle_matrix
813
814! **************************************************************************************************
815!> \brief Write structure data requested by a separate structure data input
816!> section to the output unit.
817!> input_section can be either motion_section or subsys_section.
818!>
819!> \param particle_set ...
820!> \param cell ...
821!> \param input_section ...
822!> \date 11.03.04
823!> \par History
824!> Recovered (23.03.06,MK)
825!> \author MK
826!> \version 1.0
827! **************************************************************************************************
828 SUBROUTINE write_structure_data(particle_set, cell, input_section)
829 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
830 TYPE(cell_type), POINTER :: cell
831 TYPE(section_vals_type), POINTER :: input_section
832
833 CHARACTER(LEN=*), PARAMETER :: routinen = 'write_structure_data'
834
835 CHARACTER(LEN=default_string_length) :: string, unit_str
836 INTEGER :: handle, i, i_rep, iw, n, n_rep, n_vals, &
837 natom, new_size, old_size, wrk2(2), &
838 wrk3(3), wrk4(4)
839 INTEGER, ALLOCATABLE, DIMENSION(:) :: work
840 INTEGER, DIMENSION(:), POINTER :: atomic_indices, index_list
841 LOGICAL :: unique
842 REAL(kind=dp) :: conv, dab
843 REAL(kind=dp), DIMENSION(3) :: r, rab, rbc, rcd, s
844 TYPE(cp_logger_type), POINTER :: logger
845 TYPE(section_vals_type), POINTER :: section
846
847 CALL timeset(routinen, handle)
848 NULLIFY (atomic_indices)
849 NULLIFY (index_list)
850 NULLIFY (logger)
851 NULLIFY (section)
852 string = ""
853
854 logger => cp_get_default_logger()
855 iw = cp_print_key_unit_nr(logger=logger, &
856 basis_section=input_section, &
857 print_key_path="PRINT%STRUCTURE_DATA", &
858 extension=".coordLog")
859
860 CALL section_vals_val_get(input_section, "PRINT%STRUCTURE_DATA%UNIT", c_val=unit_str)
861 conv = cp_unit_from_cp2k(1.0_dp, trim(unit_str))
862 CALL uppercase(unit_str)
863 IF (iw > 0) THEN
864 natom = SIZE(particle_set)
865 section => section_vals_get_subs_vals(section_vals=input_section, &
866 subsection_name="PRINT%STRUCTURE_DATA")
867
868 WRITE (unit=iw, fmt="(/,T2,A)") "REQUESTED STRUCTURE DATA"
869 ! Print the requested atomic position vectors
870 CALL section_vals_val_get(section_vals=section, &
871 keyword_name="POSITION", &
872 n_rep_val=n_rep)
873 IF (n_rep > 0) THEN
874 WRITE (unit=iw, fmt="(/,T3,A,/)") &
875 "Position vectors r(i) of the atoms i in "//trim(unit_str)
876 old_size = 0
877 DO i_rep = 1, n_rep
878 CALL section_vals_val_get(section_vals=section, &
879 keyword_name="POSITION", &
880 i_rep_val=i_rep, &
881 i_vals=atomic_indices)
882 n_vals = SIZE(atomic_indices)
883 new_size = old_size + n_vals
884 CALL reallocate(index_list, 1, new_size)
885 index_list(old_size + 1:new_size) = atomic_indices(1:n_vals)
886 old_size = new_size
887 END DO
888 ALLOCATE (work(new_size))
889 CALL sort(index_list, new_size, work)
890 DEALLOCATE (work)
891 DO i = 1, new_size
892 WRITE (unit=string, fmt="(A,I0,A)") "(", index_list(i), ")"
893 IF ((index_list(i) < 1) .OR. (index_list(i) > natom)) THEN
894 WRITE (unit=iw, fmt="(T3,A)") &
895 "Invalid atomic index "//trim(string)//" specified. Print request is ignored."
896 cycle
897 END IF
898 IF (i > 1) THEN
899 ! Skip redundant indices
900 IF (index_list(i) == index_list(i - 1)) cycle
901 END IF
902 WRITE (unit=iw, fmt="(T3,A,T20,A,3F13.6)") &
903 "r"//trim(string), "=", pbc(particle_set(index_list(i))%r(1:3), cell)*conv
904 END DO
905 DEALLOCATE (index_list)
906 END IF
907
908 ! Print the requested atomic position vectors in scaled coordinates
909 CALL section_vals_val_get(section_vals=section, &
910 keyword_name="POSITION_SCALED", &
911 n_rep_val=n_rep)
912 IF (n_rep > 0) THEN
913 WRITE (unit=iw, fmt="(/,T3,A,/)") &
914 "Position vectors s(i) of the atoms i in scaled coordinates"
915 old_size = 0
916 DO i_rep = 1, n_rep
917 CALL section_vals_val_get(section_vals=section, &
918 keyword_name="POSITION_SCALED", &
919 i_rep_val=i_rep, &
920 i_vals=atomic_indices)
921 n_vals = SIZE(atomic_indices)
922 new_size = old_size + n_vals
923 CALL reallocate(index_list, 1, new_size)
924 index_list(old_size + 1:new_size) = atomic_indices(1:n_vals)
925 old_size = new_size
926 END DO
927 ALLOCATE (work(new_size))
928 CALL sort(index_list, new_size, work)
929 DEALLOCATE (work)
930 DO i = 1, new_size
931 WRITE (unit=string, fmt="(A,I0,A)") "(", index_list(i), ")"
932 IF ((index_list(i) < 1) .OR. (index_list(i) > natom)) THEN
933 WRITE (unit=iw, fmt="(T3,A)") &
934 "Invalid atomic index "//trim(string)//" specified. Print request is ignored."
935 cycle
936 END IF
937 IF (i > 1) THEN
938 ! Skip redundant indices
939 IF (index_list(i) == index_list(i - 1)) cycle
940 END IF
941 r(1:3) = pbc(particle_set(index_list(i))%r(1:3), cell)
942 CALL real_to_scaled(s, r, cell)
943 WRITE (unit=iw, fmt="(T3,A,T20,A,3F13.6)") &
944 "s"//trim(string), "=", s(1:3)
945 END DO
946 DEALLOCATE (index_list)
947 END IF
948
949 ! Print the requested distances
950 CALL section_vals_val_get(section_vals=section, &
951 keyword_name="DISTANCE", &
952 n_rep_val=n)
953 IF (n > 0) THEN
954 WRITE (unit=iw, fmt="(/,T3,A,/)") &
955 "Distance vector r(i,j) between the atom i and j in "// &
956 trim(unit_str)
957 DO i = 1, n
958 CALL section_vals_val_get(section_vals=section, &
959 keyword_name="DISTANCE", &
960 i_rep_val=i, &
961 i_vals=atomic_indices)
962 string = ""
963 WRITE (unit=string, fmt="(A,2(I0,A))") &
964 "(", atomic_indices(1), ",", atomic_indices(2), ")"
965 wrk2 = atomic_indices
966 CALL sort_unique(wrk2, unique)
967 IF (((wrk2(1) >= 1) .AND. (wrk2(SIZE(wrk2)) <= natom)) .AND. unique) THEN
968 rab(:) = pbc(particle_set(atomic_indices(1))%r(:), &
969 particle_set(atomic_indices(2))%r(:), cell)
970 dab = sqrt(rab(1)*rab(1) + rab(2)*rab(2) + rab(3)*rab(3))
971 WRITE (unit=iw, fmt="(T3,A,T20,A,3F13.6,3X,A,F13.6)") &
972 "r"//trim(string), "=", rab(:)*conv, &
973 "|r| =", dab*conv
974 ELSE
975 WRITE (unit=iw, fmt="(T3,A)") &
976 "Invalid atomic indices "//trim(string)//" specified. Print request is ignored."
977 END IF
978 END DO
979 END IF
980
981 ! Print the requested angles
982 CALL section_vals_val_get(section_vals=section, &
983 keyword_name="ANGLE", &
984 n_rep_val=n)
985 IF (n > 0) THEN
986 WRITE (unit=iw, fmt="(/,T3,A,/)") &
987 "Angle a(i,j,k) between the atomic distance vectors r(j,i) and "// &
988 "r(j,k) in DEGREE"
989 DO i = 1, n
990 CALL section_vals_val_get(section_vals=section, &
991 keyword_name="ANGLE", &
992 i_rep_val=i, &
993 i_vals=atomic_indices)
994 string = ""
995 WRITE (unit=string, fmt="(A,3(I0,A))") &
996 "(", atomic_indices(1), ",", atomic_indices(2), ",", atomic_indices(3), ")"
997 wrk3 = atomic_indices
998 CALL sort_unique(wrk3, unique)
999 IF (((wrk3(1) >= 1) .AND. (wrk3(SIZE(wrk3)) <= natom)) .AND. unique) THEN
1000 rab(:) = pbc(particle_set(atomic_indices(1))%r(:), &
1001 particle_set(atomic_indices(2))%r(:), cell)
1002 rbc(:) = pbc(particle_set(atomic_indices(2))%r(:), &
1003 particle_set(atomic_indices(3))%r(:), cell)
1004 WRITE (unit=iw, fmt="(T3,A,T26,A,F9.3)") &
1005 "a"//trim(string), "=", angle(-rab, rbc)*degree
1006 ELSE
1007 WRITE (unit=iw, fmt="(T3,A)") &
1008 "Invalid atomic indices "//trim(string)//" specified. Print request is ignored."
1009 END IF
1010 END DO
1011 END IF
1012
1013 ! Print the requested dihedral angles
1014 CALL section_vals_val_get(section_vals=section, &
1015 keyword_name="DIHEDRAL_ANGLE", &
1016 n_rep_val=n)
1017 IF (n > 0) THEN
1018 WRITE (unit=iw, fmt="(/,T3,A,/)") &
1019 "Dihedral angle d(i,j,k,l) between the planes (i,j,k) and (j,k,l) "// &
1020 "in DEGREE"
1021 DO i = 1, n
1022 CALL section_vals_val_get(section_vals=section, &
1023 keyword_name="DIHEDRAL_ANGLE", &
1024 i_rep_val=i, &
1025 i_vals=atomic_indices)
1026 string = ""
1027 WRITE (unit=string, fmt="(A,4(I0,A))") &
1028 "(", atomic_indices(1), ",", atomic_indices(2), ",", &
1029 atomic_indices(3), ",", atomic_indices(4), ")"
1030 wrk4 = atomic_indices
1031 CALL sort_unique(wrk4, unique)
1032 IF (((wrk4(1) >= 1) .AND. (wrk4(SIZE(wrk4)) <= natom)) .AND. unique) THEN
1033 rab(:) = pbc(particle_set(atomic_indices(1))%r(:), &
1034 particle_set(atomic_indices(2))%r(:), cell)
1035 rbc(:) = pbc(particle_set(atomic_indices(2))%r(:), &
1036 particle_set(atomic_indices(3))%r(:), cell)
1037 rcd(:) = pbc(particle_set(atomic_indices(3))%r(:), &
1038 particle_set(atomic_indices(4))%r(:), cell)
1039 WRITE (unit=iw, fmt="(T3,A,T26,A,F9.3)") &
1040 "d"//trim(string), "=", dihedral_angle(rab, rbc, rcd)*degree
1041 ELSE
1042 WRITE (unit=iw, fmt="(T3,A)") &
1043 "Invalid atomic indices "//trim(string)//" specified. Print request is ignored."
1044 END IF
1045 END DO
1046 END IF
1047 END IF
1048 CALL cp_print_key_finished_output(iw, logger, input_section, &
1049 "PRINT%STRUCTURE_DATA")
1050
1051 CALL timestop(handle)
1052
1053 END SUBROUTINE write_structure_data
1054
1055! **************************************************************************************************
1056!> \brief Write the final geometry and cell information to files
1057!> \param particle_set pointer to particles with atm_name, element_symbol and position
1058!> \param cell pointer to cell with abc, angle_alpha, angle_beta, angle_gamma and deth
1059!> \param input_section pointer to motion_section which has PRINT%FINAL_STRUCTURE
1060!> \param conv flag for whether convergence is achieved or not in optimization
1061!> \param keep_angles flag for whether cell optimization keeps initial angles
1062!> \param keep_symmetry flag for whether cell optimization keeps initial symmetry
1063!> \param keep_volume flag for whether cell optimization keeps initial volume
1064!> \param gopt_env_label the geometry optimization label "GEO_OPT", "CELL_OPT", ...
1065!> \param constraint_label label for directions with constraint in cell optimization
1066!> \par Intended to be invoked in gopt_f_methods:write_final_info.
1067!> This implementation does not consider higher space groups even if
1068!> one is detected, and the chemical formulae are neither written in
1069!> the sorted "Hill notation" nor expressed in groups of molecules.
1070!> Other potentially useful but yet to be written information includes:
1071!> the external pressure from CELL_OPT/EXTERNAL_POTENTIAL and the
1072!> stress tensor (virial) for CELL_OPT;
1073!> the fixed atoms from MOTION/CONSTRAINT/FIXED_ATOMS for all.
1074!>
1075!> History
1076!> 04.2026 - Created
1077!> \author HE Zilong
1078!> \version 1.0
1079! **************************************************************************************************
1080 SUBROUTINE write_final_structure(particle_set, cell, input_section, conv, &
1081 keep_angles, keep_symmetry, keep_volume, &
1082 gopt_env_label, constraint_label)
1083 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
1084 TYPE(cell_type), INTENT(IN), POINTER :: cell
1085 TYPE(section_vals_type), INTENT(IN), POINTER :: input_section
1086 LOGICAL, INTENT(IN) :: conv, keep_angles, keep_symmetry, &
1087 keep_volume
1088 CHARACTER(LEN=default_string_length), INTENT(IN) :: gopt_env_label
1089 CHARACTER(LEN=4), INTENT(IN) :: constraint_label
1090
1091 CHARACTER(len=*), PARAMETER :: routinen = 'write_final_structure'
1092
1093 CHARACTER(LEN=2) :: element_symbol
1094 CHARACTER(LEN=2), ALLOCATABLE :: element_list(:)
1095 CHARACTER(LEN=5) :: perd_str
1096 CHARACTER(LEN=:), ALLOCATABLE :: formula_structural, formula_sum
1097 CHARACTER(LEN=default_path_length) :: cell_str, record
1098 CHARACTER(LEN=default_string_length) :: atm_name, f_cif, f_cif_label, &
1099 f_cif_type_symbol, f_xyz
1100 CHARACTER(LEN=default_string_length), ALLOCATABLE :: cif_label(:), cif_type_symbol(:), &
1101 xyz_label(:)
1102 CHARACTER(LEN=timestamp_length) :: timestamp
1103 INTEGER :: elem_seen, file_unit, gcd_all, handle, i, iatom, ielem, natom, output_unit, &
1104 symmetry_id, w_cif_label, w_cif_type_symbol, w_xyz_label
1105 INTEGER, ALLOCATABLE :: count_list(:)
1106 INTEGER, DIMENSION(3) :: periodic
1107 LOGICAL :: dummy, elem_in_list, orthorhombic, &
1108 write_cif, write_xyz
1109 REAL(kind=dp) :: angle_alpha, angle_beta, angle_gamma, &
1110 deth
1111 REAL(kind=dp), DIMENSION(3) :: abc, r, s
1112 REAL(kind=dp), DIMENSION(3, 3) :: hmat
1113 TYPE(cp_logger_type), POINTER :: logger
1114 TYPE(enumeration_type), POINTER :: enum
1115 TYPE(keyword_type), POINTER :: symmetry_keyword
1116 TYPE(section_type), POINTER :: tmp_cell_section
1117 TYPE(section_vals_type), POINTER :: print_key
1118
1119 CALL timeset(routinen, handle)
1120
1121 NULLIFY (enum, logger, symmetry_keyword, print_key, tmp_cell_section)
1122 logger => cp_get_default_logger()
1123 output_unit = cp_logger_get_default_io_unit(logger)
1124 print_key => section_vals_get_subs_vals(input_section, "PRINT%FINAL_STRUCTURE")
1125
1126 ! Collect cell information
1127 perd_str = "F F F"
1128 CALL get_cell(cell, alpha=angle_alpha, beta=angle_beta, gamma=angle_gamma, &
1129 deth=deth, orthorhombic=orthorhombic, abc=abc, periodic=periodic, &
1130 h=hmat, symmetry_id=symmetry_id)
1131 IF (periodic(1) == 1) perd_str(1:1) = "T"
1132 IF (periodic(2) == 1) perd_str(3:3) = "T"
1133 IF (periodic(3) == 1) perd_str(5:5) = "T"
1134 CALL create_cell_section(tmp_cell_section)
1135 symmetry_keyword => section_get_keyword(tmp_cell_section, "SYMMETRY")
1136 CALL keyword_get(symmetry_keyword, enum=enum)
1137 ! cell_str is default_path_length which is longer
1138 ! than default_string_length and should be enough
1139 WRITE (unit=cell_str, fmt="(9(1X,F19.10))") &
1140 cp_unit_from_cp2k(hmat(1, 1), "angstrom"), &
1141 cp_unit_from_cp2k(hmat(2, 1), "angstrom"), &
1142 cp_unit_from_cp2k(hmat(3, 1), "angstrom"), &
1143 cp_unit_from_cp2k(hmat(1, 2), "angstrom"), &
1144 cp_unit_from_cp2k(hmat(2, 2), "angstrom"), &
1145 cp_unit_from_cp2k(hmat(3, 2), "angstrom"), &
1146 cp_unit_from_cp2k(hmat(1, 3), "angstrom"), &
1147 cp_unit_from_cp2k(hmat(2, 3), "angstrom"), &
1148 cp_unit_from_cp2k(hmat(3, 3), "angstrom")
1149
1150 ! Collect atom information
1151 natom = SIZE(particle_set)
1152 ALLOCATE (element_list(nelem + 1), count_list(nelem + 1))
1153 count_list(:) = 0
1154 ALLOCATE (cif_label(natom), cif_type_symbol(natom), xyz_label(natom))
1155 elem_seen = 0
1156 w_cif_type_symbol = 0
1157 w_cif_label = 0
1158 w_xyz_label = 0
1159 atom_loop: DO iatom = 1, natom
1160 elem_in_list = .false.
1161 CALL get_atomic_kind(atomic_kind=particle_set(iatom)%atomic_kind, &
1162 name=atm_name, element_symbol=element_symbol)
1163 cif_type_symbol(iatom) = trim(atm_name)
1164 ! From write_particle_coordinates above it seems possible
1165 ! for some atoms to have empty element symbols; whatever
1166 ! these are, do not count them in the chemical formula
1167 IF (len_trim(element_symbol) == 0) THEN
1168 dummy = qmmm_ff_precond_only_qm(id1=atm_name)
1169 cif_label(iatom) = trim(atm_name)//trim(adjustl(cp_to_string(iatom)))
1170 xyz_label(iatom) = trim(atm_name)
1171 ELSE
1172 cif_label(iatom) = trim(element_symbol)//trim(adjustl(cp_to_string(iatom)))
1173 xyz_label(iatom) = trim(element_symbol)
1174 elem_loop: DO ielem = 1, elem_seen
1175 IF (element_list(ielem) == element_symbol) THEN
1176 elem_in_list = .true.
1177 count_list(ielem) = count_list(ielem) + 1
1178 EXIT elem_loop
1179 END IF
1180 END DO elem_loop
1181 IF (.NOT. elem_in_list) THEN
1182 elem_seen = elem_seen + 1
1183 element_list(elem_seen) = element_symbol
1184 count_list(elem_seen) = 1
1185 END IF
1186 END IF
1187 IF (len_trim(cif_type_symbol(iatom)) > w_cif_type_symbol) &
1188 w_cif_type_symbol = len_trim(cif_type_symbol(iatom))
1189 IF (len_trim(cif_label(iatom)) > w_cif_label) &
1190 w_cif_label = len_trim(cif_label(iatom))
1191 IF (len_trim(xyz_label(iatom)) > w_xyz_label) &
1192 w_xyz_label = len_trim(xyz_label(iatom))
1193 END DO atom_loop
1194
1195 ! Determine the format of each line in cif considering width of cif_type_symbol and cif_label
1196 ! The fields are, in order:
1197 ! _atom_site_type_symbol, _atom_site_label, _atom_site_symmetry_multiplicity,
1198 ! _atom_site_fract_x, _atom_site_fract_y, _atom_site_fract_z, _atom_site_occupancy
1199 ! in which:
1200 ! _atom_site_type_symbol is taken as atm_name
1201 ! _atom_site_label is taken as element_symbol//iatom
1202 ! _atom_site_symmetry_multiplicity and _atom_site_occupancy are always 1
1203 f_cif_type_symbol = "A"//trim(adjustl(cp_to_string(w_cif_type_symbol + 4)))
1204 f_cif_label = "A"//trim(adjustl(cp_to_string(w_cif_label + 4)))
1205 f_cif = "(T3,"//trim(f_cif_type_symbol)//","//trim(f_cif_label)//",I4,3F14.8,F8.2)"
1206
1207 ! Determine the format of each line in xyz
1208 f_xyz = "(T2,A"//trim(adjustl(cp_to_string(w_xyz_label + 4)))//",1X,3F20.10)"
1209
1210 ! Determine formula_sum
1211 cpassert(elem_seen > 0)
1212 cpassert(count_list(1) > 0)
1213 formula_sum = "'"
1214 DO ielem = 1, elem_seen
1215 formula_sum = formula_sum//trim(adjustl(element_list(ielem)))
1216 formula_sum = formula_sum//trim(adjustl(cp_to_string(count_list(ielem))))
1217 formula_sum = formula_sum//" "
1218 END DO
1219 formula_sum = trim(adjustl(formula_sum))//"'"
1220
1221 ! Determine formula_structural and Z
1222 gcd_all = count_list(1)
1223 DO ielem = 1, elem_seen
1224 IF (count_list(ielem) /= 0) THEN
1225 gcd_all = gcd(gcd_all, count_list(ielem))
1226 END IF
1227 END DO
1228 IF (gcd_all > 1) count_list = count_list/gcd_all
1229 formula_structural = "'"
1230 DO ielem = 1, elem_seen
1231 formula_structural = formula_structural//trim(adjustl(element_list(ielem)))
1232 formula_structural = formula_structural//trim(adjustl(cp_to_string(count_list(ielem))))
1233 formula_structural = formula_structural//" "
1234 END DO
1235 formula_structural = trim(adjustl(formula_structural))//"'"
1236
1237 ! Write XYZ
1238 CALL section_vals_val_get(print_key, "PRINT_XYZ", l_val=write_xyz)
1239 IF (write_xyz) THEN
1240 ! Print a message to log
1241 record = cp_print_key_generate_filename(logger, print_key, &
1242 extension=".xyz", &
1243 my_local=.false.)
1244 IF (output_unit > 0) THEN
1245 IF (conv) THEN
1246 WRITE (unit=output_unit, fmt="(/,T2,A)") &
1247 routinen//": Optimization converged, writing XYZ file gladly:"
1248 ELSE
1249 WRITE (unit=output_unit, fmt="(/,T2,A)") &
1250 routinen//": Optimization not yet converged, writing XYZ file anyway:"
1251 END IF
1252 WRITE (unit=output_unit, fmt="(T3,A)") trim(record)
1253 END IF
1254
1255 ! Make timestamp for the file
1256 CALL m_timestamp(timestamp)
1257
1258 ! Prepare file unit and write to it
1259 file_unit = cp_print_key_unit_nr(logger, input_section, "PRINT%FINAL_STRUCTURE", &
1260 file_status="REPLACE", extension=".xyz")
1261 IF (file_unit > 0) THEN
1262 WRITE (unit=file_unit, fmt="(I8)") &
1263 natom
1264 WRITE (unit=file_unit, fmt="(A)") &
1265 'Lattice="'//trim(adjustl(cell_str))// &
1266 '" Properties=species:S:1:pos:R:3 pbc="'// &
1267 trim(perd_str)//'"'
1268 DO iatom = 1, natom
1269 DO i = 1, 3
1270 r(i) = cp_unit_from_cp2k(particle_set(iatom)%r(i), "angstrom")
1271 END DO
1272 WRITE (unit=file_unit, fmt=trim(f_xyz)) &
1273 xyz_label(iatom), r(1:3)
1274 END DO
1275 END IF
1276 END IF
1277
1278 ! Write CIF
1279 CALL section_vals_val_get(print_key, "PRINT_CIF", l_val=write_cif)
1280 IF (write_cif) THEN
1281 ! Print a message to log
1282 record = cp_print_key_generate_filename(logger, print_key, &
1283 extension=".cif", &
1284 my_local=.false.)
1285 IF (output_unit > 0) THEN
1286 IF (conv) THEN
1287 WRITE (unit=output_unit, fmt="(/,T2,A)") &
1288 routinen//": Optimization converged, writing CIF file gladly:"
1289 ELSE
1290 WRITE (unit=output_unit, fmt="(/,T2,A)") &
1291 routinen//": Optimization not yet converged, writing CIF file anyway:"
1292 END IF
1293 WRITE (unit=output_unit, fmt="(T3,A)") trim(record)
1294 END IF
1295
1296 ! Make timestamp for the file
1297 CALL m_timestamp(timestamp)
1298
1299 ! Prepare file unit and write to it
1300 file_unit = cp_print_key_unit_nr(logger, input_section, "PRINT%FINAL_STRUCTURE", &
1301 file_status="REPLACE", extension=".cif")
1302 IF (file_unit > 0) THEN
1303 ! Generic information
1304 WRITE (unit=file_unit, fmt="(A)") &
1305 "# CIF file created by CP2K "//trim(modulen)//":"//trim(routinen)
1306 WRITE (unit=file_unit, fmt="(A)") &
1307 "data_"//trim(logger%iter_info%project_name)
1308 WRITE (unit=file_unit, fmt="(A,T39,A)") &
1309 "_audit_creation_date", timestamp(:10)
1310 WRITE (unit=file_unit, fmt="(A,/,A,/,A)") &
1311 "_audit_creation_method", ";", &
1312 trim(cp2k_version)//" (revision "//trim(compile_revision)//")"
1313 WRITE (unit=file_unit, fmt="(A,/,A,/,A,/,A)") &
1314 "Project name "//trim(logger%iter_info%project_name), &
1315 "submitted by "//trim(r_user_name)//"@"//trim(r_host_name), &
1316 "processed in "//trim(r_cwd), &
1317 "generated at "//trim(timestamp)
1318 WRITE (unit=file_unit, fmt="(T2,A)") &
1319 "- Optimization type: "//trim(gopt_env_label)
1320 IF (conv) THEN
1321 WRITE (unit=file_unit, fmt="(T2,A)") &
1322 "- Optimization converged: TRUE"
1323 ELSE
1324 WRITE (unit=file_unit, fmt="(T2,A)") &
1325 "- Optimization converged: FALSE"
1326 END IF
1327 WRITE (unit=file_unit, fmt="(T2,A)") &
1328 "- Requested initial cell symmetry: "//trim(enum_i2c(enum, symmetry_id))
1329 IF (orthorhombic) THEN
1330 WRITE (unit=file_unit, fmt="(T2,A)") &
1331 "- Cell is numerically orthorhombic: TRUE"
1332 ELSE
1333 WRITE (unit=file_unit, fmt="(T2,A)") &
1334 "- Cell is numerically orthorhombic: FALSE"
1335 END IF
1336 WRITE (unit=file_unit, fmt="(T2,A)") &
1337 "- Periodicity of cell: "//trim(perd_str)
1338 IF (gopt_env_label == "CELL_OPT") THEN
1339 WRITE (unit=file_unit, fmt="(T2,A)") &
1340 "- Cell is subject to optimization: TRUE"
1341 WRITE (unit=file_unit, fmt="(T2,A)") &
1342 "- Cell has constraint on direction: "//trim(adjustl(constraint_label))
1343 IF (keep_angles) THEN
1344 WRITE (unit=file_unit, fmt="(T2,A)") &
1345 "- Keep angles between the cell vectors during optimization: TRUE"
1346 ELSE
1347 WRITE (unit=file_unit, fmt="(T2,A)") &
1348 "- Keep angles between the cell vectors during optimization: FALSE"
1349 END IF
1350 IF (keep_symmetry) THEN
1351 WRITE (unit=file_unit, fmt="(T2,A)") &
1352 "- Keep initial cell symmetry during optimization: TRUE"
1353 ELSE
1354 WRITE (unit=file_unit, fmt="(T2,A)") &
1355 "- Keep initial cell symmetry during optimization: FALSE"
1356 END IF
1357 IF (keep_volume) THEN
1358 WRITE (unit=file_unit, fmt="(T2,A)") &
1359 "- Keep initial cell volume during optimization: TRUE"
1360 ELSE
1361 WRITE (unit=file_unit, fmt="(T2,A)") &
1362 "- Keep initial cell volume during optimization: FALSE"
1363 END IF
1364 ELSE
1365 WRITE (unit=file_unit, fmt="(T2,A)") &
1366 "- Cell is subject to optimization: FALSE"
1367 END IF
1368 WRITE (unit=file_unit, fmt="(T2,A)") &
1369 "- Final cell vectors A, B, C by rows [angstrom]:"
1370 DO i = 1, 3
1371 WRITE (unit=file_unit, fmt="(T3,3(1X,F19.10))") &
1372 cp_unit_from_cp2k(hmat(1, i), "angstrom"), &
1373 cp_unit_from_cp2k(hmat(2, i), "angstrom"), &
1374 cp_unit_from_cp2k(hmat(3, i), "angstrom")
1375 END DO
1376 WRITE (unit=file_unit, fmt="(A)") ";"
1377 ! Data of cell and geometry
1378 WRITE (unit=file_unit, fmt="(/,A,T44,A)") &
1379 "_symmetry_space_group_name_H-M", "'P 1'"
1380 WRITE (unit=file_unit, fmt="(A,T31,F18.8)") &
1381 "_cell_length_a", cp_unit_from_cp2k(abc(1), "angstrom")
1382 WRITE (unit=file_unit, fmt="(A,T31,F18.8)") &
1383 "_cell_length_b", cp_unit_from_cp2k(abc(2), "angstrom")
1384 WRITE (unit=file_unit, fmt="(A,T31,F18.8)") &
1385 "_cell_length_c", cp_unit_from_cp2k(abc(3), "angstrom")
1386 WRITE (unit=file_unit, fmt="(A,T31,F18.8)") &
1387 "_cell_angle_alpha", angle_alpha
1388 WRITE (unit=file_unit, fmt="(A,T31,F18.8)") &
1389 "_cell_angle_beta", angle_beta
1390 WRITE (unit=file_unit, fmt="(A,T31,F18.8)") &
1391 "_cell_angle_gamma", angle_gamma
1392 WRITE (unit=file_unit, fmt="(A,T48,A)") &
1393 "_symmetry_Int_Tables_number", "1"
1394 WRITE (unit=file_unit, fmt="(A,T36,A)") &
1395 "_chemical_formula_structural", formula_structural
1396 WRITE (unit=file_unit, fmt="(A,T36,A)") &
1397 "_chemical_formula_sum", formula_sum
1398 WRITE (unit=file_unit, fmt="(A,T31,F18.8)") &
1399 "_cell_volume", cp_unit_from_cp2k(abs(deth), "angstrom^3")
1400 WRITE (unit=file_unit, fmt="(A,T41,I8)") &
1401 "_cell_formula_units_Z", gcd_all
1402 WRITE (unit=file_unit, fmt="(A,/,T2,A,/,T2,A,/,T3,A)") &
1403 "loop_", "_symmetry_equiv_pos_site_id", &
1404 "_symmetry_equiv_pos_as_xyz", "1 'x, y, z'"
1405 WRITE (unit=file_unit, fmt="(A,/,T2,A,/,T2,A,/,T2,A,/,T2,A,/,T2,A,/,T2,A,/,T2,A)") &
1406 "loop_", "_atom_site_type_symbol", "_atom_site_label", &
1407 "_atom_site_symmetry_multiplicity", "_atom_site_fract_x", &
1408 "_atom_site_fract_y", "_atom_site_fract_z", "_atom_site_occupancy"
1409 DO iatom = 1, natom
1410 r(1:3) = pbc(particle_set(iatom)%r(1:3), cell)
1411 CALL real_to_scaled(s, r, cell)
1412 DO i = 1, 3
1413 s(i) = modulo(s(i), 1.0_dp)
1414 END DO
1415 WRITE (unit=file_unit, fmt=trim(f_cif)) &
1416 cif_type_symbol(iatom), cif_label(iatom), 1, s(1:3), 1.0_dp
1417 END DO
1418 END IF
1419 END IF
1420
1421 ! Finish
1422 DEALLOCATE (element_list, count_list, formula_structural, &
1423 formula_sum, cif_label, cif_type_symbol, &
1424 xyz_label)
1425 CALL section_release(tmp_cell_section)
1426 CALL cp_print_key_finished_output(file_unit, logger, input_section, &
1427 "PRINT%FINAL_STRUCTURE")
1428 IF (output_unit > 0) &
1429 WRITE (unit=output_unit, fmt='(/,T2,A)') &
1430 routinen//": Done!"
1431
1432 CALL timestop(handle)
1433
1434 END SUBROUTINE write_final_structure
1435
1436END 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, npgf_seg_sum, ccon)
...
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:535
subroutine, public cell_release(cell)
releases the given cell (see doc/ReferenceCounting.html)
Definition cell_types.F:608
subroutine, public cell_clone(cell_in, cell_out, tag)
Clone cell variable.
Definition cell_types.F:118
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:210
some minimal info about CP2K, including its version and license
Definition cp2k_info.F:20
character(len=default_string_length), public r_host_name
Definition cp2k_info.F:72
character(len= *), parameter, public compile_revision
Definition cp2k_info.F:43
character(len= *), parameter, public cp2k_version
Definition cp2k_info.F:47
character(len=default_path_length), public r_cwd
Definition cp2k_info.F:71
character(len=default_string_length), public r_user_name
Definition cp2k_info.F:72
various routines to log and control the output. The idea is that decisions about where to log should ...
integer function, public cp_logger_get_default_io_unit(logger)
returns the unit nr for the ionode (-1 on all other processors) skips as well checks if the procs cal...
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)
...
character(len=default_path_length) function, public cp_print_key_generate_filename(logger, print_key, middle_name, extension, my_local)
Utility function that returns a unit number to write the print key. Might open a file with a unique f...
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:1178
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_extxyz
integer, parameter, public dump_atomic
integer, parameter, public dump_dcd_aligned_cell
integer, parameter, public dump_dcd
builds the subsystem section of the input
subroutine, public create_cell_section(section, periodic)
creates the cell section
represents an enumeration, i.e. a mapping between integers and strings
character(len=default_string_length) function, public enum_i2c(enum, i)
maps an integer to a string
represents keywords in an input
subroutine, public keyword_get(keyword, names, usage, description, type_of_var, n_var, default_value, lone_keyword_value, repeats, enum, citations)
...
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
recursive subroutine, public section_release(section)
releases the given keyword list (see doc/ReferenceCounting.html)
recursive type(keyword_type) function, pointer, public section_get_keyword(section, keyword_name)
returns the requested keyword
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 default_path_length
Definition kinds.F:58
integer, parameter, public sp
Definition kinds.F:33
Machine interface based on Fortran 2003 and POSIX.
Definition machine.F:17
integer, parameter, public timestamp_length
Definition machine.F:58
subroutine, public m_timestamp(timestamp)
Returns a human readable timestamp.
Definition machine.F:393
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:184
elemental integer function, public gcd(a, b)
computes the greatest common divisor of two number
Definition mathlib.F:1288
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:476
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_final_structure(particle_set, cell, input_section, conv, keep_angles, keep_symmetry, keep_volume, gopt_env_label, constraint_label)
Write the final geometry and cell information to files.
subroutine, public write_particle_matrix(matrix, particle_set, iw, el_per_part, ilist, parts_per_line)
...
subroutine, public get_particle_set(particle_set, qs_kind_set, first_sgf, last_sgf, nsgf, nmao, basis, ncgf)
Get the components of a particle set.
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 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...
Periodic Table related data definitions.
integer, parameter, public nelem
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, cneo_potential, se_parameter, dftb_parameter, xtb_parameter, dftb3_param, zatom, 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_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, monovalent, floating, name, element_symbol, pao_basis_size, pao_model_file, 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:60
type of a logger, at the moment it contains just a print level starting at which level it should be l...
represent a keyword in the input
represent a section of the input file
Provides all information about a quickstep kind.