Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
pyG4ParticleTable.cc
Go to the documentation of this file.
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 // $Id: pyG4ParticleTable.cc 86749 2014-11-17 15:03:05Z gcosmo $
27 // ====================================================================
28 // pyG4ParticleTable.cc
29 //
30 // 2005 Q
31 // ====================================================================
32 #include <boost/python.hpp>
33 #include "G4Version.hh"
34 #include "G4ParticleTable.hh"
35 
36 using namespace boost::python;
37 
38 // ====================================================================
39 // thin wrappers
40 // ====================================================================
41 namespace pyG4ParticleTable {
42 
43 // contains...
46 
49 
50 // FindParticle...
53 
56 
59 
60 // FindAntiParticle...
63 
66 
69 
70 // DumpTable
71 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_DumpTable, DumpTable, 0, 1)
72 
73 
74 // --------------------------------------------------------------------
75 // GetParticleList (returning python list)
76 
77 list GetParticleList(G4ParticleTable* particleTable)
78 {
79  list particleList;
81  theParticleIterator= particleTable-> GetIterator();
82  theParticleIterator-> reset();
83  while( (*theParticleIterator)() ){
84  G4ParticleDefinition* particle= theParticleIterator-> value();
85  particleList.append(&particle);
86  }
87 
88  return particleList;
89 }
90 
91 }
92 
93 using namespace pyG4ParticleTable;
94 
95 // ====================================================================
96 // module definition
97 // ====================================================================
99 {
100  class_<G4ParticleTable, G4ParticleTable*, boost::noncopyable>
101  ("G4ParticleTable", "particle table", no_init)
102  // ---
103  .def("GetParticleTable", &G4ParticleTable::GetParticleTable,
104  return_value_policy<reference_existing_object>())
105  .staticmethod("GetParticleTable")
106  .def("contains", f1_contains)
107  .def("contains", f2_contains)
108  .def("entries", &G4ParticleTable::entries)
109  .def("size", &G4ParticleTable::size)
110  // ---
111  .def("GetParticle", &G4ParticleTable::GetParticle,
112  return_value_policy<reference_existing_object>())
113  .def("GetParticleName", &G4ParticleTable::GetParticleName,
114  return_value_policy<return_by_value>())
115  .def("FindParticle", f1_FindParticle,
116  return_value_policy<reference_existing_object>())
117  .def("FindParticle", f2_FindParticle,
118  return_value_policy<reference_existing_object>())
119  .def("FindParticle", f3_FindParticle,
120  return_value_policy<reference_existing_object>())
121  .def("FindAntiParticle", f1_FindAntiParticle,
122  return_value_policy<reference_existing_object>())
123  .def("FindAntiParticle", f2_FindAntiParticle,
124  return_value_policy<reference_existing_object>())
125  .def("FindAntiParticle", f3_FindAntiParticle,
126  return_value_policy<reference_existing_object>())
127  .def("DumpTable", &G4ParticleTable::DumpTable, f_DumpTable())
128  //.def("GetIonTable", &G4ParticleTable::GetIonTable,
129  //...)
130  //.def("GetShortLivedTable", &G4ParticleTable::GetShortLivedTable,
131  //...)
132  .def("SetVerboseLevel", &G4ParticleTable::SetVerboseLevel)
133  .def("GetVerboseLevel", &G4ParticleTable::GetVerboseLevel)
134  .def("SetReadiness", &G4ParticleTable::SetReadiness)
135  .def("GetReadiness", &G4ParticleTable::GetReadiness)
136  // ---
137  // additionals
138  .def("GetParticleList", GetParticleList)
139  ;
140 }
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
G4ParticleDefinition * FindAntiParticle(G4int PDGEncoding)
G4ParticleDefinition *(G4ParticleTable::* f2_FindAntiParticle)(const G4String &)
void DumpTable(const G4String &particle_name="ALL")
const G4String & GetParticleName(G4int index) const
BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_CreateTubeVolume, CreateTubeVolume, 4, 6) BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_CreateConeVolume
G4ParticleDefinition *(G4ParticleTable::* f1_FindParticle)(G4int)
void SetVerboseLevel(G4int value)
void SetReadiness(G4bool val=true)
G4bool GetReadiness() const
int G4int
Definition: G4Types.hh:78
G4bool(G4ParticleTable::* f2_contains)(const G4String &) const
static MCGIDI_particle * particleList
G4ParticleDefinition * GetParticle(G4int index) const
const XML_Char int const XML_Char * value
Definition: expat.h:331
bool G4bool
Definition: G4Types.hh:79
#define theParticleIterator
list GetParticleList(G4ParticleTable *particleTable)
G4bool contains(const G4ParticleDefinition *particle) const
static G4ParticleTable * GetParticleTable()
G4bool(G4ParticleTable::* f1_contains)(const G4ParticleDefinition *) const
G4ParticleDefinition *(G4ParticleTable::* f3_FindAntiParticle)(const G4ParticleDefinition *)
G4ParticleDefinition *(G4ParticleTable::* f2_FindParticle)(const G4String &)
G4int size() const
G4int GetVerboseLevel() const
G4int entries() const
G4ParticleDefinition *(G4ParticleTable::* f3_FindParticle)(const G4ParticleDefinition *)
G4ParticleDefinition *(G4ParticleTable::* f1_FindAntiParticle)(G4int)
void export_G4ParticleTable()