Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/advanced/human_phantom/src/G4MIRDLeftLegBone.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 // Previous authors: G. Guerrieri, S. Guatelli and M. G. Pia, INFN Genova, Italy
 26 // Authors (since 2007): S. Guatelli, University of Wollongong, Australia
 27 // 
 28 
 29 #include "G4MIRDLeftLegBone.hh"
 30 
 31 #include "globals.hh"
 32 #include "G4SystemOfUnits.hh"
 33 #include "G4SDManager.hh"
 34 #include "G4VisAttributes.hh"
 35 #include "G4HumanPhantomMaterial.hh"
 36 #include "G4EllipticalTube.hh"
 37 #include "G4RotationMatrix.hh"
 38 #include "G4ThreeVector.hh"
 39 #include "G4VPhysicalVolume.hh"
 40 #include "G4PVPlacement.hh"
 41 #include "G4Cons.hh"
 42 #include "G4UnionSolid.hh"
 43 #include "G4HumanPhantomColour.hh"
 44 
 45 G4VPhysicalVolume* G4MIRDLeftLegBone::Construct(const G4String& volumeName,G4VPhysicalVolume* mother, 
 46             const G4String& colourName, G4bool wireFrame,G4bool)
 47 {
 48  
 49   auto* material = new G4HumanPhantomMaterial();
 50    
 51   G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl;
 52    
 53   auto* skeleton = material -> GetMaterial("skeleton");
 54   
 55   delete material;
 56  
 57   G4double dz = 79.8 * cm;
 58   G4double rmin1 = 0.0 * cm;
 59   G4double rmin2 = 0.0 * cm;
 60   G4double rmax1 = 1. * cm;
 61   G4double rmax2 = 3.5 * cm;
 62   G4double startphi = 0. * degree;
 63   G4double deltaphi = 360. * degree;
 64 
 65   auto* leg_bone = new G4Cons("OneLeftLegBone",  
 66         rmin1, rmax1, 
 67         rmin2, rmax2, dz/2., 
 68         startphi, deltaphi);
 69 
 70   auto* logicLeftLegBone = new G4LogicalVolume(leg_bone, skeleton,"logical" + volumeName,
 71                    nullptr, nullptr, nullptr);
 72 
 73 
 74   // Define rotation and position here!
 75   G4VPhysicalVolume* physLeftLegBone = new G4PVPlacement(nullptr,
 76                G4ThreeVector(0.0 * cm, 0.0, 0.1*cm),
 77                "physicalLeftLegBone",
 78                logicLeftLegBone,
 79                mother,
 80                false,
 81                0, true);
 82 
 83 
 84   // Visualization Attributes
 85   //G4VisAttributes* LeftLegBoneVisAtt = new G4VisAttributes(G4Colour(0.46,0.53,0.6));
 86   auto* colourPointer = new G4HumanPhantomColour();
 87   G4Colour colour = colourPointer -> GetColour(colourName);
 88   auto* LeftLegBoneVisAtt = new G4VisAttributes(colour);
 89  
 90   LeftLegBoneVisAtt->SetForceSolid(wireFrame);
 91   logicLeftLegBone->SetVisAttributes(LeftLegBoneVisAtt);
 92 
 93   G4cout << "LeftLegBone created !!!!!!" << G4endl;
 94 
 95   // Testing LeftLegBone Volume
 96   G4double LeftLegBoneVol = logicLeftLegBone->GetSolid()->GetCubicVolume();
 97   G4cout << "Volume of LeftLegBone = " << LeftLegBoneVol/cm3 << " cm^3" << G4endl;
 98   
 99   // Testing LeftLegBone Material
100   G4String LeftLegBoneMat = logicLeftLegBone->GetMaterial()->GetName();
101   G4cout << "Material of LeftLegBone = " << LeftLegBoneMat << G4endl;
102   
103   // Testing Density
104   G4double LeftLegBoneDensity = logicLeftLegBone->GetMaterial()->GetDensity();
105   G4cout << "Density of Material = " << LeftLegBoneDensity*cm3/g << " g/cm^3" << G4endl;
106 
107   // Testing Mass
108   G4double LeftLegBoneMass = (LeftLegBoneVol)*LeftLegBoneDensity;
109   G4cout << "Mass of LeftLegBone = " << LeftLegBoneMass/gram << " g" << G4endl;
110 
111   
112   return physLeftLegBone;
113 }
114