All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
spacecharge::SpaceChargeSBND Class Reference

#include <SpaceChargeSBND.h>

Inheritance diagram for spacecharge::SpaceChargeSBND:
spacecharge::SpaceCharge

Public Member Functions

 SpaceChargeSBND (fhicl::ParameterSet const &pset)
 
 SpaceChargeSBND (SpaceChargeSBND const &)=delete
 
virtual ~SpaceChargeSBND ()=default
 
bool Configure (fhicl::ParameterSet const &pset)
 
bool Update (uint64_t ts=0)
 
bool EnableSimSpatialSCE () const override
 
bool EnableSimEfieldSCE () const override
 
bool EnableCalSpatialSCE () const override
 
bool EnableCalEfieldSCE () const override
 
bool EnableCorrSCE () const override
 
geo::Vector_t GetPosOffsets (geo::Point_t const &point) const override
 
geo::Vector_t GetEfieldOffsets (geo::Point_t const &point) const override
 
geo::Vector_t GetCalPosOffsets (geo::Point_t const &point, int const &TPCid=1) const override
 
geo::Vector_t GetCalEfieldOffsets (geo::Point_t const &point, int const &TPCid=1) const override
 
- Public Member Functions inherited from spacecharge::SpaceCharge
 SpaceCharge (const SpaceCharge &)=delete
 
 SpaceCharge (SpaceCharge &&)=delete
 
SpaceChargeoperator= (const SpaceCharge &)=delete
 
SpaceChargeoperator= (SpaceCharge &&)=delete
 
virtual ~SpaceCharge ()=default
 

Protected Member Functions

std::vector< double > GetPosOffsetsParametric (double xVal, double yVal, double zVal) const
 
double GetOnePosOffsetParametric (double xVal, double yVal, double zVal, std::string axis) const
 
std::vector< double > GetEfieldOffsetsParametric (double xVal, double yVal, double zVal) const
 
double GetOneEfieldOffsetParametric (double xVal, double yVal, double zVal, std::string axis) const
 
double TransformX (double xVal) const
 
double TransformY (double yVal) const
 
double TransformZ (double zVal) const
 
bool IsInsideBoundaries (double xVal, double yVal, double zVal) const
 
- Protected Member Functions inherited from spacecharge::SpaceCharge
 SpaceCharge ()=default
 

Protected Attributes

int initialSpatialFitPolN [3] = {3, 4, 3}
 
int intermediateSpatialFitPolN [3] = {4, 4, 4}
 
int initialEFieldFitPolN [3] = {3, 3, 3}
 
int intermediateEFieldFitPolN [3] = {6, 4, 4}
 
double DriftField = 500.0
 
bool fEnableSimSpatialSCE
 
bool fEnableSimEfieldSCE
 
bool fEnableCalSpatialSCE
 
bool fEnableCalEfieldSCE
 
bool fEnableCorrSCE
 
std::string fRepresentationType
 
std::string fInputFilename
 
std::vector< TH3F * > SCEhistograms = std::vector<TH3F*>(9)
 
TGraph * gSpatialGraphX [99][99]
 
TF1 * intermediateSpatialFitFunctionX [99]
 
TF1 * initialSpatialFitFunctionX
 
TGraph * gSpatialGraphY [99][99]
 
TF1 * intermediateSpatialFitFunctionY [99]
 
TF1 * initialSpatialFitFunctionY
 
TGraph * gSpatialGraphZ [99][99]
 
TF1 * intermediateSpatialFitFunctionZ [99]
 
TF1 * initialSpatialFitFunctionZ
 
TGraph * gEFieldGraphX [99][99]
 
TF1 * intermediateEFieldFitFunctionX [99]
 
TF1 * initialEFieldFitFunctionX
 
TGraph * gEFieldGraphY [99][99]
 
TF1 * intermediateEFieldFitFunctionY [99]
 
TF1 * initialEFieldFitFunctionY
 
TGraph * gEFieldGraphZ [99][99]
 
TF1 * intermediateEFieldFitFunctionZ [99]
 
TF1 * initialEFieldFitFunctionZ
 

Detailed Description

Definition at line 20 of file SpaceChargeSBND.h.

Constructor & Destructor Documentation

spacecharge::SpaceChargeSBND::SpaceChargeSBND ( fhicl::ParameterSet const &  pset)
explicit

Definition at line 19 of file SpaceChargeSBND.cxx.

20 {
21  Configure(pset);
22 }
bool Configure(fhicl::ParameterSet const &pset)
spacecharge::SpaceChargeSBND::SpaceChargeSBND ( SpaceChargeSBND const &  )
delete
virtual spacecharge::SpaceChargeSBND::~SpaceChargeSBND ( )
virtualdefault

Member Function Documentation

