gibuu is hosted by Hepforge, IPPP Durham
GiBUU

TABLE OF CONTENTS


/hist2Df90 [ Modules ]

[ Top ] [ Modules ]

NAME

module hist2Df90

PURPOSE

Encapsulate all routines and datas for 2D Histograms.

Features of Histograms provided by this module:

  • store paramaters of the (x1,x2)-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
  • provide simple histogram arithmetic
  • many others...

NOTES

Programming is fairly similar to "module histf90"

INPUTS

...(This module needs no input)


hist2Df90/histogram2D [ Types ]

[ Top ] [ hist2Df90 ] [ Types ]

NAME

type histogram2D

PURPOSE

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

SOURCE

  type histogram2D
     real             :: xMin(2)        ! smallest x-value: (x1,x2)
     real             :: xMax(2)        ! largest x-value:  (x1,x2)
     real             :: xBin(2)        ! width of x-Bin:   (x1,x2)
     real             :: xExtreme(2,2)  ! extremes of used x-values
                                        !  (x1,x2), (min,max)
     character*(NameLength) :: Name     ! name to be written
     real,allocatable :: yVal(:,:,:)    ! histogram values:
                                        !  (x1),(x2),(yy1,yy2,yy3)
  end type histogram2D

hist2Df90/AddHist2D [ Subroutines ]

[ Top ] [ hist2Df90 ] [ Subroutines ]

NAME

subroutine AddHist2D(H, x,y,y2)

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

PURPOSE

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

INPUTS

or:

  • type(histogram2D) :: H1,H2 -- Histograms to be used
  • real,dimension(2) :: x -- (x1,x2)-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.


hist2Df90/CreateHist2D [ Subroutines ]

[ Top ] [ hist2Df90 ] [ Subroutines ]

NAME

subroutine CreateHist2D(H, HName,x1,x2,bin,verbose)

PURPOSE

This is the Constructor of a 2D-Histogram! Allocate Memory for the entries and put additional variables to their default.

INPUTS

  • type(histogram2D) :: H -- Histogram to be created
  • character*(*) :: HName -- Name of Histogram
  • real,dimension(2) :: x1,x2,bin -- Minimal/maximal values for x-coordinates to be considered, bin-width
  • logical :: verbose -- switch for verbosity [OPTIONAL]

OUTPUT

H is changed


hist2Df90/WriteHist2D_SYSTEM [ Subroutines ]

[ Top ] [ hist2Df90 ] [ Subroutines ]

NAME

subroutine WriteHist2D_SYSTEM(H, iFile)

PURPOSE

Write out System information of the histogram

INPUTS

  • type(histogram2D) :: H -- Histogram to be used
  • integer :: iFile -- File number output to redirect

OUTPUT

written to file iFile


hist2Df90/RemoveHist2D [ Subroutines ]

[ Top ] [ hist2Df90 ] [ Subroutines ]

NAME

subroutine RemoveHist2D(H)

PURPOSE

Free the allocated memory

INPUTS

OUTPUT

H is changed


hist2Df90/CopyHist2D [ Subroutines ]

[ Top ] [ hist2Df90 ] [ Subroutines ]

NAME

subroutine CopyHist2D(H1,H2)

PURPOSE

Copies Histogram H1 into H2 (H2 is overwritten)

INPUTS

OUTPUT

H2 is changed


hist2Df90/WriteHist2D_Gnuplot [ Subroutines ]

[ Top ] [ hist2Df90 ] [ Subroutines ]

NAME

subroutine WriteHist2D_Gnuplot(H,iFile,add,mul, iColumn,DoAve,MaxVal,H2,file,dump,SwapXY)

PURPOSE

Write out the 2D-histogram. Format is suitable as input to gnuplots "splot" command

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

If iColumn is not given, all 3 columns are writen.

INPUTS

  • type(histogram2D) :: 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]
  • integer :: iColumn -- Column to write [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]
  • logical :: SwapXY -- write data sorted according y:x:z instead of x:y:z [OPTIONAL]

OUTPUT

