void retris(complex *data,complex *a,complex *c, complex *b, complex endl,complex endr, int nx, complex *d) { int ix; complex *e,den; complex *f; e=alloc1complex(nx); f=alloc1complex(nx); e[0]=cdiv(cneg(a[0]),endl); f[0]=cdiv(d[0],endl); for(ix=1;ix<nx-1;++ix){ den=cadd(b[ix],cmul(c[ix],e[ix-1])); e[ix]=cdiv(cneg(a[ix]),den); f[ix]=cdiv(csub(d[ix],cmul(f[ix-1],c[ix])),den); } data[nx-1]=cdiv(csub(d[nx-1],cmul(f[nx-2],c[nx-2])),cadd(endr,cmul(c[nx-2],e[nx-2]))); for(ix=nx-2;ix>-1;--ix) data[ix]=cadd(cmul(data[ix+1],e[ix]),f[ix]); free1complex(e); free1complex(f); return; }
void integ(float **mig,int nz,float dz,int nx,int m,float **migi) /* integration of a two-dimensional array input: mig[nx][nz] two-dimensional array output: migi[nx][nz+2*m] integrated array */ { int nfft, nw, ix, iz, iw; float *amp, dw, *rt; complex *ct; /* Set up FFT parameters */ nfft = npfaro(nz+m, 2 * (nz+m)); if (nfft >= SU_NFLTS || nfft >= 720720) err("Padded nt=%d -- too big", nfft); nw = nfft/2 + 1; dw = 2.0*PI/(nfft*dz); amp = ealloc1float(nw); for(iw=1; iw<nw; ++iw) amp[iw] = 0.5/(nfft*(1-cos(iw*dw*dz))); amp[0] = amp[1]; /* Allocate fft arrays */ rt = ealloc1float(nfft); ct = ealloc1complex(nw); for(ix=0; ix<nx; ++ix) { memcpy(rt, mig[ix], nz*FSIZE); memset((void *) (rt + nz), 0, (nfft-nz)*FSIZE); pfarc(1, nfft, rt, ct); /* Integrate traces */ for(iw=0; iw<nw; ++iw){ ct[iw].i = ct[iw].i*amp[iw]; ct[iw].r = ct[iw].r*amp[iw]; } pfacr(-1, nfft, ct, rt); for (iz=0; iz<m; ++iz) migi[ix][iz] = rt[nfft-m+iz]; for (iz=0; iz<nz+m; ++iz) migi[ix][iz+m] = rt[iz]; } free1float(amp); free1float(rt); free1complex(ct); }
void dipfilt(float k,float dpx, float dt, int np, int nw, int nt, float **div,complex *p,complex *q) /********************************************************************* Jacubowicz filter to apply dip-dependent divergence correction ********************************************************************** Input: k wavenumber dpx dip sampling interval dt time sampling interval np number of reflection slopes nt number of time samples nw number of frequency samples div amplitude table p array[nt] containing input p(t,k) Output: q array[nt] containing divergence corrected output q(t,k) *********************************************************************/ { int ip,iw,it,iwl,iwh; float dw,wny,pscl,fftscl,pmin,ph,pm,pl; complex *kq,*qq; /* allocate workspace */ kq = ealloc1complex(nw); qq = ealloc1complex(nw); dw = 2.0*PI/(nw*dt); wny = PI/dt; pscl = 0.5; fftscl = 1.0/nw; pmin = k/wny; /* initialize qq */ for (iw=0; iw<nw; iw++){ qq[iw].r = 0.0; qq[iw].i = 0.0; } for (ip=np-1; ip>=0; ip--){ ph=dpx*(ip+1); pm=dpx*ip; pl=dpx*(ip-1); /*if (ip==np-1) ph=pm;*/ /* define frequency range */ iwl = k/ph/dw + 1; iwh = k/pl/dw; if (pl<1.01*pmin) iwh = (nw%2 ? (nw-1)/2 : nw/2-1); if (pm<1.01*pmin) iwh = (nw%2 ? nw/2 : nw/2-1); /* sum over frequency */ if (iwh>=iwl){ for (it=0; it<nt; it++){ kq[it].r = p[it].r*div[ip][it]; kq[it].i = p[it].i*div[ip][it]; } for (it=nt; it<nw; it++){ kq[it].r = 0.0; kq[it].i = 0.0; } pfacc(1,nw,kq); /* dip filter positive frequencies */ for (iw=iwl; iw<=iwh; iw++){ qq[iw].r += kq[iw].r*pscl; qq[iw].i += kq[iw].i*pscl; } /* dip filter negative frequencies */ iwl=nw-iwl; iwh=nw-iwh; for (iw=iwh; iw<=iwl; iw++){ qq[iw].r += kq[iw].r*pscl; qq[iw].i += kq[iw].i*pscl; } } if (pm<1.01*pmin) break; } /* Fourier transform w to t */ pfacc(-1,nw,qq); for (it=0; it<nt; it++){ q[it].r = qq[it].r*fftscl; q[it].i = qq[it].i*fftscl; } /* free workspace */ free1complex(kq); free1complex(qq); }
void gazdagvt (float k, int nt, float dt, float ft, int ntau, float dtau, float ftau, float *vt, complex *p, complex *q, float qual, float gainceil) /***************************************************************************** Gazdag's phase-shift zero-offset migration for one wavenumber adapted to v(tau) velocity profile ****************************************************************************** Input: k wavenumber nt number of time samples dt time sampling interval ft first time sample ntau number of migrated time samples dtau migrated time sampling interval ftau first migrated time sample vt velocity v[tau] p array[nt] containing data to be migrated Output: q array[ntau] containing migrated data ******************************************************************************/ { int ntfft,nw,it,itau,iw; float dw,fw,tmax,w,tau,phase,coss, *cumgain, gain, alpha; complex cshift,*pp; /* determine frequency sampling */ ntfft = npfa(nt); nw = ntfft; dw = 2.0*PI/(ntfft*dt); fw = -PI/dt; /* determine maximum time */ tmax = ft+(nt-1)*dt; /* allocate workspace */ pp = alloc1complex(nw); cumgain = alloc1float(nw); for (iw=0; iw<nw; iw++) cumgain[iw] = 1.0; /* pad with zeros and Fourier transform t to w, with w centered */ for (it=0; it<nt; it++) pp[it] = (it%2 ? cneg(p[it]) : p[it]); for (it=nt; it<ntfft; it++) pp[it] = cmplx(0.0,0.0); pfacc(1,ntfft,pp); /* account for non-zero ft and non-zero ftau */ for (itau=0 ; itau < ftau ; itau++){ for (iw=0,w=fw; iw<nw; iw++,w+=dw) { if (w==0.0) w = 1e-10/dt; coss = 1.0-pow(0.5 * vt[itau] * k/w,2.0); if (coss>=pow(ftau/tmax,2.0)) { phase = w*(ft-ftau*sqrt(coss)); cshift = cmplx(cos(phase),sin(phase)); pp[iw] = cmul(pp[iw],cshift); } else { pp[iw] = cmplx(0.0,0.0); } } } /* loop over migrated times tau */ for (itau=0,tau=ftau; itau<ntau; itau++,tau+=dtau) { /* initialize migrated sample */ q[itau] = cmplx(0.0,0.0); /* loop over frequencies w */ for (iw=0,w=fw; iw<nw; iw++,w+=dw) { /* accumulate image (summed over frequency) */ q[itau] = cadd(q[itau],pp[iw]); /* compute cosine squared of propagation angle */ if (w==0.0) w = 1e-10/dt; coss = 1.0-pow(0.5 * vt[itau] * k/w,2.0); /* if wave could have been recorded in time */ if (coss>=pow(tau/tmax,2.0)) { /* extrapolate down one migrated time step */ phase = -w*dtau*sqrt(coss); cshift = cmplx(cos(phase),sin(phase)); /* apply gain until gain ceiling is reached */ if (cumgain[iw] < gainceil) { alpha = w/(2.0*vt[itau]*qual); gain = exp(fabs(0.5*vt[itau]*dtau*alpha)); pp[iw] = cmul(pp[iw],crmul(cshift,gain)); cumgain[iw] *= gain; } else { pp[iw] = cmplx(0.0,0.0); } /* else, if wave couldn't have been recorded in time */ } else { /* zero the wave */ pp[iw] = cmplx(0.0,0.0); } } /* scale accumulated image just as we would for an FFT */ q[itau] = crmul(q[itau],1.0/nw); } /* free workspace */ free1complex(pp); free1float(cumgain); }
int main (int argc, char **argv) { int nt; /* number of time samples */ int nz; /* number of migrated depth samples */ int nx; /* number of horizontal samples */ int nxshot; /* number of shots to be migrated */ /*int nxshot_orig;*/ /* first value of nxshot */ int iz,iw,ix,it; /* loop counters */ int igx; /* integerized gx value */ int ntfft; /* fft size */ int nw,truenw; /* number of wave numbers */ int dip=79; /* dip angle */ float sx,gx; /* x source and geophone location */ float gxmin=0.0,gxmax=0.0; /* x source and geophone location */ float min_sx_gx; /* min(sx,gx) */ float oldgx; /* old gx position */ /* float oldgxmin; */ /* old gx position */ /* float oldgxmax; */ /* old gx position */ float oldsx=0.0; /* old sx position */ int isx=0,nxo; /* index for source and geophone */ int oldisx=0; /* old value of source index */ int oldigx=0; /* old value of integerized gx value */ int ix1,ix2,ix3,ixshot; /* dummy index */ int lpad,rpad; /* padding on both sides of the migrated section */ float *wl=NULL,*wtmp=NULL; float fmax; float f1,f2,f3,f4; int nf1,nf2,nf3,nf4; int ntw; float dt=0.004,dz; /* time and depth sampling interval */ float dw; /* frequency sampling interval */ float fw; /* first frequency */ float w; /* frequency */ float dx; /* spatial sampling interval */ float **p=NULL; /* input data */ float **cresult=NULL; /* output result */ float v1; /* average velocity */ double kz2; float **v=NULL,**vp=NULL;/* pointers for the velocity profile */ complex cshift2; complex *wlsp=NULL; /* complex input,output */ complex **cp=NULL; /* ... */ complex **cp1=NULL; /* ... */ complex **cq=NULL; /* ... */ char *vfile=""; /* name of file containing velocities */ FILE *vfp=NULL; int verbose; /* verbose flag */ /* hook up getpar to handle the parameters */ initargs(argc,argv); requestdoc(1); /* get required parameters */ MUSTGETPARINT("nz",&nz); MUSTGETPARINT("nxo",&nxo); MUSTGETPARFLOAT("dz",&dz); MUSTGETPARSTRING("vfile",&vfile); MUSTGETPARINT("nxshot",&nxshot); /* get optional parameters */ if (!getparfloat("fmax",&fmax)) fmax = 25.0; if (!getparfloat("f1",&f1)) f1 = 10.0; if (!getparfloat("f2",&f2)) f2 = 20.0; if (!getparfloat("f3",&f3)) f3 = 40.0; if (!getparfloat("f4",&f4)) f4 = 50.0; if (!getparint("lpad",&lpad)) lpad=9999; if (!getparint("rpad",&rpad)) rpad=9999; if (!getparint("dip",&dip)) dip=79; if (!getparint("verbose",&verbose)) verbose = 0; /* allocating space */ cresult = alloc2float(nz,nxo); vp = alloc2float(nxo,nz); /* load velicoty file */ vfp=efopen(vfile,"r"); efread(vp[0],FSIZE,nz*nxo,vfp); efclose(vfp); /* zero out cresult array */ memset((void *) cresult[0], 0, nxo*nz*FSIZE); /* save value of nxshot */ /* nxshot_orig=nxshot; */ /* get info from first trace */ if (!gettr(&tr)) err("can't get first trace"); nt = tr.ns; get_sx_gx(&sx,&gx); min_sx_gx = MIN(sx,gx); sx = sx - min_sx_gx; gx = gx - min_sx_gx; /* let user give dt and/or dx from command line */ if (!getparfloat("dt", &dt)) { if (tr.dt) { /* is dt field set? */ dt = ((double) tr.dt)/1000000.0; } else { /* dt not set, assume 4 ms */ dt = 0.004; if(verbose) warn("tr.dt not set, assuming dt=0.004"); } } if (!getparfloat("dx",&dx)) { if (tr.d2) { /* is d2 field set? */ dx = tr.d2; } else { dx = 1.0; if(verbose) warn("tr.d2 not set, assuming d2=1.0"); } } checkpars(); oldisx=0; do { /* begin loop over shots */ /* determine frequency sampling interval*/ ntfft = npfar(nt); nw = ntfft/2+1; dw = 2.0*PI/(ntfft*dt); /* compute the index of the frequency to be migrated*/ fw=2.0*PI*f1; nf1=fw/dw+0.5; fw=2.0*PI*f2; nf2=fw/dw+0.5; fw=2.0*PI*f3; nf3=fw/dw+0.5; fw=2.0*PI*f4; nf4=fw/dw+0.5; /* the number of frequencies to migrated*/ truenw=nf4-nf1+1; fw=0.0+nf1*dw; if(verbose) warn("nf1=%d nf2=%d nf3=%d nf4=%d nw=%d",nf1,nf2,nf3,nf4,truenw); /* allocate space */ wl=alloc1float(ntfft); wlsp=alloc1complex(nw); /* generate the Ricker wavelet */ wtmp=ricker(fmax,dt,&ntw); /* zero out wl[] array */ memset((void *) wl, 0, ntfft*FSIZE); /* CHANGE BY CHRIS STOLK, Dec. 11, 2005 */ /* The next two lines are the old code, */ /* it is erroneous because the peak of */ /* the wavelet occurs at positive time */ /* instead of time zero. */ /* for(it=0;it<ntw;it++) wl[it]=wtmp[it]; */ /* New code: we put in the wavelet in a centered fashion */ for(it=0;it<ntw;it++) wl[(it-ntw/2+ntfft) % ntfft]=wtmp[it]; /* End of new code */ free1float(wtmp); /* fourier transform wl array */ pfarc(-1,ntfft,wl,wlsp); /* allocate space */ p = alloc2float(ntfft,nxo); cq = alloc2complex(nw,nxo); /* zero out p[][] array */ memset((void *) p[0], 0, ntfft*nxo*FSIZE); /* initialize a number of items before looping over traces */ nx = 0; igx=0; oldigx=0; oldsx=sx; oldgx=gx; /* oldgxmax=gxmax; */ /* oldgxmin=gxmin; */ do { /* begin looping over traces within a shot gather */ memcpy( (void *) p[igx], (const void *) tr.data,nt*FSIZE); /* get sx and gx */ get_sx_gx(&sx,&gx); sx = (sx - min_sx_gx); gx = (gx - min_sx_gx); igx = NINT(gx/dx); if (igx==oldigx) warn("repeated igx!!! check dx or scalco value!!!"); oldigx = igx; if(gxmin>gx)gxmin=gx; if(gxmax<gx)gxmax=gx; if(verbose) warn(" inside loop: min_sx_gx %f isx %d igx %d gx %f sx %f",min_sx_gx,isx,igx,gx,sx); /* sx, gx must increase monotonically */ if (!(oldsx <= sx) ) err("sx field must be monotonically increasing!"); if (!(oldgx <= gx) ) err("gx field must be monotonically increasing!"); ++nx; } while(gettr(&tr) && sx==oldsx); isx=NINT(oldsx/dx); ixshot=isx; if (isx==oldisx) warn("repeated isx!!! check dx or scalco value!!!"); oldisx=isx; if(verbose) { warn("sx %f, gx %f , gxmin %f gxmax %f nx %d",sx,gx,gxmin,gxmax, nx); warn("isx %d igx %d ixshot %d" ,isx,igx,ixshot); } /* transform the shot gather from time to frequency domain */ pfa2rc(1,1,ntfft,nxo,p[0],cq[0]); /* compute the most left and right index for the migrated */ /* section */ ix1=NINT(oldsx/dx); ix2=NINT(gxmin/dx); ix3=NINT(gxmax/dx); if(ix1>=ix3)ix3=ix1; if(ix1<=ix2)ix2=ix1; ix2-=lpad; ix3+=rpad; if(ix2<0)ix2=0; if(ix3>nxo-1)ix3=nxo-1; /* the total traces to be migrated */ nx=ix3-ix2+1; nw=truenw; /* allocate space for velocity profile within the aperature */ v=alloc2float(nx,nz); for(iz=0;iz<nz;iz++) for(ix=0;ix<nx;ix++) v[iz][ix]=vp[iz][ix+ix2]; /* allocate space */ cp = alloc2complex(nx,nw); cp1 = alloc2complex(nx,nw); /* transpose the frequency domain data from */ /* data[ix][iw] to data[iw][ix] and apply a */ /* Hamming at the same time */ for (ix=0; ix<nx; ++ix) { for (iw=0; iw<nw; iw++){ float tmpp=0.0,tmppp=0.0; if(iw>=(nf1-nf1)&&iw<=(nf2-nf1)){ tmpp=PI/(nf2-nf1); tmppp=tmpp*(iw-nf1)-PI; tmpp=0.54+0.46*cos(tmppp); cp[iw][ix]=crmul(cq[ix+ix2][iw+nf1],tmpp); } else { if(iw>=(nf3-nf1)&&iw<=(nf4-nf1)) { tmpp=PI/(nf4-nf3); tmppp=tmpp*(iw-nf3); tmpp=0.54+0.46*cos(tmppp); cp[iw][ix]=crmul(cq[ix+ix2][iw+nf1],tmpp); } else { cp[iw][ix]=cq[ix+ix2][iw+nf1]; } } cp1[iw][ix]=cmplx(0.0,0.0); } } for(iw=0;iw<nw;iw++){ cp1[iw][ixshot-ix2]=wlsp[iw+nf1]; } if(verbose) { warn("ixshot %d ix %d ix1 %d ix2 %d ix3 %d",ixshot,ix,ix1,ix2,ix3); warn("oldsx %f ",oldsx); } free2float(p); free2complex(cq); free1float(wl); free1complex(wlsp); /* loops over depth */ for(iz=0; iz<nz; ++iz) { /* the imaging condition */ for(ix=0; ix<nx; ++ix){ for(iw=0,w=fw;iw<nw;w+=dw,iw++){ complex tmp; float ratio=10.0; if(fabs(ix+ix2-ixshot)*dx<ratio*iz*dz) tmp=cmul(cp[iw][ix],cp1[iw][ix]); else tmp=cmplx(0.0,0.0); cresult[ix+ix2][iz]+=tmp.r/ntfft; } } /* get the average velocity */ v1=0.0; for(ix=0;ix<nx;++ix) v1+=v[iz][ix]/nx; /* compute time-invariant wavefield */ for(ix=0;ix<nx;++ix) { for(iw=0,w=fw;iw<nw;w+=dw,++iw) { kz2=-(1.0/v1)*w*dz; cshift2=cmplx(cos(kz2),sin(kz2)); cp[iw][ix]=cmul(cp[iw][ix],cshift2); cp1[iw][ix]=cmul(cp1[iw][ix],cshift2); } } /* wave-propagation using finite-difference method */ fdmig(cp, nx, nw,v[iz],fw,dw,dz,dx,dt,dip); fdmig(cp1,nx, nw,v[iz],fw,dw,dz,dx,dt,dip); /* apply thin lens term here */ for(ix=0;ix<nx;++ix) { for(iw=0,w=fw;iw<nw;w+=dw,++iw){ kz2=-(1.0/v[iz][ix]-1.0/v1)*w*dz; cshift2=cmplx(cos(kz2),sin(kz2)); cp[iw][ix]=cmul(cp[iw][ix],cshift2); cp1[iw][ix]=cmul(cp1[iw][ix],cshift2); } } } free2complex(cp); free2complex(cp1); free2float(v); --nxshot; } while(nxshot); /* restore header fields and write output */ for(ix=0; ix<nxo; ix++){ tr.ns = nz; tr.d1 = dz; tr.d2 = dx; tr.offset = 0; tr.cdp = tr.tracl = ix; memcpy( (void *) tr.data, (const void *) cresult[ix],nz*FSIZE); puttr(&tr); } return(CWP_Exit()); }
void fdmig( complex **cp, int nx, int nw, float *v,float fw,float dw,float dz,float dx,float dt,int dip) { int iw,ix,step=1; float *s1,*s2,w,coefa[5],coefb[5],v1,vn,trick=0.1; complex cp2,cp3,cpnm1,cpnm2; complex a1,a2,b1,b2; complex endl,endr; complex *data,*d,*a,*b,*c; s1=alloc1float(nx); s2=alloc1float(nx); data=alloc1complex(nx); d=alloc1complex(nx); a=alloc1complex(nx); b=alloc1complex(nx); c=alloc1complex(nx); if(dip==45){ coefa[0]=0.5;coefb[0]=0.25; step=1; } if(dip==65){ coefa[0]=0.478242060;coefb[0]=0.376369527; step=1; } if(dip==79){ coefa[0]=coefb[0]=0.4575; step=1; } if(dip==80){ coefa[1]=0.040315157;coefb[1]=0.873981642; coefa[0]=0.457289566;coefb[0]=0.222691983; step=2; } if(dip==87){ coefa[2]=0.00421042;coefb[2]=0.972926132; coefa[1]=0.081312882;coefb[1]=0.744418059; coefa[0]=0.414236605;coefb[0]=0.150843924; step=3; } if(dip==89){ coefa[3]=0.000523275;coefb[3]=0.994065088; coefa[2]=0.014853510;coefb[2]=0.919432661; coefa[1]=0.117592008;coefb[1]=0.614520676; coefa[0]=0.367013245;coefb[0]=0.105756624; step=4; } if(dip==90){ coefa[4]=0.000153427;coefb[4]=0.997370236; coefa[3]=0.004172967;coefb[3]=0.964827992; coefa[2]=0.033860918;coefb[2]=0.824918565; coefa[1]=0.143798076;coefb[1]=0.483340757; coefa[0]=0.318013812;coefb[0]=0.073588213; step=5; } v1=v[0];vn=v[nx-1]; do { step--; for(iw=0,w=fw;iw<nw;iw++,w+=dw){ if(fabs(w)<=1.0e-10)w=1.0e-10/dt; for(ix=0;ix<nx;ix++){ s1[ix]=(v[ix]*v[ix])*coefb[step]/(dx*dx*w*w)+trick; s2[ix]=-v[ix]*dz*coefa[step]/(w*dx*dx)*0.5; } for(ix=0;ix<nx;ix++){ data[ix]=cp[iw][ix]; } cp2=data[1]; cp3=data[2]; cpnm1=data[nx-2]; cpnm2=data[nx-3]; a1=crmul(cmul(cp2,conjg(cp3)),2.0); b1=cadd(cmul(cp2,conjg(cp2)),cmul(cp3,conjg(cp3))); if(b1.r==0.0 && b1.i==0.0) a1=cwp_cexp(cmplx(0.0,-w*dx*0.5/v1)); else a1=cdiv(a1,b1); if(a1.i>0.0)a1=cwp_cexp(cmplx(0.0,-w*dx*0.5/v1)); a2=crmul(cmul(cpnm1,conjg(cpnm2)),2.0); b2=cadd(cmul(cpnm1,conjg(cpnm1)),cmul(cpnm2,conjg(cpnm2))); if(b2.r==0.0 && b2.i==0.0) a2=cwp_cexp(cmplx(0.0,-w*dx*0.5/vn)); else a2=cdiv(a2,b2); if(a2.i>0.0)a2=cwp_cexp(cmplx(0.0,-w*dx*0.5/vn)); for(ix=0;ix<nx;ix++){ a[ix]=cmplx(s1[ix],s2[ix]); b[ix]=cmplx(1.0-2.0*s1[ix],-2.0*s2[ix]); } for(ix=1;ix<nx-1;ix++){ d[ix]=cadd(cadd(cmul(data[ix+1],a[ix+1]),cmul(data[ix-1],a[ix-1])), cmul(data[ix],b[ix])); } d[0]=cadd(cmul(cadd(b[0],cmul(a[0],a1)),data[0]),cmul(data[1],a[1])); d[nx-1]=cadd(cmul(cadd(b[nx-1],cmul(a[nx-1],a2)),data[nx-1]), cmul(data[nx-2],a[nx-2])); for(ix=0;ix<nx;ix++){ data[ix]=cmplx(s1[ix],-s2[ix]); b[ix]=cmplx(1.0-2.0*s1[ix],2.0*s2[ix]); } endl=cadd(b[0],cmul(data[0],a1)); endr=cadd(b[nx-1],cmul(data[nx-1],a2)); for(ix=1;ix<nx-1;ix++){ a[ix]=data[ix+1]; c[ix]=data[ix-1]; } a[0]=data[1]; c[nx-1]=data[nx-2]; retris(data,a,c,b,endl,endr,nx,d); for(ix=0;ix<nx;ix++){ cp[iw][ix]=data[ix]; } } }while(step); free1complex(data); free1complex(d); free1complex(b); free1complex(c); free1complex(a); free1float(s1); free1float(s2); return; }
int main (int argc, char **argv) { int nt; /* number of time samples */ int nz; /* number of migrated depth samples */ int nx,nxshot; /* number of midpoints,shotgathers, the folds in a shot gather */ int flag=1; /*flag to use ft or meter as the unit*/ int dip=65; /*maximum dip angle to migrate*/ int iz,iw,ix,it,oldsx; /* loop counters*/ int ntfft; /* fft size*/ int nw; /* number of wave numbers */ int mytid,tids[NNTASKS],msgtype,rc,i;/*variable for PVM function*/ int nw1,task; int lpad=9999,rpad=9999; /*zero-traces padded on left and right sides*/ float f1,f2,f3,f4; /*frequencies to build the Hamming window*/ int nf1,nf2,nf3,nf4; /*the index for above frequencies*/ int NTASKS=0; /*number of slave tasks to start*/ char cpu_name[NNTASKS][80]; /*strings to store the computers' name*/ int flag_cpu=0; /*flag to control if using NTASKS variable*/ float sx,gxmin,gxmax; /*location of geophone and receivers*/ int isx,nxo,ifx=0; /*index for geophone and receivers*/ int ix1,ix2,ix3,il,ir; /*dummy index*/ float *wl,*wtmp; /*pointers for the souce function*/ float Fmax=25; /*peak frequency to make the Ricker wavelet*/ int ntw,truenw; /*number of frequencies to be migrated*/ float dt=0.004,dz; /*time, depth sampling interval*/ float ft; /*first time sample*/ float dw; /*frequency sampling interval*/ float fw; /*first frequency*/ float dx; /*spatial sampling interval*/ float **p,**cresult,**result_tmp; /* input, output data*/ float **v; /*double pointer direct to velocity structure*/ complex *wlsp,**cp,**cq,**cq1; /*pointers for internal usage*/ char *vfile=""; /* name of file containing velocities */ char *cpufile=""; /* name of file containing CPU name */ FILE *vfp,*cpu_fp; /* hook up getpar to handle the parameters */ initargs(argc,argv); requestdoc(1); /* get optional parameters */ if (!getparfloat("ft",&ft)) ft = 0.0; if (!getparint("nz",&nz)) err("nz must be specified"); if (!getparfloat("dz",&dz)) err("dz must be specified"); if (!getparstring("vfile", &vfile)) err("vfile must be specified"); if (!getparint("nxo",&nxo)) err("nxo must be specified"); if (!getparint("nxshot",&nxshot)) err("nshot must be specified"); if (!getparfloat("Fmax",&Fmax)) err("Fmax must be specified"); if (!getparfloat("f1",&f1)) f1 = 10.0; if (!getparfloat("f2",&f2)) f2 = 20.0; if (!getparfloat("f3",&f3)) f3 = 40.0; if (!getparfloat("f4",&f4)) f4 = 50.0; if (!getparint("lpad",&lpad)) lpad=9999; if (!getparint("rpad",&rpad)) rpad=9999; if (!getparint("flag",&flag)) flag=1; if (!getparint("dip",&dip)) dip=65; if (getparstring("cpufile", &cpufile)){ cpu_fp=fopen(cpufile,"r"); NTASKS=0; while(!feof(cpu_fp)){ fscanf(cpu_fp,"%s",cpu_name[NTASKS]); NTASKS++; } NTASKS-=1; flag_cpu=1; } else /*if cpufile not specified, the use NTASKS*/ if (!getparint("NTASKS",&NTASKS)) err("No CPUfile specified, NTASKS must be specified"); /*allocate space for the velocity profile*/ tshot=nxshot; v=alloc2float(nxo,nz); /*load velicoty file*/ vfp=efopen(vfile,"r"); efread(v[0],FSIZE,nz*nxo,vfp); efclose(vfp); /*PVM communication starts here*/ mytid=pvm_mytid(); /*get my pid*/ task=NTASKS; warn("\n %d",task); rc=0; /*spawn slave processes here*/ if(!flag_cpu){ rc=pvm_spawn(child,NULL,PvmTaskDefault,"",task,tids); } else{ for(i=0;i<NTASKS;i++){ rc+=pvm_spawn(child,NULL,PvmTaskHost,cpu_name[i],1,&tids[i]); } } /*show the pid of slaves if*/ for(i=0;i<NTASKS;i++){ if(tids[i]<0)warn("\n %d",tids[i]); else warn("\nt%x\t",tids[i]); } /*if not all the slaves start, then quit*/ if(rc<NTASKS){ warn("error");pvm_exit();exit(1);} /*broadcast the global parameters nxo,nz,dip to all slaves*/ pvm_initsend(PvmDataDefault); rc=pvm_pkint(&nxo,1,1); rc=pvm_pkint(&nz,1,1); rc=pvm_pkint(&dip,1,1); msgtype=PARA_MSGTYPE; task=NTASKS; rc=pvm_mcast(tids,task,msgtype); /*broadcast the velocity profile to all slaves*/ pvm_initsend(PvmDataDefault); rc=pvm_pkfloat(v[0],nxo*nz,1); msgtype=VEL_MSGTYPE; rc=pvm_mcast(tids,task,msgtype); /*free the space for velocity profile*/ free2float(v); /*loop over shot gathers begin here*/ loop: /* get info from first trace */ if (!gettr(&tr)) err("can't get first trace"); nt = tr.ns; /* let user give dt and/or dx from command line */ if (!getparfloat("dt", &dt)) { if (tr.dt) { /* is dt field set? */ dt = ((double) tr.dt)/1000000.0; } else { /* dt not set, assume 4 ms */ dt = 0.004; warn("tr.dt not set, assuming dt=0.004"); } } if (!getparfloat("dx",&dx)) { if (tr.d2) { /* is d2 field set? */ dx = tr.d2; } else { dx = 1.0; warn("tr.d2 not set, assuming d2=1.0"); } } sx=tr.sx; isx=sx/dx; gxmin=gxmax=tr.gx; oldsx=sx; /* determine frequency sampling interval*/ ntfft = npfar(nt); nw = ntfft/2+1; dw = 2.0*PI/(ntfft*dt); /*compute the index of the frequency to be migrated*/ fw=2.0*PI*f1; nf1=fw/dw+0.5; fw=2.0*PI*f2; nf2=fw/dw+0.5; fw=2.0*PI*f3; nf3=fw/dw+0.5; fw=2.0*PI*f4; nf4=fw/dw+0.5; /*the number of frequency to migrated*/ truenw=nf4-nf1+1; fw=0.0+nf1*dw; warn("nf1=%d nf2=%d nf3=%d nf4=%d nw=%d",nf1,nf2,nf3,nf4,truenw); fw=0.0; /* allocate space */ wl=alloc1float(ntfft); wlsp=alloc1complex(nw); /*generate the Ricker wavelet*/ wtmp=ricker(Fmax,dt,&ntw); for(it=0;it<ntfft;it++) wl[it]=0.0; for(it=0;it<ntw-12;it++) wl[it]=wtmp[it+12]; free1float( wtmp); /*Fourier transform the Ricker wavelet to frequency domain*/ pfarc(-1,ntfft,wl,wlsp); /* allocate space */ p = alloc2float(ntfft,nxo); cp = alloc2complex(nw,nxo); for (ix=0; ix<nxo; ix++) for (it=0; it<ntfft; it++) p[ix][it] = 0.0; /*read in a single shot gather*/ ix=tr.gx/dx; memcpy( (void *) p[ix], (const void *) tr.data,nt*FSIZE); nx = 0; while(gettr(&tr)){ int igx; if(tr.sx!=oldsx){ fseek(stdin,(long)(-240-nt*4),SEEK_CUR); break;} igx=tr.gx/dx; memcpy( (void *) p[igx], (const void *) tr.data,nt*FSIZE); if(gxmin>tr.gx)gxmin=tr.gx; if(gxmax<tr.gx)gxmax=tr.gx; nx++; oldsx=tr.sx; } warn("\nnx= %d",nx); warn("sx %f , gxmin %f gxmax %f",sx,gxmin,gxmax); /*transform the shot gather from time to frequency domain*/ pfa2rc(1,1,ntfft,nxo,p[0],cp[0]); /*compute the most left and right index for the migrated section*/ ix1=sx/dx; ix2=gxmin/dx; ix3=gxmax/dx; if(ix1>=ix3)ix3=ix1; if(ix1<=ix2)ix2=ix1; il=ix2; ir=ix3; ix2-=lpad; ix3+=rpad; if(ix2<0)ix2=0; if(ix3>nxo-1)ix3=nxo-1; /*the total traces to be migrated*/ nx=ix3-ix2+1; /*allocate space*/ cq = alloc2complex(nx,nw); cq1 = alloc2complex(nx,nw); /*transpose the frequency domain data from data[ix][iw] to data[iw][ix] and apply a Hamming at the same time*/ for (ix=0; ix<nx; ix++) for (iw=0; iw<nw; iw++){ float tmpp=0.0,tmppp=0.0; if(iw<nf1||iw>nf4) cq[iw][ix]=cmplx(0.0,0.0); else{ if(iw>=nf1&&iw<=nf2){tmpp=PI/(nf2-nf1);tmppp=tmpp*(iw-nf1)-PI;tmpp=0.54+0.46*cos(tmppp); cq[iw][ix]=crmul(cp[ix+ix2][iw],tmpp);} else{ if(iw>=nf3&&iw<=nf4){tmpp=PI/(nf4-nf3);tmppp=tmpp*(iw-nf3);tmpp=0.54+0.46*cos(tmppp); cq[iw][ix]=crmul(cp[ix+ix2][iw],tmpp);} else {cq[iw][ix]=cp[ix+ix2][iw];} } } cq[iw][ix]=cp[ix+ix2][iw]; cq1[iw][ix]=cmplx(0.0,0.0); } ix=sx/dx-ifx; warn("ix %d",ix); for(iw=0;iw<nw;iw++){ cq1[iw][ix-ix2]=wlsp[iw]; } free2float(p); free2complex(cp); free1float(wl); free1complex(wlsp); /*if the horizontal spacing interval is in feet, convert it to meter*/ if(!flag) dx*=0.3048; /*start of the timing function*/ time(&t1); /* send local parameters to all slaves*/ pvm_initsend(PvmDataDefault); ix=15; rc=pvm_pkint(&ix,1,1); rc=pvm_pkint(&ntfft,1,1); rc=pvm_pkint(&ix2,1,1); rc=pvm_pkint(&ix3,1,1); rc=pvm_pkint(&isx,1,1); rc=pvm_pkint(&il,1,1); rc=pvm_pkint(&ir,1,1); rc=pvm_pkfloat(&dx,1,1); rc=pvm_pkfloat(&dz,1,1); rc=pvm_pkfloat(&dw,1,1); rc=pvm_pkfloat(&dt,1,1); msgtype=PARA_MSGTYPE; task=NTASKS; rc=pvm_mcast(tids,task,msgtype); /* send all the frequency to slaves*/ count=NTASKS*5; /*count is the number of frequency components in a shot gather*/ nw=truenw; nw1=nw/(count); if(nw1==0)nw1=1; total=count=ceil(nw*1.0/nw1); /* if it is the first shot gather, send equal data to all the slaves, then for the following shot gathers, only send data when slave requests*/ if(nxshot==tshot){ for(i=0;i<NTASKS;i++){ float *tmpp; float fw1; int nww,byte,nwww; pvm_initsend(PvmDataDefault); nww=nf1+i*nw1;fw1=fw+nww*dw; nwww=nw1; byte=UnDone; rc=pvm_pkint(&byte,1,1); rc=pvm_pkfloat(&fw1,1,1); rc=pvm_pkint(&nwww,1,1); rc=pvm_pkfloat((float *)cq[nww],nx*nwww*2,1); rc=pvm_pkfloat((float *)cq1[nww],nx*nwww*2,1); msgtype=DATA_MSGTYPE; pvm_send(tids[i],msgtype); } count-=NTASKS; } while(count){ int tid0,bufid; float *tmpp; float fw1; int nww,byte,nwww; int i; i=total-count; msgtype=COM_MSGTYPE; bufid=pvm_recv(-1,msgtype); rc=pvm_upkint(&tid0,1,1); pvm_freebuf(bufid); pvm_initsend(PvmDataDefault); nww=nf1+i*nw1;fw1=fw+nww*dw; if(i==total-1)nwww=nw-nw1*i; else nwww=nw1; byte=UnDone; rc=pvm_pkint(&byte,1,1); rc=pvm_pkfloat(&fw1,1,1); rc=pvm_pkint(&nwww,1,1); rc=pvm_pkfloat((float *)cq[nww],nx*nwww*2,1); rc=pvm_pkfloat((float *)cq1[nww],nx*nwww*2,1); msgtype=DATA_MSGTYPE; pvm_send(tid0,msgtype); count--; } ix=Done; pvm_initsend(PvmDataDefault); rc=pvm_pkint(&ix,1,1); msgtype=DATA_MSGTYPE; pvm_mcast(tids,task,msgtype); free2complex(cq); free2complex(cq1); time(&t2); warn("\n %d shot been finished in %f seconds, Ntask=%d",nxshot,difftime(t2,t1),NTASKS); nxshot--; if(nxshot)goto loop; /*when all the shot gathers done, send signal to all slaves to request the partial imaging*/ ix=FinalDone; pvm_initsend(PvmDataDefault); rc=pvm_pkint(&ix,1,1); msgtype=PARA_MSGTYPE; pvm_mcast(tids,task,msgtype); /*allocate space for the final image*/ cresult = alloc2float(nz,nxo); for(ix=0;ix<nxo;ix++) for(iz=0;iz<nz;iz++) { cresult[ix][iz]=0.0; } result_tmp= alloc2float(nz,nxo); /*receive partial image from all the slaves*/ msgtype=RESULT_MSGTYPE; i=0; while(i<NTASKS){ int bufid; bufid=pvm_recv(-1,msgtype); rc=pvm_upkfloat(result_tmp[0],nxo*nz,1); pvm_freebuf(bufid); for(ix=0;ix<nxo;ix++) for(iz=0;iz<nz;iz++) { cresult[ix][iz]+=result_tmp[ix][iz]; } i=i+1; warn("\n i=%d been received",i); } /*send signal to all slaves to kill themselves*/ pvm_initsend(PvmDataDefault); pvm_mcast(tids,task,COM_MSGTYPE); /*output the final image*/ for(ix=0; ix<nxo; ix++){ tr.ns = nz ; tr.dt = dz*1000000.0 ; tr.d2 = dx; tr.offset = 0; tr.cdp = tr.tracl = ix; memcpy( (void *) tr.data, (const void *) cresult[ix],nz*FSIZE); puttr(&tr); } pvm_exit(); return EXIT_SUCCESS; }
int SeisPipe2D(DsuTask *zz) { int nt,nx,nz,nw,ntpad,ntfft; int it,ix,iz,izz,iw,iw0,iw1,iw2,iw3,iwmin,iwmax; int nfreqs,verbose; float dt,dx,dy,dz,dw; float freqs[4],fw,w,scale,fftscl; float *p, **v, *wdxov,*sx; complex *cpx; float **qx; void *TabInfo; eTable *et; char msg[80]; int info, ToTid, MasterTid; int sz, pz, pei; int SeisIntPars[20]; float SeisFloPars[20]; /* Receive process control information */ MsgLog(zz, "Receiving Control info ... " ); MasterTid = RecvInt(SeisIntPars, 2, -1, MsgCntl); pei = SeisIntPars[0]; ToTid = SeisIntPars[1]; MsgLog(zz, " Ready \n"); /* Receive: efile and other pars ... */ MsgLog(zz, "Receiving parameters ..." ); TabInfo = RecvBytes(-1, MsgTable); RecvFI(SeisFloPars, 10, SeisIntPars, 10, -1, -1); MsgLog(zz, " Ready \n" ); /* get integer parameters */ nt = SeisIntPars[0]; nx = SeisIntPars[1]; nz = SeisIntPars[2]; ntpad = SeisIntPars[3]; verbose = SeisIntPars[4]; sz = SeisIntPars[5]; pz = SeisIntPars[6]; /* get Floating point parameters */ dt = SeisFloPars[0]; dx = SeisFloPars[1]; dz = SeisFloPars[2]; freqs[0] = SeisFloPars[3]; freqs[1] = SeisFloPars[4]; freqs[2] = SeisFloPars[5]; freqs[3] = SeisFloPars[6]; sz = nz / pz; if (pei == (pz - 1)) sz += nz % pz; sprintf(msg, "Receiving Velocity info (pei = %d, sz = %d) ... ", pei, sz); MsgLog(zz,msg); v = alloc2float(nx, sz); for (iz=0; iz<sz; ++iz) RecvFloat(v[iz], nx, -1, MsgVel); MsgLog(zz, " Ready \n" ); /* determine frequency w sampling */ ntfft = npfar(nt+ntpad); nw = ntfft/2+1; dw = 2.0*PI/(ntfft*dt); iwmin = MAX(0,MIN(nw-1,NINT(2.0*PI*freqs[0]/dw))); iwmax = MAX(0,MIN(nw-1,NINT(2.0*PI*freqs[3]/dw))); /* read extrapolator table */ et = ezread(TabInfo); /* pret(zz -> fp_log, et); */ /* allocate workspace */ MsgLog(zz, "Allocating space ... "); qx = alloc2float(nx,sz); sx = alloc1float(nx); wdxov = alloc1float(nx); cpx = alloc1complex(nx); MsgLog(zz, " Ready \n"); sprintf(msg, "Process (%d) starting loop on depth steps(%d,%d)\n", pei, pei*(nz/pz), pei*(nz/pz) + sz); MsgLog(zz, msg); /* Cleanup qx */ for (iz=0; iz<sz; ++iz) for (ix=0; ix<nx; ++ix) qx[iz][ix] = 0.0; /* loop over frequencies w */ for (iw=iwmin,w=iwmin*dw; iw<iwmax; ++iw,w+=dw) { if (verbose && !(iw%1)) { sprintf(msg, "\t%d/%d\n",iw,iwmax); MsgLog(zz, msg); } /* load wavefield */ RecvCplx(cpx, nx, -1, MsgSlice); /* loop over depth steps nz */ for (iz=0; iz<sz; ++iz) { /* compute 2.0*dx/v(x) and zero migrated data */ for (ix=0; ix<nx; ix++) sx[ix] = 2.0*dx/v[iz][ix]; /* make w*dx/v(x) */ for (ix=0; ix<nx; ++ix) wdxov[ix] = w*sx[ix]; /* extrapolate wavefield */ etextrap(et,nx,wdxov,cpx); /* accumulate migrated data */ for (ix=0; ix<nx; ++ix) qx[iz][ix] += cpx[ix].r; } /* Send down the wavefield */ if ( pei != (pz -1)) SendCplx(cpx, nx, ToTid, MsgSlice); } /* End of loop for iw */ for (iz=0; iz<sz; iz++) { izz = pei*(nz/pz) + iz; if (verbose) { sprintf(msg, "Sending values for iz = %d\n",izz); MsgLog(zz, msg); } SendFI(qx[iz], nx, &izz, 1, MasterTid, MsgDepth); } sprintf(msg, "End of processing for (%d)\n",pei); MsgLog(zz, msg); /* free workspace */ free1float(sx); free1float(wdxov); free2float(v); free2float(qx); free1complex(cpx); pvm_exit(); return(0); }
/* September 1995 */ #include "stratInv.h" segy tr; /* reading data */ void inputData(char* dataFile) { /* declaration of variables */ int iS, iR, iF, iF1, iF2; /* generic counters */ int ns; /* # of samples */ int wL; /* window length */ float *buffer = NULL; /* to input data */ float window; /* windowing purposes */ complex *bufferC = NULL; /* to Fourier transform the input data */ FILE *fp; /* input file */ /* memory for bufferC */ bufferC = alloc1complex(info->nSamples / 2 + 1); fp = fopen(dataFile,"r"); if (fp == NULL) err("Can't open input data file!\n"); for (iR = 0; iR < info->nR; iR++) { fgettr(fp, &tr); ns = tr.ns; /* DD fprintf(stderr, "ns %d\n", ns);*/ /* allocating memory */ if (iR == 0) buffer = alloc1float(MAX(ns, info->nSamples)); /* reseting */ for (iS = 0; iS < MAX(ns, info->nSamples); iS++) buffer[iS] = 0; memcpy(buffer, tr.data, ns * FSIZE); /* buffer -> dataObs and compensating for complex frequency */ for (iS = 0; iS < info->nSamples; iS++) { buffer[iS] *= exp(-info->tau * iS * dt); /* DD fprintf(stderr, "buffer[%d] : %f\n", iS, buffer[iS]);*/ } /* going to the Fourier domain */ pfarc(-1, info->nSamples, buffer, bufferC); /* windowing (PERC_WINDOW) spectrum */ iF1 = NINT(info->f1 / info->dF); iF2 = NINT(info->f2 / info->dF); wL = info->nF * PERC_WINDOW / 2; wL = 2 * wL + 1; for (iS = 0, iF = 0; iF < info->nSamples / 2 + 1; iF++) { window = 0; if (iF < iF1 || iF >= iF2) { bufferC[iF] = cmplx(0, 0); } else if (iF - iF1 < (wL - 1) / 2) { window = .42 - .5 * cos(2 * PI * (float) iS / ((float) (wL - 1))) + .08 * cos(4 * PI * (float) iS / ((float) (wL - 1))); bufferC[iF].r *= window; bufferC[iF].i *= window; iS++; } else if (iF - iF1 >= info->nF - (wL - 1) / 2) { iS++; window = .42 - .5 * cos(2 * PI * (float) iS / ((float) (wL - 1))) + .08 * cos(4 * PI * (float) iS / ((float) (wL - 1))); bufferC[iF].r *= window; bufferC[iF].i *= window; } } /* going back to time domain */ pfacr(1, info->nSamples, bufferC, buffer); /* copying to dataObs within target window and scaling */ for (iF = 0, iS = NINT(t1 / dt); iS <= NINT(t2 / dt); iS++, iF++) { dataObs[iR][iF] = (scaleData * buffer[iS]) / (float) info->nSamples; /* DD fprintf(stderr, "%d %d %f %f %f %f\n", iR, iF, dataObs[iR][iF], info->f1, info->f2, scaleData);*/ } } /* DD fprintf(stderr, "energy %f\n", auxm1 / (nDM * info->nR)); fwrite(&dataObs[0][0], sizeof(float), nDM * info->nR, stdout); exit(-1);*/ /* freeing memory */ free1float(buffer); free1complex(bufferC); fclose(fp); }
void do_minphdec(float *tr,int nt, float *filter,int fnl,int fnr,float prw) { float *rtr; float *rtx; complex *f; complex *w; complex a; int iamp; float amp; float ampm=-1.0e+20; float amps; float *am; float *ph; float mean=0.0; float sum=0.0; int nfftc; int nf; int i,j; /* counter */ float snfftc; /* Set up pfa fft */ nfftc = npfao(nt,LOOKFAC*nt); if (nfftc >= SU_NFLTS || nfftc >= PFA_MAX) err("Padded nt=%d--too big", nfftc); nf = nfftc/2 + 1; snfftc=1.0/nfftc; rtr = ealloc1float(nfftc); rtx = ealloc1float(nf); f = ealloc1complex(nfftc); w = ealloc1complex(nfftc); am = ealloc1float(nf); ph = ealloc1float(nf); /* clean the arrays */ memset( (void *) w, (int) '\0', nfftc*sizeof(complex)); memset( (void *) rtr, (int) '\0', nfftc*FSIZE); /* Cross correlation */ xcor(nt,0,tr,nt,0,tr,nf,0,rtr); /* FFT */ pfarc(1, nfftc,rtr,w); /* stabilize */ for(i=0;i<nf;i++) { am[i] += am[i]*prw; } /* Normalize */ for(i=0;i<nf;i++) { a=w[i]; am[i]= sqrt(a.r*a.r+a.i*a.i); sum += am[i]; if(am[i]!=0) ph[i] = atan2(a.i,a.r); else ph[i]=0; } sum *= 1.0/nf; sum = 1.0/sum; sscal(nf,sum,am,1); /* Smooth the apmlitude spectra */ if(fnl!=0) conv (fnl+fnr+1,-fnl,filter,nf,0,am,nf,0,am); fprintf(stderr," %f\n",sum); for(i=0;i<nf;i++) { w[i].r = am[i]*cos(ph[i]); w[i].i = am[i]*sin(ph[i]); } for(i=nf,j=nf-1;i<nfftc;i++,j--) { w[i].r = am[j]*cos(ph[j]); w[i].i = am[j]*sin(ph[j]); } /* log spectra */ for (i = 0; i < nfftc; ++i) w[i] = crmul(clog(cmul(w[i],conjg(w[i]))),0.5); /* Hilbert transform */ pfacc(-1,nfftc,w); for (i=0; i<nfftc; ++i) { w[i].r *=snfftc; w[i].i *=snfftc; } for(i=1;i<nfftc/2;i++) w[i] = cadd(w[i],w[i]); for(i=nfftc/2;i<nfftc;i++) w[i] = cmplx(0,0); pfacc(1,nfftc,w); /* end of Hilbert transform */ /* exponentiate */ for(i=0;i<nfftc;i++) w[i] = cexp(w[i]); /* inverse filter */ for(i=0;i<nfftc;i++) f[i] = cdiv(cmplx(1.0,0),w[i]); /* Load trace into tr (zero-padded) */ memset( (void *) w, (int) '\0',nfftc*sizeof(complex)); for(i=0;i<nt;i++) w[i].r = tr[i]; /* Trace to frequency domain */ pfacc(1,nfftc,w); /* apply filter */ for(i=0;i<nfftc;i++) w[i] = cmul(w[i],f[i]); /* Time domain */ pfacr(-1, nfftc,w,rtr); for(i=0;i<nt;i++) rtr[i] *=snfftc; memcpy( (void *) tr, (const void *) rtr, nt*FSIZE); free1float(rtr); free1float(am); free1float(ph); free1complex(f); free1complex(w); }
void fdmig( complex **cp, int nx, int nw, float *v,float fw,float dw,float dz,float dx,float dt,float vc,int dip) { int iw,ix; float *p,*s1,*s2,w,coefa,coefb,v1,vn,trick=0.1; complex cp2,cp3,cpnm1,cpnm2; complex a1,a2,b1,b2; complex endl,endr; complex *data,*d,*a,*b,*c; p=alloc1float(nx); s1=alloc1float(nx); s2=alloc1float(nx); data=alloc1complex(nx); d=alloc1complex(nx); a=alloc1complex(nx); b=alloc1complex(nx); c=alloc1complex(nx); for(ix=0;ix<nx;ix++){ p[ix]=vc/v[ix]; p[ix]=(p[ix]*p[ix]+p[ix]+1.0); } if(dip!=65){ coefa=0.5;coefb=0.25; } else { coefa=0.4784689; coefb=0.37607656; } v1=v[0]; vn=v[nx-1]; for(iw=0,w=fw;iw<nw;iw++,w+=dw){ if(fabs(w)<=1.0e-10)w=1.0e-10/dt; for(ix=0;ix<nx;ix++){ s1[ix]=(v[ix]*v[ix])*p[ix]*coefb/(dx*dx*w*w)+trick; s2[ix]=-(1-vc/v[ix])*v[ix]*dz*coefa/(w*dx*dx)*0.5; } for(ix=0;ix<nx;ix++){ data[ix]=cp[iw][ix]; } cp2=data[1]; cp3=data[2]; cpnm1=data[nx-2]; cpnm2=data[nx-3]; a1=crmul(cmul(cp2,conjg(cp3)),2.0); b1=cadd(cmul(cp2,conjg(cp2)),cmul(cp3,conjg(cp3))); if(b1.r==0.0 && b1.i==0.0) a1=cwp_cexp(cmplx(0.0,-w*dx*0.5/v1)); else a1=cdiv(a1,b1); if(a1.i>0.0) a1=cwp_cexp(cmplx(0.0,-w*dx*0.5/v1)); a2=crmul(cmul(cpnm1,conjg(cpnm2)),2.0); b2=cadd(cmul(cpnm1,conjg(cpnm1)),cmul(cpnm2,conjg(cpnm2))); if(b2.r==0.0 && b2.i==0.0) a2=cwp_cexp(cmplx(0.0,-w*dx*0.5/vn)); else a2=cdiv(a2,b2); if(a2.i>0.0) a2=cwp_cexp(cmplx(0.0,-w*dx*0.5/vn)); for(ix=0;ix<nx;ix++){ a[ix]=cmplx(s1[ix],s2[ix]); b[ix]=cmplx(1.0-2.0*s1[ix],-2.0*s2[ix]); } for(ix=1;ix<nx-1;ix++){ d[ix]=cadd(cadd(cmul(data[ix+1],a[ix+1]), cmul(data[ix-1],a[ix-1])), cmul(data[ix],b[ix])); } d[0]=cadd(cmul(cadd(b[0],cmul(a[0],a1)), data[0]),cmul(data[1],a[1])); d[nx-1]=cadd(cmul(cadd(b[nx-1], cmul(a[nx-1],a2)),data[nx-1]), cmul(data[nx-2],a[nx-2])); for(ix=0;ix<nx;ix++){ data[ix]=cmplx(s1[ix],-s2[ix]); b[ix]=cmplx(1.0-2.0*s1[ix],2.0*s2[ix]); } endl=cadd(b[0],cmul(data[0],a1)); endr=cadd(b[nx-1],cmul(data[nx-1],a2)); for(ix=1;ix<nx-1;ix++){ a[ix]=data[ix+1]; c[ix]=data[ix-1]; } a[0]=data[1]; c[nx-1]=data[nx-2]; retris(data,a,c,b,endl,endr,nx,d); for(ix=0;ix<nx;ix++){ cp[iw][ix]=data[ix]; } } free1complex(data); free1float(p); free1complex(d); free1complex(b); free1complex(c); free1complex(a); free1float(s1); free1float(s2); return; }
int main (int argc, char **argv) { int nt; /* number of time samples */ int nz; /* number of migrated depth samples */ int nx; /* number of horizontal samples */ int nxshot; /* number of shots to be migrated */ int iz,iw,ix,it,ik; /* loop counters */ int igx; /* integerized gx value */ int ntfft,nxfft; /* fft size */ int nw,truenw,nk; /* number of wave numbers */ int dip=65; /* dip angle */ int oldigx=0; /* old value of integerized gx value */ int oldisx=0; /* old value of integerized sx value */ float sx,gx; /* x source and geophone location */ float gxmin=0.0,gxmax=0.0; /* x source and geophone location */ float min_sx_gx; /* min(sx,gx) */ float oldgx; /* old gx position */ float oldgxmin; /* old gx position */ float oldgxmax; /* old gx position */ float oldsx=0.0; /* old sx position */ int isx=0,nxo; /* index for source and geophone */ int ix1,ix2,ix3,ixshot,il=0,ir=0; /* dummy index */ int lpad,rpad; /* padding on both sides of the migrated section */ float *wl=NULL,*wtmp=NULL; float fmax; float f1,f2,f3,f4; int nf1,nf2,nf3,nf4; int ntw; float dt=0.004,dz; /* time and depth sampling interval */ float dw,dk; /* wavenumber and frequency sampling interval */ float fw,fk; /* first wavenumber and frequency */ float w,k; /* wavenumber and frequency */ float dx; /* spatial sampling interval */ float **p=NULL; float **cresult=NULL; /* input, output data */ float v1,vmin; double kz1,kz2; double phase1; float **v=NULL; float **vp=NULL; complex cshift1,cshift2; complex *wlsp=NULL; complex **cp=NULL; complex **cp1=NULL; complex **cq=NULL; complex **cq1=NULL; /*complex input,output*/ char *vfile=""; /* name of file containing velocities */ FILE *vfp=NULL; int verbose; /* verbose flag */ /* hook up getpar to handle the parameters */ initargs(argc,argv); requestdoc(1); /* get optional parameters */ MUSTGETPARINT("nz",&nz); MUSTGETPARFLOAT("dz",&dz); MUSTGETPARSTRING("vfile", &vfile); MUSTGETPARINT("nxo",&nxo); MUSTGETPARINT("nxshot",&nxshot); if (!getparfloat("fmax",&fmax)) fmax = 25. ; if (!getparfloat("f1",&f1)) f1 = 10.0; if (!getparfloat("f2",&f2)) f2 = 20.0; if (!getparfloat("f3",&f3)) f3 = 40.0; if (!getparfloat("f4",&f4)) f4 = 50.0; if (!getparint("lpad",&lpad)) lpad=9999; if (!getparint("rpad",&rpad)) rpad=9999; if (!getparint("dip",&dip)) dip=65; if (!getparint("verbose",&verbose)) verbose = 0; /* allocate space */ cresult = alloc2float(nz,nxo); vp=alloc2float(nxo,nz); /* load velocity file */ vfp=efopen(vfile,"r"); efread(vp[0],FSIZE,nz*nxo,vfp); efclose(vfp); /* zero out cresult array */ memset((void *) cresult[0], 0, nxo*nz*FSIZE); if (!gettr(&tr)) err("can't get first trace"); nt = tr.ns; get_sx_gx(&sx,&gx); min_sx_gx = MIN(sx,gx); gxmin=gxmax=gx; erewind(stdin); /* sx = sx - min_sx_gx; gx = gx - min_sx_gx; */ /* let user give dt and/or dx from command line */ if (!getparfloat("dt", &dt)) { if (tr.dt) { /* is dt field set? */ dt = ((double) tr.dt)/1000000.0; } else { /* dt not set, assume 4 ms */ dt = 0.004; warn("tr.dt not set, assuming dt=0.004"); } } if (!getparfloat("dx",&dx)) { if (tr.d2) { /* is d2 field set? */ dx = tr.d2; } else { dx = 1.0; warn("tr.d2 not set, assuming d2=1.0"); } } do { /* begin loop over shots */ /* determine frequency sampling interval*/ ntfft = npfar(nt); nw = ntfft/2+1; dw = 2.0*PI/(ntfft*dt); /* compute the index of the frequency to be migrated */ fw=2.0*PI*f1; nf1=fw/dw+0.5; fw=2.0*PI*f2; nf2=fw/dw+0.5; fw=2.0*PI*f3; nf3=fw/dw+0.5; fw=2.0*PI*f4; nf4=fw/dw+0.5; /* the number of frequencies to migrated */ truenw=nf4-nf1+1; fw=0.0+nf1*dw; if (verbose) warn("nf1=%d nf2=%d nf3=%d nf4=%d nw=%d",nf1,nf2,nf3,nf4,truenw); /* allocate space */ wl=alloc1float(ntfft); wlsp=alloc1complex(nw); /* generate the Ricker wavelet */ wtmp=ricker(fmax,dt,&ntw); /* zero out wl[] array */ memset((void *) wl, 0, ntfft*FSIZE); /* CHANGE BY CHRIS STOLK, Dec. 11, 2005 */ /* The next two lines are the old code, */ /* it is erroneous because the peak of */ /* the wavelet occurs at positive time */ /* instead of time zero. */ for(it=0;it<ntw;it++) wl[it]=wtmp[it]; /* New code: we put in the wavelet in a centered fashion */ /* for(it=0;it<ntw;it++) { wl[(it-ntw/2+ntfft) % ntfft]=wtmp[it]; } */ /* warn("%12i %12f \n",(it-ntw/2+ntfft) % ntfft,wtmp[it]); */ /* End of new code */ free1float(wtmp); /* fourier transform wl array */ pfarc(-1,ntfft,wl,wlsp); /* CS TEST: this was used to output the array wlsp (the wavelet in the frequency domain) to the file CSinfo, no longer needed and commented out */ /* FILE *CSinfo; CSinfo=fopen("CSinfo","w"); fprintf(CSinfo,"ntfft=%10i\n",ntfft); fprintf(CSinfo,"ntw=%10i\n",ntw); for(iw=0;iw<ntfft/2+1;iw++) fprintf(CSinfo,"%12f %12f \n",wlsp[iw].r,wlsp[iw].i); fclose(CSinfo); */ /* conclusion from the analysis of this info: the wavelet (whose fourier transform is in wlsp) is not zero phase!!! so there is a timeshift error!!! Conclusion obtained dec 11 2005 */ /* CS */ /* allocate space */ p = alloc2float(ntfft,nxo); cq = alloc2complex(nw,nxo); /* zero out p[][] array */ memset((void *) p[0], 0, ntfft*nxo*FSIZE); /* initialize a number of items before looping over traces */ nx = 0; if (gx < 0 ) { igx=gx/dx + nxo; } else { igx=gx/dx ; } oldigx=igx; oldsx=sx; oldgx=gx; oldgxmax=gxmax; oldgxmin=gxmin; while(gettr(&tr)) { /* begin looping over traces within a shot gather */ /* get sx and gx */ get_sx_gx(&sx,&gx); /* warn("%d nx=%d", igx, nx); sx = (sx - min_sx_gx); gx = (gx - min_sx_gx); */ if (gx < 0 ) { igx=gx/dx + nxo; } else { igx=gx/dx ; } if (igx==oldigx) warn("repeated igx!!! check dx or scalco value!!!"); oldigx = igx; if(tr.sx!=oldsx){ efseeko(stdin,(off_t)(-240-nt*4),SEEK_CUR); break;} if(gxmin>gx)gxmin=gx; if(gxmax<gx)gxmax=gx; if(verbose) warn(" inside loop: min_sx_gx %f isx %d igx %d gx %f sx %f",min_sx_gx,isx,igx,gx,sx); /* sx, gx must increase monotonically */ if (!(oldsx <= sx) ) err("sx field must be monotonically increasing!"); if (!(oldgx <= gx) ) err("gx field must be monotonically increasing!"); memcpy( (void *) p[igx], (const void *) tr.data,nt*FSIZE); ++nx; } isx=oldsx/dx; if (isx==oldisx) warn("repeated isx!!! check dx or scalco value!!!"); oldisx=isx; ixshot=isx; if(verbose) { warn("sx %f, gx %f , gxmin %f gxmax %f nx %d",sx,gx,gxmin,gxmax, nx); warn("isx %d igx %d ixshot %d" ,isx,igx,ixshot); } /* transform the shot gather from time to frequency domain */ pfa2rc(1,1,ntfft,nxo,p[0],cq[0]); /* compute the most left and right index for the migrated */ /* section */ ix1=oldsx/dx; ix2=gxmin/dx; ix3=gxmax/dx; if(ix1>=ix3)ix3=ix1; if(ix1<=ix2)ix2=ix1; il=ix2; ir=ix3; ix2-=lpad; ix3+=rpad; if(ix2<0)ix2=0; if(ix3>nxo-1)ix3=nxo-1; /* the total traces to be migrated */ nx=ix3-ix2+1; nw=truenw; /* determine wavenumber sampling (for complex to complex FFT) */ nxfft = npfa(nx); nk = nxfft; dk = 2.0*PI/(nxfft*dx); fk = -PI/dx; /* allocate space for velocity profile within the aperature */ v=alloc2float(nx,nz); for(iz=0;iz<nz;iz++) for(ix=0;ix<nx;ix++) v[iz][ix]=vp[iz][ix+ix2]; /* allocate space */ cp = alloc2complex(nx,nw); cp1 = alloc2complex(nx,nw); /* transpose the frequency domain data from */ /* data[ix][iw] to data[iw][ix] and apply a */ /* Hamming at the same time */ for (ix=0; ix<nx; ix++) { for (iw=0; iw<nw; iw++){ float tmpp=0.0,tmppp=0.0; if(iw>=(nf1-nf1)&&iw<=(nf2-nf1)){ tmpp=PI/(nf2-nf1); tmppp=tmpp*(iw-nf1)-PI; tmpp=0.54+0.46*cos(tmppp); cp[iw][ix]=crmul(cq[ix+ix2][iw+nf1],tmpp); } else { if(iw>=(nf3-nf1)&&iw<=(nf4-nf1)){ tmpp=PI/(nf4-nf3); tmppp=tmpp*(iw-nf3); tmpp=0.54+0.46*cos(tmppp); cp[iw][ix]=crmul(cq[ix+ix2][iw+nf1],tmpp); } else { cp[iw][ix]=cq[ix+ix2][iw+nf1];} } cp1[iw][ix]=cmplx(0.0,0.0); } } for(iw=0;iw<nw;iw++){ cp1[iw][ixshot-ix2]=wlsp[iw+nf1]; } if(verbose) { warn("ixshot %d ix %d ix1 %d ix2 %d ix3 %d",ixshot,ix,ix1,ix2,ix3); warn("oldsx %f ",oldsx); } free2float(p); free2complex(cq); free1float(wl); free1complex(wlsp); /* allocating space */ cq=alloc2complex(nxfft,nw); cq1=alloc2complex(nxfft,nw); /* loops over depth */ for(iz=0;iz<nz;++iz){ /* the imaging condition */ for(ix=0;ix<nx;ix++){ for(iw=0,w=fw;iw<nw;w+=dw,iw++){ complex tmp; float ratio=10.0; if(fabs(ix+ix2-ixshot)*dx<ratio*iz*dz) tmp=cmul(cp[iw][ix],cp1[iw][ix]); else tmp=cmplx(0.0,0.0); cresult[ix+ix2][iz]+=tmp.r/ntfft; } } /* get the minimum velocity */ vmin=0; for(ix=il-ix2;ix<=ir-ix2;ix++){ vmin+=1.0/v[iz][ix]/(ir-il+1); } vmin=1.0/vmin; /* compute the shifted wavefield */ for (ik=0;ik<nx;++ik) { for (iw=0; iw<nw; ++iw) { cq[iw][ik] = ik%2 ? cneg(cp[iw][ik]) : cp[iw][ik]; cq1[iw][ik] = ik%2 ? cneg(cp1[iw][ik]) : cp1[iw][ik]; } } /* zero out cq[][] cq1[][] */ for (ik=nx; ik<nk; ++ik) { for (iw=0; iw<nw; ++iw) { cq[iw][ik] = cmplx(0.0,0.0); cq1[iw][ik] = cmplx(0.0,0.0); } } /* FFT to W-K domain */ pfa2cc(-1,1,nk,nw,cq[0]); pfa2cc(-1,1,nk,nw,cq1[0]); v1=vmin; for(ik=0,k=fk;ik<nk;++ik,k+=dk) { for(iw=0,w=fw;iw<nw;++iw,w+=dw){ if(w==0.0)w=1.0e-10/dt; kz1=1.0-pow(v1*k/w,2.0); if(kz1>0.15){ phase1 = -w*sqrt(kz1)*dz/v1; cshift1 = cmplx(cos(phase1), sin(phase1)); cq[iw][ik] = cmul(cq[iw][ik],cshift1); cq1[iw][ik] = cmul(cq1[iw][ik],cshift1); } else { cq[iw][ik] = cq1[iw][ik] = cmplx(0.0,0.0); } } } pfa2cc(1,1,nk,nw,cq[0]); pfa2cc(1,1,nk,nw,cq1[0]); for(ix=0;ix<nx;++ix) { for(iw=0,w=fw;iw<nw;w+=dw,++iw){ float a=0.015,g=1.0; int I=10; if(ix<=I)g=exp(-a*(I-ix)*(I-ix)); if(ix>=nx-I)g=exp(-a*(-nx+I+ix)*(-nx+I+ix)); cq[iw][ix] = crmul( cq[iw][ix],1.0/nxfft); cq[iw][ix] =ix%2 ? cneg(cq[iw][ix]) : cq[iw][ix]; kz2=(1.0/v1-1.0/v[iz][ix])*w*dz; cshift2=cmplx(cos(kz2),sin(kz2)); cp[iw][ix]=cmul(cq[iw][ix],cshift2); cq1[iw][ix] = crmul( cq1[iw][ix],1.0/nxfft); cq1[iw][ix] =ix%2 ? cneg(cq1[iw][ix]) : cq1[iw][ix]; cp1[iw][ix]=cmul(cq1[iw][ix],cshift2); } } } free2complex(cp); free2complex(cp1); free2complex(cq); free2complex(cq1); free2float(v); --nxshot; } while(nxshot); /* restore header fields and write output */ for(ix=0; ix<nxo; ix++){ tr.ns = nz; tr.d1 = dz; tr.d2 = dx; tr.offset = 0; tr.cdp = tr.tracl = ix; memcpy( (void *) tr.data, (const void *) cresult[ix],nz*FSIZE); puttr(&tr); } return(CWP_Exit()); }
float modeling() { /* declaration of variables */ FILE *fp; /* to report results */ int iF, iF1, iR, offset, iT1, iT2, iS, iProc, i, k; /* counters */ int wL; /* window length */ int die; /* die processor flag */ int FReceived; /* number of frequencies processed */ int apl_pid; /* PVM process id control */ int pid; /* process id */ int processControl; /* monitoring PVM start */ int FInfo[2]; /* frequency delimiters */ float wallcpu; /* wall clock time */ float oF; /* value of the objective function */ float residue; /* data residue */ float wdw; /* windowing purposes */ float *buffer, *bufferRCD; /* auxiliary buffers */ /* upgoing waves */ complex **dataS; /* synthethics in the frequency domain */ complex *bufferC; /* auxiliary buffer */ complex **freqPart; /* frequency arrays sent by the slaves */ /* Clean up log files */ CleanLog(); /* Reseting synchronization flags */ for (i = 0; i < nFreqPart; i++) { statusFreq[i][2] = 0; } /* allocating some memory */ dataS = alloc2complex(info->nF, info->nR); buffer = alloc1float(info->nSamples); bufferRCD = alloc1float(info->nSamples); bufferC = alloc1complex(info->nSamples / 2 + 1); freqPart = alloc2complex(info->nFreqProc, info->nR); /* reseting */ for (iF = 0; iF < info->nSamples / 2 + 1; iF++) bufferC[iF] = zeroC; for (iS = 0; iS < info->nSamples; iS++) { buffer[iS] = 0; bufferRCD[iS] = 0; } /* DD fprintf(stderr, "nF -> %d\n", info->nF);*/ fprintf(stderr, "Starting communication with PVM for modeling\n"); /* starting communication with PVM */ if ((apl_pid = pvm_mytid()) < 0) { pvm_perror("Error enrolling master process"); exit(-1); } processControl = CreateSlaves(processes, PROCESS_MODELING, nProc); if (processControl != nProc) { fprintf(stderr,"Problem starting PVM daemons\n"); exit(-1); } /* converting to velocities */ if (IMPEDANCE) { for (i = 0; i < info->nL + 1; i++) { alpha[i] /= rho[i]; beta[i] /= rho[i]; } } /* Broadcasting all processes common information */ BroadINFO(info, 1, processes, nProc, GENERAL_INFORMATION); /* sending all profiles */ BroadFloat(thick, info->nL + 1, processes, nProc, THICKNESS); BroadFloat(rho, info->nL + 1, processes, nProc, DENSITY); BroadFloat(alpha, info->nL + 1, processes, nProc, ALPHAS); BroadFloat(qP, info->nL + 1, processes, nProc, QALPHA); BroadFloat(beta, info->nL + 1, processes, nProc, BETAS); BroadFloat(qS, info->nL + 1, processes, nProc, QBETA); /* sending frequency partitions for each process */ for (iProc = 0; iProc < nProc; iProc++) { FInfo[0] = statusFreq[iProc][0]; FInfo[1] = statusFreq[iProc][1]; if (info->verbose) fprintf(stderr, "Master sending frequencies [%d, %d] out of %d to slave Modeling %d [id:%d]\n", FInfo[0], FInfo[1], info->nF, iProc, processes[iProc]); procInfo[iProc][0] = FInfo[0]; procInfo[iProc][1] = FInfo[1]; SendInt(FInfo, 2, processes[iProc], FREQUENCY_LIMITS); statusFreq[iProc][2] = 1; } /* waiting modelled frequencies */ /* master process will send more frequencies if there's more work to do */ /* measuring elapsed time */ wallcpu = walltime(); /* reseting frequency counter */ FReceived = 0; while (FOREVER) { pid = RecvCplx(freqPart[0], info->nR * info->nFreqProc, -1, FREQUENCY_PARTITION); /* finding the frequency limits of this process */ /* DD fprintf(stderr, "Master finding the frequency limits of this process\n"); */ iProc = 0; while (pid != processes[iProc]) iProc++; /* DD fprintf(stderr, "iProc %d pid %d\n", iProc, pid);*/ /* copying into proper place of the total frequency array */ for (iR = 0; iR < info->nR; iR++) { for (k = 0, i = procInfo[iProc][0]; i <= procInfo[iProc][1]; i++, k++) { dataS[iR][i - initF] = freqPart[iR][k]; } } /* summing frequencies that are done */ FReceived += procInfo[iProc][1] - procInfo[iProc][0] + 1; if (info->verbose) fprintf(stderr, "Master received %d frequencies, remaining %d\n", FReceived, info->nF - FReceived); /* defining new frequency limits */ i = 0; while (i < nFreqPart && statusFreq[i][2]) i++; /* DD fprintf(stderr, "i %d nFreqPart %d\n", i, nFreqPart);*/ if (i < nFreqPart) { /* there is still more work to be done */ /* tell this process to not die */ die = 0; SendInt(&die, 1, processes[iProc], DIE); FInfo[0] = statusFreq[i][0]; FInfo[1] = statusFreq[i][1]; if (info->verbose) fprintf(stderr, "Master sending frequencies [%d, %d] to slave %d\n", FInfo[0], FInfo[1], processes[iProc]); procInfo[iProc][0] = FInfo[0]; procInfo[iProc][1] = FInfo[1]; SendInt(FInfo, 2, processes[iProc], FREQUENCY_LIMITS); statusFreq[i][2] = 1; } else { /* tell this process to die since there is no more work to do */ if (info->verbose) fprintf(stderr, "Master ''killing'' slave %d\n", processes[iProc]); die = 1; SendInt(&die, 1, processes[iProc], DIE); } /* a check to get out the loop */ if (FReceived >= info->nF) break; } /* quitting PVM */ EndOfMaster(); /* getting elapsed time */ wallcpu = walltime() - wallcpu; fprintf(stderr, "Modeling wall clock time = %f seconds\n", wallcpu); /* back to impedances*/ if (IMPEDANCE) { for (i = 0; i < info->nL + 1; i++) { alpha[i] *= rho[i]; beta[i] *= rho[i]; } } /* computing the objective function for the time window */ for (oF = 0, residue = 0, iR = 0; iR < info->nR; iR++) { /* windowing as it was done to the input data */ iT1 = NINT(info->f1 / info->dF); iT2 = NINT(info->f2 / info->dF); wL = info->nF * PERC_WINDOW / 2; wL = 2 * wL + 1; for (iS = 0, iF = 0; iF < info->nSamples / 2 + 1; iF++) { if (iF < iT1 || iF >= iT2) { bufferC[iF] = cmplx(0, 0); } else if (iF - iT1 < (wL - 1) / 2) { wdw = .42 - .5 * cos(2 * PI * (float) iS / ((float) (wL - 1))) + .08 * cos(4 * PI * (float) iS / ((float) (wL - 1))); bufferC[iF].r = dataS[iR][iF - iT1].r * wdw; bufferC[iF].i = dataS[iR][iF - iT1].i * wdw; iS++; } else if (iF - iT1 >= info->nF - (wL - 1) / 2) { iS++; wdw = .42 - .5 * cos(2 * PI * (float) iS / ((float) (wL - 1))) + .08 * cos(4 * PI * (float) iS / ((float) (wL - 1))); bufferC[iF].r = dataS[iR][iF - iT1].r * wdw; bufferC[iF].i = dataS[iR][iF - iT1].i * wdw; } else { bufferC[iF] = dataS[iR][iF - iT1]; } } /* going to time domain */ /* DD fprintf(stderr, "going to time domain \n");*/ pfacr(1, info->nSamples, bufferC, buffer); /* muting ? */ if (MUTE) { for (iS = 0; iS <= NINT(t1Mute[iR] / dt); iS++) { buffer[iS] = 0; } } /* and computing data misfit and likelihood function */ iS = NINT(t1 / dt); for (iT1 = 0; iT1 < nDM; iT1++) { bufferRCD[iT1 + iS] = 0; for (offset = iT1, iT2 = 0; iT2 < nDM; iT2++) { bufferRCD[iT1 + iS] += (buffer[iT2 + iS] - dataObs[iR][iT2]) * CD[offset]; offset += MAX(SGN0(iT1 - iT2) * (nDM - 1 - iT2), 1); } oF += (buffer[iT1 + iS] - dataObs[iR][iT1]) * bufferRCD[iT1 + iS]; residue += (buffer[iT1 + iS] - dataObs[iR][iT1]) * (buffer[iT1 + iS] - dataObs[iR][iT1]); /* DD fprintf(stdout, "%d %f %f %f %f %f %d %f %f\n", nTotalSamples, oF, dt, auxm1, info->tau, residue, iT1, buffer[iT1], dataObs[iR][iT1 - NINT(t1 / dt)]); */ } /* windowing bufferRCD */ iT1 = NINT(t1 / dt); iT2 = NINT(t2 / dt); wL = nDM * PERC_WINDOW / 2; wL = 2 * wL + 1; for (iS = 0, iF = 0; iF < info->nSamples; iF++) { if (iF < iT1 || iF >= iT2) { bufferRCD[iF] = 0; } else if (iF - iT1 < (wL - 1) / 2) { wdw = .42 - .5 * cos(2 * PI * (float) iS / ((float) (wL - 1))) + .08 * cos(4 * PI * (float) iS / ((float) (wL - 1))); bufferRCD[iF] *= wdw; iS++; } else if (iF - iT1 >= nDM - (wL - 1) / 2) { iS++; wdw = .42 - .5 * cos(2 * PI * (float) iS / ((float) (wL - 1))) + .08 * cos(4 * PI * (float) iS / ((float) (wL - 1))); bufferRCD[iF] *= wdw; } } /* going back to Fourier domain */ pfarc(-1, info->nSamples, bufferRCD, bufferC); for (iF1 = 0, iF = NINT(info->f1 / info->dF); iF <= NINT(info->f2 / info->dF); iF++, iF1++) { resCD[iR][iF1] = bufferC[iF]; } } /* considering the .5 factor of the exponent of the Gaussian */ /* and normalizing the likelihood by the number of samples */ oF /= (2 * nTotalSamples); /* freeing some memory */ /* allocating some memory */ free2complex(dataS); free1float(buffer); free1float(bufferRCD); free1complex(bufferC); free2complex(freqPart); /* considering the regularizaton or model covariance term */ if (PRIOR) { auxm1 = 1. / (float) (numberPar * limRange); /* normalization */ for (auxm2 = 0, iF = 0; iF < limRange; iF++) { for (offset = iF, iF1 = 0; iF1 < limRange; iF1++) { if (vpFrechet) { auxm2 += (alpha[iF + lim[0]] - alphaMean[iF + lim[0]]) * CMvP[offset] * auxm1 * (alpha[iF1 + lim[0]] - alphaMean[iF1 + lim[0]]); } if (vsFrechet) { auxm2 += (beta[iF + lim[0]] - betaMean[iF + lim[0]]) * CMvS[offset] * auxm1 * (beta[iF1 + lim[0]] - betaMean[iF1 + lim[0]]); } if (rhoFrechet) { auxm2 += (rho[iF + lim[0]] - rhoMean[iF + lim[0]]) * CMrho[offset] * auxm1 * (rho[iF1 + lim[0]] - rhoMean[iF1 + lim[0]]); } offset += MAX(SGN0(iF - iF1) * (limRange - 1 - iF1), 1); } } } /* getting normalization factor */ fp = fopen("report", "a"); fprintf(fp,"-----------------------\n"); if (modCount == 0) { oFNorm = oF; fprintf(fp,">> Normalization constant for objective function: %f <<\n", oFNorm); } /* normalizing residue */ residue /= (nTotalSamples); if (!DATACOV && noiseVar == 0) noiseVar = residue / 10.; if (PRIOR) { fprintf(fp, "residue at iteration [%d] : Data residue variance %f , Noise variance %f , Likelihood %f , Prior %f\n", modCount, residue, noiseVar, oF / oFNorm, auxm2 / oFNorm); } else { fprintf(fp,"residue at iteration [%d] : Data residue variance %f , Noise variance %f , Likelihood %f , No Prior\n", modCount, residue, noiseVar, oF / oFNorm); } /* checking if we reached noise variance with the data residue */ if (residue / noiseVar <= 1) { /* DATA IS FIT, stop the procedure */ fprintf(fp, "[][][][][][][][][][][][][][][][][][][][]\n"); fprintf(fp, "DATA WAS FIT UP TO 1 VARIANCE!\n"); fprintf(fp, "[][][][][][][][][][][][][][][][][][][][]\n"); exit(0); } /* adding Likelihood and Prior */ if (PRIOR) oF += auxm2 / 2; fprintf(fp,"TOTAL residue at iteration [%d] : %f\n", modCount, oF / oFNorm); fprintf(fp,"-----------------------\n"); fclose(fp); /* returning objective function value */ return(oF / oFNorm); }
int main( int argc, char *argv[] ) { int ntr=0; /* number of traces */ int ntrv=0; /* number of traces */ int ns=0; int nsv=0; float dt; float dtv; cwp_String fs; cwp_String fv; FILE *fps; FILE *fpv; FILE *headerfp; float *data; /* data matrix of the migration volume */ float *vel; /* velocity matrix */ float *velfi; /* velocity function interpolated to ns values*/ float *velf; /* velocity function */ float *vdt; float *ddt; float *ap; /* array of apperture values in m */ float apr; /* array of apperture values in m */ int *apt=NULL; /* array of apperture time limits in mig. gath*/ float r; /* maximum radius with a given apperture */ float ir2; /* r/d2 */ float ir3; /* r/d3 */ float d2; /* spatial sampling int. in dir 2. */ float d3; /* spatial sampling int. in dir 3. */ float **mgd=NULL; /* migration gather data */ float *migt; /* migrated data trace */ int **mgdnz=NULL; /* migration gather data non zero samples*/ float dm; /* migration gather spatial sample int. */ int im; /* number of traces in migration gather */ int *mtnz; /* migrated trace data non zero smaples */ char **dummyi; /* index array that the trace contains zeros only */ float fac; /* velocity scale factor */ int sphr; /* spherical divergence flag */ int imt; /* mute time sample of trace */ float tmp; int imoff; int **igtr=NULL; int nigtr; int n2; int n3; int verbose; /* phase shift filter stuff */ float power; /* power of i omega applied to data */ float amp; /* amplitude associated with the power */ float arg; /* argument of power */ float phasefac; /* phase factor */ float phase; /* phase shift = phasefac*PI */ complex exparg; /* cexp(I arg) */ register float *rt; /* real trace */ register complex *ct; /* complex transformed trace */ complex *filt; /* complex power */ float omega; /* circular frequency */ float domega; /* circular frequency spacing (from dt) */ float sign; /* sign in front of i*omega default -1 */ int nfft; /* number of points in nfft */ int nf; /* number of frequencies (incl Nyq) */ float onfft; /* 1 / nfft */ size_t nzeros; /* number of padded zeroes in bytes */ initargs(argc, argv); requestdoc(1); MUSTGETPARSTRING("fs",&fs); MUSTGETPARSTRING("fv",&fv); MUSTGETPARINT("n2",&n2); MUSTGETPARINT("n3",&n3); MUSTGETPARFLOAT("d2",&d2); MUSTGETPARFLOAT("d3",&d3); if (!getparfloat("dm", &dm)) dm=(d2+d3)/2.0; /* open datafile */ fps = efopen(fs,"r"); fpv = efopen(fv,"r"); /* Open tmpfile for headers */ headerfp = etmpfile(); /* get information from the first data trace */ ntr = fgettra(fps,&tr,0); if(n2*n3!=ntr) err(" Number of traces in file %d not equal to n2*n3 %d \n", ntr,n2*n3); ns=tr.ns; if (!getparfloat("dt", &dt)) dt = ((float) tr.dt)/1000000.0; if (!dt) { dt = .002; warn("dt not set, assumed to be .002"); } /* get information from the first velocity trace */ ntrv = fgettra(fpv,&trv,0); if(ntrv!=ntr) err(" Number of traces in velocity file %d differ from %d \n", ntrv,ntr); nsv=trv.ns; if (!getparfloat("dtv", &dtv)) dtv = ((float) trv.dt)/1000000.0; if (!dtv) { dtv = .002; warn("dtv not set, assumed to be .002 for velocity"); } if (!getparfloat("fac", &fac)) fac=2.0; if (!getparint("verbose", &verbose)) verbose=0; if (!getparint("sphr", &sphr)) sphr=0; if (!getparfloat("apr", &apr)) apr=75; apr*=3.141592653/180; /* allocate arrays */ data = bmalloc(sizeof(float),ns,ntr); vel = bmalloc(sizeof(float),nsv,ntr); velf = ealloc1float(nsv); velfi = ealloc1float(ns); migt = ealloc1float(ns); vdt = ealloc1float(nsv); ddt = ealloc1float(ns); ap = ealloc1float(ns); mtnz = ealloc1int(ns); dummyi = (char **) ealloc2(n2,n3,sizeof(char)); /* Times to do interpolation of velocity from sparse sampling */ /* to fine sampling of the data */ { register int it; for(it=0;it<nsv;it++) vdt[it]=it*dtv; for(it=0;it<ns;it++) ddt[it]=it*dt; } /* Read traces into data */ /* Store headers in tmpfile */ ntr=0; erewind(fps); erewind(fpv); { register int i2,i3; for(i3=0;i3<n3;i3++) for(i2=0;i2<n2;i2++) { fgettr(fps,&tr); fgettr(fpv,&trv); if(tr.trid > 2) dummyi[i3][i2]=1; else dummyi[i3][i2]=0; efwrite(&tr, 1, HDRBYTES, headerfp); bmwrite(data,1,0,i3*n2+i2,ns,tr.data); bmwrite(vel,1,0,i3*n2+i2,nsv,trv.data); } erewind(headerfp); /* set up the phase filter */ power = 1.0;sign = 1.0;phasefac = 0.5; phase = phasefac * PI; /* Set up for fft */ nfft = npfaro(ns, LOOKFAC * ns); if (nfft >= SU_NFLTS || nfft >= PFA_MAX) err("Padded nt=%d -- too big", nfft); nf = nfft/2 + 1; onfft = 1.0 / nfft; nzeros = (nfft - ns) * FSIZE; domega = TWOPI * onfft / dt; /* Allocate fft arrays */ rt = ealloc1float(nfft); ct = ealloc1complex(nf); filt = ealloc1complex(nf); /* Set up args for complex power evaluation */ arg = sign * PIBY2 * power + phase; exparg = cexp(crmul(I, arg)); { register int i; for (i = 0 ; i < nf; ++i) { omega = i * domega; /* kludge to handle omega=0 case for power < 0 */ if (power < 0 && i == 0) omega = FLT_MAX; /* calculate filter */ amp = pow(omega, power) * onfft; filt[i] = crmul(exparg, amp); } } /* set up constants for migration */ if(verbose) fprintf(stderr," Setting up constants....\n"); r=0; for(i3=0;i3<n3;i3++) for(i2=0;i2<n2;i2++) { if(dummyi[i3][i2] < 1) { /* get the velocity function */ bmread(vel,1,0,i3*n2+i2,nsv,velf); /* linear interpolation from nsv to ns values */ intlin(nsv,vdt,velf,velf[0],velf[nsv-1],ns,ddt,velfi); /* Apply scale factor to velocity */ { register int it; for(it=0;it<ns;it++) velfi[it] *=fac; } /* compute maximum radius from apperture and velocity */ { register int it; for(it=0;it<ns;it++) ap[it] = ddt[it]*velfi[it]*tan(apr)/2.0; } tmp = ap[isamax(ns,ap,1)]; if(tmp>r) r=tmp; } } r=MIN(r,sqrt(SQR((n2-1)*d2)+SQR((n3-1)*d3))); ir2 = (int)(2*r/d2)+1; ir3 = (int)(2*r/d3)+1; im = (int)(r/dm)+1; /* allocate migration gather */ mgd = ealloc2float(ns,im); mgdnz = ealloc2int(ns,im); apt = ealloc1int(im); /* set up the stencil for selecting traces */ igtr = ealloc2int(ir2*ir3,2); stncl(r, d2, d3,igtr,&nigtr); if(verbose) { fprintf(stderr," Maximum radius %f\n",r); fprintf(stderr," Maximum offset %f\n", sqrt(SQR((n2-1)*d2)+SQR((n3-1)*d3))); } /* main processing loop */ for(i3=0;i3<n3;i3++) for(i2=0;i2<n2;i2++) { memset( (void *) tr.data, (int) '\0',ns*FSIZE); if(dummyi[i3][i2] < 1) { memset( (void *) mgd[0], (int) '\0',ns*im*FSIZE); memset( (void *) mgdnz[0], (int) '\0',ns*im*ISIZE); /* get the velocity function */ bmread(vel,1,0,i3*n2+i2,nsv,velf); /* linear interpolation from nsv to ns values */ intlin(nsv,vdt,velf,velf[0],velf[nsv-1],ns,ddt,velfi); /* Apply scale factor to velocity */ { register int it; for(it=0;it<ns;it++) velfi[it] *=fac; } /* create the migration gather */ { register int itr,ist2,ist3; for(itr=0;itr<nigtr;itr++) { ist2=i2+igtr[0][itr]; ist3=i3+igtr[1][itr]; if(ist2 >= 0 && ist2 <n2) if(ist3 >= 0 && ist3 <n3) { if(dummyi[ist3][ist2] <1) { imoff = (int) ( sqrt(SQR(igtr[0][itr]*d2) +SQR(igtr[1][itr]*d3))/dm+0.5); bmread(data,1,0,ist3*n2+ist2,ns,tr.data); imoff=MIN(imoff,im-1); { register int it; /* get the mute time for this offset, apperture and velocity */ xindex(ns,ap,imoff*dm,&imt); for(it=imt;it<ns;it++) if(tr.data[it]!=0) { mgd[imoff][it]+=tr.data[it]; mgdnz[imoff][it]+=1; } } } } } } /* normalize the gather */ { register int ix,it; for(ix=0;ix<im;ix++) for(it=0;it<ns;it++) if(mgdnz[ix][it] > 1) mgd[ix][it] /=(float) mgdnz[ix][it]; } memset( (void *) tr.data, (int) '\0',ns*FSIZE); memset( (void *) mtnz, (int) '\0',ns*ISIZE); /* do a knmo */ { register int ix,it; for(ix=0;ix<im;ix++) { /* get the mute time for this offset, apperture and velocity */ xindex(ns,ap,ix*dm,&imt); knmo(mgd[ix],migt,ns,velfi,0,ix*dm,dt,imt,sphr); /* stack the gather */ for(it=0;it<ns;it++) { if(migt[it]!=0.0) { tr.data[it] += migt[it]; mtnz[it]++; } /* tr.data[it] += mgd[ix][it]; */ } } } { register int it; for(it=0;it<ns;it++) if(mtnz[it]>1) tr.data[it] /=(float)mtnz[it]; } /*Do the phase filtering before the trace is released*/ /* Load trace into rt (zero-padded) */ memcpy( (void *) rt, (const void *) tr.data, ns*FSIZE); memset((void *) (rt + ns), (int) '\0', nzeros); pfarc(1, nfft, rt, ct); { register int i; for (i = 0; i < nf; ++i) ct[i] = cmul(ct[i], filt[i]); } pfacr(-1, nfft, ct, rt); memcpy( (void *) tr.data, (const void *) rt, ns*FSIZE); } /* end of dummy if */ /* spit out the gather */ efread(&tr, 1, HDRBYTES, headerfp); puttr(&tr); if(verbose) fprintf(stderr," %d %d\n",i2,i3); } /* end of i2 loop */ } /* end of i3 loop */ /* This should be the last thing */ efclose(headerfp); /* Free memory */ free2int(igtr); free2float(mgd); free2int(mgdnz); free1int(apt); bmfree(data); bmfree(vel); free1float(velfi); free1float(velf); free1float(ddt); free1float(vdt); free1float(ap); free1int(mtnz); free1float(migt); free1float(rt); free1complex(ct); free1complex(filt); free2((void **) dummyi); return EXIT_SUCCESS; }
static void dmooff (float offset, float fmax, int nx, float dx, int nt, float dt, float ft, float *vrms, float **ptx, float gamma, float *boh, float *zoh, int ntable, float s1, float s2, float sign) /***************************************************************************** perform dmo in f-k domain for one offset ****************************************************************************** Input: offset source receiver offset fmax maximum frequency s1 DMO stretch factor s2 DMO stretch factor sign sign of shift nx number of midpoints dx midpoint sampling interval nt number of time samples dt time sampling interval ft first time vrms array[nt] of rms velocities ptx array[nx][nt] for p(t,x), zero-padded for fft (see notes) Output: ptx array[nx][nt] for dmo-corrected p(t,x) ****************************************************************************** Notes: To avoid having to allocate a separate work space that is larger than the array ptx[nx][nt], ptx must be zero-padded to enable Fourier transform from x to k via prime factor FFT. nxpad (nx after zero-padding) can be estimated by nxpad = 2+npfar(nx+(int)(0.5*ABS(offset/dx))); where npfar() is a function that returns a valid length for real-to-complex prime factor FFT. ptx[nx] to ptx[nxfft-1] must point to zeros. ****************************************************************************** Author: Dave Hale, Colorado School of Mines, 08/08/91 *****************************************************************************/ { int nxfft,itmin,nu,nufft,nw,nk,ix,iu,iw,ik,it,iwn, iwmin,iwmax,nupad,ikmax,*itn,lnt; float dw,dk,tcon,wwscl,wwscl2,scale,scales,kmax,lt, amp,phase,fr,fi,pwr,pwi,cphase,sphase,os1,os2, wmin,wmax,fftscl,du,fu,w,k,*uoft,*tofu,g1,h,vmin,hk,t,*bboh=NULL; complex czero=cmplx(0.0,0.0),**ptk,*pu,*pw; /* number of cdps after padding for fft */ nxfft = npfar(nx+(int)(0.5*ABS(offset/dx))); /* get minimum time of first non-zero sample */ for (ix=0,itmin=nt; ix<nx; ++ix) { for (it=0; it<itmin && ptx[ix][it]==0.0; ++it); itmin = it; } /* if all zeros, simply return */ if (itmin>=nt) return; /* make stretch and compress functions t(u) and u(t) */ maketu(offset,itmin,fmax,nt,dt,ft,vrms,&uoft,&nu,&du,&fu,&tofu,&tcon); /* constants depending on gamma, offset, vrms[0], and nt */ g1 = 2.0*sqrt(gamma)/(1.0+gamma); h = offset/2.0; lnt = nt-1; lt = (nt-1)*dt; vmin = 0.5*MIN((1.0+1.0/gamma)*vrms[0],(1.0+gamma)*vrms[0]); /* if gamma != 1, get bboh[] and itn[] for this offset */ if(gamma!=1.0) getbbohitn (offset,itmin,nt,dt,vrms,ntable,boh,zoh, \ gamma,&bboh,&itn); /* inverse of dmo stretch/squeeze factors */ os1 = 1.0/s1; os2 = 1.0/s2; /* maximum DMO shift (in samples) for any wavenumber k */ nupad = 1.5*((s1+s2)/2.0)*tcon/du; /* frequency sampling */ nufft = npfa(nu+nupad); nw = nufft; dw = 2.0*PI/(nufft*du); /* allocate workspace */ pu = pw = ealloc1complex(nufft); /* wavenumber sampling and maximum wavenumber to apply dmo */ nk = nxfft/2+1; dk = 2.0*PI/ABS(nxfft*dx); kmax = PI/ABS(dx); ikmax = NINT(kmax/dk); /* pointers to complex p(t,k) */ ptk = (complex**)ealloc1(nk,sizeof(complex*)); for (ik=0; ik<nk; ++ik) ptk[ik] = (complex*)ptx[0]+ik*nt; /* fft scale factor */ fftscl = (float)nk/(float)(ikmax+1)/(nufft*nxfft); /* Fourier transform p(t,x) to p(t,k) */ pfa2rc(-1,2,nt,nxfft,ptx[0],ptk[0]); /* loop over wavenumbers less than maximum */ for (ik=0,k=0.0; ik<=ikmax; ++ik,k+=dk) { /* apply time vriant linear phase shift */ hk=sign*h*k; for (it=lnt,t=lt; it>=itmin; --it,t-=dt) { /* calculate phase-shift=boh*h*k, h=offset/2 */ if(gamma != 1.0) phase = bboh[it]*hk; else phase = 0.0; /* phase shift p(t,k) */ cphase=cos(phase); sphase=sin(phase); fr = ptk[ik][it].r; fi = ptk[ik][it].i; ptk[ik][it].r = fr*cphase + fi*sphase; ptk[ik][it].i = -fr*sphase + fi*cphase; } /* stretch p(t;k) to p(u) */ ints8c(nt,dt,ft,ptk[ik],czero,czero,nu,tofu,pu); /* pad with zeros and Fourier transform p(u) to p(w) */ for (iu=nu; iu<nufft; ++iu) pu[iu].r = pu[iu].i = 0.0; pfacc(1,nufft,pu); /* minimum and maximum frequencies to process */ wmin = ABS(0.5*vmin*k); wmax = ABS(PI/du); iwmin = MAX(1,MIN((nw-1)/2,NINT(wmin/dw))); iwmax = MAX(0,MIN((nw-1)/2,NINT(wmax/dw))); /* constant independent of w */ wwscl = os1*pow(g1*hk/tcon,2.0); wwscl2 = wwscl*os2/os1; /* zero dc (should be zero anyway) */ pw[0].r = pw[0].i = 0.0; /* zero frequencies below minimum */ for (iw=1,iwn=nw-iw; iw<iwmin; ++iw,--iwn) pw[iw].r = pw[iw].i = pw[iwn].r = pw[iwn].i = 0.0; /* do dmo between minimum and maximum frequencies */ for (iw=iwmin,iwn=nw-iwmin,w=iwmin*dw; iw<=iwmax; ++iw,--iwn,w+=dw) { /* positive w */ scales = 1.0+wwscl/(w*w); scale = sqrt(scales); phase = s1*w*tcon*(scale-1.0); amp = fftscl*(1.0-s1+s1/scale); fr = amp*cos(phase); fi = amp*sin(phase); pwr = pw[iw].r; pwi = pw[iw].i; pw[iw].r = pwr*fr-pwi*fi; pw[iw].i = pwr*fi+pwi*fr; /* negative w */ scales = 1.0+wwscl2/(w*w); scale = sqrt(scales); phase = s2*w*tcon*(scale-1.0); amp = fftscl*(1.0-s2+s2/scale); fr = amp*cos(phase); fi = amp*sin(phase); pwr = pw[iwn].r; pwi = pw[iwn].i; pw[iwn].r = pwr*fr+pwi*fi; pw[iwn].i = pwi*fr-pwr*fi; } /* zero frequencies above maximum to Nyquist (if present) */ for (iw=iwmax+1,iwn=nw-iw; iw<=nw/2; ++iw,--iwn) pw[iw].r = pw[iw].i = pw[iwn].r = pw[iwn].i = 0.0; /* Fourier transform p(w) to p(u) */ pfacc(-1,nufft,pu); /* compress p(u) to p(t;k) */ ints8c(nu,du,fu,pu,czero,czero,nt,uoft,ptk[ik]); } /* zero wavenumber between maximum and Nyquist */ for (; ik<nk; ++ik) for (it=0; it<nt; ++it) ptk[ik][it].r = ptk[ik][it].i = 0.0; /* Fourier transform p(t,k) to p(t,x) */ pfa2cr(1,2,nt,nxfft,ptk[0],ptx[0]); /* free workspace */ free1float(tofu); free1float(uoft); free1complex(pu); free1(ptk); }
static void dmooff (float offset, float fmax, float sdmo, int nx, float dx, int nt, float dt, float ft, float *vrms, float **ptx) /***************************************************************************** perform dmo in f-k domain for one offset ****************************************************************************** Input: offset source receiver offset fmax maximum frequency sdmo DMO stretch factor nx number of midpoints dx midpoint sampling interval nt number of time samples dt time sampling interval ft first time vrms array[nt] of rms velocities ptx array[nx][nt] for p(t,x), zero-padded for fft (see notes) Output: ptx array[nx][nt] for dmo-corrected p(t,x) ****************************************************************************** Notes: To avoid having to allocate a separate work space that is larger than the array ptx[nx][nt], ptx must be zero-padded to enable Fourier transform from x to k via prime factor FFT. nxpad (nx after zero-padding) can be estimated by nxpad = 2+npfar(nx+(int)(0.5*ABS(offset/dx))); where npfar() is a function that returns a valid length for real-to-complex prime factor FFT. ptx[nx] to ptx[nxfft-1] must point to zeros. ****************************************************************************** Author: Dave Hale, Colorado School of Mines, 08/08/91 *****************************************************************************/ { int nxfft,itmin,nu,nufft,nw,nk,ix,iu,iw,ik,it,iwn, iwmin,iwmax,nupad,ikmax; float dw,dk,tcon,wwscl,scale,scales,kmax, amp,phase,fr,fi,pwr,pwi, wmin,wmax,fftscl,du,fu,w,k,osdmo,*uoft,*tofu; complex czero=cmplx(0.0,0.0),**ptk,*pu,*pw; /* number of cdps after padding for fft */ nxfft = npfar(nx+(int)(0.5*ABS(offset/dx))); /* get minimum time of first non-zero sample */ for (ix=0,itmin=nt; ix<nx; ++ix) { for (it=0; it<itmin && ptx[ix][it]==0.0; ++it); itmin = it; } /* if all zeros, simply return */ if (itmin>=nt) return; /* make stretch and compress functions t(u) and u(t) */ maketu(offset,itmin,fmax,nt,dt,ft,vrms,&uoft,&nu,&du,&fu,&tofu,&tcon); /* adjust DMO stretch factor for nominal error in log stretch; */ /* solve sdmo*(sqrt(1-a/sdmo)-1) = 0.5*log(1-a), where a=0.5 */ sdmo *= .62; /* inverse of dmo stretch factor */ osdmo = 1.0/sdmo; /* maximum DMO shift (in samples) for any wavenumber k */ nupad = 1.5*sdmo*tcon/du; /* frequency sampling */ nufft = npfa(nu+nupad); nw = nufft; dw = 2.0*PI/(nufft*du); /* allocate workspace */ pu = pw = ealloc1complex(nufft); /* wavenumber sampling and maximum wavenumber to apply dmo */ nk = nxfft/2+1; dk = 2.0*PI/ABS(nxfft*dx); kmax = PI/ABS(dx); ikmax = NINT(kmax/dk); /* pointers to complex p(t,k) */ ptk = (complex**)ealloc1(nk,sizeof(complex*)); for (ik=0; ik<nk; ++ik) ptk[ik] = (complex*)ptx[0]+ik*nt; /* fft scale factor */ fftscl = (float)nk/(float)(ikmax+1)/(nufft*nxfft); /* Fourier transform p(t,x) to p(t,k) */ pfa2rc(-1,2,nt,nxfft,ptx[0],ptk[0]); /* loop over wavenumbers less than maximum */ for (ik=0,k=0.0; ik<=ikmax; ++ik,k+=dk) { /* stretch p(t;k) to p(u) */ ints8c(nt,dt,ft,ptk[ik],czero,czero,nu,tofu,pu); /* pad with zeros and Fourier transform p(u) to p(w) */ for (iu=nu; iu<nufft; ++iu) pu[iu].r = pu[iu].i = 0.0; pfacc(1,nufft,pu); /* minimum and maximum frequencies to process */ wmin = ABS(0.5*vrms[0]*k); wmax = ABS(PI/du); iwmin = MAX(1,MIN((nw-1)/2,NINT(wmin/dw))); iwmax = MAX(0,MIN((nw-1)/2,NINT(wmax/dw))); /* constant independent of w */ wwscl = osdmo*pow(k*0.5*offset/tcon,2.0); /* zero dc (should be zero anyway) */ pw[0].r = pw[0].i = 0.0; /* zero frequencies below minimum */ for (iw=1,iwn=nw-iw; iw<iwmin; ++iw,--iwn) pw[iw].r = pw[iw].i = pw[iwn].r = pw[iwn].i = 0.0; /* do dmo between minimum and maximum frequencies */ for (iw=iwmin,iwn=nw-iwmin,w=iwmin*dw; iw<=iwmax; ++iw,--iwn,w+=dw) { scales = 1.0+wwscl/(w*w); scale = sqrt(scales); phase = sdmo*w*tcon*(scale-1.0); amp = fftscl*(1.0-sdmo+sdmo/scale); fr = amp*cos(phase); fi = amp*sin(phase); pwr = pw[iw].r; pwi = pw[iw].i; pw[iw].r = pwr*fr-pwi*fi; pw[iw].i = pwr*fi+pwi*fr; pwr = pw[iwn].r; pwi = pw[iwn].i; pw[iwn].r = pwr*fr+pwi*fi; pw[iwn].i = pwi*fr-pwr*fi; } /* zero frequencies above maximum to Nyquist (if present) */ for (iw=iwmax+1,iwn=nw-iw; iw<=nw/2; ++iw,--iwn) pw[iw].r = pw[iw].i = pw[iwn].r = pw[iwn].i = 0.0; /* Fourier transform p(w) to p(u) */ pfacc(-1,nufft,pu); /* compress p(u) to p(t;k) */ ints8c(nu,du,fu,pu,czero,czero,nt,uoft,ptk[ik]); } /* zero wavenumber between maximum and Nyquist */ for (; ik<nk; ++ik) for (it=0; it<nt; ++it) ptk[ik][it].r = ptk[ik][it].i = 0.0; /* Fourier transform p(t,k) to p(t,x) */ pfa2cr(1,2,nt,nxfft,ptk[0],ptx[0]); /* free workspace */ free1float(tofu); free1float(uoft); free1complex(pu); free1(ptk); }