source: trunk/MagicSoft/Mars/mjobs/MJOptimize.cc@ 6949

Last change on this file since 6949 was 6949, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 28.4 KB
Line 
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, 9/2004 <mailto:tbretz@astro.uni-wuerzburg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2004
21!
22!
23\* ======================================================================== */
24
25/////////////////////////////////////////////////////////////////////////////
26//
27// MJOptimize
28//
29// Class for otimizing the parameters of the supercuts
30//
31// Minimization Control
32// ====================
33//
34// To choose the minimization algorithm use:
35// void SetOptimizer(Optimizer_t o);
36//
37// Allowed options are:
38// enum Optimizer_t
39// {
40// kMigrad, // Minimize by the method of Migrad
41// kSimplex, // Minimize by the method of Simplex
42// kMinimize, // Migrad + Simplex (if Migrad fails)
43// kMinos, // Minos error determination
44// kImprove, // Local minimum search
45// kSeek, // Minimize by the method of Monte Carlo
46// kNone // Skip optimization
47// };
48//
49// For more details on the methods see TMinuit.
50//
51//
52// You can change the behaviour of the minimization using
53//
54// void SetNumMaxCalls(UInt_t num=0);
55// void SetTolerance(Float_t tol=0);
56//
57// While NumMaxCalls is the first, Tolerance the second arguement.
58// For more details start root and type
59//
60// gMinuit->mnhelp("command")
61//
62// while command can be
63// * MIGRAD
64// * SIMPLEX
65// * MINIMIZE
66// * MINOS
67// * IMPROVE
68// * SEEK
69//
70// The default (num==0 and tol==0) should always give you the
71// corresponding defaults used in Minuit.
72//
73//
74// FIXME: Implement changing cut in hadronness...
75// FIXME: Show MHSignificance on MStatusDisplay during filling...
76// FIXME: Choose step-size percentage as static data membewr
77// FIXME: Choose minimization method
78//
79/////////////////////////////////////////////////////////////////////////////
80#include "MJOptimize.h"
81
82#include <TMinuit.h>
83#include <TVirtualFitter.h>
84
85#include <TStopwatch.h>
86
87#include <TCanvas.h>
88
89#include <TGraph.h>
90#include <TMultiGraph.h>
91
92#include "MParList.h"
93#include "MTaskList.h"
94#include "MGeomCamCT1.h"
95#include "MFEventSelector.h"
96#include "MFDataMember.h"
97#include "MReadTree.h"
98#include "MHMatrix.h"
99#include "MEnergyEstimate.h"
100#include "MMatrixLoop.h"
101#include "MChisqEval.h"
102#include "MEvtLoop.h"
103#include "MDataElement.h"
104#include "MDataMember.h"
105#include "MLog.h"
106#include "MLogManip.h"
107#include "MParameters.h"
108#include "MFillH.h"
109#include "MF.h"
110#include "MFilterList.h"
111#include "../mfilter/MFMagicCuts.h"
112#include "MContinue.h"
113#include "MGeomCamMagic.h"
114#include "MHillasSrcCalc.h"
115#include "MHMatrix.h"
116#include "MMatrixLoop.h"
117#include "../mhflux/MHAlpha.h"
118#include "MStatusDisplay.h"
119#include "../mhflux/MHEnergyEst.h"
120#include "MDirIter.h"
121#include "MSequence.h"
122
123using namespace std;
124
125//------------------------------------------------------------------------
126//
127// fcn calculates the function to be minimized (using TMinuit::Migrad)
128//
129void MJOptimize::fcn(Int_t &npar, Double_t *gin, Double_t &f, Double_t *par, Int_t iflag)
130{
131 MJOptimize *optim = (MJOptimize*)gMinuit->GetObjectFit();
132
133 // WORKAROUND --- FOR WHAT?
134 if (gMinuit->fEpsi<1e-2)
135 {
136 *optim->fLog << warn << "WARNING - For unknown reasons: fEspi<1e-100... resetting to 0.01." << endl;
137 gMinuit->fEpsi = 0.01;
138 }
139
140 TMinuit *minuit = gMinuit;
141 f = optim->Fcn(TArrayD(TMath::Min(gMinuit->fMaxpar, optim->fParameters.GetSize()), par), minuit);
142 gMinuit = minuit;
143
144}
145
146Double_t MJOptimize::Fcn(const TArrayD &par, TMinuit *minuit)
147{
148 if (fEvtLoop->GetDisplay()!=fDisplay)
149 return 0;
150 /*
151 switch(iflag)
152 {
153 case 1: // first call
154 case 2: // calc derivative
155 break;
156 case 3:
157 // last call
158 MStatusDisplay *d = new MStatusDisplay;
159 fEvtLoop->SetDisplay(d);
160 break;
161 }
162 */
163 MParList *plist = fEvtLoop->GetParList();
164
165 MParameterD *eval = (MParameterD*)plist->FindObject("MinimizationValue", "MParameterD");
166 MParContainer *pars = (MParContainer*)plist->FindObject("MParameters", "MParContainer");
167
168 MRead *read = (MRead*)plist->FindObject("MTaskList")->FindObject("MRead");
169 if (read)
170 read->Rewind();
171
172 if (fDebug>=0)
173 {
174 *fLog << inf << "New Set: ";
175 for (Int_t i=0; i<fParameters.GetSize(); i++)
176 *fLog << par[i] << " ";
177 *fLog << endl;
178 }
179
180 pars->SetVariables(par);
181 eval->SetVal(0);
182
183 if (fDebug<3)
184 gLog.SetNullOutput(kTRUE);
185
186 TStopwatch clock;
187 clock.Start();
188 fEvtLoop->Eventloop(fNumEvents);
189 clock.Stop();
190
191 if (fDebug<3)
192 gLog.SetNullOutput(kFALSE);
193
194 const Double_t f = eval->GetVal();
195
196 if (fDebug>=0)
197 *fLog << inf << "Result F=" << f << endl;
198
199 if (fDebug>=1 && minuit)
200 {
201 Double_t fmin, fedm, errdef;
202 Int_t n1, n2, istat;
203 minuit->mnstat(fmin, fedm, errdef, n1, n2, istat);
204 *fLog << inf << underline << "Minimization Status so far:" << endl;
205 *fLog << " Calls: " << minuit->fNfcn << " (max=" << gMinuit->fMaxIterations << ")" << endl;
206 *fLog << " Parameters: fixed=" << gMinuit->fNpfix << ", free=" << gMinuit->fNpar << endl;
207 *fLog << " Func min: " << fmin << " (Epsi=" << gMinuit->fEpsi << ", Apsi=" << gMinuit->fApsi << ")" << endl;
208 *fLog << " Found edm: " << fedm << endl;
209 *fLog << " ErrDef: " << errdef << endl;
210 *fLog << " Status: ";
211
212 switch (istat)
213 {
214 case 0: *fLog << "n/a" << endl; break;
215 case 1: *fLog << "approximation only, not accurate" << endl; break;
216 case 2: *fLog << "full matrix, but forced positive-definite" << endl; break;
217 case 3: *fLog << "full accurate covariance matrix" << endl; break;
218 default: *fLog << "undefined" << endl; break;
219 }
220 }
221
222 if (fDebug>=1)
223 {
224 clock.Print();
225 fEvtLoop->GetTaskList()->PrintStatistics();
226 }
227
228 return f;
229}
230
231MJOptimize::MJOptimize() : fDebug(-1), fNumEvents(0), fType(kSimplex), fNumMaxCalls(0), fTolerance(0), fTestTrain(0)
232{
233 fRules.SetOwner();
234 fFilter.SetOwner();
235
236 fNamesOn.SetOwner();
237 fNamesOff.SetOwner();
238}
239
240//------------------------------------------------------------------------
241//
242// Add seqeunces from list to reader
243//
244Bool_t MJOptimize::AddSequences(MRead &read, TList &list) const
245{
246 MDirIter files;
247
248 TIter Next(&list);
249 TObject *o=0;
250 while ((o=Next()))
251 {
252 MSequence seq(o->GetName());
253 if (!seq.IsValid())
254 return kFALSE;
255
256 seq.SetupDatRuns(files, o->GetTitle(), "I");
257 }
258
259 return read.AddFiles(files)>0;
260}
261
262//------------------------------------------------------------------------
263//
264// Add on-sequences:
265// - fname: sequence file name (with path)
266// - dir: directory were image files are stored
267//
268void MJOptimize::AddSequenceOn(const char *fname, const char *dir)
269{
270 fNamesOn.Add(new TNamed(fname, dir));
271}
272
273//------------------------------------------------------------------------
274//
275// Add off-sequences:
276// - fname: sequence file name (with path)
277// - dir: directory were image files are stored
278//
279void MJOptimize::AddSequenceOff(const char *fname, const char *dir)
280{
281 fNamesOff.Add(new TNamed(fname, dir));
282}
283
284//------------------------------------------------------------------------
285//
286// Empty list of on- and off-sequences
287//
288void MJOptimize::ResetSequences()
289{
290 fNamesOn.Delete();
291 fNamesOff.Delete();
292}
293
294//------------------------------------------------------------------------
295//
296// Add a parameter used in your filters (see AddFilter) The parameter
297// index is returned,
298//
299// Int_t idx = AddParameter("log10(MHillas.fSize)");
300//
301// The indices area starting with 0 always.
302//
303Int_t MJOptimize::AddParameter(const char *rule)
304{
305 fRules.Add(new TNamed(rule, ""));
306 return fRules.GetSize()-1;
307}
308
309//------------------------------------------------------------------------
310//
311// Add a filter which can be applied in the optimization (for deatils
312// see correspodning Run function) You can use the indices you got by
313// AddParameter, eg
314//
315// AddFilter("M[0] < 3.2");
316//
317// if used in optimization you can do
318//
319// AddFilter("M[0] < [0]");
320//
321// for more details, see SetParameter and FixParameter
322//
323void MJOptimize::AddFilter(const char *rule)
324{
325 fFilter.Add(new MF(rule));
326}
327
328//------------------------------------------------------------------------
329//
330// Add a cut which is used to fill the matrix, eg "MMcEvt.fOartId<1.5"
331// (The rule is applied, nit inverted: The matrix is filled with
332// the events fullfilling the condition)
333//
334void MJOptimize::AddPreCut(const char *rule)
335{
336 MFilter *f = new MF(rule);
337 f->SetBit(kCanDelete);
338 AddPreCut(f);
339}
340
341//------------------------------------------------------------------------
342//
343// Add a cut which is used to fill the matrix. If kCanDelete is set
344// MJOptimize takes the ownership.
345//
346void MJOptimize::AddPreCut(MFilter *f)
347{
348 fPreCuts.Add(f);
349}
350
351//------------------------------------------------------------------------
352//
353// Set the fParameters Array accoring to par.
354//
355void MJOptimize::SetParameters(const TArrayD &par)
356{
357 fParameters = par;
358}
359
360//------------------------------------------------------------------------
361//
362// Set the number of events processed by the eventloop. (Be carfull,
363// if you are doing on-off analysis and you only process the first
364// 1000 events which are on-events only the optimization may not work)
365//
366void MJOptimize::SetNumEvents(UInt_t n)
367{
368 fNumEvents=n;
369}
370
371//------------------------------------------------------------------------
372//
373// Set a debug level, which tells the optimization how much information
374// is displayed about and in the running eventloop.
375//
376void MJOptimize::SetDebug(UInt_t n)
377{
378 fDebug=n;
379}
380
381//------------------------------------------------------------------------
382//
383// Set a optimization algorithm to be used. For more information see
384// TMinuit.
385//
386// Available Algorithms are:
387// kMigrad, // Minimize by the method of Migrad
388// kSimplex, // Minimize by the method of Simplex
389// kSeek // Minimize by the method of Monte Carlo
390//
391void MJOptimize::SetOptimizer(Optimizer_t o)
392{
393 fType = o;
394}
395
396//------------------------------------------------------------------------
397//
398// If a status didplay is set, search for tab "Optimizer".
399// If not found, create it.
400// In the tab search for TMultiGraph "Parameters".
401// If not found create it.
402// If empty create TGraphs.
403// Check number of graphs vs. number of parameters.
404// return TList with graphs.
405//
406TList *MJOptimize::GetPlots() const
407{
408 if (!fDisplay)
409 return NULL;
410
411 TCanvas *c=fDisplay->GetCanvas("Optimizer");
412 if (!c)
413 c = &fDisplay->AddTab("Optimizer");
414
415 TMultiGraph *mg = dynamic_cast<TMultiGraph*>(c->FindObject("Parameters"));
416 if (!mg)
417 mg = new TMultiGraph("Parameters", "Parameters of optimization");
418
419 TList *l = mg->GetListOfGraphs();
420 if (!l)
421 {
422 const Int_t n = fParameters.GetSize();
423 for (int i=0; i<n+1; i++)
424 {
425 TGraph *g = new TGraph;
426 if (i==n)
427 g->SetLineColor(kBlue);
428 mg->Add(g, "");
429 AddPoint(mg->GetListOfGraphs(), i, i==n?1:fParameters[i]);
430 }
431 mg->SetBit(kCanDelete);
432 mg->Draw("al*");
433
434 l = mg->GetListOfGraphs();
435 }
436
437 return l->GetSize() == fParameters.GetSize()+1 ? l : NULL;
438}
439
440//------------------------------------------------------------------------
441//
442// Add a point with y=val as last point in idx-th Tgraph of list l.
443//
444void MJOptimize::AddPoint(TList *l, Int_t idx, Float_t val) const
445{
446 if (!l)
447 return;
448
449 TGraph *gr = (TGraph*)l->At(idx);
450 gr->SetPoint(gr->GetN(), gr->GetN(), val);
451}
452
453Int_t MJOptimize::Minuit(TMinuit &minuit, const char *cmd) const
454{
455 Int_t err;
456 Double_t tmp[2] = { fNumMaxCalls, fTolerance };
457 minuit.mnexcm(cmd, tmp, 2, err);
458
459 switch (err)
460 {
461 case 0:
462 *fLog << inf << GetDescriptor() << " TMinuit::mnexcm excuted normally." << endl;
463 break;
464 case 1:
465 *fLog << warn << GetDescriptor() << " TMinuit::mnexcm command is blank... ignored." << endl;
466 break;
467 case 2:
468 *fLog << warn << GetDescriptor() << " TMinuit::mnexcm command-line syntax error... ignored." << endl;
469 break;
470 case 3:
471 *fLog << warn << GetDescriptor() << " TMinuit::mnexcm unknown command... ignored." << endl;
472 break;
473 case 4:
474 *fLog << warn << GetDescriptor() << " TMinuit::mnexcm - Abnormal termination (eg Migrad not converged)" << endl;
475 break;
476 /*
477 case 5:
478 *fLog << inf << GetDescriptor() << " TMinuit::mnexcm - Parameters requested." << endl;
479 break;
480 case 6:
481 *fLog << inf << GetDescriptor() << " TMinuit::mnexcm - SET INPUT returned." << endl;
482 break;
483 case 7:
484 *fLog << inf << GetDescriptor() << " TMinuit::mnexcm - SET TITLE returned." << endl;
485 break;
486 case 8:
487 *fLog << inf << GetDescriptor() << " TMinuit::mnexcm - SET COVAR returned." << endl;
488 break;
489 case 9:
490 *fLog << inf << GetDescriptor() << " TMinuit::mnexcm - reserved." << endl;
491 break;
492 case 10:
493 *fLog << inf << GetDescriptor() << " TMinuit::mnexcm - END returned." << endl;
494 break;
495 case 11:
496 *fLog << inf << GetDescriptor() << " TMinuit::mnexcm - EXIT or STOP returned." << endl;
497 break;
498 case 12:
499 *fLog << inf << GetDescriptor() << " TMinuit::mnexcm - RETURN returned." << endl;
500 break;*/
501 }
502
503 return err;
504}
505
506Bool_t MJOptimize::Optimize(MEvtLoop &evtloop)
507{
508 if (fType==kNone)
509 return kTRUE;
510
511 gMinuit = new TMinuit(fParameters.GetSize());
512
513 gMinuit->SetFCN(fcn);
514 gMinuit->SetObjectFit(this);
515 gMinuit->SetPrintLevel(-1); // Don't print when DefineParameter
516
517 //
518 // Set starting values and step sizes for parameters
519 //
520 for (Int_t i=0; i<fParameters.GetSize(); i++)
521 {
522 TString name = "par[";
523 name += i;
524 name += "]";
525 Double_t vinit = fParameters[i];
526 Double_t step = fStep[i];
527
528 Double_t limlo = fLimLo[i];
529 Double_t limup = fLimUp[i];
530
531 Bool_t rc = gMinuit->DefineParameter(i, name, vinit, step, limlo, limup);
532 if (rc)
533 {
534 *fLog << err << dbginf << "Error in defining parameter #" << i << endl;
535 return kFALSE;
536 }
537
538 if (step==0)
539 gMinuit->FixParameter(i);
540 }
541
542 gMinuit->SetPrintLevel(1); // Switch on pritning again
543 gMinuit->mnprin(1,0); // Print all parameters
544
545 fEvtLoop = &evtloop;
546
547 TList *g=GetPlots();
548
549 // Now ready for minimization step:
550 TStopwatch clock;
551 clock.Start();
552 switch (fType)
553 {
554 case kSimplex:
555 Simplex(*gMinuit);
556 break;
557 case kMigrad:
558 Migrad(*gMinuit);
559 break;
560 case kMinimize:
561 Minimize(*gMinuit);
562 break;
563 case kMinos:
564 Minos(*gMinuit);
565 break;
566 case kImprove:
567 Improve(*gMinuit);
568 break;
569 case kSeek:
570 Seek(*gMinuit);
571 break;
572 case kNone: // Should never happen
573 return kFALSE;
574 }
575 clock.Stop();
576 clock.Print();
577
578 if (evtloop.GetDisplay()!=fDisplay)
579 {
580 *fLog << inf << "Optimization aborted by user." << endl;
581 fDisplay = 0;
582 return kFALSE;
583 }
584
585 *fLog << inf << "Resulting Chisq: " << gMinuit->fAmin << endl;
586
587 //
588 // Update values of fA, fB:
589 //
590 for (Int_t i=0; i<fParameters.GetSize(); i++)
591 {
592 Double_t x1, x2;
593 gMinuit->GetParameter(i,x1,x2);
594 fParameters[i] = x1;
595 cout << i << ": " << fParameters[i] << endl;
596
597 AddPoint(g, i, x1);
598 }
599 AddPoint(g, fParameters.GetSize(), gMinuit->fAmin);
600
601 delete gMinuit;
602
603 return kTRUE;
604}
605
606//------------------------------------------------------------------------
607//
608// Optimize allows to use the optimizing by an eventloop based on
609// some requirements.
610//
611// 1) The tasklist to be executed must have the name MTaskList and
612// be an entry in the parameterlist.
613//
614// 2) The reading task (MReadMarsFile, MMatrixLoop) must have the name
615// "MRead". If it derives from MRead Rewind() must be implemented,
616// otherwise it must start reading from scratch each time its
617// PreProcess is called.
618//
619// 3) The parameters to be optimized must be accessed through (currently)
620// a single parameter container (MParContainer) called MParameters.
621// The parameters are set through SetVariables.
622//
623// 4) The result of a single function call for minimization (eg. chisquare)
624// must be available after the eventloop in a container of type
625// MParameterD with the name "MinimizationResult".
626//
627// 5) The parameters to start with must have been set using
628// MJOptimize::SetParameter or MJOptimize::SetParameters and
629// MJOptimize::FixParameter
630//
631// The behaviour of the optimization can be changed using:
632// void SetNumEvents(UInt_t n);
633// void SetDebug(UInt_t n);
634// void SetOptimizer(Optimizer_t o);
635//
636// After optimization the resulting parameters are set and another eventloop
637// with a MStatusDisplay is set is called. The resulting parameters can be
638// accessed using: GetParameters()
639//
640// To be fixed:
641// - MStatusDisplay should show status while optimization is running
642// - write result into MStatusDisplay
643// - save result
644//
645Bool_t MJOptimize::Optimize(MParList &parlist)
646{
647 // Checks to make sure, that fcn doesn't crash
648 if (!parlist.FindCreateObj("MParameterD", "MinimizationValue"))
649 return kFALSE;
650
651 if (!parlist.FindObject("MParameters", "MParContainer"))
652 {
653 *fLog << err << "MParameters [MParContainer] not found... abort." << endl;
654 return kFALSE;
655 }
656
657 MMatrixLoop *loop = dynamic_cast<MMatrixLoop*>(parlist.FindTask("MRead"));
658
659 TString txt("Starting ");
660 switch (fType)
661 {
662 case kMigrad: txt += "Migrad"; break;
663 case kMinimize: txt += "Minimize"; break;
664 case kMinos: txt += "Minos"; break;
665 case kImprove: txt += "Improve"; break;
666 case kSimplex: txt += "Simplex"; break;
667 case kSeek: txt += "Seek"; break;
668 case kNone: txt += "no"; break;
669 }
670 txt += " optimization";
671
672 fLog->Separator(txt);
673
674 // Setup eventloop
675 MEvtLoop evtloop;
676 evtloop.SetParList(&parlist);
677 evtloop.SetDisplay(fDisplay); // set display for evtloop and all childs
678 parlist.SetDisplay(0); // reset display for all contents of parlist and tasklist
679 evtloop.SetPrivateDisplay(); // prevent display from being cascaded again in PreProcess
680
681 *fLog << inf << "Number of Parameters: " << fParameters.GetSize() << endl;
682
683 // In case the reader is the matrix loop and testrain is enabled
684 // switch on even mode...
685 if (loop && TMath::Abs(fTestTrain)>0)
686 loop->SetOperationMode(fTestTrain>0?MMatrixLoop::kEven:MMatrixLoop::kOdd);
687
688 if (!Optimize(evtloop))
689 return kFALSE;
690
691 gMinuit = 0;
692
693 fEvtLoop->SetDisplay(fDisplay);
694 if (!Fcn(fParameters))
695 return kFALSE;
696
697 // In case the reader is the matrix loop and testrain is enabled
698 // switch on odd mode...
699 if (!loop || fTestTrain==0)
700 return kTRUE;
701
702 loop->SetOperationMode(fTestTrain<0?MMatrixLoop::kEven:MMatrixLoop::kOdd);
703
704 // Done already in Fcn
705 // list.SetVariables(fParameters);
706 return Fcn(fParameters);
707}
708
709void MJOptimize::AddRulesToMatrix(MHMatrix &m) const
710{
711 TIter Next1(&fRules);
712 TObject *o1=0;
713 while ((o1=Next1()))
714 m.AddColumn(o1->GetName());
715}
716
717//------------------------------------------------------------------------
718//
719// Fill matrix m by using read. Use rules as a filter if userules.
720//
721Bool_t MJOptimize::FillMatrix(MReadTree &read, MParList &parlist, Bool_t userules)
722{
723 MHMatrix *m = (MHMatrix*)parlist.FindObject("M", "MHMatrix");
724 if (!m)
725 {
726 *fLog << err << "MJOptimize::FillMatrix - ERROR: M [MHMatrix] not found in parlist... abort." << endl;
727 return kFALSE;
728 }
729
730 m->Print("cols");
731
732 //MParList parlist;
733
734 // MGeomCamMagic cam;
735 // parlist.AddToList(&cam);
736
737 MTaskList tlist;
738 parlist.Replace(&tlist);
739
740 MFillH fillh(m);
741
742 tlist.AddToList(&read);
743
744 MFilterList list;
745 if (!list.AddToList(fPreCuts))
746 *fLog << err << "ERROR - Calling MFilterList::AddToList for fPreCuts failed!" << endl;
747 if (userules)
748 SetupFilters(list);
749 list.SetName("PreCuts"); // reset Name set by SetupFilters
750 list.SetInverted(kFALSE); // reset inversion set by SetupFilters
751 fillh.SetFilter(&list);
752 tlist.AddToList(&list);
753
754 tlist.AddToList(&fillh);
755
756 MEvtLoop fillloop;
757 fillloop.SetParList(&parlist);
758 fillloop.SetDisplay(fDisplay);
759 if (!fillloop.Eventloop(fNumEvents))
760 {
761 *fLog << err << "Filling matrix failed..." << endl;
762 return kFALSE;
763 }
764 tlist.PrintStatistics();
765
766 *fLog << inf << "Read events from file '" << read.GetFileName() << "'" << endl;
767
768 if (fillloop.GetDisplay()!=fDisplay)
769 {
770 fDisplay = 0;
771 *fLog << inf << "Optimization aborted by user." << endl;
772 return kFALSE;
773 }
774
775 m->Print("size");
776
777 return kTRUE;
778}
779
780//------------------------------------------------------------------------
781//
782// Adds all filters to MFilterList
783//
784void MJOptimize::SetupFilters(MFilterList &list, MFilter *filter) const
785{
786 list.SetName("MParameters");
787 list.SetInverted();
788
789 if (filter)
790 {
791 if (fFilter.GetSize()>0)
792 {
793 *fLog << inf;
794 *fLog << "INFORMATION - You are using an external filter and internal filters." << endl;
795 *fLog << " Please make sure that all parameters '[i]' are starting" << endl;
796 *fLog << " behind the number of parameters of the external filter." << endl;
797 }
798 list.AddToList(filter);
799 }
800
801 if (!list.AddToList(fFilter))
802 *fLog << err << "ERROR - Calling MFilterList::AddToList fFilter failed!" << endl;
803
804 *fLog << inf << "Filter: ";
805 list.Print();
806 *fLog << endl;
807}
808
809//------------------------------------------------------------------------
810//
811Bool_t MJOptimize::Run(const char *fname, MFilter *filter, MAlphaFitter *fit)
812{
813 if (fParameters.GetSize()==0)
814 {
815 *fLog << err << "Sorry, no parameters defined." << endl;
816 return kFALSE;
817 }
818
819 fLog->Separator("Preparing On-only-optimization");
820
821 MParList parlist;
822
823 MGeomCamMagic geom; // For GetConvMm2Deg
824 parlist.AddToList(&geom);
825
826 MHMatrix m("M");
827 AddRulesToMatrix(m);
828 parlist.AddToList(&m);
829
830 MHAlpha hist;
831 hist.SkipHistTime();
832 hist.SkipHistTheta();
833 hist.SkipHistEnergy();
834 hist.InitMapping(&m);
835
836 if (filter && filter->InheritsFrom(MFMagicCuts::Class()))
837 ((MFMagicCuts*)filter)->InitMapping(&m);
838
839 MReadTree read("Events");
840 read.DisableAutoScheme(); // AutoScheme doesn't seem to be faster!
841 if (fname)
842 read.AddFile(fname);
843 else
844 AddSequences(read, fNamesOn);
845 if (!FillMatrix(read, parlist))
846 return kFALSE;
847
848 MTaskList tasklist;
849 parlist.Replace(&tasklist);
850 if (fit)
851 parlist.AddToList(fit);
852
853 MFilterList list;
854 SetupFilters(list, filter);
855
856 MContinue contin(&list);
857 parlist.AddToList(&list);
858
859 MFillH fill(&hist);
860
861 MMatrixLoop loop(&m);
862
863 tasklist.AddToList(&loop);
864 tasklist.AddToList(&list);
865 tasklist.AddToList(&contin);
866 tasklist.AddToList(&fill);
867
868 if (!Optimize(parlist))
869 return kFALSE;
870
871 TObjArray cont;
872 cont.Add(&contin);
873 return WriteContainer(cont, fNameOut);
874}
875
876//------------------------------------------------------------------------
877//
878// Make sure, that filter->GetDataMember is correctly implemented!!!!
879//
880Bool_t MJOptimize::RunOnOff(const char *fname, MFilter *filter, MAlphaFitter *fit, const char *tree)
881{
882 if (fParameters.GetSize()==0)
883 {
884 *fLog << err << "Sorry, no parameters defined." << endl;
885 return kFALSE;
886 }
887
888 fLog->Separator("Preparing On/Off-optimization");
889
890 MParList parlist;
891
892 MGeomCamMagic geom; // For GetConvMm2Deg
893 parlist.AddToList(&geom);
894
895 MHMatrix m("M");
896 AddRulesToMatrix(m);
897 parlist.AddToList(&m);
898
899 const Int_t idxdatatype = m.AddColumn("DataType.fVal");
900
901 MHAlpha histon, histof("MHAlphaOff");
902 histon.SkipHistTime();
903 histon.SkipHistTheta();
904 //histon.SkipHistEnergy();
905 histof.SkipHistTime();
906 histof.SkipHistTheta();
907 //histof.SkipHistEnergy();
908 histon.InitMapping(&m, 0);
909 histof.InitMapping(&m, 0);
910
911 if (filter && filter->InheritsFrom(MFMagicCuts::Class()))
912 ((MFMagicCuts*)filter)->InitMapping(&m);
913
914 parlist.AddToList(&histon);
915 parlist.AddToList(&histof);
916
917 if (fname)
918 {
919 MReadTree read(tree);
920 read.DisableAutoScheme(); // AutoScheme doesn't seem to be faster!
921 read.AddFile(fname);
922 if (!FillMatrix(read, parlist))
923 return kFALSE;
924 }
925 else
926 {
927 MParameterI par("DataType");
928 parlist.AddToList(&par);
929
930 gLog.Separator("Reading On-Data");
931 par.SetVal(1);
932 MReadTree readon(tree);
933 readon.DisableAutoScheme(); // AutoScheme doesn't seem to be faster!
934 AddSequences(readon, fNamesOn);
935 if (!FillMatrix(readon, parlist))
936 return kFALSE;
937
938 gLog.Separator("Reading Off-Data");
939 par.SetVal(0);
940 MReadTree readoff(tree);
941 readoff.DisableAutoScheme(); // AutoScheme doesn't seem to be faster!
942 AddSequences(readoff, fNamesOff);
943 if (!FillMatrix(readoff, parlist))
944 return kFALSE;
945 }
946
947 MTaskList tasklist;
948 parlist.Replace(&tasklist);
949 if (fit)
950 parlist.AddToList(fit);
951
952 MFilterList list;
953 SetupFilters(list, filter);
954
955 MContinue contin(&list);
956 parlist.AddToList(&list);
957
958 MFillH fillon(&histon, "", "FillHistOn");
959 MFillH fillof(&histof, "", "FillHistOff");
960
961 MF f0(Form("M[%d]<0.5", idxdatatype), "FilterOffData");
962 MF f1(Form("M[%d]>0.5", idxdatatype), "FilterOnData");
963
964 fillof.SetFilter(&f0);
965 fillon.SetFilter(&f1);
966
967 MMatrixLoop loop(&m);
968
969 tasklist.AddToList(&loop);
970 tasklist.AddToList(&list);
971 tasklist.AddToList(&contin);
972 tasklist.AddToList(&f0);
973 tasklist.AddToList(&f1);
974 tasklist.AddToList(&fillof);
975 tasklist.AddToList(&fillon);
976
977 if (!Optimize(parlist))
978 return kFALSE;
979
980 TObjArray cont;
981 cont.Add(&contin);
982 return WriteContainer(cont, fNameOut);
983}
984
985//------------------------------------------------------------------------
986//
987// Read all events from file which do match rules and optimize
988// energy estimator.
989//
990Bool_t MJOptimize::RunEnergy(const char *fname, const char *rule)
991{
992 if (fParameters.GetSize()==0)
993 {
994 *fLog << err << "Sorry, no parameters defined." << endl;
995 return kFALSE;
996 }
997
998 fLog->Separator("Preparing Energy optimization");
999
1000 MParList parlist;
1001
1002 MParameterI par("DataType");
1003 par.SetVal(1);
1004 parlist.AddToList(&par);
1005
1006 MFDataMember filter("DataType.fVal", '>', 0.5);
1007 fPreCuts.Add(&filter);
1008
1009 MGeomCamMagic geom; // For GetConvMm2Deg
1010 parlist.AddToList(&geom);
1011
1012 MHMatrix m("M");
1013 AddRulesToMatrix(m);
1014 parlist.AddToList(&m);
1015
1016 MHEnergyEst hist;
1017 hist.InitMapping(&m);
1018
1019 MEnergyEstimate est("MParameters");
1020 est.SetRule(rule);
1021 parlist.AddToList(&est);
1022
1023 MReadTree read("Events");
1024 // NECESSARY BECAUSE OF MDataFormula GetRules missing
1025 read.DisableAutoScheme();
1026 if (fname)
1027 read.AddFile(fname);
1028 else
1029 AddSequences(read, fNamesOn);
1030 if (!FillMatrix(read, parlist, kTRUE))
1031 return kFALSE;
1032
1033 fPreCuts.Remove(&filter);
1034
1035 MTaskList tasklist;
1036 parlist.Replace(&tasklist);
1037
1038 MFillH fill(&hist);
1039
1040 MMatrixLoop loop(&m);
1041
1042 tasklist.AddToList(&loop);
1043 tasklist.AddToList(&est);
1044 tasklist.AddToList(&fill);
1045
1046 if (!Optimize(parlist))
1047 return kFALSE;
1048
1049 hist.Print();
1050
1051 TObjArray cont;
1052 cont.Add(&est);
1053 return WriteContainer(cont, fNameOut);
1054}
Note: See TracBrowser for help on using the repository browser.