mirror of
https://github.com/fjosw/pyerrors.git
synced 2025-03-15 06:40:24 +01:00
docs: documentation for matrix operations extended, docstrings in
input.misc fixed
This commit is contained in:
parent
9aab654256
commit
05d36dc70b
2 changed files with 23 additions and 15 deletions
|
@ -161,7 +161,7 @@ Passing arguments to the `gamma_method` still dominates over the dictionaries.
|
|||
|
||||
## Irregular Monte Carlo chains
|
||||
|
||||
Irregular Monte Carlo chains can be initialized with the parameter `idl`.
|
||||
`Obs` objects defined on irregular Monte Carlo chains can be initialized with the parameter `idl`.
|
||||
|
||||
```python
|
||||
# Observable defined on configurations 20 to 519
|
||||
|
@ -187,6 +187,8 @@ obs3.details()
|
|||
|
||||
```
|
||||
|
||||
`Obs` objects defined on regular and irregular histories of the same ensemble can be computed with each other and the correct error propagation and estimation is automatically taken care of.
|
||||
|
||||
**Warning:** Irregular Monte Carlo chains can result in odd patterns in the autocorrelation functions.
|
||||
Make sure to check the autocorrelation time with e.g. `pyerrors.obs.Obs.plot_rho` or `pyerrors.obs.Obs.plot_tauint`.
|
||||
|
||||
|
@ -241,7 +243,7 @@ my_new_corr = 0.3 * my_corr[2] * my_corr * my_corr + 12 / my_corr
|
|||
- `Corr.correlate` constructs a disconnected correlation function from the correlator and another `Corr` or `Obs` object.
|
||||
- `Corr.reweight` reweights the correlator.
|
||||
|
||||
`pyerrors` can also handle matrices of correlation functions and extract energy states from these matrices via a generalized eigenvalue problem (see `pyerrors.correlators.Corr.GEVP).
|
||||
`pyerrors` can also handle matrices of correlation functions and extract energy states from these matrices via a generalized eigenvalue problem (see `pyerrors.correlators.Corr.GEVP`).
|
||||
|
||||
For the full API see `pyerrors.correlators.Corr`.
|
||||
|
||||
|
@ -273,17 +275,26 @@ print(my_derived_cobs)
|
|||
`pyerrors.roots`
|
||||
|
||||
# Matrix operations
|
||||
`pyerrors.linalg`
|
||||
`pyerrors` provides wrappers for `Obs`-valued matrix operations based on `numpy.linalg`. The supported functions include:
|
||||
- `inv` for the matrix inverse.
|
||||
- `cholseky` for the Cholesky decomposition.
|
||||
- `det` for the matrix determinant.
|
||||
- `eigh` for eigenvalues and eigenvectors of hermitean matrices.
|
||||
- `eig` for eigenvalues of general matrices.
|
||||
- `pinv` for the Moore-Penrose pseudoinverse.
|
||||
- `svd` for the singular-value-decomposition.
|
||||
|
||||
For the full API see `pyerrors.linalg`.
|
||||
|
||||
# Export data
|
||||
The preferred exported file format within `pyerrors` is
|
||||
The preferred exported file format within `pyerrors` is json.gz
|
||||
|
||||
## Jackknife samples
|
||||
For comparison with other analysis workflows `pyerrors` can generate jackknife samples from an `Obs` object or import jackknife samples into an `Obs` object.
|
||||
See `pyerrors.obs.Obs.export_jackknife` and `pyerrors.obs.import_jackknife` for details.
|
||||
|
||||
# Input
|
||||
`pyerrors.input`
|
||||
`pyerrors` includes an `input` submodule in which input routines and parsers for the output of various numerical programs are contained. For details see `pyerrors.input`.
|
||||
'''
|
||||
from .obs import *
|
||||
from .correlators import *
|
||||
|
|
|
@ -1,6 +1,3 @@
|
|||
#!/usr/bin/env python
|
||||
# coding: utf-8
|
||||
|
||||
import os
|
||||
import fnmatch
|
||||
import re
|
||||
|
@ -12,11 +9,12 @@ from ..obs import Obs
|
|||
def read_pbp(path, prefix, **kwargs):
|
||||
"""Read pbp format from given folder structure. Returns a list of length nrw
|
||||
|
||||
Keyword arguments
|
||||
-----------------
|
||||
r_start -- list which contains the first config to be read for each replicum
|
||||
r_stop -- list which contains the last config to be read for each replicum
|
||||
|
||||
Parameters
|
||||
----------
|
||||
r_start : list
|
||||
list which contains the first config to be read for each replicum
|
||||
r_stop : list
|
||||
list which contains the last config to be read for each replicum
|
||||
"""
|
||||
|
||||
extract_nfct = 1
|
||||
|
@ -66,7 +64,6 @@ def read_pbp(path, prefix, **kwargs):
|
|||
tmp_array = []
|
||||
with open(path + '/' + ls[rep], 'rb') as fp:
|
||||
|
||||
# header
|
||||
t = fp.read(4) # number of reweighting factors
|
||||
if rep == 0:
|
||||
nrw = struct.unpack('i', t)[0]
|
||||
|
@ -74,7 +71,7 @@ def read_pbp(path, prefix, **kwargs):
|
|||
deltas.append([])
|
||||
else:
|
||||
if nrw != struct.unpack('i', t)[0]:
|
||||
raise Exception('Error: different number of reweighting factors for replicum', rep)
|
||||
raise Exception('Error: different number of factors for replicum', rep)
|
||||
|
||||
for k in range(nrw):
|
||||
tmp_array.append([])
|
||||
|
|
Loading…
Add table
Reference in a new issue