write to file number

NOTES

The Histogram Data is not affected!!!

if DoAve is given: Write out the 2D-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.


hist2Df90/WriteHist2D_Gnuplot_Bspline [ Subroutines ]

[ Top ] [ hist2Df90 ] [ Subroutines ]

NAME

subroutine WriteHist2D_Gnuplot_Bspline(H,iFile,add,mul, iColumn,nPoints,DoAve,MaxVal)

PURPOSE

Write out the 2D-histogram. Format is suitable as input to gnuplots "splot" command

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

If iColumn is not given, all 3 columns are writen.

INPUTS

  • type(histogram2D) :: H -- Histogram to be used
  • integer :: iFile -- File number output to redirect
  • real :: add -- factor to add [OPTIONAL]
  • real :: mul -- factor to multiply [OPTIONAL]
  • integer :: iColumn -- Column to write [OPTIONAL]
  • integer :: nPoints -- number of sub-points [OPTIONAL]
  • logical :: DoAve -- write also "average" (cf Notes) [OPTIONAL]
  • real :: maxVal -- value used instead "divide by zero" [OPTIONAL]

OUTPUT

write to file number

NOTES

The Histogram Data is not affected!!!

cf. WriteHist2D_Gnuplot


hist2Df90/DivideHist2D [ Subroutines ]

[ Top ] [ hist2Df90 ] [ Subroutines ]

NAME

subroutine DivideHist2D(H1,H2, maxVal, MulWithBin)

PURPOSE

calculate H1 = H1/H2 on a bin-by-bin basis

INPUTS

  • type(histogram2D) :: H1,H2 -- Histograms to be used
  • real :: maxVal -- value used instead "divide by zero" [OPTIONAL]
  • logical :: MulWithBin -- flag: multiply with bin-width [OPTIONAL]

Without explicit "MulWithBin" parameter, the default behaviour is as with "MulWithBin=.TRUE."

OUTPUT

H1 is changed

NOTES

This Dividing-Routine ignores all EXTREME values etc. !

The "MulWithBin=.TRUE." feature is provided, because during output the bin entry values are divided by the bin-width (as usual). Therefore if you divide a histogram by itself and writ it out, you would not get the bin-Values "1" but "1/bin-width". This flag here corrects for this.


hist2Df90/ReadHist2D_Gnuplot [ Subroutines ]

[ Top ] [ hist2Df90 ] [ Subroutines ]

NAME

subroutine ReadHist2D_Gnuplot(H,iFile,add,mul, iColumn, DoAve)

PURPOSE

Read in a 2D-histogram from a file, which has been written in a gnuplot format.


hist2Df90/DumpHist2D [ Subroutines ]

[ Top ] [ hist2Df90 ] [ Subroutines ]

NAME

subroutine DumpHist2D(H,file,iFile)

PURPOSE

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

INPUTS

  • type(histogram2D):: 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


hist2Df90/FetchHist2D [ Subroutines ]

[ Top ] [ hist2Df90 ] [ Subroutines ]

NAME

subroutine FetchHist2D(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(histogram2D) :: H -- Histogram to be used
  • 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!


hist2Df90/IntegrateHist2D [ Subroutines ]

[ Top ] [ hist2Df90 ] [ Subroutines ]

NAME

subroutine IntegrateHist2D(H2D,H,axis)

PURPOSE

Create a 1D histogram, which represents the integration along one axis of the 2D histogram

INPUTS

  • type(histogram2D) :: H2D -- 2D histogram to integrate
  • integer :: axis -- axis to preserve (integration is done for the other axis)

OUTPUT


hist2Df90/AverageHist2D [ Subroutines ]

[ Top ] [ hist2Df90 ] [ Subroutines ]

NAME

subroutine AverageHist2D(H2D,H,axis)

PURPOSE

Create a 1D histogram, which represents the average along one axis of the 2D histogram

INPUTS

  • type(histogram2D) :: H2D -- 2D histogram to integrate
  • integer :: axis -- axis to preserve (integration is done for the other axis)

OUTPUT