BLCMDreweightprocess Class Reference

Inheritance diagram for BLCMDreweightprocess:

BLCommand BLCallback

List of all members.


Detailed Description

class BLCMDreweightprocess will re-weight physics processes.

This command permits the user to artificially increase (or decrease) the cross-section of one or more physics processes applied to one or more particles. This can be used to reduce the variance of a simulation intended to measure a rare process. The track weights are adjusted accordingly.

For example, in a simulation with 1E6 beam particles, essentially none will get through an absorber that is 10 interaction lengths long. To determine how many make it through for a beam of 1E12 particles, the interaction cross-section must be reduced.

Public Member Functions

 BLCMDreweightprocess ()
 Constructor.
virtual G4String commandName ()
 commandName() returns "reweightprocess".
virtual int command (BLArgumentVector &argv, BLArgumentMap &namedArgs)
 command() implements the reweightprocess command.
virtual void defineNamedArgs ()
 defineNamedArgs() defines the named arguments for this command.
void callback (int type)
 callback() from BLCallback.

Private Member Functions

bool isMatchingParticle (G4ParticleDefinition *pd)
bool isMatchingProcess (const G4VProcess *proc)
bool wrapProcess (G4ProcessManager *pmgr, G4VProcess *proc)

Private Attributes

G4String particle
G4String process
G4double ratio


Constructor & Destructor Documentation

BLCMDreweightprocess::BLCMDreweightprocess (  ) 

Constructor.

References BLCMDTYPE_PHYSICS, particle, process, ratio, BLCommand::registerCommand(), BLCommand::setDescription(), and BLCommand::setSynopsis().

Referenced by command().

00134                                            : BLCommand(), BLCallback()
00135 {
00136         registerCommand(BLCMDTYPE_PHYSICS);
00137         setSynopsis("modify the cross-section of a physics process.");
00138         setDescription("This command will modify the cross-section of a "
00139                 "physics process, modifying the track weights so that "
00140                 "weighted histograms give the same statistical result as if "
00141                 "the command were not used. Used properly, this can greatly "
00142                 "reduce the variance of the result.\n\n"
00143                 "Care should be taken to ensure that all regions that should "
00144                 "be sampled actually are sampled. For instance, if ratio>1 "
00145                 "the interaction length will be reduced, and deep inside an "
00146                 "absorber there may be no sampling because no simulated tracks "
00147                 "ever get there, even though real tracks will. If ratio<1 then "
00148                 "the upstream regions of absorbers will be under sampled; this "
00149                 "is usually OK, as the desired result is to increase the "
00150                 "sampling deep inside the absorber.\n\n"
00151                 "For ratio>>1 this command can be used to examine rare "
00152                 "processes. This can induce multiple rare interactions in an "
00153                 "event when normally none would be expected; the weights will "
00154                 "still correctly correspond to the real interaction, even "
00155                 "though the event topologies don't.\n\n"
00156                 "This command cannot reweight any continuous process (e.g. "
00157                 "multiple scattering, ionization energy loss, etc.) -- it "
00158                 "will issue a fatal exception if applied to such a process.\n\n"
00159                 "This command should not be applied to other processes that "
00160                 "re-weight tracks (e.g. the neutrino command). "
00161                 "Indeed it probably won't give the correct weights if any "
00162                 "such process applies to the particle (except for itself, it "
00163                 "cannot determine the unmodified interaction length of such "
00164                 "procsses, which is needed to compute the weight).\n\n"
00165                 "The re-weighting applies to both PostStep and AtRest "
00166                 "processes, but some AtRest processes do not work with this "
00167                 "re-weighting; for instance, Decay works properly in PostStep "
00168                 "for a moving particle, but not for a stopped one AtRest. "
00169                 "This is related to the Geant4 limitation that exactly one "
00170                 "process be active AtRest, and exactly one step be taken.\n\n"
00171                 "Be sure to test your use of this commnd for a simple physical "
00172                 "situation before believing its results."
00173         );
00174 
00175         particle = "";
00176         process = "";
00177         ratio = 1.0;
00178 }


Member Function Documentation

bool BLCMDreweightprocess::isMatchingParticle ( G4ParticleDefinition *  pd  )  [inline, private]

References BLCommand::matchList(), and particle.

Referenced by callback().

00059             { return particle.size()==0 || 
00060                                 matchList(pd->GetParticleName(),particle); }

bool BLCMDreweightprocess::isMatchingProcess ( const G4VProcess *  proc  )  [inline, private]

References BLCommand::matchList(), and process.

Referenced by callback().

00062             { return matchList(proc->GetProcessName(),process); }

bool BLCMDreweightprocess::wrapProcess ( G4ProcessManager *  pmgr,
G4VProcess *  proc 
) [private]

References ratio.

Referenced by callback().

00243 {
00244         // cannot re-weight any continuous processes (don't know how to do it)
00245         if(pmgr->GetAlongStepIndex(proc) >= 0) return false;
00246 
00247         ReweightProcess *rp = new ReweightProcess(pmgr,proc,ratio);
00248 
00249         // fix the positions of all processes
00250         for(int i=idxAtRest; i<=idxPostStep; ++i) {
00251                 G4ProcessVector *pv = pmgr->GetProcessVector(
00252                                         (G4ProcessVectorDoItIndex)i,typeDoIt);
00253                 for(int k=1; k<pv->size(); ++k) { // 0 is Transportation
00254                         pmgr->SetProcessOrdering((*pv)[k],
00255                                         (G4ProcessVectorDoItIndex)i,k);
00256                 }
00257         }
00258 
00259         int orderAtRest = pmgr->GetProcessOrdering(proc,idxAtRest);
00260         int orderAlongStep = ordLast;
00261         int orderPostStep = pmgr->GetProcessOrdering(proc,idxPostStep);
00262 
00263         pmgr->RemoveProcess(proc);
00264 
00265         pmgr->AddProcess(rp,orderAtRest,orderAlongStep,orderPostStep);
00266 
00267         return true;
00268 }

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

commandName() returns "reweightprocess".

Implements BLCommand.

00069 { return "reweightprocess"; }

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

command() implements the reweightprocess command.

Implements BLCommand.

References BLCMDreweightprocess(), BLManager::getObject(), BLCommand::handleNamedArgs(), BLCommand::print(), and BLManager::registerCallback().

00181 {
00182         BLCMDreweightprocess *p = 
00183                         new BLCMDreweightprocess(defaultBLCMDreweightprocess);
00184 
00185         int retval = p->handleNamedArgs(namedArgs);
00186 
00187         // too early to modify processes here
00188         BLManager::getObject()->registerCallback(p,0);
00189 
00190         p->print("");
00191 
00192         return retval;
00193 }

void BLCMDreweightprocess::defineNamedArgs (  )  [virtual]

defineNamedArgs() defines the named arguments for this command.

Reimplemented from BLCommand.

References BLCommand::argDouble(), BLCommand::argString(), particle, process, and ratio.

00196 {
00197         argString(particle,"particle","Comma-separated list of particle "
00198                                                 "patterns ('' => all).");
00199         argString(process,"process","Comma-separated list of process patterns.");
00200         argDouble(ratio,"ratio","Ratio of artificial to real cross-section.");
00201 }

void BLCMDreweightprocess::callback ( int  type  )  [virtual]

callback() from BLCallback.

Reimplemented from BLCallback.

References BLAssert, BLManager::getObject(), isMatchingParticle(), isMatchingProcess(), BLManager::showAllExceptions(), and wrapProcess().

00204 {
00205         BLManager *blmanager = BLManager::getObject();
00206         if(type == 0) {
00207                 // wrap matching processes for matching particles
00208                 G4ParticleTable::G4PTblDicIterator *theParticleIterator =
00209                         G4ParticleTable::GetParticleTable()->GetIterator();
00210                 theParticleIterator->reset();
00211                 while((*theParticleIterator)()) {
00212                         G4ParticleDefinition *pd = theParticleIterator->value();
00213                         if(!isMatchingParticle(pd)) continue;
00214                         G4ProcessManager *pmgr = pd->GetProcessManager();
00215                         BLAssert(pmgr != 0);
00216                         int n=pmgr->GetProcessListLength();
00217                         G4ProcessVector *pv=pmgr->GetProcessList();
00218                         for(int i=0; i<n; ++i) {
00219                             G4VProcess *p=(*pv)[i];
00220                             if(!isMatchingProcess(p)) continue;
00221                             char tmp[256];
00222                             sprintf(tmp,"process=%s particle=%s",
00223                                                 p->GetProcessName().c_str(),
00224                                                 pd->GetParticleName().c_str());
00225                             if(wrapProcess(pmgr,p))  {
00226                                 bool flag=blmanager->showAllExceptions(true);
00227                                 G4Exception("reweightprocess",
00228                                                 "Process re-weighted",
00229                                                 JustWarning, tmp);
00230                                 blmanager->showAllExceptions(flag);
00231                             } else {
00232                                 G4Exception("reweightprocess",
00233                                                 "Process cannot be re-weighted",
00234                                                 FatalException, tmp);
00235                             }
00236                         }
00237                 }
00238         }
00239 }


Member Data Documentation

G4String BLCMDreweightprocess::particle [private]

G4String BLCMDreweightprocess::process [private]

G4double BLCMDreweightprocess::ratio [private]


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