KaliVeda  1.12/06
Heavy-Ion Analysis Toolkit
KVParticle.cpp
Go to the documentation of this file.
1 /***************************************************************************
2 $Id: KVParticle.cpp,v 1.50 2009/04/28 08:59:05 franklan Exp $
3  kvparticle.cpp - description
4  -------------------
5  begin : Sun May 19 2002
6  copyright : (C) 2002 by J.D. Frankland
7  email : frankland@ganil.fr
8  ***************************************************************************/
9 
10 /***************************************************************************
11  * *
12  * This program is free software; you can redistribute it and/or modify *
13  * it under the terms of the GNU General Public License as published by *
14  * the Free Software Foundation; either version 2 of the License, or *
15  * (at your option) any later version. *
16  * *
17  ***************************************************************************/
18 #include "TMath.h"
19 #include "Riostream.h"
20 #include "KVPosition.h"
21 #include "KVParticle.h"
22 #include "KVList.h"
23 #include "TRotation.h"
24 #include "TLorentzRotation.h"
25 #include "TObjString.h"
26 #include "TClass.h"
27 #include "KVKinematicalFrame.h"
28 
30 
31 using namespace std;
32 
34 
35 
36 
37 
38 
40 KVParticle::KVParticle() : fParameters("ParticleParameters", "Parameters associated with a particle in an event")
41 {
42  init();
43 }
44 
45 
46 
49 
51 {
52  //default initialisation
53  fE0 = 0;
54  SetFrameName("");
55  fGroups.SetOwner(kTRUE);
56 }
57 
58 
59 
62 
64 {
65  //copy ctor
66  init();
67 #if ROOT_VERSION_CODE >= ROOT_VERSION(3,4,0)
68  obj.Copy(*this);
69 #else
70  ((KVParticle&) obj).Copy(*this);
71 #endif
72 }
73 
74 
75 
78 
79 KVParticle::KVParticle(Double_t m, TVector3& p) : fParameters("ParticleParameters", "Parameters associated with a particle in an event")
80 {
81  //create particle with given mass and momentum vector
82  init();
83  SetMass(m);
84  SetMomentum(p);
85 }
86 
87 
88 
91 
92 KVParticle::KVParticle(Double_t m, Double_t px, Double_t py, Double_t pz) : fParameters("ParticleParameters", "Parameters associated with a particle in an event")
93 {
94  //create particle with given mass and momentum vector
95  init();
96  SetMass(m);
97  SetMomentum(px, py, pz);
98 }
99 
100 
101 
104 
105 KVParticle::~KVParticle()
106 {
107  //Info("~KVParticle","%p",this);
108  Clear();
109 }
110 
111 
112 
116 
118 {
119  //Static function.
120  //Returns speed of light in cm/ns units.
121  return kSpeedOfLight;
122 }
123 
124 
125 
135 
137  Double_t thmax, Double_t phmin,
138  Double_t phmax, Option_t* opt)
139 {
140  //Give randomly directed momentum to particle with kinetic energy T
141  //Direction will be between (thmin,thmax) [degrees] limits in polar angle,
142  //and (phmin,phmax) [degrees] limits in azimuthal angle.
143  //
144  //If opt = "" or "isotropic" (default) : direction is isotropically distributed over the solid angle
145  //If opt = "random" : direction is randomly distributed over solid angle
146  //
147  //Based on KVPosition::GetRandomDirection().
148 
149  Double_t p = (T + M()) * (T + M()) - M2();
150  if (p > 0.)
151  p = (TMath::Sqrt(p)); // calculate momentum
152  else
153  p = 0.;
154 
155  TVector3 dir;
156  KVPosition pos(thmin, thmax, phmin, phmax);
157  dir = pos.GetRandomDirection(opt); // get isotropic unit vector dir
158  if (p && dir.Mag())
159  dir.SetMag(p); // set magnitude of vector to momentum required
160  SetMomentum(dir); // set momentum 4-vector
161 }
162 
163 
164 
168 
170 {
171  //set momentum with kinetic energy t and unit direction vector d
172  //(d is normalised first in case it is not a unit vector)
173 
174  Double_t p = (T + M()) * (T + M()) - M2();
175  TVector3 pdir;
176  TVector3 unit_dir = dir.Unit();
177  if (p > 0.) {
178  p = (TMath::Sqrt(p));
179  pdir = p * unit_dir;
180  }
181  SetMomentum(pdir);
182 };
183 
184 
185 
188 
190 {
191  // recursive print out of all defined kinematical frames
192 
193  fmt += "\t";
194  if (fBoosted.GetEntries()) {
195  TIter next(&fBoosted);
196  KVKinematicalFrame* frame;
197  while ((frame = (KVKinematicalFrame*)next())) {
198  cout << fmt << " " << frame->GetName() << ": ";
199  KVParticle* part = frame->GetParticle();
200  cout << " Theta=" << part->GetTheta() << " Phi=" << part->GetPhi()
201  << " KE=" << part->GetKE() << " Vpar=" << part->GetVpar() << endl;
202  part->print_frames(fmt);
203  }
204  }
205 }
206 
207 
208 
211 
213 {
214  // print out characteristics of particle
215 
216  cout << "KVParticle mass=" << M() <<
217  " Theta=" << GetTheta() << " Phi=" << GetPhi()
218  << " KE=" << GetKE() << " Vpar=" << GetVpar() << endl;
219  print_frames();
220  GetParameters()->Print();
221 }
222 
223 
224 
229 
231 {
232  //Change particle KE, keeping momentum direction constant
233  //If momentum is zero (i.e. no direction defined) the particle will be given
234  //a velocity in the positive z-direction
235 
236  if (ecin > 0.) {
237  Double_t et = M() + ecin; // new total energy = KE + m
238  Double_t pmod = 0;
239  if (et * et > M2()) {
240  pmod = TMath::Sqrt(et * et - M2()); // p**2 = E**2 - m**2
241  TVector3 newp(Px(), Py(), Pz());
242  if (pmod && newp.Mag()) {
243  newp.SetMag(pmod);
244  } else {
245  newp.SetXYZ(0, 0, 1);
246  newp.SetMag(pmod);
247  }
248  SetMomentum(newp);
249  } else {
250  SetMomentum(0., 0., 0.);
251  }
252  } else
253  SetMomentum(0., 0., 0.);
254 }
255 
256 
257 #if ROOT_VERSION_CODE >= ROOT_VERSION(3,4,0)
258 
265 
266 void KVParticle::Copy(TObject& obj) const
267 #else
268 void KVParticle::Copy(TObject& obj)
269 #endif
270 {
271  // Copy this to obj
272  // Particle kinematics are copied using operator=(const KVParticle&)
273  // List of particle's groups is copied
274  // The particle's name is copied
275  // The list of parameters associated with the particle is copied
276 
277  ((KVParticle&) obj) = *this;
278  ((KVParticle&) obj).SetGroups(this->GetGroups());
279  ((KVParticle&) obj).SetName(this->GetName());
280  fParameters.Copy(((KVParticle&) obj).fParameters);
281 }
282 
283 
284 
285 
288 
290 {
291  //Reset particle properties i.e. before creating/reading a new event
292 
293  SetXYZM(0., 0., 0., 0.);
294  if (fE0) {
295  delete fE0;
296  fE0 = 0;
297  }
298  ResetIsOK(); //in case IsOK() status was set "by hand" in previous event
300  fParameters.Clear();
301  fGroups.Clear();
302  fBoosted.Delete();
303 }
304 
305 
306 
311 
313 {
314  //Determine whether this particle is considered "good" or not for analysis,
315  //depending on a previous call to SetIsOK(Bool_t flag).
316  //If SetIsOK has not been called, IsOK returns kTRUE by default.
317 
318  if (TestBit(kIsOKSet)) { //status set by SetIsOK()
319  return TestBit(kIsOK);
320  }
321  //Default if SetIsOK not used: accept particle
322  return kTRUE;
323 }
324 
325 
326 
327 
332 
334 {
335  //Set acceptation/rejection status for the particle.
336  //
337  //In order to 'forget' this status (accept all particles) use ResetIsOK()
338 
339  SetBit(kIsOK, flag);
340  SetBit(kIsOKSet);
341 }
342 
343 
344 
346 
347 void KVParticle::ls(Option_t* option) const
348 {
349  std::cout << option << GetName() << ":" << GetFrameName() << ":" << this << "\n";
350  if (fBoosted.GetEntries()) {
351  TString nopt = option;
352  nopt += " ";
353  TIter next(&fBoosted);
355  while ((p = (KVKinematicalFrame*)next())) p->GetParticle()->ls(nopt.Data());
356  }
357 }
358 
359 
360 
361 
364 
366 {
367  //KVParticle assignment operator.
368 
370  if (rhs.GetPInitial()) SetE0(rhs.GetPInitial());
371  return *this;
372 }
373 
374 
375 
376 
384 
386 {
387  //Used for simulated particles after passage through some absorber/detector.
388  //The passage of a particle through the different absorbers modifies its
389  //kinetic energies, indeed the particle may be stopped in the detector.
390  //Calling this method will reset the particle's momentum to its
391  //initial value i.e. before it entered the first absorber.
392  //Particles which have not encountered any absorbers/detectors are left as they are.
393 
394  if (IsDetected())
396 }
397 
398 
399 
400 
403 
404 void KVParticle::SetName(const Char_t* nom)
405 {
406  //Set Name of the particle
407  fName.Form("%s", nom);
408 
409 }
410 
411 
412 
415 
417 {
418  // return the field fName
419  return fName.Data();
420 }
421 
422 
423 
427 
429 {
430  // Dummy implementation of AddGroup(const Char_t* groupname, KVParticleCondition*)
431  // Does nothing. Real implementation is in KVNucleus::AddGroup_Withcondition.
432  Warning("AddGroup_Withcondition", "DUUUUUUUUUUUUUMYYYYYYY do nothing");
433 };
434 
435 
436 
447 
448 void KVParticle::AddGroup_Sanscondition(const Char_t* groupname, const Char_t* from)
449 {
450  // Implementation of AddGroup_Sansconditioncon(st Char_t*, const Char_t*)
451  // Can be overridden in child classes [instead of AddGroup(const Char_t*, const Char_t*),
452  // which cannot]
453  // if this method is overridde in child class
454  // the line
455  // if (!fGroups) CreateGroups();
456  // has to be included
457  // The group will be added to all 'particles' representing different kinematical frames
458  // for this particle
459 
460  TString sfrom(from);
461  sfrom.ToUpper();
462  TString sgroupname(groupname);
463  sgroupname.ToUpper();
464 
465  if (BelongsToGroup(sfrom.Data()) && !BelongsToGroup(sgroupname.Data())) {
466  fGroups.Add(new TObjString(sgroupname.Data()));
467  if (fBoosted.GetEntries()) {
468  // recursively add to all boosted particles
469  TIter it(&fBoosted);
471  while ((f = (KVKinematicalFrame*)it())) {
472  f->GetParticle()->AddGroup(sgroupname);
473  }
474  }
475  }
476 }
477 
478 
479 
484 
485 void KVParticle::AddGroup(const Char_t* groupname, const Char_t* from)
486 {
487  // Associate this particle with the given named group.
488  // Optional argument "from" allows to put a condition on the already stored
489  // group list, is set to "" by default
490 
491  AddGroup_Sanscondition(groupname, from);
492 }
493 
494 
495 
496 
497 
504 
505 void KVParticle::AddGroup(const Char_t* groupname, KVParticleCondition* cond)
506 {
507  //define and store a group name from a condition on the particle
508  //
509  // Apply the method to all particles stored in fBoosted
510  // SetParticleClassName has to be set before using this method if you use
511  // in the KVParticleCondistion a specific method of a derived KVNucleus class
512 
513  AddGroup_Withcondition(groupname, cond);
514 }
515 
516 
517 
521 
523 {
524  //Define for the particle a new list of groups
525  //if there is an existing list, it's deleted
526  fGroups.Clear();
527  AddGroups(un);
528 }
529 
530 
531 
534 
536 {
537  //list of groups added to the current one
538  TObjString* os = 0;
539  TIter no(un);
540  while ((os = (TObjString*)no.Next())) {
541  AddGroup(os->GetName());
542  }
543 
544 }
545 
546 
547 
552 
554 {
555  // Returns the total number of defined kinematical frames for this particle.
556  // This includes all frames which are defined with respect to other frames
557  // (i.e. we count recursively the contents of all fBoosted lists)
558 
559  if (!fBoosted.GetEntries()) return 0;
560  TIter it(&fBoosted);
562  Int_t nf = 0;
563  while ((p = (KVKinematicalFrame*)it())) {
564  nf += (1 + p->GetParticle()->GetNumberOfDefinedFrames());
565  }
566  return nf;
567 }
568 
569 
572 
574 {
575  //return the number of defined groups for the particle
576  return fGroups.GetEntries();
577 }
578 
579 
580 
583 
585 {
586  //return the KVUniqueNameList pointeur where list of groups are stored
587  return (KVUniqueNameList*)&fGroups;
588 }
589 
590 
591 
596 
597 Bool_t KVParticle::BelongsToGroup(const Char_t* groupname) const
598 {
599  //Check if particle belong to a given group
600  //return kTRUE if groupname="".
601  //return kFALSE if no group has be defined
602 
603  TString sgroupname(groupname);
604  sgroupname.ToUpper();
605  //Important for KVEvent::GetNextParticle()
606  if (sgroupname.IsNull()) return kTRUE;
607  //retourne kFALSE si aucun groupe n'est defini
608  if (!fGroups.GetEntries()) return kFALSE;
609  if (fGroups.FindObject(sgroupname.Data())) return kTRUE;
610  return kFALSE;
611 }
612 
613 
614 
618 
619 void KVParticle::RemoveGroup(const Char_t* groupname)
620 {
621  // Remove group from list of groups
622  // Apply the method to all particles stored in fBoosted
623  if (!fGroups.GetEntries()) return;
624  TString sgroupname(groupname);
625  sgroupname.ToUpper();
626 
627  TObjString* os = 0;
628  if ((os = (TObjString*)fGroups.FindObject(sgroupname.Data()))) {
629  delete fGroups.Remove(os);
630  if (fBoosted.GetEntries()) {
631  TIter it(&fBoosted);
633  while ((f = (KVKinematicalFrame*)it())) {
634  f->GetParticle()->RemoveGroup(sgroupname);
635  }
636  }
637  }
638 }
639 
640 
641 
645 
647 {
648  //Remove all groups
649  // Apply the method to all particles stored in fBoosted
650  fGroups.Clear();
651  if (fBoosted.GetEntries()) {
652  TIter it(&fBoosted);
654  while ((f = (KVKinematicalFrame*)it())) {
655  f->GetParticle()->RemoveAllGroups();
656  }
657  }
658 }
659 
660 
661 
662 
665 
666 void KVParticle::ListGroups(void) const
667 {
668  //List all stored groups
669  if (!fGroups.GetEntries()) {
670  cout << "Cette particle n appartient a aucun groupe" << endl;
671  return;
672  } else {
673  cout << "--------------------------------------------------" << endl;
674  cout << "Liste des groupes auxquels la particule appartient" << endl;
675  }
676  TObjString* os = 0;
677  TIter no(GetGroups());
678  while ((os = (TObjString*)no.Next())) cout << os->GetName() << endl;
679  cout << "--------------------------------------------------" << endl;
680 }
681 
682 
683 
687 
689 {
690  // Use this method to obtain 'on-the-fly' some information on particle kinematics
691  // in a different reference frame. The default kinematics of the particle are unchanged.
692 
693  KVParticle p;
694  p.Set4Mom(*this);
695  KVKinematicalFrame(&p, t);
696  return p;
697 }
698 
699 
700 
706 
708 {
709  // Permanently modify kinematics of particle according to the given transformation.
710  // You can optionally set the name of this new default kinematics.
711  // NB the current kinematics will be lost. If you want to keep it after changing the
712  // default kinematics, define the new frame with SetFrame and then use ChangeDefaultFrame.
713 
714  KVKinematicalFrame(this, t);
715  if (name != "") SetFrameName(name);
716 }
717 
718 
719 
725 
726 void KVParticle::ChangeDefaultFrame(const Char_t* newdef, const Char_t* defname)
727 {
728  // Make existing reference frame 'newdef' the new default frame for particle kinematics.
729  // The current default frame will then be accessible from the list of frames
730  // using its name (if set with SetFrameName). You can change/set the name of the previous
731  // default frame with 'defname'
732 
733  TString _defname(defname);
734  if (_defname == "") _defname = GetFrameName();
735  // get list of all parents of new default
736  TList parents;
737  TString ff = newdef;
739  do {
740  f = get_parent_frame(ff);
741  if (f) {
742  parents.Add(f);
743  ff = f->GetName();
744  }
745  } while (f);
746  // modify tree structure
747  TIter it(&parents);
748  KVKinematicalFrame* newdframe = get_frame(newdef);
749  KVKinematicalFrame* save_newdef = newdframe;
750  KVFrameTransform trans, next_trans = newdframe->GetTransform();
751  while ((f = (KVKinematicalFrame*)it())) {
752  f->GetParticle()->GetListOfFrames()->Remove(newdframe);
753  trans = next_trans;
754  next_trans = f->GetTransform();
755  newdframe->GetParticle()->GetListOfFrames()->Add(f);
756  f->SetTransform(trans.Inverse());
757  newdframe = f;
758  }
759  GetListOfFrames()->Remove(newdframe);
760  // copy current default kinematics particle momentum/energy
761  TLorentzVector old_def_p(*this);
762  // set momentum/energy of new default kinematics
763  Set4Mom(*(save_newdef->GetParticle()));
764  save_newdef->GetParticle()->Set4Mom(old_def_p);
765  save_newdef->SetTransform(next_trans.Inverse());
766  save_newdef->SetName(_defname);
767  newdframe->GetParticle()->GetListOfFrames()->Add(save_newdef);
768  // copy frame list from new default frame
769  TList frame_list;
770  frame_list.AddAll(save_newdef->GetParticle()->GetListOfFrames());
771  save_newdef->GetParticle()->GetListOfFrames()->Clear("nodelete");
772  save_newdef->GetParticle()->GetListOfFrames()->AddAll(&fBoosted);
773  fBoosted.Clear("nodelete");
774  fBoosted.AddAll(&frame_list);
775  SetFrameName(newdef);
776 }
777 
778 
779 
839 
840 void KVParticle::SetFrame(const Char_t* frame, const KVFrameTransform& ft)
841 {
842  //Define a Lorentz-boosted and/or rotated frame in which to calculate this particle's momentum and energy.
843  //
844  //The new frame will have the name given in the string "frame", which can then be used to
845  //access the kinematics of the particle in different frames using GetFrame() (frame names are case-insensitive).
846  //Calling this method with the name of an existing frame will update the kinematics of the particle
847  //in that frame using the given transform (note that kinematics in all frames defined as 'subframes' of
848  //this frame will also be updated).
849  //
850  //USING BOOSTS
851  //The boost velocity vector is that of the boosted frame with respect to the original frame of the particles in the event.
852  //The velocity vector can be given either in cm/ns units (default) or in units of 'c' (beta=kTRUE).
853  //
854  //E.g. to define a frame moving at 0.1c in the +ve z-direction with respect to the original
855  //event frame:
856  //
857  // (...supposing a valid pointer KVParticle* my_part...)
858  // TVector3 vframe(0,0,0.1);
859  // my_part->SetFrame("my_frame", KVFrameTransform(vframe, kTRUE));
860  //
861  //or with velocity in cm/ns units (default):
862  // TVector3 vframe(0,0,3);
863  // my_part->SetFrame("my_frame", KVFrameTransform(vframe));
864  // OR my_part->SetFrame("my_frame", vframe);
865  //
866  //USING ROTATIONS
867  //According to the conventions adopted for the TRotation and TLorentzRotation classes,
868  //we actually use the inverse of the TLorentzRotation to make the transformation,
869  //to get the coordinates corresponding to a rotated coordinate system, not the coordinates
870  //of a rotated vector in the same coordinate system
871  //=> you do not need to invert the transformation matrix
872  //
873  //E.g. if you want to define a new frame whose coordinate axes are rotated with respect
874  //to the original axes, you can set up a TRotation like so:
875  //
876  // TRotation rot;
877  // TVector3 newX, newY, newZ; // the new coordinate axes
878  // rot.RotateAxes(newX, newY, newZ);
879  //
880  //If you are using one of the two global variables which calculate the event tensor
881  //(KVTensP and KVTensPCM) you can obtain the transformation to the tensor frame
882  //using:
883  //
884  // TRotation rot;
885  // KVTensP* tens_gv;// pointer to tensor global variable
886  // tens_gv->GetTensor()->GetRotation(rot);// see KVTenseur3::GetRotation
887  //
888  //Then the new frame can be defined by
889  // my_part->SetFrame("my_frame", KVFrameTransform(rot));
890  // OR my_part->SetFrame("my_frame", rot);
891  //
892  //USING COMBINED BOOST AND ROTATION
893  //You can define a frame using both a boost and a rotation like so:
894  // my_part->SetFrame("my_frame", KVFrameTransform(vframe,rot,kTRUE));
895  //
896  //ACCESSING KINEMATICS IN NEW FRAMES
897  //In order to access the kinematics of the particle in the new frame:
898  //
899  // my_part->GetFrame("my_frame")->GetTransverseEnergy();// transverse energy in "my_frame"
900 
901  if (!strcmp(frame, "")) return;
902 
903  KVKinematicalFrame* tmp = get_frame(frame);
904  if (!tmp) {
905  //if this frame has not already been defined, create a new one
906  tmp = new KVKinematicalFrame(frame, this, ft);
907  fBoosted.Add(tmp);
908  } else
909  tmp->ApplyTransform(this, ft);
910 }
911 
912 
913 
951 
952 KVParticle const* KVParticle::GetFrame(const Char_t* frame, Bool_t warn_and_return_null_if_unknown) const
953 {
954  // Return the momentum of the particle in the Lorentz-boosted frame corresponding to the name
955  // "frame" given as argument (see SetFrame() for definition of different frames).
956  // If the default frame name has been set (see KVEvent::SetFrameName) and 'frame' is the
957  // name of this default frame (KVParticle::fFrameName), we return the address of the particle
958  // itself.
959  //
960  // This frame may have been defined by a direct transformation of the original kinematics of the
961  // particle (using SetFrame(newframe,...)) or by a transformation of the kinematics in another
962  // user-defined frame (using SetFrame(newframe,oldframe,...)).
963  //
964  // Note that frames are not "dynamic": if any changes are made to the original particle's kinematics
965  // after definition of a frame, if you want these changes to affect also the other frames you
966  // need to update them by hand by calling KVParticle::UpdateAllFrames().
967  //
968  // Frame names are case insensitive: "CM" or "cm" or "Cm" are all good...
969  //
970  // By default, if no frame with the given name is found, we return nullptr and print a warning.
971  // If `warn_and_return_null_if_unknown=kFALSE`, we return the address of the particle itself,
972  // i.e. the original/default kinematics.
973  // [Note that this is an inversion of the previous default behaviour]
974  //
975  // Note that the properties of the particle returned by this method can not be modified:
976  // this is deliberate, as any modifications e.g. to kinematics will have no effect
977  // in any other frames.
978  //
979  // The returned pointer corresponds to a "pseudoparticle" in the desired frame,
980  // therefore you can use any KVParticle method in order to access the kinematics of the
981  // particle in the boosted frame, e.g.
982  //
983  //~~~~~~~~~~~~~~~~~~~~
984  // (...supposing a valid pointer KVParticle* my_part...)
985  // my_part->GetFrame("cm_frame")->GetVpar();// //el velocity in "cm_frame"
986  // my_part->GetFrame("QP_frame")->GetTheta();// polar angle in "QP_frame"
987  // etc. etc.
988  //~~~~~~~~~~~~~~~~~~~~
989  //
990 
991  if (!fFrameName.CompareTo(frame, TString::kIgnoreCase)) return (KVParticle const*)this;
992  KVKinematicalFrame* f = get_frame(frame);
993  return f ? (KVParticle const*)f->GetParticle() :
994  (warn_and_return_null_if_unknown ?
995  Warning("GetFrame(const Char_t*)", "No frame \"%s\" defined for particle. 0x0 returned.",
996 #ifndef WITH_CPP11
997  frame), (KVParticle*)nullptr
998 #else
999  frame), nullptr
1000 #endif
1001  : this);
1002 }
1003 
1004 
1005 
1009 
1011 {
1012  // Call this method to update particle kinematics in all defined frames if you change
1013  // the kinematics of the particle in its original/default frame.
1014 
1015  if (fBoosted.GetEntries()) {
1016  TIter it(&fBoosted);
1018  while ((f = (KVKinematicalFrame*)it())) {
1019  f->ReapplyTransform(this);
1020  // recursively apply to all subframes
1021  f->GetParticle()->UpdateAllFrames();
1022  }
1023  }
1024 }
1025 
1026 
1027 
1032 
1034 {
1035  // PRIVATE method for internal use only
1036  // This allows to modify the returned frame, i.e. in order to define
1037  // new frames in SetFrame(newframe,oldframe,...)
1038 
1039  if (!fBoosted.GetEntries() || !strcmp(frame, "")) {
1040  // no frames defined or no frame name given
1041  return nullptr;
1042  }
1043  TString _frame(frame);
1044  TIter it(&fBoosted);
1045  KVKinematicalFrame* p(nullptr), *f(nullptr);
1046  while ((p = f = (KVKinematicalFrame*)it())) {
1047  if (!_frame.CompareTo(p->GetName(), TString::kIgnoreCase)) break;
1048  // look for subframe
1049  if ((f = p->GetParticle()->get_frame(_frame))) break;
1050  }
1051  return f;
1052 }
1053 
1054 
1055 
1059 
1061 {
1062  // PRIVATE method for internal use only
1063  // Returns pointer to parent frame of 'f'
1064  if (!fBoosted.GetEntries() || !strcmp(f, "")) {
1065  // no frames defined or no frame name given
1066  return nullptr;
1067  }
1068  TString _frame(f);
1069  TIter it(&fBoosted);
1070  KVKinematicalFrame* p(nullptr), *r(nullptr);
1071  while ((p = (KVKinematicalFrame*)it())) {
1072  if (!_frame.CompareTo(p->GetName(), TString::kIgnoreCase)) return F;
1073  // look for subframe
1074  if ((r = p->GetParticle()->get_parent_frame(_frame, p))) return r;
1075  }
1076  return nullptr;
1077 }
1078 
1079 
1080 
1081 
1082 
1087 
1088 void KVParticle::SetFrame(const Char_t* newframe, const Char_t* oldframe, const KVFrameTransform& ft)
1089 {
1090  // Define new kinematical frame by transformation from existing frame
1091  // See SetFrame(const Char_t*,const KVFrameTransform&) for details on
1092  // defining kinematically-transformed frames.
1093 
1094  KVKinematicalFrame* f = get_frame(oldframe);
1095  if (!f) {
1096  Error("SetFrame(newframe,oldframe)", "oldframe=%s does not exist!", oldframe);
1097  return;
1098  }
1099  f->GetParticle()->SetFrame(newframe, ft);
1100 }
1101 
1102 
1103 
1104 
1107 
1109 {
1110  //returns velocity vector in cm/ns units
1111  TVector3 beta;
1112  if (E()) {
1113  beta = GetMomentum() * (1. / E());
1114  } else {
1115  beta.SetXYZ(0, 0, 0);
1116  }
1117  return (kSpeedOfLight * beta);
1118 }
1119 
1120 
1121 
1122 
1126 
1128 {
1129  //returns transverse velocity in cm/ns units
1130  //sign is +ve if py>=0, -ve if py<0
1131  return (GetV().y() >= 0.0 ? GetV().Perp() : -(GetV().Perp()));
1132 }
1133 
1134 
1135 
1136 
1142 
1144  Option_t* opt)
1145 {
1146  // Set Momentum components (in MeV/c)
1147  // if option is "cart" or "cartesian" we give cartesian components (x,y,z)
1148  // if option is "spher" or "spherical" we give components (rho,theta,phi) in spherical system
1149  // with theta, phi in DEGREES
1150  if (!strcmp("spher", opt) || !strcmp("spherical", opt)) {
1151  TVector3 pvec(0., 0., 1.);
1152  pvec.SetMag(px);
1153  pvec.SetTheta(TMath::Pi() * py / 180.);
1154  pvec.SetPhi(TMath::Pi() * pz / 180.);
1155  SetVectM(pvec, M());
1156  } else {
1157  if (strcmp("cart", opt) && strcmp("cartesian", opt)) {
1158  Warning("SetMomentum(Double_t,Double_t,Double_t,Option_t*)",
1159  "Unkown coordinate system\n known system are :\n\t\"cartesian\" or \"cart\" (default)\n\t\"spherical\" or \"spher\"");
1160  Warning("SetMomentum(Double_t,Double_t,Double_t,Option_t*)",
1161  "default used.");
1162  }
1163  TVector3 pvec(px, py, pz);
1164  SetVectM(pvec, M());
1165  }
1166 }
1167 
1168 
1169 
1172 
1174 {
1175  // Set velocity of particle (in cm/ns units)
1176  Double_t gamma = 1. / kSpeedOfLight / sqrt(1 - (vel.Mag2() / pow(kSpeedOfLight, 2)));
1177  TVector3 p = GetMass() * gamma * vel;
1178  SetMomentum(p);
1179 }
1180 
1181 
1182 
1187 
1188 void KVParticle::Streamer(TBuffer& R__b)
1189 {
1190  // Stream an object of class KVParticle.
1191  // When reading: If parameter "frameName" is set, use it to set non-persistent
1192  // fFrameName member (used by GetFrame)
1193 
1194  if (R__b.IsReading()) {
1195  R__b.ReadClassBuffer(KVParticle::Class(), this);
1196  if (GetParameters()->HasStringParameter("frameName")) fFrameName = GetParameters()->GetStringValue("frameName");
1197  } else {
1198  R__b.WriteClassBuffer(KVParticle::Class(), this);
1199  }
1200 }
1201 
1202 
int Int_t
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
ROOT::R::TRInterface & r
#define f(i)
char Char_t
const Bool_t kFALSE
bool Bool_t
double Double_t
const Bool_t kTRUE
const char Option_t
double pow(double, double)
Utility class for kinematical transformations of KVParticle class.
Kinematical representation of a particle in different reference frames.
virtual void Print(Option_t *opt="") const
virtual void Clear(Option_t *opt="")
const Char_t * GetStringValue(const Char_t *name) const
Handles particle selection criteria for data analysis classes ,.
Base class for relativistic kinematics of massive particles.
Definition: KVParticle.h:398
void AddGroups(KVUniqueNameList *un)
list of groups added to the current one
Definition: KVParticle.cpp:535
void SetIsOK(Bool_t flag=kTRUE)
Definition: KVParticle.cpp:333
TVector3 * GetPInitial() const
Definition: KVParticle.h:687
virtual void SetMass(Double_t m)
Definition: KVParticle.h:527
void RemoveGroup(const Char_t *groupname)
Definition: KVParticle.cpp:619
void AddGroup(const Char_t *groupname, const Char_t *from="")
Definition: KVParticle.cpp:485
KVUniqueNameList * GetGroups() const
return the KVUniqueNameList pointeur where list of groups are stored
Definition: KVParticle.cpp:584
void UpdateAllFrames()
TVector3 GetMomentum() const
Definition: KVParticle.h:565
void ResetIsOK()
Definition: KVParticle.h:665
void RemoveAllGroups()
Definition: KVParticle.cpp:646
KVUniqueNameList fGroups
list of TObjString for manage different group name
Definition: KVParticle.h:408
KVNameValueList * GetParameters() const
Definition: KVParticle.h:735
virtual Bool_t IsOK()
Definition: KVParticle.cpp:312
Double_t GetTheta() const
Definition: KVParticle.h:641
KVParticle & operator=(const KVParticle &rhs)
KVParticle assignment operator.
Definition: KVParticle.cpp:365
void SetVectM(const TVector3 &spatial, Double_t mass)
Definition: KVParticle.h:416
const Char_t * GetFrameName(void) const
Definition: KVParticle.h:722
void ResetEnergy()
Definition: KVParticle.cpp:385
void SetMomentum(const TVector3 &v)
Definition: KVParticle.h:535
virtual void AddGroup_Withcondition(const Char_t *, KVParticleCondition *)
Definition: KVParticle.cpp:428
void SetVelocity(const TVector3 &)
Set velocity of particle (in cm/ns units)
const Char_t * GetName() const
return the field fName
Definition: KVParticle.cpp:416
Bool_t IsDetected()
Definition: KVParticle.h:695
static Double_t C()
Definition: KVParticle.cpp:117
void SetKE(Double_t ecin)
Definition: KVParticle.cpp:230
virtual void Clear(Option_t *opt="")
Reset particle properties i.e. before creating/reading a new event.
Definition: KVParticle.cpp:289
KVList fBoosted
list of momenta of the particle in different Lorentz-boosted frames
Definition: KVParticle.h:407
void ls(Option_t *option="") const
Definition: KVParticle.cpp:347
void init()
default initialisation
Definition: KVParticle.cpp:50
Double_t GetPhi() const
Definition: KVParticle.h:649
KVNameValueList fParameters
a general-purpose list of parameters associated with this particle
Definition: KVParticle.h:484
virtual void Copy(TObject &) const
Definition: KVParticle.cpp:266
KVList * GetListOfFrames(void)
Definition: KVParticle.h:670
void Set4Mom(const TLorentzVector &p)
Definition: KVParticle.h:550
KVKinematicalFrame * get_parent_frame(const Char_t *, KVKinematicalFrame *F=nullptr) const
virtual void Print(Option_t *t="") const
print out characteristics of particle
Definition: KVParticle.cpp:212
KVParticle InFrame(const KVFrameTransform &)
Definition: KVParticle.cpp:688
void SetE0(TVector3 *e=0)
Definition: KVParticle.h:676
TString fFrameName
non-persistent frame name field, sets when calling SetFrame method
Definition: KVParticle.h:406
Int_t GetNumberOfDefinedFrames(void)
Definition: KVParticle.cpp:553
void SetFrameName(const Char_t *framename)
Definition: KVParticle.h:726
static Double_t kSpeedOfLight
speed of light in cm/ns
Definition: KVParticle.h:409
Double_t GetKE() const
Definition: KVParticle.h:575
Double_t GetVpar() const
Definition: KVParticle.h:636
void SetXYZM(Double_t x, Double_t y, Double_t z, Double_t m)
Definition: KVParticle.h:472
void SetName(const Char_t *nom)
Set Name of the particle.
Definition: KVParticle.cpp:404
void SetFrame(const Char_t *frame, const KVFrameTransform &)
Definition: KVParticle.cpp:840
KVParticle const * GetFrame(const Char_t *frame, Bool_t warn_and_return_null_if_unknown=kTRUE) const
Definition: KVParticle.cpp:952
void ChangeFrame(const KVFrameTransform &, const KVString &="")
Definition: KVParticle.cpp:707
void SetGroups(KVUniqueNameList *un)
Definition: KVParticle.cpp:522
Bool_t BelongsToGroup(const Char_t *groupname) const
Definition: KVParticle.cpp:597
TVector3 * fE0
the momentum of the particle before it is slowed/stopped by an absorber
Definition: KVParticle.h:483
void ChangeDefaultFrame(const Char_t *, const Char_t *defname="")
Definition: KVParticle.cpp:726
void SetRandomMomentum(Double_t T, Double_t thmin, Double_t thmax, Double_t phmin, Double_t phmax, Option_t *opt="isotropic")
Definition: KVParticle.cpp:136
Int_t GetNumberOfDefinedGroups(void)
return the number of defined groups for the particle
Definition: KVParticle.cpp:573
TVector3 GetV() const
Definition: KVParticle.h:632
TString fName
non-persistent name field - Is useful
Definition: KVParticle.h:405
KVKinematicalFrame * get_frame(const Char_t *) const
Double_t GetMass() const
Definition: KVParticle.h:531
virtual void AddGroup_Sanscondition(const Char_t *groupname, const Char_t *from="")
Definition: KVParticle.cpp:448
void ListGroups(void) const
List all stored groups.
Definition: KVParticle.cpp:666
TVector3 GetVelocity() const
returns velocity vector in cm/ns units
Double_t GetVperp() const
void print_frames(TString fmt="") const
recursive print out of all defined kinematical frames
Definition: KVParticle.cpp:189
Base class used for handling geometry in a multidetector array.
Definition: KVPosition.h:90
virtual TVector3 GetRandomDirection(Option_t *t="isotropic")
Definition: KVPosition.cpp:242
virtual void Clear(Option_t *option="")
virtual void Add(TObject *obj)
virtual TObject * Remove(TObject *obj)
Remove object from list.
virtual void Delete(Option_t *option="")
virtual TObject * FindObject(const char *name) const
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:72
Optimised list in which named objects can only be placed once.
virtual void Add(TObject *obj)
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
Bool_t IsReading() const
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
virtual void AddAll(const TCollection *col)
virtual Int_t GetEntries() const
TObject * Next()
virtual void Add(TObject *obj)
TLorentzRotation Inverse() const
TLorentzVector & operator=(const TLorentzVector &)
Double_t Px() const
Double_t M() const
Double_t M2() const
Double_t Pz() const
Double_t Py() const
Double_t Perp() const
Double_t E() const
Double_t T() const
const char * GetName() const
void SetBit(UInt_t f)
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
virtual void Warning(const char *method, const char *msgfmt,...) const
virtual void Error(const char *method, const char *msgfmt,...) const
void ResetBit(UInt_t f)
int CompareTo(const char *cs, ECaseCompare cmp=kExact) const
void ToUpper()
const char * Data() const
Bool_t IsNull() const
void Form(const char *fmt,...)
void SetXYZ(Double_t x, Double_t y, Double_t z)
void SetPhi(Double_t)
Double_t Mag2() const
TVector3 Unit() const
Double_t Mag() const
void SetMag(Double_t)
void SetTheta(Double_t)
RooCmdArg Parameters(const RooArgSet &params)
double beta(double x, double y)
VecExpr< UnaryOp< Sqrt< T >, SVector< T, D >, T >, T, D > sqrt(const SVector< T, D > &rhs)
Double_t y[n]
#define F(x, y, z)
const long double m
Definition: KVUnits.h:70
double gamma(double x)
constexpr Double_t C()
Double_t Sqrt(Double_t x)
constexpr Double_t Pi()