G4PhysicalVolumeMassScene.hh

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 //
00027 // $Id$
00028 //
00029 // 
00030 // John Allison  12th September 2004
00031 
00032 // Class Description:
00033 //
00034 // Calculates the mass of a geometry tree taking into account daughters
00035 // up to the depth specified in the G4PhysicalVolumeModel.  Culling is
00036 // ignored so that all volumes are seen.
00037 //
00038 // Do not use this for a "parallel world" for which materials are not
00039 // defined.  Use only for the material world.
00040 //
00041 // The calculation is quite tricky, since it involves subtracting the
00042 // mass of that part of the mother that is occupied by each daughter and
00043 // then adding the mass of the daughter, and so on down the heirarchy.
00044 //
00045 // Usage for a given G4PhysicalVolumeModel* pvModel:
00046 //   G4PhysicalVolumeMassScene massScene(pvModel);
00047 //   pvModel->DescribeYourselfTo (massScene);
00048 //   G4double volume = massScene.GetVolume();
00049 //   G4double mass = massScene.GetMass();
00050 //   massScene.Reset();
00051 // See, for example, G4ASCIITreeSceneHandler::EndModeling().
00052 
00053 #ifndef G4PHYSICALVOLUMEMASSSCENE_HH
00054 #define G4PHYSICALVOLUMEMASSSCENE_HH
00055 
00056 #include "G4VGraphicsScene.hh"
00057 
00058 #include "G4Box.hh"
00059 #include "G4Cons.hh"
00060 #include "G4Tubs.hh"
00061 #include "G4Trd.hh"
00062 #include "G4Trap.hh"
00063 #include "G4Sphere.hh"
00064 #include "G4Para.hh"
00065 #include "G4Torus.hh"
00066 #include "G4Polycone.hh"
00067 #include "G4Polyhedra.hh"
00068 #include <deque>
00069 
00070 class G4VPhysicalVolume;
00071 class G4LogicalVolume;
00072 class G4PhysicalVolumeModel;
00073 class G4Material;
00074 
00075 class G4PhysicalVolumeMassScene: public G4VGraphicsScene {
00076 
00077 public:
00078   G4PhysicalVolumeMassScene (G4PhysicalVolumeModel*);
00079   virtual ~G4PhysicalVolumeMassScene ();
00080 
00081 public: // With description
00082 
00083   G4double GetVolume () const {return fVolume;}
00084   // Overall volume.
00085 
00086   G4double GetMass () const {return fMass;}
00087   // Mass of whole tree, i.e., accounting for all daughters.
00088 
00089   void Reset ();
00090   // Reset for subsequent re-use.
00091 
00092 public:
00093 
00094   // Force execution of AccrueMass for all solids...
00095   void PreAddSolid (const G4Transform3D&, const G4VisAttributes&) {}
00096   void PostAddSolid () {}
00097   void AddSolid (const G4Box& solid) {AccrueMass (solid);}
00098   void AddSolid (const G4Cons & solid) {AccrueMass (solid);}
00099   void AddSolid (const G4Tubs& solid) {AccrueMass (solid);}
00100   void AddSolid (const G4Trd& solid) {AccrueMass (solid);}
00101   void AddSolid (const G4Trap& solid) {AccrueMass (solid);}
00102   void AddSolid (const G4Sphere& solid) {AccrueMass (solid);}
00103   void AddSolid (const G4Para& solid) {AccrueMass (solid);}
00104   void AddSolid (const G4Torus& solid) {AccrueMass (solid);}
00105   void AddSolid (const G4Polycone& solid) {AccrueMass (solid);}
00106   void AddSolid (const G4Polyhedra& solid) {AccrueMass (solid);}
00107   void AddSolid (const G4VSolid& solid) {AccrueMass (solid);}
00108   void AddCompound (const G4VTrajectory&) {}
00109   void AddCompound (const G4VHit&) {}
00110   void AddCompound (const G4VDigi&) {}
00111   void AddCompound (const G4THitsMap<G4double>&) {}
00112 
00114   // Functions not used but required by the abstract interface.
00115 
00116   virtual void BeginPrimitives (const G4Transform3D&) {}
00117   virtual void EndPrimitives () {}
00118   virtual void BeginPrimitives2D (const G4Transform3D&) {}
00119   virtual void EndPrimitives2D () {}
00120   virtual void AddPrimitive (const G4Polyline&)   {}
00121   virtual void AddPrimitive (const G4Scale&)      {}
00122   virtual void AddPrimitive (const G4Text&)       {}
00123   virtual void AddPrimitive (const G4Circle&)     {}
00124   virtual void AddPrimitive (const G4Square&)     {}
00125   virtual void AddPrimitive (const G4Polymarker&) {}
00126   virtual void AddPrimitive (const G4Polyhedron&) {}
00127   virtual void AddPrimitive (const G4NURBS&)      {}
00128 
00129 private:
00130   void AccrueMass (const G4VSolid&);
00131   G4PhysicalVolumeModel* fpPVModel;
00132   G4double fVolume;
00133   G4double fMass;
00134   G4VPhysicalVolume* fpLastPV;
00135   G4int fPVPCount;
00136   G4int fLastDepth;
00137   G4double fLastDensity;
00138   std::deque<G4double> fDensityStack;
00139 };
00140 
00141 #endif

Generated on Mon May 27 17:49:19 2013 for Geant4 by  doxygen 1.4.7