[807] | 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 | |
---|
| 35 | #include "G4MIRDLeftArmBone.hh" |
---|
| 36 | #include "globals.hh" |
---|
| 37 | #include "G4SDManager.hh" |
---|
| 38 | #include "G4VisAttributes.hh" |
---|
| 39 | #include "G4HumanPhantomMaterial.hh" |
---|
| 40 | #include "G4EllipticalTube.hh" |
---|
| 41 | #include "G4RotationMatrix.hh" |
---|
| 42 | #include "G4ThreeVector.hh" |
---|
| 43 | #include "G4VPhysicalVolume.hh" |
---|
| 44 | #include "G4PVPlacement.hh" |
---|
| 45 | #include "G4UnionSolid.hh" |
---|
| 46 | #include "G4EllipticalCone.hh" |
---|
| 47 | #include "G4HumanPhantomColour.hh" |
---|
| 48 | |
---|
| 49 | G4MIRDLeftArmBone::G4MIRDLeftArmBone() |
---|
| 50 | { |
---|
| 51 | } |
---|
| 52 | |
---|
| 53 | G4MIRDLeftArmBone::~G4MIRDLeftArmBone() |
---|
| 54 | { |
---|
| 55 | } |
---|
| 56 | |
---|
| 57 | G4VPhysicalVolume* G4MIRDLeftArmBone::Construct(const G4String& volumeName,G4VPhysicalVolume* mother, |
---|
| 58 | const G4String& colourName, G4bool wireFrame, G4bool sensitivity) |
---|
| 59 | { |
---|
| 60 | // Remind! the elliptical cone gives problems! Intersections of volumes, |
---|
| 61 | // wrong calculation of the volume! |
---|
| 62 | |
---|
| 63 | G4HumanPhantomMaterial* material = new G4HumanPhantomMaterial(); |
---|
| 64 | |
---|
| 65 | G4cout << "Construct " << volumeName <<G4endl; |
---|
| 66 | |
---|
| 67 | G4Material* skeleton = material -> GetMaterial("skeleton"); |
---|
| 68 | |
---|
| 69 | delete material; |
---|
| 70 | |
---|
| 71 | G4double dx = 1.4 * cm;//a |
---|
| 72 | G4double dy = 2.7 * cm;//b |
---|
| 73 | // G4double dz= 46. * cm;//z0 |
---|
| 74 | |
---|
| 75 | //G4EllipticalCone* arm = new G4EllipticalCone("OneLeftArmBone",dx/2.,dy/2.,dz, 34.5 *cm); |
---|
| 76 | G4EllipticalTube* leftArm = new G4EllipticalTube("OneLeftArmBone",dx,dy,34.5 *cm); |
---|
| 77 | |
---|
| 78 | G4LogicalVolume* logicLeftArmBone = new G4LogicalVolume(leftArm, |
---|
| 79 | skeleton, |
---|
| 80 | "logical" + volumeName, |
---|
| 81 | 0, 0,0); |
---|
| 82 | |
---|
| 83 | G4RotationMatrix* matrix = new G4RotationMatrix(); |
---|
| 84 | matrix -> rotateX(180. * degree); |
---|
| 85 | |
---|
| 86 | G4VPhysicalVolume* physLeftArmBone = new G4PVPlacement(matrix, |
---|
| 87 | G4ThreeVector(18.4 * cm, 0.0, -0.5*cm), |
---|
| 88 | //-x0 |
---|
| 89 | "physicalLeftArmBone", |
---|
| 90 | logicLeftArmBone, |
---|
| 91 | mother, |
---|
| 92 | false,0,true); |
---|
| 93 | |
---|
| 94 | |
---|
| 95 | // Sensitive Body Part |
---|
| 96 | if (sensitivity==true) |
---|
| 97 | { |
---|
| 98 | G4SDManager* SDman = G4SDManager::GetSDMpointer(); |
---|
| 99 | logicLeftArmBone->SetSensitiveDetector( SDman->FindSensitiveDetector("BodyPartSD") ); |
---|
| 100 | } |
---|
| 101 | |
---|
| 102 | |
---|
| 103 | // Visualization Attributes |
---|
| 104 | |
---|
| 105 | G4HumanPhantomColour* colourPointer = new G4HumanPhantomColour(); |
---|
| 106 | G4Colour colour = colourPointer -> GetColour(colourName); |
---|
| 107 | G4VisAttributes* LeftArmBoneVisAtt = new G4VisAttributes(colour); |
---|
| 108 | LeftArmBoneVisAtt->SetForceSolid(wireFrame); |
---|
| 109 | logicLeftArmBone->SetVisAttributes(LeftArmBoneVisAtt); |
---|
| 110 | |
---|
| 111 | G4cout << "LeftArmBone created !!!!!!" << G4endl; |
---|
| 112 | |
---|
| 113 | // Testing LeftArmBone Volume |
---|
| 114 | G4double LeftArmBoneVol = logicLeftArmBone->GetSolid()->GetCubicVolume(); |
---|
| 115 | G4cout << "Volume of LeftArmBone = " << LeftArmBoneVol/cm3 << " cm^3" << G4endl; |
---|
| 116 | |
---|
| 117 | // Testing LeftArmBone Material |
---|
| 118 | G4String LeftArmBoneMat = logicLeftArmBone->GetMaterial()->GetName(); |
---|
| 119 | G4cout << "Material of LeftArmBone = " << LeftArmBoneMat << G4endl; |
---|
| 120 | |
---|
| 121 | // Testing Density |
---|
| 122 | G4double LeftArmBoneDensity = logicLeftArmBone->GetMaterial()->GetDensity(); |
---|
| 123 | G4cout << "Density of Material = " << LeftArmBoneDensity*cm3/g << " g/cm^3" << G4endl; |
---|
| 124 | |
---|
| 125 | // Testing Mass |
---|
| 126 | G4double LeftArmBoneMass = (LeftArmBoneVol)*LeftArmBoneDensity; |
---|
| 127 | G4cout << "Mass of LeftArmBone = " << LeftArmBoneMass/gram << " g" << G4endl; |
---|
| 128 | |
---|
| 129 | return physLeftArmBone; |
---|
| 130 | } |
---|