Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
CellParameterisation Class Reference

#include <CellParameterisation.hh>

Inheritance diagram for CellParameterisation:
Collaboration diagram for CellParameterisation:

Public Member Functions

 CellParameterisation (G4Material *nucleus1, G4Material *cytoplasm1, G4Material *nucleus2, G4Material *cytoplasm2, G4Material *nucleus3, G4Material *cytoplasm3)
 
virtual ~CellParameterisation ()
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const
 
void ComputeDimensions (G4Box &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Tubs &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Trd &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Trap &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Cons &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Sphere &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Ellipsoid &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Orb &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Torus &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Para &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Polycone &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Polyhedra &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Hype &, const G4int, const G4VPhysicalVolume *) const
 
G4int GetNoBoxes ()
 
G4MaterialComputeMaterial (const G4int copyNo, G4VPhysicalVolume *physVol, const G4VTouchable *)
 
G4int GetPhantomTotalPixels ()
 
G4int GetNucleusTotalPixels ()
 
G4int GetCytoplasmTotalPixels ()
 
G4double GetPixelSizeX ()
 
G4double GetPixelSizeY ()
 
G4double GetPixelSizeZ ()
 
G4double GetCytoplasmMass ()
 
G4double GetNucleusMass ()
 
G4ThreeVector GetVoxelThreeVector (G4int i)
 
G4double GetMaterialVector (G4int i)
 
G4double GetMassVector (G4int i)
 
G4int GetTissueType (G4int i)
 
- Public Member Functions inherited from G4VPVParameterisation
 G4VPVParameterisation ()
 
virtual ~G4VPVParameterisation ()
 
virtual G4VSolidComputeSolid (const G4int, G4VPhysicalVolume *)
 
virtual G4bool IsNested () const
 
virtual G4VVolumeMaterialScannerGetMaterialScanner ()
 

Static Public Member Functions

static CellParameterisationInstance ()
 

Detailed Description

Definition at line 43 of file CellParameterisation.hh.

Constructor & Destructor Documentation

CellParameterisation::CellParameterisation ( G4Material nucleus1,
G4Material cytoplasm1,
G4Material nucleus2,
G4Material cytoplasm2,
G4Material nucleus3,
G4Material cytoplasm3 
)

Definition at line 46 of file CellParameterisation.cc.

50 {
51  fNucleusMaterial1 = nucleus1;
52  fCytoplasmMaterial1 = cytoplasm1;
53  fNucleusMaterial2 = nucleus2;
54  fCytoplasmMaterial2 = cytoplasm2;
55  fNucleusMaterial3 = nucleus3;
56  fCytoplasmMaterial3 = cytoplasm3;
57 
58  G4int ncols,nlines;
59  G4int shiftX, shiftY, shiftZ;
60  G4double x,y,z,mat,den,tmp,density;
61  G4double denCyto1, denCyto2, denCyto3, denNucl1, denNucl2, denNucl3;
62 
63  density=0;
64  denCyto1=0;
65  denCyto2=0;
66  denCyto3=0;
67  denNucl1=0;
68  denNucl2=0;
69  denNucl3=0;
70 
71  ncols=0; nlines=0;
72 
73  // READ PHANTOM
74 
75  fNucleusMass = 0;
76  fCytoplasmMass = 0;
77 
78  FILE *fMap;
79  fMap = fopen("phantom.dat","r");
80 
81  while (1)
82  {
83  if (nlines == 0)
84  {
85  ncols = fscanf(fMap,"%i %i %i",&fPhantomTotalPixels,&fNucleusTotalPixels,&fCytoplasmTotalPixels);
86  fMapCell = new G4ThreeVector[fPhantomTotalPixels];
87  fMaterial = new G4double[fPhantomTotalPixels];
88  fMass = new G4double[fPhantomTotalPixels];
89  fTissueType = new G4int[fPhantomTotalPixels];
90  }
91 
92  if (nlines == 1)
93  {
94  ncols = fscanf(fMap,"%lf %lf %lf",&fDimCellBoxX,&fDimCellBoxY,&fDimCellBoxZ);
95 
96  fDimCellBoxX=fDimCellBoxX*micrometer;
97  fDimCellBoxY=fDimCellBoxY*micrometer;
98  fDimCellBoxZ=fDimCellBoxZ*micrometer;
99  }
100 
101  if (nlines == 2) ncols = fscanf(fMap,"%i %i %i",&shiftX,&shiftY,&shiftZ); // VOXEL SHIFT IN Z ASSUMED TO BE NEGATIVE
102 
103  if (nlines == 3) ncols = fscanf(fMap,"%lf %lf %lf",&denCyto1, &denCyto2, &denCyto3);
104 
105  if (nlines == 4) ncols = fscanf(fMap,"%lf %lf %lf",&denNucl1, &denNucl2, &denNucl3);
106 
107  if (nlines > 4) ncols = fscanf(fMap,"%lf %lf %lf %lf %lf %lf",&x,&y,&z,&mat,&den,&tmp);
108 
109  if (ncols < 0) break;
110 
111  G4ThreeVector v(x+shiftX,y+shiftY,z-1500/(fDimCellBoxZ/micrometer)-shiftZ); // VOXEL SHIFT IN ORDER TO CENTER PHANTOM
112 
113  if (nlines>4)
114  {
115 
116  fMapCell[nlines-5]=v;
117  fMaterial[nlines-5]=mat;
118  fMass[nlines-5]=den;
119 
120  // fTissueType: 1 is Cytoplasm - 2 is Nucleus
121 
122  if( fMaterial[nlines-5] == 2 ) // fMaterial 2 is nucleus
123  {
124  if( fMass[nlines-5] == 1 )
125  {
126  fTissueType[nlines-5]=2;
127  }
128  if( fMass[nlines-5] == 2 )
129  {
130  fTissueType[nlines-5]=2;
131  }
132  if( fMass[nlines-5] == 3 )
133  {
134  fTissueType[nlines-5]=2;
135  }
136  }
137 
138  else if( fMaterial[nlines-5] == 1 ) // fMaterial 1 is cytoplasm
139  {
140  if( fMass[nlines-5] == 1 )
141  {
142  fTissueType[nlines-5]=1;
143  }
144  if( fMass[nlines-5] == 2 )
145  {
146  fTissueType[nlines-5]=2;
147  }
148  if( fMass[nlines-5] == 3 )
149  {
150  fTissueType[nlines-5]=1;
151  }
152  }
153 
154  //
155 
156 
157  if (std::abs(mat-2)<1.e-30) // NUCLEUS
158  {
159  if (std::abs(den-1)<1.e-30) density = denNucl1*(g/cm3);
160  if (std::abs(den-2)<1.e-30) density = denNucl2*(g/cm3);
161  if (std::abs(den-3)<1.e-30) density = denNucl3*(g/cm3);
162  fNucleusMass = fNucleusMass + density * fDimCellBoxX * fDimCellBoxY * fDimCellBoxZ ;
163  }
164 
165  if (std::abs(mat-1)<1.e-30) // CYTOPLASM
166  {
167  if (std::abs(den-1)<1e-30) density = denCyto1*(g/cm3);
168  if (std::abs(den-2)<1e-30) density = denCyto2*(g/cm3);
169  if (std::abs(den-3)<1e-30) density = denCyto3*(g/cm3);
170  fCytoplasmMass = fCytoplasmMass + density * fDimCellBoxX * fDimCellBoxY * fDimCellBoxZ ;
171  }
172 
173  }
174 
175  nlines++;
176  }
177  fclose(fMap);
178 
179  // NUCLEUS IN GREEN
180 
181  fNucleusAttributes1 = new G4VisAttributes;
182  fNucleusAttributes1->SetColour(G4Colour(0,.8,0));
183  fNucleusAttributes1->SetForceSolid(false);
184 
185  fNucleusAttributes2 = new G4VisAttributes;
186  fNucleusAttributes2->SetColour(G4Colour(0,.9,0));
187  fNucleusAttributes2->SetForceSolid(false);
188 
189  fNucleusAttributes3 = new G4VisAttributes;
190  fNucleusAttributes3->SetColour(G4Colour(0,1,0));
191  fNucleusAttributes3->SetForceSolid(false);
192 
193  // CYTOPLASM IN RED
194 
195  fCytoplasmAttributes1 = new G4VisAttributes;
196  fCytoplasmAttributes1->SetColour(G4Colour(1,0,0));
197  fCytoplasmAttributes1->SetForceSolid(false);
198 
199  fCytoplasmAttributes2 = new G4VisAttributes; // nucleoli in yellow
200  fCytoplasmAttributes2->SetColour(G4Colour(1.,1.,0));
201  fCytoplasmAttributes2->SetForceSolid(false);
202 
203  fCytoplasmAttributes3 = new G4VisAttributes;
204  fCytoplasmAttributes3->SetColour(G4Colour(1,0,0));
205  fCytoplasmAttributes3->SetForceSolid(false);
206 
207  //
208 
209  gInstance = this;
210 
211  }
void SetColour(const G4Colour &)
tuple x
Definition: test.py:50
int G4int
Definition: G4Types.hh:78
void SetForceSolid(G4bool=true)
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
static constexpr double cm3
Definition: G4SIunits.hh:121
tuple v
Definition: test.py:18
tuple z
Definition: test.py:28
double G4double
Definition: G4Types.hh:76
static constexpr double micrometer
Definition: G4SIunits.hh:100

Here is the call graph for this function:

CellParameterisation::~CellParameterisation ( )
virtual

Definition at line 215 of file CellParameterisation.cc.

216 {
217  delete[] fMapCell;
218  delete[] fMaterial;
219  delete[] fMass;
220  delete[] fTissueType;
221 }

Member Function Documentation

void CellParameterisation::ComputeDimensions ( G4Box ,
const G4int  ,
const G4VPhysicalVolume  
) const
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 240 of file CellParameterisation.cc.

241 {}
void CellParameterisation::ComputeDimensions ( G4Tubs ,
const G4int  ,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 61 of file CellParameterisation.hh.

63  {}
void CellParameterisation::ComputeDimensions ( G4Trd ,
const G4int  ,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 65 of file CellParameterisation.hh.

67  {}
void CellParameterisation::ComputeDimensions ( G4Trap ,
const G4int  ,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 69 of file CellParameterisation.hh.

71  {}
void CellParameterisation::ComputeDimensions ( G4Cons ,
const G4int  ,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 73 of file CellParameterisation.hh.

75  {}
void CellParameterisation::ComputeDimensions ( G4Sphere ,
const G4int  ,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 77 of file CellParameterisation.hh.

79  {}
void CellParameterisation::ComputeDimensions ( G4Ellipsoid ,
const G4int  ,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 81 of file CellParameterisation.hh.

83  {}
void CellParameterisation::ComputeDimensions ( G4Orb ,
const G4int  ,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 85 of file CellParameterisation.hh.

87  {}
void CellParameterisation::ComputeDimensions ( G4Torus ,
const G4int  ,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 89 of file CellParameterisation.hh.

91  {}
void CellParameterisation::ComputeDimensions ( G4Para ,
const G4int  ,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 93 of file CellParameterisation.hh.

95  {}
void CellParameterisation::ComputeDimensions ( G4Polycone ,
const G4int  ,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 97 of file CellParameterisation.hh.

99  {}
void CellParameterisation::ComputeDimensions ( G4Polyhedra ,
const G4int  ,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 101 of file CellParameterisation.hh.

103  {}
void CellParameterisation::ComputeDimensions ( G4Hype ,
const G4int  ,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 105 of file CellParameterisation.hh.

107  {}
G4Material * CellParameterisation::ComputeMaterial ( const G4int  copyNo,
G4VPhysicalVolume physVol,
const G4VTouchable  
)
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 246 of file CellParameterisation.cc.

249 {
250  if( fMaterial[copyNo] == 2 ) // fMaterial 2 is nucleus
251  {
252  if( fMass[copyNo] == 1 )
253  {
254  physVol->GetLogicalVolume()->SetVisAttributes( fNucleusAttributes1 );
255  return fNucleusMaterial1;
256  }
257  if( fMass[copyNo] == 2 )
258  {
259  physVol->GetLogicalVolume()->SetVisAttributes( fNucleusAttributes2 );
260  return fNucleusMaterial2;
261  }
262  if( fMass[copyNo] == 3 )
263  {
264  physVol->GetLogicalVolume()->SetVisAttributes( fNucleusAttributes3 );
265  return fNucleusMaterial3;
266  }
267  }
268 
269  else if( fMaterial[copyNo] == 1 ) // fMaterial 1 is cytoplasm
270  {
271  if( fMass[copyNo] == 1 )
272  {
273  physVol->GetLogicalVolume()->SetVisAttributes( fCytoplasmAttributes1 );
274  return fCytoplasmMaterial1;
275  }
276  if( fMass[copyNo] == 2 )
277  {
278  // nucleoli so taken as nucleus !
279  physVol->GetLogicalVolume()->SetVisAttributes( fCytoplasmAttributes2 );
280  return fCytoplasmMaterial2;
281  }
282  if( fMass[copyNo] == 3 )
283  {
284  physVol->GetLogicalVolume()->SetVisAttributes( fCytoplasmAttributes3 );
285  return fCytoplasmMaterial3;
286  }
287  }
288 
289  return physVol->GetLogicalVolume()->GetMaterial();
290 }
G4Material * GetMaterial() const
G4LogicalVolume * GetLogicalVolume() const
void SetVisAttributes(const G4VisAttributes *pVA)

Here is the call graph for this function:

void CellParameterisation::ComputeTransformation ( const G4int  copyNo,
G4VPhysicalVolume physVol 
) const
virtual

Implements G4VPVParameterisation.

Definition at line 226 of file CellParameterisation.cc.

227 {
229  origin(
230  fMapCell[copyNo].x()*fDimCellBoxX,
231  fMapCell[copyNo].y()*fDimCellBoxY,
232  fMapCell[copyNo].z()*fDimCellBoxZ);
233 
234  physVol->SetTranslation(origin);
235 }
tuple x
Definition: test.py:50
void SetTranslation(const G4ThreeVector &v)
tuple z
Definition: test.py:28

Here is the call graph for this function:

G4double CellParameterisation::GetCytoplasmMass ( )
inline

Definition at line 123 of file CellParameterisation.hh.

123 {return fCytoplasmMass;}
G4int CellParameterisation::GetCytoplasmTotalPixels ( )
inline

Definition at line 119 of file CellParameterisation.hh.

119 {return fCytoplasmTotalPixels;}
G4double CellParameterisation::GetMassVector ( G4int  i)
inline

Definition at line 128 of file CellParameterisation.hh.

128 {return fMass[i];}
G4double CellParameterisation::GetMaterialVector ( G4int  i)
inline

Definition at line 127 of file CellParameterisation.hh.

127 {return fMaterial[i];}
G4int CellParameterisation::GetNoBoxes ( )
inline

Definition at line 109 of file CellParameterisation.hh.

109 {return fPhantomTotalPixels;}
G4double CellParameterisation::GetNucleusMass ( )
inline

Definition at line 124 of file CellParameterisation.hh.

124 {return fNucleusMass;}
G4int CellParameterisation::GetNucleusTotalPixels ( )
inline

Definition at line 118 of file CellParameterisation.hh.

118 {return fNucleusTotalPixels;}
G4int CellParameterisation::GetPhantomTotalPixels ( )
inline

Definition at line 117 of file CellParameterisation.hh.

117 {return fPhantomTotalPixels;}
G4double CellParameterisation::GetPixelSizeX ( )
inline

Definition at line 120 of file CellParameterisation.hh.

120 {return fDimCellBoxX;}
G4double CellParameterisation::GetPixelSizeY ( )
inline

Definition at line 121 of file CellParameterisation.hh.

121 {return fDimCellBoxY;}
G4double CellParameterisation::GetPixelSizeZ ( )
inline

Definition at line 122 of file CellParameterisation.hh.

122 {return fDimCellBoxZ;}
G4int CellParameterisation::GetTissueType ( G4int  i)
inline

Definition at line 129 of file CellParameterisation.hh.

129 {return fTissueType[i];}
G4ThreeVector CellParameterisation::GetVoxelThreeVector ( G4int  i)
inline

Definition at line 126 of file CellParameterisation.hh.

126 {return fMapCell[i];}
static CellParameterisation* CellParameterisation::Instance ( void  )
inlinestatic

Definition at line 133 of file CellParameterisation.hh.

134  {
135  return gInstance;
136  }

The documentation for this class was generated from the following files: