QwGeant4
QweakSimBeamLine Class Reference

Definition of the BeamLine. More...

#include <QweakSimBeamLine.hh>

+ Collaboration diagram for QweakSimBeamLine:

Public Member Functions

 QweakSimBeamLine (QweakSimUserInformation *myUI)
 
 ~QweakSimBeamLine ()
 
void ConstructComponent (G4VPhysicalVolume *)
 
void DestroyComponent ()
 
void SetBeamLineCenterPositionInZ (G4double)
 
G4double GetBeamLineCenterPositionInZ ()
 
void SetBeamLineMaterial (G4String)
 
G4String GetBeamLineMaterial ()
 

Private Attributes

QweakSimMaterialpMaterial
 
G4Material * BeamPipe_Material
 
G4Material * Shield_Material
 
G4LogicalVolume * R1_Pipe_Logical
 
G4LogicalVolume * R1_Flange_Logical
 
G4LogicalVolume * DS_R1_Pipe_Logical
 
G4LogicalVolume * DS_R1_Flange_Logical
 
G4LogicalVolume * DS_R1_Bellow_Logical
 
G4LogicalVolume * R2_Pipe_Logical
 
G4LogicalVolume * R2_Flange_Logical
 
G4LogicalVolume * DS_R2_Pipe_Logical
 
G4LogicalVolume * R2_RotatorPipe_Logical
 
G4LogicalVolume * R3_Pipe_Logical
 
G4LogicalVolume * R3_Flange_Logical
 
G4LogicalVolume * DS_R3_Pipe_Logical
 
G4LogicalVolume * R3_US_Wall_Pipe_Logical
 
G4LogicalVolume * Sub_LeadBox_in_Wall_Logical
 
G4LogicalVolume * R3_Wall_Pipe_Logical
 
G4LogicalVolume * Sub_LeadBox_Extent_Logical
 
G4LogicalVolume * DS_18inch_Pipe1_Logical
 
G4LogicalVolume * DS_18inch_Pipe2_Logical
 
G4LogicalVolume * DS_24inch_Pipe_Logical
 
G4LogicalVolume * DS_24inch_Pipe_Flange_Logical
 
G4double beamline_ZPos
 
QweakSimBeamLineMessengerbeamlineMessenger
 
QweakSimUserInformationmyUserInfo
 

Detailed Description

Definition of the BeamLine.

Placeholder for a long explaination

Definition at line 60 of file QweakSimBeamLine.hh.

Constructor & Destructor Documentation

QweakSimBeamLine::QweakSimBeamLine ( QweakSimUserInformation myUI)

Definition at line 36 of file QweakSimBeamLine.cc.

References beamlineMessenger, BeamPipe_Material, DS_18inch_Pipe1_Logical, DS_18inch_Pipe2_Logical, DS_24inch_Pipe_Flange_Logical, DS_24inch_Pipe_Logical, DS_R1_Bellow_Logical, DS_R1_Flange_Logical, DS_R1_Pipe_Logical, DS_R2_Pipe_Logical, DS_R3_Pipe_Logical, QweakSimMaterial::GetInstance(), QweakSimMaterial::GetMaterial(), myUserInfo, pMaterial, QweakSimBeamLineMessenger::QweakSimBeamLineMessenger(), R1_Flange_Logical, R1_Pipe_Logical, R2_Flange_Logical, R2_Pipe_Logical, R2_RotatorPipe_Logical, R3_Flange_Logical, R3_Pipe_Logical, R3_US_Wall_Pipe_Logical, R3_Wall_Pipe_Logical, Shield_Material, Sub_LeadBox_Extent_Logical, and Sub_LeadBox_in_Wall_Logical.

37 {
38  G4cout << G4endl << "###### Calling QweakSimBeamLine::QweakBeamLine() " << G4endl << G4endl;
39 
40  myUserInfo = myUI;
41 
42  //BeamLineContainer_Logical = NULL;
43  //BeamLineContainer_Physical = NULL;
44 
45  R1_Pipe_Logical= NULL;
46  R1_Flange_Logical= NULL;
47  DS_R1_Pipe_Logical= NULL;
50  R2_Pipe_Logical= NULL;
51  R2_Flange_Logical= NULL;
52  DS_R2_Pipe_Logical= NULL;
54  R3_Pipe_Logical= NULL;
55  R3_Flange_Logical= NULL;
56  DS_R3_Pipe_Logical = NULL;
65 
67 
68  //BeamLineContainer_Material = pMaterial->GetMaterial("Vacuum");
71 
73 }
G4Material * BeamPipe_Material
G4Material * Shield_Material
G4LogicalVolume * R1_Pipe_Logical
static QweakSimMaterial * GetInstance()
G4LogicalVolume * R1_Flange_Logical
QweakSimMaterial * pMaterial
G4LogicalVolume * Sub_LeadBox_in_Wall_Logical
G4LogicalVolume * R2_Pipe_Logical
G4LogicalVolume * R3_Flange_Logical
G4LogicalVolume * DS_R1_Pipe_Logical
G4LogicalVolume * Sub_LeadBox_Extent_Logical
G4LogicalVolume * DS_24inch_Pipe_Logical
G4LogicalVolume * DS_R1_Bellow_Logical
Placeholder for a long explaination.
G4Material * GetMaterial(G4String material)
QweakSimBeamLineMessenger * beamlineMessenger
G4LogicalVolume * DS_18inch_Pipe1_Logical
G4LogicalVolume * DS_R3_Pipe_Logical
G4LogicalVolume * R2_RotatorPipe_Logical
G4LogicalVolume * R2_Flange_Logical
G4LogicalVolume * R3_Pipe_Logical
G4LogicalVolume * R3_US_Wall_Pipe_Logical
G4LogicalVolume * DS_R2_Pipe_Logical
QweakSimUserInformation * myUserInfo
G4LogicalVolume * DS_18inch_Pipe2_Logical
G4LogicalVolume * DS_R1_Flange_Logical
G4LogicalVolume * DS_24inch_Pipe_Flange_Logical
G4LogicalVolume * R3_Wall_Pipe_Logical

+ Here is the call graph for this function:

QweakSimBeamLine::~QweakSimBeamLine ( )

Definition at line 76 of file QweakSimBeamLine.cc.

References beamlineMessenger.

77 {
78  delete beamlineMessenger;
79 }
QweakSimBeamLineMessenger * beamlineMessenger

Member Function Documentation

void QweakSimBeamLine::ConstructComponent ( G4VPhysicalVolume *  MotherVolume)

Definition at line 83 of file QweakSimBeamLine.cc.

References BeamPipe_Material, DS_18inch_Pipe1_Logical, DS_18inch_Pipe2_Logical, DS_24inch_Pipe_Flange_Logical, DS_24inch_Pipe_Logical, DS_R1_Bellow_Logical, DS_R1_Flange_Logical, DS_R1_Pipe_Logical, DS_R2_Pipe_Logical, DS_R3_Pipe_Logical, pSurfChk, R1_Flange_Logical, R1_Pipe_Logical, R2_Flange_Logical, R2_Pipe_Logical, R2_RotatorPipe_Logical, R3_Flange_Logical, R3_Pipe_Logical, R3_US_Wall_Pipe_Logical, R3_Wall_Pipe_Logical, Shield_Material, Sub_LeadBox_Extent_Logical, and Sub_LeadBox_in_Wall_Logical.

