AIDA GELINA BRIKEN nToF CRIB ISOLDE CIRCE nTOFCapture DESPEC DTAS EDI_PSA 179Ta CARME StellarModelling DCF K40
  BRIKEN  ELOG logo
Fields marked with * are required
Entry time:Fri Nov 11 08:07:18 2016
Author*:
Subject*:

Encoding:
        
Attachment 1: HybridBRIKENPhysicalHistogram.cpp   
//
//  HybridBRIKENPhysicalHistogram.cpp
//  
//
//  Created by Bertis Charles Rasco on 11/4/16.
//
//

#include "HybridBRIKENPhysicalHistogram.hpp"

#include <cmath>
#include <fstream>
#include <iostream>
#include <sstream> 
#include <string>

#include <stdio.h>
//#include <stdlib.h>

#include "TH1F.h"
#include "TRandom3.h"
#include "TFile.h"

// ROOT Needs this for some reason.
HybridBRIKENHistogram::HybridBRIKENHistogram(): TH2Poly()
{
}

HybridBRIKENHistogram::HybridBRIKENHistogram( const char *name, const char *title, Double_t xlow, Double_t xup, Double_t ylow, Double_t yup ):
  TH2Poly( name, title, xlow, xup, ylow, yup )
{
  this->ReadinTubes();
  
  this->FillTubeCountsRandom();
  
//  this->SetupTubes();
//  this->Draw();
}

HybridBRIKENHistogram::HybridBRIKENHistogram( const char *name, const char *title, Double_t xlow, Double_t xup, Double_t ylow, Double_t yup, TString inputFileName ):
  TH2Poly( name, title, xlow, xup, ylow, yup )
{
  m_InputFileName = inputFileName;
  
  
  this->ReadinTubes();
  this->ReadInDataFromFile();
//  this->SetupTubes();
//  this->Draw();
}


HybridBRIKENHistogram::~HybridBRIKENHistogram()
{
}

void HybridBRIKENHistogram::ReadinTubes( void )
{
  std::ifstream file("HybridBRIKEN_Tube_Info_1.txt");
  std::string str; 
  
  Int_t lineCount = 0;
  Int_t tubeCount = 0;
  Double_t rad;
  Double_t pos[2];
  Int_t tubeID;
  Int_t ringID;
  Char_t tubeType[20];
  
//  std::stringstream ss;
//  std::string::size_type sz;
  
  while ( std::getline(file, str) )
  {
    // Process str
    std::stringstream ss( str);
    if( lineCount > 0 )
    {
      tubeType[0] = 'X';
      
      ss >> tubeID >> tubeType >> rad >> pos[0] >> pos[1] >> ringID;
      rad *= 0.5;// 2R is what is in file.
      
//      tubeID = std::stoi( str, &sz );
      if( tubeID > 0 && (tubeType[0]=='R' || tubeType[0]=='U' || tubeType[0]=='O' ) )
      {
        if( tubeType[0] != 'O' ) m_Tubes.push_back( new HeTube( rad, pos[0], pos[1], tubeID, ringID, tubeType ) );
        else if( tubeType[0] == 'O' )
        {
          tubeType[1] = tubeType[4];
          m_Tubes.push_back( new HeTube( rad, pos[0], pos[1], tubeID, ringID, tubeType ) );
        }
        tubeCount++;
      }
    }
    
//    if( lineCount < 11 )
//    {
//      std::cout << str << std::endl;
//      std::cout << tubeID << ", " << tubeType[0] << tubeType[1] << ", " << rad << ", " << pos[0] << ", " << pos[1] << ", " << ringID << std::endl;
//      ss >> tubeID >> tubeType >> rad >> pos[0] >> pos[1] >> ringID;
//    }
    
    lineCount++;
  }

  Int_t binID;
  
//  Double_t minX = m_Tubes[0]->pos[0] - m_Tubes[0]->radius;
//  Double_t maxX = m_Tubes[0]->pos[0] + m_Tubes[0]->radius;
//  Double_t minY = m_Tubes[0]->pos[1] - m_Tubes[0]->radius;
//  Double_t maxY = m_Tubes[0]->pos[1] + m_Tubes[0]->radius;
// Fix   
//  binID = this->AddBin( minX, minY, maxX, maxY );
  binID = this->AddCircularBin( 0 );
  
//  std::cout << "X: " << this->GetNbinsX() << std::endl;
//  std::cout << "Y: " << this->GetNbinsY() << std::endl;
//  std::cout << "Z: " << this->GetNbinsZ() << std::endl;

//  std::cout << "BinID: " << binID << std::endl;
//  this->SetBinContent( binID, 9);

  for( Int_t i = 1; i < (Int_t) m_Tubes.size(); i++ )
  {
//    minX = m_Tubes[i]->pos[0] - m_Tubes[i]->radius;
//    maxX = m_Tubes[i]->pos[0] + m_Tubes[i]->radius;
//    minY = m_Tubes[i]->pos[1] - m_Tubes[i]->radius;
//    maxY = m_Tubes[i]->pos[1] + m_Tubes[i]->radius;
// Fix AddCircularBin
//    binID = this->AddBin( minX, minY, maxX, maxY );
    binID = this->AddCircularBin( i );
//    std::cout << "BinID: " << binID << std::endl;
//    this->SetBinContent( binID, 15 + binID * 11);
  }
  
  std::cout << "Number of Lines: " << lineCount << std::endl;// = How many tubes in total?
  std::cout << "Number of Tubes: " << binID << std::endl;// = How many tubes in total?
  
//  std::cout << "Number of Bins: " << this->GetNcells() << std::endl;// = How many tubes in total?
}

