gibuu is hosted by Hepforge, IPPP Durham
GiBUU

TABLE OF CONTENTS


/histf90 [ Modules ]

[ Top ] [ Modules ]

NAME

module histf90

PURPOSE

Encapsulate all routines and datas for 1D Histograms.

Features of Histograms provided by this module:

  • store paramaters of the x-binning
  • enable two y-values (y and y2)
  • track keeping of under-/over-score the given extreme values of x.
  • provide simple-to-understand output routines (cf. WriteHist)
  • provide simple histogram arithmetic (not yet implemented)
  • many others...

Every Histogram prints his own multicolumn output. A multicolumn output of many different histograms for the same x-value is not implemented. This is done by the module "histMPf90".

INPUTS

...(This module needs no input)


histf90/histogram [ Types ]

[ Top ] [ histf90 ] [ Types ]

NAME

type histogram

PURPOSE

Type definition to store all information for a 1D Histogram.

SOURCE

  type, public :: histogram
     real             :: xMin        ! smallest x-value
     real             :: xMax        ! largest x-value
     real             :: xBin        ! width of x-Bin
     real             :: xExtreme(2) ! extremes of used x-values
                                     ! (min,max)
     character*(NameLength) :: Name  ! name to be written
     real,allocatable :: yVal(:,:)   ! histogram values:
                                     ! (x), (yy1,yy2,yy3...)
                                     ! bin 0,-1 : extreme values !!!
     logical :: initialized = .false.! flag to indicate whether the histogram has been initialized
  end type histogram

histf90/AddHist [ Subroutines ]

[ Top ] [ histf90 ] [ Subroutines ]

NAME

subroutine AddHist(H, x,y,y2)

subroutine AddHist(H1,H2, x,y,y2)

PURPOSE

Add to the given histogram(s) at the given x-value the weight y (y2). If y2 is given, also the second column is filled.

INPUTS

  • type(histogram) :: H -- Histogram to be used

or:

  • type(histogram) :: H1,H2 -- Histograms to be used
  • real :: x -- x-value
  • real :: y -- weight to be added
  • real :: y2 -- second weight to be added [OPTIONAL]

OUTPUT

H is changed

NOTES

This routine is overloaded: Giving two histograms at the same time is like calling the routine twice.

This overloading is useful, if you want to fill some histogram in an array, but also keep the sum of all histograms of this array in the histogram at position 0: call AddHist( (/hist(0), hist(iH)/), ...)


histf90/ClearHist [ Subroutines ]

[ Top ] [ histf90 ] [ Subroutines ]

NAME

subroutine ClearHist(H)

PURPOSE

Sets the histogram to zero again

INPUTS

  • type(histogram) :: H -- Histogram to be cleared

OUTPUT

H is changed


histf90/CreateHist [ Subroutines ]

[ Top ] [ histf90 ] [ Subroutines ]

NAME

subroutine CreateHist(H, HName, xmin, xmax, bin)

PURPOSE

This is the constructor of a 1D-Histogram! It allocates memory for the entries and puts additional variables to their default.

INPUTS

  • type(histogram) :: H -- Histogram to be created
  • character*(*) :: HName -- Name of Histogram
  • real :: xmin,xmax -- Minimal/maximal value for x-coordinate to be considered
  • real :: bin -- bin width

OUTPUT

H is changed


histf90/RemoveHist [ Subroutines ]

[ Top ] [ histf90 ] [ Subroutines ]

NAME

subroutine RemoveHist(H)

PURPOSE

Free the allocated memory

INPUTS

  • type(histogram) :: H -- Histogram to be used

OUTPUT

H is changed


histf90/WriteHeader [ Subroutines ]

[ Top ] [ histf90 ] [ Subroutines ]

NAME

subroutine WriteHeader(H,iFile,mul)

PURPOSE

A header is written:

  • Underflow: sum of all calls with a x-value less than x-min
  • Entries: sum of all entries listed in section "Data" below
  • Overflow: sum of all calls with a x-value above x-max
  • Average: average value of all bins
  • Extrema: the smallest/biggest x-values which had been added
  • counted: the range of x-values, which is considered in "Entries"

Summing "Underflow"+"Entries"+"Overflow" gives the number of ALL calls, i.e. the integral from -infty upto +infty.

Output of "Underflow","Entries","Overflow" is split into 3 columns, corresponding to the columns 2)..4) of the Data-Block, see below.