bool spacecharge::SpaceChargeSBND::Configure ( fhicl::ParameterSet const &  pset)

Definition at line 24 of file SpaceChargeSBND.cxx.

25 {
26  fEnableSimSpatialSCE = pset.get<bool>("EnableSimSpatialSCE");
27  fEnableSimEfieldSCE = pset.get<bool>("EnableSimEfieldSCE");
28  //fEnableCorrSCE = pset.get<bool>("EnableCorrSCE");
29  fEnableCalSpatialSCE = pset.get<bool>("EnableCalSpatialSCE");
30  fEnableCalEfieldSCE = pset.get<bool>("EnableCalEfieldSCE");
31 
32  if((fEnableSimSpatialSCE == true) || (fEnableSimEfieldSCE == true))
33  {
34  fRepresentationType = pset.get<std::string>("RepresentationType");
35  fInputFilename = pset.get<std::string>("InputFilename");
36 
37  std::string fname;
38  cet::search_path sp("FW_SEARCH_PATH");
39  sp.find_file(fInputFilename, fname);
40 
41  std::unique_ptr<TFile> infile(new TFile(fname.c_str(), "READ"));
42  if(!infile->IsOpen())
43  {
44  throw cet::exception("SpaceChargeSBND") << "Could not find the space charge effect file '" << fname << "'!\n";
45  }
46 
47  if(fRepresentationType == "Voxelized_TH3"){
48  std::cout << "begin loading voxelized TH3s..." << std::endl;
49 
50  //Load in histograms
51  TH3F* hTrueFwdX = (TH3F*) infile->Get("TrueFwd_Displacement_X");
52  TH3F* hTrueFwdY = (TH3F*) infile->Get("TrueFwd_Displacement_Y");
53  TH3F* hTrueFwdZ = (TH3F*) infile->Get("TrueFwd_Displacement_Z");
54  TH3F* hTrueBkwdX = (TH3F*) infile->Get("TrueBkwd_Displacement_X");
55  TH3F* hTrueBkwdY = (TH3F*) infile->Get("TrueBkwd_Displacement_Y");
56  TH3F* hTrueBkwdZ = (TH3F*) infile->Get("TrueBkwd_Displacement_Z");
57  TH3F* hTrueEFieldX = (TH3F*) infile->Get("True_ElecField_X");
58  TH3F* hTrueEFieldY = (TH3F*) infile->Get("True_ElecField_Y");
59  TH3F* hTrueEFieldZ = (TH3F*) infile->Get("True_ElecField_Z");
60 
61  //https://root.cern.ch/doc/master/classTH1.html#a0367fe04ae8709fd4b82795d0a5462c3
62  //Set hist directories so they can be referenced elsewhere
63  //This needs to be done because they were read in from ext file
64  //Note this is not a property of the TH3F, so does't survive copying
65  hTrueFwdX->SetDirectory(0);
66  hTrueFwdY->SetDirectory(0);
67  hTrueFwdZ->SetDirectory(0);
68  hTrueBkwdX->SetDirectory(0);
69  hTrueBkwdY->SetDirectory(0);
70  hTrueBkwdZ->SetDirectory(0);
71  hTrueEFieldX->SetDirectory(0);
72  hTrueEFieldY->SetDirectory(0);
73  hTrueEFieldZ->SetDirectory(0);
74 
75  //SCEhistograms can be accessed globally in this script
76  SCEhistograms = {hTrueFwdX, hTrueFwdY, hTrueFwdZ,
77  hTrueBkwdX, hTrueBkwdY, hTrueBkwdZ,
78  hTrueEFieldX, hTrueEFieldY, hTrueEFieldZ};
79 
80 
81  std::cout << "...finished loading TH3s" << std::endl;
82  }else if(fRepresentationType == "Parametric")
83  {
84  for(int i = 0; i < initialSpatialFitPolN[0] + 1; i++)
85  {
86  for(int j = 0; j < intermediateSpatialFitPolN[0] + 1; j++)
87  {
88  gSpatialGraphX[i][j] = (TGraph*)infile->Get(Form("deltaX/g%i_%i", i, j));
89 
90  }
91  intermediateSpatialFitFunctionX[i] = new TF1(Form("intermediateSpatialFitFunctionX_%i", i), Form("pol%i", intermediateSpatialFitPolN[0]));
92  }
93  for(int i = 0; i < initialSpatialFitPolN[1] + 1; i++)
94  {
95  for(int j = 0; j < intermediateSpatialFitPolN[1] + 1; j++)
96  {
97  gSpatialGraphY[i][j] = (TGraph*)infile->Get(Form("deltaY/g%i_%i", i, j));
98 
99  }
100  intermediateSpatialFitFunctionY[i] = new TF1(Form("intermediateSpatialFitFunctionY_%i", i), Form("pol%i", intermediateSpatialFitPolN[1]));
101  }
102  for(int i = 0; i < initialSpatialFitPolN[2] + 1; i++)
103  {
104  for(int j = 0; j < intermediateSpatialFitPolN[2] + 1; j++)
105  {
106  gSpatialGraphZ[i][j] = (TGraph*)infile->Get(Form("deltaZ/g%i_%i", i, j));
107 
108  }
109  intermediateSpatialFitFunctionZ[i] = new TF1(Form("intermediateSpatialFitFunctionZ_%i", i), Form("pol%i", intermediateSpatialFitPolN[2]));
110  }
111 
112  initialSpatialFitFunctionX = new TF1("initialSpatialFitFunctionX", Form("pol%i", initialSpatialFitPolN[0]));
113  initialSpatialFitFunctionY = new TF1("initialSpatialFitFunctionY", Form("pol%i", initialSpatialFitPolN[1]));
114  initialSpatialFitFunctionZ = new TF1("initialSpatialFitFunctionZ", Form("pol%i", initialSpatialFitPolN[2]));
115 
116  for(int i = 0; i < initialEFieldFitPolN[0] + 1; i++)
117  {
118  for(int j = 0; j < intermediateEFieldFitPolN[0] + 1; j++)
119  {
120  gEFieldGraphX[i][j] = (TGraph*)infile->Get(Form("deltaEx/g%i_%i", i, j));
121 
122  }
123  intermediateEFieldFitFunctionX[i] = new TF1(Form("intermediateEFieldFitFunctionX_%i", i), Form("pol%i", intermediateEFieldFitPolN[0]));
124  }
125  for(int i = 0; i < initialEFieldFitPolN[1] + 1; i++)
126  {
127  for(int j = 0; j < intermediateEFieldFitPolN[1] + 1; j++)
128  {
129  gEFieldGraphY[i][j] = (TGraph*)infile->Get(Form("deltaEy/g%i_%i", i, j));
130 
131  }
132  intermediateEFieldFitFunctionY[i] = new TF1(Form("intermediateEFieldFitFunctionY_%i", i), Form("pol%i", intermediateEFieldFitPolN[1]));
133  }
134  for(int i = 0; i < initialEFieldFitPolN[2] + 1; i++)
135  {
136  for(int j = 0; j < intermediateEFieldFitPolN[2] + 1; j++)
137  {
138  gEFieldGraphZ[i][j] = (TGraph*)infile->Get(Form("deltaEz/g%i_%i", i, j));
139 
140  }
141  intermediateEFieldFitFunctionZ[i] = new TF1(Form("intermediateEFieldFitFunctionZ_%i", i), Form("pol%i", intermediateEFieldFitPolN[2]));
142  }
143 
144  initialEFieldFitFunctionX = new TF1("initialEFieldFitFunctionX", Form("pol%i", initialEFieldFitPolN[0]));
145  initialEFieldFitFunctionY = new TF1("initialEFieldFitFunctionY", Form("pol%i", initialEFieldFitPolN[1]));
146  initialEFieldFitFunctionZ = new TF1("initialEFieldFitFunctionZ", Form("pol%i", initialEFieldFitPolN[2]));
147  }else{
148  std::cout << "fRepresentationType not known!!!" << std::endl;
149  }
150  infile->Close();
151  }
152 
153  if(fEnableCorrSCE == true)
154  {
155  // Grab other parameters from pset
156  }
157  return true;
158 }
string fname
Definition: demo.py:5
TGraph * gSpatialGraphY[99][99]
TGraph * gSpatialGraphX[99][99]
TGraph * gEFieldGraphZ[99][99]
TGraph * gEFieldGraphX[99][99]
TGraph * gSpatialGraphZ[99][99]
std::vector< TH3F * > SCEhistograms
TGraph * gEFieldGraphY[99][99]
BEGIN_PROLOG could also be cout
bool spacecharge::SpaceChargeSBND::EnableCalEfieldSCE ( ) const
overridevirtual

