Geant4 Cross Reference |
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 // Previous authors: G. Guerrieri, S. Guatelli and M. G. Pia, INFN Genova, Italy 27 // Authors (since 2007): S. Guatelli, University of Wollongong, Australia 28 // 29 // 30 #include "G4MIRDLeftClavicle.hh" 31 32 #include "globals.hh" 33 #include "G4SystemOfUnits.hh" 34 #include "G4SDManager.hh" 35 #include "G4Cons.hh" 36 #include "G4VisAttributes.hh" 37 #include "G4HumanPhantomMaterial.hh" 38 #include "G4EllipticalTube.hh" 39 #include "G4ThreeVector.hh" 40 #include "G4VPhysicalVolume.hh" 41 #include "G4RotationMatrix.hh" 42 #include "G4LogicalVolume.hh" 43 #include "G4PVPlacement.hh" 44 #include "G4SubtractionSolid.hh" 45 #include "G4HumanPhantomColour.hh" 46 #include "G4Box.hh" 47 #include "G4Torus.hh" 48 49 G4VPhysicalVolume* G4MIRDLeftClavicle::Construct(const G4String& volumeName, G4VPhysicalVolume* mother, 50 const G4String& colourName, G4bool wireFrame,G4bool) 51 { 52 53 G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl; 54 55 auto* material = new G4HumanPhantomMaterial(); 56 auto* skeleton = material -> GetMaterial("skeleton"); 57 58 G4double rMin = 0*cm; 59 G4double rMax = 0.7883*cm; 60 G4double rTor = 10*cm; 61 G4double pSPhi = 298.15*degree; 62 G4double pDPhi = 0.7*rad; 63 64 auto* clavicle = new G4Torus("Clavicle",rMin,rMax,rTor,pSPhi,pDPhi); 65 66 auto* logicLeftClavicle = new G4LogicalVolume(clavicle, 67 skeleton, 68 "logical" + volumeName, 69 nullptr, nullptr, nullptr); 70 71 72 G4VPhysicalVolume* physLeftClavicle = new G4PVPlacement(nullptr, 73 G4ThreeVector(0.*cm,2. *cm,33.25*cm), 74 "physicalLeftClavicle", 75 logicLeftClavicle, 76 mother, 77 false, 78 0, true); 79 80 // Visualization Attributes 81 //G4VisAttributes* LeftClavicleVisAtt = new G4VisAttributes(G4Colour(0.94,0.5,0.5)); 82 auto* colourPointer = new G4HumanPhantomColour(); 83 G4Colour colour = colourPointer -> GetColour(colourName); 84 auto* LeftClavicleVisAtt = new G4VisAttributes(colour); 85 LeftClavicleVisAtt->SetForceSolid(wireFrame); 86 logicLeftClavicle->SetVisAttributes(LeftClavicleVisAtt); 87 88 G4cout << "LeftClavicle created !!!!!!" << G4endl; 89 90 // Testing LeftClavicle Volume 91 G4double LeftClavicleVol = logicLeftClavicle->GetSolid()->GetCubicVolume(); 92 G4cout << "Volume of LeftClavicle = " << LeftClavicleVol/cm3 << " cm^3" << G4endl; 93 94 // Testing LeftClavicle Material 95 G4String LeftClavicleMat = logicLeftClavicle->GetMaterial()->GetName(); 96 G4cout << "Material of LeftClavicle = " << LeftClavicleMat << G4endl; 97 98 // Testing Density 99 G4double LeftClavicleDensity = logicLeftClavicle->GetMaterial()->GetDensity(); 100 G4cout << "Density of Material = " << LeftClavicleDensity*cm3/g << " g/cm^3" << G4endl; 101 102 // Testing Mass 103 G4double LeftClavicleMass = (LeftClavicleVol)*LeftClavicleDensity; 104 G4cout << "Mass of LeftClavicle = " << LeftClavicleMass/gram << " g" << G4endl; 105 106 107 return physLeftClavicle; 108 } 109