QwGeant4
QweakSimTarget Class Reference

Definition of the Target. More...

#include <QweakSimTarget.hh>

+ Collaboration diagram for QweakSimTarget:

Public Member Functions

 QweakSimTarget (QweakSimUserInformation *myUI)
 
 ~QweakSimTarget ()
 
void ConstructComponent (G4VPhysicalVolume *)
 
void DestroyComponent ()
 
void SetTargetCenterPositionInZ (G4double)
 
G4double GetTargetCenterPositionInZ ()
 
void SetTargetMaterial (G4String)
 
G4String GetTargetMaterial ()
 
void SetTargetCellMaterial (G4String)
 
G4String GetTargetCellMaterial ()
 
void SetTargetEntranceWindowMaterial (G4String)
 
G4String GetTargetEntranceMaterial ()
 
void SetTargetExitWindowMaterial (G4String)
 
G4String GetTargetExitWindowMaterial ()
 
void SetTargetExitWindowNippleMaterial (G4String)
 
G4String GetTargetExitWindowNippleMaterial ()
 
G4LogicalVolume * getTargetLogicalVolume ()
 
G4VPhysicalVolume * getTargetPhysicalVolume ()
 
void SetTargetLength (G4double)
 
G4double GetTargetLength ()
 
void SetTarget (G4String)
 
void SetTargetEntranceWindowLength (G4double)
 
void SetTargetExitWindowLength (G4double)
 
void SetTargetExitWindowNippleLength (G4double)
 
G4double CalculateLuminosity (G4double mass, G4double density, G4double length, G4double pweight)
 

Private Member Functions

void CalculateTargetPositions ()
 
void ConstructTargetContainer ()
 
void ConstructScatteringChamberWindow ()
 
void ConstructTargetCell ()
 
void ConstructTargetMaterial ()
 
void ConstructTargetEntranceWindow ()
 
void ConstructTargetExitWindow ()
 
void ConstructTargetExitWindowNipple ()
 

Private Attributes

G4VPhysicalVolume * theMotherPV
 
QweakSimMaterialpMaterial
 
G4LogicalVolume * TargetCell_Logical
 
G4VPhysicalVolume * TargetCell_Physical
 
G4Material * TargetCell_Material
 
G4LogicalVolume * TargetEntranceWindow_Logical
 
G4VPhysicalVolume * TargetEntranceWindow_Physical
 
G4Material * TargetEntranceWindow_Material
 
G4LogicalVolume * TargetExitWindow_Logical
 
G4VPhysicalVolume * TargetExitWindow_Physical
 
G4Material * TargetExitWindow_Material
 
G4LogicalVolume * TargetExitWindowNipple_Logical
 
G4VPhysicalVolume * TargetExitWindowNipple_Physical
 
G4Material * TargetExitWindowNipple_Material
 
G4LogicalVolume * ScatteringChamberWindow_Logical
 
G4VPhysicalVolume * ScatteringChamberWindow_Physical
 
G4Material * ScatteringChamberWindow_Material
 
G4LogicalVolume * TargetMaterial_Logical
 
G4VPhysicalVolume * TargetMaterial_Physical
 
G4Material * TargetMaterial_Material
 
G4LogicalVolume * TargetContainer_Logical
 
G4VPhysicalVolume * TargetContainer_Physical
 
G4Material * TargetContainer_Material
 
G4double targetCellEntranceWindowThickness
 
G4double targetCellExitWindowThickness
 
G4double targetCellExitWindowNippleThickness
 
G4double targetCellWallThickness
 
G4double targetCellInnerLength
 
G4double targetCellOuterLength
 
G4double targetCellFrontRadiusMin
 
G4double targetCellFrontInnerRadiusMax
 
G4double targetCellFrontOuterRadiusMax
 
G4double targetCellBackRadiusMin
 
G4double targetCellBackInnerRadiusMax
 
G4double targetCellBackOuterRadiusMax
 
G4double targetCellExitWindowNippleRadius
 
G4double ScatteringChamberWindowRadius
 
G4double ScatteringChamberWindowThickness
 
G4double targetCellStartingPhi
 
G4double targetCellDeltaPhi
 
G4double targetZPos
 
QweakSimTargetMessengertargetMessenger
 
QweakSimUserInformationmyUserInfo
 
G4ThreeVector positionTarget
 
G4ThreeVector positionTargetEntranceWindow
 
G4ThreeVector positionTargetExitWindow
 
G4ThreeVector positionScatteringChamberWindow
 
G4VSensitiveDetector * TargetSD
 

Detailed Description

Definition of the Target.

Placeholder for a long explaination

Definition at line 65 of file QweakSimTarget.hh.

Constructor & Destructor Documentation

QweakSimTarget::QweakSimTarget ( QweakSimUserInformation myUI)

Definition at line 36 of file QweakSimTarget.cc.

References CalculateLuminosity(), QweakSimMaterial::GetInstance(), QweakSimMaterial::GetMaterial(), inch, mil, myUserInfo, pMaterial, ScatteringChamberWindow_Logical, ScatteringChamberWindow_Material, ScatteringChamberWindow_Physical, ScatteringChamberWindowRadius, ScatteringChamberWindowThickness, TargetCell_Logical, TargetCell_Material, TargetCell_Physical, targetCellBackInnerRadiusMax, targetCellBackOuterRadiusMax, targetCellBackRadiusMin, targetCellDeltaPhi, targetCellEntranceWindowThickness, targetCellExitWindowNippleRadius, targetCellExitWindowNippleThickness, targetCellExitWindowThickness, targetCellFrontInnerRadiusMax, targetCellFrontOuterRadiusMax, targetCellFrontRadiusMin, targetCellInnerLength, targetCellOuterLength, targetCellStartingPhi, targetCellWallThickness, QweakSimUserInformation::TargetCenterPositionZ, TargetContainer_Logical, TargetContainer_Material, TargetContainer_Physical, TargetEntranceWindow_Logical, TargetEntranceWindow_Material, TargetEntranceWindow_Physical, QweakSimUserInformation::TargetEntranceWindowThickness, TargetExitWindow_Logical, TargetExitWindow_Material, TargetExitWindow_Physical, TargetExitWindowNipple_Logical, TargetExitWindowNipple_Material, TargetExitWindowNipple_Physical, QweakSimUserInformation::TargetExitWindowNippleThickness, QweakSimUserInformation::TargetExitWindowThickness, QweakSimUserInformation::TargetLength, QweakSimUserInformation::TargetLuminosityDSALDummy2, QweakSimUserInformation::TargetLuminosityDSALDummy4, QweakSimUserInformation::TargetLuminosityDSALDummy4_Al, QweakSimUserInformation::TargetLuminosityDSALDummy4_Cr, QweakSimUserInformation::TargetLuminosityDSALDummy4_Cu, QweakSimUserInformation::TargetLuminosityDSALDummy4_Fe, QweakSimUserInformation::TargetLuminosityDSALDummy4_Mg, QweakSimUserInformation::TargetLuminosityDSALDummy4_Mn, QweakSimUserInformation::TargetLuminosityDSALDummy4_Si, QweakSimUserInformation::TargetLuminosityDSALDummy4_Ti, QweakSimUserInformation::TargetLuminosityDSALDummy4_Zn, QweakSimUserInformation::TargetLuminosityDSALDummy8, QweakSimUserInformation::TargetLuminosityDSALWindow, QweakSimUserInformation::TargetLuminosityDSCDummy, QweakSimUserInformation::TargetLuminosityLH2, QweakSimUserInformation::TargetLuminosityUSALDummy1, QweakSimUserInformation::TargetLuminosityUSALDummy2, QweakSimUserInformation::TargetLuminosityUSALDummy4, QweakSimUserInformation::TargetLuminosityUSALWindow, QweakSimUserInformation::TargetLuminosityUSCDummy, TargetMaterial_Logical, TargetMaterial_Material, TargetMaterial_Physical, targetMessenger, TargetSD, QweakSimUserInformation::TargetThicknessDSALDummy2, QweakSimUserInformation::TargetThicknessDSALDummy4, QweakSimUserInformation::TargetThicknessDSALDummy8, QweakSimUserInformation::TargetThicknessDSCDummy, QweakSimUserInformation::TargetThicknessUSALDummy1, QweakSimUserInformation::TargetThicknessUSALDummy2, QweakSimUserInformation::TargetThicknessUSALDummy4, QweakSimUserInformation::TargetThicknessUSCDummy, targetZPos, and theMotherPV.

