qdyn.io module

Routines for reading and writing files compatible with QDYN

Summary

Functions:

datablock

Iterator over the lines inside the file with the given filename inside the given block only.

fix_fortran_exponent

In 3-digit exponents, Fortran drops the ‘E’.

iterate_psi_amplitudes

Iterate over blocks of wave functions stored in filename

open_file

Wrapper around io.open(), allowing file to be a file handle.

print_matrix

Print a numpy complex matrix to screen, or to a file if outfile is given.

read_ascii_dump

Read a file in the format written by the QDYN dump_ascii_* routines and return its data as a nested OrderedDict.

read_cmplx_array

Read a complex array from a file.

read_cmplx_matrix

Read in complex matrix from file (as written by the Fortran QDYN print_cmplx_matrix routine).

read_complex

Convert a string to a complex number

read_indexed_matrix

Read in a matrix from the file with the given filename

read_psi_amplitudes

Read the wave function of size n from file.

read_real

Convert a string to a real number

tempinput

Context manager providing a temporary filename for a file containing the given data.

write_cmplx_array

Write a complex array to file.

write_indexed_matrix

Write the given matrix to file in indexed format (1-based indexing)

write_psi_amplitudes

Write the wave function to file in the same format as the write_psi_amplitudes Fortran routine

writetotxt

Inverse function to numpy.genfromtxt and similar to numpy.savetxt, but allowing to write multiple numpy arrays as columns to a text file.

Reference

qdyn.io.open_file(file, mode='r', **kwargs)[source]

Wrapper around io.open(), allowing file to be a file handle.

Any kwargs are passed to io.open(). The file parameter may also be given as the string ‘-‘, which indicates sys.stdin for mode='r' and sys.stdout for mode='w'. If file is an open file handle, its mode and encoding are checked against the arguments; the file handle will not be closed one exit.

In Python 2, if no kwargs beyond mode are given and file is a string giving a pathname or a file descriptor, the built-in open() function is used instead of io.open(). As a consequence, regular (unencoded) strings can be read or written directly from or to the resulting file handle, without having to open it in binary mode. This behavior provides maximal compatibility between Python 2 and 3.

This function must be used as a context manager:

>>> with open_file('-', 'w') as out_fh:
...      written_bytes = out_fh.write("Hello World")
Hello World
Raises

IOError – If the file cannot be opened (see io.open()), if file is open file handle with the incorrect mode or encoding, or if file is ‘-‘ and mode is neither ‘r’ nor ‘w’.

qdyn.io.tempinput(data, binary=False)[source]

Context manager providing a temporary filename for a file containing the given data. If binary is True, the data will be written as-is, and must be suitable for writing in binary mode. Otherwise, if encoding the given data to utf-8 is at all possible, the temporary file will be a text file with utf-8 encoding. The file is deleted on leaving the context.

>>> test_str = '''
... In the world of the very small, where particle and wave
... aspects of reality are equally significant, things do not
... behave in any way that we can understand from our experience
... of the everyday world...all pictures are false, and there is
... no physical analogy we can make to understand what goes on
... inside atoms. Atoms behave like atoms, nothing else.'''
>>> with tempinput(test_str) as filename:
...     with open_file(filename) as in_fh:
...         for line in in_fh:
...             print(line.strip())

In the world of the very small, where particle and wave
aspects of reality are equally significant, things do not
behave in any way that we can understand from our experience
of the everyday world...all pictures are false, and there is
no physical analogy we can make to understand what goes on
inside atoms. Atoms behave like atoms, nothing else.
qdyn.io.write_indexed_matrix(matrix, filename, comment=None, line_formatter=None, header=None, hermitian=False, limit=0.0)[source]

Write the given matrix to file in indexed format (1-based indexing)

Parameters
  • matrix (numpy.matrix, numpy.ndarray, qutip.Qobj or scipy.sparse.spmatrix) – Matrix to write to file. If an ndarray is given, it has to be 2-dimenstional

  • filename (str) – Name of file to write to

  • comment (str or list(str)) – Comment line, or array of comment lines to write to the top of the file. Each line that does not start with ‘#’ will have “# ” prepended.

  • line_formatter (callable) –

    Function that takes three arguments i, j, v (row index, column index, and complex value matrix[i,j]) and returns a line to be written to file. If the function returns None for any input data, no line will be written to file. If not given, defaults to

    lambda i, j, v: “%8d%8d%25.16E” % (i, j, v.real)

    if matrix is real and

    lambda i, j, v: “%8d%8d%25.16E%25.16E” % (i, j, v.real, v.imag)

    if matrix is complex.

  • header (str) – Header line to be written before any data. Must start with either ‘#’ or a space, in which case the leading space will be replaced with ‘#’. Defaults to a header line suitable for the default line_formatter

  • hermitian (bool) – If True, write only entries from the upper triangle

  • limit (float) – Only values with an absolute value greater than limit are written

