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
89
90 implicit none
91
92 private
93 public :: casida_run
94
95 integer, parameter :: &
96 CASIDA_EPS_DIFF = 1, &
100 casida_casida = 16
101
102 integer, parameter :: &
103 SOLVER_ELPA = 1, &
105
107 type casida_t
108 private
109 integer :: type
110 ! !< CASIDA_VARIATIONAL | CASIDA_CASIDA
111
112 logical :: states_are_real
113 integer, allocatable :: n_occ(:)
114 integer, allocatable :: n_unocc(:)
115 integer :: nst
116 integer :: nik
117 integer :: space_dim
118 integer :: el_per_state
119 character(len=80) :: trandens
120 character(len=80) :: print_exst
121 real(real64) :: weight_thresh
122 logical :: triplet
123 logical :: calc_forces
124 logical :: calc_forces_kernel
125 logical :: calc_forces_scf
126 logical :: herm_conj
127 type(restart_t) :: restart_load
128 type(restart_t) :: restart_dump
129
130 logical, allocatable :: is_included(:,:,:)
131 integer :: n_pairs
132 type(states_pair_t), allocatable :: pair(:)
133 integer, allocatable :: index(:,:,:)
134 integer, allocatable :: ind(:)
135
136 real(real64), allocatable :: dmat(:,:)
137 real(real64), allocatable :: dmat_save(:,:)
138 complex(real64), allocatable :: zmat(:,:)
139 complex(real64), allocatable :: zmat_save(:,:)
140 real(real64), allocatable :: w(:)
141 real(real64), allocatable :: dtm(:, :)
142 complex(real64), allocatable :: ztm(:, :)
143 real(real64), allocatable :: f(:)
144 real(real64), allocatable :: s(:)
145
146 real(real64), allocatable :: rho(:,:)
147 real(real64), allocatable :: fxc(:,:,:)
148 real(real64) :: kernel_lrc_alpha
149 logical :: gga
150 real(real64), allocatable :: fxc_grad(:,:,:,:,:)
151 real(real64), allocatable :: fxc_grad_spin(:,:,:,:)
152
153 real(real64), allocatable :: dmat2(:,:)
154 complex(real64), allocatable :: zmat2(:,:)
155 real(real64), allocatable :: dlr_hmat2(:,:)
156 complex(real64), allocatable :: zlr_hmat2(:,:)
157 real(real64), allocatable :: forces(:,:,:)
158 real(real64), allocatable :: dw2(:)
159 real(real64), allocatable :: zw2(:)
160
161 ! variables for momentum-transfer-dependent calculation
162 logical :: qcalc
163 real(real64), allocatable :: qvector(:)
164 real(real64), allocatable :: qf(:)
165 real(real64), allocatable :: qf_avg(:)
166 integer :: avg_order
167
168 logical :: parallel_in_eh_pairs
169 logical :: parallel_in_domains
170 logical :: distributed_matrix
171 logical :: write_matrix
172 integer :: parallel_solver
173 type(mpi_grp_t) :: mpi_grp
174 logical :: fromScratch
175 logical :: has_photons
176 integer :: pt_nmodes
177 type(photon_mode_t), pointer :: photon_modes => null()
178
179 integer :: n, nb_rows, nb_cols, block_size
180 type(blacs_proc_grid_t) :: proc_grid
181 integer :: desc(BLACS_DLEN)
182 type(MPI_Datatype) :: darray
183 end type casida_t
184
186 private
187 integer :: qi
188 integer :: qa
189 integer :: qk
190 real(real64), allocatable :: dpot(:)
191 complex(real64), allocatable :: zpot(:)
192 end type casida_save_pot_t
193
194contains
195
196 subroutine casida_run(system, from_scratch)
197 class(*), intent(inout) :: system
198 logical, intent(in) :: from_scratch
199
200 push_sub(casida_run)
201
202 select type (system)
203 class is (multisystem_basic_t)
204 message(1) = "CalculationMode = casida not implemented for multi-system calculations"
205 call messages_fatal(1, namespace=system%namespace)
206 type is (electrons_t)
207 call casida_run_legacy(system, from_scratch)
208 end select
210 pop_sub(casida_run)
211 end subroutine casida_run
213 ! ---------------------------------------------------------
214 subroutine casida_run_legacy(sys, fromScratch)
215 type(electrons_t), intent(inout) :: sys
216 logical, intent(in) :: fromscratch
218 type(casida_t) :: cas
219 type(block_t) :: blk
220 integer :: idir, theorylevel, iatom, ierr, default_int
221 character(len=100) :: restart_filename
222 logical :: is_frac_occ
223 type(restart_t) :: gs_restart
224
226 call profiling_in('CASIDA')
228 if (sys%hm%pcm%run_pcm) then
229 call messages_not_implemented("PCM for CalculationMode /= gs or td", namespace=sys%namespace)
230 end if
232 if (sys%space%is_periodic()) then
233 message(1) = "Casida oscillator strengths will be incorrect in periodic systems."
234 call messages_warning(1, namespace=sys%namespace)
235 end if
237 if (kpoints_number(sys%kpoints) > 1) then
238 ! Hartree matrix elements may not be correct, not tested anyway. --DAS
239 call messages_not_implemented("Casida with k-points", namespace=sys%namespace)
240 end if
241 if (family_is_mgga_with_exc(sys%hm%xc)) then
242 call messages_not_implemented("Casida with MGGA and non-local terms", namespace=sys%namespace)
243 end if
244 if (sys%hm%lda_u_level /= dft_u_none) then
245 call messages_not_implemented("Casida with DFT+U", namespace=sys%namespace)
246 end if
247 if (sys%hm%theory_level == hartree_fock) then
248 call messages_not_implemented("Casida for Hartree-Fock", namespace=sys%namespace)
249 end if
250 if (sys%hm%theory_level == generalized_kohn_sham_dft) then
251 call messages_not_implemented("Casida for generalized Kohn-Sham", namespace=sys%namespace)
252 end if
254 message(1) = 'Info: Starting Casida linear-response calculation.'
255 call messages_info(1, namespace=sys%namespace)
256
257 call gs_restart%init(sys%namespace, restart_gs, restart_type_load, sys%mc, ierr, mesh=sys%gr, exact=.true.)
258 if (ierr == 0) then
259 call states_elec_look_and_load(gs_restart, sys%namespace, sys%space, sys%st, sys%gr, sys%kpoints)
260 call gs_restart%end()
261 else
262 message(1) = "Previous gs calculation is required."
263 call messages_fatal(1, namespace=sys%namespace)
264 end if
266 cas%el_per_state = sys%st%smear%el_per_state
267 cas%nst = sys%st%nst
268 cas%nik = sys%st%nik
269 cas%space_dim = sys%space%dim
270 safe_allocate(cas%n_occ(1:sys%st%nik))
271 safe_allocate(cas%n_unocc(1:sys%st%nik))
273 call states_elec_count_pairs(sys%st, sys%namespace, cas%n_pairs, cas%n_occ, cas%n_unocc, cas%is_included, is_frac_occ)
274 if (is_frac_occ) then
275 call messages_not_implemented("Casida with partial occupations", namespace=sys%namespace)
276 ! Formulas are in Casida 1995 reference. The occupations are not used at all here currently.
277 end if
278
279 select case (sys%st%d%ispin)
281 write(message(1),'(a,i4,a)') "Info: Found", cas%n_occ(1), " occupied states."
282 write(message(2),'(a,i4,a)') "Info: Found", cas%n_unocc(1), " unoccupied states."
283 call messages_info(2, namespace=sys%namespace)
285 write(message(1),'(a,i4,a)') "Info: Found", cas%n_occ(1), " occupied states with spin up."
286 write(message(2),'(a,i4,a)') "Info: Found", cas%n_unocc(1), " unoccupied states with spin up."
287 write(message(3),'(a,i4,a)') "Info: Found", cas%n_occ(2), " occupied states with spin down."
288 write(message(4),'(a,i4,a)') "Info: Found", cas%n_unocc(2), " unoccupied states with spin down."
289 call messages_info(4, namespace=sys%namespace)
290 end select
292
293 ! setup Hamiltonian, without recalculating eigenvalues (use the ones from the restart information)
294 message(1) = 'Info: Setting up Hamiltonian.'
295 call messages_info(1, namespace=sys%namespace)
296 call v_ks_h_setup(sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st, sys%ks, &
297 sys%hm, calc_eigenval=.false.)
298
299 !%Variable CasidaTheoryLevel
300 !%Type flag
301 !%Section Linear Response::Casida
302 !%Default <tt>eps_diff + petersilka + lrtddft_casida</tt>
303 !%Description
304 !% Choose which electron-hole matrix-based theory levels to use in calculating excitation energies.
305 !% More than one may be used to take advantage of the significant commonality between the calculations.
306 !% <tt>variational</tt> and <tt>lrttdft_casida</tt> are not usable with complex wavefunctions.
307 !% Note the restart data saved by each theory level is compatible with all the others.
308 !%Option eps_diff 1
309 !% Difference of eigenvalues, <i>i.e.</i> independent-particle approximation.
310 !%Option petersilka 2
311 !% The Petersilka approximation uses only elements of the Tamm-Dancoff matrix between degenerate
312 !% transitions (if no degeneracy, this is just the diagonal elements). Also called the "single-pole" approximation.
313 !% This is acceptable if there is little mixing between single-particle transitions.
314 !% Ref: M Petersilka, UJ Gossmann, and EKU Gross, <i>Phys. Rev. Lett.</i> <b>76</b>, 1212 (1996);
315 !% T Grabo, M Petersilka,and EKU Gross, <i>Theochem</i> <b>501-502</b> 353 (2000).
316 !%Option tamm_dancoff 4
317 !% The Tamm-Dancoff approximation uses only occupied-unoccupied transitions and not
318 !% unoccupied-occupied transitions.
319 !% Ref: S Hirata and M Head-Gordon, <i>Chem. Phys. Lett.</i> <b>314</b>, 291 (1999).
320 !%Option variational 8
321 !% Second-order constrained variational theory CV(2)-DFT. Only applies to real wavefunctions.
322 !% Ref: T Ziegler, M Seth, M Krykunov, J Autschbach, and F Wang,
323 !% <i>J. Chem. Phys.</i> <b>130</b>, 154102 (2009).
324 !%Option lrtddft_casida 16
325 !% The full Casida method. Only applies to real wavefunctions.
326 !% Ref: C Jamorski, ME Casida, and DR Salahub, <i>J. Chem. Phys.</i> <b>104</b>, 5134 (1996)
327 !% and ME Casida, "Time-dependent density functional response theory for molecules,"
328 !% in <i>Recent Advances in Density Functional Methods</i>, edited by DE Chong, vol. 1
329 !% of <i>Recent Advances in Computational Chemistry</i>, pp. 155-192 (World Scientific,
330 !% Singapore, 1995).
331 !%End
332
333 call parse_variable(sys%namespace, 'CasidaTheoryLevel', casida_eps_diff + casida_petersilka + casida_casida, theorylevel)
334
335 if (states_are_complex(sys%st)) then
336 if ((bitand(theorylevel, casida_variational) /= 0 &
337 .or. bitand(theorylevel, casida_casida) /= 0)) then
338 message(1) = "Variational and full Casida theory levels do not apply to complex wavefunctions."
339 call messages_fatal(1, only_root_writes = .true., namespace=sys%namespace)
340 ! see section II.D of CV(2) paper regarding this assumption. Would be Eq. 30 with complex wfns.
341 end if
342 end if
343
344 ! This variable is documented in xc_oep_init.
345 call parse_variable(sys%namespace, 'EnablePhotons', .false., cas%has_photons)
346 cas%pt_nmodes = 0
347 if (cas%has_photons) then
348 call messages_experimental('EnablePhotons = yes', namespace=sys%namespace)
349 cas%photon_modes => sys%photons%modes
350 call photon_mode_set_n_electrons(cas%photon_modes, sys%st%qtot)
351 write(message(1), '(a,i7,a)') 'INFO: Solving Casida equation with ', &
352 cas%photon_modes%nmodes, ' photon modes.'
353 write(message(2), '(a)') 'as described in ACS Photonics 2019, 6, 11, 2757-2778.'
354 call messages_info(2, namespace=sys%namespace)
355 cas%pt_nmodes = cas%photon_modes%nmodes
356 end if
357
358 !%Variable CasidaTransitionDensities
359 !%Type string
360 !%Section Linear Response::Casida
361 !%Default write none
362 !%Description
363 !% Specifies which transition densities are to be calculated and written down. The
364 !% transition density for the many-body state <i>n</i> will be written to a file called
365 !% <tt>rho_0n</tt> prefixed by the theory level. Format is set by <tt>OutputFormat</tt>.
366 !%
367 !% This variable is a string in list form, <i>i.e.</i> expressions such as "1,2-5,8-15" are
368 !% valid.
369 !%End
370 call parse_variable(sys%namespace, 'CasidaTransitionDensities', "0", cas%trandens)
371
372 if (cas%trandens /= "0") then
373 call io_function_read_what_how_when(sys%namespace, sys%space, sys%outp%what,&
374 sys%outp%how, sys%outp%output_interval)
375 end if
376
377 !%Variable CasidaMomentumTransfer
378 !%Type block
379 !%Section Linear Response::Casida
380 !%Default 0
381 !%Description
382 !% Momentum-transfer vector for the calculation of the dynamic structure
383 !% factor. When this variable is set, the transition rates are determined
384 !% using an exponential operator instead of the normal dipole one.
385 !%End
386
387 safe_allocate(cas%qvector(1:cas%space_dim))
388 if (parse_block(sys%namespace, 'CasidaMomentumTransfer', blk) == 0) then
389 do idir = 1, cas%space_dim
390 call parse_block_float(blk, 0, idir - 1, cas%qvector(idir))
391 end do
392 call parse_block_end(blk)
393 call messages_experimental("IXS/EELS transition rate calculation", namespace=sys%namespace)
394 message(1) = "Info: Calculating IXS/EELS transition rates."
395 call messages_info(1, namespace=sys%namespace)
396 cas%qcalc = .true.
397
398 !%Variable CasidaQuadratureOrder
399 !%Type integer
400 !%Section Linear Response::Casida
401 !%Default 5
402 !%Description
403 !% Only applies if <tt>CasidaMomentumTransfer</tt> is nonzero.
404 !% Directionally averaged dynamic structure factor is calculated by
405 !% averaging over the results from a set of <math>\vec{q}</math>-vectors. The vectors
406 !% are generated using Gauss-Legendre quadrature scheme [see <i>e.g.</i>
407 !% K. Atkinson, <i>J. Austral. Math. Soc.</i> <b>23</b>, 332 (1982)], and this
408 !% variable determines the order of the scheme.
409 !%End
410 call parse_variable(sys%namespace, 'CasidaQuadratureOrder', 5, cas%avg_order)
411 else
412 cas%qvector(:) = m_zero
413 cas%qcalc = .false.
414 end if
415
416 !%Variable CasidaCalcTriplet
417 !%Type logical
418 !%Section Linear Response::Casida
419 !%Default false
420 !%Description
421 !% For a non-spin-polarized ground state, singlet or triplet excitations can be calculated
422 !% using different matrix elements. Default is to calculate singlets. This variable has no
423 !% effect for a spin-polarized calculation.
424 !%End
425 if (sys%st%d%ispin == unpolarized) then
426 call parse_variable(sys%namespace, 'CasidaCalcTriplet', .false., cas%triplet)
427 else
428 cas%triplet = .false.
429 end if
430
431 if (cas%triplet) then
432 message(1) = "Info: Using triplet kernel. Oscillator strengths will be for spin magnetic-dipole field."
433 call messages_info(1, namespace=sys%namespace)
434 end if
435
436 !%Variable CasidaHermitianConjugate
437 !%Type logical
438 !%Section Linear Response::Casida
439 !%Default false
440 !%Description
441 !% The Casida matrix is Hermitian, so it should not matter whether we calculate the upper or
442 !% lower diagonal. Numerical issues may cause small differences however. Use this variable to
443 !% calculate the Hermitian conjugate of the usual matrix, for testing.
444 !%End
445 call parse_variable(sys%namespace, 'CasidaHermitianConjugate', .false., cas%herm_conj)
446
447 !%Variable CasidaDistributedMatrix
448 !%Type logical
449 !%Section Linear Response::Casida
450 !%Default false
451 !%Description
452 !% Large matrices with more than a few thousand rows and columns usually do
453 !% not fit into the memory of one processor anymore. With this option, the
454 !% Casida matrix is distributed in block-cyclic fashion over all cores in the
455 !% ParOther group. The diagonalization is done in parallel using ScaLAPACK
456 !% or ELPA, if available. For very large matrices (>100000), only the
457 !% ParOther strategy should be used because the diagonalization dominates
458 !% the run time of the computation.
459 !%End
460 call parse_variable(sys%namespace, 'CasidaDistributedMatrix', .false., cas%distributed_matrix)
461#ifndef HAVE_SCALAPACK
462 if (cas%distributed_matrix) then
463 message(1) = "ScaLAPACK layout requested, but code not compiled with ScaLAPACK"
464 call messages_fatal(1, namespace=sys%namespace)
465 end if
466#endif
467 call messages_obsolete_variable(sys%namespace, 'CasidaUseScalapackLayout', 'CasidaDistributedMatrix')
468
469 !%Variable CasidaWriteDistributedMatrix
470 !%Type logical
471 !%Section Linear Response::Casida
472 !%Default false
473 !%Description
474 !% Set to true to write out the full distributed Casida matrix to a file
475 !% using MPI-IO.
476 !%End
477 call parse_variable(sys%namespace, 'CasidaWriteDistributedMatrix', .false., cas%write_matrix)
478 if (.not. cas%distributed_matrix .and. cas%write_matrix) then
479 message(1) = "CasidaWriteDistributedMatrix con only be used with CasidaDistributedMatrix"
480 call messages_fatal(1, namespace=sys%namespace)
481 end if
482
483 !%Variable CasidaParallelEigensolver
484 !%Type integer
485 !%Section Linear Response::Casida
486 !%Description
487 !% Choose library to use for solving the parallel eigenproblem
488 !% of the Casida problem. This options is only relevant if a
489 !% distributed matrix is used (CasidaDistributedMatrix=true).
490 !% By default, elpa is chosen if available.
491 !%Option casida_elpa 1
492 !% Use ELPA library as parallel eigensolver
493 !%Option casida_scalapack 2
494 !% Use Scalapack as parallel eigensolver
495 !%End
496#ifdef HAVE_ELPA
497 default_int = solver_elpa
498#else
499 default_int = solver_scalapack
500#endif
501 call parse_variable(sys%namespace, 'CasidaParallelEigensolver', default_int, cas%parallel_solver)
502 if (.not. varinfo_valid_option('CasidaParallelEigensolver', cas%parallel_solver)) then
503 call messages_input_error(sys%namespace, 'CasidaParallelEigensolver')
504 end if
505#ifndef HAVE_ELPA
506 if (cas%distributed_matrix .and. cas%parallel_solver == solver_elpa) then
507 message(1) = "ELPA solver requested, but code not compiled with ELPA"
508 call messages_fatal(1, namespace=sys%namespace)
509 end if
510#endif
511
512 !%Variable CasidaPrintExcitations
513 !%Type string
514 !%Section Linear Response::Casida
515 !%Default write all
516 !%Description
517 !% Specifies which excitations are written at the end of the calculation.
518 !%
519 !% This variable is a string in list form, <i>i.e.</i> expressions such as "1,2-5,8-15" are
520 !% valid.
521 !%End
522 call parse_variable(sys%namespace, 'CasidaPrintExcitations', "all", cas%print_exst)
523 if (cas%distributed_matrix) then
524 ! do not print excited states -> too many files generated!
525 cas%print_exst = "none"
526 message(1) = "Using ScaLAPACK layout, thus disabling output of excited states."
527 message(2) = "This options creates too many files for large Casida matrices."
528 call messages_info(2, namespace=sys%namespace)
529 end if
530
531 !%Variable CasidaWeightThreshold
532 !%Type float
533 !%Section Linear Response::Casida
534 !%Default -1.
535 !%Description
536 !% Specifies the threshold value for which the individual excitations are printed.
537 !% i.e. juste-h pairs with weight larger than this threshold will be printed.
538 !%
539 !% If a negative value (default) is set, all coefficients will be printed.
540 !% For many case, a 0.01 value is a valid option.
541 !%End
542 call parse_variable(sys%namespace, 'CasidaWeightThreshold', -m_one, cas%weight_thresh)
543 if (cas%weight_thresh > m_one) then
544 message(1) = 'Casida coefficients have values between 0 and 1'
545 message(2) = 'Threshold values reset to default value'
546 call messages_warning(2, namespace=sys%namespace)
547 cas%weight_thresh = -m_one
548 end if
549
550 !%Variable CasidaCalcForces
551 !%Type logical
552 !%Section Linear Response::Casida
553 !%Default false
554 !%Description
555 !% (Experimental) Enable calculation of excited-state forces. Requires previous <tt>vib_modes</tt> calculation.
556 !%End
557 call parse_variable(sys%namespace, 'CasidaCalcForces', .false., cas%calc_forces)
558 if (cas%calc_forces) then
559 call messages_experimental("Excited-state forces calculation", namespace=sys%namespace)
560
561 !%Variable CasidaCalcForcesKernel
562 !%Type logical
563 !%Section Linear Response::Casida
564 !%Default true
565 !%Description
566 !% If false, the derivative of the kernel will not be included in the excited-state force calculation.
567 !%End
568 call parse_variable(sys%namespace, 'CasidaCalcForcesKernel', .true., cas%calc_forces_kernel)
569
570 !%Variable CasidaCalcForcesSCF
571 !%Type logical
572 !%Section Linear Response::Casida
573 !%Default false
574 !%Description
575 !% If true, the ground-state forces will be included in the excited-state forces, so they are total forces.
576 !% If false, the excited-state forces that are produced are only the gradients of the excitation energy.
577 !%End
578 call parse_variable(sys%namespace, 'CasidaCalcForcesSCF', .false., cas%calc_forces_scf)
579
580 if (cas%distributed_matrix) then
581 message(1) = "Info: Forces calculation not compatible with ScaLAPACK layout."
582 message(2) = "Using normal layout."
583 call messages_info(2, namespace=sys%namespace)
584 cas%distributed_matrix = .false.
585 end if
586 end if
587
588 ! Initialize structure
589 call casida_type_init(cas, sys)
590
591 cas%fromScratch = fromscratch
592
593 if (cas%fromScratch) then ! remove old restart files
594 if (cas%triplet) then
595 call cas%restart_dump%rm('kernel_triplet')
596 else
597 call cas%restart_dump%rm('kernel')
598 end if
599
600 if (cas%calc_forces) then
601 do iatom = 1, sys%ions%natoms
602 do idir = 1, cas%space_dim
603 write(restart_filename,'(a,i6.6,a,i1)') 'lr_kernel_', iatom, '_', idir
604 if (cas%triplet) restart_filename = trim(restart_filename)//'_triplet'
605 call cas%restart_dump%rm(restart_filename)
606
607 write(restart_filename,'(a,i6.6,a,i1)') 'lr_hmat1_', iatom, '_', idir
608 call cas%restart_dump%rm(restart_filename)
609 end do
610 end do
611 end if
612 end if
613
614 ! First, print the differences between KS eigenvalues (first approximation to the excitation energies).
615 if (bitand(theorylevel, casida_eps_diff) /= 0) then
616 message(1) = "Info: Approximating resonance energies through KS eigenvalue differences"
617 call messages_info(1, namespace=sys%namespace)
618 cas%type = casida_eps_diff
619 call casida_work(sys, cas)
620 end if
621
622 if (sys%st%d%ispin /= spinors) then
623
624 if (bitand(theorylevel, casida_tamm_dancoff) /= 0) then
625 call messages_experimental("Tamm-Dancoff calculation", namespace=sys%namespace)
626 message(1) = "Info: Calculating matrix elements in the Tamm-Dancoff approximation"
627 call messages_info(1, namespace=sys%namespace)
628 cas%type = casida_tamm_dancoff
629 call casida_work(sys, cas)
630 end if
631
632 if (bitand(theorylevel, casida_variational) /= 0) then
633 call messages_experimental("CV(2)-DFT calculation", namespace=sys%namespace)
634 message(1) = "Info: Calculating matrix elements with the CV(2)-DFT theory"
635 call messages_info(1, namespace=sys%namespace)
636 cas%type = casida_variational
637 call casida_work(sys, cas)
638 end if
639
640 if (bitand(theorylevel, casida_casida) /= 0) then
641 message(1) = "Info: Calculating matrix elements with the full Casida method"
642 call messages_info(1, namespace=sys%namespace)
643 cas%type = casida_casida
644 call casida_work(sys, cas)
645 end if
646
647 ! Doing this first, if doing the others later, takes longer, because we would use
648 ! each Poisson solution for only one matrix element instead of a whole column.
649 if (bitand(theorylevel, casida_petersilka) /= 0) then
650 message(1) = "Info: Calculating resonance energies via the Petersilka approximation"
651 call messages_info(1, namespace=sys%namespace)
652 cas%type = casida_petersilka
653 call casida_work(sys, cas)
654 end if
655
656 end if
657
658 call casida_type_end(cas)
659
660 call profiling_out('CASIDA')
661 pop_sub(casida_run_legacy)
662 end subroutine casida_run_legacy
663
664 ! ---------------------------------------------------------
666 subroutine casida_type_init(cas, sys)
667 type(casida_t), intent(inout) :: cas
668 type(electrons_t), intent(in) :: sys
669
670 integer :: ist, ast, jpair, ik, ierr
671#ifdef HAVE_SCALAPACK
672 integer :: np, np_rows, np_cols, ii, info
673#endif
674
675 push_sub(casida_type_init)
676
677 cas%kernel_lrc_alpha = sys%ks%xc%kernel_lrc_alpha
678 cas%states_are_real = states_are_real(sys%st)
679 if (cas%distributed_matrix .and. .not. cas%states_are_real) then
680 call messages_not_implemented("Complex wavefunctions with ScaLAPACK layout", namespace=sys%namespace)
681 end if
682
683 write(message(1), '(a,i9)') "Number of occupied-unoccupied pairs: ", cas%n_pairs
684 call messages_info(1, namespace=sys%namespace)
685
686 if (cas%n_pairs < 1) then
687 message(1) = "No Casida pairs -- maybe there are no unoccupied states?"
688 call messages_fatal(1, only_root_writes = .true., namespace=sys%namespace)
689 end if
690
691 if (mpi_world%is_root()) write(*, "(1x)")
692
693 ! now let us take care of initializing the parallel stuff
694 cas%parallel_in_eh_pairs = multicomm_strategy_is_parallel(sys%mc, p_strategy_other)
695 if (cas%parallel_in_eh_pairs) then
696 call mpi_grp_init(cas%mpi_grp, sys%mc%group_comm(p_strategy_other))
697 else
698 call mpi_grp_init(cas%mpi_grp, mpi_comm_undefined)
699 end if
700 cas%parallel_in_domains = multicomm_strategy_is_parallel(sys%mc, p_strategy_domains)
701
702 if (cas%distributed_matrix .and. .not. cas%parallel_in_eh_pairs) then
703 message(1) = "ScaLAPACK layout requested, but 'Other' parallelization strategy not available."
704 message(2) = "Please set ParOther to use the ScaLAPACK layout."
705 message(3) = "Continuing without ScaLAPACK layout."
706 call messages_info(3, namespace=sys%namespace)
707 cas%distributed_matrix = .false.
708 end if
709
710 ! dimension of matrix
711 cas%n = cas%n_pairs + cas%pt_nmodes
712
713 ! initialize block-cyclic matrix
714 if (cas%distributed_matrix) then
715#ifdef HAVE_SCALAPACK
716 ! processor layout: always use more processors for rows, this leads to
717 ! better load balancing when computing the matrix elements
718 np = cas%mpi_grp%size
719 np_cols = 1
720 if (np > 3) then
721 do ii = floor(sqrt(real(np))), 2, -1
722 if (mod(np, ii) == 0) then
723 np_cols = ii
724 exit
725 end if
726 end do
727 end if
728 np_rows = np / np_cols
729
730 ! recommended block size: 64, take smaller value for smaller matrices for
731 ! better load balancing
732 cas%block_size = min(64, cas%n / np_rows)
733 ! limit to a minimum block size of 5 for diagonalization efficiency
734 cas%block_size = max(5, cas%block_size)
735 write(message(1), '(A,I5,A,I5,A,I5,A)') 'Parallel layout: using block size of ',&
736 cas%block_size, ' and a processor grid with ', np_rows, 'x', np_cols, &
737 ' processors (rows x cols)'
738 call messages_info(1, namespace=sys%namespace)
739
740 call blacs_proc_grid_init(cas%proc_grid, cas%mpi_grp, procdim = (/np_rows, np_cols/))
741
742 ! get size of local matrices
743 cas%nb_rows = numroc(cas%n, cas%block_size, cas%proc_grid%myrow, 0, cas%proc_grid%nprow)
744 cas%nb_cols = numroc(cas%n, cas%block_size, cas%proc_grid%mycol, 0, cas%proc_grid%npcol)
745
746 ! get ScaLAPACK descriptor
747 call descinit(cas%desc(1), cas%n, cas%n, cas%block_size, cas%block_size, 0, 0, &
748 cas%proc_grid%context, cas%nb_rows, info)
749#endif
750 else
751 ! set to full size
752 cas%nb_rows = cas%n
753 cas%nb_cols = cas%n
754 end if
755
756
757 ! allocate stuff
758 safe_allocate(cas%pair(1:cas%n))
759 if (cas%states_are_real) then
760 safe_allocate( cas%dmat(1:cas%nb_rows, 1:cas%nb_cols))
761 safe_allocate( cas%dtm(1:cas%n, 1:cas%space_dim))
762 else
763 ! caution: ScaLAPACK layout not yet tested for complex wavefunctions!
764 safe_allocate( cas%zmat(1:cas%nb_rows, 1:cas%nb_cols))
765 safe_allocate( cas%ztm(1:cas%n, 1:cas%space_dim))
766 end if
767 safe_allocate( cas%f(1:cas%n))
768 safe_allocate( cas%s(1:cas%n_pairs))
769 safe_allocate( cas%w(1:cas%n))
770 safe_allocate( cas%index(1:maxval(cas%n_occ), cas%nst - maxval(cas%n_unocc) + 1:cas%nst, 1:cas%nik))
771 safe_allocate( cas%ind(1:cas%n))
772
773 if (cas%calc_forces) then
774 if (cas%states_are_real) then
775 safe_allocate(cas%dmat_save(1:cas%n_pairs, 1:cas%n_pairs))
776 else
777 safe_allocate(cas%zmat_save(1:cas%n_pairs, 1:cas%n_pairs))
778 end if
779 safe_allocate(cas%forces(1:cas%space_dim, 1:sys%ions%natoms, 1:cas%n_pairs))
780 end if
781
782 if (cas%qcalc) then
783 safe_allocate( cas%qf (1:cas%n_pairs))
784 safe_allocate( cas%qf_avg(1:cas%n_pairs))
785 end if
786
787 cas%index(:,:,:) = 0
788
789 ! create pairs
790 jpair = 1
791 do ik = 1, cas%nik
792 do ast = cas%n_occ(ik) + 1, cas%nst
793 do ist = 1, cas%n_occ(ik)
794 if (cas%is_included(ist, ast, ik)) then
795 cas%index(ist, ast, ik) = jpair
796 cas%pair(jpair)%i = ist
797 cas%pair(jpair)%a = ast
798 cas%pair(jpair)%kk = ik
799 jpair = jpair + 1
800 end if
801 end do
802 end do
803 end do
804
805 if (cas%has_photons) then
806 ! create pairs for photon modes (negative number refers to photonic excitation)
807 do ik = 1, cas%pt_nmodes
808 cas%pair(cas%n_pairs + ik)%i = 1
809 cas%pair(cas%n_pairs + ik)%a = -ik
810 cas%pair(cas%n_pairs + ik)%kk = -ik
811 end do
812 end if
813
814 safe_deallocate_a(cas%is_included)
815
816 call cas%restart_dump%init(sys%namespace, restart_casida, restart_type_dump, sys%mc, ierr)
817 call cas%restart_load%init(sys%namespace, restart_casida, restart_type_load, sys%mc, ierr)
818
819 pop_sub(casida_type_init)
820 end subroutine casida_type_init
821
822
823 ! ---------------------------------------------------------
824 subroutine casida_type_end(cas)
825 type(casida_t), intent(inout) :: cas
826
827 push_sub(casida_type_end)
828
829 assert(allocated(cas%pair))
830 safe_deallocate_a(cas%pair)
831 safe_deallocate_a(cas%index)
832 if (cas%states_are_real) then
833 safe_deallocate_a(cas%dmat)
834 safe_deallocate_a(cas%dtm)
835 else
836 safe_deallocate_a(cas%zmat)
837 safe_deallocate_a(cas%ztm)
838 end if
839 safe_deallocate_a(cas%s)
840 safe_deallocate_a(cas%f)
841 safe_deallocate_a(cas%w)
842 safe_deallocate_a(cas%ind)
843
844 if (cas%qcalc) then
845 safe_deallocate_a(cas%qf)
846 safe_deallocate_a(cas%qf_avg)
847 end if
848
849 safe_deallocate_a(cas%n_occ)
850 safe_deallocate_a(cas%n_unocc)
851
852 if (cas%calc_forces) then
853 if (cas%states_are_real) then
854 safe_deallocate_a(cas%dmat_save)
855 else
856 safe_deallocate_a(cas%zmat_save)
857 end if
858 safe_deallocate_a(cas%forces)
859 end if
860
861 call cas%restart_dump%end()
862 call cas%restart_load%end()
863
864 if (cas%distributed_matrix) then
865#ifdef HAVE_SCALAPACK
866 call blacs_proc_grid_end(cas%proc_grid)
867#endif
868 end if
869
870 safe_deallocate_a(cas%qvector)
871
872 pop_sub(casida_type_end)
873 end subroutine casida_type_end
874
875
876 ! ---------------------------------------------------------
879 subroutine casida_work(sys, cas)
880 type(electrons_t), target, intent(inout) :: sys
881 type(casida_t), intent(inout) :: cas
882
883 type(states_elec_t), pointer :: st
884 type(grid_t), pointer :: gr
886 real(real64), allocatable :: rho_spin(:, :)
887 real(real64), allocatable :: fxc_spin(:,:,:)
888 character(len=100) :: restart_filename
889
890 push_sub(casida_work)
891
892 ! sanity checks
893 assert(cas%type >= casida_eps_diff .and. cas%type <= casida_casida)
894
895 ! some shortcuts
896 st => sys%st
897 gr => sys%gr
898
899 ! initialize stuff
900 if (cas%states_are_real) then
901 cas%dmat = m_zero
902 cas%dtm = m_zero
903 else
904 cas%zmat = m_zero
905 cas%ztm = m_zero
906 end if
907 cas%f = m_zero
908 cas%w = m_zero
909 cas%s = m_zero
910 if (cas%qcalc) then
911 cas%qf = m_zero
912 cas%qf_avg = m_zero
913 end if
914
915 if (cas%type /= casida_eps_diff .or. cas%calc_forces) then
916 ! We calculate here the kernel, since it will be needed later.
917 safe_allocate(cas%rho(1:gr%np_part, 1:st%d%nspin))
918 safe_allocate(cas%fxc(1:gr%np, 1:st%d%nspin, 1:st%d%nspin))
919 cas%gga = in_family(sys%ks%xc%family, [xc_family_gga])
920 if (cas%gga) then
921 safe_allocate(cas%fxc_grad(1:gr%np, 1:gr%der%dim, 1:gr%der%dim, 1:st%d%nspin, 1:st%d%nspin))
922 if (st%d%ispin == spin_polarized) then
923 safe_allocate(cas%fxc_grad_spin(1:gr%np, 1:gr%der%dim, 1:st%d%nspin, 1:st%d%nspin))
924 else
925 safe_allocate(cas%fxc_grad_spin(0, 0, 0, 0))
926 end if
927 else
928 safe_allocate(cas%fxc_grad(0, 0, 0, 0, 0))
929 safe_allocate(cas%fxc_grad_spin(0, 0, 0, 0))
930 end if
931
932 call states_elec_total_density(st, gr, cas%rho)
933 if (cas%triplet) then
934 safe_allocate(rho_spin(1:gr%np_part, 1:2))
935 safe_allocate(fxc_spin(1:gr%np, 1:2, 1:2))
936
937 rho_spin(:, 1) = m_half * cas%rho(:, 1)
938 rho_spin(:, 2) = m_half * cas%rho(:, 1)
939
940 call xc_get_fxc(sys%ks%xc, gr, sys%namespace, rho_spin, spin_polarized, fxc_spin)
941 cas%fxc(:, 1, 1) = m_half * (fxc_spin(:, 1, 1) - fxc_spin(:, 1, 2))
942
943 safe_deallocate_a(rho_spin)
944 safe_deallocate_a(fxc_spin)
945 else
946 if (cas%gga) then
947 call xc_get_fxc(sys%ks%xc, gr, sys%namespace, cas%rho, st%d%ispin, cas%fxc, cas%fxc_grad, &
948 cas%fxc_grad_spin)
949 else
950 call xc_get_fxc(sys%ks%xc, gr, sys%namespace, cas%rho, st%d%ispin, cas%fxc)
951 end if
952 end if
953
954 if (sys%ks%sic%level == sic_adsic) then
955 call fxc_add_adsic(sys%namespace, sys%ks, st, gr, cas)
956 end if
957
958 end if
959
960 restart_filename = 'kernel'
961 if (cas%triplet) restart_filename = trim(restart_filename)//'_triplet'
962
963 select case (cas%type)
964 case (casida_eps_diff)
965 call solve_eps_diff()
967 if (cas%states_are_real) then
968 call dcasida_get_matrix(cas, sys%namespace, sys%hm, st, sys%ks, gr, cas%dmat, cas%fxc, &
969 cas%fxc_grad, cas%fxc_grad_spin, restart_filename)
970 call dcasida_solve(cas, sys)
971 else
972 call zcasida_get_matrix(cas, sys%namespace, sys%hm, st, sys%ks, gr, cas%zmat, cas%fxc, &
973 cas%fxc_grad, cas%fxc_grad_spin, restart_filename)
974 call zcasida_solve(cas, sys)
975 end if
976 end select
977
978 ! compute oscillator strengths on all processes for the ScaLAPACK layout
979 if (cas%mpi_grp%is_root() .or. cas%distributed_matrix) then
980 if (cas%states_are_real) then
981 call doscillator_strengths(cas, gr, st)
982 else
983 call zoscillator_strengths(cas, gr, st)
984 end if
985 end if
986
987 if (cas%calc_forces) then
988 if (cas%states_are_real) then
989 call dcasida_forces(cas, sys, gr, st)
990 else
991 call zcasida_forces(cas, sys, gr, st)
992 end if
993 end if
994
995 if (cas%states_are_real) then
996 call dcasida_write(cas, sys)
997 else
998 call zcasida_write(cas, sys)
999 end if
1000
1001 ! clean up
1002 if (cas%type /= casida_eps_diff .or. cas%calc_forces) then
1003 safe_deallocate_a(cas%fxc)
1004 safe_deallocate_a(cas%fxc_grad)
1005 safe_deallocate_a(cas%fxc_grad_spin)
1006 safe_deallocate_a(cas%rho)
1007 end if
1008
1009 pop_sub(casida_work)
1010
1011 contains
1012
1013 ! ---------------------------------------------------------
1014 subroutine solve_eps_diff
1015
1016 integer :: ia
1017 real(real64), allocatable :: w(:)
1018
1019 push_sub(casida_work.solve_eps_diff)
1020
1021 ! initialize progress bar
1022 if (mpi_world%is_root()) call loct_progress_bar(-1, cas%n_pairs)
1023
1024 do ia = 1, cas%n_pairs
1025 cas%w(ia) = st%eigenval(cas%pair(ia)%a, cas%pair(ia)%kk) - &
1026 st%eigenval(cas%pair(ia)%i, cas%pair(ia)%kk)
1027 if (cas%w(ia) < -m_epsilon) then
1028 message(1) = "There is a negative unocc-occ KS eigenvalue difference for"
1029 write(message(2),'("states ",I5," and ",I5," of k-point ",I5,".")') cas%pair(ia)%i, cas%pair(ia)%a, cas%pair(ia)%kk
1030 message(3) = "This indicates an inconsistency between gs, unocc, and/or casida calculations."
1031 call messages_fatal(3, only_root_writes = .true., namespace=sys%namespace)
1032 end if
1033 if (mpi_world%is_root()) call loct_progress_bar(ia, cas%n_pairs)
1034 end do
1035
1036 safe_allocate(w(1:cas%n_pairs))
1037 w = cas%w
1038 call sort(w, cas%ind)
1039 safe_deallocate_a(w)
1040
1041 if (mpi_world%is_root()) write(*, "(1x)")
1042
1044 end subroutine solve_eps_diff
1045
1046 ! ---------------------------------------------------------
1047 subroutine fxc_add_adsic(namespace, ks, st, gr, cas)
1048 type(namespace_t), intent(in) :: namespace
1049 type(v_ks_t), intent(in) :: ks
1050 type(states_elec_t), intent(in) :: st
1051 type(grid_t), intent(in) :: gr
1052 type(casida_t), intent(inout) :: cas
1053
1054 real(real64), allocatable :: rho(:, :)
1055 real(real64), allocatable :: fxc_sic(:,:,:)
1056
1057 push_sub(casida_work.fxc_add_adsic)
1058
1059 !Check spin and triplets
1060 if (st%d%ispin /= unpolarized) then
1061 message(1) = "Casida calculation with ADSIC not implemented for spin-polarized calculations."
1062 call messages_fatal(1, namespace=sys%namespace)
1063 end if
1064 if (cas%triplet) then
1065 message(1) = "Casida calculation with ADSIC not implemented for triplet excitations."
1066 call messages_fatal(1, namespace=sys%namespace)
1067 end if
1068 if (cas%gga) then
1069 message(1) = "Casida calculation with ADSIC not implemented with GGAs."
1070 call messages_fatal(1, namespace=sys%namespace)
1071 end if
1072
1073 safe_allocate(fxc_sic(1:gr%np, 1:st%d%nspin, 1:st%d%nspin))
1074 safe_allocate(rho(1:gr%np, 1:st%d%nspin))
1075 rho = cas%rho/st%qtot
1076
1077 call xc_get_fxc(ks%xc, gr, namespace, rho, 1, fxc_sic)
1078
1079 cas%fxc = cas%fxc - fxc_sic/st%qtot
1080
1081 safe_deallocate_a(rho)
1082 safe_deallocate_a(fxc_sic)
1083
1084 pop_sub(casida_work.fxc_add_adsic)
1085 end subroutine fxc_add_adsic
1086
1087 end subroutine casida_work
1088
1089 ! ---------------------------------------------------------
1090 real(real64) function casida_matrix_factor(cas, sys)
1091 type(casida_t), intent(in) :: cas
1092 type(electrons_t), intent(in) :: sys
1093
1094 push_sub(casida_matrix_factor)
1095
1096 casida_matrix_factor = m_one
1097
1098 if (cas%type == casida_variational) then
1099 casida_matrix_factor = m_two * casida_matrix_factor
1100 end if
1101
1102 if (sys%st%d%ispin == unpolarized) then
1103 casida_matrix_factor = m_two * casida_matrix_factor
1104 end if
1105
1106 pop_sub(casida_matrix_factor)
1107
1109
1110 ! ---------------------------------------------------------
1111 subroutine qcasida_write(cas, namespace)
1112 type(casida_t), intent(in) :: cas
1113 type(namespace_t), intent(in) :: namespace
1114
1115 integer :: iunit, ia
1116
1117 if (.not. mpi_world%is_root()) return
1118
1119 push_sub(qcasida_write)
1120
1121 call io_mkdir(casida_dir, namespace)
1122 iunit = io_open(casida_dir//'q'//trim(theory_name(cas)), namespace, action='write')
1123 write(iunit, '(a1,a14,1x,a24,1x,a24,1x,a10,3es15.8,a2)') '#','E' , '|<f|exp(iq.r)|i>|^2', &
1124 '<|<f|exp(iq.r)|i>|^2>','; q = (',cas%qvector(1:cas%space_dim),')'
1125 write(iunit, '(a1,a14,1x,a24,1x,a24,1x,10x,a15)') '#', trim(units_abbrev(units_out%energy)), &
1126 trim('-'), &
1127 trim('-'), &
1128 trim('a.u.')
1129
1130 if (cas%avg_order == 0) then
1131 do ia = 1, cas%n_pairs
1132 write(iunit, '(es15.8,es15.8)') units_from_atomic(units_out%energy, cas%w(cas%ind(ia))), cas%qf(cas%ind(ia))
1133 end do
1134 else
1135 do ia = 1, cas%n_pairs
1136 write(iunit, '(3es15.8)') units_from_atomic(units_out%energy, cas%w(cas%ind(ia))), &
1137 cas%qf (cas%ind(ia)), &
1138 cas%qf_avg(cas%ind(ia))
1139 end do
1140 end if
1141
1142 call io_close(iunit)
1143
1144 pop_sub(qcasida_write)
1145
1146 end subroutine qcasida_write
1147
1148 ! ---------------------------------------------------------
1149 character(len=80) pure function theory_name(cas)
1150 type(casida_t), intent(in) :: cas
1152 select case (cas%type)
1153 case (casida_eps_diff)
1154 theory_name = "eps_diff"
1155 case (casida_petersilka)
1156 theory_name = "petersilka"
1157 case (casida_tamm_dancoff)
1158 theory_name = "tamm_dancoff"
1159 case (casida_variational)
1160 theory_name = "variational"
1161 case (casida_casida)
1162 theory_name = "casida"
1163 case default
1164 theory_name = "unknown"
1165 end select
1166
1167 end function theory_name
1168
1169 logical function isnt_degenerate(cas, st, ia, jb)
1170 type(casida_t), intent(in) :: cas
1171 type(states_elec_t), intent(in) :: st
1172 integer, intent(in) :: ia
1173 integer, intent(in) :: jb
1174
1175 push_sub(isnt_degenerate)
1176
1177 isnt_degenerate = (abs((st%eigenval(cas%pair(ia)%a, cas%pair(ia)%kk) - st%eigenval(cas%pair(ia)%i, cas%pair(ia)%kk)) &
1178 - (st%eigenval(cas%pair(jb)%a, cas%pair(jb)%kk) - st%eigenval(cas%pair(jb)%i, cas%pair(jb)%kk))) > 1e-8_real64)
1179
1180 pop_sub(isnt_degenerate)
1181 end function isnt_degenerate
1182
1183 integer function get_global_row(cas, jb_local) result(jb)
1184 implicit none
1185 type(casida_t), intent(inout) :: cas
1186 integer, intent(in) :: jb_local
1187
1188 if (.not. cas%distributed_matrix) then
1189 jb = jb_local
1190 else
1191#ifdef HAVE_SCALAPACK
1192 jb = indxl2g(jb_local, cas%block_size, cas%proc_grid%myrow, 0, cas%proc_grid%nprow)
1193#endif
1194 end if
1195 end function get_global_row
1196
1197 integer function get_global_col(cas, ia_local) result(ia)
1198 implicit none
1199 type(casida_t), intent(inout) :: cas
1200 integer, intent(in) :: ia_local
1201
1202 if (.not. cas%distributed_matrix) then
1203 ia = ia_local
1204 else
1205#ifdef HAVE_SCALAPACK
1206 ia = indxl2g(ia_local, cas%block_size, cas%proc_grid%mycol, 0, cas%proc_grid%npcol)
1207#endif
1208 end if
1209 end function get_global_col
1211 subroutine local_indices(cas, ia, jb, on_this_processor, ia_local, jb_local)
1212 implicit none
1213 type(casida_t), intent(in) :: cas
1214 integer, intent(in) :: ia, jb
1215 logical, intent(out) :: on_this_processor
1216 integer, intent(out) :: ia_local, jb_local
1217#ifdef HAVE_SCALAPACK
1218 integer :: ia_proc, jb_proc
1219#endif
1220
1221 if (.not. cas%distributed_matrix) then
1222 on_this_processor = .true.
1223 ia_local = ia
1224 jb_local = jb
1225 else
1226#ifdef HAVE_SCALAPACK
1227 ia_proc = indxg2p(ia, cas%block_size, cas%proc_grid%mycol, 0, cas%proc_grid%npcol)
1228 jb_proc = indxg2p(jb, cas%block_size, cas%proc_grid%myrow, 0, cas%proc_grid%nprow)
1229 if (cas%proc_grid%mycol == ia_proc .and. cas%proc_grid%myrow == jb_proc) then
1230 on_this_processor = .true.
1231 ia_local = indxg2l(ia, cas%block_size, cas%proc_grid%mycol, 0, cas%proc_grid%npcol)
1232 jb_local = indxg2l(jb, cas%block_size, cas%proc_grid%myrow, 0, cas%proc_grid%nprow)
1233 else
1234 on_this_processor = .false.
1235 ia_local = -1
1236 jb_local = -1
1237 end if
1238#endif
1239 end if
1240 end subroutine local_indices
1241
1242#include "undef.F90"
1243#include "real.F90"
1244#include "casida_inc.F90"
1245#include "undef.F90"
1246#include "complex.F90"
1247#include "casida_inc.F90"
1248
1249end module casida_oct_m
1250
1251!! Local Variables:
1252!! mode: f90
1253!! coding: utf-8
1254!! End:
subroutine fxc_add_adsic(namespace, ks, st, gr, cas)
Definition: casida.F90:1109
subroutine solve_eps_diff
Definition: casida.F90:1076
This is the common interface to a sorting routine. It performs the shell algorithm,...
Definition: sort.F90:151
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:1259
integer, parameter solver_scalapack
Definition: casida.F90:197
subroutine zoscillator_strengths(cas, mesh, st)
Definition: casida.F90:2966
integer, parameter casida_petersilka
Definition: casida.F90:190
integer, parameter casida_casida
Definition: casida.F90:190
subroutine casida_type_init(cas, sys)
allocates stuff, and constructs the arrays pair_i and pair_j
Definition: casida.F90:762
integer function get_global_row(cas, jb_local)
Definition: casida.F90:1245
integer, parameter casida_eps_diff
Definition: casida.F90:190
character(len=80) pure function theory_name(cas)
Definition: casida.F90:1211
subroutine dcasida_forces(cas, sys, gr, st)
Definition: casida.F90:2252
subroutine dcasida_get_matrix(cas, namespace, hm, st, ks, gr, matrix, fxc, fxc_grad, fxc_grad_spin, restart_file, is_forces)
Definition: casida.F90:1775
subroutine local_indices(cas, ia, jb, on_this_processor, ia_local, jb_local)
Definition: casida.F90:1273
subroutine zcasida_solve(cas, sys)
Definition: casida.F90:4118
subroutine casida_work(sys, cas)
this subroutine calculates electronic excitation energies using the matrix formulation of M....
Definition: casida.F90:941
subroutine zcasida_get_matrix(cas, namespace, hm, st, ks, gr, matrix, fxc, fxc_grad, fxc_grad_spin, restart_file, is_forces)
Definition: casida.F90:3380
integer, parameter casida_variational
Definition: casida.F90:190
subroutine qcasida_write(cas, namespace)
Definition: casida.F90:1173
subroutine zcasida_write(cas, sys)
Definition: casida.F90:4255
logical function isnt_degenerate(cas, st, ia, jb)
Definition: casida.F90:1231
subroutine, public casida_run(system, from_scratch)
Definition: casida.F90:292
subroutine doscillator_strengths(cas, mesh, st)
Definition: casida.F90:1361
subroutine dcasida_write(cas, sys)
Definition: casida.F90:2650
real(real64) function casida_matrix_factor(cas, sys)
Definition: casida.F90:1152
subroutine casida_run_legacy(sys, fromScratch)
Definition: casida.F90:310
integer, parameter casida_tamm_dancoff
Definition: casida.F90:190
subroutine zcasida_forces(cas, sys, gr, st)
Definition: casida.F90:3857
subroutine casida_type_end(cas)
Definition: casida.F90:886
subroutine dcasida_solve(cas, sys)
Definition: casida.F90:2513
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:853
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:192
real(real64), parameter, public m_zero
Definition: global.F90:190
real(real64), parameter, public m_epsilon
Definition: global.F90:206
real(real64), parameter, public m_half
Definition: global.F90:196
real(real64), parameter, public m_one
Definition: global.F90:191
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:1101
A module to handle KS potential, without the external potential.
integer, parameter, public hartree_fock
integer, parameter, public generalized_kohn_sham_dft
integer, parameter, public dft_u_none
Definition: lda_u.F90:203
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_not_implemented(feature, namespace)
Definition: messages.F90:1097
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:531
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1029
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:416
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:697
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1069
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:600
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:310
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:621
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:625
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:694
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:586
subroutine, public xc_get_fxc(xcs, gr, namespace, rho, ispin, fxc, fxc_grad, fxc_grad_spin)
Returns the exchange-correlation kernel.
Definition: xc.F90:1977
pure logical function, public in_family(family, xc_families)
Definition: xc.F90:614
integer, parameter, public sic_adsic
Averaged density SIC.
Definition: xc_sic.F90:150
This class contains all parameters, needed for Casida calculations.
Definition: casida.F90:202
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)