]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Add 3D function with minimal memory footprint
authorakisiel <akisiel@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 10 Mar 2011 10:19:42 +0000 (10:19 +0000)
committerakisiel <akisiel@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 10 Mar 2011 10:19:42 +0000 (10:19 +0000)
PWG2/CMakelibPWG2femtoscopy.pkg
PWG2/FEMTOSCOPY/AliFemto/AliFemtoCorrFctn3DLCMSSym.cxx [new file with mode: 0644]
PWG2/FEMTOSCOPY/AliFemto/AliFemtoCorrFctn3DLCMSSym.h [new file with mode: 0644]
PWG2/PWG2femtoscopyUserLinkDef.h

index 0bb91880ec7d86095b7211d765d9991b06f873f6..fc1e56f50ffd74312482db0dbb611da041124e98 100644 (file)
@@ -25,7 +25,7 @@
 # SHLIBS - Shared Libraries and objects for linking (Executables only)           #
 #--------------------------------------------------------------------------------#
 
-set ( SRCS  FEMTOSCOPY/AliFemto/AliFemtoSimpleAnalysis.cxx FEMTOSCOPY/AliFemto/AliFemtoLikeSignAnalysis.cxx FEMTOSCOPY/AliFemto/AliFemtoVertexAnalysis.cxx FEMTOSCOPY/AliFemto/AliFemtoVertexMultAnalysis.cxx FEMTOSCOPY/AliFemto/AliFemtoAnalysisAzimuthal.cxx FEMTOSCOPY/AliFemto/AliFemtoAnalysisReactionPlane.cxx FEMTOSCOPY/AliFemto/AliFemtoBPLCMS3DCorrFctn.cxx FEMTOSCOPY/AliFemto/AliFemtoQinvCorrFctn.cxx FEMTOSCOPY/AliFemto/AliFemtoBasicEventCut.cxx FEMTOSCOPY/AliFemto/AliFemtoBasicTrackCut.cxx FEMTOSCOPY/AliFemto/AliFemtoDummyPairCut.cxx FEMTOSCOPY/AliFemto/AliFemtoCoulomb.cxx FEMTOSCOPY/AliFemto/AliFemtoCutMonitorHandler.cxx FEMTOSCOPY/AliFemto/AliFemtoEvent.cxx FEMTOSCOPY/AliFemto/AliFemtoKink.cxx FEMTOSCOPY/AliFemto/AliFemtoManager.cxx FEMTOSCOPY/AliFemto/AliFemtoPair.cxx FEMTOSCOPY/AliFemto/AliFemtoParticle.cxx FEMTOSCOPY/AliFemto/AliFemtoPicoEvent.cxx FEMTOSCOPY/AliFemto/AliFemtoPicoEventCollectionVectorHideAway.cxx FEMTOSCOPY/AliFemto/AliFemtoTrack.cxx FEMTOSCOPY/AliFemto/AliFemtoV0.cxx FEMTOSCOPY/AliFemto/AliFemtoXi.cxx FEMTOSCOPY/AliFemto/AliFmHelix.cxx FEMTOSCOPY/AliFemto/AliFmPhysicalHelix.cxx FEMTOSCOPY/AliFemto/TpcLocalTransform.cxx FEMTOSCOPY/AliFemto/AliFemtoEventReader.cxx FEMTOSCOPY/AliFemto/AliFemtoEventReaderStandard.cxx FEMTOSCOPY/AliFemto/AliFemtoEventReaderESD.cxx FEMTOSCOPY/AliFemto/AliFemtoEventReaderESDChain.cxx FEMTOSCOPY/AliFemto/AliFemtoEventReaderESDChainKine.cxx FEMTOSCOPY/AliFemto/AliFemtoModelHiddenInfo.cxx FEMTOSCOPY/AliFemto/AliFemtoModelGlobalHiddenInfo.cxx FEMTOSCOPY/AliFemto/AliFemtoModelGausLCMSFreezeOutGenerator.cxx FEMTOSCOPY/AliFemto/AliFemtoModelWeightGenerator.cxx FEMTOSCOPY/AliFemto/AliFemtoModelWeightGeneratorBasic.cxx FEMTOSCOPY/AliFemto/AliFemtoModelManager.cxx FEMTOSCOPY/AliFemto/AliFemtoModelCorrFctn.cxx FEMTOSCOPY/AliFemto/AliFemtoModelFreezeOutGenerator.cxx FEMTOSCOPY/AliFemto/AliFemtoModelWeightGeneratorLednicky.cxx FEMTOSCOPY/AliFemto/AliFemtoCutMonitorParticleYPt.cxx FEMTOSCOPY/AliFemto/AliFemtoCutMonitorParticleVertPos.cxx FEMTOSCOPY/AliFemto/AliFemtoCutMonitorParticlePID.cxx FEMTOSCOPY/AliFemto/AliFemtoCutMonitorEventMult.cxx FEMTOSCOPY/AliFemto/AliFemtoCutMonitorEventVertex.cxx FEMTOSCOPY/AliFemto/AliFemtoKTPairCut.cxx FEMTOSCOPY/AliFemto/AliFemtoCorrFctnNonIdDR.cxx FEMTOSCOPY/AliFemto/AliFemtoCorrFctn3DSpherical.cxx FEMTOSCOPY/AliFemto/AliFemtoEventReaderAOD.cxx FEMTOSCOPY/AliFemto/AliFemtoEventReaderAODChain.cxx FEMTOSCOPY/AliFemto/AliFemtoAODTrackCut.cxx FEMTOSCOPY/AliFemto/AliFemtoCutMonitor.cxx FEMTOSCOPY/AliFemto/AliFemtoCorrFctn.cxx FEMTOSCOPY/AliFemto/AliAnalysisTaskFemto.cxx FEMTOSCOPY/AliFemto/AliTwoTrackRes.cxx)
+set ( SRCS  FEMTOSCOPY/AliFemto/AliFemtoSimpleAnalysis.cxx FEMTOSCOPY/AliFemto/AliFemtoLikeSignAnalysis.cxx FEMTOSCOPY/AliFemto/AliFemtoVertexAnalysis.cxx FEMTOSCOPY/AliFemto/AliFemtoVertexMultAnalysis.cxx FEMTOSCOPY/AliFemto/AliFemtoAnalysisAzimuthal.cxx FEMTOSCOPY/AliFemto/AliFemtoAnalysisReactionPlane.cxx FEMTOSCOPY/AliFemto/AliFemtoBPLCMS3DCorrFctn.cxx FEMTOSCOPY/AliFemto/AliFemtoCorrFctn3DLCMSSym.cxx FEMTOSCOPY/AliFemto/AliFemtoQinvCorrFctn.cxx FEMTOSCOPY/AliFemto/AliFemtoBasicEventCut.cxx FEMTOSCOPY/AliFemto/AliFemtoEventCutEstimators.cxx FEMTOSCOPY/AliFemto/AliFemtoBasicTrackCut.cxx FEMTOSCOPY/AliFemto/AliFemtoDummyPairCut.cxx FEMTOSCOPY/AliFemto/AliFemtoCoulomb.cxx FEMTOSCOPY/AliFemto/AliFemtoCutMonitorHandler.cxx FEMTOSCOPY/AliFemto/AliFemtoEvent.cxx FEMTOSCOPY/AliFemto/AliFemtoKink.cxx FEMTOSCOPY/AliFemto/AliFemtoManager.cxx FEMTOSCOPY/AliFemto/AliFemtoPair.cxx FEMTOSCOPY/AliFemto/AliFemtoParticle.cxx FEMTOSCOPY/AliFemto/AliFemtoPicoEvent.cxx FEMTOSCOPY/AliFemto/AliFemtoPicoEventCollectionVectorHideAway.cxx FEMTOSCOPY/AliFemto/AliFemtoTrack.cxx FEMTOSCOPY/AliFemto/AliFemtoV0.cxx FEMTOSCOPY/AliFemto/AliFemtoXi.cxx FEMTOSCOPY/AliFemto/AliFmHelix.cxx FEMTOSCOPY/AliFemto/AliFmPhysicalHelix.cxx FEMTOSCOPY/AliFemto/TpcLocalTransform.cxx FEMTOSCOPY/AliFemto/AliFemtoEventReader.cxx FEMTOSCOPY/AliFemto/AliFemtoEventReaderStandard.cxx FEMTOSCOPY/AliFemto/AliFemtoEventReaderESD.cxx FEMTOSCOPY/AliFemto/AliFemtoEventReaderESDChain.cxx FEMTOSCOPY/AliFemto/AliFemtoEventReaderESDChainKine.cxx FEMTOSCOPY/AliFemto/AliFemtoModelHiddenInfo.cxx FEMTOSCOPY/AliFemto/AliFemtoModelGlobalHiddenInfo.cxx FEMTOSCOPY/AliFemto/AliFemtoModelGausLCMSFreezeOutGenerator.cxx FEMTOSCOPY/AliFemto/AliFemtoModelWeightGenerator.cxx FEMTOSCOPY/AliFemto/AliFemtoModelWeightGeneratorBasic.cxx FEMTOSCOPY/AliFemto/AliFemtoModelManager.cxx FEMTOSCOPY/AliFemto/AliFemtoModelCorrFctn.cxx FEMTOSCOPY/AliFemto/AliFemtoModelFreezeOutGenerator.cxx FEMTOSCOPY/AliFemto/AliFemtoModelWeightGeneratorLednicky.cxx FEMTOSCOPY/AliFemto/AliFemtoCutMonitorParticleYPt.cxx FEMTOSCOPY/AliFemto/AliFemtoCutMonitorParticleVertPos.cxx FEMTOSCOPY/AliFemto/AliFemtoCutMonitorParticlePID.cxx FEMTOSCOPY/AliFemto/AliFemtoCutMonitorEventMult.cxx FEMTOSCOPY/AliFemto/AliFemtoCutMonitorEventVertex.cxx FEMTOSCOPY/AliFemto/AliFemtoKTPairCut.cxx FEMTOSCOPY/AliFemto/AliFemtoCorrFctnNonIdDR.cxx FEMTOSCOPY/AliFemto/AliFemtoCorrFctn3DSpherical.cxx FEMTOSCOPY/AliFemto/AliFemtoEventReaderAOD.cxx FEMTOSCOPY/AliFemto/AliFemtoEventReaderAODChain.cxx FEMTOSCOPY/AliFemto/AliFemtoAODTrackCut.cxx FEMTOSCOPY/AliFemto/AliFemtoCutMonitor.cxx FEMTOSCOPY/AliFemto/AliFemtoCorrFctn.cxx FEMTOSCOPY/AliFemto/AliAnalysisTaskFemto.cxx FEMTOSCOPY/AliFemto/AliTwoTrackRes.cxx)
 
 string ( REPLACE ".cxx" ".h" HDRS "${SRCS}" )
 
