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.
This computation uses a Poisson solver in the beam frame. The beam frame is determined from the reference particle's position and 4-momentum.
It is best to terminate the simulation with "trackcuts maxTime=...", or minActive=.... If it is terminated because particles leave the world or hit an object with kill=1, when the # particles in the bunch gets very low PoissonConvolve3D may not converge, aborting the run. Fortunately the NTuples wre written, so this is OK, it just looks ugly.
Note that fewer than minActive particles in the bunch is treated as a NORMAL way to end the simulation, not an error. But messages beginning with "***" are printed to alert the user that this happened.
This version requires the reference particle's momentum to be parallel to the z axis.
Public Member Functions | |
BLCMDspacecharge () | |
G4String | commandName () |
int | command (BLArgumentVector &argv, BLArgumentMap &namedArgs) |
void | defineNamedArgs () |
virtual void | beginCollectiveTracking (std::vector< BLTrackData > &v) |
virtual void | collectiveStep (std::vector< BLTrackData > &v) |
virtual void | endCollectiveTracking (std::vector< BLTrackData > &v) |
virtual void | PreUserTrackingAction (const G4Track *track) |
TrackingAction and SteppingAction are used only to get the points of each reference trajectory into bunch[]. PreUserTrackingAction() creates each Bunch object, seting currentBunch which is used by the others. | |
virtual void | PostUserTrackingAction (const G4Track *track) |
virtual void | UserSteppingAction (const G4Step *step) |
virtual void | addFieldValue (const G4double point[4], G4double field[6]) const |
Private Attributes | |
G4double | deltaT |
G4double | charge |
G4int | nx |
G4int | ny |
G4int | nz |
G4double | dx |
G4double | dy |
G4double | dz |
G4double | maxBeta |
G4int | verbose |
G4int | ignoreFieldWhenTracking |
G4int | useApproximationOnly |
G4int | fixedGrid |
double | percentile |
int | nxApprox |
int | nyApprox |
int | nzApprox |
int | minActive |
BLManager * | manager |
BLRunManager * | runManager |
int | nstep |
std::vector< Bunch > | bunch |
Bunch * | currentBunch |
bool | tracking |
BLCMDspacecharge::BLCMDspacecharge | ( | ) |
References BLCMDTYPE_PHYSICS, charge, currentBunch, deltaT, dx, dy, dz, fixedGrid, ignoreFieldWhenTracking, manager, maxBeta, minActive, nstep, nx, nxApprox, ny, nyApprox, nz, nzApprox, percentile, BLCommand::registerCommand(), runManager, BLCommand::setDescription(), BLCommand::setSynopsis(), tracking, useApproximationOnly, and verbose.
00331 : BLCollectiveComputation(), 00332 BLManager::SteppingAction(), BLManager::TrackingAction(), 00333 BLElementField(), bunch() 00334 { 00335 registerCommand(BLCMDTYPE_PHYSICS); 00336 setSynopsis("Beam-frame Green's function space charge computation"); 00337 setDescription("This is a space charge computation for bunched beams. " 00338 "It uses a grid in the beam frame to solve Poisson's equation " 00339 "via a Green's function with infinite boundary conditions; " 00340 "the E field is boosted back to the lab frame E and B for " 00341 "tracking.\n\n" 00342 "Macro-particles are used to enable the simulation of larger " 00343 "bunches than can be feasibly simulated as individual " 00344 "particles; the macro-particles have zero radius, " 00345 "but are pro-rated into the nearest eight grid points when " 00346 "placed into the grid. This computation can handle up to about " 00347 "a million macro-particles, but 100,000 is more sensible for " 00348 "all but the simplest physical situations.\n\n" 00349 "The bunch is created from the beam tracks before tracking " 00350 "begins. There is one bunch for each reference particle. " 00351 "Particles in the bunch must be the same particle as the " 00352 "reference, must initially be within {dx,dy,dz} of the " 00353 "reference particle, and when boosted to the reference " 00354 "particle's rest frame must initially have beta < maxBeta.\n\n" 00355 "After boosting the particles to the beam frame, they are " 00356 "placed into the grid, pro-rating to the eight nearest grid " 00357 "points. The grid is dynamically re-sized to keep the " 00358 "99th percentile of the particles between 0.5 and 0.67 of " 00359 "the grid size. This maintains a reasonable balance between " 00360 "resolution of grid points within the bunch and covering all " 00361 "of the particles. Particles located at >85% of the grid size " 00362 "do not contribute to the field computation, but are tracked " 00363 "using the field of the rest of the bunch (and other bunches)." 00364 "\n\n" 00365 "The grid has {nx,ny,nz} points; there is a small " 00366 "computational advantage to using powers of 2, but any values " 00367 ">1 can be used. For efficiency, the convolution of the " 00368 "Green's function with the charge grid is performed using " 00369 "FFTs; the grid is doubled in each dimension with the proper " 00370 "symmetry applied to the Green's function, so the cyclical " 00371 "convolution of the FFTs gives the proper potential with " 00372 "infinite boundary conditions.\n\n" 00373 "Outside the grid an approximation is used. An approximation " 00374 "grid is constructed, with the same size of the Poisson grid, " 00375 "but using {nxApprox,nyApprox,nzApprox} points. Particles are " 00376 "placed into this approximation grid, and the mean position " 00377 "is kept as well as the charge. Approximation grid points " 00378 "with less than 1% of the total charge are " 00379 "consolidated with their inner neighbors. The non-zero " 00380 "approximation grid points are treated as point charges when " 00381 "computing the potential outside the grid. For reasonably " 00382 "Gaussian bunches, {7,7,7} are reasonable values for the " 00383 "approximation grid sizes.\n\n" 00384 "The E field in the beam frame is computed via the derivatives " 00385 "of the linear interpolating function using the eight nearest " 00386 "grid points, and is boosted back to the lab " 00387 "frame E and B for tracking by the usual Geant4 routines.\n\n" 00388 "Bunch particles that get destroyed cease contributing to the " 00389 "bunch. As the bunch particles are selected during start-up, " 00390 "no additional particles are ever added to a bunch. This " 00391 "algorithm handles multiple bunches of any particle types.\n\n" 00392 "NOTE: For now, the reference MUST be parallel to the Z " 00393 "axis."); 00394 00395 // initialize class variables here 00396 deltaT = -1.0; 00397 charge = -1.0; 00398 nx = 65; 00399 ny = 65; 00400 nz = 65; 00401 dx = -1.0; 00402 dy = -1.0; 00403 dz = -1.0; 00404 maxBeta = 0.1; 00405 verbose = -1; 00406 ignoreFieldWhenTracking = 0; 00407 useApproximationOnly = 0; 00408 fixedGrid = 0; 00409 percentile = 99; 00410 nxApprox = 7; 00411 nyApprox = 7; 00412 nzApprox = 7; 00413 minActive = -95; 00414 manager = 0; // not available yet -- set in command(...) 00415 runManager = 0; // not available yet -- set in command(...) 00416 nstep = 0; 00417 currentBunch = 0; 00418 tracking = false; 00419 }
G4String BLCMDspacecharge::commandName | ( | ) | [inline, virtual] |
int BLCMDspacecharge::command | ( | BLArgumentVector & | argv, | |
BLArgumentMap & | namedArgs | |||
) | [virtual] |
Implements BLCommand.
References BLGlobalField::addElementField(), charge, deltaT, dx, dy, dz, BLParam::getInt(), BLGlobalField::getObject(), BLRunManager::getObject(), BLManager::getObject(), BLCommand::handleNamedArgs(), manager, nx, nxApprox, ny, nyApprox, nz, nzApprox, Param, BLCommand::print(), BLCommand::printError(), BLRunManager::registerCollectiveComputation(), BLManager::registerSteppingAction(), BLManager::registerTrackingAction(), runManager, BLRunManager::setCollectiveMode(), BLRunManager::setDeltaT(), Bunch::setUseApproximationOnly(), Bunch::setVerbose(), Bunch::testBfield(), Bunch::testEfield(), useApproximationOnly, and verbose.
00422 { 00423 if(manager != 0) { 00424 printError("spacecharge - can be used at most once."); 00425 return -1; 00426 } 00427 00428 handleNamedArgs(namedArgs); 00429 00430 manager = BLManager::getObject(); 00431 runManager = BLRunManager::getObject(); 00432 00433 if(verbose < 0) verbose = Param.getInt("steppingVerbose"); 00434 if(verbose < 0) verbose = 0; 00435 Bunch::setVerbose(verbose); 00436 Bunch::setUseApproximationOnly(useApproximationOnly); 00437 00438 if(nxApprox < 0) nxApprox = nx; 00439 if(nyApprox < 0) nyApprox = ny; 00440 if(nzApprox < 0) nzApprox = nz; 00441 00442 if(charge <= 0.0 || deltaT <= 0.0) 00443 printError("spacecharge: charge and deltaT must both be > 0"); 00444 if(dx <= 0.0 || dy <= 0.0 || dz <= 0.0) 00445 printError("spacecharge: dx, dy, and dz must all be > 0"); 00446 00447 runManager->registerCollectiveComputation(this); 00448 runManager->setCollectiveMode(true); 00449 runManager->setDeltaT(deltaT); 00450 00451 manager->registerSteppingAction(this); 00452 manager->registerTrackingAction(this); 00453 00454 BLGlobalField::getObject()->addElementField(this); // infinite bounding 00455 // box 00456 print(""); 00457 00458 Bunch::testEfield(); 00459 Bunch::testBfield(); 00460 00461 return 0; 00462 }
void BLCMDspacecharge::defineNamedArgs | ( | ) | [virtual] |
Reimplemented from BLCommand.
References BLCommand::argDouble(), BLCommand::argInt(), charge, deltaT, dx, dy, dz, fixedGrid, ignoreFieldWhenTracking, maxBeta, minActive, nx, nxApprox, ny, nyApprox, nz, nzApprox, percentile, useApproximationOnly, and verbose.
00465 { 00466 argDouble(deltaT,"deltaT","Time step (ns).",ns); 00467 argDouble(charge,"charge","Charge of macro particles (times particle charge)."); 00468 argInt(nx,"nx","Number of grid points in x (65)."); 00469 argInt(ny,"ny","Number of grid points in y (65)."); 00470 argInt(nz,"nz","Number of grid points in z (65)."); 00471 argDouble(dx,"dx","Max distance of particle to reference in x (mm)",mm); 00472 argDouble(dy,"dy","Max distance of particle to reference in y (mm)",mm); 00473 argDouble(dz,"dz","Max distance of particle to reference in z (mm)",mm); 00474 argInt(nxApprox,"nxApprox","# bins in x in approximation (7)."); 00475 argInt(nyApprox,"nyApprox","# bins in y in approximation (7)."); 00476 argInt(nzApprox,"nzApprox","# bins in z in approximation (7)."); 00477 argDouble(maxBeta,"maxBeta","Max beta (v/c) of particle in beam frame (0.1)."); 00478 argInt(verbose,"verbose","Non-zero for verbose prints (0)."); 00479 argInt(ignoreFieldWhenTracking,"ignoreFieldWhenTracking","For testing only (0)."); 00480 argInt(useApproximationOnly,"useApproximationOnly","For testing only (0)."); 00481 argInt(fixedGrid,"fixedGrid","Nonzero prevents re-sizing the grid (0)."); 00482 argDouble(percentile,"percentile","Percentile of charge distribution used for grid sizing (99)."); 00483 argInt(minActive,"minActive","Minimum # active tracks in bunch; if " 00484 "< 0 is % of initial bunch size (-95)."); 00485 }
void BLCMDspacecharge::beginCollectiveTracking | ( | std::vector< BLTrackData > & | v | ) | [virtual] |
Reimplemented from BLCollectiveComputation.
References BLAssert, bunch, deltaT, nstep, nx, nxApprox, ny, nyApprox, nz, nzApprox, runManager, BLRunManager::setDeltaT(), and verbose.
00488 { 00489 if(verbose >= 1) printf("BLCMDspacecharge::beginCollectiveTracking v.size=%d\n",(int)v.size()); 00490 00491 runManager->setDeltaT(deltaT); 00492 nstep = 0; 00493 00494 // add tracks to Bunch-es 00495 std::map<G4ParticleDefinition*,int> orphans; 00496 int nactive=0, norphan=0; 00497 for(unsigned j=0; j<v.size(); ++j) { 00498 G4Track *track = v[j].track; 00499 G4TrackStatus trackStatus = track->GetTrackStatus(); 00500 if(trackStatus != fAlive && trackStatus != fStopButAlive) 00501 continue; 00502 ++nactive; 00503 bool match=false; 00504 for(unsigned i=0; i<bunch.size(); ++i) { 00505 if(bunch[i].addTrack(track)) 00506 goto found; 00507 if(track->GetDefinition() == bunch[i].getDefinition()) 00508 match = true; 00509 } 00510 ++orphans[track->GetDefinition()]; 00511 if(match) ++norphan; 00512 found: continue; 00513 } 00514 00515 00516 for(unsigned i=0; i<bunch.size(); ++i) { 00517 BLAssert(bunch[i].init()); 00518 printf("Bunch %d:\n",i); 00519 printf(" Bunch has %d tracks, minActive=%d\n", 00520 bunch[i].size(),bunch[i].getMinActive()); 00521 G4ThreeVector b=bunch[i].getBound(); 00522 printf(" Grid has nx,ny,nx=%d,%d,%d\n",nx,ny,nz); 00523 printf(" Initial Grid size: %.3f,%.3f,%.3f\n", 00524 b.x(),b.y(),b.z()); 00525 printf(" Initial Approximation has nx,ny,nx=%d,%d,%d " 00526 " %d charges\n",nxApprox,nyApprox,nzApprox, 00527 bunch[i].getNapprox()); 00528 } 00529 00530 if(orphans.size() > 0) { 00531 printf("spacecharge: orphan tracks:"); 00532 std::map<G4ParticleDefinition*,int>::iterator i; 00533 char sep=' '; 00534 for(i=orphans.begin(); i!=orphans.end(); ++i) { 00535 printf("%c%d*%s",sep,i->second, 00536 i->first->GetParticleName().c_str()); 00537 sep = ','; 00538 } 00539 printf("\n"); 00540 } 00541 if(norphan > nactive/1000) 00542 G4Exception("spacecharge","Too many orphan tracks",JustWarning, 00543 "More than 0.1% of matching tracks not in any Bunch"); 00544 }
void BLCMDspacecharge::collectiveStep | ( | std::vector< BLTrackData > & | v | ) | [virtual] |
Implements BLCollectiveComputation.
References bunch, nstep, and runManager.
00547 { 00548 for(unsigned i=0; i<bunch.size(); ++i) { 00549 if(!bunch[i].computeField()) { 00550 printf("***Ending computation\n"); 00551 runManager->AbortRun(); 00552 break; 00553 } 00554 } 00555 00556 ++nstep; 00557 }
void BLCMDspacecharge::endCollectiveTracking | ( | std::vector< BLTrackData > & | v | ) | [virtual] |
Reimplemented from BLCollectiveComputation.
References nstep.
00560 { 00561 printf("BLCMDspacecharge::endCollectiveTracking after %d time steps, v.size=%d\n", 00562 nstep,(int)v.size()); 00563 }
void BLCMDspacecharge::PreUserTrackingAction | ( | const G4Track * | track | ) | [virtual] |
TrackingAction and SteppingAction are used only to get the points of each reference trajectory into bunch[]. PreUserTrackingAction() creates each Bunch object, seting currentBunch which is used by the others.
Implements BLManager::TrackingAction.
References bunch, charge, currentBunch, dx, dy, dz, fixedGrid, BLManager::getState(), manager, maxBeta, minActive, nx, nxApprox, ny, nyApprox, nz, nzApprox, percentile, REFERENCE, tracking, and verbose.
00566 { 00567 tracking = true; 00568 currentBunch = 0; 00569 if(manager->getState() == REFERENCE && 00570 track->GetDefinition()->GetPDGCharge() != 0.0) { 00571 bunch.push_back(Bunch(track->GetDefinition(),nx,ny,nz,dx,dy,dz, 00572 charge,maxBeta,fixedGrid,percentile, 00573 nxApprox,nyApprox,nzApprox,minActive)); 00574 currentBunch = &bunch[bunch.size()-1]; 00575 if(verbose >= 1) printf("spacecharge: starting Reference\n"); 00576 } 00577 }
void BLCMDspacecharge::PostUserTrackingAction | ( | const G4Track * | track | ) | [virtual] |
Implements BLManager::TrackingAction.
References currentBunch, Bunch::nTrajectoryPoints(), tracking, and verbose.
00580 { 00581 tracking = false; 00582 if(currentBunch == 0) return; 00583 00584 if(currentBunch->nTrajectoryPoints() < 2) 00585 G4Exception("BLCMDspacecharge","Too few reference steps", 00586 FatalException, "Reference particle took < 2 steps"); 00587 currentBunch = 0; 00588 00589 if(verbose >= 1) printf("spacecharge: ending Reference\n"); 00590 }
void BLCMDspacecharge::UserSteppingAction | ( | const G4Step * | step | ) | [virtual] |
require the reference be parallel to the z axis
Implements BLManager::SteppingAction.
References Bunch::addReferencePoint(), and currentBunch.
00593 { 00594 if(currentBunch == 0) return; 00595 00596 G4Track *track = step->GetTrack(); 00597 00598 /// require the reference be parallel to the z axis 00599 G4ThreeVector direction = track->GetMomentumDirection(); 00600 if(direction.perp() > 1.0e-9) 00601 G4Exception("BLCMDspacecharge","Invalid Bunch momentum", 00602 FatalException, "Reference not parallel to z axis."); 00603 00604 currentBunch->addReferencePoint(track); 00605 }
void BLCMDspacecharge::addFieldValue | ( | const G4double | point[4], | |
G4double | field[6] | |||
) | const [virtual] |
Implements BLElementField.
References BEAM, bunch, BLManager::getState(), ignoreFieldWhenTracking, manager, and tracking.
00608 { 00609 if(manager->getState() != BEAM) return; 00610 if(ignoreFieldWhenTracking != 0 && tracking) return; 00611 00612 for(unsigned i=0; i<bunch.size(); ++i) { 00613 bunch[i].addFieldValue(point,field); 00614 } 00615 }
G4double BLCMDspacecharge::deltaT [private] |
Referenced by beginCollectiveTracking(), BLCMDspacecharge(), command(), and defineNamedArgs().
G4double BLCMDspacecharge::charge [private] |
Referenced by BLCMDspacecharge(), command(), defineNamedArgs(), and PreUserTrackingAction().
G4int BLCMDspacecharge::nx [private] |
Referenced by beginCollectiveTracking(), BLCMDspacecharge(), command(), defineNamedArgs(), and PreUserTrackingAction().
G4int BLCMDspacecharge::ny [private] |
Referenced by beginCollectiveTracking(), BLCMDspacecharge(), command(), defineNamedArgs(), and PreUserTrackingAction().
G4int BLCMDspacecharge::nz [private] |
Referenced by beginCollectiveTracking(), BLCMDspacecharge(), command(), defineNamedArgs(), and PreUserTrackingAction().
G4double BLCMDspacecharge::dx [private] |
Referenced by BLCMDspacecharge(), command(), defineNamedArgs(), and PreUserTrackingAction().
G4double BLCMDspacecharge::dy [private] |
Referenced by BLCMDspacecharge(), command(), defineNamedArgs(), and PreUserTrackingAction().
G4double BLCMDspacecharge::dz [private] |
Referenced by BLCMDspacecharge(), command(), defineNamedArgs(), and PreUserTrackingAction().
G4double BLCMDspacecharge::maxBeta [private] |
Referenced by BLCMDspacecharge(), defineNamedArgs(), and PreUserTrackingAction().
G4int BLCMDspacecharge::verbose [private] |
Referenced by beginCollectiveTracking(), BLCMDspacecharge(), command(), defineNamedArgs(), PostUserTrackingAction(), and PreUserTrackingAction().
G4int BLCMDspacecharge::ignoreFieldWhenTracking [private] |
Referenced by addFieldValue(), BLCMDspacecharge(), and defineNamedArgs().
G4int BLCMDspacecharge::useApproximationOnly [private] |
Referenced by BLCMDspacecharge(), command(), and defineNamedArgs().
G4int BLCMDspacecharge::fixedGrid [private] |
Referenced by BLCMDspacecharge(), defineNamedArgs(), and PreUserTrackingAction().
double BLCMDspacecharge::percentile [private] |
Referenced by BLCMDspacecharge(), defineNamedArgs(), and PreUserTrackingAction().
int BLCMDspacecharge::nxApprox [private] |
Referenced by beginCollectiveTracking(), BLCMDspacecharge(), command(), defineNamedArgs(), and PreUserTrackingAction().
int BLCMDspacecharge::nyApprox [private] |
Referenced by beginCollectiveTracking(), BLCMDspacecharge(), command(), defineNamedArgs(), and PreUserTrackingAction().
int BLCMDspacecharge::nzApprox [private] |
Referenced by beginCollectiveTracking(), BLCMDspacecharge(), command(), defineNamedArgs(), and PreUserTrackingAction().
int BLCMDspacecharge::minActive [private] |
Referenced by BLCMDspacecharge(), defineNamedArgs(), and PreUserTrackingAction().
BLManager* BLCMDspacecharge::manager [private] |
Referenced by addFieldValue(), BLCMDspacecharge(), command(), and PreUserTrackingAction().
BLRunManager* BLCMDspacecharge::runManager [private] |
Referenced by beginCollectiveTracking(), BLCMDspacecharge(), collectiveStep(), and command().
int BLCMDspacecharge::nstep [private] |
Referenced by beginCollectiveTracking(), BLCMDspacecharge(), collectiveStep(), and endCollectiveTracking().
std::vector<Bunch> BLCMDspacecharge::bunch [private] |
Referenced by addFieldValue(), beginCollectiveTracking(), collectiveStep(), and PreUserTrackingAction().
Bunch* BLCMDspacecharge::currentBunch [private] |
Referenced by BLCMDspacecharge(), PostUserTrackingAction(), PreUserTrackingAction(), and UserSteppingAction().
bool BLCMDspacecharge::tracking [private] |
Referenced by addFieldValue(), BLCMDspacecharge(), PostUserTrackingAction(), and PreUserTrackingAction().