void synthesis(complex *Refl, complex *Fop, float *Top, float *iRN, int nx, int nt, int nxs, int nts, float dt, float *xsyn, int Nsyn, float *xrcv, float *xsrc, float fxs2, float fxs, float dxs, float dxsrc, float dx, int ixa, int ixb, int ntfft, int nw, int nw_low, int nw_high,  int mode, int reci, int nshots, int *ixpossyn, int npossyn, double *tfft, int *first, int verbose)
{
    int     nfreq, size, iox, inx;
    float   scl;
    int     i, j, l, m, iw, ix, k;
    float   *rtrace, idxs;
    complex *sum, *ctrace;
    int     npe;
	static int *ixrcv;
    static double t0, t1, t;

    size  = nxs*nts;
    nfreq = ntfft/2+1;
    /* scale factor 1/N for backward FFT,
     * scale dt for correlation/convolution along time, 
     * scale dx (or dxsrc) for integration over receiver (or shot) coordinates */
    scl   = 1.0*dt/((float)ntfft);

#ifdef _OPENMP
    npe   = omp_get_max_threads();
    /* parallelisation is over number of virtual source positions (Nsyn) */
    if (npe > Nsyn) {
        vmess("Number of OpenMP threads set to %d (was %d)", Nsyn, npe);
        omp_set_num_threads(Nsyn);
    }
#endif

    t0 = wallclock_time();

    /* reset output data to zero */
    memset(&iRN[0], 0, Nsyn*nxs*nts*sizeof(float));

	idxs = 1.0/dxs;
	if (ixrcv == NULL) {
		ixrcv = (int *)malloc(nshots*nx*sizeof(int));
	}
   	for (k=0; k<nshots; k++) {
           for (i = 0; i < nx; i++) {
               ixrcv[k*nx+i] = NINT((xrcv[k*nx+i]-fxs)*idxs);
		}
	}
	ctrace = (complex *)calloc(ntfft,sizeof(complex));
	if (!*first) {
    /* transform muted Ni (Top) to frequency domain, input for next iteration  */
    	for (l = 0; l < Nsyn; l++) {
        	/* set Fop to zero, so new operator can be defined within ixpossyn points */
            //memset(&Fop[l*nxs*nw].r, 0, nxs*nw*2*sizeof(float));
            bzero(&Fop[l*nxs*nw].r, nxs*nw*2*sizeof(float));
            for (i = 0; i < npossyn; i++) {
               	rc1fft(&Top[l*size+i*nts],ctrace,ntfft,-1);
               	ix = ixpossyn[i];
               	for (iw=0; iw<nw; iw++) {
                   	Fop[l*nxs*nw+iw*nxs+ix].r = ctrace[nw_low+iw].r;
                   	Fop[l*nxs*nw+iw*nxs+ix].i = mode*ctrace[nw_low+iw].i;
               	}
            }
		}
	}
	else { /* only for first call to synthesis */
    /* transform G_d to frequency domain, over all nxs traces */
		*first=0;
    	for (l = 0; l < Nsyn; l++) {
        	/* set Fop to zero, so new operator can be defined within all ix points */
            //memset(&Fop[l*nxs*nw].r, 0, nxs*nw*2*sizeof(float));
            bzero(&Fop[l*nxs*nw].r, nxs*nw*2*sizeof(float));
            for (i = 0; i < nxs; i++) {
               	rc1fft(&Top[l*size+i*nts],ctrace,ntfft,-1);
               	for (iw=0; iw<nw; iw++) {
                   	Fop[l*nxs*nw+iw*nxs+i].r = ctrace[nw_low+iw].r;
                   	Fop[l*nxs*nw+iw*nxs+i].i = mode*ctrace[nw_low+iw].i;
               	}
            }
		}
	}
	free(ctrace);
    t1 = wallclock_time();
	*tfft += t1 - t0;

    for (k=0; k<nshots; k++) {

/*        if (verbose>=3) {
            vmess("source position:     %.2f ixpossyn=%d", xsrc[k], ixpossyn[k]);
            vmess("receiver positions:  %.2f <--> %.2f", xrcv[k*nx+0], xrcv[k*nx+nx-1]);
        }
*/
        if ((NINT(xsrc[k]-fxs2) > 0) || (NINT(xrcv[k*nx+nx-1]-fxs2) > 0) ||
            (NINT(xrcv[k*nx+nx-1]-fxs) < 0) || (NINT(xsrc[k]-fxs) < 0) || 
            (NINT(xrcv[k*nx+0]-fxs) < 0) || (NINT(xrcv[k*nx+0]-fxs2) > 0) ) {
            vwarn("source/receiver positions are outside synthesis model");
            vwarn("integration calculation is stopped at gather %d", k);
            vmess("xsrc = %.2f xrcv_1 = %.2f xrvc_N = %.2f", xsrc[k], xrcv[k*nx+0], xrcv[k*nx+nx-1]);
            break;
        }
    

		iox = 0; inx = nx;

/*================ SYNTHESIS ================*/


#pragma omp parallel default(none) \
 shared(iRN, dx, npe, nw, verbose) \
 shared(Refl, Nsyn, reci, xrcv, xsrc, xsyn, fxs, nxs, dxs) \
 shared(nx, ixa, ixb, dxsrc, iox, inx, k, nfreq, nw_low, nw_high) \
 shared(Fop, size, nts, ntfft, scl, ixrcv, stderr) \
 private(l, ix, j, m, i, sum, rtrace)
    { /* start of parallel region */
    sum   = (complex *)malloc(nfreq*sizeof(complex));
    rtrace = (float *)calloc(ntfft,sizeof(float));

#pragma omp for schedule(guided,1)
    for (l = 0; l < Nsyn; l++) {

            ix = k; 

			/* multiply R with Fop and sum over nx */
			memset(&sum[0].r,0,nfreq*2*sizeof(float));
            //for (j = 0; j < nfreq; j++) sum[j].r = sum[j].i = 0.0;
                for (j = nw_low, m = 0; j <= nw_high; j++, m++) {
            for (i = iox; i < inx; i++) {
                    sum[j].r += Refl[k*nw*nx+m*nx+i].r*Fop[l*nw*nxs+m*nxs+ixrcv[k*nx+i]].r -
                                Refl[k*nw*nx+m*nx+i].i*Fop[l*nw*nxs+m*nxs+ixrcv[k*nx+i]].i;
                    sum[j].i += Refl[k*nw*nx+m*nx+i].i*Fop[l*nw*nxs+m*nxs+ixrcv[k*nx+i]].r +
                                Refl[k*nw*nx+m*nx+i].r*Fop[l*nw*nxs+m*nxs+ixrcv[k*nx+i]].i;
                }
            }

			/* transfrom result back to time domain */
            cr1fft(sum, rtrace, ntfft, 1);

            /* dx = receiver distance */
            for (j = 0; j < nts; j++) 
                iRN[l*size+ix*nts+j] += rtrace[j]*scl*dx;

    } /* end of parallel Nsyn loop */

    free(sum);
    free(rtrace);

#pragma omp single 
{ 
#ifdef _OPENMP
    npe   = omp_get_num_threads();
#endif
}
    } /* end of parallel region */

	if (verbose>3) vmess("*** Shot gather %d processed ***", k);

    } /* end of nshots (k) loop */

    t = wallclock_time() - t0;
    if (verbose) {
        vmess("OMP: parallel region = %f seconds (%d threads)", t, npe);
    }

    return;
}
Example #2
0
int main(int argc, char **argv)
{
	modPar mod;
	recPar rec;
	srcPar src;
	shotPar shot;
	rayPar ray;
    float *velocity, *slowness, *smooth, *trueslow, **inter;
	double t0, t1, t2, tinit, tray, tio;
	size_t size;
	int nw, n1, ix, iz, ir, ixshot, izshot;
	int nt, ntfft, nfreq, ig;
	int irec, sbox, ipos, nrx, nrz, nr;
    fcoord coordsx, coordgx, Time;
    icoord grid, isrc; 
    float Jr, *ampl, *time, *ttime, *ttime_p, cp_average, *wavelet, dw, dt;
	float dxrcv, dzrcv, rdelay, tr, dt_tmp;
    segy hdr;
    char filetime[1024], fileamp[1024], *method, *file_rcvtime, *file_src;
    size_t  nwrite, nread;
	int verbose;
    complex *cmute, *cwav;
    FILE *fpt, *fpa, *fpwav, *fprcv;

	t0= wallclock_time();
	initargs(argc,argv);
	requestdoc(0);

	if(!getparint("verbose",&verbose)) verbose=0;
    if(!getparint("sbox", &sbox)) sbox = 1;
    if(!getparstring("method", &method)) method="jesper";
	if (!getparfloat("rec_delay",&rdelay)) rdelay=0.0;

	getParameters(&mod, &rec, &src, &shot, &ray, verbose);

    /* read file_src if file_rcvtime is defined */

    if (!getparstring("file_rcvtime",&file_rcvtime)) file_rcvtime=NULL;

	if (file_rcvtime != NULL) {
    	if (!getparstring("file_src",&file_src)) file_src=NULL;
		if (!getparfloat("dt",&dt)) dt=0.004;
		if (file_src != NULL ) {
        	fpwav = fopen( file_src, "r" );
        	assert( fpwav != NULL);
        	nread = fread( &hdr, 1, TRCBYTES, fpwav );
        	assert(nread == TRCBYTES);
			ntfft = optncr(MAX(hdr.ns, rec.nt));
 			wavelet = (float *)calloc(ntfft,sizeof(float));
			/* read first trace */
        	nread = fread(wavelet, sizeof(float), hdr.ns, fpwav);
        	assert (nread == hdr.ns);
        	fclose(fpwav);
		}
		else {
			ntfft = optncr(rec.nt);
 			wavelet = (float *)calloc(ntfft,sizeof(float));
			wavelet[0] = 1.0;
		}
    	nfreq = ntfft/2+1;
    	cwav    = (complex *)calloc(nfreq,sizeof(complex));
    	cmute   = (complex *)calloc(nfreq,sizeof(complex));
        rc1fft(wavelet,cwav,ntfft,-1);
    	dw      = 2*M_PI/(ntfft*dt);
	}

	/* allocate arrays for model parameters: the different schemes use different arrays */

	n1 = mod.nz;
    if(!strcmp(method,"fd")) nw = 0;
    else nw = ray.smoothwindow;

	velocity = (float *)calloc(mod.nx*mod.nz,sizeof(float));
	slowness = (float *)calloc((mod.nx+2*nw)*(mod.nz+2*nw),sizeof(float));
    trueslow = (float *)calloc(mod.nx*mod.nz,sizeof(float));

    if(!strcmp(method,"fd")) {
		ttime = (float *)calloc(mod.nx*mod.nz,sizeof(float));
	}

	/* read velocity and density files */
	readModel(mod, velocity, slowness, nw);

	/* allocate arrays for wavefield and receiver arrays */

	size = shot.n*rec.n;
    time = (float *)calloc(size,sizeof(float));
    ampl = (float *)calloc(size,sizeof(float));

	/* Sinking source and receiver arrays: 
	   If P-velocity==0 the source and receiver 
	   postions are placed deeper until the P-velocity changes. 
	   Setting the option rec.sinkvel only sinks the receiver position 
       (not the source) and uses the velocity 
	   of the first receiver to sink through to the next layer. */

/* sink receivers to value different than sinkvel */
	for (ir=0; ir<rec.n; ir++) {
		iz = rec.z[ir];
		ix = rec.x[ir];
		while(velocity[(ix)*n1+iz] == rec.sinkvel) iz++;
		rec.z[ir]=iz+rec.sinkdepth;
		rec.zr[ir]=rec.zr[ir]+(rec.z[ir]-iz)*mod.dz;
//		rec.zr[ir]=rec.z[ir]*mod.dz;
		if (verbose>3) vmess("receiver position %d at grid[ix=%d, iz=%d] = (x=%f z=%f)", ir, ix, rec.z[ir], rec.xr[ir]+mod.x0, rec.zr[ir]+mod.z0);
	}
		vmess("   - method for ray-tracing       = %s", method);
/*
*/

/* sink sources to value different than zero */
	for (izshot=0; izshot<shot.nz; izshot++) {
		for (ixshot=0; ixshot<shot.nx; ixshot++) {
			iz = shot.z[izshot];
			ix = shot.x[ixshot];
			while(velocity[(ix)*n1+iz] == 0.0) iz++;
			shot.z[izshot]=iz+src.sinkdepth; 
		}
	}

	if (verbose>3) writeSrcRecPos(&mod, &rec, &src, &shot);

    /* smooth slowness grid */
    grid.x = mod.nx;
    grid.z = mod.nz;
    grid.y = 1;
    if ( nw != 0 ) { /* smooth slowness */ 
        smooth = (float *)calloc(grid.x*grid.z,sizeof(float));
        applyMovingAverageFilter(slowness, grid, nw, 2, smooth);
        memcpy(slowness,smooth,grid.x*grid.z*sizeof(float));
        free(smooth);
	}

    /* prepare output file and headers */
    strcpy(filetime, rec.file_rcv);
    name_ext(filetime, "_time");
    fpt = fopen(filetime, "w");
    assert(fpt != NULL);

	if (ray.geomspread) {
        strcpy(fileamp, rec.file_rcv);
        name_ext(fileamp, "_amp");
        fpa = fopen(fileamp, "w");
        assert(fpa != NULL);
	}
	if (file_rcvtime != NULL) {
        fprcv = fopen(file_rcvtime, "w");
        assert(fprcv != NULL);
	}

    memset(&hdr,0,sizeof(hdr));
    hdr.scalco = -1000;
    hdr.scalel = -1000;
    hdr.trid   = 0;

	t1=wallclock_time();
	tinit = t1-t0;
    tray=0.0;
    tio=0.0;

	/* Outer loop over number of shots */
	for (izshot=0; izshot<shot.nz; izshot++) {
		for (ixshot=0; ixshot<shot.nx; ixshot++) {

	        t2=wallclock_time();
        	if (verbose) {
            	vmess("Modeling source %d at gridpoints ix=%d iz=%d", (izshot*shot.n)+ixshot, shot.x[ixshot], shot.z[izshot]);
            	vmess(" which are actual positions x=%.2f z=%.2f", mod.x0+mod.dx*shot.x[ixshot], mod.z0+mod.dz*shot.z[izshot]);
            	vmess("Receivers at gridpoint x-range ix=%d - %d", rec.x[0], rec.x[rec.n-1]);
            	vmess(" which are actual positions x=%.2f - %.2f", mod.x0+rec.xr[0], mod.x0+rec.xr[rec.n-1]);
            	vmess("Receivers at gridpoint z-range iz=%d - %d", rec.z[0], rec.z[rec.n-1]);
            	vmess(" which are actual positions z=%.2f - %.2f", mod.z0+rec.zr[0], mod.z0+rec.zr[rec.n-1]);
        	}

        	coordsx.x = shot.x[ixshot]*mod.dx;
        	coordsx.z = shot.z[izshot]*mod.dz;
        	coordsx.y = 0;

	        t1=wallclock_time();
            tio += t1-t2;

            if (!strcmp(method,"jesper")) {
#pragma omp parallel for default(shared) \
private (coordgx,irec,Time,Jr) 
        		for (irec=0; irec<rec.n; irec++) {
            		coordgx.x=rec.xr[irec];
            		coordgx.z=rec.zr[irec];
            		coordgx.y = 0;
		
            		getWaveParameter(slowness, grid, mod.dx, coordsx, coordgx, ray, &Time, &Jr);
	
            		time[((izshot*shot.nx)+ixshot)*rec.n + irec] = Time.x + Time.y + 0.5*Time.z;
            		ampl[((izshot*shot.nx)+ixshot)*rec.n + irec] = Jr;
            		if (verbose>4) vmess("JS: shot=%f,%f receiver at %f,%f T0=%f T1=%f T2=%f Jr=%f",coordsx.x, coordsx.z, coordgx.x, coordgx.z, Time.x, Time.y, Time.z, Jr); 
        		}
			}
            else if(!strcmp(method,"fd")) {
	            int mzrcv;

                isrc.x = shot.x[ixshot];
                isrc.y = 0;
                isrc.z = shot.z[izshot];

                mzrcv = 0;
                for (irec = 0; irec < rec.n; irec++) mzrcv = MAX(rec.z[irec], mzrcv);

                vidale(ttime,slowness,&isrc,grid,mod.dx,sbox, mzrcv);
        		for (irec=0; irec<rec.n; irec++) {
            		coordgx.x=mod.x0+rec.xr[irec];
            		coordgx.z=mod.z0+rec.zr[irec];
            		coordgx.y = 0;
					ipos = ((izshot*shot.nx)+ixshot)*rec.n + irec;
	
            		time[ipos] = ttime[rec.z[irec]*mod.nx+rec.x[irec]];
					/* compute average velocity between source and receiver */
 					nrx = (rec.x[irec]-isrc.x);
 					nrz = (rec.z[irec]-isrc.z);
 					nr = abs(nrx) + abs(nrz);
					cp_average = 0.0;
        			for (ir=0; ir<nr; ir++) {
						ix = isrc.x + floor((ir*nrx)/nr);
						iz = isrc.z + floor((ir*nrz)/nr);
						//fprintf(stderr,"ir=%d ix=%d iz=%d velocity=%f\n", ir, ix, iz, velocity[ix*mod.nz+iz]);
						cp_average += velocity[ix*mod.nz+iz];
					}
					cp_average = cp_average/((float)nr);
            		ampl[ipos] = sqrt(time[ipos]*cp_average);
            		if (verbose>4) vmess("FD: shot=%f,%f receiver at %f(%d),%f(%d) T=%e V=%f Ampl=%f",coordsx.x, coordsx.z, coordgx.x, rec.x[irec], coordgx.z, rec.z[irec], time[ipos], cp_average, ampl[ipos]); 
        		}
            }
	        t2=wallclock_time();
            tray += t2-t1;

        	hdr.sx     = 1000*(mod.x0+mod.dx*shot.x[ixshot]);
        	hdr.sdepth = 1000*(mod.z0+mod.dz*shot.z[izshot]);
        	hdr.selev  = (int)(-1000.0*(mod.z0+mod.dz*shot.z[izshot]));
        	hdr.fldr   = ((izshot*shot.nx)+ixshot)+1;
        	hdr.tracl  = ((izshot*shot.nx)+ixshot)+1;
        	hdr.tracf  = ((izshot*shot.nx)+ixshot)+1;
        	hdr.ntr    = shot.n;
    		hdr.dt     = (unsigned short)1;
    		hdr.trwf   = shot.n;
    		hdr.ns     = rec.n;
        	//hdr.d1     = (rec.x[1]-rec.x[0])*mod.dx; // discrete
        	hdr.d1     = (rec.xr[1]-rec.xr[0]);
        	hdr.f1     = mod.x0+rec.x[0]*mod.dx;
        	hdr.d2     = (shot.x[MIN(shot.n-1,1)]-shot.x[0])*mod.dx;
        	hdr.f2     = mod.x0+shot.x[0]*mod.dx;
			dt_tmp = (fabs(hdr.d1*((float)hdr.scalco)));
			hdr.dt	   = (unsigned short)dt_tmp;
    
        	nwrite = fwrite( &hdr, 1, TRCBYTES, fpt);
        	assert(nwrite == TRCBYTES);
        	nwrite = fwrite( &time[((izshot*shot.nx)+ixshot)*rec.n], sizeof(float), rec.n, fpt);
        	assert(nwrite == rec.n);
	    	fflush(fpt);
	    	if (ray.geomspread) {
            	nwrite = fwrite( &hdr, 1, TRCBYTES, fpa);
            	assert(nwrite == TRCBYTES);
            	nwrite = fwrite( &ampl[((izshot*shot.nx)+ixshot)*rec.n], sizeof(float), rec.n, fpa);
            	assert(nwrite == rec.n);
	        	fflush(fpa);
        	}
			if (file_rcvtime != NULL) {
    			hdr.ns     = rec.nt;
    			hdr.trwf   = rec.n;
    			hdr.ntr    = ((izshot*shot.nx)+ixshot+1)*rec.n;
    			hdr.dt     = dt*1000000;
    			hdr.d1     = dt;
        		hdr.f1     = 0.0;
        		hdr.d2     = (rec.xr[1]-rec.xr[0]);
    			hdr.f2     = mod.x0+rec.x[0]*mod.dx;

        		for (irec=0; irec<rec.n; irec++) {
					ipos = ((izshot*shot.nx)+ixshot)*rec.n + irec;
        			hdr.tracf  = irec+1;
        			hdr.tracl  = ((izshot*shot.nx)+ixshot*shot.nz)+irec+1;
        			hdr.gx     = 1000*(mod.x0+rec.xr[irec]);
        			hdr.offset = (rec.xr[irec]-shot.x[ixshot]*mod.dx);
        			hdr.gelev  = (int)(-1000*(mod.z0+rec.zr[irec]));

					tr = time[ipos]+rdelay;
                    for (ig=0; ig<nfreq; ig++) {
                        cmute[ig].r = (cwav[ig].r*cos(ig*dw*tr-M_PI/4.0)-cwav[ig].i*sin(ig*dw*tr-M_PI/4.0))/(ntfft*ampl[ipos]);
                        cmute[ig].i = (cwav[ig].i*cos(ig*dw*tr-M_PI/4.0)+cwav[ig].r*sin(ig*dw*tr-M_PI/4.0))/(ntfft*ampl[ipos]);
                    }
                	cr1fft(cmute,wavelet,ntfft,-1);
        			nwrite = fwrite( &hdr, 1, TRCBYTES, fprcv);
        			nwrite = fwrite( wavelet, sizeof(float), rec.nt, fprcv );
				}
			}
	        t1=wallclock_time();
            tio += t1-t2;
    	} /* end of ixshot loop */
	} /* end of loop over number of shots */
	fclose(fpt);
	if (file_rcvtime != NULL) fclose(fprcv);
	if (ray.geomspread) fclose(fpa);

	t1= wallclock_time();
	if (verbose) {
		vmess("*******************************************");
		vmess("************* runtime info ****************");
		vmess("*******************************************");
		vmess("Total compute time ray-tracing    = %.2f s.", t1-t0);
		vmess("   - intializing arrays and model = %.3f", tinit);
		vmess("   - ray tracing                  = %.3f", tray);
		vmess("   - writing data to file         = %.3f", tio);
	}

	/* free arrays */

	initargs(argc,argv); /* this will free the arg arrays declared */
	free(velocity);
	free(slowness);
	
	return 0;
}
Example #3
0
void synthesis(complex *Refl, complex *Fop, float *Top, float *iRN, int nx, int nt, int nxs, int nts, float dt, float *xsyn, int
Nfoc, float *xrcv, float *xsrc, int *xnx, float fxse, float fxsb, float dxs, float dxsrc, float dx, int ntfft, int
nw, int nw_low, int nw_high,  int mode, int reci, int nshots, int *ixpos, int npos, double *tfft, int *isxcount, int
*reci_xsrc,  int *reci_xrcv, float *ixmask, int verbose)
{
    int     nfreq, size, inx;
    float   scl;
    int     i, j, l, m, iw, ix, k, ixsrc, il, ik;
    float   *rtrace, idxs;
    complex *sum, *ctrace;
    int     npe;
    static int first=1, *ixrcv;
    static double t0, t1, t;

    size  = nxs*nts;
    nfreq = ntfft/2+1;
    /* scale factor 1/N for backward FFT,
     * scale dt for correlation/convolution along time, 
     * scale dx (or dxsrc) for integration over receiver (or shot) coordinates */
    scl   = 1.0*dt/((float)ntfft);

#ifdef _OPENMP
    npe   = omp_get_max_threads();
    /* parallelisation is over number of shot positions (nshots) */
    if (npe > nshots) {
        vmess("Number of OpenMP threads set to %d (was %d)", nshots, npe);
        omp_set_num_threads(nshots);
    }
#endif

    t0 = wallclock_time();

    /* reset output data to zero */
    memset(&iRN[0], 0, Nfoc*nxs*nts*sizeof(float));
    ctrace = (complex *)calloc(ntfft,sizeof(complex));

/* this first check is done to support an acquisition geometry that has more receiver than source
 * postions. In the first iteration the int R(x_r,x_s) Fop(x_r) d x_r results in a grid on x_s. 
 * so for the next interations onlt x_s traces have to be computed on Fop */
    if (!first) {
    /* transform muted Ni (Top) to frequency domain, input for next iteration  */
        for (l = 0; l < Nfoc; l++) {
            /* set Fop to zero, so new operator can be defined within ixpos points */
            memset(&Fop[l*nxs*nw].r, 0, nxs*nw*2*sizeof(float));
            for (i = 0; i < npos; i++) {
                rc1fft(&Top[l*size+i*nts],ctrace,ntfft,-1);
                ix = ixpos[i];
                for (iw=0; iw<nw; iw++) {
                    Fop[l*nxs*nw+iw*nxs+ix].r = ctrace[nw_low+iw].r;
                    Fop[l*nxs*nw+iw*nxs+ix].i = mode*ctrace[nw_low+iw].i;
                }
            }
        }
    }
    else { /* only for first call to synthesis using all nxs traces in G_d */
    /* transform G_d to frequency domain, over all nxs traces */
        first=0;
        for (l = 0; l < Nfoc; l++) {
            /* set Fop to zero, so new operator can be defined within all ix points */
            memset(&Fop[l*nxs*nw].r, 0, nxs*nw*2*sizeof(float));
            for (i = 0; i < nxs; i++) {
                rc1fft(&Top[l*size+i*nts],ctrace,ntfft,-1);
                for (iw=0; iw<nw; iw++) {
                    Fop[l*nxs*nw+iw*nxs+i].r = ctrace[nw_low+iw].r;
                    Fop[l*nxs*nw+iw*nxs+i].i = mode*ctrace[nw_low+iw].i;
                }
            }
        }
        idxs = 1.0/dxs;
        ixrcv = (int *)malloc(nshots*nx*sizeof(int));
        for (k=0; k<nshots; k++) {
            for (i = 0; i < nx; i++) {
                ixrcv[k*nx+i] = NINT((xrcv[k*nx+i]-fxsb)*idxs);
            }
        }
    }
    free(ctrace);
    t1 = wallclock_time();
    *tfft += t1 - t0;

    if (reci == 0 || reci == 1) {

/*================ SYNTHESIS ================*/

#pragma omp parallel default(none) \
 shared(iRN, dx, npe, nw, verbose, nshots, xnx) \
 shared(Refl, Nfoc, reci, xsrc, xsyn, fxsb, fxse, nxs, dxs) \
 shared(nx, dxsrc, nfreq, nw_low, nw_high) \
 shared(Fop, size, nts, ntfft, scl, ixrcv) \
 private(l, ix, j, m, i, sum, rtrace, k, ixsrc, inx)
{ /* start of parallel region */
        sum   = (complex *)malloc(nfreq*sizeof(complex));
        rtrace = (float *)calloc(ntfft,sizeof(float));

/* Loop over total number of shots */
#pragma omp for schedule(guided,1)
        for (k=0; k<nshots; k++) {
            if ((xsrc[k] < 0.999*fxsb) || (xsrc[k] > 1.001*fxse)) continue;
            ixsrc = NINT((xsrc[k] - fxsb)/dxs);
            inx = xnx[k]; /* number of traces per shot */

            for (l = 0; l < Nfoc; l++) {
                /* compute integral over receiver positions */
                /* multiply R with Fop and sum over nx */
                memset(&sum[0].r,0,nfreq*2*sizeof(float));
                for (j = nw_low, m = 0; j <= nw_high; j++, m++) {
                    for (i = 0; i < inx; i++) {
                        ix = ixrcv[k*nx+i];
                        sum[j].r += Refl[k*nw*nx+m*nx+i].r*Fop[l*nw*nxs+m*nxs+ix].r -
                                    Refl[k*nw*nx+m*nx+i].i*Fop[l*nw*nxs+m*nxs+ix].i;
                        sum[j].i += Refl[k*nw*nx+m*nx+i].i*Fop[l*nw*nxs+m*nxs+ix].r +
                                    Refl[k*nw*nx+m*nx+i].r*Fop[l*nw*nxs+m*nxs+ix].i;
                    }
                }

                /* transfrom result back to time domain */
                cr1fft(sum, rtrace, ntfft, 1);

                /* place result at source position ixsrc; dx = receiver distance */
                for (j = 0; j < nts; j++) 
                    iRN[l*size+ixsrc*nts+j] += rtrace[j]*scl*dx;
            
            } /* end of parallel Nfoc loop */

            if (verbose>4) vmess("*** Shot gather %d processed ***", k);

        } /* end of nshots (k) loop */
        free(sum);
        free(rtrace);

} /* end of parallel region */


    }     /* end of if reci */

/* if reciprocal traces are enabled start a new loop over reciprocal shot positions */
    if (reci != 0) {

#pragma omp parallel default(none) \
 shared(iRN, dx, nw, verbose) \
 shared(Refl, Nfoc, reci, xsrc, xsyn, fxsb, fxse, nxs, dxs) \
 shared(nx, dxsrc, nfreq, nw_low, nw_high) \
 shared(reci_xrcv, reci_xsrc, ixmask, isxcount) \
 shared(Fop, size, nts, ntfft, scl, ixrcv) \
 private(l, ix, j, m, i, k, sum, rtrace, ik, il, ixsrc, inx)
{ /* start of parallel region */
        sum   = (complex *)malloc(nfreq*sizeof(complex));
        rtrace = (float *)calloc(ntfft,sizeof(float));

#pragma omp for schedule(guided,1)
        for (k=0; k<nxs; k++) {
            if (isxcount[k] == 0) continue;
            ixsrc = k;
            inx = isxcount[ixsrc]; /* number of traces per reciprocal shot */

            for (l = 0; l < Nfoc; l++) {
                /* compute integral over (reciprocal) source positions */
                /* multiply R with Fop and sum over nx */
                memset(&sum[0].r,0,nfreq*2*sizeof(float));
                for (j = nw_low, m = 0; j <= nw_high; j++, m++) {
                    for (i = 0; i < inx; i++) {
                        il = reci_xrcv[ixsrc*nxs+i];
                        ik = reci_xsrc[ixsrc*nxs+i];
                        ix = NINT((xsrc[il] - fxsb)/dxs);
                        sum[j].r += Refl[il*nw*nx+m*nx+ik].r*Fop[l*nw*nxs+m*nxs+ix].r -
                                    Refl[il*nw*nx+m*nx+ik].i*Fop[l*nw*nxs+m*nxs+ix].i;
                        sum[j].i += Refl[il*nw*nx+m*nx+ik].i*Fop[l*nw*nxs+m*nxs+ix].r +
                                    Refl[il*nw*nx+m*nx+ik].r*Fop[l*nw*nxs+m*nxs+ix].i;
                    }
                }

                /* transfrom result back to time domain */
                cr1fft(sum, rtrace, ntfft, 1);

                /* place result at source position ixsrc; dxsrc = shot distance */
                for (j = 0; j < nts; j++) 
                    iRN[l*size+ixsrc*nts+j] = ixmask[ixsrc]*(iRN[l*size+ixsrc*nts+j]+rtrace[j]*scl*dxsrc);
                
            } /* end of Nfoc loop */

        } /* end of parallel reciprocal shots (k) loop */
        free(sum);
        free(rtrace);

 } /* end of parallel region */

    } /* end of if reci */

    t = wallclock_time() - t0;
    if (verbose>2) {
        vmess("OMP: parallel region = %f seconds (%d threads)", t, npe);
    }

    return;
}
Example #4
0
int getWaveletInfo(char *file_src, int *n1, int *n2, float *d1, float *d2, float *f1, float *f2, float *fmax, int *nxm, int verbose)
{
    FILE    *fp;
    size_t  nread, data_sz;
	off_t bytes, ret, trace_sz, ntraces;
    int sx_shot, gx_start, one_shot;
    int optn, nfreq, i, iwmax;
    float *trace;
	float ampl, amplmax, tampl, tamplmax; 
    complex *ctrace;
    segy hdr;
    
    if (file_src == NULL) return 0; /* Input pipe can not be handled */
    else fp = fopen( file_src, "r" );
    assert( fp != NULL);
    nread = fread( &hdr, 1, TRCBYTES, fp );
    assert(nread == TRCBYTES);
    ret = fseeko( fp, 0, SEEK_END );
	if (ret<0) perror("fseeko");
    bytes = ftello( fp );

    *n1 = hdr.ns;
    if (hdr.trid == 1 || hdr.dt != 0) {
        *d1 = ((float) hdr.dt)*1.e-6;
        *f1 = ((float) hdr.delrt)/1000.;
		if (*d1 == 0.0) *d1 = hdr.d1;
    }
    else {
        *d1 = hdr.d1;
        *f1 = hdr.f1;
    }
    *f2 = hdr.f2;

    data_sz = sizeof(float)*(*n1);
    trace_sz = sizeof(float)*(*n1)+TRCBYTES;
    ntraces  = (int) (bytes/trace_sz);
	*n2 = ntraces;

    /* check to find out number of traces in shot gather */

	optn = optncr(hdr.ns);
	nfreq = optn/2 + 1;
	ctrace = (complex *)malloc(nfreq*sizeof(complex));
    one_shot = 1;
    sx_shot  = hdr.sx;
    gx_start = hdr.gx;
    trace = (float *)malloc(optn*sizeof(float));
    fseeko( fp, TRCBYTES, SEEK_SET );

    while (one_shot) {
		memset(trace,0,optn*sizeof(float));
        nread = fread( trace, sizeof(float), hdr.ns, fp );
        assert (nread == hdr.ns);
		tamplmax = 0.0;
		for (i=0;i<hdr.ns;i++) {
			tampl = fabsf(trace[i]);
			if (tampl > tamplmax) tamplmax = tampl;
		}
		if (trace[0]*1e-3 > tamplmax) {
			fprintf(stderr,"WARNING: file_src has a large amplitude %f at t=0\n", trace[0]);
			fprintf(stderr,"This will introduce high frequencies and can cause dispersion.\n");
		}

		/* estimate maximum frequency assuming amplitude spectrum is smooth */
		rc1fft(trace,ctrace,optn,1);

		/* find maximum amplitude */
		amplmax = 0.0;
		for (i=0;i<nfreq;i++) {
			ampl = sqrt(ctrace[i].r*ctrace[i].r+ctrace[i].i*ctrace[i].i);
			if (ampl > amplmax) {
				amplmax = ampl;
				iwmax = i;
			}
		}
		/* from the maximum amplitude position look for the largest frequency
         * which has an amplitude 400 times weaker than the maximum amplitude */
		for (i=iwmax;i<nfreq;i++) {
			ampl = sqrt(ctrace[i].r*ctrace[i].r+ctrace[i].i*ctrace[i].i);
			if (400*ampl < amplmax) {
				*fmax = (i-1)*(1.0/(optn*(*d1)));
				break;
			}
		}

        nread = fread( &hdr, 1, TRCBYTES, fp );
        if (nread==0) break;
    }
	*nxm = (int)ntraces;

	if (verbose>2) {
		vmess("For file %s", file_src);
		vmess("nt=%d nx=%d", *n1, *n2);
		vmess("dt=%f dx=%f", *d1, *d2);
		vmess("fmax=%f", *fmax);
		vmess("tstart=%f", *f1);
	}

    fclose(fp);
    free(trace);
    free(ctrace);

    return 0;
}
Example #5
0
int read_FFT_DataFile(FILE *fp, complex *data, Area vel_area, int nfft, int nw, int nw_low, int *tr_read_in, int *tr_shot, int *ixmin, int *ixmax, int *iymin, int *iymax, int *sx, int *sy, int conjg, int verbose)
{
	int ix, iy, iw, i, nt, nx, ny, sxy, nfreq, pos, sign;
	float xvmin, yvmin;
	int err=0;
	int one_shot, traces_shot, fldr_s, traces_read_in;
	size_t  nread;
	float dx, dy, *trace, scl;
	complex *ctrace;
	segy *hdr;

	nx  = vel_area.nx;
	ny  = vel_area.ny;
	sxy = vel_area.sxy;
	dx  = vel_area.dx;
	dy  = vel_area.dy;
	xvmin = vel_area.xmin;
	yvmin = vel_area.ymin;

	sign   = -1;
	nfreq  = nfft/2 + 1;
	hdr    = (segy *)malloc(TRCBYTES);
    trace  =   (float *)calloc(nfft,sizeof(float));
    ctrace = (complex *)malloc(nfreq*sizeof(complex));

	one_shot=1;
	traces_shot = 0;
	traces_read_in = *tr_read_in;
	while (one_shot) {
		nread = fread( hdr, 1, TRCBYTES, fp );
		if (nread==0) {
			err = -1;
			break;
		}
		nt = hdr[0].ns;
		if (traces_shot==0) fldr_s = hdr[0].fldr;
		if ((fldr_s != hdr[0].fldr)) {
			fseek(fp, -TRCBYTES, SEEK_CUR);
			break;
		}

		if (traces_shot==0) {
			if (hdr[0].scalco < 0) scl = 1.0/fabs(hdr[0].scalco);
			else if (hdr[0].scalco == 0) scl = 1.0;
			else scl = hdr[0].scalco;
			*sx = hdr[0].sx;
			*sy = hdr[0].sy;
		}

		nread = fread( trace, sizeof(float), nt, fp );
		assert (nread == nt);
		traces_shot++;
		traces_read_in++;
	
		if (nfft > nt) memset( &trace[nt], 0, sizeof(float)*(nfft-nt) ); 
		rc1fft(trace,ctrace,nfft,sign);

		ix = (hdr[0].gx*scl-xvmin)/dx;
		iy = (hdr[0].gy*scl-yvmin)/dy;
		if (ix >=0 && ix<nx && iy>=0 && iy<ny) {
			for (iw=0; iw<nw; iw++) {
				data[iw*sxy+iy*nx+ix].r = ctrace[nw_low+iw].r;
				data[iw*sxy+iy*nx+ix].i = conjg*ctrace[nw_low+iw].i;
			}
			*ixmin = MIN(ix,*ixmin);
			*iymin = MIN(iy,*iymin);
			*ixmax = MAX(ix,*ixmax);
			*iymax = MAX(iy,*iymax);
		}
		else {
			fprintf(stderr,"*** trace at %.2f (%d), %.2f (%d) outside model\n",
				hdr[0].gx*scl, ix, hdr[0].gy*scl, iy);
		}

		if (verbose>2) {
			fprintf(stderr," trace %d: gx = %.2f(%d) gy = %.2f(%d) \n", 
				traces_shot, hdr[0].gx*scl, ix, hdr[0].gy*scl, iy);
		}

	} /* end of receiver gather */
	*tr_read_in = traces_read_in;
	*tr_shot    = traces_shot;

	free(hdr);
	free(trace);
	free(ctrace);

	return err;
}
Example #6
0
void kwZoMigr(float *data, int nx, int nt, float dt, float *velmod, int nxm, int nzm, int ixa, int ixb, float fmin, float fmax, float *xrcv, int izrcv, float ox, float dxm, float dz, int ntap, int conjg, int ndepth, float *image, int verbose, float *exrcv, int ndepthex)
{
    int     iomin, iomax, iom, ix, d, i, j;
    int     nfreq, optn, nkx, sign, endkx;
    int     ixrcv, ixmin, ixmax, ixo, ixn, ikx;
    float   k, k2, kz2, kx, kx2;
    float   dom, om, c, dkx, df, sr;
    float   *taper, scl, scl2, *pdata;
    float 	*trace;
    complex *ctrace;
    float   t0, t1;
    complex *cdata, tmp, ez;
    complex wa, *ctmp, da, *locdat;
    complex *cexrcv=(complex *) exrcv;

    /* define some constants */

    optn  = optncr(nt);
    nfreq = optn/2 + 1;
    df    = 1.0/(optn*dt);
    dom   = 2.0*M_PI*df;
    iomin = (int)MIN((fmin*dt*optn), (nfreq-1));
    iomin = MAX(iomin, 1);
    iomax = MIN((int)(fmax*dt*optn), (nfreq-1));

    /* transformation of shot record to frequency domain  */

    trace  = (float *)calloc(optn,sizeof(float));
    ctrace = (complex *)malloc(optn*sizeof(complex));
    cdata  = (complex *)calloc(nxm*nfreq, sizeof(complex));

    if (conjg) scl = -1.0;
    else scl = 1.0;

    sign  = -1;
    for (ix = 0; ix < nx; ix++) {

        memcpy(trace,&data[ix*nt],nt*sizeof(float));
        if (optn > nt)
            memset( &trace[nt], 0, sizeof(float)*(optn-nt) );

        rc1fft(trace,ctrace,optn,sign);

        ixrcv = NINT((xrcv[ix]-ox)/dxm);
        if (ixrcv < 0 || ixrcv > nxm-1) {
            fprintf(stderr,"kwZoMigr: ixrcv %f (%d) outside model\n", xrcv[ix], ixrcv);
            continue;
        }
        for (iom=0; iom<nfreq; iom++) {
            /* positioning of shot record into velocity model */
            cdata[iom*nxm+ixrcv].r = ctrace[iom].r;
            cdata[iom*nxm+ixrcv].i = ctrace[iom].i*scl;
        }
    }

    /* determine aperture to be calculated */

    ixo = nxm;
    ixn = 0;
    for (ix = 0; ix < nx; ix++) {
        ixrcv = NINT((xrcv[ix]-ox)/dxm);
        if (ixrcv < ixo) ixo = ixrcv;
        if (ixrcv > ixn) ixn = ixrcv;
    }
    ixmin = MAX(0, ixo-ixb-1);
    ixmax = MIN(ixn+ixa+1, nxm-1);
    nx    = (ixmax-ixmin)+1;

    if (verbose>=2) {
        vmess("kwZoMigr: calculation aperture: %.2f (%d) <--> %.2f (%d) (%d positions)", ixmin*dxm+ox, ixmin, ixmax*dxm+ox, ixmax, nx);
    }

    /* define some constants */

    scl   = 2.0/nfreq;
    scl   = 1.0/(dt*dxm*dxm);
    nkx   = optncc(2*ntap+nxm);
    ntap  = (nkx-nxm)/2;
    scl2  = 1.0/nkx;
    dkx   = 2.0*M_PI/(nkx*dxm);

    taper = (float *)malloc(ntap*sizeof(float));
    for (ix = 0; ix < ntap; ix++) {
        taper[ix] = exp(-1.0*(pow((0.4*(ntap-ix)/ntap), 2)));
    }

    /* calculate image at depth = 0 */

    if(izrcv==0 ) {
        for (ix = ixmin; ix <= ixmax; ix++) {
            for (iom = iomin; iom <= iomax; iom++) {
                image[ix*nzm+0] += scl*cdata[iom*nxm+ix].r;
            }
        }
    }

    t0 = wallclock_time();

    locdat  = (complex *)malloc(nkx*sizeof(complex));

    /* start extrapolation for all frequencies, depths and x-positions */

    for (iom = iomin; iom <= iomax; iom++) {
        memset(locdat,0,nkx*sizeof(complex));

        for (ix = ixmin; ix <= ixmax; ix++) {
            locdat[ntap+ix] = cdata[iom*nxm+ix];
        }
        om = iom*dom;
        d  = izrcv;

        /* start extrapolation of receiver arrays */

        for (; d < ndepth; d++) {

            /* transform to wavenumber domain */

            cc1fft(locdat, nkx, 1);

            /* Extrapolation of data */

            c = 0.0;
            for (ix = ixmin; ix <= ixmax; ix++) c += velmod[d*nxm+ix];
            k = nx*om/c;
            k2 = k*k;

            /* kx = 0 */
            ez.r = cos(k*dz);
            ez.i = -sin(k*dz);

            tmp.r  = ez.r*locdat[0].r;
            tmp.r += ez.i*locdat[0].i;
            tmp.i  = ez.r*locdat[0].i;
            tmp.i -= ez.i*locdat[0].r;
            locdat[0] = tmp;

            /* kx != 0 */
            endkx = MIN((int)(k/dkx),nkx/2);
            for (ikx = 1; ikx <= endkx; ikx++) {
                kx  = ikx*dkx;
                kx2 = kx*kx;
                kz2 = k2 - kx2;

                ez.r = cos(sqrt(kz2)*dz);
                ez.i = -sin(sqrt(kz2)*dz);

                tmp.r  = ez.r*locdat[ikx].r;
                tmp.r += ez.i*locdat[ikx].i;
                tmp.i  = ez.r*locdat[ikx].i;
                tmp.i -= ez.i*locdat[ikx].r;
                locdat[ikx] = tmp;

                tmp.r  = ez.r*locdat[nkx-ikx].r;
                tmp.r += ez.i*locdat[nkx-ikx].i;
                tmp.i  = ez.r*locdat[nkx-ikx].i;
                tmp.i -= ez.i*locdat[nkx-ikx].r;
                locdat[nkx-ikx] = tmp;
            }

            /* transform data back to space domain */

            cc1fft(locdat, nkx, -1);

            for (j = 0; j < nkx; j++) {
                locdat[j].r *= scl2;
                locdat[j].i *= scl2;
            }

            /* imaging condition */

            for (ix = ixmin; ix <= ixmax; ix++) {
                image[ix*nzm+d+1]+= scl*locdat[ntap+ix].r;
            }

            /* save extrapolated field at requested depth */

            if	(d==ndepthex-1) {
                if ( cexrcv != NULL ) {
                    for (ix = 0; ix < nxm; ix++) {
                        cexrcv[iom*nxm+ix] = locdat[ntap+ix];
                    }
                }
            }

            /* taper extrapolated data at edges */

            for (j = 0; j < ntap; j++) {
                locdat[j].r *= taper[j];
                locdat[j].i *= taper[j];
                locdat[nkx-j-1].r *= taper[j];
                locdat[nkx-j-1].i *= taper[j];
            }

        } /* end of depth loop */

    } /* end of iom loop */


    if(exrcv) wx2xt(cexrcv, exrcv, optn, nxm, nxm, optn);

    free(locdat);
    free(cdata);
    free(taper);

    t1 = wallclock_time();
    vmess("kwZoMigr took: %f seconds", t1-t0);

    return;
}
Example #7
0
int defineSource(wavPar wav, srcPar src, float **src_nwav, int reverse, int verbose)
{
    FILE    *fp;
    size_t  nread;
    int optn, nfreq, i, j, k, iwmax, tracesToDo;
	int iw, n1, namp;
    float scl, d1, df, deltom, om, tshift;
	float amp1, amp2, amp3;
	float *trace, maxampl;
    complex *ctrace, tmp;
    segy hdr;
    
	n1 = wav.nt;
	if (wav.random) { /* initialize random sequence */
		srand48(wav.seed+1);
		seedCMWC4096();
		for (i=0; i<8192; i++) {
			amp1 = dcmwc4096();
		}

	}
	else {

/* read first header and last byte to get file size */

    	fp = fopen( wav.file_src, "r" );
    	assert( fp != NULL);
    	nread = fread( &hdr, 1, TRCBYTES, fp );
    	assert(nread == TRCBYTES);
	
/* read all traces */

		tracesToDo = wav.nx;
		i = 0;
		while (tracesToDo) {
			memset(&src_nwav[i][0],0,n1*sizeof(float));
        	nread = fread(&src_nwav[i][0], sizeof(float), hdr.ns, fp);
        	assert (nread == hdr.ns);
	
        	nread = fread( &hdr, 1, TRCBYTES, fp );
        	if (nread==0) break;
			tracesToDo--;
			i++;
		}
    	fclose(fp);
	}


	optn = optncr(2*n1);
	nfreq = optn/2 + 1;
	ctrace = (complex *)calloc(nfreq,sizeof(complex));
	trace = (float *)calloc(optn,sizeof(float));

	df     = 1.0/(optn*wav.dt);
    deltom = 2.*M_PI*df;
	scl    = 1.0/optn;
	maxampl=0.0;
	iwmax = nfreq;

	for (i=0; i<wav.nx; i++) {
		if (wav.random) {
			randomWavelet(wav, src, &src_nwav[i][0], src.tbeg[i], src.tend[i], verbose);
		}
		else {
			memset(&ctrace[0].r,0,nfreq*sizeof(complex));
			memset(&trace[0],0,optn*sizeof(float));
			memcpy(&trace[0],&src_nwav[i][0],n1*sizeof(float));
			rc1fft(trace,ctrace,optn,-1);
			/* Scale source from file with -j/w (=1/(jw)) for volume source injections
         	   no scaling is applied for volume source injection rates */
            if (src.injectionrate==0) {
                for (iw=1;iw<iwmax;iw++) {
                    om = 1.0/(deltom*iw);
                    tmp.r = om*ctrace[iw].i;
                    tmp.i = -om*ctrace[iw].r;
                    ctrace[iw].r = tmp.r;
                    ctrace[iw].i = tmp.i;
                }
            }
/*
*/
            if (src.type < 6) { // shift wavelet with +1/2 DeltaT due to staggered in time 
				tshift=0.5*wav.dt;
                for (iw=1;iw<iwmax;iw++) {
            		om = deltom*iw*tshift;
            		tmp.r = ctrace[iw].r*cos(-om) - ctrace[iw].i*sin(-om);
            		tmp.i = ctrace[iw].i*cos(-om) + ctrace[iw].r*sin(-om);
                    ctrace[iw].r = tmp.r;
                    ctrace[iw].i = tmp.i;
                }
			}

			/* zero frequency iw=0 set to 0 if the next sample has amplitude==0*/
			amp1 = sqrt(ctrace[1].r*ctrace[1].r+ctrace[1].i*ctrace[1].i);
			if (amp1 == 0.0) {
				ctrace[0].r = ctrace[0].i = 0.0;
			}
			else { /* stabilization for w=0: extrapolate amplitudes to 0 */
				amp2 = sqrt(ctrace[2].r*ctrace[2].r+ctrace[2].i*ctrace[2].i);
				amp3 = sqrt(ctrace[3].r*ctrace[3].r+ctrace[3].i*ctrace[3].i);
				ctrace[0].r = amp1+(2.0*(amp1-amp2)-(amp2-amp3));
				ctrace[0].i = 0.0;
				if (ctrace[1].r < 0.0) {
					ctrace[0].r *= -1.0;
				}
			}
			for (iw=iwmax;iw<nfreq;iw++) {
				ctrace[iw].r = 0.0;
				ctrace[iw].i = 0.0;
			}
			cr1fft(ctrace,trace,optn,1);
			/* avoid a (small) spike in the last sample 
			   this is done to avoid diffraction from last wavelet sample
			   which will act as a pulse */
			if (reverse) {
				for (j=0; j<n1; j++) src_nwav[i][j] = scl*(trace[n1-j-1]-trace[0]);
//				for (j=0; j<n1; j++) src_nwav[i][j] = scl*(trace[j]-trace[optn-1]);
			}
			else {
				for (j=0; j<n1; j++) src_nwav[i][j] = scl*(trace[j]-trace[optn-1]);
			}
		}

    }
    free(ctrace);
    free(trace);

/* use random amplitude gain factor for each source */
	if (src.amplitude > 0.0) {
		namp=wav.nx*10;
		trace = (float *)calloc(2*namp,sizeof(float));
		for (i=0; i<wav.nx; i++) {
			if (src.distribution) {
				scl = gaussGen()*src.amplitude;
				k = (int)MAX(MIN(namp*(scl+5*src.amplitude)/(10*src.amplitude),namp-1),0);
				d1 = 10.0*src.amplitude/(namp-1);
			}
			else {
				scl = (float)(drand48()-0.5)*src.amplitude;
				k = (int)MAX(MIN(namp*(scl+1*src.amplitude)/(2*src.amplitude),namp-1),0);
				d1 = 2.0*src.amplitude/(namp-1);
			}

			trace[k] += 1.0;
/*			trace[i] = scl; */
			if (wav.random) n1 = wav.nsamp[i];
			else n1 = wav.nt;
			for (j=0; j<n1; j++) {
				src_nwav[i][j] *= scl;
			}
    	}
		if (verbose>2) writesufile("src_ampl.su", trace, namp, 1, -5*src.amplitude, 0.0, d1, 1);
/*
		qsort(trace,wav.nx,sizeof(float), comp);
		for (i=0; i<wav.nx; i++) {
			scl = trace[i];
			trace[i] = normal(scl, 0.0, src.amplitude);
		}
		if (verbose>2) writesufile("src_ampl.su", trace, wav.nx, 1, -5*src.amplitude, 0.0, d1, 1);
*/

		free(trace);
	}

	if (verbose>3) writesufilesrcnwav("src_nwav.su", src_nwav, wav, wav.nt, wav.nx, 0.0, 0.0, wav.dt, 1);

/* set maximum amplitude in source file to 1.0 */
/*
	assert(maxampl != 0.0);
	scl = wav.dt/(maxampl);
	scl = 1.0/(maxampl);
	for (i=0; i<wav.nx; i++) {
		for (j=0; j<n1; j++) {
			src_nwav[i*n1+j] *= scl;
		}
    }
*/

    return 0;
}