KaliVeda  1.12/06
Heavy-Ion Analysis Toolkit
KVElasticScatterEvent.cpp
Go to the documentation of this file.
1 /*
2 $Id: KVElasticScatterEvent.cpp,v 1.5 2009/01/14 15:35:50 ebonnet Exp $
3 $Revision: 1.5 $
4 $Date: 2009/01/14 15:35:50 $
5 */
6 
7 //Created by KVClassFactory on Thu Dec 11 14:45:29 2008
8 //Author: eric bonnet,,,
9 
10 #include "KVElasticScatterEvent.h"
11 #include "KVPosition.h"
12 #include "KVASMultiDetArray.h"
13 #include "KVASGroup.h"
14 #include "KVDetector.h"
15 #include "KVTelescope.h"
16 #include "TH2F.h"
17 #include "TH1F.h"
18 #include "KVNamedParameter.h"
19 #include "KVUnits.h"
20 
21 using namespace std;
22 
24 
25 
26 
27 //_______________________________________________________________________
28 
29 
33 {
34  // Default constructor
35  init();
36 
37 }
38 
39 
40 
43 
45 {
46  // Destructor
47  delete proj;
48  delete targ;
49  if (kb2) delete kb2;
50 
51  if (sim_evt) delete sim_evt;
52  if (rec_evt) delete rec_evt;
53 
54  ClearHistos();
55  ClearTrees();
56 }
57 
58 
59 
66 
68 {
69  //Initialisation des variables
70  //par defaut
71  //theta 0,180 et phi 0,360
72  //tirage isotropique
73  //noyau de reference pour la direction de diffusion Projectile
74 
75  kIPPVector.SetXYZ(0, 0, 0);
76 
77  kchoix_layer = -1;
78  kXruth_evt = 0;
79  kTreatedNevts = 0;
80 
81  th_min = 1e-3;
82  th_max = 180;
83  phi_min = 0;
84  phi_max = 360;
85 
86  lhisto = 0;
87  ltree = 0;
88 
89  targ = new KVNucleus();
90  proj = new KVNucleus();
91 
92  rec_evt = 0;
93  sim_evt = 0;
94  ktarget = 0;
95  kb2 = 0;
96 
97  ResetBit(kProjIsSet);
98  ResetBit(kTargIsSet);
99  ResetBit(kHasTarget);
100  ResetBit(kIsUpdated);
101  ResetBit(kIsDetectionOn);
102 
103  SetDiffNucleus("PROJ");
104  SetRandomOption("isotropic");
105 
106  SetDetectionOn(kFALSE);
107  ChooseKinSol(1);
108 
109 }
110 
111 
112 
115 
117 {
118  //Set contents/entries to zero for predifined histograms, trees
119  kTreatedNevts = 0;
120  ResetHistos();
121  ResetTrees();
122 
123 }
124 
125 
126 
131 
133 {
134  //Define the entrance channel using KVDBSystem object
135  //Get projectile and target via KVDBSystem::GetKinematics()
136  //Get target material using KVDBSystem::GetTarget()
137  if (sys->GetKinematics()) {
138  SetProjNucleus(sys->GetKinematics()->GetNucleus(1));
139  SetTargNucleus(sys->GetKinematics()->GetNucleus(2));
140  SetTargetMaterial(sys->GetTarget());
141  }
142  else {
143  Error("SetSystem", "KVDBSystem pointer is not valid");
144  return;
145  }
146 }
147 
148 
149 
154 
156 {
157  //Define the entrance channel
158  //zp, ap, ekin, atomic number, mass number and kinetic energy (MeV) of the projectile
159  //zt, at, atomic number, mass number of the target
160  KVNucleus nn(zp, ap, ekin);
161  SetProjNucleus(&nn);
162  nn.SetZAandE(zt, at, 0.0);
163  SetTargNucleus(&nn);
164 
165 }
166 
167 
168 
171 
173 {
174  //Define new target nucleus
175  if (!nuc) return;
176  nuc->Copy(*targ);
177  targ->SetName("TARG");
178  SetBit(kTargIsSet);
179  ResetBit(kIsUpdated);
180 
181 }
182 
183 
184 
189 
191 {
192  //Defini le noyau auquel se réfère la direction de diffusion (theta, phi)
193  //name="PROJ" (default) diffusion du projectile
194  //name="TARG" diffusion de la cible
195 
196  if (name == "TARG") {
197  kDiffNuc = 4;
198  }
199  else {
200  kDiffNuc = 3;
201  if (name != "PROJ") {
202  Warning("SetDiffNucleus", "%s Le nom en argument doit etre PROJ ou TARG, par defaut on choisit le projectile", name.Data());
203  }
204  }
205 }
206 
207 
208 
213 
215 {
216  //Defini le mode de tirage aleatoire pour l'angle polaire
217  //opt="isotropic" ou "" (defaut) ou "random"
218  //voir KVPosition::GetRandomDirection()
219  kRandomOpt = opt;
220 
221 }
222 
223 
224 
227 
229 {
230 
231  //retoune kTRUE si l'option choisi est isotrope
232  return strcmp(kRandomOpt, "random");
233 
234 }
235 
236 
237 
245 
247 {
248  //Dans le cas d'une cinematique inverse (Zproj>Ztarget)
249  //deux solutions cinetiques sont possibles pour un meme angle de
250  //diffusion du projectile
251  //choix=1, on traite seulement la premiere solution ie diffusion a l arriere de la cible
252  //choix=2, on traite seulement la deuxieme solution ie diffusion a l avant de la cible
253  //choix=0, les deux solution sont traitees avec la meme probabilite
254  if (choix >= 0 && choix <= 2)
255  kChoixSol = choix;
256 
257 }
258 
259 
260 
263 
265 {
266  //Define new projectile nucleus
267  if (!nuc) return;
268  nuc->Copy(*proj);
269  proj->SetName("PROJ");
270  SetBit(kProjIsSet);
271  ResetBit(kIsUpdated);
272  proj->SetE0();
273 
274 }
275 
276 
277 
280 
282 {
283  //return the current projectile ("PROJ") or the target ("TARG") nucleus
284  KVString sname(name);
285  if (sname == "PROJ") return proj;
286  else if (sname == "TARG") return targ;
287  else return 0;
288 
289 }
290 
291 
292 
295 
297 {
298  //return the current projectile (ii=1) or the target (ii==2) nucleus
299  if (ii == 1) return proj;
300  else if (ii == 2) return targ;
301  else return 0;
302 
303 }
304 
305 
306 
310 
312 {
313  //Define a new target material where the nuclei will be propagated
314  // if IsRandomized=kTRUE, the interaction point are randomly determined
315  if (!targ) return;
316  if (ktarget) delete ktarget;
317 
318  ktarget = new KVTarget(*targ);
319  ktarget->SetRandomized(IsRandomized);
320  SetBit(kHasTarget);
321  ResetBit(kIsUpdated);
322 }
323 
324 
325 
327 
329 {
330 
331  if (gMultiDetArray) {
334  SetBit(kIsDetectionOn, On);
335  ResetBit(kIsUpdated);
336 
337  }
338  else {
339  if (On)
340  Warning("MakeDetection", "Detection asked but no multiDetArray defined");
341  }
342 
343 }
344 
345 
346 
349 
351 {
352  //return the current target material
353  return ktarget;
354 
355 }
356 
357 
358 
373 
375 {
376  //Protected Method called by ValidateEntrance() method
377  //Genere the KV2Body object which manage the 2 body kinematics
378  //for the elastik scatter
379  //
380  //Store the original momentum of the projectile nuclei
381  //
382  //
383  //Define the KVEvent and KVReconstructedEvent pointer
384  //where are stored the projectile/target nuclei couple after diffusion / detection
385  //StartEvents() methods
386  //
387  //Make a copy of projectile and target nuclei for the KVEvent
388  //
389 
390  Info("GenereKV2Body", "On genere le KV2Body");
391  if (kb2) delete kb2;
392  kb2 = new KV2Body(*proj, *targ);
393  kb2->CalculateKinematics();
394 
395  //GetNucleus("PROJ")->SetE0();
396 
397  //Creer ou clear les deux pointeurs associes aux evts simules et reconstruits
398  StartEvents();
399 
400  GetNucleus("PROJ")->Copy(*sim_evt->AddParticle());
401  GetNucleus("TARG")->Copy(*sim_evt->AddParticle());
402 }
403 
404 
405 
409 
411 {
412  //Define the KVEvent and KVReconstructedEvent pointer
413  //where are stored the projectile/target nuclei couple after diffusion / detection
414  if (!sim_evt) sim_evt = new KVEvent();
415  else sim_evt->Clear();
416 
417  if (!rec_evt) rec_evt = new KVReconstructedEvent();
418  else rec_evt->Clear();
419 
420 }
421 
422 
423 
428 
430 {
431 
432  //Define the nucleus target from the type
433  //of the given layer which belongs to the predefined target
434  //layer_name has to be the chemical symbol of the material
435 
436  if (!IsTargMatSet())
437  return kFALSE;
438 
439  if (GetTarget()->GetLayers()->GetEntries() == 0) return kFALSE;
440  KVMaterial* mat = 0;
441  if (layer_name == "") {
442  kchoix_layer = 1;
443  mat = GetTarget()->GetLayerByIndex(kchoix_layer);
444  }
445  else {
446  if (!(mat = GetTarget()->GetLayer(layer_name))) {
447  printf("le nom du layer %s ne correspond a aucun present dans le cible\n", layer_name.Data());
448  printf("Attention le layer doit etre appele par son symbol (ie Calcium -> Ca)");
449  ktarget->GetLayers()->Print();
450  return kFALSE;
451  }
452  kchoix_layer = GetTarget()->GetLayerIndex(layer_name);
453  }
454 
455  KVNucleus* nuc = new KVNucleus();
456  nuc->SetZandA(TMath::Nint(mat->GetZ()), TMath::Nint(mat->GetMass()));
457 
458  SetTargNucleus(nuc);
459 
460  return kTRUE;
461 
462 }
463 
464 
465 
491 
493 {
494  //Check if there is :
495  // - one define projectile nuclei : SetProjNucleus()
496  // - one define target nuclei : SetTargNucleus()
497  // - one define material target : SetTargetMaterial() [Optionnel]
498  //Si pas de noyau cible défini mais une cible est definie le noyau cible est definie a partir de celle-ci
499  //Si la cible comporte plusieurs layers, le premier est choisi pour definir le noyau cible
500  // see DefineTargetNucleusFromLayer() method
501  //
502  //
503  //Check if the gMultiDetArray is valid, put it in SimMode in order to able the detection and reconstruction
504  //If GetTarget() return kTRUE, put it in the gMultiDetArray
505  //If not, check if there is already one in the gMultiDetArray
506  //If there is one, use it for the following
507  //If there is no target material at all make diffusion without
508  //
509  //Generate the KV2Body object to calculate kinematics of the elastic scatter
510  //
511  //if histograms and trees is defined do nothing for this objects
512  //if not DefineTrees() and DefineHistos() are called.
513  //if you want to regenerate histograms and/or trees
514  //call ClearHistos() and/or ClearTrees() before using ValidateEntrance()
515  //
516  //Return kTRUE if everything is ready
517  //
518 
519  if (!IsProjNucSet()) {
520  Error("ValidateEntrance", "Il n'y a pas de noyau projectile -> use SetProjNucleus()");
521  return kFALSE;
522  }
523 
524  if (!IsTargNucSet()) {
525  Info("ValidateEntrance", "Il n'y a pas de noyau cible");
526  if (!IsTargMatSet()) {
527  Error("ValidateEntrance", "Il n'y a pas de noyau cible -> use SetTargNucleus() ou SetTargetMaterial");
528  return kFALSE;
529  }
530  else {
531  if (DefineTargetNucleusFromLayer())
532  Info("ValidateEntrance", "Definition du noyau cible via DefineTargetNucleusFromLayer()");
533  }
534  }
535 
536  if (IsDetectionOn()) {
537  if (GetTarget()) {
538  if (gMultiDetArray->GetTarget() && gMultiDetArray->GetTarget() == ktarget) {
539 
540  }
541  else {
542  //Fait un clone
543  gMultiDetArray->SetTarget(GetTarget());
544  delete ktarget;
545  ktarget = gMultiDetArray->GetTarget();
546  }
547  }
548  else if (gMultiDetArray->GetTarget()) {
549  ktarget = gMultiDetArray->GetTarget();
550  }
551  else {
552  Warning("ValidateEntrance", "Pas de calcul de perte dans la cible ... ");
553  }
554  //DefineAngularRange(gMultiDetArray);
555  }
556  else {
557  Info("ValidateEntrance", "The elastic scatter events will not be detected/filtered");
558  }
559 
560  GenereKV2Body();
561  Info("ValidateEntrance", "Fin de GenereKV2Body");
562  //Define histograms/trees only at the first call
563  //of validate entrance
564 
565  if (!ltree) DefineTrees();
566 
567  ClearHistos();
568  DefineHistos();
569 
570  SetBit(kIsUpdated);
571 
572  return kTRUE;
573 
574 }
575 
576 
577 
578 
582 
584 {
585  //process ntimes elastic scatter
586  //if reset=kTRUE, reset histograms, trees and counter before
587  Info("Process", "Je Suis dans Process ... Youpee");
588 
589  if (!IsUpdated())
590  if (!ValidateEntrance()) return;
591 
592  if (reset) Reset();
593  Int_t nn = 0;
594  while (nn < ntimes) {
595  MakeDiffusion();
596  nn += 1;
597  if ((nn % 1000) == 0) {
598  Info("Process", "%d/%d diffusions treated", nn, ntimes);
599  }
600  }
601  Info("Process", "End of process : %d diffusions performed", kTreatedNevts);
602 }
603 
604 
605 
617 
619 {
620  //
621  //Propagation dans la cible du projectile jusqu'au point d interaction
622  // PropagateInTargetLayer();
623  //Tirage aleatoire d'un couple theta phi pour la direction de diffusion du projectile
624  //Determination de la cinematique de la voie de sortie
625  // SetAnglesForDiffusion(the,phi);
626  //Filtre si un multidetecteur est defini
627  // Filter
628  //Traitement de l evt (remplissage d'arbre ou d'histo
629  // TreateEvent();
630 
631  KVNucleus* knuc = 0;
632  sim_evt->GetParameters()->Clear();
633  while ((knuc = sim_evt->GetNextParticle())) {
634  knuc->GetParameters()->Clear();
635  knuc->RemoveAllGroups();
636  }
637 
638  //-------------------------
639  if (IsTargMatSet()) {
640  NewInteractionPointInTargetLayer();
641  PropagateInTargetLayer();
642  }
643  //-------------------------
644 
645 
646  Double_t tmin = GetTheta("min");
647  if (tmin >= kb2->GetMaxAngleLab(kDiffNuc)) {
648  GetNucleus("PROJ")->SetMomentum(*GetNucleus("PROJ")->GetPInitial());
649  return;
650  }
651  Double_t tmax = GetTheta("max");
652 
653  if (tmax <= kb2->GetMinAngleLab(kDiffNuc)) {
654  GetNucleus("PROJ")->SetMomentum(*GetNucleus("PROJ")->GetPInitial());
655  return;
656  }
657 
658  if (tmin < kb2->GetMinAngleLab(kDiffNuc)) tmin = kb2->GetMinAngleLab(kDiffNuc);
659  if (tmax > kb2->GetMaxAngleLab(kDiffNuc)) tmax = kb2->GetMaxAngleLab(kDiffNuc);
660 
661  Double_t pmin = GetPhi("min");
662  Double_t pmax = GetPhi("max");
663 
664  kposalea.SetPolarMinMax(tmin, tmax);
665  kposalea.SetAzimuthalMinMax(pmin, pmax);
666 
667  Double_t th_deg, ph_deg;
668  kposalea.GetRandomAngles(th_deg, ph_deg, kRandomOpt);
669 
670  SetAnglesForDiffusion(th_deg, ph_deg);
671 
672  ((TH2F*)lhisto->FindObject("phi_theta"))->Fill(th_deg, ph_deg);
673  ((TH1F*)lhisto->FindObject("costheta"))->Fill(TMath::Cos(TMath::DegToRad()*th_deg));
674 
675 
676  if (IsDetectionOn()) {
677  Filter();
678  }
679  else {
680  if (IsTargMatSet())
681  SortieDeCible();
682  }
683 
684  TreateEvent();
685 
686  kTreatedNevts += 1;
687 }
688 
689 
690 
695 
697 {
698  //Choose a new interaction point in the current target layer
699  //This position can be read via the GetInteractionPointInTargetLayer()
700  //method
701  if (kchoix_layer != -1) {
702  TVector3 dir = GetNucleus("PROJ")->GetMomentum();
703  ktarget->SetInteractionLayer(kchoix_layer, dir);
704  //kIPPVector = ktarget->GetInteractionPoint();
705  }
706  kIPPVector = ktarget->GetInteractionPoint(GetNucleus("PROJ"));
707  ((TH1F*)lhisto->FindObject("target_layer_depth"))->Fill(kIPPVector.Z());
708 }
709 
710 
711 
720 
722 {
723  //Apply Energy loss calculation to the entering projectile
724  //along its path in the target layer to the interation point
725  //
726  //if a gMultiDetArray is defined the outgoing (after diffusion) pathes are not treated here but
727  //in the Filter() method
728  //
729  //if not is treated in the SortieDeCible method
730 
731  Double_t eLostInTarget = GetNucleus("PROJ")->GetKE();
732  ktarget->SetIncoming(kTRUE);
733  ktarget->DetectParticle(GetNucleus("PROJ"), 0);
734  eLostInTarget -= GetNucleus("PROJ")->GetKE();
735 
736  //Energie perdue jusqu'au point d interaction
737  sim_evt->GetParticleWithName("PROJ")->GetParameters()->SetValue("TARGET In", eLostInTarget);
738  //On modifie l'energie du projectile dans KV2Body
739  //pour prendre en compte l energie deposee dans la cible
740  //avant de faire le calcul de la cinematique
741  if (GetNucleus("PROJ")->GetKE() == 0) {
742  GetNucleus("PROJ")->Print();
743  printf("%lf / %lf\n", eLostInTarget, proj->GetKE());
744  }
745  kb2->GetNucleus(1)->SetKE(GetNucleus("PROJ")->GetKE());
746  kb2->CalculateKinematics();
747 
748  ktarget->SetIncoming(kFALSE);
749 
750 }
751 
752 
753 
756 
758 {
759  //return the last interaction point in the target
760  return kIPPVector;
761 
762 }
763 
764 
765 
766 
771 
773 {
774  //Apply Energy loss calculation in target material for the outgoing projectile
775  //and target
776  //
777 
778  ktarget->SetOutgoing(kTRUE);
779  KVNucleus* knuc = 0;
780  while ((knuc = sim_evt->GetNextParticle())) {
781  knuc->SetE0();
782  Double_t eLostInTarget = knuc->GetKE();
783  ktarget->DetectParticle(knuc, 0);
784  eLostInTarget -= knuc->GetKE();
785  knuc->GetParameters()->SetValue("TARGET Out", eLostInTarget);
786  knuc->SetMomentum(*knuc->GetPInitial());
787  }
788 
789  ktarget->SetOutgoing(kFALSE);
790 
791 }
792 
793 
794 
795 
804 
806 {
807  //Determination a partir du theta choisi
808  //de l'energie cinetique du projectile diffuse
809  //All kinematics properties calculated in the KV2Body class
810  //are accessible via the KV2Body& GetKinematics() method
811  //
812  // WARNING: in inverse kinematics, there are two projectile energies
813  // for each lab angle.
814  Double_t ediff;
815  Double_t ediff1, ediff2;
816  Int_t nsol_kin;
817 
818  nsol_kin = kb2->GetELab(kDiffNuc, theta, kDiffNuc, ediff1, ediff2);
819 
820  if (nsol_kin == 1) {
821  ediff = ediff1;
822  }
823  else {
824  if (kChoixSol == 1) ediff = ediff1;
825  else if (kChoixSol == 2) {
826  ediff = ediff2;
827  }
828  else {
829  Int_t choix = TMath::Nint(gRandom->Uniform(0.5, 2.5));
830  if (choix == 2) ediff = ediff2;
831  else ediff = ediff1;
832  }
833  }
834  Bool_t sol2 = (ediff == ediff2);
835 
836  kXruth_evt = kb2->GetXSecRuthLab(theta, kDiffNuc);
837 
838  //On modifie l energie et les angles du projectile ou cible diffusé(e)
839  //puis par conservation, on deduit ceux du noyau complementaire
840  KVNucleus* knuc = 0;
841 
842  if (kDiffNuc == 3)
843  knuc = sim_evt->GetParticleWithName("PROJ");
844  else
845  knuc = sim_evt->GetParticleWithName("TARG");
846 
847  knuc->SetKE(ediff);
848  knuc->SetTheta(theta);
849  knuc->SetPhi(phi);
850  ((TH2F*)lhisto->FindObject("ek_theta"))->Fill(knuc->GetTheta(), knuc->GetKE());
851 
852  //Conservation de l impulsion
853  //Conservation de l energie tot
854  TVector3 ptot = proj->Vect() + targ->Vect();
855  Double_t etot = proj->E() + targ->E();
856  //on retire la contribution du noyau diffusé
857  ptot -= knuc->Vect();
858  etot -= knuc->E();
859  //on met a jour les pptés la cible ou projectile diffusé(e)
860 
861  if (kDiffNuc == 3)
862  knuc = sim_evt->GetParticleWithName("TARG");
863  else
864  knuc = sim_evt->GetParticleWithName("PROJ");
865 
866  TLorentzVector q(ptot, etot);
867  knuc->Set4Mom(q);
868 
869  ((TH2F*)lhisto->FindObject("ek_theta"))->Fill(knuc->GetTheta(), knuc->GetKE());
870 
871  sim_evt->SetNumber(kTreatedNevts);
872 
873  sim_evt->GetParameters()->SetValue("XRuth", kXruth_evt);
874  sim_evt->GetParameters()->SetValue("ThDiff", theta);
875  sim_evt->GetParameters()->SetValue("EkDiff", ediff1);
876  sim_evt->GetParameters()->SetValue("IPz", kIPPVector.Z());
877 
878  if (sol2)
879  sim_evt->GetParameters()->SetValue("Sol2", 1);
880 
881  //L' energie cinetique du projectile est reinitialisee
882  //pour la prochaine diffusion
883  GetNucleus("PROJ")->SetMomentum(*GetNucleus("PROJ")->GetPInitial());
884 
885 }
886 
887 
888 
893 
895 {
896  //Simulate passage of the projectile/target couple
897  //through the multidetector refered by the gMultiDetArray pointer
898  //if it is not valid do nothing
899 
900  gMultiDetArray->DetectEvent(sim_evt, rec_evt);
901 
902 
903 }
904 
905 
906 
926 
928 {
929  //Rempli l'arbre ElasticScatter
930  //Boucle sur tous les parametres associés a l evt (KVEvent::GetParameters() et au projectiles et cible
931  //qui le constituent GetParticle(1)->GetParameters()
932  // Chaque parametre devient un alias de l'arbre ElasticScatter
933  // pour une utilisation a posteriori plus facile.
934  // - pour les parametres de l'evt, on donne directement le nom du parametre
935  // - pour les particule :
936  // N1_[nom_du parametre] pour les projectiles et N2_[nom_du parametre] pour les cibles diffusés
937  //
938  // Exemple avec l'utilisation de TTree::Draw
939  // Si on veut voir le spectre en energie laissé par les projectiles diffuses dans la "CI_0601"
940  // au lieu de faire
941  // GetTree("Simulated_evts")->Draw("Simulated_evts->GetParticle(1)->GetParameters()->GetValue(\"CI_0601\")")
942  // on fera
943  // GetTree("Simulated_evts")->Draw("N1_CI_0601")
944  //
945  // Generation des correlation Energie Cinetique (Ek) vs Angle de diffusion (theta)
946  // pour tous les cas de détection
947 
948  if (rec_evt->GetMult() > 0) {
949  TTree* tt = (TTree*)ltree->FindObject("ElasticScatter");
950  tt->Fill();
951  }
952  /*
953  TString snom,stit;
954  KVNamedParameter* nm = 0;
955 
956  TIter it(sim_evt->GetParameters()->GetList());
957  while ( (nm = (KVNamedParameter* )it.Next()) ){
958  snom.Form("%s",nm->GetName());
959  if (snom.Contains(" ")) snom.ReplaceAll(" ","_");
960  stit.Form("Simulated_evts->GetParameters()->GetDoubleValue(\"%s\")",nm->GetName());
961  tt->SetAlias(snom.Data(),stit.Data());
962  }
963  TIter it1(sim_evt->GetParticle(1)->GetParameters()->GetList());
964  while ( (nm = (KVNamedParameter* )it1.Next()) ){
965  snom.Form("N1_%s",nm->GetName());
966  if (snom.Contains(" ")) snom.ReplaceAll(" ","_");
967  stit.Form("Simulated_evts->GetParticle(1)->GetParameters()->GetDoubleValue(\"%s\")",nm->GetName());
968  tt->SetAlias(snom.Data(),stit.Data());
969  }
970  TIter it2(sim_evt->GetParticle(2)->GetParameters()->GetList());
971  while ( (nm = (KVNamedParameter* )it2.Next()) ){
972  snom.Form("N2_%s",nm->GetName());
973  if (snom.Contains(" ")) snom.ReplaceAll(" ","_");
974  stit.Form("Simulated_evts->GetParticle(2)->GetParameters()->GetDoubleValue(\"%s\")",nm->GetName());
975  tt->SetAlias(snom.Data(),stit.Data());
976  }
977 
978  if (IsDetectionOn()){
979  TH2F* hh = 0;
980  KVNucleus* knuc = 0;
981  while ( (knuc = sim_evt->GetNextParticle()) ){
982 
983  if ( knuc->GetParameters()->HasParameter("DETECTED") ){
984  if ( !strcmp(knuc->GetParameters()->GetStringValue("DETECTED"),"") ){
985  snom.Form("ek_theta_DETECTED");
986  if ( !(hh = (TH2F* )lhisto->FindObject(snom.Data())) ){
987  Double_t totalE = GetNucleus(1)->GetKE();
988  lhisto->Add(new TH2F(snom.Data(),"DETECTED",180,0,180,TMath::Nint(totalE*11),0,totalE*1.1));
989  hh = (TH2F* )lhisto->Last(); hh->SetXTitle("#theta_{ lab}"); hh->SetYTitle("Ek_{ lab}");
990  }
991  }
992  else {
993  snom.Form("ek_theta_%s",knuc->GetParameters()->GetStringValue("DETECTED"));
994  if ( !(hh = (TH2F* )lhisto->FindObject(snom.Data())) ){
995  Double_t totalE = GetNucleus(1)->GetKE();
996  lhisto->Add(new TH2F(snom.Data(),knuc->GetParameters()->GetStringValue("DETECTED"),180,0,180,TMath::Nint(totalE*11),0,totalE*1.1));
997  hh = (TH2F* )lhisto->Last(); hh->SetXTitle("#theta_{ lab}"); hh->SetYTitle("Ek_{ lab}");
998  }
999  }
1000  }
1001  else if ( knuc->GetParameters()->HasParameter("UNDETECTED") ){
1002  if ( !strcmp(knuc->GetParameters()->GetStringValue("UNDETECTED"),"") ){
1003  snom.Form("ek_theta_UNDETECTED");
1004  if ( !(hh = (TH2F* )lhisto->FindObject(snom.Data())) ){
1005  Double_t totalE = GetNucleus(1)->GetKE();
1006  lhisto->Add(new TH2F(snom.Data(),"UNDETECTED",180,0,180,TMath::Nint(totalE*11),0,totalE*1.1));
1007  hh = (TH2F* )lhisto->Last(); hh->SetXTitle("#theta_{ lab}"); hh->SetYTitle("Ek_{ lab}");
1008  }
1009  }
1010  else {
1011  snom.Form("ek_theta_%s",knuc->GetParameters()->GetStringValue("UNDETECTED"));
1012  if ( !(hh = (TH2F* )lhisto->FindObject(snom.Data())) ){
1013  Double_t totalE = GetNucleus(1)->GetKE();
1014  lhisto->Add(new TH2F(snom.Data(),knuc->GetParameters()->GetStringValue("UNDETECTED"),180,0,180,TMath::Nint(totalE*11),0,totalE*1.1));
1015  hh = (TH2F* )lhisto->Last(); hh->SetXTitle("#theta_{ lab}"); hh->SetYTitle("Ek_{ lab}");
1016  }
1017  }
1018  }
1019  else {
1020  knuc->Print();
1021  }
1022  if (hh)
1023  hh->Fill(knuc->GetTheta(),knuc->GetKE());
1024  }
1025  }
1026  */
1027 
1028 }
1029 
1030 
1031 
1033 
1035 {
1036 
1037  kb2->Print();
1038 
1039  printf("#####################\n");
1040  printf("## KVElasticScatterEvent::Print() ##\n");
1041  printf("# Diffusion elastique traitee :\n");
1042  printf("# %s+%s@%1.2lf MeV/A\n",
1043  GetNucleus(1)->GetSymbol(),
1044  GetNucleus(2)->GetSymbol(),
1045  GetNucleus(1)->GetKE() / GetNucleus(1)->GetA()
1046  );
1047  if (IsTargMatSet()) {
1048  printf("-------------------------\n");
1049  printf("# Propagation dans une cible de:\n");
1050  for (Int_t nn = 0; nn < GetTarget()->GetLayers()->GetEntries(); nn += 1) {
1051  Double_t epaiss = GetTarget()->GetLayerByIndex(nn + 1)->GetAreaDensity() / (KVUnits::mg / pow(KVUnits::cm, 2));
1052  printf("#\ttype:%s epaisseur:%lf (mg/cm**2) / %lf\n",
1053  GetTarget()->GetLayerByIndex(nn + 1)->GetType(),
1054  epaiss,
1055  GetTarget()->GetLayerByIndex(nn + 1)->GetThickness()
1056  );
1057  }
1058  }
1059  printf("-------------------------\n");
1060  if (IsDetectionOn()) {
1061  printf("# Detection par %s\n", gMultiDetArray->GetName());
1062  }
1063  printf("#####################\n");
1064 
1065 }
1066 
1067 
1074 
1076 {
1077  //Definition of control histograms
1078  //- phi_theta : filled with angles choosen to determine the direction of the diffused projectile
1079  //- target_layer_depth : interaction point position in the target
1080  //- ek_theta : filled with energies and polar angles of projectile and target nuclei after diffusion
1081  //they are detected by the multidetarray
1082  Info("DefineHistos", "DefineHistos");
1083  lhisto = new KVHashList();
1084  lhisto->SetOwner(kTRUE);
1085  TH2F* h2 = 0;
1086  TH1F* h1 = 0;
1087 
1088  lhisto->Add(new TH2F("phi_theta", "phi_theta", 180, 0, 180, 360, 0, 360));
1089  h2 = (TH2F*)lhisto->Last();
1090  h2->SetXTitle("#theta_{ lab}");
1091  h2->SetYTitle("#phi_{ lab}");
1092  lhisto->Add(new TH1F("costheta", "costheta", 200, -1, 1));
1093  h1 = (TH1F*)lhisto->Last();
1094  h1->SetXTitle("cos #theta_{ lab}");
1095  if (IsTargMatSet()) {
1096  Info("DefineHistos", "Creation de l histo interaction dans la cible");
1097  Float_t thickness = GetTarget()->GetThickness();
1098  lhisto->Add(new TH1F("target_layer_depth", "target_layer_depth", TMath::Nint(thickness * 110), 0, thickness * 1.1));
1099  h1 = (TH1F*)lhisto->Last();
1100  h1->SetXTitle("IP position (mg / cm**2)");
1101  }
1102  Float_t totalE = GetNucleus(1)->GetKE();
1103  lhisto->Add(new TH2F("ek_theta", "ek_theta", 180, 0, 180, TMath::Nint(totalE * 11), 0, totalE * 1.1));
1104  h2 = (TH2F*)lhisto->Last();
1105  h2->SetXTitle("#theta_{ lab}");
1106  h2->SetYTitle("Ek_{ lab}");
1107 
1108 }
1109 
1110 
1111 
1114 
1116 {
1117  //return the list where histo are stored
1118  return lhisto;
1119 
1120 }
1121 
1122 
1125 
1127 {
1128 
1129  //Reset histo in the list
1130  lhisto->Execute("Reset", "");
1131 
1132 }
1133 
1134 
1137 
1139 {
1140  //Efface la liste des histo et leur contenu et met le pointeur a zero
1141  if (lhisto) delete lhisto;
1142  lhisto = 0;
1143 
1144 }
1145 
1146 
1147 
1158 
1160 {
1161  //Definition of trees
1162  //
1163  //Par defaut un seul arbre : ElasticScatter where simulated events are stored
1164  //sous forme de KVEvent
1165  //Lors du remplissage de l arbre ( methode TreateEvent)
1166  //les parametres associes au KVEvent::GetParameters et au KVNucleus::GetParameters (projectile et cible)
1167  //sont scannés et leur nom et le moyen d'y accéder ajouté aux alias de l arbre pour une utilisation
1168  //plus aisé de celui_ci
1169  //
1170 
1171  Info("DefineTrees", "DefineTrees");
1172  ltree = new KVHashList();
1173  ltree->SetOwner(kTRUE);
1174  TTree* tt = 0;
1175 
1176  tt = new TTree("ElasticScatter", IsA()->GetName());
1177  KVEvent::MakeEventBranch(tt, "Simulated_evts", "KVEvent", sim_evt);
1178  ltree->Add(tt);
1179 }
1180 
1181 
1182 
1185 
1187 {
1188  //return the list where histo are stored
1189  return ltree;
1190 
1191 }
1192 
1193 
1194 
1197 
1199 {
1200  //Efface la liste des arbres et leur contenu et met le pointeur a zero
1201  if (ltree) delete ltree;
1202  ltree = 0;
1203 
1204 }
1205 
1206 
1207 
1211 
1213 {
1214  //Reset the tree contents
1215  //and aliases for the "ElasticScatter" tree
1216  ltree->Execute("Reset", "");
1217  if (((TTree*)ltree->FindObject("ElasticScatter"))->GetListOfAliases())
1218  ((TTree*)ltree->FindObject("ElasticScatter"))->GetListOfAliases()->Clear();
1219 
1220 }
1221 
1222 
1223 
1229 
1231 {
1232  //Define in which angular (polar and azimuthal) range
1233  //The projectile diffusion direction will be randomized
1234  //If this method is not used
1235  //Default range is \theta [0,180] and \phi [0,360]
1236  if (tmin != -1) th_min = tmin;
1237  if (tmax != -1) th_max = tmax;
1238  if (pmin != -1) phi_min = pmin;
1239  if (pmax != -1) phi_max = pmax;
1240 
1241 }
1242 
1243 
1244 
1253 
1255 {
1256  //Define in which angular (polar and azimuthal) range
1257  //The projectile diffusion direction will be randomized
1258  //From the geometry of the given object obj
1259  //This method is defined for object derived from
1260  // - KVPosition (ie KVTelescope KVRing KVGroup etc ...)
1261  // - KVDetector (in this case, the KVTelescope it belongs to is used)
1262  // - KVMultiDetArray
1263 
1264  Double_t tmin = -1, tmax = -1, pmin = -1, pmax = -1;
1265  if (obj->InheritsFrom("KVPosition")) {
1266  KVPosition* pos_obj = (KVPosition*)obj;
1267  tmin = pos_obj->GetThetaMin();
1268  tmax = pos_obj->GetThetaMax();
1269  pmin = pos_obj->GetPhiMin();
1270  pmax = pos_obj->GetPhiMax();
1271  }
1272  else if (obj->InheritsFrom("KVDetector")) {
1273  KVTelescope* pos_obj = (KVTelescope*)((KVDetector*)obj)->GetParentStructure("TELESCOPE");
1274  tmin = pos_obj->GetThetaMin();
1275  tmax = pos_obj->GetThetaMax();
1276  pmin = pos_obj->GetPhiMin();
1277  pmax = pos_obj->GetPhiMax();
1278  }
1279  else if (obj->InheritsFrom("KVASMultiDetArray")) {
1280  Warning("DefineAngularRange(KVASMultiDetArray*)", "Needs reimplementing");
1281  /*KVASMultiDetArray* pos_obj=(KVASMultiDetArray* )obj;
1282  KVSeqCollection* ll = pos_obj->GetGroups();
1283  KVASGroup* gr=0;
1284  Double_t tmin2=180,tmax2=0;
1285  Double_t pmin2=360,pmax2=0;
1286  for (Int_t nl=0;nl<ll->GetEntries();nl+=1){
1287  gr = (KVASGroup* )ll->At(nl);
1288  if (gr->GetThetaMin()<tmin2) tmin2 = gr->GetThetaMin();
1289  if (gr->GetThetaMax()>tmax2) tmax2 = gr->GetThetaMax();
1290  if (gr->GetPhiMin()<pmin2) pmin2 = gr->GetPhiMin();
1291  if (gr->GetPhiMax()>pmax2) pmax2 = gr->GetPhiMax();
1292  }
1293  tmin = tmin2;
1294  tmax = tmax2;
1295  pmin = pmin2;
1296  pmax = pmax2;*/
1297  }
1298  else {
1299  printf("les objects de type %s ne sont pas implemente dans KVElasticScatterEvent::DefineAngularRange\n", obj->IsA()->GetName());
1300  return;
1301  }
1302 
1303  DefineAngularRange(tmin, tmax, pmin, pmax);
1304 
1305 }
1306 
1307 
1308 
1312 
1314 {
1315  //Return the limite in theta range (polar angle)
1316  //opt has to be "min" or "max"
1317  if (opt == "min") return th_min;
1318  else if (opt == "max") return th_max;
1319  else return -1;
1320 
1321 }
1322 
1323 
1324 
1328 
1330 {
1331 
1332  //Return the limite in phi range (azimuthal angle)
1333  //opt has to be "min" or "max"
1334  if (opt == "min") return phi_min;
1335  else if (opt == "max") return phi_max;
1336  else return -1;
1337 
1338 }
1339 
1340 
int Int_t
KVMultiDetArray * gMultiDetArray
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
#define e(i)
char Char_t
const Bool_t kFALSE
bool Bool_t
double Double_t
float Float_t
const Bool_t kTRUE
const char Option_t
float * q
double pow(double, double)
R__EXTERN TRandom * gRandom
Relativistic binary kinematics calculator.
Definition: KV2Body.h:165
KVNucleus * GetNucleus(Int_t i) const
Definition: KV2Body.cpp:457
Database class used to store information on different colliding systems studied during an experiment.
Definition: KVDBSystem.h:51
KVTarget * GetTarget() const
Definition: KVDBSystem.h:78
KV2Body * GetKinematics()
Definition: KVDBSystem.cpp:79
Base class for detector geometry description, interface to energy-loss calculations.
Definition: KVDetector.h:121
simulate ElasticScatterEvent and answer of a given (multi-)detector : A + B -> A + B
virtual void SetProjNucleus(KVNucleus *nuc)
Define new projectile nucleus.
virtual void SetSystem(KVDBSystem *sys)
virtual TVector3 & GetInteractionPointInTargetLayer()
return the last interaction point in the target
void Print(Option_t *="") const
virtual void ClearTrees()
Efface la liste des arbres et leur contenu et met le pointeur a zero.
virtual void DefineAngularRange(TObject *)
void SetDiffNucleus(KVString name="PROJ")
virtual void SetTargNucleus(KVNucleus *nuc)
Define new target nucleus.
Bool_t IsIsotropic() const
retoune kTRUE si l'option choisi est isotrope
virtual void SetAnglesForDiffusion(Double_t theta, Double_t phi)
void SetRandomOption(Option_t *opt="isotropic")
KVHashList * GetTrees() const
return the list where histo are stored
KVTarget * GetTarget() const
return the current target material
virtual void Process(Int_t ntimes=1, Bool_t reset=kTRUE)
virtual void ClearHistos()
Efface la liste des histo et leur contenu et met le pointeur a zero.
virtual void ResetHistos()
Reset histo in the list.
virtual void SetTargetMaterial(KVTarget *targ, Bool_t IsRandomized=kTRUE)
void SetDetectionOn(Bool_t On=kTRUE)
KVHashList * GetHistos() const
return the list where histo are stored
virtual ~KVElasticScatterEvent()
Destructor.
virtual void Reset()
Set contents/entries to zero for predifined histograms, trees.
KVNucleus * GetNucleus(const Char_t *name) const
return the current projectile ("PROJ") or the target ("TARG") nucleus
void ChooseKinSol(Int_t choix=1)
Bool_t DefineTargetNucleusFromLayer(KVString layer_name="")
Double_t GetPhi(KVString opt) const
Double_t GetTheta(KVString opt) const
Base class container for multi-particle events.
Definition: KVEvent.h:176
static void MakeEventBranch(TTree *tree, const TString &branchname, const TString &classname, T &event, Int_t bufsize=10000000)
Definition: KVEvent.h:468
Extended version of ROOT THashList.
Definition: KVHashList.h:28
Description of physical materials used to construct detectors; interface to range tables.
Definition: KVMaterial.h:41
Double_t GetZ() const
Returns atomic number of material.
Definition: KVMaterial.cpp:320
Double_t GetMass() const
Returns atomic mass of material. Will be isotopic mass if set.
Definition: KVMaterial.cpp:252
virtual void DetectEvent(KVEvent *event, KVReconstructedEvent *rec_event, const Char_t *detection_frame="")
KVTarget * GetTarget()
void SetFilterType(Int_t t)
void SetTarget(const Char_t *material, const Float_t thickness)
virtual void SetSimMode(Bool_t on=kTRUE)
void SetValue(const Char_t *name, value_type value)
virtual void Clear(Option_t *opt="")
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:125
virtual void Copy(TObject &) const
Copy this KVNucleus into the KVNucleus object referenced by "obj".
Definition: KVNucleus.cpp:834
void SetZandA(Int_t z, Int_t a)
Set atomic number and mass number.
Definition: KVNucleus.cpp:721
void SetZAandE(Int_t z, Int_t a, Double_t ekin)
Set atomic number, mass number, and kinetic energy in MeV.
Definition: KVNucleus.cpp:733
void SetTheta(Double_t theta)
Definition: KVParticle.h:654
TVector3 * GetPInitial() const
Definition: KVParticle.h:687
void SetPhi(Double_t phi)
Definition: KVParticle.h:658
void RemoveAllGroups()
Definition: KVParticle.cpp:646
KVNameValueList * GetParameters() const
Definition: KVParticle.h:735
Double_t GetTheta() const
Definition: KVParticle.h:641
void SetMomentum(const TVector3 &v)
Definition: KVParticle.h:535
void SetKE(Double_t ecin)
Definition: KVParticle.cpp:230
void Set4Mom(const TLorentzVector &p)
Definition: KVParticle.h:550
void SetE0(TVector3 *e=0)
Definition: KVParticle.h:676
Double_t GetKE() const
Definition: KVParticle.h:575
Base class used for handling geometry in a multidetector array.
Definition: KVPosition.h:90
Double_t GetPhiMax() const
Definition: KVPosition.h:145
Double_t GetPhiMin() const
Definition: KVPosition.h:141
Double_t GetThetaMin() const
Definition: KVPosition.h:149
Double_t GetThetaMax() const
Definition: KVPosition.h:153
Physical event reconstructed from data measured with a detector array using implemented identificatio...
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:72
Calculation/correction of energy losses of particles through an experimental target.
Definition: KVTarget.h:126
Associates two detectors placed one behind the other.
Definition: KVTelescope.h:35
virtual Bool_t Add(const TH1 *h, const TH1 *h2, Double_t c1=1, Double_t c2=1)
virtual void SetXTitle(const char *title)
virtual void SetYTitle(const char *title)
TVector3 Vect() const
Double_t E() const
virtual const char * GetName() const
virtual const char * GetName() const
virtual Bool_t InheritsFrom(const char *classname) const
virtual Double_t Uniform(Double_t x1, Double_t x2)
const char * Data() const
RVec< T > Filter(const RVec< T > &v, F &&f)
TH1F * h1
const long double mg
Definition: KVUnits.h:74
const long double cm
Definition: KVUnits.h:66
void Info(const char *location, const char *va_(fmt),...)
void Error(const char *location, const char *va_(fmt),...)
void Warning(const char *location, const char *va_(fmt),...)
Type GetType(const std::string &Name)
Int_t Nint(T x)
constexpr Double_t DegToRad()
Double_t Cos(Double_t)
auto * tt