84 {
85  G4cout << G4endl << "###### Calling QweakSimBeamLine::ConstructComponent() " << G4endl << G4endl;
86 
87 
88  // define BeamLine position
89  G4ThreeVector positionBeamPipe = G4ThreeVector(0,0,0);
90 
91 //--------------------------------------
92  // define beampipe solid volume
93  G4cout << G4endl << "###### QweakSimBeamLine: Define Beam Pipe" << G4endl << G4endl;
94 
95  // define R1 - GEM Region pipe
96  G4double collimator_center = -422.0195*cm;
97  G4Tubs* R1_Pipe_Solid = new G4Tubs("R1_Pipe_Solid",
98  2.0447*cm,
99  2.413*cm,
100  24.9174*cm,
101  0,
102  360*degree);
103 
104  R1_Pipe_Logical = new G4LogicalVolume(R1_Pipe_Solid,
106  "R1_Pipe_Log",
107  0,0,0);
108 
109  new G4PVPlacement(0,
110  G4ThreeVector(0,0,collimator_center-115.4186*cm),
111  "R1_Pipe",
113  MotherVolume,
114  false,
115  0,
116  pSurfChk);
117 
118  // define R1 - GEM Region flange
119  G4Tubs* R1_Flange_Solid = new G4Tubs("R1_Flange_Solid",
120  2.413*cm,
121  4.445*cm,
122  0.8*cm,
123  0,
124  360*degree);
125 
126  R1_Flange_Logical = new G4LogicalVolume(R1_Flange_Solid,
128  "R1_Flange_Log",
129  0,0,0);
130 
131  new G4PVPlacement(0,
132  G4ThreeVector(0,0,collimator_center-139.636*cm),
133  "R1_Flange",
135  MotherVolume,
136  false,
137  0,
138  pSurfChk);
139 
140  // define Downstream R1 - GEM Region pipe
141  G4Tubs* DS_R1_Pipe_Solid = new G4Tubs("DS_R1_Pipe_Solid",
142  3.8862*cm,
143  4.445*cm,
144  63.45085*cm,
145  0,
146  360*degree);
147 
148  DS_R1_Pipe_Logical = new G4LogicalVolume(DS_R1_Pipe_Solid,
150  "DS_R1_Pipe_Log",
151  0,0,0);
152 
153  new G4PVPlacement(0,
154  G4ThreeVector(0,0,collimator_center-27.15035*cm),
155  "DS_R1_Pipe",
157  MotherVolume,
158  false,
159  0,
160  pSurfChk);
161 
162  // define Downstream R1 - GEM Region flange
163  G4Tubs* DS_R1_Flange_Solid = new G4Tubs("DS_R1_Flange_Solid",
164  2.413*cm,
165  4.445*cm,
166  0.5588*cm,
167  0,
168  360*degree);
169 
170  DS_R1_Flange_Logical = new G4LogicalVolume(DS_R1_Flange_Solid,
172  "DS_R1_Flange_Log",
173  0,0,0);
174 
175  new G4PVPlacement(0,
176  G4ThreeVector(0,0,collimator_center-90.0424*cm),
177  "DS_R1_Flange",
179  MotherVolume,
180  false,
181  0,
182  pSurfChk);
183 
184  // define Downstream R1 - GEM Region bellow
185  G4Tubs* DS_R1_Bellow_Solid = new G4Tubs("DS_R1_Bellow_Solid",
186  4.445*cm,
187  13.5*cm,
188  7.5*cm,
189  0,
190  360*degree);
191 
192  DS_R1_Bellow_Logical = new G4LogicalVolume(DS_R1_Bellow_Solid,
194  "DS_R1_Bellow_Log",
195  0,0,0);
196 
197  new G4PVPlacement(0,
198  G4ThreeVector(0,0,collimator_center+28.8*cm),
199  "DS_R1_Bellow",
201  MotherVolume,
202  false,
203  0,
204  pSurfChk);
205 
206  // define R2 - inside primary collimator pipe
207  G4Tubs* R2_Pipe_Solid = new G4Tubs("R2_Pipe_Solid",
208  10.0711*cm,
209  10.9479*cm,
210  7.5*cm,
211  0,
212  360*degree);
213 
214  R2_Pipe_Logical = new G4LogicalVolume(R2_Pipe_Solid,
216  "R2_Pipe_Log",
217  0,0,0);
218 
219  new G4PVPlacement(0,
220  G4ThreeVector(0,0,collimator_center+43.8005*cm),
221  "R2_Pipe",
223  MotherVolume,
224  false,
225  0,
226  pSurfChk);
227 
228  // define R2 - inside primary collimator flange
229  G4Tubs* R2_Flange_Solid = new G4Tubs("R2_Flange_Solid",
230  4.445*cm,
231  10.9479*cm,
232  0.4384*cm,
233  0,
234  360*degree);
235 
236  R2_Flange_Logical = new G4LogicalVolume(R2_Flange_Solid,
238  "R2_Flange_Log",
239  0,0,0);
240 
241  new G4PVPlacement(0,
242  G4ThreeVector(0,0,collimator_center+43.8005*cm-7.0616*cm),
243  "R2_Flange",
245  MotherVolume,
246  false,
247  0,
248  pSurfChk);
249 
250  // define downstream R2 pipe
251  G4Tubs* DS_R2_Pipe_Solid = new G4Tubs("DS_R2_Pipe_Solid",
252  10.0711*cm,
253  10.9479*cm,
254  49.43*cm,
255  0,
256  360*degree);
257 
258  DS_R2_Pipe_Logical = new G4LogicalVolume(DS_R2_Pipe_Solid,
260  "DS_R2_Pipe_Log",
261  0,0,0);
262 
263  new G4PVPlacement(0,
264  G4ThreeVector(0,0,collimator_center+100.7305*cm),
265  "DS_R2_Pipe",
267  MotherVolume,
268  false,
269  0,
270  pSurfChk);
271 
272  // define R2 - rotator pipe
273  G4Tubs* R2_RotatorPipe_Solid = new G4Tubs("R2_RotatorPipe_Solid",
274  15.1613*cm,
275  16.1925*cm,
276  42.4307*cm,
277  0,
278  360*degree);
279 
280  R2_RotatorPipe_Logical = new G4LogicalVolume(R2_RotatorPipe_Solid,
282  "R2_RotatorPipe_Log",
283  0,0,0);
284 
285  new G4PVPlacement(0,
286  G4ThreeVector(0,0,collimator_center+107.7298*cm),
287  "R2_RotatorPipe",
289  MotherVolume,
290  false,
291  0,
292  pSurfChk);
293 
294 
295  // define R3 - inside collimator pipe
296  G4Tubs* R3_Pipe_Solid = new G4Tubs("R3_Pipe_Solid",
297  10.0711*cm,
298  10.9479*cm,
299  5.615*cm,
300  0,
301  360*degree);
302 
303  R3_Pipe_Logical = new G4LogicalVolume(R3_Pipe_Solid,
305  "R3_Pipe_Log",
306  0,0,0);
307 
308  new G4PVPlacement(0,
309  G4ThreeVector(0,0,collimator_center+155.7755*cm),
310  "R3_Pipe",
312  MotherVolume,
313  false,
314  0,
315  pSurfChk);
316 
317  // define R3 - after collimator pipe
318  G4Tubs* DS_R3_Pipe_Solid = new G4Tubs("DS_R3_Pipe_Solid",
319  12.7381*cm,
320  13.6525*cm,
321  256.8575*cm,
322  0,
323  360*degree);
324 
325  DS_R3_Pipe_Logical = new G4LogicalVolume(DS_R3_Pipe_Solid,
327  "DS_R3_Pipe_Log",
328  0,0,0);
329 
330  new G4PVPlacement(0,
331  G4ThreeVector(0,0,-3.9065*cm),
332  "DS_R3_Pipe",
334  MotherVolume,
335  false,
336  0,
337  pSurfChk);
338 
339  // define R3 - flange
340  G4Tubs* R3_Flange_Solid = new G4Tubs("R3_Flange_Solid",
341  13.6525*cm,
342  22.86*cm,
343  1.4224*cm,
344  0,
345  360*degree);
346 
347  R3_Flange_Logical = new G4LogicalVolume(R3_Flange_Solid,
349  "R3_Flange_Log",
350  0,0,0);
351 
352  new G4PVPlacement(0,
353  G4ThreeVector(0,0,251.5286*cm),
354  "R3_Flange",
356  MotherVolume,
357  false,
358  0,
359  pSurfChk);
360 
361  // define R3 - before shield wall
362  G4Tubs* R3_US_Wall_Pipe_Solid = new G4Tubs("R3_US_Wall_Pipe_Solid",
363  21.9075*cm,
364  22.86*cm,
365  23.5245*cm,
366  0,
367  360*degree);
368 
369  R3_US_Wall_Pipe_Logical = new G4LogicalVolume(R3_US_Wall_Pipe_Solid,
371  "R3_US_Wall_Pipe_Log",
372  0,0,0);
373 
374  new G4PVPlacement(0,
375  G4ThreeVector(0,0,276.4755*cm),
376  "R3_US_Wall_Pipe",
378  MotherVolume,
379  false,
380  0,
381  pSurfChk);
382 
383  // define 2" lead box in the shield wall
384  G4Box* LeadBox_in_Wall_Solid = new G4Box("LeadBox_in_Wall_Solid",
385  27.94*cm,
386  27.94*cm,
387  40*cm);
388 
389  G4Tubs* R3_Wall_Pipe_Vacuum_Solid = new G4Tubs("R3_Wall_Pipe_Vacuum_Solid",
390  0.0*cm,
391  (22.86+0.1)*cm,
392  (40.0+0.1)*cm,
393  0,
394  360*degree);
395 
396  G4SubtractionSolid* Sub_LeadBox_in_Wall_Solid = new G4SubtractionSolid("LeadBox_in_Wall_Solid - R3_Wall_Pipe_Vacuum_Solid",
397  LeadBox_in_Wall_Solid,
398  R3_Wall_Pipe_Vacuum_Solid);
399 
400  Sub_LeadBox_in_Wall_Logical = new G4LogicalVolume(Sub_LeadBox_in_Wall_Solid,
402  "Sub_LeadBox_in_Wall_Log",
403  0,0,0);
404 
405  new G4PVPlacement(0,
406  G4ThreeVector(0,0,340.0*cm),
407  "Sub_LeadBox_in_Wall",
409  MotherVolume,
410  false,
411  0,
412  pSurfChk);
413 
414  // define R3 - pipe in shield wall
415  G4Tubs* R3_Wall_Pipe_Solid = new G4Tubs("R3_Wall_Pipe_Solid",
416  21.9075*cm,
417  22.86*cm,
418  40.0*cm,
419  0,
420  360*degree);
421 
422  R3_Wall_Pipe_Logical = new G4LogicalVolume(R3_Wall_Pipe_Solid,
424  "R3_Wall_Pipe_Log",
425  0,0,0);
426 
427  new G4PVPlacement(0,
428  G4ThreeVector(0,0,340.0*cm),
429  "R3_Wall_Pipe",
431  MotherVolume,
432  false,
433  0,
434  pSurfChk);
435 
436  // define 2" lead box, extends only to Z=622.625
437  G4Box* LeadBox_Extent_Solid = new G4Box("LeadBox_Extent_Solid",
438  27.94*cm,
439  27.94*cm,
440  121.3125*cm);
441 
442  G4Tubs* DS_18inch_Pipe1_Vacuum_Solid = new G4Tubs("DS_18inch_Pipe1_Vacuum_Solid",
443  0.0*cm,
444  (22.86+0.1)*cm,
445  (121.3125+0.1)*cm,
446  0,
447  360*degree);
448 
449  G4SubtractionSolid* Sub_LeadBox_Extent_Solid = new G4SubtractionSolid("LeadBox_Extent_Solid - DS_18inch_Pipe1_Vacuum_Solid",
450  LeadBox_Extent_Solid,
451  DS_18inch_Pipe1_Vacuum_Solid);
452 
453  Sub_LeadBox_Extent_Logical = new G4LogicalVolume(Sub_LeadBox_Extent_Solid,
455  "Sub_LeadBox_Extent_Log",
456  0,0,0);
457 
458  new G4PVPlacement(0,
459  G4ThreeVector(0,0,501.3125*cm),
460  "Sub_LeadBox_Extent",
462  MotherVolume,
463  false,
464  0,
465  pSurfChk);
466 
467 
468  // define Downstream 18" pipe within lead box
469  G4Tubs* DS_18inch_Pipe1_Solid = new G4Tubs("DS_18inch_Pipe1_Solid",
470  21.9075*cm,
471  22.86*cm,
472  121.3125*cm,
473  0,
474  360*degree);
475 
476  DS_18inch_Pipe1_Logical = new G4LogicalVolume(DS_18inch_Pipe1_Solid,
478  "DS_18inch_Pipe1_Log",
479  0,0,0);
480 
481  new G4PVPlacement(0,
482  G4ThreeVector(0,0,501.3125*cm),
483  "DS_18inch_Pipe1",
485  MotherVolume,
486  false,
487  0,
488  pSurfChk);
489 
490  // define Downstream 18" pipe after lead box
491  G4Tubs* DS_18inch_Pipe2_Solid = new G4Tubs("DS_18inch_Pipe2_Solid",
492  21.9075*cm,
493  22.86*cm,
494  7.568*cm,
495  0,
496  360*degree);
497 
498  DS_18inch_Pipe2_Logical = new G4LogicalVolume(DS_18inch_Pipe2_Solid,
500  "DS_18inch_Pipe2_Log",
501  0,0,0);
502 
503  new G4PVPlacement(0,
504  G4ThreeVector(0,0,630.193*cm),
505  "DS_18inch_Pipe2",
507  MotherVolume,
508  false,
509  0,
510  pSurfChk);
511 
512  // define Downstream 24" pipe
513  G4Tubs* DS_24inch_Pipe_Solid = new G4Tubs("DS_24inch_Pipe_Solid",
514  29.5275*cm,
515  30.48*cm,
516  31.1195*cm,
517  0,
518  360*degree);
519 
520  DS_24inch_Pipe_Logical = new G4LogicalVolume(DS_24inch_Pipe_Solid,
522  "DS_24inch_Pipe_Log",
523  0,0,0);
524 
525  new G4PVPlacement(0,
526  G4ThreeVector(0,0,668.8805*cm),
527  "DS_24inch_Pipe",
529  MotherVolume,
530  false,
531  0,
532  pSurfChk);
533 
534  // define Downstream 24" pipe flange
535  G4Tubs* DS_24inch_Pipe_Flange_Solid = new G4Tubs("DS_24inch_Pipe_Flange_Solid",
536  22.86*cm,
537  29.5275*cm,
538  1.4224*cm,
539  0,
540  360*degree);
541 
542  DS_24inch_Pipe_Flange_Logical = new G4LogicalVolume(DS_24inch_Pipe_Flange_Solid,
544  "DS_24inch_Pipe_Flange_Log",
545  0,0,0);
546 
547  new G4PVPlacement(0,
548  G4ThreeVector(0,0,(637.761+1.4224)*cm),
549  "DS_24inch_Pipe",
551  MotherVolume,
552  false,
553  0,
554  pSurfChk);
555 
556 //--------------------------------------
557 
558  G4cout << G4endl << "###### QweakSimBeamLine: Setting Attributes " << G4endl << G4endl;
559 
560  G4Colour red (1.,0.,0.);
561  G4Colour orange ( 255/255., 127/255., 0/255.);
562  G4Colour blue ( 0/255., 0/255., 255/255.);
563  G4Colour magenta ( 255/255., 0/255., 255/255.);
564  G4Colour grey ( 127/255., 127/255., 127/255.);
565  G4Colour lightgrey ( 220/255., 220/255., 220/255.);
566  G4Colour lightblue ( 139/255., 208/255., 255/255.);
567  G4Colour lightorange ( 255/255., 189/255., 165/255.);
568  G4Colour khaki3 ( 205/255., 198/255., 115/255.);
569  G4Colour brown (178/255., 102/255., 26/255.);
570 
571 
572 
573  G4VisAttributes* Shield_VisAtt = new G4VisAttributes(brown);
574  Shield_VisAtt -> SetVisibility(true);
575  //BeamPipe_VisAtt -> SetForceWireframe(true);
576  DS_R1_Bellow_Logical -> SetVisAttributes(Shield_VisAtt);
577  Sub_LeadBox_in_Wall_Logical->SetVisAttributes(Shield_VisAtt);
578  Sub_LeadBox_Extent_Logical->SetVisAttributes(Shield_VisAtt);
579  G4cout << G4endl << "###### Leaving QweakSimBeamLine::ConstructComponent() " << G4endl << G4endl;
580 
581  G4VisAttributes* Pipe_VisAtt = new G4VisAttributes(lightgrey);
582  Pipe_VisAtt -> SetVisibility(true);
583  R1_Pipe_Logical->SetVisAttributes(Pipe_VisAtt);
584  R1_Flange_Logical->SetVisAttributes(Pipe_VisAtt);
585  DS_R1_Pipe_Logical->SetVisAttributes(Pipe_VisAtt);
586  DS_R1_Flange_Logical->SetVisAttributes(Pipe_VisAtt);
587  R2_Pipe_Logical->SetVisAttributes(Pipe_VisAtt);
588  R2_Flange_Logical->SetVisAttributes(Pipe_VisAtt);
589  DS_R2_Pipe_Logical->SetVisAttributes(Pipe_VisAtt);
590  R2_RotatorPipe_Logical->SetVisAttributes(Pipe_VisAtt);
591  R3_Pipe_Logical->SetVisAttributes(Pipe_VisAtt);
592  R3_Flange_Logical->SetVisAttributes(Pipe_VisAtt);
593  DS_R3_Pipe_Logical->SetVisAttributes(Pipe_VisAtt);
594  R3_US_Wall_Pipe_Logical->SetVisAttributes(Pipe_VisAtt);
595  R3_Wall_Pipe_Logical->SetVisAttributes(Pipe_VisAtt);
596  DS_18inch_Pipe1_Logical->SetVisAttributes(Pipe_VisAtt);
597  DS_18inch_Pipe2_Logical->SetVisAttributes(Pipe_VisAtt);
598  DS_24inch_Pipe_Logical->SetVisAttributes(Pipe_VisAtt);
599  DS_24inch_Pipe_Flange_Logical->SetVisAttributes(Pipe_VisAtt);
600 
601 } // end of QweakSimBeamLine::ConstructComponent()
G4Material * BeamPipe_Material
G4Material * Shield_Material
G4LogicalVolume * R1_Pipe_Logical
G4LogicalVolume * R1_Flange_Logical
G4LogicalVolume * Sub_LeadBox_in_Wall_Logical
G4LogicalVolume * R2_Pipe_Logical
G4LogicalVolume * R3_Flange_Logical
static const G4bool pSurfChk
G4LogicalVolume * DS_R1_Pipe_Logical
G4LogicalVolume * Sub_LeadBox_Extent_Logical
G4LogicalVolume * DS_24inch_Pipe_Logical
G4LogicalVolume * DS_R1_Bellow_Logical
G4LogicalVolume * DS_18inch_Pipe1_Logical
G4LogicalVolume * DS_R3_Pipe_Logical
G4LogicalVolume * R2_RotatorPipe_Logical
G4LogicalVolume * R2_Flange_Logical
G4LogicalVolume * R3_Pipe_Logical
G4LogicalVolume * R3_US_Wall_Pipe_Logical
G4LogicalVolume * DS_R2_Pipe_Logical
G4LogicalVolume * DS_18inch_Pipe2_Logical
G4LogicalVolume * DS_R1_Flange_Logical
G4LogicalVolume * DS_24inch_Pipe_Flange_Logical
G4LogicalVolume * R3_Wall_Pipe_Logical
void QweakSimBeamLine::DestroyComponent ( )

