source: trunk/examples/advanced/human_phantom/src/G4MIRDLeftKidney.cc @ 1317

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

update

File size: 5.2 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// Authors: S. Guatelli and M. G. Pia, INFN Genova, Italy
27//
28// Based on code developed by the undergraduate student G. Guerrieri
29// Note: this is a preliminary beta-version of the code; an improved
30// version will be distributed in the next Geant4 public release, compliant
31// with the design in a forthcoming publication, and subject to a
32// design and code review.
33//
34
35#include "G4MIRDLeftKidney.hh"
36
37#include "globals.hh"
38#include "G4SDManager.hh"
39#include "G4VisAttributes.hh"
40#include "G4HumanPhantomMaterial.hh"
41#include "G4SDManager.hh"
42#include "G4PVPlacement.hh"
43#include "G4SubtractionSolid.hh"
44#include "G4Ellipsoid.hh"
45#include "G4ThreeVector.hh"
46#include "G4VPhysicalVolume.hh"
47#include "G4RotationMatrix.hh"
48#include "G4Material.hh"
49#include "G4EllipticalTube.hh"
50#include "G4Box.hh"
51#include "G4UnionSolid.hh"
52#include "G4HumanPhantomColour.hh"
53G4MIRDLeftKidney::G4MIRDLeftKidney()
54{
55}
56
57G4MIRDLeftKidney::~G4MIRDLeftKidney()
58{
59}
60
61G4VPhysicalVolume* G4MIRDLeftKidney::Construct(const G4String& volumeName,G4VPhysicalVolume* mother,
62                                                    const G4String& colourName, G4bool wireFrame, G4bool sensitivity)
63{
64  G4cout << "Construct "<< volumeName << G4endl;
65 
66 G4HumanPhantomMaterial* material = new G4HumanPhantomMaterial();
67 G4Material* soft = material -> GetMaterial("soft_tissue");
68 delete material;
69 
70 G4double ax= 4.5 *cm; //a
71 G4double by= 1.5 *cm; //b
72 G4double cz= 5.5 *cm; //c
73 
74 G4VSolid* oneLeftKidney = new G4Ellipsoid("OneLeftKidney",ax, by, cz); 
75 
76 G4double xx = 6. * cm; 
77 G4double yy = 12.00*cm; 
78 G4double zz = 12.00*cm;
79 G4VSolid* subtrLeftKidney = new G4Box("SubtrLeftKidney",xx/2., yy/2., zz/2.);
80 
81 G4SubtractionSolid* leftKidney = new G4SubtractionSolid("LeftKidney",
82                                                     oneLeftKidney,
83                                                     subtrLeftKidney,
84                                                     0, 
85                                                     G4ThreeVector(-6. *cm, // x0
86                                                                   0.0 *cm,
87                                                                   0.0 * cm));
88
89  G4LogicalVolume* logicLeftKidney = new G4LogicalVolume(leftKidney,
90                                                     soft,
91                                                     "logical" + volumeName,
92                                                     0, 0, 0);
93
94  G4VPhysicalVolume* physLeftKidney = new G4PVPlacement(0 ,G4ThreeVector(6.*cm,  // xo
95                                                                     6. *cm, //yo
96                                                                     -2.50 *cm),//zo
97                               "physicalLeftKidney", logicLeftKidney,
98                               mother,
99                               false,
100                               0, true);
101
102  // Sensitive Body Part
103  if (sensitivity==true)
104  { 
105    G4SDManager* SDman = G4SDManager::GetSDMpointer();
106    logicLeftKidney->SetSensitiveDetector( SDman->FindSensitiveDetector("BodyPartSD") );
107  }
108
109  // Visualization Attributes
110  //  G4VisAttributes* LeftKidneyVisAtt = new G4VisAttributes(G4Colour(0.72,0.52,0.04));
111  G4HumanPhantomColour* colourPointer = new G4HumanPhantomColour();
112  G4Colour colour = colourPointer -> GetColour(colourName);
113  G4VisAttributes* LeftKidneyVisAtt = new G4VisAttributes(colour);
114  LeftKidneyVisAtt->SetForceSolid(wireFrame);
115  logicLeftKidney->SetVisAttributes(LeftKidneyVisAtt);
116
117  G4cout << "Left LeftKidney created !!!!!!" << G4endl;
118
119  // Testing LeftKidney Volume
120  G4double LeftKidneyVol = logicLeftKidney->GetSolid()->GetCubicVolume();
121  G4cout << "Volume of LeftKidney = " << LeftKidneyVol/cm3 << " cm^3" << G4endl;
122 
123  // Testing LeftKidney Material
124  G4String LeftKidneyMat = logicLeftKidney->GetMaterial()->GetName();
125  G4cout << "Material of LeftKidney = " << LeftKidneyMat << G4endl;
126 
127  // Testing Density
128  G4double LeftKidneyDensity = logicLeftKidney->GetMaterial()->GetDensity();
129  G4cout << "Density of Material = " << LeftKidneyDensity*cm3/g << " g/cm^3" << G4endl;
130
131  // Testing Mass
132  G4double LeftKidneyMass = (LeftKidneyVol)*LeftKidneyDensity;
133  G4cout << "Mass of LeftKidney = " << LeftKidneyMass/gram << " g" << G4endl;
134
135 
136  return physLeftKidney;
137}
Note: See TracBrowser for help on using the repository browser.