checkfd_ssg.c 29.5 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 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 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 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325
/*------------------------------------------------------------------------
 * 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>.
--------------------------------------------------------------------------*/

/*--------------------------------------------------------------------------*/
/* $Id: checkfd_ssg.c,v 1.1.1.1 2009/12/09 14:29:22 sjetschny Exp $ */
/*-------------------------------------------------------------
 * - Check SSG O(2,4) for stability and grid dispersion.
 *   If the stability criterion is not fullfilled the program will
 *   terminate.                   
 * - Check data output directories and files for accessibility.
 +   If any is not accessable form any PE the program will terminate.
 *  ----------------------------------------------------------*/

#include <libgen.h>
#include <unistd.h>
#include "fd.h"

void err2(char errformat[],char errfilename[]){
	char outtxt[STRING_SIZE];
	sprintf(outtxt,errformat,errfilename);
	err(outtxt);
}

void checkfd(FILE *fp, float *** prho, float *** ppi, float *** pu,
float *** ptaus, float *** ptaup, float *peta, float **srcpos, int nsrc, int **recpos, int ntr){

	extern float DX, DY, DZ, DT, TS, TIME, TSNAP2;
	extern float XREC1, XREC2, YREC1, YREC2, ZREC1, ZREC2;
        extern int NX, NY, NZ, L, MYID, IDX, IDY, IDZ, FW, POS[4], NT, NDT, NDTSHIFT;
	extern int FDCOEFF;
	extern int READREC, NPROCX,NPROCY,NPROCZ, FW, ABS_TYPE, SRCREC, FREE_SURF;
	extern int SNAP, SEISMO, SEIS_FORMAT[6], SNAP_FORMAT;
	/*extern int RUN_MULTIPLE_SHOTS; no determination is done for the output check whether the simulation runs with one or multiple shot
		-> directorys specified in input file should work in both cases */
	extern char SEIS_FILE[STRING_SIZE], SNAP_FILE[STRING_SIZE];
	extern char SOURCE_FILE[STRING_SIZE], REC_FILE[STRING_SIZE];

	/* local variables */

	float  c=0.0, cmax_p=0.0, cmin_p=1e9, cmax_s=0.0, cmin_s=1e9, fmax, cwater=1.0e-1;
	float snapoutx=0.0, snapouty=0.0, snapoutz=0.0, dhmax, dhmin;
	float  cmax=0.0, cmin=1e9, sum, dtstab, ts, cmax_r, cmin_r, g=0.0, gamma=0.0;
	float srec_minx=DX*NX*NPROCX+1, srec_miny=DY*NY*NPROCY+1, srec_minz=DZ*NZ*NPROCZ+1;
	float srec_maxx=-1.0, srec_maxy=-1.0, srec_maxz=-1.0;
	const float w=2.0*PI/TS; /*center frequency of source*/
	int nfw=FW;
	int i, j, k, l, ny1=1, nx, ny, nz;
	extern int FDORDER;
	char xfile[STRING_SIZE], errormessage[STRING_SIZE], xmod[4], file_ext[8];
	FILE *fpcheck;

	nx=NX; ny=NY; nz=NZ;

	/* low Q frame not yet applied as a absorbing boundary */
	/* if (!FREE_SURF) ny1=1+nfw;*/
	nfw=0;
	

	/*printf(" checkfd: TS= %f \n",TS);*/
	
	/* find maximum model phase velocity of shear waves at infinite
	      frequency within the whole model */
	for (k=1+nfw;k<=(nz-nfw);k++){
		for (i=1+nfw;i<=(nx-nfw);i++){
			for (j=ny1;j<=(ny-nfw);j++){
				if (prho[j][i][k]==0.0)
					printf("prho=0 detected from MYID=%d at i=%d,j=%d,k=%d\n",MYID,i,j,k);
					else{
						if(L){
						c=sqrt(pu[j][i][k]*(1.0+L*ptaus[j][i][k])/prho[j][i][k]);}
						else c=sqrt(pu[j][i][k]/prho[j][i][k]);
					}
					
				if (cmax_s<c) cmax_s=c;
				/* find minimum model phase velocity of shear waves at center
					       frequency of the source */
				sum=0.0;
				for (l=1;l<=L;l++){
					ts=DT/peta[l];
					sum=sum+((w*w*ptaus[j][i][k]*ts*ts)/(1.0+w*w*ts*ts));
				}
				if (prho[j][i][k]==0.0)
					printf("prho=0 detected from MYID=%d at i=%d,j=%d,k=%d\n",MYID,i,j,k);
					else
					c=sqrt(pu[j][i][k]*(1.0+sum)/prho[j][i][k]);
				if ((c>cwater)&&(cmin_s>c)) cmin_s=c;
			}
		}
	}



	/* find maximum model phase velocity of P-waves at infinite
		 frequency within the whole model */
	for (k=1+nfw;k<=(nz-nfw);k++){
		for (i=1+nfw;i<=(nx-nfw);i++){
			for (j=ny1;j<=(ny-nfw);j++){
				if (prho[j][i][k]==0.0)
					printf("prho=0 detected from MYID=%d at i=%d,j=%d,k=%d\n",MYID,i,j,k);
					else{
						if(L)	c=sqrt(ppi[j][i][k]*(1.0+L*ptaup[j][i][k])/prho[j][i][k]);
						else    c=sqrt(ppi[j][i][k]/prho[j][i][k]);
					  }
				if (cmax_p<c) cmax_p=c;
				/* find minimum model phase velocity of P-waves at center
					       frequency of the source */
				sum=0.0;
				for (l=1;l<=L;l++){
					ts=DT/peta[l];
					sum=sum+((w*w*ptaup[j][i][k]*ts*ts)/(1.0+w*w*ts*ts));
				}
				if (prho[j][i][k]==0.0)
					printf("prho=0 detected from MYID=%d at i=%d,j=%d,k=%d\n",MYID,i,j,k);
					else
					c=sqrt(ppi[j][i][k]*(1.0+sum)/prho[j][i][k]);
				if (cmin_p>c) cmin_p=c;
			}
		}
	}


	fprintf(fp,"\n\n\n **Message from checkfd (printed by PE %d):\n",MYID);
	fprintf(fp," Minimum and maximum P-wave and S-wave velocities within subvolume: \n ");
	fprintf(fp," MYID\t Vp_min(f=fc) \t Vp_max(f=inf) \t Vs_min(f=fc) \t Vsmax(f=inf) \n");
	fprintf(fp," %d \t %8.2f \t %8.2f \t %8.2f \t %8.2f \n\n\n", MYID, cmin_p, cmax_p, cmin_s, cmax_s);

	if (cmax_s>cmax_p) cmax=cmax_s; 
	else cmax=cmax_p;
	if (cmin_s<cmin_p) cmin=cmin_s; 
	else cmin=cmin_p;

	/* find global maximum for Vp and global minimum for Vs*/
	MPI_Allreduce(&cmax,&cmax_r,1,MPI_FLOAT,MPI_MAX,MPI_COMM_WORLD);
	MPI_Allreduce(&cmin,&cmin_r,1,MPI_FLOAT,MPI_MIN,MPI_COMM_WORLD);
	cmax=cmax_r;
	cmin=cmin_r;	

	/* if (MYID==0){		checkfd is performed by MID=0 only */
	
	fprintf(fp," Global values for entire model: \n");
	fprintf(fp," Vp_max= %8.2f m/s \t Vs_min=%8.2f m/s \n", cmax,cmin);
		
	
	/* estimate number of grid points per minimum wavelength to avoid numerical dispersion */
	/* Taylor */
	if (FDCOEFF==1) {
		switch (FDORDER){
		case 2: g=12.0;
			break;
		case 4: g=8.0;
			break;
		case 6: g=7.0;
			break;
		case 8: g=6.0;
			break;
		case 10: g=5.0;
			break;
		case 12: g=4.0;
			break;
		default: err(" invalid FDORDER in checkfd_ssg() !");
		        break;
		}
	}
	
	/* Holberg */
	if (FDCOEFF==2) {
	        switch (FDORDER){
		case 2: g=12.0;
			break;
		case 4: g=8.32;
			break;
		case 6: g=4.77;
			break;
		case 8: g=3.69;
			break;
		case 10: g=3.19;	
			break;
		case 12: g=2.91;
			break;
		default: err(" invalid FDORDER in checkfd_ssg() !");
		        break;
		}
	}

	MPI_Barrier(MPI_COMM_WORLD);
	
	
	fprintf(fp,"\n\n ------------------ CHECK OUTPUT FILES --------------------------\n");	
	
	/* The original checks might delete files accidentally that would not be overwritten anyway. 
	   and did not test accessibility from all CPUs which may be vary, especially in distributed clusters */
			
	/*Checking SNAP Output */
	/*-------------------- */
	
	if (SNAP>0){
		
		switch(SNAP_FORMAT){
			case 1: 
				sprintf(file_ext,".su");
				strcpy(xmod,"ab");
				break;
			case 2: 
				sprintf(file_ext,".asc");
				strcpy(xmod,"a");
				break;
			case 3: 
				sprintf(file_ext,".bin");
				strcpy(xmod,"ab");
				break;
			default: err(" Sorry. Snapshot format (SNAP_FORMAT) unknown. \n");
		}		
		
		
		fprintf(fp," Check accessibility for snapshot files ... \n");
		
		switch(SNAP){
			case 1 :
				sprintf(xfile,"%s%s.x.%i%i%i",SNAP_FILE,file_ext,POS[1],POS[2],POS[3]);
				fprintf(fp,"    Check accessibility for snapshot file %s... \n",xfile);
				sprintf(errormessage,"PE %i cannot write snapshot %s!",MYID,xfile);
				if (((fpcheck=fopen(xfile,xmod))==NULL))err(errormessage);
				else fclose(fpcheck); 
				remove(xfile);
				
				sprintf(xfile,"%s%s.y.%i%i%i",SNAP_FILE,file_ext,POS[1],POS[2],POS[3]);
				fprintf(fp,"    Check accessibility for snapshot file %s... \n",xfile);
				sprintf(errormessage,"PE %i cannot write snapshot %s!",MYID,xfile);
				if (((fpcheck=fopen(xfile,xmod))==NULL))err(errormessage);
				else fclose(fpcheck); 
				remove(xfile);
							
				sprintf(xfile,"%s%s.z.%i%i%i",SNAP_FILE,file_ext,POS[1],POS[2],POS[3]);
				fprintf(fp,"    Check accessibility for snapshot file %s... \n",xfile);
				sprintf(errormessage,"PE %i cannot write snapshot %s!",MYID,xfile);
				if (((fpcheck=fopen(xfile,xmod))==NULL))err(errormessage);
				else fclose(fpcheck); 
				remove(xfile);
				
				break;			
			case 2 :
				sprintf(xfile,"%s%s.p.%i%i%i",SNAP_FILE,file_ext,POS[1],POS[2],POS[3]);
				//fprintf(fp,"    Check accessibility for snapshot file %s... \n",xfile);
				sprintf(errormessage,"PE %i cannot write snapshot %s!",MYID,xfile);
				if (((fpcheck=fopen(xfile,xmod))==NULL))err(errormessage);
				else fclose(fpcheck); 
				remove(xfile);
				break;
			case 4 :
				sprintf(xfile,"%s%s.x.%i%i%i",SNAP_FILE,file_ext,POS[1],POS[2],POS[3]);
				fprintf(fp,"    Check accessibility for snapshot file %s... \n",xfile);
				sprintf(errormessage,"PE %i cannot write snapshot %s!",MYID,xfile);
				if (((fpcheck=fopen(xfile,xmod))==NULL))err(errormessage);
				else fclose(fpcheck); 
				remove(xfile);
				
				sprintf(xfile,"%s%s.y.%i%i%i",SNAP_FILE,file_ext,POS[1],POS[2],POS[3]);
				fprintf(fp,"    Check accessibility for snapshot file %s... \n",xfile);
				sprintf(errormessage,"PE %i cannot write snapshot %s!",MYID,xfile);
				if (((fpcheck=fopen(xfile,xmod))==NULL))err(errormessage);
				else fclose(fpcheck); 
				remove(xfile);
				
				sprintf(xfile,"%s%s.z.%i%i%i",SNAP_FILE,file_ext,POS[1],POS[2],POS[3]);
				fprintf(fp,"    Check accessibility for snapshot file %s... \n",xfile);
				sprintf(errormessage,"PE %i cannot write snapshot %s!",MYID,xfile);
				if (((fpcheck=fopen(xfile,xmod))==NULL))err(errormessage);
				else fclose(fpcheck); 
				remove(xfile);
				
			case 3 :	
				sprintf(xfile,"%s%s.div.%i%i%i",SNAP_FILE,file_ext,POS[1],POS[2],POS[3]);
				fprintf(fp,"    Check accessibility for snapshot file %s... \n",xfile);
				sprintf(errormessage,"PE %i cannot write snapshot %s!",MYID,xfile);
				if (((fpcheck=fopen(xfile,xmod))==NULL))err(errormessage);
				else fclose(fpcheck); 
				remove(xfile);
				
				sprintf(xfile,"%s%s.rot.%i%i%i",SNAP_FILE,file_ext,POS[1],POS[2],POS[3]);		   
				fprintf(fp,"    Check accessibility for snapshot file %s... \n",xfile);
				sprintf(errormessage,"PE %i cannot write snapshot %s!",MYID,xfile);
				if (((fpcheck=fopen(xfile,xmod))==NULL))err(errormessage);
				else fclose(fpcheck); 
				remove(xfile);
				
				break;		   
		}
	}
	
		
	/*Checking SEISMOGRAM Output Particle velocities */
	/*-------------------------------------- */
	if (SEISMO>0) {	
		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;
		}
	
		/*if ((RUN_MULTIPLE_SHOTS)||(SEIS_FORMAT[4])) {
		 possibly many files ... -> perform check of write and execute permission for directories, only 
		
			fprintf(fp," Check accessibility for seismogram files (multiple shots) ... \n");
			
			switch (SEISMO){
			case 1:  particle velocities only 
326
				sprintf(xfile,"%s.ifos.tmp",dirname(SEIS_FILE_VX));
Simone Butzer's avatar
Simone Butzer committed
327
				if (access(xfile,W_OK|X_OK)==-1) err2(" PE %d cannot write seismograms to %s!",xfile);
328
				sprintf(xfile,"%s.ifos.tmp",dirname(SEIS_FILE_VY));
Simone Butzer's avatar
Simone Butzer committed
329
				if (access(xfile,W_OK|X_OK)==-1) err2(" PE %d cannot write seismograms to %s!",xfile);
330
				sprintf(xfile,"%s.ifos.tmp",dirname(SEIS_FILE_VZ));
Simone Butzer's avatar
Simone Butzer committed
331 332 333
				if (access(xfile,W_OK|X_OK)==-1) err2(" PE %d cannot write seismograms to %s!",xfile);
				break;
			case 2 :  pressure only 
334
				sprintf(xfile,"%s.ifos.tmp",dirname(SEIS_FILE_P));
Simone Butzer's avatar
Simone Butzer committed
335 336 337
				if (access(xfile,W_OK|X_OK)==-1) err2(" PE %d cannot write seismograms to %s!",xfile);
				break;
			case 4 :
338
				sprintf(xfile,"%s.ifos.tmp",dirname(SEIS_FILE_VX));
Simone Butzer's avatar
Simone Butzer committed
339
				if (access(xfile,W_OK|X_OK)==-1) err2(" PE %d cannot write seismograms to %s!",xfile);	
340
				sprintf(xfile,"%s.ifos.tmp",dirname(SEIS_FILE_VY));
Simone Butzer's avatar
Simone Butzer committed
341
				if (access(xfile,W_OK|X_OK)==-1) err2(" PE %d cannot write seismograms to %s!",xfile);
342
				sprintf(xfile,"%s.ifos.tmp",dirname(SEIS_FILE_VZ));
Simone Butzer's avatar
Simone Butzer committed
343 344
				if (access(xfile,W_OK|X_OK)==-1) err2(" PE %d cannot write seismograms to %s!",xfile);							
			case 3 : curl and div only 
345
				sprintf(xfile,"%s.ifos.tmp",dirname(SEIS_FILE_DIV));
Simone Butzer's avatar
Simone Butzer committed
346
				if (access(xfile,W_OK|X_OK)==-1) err2(" PE %d cannot write seismograms to %s!",xfile);
347
				sprintf(xfile,"%s.ifos.tmp",dirname(SEIS_FILE_CURL));
Simone Butzer's avatar
Simone Butzer committed
348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719
				if (access(xfile,W_OK|X_OK)==-1) err2(" PE %d cannot write seismograms to %s!",xfile);	
				break;
			}		
		}
		else{*/

			fprintf(fp," Check accessibility for seismogram files  ... \n");
		
			if (SEIS_FORMAT[0]==2) strcpy(xmod,"a");
			else strcpy(xmod,"ab");			
			switch (SEISMO){
			case 1: /* particle velocities only */
				sprintf(xfile,"%s_vx.%s.%d.tmp",SEIS_FILE,file_ext,MYID);
				fprintf(fp,"    Check accessibility for seismogram file %s... \n",xfile); /*in case of number of PE's=500, there will be 500 messages, too many to be displayed! */
				sprintf(errormessage,"PE %i cannot write seismogram file %s!",MYID,xfile);
				//if (access(xfile,W_OK|X_OK)==-1) err(errormessage);
				if (((fpcheck=fopen(xfile,xmod))==NULL) && (MYID==0))err(errormessage);
				else fclose(fpcheck); 
				remove(xfile);
				
				sprintf(xfile,"%s_vy.%s.%d.tmp",SEIS_FILE,file_ext,MYID);
				fprintf(fp,"    Check accessibility for seismogram file %s... \n",xfile);
				if (((fpcheck=fopen(xfile,xmod))==NULL) && (MYID==0))err("PE 0 cannot write seismograms!");
				else fclose(fpcheck);
				remove(xfile);
				
				sprintf(xfile,"%s_vz.%s.%d.tmp",SEIS_FILE,file_ext,MYID);
				fprintf(fp,"    Check accessibility for seismogram file %s... \n",xfile);
				if (((fpcheck=fopen(xfile,xmod))==NULL) && (MYID==0))err("PE 0 cannot write seismograms!");
				else fclose(fpcheck); 
				remove(xfile);
				
				break;
			case 2 : /* pressure only */
				sprintf(xfile,"%s_p.%s.%d.tmp",SEIS_FILE,file_ext,MYID);
				fprintf(fp,"    Check accessibility for seismogram file %s... \n",xfile);
				if (((fpcheck=fopen(xfile,xmod))==NULL) && (MYID==0))err("PE 0 cannot write seismograms!");
				else fclose(fpcheck); 
				remove(xfile);
								
				break;
				
			case 4 : /* everything */
				sprintf(xfile,"%s_vx.%s.%d.tmp",SEIS_FILE,file_ext,MYID);
				fprintf(fp,"    Check accessibility for seismogram file %s... \n",xfile); /*in case of number of PE's=500, there will be 500 messages, too many to be displayed! */
				sprintf(errormessage,"PE %i cannot write seismogram file %s!",MYID,xfile);
				//if (access(xfile,W_OK|X_OK)==-1) err(errormessage);
				if (((fpcheck=fopen(xfile,xmod))==NULL) && (MYID==0))err(errormessage);
				else fclose(fpcheck); 
				remove(xfile);
				
				sprintf(xfile,"%s_vy.%s.%d.tmp",SEIS_FILE,file_ext,MYID);
				fprintf(fp,"    Check accessibility for seismogram file %s... \n",xfile);
				if (((fpcheck=fopen(xfile,xmod))==NULL) && (MYID==0))err("PE 0 cannot write seismograms!");
				else fclose(fpcheck);
				remove(xfile);
				
				sprintf(xfile,"%s_vz.%s.%d.tmp",SEIS_FILE,file_ext,MYID);
				fprintf(fp,"    Check accessibility for seismogram file %s... \n",xfile);
				if (((fpcheck=fopen(xfile,xmod))==NULL) && (MYID==0))err("PE 0 cannot write seismograms!");
				else fclose(fpcheck); 
				remove(xfile);


			case 3 : /* curl and div only */
				sprintf(xfile,"%s_div.%s.%d.tmp",SEIS_FILE,file_ext,MYID);
				fprintf(fp,"    Check accessibility for seismogram file %s... \n",xfile);
				if (((fpcheck=fopen(xfile,xmod))==NULL) && (MYID==0))err("PE 0 cannot write seismograms!");
				else fclose(fpcheck);  
				remove(xfile);
				
				sprintf(xfile,"%s_rot.%s.%d.tmp",SEIS_FILE,file_ext,MYID);
				fprintf(fp,"    Check accessibility for seismogram file %s... \n",xfile);
				if (((fpcheck=fopen(xfile,xmod))==NULL) && (MYID==0))err("PE 0 cannot write seismograms!");
				else fclose(fpcheck); 
				remove(xfile);
						
				break;
		
			
		}		
	}

	
	
	/******************************************************************************************/
	
        /* calculate maximum grid point distance */
	dhmax=DX;
	if(DY>dhmax){dhmax=DY;}
	if(DZ>dhmax){dhmax=DZ;}	
	
	fprintf(fp,"\n\n ------------------ CHECK FOR GRID DISPERSION --------------------\n");
	fprintf(fp," To satisfactorily limit grid dispersion the number of gridpoints \n");
	fprintf(fp,"    per minimum wavelength (of S-waves) should be 6 (better more).\n");
	fprintf(fp," Here the minimum wavelength is assumed to be minimum model phase velocity \n");
	fprintf(fp,"    (of S-waves) at maximum frequency of the source\n");
	fprintf(fp,"    devided by maximum frequency of the source.\n");
	fmax=2.0/TS;
	fprintf(fp," Maximum frequency of the source is approximately %8.2f Hz\n",fmax);
	fprintf(fp,"    (= two times FC specified in the sources file)\n");
	fprintf(fp," The minimum wavelength (of S-waves) in the following simulation will\n");
	fprintf(fp,"    be %4.4f meter.\n", cmin/fmax);
	fprintf(fp," Thus, at the chosen FD_ORDER of %i the recommended value for DH is %4.4f meter \n    (%2.2f gridpoints per minimum wavelength).\n", FDORDER,cmin/fmax/g,g);
	fprintf(fp," You have specified DH= %4.4f meter.\n\n", dhmax);
	if ((dhmax>(cmin/fmax/g))&&(MYID==0))
	warning(" Grid dispersion will influence wave propagation, choose smaller grid spacing (DH).");
        
	/*  gamma for stability estimation (Taylor) */
	/*switch (FDORDER){
	case 2: gamma=1.0;
		break;
	case 4: gamma=7.0/6.0;
		break;
	case 6: gamma=149.0/120.0;
		break;
	case 8: gamma=2161.0/1680.0;
		break;
	case 10: gamma=53089.0/40320.0;
		break;
	case 12: gamma=1187803.0/887040.0;
		break;
	default: err(" invalid FDORDER in checkfd_acoustic() !");
	        break;
	}*/
	
	/*  gamma for stability estimation (Holberg) */
	switch (FDORDER){
	case 2: gamma=1.0;
		break;
	case 4: gamma=1.184614;
		break;
	case 6: gamma=1.283482;
		break;
	case 8: gamma=1.345927;
		break;
	case 10: gamma=1.38766;
		break;
	case 12: gamma=1.417065;
		break;
	default: err(" invalid FDORDER in checkfd_acoustic() !");
	        break;
	}


        /* calculate minimum grid point distance */
	dhmin=DX;
	if(DY<dhmin){dhmin=DY;}
	if(DZ<dhmin){dhmin=DZ;}
	
	fprintf(fp," \n\n ----------------------- CHECK FOR STABILITY ---------------------\n");
	fprintf(fp," The following simulation is stable provided that\n\n");
	fprintf(fp," \t p=cmax*DT/DH < 1/(gamma*sqrt(3)) = %2.3f,\n\n",1/(gamma*sqrt(3)));
	fprintf(fp," where cmax is the maximum phase velocity at infinite frequency,\n");
	fprintf(fp," In the current simulation cmax is %8.2f m/s .\n\n",cmax);
	fprintf(fp," DT is the timestep and DH is the grid size.\n\n");
	fprintf(fp," gamma = %2.3f for spatial FD operators of order %d .\n\n", gamma, FDORDER);
	
	dtstab=dhmin/(gamma*sqrt(3.0)*cmax);
	fprintf(fp," In this simulation the stability limit for timestep DT is %e seconds .\n",dtstab);
	fprintf(fp," You have specified DT= %e s.\n", DT);
	
	if (DT>dtstab)
		err(" The simulation will get unstable, choose smaller DT. ");
	else fprintf(fp," The simulation will be stable.\n");
	
	fprintf(fp," \n\n --------------------- CHECK FOR INPUT ERRORS ---------------------\n");

	fprintf(fp," Checking the time of wave propagation. \n");
	
	if (SNAP){
	   fprintf(fp," Checking the snapshot parameters. \n");
	   if ((TSNAP2>TIME) && (MYID==0)){
		   sprintf(errormessage,"TSNAP2 = %e (last snapshot) > Time of wave propagation %e. Choose smaller TSNAP2.",TSNAP2, TIME);	   	
		   err(errormessage); /* if TSNAP2>simulation TIME, snapmerge will generate "additional" snapshots out of nowhere, thus, snapshot files size blow up */
		   }

	   snapoutx=NX/(float)IDX;
	   snapouty=NY/(float)IDY;
	   snapoutz=NZ/(float)IDZ;
	   fprintf(fp," Output of snapshot gridpoints per node (NX/NPROCX/IDX) %8.2f .\n", snapoutx);
	   fprintf(fp," Output of snapshot gridpoints per node (NY/NPROCY/IDY) %8.2f .\n", snapoutz);
	   fprintf(fp," Output of snapshot gridpoints per node (NZ/NPROCZ/IDZ) %8.2f .\n", snapouty);

	   if (snapoutx-(int)snapoutx>0)
		   err("\n\n Ratio NX-NPROCX-IDX must be whole-numbered \n\n");
	   if (snapouty-(int)snapouty>0)
		   err("\n\n Ratio NY-NPROCY-IDY must be whole-numbered \n\n");
	   if (snapoutz-(int)snapoutz>0)
		   err("\n\n Ratio NZ-NPROCZ-IDZ must be whole-numbered \n\n");
	}
	

	if ((SEISMO)&& (MYID==0)){
		fprintf(fp," Checking the number of seismogram samples. \n");
		fprintf(fp,"    Number of timesteps %d.\n", NT);
		fprintf(fp,"    Seismogram sampling rate in timesteps %d.\n", NDT);
		fprintf(fp,"    Number of seismogram output samples %d.\n", NT/NDT-NDTSHIFT);
	
	} 
	 
	if ((SEISMO>0) && (MYID==0)) {
		srec_minx=DX*NX*NPROCX+1, srec_miny=DY*NY*NPROCY+1, srec_minz=DZ*NZ*NPROCZ+1;
		srec_maxx=-1.0, srec_maxy=-1.0, srec_maxz=-1.0;
		fprintf(fp," Checking for receiver position(s) specified in input file.\n");
		fprintf(fp,"    Global grid size in m: %5.2f (x) : %5.2f (y) : %5.2f (z) :.\n",NX*DX*NPROCX,NZ*DZ*NPROCZ,NY*DY*NPROCY);
		if (ABS_TYPE==2){
			if (FREE_SURF==0) fprintf(fp,"    Global grid size in m(-width of abs.boundary) : \n        %5.2f-%5.2f (x in m) : %5.2f-%5.2f (y in m) : %5.2f-%5.2f (z in	m).\n",FW*DX,NX*DX*NPROCX-FW*DX,FW*DZ,NZ*DZ*NPROCZ-FW*DZ,FW*DY,NY*DY*NPROCY-FW*DY);
			if (FREE_SURF==1) fprintf(fp,"    Global grid size in m(-width of abs.boundary) : \n        %5.2f-%5.2f (x in m) : %5.2f-%5.2f (y in m) : %5.2f-%5.2f (z in m).\n",FW*DX,NX*DX*NPROCX-FW*DX,FW*DZ,NZ*DZ*NPROCZ-FW*DZ,DY,NY*DY*NPROCY-FW*DY);
		}

		
		/* find maximum and minimum source positions coordinate ---- from input file*/
		/*for usability reasons, "z" - as commonly used - denotes the depth (vertical direction),
                  however, internally "y" is used for the vertical coordinate,
                  we simply switch the "y" and "z" coordinate as read in the input file,
                  therefore we determine the minimum/maximum position in y-direction by the ZREC1 variable and vice versa.
		  this has to be considered for the receiver line coordinates specified in both the input file and separate source/receiver files*/

		if (READREC==0) {
			if (XREC1>XREC2) {
				srec_maxx=XREC1;
				srec_minx=XREC2;
			}
			else {
				srec_maxx=XREC2;
				srec_minx=XREC1;
			}
			if (YREC1>YREC2) {
				srec_maxz=YREC1;
				srec_minz=YREC2;
			}			
			else {
				srec_maxz=YREC2;
				srec_minz=YREC1;
			}			
			if (ZREC1>ZREC2) {
				srec_maxy=ZREC1;
				srec_miny=ZREC2;
			}
			else {
				srec_maxy=ZREC2;
				srec_miny=ZREC1;
			}
		}
		if (READREC==1) {
			/* find maximum and minimum source positions coordinate ---- from receiver file*/
			for (k=1;k<=ntr;k++){
				/* find maximum source positions coordinate*/
				if ((recpos[1][k]*DX)>srec_maxx) srec_maxx=recpos[1][k]*DX;
				if ((recpos[3][k]*DZ)>srec_maxy) srec_maxy=recpos[3][k]*DZ;
				if ((recpos[2][k]*DY)>srec_maxz) srec_maxz=recpos[2][k]*DY;
				/* find minimum source positions coordinate*/
				if ((recpos[1][k]*DX)<srec_minx) srec_minx=recpos[1][k]*DX;
				if ((recpos[3][k]*DY)<srec_miny) srec_miny=recpos[3][k]*DZ;
				if ((recpos[2][k]*DY)<srec_minz) srec_minz=recpos[2][k]*DY;
			}
			fprintf(fp,"    Number of receiver positions in receiver file %s : %i \n", REC_FILE, ntr);
		}
		
		
		fprintf(fp,"    Minimum receiver position coordinates : %5.2f (x) : %5.2f (y) : %5.2f (z).\n",srec_minx,srec_miny, srec_minz);
		fprintf(fp,"    Maximum receiver position coordinates : %5.2f (x) : %5.2f (y) : %5.2f (z).\n",srec_maxx,srec_maxy, srec_maxz);
		/* checking if receiver coordinate of first receiver in line specified in input-file is inside the global grid */
		if ((((srec_maxx<0.0) || (srec_maxy<0.0)) || ((srec_maxz<0.0))) || (((srec_minx<0.0) || (srec_miny<0.0)) || ((srec_minz<0.0)))) {
			err("\n\n Coordinate of at least one receiver location is outside the global grid. \n\n");
		}
		if (((srec_maxx>NX*DX*NPROCX) || (srec_maxy>NZ*DZ*NPROCZ)) || ((srec_maxz>NY*DY*NPROCY))) {
			err("\n\n Coordinate of at least one receiver location is outside the global grid. \n\n");
		}
		/* checking if receiver coordinate of first receiver in line specified in input-file is outside the Absorbing Boundary  */
		if (ABS_TYPE==2){
			if ((((srec_maxx<(FW*DX)) || (srec_maxy<(FW*DZ))) || ((srec_maxx<(FW*DY))&& !(FREE_SURF==1))) || (((srec_minx<(FW*DX)) || (srec_miny<(FW*DZ)))) || ((srec_minz<(FW*DY)) && !(FREE_SURF==1))) {
				/* this warning appears, when at least a single receiver is located in AB between 0 - FW+DX/DX/DZ ("inner boundary")*/
				warning("\n\n Coordinate of at least one receiver location is inside the Absorbing Boundary (warning 1). \n\n");
			}
			if (((srec_maxx>(NX*DX*NPROCX-FW*DX)) || (srec_maxy>(NZ*DZ*NPROCZ-FW*DZ)) || ((srec_maxz>(NY*DY*NPROCY-FW*DY))))) {
				/* this warning appears, when at least a single receiver is located in AB between NX/NY/NZ-FW+DX/DX/DZ and NX/NY/NZ ("outer boundary")*/
				warning("\n\n Coordinate of at least one receiver location is inside the Absorbing Boundary (warning 2). \n\n");
			}
		}
		
		fprintf(fp," ... complete, receiver position specified in input file are located within the global grid.\n");
	
	}
	
	if ((SRCREC==1)&& (MYID==0)){
		srec_minx=DX*NX*NPROCX+1, srec_miny=DY*NY*NPROCY+1, srec_minz=DZ*NZ*NPROCZ+1;
		srec_maxx=-1.0, srec_maxy=-1.0, srec_maxz=-1.0;
		fprintf(fp," Checking for source position(s) specified in source file. \n");
		fprintf(fp,"    Global grid size in m: %5.2f (x) : %5.2f (y) : %5.2f (z) :.\n",NX*DX*NPROCX,NZ*DZ*NPROCZ,NY*DY*NPROCY);
		if (ABS_TYPE==2){
			if (FREE_SURF==0) fprintf(fp,"    Global grid size in m (-width of abs.boundary) : \n        %5.2f-%5.2f (x in m) : %5.2f-%5.2f (y in m) : %5.2f-%5.2f (z in	m).\n",FW*DX,NX*DX*NPROCX-FW*DX,FW*DZ,NZ*DZ*NPROCZ-FW*DZ,FW*DY,NY*DY*NPROCY-FW*DY);
			if (FREE_SURF==1) fprintf(fp,"    Global grid size in m(-width of abs.boundary) : \n        %5.2f-%5.2f (x in m) : %5.2f-%5.2f (y in m) : %5.2f-%5.2f (z in m).\n",FW*DX,NX*DX*NPROCX-FW*DX,FW*DZ,NZ*DZ*NPROCZ-FW*DZ,DY,NY*DY*NPROCY-FW*DY);
		}
		
		/*only for testing
		fprintf(fp," initial min : %5.2f (x) %5.2f (y) %5.2f (z). \n",srec_minx,srec_miny,srec_minz);
		fprintf(fp," initial max : %5.2f (x) %5.2f (y) %5.2f (z). \n",srec_maxx,srec_maxy,srec_maxz);
		
		fprintf(fp," initial model size  : %i (Nx) %i (Ny) %i (Nz). \n",NX,NY,NZ);
		fprintf(fp," initial model size  : %5.2f (dx) %5.2f (dy) %5.2f (dz). \n",DX,DY,DZ);*/
		for (k=1;k<=nsrc;k++){
			/* find maximum source positions coordinate*/
			if (srcpos[1][k]>srec_maxx) srec_maxx=srcpos[1][k];
			if (srcpos[3][k]>srec_maxy) srec_maxy=srcpos[3][k];
			if (srcpos[2][k]>srec_maxz) srec_maxz=srcpos[2][k];
			/* find minimum source positions coordinate*/
			if (srcpos[1][k]<srec_minx) srec_minx=srcpos[1][k];
			if (srcpos[3][k]<srec_miny) srec_miny=srcpos[3][k];
			if (srcpos[2][k]<srec_minz) srec_minz=srcpos[2][k];
		}
		/*only for testing
		fprintf(fp," new min : %5.2f (x) %5.2f (y) %5.2f (z). \n",srec_minx,srec_miny,srec_minz);
		fprintf(fp," new max : %5.2f (x) %5.2f (y) %5.2f (z). \n",srec_maxx,srec_maxy,srec_maxz);*/
		
		fprintf(fp,"    Number of source positions in source file %s. : %i.\n", SOURCE_FILE, nsrc);
		fprintf(fp,"    Minimum source position coordinates : %5.2f (x) : %5.2f (y) : %5.2f (z).\n",srec_minx,srec_miny, srec_minz);
		fprintf(fp,"    Maximum source position coordinates : %5.2f (x) : %5.2f (y) : %5.2f (z).\n",srec_maxx,srec_maxy, srec_maxz);
		/* checking if receiver coordinate of first receiver in line specified in input-file is inside the global grid */
		if ((((srec_maxx<0.0) || (srec_maxy<0.0)) || ((srec_maxz<0.0))) || (((srec_minx<0.0) || (srec_miny<0.0)) || ((srec_minz<0.0)))) {
			err("\n\n Coordinate of at least one source location is outside the global grid. \n\n");
		}
		if (((srec_maxx>NX*DX*NPROCX) || (srec_maxy>NZ*DZ*NPROCZ)) || ((srec_maxz>NY*DY*NPROCY))) {
			err("\n\n Coordinate of at least one source location is outside the global grid. \n\n");
		}
		/* checking if receiver coordinate of first receiver in line specified in input-file is outside the Absorbing Boundary  */
		if (ABS_TYPE==2){
			if ((((srec_maxx<(FW*DX)) || (srec_maxy<(FW*DZ))) || ((srec_maxx<(FW*DY))&& !(FREE_SURF==1))) || (((srec_minx<(FW*DX)) || (srec_miny<(FW*DZ)))) || ((srec_minz<(FW*DY)) && !(FREE_SURF==1))) {
				/* this warning appears, when at least a single receiver is located in AB between 0 - FW+DX/DX/DZ ("inner boundary")*/
				warning("\n\n Coordinate of at least one source location is inside the Absorbing Boundary (warning 1). \n\n");
			}
			if (((srec_maxx>(NX*DX*NPROCX-FW*DX)) || (srec_maxy>(NZ*DZ*NPROCZ-FW*DZ)) || ((srec_maxz>(NY*DY*NPROCY-FW*DY))))) {
				/* this warning appears, when at least a single receiver is located in AB between NX/NY/NZ-FW+DX/DX/DZ and NX/NY/NZ ("outer boundary")*/
				warning("\n\n Coordinate of at least one source location is inside the Absorbing Boundary (warning 2). \n\n");
			}
		}
			
			fprintf(fp," ...complete, all source position(s) specified in source file are located within the global grid.\n");	
		
		
	}

	
	/* 		if ((NT/NDT-NDTSHIFT)>32000)
	err("\n\n Seismogram samples must be less than 32.0000 \n\n"); */ /*moved here from inside the previous statement for reason below */
	/* checked in writepar.c now. SU and SEG-Y allow 32767 samples, furthermore
	the exist programs allowing 65535 samples and pseudo-SEG-Y formats allowing
	almost arbitrarily long traces. And for binary and textual output the limit is
	arbitrary. */


	
	fprintf(fp,"\n\n ----------------------- ABSORBING BOUNDARY ------------------------\n");
	
	if (ABS_TYPE==2){
	  
	fprintf(fp," Width (FW) of absorbing frame should be at least 30 gridpoints.\n");
	fprintf(fp," You have specified a width of %i gridpoints.\n",FW);
	if ((FW<30)&&(MYID==0)) 
		warning(" Be aware of strong artificial reflections from grid boundaries ! \n");}
		
	if (ABS_TYPE==2){
	fprintf(fp," Width (FW) of CPML-frame should be at least 20 gridpoints.\n");
	fprintf(fp," You have specified a width of %i gridpoints.\n",FW);
	if ((FW<20)&&(MYID==0)) 
		warning(" Be aware of strong artificial reflections from grid boundaries ! \n");}



}