Definition at line 642 of file QweakSimBeamLine.cc.

643 {
644 }
G4double QweakSimBeamLine::GetBeamLineCenterPositionInZ ( )

Definition at line 656 of file QweakSimBeamLine.cc.

References beamline_ZPos.

657 {
658  G4cout << G4endl << "###### Calling QweakSimBeamLine::GetBeamLineCenterPositionInZ() " << G4endl << G4endl;
659  return beamline_ZPos;
660 }
G4String QweakSimBeamLine::GetBeamLineMaterial ( )
void QweakSimBeamLine::SetBeamLineCenterPositionInZ ( G4double  zPos)

Definition at line 648 of file QweakSimBeamLine.cc.

References beamline_ZPos.

Referenced by QweakSimBeamLineMessenger::SetNewValue().

649 {
650  G4cout << G4endl << "###### Calling QweakSimBeamLine::SetBeamLineCenterPositionInZ() " << G4endl << G4endl;
651 
652  beamline_ZPos = zPos;
653  //BeamLineContainer_Physical->SetTranslation(G4ThreeVector(0.,0., zPos));
654 }

+ Here is the caller graph for this function:

void QweakSimBeamLine::SetBeamLineMaterial ( G4String  materialName)

Definition at line 605 of file QweakSimBeamLine.cc.

