| 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, 12/2003 <mailto:tbretz@astro.uni-wuerzburg.de>
|
|---|
| 19 | !
|
|---|
| 20 | ! Copyright: MAGIC Software Development, 2000-2005
|
|---|
| 21 | !
|
|---|
| 22 | !
|
|---|
| 23 | \* ======================================================================== */
|
|---|
| 24 |
|
|---|
| 25 | /////////////////////////////////////////////////////////////////////////////
|
|---|
| 26 | //
|
|---|
| 27 | // MTFillMatrix
|
|---|
| 28 | //
|
|---|
| 29 | // Use this tool to fill eg trainings and test-matrices, while the matrix
|
|---|
| 30 | // to be filled can be a real matrix (MHMatrix) or a file (MWriteRootFile)
|
|---|
| 31 | // or both.
|
|---|
| 32 | //
|
|---|
| 33 | // First create a reference histogram (MH3). For more details see
|
|---|
| 34 | // MFEventSelector2 which is used to select events according to the
|
|---|
| 35 | // reference histogram.
|
|---|
| 36 | //
|
|---|
| 37 | // If no reference histogram is available the number of events are
|
|---|
| 38 | // randomly choosen from the sample with a probability which fits
|
|---|
| 39 | // the total destination number of events.
|
|---|
| 40 | //
|
|---|
| 41 | // Here is an example of how to choose 1000 events somehow distributed in
|
|---|
| 42 | // size from a file.
|
|---|
| 43 | // -----------------------------------------------------------------------
|
|---|
| 44 | //
|
|---|
| 45 | // MH3 ref("MHillas.fSize"); // choose a predefined size distribution
|
|---|
| 46 | // // Now setup the distribution
|
|---|
| 47 | //
|
|---|
| 48 | // MHMatrix matrix1; // create matrix to fill
|
|---|
| 49 | // matrix.AddColumn("MHillas.fWidth"); // setup columns of your matrix
|
|---|
| 50 | //
|
|---|
| 51 | // MReadMarsFile read("myfile.root"); // Setup your 'reader'
|
|---|
| 52 | // read.DisableAutoScheme(); // make sure everything is read
|
|---|
| 53 | //
|
|---|
| 54 | // MTFillMatrix fill(&ref); // setup MTFillMatrix
|
|---|
| 55 | // fill.SetNumDestEvents1(1000); // setup number of events to select
|
|---|
| 56 | // fill.SetDestMatrix1(&matrix1); // setup destination matrix
|
|---|
| 57 | // if (!fill.Process()) // check if everything worked
|
|---|
| 58 | // return;
|
|---|
| 59 | // fill.WriteMatrix1("myoutput.root"); // write matrix to file
|
|---|
| 60 | //
|
|---|
| 61 | //
|
|---|
| 62 | // To get two matrices instead of one (splitted randomly) you can add
|
|---|
| 63 | // the following before calling Process():
|
|---|
| 64 | // ------------------------------------------------------------------------
|
|---|
| 65 | //
|
|---|
| 66 | // MHMatrix matrix2;
|
|---|
| 67 | // fill.SetNumDestEvents2(500); // setup number of events to select
|
|---|
| 68 | // fill.SetDestMatrix2(&matrix2); // setup destination matrix
|
|---|
| 69 | // [...]
|
|---|
| 70 | // fill.WriteMatrix2("myoutput2.root");
|
|---|
| 71 | //
|
|---|
| 72 | //
|
|---|
| 73 | // To write both matrices into a single file use:
|
|---|
| 74 | // ----------------------------------------------
|
|---|
| 75 | //
|
|---|
| 76 | // fill.WriteMatrices("myfile.root");
|
|---|
| 77 | //
|
|---|
| 78 | //
|
|---|
| 79 | // Task List execution
|
|---|
| 80 | // ----------------------------------------------
|
|---|
| 81 | //
|
|---|
| 82 | // The tasklist is excuted in the follwowing order:
|
|---|
| 83 | // - fReader set by SetReader
|
|---|
| 84 | // - fPreTasks in the order as set with the SetPreTask functions
|
|---|
| 85 | // - fPreCuts in the order as set with the SetPreCut functions
|
|---|
| 86 | // - Selection the event selection
|
|---|
| 87 | // - fPostTasks in the order as set with the SetPostTask functions
|
|---|
| 88 | // - Fill Matrix
|
|---|
| 89 | // - Write output
|
|---|
| 90 | //
|
|---|
| 91 | /////////////////////////////////////////////////////////////////////////////
|
|---|
| 92 | #include "MTFillMatrix.h"
|
|---|
| 93 |
|
|---|
| 94 | #include <TFile.h>
|
|---|
| 95 | #include <TMath.h>
|
|---|
| 96 |
|
|---|
| 97 | // environment
|
|---|
| 98 | #include "MHMatrix.h"
|
|---|
| 99 |
|
|---|
| 100 | #include "MLog.h"
|
|---|
| 101 | #include "MLogManip.h"
|
|---|
| 102 |
|
|---|
| 103 | // eventloop
|
|---|
| 104 | #include "MParList.h"
|
|---|
| 105 | #include "MTaskList.h"
|
|---|
| 106 | #include "MEvtLoop.h"
|
|---|
| 107 |
|
|---|
| 108 | // tasks
|
|---|
| 109 | #include "MRead.h"
|
|---|
| 110 | #include "MFillH.h"
|
|---|
| 111 | #include "MContinue.h"
|
|---|
| 112 |
|
|---|
| 113 | // filters
|
|---|
| 114 | #include "MFDataPhrase.h"
|
|---|
| 115 | #include "MFilterList.h"
|
|---|
| 116 | #include "MFEventSelector.h"
|
|---|
| 117 | #include "MFEventSelector2.h"
|
|---|
| 118 |
|
|---|
| 119 | ClassImp(MTFillMatrix);
|
|---|
| 120 |
|
|---|
| 121 | using namespace std;
|
|---|
| 122 |
|
|---|
| 123 | // --------------------------------------------------------------------------
|
|---|
| 124 | //
|
|---|
| 125 | // Print Size and contained columns.
|
|---|
| 126 | // Check whether the number of generated events is compatible with
|
|---|
| 127 | // the number of requested events.
|
|---|
| 128 | //
|
|---|
| 129 | Bool_t MTFillMatrix::CheckResult(MHMatrix *m, Int_t num) const
|
|---|
| 130 | {
|
|---|
| 131 | if (!m || num==0)
|
|---|
| 132 | return kTRUE;
|
|---|
| 133 |
|
|---|
| 134 | m->Print("SizeCols");
|
|---|
| 135 |
|
|---|
| 136 | const Int_t generated = m->GetM().GetNrows();
|
|---|
| 137 | if (TMath::Abs(generated-num) <= TMath::Sqrt(9.0*num))
|
|---|
| 138 | return kTRUE;
|
|---|
| 139 |
|
|---|
| 140 | *fLog << warn << "WARNING - No. of generated events (";
|
|---|
| 141 | *fLog << generated << ") incompatible with requested events (";
|
|---|
| 142 | *fLog << num << ") for " << m->GetDescriptor() << endl;
|
|---|
| 143 |
|
|---|
| 144 | return kFALSE;
|
|---|
| 145 | }
|
|---|
| 146 |
|
|---|
| 147 | // --------------------------------------------------------------------------
|
|---|
| 148 | //
|
|---|
| 149 | // Write the given MHMatrix with its name as default name to a
|
|---|
| 150 | // file fname.
|
|---|
| 151 | //
|
|---|
| 152 | Bool_t MTFillMatrix::WriteMatrix(MHMatrix *m, const TString &fname, Int_t i) const
|
|---|
| 153 | {
|
|---|
| 154 | if (!m)
|
|---|
| 155 | {
|
|---|
| 156 | *fLog << "ERROR - Unable to write matrix #" << i << " (=NULL)... ignored." << endl;
|
|---|
| 157 | return kFALSE;
|
|---|
| 158 | }
|
|---|
| 159 | if (fname.IsNull())
|
|---|
| 160 | {
|
|---|
| 161 | *fLog << "ERROR - Unable to write matrix, file name empty." << endl;
|
|---|
| 162 | return kFALSE;
|
|---|
| 163 | }
|
|---|
| 164 |
|
|---|
| 165 | TFile file(fname, "RECREATE", m->GetTitle());
|
|---|
| 166 | m->Write();
|
|---|
| 167 | return kTRUE;
|
|---|
| 168 | }
|
|---|
| 169 |
|
|---|
| 170 | void MTFillMatrix::Init(const char *name, const char *title)
|
|---|
| 171 | {
|
|---|
| 172 | fName = name ? name : "MFillMatrix";
|
|---|
| 173 | fTitle = title ? title : "Tool to fill MHMatrix from file";
|
|---|
| 174 | }
|
|---|
| 175 |
|
|---|
| 176 | MTFillMatrix::MTFillMatrix(const char *name, const char *title)
|
|---|
| 177 | : fReference(0), fReader(0), fDestMatrix1(0), fDestMatrix2(0),
|
|---|
| 178 | fNumDestEvents1(0), fNumDestEvents2(0), fNumMaxEvents(0),
|
|---|
| 179 | fWriteFile1(0), fWriteFile2(0)
|
|---|
| 180 | {
|
|---|
| 181 | Init(name, title);
|
|---|
| 182 | }
|
|---|
| 183 |
|
|---|
| 184 | // --------------------------------------------------------------------------
|
|---|
| 185 | //
|
|---|
| 186 | // Constructor. Takes an MH3 as argument. This MH3 is the reference
|
|---|
| 187 | // distribution to fill the matrix. More information can be found
|
|---|
| 188 | // at MFEventSelector2 which is used to select the events.
|
|---|
| 189 | //
|
|---|
| 190 | // If no MH3 *ref is given the events are randomly selected from the
|
|---|
| 191 | // total sample - this may result in samples which don't have exactly
|
|---|
| 192 | // the predefined size, but it is much faster.
|
|---|
| 193 | //
|
|---|
| 194 | MTFillMatrix::MTFillMatrix(const MH3 *ref, const char *name, const char *title)
|
|---|
| 195 | : fReference(0), fReader(0), fDestMatrix1(0), fDestMatrix2(0),
|
|---|
| 196 | fNumDestEvents1(0), fNumDestEvents2(0), fNumMaxEvents(0),
|
|---|
| 197 | fWriteFile1(0), fWriteFile2(0)
|
|---|
| 198 | {
|
|---|
| 199 | Init(name, title);
|
|---|
| 200 | if (ref)
|
|---|
| 201 | fReference = (MH3*)ref->Clone();
|
|---|
| 202 | }
|
|---|
| 203 |
|
|---|
| 204 | MTFillMatrix::~MTFillMatrix()
|
|---|
| 205 | {
|
|---|
| 206 | if (fReference)
|
|---|
| 207 | delete fReference;
|
|---|
| 208 | }
|
|---|
| 209 |
|
|---|
| 210 | //------------------------------------------------------------------------
|
|---|
| 211 | //
|
|---|
| 212 | // Function serving AddPreCuts, AddPreTasks and AddPostTasks
|
|---|
| 213 | //
|
|---|
| 214 | void MTFillMatrix::Add(const TList &src, const TClass *cls, TList &dest)
|
|---|
| 215 | {
|
|---|
| 216 | TIter Next(&src);
|
|---|
| 217 | TObject *obj=0;
|
|---|
| 218 | while ((obj=Next()))
|
|---|
| 219 | if (obj->InheritsFrom(cls))
|
|---|
| 220 | dest.Add(obj);
|
|---|
| 221 | }
|
|---|
| 222 |
|
|---|
| 223 | //------------------------------------------------------------------------
|
|---|
| 224 | //
|
|---|
| 225 | // Add a cut which is used to fill the matrix, eg "MMcEvt.fOartId<1.5"
|
|---|
| 226 | // (The rule is applied, nit inverted: The matrix is filled with
|
|---|
| 227 | // the events fullfilling the condition)
|
|---|
| 228 | //
|
|---|
| 229 | void MTFillMatrix::AddPreCut(const char *rule)
|
|---|
| 230 | {
|
|---|
| 231 | MFilter *f = new MFDataPhrase(rule);
|
|---|
| 232 | f->SetBit(kCanDelete);
|
|---|
| 233 | AddPreCut(f);
|
|---|
| 234 | }
|
|---|
| 235 |
|
|---|
| 236 | //------------------------------------------------------------------------
|
|---|
| 237 | //
|
|---|
| 238 | // Add a cut which is used to fill the matrix. If kCanDelete is set
|
|---|
| 239 | // MJOptimize takes the ownership.
|
|---|
| 240 | //
|
|---|
| 241 | void MTFillMatrix::AddPreCut(MFilter *f)
|
|---|
| 242 | {
|
|---|
| 243 | fPreCuts.Add(f);
|
|---|
| 244 | }
|
|---|
| 245 |
|
|---|
| 246 | //------------------------------------------------------------------------
|
|---|
| 247 | //
|
|---|
| 248 | // Add all entries deriving from MFilter from list to PreCuts.
|
|---|
| 249 | // The ownership is not affected.
|
|---|
| 250 | //
|
|---|
| 251 | void MTFillMatrix::AddPreCuts(const TList &list)
|
|---|
| 252 | {
|
|---|
| 253 | Add(list, MFilter::Class(), fPreCuts);
|
|---|
| 254 | }
|
|---|
| 255 |
|
|---|
| 256 | //------------------------------------------------------------------------
|
|---|
| 257 | //
|
|---|
| 258 | // Add a task which is executed before the precuts. If kCanDelete is set
|
|---|
| 259 | // MJOptimize takes the ownership.
|
|---|
| 260 | //
|
|---|
| 261 | void MTFillMatrix::AddPreTask(MTask *t)
|
|---|
| 262 | {
|
|---|
| 263 | fPreTasks.Add(t);
|
|---|
| 264 | }
|
|---|
| 265 |
|
|---|
| 266 | //------------------------------------------------------------------------
|
|---|
| 267 | //
|
|---|
| 268 | // Add all entries deriving from MTask from list to PreTasks.
|
|---|
| 269 | // The ownership is not affected.
|
|---|
| 270 | //
|
|---|
| 271 | void MTFillMatrix::AddPreTasks(const TList &list)
|
|---|
| 272 | {
|
|---|
| 273 | Add(list, MTask::Class(), fPreTasks);
|
|---|
| 274 | }
|
|---|
| 275 |
|
|---|
| 276 | //------------------------------------------------------------------------
|
|---|
| 277 | //
|
|---|
| 278 | // Add a task which is executed after the precuts. If kCanDelete is set
|
|---|
| 279 | // MJOptimize takes the ownership.
|
|---|
| 280 | //
|
|---|
| 281 | void MTFillMatrix::AddPostTask(MTask *t)
|
|---|
| 282 | {
|
|---|
| 283 | fPostTasks.Add(t);
|
|---|
| 284 | }
|
|---|
| 285 |
|
|---|
| 286 | //------------------------------------------------------------------------
|
|---|
| 287 | //
|
|---|
| 288 | // Add all entries deriving from MTask from list to PostTasks.
|
|---|
| 289 | // The ownership is not affected.
|
|---|
| 290 | //
|
|---|
| 291 | void MTFillMatrix::AddPostTasks(const TList &list)
|
|---|
| 292 | {
|
|---|
| 293 | Add(list, MTask::Class(), fPostTasks);
|
|---|
| 294 | }
|
|---|
| 295 |
|
|---|
| 296 | // --------------------------------------------------------------------------
|
|---|
| 297 | //
|
|---|
| 298 | // Fill the matrix (FIXME: Flow diagram missing)
|
|---|
| 299 | //
|
|---|
| 300 | Bool_t MTFillMatrix::Process(const MParList &parlist)
|
|---|
| 301 | {
|
|---|
| 302 | if (!fReader)
|
|---|
| 303 | {
|
|---|
| 304 | *fLog << err << "ERROR - No task to read data was given... abort." << endl;
|
|---|
| 305 | return kFALSE;
|
|---|
| 306 | }
|
|---|
| 307 |
|
|---|
| 308 | *fLog << inf;
|
|---|
| 309 | fLog->Separator(GetDescriptor());
|
|---|
| 310 | *fLog << "Fill " << fDestMatrix1->GetDescriptor() << " with " << fNumDestEvents1 << endl;
|
|---|
| 311 | if (fDestMatrix2)
|
|---|
| 312 | *fLog << "Fill " << fDestMatrix2->GetDescriptor() << " with " << fNumDestEvents2 << endl;
|
|---|
| 313 | *fLog << "Distribution choosen ";
|
|---|
| 314 | if (fReference && fReference->GetHist().GetEntries()>0)
|
|---|
| 315 | *fLog << "from " << fReference->GetDescriptor();
|
|---|
| 316 | else
|
|---|
| 317 | *fLog << "randomly";
|
|---|
| 318 | *fLog << endl;
|
|---|
| 319 |
|
|---|
| 320 | //
|
|---|
| 321 | // Create parameter list and task list, add tasklist to parlist
|
|---|
| 322 | //
|
|---|
| 323 | parlist.Print();
|
|---|
| 324 | MParList plist(parlist);
|
|---|
| 325 | MTaskList tlist;
|
|---|
| 326 | plist.AddToList(&tlist);
|
|---|
| 327 |
|
|---|
| 328 | //
|
|---|
| 329 | // A selector to select a given number of events from a sample
|
|---|
| 330 | //
|
|---|
| 331 | // FIXME: Merge MFEventSelector and MFEventSelector2
|
|---|
| 332 | MFilter *selector=0;
|
|---|
| 333 | if (fNumDestEvents1>0 || fNumDestEvents2>0)
|
|---|
| 334 | {
|
|---|
| 335 | if (fReference)
|
|---|
| 336 | {
|
|---|
| 337 | // Case of a reference/nominal distribution
|
|---|
| 338 | // The events must be read before selection
|
|---|
| 339 | MFEventSelector2 *sel = new MFEventSelector2(*fReference);
|
|---|
| 340 | sel->SetNumMax(fNumDestEvents1+fNumDestEvents2);
|
|---|
| 341 | sel->SetInverted();
|
|---|
| 342 |
|
|---|
| 343 | selector = sel;
|
|---|
| 344 | }
|
|---|
| 345 | else
|
|---|
| 346 | {
|
|---|
| 347 | // Case of a random distribution
|
|---|
| 348 | // The events can be selected before reading
|
|---|
| 349 | MFEventSelector *sel = new MFEventSelector;
|
|---|
| 350 | sel->SetNumSelectEvts(fNumDestEvents1+fNumDestEvents2);
|
|---|
| 351 | fReader->SetSelector(sel);
|
|---|
| 352 |
|
|---|
| 353 | selector = sel;
|
|---|
| 354 | }
|
|---|
| 355 | }
|
|---|
| 356 |
|
|---|
| 357 | //
|
|---|
| 358 | // Continue for PreCuts
|
|---|
| 359 | //
|
|---|
| 360 | MFilterList list;
|
|---|
| 361 | list.SetName("PreCuts");
|
|---|
| 362 | if (!list.AddToList(fPreCuts))
|
|---|
| 363 | *fLog << err << "ERROR - Calling MFilterList::AddToList for fPreCuts failed!" << endl;
|
|---|
| 364 |
|
|---|
| 365 | MContinue cont0(&list);
|
|---|
| 366 | cont0.SetInverted();
|
|---|
| 367 |
|
|---|
| 368 | //
|
|---|
| 369 | // Continue for all events which are not (SetInverted())
|
|---|
| 370 | // selected by the 'selector'
|
|---|
| 371 | //
|
|---|
| 372 | MContinue cont(selector);
|
|---|
| 373 |
|
|---|
| 374 | //
|
|---|
| 375 | // Create a filter doing a random split
|
|---|
| 376 | //
|
|---|
| 377 | const Double_t prob = (Double_t)fNumDestEvents1/(fNumDestEvents1+fNumDestEvents2);
|
|---|
| 378 | MFEventSelector split;
|
|---|
| 379 | split.SetSelectionRatio(prob);
|
|---|
| 380 |
|
|---|
| 381 | //
|
|---|
| 382 | // Create the logical inverted filter for 'split'
|
|---|
| 383 | //
|
|---|
| 384 | MFilterList invsplit;
|
|---|
| 385 | invsplit.AddToList(&split);
|
|---|
| 386 | invsplit.SetInverted();
|
|---|
| 387 |
|
|---|
| 388 | //
|
|---|
| 389 | // The two tasks filling the two matrices
|
|---|
| 390 | //
|
|---|
| 391 | MFillH fill1(fDestMatrix1);
|
|---|
| 392 | MFillH fill2(fDestMatrix2);
|
|---|
| 393 | fill1.SetFilter(&split);
|
|---|
| 394 | fill2.SetFilter(&invsplit);
|
|---|
| 395 |
|
|---|
| 396 | // entries in MTaskList
|
|---|
| 397 | tlist.AddToList(fReader); // Read events
|
|---|
| 398 | tlist.AddToList(fPreTasks, 0); // PreTasks
|
|---|
| 399 | if (fPreCuts.GetEntries()>0)
|
|---|
| 400 | tlist.AddToList(&cont0); // PreCuts
|
|---|
| 401 | if (fReference && selector)
|
|---|
| 402 | tlist.AddToList(&cont); // select a sample of events
|
|---|
| 403 | tlist.AddToList(&invsplit); // process invsplit (which implicitly processes split)
|
|---|
| 404 | tlist.AddToList(fPostTasks, 0); // PostTasks
|
|---|
| 405 | if (fDestMatrix1)
|
|---|
| 406 | tlist.AddToList(&fill1); // fill matrix 1
|
|---|
| 407 | if (fDestMatrix2)
|
|---|
| 408 | tlist.AddToList(&fill2); // fill matrix 2
|
|---|
| 409 | if (fWriteFile1)
|
|---|
| 410 | {
|
|---|
| 411 | fWriteFile1->SetFilter(&split);
|
|---|
| 412 | tlist.AddToList(fWriteFile1);
|
|---|
| 413 | }
|
|---|
| 414 | if (fWriteFile2)
|
|---|
| 415 | {
|
|---|
| 416 | fWriteFile2->SetFilter(&invsplit);
|
|---|
| 417 | tlist.AddToList(fWriteFile2);
|
|---|
| 418 | }
|
|---|
| 419 |
|
|---|
| 420 | tlist.SetAccelerator(MTask::kAccDontReset|MTask::kAccDontTime);
|
|---|
| 421 |
|
|---|
| 422 | //
|
|---|
| 423 | // Execute the eventloop
|
|---|
| 424 | //
|
|---|
| 425 | MEvtLoop evtloop(fName);
|
|---|
| 426 | evtloop.SetParList(&plist);
|
|---|
| 427 | evtloop.SetDisplay(fDisplay);
|
|---|
| 428 | evtloop.SetLogStream(fLog);
|
|---|
| 429 |
|
|---|
| 430 | const Bool_t rc = evtloop.Eventloop(fNumMaxEvents);
|
|---|
| 431 |
|
|---|
| 432 | if (selector)
|
|---|
| 433 | delete selector;
|
|---|
| 434 |
|
|---|
| 435 | if (!rc)
|
|---|
| 436 | {
|
|---|
| 437 | *fLog << err << GetDescriptor() << ": Failed." << endl;
|
|---|
| 438 | return kFALSE;
|
|---|
| 439 | }
|
|---|
| 440 |
|
|---|
| 441 | // Check the result of filling...
|
|---|
| 442 | CheckResult(fDestMatrix1, fNumDestEvents1);
|
|---|
| 443 | CheckResult(fDestMatrix2, fNumDestEvents2);
|
|---|
| 444 |
|
|---|
| 445 | *fLog << inf << GetDescriptor() << ": Done." << endl;
|
|---|
| 446 | return kTRUE;
|
|---|
| 447 | }
|
|---|
| 448 |
|
|---|
| 449 | // --------------------------------------------------------------------------
|
|---|
| 450 | //
|
|---|
| 451 | // Write both matrices to a file. Return kFALSE if writing one of them
|
|---|
| 452 | // failed or both have the same name.
|
|---|
| 453 | //
|
|---|
| 454 | Bool_t MTFillMatrix::WriteMatrices(const TString &fname) const
|
|---|
| 455 | {
|
|---|
| 456 | if (fDestMatrix1 && fDestMatrix2 &&
|
|---|
| 457 | (TString)fDestMatrix1->GetName()==(TString)fDestMatrix2->GetName())
|
|---|
| 458 | {
|
|---|
| 459 | *fLog << "ERROR - Cannot write matrices (both have the same name)... ignored." << endl;
|
|---|
| 460 | return kFALSE;
|
|---|
| 461 | }
|
|---|
| 462 |
|
|---|
| 463 | return WriteMatrix1(fname) && WriteMatrix2(fname);
|
|---|
| 464 | }
|
|---|
| 465 |
|
|---|
| 466 | Int_t MTFillMatrix::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
|
|---|
| 467 | {
|
|---|
| 468 | Bool_t rc = kFALSE;
|
|---|
| 469 | if (IsEnvDefined(env, prefix, "NumDestEvents1", print))
|
|---|
| 470 | {
|
|---|
| 471 | rc = kTRUE;
|
|---|
| 472 | SetNumDestEvents1(GetEnvValue(env, prefix, "NumDestEvents1", fNumDestEvents1));
|
|---|
| 473 | }
|
|---|
| 474 | if (IsEnvDefined(env, prefix, "NumDestEvents2", print))
|
|---|
| 475 | {
|
|---|
| 476 | rc = kTRUE;
|
|---|
| 477 | SetNumDestEvents2(GetEnvValue(env, prefix, "NumDestEvents2", fNumDestEvents2));
|
|---|
| 478 | }
|
|---|
| 479 | return rc;
|
|---|
| 480 | }
|
|---|