Index: trunk/MagicSoft/Mars/Changelog
===================================================================
--- trunk/MagicSoft/Mars/Changelog	(revision 2161)
+++ trunk/MagicSoft/Mars/Changelog	(revision 2162)
@@ -1,3 +1,15 @@
                                                  -*-*- END OF LINE -*-*-
+
+ 2003/06/05: Abelardo Moralejo
+
+   * mhist/MHTrigLvl0.[h,cc]:
+     - added. This is intended to find "hot" pixels firing too often
+       or pixels firing too rarely. Very preliminar!
+
+   * macros/pixfirerate.C:
+     - added. An example on 
+
+   * mhist/Makefile, HistLinkDef.h :
+     added new class.
 
  2003/06/05: Thomas Bretz
Index: trunk/MagicSoft/Mars/macros/pixfirerate.C
===================================================================
--- trunk/MagicSoft/Mars/macros/pixfirerate.C	(revision 2162)
+++ trunk/MagicSoft/Mars/macros/pixfirerate.C	(revision 2162)
@@ -0,0 +1,133 @@
+/* ======================================================================== *\
+!
+! *
+! * This file is part of MARS, the MAGIC Analysis and Reconstruction
+! * Software. It is distributed to you in the hope that it can be a useful
+! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
+! * It is distributed WITHOUT ANY WARRANTY.
+! *
+! * Permission to use, copy, modify and distribute this software and its
+! * documentation for any purpose is hereby granted without fee,
+! * provided that the above copyright notice appear in all copies and
+! * that both that copyright notice and this permission notice appear
+! * in supporting documentation. It is provided "as is" without express
+! * or implied warranty.
+! *
+!
+!
+!   Author(s): Abelardo Moralejo <mailto:moralejo@pd.infn.it>
+!
+!   Copyright: MAGIC Software Development, 2000-2003
+!
+!
+\* ======================================================================== */
+
+/////////////////////////////////////////////////////////////////////////////
+//
+// pixfirerate.C
+//
+// This macro can help to  find "hot" pixels firing too often
+// or pixels firing too rarely. It plots camera displays showing how many 
+// times the FADC integral of each pixel has been found to be above pedestal
+// by 10, 20, 50, 100 or 200 ADC counts. As "pedestal"  we now use simply 
+// the signal in the first ADC slice, which seemed reasonable from a first
+// look at the available test data.
+//
+/////////////////////////////////////////////////////////////////////////////
+
+void pixfirerate(TString filename="20030603_01731_D_cosmics_E.root")
+{
+    //
+    // Update frequency by default = 1Hz
+    //
+    MStatusDisplay *d = new MStatusDisplay;
+
+    // Set update time to 5s
+    // d->SetUpdateTime(5000);
+
+    // Disable online update
+    // d->SetUpdateTime(-1);
+
+    d->SetLogStream(&gLog, kTRUE); // Disables output to stdout
+    gLog.SetOutputFile("status.log", kTRUE); // Enable output to file
+    //gLog.EnableOutputDevice(MLog::eStdout); // Enable output to stdout again
+
+    //
+    // Create a empty Parameter List and an empty Task List
+    // The tasklist is identified in the eventloop by its name
+    //
+    MParList  plist;
+    MTaskList tlist;
+    plist.AddToList(&tlist);
+
+    // The geometry container must be created by yourself to make sure
+    // that you don't choose a wrong geometry by mistake
+    //
+    MGeomCamMagic geomcam;
+    plist.AddToList(&geomcam);
+
+    //
+    // Now setup the tasks and tasklist:
+    // ---------------------------------
+    //
+    MReadMarsFile read("Events");
+    read.DisableAutoScheme();
+
+    MHTrigLvl0 trigmap1(10.,"Above 10");
+    MHTrigLvl0 trigmap2(20.,"Above 20");
+    MHTrigLvl0 trigmap3(50.,"Above 50");
+    MHTrigLvl0 trigmap4(100.,"Above 100");
+    MHTrigLvl0 trigmap5(200.,"Above 200");
+
+    plist.AddToList(&trigmap1);
+    plist.AddToList(&trigmap2);
+    plist.AddToList(&trigmap3);
+    plist.AddToList(&trigmap4);
+    plist.AddToList(&trigmap5);
+
+    MFillH fill1("Above 10","MRawEvtData");
+    MFillH fill2("Above 20","MRawEvtData");
+    MFillH fill3("Above 50","MRawEvtData");
+    MFillH fill4("Above 100","MRawEvtData");
+    MFillH fill5("Above 200","MRawEvtData");
+
+    read.AddFile(filename);
+
+    tlist.AddToList(&read);
+    tlist.AddToList(&fill1);
+    tlist.AddToList(&fill2);
+    tlist.AddToList(&fill3);
+    tlist.AddToList(&fill4);
+    tlist.AddToList(&fill5);
+
+    MEvtLoop evtloop;
+    evtloop.SetParList(&plist);
+    evtloop.SetDisplay(d);
+
+    //
+    // Execute your analysis
+    //
+    if (!evtloop.Eventloop())
+        return;
+
+    tlist.PrintStatistics();
+
+    //
+    // Make sure the display hasn't been deleted by the user while the
+    // eventloop was running.
+    //
+    if ((d = evtloop.GetDisplay()))
+    {
+        // Save data in a postscriptfile (status.ps)
+        d->SaveAsPS();
+        /*
+         * ----------- Write status to a root file ------------
+         *
+         TFile file("status.root", "RECREATE");
+         d->Write();
+         file.Close();
+         delete d;
+         */
+    }
+
+}
Index: trunk/MagicSoft/Mars/mhist/HistLinkDef.h
===================================================================
--- trunk/MagicSoft/Mars/mhist/HistLinkDef.h	(revision 2161)
+++ trunk/MagicSoft/Mars/mhist/HistLinkDef.h	(revision 2162)
@@ -41,5 +41,6 @@
 
 #pragma link C++ class MHCompProb+;
