KaliVeda
Toolkit for HIC analysis
KVIDGridManager.cpp
1 /***************************************************************************
2  KVIDGridManager.cpp - description
3  -------------------
4  begin : Jan 27 2005
5  copyright : (C) 2005 by J.D. Frankland
6  email : frankland@ganil.fr
7 
8 $Id: KVIDGridManager.cpp,v 1.13 2009/03/03 14:27:15 franklan Exp $
9 ***************************************************************************/
10 
11 /***************************************************************************
12  * *
13  * This program is free software; you can redistribute it and/or modify *
14  * it under the terms of the GNU General Public License as published by *
15  * the Free Software Foundation; either version 2 of the License, or *
16  * (at your option) any later version. *
17  * *
18  ***************************************************************************/
19 
20 #include "KVIDGridManager.h"
21 #include "KVString.h"
22 #include "TClass.h"
23 #include "TROOT.h"
24 
25 using namespace std;
26 
28 KVIDGridManager* gIDGridManager = nullptr;
29 
30 
38 
40  : fGrids{"TList"}
41 {
42  //Default constructor
43  //Initialise global pointer gIDGridManager
44  //
45  //Create list for ID grids
46  //
47  //Will delete grids if owns_grids=true
48  gIDGridManager = this;
50  fGrids.Connect("Modified()", "KVIDGridManager", this, "Modified()");
51  fGrids.SetOwner(owns_grids);
52 }
53 
54 
55 
59 
60 KVIDGridManager::~KVIDGridManager()
61 {
62  //Destructor
63  //Reset global pointer gIDGridManager
64 
65  if (gIDGridManager == this)
66  gIDGridManager = 0;
67  fGrids.Disconnect("Modified()", this, "Modified()");
69 }
70 
71 
72 
75 
77 {
78  // Add a grid to the collection
79 
80  fGrids.Add(grid);
81 }
82 
83 
84 
89 
91 {
92  //Remove grid from manager's list and delete it
93  //update flag allows to disable the emission of the 'Modified' signal in case the GUI
94  //is deleting a list of grids - in this case we don't want to update until the end
95 
96  if (!update) fGrids.Disconnect("Modified()", this, "Modified()");
97  fGrids.Remove(grid);
98  delete grid;
99  if (!update) fGrids.Connect("Modified()", "KVIDGridManager", this, "Modified()");
100 }
101 
102 
103 
106 
108 {
109  //Delete all grids and empty list, ready to start anew
110 
111  fGrids.Disconnect("Modified()", this, "Modified()");
113  fGrids.Clear();
114  Modified(); // emit signal to say something changed
115  fGrids.Connect("Modified()", "KVIDGridManager", this, "Modified()");
116 }
117 
118 
119 
131 
133 {
134  //read file, create grids corresponding to information in file.
135  //note: any existing grids are not destroyed, use Clear() beforehand if you want to
136  //start afresh and anew (ais athat aOK?)
137  //
138  // N.B. the links between each grid and the IDtelescope(s) for which it is to be used
139  // are not set up by this method.
140  //
141  //the list of grids created by reading the file can be accessed with method
142  //GetLastReadGrids() after calling this method
143 
144  // clear list of read grids
146 
147  Bool_t is_it_ok = kFALSE;
148  ifstream gridfile(filename);
149  if (!gridfile.good()) {
150  Error("ReadAsciiFile", "File %s cannot be opened", filename);
151  return is_it_ok;
152  }
153 
154  KVString s;
155 
156  fGrids.Disconnect("Modified()", this, "Modified()");
157  while (gridfile.good()) {
158  //read a line
159  s.ReadLine(gridfile);
160  if (s.BeginsWith("++")) {
161  //New grid
162  //Get name of class by stripping off the '+' at the start of the line
163  s.Remove(0, 2);
164  /************ BACKWARDS COMPATIBILITY FIX *************
165  Old grid files may contain obsolete KVIDZGrid class
166  We replace by KVIDZAGrid with SetOnlyZId(kTRUE)
167  */
168  //Make new grid using this class
169  KVIDGraph* grid = 0;
170  Bool_t onlyz = kFALSE;
171  if (s == "KVIDZGrid") {
172  s = "KVIDZAGrid";
173  onlyz = kTRUE;
174  }
175  TClass* clas = TClass::GetClass(s.Data());
176  if (!clas) {
177  Fatal("ReadAsciiFile",
178  "Cannot load TClass information for %s", s.Data());
179  }
180  grid = (KVIDGraph*) clas->New();
181  fLastReadGrids.Add(grid);
182  //read grid
183  grid->ReadFromAsciiFile(gridfile);
184  if (onlyz) grid->SetOnlyZId(kTRUE);
185  }
186  }
187 
188  gridfile.close();
189  is_it_ok = kTRUE;
190  Modified(); // emit signal to say something changed
191  fGrids.Connect("Modified()", "KVIDGridManager", this, "Modified()");
192  return is_it_ok;
193 }
194 
195 
196 
199 
201 {
202  // Add all grids read the last time ReadAsciiFile() was called to the grid manager
203 
204  TIter it(&fLastReadGrids);KVIDGraph* igr;
205  while( (igr = (KVIDGraph*)it()) )
206  fGrids.Add(igr);
207 }
208 
209 
210 
216 
217 Int_t KVIDGridManager::WriteAsciiFile(const Char_t* filename, const TCollection* selection)
218 {
219  // Write grids in file 'filename'.
220  // If selection=0 (default), write all grids.
221  // If selection!=0, write only grids in list.
222  // Returns number of grids written in file.
223 
224  ofstream gridfile(filename);
225 
226  const TCollection* list_of_grids = (selection ? selection : &fGrids);
227  TIter next(list_of_grids);
228  KVIDGraph* grid = 0;
229  Int_t n_saved = 0;
230  while ((grid = (KVIDGraph*) next())) {
231 
232  grid->WriteToAsciiFile(gridfile);
233  Info("WriteAsciiFile", "%s saved", grid->GetName());
234  n_saved++;
235 
236  }
237 
238  gridfile.close();
239  return n_saved;
240 }
241 
242 
243 
246 
248 {
249  //Return pointer to grid with name "name"
250  return (KVIDGraph*) GetGrids()->FindObjectByName(name);
251 }
252 
253 
254 
257 
259 {
260  //Opens GUI for managing grids
261  if (gROOT->IsBatch()) {
262  Warning("StartViewer", "To launch graphical interface, you should not use ROOT in batch mode");
263  return;
264  }
265  TClass* cl = TClass::GetClass("KVIDGridManagerGUI");
266  cl->New();
267 }
268 
269 
270 
274 
276 {
277  // Replace contents of KVString with a comma-separated list of all
278  // different labels of ID telescopes associated with current list of ID grids.
279 
280  list = "";
281  TIter next(&fGrids);
282  KVIDGraph* grid = 0;
283  KVString lab;
284  while ((grid = (KVIDGraph*) next())) {
285  //cout << "grid=" << grid->GetName() << " label=" << grid->GetIDTelescopeLabel() << endl;
286  lab.Form("/%s/", grid->GetIDTelescopeLabel());
287  if (!list.Contains(lab)) list.Append(lab);
288  }
289  //cout << "list=" << list << endl;
290  list.ReplaceAll("//", ",");
291  list.ReplaceAll("/", "");
292  if (list.EqualTo(",")) list = "";
293 }
294 
295 
296 
300 
302 {
303  // Initialize all grids in ID grid manager's list, i.e. we call the Initialize() method
304  // of every grid/graph.
305  TIter next(&fGrids);
306  KVIDGraph* gr = 0;
307  while ((gr = (KVIDGraph*) next())) gr->Initialize();
308 }
309 
310 
311 
int Int_t
bool Bool_t
char Char_t
constexpr Bool_t kFALSE
constexpr Bool_t kTRUE
const char Option_t
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char filename
char name[80]
#define gROOT
void Error(const char *method, const char *msgfmt,...) const override
Definition: KVBase.cpp:1666
void Warning(const char *method, const char *msgfmt,...) const override
Definition: KVBase.cpp:1653
Base class for particle identification in a 2D map.
Definition: KVIDGraph.h:32
virtual void WriteToAsciiFile(std::ofstream &gridfile)
Definition: KVIDGraph.cpp:442
virtual void ReadFromAsciiFile(std::ifstream &gridfile)
Definition: KVIDGraph.cpp:601
const Char_t * GetName() const override
Definition: KVIDGraph.cpp:1344
const Char_t * GetIDTelescopeLabel() const
Definition: KVIDGraph.h:456
virtual void SetOnlyZId(Bool_t yes=kTRUE)
Definition: KVIDGraph.cpp:1508
Handles a stock of identification grids to be used by one or more identification telescopes.
void Clear(Option_t *opt="") override
Delete all grids and empty list, ready to start anew.
void GetListOfIDTelescopeLabels(KVString &)
void DeleteGrid(KVIDGraph *, Bool_t update=kTRUE)
KVSeqCollection fGrids
collection of all ID graphs handled by manager
KVIDGridManager(Bool_t owns_grids=kFALSE)
Int_t WriteAsciiFile(const Char_t *filename, const TCollection *selection=0)
KVSeqCollection * GetGrids()
void AddLastReadGrids()
Add all grids read the last time ReadAsciiFile() was called to the grid manager.
Bool_t ReadAsciiFile(const Char_t *filename)
void StartViewer() const
Opens GUI for managing grids.
TList fLastReadGrids
list of grids created by last call to ReadAsciiFile
KVIDGraph * GetGrid(const Char_t *name)
Return pointer to grid with name "name".
void Initialize(Option_t *="")
void AddGrid(KVIDGraph *)
Add a grid to the collection.
virtual void SendModifiedSignals(Bool_t yes=kTRUE)
TObject * Remove(TObject *obj) override
Remove object from list.
void Add(TObject *obj) override
void Clear(Option_t *option="") override
void SetOwner(Bool_t enable=kTRUE) override
virtual TObject * FindObjectByName(const Char_t *name) const
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:73
void * New(ENewType defConstructor=kClassNew, Bool_t quiet=kFALSE) const
static TClass * GetClass(Bool_t load=kTRUE, Bool_t silent=kFALSE)
void Clear(Option_t *option="") override
void Add(TObject *obj) override
virtual void Fatal(const char *method, const char *msgfmt,...) const
virtual void Info(const char *method, const char *msgfmt,...) const
Bool_t EqualTo(const char *cs, ECaseCompare cmp=kExact) const
TString & Append(char c, Ssiz_t rep=1)
void Form(const char *fmt,...)
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
TString & ReplaceAll(const char *s1, const char *s2)
TGraphErrors * gr
void update(const LAYERDATA &prevLayerData, LAYERDATA &currLayerData, double factorWeightDecay, EnumRegularization regularization)
ClassImp(TPyArg)