60#include "./base/base_uses.f90"
76 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'particle_methods'
100 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
101 INTEGER,
DIMENSION(:),
INTENT(INOUT),
OPTIONAL :: first_sgf, last_sgf, nsgf, nmao
104 INTEGER :: ikind, iparticle, isgf, nparticle, ns
106 cpassert(
ASSOCIATED(particle_set))
108 nparticle =
SIZE(particle_set)
109 IF (
PRESENT(first_sgf))
THEN
110 cpassert(
SIZE(first_sgf) >= nparticle)
112 IF (
PRESENT(last_sgf))
THEN
113 cpassert(
SIZE(last_sgf) >= nparticle)
115 IF (
PRESENT(nsgf))
THEN
116 cpassert(
SIZE(nsgf) >= nparticle)
118 IF (
PRESENT(nmao))
THEN
119 cpassert(
SIZE(nmao) >= nparticle)
122 IF (
PRESENT(first_sgf) .OR.
PRESENT(last_sgf) .OR.
PRESENT(nsgf))
THEN
124 DO iparticle = 1, nparticle
125 CALL get_atomic_kind(particle_set(iparticle)%atomic_kind, kind_number=ikind)
126 IF (
PRESENT(basis))
THEN
127 IF (
ASSOCIATED(basis(ikind)%gto_basis_set))
THEN
135 IF (
PRESENT(nsgf)) nsgf(iparticle) = ns
136 IF (
PRESENT(first_sgf)) first_sgf(iparticle) = isgf + 1
138 IF (
PRESENT(last_sgf)) last_sgf(iparticle) = isgf
142 IF (
PRESENT(first_sgf))
THEN
143 IF (
SIZE(first_sgf) > nparticle) first_sgf(nparticle + 1) = isgf + 1
146 IF (
PRESENT(nmao))
THEN
147 DO iparticle = 1, nparticle
148 CALL get_atomic_kind(particle_set(iparticle)%atomic_kind, kind_number=ikind)
180 content, title, cell, array, unit_conv, &
181 charge_occup, charge_beta, &
182 charge_extended, print_kind)
185 INTEGER :: iunit, output_format
186 CHARACTER(LEN=*) :: content, title
187 TYPE(
cell_type),
OPTIONAL,
POINTER :: cell
188 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN),
OPTIONAL :: array
189 REAL(kind=
dp),
INTENT(IN),
OPTIONAL :: unit_conv
190 LOGICAL,
INTENT(IN),
OPTIONAL :: charge_occup, charge_beta, &
191 charge_extended, print_kind
193 CHARACTER(len=*),
PARAMETER :: routinen =
'write_particle_coordinates'
195 CHARACTER(LEN=120) :: line
196 CHARACTER(LEN=2) :: element_symbol
197 CHARACTER(LEN=4) :: name
198 CHARACTER(LEN=default_string_length) :: atm_name, my_format
199 INTEGER :: handle, iatom, natom
200 LOGICAL :: dummy, my_charge_beta, &
201 my_charge_extended, my_charge_occup, &
203 REAL(kind=
dp) :: angle_alpha, angle_beta, angle_gamma, &
205 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: arr
206 REAL(kind=
dp),
DIMENSION(3) :: abc, angles, f, r, v
207 REAL(kind=
dp),
DIMENSION(3, 3) :: h
208 REAL(kind=
sp),
ALLOCATABLE,
DIMENSION(:) :: x4, y4, z4
213 CALL timeset(routinen, handle)
215 natom =
SIZE(particle_set)
216 IF (
PRESENT(array))
THEN
217 SELECT CASE (trim(content))
218 CASE (
"POS_VEL",
"POS_VEL_FORCE")
219 cpabort(
"Illegal usage")
223 IF (
PRESENT(unit_conv))
THEN
226 SELECT CASE (output_format)
228 my_print_kind = .false.
229 IF (
PRESENT(print_kind)) my_print_kind = print_kind
230 WRITE (iunit,
"(I8)") natom
231 WRITE (iunit,
"(A)") trim(title)
234 element_symbol=element_symbol)
235 IF (len_trim(element_symbol) == 0 .OR. my_print_kind)
THEN
240 atm_name = trim(atm_name)
242 my_format =
"(T2,A2,"
243 atm_name = trim(element_symbol)
245 SELECT CASE (trim(content))
247 IF (
PRESENT(array))
THEN
250 r(:) = particle_set(iatom)%r(:)
252 WRITE (iunit, trim(my_format)//
"1X,3F20.10)") trim(atm_name), r(1:3)*factor
254 IF (
PRESENT(array))
THEN
257 v(:) = particle_set(iatom)%v(:)
259 WRITE (iunit, trim(my_format)//
"1X,3F20.10)") trim(atm_name), v(1:3)*factor
261 IF (
PRESENT(array))
THEN
262 f(:) = array((iatom - 1)*3 + 1:(iatom - 1)*3 + 3)
264 f(:) = particle_set(iatom)%f(:)
266 WRITE (iunit, trim(my_format)//
"1X,3F20.10)") trim(atm_name), f(1:3)*factor
267 CASE (
"FORCE_MIXING_LABELS")
268 IF (
PRESENT(array))
THEN
269 f(:) = array((iatom - 1)*3 + 1:(iatom - 1)*3 + 3)
271 f(:) = particle_set(iatom)%f(:)
273 WRITE (iunit, trim(my_format)//
"1X,3F20.10)") trim(atm_name), f(1:3)*factor
278 SELECT CASE (trim(content))
280 IF (
PRESENT(array))
THEN
283 r(:) = particle_set(iatom)%r(:)
285 WRITE (iunit,
"(3F20.10)") r(1:3)*factor
287 IF (
PRESENT(array))
THEN
290 v(:) = particle_set(iatom)%v(:)
292 WRITE (iunit,
"(3F20.10)") v(1:3)*factor
294 IF (
PRESENT(array))
THEN
295 f(:) = array((iatom - 1)*3 + 1:(iatom - 1)*3 + 3)
297 f(:) = particle_set(iatom)%f(:)
299 WRITE (iunit,
"(3F20.10)") f(1:3)*factor
300 CASE (
"FORCE_MIXING_LABELS")
301 IF (
PRESENT(array))
THEN
302 f(:) = array((iatom - 1)*3 + 1:(iatom - 1)*3 + 3)
304 f(:) = particle_set(iatom)%f(:)
306 WRITE (iunit,
"(3F20.10)") f(1:3)*factor
310 IF (.NOT. (
PRESENT(cell)))
THEN
311 cpabort(
"Cell is not present! Report this bug!")
313 CALL get_cell(cell, alpha=angle_alpha, beta=angle_beta,
gamma=angle_gamma, &
322 CALL cell_clone(cell, cell_dcd, tag=
"CELL_DCD")
323 angles(1) = angle_alpha/
degree
324 angles(2) = angle_beta/
degree
325 angles(3) = angle_gamma/
degree
328 h(1:3, 1:3) = matmul(cell_dcd%hmat(1:3, 1:3), cell%h_inv(1:3, 1:3))
331 ALLOCATE (arr(3, natom))
332 IF (
PRESENT(array))
THEN
333 arr(1:3, 1:natom) = reshape(array, (/3, natom/))
335 SELECT CASE (trim(content))
338 arr(1:3, iatom) = particle_set(iatom)%r(1:3)
342 arr(1:3, iatom) = particle_set(iatom)%v(1:3)
346 arr(1:3, iatom) = particle_set(iatom)%f(1:3)
349 cpabort(
"Illegal DCD dump type")
356 x4(1:natom) = real(matmul(h(1, 1:3), arr(1:3, 1:natom)), kind=
sp)
357 y4(1:natom) = real(matmul(h(2, 1:3), arr(1:3, 1:natom)), kind=
sp)
358 z4(1:natom) = real(matmul(h(3, 1:3), arr(1:3, 1:natom)), kind=
sp)
360 x4(1:natom) = real(arr(1, 1:natom), kind=
sp)
361 y4(1:natom) = real(arr(2, 1:natom), kind=
sp)
362 z4(1:natom) = real(arr(3, 1:natom), kind=
sp)
364 WRITE (iunit) abc(1)*factor, angle_gamma, abc(2)*factor, &
365 angle_beta, angle_alpha, abc(3)*factor
366 WRITE (iunit) x4*real(factor, kind=
sp)
367 WRITE (iunit) y4*real(factor, kind=
sp)
368 WRITE (iunit) z4*real(factor, kind=
sp)
375 my_charge_occup = .false.
376 IF (
PRESENT(charge_occup)) my_charge_occup = charge_occup
377 my_charge_beta = .false.
378 IF (
PRESENT(charge_beta)) my_charge_beta = charge_beta
379 my_charge_extended = .false.
380 IF (
PRESENT(charge_extended)) my_charge_extended = charge_extended
381 IF (len_trim(title) > 0)
THEN
382 WRITE (unit=iunit, fmt=
"(A6,T11,A)") &
383 "REMARK", trim(title)
385 CALL get_cell(cell, alpha=angle_alpha, beta=angle_beta,
gamma=angle_gamma, abc=abc)
397 WRITE (unit=iunit, fmt=
"(A6,3F9.3,3F7.2)") &
398 "CRYST1", abc(1:3)*factor, angle_alpha, angle_beta, angle_gamma
399 WRITE (unit=line(1:6), fmt=
"(A6)")
"ATOM "
420 element_symbol=element_symbol, name=atm_name, &
421 fist_potential=fist_potential, shell=shell)
422 IF (len_trim(element_symbol) == 0)
THEN
425 name = trim(atm_name)
426 IF (
ASSOCIATED(fist_potential))
THEN
431 IF (
ASSOCIATED(shell))
CALL get_shell(shell=shell, charge=qeff)
432 WRITE (unit=line(1:6), fmt=
"(A6)")
"ATOM "
433 WRITE (unit=line(7:11), fmt=
"(I5)")
modulo(iatom, 100000)
434 WRITE (unit=line(13:16), fmt=
"(A4)") adjustl(name)
437 SELECT CASE (trim(content))
439 IF (
PRESENT(array))
THEN
442 r(:) = particle_set(iatom)%r(:)
444 WRITE (unit=line(31:54), fmt=
"(3F8.3)") r(1:3)*factor
446 cpabort(
"PDB dump only for trajectory available")
448 IF (my_charge_occup)
THEN
449 WRITE (unit=line(55:60), fmt=
"(F6.2)") qeff
451 WRITE (unit=line(55:60), fmt=
"(F6.2)") 0.0_dp
453 IF (my_charge_beta)
THEN
454 WRITE (unit=line(61:66), fmt=
"(F6.2)") qeff
456 WRITE (unit=line(61:66), fmt=
"(F6.2)") 0.0_dp
459 WRITE (unit=line(77:78), fmt=
"(A2)") adjustr(trim(element_symbol))
460 IF (my_charge_extended)
THEN
461 WRITE (unit=line(81:), fmt=
"(SP,F0.8)") qeff
463 WRITE (unit=iunit, fmt=
"(A)") trim(line)
465 WRITE (unit=iunit, fmt=
"(A)")
"END"
467 cpabort(
"Illegal dump type")
470 CALL timestop(handle)
486 REAL(kind=
dp),
DIMENSION(:),
OPTIONAL :: charges
488 CHARACTER(LEN=default_string_length) :: name, unit_str
489 INTEGER :: iatom, ikind, iw, natom
490 REAL(kind=
dp) :: conv, mass, qcore, qeff, qshell
499 "PRINT%ATOMIC_COORDINATES", extension=
".coordLog")
505 WRITE (unit=iw, fmt=
"(/,/,T2,A)") &
506 "MODULE FIST: ATOMIC COORDINATES IN "//trim(unit_str)
507 WRITE (unit=iw, fmt=
"(/,T4,A,T30,A,T44,A,T58,A,T66,A,T77,A)") &
508 "Atom Kind Name",
"X",
"Y",
"Z",
"q(eff)",
"Mass"
509 natom =
SIZE(particle_set)
517 IF (
PRESENT(charges)) qeff = charges(iatom)
518 IF (
ASSOCIATED(shell_kind))
THEN
522 qeff = qcore + qshell
524 WRITE (unit=iw, fmt=
"(T2,I6,1X,I4,1X,A7,3(1X,F13.6),2(1X,F8.4))") &
525 iatom, ikind, name, particle_set(iatom)%r(1:3)*conv, qeff, mass/
massunit
531 "PRINT%ATOMIC_COORDINATES")
548 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
550 CHARACTER(LEN=*),
INTENT(IN) :: label
552 CHARACTER(len=*),
PARAMETER :: routinen =
'write_qs_particle_coordinates'
554 CHARACTER(LEN=2) :: element_symbol
555 CHARACTER(LEN=default_string_length) :: unit_str
556 INTEGER :: handle, iatom, ikind, iw, natom, z
557 REAL(kind=
dp) :: conv, mass, zeff
560 CALL timeset(routinen, handle)
565 "PRINT%ATOMIC_COORDINATES", extension=
".coordLog")
571 WRITE (unit=iw, fmt=
"(/,/,T2,A)") &
572 "MODULE "//trim(label)//
": ATOMIC COORDINATES IN "//trim(unit_str)
573 WRITE (unit=iw, fmt=
"(/,T4,A,T30,A,T44,A,T58,A,T66,A,T77,A)") &
574 "Atom Kind Element",
"X",
"Y",
"Z",
"Z(eff)",
"Mass"
575 natom =
SIZE(particle_set)
579 element_symbol=element_symbol, &
583 WRITE (unit=iw, fmt=
"(T2,I6,1X,I4,1X,A2,1X,I4,3(1X,F13.6),2(1X,F8.4))") &
584 iatom, ikind, element_symbol, z, particle_set(iatom)%r(1:3)*conv, zeff, mass/
massunit
590 "PRINT%ATOMIC_COORDINATES")
592 CALL timestop(handle)
611 CHARACTER(len=*),
PARAMETER :: routinen =
'write_particle_distances'
613 CHARACTER(LEN=default_string_length) :: unit_str
614 INTEGER :: handle, iatom, iw, jatom, natom
615 INTEGER,
DIMENSION(3) :: periodic
617 REAL(kind=
dp) :: conv, dab, dab_abort, dab_min, dab_warn
618 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: distance_matrix
619 REAL(kind=
dp),
DIMENSION(3) :: rab
622 CALL timeset(routinen, handle)
624 cpassert(
ASSOCIATED(particle_set))
625 cpassert(
ASSOCIATED(cell))
626 cpassert(
ASSOCIATED(subsys_section))
631 "PRINT%INTERATOMIC_DISTANCES", extension=
".distLog")
635 CALL section_vals_val_get(subsys_section,
"PRINT%INTERATOMIC_DISTANCES%CHECK_INTERATOMIC_DISTANCES", &
636 r_val=dab_min, explicit=explicit)
640 natom =
SIZE(particle_set)
644 IF (explicit .OR. (iw > 0) .OR. (natom <= 2000))
THEN
645 IF (dab_min > 0.0_dp)
THEN
646 dab_warn = dab_min*conv
647 ELSE IF (dab_min < 0.0_dp)
THEN
648 dab_abort = abs(dab_min)*conv
652 IF ((iw > 0) .OR. (dab_abort > 0.0_dp) .OR. (dab_warn > 0.0_dp))
THEN
653 CALL get_cell(cell=cell, periodic=periodic)
655 ALLOCATE (distance_matrix(natom, natom))
656 distance_matrix(:, :) = 0.0_dp
659 DO jatom = iatom + 1, natom
660 rab(:) =
pbc(particle_set(iatom)%r(:), &
661 particle_set(jatom)%r(:), cell)
662 dab = sqrt(rab(1)*rab(1) + rab(2)*rab(2) + rab(3)*rab(3))*conv
663 IF (dab_abort > 0.0_dp)
THEN
665 IF (dab < dab_abort)
THEN
666 CALL cp_abort(__location__,
"The distance between the atoms "// &
667 trim(adjustl(
cp_to_string(iatom, fmt=
"(I8)")))//
" and "// &
668 trim(adjustl(
cp_to_string(jatom, fmt=
"(I8)")))//
" is only "// &
670 trim(adjustl(unit_str))//
" and thus smaller than the requested threshold of "// &
671 trim(adjustl(
cp_to_string(dab_abort, fmt=
"(F6.3)")))//
" "// &
672 trim(adjustl(unit_str)))
675 IF (dab < dab_warn)
THEN
677 CALL cp_warn(__location__,
"The distance between the atoms "// &
678 trim(adjustl(
cp_to_string(iatom, fmt=
"(I8)")))//
" and "// &
679 trim(adjustl(
cp_to_string(jatom, fmt=
"(I8)")))//
" is only "// &
681 trim(adjustl(unit_str))//
" and thus smaller than the threshold of "// &
682 trim(adjustl(
cp_to_string(dab_warn, fmt=
"(F6.3)")))//
" "// &
683 trim(adjustl(unit_str)))
686 distance_matrix(iatom, jatom) = dab
687 distance_matrix(jatom, iatom) = distance_matrix(iatom, jatom)
693 WRITE (unit=iw, fmt=
"(/,/,T2,A)") &
694 "INTERATOMIC DISTANCES IN "//trim(unit_str)
696 IF (
ALLOCATED(distance_matrix))
DEALLOCATE (distance_matrix)
699 "PRINT%INTERATOMIC_DISTANCES")
702 CALL timestop(handle)
716 REAL(kind=
dp),
DIMENSION(:, :) :: matrix
718 INTEGER,
INTENT(IN) :: iw
719 INTEGER,
INTENT(IN),
OPTIONAL :: el_per_part
720 INTEGER,
DIMENSION(:),
OPTIONAL,
POINTER :: ilist
721 INTEGER,
INTENT(IN),
OPTIONAL :: parts_per_line
723 CHARACTER(LEN=2) :: element_symbol
724 CHARACTER(LEN=default_string_length) :: fmt_string1, fmt_string2
725 INTEGER :: from, i, iatom, icol, jatom, katom, &
726 my_el_per_part, my_parts_per_line, &
728 INTEGER,
DIMENSION(:),
POINTER :: my_list
731 IF (
PRESENT(el_per_part)) my_el_per_part = el_per_part
732 my_parts_per_line = 5
733 IF (
PRESENT(parts_per_line)) my_parts_per_line = max(parts_per_line, 1)
734 WRITE (fmt_string1, fmt=
'(A,I0,A)') &
735 "(/,T2,11X,", my_parts_per_line,
"(4X,I5,4X))"
736 WRITE (fmt_string2, fmt=
'(A,I0,A)') &
737 "(T2,I5,2X,A2,2X,", my_parts_per_line,
"(1X,F12.6))"
738 IF (
PRESENT(ilist))
THEN
741 natom =
SIZE(particle_set)
743 ALLOCATE (my_list(natom))
744 IF (
PRESENT(ilist))
THEN
751 natom = natom*my_el_per_part
752 DO jatom = 1, natom, my_parts_per_line
754 to = min(from + my_parts_per_line - 1, natom)
755 WRITE (unit=iw, fmt=trim(fmt_string1)) (icol, icol=from, to)
757 katom = iatom/my_el_per_part
758 IF (mod(iatom, my_el_per_part) /= 0) katom = katom + 1
759 CALL get_atomic_kind(atomic_kind=particle_set(my_list(katom))%atomic_kind, &
760 element_symbol=element_symbol)
761 WRITE (unit=iw, fmt=trim(fmt_string2)) &
762 iatom, element_symbol, &
763 (matrix(iatom, icol), icol=from, to)
790 CHARACTER(LEN=*),
PARAMETER :: routinen =
'write_structure_data'
792 CHARACTER(LEN=default_string_length) :: string, unit_str
793 INTEGER :: handle, i, i_rep, iw, n, n_rep, n_vals, &
794 natom, new_size, old_size, wrk2(2), &
796 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: work
797 INTEGER,
DIMENSION(:),
POINTER :: atomic_indices, index_list
799 REAL(kind=
dp) :: conv, dab
800 REAL(kind=
dp),
DIMENSION(3) :: r, rab, rbc, rcd, s
804 CALL timeset(routinen, handle)
805 NULLIFY (atomic_indices)
813 basis_section=input_section, &
814 print_key_path=
"PRINT%STRUCTURE_DATA", &
815 extension=
".coordLog")
821 natom =
SIZE(particle_set)
823 subsection_name=
"PRINT%STRUCTURE_DATA")
825 WRITE (unit=iw, fmt=
"(/,T2,A)")
"REQUESTED STRUCTURE DATA"
828 keyword_name=
"POSITION", &
831 WRITE (unit=iw, fmt=
"(/,T3,A,/)") &
832 "Position vectors r(i) of the atoms i in "//trim(unit_str)
836 keyword_name=
"POSITION", &
838 i_vals=atomic_indices)
839 n_vals =
SIZE(atomic_indices)
840 new_size = old_size + n_vals
842 index_list(old_size + 1:new_size) = atomic_indices(1:n_vals)
845 ALLOCATE (work(new_size))
846 CALL sort(index_list, new_size, work)
849 WRITE (unit=string, fmt=
"(A,I0,A)")
"(", index_list(i),
")"
850 IF ((index_list(i) < 1) .OR. (index_list(i) > natom))
THEN
851 WRITE (unit=iw, fmt=
"(T3,A)") &
852 "Invalid atomic index "//trim(string)//
" specified. Print request is ignored."
857 IF (index_list(i) == index_list(i - 1)) cycle
859 WRITE (unit=iw, fmt=
"(T3,A,T20,A,3F13.6)") &
860 "r"//trim(string),
"=",
pbc(particle_set(index_list(i))%r(1:3), cell)*conv
862 DEALLOCATE (index_list)
867 keyword_name=
"POSITION_SCALED", &
870 WRITE (unit=iw, fmt=
"(/,T3,A,/)") &
871 "Position vectors s(i) of the atoms i in scaled coordinates"
875 keyword_name=
"POSITION_SCALED", &
877 i_vals=atomic_indices)
878 n_vals =
SIZE(atomic_indices)
879 new_size = old_size + n_vals
881 index_list(old_size + 1:new_size) = atomic_indices(1:n_vals)
884 ALLOCATE (work(new_size))
885 CALL sort(index_list, new_size, work)
888 WRITE (unit=string, fmt=
"(A,I0,A)")
"(", index_list(i),
")"
889 IF ((index_list(i) < 1) .OR. (index_list(i) > natom))
THEN
890 WRITE (unit=iw, fmt=
"(T3,A)") &
891 "Invalid atomic index "//trim(string)//
" specified. Print request is ignored."
896 IF (index_list(i) == index_list(i - 1)) cycle
898 r(1:3) =
pbc(particle_set(index_list(i))%r(1:3), cell)
900 WRITE (unit=iw, fmt=
"(T3,A,T20,A,3F13.6)") &
901 "s"//trim(string),
"=", s(1:3)
903 DEALLOCATE (index_list)
908 keyword_name=
"DISTANCE", &
911 WRITE (unit=iw, fmt=
"(/,T3,A,/)") &
912 "Distance vector r(i,j) between the atom i and j in "// &
916 keyword_name=
"DISTANCE", &
918 i_vals=atomic_indices)
920 WRITE (unit=string, fmt=
"(A,2(I0,A))") &
921 "(", atomic_indices(1),
",", atomic_indices(2),
")"
922 wrk2 = atomic_indices
924 IF (((wrk2(1) >= 1) .AND. (wrk2(
SIZE(wrk2)) <= natom)) .AND. unique)
THEN
925 rab(:) =
pbc(particle_set(atomic_indices(1))%r(:), &
926 particle_set(atomic_indices(2))%r(:), cell)
927 dab = sqrt(rab(1)*rab(1) + rab(2)*rab(2) + rab(3)*rab(3))
928 WRITE (unit=iw, fmt=
"(T3,A,T20,A,3F13.6,3X,A,F13.6)") &
929 "r"//trim(string),
"=", rab(:)*conv, &
932 WRITE (unit=iw, fmt=
"(T3,A)") &
933 "Invalid atomic indices "//trim(string)//
" specified. Print request is ignored."
940 keyword_name=
"ANGLE", &
943 WRITE (unit=iw, fmt=
"(/,T3,A,/)") &
944 "Angle a(i,j,k) between the atomic distance vectors r(j,i) and "// &
948 keyword_name=
"ANGLE", &
950 i_vals=atomic_indices)
952 WRITE (unit=string, fmt=
"(A,3(I0,A))") &
953 "(", atomic_indices(1),
",", atomic_indices(2),
",", atomic_indices(3),
")"
954 wrk3 = atomic_indices
956 IF (((wrk3(1) >= 1) .AND. (wrk3(
SIZE(wrk3)) <= natom)) .AND. unique)
THEN
957 rab(:) =
pbc(particle_set(atomic_indices(1))%r(:), &
958 particle_set(atomic_indices(2))%r(:), cell)
959 rbc(:) =
pbc(particle_set(atomic_indices(2))%r(:), &
960 particle_set(atomic_indices(3))%r(:), cell)
961 WRITE (unit=iw, fmt=
"(T3,A,T26,A,F9.3)") &
964 WRITE (unit=iw, fmt=
"(T3,A)") &
965 "Invalid atomic indices "//trim(string)//
" specified. Print request is ignored."
972 keyword_name=
"DIHEDRAL_ANGLE", &
975 WRITE (unit=iw, fmt=
"(/,T3,A,/)") &
976 "Dihedral angle d(i,j,k,l) between the planes (i,j,k) and (j,k,l) "// &
980 keyword_name=
"DIHEDRAL_ANGLE", &
982 i_vals=atomic_indices)
984 WRITE (unit=string, fmt=
"(A,4(I0,A))") &
985 "(", atomic_indices(1),
",", atomic_indices(2),
",", &
986 atomic_indices(3),
",", atomic_indices(4),
")"
987 wrk4 = atomic_indices
989 IF (((wrk4(1) >= 1) .AND. (wrk4(
SIZE(wrk4)) <= natom)) .AND. unique)
THEN
990 rab(:) =
pbc(particle_set(atomic_indices(1))%r(:), &
991 particle_set(atomic_indices(2))%r(:), cell)
992 rbc(:) =
pbc(particle_set(atomic_indices(2))%r(:), &
993 particle_set(atomic_indices(3))%r(:), cell)
994 rcd(:) =
pbc(particle_set(atomic_indices(3))%r(:), &
995 particle_set(atomic_indices(4))%r(:), cell)
996 WRITE (unit=iw, fmt=
"(T3,A,T26,A,F9.3)") &
999 WRITE (unit=iw, fmt=
"(T3,A)") &
1000 "Invalid atomic indices "//trim(string)//
" specified. Print request is ignored."
1006 "PRINT%STRUCTURE_DATA")
1008 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)
...
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.
various routines to log and control the output. The idea is that decisions about where to log should ...
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer function, public cp_print_key_unit_nr(logger, basis_section, print_key_path, extension, middle_name, local, log_filename, ignore_should_output, file_form, file_position, file_action, file_status, do_backup, on_file, is_new_file, mpi_io, fout)
...
subroutine, public cp_print_key_finished_output(unit_nr, logger, basis_section, print_key_path, local, ignore_should_output, on_file, mpi_io)
should be called after you finish working with a unit obtained with cp_print_key_unit_nr,...
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 sp
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.
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_particle_matrix(matrix, particle_set, iw, el_per_part, ilist, parts_per_line)
...
subroutine, public write_structure_data(particle_set, cell, input_section)
Write structure data requested by a separate structure data input section to the output unit....
subroutine, public write_particle_distances(particle_set, cell, subsys_section)
Write the matrix of the particle distances to the output unit.
subroutine, public get_particle_set(particle_set, qs_kind_set, first_sgf, last_sgf, nsgf, nmao, basis)
Get the components of a particle set.
subroutine, public write_particle_coordinates(particle_set, iunit, output_format, content, title, cell, array, unit_conv, charge_occup, charge_beta, charge_extended, print_kind)
Should be able to write a few formats e.g. xmol, and some binary format (dcd) some format can be used...
Define the data structure for the particle information.
pure real(kind=dp) function, dimension(3), public get_particle_pos_or_vel(iatom, particle_set, vector)
Return the atomic position or velocity of atom iatom in x from a packed vector even if core-shell par...
Definition of physical constants:
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, 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, 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.