Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4MIRDSpleen.cc
Go to the documentation of this file.
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 "G4MIRDSpleen.hh"
35 
36 #include "globals.hh"
37 #include "G4SystemOfUnits.hh"
38 #include "G4SDManager.hh"
39 #include "G4VisAttributes.hh"
40 #include "G4ThreeVector.hh"
41 #include "G4VPhysicalVolume.hh"
42 #include "G4RotationMatrix.hh"
43 #include "G4Material.hh"
44 #include "G4LogicalVolume.hh"
46 #include "G4VPhysicalVolume.hh"
47 #include "G4PVPlacement.hh"
48 #include "G4Ellipsoid.hh"
49 #include "G4HumanPhantomColour.hh"
50 
52 {
53 }
54 
56 {
57 
58 }
59 
61  const G4String& colourName, G4bool wireFrame, G4bool sensitivity)
62 {
63 
64  G4cout << "Construct "<< volumeName << G4endl;
66  G4Material* soft = material -> GetMaterial("soft_tissue");
67  delete material;
68 
69  G4double ax= 3.5 *cm;
70  G4double by= 2. *cm;
71  G4double cz= 6. * cm;
72 
73  G4Ellipsoid* spleen = new G4Ellipsoid("spleen", ax, by, cz);
74 
75 
76  G4LogicalVolume* logicSpleen = new G4LogicalVolume(spleen, soft,
77  "logical" + volumeName,
78  0, 0, 0);
79 
80  // Define rotation and position here!
81  G4VPhysicalVolume* physSpleen = new G4PVPlacement(0,
82  G4ThreeVector(11. *cm, 3. *cm, 2.*cm), // ztrans = half trunk lenght - z0
83  "physicalSpleen",
84  logicSpleen,
85  mother,
86  false,
87  0, true);
88 
89  // Sensitive Body Part
90  if (sensitivity==true)
91  {
93  logicSpleen->SetSensitiveDetector( SDman->FindSensitiveDetector("BodyPartSD") );
94  }
95 
96  // Visualization Attributes
97  // G4VisAttributes* SpleenVisAtt = new G4VisAttributes(G4Colour(0.41,0.41,0.41));
98  G4HumanPhantomColour* colourPointer = new G4HumanPhantomColour();
99  G4Colour colour = colourPointer -> GetColour(colourName);
100  G4VisAttributes* SpleenVisAtt = new G4VisAttributes(colour);
101  SpleenVisAtt->SetForceSolid(wireFrame);
102  logicSpleen->SetVisAttributes(SpleenVisAtt);
103 
104  G4cout << "Spleen created !!!!!!" << G4endl;
105 
106  // Testing Spleen Volume
107  G4double SpleenVol = logicSpleen->GetSolid()->GetCubicVolume();
108  G4cout << "Volume of Spleen = " << SpleenVol/cm3 << " cm^3" << G4endl;
109 
110  // Testing Spleen Material
111  G4String SpleenMat = logicSpleen->GetMaterial()->GetName();
112  G4cout << "Material of Spleen = " << SpleenMat << G4endl;
113 
114  // Testing Density
115  G4double SpleenDensity = logicSpleen->GetMaterial()->GetDensity();
116  G4cout << "Density of Material = " << SpleenDensity*cm3/g << " g/cm^3" << G4endl;
117 
118  // Testing Mass
119  G4double SpleenMass = (SpleenVol)*SpleenDensity;
120  G4cout << "Mass of Spleen = " << SpleenMass/gram << " g" << G4endl;
121 
122 
123 
124  return physSpleen;
125 }