File: timepropagation.rst

package info (click to toggle)
gpaw 21.1.0-1
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 14,492 kB
  • sloc: python: 121,997; ansic: 14,138; sh: 1,125; csh: 139; makefile: 43
file content (220 lines) | stat: -rw-r--r-- 9,270 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
.. module:: gpaw.tddft
.. _timepropagation:

======================
Time-propagation TDDFT
======================

Optical photoabsorption spectrum as well as nonlinear effects can be
studied using time propagation TDDFT. This approach
scales better than linear response, but the prefactor is so large that
for small and moderate systems linear response is significantly
faster.


------------
Ground state
------------

To obtain the ground state for TDDFT, one has to just do a standard ground state
with a larger simulation box. A proper distance from any atom to edge of the
simulation box is problem dependent, but a minimum reasonable value is around
6 Ångströms and recommended between 8-10 Ång. In TDDFT, one can use larger
grid spacing than for geometry optimization. For example, if you use h=0.25
for geometry optimization, try h=0.3 for TDDFT. This saves a lot of time.

A good way to start is to use too small box (vacuum=6.0), too large grid
spacing (h=0.35), and too large time step (dt=16.0). Then repeat the simulation
with better parameter values and compare. Probably lowest peaks are already
pretty good, and far beyond the ionization limit, in the continuum, the spectrum
is not going to converge anyway. The first run takes only fraction of
the time of the second run.

For a parallel-over-states TDDFT calculation, you must choose the number
of states so, that these can be distributed equally to processors. For
example, if you have 79 occupied states and you want to use 8 processes
in parallelization over states, add one unoccupied state to get 80 states
in total.


Ground state example::

  # Standard magic
  from ase import Atoms
  from gpaw import GPAW
  
  # Beryllium atom
  atoms = Atoms(symbols='Be',
                positions=[(0, 0, 0)],
                pbc=False)
  
  # Add 6.0 ang vacuum around the atom
  atoms.center(vacuum=6.0)
  
  # Create GPAW calculator
  calc = GPAW(nbands=1, h=0.3)
  # Attach calculator to atoms
  atoms.calc = calc
  
  # Calculate the ground state
  energy = atoms.get_potential_energy()
  
  # Save the ground state
  calc.write('be_gs.gpw', 'all')


Time-propagation TDDFT is also available in :ref:`LCAO mode <lcaotddft>`.

--------------------------------
Optical photoabsorption spectrum
--------------------------------

Optical photoabsorption spectrum can be obtained by applying a weak
delta pulse of dipole electric field, and then letting the system evolve
freely while recording the dipole moment. A time-step around 4.0-8.0
attoseconds is reasonable. The total simulation time should be few tens
of femtoseconds depending on the desired resolution.


Example::

  from gpaw.tddft import *
  
  time_step = 8.0                  # 1 attoseconds = 0.041341 autime
  iterations = 2500                # 2500 x 8 as => 20 fs
  kick_strength = [0.0,0.0,1e-3]   # Kick to z-direction
  
  # Read ground state
  td_calc = TDDFT('be_gs.gpw')
  
  # Kick with a delta pulse to z-direction
  td_calc.absorption_kick(kick_strength=kick_strength)
  
  # Propagate, save the time-dependent dipole moment to 'be_dm.dat',
  # and use 'be_td.gpw' as restart file
  td_calc.propagate(time_step, iterations, 'be_dm.dat', 'be_td.gpw')

  # Calculate photoabsorption spectrum and write it to 'be_spectrum_z.dat'
  photoabsorption_spectrum('be_dm.dat', 'be_spectrum_z.dat')

.. note::

  Make sure to number of iterations is divisible by the dump interval
  such that the last iteration will be stored in the restart file.
  Otherwise append td_calc.write('be_td.gpw', mode='all') to the script.

When propagating after an absorption kick has been applied, it is a good
idea to periodically write the time-evolution state to a restart file.
This ensures that you can resume adding data to the dipole moment file
if you experience artificial oscillations in the spectrum because the total
simulation time was too short.