INPUTS

  • type(histogramMC) :: H -- Histogram
  • integer :: iFile -- file number to write to
  • real, optional :: mul

OUTPUT

write to file


histf90/WriteHist [ Subroutines ]

[ Top ] [ histf90 ] [ Subroutines ]

NAME

subroutine WriteHist(H, iFile, add, mul, DoAve, maxVal, H2, file, dump)

PURPOSE

Write out the histogram.

The entries are multiplied by 'mul' and 'add' is added.

INPUTS

  • type(histogram) :: H -- Histogram to be used
  • integer :: iFile -- File number output to redirect [OPTIONAL]
  • real :: add -- factor to add [OPTIONAL]
  • real :: mul -- factor to multiply [OPTIONAL]
  • logical :: DoAve -- write also "average" (cf Notes) [OPTIONAL]
  • real :: maxVal -- value used instead "divide by zero" [OPTIONAL]
  • type(histogram) :: H2 -- Histogram to divide by [OPTIONAL]
  • character*(*) :: file -- name of file to open and close [OPTIONAL]
  • logical :: dump -- if true, also dump it binary to file [OPTIONAL]

OUTPUT

Write to file number 'iFile'.

4 columns are written in the data section:

  • 1) x-value (i.e. middle of bin)
  • 2) y-value
  • 3) number of entries that lead to y-value
  • 4) y2-value (if used, otherwise 0)

Columns 2)..4) are divided by the bin-width.

If the (optional) parameter "add" was given, this value is added to the written values of columns 2)...4). (E.g. this is used, if one wants to create logplots with xmgr/grace: the value "0.0" is not allowed as input and destroys everything. Therefore calling this routine with the argument "add=1e-20" prohibits writing of "0.0"-values and everything is fine.)

If the (optional) parameter "mul" was given, all written values of columns 2)...4) are multiplied by this value. (E.g. this is used, if one wants to divide the output by a "number of runs" value: "mul=1./NumberOfRuns".)

Column 4) "y2" provides you with a simple way, to not calculate only one histogram "y", but simultanousely a very closed connected one. (E.g. if you want to calculate ratios "y2"/"y")

If DoAve is given: Write out the histogram, but now also divide column 3 by column 1. If column 2 is zero, i.e. no entries, then parameter maxVal is written instead. This column output is NOT divided by bin widths.

Plotting programs normally neglect the header lines and columns 3) and 4) -- leaving us with the use of these routines as expected.

If the (optional) parameter H2 is given, the first all entries of the histogram are divided by the entries of histogram H2, column 1.

If the (optional) parameter 'file' is given, this routine first opens the file with this name (using stream number iFile), rewinds it. The file is closed at the end of the routine.

NOTES

The histogram data is not affected!!!


histf90/WriteHist_Integrated [ Subroutines ]

[ Top ] [ histf90 ] [ Subroutines ]

NAME

subroutine WriteHist_Integrated(H, file, backward, normalize)

PURPOSE

Write out the histogram, integrating the data over x.

INPUTS

  • type(histogram) :: H -- Histogram to be used
  • character*(*) :: file -- File name to write to
  • logical :: backward -- integrate from front or back?
  • logical :: normalize -- normalize by total integral

OUTPUT

write to file

NOTES

The histogram data is not affected!!!


histf90/WriteHist_Spline [ Subroutines ]

[ Top ] [ histf90 ] [ Subroutines ]

NAME

subroutine WriteHist_Spline(H,file,add,mul,nPoints)

PURPOSE

Write out the histogram.

The entries are multiplied by 'mul' and 'add' is added.

The binning is subdivided into "nPoints" points per bin, where every point is calculated via cubic spline interpolation. (Default: nPoints=5)

INPUTS

  • type(histogram) :: H -- Histogram to be used
  • character*(*) :: file -- File name to write to
  • real :: add -- factor to add [OPTIONAL]
  • real :: mul -- factor to multiply [OPTIONAL]
  • integer :: nPoints -- number of sub-points [OPTIONAL]

OUTPUT

write to file number

cf. WriteHist

NOTES

The histogram data is not affected!!!


histf90/WriteHist_BSpline [ Subroutines ]

[ Top ] [ histf90 ] [ Subroutines ]

NAME

subroutine WriteHist_BSpline(H,file,add,mul,nPoints)

PURPOSE

Write out the histogram.