void HybridBRIKENHistogram::SetupTubes( void )
{
// 1" tubes
  const Int_t maxN1 = 8;
  for( Int_t i = 0; i < maxN1; i++ )
  {
    m_Tubes.push_back( new HeTube(0.5, 2.0 * cos( i * 6.28 / maxN1 ), 2.0 * sin(i * 6.28 / maxN1 ), i, 1, "RI" ) );
  }

// 2" tube at 14" radius
  const Int_t maxN = 3;
  for( Int_t i = 0; i < maxN; i++ )
  {
    m_Tubes.push_back( new HeTube(1.0, 14.0 * cos( i * 6.28 / maxN ), 14.0 * sin(i * 6.28 / maxN ), i, 2, "O2" ) );
  }
  
  Int_t binID;
  
  Double_t minX = m_Tubes[0]->pos[0] - m_Tubes[0]->radius;
  Double_t maxX = m_Tubes[0]->pos[0] + m_Tubes[0]->radius;
  Double_t minY = m_Tubes[0]->pos[1] - m_Tubes[0]->radius;
  Double_t maxY = m_Tubes[0]->pos[1] + m_Tubes[0]->radius;

// Fix this  
  binID = this->AddBin( minX, minY, maxX, maxY );
  
//  std::cout << binID << std::endl;
//  this->SetBinContent( binID, 9);
  
  for( Int_t i = 0; i < (Int_t) m_Tubes.size(); i++ )
  {
    minX = m_Tubes[i]->pos[0] - m_Tubes[i]->radius;
    maxX = m_Tubes[i]->pos[0] + m_Tubes[i]->radius;
    minY = m_Tubes[i]->pos[1] - m_Tubes[i]->radius;
    maxY = m_Tubes[i]->pos[1] + m_Tubes[i]->radius;
    
// Fix this  
    binID = this->AddBin( minX, minY, maxX, maxY );
  
//    std::cout << binID << std::endl;
//    this->SetBinContent( binID, 15 + binID * 11);
  }
}

void HybridBRIKENHistogram::ReadInDataFromFile( void )
{
  const Int_t lowHeCut = 600;//20;// Cut is by bin not by energy. But for 1 keV per bin they are the same.
  const Int_t highHeCut = 900;//200;// Cut is by bin not by energy. But for 1 keV per bin they are the same.
  
  m_InputFile = new TFile( m_InputFileName );
  
  Int_t binID;
  Double_t counts;

  char tempCh[20];
  for( Int_t i = 0; i < m_Tubes.size(); i++ )
  {
    binID = i + 1;
    TString histName("He");
    if( binID < 10 )
    {
      histName.Append("00");
    }
    else if( binID < 100 )
    {
      histName.Append("0");
    }
//    else if( binID < 1000 )
//    {
//    }
    
//    histName.Append( TString::Itoa( (Int_t) binID, 10 ) );// Works for ROOT 5.34, but not for 5.32
//      itoa(binID,tempCh,10);// Does not work for 5.32 unkown if works for 5.34
      sprintf (tempCh, "%d", binID);// works for 5.32
      histName.Append( TString( tempCh ) );  
//    std::cout << histName << std::endl;
    
    TH1F *histogram = (TH1F*)m_InputFile->Get( histName.Data() );
    counts = histogram->Integral( lowHeCut, highHeCut );
    
    this->SetBinContent( binID, counts);
  }
}

void HybridBRIKENHistogram::FillTubeCountsRandom( void )
{
  TRandom3 rand;
//  for( Int_t i = 0; i < this->GetNbins() ; i++)
//  {
//    hist->Fill( rand.Uniform(1.0,140.0) );
//  }
  for( Int_t i = 1; i <= (Int_t) m_Tubes.size(); i++ )
  {
    this->SetBinContent(i, rand.Gaus(1000., 50.0) );
  }
  
}

Int_t HybridBRIKENHistogram::AddCircularBin( Int_t tubeID )
{
  Int_t binID;
  Double_t r = m_Tubes[tubeID]->radius;
  Double_t theta;
  
  const Int_t N = 16;
  Double_t* x;
  Double_t* y;
  
  x = new Double_t[N];
  y = new Double_t[N];
  
  
  for( Int_t i = 0; i < N; i++ )
  {
    theta = i * 2.0 * 3.1415 / N;
    x[i] = m_Tubes[tubeID]->pos[0] + r * cos( theta );
    y[i] = m_Tubes[tubeID]->pos[1] + r * sin( theta );
  }
  
  binID = this->AddBin( N, x, y );
  
//  Double_t minX = m_Tubes[tubeID]->pos[0] - m_Tubes[tubeID]->radius;
//  Double_t maxX = m_Tubes[tubeID]->pos[0] + m_Tubes[tubeID]->radius;
//  Double_t minY = m_Tubes[tubeID]->pos[1] - m_Tubes[tubeID]->radius;
//  Double_t maxY = m_Tubes[tubeID]->pos[1] + m_Tubes[tubeID]->radius;
//
////  
//  binID = this->AddBin( minX, minY, maxX, maxY );
  
  return binID;
}

//***************************************************************************************
// Main Program...
void HybridBRIKENPhysicalHistogram( void )
{
  Double_t halfHDPEWidth = 500.0;// mm
//  TString inputFileName( "161103_1426_152Eu.root" );
//  TString inputFileName( "161103_1319_252Cf.root" );
//  TString inputFileName( "161107_0833_Overnight.root" );
//  TString inputFileName( "161108_1609_ShortBackground.root" );
//  TString inputFileName( "161013_1518_137Cs60Co.root" );
  TString inputFileName( "161108_2323_717273Ni.root" );

  
  HybridBRIKENHistogram *hist = new HybridBRIKENHistogram("hBR", "HybridBRIKEN", -halfHDPEWidth, halfHDPEWidth, -halfHDPEWidth, halfHDPEWidth, inputFileName );
//  HybridBRIKENHistogram *hist = new HybridBRIKENHistogram("hBR", "HybridBRIKEN", -halfHDPEWidth, halfHDPEWidth, -halfHDPEWidth, halfHDPEWidth );
    
  hist->Draw("COLZ");  
}
Attachment 2: HybridBRIKENPhysicalHistogram.hpp   
//
//  HybridBRIKENPhysicalHistogram.hpp
//  
//
//  Created by Bertis Charles Rasco on 11/4/16.
//
//

#ifndef HybridBRIKENPhysicalHistogram_hpp
#define HybridBRIKENPhysicalHistogram_hpp

#include <vector>

#include "TFile.h"
#include "TH2Poly.h"
#include "TObject.h"
#include "TString.h"

//class TFile;

class HeTube
{
public:
 HeTube( Double_t r, Double_t x, Double_t y, Int_t tID, Int_t rID, char tubeTypeChars[2] )
  { radius = r; pos[0] = x; pos[1] = y; tubeID = tID; ringID = rID; 
    tubeType[0] = tubeTypeChars[0]; tubeType[1] = tubeTypeChars[1]; }

//private:
  Double_t radius;
  Double_t pos[2];
  Int_t tubeID;
  Int_t ringID;
  Char_t tubeType[2]; // RI (RIKEN), UP (UPC), O1 (ORNL1), O2 (ORNL2)
};


class HybridBRIKENHistogram : public TH2Poly
{
public:
  HybridBRIKENHistogram();// ROOT needs this to be defined if it is inherited from a ROOT class.
  HybridBRIKENHistogram( const char *name, const char *title, Double_t xlow, Double_t xup, Double_t ylow, Double_t yup );
  HybridBRIKENHistogram( const char *name, const char *title, Double_t xlow, Double_t xup, Double_t ylow, Double_t yup,
    TString inputFileName );

  ~HybridBRIKENHistogram();
  
  ClassDef(HybridBRIKENHistogram,1);// ROOT also needs this if it is inherited from a ROOT class. Not sure what the 1 means.
private:

  void ReadinTubes( void );
  void SetupTubes( void );
  void ReadInDataFromFile( void );

  void FillTubeCountsRandom( void );

  Int_t AddCircularBin( Int_t tubeID );
  
  TString m_InputFileName;
  TFile *m_InputFile;
  
  std::vector<HeTube*> m_Tubes;
};

#endif /* HybridBRIKENPhysicalHistogram_hpp */
Attachment 3: HybridBRIKEN_Tube_Info_1.txt   
ID Type 2R(mm) X(mm) Y(mm) RING
1 RIKEN 27.5 -77.45 -39.25 1
2 RIKEN 27.5 -77.45 0 1
3 RIKEN 27.5 -77.45 39.25 1
4 RIKEN 27.5 -106.052 -19.625 2
5 RIKEN 27.5 -106.052 19.625 2
6 RIKEN 27.5 -134.655 0 3
7 RIKEN 27.5  77.45 -39.25 1
8 RIKEN 27.5  77.45 0 1
9 RIKEN 27.5  77.45 39.25 1
10 RIKEN 27.5  106.052 -19.625 2
11 RIKEN 27.5  106.052 19.625 2
12 RIKEN 27.5  134.655 0 3
13 RIKEN 27.5 -177.45 -39.25 1
14 RIKEN 27.5 -177.45 0 1
15 RIKEN 27.5 -177.45 39.25 1
16 RIKEN 27.5 -206.052 -19.625 2
17 RIKEN 27.5 -206.052 19.625 2
18 RIKEN 27.5 -234.655 0 3
19 RIKEN 27.5  177.45 -39.25 1
20 RIKEN 27.5  177.45 0 1
21 RIKEN 27.5  177.45 39.25 1
22 RIKEN 27.5  206.052 -19.625 2
23 RIKEN 27.5  206.052 19.625 2
24 RIKEN 27.5  234.655 0 3
25 UPC 27.5 -118.3 77.45 3
26 UPC 27.5 -69.7 77.45 1
27 UPC 27.5 -34.85 77.45 1
28 UPC 27.5 0 77.45 1
29 UPC 27.5 34.85 77.45 1
30 UPC 27.5 69.7 77.45 1
31 UPC 27.5 118.3 77.45 3
32 UPC 27.5 -97.2 116.7 3
33 UPC 27.5 -60.362 112.3 2
34 UPC 27.5 -17.425 107.631 2
35 UPC 27.5 17.425 107.631 2
36 UPC 27.5 60.362 112.3 2
37 UPC 27.5 97.2 116.7 3
38 UPC 27.5 -73.928 151.378 3
39 UPC 27.5 -34.85 137.812 3
40 UPC 27.5 0 147.15 3
41 UPC 27.5 34.85 137.812 3
42 UPC 27.5 73.928 151.378 3
43 UPC 27.5 -27.0595 178.438 3
44 UPC 27.5 27.0595 178.438 3
45 UPC 27.5 -118.3 -77.45 3
46 UPC 27.5 -69.7 -77.45 1
47 UPC 27.5 -34.85 -77.45 1
48 UPC 27.5 0 -77.45 1
49 UPC 27.5 34.85 -77.45 1
50 UPC 27.5 69.7 -77.45 1
51 UPC 27.5 118.3 -77.45 3
52 UPC 27.5 -97.2 -116.7 3
53 UPC 27.5 -60.362 -112.3 2
54 UPC 27.5 -17.425 -107.631 2
55 UPC 27.5 17.425 -107.631 2
56 UPC 27.5 60.362 -112.3 2
57 UPC 27.5 97.2 -116.7 3
58 UPC 27.5 -73.928 -151.378 3
59 UPC 27.5 -34.85 -137.812 3
60 UPC 27.5 0 -147.15 3
61 UPC 27.5 34.85 -137.812 3
62 UPC 27.5 73.928 -151.378 3
63 UPC 27.5 -27.0595 -178.438 3
64 UPC 27.5 27.0595 -178.438 3
65 ORNL1 27.5 -215.5 77.45 5
66 ORNL1 27.5 -166.9 77.45 4
67 ORNL1 27.5 166.9 77.45 4
68 ORNL1 27.5 215.5 77.45 5
69 ORNL1 27.5 -301.012 214.133 7
70 ORNL1 27.5 -252.393 235.169 7
71 ORNL1 27.5 252.393 235.169 7
72 ORNL1 27.5 301.012 214.133 7
73 ORNL1 27.5 -215.5 -77.45 5
74 ORNL1 27.5 -166.9 -77.45 4
75 ORNL1 27.5 166.9 -77.45 4
76 ORNL1 27.5 215.5 -77.45 5
77 ORNL1 27.5 -301.012 -214.133 7
78 ORNL1 27.5 -252.393 -235.169 7
79 ORNL1 27.5 252.393 -235.169 7
80 ORNL1 27.5 301.012 -214.133 7
81 ORNL2 53 -139.932 142.007 4
82 ORNL2 53 -70.014 198.718 4
83 ORNL2 53 0 217.478 4
84 ORNL2 53 70.014 198.718 4
85 ORNL2 53 139.932 142.007 4
86 ORNL2 53 -214.903 127.3 5
87 ORNL2 53 -133.044 210.494 5
88 ORNL2 53 -48.6976 259.192 5
89 ORNL2 53 48.6976 259.192 5
90 ORNL2 53 133.044 210.494 5
91 ORNL2 53 214.903 127.3 5
92 ORNL2 53 -273.203 102.2 6
93 ORNL2 53 -263.894 172.91 6
94 ORNL2 53 -204.385 193.709 6
95 ORNL2 53 -118.139 282.073 6
96 ORNL2 53 0 313.728 6
97 ORNL2 53 118.139 282.073 6
98 ORNL2 53 204.385 193.709 6
99 ORNL2 53 263.894 172.91 6
100 ORNL2 53 273.203 102.2 6
101 ORNL2 53 -339.453 102.2 7
102 ORNL2 53 -323.857 159.305 7
103 ORNL2 53 -201.103 278.553 7
104 ORNL2 53 -166.264 365.428 7
105 ORNL2 53 -73.6089 352.162 7
106 ORNL2 53 73.6089 352.162 7
107 ORNL2 53 166.264 365.428 7
108 ORNL2 53 201.103 278.553 7
109 ORNL2 53 323.857 159.305 7
110 ORNL2 53 339.453 102.2 7
111 ORNL2 53 -139.932 -142.007 4
112 ORNL2 53 -70.014 -198.718 4
113 ORNL2 53 0 -217.478 4
114 ORNL2 53 70.014 -198.718 4
115 ORNL2 53 139.932 -142.007 4
116 ORNL2 53 -214.903 -127.3 5
117 ORNL2 53 -133.044 -210.494 5
118 ORNL2 53 -48.6976 -259.192 5
119 ORNL2 53 48.6976 -259.192 5
120 ORNL2 53 133.044 -210.494 5
121 ORNL2 53 214.903 -127.3 5
122 ORNL2 53 -273.203 -102.2 6
123 ORNL2 53 -263.894 -172.91 6
124 ORNL2 53 -204.385 -193.709 6
125 ORNL2 53 -118.139 -282.073 6
126 ORNL2 53 0 -313.728 6
127 ORNL2 53 118.139 -282.073 6
128 ORNL2 53 204.385 -193.709 6
129 ORNL2 53 263.894 -172.91 6
130 ORNL2 53 273.203 -102.2 6
131 ORNL2 53 -339.453 -102.2 7
132 ORNL2 53 -323.857 -159.305 7
133 ORNL2 53 -201.103 -278.553 7
134 ORNL2 53 -166.264 -365.428 7
135 ORNL2 53 -73.6089 -352.162 7
136 ORNL2 53 73.6089 -352.162 7
137 ORNL2 53 166.264 -365.428 7
138 ORNL2 53 201.103 -278.553 7
139 ORNL2 53 323.857 -159.305 7
140 ORNL2 53 339.453 -102.2 7
Attachment 4:   
Drop attachments here...
ELOG V3.1.4-unknown