00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033 #include "G4LogicalVolumeModel.hh"
00034
00035 #include "G4VSolid.hh"
00036 #include "G4LogicalVolume.hh"
00037 #include "G4PVPlacement.hh"
00038 #include "G4ModelingParameters.hh"
00039 #include "G4VGraphicsScene.hh"
00040 #include "G4DrawVoxels.hh"
00041 #include "G4VSensitiveDetector.hh"
00042 #include "G4VReadOutGeometry.hh"
00043
00044 G4LogicalVolumeModel::G4LogicalVolumeModel
00045 (G4LogicalVolume* pLV,
00046 G4int soughtDepth,
00047 G4bool booleans,
00048 G4bool voxels,
00049 G4bool readout,
00050 const G4Transform3D& modelTransformation,
00051 const G4ModelingParameters* pMP):
00052
00053
00054
00055
00056
00057
00058 G4PhysicalVolumeModel
00059 (new G4PVPlacement (0,
00060 G4ThreeVector(),
00061 "PhysVol representation of LogVol " + pLV -> GetName (),
00062 pLV,
00063 0,
00064 false,
00065 0),
00066 soughtDepth,
00067 modelTransformation,
00068 pMP,
00069 true),
00070 fpLV (pLV),
00071 fBooleans (booleans),
00072 fVoxels (voxels),
00073 fReadout (readout)
00074 {
00075 fType = "G4LogicalVolumeModel";
00076 fGlobalTag = fpLV -> GetName ();
00077 fGlobalDescription = "G4LogicalVolumeModel " + fGlobalTag;
00078 }
00079
00080 G4LogicalVolumeModel::~G4LogicalVolumeModel () {}
00081
00082 void G4LogicalVolumeModel::DescribeYourselfTo
00083 (G4VGraphicsScene& sceneHandler) {
00084
00085
00086 const G4ModelingParameters* tmpMP = fpMP;
00087 G4ModelingParameters nonCulledMP;
00088 if (fpMP) nonCulledMP = *fpMP;
00089 nonCulledMP.SetCulling (false);
00090 fpMP = &nonCulledMP;
00091 G4PhysicalVolumeModel::DescribeYourselfTo (sceneHandler);
00092 fpMP = tmpMP;
00093
00094 if (fVoxels) {
00095 if (fpTopPV->GetLogicalVolume()->GetVoxelHeader()) {
00096
00097 G4DrawVoxels dv;
00098 G4PlacedPolyhedronList* pPPL =
00099 dv.CreatePlacedPolyhedra (fpTopPV -> GetLogicalVolume ());
00100 for (size_t i = 0; i < pPPL -> size (); i++) {
00101 const G4Transform3D& transform = (*pPPL)[i].GetTransform ();
00102 const G4Polyhedron& polyhedron = (*pPPL)[i].GetPolyhedron ();
00103 sceneHandler.BeginPrimitives (transform);
00104 sceneHandler.AddPrimitive (polyhedron);
00105 sceneHandler.EndPrimitives ();
00106 }
00107 delete pPPL;
00108 }
00109 }
00110
00111 if (fReadout) {
00112
00113 G4VSensitiveDetector* sd = fpLV->GetSensitiveDetector();
00114 if (sd) {
00115 G4VReadOutGeometry* roGeom = sd->GetROgeometry();
00116 if (roGeom) {
00117 G4VPhysicalVolume* roWorld = roGeom->GetROWorld();
00118 G4cout << "Readout geometry \"" << roGeom->GetName()
00119 << "\" with top physical volume \""
00120 << roWorld->GetName()
00121 << "\"" << G4endl;
00122 G4PhysicalVolumeModel pvModel(roWorld);
00123 pvModel.SetModelingParameters(fpMP);
00124 pvModel.DescribeYourselfTo(sceneHandler);
00125 }
00126 }
00127 }
00128 }
00129
00130
00131
00132 void G4LogicalVolumeModel::DescribeSolid
00133 (const G4Transform3D& theAT,
00134 G4VSolid* pSol,
00135 const G4VisAttributes* pVisAttribs,
00136 G4VGraphicsScene& sceneHandler) {
00137
00138 if (fBooleans) {
00139
00140 G4VSolid* pSol0 = pSol -> GetConstituentSolid (0);
00141 if (pSol0) {
00142 G4VSolid* pSol1 = pSol -> GetConstituentSolid (1);
00143 if (!pSol1) {
00144 G4Exception
00145 ("G4PhysicalVolumeModel::DescribeSolid",
00146 "modeling0001", FatalException,
00147 "2nd component solid in Boolean is missing.");
00148 }
00149
00150 G4VisAttributes constituentAttributes;
00151 constituentAttributes.SetForceWireframe(true);
00152 DescribeSolid (theAT, pSol0, &constituentAttributes, sceneHandler);
00153 DescribeSolid (theAT, pSol1, &constituentAttributes, sceneHandler);
00154 }
00155 }
00156
00157
00158 sceneHandler.PreAddSolid (theAT, *pVisAttribs);
00159 pSol -> DescribeYourselfTo (sceneHandler);
00160 sceneHandler.PostAddSolid ();
00161 }