BLCMDspacechargelw Class Reference

Inheritance diagram for BLCMDspacechargelw:

BLCommand BLCollectiveComputation BLManager::SteppingAction BLManager::TrackingAction BLElementField

List of all members.


Detailed Description

class BLCMDspacechargelw - Lienard-Wiechert space-charge computation This computation works within the Geant4 tracking framework, and does the computation in the various callbacks.

The algorithm tracks individual particles, and computes their field using macro-particles that are co-located with the tracked particles. So it computes E and B fields from the macro-particles and adds this in to the BLGlobalField, which Geant4 uses for tracking. The field from each macro-particle is computed via Lienard-Wiechert potentials. When necessary, the trajectories of other macro-particles are linearly extrapolated from their current positions (because of the retarded time, this is needed only needed for macro-particles that are close to the particle currently being tracked). To handle particle creation and destruction, each macro-particle has time limits, and does not contribute any field if its trajectory is evaluated outside those limits. The expensive N^2 computation is done in the beginTrack() callback, not in the much more frequently called addElementField().

For a given particle being tracked, this version first evaluates the other MacroParticles' field at its position at the start of the time step; it then linearly extrapolates the particle's trajectory to the next time step, and evaluates the other MacroParticles' field there; addElementField() then interpolates between those values linearly in global time. As the interpolation is linear, it gives a reasonable value even a short time outside the time step.

Public Member Functions

 BLCMDspacechargelw ()
G4String commandName ()
int command (BLArgumentVector &argv, BLArgumentMap &namedArgs)
void defineNamedArgs ()
void setStochasticsOff (G4Track *track)
void restoreStochastics (G4Track *track)
virtual void beginCollectiveTracking (std::vector< BLTrackData > &v)
virtual void beginTrack (std::vector< BLTrackData > &v, int _index)
virtual void collectiveStep (std::vector< BLTrackData > &v)
virtual void endCollectiveTracking (std::vector< BLTrackData > &v)
virtual void UserSteppingAction (const G4Step *step)
virtual void PreUserTrackingAction (const G4Track *track)
virtual void PostUserTrackingAction (const G4Track *track)
virtual void addFieldValue (const G4double point[4], G4double field[6]) const

Private Attributes

G4double deltaT
G4double radius
G4double charge
G4int trackTwice
G4int verbose
G4int K
G4int ignoreFieldWhenTracking
bool tracking
std::vector< MacroParticlemp
int index
G4double pos0 [4]
G4double pos1 [4]
G4double field0 [6]
G4double field1 [6]
std::vector< Monitormon
BLNTuplentuple
BLManagermanager
BLRunManagerrunManager
BLPhysicsphysics
int nstep
std::vector< int > restoreList
G4double rangeCut

Classes

struct  Monitor


Constructor & Destructor Documentation

BLCMDspacechargelw::BLCMDspacechargelw (  ) 

References BLCMDTYPE_PHYSICS, charge, deltaT, ignoreFieldWhenTracking, index, K, manager, nstep, ntuple, physics, radius, rangeCut, BLCommand::registerCommand(), runManager, BLCommand::setDescription(), BLCommand::setSynopsis(), tracking, trackTwice, and verbose.

