REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestDetectorHitsRotationProcess.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 
66 #include "TRestDetectorHitsRotationProcess.h"
67 
68 using namespace std;
69 
71 
72 TRestDetectorHitsRotationProcess::TRestDetectorHitsRotationProcess() { Initialize(); }
73 
86 TRestDetectorHitsRotationProcess::TRestDetectorHitsRotationProcess(const char* configFilename) {
87  Initialize();
88  LoadConfigFromFile(configFilename);
89 }
90 
91 TRestDetectorHitsRotationProcess::~TRestDetectorHitsRotationProcess() {}
92 
94  SetSectionName(this->ClassName());
95  SetLibraryVersion(LIBRARY_VERSION);
96 
97  fInputEvent = nullptr;
98  fOutputEvent = new TRestDetectorHitsEvent();
99 }
100 
102  fInputEvent = (TRestDetectorHitsEvent*)inputEvent;
103  fOutputEvent->SetEventInfo(fInputEvent);
104 
105  for (unsigned int hit = 0; hit < fInputEvent->GetNumberOfHits(); hit++) {
106  TVector3 position = {fInputEvent->GetX(hit), fInputEvent->GetY(hit), fInputEvent->GetZ(hit)};
107  const auto type = fInputEvent->GetType(hit);
108  const auto energy = fInputEvent->GetEnergy(hit);
109  const auto time = fInputEvent->GetTime(hit);
110 
111  if (type != XYZ) {
112  fOutputEvent->AddHit(position, energy, time, type);
113  continue;
114  }
115 
116  // veto hits won't be rotated because they are not XYZ
117 
118  position -= fCenter;
119  position.Rotate(fAngle, fAxis);
120  position += fCenter;
121 
122  fOutputEvent->AddHit(position, energy, time, type);
123  }
124 
125  return fOutputEvent;
126 }
127 
130 
131  fAxis = Get3DVectorParameterWithUnits("axis", fAxis);
132  fAxis.Unit();
133  fCenter = Get3DVectorParameterWithUnits("center", fCenter);
134 }
135 
137  BeginPrintProcess();
138 
139  RESTMetadata << " - Rotation center : ( " << fCenter.X() << ", " << fCenter.Y() << ", " << fCenter.Z()
140  << ") mm" << RESTendl;
141  RESTMetadata << " - Rotation axis : ( " << fAxis.X() << ", " << fAxis.Y() << ", " << fAxis.Z() << ")"
142  << RESTendl;
143  RESTMetadata << " - Rotation angle : " << fAngle * units("degrees") << " degrees" << RESTendl;
144 
145  EndPrintProcess();
146 }
A process to rotate hits from a given center and axis.
void InitFromConfigFile() override
To make settings from rml file. This method must be implemented in the derived class.
void Initialize() override
Making default settings.
void PrintMetadata() override
Implemented it in the derived metadata class to print out specific metadata information.
TRestEvent * ProcessEvent(TRestEvent *inputEvent) override
Process one event.
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