Implements spacecharge::SpaceCharge.

Definition at line 195 of file SpaceChargeSBND.cxx.

196 {
197  return fEnableCalEfieldSCE;
198 }
bool spacecharge::SpaceChargeSBND::EnableCalSpatialSCE ( ) const
overridevirtual

Implements spacecharge::SpaceCharge.

Definition at line 189 of file SpaceChargeSBND.cxx.

190 {
191  return fEnableCalSpatialSCE;
192 }
bool spacecharge::SpaceChargeSBND::EnableCorrSCE ( ) const
inlineoverridevirtual

Implements spacecharge::SpaceCharge.

Definition at line 36 of file SpaceChargeSBND.h.

bool EnableCalEfieldSCE() const override
bool EnableCalSpatialSCE() const override
bool spacecharge::SpaceChargeSBND::EnableSimEfieldSCE ( ) const
overridevirtual

Implements spacecharge::SpaceCharge.

Definition at line 177 of file SpaceChargeSBND.cxx.

178 {
179  return fEnableSimEfieldSCE;
180 }
bool spacecharge::SpaceChargeSBND::EnableSimSpatialSCE ( ) const
overridevirtual

Implements spacecharge::SpaceCharge.

Definition at line 171 of file SpaceChargeSBND.cxx.

172 {
173  return fEnableSimSpatialSCE;
174 }
geo::Vector_t spacecharge::SpaceChargeSBND::GetCalEfieldOffsets ( geo::Point_t const &  point,
int const &  TPCid = 1 
) const
inlineoverridevirtual

