61#include "./base/base_uses.f90"
77 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'particle_methods'
101 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
102 INTEGER,
DIMENSION(:),
INTENT(INOUT),
OPTIONAL :: first_sgf, last_sgf, nsgf, nmao
105 INTEGER :: ikind, iparticle, isgf, nparticle, ns
107 cpassert(
ASSOCIATED(particle_set))
109 nparticle =
SIZE(particle_set)
110 IF (
PRESENT(first_sgf))
THEN
111 cpassert(
SIZE(first_sgf) >= nparticle)
113 IF (
PRESENT(last_sgf))
THEN
114 cpassert(
SIZE(last_sgf) >= nparticle)
116 IF (
PRESENT(nsgf))
THEN
117 cpassert(
SIZE(nsgf) >= nparticle)
119 IF (
PRESENT(nmao))
THEN
120 cpassert(
SIZE(nmao) >= nparticle)
123 IF (
PRESENT(first_sgf) .OR.
PRESENT(last_sgf) .OR.
PRESENT(nsgf))
THEN
125 DO iparticle = 1, nparticle
126 CALL get_atomic_kind(particle_set(iparticle)%atomic_kind, kind_number=ikind)
127 IF (
PRESENT(basis))
THEN
128 IF (
ASSOCIATED(basis(ikind)%gto_basis_set))
THEN
136 IF (
PRESENT(nsgf)) nsgf(iparticle) = ns
137 IF (
PRESENT(first_sgf)) first_sgf(iparticle) = isgf + 1
139 IF (
PRESENT(last_sgf)) last_sgf(iparticle) = isgf
143 IF (
PRESENT(first_sgf))
THEN
144 IF (
SIZE(first_sgf) > nparticle) first_sgf(nparticle + 1) = isgf + 1
147 IF (
PRESENT(nmao))
THEN
148 DO iparticle = 1, nparticle
149 CALL get_atomic_kind(particle_set(iparticle)%atomic_kind, kind_number=ikind)
181 content, title, cell, array, unit_conv, &
182 charge_occup, charge_beta, &
183 charge_extended, print_kind)
186 INTEGER :: iunit, output_format
187 CHARACTER(LEN=*) :: content, title
188 TYPE(
cell_type),
OPTIONAL,
POINTER :: cell
189 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN),
OPTIONAL :: array
190 REAL(kind=
dp),
INTENT(IN),
OPTIONAL :: unit_conv
191 LOGICAL,
INTENT(IN),
OPTIONAL :: charge_occup, charge_beta, &
192 charge_extended, print_kind
194 CHARACTER(len=*),
PARAMETER :: routinen =
'write_particle_coordinates'
196 CHARACTER(LEN=120) :: line
197 CHARACTER(LEN=2) :: element_symbol
198 CHARACTER(LEN=4) :: name
199 CHARACTER(LEN=default_string_length) :: atm_name, my_format
200 INTEGER :: handle, iatom, natom
201 LOGICAL :: dummy, my_charge_beta, &
202 my_charge_extended, my_charge_occup, &
204 REAL(kind=
dp) :: angle_alpha, angle_beta, angle_gamma, &
206 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: arr
207 REAL(kind=
dp),
DIMENSION(3) :: abc, angles, f, r, v
208 REAL(kind=
dp),
DIMENSION(3, 3) :: h
209 REAL(kind=
sp),
ALLOCATABLE,
DIMENSION(:) :: x4, y4, z4
214 CALL timeset(routinen, handle)
216 natom =
SIZE(particle_set)
217 IF (
PRESENT(array))
THEN
218 SELECT CASE (trim(content))
219 CASE (
"POS_VEL",
"POS_VEL_FORCE")
220 cpabort(
"Illegal usage")
224 IF (
PRESENT(unit_conv))
THEN
227 SELECT CASE (output_format)
229 my_print_kind = .false.
230 IF (
PRESENT(print_kind)) my_print_kind = print_kind
231 WRITE (iunit,
"(I8)") natom
232 WRITE (iunit,
"(A)") trim(title)
235 element_symbol=element_symbol)
236 IF (len_trim(element_symbol) == 0 .OR. my_print_kind)
THEN
241 atm_name = trim(atm_name)
243 my_format =
"(T2,A2,"
244 atm_name = trim(element_symbol)
246 SELECT CASE (trim(content))
248 IF (
PRESENT(array))
THEN
251 r(:) = particle_set(iatom)%r(:)
253 WRITE (iunit, trim(my_format)//
"1X,3F20.10)") trim(atm_name), r(1:3)*factor
255 IF (
PRESENT(array))
THEN
258 v(:) = particle_set(iatom)%v(:)
260 WRITE (iunit, trim(my_format)//
"1X,3F20.10)") trim(atm_name), v(1:3)*factor
262 IF (
PRESENT(array))
THEN
263 f(:) = array((iatom - 1)*3 + 1:(iatom - 1)*3 + 3)
265 f(:) = particle_set(iatom)%f(:)
267 WRITE (iunit, trim(my_format)//
"1X,3F20.10)") trim(atm_name), f(1:3)*factor
268 CASE (
"FORCE_MIXING_LABELS")
269 IF (
PRESENT(array))
THEN
270 f(:) = array((iatom - 1)*3 + 1:(iatom - 1)*3 + 3)
272 f(:) = particle_set(iatom)%f(:)
274 WRITE (iunit, trim(my_format)//
"1X,3F20.10)") trim(atm_name), f(1:3)*factor
279 SELECT CASE (trim(content))
281 IF (
PRESENT(array))
THEN
284 r(:) = particle_set(iatom)%r(:)
286 WRITE (iunit,
"(3F20.10)") r(1:3)*factor
288 IF (
PRESENT(array))
THEN
291 v(:) = particle_set(iatom)%v(:)
293 WRITE (iunit,
"(3F20.10)") v(1:3)*factor
295 IF (
PRESENT(array))
THEN
296 f(:) = array((iatom - 1)*3 + 1:(iatom - 1)*3 + 3)
298 f(:) = particle_set(iatom)%f(:)
300 WRITE (iunit,
"(3F20.10)") f(1:3)*factor
301 CASE (
"FORCE_MIXING_LABELS")
302 IF (
PRESENT(array))
THEN
303 f(:) = array((iatom - 1)*3 + 1:(iatom - 1)*3 + 3)
305 f(:) = particle_set(iatom)%f(:)
307 WRITE (iunit,
"(3F20.10)") f(1:3)*factor
311 IF (.NOT. (
PRESENT(cell)))
THEN
312 cpabort(
"Cell is not present! Report this bug!")
314 CALL get_cell(cell, alpha=angle_alpha, beta=angle_beta,
gamma=angle_gamma, &
323 CALL cell_clone(cell, cell_dcd, tag=
"CELL_DCD")
324 angles(1) = angle_alpha/
degree
325 angles(2) = angle_beta/
degree
326 angles(3) = angle_gamma/
degree
329 h(1:3, 1:3) = matmul(cell_dcd%hmat(1:3, 1:3), cell%h_inv(1:3, 1:3))
332 ALLOCATE (arr(3, natom))
333 IF (
PRESENT(array))
THEN
334 arr(1:3, 1:natom) = reshape(array, [3, natom])
336 SELECT CASE (trim(content))
339 arr(1:3, iatom) = particle_set(iatom)%r(1:3)
343 arr(1:3, iatom) = particle_set(iatom)%v(1:3)
347 arr(1:3, iatom) = particle_set(iatom)%f(1:3)
350 cpabort(
"Illegal DCD dump type")
357 x4(1:natom) = real(matmul(h(1, 1:3), arr(1:3, 1:natom)), kind=
sp)
358 y4(1:natom) = real(matmul(h(2, 1:3), arr(1:3, 1:natom)), kind=
sp)
359 z4(1:natom) = real(matmul(h(3, 1:3), arr(1:3, 1:natom)), kind=
sp)
361 x4(1:natom) = real(arr(1, 1:natom), kind=
sp)
362 y4(1:natom) = real(arr(2, 1:natom), kind=
sp)
363 z4(1:natom) = real(arr(3, 1:natom), kind=
sp)
365 WRITE (iunit) abc(1)*factor, angle_gamma, abc(2)*factor, &
366 angle_beta, angle_alpha, abc(3)*factor
367 WRITE (iunit) x4*real(factor, kind=
sp)
368 WRITE (iunit) y4*real(factor, kind=
sp)
369 WRITE (iunit) z4*real(factor, kind=
sp)
376 my_charge_occup = .false.
377 IF (
PRESENT(charge_occup)) my_charge_occup = charge_occup
378 my_charge_beta = .false.
379 IF (
PRESENT(charge_beta)) my_charge_beta = charge_beta
380 my_charge_extended = .false.
381 IF (
PRESENT(charge_extended)) my_charge_extended = charge_extended
382 IF (len_trim(title) > 0)
THEN
383 WRITE (unit=iunit, fmt=
"(A6,T11,A)") &
384 "REMARK", trim(title)
386 CALL get_cell(cell, alpha=angle_alpha, beta=angle_beta,
gamma=angle_gamma, abc=abc)
398 WRITE (unit=iunit, fmt=
"(A6,3F9.3,3F7.2)") &
399 "CRYST1", abc(1:3)*factor, angle_alpha, angle_beta, angle_gamma
400 WRITE (unit=line(1:6), fmt=
"(A6)")
"ATOM "
421 element_symbol=element_symbol, name=atm_name, &
422 fist_potential=fist_potential, shell=shell)
423 IF (len_trim(element_symbol) == 0)
THEN
426 name = trim(atm_name)
427 IF (
ASSOCIATED(fist_potential))
THEN
432 IF (
ASSOCIATED(shell))
CALL get_shell(shell=shell, charge=qeff)
433 WRITE (unit=line(1:6), fmt=
"(A6)")
"ATOM "
434 WRITE (unit=line(7:11), fmt=
"(I5)")
modulo(iatom, 100000)
435 WRITE (unit=line(13:16), fmt=
"(A4)") adjustl(name)
438 SELECT CASE (trim(content))
440 IF (
PRESENT(array))
THEN
443 r(:) = particle_set(iatom)%r(:)
445 WRITE (unit=line(31:54), fmt=
"(3F8.3)") r(1:3)*factor
447 cpabort(
"PDB dump only for trajectory available")
449 IF (my_charge_occup)
THEN
450 WRITE (unit=line(55:60), fmt=
"(F6.2)") qeff
452 WRITE (unit=line(55:60), fmt=
"(F6.2)") 0.0_dp
454 IF (my_charge_beta)
THEN
455 WRITE (unit=line(61:66), fmt=
"(F6.2)") qeff
457 WRITE (unit=line(61:66), fmt=
"(F6.2)") 0.0_dp
460 WRITE (unit=line(77:78), fmt=
"(A2)") adjustr(trim(element_symbol))
461 IF (my_charge_extended)
THEN
462 WRITE (unit=line(81:), fmt=
"(SP,F0.8)") qeff
464 WRITE (unit=iunit, fmt=
"(A)") trim(line)
466 WRITE (unit=iunit, fmt=
"(A)")
"END"
468 cpabort(
"Illegal dump type")
471 CALL timestop(handle)
487 REAL(kind=
dp),
DIMENSION(:),
OPTIONAL :: charges
489 CHARACTER(LEN=default_string_length) :: name, unit_str
490 INTEGER :: iatom, ikind, iw, natom
491 REAL(kind=
dp) :: conv, mass, qcore, qeff, qshell
500 "PRINT%ATOMIC_COORDINATES", extension=
".coordLog")
506 WRITE (unit=iw, fmt=
"(/,/,T2,A)") &
507 "MODULE FIST: ATOMIC COORDINATES IN "//trim(unit_str)
508 WRITE (unit=iw, fmt=
"(/,T4,A,T30,A,T44,A,T58,A,T66,A,T77,A)") &
509 "Atom Kind Name",
"X",
"Y",
"Z",
"q(eff)",
"Mass"
510 natom =
SIZE(particle_set)
518 IF (
PRESENT(charges)) qeff = charges(iatom)
519 IF (
ASSOCIATED(shell_kind))
THEN
523 qeff = qcore + qshell
525 WRITE (unit=iw, fmt=
"(T2,I6,1X,I4,1X,A7,3(1X,F13.6),2(1X,F8.4))") &
526 iatom, ikind, name, particle_set(iatom)%r(1:3)*conv, qeff, mass/
massunit
532 "PRINT%ATOMIC_COORDINATES")
549 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
551 CHARACTER(LEN=*),
INTENT(IN) :: label
553 CHARACTER(len=*),
PARAMETER :: routinen =
'write_qs_particle_coordinates'
555 CHARACTER(LEN=2) :: element_symbol
556 CHARACTER(LEN=default_string_length) :: unit_str
557 INTEGER :: handle, iatom, ikind, iw, natom, z
558 REAL(kind=
dp) :: conv, mass, zeff
561 CALL timeset(routinen, handle)
566 "PRINT%ATOMIC_COORDINATES", extension=
".coordLog")
572 WRITE (unit=iw, fmt=
"(/,/,T2,A)") &
573 "MODULE "//trim(label)//
": ATOMIC COORDINATES IN "//trim(unit_str)
574 WRITE (unit=iw, fmt=
"(/,T4,A,T30,A,T44,A,T58,A,T66,A,T77,A)") &
575 "Atom Kind Element",
"X",
"Y",
"Z",
"Z(eff)",
"Mass"
576 natom =
SIZE(particle_set)
580 element_symbol=element_symbol, &
584 WRITE (unit=iw, fmt=
"(T2,I6,1X,I4,1X,A2,1X,I4,3(1X,F13.6),2(1X,F8.4))") &
585 iatom, ikind, element_symbol, z, particle_set(iatom)%r(1:3)*conv, zeff, mass/
massunit
591 "PRINT%ATOMIC_COORDINATES")
593 CALL timestop(handle)
612 CHARACTER(len=*),
PARAMETER :: routinen =
'write_particle_distances'
614 CHARACTER(LEN=default_string_length) :: unit_str
615 INTEGER :: handle, iatom, iw, jatom, natom
616 INTEGER,
DIMENSION(3) :: periodic
618 REAL(kind=
dp) :: conv, dab, dab_abort, dab_min, dab_warn
619 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: distance_matrix
620 REAL(kind=
dp),
DIMENSION(3) :: rab
623 CALL timeset(routinen, handle)
625 cpassert(
ASSOCIATED(particle_set))
626 cpassert(
ASSOCIATED(cell))
627 cpassert(
ASSOCIATED(subsys_section))
632 "PRINT%INTERATOMIC_DISTANCES", extension=
".distLog")
636 CALL section_vals_val_get(subsys_section,
"PRINT%INTERATOMIC_DISTANCES%CHECK_INTERATOMIC_DISTANCES", &
637 r_val=dab_min, explicit=explicit)
641 natom =
SIZE(particle_set)
645 IF (explicit .OR. (iw > 0) .OR. (natom <= 2000))
THEN
646 IF (dab_min > 0.0_dp)
THEN
647 dab_warn = dab_min*conv
648 ELSE IF (dab_min < 0.0_dp)
THEN
649 dab_abort = abs(dab_min)*conv
653 IF ((iw > 0) .OR. (dab_abort > 0.0_dp) .OR. (dab_warn > 0.0_dp))
THEN
654 CALL get_cell(cell=cell, periodic=periodic)
656 ALLOCATE (distance_matrix(natom, natom))
657 distance_matrix(:, :) = 0.0_dp
660 DO jatom = iatom + 1, natom
661 rab(:) =
pbc(particle_set(iatom)%r(:), &
662 particle_set(jatom)%r(:), cell)
663 dab = sqrt(rab(1)*rab(1) + rab(2)*rab(2) + rab(3)*rab(3))*conv
664 IF (dab_abort > 0.0_dp)
THEN
666 IF (dab < dab_abort)
THEN
667 CALL cp_abort(__location__,
"The distance between the atoms "// &
668 trim(adjustl(
cp_to_string(iatom, fmt=
"(I8)")))//
" and "// &
669 trim(adjustl(
cp_to_string(jatom, fmt=
"(I8)")))//
" is only "// &
671 trim(adjustl(unit_str))//
" and thus smaller than the requested threshold of "// &
672 trim(adjustl(
cp_to_string(dab_abort, fmt=
"(F6.3)")))//
" "// &
673 trim(adjustl(unit_str)))
676 IF (dab < dab_warn)
THEN
678 CALL cp_warn(__location__,
"The distance between the atoms "// &
679 trim(adjustl(
cp_to_string(iatom, fmt=
"(I8)")))//
" and "// &
680 trim(adjustl(
cp_to_string(jatom, fmt=
"(I8)")))//
" is only "// &
682 trim(adjustl(unit_str))//
" and thus smaller than the threshold of "// &
683 trim(adjustl(
cp_to_string(dab_warn, fmt=
"(F6.3)")))//
" "// &
684 trim(adjustl(unit_str)))
687 distance_matrix(iatom, jatom) = dab
688 distance_matrix(jatom, iatom) = distance_matrix(iatom, jatom)
694 WRITE (unit=iw, fmt=
"(/,/,T2,A)") &
695 "INTERATOMIC DISTANCES IN "//trim(unit_str)
697 IF (
ALLOCATED(distance_matrix))
DEALLOCATE (distance_matrix)
700 "PRINT%INTERATOMIC_DISTANCES")
703 CALL timestop(handle)
717 REAL(kind=
dp),
DIMENSION(:, :) :: matrix
719 INTEGER,
INTENT(IN) :: iw
720 INTEGER,
INTENT(IN),
OPTIONAL :: el_per_part
721 INTEGER,
DIMENSION(:),
OPTIONAL,
POINTER :: ilist
722 INTEGER,
INTENT(IN),
OPTIONAL :: parts_per_line
724 CHARACTER(LEN=2) :: element_symbol
725 CHARACTER(LEN=default_string_length) :: fmt_string1, fmt_string2
726 INTEGER :: from, i, iatom, icol, jatom, katom, &
727 my_el_per_part, my_parts_per_line, &
729 INTEGER,
DIMENSION(:),
POINTER :: my_list
732 IF (
PRESENT(el_per_part)) my_el_per_part = el_per_part
733 my_parts_per_line = 5
734 IF (
PRESENT(parts_per_line)) my_parts_per_line = max(parts_per_line, 1)
735 WRITE (fmt_string1, fmt=
'(A,I0,A)') &
736 "(/,T2,11X,", my_parts_per_line,
"(4X,I5,4X))"
737 WRITE (fmt_string2, fmt=
'(A,I0,A)') &
738 "(T2,I5,2X,A2,2X,", my_parts_per_line,
"(1X,F12.6))"
739 IF (
PRESENT(ilist))
THEN
742 natom =
SIZE(particle_set)
744 ALLOCATE (my_list(natom))
745 IF (
PRESENT(ilist))
THEN
752 natom = natom*my_el_per_part
753 DO jatom = 1, natom, my_parts_per_line
755 to = min(from + my_parts_per_line - 1, natom)
756 WRITE (unit=iw, fmt=trim(fmt_string1)) (icol, icol=from, to)
758 katom = iatom/my_el_per_part
759 IF (mod(iatom, my_el_per_part) /= 0) katom = katom + 1
760 CALL get_atomic_kind(atomic_kind=particle_set(my_list(katom))%atomic_kind, &
761 element_symbol=element_symbol)
762 WRITE (unit=iw, fmt=trim(fmt_string2)) &
763 iatom, element_symbol, &
764 (matrix(iatom, icol), icol=from, to)
791 CHARACTER(LEN=*),
PARAMETER :: routinen =
'write_structure_data'
793 CHARACTER(LEN=default_string_length) :: string, unit_str
794 INTEGER :: handle, i, i_rep, iw, n, n_rep, n_vals, &
795 natom, new_size, old_size, wrk2(2), &
797 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: work
798 INTEGER,
DIMENSION(:),
POINTER :: atomic_indices, index_list
800 REAL(kind=
dp) :: conv, dab
801 REAL(kind=
dp),
DIMENSION(3) :: r, rab, rbc, rcd, s
805 CALL timeset(routinen, handle)
806 NULLIFY (atomic_indices)
814 basis_section=input_section, &
815 print_key_path=
"PRINT%STRUCTURE_DATA", &
816 extension=
".coordLog")
822 natom =
SIZE(particle_set)
824 subsection_name=
"PRINT%STRUCTURE_DATA")
826 WRITE (unit=iw, fmt=
"(/,T2,A)")
"REQUESTED STRUCTURE DATA"
829 keyword_name=
"POSITION", &
832 WRITE (unit=iw, fmt=
"(/,T3,A,/)") &
833 "Position vectors r(i) of the atoms i in "//trim(unit_str)
837 keyword_name=
"POSITION", &
839 i_vals=atomic_indices)
840 n_vals =
SIZE(atomic_indices)
841 new_size = old_size + n_vals
843 index_list(old_size + 1:new_size) = atomic_indices(1:n_vals)
846 ALLOCATE (work(new_size))
847 CALL sort(index_list, new_size, work)
850 WRITE (unit=string, fmt=
"(A,I0,A)")
"(", index_list(i),
")"
851 IF ((index_list(i) < 1) .OR. (index_list(i) > natom))
THEN
852 WRITE (unit=iw, fmt=
"(T3,A)") &
853 "Invalid atomic index "//trim(string)//
" specified. Print request is ignored."
858 IF (index_list(i) == index_list(i - 1)) cycle
860 WRITE (unit=iw, fmt=
"(T3,A,T20,A,3F13.6)") &
861 "r"//trim(string),
"=",
pbc(particle_set(index_list(i))%r(1:3), cell)*conv
863 DEALLOCATE (index_list)
868 keyword_name=
"POSITION_SCALED", &
871 WRITE (unit=iw, fmt=
"(/,T3,A,/)") &
872 "Position vectors s(i) of the atoms i in scaled coordinates"
876 keyword_name=
"POSITION_SCALED", &
878 i_vals=atomic_indices)
879 n_vals =
SIZE(atomic_indices)
880 new_size = old_size + n_vals
882 index_list(old_size + 1:new_size) = atomic_indices(1:n_vals)
885 ALLOCATE (work(new_size))
886 CALL sort(index_list, new_size, work)
889 WRITE (unit=string, fmt=
"(A,I0,A)")
"(", index_list(i),
")"
890 IF ((index_list(i) < 1) .OR. (index_list(i) > natom))
THEN
891 WRITE (unit=iw, fmt=
"(T3,A)") &
892 "Invalid atomic index "//trim(string)//
" specified. Print request is ignored."
897 IF (index_list(i) == index_list(i - 1)) cycle
899 r(1:3) =
pbc(particle_set(index_list(i))%r(1:3), cell)
901 WRITE (unit=iw, fmt=
"(T3,A,T20,A,3F13.6)") &
902 "s"//trim(string),
"=", s(1:3)
904 DEALLOCATE (index_list)
909 keyword_name=
"DISTANCE", &
912 WRITE (unit=iw, fmt=
"(/,T3,A,/)") &
913 "Distance vector r(i,j) between the atom i and j in "// &
917 keyword_name=
"DISTANCE", &
919 i_vals=atomic_indices)
921 WRITE (unit=string, fmt=
"(A,2(I0,A))") &
922 "(", atomic_indices(1),
",", atomic_indices(2),
")"
923 wrk2 = atomic_indices
925 IF (((wrk2(1) >= 1) .AND. (wrk2(
SIZE(wrk2)) <= natom)) .AND. unique)
THEN
926 rab(:) =
pbc(particle_set(atomic_indices(1))%r(:), &
927 particle_set(atomic_indices(2))%r(:), cell)
928 dab = sqrt(rab(1)*rab(1) + rab(2)*rab(2) + rab(3)*rab(3))
929 WRITE (unit=iw, fmt=
"(T3,A,T20,A,3F13.6,3X,A,F13.6)") &
930 "r"//trim(string),
"=", rab(:)*conv, &
933 WRITE (unit=iw, fmt=
"(T3,A)") &
934 "Invalid atomic indices "//trim(string)//
" specified. Print request is ignored."
941 keyword_name=
"ANGLE", &
944 WRITE (unit=iw, fmt=
"(/,T3,A,/)") &
945 "Angle a(i,j,k) between the atomic distance vectors r(j,i) and "// &
949 keyword_name=
"ANGLE", &
951 i_vals=atomic_indices)
953 WRITE (unit=string, fmt=
"(A,3(I0,A))") &
954 "(", atomic_indices(1),
",", atomic_indices(2),
",", atomic_indices(3),
")"
955 wrk3 = atomic_indices
957 IF (((wrk3(1) >= 1) .AND. (wrk3(
SIZE(wrk3)) <= natom)) .AND. unique)
THEN
958 rab(:) =
pbc(particle_set(atomic_indices(1))%r(:), &
959 particle_set(atomic_indices(2))%r(:), cell)
960 rbc(:) =
pbc(particle_set(atomic_indices(2))%r(:), &
961 particle_set(atomic_indices(3))%r(:), cell)
962 WRITE (unit=iw, fmt=
"(T3,A,T26,A,F9.3)") &
965 WRITE (unit=iw, fmt=
"(T3,A)") &
966 "Invalid atomic indices "//trim(string)//
" specified. Print request is ignored."
973 keyword_name=
"DIHEDRAL_ANGLE", &
976 WRITE (unit=iw, fmt=
"(/,T3,A,/)") &
977 "Dihedral angle d(i,j,k,l) between the planes (i,j,k) and (j,k,l) "// &
981 keyword_name=
"DIHEDRAL_ANGLE", &
983 i_vals=atomic_indices)
985 WRITE (unit=string, fmt=
"(A,4(I0,A))") &
986 "(", atomic_indices(1),
",", atomic_indices(2),
",", &
987 atomic_indices(3),
",", atomic_indices(4),
")"
988 wrk4 = atomic_indices
990 IF (((wrk4(1) >= 1) .AND. (wrk4(
SIZE(wrk4)) <= natom)) .AND. unique)
THEN
991 rab(:) =
pbc(particle_set(atomic_indices(1))%r(:), &
992 particle_set(atomic_indices(2))%r(:), cell)
993 rbc(:) =
pbc(particle_set(atomic_indices(2))%r(:), &
994 particle_set(atomic_indices(3))%r(:), cell)
995 rcd(:) =
pbc(particle_set(atomic_indices(3))%r(:), &
996 particle_set(atomic_indices(4))%r(:), cell)
997 WRITE (unit=iw, fmt=
"(T3,A,T26,A,F9.3)") &
1000 WRITE (unit=iw, fmt=
"(T3,A)") &
1001 "Invalid atomic indices "//trim(string)//
" specified. Print request is ignored."
1007 "PRINT%STRUCTURE_DATA")
1009 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, 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.