qdyn.io.read_indexed_matrix(filename, format='coo', shape=None, expand_hermitian=False, val_real=False)[source]

Read in a matrix from the file with the given filename

The file must contain a description in indexed format, like this:

# row col re(val) im(val)

0 1 1.0 0.0 1 0 0.0 1.0

The fourth column is optional, if not present, the result will be real.

Return a matrix in any of the numpy/scipy sparse (or non-sparse) formats. See the documentation of scipy.sparse for information about the different sparse formats

Parameters
  • filename (str) – Name of file from which to read the matrix

  • format (str) – Result type: * ‘coo’ (default): scipy.sparse.coo.coo_matrix * ‘array’: numpy.ndarray * ‘dense’: numpy.matrixlib.defmatrix.matrix * ‘bsr’: scipy.sparse.bsr.bsr_matrix * ‘csc’: scipy.sparse.csc.csc_matrix * ‘csr’: scipy.sparse.csr.csr_matrix * ‘dia’: scipy.sparse.dia.dia_matrix * ‘dok’: scipy.sparse.dok.dok_matrix * ‘lil’: scipy.sparse.lil.lil_matrix

  • shape (int or (int,int)) – If given, shape of the resulting matrix. If not given, will be determined from largest occurring index in the data from the input file

  • expand_hermitian (bool) – If True, the input file must contain data only for the upper or lower triangle. The oterh triangle will be set with the complex conjugate values.

  • val_real (bool) – If True, only read 3 columns from the input file (i, j, value), even if more columns are present in the file, and return a real matrix.

qdyn.io.print_matrix(M, matrix_name=None, limit=1e-14, fmt='%9.2E', compress=False, zero_as_blank=False, out=None)[source]

Print a numpy complex matrix to screen, or to a file if outfile is given. Values below the given limit are printed as zero

Parameters
  • M (numpy.ndarray, scipy.sparse.spmatrix) – Matrix to print. In addition to a standard dense matrix, may also be any scipy sparse matrix in a format where M[i,j] is defined.

  • matrix_name (str) – Name of matrix

  • limit (float) – Any number (real or imaginary part) whose absolute value is smaller than this limit will be printed as 0.0.

  • fmt (str or callable) – Format of each entry (both for real and imaginary part). If str, must be an “old-style” format string the formats a single floating value. If a callable, the callable must return a string of fixed length when passed a number. The string returned by fmt(0) will be used for values that are exactly zero, whereas the string returned by fmt(0.0) will be used for values that are below limit.

  • compress (bool) – If True, remove spaces to compress the output to be narrower. Real and imaginary parts will no longer be aligned.

  • zero_as_blank (bool) – If True, represent entries that are exactly zero as blank strings

  • out (file) – open filehandle. If None, print to stdout

Examples