Implements spacecharge::SpaceCharge.

Definition at line 41 of file SpaceChargeSBND.h.

41 { return {0.,0.,0.}; }
geo::Vector_t spacecharge::SpaceChargeSBND::GetCalPosOffsets ( geo::Point_t const &  point,
int const &  TPCid = 1 
) const
overridevirtual

Implements spacecharge::SpaceCharge.

Definition at line 240 of file SpaceChargeSBND.cxx.

241 {
242  std::vector<double> theCalPosOffsets;
243  double xx=point.X(), yy=point.Y(), zz=point.Z();
244 
245  if(fRepresentationType == "Voxelized_TH3"){
246  //handle OOAV by projecting edge cases
247  if(xx<-199.999){xx=-199.999;}
248  else if(xx>199.999){xx=199.999;}
249  if(yy<-199.999){yy=-199.999;}
250  else if(yy>199.999){yy=199.999;}
251  if(zz<0.001){zz=0.001;}
252  else if(zz>499.999){zz=499.999;}
253  //correct for charge drifted across cathode
254  if ((TPCid == 0) and (xx > -2.5)) { xx = -2.5; }
255  if ((TPCid == 1) and (xx < 2.5)) { xx = 2.5; }
256  double offset_x=0., offset_y=0., offset_z=0.;
257  offset_x = SCEhistograms.at(3)->Interpolate(xx,yy,zz);
258  offset_y = SCEhistograms.at(4)->Interpolate(xx,yy,zz);
259  offset_z = SCEhistograms.at(5)->Interpolate(xx,yy,zz);
260  theCalPosOffsets = {offset_x, offset_y, offset_z};
261 
262  }else if(fRepresentationType == "Parametric"){
263  //this is not supported for parametric
264  std::cout << "Change Representation Type to Voxelized TH3 if you want to use the backward offset function" << std::endl;
265  theCalPosOffsets.resize(3, 0.0);
266  }
267 
268  return { theCalPosOffsets[0], theCalPosOffsets[1], theCalPosOffsets[2] };
269 }
return match has_match and(match.match_pdg==11 or match.match_pdg==-11)
std::vector< TH3F * > SCEhistograms
BEGIN_PROLOG could also be cout
geo::Vector_t spacecharge::SpaceChargeSBND::GetEfieldOffsets ( geo::Point_t const &  point) const
overridevirtual

Implements spacecharge::SpaceCharge.

Definition at line 383 of file SpaceChargeSBND.cxx.

384 {
385  std::vector<double> theEfieldOffsets;
386  double xx=point.X(), yy=point.Y(), zz=point.Z();
387  double offset_x=0., offset_y=0., offset_z=0.;
388 
389  if(fRepresentationType == "Voxelized_TH3"){
390  //handle OOAV by projecting edge cases
391  if(xx<-199.999){xx=-199.999;}
392  else if(xx>199.999){xx=199.999;}
393  if(yy<-199.999){yy=-199.999;}
394  else if(yy>199.999){yy=199.999;}
395  if(zz<0.001){zz=0.001;}
396  else if(zz>499.999){zz=499.999;}
397  offset_x = SCEhistograms.at(6)->Interpolate(xx, yy, zz);
398  offset_y = SCEhistograms.at(7)->Interpolate(xx, yy, zz);
399  offset_z = SCEhistograms.at(8)->Interpolate(xx, yy, zz);
400 
401  theEfieldOffsets = {offset_x, offset_y, offset_z};
402 
403  }else if(fRepresentationType == "Parametric"){
404 
405  if(IsInsideBoundaries(point.X(), point.Y(), point.Z()) == false){
406  theEfieldOffsets.resize(3, 0.0);
407  }
408  else
409  {
410  theEfieldOffsets = GetEfieldOffsetsParametric(point.X(), point.Y(), point.Z());
411 
412  // GetOneEfieldOffsetParametric returns V/m
413  // The E-field offsets are returned as -dEx/|E_nominal|, -dEy/|E_nominal|, and -dEz/|E_nominal| where |E_nominal| is DriftField
414  theEfieldOffsets[0] = -1.0 * theEfieldOffsets[0] / (100.0 * DriftField);
415  theEfieldOffsets[1] = -1.0 * theEfieldOffsets[1] / (100.0 * DriftField);
416  theEfieldOffsets[2] = -1.0 * theEfieldOffsets[2] / (100.0 * DriftField);
417  }
418  }
419 
420  return { theEfieldOffsets[0], theEfieldOffsets[1], theEfieldOffsets[2] };
421 }
bool IsInsideBoundaries(double xVal, double yVal, double zVal) const
std::vector< double > GetEfieldOffsetsParametric(double xVal, double yVal, double zVal) const
std::vector< TH3F * > SCEhistograms
std::vector< double > spacecharge::SpaceChargeSBND::GetEfieldOffsetsParametric ( double  xVal,
double  yVal,
double  zVal 
) const
protected

