#include <BLWindowShape.hh>
This is intended to be used to create a G4PolyCone for the window, a G4PolyCone for its contents, and a G4Tubs for the flange and possibly the pipe to which the flange is bolted (bolts are not modeled).
This class may be subclassed to implement other methods of constructing a window shape.
The standard constructor reads data from filename. The file format is a series of lines: First character # means comment, * means printed comment. Comment and blank lines are ignored. Units are mm. The first line contains the 4 flange variables: innerR outerR insideZ outsideZ The remaining lines contain 3 values for a given radius: r z t The first line must have r=0.0 and have the largest z value and the largest z+t value. z is the inside of the window, z+t is the outside of the window. All values must be positive. Successive r values must increase by at least 0.010 mm. Successive z values and z+t values must decrease by at least 0.010 mm. flangeInnerRadius should equal the last r value. Any z origin may be used (it will be subtracted away). This is intended to be easy to interface to window design spreadsheets (export a list of 3 columns delimited by spaces, then add the appropriate comments and flange values at the top).
Public Member Functions | |
BLWindowShape () | |
default constructor. | |
BLWindowShape (G4String filename) | |
Constructor. Reads directory/filename for the data. | |
virtual | ~BLWindowShape () |
destructor. | |
Public Attributes | |
G4double | flangeInnerRadius |
G4double | flangeOuterRadius |
G4double | flangeInsideZ |
G4double | flangeOutsideZ |
std::vector< G4double > | r |
std::vector< G4double > | z |
std::vector< G4double > | t |
BLWindowShape::BLWindowShape | ( | ) | [inline] |
default constructor.
References flangeInnerRadius, flangeInsideZ, flangeOuterRadius, and flangeOutsideZ.
00066 : r(), z(), t() { flangeInnerRadius=0.0; 00067 flangeOuterRadius=0.0; flangeOutsideZ=0.0; flangeInsideZ=0.0; }
BLWindowShape::BLWindowShape | ( | G4String | filename | ) |
Constructor. Reads directory/filename for the data.
References flangeInnerRadius, flangeInsideZ, flangeOuterRadius, flangeOutsideZ, BLCommand::printError(), r, t, and z.
00025 : r(), z(), t() 00026 { 00027 flangeInnerRadius = 0.0; 00028 flangeOuterRadius = 0.0; 00029 flangeInsideZ = 0.0; 00030 flangeOutsideZ = 0.0; 00031 00032 G4String file = filename; 00033 00034 std::ifstream in; 00035 if(filename != "") 00036 in.open(file.c_str()); 00037 if(!in.good()) { 00038 BLCommand::printError("BLWindowShape Cannot open file '%s'", 00039 file.c_str()); 00040 return; 00041 } 00042 00043 bool flangeLine = true; 00044 int lineno = 0; 00045 while(in.good()) { 00046 char line[1024+1]; 00047 line[0] = '\0'; 00048 in.getline(line,sizeof(line)); 00049 ++lineno; 00050 if(line[0] == '*') 00051 printf("%s\n",line); 00052 if(line[0] == '\0' || line[0] == '#' || line[0] == '*') 00053 continue; 00054 int n = 0; 00055 if(flangeLine) { 00056 flangeLine = false; 00057 n = sscanf(line,"%lf %lf %lf %lf",&flangeInnerRadius, 00058 &flangeOuterRadius,&flangeInsideZ,&flangeOutsideZ); 00059 n -= 1; // 4 items here, check is for 3 00060 } else { 00061 G4double rr=0.0,zz=0.0,tt=0.0; 00062 n = sscanf(line,"%lf %lf %lf",&rr,&zz,&tt); 00063 r.push_back(rr); 00064 z.push_back(zz); 00065 t.push_back(tt); 00066 } 00067 if(n != 3) { 00068 BLCommand::printError("BLWindowShape Invalid syntax" 00069 " '%s' line %d\n",file.c_str(),lineno); 00070 break; 00071 } 00072 } 00073 in.close(); 00074 }
G4double BLWindowShape::flangeInnerRadius |
Referenced by BLWindowShape(), and BLCMDabsorber::constructSolids().
G4double BLWindowShape::flangeOuterRadius |
Referenced by BLWindowShape(), and BLCMDabsorber::constructSolids().
G4double BLWindowShape::flangeInsideZ |
Referenced by BLWindowShape().
G4double BLWindowShape::flangeOutsideZ |
Referenced by BLWindowShape(), and BLCMDabsorber::constructSolids().
std::vector<G4double> BLWindowShape::r |
Referenced by BLWindowShape(), and BLCMDabsorber::constructSolids().
std::vector<G4double> BLWindowShape::z |
Referenced by BLWindowShape(), and BLCMDabsorber::constructSolids().
std::vector<G4double> BLWindowShape::t |
Referenced by BLWindowShape(), and BLCMDabsorber::constructSolids().