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 xkolmog(sf_complex *trace1, sf_complex *trace2) /*< convert Fourier-domain cross-correlation to minimum-phase >*/ { int i1; const double eps=1.e-32; for (i1=0; i1 < nk; i1++) { #ifdef SF_HAS_COMPLEX_H fft1[i1] = clogf(trace1[i1]+eps)/nk; #else fft1[i1] = sf_crmul(clogf(sf_cadd(trace1[i1],sf_cmplx(eps,0.))), 1.0/nk); #endif } /* Inverse transform */ kiss_fft(invs,(const kiss_fft_cpx *) fft1, (kiss_fft_cpx *) trace1); #ifdef SF_HAS_COMPLEX_H trace1[0] *= 0.5; trace2[0] = trace1[0]; trace1[nk/2] *= 0.5; trace2[nk/2] = trace1[nk/2]; #else trace1[0] = sf_crmul(trace1[0], 0.5); trace2[0] = trace1[0]; trace1[nk/2] = sf_crmul(trace1[nk/2],0.5); trace2[nk/2] = trace1[nk/2]; #endif for (i1=1+nk/2; i1 < nk; i1++) { trace2[nk-i1] = trace1[i1]; trace1[i1] = sf_cmplx(0.,0.); trace2[i1] = sf_cmplx(0.,0.); } /* Fourier transform */ kiss_fft(forw,(const kiss_fft_cpx *) trace1, (kiss_fft_cpx *) fft1); kiss_fft(forw,(const kiss_fft_cpx *) trace2, (kiss_fft_cpx *) fft2); for (i1=0; i1 < nk; i1++) { #ifdef SF_HAS_COMPLEX_H fft1[i1] = cexpf(fft1[i1])/nk; fft2[i1] = cexpf(fft2[i1])/nk; #else fft1[i1] = sf_crmul(cexpf(fft1[i1]),1./nk); fft2[i1] = sf_crmul(cexpf(fft2[i1]),1./nk); #endif } /* Inverse transform */ kiss_fft(invs,(const kiss_fft_cpx *) fft1, (kiss_fft_cpx *) trace1); kiss_fft(invs,(const kiss_fft_cpx *) fft2, (kiss_fft_cpx *) trace2); for (i1=0; i1 < nk; i1++) { trace2[i1] = conjf(trace2[i1]); } }
/*------------------------------------------------------------*/ void sf_ompfft3a1(bool inv /* inverse/forward flag */, kiss_fft_cpx ***pp /* [n1][n2][n3] */, ompfft3d fft, int ompith) /*< apply FFT on axis 1 >*/ { int i1, i2, i3; if (inv) { /* IFT 1 */ for (i3=0; i3 < fft->n3; i3++) { for(i2=0; i2 < fft->n2; i2++) { #ifdef _OPENMP #pragma omp critical #endif kiss_fft(fft->invs[ompith],pp[i3][i2],pp[i3][i2]); } } /* scaling */ for (i3=0; i3 < fft->n3; i3++) { for (i2=0; i2 < fft->n2; i2++) { for(i1=0; i1 < fft->n1; i1++) { pp[i3][i2][i1] = sf_crmul(pp[i3][i2][i1],fft->scale); } } } } else { /* scaling */ for (i3=0; i3 < fft->n3; i3++) { for (i2=0; i2 < fft->n2; i2++) { for(i1=0; i1 < fft->n1; i1++) { pp[i3][i2][i1] = sf_crmul(pp[i3][i2][i1],fft->scale); } } } /* FFT 1 */ for (i3=0; i3 < fft->n3; i3++) { for(i2=0; i2 < fft->n2; i2++) { #ifdef _OPENMP #pragma omp critical #endif kiss_fft(fft->forw[ompith],pp[i3][i2],pp[i3][i2]); } } } }
/*------------------------------------------------------------*/ void fft2(bool inv /* inverse/forward flag */, kiss_fft_cpx **pp /* [1...n2][1...n1] */) /*< Apply 2-D FFT >*/ { int i1,i2; if (inv) { /* IFT 1 */ for(i2=0; i2 < n2; i2++) { kiss_fft(invs1,pp[i2],pp[i2]); } /* IFT 2 */ for(i1=0; i1 < n1; i1++) { kiss_fft_stride(invs2,pp[0]+i1,trace2,n1); for(i2=0; i2<n2; i2++) { pp[i2][i1] = trace2[i2]; } } /* scaling */ for (i2=0; i2<n2; i2++) { for(i1=0; i1 < n1; i1++) { pp[i2][i1] = sf_crmul(pp[i2][i1],fftscale); } } } else { /* scaling */ for (i2=0; i2<n2; i2++) { for(i1=0; i1 < n1; i1++) { pp[i2][i1] = sf_crmul(pp[i2][i1],fftscale); } } /* FFT 2 */ for(i1=0; i1 < n1; i1++) { kiss_fft_stride(forw2,pp[0]+i1,trace2,n1); for(i2=0; i2<n2; i2++) { pp[i2][i1] = trace2[i2]; } } /* FFT 1 */ for(i2=0; i2 < n2; i2++) { kiss_fft(forw1,pp[i2],pp[i2]); } } }
/*------------------------------------------------------------*/ void sf_fft3a3(bool inv /* inverse/forward flag */, kiss_fft_cpx ***pp /* [n1][n2][n3] */, sf_fft3d fft) /*< apply FFT on axis 3 >*/ { int i1, i2, i3; if (inv) { /* IFT 3 */ for (i2=0; i2 < fft->n2; i2++) { for(i1=0; i1 < fft->n1; i1++) { kiss_fft_stride(fft->invs,pp[0][0]+i1+i2*fft->n1,fft->trace,fft->n1*fft->n2); for(i3=0; i3 < fft->n3; i3++) { pp[i3][i2][i1] = fft->trace[i3]; } } } /* scaling */ for (i3=0; i3 < fft->n3; i3++) { for (i2=0; i2 < fft->n2; i2++) { for(i1=0; i1 < fft->n1; i1++) { pp[i3][i2][i1] = sf_crmul(pp[i3][i2][i1],fft->scale); } } } } else { /* scaling */ for (i3=0; i3 < fft->n3; i3++) { for (i2=0; i2 < fft->n2; i2++) { for(i1=0; i1 < fft->n1; i1++) { pp[i3][i2][i1] = sf_crmul(pp[i3][i2][i1],fft->scale); } } } /* FFT 3 */ for (i2=0; i2 < fft->n2; i2++) { for(i1=0; i1 < fft->n1; i1++) { kiss_fft_stride(fft->forw,pp[0][0]+i1+i2*fft->n1,fft->trace,fft->n1*fft->n2); for(i3=0; i3 < fft->n3; i3++) { pp[i3][i2][i1] = fft->trace[i3]; } } } } }
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 icfft2(sf_complex *out /* [n1*n2] */, sf_complex *inp /* [nk*n2] */) /*< 2-D inverse FFT >*/ { int i1, i2; #pragma omp parallel for private(i2,i1) default(shared) for (i2=0; i2<local_n0; i2++) { for (i1=0; i1<nk; i1++) { dd[i2*nk+i1]=inp[i2*nk+i1]; } } fftwf_execute(icfg); /* FFT centering and normalization*/ #pragma omp parallel for private(i2,i1) default(shared) for (i2=0; i2<local_n0; i2++) { for (i1=0; i1<n1; i1++) { #ifdef SF_HAS_COMPLEX_H out[i2*n1+i1] = ((((i2+local_0_start)%2==0)==(i1%2==0))? wt:-wt) * cc[i2*n1+i1]; #else out[i2*n1+i1] = sf_crmul(cc[i2*n1+i1],((((i2+local_0_start)%2==0)==(i1%2==0))? wt:-wt)); #endif } } }
void explsourcet(sf_complex *curr/*@out@*/, sf_complex *vwavlet, int vit, int vsx, int vsz, int nx2, int nz2, srcpar vps/*decay parameters*/) /*< explosive source >*/ { float phi = 0.0; int cent = (int)vps->range/2; int ix, iz; if (vps->decay ==1){ #ifdef _OPENMP #pragma omp parallel for private(ix,iz,phi) #endif for (ix=0; ix<2*cent; ix++) for (iz=0; iz<2*cent; iz++) { phi = exp( -1*vps->alpha*vps->alpha*((ix-cent)*(ix-cent)+(iz-cent)*(iz-cent)) ); #ifdef SF_HAS_COMPLEX_H curr[(vsx-cent+ix)*nz2+(vsz-cent+iz)] += vwavlet[vit]*phi; #else curr[(vsx-cent+ix)*nz2+(vsz-cent+iz)] += sf_crmul(vwavlet[vit],phi); #endif } } else { curr[vsx*nz2+vsz] += vwavlet[vit]; } }
void icfft2(sf_complex *out /* [n1*n2] */, sf_complex *inp /* [nkk*n2] */) /*< 2-D inverse FFT >*/ { int i1, i2; #ifdef SF_HAS_FFTW fftwf_execute(icfg); #else for (i1=0; i1 < nkk; i1++) { kiss_fft_stride(icfg2,(kiss_fft_cpx *) (inp+i1),ctrace2,nkk); for (i2=0; i2<n2; i2++) { temp[i2][i1] = ctrace2[i2]; } } for (i2=0; i2 < n2; i2++) { kiss_fft_stride(icfg1,temp[i2],(kiss_fft_cpx *) cc[i2],1); } #endif /* FFT centering and normalization*/ for (i2=0; i2<n2; i2++) { for (i1=0; i1<n1; i1++) { #ifdef SF_HAS_COMPLEX_H out[i2*n1+i1] = (((i2%2==0)==(i1%2==0))? wt:-wt) * cc[i2][i1]; #else out[i2*n1+i1] = sf_crmul(cc[i2][i1],(((i2%2==0)==(i1%2==0))? wt:-wt)); #endif } } }
static void triple (int o, int d, int nx, int nb, const sf_complex* x, sf_complex* tmp, bool box) { int i; float wt; sf_complex xi; for (i=0; i < nx + 2*nb; i++) { tmp[i] = sf_cmplx(0.,0.); } if (box) { wt = 1.0/(2*nb-1); for (i=0; i < nx; i++) { #ifdef SF_HAS_COMPLEX_H xi = wt*x[o+i*d]; tmp[i+1] += xi; tmp[i+2*nb] -= xi; #else xi = sf_crmul(x[o+i*d],wt); tmp[i+1] = sf_cadd(tmp[i+1],xi); tmp[i+2*nb] = sf_cadd(tmp[i+2*nb],sf_cneg(xi)); #endif } } else { wt = 1.0/(nb*nb); for (i=0; i < nx; i++) { #ifdef SF_HAS_COMPLEX_H xi = wt*x[o+i*d]; tmp[i] -= xi; tmp[i+nb] += 2*xi; tmp[i+2*nb] -= xi; #else xi = sf_crmul(x[o+i*d],wt); tmp[i] = sf_cadd(tmp[i],sf_cneg(xi)); tmp[i+nb] = sf_cadd(tmp[i+nb],sf_crmul(xi,2.)); tmp[i+2*nb] = sf_cadd(tmp[i+2*nb],sf_cneg(xi)); #endif } } }
/*------------------------------------------------------------*/ void sf_fft3a1(bool inv /* inverse/forward flag */, kiss_fft_cpx ***pp /* [n1][n2][n3] */, sf_fft3d fft) /*< apply FFT on axis 1 >*/ { int i1, i2, i3; if (inv) { /* IFT 1 */ for (i3=0; i3 < fft->n3; i3++) { for(i2=0; i2 < fft->n2; i2++) { kiss_fft(fft->invs,pp[i3][i2],pp[i3][i2]); } } /* scaling */ for (i3=0; i3 < fft->n3; i3++) { for (i2=0; i2 < fft->n2; i2++) { for(i1=0; i1 < fft->n1; i1++) { pp[i3][i2][i1] = sf_crmul(pp[i3][i2][i1],fft->scale); } } } } else { /* scaling */ for (i3=0; i3 < fft->n3; i3++) { for (i2=0; i2 < fft->n2; i2++) { for(i1=0; i1 < fft->n1; i1++) { pp[i3][i2][i1] = sf_crmul(pp[i3][i2][i1],fft->scale); } } } /* FFT 1 */ for (i3=0; i3 < fft->n3; i3++) { for(i2=0; i2 < fft->n2; i2++) { kiss_fft(fft->forw,pp[i3][i2],pp[i3][i2]); } } } }
static void init_wave(int init, int nx, float dx, int nz, float dz, sf_complex *pp /* [nx] */, float wov, int nw, int iw) { int ix; float x,x0,z0,phase,amp; x0 = nx*dx/3; z0 = nz*dz/3; switch(init) { case 1: /* planar wave @ 15deg */ for (ix=0; ix < nx; ix++) { x = (ix+1)*dx - x0; phase = wov*x*sinf(15*SF_PI/180.); pp[ix] = cexpf(sf_cmplx(0.,phase)); } break; case 2: /* expanding spherical wave */ for (ix=0; ix < nx; ix++) { x = (ix+1)*dx - x0; phase = wov*hypotf(z0,x); pp[ix] = cexpf(sf_cmplx(0.,phase)); } break; case 3: /* point source */ for (ix=0; ix < nx; ix++) { pp[ix]=sf_cmplx(0.,0.); } pp[nx/3-1] = sf_cmplx(1.,0.); break; case 4: /* collapsing spherical wave */ for (ix=0; ix < nx; ix++) { x = (ix+1)*dx - x0; phase = -wov*hypotf(z0,x); pp[ix] = cexpf(sf_cmplx(0.,phase)); } break; default: sf_error("Unknown init=%d",init); } amp = (nw-iw+1.0)/nw; amp = cosf((1-amp)*(0.5*SF_PI)); amp *= amp; for (ix=0; ix < nx; ix++) { #ifdef SF_HAS_COMPLEX_H pp[ix] *= amp; #else pp[ix] = sf_crmul(pp[ix],amp); #endif } }
void rweone_fft( bool inv, kiss_fft_cpx *d) /*< apply FFT >*/ { int ig; if(inv) { kiss_fft(invs,d,d); for(ig=0;ig<ag.n;ig++) { d[ig] = sf_crmul(d[ig],ffts); } } else { for(ig=0;ig<ag.n;ig++) { d[ig] = sf_crmul(d[ig],ffts); } kiss_fft(forw,d,d); } }
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_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 rweone_tap(sf_complex *v) /*< apply taper >*/ { int ig; for(ig=0;ig<ag.n;ig++) { #ifdef SF_HAS_COMPLEX_H v[ig] *= tap[ig]; #else v[ig] = sf_crmul(v[ig],tap[ig]); #endif } }
void icfft3(sf_complex *out /* [n1*n2*n3] */, sf_complex *inp /* [nk*n2*n3] */) /*< 3-D inverse FFT >*/ { int i1, i2, i3; #ifdef SF_HAS_FFTW fftwf_execute(icfg); #else /* IFFT over third axis */ for (i2=0; i2 < n2; i2++) { for (i1=0; i1 < nk; i1++) { kiss_fft_stride(icfg3,(kiss_fft_cpx *) (inp+i2*nk+i1),ctrace3,nk*n2); for (i3=0; i3<n3; i3++) { tmp[i3][i2][i1] = ctrace3[i3]; } } } /* IFFT over second axis */ for (i3=0; i3 < n3; i3++) { for (i1=0; i1 < nk; i1++) { kiss_fft_stride(icfg2,tmp[i3][0]+i1,ctrace2,nk); for (i2=0; i2<n2; i2++) { tmp[i3][i2][i1] = ctrace2[i2]; } } } /* IFFT over first axis */ for (i3=0; i3 < n3; i3++) { for (i2=0; i2 < n2; i2++) { kiss_fft_stride(icfg1,tmp[i3][i2],(kiss_fft_cpx *) cc[i3][i2],1); } } #endif /* FFT centering and normalization */ for (i3=0; i3<n3; i3++) { for (i2=0; i2<n2; i2++) { for (i1=0; i1<n1; i1++) { #ifdef SF_HAS_COMPLEX_H out[(i3*n2+i2)*n1+i1] = ((((i3%2==0)==(i2%2==0))==(i1%2==0))? wt:-wt)*cc[i3][i2][i1]; #else out[(i3*n2+i2)*n1+i1] = sf_crmul(cc[i3][i2][i1],((((i3%2==0)==(i2%2==0))==(i1%2==0))? wt:-wt)); #endif } } } }
sf_complex Green(float r1,float r2,float r3,float s1,float s2,float s3,float omega) /*< Green's function >*/ { double tt,amp; sf_complex val; GreenTtAmp(r1,r2,r3,s1,s2,s3,&tt,&); #ifdef SF_HAS_COMPLEX_H val=amp*cexpf(sf_cmplx(0.,omega*tt)); #else val=sf_crmul(cexpf(sf_cmplx(0.,omega*tt)),amp); #endif return (val); }
void taper2(sf_complex** tt /* [n2][n1] tapered array (in and out) */) /*< 2-D taper >*/ { int it,i2,i1; float gain; for (it=0; it < nt2; it++) { gain = tap2[it]; for (i1=0; i1 < n1; i1++) { #ifdef SF_HAS_COMPLEX_H if (b2) tt[ it ][i1] *= gain; ; tt[n2-it-1][i1] *= gain; #else if (b2) tt[ it ][i1] = sf_crmul(tt[ it ][i1],gain); ; tt[n2-it-1][i1] = sf_crmul(tt[n2-it-1][i1],gain); #endif } } for (it=0; it < nt1; it++) { gain = tap1[it]; for (i2=0; i2 < n2; i2++) { #ifdef SF_HAS_COMPLEX_H if (b1) tt[i2][ it ] *= gain; ; tt[i2][n1-it-1] *= gain; #else if (b1) tt[i2][ it ] = sf_crmul(tt[i2][ it ],gain); ; tt[i2][n1-it-1] = sf_crmul(tt[i2][n1-it-1],gain); #endif } } }
void cblas_csscal(int n, float alpha, void *x, int sx) /*< x = alpha*x >*/ { int i, ix; sf_complex* c; c = (sf_complex*) x; for (i=0; i < n; i++) { ix = i*sx; #ifdef SF_HAS_COMPLEX_H c[ix] *= alpha; #else c[ix] = sf_crmul(c[ix],alpha); #endif } }
void icfft2(sf_complex *out /* [n1*local_n0] */, sf_complex *inp /* [nk*local_n0] */) /*< 2-D inverse FFT >*/ { int i1, i2, ith=0; #ifdef _OPENMP #pragma omp parallel for private(i2,ith) default(shared) #endif for (i2=0; i2 < local_n0; i2++) { #ifdef _OPENMP ith = omp_get_thread_num(); #endif kiss_fft_stride(icfg1[ith],(kiss_fft_cpx *) inp+i2*nk,tmp+i2*nk,1); } fftwf_execute(cfg); #ifdef _OPENMP #pragma omp parallel for private(i1,i2,ith) default(shared) #endif for (i1=0; i1 < local_n1; i1++) { #ifdef _OPENMP ith = omp_get_thread_num(); #endif kiss_fft_stride(icfg2[ith],tmp+i1*n2,ctrace2[ith],1); for (i2=0; i2<n2; i2++) { tmp[i1*n2+i2] = ctrace2[ith][i2]; } } fftwf_execute(icfg); /* FFT centering and normalization*/ #pragma omp parallel for private(i2,i1) default(shared) for (i2=0; i2<local_n0; i2++) { for (i1=0; i1<n1; i1++) { #ifdef SF_HAS_COMPLEX_H out[i2*n1+i1] = ((((i2+local_0_start)%2==0)==(i1%2==0))? wt:-wt) * tmp2[i2*n1+i1]; #else out[i2*n1+i1] = sf_crmul(tmp2[i2*n1+i1],((((i2+local_0_start)%2==0)==(i1%2==0))? wt:-wt)); #endif } } }
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_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 icfft2(sf_complex *out /* [n1*n2] */, sf_complex *inp /* [nk*n2] */) /*< 2-D inverse FFT >*/ { int i1, i2; #ifdef SF_HAS_FFTW #ifdef _OPENMP #pragma omp parallel for private(i2,i1) default(shared) #endif for (i2=0; i2<n2; i2++) { for (i1=0; i1<nk; i1++) { dd[i2][i1]=inp[i2*nk+i1]; } } fftwf_execute(icfg); #else for (i1=0; i1 < nk; i1++) { kiss_fft_stride(icfg2,(kiss_fft_cpx *) (inp+i1),ctrace2,nk); for (i2=0; i2<n2; i2++) { tmp[i2][i1] = ctrace2[i2]; } } for (i2=0; i2 < n2; i2++) { kiss_fft_stride(icfg1,tmp[i2],(kiss_fft_cpx *) cc[i2],1); } #endif /* FFT centering and normalization*/ #ifdef _OPENMP #pragma omp parallel for private(i2,i1) default(shared) #endif for (i2=0; i2<n2; i2++) { for (i1=0; i1<n1; i1++) { #ifdef SF_HAS_COMPLEX_H out[i2*n1+i1] = (((i2%2==0)==(i1%2==0))? wt:-wt) * cc[i2][i1]; #else out[i2*n1+i1] = sf_crmul(cc[i2][i1],(((i2%2==0)==(i1%2==0))? wt:-wt)); #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 rweone_mrs( sf_complex *v, sf_complex *d, float *m, int ir) /*< combine MRS >*/ { int ig, iloop; for(ig=0;ig<ag.n;ig++) { if(m[ig]-1 == ir) { mtt[ig] = 1.; } else { mtt[ig] = 0.; } } for(iloop=0;iloop<nloop;iloop++) { msk[0] = mtt[1]; for(ig=1;ig<ag.n-1;ig++) { msk[ig] = 0.5 * (mtt[ig-1] + mtt[ig+1] ); } msk[ag.n-1] = mtt[ag.n-2]; for(ig=0;ig<ag.n;ig++) { mtt[ig] = msk[ig]; } } for(ig=0;ig<ag.n;ig++) { #ifdef SF_HAS_COMPLEX_H d[ig] += msk[ig] * v[ig]; #else d[ig] = sf_cadd(d[ig],sf_crmul(v[ig],msk[ig])); #endif } }
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 main(int argc, char* argv[]) { int i, n, n1; float *dat=NULL, *adat=NULL, t, pclip, d; sf_complex *cdat=NULL; sf_file in=NULL, out=NULL; sf_init(argc,argv); in = sf_input("in"); out = sf_output("out"); n = sf_filesize(in); adat = sf_floatalloc(n); if (!sf_getfloat("pclip",&pclip)) sf_error("Need pclip="); /* percentage to clip */ n1 = 0.5+n*(1.-0.01*pclip); if (n1 < 0) n1=0; if (n1 >= n) n1=n-1; if (SF_FLOAT == sf_gettype(in)) { dat = sf_floatalloc(n); sf_floatread(dat,n,in); for (i=0; i < n; i++) { adat[i] = fabsf(dat[i]); } } else if (SF_COMPLEX == sf_gettype(in)) { cdat = sf_complexalloc(n); sf_complexread(cdat,n,in); for (i=0; i < n; i++) { adat[i] = cabsf(cdat[i]); } } else { sf_error("Need float or complex input"); } t = sf_quantile(n1,n,adat); if (NULL != dat) { for (i=0; i < n; i++) { d = dat[i]; if (d < -t) { dat[i] = d+t; } else if (d > t) { dat[i] = d-t; } else { dat[i] = 0.; } } sf_floatwrite(dat,n,out); } else { for (i=0; i < n; i++) { d = cabsf(cdat[i]); if (d < -t) { #ifdef SF_HAS_COMPLEX_H cdat[i] *= (d+t)/d; #else cdat[i] = sf_crmul(cdat[i],(d+t)/d); #endif } else if (d > t) { #ifdef SF_HAS_COMPLEX_H cdat[i] *= (d-t)/d; #else cdat[i] = sf_crmul(cdat[i],(d-t)/d); #endif } else { cdat[i] = sf_cmplx(0.,0.); } } sf_complexwrite(cdat,n,out); } exit(0); }
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; }