#include <G4ASCIITreeSceneHandler.hh>
Inheritance diagram for G4ASCIITreeSceneHandler:
Public Member Functions | |
G4ASCIITreeSceneHandler (G4VGraphicsSystem &system, const G4String &name) | |
virtual | ~G4ASCIITreeSceneHandler () |
virtual void | BeginModeling () |
virtual void | EndModeling () |
Protected Types | |
typedef std::set< G4LogicalVolume * >::iterator | LVSetIterator |
typedef G4PhysicalVolumeModel::G4PhysicalVolumeNodeID | PVNodeID |
typedef std::vector< PVNodeID > | PVPath |
typedef std::set< PVPath >::iterator | ReplicaSetIterator |
typedef std::set< PVPath >::reverse_iterator | ReplicaSetReverseIterator |
Protected Member Functions | |
virtual void | RequestPrimitives (const G4VSolid &) |
void | WriteHeader (std::ostream &) |
Protected Attributes | |
std::ostream * | fpOutFile |
std::ofstream | fOutFile |
std::ostringstream | fRestOfLine |
const G4VPhysicalVolume * | fpLastPV |
G4String | fLastPVName |
G4int | fLastCopyNo |
G4int | fLastNonSequentialCopyNo |
std::set< G4LogicalVolume * > | fLVSet |
std::set< PVPath > | fReplicaSet |
Definition at line 48 of file G4ASCIITreeSceneHandler.hh.
typedef std::set<G4LogicalVolume*>::iterator G4ASCIITreeSceneHandler::LVSetIterator [protected] |
Definition at line 73 of file G4ASCIITreeSceneHandler.hh.
typedef G4PhysicalVolumeModel::G4PhysicalVolumeNodeID G4ASCIITreeSceneHandler::PVNodeID [protected] |
Definition at line 74 of file G4ASCIITreeSceneHandler.hh.
typedef std::vector<PVNodeID> G4ASCIITreeSceneHandler::PVPath [protected] |
Definition at line 75 of file G4ASCIITreeSceneHandler.hh.
typedef std::set<PVPath>::iterator G4ASCIITreeSceneHandler::ReplicaSetIterator [protected] |
Definition at line 77 of file G4ASCIITreeSceneHandler.hh.
typedef std::set<PVPath>::reverse_iterator G4ASCIITreeSceneHandler::ReplicaSetReverseIterator [protected] |
Definition at line 78 of file G4ASCIITreeSceneHandler.hh.
G4ASCIITreeSceneHandler::G4ASCIITreeSceneHandler | ( | G4VGraphicsSystem & | system, | |
const G4String & | name | |||
) |
Definition at line 54 of file G4ASCIITreeSceneHandler.cc.
00055 : 00056 G4VTreeSceneHandler(system, name), 00057 fpOutFile(0), 00058 fpLastPV(0), 00059 fLastCopyNo(-99), 00060 fLastNonSequentialCopyNo(-99) 00061 {}
G4ASCIITreeSceneHandler::~G4ASCIITreeSceneHandler | ( | ) | [virtual] |
void G4ASCIITreeSceneHandler::BeginModeling | ( | ) | [virtual] |
Reimplemented from G4VTreeSceneHandler.
Definition at line 65 of file G4ASCIITreeSceneHandler.cc.
References G4VTreeSceneHandler::BeginModeling(), fOutFile, fpOutFile, G4cout, G4endl, G4VSceneHandler::GetGraphicsSystem(), and WriteHeader().
00065 { 00066 00067 G4VTreeSceneHandler::BeginModeling (); // To re-use "culling off" code. 00068 00069 const G4ASCIITree* pSystem = (G4ASCIITree*)GetGraphicsSystem(); 00070 const G4String& outFileName = pSystem -> GetOutFileName(); 00071 if (outFileName == "G4cout") { 00072 fpOutFile = &G4cout; 00073 } else { 00074 fOutFile.open (outFileName); 00075 fpOutFile = &fOutFile; 00076 } 00077 00078 G4cout << "G4ASCIITreeSceneHandler::BeginModeling: writing to "; 00079 if (outFileName == "G4cout") { 00080 G4cout << "G4 standard output (G4cout)"; 00081 } else { 00082 G4cout << "file \"" << outFileName << "\""; 00083 } 00084 G4cout << G4endl; 00085 00086 WriteHeader (G4cout); G4cout << G4endl; 00087 if (outFileName != "G4cout") { 00088 WriteHeader (fOutFile); fOutFile << std::endl; 00089 } 00090 }
void G4ASCIITreeSceneHandler::EndModeling | ( | ) | [virtual] |
Reimplemented from G4VTreeSceneHandler.
Definition at line 113 of file G4ASCIITreeSceneHandler.cc.
References G4VTreeSceneHandler::EndModeling(), fLastCopyNo, fLastNonSequentialCopyNo, fLastPVName, fLVSet, fOutFile, fpLastPV, fpOutFile, G4VSceneHandler::fpScene, fReplicaSet, fRestOfLine, G4BestUnit, G4cout, G4endl, G4VSceneHandler::GetGraphicsSystem(), G4PhysicalVolumeMassScene::GetMass(), G4TransportationManager::GetNavigatorForTracking(), G4Scene::GetRunDurationModelList(), G4TransportationManager::GetTransportationManager(), G4ASCIITree::GetVerbosity(), G4PhysicalVolumeMassScene::GetVolume(), G4Navigator::GetWorldVolume(), and G4PhysicalVolumeModel::UNLIMITED.
00113 { 00114 const G4ASCIITree* pSystem = (G4ASCIITree*) GetGraphicsSystem(); 00115 const G4int verbosity = pSystem->GetVerbosity(); 00116 const G4int detail = verbosity % 10; 00117 const G4String& outFileName = pSystem -> GetOutFileName(); 00118 00119 // Output left over copy number, if any... 00120 if (fLastCopyNo != fLastNonSequentialCopyNo) { 00121 if (fLastCopyNo == fLastNonSequentialCopyNo + 1) *fpOutFile << ','; 00122 else *fpOutFile << '-'; 00123 *fpOutFile << fLastCopyNo; 00124 } 00125 // Output outstanding rest of line, if any... 00126 if (!fRestOfLine.str().empty()) *fpOutFile << fRestOfLine.str(); 00127 fRestOfLine.str(""); 00128 fpLastPV = 0; 00129 fLastPVName.clear(); 00130 fLastCopyNo = -99; 00131 fLastNonSequentialCopyNo = -99; 00132 00133 // This detail to G4cout regardless of outFileName... 00134 if (detail >= 4) { 00135 G4cout << "Calculating mass(es)..." << G4endl; 00136 const std::vector<G4Scene::Model>& models = fpScene->GetRunDurationModelList(); 00137 std::vector<G4Scene::Model>::const_iterator i; 00138 for (i = models.begin(); i != models.end(); ++i) { 00139 G4PhysicalVolumeModel* pvModel = 00140 dynamic_cast<G4PhysicalVolumeModel*>(i->fpModel); 00141 if (pvModel) { 00142 if (pvModel->GetTopPhysicalVolume() == 00143 G4TransportationManager::GetTransportationManager() 00144 ->GetNavigatorForTracking()->GetWorldVolume()) { 00145 const G4ModelingParameters* tempMP = 00146 pvModel->GetModelingParameters(); 00147 G4ModelingParameters mp; // Default - no culling. 00148 pvModel->SetModelingParameters (&mp); 00149 G4PhysicalVolumeMassScene massScene(pvModel); 00150 pvModel->DescribeYourselfTo (massScene); 00151 G4double volume = massScene.GetVolume(); 00152 G4double mass = massScene.GetMass(); 00153 00154 G4cout << "Overall volume of \"" 00155 << pvModel->GetTopPhysicalVolume()->GetName() 00156 << "\":" 00157 << pvModel->GetTopPhysicalVolume()->GetCopyNo() 00158 << ", is " 00159 << G4BestUnit (volume, "Volume") 00160 << " and the daughter-included mass"; 00161 G4int requestedDepth = pvModel->GetRequestedDepth(); 00162 if (requestedDepth == G4PhysicalVolumeModel::UNLIMITED) { 00163 G4cout << " to unlimited depth"; 00164 } else { 00165 G4cout << ", ignoring daughters at depth " 00166 << requestedDepth 00167 << " and below,"; 00168 } 00169 G4cout << " is " << G4BestUnit (mass, "Mass") 00170 << G4endl; 00171 00172 pvModel->SetModelingParameters (tempMP); 00173 } 00174 } 00175 } 00176 } 00177 00178 if (outFileName != "G4cout") { 00179 fOutFile.close(); 00180 G4cout << "Output file \"" << outFileName << "\" closed." << G4endl; 00181 } 00182 fLVSet.clear(); 00183 fReplicaSet.clear(); 00184 G4cout << "G4ASCIITreeSceneHandler::EndModeling" << G4endl; 00185 G4VTreeSceneHandler::EndModeling (); // To re-use "culling off" code. 00186 }
void G4ASCIITreeSceneHandler::RequestPrimitives | ( | const G4VSolid & | ) | [protected, virtual] |
Reimplemented from G4VSceneHandler.
Definition at line 188 of file G4ASCIITreeSceneHandler.cc.
References G4PhysicalVolumeModel::CurtailDescent(), fLastCopyNo, fLastNonSequentialCopyNo, fLastPVName, fLVSet, fpLastPV, G4VSceneHandler::fpModel, fpOutFile, fReplicaSet, fRestOfLine, G4BestUnit, G4Exception(), G4VPhysicalVolume::GetCopyNo(), G4PhysicalVolumeModel::GetCurrentLV(), G4PhysicalVolumeModel::GetCurrentMaterial(), G4PhysicalVolumeModel::GetCurrentPV(), G4Material::GetDensity(), G4PhysicalVolumeModel::GetDrawnPVPath(), G4VSolid::GetEntityType(), G4VSceneHandler::GetGraphicsSystem(), G4VPhysicalVolume::GetLogicalVolume(), G4LogicalVolume::GetMass(), G4Material::GetName(), G4VSolid::GetName(), G4VReadOutGeometry::GetName(), G4VSensitiveDetector::GetName(), G4LogicalVolume::GetName(), G4VPhysicalVolume::GetName(), G4LogicalVolume::GetNoDaughters(), G4VPhysicalVolume::GetParameterisation(), G4VPhysicalVolume::GetReplicationData(), G4VSensitiveDetector::GetROgeometry(), G4LogicalVolume::GetSensitiveDetector(), G4ASCIITree::GetVerbosity(), G4VPhysicalVolume::IsParameterised(), G4VPhysicalVolume::IsReplicated(), and JustWarning.
00188 { 00189 00190 G4PhysicalVolumeModel* pPVModel = 00191 dynamic_cast<G4PhysicalVolumeModel*>(fpModel); 00192 if (!pPVModel) return; // Not from a G4PhysicalVolumeModel. 00193 00194 // This call comes from a G4PhysicalVolumeModel. drawnPVPath is 00195 // the path of the current drawn (non-culled) volume in terms of 00196 // drawn (non-culled) ancesters. Each node is identified by a 00197 // PVNodeID object, which is a physical volume and copy number. It 00198 // is a vector of PVNodeIDs corresponding to the geometry hierarchy 00199 // actually selected, i.e., not culled. 00200 // The following typedef's already set in header file... 00201 //typedef G4PhysicalVolumeModel::G4PhysicalVolumeNodeID PVNodeID; 00202 //typedef std::vector<PVNodeID> PVPath; 00203 const PVPath& drawnPVPath = pPVModel->GetDrawnPVPath(); 00204 //G4int currentDepth = pPVModel->GetCurrentDepth(); 00205 G4VPhysicalVolume* pCurrentPV = pPVModel->GetCurrentPV(); 00206 G4LogicalVolume* pCurrentLV = pPVModel->GetCurrentLV(); 00207 G4Material* pCurrentMaterial = pPVModel->GetCurrentMaterial(); 00208 00209 G4ASCIITree* pSystem = (G4ASCIITree*)GetGraphicsSystem(); 00210 G4int verbosity = pSystem->GetVerbosity(); 00211 G4int detail = verbosity % 10; 00212 00213 if (verbosity < 10 && pCurrentPV->IsReplicated()) { 00214 // See if this has been a replica found before with same mother LV... 00215 PVPath::const_reverse_iterator thisID = drawnPVPath.rbegin(); 00216 PVPath::const_reverse_iterator motherID = ++drawnPVPath.rbegin(); 00217 G4bool ignore = false; 00218 for (ReplicaSetIterator i = fReplicaSet.begin(); i != fReplicaSet.end(); 00219 ++i) { 00220 if (i->back().GetPhysicalVolume()->GetLogicalVolume() == 00221 thisID->GetPhysicalVolume()->GetLogicalVolume()) { 00222 // For each one previously found (if more than one, they must 00223 // have different mothers)... 00224 // To avoid compilation errors on VC++ .Net 7.1... 00225 // Previously: 00226 // PVNodeID previousMotherID = ++(i->rbegin()); 00227 // (Should that have been: PVNodeID::const_iterator previousMotherID?) 00228 // Replace 00229 // previousMotherID == i->rend() 00230 // by 00231 // i->size() <= 1 00232 // Replace 00233 // previousMotherID != i->rend() 00234 // by 00235 // i->size() > 1 00236 // Replace 00237 // previousMotherID-> 00238 // by 00239 // (*i)[i->size() - 2]. 00240 if (motherID == drawnPVPath.rend() && 00241 i->size() <= 1) 00242 ignore = true; // Both have no mother. 00243 if (motherID != drawnPVPath.rend() && 00244 i->size() > 1 && 00245 motherID->GetPhysicalVolume()->GetLogicalVolume() == 00246 (*i)[i->size() - 2].GetPhysicalVolume()->GetLogicalVolume()) 00247 ignore = true; // Same mother LV... 00248 } 00249 } 00250 if (ignore) { 00251 pPVModel->CurtailDescent(); 00252 return; 00253 } 00254 } 00255 00256 const G4String& currentPVName = pCurrentPV->GetName(); 00257 const G4int currentCopyNo = pCurrentPV->GetCopyNo(); 00258 00259 if (verbosity < 10 && 00260 currentPVName == fLastPVName && 00261 currentCopyNo != fLastCopyNo) { 00262 // For low verbosity, trap series of G4PVPlacements with the same name but 00263 // different copy number (replicas should have been taken out above)... 00264 // Check... 00265 if (pCurrentPV->IsReplicated()) { 00266 G4Exception("G4ASCIITreeSceneHandler::RequestPrimitives", 00267 "vistree0001", 00268 JustWarning, 00269 "Replica unexpected"); 00270 } 00271 // Check mothers are identical... 00272 else if (pCurrentLV == (fpLastPV? fpLastPV->GetLogicalVolume(): 0)) { 00273 if (currentCopyNo != fLastCopyNo + 1) { 00274 // Non-sequential copy number... 00275 *fpOutFile << ',' << currentCopyNo; 00276 fLastNonSequentialCopyNo = currentCopyNo; 00277 } 00278 fLastCopyNo = currentCopyNo; 00279 pPVModel->CurtailDescent(); 00280 return; 00281 } 00282 } 00283 fpLastPV = pCurrentPV; 00284 00285 // High verbosity or a new G4VPhysicalVolume... 00286 // Output copy number, if any, from previous invocation... 00287 if (fLastCopyNo != fLastNonSequentialCopyNo) { 00288 if (fLastCopyNo == fLastNonSequentialCopyNo + 1) *fpOutFile << ','; 00289 else *fpOutFile << '-'; 00290 *fpOutFile << fLastCopyNo; 00291 } 00292 // Output rest of line, if any, from previous invocation... 00293 if (!fRestOfLine.str().empty()) *fpOutFile << fRestOfLine.str(); 00294 fRestOfLine.str(""); 00295 fLastPVName = currentPVName; 00296 fLastCopyNo = currentCopyNo; 00297 fLastNonSequentialCopyNo = currentCopyNo; 00298 // Indent according to drawn path depth... 00299 for (size_t i = 0; i < drawnPVPath.size(); i++ ) *fpOutFile << " "; 00300 // Start next line... 00301 *fpOutFile << "\"" << currentPVName 00302 << "\":" << currentCopyNo; 00303 00304 if (pCurrentPV->IsReplicated()) { 00305 if (verbosity < 10) { 00306 // Add printing for replicas (when replicas are ignored)... 00307 EAxis axis; 00308 G4int nReplicas; 00309 G4double width; 00310 G4double offset; 00311 G4bool consuming; 00312 pCurrentPV->GetReplicationData(axis,nReplicas,width,offset,consuming); 00313 G4VPVParameterisation* pP = pCurrentPV->GetParameterisation(); 00314 if (pP) { 00315 if (detail < 3) { 00316 fReplicaSet.insert(drawnPVPath); 00317 if (nReplicas > 2) fRestOfLine << '-'; 00318 else fRestOfLine << ','; 00319 fRestOfLine << nReplicas - 1 00320 << " (" << nReplicas << " parametrised volumes)"; 00321 } 00322 } 00323 else { 00324 fReplicaSet.insert(drawnPVPath); 00325 if (nReplicas > 2) fRestOfLine << '-'; 00326 else fRestOfLine << ','; 00327 fRestOfLine << nReplicas - 1 00328 << " (" << nReplicas << " replicas)"; 00329 } 00330 } 00331 } else { 00332 if (fLVSet.find(pCurrentLV) != fLVSet.end()) { 00333 if (verbosity < 10) { 00334 // Add printing for repeated LV (if it has daughters)... 00335 if (pCurrentLV->GetNoDaughters()) fRestOfLine << " (repeated LV)"; 00336 // ...and curtail descent. 00337 pPVModel->CurtailDescent(); 00338 } 00339 } 00340 } 00341 00342 if (detail >= 1) { 00343 fRestOfLine << " / \"" 00344 << pCurrentLV->GetName() << "\""; 00345 G4VSensitiveDetector* sd = pCurrentLV->GetSensitiveDetector(); 00346 if (sd) { 00347 fRestOfLine << " (SD=\"" << sd->GetName() << "\""; 00348 G4VReadOutGeometry* roGeom = sd->GetROgeometry(); 00349 if (roGeom) { 00350 fRestOfLine << ",RO=\"" << roGeom->GetName() << "\""; 00351 } 00352 fRestOfLine << ")"; 00353 } 00354 } 00355 00356 if (detail >= 2) { 00357 fRestOfLine << " / \"" 00358 << solid.GetName() 00359 << "\"(" 00360 << solid.GetEntityType() << ")"; 00361 } 00362 00363 if (detail >= 3) { 00364 fRestOfLine << ", " 00365 << G4BestUnit(((G4VSolid&)solid).GetCubicVolume(),"Volume") 00366 << ", "; 00367 if (pCurrentMaterial) { 00368 fRestOfLine 00369 << G4BestUnit(pCurrentMaterial->GetDensity(), "Volumic Mass") 00370 << " (" << pCurrentMaterial->GetName() << ")"; 00371 } else { 00372 fRestOfLine << "(No material)"; 00373 } 00374 } 00375 00376 if (detail >= 5) { 00377 if (pCurrentMaterial) { 00378 G4Material* pMaterial = const_cast<G4Material*>(pCurrentMaterial); 00379 // ...and find daughter-subtracted mass... 00380 G4double daughter_subtracted_mass = pCurrentLV->GetMass 00381 (pCurrentPV->IsParameterised(), // Force if parametrised. 00382 false, // Do not propagate - calculate for this volume minus 00383 // volume of daughters. 00384 pMaterial); 00385 G4double daughter_subtracted_volume = 00386 daughter_subtracted_mass / pCurrentMaterial->GetDensity(); 00387 fRestOfLine << ", " 00388 << G4BestUnit(daughter_subtracted_volume,"Volume") 00389 << ", " 00390 << G4BestUnit(daughter_subtracted_mass,"Mass"); 00391 } 00392 } 00393 00394 if (fLVSet.find(pCurrentLV) == fLVSet.end()) { 00395 fLVSet.insert(pCurrentLV); // Record new logical volume. 00396 } 00397 00398 fRestOfLine << std::endl; 00399 00400 return; 00401 }
void G4ASCIITreeSceneHandler::WriteHeader | ( | std::ostream & | ) | [protected] |
Definition at line 92 of file G4ASCIITreeSceneHandler.cc.
References G4ASCIITreeMessenger::fVerbosityGuidance, G4VSceneHandler::GetGraphicsSystem(), and G4ASCIITree::GetVerbosity().
Referenced by BeginModeling().
00093 { 00094 const G4ASCIITree* pSystem = (G4ASCIITree*)GetGraphicsSystem(); 00095 const G4int verbosity = pSystem->GetVerbosity(); 00096 const G4int detail = verbosity % 10; 00097 os << "# Set verbosity with \"/vis/ASCIITree/verbose <verbosity>\":"; 00098 for (size_t i = 0; 00099 i < G4ASCIITreeMessenger::fVerbosityGuidance.size(); ++i) { 00100 os << "\n# " << G4ASCIITreeMessenger::fVerbosityGuidance[i]; 00101 } 00102 os << "\n# Now printing with verbosity " << verbosity; 00103 os << "\n# Format is: PV:n"; 00104 if (detail >= 1) os << " / LV (SD,RO)"; 00105 if (detail >= 2) os << " / Solid(type)"; 00106 if (detail >= 3) os << ", volume, density"; 00107 if (detail >= 5) os << ", daughter-subtracted volume and mass"; 00108 os << 00109 "\n# Abbreviations: PV = Physical Volume, LV = Logical Volume," 00110 "\n# SD = Sensitive Detector, RO = Read Out Geometry."; 00111 }
G4int G4ASCIITreeSceneHandler::fLastCopyNo [protected] |
Definition at line 69 of file G4ASCIITreeSceneHandler.hh.
Referenced by EndModeling(), and RequestPrimitives().
Definition at line 70 of file G4ASCIITreeSceneHandler.hh.
Referenced by EndModeling(), and RequestPrimitives().
G4String G4ASCIITreeSceneHandler::fLastPVName [protected] |
Definition at line 68 of file G4ASCIITreeSceneHandler.hh.
Referenced by EndModeling(), and RequestPrimitives().
std::set<G4LogicalVolume*> G4ASCIITreeSceneHandler::fLVSet [protected] |
Definition at line 72 of file G4ASCIITreeSceneHandler.hh.
Referenced by EndModeling(), and RequestPrimitives().
std::ofstream G4ASCIITreeSceneHandler::fOutFile [protected] |
Definition at line 63 of file G4ASCIITreeSceneHandler.hh.
Referenced by BeginModeling(), and EndModeling().
const G4VPhysicalVolume* G4ASCIITreeSceneHandler::fpLastPV [protected] |
Definition at line 67 of file G4ASCIITreeSceneHandler.hh.
Referenced by EndModeling(), and RequestPrimitives().
std::ostream* G4ASCIITreeSceneHandler::fpOutFile [protected] |
Definition at line 62 of file G4ASCIITreeSceneHandler.hh.
Referenced by BeginModeling(), EndModeling(), and RequestPrimitives().
std::set<PVPath> G4ASCIITreeSceneHandler::fReplicaSet [protected] |
Definition at line 76 of file G4ASCIITreeSceneHandler.hh.
Referenced by EndModeling(), and RequestPrimitives().
std::ostringstream G4ASCIITreeSceneHandler::fRestOfLine [protected] |
Definition at line 64 of file G4ASCIITreeSceneHandler.hh.
Referenced by EndModeling(), and RequestPrimitives().