>>> import numpy as np
>>> M = np.array([[1.0, 2.0, 0.0], [-1.0j, 2.0, 1.0e-20],
... [1+1j, 1.0e-9, -1.0]])
>>> print_matrix(M)
{ 1.00E+00,      0.0}( 2.00E+00,      0.0)(        0,        0)
(      0.0,-1.00E+00){ 2.00E+00,      0.0}(      0.0,      0.0)
( 1.00E+00, 1.00E+00)( 1.00E-09,      0.0){-1.00E+00,      0.0}
>>> print_matrix(M, limit=1.0e-5)
{ 1.00E+00,      0.0}( 2.00E+00,      0.0)(        0,        0)
(      0.0,-1.00E+00){ 2.00E+00,      0.0}(      0.0,      0.0)
( 1.00E+00, 1.00E+00)(      0.0,      0.0){-1.00E+00,      0.0}
>>> print_matrix(M, compress=True)
{1.00E+00,     0.0}(2.00E+00,     0.0)(       0,       0)
(    0.0,-1.00E+00){2.00E+00,     0.0}(     0.0,     0.0)
(1.00E+00,1.00E+00)(1.00E-09,     0.0){-1.00E+00,    0.0}
>>> print_matrix(M, compress=True, zero_as_blank=True)
{1.00E+00,     0.0}(2.00E+00,     0.0)(                 )
(    0.0,-1.00E+00){2.00E+00,     0.0}(     0.0,     0.0)
(1.00E+00,1.00E+00)(1.00E-09,     0.0){-1.00E+00,    0.0}
>>> M[2,1] = 1.0
>>> print_matrix(M, fmt="%5.1f")
{  1.0,  0.0}(  2.0,  0.0)(    0,    0)
(  0.0, -1.0){  2.0,  0.0}(  0.0,  0.0)
(  1.0,  1.0)(  1.0,  0.0){ -1.0,  0.0}
>>> def compact_exp_fmt(x):
...     if x == 0:
...         return '%7d' % 0
...     else:  # single-digit exponent
...         s = "%8.1e" % x
...         base, exp = s.split('e')
...         return base + 'e%+d' % int(exp)
>>> print_matrix(M, compress=True, zero_as_blank=True, fmt=compact_exp_fmt)
{1.0e+0,     0}(2.0e+0,     0)(             )
(    0,-1.0e+0){2.0e+0,     0}(     0,     0)
(1.0e+0,1.0e+0)(1.0e+0,     0){-1.0e+0,    0}
>>> print_matrix(M, matrix_name="M", fmt="%5.1f")
M = [
{  1.0,  0.0}(  2.0,  0.0)(    0,    0)
(  0.0, -1.0){  2.0,  0.0}(  0.0,  0.0)
(  1.0,  1.0)(  1.0,  0.0){ -1.0,  0.0}
]
>>> import scipy.sparse
>>> print_matrix(scipy.sparse.csr_matrix(M), matrix_name="M", fmt="%5.1f")
M = [
{  1.0,  0.0}(  2.0,  0.0)(    0,    0)
(  0.0, -1.0){  2.0,  0.0}(  0.0,  0.0)
(  1.0,  1.0)(  1.0,  0.0){ -1.0,  0.0}
]
qdyn.io.fix_fortran_exponent(num_str)[source]

In 3-digit exponents, Fortran drops the ‘E’. Return a string with the ‘E’ restored.

>>> print(fix_fortran_exponent("1.0-100"))
1.0E-100
>>> print(fix_fortran_exponent("1.0E-99"))
1.0E-99
qdyn.io.read_complex(s)[source]

Convert a string to a complex number

>>> read_complex("1.0 -2.0-100")
(1-2e-100j)
qdyn.io.read_real(str)[source]

Convert a string to a real number

This works for Fortran-formatted numbers with a missing ‘E’ sign

>>> read_real("-2.0-100")
-2e-100
qdyn.io.read_cmplx_array(filename, **kwargs)[source]

Read a complex array from a file. The file must contain two columns (real and imaginary part). This routine is equivalent to the Fortran QDYN read_cmplx_array routine

Parameters

Notes

You may use datablock() as a wrapper for fileanme in order to read a specific block from a file that contains multiple blocks

qdyn.io.read_cmplx_matrix(filename)[source]

Read in complex matrix from file (as written by the Fortran QDYN print_cmplx_matrix routine).

Return a two-dimensional, double precision complex Numpy array

Parameters

filename – str or file-like object from which to read gate, or file-like object with equivalent content

qdyn.io.datablock(filename, block, decoding=None)[source]

Iterator over the lines inside the file with the given filename inside the given block only. Blocks must be separated by two empty lines.

Parameters
  • filename (str) – Name of file from which to read

  • block (int) – One-based Index of the desired block. A value of less than 1 or more than the number of block available in the file will not cause an error, but result in an empty iterator

  • decoding (None or str) – By default, the resulting lines are byte strings (so that datablock() can wrap fname in numpy.genfromtxt()). If decoding is given different from None, the resulting strings will be decoded instead.

qdyn.io.write_cmplx_array(carray, filename, header=None, fmtstr='%25.16E', append=False, comment=None)[source]

Write a complex array to file. Equivalent to the Fortran QDYN write_cmplx_array routine. Two columns will be written to the output file (real and imaginary part of carray)

Parameters
  • filename (str) – Name of file to which to write the array

  • header (str or None) – A header line, to be written immediately before the data. Should start with a ‘#’

  • fmtstr (str or None) – The format to use for reach of the two columns

  • append (bool) – If True, append to existing files, with a separator of two blank lines

  • comment (str or None) – A comment line, to be written at the top of the file (before the header). Should start with a ‘#’

qdyn.io.writetotxt(fname, *args, **kwargs)[source]