diff --git a/PWG2/FEMTOSCOPY/AliFemto/AliFemtoCorrFctn3DLCMSSym.cxx b/PWG2/FEMTOSCOPY/AliFemto/AliFemtoCorrFctn3DLCMSSym.cxx
new file mode 100644 (file)
index 0000000..465567b
--- /dev/null
@@ -0,0 +1,153 @@
+///////////////////////////////////////////////////////////////////////////
+//                                                                       //
+// AliFemtoCorrFctn3DLCMSSym: a class to calculate 3D correlation        //
+// for pairs of identical particles.                                     //
+// In analysis the function should be first created in a macro, then     //
+// added to the analysis, and at the end of the macro the procedure to   //
+// write out histograms should be called.                                //
+//                                                                       //
+///////////////////////////////////////////////////////////////////////////
+
+#include "AliFemtoCorrFctn3DLCMSSym.h"
+#include <cstdio>
+
+#ifdef __ROOT__ 
+ClassImp(AliFemtoCorrFctn3DLCMSSym)
+#endif
+
+//____________________________
+AliFemtoCorrFctn3DLCMSSym::AliFemtoCorrFctn3DLCMSSym(char* title, const int& nbins, const float& QHi)
+  :
+  AliFemtoCorrFctn(),
+  fNumerator(0),
+  fDenominator(0)
+{
+  // Basic constructor
+
+  // set up numerator
+  char tTitNum[100] = "Num";
+  strncat(tTitNum,title, 100);
+  fNumerator = new TH3F(tTitNum,title,nbins,-QHi,QHi,nbins,-QHi,QHi,nbins/2,0.0,QHi);
+  // set up denominator
+  char tTitDen[100] = "Den";
+  strncat(tTitDen,title, 100);
+  fDenominator = new TH3F(tTitDen,title,nbins,-QHi,QHi,nbins,-QHi,QHi,nbins/2,0.0,QHi);
+
+  // to enable error bar calculation...
+  fNumerator->Sumw2();
+  fDenominator->Sumw2();
+}
+
+AliFemtoCorrFctn3DLCMSSym::AliFemtoCorrFctn3DLCMSSym(const AliFemtoCorrFctn3DLCMSSym& aCorrFctn) :
+  AliFemtoCorrFctn(aCorrFctn),
+  fNumerator(0),
+  fDenominator(0)
+{
+  // Copy constructor
+  fNumerator = new TH3F(*aCorrFctn.fNumerator);
+  fDenominator = new TH3F(*aCorrFctn.fDenominator);
+}
+//____________________________
+AliFemtoCorrFctn3DLCMSSym::~AliFemtoCorrFctn3DLCMSSym(){
+  // Destructor
+  delete fNumerator;
+  delete fDenominator;
+}
+//_________________________
+AliFemtoCorrFctn3DLCMSSym& AliFemtoCorrFctn3DLCMSSym::operator=(const AliFemtoCorrFctn3DLCMSSym& aCorrFctn)
+{
+  // assignment operator
+  if (this == &aCorrFctn)
+    return *this;
+
+  if (fNumerator) delete fNumerator;
+  fNumerator = new TH3F(*aCorrFctn.fNumerator);
+  if (fDenominator) delete fDenominator;
+  fDenominator = new TH3F(*aCorrFctn.fDenominator);
+
+  return *this;
+}
+
+//_________________________
+void AliFemtoCorrFctn3DLCMSSym::WriteOutHistos(){
+  // Write out all histograms to file
+  fNumerator->Write();
+  fDenominator->Write();
+}
+//______________________________
+TList* AliFemtoCorrFctn3DLCMSSym::GetOutputList()
+{
+  // Prepare the list of objects to be written to the output
+  TList *tOutputList = new TList();
+
+  tOutputList->Add(fNumerator); 
+  tOutputList->Add(fDenominator);  
+
+  return tOutputList;
+}
+
+//_________________________
+void AliFemtoCorrFctn3DLCMSSym::Finish(){
+  // here is where we should normalize, fit, etc...
+
+}
+
+//____________________________
+AliFemtoString AliFemtoCorrFctn3DLCMSSym::Report(){
+  // Construct the report
+  string stemp = "LCMS Frame Bertsch-Pratt 3D Correlation Function Report:\n";
+  char ctemp[100];
+  snprintf(ctemp , 100, "Number of entries in numerator:\t%E\n",fNumerator->GetEntries());
+  stemp += ctemp;
+  snprintf(ctemp , 100, "Number of entries in denominator:\t%E\n",fDenominator->GetEntries());
+  stemp += ctemp;
+
+  if (fPairCut){
+    snprintf(ctemp , 100, "Here is the PairCut specific to this CorrFctn\n");
+    stemp += ctemp;
+    stemp += fPairCut->Report();
+  }
+  else{
+    snprintf(ctemp , 100, "No PairCut specific to this CorrFctn\n");
+    stemp += ctemp;
+  }
+
+  //  
+  AliFemtoString returnThis = stemp;
+  return returnThis;
+}
+//____________________________
+void AliFemtoCorrFctn3DLCMSSym::AddRealPair( AliFemtoPair* pair){
+  // perform operations on real pairs
+  if (fPairCut){
+    if (!(fPairCut->Pass(pair))) return;
+  }
+
+  double qOut = (pair->QOutCMS());
+  double qSide = (pair->QSideCMS());
+  double qLong = (pair->QLongCMS());
+
+  if (qLong > 0.0)
+    fNumerator->Fill(qOut,qSide,qLong);
+  else
+    fNumerator->Fill(-qOut,-qSide,-qLong);
+    
+}
+//____________________________
+void AliFemtoCorrFctn3DLCMSSym::AddMixedPair( AliFemtoPair* pair){
+  // perform operations on mixed pairs
+  if (fPairCut){
+    if (!(fPairCut->Pass(pair))) return;
+  }
+
+  double qOut = (pair->QOutCMS());
+  double qSide = (pair->QSideCMS());
+  double qLong = (pair->QLongCMS());
+
+  if (qLong > 0.0)
+    fDenominator->Fill(qOut,qSide,qLong,1.0);
+  else
+    fDenominator->Fill(-qOut,-qSide,-qLong,1.0);
+}
+
+
diff --git a/PWG2/FEMTOSCOPY/AliFemto/AliFemtoCorrFctn3DLCMSSym.h b/PWG2/FEMTOSCOPY/AliFemto/AliFemtoCorrFctn3DLCMSSym.h
new file mode 100644 (file)
index 0000000..561d68c
--- /dev/null
@@ -0,0 +1,49 @@
+///////////////////////////////////////////////////////////////////////////
+//                                                                       //
+// AliFemtoCorrFctn3DLCMSSym: a class to calculate 3D correlation        //
+// for pairs of identical particles vs. Bertsh-Pratt coordinates.        //
+//                                                                       //
+///////////////////////////////////////////////////////////////////////////
+
+#ifndef ALIFEMTOCORRFCTN3DLCMS_H
+#define ALIFEMTOCORRFCTN3DLCMS_H
+
+#include "AliFemtoCorrFctn.h"
+#include "AliFemtoPairCut.h"
+#include "TH3F.h"
+
+class AliFemtoCorrFctn3DLCMSSym : public AliFemtoCorrFctn {
+public:
+  AliFemtoCorrFctn3DLCMSSym(char* title, const int& nbins, const float& QHi);
+  AliFemtoCorrFctn3DLCMSSym(const AliFemtoCorrFctn3DLCMSSym& aCorrFctn);
+  virtual ~AliFemtoCorrFctn3DLCMSSym();
+
+  AliFemtoCorrFctn3DLCMSSym& operator=(const AliFemtoCorrFctn3DLCMSSym& aCorrFctn);
+
+  virtual AliFemtoString Report();
+  virtual void AddRealPair( AliFemtoPair* aPair);
+  virtual void AddMixedPair( AliFemtoPair* aPair);
+
+  virtual void Finish();
+
+  TH3F* Numerator();
+  TH3F* Denominator();
+
+  void WriteOutHistos();
+  virtual TList* GetOutputList();
+
+private:
+
+  TH3F* fNumerator;         // numerator
+  TH3F* fDenominator;       // denominator
+
+#ifdef __ROOT__
+  ClassDef(AliFemtoCorrFctn3DLCMSSym, 1)
+#endif
+};
+
+inline  TH3F* AliFemtoCorrFctn3DLCMSSym::Numerator(){return fNumerator;}
+inline  TH3F* AliFemtoCorrFctn3DLCMSSym::Denominator(){return fDenominator;}
+
+#endif
+
index 2267a4407df589769c9c5c4f95757acaa82c5512..14f500e36fc836cdc30e1c06ef7bc8b84f016967 100644 (file)
@@ -5,6 +5,7 @@
 #pragma link C++ class AliFemtoShareQualityPairCut;
 #pragma link C++ class AliFemtoShareQualityKTPairCut;
 #pragma link C++ class AliFemtoShareQualityTPCEntranceSepPairCut;
+#pragma link C++ class AliFemtoPairCutRadialDistance;
 #pragma link C++ class AliFemtoESDTrackCut;
 #pragma link C++ class AliFemtoShareQualityCorrFctn;
 #pragma link C++ class AliFemtoTPCInnerCorrFctn;