48 use,
intrinsic :: iso_fortran_env
96 integer,
parameter :: &
97 CASIDA_EPS_DIFF = 1, &
103 integer,
parameter :: &
113 logical :: states_are_real
114 integer,
allocatable :: n_occ(:)
115 integer,
allocatable :: n_unocc(:)
119 integer :: el_per_state
120 character(len=80) :: trandens
121 character(len=80) :: print_exst
122 real(real64) :: weight_thresh
124 logical :: calc_forces
125 logical :: calc_forces_kernel
126 logical :: calc_forces_scf
128 type(restart_t) :: restart_load
129 type(restart_t) :: restart_dump
131 logical,
allocatable :: is_included(:,:,:)
133 type(states_pair_t),
allocatable :: pair(:)
134 integer,
allocatable :: index(:,:,:)
135 integer,
allocatable :: ind(:)
137 real(real64),
allocatable :: dmat(:,:)
138 real(real64),
allocatable :: dmat_save(:,:)
139 complex(real64),
allocatable :: zmat(:,:)
140 complex(real64),
allocatable :: zmat_save(:,:)
141 real(real64),
allocatable :: w(:)
142 real(real64),
allocatable :: dtm(:, :)
143 complex(real64),
allocatable :: ztm(:, :)
144 real(real64),
allocatable :: f(:)
145 real(real64),
allocatable :: s(:)
147 real(real64),
allocatable :: rho(:,:)
148 real(real64),
allocatable :: fxc(:,:,:)
149 real(real64) :: kernel_lrc_alpha
151 real(real64),
allocatable :: fxc_grad(:,:,:,:,:)
152 real(real64),
allocatable :: fxc_grad_spin(:,:,:,:)
154 real(real64),
allocatable :: dmat2(:,:)
155 complex(real64),
allocatable :: zmat2(:,:)
156 real(real64),
allocatable :: dlr_hmat2(:,:)
157 complex(real64),
allocatable :: zlr_hmat2(:,:)
158 real(real64),
allocatable :: forces(:,:,:)
159 real(real64),
allocatable :: dw2(:)
160 real(real64),
allocatable :: zw2(:)
164 real(real64),
allocatable :: qvector(:)
165 real(real64),
allocatable :: qf(:)
166 real(real64),
allocatable :: qf_avg(:)
169 logical :: parallel_in_eh_pairs
170 logical :: parallel_in_domains
171 logical :: distributed_matrix
172 logical :: write_matrix
173 integer :: parallel_solver
174 type(mpi_grp_t) :: mpi_grp
175 logical :: fromScratch
176 logical :: has_photons
178 type(photon_mode_t),
pointer :: photon_modes => null()
180 integer :: n, nb_rows, nb_cols, block_size
181 type(blacs_proc_grid_t) :: proc_grid
182 integer :: desc(BLACS_DLEN)
183 type(MPI_Datatype) :: darray
191 real(real64),
allocatable :: dpot(:)
192 complex(real64),
allocatable :: zpot(:)
198 class(*),
intent(inout) :: system
199 logical,
intent(in) :: from_scratch
205 message(1) =
"CalculationMode = casida not implemented for multi-system calculations"
217 logical,
intent(in) :: fromscratch
221 integer :: idir, theorylevel, iatom, ierr, default_int
222 character(len=100) :: restart_filename
223 logical :: is_frac_occ
229 if (sys%hm%pcm%run_pcm)
then
233 if (sys%space%is_periodic())
then
234 message(1) =
"Casida oscillator strengths will be incorrect in periodic systems."
255 message(1) =
'Info: Starting Casida linear-response calculation.'
267 call gs_restart%end()
269 message(1) =
"Previous gs calculation is required."
273 cas%el_per_state = sys%st%smear%el_per_state
276 cas%space_dim = sys%space%dim
277 safe_allocate(cas%n_occ(1:sys%st%nik))
278 safe_allocate(cas%n_unocc(1:sys%st%nik))
280 call states_elec_count_pairs(sys%st, sys%namespace, cas%n_pairs, cas%n_occ, cas%n_unocc, cas%is_included, is_frac_occ)
281 if (is_frac_occ)
then
286 select case (sys%st%d%ispin)
288 write(
message(1),
'(a,i4,a)')
"Info: Found", cas%n_occ(1),
" occupied states."
289 write(
message(2),
'(a,i4,a)')
"Info: Found", cas%n_unocc(1),
" unoccupied states."
292 write(
message(1),
'(a,i4,a)')
"Info: Found", cas%n_occ(1),
" occupied states with spin up."
293 write(
message(2),
'(a,i4,a)')
"Info: Found", cas%n_unocc(1),
" unoccupied states with spin up."
294 write(
message(3),
'(a,i4,a)')
"Info: Found", cas%n_occ(2),
" occupied states with spin down."
295 write(
message(4),
'(a,i4,a)')
"Info: Found", cas%n_unocc(2),
" unoccupied states with spin down."
301 message(1) =
'Info: Setting up Hamiltonian.'
303 call v_ks_h_setup(sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st, sys%ks, &
304 sys%hm, calc_eigenval=.false.)
345 message(1) =
"Variational and full Casida theory levels do not apply to complex wavefunctions."
352 call parse_variable(sys%namespace,
'EnablePhotons', .false., cas%has_photons)
354 if (cas%has_photons)
then
356 cas%photon_modes => sys%photons%modes
358 write(
message(1),
'(a,i7,a)')
'INFO: Solving Casida equation with ', &
359 cas%photon_modes%nmodes,
' photon modes.'
360 write(
message(2),
'(a)')
'as described in ACS Photonics 2019, 6, 11, 2757-2778.'
362 cas%pt_nmodes = cas%photon_modes%nmodes
377 call parse_variable(sys%namespace,
'CasidaTransitionDensities',
"0", cas%trandens)
379 if (cas%trandens /=
"0")
then
381 sys%outp%how, sys%outp%output_interval)
394 safe_allocate(cas%qvector(1:cas%space_dim))
395 if (
parse_block(sys%namespace,
'CasidaMomentumTransfer', blk) == 0)
then
396 do idir = 1, cas%space_dim
401 message(1) =
"Info: Calculating IXS/EELS transition rates."
417 call parse_variable(sys%namespace,
'CasidaQuadratureOrder', 5, cas%avg_order)
433 call parse_variable(sys%namespace,
'CasidaCalcTriplet', .false., cas%triplet)
435 cas%triplet = .false.
438 if (cas%triplet)
then
439 message(1) =
"Info: Using triplet kernel. Oscillator strengths will be for spin magnetic-dipole field."
452 call parse_variable(sys%namespace,
'CasidaHermitianConjugate', .false., cas%herm_conj)
467 call parse_variable(sys%namespace,
'CasidaDistributedMatrix', .false., cas%distributed_matrix)
468#ifndef HAVE_SCALAPACK
469 if (cas%distributed_matrix)
then
470 message(1) =
"ScaLAPACK layout requested, but code not compiled with ScaLAPACK"
484 call parse_variable(sys%namespace,
'CasidaWriteDistributedMatrix', .false., cas%write_matrix)
485 if (.not. cas%distributed_matrix .and. cas%write_matrix)
then
486 message(1) =
"CasidaWriteDistributedMatrix con only be used with CasidaDistributedMatrix"
504 default_int = solver_elpa
508 call parse_variable(sys%namespace,
'CasidaParallelEigensolver', default_int, cas%parallel_solver)
513 if (cas%distributed_matrix .and. cas%parallel_solver == solver_elpa)
then
514 message(1) =
"ELPA solver requested, but code not compiled with ELPA"
529 call parse_variable(sys%namespace,
'CasidaPrintExcitations',
"all", cas%print_exst)
530 if (cas%distributed_matrix)
then
532 cas%print_exst =
"none"
533 message(1) =
"Using ScaLAPACK layout, thus disabling output of excited states."
534 message(2) =
"This options creates too many files for large Casida matrices."
550 if (cas%weight_thresh >
m_one)
then
551 message(1) =
'Casida coefficients have values between 0 and 1'
552 message(2) =
'Threshold values reset to default value'
554 cas%weight_thresh = -
m_one
564 call parse_variable(sys%namespace,
'CasidaCalcForces', .false., cas%calc_forces)
565 if (cas%calc_forces)
then
575 call parse_variable(sys%namespace,
'CasidaCalcForcesKernel', .
true., cas%calc_forces_kernel)
585 call parse_variable(sys%namespace,
'CasidaCalcForcesSCF', .false., cas%calc_forces_scf)
587 if (cas%distributed_matrix)
then
588 message(1) =
"Info: Forces calculation not compatible with ScaLAPACK layout."
589 message(2) =
"Using normal layout."
591 cas%distributed_matrix = .false.
598 cas%fromScratch = fromscratch
600 if (cas%fromScratch)
then
601 if (cas%triplet)
then
602 call cas%restart_dump%rm(
'kernel_triplet')
604 call cas%restart_dump%rm(
'kernel')
607 if (cas%calc_forces)
then
608 do iatom = 1, sys%ions%natoms
609 do idir = 1, cas%space_dim
610 write(restart_filename,
'(a,i6.6,a,i1)')
'lr_kernel_', iatom,
'_', idir
611 if (cas%triplet) restart_filename = trim(restart_filename)//
'_triplet'
612 call cas%restart_dump%rm(restart_filename)
614 write(restart_filename,
'(a,i6.6,a,i1)')
'lr_hmat1_', iatom,
'_', idir
615 call cas%restart_dump%rm(restart_filename)
622 if (
bitand(theorylevel, casida_eps_diff) /= 0)
then
623 message(1) =
"Info: Approximating resonance energies through KS eigenvalue differences"
625 cas%type = casida_eps_diff
629 if (sys%st%d%ispin /=
spinors)
then
633 message(1) =
"Info: Calculating matrix elements in the Tamm-Dancoff approximation"
641 message(1) =
"Info: Calculating matrix elements with the CV(2)-DFT theory"
648 message(1) =
"Info: Calculating matrix elements with the full Casida method"
657 message(1) =
"Info: Calculating resonance energies via the Petersilka approximation"
674 type(
casida_t),
intent(inout) :: cas
677 integer :: ist, ast, jpair, ik, ierr
679 integer :: np, np_rows, np_cols, ii, info
684 cas%kernel_lrc_alpha = sys%ks%xc%lrc%alpha
686 if (cas%distributed_matrix .and. .not. cas%states_are_real)
then
690 write(
message(1),
'(a,i9)')
"Number of occupied-unoccupied pairs: ", cas%n_pairs
693 if (cas%n_pairs < 1)
then
694 message(1) =
"No Casida pairs -- maybe there are no unoccupied states?"
698 if (
mpi_world%is_root())
write(*,
"(1x)")
702 if (cas%parallel_in_eh_pairs)
then
709 if (cas%distributed_matrix .and. .not. cas%parallel_in_eh_pairs)
then
710 message(1) =
"ScaLAPACK layout requested, but 'Other' parallelization strategy not available."
711 message(2) =
"Please set ParOther to use the ScaLAPACK layout."
712 message(3) =
"Continuing without ScaLAPACK layout."
714 cas%distributed_matrix = .false.
718 cas%n = cas%n_pairs + cas%pt_nmodes
721 if (cas%distributed_matrix)
then
725 np = cas%mpi_grp%size
729 if (mod(np, ii) == 0)
then
735 np_rows = np / np_cols
739 cas%block_size = min(64, cas%n / np_rows)
741 cas%block_size = max(5, cas%block_size)
742 write(
message(1),
'(A,I5,A,I5,A,I5,A)')
'Parallel layout: using block size of ',&
743 cas%block_size,
' and a processor grid with ', np_rows,
'x', np_cols, &
744 ' processors (rows x cols)'
750 cas%nb_rows = numroc(cas%n, cas%block_size, cas%proc_grid%myrow, 0, cas%proc_grid%nprow)
751 cas%nb_cols = numroc(cas%n, cas%block_size, cas%proc_grid%mycol, 0, cas%proc_grid%npcol)
754 call descinit(cas%desc(1), cas%n, cas%n, cas%block_size, cas%block_size, 0, 0, &
755 cas%proc_grid%context, cas%nb_rows, info)
765 safe_allocate(cas%pair(1:cas%n))
766 if (cas%states_are_real)
then
767 safe_allocate( cas%dmat(1:cas%nb_rows, 1:cas%nb_cols))
768 safe_allocate( cas%dtm(1:cas%n, 1:cas%space_dim))
771 safe_allocate( cas%zmat(1:cas%nb_rows, 1:cas%nb_cols))
772 safe_allocate( cas%ztm(1:cas%n, 1:cas%space_dim))
774 safe_allocate( cas%f(1:cas%n))
775 safe_allocate( cas%s(1:cas%n_pairs))
776 safe_allocate( cas%w(1:cas%n))
777 safe_allocate( cas%index(1:maxval(cas%n_occ), cas%nst - maxval(cas%n_unocc) + 1:cas%nst, 1:cas%nik))
778 safe_allocate( cas%ind(1:cas%n))
780 if (cas%calc_forces)
then
781 if (cas%states_are_real)
then
782 safe_allocate(cas%dmat_save(1:cas%n_pairs, 1:cas%n_pairs))
784 safe_allocate(cas%zmat_save(1:cas%n_pairs, 1:cas%n_pairs))
786 safe_allocate(cas%forces(1:cas%space_dim, 1:sys%ions%natoms, 1:cas%n_pairs))
790 safe_allocate( cas%qf (1:cas%n_pairs))
791 safe_allocate( cas%qf_avg(1:cas%n_pairs))
799 do ast = cas%n_occ(ik) + 1, cas%nst
800 do ist = 1, cas%n_occ(ik)
801 if (cas%is_included(ist, ast, ik))
then
802 cas%index(ist, ast, ik) = jpair
803 cas%pair(jpair)%i = ist
804 cas%pair(jpair)%a = ast
805 cas%pair(jpair)%kk = ik
812 if (cas%has_photons)
then
814 do ik = 1, cas%pt_nmodes
815 cas%pair(cas%n_pairs + ik)%i = 1
816 cas%pair(cas%n_pairs + ik)%a = -ik
817 cas%pair(cas%n_pairs + ik)%kk = -ik
821 safe_deallocate_a(cas%is_included)
832 type(
casida_t),
intent(inout) :: cas
836 assert(
allocated(cas%pair))
837 safe_deallocate_a(cas%pair)
838 safe_deallocate_a(cas%index)
839 if (cas%states_are_real)
then
840 safe_deallocate_a(cas%dmat)
841 safe_deallocate_a(cas%dtm)
843 safe_deallocate_a(cas%zmat)
844 safe_deallocate_a(cas%ztm)
846 safe_deallocate_a(cas%s)
847 safe_deallocate_a(cas%f)
848 safe_deallocate_a(cas%w)
849 safe_deallocate_a(cas%ind)
852 safe_deallocate_a(cas%qf)
853 safe_deallocate_a(cas%qf_avg)
856 safe_deallocate_a(cas%n_occ)
857 safe_deallocate_a(cas%n_unocc)
859 if (cas%calc_forces)
then
860 if (cas%states_are_real)
then
861 safe_deallocate_a(cas%dmat_save)
863 safe_deallocate_a(cas%zmat_save)
865 safe_deallocate_a(cas%forces)
868 call cas%restart_dump%end()
869 call cas%restart_load%end()
871 if (cas%distributed_matrix)
then
877 safe_deallocate_a(cas%qvector)
888 type(
casida_t),
intent(inout) :: cas
891 type(
grid_t),
pointer :: gr
893 real(real64),
allocatable :: rho_spin(:, :)
894 real(real64),
allocatable :: fxc_spin(:,:,:)
895 character(len=100) :: restart_filename
900 assert(cas%type >= casida_eps_diff .and. cas%type <=
casida_casida)
907 if (cas%states_are_real)
then
922 if (cas%type /= casida_eps_diff .or. cas%calc_forces)
then
924 safe_allocate(cas%rho(1:gr%np_part, 1:st%d%nspin))
925 safe_allocate(cas%fxc(1:gr%np, 1:st%d%nspin, 1:st%d%nspin))
926 cas%gga =
in_family(sys%ks%xc%kernel_family, [xc_family_gga])
928 safe_allocate(cas%fxc_grad(1:gr%np, 1:gr%der%dim, 1:gr%der%dim, 1:st%d%nspin, 1:st%d%nspin))
930 safe_allocate(cas%fxc_grad_spin(1:gr%np, 1:gr%der%dim, 1:st%d%nspin, 1:st%d%nspin))
932 safe_allocate(cas%fxc_grad_spin(0, 0, 0, 0))
935 safe_allocate(cas%fxc_grad(0, 0, 0, 0, 0))
936 safe_allocate(cas%fxc_grad_spin(0, 0, 0, 0))
940 if (cas%triplet)
then
941 safe_allocate(rho_spin(1:gr%np_part, 1:2))
942 safe_allocate(fxc_spin(1:gr%np, 1:2, 1:2))
944 rho_spin(:, 1) =
m_half * cas%rho(:, 1)
945 rho_spin(:, 2) =
m_half * cas%rho(:, 1)
948 cas%fxc(:, 1, 1) =
m_half * (fxc_spin(:, 1, 1) - fxc_spin(:, 1, 2))
950 safe_deallocate_a(rho_spin)
951 safe_deallocate_a(fxc_spin)
954 call xc_get_fxc(sys%ks%xc, gr, sys%namespace, cas%rho, st%d%ispin, cas%fxc, cas%fxc_grad, &
957 call xc_get_fxc(sys%ks%xc, gr, sys%namespace, cas%rho, st%d%ispin, cas%fxc)
962 if (cas%triplet)
then
963 message(1) =
"Casida calculation with ADSIC not implemented for triplet excitations."
971 restart_filename =
'kernel'
972 if (cas%triplet) restart_filename = trim(restart_filename)//
'_triplet'
974 select case (cas%type)
975 case (casida_eps_diff)
978 if (cas%states_are_real)
then
979 call dcasida_get_matrix(cas, sys%namespace, sys%hm, st, sys%ks, gr, cas%dmat, cas%fxc, &
980 cas%fxc_grad, cas%fxc_grad_spin, restart_filename)
983 call zcasida_get_matrix(cas, sys%namespace, sys%hm, st, sys%ks, gr, cas%zmat, cas%fxc, &
984 cas%fxc_grad, cas%fxc_grad_spin, restart_filename)
990 if (cas%mpi_grp%is_root() .or. cas%distributed_matrix)
then
991 if (cas%states_are_real)
then
998 if (cas%calc_forces)
then
999 if (cas%states_are_real)
then
1006 if (cas%states_are_real)
then
1013 if (cas%type /= casida_eps_diff .or. cas%calc_forces)
then
1014 safe_deallocate_a(cas%fxc)
1015 safe_deallocate_a(cas%fxc_grad)
1016 safe_deallocate_a(cas%fxc_grad_spin)
1017 safe_deallocate_a(cas%rho)
1028 real(real64),
allocatable :: w(:)
1035 do ia = 1, cas%n_pairs
1036 cas%w(ia) = st%eigenval(cas%pair(ia)%a, cas%pair(ia)%kk) - &
1037 st%eigenval(cas%pair(ia)%i, cas%pair(ia)%kk)
1039 message(1) =
"There is a negative unocc-occ KS eigenvalue difference for"
1040 write(
message(2),
'("states ",I5," and ",I5," of k-point ",I5,".")') cas%pair(ia)%i, cas%pair(ia)%a, cas%pair(ia)%kk
1041 message(3) =
"This indicates an inconsistency between gs, unocc, and/or casida calculations."
1047 safe_allocate(w(1:cas%n_pairs))
1049 call sort(w, cas%ind)
1050 safe_deallocate_a(w)
1052 if (
mpi_world%is_root())
write(*,
"(1x)")
1060 real(real64) function casida_matrix_factor(cas, sys)
1064 push_sub(casida_matrix_factor)
1066 casida_matrix_factor =
m_one
1069 casida_matrix_factor =
m_two * casida_matrix_factor
1073 casida_matrix_factor =
m_two * casida_matrix_factor
1076 pop_sub(casida_matrix_factor)
1083 type(namespace_t),
intent(in) :: namespace
1085 integer :: iunit, ia
1087 if (.not. mpi_world%is_root())
return
1091 call io_mkdir(casida_dir, namespace)
1092 iunit = io_open(casida_dir//
'q'//trim(
theory_name(cas)), namespace, action=
'write')
1093 write(iunit,
'(a1,a14,1x,a24,1x,a24,1x,a10,3es15.8,a2)')
'#',
'E' ,
'|<f|exp(iq.r)|i>|^2', &
1094 '<|<f|exp(iq.r)|i>|^2>',
'; q = (',cas%qvector(1:cas%space_dim),
')'
1095 write(iunit,
'(a1,a14,1x,a24,1x,a24,1x,10x,a15)')
'#', trim(units_abbrev(units_out%energy)), &
1100 if (cas%avg_order == 0)
then
1101 do ia = 1, cas%n_pairs
1102 write(iunit,
'(es15.8,es15.8)') units_from_atomic(units_out%energy, cas%w(cas%ind(ia))), cas%qf(cas%ind(ia))
1105 do ia = 1, cas%n_pairs
1106 write(iunit,
'(3es15.8)') units_from_atomic(units_out%energy, cas%w(cas%ind(ia))), &
1107 cas%qf (cas%ind(ia)), &
1108 cas%qf_avg(cas%ind(ia))
1112 call io_close(iunit)
1119 character(len=80) pure function theory_name(cas)
1122 select case (cas%type)
1141 type(states_elec_t),
intent(in) :: st
1142 integer,
intent(in) :: ia
1143 integer,
intent(in) :: jb
1147 isnt_degenerate = (abs((st%eigenval(cas%pair(ia)%a, cas%pair(ia)%kk) - st%eigenval(cas%pair(ia)%i, cas%pair(ia)%kk)) &
1148 - (st%eigenval(cas%pair(jb)%a, cas%pair(jb)%kk) - st%eigenval(cas%pair(jb)%i, cas%pair(jb)%kk))) > 1e-8_real64)
1155 type(
casida_t),
intent(inout) :: cas
1156 integer,
intent(in) :: jb_local
1158 if (.not. cas%distributed_matrix)
then
1161#ifdef HAVE_SCALAPACK
1162 jb = indxl2g(jb_local, cas%block_size, cas%proc_grid%myrow, 0, cas%proc_grid%nprow)
1169 type(
casida_t),
intent(inout) :: cas
1170 integer,
intent(in) :: ia_local
1172 if (.not. cas%distributed_matrix)
then
1175#ifdef HAVE_SCALAPACK
1176 ia = indxl2g(ia_local, cas%block_size, cas%proc_grid%mycol, 0, cas%proc_grid%npcol)
1181 subroutine local_indices(cas, ia, jb, on_this_processor, ia_local, jb_local)
1184 integer,
intent(in) :: ia, jb
1185 logical,
intent(out) :: on_this_processor
1186 integer,
intent(out) :: ia_local, jb_local
1187#ifdef HAVE_SCALAPACK
1188 integer :: ia_proc, jb_proc
1191 if (.not. cas%distributed_matrix)
then
1192 on_this_processor = .
true.
1196#ifdef HAVE_SCALAPACK
1197 ia_proc = indxg2p(ia, cas%block_size, cas%proc_grid%mycol, 0, cas%proc_grid%npcol)
1198 jb_proc = indxg2p(jb, cas%block_size, cas%proc_grid%myrow, 0, cas%proc_grid%nprow)
1199 if (cas%proc_grid%mycol == ia_proc .and. cas%proc_grid%myrow == jb_proc)
then
1201 ia_local = indxg2l(ia, cas%block_size, cas%proc_grid%mycol, 0, cas%proc_grid%npcol)
1202 jb_local = indxg2l(jb, cas%block_size, cas%proc_grid%myrow, 0, cas%proc_grid%nprow)
1204 on_this_processor = .false.
1214#include "casida_inc.F90"
1216#include "complex.F90"
1217#include "casida_inc.F90"
subroutine solve_eps_diff
This is the common interface to a sorting routine. It performs the shell algorithm,...
double floor(double __x) __attribute__((__nothrow__
This module implements batches of mesh functions.
This module provides the BLACS processor grid.
subroutine, public blacs_proc_grid_init(this, mpi_grp, procdim)
Initializes a blacs context from an MPI communicator with topological information.
subroutine, public blacs_proc_grid_end(this)
This module handles the calculation mode.
integer, parameter, public p_strategy_other
something else like e-h pairs
integer, parameter, public p_strategy_domains
parallelization in domains
This module implements the Casida equations for excited states.
integer function get_global_col(cas, ia_local)
integer, parameter solver_scalapack
subroutine zoscillator_strengths(cas, mesh, st)
integer, parameter casida_petersilka
integer, parameter casida_casida
subroutine casida_type_init(cas, sys)
allocates stuff, and constructs the arrays pair_i and pair_j
integer function get_global_row(cas, jb_local)
integer, parameter casida_eps_diff
character(len=80) pure function theory_name(cas)
subroutine dcasida_forces(cas, sys, gr, st)
subroutine dcasida_get_matrix(cas, namespace, hm, st, ks, gr, matrix, fxc, fxc_grad, fxc_grad_spin, restart_file, is_forces)
subroutine local_indices(cas, ia, jb, on_this_processor, ia_local, jb_local)
subroutine zcasida_solve(cas, sys)
subroutine casida_work(sys, cas)
this subroutine calculates electronic excitation energies using the matrix formulation of M....
subroutine zcasida_get_matrix(cas, namespace, hm, st, ks, gr, matrix, fxc, fxc_grad, fxc_grad_spin, restart_file, is_forces)
integer, parameter casida_variational
subroutine qcasida_write(cas, namespace)
subroutine zcasida_write(cas, sys)
logical function isnt_degenerate(cas, st, ia, jb)
subroutine, public casida_run(system, from_scratch)
subroutine doscillator_strengths(cas, mesh, st)
subroutine dcasida_write(cas, sys)
real(real64) function casida_matrix_factor(cas, sys)
subroutine casida_run_legacy(sys, fromScratch)
integer, parameter casida_tamm_dancoff
subroutine zcasida_forces(cas, sys, gr, st)
subroutine casida_type_end(cas)
subroutine dcasida_solve(cas, sys)
This module implements a calculator for the density and defines related functions.
subroutine, public states_elec_total_density(st, mesh, total_rho)
This routine calculates the total electronic density.
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
integer, parameter, public unpolarized
Parameters...
integer, parameter, public spinors
integer, parameter, public spin_polarized
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
integer, parameter, public hartree_fock
integer, parameter, public generalized_kohn_sham_dft
real(real64), parameter, public m_epsilon
real(real64), parameter, public m_half
real(real64), parameter, public m_one
This module implements the underlying real-space grid.
subroutine, public io_function_read_what_how_when(namespace, space, what, how, output_interval, what_tag_in, how_tag_in, output_interval_tag_in, ignore_error)
integer pure function, public kpoints_number(this)
A module to handle KS potential, without the external potential.
integer, parameter, public dft_u_none
System information (time, memory, sysname)
subroutine, public loct_progress_bar(a, maxcount)
A wrapper around the progress bar, such that it can be silenced without needing to dress the call wit...
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
subroutine, public messages_not_implemented(feature, namespace)
character(len=512), private msg
subroutine, public messages_warning(no_lines, all_nodes, namespace)
subroutine, public messages_obsolete_variable(namespace, name, rep)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_input_error(namespace, var, details, row, column)
subroutine, public messages_experimental(name, namespace)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
type(mpi_comm), parameter, public mpi_comm_undefined
used to indicate a communicator has not been initialized
type(mpi_grp_t), public mpi_world
subroutine mpi_grp_init(grp, comm)
Initialize MPI group instance.
This module handles the communicators for the various parallelization strategies.
logical pure function, public multicomm_strategy_is_parallel(mc, level)
This module implements the basic mulsisystem class, a container system for other systems.
integer function, public parse_block(namespace, name, blk, check_varinfo_)
subroutine, public photon_mode_set_n_electrons(this, qtot)
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
integer, parameter, public restart_casida
integer, parameter, public restart_gs
integer, parameter, public restart_type_dump
integer, parameter, public restart_type_load
This module contains interfaces for ScaLAPACK routines Interfaces are from http:
This module is intended to contain "only mathematical" functions and procedures.
pure logical function, public states_are_complex(st)
pure logical function, public states_are_real(st)
This module handles spin dimensions of the states and the k-point distribution.
subroutine, public states_elec_count_pairs(st, namespace, n_pairs, n_occ, n_unocc, is_included, is_frac_occ)
number of occupied-unoccupied pairs for Casida
This module handles reading and writing restart information for the states_elec_t.
subroutine, public states_elec_look_and_load(restart, namespace, space, st, mesh, kpoints, is_complex, packed)
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
This module defines the unit system, used for input and output.
This module is intended to contain simple general-purpose utility functions and procedures.
subroutine, public v_ks_h_setup(namespace, space, gr, ions, ext_partners, st, ks, hm, calc_eigenval, calc_current)
subroutine, public xc_get_fxc(xcs, gr, namespace, rho, ispin, fxc, fxc_grad, fxc_grad_spin)
Returns the exchange-correlation kernel.
logical pure function, public family_is_mgga_with_exc(xcs)
Is the xc function part of the mGGA family with an energy functional.
pure logical function, public in_family(family, xc_families)
subroutine, public xc_write_fxc_info(xcs, iunit, namespace)
subroutine, public xc_sic_write_info(sic, iunit, namespace)
integer, parameter, public sic_adsic
Averaged density SIC.
subroutine, public xc_sic_add_fxc_adsic(namespace, xc, st, gr, rho, fxc)
Adds to fxc the ADSIC contribution.
This class contains all parameters, needed for Casida calculations.
Class describing the electron system.
Description of the grid, containing information on derivatives, stencil, and symmetries.
Container class for lists of system_oct_m::system_t.
The states_elec_t class contains all electronic wave functions.