81#include "./base/base_uses.f90"
98 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'particle_methods'
123 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
124 INTEGER,
DIMENSION(:),
INTENT(INOUT),
OPTIONAL :: first_sgf, last_sgf, nsgf, nmao
126 INTEGER,
DIMENSION(:),
INTENT(INOUT),
OPTIONAL :: ncgf
128 INTEGER :: ikind, iparticle, isgf, nparticle, ns
130 cpassert(
ASSOCIATED(particle_set))
132 nparticle =
SIZE(particle_set)
133 IF (
PRESENT(first_sgf))
THEN
134 cpassert(
SIZE(first_sgf) >= nparticle)
136 IF (
PRESENT(last_sgf))
THEN
137 cpassert(
SIZE(last_sgf) >= nparticle)
139 IF (
PRESENT(nsgf))
THEN
140 cpassert(
SIZE(nsgf) >= nparticle)
142 IF (
PRESENT(nmao))
THEN
143 cpassert(
SIZE(nmao) >= nparticle)
145 IF (
PRESENT(ncgf))
THEN
146 cpassert(
SIZE(ncgf) >= nparticle)
149 IF (
PRESENT(first_sgf) .OR.
PRESENT(last_sgf) .OR.
PRESENT(nsgf))
THEN
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
162 IF (
PRESENT(nsgf)) nsgf(iparticle) = ns
163 IF (
PRESENT(first_sgf)) first_sgf(iparticle) = isgf + 1
165 IF (
PRESENT(last_sgf)) last_sgf(iparticle) = isgf
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
185 IF (
PRESENT(first_sgf))
THEN
186 IF (
SIZE(first_sgf) > nparticle) first_sgf(nparticle + 1) = isgf + 1
189 IF (
PRESENT(nmao))
THEN
190 DO iparticle = 1, nparticle
191 CALL get_atomic_kind(particle_set(iparticle)%atomic_kind, kind_number=ikind)
223 content, title, cell, array, unit_conv, &
224 charge_occup, charge_beta, &
225 charge_extended, print_kind)
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
236 CHARACTER(len=*),
PARAMETER :: routinen =
'write_particle_coordinates'
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, &
246 REAL(kind=
dp) :: angle_alpha, angle_beta, angle_gamma, &
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
256 CALL timeset(routinen, handle)
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")
266 IF (
PRESENT(unit_conv))
THEN
269 SELECT CASE (output_format)
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)
277 element_symbol=element_symbol)
278 IF (len_trim(element_symbol) == 0 .OR. my_print_kind)
THEN
283 atm_name = trim(atm_name)
285 my_format =
"(T2,A2,"
286 atm_name = trim(element_symbol)
288 SELECT CASE (trim(content))
290 IF (
PRESENT(array))
THEN
293 r(:) = particle_set(iatom)%r(:)
295 WRITE (iunit, trim(my_format)//
"1X,3F20.10)") trim(atm_name), r(1:3)*factor
297 IF (
PRESENT(array))
THEN
300 v(:) = particle_set(iatom)%v(:)
302 WRITE (iunit, trim(my_format)//
"1X,3F20.10)") trim(atm_name), v(1:3)*factor
304 IF (
PRESENT(array))
THEN
305 f(:) = array((iatom - 1)*3 + 1:(iatom - 1)*3 + 3)
307 f(:) = particle_set(iatom)%f(:)
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)
314 f(:) = particle_set(iatom)%f(:)
316 WRITE (iunit, trim(my_format)//
"1X,3F20.10)") trim(atm_name), f(1:3)*factor
321 SELECT CASE (trim(content))
323 IF (
PRESENT(array))
THEN
326 r(:) = particle_set(iatom)%r(:)
328 WRITE (iunit,
"(3F20.10)") r(1:3)*factor
330 IF (
PRESENT(array))
THEN
333 v(:) = particle_set(iatom)%v(:)
335 WRITE (iunit,
"(3F20.10)") v(1:3)*factor
337 IF (
PRESENT(array))
THEN
338 f(:) = array((iatom - 1)*3 + 1:(iatom - 1)*3 + 3)
340 f(:) = particle_set(iatom)%f(:)
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)
347 f(:) = particle_set(iatom)%f(:)
349 WRITE (iunit,
"(3F20.10)") f(1:3)*factor
353 IF (.NOT. (
PRESENT(cell)))
THEN
354 cpabort(
"Cell is not present! Report this bug!")
356 CALL get_cell(cell, alpha=angle_alpha, beta=angle_beta,
gamma=angle_gamma, &
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
371 h(1:3, 1:3) = matmul(cell_dcd%hmat(1:3, 1:3), cell%h_inv(1:3, 1:3))
374 ALLOCATE (arr(3, natom))
375 IF (
PRESENT(array))
THEN
376 arr(1:3, 1:natom) = reshape(array, [3, natom])
378 SELECT CASE (trim(content))
381 arr(1:3, iatom) = particle_set(iatom)%r(1:3)
385 arr(1:3, iatom) = particle_set(iatom)%v(1:3)
389 arr(1:3, iatom) = particle_set(iatom)%f(1:3)
392 cpabort(
"Illegal DCD dump type")
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)
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)
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)
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)
428 CALL get_cell(cell, alpha=angle_alpha, beta=angle_beta,
gamma=angle_gamma, abc=abc)
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 "
463 element_symbol=element_symbol, name=atm_name, &
464 fist_potential=fist_potential, shell=shell)
465 IF (len_trim(element_symbol) == 0)
THEN
468 name = trim(atm_name)
469 IF (
ASSOCIATED(fist_potential))
THEN
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)
480 SELECT CASE (trim(content))
482 IF (
PRESENT(array))
THEN
485 r(:) = particle_set(iatom)%r(:)
487 WRITE (unit=line(31:54), fmt=
"(3F8.3)") r(1:3)*factor
489 cpabort(
"PDB dump only for trajectory available")
491 IF (my_charge_occup)
THEN
492 WRITE (unit=line(55:60), fmt=
"(F6.2)") qeff
494 WRITE (unit=line(55:60), fmt=
"(F6.2)") 0.0_dp
496 IF (my_charge_beta)
THEN
497 WRITE (unit=line(61:66), fmt=
"(F6.2)") qeff
499 WRITE (unit=line(61:66), fmt=
"(F6.2)") 0.0_dp
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
506 WRITE (unit=iunit, fmt=
"(A)") trim(line)
508 WRITE (unit=iunit, fmt=
"(A)")
"END"
510 cpabort(
"Illegal dump type")
513 CALL timestop(handle)
529 REAL(kind=
dp),
DIMENSION(:),
OPTIONAL :: charges
531 CHARACTER(LEN=default_string_length) :: name, unit_str
532 INTEGER :: iatom, ikind, iw, natom
533 REAL(kind=
dp) :: conv, mass, qcore, qeff, qshell
542 "PRINT%ATOMIC_COORDINATES", extension=
".coordLog")
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)
560 IF (
PRESENT(charges)) qeff = charges(iatom)
561 IF (
ASSOCIATED(shell_kind))
THEN
565 qeff = qcore + qshell
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
574 "PRINT%ATOMIC_COORDINATES")
591 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
593 CHARACTER(LEN=*),
INTENT(IN) :: label
595 CHARACTER(len=*),
PARAMETER :: routinen =
'write_qs_particle_coordinates'
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
603 CALL timeset(routinen, handle)
608 "PRINT%ATOMIC_COORDINATES", extension=
".coordLog")
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)
622 element_symbol=element_symbol, &
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
633 "PRINT%ATOMIC_COORDINATES")
635 CALL timestop(handle)
654 CHARACTER(len=*),
PARAMETER :: routinen =
'write_particle_distances'
656 CHARACTER(LEN=default_string_length) :: unit_str
657 INTEGER :: handle, iatom, iw, jatom, natom
658 INTEGER,
DIMENSION(3) :: periodic
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
665 CALL timeset(routinen, handle)
667 cpassert(
ASSOCIATED(particle_set))
668 cpassert(
ASSOCIATED(cell))
669 cpassert(
ASSOCIATED(subsys_section))
674 "PRINT%INTERATOMIC_DISTANCES", extension=
".distLog")
678 CALL section_vals_val_get(subsys_section,
"PRINT%INTERATOMIC_DISTANCES%CHECK_INTERATOMIC_DISTANCES", &
679 r_val=dab_min, explicit=explicit)
683 natom =
SIZE(particle_set)
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
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)
698 ALLOCATE (distance_matrix(natom, natom))
699 distance_matrix(:, :) = 0.0_dp
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
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 "// &
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)))
718 IF (dab < dab_warn)
THEN
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 "// &
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)))
729 distance_matrix(iatom, jatom) = dab
730 distance_matrix(jatom, iatom) = distance_matrix(iatom, jatom)
736 WRITE (unit=iw, fmt=
"(/,/,T2,A)") &
737 "INTERATOMIC DISTANCES IN "//trim(unit_str)
739 IF (
ALLOCATED(distance_matrix))
DEALLOCATE (distance_matrix)
742 "PRINT%INTERATOMIC_DISTANCES")
745 CALL timestop(handle)
759 REAL(kind=
dp),
DIMENSION(:, :) :: matrix
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
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, &
771 INTEGER,
DIMENSION(:),
POINTER :: my_list
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
784 natom =
SIZE(particle_set)
786 ALLOCATE (my_list(natom))
787 IF (
PRESENT(ilist))
THEN
794 natom = natom*my_el_per_part
795 DO jatom = 1, natom, my_parts_per_line
797 to = min(from + my_parts_per_line - 1, natom)
798 WRITE (unit=iw, fmt=trim(fmt_string1)) (icol, icol=from, to)
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)
833 CHARACTER(LEN=*),
PARAMETER :: routinen =
'write_structure_data'
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), &
839 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: work
840 INTEGER,
DIMENSION(:),
POINTER :: atomic_indices, index_list
842 REAL(kind=
dp) :: conv, dab
843 REAL(kind=
dp),
DIMENSION(3) :: r, rab, rbc, rcd, s
847 CALL timeset(routinen, handle)
848 NULLIFY (atomic_indices)
856 basis_section=input_section, &
857 print_key_path=
"PRINT%STRUCTURE_DATA", &
858 extension=
".coordLog")
864 natom =
SIZE(particle_set)
866 subsection_name=
"PRINT%STRUCTURE_DATA")
868 WRITE (unit=iw, fmt=
"(/,T2,A)")
"REQUESTED STRUCTURE DATA"
871 keyword_name=
"POSITION", &
874 WRITE (unit=iw, fmt=
"(/,T3,A,/)") &
875 "Position vectors r(i) of the atoms i in "//trim(unit_str)
879 keyword_name=
"POSITION", &
881 i_vals=atomic_indices)
882 n_vals =
SIZE(atomic_indices)
883 new_size = old_size + n_vals
885 index_list(old_size + 1:new_size) = atomic_indices(1:n_vals)
888 ALLOCATE (work(new_size))
889 CALL sort(index_list, new_size, work)
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."
900 IF (index_list(i) == index_list(i - 1)) cycle
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
905 DEALLOCATE (index_list)
910 keyword_name=
"POSITION_SCALED", &
913 WRITE (unit=iw, fmt=
"(/,T3,A,/)") &
914 "Position vectors s(i) of the atoms i in scaled coordinates"
918 keyword_name=
"POSITION_SCALED", &
920 i_vals=atomic_indices)
921 n_vals =
SIZE(atomic_indices)
922 new_size = old_size + n_vals
924 index_list(old_size + 1:new_size) = atomic_indices(1:n_vals)
927 ALLOCATE (work(new_size))
928 CALL sort(index_list, new_size, work)
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."
939 IF (index_list(i) == index_list(i - 1)) cycle
941 r(1:3) =
pbc(particle_set(index_list(i))%r(1:3), cell)
943 WRITE (unit=iw, fmt=
"(T3,A,T20,A,3F13.6)") &
944 "s"//trim(string),
"=", s(1:3)
946 DEALLOCATE (index_list)
951 keyword_name=
"DISTANCE", &
954 WRITE (unit=iw, fmt=
"(/,T3,A,/)") &
955 "Distance vector r(i,j) between the atom i and j in "// &
959 keyword_name=
"DISTANCE", &
961 i_vals=atomic_indices)
963 WRITE (unit=string, fmt=
"(A,2(I0,A))") &
964 "(", atomic_indices(1),
",", atomic_indices(2),
")"
965 wrk2 = atomic_indices
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, &
975 WRITE (unit=iw, fmt=
"(T3,A)") &
976 "Invalid atomic indices "//trim(string)//
" specified. Print request is ignored."
983 keyword_name=
"ANGLE", &
986 WRITE (unit=iw, fmt=
"(/,T3,A,/)") &
987 "Angle a(i,j,k) between the atomic distance vectors r(j,i) and "// &
991 keyword_name=
"ANGLE", &
993 i_vals=atomic_indices)
995 WRITE (unit=string, fmt=
"(A,3(I0,A))") &
996 "(", atomic_indices(1),
",", atomic_indices(2),
",", atomic_indices(3),
")"
997 wrk3 = atomic_indices
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)") &
1007 WRITE (unit=iw, fmt=
"(T3,A)") &
1008 "Invalid atomic indices "//trim(string)//
" specified. Print request is ignored."
1015 keyword_name=
"DIHEDRAL_ANGLE", &
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) "// &
1023 keyword_name=
"DIHEDRAL_ANGLE", &
1025 i_vals=atomic_indices)
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
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)") &
1042 WRITE (unit=iw, fmt=
"(T3,A)") &
1043 "Invalid atomic indices "//trim(string)//
" specified. Print request is ignored."
1049 "PRINT%STRUCTURE_DATA")
1051 CALL timestop(handle)
1081 keep_angles, keep_symmetry, keep_volume, &
1082 gopt_env_label, constraint_label)
1084 TYPE(
cell_type),
INTENT(IN),
POINTER :: cell
1086 LOGICAL,
INTENT(IN) :: conv, keep_angles, keep_symmetry, &
1088 CHARACTER(LEN=default_string_length),
INTENT(IN) :: gopt_env_label
1089 CHARACTER(LEN=4),
INTENT(IN) :: constraint_label
1091 CHARACTER(len=*),
PARAMETER :: routinen =
'write_final_structure'
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(:), &
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, &
1111 REAL(kind=
dp),
DIMENSION(3) :: abc, r, s
1112 REAL(kind=
dp),
DIMENSION(3, 3) :: hmat
1119 CALL timeset(routinen, handle)
1121 NULLIFY (enum, logger, symmetry_keyword, print_key, tmp_cell_section)
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"
1139 WRITE (unit=cell_str, fmt=
"(9(1X,F19.10))") &
1151 natom =
SIZE(particle_set)
1152 ALLOCATE (element_list(
nelem + 1), count_list(
nelem + 1))
1154 ALLOCATE (cif_label(natom), cif_type_symbol(natom), xyz_label(natom))
1156 w_cif_type_symbol = 0
1159 atom_loop:
DO iatom = 1, natom
1160 elem_in_list = .false.
1162 name=atm_name, element_symbol=element_symbol)
1163 cif_type_symbol(iatom) = trim(atm_name)
1167 IF (len_trim(element_symbol) == 0)
THEN
1169 cif_label(iatom) = trim(atm_name)//trim(adjustl(
cp_to_string(iatom)))
1170 xyz_label(iatom) = trim(atm_name)
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
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
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))
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)"
1208 f_xyz =
"(T2,A"//trim(adjustl(
cp_to_string(w_xyz_label + 4)))//
",1X,3F20.10)"
1211 cpassert(elem_seen > 0)
1212 cpassert(count_list(1) > 0)
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//
" "
1219 formula_sum = trim(adjustl(formula_sum))//
"'"
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))
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//
" "
1235 formula_structural = trim(adjustl(formula_structural))//
"'"
1244 IF (output_unit > 0)
THEN
1246 WRITE (unit=output_unit, fmt=
"(/,T2,A)") &
1247 routinen//
": Optimization converged, writing XYZ file gladly:"
1249 WRITE (unit=output_unit, fmt=
"(/,T2,A)") &
1250 routinen//
": Optimization not yet converged, writing XYZ file anyway:"
1252 WRITE (unit=output_unit, fmt=
"(T3,A)") trim(record)
1260 file_status=
"REPLACE", extension=
".xyz")
1261 IF (file_unit > 0)
THEN
1262 WRITE (unit=file_unit, fmt=
"(I8)") &
1264 WRITE (unit=file_unit, fmt=
"(A)") &
1265 'Lattice="'//trim(adjustl(cell_str))// &
1266 '" Properties=species:S:1:pos:R:3 pbc="'// &
1272 WRITE (unit=file_unit, fmt=trim(f_xyz)) &
1273 xyz_label(iatom), r(1:3)
1285 IF (output_unit > 0)
THEN
1287 WRITE (unit=output_unit, fmt=
"(/,T2,A)") &
1288 routinen//
": Optimization converged, writing CIF file gladly:"
1290 WRITE (unit=output_unit, fmt=
"(/,T2,A)") &
1291 routinen//
": Optimization not yet converged, writing CIF file anyway:"
1293 WRITE (unit=output_unit, fmt=
"(T3,A)") trim(record)
1301 file_status=
"REPLACE", extension=
".cif")
1302 IF (file_unit > 0)
THEN
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",
";", &
1313 WRITE (unit=file_unit, fmt=
"(A,/,A,/,A,/,A)") &
1314 "Project name "//trim(logger%iter_info%project_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)
1321 WRITE (unit=file_unit, fmt=
"(T2,A)") &
1322 "- Optimization converged: TRUE"
1324 WRITE (unit=file_unit, fmt=
"(T2,A)") &
1325 "- Optimization converged: FALSE"
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"
1333 WRITE (unit=file_unit, fmt=
"(T2,A)") &
1334 "- Cell is numerically orthorhombic: FALSE"
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"
1347 WRITE (unit=file_unit, fmt=
"(T2,A)") &
1348 "- Keep angles between the cell vectors during optimization: FALSE"
1350 IF (keep_symmetry)
THEN
1351 WRITE (unit=file_unit, fmt=
"(T2,A)") &
1352 "- Keep initial cell symmetry during optimization: TRUE"
1354 WRITE (unit=file_unit, fmt=
"(T2,A)") &
1355 "- Keep initial cell symmetry during optimization: FALSE"
1357 IF (keep_volume)
THEN
1358 WRITE (unit=file_unit, fmt=
"(T2,A)") &
1359 "- Keep initial cell volume during optimization: TRUE"
1361 WRITE (unit=file_unit, fmt=
"(T2,A)") &
1362 "- Keep initial cell volume during optimization: FALSE"
1365 WRITE (unit=file_unit, fmt=
"(T2,A)") &
1366 "- Cell is subject to optimization: FALSE"
1368 WRITE (unit=file_unit, fmt=
"(T2,A)") &
1369 "- Final cell vectors A, B, C by rows [angstrom]:"
1371 WRITE (unit=file_unit, fmt=
"(T3,3(1X,F19.10))") &
1376 WRITE (unit=file_unit, fmt=
"(A)")
";"
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)") &
1382 WRITE (unit=file_unit, fmt=
"(A,T31,F18.8)") &
1384 WRITE (unit=file_unit, fmt=
"(A,T31,F18.8)") &
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)") &
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"
1410 r(1:3) =
pbc(particle_set(iatom)%r(1:3), cell)
1413 s(i) =
modulo(s(i), 1.0_dp)
1415 WRITE (unit=file_unit, fmt=trim(f_cif)) &
1416 cif_type_symbol(iatom), cif_label(iatom), 1, s(1:3), 1.0_dp
1422 DEALLOCATE (element_list, count_list, formula_structural, &
1423 formula_sum, cif_label, cif_type_symbol, &
1427 "PRINT%FINAL_STRUCTURE")
1428 IF (output_unit > 0) &
1429 WRITE (unit=output_unit, fmt=
'(/,T2,A)') &
1432 CALL timestop(handle)
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.
subroutine, public real_to_scaled(s, r, cell)
Transform real to scaled cell coordinates. s=h_inv*r.
subroutine, public cell_release(cell)
releases the given cell (see doc/ReferenceCounting.html)
subroutine, public cell_clone(cell_in, cell_out, tag)
Clone cell variable.
subroutine, public get_cell(cell, alpha, beta, gamma, deth, orthorhombic, abc, periodic, h, h_inv, symmetry_id, tag)
Get informations about a simulation cell.
some minimal info about CP2K, including its version and license
character(len=default_string_length), public r_host_name
character(len= *), parameter, public compile_revision
character(len= *), parameter, public cp2k_version
character(len=default_path_length), public r_cwd
character(len=default_string_length), public r_user_name
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,...
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 of the atomic potential types.
Calculation of the incomplete Gamma function F_n(t) for multi-center integrals over Cartesian Gaussia...
Defines the basic variable types.
integer, parameter, public dp
integer, parameter, public default_string_length
integer, parameter, public default_path_length
integer, parameter, public sp
Machine interface based on Fortran 2003 and POSIX.
integer, parameter, public timestamp_length
subroutine, public m_timestamp(timestamp)
Returns a human readable timestamp.
Definition of mathematical constants and functions.
real(kind=dp), parameter, public degree
Collection of simple mathematical functions and subroutines.
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.
elemental integer function, public gcd(a, b)
computes the greatest common divisor of two number
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,...
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:
real(kind=dp), parameter, public massunit
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.
Type defining parameters related to the simulation cell.
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.