]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Added charge correlation between cathods.
authoregangler <egangler@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 20 Dec 2000 13:00:22 +0000 (13:00 +0000)
committeregangler <egangler@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 20 Dec 2000 13:00:22 +0000 (13:00 +0000)
In Config_slat.C, use
 MUON->Chamber(chamber-1).SetChargeCorrel(0.11); to set the RMS of
 q1/q2 to 11 % (number from Alberto)
 This is stored in AliMUONChamber fChargeCorrel member.
 At generation time, when a tracks enters the volume,
 AliMUONv1::StepManager calls
 AliMUONChamber::ChargeCorrelationInit() to set the current value of
 fCurrentCorrel which is then used at Disintegration level to scale
 appropriately the PadHit charges.

MUON/AliMUONChamber.cxx
MUON/AliMUONChamber.h
MUON/AliMUONv1.cxx
MUON/Config_slat.C

index 91eec1f45e8c8620d13019f7a67c15d3d4ad7b3b..f5b79d031dc55415b8ce73d7bc34f3a7b9d324c5 100644 (file)
@@ -14,6 +14,9 @@
  **************************************************************************/
 /*
 $Log$
+Revision 1.6  2000/10/09 14:01:12  morsch
+Double inclusion of AliResponse removed.
+
 Revision 1.5  2000/07/03 11:54:57  morsch
 AliMUONSegmentation and AliMUONHitMap have been replaced by AliSegmentation and AliHitMap in STEER
 The methods GetPadIxy and GetPadXxy of AliMUONSegmentation have changed name to GetPadI and GetPadC.
@@ -44,9 +47,14 @@ Revision 1.1.2.3  2000/05/05 10:09:52  morsch
 Log messages included
 */
 
+// --- MUON includes ---
 #include "AliMUONChamber.h"
 
+// --- ROOT includes ---
+
+#include "TRandom.h"
 #include "TMath.h"
+
 ClassImp(AliMUONChamber)       
 
     AliMUONChamber::AliMUONChamber()
@@ -59,6 +67,9 @@ ClassImp(AliMUONChamber)
     fnsec=1;
     fReconstruction=0;
     fId=0;
+    // to avoid mistakes if ChargeCorrelInit is not called
+    fChargeCorrel = 0;
+    fCurrentCorrel =1;
 }
 
     AliMUONChamber::AliMUONChamber(Int_t id) 
@@ -71,6 +82,9 @@ ClassImp(AliMUONChamber)
     fnsec=1;
     fReconstruction=0;
     fId=id;
+    // to avoid mistakes if ChargeCorrelInit is not called
+    fChargeCorrel = 0;
+    fCurrentCorrel =1;
 }
 
 AliMUONChamber::~AliMUONChamber() 
@@ -128,6 +142,16 @@ void    AliMUONChamber::SigGenInit(Float_t x, Float_t y, Float_t z)
     }
 }
 
+void AliMUONChamber::ChargeCorrelationInit() {
+// Initialisation of charge correlation for current hit
+// the value is stored, and then used by Disintegration
+if (fnsec==1) 
+    fCurrentCorrel =1;
+else 
+    // exponential is here to avoid eventual problems in 0 
+    fCurrentCorrel = TMath::Exp(gRandom->Gaus(0,fChargeCorrel/2));
+}
+
 void AliMUONChamber::DisIntegration(Float_t eloss, Float_t tof, 
                                    Float_t xhit, Float_t yhit, Float_t zhit,
                                    Int_t& nnew,Float_t newclust[6][500]) 
