Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
NeuronLoadDataFile.cc
이 파일의 문서화 페이지로 가기
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 // This example is provided by the Geant4-DNA collaboration
27 // Any report or published results obtained using the Geant4-DNA software
28 // shall cite the following Geant4-DNA collaboration publication:
29 // Med. Phys. 37 (2010) 4692-4708
30 // and papers
31 // M. Batmunkh et al. J Radiat Res Appl Sci 8 (2015) 498-507
32 // O. Belov et al. Physica Medica 32 (2016) 1510-1520
33 // The Geant4-DNA web site is available at http://geant4-dna.org
34 //
35 // -------------------------------------------------------------------
36 // November 2016
37 // -------------------------------------------------------------------
38 //
39 // $Id:
40 //
43 
44 #include "NeuronLoadDataFile.hh"
45 //#include "NeuronLoadMessenger.hh"
46 #include "G4VPhysicalVolume.hh"
47 #include "G4LogicalVolume.hh"
48 #include "G4SystemOfUnits.hh"
49 #include "G4UnitsTable.hh"
50 #include "G4PhysicalConstants.hh"
51 #include "G4Colour.hh"
52 #include "G4VisAttributes.hh"
53 #include "G4RotationMatrix.hh"
54 #include "G4ios.hh"
55 #include <algorithm>
56 #include <fstream>
57 #include <iostream>
58 #include <limits>
59 #include <cmath>
60 #include <sstream>
61 #include <string>
62 #include <stdlib.h>
63 //define if the program is running with Geant4
64 #define GEANT4
65 #ifdef GEANT4
66 //Specific to Geant4, globals.hh is used for G4cout
67 #include "globals.hh"
68 #endif
69 #include "CommandLineParser.hh"
70 #include "G4UImanager.hh"
71 
72 using namespace std;
73 using namespace G4DNAPARSER;
74 
75 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
76 
78 {
79  //CommandLineParser* parser = CommandLineParser::GetParser();
80  Command* commandLine(0);
81 
82  // 1. Load single neuron morphology and obtain parameters.
83  // Default SWC file name of neuron
84  fNeuronFileNameSWC = G4String("GranuleCell-Nr2.CNG.swc");
85 
86  // 2. Load neural network and obtain parameters.
87  // Default prepared data filename of neural network with single/multi-layer.
88  // Small network of 10 pyramidal neurons with single layer
89  fNeuronFileNameDATA = G4String("NeuralNETWORK.dat");
90 
91  // Load/change SWC or DAT as "CommandLineParser" class
92  if((commandLine=CommandLineParser::GetParser()->GetCommandIfActive("-swc")))
93  {
94  fNeuronFileNameSWC = G4String(commandLine->GetOption());
95  SingleNeuronSWCfile(fNeuronFileNameSWC);
96  }
97  //if (CommandLineParser::GetParser()->GetCommandIfActive("-network"))
98  // {
99  // NeuralNetworkDATAfile(fNeuronFileNameDATA);
100  // }
101  if ((commandLine =CommandLineParser::GetParser()->
102  GetCommandIfActive("-network")))
103  {
104  fNeuronFileNameDATA = G4String(commandLine->GetOption());
105  NeuralNetworkDATAfile(fNeuronFileNameDATA);
106  }
107  else
108  {
109  SingleNeuronSWCfile(fNeuronFileNameSWC);
110  }
111 
112 }
113 
114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
115 
117 {
118 
119  // -----------
120  // 12 November 2012 - code created
121  // -------------------------------------------------------------------
122  // November 2012: First model of neuron[*] adapted into Geant4 microdosimetry
123  // from Claiborne`s database[**] by M. Batmunkh.
124  // February 2013: Loading SWC file from NeuronMorpho.Org[***]
125  // suggested by L. Bayarchimeg.
126  // [*] http://lt-jds.jinr.ru/record/62124/files/lrb_e_2012.pdf
127  // [**] http://www.utsa.edu/claibornelab/
128  // [***] http://neuromorpho.org
129  // -------------------------------------------------------------------
130 
131  G4String sLine = "";
132  std::ifstream infile;
133  infile.open(filename.c_str());
134  if (!infile)
135  {
136 #ifdef GEANT4
137  G4cout<<"\n NeuronLoadDataFile::SingleNeuronSWCfile >> datafile "
138  <<filename<<" not found !!!!"<<G4endl;
139  exit(0);
140 #endif
141  }
142  else
143  {
144 #ifdef G4VERBOSE
145  G4cout<<"\n NeuronLoadDataFile::SingleNeuronSWCfile >> opening filename: "
146  << "\n" <<'\t'<<'\t'<<'\t'<<'\t'<<'\t'<<"' "<<filename
147  << " ' \n"<< G4endl;
148 #endif
149 
150  G4int nrows,nlines;
151  nrows=0; nlines=0;
152  while (getline(infile, sLine))
153  {
154  fnNn = new G4int[nrows];
155  fpNn = new G4int[nrows];
156  fnNd = new G4int[nrows];
157  fpNd = new G4int[nrows];
158  fnNa = new G4int[nrows];
159  fpNa = new G4int[nrows];
160  nrows++;
161  }
162  infile.close();
163  //G4cout << " number of tracing points: "<< nrows <<G4endl;
164  infile.open(filename.c_str());
165 
166  fnbSomacomp = 0 ; // total number of compartment into Soma
167  fnbDendritecomp = 0 ; // total number of compartment into Dendrites
168  fnbAxoncomp = 0 ; // total number of compartment into Axon
169  fnbSpinecomp = 0 ; // total number of compartment into Spines
170  G4double TotVolSoma, TotVolDend, TotVolAxon, TotVolSpine;
171  TotVolSoma=TotVolDend=TotVolAxon=TotVolSpine=0.;
172  G4double TotSurfSoma, TotSurfDend, TotSurfAxon, TotSurfSpine;
173  TotSurfSoma=TotSurfDend=TotSurfAxon=TotSurfSpine=0.;
174  G4int nNcomp; // current index of neuronal compartment
175  G4int typeNcomp; // type of neuron structures: soma, axon, dendrite, etc.
176  G4double x,y,z; // cartesian coordinates of each compartment in micrometer
177  G4double radius; // radius of each compartment in micrometer
178  G4int pNcomp; // linked compartment, indicates branch points of dendrites
179  G4double minX,minY,minZ;
180  G4double maxX,maxY,maxZ;
181  G4double maxRad = -1e+09;
182  minX=minY=minZ=1e+09;
183  maxX=maxY=maxZ=-1e+09;
184  G4double density = 1.0 * (g/cm3) ; // water medium
185  G4double Piconst = (4.0/3.0)*pi ;
186 
187  fMassSomacomp = new G4double[nrows];
188  fMassSomaTot = 0.0 ;
189  fPosSomacomp = new G4ThreeVector[nrows];
190  fRadSomacomp = new G4double[nrows];
191  G4ThreeVector * PosDendcomp= new G4ThreeVector[nrows];
192  fRadDendcomp = new G4double[nrows];
193  fHeightDendcomp= new G4double[nrows];
194  fMassDendcomp = new G4double[nrows];
195  fMassDendTot = 0.0 ;
196  fDistADendSoma = new G4double[nrows];
197  fDistBDendSoma = new G4double[nrows];
198  fPosDendcomp = new G4ThreeVector[nrows];
199  fRotDendcomp = new G4RotationMatrix[nrows];
200  G4ThreeVector * PosAxoncomp= new G4ThreeVector[nrows];
201  fRadAxoncomp = new G4double[nrows];
202  fHeightAxoncomp= new G4double[nrows];
203  fMassAxoncomp = new G4double[nrows];
204  fMassAxonTot = 0.0 ;
205  fDistAxonsoma = new G4double[nrows];
206  fPosAxoncomp = new G4ThreeVector[nrows];
207  fRotAxoncomp = new G4RotationMatrix[nrows];
208  fMassSpinecomp = new G4double[nrows];
209  fMassSpineTot = 0.0 ;
210  fPosSpinecomp = new G4ThreeVector[nrows];
211  fRadSpinecomp = new G4double[nrows];
212  G4ThreeVector * PosNeuroncomp= new G4ThreeVector[nrows];
213  fRadNeuroncomp = new G4double[nrows];
214  fHeightNeuroncomp = new G4double[nrows];
215  fDistNeuronsoma = new G4double[nrows];
216  fPosNeuroncomp = new G4ThreeVector[nrows];
217  fRotNeuroncomp = new G4RotationMatrix[nrows];
218  fPosNeuroncomp = new G4ThreeVector[nrows];
219  fRadNeuroncomp = new G4double[nrows];
220  fTypeN = new G4int[nrows];
221 
222  // to read datafile containing numbers, alphabets and symbols..,
223  while (getline(infile, sLine))
224  {
225  std::istringstream form(sLine);
226  G4String token;
227  while (getline(form, token, ':'))
228  {
229  std::istringstream found(token);
230  while (found >> nNcomp >> typeNcomp >> x >> y >> z >> radius >> pNcomp)
231  {
232  // =======================================================================
233  // to find the largest and the smallest values of compartment positions
234  // for parameters of bounding slice, sphere medium and shift of neuron.
235  if (minX > x) minX = x;
236  if (minY > y) minY = y;
237  if (minZ > z) minZ = z;
238  if (maxX < x) maxX = x;
239  if (maxY < y) maxY = y;
240  if (maxZ < z) maxZ = z;
241  // max diameter of compartments
242  if (maxRad < radius) maxRad = radius;
243 
244  // =======================================================================
245  // Soma compartments represented as Sphere or Ellipsoid solid
246  if (typeNcomp == 1)
247  {
248  // Sphere volume and surface area
249  G4double VolSomacomp = Piconst*pow(radius*um,3.) ;
250  TotVolSoma = TotVolSoma + VolSomacomp;
251  G4double SurSomacomp = 3.*Piconst*pow(radius*um,2.) ;
252  TotSurfSoma = TotSurfSoma + SurSomacomp;
253  // OR
254  // Ellipsoid volume and Approximate formula of surface area
255  //G4double VolSomacomp = Piconst*(Ra*um)*(Rb*um)*(Rc*um);
256  //G4double SurSomacomp = 3.*Piconst*pow((pow(Ra,1.6075)*pow(Rb,1.6075)+
257  //pow(Ra,1.6075)*pow(Rc,1.6075)+pow(Rb,1.6075)*pow(Rc,1.6075))/3.,0.622084);
258  fMassSomacomp[fnbSomacomp] = density*VolSomacomp;
259  fMassSomaTot = fMassSomaTot + fMassSomacomp[fnbSomacomp];
260  G4ThreeVector vSoma (x ,y ,z);
261  fPosSomacomp [fnbSomacomp] = vSoma;
262  fRadSomacomp [fnbSomacomp]= radius;
263  // no rotate
264  // OR
265  // RotationMatrix for Ellipsoid solid
266  // ....
267  fnbSomacomp++ ;
268  }
269  // =======================================================================
270  // Apical and basal dendritic compartments represented as cylinderical solid
271  if (typeNcomp == 3 || typeNcomp == 4)
272  {
273  G4ThreeVector vDend (x ,y ,z);
274  // Position and Radius of compartments
275  PosDendcomp [fnbDendritecomp] = vDend;
276  fRadDendcomp [fnbDendritecomp]= radius;
277  fnNd[fnbDendritecomp]= nNcomp-(fnbSomacomp+fnbAxoncomp)-1;
278  fpNd[fnbDendritecomp]= pNcomp-(fnbSomacomp+fnbAxoncomp)-1;
279  // To join two tracing points along the dendritic branches.
280  // To calculate length, center and rotation angles of each cylinder
281 
282  // Center-position of each cylinder
283  G4double Dendxx= PosDendcomp[fnNd[fnbDendritecomp]].x()+
284  PosDendcomp[fpNd[fnbDendritecomp]].x();
285  G4double Dendyy= PosDendcomp[fnNd[fnbDendritecomp]].y()+
286  PosDendcomp[fpNd[fnbDendritecomp]].y();
287  G4double Dendzz= PosDendcomp[fnNd[fnbDendritecomp]].z()+
288  PosDendcomp[fpNd[fnbDendritecomp]].z();
289  G4ThreeVector translmDend = G4ThreeVector(Dendxx/2. ,
290  Dendyy/2. , Dendzz/2.) ;
291  fPosDendcomp [fnbDendritecomp] = translmDend;
292  // delta of position A and position B of cylinder
293  G4double Dendx, Dendy, Dendz;
294  //primary dendritic branch should be connect with Soma
295  if (fpNd[fnbDendritecomp] == -fnbSomacomp)
296  {
297  Dendx= PosDendcomp[fnNd[fnbDendritecomp]].x()-
298  (fPosSomacomp[0].x()+fRadSomacomp[0]);
299  Dendy= PosDendcomp[fnNd[fnbDendritecomp]].y()-
300  (fPosSomacomp[0].y()+fRadSomacomp[0]);
301  Dendz= PosDendcomp[fnNd[fnbDendritecomp]].z()-
302  (fPosSomacomp[0].z()+fRadSomacomp[0]);
303  }
304  else
305  {
306  Dendx= PosDendcomp[fnNd[fnbDendritecomp]].x()-
307  PosDendcomp[fpNd[fnbDendritecomp]].x();
308  Dendy= PosDendcomp[fnNd[fnbDendritecomp]].y()-
309  PosDendcomp[fpNd[fnbDendritecomp]].y();
310  Dendz= PosDendcomp[fnNd[fnbDendritecomp]].z()-
311  PosDendcomp[fpNd[fnbDendritecomp]].z();
312  }
313  G4double lengthDendcomp = std::sqrt(Dendx*Dendx+Dendy*Dendy+Dendz*Dendz);
314  // Height of compartment
315  fHeightDendcomp [fnbDendritecomp]= lengthDendcomp;
316 
317  // Distance from Soma
318  G4double DendDisx= fPosSomacomp[0].x()-
319  fPosDendcomp [fnbDendritecomp].x();
320  G4double DendDisy= fPosSomacomp[0].y()-
321  fPosDendcomp [fnbDendritecomp].y();
322  G4double DendDisz= fPosSomacomp[0].z()-
323  fPosDendcomp [fnbDendritecomp].z();
324  if (typeNcomp == 3) fDistADendSoma[fnbDendritecomp] =
325  std::sqrt(DendDisx*DendDisx + DendDisy*DendDisy + DendDisz*DendDisz);
326  if (typeNcomp == 4) fDistBDendSoma[fnbDendritecomp] =
327  std::sqrt(DendDisx*DendDisx + DendDisy*DendDisy + DendDisz*DendDisz);
328 
329  // Cylinder volume and surface area
330  G4double VolDendcomp = pi*pow(radius*um,2)*(lengthDendcomp*um);
331  TotVolDend = TotVolDend + VolDendcomp;
332  G4double SurDendcomp = 2.*pi*radius*um*(radius+lengthDendcomp)*um;
333  TotSurfDend = TotSurfDend + SurDendcomp;
334  fMassDendcomp[fnbDendritecomp] = density*VolDendcomp;
335  fMassDendTot = fMassDendTot + fMassDendcomp[fnbDendritecomp];
336 
337  Dendx=Dendx/lengthDendcomp;
338  Dendy=Dendy/lengthDendcomp;
339  Dendz=Dendz/lengthDendcomp;
340 
341  // Euler angles of each compartment
342  G4ThreeVector directionDend = G4ThreeVector(Dendx,Dendy,Dendz);
343  G4double theta_eulerDend = directionDend.theta();
344  G4double phi_eulerDend = directionDend.phi();
345  G4double psi_eulerDend = 0;
346 
347  //Rotation Matrix, Euler constructor build inverse matrix.
348  G4RotationMatrix rotmDendInv = G4RotationMatrix(
349  phi_eulerDend+pi/2,
350  theta_eulerDend,
351  psi_eulerDend);
352  G4RotationMatrix rotmDend = rotmDendInv.inverse();
353  /*
354  // To convert from Rotation Matrix after inverse to Euler angles
355  G4double cosX = std::sqrt (rotmDend.xx()*rotmDend.xx() +
356  rotmDend.yx()*rotmDend.yx()) ;
357  G4double euX, euY, euZ;
358  if (cosX > 16*FLT_EPSILON)
359  {
360  euX = std::atan2 (rotmDend.zy(),rotmDend.zz());
361  euY = std::atan2 (-rotmDend.zx(),cosX);
362  euZ = std::atan2 (rotmDend.yx(),rotmDend.xx());
363  }
364  else
365  {
366  euX = std::atan2 (-rotmDend.yz(),rotmDend.yy());
367  euY = std::atan2 (-rotmDend.zx(),cosX);
368  euZ = 0. ;
369  }
370  G4RotationMatrix * rot = new G4RotationMatrix();
371  rot->rotateX(euX);
372  rot->rotateY(euY);
373  rot->rotateZ(euZ);
374  */
375  fRotDendcomp [fnbDendritecomp]= rotmDend ;
376 
377  fnbDendritecomp++ ;
378  }
379 
380  // =======================================================================
381  // Axon compartments represented as cylinderical solid
382  if (typeNcomp == 2)
383  {
384  G4ThreeVector vAxon (x ,y ,z);
385  // Position and Radius of compartments
386  PosAxoncomp [fnbAxoncomp] = vAxon;
387  fRadAxoncomp [fnbAxoncomp]= radius;
388  fnNa[fnbAxoncomp]= nNcomp-(fnbSomacomp+fnbDendritecomp)-1;
389  fpNa[fnbAxoncomp]= pNcomp-(fnbSomacomp+fnbDendritecomp)-1;
390  // To join two tracing points in loaded SWC data file.
391  // To calculate length, center and rotation angles of each cylinder
392 
393  // Center-position of each cylinder
394  G4double Axonxx= PosAxoncomp[fnNa[fnbAxoncomp]].x()+
395  PosAxoncomp[fpNa[fnbAxoncomp]].x();
396  G4double Axonyy= PosAxoncomp[fnNa[fnbAxoncomp]].y()+
397  PosAxoncomp[fpNa[fnbAxoncomp]].y();
398  G4double Axonzz= PosAxoncomp[fnNa[fnbAxoncomp]].z()+
399  PosAxoncomp[fpNa[fnbAxoncomp]].z();
400  G4ThreeVector translmAxon = G4ThreeVector(Axonxx/2. ,
401  Axonyy/2. , Axonzz/2.) ;
402  fPosAxoncomp [fnbAxoncomp] = translmAxon;
403  // delta of position A and position B of cylinder
404  G4double Axonx, Axony, Axonz;
405  //primary axon point should be connect with Soma
406  if (fpNa[fnbAxoncomp] == -(fnbSomacomp+fnbDendritecomp))
407  {
408  Axonx= PosAxoncomp[fnNa[fnbAxoncomp]].x()-
409  (fPosSomacomp[0].x()+fRadSomacomp[0]);
410  Axony= PosAxoncomp[fnNa[fnbAxoncomp]].y()-
411  (fPosSomacomp[0].y()+fRadSomacomp[0]);
412  Axonz= PosAxoncomp[fnNa[fnbAxoncomp]].z()-
413  (fPosSomacomp[0].z()+fRadSomacomp[0]);
414  }
415  else
416  {
417  Axonx= PosAxoncomp[fnNa[fnbAxoncomp]].x()-
418  PosAxoncomp[fpNa[fnbAxoncomp]].x();
419  Axony= PosAxoncomp[fnNa[fnbAxoncomp]].y()-
420  PosAxoncomp[fpNa[fnbAxoncomp]].y();
421  Axonz= PosAxoncomp[fnNa[fnbAxoncomp]].z()-
422  PosAxoncomp[fpNa[fnbAxoncomp]].z();
423  }
424  G4double lengthAxoncomp = std::sqrt(Axonx*Axonx+Axony*Axony+Axonz*Axonz);
425  // Height of compartment
426  fHeightAxoncomp [fnbAxoncomp]= lengthAxoncomp;
427 
428  // Distance from Soma
429  G4double AxonDisx= fPosSomacomp[0].x()-
430  fPosAxoncomp [fnbAxoncomp].x();
431  G4double AxonDisy= fPosSomacomp[0].y()-
432  fPosAxoncomp [fnbAxoncomp].y();
433  G4double AxonDisz= fPosSomacomp[0].z()-
434  fPosAxoncomp [fnbAxoncomp].z();
435  fDistAxonsoma[fnbAxoncomp] = std::sqrt(AxonDisx*AxonDisx +
436  AxonDisy*AxonDisy + AxonDisz*AxonDisz);
437 
438  // Cylinder volume and surface area
439  G4double VolAxoncomp = pi*pow(radius*um,2)*(lengthAxoncomp*um);
440  TotVolAxon = TotVolAxon + VolAxoncomp;
441  G4double SurAxoncomp = 2.*pi*radius*um*(radius+lengthAxoncomp)*um;
442  TotSurfAxon = TotSurfAxon + SurAxoncomp;
443  fMassAxoncomp[fnbAxoncomp] = density*VolAxoncomp;
444  fMassAxonTot = fMassAxonTot + fMassAxoncomp[fnbAxoncomp];
445  Axonx=Axonx/lengthAxoncomp;
446  Axony=Axony/lengthAxoncomp;
447  Axonz=Axonz/lengthAxoncomp;
448 
449  // Euler angles of each compartment
450  G4ThreeVector directionAxon = G4ThreeVector(Axonx,Axony,Axonz);
451  G4double theta_eulerAxon = directionAxon.theta();
452  G4double phi_eulerAxon = directionAxon.phi();
453  G4double psi_eulerAxon = 0;
454 
455  //Rotation Matrix, Euler constructor build inverse matrix.
456  G4RotationMatrix rotmAxonInv = G4RotationMatrix(
457  phi_eulerAxon+pi/2,
458  theta_eulerAxon,
459  psi_eulerAxon);
460  G4RotationMatrix rotmAxon = rotmAxonInv.inverse();
461  fRotAxoncomp [fnbAxoncomp]= rotmAxon ;
462 
463  fnbAxoncomp++ ;
464  }
465  // =======================================================================
466  // checking additional types
467  if (typeNcomp != 1 && typeNcomp != 2 && typeNcomp != 3 && typeNcomp != 4)
468  {
469  G4cout << " Additional types:--> "<< typeNcomp <<G4endl;
470  }
471 
472  // If tracing points including spines, user can be define spine morphology
473  // including stubby, mushroom, thin, long thin, filopodia and
474  // branched with heads and necks!
475 
476  if (typeNcomp == 5)
477  {
478  // Sphere volume and surface area
479  G4double VolSpinecomp = Piconst*pow(radius*um,3.) ;
480  TotVolSpine = TotVolSpine + VolSpinecomp;
481  G4double SurSpinecomp = 3.*Piconst*pow(radius*um,2.) ;
482  TotSurfSpine = TotSurfSpine + SurSpinecomp;
483  fMassSpinecomp[fnbSpinecomp] = density*VolSpinecomp;
484  fMassSpineTot = fMassSpineTot + fMassSpinecomp[fnbSpinecomp];
485  // OR
486  // Ellipsoid volume and Approximate formula of surface area
487  // ...
488  G4ThreeVector vSpine (x ,y ,z);
489  fPosSpinecomp [fnbSpinecomp] = vSpine;
490  fRadSpinecomp [fnbSpinecomp] = radius;
491  // no rotate
492  // OR
493  // RotationMatrix for Ellipsoid solid
494  // ....
495  fnbSpinecomp++ ;
496  }
497 
498  // =======================================================================
499  // Neuron- all compartments allocate in ThreeVector
500  fTypeN[nlines] = typeNcomp;
501  G4ThreeVector vNeuron (x ,y ,z);
502  // Position and Radius of compartments
503  PosNeuroncomp [nlines] = vNeuron;
504  fRadNeuroncomp [nlines]= radius;
505  fnNn[nlines]= nNcomp-1;
506  fpNn[nlines]= pNcomp-1;
507  // To join two tracing points in loaded SWC data file.
508  // To calculate length, center and rotation angles of each cylinder
509 
510  // Center-position of each cylinder
511  G4double Neuronxx= PosNeuroncomp[fnNn[nlines]].x()+
512  PosNeuroncomp[fpNn[nlines]].x();
513  G4double Neuronyy= PosNeuroncomp[fnNn[nlines]].y()+
514  PosNeuroncomp[fpNn[nlines]].y();
515  G4double Neuronzz= PosNeuroncomp[fnNn[nlines]].z()+
516  PosNeuroncomp[fpNn[nlines]].z();
517  G4ThreeVector translmNeuron = G4ThreeVector(Neuronxx/2. ,
518  Neuronyy/2. , Neuronzz/2.) ;
519  fPosNeuroncomp [nlines] = translmNeuron;
520  // delta of position A and position B of cylinder
521  G4double Neuronx, Neurony, Neuronz;
522  //primary point
523  if (fpNn[nlines] == -2)
524  {
525  Neuronx= PosNeuroncomp[fnNn[nlines]].x()-
526  fPosNeuroncomp[0].x();
527  Neurony= PosNeuroncomp[fnNn[nlines]].y()-
528  fPosNeuroncomp[0].y();
529  Neuronz= PosNeuroncomp[fnNn[nlines]].z()-
530  fPosNeuroncomp[0].z();
531  }
532  else
533  {
534  Neuronx= PosNeuroncomp[fnNn[nlines]].x()-
535  PosNeuroncomp[fpNn[nlines]].x();
536  Neurony= PosNeuroncomp[fnNn[nlines]].y()-
537  PosNeuroncomp[fpNn[nlines]].y();
538  Neuronz= PosNeuroncomp[fnNn[nlines]].z()-
539  PosNeuroncomp[fpNn[nlines]].z();
540  }
541  G4double lengthNeuroncomp = std::sqrt(Neuronx*Neuronx+
542  Neurony*Neurony+Neuronz*Neuronz);
543  // Height of compartment
544  fHeightNeuroncomp [nlines]= lengthNeuroncomp;
545  // Distance from Soma
546  G4double NeuronDisx= fPosNeuroncomp[0].x()-
547  fPosNeuroncomp [nlines].x();
548  G4double NeuronDisy= fPosNeuroncomp[0].y()-
549  fPosNeuroncomp [nlines].y();
550  G4double NeuronDisz= fPosNeuroncomp[0].z()-
551  fPosNeuroncomp [nlines].z();
552  fDistNeuronsoma[nlines] = std::sqrt(NeuronDisx*NeuronDisx +
553  NeuronDisy*NeuronDisy + NeuronDisz*NeuronDisz);
554  /*
555  // Cylinder volume and surface area
556  G4double VolNeuroncomp = pi*pow(radius*um,2)*(lengthNeuroncomp*um);
557  fTotVolNeuron = fTotVolNeuron + VolNeuroncomp;
558  G4double SurNeuroncomp = 2.*pi*radius*um*
559  (radius+lengthNeuroncomp)*um;
560  fTotSurfNeuron = TotSurfNeuron + SurNeuroncomp;
561  fMassNeuroncomp[nlines] = density*VolNeuroncomp;
562  MassNeuronTot = MassNeuronTot + fMassNeuroncomp[nlines]; */
563  Neuronx=Neuronx/lengthNeuroncomp;
564  Neurony=Neurony/lengthNeuroncomp;
565  Neuronz=Neuronz/lengthNeuroncomp;
566 
567  // Euler angles of each compartment
568  G4ThreeVector directionNeuron = G4ThreeVector(Neuronx,Neurony,Neuronz);
569  G4double theta_eulerNeuron = directionNeuron.theta();
570  G4double phi_eulerNeuron = directionNeuron.phi();
571  G4double psi_eulerNeuron = 0;
572 
573  //Rotation Matrix, Euler constructor build inverse matrix.
574  G4RotationMatrix rotmNeuronInv = G4RotationMatrix(
575  phi_eulerNeuron+pi/2,
576  theta_eulerNeuron,
577  psi_eulerNeuron);
578  G4RotationMatrix rotmNeuron = rotmNeuronInv.inverse();
579  fRotNeuroncomp [nlines]= rotmNeuron ;
580 
581  nlines++;
582  }
583  }
584  }
585  infile.close();
586  // =======================================================================
587 
588  fnbNeuroncomp = nlines ;
589  G4cout << " Total number of compartments into Neuron : "
590  << fnbNeuroncomp<<G4endl;
591  G4cout << "\n"<<G4endl;
592 
593  // to calculate SHIFT value for neuron translation
594  fshiftX = (minX + maxX)/2. ;
595  fshiftY = (minY + maxY)/2. ;
596  fshiftZ = (minZ + maxZ)/2. ;
597 
598  // width, height, depth of bounding slice volume
599  //maxRad = 0.0 ;
600  fwidthB = std::fabs(minX - maxX) + maxRad;
601  fheightB = std::fabs(minY - maxY) + maxRad;
602  fdepthB = std::fabs(minZ - maxZ) + maxRad;
603 
604  // diagonal length of bounding slice, that give diameter of sphere
605  // for particle direction and fluence!
606  fdiagnlLength = std::sqrt(fwidthB*fwidthB + fheightB*fheightB
607  + fdepthB*fdepthB);
608 
609  fTotVolNeuron = TotVolSoma+TotVolDend+TotVolAxon;
610  fTotSurfNeuron = TotSurfSoma+TotSurfDend+TotSurfAxon;
611  fTotMassNeuron = fMassSomaTot+fMassDendTot+fMassAxonTot;
612 
613  fTotVolSlice = fwidthB*um*fheightB*um*fdepthB*um;
614  fTotSurfSlice = 2*(fwidthB*um*fheightB*um+fheightB*um*fdepthB*um+
615  fwidthB*um*fdepthB*um);
616  fTotMassSlice = 1.0 * (g/cm3) *fTotVolSlice;
617 
618  fTotVolMedium = Piconst*pow(fdiagnlLength*um/2.,3.) ;
619  fTotSurfMedium = 3.*Piconst*pow(fdiagnlLength*um/2.,2);
620  fTotMassMedium = 1.0 * (g/cm3) *fTotVolMedium;
621 
622  // Soma in Violet with opacity
623  fSomaColour = new G4VisAttributes;
624  fSomaColour->SetColour(G4Colour(G4Colour(0.85,0.44,0.84))); // ,1.0
625  fSomaColour->SetForceSolid(true); // true
626  fSomaColour->SetVisibility(true);
627 
628  // Dendrites in Dark-Blue
629  fDendColour = new G4VisAttributes;
630  fDendColour->SetColour(G4Colour(G4Colour(0.0, 0.0, 0.5)));
631  fDendColour->SetForceSolid(true);
632  //fDendColour->SetVisibility(true);
633 
634  // Axon in Maroon
635  fAxonColour = new G4VisAttributes;
636  fAxonColour->SetColour(G4Colour(G4Colour(0.5, 0.0, 0.0)));
637  fAxonColour->SetForceSolid(true);
638  fAxonColour->SetVisibility(true);
639 
640  // Spines in Dark-Green
641  fSpineColour = new G4VisAttributes;
642  fSpineColour->SetColour(G4Colour(G4Colour(0.0 , 100/255. , 0.0)));
643  fSpineColour->SetForceSolid(true);
644  fSpineColour->SetVisibility(true);
645 
646  // Whole neuron in semitransparent navy blue
647  fNeuronColour = new G4VisAttributes;
648  fNeuronColour->SetColour(G4Colour(G4Colour(0.0,0.4,0.8,0.5)));
649  fNeuronColour->SetForceSolid(true);
650  fNeuronColour->SetVisibility(true);
651 
652  }
653 }
654 
655 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
656 
657 // Load prepared data file of neural network with single and multiple layers
658 
660  (const G4String& filename)
661 {
662 
663  G4String sLine = "";
664  std::ifstream infile;
665  infile.open(filename.c_str());
666  if (!infile)
667  {
668 #ifdef GEANT4
669  G4cout<<" \n NeuronLoadDataFile::NeuralNetworkDATAfile >> datafile "
670  <<filename<<" not found !!!!"<<G4endl;
671  exit(0);
672 #endif
673  }
674  else
675  {
676 #ifdef G4VERBOSE
677  G4cout<< " NeuronLoadDataFile::NeuralNetworkDATAfile >> opening filename: "
678  << "\n" <<'\t'<<'\t'<<'\t'<<'\t'<<'\t'<<"' "<<filename
679  << " ' \n"<< G4endl;
680 #endif
681 
682  G4int nlines, nbSoma, nbDendrite;
683  nlines=0;
684  fnbSomacomp = 0 ; // total number of compartment into Soma
685  fnbDendritecomp = 0 ; // total number of compartment into Dendrites
686  fnbAxoncomp = 0 ; // total number of compartment into Axon
687  fnbSpinecomp = 0 ; // total number of compartment into Spines
688  G4double TotVolSoma, TotVolDend, TotVolAxon;
689  TotVolSoma=TotVolDend=TotVolAxon=0.;
690  G4double TotSurfSoma, TotSurfDend, TotSurfAxon;
691  TotSurfSoma=TotSurfDend=TotSurfAxon=0.;
692  //G4int nNmorph; // current index of neuronal morphology
693  G4int typeNcomp; // types of structure: soma, axon, apical dendrite, etc.
694  G4double x1,y1,z1,x2,y2,z2; // cartesian coordinates of each compartment
695  G4double radius; // radius of each compartment in micrometer
696  G4double height; // height of each compartment in micrometer
697  //G4double minX,minY,minZ; //minimum
698  //G4double maxX,maxY,maxZ; //maximum
699  G4double maxRad = -1e+09;
700  //minX=minY=minZ=1e+09;
701  //maxX=maxY=maxZ=-1e+09;
702  G4double density = 1.0 * (g/cm3) ; // water medium
703  G4double Piconst = (4.0/3.0)*pi ;
704 
705  while (getline(infile, sLine))
706  {
707  std::istringstream form(sLine);
708  if (nlines == 0) {
709  // to read total number of compartments
710  form >> fnbNeuroncomp >> nbSoma >> nbDendrite ;
711  fMassSomacomp = new G4double[nbSoma];
712  fMassSomaTot = 0.0 ;
713  fPosSomacomp = new G4ThreeVector[nbSoma];
714  fRadSomacomp = new G4double[nbSoma];
715  fRadDendcomp = new G4double[nbDendrite];
716  fHeightDendcomp = new G4double[nbDendrite];
717  fMassDendcomp = new G4double[nbDendrite];
718  fMassDendTot = 0.0 ;
719  fDistADendSoma = new G4double[nbDendrite];
720  fDistBDendSoma = new G4double[nbDendrite];
721  fPosDendcomp = new G4ThreeVector[nbDendrite];
722  fRotDendcomp = new G4RotationMatrix[nbDendrite];
723  }
724  // =======================================================================
725  // Soma compartments represented as Sphere or Ellipsoid solid
726  if (nlines > 0 && nlines <= nbSoma) // Total number of Soma compartments
727  {
728  form >> typeNcomp >> x1 >> y1 >> z1 >> radius ;
729  if (typeNcomp !=1) break;
730  // max diameter of compartments
731  if (maxRad < radius) maxRad = radius;
732  // Sphere volume and surface area
733  G4double VolSomacomp = Piconst*pow(radius*um,3.) ;
734  TotVolSoma = TotVolSoma + VolSomacomp;
735  G4double SurSomacomp = 3.*Piconst*pow(radius*um,2.) ;
736  TotSurfSoma = TotSurfSoma + SurSomacomp;
737  fMassSomacomp[fnbSomacomp] = density*VolSomacomp;
738  fMassSomaTot = fMassSomaTot + fMassSomacomp[fnbSomacomp];
739  // OR
740  // Ellipsoid volume and Approximate formula of surface area
741  //G4double VolSomacomp = Piconst*(Ra*um)*(Rb*um)*(Rc*um);
742  //G4double SurSomacomp = 3.*Piconst*pow((pow(Ra,1.6075)*pow(Rb,1.6075)+
743  //pow(Ra,1.6075)*pow(Rc,1.6075)+pow(Rb,1.6075)*pow(Rc,1.6075))/3.,0.622084);
744 
745  G4ThreeVector vSoma (x1 ,y1 ,z1);
746  fPosSomacomp [fnbSomacomp] = vSoma;
747  fRadSomacomp [fnbSomacomp]= radius;
748 
749  // RotationMatrix for Ellipsoid solid
750  // ....
751  fnbSomacomp++ ;
752  }
753  // =======================================================================
754  // Apical and basal dendritic compartments represented as cylinderical solid
755  if (nlines > nbSoma && nlines <= fnbNeuroncomp)
756  {
757  form >> typeNcomp >> x1 >> y1 >> z1 >> x2 >> y2 >> z2 >> radius >> height;
758  if (typeNcomp != 3 ) break; // || typeNcomp != 4
759 
760  // To calculate length, center and rotation angles of each cylinder
761  // Center-position of each cylinder
762  G4double Dendxx= x1 + x2;
763  G4double Dendyy= y1 + y2;
764  G4double Dendzz= z1 + z2;
765  G4ThreeVector translmDend = G4ThreeVector(Dendxx/2. ,
766  Dendyy/2. , Dendzz/2.) ;
767  fPosDendcomp [fnbDendritecomp] = translmDend;
768  fRadDendcomp [fnbDendritecomp]= radius;
769  G4double lengthDendcomp = height;
770  // Height of compartment
771  fHeightDendcomp [fnbDendritecomp]= lengthDendcomp;
772  // Distance from Soma
773 
774  // Cylinder volume and surface area
775  G4double VolDendcomp = pi*pow(radius*um,2)*(lengthDendcomp*um);
776  TotVolDend = TotVolDend + VolDendcomp;
777  G4double SurDendcomp = 2.*pi*radius*um*(radius+lengthDendcomp)*um;
778  TotSurfDend = TotSurfDend + SurDendcomp;
779  fMassDendcomp[fnbDendritecomp] = density*VolDendcomp;
780  fMassDendTot = fMassDendTot + fMassDendcomp[fnbDendritecomp];
781 
782  G4double Dendx= x1 - x2;
783  G4double Dendy= y1 - y2;
784  G4double Dendz= z1 - z2;
785  Dendx=Dendx/lengthDendcomp;
786  Dendy=Dendy/lengthDendcomp;
787  Dendz=Dendz/lengthDendcomp;
788 
789  // Euler angles of each compartment
790  G4ThreeVector directionDend = G4ThreeVector(Dendx,Dendy,Dendz);
791  G4double theta_eulerDend = directionDend.theta();
792  G4double phi_eulerDend = directionDend.phi();
793  G4double psi_eulerDend = 0;
794 
795  //Rotation Matrix, Euler constructor build inverse matrix.
796  G4RotationMatrix rotmDendInv = G4RotationMatrix(
797  phi_eulerDend+pi/2,
798  theta_eulerDend,
799  psi_eulerDend);
800  G4RotationMatrix rotmDend = rotmDendInv.inverse();
801 
802  fRotDendcomp [fnbDendritecomp]= rotmDend ;
803  //G4Transform3D transformDend = G4Transform3D(rotmDend,translmDend);
804  //fRotTransDendPos [fnbDendritecomp]= transformDend ;
805  fnbDendritecomp++ ;
806 
807  }
808 
809  nlines++;
810  }
811 
812  // =======================================================================
813 
814  G4cout << " Total number of compartments into Neuron : "<<
815  fnbNeuroncomp <<G4endl;
816  G4cout << "\n"<<G4endl;
817 
818  // to calculate SHIFT value for neuron translation
819  fshiftX = 0.; //(minX + maxX)/2. ;
820  fshiftY = 0.; //(minY + maxY)/2. ;
821  fshiftZ = 0.; //(minZ + maxZ)/2. ;
822 
823  // width, height, depth of bounding slice volume
824  //maxRad = 0.0 ;
825  fwidthB = 640.;
826  fheightB = 280.;
827  fdepthB = 25.;
828  // diagonal length of bounding slice, that give diameter of sphere
829  // for particle direction and fluence!
830  fdiagnlLength = std::sqrt(fwidthB*fwidthB + fheightB*fheightB
831  + fdepthB*fdepthB);
832 
833  fTotVolNeuron = TotVolSoma+TotVolDend+TotVolAxon;
834  fTotSurfNeuron = TotSurfSoma+TotSurfDend+TotSurfAxon;
835  fTotMassNeuron = fMassSomaTot+fMassDendTot+fMassAxonTot;
836 
837  fTotVolSlice = fwidthB*um*fheightB*um*fdepthB*um;
838  fTotSurfSlice = 2*(fwidthB*um*fheightB*um+fheightB*um*fdepthB*um+
839  fwidthB*um*fdepthB*um);
840  fTotMassSlice = 1.0 * (g/cm3) *fTotVolSlice;
841 
842  fTotVolMedium = Piconst*pow(fdiagnlLength*um/2.,3.) ;
843  fTotSurfMedium = 3.*Piconst*pow(fdiagnlLength*um/2.,2);
844  fTotMassMedium = 1.0 * (g/cm3) *fTotVolMedium;
845 
846  }
847  infile.close();
848 }
849 
850 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
851 
853 {
854 delete[] fMassSomacomp ;
855 delete[] fPosSomacomp ;
856 delete[] fRadSomacomp ;
857 delete[] fRadDendcomp ;
858 delete[] fHeightDendcomp;
859 delete[] fMassDendcomp ;
860 delete[] fDistADendSoma ;
861 delete[] fDistBDendSoma ;
862 delete[] fPosDendcomp ;
863 delete[] fRotDendcomp ;
864 delete[] fRadAxoncomp ;
865 delete[] fHeightAxoncomp;
866 delete[] fMassAxoncomp ;
867 delete[] fDistAxonsoma ;
868 delete[] fPosAxoncomp ;
869 delete[] fRotAxoncomp ;
870 delete[] fRadNeuroncomp ;
871 delete[] fHeightNeuroncomp;
872 delete[] fMassNeuroncomp ;
873 delete[] fDistNeuronsoma ;
874 delete[] fPosNeuroncomp ;
875 delete[] fRotNeuroncomp ;
876 }
877 
878 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
879 
881 (const G4int copyNo, G4VPhysicalVolume* physVol) const
882 {
883 /*
884 // sphere transformation
885  G4ThreeVector
886  originSoma(
887  (fPosSomacomp[copyNo].x()-fshiftX) * um,
888  (fPosSomacomp[copyNo].y()-fshiftY) * um,
889  (fPosSomacomp[copyNo].z()-fshiftZ) * um
890  );
891  physVol->SetRotation(0);
892  physVol->SetTranslation(originSoma);
893  */
894 
895 // cylinder rotation and transformation
896 
897 // to calculate Euler angles from Rotation Matrix after Inverse!
898 //
899  G4RotationMatrix rotmNeuron = G4RotationMatrix(fRotNeuroncomp[copyNo]);
900  G4double cosX = std::sqrt (rotmNeuron.xx()*rotmNeuron.xx() +
901  rotmNeuron.yx()*rotmNeuron.yx()) ;
902  G4double euX, euY, euZ;
903  if (cosX > 16*FLT_EPSILON)
904  {
905  euX = std::atan2 (rotmNeuron.zy(),rotmNeuron.zz());
906  euY = std::atan2 (-rotmNeuron.zx(),cosX);
907  euZ = std::atan2 (rotmNeuron.yx(),rotmNeuron.xx());
908  }
909  else
910  {
911  euX = std::atan2 (-rotmNeuron.yz(),rotmNeuron.yy());
912  euY = std::atan2 (-rotmNeuron.zx(),cosX);
913  euZ = 0. ;
914  }
915  G4RotationMatrix* rot = new G4RotationMatrix();
916  rot->rotateX(euX);
917  rot->rotateY(euY);
918  rot->rotateZ(euZ);
919 
920  physVol->SetRotation(rot);
921 
922 // shift of cylinder compartments
924  originNeuron(
925  (fPosNeuroncomp[copyNo].x()-fshiftX) * um,
926  (fPosNeuroncomp[copyNo].y()-fshiftY) * um,
927  (fPosNeuroncomp[copyNo].z()-fshiftZ) * um
928  );
929  physVol->SetTranslation(originNeuron);
930 
931 }
932 
933 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
934 /*
935 G4VSolid* NeuronLoadDataFile::ComputeSolid(const G4int copyNo,
936  G4VPhysicalVolume* physVol )
937 {
938  G4VSolid* solid;
939  if( typeNcomp[copyNo] == 1 )
940  {
941  solid = sphereComp ;
942  }
943  else if( typeNcomp[copyNo] == 3 || typeNcomp[copyNo] == 4 )
944  {
945  solid = cylinderComp ;
946  }
947  return solid;
948 }
949 */
950 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
951 
953 (G4Tubs& fcylinderComp, const G4int copyNo, const G4VPhysicalVolume*) const
954 {
955  fcylinderComp.SetInnerRadius(0*um);
956  fcylinderComp.SetOuterRadius(fRadNeuroncomp[copyNo]*um);
957  fcylinderComp.SetZHalfLength(fHeightNeuroncomp[copyNo]*um /2.);
958  fcylinderComp.SetStartPhiAngle(0.*deg);
959  fcylinderComp.SetDeltaPhiAngle(360.*deg);
960 }
961 
962 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
963 /*
964 void NeuronLoadDataFile::ComputeDimensions
965 (G4Sphere& sphereComp, const G4int copyNo, const G4VPhysicalVolume*) const
966 {
967  fsphereComp.SetInnerRadius(0);
968  fsphereComp.SetOuterRadius(fRadSomacomp[copyNo] * um);
969  fsphereComp.SetStartPhiAngle(0.*deg);
970  fsphereComp.SetDeltaPhiAngle(360.*deg);
971  fsphereComp.SetStartThetaAngle(0.*deg);
972  fsphereComp.SetDeltaThetaAngle(180.*deg);
973 }
974 */
975 /*
976 #if 1
977 (G4Sphere& somaS, const G4int copyNo, const G4VPhysicalVolume* physVol) const
978 #else
979 (G4Tubs& dendritesS, const G4int copyNo, const G4VPhysicalVolume* ) const
980 #endif
981 {
982  G4LogicalVolume* LogicalVolume = physVol->GetLogicalVolume();
983 
984  G4PhysicalVolume* somaPV = LogicalVolume->GetDaughter(0);
985  G4LogicalVolume* somaLV = somaPV->GetLogicalVolume();
986  G4Sphere* somaS = (G4Sphere*)somaLV->GetSolid();
987 
988  G4PhysicalVolume* dendritesPV = LogicalVolume->GetDaughter(0);
989  G4LogicalVolume* dendritesLV = dendritesPV->GetLogicalVolume();
990  G4Tubs* dendritesS = (G4Tubs*)dendritesLV->GetSolid();
991 }
992 */
993 
994 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Float_t x
Definition: compare.C:6
virtual const G4String & GetOption()
CLHEP::Hep3Vector G4ThreeVector
double minY
Definition: plot_hist.C:9
void SetInnerRadius(G4double newRMin)
void SetOuterRadius(G4double newRMax)
Float_t y1[n_points_granero]
Definition: compare.C:5
double yx() const
Implementation of the NeuronLoadDataFile class.
#define FLT_EPSILON
Definition: templates.hh:103
Definition: G4Tubs.hh:85
double zz() const
Float_t x1[n_points_granero]
Definition: compare.C:5
CLHEP::HepRotation G4RotationMatrix
#define G4endl
Definition: G4ios.hh:61
Float_t y
Definition: compare.C:6
Double_t z
double yz() const
double z() const
void SetTranslation(const G4ThreeVector &v)
double maxY
Definition: plot_hist.C:10
void SingleNeuronSWCfile(const G4String &filename)
static constexpr double um
Definition: G4SIunits.hh:113
double theta() const
Float_t y2[n_points_geant4]
Definition: compare.C:26
static const G4int maxZ
nlines
static constexpr double g
Definition: G4SIunits.hh:183
void ComputeDimensions(G4Tubs &cylinderComp, const G4int copyNo, const G4VPhysicalVolume *) const
double zx() const
double G4double
Definition: G4Types.hh:76
double xx() const
void SetStartPhiAngle(G4double newSPhi, G4bool trig=true)
static constexpr double deg
Definition: G4SIunits.hh:152
void SetZHalfLength(G4double newDz)
Double_t radius
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
void NeuralNetworkDATAfile(const G4String &filename)
void SetDeltaPhiAngle(G4double newDPhi)
HepRotation & rotateY(double delta)
Definition: Rotation.cc:79
double zy() const
int G4int
Definition: G4Types.hh:78
double phi() const
double yy() const
HepRotation & rotateZ(double delta)
Definition: Rotation.cc:92
G4GLOB_DLL std::ostream G4cout
double x() const
void SetColour(const G4Colour &)
void SetRotation(G4RotationMatrix *)
static constexpr double pi
Definition: G4SIunits.hh:75
double y() const
Float_t x2[n_points_geant4]
Definition: compare.C:26
static constexpr double cm3
Definition: G4SIunits.hh:121
HepRotation & rotateX(double delta)
Definition: Rotation.cc:66
HepRotation inverse() const