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

Last change on this file since 6921 was 6907, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 27.3 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)
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 fPreCuts.Add(new MF(rule));
336}
337
338//------------------------------------------------------------------------
339//
340// Set the fParameters Array accoring to par.
341//
342void MJOptimize::SetParameters(const TArrayD &par)
343{
344 fParameters = par;
345}
346
347//------------------------------------------------------------------------
348//
349// Set the number of events processed by the eventloop. (Be carfull,
350// if you are doing on-off analysis and you only process the first
351// 1000 events which are on-events only the optimization may not work)
352//
353void MJOptimize::SetNumEvents(UInt_t n)
354{
355 fNumEvents=n;
356}
357
358//------------------------------------------------------------------------
359//
360// Set a debug level, which tells the optimization how much information
361// is displayed about and in the running eventloop.
362//
363void MJOptimize::SetDebug(UInt_t n)
364{
365 fDebug=n;
366}
367
368//------------------------------------------------------------------------
369//
370// Set a optimization algorithm to be used. For more information see
371// TMinuit.
372//
373// Available Algorithms are:
374// kMigrad, // Minimize by the method of Migrad
375// kSimplex, // Minimize by the method of Simplex
376// kSeek // Minimize by the method of Monte Carlo
377//
378void MJOptimize::SetOptimizer(Optimizer_t o)
379{
380 fType = o;
381}
382
383//------------------------------------------------------------------------
384//
385// If a status didplay is set, search for tab "Optimizer".
386// If not found, create it.
387// In the tab search for TMultiGraph "Parameters".
388// If not found create it.
389// If empty create TGraphs.
390// Check number of graphs vs. number of parameters.
391// return TList with graphs.
392//
393TList *MJOptimize::GetPlots() const
394{
395 if (!fDisplay)
396 return NULL;
397
398 TCanvas *c=fDisplay->GetCanvas("Optimizer");
399 if (!c)
400 c = &fDisplay->AddTab("Optimizer");
401
402 TMultiGraph *mg = dynamic_cast<TMultiGraph*>(c->FindObject("Parameters"));
403 if (!mg)
404 mg = new TMultiGraph("Parameters", "Parameters of optimization");
405
406 TList *l = mg->GetListOfGraphs();
407 if (!l)
408 {
409 const Int_t n = fParameters.GetSize();
410 for (int i=0; i<n+1; i++)
411 {
412 TGraph *g = new TGraph;
413 if (i==n)
414 g->SetLineColor(kBlue);
415 mg->Add(g, "");
416 AddPoint(mg->GetListOfGraphs(), i, i==n?1:fParameters[i]);
417 }
418 mg->SetBit(kCanDelete);
419 mg->Draw("al*");
420
421 l = mg->GetListOfGraphs();
422 }
423
424 return l->GetSize() == fParameters.GetSize()+1 ? l : NULL;
425}
426
427//------------------------------------------------------------------------
428//
429// Add a point with y=val as last point in idx-th Tgraph of list l.
430//
431void MJOptimize::AddPoint(TList *l, Int_t idx, Float_t val) const
432{
433 if (!l)
434 return;
435
436 TGraph *gr = (TGraph*)l->At(idx);
437 gr->SetPoint(gr->GetN(), gr->GetN(), val);
438}
439
440Int_t MJOptimize::Minuit(TMinuit &minuit, const char *cmd) const
441{
442 Int_t err;
443 Double_t tmp[2] = { fNumMaxCalls, fTolerance };
444 minuit.mnexcm(cmd, tmp, 2, err);
445
446 switch (err)
447 {
448 case 0:
449 *fLog << inf << GetDescriptor() << " TMinuit::mnexcm excuted normally." << endl;
450 break;
451 case 1:
452 *fLog << warn << GetDescriptor() << " TMinuit::mnexcm command is blank... ignored." << endl;
453 break;
454 case 2:
455 *fLog << warn << GetDescriptor() << " TMinuit::mnexcm command-line syntax error... ignored." << endl;
456 break;
457 case 3:
458 *fLog << warn << GetDescriptor() << " TMinuit::mnexcm unknown command... ignored." << endl;
459 break;
460 case 4:
461 *fLog << warn << GetDescriptor() << " TMinuit::mnexcm - Abnormal termination (eg Migrad not converged)" << endl;
462 break;
463 /*
464 case 5:
465 *fLog << inf << GetDescriptor() << " TMinuit::mnexcm - Parameters requested." << endl;
466 break;
467 case 6:
468 *fLog << inf << GetDescriptor() << " TMinuit::mnexcm - SET INPUT returned." << endl;
469 break;
470 case 7:
471 *fLog << inf << GetDescriptor() << " TMinuit::mnexcm - SET TITLE returned." << endl;
472 break;
473 case 8:
474 *fLog << inf << GetDescriptor() << " TMinuit::mnexcm - SET COVAR returned." << endl;
475 break;
476 case 9:
477 *fLog << inf << GetDescriptor() << " TMinuit::mnexcm - reserved." << endl;
478 break;
479 case 10:
480 *fLog << inf << GetDescriptor() << " TMinuit::mnexcm - END returned." << endl;
481 break;
482 case 11:
483 *fLog << inf << GetDescriptor() << " TMinuit::mnexcm - EXIT or STOP returned." << endl;
484 break;
485 case 12:
486 *fLog << inf << GetDescriptor() << " TMinuit::mnexcm - RETURN returned." << endl;
487 break;*/
488 }
489
490 return err;
491}
492
493Bool_t MJOptimize::Optimize(MEvtLoop &evtloop)
494{
495 if (fType==kNone)
496 return kTRUE;
497
498 gMinuit = new TMinuit(fParameters.GetSize());
499
500 gMinuit->SetFCN(fcn);
501 gMinuit->SetObjectFit(this);
502 gMinuit->SetPrintLevel(-1); // Don't print when DefineParameter
503
504 //
505 // Set starting values and step sizes for parameters
506 //
507 for (Int_t i=0; i<fParameters.GetSize(); i++)
508 {
509 TString name = "par[";
510 name += i;
511 name += "]";
512 Double_t vinit = fParameters[i];
513 Double_t step = fStep[i];
514
515 Double_t limlo = fLimLo[i];
516 Double_t limup = fLimUp[i];
517
518 Bool_t rc = gMinuit->DefineParameter(i, name, vinit, step, limlo, limup);
519 if (rc)
520 {
521 *fLog << err << dbginf << "Error in defining parameter #" << i << endl;
522 return kFALSE;
523 }
524
525 if (step==0)
526 gMinuit->FixParameter(i);
527 }
528
529 gMinuit->SetPrintLevel(1); // Switch on pritning again
530 gMinuit->mnprin(1,0); // Print all parameters
531
532 fEvtLoop = &evtloop;
533
534 TList *g=GetPlots();
535
536 // Now ready for minimization step:
537 TStopwatch clock;
538 clock.Start();
539 switch (fType)
540 {
541 case kSimplex:
542 Simplex(*gMinuit);
543 break;
544 case kMigrad:
545 Migrad(*gMinuit);
546 break;
547 case kMinimize:
548 Minimize(*gMinuit);
549 break;
550 case kMinos:
551 Minos(*gMinuit);
552 break;
553 case kImprove:
554 Improve(*gMinuit);
555 break;
556 case kSeek:
557 Seek(*gMinuit);
558 break;
559 case kNone: // Should never happen
560 return kFALSE;
561 }
562 clock.Stop();
563 clock.Print();
564
565 if (evtloop.GetDisplay()!=fDisplay)
566 {
567 *fLog << inf << "Optimization aborted by user." << endl;
568 fDisplay = 0;
569 return kFALSE;
570 }
571
572 *fLog << inf << "Resulting Chisq: " << gMinuit->fAmin << endl;
573
574 //
575 // Update values of fA, fB:
576 //
577 for (Int_t i=0; i<fParameters.GetSize(); i++)
578 {
579 Double_t x1, x2;
580 gMinuit->GetParameter(i,x1,x2);
581 fParameters[i] = x1;
582 cout << i << ": " << fParameters[i] << endl;
583
584 AddPoint(g, i, x1);
585 }
586 AddPoint(g, fParameters.GetSize(), gMinuit->fAmin);
587
588 delete gMinuit;
589
590 return kTRUE;
591}
592
593//------------------------------------------------------------------------
594//
595// Optimize allows to use the optimizing by an eventloop based on
596// some requirements.
597//
598// 1) The tasklist to be executed must have the name MTaskList and
599// be an entry in the parameterlist.
600//
601// 2) The reading task (MReadMarsFile, MMatrixLoop) must have the name
602// "MRead". If it derives from MRead Rewind() must be implemented,
603// otherwise it must start reading from scratch each time its
604// PreProcess is called.
605//
606// 3) The parameters to be optimized must be accessed through (currently)
607// a single parameter container (MParContainer) called MParameters.
608// The parameters are set through SetVariables.
609//
610// 4) The result of a single function call for minimization (eg. chisquare)
611// must be available after the eventloop in a container of type
612// MParameterD with the name "MinimizationResult".
613//
614// 5) The parameters to start with must have been set using
615// MJOptimize::SetParameter or MJOptimize::SetParameters and
616// MJOptimize::FixParameter
617//
618// The behaviour of the optimization can be changed using:
619// void SetNumEvents(UInt_t n);
620// void SetDebug(UInt_t n);
621// void SetOptimizer(Optimizer_t o);
622//
623// After optimization the resulting parameters are set and another eventloop
624// with a MStatusDisplay is set is called. The resulting parameters can be
625// accessed using: GetParameters()
626//
627// To be fixed:
628// - MStatusDisplay should show status while optimization is running
629// - write result into MStatusDisplay
630// - save result
631//
632Bool_t MJOptimize::Optimize(MParList &parlist)
633{
634 // Checks to make sure, that fcn doesn't crash
635 if (!parlist.FindCreateObj("MParameterD", "MinimizationValue"))
636 return kFALSE;
637
638 if (!parlist.FindObject("MParameters", "MParContainer"))
639 {
640 *fLog << err << "MParameters [MParContainer] not found... abort." << endl;
641 return kFALSE;
642 }
643
644 TString txt("Starting ");
645 switch (fType)
646 {
647 case kMigrad: txt += "Migrad"; break;
648 case kMinimize: txt += "Minimize"; break;
649 case kMinos: txt += "Minot"; break;
650 case kImprove: txt += "Improve"; break;
651 case kSimplex: txt += "Simplex"; break;
652 case kSeek: txt += "Seek"; break;
653 case kNone: txt += "no"; break;
654 }
655 txt += " optimization";
656
657 fLog->Separator(txt);
658
659 // Setup eventloop
660 MEvtLoop evtloop;
661 evtloop.SetParList(&parlist);
662 evtloop.SetDisplay(fDisplay); // set display for evtloop and all childs
663 parlist.SetDisplay(0); // reset display for all contents of parlist and tasklist
664 evtloop.SetPrivateDisplay(); // prevent display from being cascaded again in PreProcess
665
666 *fLog << inf << "Number of Parameters: " << fParameters.GetSize() << endl;
667
668 if (!Optimize(evtloop))
669 return kFALSE;
670
671 gMinuit = 0;
672
673 fEvtLoop->SetDisplay(fDisplay);
674 return Fcn(fParameters);
675 // Done already in Fcn
676 // list.SetVariables(fParameters);
677}
678
679void MJOptimize::AddRulesToMatrix(MHMatrix &m) const
680{
681 TIter Next1(&fRules);
682 TObject *o1=0;
683 while ((o1=Next1()))
684 m.AddColumn(o1->GetName());
685}
686
687//------------------------------------------------------------------------
688//
689// Fill matrix m by using read. Use rules as a filter if userules.
690//
691Bool_t MJOptimize::FillMatrix(MReadTree &read, MParList &parlist, Bool_t userules)
692{
693 MHMatrix *m = (MHMatrix*)parlist.FindObject("M", "MHMatrix");
694 if (!m)
695 {
696 *fLog << err << "MJOptimize::FillMatrix - ERROR: M [MHMatrix] not found in parlist... abort." << endl;
697 return kFALSE;
698 }
699
700 m->Print("cols");
701
702 //MParList parlist;
703
704 // MGeomCamMagic cam;
705 // parlist.AddToList(&cam);
706
707 MTaskList tlist;
708 parlist.Replace(&tlist);
709
710 MFillH fillh(m);
711
712 tlist.AddToList(&read);
713
714 MFilterList list;
715 if (!list.AddToList(fPreCuts))
716 *fLog << err << "ERROR - Calling MFilterList::AddToList for fPreCuts failed!" << endl;
717 if (userules)
718 SetupFilters(list);
719 list.SetName("PreCuts"); // reset Name set by SetupFilters
720 list.SetInverted(kFALSE); // reset inversion set by SetupFilters
721 fillh.SetFilter(&list);
722 tlist.AddToList(&list);
723
724 tlist.AddToList(&fillh);
725
726 MEvtLoop fillloop;
727 fillloop.SetParList(&parlist);
728 fillloop.SetDisplay(fDisplay);
729 if (!fillloop.Eventloop(fNumEvents))
730 {
731 *fLog << err << "Filling matrix failed..." << endl;
732 return kFALSE;
733 }
734 tlist.PrintStatistics();
735
736 *fLog << inf << "Read events from file '" << read.GetFileName() << "'" << endl;
737
738 if (fillloop.GetDisplay()!=fDisplay)
739 {
740 fDisplay = 0;
741 *fLog << inf << "Optimization aborted by user." << endl;
742 return kFALSE;
743 }
744
745 m->Print("size");
746
747 return kTRUE;
748}
749
750//------------------------------------------------------------------------
751//
752// Adds all filters to MFilterList
753//
754void MJOptimize::SetupFilters(MFilterList &list, MFilter *filter) const
755{
756 list.SetName("MParameters");
757 list.SetInverted();
758
759 if (filter)
760 {
761 if (fFilter.GetSize()>0)
762 {
763 *fLog << inf;
764 *fLog << "INFORMATION - You are using an external filter and internal filters." << endl;
765 *fLog << " Please make sure that all parameters '[i]' are starting" << endl;
766 *fLog << " behind the number of parameters of the external filter." << endl;
767 }
768 list.AddToList(filter);
769 }
770
771 if (!list.AddToList(fFilter))
772 *fLog << err << "ERROR - Calling MFilterList::AddToList fFilter failed!" << endl;
773
774 *fLog << inf << "Filter: ";
775 list.Print();
776 *fLog << endl;
777}
778
779//------------------------------------------------------------------------
780//
781Bool_t MJOptimize::Run(const char *fname, MFilter *filter, MAlphaFitter *fit)
782{
783 if (fParameters.GetSize()==0)
784 {
785 *fLog << err << "Sorry, no parameters defined." << endl;
786 return kFALSE;
787 }
788
789 fLog->Separator("Preparing On-only-optimization");
790
791 MParList parlist;
792
793 MGeomCamMagic geom; // For GetConvMm2Deg
794 parlist.AddToList(&geom);
795
796 MHMatrix m("M");
797 AddRulesToMatrix(m);
798 parlist.AddToList(&m);
799
800 MHAlpha hist;
801 hist.SkipHistTime();
802 hist.SkipHistTheta();
803 hist.SkipHistEnergy();
804 hist.InitMapping(&m);
805
806 if (filter && filter->InheritsFrom(MFMagicCuts::Class()))
807 ((MFMagicCuts*)filter)->InitMapping(&m);
808
809 MReadTree read("Events");
810 read.DisableAutoScheme(); // AutoScheme doesn't seem to be faster!
811 if (fname)
812 read.AddFile(fname);
813 else
814 AddSequences(read, fNamesOn);
815 if (!FillMatrix(read, parlist))
816 return kFALSE;
817
818 MTaskList tasklist;
819 parlist.Replace(&tasklist);
820 if (fit)
821 parlist.AddToList(fit);
822
823 MFilterList list;
824 SetupFilters(list, filter);
825
826 MContinue contin(&list);
827 parlist.AddToList(&list);
828
829 MFillH fill(&hist);
830
831 MMatrixLoop loop(&m);
832
833 tasklist.AddToList(&loop);
834 tasklist.AddToList(&list);
835 tasklist.AddToList(&contin);
836 tasklist.AddToList(&fill);
837
838 if (!Optimize(parlist))
839 return kFALSE;
840
841 TObjArray cont;
842 cont.Add(&contin);
843 return WriteContainer(cont, fNameOut);
844}
845
846//------------------------------------------------------------------------
847//
848// Make sure, that filter->GetDataMember is correctly implemented!!!!
849//
850Bool_t MJOptimize::RunOnOff(const char *fname, MFilter *filter, MAlphaFitter *fit, const char *tree)
851{
852 if (fParameters.GetSize()==0)
853 {
854 *fLog << err << "Sorry, no parameters defined." << endl;
855 return kFALSE;
856 }
857
858 fLog->Separator("Preparing On/Off-optimization");
859
860 MParList parlist;
861
862 MGeomCamMagic geom; // For GetConvMm2Deg
863 parlist.AddToList(&geom);
864
865 MHMatrix m("M");
866 AddRulesToMatrix(m);
867 parlist.AddToList(&m);
868
869 const Int_t idxdatatype = m.AddColumn("DataType.fVal");
870
871 MHAlpha histon, histof("MHAlphaOff");
872 histon.SkipHistTime();
873 histon.SkipHistTheta();
874 //histon.SkipHistEnergy();
875 histof.SkipHistTime();
876 histof.SkipHistTheta();
877 //histof.SkipHistEnergy();
878 histon.InitMapping(&m, 0);
879 histof.InitMapping(&m, 0);
880
881 if (filter && filter->InheritsFrom(MFMagicCuts::Class()))
882 ((MFMagicCuts*)filter)->InitMapping(&m);
883
884 parlist.AddToList(&histon);
885 parlist.AddToList(&histof);
886
887 if (fname)
888 {
889 MReadTree read(tree);
890 read.DisableAutoScheme(); // AutoScheme doesn't seem to be faster!
891 read.AddFile(fname);
892 if (!FillMatrix(read, parlist))
893 return kFALSE;
894 }
895 else
896 {
897 MParameterI par("DataType");
898 parlist.AddToList(&par);
899
900 gLog.Separator("Reading On-Data");
901 par.SetVal(1);
902 MReadTree readon(tree);
903 readon.DisableAutoScheme(); // AutoScheme doesn't seem to be faster!
904 AddSequences(readon, fNamesOn);
905 if (!FillMatrix(readon, parlist))
906 return kFALSE;
907
908 gLog.Separator("Reading Off-Data");
909 par.SetVal(0);
910 MReadTree readoff(tree);
911 readoff.DisableAutoScheme(); // AutoScheme doesn't seem to be faster!
912 AddSequences(readoff, fNamesOff);
913 if (!FillMatrix(readoff, parlist))
914 return kFALSE;
915 }
916
917 MTaskList tasklist;
918 parlist.Replace(&tasklist);
919 if (fit)
920 parlist.AddToList(fit);
921
922 MFilterList list;
923 SetupFilters(list, filter);
924
925 MContinue contin(&list);
926 parlist.AddToList(&list);
927
928 MFillH fillon(&histon, "", "FillHistOn");
929 MFillH fillof(&histof, "", "FillHistOff");
930
931 MF f0(Form("M[%d]<0.5", idxdatatype), "FilterOffData");
932 MF f1(Form("M[%d]>0.5", idxdatatype), "FilterOnData");
933
934 fillof.SetFilter(&f0);
935 fillon.SetFilter(&f1);
936
937 MMatrixLoop loop(&m);
938
939 tasklist.AddToList(&loop);
940 tasklist.AddToList(&list);
941 tasklist.AddToList(&contin);
942 tasklist.AddToList(&f0);
943 tasklist.AddToList(&f1);
944 tasklist.AddToList(&fillof);
945 tasklist.AddToList(&fillon);
946
947 if (!Optimize(parlist))
948 return kFALSE;
949
950 TObjArray cont;
951 cont.Add(&contin);
952 return WriteContainer(cont, fNameOut);
953}
954
955//------------------------------------------------------------------------
956//
957// Read all events from file which do match rules and optimize
958// energy estimator.
959//
960Bool_t MJOptimize::RunEnergy(const char *fname, const char *rule)
961{
962 if (fParameters.GetSize()==0)
963 {
964 *fLog << err << "Sorry, no parameters defined." << endl;
965 return kFALSE;
966 }
967
968 fLog->Separator("Preparing Energy optimization");
969
970 MParList parlist;
971
972 MGeomCamMagic geom; // For GetConvMm2Deg
973 parlist.AddToList(&geom);
974
975 MHMatrix m("M");
976 AddRulesToMatrix(m);
977 parlist.AddToList(&m);
978
979 MHEnergyEst hist;
980 hist.InitMapping(&m);
981
982 MEnergyEstimate est("MParameters");
983 est.SetRule(rule);
984 parlist.AddToList(&est);
985
986 MReadTree read("Events");
987 // NECESSARY BECAUSE OF MDataFormula GetRules missing
988 read.DisableAutoScheme();
989 if (fname)
990 read.AddFile(fname);
991 else
992 AddSequences(read, fNamesOn);
993 if (!FillMatrix(read, parlist, kTRUE))
994 return kFALSE;
995
996 MTaskList tasklist;
997 parlist.Replace(&tasklist);
998
999 MFillH fill(&hist);
1000
1001 MMatrixLoop loop(&m);
1002
1003 tasklist.AddToList(&loop);
1004 tasklist.AddToList(&est);
1005 tasklist.AddToList(&fill);
1006
1007 if (!Optimize(parlist))
1008 return kFALSE;
1009
1010 TObjArray cont;
1011 cont.Add(&est);
1012 return WriteContainer(cont, fNameOut);
1013}
Note: See TracBrowser for help on using the repository browser.