References DS_18inch_Pipe1_Logical, DS_18inch_Pipe2_Logical, DS_24inch_Pipe_Flange_Logical, DS_24inch_Pipe_Logical, DS_R1_Flange_Logical, DS_R1_Pipe_Logical, DS_R2_Pipe_Logical, DS_R3_Pipe_Logical, R1_Flange_Logical, R1_Pipe_Logical, R2_Flange_Logical, R2_Pipe_Logical, R2_RotatorPipe_Logical, R3_Flange_Logical, R3_Pipe_Logical, R3_US_Wall_Pipe_Logical, and R3_Wall_Pipe_Logical.

Referenced by QweakSimBeamLineMessenger::SetNewValue().

606 {
607  // search the material by its name
608  G4Material* pttoMaterial = G4Material::GetMaterial(materialName);
609  if (pttoMaterial)
610  {
611  G4cout << "==== Changing Beam Pipe Material: Looking up Material " << G4endl;
612  R1_Pipe_Logical->SetMaterial(pttoMaterial);
613  R1_Flange_Logical->SetMaterial(pttoMaterial);
614  DS_R1_Pipe_Logical->SetMaterial(pttoMaterial);
615  DS_R1_Flange_Logical->SetMaterial(pttoMaterial);
616  //DS_R1_Bellow_Logical->SetMaterial(pttoMaterial);
617  R2_Pipe_Logical->SetMaterial(pttoMaterial);
618  R2_Flange_Logical->SetMaterial(pttoMaterial);
619  DS_R2_Pipe_Logical->SetMaterial(pttoMaterial);
620  R2_RotatorPipe_Logical->SetMaterial(pttoMaterial);
621  R3_Pipe_Logical->SetMaterial(pttoMaterial);
622  R3_Flange_Logical->SetMaterial(pttoMaterial);
623  DS_R3_Pipe_Logical->SetMaterial(pttoMaterial);
624  R3_US_Wall_Pipe_Logical->SetMaterial(pttoMaterial);
625  //Sub_LeadBox_in_Wall_Logical->SetMaterial(pttoMaterial);
626  R3_Wall_Pipe_Logical->SetMaterial(pttoMaterial);
627  //Sub_LeadBox_Extent_Logical->SetMaterial(pttoMaterial);
628  DS_18inch_Pipe1_Logical->SetMaterial(pttoMaterial);
629  DS_18inch_Pipe2_Logical->SetMaterial(pttoMaterial);
630  DS_24inch_Pipe_Logical->SetMaterial(pttoMaterial);
631  DS_24inch_Pipe_Flange_Logical->SetMaterial(pttoMaterial);
632  G4cout << "==== Changing Beam Pipe Material: Now the BeamLine is made of " << materialName << G4endl;
633  }
634  else {
635  G4cerr << "==== ERROR: Changing Beam Pipe Material failed" << G4endl;
636  }
637 
638 }
G4LogicalVolume * R1_Pipe_Logical
G4LogicalVolume * R1_Flange_Logical
G4LogicalVolume * R2_Pipe_Logical
G4LogicalVolume * R3_Flange_Logical
G4LogicalVolume * DS_R1_Pipe_Logical
G4LogicalVolume * DS_24inch_Pipe_Logical
G4LogicalVolume * DS_18inch_Pipe1_Logical
G4LogicalVolume * DS_R3_Pipe_Logical
G4LogicalVolume * R2_RotatorPipe_Logical
G4LogicalVolume * R2_Flange_Logical
G4LogicalVolume * R3_Pipe_Logical
G4LogicalVolume * R3_US_Wall_Pipe_Logical
G4LogicalVolume * DS_R2_Pipe_Logical
G4LogicalVolume * DS_18inch_Pipe2_Logical
G4LogicalVolume * DS_R1_Flange_Logical
G4LogicalVolume * DS_24inch_Pipe_Flange_Logical
G4LogicalVolume * R3_Wall_Pipe_Logical

