1
C> @file rt_tddft_cs_tdfock.F
4
C--------------------------------------------------------------------
5
C> Builds time-dependent closed-shell Fock matrix.
7
C> Constructs the time-dependent Fock matrix including
8
C> building Fock matrix from density matrix in AO basis,
9
C> and calculating and adding dipole interation with external uniform E-field.
10
C--------------------------------------------------------------------
11
logical function rt_tddft_cs_tdfock (params, tt, g_zdens_ao,
12
$ energies, g_zfock_ao)
17
#include "mafdecls.fh"
22
#include "matutils.fh"
23
#include "rt_tddft.fh"
27
type(rt_params_t), intent(in) :: params !< struct containing parameters
28
double precision, intent(in) :: tt !< current time
29
integer, intent(in) :: g_zdens_ao !< complex dens mat, ns_ao x ns_ao
33
type(rt_energies_t), intent(out) :: energies !< time-dependent energies
34
integer, intent(in) :: g_zfock_ao !< complex fock mat, ns_ao x ns_ao
38
character(*), parameter :: pname = "rt_tddft_cs_tdfock: "
42
type(rt_vector_t) field
43
type(rt_quad_t) field_grad
44
integer g_zscr(2) !scratch, ns_ao x ns_ao; alpha beta
45
integer g_zscr_mo ! ns_mo x ns_mo dcpl
48
double complex zval1, zval2, zscale
50
integer lvals, ivals, g_zevecs, i, g_ztmp_mo
51
integer lvre, ivre, lvim, ivim
54
call rt_tddft_cs_confirm (params)
58
n = params%ns_ao ! alias
61
if (.not.ga_create(mt_dcpl,n ,n ,"zscr1", 0, 0, g_zscr(1)))
62
$ call errquit ("failed to create zscr1", 0, GA_ERR)
63
if (.not.ga_create(mt_dcpl,n ,n ,"zscr2", 0, 0, g_zscr(2)))
64
$ call errquit ("failed to create zscr2", 0, GA_ERR)
68
C Build new complex fock mat from complex dens mat; this also
69
C calculates energies. Note, the input g_zdens_ao is in AO basis,
70
C and the output g_zfock_ao is also in AO basis.
72
call zfock_cs_build (params, g_zdens_ao, energies, g_zfock_ao)
73
call ga_sync () !XXX needed?
77
C Compute dipole interaction (updates field values inside), and add
78
C to Fock matrix. Since this is closedshell, we have already
79
C checked that all fields only act on the full spin. Thus we just
80
C add the alpha part to the fock matrix.
82
call ga_zero (g_zscr(1))
83
call ga_zero (g_zscr(2))
84
call rt_tddft_calc_excite (params, tt, g_zscr)
85
call ga_add (z1, g_zfock_ao, zn1, g_zscr(1), g_zfock_ao)
90
if (.not.ga_destroy(g_zscr(1)))
91
$ call errquit ("failed to destroy zscr1", 0, GA_ERR)
92
if (.not.ga_destroy(g_zscr(2)))
93
$ call errquit ("failed to destroy zscr2", 0, GA_ERR)
96
CXXX [KAL]: have an option to return false?
97
rt_tddft_cs_tdfock= .true.