GeFiCa
Germanium detector Field Calculator
Segmented.cc
Go to the documentation of this file.
1 #include "Units.h"
2 #include "Segmented.h"
3 using namespace GeFiCa;
4 
5 Segmented::Segmented(const char *name, const char *title) :
6  Detector(name, title), Radius(3.5*cm), BoreR(0.5*cm),
7  Nphi(6), Nz(3), SegmentId(1)
8 { Height=5*cm; Bias.push_back(2*kV); }
9 //______________________________________________________________________________
10 //
12 {
13  if (Radius<=0) {
14  Error("CheckConfigurations", "Radius==%.1f!", Radius);
15  abort();
16  }
17  if (BoreR<=0) {
18  Error("CheckConfigurations", "BoreR==%.1f!", BoreR);
19  abort();
20  }
21  if (BoreR>=Radius) {
22  Error("CheckConfigurations",
23  "BoreR (%.1f) >= Radius (%.1f)!", BoreR, Radius);
24  abort();
25  }
26  if (Nphi==0) {
27  Error("CheckConfigurations",
28  "Total number of segments in phi cannot be zero!");
29  abort();
30  }
31  if (Nz==0) {
32  Error("CheckConfigurations",
33  "Total number of segments in z cannot be zero!");
34  abort();
35  }
36  if (SegmentId>Nphi*Nz) {
37  Error("CheckConfigurations",
38  "SegmentId(%zu)>Nphi(%zu)*Nz(%zu)!",SegmentId, Nphi, Nz);
39  abort();
40  }
41  if (SegmentId==0) {
42  Info("CheckConfigurations", "SegmentId==0, "
43  "please use TrueCoaxial for the core electrode.");
44  abort();
45  }
46 }
47 //______________________________________________________________________________
48 //
49 #include <TLine.h>
50 #include <TEllipse.h>
51 void Segmented::Draw(Option_t* option)
52 {
53  TString pointOfView(option); pointOfView.ToLower();
54  if (pointOfView.Contains("top")) {
55  double x=Radius*cos(2*Pi/Nphi), y=Radius*sin(2*Pi/Nphi);
56 
57  TLine *l1 = new TLine(-Radius,0,Radius,0);
58  l1->SetLineColor(kBlack); l1->SetLineStyle(kDashed); l1->Draw();
59  TLine *l2 = new TLine(-x,-y,x,y);
60  l2->SetLineColor(kBlack); l2->SetLineStyle(kDashed); l2->Draw();
61  TLine *l3 = new TLine(-x,y,x,-y);
62  l3->SetLineColor(kBlack); l3->SetLineStyle(kDashed); l3->Draw();
63 
64  TEllipse *e1 = new TEllipse(0,0,BoreR,BoreR); e1->Draw();
65  TEllipse *e2 = new TEllipse(0,0,Radius,Radius); e2->SetFillStyle(0);
66  e2->Draw();
67  }
68 }
std::vector< double > Bias
bias on electrodes
Definition: Detector.h:35
size_t Nphi
total number of segments in phi
Definition: Segmented.h:16
static const double Pi
Definition: Units.h:29
double Radius
radius of the crystal
Definition: Segmented.h:14
static const double cm
centimeter
Definition: Units.h:12
static const double kV
kilo volt
Definition: Units.h:21
void Draw(Option_t *option="")
Definition: Segmented.cc:51
size_t SegmentId
segment Id in [0, Nphi*Nz]
Definition: Segmented.h:18
Detector & crystal properties.
Definition: Detector.h:32
double Height
height of crystal
Definition: Detector.h:13
double BoreR
radius of the bore hole
Definition: Segmented.h:15
A file defining commonly used units & constants.
void CheckConfigurations()
Definition: Segmented.cc:11
size_t Nz
total number of segments in z
Definition: Segmented.h:17
The only namespace in GeFiCa.
Definition: Detector.h:6
Segmented(const char *name="sip", const char *title="segmented detector")
Definition: Segmented.cc:5