111#include "./base/base_uses.f90"
117 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'qs_tddfpt2_properties'
120 INTEGER,
PARAMETER,
PRIVATE :: nderivs = 3
121 INTEGER,
PARAMETER,
PRIVATE :: maxspins = 2
154 TYPE(
cp_fm_type),
ALLOCATABLE,
DIMENSION(:, :), &
155 INTENT(inout) :: dipole_op_mos_occ
161 CHARACTER(LEN=*),
PARAMETER :: routinen =
'tddfpt_dipole_operator'
163 INTEGER :: handle, i_cos_sin, icol, ideriv, irow, &
164 ispin, jderiv, nao, ncols_local, &
165 ndim_periodic, nrows_local, nspins
166 INTEGER,
DIMENSION(:),
POINTER :: col_indices, row_indices
167 INTEGER,
DIMENSION(maxspins) :: nmo_occ, nmo_virt
168 REAL(kind=
dp) :: eval_occ
169 REAL(kind=
dp),
CONTIGUOUS,
DIMENSION(:, :), &
170 POINTER :: local_data_ediff, local_data_wfm
171 REAL(kind=
dp),
DIMENSION(3) :: kvec, reference_point
174 TYPE(
cp_cfm_type),
ALLOCATABLE,
DIMENSION(:) :: gamma_00, gamma_inv_00
176 TYPE(
cp_fm_type) :: ediff_inv, wfm_ao_ao, wfm_mo_virt_mo_occ
177 TYPE(
cp_fm_type),
ALLOCATABLE,
DIMENSION(:) :: s_mos_virt
178 TYPE(
cp_fm_type),
ALLOCATABLE,
DIMENSION(:, :) :: dberry_mos_occ, gamma_real_imag, opvec
179 TYPE(
dbcsr_p_type),
DIMENSION(:),
POINTER :: berry_cossin_xyz, matrix_s, rrc_xyz, scrm
187 CALL timeset(routinen, handle)
189 NULLIFY (blacs_env, cell, matrix_s, pw_env)
190 CALL get_qs_env(qs_env, blacs_env=blacs_env, cell=cell, matrix_s=matrix_s, pw_env=pw_env)
192 nspins =
SIZE(gs_mos)
195 nmo_occ(ispin) =
SIZE(gs_mos(ispin)%evals_occ)
196 nmo_virt(ispin) =
SIZE(gs_mos(ispin)%evals_virt)
200 ALLOCATE (dipole_op_mos_occ(nderivs, nspins))
202 CALL cp_fm_get_info(gs_mos(ispin)%mos_occ, matrix_struct=fm_struct)
204 DO ideriv = 1, nderivs
205 CALL cp_fm_create(dipole_op_mos_occ(ideriv, ispin), fm_struct)
210 ALLOCATE (s_mos_virt(nspins))
212 CALL cp_fm_get_info(gs_mos(ispin)%mos_virt, matrix_struct=fm_struct)
215 gs_mos(ispin)%mos_virt, &
217 ncol=nmo_virt(ispin), alpha=1.0_dp, beta=0.0_dp)
221 CALL pw_env_get(pw_env, poisson_env=poisson_env)
222 ndim_periodic = count(poisson_env%parameters%periodic == 1)
225 IF (tddfpt_control%dipole_form == 0)
THEN
226 CALL get_qs_env(qs_env, dft_control=dft_control)
227 IF (dft_control%qs_control%xtb)
THEN
228 IF (ndim_periodic == 0)
THEN
238 SELECT CASE (tddfpt_control%dipole_form)
240 IF (ndim_periodic /= 3)
THEN
241 CALL cp_warn(__location__, &
242 "Fully periodic Poisson solver (PERIODIC xyz) "// &
243 "or a large supercell in non-periodic directions is needed "// &
244 "for oscillator strengths based on the Berry phase formula")
247 NULLIFY (berry_cossin_xyz)
255 CALL dbcsr_copy(berry_cossin_xyz(i_cos_sin)%matrix, matrix_s(1)%matrix)
259 ALLOCATE (gamma_00(nspins), gamma_inv_00(nspins), gamma_real_imag(2, nspins), opvec(2, nspins))
260 ALLOCATE (dberry_mos_occ(nderivs, nspins))
264 ncol_global=nmo_occ(ispin), context=blacs_env)
270 CALL cp_fm_create(gamma_real_imag(i_cos_sin, ispin), fm_struct)
275 CALL cp_fm_get_info(gs_mos(ispin)%mos_occ, matrix_struct=fm_struct)
281 DO ideriv = 1, nderivs
282 CALL cp_fm_create(dberry_mos_occ(ideriv, ispin), fm_struct)
287 DO ideriv = 1, nderivs
288 kvec(:) =
twopi*cell%h_inv(ideriv, :)
290 CALL dbcsr_set(berry_cossin_xyz(i_cos_sin)%matrix, 0.0_dp)
293 berry_cossin_xyz(2)%matrix, kvec)
300 gs_mos(ispin)%mos_occ, &
301 opvec(i_cos_sin, ispin), &
302 ncol=nmo_occ(ispin), alpha=1.0_dp, beta=0.0_dp)
305 CALL parallel_gemm(
'T',
'N', nmo_occ(ispin), nmo_occ(ispin), nao, &
306 1.0_dp, gs_mos(ispin)%mos_occ, opvec(1, ispin), &
307 0.0_dp, gamma_real_imag(1, ispin))
309 CALL parallel_gemm(
'T',
'N', nmo_occ(ispin), nmo_occ(ispin), nao, &
310 -1.0_dp, gs_mos(ispin)%mos_occ, opvec(2, ispin), &
311 0.0_dp, gamma_real_imag(2, ispin))
314 msourcei=gamma_real_imag(2, ispin), &
315 mtarget=gamma_00(ispin))
322 mtargetr=gamma_real_imag(1, ispin), &
323 mtargeti=gamma_real_imag(2, ispin))
326 CALL parallel_gemm(
"N",
"N", nao, nmo_occ(ispin), nmo_occ(ispin), &
327 1.0_dp, opvec(1, ispin), gamma_real_imag(2, ispin), &
328 0.0_dp, dipole_op_mos_occ(1, ispin))
329 CALL parallel_gemm(
"N",
"N", nao, nmo_occ(ispin), nmo_occ(ispin), &
330 -1.0_dp, opvec(2, ispin), gamma_real_imag(1, ispin), &
331 1.0_dp, dipole_op_mos_occ(1, ispin))
333 DO jderiv = 1, nderivs
335 cell%hmat(jderiv, ideriv), dipole_op_mos_occ(1, ispin))
343 DO ispin = nspins, 1, -1
347 DEALLOCATE (gamma_00, gamma_inv_00)
360 DO ideriv = 1, nderivs
361 CALL cp_fm_to_fm(dberry_mos_occ(ideriv, ispin), dipole_op_mos_occ(ideriv, ispin))
369 IF (ndim_periodic /= 0)
THEN
370 CALL cp_warn(__location__, &
371 "Non-periodic Poisson solver (PERIODIC none) "// &
372 "or a large supercell approach is needed "// &
373 "for oscillator strengths based on the length operator")
380 DO ideriv = 1, nderivs
382 CALL dbcsr_copy(rrc_xyz(ideriv)%matrix, matrix_s(1)%matrix)
386 reference=tddfpt_control%dipole_reference, &
387 ref_point=tddfpt_control%dipole_ref_point)
389 CALL rrc_xyz_ao(op=rrc_xyz, qs_env=qs_env, rc=reference_point, order=1, &
390 minimum_image=.false., soft=.false.)
394 DO ideriv = 1, nderivs
396 gs_mos(ispin)%mos_occ, &
397 dipole_op_mos_occ(ideriv, ispin), &
398 ncol=nmo_occ(ispin), alpha=1.0_dp, beta=0.0_dp)
407 CALL get_qs_env(qs_env, ks_env=ks_env, sab_orb=sab_orb)
410 basis_type_a=
"ORB", basis_type_b=
"ORB", &
414 DO ideriv = 1, nderivs
416 gs_mos(ispin)%mos_occ, &
417 dipole_op_mos_occ(ideriv, ispin), &
418 ncol=nmo_occ(ispin), alpha=1.0_dp, beta=0.0_dp)
427 CALL get_qs_env(qs_env, ks_env=ks_env, sab_orb=sab_orb)
430 basis_type_a=
"ORB", basis_type_b=
"ORB", &
436 ncol_global=nmo_occ(ispin), context=blacs_env)
441 CALL cp_fm_get_info(ediff_inv, nrow_local=nrows_local, ncol_local=ncols_local, &
442 row_indices=row_indices, col_indices=col_indices, local_data=local_data_ediff)
443 CALL cp_fm_get_info(wfm_mo_virt_mo_occ, local_data=local_data_wfm)
448 DO icol = 1, ncols_local
450 eval_occ = gs_mos(ispin)%evals_occ(col_indices(icol))
452 DO irow = 1, nrows_local
455 local_data_ediff(irow, icol) = 1.0_dp/(gs_mos(ispin)%evals_virt(row_indices(irow)) - eval_occ)
460 DO ideriv = 1, nderivs
462 gs_mos(ispin)%mos_occ, &
463 dipole_op_mos_occ(ideriv, ispin), &
464 ncol=nmo_occ(ispin), alpha=1.0_dp, beta=0.0_dp)
466 CALL parallel_gemm(
'T',
'N', nmo_virt(ispin), nmo_occ(ispin), nao, &
467 1.0_dp, gs_mos(ispin)%mos_virt, dipole_op_mos_occ(ideriv, ispin), &
468 0.0_dp, wfm_mo_virt_mo_occ)
477 DO icol = 1, ncols_local
478 DO irow = 1, nrows_local
479 local_data_wfm(irow, icol) = local_data_wfm(irow, icol)*local_data_ediff(irow, icol)
484 CALL parallel_gemm(
'N',
'N', nao, nmo_occ(ispin), nmo_virt(ispin), &
485 1.0_dp, s_mos_virt(ispin), wfm_mo_virt_mo_occ, &
486 0.0_dp, dipole_op_mos_occ(ideriv, ispin))
495 cpabort(
"Unimplemented form of the dipole operator")
501 CALL timestop(handle)
530 dipole_op_mos_occ, dipole_form)
531 INTEGER,
INTENT(in) :: log_unit
532 TYPE(
cp_fm_type),
DIMENSION(:, :),
INTENT(in) :: evects
533 REAL(kind=
dp),
DIMENSION(:),
INTENT(in) :: evals
536 REAL(kind=
dp),
DIMENSION(:),
INTENT(inout) :: ostrength
537 INTEGER,
INTENT(in) :: mult
538 TYPE(
cp_fm_type),
DIMENSION(:, :),
INTENT(in) :: dipole_op_mos_occ
539 INTEGER,
INTENT(in) :: dipole_form
541 CHARACTER(LEN=*),
PARAMETER :: routinen =
'tddfpt_print_summary'
543 CHARACTER(len=1) :: lsd_str
544 CHARACTER(len=20) :: mult_str
545 INTEGER :: handle, i, ideriv, ispin, istate, j, &
546 nactive, nao, nocc, nspins, nstates
547 REAL(kind=
dp) :: osc_strength
548 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :, :) :: trans_dipoles
552 CALL timeset(routinen, handle)
554 nspins =
SIZE(evects, 1)
555 nstates =
SIZE(evects, 2)
564 IF (log_unit > 0)
THEN
566 WRITE (log_unit,
'(/,1X,A1,A,1X,A)') lsd_str,
"-TDDFPT states of multiplicity", trim(mult_str)
567 SELECT CASE (dipole_form)
569 WRITE (log_unit,
'(1X,A,/)')
"Transition dipoles calculated using Berry operator formulation"
571 WRITE (log_unit,
'(1X,A,/)')
"Transition dipoles calculated using length formulation"
573 WRITE (log_unit,
'(1X,A,/)')
"Transition dipoles calculated using velocity formulation"
575 WRITE (log_unit,
'(1X,A,/)')
"Transition dipoles calculated using old velocity formulation"
577 WRITE (log_unit,
'(1X,A,/)')
"Transition dipoles calculated using SCF-MO moment formulation"
579 cpabort(
"Unimplemented form of the dipole operator")
582 WRITE (log_unit,
'(T10,A,T19,A,T37,A,T69,A)')
"State",
"Excitation", &
583 "Transition dipole (a.u.)",
"Oscillator"
584 WRITE (log_unit,
'(T10,A,T19,A,T37,A,T49,A,T61,A,T67,A)')
"number",
"energy (eV)", &
585 "x",
"y",
"z",
"strength (a.u.)"
586 WRITE (log_unit,
'(T10,72("-"))')
590 ALLOCATE (trans_dipoles(nstates, nderivs, nspins))
591 trans_dipoles(:, :, :) = 0.0_dp
594 IF (nspins > 1 .OR. mult == 1)
THEN
596 CALL cp_fm_get_info(dipole_op_mos_occ(1, ispin), nrow_global=nao, ncol_global=nocc)
598 IF (nocc == nactive)
THEN
599 DO ideriv = 1, nderivs
600 CALL cp_fm_trace(evects(ispin, :), dipole_op_mos_occ(ideriv, ispin), &
601 trans_dipoles(:, ideriv, ispin))
604 CALL cp_fm_get_info(evects(ispin, 1), matrix_struct=matrix_struct)
606 DO ideriv = 1, nderivs
608 j = gs_mos(ispin)%index_active(i)
609 CALL cp_fm_to_fm(dipole_op_mos_occ(ideriv, ispin), dipact, &
610 ncol=1, source_start=j, target_start=i)
612 CALL cp_fm_trace(evects(ispin, :), dipact, trans_dipoles(:, ideriv, ispin))
618 IF (nspins == 1)
THEN
619 trans_dipoles(:, :, 1) = sqrt(2.0_dp)*trans_dipoles(:, :, 1)
621 trans_dipoles(:, :, 1) = trans_dipoles(:, :, 1) + trans_dipoles(:, :, 2)
626 DO istate = 1, nstates
628 SELECT CASE (dipole_form)
630 osc_strength = 2.0_dp/3.0_dp*evals(istate)*sum(trans_dipoles(istate, :, 1)**2)
632 osc_strength = 2.0_dp/3.0_dp*evals(istate)*sum(trans_dipoles(istate, :, 1)**2)
634 osc_strength = 2.0_dp/3.0_dp/evals(istate)*sum(trans_dipoles(istate, :, 1)**2)
636 osc_strength = 2.0_dp/3.0_dp*evals(istate)*sum(trans_dipoles(istate, :, 1)**2)
638 cpabort(
"Unimplemented form of the dipole operator")
641 ostrength(istate) = osc_strength
642 IF (log_unit > 0)
THEN
643 WRITE (log_unit,
'(1X,A,T9,I7,T19,F11.5,T31,3(1X,ES11.4E2),T69,ES12.5E2)') &
644 "TDDFPT|", istate, evals(istate)*
evolt, trans_dipoles(istate, 1:nderivs, 1), osc_strength
649 IF (log_unit > 0)
THEN
650 WRITE (log_unit,
'(/,T2,A,E16.8)')
'TDDFPT : CheckSum E = ', sqrt(sum(evals**2))
651 WRITE (log_unit,
'(/,T2,A,E16.8)')
'TDDFPT : CheckSum F = ', sqrt(sum(ostrength**2))
654 DEALLOCATE (trans_dipoles)
656 CALL timestop(handle)
675 INTEGER,
INTENT(in) :: log_unit
676 TYPE(
cp_fm_type),
DIMENSION(:, :),
INTENT(in) :: evects
677 REAL(kind=
dp),
DIMENSION(:),
INTENT(in) :: evals
682 REAL(kind=
dp),
INTENT(in) :: min_amplitude
684 CHARACTER(LEN=*),
PARAMETER :: routinen =
'tddfpt_print_excitation_analysis'
686 CHARACTER(len=5) :: spin_label, spin_label2
687 INTEGER :: handle, icol, iproc, irow, ispin, &
688 istate, nao, ncols_local, nrows_local, &
689 nspins, nstates, spin2, state_spin, &
691 INTEGER(kind=int_8) :: iexc, imo_act, imo_occ, imo_virt, ind, &
692 nexcs, nexcs_local, nexcs_max_local, &
693 nmo_virt_occ, nmo_virt_occ_alpha
694 INTEGER(kind=int_8),
ALLOCATABLE,
DIMENSION(:) :: inds_local, inds_recv, nexcs_recv
695 INTEGER(kind=int_8),
DIMENSION(1) :: nexcs_send
696 INTEGER(kind=int_8),
DIMENSION(maxspins) :: nactive8, nmo_occ8, nmo_virt8
697 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: inds
698 INTEGER,
DIMENSION(:),
POINTER :: col_indices, row_indices
699 INTEGER,
DIMENSION(maxspins) :: nactive, nmo_occ, nmo_virt
700 LOGICAL :: do_exc_analysis
701 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: weights_local, weights_neg_abs_recv, &
703 REAL(kind=
dp),
CONTIGUOUS,
DIMENSION(:, :), &
704 POINTER :: local_data
707 TYPE(
cp_fm_type),
ALLOCATABLE,
DIMENSION(:) :: s_mos_virt, weights_fm
710 TYPE(
mp_request_type),
ALLOCATABLE,
DIMENSION(:) :: recv_handlers, recv_handlers2
712 CALL timeset(routinen, handle)
714 nspins =
SIZE(gs_mos, 1)
715 nstates =
SIZE(evects, 2)
716 do_exc_analysis = min_amplitude < 1.0_dp
718 CALL cp_fm_get_info(gs_mos(1)%mos_occ, context=blacs_env, para_env=para_env)
722 nactive(ispin) = gs_mos(ispin)%nmo_active
723 nactive8(ispin) = int(nactive(ispin), kind=
int_8)
724 nmo_occ(ispin) =
SIZE(gs_mos(ispin)%evals_occ)
725 nmo_virt(ispin) =
SIZE(gs_mos(ispin)%evals_virt)
726 nmo_occ8(ispin) =
SIZE(gs_mos(ispin)%evals_occ, kind=
int_8)
727 nmo_virt8(ispin) =
SIZE(gs_mos(ispin)%evals_virt, kind=
int_8)
731 IF (do_exc_analysis)
THEN
732 cpassert(log_unit <= 0 .OR. para_env%is_source())
733 nmo_virt_occ_alpha = int(nmo_virt(1),
int_8)*int(nmo_occ(1),
int_8)
735 IF (log_unit > 0)
THEN
736 WRITE (log_unit,
"(1X,A)")
"", &
737 "-------------------------------------------------------------------------------", &
738 "- Excitation analysis -", &
739 "-------------------------------------------------------------------------------"
740 WRITE (log_unit,
'(8X,A,T27,A,T49,A,T69,A)')
"State",
"Occupied",
"Virtual",
"Excitation"
741 WRITE (log_unit,
'(8X,A,T28,A,T49,A,T69,A)')
"number",
"orbital",
"orbital",
"amplitude"
742 WRITE (log_unit,
'(1X,79("-"))')
744 IF (nspins == 1)
THEN
753 spin_label2 =
'(bet)'
757 ALLOCATE (s_mos_virt(
SIZE(evects, 1)), weights_fm(
SIZE(evects, 1)))
758 DO ispin = 1,
SIZE(evects, 1)
764 CALL cp_fm_get_info(gs_mos(spin2)%mos_virt, matrix_struct=fm_struct)
767 gs_mos(spin2)%mos_virt, &
769 ncol=nmo_virt(spin2), alpha=1.0_dp, beta=0.0_dp)
772 CALL cp_fm_struct_create(fm_struct, nrow_global=nmo_virt(spin2), ncol_global=nactive(ispin), &
780 DO ispin = 1,
SIZE(evects, 1)
781 CALL cp_fm_get_info(weights_fm(ispin), nrow_local=nrows_local, ncol_local=ncols_local)
782 nexcs_max_local = nexcs_max_local + int(nrows_local,
int_8)*int(ncols_local,
int_8)
785 ALLOCATE (weights_local(nexcs_max_local), inds_local(nexcs_max_local))
787 DO istate = 1, nstates
793 DO ispin = 1,
SIZE(evects, 1)
800 CALL parallel_gemm(
'T',
'N', nmo_virt(spin2), nactive(ispin), nao, 1.0_dp, s_mos_virt(ispin), &
801 evects(ispin, istate), 0.0_dp, weights_fm(ispin))
803 CALL cp_fm_get_info(weights_fm(ispin), nrow_local=nrows_local, ncol_local=ncols_local, &
804 row_indices=row_indices, col_indices=col_indices, local_data=local_data)
807 DO icol = 1, ncols_local
808 DO irow = 1, nrows_local
809 IF (abs(local_data(irow, icol)) >= min_amplitude)
THEN
811 nexcs_local = nexcs_local + 1
813 weights_local(nexcs_local) = local_data(irow, icol)
815 inds_local(nexcs_local) = nmo_virt_occ + int(row_indices(irow),
int_8) + &
816 int(col_indices(icol) - 1,
int_8)*nmo_virt8(spin2)
821 nmo_virt_occ = nmo_virt_occ + nmo_virt8(spin2)*nmo_occ8(ispin)
824 IF (para_env%is_source())
THEN
826 ALLOCATE (nexcs_recv(para_env%num_pe), recv_handlers(para_env%num_pe), recv_handlers2(para_env%num_pe))
829 DO iproc = 1, para_env%num_pe
830 IF (iproc - 1 /= para_env%mepos)
THEN
831 CALL para_env%irecv(nexcs_recv(iproc:iproc), iproc - 1, recv_handlers(iproc), 0)
833 nexcs_recv(iproc) = nexcs_local
837 DO iproc = 1, para_env%num_pe
838 IF (iproc - 1 /= para_env%mepos) &
839 CALL recv_handlers(iproc)%wait()
844 DO iproc = 1, para_env%num_pe
845 nexcs = nexcs + nexcs_recv(iproc)
849 ALLOCATE (weights_recv(nexcs), weights_neg_abs_recv(nexcs))
850 ALLOCATE (inds_recv(nexcs), inds(nexcs))
853 DO iproc = 1, para_env%num_pe
854 IF (nexcs_recv(iproc) > 0)
THEN
855 IF (iproc - 1 /= para_env%mepos)
THEN
857 CALL para_env%irecv(weights_recv(nmo_virt_occ + 1:nmo_virt_occ + nexcs_recv(iproc)), &
858 iproc - 1, recv_handlers(iproc), 1)
860 CALL para_env%irecv(inds_recv(nmo_virt_occ + 1:nmo_virt_occ + nexcs_recv(iproc)), &
861 iproc - 1, recv_handlers2(iproc), 2)
864 weights_recv(nmo_virt_occ + 1:nmo_virt_occ + nexcs_recv(iproc)) = weights_local(1:nexcs_recv(iproc))
865 inds_recv(nmo_virt_occ + 1:nmo_virt_occ + nexcs_recv(iproc)) = inds_local(1:nexcs_recv(iproc))
868 nmo_virt_occ = nmo_virt_occ + nexcs_recv(iproc)
872 DO iproc = 1, para_env%num_pe
873 IF (iproc - 1 /= para_env%mepos .AND. nexcs_recv(iproc) > 0)
THEN
874 CALL recv_handlers(iproc)%wait()
875 CALL recv_handlers2(iproc)%wait()
879 DEALLOCATE (nexcs_recv, recv_handlers, recv_handlers2)
882 nexcs_send(1) = nexcs_local
883 CALL para_env%isend(nexcs_send, para_env%source, send_handler, 0)
884 CALL send_handler%wait()
886 IF (nexcs_local > 0)
THEN
888 CALL para_env%isend(weights_local(1:nexcs_local), para_env%source, send_handler, 1)
890 CALL para_env%isend(inds_local(1:nexcs_local), para_env%source, send_handler2, 2)
892 CALL send_handler%wait()
893 CALL send_handler2%wait()
899 IF (para_env%is_source() .AND. log_unit > 0)
THEN
900 weights_neg_abs_recv(:) = -abs(weights_recv)
901 CALL sort(weights_neg_abs_recv, int(nexcs), inds)
903 WRITE (log_unit,
'(T7,I8,F10.5,A)') istate, evals(istate)*
evolt,
" eV"
914 ind = inds_recv(inds(iexc)) - 1
916 IF (ind < nmo_virt_occ_alpha)
THEN
920 spin_label2 =
'(alp)'
924 ind = ind - nmo_virt_occ_alpha
926 spin_label2 =
'(bet)'
929 imo_act = ind/nmo_virt8(state_spin2) + 1
930 imo_occ = gs_mos(state_spin)%index_active(imo_act)
931 imo_virt = mod(ind, nmo_virt8(state_spin2)) + 1
933 WRITE (log_unit,
'(T27,I8,1X,A5,T48,I8,1X,A5,T70,F9.6)') imo_occ, spin_label, &
934 nmo_occ8(state_spin2) + imo_virt, spin_label2, weights_recv(inds(iexc))
939 IF (para_env%is_source()) &
940 DEALLOCATE (weights_recv, weights_neg_abs_recv, inds_recv, inds)
943 DEALLOCATE (weights_local, inds_local)
944 IF (log_unit > 0)
THEN
945 WRITE (log_unit,
"(1X,A)") &
946 "-------------------------------------------------------------------------------"
953 CALL timestop(handle)
972 TYPE(
cp_fm_type),
DIMENSION(:, :),
INTENT(in) :: evects
973 REAL(kind=
dp),
DIMENSION(:),
INTENT(in) :: evals, ostrength
979 CHARACTER(LEN=*),
PARAMETER :: routinen =
'tddfpt_print_nto_analysis'
980 INTEGER,
PARAMETER :: ntomax = 10
982 CHARACTER(LEN=20),
DIMENSION(2) :: nto_name
983 INTEGER :: handle, i, ia, icg, iounit, ispin, &
984 istate, j, nao, nlist, nmax, nmo, &
985 nnto, nspins, nstates
986 INTEGER,
DIMENSION(2) :: iv
987 INTEGER,
DIMENSION(2, ntomax) :: ia_index
988 INTEGER,
DIMENSION(:),
POINTER :: slist, stride
989 LOGICAL :: append_cube, cube_file, explicit
990 REAL(kind=
dp) :: os_threshold, sume, threshold
991 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: eigvals
992 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: eigenvalues
993 REAL(kind=
dp),
DIMENSION(ntomax) :: ia_eval
996 TYPE(
cp_fm_type) :: sev, smat, tmat, wmat, work, wvec
997 TYPE(
cp_fm_type),
ALLOCATABLE,
DIMENSION(:) :: teig
999 TYPE(
mo_set_type),
ALLOCATABLE,
DIMENSION(:) :: nto_set
1001 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
1004 CALL timeset(routinen, handle)
1015 CALL section_vals_val_get(print_section,
"NTO_ANALYSIS%INTENSITY_THRESHOLD", r_val=os_threshold)
1025 IF (iounit > 0)
THEN
1026 WRITE (iounit,
"(1X,A)")
"", &
1027 "-------------------------------------------------------------------------------", &
1028 "- Natural Orbital analysis -", &
1029 "-------------------------------------------------------------------------------"
1032 nspins =
SIZE(evects, 1)
1033 nstates =
SIZE(evects, 2)
1036 DO istate = 1, nstates
1037 IF (os_threshold > ostrength(istate))
THEN
1038 IF (iounit > 0)
THEN
1039 WRITE (iounit,
"(1X,A,I6)")
" Skipping state ", istate
1044 IF (.NOT. any(slist == istate))
THEN
1045 IF (iounit > 0)
THEN
1046 WRITE (iounit,
"(1X,A,I6)")
" Skipping state ", istate
1051 IF (iounit > 0)
THEN
1052 WRITE (iounit,
"(1X,A,I6,T30,F10.5,A)")
" STATE NR. ", istate, evals(istate)*
evolt,
" eV"
1055 DO ispin = 1, nspins
1056 CALL cp_fm_get_info(evects(ispin, istate), matrix_struct=fm_struct, ncol_global=nmo)
1057 nmax = max(nmax, nmo)
1059 ALLOCATE (eigenvalues(nmax, nspins))
1060 eigenvalues = 0.0_dp
1063 nto_name(1) =
'Hole_states'
1064 nto_name(2) =
'Particle_states'
1065 ALLOCATE (nto_set(2))
1067 CALL allocate_mo_set(nto_set(i), nao, ntomax, 0, 0.0_dp, 1.0_dp, 0.0_dp)
1072 CALL init_mo_set(nto_set(i), fm_ref=tmat, name=nto_name(i))
1077 ALLOCATE (teig(nspins))
1080 DO ispin = 1, nspins
1081 associate(ev => evects(ispin, istate))
1082 CALL cp_fm_get_info(ev, matrix_struct=fm_struct, ncol_global=nmo)
1085 nrow_global=nmo, ncol_global=nmo)
1089 CALL parallel_gemm(
'T',
'N', nmo, nmo, nao, 1.0_dp, ev, sev, 0.0_dp, tmat)
1103 iv = maxloc(eigenvalues)
1104 ia_eval(i) = eigenvalues(iv(1), iv(2))
1105 ia_index(1:2, i) = iv(1:2)
1106 sume = sume + ia_eval(i)
1107 eigenvalues(iv(1), iv(2)) = 0.0_dp
1109 IF (sume > threshold)
EXIT
1115 ispin = ia_index(2, i)
1119 nrow_global=nmo, ncol_global=1)
1128 CALL parallel_gemm(
'N',
'N', nao, 1, nmo, 1.0_dp, gs_mos(ispin)%mos_occ, &
1130 CALL cp_fm_to_fm(wvec, nto_set(1)%mo_coeff, 1, 1, i)
1137 DO ispin = 1, nspins
1138 associate(ev => evects(ispin, istate))
1139 CALL cp_fm_get_info(ev, matrix_struct=fm_struct, nrow_global=nao, ncol_global=nmo)
1140 ALLOCATE (eigvals(nao))
1145 CALL cp_fm_struct_create(fmstruct=fm_mo_struct, template_fmstruct=fm_struct, &
1146 nrow_global=nao, ncol_global=nao)
1147 CALL cp_fm_create(tmat, fm_mo_struct)
1148 CALL cp_fm_create(smat, fm_mo_struct)
1149 CALL cp_fm_create(wmat, fm_mo_struct)
1150 CALL cp_fm_create(work, fm_mo_struct)
1151 CALL cp_fm_struct_release(fm_mo_struct)
1152 CALL copy_dbcsr_to_fm(matrix_s, smat)
1153 CALL parallel_gemm(
'N',
'T', nao, nao, nmo, 1.0_dp, sev, sev, 0.0_dp, tmat)
1154 CALL cp_fm_geeig(tmat, smat, wmat, eigvals, work)
1156 IF (ispin == ia_index(2, i))
THEN
1159 IF (abs(eigvals(j) - ia_eval(i)) < 1.e-6_dp)
THEN
1165 CALL cp_warn(__location__, &
1166 "Could not locate particle state associated with hole state.")
1168 CALL cp_fm_to_fm(wmat, nto_set(2)%mo_coeff, 1, icg, i)
1172 DEALLOCATE (eigvals)
1173 CALL cp_fm_release(sev)
1174 CALL cp_fm_release(tmat)
1175 CALL cp_fm_release(smat)
1176 CALL cp_fm_release(wmat)
1177 CALL cp_fm_release(work)
1180 IF (iounit > 0)
THEN
1183 sume = sume + ia_eval(i)
1184 WRITE (iounit,
"(T6,A,i2,T30,A,i1,T42,A,F8.5,T63,A,F8.5)") &
1185 "Particle-Hole state:", i,
" Spin:", ia_index(2, i), &
1186 "Eigenvalue:", ia_eval(i),
" Sum Eigv:", sume
1190 nto_section => section_vals_get_subs_vals(print_section,
"NTO_ANALYSIS")
1191 CALL section_vals_val_get(nto_section,
"CUBE_FILES", l_val=cube_file)
1192 CALL section_vals_val_get(nto_section,
"STRIDE", i_vals=stride)
1193 CALL section_vals_val_get(nto_section,
"APPEND", l_val=append_cube)
1195 CALL print_nto_cubes(qs_env, nto_set, istate, stride, append_cube, nto_section)
1197 CALL get_qs_env(qs_env, qs_kind_set=qs_kind_set, particle_set=particle_set, cell=cell)
1198 molden_section => section_vals_get_subs_vals(print_section,
"MOS_MOLDEN")
1199 CALL write_mos_molden(nto_set, qs_kind_set, particle_set, molden_section, cell=cell, qs_env=qs_env)
1201 DEALLOCATE (eigenvalues)
1202 CALL cp_fm_release(teig)
1205 CALL deallocate_mo_set(nto_set(i))
1207 DEALLOCATE (nto_set)
1210 IF (iounit > 0)
THEN
1211 WRITE (iounit,
"(1X,A)") &
1212 "-------------------------------------------------------------------------------"
1217 CALL timestop(handle)
1234 do_directional_exciton_descriptors, qs_env)
1235 INTEGER,
INTENT(in) :: log_unit
1236 TYPE(cp_fm_type),
DIMENSION(:, :),
INTENT(in) :: evects
1237 TYPE(tddfpt_ground_state_mos),
DIMENSION(:), &
1238 INTENT(in) :: gs_mos
1239 TYPE(dbcsr_type),
POINTER :: matrix_s
1240 LOGICAL,
INTENT(IN) :: do_directional_exciton_descriptors
1241 TYPE(qs_environment_type),
INTENT(IN),
POINTER :: qs_env
1243 CHARACTER(LEN=*),
PARAMETER :: routinen =
'tddfpt_print_exciton_descriptors'
1245 CHARACTER(LEN=4) :: prefix_output
1246 INTEGER :: handle, ispin, istate, n_moments_quad, &
1247 nactive, nao, nspins, nstates
1248 INTEGER,
DIMENSION(maxspins) :: nmo_occ, nmo_virt
1249 LOGICAL :: print_checkvalue
1250 REAL(kind=dp),
ALLOCATABLE,
DIMENSION(:) :: ref_point_multipole
1251 TYPE(cp_blacs_env_type),
POINTER :: blacs_env
1252 TYPE(cp_fm_struct_type),
POINTER :: fm_struct_mo_coeff, &
1253 fm_struct_s_mos_virt, fm_struct_x_ia_n
1254 TYPE(cp_fm_type),
ALLOCATABLE,
DIMENSION(:) :: eigvec_x_ia_n, fm_multipole_ab, &
1255 fm_multipole_ai, fm_multipole_ij, &
1257 TYPE(cp_fm_type),
DIMENSION(:),
POINTER :: mo_coeff
1258 TYPE(exciton_descr_type),
ALLOCATABLE, &
1259 DIMENSION(:) :: exc_descr
1261 CALL timeset(routinen, handle)
1263 nspins =
SIZE(evects, 1)
1264 nstates =
SIZE(evects, 2)
1266 cpassert(nspins == 1)
1268 CALL cp_fm_get_info(gs_mos(1)%mos_occ, context=blacs_env)
1269 CALL dbcsr_get_info(matrix_s, nfullrows_total=nao)
1271 DO ispin = 1, nspins
1272 nmo_occ(ispin) =
SIZE(gs_mos(ispin)%evals_occ)
1273 nmo_virt(ispin) =
SIZE(gs_mos(ispin)%evals_virt)
1277 ALLOCATE (mo_coeff(nspins))
1278 CALL cp_fm_struct_create(fm_struct_mo_coeff, nrow_global=nao, ncol_global=nao, &
1280 DO ispin = 1, nspins
1281 CALL cp_fm_create(mo_coeff(ispin), fm_struct_mo_coeff)
1282 CALL cp_fm_to_fm_submat_general(gs_mos(ispin)%mos_occ, &
1291 CALL cp_fm_to_fm_submat_general(gs_mos(ispin)%mos_virt, &
1298 nmo_occ(ispin) + 1, &
1301 CALL cp_fm_struct_release(fm_struct_mo_coeff)
1306 ALLOCATE (ref_point_multipole(3))
1307 ALLOCATE (fm_multipole_ij(n_moments_quad))
1308 ALLOCATE (fm_multipole_ab(n_moments_quad))
1309 ALLOCATE (fm_multipole_ai(n_moments_quad))
1311 CALL get_multipoles_mo(fm_multipole_ai, fm_multipole_ij, fm_multipole_ab, &
1312 qs_env, mo_coeff, ref_point_multipole, 2, &
1313 nmo_occ(1), nmo_virt(1), blacs_env)
1315 CALL cp_fm_release(mo_coeff)
1318 ALLOCATE (s_mos_virt(nspins), eigvec_x_ia_n(nspins))
1319 DO ispin = 1, nspins
1320 CALL cp_fm_get_info(gs_mos(ispin)%mos_virt, matrix_struct=fm_struct_s_mos_virt)
1321 CALL cp_fm_create(s_mos_virt(ispin), fm_struct_s_mos_virt)
1322 NULLIFY (fm_struct_s_mos_virt)
1323 CALL cp_dbcsr_sm_fm_multiply(matrix_s, &
1324 gs_mos(ispin)%mos_virt, &
1325 s_mos_virt(ispin), &
1326 ncol=nmo_virt(ispin), alpha=1.0_dp, beta=0.0_dp)
1328 CALL cp_fm_struct_create(fm_struct_x_ia_n, nrow_global=nmo_occ(ispin), ncol_global=nmo_virt(ispin), &
1330 CALL cp_fm_create(eigvec_x_ia_n(ispin), fm_struct_x_ia_n)
1331 CALL cp_fm_struct_release(fm_struct_x_ia_n)
1333 ALLOCATE (exc_descr(nstates))
1334 DO istate = 1, nstates
1335 DO ispin = 1, nspins
1336 CALL cp_fm_set_all(eigvec_x_ia_n(ispin), 0.0_dp)
1343 CALL cp_fm_get_info(evects(ispin, istate), ncol_global=nactive)
1344 IF (nactive /= nmo_occ(ispin))
THEN
1345 CALL cp_abort(__location__, &
1346 "Reduced active space excitations not implemented")
1348 CALL parallel_gemm(
'T',
'N', nmo_occ(ispin), nmo_virt(ispin), nao, 1.0_dp, &
1349 evects(ispin, istate), s_mos_virt(ispin), 0.0_dp, eigvec_x_ia_n(ispin))
1351 CALL get_exciton_descriptors(exc_descr, eigvec_x_ia_n(ispin), &
1352 fm_multipole_ij, fm_multipole_ab, &
1354 istate, nmo_occ(ispin), nmo_virt(ispin))
1357 CALL cp_fm_release(eigvec_x_ia_n)
1358 CALL cp_fm_release(s_mos_virt)
1359 CALL cp_fm_release(fm_multipole_ai)
1360 CALL cp_fm_release(fm_multipole_ij)
1361 CALL cp_fm_release(fm_multipole_ab)
1364 print_checkvalue = .true.
1366 CALL print_exciton_descriptors(exc_descr, ref_point_multipole, log_unit, &
1367 nstates, print_checkvalue, do_directional_exciton_descriptors, &
1368 prefix_output, qs_env)
1370 DEALLOCATE (ref_point_multipole)
1371 DEALLOCATE (exc_descr)
1373 CALL timestop(handle)
1384 SUBROUTINE project_vector(vin, vout, mos_occ, matrix_s)
1385 TYPE(dbcsr_type) :: vin, vout
1386 TYPE(cp_fm_type),
INTENT(IN) :: mos_occ
1387 TYPE(dbcsr_type),
POINTER :: matrix_s
1389 CHARACTER(LEN=*),
PARAMETER :: routinen =
'project_vector'
1391 INTEGER :: handle, nao, nmo
1392 REAL(kind=dp) :: norm(1)
1393 TYPE(cp_fm_struct_type),
POINTER :: fm_struct, fm_vec_struct
1394 TYPE(cp_fm_type) :: csvec, svec, vec
1396 CALL timeset(routinen, handle)
1398 CALL cp_fm_get_info(mos_occ, matrix_struct=fm_struct, nrow_global=nao, ncol_global=nmo)
1399 CALL cp_fm_struct_create(fmstruct=fm_vec_struct, template_fmstruct=fm_struct, &
1400 nrow_global=nao, ncol_global=1)
1401 CALL cp_fm_create(vec, fm_vec_struct)
1402 CALL cp_fm_create(svec, fm_vec_struct)
1403 CALL cp_fm_struct_release(fm_vec_struct)
1404 CALL cp_fm_struct_create(fmstruct=fm_vec_struct, template_fmstruct=fm_struct, &
1405 nrow_global=nmo, ncol_global=1)
1406 CALL cp_fm_create(csvec, fm_vec_struct)
1407 CALL cp_fm_struct_release(fm_vec_struct)
1409 CALL copy_dbcsr_to_fm(vin, vec)
1410 CALL cp_dbcsr_sm_fm_multiply(matrix_s, vec, svec, ncol=1, alpha=1.0_dp, beta=0.0_dp)
1411 CALL parallel_gemm(
'T',
'N', nmo, 1, nao, 1.0_dp, mos_occ, svec, 0.0_dp, csvec)
1412 CALL parallel_gemm(
'N',
'N', nao, 1, nmo, -1.0_dp, mos_occ, csvec, 1.0_dp, vec)
1413 CALL cp_fm_vectorsnorm(vec, norm)
1414 cpassert(norm(1) > 1.e-14_dp)
1415 norm(1) = sqrt(1._dp/norm(1))
1416 CALL cp_fm_scale(norm(1), vec)
1417 CALL copy_fm_to_dbcsr(vec, vout, keep_sparsity=.false.)
1419 CALL cp_fm_release(csvec)
1420 CALL cp_fm_release(svec)
1421 CALL cp_fm_release(vec)
1423 CALL timestop(handle)
1425 END SUBROUTINE project_vector
1433 SUBROUTINE vec_product(va, vb, res)
1434 TYPE(dbcsr_type) :: va, vb
1435 REAL(kind=dp),
INTENT(OUT) :: res
1437 CHARACTER(LEN=*),
PARAMETER :: routinen =
'vec_product'
1439 INTEGER :: handle, icol, irow
1441 REAL(kind=dp),
DIMENSION(:, :),
POINTER :: vba, vbb
1442 TYPE(dbcsr_iterator_type) :: iter
1443 TYPE(mp_comm_type) :: group
1445 CALL timeset(routinen, handle)
1449 CALL dbcsr_get_info(va, group=group)
1450 CALL dbcsr_iterator_start(iter, va)
1451 DO WHILE (dbcsr_iterator_blocks_left(iter))
1452 CALL dbcsr_iterator_next_block(iter, irow, icol, vba)
1453 CALL dbcsr_get_block_p(vb, row=irow, col=icol, block=vbb, found=found)
1454 res = res + sum(vba*vbb)
1457 CALL dbcsr_iterator_stop(iter)
1460 CALL timestop(handle)
1462 END SUBROUTINE vec_product
1473 SUBROUTINE print_nto_cubes(qs_env, mos, istate, stride, append_cube, print_section)
1475 TYPE(qs_environment_type),
POINTER :: qs_env
1476 TYPE(mo_set_type),
DIMENSION(:),
INTENT(IN) :: mos
1477 INTEGER,
INTENT(IN) :: istate
1478 INTEGER,
DIMENSION(:),
POINTER :: stride
1479 LOGICAL,
INTENT(IN) :: append_cube
1480 TYPE(section_vals_type),
POINTER :: print_section
1482 CHARACTER(LEN=default_path_length) :: filename, my_pos_cube, title
1483 INTEGER :: i, iset, nmo, unit_nr
1485 TYPE(atomic_kind_type),
DIMENSION(:),
POINTER :: atomic_kind_set
1486 TYPE(cell_type),
POINTER :: cell
1487 TYPE(cp_fm_type),
POINTER :: mo_coeff
1488 TYPE(cp_logger_type),
POINTER :: logger
1489 TYPE(dft_control_type),
POINTER :: dft_control
1490 TYPE(particle_list_type),
POINTER :: particles
1491 TYPE(particle_type),
DIMENSION(:),
POINTER :: particle_set
1492 TYPE(pw_c1d_gs_type) :: wf_g
1493 TYPE(pw_env_type),
POINTER :: pw_env
1494 TYPE(pw_pool_p_type),
DIMENSION(:),
POINTER :: pw_pools
1495 TYPE(pw_pool_type),
POINTER :: auxbas_pw_pool
1496 TYPE(pw_r3d_rs_type) :: wf_r
1497 TYPE(qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
1498 TYPE(qs_subsys_type),
POINTER :: subsys
1500 logger => cp_get_default_logger()
1502 CALL get_qs_env(qs_env=qs_env, dft_control=dft_control, pw_env=pw_env)
1503 CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool, pw_pools=pw_pools)
1504 CALL auxbas_pw_pool%create_pw(wf_r)
1505 CALL auxbas_pw_pool%create_pw(wf_g)
1507 CALL get_qs_env(qs_env, subsys=subsys)
1508 CALL qs_subsys_get(subsys, particles=particles)
1510 my_pos_cube =
"REWIND"
1511 IF (append_cube)
THEN
1512 my_pos_cube =
"APPEND"
1515 CALL get_qs_env(qs_env=qs_env, &
1516 atomic_kind_set=atomic_kind_set, &
1517 qs_kind_set=qs_kind_set, &
1519 particle_set=particle_set)
1522 CALL get_mo_set(mo_set=mos(iset), mo_coeff=mo_coeff, nmo=nmo)
1524 CALL calculate_wavefunction(mo_coeff, i, wf_r, wf_g, atomic_kind_set, qs_kind_set, &
1525 cell, dft_control, particle_set, pw_env)
1527 WRITE (filename,
'(a4,I3.3,I2.2,a11)')
"NTO_STATE", istate, i,
"_Hole_State"
1528 ELSEIF (iset == 2)
THEN
1529 WRITE (filename,
'(a4,I3.3,I2.2,a15)')
"NTO_STATE", istate, i,
"_Particle_State"
1532 unit_nr = cp_print_key_unit_nr(logger, print_section,
'', extension=
".cube", &
1533 middle_name=trim(filename), file_position=my_pos_cube, &
1534 log_filename=.false., ignore_should_output=.true., mpi_io=mpi_io)
1536 WRITE (title, *)
"Natural Transition Orbital Hole State", i
1537 ELSEIF (iset == 2)
THEN
1538 WRITE (title, *)
"Natural Transition Orbital Particle State", i
1540 CALL cp_pw_to_cube(wf_r, unit_nr, title, particles=particles, stride=stride, mpi_io=mpi_io)
1541 CALL cp_print_key_finished_output(unit_nr, logger, print_section,
'', &
1542 ignore_should_output=.true., mpi_io=mpi_io)
1546 CALL auxbas_pw_pool%give_back_pw(wf_g)
1547 CALL auxbas_pw_pool%give_back_pw(wf_r)
1549 END SUBROUTINE print_nto_cubes
Define the atomic kind types and their sub types.
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public martin2003
Routines for printing information in context of the BSE calculation.
subroutine, public print_exciton_descriptors(exc_descr, ref_point_multipole, unit_nr, num_print_exc_descr, print_checkvalue, print_directional_exc_descr, prefix_output, qs_env)
Prints exciton descriptors, cf. Mewes et al., JCTC 14, 710-725 (2018)
Routines for computing excitonic properties, e.g. exciton diameter, from the BSE.
subroutine, public get_exciton_descriptors(exc_descr, fm_x_ia, fm_multipole_ij_trunc, fm_multipole_ab_trunc, fm_multipole_ai_trunc, i_exc, homo, virtual, fm_y_ia)
...
Auxiliary routines for GW + Bethe-Salpeter for computing electronic excitations.
subroutine, public get_multipoles_mo(fm_multipole_ai_trunc, fm_multipole_ij_trunc, fm_multipole_ab_trunc, qs_env, mo_coeff, rpoint, n_moments, homo_red, virtual_red, context_bse)
...
Handles all functions related to the CELL.
methods related to the blacs parallel environment
Basic linear algebra operations for complex full matrices.
subroutine, public cp_cfm_solve(matrix_a, general_a, determinant)
Solve the system of linear equations A*b=A_general using LU decomposition. Pay attention that both ma...
Represents a complex full matrix distributed on many processors.
subroutine, public cp_cfm_release(matrix)
Releases a full matrix.
subroutine, public cp_fm_to_cfm(msourcer, msourcei, mtarget)
Construct a complex full matrix by taking its real and imaginary parts from two separate real-value f...
subroutine, public cp_cfm_create(matrix, matrix_struct, name, nrow, ncol, set_zero)
Creates a new full matrix with the given structure.
subroutine, public cp_cfm_set_all(matrix, alpha, beta)
Set all elements of the full matrix to alpha. Besides, set all diagonal matrix elements to beta (if g...
subroutine, public cp_cfm_to_fm(msource, mtargetr, mtargeti)
Copy real and imaginary parts of a complex full matrix into separate real-value full matrices.
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
logical function, public dbcsr_iterator_blocks_left(iterator)
...
subroutine, public dbcsr_iterator_stop(iterator)
...
subroutine, public dbcsr_copy(matrix_b, matrix_a, name, keep_sparsity, keep_imaginary)
...
subroutine, public dbcsr_get_block_p(matrix, row, col, block, found, row_size, col_size)
...
subroutine, public dbcsr_get_info(matrix, nblkrows_total, nblkcols_total, nfullrows_total, nfullcols_total, nblkrows_local, nblkcols_local, nfullrows_local, nfullcols_local, my_prow, my_pcol, local_rows, local_cols, proc_row_dist, proc_col_dist, row_blk_size, col_blk_size, row_blk_offset, col_blk_offset, distribution, name, matrix_type, group)
...
subroutine, public dbcsr_init_p(matrix)
...
subroutine, public dbcsr_iterator_next_block(iterator, row, column, block, block_number_argument_has_been_removed, row_size, col_size, row_offset, col_offset, transposed)
...
subroutine, public dbcsr_iterator_start(iterator, matrix, shared, dynamic, dynamic_byrows)
...
subroutine, public dbcsr_set(matrix, alpha)
...
DBCSR operations in CP2K.
subroutine, public cp_dbcsr_sm_fm_multiply(matrix, fm_in, fm_out, ncol, alpha, beta)
multiply a dbcsr with a fm matrix
subroutine, public copy_dbcsr_to_fm(matrix, fm)
Copy a DBCSR matrix to a BLACS matrix.
subroutine, public copy_fm_to_dbcsr(fm, matrix, keep_sparsity)
Copy a BLACS matrix to a dbcsr matrix.
Basic linear algebra operations for full matrices.
subroutine, public cp_fm_scale_and_add(alpha, matrix_a, beta, matrix_b)
calc A <- alpha*A + beta*B optimized for alpha == 1.0 (just add beta*B) and beta == 0....
subroutine, public cp_fm_scale(alpha, matrix_a)
scales a matrix matrix_a = alpha * matrix_b
used for collecting some of the diagonalization schemes available for cp_fm_type. cp_fm_power also mo...
subroutine, public cp_fm_geeig(amatrix, bmatrix, eigenvectors, eigenvalues, work)
General Eigenvalue Problem AX = BXE. Use cuSOLVERMp directly when requested and large enough; otherwi...
subroutine, public choose_eigv_solver(matrix, eigenvectors, eigenvalues, info)
Choose the Eigensolver depending on which library is available ELPA seems to be unstable for small sy...
represent the structure of a full matrix
subroutine, public cp_fm_struct_create(fmstruct, para_env, context, nrow_global, ncol_global, nrow_block, ncol_block, descriptor, first_p_pos, local_leading_dimension, template_fmstruct, square_blocks, force_block)
allocates and initializes a full matrix structure
subroutine, public cp_fm_struct_release(fmstruct)
releases a full matrix structure
represent a full matrix distributed on many processors
subroutine, public cp_fm_get_info(matrix, name, nrow_global, ncol_global, nrow_block, ncol_block, nrow_local, ncol_local, row_indices, col_indices, local_data, context, nrow_locals, ncol_locals, matrix_struct, para_env)
returns all kind of information about the full matrix
subroutine, public cp_fm_create(matrix, matrix_struct, name, use_sp, nrow, ncol, set_zero)
creates a new full matrix with the given structure
subroutine, public cp_fm_vectorsnorm(matrix, norm_array)
find the inorm of each column norm_{j}= sqrt( \sum_{i} A_{ij}*A_{ij} )
subroutine, public cp_fm_to_fm_submat_general(source, destination, nrows, ncols, s_firstrow, s_firstcol, d_firstrow, d_firstcol, global_context)
General copy of a submatrix of fm matrix to a submatrix of another fm matrix. The two matrices can ha...
subroutine, public cp_fm_set_all(matrix, alpha, beta)
set all elements of a matrix to the same value, and optionally the diagonal to a different one
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)
...
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,...
integer, parameter, public cp_p_file
integer function, public cp_print_key_should_output(iteration_info, basis_section, print_key_path, used_print_key, first_time)
returns what should be done with the given property if btest(res,cp_p_store) then the property should...
A wrapper around pw_to_cube() which accepts particle_list_type.
subroutine, public cp_pw_to_cube(pw, unit_nr, title, particles, zeff, stride, max_file_size_mb, zero_tails, silent, mpi_io)
...
Defines the basic variable types.
integer, parameter, public int_8
integer, parameter, public dp
integer, parameter, public default_path_length
Definition of mathematical constants and functions.
complex(kind=dp), parameter, public z_one
real(kind=dp), parameter, public twopi
complex(kind=dp), parameter, public z_zero
Interface to the message passing library MPI.
Functions handling the MOLDEN format. Split from mode_selective.
subroutine, public write_mos_molden(mos, qs_kind_set, particle_set, print_section, cell, unoccupied_orbs, unoccupied_evals, qs_env, calc_energies)
Write out the MOs in molden format for visualisation.
Calculates the moment integrals <a|r^m|b>
subroutine, public get_reference_point(rpoint, drpoint, qs_env, fist_env, reference, ref_point, ifirst, ilast)
...
basic linear algebra operations for full matrixes
represent a simple array based list of the given type
Define the data structure for the particle information.
Definition of physical constants:
real(kind=dp), parameter, public evolt
container for various plainwaves related things
subroutine, public pw_env_get(pw_env, pw_pools, cube_info, gridlevel_info, auxbas_pw_pool, auxbas_grid, auxbas_rs_desc, auxbas_rs_grid, rs_descs, rs_grids, xc_pw_pool, vdw_pw_pool, poisson_env, interp_section)
returns the various attributes of the pw env
functions related to the poisson solver on regular grids
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
Calculate the plane wave density by collocating the primitive Gaussian functions (pgf).
subroutine, public calculate_wavefunction(mo_vectors, ivector, rho, rho_gspace, atomic_kind_set, qs_kind_set, cell, dft_control, particle_set, pw_env, basis_type)
maps a given wavefunction on the grid
subroutine, public get_qs_env(qs_env, atomic_kind_set, qs_kind_set, cell, super_cell, cell_ref, use_ref_cell, kpoints, dft_control, mos, sab_orb, sab_all, qmmm, qmmm_periodic, mimic, sac_ae, sac_ppl, sac_lri, sap_ppnl, sab_vdw, sab_scp, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_pp, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, sab_cneo, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, xcint_weights, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, rhoz_cneo_set, ecoul_1c, rho0_s_rs, rho0_s_gs, rhoz_cneo_s_rs, rhoz_cneo_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, harris_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, wanniercentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs, do_rixs, tb_tblite)
Get the QUICKSTEP environment.
Define the quickstep kind type and their sub types.
Definition and initialisation of the mo data type.
subroutine, public set_mo_set(mo_set, maxocc, homo, lfomo, nao, nelectron, n_el_f, nmo, eigenvalues, occupation_numbers, uniform_occupation, kts, mu, flexible_electron_count)
Set the components of a MO set data structure.
subroutine, public init_mo_set(mo_set, fm_pool, fm_ref, fm_struct, name, counter)
initializes an allocated mo_set. eigenvalues, mo_coeff, occupation_numbers are valid only after this ...
subroutine, public allocate_mo_set(mo_set, nao, nmo, nelectron, n_el_f, maxocc, flexible_electron_count)
Allocates a mo set and partially initializes it (nao,nmo,nelectron, and flexible_electron_count are v...
subroutine, public deallocate_mo_set(mo_set)
Deallocate a wavefunction data structure.
subroutine, public get_mo_set(mo_set, maxocc, homo, lfomo, nao, nelectron, n_el_f, nmo, eigenvalues, occupation_numbers, mo_coeff, mo_coeff_b, uniform_occupation, kts, mu, flexible_electron_count)
Get the components of a MO set data structure.
Calculates the moment integrals <a|r^m|b> and <a|r x d/dr|b>
subroutine, public build_berry_moment_matrix(qs_env, cosmat, sinmat, kvec, sab_orb_external, basis_type)
...
Define the neighbor list data types and the corresponding functionality.
subroutine, public rrc_xyz_ao(op, qs_env, rc, order, minimum_image, soft)
Calculation of the components of the dipole operator in the length form by taking the relative positi...
Calculation of overlap matrix, its derivatives and forces.
subroutine, public build_overlap_matrix(ks_env, matrix_s, matrixkp_s, matrix_name, nderivative, basis_type_a, basis_type_b, sab_nl, calculate_forces, matrix_p, matrixkp_p, ext_kpoints)
Calculation of the overlap matrix over Cartesian Gaussian functions.
types that represent a quickstep subsys
subroutine, public qs_subsys_get(subsys, atomic_kinds, atomic_kind_set, particles, particle_set, local_particles, molecules, molecule_set, molecule_kinds, molecule_kind_set, local_molecules, para_env, colvar_p, shell_particles, core_particles, gci, multipoles, natom, nparticle, ncore, nshell, nkind, atprop, virial, results, cell, cell_ref, use_ref_cell, energy, force, qs_kind_set, cp_subsys, nelectron_total, nelectron_spin)
...
subroutine, public tddfpt_dipole_operator(dipole_op_mos_occ, tddfpt_control, gs_mos, qs_env)
Compute the action of the dipole operator on the ground state wave function.
subroutine, public tddfpt_print_nto_analysis(qs_env, evects, evals, ostrength, gs_mos, matrix_s, print_section)
Print natural transition orbital analysis.
subroutine, public tddfpt_print_excitation_analysis(log_unit, evects, evals, gs_mos, matrix_s, spinflip, min_amplitude)
Print excitation analysis.
subroutine, public tddfpt_print_exciton_descriptors(log_unit, evects, gs_mos, matrix_s, do_directional_exciton_descriptors, qs_env)
Print exciton descriptors, cf. Mewes et al., JCTC 14, 710-725 (2018)
subroutine, public tddfpt_print_summary(log_unit, evects, evals, gs_mos, ostrength, mult, dipole_op_mos_occ, dipole_form)
Print final TDDFPT excitation energies and oscillator strengths.
Utilities for string manipulations.
subroutine, public integer_to_string(inumber, string)
Converts an integer number to a string. The WRITE statement will return an error message,...
All kind of helpful little routines.
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
Represent a complex full matrix.
keeps the information about the structure of a full matrix
type of a logger, at the moment it contains just a print level starting at which level it should be l...
stores all the informations relevant to an mpi environment
represent a list of objects
contained for different pw related things
environment for the poisson solver
to create arrays of pools
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
Provides all information about a quickstep kind.
calculation environment to calculate the ks matrix, holds all the needed vars. assumes that the core ...
Ground state molecular orbitals.