| 1 | /* ======================================================================== *\
|
|---|
| 2 | !
|
|---|
| 3 | ! *
|
|---|
| 4 | ! * This file is part of MARS, the MAGIC Analysis and Reconstruction
|
|---|
| 5 | ! * Software. It is distributed to you in the hope that it can be a useful
|
|---|
| 6 | ! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
|
|---|
| 7 | ! * It is distributed WITHOUT ANY WARRANTY.
|
|---|
| 8 | ! *
|
|---|
| 9 | ! * Permission to use, copy, modify and distribute this software and its
|
|---|
| 10 | ! * documentation for any purpose is hereby granted without fee,
|
|---|
| 11 | ! * provided that the above copyright notice appear in all copies and
|
|---|
| 12 | ! * that both that copyright notice and this permission notice appear
|
|---|
| 13 | ! * in supporting documentation. It is provided "as is" without express
|
|---|
| 14 | ! * or implied warranty.
|
|---|
| 15 | ! *
|
|---|
| 16 | !
|
|---|
| 17 | !
|
|---|
| 18 | ! Author(s): Thomas Bretz, 1/2004 <mailto:tbretz@astro.uni-wuerzburg.de>
|
|---|
| 19 | !
|
|---|
| 20 | ! Copyright: MAGIC Software Development, 2000-2004
|
|---|
| 21 | !
|
|---|
| 22 | !
|
|---|
| 23 | \* ======================================================================== */
|
|---|
| 24 |
|
|---|
| 25 | /////////////////////////////////////////////////////////////////////////////
|
|---|
| 26 | //
|
|---|
| 27 | // MJStar
|
|---|
| 28 | //
|
|---|
| 29 | /////////////////////////////////////////////////////////////////////////////
|
|---|
| 30 | #include "MJStar.h"
|
|---|
| 31 |
|
|---|
| 32 | #include <TEnv.h>
|
|---|
| 33 | #include <TFile.h>
|
|---|
| 34 |
|
|---|
| 35 | #include "MLog.h"
|
|---|
| 36 | #include "MLogManip.h"
|
|---|
| 37 |
|
|---|
| 38 | #include "MDirIter.h"
|
|---|
| 39 | #include "MParList.h"
|
|---|
| 40 | #include "MTaskList.h"
|
|---|
| 41 | #include "MEvtLoop.h"
|
|---|
| 42 |
|
|---|
| 43 | #include "MStatusDisplay.h"
|
|---|
| 44 |
|
|---|
| 45 | #include "MH3.h"
|
|---|
| 46 | #include "MHVsTime.h"
|
|---|
| 47 | #include "MHCamEvent.h"
|
|---|
| 48 | #include "MHCamEventRot.h"
|
|---|
| 49 | #include "MBinning.h"
|
|---|
| 50 |
|
|---|
| 51 | #include "MReadReports.h"
|
|---|
| 52 | #include "MReadMarsFile.h"
|
|---|
| 53 | #include "MFDeltaT.h"
|
|---|
| 54 | #include "MContinue.h"
|
|---|
| 55 | #include "MGeomApply.h"
|
|---|
| 56 | #include "MEventRateCalc.h"
|
|---|
| 57 | #include "MImgCleanStd.h"
|
|---|
| 58 | #include "MHillasCalc.h"
|
|---|
| 59 | #include "MFillH.h"
|
|---|
| 60 | #include "MWriteRootFile.h"
|
|---|
| 61 |
|
|---|
| 62 | #include "MObservatory.h"
|
|---|
| 63 | #include "MPointingPosCalc.h"
|
|---|
| 64 |
|
|---|
| 65 | ClassImp(MJStar);
|
|---|
| 66 |
|
|---|
| 67 | using namespace std;
|
|---|
| 68 |
|
|---|
| 69 | // --------------------------------------------------------------------------
|
|---|
| 70 | //
|
|---|
| 71 | // Default constructor.
|
|---|
| 72 | //
|
|---|
| 73 | // Sets fRuns to 0, fExtractor to NULL, fDataCheck to kFALSE
|
|---|
| 74 | //
|
|---|
| 75 | MJStar::MJStar(const char *name, const char *title)
|
|---|
| 76 | {
|
|---|
| 77 | fName = name ? name : "MJStar";
|
|---|
| 78 | fTitle = title ? title : "Standard analysis and reconstruction";
|
|---|
| 79 | }
|
|---|
| 80 |
|
|---|
| 81 | Bool_t MJStar::WriteResult()
|
|---|
| 82 | {
|
|---|
| 83 | if (fPathOut.IsNull())
|
|---|
| 84 | {
|
|---|
| 85 | *fLog << inf << "No output path specified via SetPathOut - no output written." << endl;
|
|---|
| 86 | return kTRUE;
|
|---|
| 87 | }
|
|---|
| 88 |
|
|---|
| 89 | const TString oname = Form("%s/star%06d.root", (const char*)fPathOut, fSequence.GetSequence());
|
|---|
| 90 |
|
|---|
| 91 | *fLog << inf << "Writing to file: " << oname << endl;
|
|---|
| 92 |
|
|---|
| 93 | TFile file(oname, "RECREATE");
|
|---|
| 94 |
|
|---|
| 95 | *fLog << inf << " - MStatusDisplay..." << flush;
|
|---|
| 96 | if (fDisplay && fDisplay->Write()<=0)
|
|---|
| 97 | {
|
|---|
| 98 | *fLog << err << "Unable to write MStatusDisplay to " << oname << endl;
|
|---|
| 99 | return kFALSE;
|
|---|
| 100 | }
|
|---|
| 101 | *fLog << inf << "ok." << endl;
|
|---|
| 102 |
|
|---|
| 103 | return kTRUE;
|
|---|
| 104 | }
|
|---|
| 105 |
|
|---|
| 106 | Bool_t MJStar::ProcessFile(Bool_t ismc)
|
|---|
| 107 | {
|
|---|
| 108 | if (!fSequence.IsValid())
|
|---|
| 109 | {
|
|---|
| 110 | *fLog << err << "ERROR - Sequence invalid!" << endl;
|
|---|
| 111 | return kFALSE;
|
|---|
| 112 | }
|
|---|
| 113 |
|
|---|
| 114 | //if (!CheckEnv())
|
|---|
| 115 | // return kFALSE;
|
|---|
| 116 |
|
|---|
| 117 | CheckEnv();
|
|---|
| 118 |
|
|---|
| 119 | // --------------------------------------------------------------------------------
|
|---|
| 120 |
|
|---|
| 121 | *fLog << inf;
|
|---|
| 122 | fLog->Separator(GetDescriptor());
|
|---|
| 123 | *fLog << "Calculate image parameters from sequence ";
|
|---|
| 124 | *fLog << fSequence.GetName() << endl;
|
|---|
| 125 | *fLog << endl;
|
|---|
| 126 |
|
|---|
| 127 | // --------------------------------------------------------------------------------
|
|---|
| 128 |
|
|---|
| 129 | MDirIter iter;
|
|---|
| 130 | const Int_t n0 = fSequence.SetupDatRuns(iter, fPathData, "Y");
|
|---|
| 131 | const Int_t n1 = fSequence.GetNumDatRuns();
|
|---|
| 132 | if (n0==0)
|
|---|
| 133 | {
|
|---|
| 134 | *fLog << err << "ERROR - No input files of sequence found!" << endl;
|
|---|
| 135 | return kFALSE;
|
|---|
| 136 | }
|
|---|
| 137 | if (n0!=n1)
|
|---|
| 138 | {
|
|---|
| 139 | *fLog << err << "ERROR - Number of files found (" << n0 << ") doesn't match number of files in sequence (" << n1 << ")" << endl;
|
|---|
| 140 | return kFALSE;
|
|---|
| 141 | }
|
|---|
| 142 |
|
|---|
| 143 | // Setup Parlist
|
|---|
| 144 | MParList plist;
|
|---|
| 145 | plist.AddToList(this); // take care of fDisplay!
|
|---|
| 146 |
|
|---|
| 147 | MObservatory obs;
|
|---|
| 148 | plist.AddToList(&obs);
|
|---|
| 149 |
|
|---|
| 150 | // Setup Tasklist
|
|---|
| 151 | MTaskList tlist;
|
|---|
| 152 | plist.AddToList(&tlist);
|
|---|
| 153 |
|
|---|
| 154 | MReadReports readreal;
|
|---|
| 155 | readreal.AddTree("Events", "MTime.", kTRUE);
|
|---|
| 156 | readreal.AddTree("Drive");
|
|---|
| 157 | //read.AddTree("Trigger");
|
|---|
| 158 | //read.AddTree("Camera");
|
|---|
| 159 | //read.AddTree("CC");
|
|---|
| 160 | //read.AddTree("Currents");
|
|---|
| 161 | readreal.AddFiles(iter);
|
|---|
| 162 |
|
|---|
| 163 | MReadMarsFile readmc("Events");
|
|---|
| 164 | readmc.DisableAutoScheme();
|
|---|
| 165 | readmc.AddFiles(iter);
|
|---|
| 166 |
|
|---|
| 167 | // ------------------ Setup general tasks ----------------
|
|---|
| 168 |
|
|---|
| 169 | MFDeltaT fdeltat;
|
|---|
| 170 | MContinue cont(&fdeltat, "FilterDeltaT", "Filter events with wrong timing");
|
|---|
| 171 | cont.SetInverted();
|
|---|
| 172 |
|
|---|
| 173 | MGeomApply apply; // Only necessary to craete geometry
|
|---|
| 174 | MEventRateCalc rate;
|
|---|
| 175 | //MEventRateCalc rate1; // 5min
|
|---|
| 176 | rate.SetNumEvents(1200);
|
|---|
| 177 | //rate1.SetNumEvents(60000);
|
|---|
| 178 | //rate1.SetNameEventRate("MEventRate2");
|
|---|
| 179 | //rate1.SetNameTimeRate("MTimeRate2");
|
|---|
| 180 |
|
|---|
| 181 | /*
|
|---|
| 182 | MEventRateCalc rate10000;
|
|---|
| 183 | rate10000.SetNameEventRate("MEventRate10000");
|
|---|
| 184 | rate10000.SetNumEvents(10000);
|
|---|
| 185 | */
|
|---|
| 186 | //MBadPixelsMerge merge(&badpix);
|
|---|
| 187 | MImgCleanStd clean;
|
|---|
| 188 | MHillasCalc hcalc;
|
|---|
| 189 |
|
|---|
| 190 | // ------------------ Setup histograms and fill tasks ----------------
|
|---|
| 191 | MHCamEvent evt0a(0, "Cleaned", "Signal after Cleaning;;S [\\gamma]");
|
|---|
| 192 | MHCamEvent evt0b(0, "UsedPix", "Pixels marked Used;;Used [%]");
|
|---|
| 193 | evt0b.SetThreshold(0);
|
|---|
| 194 |
|
|---|
| 195 | //MHCamEventRot evt0r("UsedRot", "Pixels marked Used (derotated)");
|
|---|
| 196 | //evt0r.SetThreshold(0);
|
|---|
| 197 |
|
|---|
| 198 | MH3 h1("MEventRate.fRate");
|
|---|
| 199 | h1.SetName("MHEventRate");
|
|---|
| 200 | h1.SetTitle("Event Rate distribution;R [Hz];Counts");
|
|---|
| 201 | h1.SetLogy();
|
|---|
| 202 | /*
|
|---|
| 203 | MH3 h12("MEventRate10000.fRate");
|
|---|
| 204 | h12.SetName("MHEventRate");
|
|---|
| 205 | h12.SetLogy();
|
|---|
| 206 | */
|
|---|
| 207 | MBinning b1("BinningMHEventRate");
|
|---|
| 208 | b1.SetEdges(150, 0, 1500);
|
|---|
| 209 | plist.AddToList(&b1);
|
|---|
| 210 |
|
|---|
| 211 | MHVsTime hvs("MEventRate.fRate");
|
|---|
| 212 | hvs.SetTitle("Rate per 500 events;R [Hz];Counts");
|
|---|
| 213 | hvs.SetNumEvents(500);
|
|---|
| 214 |
|
|---|
| 215 | //MContinue cont1("MEventRate2.fRate/MEventRate.fRate>1.1");
|
|---|
| 216 | //MContinue cont2("MEventRate.fRate/MEventRate2.fRate>1.1");
|
|---|
| 217 |
|
|---|
| 218 | MFillH fillvs(&hvs, "MTime", "FillEventRate10s");
|
|---|
| 219 |
|
|---|
| 220 | MFillH fill0a(&evt0a, "MCerPhotEvt", "FillCerPhotEvt");
|
|---|
| 221 | MFillH fill0b(&evt0b, "MCerPhotEvt", "FillCntUsedPixels");
|
|---|
| 222 | //MFillH fill0r(&evt0r, "MCerPhotEvt", "FillCntUsedRotated");
|
|---|
| 223 | MFillH fill1("MHHillas", "MHillas", "FillHillas");
|
|---|
| 224 | MFillH fill2("MHHillasExt", "", "FillHillasExt");
|
|---|
| 225 | MFillH fill3("MHHillasSrc", "MHillasSrc", "FillHillasSrc");
|
|---|
| 226 | MFillH fill4("MHImagePar", "MImagePar", "FillImagePar");
|
|---|
| 227 | MFillH fill5("MHNewImagePar", "MNewImagePar", "FillNewImagePar");
|
|---|
| 228 | //MFillH fill6("MHImageParTime","MImageParTime","FillImageParTime");
|
|---|
| 229 | //MFillH fill7("MHNewImagePar2","MNewImagePar2","FillNewImagePar2");
|
|---|
| 230 | MFillH fill8(&h1, "", "FillEventRate");
|
|---|
| 231 | MFillH fill9("MHEffectiveOnTime", "MTime", "FillEffOnTime");
|
|---|
| 232 | //MFillH fillb(&h12, "", "FillEvtRate2");
|
|---|
| 233 | //MFillH fill9("MHCerPhot");
|
|---|
| 234 |
|
|---|
| 235 | //fill0r.SetDrawOption("colz");
|
|---|
| 236 | fill8.SetNameTab("EvtRate");
|
|---|
| 237 | fill9.SetNameTab("EffOnTime");
|
|---|
| 238 |
|
|---|
| 239 | // ------------------ Setup write task ----------------
|
|---|
| 240 |
|
|---|
| 241 | MWriteRootFile write(2, Form("%s{s/_Y_/_I_}", fPathOut.Data()), fOverwrite);
|
|---|
| 242 | // Data
|
|---|
| 243 | write.AddContainer("MHillas", "Events");
|
|---|
| 244 | write.AddContainer("MHillasExt", "Events");
|
|---|
| 245 | write.AddContainer("MHillasSrc", "Events");
|
|---|
| 246 | write.AddContainer("MImagePar", "Events");
|
|---|
| 247 | write.AddContainer("MNewImagePar", "Events");
|
|---|
| 248 | //write.AddContainer("MNewImagePar2", "Events");
|
|---|
| 249 | //write.AddContainer("MImageParTime", "Events");
|
|---|
| 250 | write.AddContainer("MRawEvtHeader", "Events");
|
|---|
| 251 | if (ismc)
|
|---|
| 252 | {
|
|---|
| 253 | write.AddContainer("MPointingPos", "Events");
|
|---|
| 254 | // Monte Carlo
|
|---|
| 255 | write.AddContainer("MMcEvt", "Events");
|
|---|
| 256 | write.AddContainer("MMcTrig", "Events");
|
|---|
| 257 | // Monte Carlo Headers
|
|---|
| 258 | write.AddContainer("MMcTrigHeader", "RunHeaders");
|
|---|
| 259 | write.AddContainer("MMcConfigRunHeader", "RunHeaders");
|
|---|
| 260 | write.AddContainer("MMcCorsikaRunHeader", "RunHeaders");
|
|---|
| 261 | }
|
|---|
| 262 | else
|
|---|
| 263 | {
|
|---|
| 264 | write.AddContainer("MTime", "Events");
|
|---|
| 265 | // Run Header
|
|---|
| 266 | write.AddContainer("MRawRunHeader", "RunHeaders");
|
|---|
| 267 | write.AddContainer("MBadPixelsCam", "RunHeaders");
|
|---|
| 268 | write.AddContainer("MGeomCam", "RunHeaders");
|
|---|
| 269 | //write.AddContainer("MObservatory", "RunHeaders");
|
|---|
| 270 | // Drive
|
|---|
| 271 | //write.AddContainer("MSrcPosCam", "Drive");
|
|---|
| 272 | write.AddContainer("MPointingPos", "Drive");
|
|---|
| 273 | write.AddContainer("MReportDrive", "Drive");
|
|---|
| 274 | write.AddContainer("MTimeDrive", "Drive");
|
|---|
| 275 | // Effective On Time
|
|---|
| 276 | write.AddContainer("MEffectiveOnTime", "EffectiveOnTime");
|
|---|
| 277 | write.AddContainer("MTimeEffectiveOnTime", "EffectiveOnTime");
|
|---|
| 278 | }
|
|---|
| 279 |
|
|---|
| 280 | MTaskList tlist2;
|
|---|
| 281 | tlist2.AddToList(&apply);
|
|---|
| 282 | if (!ismc)
|
|---|
| 283 | {
|
|---|
| 284 | tlist2.AddToList(&cont);
|
|---|
| 285 | tlist2.AddToList(&rate);
|
|---|
| 286 | //tlist2.AddToList(&rate1);
|
|---|
| 287 | tlist2.AddToList(&fillvs);
|
|---|
| 288 | //tlist2.AddToList(&cont1);
|
|---|
| 289 | //tlist2.AddToList(&cont2);
|
|---|
| 290 | tlist2.AddToList(&fill8);
|
|---|
| 291 | tlist2.AddToList(&fill9);
|
|---|
| 292 | }
|
|---|
| 293 | //tlist2.AddToList(&fillb);
|
|---|
| 294 | tlist2.AddToList(&clean);
|
|---|
| 295 | tlist2.AddToList(&fill0a);
|
|---|
| 296 | tlist2.AddToList(&fill0b);
|
|---|
| 297 | //tlist2.AddToList(&fill0r);
|
|---|
| 298 | tlist2.AddToList(&hcalc);
|
|---|
| 299 | tlist2.AddToList(&fill1);
|
|---|
| 300 | tlist2.AddToList(&fill2);
|
|---|
| 301 | tlist2.AddToList(&fill3);
|
|---|
| 302 | tlist2.AddToList(&fill4);
|
|---|
| 303 | tlist2.AddToList(&fill5);
|
|---|
| 304 | //tlist2.AddToList(&fill6);
|
|---|
| 305 | //tlist2.AddToList(&fill7);
|
|---|
| 306 | //tlist2.AddToList(&fill9);
|
|---|
| 307 |
|
|---|
| 308 | MPointingPosCalc pcalc;
|
|---|
| 309 | //MSrcPosFromModel srcpos;
|
|---|
| 310 |
|
|---|
| 311 | MTaskList tlist3;
|
|---|
| 312 | tlist3.AddToList(&pcalc);
|
|---|
| 313 | //tlist3.AddToList(&srcpos);
|
|---|
| 314 |
|
|---|
| 315 | tlist.AddToList(ismc ? (MTask*)&readmc : (MTask*)&readreal);
|
|---|
| 316 | tlist.AddToList(&tlist3, "Drive");
|
|---|
| 317 | tlist.AddToList(&tlist2, "Events");
|
|---|
| 318 | tlist.AddToList(&write);
|
|---|
| 319 |
|
|---|
| 320 | // Create and setup the eventloop
|
|---|
| 321 | MEvtLoop evtloop(fName);
|
|---|
| 322 | evtloop.SetParList(&plist);
|
|---|
| 323 | evtloop.SetDisplay(fDisplay);
|
|---|
| 324 | evtloop.SetLogStream(fLog);
|
|---|
| 325 | if (!SetupEnv(evtloop))
|
|---|
| 326 | return kFALSE;
|
|---|
| 327 |
|
|---|
| 328 | // Execute first analysis
|
|---|
| 329 | if (!evtloop.Eventloop(fMaxEvents))
|
|---|
| 330 | {
|
|---|
| 331 | *fLog << err << GetDescriptor() << ": Failed." << endl;
|
|---|
| 332 | return kFALSE;
|
|---|
| 333 | }
|
|---|
| 334 |
|
|---|
| 335 | tlist.PrintStatistics();
|
|---|
| 336 |
|
|---|
| 337 | if (!WriteResult())
|
|---|
| 338 | return kFALSE;
|
|---|
| 339 |
|
|---|
| 340 | *fLog << all << GetDescriptor() << ": Done." << endl;
|
|---|
| 341 | *fLog << endl << endl;
|
|---|
| 342 |
|
|---|
| 343 | return kTRUE;
|
|---|
| 344 | }
|
|---|