Geant4 Cross Reference |
1 // 2 // ******************************************************************** 3 // * License and Disclaimer * 4 // * * 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. * 10 // * * 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitation of liability. * 17 // * * 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************************************** 25 // 26 // 27 // 28 /// \file ScanDamage.hh 29 /// \brief Definition of the ScanDamage class 30 31 #ifndef ScanDamage_h 32 #define ScanDamage_h 33 34 #include <map> 35 #include <vector> 36 #include <string> 37 #include <tuple> 38 #include <set> 39 #include "Damage.hh" 40 #include <filesystem> 41 namespace fs = std::filesystem; 42 class TFile; 43 using ullint = unsigned long long int; 44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 45 46 struct VoxelData 47 { 48 VoxelData(int chromo, int domain, ullint firstBpCN) 49 { 50 fChromosome = chromo; 51 fDomain = domain; 52 fFirstBpCopyNum = firstBpCN; 53 } 54 55 ~VoxelData() {} 56 57 int fChromosome{0}; 58 int fDomain{0}; 59 ullint fFirstBpCopyNum{0}; 60 }; 61 62 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 63 64 typedef std::vector<std::vector<ullint> > Table; 65 66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 67 68 class Trier { 69 public: 70 bool operator()(const std::vector<ullint>& a, const std::vector<ullint>& b) 71 { 72 bool bb = false; 73 if(a[0] < b[0]) bb = true; 74 return bb; 75 } 76 }; 77 78 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 79 80 class ScanDamage 81 { 82 public: 83 ScanDamage(); 84 ~ScanDamage() = default; 85 std::map<unsigned int,std::map<unsigned int,std::vector<Damage> > > ExtractDamage(); 86 void SetThresholdEnergy(double e) {fThresholdEnergy = e;} 87 void SetProbabilityForIndirectSBSelection(double p) {fProbabilityForIndirectSB = p;} 88 double GetThresholdEnergy() {return fThresholdEnergy;} 89 double GetProbabilityForIndirectSBSelection() {return fProbabilityForIndirectSB;} 90 std::map<int, Table> GetMergedSBData() {return fMergedTables;} 91 double GetEdepSumInNucleus() {return fEdepSumInNucleus;} //eV 92 double GetTotalNbBpPlacedInGeo() {return fTotalNbBpPlacedInGeo;} 93 double GetTotalNbHistonePlacedInGeo() {return fTotalNbHistonePlacedInGeo;} 94 double GetNucleusVolume() {return fNucleusVolume;} 95 double GetNucleusMassDensity() {return fNucleusMassDensity;} 96 double GetNucleusMass() {return fNucleusMass;} 97 std::map<int,ullint> GetChromosomeBpSizesMap() {return fChromosomeBpMap;} 98 void SkipScanningIndirectDamage() {fSkipScanningIndirectDamage = true;} 99 bool SkippedScanningIndirectDamage() {return fSkipScanningIndirectDamage;} 100 private: 101 void ScanDamageFromPhys(); 102 void ScanDamageFromChem(); 103 void RetrieveVoxelBp(); 104 void FillVoxelData(); 105 void AnaPhysRootFile(const std::string fileName); 106 void AnaChemRootFile(fs::directory_entry entry); 107 void AnaPhysRootTree1(TFile*); 108 void AnaPhysRootTree2(TFile*); 109 void SortPhysTableWithSelection(); 110 void SortChemTableWithSelection(); 111 void ReadCellandVoxelDefFilePaths(); 112 void MergeDamageFromPhysChem(); 113 std::tuple<unsigned int, unsigned int> GetEventNberAndVoxelNberFromChemRoot(const std::string fileNam); 114 double fThresholdEnergy{17.5};//eV 115 std::string fCellDefFilePath{""}; 116 std::set<std::string> fVoxelDefFilesList; 117 std::map<std::string, int> fBpPerVoxel; 118 std::vector<VoxelData> fVoxels; 119 std::map<int, Table> fphysTables, fphysSlectedTables, fchemTables, fchemSlectedTables,fMergedTables; 120 std::map<unsigned int,std::map<unsigned int,std::vector<Damage> > > fDamage; 121 double fEdepSumInNucleus{0}; //eV 122 int corruptedFiles = 0; 123 double fTotalNbBpPlacedInGeo{0}; 124 double fTotalNbHistonePlacedInGeo{0}; 125 double fNucleusVolume{0}; 126 double fNucleusMassDensity{0}; 127 double fNucleusMass{0}; 128 double fProbabilityForIndirectSB{0.4}; 129 std::map<int,ullint> fChromosomeBpMap; //Store number of Bp in each Chomosomes; 130 bool fSkipScanningIndirectDamage{false}; 131 }; 132 133 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 134 135 #endif 136 137 138