Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-10-14 08:06:13

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 //
0027 /// \file ScanDamage.cc
0028 /// \brief Implementation of the ScanDamage class
0029 
0030 #include "ScanDamage.hh"
0031 #include "ParametersParser.hh"
0032 
0033 #include "TSystemDirectory.h"
0034 #include "TFile.h"
0035 #include "TTree.h"
0036 #include "TRandom.h"
0037 
0038 
0039 #include <iostream>
0040 #include <fstream>
0041 #include <sstream>
0042 #include <iostream>
0043 #include <filesystem>
0044 namespace fs = std::filesystem;
0045 TRandom gRandomGen; 
0046 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0047 
0048 ScanDamage::ScanDamage(): fSkipScanningIndirectDamage(false)
0049 {
0050     fThresholdEnergy = 17.5; //eV
0051     fEdepSumInNucleus = 0;
0052     fProbabilityForIndirectSB = 0.40; // 40%
0053 }
0054 
0055 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0056 
0057 std::map<unsigned int,std::map<unsigned int,std::vector<Damage> > > ScanDamage::ExtractDamage(){
0058     ReadCellandVoxelDefFilePaths();
0059     fMergedTables.clear();
0060     fDamage.clear();
0061     RetrieveVoxelBp();
0062     FillVoxelData();
0063     ScanDamageFromPhys();
0064     if (!fSkipScanningIndirectDamage) ScanDamageFromChem();
0065     MergeDamageFromPhysChem();
0066     for (const auto& [pChrom,table] : fMergedTables) {
0067         unsigned int evt;
0068         unsigned int strand;
0069         ullint cpyNb;
0070         unsigned int isBase;
0071         double time;
0072         double edep;
0073         double origin;
0074         Damage::DamageType pType;
0075         Damage::DamageCause pOrigin;
0076         std::map<unsigned int,std::vector<Damage> > perChromoDamage;
0077         for (const auto& v : table) {
0078             evt = v.at(0);
0079             strand = v.at(1);
0080             cpyNb = v.at(2);
0081             isBase = v.at(3);
0082             time = v.at(4);
0083             edep = v.at(5);
0084             origin = v.at(6);
0085             if(isBase==0)
0086                 pType=Damage::DamageType::fBackbone;
0087             else
0088                 pType=Damage::DamageType::fBase;
0089             if(origin==0)
0090                 pOrigin=Damage::DamageCause::fDirect;
0091             else
0092                 pOrigin=Damage::DamageCause::fIndirect;
0093             if (perChromoDamage.find(evt) == perChromoDamage.end()) {
0094                 std::vector<Damage> dm{Damage(pType,pChrom,evt,strand,cpyNb,Position(0,0,0),
0095                                             pOrigin,Damage::DamageChromatin::fUnspecified)};
0096                 perChromoDamage.insert({evt,dm});
0097             } else perChromoDamage[evt].push_back(Damage(pType,pChrom,evt,strand,cpyNb,Position(0,0,0),
0098                                             pOrigin,Damage::DamageChromatin::fUnspecified));
0099         }
0100         fDamage.insert({pChrom,perChromoDamage});
0101     }
0102     return fDamage;
0103 }
0104 
0105 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0106 
0107 void ScanDamage::RetrieveVoxelBp()
0108 {
0109     fBpPerVoxel.clear();
0110     
0111     for (const auto & entry : fVoxelDefFilesList) {
0112         std::ifstream file(entry);
0113         if(!file.good() )
0114         {
0115             std::cerr<<"**** Fatal Error *****"<<std::endl;
0116             std::cerr<<"ScanDamage::RetrieveVoxelBp: No file named "<<entry<<std::endl;
0117             std::cerr<<"*************** *****"<<std::endl;
0118             exit(EXIT_FAILURE);
0119         }
0120         std::string voxelName = "noName";
0121         std::string line;
0122         bool foundName = false;
0123         bool foundNumOfBp = false;
0124         while(std::getline(file, line)
0125             && !foundNumOfBp)
0126         {
0127             std::istringstream iss(line);
0128             std::string flag;
0129             iss >> flag;
0130             std::string charac;
0131             iss >> charac;
0132             // Look for the name of the voxel
0133             if(flag=="_Name")
0134             {
0135                 voxelName = charac;
0136                 foundName = true;
0137             }
0138             // Look for the flag "_Number"
0139             // And the characteristic "voxelBasePair"
0140             if(flag=="_Number" && charac=="voxelBasePair")
0141             {
0142                 int numOfBp;
0143                 iss >> numOfBp;
0144                 if(!foundName)
0145                 {
0146                     std::cerr<<"*** Fatal Error ***"<<std::endl;
0147                     std::cerr<<"ScanDamage::RetrieveVoxelBp: The number of bp was found before the name "
0148                     <<"of the voxel... This is an unexpected case."<<std::endl;
0149                     std::cerr<<"******"<<std::endl;
0150                     exit(EXIT_FAILURE);
0151                 }
0152                 else
0153                 {
0154                     fBpPerVoxel[voxelName] = numOfBp;
0155                     std::cout<<voxelName<<" has "<<numOfBp<<" bp"<<std::endl;
0156                     foundNumOfBp = true;
0157                 }
0158             }
0159         }
0160         file.close();
0161     }
0162     
0163     if (fBpPerVoxel.size() == 0) {
0164         std::cerr<<"**** Fatal Error *****"<<std::endl;
0165         std::cerr<<"ScanDamage::RetrieveVoxelBp: No Bp found in voxel definition files. \n Or make"
0166         <<" sure that file imp.info exists in working directory!"<<std::endl;
0167         std::cerr<<"*************** *****"<<std::endl;
0168         exit(EXIT_FAILURE);
0169     }
0170 }
0171 
0172 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0173 
0174 void ScanDamage::FillVoxelData()
0175 {   
0176     std::ifstream file(fCellDefFilePath);
0177     if(!file.good() )
0178     {
0179         std::cerr<<"**** Fatal Error *****"<<std::endl;
0180         std::cerr<<"FillVoxelData: No file named "<<fCellDefFilePath<<std::endl;
0181         std::cerr<<"*************** *****"<<std::endl;
0182         exit(EXIT_FAILURE);
0183     }
0184     ullint bpCount = 0;
0185     unsigned int voxelCount = 0;
0186     int chromo_previous = 0;
0187     // Read the file line by line
0188     std::string line;
0189     while(std::getline(file, line) )
0190     {
0191         std::istringstream iss(line);
0192         std::string flag;
0193         iss >> flag;
0194         // If the flag correspond to the placement of a voxel
0195         if(flag == "_pl")
0196         {
0197             std::string voxelName;
0198             iss >> voxelName;
0199             int chromo;
0200             iss >> chromo;
0201             int domain;
0202             iss >> domain;
0203             // If we change of chromosome then reset the number of bp.
0204             // Each chromosome starts at 0 bp.
0205             if(chromo != chromo_previous)
0206             {
0207                 bpCount = 0;
0208                 chromo_previous = chromo;
0209             }
0210             // Fill the data structure
0211             fVoxels.push_back( VoxelData(chromo, domain, bpCount) );
0212             int numBpinthisvoxel = int(fBpPerVoxel[voxelName]);
0213             bpCount += numBpinthisvoxel;
0214             if (fChromosomeBpMap.find(chromo) == fChromosomeBpMap.end()) {
0215                 fChromosomeBpMap.insert({chromo,numBpinthisvoxel});
0216             } else {
0217                 fChromosomeBpMap[chromo] += numBpinthisvoxel;
0218             }
0219             voxelCount++;
0220         }
0221     }
0222     file.close();
0223     
0224 
0225     if (fVoxels.size() == 0) {
0226         std::cerr<<"**** Fatal Error *****"<<std::endl;
0227         std::cerr<<"ScanDamage::FillVoxelData: NofVoxels info found in files "<<fCellDefFilePath<<std::endl;
0228         std::cerr<<"*************** *****"<<std::endl;
0229         exit(EXIT_FAILURE);
0230     }
0231     std::cout<<"Num of voxels: "<< fVoxels.size()<<" placed in Cell Nucleus."<<std::endl;
0232     std::cout<<"=====Choromosome sizes====="<<std::endl;
0233     std::cout<<"Chromosome ID\t Number of Bp"<<std::endl;
0234     for (auto const& [chrom, nBp] : fChromosomeBpMap) {
0235          std::cout<<chrom<< "\t"<<nBp<<std::endl;
0236     }
0237     std::cout<<"==========================="<<std::endl;
0238 }
0239 
0240 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0241 
0242 void ScanDamage::ScanDamageFromPhys()
0243 {
0244     std::cout<<"===== Start Scanning Damages From Phys =====\n";
0245     fEdepSumInNucleus = 0;
0246     fphysTables.clear();
0247     fphysSlectedTables.clear();
0248     fs::path currentP{"phys_output"};
0249     fs::file_status s = fs::file_status{};
0250     auto isExist = fs::status_known(s) ? fs::exists(s) : fs::exists(currentP);
0251     if (isExist) {
0252         bool isFoundRootFiles = false;
0253         for (const auto entry : fs::directory_iterator(currentP)) {
0254             if (entry.path().extension() == ".root") {
0255                 std::cout <<"ScanDamageFromPhys(): Processing file: "<< entry.path().filename()<< std::endl;
0256                 AnaPhysRootFile(entry.path());
0257                 if (!isFoundRootFiles) isFoundRootFiles=true;
0258             }
0259         }
0260         if (!isFoundRootFiles) {
0261             std::cout<<"=====>> No root files found in folder \"phys_ouput\"!!! Skip Scanning Damages From Phys =====\n";
0262         }
0263         if (fphysTables.size() > 0) {
0264             SortPhysTableWithSelection();
0265         }
0266     } else {
0267         std::cout<<"=====>> Cannot find folder \"phys_ouput\"!!! Skip Scanning Damages From Phys =====\n";
0268     }
0269 
0270     std::cout<<"===== End Scanning Damages From Phys =====\n";
0271 }
0272 
0273 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0274 
0275 void ScanDamage::ScanDamageFromChem()
0276 {
0277     std::cout<<"===== Start Scanning Damages From Chem =====\n";
0278     std::string fChemOutFolderName = ParametersParser::Instance()->GetChemOutFolderName();
0279     if (fChemOutFolderName == "") fChemOutFolderName = "chem_output";
0280     fs::path currentP{fChemOutFolderName};
0281     fs::file_status s = fs::file_status{};
0282     auto isExist = fs::status_known(s) ? fs::exists(s) : fs::exists(currentP);
0283     if (isExist) {
0284         bool isFoundRootFiles = false;
0285         for (const auto entry : fs::directory_iterator(currentP)) {
0286             if (entry.path().extension() == ".root") {
0287                 AnaChemRootFile(entry);
0288                 if (!isFoundRootFiles) isFoundRootFiles=true;
0289             }
0290         }
0291         if (!isFoundRootFiles) {
0292             std::cout<<"=====>> No root files found in folder \""<<fChemOutFolderName<<"\"!!! Skip Scanning Damages From Chem =====\n";
0293             fSkipScanningIndirectDamage = true;
0294         }
0295         if (fchemTables.size() > 0) {
0296             SortChemTableWithSelection();
0297         }
0298     } else {
0299         std::cout<<"=====>> Cannot find folder \""<<fChemOutFolderName<<"\"!!! Skip Scanning Damages From Chem =====\n";
0300         fSkipScanningIndirectDamage = true;
0301     }
0302     
0303     std::cout<<"===== End Scanning Damages From Chem =====\n";
0304 }
0305 
0306 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0307 
0308 void ScanDamage::AnaPhysRootFile(const std::string fileName)
0309 {
0310     TFile* f = new TFile(fileName.c_str());
0311     if(f->IsZombie() ){
0312         // File is corrupted
0313         std::cerr<<"*********** Warning *************"<<std::endl;
0314         std::cerr<<"The file "<<fileName<<" seems to be corrupted..."<<std::endl;
0315         std::cerr<<"We will skip it."<<std::endl;
0316         std::cerr<<"**********************************"<<std::endl;
0317         return;
0318     }
0319     AnaPhysRootTree1(f);
0320     AnaPhysRootTree2(f);
0321     f->Close();
0322 }
0323 
0324 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0325 
0326 void ScanDamage::AnaChemRootFile(fs::directory_entry entry)
0327 {
0328     std::string fnameWithoutExtension = entry.path().stem().string();
0329     auto [eventNumber, voxelNumber] = GetEventNberAndVoxelNberFromChemRoot(fnameWithoutExtension);
0330     
0331     // Voxel data retrieval
0332     VoxelData& voxelData = fVoxels.at(voxelNumber);
0333 
0334     int chromo = voxelData.fChromosome;
0335     ullint firstBpNum = voxelData.fFirstBpCopyNum;
0336     // *******************
0337     // Analyse of the ntuple to detect SB and associate them with a bpNumCorrected
0338     // *******************
0339 
0340     // Load the file, the directory and the ntuple
0341     TFile f(entry.path().c_str());
0342     if(f.IsZombie() ){
0343         corruptedFiles++;
0344         // File is corrupted
0345         std::cerr<<"*********** Warning *************"<<std::endl;
0346         std::cerr<<"The file "<<entry.path().string()<<" seems to be corrupted..."<<std::endl;
0347         std::cerr<<"We will skip it."<<std::endl;
0348         std::cerr<<"Number of corrupted files: "<< corruptedFiles<<std::endl;
0349         std::cerr<<"**********************************"<<std::endl;
0350     } else {
0351         TDirectoryFile *d = dynamic_cast<TDirectoryFile*> (f.Get("ntuple") );
0352         TTree* chemTree = (TTree*) d->Get("ntuple_2");
0353 
0354         if( (int) chemTree->GetEntries() >0)
0355         {
0356             int strand;
0357             int copyNumber;
0358             double xp;
0359             double yp;
0360             double zp;
0361             double time;
0362             int base;
0363             chemTree->SetBranchAddress("strand", &strand);
0364             chemTree->SetBranchAddress("copyNumber", &copyNumber);
0365             chemTree->SetBranchAddress("xp", &xp);
0366             chemTree->SetBranchAddress("yp", &yp);
0367             chemTree->SetBranchAddress("zp", &zp);
0368             chemTree->SetBranchAddress("time", &time);
0369             chemTree->SetBranchAddress("base", &base);
0370             unsigned int entryNumber = (int) chemTree->GetEntries();
0371             for (unsigned int e=0;e<entryNumber;e++)
0372             {
0373                 chemTree->GetEntry(e);
0374                 ullint cpNumCorrected = firstBpNum+int(copyNumber);
0375                 std::vector<ullint> newLine{ 
0376                     (ullint)eventNumber,(ullint)strand,cpNumCorrected,
0377                     (ullint)base,(ullint)(time*1000000)};
0378                 auto itr = fchemTables.find(chromo);
0379                 if ( itr == fchemTables.end()) {
0380                     Table tbforThisChro{newLine};
0381                     fchemTables.insert({chromo,tbforThisChro});
0382                 } else (itr->second).push_back(newLine);
0383             }
0384         }
0385     }//
0386 }
0387 
0388 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0389 
0390 std::tuple<unsigned int, unsigned int> ScanDamage::GetEventNberAndVoxelNberFromChemRoot(
0391     const std::string fileNameWithoutExtension)
0392 {
0393     unsigned int evnN, volxelN;
0394     auto fristPos = fileNameWithoutExtension.find_first_of("_");
0395     auto secondPos = fileNameWithoutExtension.substr(fristPos+1).find_first_of("_");
0396     auto lastPos = fileNameWithoutExtension.find_last_of("_");
0397     evnN = std::stoul(fileNameWithoutExtension.substr(fristPos+1,secondPos));
0398     volxelN = std::stoul(fileNameWithoutExtension.substr(lastPos+1));
0399     return  {evnN, volxelN};
0400 }
0401 
0402 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0403 
0404 void ScanDamage::AnaPhysRootTree1(TFile* f)
0405 {
0406     TDirectoryFile* d = dynamic_cast<TDirectoryFile*> (f->Get("ntuple") );
0407     TTree* tPhys = dynamic_cast<TTree*> (d->Get("ntuple_1") );
0408 
0409     if( tPhys->GetEntries()  > 0)
0410         {
0411         int flagParticle;
0412         int flagParentID;
0413         int flagProcess;
0414         double x;
0415         double y;
0416         double z;
0417         double edep;
0418         int eventNumber;
0419         int volumeName;
0420         int copyNumber;
0421         int lastMetVoxelCopyNum;
0422 
0423         tPhys->SetBranchAddress("flagParticle", &flagParticle);
0424         tPhys->SetBranchAddress("flagParentID", &flagParentID);
0425         tPhys->SetBranchAddress("flagProcess", &flagProcess);
0426         tPhys->SetBranchAddress("x", &x);
0427         tPhys->SetBranchAddress("y", &y);
0428         tPhys->SetBranchAddress("z", &z);
0429         tPhys->SetBranchAddress("edep", &edep);
0430         tPhys->SetBranchAddress("eventNumber", &eventNumber);
0431         tPhys->SetBranchAddress("volumeName", &volumeName);
0432         tPhys->SetBranchAddress("copyNumber", &copyNumber);
0433         tPhys->SetBranchAddress("lastMetVoxelCopyNum", &lastMetVoxelCopyNum);
0434 
0435         unsigned int entryNumber =  tPhys->GetEntries() ;
0436 
0437         // Loop on all the "lines" (ie entry) of the ntuple
0438         for(unsigned int e=0; e<entryNumber; e++)
0439         {
0440             // Set all the variables to the values corresponding to the entry number
0441             tPhys->GetEntry(e);
0442             // Check if the process is an ionisation
0443             // Only ionisation should trigger the removal of a DNA molecule from the chemical step
0444             if(flagProcess == 13 // e-_DNAIonisation
0445                     || flagProcess == 113 // e-_DNAPTBIonisation
0446                     || flagProcess == 18 // proton_DNAIonisation
0447                     || flagProcess == 21 // hydrogen_DNAIonisation
0448                     || flagProcess == 24 // alpha_DNAIonisation
0449                     || flagProcess == 27 // alpha+_DNAIonisation
0450                     || flagProcess == 31 // helium_DNAIonisation
0451                     || flagProcess == 12 // e-_DNAExcitation
0452                     || flagProcess == 112 // e-_DNAPTBExcitation
0453                     || flagProcess == 15 // e-_DNAVibExcitation
0454                     || flagProcess == 17 // proton_DNAExcitation
0455                     || flagProcess == 20 // hydrogen_DNAExcitation
0456                     || flagProcess == 23 // alpha_DNAExcitation
0457                     || flagProcess == 26 // alpha+_DNAExcitation
0458                     || flagProcess == 30 // helium_DNAExcitation
0459                     ) {
0460                 // Check the interaction happened in a dna molecule or its hydration shell
0461                 if(volumeName == 1 // d1
0462                         || volumeName == 11 // p1
0463                         || volumeName == 2 // d2
0464                         || volumeName == 22 // p2
0465                         || volumeName == 7 // d1_w
0466                         || volumeName == 71 // p1_w
0467                         || volumeName == 8 // d2_w
0468                         || volumeName == 81 // p2_w
0469                         )
0470                 {
0471                     // *************
0472 
0473                     // Retrieve the voxel copy number
0474                     double voxelCopyNumber  = lastMetVoxelCopyNum;
0475                     if (voxelCopyNumber >= 0 && voxelCopyNumber<fVoxels.size()){
0476                         // Chromosome, domain and firstNucleotideNum
0477                         const VoxelData& voxelData = fVoxels.at(size_t(voxelCopyNumber) );
0478                         int chromo = voxelData.fChromosome;
0479                         int domain = voxelData.fDomain;
0480                         ullint firstBpCN = voxelData.fFirstBpCopyNum;
0481                         ullint cpNumCorrected = firstBpCN+int(copyNumber);
0482 
0483                         // Get the event number
0484                         double eventNum = eventNumber;
0485 
0486                         // Determine the strand
0487                         double strand (-1);
0488                         if(volumeName==1
0489                                 || volumeName==11
0490                                 || volumeName==7
0491                                 || volumeName==71
0492                                 || volumeName==6 // ade
0493                                 || volumeName==9 // ade
0494                                 || volumeName==4 // gua
0495                                 || volumeName==10) // gua
0496                             strand = 1;
0497                         else if(volumeName==2
0498                                 || volumeName==22
0499                                 || volumeName==8
0500                                 || volumeName==81
0501                                 || volumeName==5 // thy
0502                                 || volumeName==12 // thy
0503                                 || volumeName==3 // cyto
0504                                 || volumeName==13) // cyto
0505                             strand = 2;
0506                         // Check if the chromo has already been registered
0507                         std::vector<ullint> newLine{ 
0508                             (ullint)eventNum,(ullint)strand,cpNumCorrected,
0509                             (ullint)volumeName,(ullint)flagProcess,(ullint)(edep*1000000)};
0510                             // *1000000 in edep is taken back after. This is done 
0511                             // because the number is "unsigned long long int" instead of "double".
0512                         auto itr = fphysTables.find(chromo);
0513                         if ( itr == fphysTables.end()) {
0514                             Table tbforThisChro{newLine};
0515                             fphysTables.insert({chromo,tbforThisChro});
0516                         } else (itr->second).push_back(newLine);
0517                     }
0518                 }
0519             }
0520         }
0521     }
0522 }
0523 
0524 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0525 
0526 void ScanDamage::AnaPhysRootTree2(TFile* f)
0527 {
0528     TDirectoryFile* d2 = dynamic_cast<TDirectoryFile*> (f->Get("ntuple") );
0529     TTree* tPhys2 = dynamic_cast<TTree*> (d2->Get("ntuple_3") );
0530 
0531     if( int(tPhys2->GetEntries() ) > 0)
0532     {
0533         double edep;
0534         int eventNumber;
0535 
0536         tPhys2->SetBranchAddress("edep", &edep);
0537         tPhys2->SetBranchAddress("eventNumber", &eventNumber);
0538 
0539         unsigned int entryNumber = int( tPhys2->GetEntries() );
0540 
0541         // Loop on all the "lines" (ie entry) of the ntuple
0542         for(unsigned int e=0; e<entryNumber; e++)
0543         {
0544             tPhys2->GetEntry(e);
0545             fEdepSumInNucleus += edep;
0546         }
0547     }
0548 }
0549 
0550 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0551 
0552 void ScanDamage::SortPhysTableWithSelection()
0553 {
0554     for(const auto [chrom, physTable] : fphysTables)
0555     {
0556         // Final table
0557         Table physTableWithSelection;
0558 
0559         std::map<ullint,std::map<ullint,std::map<ullint, ullint > > > energyMap;
0560         
0561         // Loop on all the lines of the table
0562         for(unsigned int line=0, eline=physTable.size(); line<eline; ++line)
0563         {
0564             ullint eventNum = physTable[line][0];
0565             ullint strand = physTable[line][1];
0566             ullint copyNumber = physTable[line][2];
0567             ullint volumeFlag = physTable[line][3];
0568             ullint processFlagr = physTable[line][4];
0569             ullint energy = physTable[line][5];
0570 
0571             // Cumulate the energy value
0572             energyMap[eventNum][strand][copyNumber] += energy;
0573         }
0574 
0575         int notDuplicatedLine = 0;
0576 
0577         // Loop on all the events
0578         std::map<ullint, std::map<ullint, std::map<ullint, ullint> > >::iterator iit = energyMap.begin();
0579         std::map<ullint, std::map<ullint, std::map<ullint, ullint> > >::iterator iite = energyMap.end();
0580         for(; iit!=iite;++iit)
0581         {
0582             ullint eventNum = iit->first;
0583 
0584             // Loop on all the strands
0585             std::map<ullint, std::map<ullint, ullint> >::iterator itt = iit->second.begin();
0586             std::map<ullint, std::map<ullint, ullint> >::iterator itte = iit->second.end();
0587             for(; itt!=itte;++itt)
0588             {
0589                 ullint strand = itt->first;
0590                 // Loop on all the copy numbers
0591                 std::map<ullint, ullint>::iterator ittt = itt->second.begin();
0592                 std::map<ullint, ullint>::iterator ittte = itt->second.end();
0593                 for(; ittt!=ittte;++ittt)
0594                 {
0595                     ullint copyNumber = ittt->first;
0596                     double currentE = double(ittt->second) / 1000000; // eV
0597                     // Energy condition(s) are set here
0598                     bool fill = false;                  
0599                      // Threshold condition
0600                      if(currentE < fThresholdEnergy)
0601                           fill=false;
0602                      else
0603                           fill=true;
0604 
0605                     if(fill)
0606                     {
0607                         // Add a line
0608                         physTableWithSelection.push_back(std::vector<ullint>());
0609                         // Fill the line
0610                         physTableWithSelection[notDuplicatedLine].push_back(eventNum);
0611                         physTableWithSelection[notDuplicatedLine].push_back(strand);
0612                         physTableWithSelection[notDuplicatedLine].push_back(copyNumber);
0613                         physTableWithSelection[notDuplicatedLine].push_back(0);
0614                         physTableWithSelection[notDuplicatedLine].push_back(0);
0615                         physTableWithSelection[notDuplicatedLine].push_back((ullint)(currentE*1000000) );
0616                         physTableWithSelection[notDuplicatedLine].push_back(0);
0617                         ++notDuplicatedLine;
0618                     }
0619                 }
0620             }
0621         }
0622         // *******************************************
0623         // Print the "physTableWithSelection" table for the current chromosome
0624         // *******************************************
0625         std::cout << "### Phys SB for chromosome "<<chrom<<" : " << physTableWithSelection.size() << " ###" << std::endl;
0626         if (physTableWithSelection.size()>0) fphysSlectedTables.insert({chrom,physTableWithSelection});
0627     }
0628     fphysTables.clear();//Free memory
0629 }
0630 
0631 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0632 
0633 void ScanDamage::SortChemTableWithSelection()
0634 {
0635     for (const auto& [chromo,chemTable] : fchemTables)
0636     {
0637         // Final table
0638         Table chemTableWithSelection;
0639         
0640         int notDuplicatedLine = 0;
0641 
0642         for(unsigned int line=0, eline=chemTable.size(); line<eline; ++line)
0643         {
0644             ullint eventNum = chemTable[line][0];
0645             ullint strand = chemTable[line][1];
0646             ullint copyNumber = chemTable[line][2];
0647             ullint base = chemTable[line][3];
0648             ullint time = chemTable[line][4];
0649         
0650             // Random number between 0 and <1
0651             if (base==1) {
0652                 // Add a line
0653                 chemTableWithSelection.push_back(std::vector<ullint>());
0654                 // Fill the line
0655                 chemTableWithSelection[notDuplicatedLine].push_back(eventNum);
0656                 chemTableWithSelection[notDuplicatedLine].push_back(strand);
0657                 chemTableWithSelection[notDuplicatedLine].push_back(copyNumber);
0658                 chemTableWithSelection[notDuplicatedLine].push_back(base);
0659                 chemTableWithSelection[notDuplicatedLine].push_back(time);
0660                 chemTableWithSelection[notDuplicatedLine].push_back(0);
0661                 chemTableWithSelection[notDuplicatedLine].push_back(1);              
0662                 ++notDuplicatedLine;
0663             }
0664             else { 
0665                 //double r = double(std::rand() ) / RAND_MAX;
0666                 double r = gRandomGen.Rndm(); 
0667                 //std::cout<<r<<std::endl;
0668                 if(r <= fProbabilityForIndirectSB){
0669                     // Add a line
0670                     chemTableWithSelection.push_back(std::vector<ullint>());
0671                     // Fill the line
0672                     chemTableWithSelection[notDuplicatedLine].push_back(eventNum);
0673                     chemTableWithSelection[notDuplicatedLine].push_back(strand);
0674                     chemTableWithSelection[notDuplicatedLine].push_back(copyNumber);
0675                     chemTableWithSelection[notDuplicatedLine].push_back(base);
0676                     chemTableWithSelection[notDuplicatedLine].push_back(time);
0677                     chemTableWithSelection[notDuplicatedLine].push_back(0);
0678                     chemTableWithSelection[notDuplicatedLine].push_back(1);                          
0679                     ++notDuplicatedLine;
0680                 }
0681             }
0682         }       
0683         if (chemTableWithSelection.size() > 0) fchemSlectedTables.insert({chromo,chemTableWithSelection});
0684     }
0685     fchemTables.clear();//free memory
0686 }
0687 
0688 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0689 
0690 void ScanDamage::MergeDamageFromPhysChem()
0691 {
0692     for(int  chromo=0; chromo<46; chromo++)
0693     {
0694         // *****************************
0695         // Add one table after the other
0696         // *****************************
0697 
0698         // MergedTable to be built
0699         Table mergedTable;
0700 
0701         auto itr = fphysSlectedTables.find(chromo);
0702         if(itr != fphysSlectedTables.end())
0703         {
0704             Table physTable= itr->second;
0705             for(auto const& vec : physTable) mergedTable.push_back(vec);
0706         }
0707         itr = fchemSlectedTables.find(chromo);
0708         if(itr != fchemSlectedTables.end())
0709         {
0710             Table chemTable = itr->second;
0711             for(auto const& vec : chemTable) mergedTable.push_back(vec);
0712         }
0713 
0714         // *******************************************
0715         // Sort the merged table to put the event in the correct order
0716         // *******************************************
0717         Trier tri;
0718         Table::iterator it = mergedTable.begin();
0719         Table::iterator ite = mergedTable.end();
0720         std::sort(it, ite, tri);
0721 
0722         // *******************************************
0723         // Delete duplicate SB
0724         // *******************************************
0725 
0726         std::map<ullint,std::map<ullint,std::map<ullint,std::map<ullint,std::vector<ullint>>>>> removeDuplicateValueMap;
0727 
0728         // Put all the values of the mergedTable in the map created just above.
0729         // Loop on all the lines of the mergedTable
0730         for(unsigned int line=0; line<mergedTable.size(); ++line)
0731         {
0732             ullint eventNum = mergedTable[line][0];
0733             ullint strand = mergedTable[line][1];
0734             ullint copyNumber = mergedTable[line][2];
0735             ullint base = mergedTable[line][3];
0736             std::vector<ullint> lineV;
0737             // If more elements are presents, add them here
0738             int lineSize = mergedTable[line].size();
0739             if(lineSize > 4)
0740             {
0741                 for(int i=4; i<lineSize; i++)
0742                 {
0743                     lineV.push_back(mergedTable[line][i]);
0744                 }
0745                 removeDuplicateValueMap[eventNum][strand][copyNumber][base]=lineV;
0746             }
0747             else
0748             {
0749                 lineV.push_back(0);
0750                 removeDuplicateValueMap[eventNum][strand][copyNumber][base]=lineV;
0751             }
0752         }
0753 
0754         // *******************************************
0755         // Create the "mergedTableWithoutDuplicatedSB" table
0756         // *******************************************
0757 
0758         // At this point, we have created a map named "removeDuplicateValueMap" that organized all the mergedTable content
0759         // AND that does not contain any duplicate because of the override caracteristic of a map.
0760         // Indeed, doing map[2] = "hello" followed by map[2]  = "bye" will put map[2] value to "bye" because it overrided the first "hello".
0761         // This is a cheap way to remove duplicates by overriding them.
0762 
0763         // The next part is dedicated to the creation of the "mergedTableWithoutDuplicatedSB" table from the "removeDuplicateValueMap" map.
0764         // Only the event, strand and copynumber and  will be put in this final table.
0765         Table mergedTableWithoutDuplicatedSB;
0766         Table mergedTableWithoutDuplicatedSBandbases;
0767         int notDuplicatedLine = 0;
0768         int notDuplicatedLine2= 0;
0769         // Loop on all the events
0770         std::map<ullint,std::map<ullint,std::map<ullint,std::map<ullint,std::vector<ullint>>>>>::iterator 
0771             iit = removeDuplicateValueMap.begin();
0772         std::map<ullint,std::map<ullint,std::map<ullint,std::map<ullint,std::vector<ullint>>>>>::iterator 
0773             iite = removeDuplicateValueMap.end();
0774         for(; iit!=iite;++iit)
0775         {
0776             ullint eventNum = iit->first;
0777             // Loop on all the strands
0778             std::map<ullint, std::map<ullint, std::map<ullint, std::vector<ullint > > > >::iterator 
0779                 itt = iit->second.begin();
0780             std::map<ullint, std::map<ullint, std::map<ullint, std::vector<ullint > > > >::iterator 
0781                 itte = iit->second.end();
0782             for(; itt!=itte;++itt)
0783             {
0784                 ullint strand = itt->first;
0785                 // Loop on all the copy numbers
0786                 std::map<ullint, std::map<ullint, std::vector<ullint > > >::iterator ittt = itt->second.begin();
0787                 std::map<ullint, std::map<ullint, std::vector<ullint > > >::iterator ittte = itt->second.end();
0788                 for(; ittt!=ittte;++ittt)
0789                 {
0790                     ullint copyNumber = ittt->first;                    
0791                     // Loop on all the base flag
0792                     std::map<ullint, std::vector<ullint > >::iterator itttt = ittt->second.begin();
0793                     std::map<ullint, std::vector<ullint > >::iterator itttte = ittt->second.end();
0794                     for(; itttt!=itttte;++itttt)
0795                     {
0796                         ullint base = itttt->first;
0797                         // Fill the table
0798                         if (strand>0 && strand<3)
0799                         {
0800                             // Add a line
0801                             mergedTableWithoutDuplicatedSB.push_back(std::vector<ullint>());                           
0802                             // Fill the line
0803                             mergedTableWithoutDuplicatedSB[notDuplicatedLine].push_back(eventNum);
0804                             mergedTableWithoutDuplicatedSB[notDuplicatedLine].push_back(strand);
0805                             mergedTableWithoutDuplicatedSB[notDuplicatedLine].push_back(copyNumber);
0806                             mergedTableWithoutDuplicatedSB[notDuplicatedLine].push_back(base);
0807                             mergedTableWithoutDuplicatedSB[notDuplicatedLine].push_back(
0808                                 removeDuplicateValueMap[eventNum][strand][copyNumber][base][0]);
0809                             mergedTableWithoutDuplicatedSB[notDuplicatedLine].push_back(
0810                                 removeDuplicateValueMap[eventNum][strand][copyNumber][base][1]);
0811                             mergedTableWithoutDuplicatedSB[notDuplicatedLine].push_back(
0812                                 removeDuplicateValueMap[eventNum][strand][copyNumber][base][2]);
0813                             ++notDuplicatedLine;                           
0814                             if (base==0)
0815                             {
0816                                 // Add a line
0817                                 mergedTableWithoutDuplicatedSBandbases.push_back(std::vector<ullint>());                              
0818                                 // Fill the line
0819                                 mergedTableWithoutDuplicatedSBandbases[notDuplicatedLine2].push_back(eventNum);
0820                                 mergedTableWithoutDuplicatedSBandbases[notDuplicatedLine2].push_back(strand);
0821                                 mergedTableWithoutDuplicatedSBandbases[notDuplicatedLine2].push_back(copyNumber);                               
0822                                 ++notDuplicatedLine2;
0823                             }
0824                         }
0825                     }
0826                 }
0827             }
0828         }
0829 
0830         // *******************************************
0831         // Print the "mergedTableWithoutDuplicatedSB" table for the current chromosome
0832         // *******************************************
0833         if (mergedTableWithoutDuplicatedSB.size() > 0) {
0834             //PrintTable(fMergeFolder + "/chromo_"+std::to_string(chromo)+".dat", 
0835             //mergedTableWithoutDuplicatedSB, "eventNum, strand, copyNumber, isbase, 
0836             //time(ns*1000000), edep, phy:0 chem:1");
0837             fMergedTables.insert({chromo,mergedTableWithoutDuplicatedSB});
0838         }
0839         mergedTable.clear();
0840         mergedTableWithoutDuplicatedSBandbases.clear();
0841         mergedTableWithoutDuplicatedSB.clear();
0842     }
0843     //free memory:
0844     fphysSlectedTables.clear();
0845     fchemSlectedTables.clear();
0846 }
0847 
0848 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0849 
0850 void ScanDamage::ReadCellandVoxelDefFilePaths()
0851 {
0852     fs::path thisP = fs::current_path();
0853     for (const auto entry : fs::directory_iterator(thisP)){
0854         if (entry.path().filename() == "imp.info") {
0855             std::ifstream file(entry.path().c_str());
0856             if(!file.good() ){
0857                 std::cerr<<"**** Fatal Error *****"<<std::endl;
0858                 std::cerr<<"ScanDamage::ReadCellandVoxelDefFilePaths(): File corupted: "
0859                 <<entry.path()<<std::endl;
0860                 std::cerr<<"*************** *****"<<std::endl;
0861                 exit(EXIT_FAILURE);
0862             }
0863 
0864             std::string line;
0865             while(std::getline(file, line) ){
0866                 std::istringstream iss(line);
0867                 std::string flag;
0868                 iss >> flag;
0869                 if ( flag == "_geovolxelpath") {
0870                     std::string voxname;
0871                     iss >> voxname;
0872                     fVoxelDefFilesList.insert(voxname);
0873                 }
0874                 if ( flag == "_geocellpath") {
0875                     std::string cellpname;
0876                     iss >> cellpname;
0877                     fCellDefFilePath = cellpname;
0878                 }
0879                 if ( flag == "_numberOfBasepairs") {
0880                     iss >> fTotalNbBpPlacedInGeo;
0881                 }
0882                 if ( flag == "_numberOfHistones") {
0883                     iss >> fTotalNbHistonePlacedInGeo;
0884                 }
0885                 if ( flag == "_nucleusVolume") {
0886                     iss >> fNucleusVolume;
0887                 }
0888                 if ( flag == "_nucleusMassDensity") {
0889                     iss >> fNucleusMassDensity;
0890                 }
0891                 if ( flag == "_nucleusMass") {
0892                     iss >> fNucleusMass;
0893                 }
0894             }
0895             file.close();
0896         }
0897     }
0898 }
0899 
0900 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....