+ Here is the caller graph for this function:

Field Documentation

G4double QweakSimBeamLine::beamline_ZPos
private
QweakSimBeamLineMessenger* QweakSimBeamLine::beamlineMessenger
private

Definition at line 112 of file QweakSimBeamLine.hh.

Referenced by QweakSimBeamLine(), and ~QweakSimBeamLine().

G4Material* QweakSimBeamLine::BeamPipe_Material
private

Definition at line 86 of file QweakSimBeamLine.hh.

Referenced by ConstructComponent(), and QweakSimBeamLine().

G4LogicalVolume* QweakSimBeamLine::DS_18inch_Pipe1_Logical
private

Definition at line 105 of file QweakSimBeamLine.hh.

Referenced by ConstructComponent(), QweakSimBeamLine(), and SetBeamLineMaterial().

G4LogicalVolume* QweakSimBeamLine::DS_18inch_Pipe2_Logical
private

Definition at line 106 of file QweakSimBeamLine.hh.

Referenced by ConstructComponent(), QweakSimBeamLine(), and SetBeamLineMaterial().

G4LogicalVolume* QweakSimBeamLine::DS_24inch_Pipe_Flange_Logical
private

Definition at line 108 of file QweakSimBeamLine.hh.

Referenced by ConstructComponent(), QweakSimBeamLine(), and SetBeamLineMaterial().

