-
Notifications
You must be signed in to change notification settings - Fork 1
QMol_TDDFT_SSO_6BM
Sixth-order optimized Blanes and Moan [Blanes 2002] symplectic split-operator scheme for TDDFT simulations.
Use QMol_TDDFT_SSO_6BM
to propagate the TDDFT dynamics of a spin-polarized or restricted DFT model with a Cartesian-grid discretization domain using a sixth-order optimized Blanes and Moan [Blanes 2002] symplectic split-operator scheme. The time propagation can be performed field free or with an external driving field in the dipole approximation and either length or velocity gauge. The derivation and details of the implementation of the symplectic split-operator can be found in [Mauger 2024]. During the TDDFT propagation, several observables can be computed on-the-fly and stored in output structures in the QMol_TDDFT_SSO_6BM
object. Each of these can be activated independently and may define their own sampling times:
- Save the DFT object into individual files
- Save the dipole, dipole velocity, and dipole acceleration signals
- Save DFT and orbital energies
- Save the external field information
- Save the ionization signal
- Save the Kohn-Sham orbitals (projection) and one-body density
- Save output functions of the one-body density and Kohn-Sham orbitals
- Save restart data file
See the TDDFT documentation page for examples of using most of these features. Note that QMol_TDDFT_SSO_6BM
can only propagate DFT models with local exchange-correlation potential (i.e., LDA and GGA type). Any implicit potential is ignored and triggers warning.
QMol_TDDFT_SSO_6BM
is a handle class overloading QMol_TDDFT_SSO
(and thus QMol_TDDFT_sympSplitOp
and QMol_TDDFT
).
Whether to display progress of the TDDFT calculations as they go on [ true (default) | false ]
Time propagation vector [ vector (default []) ]
- In all TDDFT propagation simulations
T(1)
specifies the starting time andT(end)
the ending time. - In forward time propagation simulations (
timeStep > 0
), time values inT
must be all increasing, while for backward time propagation (timeStep < 0
), they must be all decreasing. Time increments in the vectorT
need not be equally spaced. - When
display
is activated,T
specifies the intermediate times for the time-propagation progress display. -
T
is also the default time sampling for saved results that do not define their own. - Note that the time step used in TDDFT propagations is set independently of
T
.
Time step used for the time propagation [ scalar (default 0.01) ]
- Positive (resp. negative)
timeStep
defines forward (resp. backward) time propagation.
Splitting motif [ 'VTV' (default) | 'TVT' ]
- Whether the splitting starts with applying the potential (
'VTV'
) or kinetic ('TVT'
) part of the DFT Hamiltonian.
Absorber at the edges of the domain [ [] (default) | mask object | CAP object ]
- The boundary absorber aims at eliminating outgoing wave packets from the discretization domain and avoid spurious reflections at the edges of the domain.
- Empty
absorbingBoundary
does not implement any absorbing boundary method and any part of the wave packet reaching the edges of the domain will be reflected (or artificially reappears on the other side of the domain through periodic effects). - Mask absorbers are applied at the end of each propagation step. They are the easiest to use and implement but are first order, irrespective of the order of the time propagation scheme.
- Complex absorbing potentials (CAPs) are integrated within the propagation scheme, as an imaginary potential term that causes exponential decay of the wave packet at the edges. Under the right circumstances, CAPs can preserve the order of the propagation schemes.
External driving field [ [] (default) |
QMol_extField_dipole
object ]
- When performing the TDDFT propagation, any
electricField
,electricFieldDerivative
, orpotentialVector
are repackaged intoexternalField
. - If no electric field or potential vector is defined, either through
externalField
,electricField
, orpotentialVector
, the TDDFT propagation is performed field free.
Gauge in which the external driving field is described [ [] (default) | 'none' | 'length' | 'velocity' ]
-
'none'
ignores any input external field and propagates the field-free TDDFT dynamics - The decision tree for the selection of the default gauge is detailed below.
Electric field of the external driving field [ [] (default) | function handle | griddedInterpolant ]
- For developers: When initializing the class, non-empty
electricField
is moved into theexternalField
property object (if needed creating the object first). During run time,FE
contains the value of the electric field used in the propagation (where relevant).
Derivative of the electric field of the external driving field [ [] (default) | function handle | griddedInterpolant ]
-
griddedInterpolant
electricFieldDerivative
is provided for general support but note that non-function-handle electric field derivatives are computed self-consistently from the electric-field or potential-vector components. - For developers: When initializing the class, non-empty
electricFieldDerivative
is moved into theexternalField
property object (if needed creating the object first). During run time,FDE
contains the value of the derivative of the electric field used in the propagation (where relevant).
Potential vector of the external driving field [ [] (default) | function handle | griddedInterpolant ]
- For developers: When initializing the class, non-empty
potentialVector
is moved into theexternalField
property object (if needed creating the object first). During run time,FA
contains the value of the potential vector used in the propagation (where relevant).
Time step used for computing time derivatives [ nonnegative scalar (default 1e-5) ]
During the time propagation, copies of the DFT object can be saved in separate MATLAB files.
Activate saving the DFT object into separate files [ true | false (default) ]
Name for the files in which the DFT objects are saved [ character array (default 'QMolGrid--TDDFT--DFT') ]
- The DFT objects saved at different times are put in separate MATLAB (.mat) file, with names starting with
saveDFTFileName
to which the iteration index is appended. Each file produced also contains a scalar variablet
with the time information. - One may specify a location where to create the files by indicating the folder path in
saveDFTFileName
-
saveDFTFileName
is irrelevant whensaveDFT == false
Times at which to save the DFT object into a file [ [] (default) | positive scalar | negative integer | vector | 'all' ]
- Empty
saveDFTTime
uses the same values as intime
for when to save the DFT object into MATLAB files. - A positive scalar specifies the sampling time step when to save the DFT object into MATLAB files. For forward time propagation, it is equivalent to
time(1):saveDFTTime:time(end)
, and for backward ones totime(1):-saveDFTTime:time(end)
. - A negative integer specifies the number of propagation time steps between saves. For forward time propagation, it is equivalent to
time(1):abs(saveDFTTime)*timeStep:time(end)
, and for backward ones totime(1):-abs(saveDFTTime)*timeStep:time(end)
. - A vector provides user-defined times at which to save the DFT object to a MATLAB file.
-
'all'
saves the DFT object into a MATLAB file after every time step. Warning: this is slow and may result in a very large number of file and/or large disk usage. - See the note on how intermediate time saving is performed during the TDDFT propagation below.
Whether to calculate and save the dipole signal during the TDDFT propagation [ true | false (default) ]
Index of the Kohn-Sham orbitals for which to compute the orbital-resolved dipole signal [ [] (default) | index vector | cell | 'all' ]
- This is irrelevant if
saveDipole
isfalse
- Empty
saveDipoleOrbitalIndex
disables the orbital-resolved dipole signal calculation - For spin-restricted models, specify the indexes of the orbitals for which to compute the orbital-resolved dipole signal
- For spin-polarized models, combine the up- and down-spin channel orbital indexes in a cell of the form
{ind_up,ind_down}
-
'all'
computes the orbital-resolved dipole for all orbitals
Times at which to compute and save the dipole signal [ [] (default) | positive scalar | negative integer | vector | 'all' ]
- Empty
saveDipoleTime
uses the same values as intime
for when to compute and save the dipole signal. - A positive scalar specifies the sampling time step between successive computations of the dipole signal. For forward time propagation, it is equivalent to
time(1):saveDipoleTime:time(end)
, and for backward ones totime(1):-saveDipoleTime:time(end)
. - A negative integer specifies the number of propagation time steps between dipole-signal computations. For forward time propagation, it is equivalent to
time(1):abs(saveDipoleTime)*timeStep:time(end)
, and for backward ones totime(1):-abs(saveDipoleTime)*timeStep:time(end)
. - A vector provides user-defined times at which to compute and save the dipole signal.
-
'all'
computes and saves the dipole signal after every time step. Warning: this may be slow. - See the note on how intermediate time saving is performed during the TDDFT propagation below.
Whether to calculate and save the dipole velocity signal during the TDDFT propagation [ true | false (default) ]
Index of the Kohn-Sham orbitals for which to compute the orbital-resolved dipole velocity signal [ [] (default) | index vector | cell | 'all' ]
- This is irrelevant if
saveDipoleVelocity
isfalse
- Empty
saveDipoleVelocityOrbitalIndex
disables the orbital-resolved dipole velocity signal calculation - For spin-restricted models, specify the indexes of the orbitals for which to compute the orbital-resolved dipole velocity signal
- For spin-polarized models, combine the up- and down-spin channel orbital indexes in a cell of the form
{ind_up,ind_down}
-
'all'
computes the orbital-resolved dipole velocity for all orbitals
Times at which to compute and save the dipole velocity signal [ 'dipole' (default) | positive scalar | negative integer | vector | 'all' ]
-
'dipole'
saveDipoleVelocityTime
uses the same times as forsaveDipoleTime
- A positive scalar specifies the sampling time step between successive computations of the dipole velocity signal. For forward time propagation, it is equivalent to
time(1):saveDipoleVelocityTime:time(end)
, and for backward ones totime(1):-saveDipoleVelocityTime:time(end)
. - A negative integer specifies the number of propagation time steps between dipole-velocity-signal computations. For forward time propagation, it is equivalent to
time(1):abs(saveDipoleVelocityTime)*timeStep:time(end)
, and for backward ones totime(1):-abs(saveDipoleVelocityTime)*timeStep:time(end)
. - A vector provides user-defined times at which to compute and save the dipole velocity signal.
-
'all'
computes and saves the dipole velocity signal after every time step. Warning: this may be slow. - See the note on how intermediate time saving is performed during the TDDFT propagation below.
Whether to calculate and save the dipole acceleration signal during the TDDFT propagation [ true | false (default) ]
Index of the Kohn-Sham orbitals for which to compute the orbital-resolved dipole acceleration signal [ [] (default) | index vector | cell | 'all' ]
- This is irrelevant if
saveDipoleAcceleration
isfalse
- Empty
saveDipoleAccelerationOrbitalIndex
disables the orbital-resolved dipole acceleration signal calculation - For spin-restricted models, specify the indexes of the orbitals for which to compute the orbital-resolved dipole acceleration signal
- For spin-polarized models, combine the up- and down-spin channel orbital indexes in a cell of the form
{ind_up,ind_down}
-
'all'
computes the orbital-resolved dipole acceleration for all orbitals
Times at which to compute and save the dipole acceleration signal [ 'dipole' (default) | positive scalar | negative integer | vector | 'all' ]
-
'dipole'
saveDipoleAccelerationTime
uses the same times as forsaveDipoleTime
- A positive scalar specifies the sampling time step between successive computations of the dipole acceleration signal. For forward time propagation, it is equivalent to
time(1):saveDipoleAccelerationTime:time(end)
, and for backward ones totime(1):-saveDipoleAccelerationTime:time(end)
. - A negative integer specifies the number of propagation time steps between dipole-acceleration-signal computations. For forward time propagation, it is equivalent to
time(1):abs(saveDipoleAccelerationTime)*timeStep:time(end)
, and for backward ones totime(1):-abs(saveDipoleAccelerationTime)*timeStep:time(end)
. - A vector provides user-defined times at which to compute and save the dipole acceleration signal.
-
'all'
computes and saves the dipole acceleration signal after every time step. Warning: this may be slow. - See the note on how intermediate time saving is performed during the TDDFT propagation below.
Whether to track the DFT energy during the TDDFT propagation [ true | false (default) ]
Times at which to compute and save the DFT energy [ [] (default) | positive scalar | negative integer | vector | 'all' ]
- Empty
saveEnergyDFTTime
uses the same values as intime
for when to compute and save the DFT energy - A positive scalar specifies the sampling time step between successive computations of the DFT energy. For forward time propagation, it is equivalent to
time(1):saveEnergyDFTTime:time(end)
, and for backward ones totime(1):-saveEnergyDFTTime:time(end)
. - A negative integer specifies the number of propagation time steps between DFT-energy computations. For forward time propagation, it is equivalent to
time(1):abs(saveEnergyDFTTime)*timeStep:time(end)
, and for backward ones totime(1):-abs(saveEnergyDFTTime)*timeStep:time(end)
. - A vector provides user-defined times at which to compute and save the DFT energy.
-
'all'
computes and saves the DFT energy after every time step. Warning: this may be slow. - See the note on how intermediate time saving is performed during the TDDFT propagation below.
Whether to track the orbital energies during the TDDFT propagation [ true | false (default) ]
- For a Kohn-Sham orbital
$|\phi \rangle$ , the orbital energy is defined as$\langle \phi |{\hat{\mathcal{H}} }_{{\mathrm{D}\mathrm{F}\mathrm{T}}} |\phi \rangle$ with${\hat{\mathcal{H}} }_{{\mathrm{D}\mathrm{F}\mathrm{T}}}$ the DFT Hamiltonian operator. - In the velocity gauge, the external potential vector is accounted for in the kinetic-operator of
${\hat{\mathcal{H}} }_{{\mathrm{D}\mathrm{F}\mathrm{T}}}$ while, in the length gauge, the external electric field is ignored. - Note that the energy of all the Kohn-Sham orbitals are computed and saved. One may use an installable output function of the Kohn-Sham orbitals to track the energy of a subset of orbitals.
Times at which to compute and save the Kohn-Sham orbital energies [ [] (default) | positive scalar | negative integer | vector | 'all' ]
- Empty
saveEnergyOrbitalTime
uses the same values as intime
for when to compute and save the Kohn-Sham orbital energies - A positive scalar specifies the sampling time step between successive computations of the Kohn-Sham orbital energies. For forward time propagation, it is equivalent to
time(1):saveEnergyOrbitalTime:time(end)
, and for backward ones totime(1):-saveEnergyOrbitalTime:time(end)
. - A negative integer specifies the number of propagation time steps between Kohn-Sham-orbital-energy computations. For forward time propagation, it is equivalent to
time(1):abs(saveEnergyOrbitalTime)*timeStep:time(end)
, and for backward ones totime(1):-abs(saveEnergyOrbitalTime)*timeStep:time(end)
. - A vector provides user-defined times at which to compute and save the Kohn-Sham-orbital energies.
-
'all'
computes and saves the Kohn-Sham orbital energies after every time step. Warning: this may be slow. - See the note on how intermediate time saving is performed during the TDDFT propagation below.
For TDDFT simulation with external driving field, save the external field in the output structures [ true | false (default) ]
- When activated (
saveExternalField = true
), the information about the external driving field at the sampled times is added to each of the output structures. - For practical reasons, the values for the external driving field may be slightly different from that of
externalField
. This features enables keeping the actual external-field values used throughout the propagation.
The ionization signal tracks how much density leaves the simulation domain by being absorbed at the boundaries.
Whether to calculate and save the ionization signal during the TDDFT propagation [ true | false (default) ]
Index of the Kohn-Sham orbitals for which to compute the orbital-resolved ionization signal [ [] (default) | index vector | cell | 'all' ]
- This is irrelevant if
saveIonization
isfalse
- Empty
saveIonizationOrbitalIndex
disables the orbital-resolved ionization signal calculation - For spin-restricted models, specify the indexes of the orbitals for which to compute the orbital-resolved ionization signal
- For spin-polarized models, combine the up- and down-spin channel orbital indexes in a cell of the form
{ind_up,ind_down}
-
'all'
computes the orbital-resolved ionization for all orbitals
Times at which to compute and save the ionization signal [ [] (default) | positive scalar | negative integer | vector | 'all' ]
- Empty
saveIonizationTime
uses the same values as intime
for when to compute and save the ionization signal. - A positive scalar specifies the sampling time step between successive computations of the ionization signal. For forward time propagation, it is equivalent to
time(1):saveIonizationTime:time(end)
, and for backward ones totime(1):-saveIonizationTime:time(end)
. - A negative integer specifies the number of propagation time steps between ionization-signal computations. For forward time propagation, it is equivalent to
time(1):abs(saveIonizationTime)*timeStep:time(end)
, and for backward ones totime(1):-abs(saveIonizationTime)*timeStep:time(end)
. - A vector provides user-defined times at which to compute and save the ionization signal.
-
'all'
computes and saves the ionization signal after every time step. Warning: this may be slow. - See the note on how intermediate time saving is performed during the TDDFT propagation below.
Warning: The memory requirements for saving the orbitals or density throughout the TDDFT propagation may grow very fast and trigger an error (if MATLAB runs out of memory). Instead, if the full orbital/one-body density is required consider saving the DFT object in separate files, or if only the result of a specify calculation on the orbital/one-body density is required consider consider using installable output functions.
Whether to save the one-body density during the TDDFT propagation [ true | false (default) ]
Times at which to save the one-body density [ [] (default) | positive scalar | negative integer | vector | 'all' ]
- Empty
saveDensityTime
uses the same values as intime
for when to save the one-body density. - A positive scalar specifies the sampling time step between successive saving of the one-body density. For forward time propagation, it is equivalent to
time(1):saveDensityTime:time(end)
, and for backward ones totime(1):-saveDensityTime:time(end)
. - A negative integer specifies the number of propagation time steps between saving the one-body density. For forward time propagation, it is equivalent to
time(1):abs(saveDensityTime)*timeStep:time(end)
, and for backward ones totime(1):-abs(saveDensityTime)*timeStep:time(end)
. - A vector provides user-defined times at which to save the one-body density.
-
'all'
saves the one-body density after every time step. Warning: generally discouraged as it likely results in very large memory requirements or causes an out-of-memory error -- see the warning above. - See the note on how intermediate time saving is performed during the TDDFT propagation below.
Whether to save the Kohn-Sham orbitals during the TDDFT propagation [ true | false (default) ]
Index of the Kohn-Sham orbitals to save [ [] (default) | index vector | cell | 'all' ]
- This is irrelevant if
saveOrbital
isfalse
- Empty
saveOrbitalIndex
disables saving of the Kohn-Sham orbitals and is technically equivalent tosaveOrbital = false
- For spin-restricted models, specify the indexes of the orbitals to save
- For spin-polarized models, combine the up- and down-spin channel orbital indexes in a cell of the form
{ind_up,ind_down}
-
'all'
saves all the Kohn-Sham orbitals
Times at which to save the Kohn-Sham orbitals [ [] (default) | positive scalar | negative integer | vector | 'all' ]
- Empty
saveOrbitalTime
uses the same values as intime
for when to save the Kohn-Sham orbitals. - A positive scalar specifies the sampling time step between successive saving of the Kohn-Sham orbitals. For forward time propagation, it is equivalent to
time(1):saveOrbitalTime:time(end)
, and for backward ones totime(1):-saveOrbitalTime:time(end)
. - A negative integer specifies the number of propagation time steps between saving Kohn-Sham orbitals. For forward time propagation, it is equivalent to
time(1):abs(saveOrbitalTime)*timeStep:time(end)
, and for backward ones totime(1):-abs(saveOrbitalTime)*timeStep:time(end)
. - A vector provides user-defined times at which to save the Kohn-Sham orbitals.
-
'all'
saves the Kohn-Sham orbitals after every time step. Warning: generally discouraged as it likely results in very large memory requirements or causes an out-of-memory error -- see the warning above. - See the note on how intermediate time saving is performed during the TDDFT propagation below.
Whether to save the projection of the Kohn-Sham orbitals onto a specific basis during the TDDFT propagation [ true | false (default) ]
Projection basis [ [] (default) | matrix basis |
QMol_disc_basis
]
- Empty
saveOrbitalProjectionBasis
uses the initial Kohn-Sham orbitals as the projection basis - matrix
saveOrbitalProjectionBasis
specifies the projection basis to be used, with the projection vectors specified in each column of the matrix. For spin-polarized models, the same projection basis used for both the projection of the up- and down-spin Kohn-Sham orbitals. Alternatively, one may specify separate basis sets by collecting them into a cell{matrix_basis_up, matrix_basis_down}
. -
QMol_disc_basis
saveOrbitalProjectionBasis
uses the basis set of the associated domain definition. For spin-polarized models, the same projection basis used for both the projection of the up- and down-spin Kohn-Sham orbitals. Alternatively, one may specify separate basis sets by collecting them into a cell{QMol_disc_basis_up, QMol_disc_basis_down}
. - For both matrix and
QMol_disc_basis
, the user-defined basis set is assumed to be a proper orthonormal family defined over the same domain grid as the DFT model being propagated. The TDDFT propagator does not check for this and will produce erroneous results, or trigger an error, otherwise.
Index of the Kohn-Sham orbitals for which to perform the projection [ [] (default) | index vector | cell | 'all' ]
- This is irrelevant if
saveOrbitalProjectionBasis
isfalse
- Empty
saveOrbitalProjectionIndex
disables saving of the Kohn-Sham orbitals and is technically equivalent tosaveOrbitalProjectionBasis = false
- For spin-restricted models, specify the indexes of the orbitals for which to perform the projection
- For spin-polarized models, combine the up- and down-spin channel orbital indexes in a cell of the form
{ind_up,ind_down}
-
'all'
saves the projection for all the Kohn-Sham orbitals
Times at which to save the projection of the Kohn-Sham orbitals [ [] (default) | positive scalar | negative integer | vector | 'all' ]
- Empty
saveOrbitalProjectionTime
uses the same values as intime
for when to save the projection of the Kohn-Sham orbitals. - A positive scalar specifies the sampling time step between successive savings of the Kohn-Sham orbitals' projection. For forward time propagation, it is equivalent to
time(1):saveOrbitalProjectionTime:time(end)
, and for backward ones totime(1):-saveOrbitalProjectionTime:time(end)
. - A negative integer specifies the number of propagation time steps between saving Kohn-Sham orbitals' projection. For forward time propagation, it is equivalent to
time(1):abs(saveOrbitalProjectionTime)*timeStep:time(end)
, and for backward ones totime(1):-abs(saveOrbitalProjectionTime)*timeStep:time(end)
. - A vector provides user-defined times at which to save the Kohn-Sham orbitals' projection.
-
'all'
saves the Kohn-Sham orbitals' projection after every time step. Warning: this may be slow. - See the note on how intermediate time saving is performed during the TDDFT propagation below.
These enable on-the-fly computation and saving of user-defined observables without having to save the DFT, Kohn-Sham orbitals, or one-body density.
Installable output function of the one-body density [ [] (default) | function handle ]
- Leave empty to disable the feature
- Provide a function handle to enable the feature. The signature for the function is
fun(rho,t)
whererho
is a one-body density object andt
is the time (scalar). The function handle may return an array of arbitrary size and shape, but must return at least a scalar and the shape of the output must remain constant throughout the TDDFT propagation.
Times at which to evaluate and save the installable output function of the one-body density [ [] (default) | positive scalar | negative integer | vector | 'all' ]
- Empty
saveOutputFunctionDensityTime
uses the same values as intime
for when to evaluate and save the output function. - A positive scalar specifies the sampling time step between successive evaluation and saving of the output function. For forward time propagation, it is equivalent to
time(1):saveOutputFunctionDensityTime:time(end)
, and for backward ones totime(1):-saveOutputFunctionDensityTime:time(end)
. - A negative integer specifies the number of propagation time steps between evaluations and saving of the output function. For forward time propagation, it is equivalent to
time(1):abs(saveOutputFunctionDensityTime)*timeStep:time(end)
, and for backward ones totime(1):-abs(saveOutputFunctionDensityTime)*timeStep:time(end)
. - A vector provides user-defined times at which to evaluate and save the output function.
-
'all'
evaluates and saves the output function after every time step. Warning: this may be slow. - See the note on how intermediate time saving is performed during the TDDFT propagation below.
Installable output function of the Kohn-Sham orbitals [ [] (default) | function handle ]
- Leave empty to disable the feature
- Provide a function handle to enable the feature. The signature for the function is
fun(KSO,t)
whereKSO
is a Kohn-Sham orbital object andt
is the time (scalar). The functional handle may return an array of arbitrary size and shape, but must return at least a scalar and the shape of the output must remain constant throughout the TDDFT propagation. - Warning: The DFT model Kohn-Sham orbitals are passed by reference to the output function. Thus modifying the Kohn-Sham orbitals in the output function will also modify them for the DFT model (and thus TDDFT propagation) and will likely result in erroneous results or produce an error.
Times at which to evaluate and save the installable output function of the Kohn-Sham orbitals [ [] (default) | positive scalar | negative integer | vector | 'all' ]
- Empty
saveOutputFunctionOrbitalTime
uses the same values as intime
for when to evaluate and save the output function. - A positive scalar specifies the sampling time step between successive evaluation and saving of the output function. For forward time propagation, it is equivalent to
time(1):saveOutputFunctionOrbitalTime:time(end)
, and for backward ones totime(1):-saveOutputFunctionOrbitalTime:time(end)
. - A negative integer specifies the number of propagation time steps between evaluations and saving of the output function. For forward time propagation, it is equivalent to
time(1):abs(saveOutputFunctionOrbitalTime)*timeStep:time(end)
, and for backward ones totime(1):-abs(saveOutputFunctionOrbitalTime)*timeStep:time(end)
. - A vector provides user-defined times at which to evaluate and save the output function.
-
'all'
evaluates and saves the output function after every time step. Warning: this may be slow. - See the note on how intermediate time saving is performed during the TDDFT propagation below.
During the time propagation, a restart file can be generated to enable resuming the simulation if it is cut short. Stopping the TDDFT simulation while it writes the restart info may lead to a corrupted file from which restart will not be possible. The restart file contains a copy of the TDDFT-propagator and DFT-model objects, respectively called TDDFT
and DFT
.
Activate the generation of a restart file [ true | false (default) ]
Name for the restart file [ character array (default 'QMolGrid--TDDFT--Restart.mat') ]
-
saveRestartFileName
must be a valid MATLAB file name ('.mat' extension). -
saveRestartFileName
is irrelevant whensaveRestart == false
Times at which to generate or update the restart file [ [] (default) | positive scalar | negative integer | vector | 'all' ]
- Empty
saveRestartTime
uses the same values as intime
for when to generate and update the restart file. - A positive scalar specifies the sampling time step when to generate and update the restart file. For forward time propagation, it is equivalent to
time(1):saveRestartTime:time(end)
, and for backward ones totime(1):-saveRestartTime:time(end)
. - A negative integer specifies the number of propagation time steps between the generation and update of the restart file. For forward time propagation, it is equivalent to
time(1):abs(saveRestartTime)*timeStep:time(end)
, and for backward ones totime(1):-abs(saveRestartTime)*timeStep:time(end)
. - A vector provides user-defined times at which to generate and update the restart file.
-
'all'
generates and updates the restart file after each time step. Warning: this is very slow and is strongly discouraged. - Note: In all cases no restart file is generated for the initial time and final propagation time-step.
- See the note on how intermediate time saving is performed during the TDDFT propagation below.
During a TDDFT propagation, the results of on-the-fly calculations are stored in structures in the QMol_TDDFT
object. Note that QMol_TDDFT
does not interpolate its time propagation to fit user-supplied sample times, Instead, the results are saved at the closest propagation time steps, excluding duplicate times. This may result in sampled times that are different, or have a different (smaller) number of elements, from the ones specified above. Notably, using a sampling time that is not a multiple of the propagation time step may result in uneven saved time sampling. The actual times at with output are saved is included in each of the output structure.
When saveExternalField
is true
, each output-result structure includes the following fields:
- In the length gauge,
outStructure.electricField
contains the specific values for the electric field used in the TDSE propagation at the saved times. - In the velocity gauge,
outStructure.potentialVector
contains the specific values for the potential vector used in the TDSE propagation at the saved times. - In all cases, if the potential vector, electric field, or derivative of the electric field are calculated by the propagator during the simulations (see the decision tree below), those are included in
outStructure.potentialVector
,outStructure.electricField
, andoutStructure.electricFieldDerivative
, respectively.
In restart mode, each of the output structure contains a handful of additional fields to the ones listed below that are relevant for run-time calculations. These fields are saved in the MATALB restart file but removed upon completion of the TDDFT propagation.
Result of dipole-signal calculations during the TDDFT propagation [ [] | structure ]
- TDDFT propagation that does not compute the dipole signal (
saveDipole = false
) leavesoutDipole
empty. Otherwise: -
outDipole.time
defines the precise times at which the dipole signal is computed. -
outDipole.total
contains the corresponding total dipole signal and, for spin-polarized DFT models,outDipole.totalUp
andoutDipole.totalDown
contain the total dipole signal in the up- and down-spin channels, respectively. - Orbital-resolved dipole signals are stored in
outDipole.orbital_x
,outDipole.orbitalUp_x
, andoutDipole.orbitalDown_x
. The indices of the orbitals associated with each of these are included inoutDipole.indexOrbital
,outDipole.indexOrbitalUp
, andoutDipole.indexOrbitalDown
.
Result of dipole-velocity-signal calculations during the TDDFT propagation [ [] | structure ]
- TDDFT propagation that does not compute the dipole signal (
saveDipoleVelocity = false
) leavesoutDipoleVelocity
empty. Otherwise: -
outDipoleVelocity.time
defines the precise times at which the dipole velocity signal is computed. -
outDipoleVelocity.total
contains the corresponding total dipole velocity signal and, for spin-polarized DFT models,outDipoleVelocity.totalUp
andoutDipoleVelocity.totalDown
contain the total dipole velocity signal in the up- and down-spin channels, respectively. - Orbital-resolved dipole velocity signals are stored in
outDipoleVelocity.orbital_x
,outDipoleVelocity.orbitalUp_x
, andoutDipoleVelocity.orbitalDown_x
. The indices of the orbitals associated with each of these are included inoutDipoleVelocity.indexOrbital
,outDipoleVelocity.indexOrbitalUp
, andoutDipoleVelocity.indexOrbitalDown
.
Result of dipole-acceleration-signal calculations during the TDDFT propagation [ [] | structure ]
- TDDFT propagation that does not compute the dipole signal (
saveDipoleAcceleration = false
) leavesoutDipoleAcceleration
empty. Otherwise: -
outDipoleAcceleration.time
defines the precise times at which the dipole acceleration signal is computed. -
outDipoleAcceleration.total
contains the corresponding total dipole acceleration signal and, for spin-polarized DFT models,outDipoleAcceleration.totalUp
andoutDipoleAcceleration.totalDown
contain the total dipole acceleration signal in the up- and down-spin channels, respectively. - Orbital-resolved dipole velocity signals are stored in
outDipoleAcceleration.orbital_x
,outDipoleAcceleration.orbitalUp_x
, andoutDipoleAcceleration.orbitalDown_x
. The indices of the orbitals associated with each of these are included inoutDipoleAcceleration.indexOrbital
,outDipoleAcceleration.indexOrbitalUp
, andoutDipoleAcceleration.indexOrbitalDown
.
Result of the DFT-energy calculations during the TDDFT propagation [ [] | structure ]
- TDDFT propagation that does track the DFT energy (
saveEnergyDFT = false
) leavesoutEnergyDFT
empty. Otherwise: -
outEnergyDFT.time
defines the precise times at which the DFT energies are computed. -
outEnergyDFT.total
contains the total DFT energy, which should be constant (within the precision of the propagation scheme). -
outEnergyDFT.kinetic
contains the kinetic-energy component, including the driving-field potential vector in the velocity gauge. For spin-polarized models the two rows contain the up- and down-spin kinetic energy components, respectively. -
outEnergyDFT.external
contains the external-energy component. For spin-polarized models the two rows contain the up- and down-spin external energy components, respectively -
outEnergyDFT.Hartree
contains the Hartree-energy component. -
outEnergyDFT.exchangeCorrelation
contains the exchange-correlation-energy component. -
outEnergyDFT.externalField
contains the energy contribution from the driving electric field in the length gauge. For spin-polarized models the two rows contain the up- and down-spin energy components, respectively. -
outEnergyDFT.autonomization
contains the energy brought in and out of the system by the external driving field (if any).
Result of the Kohn-Sham-orbital-energy calculations during the TDDFT propagation [ [] | structure ]
- TDDFT propagation that does track the DFT energy (
saveEnergyOrbital = false
) leavesoutEnergyDFT
empty. Otherwise: - For a Kohn-Sham orbital
$|\phi \rangle$ , the orbital energy is defined as$\langle \phi |{\hat{\mathcal{H}} }_{{\mathrm{D}\mathrm{F}\mathrm{T}}} |\phi \rangle$ with${\hat{\mathcal{H}} }_{{\mathrm{D}\mathrm{F}\mathrm{T}}}$ the DFT Hamiltonian operator. In the velocity gauge, the external potential vector is accounted for in the kinetic-operator of${\hat{\mathcal{H}} }_{{\mathrm{D}\mathrm{F}\mathrm{T}}}$ while, in the length gauge, the external electric field is ignored. -
outEnergyOrbital.time
defines the precise times at which the Kohn-Sham-orbital energies are computed - For spin-restricted models, each row of
outEnergyOrbital.orbital
contains the energies of the corresponding Kohn-Sham orbital - For spin-polarized models, each row of
outEnergyOrbital.orbitalUp
andoutEnergyOrbital.orbitalDown
contains the energies of the corresponding Kohn-Sham orbital for the up- and down-spin channels, respectively.
Result of the ionization calculations during the TDDFT propagation [ [] | structure ]
- TDDFT propagation that does track ionization (
saveIonization = false
) leavesoutIonization
empty. Otherwise: -
outIonization.time
defines the precise times at which the ionization is computed. -
outIonization.total
contains the total ionization. - For spin-polarized models,
outIonization.totalUp
andoutIonization.totalDown
contain the total ionization for the up- and down-spin channels, respectively. - Orbital-resolved ionization signals are stored in
outIonization.orbital
,outIonization.orbitalUp
, andoutIonization.orbitalDown
. The indices of the orbitals associated with each of these are included inoutIonization.indexOrbital
,outIonization.indexOrbitalUp
, andoutIonization.indexOrbitalDown
.
One-body densities during the TDDFT propagation [ [] | structure ]
- TDDFT propagation that does not save the one-body density (
saveDensity = false
) leavesoutDensity
empty. Otherwise: - The specific shape of the saved one-body densities is implementation dependent -- see the corresponding documentation for details. At minimum, each implementation defines:
-
outDensity.time
defines the precise times at which the one-body density is saved. -
outDensity.total
contains the corresponding one-body densities and, for spin-polarized DFT models,outDensity.totalUp
andoutDensity.totalDown
contain the one-body densities in the up- and down-spin channels respectively.
Kohn-Sham orbitals during the TDDFT propagation [ [] | structure ]
- TDDFT propagation that does not save the Kohn-Sham orbitals (
saveOrbital = false
) leavesoutOrbital
empty. Otherwise: - The specific shape of the saved Kohn-Sham orbitals is implementation dependent -- see the corresponding documentation for details. At minimum, each implementation defines:
-
outOrbital.time
defines the precise times at which the Kohn-Sham orbitals are saved. - For spin-restricted models,
outOrbital.orbital
contains the saved Kohn-Sham orbitals andoutOrbital.indexOrbital
the corresponding orbital indices. - For spin-polarized models,
outOrbital.orbitalUp
andoutOrbital.orbitalDown
contains the saved Kohn-Sham orbitals for the up- and down-spin channels, respectively, andoutOrbital.indexOrbitalUp
andoutOrbital.indexOrbitalDown
the corresponding orbital indices.
Projection of the Kohn-Sham orbitals during the TDDFT propagation [ [] | structure ]
- TDDFT propagation that does not save the projection of the Kohn-Sham orbitals (
saveOrbitalProjection = false
) leavesoutOrbitalProjection
empty. Otherwise: - The specific shape of the saved Kohn-Sham orbitals' projection is implementation dependent -- see the corresponding documentation for details. At minimum, each implementation defines:
-
outOrbitalProjection.time
defines the precise times at which the Kohn-Sham orbitals' projection is saved. - For spin-restricted models,
outOrbitalProjection.orbital
contains the saved Kohn-Sham orbitals' projection andoutOrbitalProjection.indexOrbital
the corresponding orbital indices. - For spin-polarized models,
outOrbitalProjection.orbitalUp
andoutOrbitalProjection.orbitalDown
contains the saved Kohn-Sham orbitals' projection for the up- and down-spin channels, respectively, andoutOrbitalProjection.indexOrbitalUp
andoutOrbitalProjection.indexOrbitalDown
the corresponding orbital indices. - The structure also retains a copy of the projection basis as
QMol_disc_basis
object(s).
Result of the installable output function of the one-body density [ [] | structure ]
- TDDFT propagation that does not define an installable output function of the density (
saveOutputFunctionDensity = []
) leavesoutOutputFunctionDensity
empty. Otherwise: -
outOutputFunctionDensity.time
defines the precise times at which the result of the output function is saved. -
outOutputFunctionDensity.result
contains the results of the output functions. If the installable function returns a scalar or a column vector,outOutputFunctionDensity.result
is aN-by-numel(outOutputFunctionDensity.time)
matrix withN
the number of elements in the output. Otherwise,outOutputFunctionDensity.result
is anM-by-numel(outOutputFunctionDensity.time)
array withM
the size of the output. -
outOutputFunctionDensity.shape
contains the shape of the output function (N
orM
, as defined in the previous bullet point).
Result of the installable output function of the Kohn-Sham orbitals [ [] | structure ]
- TDDFT propagation that does not define an installable output function of the orbitals (
saveOutputFunctionOrbital = []
) leavesoutOutputFunctionOrbital
empty. Otherwise: -
outOutputFunctionOrbital.time
defines the precise times at which the result of the output function is saved. -
outOutputFunctionOrbital.result
contains the results of the output functions. If the installable function returns a scalar or a column vector,outOutputFunctionOrbital.result
is aN-by-numel(outOutputFunctionOrbital.time)
matrix withN
the number of elements in the output. Otherwise,outOutputFunctionOrbital.result
is anM-by-numel(outOutputFunctionOrbital.time)
array withM
the size of the output. -
outOutputFunctionOrbital.shape
contains the shape of the output function (N
orM
, as defined in the previous bullet point).
Create a sixth-order optimized Blanes and Moan [Blanes 2002] symplectic TDDFT-propagator object with empty class properties.
obj = QMol_TDDFT_SSO_6BM;
Create a TDDFT-propagator object with the name
properties set to the specified value
. Several name-value
pairs can be specified consecutively. Suitable name
is any of the TDDFT class properties and is case insensitive.
obj = QMol_TDDFT_SSO_6BM(name1,value1);
obj = QMol_TDDFT_SSO_6BM(name1,value1,name2,value2,___);
Update the name
properties of a TDDFT-model object to the specified value
. Several name-value
pairs can be specified consecutively. Suitable name
is any of the TDDFT class properties and is case insensitive. In restart mode, output result structures can also be edited with set
. QMol_TDDFT
does not check the integrity of its input/output component during a restart and this feature should only be considered by advanced users.
obj.set(name1,value1);
obj.set(name1,value1,name2,value2,___);
This is the common name-value pair assignment method used throughout the QMol-grid package.
Clear all class properties
obj.clear;
Clear a specific set of the class properties. Suitable name
is any of the TDDFT class properties and is case insensitive.
obj.clear(name1,name2,___);
This is the common clear
method available to all classes throughout the QMol-grid package. The clear
method can be used to delete specific properties before saving an instance of the TDDFT propagator class and propagation results.
Initialize a sixth-order optimized Blanes and Moan [Blanes 2002] symplectic TDDFT-propagator object without allocating the output-result structures
obj.initialize(DFT);
-
DFT
is the DFT-model handle object, i.e.,QMol_DFT_spinPol
orQMol_DFT_spinRes
, that describes the DFT system to propagate. - For developers:
initialize
has also a specific interface when called from thepropagate
method that initializes the absorbing boundaries (if any) and determines whether the output-result structures should be initialized too. Overloading classes should avoid redefininginitialize
; If the overloading class needs to perform some initialization (at the beginning of a TDDFT propagation or upon restart), overload theinitializeChildren
method instead.
Get an estimate of the memory help by a QMol_TDDFT_SSO_6BM
object with either
mem = obj.getMemoryProfile;
mem = obj.getMemoryProfile(false);
- The object must be
initialize
d for the memory footprint evaluation. If not already,getMemoryProfile
initializes theDFT
model object. - The estimate includes the (1)
DFT
model, (2) TDDFT propagator, and (3) output results. Note that all these components may not be used in actual simulations and the memory estimate tries to be conservative. On the other hand, it only includes the discretization of member components on the domain grid and ignores other (expectedly small) properties. - The output
mem
is the estimated size in bytes.
Additionally display the detail of the memory footprint with
mem = obj.getMemoryProfile(true);
Display the run-time documentation for the specific configuration of a QMol_TDDFT_SSO_6BM
object, which must have been initialize
d beforehand
obj.showDocumentation;
Propagate the TDDFT dynamics starting from a DFT model object (from scratch)
obj.propagate(DFT);
-
DFT
is the DFT-model handle object, i.e.,QMol_DFT_spinPol
orQMol_DFT_spinRes
, that describes the DFT system to propagate. - Neither the TDDFT-propagation
obj
nor the DFT-modelDFT
objects need beinitialize
d. In all cases both are (re)initialized at the beginning of the TDDFT propagation to ensure proper and consistent linkage and setup between them.
Restart a TDDFT propagation
obj.propagate('restart');
- This can be performed after loading the corresponding restart MATLAB file into the workspace.
If no specific gauge is specified with externalFieldGauge
, the gauge is selected following the decision tree:
- Length gauge if
externalField.electricField
orelectricField
is a function handle, otherwise: - Velocity gauge if
externalField.potentialVector
orpotentialVector
is a function handle, otherwise: - Length gauge if
externalField.electricField
orelectricField
is agriddedInterpolant
, otherwise: - Velocity gauge if
externalField.potentialVector
orpotentialVector
is agriddedInterpolant
, otherwise: - Field free (ignore any input field)
Depending on the types of input, not all provided external-field may be used in the simulations. For the potential vector and electric field -- recall that electricField
and potentialVector
are repackaged into externalField
at the beginning of the simulation:
- If both are provided as function handles, then
externalField.electricField
andexternalField.potentialVector
are called whenever the values for the potential vector or electric field are required. Note that the TDDFT propagator does not check whether the provided potential vector and electric field are consistent with each other ($E(t)=-\partial_t A(t)$ ) and providing non-matching components will likely result in erroneous results. Otherwise, - If only one of them is provided as a function handle, then both the potential vector and electric field are computed from that same function handle (see below for details on how this is done). Otherwise,
- In the length gauge, if
externalField.electricField
is agriddedInterpolant
then it is used to compute both the electric field and potential vector, otherwiseexternalField.potentialVector
is used to compute them both. - In the velocity gauge, if
externalField.potentialVector
is agriddedInterpolant
then it is used to compute both the electric field and potential vector, otherwiseexternalField.electricField
is used to compute them both. - If none of the
externalField.electricField
andexternalField.potentialVector
are defined, the TDDFT propagation is performed field free irrespective of the choice ofexternalFieldGauge
and even if anexternalField.electricFieldDerivative
is defined.
The choice of ignoring some (interpolated) input fields is made to ensure the self-consistency of the TDDFT dynamics, and associated observables, during the numerical propagation. To force using an interpolated field input, one can encapsulate it into a function handle @(t) interpolatedField(t)
.
For the electric-field derivative:
- If
externalField.electricFieldDerivative
is defined as a function handle, thenexternalField.electricFieldDerivative
is called when the derivative of the electric field is required. - Otherwise, the electric-field derivative is computed numerically from
externalField.electricField
orexternalField.potentialVector
, depending on which field component(s) are being used in the TDDFT propagation.
Numerical computation of missing or ignored field components:
- Where required, the potential vector is computed from the electric field using a Simpson's 1/3 quadrature rule.
- Where required, the electric field is computed from the potential vector using a second-order centered finite difference with
diffDT
time step. - Where required, the electric field derivative is computed from the electric field or the potential vector using a second-order centered finite difference with
diffDT
time step.
See the main documentation page for examples of TDDFT simulations with various output computations.
For consistency with the rest of the QMol-grid package, QMol_TDDFT_SSO_6BM
defines an associated test suite. Run the test suite for the class in normal or summary mode respectively with
QMol_test.test('TDDFT_SSO_6BM');
QMol_test.test('-summary','TDDFT_SSO_6BM');
[Blanes 2002] S. Blanes and P.C. Moan, "Practical symplectic partitioned Runge-Kutta and Runge-Kutta-Nystrom methods," Journal of Computational and Applied Mathematics 142, 313 (2002).
[Mauger 2024] F. Mauger, C. Chandre, M.B. Gaarde, K. Lopata, and K.J. Schafer, "Hamiltonian formulation and symplectic split-operator schemes for time-dependent density-functional-theory equations of electron dynamics in molecules," Communications in Nonlinear Science and Numerical Simulation 129, 107685 (2024).
-
QMol_TDDFT_SSO_6BM
was introduced in version 01.10
This wiki is a copy of the documentation provided with the QMol-grid package (accessible in MATLAB documentation, via the "Supplemental Software" section).
Copyright © 2024, Francois Mauger, all right reserved.
Density-functional theory (DFT)
QMol_DFT_density
QMol_DFT_eigs
QMol_DFT_eig_basis
QMol_DFT_orbital
QMol_DFT_orbital_basis
QMol_DFT_profiler
QMol_DFT_SCF_Anderson
QMol_DFT_spinPol
QMol_DFT_spinRes
QMol_DFT_Vc_LDA_soft
QMol_DFT_Vext
QMol_DFT_Vh_conv
QMol_DFT_Vh_fft
QMol_DFT_Vks
QMol_DFT_Vks_basis
QMol_DFT_Vks_grad
QMol_DFT_Vx_LDA_exp
QMol_DFT_Vx_LDA_soft
QMol_DFT_Vx_XX_conv
QMol_DFT_Vx_XX_fft
Tutorials
- Tutorial 1: Schrödinger-equation ground state
- Tutorial 2: Schrödinger-equation input and output
- Tutorial 3: Time-dependent Schrödinger equation
- Tutorial 4: Time-dependent Schrödinger-equation input and output
- Tutorial 5: Density-functional theory ground state
- Tutorial 6: Time-dependent density-functional theory
- Tutorial 7: Time-dependent density-functional theory input and output