Inverse function to numpy.genfromtxt and similar to numpy.savetxt, but allowing to write multiple numpy arrays as columns to a text file. Also, handle headers/footers more intelligently.

The first argument is the filename or handler to which to write, followed by an arbitrary number of numpy arrays, to be written as columns in the file (real arrays will produce once column, complex arrays two). The remaining keyword arguments are passed directly to numpy.savetxt (with fixes to the header/footer lines, as described below)

Parameters
  • fname (str) – filename or file handle

  • args (numpy.ndarray) – Numpy arrays to write to fname. All arrays must have the same length

  • fmt (str, list(str)) – A single format (e.g. ‘%10.5f’), a sequence of formats, or a multi-format string, e.g. ‘Iteration %d – %10.5f’, in which case delimiter is ignored. For a complex array in *args, a format for the real and imaginary parts must be given. Defaults to ‘%25.16E’.

  • delimiter (str) – Character separating columns. Defaults to ‘’

  • header (str, list(str)) – String that will be written at the beginning of the file. If sequence of strings, multiple lines will be written.

  • footer (str, list(str)) – String that will be written at the end of the file. If sequence of strings, multiple lines will be written

  • comments (str) – String that will be prepended to the each line of the header and footer strings, to mark them as comments. Defaults to ‘# ‘

Note

The header and footer lines are handled more intelligently than by the numpy.savetxt routine. First, header and footer may be an array of lines instead of just a (multiline) string. Second, each line in the header may or may not already include the comments prefix. If a line does not include the comments prefix yet, but starts with a sufficient number of spaces, the comments prefix will not be prepended to the line in output, but will overwrite the beginning of the line, so as not to change the line length. E.g. giving header=” time [ns]” will result in a header line of # time [ns] in the output, not # time [ns].

Further explanation of the fmt parameter (%[flag]width[.precision]specifier):

flags:

- : left justify

+ : Forces to precede result with + or -.

0Left pad the number with zeros instead of space (see

width).

width:

Minimum number of characters to be printed. The value is not truncated if it has more characters.

precision:
  • For integer specifiers (eg. d,i), the minimum number of digits.

  • For e, E and f specifiers, the number of digits to print after the decimal point.

  • For g and G, the maximum number of significant digits.

specifiers (partial list):

c : character

d or i : signed decimal integer

e or E : scientific notation with e or E.

f : decimal floating point

g,G : use the shorter of e,E or f

For more details, see numpy.savetxt

qdyn.io.read_ascii_dump(filename, convert_boolean=True, flatten=False)[source]

Read a file in the format written by the QDYN dump_ascii_* routines and return its data as a nested OrderedDict. The QDYN type of the dumped data structure is recorded in the typename attribute of the result.

Parameters
  • filename (str) – name of file from which to read data

  • convert_boolean (bool) – Convert strings ‘T’, ‘F’ to Boolean values True and False

  • flatten (bool) – If True, numerical array data is returned flattend (one-dimensional). If False, it is reshaped to the shape defined in the dump file.

qdyn.io.write_psi_amplitudes(psi, filename)[source]

Write the wave function to file in the same format as the write_psi_amplitudes Fortran routine

Parameters
  • psi (numpy.ndarray) – Array of complex probability amplitudes

  • filename (str) – Name of file to which to write

qdyn.io.read_psi_amplitudes(filename, n, block=1, normalize=True)[source]

Read the wave function of size n from file. For ‘block=1’, inverse to write_psi_amplitudes. Returns complex or real numpy array

By specifying blocks, data may be read from a file that contains multiple wave functions, in the format generated e.g. by the qdyn_prop_traj --write-all-states utility

Parameters
  • filename (str) – Name of file from which to read data

  • n (int) – dimension of the Hilbert space (i.e. size of returned vector)

  • block (int) – One-based index of the block to read from filename, if the file contains multiple block. Blocks must be separated by exactly two empty lines

  • normalize (bool) – Whether or not to normalize the wave function

qdyn.io.iterate_psi_amplitudes(filename, n, start_from_block=1, normalize=True)[source]

Iterate over blocks of wave functions stored in filename

This is equivalent to the following pseudocode:

iter([read_psi_amplitudes(filename, n, block, normalize)
      for block >= start_from_block])

Hoever, the data file is only traversed once, and memory is only required for one block.

Parameters
  • filename (str) – Name of file from which to read data

  • n (int) – dimension of the Hilbert space (i.e. size of returned vector)

  • start_from_block (int) – One-based index of block from which to start the iterator.

  • normalize (bool) – Whether or not to normalize the wave function