@@ -150,8 +174,10 @@ void AliMUONChamber::DisIntegration(Float_t eloss, Float_t tof,
     Float_t qcheck=0, qp;
     nnew=0;
     
+    // Cathode plane loop
     for (Int_t i=1; i<=fnsec; i++) {
        qcheck=0;
+       Float_t qcath = qtot * (i==1? fCurrentCorrel : 1/fCurrentCorrel);
        AliSegmentation * segmentation=
            (AliSegmentation *) (*fSegmentation)[i-1];
        for (segmentation->FirstPad(xhit, yhit, zhit, dx, dy); 
@@ -163,17 +189,17 @@ void AliMUONChamber::DisIntegration(Float_t eloss, Float_t tof,
 //
 //
            if (qp > 1.e-4) {
-               qcheck+=qp*qtot;
+               qcheck+=qp*qcath;
            //
            // --- store signal information
-               newclust[0][nnew]=qtot;                     // total charge
+               newclust[0][nnew]=qcath;                     // total charge
                newclust[1][nnew]=segmentation->Ix();       // ix-position of pad
                newclust[2][nnew]=segmentation->Iy();       // iy-position of pad
-               newclust[3][nnew]=qp * qtot;                // charge on pad
+               newclust[3][nnew]=qp * qcath;                // charge on pad
                newclust[4][nnew]=segmentation->ISector();  // sector id
                newclust[5][nnew]=(Float_t) i;              // counter
                nnew++;
-//             if (i==2) printf("\n i, nnew, q %d %d %f", i, nnew, qp*qtot);
+//             if (i==2) printf("\n i, nnew, q %d %d %f", i, nnew, qp*qcath);
                
            }
        } // Pad loop
index 61639d1bec2f5bb37fa99b7b5c28a813104306df..835b7d7078d214a35f13add99550a07446b81a33 100644 (file)
@@ -81,6 +81,9 @@ public TObject
 //
 // Initialisation of segmentation for hit  
   virtual void    SigGenInit(Float_t x, Float_t y, Float_t z);
+// Initialisation of charge fluctuation for given hit
+  virtual void    ChargeCorrelationInit();
+
 // Configuration forwarding
 //
 // Define signal distribution region
@@ -109,6 +112,7 @@ public TObject
   virtual Float_t DAlu() {return fdAlu;}  
   virtual void SetDGas(Float_t DGas) {fdGas = DGas;}
   virtual void SetDAlu(Float_t DAlu) {fdAlu = DAlu;}  
+  virtual void SetChargeCorrel(Float_t correl) {fChargeCorrel = correl;}
 // assignment operator  
   AliMUONChamber& operator =(const AliMUONChamber& rhs);
   
@@ -121,6 +125,9 @@ public TObject
   Int_t   fnsec; // number of semented cathode planes
   Float_t frMin; // innermost sensitive radius
   Float_t frMax; // outermost sensitive radius
+  Float_t fChargeCorrel; // amplitude of charge correlation on 2 cathods
+                         // is RMS of ln(q1/q2)
+  Float_t fCurrentCorrel; //! charge correlation for current hit.
 
   TObjArray              *fSegmentation;    // pointer to segmentation
   AliMUONClusterFinderVS *fReconstruction;  // pointer to reconstruction
index f2534545f45076d4a9c5b47f6d7a94e3fb682c06..21845a8fa9d543890ee001d5c7d43d346bccbc72 100644 (file)
 
 /*
 $Log$
+Revision 1.20  2000/12/04 17:48:23  gosset
+Modifications for stations 1 et 2 mainly:
+* station 1 with 4 mm gas gap and smaller cathode segmentation...
+* stations 1 and 2 with "grey" frame crosses
+* mean noise at 1.5 ADC channel
+* Ar-CO2 gas (80%+20%)
+
 Revision 1.19  2000/12/02 17:15:46  morsch
 Correction of dead zones in inner regions of stations 3-5
 Correction of length of slats 3 and 9 of station 4.
@@ -1922,7 +1929,8 @@ void AliMUONv1::StepManager()
   TLorentzVector mom;
   Float_t        theta,phi;
   Float_t        destep, step;
-  
+  Float_t        fCharge1=1;
+
   static Float_t eloss, eloss2, xhit, yhit, zhit, tof, tlength;
   const  Float_t kBig=1.e10;
   //  modifs perso
@@ -2002,6 +2010,7 @@ void AliMUONv1::StepManager()
       xhit    = pos[0];
       yhit    = pos[1];      
       zhit    = pos[2];      
+      Chamber(idvol).ChargeCorrelationInit();
       // Only if not trigger chamber
 
       
index c34211422bf1fc09863800bb534e632a8b9626d4..7882b1cd05c4d7e90cf0db89720313e3e5c524ef 100644 (file)
@@ -412,6 +412,7 @@ AliMUON *MUON  = new AliMUONv1("MUON","normal MUON");
 
 //  MUON->SetResponseModel(chamber-1, response0);          
  MUON->SetResponseModel(chamber-1, responseSt1); // special response       
+ MUON->Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
 
  chamber=2;
 //^^^^^^^^^
@@ -438,6 +439,7 @@ AliMUON *MUON  = new AliMUONv1("MUON","normal MUON");
 
 //  MUON->SetResponseModel(chamber-1, response0);          
  MUON->SetResponseModel(chamber-1, responseSt1); // special response
+ MUON->Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
            
 //
 //--------------------------------------------------------
@@ -467,6 +469,7 @@ AliMUON *MUON  = new AliMUONv1("MUON","normal MUON");
  MUON->SetSegmentationModel(chamber-1, 2, seg32);
 
  MUON->SetResponseModel(chamber-1, response0);     
+ MUON->Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
 
  chamber=4;
 //^^^^^^^^^
@@ -489,6 +492,7 @@ AliMUON *MUON  = new AliMUONv1("MUON","normal MUON");
  MUON->SetSegmentationModel(chamber-1, 2, seg42);
 
  MUON->SetResponseModel(chamber-1, response0);     
+ MUON->Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
 
 
 //--------------------------------------------------------
@@ -534,6 +538,7 @@ AliMUON *MUON  = new AliMUONv1("MUON","normal MUON");
  seg52->SetPadDivision(nseg3);
  MUON->SetSegmentationModel(chamber-1, 2, seg52);
  MUON->SetResponseModel(chamber-1, response0);     
+ MUON->Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
 
  chamber=6;
  MUON->SetNsec(chamber-1,2);
@@ -559,6 +564,7 @@ AliMUON *MUON  = new AliMUONv1("MUON","normal MUON");
  seg62->SetPadDivision(nseg3);
  MUON->SetSegmentationModel(chamber-1, 2, seg62);
  MUON->SetResponseModel(chamber-1, response0);     
+ MUON->Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
 
 //--------------------------------------------------------
 // Configuration for Chamber TC7/8  (Station 4) ----------           
@@ -610,6 +616,7 @@ AliMUON *MUON  = new AliMUONv1("MUON","normal MUON");
  seg72->SetPadDivision(nseg4);
 
  MUON->SetResponseModel(chamber-1, response0);     
+ MUON->Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
 
  chamber=8;
 //^^^^^^^^^
@@ -640,6 +647,7 @@ AliMUON *MUON  = new AliMUONv1("MUON","normal MUON");
  seg82->SetPadDivision(nseg4);
 
  MUON->SetResponseModel(chamber-1, response0);     
+ MUON->Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
 
 
 //--------------------------------------------------------
@@ -692,6 +700,7 @@ AliMUON *MUON  = new AliMUONv1("MUON","normal MUON");
  seg92->SetPadDivision(nseg4);
 
  MUON->SetResponseModel(chamber-1, response0);     
+ MUON->Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
 
  chamber=10;
 //^^^^^^^^^
@@ -722,6 +731,7 @@ AliMUON *MUON  = new AliMUONv1("MUON","normal MUON");
  seg102->SetPadDivision(nseg4);
 
  MUON->SetResponseModel(chamber-1, response0);     
+ MUON->Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
 
 //--------------------------------------------------------
 // Configuration for Trigger Stations -------------------- 
@@ -739,6 +749,7 @@ AliMUON *MUON  = new AliMUONv1("MUON","normal MUON");
  MUON->SetSegmentationModel(chamber-1, 2, seg112);
 
  MUON->SetResponseModel(chamber-1, responseTrigger0);      
+ MUON->Chamber(chamber-1).SetChargeCorrel(0); // same charge on cathodes
 
  
  chamber=12;
@@ -748,7 +759,9 @@ AliMUON *MUON  = new AliMUONv1("MUON","normal MUON");
  AliMUONSegmentationTriggerY *seg122=new AliMUONSegmentationTriggerY;
  MUON->SetSegmentationModel(chamber-1, 2, seg122);
 
- MUON->SetResponseModel(chamber-1, responseTrigger0);      
+ MUON->SetResponseModel(chamber-1, responseTrigger0);
+ MUON->Chamber(chamber-1).SetChargeCorrel(0); // same charge on cathodes
+      
   chamber=13;
  MUON->SetNsec(chamber-1,2);
  AliMUONSegmentationTriggerX *seg131=new AliMUONSegmentationTriggerX;
@@ -756,6 +769,7 @@ AliMUON *MUON  = new AliMUONv1("MUON","normal MUON");
  AliMUONSegmentationTriggerY *seg132=new AliMUONSegmentationTriggerY;
  MUON->SetSegmentationModel(chamber-1, 2, seg132);
  MUON->SetResponseModel(chamber-1, responseTrigger0);      
+ MUON->Chamber(chamber-1).SetChargeCorrel(0); // same charge on cathodes
 
  chamber=14;
  MUON->SetNsec(chamber-1,2);
@@ -765,6 +779,7 @@ AliMUON *MUON  = new AliMUONv1("MUON","normal MUON");
  MUON->SetSegmentationModel(chamber-1, 2, seg142);
 
  MUON->SetResponseModel(chamber-1, responseTrigger0); 
+ MUON->Chamber(chamber-1).SetChargeCorrel(0); // same charge on cathodes
 
 }