source: branches/removing_cpp11_features/mtemp/mmpi/macros/findstars.C@ 18730

Last change on this file since 18730 was 4440, checked in by rwagner, 20 years ago
*** empty log message ***
File size: 5.0 KB
Line 
1void ReadSetup(TString fname, MAstroCamera &cam)
2{
3 MMcConfigRunHeader *config=0;
4 MGeomCam *geom=0;
5
6 TFile file(fname);
7 TTree *tree = (TTree*)file.Get("RunHeaders");
8 tree->SetBranchAddress("MMcConfigRunHeader", &config);
9 if (tree->GetBranch("MGeomCam"))
10 tree->SetBranchAddress("MGeomCam", &geom);
11 tree->GetEntry(0);
12
13 cam.SetMirrors(*config->GetMirrors());
14 cam.SetGeom(*geom);
15}
16
17void findstars(const TString filename="20040422_23213_D_Mrk421_E.root", const TString directory="/data/MAGIC/Period016/rootdata2/2004_04_22/", const UInt_t numEvents = 0)
18{
19
20 MParList plist;
21 MTaskList tlist;
22 plist.AddToList(&tlist);
23
24 MGeomCamMagic geomcam;
25 MCameraDC dccam;
26 MStarLocalCam starcam;
27
28 plist.AddToList(&geomcam);
29 plist.AddToList(&dccam);
30 plist.AddToList(&starcam);
31
32 // Reads the trees of the root file and the analysed branches
33 MReadReports read;
34 read.AddTree("Currents");
35 read.AddTree("Drive"); // If you do not include Drive info, the MFindStars class
36 // will not use the star catalog method
37 read.AddFile(directory+filename); // after the reading of the trees!!!
38 read.AddToBranchList("MReportCurrents.*");
39 read.AddToBranchList("MReportDrive.*");
40
41 MGeomApply geomapl;
42 TString continuoslightfile =
43 "/data/MAGIC/Period016/rootdata/2004_04_16/20040416_22368_P_Off3c279-2CL100_E.root";
44 Float_t mindc = 0.9; //[uA]
45
46 MCalibrateDC dccal;
47 dccal.SetFileName(continuoslightfile);
48 dccal.SetMinDCAllowed(mindc);
49
50 const Int_t numblind = 5;
51 const Short_t x[numblind] = { 47, 124, 470, 475, 571};
52 const TArrayS blindpixels(numblind,(Short_t*)x);
53 Float_t ringinterest = 100; //[mm]
54 Float_t tailcut = 2.5;
55 UInt_t integratedevents = 1;
56
57 // We need the MAGIC mirror geometry from a MC header:
58 TString geometryfile = "/mcdata/standard/camera/NSB_013/Gamma/Gamma_zbin9_90_7_1480to1489_w0.root";
59
60 // We need the Bright Star Catalog:
61 TString catalogfile = "/home/rwagner/bsc5.dat";
62
63 MFindStars findstars;
64 findstars.SetBlindPixels(blindpixels);
65 findstars.SetRingInterest(ringinterest);
66 findstars.SetDCTailCut(tailcut);
67 findstars.SetNumIntegratedEvents(integratedevents);
68 findstars.SetMinuitPrintOutLevel(-1);
69 findstars.SetGeometryFile(geometryfile);
70 findstars.SetBSCFile(catalogfile);
71 const Double_t ra = MAstro::Hms2Rad(11, 4, 26);
72 const Double_t dec = MAstro::Dms2Rad(38, 12, 36);
73 findstars.SetRaDec(ra,dec);
74 findstars.SetLimMag(8);
75 findstars.SetRadiusFOV(1.5);
76
77 tlist.AddToList(&geomapl);
78 tlist.AddToList(&read);
79 tlist.AddToList(&dccal);
80 tlist.AddToList(&findstars, "Currents");
81
82 // The following lines you only need if in addition you want to display
83 // independent MAstroCamera output
84 //
85 // TString fname = "/mcdata/standard/camera/NSB_013/Gamma/Gamma_zbin9_90_7_1480to1489_w0.root";
86 // MObservatory magic1;
87 // const Double_t ra = MAstro::Hms2Rad(11, 4, 26); //Mkn421
88 // const Double_t dec = MAstro::Dms2Rad(38, 12, 36);
89 //
90 // MAstroCamera stars;
91 // ReadSetup(fname, stars);
92 // stars.SetLimMag(9);
93 // stars.SetRadiusFOV(3);
94 // stars.SetRaDec(ra, dec);
95 // stars.ReadBSC("/home/rwagner/bsc5.dat");
96 // stars.SetObservatory(magic1);
97
98 MEvtLoop evtloop;
99 evtloop.SetParList(&plist);
100
101 if (!evtloop.PreProcess())
102 return;
103
104 MHCamera display(geomcam);
105 display.SetPrettyPalette();
106 display.Draw();
107 gPad->cd(1);
108 starcam.Draw();
109
110 UInt_t numevents=0;
111
112 while (tlist.Process())
113 {
114 numevents++;
115 if (numevents%integratedevents==0)
116 {
117 display.SetCamContent(findstars.GetDisplay());
118 gPad->Modified();
119 gPad->Update();
120 // This line prints the results:
121 // starcam.Print();
122 // This is how to access the TList of stars:
123 // TList* starlist = starcam.GetList();
124
125 // This is how to iterate over stars found:
126 // TIter Next(starlist);
127 // MStarLocalPos* star;
128 // UInt_t starnum = 0;
129 // cout << filename << " ";
130 // cout << "Iterating over list" << endl;
131 // while ((star=(MStarLocalPos*)Next()))
132 // {
133 // cout << "star[" << starnum << "] ";
134 // cout << star->GetMeanX() << " "
135 // << star->GetMeanY() << " ";
136 // starnum++;
137 // }
138 // cout << endl;
139
140 }//integratedevents
141
142 MTime time;
143 time.Set(2004, 4, 22, 21, 51, 15);
144
145 //superimpose star picture
146 // stars.SetTime(time);
147 // TObject *o = stars.Clone();
148 // o->SetBit(kCanDelete);
149 // o->Draw();
150
151 if (!HandleInput())
152 break;
153
154 }
155
156 evtloop.PostProcess();
157
158 tlist.PrintStatistics();
159
160}
161
162
163Bool_t HandleInput()
164{
165 TTimer timer("gSystem->ProcessEvents();", 50, kFALSE);
166
167 while (1)
168 {
169 //
170 // While reading the input process gui events asynchronously
171 //
172 timer.TurnOn();
173 TString input = Getline("Type 'q' to exit, <return> to go on: ");
174 timer.TurnOff();
175
176 if (input=="q\n")
177 return kFALSE;
178
179 if (input=="\n")
180 return kTRUE;
181 };
182
183 return kFALSE;
184}
Note: See TracBrowser for help on using the repository browser.