|
||||
File indexing completed on 2025-01-18 09:58:16
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 // 0027 // 0028 // 0029 //--------------------------------------------------------------- 0030 // 0031 // G4FastSimulationManager.hh 0032 // 0033 // Description: 0034 // Manages the Fast Simulation models attached to a envelope. 0035 // 0036 // History: 0037 // Oct 97: Verderi && MoraDeFreitas - First Implementation. 0038 // 0039 //--------------------------------------------------------------- 0040 0041 #ifndef G4FastSimulationManager_h 0042 #define G4FastSimulationManager_h 1 0043 0044 #include "G4FastSimulationVector.hh" 0045 #include "G4FastStep.hh" 0046 #include "G4FastTrack.hh" 0047 #include "G4LogicalVolume.hh" 0048 #include "G4ParticleDefinition.hh" 0049 #include "G4ParticleTable.hh" 0050 #include "G4Region.hh" 0051 #include "G4RotationMatrix.hh" 0052 #include "G4ThreeVector.hh" 0053 #include "G4Transform3D.hh" 0054 #include "G4VFastSimulationModel.hh" 0055 #include "G4VParticleChange.hh" 0056 #include "G4VPhysicalVolume.hh" 0057 #include "G4ios.hh" 0058 #include "globals.hh" 0059 0060 //------------------------------------------- 0061 // 0062 // G4FastSimulationManager class 0063 // 0064 //------------------------------------------- 0065 0066 // Class Description: 0067 // The G4VFastSimulationModel objects are attached to the envelope 0068 // through a G4FastSimulationManager. 0069 // This object will manage the list of models and will message them 0070 // at tracking time. 0071 // 0072 0073 class G4FastSimulationManager 0074 { 0075 public: // with description 0076 //------------------------ 0077 // Constructor/Destructor 0078 //------------------------ 0079 // Only one Constructor. By default the envelope can 0080 // be placed n-Times. 0081 // If the user is sure that it is placed just one time, 0082 // the IsUnique flag should be set TRUE to avoid the 0083 // G4AffineTransform re-calculations each time we reach 0084 // the envelope. 0085 0086 G4FastSimulationManager(G4Envelope* anEnvelope, G4bool IsUnique = FALSE); 0087 // This is the only constructor. In this constructor you specify 0088 // the envelope by giving the G4Region (typedef-ed as G4Envelope) 0089 // pointer. The G4FastSimulationManager object will bind itself to 0090 // this envelope and will notify this G4Region to become an envelope. 0091 // If you know that this region is used for only one logical volume, 0092 // you can turn the IsUnique boolean to "true" to allow some optimization. 0093 // 0094 // Note that if you choose to use the G4VFastSimulationModel(const G4String&, 0095 // G4Region*, G4bool) constructor for you model, the G4FastSimulationManager 0096 // will be constructed using the given G4Region* and G4bool values of the 0097 // model constructor. 0098 // 0099 0100 // Destructor 0101 ~G4FastSimulationManager(); 0102 0103 // Add a model to the Model List. 0104 void AddFastSimulationModel(G4VFastSimulationModel*); 0105 0106 // Remove a model from the Model List. 0107 void RemoveFastSimulationModel(G4VFastSimulationModel*); 0108 0109 // Activate a model in the Model List. 0110 G4bool ActivateFastSimulationModel(const G4String&); 0111 0112 // Inactivate a model in the Model List. 0113 G4bool InActivateFastSimulationModel(const G4String&); 0114 0115 // Methods for print/control commands 0116 void ListTitle() const; 0117 void ListModels() const; 0118 void ListModels(const G4ParticleDefinition*) const; 0119 void ListModels(const G4String& aName) const; 0120 const G4Envelope* GetEnvelope() const; 0121 0122 G4VFastSimulationModel* GetFastSimulationModel(const G4String& modelName, 0123 const G4VFastSimulationModel* previousFound, 0124 G4bool& foundPrevious) const; 0125 0126 const std::vector<G4VFastSimulationModel*>& GetFastSimulationModelList() const 0127 { 0128 return ModelList; 0129 } 0130 0131 void FlushModels(); 0132 0133 //---------------------------------------------- 0134 // Interface methods for the 0135 // G4FastSimulationManagerProcess process. 0136 //---------------------------------------------- 0137 // Trigger 0138 G4bool PostStepGetFastSimulationManagerTrigger(const G4Track&, const G4Navigator* a = nullptr); 0139 // DoIt 0140 G4VParticleChange* InvokePostStepDoIt(); 0141 0142 // AtRest methods: 0143 G4bool AtRestGetFastSimulationManagerTrigger(const G4Track&, const G4Navigator* a = nullptr); 0144 G4VParticleChange* InvokeAtRestDoIt(); 0145 0146 // For management 0147 G4bool operator==(const G4FastSimulationManager&) const; 0148 0149 private: 0150 // Private members : 0151 G4FastTrack fFastTrack; 0152 G4FastStep fFastStep; 0153 G4VFastSimulationModel* fTriggedFastSimulationModel{nullptr}; 0154 G4FastSimulationVector<G4VFastSimulationModel> ModelList; 0155 G4FastSimulationVector<G4VFastSimulationModel> fInactivatedModels; 0156 0157 G4ParticleDefinition* fLastCrossedParticle{nullptr}; 0158 G4FastSimulationVector<G4VFastSimulationModel> fApplicableModelList; 0159 0160 // -- *** depracating, to be dropped @ next major release: 0161 G4FastSimulationVector<G4Transform3D> GhostPlacements; 0162 }; 0163 0164 inline void G4FastSimulationManager::AddFastSimulationModel(G4VFastSimulationModel* fsm) 0165 { 0166 ModelList.push_back(fsm); 0167 // forces the fApplicableModelList to be rebuild 0168 fLastCrossedParticle = nullptr; 0169 } 0170 0171 inline void G4FastSimulationManager::RemoveFastSimulationModel(G4VFastSimulationModel* fsm) 0172 { 0173 if (ModelList.remove(fsm) == nullptr) fInactivatedModels.remove(fsm); 0174 // forces the fApplicableModelList to be rebuild 0175 fLastCrossedParticle = nullptr; 0176 } 0177 0178 inline G4bool G4FastSimulationManager::operator==(const G4FastSimulationManager& fsm) const 0179 { 0180 return this == &fsm; 0181 } 0182 0183 inline const G4Envelope* G4FastSimulationManager::GetEnvelope() const 0184 { 0185 return fFastTrack.GetEnvelope(); 0186 } 0187 0188 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |