README.md 8.99 KB
Newer Older
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

# aderrors - Error analysis of Monte Carlo data with Automatic Differentiation

`aderrors` is a fortran implementation of the $`\Gamma`$-method for
analysis of Monte Carlo data. It uses Automatic Differentiation to
perform **exact linear error propagation**. It preforms the computation
of gradients and Hessians of arbitrary functions, which allows a
robust and **exact error propagation even in iterative algorithms**. 

- [Features](#features)
- [Examples](#examples)
  - [Simple analysis of MC data](#simple-analysis-of-mc-data)
  - [A complete example](#a-complete-example)
  - [A calculator with uncertainties](#a-calculator-with-uncertainties)
- [Installation](#installation)
  - [Requirements](#requirements)
  - [Instructions](#instructions)
  - [Using the library](#using-the-library)
- [Full documentation](#full-documentation)
- [How to cite](#how-to-cite)

## Features

- **Exact** linear error propagation, even in iterative algorithms
  (i.e. error propagation in fit parameters).
- Handles data from **any number of ensembles** (i.e. simulations with
  different parameters).
- Support for **replicas** (i.e. several runs with the same simulation
  parameters). 
30
- Irrelagular MC measurements are handled transparently.
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 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265
- Standalone **portable** implementation without any external
  dependencies.
- **Fast** computation of autocorrelation functions with the [FFT
  package](http://www.kurims.kyoto-u.ac.jp/~ooura/fft.html) (included
  in the distribution).
- **Exact** determination of gradients and Hessians of arbitrary
  functions. 

## Examples 

This is just a small collection of examples. Note that the basic data type is
`uwreal`, that is able to handle MC histories and data with
errors. The error is determined by calling the method `uwerr`
on the data type. More examples can be found in the `test` directory
of the distribution. They are explained in the documentation `doc/aderrors.pdf`.

All these examples use the `module simulator` (available in
`test/simulator.f90`) to generate autocorrelated data with
autocorrelation function $`\Gamma(t) = \sum_k \lambda_k
e^{-|t|/\tau_k}`$. 

### Simple analysis of MC data

This code performs a simple error analysis of the data with the default
parameters ($`S_{\tau}=4`$ for automatic windowing, no tail added to the
autocorrelation function). The result of the analysis with the
$`\Gamma`$-method is compared with the exact values of the error and
$`\tau_{\rm int}`$ returned by the `module simulator`. This example is
included in the distribution in the file `test/simple.f90`.

```fortran
program simple

  use ISO_FORTRAN_ENV, Only : error_unit, output_unit
  use numtypes
  use aderrors
  use simulator
  
  implicit none

  integer, parameter  :: nd = 20000
  type (uwreal)       :: x
  real (kind=DP)      :: data_x(nd), err, ti, tau(5), lam(5)


  ! Fill arrays data_x(:) with autocorrelated
  ! data from the module simulator.
  tau = (/1.0_DP, 3.0_DP, 4.0_DP, 5.0_DP, 7.34_DP/)
  lam = (/1.00_DP, 0.87_DP, 1.23_DP, 0.56_DP, 0.87_DP/)

  call gen_series(data_x,  err, ti, tau, lam, 0.3_DP)

  ! Load data_x(:) measurements in variable x. Use
  ! default settings (Stau=4, texp=0, 1 replica)
  x = data_x
  
  ! Perform error analysis (optimal window)
  call x%uwerr() 

  ! Print results
  write(*,'(1A,1I6,1A)')'** Measurements: ', nd, ' **   '
  write(*,100)'  - Gamma-method:  ', x%value(),  " +/- ", x%error(), '( tauint: ', &
       x%taui(1), " +/- ", x%dtaui(1), ')'
  write(*,100)'  - Exact:         ', 0.3_DP,     " +/- ", err, '( tauint: ', ti, " )"

100 FORMAT((2X,1A,1F8.5,1A,1F7.5,5X,1A,1F0.2,1A,1F7.5,1A))
  
  
  stop
end program simple
```

Running this code gives as output
```
** Measurements:  20000 **   
    - Gamma-method:   0.29883 +/- 0.04271     ( tauint: 4.13 +/- 0.48344)
    - Exact:          0.30000 +/- 0.04074     ( tauint: 3.82 )
```

### A complete example

This example is included in the distribution in the file
`test/complete.f90`. The analysis is performed
for a complicated non-linear function of the two primary
observables. It also uses replica for the MC ensemble labeled 1. 

```fortran
program complete

  use ISO_FORTRAN_ENV, Only : error_unit, output_unit
  use numtypes
  use constants
  use aderrors
  use simulator
  
  implicit none

  integer, parameter  :: nd = 5000, nrep=4
  type (uwreal)       :: x, y, z
  integer             :: iflog, ivrep(nrep)=(/1000,30,3070,900/), i, is, ie
  real (kind=DP)      :: data_x(nd), data_y(nd/2), err, ti, texp
  real (kind=DP)      :: tau(4), &
       lam_x(4)=(/1.0_DP, 0.70_DP, 0.40_DP, 0.40_DP/), &
       lam_y(4)=(/2.3_DP, 0.40_DP, 0.20_DP, 0.90_DP/)
  character (len=200) :: flog='history_z.log'


  ! Fill arrays data_x(:) with autocorrelated
  ! data from the module simulator. Use nrep replica
  tau = (/1.0_DP, 3.0_DP, 12.0_DP, 75.0_DP/)
  texp = maxval(tau)

  is = 1
  do i = 1, nrep
     ie = is + ivrep(i) - 1
     call gen_series(data_x(is:ie),  err, ti, tau, lam_x, 0.3_DP)
     is = ie + 1
  end do

  ! Fill data_y(:) with different values of tau also using
  ! module simulator
  forall (i=1:4) tau(i) = real(2*i,kind=DP)
  call gen_series(data_y, err, ti, tau, lam_y, 1.3_DP)


  ! Load data_x(:) measurements in variable x.
  ! Set replica vector, exponential autocorrelation time
  ! and ensemble ID.
  x = data_x
  call x%set_id(1)
  call x%set_replica(ivrep)
  call x%set_texp(texp)
  
  ! Load data_y(:) measurements in variable y
  y = data_y
  call y%set_id(2)

  ! Exact, transparent error propagation
  z = sin(x)/(cos(y) + 1.0_DP)

  ! Attach tail in ensemble with ID 1 when signal in the
  ! normalized auto-correlation function equals its error
  call z%set_dsig(1.0_DP,1) 
  ! Set Stau=3 for automatic window in ensemble with ID 2
  call z%set_stau(3.0_DP,2)
  ! Perform error analysis (tails, optimal window,...)
  call z%uwerr() 

  ! Print results and output details to flog  
  write(*,'(1A,1F8.5,1A,1F8.5)')'** Observable z:   ', z%value(),  " +/- ", z%error()
  do i = 1, z%neid()
     write(*,'(3X,1A,1I3,3X,1F5.2,"%")',advance="no")'Contribution to error from ensemble ID', &
          z%eid(i), 100.0_DP*z%error_src(i)
     write(*,'(2X,1A,1F0.4,1A,1F8.4,1A)')'(tau int: ', z%taui(i), " +/- ", z%dtaui(i), ")"
  end do

  open(newunit=iflog, file=trim(flog))
  call z%print_hist(iflog)
  close(iflog)
  
  stop
end program complete
```

Running this code gives as output
```
** Observable z:    0.24426 +/-  0.05374
   Contribution to error from ensemble ID  1   83.93%  (tau int: 5.8333 +/-   2.0772)
   Contribution to error from ensemble ID  2   16.07%  (tau int: 2.5724 +/-   0.5268)
```


The file `history_z.log` contains the MC histories, normalized
autocorrelation functions and $`\tau_{\rm int}`$. The format is
a simple text file: it can be processed by shell scripts
(see `tools/plot/plot_hist.sh`) to produce the following graphics. 

<img src="plots/ID1-hist.png" alt="ID1-hist" width="800"/>
</br>
<img src="plots/ID1-auto.png" alt="ID1-hist" width="400"/>
<img src="plots/ID1-taui.png" alt="ID1-hist" width="400"/>

### A calculator with uncertainties

The module also handles simple data with errors (i.e. **not** a MC
history). This is extremely useful since in many occasions we have to
use data from statistically independent sources where only the central
value and the error is available. 

With the module `aderrors` these kind of data in handled
transparently, just as if it were another ensemble (see `test/calculator.f90`).

```fortran
program calculator

  use aderrors

  implicit none

  type (uwreal) :: x, y, z, t
  
  x = (/1.223_8, 0.012_8/) ! x = 1.223 +/- 0.012
  y = sin(2.0_8*x)
  z = 1.0_8 + 2.0_8 * sin(x)*cos(x)
  t = z - y
  call t%uwerr()
  write(*,'(1A,1F18.16,1A,1F18.16)')'Exactly one: ', t%value(), " +/- ", t%error() ! 1.0 +/- 0.0
  
  stop
end program calculator
```

Running this code gives as output
```
Exactly one: 1.0000000000000000 +/- 0.0000000000000000
```

## Installation

### Requirements

The code is strict `fortran 2008` compliant. Any compiler that
supports this standard can be used. The code has been tested with
`gfortran 6.X`, `gfortran 7.X`, `gfortran 8.X`, `intel fortran 17`,
`intel fortran 18`. 

Note that `gfortran 4.X` and `gfortran 5.X` do not support
`submodules` (part of the `fortran 2008` standard). This code will not
work whith these versions. 

### Instructions

1. Download or clone the repository.
1. Edit the `Makefile` in the `build` directory. Change the compiler
   command/options (variables `FC` and `FOPT`). 
266 267
1. Compile the library with `GNU make`.
1. Optionally build/run the test codes with `make test`. Executabes will
268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283
   be placed in the `test` directory.
1. If preferred, move the contents of the `include` and `lib`
   directories somewhere else.

### Using the library

1. Compile your programs with `-I <dir>/include`.
1. Link your programs with the `-L <dir>/lib` and `-laderr` options. 

## Full documentation

Look into the `doc/aderrors.pdf` file.

## How to cite

If you use this package for a scientific publication, please cite the
284 285 286 287 288
original work:

"Automatic differentiation for error analysis of Monte Carlo data"
Alberto Ramos. 
[arXiv:1809.01289](https://arxiv.org/abs/1809.01289)