void preparation_s(sf_file Fv, sf_file Fw, sf_acqui acpar, sf_sou soupar, sf_vec_s array, sf_seis seispar) /*< read data, initialize variables and prepare acquisition geometry >*/ { int i, nb, nzx; float sx, xend, rbegin, rend, tmp; int nplo=3, nphi=3, nt; float eps=0.0001; sf_butter blo=NULL, bhi=NULL; /* absorbing boundary coefficients */ nb=acpar->nb; acpar->bc=sf_floatalloc(nb); for(i=0; i<nb; i++){ tmp=acpar->coef*(nb-i); acpar->bc[i]=expf(-tmp*tmp); } /* padding variables */ acpar->padnx=acpar->nx+2*nb; acpar->padnz=acpar->nz+2*nb; acpar->padx0=acpar->x0-nb*acpar->dx; acpar->padz0=acpar->z0-nb*acpar->dz; /* acquisition parameters */ acpar->ds_v=acpar->ds/acpar->dx+0.5; acpar->s0_v=(acpar->s0-acpar->x0)/acpar->dx+0.5+nb; acpar->sz += nb; acpar->dr_v=acpar->dr/acpar->dx+0.5; acpar->r0_v=sf_intalloc(acpar->ns); acpar->r02=sf_intalloc(acpar->ns); acpar->nr2=sf_intalloc(acpar->ns); acpar->rz += nb; xend=acpar->x0+(acpar->nx-1)*acpar->dx; if(acpar->acqui_type==1){ for(i=0; i<acpar->ns; i++){ acpar->r0_v[i]=acpar->r0/acpar->dx+0.5+nb; acpar->r02[i]=0; acpar->nr2[i]=acpar->nr; } }else{ for(i=0; i<acpar->ns; i++){ sx=acpar->s0+acpar->ds*i; rbegin=(sx+acpar->r0 <acpar->x0)? acpar->x0 : sx+acpar->r0; rend=sx+acpar->r0 +(acpar->nr-1)*acpar->dr; rend=(rend < xend)? rend : xend; acpar->r0_v[i]=(rbegin-acpar->x0)/acpar->dx+0.5+nb; acpar->r02[i]=(rbegin-sx-acpar->r0)/acpar->dx+0.5; acpar->nr2[i]=(rend-rbegin)/acpar->dr+1.5; } } /* read model parameters */ nzx=acpar->nz*acpar->nx; nt=acpar->nt; array->vv=sf_floatalloc(nzx); array->ww=sf_floatalloc(nt); sf_floatread(array->vv, nzx, Fv); sf_floatread(array->ww, nt, Fw); /* bandpass the wavelet */ soupar->flo *= acpar->dt; soupar->fhi *= acpar->dt; if(soupar->flo > eps) blo=sf_butter_init(false, soupar->flo, nplo); if(soupar->fhi < 0.5-eps) bhi=sf_butter_init(true, soupar->fhi, nphi); if(NULL != blo){ sf_butter_apply(blo, nt, array->ww); sf_reverse(nt, array->ww); sf_butter_apply(blo, nt, array->ww); sf_reverse(nt, array->ww); sf_butter_close(blo); } if(NULL != bhi){ sf_butter_apply(bhi, nt, array->ww); sf_reverse(nt, array->ww); sf_butter_apply(bhi, nt, array->ww); sf_reverse(nt, array->ww); sf_butter_close(bhi); } /* simultaneous source */ nsource=soupar->nsource; dsource=soupar->dsource; /* seislet regularization */ if(seispar!=NULL){ seislet=true; seislet_init(acpar->nz, acpar->nx, true, true, seispar->eps, seispar->order, seispar->type[0]); seislet_set(seispar->dip); pclip=seispar->pclip; ss=sf_floatalloc(nzx); } }
void preparation_q(sf_file Fv, sf_file Fq, sf_file Ftau, sf_file Fw, sf_acqui acpar, sf_sou soupar, sf_vec_q array, float f0) /*< read data, initialize variables and prepare acquisition geometry >*/ { int i, nb, nzx; float sx, xend, rbegin, rend, tmp; float *qq; int nplo=3, nphi=3, nt; float eps=0.0001; sf_butter blo=NULL, bhi=NULL; /* absorbing boundary coefficients */ nb=acpar->nb; acpar->bc=sf_floatalloc(nb); for(i=0; i<nb; i++){ tmp=acpar->coef*(nb-i); acpar->bc[i]=expf(-tmp*tmp); } /* padding variables */ acpar->padnx=acpar->nx+2*nb; acpar->padnz=acpar->nz+2*nb; acpar->padx0=acpar->x0-nb*acpar->dx; acpar->padz0=acpar->z0-nb*acpar->dz; /* acquisition parameters */ acpar->ds_v=acpar->ds/acpar->dx+0.5; acpar->s0_v=(acpar->s0-acpar->x0)/acpar->dx+0.5+nb; acpar->sz += nb; acpar->dr_v=acpar->dr/acpar->dx+0.5; acpar->r0_v=sf_intalloc(acpar->ns); acpar->r02=sf_intalloc(acpar->ns); acpar->nr2=sf_intalloc(acpar->ns); acpar->rz += nb; xend=acpar->x0+(acpar->nx-1)*acpar->dx; if(acpar->acqui_type==1){ for(i=0; i<acpar->ns; i++){ acpar->r0_v[i]=acpar->r0/acpar->dx+0.5+nb; acpar->r02[i]=0; acpar->nr2[i]=acpar->nr; } }else{ for(i=0; i<acpar->ns; i++){ sx=acpar->s0+acpar->ds*i; rbegin=(sx+acpar->r0 <acpar->x0)? acpar->x0 : sx+acpar->r0; rend=sx+acpar->r0 +(acpar->nr-1)*acpar->dr; rend=(rend < xend)? rend : xend; acpar->r0_v[i]=(rbegin-acpar->x0)/acpar->dx+0.5+nb; acpar->r02[i]=(rbegin-sx-acpar->r0)/acpar->dx+0.5; acpar->nr2[i]=(rend-rbegin)/acpar->dr+1.5; } } /* read model parameters */ nzx=acpar->nz*acpar->nx; nt=acpar->nt; array->vv=sf_floatalloc(nzx); qq=sf_floatalloc(nzx); array->tau=sf_floatalloc(nzx); array->taus=sf_floatalloc(nzx); array->ww=sf_floatalloc(nt); sf_floatread(array->vv, nzx, Fv); sf_floatread(qq, nzx, Fq); sf_floatread(array->tau, nzx, Ftau); sf_floatread(array->ww, nt, Fw); /* calculate taus */ for(i=0; i<nzx; i++){ array->taus[i]=(sqrtf(qq[i]*qq[i]+1)-1.)/(2.*SF_PI*f0*qq[i]); } /* bandpass the wavelet */ soupar->flo *= acpar->dt; soupar->fhi *= acpar->dt; if(soupar->flo > eps) blo=sf_butter_init(false, soupar->flo, nplo); if(soupar->fhi < 0.5-eps) bhi=sf_butter_init(true, soupar->fhi, nphi); if(NULL != blo){ sf_butter_apply(blo, nt, array->ww); sf_reverse(nt, array->ww); sf_butter_apply(blo, nt, array->ww); sf_reverse(nt, array->ww); sf_butter_close(blo); } if(NULL != bhi){ sf_butter_apply(bhi, nt, array->ww); sf_reverse(nt, array->ww); sf_butter_apply(bhi, nt, array->ww); sf_reverse(nt, array->ww); sf_butter_close(bhi); } free(qq); }
//JS //static int counter=0; void gradient_pas_av(float *x, float *fcost, float *grad) /*< acoustic velocity gradient >*/ { int ix, iz, is, ir, it, wit, iturn; int rx; float temp, dmax; float **p0, **p1, **p2, **term, **tmparray, ***wave, **pp; float *sendbuf, *recvbuf; /* //JS sf_file Fwfl1, Fwfl2, Fres; counter++; if (cpuid==0 && counter==3) { Fwfl1=sf_output("Fwfl1"); Fwfl2=sf_output("Fwfl2"); Fres=sf_output("Fres"); sf_putint(Fwfl1,"n1",padnz); sf_putint(Fwfl1,"n2",padnx); sf_putint(Fwfl1,"n3",(nt-1)/50+1); sf_putint(Fwfl2,"n1",padnz); sf_putint(Fwfl2,"n2",padnx); sf_putint(Fwfl2,"n3",(nt-1)/50+1); sf_putint(Fres,"n1",nt); sf_putint(Fres,"n2",nr); } */ /* initialize fcost */ *fcost=0.; /* update velocity */ pad2d(x, vv, nz, nx, nb); /* initialize gradient */ memset(grad, 0., nzx*sizeof(float)); /* memory allocation */ p0=sf_floatalloc2(padnz, padnx); p1=sf_floatalloc2(padnz, padnx); p2=sf_floatalloc2(padnz, padnx); term=sf_floatalloc2(padnz, padnx); wave=sf_floatalloc3(nz, nx, wnt); pp=sf_floatalloc2(nt, nr); iturn=0; for(is=cpuid; is<ns; is+=numprocs){ if(cpuid==0) sf_warning("###### is=%d ######", is+1); memset(p0[0], 0., padnzx*sizeof(float)); memset(p1[0], 0., padnzx*sizeof(float)); memset(p2[0], 0., padnzx*sizeof(float)); memset(pp[0], 0., nr*nt*sizeof(float)); wit=0; /* forward propagation */ for(it=0; it<nt; it++){ if(verb) sf_warning("Forward propagation it=%d;", it); /* output predicted data */ for(ir=0; ir<nr; ir++){ rx=r0_v[0]+ir*dr_v; pp[ir][it]=p1[rx][rz]; } /* save wavefield */ if(it%interval==0){ #ifdef _OPENMP #pragma omp parallel for \ private(ix,iz) \ shared(wave,p1,wit,nb,nx,nz) #endif for(ix=0; ix<nx; ix++) for(iz=0; iz<nz; iz++) wave[wit][ix][iz]=p1[ix+nb][iz+nb]; wit++; } /* //JS if(is==0 && counter==3 && it%50==0) sf_floatwrite(p1[0],padnzx,Fwfl1); */ /* laplacian operator */ laplace(p1, term, padnx, padnz, dx2, dz2); /* load source */ #ifdef _OPENMP #pragma omp parallel for \ private(ix,iz) \ shared(term,nb,ww3,it) #endif for(ix=0; ix<nx; ix++){ for(iz=0; iz<nz; iz++){ term[ix+nb][iz+nb] += ww3[iturn][it][ix][iz]; } } /* update */ #ifdef _OPENMP #pragma omp parallel for \ private(ix,iz) \ shared(p0,p1,p2,vv,term,padnx,padnz,dt2) #endif for(ix=4; ix<padnx-4; ix++){ for(iz=4; iz<padnz-4; iz++){ p2[ix][iz]=2*p1[ix][iz]-p0[ix][iz]+vv[ix][iz]*vv[ix][iz]*dt2*term[ix][iz]; } } /* swap wavefield pointer of different time steps */ tmparray=p0; p0=p1; p1=p2; p2=tmparray; /* boundary condition */ apply_sponge(p0, bc, padnx, padnz, nb); apply_sponge(p1, bc, padnx, padnz, nb); } // end of time loop /* check */ if(wit != wnt) sf_error("Incorrect number of wavefield snapshots"); wit--; /* calculate data residual and data misfit */ for(ir=0; ir<nr; ir++){ for(it=0; it<nt; it++){ pp[ir][it]=dd[iturn][ir][it]-pp[ir][it]; *fcost += 0.5*pp[ir][it]*pp[ir][it]; } } /* window the data residual */ for(ir=0; ir<nr; ir++){ /* multiscale */ if(NULL != blo){ sf_butter_apply(blo, nt, pp[ir]); sf_reverse(nt, pp[ir]); sf_butter_apply(blo, nt, pp[ir]); sf_reverse(nt, pp[ir]); } if(NULL != bhi){ sf_butter_apply(bhi, nt, pp[ir]); sf_reverse(nt, pp[ir]); sf_butter_apply(bhi, nt, pp[ir]); sf_reverse(nt, pp[ir]); } for(it=0; it<nt; it++){ pp[ir][it] *= weight[ir][it]; } } /* // JS if(is==0 && counter==3) sf_floatwrite(pp[0], nr*nt, Fres); */ /* initialization */ memset(p0[0], 0., padnzx*sizeof(float)); memset(p1[0], 0., padnzx*sizeof(float)); memset(p2[0], 0., padnzx*sizeof(float)); memset(term[0], 0., padnzx*sizeof(float)); /* backward propagation */ for(it=nt-1; it>=0; it--){ if(verb) sf_warning("Backward propagation it=%d;", it); /* laplacian operator */ laplace(p1, term, padnx, padnz, dx2, dz2); /* load data residual*/ for(ir=0; ir<nr; ir++){ rx=r0_v[0]+ir*dr_v; term[rx][rz] += pp[ir][it]; } /* update */ #ifdef _OPENMP #pragma omp parallel for \ private(ix,iz) \ shared(p0,p1,p2,vv,term,padnx,padnz,dt2) #endif for(ix=4; ix<padnx-4; ix++){ for(iz=4; iz<padnz-4; iz++){ p2[ix][iz]=2*p1[ix][iz]-p0[ix][iz]+vv[ix][iz]*vv[ix][iz]*dt2*term[ix][iz]; } } /* // JS if(is==0 && counter==3 && it%50==0) sf_floatwrite(p1[0],padnzx,Fwfl2); */ /* calculate gradient */ if(it%interval==0){ if(wit != wnt-1 && wit != 0){ // avoid the first and last time step #ifdef _OPENMP #pragma omp parallel for \ private(ix,iz,temp) \ shared(nx,nz,vv,wave,p1,wit,wdt2,grad) #endif for(ix=0; ix<nx; ix++){ for(iz=0; iz<nz; iz++){ temp=vv[ix+nb][iz+nb]; temp=temp*temp*temp; temp=-2./temp; grad[ix*nz+iz] += gwt[iturn][it][ix][iz]*(wave[wit+1][ix][iz]-2.*wave[wit][ix][iz]+wave[wit-1][ix][iz])/wdt2*p1[ix+nb][iz+nb]*temp; } } } wit--; } /* swap wavefield pointer of different time steps */ tmparray=p0; p0=p1; p1=p2; p2=tmparray; /* boundary condition */ apply_sponge(p0, bc, padnx, padnz, nb); apply_sponge(p1, bc, padnx, padnz, nb); } // end of time loop iturn++; /* // JS if(is==0 && counter==3){ sf_fileclose(Fwfl1); sf_fileclose(Fwfl2); sf_fileclose(Fres); } sf_warning("---counter=%d fcost=%3.3e---",counter, *fcost); */ } // end of shot loop MPI_Barrier(comm); /* misfit reduction */ if(cpuid==0){ sendbuf=MPI_IN_PLACE; recvbuf=fcost; }else{ sendbuf=fcost; recvbuf=NULL; } MPI_Reduce(sendbuf, recvbuf, 1, MPI_FLOAT, MPI_SUM, 0, comm); MPI_Bcast(fcost, 1, MPI_FLOAT, 0, comm); /* gradient reduction */ if(cpuid==0){ sendbuf=MPI_IN_PLACE; recvbuf=grad; }else{ sendbuf=grad; recvbuf=NULL; } MPI_Reduce(sendbuf, recvbuf, nzx, MPI_FLOAT, MPI_SUM, 0, comm); MPI_Bcast(grad, nzx, MPI_FLOAT, 0, comm); /* scaling gradient */ if(first){ dmax=0.; for(ix=0; ix<nzx; ix++) if(fabsf(grad[ix])>dmax) dmax=fabsf(grad[ix]); scaling=0.1/dmax; first=false; } /* smooth gradient */ gradient_smooth2b(grectx, grectz, nx, nz, waterz, waterzb, scaling, grad); /* free allocated memory */ free(*p0); free(p0); free(*p1); free(p1); free(*p2); free(p2); free(*pp); free(pp); free(**wave); free(*wave); free(wave); free(*term); free(term); }