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 "G4MIRDMiddleLowerSpine.hh" 31 32 #include "globals.hh" 33 #include "G4SystemOfUnits.hh" 34 #include "G4SDManager.hh" 35 #include "G4VisAttributes.hh" 36 #include "G4HumanPhantomMaterial.hh" 37 #include "G4EllipticalTube.hh" 38 #include "G4RotationMatrix.hh" 39 #include "G4ThreeVector.hh" 40 #include "G4VPhysicalVolume.hh" 41 #include "G4PVPlacement.hh" 42 #include "G4UnionSolid.hh" 43 #include "G4HumanPhantomColour.hh" 44 45 G4VPhysicalVolume* G4MIRDMiddleLowerSpine::Construct(const G4String& volumeName, 46 G4VPhysicalVolume* mother, 47 const G4String& colourName 48 , G4bool wireFrame, G4bool ) 49 { 50 auto* material = new G4HumanPhantomMaterial(); 51 52 G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl; 53 54 auto* skeleton = material -> GetMaterial("skeleton"); 55 56 delete material; 57 58 G4double dx = 2. *cm; 59 G4double dy = 2.5 *cm; 60 G4double dz = 24. *cm; 61 62 auto* middleLowerSpine = new G4EllipticalTube("MiddleLowerSpine",dx, dy, dz); 63 64 auto* logicMiddleLowerSpine = new G4LogicalVolume( middleLowerSpine, skeleton, 65 "logical" + volumeName, 66 nullptr, nullptr, nullptr); 67 // Define rotation and position here! 68 G4VPhysicalVolume* physMiddleLowerSpine = new G4PVPlacement(nullptr,G4ThreeVector(0.0 *cm, 5.5 * cm,11. * cm), 69 "physicalMiddleLowerSpine", 70 logicMiddleLowerSpine, 71 mother, 72 false, 73 0, true); 74 75 76 // Visualization Attributes 77 // G4VisAttributes* MiddleLowerSpineVisAtt = new G4VisAttributes(G4Colour(0.46,0.53,0.6)); 78 79 auto* colourPointer = new G4HumanPhantomColour(); 80 G4Colour colour = colourPointer -> GetColour(colourName); 81 auto* MiddleLowerSpineVisAtt = new G4VisAttributes(colour); 82 MiddleLowerSpineVisAtt->SetForceSolid(wireFrame); 83 logicMiddleLowerSpine->SetVisAttributes(MiddleLowerSpineVisAtt); 84 85 G4cout << "MiddleLowerSpine created !!!!!!" << G4endl; 86 87 // Testing MiddleLowerSpine Volume 88 G4double MiddleLowerSpineVol = logicMiddleLowerSpine->GetSolid()->GetCubicVolume(); 89 G4cout << "Volume of MiddleLowerSpine = " << MiddleLowerSpineVol/cm3 << " cm^3" << G4endl; 90 91 // Testing MiddleLowerSpine Material 92 G4String MiddleLowerSpineMat = logicMiddleLowerSpine->GetMaterial()->GetName(); 93 G4cout << "Material of MiddleLowerSpine = " << MiddleLowerSpineMat << G4endl; 94 95 // Testing Density 96 G4double MiddleLowerSpineDensity = logicMiddleLowerSpine->GetMaterial()->GetDensity(); 97 G4cout << "Density of Material = " << MiddleLowerSpineDensity*cm3/g << " g/cm^3" << G4endl; 98 99 // Testing Mass 100 G4double MiddleLowerSpineMass = (MiddleLowerSpineVol)*MiddleLowerSpineDensity; 101 G4cout << "Mass of MiddleLowerSpine = " << MiddleLowerSpineMass/gram << " g" << G4endl; 102 103 return physMiddleLowerSpine; 104 } 105