Octopus
casida.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2012-2013 D. Strubbe
3!! Copyright (C) 2017-2018 J. Flick, S. Ohlmann
4!!
5!! This program is free software; you can redistribute it and/or modify
6!! it under the terms of the GNU General Public License as published by
7!! the Free Software Foundation; either version 2, or (at your option)
8!! any later version.
9!!
10!! This program is distributed in the hope that it will be useful,
11!! but WITHOUT ANY WARRANTY; without even the implied warranty of
12!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13!! GNU General Public License for more details.
14!!
15!! You should have received a copy of the GNU General Public License
16!! along with this program; if not, write to the Free Software
17!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
18!! 02110-1301, USA.
19!!
20
21#include "global.h"
22
24
25module casida_oct_m
26 use batch_oct_m
29 use comm_oct_m
30 use debug_oct_m
35#ifdef HAVE_ELPA
36 use elpa
37#endif
40 use forces_oct_m
42 use global_oct_m
43 use grid_oct_m
45 use io_oct_m
47 use iso_c_binding
48 use, intrinsic :: iso_fortran_env
53 use lda_u_oct_m
54 use loct_oct_m
56 use mesh_oct_m
59 use mpi_oct_m
63 use parser_oct_m
64 use pblas_oct_m
65 use pcm_oct_m
74 use sort_oct_m
75 use space_oct_m
81 use unit_oct_m
83 use utils_oct_m
85 use v_ks_oct_m
86 use xc_oct_m
87 use xc_sic_oct_m
88 use xc_f03_lib_m
90
91 implicit none
92
93 private
94 public :: casida_run
95
96 integer, parameter :: &
97 CASIDA_EPS_DIFF = 1, &
100 casida_variational = 8, &
101 casida_casida = 16
102
103 integer, parameter :: &
104 SOLVER_ELPA = 1, &
106
108 type casida_t
109 private
110 integer :: type
111 ! !< CASIDA_VARIATIONAL | CASIDA_CASIDA
112
113 logical :: states_are_real
114 integer, allocatable :: n_occ(:)
115 integer, allocatable :: n_unocc(:)
116 integer :: nst
117 integer :: nik
118 integer :: space_dim
119 integer :: el_per_state
120 character(len=80) :: trandens
121 character(len=80) :: print_exst
122 real(real64) :: weight_thresh
123 logical :: triplet
124 logical :: calc_forces
125 logical :: calc_forces_kernel
126 logical :: calc_forces_scf
127 logical :: herm_conj
128 type(restart_t) :: restart_load
129 type(restart_t) :: restart_dump
130
131 logical, allocatable :: is_included(:,:,:)
132 integer :: n_pairs
133 type(states_pair_t), allocatable :: pair(:)
134 integer, allocatable :: index(:,:,:)
135 integer, allocatable :: ind(:)
136
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(:)
146
147 real(real64), allocatable :: rho(:,:)
148 real(real64), allocatable :: fxc(:,:,:)
149 real(real64) :: kernel_lrc_alpha
150 logical :: gga
151 real(real64), allocatable :: fxc_grad(:,:,:,:,:)
152 real(real64), allocatable :: fxc_grad_spin(:,:,:,:)
153
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(:)
161
162 ! variables for momentum-transfer-dependent calculation
163 logical :: qcalc
164 real(real64), allocatable :: qvector(:)
165 real(real64), allocatable :: qf(:)
166 real(real64), allocatable :: qf_avg(:)
167 integer :: avg_order
168
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
177 integer :: pt_nmodes
178 type(photon_mode_t), pointer :: photon_modes => null()
179
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
184 end type casida_t
185
187 private
188 integer :: qi
189 integer :: qa
190 integer :: qk
191 real(real64), allocatable :: dpot(:)
192 complex(real64), allocatable :: zpot(:)
193 end type casida_save_pot_t
194
195contains
196
197 subroutine casida_run(system, from_scratch)
198 class(*), intent(inout) :: system
199 logical, intent(in) :: from_scratch
200
201 push_sub(casida_run)
202
203 select type (system)
204 class is (multisystem_basic_t)
205 message(1) = "CalculationMode = casida not implemented for multi-system calculations"
206 call messages_fatal(1, namespace=system%namespace)
207 type is (electrons_t)
208 call casida_run_legacy(system, from_scratch)
209 end select
211 pop_sub(casida_run)
212 end subroutine casida_run
214 ! ---------------------------------------------------------
215 subroutine casida_run_legacy(sys, fromScratch)
216 type(electrons_t), intent(inout) :: sys
217 logical, intent(in) :: fromscratch
219 type(casida_t) :: cas
220 type(block_t) :: blk
221 integer :: idir, theorylevel, iatom, ierr, default_int
222 character(len=100) :: restart_filename
223 logical :: is_frac_occ
224 type(restart_t) :: gs_restart
225
227 call profiling_in('CASIDA')
229 if (sys%hm%pcm%run_pcm) then
230 call messages_not_implemented("PCM for CalculationMode /= gs or td", namespace=sys%namespace)
231 end if
233 if (sys%space%is_periodic()) then
234 message(1) = "Casida oscillator strengths will be incorrect in periodic systems."
235 call messages_warning(1, namespace=sys%namespace)
236 end if
238 if (kpoints_number(sys%kpoints) > 1) then
239 ! Hartree matrix elements may not be correct, not tested anyway. --DAS
240 call messages_not_implemented("Casida with k-points", namespace=sys%namespace)
241 end if
242 if (family_is_mgga_with_exc(sys%hm%xc)) then
243 call messages_not_implemented("Casida with MGGA and non-local terms", namespace=sys%namespace)
244 end if
245 if (sys%hm%lda_u_level /= dft_u_none) then
246 call messages_not_implemented("Casida with DFT+U", namespace=sys%namespace)
247 end if
248 if (sys%hm%theory_level == hartree_fock) then
249 call messages_not_implemented("Casida for Hartree-Fock", namespace=sys%namespace)
250 end if
251 if (sys%hm%theory_level == generalized_kohn_sham_dft) then
252 call messages_not_implemented("Casida for generalized Kohn-Sham", namespace=sys%namespace)
253 end if
255 message(1) = 'Info: Starting Casida linear-response calculation.'
256 call messages_info(1, namespace=sys%namespace)
257
258 call messages_print_with_emphasis(msg="XC Kernel level", namespace=sys%namespace)
259 call xc_write_fxc_info(sys%ks%xc, namespace=sys%namespace)
260 call xc_sic_write_info(sys%ks%sic, namespace=sys%namespace)
261 call messages_print_with_emphasis(namespace=sys%namespace)
263
264 call gs_restart%init(sys%namespace, restart_gs, restart_type_load, sys%mc, ierr, mesh=sys%gr, exact=.true.)
265 if (ierr == 0) then
266 call states_elec_look_and_load(gs_restart, sys%namespace, sys%space, sys%st, sys%gr, sys%kpoints)
267 call gs_restart%end()
268 else
269 message(1) = "Previous gs calculation is required."
270 call messages_fatal(1, namespace=sys%namespace)
271 end if
273 cas%el_per_state = sys%st%smear%el_per_state
274 cas%nst = sys%st%nst
275 cas%nik = sys%st%nik
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))
279
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
282 call messages_not_implemented("Casida with partial occupations", namespace=sys%namespace)
283 ! Formulas are in Casida 1995 reference. The occupations are not used at all here currently.
284 end if
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."
290 call messages_info(2, namespace=sys%namespace)
291 case (spin_polarized)
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."
296 call messages_info(4, namespace=sys%namespace)
297 end select
298
299
300 ! setup Hamiltonian, without recalculating eigenvalues (use the ones from the restart information)
301 message(1) = 'Info: Setting up Hamiltonian.'
302 call messages_info(1, namespace=sys%namespace)
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.)
305
306 !%Variable CasidaTheoryLevel
307 !%Type flag
308 !%Section Linear Response::Casida
309 !%Default <tt>eps_diff + petersilka + lrtddft_casida</tt>
310 !%Description
311 !% Choose which electron-hole matrix-based theory levels to use in calculating excitation energies.
312 !% More than one may be used to take advantage of the significant commonality between the calculations.
313 !% <tt>variational</tt> and <tt>lrttdft_casida</tt> are not usable with complex wavefunctions.
314 !% Note the restart data saved by each theory level is compatible with all the others.
315 !%Option eps_diff 1
316 !% Difference of eigenvalues, <i>i.e.</i> independent-particle approximation.
317 !%Option petersilka 2
318 !% The Petersilka approximation uses only elements of the Tamm-Dancoff matrix between degenerate
319 !% transitions (if no degeneracy, this is just the diagonal elements). Also called the "single-pole" approximation.
320 !% This is acceptable if there is little mixing between single-particle transitions.
321 !% Ref: M Petersilka, UJ Gossmann, and EKU Gross, <i>Phys. Rev. Lett.</i> <b>76</b>, 1212 (1996);
322 !% T Grabo, M Petersilka,and EKU Gross, <i>Theochem</i> <b>501-502</b> 353 (2000).
323 !%Option tamm_dancoff 4
324 !% The Tamm-Dancoff approximation uses only occupied-unoccupied transitions and not
325 !% unoccupied-occupied transitions.
326 !% Ref: S Hirata and M Head-Gordon, <i>Chem. Phys. Lett.</i> <b>314</b>, 291 (1999).
327 !%Option variational 8
328 !% Second-order constrained variational theory CV(2)-DFT. Only applies to real wavefunctions.
329 !% Ref: T Ziegler, M Seth, M Krykunov, J Autschbach, and F Wang,
330 !% <i>J. Chem. Phys.</i> <b>130</b>, 154102 (2009).
331 !%Option lrtddft_casida 16
332 !% The full Casida method. Only applies to real wavefunctions.
333 !% Ref: C Jamorski, ME Casida, and DR Salahub, <i>J. Chem. Phys.</i> <b>104</b>, 5134 (1996)
334 !% and ME Casida, "Time-dependent density functional response theory for molecules,"
335 !% in <i>Recent Advances in Density Functional Methods</i>, edited by DE Chong, vol. 1
336 !% of <i>Recent Advances in Computational Chemistry</i>, pp. 155-192 (World Scientific,
337 !% Singapore, 1995).
338 !%End
339
340 call parse_variable(sys%namespace, 'CasidaTheoryLevel', casida_eps_diff + casida_petersilka + casida_casida, theorylevel)
341
342 if (states_are_complex(sys%st)) then
343 if ((bitand(theorylevel, casida_variational) /= 0 &
344 .or. bitand(theorylevel, casida_casida) /= 0)) then
345 message(1) = "Variational and full Casida theory levels do not apply to complex wavefunctions."
346 call messages_fatal(1, only_root_writes = .true., namespace=sys%namespace)
347 ! see section II.D of CV(2) paper regarding this assumption. Would be Eq. 30 with complex wfns.
348 end if
349 end if
350
351 ! This variable is documented in xc_oep_init.
352 call parse_variable(sys%namespace, 'EnablePhotons', .false., cas%has_photons)
353 cas%pt_nmodes = 0
354 if (cas%has_photons) then
355 call messages_experimental('EnablePhotons = yes', namespace=sys%namespace)
356 cas%photon_modes => sys%photons%modes
357 call photon_mode_set_n_electrons(cas%photon_modes, sys%st%qtot)
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.'
361 call messages_info(2, namespace=sys%namespace)
362 cas%pt_nmodes = cas%photon_modes%nmodes
363 end if
364
365 !%Variable CasidaTransitionDensities
366 !%Type string
367 !%Section Linear Response::Casida
368 !%Default write none
369 !%Description
370 !% Specifies which transition densities are to be calculated and written down. The
371 !% transition density for the many-body state <i>n</i> will be written to a file called
372 !% <tt>rho_0n</tt> prefixed by the theory level. Format is set by <tt>OutputFormat</tt>.
373 !%
374 !% This variable is a string in list form, <i>i.e.</i> expressions such as "1,2-5,8-15" are
375 !% valid.
376 !%End
377 call parse_variable(sys%namespace, 'CasidaTransitionDensities', "0", cas%trandens)
378
379 if (cas%trandens /= "0") then
380 call io_function_read_what_how_when(sys%namespace, sys%space, sys%outp%what,&
381 sys%outp%how, sys%outp%output_interval)
382 end if
383
384 !%Variable CasidaMomentumTransfer
385 !%Type block
386 !%Section Linear Response::Casida
387 !%Default 0
388 !%Description
389 !% Momentum-transfer vector for the calculation of the dynamic structure
390 !% factor. When this variable is set, the transition rates are determined
391 !% using an exponential operator instead of the normal dipole one.
392 !%End
393
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
397 call parse_block_float(blk, 0, idir - 1, cas%qvector(idir))
398 end do
399 call parse_block_end(blk)
400 call messages_experimental("IXS/EELS transition rate calculation", namespace=sys%namespace)
401 message(1) = "Info: Calculating IXS/EELS transition rates."
402 call messages_info(1, namespace=sys%namespace)
403 cas%qcalc = .true.
404
405 !%Variable CasidaQuadratureOrder
406 !%Type integer
407 !%Section Linear Response::Casida
408 !%Default 5
409 !%Description
410 !% Only applies if <tt>CasidaMomentumTransfer</tt> is nonzero.
411 !% Directionally averaged dynamic structure factor is calculated by
412 !% averaging over the results from a set of <math>\vec{q}</math>-vectors. The vectors
413 !% are generated using Gauss-Legendre quadrature scheme [see <i>e.g.</i>
414 !% K. Atkinson, <i>J. Austral. Math. Soc.</i> <b>23</b>, 332 (1982)], and this
415 !% variable determines the order of the scheme.
416 !%End
417 call parse_variable(sys%namespace, 'CasidaQuadratureOrder', 5, cas%avg_order)
418 else
419 cas%qvector(:) = m_zero
420 cas%qcalc = .false.
421 end if
422
423 !%Variable CasidaCalcTriplet
424 !%Type logical
425 !%Section Linear Response::Casida
426 !%Default false
427 !%Description
428 !% For a non-spin-polarized ground state, singlet or triplet excitations can be calculated
429 !% using different matrix elements. Default is to calculate singlets. This variable has no
430 !% effect for a spin-polarized calculation.
431 !%End
432 if (sys%st%d%ispin == unpolarized) then
433 call parse_variable(sys%namespace, 'CasidaCalcTriplet', .false., cas%triplet)
434 else
435 cas%triplet = .false.
436 end if
437
438 if (cas%triplet) then
439 message(1) = "Info: Using triplet kernel. Oscillator strengths will be for spin magnetic-dipole field."
440 call messages_info(1, namespace=sys%namespace)
441 end if
442
443 !%Variable CasidaHermitianConjugate
444 !%Type logical
445 !%Section Linear Response::Casida
446 !%Default false
447 !%Description
448 !% The Casida matrix is Hermitian, so it should not matter whether we calculate the upper or
449 !% lower diagonal. Numerical issues may cause small differences however. Use this variable to
450 !% calculate the Hermitian conjugate of the usual matrix, for testing.
451 !%End
452 call parse_variable(sys%namespace, 'CasidaHermitianConjugate', .false., cas%herm_conj)
453
454 !%Variable CasidaDistributedMatrix
455 !%Type logical
456 !%Section Linear Response::Casida
457 !%Default false
458 !%Description
459 !% Large matrices with more than a few thousand rows and columns usually do
460 !% not fit into the memory of one processor anymore. With this option, the
461 !% Casida matrix is distributed in block-cyclic fashion over all cores in the
462 !% ParOther group. The diagonalization is done in parallel using ScaLAPACK
463 !% or ELPA, if available. For very large matrices (>100000), only the
464 !% ParOther strategy should be used because the diagonalization dominates
465 !% the run time of the computation.
466 !%End
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"
471 call messages_fatal(1, namespace=sys%namespace)
472 end if
473#endif
474 call messages_obsolete_variable(sys%namespace, 'CasidaUseScalapackLayout', 'CasidaDistributedMatrix')
475
476 !%Variable CasidaWriteDistributedMatrix
477 !%Type logical
478 !%Section Linear Response::Casida
479 !%Default false
480 !%Description
481 !% Set to true to write out the full distributed Casida matrix to a file
482 !% using MPI-IO.
483 !%End
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"
487 call messages_fatal(1, namespace=sys%namespace)
488 end if
489
490 !%Variable CasidaParallelEigensolver
491 !%Type integer
492 !%Section Linear Response::Casida
493 !%Description
494 !% Choose library to use for solving the parallel eigenproblem
495 !% of the Casida problem. This options is only relevant if a
496 !% distributed matrix is used (CasidaDistributedMatrix=true).
497 !% By default, elpa is chosen if available.
498 !%Option casida_elpa 1
499 !% Use ELPA library as parallel eigensolver
500 !%Option casida_scalapack 2
501 !% Use Scalapack as parallel eigensolver
502 !%End
503#ifdef HAVE_ELPA
504 default_int = solver_elpa
505#else
506 default_int = solver_scalapack
507#endif
508 call parse_variable(sys%namespace, 'CasidaParallelEigensolver', default_int, cas%parallel_solver)
509 if (.not. varinfo_valid_option('CasidaParallelEigensolver', cas%parallel_solver)) then
510 call messages_input_error(sys%namespace, 'CasidaParallelEigensolver')
511 end if
512#ifndef HAVE_ELPA
513 if (cas%distributed_matrix .and. cas%parallel_solver == solver_elpa) then
514 message(1) = "ELPA solver requested, but code not compiled with ELPA"
515 call messages_fatal(1, namespace=sys%namespace)
516 end if
517#endif
518
519 !%Variable CasidaPrintExcitations
520 !%Type string
521 !%Section Linear Response::Casida
522 !%Default write all
523 !%Description
524 !% Specifies which excitations are written at the end of the calculation.
525 !%
526 !% This variable is a string in list form, <i>i.e.</i> expressions such as "1,2-5,8-15" are
527 !% valid.
528 !%End
529 call parse_variable(sys%namespace, 'CasidaPrintExcitations', "all", cas%print_exst)
530 if (cas%distributed_matrix) then
531 ! do not print excited states -> too many files generated!
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."
535 call messages_info(2, namespace=sys%namespace)
536 end if
537
538 !%Variable CasidaWeightThreshold
539 !%Type float
540 !%Section Linear Response::Casida
541 !%Default -1.
542 !%Description
543 !% Specifies the threshold value for which the individual excitations are printed.
544 !% i.e. juste-h pairs with weight larger than this threshold will be printed.
545 !%
546 !% If a negative value (default) is set, all coefficients will be printed.
547 !% For many case, a 0.01 value is a valid option.
548 !%End
549 call parse_variable(sys%namespace, 'CasidaWeightThreshold', -m_one, cas%weight_thresh)
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'
553 call messages_warning(2, namespace=sys%namespace)
554 cas%weight_thresh = -m_one
555 end if
556
557 !%Variable CasidaCalcForces
558 !%Type logical
559 !%Section Linear Response::Casida
560 !%Default false
561 !%Description
562 !% (Experimental) Enable calculation of excited-state forces. Requires previous <tt>vib_modes</tt> calculation.
563 !%End
564 call parse_variable(sys%namespace, 'CasidaCalcForces', .false., cas%calc_forces)
565 if (cas%calc_forces) then
566 call messages_experimental("Excited-state forces calculation", namespace=sys%namespace)
567
568 !%Variable CasidaCalcForcesKernel
569 !%Type logical
570 !%Section Linear Response::Casida
571 !%Default true
572 !%Description
573 !% If false, the derivative of the kernel will not be included in the excited-state force calculation.
574 !%End
575 call parse_variable(sys%namespace, 'CasidaCalcForcesKernel', .true., cas%calc_forces_kernel)
576
577 !%Variable CasidaCalcForcesSCF
578 !%Type logical
579 !%Section Linear Response::Casida
580 !%Default false
581 !%Description
582 !% If true, the ground-state forces will be included in the excited-state forces, so they are total forces.
583 !% If false, the excited-state forces that are produced are only the gradients of the excitation energy.
584 !%End
585 call parse_variable(sys%namespace, 'CasidaCalcForcesSCF', .false., cas%calc_forces_scf)
586
587 if (cas%distributed_matrix) then
588 message(1) = "Info: Forces calculation not compatible with ScaLAPACK layout."
589 message(2) = "Using normal layout."
590 call messages_info(2, namespace=sys%namespace)
591 cas%distributed_matrix = .false.
592 end if
593 end if
594
595 ! Initialize structure
596 call casida_type_init(cas, sys)
597
598 cas%fromScratch = fromscratch
599
600 if (cas%fromScratch) then ! remove old restart files
601 if (cas%triplet) then
602 call cas%restart_dump%rm('kernel_triplet')
603 else
604 call cas%restart_dump%rm('kernel')
605 end if
606
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)
613
614 write(restart_filename,'(a,i6.6,a,i1)') 'lr_hmat1_', iatom, '_', idir
615 call cas%restart_dump%rm(restart_filename)
616 end do
617 end do
618 end if
619 end if
620
621 ! First, print the differences between KS eigenvalues (first approximation to the excitation energies).
622 if (bitand(theorylevel, casida_eps_diff) /= 0) then
623 message(1) = "Info: Approximating resonance energies through KS eigenvalue differences"
624 call messages_info(1, namespace=sys%namespace)
625 cas%type = casida_eps_diff
626 call casida_work(sys, cas)
627 end if
628
629 if (sys%st%d%ispin /= spinors) then
630
631 if (bitand(theorylevel, casida_tamm_dancoff) /= 0) then
632 call messages_experimental("Tamm-Dancoff calculation", namespace=sys%namespace)
633 message(1) = "Info: Calculating matrix elements in the Tamm-Dancoff approximation"
634 call messages_info(1, namespace=sys%namespace)
635 cas%type = casida_tamm_dancoff
636 call casida_work(sys, cas)
637 end if
638
639 if (bitand(theorylevel, casida_variational) /= 0) then
640 call messages_experimental("CV(2)-DFT calculation", namespace=sys%namespace)
641 message(1) = "Info: Calculating matrix elements with the CV(2)-DFT theory"
642 call messages_info(1, namespace=sys%namespace)
643 cas%type = casida_variational
644 call casida_work(sys, cas)
645 end if
646
647 if (bitand(theorylevel, casida_casida) /= 0) then
648 message(1) = "Info: Calculating matrix elements with the full Casida method"
649 call messages_info(1, namespace=sys%namespace)
650 cas%type = casida_casida
651 call casida_work(sys, cas)
652 end if
653
654 ! Doing this first, if doing the others later, takes longer, because we would use
655 ! each Poisson solution for only one matrix element instead of a whole column.
656 if (bitand(theorylevel, casida_petersilka) /= 0) then
657 message(1) = "Info: Calculating resonance energies via the Petersilka approximation"
658 call messages_info(1, namespace=sys%namespace)
659 cas%type = casida_petersilka
660 call casida_work(sys, cas)
661 end if
662
663 end if
664
665 call casida_type_end(cas)
666
667 call profiling_out('CASIDA')
668 pop_sub(casida_run_legacy)
669 end subroutine casida_run_legacy
670
671 ! ---------------------------------------------------------
673 subroutine casida_type_init(cas, sys)
674 type(casida_t), intent(inout) :: cas
675 type(electrons_t), intent(in) :: sys
676
677 integer :: ist, ast, jpair, ik, ierr
678#ifdef HAVE_SCALAPACK
679 integer :: np, np_rows, np_cols, ii, info
680#endif
681
682 push_sub(casida_type_init)
683
684 cas%kernel_lrc_alpha = sys%ks%xc%lrc%alpha
685 cas%states_are_real = states_are_real(sys%st)
686 if (cas%distributed_matrix .and. .not. cas%states_are_real) then
687 call messages_not_implemented("Complex wavefunctions with ScaLAPACK layout", namespace=sys%namespace)
688 end if
689
690 write(message(1), '(a,i9)') "Number of occupied-unoccupied pairs: ", cas%n_pairs
691 call messages_info(1, namespace=sys%namespace)
692
693 if (cas%n_pairs < 1) then
694 message(1) = "No Casida pairs -- maybe there are no unoccupied states?"
695 call messages_fatal(1, only_root_writes = .true., namespace=sys%namespace)
696 end if
697
698 if (mpi_world%is_root()) write(*, "(1x)")
699
700 ! now let us take care of initializing the parallel stuff
701 cas%parallel_in_eh_pairs = multicomm_strategy_is_parallel(sys%mc, p_strategy_other)
702 if (cas%parallel_in_eh_pairs) then
703 call mpi_grp_init(cas%mpi_grp, sys%mc%group_comm(p_strategy_other))
704 else
705 call mpi_grp_init(cas%mpi_grp, mpi_comm_undefined)
706 end if
707 cas%parallel_in_domains = multicomm_strategy_is_parallel(sys%mc, p_strategy_domains)
708
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."
713 call messages_info(3, namespace=sys%namespace)
714 cas%distributed_matrix = .false.
715 end if
716
717 ! dimension of matrix
718 cas%n = cas%n_pairs + cas%pt_nmodes
719
720 ! initialize block-cyclic matrix
721 if (cas%distributed_matrix) then
722#ifdef HAVE_SCALAPACK
723 ! processor layout: always use more processors for rows, this leads to
724 ! better load balancing when computing the matrix elements
725 np = cas%mpi_grp%size
726 np_cols = 1
727 if (np > 3) then
728 do ii = floor(sqrt(real(np))), 2, -1
729 if (mod(np, ii) == 0) then
730 np_cols = ii
731 exit
732 end if
733 end do
734 end if
735 np_rows = np / np_cols
736
737 ! recommended block size: 64, take smaller value for smaller matrices for
738 ! better load balancing
739 cas%block_size = min(64, cas%n / np_rows)
740 ! limit to a minimum block size of 5 for diagonalization efficiency
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)'
745 call messages_info(1, namespace=sys%namespace)
746
747 call blacs_proc_grid_init(cas%proc_grid, cas%mpi_grp, procdim = (/np_rows, np_cols/))
748
749 ! get size of local matrices
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)
752
753 ! get ScaLAPACK descriptor
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)
756#endif
757 else
758 ! set to full size
759 cas%nb_rows = cas%n
760 cas%nb_cols = cas%n
761 end if
762
763
764 ! allocate stuff
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))
769 else
770 ! caution: ScaLAPACK layout not yet tested for complex wavefunctions!
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))
773 end if
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))
779
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))
783 else
784 safe_allocate(cas%zmat_save(1:cas%n_pairs, 1:cas%n_pairs))
785 end if
786 safe_allocate(cas%forces(1:cas%space_dim, 1:sys%ions%natoms, 1:cas%n_pairs))
787 end if
788
789 if (cas%qcalc) then
790 safe_allocate( cas%qf (1:cas%n_pairs))
791 safe_allocate( cas%qf_avg(1:cas%n_pairs))
792 end if
793
794 cas%index(:,:,:) = 0
795
796 ! create pairs
797 jpair = 1
798 do ik = 1, cas%nik
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
806 jpair = jpair + 1
807 end if
808 end do
809 end do
810 end do
811
812 if (cas%has_photons) then
813 ! create pairs for photon modes (negative number refers to photonic excitation)
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
818 end do
819 end if
820
821 safe_deallocate_a(cas%is_included)
822
823 call cas%restart_dump%init(sys%namespace, restart_casida, restart_type_dump, sys%mc, ierr)
824 call cas%restart_load%init(sys%namespace, restart_casida, restart_type_load, sys%mc, ierr)
825
826 pop_sub(casida_type_init)
827 end subroutine casida_type_init
828
829
830 ! ---------------------------------------------------------
831 subroutine casida_type_end(cas)
832 type(casida_t), intent(inout) :: cas
833
834 push_sub(casida_type_end)
835
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)
842 else
843 safe_deallocate_a(cas%zmat)
844 safe_deallocate_a(cas%ztm)
845 end if
846 safe_deallocate_a(cas%s)
847 safe_deallocate_a(cas%f)
848 safe_deallocate_a(cas%w)
849 safe_deallocate_a(cas%ind)
850
851 if (cas%qcalc) then
852 safe_deallocate_a(cas%qf)
853 safe_deallocate_a(cas%qf_avg)
854 end if
855
856 safe_deallocate_a(cas%n_occ)
857 safe_deallocate_a(cas%n_unocc)
858
859 if (cas%calc_forces) then
860 if (cas%states_are_real) then
861 safe_deallocate_a(cas%dmat_save)
862 else
863 safe_deallocate_a(cas%zmat_save)
864 end if
865 safe_deallocate_a(cas%forces)
866 end if
867
868 call cas%restart_dump%end()
869 call cas%restart_load%end()
870
871 if (cas%distributed_matrix) then
872#ifdef HAVE_SCALAPACK
873 call blacs_proc_grid_end(cas%proc_grid)
874#endif
875 end if
876
877 safe_deallocate_a(cas%qvector)
878
879 pop_sub(casida_type_end)
880 end subroutine casida_type_end
881
882
883 ! ---------------------------------------------------------
886 subroutine casida_work(sys, cas)
887 type(electrons_t), target, intent(inout) :: sys
888 type(casida_t), intent(inout) :: cas
889
890 type(states_elec_t), pointer :: st
891 type(grid_t), pointer :: gr
893 real(real64), allocatable :: rho_spin(:, :)
894 real(real64), allocatable :: fxc_spin(:,:,:)
895 character(len=100) :: restart_filename
896
897 push_sub(casida_work)
898
899 ! sanity checks
900 assert(cas%type >= casida_eps_diff .and. cas%type <= casida_casida)
901
902 ! some shortcuts
903 st => sys%st
904 gr => sys%gr
905
906 ! initialize stuff
907 if (cas%states_are_real) then
908 cas%dmat = m_zero
909 cas%dtm = m_zero
910 else
911 cas%zmat = m_zero
912 cas%ztm = m_zero
913 end if
914 cas%f = m_zero
915 cas%w = m_zero
916 cas%s = m_zero
917 if (cas%qcalc) then
918 cas%qf = m_zero
919 cas%qf_avg = m_zero
920 end if
921
922 if (cas%type /= casida_eps_diff .or. cas%calc_forces) then
923 ! We calculate here the kernel, since it will be needed later.
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])
927 if (cas%gga) then
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))
929 if (st%d%ispin == spin_polarized) then
930 safe_allocate(cas%fxc_grad_spin(1:gr%np, 1:gr%der%dim, 1:st%d%nspin, 1:st%d%nspin))
931 else
932 safe_allocate(cas%fxc_grad_spin(0, 0, 0, 0))
933 end if
934 else
935 safe_allocate(cas%fxc_grad(0, 0, 0, 0, 0))
936 safe_allocate(cas%fxc_grad_spin(0, 0, 0, 0))
937 end if
938
939 call states_elec_total_density(st, gr, cas%rho)
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))
943
944 rho_spin(:, 1) = m_half * cas%rho(:, 1)
945 rho_spin(:, 2) = m_half * cas%rho(:, 1)
946
947 call xc_get_fxc(sys%ks%xc, gr, sys%namespace, rho_spin, spin_polarized, fxc_spin)
948 cas%fxc(:, 1, 1) = m_half * (fxc_spin(:, 1, 1) - fxc_spin(:, 1, 2))
949
950 safe_deallocate_a(rho_spin)
951 safe_deallocate_a(fxc_spin)
952 else
953 if (cas%gga) then
954 call xc_get_fxc(sys%ks%xc, gr, sys%namespace, cas%rho, st%d%ispin, cas%fxc, cas%fxc_grad, &
955 cas%fxc_grad_spin)
956 else
957 call xc_get_fxc(sys%ks%xc, gr, sys%namespace, cas%rho, st%d%ispin, cas%fxc)
958 end if
959 end if
960
961 if (sys%ks%sic%level == sic_adsic) then
962 if (cas%triplet) then
963 message(1) = "Casida calculation with ADSIC not implemented for triplet excitations."
964 call messages_fatal(1, namespace=sys%namespace)
965 end if
966 call xc_sic_add_fxc_adsic(sys%namespace, sys%ks%xc, st, gr, cas%rho, cas%fxc)
967 end if
968
969 end if
970
971 restart_filename = 'kernel'
972 if (cas%triplet) restart_filename = trim(restart_filename)//'_triplet'
973
974 select case (cas%type)
975 case (casida_eps_diff)
976 call solve_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)
981 call dcasida_solve(cas, sys)
982 else
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)
985 call zcasida_solve(cas, sys)
986 end if
987 end select
988
989 ! compute oscillator strengths on all processes for the ScaLAPACK layout
990 if (cas%mpi_grp%is_root() .or. cas%distributed_matrix) then
991 if (cas%states_are_real) then
992 call doscillator_strengths(cas, gr, st)
993 else
994 call zoscillator_strengths(cas, gr, st)
995 end if
996 end if
997
998 if (cas%calc_forces) then
999 if (cas%states_are_real) then
1000 call dcasida_forces(cas, sys, gr, st)
1001 else
1002 call zcasida_forces(cas, sys, gr, st)
1003 end if
1004 end if
1005
1006 if (cas%states_are_real) then
1007 call dcasida_write(cas, sys)
1008 else
1009 call zcasida_write(cas, sys)
1010 end if
1011
1012 ! clean up
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)
1018 end if
1019
1020 pop_sub(casida_work)
1021
1022 contains
1023
1024 ! ---------------------------------------------------------
1025 subroutine solve_eps_diff
1026
1027 integer :: ia
1028 real(real64), allocatable :: w(:)
1029
1030 push_sub(casida_work.solve_eps_diff)
1031
1032 ! initialize progress bar
1033 if (mpi_world%is_root()) call loct_progress_bar(-1, cas%n_pairs)
1034
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)
1038 if (cas%w(ia) < -m_epsilon) then
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."
1042 call messages_fatal(3, only_root_writes = .true., namespace=sys%namespace)
1043 end if
1044 if (mpi_world%is_root()) call loct_progress_bar(ia, cas%n_pairs)
1045 end do
1046
1047 safe_allocate(w(1:cas%n_pairs))
1048 w = cas%w
1049 call sort(w, cas%ind)
1050 safe_deallocate_a(w)
1051
1052 if (mpi_world%is_root()) write(*, "(1x)")
1053
1055 end subroutine solve_eps_diff
1056
1057 end subroutine casida_work
1058
1059 ! ---------------------------------------------------------
1060 real(real64) function casida_matrix_factor(cas, sys)
1061 type(casida_t), intent(in) :: cas
1062 type(electrons_t), intent(in) :: sys
1063
1064 push_sub(casida_matrix_factor)
1065
1066 casida_matrix_factor = m_one
1067
1068 if (cas%type == casida_variational) then
1069 casida_matrix_factor = m_two * casida_matrix_factor
1070 end if
1071
1072 if (sys%st%d%ispin == unpolarized) then
1073 casida_matrix_factor = m_two * casida_matrix_factor
1074 end if
1075
1076 pop_sub(casida_matrix_factor)
1077
1078 end function casida_matrix_factor
1079
1080 ! ---------------------------------------------------------
1081 subroutine qcasida_write(cas, namespace)
1082 type(casida_t), intent(in) :: cas
1083 type(namespace_t), intent(in) :: namespace
1084
1085 integer :: iunit, ia
1087 if (.not. mpi_world%is_root()) return
1088
1089 push_sub(qcasida_write)
1090
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)), &
1096 trim('-'), &
1097 trim('-'), &
1098 trim('a.u.')
1099
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))
1103 end do
1104 else
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))
1109 end do
1110 end if
1111
1112 call io_close(iunit)
1113
1114 pop_sub(qcasida_write)
1115
1116 end subroutine qcasida_write
1117
1118 ! ---------------------------------------------------------
1119 character(len=80) pure function theory_name(cas)
1120 type(casida_t), intent(in) :: cas
1122 select case (cas%type)
1123 case (casida_eps_diff)
1124 theory_name = "eps_diff"
1125 case (casida_petersilka)
1126 theory_name = "petersilka"
1127 case (casida_tamm_dancoff)
1128 theory_name = "tamm_dancoff"
1129 case (casida_variational)
1130 theory_name = "variational"
1131 case (casida_casida)
1132 theory_name = "casida"
1133 case default
1134 theory_name = "unknown"
1135 end select
1136
1137 end function theory_name
1138
1139 logical function isnt_degenerate(cas, st, ia, jb)
1140 type(casida_t), intent(in) :: cas
1141 type(states_elec_t), intent(in) :: st
1142 integer, intent(in) :: ia
1143 integer, intent(in) :: jb
1144
1145 push_sub(isnt_degenerate)
1146
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)
1149
1150 pop_sub(isnt_degenerate)
1151 end function isnt_degenerate
1152
1153 integer function get_global_row(cas, jb_local) result(jb)
1154 implicit none
1155 type(casida_t), intent(inout) :: cas
1156 integer, intent(in) :: jb_local
1157
1158 if (.not. cas%distributed_matrix) then
1159 jb = jb_local
1160 else
1161#ifdef HAVE_SCALAPACK
1162 jb = indxl2g(jb_local, cas%block_size, cas%proc_grid%myrow, 0, cas%proc_grid%nprow)
1163#endif
1164 end if
1165 end function get_global_row
1166
1167 integer function get_global_col(cas, ia_local) result(ia)
1168 implicit none
1169 type(casida_t), intent(inout) :: cas
1170 integer, intent(in) :: ia_local
1171
1172 if (.not. cas%distributed_matrix) then
1173 ia = ia_local
1174 else
1175#ifdef HAVE_SCALAPACK
1176 ia = indxl2g(ia_local, cas%block_size, cas%proc_grid%mycol, 0, cas%proc_grid%npcol)
1177#endif
1178 end if
1179 end function get_global_col
1181 subroutine local_indices(cas, ia, jb, on_this_processor, ia_local, jb_local)
1182 implicit none
1183 type(casida_t), intent(in) :: cas
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
1189#endif
1190
1191 if (.not. cas%distributed_matrix) then
1192 on_this_processor = .true.
1193 ia_local = ia
1194 jb_local = jb
1195 else
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
1200 on_this_processor = .true.
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)
1203 else
1204 on_this_processor = .false.
1205 ia_local = -1
1206 jb_local = -1
1207 end if
1208#endif
1209 end if
1210 end subroutine local_indices
1211
1212#include "undef.F90"
1213#include "real.F90"
1214#include "casida_inc.F90"
1215#include "undef.F90"
1216#include "complex.F90"
1217#include "casida_inc.F90"
1218
1219end module casida_oct_m
1220
1221!! Local Variables:
1222!! mode: f90
1223!! coding: utf-8
1224!! End:
subroutine solve_eps_diff
Definition: casida.F90:1087
This is the common interface to a sorting routine. It performs the shell algorithm,...
Definition: sort.F90:156
double floor(double __x) __attribute__((__nothrow__
This module implements batches of mesh functions.
Definition: batch.F90:135
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.
Definition: casida.F90:120
integer function get_global_col(cas, ia_local)
Definition: casida.F90:1229
integer, parameter solver_scalapack
Definition: casida.F90:198
subroutine zoscillator_strengths(cas, mesh, st)
Definition: casida.F90:2944
integer, parameter casida_petersilka
Definition: casida.F90:191
integer, parameter casida_casida
Definition: casida.F90:191
subroutine casida_type_init(cas, sys)
allocates stuff, and constructs the arrays pair_i and pair_j
Definition: casida.F90:769
integer function get_global_row(cas, jb_local)
Definition: casida.F90:1215
integer, parameter casida_eps_diff
Definition: casida.F90:191
character(len=80) pure function theory_name(cas)
Definition: casida.F90:1181
subroutine dcasida_forces(cas, sys, gr, st)
Definition: casida.F90:2228
subroutine dcasida_get_matrix(cas, namespace, hm, st, ks, gr, matrix, fxc, fxc_grad, fxc_grad_spin, restart_file, is_forces)
Definition: casida.F90:1746
subroutine local_indices(cas, ia, jb, on_this_processor, ia_local, jb_local)
Definition: casida.F90:1243
subroutine zcasida_solve(cas, sys)
Definition: casida.F90:4104
subroutine casida_work(sys, cas)
this subroutine calculates electronic excitation energies using the matrix formulation of M....
Definition: casida.F90:948
subroutine zcasida_get_matrix(cas, namespace, hm, st, ks, gr, matrix, fxc, fxc_grad, fxc_grad_spin, restart_file, is_forces)
Definition: casida.F90:3359
integer, parameter casida_variational
Definition: casida.F90:191
subroutine qcasida_write(cas, namespace)
Definition: casida.F90:1143
subroutine zcasida_write(cas, sys)
Definition: casida.F90:4241
logical function isnt_degenerate(cas, st, ia, jb)
Definition: casida.F90:1201
subroutine, public casida_run(system, from_scratch)
Definition: casida.F90:293
subroutine doscillator_strengths(cas, mesh, st)
Definition: casida.F90:1331
subroutine dcasida_write(cas, sys)
Definition: casida.F90:2628
real(real64) function casida_matrix_factor(cas, sys)
Definition: casida.F90:1122
subroutine casida_run_legacy(sys, fromScratch)
Definition: casida.F90:311
integer, parameter casida_tamm_dancoff
Definition: casida.F90:191
subroutine zcasida_forces(cas, sys, gr, st)
Definition: casida.F90:3841
subroutine casida_type_end(cas)
Definition: casida.F90:893
subroutine dcasida_solve(cas, sys)
Definition: casida.F90:2491
This module implements a calculator for the density and defines related functions.
Definition: density.F90:122
subroutine, public states_elec_total_density(st, mesh, total_rho)
This routine calculates the total electronic density.
Definition: density.F90:857
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
Definition: global.F90:193
real(real64), parameter, public m_zero
Definition: global.F90:191
integer, parameter, public hartree_fock
Definition: global.F90:241
integer, parameter, public generalized_kohn_sham_dft
Definition: global.F90:241
real(real64), parameter, public m_epsilon
Definition: global.F90:207
real(real64), parameter, public m_half
Definition: global.F90:197
real(real64), parameter, public m_one
Definition: global.F90:192
This module implements the underlying real-space grid.
Definition: grid.F90:119
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)
Definition: io.F90:116
integer pure function, public kpoints_number(this)
Definition: kpoints.F90:1130
A module to handle KS potential, without the external potential.
integer, parameter, public dft_u_none
Definition: lda_u.F90:205
System information (time, memory, sysname)
Definition: loct.F90:117
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...
Definition: loct.F90:276
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:897
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1067
character(len=512), private msg
Definition: messages.F90:166
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:524
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:999
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:409
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:690
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1039
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:593
type(mpi_comm), parameter, public mpi_comm_undefined
used to indicate a communicator has not been initialized
Definition: mpi.F90:138
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:272
subroutine mpi_grp_init(grp, comm)
Initialize MPI group instance.
Definition: mpi.F90:341
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
logical pure function, public multicomm_strategy_is_parallel(mc, level)
Definition: multicomm.F90:719
This module implements the basic mulsisystem class, a container system for other systems.
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:615
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.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
integer, parameter, public restart_casida
Definition: restart.F90:156
integer, parameter, public restart_gs
Definition: restart.F90:156
integer, parameter, public restart_type_dump
Definition: restart.F90:183
integer, parameter, public restart_type_load
Definition: restart.F90:183
This module contains interfaces for ScaLAPACK routines Interfaces are from http:
Definition: scalapack.F90:133
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:119
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.
Definition: unit.F90:134
This module defines the unit system, used for input and output.
This module is intended to contain simple general-purpose utility functions and procedures.
Definition: utils.F90:120
subroutine, public v_ks_h_setup(namespace, space, gr, ions, ext_partners, st, ks, hm, calc_eigenval, calc_current)
Definition: v_ks.F90:700
subroutine, public xc_get_fxc(xcs, gr, namespace, rho, ispin, fxc, fxc_grad, fxc_grad_spin)
Returns the exchange-correlation kernel.
Definition: xc_kernel.F90:172
Definition: xc.F90:116
logical pure function, public family_is_mgga_with_exc(xcs)
Is the xc function part of the mGGA family with an energy functional.
Definition: xc.F90:599
pure logical function, public in_family(family, xc_families)
Definition: xc.F90:627
subroutine, public xc_write_fxc_info(xcs, iunit, namespace)
Definition: xc.F90:255
subroutine, public xc_sic_write_info(sic, iunit, namespace)
Definition: xc_sic.F90:259
integer, parameter, public sic_adsic
Averaged density SIC.
Definition: xc_sic.F90:153
subroutine, public xc_sic_add_fxc_adsic(namespace, xc, st, gr, rho, fxc)
Adds to fxc the ADSIC contribution.
Definition: xc_sic.F90:460
This class contains all parameters, needed for Casida calculations.
Definition: casida.F90:203
Class describing the electron system.
Definition: electrons.F90:220
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:171
Container class for lists of system_oct_m::system_t.
The states_elec_t class contains all electronic wave functions.
int true(void)