TABLE OF CONTENTS


ABINIT/m_rttddft_driver [ Modules ]

[ Top ] [ Modules ]

NAME

  m_rttddft_driver

FUNCTION

  Real-Time Time Dependent DFT (RT-TDDFT)

COPYRIGHT

  Copyright (C) 2021-2024 ABINIT group (FB)
  This file is distributed under the terms of the
  GNU General Public License, see ~abinit/COPYING
  or http://www.gnu.org/copyleft/gpl.txt .

SOURCE

16 #if defined HAVE_CONFIG_H
17 #include "config.h"
18 #endif
19 
20 #include "abi_common.h"
21 
22 module m_rttddft_driver
23 
24  use defs_basis
25  use defs_abitypes,        only: MPI_type
26  use defs_datatypes,       only: pseudopotential_type
27 
28  use m_dtfil,              only: datafiles_type
29  use m_dtset,              only: dataset_type
30  use m_errors,             only: msg_hndl
31  use m_pawang,             only: pawang_type
32  use m_pawrad,             only: pawrad_type
33  use m_pawtab,             only: pawtab_type
34  use m_rttddft_output,     only: rttddft_output
35  use m_rttddft_tdks,       only: tdks_type
36  use m_rttddft_propagate,  only: rttddft_propagate_ele
37  use m_rttddft_properties, only: rttddft_calc_density, &
38                                & rttddft_calc_etot
39  use m_specialmsg,         only: wrtout
40  use m_time,               only: cwtime
41 
42  implicit none
43 
44  private

m_rttddft_driver/rttddft [ Functions ]

[ Top ] [ m_rttddft_driver ] [ Functions ]

NAME

  rttddft

FUNCTION

  Primary routine for real-time TDDFT calculations.
  1) Initialization: create main tdks (Time-Dependent Kohn-Sham) object
    - intialize various important parameters
    - read intial KS orbitals in WFK file
    - compute initial density from KS orbitals
  2) Propagation loop (rttddft_propagate):
    for i = 1, ntime
       - propagate KS orbitals
       - propagate nuclei if requested (Erhenfest dynamics)
       - compute new density
       - Compute and print requested properties
  3) Final printout and finalize

INPUTS

  codvsn = code version
  dtfil <type datafiles_type> = infos about file names
  dtset <type(dataset_type)> = all input variables for this dataset
  mpi_enreg <MPI_type> = MPI-parallelisation information
  pawang <type(pawang_type)> = paw angular mesh and related data
  pawrad(ntypat*usepaw) <type(pawrad_type)> = paw radial mesh and related data
  pawtab(ntypat*usepaw) <type(pawtab_type)> = paw tabulated starting data
  psps <type(pseudopotential_type)> = variables related to pseudopotentials

OUTPUT

SOURCE

 84 subroutine rttddft(codvsn,dtfil,dtset,mpi_enreg,pawang,pawrad,pawtab,psps)
 85 
 86  implicit none
 87 
 88  !Arguments ------------------------------------
 89  !scalars
 90  character(len=8),           intent(in)    :: codvsn
 91  type(dataset_type),         intent(inout) :: dtset
 92  type(datafiles_type),       intent(inout) :: dtfil
 93  type(MPI_type),             intent(inout) :: mpi_enreg
 94  type(pawang_type),          intent(inout) :: pawang
 95  type(pseudopotential_type), intent(inout) :: psps
 96  !arrays
 97  type(pawrad_type), target,  intent(inout) :: pawrad(psps%ntypat*psps%usepaw)
 98  type(pawtab_type), target,  intent(inout) :: pawtab(psps%ntypat*psps%usepaw)
 99 
100  !Local variables-------------------------------
101  !scalars
102  character(len=500)   :: msg
103  integer              :: istep
104  type(tdks_type)      :: tdks
105  !arrays
106  real(dp)             :: cpu, wall, gflops
107 
108 ! ***********************************************************************
109 
110  write(msg,'(7a)') ch10,'---------------------------------------------------------------------------',&
111                  & ch10,'----------------   Starting real-time time dependent DFT   ----------------',&
112                  & ch10,'---------------------------------------------------------------------------',ch10
113  call wrtout(ab_out,msg)
114  if (do_write_log) call wrtout(std_out,msg)
115 
116  !FB: Do not allow RT-TDDFT for now. Still under development
117  write(msg,'(3a)') ch10,'Real-time time dependent DFT is not yet available. Still under development..',ch10
118  call wrtout(ab_out,msg)
119  if (do_write_log) call wrtout(std_out,msg)
120  ABI_ERROR(msg)
121 
122  !** 1) Initialization: create main tdks (Time-Dependent Kohn-Sham) object
123  write(msg,'(3a)') ch10,'---------------------------   Initialization   ----------------------------',ch10
124  call wrtout(ab_out,msg)
125  if (do_write_log) call wrtout(std_out,msg)
126 
127  call tdks%init(codvsn,dtfil,dtset,mpi_enreg,pawang,pawrad,pawtab,psps)
128 
129  !Compute initial electronic density
130  call rttddft_calc_density(dtset,mpi_enreg,psps,tdks)
131 
132  !** 2) Propagation loop
133  write(msg,'(3a)') ch10,'-------------------------   Starting propagation   ------------------------',ch10
134  call wrtout(ab_out,msg)
135  if (do_write_log) call wrtout(std_out,msg)
136  
137  do istep = tdks%first_step, tdks%first_step+tdks%ntime-1
138 
139    call cwtime(cpu, wall, gflops, "start")
140  
141    !Perform electronic step
142    !Compute new WF at time t and energy contribution at time t-dt
143    call rttddft_propagate_ele(dtset,istep,mpi_enreg,psps,tdks)
144 
145    !Compute total energy at time t-dt
146    call rttddft_calc_etot(dtset,tdks%energies,tdks%etot,tdks%occ)
147 
148    !Compute new electronic density at t
149    call rttddft_calc_density(dtset,mpi_enreg,psps,tdks)
150 
151    !Compute and output useful electronic values
152    call rttddft_output(dtfil,dtset,istep,mpi_enreg,psps,tdks)
153 
154    !TODO: If Ehrenfest dynamics perform nuclear step
155    !call rttddft_propagate_nuc(dtset,istep,mpi_enreg,psps,tdks)
156 
157    call cwtime(cpu,wall,gflops,"stop")
158    write(msg,'(a,2f8.2,a)') 'Time - cpu, wall (sec):', cpu, wall, ch10
159    call wrtout(ab_out,msg)
160    if (do_write_log) call wrtout(std_out,msg)
161 
162  end do
163 
164  !** 3) Final Output and free memory
165  write(msg,'(7a)') ch10,'---------------------------------------------------------------------------',&
166                  & ch10,'----------------   Finished real-time time dependent DFT   ----------------',&
167                  & ch10,'---------------------------------------------------------------------------',ch10
168  call wrtout(ab_out,msg)
169  if (do_write_log) call wrtout(std_out,msg)
170 
171  call tdks%free(dtset,mpi_enreg,psps)
172 
173 end subroutine rttddft