37 {
38  G4cout << G4endl << "###### Calling QweakSimTarget::QweakTarget() " << G4endl << G4endl;
39 
40  myUserInfo = myUI;
41 
42  theMotherPV = NULL;
43 
47 
48  TargetCell_Logical = NULL;
49  TargetCell_Physical = NULL;
50  TargetCell_Material = NULL;
51 
55 
59 
63 
67 
71 
72 
74 
82 
83  // define target geometry values
84 
85  targetCellEntranceWindowThickness = 4.105*mil; //Avg of Run 1 & 2 (Greg's Target Elog 19)
86  targetCellExitWindowThickness = 29.06*mil; //Greg's Target Elog 19
87  targetCellExitWindowNippleThickness = 5.165*mil; //Avg of Run 1 & 2 (Greg's Target Elog 19)
89 
91 
95  + targetCellExitWindowThickness; // Full length of Target
96 
97  targetCellFrontRadiusMin = 0.0*cm;
100 
104 
105  targetCellStartingPhi = 0.0*deg;
106  targetCellDeltaPhi = 360*deg;
107 
109 
112 
116 
117  //Thicknesses taken from 2012 target survey (https://qweak.jlab.org/elog/Target/21)
124  myUserInfo->TargetThicknessUSCDummy = 0.9973*mm;
125  myUserInfo->TargetThicknessDSCDummy = 3.1876*mm;
126 
127  G4double densityLH2 = 0.0713 /(cm*cm*cm); // [g/cm^3]
128  G4double densityAL = 2.804 /(cm*cm*cm); // Avg density, Greg's Target Elog 21
129  G4double densityUSC = 1.70 /(cm*cm*cm); //
130  G4double densityDSC = 2.205 /(cm*cm*cm); //
131 
132  // Molar masses taken from PDG:
133  // J. Beringer et al. (Particle Data Group), Phys. Rev. D86, 010001 (2012).
134  G4double massLH2 = 1.00794; // [g/mol]
135  G4double massAL = 26.9815386; // [g/mol]
136  G4double massC = 12.0107; // [g/mol]
137  //Added mass for the elements in the alloy (PDG - 2015)
138  G4double massZn = 65.38; //[g/mol]
139  G4double massMg = 24.3050; //[g/mol]
140  G4double massCu = 63.546; //[g/mol]
141  G4double massCr = 51.9961; //[g/mol]
142  G4double massFe = 55.845; //[g/mol]
143  G4double massSi = 28.0855; //[g/mol]
144  G4double massMn = 54.938044; //[g/mol]
145  G4double massTi = 47.867; //[g/mol]
146 
147  //Weighting factors for calculation of luminosity of elements in the aluminum alloy
148  //Only implementing for DS 4% Al dummy target. Weighting factors take from Greg's
149  //Target ELOG 30 - K.Bartlett 2017/7/13
150  G4double weight_std = 1.0;//Standard weight for hydrogen and other solid targets.
151  G4double weight_Al = 0.8870;//Aluminum's % by weight[unitless]
152  G4double weight_Zn = 0.0630;//Zinc's % by weight[unitless]
153  G4double weight_Mg = 0.0270;//Magnesium's % by weight
154  G4double weight_Cu = 0.0180;//Copper's % by weight
155  G4double weight_Cr = 0.0021;//Chromium's % by weight
156  G4double weight_Fe = 0.0012;//Iron's % by weight
157  G4double weight_Si = 0.0010;//Silicon's % by weight
158  G4double weight_Mn = 0.0004;//Manganese's % by weight
159  G4double weight_Ti = 0.0003;//Titanium's % by weight
160 
161 
162  myUserInfo->TargetLuminosityLH2 = CalculateLuminosity(massLH2, densityLH2, myUserInfo->TargetLength, weight_std);
173 
174  //Calculate luminosities for elements in alloy for DS 4%.
184 
185 
186 
187  TargetSD = NULL;
188 
190 }
G4Material * ScatteringChamberWindow_Material
G4VPhysicalVolume * TargetEntranceWindow_Physical
G4double targetCellOuterLength
G4double ScatteringChamberWindowRadius
G4double CalculateLuminosity(G4double mass, G4double density, G4double length, G4double pweight)
static QweakSimMaterial * GetInstance()
G4double targetCellWallThickness
G4VSensitiveDetector * TargetSD
G4LogicalVolume * TargetExitWindow_Logical
Scans the input file for /Target/xyz commands.
G4double targetCellFrontInnerRadiusMax
G4VPhysicalVolume * TargetMaterial_Physical
QweakSimMaterial * pMaterial
G4Material * TargetMaterial_Material
G4Material * TargetExitWindow_Material
G4double targetCellExitWindowNippleRadius
G4double targetCellStartingPhi
G4LogicalVolume * TargetMaterial_Logical
G4VPhysicalVolume * TargetCell_Physical
G4VPhysicalVolume * TargetExitWindow_Physical
static const G4double inch
static const G4double mil
G4double targetCellFrontRadiusMin
G4double targetCellEntranceWindowThickness
G4Material * GetMaterial(G4String material)
G4double targetCellInnerLength
G4Material * TargetExitWindowNipple_Material
G4LogicalVolume * ScatteringChamberWindow_Logical
G4LogicalVolume * TargetContainer_Logical
G4VPhysicalVolume * ScatteringChamberWindow_Physical
G4LogicalVolume * TargetExitWindowNipple_Logical
G4double targetCellBackInnerRadiusMax
G4VPhysicalVolume * TargetContainer_Physical
G4LogicalVolume * TargetEntranceWindow_Logical
G4double targetCellDeltaPhi
G4double targetCellBackOuterRadiusMax
G4Material * TargetContainer_Material
G4double ScatteringChamberWindowThickness
G4double targetCellBackRadiusMin
G4Material * TargetEntranceWindow_Material
G4double targetCellExitWindowThickness
G4double targetCellExitWindowNippleThickness
G4Material * TargetCell_Material
G4VPhysicalVolume * theMotherPV
QweakSimTargetMessenger * targetMessenger
G4VPhysicalVolume * TargetExitWindowNipple_Physical
G4LogicalVolume * TargetCell_Logical
G4double targetCellFrontOuterRadiusMax
QweakSimUserInformation * myUserInfo

+ Here is the call graph for this function:

QweakSimTarget::~QweakSimTarget ( )

Definition at line 193 of file QweakSimTarget.cc.

References targetMessenger.

194 {
195  delete targetMessenger;
196 }
QweakSimTargetMessenger * targetMessenger

Member Function Documentation

G4double QweakSimTarget::CalculateLuminosity ( G4double  mass,
G4double  density,
G4double  length,
G4double  pweight 
)

Definition at line 923 of file QweakSimTarget.cc.

Referenced by QweakSimTarget().

924 {
925  G4double Lum = 0.0; // Luminosity
926  G4double N_b = 6.241e12; // [Hz/uA] # of particles in the beam (from definition of ampere)
927  G4double N_A = 6.02214129e23; // Avagadro's number
928 
929  Lum = N_b*length*pweight*density*N_A/mass; // units [Hz/(uA*mm^2)]
930 
931  Lum *= 1e-31; // Conversion from [Hz/(uA*mm^2)] -> [kHz/(uA*ub)]
932 
933  return Lum;
934 }

+ Here is the caller graph for this function:

void QweakSimTarget::CalculateTargetPositions ( )
private

Definition at line 266 of file QweakSimTarget.cc.

References positionScatteringChamberWindow, positionTarget, positionTargetEntranceWindow, positionTargetExitWindow, targetCellEntranceWindowThickness, targetCellExitWindowThickness, targetCellInnerLength, targetCellOuterLength, and TargetContainer_Physical.

Referenced by ConstructComponent(), SetTargetEntranceWindowLength(), SetTargetExitWindowLength(), SetTargetExitWindowNippleLength(), and SetTargetLength().