Definition at line 424 of file SpaceChargeSBND.cxx.

425 {
426  std::vector<double> theEfieldOffsetsParametric;
427 
428  double xValNew = TransformX(xVal);
429  double yValNew = TransformY(yVal);
430  double zValNew = TransformZ(zVal);
431 
432  theEfieldOffsetsParametric.push_back(GetOneEfieldOffsetParametric(xValNew, yValNew, zValNew, "X"));
433  theEfieldOffsetsParametric.push_back(GetOneEfieldOffsetParametric(xValNew, yValNew, zValNew, "Y"));
434  theEfieldOffsetsParametric.push_back(GetOneEfieldOffsetParametric(xValNew, yValNew, zValNew, "Z"));
435 
436  return theEfieldOffsetsParametric;
437 }
double TransformZ(double zVal) const
double TransformY(double yVal) const
double GetOneEfieldOffsetParametric(double xVal, double yVal, double zVal, std::string axis) const
double TransformX(double xVal) const
double spacecharge::SpaceChargeSBND::GetOneEfieldOffsetParametric ( double  xVal,
double  yVal,
double  zVal,
std::string  axis 
) const
protected

Definition at line 440 of file SpaceChargeSBND.cxx.

441 {
442  double parA[99][99];
443  double parB[99];
444 
445  for(int i = 0; i < 99; i++)
446  {
447  for(int j = 0; j < 99; j++)
448  {
449  parA[i][j] = 0.0;
450  }
451  parB[i] = 0.0;
452  }
453  if(axis == "X")
454  {
455  for(int i = 0; i < initialEFieldFitPolN[0] + 1; i++)
456  {
457  for(int j = 0; j < intermediateEFieldFitPolN[0] + 1; j++)
458  {
459  parA[i][j] = gEFieldGraphX[i][j]->Eval(zValNew);
460  }
461  intermediateEFieldFitFunctionX[i]->SetParameters(parA[i]);
462  }
463  }
464  else if(axis == "Y")
465  {
466  for(int i = 0; i < initialEFieldFitPolN[1] + 1; i++)
467  {
468  for(int j = 0; j < intermediateEFieldFitPolN[1] + 1; j++)
469  {
470  parA[i][j] = gEFieldGraphY[i][j]->Eval(zValNew);
471  }
472  intermediateEFieldFitFunctionY[i]->SetParameters(parA[i]);
473  }
474  }
475  else if(axis == "Z")
476  {
477  for(int i = 0; i < initialEFieldFitPolN[2] + 1; i++)
478  {
479  for(int j = 0; j < intermediateEFieldFitPolN[2] + 1; j++)
480  {
481  parA[i][j] = gEFieldGraphZ[i][j]->Eval(zValNew);
482  }
483  intermediateEFieldFitFunctionZ[i]->SetParameters(parA[i]);
484  }
485  }
486 
487  double aValNew;
488  double bValNew;
489 
490  if(axis == "Y")
491  {
492  aValNew = xValNew;
493  bValNew = yValNew;
494  }
495  else
496  {
497  aValNew = yValNew;
498  bValNew = xValNew;
499  }
500 
501  double offsetValNew = 0.0;
502  if(axis == "X")
503  {
504  for(int i = 0; i < initialEFieldFitPolN[0] + 1; i++)
505  {
506  parB[i] = intermediateEFieldFitFunctionX[i]->Eval(aValNew);
507  }
508  initialEFieldFitFunctionX->SetParameters(parB);
509  offsetValNew = initialEFieldFitFunctionX->Eval(bValNew);
510  }
511  else if(axis == "Y")
512  {
513  for(int i = 0; i < initialEFieldFitPolN[1] + 1; i++)
514  {
515  parB[i] = intermediateEFieldFitFunctionY[i]->Eval(aValNew);
516  }
517  initialEFieldFitFunctionY->SetParameters(parB);
518  offsetValNew = initialEFieldFitFunctionY->Eval(bValNew);
519  }
520  else if(axis == "Z")
521  {
522  for(int i = 0; i < initialEFieldFitPolN[2] + 1; i++)
523  {
524  parB[i] = intermediateEFieldFitFunctionZ[i]->Eval(aValNew);
525  }
526  initialEFieldFitFunctionZ->SetParameters(parB);
527  offsetValNew = initialEFieldFitFunctionZ->Eval(bValNew);
528  }
529 
530  return offsetValNew;
531 }
TGraph * gEFieldGraphZ[99][99]
TGraph * gEFieldGraphX[99][99]
TGraph * gEFieldGraphY[99][99]
double spacecharge::SpaceChargeSBND::GetOnePosOffsetParametric ( double  xVal,
double  yVal,
double  zVal,
std::string  axis 
) const
protected

Definition at line 290 of file SpaceChargeSBND.cxx.

291 {
292  double parA[99][99];
293  double parB[99];
294 
295  for(int i = 0; i < 99; i++)
296  {
297  for(int j = 0; j < 99; j++)
298  {
299  parA[i][j] = 0.0;
300  }
301  parB[i] = 0.0;
302  }
303  if(axis == "X")
304  {
305  for(int i = 0; i < initialSpatialFitPolN[0] + 1; i++)
306  {
307  for(int j = 0; j < intermediateSpatialFitPolN[0] + 1; j++)
308  {
309  parA[i][j] = gSpatialGraphX[i][j]->Eval(zValNew);
310  }
311  intermediateSpatialFitFunctionX[i]->SetParameters(parA[i]);
312  }
313  }
314  else if(axis == "Y")
315  {
316  for(int i = 0; i < initialSpatialFitPolN[1] + 1; i++)
317  {
318  for(int j = 0; j < intermediateSpatialFitPolN[1] + 1; j++)
319  {
320  parA[i][j] = gSpatialGraphY[i][j]->Eval(zValNew);
321  }
322  intermediateSpatialFitFunctionY[i]->SetParameters(parA[i]);
323  }
324  }
325  else if(axis == "Z")
326  {
327  for(int i = 0; i < initialSpatialFitPolN[2] + 1; i++)
328  {
329  for(int j = 0; j < intermediateSpatialFitPolN[2] + 1; j++)
330  {
331  parA[i][j] = gSpatialGraphZ[i][j]->Eval(zValNew);
332  }
333  intermediateSpatialFitFunctionZ[i]->SetParameters(parA[i]);
334  }
335  }
336 
337  double aValNew;
338  double bValNew;
339 
340  if(axis == "Y")
341  {
342  aValNew = xValNew;
343  bValNew = yValNew;
344  }
345  else
346  {
347  aValNew = yValNew;
348  bValNew = xValNew;
349  }
350  double offsetValNew = 0.0;
351  if(axis == "X")
352  {
353  for(int i = 0; i < initialSpatialFitPolN[0] + 1; i++)
354  {
355  parB[i] = intermediateSpatialFitFunctionX[i]->Eval(aValNew);
356  }
357  initialSpatialFitFunctionX->SetParameters(parB);
358  offsetValNew = initialSpatialFitFunctionX->Eval(bValNew);
359  }
360  else if(axis == "Y")
361  {
362  for(int i = 0; i < initialSpatialFitPolN[1] + 1; i++)
363  {
364  parB[i] = intermediateSpatialFitFunctionY[i]->Eval(aValNew);
365  }
366  initialSpatialFitFunctionY->SetParameters(parB);
367  offsetValNew = initialSpatialFitFunctionY->Eval(bValNew);
368  }
369  else if(axis == "Z")
370  {
371  for(int i = 0; i < initialSpatialFitPolN[2] + 1; i++)
372  {
373  parB[i] = intermediateSpatialFitFunctionZ[i]->Eval(aValNew);
374  }
375  initialSpatialFitFunctionZ->SetParameters(parB);
376  offsetValNew = initialSpatialFitFunctionZ->Eval(bValNew);
377  }
378 
379  return offsetValNew;
380 }
TGraph * gSpatialGraphY[99][99]
TGraph * gSpatialGraphX[99][99]
TGraph * gSpatialGraphZ[99][99]
geo::Vector_t spacecharge::SpaceChargeSBND::GetPosOffsets ( geo::Point_t const &  point) const
overridevirtual

Implements spacecharge::SpaceCharge.

Definition at line 201 of file SpaceChargeSBND.cxx.

202 {
203  std::vector<double> thePosOffsets;
204  double xx=point.X(), yy=point.Y(), zz=point.Z();
205 
206  if(fRepresentationType == "Voxelized_TH3"){
207  //handle OOAV by projecting edge cases
208  if(xx<-199.999){xx=-199.999;}
209  else if(xx>199.999){xx=199.999;}
210  if(yy<-199.999){yy=-199.999;}
211  else if(yy>199.999){yy=199.999;}
212  if(zz<0.001){zz=0.001;}
213  else if(zz>499.999){zz=499.999;}
214  //larsim requires negative sign in TPC 0
215  int corr = 1;
216  if (xx < 0) { corr = -1; }
217  double offset_x=0., offset_y=0., offset_z=0.;
218  offset_x = corr*SCEhistograms.at(0)->Interpolate(xx,yy,zz);
219  offset_y = SCEhistograms.at(1)->Interpolate(xx,yy,zz);
220  offset_z = SCEhistograms.at(2)->Interpolate(xx,yy,zz);
221  thePosOffsets = {offset_x, offset_y, offset_z};
222 
223  }else if(fRepresentationType == "Parametric"){
224  if(IsInsideBoundaries(point.X(), point.Y(), point.Z()) == false){
225  thePosOffsets.resize(3, 0.0);
226  }else{
227  // GetPosOffsetsParametric returns m; the PosOffsets should be in cm
228  thePosOffsets = GetPosOffsetsParametric(xx, yy, zz);
229  for(int i=0; i<3; i++){thePosOffsets[i]=100.*thePosOffsets[i];}
230  }
231  }
232  else{
233  thePosOffsets.resize(3, 0.0);
234  }
235 
236  return { thePosOffsets[0], thePosOffsets[1], thePosOffsets[2] };
237 }
std::vector< double > GetPosOffsetsParametric(double xVal, double yVal, double zVal) const
bool IsInsideBoundaries(double xVal, double yVal, double zVal) const
std::vector< TH3F * > SCEhistograms
std::vector< double > spacecharge::SpaceChargeSBND::GetPosOffsetsParametric ( double  xVal,
double  yVal,
double  zVal 
) const
protected

Definition at line 274 of file SpaceChargeSBND.cxx.

275 {
276  std::vector<double> thePosOffsetsParametric;
277 
278  double xValNew = TransformX(xVal);
279  double yValNew = TransformY(yVal);
280  double zValNew = TransformZ(zVal);
281 
282  thePosOffsetsParametric.push_back(GetOnePosOffsetParametric(xValNew, yValNew, zValNew, "X"));
283  thePosOffsetsParametric.push_back(GetOnePosOffsetParametric(xValNew, yValNew, zValNew, "Y"));
284  thePosOffsetsParametric.push_back(GetOnePosOffsetParametric(xValNew, yValNew, zValNew, "Z"));
285 
286  return thePosOffsetsParametric;
287 }
double TransformZ(double zVal) const
double TransformY(double yVal) const
double TransformX(double xVal) const
double GetOnePosOffsetParametric(double xVal, double yVal, double zVal, std::string axis) const
bool spacecharge::SpaceChargeSBND::IsInsideBoundaries ( double  xVal,
double  yVal,
double  zVal 
) const
protected

Definition at line 565 of file SpaceChargeSBND.cxx.

566 {
567  bool isInside = true;
568 
569  if((xVal < -196.5) || (xVal > 196.5) || (yVal < -200.0) || (yVal > 200.0) || (zVal < 0) || (zVal > 500.0))
570  {
571  isInside = false;
572  }
573 
574  return isInside;
575 }
double spacecharge::SpaceChargeSBND::TransformX ( double  xVal) const
protected

Definition at line 535 of file SpaceChargeSBND.cxx.

536 {
537  xVal = xVal / 100.0;
538 
539  // We use the same map twice; 0 is cathod
540  // Map [0, 196.5] to [0, 2.0]
541  double xValNew = (2.0 / 1.965) * fabs(xVal);
542 
543  return xValNew;
544 }
double spacecharge::SpaceChargeSBND::TransformY ( double  yVal) const
protected

Definition at line 548 of file SpaceChargeSBND.cxx.

549 {
550  yVal = yVal / 100.0;
551 
552  return (yVal + 2.0);
553 }
double spacecharge::SpaceChargeSBND::TransformZ ( double  zVal) const
protected

Definition at line 557 of file SpaceChargeSBND.cxx.

558 {
559  return (zVal / 100.0);
560 }
bool spacecharge::SpaceChargeSBND::Update ( uint64_t  ts = 0)

Definition at line 160 of file SpaceChargeSBND.cxx.

161 {
162  if (ts == 0)
163  {
164  return false;
165  }
166 
167  return true;
168 }

Member Data Documentation

double spacecharge::SpaceChargeSBND::DriftField = 500.0
protected

Definition at line 50 of file SpaceChargeSBND.h.

bool spacecharge::SpaceChargeSBND::fEnableCalEfieldSCE
protected

Definition at line 55 of file SpaceChargeSBND.h.