00261                                        : BLCollectiveComputation(), BLManager::SteppingAction(),
00262                 BLManager::TrackingAction(), BLElementField(), mp(), mon(),
00263                 restoreList()
00264 {
00265         registerCommand(BLCMDTYPE_PHYSICS);
00266         setSynopsis("Lienard-Wiechert space charge computation");
00267         setDescription("This is a space charge computation that uses "
00268                 "macro-particles to simulate more particles than is feasible "
00269                 "to track individually. Each macro-particle is tracked as a "
00270                 "single particle, but its charge is multiplied by the "
00271                 "macro-particle charge when computing the field. The radius "
00272                 "of the macro-particle is used to avoid the singularity from "
00273                 "a point charge; outside the radius the macro-particle is "
00274                 "treated as a point charge; inside the radius the point-charge "
00275                 "field is multiplied by (r/radius)^K (radius and K are "
00276                 "parameters).\n\n"
00277                 "The trajectory of every particle is kept, and when computing "
00278                 "the field at a point, the intersection of the point's past "
00279                 "lightcone with the trajectory is used to determine the field "
00280                 "from the macro-particle; there is a loop over all particles "
00281                 "except the one currently being tracked. This computation "
00282                 "scales as N^2, where N is the number of macro-particles; that "
00283                 "makes it computationally infeasible for more than a few "
00284                 "hundred macro-particles. But for the particles used, it is "
00285                 "correct to within the following approximations: a) using "
00286                 "macro-particles, b) linearly interpolating between steps, "
00287                 "c) omitting the radiation term in the L-W potential.\n\n"
00288                 "The fields are computed using eq. 63.8-9 (p 162) of Landau "
00289                 "and Lifshitz, _Classical_Theory_of_Fields_, ignoring the "
00290                 "radiation term. The computed fields are used in the usual "
00291                 "Geant4 tracking. Particle creation and destruction are "
00292                 "handled properly.\n\n"
00293                 "This algorithm is primarily intended to test other space "
00294                 "charge algorithms.");
00295 
00296         // initialize class variables here
00297         deltaT = -1.0;
00298         radius = -1.0;
00299         charge = -1.0;
00300         trackTwice = 0;
00301         verbose = -1;
00302         K = 1;
00303         ignoreFieldWhenTracking = 0;
00304         tracking = false;
00305         index = -1;
00306         ntuple = 0;
00307         manager = 0;    // not available yet -- set in command(...)
00308         runManager = 0; // not available yet -- set in command(...)
00309         physics = 0;    // not available yet -- set in beginCollectiveTracking
00310         nstep = 0;
00311         rangeCut = 0.0;
00312 }


Member Function Documentation

G4String BLCMDspacechargelw::commandName (  )  [inline, virtual]

Implements BLCommand.

00233 { return "spacechargelw"; }

int BLCMDspacechargelw::command ( BLArgumentVector argv,
BLArgumentMap namedArgs 
) [virtual]

Implements BLCommand.

References BLGlobalField::addElementField(), charge, deltaT, BLParam::getInt(), BLGlobalField::getObject(), BLRunManager::getObject(), BLManager::getObject(), BLCommand::handleNamedArgs(), manager, Param, BLCommand::print(), BLCommand::printError(), radius, BLRunManager::registerCollectiveComputation(), BLManager::registerSteppingAction(), BLManager::registerTrackingAction(), runManager, BLRunManager::setCollectiveMode(), BLRunManager::setDeltaT(), MacroParticle::testUnits(), MacroParticle::verbose, and verbose.

00315 {
00316         manager = BLManager::getObject();
00317         runManager = BLRunManager::getObject();
00318 
00319         handleNamedArgs(namedArgs);
00320 
00321         if(verbose < 0) verbose=Param.getInt("steppingVerbose");
00322         if(verbose < 0) verbose = 0;
00323         MacroParticle::verbose = verbose;
00324 
00325         if(radius <= 0.0 || charge <= 0.0 || deltaT <= 0.0)
00326                 printError("spacechargelw: radius, charge, deltaT must all be > 0");
00327 
00328         if(K < 0 || K > 4)
00329                 printError("spacechargelw: K must be one of {0,1,2,3,4}");
00330 
00331         runManager->registerCollectiveComputation(this);
00332         runManager->setCollectiveMode(true);
00333         runManager->setDeltaT(deltaT);
00334 
00335         manager->registerSteppingAction(this);
00336         manager->registerTrackingAction(this);
00337 
00338         BLGlobalField::getObject()->addElementField(this); // infinite bounding
00339                                                            // box 
00340         print("");
00341 
00342         MacroParticle::testUnits();
00343 
00344         return 0;
00345 }