G4LogicalVolume* QweakSimBeamLine::DS_24inch_Pipe_Logical
private

Definition at line 107 of file QweakSimBeamLine.hh.

Referenced by ConstructComponent(), QweakSimBeamLine(), and SetBeamLineMaterial().

G4LogicalVolume* QweakSimBeamLine::DS_R1_Bellow_Logical
private

Definition at line 93 of file QweakSimBeamLine.hh.

Referenced by ConstructComponent(), and QweakSimBeamLine().

G4LogicalVolume* QweakSimBeamLine::DS_R1_Flange_Logical
private

Definition at line 92 of file QweakSimBeamLine.hh.

Referenced by ConstructComponent(), QweakSimBeamLine(), and SetBeamLineMaterial().

G4LogicalVolume* QweakSimBeamLine::DS_R1_Pipe_Logical
private

Definition at line 91 of file QweakSimBeamLine.hh.

Referenced by ConstructComponent(), QweakSimBeamLine(), and SetBeamLineMaterial().

G4LogicalVolume* QweakSimBeamLine::DS_R2_Pipe_Logical
private

Definition at line 96 of file QweakSimBeamLine.hh.

Referenced by ConstructComponent(), QweakSimBeamLine(), and SetBeamLineMaterial().

G4LogicalVolume* QweakSimBeamLine::DS_R3_Pipe_Logical
private

