source: trunk/examples/advanced/human_phantom/src/G4MIRDLeftAdrenal.cc @ 1333

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

update

  • Property svn:executable set to *
File size: 4.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// 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 "G4MIRDLeftAdrenal.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"
53G4MIRDLeftAdrenal::G4MIRDLeftAdrenal()
54{
55}
56
57G4MIRDLeftAdrenal::~G4MIRDLeftAdrenal()
58{
59}
60
61G4VPhysicalVolume* G4MIRDLeftAdrenal::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= 1.5 *cm; //a
71 G4double by= 0.5 *cm; //b
72 G4double cz= 5.0 *cm; //c
73 
74 G4VSolid* leftAdrenal = new G4Ellipsoid("OneLeftAdrenal",ax, by, cz, 0. *cm, cz); 
75 
76 
77 G4LogicalVolume* logicLeftAdrenal = new G4LogicalVolume(leftAdrenal,
78                                                     soft,
79                                                     "logical" + volumeName,
80                                                     0, 0, 0);
81
82  G4VPhysicalVolume* physLeftAdrenal = new G4PVPlacement(0 ,G4ThreeVector(4.5*cm,  // xo
83                                                                     6.5 *cm, //yo
84                                                                     3. *cm),//zo
85                               "physicalLeftAdrenal", logicLeftAdrenal,
86                               mother,
87                               false,
88                               0, true);
89
90  // Sensitive Body Part
91  if (sensitivity==true)
92  { 
93    G4SDManager* SDman = G4SDManager::GetSDMpointer();
94    logicLeftAdrenal->SetSensitiveDetector( SDman->FindSensitiveDetector("BodyPartSD") );
95  }
96
97  // Visualization Attributes
98  //  G4VisAttributes* LeftAdrenalVisAtt = new G4VisAttributes(G4Colour(0.72,0.52,0.04));
99  G4HumanPhantomColour* colourPointer = new G4HumanPhantomColour();
100  G4Colour colour = colourPointer -> GetColour(colourName);
101  G4VisAttributes* LeftAdrenalVisAtt = new G4VisAttributes(colour);
102  LeftAdrenalVisAtt->SetForceSolid(wireFrame);
103  logicLeftAdrenal->SetVisAttributes(LeftAdrenalVisAtt);
104
105  G4cout << "Left LeftAdrenal created !!!!!!" << G4endl;
106
107  // Testing LeftAdrenal Volume
108  G4double LeftAdrenalVol = logicLeftAdrenal->GetSolid()->GetCubicVolume();
109  G4cout << "Volume of LeftAdrenal = " << LeftAdrenalVol/cm3 << " cm^3" << G4endl;
110 
111  // Testing LeftAdrenal Material
112  G4String LeftAdrenalMat = logicLeftAdrenal->GetMaterial()->GetName();
113  G4cout << "Material of LeftAdrenal = " << LeftAdrenalMat << G4endl;
114 
115  // Testing Density
116  G4double LeftAdrenalDensity = logicLeftAdrenal->GetMaterial()->GetDensity();
117  G4cout << "Density of Material = " << LeftAdrenalDensity*cm3/g << " g/cm^3" << G4endl;
118
119  // Testing Mass
120  G4double LeftAdrenalMass = (LeftAdrenalVol)*LeftAdrenalDensity;
121  G4cout << "Mass of LeftAdrenal = " << LeftAdrenalMass/gram << " g" << G4endl;
122
123 
124  return physLeftAdrenal;
125}
Note: See TracBrowser for help on using the repository browser.