void BLCMDspacechargelw::defineNamedArgs (  )  [virtual]

Reimplemented from BLCommand.

References BLCommand::argDouble(), BLCommand::argInt(), charge, deltaT, ignoreFieldWhenTracking, K, radius, trackTwice, and verbose.

00348 {
00349         argDouble(deltaT,"deltaT","Time step (ns).",ns);
00350         argDouble(radius,"radius","Radius of macro-particles (mm).",mm);
00351         argDouble(charge,"charge","Charge of macro-particles (times particle charge).");
00352         argInt(trackTwice,"trackTwice","0=linear extrapolation, 1=track (0)");
00353         argInt(verbose,"verbose","Non-zero for verbose prints (0).");
00354         argInt(K,"K","Exponent for macro-particle density (1).");
00355         argInt(ignoreFieldWhenTracking,"ignoreFieldWhenTracking","For testing only (0).");
00356 }       

void BLCMDspacechargelw::setStochasticsOff ( G4Track *  track  ) 

References BLPhysics::getRangeCut(), BLPhysics::getStochasticsEnabled(), BLPhysics::isStochasticProcess(), physics, rangeCut, and restoreList.

Referenced by beginTrack().

00359 {
00360         restoreList.clear();
00361 
00362         if(!physics->getStochasticsEnabled()) return;
00363 
00364         rangeCut = physics->getRangeCut();
00365 
00366         G4UImanager* UI = G4UImanager::GetUIpointer();
00367         // turn off delta rays (and other particle production) by
00368         // setting production cut very high in energy by setting
00369         // the required range very large.
00370         UI->ApplyCommand("/run/setCut 1 parsec");
00371         // Turn off enegry-loss straggling
00372         UI->ApplyCommand("/process/eLoss/fluct false");
00373 
00374         // run through list of processes, disabling stochastic ones
00375         G4ProcessManager *pm = track->GetDefinition()->GetProcessManager();
00376         G4ProcessVector *pv=pm->GetProcessList();
00377         for(int i=0; i<pv->size(); ++i) {
00378                 if(physics->isStochasticProcess((*pv)[i])) {
00379                         if(pm->GetProcessActivation(i))
00380                                 restoreList.push_back(i);
00381                         pm->SetProcessActivation(i,false);
00382                 }
00383         }
00384 }

void BLCMDspacechargelw::restoreStochastics ( G4Track *  track  ) 

References BLPhysics::getStochasticsEnabled(), physics, rangeCut, and restoreList.

Referenced by beginTrack().

00387 {
00388         if(!physics->getStochasticsEnabled()) return;
00389 
00390         G4UImanager* UI = G4UImanager::GetUIpointer();
00391         static char setCmd[64] = { 0 };
00392         if(setCmd[0] == '\0') sprintf(setCmd,"/run/setCut %g mm",rangeCut/mm);
00393         UI->ApplyCommand(setCmd);
00394         UI->ApplyCommand("/process/eLoss/fluct true");
00395 
00396         // run through list of disabled processes, restoring them
00397         G4ProcessManager *pm = track->GetDefinition()->GetProcessManager();
00398         for(unsigned int j=0; j<restoreList.size(); ++j) {
00399                 pm->SetProcessActivation(restoreList[j],true);
00400         }
00401 }

void BLCMDspacechargelw::beginCollectiveTracking ( std::vector< BLTrackData > &  v  )  [virtual]

Reimplemented from BLCollectiveComputation.

References charge, BLNTuple::create(), deltaT, FIELDS, BLManager::getPhysics(), K, manager, mon, mp, nstep, ntuple, physics, radius, runManager, BLRunManager::setDeltaT(), verbose, and MacroParticle::verbose.

00404 {
00405         physics = manager->getPhysics();
00406 
00407         // initialize mp[] and mon[]
00408         mp.clear();
00409         mon.clear();
00410         for(unsigned i=0; i<v.size(); ++i) {
00411                 mp.push_back(MacroParticle(v[i].track,charge,radius,K));
00412                 mon.push_back(Monitor());
00413         }
00414 
00415         runManager->setDeltaT(deltaT);
00416 
00417         nstep = 0;
00418         MacroParticle::verbose = verbose;
00419 
00420         ntuple = BLNTuple::create("","NTuples","CollectiveMonitor",FIELDS,"");
00421 }

void BLCMDspacechargelw::beginTrack ( std::vector< BLTrackData > &  v,
int  _index 
) [virtual]

Reimplemented from BLCollectiveComputation.

References BLNTuple::appendRow(), BEAM, BLAssert, charge, deltaT, E, field0, field1, BLManager::getState(), BLManager::getSteppingVerbose(), index, K, manager, mon, mp, nstep, ntuple, pos0, pos1, BLRunManager::processOneTrack(), radius, restoreStochastics(), runManager, BLManager::setState(), setStochasticsOff(), SPECIAL, trackTwice, UNSET, and verbose.

00424 {
00425         if(manager->getState() != BEAM) return;
00426 
00427         index = _index;
00428         BLAssert(index >= 0 && (unsigned)index<v.size());
00429 
00430         G4Track *track = v[index].track;
00431         if(verbose==1 || verbose == 2) 
00432                 printf("spacechargelw::beginTrack(%d)\n",index);
00433 
00434         // if this is a new (secondary) particle, add it to mp[] and mon[]
00435         if((unsigned)index >= mp.size()) {
00436                 BLAssert((unsigned)index == mp.size());
00437                 mp.push_back(MacroParticle(track,charge,radius,K));
00438                 mp[index].setTmin(track->GetGlobalTime());
00439                 mon.push_back(Monitor());
00440                 if(verbose >= 3)
00441                         printf("spacechargelw::beginTrack(%d) added to mp[]\n",index);
00442         }
00443 
00444         // compute field of all other MacroParticle-s at this track's position
00445         G4ThreeVector p=track->GetPosition();
00446         pos0[0] = p[0];
00447         pos0[1] = p[1];
00448         pos0[2] = p[2];
00449         pos0[3] = track->GetGlobalTime();
00450         field0[0]=field0[1]=field0[2]=field0[3]=field0[4]=field0[5]=0.0;
00451         for(unsigned i=0; i<mp.size(); ++i) {
00452                 if(i == (unsigned)index) continue;
00453                 mp[i].addField(pos0,field0);
00454         }
00455         if(verbose >= 3)
00456                 printf("spacechargelw::beginTrack(%4d) pos0:%.3f,%.3f,%.3f,%.3f field0:%.3f,%.3f,%.3f,%.3f,%.3f,%.3f\n",
00457                         index,pos0[0],pos0[1],pos0[2],pos0[3],
00458                         field0[0]/tesla,field0[1]/tesla,field0[2]/tesla,
00459                         field0[3]/(megavolt/meter),field0[4]/(megavolt/meter),
00460                         field0[5]/(megavolt/meter));
00461 
00462         if(trackTwice) {
00463                 // track with stochastics off, full external field, and
00464                 // constant space charge field (at track's current position).
00465                 int isave = index;
00466                 pos1[0] = pos0[0]; pos1[1] = pos0[1]; pos1[2] = pos0[2];
00467                 pos1[3] = pos0[3] + deltaT; // values for a constant field
00468                 G4Track *temp = new G4Track(*track);
00469                 manager->setState(SPECIAL);
00470                 setStochasticsOff(temp);
00471                 if(manager->getSteppingVerbose() != 0)
00472                         printf("spacechargelw: pre-tracking begun\n");
00473                 runManager->processOneTrack(temp); // will be suspended
00474                 if(manager->getSteppingVerbose() != 0)
00475                         printf("spacechargelw: pre-tracking ended\n");
00476                 restoreStochastics(temp);
00477                 manager->setState(BEAM);
00478                 p = temp->GetPosition();
00479                 pos1[0] = p[0];
00480                 pos1[1] = p[1];
00481                 pos1[2] = p[2];
00482                 pos1[3] = temp->GetGlobalTime();
00483                 delete temp;
00484                 index = isave;
00485                 // protect against immediate track death
00486                 if(pos1[3]-pos0[3] < 1.0E-4*ns) pos1[3] += 1.0E-4*ns;
00487         } else {
00488                 // Linearly extrapolate this track to the next time step.
00489                 p = track->GetMomentumDirection();
00490                 G4double speed = track->GetVelocity();
00491                 pos1[0] = pos0[0] + p[0]*speed*deltaT;
00492                 pos1[1] = pos0[1] + p[1]*speed*deltaT;
00493                 pos1[2] = pos0[2] + p[2]*speed*deltaT;
00494                 pos1[3] = pos0[3] + deltaT;
00495         }
00496 
00497         // compute field of all other MacroParticle-s at extrapolated position
00498         // (this can be at a global time beyond the known trajectory of some of 
00499         //  the other macro-particles; they are linearly extrapolated inside
00500         //  their addfield()).
00501         field1[0]=field1[1]=field1[2]=field1[3]=field1[4]=field1[5]=0.0;
00502         for(unsigned i=0; i<mp.size(); ++i) {
00503                 if(i == (unsigned)index) continue;
00504                 mp[i].addField(pos1,field1);
00505         }
00506         if(verbose >= 3)
00507                 printf("                        pos1:%.3f,%.3f,%.3f,%.3f field1:%.3f,%.3f,%.3f,%.3f,%.3f,%.3f\n",
00508                         pos1[0],pos1[1],pos1[2],pos1[3],
00509                         field1[0]/tesla,field1[1]/tesla,field1[2]/tesla,
00510                         field1[3]/(megavolt/meter),field1[4]/(megavolt/meter),
00511                         field1[5]/(megavolt/meter));
00512 
00513         // monitor the differences
00514         // (mon[index] contains the previous extrapolated values)
00515         if(mon[index].pos[0] != UNSET && fabs(pos0[3]-mon[index].pos[3]) < 0.002*ns) {
00516                 G4ThreeVector deltaPos, deltaB, deltaE;
00517                 for(int i=0; i<3; ++i) {
00518                         deltaPos[i] = pos0[i] - mon[index].pos[i];
00519                         deltaB[i] = field0[i] - mon[index].field[i];
00520                         deltaE[i] = field0[i+3] - mon[index].field[i+3];
00521                 }
00522                 if(verbose >= 3)
00523                     printf("deltaPos=%.6f,%.6f,%.6f deltaB=%.6f,%.6f,%.6f "
00524                         "deltaE=%.5f,%.6f,%.6f\n",
00525                         deltaPos[0]/mm,deltaPos[1]/mm,deltaPos[2]/mm,
00526                         deltaB[0]/tesla,deltaB[1]/tesla,deltaB[2]/tesla,
00527                         deltaE[0]/(megavolt/meter),deltaE[1]/(megavolt/meter),
00528                         deltaE[2]/(megavolt/meter));
00529                 double data[11];
00530                 data[0] = pos0[0];      // x
00531                 data[1] = pos0[1];      // y
00532                 data[2] = pos0[2];      // z
00533                 data[3] = pos0[3];      // t
00534                 data[4] = nstep;        // nStep
00535                 data[5] = mon[index].deltaT;    // deltaT
00536                 data[6] = sqrt(field0[0]*field0[0]+field0[1]*field0[1]+
00537                                         field0[2]*field0[2]); // Btot
00538                 data[7] = sqrt(field0[3]*field0[3]+field0[4]*field0[4]+
00539                                         field0[5]*field0[5]); // Etot
00540                 data[8] = deltaPos.mag(); // deltaPos
00541                 data[9] = deltaB.mag(); // deltaB
00542                 data[10] = deltaE.mag();        // deltaE
00543                 ntuple->appendRow(data,11);
00544         }
00545         // put the new extrapolated values into mon[]
00546         mon[index] = Monitor(pos1,field1,deltaT);
00547 
00548         // Now that we have done the expensive N^2 computation, BLCMDspacechargelw::addField
00549         // will interpolate between pos0[] and pos1[] to approximate the field
00550         // at the trajectory of this particle with the corresponding value on 
00551         // the linearly-extrapolated trajectory.
00552 }

void BLCMDspacechargelw::collectiveStep ( std::vector< BLTrackData > &  v  )  [virtual]

Implements BLCollectiveComputation.

References BEAM, BLManager::getState(), manager, and mp.

00555 {
00556         if(manager->getState() != BEAM) return;
00557 
00558         // TODO: update radius here
00559         // TODO: update deltaT here
00560 
00561         // prune MacroParticle trajectories
00562         for(unsigned i=0; i<mp.size(); ++i)
00563                 mp[i].pruneTrajectory();
00564 }

void BLCMDspacechargelw::endCollectiveTracking ( std::vector< BLTrackData > &  v  )  [virtual]

Reimplemented from BLCollectiveComputation.

00567 {
00568 }

void BLCMDspacechargelw::UserSteppingAction ( const G4Step *  step  )  [virtual]

Implements BLManager::SteppingAction.

References BEAM, BLAssert, deltaT, BLManager::getState(), index, manager, mp, and verbose.

00571 {
00572         if(manager->getState() != BEAM) return; // includes SPECIAL
00573 
00574         BLAssert(index>=0 && (unsigned)index<mp.size());
00575 
00576         // add trajectory point unless tiny increment in time
00577         if(step->GetTrack()->GetGlobalTime()-mp[index].getMostRecentTime() >=
00578                                                                 deltaT/10.0) {
00579                 mp[index].addPoint(step->GetTrack());
00580                 if(verbose >= 3)
00581                         printf("spacechargelw::UserSteppingAction index=%d point added\n",index);
00582         } else {
00583                 if(verbose >= 3)
00584                         printf("spacechargelw::UserSteppingAction index=%d NO POINT - dt too small\n",index);
00585         }
00586 }

void BLCMDspacechargelw::PreUserTrackingAction ( const G4Track *  track  )  [virtual]

Implements BLManager::TrackingAction.

References BEAM, BLAssert, BLManager::getState(), index, manager, mp, tracking, and verbose.

00589 {
00590         if(manager->getState() != BEAM) return;
00591 
00592         if(verbose >= 3)
00593                 printf("spacechargelw::PreUserTrackingAction index=%d\n",index);
00594 
00595         BLAssert(index>=0 && (unsigned)index<mp.size());
00596 
00597         tracking = true;
00598 }

void BLCMDspacechargelw::PostUserTrackingAction ( const G4Track *  track  )  [virtual]

Implements BLManager::TrackingAction.

References BEAM, BLAssert, deltaT, BLManager::getState(), index, manager, mp, tracking, and verbose.

00601 {
00602         if(manager->getState() != BEAM) return;
00603 
00604         BLAssert(index>=0 && (unsigned)index<mp.size());
00605 
00606         tracking = false;
00607 
00608         // add this point to the trajectory (no call to UserSteppingAction
00609         // when track is suspended).
00610         if(track->GetGlobalTime()-mp[index].getMostRecentTime() > deltaT/100.0)
00611                 mp[index].addPoint(track);
00612 
00613         // set Tmax for tracks that are killed
00614         G4TrackStatus trackStatus = track->GetTrackStatus();
00615         if(trackStatus != fAlive && trackStatus != fStopButAlive && trackStatus != fSuspend) {
00616                 mp[index].setTmax(track->GetGlobalTime());
00617                 if(verbose >= 3)
00618                         printf("spacechargelw::PostUserTrackingAction index=%d track killed\n",index);
00619         } else {
00620                 if(verbose >= 3)
00621                         printf("spacechargelw::PostUserTrackingAction index=%d\n",index);
00622         }
00623 
00624         index = -1;
00625 }

void BLCMDspacechargelw::addFieldValue ( const G4double  point[4],
G4double  field[6] 
) const [virtual]

Implements BLElementField.

References BEAM, BLAssert, field0, field1, BLManager::getState(), ignoreFieldWhenTracking, index, manager, mp, pos0, pos1, tracking, and verbose.

00628 {
00629         if(manager->getState() != BEAM) return;
00630         if(ignoreFieldWhenTracking != 0 && tracking) return;
00631 
00632         double fthis[6];
00633         fthis[0]=fthis[1]=fthis[2]=fthis[3]=fthis[4]=fthis[5]=0.0;
00634 
00635         if(index >= 0) {
00636                 // interpolate field between pos0 and pos1 for this track
00637                 BLAssert(pos1[3]-pos0[3] > 1.0e-8*ns);
00638                 G4double f = (point[3]-pos0[3])/(pos1[3]-pos0[3]);
00639                 fthis[0] = (1.0-f)*field0[0] + f*field1[0];
00640                 fthis[1] = (1.0-f)*field0[1] + f*field1[1];
00641                 fthis[2] = (1.0-f)*field0[2] + f*field1[2];
00642                 fthis[3] = (1.0-f)*field0[3] + f*field1[3];
00643                 fthis[4] = (1.0-f)*field0[4] + f*field1[4];
00644                 fthis[5] = (1.0-f)*field0[5] + f*field1[5];
00645                 if(verbose >= 3)
00646                         printf("spacechargelw::addFieldValue index=%d f=%.4f point=%.3f,%.3f,%.3f,%.3f\n",
00647                                 index,f,point[0],point[1],point[2],point[3]);
00648         } else {
00649                 // not during tracking, sum all macro-particles
00650                 for(unsigned i=0; i<mp.size(); ++i) {
00651                     ((BLCMDspacechargelw*)this)->mp[i].addField(point,fthis);
00652                 }
00653         }
00654         if(verbose >= 3)
00655                 printf("spacechargelw::addFieldValue pos:%.3f,%.3f,%.3f,%.3f field:%.3f,%.3f,%.3f,%.3f,%.3f,%.3f\n",
00656                         point[0],point[1],point[2],point[3],
00657                         fthis[0]/tesla,fthis[1]/tesla,fthis[2]/tesla,
00658                         fthis[3]/(megavolt/meter),fthis[4]/(megavolt/meter),
00659                         fthis[5]/(megavolt/meter));
00660 
00661         field[0] += fthis[0];
00662         field[1] += fthis[1];
00663         field[2] += fthis[2];
00664         field[3] += fthis[3];
00665         field[4] += fthis[4];
00666         field[5] += fthis[5];
00667 }


Member Data Documentation

G4double BLCMDspacechargelw::deltaT [private]

G4double BLCMDspacechargelw::radius [private]

G4double BLCMDspacechargelw::charge [private]

G4int BLCMDspacechargelw::verbose [private]

G4int BLCMDspacechargelw::K [private]

std::vector<MacroParticle> BLCMDspacechargelw::mp [private]

G4double BLCMDspacechargelw::pos0[4] [private]

Referenced by addFieldValue(), and beginTrack().

G4double BLCMDspacechargelw::pos1[4] [private]

Referenced by addFieldValue(), and beginTrack().

G4double BLCMDspacechargelw::field0[6] [private]

Referenced by addFieldValue(), and beginTrack().

G4double BLCMDspacechargelw::field1[6] [private]

Referenced by addFieldValue(), and beginTrack().

std::vector<Monitor> BLCMDspacechargelw::mon [private]

std::vector<int> BLCMDspacechargelw::restoreList [private]

G4double BLCMDspacechargelw::rangeCut [private]


The documentation for this class was generated from the following file:
g4beamline