Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4MIRDLeftOvary.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 "G4MIRDLeftOvary.hh"
35 
36 #include "globals.hh"
37 #include "G4SystemOfUnits.hh"
38 #include "G4SDManager.hh"
39 #include "G4VisAttributes.hh"
40 #include "G4Ellipsoid.hh"
41 #include "G4ThreeVector.hh"
42 #include "G4VPhysicalVolume.hh"
43 #include "G4RotationMatrix.hh"
44 #include "G4Material.hh"
45 #include "G4LogicalVolume.hh"
47 #include "G4VPhysicalVolume.hh"
48 #include "G4PVPlacement.hh"
49 #include "G4UnionSolid.hh"
50 #include "G4HumanPhantomColour.hh"
51 
53 {
54 }
55 
57 {
58 
59 }
60 
62  const G4String& colourName, G4bool wireFrame, G4bool sensitivity)
63 {
64  G4cout << "Construct "<< volumeName << G4endl;
65 
67  G4Material* soft = material -> GetMaterial("soft_tissue");
68  delete material;
69 
70  G4double ax= 1. *cm;
71  G4double by= 0.5*cm;
72  G4double cz= 2.*cm;
73 
74  G4Ellipsoid* OneOvary = new G4Ellipsoid("OneOvary",
75  ax, by, cz);
76 
77 
78  G4LogicalVolume* logicLeftOvary = new G4LogicalVolume(OneOvary,
79  soft,
80  "logical" + volumeName,
81  0, 0, 0);
82 
83  // Define rotation and position here!
84  G4VPhysicalVolume* physLeftOvary = new G4PVPlacement(0,
85  G4ThreeVector(-6. *cm,0.0*cm, -20*cm),
86  "physicalLeftOvary",
87  logicLeftOvary,
88  mother,
89  false,
90  0, true);
91 
92  // Sensitive Body Part
93  if (sensitivity==true)
94  {
96  logicLeftOvary->SetSensitiveDetector( SDman->FindSensitiveDetector("BodyPartSD") );
97  G4cout<< SDman->FindSensitiveDetector("BodyPartSD")->GetName()<< G4endl;
98  SDman->FindSensitiveDetector("BodyPartSD")->SetVerboseLevel(1);
99  }
100 
101  // Visualization Attributes
102  //G4VisAttributes* LeftOvaryVisAtt = new G4VisAttributes(G4Colour(0.85,0.44,0.84));
103  G4HumanPhantomColour* colourPointer = new G4HumanPhantomColour();
104  G4Colour colour = colourPointer -> GetColour(colourName);
105  delete colourPointer;
106 
107  G4VisAttributes* LeftOvaryVisAtt = new G4VisAttributes(colour);
108  LeftOvaryVisAtt->SetForceSolid(wireFrame);
109  logicLeftOvary->SetVisAttributes(LeftOvaryVisAtt);
110 
111  G4cout << "LeftOvary created !!!!!!" << G4endl;
112 
113  // Testing LeftOvary Volume
114  G4double LeftOvaryVol = logicLeftOvary->GetSolid()->GetCubicVolume();
115  G4cout << "Volume of LeftOvary = " << LeftOvaryVol/cm3 << " cm^3" << G4endl;
116 
117  // Testing LeftOvary Material
118  G4String LeftOvaryMat = logicLeftOvary->GetMaterial()->GetName();
119  G4cout << "Material of LeftOvary = " << LeftOvaryMat << G4endl;
120 
121  // Testing Density
122  G4double LeftOvaryDensity = logicLeftOvary->GetMaterial()->GetDensity();
123  G4cout << "Density of Material = " << LeftOvaryDensity*cm3/g << " g/cm^3" << G4endl;
124 
125  // Testing Mass
126  G4double LeftOvaryMass = (LeftOvaryVol)*LeftOvaryDensity;
127  G4cout << "Mass of LeftOvary = " << LeftOvaryMass/gram << " g" << G4endl;
128 
129  return physLeftOvary;
130 }