Esempio n. 1
0
int main(int argc, char* argv[])
{
 
    bool hermite_false, hermite_true;
    int n1, n2, npml, pad1, pad2, ns, nw;
    float d1, d2, **v, ds, os, dw, ow;
    sf_complex ****f,  ****obs; 
    sf_file in, out, misfit, source, receiver, record;
    char *order;
    int uts, mts, is, i, j, iw, iter, niter;
    float **recloc;
    float **m_old, **m_new;
    float **d_new, **d_old;
    float **g_old, **g_new; 
    sf_complex ****r_new, ****r_old;
    sf_complex ****Fg; 
    float alpha, beta, gnorm, rnorm; 
    float *datamisfit;

    sf_init(argc, argv);

    in = sf_input("in");
    out = sf_output("out");
    misfit = sf_output("misfit");

    if (!sf_getint("uts",&uts)) uts=0;
//#ifdef _OPENMP
//    mts = omp_get_max_threads();
//#else
    mts = 1;
//#endif
    uts = (uts < 1)? mts: uts;

    hermite_false=false;
    hermite_true=true;
    /* Hermite operator */
    
    if (!sf_getint("npml",&npml)) npml=20;
    /* PML width */
    
    if (!sf_getint("niter",&niter)) niter=0; 
    /* Number of iterations */

    if (NULL == (order = sf_getstring("order"))) order="c";
    /* discretization scheme (default optimal 9-point) */

    fdprep_order(order);

    /* read input dimension */
    if (!sf_histint(in,"n1",&n1)) sf_error("No n1= in input.");
    if (!sf_histint(in,"n2",&n2)) sf_error("No n2= in input.");

    if (!sf_histfloat(in,"d1",&d1)) sf_error("No d1= in input.");
    if (!sf_histfloat(in,"d2",&d2)) sf_error("No d2= in input.");

    v = sf_floatalloc2(n1,n2);
    sf_floatread(v[0],n1*n2,in);
    
	/* PML padding */
	pad1 = n1+2*npml;
	pad2 = n2+2*npml;

    /* read receiver */ 
    if (NULL == sf_getstring("receiver")) sf_error("Need receiver="); 
    receiver = sf_input("receiver"); 
    recloc=sf_floatalloc2(n1,n2);
    sf_floatread(recloc[0],n1*n2,receiver);

    /* read source */
    if (NULL == sf_getstring("source")) sf_error("Need source=");
    source = sf_input("source");

    if (!sf_histint(source,"n3",&ns)) sf_error("No ns=.");
    if (!sf_histfloat(source,"d3",&ds)) ds=d2;
    if (!sf_histfloat(source,"o3",&os)) os=0.;

    /* read observed data */
    if (NULL == sf_getstring("record")) sf_error("Need record=");
    record = sf_input("record");

    if (!sf_histint(record,"n4",&nw)) sf_error("No nw=.");
    if (!sf_histfloat(record,"d4",&dw)) sf_error("No dw=.");
    if (!sf_histfloat(record,"o4",&ow)) sf_error("No ow=."); 

    f = sf_complexalloc4(n1,n2,ns,nw);
    obs = sf_complexalloc4(n1,n2,ns,nw);

    sf_complexread(f[0][0][0],n1*n2*ns*nw,source);
    sf_complexread(obs[0][0][0],n1*n2*ns*nw,record);
   
    /* allocate variables */
    m_old = sf_floatalloc2(n1,n2);
    m_new = sf_floatalloc2(n1,n2);
    d_old = sf_floatalloc2(n1,n2);
    d_new = sf_floatalloc2(n1,n2);
    g_old = sf_floatalloc2(n1,n2);
    g_new = sf_floatalloc2(n1,n2);
    
    r_old = sf_complexalloc4(n1,n2,ns,nw);
    r_new = sf_complexalloc4(n1,n2,ns,nw);
    Fg = sf_complexalloc4(n1,n2,ns,nw);

    /* set output dimension */
    sf_putint(out,"n1",n1);
    sf_putint(out,"n2",n2);
    sf_putint(out,"n3",niter);
    
    sf_putint(misfit,"n1",niter);
    sf_putfloat(misfit,"d1",1);
    sf_putint(misfit,"n2",1);
    datamisfit = sf_floatalloc(niter);

    rnorm = 0.0; 
    for ( iw = 0; iw < nw; iw ++ ) { 
    for ( is = 0; is < ns; is ++) { 
    for ( j = 0; j < n2 ; j++ ) { 
    for ( i = 0; i < n1; i++ ) { 
        r_old[iw][is][j][i] = obs[iw][is][j][i];
        if ( recloc[j][i] > 0.0 ) { 
            rnorm += cabsf( r_old[iw][is][j][i] * r_old[iw][is][j][i] ); 
        }
    }
    }
    }
    }
    sf_warning("rnorm = %g.",rnorm);


    sf_warning("Adjoint calculation for the first iteration.");
    /* adjoint calculation */
    adjlsm_operator(nw, ow, dw, ns, n1, n2, 
                     uts, pad1, pad2, npml, d1, d2,
                     hermite_false, hermite_true, v, f, recloc, 
                     r_old, g_old);

    /* set starting valuables */
    for (j = 0; j < n2; j++) { 
    for (i = 0; i < n2; i++ ) { 
        d_old[j][i] = g_old[j][i];
        m_old[j][i] = 0.0;
    }
    }
    
    for ( iter = 0; iter < niter; iter ++ ) { 

        sf_warning("Calculating iteration %d out of %d.",iter,niter);

        /* born forward operator */
        bornsyn_operator(nw, ow, dw, ns, n1, n2, 
                         uts, pad1, pad2, npml, d1, d2,
                         hermite_false, v, f, recloc, 
                         d_old, Fg); 

        /* calculate alpha value */
        alpha = calc_alpha(g_old, Fg, recloc, n1, n2, ns, nw);
        
        /* update model */ 
        update_model_lsm(m_old, m_new, d_old, alpha, n1, n2); 

        /* update residual */ 
        update_residual(r_old, r_new, Fg, recloc, alpha, n1, n2, ns, nw);  

        /* adjoint operator */ 
        adjlsm_operator(nw, ow, dw, ns, n1, n2, 
                         uts, pad1, pad2, npml, d1, d2,
                         hermite_false, hermite_true, v, f, recloc, 
                         r_new, g_new);

        /* update direction */
        beta = direction_cg_fletcher(g_old, d_old, g_new, d_new, n1, n2); 

        sf_warning("alpha = %g, beta = %g.",alpha, beta);

        /* update vectors */
        gnorm = 0.0 ;  
        for (j = 0; j < n2; j++ ) { 
        for (i = 0; i < n1; i++ ) { 
            d_old[j][i] = d_new[j][i]; 
            g_old[j][i] = g_new[j][i];
            m_old[j][i] = m_new[j][i]; 
            gnorm += g_old[j][i] * g_old[j][i]; 
        }
        }

        rnorm = 0.0 ; 
        for (iw = 0; iw < nw; iw ++ ) { 
        for (is = 0; is < ns; is ++ ) { 
        for (j = 0; j < n2; j ++ ) { 
        for (i = 0; i < n1; i ++ ) { 
            r_old[iw][is][j][i] = r_new[iw][is][j][i]; 
            if ( recloc[j][i] > 0.0 ) { 
                rnorm += cabsf( r_old[iw][is][j][i] * r_old[iw][is][j][i] ); 
            }
        }
        }
        }
        }
        sf_warning("gnorm = %g; rnorm = %g.",gnorm, rnorm);

        datamisfit[iter] = rnorm; 

        sf_floatwrite(m_old[0],n1*n2,out);

    }  /* end iteration */

    sf_floatwrite(datamisfit, niter, misfit); 

    exit(0);
}
Esempio n. 2
0
int main(int argc, char* argv[])
{
    bool verb, save, load;
    int npml, pad1, pad2, n1, n2; 
    int ih, nh, is, ns, iw, nw, i, j;
    SuiteSparse_long n, nz, *Ti, *Tj;
    float d1, d2, **vel, ****image, ****timage, dw, ow;
    double omega, *Tx, *Tz;
    SuiteSparse_long *Ap, *Ai, *Map;
    double *Ax, *Az, **Xx, **Xz, **Bx, **Bz;
    void *Symbolic, **Numeric;
    double Control[UMFPACK_CONTROL];
    sf_complex ***srce, ***recv;
    char *datapath, *insert, *append;
    size_t srclen, inslen;
    sf_file in, out, source, data, us, ur, timg;
    int uts, its, mts;
    sf_timer timer;
    char *order;

    sf_init(argc,argv);
    in  = sf_input("in");
    out = sf_output("out");

    if (!sf_getbool("verb",&verb)) verb=false;
    /* verbosity flag */

    if (verb)
	timer = sf_timer_init();
    else
	timer = NULL;

    if (!sf_getbool("save",&save)) save=false;
    /* save LU */

    if (!sf_getbool("load",&load)) load=false;
    /* load LU */

    if (save || load) {
	datapath = sf_histstring(in,"in");
	srclen = strlen(datapath);
	insert = sf_charalloc(6);
    } else {
	datapath = NULL;
	srclen = 0;
	insert = NULL;
	append = NULL;
    }

    if (!sf_getint("uts",&uts)) uts=0;
    /* number of OMP threads */

#ifdef _OPENMP
    mts = omp_get_max_threads();
#else
    mts = 1;
#endif

    uts = (uts < 1)? mts: uts;
    if (verb) sf_warning("Using %d out of %d threads.",uts,mts);

    if (!sf_getint("nh",&nh)) nh=0;
    /* horizontal space-lag */

    if (!sf_getint("npml",&npml)) npml=10;
    /* PML width */

    if (NULL == (order = sf_getstring("order"))) order="j";
    /* discretization scheme (default optimal 9-point) */

    fdprep_order(order);

    /* read model */
    if (!sf_histint(in,"n1",&n1)) sf_error("No n1= in input.");
    if (!sf_histint(in,"n2",&n2)) sf_error("No n2= in input.");

    if (!sf_histfloat(in,"d1",&d1)) sf_error("No d1= in input.");
    if (!sf_histfloat(in,"d2",&d2)) sf_error("No d2= in input.");

    vel = sf_floatalloc2(n1,n2);
    sf_floatread(vel[0],n1*n2,in);

    /* read source */
    if (NULL == sf_getstring("source"))
	sf_error("Need source=");
    source = sf_input("source");

    if (!sf_histint(source,"n3",&ns)) sf_error("No ns=.");
    if (!sf_histint(source,"n4",&nw)) sf_error("No nw=.");
    if (!sf_histfloat(source,"d4",&dw)) sf_error("No dw=.");
    if (!sf_histfloat(source,"o4",&ow)) sf_error("No ow=.");

    srce = sf_complexalloc3(n1,n2,ns);

    /* read receiver */
    if (NULL == sf_getstring("data"))
	sf_error("Need data=");
    data = sf_input("data");

    recv = sf_complexalloc3(n1,n2,ns);

    /* write output header */
    sf_putint(out,"n3",2*nh+1);
    sf_putfloat(out,"d3",d2);
    sf_putfloat(out,"o3",(float) -nh*d2);

    /* output source wavefield */
    if (NULL != sf_getstring("us")) {
	us = sf_output("us");
	
	sf_settype(us,SF_COMPLEX);
	sf_putint(us,"n3",ns);
	sf_putstring(us,"label3","Shot");
	sf_putstring(us,"unit3","");
	sf_putint(us,"n4",nw);
	sf_putfloat(us,"d4",dw);
	sf_putfloat(us,"o4",ow);
	sf_putstring(us,"label4","Frequency");
	sf_putstring(us,"unit4","Hz");
    } else {
	us = NULL;
    }

    /* output receiver wavefield */
    if (NULL != sf_getstring("ur")) {
	ur = sf_output("ur");
	
	sf_settype(ur,SF_COMPLEX);
	sf_putint(ur,"n3",ns);
	sf_putstring(ur,"label3","Shot");
	sf_putstring(ur,"unit3","");
	sf_putint(ur,"n4",nw);
	sf_putfloat(ur,"d4",dw);
	sf_putfloat(ur,"o4",ow);
	sf_putstring(ur,"label4","Frequency");
	sf_putstring(ur,"unit4","Hz");
    } else {
	ur = NULL;
    }

    /* output time-shift image derivative */
    if (NULL != sf_getstring("timg")) {
	timg = sf_output("timg");

	sf_putint(timg,"n3",2*nh+1);
	sf_putfloat(timg,"d3",d2);
	sf_putfloat(timg,"o3",(float) -nh*d2);

	timage = (float****) sf_alloc(uts,sizeof(float***));
	for (its=0; its < uts; its++) {
	    timage[its] = sf_floatalloc3(n1,n2,2*nh+1);
	}
    } else {
	timg = NULL;
	timage = NULL;
    }

    /* allocate temporary memory */    
    if (load) {
	Ti = NULL; Tj = NULL; Tx = NULL; Tz = NULL; 
	Ap = NULL; Ai = NULL; Map = NULL; Ax = NULL; Az = NULL;
    }
    
    Bx = (double**) sf_alloc(uts,sizeof(double*));
    Bz = (double**) sf_alloc(uts,sizeof(double*));
    Xx = (double**) sf_alloc(uts,sizeof(double*));
    Xz = (double**) sf_alloc(uts,sizeof(double*));

    image = (float****) sf_alloc(uts,sizeof(float***));
    for (its=0; its < uts; its++) {
	image[its] = sf_floatalloc3(n1,n2,2*nh+1);
    }
    
    Numeric = (void**) sf_alloc(uts,sizeof(void*));

    /* LU control */
    umfpack_zl_defaults (Control);
    Control [UMFPACK_IRSTEP] = 0;

    /* loop over frequency */
    for (iw=0; iw < nw; iw++) {
	omega = (double) 2.*SF_PI*(ow+iw*dw);

	/* PML padding */
	pad1 = n1+2*npml;
	pad2 = n2+2*npml;

	n  = fdprep_n (pad1,pad2);
	nz = fdprep_nz(pad1,pad2);

	if (!load) {
	    Ti = (SuiteSparse_long*) sf_alloc(nz,sizeof(SuiteSparse_long));
	    Tj = (SuiteSparse_long*) sf_alloc(nz,sizeof(SuiteSparse_long));
	    Tx = (double*) sf_alloc(nz,sizeof(double));
	    Tz = (double*) sf_alloc(nz,sizeof(double));
	    
	    Ap = (SuiteSparse_long*) sf_alloc(n+1,sizeof(SuiteSparse_long));
	    Ai = (SuiteSparse_long*) sf_alloc(nz,sizeof(SuiteSparse_long));
	    Map = (SuiteSparse_long*) sf_alloc(nz,sizeof(SuiteSparse_long));
	    
	    Ax = (double*) sf_alloc(nz,sizeof(double));
	    Az = (double*) sf_alloc(nz,sizeof(double));
	}
	for (its=0; its < uts; its++) {
	    Bx[its] = (double*) sf_alloc(n,sizeof(double));
	    Bz[its] = (double*) sf_alloc(n,sizeof(double));
	    Xx[its] = (double*) sf_alloc(n,sizeof(double));
	    Xz[its] = (double*) sf_alloc(n,sizeof(double));
	}

	if (verb) {
	    sf_warning("Frequency %d of %d.",iw+1,nw);
	    sf_timer_start(timer);
	}

	/* LU file (append _lu* after velocity file) */
	if (save || load) {
	    sprintf(insert,"_lu%d",iw);
	    inslen = strlen(insert);
	    
	    append = malloc(srclen+inslen+1);

	    memcpy(append,datapath,srclen-5);
	    memcpy(append+srclen-5,insert,inslen);
	    memcpy(append+srclen-5+inslen,datapath+srclen-5,5+1);
	}

	if (!load) {
	    /* assemble matrix */
	    fdprep(omega,
		   n1, n2, d1, d2, vel,
		   npml, pad1, pad2,
		   Ti, Tj, Tx, Tz);	    
	    
	    (void) umfpack_zl_triplet_to_col (n, n, nz, 
					      Ti, Tj, Tx, Tz, 
					      Ap, Ai, Ax, Az, Map);	    

	    /* LU */
	    (void) umfpack_zl_symbolic (n, n, 
					Ap, Ai, Ax, Az, 
					&Symbolic, Control, NULL);
	    
	    (void) umfpack_zl_numeric (Ap, Ai, Ax, Az, 
				       Symbolic, &Numeric[0], 
				       Control, NULL);
	    
	    /* save Numeric */
#ifdef _OPENMP
	    (void) umfpack_zl_save_numeric (Numeric[0], append);
	    
	    for (its=1; its < uts; its++) {
		(void) umfpack_zl_load_numeric (&Numeric[its], append);
	    }
	    
	    if (!save) {
		(void) remove (append);
		(void) remove ("numeric.umf");
	    }
#else
	    if (save) (void) umfpack_zl_save_numeric (Numeric[0], append);
#endif
	} else {
	    /* load Numeric */
	    for (its=0; its < uts; its++) {
		(void) umfpack_zl_load_numeric (&Numeric[its], append);
	    }
	}
	
	if (save || load) free(append);

	/* read source and data */
	sf_complexread(srce[0][0],n1*n2*ns,source);
	sf_complexread(recv[0][0],n1*n2*ns,data);

	/* loop over shots */
#ifdef _OPENMP
#pragma omp parallel for num_threads(uts) private(its,ih,j,i)
#endif
	for (is=0; is < ns; is++) {
#ifdef _OPENMP
	    its = omp_get_thread_num();
#else
	    its = 0;
#endif

	    /* source wavefield */
	    fdpad(npml,pad1,pad2, srce[is],Bx[its],Bz[its]);	    

	    (void) umfpack_zl_solve (UMFPACK_A, 
				     NULL, NULL, NULL, NULL, 
				     Xx[its], Xz[its], Bx[its], Bz[its], 
				     Numeric[its], Control, NULL);
	    
	    fdcut(npml,pad1,pad2, srce[is],Xx[its],Xz[its]);

	    /* receiver wavefield */
	    fdpad(npml,pad1,pad2, recv[is],Bx[its],Bz[its]);	    
	    
	    (void) umfpack_zl_solve (UMFPACK_At, 
				     NULL, NULL, NULL, NULL, 
				     Xx[its], Xz[its], Bx[its], Bz[its], 
				     Numeric[its], Control, NULL);
	    	    
	    fdcut(npml,pad1,pad2, recv[is],Xx[its],Xz[its]);

	    /* imaging condition */
	    for (ih=-nh; ih < nh+1; ih++) {
		for (j=0; j < n2; j++) {
		    for (i=0; i < n1; i++) {
			if (j-abs(ih) >= 0 && j+abs(ih) < n2) {
			    image[its][ih+nh][j][i] += crealf(conjf(srce[is][j-ih][i])*recv[is][j+ih][i]);
			    if (timg != NULL) timage[its][ih+nh][j][i] 
						  += crealf(2.*I*omega*conjf(srce[is][j-ih][i])*recv[is][j+ih][i]);
			}
		    }
		}
	    }
	}

	if (verb) {
	    sf_timer_stop (timer);
	    sf_warning("Finished in %g seconds.",sf_timer_get_diff_time(timer)/1.e3);
	}
	
	if (!load) (void) umfpack_zl_free_symbolic (&Symbolic);
	for (its=0; its < uts; its++) {
	    (void) umfpack_zl_free_numeric (&Numeric[its]);
	}

	if (!load) {
	    free(Ti); free(Tj); free(Tx); free(Tz);
	    free(Ap); free(Ai); free(Map);
	    free(Ax); free(Az);
	}
	for (its=0; its < uts; its++) {
	    free(Bx[its]); free(Bz[its]); free(Xx[its]); free(Xz[its]);
	}

	if (us != NULL) sf_complexwrite(srce[0][0],n1*n2*ns,us);
	if (ur != NULL) sf_complexwrite(recv[0][0],n1*n2*ns,ur);	
    }

#ifdef _OPENMP
#pragma omp parallel for num_threads(uts) private(j,i,its)
    for (ih=-nh; ih < nh+1; ih++) {
	for (j=0; j < n2; j++) {
	    for (i=0; i < n1; i++) {
		for (its=1; its < uts; its++) {
		    image[0][ih+nh][j][i] += image[its][ih+nh][j][i];
		    if (timg != NULL) timage[0][ih+nh][j][i] 
					  += timage[its][ih+nh][j][i];
		}
	    }
	}
    }
#endif
    
    sf_floatwrite(image[0][0][0],n1*n2*(2*nh+1),out);
    if (timg != NULL) sf_floatwrite(timage[0][0][0],n1*n2*(2*nh+1),timg);

    exit(0);
}
int main(int argc, char* argv[])
{

    bool hermite_false, hermite_true;
    int n1, n2, npml, pad1, pad2, ns, nw, nh;
    float d1, d2, **v, ds, os, dw, ow;
    double omega;
    sf_complex ***f, ***srcw, ***recw, ***obs, ***obs_cut;
    sf_file in, out, source, receiver, record;
    int uts, mts;
    char *order;
    int is, i, j, iw, ih;
    float ***image, **recloc;

    sf_init(argc, argv);

    in = sf_input("in");
    out = sf_output("out");

    if (!sf_getint("nh",&nh)) nh=0;

    if (!sf_getint("uts",&uts)) uts=0;
//#ifdef _OPENMP
//    mts = omp_get_max_threads();
//#else
    mts = 1;
//#endif
    uts = (uts < 1)? mts: uts;

    hermite_false=false;
    hermite_true=true;
    /* Hermite operator */

    if (!sf_getint("npml",&npml)) npml=20;
    /* PML width */

    if (NULL == (order = sf_getstring("order"))) order="j";
    /* discretization scheme (default optimal 9-point) */

    fdprep_order(order);

    /* read input dimension */
    if (!sf_histint(in,"n1",&n1)) sf_error("No n1= in input.");
    if (!sf_histint(in,"n2",&n2)) sf_error("No n2= in input.");

    if (!sf_histfloat(in,"d1",&d1)) sf_error("No d1= in input.");
    if (!sf_histfloat(in,"d2",&d2)) sf_error("No d2= in input.");

    v = sf_floatalloc2(n1,n2);
    sf_floatread(v[0],n1*n2,in);

	/* PML padding */
	pad1 = n1+2*npml;
	pad2 = n2+2*npml;

    /* read receiver */
    if (NULL == sf_getstring("receiver")) sf_error("Need receiver=");
    receiver = sf_input("receiver");
    recloc=sf_floatalloc2(n1,n2);
    sf_floatread(recloc[0],n1*n2,receiver);

    /* read source */
    if (NULL == sf_getstring("source")) sf_error("Need source=");
    source = sf_input("source");

    if (!sf_histint(source,"n3",&ns)) sf_error("No ns=.");
    if (!sf_histfloat(source,"d3",&ds)) ds=d2;
    if (!sf_histfloat(source,"o3",&os)) os=0.;

    f = sf_complexalloc3(n1,n2,ns);

    /* read observed data */
    if (NULL == sf_getstring("record")) sf_error("Need record=");
    record = sf_input("record");

    if (!sf_histint(record,"n4",&nw)) sf_error("No nw=.");
    if (!sf_histfloat(record,"d4",&dw)) sf_error("No dw=.");
    if (!sf_histfloat(record,"o4",&ow)) sf_error("No ow=.");

    obs = sf_complexalloc3(n1,n2,ns);
    obs_cut = sf_complexalloc3(n1,n2,ns);

    srcw = sf_complexalloc3(n1,n2,ns);
    recw = sf_complexalloc3(n1,n2,ns);
    image = sf_floatalloc3(n1,n2,2*nh+1);

    /* Loop over frequency */
    for (iw=0; iw<nw; iw++ ) {
        omega=(double) 2.*SF_PI*(ow+iw*dw);

        sf_warning("Calculating frequency %d out of %d for %f HZ.",iw+1,nw,ow+iw*dw);

    	sf_complexread(f[0][0],n1*n2*ns,source);
        sf_complexread(obs[0][0],n1*n2*ns,record);

        /* generate adjoint source for reverse time migration */
        genadjsrc_rtm(obs, obs_cut, recloc, n1, n2, ns);

        /* initialize sparse solver */
        sparse_init(uts, pad1, pad2);

        /* factorize matrix, change according to different frequencies and models */
        sparse_factor(omega,n1,n2,d1,d2,v,npml,pad1,pad2,uts);

        for (is=0; is < ns; is++ ) {
        for (j=0; j < n2; j++ ) {
        for (i=0; i < n1; i++ ) {
            srcw[is][j][i]=f[is][j][i];
            recw[is][j][i]=obs_cut[is][j][i];
        }
        }
        }

        /* sparse solver for source wavefield */
        sparse_solve(npml, pad1, pad2, srcw, hermite_false, ns, uts);

        /* sparse solver for receiver wavefield */
        sparse_solve(npml, pad1, pad2, recw, hermite_true, ns, uts);

        /* imaging condition */
        for (ih=-nh; ih < nh+1; ih++ ) {
        for (j=0; j<n2; j++ ) {
        for (i=0; i< n1; i++ ) {
        for (is=0; is < ns; is++ ) {
        if (j-abs(ih) >= 0 && j+abs(ih) < n2) {
            image[ih+nh][j][i] += crealf(omega*omega*conjf(srcw[is][j-ih][i])*recw[is][j+ih][i]/(v[j][i]*v[j][i]));
        }
        }
        }
        }
        }

        /* free memory */
        sparse_free(uts);
    } /* end frequency */

    sf_putint(out,"n1",n1);
    sf_putint(out,"n2",n2);
    sf_putint(out,"n3",2*nh+1);
    sf_putfloat(out,"d3",d2);
    sf_putfloat(out,"o3", (float) -nh*d2);

    sf_floatwrite(image[0][0],n1*n2*(2*nh+1),out);

    exit(0);
}
Esempio n. 4
0
int main(int argc, char* argv[])
{
 
    bool hermite_false, hermite_true, precond;
    int n1, n2, npml, pad1, pad2, ns, nw, radius;
    float d1, d2, ds, os, dw, ow;
    double omega;
    sf_complex ***f, ***obs;
    sf_file in, out, source, receiver, record, dip;
    char *order;
    int uts, mts, i, j, k, iw, niter, iter;
    float **cur_grad, **pre_grad, **cur_dir, **pre_dir;
    float *cur_grad_input, *cur_grad_smooth, *cur_grad_smooth2;
    float misfit0, misfit1, misfit2, misfitold, beta, alpha;
    float **v, **vnew, **slope, **recloc;

    sf_init(argc, argv);

    in = sf_input("in");
    out = sf_output("out");

    if (!sf_getint("niter",&niter)) niter=0;

    if (!sf_getint("uts",&uts)) uts=0;
//#ifdef _OPENMP
//    mts = omp_get_max_threads();
//#else
    mts = 1;
//#endif
    uts = (uts < 1)? mts: uts;

    hermite_false=false;
    hermite_true=true;
    /* Hermite operator */
    
    if (!sf_getint("npml",&npml)) npml=20;
    /* PML width */

    if (NULL == (order = sf_getstring("order"))) order="j";
    /* discretization scheme (default optimal 9-point) */

    fdprep_order(order);

    /* read input dimension */
    if (!sf_histint(in,"n1",&n1)) sf_error("No n1= in input.");
    if (!sf_histint(in,"n2",&n2)) sf_error("No n2= in input.");

    if (!sf_histfloat(in,"d1",&d1)) sf_error("No d1= in input.");
    if (!sf_histfloat(in,"d2",&d2)) sf_error("No d2= in input.");

    v = sf_floatalloc2(n1,n2);
    vnew = sf_floatalloc2(n1,n2);
    sf_floatread(v[0],n1*n2,in);
    
	/* PML padding */
	pad1 = n1+2*npml;
	pad2 = n2+2*npml;

    /* read receiver */ 
    if (NULL == sf_getstring("receiver")) sf_error("Need receiver="); 
    receiver = sf_input("receiver"); 
    recloc = sf_floatalloc2(n1,n2);
    sf_floatread(recloc[0],n1*n2,receiver);

    /* read source */
    if (NULL == sf_getstring("source")) sf_error("Need source=");
    source = sf_input("source");

    if (!sf_histint(source,"n3",&ns)) sf_error("No ns=.");
    if (!sf_histfloat(source,"d3",&ds)) ds=d2;
    if (!sf_histfloat(source,"o3",&os)) os=0.;

    f = sf_complexalloc3(n1,n2,ns);

    if (NULL == sf_getstring("record")) sf_error("Need record=");
    record = sf_input("record");

    if (!sf_histint(record,"n4",&nw)) sf_error("No nw=.");
    if (!sf_histfloat(record,"d4",&dw)) sf_error("No dw=.");
    if (!sf_histfloat(record,"o4",&ow)) sf_error("No ow=."); 

    /* read in slope information */
    if (!sf_getbool("precond",&precond)) precond=false;
    if (precond) {
        if (NULL == sf_getstring("dip")) sf_error("Need dip=");
        dip = sf_input("dip");
        slope = sf_floatalloc2(n1,n2);
        sf_floatread(slope[0],n1*n2,dip);

        if (!sf_getint("radius",&radius)) sf_error("Need radius=");
        cur_grad_input = sf_floatalloc(n1*n2);
        cur_grad_smooth = sf_floatalloc(n1*n2);
        cur_grad_smooth2 = sf_floatalloc(n1*n2);
    }

    sf_putint(out,"n1",n1);
    sf_putint(out,"n2",n2);
    sf_putint(out,"n3",niter*nw);

    obs = sf_complexalloc3(n1,n2,ns);
    cur_grad = sf_floatalloc2(n1,n2);
    cur_dir = sf_floatalloc2(n1,n2);
    pre_grad = sf_floatalloc2(n1,n2);
    pre_dir = sf_floatalloc2(n1,n2);

    /* Loop over frequency */
    for ( iw = 0; iw < nw; iw ++ ) { 
        omega=(double) 2.*SF_PI*(ow+iw*dw); 
            
        sf_warning("Calculating frequency %d out of %d for %f HZ.\n",iw+1,nw,ow+iw*dw);

        sf_complexread(f[0][0],n1*n2*ns,source);
        sf_complexread(obs[0][0],n1*n2*ns,record);

        misfitold = 100000000.0;
        iter = 0;

	pwsmooth_init(radius, n1, n2, 5, 0.01);
	
        while (iter < niter) { 
			sf_warning("Calculating %d out of %d iteration", iter+1, niter);

            misfit0 = adjfwi_operator(uts, pad1, pad2, omega, n1, n2, d1, d2,
                                      npml, ns, f, obs, hermite_false, hermite_true, recloc,
                                      v, cur_grad);  

            /* Preconditioning or regularization */
            if (precond) {
                sf_warning("Preconditioning gradient.");

                /* change from 2D array to 1D */
                k=0;
                for (j=0;j<n2;j++) {
                for (i=0;i<n1;i++) { 
                    cur_grad_input[k]=cur_grad[j][i];
                    k=k+1;
                }
                }
                pwsmooth_set(slope);
                pwsmooth_lop(true,false,n1*n2,n1*n2,cur_grad_smooth,cur_grad_input);  /* adjoint of pwspray */
                pwsmooth_lop(false,false,n1*n2,n1*n2,cur_grad_smooth,cur_grad_smooth2);  /* forward of pwspray */
                /* reset cur_grad to smooth version */
                k=0;
                for (j=0; j<n2;j++) { 
                for (i=0; i<n1; i++) { 
                    cur_grad[j][i]=cur_grad_smooth2[k];
                    /* cur_grad[j][i]=cur_grad_smooth[k]; */
                    k=k+1;
                }
                }
            }

            /* calculate direction */
            if (iter == 0) { 
                direction_sd(cur_grad,cur_dir,n1,n2);
                beta = 0.0; 
            } else { 
                beta = direction_cg_polak(pre_grad, pre_dir, cur_grad, cur_dir, n1, n2);
            }
            sf_warning("Finish direction calculation.");


            /* test model 1 */
            update_model_fwi(v,vnew,cur_dir,0.1,n1,n2);             
            
            misfit1 = forward_operator(uts, pad1, pad2, omega, n1, n2, d1, d2,
                                       npml, ns, f, obs, hermite_false, recloc,
                                       vnew); 

            sf_warning("Finish test1 calculation.");

            /* test model 2 */
            update_model_fwi(v,vnew,cur_dir,0.2,n1,n2);             

            misfit2 = forward_operator(uts, pad1, pad2, omega, n1, n2, d1, d2,
                                       npml, ns, f, obs, hermite_false, recloc,
                                       vnew); 
            
            sf_warning("Finish test2 calculation.");

            /* update model, quaduartic fit */ 
            alpha = (misfit2-4.0*misfit1+3.0*misfit0)/(20.0*(misfit2-2.0*misfit1+misfit0));
            if ( alpha < 0.0 ) {
                alpha = 0.001;
                sf_warning("alpha is smaller than 0.0");
            }
           
            sf_warning("In iteration %d, alpha = %f, beta = %f, misfit0 = %f.",iter+1,alpha,beta,misfit0);
            sf_warning("Test model alpha=0.1, misfit1 = %f, alpha=0.2, misfit2 = %f.\n",misfit1,misfit2);


            if ( misfit0 > misfitold ) { 
                sf_warning("Terminate at iteration %d.",iter);
                break;
            } else {
                iter++;
                update_model_fwi(v,vnew,cur_dir,alpha,n1,n2);
                for (j=0; j<n2; j++ ) {
                for (i=0; i<n1; i++ ) { 
                    v[j][i]=vnew[j][i];
                    pre_grad[j][i]=cur_grad[j][i];
                    pre_dir[j][i]=cur_dir[j][i];
                }
                }
                misfitold=misfit0;
            }
    
            sf_floatwrite(v[0],n1*n2,out);

        }  /* end iteration */ 
        sf_warning("Ending frequency %d out of %d for %f HZ.\n\n",iw+1,nw,ow+iw*dw);

    } /* end frequency */
       
    exit(0);
}
Esempio n. 5
0
void iwigrad_init(char *order,
                  int npml0,
                  int nn1, int nn2,
                  float dd1, float dd2,
                  int nh0, int ns0,
                  float ow0, float dw0, int nw0,
                  sf_fslice sfile0, sf_fslice rfile0,
                  bool load0, char *datapath0,
                  int uts0)