Definition at line 100 of file QweakSimBeamLine.hh.

Referenced by ConstructComponent(), QweakSimBeamLine(), and SetBeamLineMaterial().

QweakSimUserInformation* QweakSimBeamLine::myUserInfo
private

Definition at line 113 of file QweakSimBeamLine.hh.

Referenced by QweakSimBeamLine().

QweakSimMaterial* QweakSimBeamLine::pMaterial
private

Definition at line 80 of file QweakSimBeamLine.hh.

Referenced by QweakSimBeamLine().

G4LogicalVolume* QweakSimBeamLine::R1_Flange_Logical
private

Definition at line 90 of file QweakSimBeamLine.hh.

Referenced by ConstructComponent(), QweakSimBeamLine(), and SetBeamLineMaterial().

G4LogicalVolume* QweakSimBeamLine::R1_Pipe_Logical
private

Definition at line 89 of file QweakSimBeamLine.hh.

Referenced by ConstructComponent(), QweakSimBeamLine(), and SetBeamLineMaterial().

G4LogicalVolume* QweakSimBeamLine::R2_Flange_Logical
private

Definition at line 95 of file QweakSimBeamLine.hh.

Referenced by ConstructComponent(), QweakSimBeamLine(), and SetBeamLineMaterial().

G4LogicalVolume* QweakSimBeamLine::R2_Pipe_Logical
private

Definition at line 94 of file QweakSimBeamLine.hh.

Referenced by ConstructComponent(), QweakSimBeamLine(), and SetBeamLineMaterial().

G4LogicalVolume* QweakSimBeamLine::R2_RotatorPipe_Logical
private

Definition at line 97 of file QweakSimBeamLine.hh.

Referenced by ConstructComponent(), QweakSimBeamLine(), and SetBeamLineMaterial().

G4LogicalVolume* QweakSimBeamLine::R3_Flange_Logical
private

Definition at line 99 of file QweakSimBeamLine.hh.

Referenced by ConstructComponent(), QweakSimBeamLine(), and SetBeamLineMaterial().

G4LogicalVolume* QweakSimBeamLine::R3_Pipe_Logical
private

Definition at line 98 of file QweakSimBeamLine.hh.

Referenced by ConstructComponent(), QweakSimBeamLine(), and SetBeamLineMaterial().

G4LogicalVolume* QweakSimBeamLine::R3_US_Wall_Pipe_Logical
private

Definition at line 101 of file QweakSimBeamLine.hh.

Referenced by ConstructComponent(), QweakSimBeamLine(), and SetBeamLineMaterial().

G4LogicalVolume* QweakSimBeamLine::R3_Wall_Pipe_Logical
private

Definition at line 103 of file QweakSimBeamLine.hh.

Referenced by ConstructComponent(), QweakSimBeamLine(), and SetBeamLineMaterial().

G4Material* QweakSimBeamLine::Shield_Material
private

Definition at line 87 of file QweakSimBeamLine.hh.

Referenced by ConstructComponent(), and QweakSimBeamLine().

G4LogicalVolume* QweakSimBeamLine::Sub_LeadBox_Extent_Logical
private

Definition at line 104 of file QweakSimBeamLine.hh.

Referenced by ConstructComponent(), and QweakSimBeamLine().

G4LogicalVolume* QweakSimBeamLine::Sub_LeadBox_in_Wall_Logical
private

Definition at line 102 of file QweakSimBeamLine.hh.

Referenced by ConstructComponent(), and QweakSimBeamLine().


The documentation for this class was generated from the following files: