Public Member Functions | |
BLCMDprofile () | |
Constructor. | |
BLCMDprofile (BLCMDprofile &r) | |
virtual G4String | commandName () |
commandName() returns "profile". | |
virtual int | command (BLArgumentVector &argv, BLArgumentMap &namedArgs) |
command() implements the profile command. | |
virtual void | defineNamedArgs () |
defineNamedArgs() defines the named arguments for this command. | |
virtual void | callback (int type) |
callback from BLCallback. | |
Private Attributes | |
G4String | z |
G4String | zloop |
G4String | require |
G4String | file |
G4String | particle |
G4String | coordinates |
BLCoordinateType | coordinateType |
std::vector< Entry * > | entries |
unsigned | place |
FILE * | fd |
G4ParticleDefinition * | pd |
Classes | |
struct | Entry |
BLCMDprofile::BLCMDprofile | ( | ) |
Constructor.
References BLCMDTYPE_DATA, BLCOORD_CENTERLINE, coordinates, coordinateType, fd, file, particle, pd, place, BLCommand::registerCommand(), require, BLCommand::setDescription(), BLCommand::setSynopsis(), z, and zloop.
Referenced by command().
00090 : BLCommand(), BLCallback(), entries() 00091 { 00092 registerCommand(BLCMDTYPE_DATA); 00093 setSynopsis("write beam profile information to a file"); 00094 setDescription("This command accumulates the moments of the track\n" 00095 "distributions during the run, and at the end of\n" 00096 "run prints the mean, sigma, emittance, alpha, and\n" 00097 "beta (Twiss parameters) for the tracks. Each z\n" 00098 "position generates a line in the output file.\n\n" 00099 "Each value in z and zloop can be an expression using " 00100 "double constants and the usual C operators and " 00101 "functions."); 00102 00103 z = ""; 00104 zloop = ""; 00105 require = ""; 00106 file = "-"; 00107 particle = "mu+"; 00108 coordinates = "Centerline"; 00109 coordinateType = BLCOORD_CENTERLINE; 00110 place = 0; 00111 fd = 0; 00112 pd = 0; 00113 }
BLCMDprofile::BLCMDprofile | ( | BLCMDprofile & | r | ) |
References coordinates, coordinateType, fd, file, particle, pd, place, require, z, and zloop.
00115 : BLCommand(r), 00116 BLCallback(), entries() 00117 { 00118 z = r.z; 00119 zloop = r.zloop; 00120 require = r.require; 00121 file = r.file; 00122 particle = r.particle; 00123 coordinates = r.coordinates; 00124 coordinateType = r.coordinateType; 00125 place = r.place; 00126 fd = 0; 00127 pd = 0; 00128 }
virtual G4String BLCMDprofile::commandName | ( | ) | [inline, virtual] |
int BLCMDprofile::command | ( | BLArgumentVector & | argv, | |
BLArgumentMap & | namedArgs | |||
) | [virtual] |
command() implements the profile command.
Implements BLCommand.
References BLCMDprofile(), coordinates, coordinateType, entries, fd, file, BLAsciiFile::fopen(), BLCoordinates::getCoordinateType(), BLCommand::getList(), BLManager::getObject(), BLCommand::handleNamedArgs(), particle, pd, BLCommand::print(), BLCommand::printError(), BLManager::registerCallback(), BLManager::registerZStep(), require, z, and zloop.
00131 { 00132 BLCMDprofile *p = new BLCMDprofile(defaultProfile); 00133 00134 int retval = p->handleNamedArgs(namedArgs); 00135 00136 p->coordinateType = BLCoordinates::getCoordinateType(p->coordinates); 00137 00138 if(p->file == "-") { 00139 p->fd = stdout; 00140 } else { 00141 p->fd = BLAsciiFile::fopen(p->file); 00142 if(!p->fd) printError("profile: Cannot write file '%s'", 00143 p->file.c_str()); 00144 } 00145 00146 p->pd = G4ParticleTable::GetParticleTable()->FindParticle(p->particle); 00147 if(!p->pd) printError("profile: unknown particle '%s'\n",p->particle.c_str()); 00148 00149 BLEvaluator eval; 00150 00151 // register the z position(s) in the string z 00152 std::vector<G4double> vect = getList(p->z,","); 00153 if(vect.size() > 0) { 00154 for(unsigned i=0; i<vect.size(); ++i) { 00155 Entry *e = new Entry(vect[i],p->require,p->pd, 00156 p->coordinateType); 00157 BLManager::getObject()->registerZStep(vect[i],e,4); 00158 p->entries.push_back(e); 00159 } 00160 } else if(p->z.size() > 0) { 00161 printError("printf: Syntax error in z"); 00162 } 00163 00164 // register the z position(s) in the string zloop 00165 vect = getList(p->zloop,",:"); 00166 if(vect.size() == 3) { 00167 if(vect[2] < 1.0*mm || vect[0] > vect[1]) { 00168 printError("printf: invalid zloop '%s'", 00169 p->zloop.c_str()); 00170 } else { 00171 while(vect[0] <= vect[1]) { 00172 Entry *e = new Entry(vect[0],p->require,p->pd, 00173 p->coordinateType); 00174 BLManager::getObject()->registerZStep(vect[0], 00175 e,4); 00176 p->entries.push_back(e); 00177 vect[0] += vect[2]; 00178 } 00179 } 00180 } else if(p->zloop.size() > 0) { 00181 printError("printf: invalid zloop '%s'",p->zloop.c_str()); 00182 } 00183 00184 BLManager::getObject()->registerCallback(p,2); 00185 00186 p->print(""); 00187 00188 return retval; 00189 }
void BLCMDprofile::defineNamedArgs | ( | ) | [virtual] |
defineNamedArgs() defines the named arguments for this command.
Reimplemented from BLCommand.
References BLCommand::argString(), coordinates, file, particle, require, z, and zloop.
00192 { 00193 argString(z,"z","Comma-separated list of Z positions for profiling (mm)"); 00194 argString(zloop,"zloop","Loop in z, first:last:incr (mm)"); 00195 argString(require,"require","logical expression for cutting (default=true)"); 00196 argString(particle,"particle","Name of particle to profile (default=mu+)"); 00197 argString(file,"file","Output filename (default=stdout)"); 00198 argString(file,"filename","Synonym for file"); 00199 argString(coordinates,"coordinates","Coordinates: centerline or reference (default=c)."); 00200 }
void BLCMDprofile::callback | ( | int | type | ) | [virtual] |
callback from BLCallback.
Reimplemented from BLCallback.
References determinant(), entries, BLAsciiFile::fclose(), fd, file, pd, and z.
00259 { 00260 // after beam tracking is complete, write the beam profile data 00261 if(type != 2 || !fd) return; 00262 fprintf(fd,"# g4beamline profile\n"); 00263 fprintf(fd,"#Z N meanX sigmaX meanY sigmaY emitX emitY emitTrans betaX betaY betaTrans alphaX alphaY alphaTrans meanP\n"); 00264 fprintf(fd,"#mm - mm mm mm mm mm-rad mm-rad mm-rad mm mm mm - - - MeV/c\n"); 00265 for(unsigned i=0; i<entries.size(); ++i) { 00266 int n = entries[i]->n; 00267 G4double sumW = entries[i]->sumW; 00268 if(n <= 0 || fabs(sumW) < 1.0e-12) continue; 00269 G4double meanX = entries[i]->sumX/sumW; 00270 G4double meanY = entries[i]->sumY/sumW; 00271 G4double meanPx = entries[i]->sumPx/sumW; 00272 G4double meanPy = entries[i]->sumPy/sumW; 00273 G4double meanPz = entries[i]->sumPz/sumW; 00274 G4double meanPtot = entries[i]->sumPtot/sumW; 00275 G4double **M = entries[i]->M; 00276 M[0][0] = M[0][0]/sumW - meanX*meanX; 00277 M[0][1] = M[0][1]/sumW - meanX*meanPx; 00278 M[0][2] = M[0][2]/sumW - meanX*meanY; 00279 M[0][3] = M[0][3]/sumW - meanX*meanPy; 00280 M[1][0] = M[1][0]/sumW - meanPx*meanX; 00281 M[1][1] = M[1][1]/sumW - meanPx*meanPx; 00282 M[1][2] = M[1][2]/sumW - meanPx*meanY; 00283 M[1][3] = M[1][3]/sumW - meanPx*meanPy; 00284 M[2][0] = M[2][0]/sumW - meanY*meanX; 00285 M[2][1] = M[2][1]/sumW - meanY*meanPx; 00286 M[2][2] = M[2][2]/sumW - meanY*meanY; 00287 M[2][3] = M[2][3]/sumW - meanY*meanPy; 00288 M[3][0] = M[3][0]/sumW - meanPy*meanX; 00289 M[3][1] = M[3][1]/sumW - meanPy*meanPx; 00290 M[3][2] = M[3][2]/sumW - meanPy*meanY; 00291 M[3][3] = M[3][3]/sumW - meanPy*meanPy; 00292 G4double sigmaX = sqrt(M[0][0]); 00293 G4double sigmaY = sqrt(M[2][2]); 00294 G4double mass = pd->GetPDGMass(); 00295 G4double emitX = sqrt(M[0][0]*M[1][1]-M[1][0]*M[0][1])/mass; 00296 G4double emitY = sqrt(M[2][2]*M[3][3]-M[3][2]*M[2][3])/mass; 00297 G4double emitTrans = sqrt(sqrt(determinant(M,4)))/mass; 00298 G4double alphaX = -M[0][1]/mass/emitTrans; 00299 G4double alphaY = -M[2][3]/mass/emitTrans; 00300 G4double alphaTrans = (alphaX+alphaY)/2.0; 00301 G4double betaX = M[0][0]*meanPz/mass/emitTrans; 00302 G4double betaY= M[2][2]*meanPz/mass/emitTrans; 00303 G4double betaTrans = (betaX+betaY)/2.0; 00304 fprintf(fd,"%.1f %d %.3f %.3f %.3f %.3f %.4g %.4g %.4g %.4g %.4g %.4g %.4g %.4g %.4g %.4g\n", 00305 entries[i]->z,n,meanX,sigmaX,meanY,sigmaY,emitX,emitY,emitTrans, 00306 betaX,betaY,betaTrans,alphaX,alphaY,alphaTrans,meanPtot/MeV); 00307 } 00308 if(fd != stdout) 00309 BLAsciiFile::fclose(fd); 00310 printf("profile: wrote '%s'\n",file.c_str()); 00311 }
G4String BLCMDprofile::z [private] |
Referenced by BLCMDprofile(), callback(), command(), and defineNamedArgs().
G4String BLCMDprofile::zloop [private] |
Referenced by BLCMDprofile(), command(), and defineNamedArgs().
G4String BLCMDprofile::require [private] |
Referenced by BLCMDprofile(), command(), and defineNamedArgs().
G4String BLCMDprofile::file [private] |
Referenced by BLCMDprofile(), callback(), command(), and defineNamedArgs().
G4String BLCMDprofile::particle [private] |
Referenced by BLCMDprofile(), command(), and defineNamedArgs().
G4String BLCMDprofile::coordinates [private] |
Referenced by BLCMDprofile(), command(), and defineNamedArgs().
BLCoordinateType BLCMDprofile::coordinateType [private] |
Referenced by BLCMDprofile(), and command().
std::vector<Entry*> BLCMDprofile::entries [private] |
Referenced by callback(), and command().
unsigned BLCMDprofile::place [private] |
Referenced by BLCMDprofile().
FILE* BLCMDprofile::fd [private] |
Referenced by BLCMDprofile(), callback(), and command().
G4ParticleDefinition* BLCMDprofile::pd [private] |
Referenced by BLCMDprofile(), callback(), and command().