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 // 31 32 #include "G4MIRDRightAdrenal.hh" 33 34 #include "globals.hh" 35 #include "G4SystemOfUnits.hh" 36 #include "G4SDManager.hh" 37 #include "G4VisAttributes.hh" 38 #include "G4HumanPhantomMaterial.hh" 39 #include "G4SDManager.hh" 40 #include "G4PVPlacement.hh" 41 #include "G4SubtractionSolid.hh" 42 #include "G4Ellipsoid.hh" 43 #include "G4ThreeVector.hh" 44 #include "G4VPhysicalVolume.hh" 45 #include "G4RotationMatrix.hh" 46 #include "G4Material.hh" 47 #include "G4EllipticalTube.hh" 48 #include "G4Box.hh" 49 #include "G4UnionSolid.hh" 50 #include "G4HumanPhantomColour.hh" 51 52 53 G4VPhysicalVolume* G4MIRDRightAdrenal::Construct(const G4String& volumeName,G4VPhysicalVolume* mother, 54 const G4String& colourName, G4bool wireFrame, G4bool) 55 { 56 G4cout<<"Construct "<<volumeName<<" with mother volume "<<mother->GetName()<<G4endl; 57 58 auto* material = new G4HumanPhantomMaterial(); 59 auto* soft = material -> GetMaterial("soft_tissue"); 60 delete material; 61 62 G4double ax= 1.5 *cm; //a 63 G4double by= 0.5 *cm; //b 64 G4double cz= 5.0 *cm; //c 65 66 auto* rightAdrenal = new G4Ellipsoid("OneRightAdrenal",ax, by, cz, 0. *cm, cz); 67 68 69 auto* logicRightAdrenal = new G4LogicalVolume(rightAdrenal, 70 soft, 71 "logical" + volumeName, 72 nullptr, nullptr, nullptr); 73 74 G4VPhysicalVolume* physRightAdrenal = new G4PVPlacement(nullptr, G4ThreeVector(-4.5*cm, // xo 75 6.5 *cm, //yo 76 3. *cm),//zo 77 "physicalRightAdrenal", logicRightAdrenal, 78 mother, 79 false, 80 0, true); 81 82 // Visualization Attributes 83 // G4VisAttributes* RightAdrenalVisAtt = new G4VisAttributes(G4Colour(0.72,0.52,0.04)); 84 auto* colourPointer = new G4HumanPhantomColour(); 85 G4Colour colour = colourPointer -> GetColour(colourName); 86 auto* RightAdrenalVisAtt = new G4VisAttributes(colour); 87 RightAdrenalVisAtt->SetForceSolid(wireFrame); 88 logicRightAdrenal->SetVisAttributes(RightAdrenalVisAtt); 89 90 G4cout << "Right RightAdrenal created !!!!!!" << G4endl; 91 92 // Testing RightAdrenal Volume 93 G4double RightAdrenalVol = logicRightAdrenal->GetSolid()->GetCubicVolume(); 94 G4cout << "Volume of RightAdrenal = " << RightAdrenalVol/cm3 << " cm^3" << G4endl; 95 96 // Testing RightAdrenal Material 97 G4String RightAdrenalMat = logicRightAdrenal->GetMaterial()->GetName(); 98 G4cout << "Material of RightAdrenal = " << RightAdrenalMat << G4endl; 99 100 // Testing Density 101 G4double RightAdrenalDensity = logicRightAdrenal->GetMaterial()->GetDensity(); 102 G4cout << "Density of Material = " << RightAdrenalDensity*cm3/g << " g/cm^3" << G4endl; 103 104 // Testing Mass 105 G4double RightAdrenalMass = (RightAdrenalVol)*RightAdrenalDensity; 106 G4cout << "Mass of RightAdrenal = " << RightAdrenalMass/gram << " g" << G4endl; 107 108 return physRightAdrenal; 109 } 110