|
Multidip 1.0
Multi-photon matrix elements
|
Go to the source code of this file.
Data Types | |
| type | string |
Functions/Subroutines | |
| program | multidip_rabitt_delays |
| Calculate RABITT delays from raw multipole elements. | |
| subroutine | print_help |
| Print program usage information. | |
| subroutine | rabitt_main |
| Main program. | |
| logical function | process_command_line (ord1, ord2, ext1, ext2, sph, maxlab, polar, polarir, theta, phi, ek1, ek2, wir, properties, xuv_plus_ir_filenames, xuv_minus_ir_filenames, xuv_plus_ir_outdir, xuv_minus_ir_outdir) |
| Read options from the command line. | |
| subroutine | read_raw_multipoles (order1, order2, sph, ntarg, maxl, nesc, target_states, chains1, chains2, xuv_plus_ir_filenames, xuv_minus_ir_filenames, m_xuv_plus_ir, m_xuv_minus_ir) |
| Read multipoles from disk. | |
| subroutine | setup_chains (order, nchain, chains) |
| Assemble absorption chains. | |
| subroutine | parse_file_name (filename, order, sph, ichain, i, l, m) |
| Extract quantum numbers from raw multipole file name. | |
| subroutine | read_file_data (filename, n, buffer) |
| Read raw multipole from file. | |
| subroutine | extend_order (ord, mm, chains, ext, ek, wir, target_states, echl, prop, maxl, outdir) |
| Calculate multi-photon matrix elements by asymptotic approximation. | |
| subroutine | write_quadratic_dipoles (l, mm, prefix, suffix) |
| Write results to file. | |
| subroutine | calculate_oriented_dipoles (maxl, chains1, chains2, ntarg, nesc, mp, mm, axis, axisir, theta, phi) |
| Molecular- and recoil-frame RABITT delays. | |
| subroutine | read_real_array (filename, array) |
| Read a column from the text file. | |
| subroutine | write_to_file (stem, q) |
| Write the interference term to file. | |
|
private |
Molecular- and recoil-frame RABITT delays.
Calculate the interference term of the photoelectron ionization probability for the following configurations:
Colinear polarization of XUV and IR fields is assumed.
| [in] | maxl | Highest angular momentum |
| [in] | chains1 | List of available absorption paths for first (conjugated) pathway |
| [in] | chains2 | List of available absorption paths for second (non-conugated) pathway |
| [in] | ntarg | Number of targets to consider |
| [in] | nesc | Number of scattering energies |
| [in] | Mp | Raw ionization matrix elements for the XUV+IR pathway (Cartesian basis) |
| [in] | Mm | Raw ionization matrix elements for the XUV-IR pathway (Cartesian basis) |
| [in] | axis | Unit vector in direction of the XUV polarization |
| [in] | axisIR | Unit vector in direction of the IR polarization |
| [in] | theta | List of photoelectron emission polar angles (in degrees) |
| [in] | phi | List of photoelectron emission azimuthal angles (in degrees) |
Definition at line 1005 of file multidip_rabitt_delays.F90.
|
private |
Calculate multi-photon matrix elements by asymptotic approximation.
Use the provided matrix elements of order ord to calculate matrix elements of order ord+ext. The results are returned in the same (reallocated) array. If outdir is given, write out the resulting matrix elements to that directory.
Definition at line 809 of file multidip_rabitt_delays.F90.
| program multidip_rabitt_delays |
Calculate RABITT delays from raw multipole elements.
Having (spherical or Cartesian) raw multipole elements calculated by MULTIDIP for the two different absorption/emission pathways, evaluate the corresponding discrete time delays, unresolved (i.e. averaged) in emission direction as well as in molecular orientation. The raw multipole elements filenames are provided on the command line and they are expected to bear either the original file name assigned by MULTIDIP or additional stem suffixes, for instance "-smooth" as added by the utility program multidip_smooth_rawdips. Example usage:
See the program help for details. Note that the program can be used to calculate one-photon delays, too. For this to work, one needs to produce two sets of the one-photon raw transition dipoles that are offset in energies by 2*omega_IR. These two sets are then used with --xuv-plus-ir and --xuv-minus-ir (with smaller and larger energies, respectively). The results are then for energies in between.
When the parameters --polar, --theta and --phi are specified, molecular frame delays will be calculated as well. The parameter --polar specifies the direction of the polarization axis the polar and azimuthal angle in degrees. The parameters --theta and --phi specify the emission directions in the molecular frame; arbitrary number of angles in degrees can be given. If one of these is omitted, it is assumed to be zero. If the two arrays of angles have unequal lengths, the shorter one is cycled over as long as necessary. When --theta and/or --phi are given, the program will produce filws with the RABITT interference term called "rabitt-signal-mol-?.txt", where "?" is the index of the emission direction. Then there weill be two more sets of results:
Definition at line 57 of file multidip_rabitt_delays.F90.
|
private |
Extract quantum numbers from raw multipole file name.
Expect the multipole to be named "XXX-ab-(u,v,w)YYY.txt" and extract the characters 'a' and 'b' as well as the numbers 'u', 'v' and 'w'.
Definition at line 674 of file multidip_rabitt_delays.F90.
| subroutine multidip_rabitt_delays::print_help |
Print program usage information.
Definition at line 78 of file multidip_rabitt_delays.F90.
|
private |
Read options from the command line.
Obtain command line options, in particular the list of transition multipole files written by MULTIDIP.
| [out] | ord1 | Photon absorption order for conjugated pathway |
| [out] | ord2 | Photon absorption order for non-conjugated pathway |
| [out] | ext1 | Requested asymptotic extension of the absorption order in conjugated pathway |
| [out] | ext2 | Requested asymptotic extension of the absorption order in non-conjugated pathway |
| [out] | sph | Whether the spherical basis is used |
| [out] | maxlab | Highest angular momentum to consider in the laboratory frame |
| [out] | polar | Polarization vector of the XUV field in the molecular frame |
| [out] | polarIR | Polarization vector of the IR field in the molecular frame |
| [in,out] | theta | List of polar emission angles |
| [in,out] | phi | List of azimuthal emission angles |
| [in,out] | Ek1 | List of photoelectron kinetic energies on input (needed for extension to higher orders) |
| [in,out] | Ek2 | List of photoelectron kinetic energies on input (needed for extension to higher orders) |
| [in,out] | properties | Target properties file path |
| [out] | wIR | Energy of the IR quantum (needed for extension to higher orders) |
| [in,out] | XUV_plus_IR_filenames | List of files with transition multipoles for XUV+IR pathway |
| [in,out] | XUV_minus_IR_filenames | List of files with transition multipoles for XUV-IR pathway |
| [in,out] | XUV_plus_IR_outdir | Output directory for extended XUV+IR amplitudes |
| [in,out] | XUV_minus_IR_outdir | Output directory for extended XUV-IR amplitudes |
Definition at line 300 of file multidip_rabitt_delays.F90.
|
private |
Main program.
Read the intermediate data from disk and produces the (orientation averaged) delays.
Definition at line 135 of file multidip_rabitt_delays.F90.
|
private |
Read raw multipole from file.
Read real and imaginary part of the raw multipole for all energies from the given file. The number of read elements is returned via the parameter "n", the values themselves via "buffer". If "buffer" is not allocated or is too short, it will be reallocated by this subroutine.
Definition at line 748 of file multidip_rabitt_delays.F90.
|
private |
Read multipoles from disk.
Read all raw multipole files specified on the command line and store the data in M_XUV_plus_IR and M_XUV_minus_IR.
Definition at line 526 of file multidip_rabitt_delays.F90.
|
private |
Read a column from the text file.
Definition at line 1146 of file multidip_rabitt_delays.F90.
|
private |
Assemble absorption chains.
Construct database of all possible dipole component absorption chains. This is a two-dimensional array, where the second index numbers the individual pathway, while the first index numbers the absorbed photons. Polarizations of the photons are denoted by labels -1, 0, +1. The pathways are organized in lexicographical order. For example, the two-photon pathways look like this:
!> 1: -1, -1 !> 2: -1, 0 !> 3: -1, +1 !> 4: 0, -1 !> 5: 0, 0 !> 6: 0, +1 !> 7: +1, -1 !> 8: +1, 0 !> 9: +1, +1 !>
Definition at line 629 of file multidip_rabitt_delays.F90.
|
private |
Write results to file.
Definition at line 947 of file multidip_rabitt_delays.F90.
|
private |
Write the interference term to file.
Produce one text file per angular combination. Write the real and imaginary part of the interference term, one line per scattering energy.
| [in] | stem | Initial string to be used in the name of the file to write |
| [in] | Q | Interference terms (number of energies times number fo angles) |
Definition at line 1188 of file multidip_rabitt_delays.F90.