Skip to content
Snippets Groups Projects
Commit 83cf4f85 authored by Adrien Matta's avatar Adrien Matta :skull_crossbones:
Browse files

* npsimulation now use correctly geometry from XML file for Nebula

parent 4677314a
No related branches found
No related tags found
No related merge requests found
Pipeline #232967 passed
......@@ -116,7 +116,6 @@ void TNebulaPhysics::BuildPhysicalEvent() {
static double threshold;
// loop on Qup
for (unsigned int qup = 0; qup < QUsize ; qup++) {
rawQup = m_EventData->GetChargeUp(qup);
rawTup=-1;
rawQdown=-1;
......@@ -158,6 +157,7 @@ void TNebulaPhysics::BuildPhysicalEvent() {
}
// Got everything, do the math
if(rawTup>0){
//std::cout << "Hello 2" << std::endl;
// cal Q Up and Down
calQup=aQu[ID]*(rawQup-bQu[ID]);
calQdown=aQd[ID]*(rawQdown-bQd[ID]);
......@@ -176,7 +176,10 @@ void TNebulaPhysics::BuildPhysicalEvent() {
calTdown -= slwTd[ID]/sqrt(rawQdown-bQd[ID]);
//std::cout << calQ << " " << threshold << std::endl;
if(calQ>threshold){
//std::cout << "Hello 3" << std::endl;
calT= (calTdown+calTup)*0.5+avgT0[ID]+Cal->GetPedestal("NEBULA_T_ID"+NPL::itoa(ID));
Y=(calTdown-calTup)*DTa[ID]+DTb[ID]+Cal->GetPedestal("NEBULA_Y_ID"+NPL::itoa(ID));
......
......@@ -151,10 +151,14 @@ void Nebula::ReadConfiguration(NPL::InputParser parser){
if(NPOptionManager::getInstance()->GetVerboseLevel())
cout << "//// " << blocks.size() << " detectors found " << endl;
vector<string> cart = {"Pos","NumberOfModule","Veto","Frame"};
// define an entire wall
vector<string> wall = {"Pos","NumberOfModule","Veto","Frame"};
// use an experiment xml file to position bars individually
vector<string> xml= {"XML","Offset","InvertX","InvertY"};
for(unsigned int i = 0 ; i < blocks.size() ; i++){
if(blocks[i]->HasTokenList(cart)){
if(blocks[i]->HasTokenList(wall)){
if(NPOptionManager::getInstance()->GetVerboseLevel())
cout << endl << "//// Nebula " << i+1 << endl;
......@@ -163,6 +167,15 @@ void Nebula::ReadConfiguration(NPL::InputParser parser){
bool Veto = blocks[i]->GetInt("Veto");
bool Frame= blocks[i]->GetInt("Frame");
AddWall(Pos,NbrModule,Veto,Frame);
}
else if(blocks[i]->HasTokenList(xml)){
if(NPOptionManager::getInstance()->GetVerboseLevel())
cout << endl << "//// Nebula XML file" << i+1 << endl;
std::string xml_file = blocks[i]->GetString("XML");
G4ThreeVector Offset = NPS::ConvertVector(blocks[i]->GetTVector3("Offset","mm"));
bool InvertX = blocks[i]->GetInt("InvertX");
bool InvertY = blocks[i]->GetInt("InvertY");
ReadXML(xml_file,Offset,InvertX,InvertY);
}
else{
cout << "ERROR: check your input file formatting " << endl;
......@@ -176,15 +189,60 @@ void Nebula::ReadConfiguration(NPL::InputParser parser){
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
void Nebula::ReadXML(std::string xml_file,G4ThreeVector offset, bool InvertX,bool InvertY){
NPL::XmlParser xml;
xml.LoadFile(xml_file);
std::vector<NPL::XML::block*> b = xml.GetAllBlocksWithName("NEBULA");
int NumberOfBars=0;
for(unsigned int i = 0 ; i < b.size() ; i++){
NumberOfBars++;
unsigned int id = b[i]->AsInt("ID");
// position
auto PositionX = b[i]->AsDouble("PosX");
auto PositionY = b[i]->AsDouble("PosY");
auto PositionZ = b[i]->AsDouble("PosZ");
// SubLayer 0 is use for Veto
auto SubLayer = b[i]->AsInt("SubLayer");
// Name "NoUseX" is used to silence bars
auto nousestr = b[i]->AsString("NAME");
// Remove unused bar
if(nousestr.find("NoUse")==std::string::npos && PositionX!=PositionZ){
if(InvertX)
PositionX*=-1;
if(InvertY)
PositionY*=-1;
m_PositionBar[id]= G4ThreeVector(PositionX,PositionY,PositionZ)+offset;
if(SubLayer)
m_IsVetoBar[id]= false;
else
m_IsVetoBar[id]=true;
}
}
cout << " -> " << NumberOfBars << " bars found" << endl;
}
// Construct detector and inialise sensitive part.
// Called After DetecorConstruction::AddDetector Method
void Nebula::ConstructDetector(G4LogicalVolume* world){
// Start with XML case
G4RotationMatrix* Rot = new G4RotationMatrix();
for(auto pos : m_PositionBar){
if(!m_IsVetoBar[pos.first]){
new G4PVPlacement(G4Transform3D(*Rot,pos.second),
BuildModule(),
"NebulaModule",world,false,pos.first);
}
else{
new G4PVPlacement(G4Transform3D(*Rot,pos.second),
BuildVeto(),
"NebulaModule",world,false,pos.first);
}
}
unsigned int nbrM = 1 ;
unsigned int nbrV = 1 ;
for (unsigned short i = 0 ; i < m_Pos.size() ; i++) {
for (unsigned short m = 0 ; m < m_NbrModule[i] ; m++) {
G4RotationMatrix* Rot = new G4RotationMatrix();
double offset = (Nebula_NS::ModuleWidth+Nebula_NS::InterModule)*(-m_NbrModule[i]*0.5+m)+Nebula_NS::ModuleWidth*0.5;
G4ThreeVector Offset(offset,0,0);
new G4PVPlacement(G4Transform3D(*Rot,m_Pos[i]+Offset),
......@@ -195,7 +253,6 @@ void Nebula::ConstructDetector(G4LogicalVolume* world){
if(m_HasVeto[i]){
if(m_NbrModule[i] > 15){
for (unsigned short m = 0 ; m < Nebula_NS::VetoPerWall ; m++) {
G4RotationMatrix* Rot = new G4RotationMatrix();
double offset = (Nebula_NS::VetoWidth+Nebula_NS::InterVeto)*(-Nebula_NS::VetoPerWall*0.5+m)+Nebula_NS::VetoWidth*0.5;
G4ThreeVector Offset(offset,0,-Nebula_NS::WallToVeto);
new G4PVPlacement(G4Transform3D(*Rot,m_Pos[i]+Offset),
......@@ -205,7 +262,6 @@ void Nebula::ConstructDetector(G4LogicalVolume* world){
}
else{
for (unsigned short m = 0 ; m < Nebula_NS::VetoPerExpand ; m++) {
G4RotationMatrix* Rot = new G4RotationMatrix();
double offset = (Nebula_NS::VetoWidth+Nebula_NS::InterVeto)*(-Nebula_NS::VetoPerExpand*0.5+m)+Nebula_NS::VetoWidth*0.5;
G4ThreeVector Offset(offset,0,-Nebula_NS::WallToVeto);
new G4PVPlacement(G4Transform3D(*Rot,m_Pos[i]+Offset),
......
......@@ -24,6 +24,7 @@
// C++ header
#include <string>
#include <vector>
#include <map>
using namespace std;
// G4 headers
......@@ -36,6 +37,7 @@ using namespace std;
#include "NPSVDetector.hh"
#include "TNebulaData.h"
#include "NPInputParser.h"
#include "NPXmlParser.h"
class Nebula : public NPS::VDetector{
////////////////////////////////////////////////////
......@@ -51,11 +53,14 @@ class Nebula : public NPS::VDetector{
public:
// Cartesian
void AddWall(G4ThreeVector POS, int NbrModule,bool veto,bool frame);
void ReadXML(std::string xml_file,G4ThreeVector offset, bool InvertX,bool InvertY);
G4LogicalVolume* BuildModule();
G4LogicalVolume* BuildVeto();
private:
std::map<unsigned int , G4ThreeVector> m_PositionBar;
std::map<unsigned int , bool> m_IsVetoBar;
G4LogicalVolume* m_Module;
G4LogicalVolume* m_Veto;
double Energy;
......
......@@ -15,6 +15,7 @@ NEBULA
NumberOfModule= 15
Veto= 1
Frame= 1
XML= ./db/NEBULA.xml
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
NEBULA
POS= 0 0 12.15 m
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment