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

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