-#pragma link C++ class MHHadronness+;                                          
+#pragma link C++ class MHHadronness+;    
+#pragma link C++ class MHTrigLvl0+;                                
 
 #endif
Index: trunk/MagicSoft/Mars/mhist/MHTrigLvl0.cc
===================================================================
--- trunk/MagicSoft/Mars/mhist/MHTrigLvl0.cc	(revision 2162)
+++ trunk/MagicSoft/Mars/mhist/MHTrigLvl0.cc	(revision 2162)
@@ -0,0 +1,232 @@
+/* ======================================================================== *\
+!
+! *
+! * This file is part of MARS, the MAGIC Analysis and Reconstruction
+! * Software. It is distributed to you in the hope that it can be a useful
+! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
+! * It is distributed WITHOUT ANY WARRANTY.
+! *
+! * Permission to use, copy, modify and distribute this software and its
+! * documentation for any purpose is hereby granted without fee,
+! * provided that the above copyright notice appear in all copies and
+! * that both that copyright notice and this permission notice appear
+! * in supporting documentation. It is provided "as is" without express
+! * or implied warranty.
+! *
+!
+!
+!   Author(s): Abelardo Moralejo 06/2003 <mailto:moralejo@pd.infn.it>
+!
+!   Copyright: MAGIC Software Development, 2000-2003
+!
+!
+\* ======================================================================== */
+
+/////////////////////////////////////////////////////////////////////////////
+//
+// MHTrigLvl0
+//
+// This is intended to be a sort of "level 0 trigger display". What it really
+// does is to store the number of events of a data file in which each pixel 
+// has gone above a given threshold (fPixelThreshold) which is chosen when
+// calling the constructor. Displaying a camera view with these values can
+// help identify noisy pixels. See the macro pixfixrate.C to see an example
+// of its use. Many things are to be fixed. Only inner pixels are shown now
+// (which are anyhow those involved in the trigger), and the camera geometry
+// (number of pixels, and how many inner ones) is not yet read from the input 
+// file. 
+// The "pedestal" we are using is just the signal in the first ADC slice 
+// (seems reasonable from the inspection of the available test data files).
+//
+//
+/////////////////////////////////////////////////////////////////////////////
+#include "MHTrigLvl0.h"
+
+#include <TCanvas.h>
+
+#include "MLog.h"
+#include "MLogManip.h"
+
+#include "MParList.h"
+#include "MRawEvtData.h"
+#include "MRawEvtPixelIter.h"
+#include "MCamDisplay.h"
+
+#include "MGeomCam.h"
+#include "MGeomPix.h"
+
+ClassImp(MHTrigLvl0);
+
+// --------------------------------------------------------------------------
+//
+// Reset all pixels to 0 and reset fEntries to 0.
+//
+void MHTrigLvl0::Clear()
+{
+    fSum.Reset();
+    for (int i=0; i<577; i++)
+    {
+        fSum.AddPixel(i, 0, 0);
+        fSum[i].SetPixelUnused();
+    }
+    fSum.FixSize();
+    fEntries = 0;
+}
+
+// --------------------------------------------------------------------------
+//
+// Initialize the name and title of the task.
+// Resets the sum histogram
+//
+MHTrigLvl0::MHTrigLvl0(const Float_t pixelthreshold, 
+		       const char *name, const char *title)
+    : fCam(NULL), fRawEvt(NULL), fDispl(NULL)
+{
+    //
+    //   set the name and title of this object
+    //
+    fName  = name  ? name  : "MHTrigLvl0";
+    fTitle = title ? title : "Number of hits above per pixel";
+    fPixelThreshold = pixelthreshold;
+
+    Clear();
+}
+
+// --------------------------------------------------------------------------
+//
+// Delete the corresponding camera display if available
+//
+MHTrigLvl0::~MHTrigLvl0()
+{
+    if (fDispl)
+        delete fDispl;
+}
+
+// --------------------------------------------------------------------------
+//
+// Get the event (MRawEvtData) the histogram might be filled with. If
+// it is not given, it is assumed, that it is filled with the argument
+// of the Fill function.
+// Looks for the camera geometry MGeomCam and resets the sum histogram.
+//
+Bool_t MHTrigLvl0::SetupFill(const MParList *plist)
+{
+    fRawEvt = (MRawEvtData*)plist->FindObject("MRawEvtData");
+    if (!fRawEvt)
+    {
+        *fLog << dbginf << "MRawEvtData not found... aborting." << endl;
+        return kFALSE;
+    }
+
+    fCam = (MGeomCam*)plist->FindObject("MGeomCam");
+    if (!fCam)
+        *fLog << warn << GetDescriptor() << ": No MGeomCam found." << endl;
+
+    Clear();
+
+    return kTRUE;
+}
+
+// --------------------------------------------------------------------------
+//
+// Fill the histograms with data from a MCerPhotEvt-Container.
+//
+Bool_t MHTrigLvl0::Fill(const MParContainer *par, const Stat_t w)
+{
+    MRawEvtData *rawevt = par ? (MRawEvtData*) par : fRawEvt;
+    if (!rawevt)
+    {
+        *fLog << err << dbginf << "No MRawEvtData found..." << endl;
+        return kFALSE;
+    }
+
+    MRawEvtPixelIter pixel(rawevt);
+
+    while(pixel.Next())
+    {
+      const UInt_t pixid = pixel.GetPixelId();
+
+      // FIXME: number of inner pixels should be read from file
+      if (pixid > 396)
+	break;
+      if (pixid == 0)
+	continue;
+
+      fSum[pixid].SetPixelUsed();
+
+      //
+      // FIXME: we now use as "pedestal" the value of the first ADC slice!
+      //
+      Float_t baseline = rawevt->GetNumHiGainSamples() *
+	pixel.GetHiGainSamples()[0];
+
+      Float_t pixel_signal = pixel.GetSumHiGainSamples() - baseline;
+
+      Float_t pixel_is_on = ( pixel_signal > fPixelThreshold)? 1. : 0.;
+
+      fSum[pixid].AddNumPhotons(pixel_is_on);
+    }
+
+    fEntries++;
+
+    return kTRUE;
+}
+
+// --------------------------------------------------------------------------
+//
+// Set to Unused outer pixels.
+//
+Bool_t MHTrigLvl0::Finalize()
+{
+  //
+  // Show only pixels in the inner region of the camera:
+  // (otherwise, problem with the too different ranges) 
+  //
+  for (Int_t i=0; i<577; i++)
+    // FIXME: read number of total and inner pixels from file
+    {
+      if (i > 396)
+	fSum[i].SetPixelUnused();
+    }
+
+  //  fSum.Scale(fEntries); Now disabled, scale was not readable otherwise.
+
+  return kTRUE;
+}
+
+// --------------------------------------------------------------------------
+//
+// Draw the present 'fill status'
+//
+void MHTrigLvl0::Draw(Option_t *)
+{
+    if (!fCam)
+    {
+        *fLog << warn << "WARNING - Cannot draw " << GetDescriptor() << ": No Camera Geometry available." << endl;
+        return;
+    }
+
+    TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this, 750, 600);
+    pad->SetBorderMode(0);
+
+    AppendPad("");
+}
+
+// --------------------------------------------------------------------------
+//
+// If a camera display is not yet assigned, assign a new one.
+//
+void MHTrigLvl0::Paint(Option_t *option)
+{
+    if (!fCam)
+    {
+        *fLog << warn << "WARNING - Cannot paint " << GetDescriptor() << ": No Camera Geometry available." << endl;
+        return;
+    }
+
+    if (!fDispl)
+        fDispl = new MCamDisplay(fCam);
+
+    fDispl->FillPhotNum(fSum);
+    fDispl->Paint();
+}
Index: trunk/MagicSoft/Mars/mhist/MHTrigLvl0.h
===================================================================
--- trunk/MagicSoft/Mars/mhist/MHTrigLvl0.h	(revision 2162)
+++ trunk/MagicSoft/Mars/mhist/MHTrigLvl0.h	(revision 2162)
@@ -0,0 +1,46 @@
+#ifndef MARS_MHTrigLvl0
+#define MARS_MHTrigLvl0
+
+#ifndef MARS_MH
+#include "MH.h"
+#endif
+
+#ifndef MARS_MCerPhotEvt
+#include "MCerPhotEvt.h"
+#endif
+
+class TH1D;
+class MCamDisplay;
+class MRawEvtData;
+
+class MHTrigLvl0 : public MH
+{
+private:
+    MCerPhotEvt  fSum;            // storing the sum of triggers
+    Int_t        fEntries;        // number of entries in the histogram
+    MGeomCam    *fCam;            // the present geometry
+    MRawEvtData *fRawEvt;         //! ADC info of the current event
+    MCamDisplay *fDispl;          //! the camera display
+    Float_t      fPixelThreshold; // Threshold (ADC counts a.p.) to consider
+                                  // a pixel as "fired".
+
+public:
+    MHTrigLvl0(const Float_t pixelthreshold = 0, 
+	       const char *name=NULL, const char *title=NULL);
+    ~MHTrigLvl0();
+
+    void Clear();
+
+    Bool_t SetupFill(const MParList *pList);
+    Bool_t Fill(const MParContainer *par, const Stat_t w=1);
+    Bool_t Finalize();
+
+    const MCerPhotEvt &GetSum() const { return fSum; }
+
+    void Draw(Option_t *opt="");
+    void Paint(Option_t *option="");
+
+    ClassDef(MHTrigLvl0, 1) // Histogram to sum level 0 triggers in all pixels
+};
+
+#endif
Index: trunk/MagicSoft/Mars/mhist/Makefile
===================================================================
--- trunk/MagicSoft/Mars/mhist/Makefile	(revision 2161)
+++ trunk/MagicSoft/Mars/mhist/Makefile	(revision 2162)
@@ -58,5 +58,6 @@
 	   MHSigmaPixel.cc \
 	   MHSigmabarTheta.cc \
-	   MHSigmaTheta.cc
+	   MHSigmaTheta.cc \
+	   MHTrigLvl0.cc
 
 SRCS    = $(SRCFILES)