266  {
267 
268  // define Target position
269  positionTarget = G4ThreeVector(0,0,0);
272  positionScatteringChamberWindow = G4ThreeVector(0,0,0.5*targetCellInnerLength + 0.5*targetCellExitWindowThickness + 45.0*cm); // Peiqing: should be located at z=-583.41 cm
273 
276  + targetCellExitWindowThickness; // Full length of Target
277 
279  G4cout << "TargetContainer_Physical = " << TargetContainer_Physical->GetTranslation()/cm << " cm" << G4endl;
280  }
281 
282  G4cout << "targetCellEntranceWindowThickness: " << targetCellEntranceWindowThickness/cm << " cm" << G4endl;
283  G4cout << "positionTargetEntranceWindow: " << positionTargetEntranceWindow/cm << " cm" << G4endl;
284  G4cout << "positionTarget: " << positionTarget/cm << " cm" << G4endl;
285  G4cout << "positionTargetExitWindow: " << positionTargetExitWindow/cm << " cm" << G4endl;
286  G4cout << "targetCellExitWindowThickness: " << targetCellExitWindowThickness/cm << " cm" << G4endl;
287  G4cout << "positionScatteringChamberWindow: " << positionScatteringChamberWindow/cm << " cm" << G4endl;
288  G4cout << "targetCellInnerLength = " << targetCellInnerLength/cm << " cm" << G4endl;
289  G4cout << "targetCellOuterLength = " << targetCellOuterLength/cm << " cm" << G4endl;
290 }
G4double targetCellOuterLength
G4ThreeVector positionTarget
G4ThreeVector positionScatteringChamberWindow
G4double targetCellEntranceWindowThickness
G4double targetCellInnerLength
G4ThreeVector positionTargetExitWindow
G4VPhysicalVolume * TargetContainer_Physical
G4ThreeVector positionTargetEntranceWindow
G4double targetCellExitWindowThickness

+ Here is the caller graph for this function:

void QweakSimTarget::ConstructComponent ( G4VPhysicalVolume *  MotherVolume)

Definition at line 200 of file QweakSimTarget.cc.

References CalculateTargetPositions(), ConstructScatteringChamberWindow(), ConstructTargetCell(), ConstructTargetContainer(), ConstructTargetEntranceWindow(), ConstructTargetExitWindow(), ConstructTargetExitWindowNipple(), ConstructTargetMaterial(), ScatteringChamberWindow_Logical, TargetCell_Logical, TargetContainer_Logical, TargetEntranceWindow_Logical, TargetExitWindow_Logical, TargetMaterial_Logical, TargetSD, and theMotherPV.

201 {
202  G4cout << G4endl << "###### Calling QweakSimTarget::ConstructComponent() " << G4endl << G4endl;
203 
204  theMotherPV = MotherVolume;
205 
207 
208  ConstructTargetContainer(); // scattering chamber
210 
211  ConstructTargetCell(); // Al cell without end caps
212  ConstructTargetMaterial(); // LH2 for production target
213 
217 
218 //--------------------------------------
219 
220  G4cout << G4endl << "###### QweakSimTarget: Setting Attributes " << G4endl << G4endl;
221 
222  G4Colour blue (0.,0.,1.);
223  G4Colour red (1.,0.,0.);
224 
225  G4VisAttributes* TargetContainer_VisAtt = new G4VisAttributes(red);
226  TargetContainer_VisAtt -> SetVisibility(false);
227  //TargetContainer_VisAtt -> SetForceWireframe(true);
228  TargetContainer_Logical -> SetVisAttributes(TargetContainer_VisAtt);
229 
230  G4VisAttributes* TargetCell_VisAtt = new G4VisAttributes(blue);
231  TargetCell_VisAtt -> SetVisibility(true);
232  //TargetCell_VisAtt -> SetForceWireframe(true);
233  TargetCell_Logical -> SetVisAttributes(TargetCell_VisAtt);
234 
235  G4VisAttributes* TargetWindow_VisAtt = new G4VisAttributes(blue);
236  TargetWindow_VisAtt -> SetVisibility(true);
237  //TargetWindow_VisAtt -> SetForceWireframe(true);
238  TargetEntranceWindow_Logical -> SetVisAttributes(TargetWindow_VisAtt);
239  TargetExitWindow_Logical -> SetVisAttributes(TargetWindow_VisAtt);
240 
241  G4VisAttributes* ScatteringChamberWindow_VisAtt = new G4VisAttributes(red);
242  ScatteringChamberWindow_VisAtt -> SetVisibility(true);
243  //ScatteringChamberWindow_VisAtt -> SetForceWireframe(true);
244  ScatteringChamberWindow_Logical -> SetVisAttributes(ScatteringChamberWindow_VisAtt);
245 
246  G4VisAttributes* TargetMaterial_VisAtt = new G4VisAttributes(red);
247  TargetMaterial_VisAtt -> SetVisibility(true);
248  //TargetVisAtt -> SetForceWireframe(true);
249  TargetMaterial_Logical -> SetVisAttributes(TargetMaterial_VisAtt);
250 
251  ///////////////////////////////////////////
252  // Define Sensitive Detectors to Target //
253  ///////////////////////////////////////////
254 
255  G4SDManager* SDman = G4SDManager::GetSDMpointer();
256 
257  TargetSD = new QweakSimTarget_DetectorSD("TargetSD");
258  SDman->AddNewDetector(TargetSD);
259  ScatteringChamberWindow_Logical->SetSensitiveDetector(TargetSD);
260 
261 
262  G4cout << G4endl << "###### Leaving QweakSimTarget::ConstructComponent() " << G4endl << G4endl;
263 
264 } // end of QweakSimTarget::ConstructComponent()
G4VSensitiveDetector * TargetSD
void ConstructTargetExitWindow()
void ConstructTargetContainer()
G4LogicalVolume * TargetExitWindow_Logical
void ConstructTargetExitWindowNipple()
void ConstructTargetMaterial()
void CalculateTargetPositions()
G4LogicalVolume * TargetMaterial_Logical
void ConstructTargetEntranceWindow()
void ConstructTargetCell()
G4LogicalVolume * ScatteringChamberWindow_Logical
G4LogicalVolume * TargetContainer_Logical
void ConstructScatteringChamberWindow()
G4LogicalVolume * TargetEntranceWindow_Logical
G4VPhysicalVolume * theMotherPV
G4LogicalVolume * TargetCell_Logical

+ Here is the call graph for this function:

void QweakSimTarget::ConstructScatteringChamberWindow ( )
private

Definition at line 509 of file QweakSimTarget.cc.

References positionScatteringChamberWindow, pSurfChk, ScatteringChamberWindow_Logical, ScatteringChamberWindow_Material, ScatteringChamberWindow_Physical, ScatteringChamberWindowRadius, ScatteringChamberWindowThickness, targetCellDeltaPhi, targetCellStartingPhi, and TargetContainer_Physical.

Referenced by ConstructComponent().

