REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestDetectorHitsSpecularProcess.cxx
1 /*************************************************************************
2  * This file is part of the REST software framework. *
3  * *
4  * Copyright (C) 2016 GIFNA/TREX (University of Zaragoza) *
5  * For more information see http://gifna.unizar.es/trex *
6  * *
7  * REST is free software: you can redistribute it and/or modify *
8  * it under the terms of the GNU General Public License as published by *
9  * the Free Software Foundation, either version 3 of the License, or *
10  * (at your option) any later version. *
11  * *
12  * REST is distributed in the hope that it will be useful, *
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of *
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
15  * GNU General Public License for more details. *
16  * *
17  * You should have a copy of the GNU General Public License along with *
18  * REST in $REST_PATH/LICENSE. *
19  * If not, see http://www.gnu.org/licenses/. *
20  * For the list of contributors see $REST_PATH/CREDITS. *
21  *************************************************************************/
22 
67 #include "TRestDetectorHitsSpecularProcess.h"
68 
69 using namespace std;
70 
72 
73 TRestDetectorHitsSpecularProcess::TRestDetectorHitsSpecularProcess() { Initialize(); }
74 
87 TRestDetectorHitsSpecularProcess::TRestDetectorHitsSpecularProcess(const char* configFilename) {
88  Initialize();
89  LoadConfigFromFile(configFilename);
90 }
91 
92 TRestDetectorHitsSpecularProcess::~TRestDetectorHitsSpecularProcess() {}
93 
95  SetSectionName(this->ClassName());
96  SetLibraryVersion(LIBRARY_VERSION);
97 
98  fInputEvent = nullptr;
99  fOutputEvent = new TRestDetectorHitsEvent();
100 }
101 
103  fInputEvent = (TRestDetectorHitsEvent*)inputEvent;
104  fOutputEvent->SetEventInfo(fInputEvent);
105 
106  for (unsigned int hit = 0; hit < fInputEvent->GetNumberOfHits(); hit++) {
107  TVector3 position(fInputEvent->GetX(hit), fInputEvent->GetY(hit), fInputEvent->GetZ(hit));
108  const auto type = fInputEvent->GetType(hit);
109  const auto energy = fInputEvent->GetEnergy(hit);
110  const auto time = fInputEvent->GetTime(hit);
111 
112  if (type != XYZ) {
113  fOutputEvent->AddHit(position, energy, time, type);
114  continue;
115  }
116 
117  const TVector3 relativePosition = position - fPosition;
118  const TVector3 reflectionVector = relativePosition - 2 * relativePosition.Dot(fNormal) * fNormal;
119 
120  position = fPosition + reflectionVector;
121 
122  fOutputEvent->AddHit(position, energy, time, type);
123  }
124  return fOutputEvent;
125 }
126 
129 
130  fNormal = Get3DVectorParameterWithUnits("normal", {0, 0, 1});
131  fNormal = fNormal.Unit();
132  fPosition = Get3DVectorParameterWithUnits("position", {0, 0, 1});
133 }
134 
136  BeginPrintProcess();
137 
138  RESTMetadata << " - Plane normal : ( " << fNormal.X() << ", " << fNormal.Y() << ", " << fNormal.Z()
139  << ") mm" << RESTendl;
140  RESTMetadata << " - Plane position : ( " << fPosition.X() << ", " << fPosition.Y() << ", "
141  << fPosition.Z() << ")" << RESTendl;
142 
143  EndPrintProcess();
144 }
A process to create specular hit images using a plane definition.
TRestEvent * ProcessEvent(TRestEvent *inputEvent) override
Process one event.
void InitFromConfigFile() override
To make settings from rml file. This method must be implemented in the derived class.
void PrintMetadata() override
Implemented it in the derived metadata class to print out specific metadata information.
void Initialize() override
Making default settings.
virtual void InitFromConfigFile() override
To make settings from rml file. This method must be implemented in the derived class.
A base class for any REST event.
Definition: TRestEvent.h:38
void SetEventInfo(TRestEvent *eve)
Definition: TRestEvent.cxx:137