The entries are multiplied by 'mul' and 'add' is added.

The binning is subdivided into "nPoints" points per bin, where every point is calculated via cubic B-spline interpolation. (Default: nPoints=5)

INPUTS

  • type(histogram) :: H -- Histogram to be used
  • integer :: file -- File name to write to
  • real :: add -- factor to add [OPTIONAL]
  • real :: mul -- factor to multiply [OPTIONAL]
  • integer :: nPoints -- number of sub-points [OPTIONAL]

OUTPUT

write to file number

cf. WriteHist

NOTES

The histogram data is not affected!!!


histf90/WriteHist_Gauss [ Subroutines ]

[ Top ] [ histf90 ] [ Subroutines ]

NAME

subroutine WriteHist_Gauss(H,file,width_in)

PURPOSE

Write out the histogram smeared with a gaussian distribution

INPUTS

  • type(histogram) :: H -- Histogram to be used
  • character*(*) :: file -- name of file to write to
  • real, OPTIONAL :: width_in -- width of gaussian

OUTPUT

write to file number

NOTES

The histogram data is not affected!!!


histf90/WriteHist_Novo [ Subroutines ]

[ Top ] [ histf90 ] [ Subroutines ]

NAME

subroutine WriteHist_Novo(H,file,w,d)

PURPOSE

Write out the histogram smeared with a Novosibirsk distribution

INPUTS

  • type(histogram) :: H -- Histogram to be used
  • character*(*) :: file -- name of file to write to
  • real :: w,d -- width and skewness of Novosibirsk function

OUTPUT

write to file number

NOTES

The histogram data is not affected!!!


histf90/ReadHist [ Functions ]

[ Top ] [ histf90 ] [ Functions ]

NAME

logical function ReadHist(H, file)

PURPOSE

Read in a histogram from a file.

INPUTS

  • type(histogram) :: H -- Histogram to be used
  • character*(*) :: file -- name of file to read

OUTPUT

  • H is changed
  • return value indicates success or failure

NOTES

Only data points are read (no under-/overflow etc).


histf90/SumHist [ Subroutines ]

[ Top ] [ histf90 ] [ Subroutines ]

NAME

subroutine SumHist(A, B, w, doCheck)

PURPOSE

Perform a summation of two histograms: A = A + B. If a weight 'w' is given, do a weighted summation: A = A + w*B

INPUTS

  • type(histogram) :: A, B -- Histograms to be used
  • real, optional :: w -- optional weighting factor
  • logical, optional :: doCheck -- check array sizes or not

OUTPUT

A is changed.


histf90/CopyHist [ Subroutines ]

[ Top ] [ histf90 ] [ Subroutines ]

NAME

subroutine CopyHist(A, B)

PURPOSE

Perform a deep copy of two histograms: A = B.

INPUTS

  • type(histogram) :: A, B -- Histograms to be used

OUTPUT

A is changed.


histf90/DumpHist [ Subroutines ]

[ Top ] [ histf90 ] [ Subroutines ]

NAME

subroutine DumpHist(H,file,iFile, add,mul)

PURPOSE

Write all the histogram information unformatted (i.e. binary) to a file

INPUTS

  • type(histogram) :: H -- Histogram to be used
  • character*(*) :: file -- name of file to open and close
  • integer,OPTIONAL:: iFile -- File number output to redirect [OPTIONAL]
  • real :: add -- factor to add [OPTIONAL]
  • real :: mul -- factor to multiply [OPTIONAL]

OUTPUT

H is written UNFORMATTED to the given file


histf90/FetchHist [ Subroutines ]

[ Top ] [ histf90 ] [ Subroutines ]

NAME

subroutine FetchHist(H,file,iFile, add,mul,flagOK)

PURPOSE

Read in all the histogram information previously dumped unformatted (i.e. binary) to a file

INPUTS

  • character*(*) :: file -- name of file to open and close
  • integer,OPTIONAL:: iFile -- File number input to redirect [OPTIONAL]

OUTPUT

  • type(histogram) :: H -- Histogram to be used
  • real :: add -- factor to add [OPTIONAL]
  • real :: mul -- factor to multiply [OPTIONAL]
  • logical :: flagOK -- flag, if reading was okay [OPTIONAL]

H is read UNFORMATTED from the given file. Sizes are calculated as in CreateHist, also memory is allocated.

NOTES

No checks about input are performed!