From 04df85d32488074f211d1381f7c51c72871f1a2e Mon Sep 17 00:00:00 2001 From: meyanne <anneclairemeyer@yahoo.fr> Date: Fri, 26 Feb 2021 21:54:52 +0100 Subject: [PATCH] ComptonTelescope (DSSSD): - Change ComptonCAM.detector file: separated blocks for dsssd / calorimeter / telescope add cartesian positions for dsssd - Keep initial ComptonCAM.detector file = ComptonCAM_ini.detector - Add position of interaction (x,y) for dsssd: in Physics and Spectra changes in ReadConfiguration to read the separated detectors blocks add AddDetectorDSSSD method to get (x,y) strips positions add impact matrix control spectra --- .../TComptonTelescopePhysics.cxx | 132 ++++++++++++++++-- .../TComptonTelescopePhysics.h | 9 +- .../TComptonTelescopeSpectra.cxx | 33 ++--- .../online/ComptonCAM.detector | 32 +++-- .../online/ComptonCAM_ini.detector | 17 +++ 5 files changed, 170 insertions(+), 53 deletions(-) create mode 100644 Projects/ComptonTelescope/online/ComptonCAM_ini.detector diff --git a/NPLib/Detectors/ComptonTelescope/TComptonTelescopePhysics.cxx b/NPLib/Detectors/ComptonTelescope/TComptonTelescopePhysics.cxx index 012fc84f2..5cfec3ed0 100644 --- a/NPLib/Detectors/ComptonTelescope/TComptonTelescopePhysics.cxx +++ b/NPLib/Detectors/ComptonTelescope/TComptonTelescopePhysics.cxx @@ -671,7 +671,62 @@ void TComptonTelescopePhysics::Clear() /////////////////////////////////////////////////////////////////////////// void TComptonTelescopePhysics::ReadConfiguration(NPL::InputParser parser){ - vector<NPL::InputBlock*> blocks = parser.GetAllBlocksWithToken("ComptonTelescope"); + // DSSSD blocks + vector<NPL::InputBlock*> blocksDSSSD = parser.GetAllBlocksWithToken("DSSSD"); + if(NPOptionManager::getInstance()->GetVerboseLevel()) + cout << "//// " << blocksDSSSD.size() << " DSSSD detectors found " << endl; + + for(unsigned int i = 0 ; i < blocksDSSSD.size() ; i++){ + + vector<string> cartDSSSD = {"X0_Y0", "X31_Y0", "X0_Y31", "X31_Y31"}; + + if(blocksDSSSD[i]->HasTokenList(cartDSSSD)){ + cout << endl << "//// DSSSD " << i+1 << endl; + TVector3 A = blocksDSSSD[i]->GetTVector3("X0_Y0", "mm"); + TVector3 B = blocksDSSSD[i]->GetTVector3("X31_Y0", "mm"); + TVector3 C = blocksDSSSD[i]->GetTVector3("X0_Y31", "mm"); + TVector3 D = blocksDSSSD[i]->GetTVector3("X31_Y31", "mm"); + + int nbr_detDSSSD = blocksDSSSD[i]->GetInt("NUMBER_DSSSD"); + double sizeDSSSD = blocksDSSSD[i]->GetDouble("SIZE_DSSSD","mm"); + double thicknessDSSSD = blocksDSSSD[i]->GetDouble("THICKNESS_DSSSD","mm"); + int nbr_stripDSSSD = blocksDSSSD[i]->GetInt("NUMBER_STRIPS"); + + AddDetectorDSSSD(A,B,C,D,sizeDSSSD,nbr_stripDSSSD); + } + else { + cout << "ERROR: Missing token for DSSSD blocks, check your input file" << endl; + exit(1); + } + } // end DSSSD blocks + + // Calorimeter blocks + vector<NPL::InputBlock*> blocksCalo = parser.GetAllBlocksWithToken("CALORIMETER"); + if(NPOptionManager::getInstance()->GetVerboseLevel()) + cout << "//// " << blocksCalo.size() << " calorimeter detectors found " << endl; + + for (unsigned int i = 0 ; i < blocksCalo.size() ; i++){ + int nbr_detCalorimeter = blocksCalo[i]->GetInt("NUMBER_CALORIMETER"); + double thickness_cal = blocksCalo[i]->GetDouble("THICKNESS_CALORIMETER","mm"); + int npixels_cal = blocksCalo[i]->GetInt("NPIXELS_CALORIMETER"); + } // end calorimeter blocks + + // ComptonTelescope blocks + vector<NPL::InputBlock*> blocksCompton = parser.GetAllBlocksWithToken("ComptonTelescope"); + if(NPOptionManager::getInstance()->GetVerboseLevel()) + cout << "//// " << blocksCompton.size() << " ComptonTelescope found " << endl; + + for (unsigned int i = 0 ; i < blocksCompton.size() ; i++){ + int nb_tracker = blocksCompton[i]->GetInt("TRACKER"); + double inter = blocksCompton[i]->GetDouble("DISTANCE_INTER_DSSSD","mm"); + double distance_cal = blocksCompton[i]->GetDouble("DISTANCE_TRACKER_CALORIMETER","mm"); + int vis= blocksCompton[i]->GetInt("VIS"); + } + + InitializeStandardParameter(); + ReadAnalysisConfig(); +} +/* vector<NPL::InputBlock*> blocks = parser.GetAllBlocksWithToken("ComptonTelescope"); if(NPOptionManager::getInstance()->GetVerboseLevel()) cout << "//// " << blocks.size() << " detectors found " << endl; @@ -696,17 +751,7 @@ void TComptonTelescopePhysics::ReadConfiguration(NPL::InputParser parser){ int vis= blocks[i]->GetInt("VIS"); AddComptonTelescope(R); } - - else{ - cout << "ERROR: check your input file formatting " << endl; - exit(1); - } - } - - InitializeStandardParameter(); - ReadAnalysisConfig(); -} - +*/ /////////////////////////////////////////////////////////////////////////// void TComptonTelescopePhysics::AddParameterToCalibrationManager() @@ -780,6 +825,64 @@ void TComptonTelescopePhysics::InitializeRootOutput() ///// Specific to ComptonTelescopeArray //// +void TComptonTelescopePhysics::AddDetectorDSSSD(TVector3 C_X0_Y0, TVector3 C_X31_Y0, TVector3 C_X0_Y31, TVector3 C_X31_Y31, double size_dsssd, int nb_strip) +{ + m_NumberOfDetectors++; + + // remove warning using C_X31_Y31 + C_X31_Y31.Unit(); + + // Vector U on Module Face (parallele to Y/Back Strip) + // NB: Y strips are allong X axis + TVector3 U = C_X31_Y0 - C_X0_Y0; + U = U.Unit(); + + // Vector V on Module Face (parallele to X Strip) + TVector3 V = C_X0_Y31 - C_X0_Y0; + V = V.Unit(); + + // Position Vector of Strip Center + TVector3 StripCenter = TVector3(0,0,0); + + // Position Vector of X=0 Y=0 strip + TVector3 Strip_0_0; + + // Buffer object to fill Position Array + vector<double> lineX; + vector<double> lineY; + vector<double> lineZ; + + vector< vector< double > > OneModuleStripPositionX; + vector< vector< double > > OneModuleStripPositionY; + vector< vector< double > > OneModuleStripPositionZ; + + // strip pitch + double stripPitch = size_dsssd/(double)nb_strip; + // Moving StripCenter to strip center of 0.0 corner + Strip_0_0 = C_X0_Y0 + (U+V) * (stripPitch/2.); + + for (int i = 0; i < nb_strip; i++) { + lineX.clear(); + lineY.clear(); + lineZ.clear(); + + for (int j = 0; j < nb_strip; j++) { + StripCenter = Strip_0_0 + stripPitch*(i*U + j*V); + lineX.push_back( StripCenter.X() ); + lineY.push_back( StripCenter.Y() ); + lineZ.push_back( StripCenter.Z() ); + } + + OneModuleStripPositionX.push_back(lineX); + OneModuleStripPositionY.push_back(lineY); + OneModuleStripPositionZ.push_back(lineZ); + } + + m_StripPositionX.push_back(OneModuleStripPositionX); + m_StripPositionY.push_back(OneModuleStripPositionY); + m_StripPositionZ.push_back(OneModuleStripPositionZ); +} + void TComptonTelescopePhysics::AddComptonTelescope(double Z) { m_NumberOfDetectors++; @@ -787,8 +890,6 @@ void TComptonTelescopePhysics::AddComptonTelescope(double Z) // needed if solid angle analysis are needed } - - TVector3 TComptonTelescopePhysics::GetDetectorNormal( const int i) const{ /* TVector3 U = TVector3 ( GetStripPositionX( DetectorNumber[i] , 24 , 1 ) , GetStripPositionY( DetectorNumber[i] , 24 , 1 ) , @@ -814,10 +915,11 @@ TVector3 TComptonTelescopePhysics::GetDetectorNormal( const int i) const{ } -TVector3 TComptonTelescopePhysics::GetPositionOfInteraction(const int i) const{ +TVector3 TComptonTelescopePhysics::GetPositionOfInteractionDSSSD(const int i) const{ TVector3 Position = TVector3 ( GetStripPositionX( DetectorNumber[i] , Strip_Front[i] , Strip_Back[i] ) , GetStripPositionY( DetectorNumber[i] , Strip_Front[i] , Strip_Back[i] ) , GetStripPositionZ( DetectorNumber[i] , Strip_Front[i] , Strip_Back[i] ) ) ; + return(Position) ; } diff --git a/NPLib/Detectors/ComptonTelescope/TComptonTelescopePhysics.h b/NPLib/Detectors/ComptonTelescope/TComptonTelescopePhysics.h index b7bb09b50..3c90dd1f9 100644 --- a/NPLib/Detectors/ComptonTelescope/TComptonTelescopePhysics.h +++ b/NPLib/Detectors/ComptonTelescope/TComptonTelescopePhysics.h @@ -140,6 +140,7 @@ class TComptonTelescopePhysics : public TObject, public NPL::VDetector void ReadAnalysisConfig(); // Add a Detector + void AddDetectorDSSSD(TVector3 C_X0_Y0, TVector3 C_X31_Y0, TVector3 C_X0_Y31, TVector3 C_X31_Y31, double size_dsssd, int nb_strip); void AddComptonTelescope(double Z); // Give and external TComptonTelescopeData object to TComptonTelescopePhysics @@ -151,9 +152,9 @@ class TComptonTelescopePhysics : public TObject, public NPL::VDetector TComptonTelescopeData* GetPreTreatedData() const {return m_PreTreatedData;} // Use to access the strip position - double GetStripPositionX(const int N, const int Front, const int Back) const {return m_StripPositionX[N-1][Front-1][Back-1];}; - double GetStripPositionY(const int N, const int Front, const int Back) const {return m_StripPositionY[N-1][Front-1][Back-1];}; - double GetStripPositionZ(const int N, const int Front, const int Back) const {return m_StripPositionZ[N-1][Front-1][Back-1];}; + double GetStripPositionX(const int N, const int Front, const int Back) const {return m_StripPositionX[N-1][Front][Back];}; + double GetStripPositionY(const int N, const int Front, const int Back) const {return m_StripPositionY[N-1][Front][Back];}; + double GetStripPositionZ(const int N, const int Front, const int Back) const {return m_StripPositionZ[N-1][Front][Back];}; double GetNumberOfDetectors() const {return m_NumberOfDetectors;}; @@ -169,7 +170,7 @@ class TComptonTelescopePhysics : public TObject, public NPL::VDetector double GetFrontTime(const int i) {return Front_Time[i];}; double GetBackTime(const int i) {return Back_Time[i];}; - TVector3 GetPositionOfInteraction(const int i) const; + TVector3 GetPositionOfInteractionDSSSD(const int i) const; TVector3 GetDetectorNormal(const int i) const; private: // Parameter used in the analysis diff --git a/NPLib/Detectors/ComptonTelescope/TComptonTelescopeSpectra.cxx b/NPLib/Detectors/ComptonTelescope/TComptonTelescopeSpectra.cxx index bf6693dc9..a44464a37 100644 --- a/NPLib/Detectors/ComptonTelescope/TComptonTelescopeSpectra.cxx +++ b/NPLib/Detectors/ComptonTelescope/TComptonTelescopeSpectra.cxx @@ -300,38 +300,30 @@ void TComptonTelescopeSpectra::InitPhysicsSpectra() AddHisto1D(name, name, 5000, -2000, 5e9, "COMPTONTELESCOPE/PHY/TIME"); // Hit pattern - name = "CT"+NPL::itoa(i+1)+"_DSSSD_HIT_PHY"; + name = "CT"+NPL::itoa(i+1)+"_DSSSD"+NPL::itoa(j+1)+"_HIT_PHY"; int ntotF = fNumberOfStripsFront * fNumberOfDetectors; int ntotB = fNumberOfStripsBack * fNumberOfDetectors; AddHisto2D(name, name, ntotF, 0, ntotF, ntotB, 0, ntotB, "COMPTONTELESCOPE/PHY/HITPATTERN"); - -/* + // X-Y Impact Matrix name = "CT"+NPL::itoa(i+1)+"_DSSSD"+NPL::itoa(j+1)+"_IMPACT_MATRIX"; - AddHisto2D(name, name,500,-150,150,500,-150,150, "COMPTONTELESCOPE/PHY/IMPACTMATRIX"); -*/ + AddHisto2D(name, name, 80, -40, 40, 80, -40, 40, "COMPTONTELESCOPE/PHY/IMPACTMATRIX"); + } - // X-Y Energy Correlation -// for (unsigned int i = 0 ; i < fNumberOfTelescope ; i++) { // loop on number of detectors + // X-Y Energy Correlation name = "CT"+NPL::itoa(i+1)+"_XY_COR"; AddHisto2D(name, name,500,0,50,500,0,50, "COMPTONTELESCOPE/PHY"); -// } - // Calorimeter energy spectrum -// for (unsigned int i = 0 ; i < fNumberOfTelescope ; i++) { // loop on number of detectors + // Calorimeter energy spectrum name = "CT"+NPL::itoa(i+1)+"_CALOR_SPECTRUM"; AddHisto1D(name, name, 1000, 1, 2000, "COMPTONTELESCOPE/PHY/CALOR"); -// } - // Position on calorimeter -// for (unsigned int i = 0 ; i < fNumberOfTelescope ; i++) { // loop on number of detectors + // Position on calorimeter name = "CT"+NPL::itoa(i+1)+"_CALOR_POS"; AddHisto2D(name, name, 8, -24, 24, 8, -24, 24, "COMPTONTELESCOPE/PHY/CALOR_POS"); -// } - // Sum spectrum -// for (unsigned int i = 0; i < fNumberOfTelescope; i++) { // loop on number of detectors + // Sum spectrum name = "CT"+NPL::itoa(i+1)+"_SUM_SPECTRUM"; AddHisto1D(name, name, 1000, 1, 2000, "COMPTONTELESCOPE/PHY/CALOR"); @@ -844,20 +836,19 @@ void TComptonTelescopeSpectra::FillPhysicsSpectra(TComptonTelescopePhysics* Phys // Hit pattern for (unsigned int i = 0; i < Physics->GetEventMultiplicity(); i++) { - name = "CT"+NPL::itoa(Physics->GetTowerNumber(i))+"_DSSSD_HIT_PHY"; + name = "CT"+NPL::itoa(Physics->GetTowerNumber(i))+"_DSSSD"+NPL::itoa(Physics->GetDetectorNumber(i))+"_HIT_PHY"; family= "COMPTONTELESCOPE/PHY/HITPATTERN"; FillSpectra(family,name,Physics->GetFrontStrip(i) + fNumberOfStripsFront*(Physics->GetDetectorNumber(i)-1), Physics->GetBackStrip(i) + fNumberOfStripsBack*(Physics->GetDetectorNumber(i)-1)); } // X-Y Impact Matrix - /* for(unsigned int i = 0 ; i < Physics->GetEventMultiplicity(); i++){ + for(unsigned int i = 0 ; i < Physics->GetEventMultiplicity(); i++){ name = "CT"+NPL::itoa(Physics->GetTowerNumber(i))+"_DSSSD"+NPL::itoa(Physics->GetDetectorNumber(i))+"_IMPACT_MATRIX"; family= "COMPTONTELESCOPE/PHY/IMPACTMATRIX"; - double x = Physics->GetPositionOfInteraction(i).x(); - double y = Physics->GetPositionOfInteraction(i).y(); + double x = Physics->GetPositionOfInteractionDSSSD(i).x(); + double y = Physics->GetPositionOfInteractionDSSSD(i).y(); FillSpectra(family,name,x,y); } -*/ //// Calorimeter diff --git a/Projects/ComptonTelescope/online/ComptonCAM.detector b/Projects/ComptonTelescope/online/ComptonCAM.detector index cdfcea260..1d61c2e51 100644 --- a/Projects/ComptonTelescope/online/ComptonCAM.detector +++ b/Projects/ComptonTelescope/online/ComptonCAM.detector @@ -1,17 +1,23 @@ -%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% +%%%%%% DSSSD 1 %%%%%% +DSSSD + NUMBER_DSSSD= 1 + X0_Y0= 32.0 32.0 0.0 + X31_Y0= -32.0 32.0 0.0 + X0_Y31= 32.0 -32.0 0.0 + X31_Y31= -32.0 -32.0 0.0 + SIZE_DSSSD= 64 mm + THICKNESS_DSSSD= 1.5 mm + NUMBER_STRIPS= 32 + +%%%%%% CALORIMETER 1 %%%%%% +CALORIMETER + NUMBER_CALORIMETER= 4 + THICKNESS_CALORIMETER= 10 mm + NPIXELS_CALORIMETER= 64 + +%%%%%% ComptonTelescope %%%%%% ComptonTelescope - R= 150 mm - THETA= 0 deg - PHI= 0 deg - BETA= 0 0 0 deg - SIZE_DSSSD= 64 mm - NUMBER_DSSSD= 1 + TRACKER= 1 DISTANCE_INTER_DSSSD= 10 mm - THICKNESS_DSSSD= 1.5 mm - NUMBER_STRIPS= 32 DISTANCE_TRACKER_CALORIMETER= 20 mm - THICKNESS_CALORIMETER= 10 mm - NPIXELS_CALORIMETER= 64 - TRACKER= 1 - CALORIMETER= 4 VIS= 1 diff --git a/Projects/ComptonTelescope/online/ComptonCAM_ini.detector b/Projects/ComptonTelescope/online/ComptonCAM_ini.detector new file mode 100644 index 000000000..cdfcea260 --- /dev/null +++ b/Projects/ComptonTelescope/online/ComptonCAM_ini.detector @@ -0,0 +1,17 @@ +%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% +ComptonTelescope + R= 150 mm + THETA= 0 deg + PHI= 0 deg + BETA= 0 0 0 deg + SIZE_DSSSD= 64 mm + NUMBER_DSSSD= 1 + DISTANCE_INTER_DSSSD= 10 mm + THICKNESS_DSSSD= 1.5 mm + NUMBER_STRIPS= 32 + DISTANCE_TRACKER_CALORIMETER= 20 mm + THICKNESS_CALORIMETER= 10 mm + NPIXELS_CALORIMETER= 64 + TRACKER= 1 + CALORIMETER= 4 + VIS= 1 -- GitLab