bool spacecharge::SpaceChargeSBND::fEnableCalSpatialSCE
protected

Definition at line 54 of file SpaceChargeSBND.h.

bool spacecharge::SpaceChargeSBND::fEnableCorrSCE
protected

Definition at line 56 of file SpaceChargeSBND.h.

bool spacecharge::SpaceChargeSBND::fEnableSimEfieldSCE
protected

Definition at line 53 of file SpaceChargeSBND.h.

bool spacecharge::SpaceChargeSBND::fEnableSimSpatialSCE
protected

Definition at line 52 of file SpaceChargeSBND.h.

std::string spacecharge::SpaceChargeSBND::fInputFilename
protected

Definition at line 59 of file SpaceChargeSBND.h.

std::string spacecharge::SpaceChargeSBND::fRepresentationType
protected

Definition at line 58 of file SpaceChargeSBND.h.

TGraph* spacecharge::SpaceChargeSBND::gEFieldGraphX[99][99]
protected

Definition at line 83 of file SpaceChargeSBND.h.

TGraph* spacecharge::SpaceChargeSBND::gEFieldGraphY[99][99]
protected

Definition at line 86 of file SpaceChargeSBND.h.

TGraph* spacecharge::SpaceChargeSBND::gEFieldGraphZ[99][99]
protected

Definition at line 89 of file SpaceChargeSBND.h.

TGraph* spacecharge::SpaceChargeSBND::gSpatialGraphX[99][99]
protected

Definition at line 73 of file SpaceChargeSBND.h.

TGraph* spacecharge::SpaceChargeSBND::gSpatialGraphY[99][99]
protected

Definition at line 76 of file SpaceChargeSBND.h.

TGraph* spacecharge::SpaceChargeSBND::gSpatialGraphZ[99][99]
protected

Definition at line 79 of file SpaceChargeSBND.h.

TF1* spacecharge::SpaceChargeSBND::initialEFieldFitFunctionX
protected

Definition at line 85 of file SpaceChargeSBND.h.

TF1* spacecharge::SpaceChargeSBND::initialEFieldFitFunctionY
protected

Definition at line 88 of file SpaceChargeSBND.h.

TF1* spacecharge::SpaceChargeSBND::initialEFieldFitFunctionZ
protected

Definition at line 91 of file SpaceChargeSBND.h.

int spacecharge::SpaceChargeSBND::initialEFieldFitPolN[3] = {3, 3, 3}
protected

Definition at line 48 of file SpaceChargeSBND.h.

TF1* spacecharge::SpaceChargeSBND::initialSpatialFitFunctionX
protected

Definition at line 75 of file SpaceChargeSBND.h.

TF1* spacecharge::SpaceChargeSBND::initialSpatialFitFunctionY
protected

Definition at line 78 of file SpaceChargeSBND.h.

TF1* spacecharge::SpaceChargeSBND::initialSpatialFitFunctionZ
protected

Definition at line 81 of file SpaceChargeSBND.h.

int spacecharge::SpaceChargeSBND::initialSpatialFitPolN[3] = {3, 4, 3}
protected

Definition at line 46 of file SpaceChargeSBND.h.

TF1* spacecharge::SpaceChargeSBND::intermediateEFieldFitFunctionX[99]
protected

Definition at line 84 of file SpaceChargeSBND.h.

TF1* spacecharge::SpaceChargeSBND::intermediateEFieldFitFunctionY[99]
protected

Definition at line 87 of file SpaceChargeSBND.h.

TF1* spacecharge::SpaceChargeSBND::intermediateEFieldFitFunctionZ[99]
protected

Definition at line 90 of file SpaceChargeSBND.h.

int spacecharge::SpaceChargeSBND::intermediateEFieldFitPolN[3] = {6, 4, 4}
protected

Definition at line 49 of file SpaceChargeSBND.h.

TF1* spacecharge::SpaceChargeSBND::intermediateSpatialFitFunctionX[99]
protected

Definition at line 74 of file SpaceChargeSBND.h.

TF1* spacecharge::SpaceChargeSBND::intermediateSpatialFitFunctionY[99]
protected

Definition at line 77 of file SpaceChargeSBND.h.

TF1* spacecharge::SpaceChargeSBND::intermediateSpatialFitFunctionZ[99]
protected

Definition at line 80 of file SpaceChargeSBND.h.

int spacecharge::SpaceChargeSBND::intermediateSpatialFitPolN[3] = {4, 4, 4}
protected

Definition at line 47 of file SpaceChargeSBND.h.

std::vector<TH3F*> spacecharge::SpaceChargeSBND::SCEhistograms = std::vector<TH3F*>(9)
protected

Definition at line 71 of file SpaceChargeSBND.h.


The documentation for this class was generated from the following files: