source: trunk/examples/advanced/microbeam/src/MicrobeamCellParameterisation.cc @ 1282

Last change on this file since 1282 was 807, checked in by garnier, 16 years ago

update

  • Property svn:executable set to *
File size: 6.9 KB
Line 
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// $Id: MicrobeamCellParameterisation.cc,v 1.5 2006/06/29 16:05:23 gunter Exp $
28// -------------------------------------------------------------------
29
30#include "G4VPhysicalVolume.hh"
31#include "G4LogicalVolume.hh"
32
33#include "MicrobeamCellParameterisation.hh"
34
35//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
36
37MicrobeamCellParameterisation::MicrobeamCellParameterisation
38(G4int NoBoxes, G4float DimBoxX, G4float DimBoxY, G4float DimBoxZ, 
39 G4Material * nucleus1,  G4Material * cytoplasm1,
40 G4Material * nucleus2,  G4Material * cytoplasm2,
41 G4Material * nucleus3,  G4Material * cytoplasm3
42 )
43{
44   nucleusMaterial1 = nucleus1;
45   cytoplasmMaterial1 = cytoplasm1;
46   nucleusMaterial2 = nucleus2;
47   cytoplasmMaterial2 = cytoplasm2;
48   nucleusMaterial3 = nucleus3;
49   cytoplasmMaterial3 = cytoplasm3;
50
51   NoCellBoxes = NoBoxes;
52   DimCellBoxX = DimBoxX;
53   DimCellBoxY = DimBoxY;
54   DimCellBoxZ = DimBoxZ;
55 
56   mapCell  = new G4ThreeVector[NoCellBoxes];
57   material = new G4float[NoCellBoxes];
58   mass     = new G4float[NoCellBoxes];
59 
60   G4int ncols,nlines;
61   G4int shiftX, shiftY, shiftZ;
62   G4float x,y,z,mat,den,tmp,sizeZ; 
63   ncols=0; nlines=0;
64   
65   // READ PHANTOM
66   FILE *fMap;
67   fMap = fopen("phantom.dat","r");
68   while (1) 
69   { 
70      if (nlines >= 0  && nlines <=1 ) ncols = fscanf(fMap,"%f %f %f",&tmp,&tmp,&sizeZ);
71      if (nlines == 2) ncols = fscanf(fMap,"%i %i %i",&shiftX,&shiftY,&shiftZ); // VOXEL SHIFT IN Z ASSUMED TO BE NEGATIVE
72      if (nlines == 3) ncols = fscanf(fMap,"%f %f %f",&tmp,&tmp,&tmp);
73      if (nlines == 4) ncols = fscanf(fMap,"%f %f %f",&tmp,&tmp,&tmp);
74      if (nlines >  4) ncols = fscanf(fMap,"%f %f %f %f %f %f",&x,&y,&z,&mat,&den,&tmp);
75      if (ncols  <  0) break;
76
77      G4ThreeVector v(x+shiftX,y+shiftY,z-1500/sizeZ-shiftZ); // VOXEL SHIFT TO CENTER PHANTOM
78     
79      if (nlines>4) 
80      {
81          mapCell[nlines-5]=v; 
82          material[nlines-5]=mat;
83          mass[nlines-5]=den;
84      }   
85
86      nlines++;   
87   }
88   fclose(fMap);
89   
90  // NUCLEUS IN GREEN
91  nucleusAttributes1 = new G4VisAttributes;
92  nucleusAttributes1->SetColour(G4Colour(0,.8,0));
93  nucleusAttributes1->SetForceSolid(false);
94
95  nucleusAttributes2 = new G4VisAttributes;
96  nucleusAttributes2->SetColour(G4Colour(0,.9,0));
97  nucleusAttributes2->SetForceSolid(false);
98
99  nucleusAttributes3 = new G4VisAttributes;
100  nucleusAttributes3->SetColour(G4Colour(0,1,0));
101  nucleusAttributes3->SetForceSolid(false);
102
103// CYTOPLASM IN RED
104  cytoplasmAttributes1 = new G4VisAttributes;
105  cytoplasmAttributes1->SetColour(G4Colour(1,0,0));
106  cytoplasmAttributes1->SetForceSolid(false);
107
108  cytoplasmAttributes2 = new G4VisAttributes; // nucleoli in yellow
109  cytoplasmAttributes2->SetColour(G4Colour(1.,1.,0));
110  cytoplasmAttributes2->SetForceSolid(false);
111
112  cytoplasmAttributes3 = new G4VisAttributes;
113  cytoplasmAttributes3->SetColour(G4Colour(1,0,0));
114  cytoplasmAttributes3->SetForceSolid(false);
115}
116
117//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
118
119MicrobeamCellParameterisation::~MicrobeamCellParameterisation()
120{
121delete[] mapCell;
122delete[] material;
123delete[] mass;
124}
125
126//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
127
128void MicrobeamCellParameterisation::ComputeTransformation
129(const G4int copyNo, G4VPhysicalVolume* physVol) const
130{
131  G4ThreeVector
132    origin(
133      mapCell[copyNo].x()*DimCellBoxX*2,
134      mapCell[copyNo].y()*DimCellBoxY*2,
135      mapCell[copyNo].z()*DimCellBoxZ*2);
136
137  physVol->SetTranslation(origin);   
138}
139
140//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
141
142void MicrobeamCellParameterisation::ComputeDimensions
143(G4Box& /*trackerChamber*/, const G4int /*copyNo*/, const G4VPhysicalVolume*) const
144{
145}
146
147//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
148
149G4Material*
150MicrobeamCellParameterisation::ComputeMaterial(const G4int copyNo,
151                                               G4VPhysicalVolume* physVol,
152                                               const G4VTouchable*)
153{
154    if( material[copyNo] == 2 ) // material 2 is nucleus
155        {
156         if( mass[copyNo] == 1 ) 
157                {
158                physVol->SetName("physicalNucleus");
159                physVol->GetLogicalVolume()->SetVisAttributes( nucleusAttributes1 );
160                return nucleusMaterial1;
161                }
162         if( mass[copyNo] == 2 ) 
163                {
164                physVol->SetName("physicalNucleus");
165                physVol->GetLogicalVolume()->SetVisAttributes( nucleusAttributes2 );
166                return nucleusMaterial2;
167                }
168         if( mass[copyNo] == 3 ) 
169                {
170                physVol->SetName("physicalNucleus");
171                physVol->GetLogicalVolume()->SetVisAttributes( nucleusAttributes3 );
172                return nucleusMaterial3;
173                }
174        } 
175       
176    else if( material[copyNo] == 1 ) // material 1 is cytoplasm
177        {
178         if( mass[copyNo] == 1 ) 
179                {
180                physVol->SetName("physicalCytoplasm");
181                physVol->GetLogicalVolume()->SetVisAttributes( cytoplasmAttributes1 );
182                return cytoplasmMaterial1;
183                }
184         if( mass[copyNo] == 2 ) 
185                {
186                physVol->SetName("physicalNucleus"); // nucleoli
187                physVol->GetLogicalVolume()->SetVisAttributes( cytoplasmAttributes2 );
188                return cytoplasmMaterial2;
189                }
190         if( mass[copyNo] == 3 ) 
191                {
192                physVol->SetName("physicalCytoplasm");
193                physVol->GetLogicalVolume()->SetVisAttributes( cytoplasmAttributes3 );
194                return cytoplasmMaterial3;
195                }
196        } 
197
198    return physVol->GetLogicalVolume()->GetMaterial();
199}
200
201
Note: See TracBrowser for help on using the repository browser.