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

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

update

File size: 5.5 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#include "G4MIRDRightLung.hh"
35#include "globals.hh"
36#include "G4SDManager.hh"
37#include "G4VisAttributes.hh"
38#include "G4Ellipsoid.hh"
39#include "G4ThreeVector.hh"
40#include "G4VPhysicalVolume.hh"
41#include "G4RotationMatrix.hh"
42#include "G4Material.hh"
43#include "G4LogicalVolume.hh"
44#include "G4HumanPhantomMaterial.hh"
45#include "G4VPhysicalVolume.hh"
46#include "G4PVPlacement.hh"
47#include "G4UnionSolid.hh"
48#include "G4SubtractionSolid.hh"
49#include "G4Box.hh"
50#include "G4HumanPhantomColour.hh"
51
52G4MIRDRightLung::G4MIRDRightLung()
53{
54}
55
56G4MIRDRightLung::~G4MIRDRightLung()
57{
58
59}
60
61G4VPhysicalVolume* G4MIRDRightLung::Construct(const G4String& volumeName,G4VPhysicalVolume* mother, 
62                                                   const G4String& colourName, G4bool wireFrame,G4bool sensitivity)
63{
64
65  G4cout << "Construct " << volumeName <<G4endl;
66 
67 G4HumanPhantomMaterial* material = new G4HumanPhantomMaterial();
68 G4Material* lung_material = material -> GetMaterial("lung_material");
69 delete material;
70
71 G4double ax = 5. *cm; //a
72 G4double by = 7.5 *cm; //b
73 G4double cz = 24.*cm; //c
74 G4double zcut1 = 0.0 *cm; 
75 G4double zcut2=24. *cm;
76 
77 G4Ellipsoid* oneLung = new G4Ellipsoid("OneLung",ax, by, cz, zcut1,zcut2);
78
79 ax= 5.*cm; 
80 by= 7.5*cm; 
81 cz= 24.*cm;
82
83
84 G4Ellipsoid* subtrLung = new G4Ellipsoid("subtrLung",ax, by, cz);
85
86 // y<0
87
88 G4double dx = 5.5* cm;
89 G4double dy = 8.5 * cm;
90 G4double dz = 24. * cm;
91
92 G4Box* box = new G4Box("Box", dx, dy, dz);
93 
94
95 G4SubtractionSolid* section = new G4SubtractionSolid("BoxSub", subtrLung, box, 0, G4ThreeVector(0.*cm, 8.5* cm, 0.*cm)); 
96 //G4SubtractionSolid* section2 = new G4SubtractionSolid("BoxSub2", subtrLung, box, 0, G4ThreeVector(0.*cm, -8.5* cm, 0.*cm));
97
98 G4SubtractionSolid* lung1 =  new G4SubtractionSolid("Lung1", oneLung,
99                                               section,
100                                               0, G4ThreeVector(6.*cm,0*cm,0.0*cm));
101
102
103 G4LogicalVolume* logicRightLung = new G4LogicalVolume(lung1,lung_material,
104                                                  "logical" + volumeName, 0, 0, 0); 
105 
106 G4RotationMatrix* matrix = new G4RotationMatrix();
107 matrix -> rotateZ(-360. * deg);
108  G4VPhysicalVolume* physRightLung = new G4PVPlacement(0,G4ThreeVector(-8.50 *cm, 0.0*cm, 8.5*cm),
109                                                  "physicalRightLung",                   
110                               logicRightLung,
111                               mother,
112                               false,
113                               0, true);
114
115
116  // Sensitive Body Part
117  if (sensitivity==true)
118  { 
119    G4SDManager* SDman = G4SDManager::GetSDMpointer();
120    logicRightLung->SetSensitiveDetector( SDman->FindSensitiveDetector("BodyPartSD") );
121 }
122
123  // Visualization Attributes
124  // G4VisAttributes* RightLungVisAtt = new G4VisAttributes(G4Colour(0.25,0.41,0.88));
125  G4HumanPhantomColour* colourPointer = new G4HumanPhantomColour();
126  G4Colour colour = colourPointer -> GetColour(colourName);
127  G4VisAttributes* RightLungVisAtt = new G4VisAttributes(colour);
128   RightLungVisAtt->SetForceSolid(wireFrame);
129  logicRightLung->SetVisAttributes(RightLungVisAtt); 
130
131  G4cout << "RightLung created !!!!!!" << G4endl;
132
133  // Testing RightLung Volume
134  G4double RightLungVol = logicRightLung->GetSolid()->GetCubicVolume();
135 
136 G4cout << "Volume of RightLung = " << (RightLungVol)/cm3 << " cm^3" << G4endl;
137 
138  // Testing RightLung Material
139  G4String RightLungMat = logicRightLung->GetMaterial()->GetName();
140  G4cout << "Material of RightLung = " << RightLungMat << G4endl;
141 
142  // Testing Density
143  G4double RightLungDensity = logicRightLung->GetMaterial()->GetDensity();
144  G4cout << "Density of Material = " << RightLungDensity*cm3/g << " g/cm^3" << G4endl;
145
146  // Testing Mass
147  G4double RightLungMass = (RightLungVol)*RightLungDensity;
148  G4cout << "Mass of RightLung = " << RightLungMass/gram << " g" << G4endl;
149 
150  return physRightLung;
151}
Note: See TracBrowser for help on using the repository browser.