Example::

  from gpaw.tddft import *
  
  time_step = 8.0                  # 1 attoseconds = 0.041341 autime
  iterations = 2500                # 2500 x 8 as => 20 fs

  # Read restart file with result of previous propagation
  td_calc = TDDFT('be_td.gpw')

  # Propagate more, appending the time-dependent dipole moment to the
  # already existing 'be_dm.dat' and use 'be_td2.gpw' as restart file
  td_calc.propagate(time_step, iterations, 'be_dm.dat', 'be_td2.gpw')

  # Recalculate photoabsorption spectrum and write it to 'be_spectrum_z2.dat'
  photoabsorption_spectrum('be_dm.dat', 'be_spectrum_z2.dat')


Typically in experiments, the spherically averaged spectrum is measured.
To obtain this, one must repeat the time-propagation to each Cartesian
direction and average over the Fourier transformed dipole moments.


---------------------------------------------------
Fourier transformed density and electric near field
---------------------------------------------------

To analyze the resonances seen in the photoabsorption spectrum,
see :ref:`inducedfield`.


--------------------------------
Time propagation
--------------------------------

Since the total CPU time also depends on the number of iterations performed
by the linear solvers in each time-step, smaller time-steps around 2.0-4.0
attoseconds might prove to be faster with the ``ECN`` and ``SICN``
propagators because they have an embedded Euler step in each predictor step:

.. math::

  \tilde{\psi}_n(t+\Delta t) \approx (1 - i \hat{S}^{\;-1}_\mathrm{approx.}(t) \tilde{H}(t) \Delta t)\tilde{\psi}_n(t)

, where `\hat{S}^{\;-1}_\mathrm{approx.}` is an inexpensive operation
which approximates the inverse of the overlap operator `\hat{S}`. See
the :ref:`Developers Guide <overlaps>` for details.


Therefore, as a rule-of-thumb, choose a time-step small enough to minimize the
number of iterations performed by the linear solvers in each time-step, but
large enough to minimize the number of time-steps required to arrive at the
desired total simulation time.


--------------------------------
TDDFT reference manual
--------------------------------

The :class:`~gpaw.tddft.TDDFT` class and keywords:

===================== =============== ============== =====================================
Keyword               Type            Default        Description
===================== =============== ============== =====================================
``ground_state_file`` ``string``                     Name of the ground state file
``td_potential``      ``TDPotential`` ``None``       Time-dependent external potential
``propagator``        ``string``      ``'SICN'``     Time-propagator (``'ECN'``/``'SICN'``/``'SITE'``/``'SIKE'``)
``solver``            ``string``      ``'CSCG'``     Linear equation solver (``'CSCG'``/``'BiCGStab'``)
``tolerance``         ``float``       ``1e-8``       Tolerance for linear solver
===================== =============== ============== =====================================

Keywords for :func:`~gpaw.tddft.TDDFT.absorption_kick`:

================== =============== ================== =====================================
Keyword            Type            Default            Description
================== =============== ================== =====================================
``kick_strength``  ``float[3]``    ``[0,0,1e-3]``     Kick strength
================== =============== ================== =====================================

Keywords for :func:`~gpaw.tddft.TDDFT.propagate`:

====================== =========== =========== ================================================
Keyword                Type        Default     Description
====================== =========== =========== ================================================
``time_step``          ``float``               Time step in attoseconds (``1 autime = 24.188 as``)
``iterations``         ``integer``             Iterations
``dipole_moment_file`` ``string``  ``None``    Name of the dipole moment file
``restart_file``       ``string``  ``None``    Name of the restart file
``dump_interal``       ``integer`` ``500``     How often restart file is written
====================== =========== =========== ================================================

Keywords for :func:`gpaw.tddft.photoabsorption_spectrum`:

====================== ============ ============== ===============================================
Keyword                Type         Default        Description
====================== ============ ============== ===============================================
``dipole_moment_file`` ``string``                  Name of the dipole moment file
``spectrum_file``      ``string``                  Name of the spectrum file
``folding``            ``string``   ``Gauss``      Gaussian folding (or Lorentzian in future)
``width``              ``float``    ``0.2123``     Width of the Gaussian/Lorentzian (in eV)
``e_min``              ``float``    ``0.0``        Lowest energy shown in spectrum (in eV)
``e_max``              ``float``    ``30.0``       Highest energy shown in spectrum (in eV)
``delta_e``            ``float``    ``0.05``       Resolution of energy in spectrum (in eV)
====================== ============ ============== ===============================================

.. autoclass:: gpaw.tddft.TDDFT
   :members:

.. autofunction:: gpaw.tddft.photoabsorption_spectrum