saveseis.c 9.16 KB
Newer Older
Simone Butzer's avatar
Simone Butzer committed
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80
/*------------------------------------------------------------------------
 * Copyright (C) 2015 For the list of authors, see file AUTHORS.
 *
 * This file is part of IFOS3D.
 * 
 * IFOS3D is free software: you can redistribute it and/or modify
 * it under the terms of the GNU General Public License as published by
 * the Free Software Foundation, version 2.0 of the License only.
 * 
 * IFOS3D is distributed in the hope that it will be useful,
 * but WITHOUT ANY WARRANTY; without even the implied warranty of
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 * GNU General Public License for more details.
 * 
 * You should have received a copy of the GNU General Public License
 * along with IFOS3D. See file COPYING and/or 
 * <http://www.gnu.org/licenses/gpl-2.0.html>.
--------------------------------------------------------------------------*/

/*------------------------------------------------------------------------
 *   write seismograms to files 
 *  ----------------------------------------------------------------------*/

#include "fd.h"

void saveseis(FILE *fp, float **sectionvx, float **sectionvy,float **sectionvz,
float **sectionp, float **sectioncurl, float **sectiondiv,
int  **recpos, int  **recpos_loc, int ntr, float ** srcpos, int ishot,int ns, int obs, int iteration){ 
		
	extern int SEISMO, SEIS_FORMAT[6], MYID, RUN_MULTIPLE_SHOTS,NPROC,METHOD;	
	extern char  SEIS_FILE[STRING_SIZE], SEIS_OBS_FILE[STRING_SIZE];
	extern FILE *FP;

	char vxf[STRING_SIZE], vyf[STRING_SIZE], vzf[STRING_SIZE], curlf[STRING_SIZE], divf[STRING_SIZE], pf[STRING_SIZE],file_ext[5];
	char outfile[STRING_SIZE],infile[STRING_SIZE];
	char seisfile[STRING_SIZE];
	float * buf2,* buf1, **srcpos1=NULL;
	int nsrc=1,i,k,m,l,nt;
	FILE *fin, *fout;
	
		srcpos1=fmatrix(1,7,1,1);
		for (nt=1;nt<=7;nt++) srcpos1[nt][1]=srcpos[nt][ishot];
		
		switch (SEIS_FORMAT[0]){
		case 0: sprintf(file_ext,"sgy"); break;
		case 1: sprintf(file_ext,"su");  break;
		case 2: sprintf(file_ext,"txt"); break;
		case 3: sprintf(file_ext,"bin"); break;
		case 4: sprintf(file_ext,"sgy"); break;
		case 5: sprintf(file_ext,"sgy"); break;
		}
		/*note that internally "y" is used for the vertical coordinate,
		for usability reasons, we switch the "y" and "z" coordinate 
		so that "z" - as commonly used - denotes the depth (vertical direction)*/
		
		if(obs==0) sprintf(seisfile,"%s",SEIS_FILE);
		if(obs==1) sprintf(seisfile,"%s",SEIS_OBS_FILE);
		
	if(ntr>0){	
		if (RUN_MULTIPLE_SHOTS){
			sprintf(vxf,"%s_vx_it%d.%s.shot%d.%d",seisfile,iteration,file_ext,ishot,MYID);
			sprintf(vzf,"%s_vy_it%d.%s.shot%d.%d",seisfile,iteration,file_ext,ishot,MYID);
			sprintf(vyf,"%s_vz_it%d.%s.shot%d.%d",seisfile,iteration,file_ext,ishot,MYID);
			sprintf(curlf,"%s_rot_it%d.%s.shot%d.%d",seisfile,iteration,file_ext,ishot,MYID);
			sprintf(divf,"%s_div_it%d.%s.shot%d.%d",seisfile,iteration,file_ext,ishot,MYID);
			sprintf(pf,"%s_p_it%d.%s.shot%d.%d",seisfile,iteration,file_ext,ishot,MYID);
		}
		else{
			sprintf(vxf,"%s_vx_it%d.%s.%d",seisfile,iteration,file_ext,MYID);
			sprintf(vzf,"%s_vy_it%d.%s.%d",seisfile,iteration,file_ext,MYID);
			sprintf(vyf,"%s_vz_it%d.%s.%d",seisfile,iteration,file_ext,MYID);
			sprintf(curlf,"%s_rot_it%d.%s.%d",seisfile,iteration,file_ext,MYID);
			sprintf(divf,"%s_div_it%d.%s.%d",seisfile,iteration,file_ext,MYID);
			sprintf(pf,"%s_p_it%d.%s.%d",seisfile,iteration,file_ext,MYID);
			
		}
		

		switch (SEISMO){
		case 1 : /* particle velocities only */
tilman.metz's avatar
tilman.metz committed
81 82
			
			fprintf(fp,"\n PE %d is writing %d seismogramtraces (vx)   to \t %s \n",MYID,ntr,vxf);
Simone Butzer's avatar
Simone Butzer committed
83
			outseis(fp,fopen(vxf,"w"),1,sectionvx,recpos,recpos_loc, ntr,srcpos1,nsrc,ns,SEIS_FORMAT);
tilman.metz's avatar
tilman.metz committed
84
			fprintf(fp," PE %d is writing %d seismogramtraces (vy)   to \t %s \n",MYID,ntr,vyf);
Simone Butzer's avatar
Simone Butzer committed
85
			outseis(fp,fopen(vyf,"w"),2,sectionvy,recpos,recpos_loc, ntr,srcpos1,nsrc,ns,SEIS_FORMAT);
tilman.metz's avatar
tilman.metz committed
86
			fprintf(fp," PE %d is writing %d seismogramtraces (vz)   to \t %s \n",MYID,ntr,vzf);
Simone Butzer's avatar
Simone Butzer committed
87 88 89 90
			outseis(fp,fopen(vzf,"w"),3,sectionvz,recpos,recpos_loc, ntr,srcpos1,nsrc,ns,SEIS_FORMAT);
			
			break;
		case 2 : /* pressure only */
tilman.metz's avatar
tilman.metz committed
91
			fprintf(fp," PE %d is writing %d seismogramtraces (p)   to \t %s \n",MYID,ntr,pf);
Simone Butzer's avatar
Simone Butzer committed
92 93 94 95 96
			outseis(fp,fopen(pf,"w"), 0,sectionp, recpos,recpos_loc, ntr,srcpos1,nsrc,ns,SEIS_FORMAT);
			
			break;
		case 3 : /* curl and div only */
			
tilman.metz's avatar
tilman.metz committed
97
			fprintf(fp," PE %d is writing %d seismogramtraces (div)   to \t %s \n",MYID,ntr,divf);
Simone Butzer's avatar
Simone Butzer committed
98
			outseis(fp,fopen(divf,"w"),0,sectiondiv,recpos,recpos_loc,ntr,srcpos1,nsrc,ns,SEIS_FORMAT);
tilman.metz's avatar
tilman.metz committed
99
			fprintf(fp," PE %d is writing %d seismogramtraces (curl)  to \t %s \n",MYID,ntr,curlf);
Simone Butzer's avatar
Simone Butzer committed
100 101 102 103
			outseis(fp,fopen(curlf,"w"),0,sectioncurl,recpos,recpos_loc,ntr,srcpos1,nsrc,ns,SEIS_FORMAT);	
			
			break;	
		case 4 : /* everything */
tilman.metz's avatar
tilman.metz committed
104
			fprintf(fp,"\n PE %d is writing %d seismogramtraces (vx)   to \t %s \n",MYID,ntr,vxf);
Simone Butzer's avatar
Simone Butzer committed
105
			outseis(fp,fopen(vxf,"w"),1,sectionvx,recpos,recpos_loc, ntr,srcpos1,nsrc,ns,SEIS_FORMAT);
tilman.metz's avatar
tilman.metz committed
106
			fprintf(fp," PE %d is writing %d seismogramtraces (vy)   to \t %s \n",MYID,ntr,vyf);
Simone Butzer's avatar
Simone Butzer committed
107
			outseis(fp,fopen(vyf,"w"),2,sectionvy,recpos,recpos_loc, ntr,srcpos1,nsrc,ns,SEIS_FORMAT);
tilman.metz's avatar
tilman.metz committed
108
			fprintf(fp," PE %d is writing %d seismogramtraces (vz)   to \t %s \n",MYID,ntr,vzf);
Simone Butzer's avatar
Simone Butzer committed
109 110
			outseis(fp,fopen(vzf,"w"),3,sectionvz,recpos,recpos_loc, ntr,srcpos1,nsrc,ns,SEIS_FORMAT);
			
tilman.metz's avatar
tilman.metz committed
111
			fprintf(fp," PE %d is writing %d seismogramtraces (p)    to \t %s \n",MYID,ntr,pf);
Simone Butzer's avatar
Simone Butzer committed
112 113
			outseis(fp,fopen(pf,"w"), 0,sectionp, recpos,recpos_loc, ntr,srcpos1,nsrc,ns,SEIS_FORMAT);

tilman.metz's avatar
tilman.metz committed
114
			fprintf(fp," PE %d is writing %d seismogramtraces (div)  to \t %s \n",MYID,ntr,divf);
Simone Butzer's avatar
Simone Butzer committed
115
			outseis(fp,fopen(divf,"w"),0,sectiondiv,recpos,recpos_loc,ntr,srcpos1,nsrc,ns,SEIS_FORMAT);
tilman.metz's avatar
tilman.metz committed
116
			fprintf(fp," PE %d is writing %d seismogramtraces (curl) to \t %s \n",MYID,ntr,curlf);
Simone Butzer's avatar
Simone Butzer committed
117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252
			outseis(fp,fopen(curlf,"w"),0,sectioncurl,recpos,recpos_loc,ntr,srcpos1,nsrc,ns,SEIS_FORMAT);	
			break;
			
	      }     	
	}
	
	MPI_Barrier(MPI_COMM_WORLD);
	if(MYID==0){
		sprintf(outfile,"%s_vx_it%d.%s.shot%d",seisfile,iteration,file_ext,ishot);
		/*fprintf(FP,"outfile:%s",outfile);*/
		fout=fopen(outfile,"w");
		buf2=vector(1,ns);
		buf1=vector(1,240/sizeof(float));	
			for(i=0;i<=NPROC-1;i++){
				k=0;
				sprintf(infile,"%s_vx_it%d.%s.shot%d.%d",seisfile,iteration,file_ext,ishot,i);
				fin=fopen(infile,"r");
				if(fin!=0){
					fseek(fin,0,SEEK_END);
					k=ftell(fin);
					if(k%(240+4*ns)==0)l=k/(240+4*ns);
					else{ fprintf(FP,"attention: Error in seisfile");
						l=0;
					}
					/*fprintf(FP,"l=%d",l);
					fprintf(FP,"infile:%s,MYID=%d\n",infile,i);*/
					fseek(fin,0,SEEK_SET);
					
					for(m=1;m<=l;m++){
					  /*fprintf(FP,"m=%d",m);*/
					  
					  
					fread(buf1,sizeof(float),240/sizeof(float),fin);
					fread(buf2,sizeof(float),ns,fin);
							
					fwrite(buf1,sizeof(float),240/sizeof(float),fout);
					fwrite(buf2,sizeof(float),ns,fout);
					}
					fclose(fin);
					
				}
			}
		fclose(fout);
		free_vector(buf2,1,240/sizeof(float));
		free_vector(buf1,1,ns);
	}

	if(MYID==0){
		sprintf(outfile,"%s_vy_it%d.%s.shot%d",seisfile,iteration,file_ext,ishot);
		/*fprintf(FP,"outfile:%s",outfile);*/
		fout=fopen(outfile,"w");
		buf2=vector(1,ns);
		buf1=vector(1,240/sizeof(float));
		k=0;
		
		for(i=0;i<=NPROC-1;i++){
			sprintf(infile,"%s_vy_it%d.%s.shot%d.%d",seisfile,iteration,file_ext,ishot,i);
			fin=fopen(infile,"r");
			if(fin!=0){
				fseek(fin,0,SEEK_END);
				k=ftell(fin);
				if(k%(240+4*ns)==0)l=k/(240+4*ns);
				else{ fprintf(FP,"attention: Error in seisfile");
					l=0;
				}
				/*fprintf(FP,"l=%d",l);
				fprintf(FP,"infile:%s,MYID=%d",infile,MYID);*/
				fseek(fin,0,SEEK_SET);
				
				for(m=1;m<=l;m++){
				  /*fprintf(FP,"m=%d",m);*/
				  
				  
				fread(buf1,sizeof(float),240/sizeof(float),fin);
				fread(buf2,sizeof(float),ns,fin);
						
				fwrite(buf1,sizeof(float),240/sizeof(float),fout);
				fwrite(buf2,sizeof(float),ns,fout);
				}
				fclose(fin);
			}
		}
		fclose(fout);
		free_vector(buf2,1,240/sizeof(float));
		free_vector(buf1,1,ns);
	}
	
	if(MYID==0){
		sprintf(outfile,"%s_vz_it%d.%s.shot%d",seisfile,iteration,file_ext,ishot);
		/*fprintf(FP,"outfile:%s",outfile);*/
		fout=fopen(outfile,"w");
		buf2=vector(1,ns);
		buf1=vector(1,240/sizeof(float));
		k=0;
		for(i=0;i<=NPROC-1;i++){
			sprintf(infile,"%s_vz_it%d.%s.shot%d.%d",seisfile,iteration,file_ext,ishot,i);
			fin=fopen(infile,"r");
			if(fin!=0){
				fseek(fin,0,SEEK_END);
				k=ftell(fin);
				if(k%(240+4*ns)==0)l=k/(240+4*ns);
				else{ fprintf(FP,"attention: Error in seisfile");
					l=0;
				}
				/*fprintf(FP,"l=%d",l);
				fprintf(FP,"infile:%s,MYID=%d",infile,MYID);*/
				fseek(fin,0,SEEK_SET);
				for(m=1;m<=l;m++){
				  				  
				  
				fread(buf1,sizeof(float),240/sizeof(float),fin);
				fread(buf2,sizeof(float),ns,fin);
						
				fwrite(buf1,sizeof(float),240/sizeof(float),fout);
				fwrite(buf2,sizeof(float),ns,fout);
				}
				fclose(fin);
				
			}
		}
		fclose(fout);
		free_vector(buf2,1,240/sizeof(float));
		free_vector(buf1,1,ns);
	}
	free_matrix(srcpos1,1,7,1,1);
	MPI_Barrier(MPI_COMM_WORLD);
	
	if(iteration>0&&ntr>0&&METHOD&&!obs){
		sprintf(infile,"%s_vx_it%d.%s.shot%d.%d",seisfile,iteration,file_ext,ishot,MYID);
		remove(infile);
		sprintf(infile,"%s_vy_it%d.%s.shot%d.%d",seisfile,iteration,file_ext,ishot,MYID);
		remove(infile);
		sprintf(infile,"%s_vz_it%d.%s.shot%d.%d",seisfile,iteration,file_ext,ishot,MYID);
		remove(infile);
	}
}