4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
12 * Copyright (c) 1991-1999
13 * BIOSON Research Institute, Dept. of Biophysical Chemistry
14 * University of Groningen, The Netherlands
17 * GROMACS: A message-passing parallel molecular dynamics implementation
18 * H.J.C. Berendsen, D. van der Spoel and R. van Drunen
19 * Comp. Phys. Comm. 91, 43-56 (1995)
21 * Also check out our WWW page:
22 * http://md.chem.rug.nl/~gmx
27 * Great Red Oystrich Makes All Chemists Sane
29 static char *SRCID_ionize_c = "$Id$";
47 #define PREFIX "IONIZE: "
49 enum { eionCYL, eionSURF, eionGAUSS, eionNR };
51 enum { ecollPHOTO, ecollINELASTIC, ecollNR };
55 real fj,sigPh,sigIn,vAuger;
59 int nelec,maxelec,elec0,elmin_type;
62 /* BEGIN GLOBAL VARIABLES */
63 static int Energies[] = { 6, 8, 10, 12, 15, 20 };
64 static int ionize_seed = 1993;
65 #define NENER asize(Energies)
67 static t_electron_db edb;
69 /* END GLOBAL VARIABLES */
71 void dump_ca(FILE *fp,t_cross_atom *ca,int i,char *file,int line)
73 fprintf(fp,PREFIX"(line %d) atom %d, z = %d, n = %d, k = %d\n",
74 line,i,ca->z,ca->n,ca->k);
77 t_cross_atom *mk_cross_atom(FILE *log,t_mdatoms *md,
78 char **atomname[],int Eindex)
80 int elem_index[] = { 0, 0, 0, 0, 0, 0, 1, 2, 3, 0, 0, 0, 0, 0, 0, 0, 4, 0, 0, 0, 5, 6 };
86 fprintf(log,PREFIX"Filling data structure for ionization\n");
87 fprintf(log,PREFIX"Warning: all fj values set to 0.95 for now\n");
89 snew(elemcnt,NELEM+1);
90 for(i=0; (i<md->nr); i++) {
94 for(j=0; (j<NELEM); j++)
95 if (strncmp(cc,element[j].name,strlen(element[j].name)) == 0) {
96 ca[i].z = element[j].nel;
100 fatal_error(0,PREFIX"Don't know number of electrons for %s",
104 ca[i].sigPh = element[elem_index[ca[i].z]].cross[Eindex].photo;
105 ca[i].sigIn = element[elem_index[ca[i].z]].cross[Eindex].incoh;
106 ca[i].fj = recoil[ca[i].z].Prob_K;
109 ca[i].vAuger = 0.904;
112 ca[i].vAuger = 0.920;
115 ca[i].vAuger = 0.929;
126 fprintf(log,PREFIX"You have the following elements in your system (%d atoms):\n"PREFIX,md->nr);
127 for(j=0; (j<NELEM); j++)
129 fprintf(log," %s: %d",element[j].name,elemcnt[j]);
130 fprintf(log," atoms\n");
137 int number_K(t_cross_atom *ca)
145 int number_L(t_cross_atom *ca)
147 return ca->k-2+ca->z-ca->n;
150 real xray_cross_section(int eColl,t_cross_atom *ca)
164 c = (nK*0.5*ca->fj + nL/(ca->z-2)*(1-ca->fj))*ca->sigPh;
167 c = (ca->z-ca->n)*ca->sigIn/ca->z;
170 fatal_error(0,"No such collision type %d\n",eColl);
175 real prob_K(int eColl,t_cross_atom *ca)
179 if ((ca->z <= 2) || (ca->z == ca->n))
184 Pl = (ca->k-2+ca->z-ca->n)*(1-ca->fj)/(ca->z-2);
185 Pk = (2-ca->k)*ca->fj*0.5;
189 P = (2-ca->k)/(ca->z-ca->n);
192 fatal_error(0,"No such collision type %d\n",eColl);
197 double myexp(double x)
205 real ptheta_incoh(int Eindex,real theta)
206 /* theta should be in degrees */
208 /* These numbers generated by fitting 5 gaussians to the real function
209 * that describes the probability for theta.
210 * We use symmetry in the gaussian (see 180-angle) therefore there
211 * are fewer parameters (only 8 per energylevel).
213 static double ppp[NENER][8] = {
214 { -0.00295169, 10.4847, 0.0341099, /*-*/43.1963,
215 -0.0164054, 30.2452, 71.0311, 2.50282 },
216 { -0.00370852, 9.02037, 0.100559, /*-*/42.9962,
217 -0.0537891, 35.5077, 71.4305, 1.05515 },
218 { -0.00427039, 7.86831, 0.118042, /*-*/45.9846,
219 -0.0634505, 38.6134, 70.3857, 0.240082 },
220 { -0.004514, 7.0728, 0.13464, /*-*/48.213,
221 -0.0723, 41.06, 69.38, -0.02 },
222 { -0.00488796, 5.87988, 0.159574, /*-*/51.5556,
223 -0.0855767, 44.7307, 69.0251, -0.414604 },
224 { -0.00504604, 4.56299, 0.201064, /*-*/54.8599,
225 -0.107153, 48.7016, 68.8212, -0.487699 }
227 double g1,g2,g3,g4,g5,ptheta;
229 g1 = myexp(-0.5*sqr((theta-ppp[Eindex][7])/ppp[Eindex][1]));
230 g2 = myexp(-0.5*sqr((theta-180+ppp[Eindex][7])/ppp[Eindex][1]));
231 g3 = myexp(-0.5*sqr((theta-90)/ppp[Eindex][3]));
232 g4 = myexp(-0.5*sqr((theta-ppp[Eindex][6])/ppp[Eindex][5]));
233 g5 = myexp(-0.5*sqr((theta-180+ppp[Eindex][6])/ppp[Eindex][5]));
235 ptheta = ppp[Eindex][0]*(g1+g2) + ppp[Eindex][2]*g3 + ppp[Eindex][4]*(g4+g5);
240 real rand_theta_incoh(int Eindex,int *seed)
244 static bool bFirst = TRUE;
246 static int i,j,cur=1;
247 real theta,sum,rrr,dx;
250 dx = 90.0/(real)NINTP;
252 /* Compute cumulative integrals of all probability distributions */
254 for(i=0; (i<NENER); i++) {
255 snew(intp[i],NINTP+1);
256 y[prev] = ptheta_incoh(i,0.0);
258 for(j=1; (j<=NINTP); j++) {
259 y[cur] = ptheta_incoh(i,j*dx);
261 intp[i][j] = intp[i][j-1] + (y[cur]+y[prev])*dx;
266 fprintf(debug,"Integrated probability functions for theta incoherent\n");
267 for(j=0; (j<NINTP); j++) {
268 fprintf(debug,"%10f",dx*j);
269 for(i=0; (i<NENER); i++)
270 fprintf(debug," %10f",intp[i][j]);
278 for(j=0; (j<NINTP) && (rrr > intp[Eindex][j]); j++)
281 return (j-1+(rrr-intp[Eindex][j-1])/(intp[Eindex][j]-intp[Eindex][j-1]))*dx;
284 static void polar2cart(real phi,real theta,rvec v)
286 v[XX] = cos(phi)*sin(theta);
287 v[YY] = sin(phi)*sin(theta);
291 void rand_vector(rvec v,int *seed)
295 theta = 180.0*rando(seed);
296 phi = 360.0*rando(seed);
297 polar2cart(phi,theta,v);
300 real electron_cross_section(FILE *fp,rvec v,real mass,int nelec)
302 /* Compute cross section for electrons */
303 real T,B,U,S,Q,R,N,t,u,lnt,sigma;
304 real a0 = 0.05292; /* nm */
306 /* Have to determine T (kinetic energy of electron) */
307 T = 0.5*mass*iprod(v,v);
309 /* R is the binding energy of the electron in hydrogen */
310 R = 13.61*ELECTRONVOLT;
312 /* Have to determine the binding energy B, differs per orbital of course */
315 /* Have to determine the orbital kinetic energy U */
318 /* Have to know number of electrons */
321 /* Magic constant Q */
324 /* Some help variables */
327 S = 4*M_PI*sqr(a0)*N*sqr(R/B);
330 /* Resulting variable */
331 sigma = (S/(t+u+1))*( 0.5*Q*lnt*(1-1/sqr(t)) + (2-Q)*(1-1/t-lnt/(t+1)) );
336 static bool analyze_electrons(FILE *fp,t_electron_db *edb,
337 int natom,char **atomname[])
342 cc = getenv("GENERATE_ELECTRONS");
343 if (sscanf(cc,"%d",&etp) == 1) {
344 for(i=0; (i<natom); i++) {
345 if (strcmp(*atomname[i],"EL") == 0)
350 edb->maxelec = natom-i;
351 edb->elmin_type = etp;
352 fprintf(fp,PREFIX"There are %d possible electrons\n",edb->maxelec);
357 fprintf(fp,PREFIX"No electron features today.\n");
362 void add_electron(FILE *fp,t_mdatoms *md,t_electron_db *edb,int ion,
363 rvec x[],rvec v[],rvec dv,real dt)
368 if (edb->nelec < edb->maxelec) {
369 ee = edb->elec0+edb->nelec++;
370 md->chargeA[ee] = md->chargeB[ee] = md->chargeT[ee] = -1;
371 md->typeA[ee] = md->typeB[ee] = edb->elmin_type;
374 svmul(-md->massA[ion]*md->invmass[ee],dv,v[ee]);
375 /* Do a first step to prevent the electron from being on top of the
376 * nucleus, move it 0.05 nm from the nucleus
378 nv = 1.0/norm(v[ee]);
379 for(m=0; (m<DIM); m++)
380 x[ee][m] = x[ion][m] + v[ee][m]*nv*0.05;
383 fatal_error(0,PREFIX"No more particles to turn into electrons\n");
386 bool khole_decay(FILE *fp,t_cross_atom *ca,rvec x[],rvec v[],int ion,
387 int *seed,real dt,bool bElectron,
388 t_mdatoms *md,t_electron_db *edb)
394 if ((ca->vAuger < 0) || (recoil[ca->z].tau == 0)) {
395 dump_ca(stderr,ca,ion,__FILE__,__LINE__);
398 if (rando(seed) < dt/recoil[ca->z].tau) {
400 fprintf(debug,"DECAY: Going to decay a k hole\n");
403 /* Generate random vector */
404 rand_vector(dv,seed);
408 fprintf(debug,"DECAY: factor=%10g, dv = (%8.3f, %8.3f, %8.3f)\n",
409 factor,dv[XX],dv[YY],dv[ZZ]);
413 /* Now put the electron in place */
415 add_electron(fp,md,edb,ion,x,v,dv,dt);
423 real electron_atom_interactions(FILE *fp,t_mdatoms *md,t_inputrec *ir,
425 rvec x[],rvec v[],rvec f[],matrix box)
427 /* Compute what the name says... */
430 real mindist2,vc,vtot,fscal,fc,dx2,dx_1,qi,*q;
432 mindist2 = sqr(0.05);
438 if (ir->ePBC != epbcNONE)
440 /* the end variable usually includes electrons */
441 e1 = min(end,edb.elec0);
442 for(i=start; (i<e1); i++) {
443 elec1 = edb.elec0 + edb.nelec;
444 qi = q[i]*ONE_4PI_EPS0;
445 for(j=edb.elec0; (j<elec1); j++) {
446 if (ir->ePBC == epbcNONE)
447 rvec_sub(x[i],x[j],dx);
449 pbc_dx(x[i],x[j],dx);
451 if (dx2 < mindist2) {
452 /* Model collision */
455 /* Do normal coulomb */
459 fscal = vc*dx_1*dx_1;
460 for(m=0; (m<DIM); m++) {
472 void ionize(FILE *fp,t_mdatoms *md,char **atomname[],real t,t_inputrec *ir,
473 rvec x[],rvec v[],int start,int end,matrix box,t_commrec *cr)
475 static FILE *xvg,*ion;
476 static char *leg[] = { "Probability", "Primary Ionization", "Integral over PI", "KHole-Decay", "Integral over KD" };
477 static bool bFirst = TRUE,bElectron = FALSE;
478 static real t0,imax,width,inv_nratoms,rho,nphot;
479 static real interval;
480 static int dq_tot,nkd_tot,ephot,mode;
481 static t_cross_atom *ca;
482 static int Eindex=-1;
484 real r,factor,ndv,E_lost=0,cross_atom,dvz,rrc;
485 real pt,ptot,pphot,pcoll[ecollNR],tmax;
486 real incoh,incoh_abs,sigmaPincoh,hboxx,hboxy,rho2;
488 bool bIonize=FALSE,bKHole,bL,bDOIT;
490 int i,j,k,kk,m,nK,nL,dq,nkh,nkdecay,elmin_type;
491 int *nionize,*nkhole,*ndecay,nbuf[2];
494 /* Get parameters for gaussian photon pulse from inputrec */
495 t0 = ir->userreal1; /* Peak of the gaussian pulse */
496 nphot = ir->userreal2; /* Intensity */
497 width = ir->userreal3; /* Width of the peak (in time) */
498 rho = ir->userreal4; /* Diameter of the focal spot (nm) */
499 ionize_seed = ir->userint1; /* Random seed for stochastic ionization */
500 ephot = ir->userint2; /* Energy of the photons */
501 mode = ir->userint3; /* Mode of ionizing */
502 interval = 0.001*ir->userint4; /* Interval between pulses (ps) */
504 if ((width <= 0) || (nphot <= 0))
505 fatal_error(0,"Your parameters for ionization are not set properly\n"
506 "width (userreal3) = %f, nphot (userreal2) = %f",
509 if ((mode < 0) || (mode >= eionNR))
510 fatal_error(0,"Ionization mode (userint3)"
511 " should be in the range 0 .. %d",eionNR-1);
515 imax = (nphot/(M_PI*sqr(rho/2)))*1e-10*1.0/(width*sqrt(2.0*M_PI));
518 imax = (nphot/(M_PI*sqr(rho/2)))*1e-10*1.0/(width*sqrt(2.0*M_PI));
521 if (ionize_seed == 0)
522 ionize_seed = make_seed();
524 for(i=0; (i<cr->nodeid); i++)
525 ionize_seed = INT_MAX*rando(&ionize_seed);
526 fprintf(fp,PREFIX"Modifying seed on parallel processor to %d\n",
530 for(Eindex=0; (Eindex < NENER) && (Energies[Eindex] != ephot); Eindex++)
533 fatal_error(0,PREFIX"Energy level of %d keV not supported",ephot);
535 /* Initiate cross section data etc. */
536 ca = mk_cross_atom(fp,md,atomname,Eindex);
540 inv_nratoms = 1.0/md->nr;
542 xvg = xvgropen("ionize.xvg","Ionization Events","Time (ps)","()");
543 xvgr_legend(xvg,asize(leg),leg);
544 ion = ffopen("ionize.log","w");
546 fprintf(fp,PREFIX"Parameters for ionization events:\n");
547 fprintf(fp,PREFIX"Imax = %g, t0 = %g, width = %g, seed = %d\n"
548 PREFIX"# Photons = %g, rho = %g, ephot = %d (keV), Electrons = %s\n",
549 imax,t0,width,ionize_seed,nphot,rho,ephot,yesno_names[bElectron]);
550 fprintf(fp,PREFIX"Electron_mass: %10.3e(keV) Atomic_mass: %10.3e(keV)\n"
551 PREFIX"Speed_of_light: %10.3e(nm/ps)\n",
552 ELECTRONMASS_keV,ATOMICMASS_keV,SPEED_OF_LIGHT);
553 fprintf(fp,PREFIX"Interval between shots: %g ps\n",interval);
554 fprintf(fp,PREFIX"Eindex = %d\n",Eindex);
555 fprintf(fp,PREFIX"Doing ionizations for atoms %d - %d\n",start,end);
557 bElectron = analyze_electrons(fp,&edb,md->nr,atomname);
564 /******************************************************
566 * H E R E S T A R T S I O N I Z A T I O N
568 ******************************************************/
570 /* Calculate probability */
573 while (t > (tmax+interval*0.5))
575 /* End when t <= t0 + (N+0.5) interval */
577 pt = imax*ir->delta_t*exp(-0.5*sqr((t-tmax)/width));
581 hboxx = 0.5*box[XX][XX];
582 hboxy = 0.5*box[YY][YY];
585 /* Width of gaussian for probability of incoherent scattering */
586 sigmaPincoh = 1/sqrt(44.0);
588 /* Arrays for ionization statistics */
589 snew(nionize,md->nr);
593 /* Loop over atoms */
594 for(i=start; (i<end); i++) {
595 /* Loop over collision types */
597 for(k=0; (k<ecollNR); k++)
598 /* Determine cross section for this collision type */
599 pcoll[k]= pt*xray_cross_section(k,&(ca[i]));
601 /* Total probability of ionisation */
602 ptot = 1 - (1-pcoll[ecollPHOTO])*(1-pcoll[ecollINELASTIC]);
604 fprintf(debug,PREFIX"Ptot = %g, t = %g\n",ptot,t);
606 /* Check whether to ionize this guy */
610 bDOIT = (((rando(&ionize_seed) < ptot) && (ca[i].n < ca[i].z)) &&
611 ((sqr(x[i][XX] - hboxx) + sqr(x[i][YY] - hboxy)) < rho2));
617 fatal_error(0,"Unknown ionization mode %d (%s, line %d)",mode,
624 /* The relative probability for a photoellastic event is given by: */
625 pphot = pcoll[ecollPHOTO]/(pcoll[ecollPHOTO]+pcoll[ecollINELASTIC]);
627 if (rando(&ionize_seed) < pphot)
632 /* If a random number is smaller than the probability for
633 * an L ionization than do that. Note that the probability
634 * may be zero (H, He), but the < instead of <= covers that.
636 nK = number_K(&ca[i]);
637 nL = number_L(&ca[i]);
638 bL = (nK == 0) || ( (nL > 0) && (rando(&ionize_seed) > prob_K(k,&(ca[i]))));
642 /* Select which one to take by yet another random numer */
645 /* Get parameters for photoelestic effect */
646 /* Note that in the article this is called 2 theta */
647 theta = DEG2RAD*gauss(70.0,26.0,&ionize_seed);
648 phi = 2*M_PI*rando(&ionize_seed);
651 E_lost = ephot-recoil[ca[i].z].E_L*(ca[i].n+1);
653 E_lost = ephot-recoil[ca[i].z].E_K;
654 if ((ca[i].z > 2) && (nL > 0))
658 fprintf(debug,"i = %d, nK = %d, nL = %d, bL = %s, bKHole = %s\n",
659 i,nK,nL,BOOL(bL),BOOL(bKHole));
666 /* Compute the components of the velocity vector */
667 factor = ((ELECTRONMASS_keV/(ATOMICMASS_keV*md->massT[i]))*
668 (SPEED_OF_LIGHT*sqrt(2*E_lost/ELECTRONMASS_keV)));
670 /* Subtract momentum of recoiling electron */
671 polar2cart(phi,theta,ddv);
672 for(m=0; (m<DIM); m++)
673 dv[m] -= factor*ddv[m];
675 add_electron(fp,md,&edb,i,x,v,dv,ir->delta_t);
678 pr_rvec(debug,0,"ELL",dv,DIM);
684 case ecollINELASTIC: {
685 real theta,phi,Ebind,Eelec;
688 Ebind = (ca[i].n+1)*recoil[ca[i].z].E_L;
690 Ebind = recoil[ca[i].z].E_K;
691 if ((ca[i].z > 2) && (nL > 0))
694 theta = DEG2RAD*rand_theta_incoh(Eindex,&ionize_seed);
695 Eelec = (sqr(ephot)/512)*(1-cos(2*theta));
696 bIonize = (Ebind <= Eelec);
697 bKHole = bKHole && bIonize;
699 fprintf(debug,PREFIX"Ebind: %g, Eelectron: %g\n",Ebind,Eelec);
701 /* Subtract momentum of recoiling photon */
702 /*phi = 2*M_PI*rando(&ionize_seed);
705 dv[XX] -= factor*cos(phi)*sin(theta);
706 dv[YY] -= factor*sin(phi)*sin(theta);
707 dv[ZZ] -= factor*cos(theta);
710 pr_rvec(debug,0,"INELL",dv,DIM);
715 fatal_error(0,"Ga direct naar de gevangenis. Ga niet langs start");
718 /* First increase the charge */
719 if (ca[i].n < ca[i].z) {
720 md->chargeA[i] += 1.0;
721 md->chargeB[i] += 1.0;
726 fprintf(debug,"Random-dv[%3d] = %10.3e,%10.3e,%10.3e,"
727 " ephot = %d, Elost=%10.3e\n",
728 i,dv[XX],dv[YY],dv[ZZ],ephot,E_lost);
731 /* Now actually add the impulse to the velocities */
732 for(m=0; (m<DIM); m++)
742 /* Now check old event: Loop over k holes! */
744 for (kk = 0; (kk < nkh); kk++)
745 if (khole_decay(fp,&(ca[i]),x,v,i,&ionize_seed,ir->delta_t,
746 bElectron,md,&edb)) {
751 if (debug && (ca[i].n > 0))
752 dump_ca(debug,&(ca[i]),i,__FILE__,__LINE__);
755 /* Sum events for statistics if necessary */
757 gmx_sumi(md->nr,nionize,cr);
758 gmx_sumi(md->nr,nkhole,cr);
759 gmx_sumi(md->nr,ndecay,cr);
760 nbuf[0] = dq; nbuf[1] = nkdecay;
762 dq = nbuf[0]; nkdecay = nbuf[1];
764 /* Now sum global events on this timestep to cumulative numbers */
770 /* Print data to the file that holds ionization events per atom */
771 fprintf(ion,"%12.8f",t);
772 for(i=0; (i<md->nr); i++) {
774 fprintf(ion," I:%d",i+1);
776 fprintf(ion," K:%d",i+1);
778 fprintf(ion," D:%d",i+1);
784 /* Print statictics to file */
785 fprintf(xvg,"%10.5f %10.3e %6d %6d %6d %6d",
786 t,pt,dq,dq_tot,nkdecay,nkd_tot);