update_s_visc_PML.c 25.3 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 33 34 35 36 37 38 39 40 41 42 43 44 45
-----------------------------------------------------------------------------------------*/

/*------------------------------------------------------------------------
 *   updating stress components at gridpoints [nx1...nx2][ny1...ny2]
 *   by a staggered grid finite difference scheme of arbitrary (FDORDER) order accuracy in space
 *   and second order accuracy in time
 *
 *  ----------------------------------------------------------------------*/

#include "fd.h"

void update_s_visc_PML(int nx1, int nx2, int ny1, int ny2,
	float **  vx, float **   vy, float **  ux, float **   uy, float **  uxy, float **   uyx, float **   sxx, float **   syy,
	float **   sxy, float ** pi, float ** u, float ** uipjp, float **rho, float *hc, int infoout,
	float ***r, float ***p, float ***q, float **fipjp, float **f, float **g, float *bip, float *bjm, float *cip, float *cjm, float ***d, float ***e, float ***dip, 
      float * K_x, float * a_x, float * b_x, float * K_x_half, float * a_x_half, float * b_x_half,
      float * K_y, float * a_y, float * b_y, float * K_y_half, float * a_y_half, float * b_y_half,
      float ** psi_vxx, float ** psi_vyy, float ** psi_vxy, float ** psi_vyx){



	int i,j, m, fdoh, h, h1, l;
	float  vxx, vyy, vxy, vyx;
	float  dhi, dthalbe;	
	extern float DT, DH;
	extern int MYID, FDORDER, FW, L;
        extern int FREE_SURF, BOUNDARY;
	extern int NPROCX, NPROCY, POS[3];
	extern FILE *FP;
46
	double time1 = 0.0, time2;
tilman.metz's avatar
tilman.metz committed
47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68
	
	float sumr=0.0, sump=0.0, sumq=0.0;
	
	

	/*dhi = DT/DH;*/
	dhi=1.0/DH;
	fdoh = FDORDER/2;
	dthalbe = DT/2.0;

	
	if (infoout && (MYID==0)){
		time1=MPI_Wtime();
		fprintf(FP,"\n **Message from update_s (printed by PE %d):\n",MYID);
		fprintf(FP," Updating stress components ...");
	}
	


	switch (FDORDER){

	case 2:
69 70 71 72 73 74 75 76 77 78 79 80
            for (j=ny1;j<=ny2;j++){
                for (i=nx1;i<=nx2;i++){
                    vxx = (  hc[1]*(vx[j][i]  -vx[j][i-1]))*dhi;
                    
                    vyx = (  hc[1]*(vy[j][i+1]-vy[j][i]))*dhi;
                    
                    vxy = (  hc[1]*(vx[j+1][i]-vx[j][i]))*dhi;
                    
                    vyy = (  hc[1]*(vy[j][i]  -vy[j-1][i]))*dhi;
                    
                    /* left boundary */
                    if((!BOUNDARY) && (POS[1]==0) && (i<=FW)){
tilman.metz's avatar
tilman.metz committed
81 82 83
                        
                        psi_vxx[j][i] = b_x[i] * psi_vxx[j][i] + a_x[i] * vxx;
                        vxx = vxx / K_x[i] + psi_vxx[j][i];
84
                        
tilman.metz's avatar
tilman.metz committed
85
                        psi_vyx[j][i] = b_x_half[i] * psi_vyx[j][i] + a_x_half[i] * vyx;
86 87 88 89 90 91
                        vyx = vyx / K_x_half[i] + psi_vyx[j][i];
                    }
                    
                    /* right boundary */
                    if((!BOUNDARY) && (POS[1]==NPROCX-1) && (i>=nx2-FW+1)){
                        
tilman.metz's avatar
tilman.metz committed
92 93 94 95
                        h1 = (i-nx2+2*FW);
                        h = i;
                        
                        psi_vxx[j][h1] = b_x[h1] * psi_vxx[j][h1] + a_x[h1] * vxx;
96 97
                        vxx = vxx / K_x[h1] + psi_vxx[j][h1];
                        
tilman.metz's avatar
tilman.metz committed
98
                        /*psi_vyx[j][h] = b_x_half[h] * psi_vyx[j][h] + a_x_half[h] * vyx;
99
                         vyx = vyx / K_x_half[h] + psi_vyx[j][h];*/
tilman.metz's avatar
tilman.metz committed
100 101
                        
                        psi_vyx[j][h1] = b_x_half[h1] * psi_vyx[j][h1] + a_x_half[h1] * vyx;
102 103 104 105 106 107 108 109
                        vyx = vyx / K_x_half[h1] + psi_vyx[j][h1];
                        
                    }
                    
                    /* top boundary */
                    if((POS[2]==0) && (!(FREE_SURF)) && (j<=FW)){
                        
                        psi_vyy[j][i] = b_y[j] * psi_vyy[j][i] + a_y[j] * vyy;
tilman.metz's avatar
tilman.metz committed
110
                        psi_vxy[j][i] = b_y_half[j] * psi_vxy[j][i] + a_y_half[j] * vxy;
111
                        
tilman.metz's avatar
tilman.metz committed
112 113
                        vyy = vyy / K_y[j] + psi_vyy[j][i];
                        vxy = vxy / K_y_half[j] + psi_vxy[j][i];
114 115 116 117 118 119 120
                        
                    }
                    
                    /* bottom boundary */
                    if((POS[2]==NPROCY-1) && (j>=ny2-FW+1)){
                        
                        h1 = (j-ny2+2*FW);
tilman.metz's avatar
tilman.metz committed
121
                        h = j;
122 123
                        
                        psi_vyy[h1][i] = b_y[h1] * psi_vyy[h1][i] + a_y[h1] * vyy;
tilman.metz's avatar
tilman.metz committed
124 125 126
                        vyy = vyy / K_y[h1] + psi_vyy[h1][i];
                        
                        /*psi_vxy[j][i] = b_y_half[j] * psi_vxy[j][i] + a_y_half[j] * vxy;
127
                         vxy = vxy / K_y_half[j] + psi_vxy[j][i];*/
tilman.metz's avatar
tilman.metz committed
128 129
                        
                        psi_vxy[h1][i] = b_y_half[h1] * psi_vxy[h1][i] + a_y_half[h1] * vxy;
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
                        vxy = vxy / K_y_half[h1] + psi_vxy[h1][i];
                        
                    }
                    
                    /* computing sums of the old memory variables */
                    sumr=sump=sumq=0.0;
                    for (l=1;l<=L;l++){
                        sumr+=r[j][i][l];
                        sump+=p[j][i][l];
                        sumq+=q[j][i][l];
                    }
                    
                    
                    /* updating components of the stress tensor, partially */
                    sxy[j][i] += (fipjp[j][i]*(vxy+vyx))+(dthalbe*sumr);
                    sxx[j][i] += (g[j][i]*(vxx+vyy))-(2.0*f[j][i]*vyy)+(dthalbe*sump);
                    syy[j][i] += (g[j][i]*(vxx+vyy))-(2.0*f[j][i]*vxx)+(dthalbe*sumq);
                    
                    uxy[j][i] = ((fipjp[j][i]/DT)*(vxy+vyx))+(0.5*sumr);
                    ux[j][i] = ((g[j][i]/DT)*(vxx+vyy))-((2.0*f[j][i]/DT)*vyy)+(0.5*sump);
                    uy[j][i] = ((g[j][i]/DT)*(vxx+vyy))-((2.0*f[j][i]/DT)*vxx)+(0.5*sumq);
                    
                    
                    /* now updating the memory-variables and sum them up*/
                    sumr=sump=sumq=0.0;
                    for (l=1;l<=L;l++){
                        r[j][i][l] = bip[l]*(r[j][i][l]*cip[l]-(dip[j][i][l]*(vxy+vyx)));
                        p[j][i][l] = bjm[l]*(p[j][i][l]*cjm[l]-(e[j][i][l]*(vxx+vyy))+(2.0*d[j][i][l]*vyy));
                        q[j][i][l] = bjm[l]*(q[j][i][l]*cjm[l]-(e[j][i][l]*(vxx+vyy))+(2.0*d[j][i][l]*vxx));
                        sumr += r[j][i][l];
                        sump += p[j][i][l];
                        sumq += q[j][i][l];
                    }
                    
                    
                    /* and now the components of the stress tensor are
                     completely updated */
                    sxy[j][i]+=(dthalbe*sumr);
                    sxx[j][i]+=(dthalbe*sump);
                    syy[j][i]+=(dthalbe*sumq);
                    
                    uxy[j][i]+=(0.5*sumr);
                    ux[j][i]+=(0.5*sump);
                    uy[j][i]+=(0.5*sumq);
                }
            }
            break;
            
        case 4:
tilman.metz's avatar
tilman.metz committed
179 180
		for (j=ny1;j<=ny2;j++){
			for (i=nx1;i<=nx2;i++){
181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199
                vxx = (  hc[1]*(vx[j][i]  -vx[j][i-1])
                       + hc[2]*(vx[j][i+1]-vx[j][i-2]))*dhi;
                
                vyx = (  hc[1]*(vy[j][i+1]-vy[j][i])
                       + hc[2]*(vy[j][i+2]-vy[j][i-1]))*dhi;
                
                vxy = (  hc[1]*(vx[j+1][i]-vx[j][i])
                       + hc[2]*(vx[j+2][i]-vx[j-1][i]))*dhi;
                
                vyy = (  hc[1]*(vy[j][i]  -vy[j-1][i])
                       + hc[2]*(vy[j+1][i]-vy[j-2][i]))*dhi;
                
                /* left boundary */
                if((!BOUNDARY) && (POS[1]==0) && (i<=FW)){
                    
                    psi_vxx[j][i] = b_x[i] * psi_vxx[j][i] + a_x[i] * vxx;
                    vxx = vxx / K_x[i] + psi_vxx[j][i];
                    
                    psi_vyx[j][i] = b_x_half[i] * psi_vyx[j][i] + a_x_half[i] * vyx;
tilman.metz's avatar
tilman.metz committed
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
                        vyx = vyx / K_x_half[i] + psi_vyx[j][i];                 
         }

        /* right boundary */                                         
        if((!BOUNDARY) && (POS[1]==NPROCX-1) && (i>=nx2-FW+1)){
		
                        h1 = (i-nx2+2*FW);
                        h = i;
                        
                        psi_vxx[j][h1] = b_x[h1] * psi_vxx[j][h1] + a_x[h1] * vxx;
                        vxx = vxx / K_x[h1] + psi_vxx[j][h1]; 

                        /*psi_vyx[j][h] = b_x_half[h] * psi_vyx[j][h] + a_x_half[h] * vyx;
                        vyx = vyx / K_x_half[h] + psi_vyx[j][h];*/
                        
                        psi_vyx[j][h1] = b_x_half[h1] * psi_vyx[j][h1] + a_x_half[h1] * vyx;
			vyx = vyx / K_x_half[h1] + psi_vyx[j][h1];
                                           
         }

	  /* top boundary */                                         
        if((POS[2]==0) && (!(FREE_SURF)) && (j<=FW)){
                                                
                        psi_vyy[j][i] = b_y[j] * psi_vyy[j][i] + a_y[j] * vyy;                                            
                        psi_vxy[j][i] = b_y_half[j] * psi_vxy[j][i] + a_y_half[j] * vxy;
                     
                        vyy = vyy / K_y[j] + psi_vyy[j][i];
                        vxy = vxy / K_y_half[j] + psi_vxy[j][i];

        }
	
	  /* bottom boundary */                                         
        if((POS[2]==NPROCY-1) && (j>=ny2-FW+1)){

                        h1 = (j-ny2+2*FW);                                        
                        h = j;
                                                
                        psi_vyy[h1][i] = b_y[h1] * psi_vyy[h1][i] + a_y[h1] * vyy;                                            
                        vyy = vyy / K_y[h1] + psi_vyy[h1][i];
                        
                        /*psi_vxy[j][i] = b_y_half[j] * psi_vxy[j][i] + a_y_half[j] * vxy;
                        vxy = vxy / K_y_half[j] + psi_vxy[j][i];*/
                        
                        psi_vxy[h1][i] = b_y_half[h1] * psi_vxy[h1][i] + a_y_half[h1] * vxy;
			vxy = vxy / K_y_half[h1] + psi_vxy[h1][i];
        
        }

	/* computing sums of the old memory variables */
			sumr=sump=sumq=0.0;
			for (l=1;l<=L;l++){
				sumr+=r[j][i][l];
				sump+=p[j][i][l];
				sumq+=q[j][i][l];
			}


257
            /* updating components of the stress tensor, partially */
tilman.metz's avatar
tilman.metz committed
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 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 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
			sxy[j][i] += (fipjp[j][i]*(vxy+vyx))+(dthalbe*sumr);
			sxx[j][i] += (g[j][i]*(vxx+vyy))-(2.0*f[j][i]*vyy)+(dthalbe*sump);
			syy[j][i] += (g[j][i]*(vxx+vyy))-(2.0*f[j][i]*vxx)+(dthalbe*sumq);
			
			uxy[j][i] = ((fipjp[j][i]/DT)*(vxy+vyx))+(0.5*sumr);
			ux[j][i] = ((g[j][i]/DT)*(vxx+vyy))-((2.0*f[j][i]/DT)*vyy)+(0.5*sump);
			uy[j][i] = ((g[j][i]/DT)*(vxx+vyy))-((2.0*f[j][i]/DT)*vxx)+(0.5*sumq);
			
			
			/* now updating the memory-variables and sum them up*/
			sumr=sump=sumq=0.0;
			for (l=1;l<=L;l++){
				r[j][i][l] = bip[l]*(r[j][i][l]*cip[l]-(dip[j][i][l]*(vxy+vyx)));
				p[j][i][l] = bjm[l]*(p[j][i][l]*cjm[l]-(e[j][i][l]*(vxx+vyy))+(2.0*d[j][i][l]*vyy));
				q[j][i][l] = bjm[l]*(q[j][i][l]*cjm[l]-(e[j][i][l]*(vxx+vyy))+(2.0*d[j][i][l]*vxx));
				sumr += r[j][i][l];
				sump += p[j][i][l];
				sumq += q[j][i][l];
			}

			/* and now the components of the stress tensor are
			   completely updated */
			sxy[j][i]+=(dthalbe*sumr);
			sxx[j][i]+=(dthalbe*sump);
			syy[j][i]+=(dthalbe*sumq);
			
			uxy[j][i]+=(0.5*sumr);
			ux[j][i]+=(0.5*sump);
			uy[j][i]+=(0.5*sumq);
			}
		}
		break;

	case 6:
		for (j=ny1;j<=ny2;j++){
			for (i=nx1;i<=nx2;i++){
			vxx = (  hc[1]*(vx[j][i]  -vx[j][i-1])
				       + hc[2]*(vx[j][i+1]-vx[j][i-2])
				       + hc[3]*(vx[j][i+2]-vx[j][i-3]))*dhi;
			
			vyx = (  hc[1]*(vy[j][i+1]-vy[j][i])
				       + hc[2]*(vy[j][i+2]-vy[j][i-1])
				       + hc[3]*(vy[j][i+3]-vy[j][i-2]))*dhi;

                        vxy = (  hc[1]*(vx[j+1][i]-vx[j][i])
				       + hc[2]*(vx[j+2][i]-vx[j-1][i])
				       + hc[3]*(vx[j+3][i]-vx[j-2][i]))*dhi;

                        vyy = (  hc[1]*(vy[j][i]  -vy[j-1][i])
				       + hc[2]*(vy[j+1][i]-vy[j-2][i])
				       + hc[3]*(vy[j+2][i]-vy[j-3][i]))*dhi; 

        /* left boundary */                                         
        if((!BOUNDARY) && (POS[1]==0) && (i<=FW)){
                        
                        psi_vxx[j][i] = b_x[i] * psi_vxx[j][i] + a_x[i] * vxx;
                        vxx = vxx / K_x[i] + psi_vxx[j][i];

                        psi_vyx[j][i] = b_x_half[i] * psi_vyx[j][i] + a_x_half[i] * vyx;
                        vyx = vyx / K_x_half[i] + psi_vyx[j][i];                 
         }

        /* right boundary */                                         
        if((!BOUNDARY) && (POS[1]==NPROCX-1) && (i>=nx2-FW+1)){
		
                        h1 = (i-nx2+2*FW);
                        h = i;
                        
                        psi_vxx[j][h1] = b_x[h1] * psi_vxx[j][h1] + a_x[h1] * vxx;
                        vxx = vxx / K_x[h1] + psi_vxx[j][h1]; 

                        /*psi_vyx[j][h] = b_x_half[h] * psi_vyx[j][h] + a_x_half[h] * vyx;
                        vyx = vyx / K_x_half[h] + psi_vyx[j][h];*/
                        
                        psi_vyx[j][h1] = b_x_half[h1] * psi_vyx[j][h1] + a_x_half[h1] * vyx;
			vyx = vyx / K_x_half[h1] + psi_vyx[j][h1];
                                           
         }

	  /* top boundary */                                         
        if((POS[2]==0) && (!(FREE_SURF)) && (j<=FW)){
                                                
                        psi_vyy[j][i] = b_y[j] * psi_vyy[j][i] + a_y[j] * vyy;                                            
                        psi_vxy[j][i] = b_y_half[j] * psi_vxy[j][i] + a_y_half[j] * vxy;
                     
                        vyy = vyy / K_y[j] + psi_vyy[j][i];
                        vxy = vxy / K_y_half[j] + psi_vxy[j][i];

        }
	
	  /* bottom boundary */                                         
        if((POS[2]==NPROCY-1) && (j>=ny2-FW+1)){

                        h1 = (j-ny2+2*FW);                                        
                        h = j;
                                                
                        psi_vyy[h1][i] = b_y[h1] * psi_vyy[h1][i] + a_y[h1] * vyy;                                            
                        vyy = vyy / K_y[h1] + psi_vyy[h1][i];
                        
                        /*psi_vxy[j][i] = b_y_half[j] * psi_vxy[j][i] + a_y_half[j] * vxy;
                        vxy = vxy / K_y_half[j] + psi_vxy[j][i];*/
                        
                        psi_vxy[h1][i] = b_y_half[h1] * psi_vxy[h1][i] + a_y_half[h1] * vxy;
			vxy = vxy / K_y_half[h1] + psi_vxy[h1][i];
        
        }

	/* computing sums of the old memory variables */
			sumr=sump=sumq=0.0;
			for (l=1;l<=L;l++){
				sumr+=r[j][i][l];
				sump+=p[j][i][l];
				sumq+=q[j][i][l];
			}


                        /* updating components of the stress tensor, partially */
			sxy[j][i] += (fipjp[j][i]*(vxy+vyx))+(dthalbe*sumr);
			sxx[j][i] += (g[j][i]*(vxx+vyy))-(2.0*f[j][i]*vyy)+(dthalbe*sump);
			syy[j][i] += (g[j][i]*(vxx+vyy))-(2.0*f[j][i]*vxx)+(dthalbe*sumq);
			
			uxy[j][i] = ((fipjp[j][i]/DT)*(vxy+vyx))+(0.5*sumr);
			ux[j][i] = ((g[j][i]/DT)*(vxx+vyy))-((2.0*f[j][i]/DT)*vyy)+(0.5*sump);
			uy[j][i] = ((g[j][i]/DT)*(vxx+vyy))-((2.0*f[j][i]/DT)*vxx)+(0.5*sumq);
			
			
			/* now updating the memory-variables and sum them up*/
			sumr=sump=sumq=0.0;
			for (l=1;l<=L;l++){
				r[j][i][l] = bip[l]*(r[j][i][l]*cip[l]-(dip[j][i][l]*(vxy+vyx)));
				p[j][i][l] = bjm[l]*(p[j][i][l]*cjm[l]-(e[j][i][l]*(vxx+vyy))+(2.0*d[j][i][l]*vyy));
				q[j][i][l] = bjm[l]*(q[j][i][l]*cjm[l]-(e[j][i][l]*(vxx+vyy))+(2.0*d[j][i][l]*vxx));
				sumr += r[j][i][l];
				sump += p[j][i][l];
				sumq += q[j][i][l];
			}

			/* and now the components of the stress tensor are
			   completely updated */
			sxy[j][i]+=(dthalbe*sumr);
			sxx[j][i]+=(dthalbe*sump);
			syy[j][i]+=(dthalbe*sumq);
			
			uxy[j][i]+=(0.5*sumr);
			ux[j][i]+=(0.5*sump);
			uy[j][i]+=(0.5*sumq);
			}
		}
		break;

	case 8:

    for (j=ny1;j<=ny2;j++){
	for (i=nx1;i<=nx2;i++){

			vxx = (  hc[1]*(vx[j][i]  -vx[j][i-1])
				       + hc[2]*(vx[j][i+1]-vx[j][i-2])
				       + hc[3]*(vx[j][i+2]-vx[j][i-3])
				       + hc[4]*(vx[j][i+3]-vx[j][i-4]))*dhi;
			
			vyx = (  hc[1]*(vy[j][i+1]-vy[j][i])
				       + hc[2]*(vy[j][i+2]-vy[j][i-1])
				       + hc[3]*(vy[j][i+3]-vy[j][i-2])
				       + hc[4]*(vy[j][i+4]-vy[j][i-3]))*dhi;

                        vxy = (  hc[1]*(vx[j+1][i]-vx[j][i])
				       + hc[2]*(vx[j+2][i]-vx[j-1][i])
				       + hc[3]*(vx[j+3][i]-vx[j-2][i])
				       + hc[4]*(vx[j+4][i]-vx[j-3][i]))*dhi;

                        vyy = (  hc[1]*(vy[j][i]  -vy[j-1][i])
				       + hc[2]*(vy[j+1][i]-vy[j-2][i])
				       + hc[3]*(vy[j+2][i]-vy[j-3][i])
				       + hc[4]*(vy[j+3][i]-vy[j-4][i]))*dhi; 

        /* left boundary */                                         
        if((!BOUNDARY) && (POS[1]==0) && (i<=FW)){
                        
                        psi_vxx[j][i] = b_x[i] * psi_vxx[j][i] + a_x[i] * vxx;
                        vxx = vxx / K_x[i] + psi_vxx[j][i];

                        psi_vyx[j][i] = b_x_half[i] * psi_vyx[j][i] + a_x_half[i] * vyx;
                        vyx = vyx / K_x_half[i] + psi_vyx[j][i];                 
         }

        /* right boundary */                                         
        if((!BOUNDARY) && (POS[1]==NPROCX-1) && (i>=nx2-FW+1)){
		
                        h1 = (i-nx2+2*FW);
                        h = i;
                        
                        psi_vxx[j][h1] = b_x[h1] * psi_vxx[j][h1] + a_x[h1] * vxx;
                        vxx = vxx / K_x[h1] + psi_vxx[j][h1]; 

                        /*psi_vyx[j][h] = b_x_half[h] * psi_vyx[j][h] + a_x_half[h] * vyx;
                        vyx = vyx / K_x_half[h] + psi_vyx[j][h];*/
                        
                        psi_vyx[j][h1] = b_x_half[h1] * psi_vyx[j][h1] + a_x_half[h1] * vyx;
			vyx = vyx / K_x_half[h1] + psi_vyx[j][h1];
                                           
         }

	  /* top boundary */                                         
        if((POS[2]==0) && (!(FREE_SURF)) && (j<=FW)){
                                                
                        psi_vyy[j][i] = b_y[j] * psi_vyy[j][i] + a_y[j] * vyy;                                            
                        psi_vxy[j][i] = b_y_half[j] * psi_vxy[j][i] + a_y_half[j] * vxy;
                     
                        vyy = vyy / K_y[j] + psi_vyy[j][i];
                        vxy = vxy / K_y_half[j] + psi_vxy[j][i];

        }
	
	  /* bottom boundary */                                         
        if((POS[2]==NPROCY-1) && (j>=ny2-FW+1)){

                        h1 = (j-ny2+2*FW);                                        
                        h = j;
                                                
                        psi_vyy[h1][i] = b_y[h1] * psi_vyy[h1][i] + a_y[h1] * vyy;                                            
                        vyy = vyy / K_y[h1] + psi_vyy[h1][i];
                        
                        /*psi_vxy[j][i] = b_y_half[j] * psi_vxy[j][i] + a_y_half[j] * vxy;
                        vxy = vxy / K_y_half[j] + psi_vxy[j][i];*/
                        
                        psi_vxy[h1][i] = b_y_half[h1] * psi_vxy[h1][i] + a_y_half[h1] * vxy;
			vxy = vxy / K_y_half[h1] + psi_vxy[h1][i];
        
        }

	/* computing sums of the old memory variables */
			sumr=sump=sumq=0.0;
			for (l=1;l<=L;l++){
				sumr+=r[j][i][l];
				sump+=p[j][i][l];
				sumq+=q[j][i][l];
			}


                        /* updating components of the stress tensor, partially */
			sxy[j][i] += (fipjp[j][i]*(vxy+vyx))+(dthalbe*sumr);
			sxx[j][i] += (g[j][i]*(vxx+vyy))-(2.0*f[j][i]*vyy)+(dthalbe*sump);
			syy[j][i] += (g[j][i]*(vxx+vyy))-(2.0*f[j][i]*vxx)+(dthalbe*sumq);
			
			uxy[j][i] = ((fipjp[j][i]/DT)*(vxy+vyx))+(0.5*sumr);
			ux[j][i] = ((g[j][i]/DT)*(vxx+vyy))-((2.0*f[j][i]/DT)*vyy)+(0.5*sump);
			uy[j][i] = ((g[j][i]/DT)*(vxx+vyy))-((2.0*f[j][i]/DT)*vxx)+(0.5*sumq);
			
			
			/* now updating the memory-variables and sum them up*/
			sumr=sump=sumq=0.0;
			for (l=1;l<=L;l++){
				r[j][i][l] = bip[l]*(r[j][i][l]*cip[l]-(dip[j][i][l]*(vxy+vyx)));
				p[j][i][l] = bjm[l]*(p[j][i][l]*cjm[l]-(e[j][i][l]*(vxx+vyy))+(2.0*d[j][i][l]*vyy));
				q[j][i][l] = bjm[l]*(q[j][i][l]*cjm[l]-(e[j][i][l]*(vxx+vyy))+(2.0*d[j][i][l]*vxx));
				sumr += r[j][i][l];
				sump += p[j][i][l];
				sumq += q[j][i][l];
			}

			/* and now the components of the stress tensor are
			   completely updated */
			sxy[j][i]+=(dthalbe*sumr);
			sxx[j][i]+=(dthalbe*sump);
			syy[j][i]+=(dthalbe*sumq);
			
			uxy[j][i]+=(0.5*sumr);
			ux[j][i]+=(0.5*sump);
			uy[j][i]+=(0.5*sumq);

   }}
		break;

/*	case 10:
		for (j=ny1;j<=ny2;j++){
			for (i=nx1;i<=nx2;i++){
				vxx = (  hc[1]*(vx[j][i]  -vx[j][i-1])
				       + hc[2]*(vx[j][i+1]-vx[j][i-2])
				       + hc[3]*(vx[j][i+2]-vx[j][i-3])
				       + hc[4]*(vx[j][i+3]-vx[j][i-4])
				       + hc[5]*(vx[j][i+4]-vx[j][i-5])
				      )*dhi;
				vyy = (  hc[1]*(vy[j][i]  -vy[j-1][i])
				       + hc[2]*(vy[j+1][i]-vy[j-2][i])
				       + hc[3]*(vy[j+2][i]-vy[j-3][i])
				       + hc[4]*(vy[j+3][i]-vy[j-4][i])
				       + hc[5]*(vy[j+4][i]-vy[j-5][i])
				      )*dhi;
				vyx = (  hc[1]*(vy[j][i+1]-vy[j][i])
				       + hc[2]*(vy[j][i+2]-vy[j][i-1])
				       + hc[3]*(vy[j][i+3]-vy[j][i-2])
				       + hc[4]*(vy[j][i+4]-vy[j][i-3])
				       + hc[5]*(vy[j][i+5]-vy[j][i-4])
				      )*dhi;
				vxy = (  hc[1]*(vx[j+1][i]-vx[j][i])
				       + hc[2]*(vx[j+2][i]-vx[j-1][i])
				       + hc[3]*(vx[j+3][i]-vx[j-2][i])
				       + hc[4]*(vx[j+4][i]-vx[j-3][i])
				       + hc[5]*(vx[j+5][i]-vx[j-4][i])
				      )*dhi;
				      
				ux[j][i] += DT*vxx;
				uy[j][i] += DT*vyy;      
		       		
				fipjp=uipjp[j][i]*DT;
				f=u[j][i]*DT;
				g=pi[j][i]*DT;
				
				sxy[j][i]+=(fipjp*(vxy+vyx));
				sxx[j][i]+=(g*(vxx+vyy))-(2.0*f*vyy);
				syy[j][i]+=(g*(vxx+vyy))-(2.0*f*vxx);*/
				
				/*
				piv = pi[j][i]*(vxx+vyy);
				sxy[j][i]+=(uipjp[j][i]*(vxy+vyx));
				sxx[j][i]+=piv-(2.0*u[j][i]*vyy);
				syy[j][i]+=piv-(2.0*u[j][i]*vxx);
				*/
/*			}
		}
		break;
		
	case 12:
		for (j=ny1;j<=ny2;j++){
			for (i=nx1;i<=nx2;i++){
				vxx = (  hc[1]*(vx[j][i]  -vx[j][i-1])
				       + hc[2]*(vx[j][i+1]-vx[j][i-2])
				       + hc[3]*(vx[j][i+2]-vx[j][i-3])
				       + hc[4]*(vx[j][i+3]-vx[j][i-4])
				       + hc[5]*(vx[j][i+4]-vx[j][i-5])
				       + hc[6]*(vx[j][i+5]-vx[j][i-6])
				      )*dhi;
				vyy = (  hc[1]*(vy[j][i]  -vy[j-1][i])
				       + hc[2]*(vy[j+1][i]-vy[j-2][i])
				       + hc[3]*(vy[j+2][i]-vy[j-3][i])
				       + hc[4]*(vy[j+3][i]-vy[j-4][i])
				       + hc[5]*(vy[j+4][i]-vy[j-5][i])
				       + hc[6]*(vy[j+5][i]-vy[j-6][i])
				      )*dhi;
				vyx = (  hc[1]*(vy[j][i+1]-vy[j][i])
				       + hc[2]*(vy[j][i+2]-vy[j][i-1])
				       + hc[3]*(vy[j][i+3]-vy[j][i-2])
				       + hc[4]*(vy[j][i+4]-vy[j][i-3])
				       + hc[5]*(vy[j][i+5]-vy[j][i-4])
				       + hc[6]*(vy[j][i+6]-vy[j][i-5])
				      )*dhi;
				vxy = (  hc[1]*(vx[j+1][i]-vx[j][i])
				       + hc[2]*(vx[j+2][i]-vx[j-1][i])
				       + hc[3]*(vx[j+3][i]-vx[j-2][i])
				       + hc[4]*(vx[j+4][i]-vx[j-3][i])
				       + hc[5]*(vx[j+5][i]-vx[j-4][i])
				       + hc[6]*(vx[j+6][i]-vx[j-5][i])
				      )*dhi;
				      
				ux[j][i] += DT*vxx;
				uy[j][i] += DT*vyy;
		       
				fipjp=uipjp[j][i]*DT;
				f=u[j][i]*DT;
				g=pi[j][i]*DT;

				sxy[j][i]+=(fipjp*(vxy+vyx));
				sxx[j][i]+=(g*(vxx+vyy))-(2.0*f*vyy);
				syy[j][i]+=(g*(vxx+vyy))-(2.0*f*vxx);
			}
		}
		break;
*/
		
	default:
		for (j=ny1;j<=ny2;j++){
			for (i=nx1;i<=nx2;i++){
				vxx = 0.0;
				vyy = 0.0;
				vyx = 0.0;
				vxy = 0.0;
				for (m=1; m<=fdoh; m++) {
					vxx += hc[m]*(vx[j][i+m-1] -vx[j][i-m]  );
					vyy += hc[m]*(vy[j+m-1][i] -vy[j-m][i]  );
					vyx += hc[m]*(vy[j][i+m]   -vy[j][i-m+1]);
					vxy += hc[m]*(vx[j+m][i]   -vx[j-m+1][i]);
				}
				vxx *= dhi;
				vyy *= dhi;
				vyx *= dhi;
				vxy *= dhi;

				sumr=sump=sumq=0.0;
				for (l=1;l<=L;l++){
					sumr+=r[j][i][l];
					sump+=p[j][i][l];
					sumq+=q[j][i][l];
				}

				sxy[j][i] += (fipjp[j][i]*(vxy+vyx))+(dthalbe*sumr);
				sxx[j][i] += (g[j][i]*(vxx+vyy))-(2.0*f[j][i]*vyy)+(dthalbe*sump);
				syy[j][i] += (g[j][i]*(vxx+vyy))-(2.0*f[j][i]*vxx)+(dthalbe*sumq);

				sumr=sump=sumq=0.0;
				for (l=1;l<=L;l++){
					r[j][i][l] = bip[l]*(r[j][i][l]*cip[l]-(dip[j][i][l]*(vxy+vyx)));
					p[j][i][l] = bjm[l]*(p[j][i][l]*cjm[l]-(e[j][i][l]*(vxx+vyy))+(2.0*d[j][i][l]*vyy));
					q[j][i][l] = bjm[l]*(q[j][i][l]*cjm[l]-(e[j][i][l]*(vxx+vyy))+(2.0*d[j][i][l]*vxx));
					sumr += r[j][i][l];
					sump += p[j][i][l];
					sumq += q[j][i][l];
				}

				sxy[j][i]+=(dthalbe*sumr);
				sxx[j][i]+=(dthalbe*sump);
				syy[j][i]+=(dthalbe*sumq);
			}
		}
		break;
		
	} /* end of switch(FDORDER) */


	if (infoout && (MYID==0)){
		time2=MPI_Wtime();
		fprintf(FP," finished (real time: %4.2f s).\n",time2-time1);
	}
}