source: trunk/environments/g4py/site-modules/primaries/MedicalBeam/pyMedicalBeam.cc @ 1358

Last change on this file since 1358 was 1337, checked in by garnier, 14 years ago

tag geant4.9.4 beta 1 + modifs locales

File size: 6.0 KB
Line 
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: pyMedicalBeam.cc,v 1.4 2006/06/29 15:30:14 gunter Exp $
27// $Name: geant4-09-04-beta-01 $
28// ====================================================================
29//   pyMedicalBeam.cc
30//
31//   [MedicalBeam]
32//   a site-module of Geant4Py
33//
34//   primary generator action for medical beam
35//
36//                                         2005 Q
37// ====================================================================
38#include <boost/python.hpp>
39#include "MedicalBeam.hh"
40#include "G4ParticleTable.hh"
41#include "G4RunManager.hh"
42
43using namespace boost::python;
44
45// ====================================================================
46// thin wrappers
47// ====================================================================
48namespace pyMedicalBeam {
49
50MedicalBeam* Construct()
51{
52  G4RunManager* runMgr= G4RunManager::GetRunManager();
53
54  MedicalBeam* medicalbeam= new MedicalBeam();
55  runMgr-> SetUserAction(medicalbeam);
56
57  return medicalbeam;
58}
59
60
61///////////////////////////////////////////////////////////////////
62void SetParticleByName(MedicalBeam* beam, const std::string& pname)
63///////////////////////////////////////////////////////////////////
64{
65  G4ParticleTable* particleTable= G4ParticleTable::GetParticleTable();
66  G4ParticleDefinition* pd= particleTable-> FindParticle(pname);
67  if (pd != 0) {
68    beam-> SetParticleDefinition(pd);
69  } else {
70    G4cout << "*** \"" << pname << "\" is not registered "
71           << "in available particle list" << G4endl;
72  }
73}
74
75////////////////////////////////////////////////
76std::string GetParticleByName(MedicalBeam* beam)
77////////////////////////////////////////////////
78{
79  const G4ParticleDefinition* pd= beam-> GetParticleDefinition();
80
81  if(pd==0) return std::string("None");
82  else return (pd-> GetParticleName()).c_str();
83}
84
85////////////////////////////////////////////////////////
86void f_SetFieldXY(MedicalBeam* beam, const list& listXY)
87////////////////////////////////////////////////////////
88{
89  G4double fx= extract<double>(listXY[0]);
90  G4double fy= extract<double>(listXY[1]);
91  beam-> SetFieldXY(fx, fy);
92}
93
94
95////////////////////////////////////
96list f_GetFieldXY(MedicalBeam* beam)
97////////////////////////////////////
98{
99  list listFieldXY;
100
101  listFieldXY.append(beam-> GetFieldX());
102  listFieldXY.append(beam-> GetFieldY());
103
104  return listFieldXY;
105}
106
107};
108
109using namespace pyMedicalBeam;
110
111// ====================================================================
112//   Expose to Python
113// ====================================================================
114
115BOOST_PYTHON_MODULE(MedicalBeam) {
116  class_<MedicalBeam, MedicalBeam*,
117    bases<G4VUserPrimaryGeneratorAction> >
118    ("MedicalBeam", "primary generator action with medical beam")
119    // ---
120    .add_property("particle", GetParticleByName, SetParticleByName)
121    .def("SetParticleByName", SetParticleByName)
122    .def("GetParticleByName", GetParticleByName)
123    // ---
124    .add_property("kineticE", &MedicalBeam::GetKineticE, 
125                              &MedicalBeam::SetKineticE)
126    .def("SetKineticE",       &MedicalBeam::SetKineticE)
127    .def("GetKineticE",       &MedicalBeam::GetKineticE)
128    // ---
129    .add_property("sourcePosition", &MedicalBeam::GetSourcePosition,
130                                    &MedicalBeam::SetSourcePosition)
131    .def("SetSourcePosition",       &MedicalBeam::SetSourcePosition)
132    .def("GetSourcePosition",       &MedicalBeam::GetSourcePosition)
133    // ---
134    .add_property("fieldShape", &MedicalBeam::GetFieldShape,
135                                &MedicalBeam::SetFieldShape)
136    .def("SetFieldShape",       &MedicalBeam::SetFieldShape)
137    .def("GetFieldShape",       &MedicalBeam::GetFieldShape)
138    // ---
139    .add_property("SSD", &MedicalBeam::GetSSD, &MedicalBeam::SetSSD)
140    .def("SetSSD",       &MedicalBeam::SetSSD)
141    .def("GetSSD",       &MedicalBeam::GetSSD)
142    // ----
143    .add_property("fieldXY", f_GetFieldXY, f_SetFieldXY)
144    .def("SetFieldXY",       f_SetFieldXY)
145    .def("GetFieldXY",       f_GetFieldXY)
146    .def("GetFieldX",        &MedicalBeam::GetFieldX)
147    .def("GetFieldY",        &MedicalBeam::GetFieldY)
148    // ---
149    .add_property("fieldR", &MedicalBeam::GetFieldR, &MedicalBeam::SetFieldR)
150    .def("SetFieldR",       &MedicalBeam::SetFieldR)
151    .def("GetFieldR",       &MedicalBeam::GetFieldR)
152    ;
153 
154  // enums...
155  enum_<MedicalBeam::FieldShape>("FieldShape")
156    .value("SQUARE", MedicalBeam::SQUARE)
157    .value("CIRCLE", MedicalBeam::CIRCLE)
158    ;
159 
160  // ---
161  def("Construct",  Construct,
162      return_value_policy<reference_existing_object>());
163}
164
Note: See TracBrowser for help on using the repository browser.