checkfd.c 17 KB
Newer Older
tilman.metz's avatar
tilman.metz committed
1
/*-----------------------------------------------------------------------------------------
2
 * Copyright (C) 2016  For the list of authors, see file AUTHORS.
tilman.metz's avatar
tilman.metz committed
3
 *
Florian Wittkamp's avatar
Florian Wittkamp committed
4
 * This file is part of IFOS.
tilman.metz's avatar
tilman.metz committed
5
 * 
Florian Wittkamp's avatar
Florian Wittkamp committed
6
 * IFOS is free software: you can redistribute it and/or modify
tilman.metz's avatar
tilman.metz committed
7 8 9
 * it under the terms of the GNU General Public License as published by
 * the Free Software Foundation, version 2.0 of the License only.
 * 
Florian Wittkamp's avatar
Florian Wittkamp committed
10
 * IFOS is distributed in the hope that it will be useful,
tilman.metz's avatar
tilman.metz committed
11 12 13 14 15
 * 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
Florian Wittkamp's avatar
Florian Wittkamp committed
16
 * along with IFOS. See file COPYING and/or <http://www.gnu.org/licenses/gpl-2.0.html>.
tilman.metz's avatar
tilman.metz committed
17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32
-----------------------------------------------------------------------------------------*/

/*-------------------------------------------------------------
 *  Check FD-Grid for stability and grid dispersion.
 *  If the stability criterion is not fullfilled the program will
 *  terminate.                   
 *
 *  ----------------------------------------------------------*/

#include <limits.h>
#include "fd.h"

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

	extern float DH, DT, TS;
	extern float XREC1, XREC2, YREC1, YREC2;
33
        extern int NX, NY, MYID, PARAMETERIZATION, FW, L, NT, NDT, ACOUSTIC;
tilman.metz's avatar
tilman.metz committed
34 35 36 37 38 39
	extern int READREC, NPROCX,NPROCY, SRCREC, FREE_SURF;
        extern int SEISMO, SEIS_FORMAT[6];
	extern char SOURCE_FILE[STRING_SIZE], REC_FILE[STRING_SIZE];

	/* local variables */

40
	float  c = 0.0, cmax_p=0.0, cmin_p=1e9, cmax_s=0.0, cmin_s=1e9, fmax, gamma;
tilman.metz's avatar
tilman.metz committed
41
	float  cmax=0.0, cmin=1e9, dtstab, dhstab, cmax_r, cmin_r;
42
	float sumu, sumpi, ws, ts, ppi_ref = 0.0, pu_ref;
tilman.metz's avatar
tilman.metz committed
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
	int nfw=iround(FW/DH);
	int i, j, k, l, ny1=1, nx, ny, nx_min, ny_min;
	float srec_minx=DH*NX*NPROCX+1, srec_miny=DH*NY*NPROCY+1;
	float srec_maxx=-1.0, srec_maxy=-1.0;
	
	
	ws=2.0*PI/TS; /*center frequency of source*/
	
	nx=NX; ny=NY; 

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

	/* find maximum model phase velocity of shear waves at infinite
	      frequency within the whole model */
	if(!ACOUSTIC){
		if (L){  /*viscoelastic*/
			for (i=1+nfw;i<=(nx-nfw);i++){
				for (j=ny1;j<=(ny-nfw);j++){
				
					sumu=0.0;
					for (l=1;l<=L;l++){
						ts=DT/peta[l];
						sumu=sumu+((ws*ws*ts*ts*ptaus[j][i])/(1.0+ws*ws*ts*ts));
					}
					
72
					if(PARAMETERIZATION==1){
73 74 75 76
						pu_ref=prho[j][i]*pu[j][i]*pu[j][i];
                    } else {
						pu_ref=pu[j][i];
                    }
tilman.metz's avatar
tilman.metz committed
77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94
						
					
					/* minimum phase velocity of shear waves */
					c=sqrt(pu_ref/(prho[j][i]*(1.0+sumu)));
									
					if (cmin_s>c) cmin_s=c;
					
					
					/* maximum phase velocity of shear waves */
					c=sqrt(pu_ref*(1.0+L*ptaus[j][i])/(prho[j][i]*(1.0+sumu)));
					
					if (cmax_s<c) cmax_s=c;
				}
			}
		} else { /*L=0, elastic*/
			for (i=1+nfw;i<=(nx-nfw);i++){
				for (j=ny1;j<=(ny-nfw);j++){

95
					if(PARAMETERIZATION==3){
96 97
						c=sqrt(pu[j][i]/prho[j][i]);
                    } else {
tilman.metz's avatar
tilman.metz committed
98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120
						c=pu[j][i];
					}

					if (cmax_s<c) cmax_s=c;
					if (cmin_s>c) cmin_s=c;
				}
			}
		}
	}
	
	/* find maximum model phase velocity of P-waves at infinite
		 frequency within the whole model */
	if (L){  /*viscoelastic*/
		for (i=1+nfw;i<=(nx-nfw);i++){
			for (j=ny1;j<=(ny-nfw);j++){
			        
				
				sumpi=0.0;
				for (l=1;l<=L;l++){
					ts=DT/peta[l];
					sumpi=sumpi+((ws*ws*ts*ts*ptaup[j][i])/(1.0+ws*ws*ts*ts));
				}
				
121
				if(PARAMETERIZATION==1){
tilman.metz's avatar
tilman.metz committed
122
					ppi_ref=prho[j][i]*ppi[j][i]*ppi[j][i];}
123
				if(PARAMETERIZATION==3){
tilman.metz's avatar
tilman.metz committed
124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142
					ppi_ref=ppi[j][i]+2*pu[j][i];}
					
				/* minimum phase velocity of P waves */
				c=sqrt(ppi_ref/(prho[j][i]*(1.0+sumpi)));
				
				if (cmin_p>c) cmin_p=c;
				
				
				/* maximum phase velocity of shear waves */
				c=sqrt(ppi_ref*(1.0+L*ptaup[j][i])/(prho[j][i]*(1.0+sumpi)));
				
				if (cmax_p<c) cmax_p=c;
			}
		}
	}
	else{ /*L=0, elastic*/
		for (i=1+nfw;i<=(nx-nfw);i++){
			for (j=ny1;j<=(ny-nfw);j++){

143
				if(PARAMETERIZATION==3){
tilman.metz's avatar
tilman.metz committed
144 145
					c=sqrt((ppi[j][i]+2.0*pu[j][i])/prho[j][i]);}

146
				if(PARAMETERIZATION==1){
tilman.metz's avatar
tilman.metz committed
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
					c=ppi[j][i];}

				if (cmax_p<c) cmax_p=c;
				if (cmin_p>c) cmin_p=c;
			}
		}

	}


	if (MYID==0){
		fprintf(fp,"\n\n\n **Message from checkfd (printed by PE %d):\n",MYID);
		if(!ACOUSTIC){
		fprintf(fp," Minimum and maximum P-wave and S-wave velocities within subvolumes: \n ");
		fprintf(fp," MYID\t Vp_min(f=fc) \t Vp_max(f=inf) \t Vs_min(f=fc) \t Vsmax(f=inf) \n");}
		else{
		fprintf(fp," Minimum and maximum P-wave velocities within subvolumes: \n ");
		fprintf(fp," MYID\t Vp_min(f=fc) \t Vp_max(f=inf) \n");}
	}
	MPI_Barrier(MPI_COMM_WORLD);
	if(!ACOUSTIC)
	fprintf(fp," %d \t %e \t %e \t %e \t %e \n", MYID, cmin_p, cmax_p, cmin_s, cmax_s);
	else
	fprintf(fp," %d \t %e \t %e \n", MYID, cmin_p, cmax_p);

	if(ACOUSTIC==0){
		if (cmax_s>cmax_p) cmax=cmax_s; 
		else cmax=cmax_p;
		if (cmin_s<cmin_p) cmin=cmin_s; 
		else cmin=cmin_p;}
	else if(ACOUSTIC==1){
		cmax=cmax_p;
		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;	

	fmax=2.0/TS;
	dhstab = (cmin/(hc[0]*fmax));
	gamma = fabs(hc[1]) + fabs(hc[2]) + fabs(hc[3]) + fabs(hc[4]) + fabs(hc[5]) + fabs(hc[6]);
	dtstab = DH/(sqrt(2)*gamma*cmax);
	/*dtstab=DH/(sqrt(2.0)*cmax);*/

	/* find global minimum for NX and NY */
	MPI_Allreduce(&NX,&nx_min,1,MPI_INT,MPI_MIN,MPI_COMM_WORLD);
	MPI_Allreduce(&NY,&ny_min,1,MPI_INT,MPI_MIN,MPI_COMM_WORLD);
	

	if (MYID == 0) {

	fprintf(fp," Global values for entire model: \n");
	if(!ACOUSTIC)
		fprintf(fp," Vp_max= %e m/s \t Vs_min=%e m/s \n\n", cmax,cmin);
	else
		fprintf(fp," Vp_max= %e m/s \t Vp_min=%e m/s \n\n", cmax,cmin);
	fprintf(fp,"\n\n ------------------ CHECK FOR GRID DISPERSION --------------------\n");
	fprintf(fp," To satisfactorily limit grid dispersion the number of gridpoints \n");
	if(!ACOUSTIC)
		fprintf(fp," per minimum wavelength (of S-waves) should be 6 (better more).\n");
	else
		fprintf(fp," per minimum wavelength (of P-waves) should be 6 (better more).\n");
	fprintf(fp," Here the minimum wavelength is assumed to be minimum model phase velocity \n");
	if(!ACOUSTIC)
		fprintf(fp," (of S-waves) at maximum frequency of the source\n");
	else
		fprintf(fp," (of P-waves) at maximum frequency of the source\n");
	fprintf(fp," devided by maximum frequency of the source.\n");
	fprintf(fp," Maximum frequency of the source is approximately %8.2f Hz\n",2.0/TS);
	if(!ACOUSTIC)
		fprintf(fp," The minimum wavelength (of S-waves) in the following simulation will\n");
	else
		fprintf(fp," The minimum wavelength (of P-waves) in the following simulation will\n");
	fprintf(fp," be %e meter.\n", cmin/fmax);
	fprintf(fp," Thus, the recommended value for DH is %e meter.\n", dhstab);
	fprintf(fp," You have specified DH= %e meter.\n\n", DH);
	if (DH>dhstab)
		warning(" Grid dispersion will influence wave propagation, choose smaller grid spacing (DH).");
	
	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/(sqrt(2)*gamma),\n\n");
	fprintf(fp," where cmax is the maximum phase velocity at infinite frequency\n");
	fprintf(fp," and gamma = sum(|FD coeff.|)\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," 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)
241
		declare_error(" The simulation will get unstable, choose smaller DT. ");
tilman.metz's avatar
tilman.metz committed
242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258
	else fprintf(fp," The simulation will be stable.\n");


        fprintf(fp," \n\n --------------------- CHECK FOR INPUT ERRORS ---------------------\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);

		/* SU and SEG-Y allow 32767 samples, furthermore the exist programs allow for 65535 samples
		and pseudo-SEG-Y formats allow foralmost arbitrarily long traces. 
		For binary and textual output the limit is arbitrary. USHRT_MAX is the limit of an unsigned short specified in limits.h */

		if ((SEIS_FORMAT[0]==1) && (NT/NDT)>( USHRT_MAX )) {
			fprintf(fp," Maximum allowed number of samples per trace in SU format: %d \n", USHRT_MAX );
259
			declare_error(" Sorry. Too many samples per receiver! \n");
tilman.metz's avatar
tilman.metz committed
260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279
		}
	}



	if ((SEISMO>0) && (MYID==0)) {
		srec_minx=DH*NX*NPROCX+1, srec_miny=DH*NY*NPROCY+1;
		srec_maxx=-1.0, srec_maxy=-1.0;
		fprintf(fp,"\n Checking for receiver position(s) specified in input file.\n");
		fprintf(fp,"    Global grid size in m: %5.2f (x) : %5.2f (y) \n",NX*DH*NPROCX,NY*DH*NPROCY);
		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) \n",FW*DH,NX*DH*NPROCX-FW*DH,FW*DH,NY*DH*NPROCY-FW*DH);
		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) \n",FW*DH,NX*DH*NPROCX-FW*DH,DH,NY*DH*NPROCY-FW*DH);

                        /* 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*/

280
		if (READREC==0 || READREC==2) {
tilman.metz's avatar
tilman.metz committed
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
			if (XREC1>XREC2) {
				srec_maxx=XREC1;
				srec_minx=XREC2;
			}
			else {
				srec_maxx=XREC2;
				srec_minx=XREC1;
			}
			if (YREC1>YREC2) {
				srec_maxy=YREC1;
				srec_miny=YREC2;
			}
			else {
				srec_maxy=YREC2;
				srec_miny=YREC1;
			}
			fprintf(fp,"    Number of receiver positions in input file : %i \n", ntr);
		}
		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]*DH)>srec_maxx) srec_maxx=recpos[1][k]*DH;
				if ((recpos[2][k]*DH)>srec_maxy) srec_maxy=recpos[2][k]*DH;
				/* find minimum source positions coordinate*/
				if ((recpos[1][k]*DH)<srec_minx) srec_minx=recpos[1][k]*DH;
				if ((recpos[2][k]*DH)<srec_miny) srec_miny=recpos[2][k]*DH;
			}
			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) \n",srec_minx,srec_miny);
		fprintf(fp,"    Maximum receiver position coordinates : %5.2f (x) : %5.2f (y) \n",srec_maxx,srec_maxy);
		/* 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_minx<0.0) || (srec_miny<0.0))) {
318
			declare_error("\n\n Coordinate of at least one receiver location is outside the global grid. \n\n");
tilman.metz's avatar
tilman.metz committed
319 320
		}
		if ((srec_maxx>NX*DH*NPROCX) || (srec_maxy>NY*DH*NPROCY)) {
321
			declare_error("\n\n Coordinate of at least one receiver location is outside the global grid. \n\n");
tilman.metz's avatar
tilman.metz committed
322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371
		}
		/* checking if receiver coordinate of first receiver in line specified in input-file is outside the Absorbing Boundary  */
		if ((srec_maxx<(FW*DH)) || (srec_minx<(FW*DH))) {
			/* 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, left boundary). \n\n");
		}
		if (srec_maxx>(NX*DH*NPROCX-FW*DH)) {
			/* 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, right boundary). \n\n");
		}
		if (srec_maxy>(NY*DH*NPROCY-FW*DH)) {
			/* 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 3, lower boundary). \n\n");
		}
		if ((srec_miny<(FW*DH)) && !(FREE_SURF)) {
			/* 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 4, top boundary). \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=DH*NX*NPROCX+1, srec_miny=DH*NY*NPROCY+1;
		srec_maxx=-1.0, srec_maxy=-1.0;
		fprintf(fp,"\n Checking for source position(s) specified in source file. \n");
		fprintf(fp,"    Global grid size in m: %5.2f (x) : %5.2f (y) \n",NX*DH*NPROCX,NY*DH*NPROCY);
		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) \n",FW*DH,NX*DH*NPROCX-FW*DH,FW*DH,NY*DH*NPROCY-FW*DH);
		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) \n",FW*DH,NX*DH*NPROCX-FW*DH,DH,NY*DH*NPROCY-FW*DH);


		for (k=1;k<=nsrc;k++){
			/* find maximum source positions coordinate*/
			if (srcpos[1][k]>srec_maxx) srec_maxx=srcpos[1][k];
			if (srcpos[2][k]>srec_maxy) srec_maxy=srcpos[2][k];
			/* find minimum source positions coordinate*/
			if (srcpos[1][k]<srec_minx) srec_minx=srcpos[1][k];
			if (srcpos[2][k]<srec_miny) srec_miny=srcpos[2][k];
		}

		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) \n",srec_minx,srec_miny);
		fprintf(fp,"    Maximum source position coordinates : %5.2f (x) : %5.2f (y) \n",srec_maxx,srec_maxy);
		/* 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_minx<0.0) || (srec_miny<0.0))) {
372
			declare_error("\n\n Coordinate of at least one source location is outside the global grid. \n\n");
tilman.metz's avatar
tilman.metz committed
373 374
		}
		if ((srec_maxx>NX*DH*NPROCX) || (srec_maxy>NY*DH*NPROCY)) {
375
			declare_error("\n\n Coordinate of at least one source location is outside the global grid. \n\n");
tilman.metz's avatar
tilman.metz committed
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
		}
		/* checking if receiver coordinate of first receiver in line specified in input-file is outside the Absorbing Boundary  */
		if ((srec_maxx<(FW*DH)) || (srec_minx<(FW*DH))) {
			/* 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, left boundary). \n\n");
		}
		if (srec_maxx>(NX*DH*NPROCX-FW*DH)) {
			/* 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, right boundary). \n\n");
		}
		if (srec_maxy>(NY*DH*NPROCY-FW*DH)) {
			/* 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 3, lower boundary). \n\n");
		}
		if ((srec_miny<(FW*DH)) && !(FREE_SURF)) {
			/* 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 4, top boundary). \n\n");
		}
		fprintf(fp," ...complete, all source position(s) specified in source file are located within the global grid.\n");

	}





	fprintf(fp,"\n\n ----------------------- ABSORBING BOUNDARY ------------------------\n");
        if((FW>nx_min)||(FW>ny_min)){
404
	  declare_error(" The width of the absorbing boundary is larger than one computational domain. Choose smaller FW or use less CPUs.");
tilman.metz's avatar
tilman.metz committed
405 406 407 408 409 410 411 412 413 414 415 416 417
	}

	fprintf(fp," Width (FW) of absorbing frame should be at least 10 gridpoints.\n");
	fprintf(fp," You have specified a width of %d gridpoints.\n",FW);
	if (FW<10) 
		warning(" Be aware of artificial reflections from grid boundaries ! \n");

	}
	
	MPI_Barrier(MPI_COMM_WORLD);
	
}