/////////////////////////////////////////////////////////////////////////////////////// /// /// FILENAME: PlotPMTShadowing /// /// BRIEF: ROOT macro to demonstrate the capability of the /// PMT shadowing utility; RAT::DU::ShadowingCalculator /// /// AUTHOR: Rob Stainforth [RPFS] /// /// REVISION HISTORY: /// 06/2015 : RPFS - First Revision, new file. /// /////////////////////////////////////////////////////////////////////////////////////// #include "RAT/DU/Utility.hh" #include "RAT/DU/LightPathCalculator.hh" #include "RAT/DU/ShadowingCalculator.hh" #include "RAT/DU/PMTInfo.hh" #include "RAT/Log.hh" #include #include "TVector3.h" #include "TCanvas.h" #include "TGraph.h" #include "TAxis.h" #include "TStyle.h" #include "TMath.h" #include "TLegend.h" #include "TPaveText.h" #include #include #include #include #include #include using namespace std; using namespace RAT; using namespace RAT::DU; /// Plot the PMTs shadowed. /// /// @param[in] xPos - x-coordinate (mm) of the light path origin /// @param[in] yPos - y-coordinate (mm) of the light path origin /// @param[in] zPos - z-coordinate (mm) of the light path origin /// @param[in] wavelength - Wavelength in nm of the light /// @param[in] localityVal - Proximity required for the light path calculator to calculate /// a refracted path to each PMT (set to = 0.0 for straight paths). TGraph* PlotShadowedPMTs( const Double_t xPos, const Double_t yPos, const Double_t zPos, const Double_t wavelength, const Double_t localityVal ) { // Required number of arguments is 6 xPos yPos zPos wavelength localityValue-mm // Print out the passed information. info << dformat( "Light Source Position: ( %.2f, %.2f, %.2f ) mm\nWavelength: %.2f nm\nLocality Threshold: %.2f mm\n", xPos, yPos, zPos, wavelength, localityVal ); // If this is being done on data that does not require remote database connection // eg.: a simple simulation with default run number (0) // We can disable the remote connections: // // NOTE: Don't do this if you are using real data!!! RAT::DB::Get()->SetAirplaneModeStatus(true); // First load all the utilities: PMT information, lightpath calculator // and the shadowing calculator. RAT::DU::Utility::Get()->LoadDBAndBeginRun(); PMTInfo pmtInfo = Utility::Get()->GetPMTInfo(); LightPathCalculator lightPath = Utility::Get()->GetLightPathCalculator(); ShadowingCalculator shadowCalc = Utility::Get()->GetShadowingCalculator(); // Declare the vectors for the start and ending positions of the light path // Start position as declared by the passed (x,y,z) arguments TVector3 startPos( xPos, yPos, zPos ); // End position is initialied below when reading through the PMT database TVector3 pmtPos( 0.0, 0.0, 0.0 ); // The equivelent energy (MeV) of the wavelength (nm) provided at the command line // The 1.0e-6 is the conversion for nm -> mm (mm being the offical RAT units of length) Double_t laserE = RAT::util::WavelengthToEnergy( wavelength * 1.0e-6 ); // Define the graph to plot the PMT positions shadowed. TGraph* shadowedPMTs = new TGraph(); Int_t graphPoint = 0; // Counters to see how many PMTs there are in total, // how many are normal type PMTs // and how many are shadowed PMTs Int_t nPMTs = 0; Int_t nNormPMTs = 0; Int_t nNormShPMTs = 0; // Loop over all the PMTs, check to see if they are of a 'normal' type, // if so: check for shadowing // Set all of the tolerances for shadowing by the detector geometry to 150 mm equivalent. shadowCalc.SetAllGeometryTolerances( 150.0 ); for ( Int_t iPMT = 0; iPMT < (Int_t)pmtInfo.GetCount(); iPMT++ ){ nPMTs++; // Check that the PMT is of normal type if ( pmtInfo.GetType( iPMT ) != 1 ){ continue; } else{ nNormPMTs++; // Get the PMT position pmtPos = pmtInfo.GetPosition( iPMT ); // Calculate the light path from the given start location to the PMT location lightPath.CalcByPosition( startPos, pmtPos, laserE, localityVal ); // Check for shadowing if ( shadowCalc.CheckForShadowing( lightPath ) == true ){ shadowedPMTs->SetPoint( graphPoint++, ( pmtPos.Phi() * ( 180.0 / TMath::Pi() ) ), pmtPos.CosTheta() ); } } lightPath.ResetValues(); } TCanvas* c1 = new TCanvas( "pmtShadowing", "PMT Shadowing", 900, 600 ); shadowedPMTs->GetXaxis()->SetTitle( "PMT Position #phi [degrees]" ); shadowedPMTs->GetYaxis()->SetTitle( "PMT Position Cos(#theta)" ); char str[1000]; sprintf( str, "Shadowed PMTs from Light Source Position ( %.1f, %.1f, %.1f ) cm", startPos.X() / 10.0, startPos.Y() / 10.0, startPos.Z() / 10.0 ); shadowedPMTs->SetTitle( str ); gStyle->SetTitleFontSize( 0.05 ); shadowedPMTs->SetMarkerColor( 1 ); shadowedPMTs->SetMarkerStyle( 7 ); shadowedPMTs->Draw( "AP" ); c1->SetMargin( 0.1, 0.1, 0.1, 0.2 ); info << dformat( "----------------------\n" ); info << dformat( "Number of PMTs: %i\n", nPMTs ); info << dformat( "Number of Normal PMTs: %i\n", nNormPMTs ); info << dformat( "Number of Shadowed Normal PMTs: %i\n", nNormShPMTs ); info << dformat( "----------------------\n" ); return shadowedPMTs; } /////////////////////////////////////// /////////////////////////////////////// /// Plot the PMTs shadowed based on geometry type. /// /// @param[in] xPos - x-coordinate (mm) of the light path origin /// @param[in] yPos - y-coordinate (mm) of the light path origin /// @param[in] zPos - z-coordinate (mm) of the light path origin /// @param[in] wavelength - Wavelength in nm of the light /// @param[in] localityVal - Proximity required for the light path calculator to calculate /// a refracted path to each PMT (set to = 0.0 for straight paths). TGraph* PlotPMTShadowingByGeometry( const Double_t xPos, const Double_t yPos, const Double_t zPos, const Double_t wavelength, const Double_t localityVal ) { // Required number of arguments is 6 xPos yPos zPos wavelength localityValue-mm // Print out the passed information. info << dformat( "Light Source Position: ( %.2f, %.2f, %.2f ) mm\nWavelength: %.2f nm\nLocality Threshold: %.2f mm\n", xPos, yPos, zPos, wavelength, localityVal ); // If this is being done on data that does not require remote database connection // eg.: a simple simulation with default run number (0) // We can disable the remote connections: // // NOTE: Don't do this if you are using real data!!! RAT::DB::Get()->SetAirplaneModeStatus(true); // First load all the utilities: PMT information, lightpath calculator // and the shadowing calculator. RAT::DU::Utility::Get()->LoadDBAndBeginRun(); PMTInfo pmtInfo = Utility::Get()->GetPMTInfo(); LightPathCalculator lightPath = Utility::Get()->GetLightPathCalculator(); ShadowingCalculator shadowCalc = Utility::Get()->GetShadowingCalculator(); // Set the tolerances for shadowing by the individual detector geometries... shadowCalc.SetAVHDRopeTolerance( 150.0 ); // AV hold-down rope tolerance: 150.0 mm shadowCalc.SetSupportRopeTolerance( 200.0 ); // Support hold-up rope tolerance: 200.0 mm shadowCalc.SetBellyPlateTolerance( 100.0 ); // Belly plate tolerance: 100.0 mm shadowCalc.SetNCDAnchorTolerance( 300.0 ); // NCD anchor tolerance: 300.0 mm shadowCalc.SetAVPipeTolerance( 200.0 ); // AV pipe tolerance: 200.0 mm shadowCalc.SetNeckBossTolerance( 200.0 ); // Neck boss tolerance: 200.0 mm // Declare the vectors for the start and ending positions of the light path // Start position as declared by the passed (x,y,z) arguments TVector3 startPos( xPos, yPos, zPos ); // End position is initialied below when reading through the PMT database TVector3 pmtPos( 0.0, 0.0, 0.0 ); // The equivelent energy (MeV) of the wavelength (nm) provided at the command line // The 1.0e-6 is the conversion for nm -> mm (mm being the offical RAT units of length) Double_t laserE = RAT::util::WavelengthToEnergy( wavelength * 1.0e-6 ); // Define the graphs to plot the PMT positions shadowed due to // various detector geometry. TGraph* zPMTs = new TGraph(); // Neck PMTs (high-z) TGraph* bpPMTs = new TGraph(); // Belly Plate PMTs TGraph* avhdPMTs = new TGraph(); // AV Hold Down Rope PMTs TGraph* srPMTs = new TGraph(); // Support Rope PMTs TGraph* ncdPMTs = new TGraph(); // NCD Anchor PMTs TGraph* pipePMTs = new TGraph(); // Pipe PMTs TGraph* nbossPMTs = new TGraph(); // Neck Boss PMTs // TGraph 'SetPoint' counters Int_t zP = 0; // Neck PMTs (high-z) Int_t bpP = 0; // Belly Plate PMTs Int_t avhdP = 0; // AV Hold Down Rope PMTs Int_t srP = 0; // Support Rope PMTs Int_t ncdP = 0; // NCD Anchor PMTs Int_t pipeP = 0; // Pipe PMTs Int_t nbossP = 0; // Neck Boss PMTs // Set these such that if only 'zPMTs' are plotted, they occupy the whole // Cos(theta), Phi region zPMTs->SetPoint( zP++, 0.0, 1.0 ); zPMTs->SetPoint( zP++, 0.0, -1.0 ); zPMTs->SetPoint( zP++, -180.0, 0.0 ); zPMTs->SetPoint( zP++, 180.0, 0.0 ); // Counters to see how many PMTs there are in total, // how many are normal type PMTs // and how many are shadowed PMTs Int_t nPMTs = 0; Int_t nNormPMTs = 0; Int_t nNormShPMTs = 0; // Loop over all the PMTs, check to see if they are of a 'normal' type, // if so: check for shadowing for ( Int_t iPMT = 0; iPMT < (Int_t)pmtInfo.GetCount(); iPMT++ ){ nPMTs++; // Check that the PMT is of normal type if ( pmtInfo.GetType( iPMT ) != 1 ){ continue; } else{ nNormPMTs++; // Get the PMT position pmtPos = pmtInfo.GetPosition( iPMT ); // Calculate the light path from the given start location to the PMT location lightPath.CalcByPosition( startPos, pmtPos, laserE, localityVal ); // Check for shadowing if ( shadowCalc.CheckForShadowing( lightPath ) == true ){ nNormShPMTs++; // Check for belly plate shadowing if ( shadowCalc.GetBellyPlateShadowFlag() == true ){ bpPMTs->SetPoint( bpP++, ( pmtPos.Phi() * ( 180.0 / TMath::Pi() ) ), pmtPos.CosTheta() ); } // Check for AV hold-down rope shadowing if ( shadowCalc.GetAVHDRopeShadowFlag() == true ){ avhdPMTs->SetPoint( avhdP++, ( pmtPos.Phi() * ( 180.0 / TMath::Pi() ) ), pmtPos.CosTheta() ); } // Check for support rope shadowing if ( shadowCalc.GetSupportRopeShadowFlag() == true ){ srPMTs->SetPoint( srP++, ( pmtPos.Phi() * ( 180.0 / TMath::Pi() ) ), pmtPos.CosTheta() ); } // Check for NCD anchor shadowing if ( shadowCalc.GetNCDAnchorShadowFlag() == true ){ ncdPMTs->SetPoint( ncdP++, ( pmtPos.Phi() * ( 180.0 / TMath::Pi() ) ), pmtPos.CosTheta() ); } // Check for AV pipe shadowing if ( shadowCalc.GetAVPipeShadowFlag() == true ){ pipePMTs->SetPoint( pipeP++, ( pmtPos.Phi() * ( 180.0 / TMath::Pi() ) ), pmtPos.CosTheta() ); } // Check for neck boss shadowing if ( shadowCalc.GetNeckBossShadowFlag() == true ){ nbossPMTs->SetPoint( nbossP++, ( pmtPos.Phi() * ( 180.0 / TMath::Pi() ) ), pmtPos.CosTheta() ); } } } lightPath.ResetValues(); } // What follows below is simply plotting etiquette TCanvas* c1 = new TCanvas( "pmtShadowing", "PMT Shadowing", 900, 600 ); zPMTs->GetXaxis()->SetTitle( "PMT Position #phi [degrees]" ); zPMTs->GetYaxis()->SetTitle( "PMT Position Cos(#theta)" ); zPMTs->SetMarkerColor( 0 ); zPMTs->SetMarkerStyle( 1 ); zPMTs->SetLineColor( 0 ); zPMTs->SetLineWidth( 2 ); zPMTs->SetFillColor( 0 ); bpPMTs->SetMarkerColor( 1 ); bpPMTs->SetMarkerStyle( 7 ); bpPMTs->SetLineColor( 1 ); bpPMTs->SetLineWidth( 2 ); bpPMTs->SetFillColor( 1 ); avhdPMTs->SetMarkerColor( 2 ); avhdPMTs->SetMarkerStyle( 7 ); avhdPMTs->SetLineColor( 2 ); avhdPMTs->SetLineWidth( 2 ); avhdPMTs->SetFillColor( 2 ); srPMTs->SetMarkerColor( 3 ); srPMTs->SetMarkerStyle( 7 ); srPMTs->SetLineColor( 3 ); srPMTs->SetLineWidth( 2 ); srPMTs->SetFillColor( 3 ); ncdPMTs->SetMarkerColor( 4 ); ncdPMTs->SetMarkerStyle( 7 ); ncdPMTs->SetLineColor( 4 ); ncdPMTs->SetLineWidth( 2 ); ncdPMTs->SetFillColor( 4 ); pipePMTs->SetMarkerColor( 6 ); pipePMTs->SetMarkerStyle( 7 ); pipePMTs->SetLineColor( 6 ); pipePMTs->SetLineWidth( 2 ); pipePMTs->SetFillColor( 6 ); nbossPMTs->SetMarkerColor( 7 ); nbossPMTs->SetMarkerStyle( 7 ); nbossPMTs->SetLineColor( 7 ); nbossPMTs->SetLineWidth( 2 ); nbossPMTs->SetFillColor( 7 ); zPMTs->Draw( "AP" ); bpPMTs->Draw( "same,P" ); avhdPMTs->Draw( "same,P" ); srPMTs->Draw( "same,P" ); ncdPMTs->Draw( "same,P" ); pipePMTs->Draw( "same,P" ); nbossPMTs->Draw( "same,P" ); char str[1000]; sprintf( str, "Shadowed PMTs from Light Source Position ( %.1f, %.1f, %.1f ) cm", startPos.X() / 10.0, startPos.Y() / 10.0, startPos.Z() / 10.0 ); zPMTs->SetTitle( str ); gStyle->SetTitleFontSize( 0.05 ); stringstream myStream; string myString = ""; /////////////////////////////////////////////////// /////////////////////////////////////////////////// TLegend* leg1 = new TLegend( 0.1, 0.80, 0.366666666, 0.88 ); myStream << setprecision( 2 ) << (double)( 100.0 * (double)bpP / nNormPMTs ); myStream >> myString; leg1->AddEntry( bpPMTs, ((string)"Belly Plates [" + myString + " %]").c_str(), "f" ); myStream.clear(); myString.clear(); myStream << setprecision( 2 ) << (double)( 100.0 * (double)avhdP / nNormPMTs ); myStream >> myString; leg1->AddEntry( avhdPMTs, ((string)"Hold-Down Ropes [" + myString + " %]").c_str(), "f" ); myStream.clear(); myString.clear(); leg1->Draw(); leg1->SetTextSize( 0.025 ); leg1->SetTextAlign( 22 ); /////////////////////////////////////////////////// /////////////////////////////////////////////////// TLegend* leg2 = new TLegend( 0.366666666, 0.80, 0.633333333, 0.88 ); myStream << setprecision( 2 ) << (double)( 100.0 * (double)srP / nNormPMTs ); myStream >> myString; leg2->AddEntry( srPMTs, ((string)"Support Ropes [" + myString + " %]").c_str(), "f" ); myStream.clear(); myString.clear(); myStream << setprecision( 2 ) << (double)( 100.0 * (double)ncdP / nNormPMTs ); myStream >> myString; leg2->AddEntry( ncdPMTs, ((string)"NCD Anchors [" + myString + " %]").c_str(), "f" ); myStream.clear(); myString.clear(); leg2->Draw(); leg2->SetTextSize( 0.025 ); leg2->SetTextAlign( 22 ); /////////////////////////////////////////////////// /////////////////////////////////////////////////// TLegend* leg3 = new TLegend( 0.633333333, 0.80, 0.90, 0.88 ); myStream << setprecision( 2 ) << (double)( 100.0 * (double)pipeP / nNormPMTs ); myStream >> myString; leg3->AddEntry( pipePMTs, ((string)"AV Pipes [" + myString + " %]").c_str(), "f" ); myStream.clear(); myString.clear(); myStream << setprecision( 2 ) << (double)( 100.0 * (double)nbossP / nNormPMTs ); myStream >> myString; leg3->AddEntry( nbossPMTs, ((string)"Neck Boss [" + myString + " %]").c_str(), "f" ); myStream.clear(); myString.clear(); leg3->Draw(); leg3->SetTextSize( 0.025 ); leg3->SetTextAlign( 22 ); /////////////////////////////////////////////////// /////////////////////////////////////////////////// myStream << setprecision( 2 ) << (double)( 100.0 * (double)nNormShPMTs / nNormPMTs ); myStream >> myString; TPaveText* totSh = new TPaveText( 0.366666666, 1.55, 0.633333333, 1.59 ); totSh->AddText( ("Total PMT Shadowing: " + myString + " %").c_str() ); totSh->Draw( "same" ); totSh->SetTextSize( 0.025 ); totSh->SetTextAlign( 22 ); myStream.clear(); myString.clear(); /////////////////////////////////////////////////// /////////////////////////////////////////////////// c1->SetMargin( 0.1, 0.1, 0.1, 0.2 ); info << dformat( "Number of Belly Plate PMTs: %i\n", bpP ); info << dformat( "Number of Support Rope PMTs: %i\n", srP ); info << dformat( "Number of AVHD Rope PMTs: %i\n", avhdP ); info << dformat( "Number of AV Pipe PMTs: %i\n", pipeP ); info << dformat( "Number of NCD Anchor PMTs: %i\n", ncdP ); info << dformat( "----------------------\n" ); info << dformat( "Number of PMTs: %i\n", nPMTs ); info << dformat( "Number of Normal PMTs: %i\n", nNormPMTs ); info << dformat( "Number of Shadowed Normal PMTs: %i\n", nNormShPMTs ); info << dformat( "----------------------\n" ); info << dformat( "Number of Shadowed PMTs: %i (%.2f %%)\n", nNormShPMTs, 100.0 * (Double_t)nNormShPMTs / nNormPMTs ); info << dformat( "Belly Plate PMTs: %.2f %%\n", 100.0 * (Double_t)bpP / nNormPMTs ); info << dformat( "Support Rope PMTs: %.2f %%\n", 100.0 * (Double_t)srP / nNormPMTs ); info << dformat( "AVHD Rope PMTs: %.2f %%\n", 100.0 * (Double_t)avhdP / nNormPMTs ); info << dformat( "AV Pipe PMTs: %.2f %%\n", 100.0 * (Double_t)pipeP / nNormPMTs ); info << dformat( "NCD Anchor PMTs: %.2f %%\n", 100.0 * (Double_t)ncdP / nNormPMTs ); info << dformat( "Neck Boss PMTs: %.2f %%\n", 100.0 * (Double_t)nbossP / nNormPMTs ); return zPMTs; }