ATLAS Offline Software
Loading...
Searching...
No Matches
ITileCalculator.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2019 CERN for the benefit of the ATLAS collaboration
3*/
4
5//************************************************************
6// Filename : TileGeoG4SDCalc.hh
7// Author : Sergey Karpov <Sergey.Karpov@cern.ch>
8// Created : July, 2013
9//
10// DESCRIPTION:
11// Sensitive Detector initialisation for TileCal G4 simulations
12// of both ordinary Hits and CalibHits
13//
14// HISTORY:
15// Nov 2013 - Work with U-shape was added (Sasha Solodkov)
16//
17//************************************************************
18
19#ifndef TILEG4INTERFACES_ITILECALCULATOR_H
20#define TILEG4INTERFACES_ITILECALCULATOR_H
21
22// Geant4 headers
23#include "G4Types.hh"
24
25#include "GaudiKernel/IService.h"
26#include "Identifier/Identifier.h"
27
28class G4Step;
29class TileGeoG4Section;
30class TileGeoG4Cell;
31class TileGeoG4LookupBuilder;
32class TileSDOptions;
33
37 G4double e_up;
38 G4double e_down;
39 double time_up;
40 double time_down;
41 //int period, tilerow; // prepared for future use
42};
43
46 TileGeoG4Section *section = nullptr;
47 TileGeoG4Cell *cell = nullptr;
48 int nModule = 0;
49 int nDetector = 0;
50 int nTower = 0;
51 int nSample = 0;
52 int nSide = 0;
53 int nrOfPMT = 0;
54 int tileSize = 0;
55 int tilePeriod = 0;
56 double totalTimeUp = 0.0;
57 double totalTimeDown = 0.0;
58 bool isNegative = false;
61 G4double edep_up = 0.0;
62 G4double edep_down = 0.0;
63 double scin_Time_up = 0.0;
64 double scin_Time_down = 0.0;
65};
66
67class ITileCalculator : virtual public IService {
68 public:
72 virtual ~ITileCalculator() {}
75
77 virtual G4bool FindTileScinSection(const G4Step*, TileHitData& hitData) const = 0;
79 virtual G4bool MakePmtEdepTime(const G4Step*, TileHitData& hitData, double& deltaTime) const = 0;
81 virtual G4bool ManageScintHit(TileHitData& hitData, double deltaTime) const = 0;
83 virtual TileMicroHit GetTileMicroHit(const G4Step*, TileHitData& hitData) const = 0;
85 virtual std::unique_ptr<TileGeoG4LookupBuilder> GetLookupBuilder() const = 0;
87 virtual const TileSDOptions* GetOptions() const = 0;
89 virtual void pmtEdepFromFCS_StepInfo(TileHitData& hitData, double ene, double yLocal, double halfYLocal, double zLocal, int Ushape) const = 0;
90
91};
92//class ITileCalculator
93
94#endif // not TILEG4INTERFACES_ITILECALCULATOR_H
virtual const TileSDOptions * GetOptions() const =0
pointer to class with all options
virtual void pmtEdepFromFCS_StepInfo(TileHitData &hitData, double ene, double yLocal, double halfYLocal, double zLocal, int Ushape) const =0
Method used by TileFastCaloSim/TileFCSmStepToTileHitVec.
virtual G4bool MakePmtEdepTime(const G4Step *, TileHitData &hitData, double &deltaTime) const =0
Calculation of pmtID, edep and scin_Time with aStep (Sergey)
virtual G4bool FindTileScinSection(const G4Step *, TileHitData &hitData) const =0
Search for the tilecal sub-section, its module and some identifiers.
virtual TileMicroHit GetTileMicroHit(const G4Step *, TileHitData &hitData) const =0
Used by FastCaloSimParamAction.
virtual std::unique_ptr< TileGeoG4LookupBuilder > GetLookupBuilder() const =0
virtual ~ITileCalculator()
virtual G4bool ManageScintHit(TileHitData &hitData, double deltaTime) const =0
Calculation of pmtID, edep and scin_Time with aStep (Sergey)
DeclareInterfaceID(ITileCalculator, 1, 0)
Variables to identify Hit objects.
double totalTimeDown
TileGeoG4Section * section
Identifier pmtID_down
G4double edep_up
TileGeoG4Cell * cell
Identifier pmtID_up
double scin_Time_down
double scin_Time_up
G4double edep_down
Identifier pmt_down
Identifier pmt_up