/*< initialization >*/
{
    fdprep_order(order);

    npml = npml0;

    n1 = nn1;
    n2 = nn2;
    d1 = dd1;
    d2 = dd2;

    nh = nh0;
    ns = ns0;
    ow = ow0;
    dw = dw0;
    nw = nw0;

    sfile = sfile0;
    rfile = rfile0;

    load = load0;
    datapath = datapath0;

    uts = uts0;

    ss[0] = 1;
    ss[1] = n1;
    ss[2] = n1*n2;

    /* LU file */
    if (load) {
        srclen = strlen(datapath);
        insert = sf_charalloc(6);
    } else {
        srclen = 0;
        insert = NULL;
        append = NULL;
    }

    /* allocate temporary space */
    us = sf_complexalloc3(n1,n2,ns);
    ur = sf_complexalloc3(n1,n2,ns);
    as = sf_complexalloc3(n1,n2,ns);
    ar = sf_complexalloc3(n1,n2,ns);

    if (load) {
        Ti = NULL;
        Tj = NULL;
        Tx = NULL;
        Tz = NULL;
        Ap = NULL;
        Ai = NULL;
        Map = NULL;
        Ax = NULL;
        Az = NULL;
    }

    Bx = (double**) sf_alloc(uts,sizeof(double*));
    Bz = (double**) sf_alloc(uts,sizeof(double*));
    Xx = (double**) sf_alloc(uts,sizeof(double*));
    Xz = (double**) sf_alloc(uts,sizeof(double*));

    Numeric = (void**) sf_alloc(uts,sizeof(void*));

    tempx = sf_floatalloc2(n1*n2,uts);
    tempr = sf_floatalloc2(n1*n2*(2*nh+1),uts);

    /* LU control */
    umfpack_zl_defaults (Control);
    Control [UMFPACK_IRSTEP] = 0;
}
Esempio n. 6
0
int main(int argc, char* argv[])
{
 
    bool hermite;
    int n1, n2, npml, pad1, pad2, ns, nw, iw;
    float d1, d2, **v, ds, os, ow, dw;
    double omega;
    sf_complex ***f;
    sf_file in, out, source;
    int uts, mts;
    char *order;
    
    sf_init(argc, argv);

    in = sf_input("in");
    out = sf_output("out");

    if (!sf_getint("uts",&uts)) uts=0;
//#ifdef _OPENMP
//    mts = omp_get_max_threads();
//#else
    mts = 1;
//#endif
    uts = (uts < 1)? mts: uts;

	sf_warning("Using %d out of %d threads!", uts, mts);

    if (!sf_getbool("hermite",&hermite)) hermite=false;
    /* Hermite operator */
    
    if (!sf_getint("npml",&npml)) npml=20;
    /* PML width */

    if (NULL == (order = sf_getstring("order"))) order="j";
    /* discretization scheme (default optimal 9-point) */

    fdprep_order(order);

    /* read input dimension */
    if (!sf_histint(in,"n1",&n1)) sf_error("No n1= in input.");
    if (!sf_histint(in,"n2",&n2)) sf_error("No n2= in input.");

    if (!sf_histfloat(in,"d1",&d1)) sf_error("No d1= in input.");
    if (!sf_histfloat(in,"d2",&d2)) sf_error("No d2= in input.");

    v = sf_floatalloc2(n1,n2);
    sf_floatread(v[0],n1*n2,in);
    
	/* PML padding */
	pad1 = n1+2*npml;
	pad2 = n2+2*npml;

    /* read source */
    if (NULL == sf_getstring("source"))
	sf_error("Need source=");
    source = sf_input("source");

    if (!sf_histint(source,"n3",&ns)) sf_error("No ns=.");
    if (!sf_histfloat(source,"d3",&ds)) ds=d2;
    if (!sf_histfloat(source,"o3",&os)) os=0.;
    if (!sf_histint(source,"n4",&nw)) sf_error("No nw=.");
    if (!sf_histfloat(source,"d4",&dw)) sf_error("No dw=.");
    if (!sf_histfloat(source,"o4",&ow)) sf_error("No ow=.");

    f = sf_complexalloc3(n1,n2,ns);

    /* write out forward simulation */
    sf_settype(out,SF_COMPLEX);

    sf_putint(out,"n3",ns);
    sf_putfloat(out,"d3",ds);
    sf_putfloat(out,"o3",os);
    sf_putstring(out,"label3","Shot");
    sf_putstring(out,"unit3","");
    sf_putint(out,"n4",nw);
    sf_putfloat(out,"d4",dw);
    sf_putfloat(out,"o4",ow);
    sf_putstring(out,"label4","Frequency");
    sf_putstring(out,"unit4","Hz");

    /* Loop over frequency */
    for (iw=0; iw<nw; iw++ ) { 
        omega=(double) 2.*SF_PI*(ow+iw*dw); 

        sf_warning("Calculating frequency %d out of %d for %f HZ.",iw+1,nw,ow+iw*dw);

        /* read in source */
        sf_complexread(f[0][0],n1*n2*ns,source);

        /* initialize sparse solver */ 
        sparse_init(uts, pad1, pad2);

        /* factorize matrix, change according to different frequencies and models */
        sparse_factor(omega, n1, n2, d1, d2, v, npml, pad1, pad2);
    
        /* sparse solver */
        sparse_solve(npml, pad1, pad2, f, hermite, ns, uts);

        /* write out wavefield */
        sf_complexwrite(f[0][0],n1*n2*ns,out);

		sparse_free(uts);
    }

    exit(0);
}