510 {
511  G4cout << G4endl << "###### Calling QweakSimTarget::ConstructScatteringChamberWindow() " << G4endl << G4endl;
512 
513  //*********************** Define Target scattering chamber Vacuum Window ************************
514  G4cout << G4endl << "###### QweakSimTarget: Define ScatteringChamberWindow_Solid " << G4endl << G4endl;
515  G4Tubs* ScatteringChamberWindow_Solid = new G4Tubs("ScatteringChamberWindow_Sol",
516  0.,
521 
522  // define ScatteringChamberWindow logical volume
523  G4cout << G4endl << "###### QweakSimTarget: Define ScatteringChamberWindow_Logical " << G4endl << G4endl;
524 
525  ScatteringChamberWindow_Logical = new G4LogicalVolume(ScatteringChamberWindow_Solid,
527  "ScatteringChamberWindow_Log",
528  0,0,0);
529 
530  // define Target vacuum window physical volume
531  G4cout << G4endl << "###### QweakSimTarget: Define ScatteringChamberWindow_Physical " << G4endl << G4endl;
532  ScatteringChamberWindow_Physical = new G4PVPlacement(0,
534  "ScatteringChamberWindow",
536  TargetContainer_Physical, //MotherVolume,
537  false,
538  0,
539  pSurfChk);
540  G4cout << G4endl << "###### Leaving QweakSimTarget::ConstructScatteringChamberWindow() " << G4endl << G4endl;
541 
542 //*************************End of Target Vacuum Window *****************************
543 }
G4Material * ScatteringChamberWindow_Material
G4double ScatteringChamberWindowRadius
static const G4bool pSurfChk
G4double targetCellStartingPhi
G4ThreeVector positionScatteringChamberWindow
G4LogicalVolume * ScatteringChamberWindow_Logical
G4VPhysicalVolume * ScatteringChamberWindow_Physical
G4VPhysicalVolume * TargetContainer_Physical
G4double targetCellDeltaPhi
G4double ScatteringChamberWindowThickness

+ Here is the caller graph for this function:

void QweakSimTarget::ConstructTargetCell ( )
private

Definition at line 321 of file QweakSimTarget.cc.

References positionTarget, pSurfChk, TargetCell_Logical, TargetCell_Material, TargetCell_Physical, targetCellBackInnerRadiusMax, targetCellBackOuterRadiusMax, targetCellDeltaPhi, targetCellFrontInnerRadiusMax, targetCellFrontOuterRadiusMax, targetCellInnerLength, targetCellOuterLength, targetCellStartingPhi, TargetContainer_Logical, and TargetContainer_Physical.

Referenced by ConstructComponent(), and SetTargetLength().

322 {
323  G4cout << G4endl << "###### Calling QweakSimTarget::ConstructTargetCell() " << G4endl << G4endl;
324 
327 
328  delete TargetCell_Physical;
329  TargetCell_Physical = NULL;
330 
331  // define target solid volume
332  G4cout << G4endl << "###### QweakSimTarget: Define TargetCell_Solid " << G4endl << G4endl;
333  // front and back inner radius are at targetCellInnerLength, not targetCellOuterLength
334  G4double correction
338  G4Cons* TargetCell_Solid = new G4Cons("QweakTargetCell_Sol",
339  targetCellFrontInnerRadiusMax+correction,// G4double pRmin1,
340  targetCellFrontOuterRadiusMax, // G4double pRmax1,
341  targetCellBackInnerRadiusMax+correction, //G4double pRmin2,
342  targetCellBackOuterRadiusMax, //G4double pRmax2,
343  0.5*targetCellOuterLength, //G4double pDz,
344  targetCellStartingPhi, //G4double pSPhi,
345  targetCellDeltaPhi); //G4double pDPhi)
346 
347  // define Target logical volume
348  G4cout << G4endl << "###### QweakSimTarget: Define TargetCell_Logical " << G4endl << G4endl;
349 
350  TargetCell_Logical = new G4LogicalVolume(TargetCell_Solid,
352  "QweakTargetCell_Log",
353  0,0,0);
354  // define Target physical volume
355  G4cout << G4endl << "###### QweakSimTarget: Define TargetCell_Physical " << G4endl << G4endl;
356  TargetCell_Physical = new G4PVPlacement(0,
358  "QweakTargetCell",
360  TargetContainer_Physical, //MotherVolume,
361  false,
362  0,
363  pSurfChk);
364  G4cout << G4endl << "###### Leaving QweakSimTarget::ConstructTargetCell() " << G4endl << G4endl;
365 
366 }
G4double targetCellOuterLength
static const G4bool pSurfChk
G4double targetCellFrontInnerRadiusMax
G4ThreeVector positionTarget
G4double targetCellStartingPhi
G4VPhysicalVolume * TargetCell_Physical
G4double targetCellInnerLength
G4LogicalVolume * TargetContainer_Logical
G4double targetCellBackInnerRadiusMax
G4VPhysicalVolume * TargetContainer_Physical
G4double targetCellDeltaPhi
G4double targetCellBackOuterRadiusMax
G4Material * TargetCell_Material
G4LogicalVolume * TargetCell_Logical
G4double targetCellFrontOuterRadiusMax

+ Here is the caller graph for this function:

void QweakSimTarget::ConstructTargetContainer ( )
private

Definition at line 293 of file QweakSimTarget.cc.

References positionTarget, pSurfChk, ScatteringChamberWindowRadius, targetCellDeltaPhi, targetCellOuterLength, targetCellStartingPhi, TargetContainer_Logical, TargetContainer_Material, TargetContainer_Physical, and theMotherPV.

Referenced by ConstructComponent().

294 {
295  G4cout << G4endl << "###### Calling QweakSimTarget::ConstructTargetContainer() " << G4endl << G4endl;
296 
297  G4Tubs* TargetContainer_Solid = new G4Tubs("TargetContainer_Sol",
298  0, //targetCellRadiusMin, jpan@nuclear.uwinnipeg.ca
300  0.5*targetCellOuterLength+50.0*cm,
303 
304  TargetContainer_Logical = new G4LogicalVolume(TargetContainer_Solid,
306  "TargetContainer_Log",
307  0,0,0);
308 
309  TargetContainer_Physical = new G4PVPlacement(0,
311  "TargetContainer",
313  theMotherPV,
314  false,
315  0,
316  pSurfChk);
317  G4cout << G4endl << "###### Leaving QweakSimTarget::ConstructTargetContainer() " << G4endl << G4endl;
318 
319 }
G4double targetCellOuterLength
G4double ScatteringChamberWindowRadius
static const G4bool pSurfChk
G4ThreeVector positionTarget
G4double targetCellStartingPhi
G4LogicalVolume * TargetContainer_Logical
G4VPhysicalVolume * TargetContainer_Physical
G4double targetCellDeltaPhi
G4Material * TargetContainer_Material
G4VPhysicalVolume * theMotherPV

+ Here is the caller graph for this function:

void QweakSimTarget::ConstructTargetEntranceWindow ( )
private

Definition at line 368 of file QweakSimTarget.cc.

References positionTargetEntranceWindow, pSurfChk, targetCellDeltaPhi, targetCellEntranceWindowThickness, targetCellFrontInnerRadiusMax, targetCellFrontRadiusMin, targetCellStartingPhi, TargetContainer_Logical, TargetContainer_Physical, TargetEntranceWindow_Logical, TargetEntranceWindow_Material, and TargetEntranceWindow_Physical.

Referenced by ConstructComponent(), SetTargetEntranceWindowLength(), and SetTargetLength().

369 {
370  G4cout << G4endl << "###### Calling QweakSimTarget::ConstructTargetEntranceWindow() " << G4endl << G4endl;
371 
374 
377 
378 //--------------------------------------
379  // define target window solid volume (front, upstream)
380  G4cout << G4endl << "###### QweakSimTarget: Define TargetEntranceWindow_Solid " << G4endl << G4endl;
381 
382  G4Tubs* TargetEntranceWindow_Solid = new G4Tubs("TargetEntranceWindow_Sol",
388 
389  // define Target window logical volume (front, upstream)
390  G4cout << G4endl << "###### QweakSimTarget: Define TargetEntranceWindow_Logical " << G4endl << G4endl;
391 
392  TargetEntranceWindow_Logical = new G4LogicalVolume(TargetEntranceWindow_Solid,
394  "QweakTargetEntranceWindow_Log",
395  0,0,0);
396 
397  //G4double MaxStepInEntranceWindow = 0.1*targetCellEntranceWindowThickness;
398  //TargetEntranceWindow_Logical->SetUserLimits(new G4UserLimits(MaxStepInEntranceWindow));
399 
400  // define Target window physical volume (front, upstream)
401  G4cout << G4endl << "###### QweakSimTarget: Define TargetEntranceWindow_Physical " << G4endl << G4endl;
402  TargetEntranceWindow_Physical = new G4PVPlacement(0,
404  "QweakTargetEntranceWindow",
406  TargetContainer_Physical, //MotherVolume,
407  false,
408  0,
409  pSurfChk);
410  G4cout << G4endl << "###### Leaving QweakSimTarget::ConstructTargetEntranceWindow() " << G4endl << G4endl;
411 
412 }
G4VPhysicalVolume * TargetEntranceWindow_Physical
static const G4bool pSurfChk
G4double targetCellFrontInnerRadiusMax
G4double targetCellStartingPhi
G4double targetCellFrontRadiusMin
G4double targetCellEntranceWindowThickness
G4LogicalVolume * TargetContainer_Logical
G4VPhysicalVolume * TargetContainer_Physical
G4LogicalVolume * TargetEntranceWindow_Logical
G4double targetCellDeltaPhi
G4ThreeVector positionTargetEntranceWindow
G4Material * TargetEntranceWindow_Material

+ Here is the caller graph for this function:

void QweakSimTarget::ConstructTargetExitWindow ( )
private

Definition at line 414 of file QweakSimTarget.cc.

References positionTargetExitWindow, pSurfChk, targetCellBackInnerRadiusMax, targetCellBackRadiusMin, targetCellDeltaPhi, targetCellExitWindowThickness, targetCellStartingPhi, TargetContainer_Logical, TargetContainer_Physical, TargetExitWindow_Logical, TargetExitWindow_Material, and TargetExitWindow_Physical.

Referenced by ConstructComponent(), SetTargetExitWindowLength(), and SetTargetLength().

415 {
416  G4cout << G4endl << "###### Calling QweakSimTarget::ConstructTargetExitWindow() " << G4endl << G4endl;
417 
420 
423 
424  // define target window solid volume (back, downstream)
425  G4cout << G4endl << "###### QweakSimTarget: Define TargetExitWindow_Solid " << G4endl << G4endl;
426 
427  G4Tubs* TargetExitWindow_Solid = new G4Tubs("TargetExitWindow_Sol",
433 
434  // define Target window logical volume (back, downstream)
435  G4cout << G4endl << "###### QweakSimTarget: Define TargetExitWindow_Logical " << G4endl << G4endl;
436 
437  TargetExitWindow_Logical = new G4LogicalVolume(TargetExitWindow_Solid,
439  "QweakTargetExitWindow_Log",
440  0,0,0);
441 
442  // Set max step size for a certain volume, see
443 // http://geant4.web.cern.ch/geant4/G4UsersDocuments/UsersGuides/
444 // ForApplicationDeveloper/html/TrackingAndPhysics/thresholdVScut.html
445  //G4double MaxStepInExitWindow = 0.1*targetCellExitWindowThickness;
446  //TargetExitWindow_Logical->SetUserLimits(new G4UserLimits(MaxStepInExitWindow));
447 
448 
449  // define Target window physical volume (back, downstream)
450  G4cout << G4endl << "###### QweakSimTarget: Define TargetExitWindow_Physical " << G4endl << G4endl;
451  TargetExitWindow_Physical = new G4PVPlacement(0,
453  "QweakTargetExitWindow",
455  TargetContainer_Physical, //MotherVolume,
456  false,
457  0,
458  pSurfChk);
459  G4cout << G4endl << "###### Leaving QweakSimTarget::ConstructTargetExitWindow() " << G4endl << G4endl;
460 
461 }
G4LogicalVolume * TargetExitWindow_Logical
static const G4bool pSurfChk
G4Material * TargetExitWindow_Material
G4double targetCellStartingPhi
G4VPhysicalVolume * TargetExitWindow_Physical
G4ThreeVector positionTargetExitWindow
G4LogicalVolume * TargetContainer_Logical
G4double targetCellBackInnerRadiusMax
G4VPhysicalVolume * TargetContainer_Physical
G4double targetCellDeltaPhi
G4double targetCellBackRadiusMin
G4double targetCellExitWindowThickness

+ Here is the caller graph for this function:

void QweakSimTarget::ConstructTargetExitWindowNipple ( )
private

Definition at line 463 of file QweakSimTarget.cc.

References positionTargetExitWindow, pSurfChk, targetCellDeltaPhi, targetCellExitWindowNippleRadius, targetCellExitWindowNippleThickness, targetCellStartingPhi, TargetContainer_Logical, TargetContainer_Physical, TargetExitWindowNipple_Logical, TargetExitWindowNipple_Material, and TargetExitWindowNipple_Physical.

Referenced by ConstructComponent(), SetTargetExitWindowNippleLength(), and SetTargetLength().

464 {
465  G4cout << G4endl << "###### Calling QweakSimTarget::ConstructTargetExitWindowNipple() " << G4endl << G4endl;
466 
469 
472 
473  // define target window Nipple solid volume (back, downstream)
474  G4cout << G4endl << "###### QweakSimTarget: Define TargetExitWindowNipple_Solid " << G4endl << G4endl;
475 
476  G4Tubs* TargetExitWindowNipple_Solid = new G4Tubs("TargetExitWindowNipple_Sol",
477  0.,
479  // targetCellExitWindowNippleRadius-0.00001*mm,
483 
484  // define Target window logical volume (back, downstream)
485  G4cout << G4endl << "###### QweakSimTarget: Define TargetExitWindowNipple_Logical " << G4endl << G4endl;
486 
487  TargetExitWindowNipple_Logical = new G4LogicalVolume(TargetExitWindowNipple_Solid,
489  "QweakTargetExitWindowNipple_Log",
490  0,0,0);
491 
492  //G4double MaxStepInExitWindowNipple = 0.1*targetCellExitWindowNippleThickness;
493  //TargetExitWindowNipple_Logical->SetUserLimits(new G4UserLimits(MaxStepInExitWindowNipple));
494 
495  // define Target window nipple physical volume (back, downstream)
496  G4cout << G4endl << "###### QweakSimTarget: Define TargetExitWindowNipple_Physical " << G4endl << G4endl;
497  TargetExitWindowNipple_Physical = new G4PVPlacement(0,
499  "QweakTargetExitWindowNipple_Physical",
501  TargetContainer_Physical, //MotherVolume,
502  false,
503  0,
504  pSurfChk);
505  G4cout << G4endl << "###### Leaving QweakSimTarget::ConstructTargetExitWindowNipple() " << G4endl << G4endl;
506 
507 }
static const G4bool pSurfChk
G4double targetCellExitWindowNippleRadius
G4double targetCellStartingPhi
G4Material * TargetExitWindowNipple_Material
G4ThreeVector positionTargetExitWindow
G4LogicalVolume * TargetContainer_Logical
G4LogicalVolume * TargetExitWindowNipple_Logical
G4VPhysicalVolume * TargetContainer_Physical
G4double targetCellDeltaPhi
G4double targetCellExitWindowNippleThickness
G4VPhysicalVolume * TargetExitWindowNipple_Physical

+ Here is the caller graph for this function:

void QweakSimTarget::ConstructTargetMaterial ( )
private

Definition at line 547 of file QweakSimTarget.cc.

References positionTarget, pSurfChk, targetCellBackInnerRadiusMax, targetCellDeltaPhi, targetCellFrontInnerRadiusMax, targetCellInnerLength, targetCellStartingPhi, TargetContainer_Logical, TargetContainer_Physical, TargetMaterial_Logical, TargetMaterial_Material, and TargetMaterial_Physical.

Referenced by ConstructComponent(), and SetTargetLength().

548 {
549  // define target material solid volume
550  G4cout << G4endl << "###### Calling QweakSimTarget::ConstructTargetMaterial() " << G4endl << G4endl;
551 
554 
557 
558  G4Cons* TargetMaterial_Solid = new G4Cons("QweakTargetMaterial_Sol",
559  0., //targetCellFrontRadiusMin, // G4double pRmin1,
560  targetCellFrontInnerRadiusMax, // G4double pRmax1,
561  0., //targetCellBackRadiusMin, //G4double pRmin2,
562  targetCellBackInnerRadiusMax, //G4double pRmax2,
563  0.5*targetCellInnerLength, //G4double pDz,
564  targetCellStartingPhi, //G4double pSPhi,
565  targetCellDeltaPhi); //G4double pDPhi)
566 
567 
568  // define Target logical volume
569  G4cout << G4endl << "###### QweakSimTarget: Define Target_Logical " << G4endl << G4endl;
570 
571  TargetMaterial_Logical = new G4LogicalVolume(TargetMaterial_Solid,
573  "QweakTargetMaterial_Log",
574  0,0,0);
575 
576  // set max step size for LH2 target
577  //G4double MaxStepInTarget = 0.05*targetCellInnerLength; //step size < 20% of target length
578  //TargetMaterial_Logical->SetUserLimits(new G4UserLimits(MaxStepInTarget));
579 
580  // define Target material physical volume
581  G4cout << G4endl << "###### QweakSimTarget: Define TargetMaterial_Physical " << G4endl << G4endl;
582  TargetMaterial_Physical = new G4PVPlacement(0,
584  "QweakTargetMaterial",
586  TargetContainer_Physical, //MotherVolume,
587  false,
588  0,
589  pSurfChk);
590  G4cout << G4endl << "###### Leaving QweakSimTarget::ConstructTargetMaterial() " << G4endl << G4endl;
591 
592 }
static const G4bool pSurfChk
G4double targetCellFrontInnerRadiusMax
G4ThreeVector positionTarget
G4VPhysicalVolume * TargetMaterial_Physical
G4Material * TargetMaterial_Material
G4double targetCellStartingPhi
G4LogicalVolume * TargetMaterial_Logical
G4double targetCellInnerLength
G4LogicalVolume * TargetContainer_Logical
G4double targetCellBackInnerRadiusMax
G4VPhysicalVolume * TargetContainer_Physical
G4double targetCellDeltaPhi

+ Here is the caller graph for this function:

void QweakSimTarget::DestroyComponent ( )

Definition at line 824 of file QweakSimTarget.cc.

825 {
826 }
G4String QweakSimTarget::GetTargetCellMaterial ( )
G4double QweakSimTarget::GetTargetCenterPositionInZ ( )

Definition at line 840 of file QweakSimTarget.cc.

References targetZPos.

841 {
842  G4cout << G4endl << "###### Calling QweakSimTarget::GetTargetCenterPositionInZ() " << G4endl << G4endl;
843  G4cout << "==== Getting Target CenterPositionZ: Now the Target Center Position in Z is " << targetZPos/cm << " cm" << G4endl;
844  return targetZPos;
845 }
G4String QweakSimTarget::GetTargetEntranceMaterial ( )
G4String QweakSimTarget::GetTargetExitWindowMaterial ( )
G4String QweakSimTarget::GetTargetExitWindowNippleMaterial ( )
G4double QweakSimTarget::GetTargetLength ( )

Definition at line 870 of file QweakSimTarget.cc.

References targetCellInnerLength.

871 {
872  G4cout << G4endl << "###### Calling QweakSimTarget::GetTargetLength() " << G4endl << G4endl;
873  G4cout << "==== Getting Target Length: Now the Target Length is " << targetCellInnerLength/cm << " cm" << G4endl;
874  return targetCellInnerLength;
875 }
G4double targetCellInnerLength
G4LogicalVolume* QweakSimTarget::getTargetLogicalVolume ( )
inline

Definition at line 92 of file QweakSimTarget.hh.

References TargetMaterial_Logical.

92 {return TargetMaterial_Logical;}
G4LogicalVolume * TargetMaterial_Logical
G4String QweakSimTarget::GetTargetMaterial ( )
G4VPhysicalVolume* QweakSimTarget::getTargetPhysicalVolume ( )
inline

Definition at line 93 of file QweakSimTarget.hh.

References TargetMaterial_Physical.

Referenced by QweakSimDetectorConstruction::ConstructQweak().

93 {return TargetMaterial_Physical;}
G4VPhysicalVolume * TargetMaterial_Physical

+ Here is the caller graph for this function:

void QweakSimTarget::SetTarget ( G4String  targName)

Definition at line 596 of file QweakSimTarget.cc.

References myUserInfo, SetTargetEntranceWindowLength(), SetTargetEntranceWindowMaterial(), SetTargetExitWindowLength(), SetTargetExitWindowMaterial(), SetTargetExitWindowNippleLength(), SetTargetExitWindowNippleMaterial(), SetTargetLength(), SetTargetMaterial(), targetCellInnerLength, QweakSimUserInformation::TargetThicknessDSALDummy2, QweakSimUserInformation::TargetThicknessDSALDummy4, QweakSimUserInformation::TargetThicknessDSALDummy8, QweakSimUserInformation::TargetThicknessDSCDummy, QweakSimUserInformation::TargetThicknessUSALDummy1, QweakSimUserInformation::TargetThicknessUSALDummy2, QweakSimUserInformation::TargetThicknessUSALDummy4, and QweakSimUserInformation::TargetThicknessUSCDummy.

Referenced by QweakSimTargetMessenger::SetNewValue().

597 { // NOTE: Can only make one call at a time. To change targets you must
598  // make a new run, otherwise the positions of the dummy target will be wrong.
599 
600  double zlength = 0.0; // Use for setting dummy position such that
601  // US face's are aligned (see Target ELOG 18).
602  // Technically the US targets are on average 1.39 mm upstream of the entrance window
603  // and the DS targets are on average 2.97 mm upstream of the exit window. For simplicity
604  // just align the US faces at -(+)17.173 cm for the US(DS) targets.
605  // K.E. Mesick 23/September/2014
606 
607  if (strcmp(targName,"LH2")==0)
608  {
609  G4cout << "==== Target is LH2 by default!! " << G4endl;
610  }
611  else if (strcmp(targName,"USAl1p")==0)
612  {
613  G4cout << "==== Changing Target to " << targName << G4endl;
617  SetTargetLength(zlength);
618  SetTargetMaterial("Vacuum");
619  SetTargetExitWindowMaterial("Vacuum");
621  G4cout << "==== Changing Target : Now the Target is " << targName << G4endl;
622  }
623  else if (strcmp(targName,"USAl2p")==0)
624  {
625  G4cout << "==== Changing Target to " << targName << G4endl;
629  SetTargetLength(zlength);
630  SetTargetMaterial("Vacuum");
631  SetTargetExitWindowMaterial("Vacuum");
633  G4cout << "==== Changing Target : Now the Target is " << targName << G4endl;
634  }
635  else if (strcmp(targName,"USAl4p")==0)
636  {
637  G4cout << "==== Changing Target to " << targName << G4endl;
641  SetTargetLength(zlength);
642  SetTargetMaterial("Vacuum");
643  SetTargetExitWindowMaterial("Vacuum");
645  G4cout << "==== Changing Target : Now the Target is " << targName << G4endl;
646  }
647  else if (strcmp(targName,"USC")==0)
648  {
649  G4cout << "==== Changing Target to " << targName << G4endl;
653  SetTargetLength(zlength);
654  SetTargetMaterial("Vacuum");
655  SetTargetExitWindowMaterial("Vacuum");
657  G4cout << "==== Changing Target : Now the Target is " << targName << G4endl;
658  }
659  else if (strcmp(targName,"DSAl2p")==0)
660  {
661  G4cout << "==== Changing Target to " << targName << G4endl;
662  zlength = targetCellInnerLength;
664  SetTargetLength(zlength);
665  SetTargetMaterial("Vacuum");
668  SetTargetExitWindowMaterial("Aluminum");
670  G4cout << "==== Changing Target : Now the Target is " << targName << G4endl;
671  }
672  else if (strcmp(targName,"DSAl4p")==0)
673  {
674  G4cout << "==== Changing Target to " << targName << G4endl;
675  zlength = targetCellInnerLength;
677  SetTargetLength(zlength);
678  SetTargetMaterial("Vacuum");
681  //Changed to AlAlloy for alloy simulations. K. Bartlett 2016/9/13
682  SetTargetExitWindowMaterial("AlAlloy");
684  G4cout << "==== Changing Target : Now the Target is " << targName << G4endl;
685  }
686  else if (strcmp(targName,"DSAl4pVacuum")==0)
687  {
688  //Special vacuum DS 4% aluminum target for radiative correction studies.
689  G4cout << "==== Changing Target to " << targName << G4endl;
690  zlength = targetCellInnerLength;
692  SetTargetLength(zlength);
693  SetTargetMaterial("Vacuum");
696  SetTargetExitWindowMaterial("Vacuum");
698  G4cout << "==== Changing Target : Now the Target is " << targName << G4endl;
699  }
700  else if (strcmp(targName,"DSAl8p")==0)
701  {
702  G4cout << "==== Changing Target to " << targName << G4endl;
703  zlength = targetCellInnerLength;
705  SetTargetLength(zlength);
706  SetTargetMaterial("Vacuum");
709  SetTargetExitWindowMaterial("Aluminum");
711  G4cout << "==== Changing Target : Now the Target is " << targName << G4endl;
712  }
713  else if (strcmp(targName,"DSC")==0)
714  {
715  G4cout << "==== Changing Target to " << targName << G4endl;
716  zlength = targetCellInnerLength;
718  SetTargetLength(zlength);
719  SetTargetMaterial("Vacuum");
722  SetTargetExitWindowMaterial("DSCarbon");
724  G4cout << "==== Changing Target : Now the Target is " << targName << G4endl;
725  }
726  else {
727  G4cerr << "==== ERROR: Changing Target failed" << G4endl;
728  }
729 
730 }
void SetTargetEntranceWindowLength(G4double)
void SetTargetEntranceWindowMaterial(G4String)
G4double targetCellInnerLength
void SetTargetExitWindowNippleLength(G4double)
void SetTargetMaterial(G4String)
void SetTargetExitWindowLength(G4double)
void SetTargetExitWindowMaterial(G4String)
void SetTargetLength(G4double)
void SetTargetExitWindowNippleMaterial(G4String)
QweakSimUserInformation * myUserInfo

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QweakSimTarget::SetTargetCellMaterial ( G4String  materialName)

Definition at line 752 of file QweakSimTarget.cc.

References TargetCell_Logical.

Referenced by QweakSimTargetMessenger::SetNewValue().

753 {
754  // search the material by its name
755  G4Material* pttoMaterial = G4Material::GetMaterial(materialName);
756  if (pttoMaterial)
757  {
758  G4cout << "==== Changing Target Cell Material: Looking up Material " << G4endl;
759  TargetCell_Logical->SetMaterial(pttoMaterial);
760  G4cout << "==== Changing Target Cell Material: Now the Target Cell is made of " << materialName << G4endl;
761  }
762  else {
763  G4cerr << "==== ERROR: Changing Target Cell Material failed" << G4endl;
764  }
765 
766 }
G4LogicalVolume * TargetCell_Logical

+ Here is the caller graph for this function:

void QweakSimTarget::SetTargetCenterPositionInZ ( G4double  zPos)

Definition at line 830 of file QweakSimTarget.cc.

References myUserInfo, QweakSimUserInformation::TargetCenterPositionZ, TargetContainer_Physical, and targetZPos.

Referenced by QweakSimTargetMessenger::SetNewValue().

831 {
832  G4cout << G4endl << "###### Calling QweakSimTarget::SetTargetCenterPositionInZ() " << G4endl << G4endl;
833 
834  targetZPos = zPos;
836  TargetContainer_Physical->SetTranslation(G4ThreeVector(0.,0.,zPos));
837  G4cout << "==== Changing Target CenterPositionZ: Now the Target Center Position in Z is " << zPos/cm << " cm" << G4endl;
838 }
G4VPhysicalVolume * TargetContainer_Physical
QweakSimUserInformation * myUserInfo

+ Here is the caller graph for this function:

void QweakSimTarget::SetTargetEntranceWindowLength ( G4double  len)

Definition at line 878 of file QweakSimTarget.cc.

References CalculateTargetPositions(), ConstructTargetEntranceWindow(), myUserInfo, targetCellEntranceWindowThickness, and QweakSimUserInformation::TargetEntranceWindowThickness.

Referenced by SetTarget().

879 {
880  G4cout << G4endl << "###### Calling QweakSimTarget::SetTargetEntranceWindowLength() " << G4endl << G4endl;
881 
884 
885  G4cout << "==== Changing Target Entrance Window Length: Now the Target Entrance Window Length is " << len/mm << " mm" << G4endl;
886 
889 
890  G4cout << G4endl << "###### Leaving QweakSimTarget::SetTargetEntranceWindowLength() " << G4endl << G4endl;
891 }
void CalculateTargetPositions()
void ConstructTargetEntranceWindow()
G4double targetCellEntranceWindowThickness
QweakSimUserInformation * myUserInfo

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QweakSimTarget::SetTargetEntranceWindowMaterial ( G4String  materialName)

Definition at line 770 of file QweakSimTarget.cc.

References TargetEntranceWindow_Logical.

Referenced by QweakSimTargetMessenger::SetNewValue(), and SetTarget().

771 {
772  // search the material by its name
773  G4Material* pttoMaterial = G4Material::GetMaterial(materialName);
774  if (pttoMaterial)
775  {
776  G4cout << "==== Changing Target Entrance Window Material: Looking up Material " << G4endl;
777  TargetEntranceWindow_Logical->SetMaterial(pttoMaterial);
778  G4cout << "==== Changing Target Entrance Window Material: Now the Target Entrance Window is made of " << materialName << G4endl;
779  }
780  else {
781  G4cerr << "==== ERROR: Changing Target Entrance Window Material failed" << G4endl;
782  }
783 
784 }
G4LogicalVolume * TargetEntranceWindow_Logical

+ Here is the caller graph for this function:

void QweakSimTarget::SetTargetExitWindowLength ( G4double  len)

Definition at line 893 of file QweakSimTarget.cc.

References CalculateTargetPositions(), ConstructTargetExitWindow(), myUserInfo, targetCellExitWindowThickness, and QweakSimUserInformation::TargetExitWindowThickness.

Referenced by SetTarget().

894 {
895  G4cout << G4endl << "###### Calling QweakSimTarget::SetTargetExitWindowLength() " << G4endl << G4endl;
896 
899  G4cout << "==== Changing Target Exit Window Length: Now the Target Exit Window Length is " << len/mm << " mm" << G4endl;
900 
903 
904  G4cout << G4endl << "###### Leaving QweakSimTarget::SetTargetExitWindowLength() " << G4endl << G4endl;
905 }
void ConstructTargetExitWindow()
void CalculateTargetPositions()
G4double targetCellExitWindowThickness
QweakSimUserInformation * myUserInfo

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QweakSimTarget::SetTargetExitWindowMaterial ( G4String  materialName)

Definition at line 788 of file QweakSimTarget.cc.

References TargetExitWindow_Logical.

Referenced by QweakSimTargetMessenger::SetNewValue(), and SetTarget().

789 {
790  // search the material by its name
791  G4Material* pttoMaterial = G4Material::GetMaterial(materialName);
792  if (pttoMaterial)
793  {
794  G4cout << "==== Changing Target Exit Window Material: Looking up Material " << G4endl;
795  TargetExitWindow_Logical->SetMaterial(pttoMaterial);
796  G4cout << "==== Changing Target Exit Window Material: Now the Target Exit Window is made of " << materialName << G4endl;
797  }
798  else {
799  G4cerr << "==== ERROR: Changing Target Exit Window Material failed" << G4endl;
800  }
801 
802 }
G4LogicalVolume * TargetExitWindow_Logical

+ Here is the caller graph for this function:

void QweakSimTarget::SetTargetExitWindowNippleLength ( G4double  len)

Definition at line 907 of file QweakSimTarget.cc.

References CalculateTargetPositions(), ConstructTargetExitWindowNipple(), myUserInfo, targetCellExitWindowNippleThickness, and QweakSimUserInformation::TargetExitWindowNippleThickness.

Referenced by SetTarget().

908 {
909  G4cout << G4endl << "###### Calling QweakSimTarget::SetTargetExitWindowNippleLength() " << G4endl << G4endl;
910 
913  G4cout << "==== Changing Target Exit Window Nipple Length: Now the Target Exit Window Nipple Length is " << len/mm << " mm" << G4endl;
914 
917 
918  G4cout << G4endl << "###### Leaving QweakSimTarget::SetTargetExitWindowNippleLength() " << G4endl << G4endl;
919 }
void ConstructTargetExitWindowNipple()
void CalculateTargetPositions()
G4double targetCellExitWindowNippleThickness
QweakSimUserInformation * myUserInfo

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QweakSimTarget::SetTargetExitWindowNippleMaterial ( G4String  materialName)

Definition at line 806 of file QweakSimTarget.cc.

References TargetExitWindowNipple_Logical.

Referenced by QweakSimTargetMessenger::SetNewValue(), and SetTarget().

807 {
808  // search the material by its name
809  G4Material* pttoMaterial = G4Material::GetMaterial(materialName);
810  if (pttoMaterial)
811  {
812  G4cout << "==== Changing Target Exit Window Nipple Material: Looking up Material " << G4endl;
813  TargetExitWindowNipple_Logical->SetMaterial(pttoMaterial);
814  G4cout << "==== Changing Target Exit Window Nipple Material: Now the Target Exit Window Nipple is made of " << materialName << G4endl;
815  }
816  else {
817  G4cerr << "==== ERROR: Changing Target Exit Window Nipple Material failed" << G4endl;
818  }
819 
820 }
G4LogicalVolume * TargetExitWindowNipple_Logical

+ Here is the caller graph for this function:

void QweakSimTarget::SetTargetLength ( G4double  len)

Definition at line 848 of file QweakSimTarget.cc.

References CalculateTargetPositions(), ConstructTargetCell(), ConstructTargetEntranceWindow(), ConstructTargetExitWindow(), ConstructTargetExitWindowNipple(), ConstructTargetMaterial(), myUserInfo, targetCellInnerLength, and QweakSimUserInformation::TargetLength.

Referenced by QweakSimTargetMessenger::SetNewValue(), and SetTarget().

849 {
850  G4cout << G4endl << "###### Calling QweakSimTarget::SetTargetLength() " << G4endl << G4endl;
851 
852  targetCellInnerLength = len;
853  myUserInfo->TargetLength = len;
854  G4cout << "==== Changing Target Length: Now the Target Length is " << len/cm << " cm" << G4endl;
855 
857 
858  // construct target cell/material of appropriate length
859  ConstructTargetCell(); // Al cell without end caps
860  ConstructTargetMaterial(); // LH2 for production target
861 
862  // construct the end caps at the correct location
866 
867  G4cout << G4endl << "###### Leaving QweakSimTarget::SetTargetLength() " << G4endl << G4endl;
868 }
void ConstructTargetExitWindow()
void ConstructTargetExitWindowNipple()
void ConstructTargetMaterial()
void CalculateTargetPositions()
void ConstructTargetEntranceWindow()
G4double targetCellInnerLength
void ConstructTargetCell()
QweakSimUserInformation * myUserInfo

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QweakSimTarget::SetTargetMaterial ( G4String  materialName)

Definition at line 734 of file QweakSimTarget.cc.

References TargetMaterial_Logical.

Referenced by QweakSimTargetMessenger::SetNewValue(), and SetTarget().

735 {
736  // search the material by its name
737  G4Material* pttoMaterial = G4Material::GetMaterial(materialName);
738  if (pttoMaterial)
739  {
740  G4cout << "==== Changing Target Material: Looking up Material " << G4endl;
741  TargetMaterial_Logical->SetMaterial(pttoMaterial);
742  G4cout << "==== Changing Target Material: Now the Target is made of " << materialName << G4endl;
743  }
744  else {
745  G4cerr << "==== ERROR: Changing Target Material failed" << G4endl;
746  }
747 
748 }
G4LogicalVolume * TargetMaterial_Logical

+ Here is the caller graph for this function:

Field Documentation

QweakSimMaterial* QweakSimTarget::pMaterial
private

Definition at line 109 of file QweakSimTarget.hh.

Referenced by QweakSimTarget().

G4ThreeVector QweakSimTarget::positionScatteringChamberWindow
private
G4ThreeVector QweakSimTarget::positionTarget
private
G4ThreeVector QweakSimTarget::positionTargetEntranceWindow
private

Definition at line 168 of file QweakSimTarget.hh.

Referenced by CalculateTargetPositions(), and ConstructTargetEntranceWindow().

G4ThreeVector QweakSimTarget::positionTargetExitWindow
private
G4LogicalVolume* QweakSimTarget::ScatteringChamberWindow_Logical
private
G4Material* QweakSimTarget::ScatteringChamberWindow_Material
private

Definition at line 129 of file QweakSimTarget.hh.

Referenced by ConstructScatteringChamberWindow(), and QweakSimTarget().

G4VPhysicalVolume* QweakSimTarget::ScatteringChamberWindow_Physical
private

Definition at line 128 of file QweakSimTarget.hh.

Referenced by ConstructScatteringChamberWindow(), and QweakSimTarget().

G4double QweakSimTarget::ScatteringChamberWindowRadius
private
G4double QweakSimTarget::ScatteringChamberWindowThickness
private

Definition at line 158 of file QweakSimTarget.hh.

Referenced by ConstructScatteringChamberWindow(), and QweakSimTarget().

G4LogicalVolume* QweakSimTarget::TargetCell_Logical
private
G4Material* QweakSimTarget::TargetCell_Material
private

Definition at line 113 of file QweakSimTarget.hh.

Referenced by ConstructTargetCell(), and QweakSimTarget().

G4VPhysicalVolume* QweakSimTarget::TargetCell_Physical
private

Definition at line 112 of file QweakSimTarget.hh.

Referenced by ConstructTargetCell(), and QweakSimTarget().

G4double QweakSimTarget::targetCellBackInnerRadiusMax
private
G4double QweakSimTarget::targetCellBackOuterRadiusMax
private

Definition at line 153 of file QweakSimTarget.hh.

Referenced by ConstructTargetCell(), and QweakSimTarget().

G4double QweakSimTarget::targetCellBackRadiusMin
private

Definition at line 151 of file QweakSimTarget.hh.

Referenced by ConstructTargetExitWindow(), and QweakSimTarget().

G4double QweakSimTarget::targetCellEntranceWindowThickness
private
G4double QweakSimTarget::targetCellExitWindowNippleRadius
private

Definition at line 155 of file QweakSimTarget.hh.

Referenced by ConstructTargetExitWindowNipple(), and QweakSimTarget().

G4double QweakSimTarget::targetCellExitWindowNippleThickness
private
G4double QweakSimTarget::targetCellExitWindowThickness
private
G4double QweakSimTarget::targetCellFrontInnerRadiusMax
private
G4double QweakSimTarget::targetCellFrontOuterRadiusMax
private

Definition at line 149 of file QweakSimTarget.hh.

Referenced by ConstructTargetCell(), and QweakSimTarget().

G4double QweakSimTarget::targetCellFrontRadiusMin
private

Definition at line 147 of file QweakSimTarget.hh.

Referenced by ConstructTargetEntranceWindow(), and QweakSimTarget().

G4double QweakSimTarget::targetCellInnerLength
private
G4double QweakSimTarget::targetCellOuterLength
private
G4double QweakSimTarget::targetCellWallThickness
private

Definition at line 142 of file QweakSimTarget.hh.

Referenced by QweakSimTarget().

G4Material* QweakSimTarget::TargetContainer_Material
private

Definition at line 137 of file QweakSimTarget.hh.

Referenced by ConstructTargetContainer(), and QweakSimTarget().

G4LogicalVolume* QweakSimTarget::TargetEntranceWindow_Logical
private
G4Material* QweakSimTarget::TargetEntranceWindow_Material
private

Definition at line 117 of file QweakSimTarget.hh.

Referenced by ConstructTargetEntranceWindow(), and QweakSimTarget().

G4VPhysicalVolume* QweakSimTarget::TargetEntranceWindow_Physical
private

Definition at line 116 of file QweakSimTarget.hh.

Referenced by ConstructTargetEntranceWindow(), and QweakSimTarget().

G4LogicalVolume* QweakSimTarget::TargetExitWindow_Logical
private
G4Material* QweakSimTarget::TargetExitWindow_Material
private

Definition at line 121 of file QweakSimTarget.hh.

Referenced by ConstructTargetExitWindow(), and QweakSimTarget().

G4VPhysicalVolume* QweakSimTarget::TargetExitWindow_Physical
private

Definition at line 120 of file QweakSimTarget.hh.

Referenced by ConstructTargetExitWindow(), and QweakSimTarget().

G4LogicalVolume* QweakSimTarget::TargetExitWindowNipple_Logical
private
G4Material* QweakSimTarget::TargetExitWindowNipple_Material
private

Definition at line 125 of file QweakSimTarget.hh.

Referenced by ConstructTargetExitWindowNipple(), and QweakSimTarget().

G4VPhysicalVolume* QweakSimTarget::TargetExitWindowNipple_Physical
private

Definition at line 124 of file QweakSimTarget.hh.

Referenced by ConstructTargetExitWindowNipple(), and QweakSimTarget().

G4LogicalVolume* QweakSimTarget::TargetMaterial_Logical
private
G4Material* QweakSimTarget::TargetMaterial_Material
private

Definition at line 133 of file QweakSimTarget.hh.

Referenced by ConstructTargetMaterial(), and QweakSimTarget().

G4VPhysicalVolume* QweakSimTarget::TargetMaterial_Physical
private
QweakSimTargetMessenger* QweakSimTarget::targetMessenger
private

Definition at line 164 of file QweakSimTarget.hh.

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

G4VSensitiveDetector* QweakSimTarget::TargetSD
private

Definition at line 172 of file QweakSimTarget.hh.

Referenced by ConstructComponent(), and QweakSimTarget().

G4double QweakSimTarget::targetZPos
private
G4VPhysicalVolume* QweakSimTarget::theMotherPV
private

Definition at line 107 of file QweakSimTarget.hh.

Referenced by ConstructComponent(), ConstructTargetContainer(), and QweakSimTarget().


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