void cbanded_solve (sf_complex *b) /*< multiply by inverse (in place) >*/ { int k, m; sf_complex t; for (k = 1; k < band; k++) { t = b[k]; for (m = 1; m <= k; m++) { #ifdef SF_HAS_COMPLEX_H t -= o[m-1][k-m] * b[k-m]; #else t = sf_csub(t,sf_cmul(o[m-1][k-m],b[k-m])); #endif } b[k] = t; } for (k = band; k < n; k++) { t = b[k]; for (m = 1; m <= band; m++) { #ifdef SF_HAS_COMPLEX_H t -= o[m-1][k-m] * b[k-m]; #else t = sf_csub(t,sf_cmul(o[m-1][k-m],b[k-m])); #endif } b[k] = t; } for (k = n-1; k >= n - band; k--) { #ifdef SF_HAS_COMPLEX_H t = b[k]/d[k]; #else t = sf_crmul(b[k],1./d[k]); #endif for (m = 0; m < n - k - 1; m++) { #ifdef SF_HAS_COMPLEX_H t -= conjf(o[m][k]) * b[k+m+1]; #else t = sf_csub(t,sf_cmul(conjf(o[m][k]),b[k+m+1])); #endif } b[k] = t; } for (k = n - band - 1; k >= 0; k--) { #ifdef SF_HAS_COMPLEX_H t = b[k]/d[k]; #else t = sf_crmul(b[k],1./d[k]); #endif for (m = 0; m < band; m++) { #ifdef SF_HAS_COMPLEX_H t -= conjf(o[m][k]) * b[k+m+1]; #else t = sf_csub(t,sf_cmul(conjf(o[m][k]),b[k+m+1])); #endif } b[k] = t; } }
void ctoeplitz_solve (const sf_complex *r /* top row of the matrix */, sf_complex *f /* inverted in place */) /*< apply the solver >*/ { int i,j; sf_complex e,c,w, bot; float v; v=crealf(r[0]); #ifdef SF_HAS_COMPLEX_H f[0] /= v; #else f[0] = sf_crmul(f[0],1./v); #endif for (j=1; j < n; j++) { e = cdprod(j,a,r); #ifdef SF_HAS_COMPLEX_H c = -e/v; #else c = sf_crmul(e,-1./v); #endif v += crealf(c)*crealf(e) + cimagf(c)*cimagf(e); for (i=1; i <= j/2; i++) { #ifdef SF_HAS_COMPLEX_H bot = a[j-i] + c*conjf(a[i]); a[i] += c*conjf(a[j-i]); #else bot = sf_cadd(a[j-i],sf_cmul(c,conjf(a[i]))); a[i] = sf_cadd(a[i],sf_cmul(c,conjf(a[j-i]))); #endif a[j-i] = bot; } a[j] = c; w = cdprod(j,f,r); #ifdef SF_HAS_COMPLEX_H c = (f[j]-w)/v; #else c = sf_crmul(sf_csub(f[j],w),1./v); #endif for (i=0; i < j; i++) { #ifdef SF_HAS_COMPLEX_H f[i] += c*conjf(a[j-i]); #else f[i] = sf_cadd(f[i],sf_cmul(c,conjf(a[j-i]))); #endif } f[j] = c; } }
void radon_lop (bool adj, bool add, int nm, int nd, sf_complex *mm, sf_complex *dd) /*< linear operator >*/ { int ix, ip; sf_complex c, d; if (nm != np || nd != nx) sf_error("%s: mismatched data sizes",__FILE__); sf_cadjnull(adj, add, nm, nd, mm, dd); for (ix=0; ix < nx; ix++) { if (adj == true) { c = dc[ix]; #ifdef SF_HAS_COMPLEX_H d = dd[ix]*c0[ix]; #else d = sf_cmul(dd[ix],c0[ix]); #endif for (ip=0; ip < np-1; ip++) { #ifdef SF_HAS_COMPLEX_H mm[ip] += d; d *= c; #else mm[ip] = sf_cadd(mm[ip],d); d = sf_cmul(d,c); #endif } #ifdef SF_HAS_COMPLEX_H mm[np-1] += d; #else mm[np-1] = sf_cadd(mm[np-1],d); #endif } else { c = conjf(dc[ix]); d = mm[np-1]; for (ip=np-2; ip >= 0; ip--) { #ifdef SF_HAS_COMPLEX_H d = d*c + mm[ip]; #else d = sf_cadd(sf_cmul(d,c),mm[ip]); #endif } #ifdef SF_HAS_COMPLEX_H dd[ix] += d*conjf(c0[ix]); #else dd[ix] = sf_cadd(dd[ix],sf_cmul(d,conjf(c0[ix]))); #endif } } }
void BornScatteredField(float xx,float xy,float xz,float *output) /*< Born scattering field >*/ { int iw; double omega,scale; sf_complex U,dU; sf_complex val,fkern; ZeroArray(u, nt); ZeroArray(du, nt); *output=0.; for (iw=0; iw<nw; iw++) { omega=ow+dw*iw; scale =cos((SF_PI/2)*(((double) iw+1)/((double) nw+1))); scale*=scale*omega; #ifdef SF_HAS_COMPLEX_H /* background field */ U=Green(xx,xy,xz,sx,sy,sz,omega)*scale; /* scattered field */ val=Green(px,py,pz,sx,sy,sz,omega)*scale; val *= Green(xx,xy,xz,px,py,pz,omega); dU = val*(-omega*omega*dv); U += dU; fkern=cexpf(sf_cmplx(0.,-omega*0.6)); val=fkern*U; #else /* background field */ U=sf_crmul(Green(xx,xy,xz,sx,sy,sz,omega),scale); /* scattered field */ val=sf_crmul(Green(px,py,pz,sx,sy,sz,omega),scale); val=sf_cmul(val,Green(xx,xy,xz,px,py,pz,omega)); dU=sf_crmul(val,-omega*omega*dv); U=sf_cadd(U,dU); fkern=cexpf(sf_cmplx(0.,-omega*0.6)); val=sf_cmul(fkern,U); #endif *output+=crealf(val); } return; }
void radon_toep (sf_complex *qq /* Toeplitz row */, float eps /* regularization */) /*< fill Toeplitz matrix for inversion >*/ { int ix, ip; sf_complex c, q; qq[0] = sf_cmplx(eps*eps,0.); for (ip=1; ip < np; ip++) { qq[ip] = czero; } for (ix=0; ix < nx; ix++) { c = conjf(dc[ix]); q = sf_cmplx(1.,0.); for (ip=0; ip < np-1; ip++) { #ifdef SF_HAS_COMPLEX_H qq[ip] += q; q *= c; #else qq[ip] = sf_cadd(qq[ip],q); q = sf_cmul(q,c); #endif } #ifdef SF_HAS_COMPLEX_H qq[np-1] += q; #else qq[np-1] = sf_cadd(qq[np-1],q); #endif } }
/*------------------------------------------------------------*/ void cam3_ssf( sf_complex w /* frequency */, sf_complex ***wx /* wavefield */, cub3d cub, cam3d cam, tap3d tap, slo3d slo, int imz, int ompith ) /*< Wavefield extrapolation by SSF >*/ { sf_complex co,cs,cr,w2,khy,kss,krr; float s, kmy, d,dsc,drc; int imy,imx,ihx,jmy, js,jr; co = sf_cmplx(0,0); #ifdef SF_HAS_COMPLEX_H w2 = w*w; #else w2 = sf_cmul(w,w); #endif #ifdef SF_HAS_COMPLEX_H LOOP( s = slo->so[ompith][ cam->jy[imy] ][ cam->is[ihx][imx] ] + slo->so[ompith][ cam->jy[imy] ][ cam->ir[ihx][imx] ]; wx[ihx][imy][imx] *= cexpf(-w*s* cub->amz.d/2); );
/*------------------------------------------------------------*/ void sft2(sf_complex **pp) /*< apply shift >*/ { int i1,i2; for (i2=0; i2<n2; i2++) { for(i1=0; i1<n1; i1++) { #ifdef SF_HAS_COMPLEX_H pp[i2][i1] *= shf1[i1]*shf2[i2]; #else pp[i2][i1] = sf_cmul(pp[i2][i1], sf_cmul(shf1[i1],shf2[i2]) ); #endif } } }
/*------------------------------------------------------------*/ void ompsft2(sf_complex **pp, fft2d fft) /*< apply shift >*/ { int i1,i2; for (i2=0; i2<fft->n2; i2++) { for (i1=0; i1<fft->n1; i1++) { #ifdef SF_HAS_COMPLEX_H pp[i2][i1] *= fft->shf1[i1]*fft->shf2[i2]; #else pp[i2][i1] = sf_cmul(pp[i2][i1], sf_cmul(fft->shf1[i1],fft->shf2[i2])); #endif } } }
void fft1_2D_fwd (float *input, float *output, int n2) /*< Fast Fourier Transform along the first axis forward>*/ { /*bool inv, sym, opt;*/ int n1, i1, i2; float *p, wt, shift; kiss_fft_cpx *pp, ce; char *label; sf_file out=NULL; kiss_fftr_cfg cfg; bool verb; verb=1; wt = sym? 1./sqrtf((float) ntopt): 1.0/ntopt; p = sf_floatalloc(ntopt); pp = (kiss_fft_cpx*) sf_complexalloc(nw); cfg = kiss_fftr_alloc(ntopt,0,NULL,NULL); for (i2=0; i2 < n2; i2++) { /*sf_floatread (p,n1,in);*/ memcpy(p,&input[i2*nt],nt*sizeof(float)); if (sym) { for (i1=0; i1 < nt; i1++) { p[i1] *= wt; } } for (i1=nt; i1 < ntopt; i1++) { p[i1]=0.0; // padding } kiss_fftr (cfg,p,pp); if (0. != ot) { for (i1=0; i1 < nw; i1++) { shift = -2.0*SF_PI*i1*dw*ot; ce.r = cosf(shift); ce.i = sinf(shift); pp[i1]=sf_cmul(pp[i1],ce); } } memcpy(&output[i2*2*nw],pp,2*nw*sizeof(float)); } free(p); free(pp); /*exit (0);*/ }
void rweone_phs( sf_complex *v, float w, float a0, float b0 ) /*< Fourier-domain phase shift >*/ { int ig,ikg; float kg; float a2,b2,k2; sf_complex iw,ikt,w2; a2 = a0*a0; b2 = b0*b0; iw = sf_cmplx(2e-3,-w); #ifdef SF_HAS_COMPLEX_H w2 = iw*iw; #else w2 = sf_cmul(iw,iw); #endif rweone_fft(false,(kiss_fft_cpx*) v); for(ig=0;ig<ag.n;ig++) { ikg = KMAP(ig,ag.n); kg = okg + ikg * dkg; k2 = kg*kg; #ifdef SF_HAS_COMPLEX_H ikt = csqrtf( w2*a2 + k2*b2 ); v[ig] *= cexpf(-ikt*at.d); #else ikt = csqrtf(sf_cadd(sf_crmul(w2,a2),sf_cmplx(k2*b2,0.))); v[ig] = sf_cmul(v[ig],cexpf(sf_crmul(ikt,-at.d))); #endif } rweone_fft( true,(kiss_fft_cpx*) v); }
void cmatmul(kiss_fft_cpx *a,kiss_fft_cpx *b,int m,int n,int k,kiss_fft_cpx *c) /*< complex matrix multiplication: a->m*n, b->n*k >*/ { int i,j,l,u; for (i=0; i<=m-1; i++) for (j=0; j<=k-1; j++) { u=i*k+j; c[u].r=0.0; c[u].i=0.0; for (l=0; l<=n-1; l++) c[u]=sf_cadd(c[u],sf_cmul(a[i*n+l],b[l*k+j])); } }
void cicai1_lop (bool adj, bool add, int nx, int ny, sf_complex* xx, sf_complex* yy) /*< linear operator >*/ { int b, x, y; if(ny != nx) sf_error("%s: size problem: %d != %d",__FILE__,ny,nx); sf_cadjnull (adj, add, nx, ny, xx, yy); for( b=0; b < nb; b++) { for( y = nb - lg; y <= ny - lg; y++) { x = y - b + lg - 1; #ifdef SF_HAS_COMPLEX_H if( adj) xx[x] += yy[y] * conjf(bb[b]); else yy[y] += xx[x] * bb[b]; #else if( adj) xx[x] = sf_cadd(xx[x],sf_cmul(yy[y],conjf(bb[b]))); else yy[y] = sf_cadd(yy[y],sf_cmul(xx[x],bb[b])); #endif } } }
void rweone_thr( sf_complex *a, sf_complex *b, sf_complex *c, sf_complex *v, int n) /*< tridiagonal solver >*/ { int i; #ifdef SF_HAS_COMPLEX_H b[0]/=a[0]; v[0]/=a[0]; #else b[0] = sf_cdiv(b[0],a[0]); v[0] = sf_cdiv(v[0],a[0]); #endif for(i=1;i<n;i++) { #ifdef SF_HAS_COMPLEX_H a[i] -= c[i-1]*b[i-1]; if(i<n-1) b[i]/=a[i]; v[i]=( v[i]-c[i-1]*v[i-1] ) / a[i]; #else a[i] = sf_csub(a[i],sf_cmul(c[i-1],b[i-1])); if(i<n-1) b[i] = sf_cdiv(b[i],a[i]); v[i]=sf_cdiv(sf_csub(v[i],sf_cmul(c[i-1],v[i-1])),a[i]); #endif } for(i=n-2;i>=0;i--) { #ifdef SF_HAS_COMPLEX_H v[i] -= b[i]*v[i+1]; #else v[i] = sf_csub(v[i],sf_cmul(b[i],v[i+1])); #endif } }
static sf_complex cdprod (int j, const sf_complex* a, const sf_complex* b) /* complex dot product */ { int i; sf_complex c; c = sf_cmplx(0.,0.); for (i=1; i <= j; i++) { #ifdef SF_HAS_COMPLEX_H c += a[j-i]*conjf(b[i]); #else c = sf_cadd(c,sf_cmul(a[j-i],conjf(b[i]))); #endif } return c; }
sf_complex plyr_val(int n, float *r, sf_complex x) /*< polynomial value >*/ { int i; sf_complex val; val = sf_cmplx(r[n],0.0); for(i=0; i<n; i++) { #ifdef SF_HAS_COMPLEX_H val = val*(x-r[i]); #else val = sf_cmul(val,sf_cadd(x,sf_cmplx(-r[i],0.0))); #endif } return val; }
sf_complex poly_val(int n, float *c, sf_complex x ) /*< polynomial value >*/ { int i; sf_complex val; val = sf_cmplx(0.0,0.0); for(i=n; i>=0; i--) { #ifdef SF_HAS_COMPLEX_H val = val*x + c[i]; #else val = sf_cadd(sf_cmul(val,x),sf_cmplx(c[i],0.0)); #endif } return val; }
void rweone_ssh( sf_complex *v, float w, float *aa) /*< space-domain phase shift >*/ { int ig; sf_complex ikz; for(ig=0;ig<ag.n;ig++) { ikz = sf_cmplx(0.,w * aa[ig]); #ifdef SF_HAS_COMPLEX_H v[ig] *= cexpf( ikz * at.d ); #else v[ig] = sf_cmul(v[ig],cexpf(sf_crmul( ikz, at.d ))); #endif } }
void rweone_spi( sf_complex *swf, sf_complex *rwf, float *iii) /*< shot-record imaging condition >*/ { int ig; rweone_tap(swf); rweone_tap(rwf); for(ig=0;ig<ag.n;ig++) { #ifdef SF_HAS_COMPLEX_H iii[ig] += crealf( conjf(swf[ig]) * rwf[ig] ); #else iii[ig] += crealf(sf_cmul(conjf(swf[ig]), rwf[ig]) ); #endif } }
/*------------------------------------------------------------*/ void sf_sft3a1(sf_complex ***pp, sft3d sft, sf_fft3d fft) /*< apply shift on axis 1 >*/ { int i1,i2,i3; for (i3=0; i3<fft->n3; i3++) { for (i2=0; i2<fft->n2; i2++) { for(i1=0; i1<fft->n1; i1++) { #ifdef SF_HAS_COMPLEX_H pp[i3][i2][i1] *= sft->www[i1]; #else pp[i3][i2][i1] = sf_cmul(pp[i3][i2][i1],sft->www[i1]); #endif } } } }
void RytovSensitivity(float xx,float xy,float xz,float *output) /*< Rytov sensitivity >*/ { int iw; double omega,scale; sf_complex U,dU; sf_complex val; *output=0.; for (iw=0; iw<nw; iw++) { omega=ow+dw*iw; scale =cos((SF_PI/2)*(((double) iw+1)/((double) nw+1))); scale*=scale; #ifdef SF_HAS_COMPLEX_H /* background field */ U=Green(rx,ry,rz,sx,sy,sz,omega)*scale; /* scattered field */ val=Green(xx,xy,xz,sx,sy,sz,omega)*scale; val*=Green(rx,ry,rz,xx,xy,xz,omega); dU=val*(-omega*omega*dv); val=U*omega; *output -= scale*cimagf(dU/U); #else /* background field */ U=sf_crmul(Green(rx,ry,rz,sx,sy,sz,omega),scale); /* scattered field */ val=sf_crmul(Green(xx,xy,xz,sx,sy,sz,omega),scale); val=sf_cmul(val,Green(rx,ry,rz,xx,xy,xz,omega)); dU=sf_crmul(val,-omega*omega*dv); val=sf_crmul(U,omega); *output -= scale*cimagf(sf_cdiv(dU,U)); #endif } return; }
void rweone_ssf( sf_complex *v, float w, float *aa, float a0) /*< split-step Fourier correction >*/ { int ig; sf_complex ikz; for(ig=0; ig<ag.n; ig++) { ikz = sf_cmplx(0.,w * (aa[ig] - a0)); #ifdef SF_HAS_COMPLEX_H v[ig] *= cexpf( ikz * ( at.d) ); #else v[ig] = sf_cmul(v[ig],cexpf(sf_crmul(ikz, at.d))); #endif } }
/*------------------------------------------------------------*/ void sf_ompsft3a2(sf_complex ***pp, ompsft3d sft, ompfft3d fft, int ompith) /*< apply shift on axis 2 >*/ { int i1,i2,i3; for (i3=0; i3<fft->n3; i3++) { for (i2=0; i2<fft->n2; i2++) { for(i1=0; i1<fft->n1; i1++) { #ifdef SF_HAS_COMPLEX_H pp[i3][i2][i1] *= sft->www[ompith][i2]; #else pp[i3][i2][i1] = sf_cmul(pp[i3][i2][i1],sft->www[ompith][i2]); #endif } } } }
float cblas_scnrm2 (int n, const void* x, int sx) /*< sum |x_i|^2 >*/ { int i, ix; float xn; const sf_complex* c; c = (const sf_complex*) x; xn = 0.0; for (i=0; i < n; i++) { ix = i*sx; #ifdef SF_HAS_COMPLEX_H xn += crealf(c[ix]*conjf(c[ix])); #else xn += crealf(sf_cmul(c[ix],conjf(c[ix]))); #endif } return xn; }
sf_complex plyr_dval(int n, float *r, sf_complex x) /*< polynomial derivative value >*/ { int i, k; sf_complex v1, v2; v1 = sf_cmplx(1.0,0.0); k = 0; for(i=0; i<n; i++) { #ifdef SF_HAS_COMPLEX_H v2 = x - r[i]; #else v2 = sf_cadd(x,sf_cmplx(-r[i],0.0)); #endif if(crealf(v2) == 0.0f && cimagf(v2) == 0.0f) k++; #ifdef SF_HAS_COMPLEX_H else v1 *= v2; #else else v1 = sf_cmul(v1,v2); #endif }
void rweone_phs_old( sf_complex *v, float w, float a0, float b0 ) /*< Fourier-domain phase shift >*/ { int ig,ikg; float kg,arg; float ta,tb,tt; sf_complex ikz; rweone_fft(false,(kiss_fft_cpx*) v); for(ig=0;ig<ag.n;ig++) { ikg = KMAP(ig,ag.n); kg = okg + ikg * dkg; ta = a0*w; tb = b0*kg; tt = tb/ta; arg = 1.0 - tt*tt; if(arg<0.) { ikz = sf_cmplx(SF_ABS(ta) * sqrtf(-arg),0.); } else { ikz = sf_cmplx(0.,ta * sqrtf(+arg)); } #ifdef SF_HAS_COMPLEX_H v[ig] *= cexpf( ikz * (-at.d)); #else v[ig] = sf_cmul(v[ig],cexpf(sf_crmul(ikz,-at.d))); #endif } rweone_fft( true,(kiss_fft_cpx*) v); }
void xkolmog_helix(cfilter cross, cfilter fac1, cfilter fac2) /*< Helix filter factorization >*/ { int ik, ih; float dw, w; dw=2.0*SF_PI/nk; for (ik=0; ik < nk; ik++) { fft1[ik]=sf_cmplx(0.,0.); w = dw*ik; for (ih=0; ih < cross->nh; ih++) { #ifdef SF_HAS_COMPLEX_H fft1[ik] += cross->flt[ih]*cexpf(sf_cmplx(0.,cross->lag[ih]*w)); #else fft1[ik] = sf_cadd(fft1[ik], sf_cmul(cross->flt[ih], cexpf(sf_cmplx(0.,cross->lag[ih]*w)))); #endif } #ifdef SF_HAS_COMPLEX_H fft1[ik] /= nk; #else fft1[ik] = sf_crmul(fft1[ik],1.0/nk); #endif } xkolmog(fft1,fft2); for (ih=0; ih < fac1->nh; ih++) { fac1->flt[ih] = fft1[fac1->lag[ih]]; fac2->flt[ih] = fft2[fac2->lag[ih]]; } }
int propnewc(sf_complex **ini, sf_complex **lt, sf_complex **rt, int nz, int nx, int nt, int m2, int nkzx, char *mode, int pad1, int snap, sf_complex **cc, sf_complex ***wvfld, bool verb, bool correct, sf_complex *alpha, sf_complex *beta) /*^*/ { /* index variables */ int it,iz,ix,im,ik,i,j,wfit; int nz2,nx2,nk,nzx2; sf_complex c; /* wavefield */ sf_complex **wave,**wave2, *curr, *currm, *cwave, *cwavem, *curr1, *curr2; nk = cfft2_init(pad1,nz,nx,&nz2,&nx2); nzx2 = nz2*nx2; if (nk!=nkzx) sf_error("nk discrepancy!"); curr = sf_complexalloc(nzx2); if (correct) { curr1 = sf_complexalloc(nzx2); curr2 = sf_complexalloc(nzx2); } currm = sf_complexalloc(nzx2); cwave = sf_complexalloc(nk); cwavem = sf_complexalloc(nk); wave = sf_complexalloc2(nk,m2); wave2 = sf_complexalloc2(nzx2,m2); icfft2_allocate(cwavem); /* initialization */ for (ix = 0; ix < nx2; ix++) { for (iz=0; iz < nz2; iz++) { j = iz+ix*nz2; if (ix<nx && iz<nz) curr[j] = ini[ix][iz]; else curr[j] = sf_cmplx(0.,0.); } } wfit = 0; /* MAIN LOOP */ for (it=0; it<nt; it++) { if(verb) sf_warning("it=%d;",it); /* outout wavefield */ if(snap>0) { if(it%snap==0 && wfit<=(int)(nt-1)/snap) { for (ix=0; ix<nx; ix++) for (iz=0; iz<nz; iz++) wvfld[wfit][ix][iz] = curr[iz+ix*nz2]; wfit++; } } if (mode[0]=='m') { /* matrix multiplication */ for (im = 0; im < m2; im++) { for (ix = 0; ix < nx; ix++) { for (iz=0; iz < nz; iz++) { i = iz+ix*nz; /* original grid */ j = iz+ix*nz2; /* padded grid */ #ifdef SF_HAS_COMPLEX_H currm[j] = lt[im][i]*curr[j]; #else currm[j] = sf_cmul(lt[im][i], curr[j]); #endif } } cfft2(currm,wave[im]); } for (ik = 0; ik < nk; ik++) { c = sf_cmplx(0.,0.); for (im = 0; im < m2; im++) { #ifdef SF_HAS_COMPLEX_H c += wave[im][ik]*rt[ik][im]; #else c += sf_cmul(wave[im][ik],rt[ik][im]); //complex multiplies complex #endif } cwave[ik] = c; } /* matrix multiplication */ for (im = 0; im < m2; im++) { for (ik = 0; ik < nk; ik++) { #ifdef SF_HAS_COMPLEX_H cwavem[ik] = cwave[ik]*rt[ik][im]; #else cwavem[ik] = sf_cmul(cwave[ik],rt[ik][im]); //complex multiplies complex #endif } icfft2(wave2[im],cwavem); } for (ix = 0; ix < nx; ix++) { for (iz=0; iz < nz; iz++) { i = iz+ix*nz; /* original grid */ j = iz+ix*nz2; /* padded grid */ c = sf_cmplx(0.,0.); for (im = 0; im < m2; im++) { #ifdef SF_HAS_COMPLEX_H c += lt[im][i]*wave2[im][j]; #else c += sf_cmul(lt[im][i], wave2[im][j]); #endif } curr[j] = c; } } if (correct) { for (ix = 0; ix < nx2; ix++) { for (iz=0; iz < nz2; iz++) { i = iz+ix*nz; /* original grid */ j = iz+ix*nz2; /* padded grid */ if (ix<nx && iz<nz) { #ifdef SF_HAS_COMPLEX_H currm[j] = curr[j]/alpha[i]; #else currm[j] = sf_cdiv(curr[j],alpha[i]); #endif } else { currm[j] = sf_cmplx(0.,0.); } } } cfft2(currm,cwave); for (ik = 0; ik < nk; ik++) { #ifdef SF_HAS_COMPLEX_H cwavem[ik] = cwave[ik]/beta[ik]; #else cwavem[ik] = sf_cdiv(cwave[ik],beta[ik]); #endif } icfft2(curr1,cwavem); for (ix = nx; ix < nx2; ix++) { for (iz=nz; iz < nz2; iz++) { j = iz+ix*nz2; /* padded grid */ curr1[j] = sf_cmplx(0.,0.); } } /**/ cfft2(curr,cwave); for (ik = 0; ik < nk; ik++) { #ifdef SF_HAS_COMPLEX_H cwavem[ik] = cwave[ik]/conjf(beta[ik]); #else cwavem[ik] = sf_cdiv(cwave[ik],conjf(beta[ik])); #endif } icfft2(curr,cwavem); for (ix = 0; ix < nx2; ix++) { for (iz=0; iz < nz2; iz++) { i = iz+ix*nz; /* original grid */ j = iz+ix*nz2; /* padded grid */ if (ix<nx && iz<nz) { #ifdef SF_HAS_COMPLEX_H curr2[j] = curr[j]/conjf(alpha[i]); #else curr2[j] = sf_cdiv(curr[j],conjf(alpha[i])); #endif } else { curr2[j] = sf_cmplx(0.,0.); } } } for (ix = 0; ix < nx2; ix++) { for (iz=0; iz < nz2; iz++) { j = iz+ix*nz2; /* padded grid */ #ifdef SF_HAS_COMPLEX_H curr[j] = (curr1[j] + curr2[j])/2.; #else curr[j] = sf_crmul(curr1[j]+curr2[j],0.5); #endif } } } } else if (mode[0]=='x') { cfft2(curr,cwave); /* matrix multiplication */ for (im = 0; im < m2; im++) { for (ik = 0; ik < nk; ik++) { #ifdef SF_HAS_COMPLEX_H cwavem[ik] = cwave[ik]*rt[ik][im]; #else cwavem[ik] = sf_cmul(cwave[ik],rt[ik][im]); //complex multiplies complex #endif } icfft2(wave2[im],cwavem); } for (ix = 0; ix < nx; ix++) { for (iz=0; iz < nz; iz++) { i = iz+ix*nz; /* original grid */ j = iz+ix*nz2; /* padded grid */ c = sf_cmplx(0.,0.); for (im = 0; im < m2; im++) { #ifdef SF_HAS_COMPLEX_H c += lt[im][i]*wave2[im][j]; #else c += sf_cmul(lt[im][i], wave2[im][j]); #endif } curr[j] = c; } } /* matrix multiplication */ for (im = 0; im < m2; im++) { for (ix = 0; ix < nx; ix++) { for (iz=0; iz < nz; iz++) { i = iz+ix*nz; /* original grid */ j = iz+ix*nz2; /* padded grid */ #ifdef SF_HAS_COMPLEX_H currm[j] = lt[im][i]*curr[j]; #else currm[j] = sf_cmul(lt[im][i], curr[j]); #endif } } cfft2(currm,wave[im]); } for (ik = 0; ik < nk; ik++) { c = sf_cmplx(0.,0.); for (im = 0; im < m2; im++) { #ifdef SF_HAS_COMPLEX_H c += wave[im][ik]*rt[ik][im]; #else c += sf_cmul(wave[im][ik],rt[ik][im]); //complex multiplies complex #endif } cwavem[ik] = c; } icfft2(curr,cwavem); if (correct) { for (ix = 0; ix < nx2; ix++) { for (iz=0; iz < nz2; iz++) { i = iz+ix*nz; /* original grid */ j = iz+ix*nz2; /* padded grid */ if (ix<nx && iz<nz) { #ifdef SF_HAS_COMPLEX_H currm[j] = curr[j]/alpha[i]; #else currm[j] = sf_cdiv(curr[j],alpha[i]); #endif } else { currm[j] = sf_cmplx(0.,0.); } } } cfft2(currm,cwave); for (ik = 0; ik < nk; ik++) { #ifdef SF_HAS_COMPLEX_H cwavem[ik] = cwave[ik]/beta[ik]; #else cwavem[ik] = sf_cdiv(cwave[ik],beta[ik]); #endif } icfft2(curr,cwavem); for (ix = nx; ix < nx2; ix++) { for (iz=nz; iz < nz2; iz++) { j = iz+ix*nz2; /* padded grid */ curr[j] = sf_cmplx(0.,0.); } } } } else if (mode[0]=='n') { /* matrix multiplication */ for (im = 0; im < m2; im++) { for (ix = 0; ix < nx; ix++) { for (iz=0; iz < nz; iz++) { i = iz+ix*nz; /* original grid */ j = iz+ix*nz2; /* padded grid */ #ifdef SF_HAS_COMPLEX_H currm[j] = lt[im][i]*curr[j]; #else currm[j] = sf_cmul(lt[im][i], curr[j]); #endif } } cfft2(currm,wave[im]); } for (ik = 0; ik < nk; ik++) { c = sf_cmplx(0.,0.); for (im = 0; im < m2; im++) { #ifdef SF_HAS_COMPLEX_H c += wave[im][ik]*rt[ik][im]; #else c += sf_cmul(wave[im][ik],rt[ik][im]); //complex multiplies complex #endif } cwavem[ik] = c; } icfft2(curr,cwavem); /* matrix multiplication */ for (im = 0; im < m2; im++) { for (ix = 0; ix < nx; ix++) { for (iz=0; iz < nz; iz++) { i = iz+ix*nz; /* original grid */ j = iz+ix*nz2; /* padded grid */ #ifdef SF_HAS_COMPLEX_H currm[j] = lt[im][i]*curr[j]; #else currm[j] = sf_cmul(lt[im][i], curr[j]); #endif } } cfft2(currm,wave[im]); } for (ik = 0; ik < nk; ik++) { c = sf_cmplx(0.,0.); for (im = 0; im < m2; im++) { #ifdef SF_HAS_COMPLEX_H c += wave[im][ik]*rt[ik][im]; #else c += sf_cmul(wave[im][ik],rt[ik][im]); //complex multiplies complex #endif } cwavem[ik] = c; } icfft2(curr,cwavem); } else if (mode[0]=='p') { cfft2(curr,cwave); /* matrix multiplication */ for (im = 0; im < m2; im++) { for (ik = 0; ik < nk; ik++) { #ifdef SF_HAS_COMPLEX_H cwavem[ik] = cwave[ik]*rt[ik][im]; #else cwavem[ik] = sf_cmul(cwave[ik],rt[ik][im]); //complex multiplies complex #endif } icfft2(wave2[im],cwavem); } for (ix = 0; ix < nx; ix++) { for (iz=0; iz < nz; iz++) { i = iz+ix*nz; /* original grid */ j = iz+ix*nz2; /* padded grid */ c = sf_cmplx(0.,0.); for (im = 0; im < m2; im++) { #ifdef SF_HAS_COMPLEX_H c += lt[im][i]*wave2[im][j]; #else c += sf_cmul(lt[im][i], wave2[im][j]); #endif } curr[j] = c; } } cfft2(curr,cwave); /* matrix multiplication */ for (im = 0; im < m2; im++) { for (ik = 0; ik < nk; ik++) { #ifdef SF_HAS_COMPLEX_H cwavem[ik] = cwave[ik]*rt[ik][im]; #else cwavem[ik] = sf_cmul(cwave[ik],rt[ik][im]); //complex multiplies complex #endif } icfft2(wave2[im],cwavem); } for (ix = 0; ix < nx; ix++) { for (iz=0; iz < nz; iz++) { i = iz+ix*nz; /* original grid */ j = iz+ix*nz2; /* padded grid */ c = sf_cmplx(0.,0.); for (im = 0; im < m2; im++) { #ifdef SF_HAS_COMPLEX_H c += lt[im][i]*wave2[im][j]; #else c += sf_cmul(lt[im][i], wave2[im][j]); #endif } curr[j] = c; } } } else sf_error("Check mode parameter!"); } /* time stepping */ if(verb) sf_warning("."); /* output final result*/ for (ix=0; ix<nx; ix++) for (iz=0; iz<nz; iz++) cc[ix][iz] = curr[iz+ix*nz2]; cfft2_finalize(); return 0; }
int main(int argc, char* argv[]) { clock_t tstart, tend; double duration; /*Flag*/ bool verb, cmplx; /*I/O*/ sf_file Fsrc,Fo,Frec; /* I/O files */ sf_file left, right; /*left/right matrix*/ sf_file Fvel, Fden, Ffft; /*Model*/ sf_axis at,az,ax; /* cube axes */ /* I/O arrays*/ float *src; /*point source, distributed source*/ float **lt, **rt; float **vel, **den, **c11; /* Grid index variables */ int it,iz,im,ik,ix,i,j; int nt,nz,nx, m2, nk, nkx, nkz, nzx, nz2, nx2, nzx2, n1, n2, pad1; float cx, cz; float kx, kz, dkx, dkz, kx0, kz0; float dx, dz, dt, d1, d2; float ox, oz; sf_complex *cwavex, *cwavez, *cwavemx, *cwavemz; float **record; float **wavex, **wavez; float *curtxx, *pretxx; float *curvx, *prevx, *curvz, *prevz; /*source*/ spara sp={0}; int srcrange; float srctrunc; bool srcdecay; float slx, slz; int spx, spz; /*options*/ float gdep; int gp; tstart = clock(); sf_init(argc,argv); if(!sf_getbool("verb",&verb)) verb=false; /* verbosity */ Fvel = sf_input("vel"); Fden = sf_input("den"); /* setup I/O files */ Fsrc = sf_input ("in" ); Fo = sf_output("out"); Frec = sf_output("rec"); /*record*/ /* Read/Write axes */ at = sf_iaxa(Fsrc,1); nt = sf_n(at); dt = sf_d(at); ax = sf_iaxa(Fvel,2); nx = sf_n(ax); dx = sf_d(ax); ox=sf_o(ax); az = sf_iaxa(Fvel,1); nz = sf_n(az); dz = sf_d(az); oz=sf_o(az); sf_oaxa(Fo,az,1); sf_oaxa(Fo,ax,2); sf_oaxa(Fo,at,3); /*set for record*/ sf_oaxa(Frec, at, 1); sf_oaxa(Frec, ax, 2); if (!sf_getbool("cmplx",&cmplx)) cmplx=false; /* use complex FFT */ if (!sf_getint("pad1",&pad1)) pad1=1; /* padding factor on the first axis */ nk = fft2_init(cmplx,pad1,nz,nx,&nz2,&nx2); nzx = nz*nx; nzx2 = nz2*nx2; /* propagator matrices */ left = sf_input("left"); right = sf_input("right"); if (!sf_histint(left,"n1",&n2) || n2 != nzx) sf_error("Need n1=%d in left",nzx); if (!sf_histint(left,"n2",&m2)) sf_error("Need n2= in left"); if (!sf_histint(right,"n1",&n2) || n2 != m2) sf_error("Need n1=%d in right",m2); if (!sf_histint(right,"n2",&n2) || n2 != nk) sf_error("Need n2=%d in right",nk); lt = sf_floatalloc2(nzx,m2); rt = sf_floatalloc2(m2,nk); sf_floatread(lt[0],nzx*m2,left); sf_floatread(rt[0],m2*nk,right); /*model veloctiy & density*/ if (!sf_histint(Fvel,"n1", &n1) || n1 != nz) sf_error("Need n1=%d in vel", nz); if (!sf_histfloat(Fvel,"d1", &d1) || d1 != dz) sf_error("Need d1=%d in vel", dz); if (!sf_histint(Fvel,"n2", &n2) || n2 != nx) sf_error("Need n2=%d in vel", nx); if (!sf_histfloat(Fvel,"d2", &d2) || d2 != dx) sf_error("Need d2=%d in vel", dx); if (!sf_histint(Fden,"n1", &n1) || n1 != nz) sf_error("Need n1=%d in den", nz); if (!sf_histfloat(Fden,"d1", &d1) || d1 != dz) sf_error("Need d1=%d in den", dz); if (!sf_histint(Fden,"n2", &n2) || n2 != nx) sf_error("Need n2=%d in den", nx); if (!sf_histfloat(Fden,"d2", &d2) || d2 != dx) sf_error("Need d2=%d in den", dx); vel = sf_floatalloc2(nz, nx); den = sf_floatalloc2(nz, nx); c11 = sf_floatalloc2(nz, nx); sf_floatread(vel[0], nzx, Fvel); sf_floatread(den[0], nzx, Fden); for (ix = 0; ix < nx; ix++) { for (iz = 0; iz < nz; iz++) { c11[ix][iz] = den[ix][iz]*vel[ix][iz]*vel[ix][iz]; } } /*parameters of fft*/ Ffft = sf_input("fft"); if (!sf_histint(Ffft,"n1", &nkz)) sf_error("Need n1 in fft"); if (!sf_histint(Ffft,"n2", &nkx)) sf_error("Need n2 in fft"); if ( nkx*nkz != nk ) sf_error("Need nk=nkx*nkz, nk=%d, nkx=%d, nkz=%d", nk, nkx, nkz); if (!sf_histfloat(Ffft,"d1", &dkz)) sf_error("Need d1 in fft"); if (!sf_histfloat(Ffft,"d2", &dkx)) sf_error("Need d2 in fft"); if (!sf_histfloat(Ffft,"o1", &kz0)) sf_error("Need o1 in fft"); if (!sf_histfloat(Ffft,"o2", &kx0)) sf_error("Need o2 in fft"); /*parameters of geometry*/ if (!sf_getfloat("gdep", &gdep)) gdep = 0.0; /*depth of geophone (meter)*/ if (gdep <0.0) sf_error("gdep need to be >=0.0"); /*source and receiver location*/ if (!sf_getfloat("slx", &slx)) slx=-1.0; /*source location x */ if (!sf_getint("spx", &spx)) spx = -1; /*source location x (index)*/ if((slx<0 && spx <0) || (slx>=0 && spx >=0 )) sf_error("Need src location"); if (slx >= 0 ) spx = (int)((slx-ox)/dx+0.5); if (!sf_getfloat("slz", &slz)) slz = -1.0; /* source location z */ if (!sf_getint("spz", &spz)) spz=-1; /*source location z (index)*/ if((slz<0 && spz <0) || (slz>=0 && spz >=0 )) sf_error("Need src location"); if (slz >= 0 ) spz = (int)((slz-ox)/dz+0.5); if (!sf_getfloat("gdep", &gdep)) gdep = -1.0; /* recorder depth on grid*/ if (!sf_getint("gp", &gp)) gp=0; /* recorder depth on index*/ if ( gdep>=oz) { gp = (int)((gdep-oz)/dz+0.5);} if (gp < 0.0) sf_error("gdep need to be >=oz"); /*source and receiver location*/ if (!sf_getbool("srcdecay", &srcdecay)) srcdecay=false; /*source decay*/ if (!sf_getint("srcrange", &srcrange)) srcrange=10; /*source decay range*/ if (!sf_getfloat("srctrunc", &srctrunc)) srctrunc=100; /*trunc source after srctrunc time (s)*/ /* read wavelet & reflectivity */ src = sf_floatalloc(nt); sf_floatread(src,nt,Fsrc); curtxx = sf_floatalloc(nzx2); curvx = sf_floatalloc(nzx2); curvz = sf_floatalloc(nzx2); pretxx = sf_floatalloc(nzx); prevx = sf_floatalloc(nzx); prevz = sf_floatalloc(nzx); cwavex = sf_complexalloc(nk); cwavez = sf_complexalloc(nk); cwavemx = sf_complexalloc(nk); cwavemz = sf_complexalloc(nk); wavex = sf_floatalloc2(nzx2,m2); wavez = sf_floatalloc2(nzx2,m2); record = sf_floatalloc2(nt,nx); ifft2_allocate(cwavemx); ifft2_allocate(cwavemz); for (iz=0; iz < nzx; iz++) { pretxx[iz]=0.; prevx[iz] =0.; prevz[iz] =0.; } for (iz=0; iz < nzx2; iz++) { curtxx[iz]=0.; curvx[iz]=0.; curvz[iz]=0.; } /* Check parameters*/ if(verb) { sf_warning("======================================"); #ifdef SF_HAS_FFTW sf_warning("FFTW is defined"); #endif #ifdef SF_HAS_COMPLEX_H sf_warning("Complex is defined"); #endif sf_warning("nx=%d nz=%d nzx=%d dx=%f dz=%f", nx, nz, nzx, dx, dz); sf_warning("nkx=%d nkz=%d dkx=%f dkz=%f nk=%d", nkx, nkz, dkx, dkz, nk); sf_warning("nx2=%d nz2=%d nzx2=%d", nx2, nz2, nzx2); sf_warning("======================================"); } /*set source*/ sp.trunc=srctrunc; sp.srange=srcrange; sp.alpha=0.5; sp.decay=srcdecay?1:0; /* MAIN LOOP */ for (it=0; it<nt; it++) { if(verb) sf_warning("it=%d/%d;",it,nt-1); /*vx, vz--- matrix multiplication */ fft2(curtxx,cwavex); /* P^(k,t) */ for (im = 0; im < m2; im++) { for (ik = 0; ik < nk; ik++) { kx = kx0+dkx*(ik/nkz); kz = kz0+dkz*(ik%nkz); #ifdef SF_HAS_COMPLEX_H cwavemz[ik] = cwavex[ik]*rt[ik][im]; cwavemx[ik] = fplus(kx,dx)*cwavemz[ik]; cwavemz[ik] = fplus(kz,dz)*cwavemz[ik]; #else cwavemz[ik] = sf_crmul(cwavex[ik],rt[ik][im]); cwavemx[ik] = sf_cmul(fplus(kx,dx), cwavemz[ik]); cwavemz[ik] = sf_cmul(fplus(kz,dz), cwavemz[ik]); #endif } ifft2(wavex[im], cwavemx); /* dp/dx */ ifft2(wavez[im], cwavemz); /* dp/dz */ } for (ix = 0; ix < nx; ix++) { for (iz = 0; iz < nz; iz++) { i = iz+ix*nz; /* original grid */ j = iz+ix*nz2; /* padded grid */ cx = 0.0; cz = 0.0; for (im=0; im<m2; im++) { cx += lt[im][i]*wavex[im][j]; cz += lt[im][i]*wavez[im][j]; } curvx[j] = -1*dt/den[ix][iz]*cx + prevx[i]; /*vx(t+dt/2) = -dt/rho*dp/dx(t) + vx(t-dt/2) */ prevx[i] = curvx[j]; curvz[j] = -1*dt/den[ix][iz]*cz + prevz[i]; prevz[i] = curvz[j]; } } /*txx--- matrix multiplication */ fft2(curvx, cwavex); fft2(curvz, cwavez); for (im = 0; im < m2; im++) { for (ik = 0; ik < nk; ik++ ) { kx = kx0 + dkx*(ik/nkz); kz = kz0 + dkz*(ik%nkz); #ifdef SF_HAS_COMPLEX_H cwavemz[ik] = cwavez[ik]*rt[ik][im]; cwavemx[ik] = cwavex[ik]*rt[ik][im]; cwavemx[ik] = fminu(kx,dx)*cwavemx[ik]; cwavemz[ik] = fminu(kz,dz)*cwavemz[ik]; #else cwavemz[ik] = sf_crmul(cwavez[ik],rt[ik][im]); cwavemx[ik] = sf_crmul(cwavex[ik],rt[ik][im]); cwavemx[ik] = sf_cmul(fplus(kx,dx), cwavemx[ik]); cwavemz[ik] = sf_cmul(fplus(kz,dz), cwavemz[ik]); #endif } ifft2(wavex[im], cwavemx); /* dux/dx */ ifft2(wavez[im], cwavemz); /* duz/dz */ } for (ix = 0; ix < nx; ix++) { for (iz = 0; iz < nz; iz++) { i = iz+ix*nz; /* original grid */ j = iz+ix*nz2; /* padded grid */ cx = 0.0; cz = 0.0; for (im=0; im<m2; im++) { cx += lt[im][i]*wavex[im][j]; cz += lt[im][i]*wavez[im][j]; } curtxx[j] = -1*dt*c11[ix][iz]*(cx+cz) + pretxx[i]; } } if ((it*dt)<=sp.trunc ) { curtxx[spz+spx*nz2] += src[it]*dt; } for (ix = 0; ix < nx; ix++) { /* write wavefield to output */ sf_floatwrite(pretxx+ix*nz,nz,Fo); } /*record*/ for (ix = 0; ix < nx; ix++){ record[ix][it] = pretxx[ix*nz+gp]; } for (ix = 0; ix < nx; ix++) { for (iz = 0; iz < nz; iz++) { i = iz+ix*nz; /* original grid */ j = iz+ix*nz2; /* padded grid */ pretxx[i] = curtxx[j]; } } }/*End of MAIN LOOP*/ if(verb) sf_warning("."); for ( ix = 0; ix < nx; ix++) { sf_floatwrite(record[ix], nt, Frec); } tend = clock(); duration=(double)(tend-tstart)/CLOCKS_PER_SEC; sf_warning(">> The CPU time of sfsglr is: %f seconds << ", duration); exit (0); }
int main(int argc, char* argv[]) { fint1 str, istr; int i1,i2, n1,n2,n3, ix,iv,ih, ib, ie, nb, nx,nv,nh, nw, next; float d1,o1,d2,o2, eps, w,x,k, v0,v2,v,v1,dv, dx, h0,dh,h, num, den, t, dw; float *trace=NULL, *strace=NULL, ***stack=NULL, ***stack2=NULL, ***cont=NULL, **image=NULL; sf_complex *ctrace=NULL, *ctrace0=NULL, shift; char *time=NULL, *space=NULL, *unit=NULL; size_t len; static kiss_fftr_cfg forw, invs; sf_file in=NULL, out=NULL; bool sembl; sf_init (argc,argv); in = sf_input("in"); out = sf_output("out"); if (!sf_histint(in,"n1",&n1)) sf_error("No n1= in input"); if (!sf_histint(in,"n2",&nx)) sf_error("No n2= in input"); if (!sf_histint(in,"n3",&nh)) sf_error("No n3= in input"); if (!sf_getint("nb",&nb)) nb=2; if (!sf_getfloat("eps",&eps)) eps=0.01; if (!sf_getint("pad",&n2)) n2=n1; if (!sf_getint("pad2",&n3)) n3=2*kiss_fft_next_fast_size((n2+1)/2); nw = n3/2+1; forw = kiss_fftr_alloc(n3,0,NULL,NULL); invs = kiss_fftr_alloc(n3,1,NULL,NULL); if (NULL == forw || NULL == invs) sf_error("KISS FFT allocation error"); if (!sf_histfloat(in,"o1",&o1)) o1=0.; o2 = o1*o1; if(!sf_histfloat(in,"d1",&d1)) sf_error("No d1= in input"); d2 = o1+(n1-1)*d1; d2 = (d2*d2 - o2)/(n2-1); if (!sf_getint("nv",&nv)) sf_error("Need nv="); if (!sf_getfloat("dv",&dv)) sf_error("Need dv="); if (!sf_getfloat("v0",&v0) && !sf_histfloat(in,"v0",&v0)) sf_error("Need v0="); if (!sf_getbool("semblance",&sembl)) sembl=true; /* if y, compute semblance; if n, stack */ if(!sf_histfloat(in,"o3",&h0)) sf_error("No o2= in input"); if(!sf_histfloat(in,"d3",&dh)) sf_error("No d2= in input"); if(!sf_histfloat(in,"d2",&dx)) sf_error("No d3= in input"); sf_putfloat(out,"o3",v0+dv); sf_putfloat(out,"d3",dv); sf_putint(out,"n3",nv); sf_putstring(out,"label3","Velocity"); if (NULL != (time = sf_histstring(in,"label1")) && NULL != (space = sf_histstring(in,"label2"))) { len = strlen(time)+strlen(space)+2; unit = sf_charalloc(len); snprintf(unit,len,"%s/%s",space,time); sf_putstring(out,"unit3",unit); free(time); free(space); } dx = 2*SF_PI/(2*kiss_fft_next_fast_size(nx-1)*dx); dw = 16*SF_PI/(d2*n3); /* 2pi * 8 */ stack = sf_floatalloc3(n1,nx,nv); stack2 = sf_floatalloc3(n1,nx,nv); cont = sf_floatalloc3(n1,nx,nv); image = sf_floatalloc2(n1,nx); trace = sf_floatalloc(n1); strace = sf_floatalloc(n3); ctrace = sf_complexalloc(nw); ctrace0 = sf_complexalloc(nw); if (!sf_getint("extend",&next)) next=4; /* trace extension */ str = fint1_init(next,n1,0); istr = fint1_init(next,n2,0); for (i1=0; i1 < n1*nx*nv; i1++) { stack[0][0][i1] = 0.; stack2[0][0][i1] = 0.; } sf_cosft_init(nx); for (ih=0; ih < nh; ih++) { sf_warning("offset %d of %d;",ih+1,nh); h = h0 + ih*dh; h *= h * 0.5; sf_floatread(image[0],n1*nx,in); for (i1=0; i1 < n1; i1++) { sf_cosft_frw(image[0],i1,n1); } for (ix=0; ix < nx; ix++) { x = ix*dx; x *= x; k = x * 0.5; fint1_set(str,image[ix]); for (i2=0; i2 < n2; i2++) { t = o2+i2*d2; t = sqrtf(t); t = (t-o1)/d1; i1 = t; if (i1 >= 0 && i1 < n1) { strace[i2] = fint1_apply(str,i1,t-i1,false); } else { strace[i2] = 0.; } } for (i2=n2; i2 < n3; i2++) { strace[i2] = 0.; } kiss_fftr(forw,strace, (kiss_fft_cpx *) ctrace0); for (iv=0; iv < nv; iv++) { v = v0 + (iv+1)* dv; v1 = h * (1./(v*v) - 1./(v0*v0)); v2 = k * ((v0*v0) - (v*v)); ctrace[0]=sf_cmplx(0.,0.); /* dc */ for (i2=1; i2 < nw; i2++) { w = i2*dw; w = v2/w+(v1-0.125*o2)*w; shift = sf_cmplx(cosf(w),sinf(w)); #ifdef SF_HAS_COMPLEX_H ctrace[i2] = ctrace0[i2] * shift; #else ctrace[i2] = sf_cmul(ctrace0[i2],shift); #endif } /* w */ kiss_fftri(invs,(const kiss_fft_cpx *) ctrace, strace); fint1_set(istr,strace); for (i1=0; i1 < n1; i1++) { t = o1+i1*d1; t = t*t; t = (t-o2)/d2; i2 = t; if (i2 >= 0 && i2 < n2) { cont[iv][ix][i1] = fint1_apply(istr,i2,t-i2,false); } else { cont[iv][ix][i1] = 0.; } } } /* v */ } /* x */ for (iv=0; iv < nv; iv++) { for (i1=0; i1 < n1; i1++) { sf_cosft_inv(cont[0][0],i1+iv*nx*n1,n1); } } for (iv=0; iv < nv; iv++) { for (ix=0; ix < nx; ix++) { for (i1=0; i1 < n1; i1++) { t = cont[iv][ix][i1]; stack[iv][ix][i1] += t; stack2[iv][ix][i1] += t*t; } /* i1 */ } /* x */ } /* v */ } /* h */ sf_warning("."); for (iv=0; iv < nv; iv++) { for (ix=0; ix < nx; ix++) { for (i1=0; i1 < n1; i1++) { ib = i1-nb > 0? i1-nb: 0; ie = i1+nb+1 < n1? i1+nb+1: n1; num = 0.; den = 0.; if (sembl) { for (i2=ib; i2 < ie; i2++) { t = stack[iv][ix][i2]; num += t*t; den += stack2[iv][ix][i2]; } den *= nh; trace[i1] = den > 0.? num/den: 0.; } else { for (i2=ib; i2 < ie; i2++) { t = stack[iv][ix][i2]; num += t; } den = nh; trace[i1] = num/(den+ FLT_EPSILON); } } sf_floatwrite (trace,n1,out); } } exit(0); }
void cbanded_const_define (float diag /* diagonal */, const sf_complex *offd /* lower off-diagonal */) /*< set matrix coefficients (constant along diagonals) >*/ { int k, m, j; sf_complex ct; float rt; d[0] = diag; for (k = 0; k < band-1; k++) { for (m = k; m >= 0; m--) { ct = offd[m]; for (j = m+1; j < k-1; j++) { #ifdef SF_HAS_COMPLEX_H ct -= o[j][k-j]*conjf(o[j-m-1][k-j])*d[k-j]; #else ct = sf_csub(ct,sf_cmul(o[j][k-j], sf_crmul(conjf(o[j-m-1][k-j]), d[k-j]))); #endif } #ifdef SF_HAS_COMPLEX_H o[m][k-m] = ct/d[k-m]; #else o[m][k-m] = sf_crmul(ct,1./d[k-m]); #endif } rt = diag; for (m = 0; m <= k; m++) { #ifdef SF_HAS_COMPLEX_H rt -= crealf(o[m][k-m]*conjf(o[m][k-m]))*d[k-m]; #else rt -= crealf(sf_cmul(o[m][k-m],conjf(o[m][k-m])))*d[k-m]; #endif } d[k+1] = rt; } for (k = band-1; k < n-1; k++) { for (m = band-1; m >= 0; m--) { ct = offd[m]; for (j = m+1; j < band; j++) { #ifdef SF_HAS_COMPLEX_H ct -= o[j][k-j]*conjf(o[j-m-1][k-j])*d[k-j]; #else ct = sf_csub(ct,sf_cmul(o[j][k-j], sf_crmul(conjf(o[j-m-1][k-j]), d[k-j]))); #endif } #ifdef SF_HAS_COMPLEX_H o[m][k-m] = ct/d[k-m]; #else o[m][k-m] = sf_crmul(ct,1./d[k-m]); #endif } rt = diag; for (m = 0; m < band; m++) { #ifdef SF_HAS_COMPLEX_H rt -= crealf(o[m][k-m]*conjf(o[m][k-m]))*d[k-m]; #else rt -= crealf(sf_cmul(o[m][k-m],conjf(o[m][k-m])))*d[k-m]; #endif } d[k+1] = rt; } }