update_v_ssg.c 86.7 KB
Newer Older
Tilman Steinweg's avatar
Tilman Steinweg 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 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 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 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208 1209 1210 1211 1212 1213 1214 1215 1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 1234 1235 1236 1237 1238 1239 1240 1241 1242 1243 1244 1245 1246 1247 1248 1249 1250 1251 1252 1253 1254 1255 1256 1257 1258 1259 1260 1261 1262 1263 1264 1265 1266 1267 1268 1269 1270 1271 1272 1273 1274 1275 1276 1277 1278 1279 1280 1281 1282 1283 1284 1285 1286 1287 1288 1289 1290 1291 1292 1293 1294 1295 1296 1297 1298 1299 1300 1301 1302 1303 1304 1305 1306 1307 1308 1309 1310 1311 1312 1313 1314 1315 1316 1317 1318 1319 1320 1321 1322 1323 1324 1325 1326 1327 1328 1329 1330 1331 1332 1333 1334 1335 1336 1337 1338 1339 1340 1341 1342 1343 1344 1345 1346 1347 1348 1349 1350 1351 1352 1353 1354 1355 1356 1357 1358 1359 1360 1361 1362 1363 1364 1365 1366 1367 1368 1369 1370 1371 1372 1373 1374 1375 1376 1377 1378 1379 1380 1381 1382 1383 1384 1385 1386 1387 1388 1389 1390 1391 1392 1393 1394 1395 1396 1397 1398 1399 1400 1401 1402 1403 1404 1405 1406 1407 1408 1409 1410 1411 1412 1413 1414 1415 1416 1417 1418 1419 1420 1421 1422 1423 1424 1425 1426 1427 1428 1429 1430 1431 1432 1433 1434 1435 1436 1437 1438 1439 1440 1441 1442 1443 1444 1445 1446 1447 1448 1449 1450 1451 1452 1453 1454 1455 1456 1457 1458 1459 1460 1461 1462 1463 1464 1465 1466 1467 1468 1469 1470 1471 1472 1473 1474 1475 1476 1477 1478 1479 1480 1481 1482 1483 1484 1485 1486 1487 1488 1489 1490 1491 1492 1493 1494 1495 1496 1497 1498 1499 1500 1501 1502 1503 1504 1505
/*------------------------------------------------------------------------
 * Copyright (C) 2011 For the list of authors, see file AUTHORS.
 *
 * This file is part of SOFI3D.
 * 
 * SOFI3D 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.
 * 
 * SOFI3D 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 SOFI3D. See file COPYING and/or 
  * <http://www.gnu.org/licenses/gpl-2.0.html>.
--------------------------------------------------------------------------*/
/*------------------------------------------------------------------------
 *   updating particle velocities at gridpoints [nx1...nx2][ny1...ny2][nz1...nz2]
 *   by a staggered grid finite difference scheme of 4th order accuracy in space
 *   and second order accuracy in time
 *  ----------------------------------------------------------------------*/

#include "fd.h"

double update_v(int nx1, int nx2, int ny1, int ny2, int nz1, int nz2,
		int nt, float *** vx, float *** vy, float *** vz,
		float *** sxx, float *** syy, float *** szz, float *** sxy,
		float *** syz, float *** sxz, float  ***  rho,  float  *** rjp, float  *** rkp, float  *** rip,
		float **  srcpos_loc, float ** signals, int nsrc, float *** absorb_coeff, int * stype, float *** svx, float *** svy, float *** svz,
        float *** svx_2, float *** svy_2, float *** svz_2, float *** svx_3, float *** svy_3, float *** svz_3,
        float *** svx_4, float *** svy_4, float *** svz_4){


	extern float DT, DX, DY, DZ, SOURCE_ALPHA, SOURCE_BETA;
	double time=0.0, time1=0.0, time2=0.0;
	extern int MYID, FDORDER, FDORDER_TIME, LOG, ABS_TYPE, FDCOEFF; /* variable "FW" removed, not in use */
	extern FILE *FP;
	extern int OUTNTIMESTEPINFO;

	int i, j, k, l;
	float  amp, alpha_rad, beta_rad;
	float b1, b2, b3, b4, b5, b6, dx, dy, dz;
	float sxx_x, sxy_y, sxz_z, syy_y, sxy_x, syz_z;
	float szz_z, sxz_x, syz_y;
    float c1, c2, c3, c4; /* Coefficients for Adam Bashforth */

    float **svx_j,**svy_j,**svz_j;
    float **svx_j_2,**svy_j_2,**svz_j_2;
    float **svx_j_3,**svy_j_3,**svz_j_3;
    float **svx_j_4,**svy_j_4,**svz_j_4;
    
    float *svx_j_i,*svy_j_i,*svz_j_i;
    float *svx_j_i_2,*svy_j_i_2,*svz_j_i_2;
    float *svx_j_i_3,*svy_j_i_3,*svz_j_i_3;
    float *svx_j_i_4,*svy_j_i_4,*svz_j_i_4;
    
    
	if (LOG)
		if ((MYID==0) && ((nt+(OUTNTIMESTEPINFO-1))%OUTNTIMESTEPINFO)==0) time1=MPI_Wtime();

    
    
    switch (FDORDER_TIME){
            
        case 2:
            
            switch (FDORDER){
                    
                case 2 :
                    
                    dx=DT/DX;
                    dy=DT/DY;
                    dz=DT/DZ;
                    
                    b1=1.0; /* Taylor coefficients*/
                    if(FDCOEFF==2){
                        b1=1.00100; } /* Holberg coefficients E=0.1 %*/
                    
                    for (j=ny1;j<=ny2;j++){
                        for (i=nx1;i<=nx2;i++){
                            for (k=nz1;k<=nz2;k++){
                                
                                sxx_x = dx*b1*(sxx[j][i+1][k]-sxx[j][i][k]);
                                sxy_y = dy*b1*(sxy[j][i][k]-sxy[j-1][i][k]);
                                sxz_z = dz*b1*(sxz[j][i][k]-sxz[j][i][k-1]); /* backward operator */
                                
                                /* updating components of particle velocities */
                                vx[j][i][k]+= ((sxx_x + sxy_y +sxz_z)/rip[j][i][k]);
                                
                                syy_y = dy*b1*(syy[j+1][i][k]-syy[j][i][k]);
                                sxy_x = dx*b1*(sxy[j][i][k]-sxy[j][i-1][k]);
                                syz_z = dz*b1*(syz[j][i][k]-syz[j][i][k-1]);
                                
                                
                                vy[j][i][k]+= ((syy_y + sxy_x + syz_z)/rjp[j][i][k]);
                                
                                szz_z = dz*b1*(szz[j][i][k+1]-szz[j][i][k]);
                                sxz_x = dx*b1*(sxz[j][i][k]-sxz[j][i-1][k]);
                                syz_y = dy*b1*(syz[j][i][k]-syz[j-1][i][k]);
                                
                                
                                vz[j][i][k]+= ((szz_z + sxz_x + syz_y)/rkp[j][i][k]);
                                
                            }
                        }
                    }
                    
                    break;
                    
                case 4 :
                    
                    dx=DT/DX;
                    dy=DT/DY;
                    dz=DT/DZ;
                    
                    b1=9.0/8.0; b2=-1.0/24.0; /* Taylor coefficients*/
                    if(FDCOEFF==2){
                        b1=1.1382; b2=-0.046414;} /* Holberg coefficients E=0.1 %*/
                    
                    for (j=ny1;j<=ny2;j++){
                        for (i=nx1;i<=nx2;i++){
                            for (k=nz1;k<=nz2;k++){
                                
                                sxx_x = dx*(b1*(sxx[j][i+1][k]-sxx[j][i][k])+b2*(sxx[j][i+2][k]-sxx[j][i-1][k]));
                                sxy_y = dy*(b1*(sxy[j][i][k]-sxy[j-1][i][k])+b2*(sxy[j+1][i][k]-sxy[j-2][i][k]));
                                sxz_z = dz*(b1*(sxz[j][i][k]-sxz[j][i][k-1])+b2*(sxz[j][i][k+1]-sxz[j][i][k-2]));
                                
                                /* updating components of particle velocities */
                                vx[j][i][k]+= (sxx_x + sxy_y +sxz_z)/rip[j][i][k];
                                
                                syy_y = dy*(b1*(syy[j+1][i][k]-syy[j][i][k])+b2*(syy[j+2][i][k]-syy[j-1][i][k]));
                                sxy_x = dx*(b1*(sxy[j][i][k]-sxy[j][i-1][k])+b2*(sxy[j][i+1][k]-sxy[j][i-2][k]));
                                syz_z = dz*(b1*(syz[j][i][k]-syz[j][i][k-1])+b2*(syz[j][i][k+1]-syz[j][i][k-2]));
                                
                                
                                vy[j][i][k]+= (syy_y + sxy_x + syz_z)/rjp[j][i][k];
                                
                                szz_z = dz*(b1*(szz[j][i][k+1]-szz[j][i][k])+b2*(szz[j][i][k+2]-szz[j][i][k-1]));
                                sxz_x = dx*(b1*(sxz[j][i][k]-sxz[j][i-1][k])+b2*(sxz[j][i+1][k]-sxz[j][i-2][k]));
                                syz_y = dy*(b1*(syz[j][i][k]-syz[j-1][i][k])+b2*(syz[j+1][i][k]-syz[j-2][i][k]));
                                
                                
                                vz[j][i][k]+= (szz_z + sxz_x + syz_y)/rkp[j][i][k];
                                
                            }
                        }
                    }
                    
                    break;
                    
                case 6 :
                    
                    dx=DT/DX;
                    dy=DT/DY;
                    dz=DT/DZ;
                    
                    b1=75.0/64.0; b2=-25.0/384.0; b3=3.0/640.0; /* Taylor coefficients*/
                    if(FDCOEFF==2){
                        b1=1.1965; b2=-0.078804; b3=0.0081781;}   /* Holberg coefficients E=0.1 %*/
                    
                    for (j=ny1;j<=ny2;j++){
                        for (i=nx1;i<=nx2;i++){
                            for (k=nz1;k<=nz2;k++){
                                
                                sxx_x = dx*(b1*(sxx[j][i+1][k]-sxx[j][i][k])+
                                            b2*(sxx[j][i+2][k]-sxx[j][i-1][k])+
                                            b3*(sxx[j][i+3][k]-sxx[j][i-2][k]));
                                
                                sxy_y = dy*(b1*(sxy[j][i][k]-sxy[j-1][i][k])+
                                            b2*(sxy[j+1][i][k]-sxy[j-2][i][k])+
                                            b3*(sxy[j+2][i][k]-sxy[j-3][i][k]));
                                
                                sxz_z = dz*(b1*(sxz[j][i][k]-sxz[j][i][k-1])+
                                            b2*(sxz[j][i][k+1]-sxz[j][i][k-2])+
                                            b3*(sxz[j][i][k+2]-sxz[j][i][k-3]));
                                
                                
                                /* updating components of particle velocities */
                                vx[j][i][k]+= (sxx_x + sxy_y +sxz_z)/rip[j][i][k];
                                
                                syy_y = dy*(b1*(syy[j+1][i][k]-syy[j][i][k])+
                                            b2*(syy[j+2][i][k]-syy[j-1][i][k])+
                                            b3*(syy[j+3][i][k]-syy[j-2][i][k]));
                                
                                sxy_x = dx*(b1*(sxy[j][i][k]-sxy[j][i-1][k])+
                                            b2*(sxy[j][i+1][k]-sxy[j][i-2][k])+
                                            b3*(sxy[j][i+2][k]-sxy[j][i-3][k]));
                                
                                syz_z = dz*(b1*(syz[j][i][k]-syz[j][i][k-1])+
                                            b2*(syz[j][i][k+1]-syz[j][i][k-2])+
                                            b3*(syz[j][i][k+2]-syz[j][i][k-3]));
                                
                                
                                vy[j][i][k]+= (syy_y + sxy_x + syz_z)/rjp[j][i][k];
                                
                                szz_z = dz*(b1*(szz[j][i][k+1]-szz[j][i][k])+
                                            b2*(szz[j][i][k+2]-szz[j][i][k-1])+
                                            b3*(szz[j][i][k+3]-szz[j][i][k-2]));
                                
                                sxz_x = dx*(b1*(sxz[j][i][k]-sxz[j][i-1][k])+
                                            b2*(sxz[j][i+1][k]-sxz[j][i-2][k])+
                                            b3*(sxz[j][i+2][k]-sxz[j][i-3][k]));
                                
                                
                                syz_y = dy*(b1*(syz[j][i][k]-syz[j-1][i][k])+
                                            b2*(syz[j+1][i][k]-syz[j-2][i][k])+
                                            b3*(syz[j+2][i][k]-syz[j-3][i][k]));
                                
                                
                                vz[j][i][k]+= (szz_z + sxz_x + syz_y)/rkp[j][i][k];
                                
                            }
                        }
                    }
                    
                    break;
                    
                case 8 :
                    
                    dx=DT/DX;
                    dy=DT/DY;
                    dz=DT/DZ;
                    
                    b1=1225.0/1024.0; b2=-245.0/3072.0; b3=49.0/5120.0; b4=-5.0/7168.0; /* Taylor coefficients*/
                    if(FDCOEFF==2){
                        b1=1.2257; b2=-0.099537; b3=0.018063; b4=-0.0026274;} /* Holberg coefficients E=0.1 %*/
                    
                    for (j=ny1;j<=ny2;j++){
                        for (i=nx1;i<=nx2;i++){
                            for (k=nz1;k<=nz2;k++){
                                
                                sxx_x = dx*(b1*(sxx[j][i+1][k]-sxx[j][i][k])+
                                            b2*(sxx[j][i+2][k]-sxx[j][i-1][k])+
                                            b3*(sxx[j][i+3][k]-sxx[j][i-2][k])+
                                            b4*(sxx[j][i+4][k]-sxx[j][i-3][k]));
                                
                                sxy_y = dy*(b1*(sxy[j][i][k]-sxy[j-1][i][k])+
                                            b2*(sxy[j+1][i][k]-sxy[j-2][i][k])+
                                            b3*(sxy[j+2][i][k]-sxy[j-3][i][k])+
                                            b4*(sxy[j+3][i][k]-sxy[j-4][i][k]));
                                
                                sxz_z = dz*(b1*(sxz[j][i][k]-sxz[j][i][k-1])+
                                            b2*(sxz[j][i][k+1]-sxz[j][i][k-2])+
                                            b3*(sxz[j][i][k+2]-sxz[j][i][k-3])+
                                            b4*(sxz[j][i][k+3]-sxz[j][i][k-4]));
                                
                                /* updating components of particle velocities */
                                vx[j][i][k]+= (sxx_x + sxy_y +sxz_z)/rip[j][i][k];
                                
                                syy_y = dy*(b1*(syy[j+1][i][k]-syy[j][i][k])+
                                            b2*(syy[j+2][i][k]-syy[j-1][i][k])+
                                            b3*(syy[j+3][i][k]-syy[j-2][i][k])+
                                            b4*(syy[j+4][i][k]-syy[j-3][i][k]));
                                
                                sxy_x = dx*(b1*(sxy[j][i][k]-sxy[j][i-1][k])+
                                            b2*(sxy[j][i+1][k]-sxy[j][i-2][k])+
                                            b3*(sxy[j][i+2][k]-sxy[j][i-3][k])+
                                            b4*(sxy[j][i+3][k]-sxy[j][i-4][k]));
                                
                                syz_z = dz*(b1*(syz[j][i][k]-syz[j][i][k-1])+
                                            b2*(syz[j][i][k+1]-syz[j][i][k-2])+
                                            b3*(syz[j][i][k+2]-syz[j][i][k-3])+
                                            b4*(syz[j][i][k+3]-syz[j][i][k-4]));
                                
                                
                                vy[j][i][k]+= (syy_y + sxy_x + syz_z)/rjp[j][i][k];
                                
                                szz_z = dz*(b1*(szz[j][i][k+1]-szz[j][i][k])+
                                            b2*(szz[j][i][k+2]-szz[j][i][k-1])+
                                            b3*(szz[j][i][k+3]-szz[j][i][k-2])+
                                            b4*(szz[j][i][k+4]-szz[j][i][k-3]));
                                
                                sxz_x = dx*(b1*(sxz[j][i][k]-sxz[j][i-1][k])+
                                            b2*(sxz[j][i+1][k]-sxz[j][i-2][k])+
                                            b3*(sxz[j][i+2][k]-sxz[j][i-3][k])+
                                            b4*(sxz[j][i+3][k]-sxz[j][i-4][k]));
                                
                                
                                syz_y = dy*(b1*(syz[j][i][k]-syz[j-1][i][k])+
                                            b2*(syz[j+1][i][k]-syz[j-2][i][k])+
                                            b3*(syz[j+2][i][k]-syz[j-3][i][k])+
                                            b4*(syz[j+3][i][k]-syz[j-4][i][k]));
                                
                                
                                vz[j][i][k]+= (szz_z + sxz_x + syz_y)/rkp[j][i][k];
                                
                            }
                        }
                    }
                    
                    break;
                    
                case 10 :
                    
                    dx=DT/DX;
                    dy=DT/DY;
                    dz=DT/DZ;
                    
                    
                    b1=19845.0/16384.0; b2=-735.0/8192.0; b3=567.0/40960.0; b4=-405.0/229376.0; b5=35.0/294912.0; /* Taylor Coefficients*/
                    if(FDCOEFF==2){
                        b1=1.2415; b2=-0.11231; b3=0.026191; b4=-0.0064682; b5=0.001191;} /* Holberg coefficients E=0.1 %*/
                    
                    for (j=ny1;j<=ny2;j++){
                        for (i=nx1;i<=nx2;i++){
                            for (k=nz1;k<=nz2;k++){
                                
                                sxx_x = dx*(b1*(sxx[j][i+1][k]-sxx[j][i][k])+
                                            b2*(sxx[j][i+2][k]-sxx[j][i-1][k])+
                                            b3*(sxx[j][i+3][k]-sxx[j][i-2][k])+
                                            b4*(sxx[j][i+4][k]-sxx[j][i-3][k])+
                                            b5*(sxx[j][i+5][k]-sxx[j][i-4][k]));
                                
                                sxy_y = dy*(b1*(sxy[j][i][k]-sxy[j-1][i][k])+
                                            b2*(sxy[j+1][i][k]-sxy[j-2][i][k])+
                                            b3*(sxy[j+2][i][k]-sxy[j-3][i][k])+
                                            b4*(sxy[j+3][i][k]-sxy[j-4][i][k])+
                                            b5*(sxy[j+4][i][k]-sxy[j-5][i][k]));
                                
                                sxz_z = dz*(b1*(sxz[j][i][k]-sxz[j][i][k-1])+
                                            b2*(sxz[j][i][k+1]-sxz[j][i][k-2])+
                                            b3*(sxz[j][i][k+2]-sxz[j][i][k-3])+
                                            b4*(sxz[j][i][k+3]-sxz[j][i][k-4])+
                                            b5*(sxz[j][i][k+4]-sxz[j][i][k-5]));
                                
                                
                                /* updating components of particle velocities */
                                vx[j][i][k]+= (sxx_x + sxy_y +sxz_z)/rip[j][i][k];
                                
                                syy_y = dy*(b1*(syy[j+1][i][k]-syy[j][i][k])+
                                            b2*(syy[j+2][i][k]-syy[j-1][i][k])+
                                            b3*(syy[j+3][i][k]-syy[j-2][i][k])+
                                            b4*(syy[j+4][i][k]-syy[j-3][i][k])+
                                            b5*(syy[j+5][i][k]-syy[j-4][i][k]));
                                
                                sxy_x = dx*(b1*(sxy[j][i][k]-sxy[j][i-1][k])+
                                            b2*(sxy[j][i+1][k]-sxy[j][i-2][k])+
                                            b3*(sxy[j][i+2][k]-sxy[j][i-3][k])+
                                            b4*(sxy[j][i+3][k]-sxy[j][i-4][k])+
                                            b5*(sxy[j][i+4][k]-sxy[j][i-5][k]));
                                
                                syz_z = dz*(b1*(syz[j][i][k]-syz[j][i][k-1])+
                                            b2*(syz[j][i][k+1]-syz[j][i][k-2])+
                                            b3*(syz[j][i][k+2]-syz[j][i][k-3])+
                                            b4*(syz[j][i][k+3]-syz[j][i][k-4])+
                                            b5*(syz[j][i][k+4]-syz[j][i][k-5]));
                                
                                
                                vy[j][i][k]+= (syy_y + sxy_x + syz_z)/rjp[j][i][k];
                                
                                szz_z = dz*(b1*(szz[j][i][k+1]-szz[j][i][k])+
                                            b2*(szz[j][i][k+2]-szz[j][i][k-1])+
                                            b3*(szz[j][i][k+3]-szz[j][i][k-2])+
                                            b4*(szz[j][i][k+4]-szz[j][i][k-3])+
                                            b5*(szz[j][i][k+5]-szz[j][i][k-4]));
                                
                                sxz_x = dx*(b1*(sxz[j][i][k]-sxz[j][i-1][k])+
                                            b2*(sxz[j][i+1][k]-sxz[j][i-2][k])+
                                            b3*(sxz[j][i+2][k]-sxz[j][i-3][k])+
                                            b4*(sxz[j][i+3][k]-sxz[j][i-4][k])+
                                            b5*(sxz[j][i+4][k]-sxz[j][i-5][k]));
                                
                                
                                syz_y = dy*(b1*(syz[j][i][k]-syz[j-1][i][k])+
                                            b2*(syz[j+1][i][k]-syz[j-2][i][k])+
                                            b3*(syz[j+2][i][k]-syz[j-3][i][k])+
                                            b4*(syz[j+3][i][k]-syz[j-4][i][k])+
                                            b5*(syz[j+4][i][k]-syz[j-5][i][k]));
                                
                                
                                vz[j][i][k]+= (szz_z + sxz_x + syz_y)/rkp[j][i][k];
                                
                            }
                        }
                    }
                    
                    break;
                    
                case 12 :
                    dx=DT/DX;
                    dy=DT/DY;
                    dz=DT/DZ;
                    
                    /* Taylor coefficients */
                    b1=160083.0/131072.0; b2=-12705.0/131072.0; b3=22869.0/1310720.0;
                    b4=-5445.0/1835008.0; b5=847.0/2359296.0; b6=-63.0/2883584;
                    
                    /* Holberg coefficients E=0.1 %*/
                    if(FDCOEFF==2){
                        b1=1.2508; b2=-0.12034; b3=0.032131; b4=-0.010142; b5=0.0029857; b6=-0.00066667;}
                    
                    
                    for (j=ny1;j<=ny2;j++){
                        for (i=nx1;i<=nx2;i++){
                            for (k=nz1;k<=nz2;k++){
                                
                                sxx_x = dx*(b1*(sxx[j][i+1][k]-sxx[j][i][k])+
                                            b2*(sxx[j][i+2][k]-sxx[j][i-1][k])+
                                            b3*(sxx[j][i+3][k]-sxx[j][i-2][k])+
                                            b4*(sxx[j][i+4][k]-sxx[j][i-3][k])+
                                            b5*(sxx[j][i+5][k]-sxx[j][i-4][k])+
                                            b6*(sxx[j][i+6][k]-sxx[j][i-5][k]));
                                
                                sxy_y = dy*(b1*(sxy[j][i][k]-sxy[j-1][i][k])+
                                            b2*(sxy[j+1][i][k]-sxy[j-2][i][k])+
                                            b3*(sxy[j+2][i][k]-sxy[j-3][i][k])+
                                            b4*(sxy[j+3][i][k]-sxy[j-4][i][k])+
                                            b5*(sxy[j+4][i][k]-sxy[j-5][i][k])+
                                            b6*(sxy[j+5][i][k]-sxy[j-6][i][k]));
                                
                                sxz_z = dz*(b1*(sxy[j][i][k]-sxy[j][i][k-1])+
                                            b2*(sxy[j][i][k+1]-sxy[j][i][k-2])+
                                            b3*(sxy[j][i][k+2]-sxy[j][i][k-3])+
                                            b4*(sxy[j][i][k+3]-sxy[j][i][k-4])+
                                            b5*(sxy[j][i][k+4]-sxy[j][i][k-5])+
                                            b6*(sxy[j][i][k+5]-sxy[j][i][k-6]));
                                
                                
                                /* updating components of particle velocities */
                                vx[j][i][k]+= (sxx_x + sxy_y +sxz_z)/rip[j][i][k];
                                
                                syy_y = dy*(b1*(syy[j+1][i][k]-syy[j][i][k])+
                                            b2*(syy[j+2][i][k]-syy[j-1][i][k])+
                                            b3*(syy[j+3][i][k]-syy[j-2][i][k])+
                                            b4*(syy[j+4][i][k]-syy[j-3][i][k])+
                                            b5*(syy[j+5][i][k]-syy[j-4][i][k])+
                                            b6*(syy[j+6][i][k]-syy[j-5][i][k]));
                                
                                sxy_x = dx*(b1*(sxy[j][i][k]-sxy[j][i-1][k])+
                                            b2*(sxy[j][i+1][k]-sxy[j][i-2][k])+
                                            b3*(sxy[j][i+2][k]-sxy[j][i-3][k])+
                                            b4*(sxy[j][i+3][k]-sxy[j][i-4][k])+
                                            b5*(sxy[j][i+4][k]-sxy[j][i-5][k])+
                                            b6*(sxy[j][i+5][k]-sxy[j][i-6][k]));
                                
                                syz_z = dz*(b1*(syz[j][i][k]-syz[j][i][k-1])+
                                            b2*(syz[j][i][k+1]-syz[j][i][k-2])+
                                            b3*(syz[j][i][k+2]-syz[j][i][k-3])+
                                            b4*(syz[j][i][k+3]-syz[j][i][k-4])+
                                            b5*(syz[j][i][k+4]-syz[j][i][k-5])+
                                            b6*(syz[j][i][k+5]-syz[j][i][k-6]));
                                
                                
                                
                                vy[j][i][k]+= (syy_y + sxy_x + syz_z)/rjp[j][i][k];
                                
                                szz_z = dz*(b1*(szz[j][i][k+1]-szz[j][i][k])+
                                            b2*(szz[j][i][k+2]-szz[j][i][k-1])+
                                            b3*(szz[j][i][k+3]-szz[j][i][k-2])+
                                            b4*(szz[j][i][k+4]-szz[j][i][k-3])+
                                            b5*(szz[j][i][k+5]-szz[j][i][k-4])+
                                            b6*(szz[j][i][k+6]-szz[j][i][k-5]));
                                
                                sxz_x = dx*(b1*(sxz[j][i][k]-sxz[j][i-1][k])+
                                            b2*(sxz[j][i+1][k]-sxz[j][i-2][k])+
                                            b3*(sxz[j][i+2][k]-sxz[j][i-3][k])+
                                            b4*(sxz[j][i+3][k]-sxz[j][i-4][k])+
                                            b5*(sxz[j][i+4][k]-sxz[j][i-5][k])+
                                            b6*(sxz[j][i+5][k]-sxz[j][i-6][k]));
                                
                                
                                syz_y = dy*(b1*(syz[j][i][k]-syz[j-1][i][k])+
                                            b2*(syz[j+1][i][k]-syz[j-2][i][k])+
                                            b3*(syz[j+2][i][k]-syz[j-3][i][k])+
                                            b4*(syz[j+3][i][k]-syz[j-4][i][k])+
                                            b5*(syz[j+4][i][k]-syz[j-5][i][k])+
                                            b6*(syz[j+5][i][k]-syz[j-6][i][k]));
                                
                                
                                vz[j][i][k]+= (szz_z + sxz_x + syz_y)/rkp[j][i][k];
                                
                            }
                        }
                    }
                    
                    break;
                    
            }
            break; /* break for FDORDER_TIME=2 */
            
        case 3:
            
            c1=25.0/24.0; c2=-1.0/12.0; c3=1.0/24.0;
            
            switch (FDORDER){
                    
                case 2 :
                    
                    dx=DT/DX;
                    dy=DT/DY;
                    dz=DT/DZ;
                    
                    b1=1.0; /* Taylor coefficients*/
                    if(FDCOEFF==2){
                        b1=1.00100; } /* Holberg coefficients E=0.1 %*/
                    
                    for (j=ny1;j<=ny2;j++){
                        svx_j=*(svx+j);svy_j=*(svy+j);svz_j=*(svz+j);
                        svx_j_2=*(svx_2+j);svy_j_2=*(svy_2+j);svz_j_2=*(svz_2+j);
                        svx_j_3=*(svx_3+j);svy_j_3=*(svy_3+j);svz_j_3=*(svz_3+j);
                        
                        
                        for (i=nx1;i<=nx2;i++){
                            svx_j_i=*(svx_j+i);svy_j_i=*(svy_j+i);svz_j_i=*(svz_j+i);
                            svx_j_i_2=*(svx_j_2+i);svy_j_i_2=*(svy_j_2+i);svz_j_i_2=*(svz_j_2+i);
                            svx_j_i_3=*(svx_j_3+i);svy_j_i_3=*(svy_j_3+i);svz_j_i_3=*(svz_j_3+i);
                            
                            for (k=nz1;k<=nz2;k++){
                                
                                sxx_x = dx*b1*(sxx[j][i+1][k]-sxx[j][i][k]);
                                sxy_y = dy*b1*(sxy[j][i][k]-sxy[j-1][i][k]);
                                sxz_z = dz*b1*(sxz[j][i][k]-sxz[j][i][k-1]); /* backward operator */
                                
                                syy_y = dy*b1*(syy[j+1][i][k]-syy[j][i][k]);
                                sxy_x = dx*b1*(sxy[j][i][k]-sxy[j][i-1][k]);
                                syz_z = dz*b1*(syz[j][i][k]-syz[j][i][k-1]);
                                
                                szz_z = dz*b1*(szz[j][i][k+1]-szz[j][i][k]);
                                sxz_x = dx*b1*(sxz[j][i][k]-sxz[j][i-1][k]);
                                syz_y = dy*b1*(syz[j][i][k]-syz[j-1][i][k]);
                                
                                *(svx_j_i+k)=sxx_x + sxy_y +sxz_z;
                                *(svy_j_i+k)=syy_y + sxy_x + syz_z;
                                *(svz_j_i+k)=szz_z + sxz_x + syz_y;
                                
                                
                                /* updating components of particle velocities */
                                vx[j][i][k]+= ((c1*(*(svx_j_i+k))+c2*(*(svx_j_i_2+k))+c3*(*(svx_j_i_3+k)))/rip[j][i][k]);
                                vy[j][i][k]+= ((c1*(*(svy_j_i+k))+c2*(*(svy_j_i_2+k))+c3*(*(svy_j_i_3+k)))/rjp[j][i][k]);
                                vz[j][i][k]+= ((c1*(*(svz_j_i+k))+c2*(*(svz_j_i_2+k))+c3*(*(svz_j_i_3+k)))/rkp[j][i][k]);
                                
                            }
                        }
                    }
                    
                    break;
                    
                case 4 :
                    
                    dx=DT/DX;
                    dy=DT/DY;
                    dz=DT/DZ;
                    
                    b1=9.0/8.0; b2=-1.0/24.0; /* Taylor coefficients*/
                    if(FDCOEFF==2){
                        b1=1.1382; b2=-0.046414;} /* Holberg coefficients E=0.1 %*/
                    
                    for (j=ny1;j<=ny2;j++){
                        svx_j=*(svx+j);svy_j=*(svy+j);svz_j=*(svz+j);
                        svx_j_2=*(svx_2+j);svy_j_2=*(svy_2+j);svz_j_2=*(svz_2+j);
                        svx_j_3=*(svx_3+j);svy_j_3=*(svy_3+j);svz_j_3=*(svz_3+j);
                        
                        
                        for (i=nx1;i<=nx2;i++){
                            svx_j_i=*(svx_j+i);svy_j_i=*(svy_j+i);svz_j_i=*(svz_j+i);
                            svx_j_i_2=*(svx_j_2+i);svy_j_i_2=*(svy_j_2+i);svz_j_i_2=*(svz_j_2+i);
                            svx_j_i_3=*(svx_j_3+i);svy_j_i_3=*(svy_j_3+i);svz_j_i_3=*(svz_j_3+i);
                            
                            for (k=nz1;k<=nz2;k++){
                                
                                sxx_x = dx*(b1*(sxx[j][i+1][k]-sxx[j][i][k])+b2*(sxx[j][i+2][k]-sxx[j][i-1][k]));
                                sxy_y = dy*(b1*(sxy[j][i][k]-sxy[j-1][i][k])+b2*(sxy[j+1][i][k]-sxy[j-2][i][k]));
                                sxz_z = dz*(b1*(sxz[j][i][k]-sxz[j][i][k-1])+b2*(sxz[j][i][k+1]-sxz[j][i][k-2]));
                                
                                syy_y = dy*(b1*(syy[j+1][i][k]-syy[j][i][k])+b2*(syy[j+2][i][k]-syy[j-1][i][k]));
                                sxy_x = dx*(b1*(sxy[j][i][k]-sxy[j][i-1][k])+b2*(sxy[j][i+1][k]-sxy[j][i-2][k]));
                                syz_z = dz*(b1*(syz[j][i][k]-syz[j][i][k-1])+b2*(syz[j][i][k+1]-syz[j][i][k-2]));
                                
                                szz_z = dz*(b1*(szz[j][i][k+1]-szz[j][i][k])+b2*(szz[j][i][k+2]-szz[j][i][k-1]));
                                sxz_x = dx*(b1*(sxz[j][i][k]-sxz[j][i-1][k])+b2*(sxz[j][i+1][k]-sxz[j][i-2][k]));
                                syz_y = dy*(b1*(syz[j][i][k]-syz[j-1][i][k])+b2*(syz[j+1][i][k]-syz[j-2][i][k]));
                                
                                
                                *(svx_j_i+k)=sxx_x + sxy_y +sxz_z;
                                *(svy_j_i+k)=syy_y + sxy_x + syz_z;
                                *(svz_j_i+k)=szz_z + sxz_x + syz_y;
                                
                                
                                /* updating components of particle velocities */
                                vx[j][i][k]+= ((c1*(*(svx_j_i+k))+c2*(*(svx_j_i_2+k))+c3*(*(svx_j_i_3+k)))/rip[j][i][k]);
                                vy[j][i][k]+= ((c1*(*(svy_j_i+k))+c2*(*(svy_j_i_2+k))+c3*(*(svy_j_i_3+k)))/rjp[j][i][k]);
                                vz[j][i][k]+= ((c1*(*(svz_j_i+k))+c2*(*(svz_j_i_2+k))+c3*(*(svz_j_i_3+k)))/rkp[j][i][k]);
                                
                                
                            }
                        }
                    }
                    
                    break;
                    
                case 6 :
                    
                    dx=DT/DX;
                    dy=DT/DY;
                    dz=DT/DZ;
                    
                    b1=75.0/64.0; b2=-25.0/384.0; b3=3.0/640.0; /* Taylor coefficients*/
                    if(FDCOEFF==2){
                        b1=1.1965; b2=-0.078804; b3=0.0081781;}   /* Holberg coefficients E=0.1 %*/
                    
                    for (j=ny1;j<=ny2;j++){
                        svx_j=*(svx+j);svy_j=*(svy+j);svz_j=*(svz+j);
                        svx_j_2=*(svx_2+j);svy_j_2=*(svy_2+j);svz_j_2=*(svz_2+j);
                        svx_j_3=*(svx_3+j);svy_j_3=*(svy_3+j);svz_j_3=*(svz_3+j);
                        
                        
                        for (i=nx1;i<=nx2;i++){
                            svx_j_i=*(svx_j+i);svy_j_i=*(svy_j+i);svz_j_i=*(svz_j+i);
                            svx_j_i_2=*(svx_j_2+i);svy_j_i_2=*(svy_j_2+i);svz_j_i_2=*(svz_j_2+i);
                            svx_j_i_3=*(svx_j_3+i);svy_j_i_3=*(svy_j_3+i);svz_j_i_3=*(svz_j_3+i);
                            
                            for (k=nz1;k<=nz2;k++){
                                
                                sxx_x = dx*(b1*(sxx[j][i+1][k]-sxx[j][i][k])+
                                            b2*(sxx[j][i+2][k]-sxx[j][i-1][k])+
                                            b3*(sxx[j][i+3][k]-sxx[j][i-2][k]));
                                
                                sxy_y = dy*(b1*(sxy[j][i][k]-sxy[j-1][i][k])+
                                            b2*(sxy[j+1][i][k]-sxy[j-2][i][k])+
                                            b3*(sxy[j+2][i][k]-sxy[j-3][i][k]));
                                
                                sxz_z = dz*(b1*(sxz[j][i][k]-sxz[j][i][k-1])+
                                            b2*(sxz[j][i][k+1]-sxz[j][i][k-2])+
                                            b3*(sxz[j][i][k+2]-sxz[j][i][k-3]));
                                
                                syy_y = dy*(b1*(syy[j+1][i][k]-syy[j][i][k])+
                                            b2*(syy[j+2][i][k]-syy[j-1][i][k])+
                                            b3*(syy[j+3][i][k]-syy[j-2][i][k]));
                                
                                sxy_x = dx*(b1*(sxy[j][i][k]-sxy[j][i-1][k])+
                                            b2*(sxy[j][i+1][k]-sxy[j][i-2][k])+
                                            b3*(sxy[j][i+2][k]-sxy[j][i-3][k]));
                                
                                syz_z = dz*(b1*(syz[j][i][k]-syz[j][i][k-1])+
                                            b2*(syz[j][i][k+1]-syz[j][i][k-2])+
                                            b3*(syz[j][i][k+2]-syz[j][i][k-3]));
                                
                                szz_z = dz*(b1*(szz[j][i][k+1]-szz[j][i][k])+
                                            b2*(szz[j][i][k+2]-szz[j][i][k-1])+
                                            b3*(szz[j][i][k+3]-szz[j][i][k-2]));
                                
                                sxz_x = dx*(b1*(sxz[j][i][k]-sxz[j][i-1][k])+
                                            b2*(sxz[j][i+1][k]-sxz[j][i-2][k])+
                                            b3*(sxz[j][i+2][k]-sxz[j][i-3][k]));
                                
                                
                                syz_y = dy*(b1*(syz[j][i][k]-syz[j-1][i][k])+
                                            b2*(syz[j+1][i][k]-syz[j-2][i][k])+
                                            b3*(syz[j+2][i][k]-syz[j-3][i][k]));
                                
                                
                                *(svx_j_i+k)=sxx_x + sxy_y +sxz_z;
                                *(svy_j_i+k)=syy_y + sxy_x + syz_z;
                                *(svz_j_i+k)=szz_z + sxz_x + syz_y;
                                
                                
                                /* updating components of particle velocities */
                                vx[j][i][k]+= ((c1*(*(svx_j_i+k))+c2*(*(svx_j_i_2+k))+c3*(*(svx_j_i_3+k)))/rip[j][i][k]);
                                vy[j][i][k]+= ((c1*(*(svy_j_i+k))+c2*(*(svy_j_i_2+k))+c3*(*(svy_j_i_3+k)))/rjp[j][i][k]);
                                vz[j][i][k]+= ((c1*(*(svz_j_i+k))+c2*(*(svz_j_i_2+k))+c3*(*(svz_j_i_3+k)))/rkp[j][i][k]);
                                
                            }
                        }
                    }
                    
                    break;
                    
                case 8 :
                    
                    dx=DT/DX;
                    dy=DT/DY;
                    dz=DT/DZ;
                    
                    b1=1225.0/1024.0; b2=-245.0/3072.0; b3=49.0/5120.0; b4=-5.0/7168.0; /* Taylor coefficients*/
                    if(FDCOEFF==2){
                        b1=1.2257; b2=-0.099537; b3=0.018063; b4=-0.0026274;} /* Holberg coefficients E=0.1 %*/
                    
                    for (j=ny1;j<=ny2;j++){
                        svx_j=*(svx+j);svy_j=*(svy+j);svz_j=*(svz+j);
                        svx_j_2=*(svx_2+j);svy_j_2=*(svy_2+j);svz_j_2=*(svz_2+j);
                        svx_j_3=*(svx_3+j);svy_j_3=*(svy_3+j);svz_j_3=*(svz_3+j);
                        
                        
                        for (i=nx1;i<=nx2;i++){
                            svx_j_i=*(svx_j+i);svy_j_i=*(svy_j+i);svz_j_i=*(svz_j+i);
                            svx_j_i_2=*(svx_j_2+i);svy_j_i_2=*(svy_j_2+i);svz_j_i_2=*(svz_j_2+i);
                            svx_j_i_3=*(svx_j_3+i);svy_j_i_3=*(svy_j_3+i);svz_j_i_3=*(svz_j_3+i);
                            
                            for (k=nz1;k<=nz2;k++){
                                sxx_x = dx*(b1*(sxx[j][i+1][k]-sxx[j][i][k])+
                                            b2*(sxx[j][i+2][k]-sxx[j][i-1][k])+
                                            b3*(sxx[j][i+3][k]-sxx[j][i-2][k])+
                                            b4*(sxx[j][i+4][k]-sxx[j][i-3][k]));
                                
                                sxy_y = dy*(b1*(sxy[j][i][k]-sxy[j-1][i][k])+
                                            b2*(sxy[j+1][i][k]-sxy[j-2][i][k])+
                                            b3*(sxy[j+2][i][k]-sxy[j-3][i][k])+
                                            b4*(sxy[j+3][i][k]-sxy[j-4][i][k]));
                                
                                sxz_z = dz*(b1*(sxz[j][i][k]-sxz[j][i][k-1])+
                                            b2*(sxz[j][i][k+1]-sxz[j][i][k-2])+
                                            b3*(sxz[j][i][k+2]-sxz[j][i][k-3])+
                                            b4*(sxz[j][i][k+3]-sxz[j][i][k-4]));
                                
                                syy_y = dy*(b1*(syy[j+1][i][k]-syy[j][i][k])+
                                            b2*(syy[j+2][i][k]-syy[j-1][i][k])+
                                            b3*(syy[j+3][i][k]-syy[j-2][i][k])+
                                            b4*(syy[j+4][i][k]-syy[j-3][i][k]));
                                
                                sxy_x = dx*(b1*(sxy[j][i][k]-sxy[j][i-1][k])+
                                            b2*(sxy[j][i+1][k]-sxy[j][i-2][k])+
                                            b3*(sxy[j][i+2][k]-sxy[j][i-3][k])+
                                            b4*(sxy[j][i+3][k]-sxy[j][i-4][k]));
                                
                                syz_z = dz*(b1*(syz[j][i][k]-syz[j][i][k-1])+
                                            b2*(syz[j][i][k+1]-syz[j][i][k-2])+
                                            b3*(syz[j][i][k+2]-syz[j][i][k-3])+
                                            b4*(syz[j][i][k+3]-syz[j][i][k-4]));
                                
                                szz_z = dz*(b1*(szz[j][i][k+1]-szz[j][i][k])+
                                            b2*(szz[j][i][k+2]-szz[j][i][k-1])+
                                            b3*(szz[j][i][k+3]-szz[j][i][k-2])+
                                            b4*(szz[j][i][k+4]-szz[j][i][k-3]));
                                
                                sxz_x = dx*(b1*(sxz[j][i][k]-sxz[j][i-1][k])+
                                            b2*(sxz[j][i+1][k]-sxz[j][i-2][k])+
                                            b3*(sxz[j][i+2][k]-sxz[j][i-3][k])+
                                            b4*(sxz[j][i+3][k]-sxz[j][i-4][k]));
                                
                                
                                syz_y = dy*(b1*(syz[j][i][k]-syz[j-1][i][k])+
                                            b2*(syz[j+1][i][k]-syz[j-2][i][k])+
                                            b3*(syz[j+2][i][k]-syz[j-3][i][k])+
                                            b4*(syz[j+3][i][k]-syz[j-4][i][k]));
                                
                                
                                *(svx_j_i+k)=sxx_x + sxy_y +sxz_z;
                                *(svy_j_i+k)=syy_y + sxy_x + syz_z;
                                *(svz_j_i+k)=szz_z + sxz_x + syz_y;
                                
                                
                                /* updating components of particle velocities */
                                vx[j][i][k]+= ((c1*(*(svx_j_i+k))+c2*(*(svx_j_i_2+k))+c3*(*(svx_j_i_3+k)))/rip[j][i][k]);
                                vy[j][i][k]+= ((c1*(*(svy_j_i+k))+c2*(*(svy_j_i_2+k))+c3*(*(svy_j_i_3+k)))/rjp[j][i][k]);
                                vz[j][i][k]+= ((c1*(*(svz_j_i+k))+c2*(*(svz_j_i_2+k))+c3*(*(svz_j_i_3+k)))/rkp[j][i][k]);
                                
                            }
                        }
                    }
                    
                    break;
                    
                case 10 :
                    
                    dx=DT/DX;
                    dy=DT/DY;
                    dz=DT/DZ;
                    
                    
                    b1=19845.0/16384.0; b2=-735.0/8192.0; b3=567.0/40960.0; b4=-405.0/229376.0; b5=35.0/294912.0; /* Taylor Coefficients*/
                    if(FDCOEFF==2){
                        b1=1.2415; b2=-0.11231; b3=0.026191; b4=-0.0064682; b5=0.001191;} /* Holberg coefficients E=0.1 %*/
                    
                    for (j=ny1;j<=ny2;j++){
                        svx_j=*(svx+j);svy_j=*(svy+j);svz_j=*(svz+j);
                        svx_j_2=*(svx_2+j);svy_j_2=*(svy_2+j);svz_j_2=*(svz_2+j);
                        svx_j_3=*(svx_3+j);svy_j_3=*(svy_3+j);svz_j_3=*(svz_3+j);
                        
                        
                        for (i=nx1;i<=nx2;i++){
                            svx_j_i=*(svx_j+i);svy_j_i=*(svy_j+i);svz_j_i=*(svz_j+i);
                            svx_j_i_2=*(svx_j_2+i);svy_j_i_2=*(svy_j_2+i);svz_j_i_2=*(svz_j_2+i);
                            svx_j_i_3=*(svx_j_3+i);svy_j_i_3=*(svy_j_3+i);svz_j_i_3=*(svz_j_3+i);
                            
                            for (k=nz1;k<=nz2;k++){
                                sxx_x = dx*(b1*(sxx[j][i+1][k]-sxx[j][i][k])+
                                            b2*(sxx[j][i+2][k]-sxx[j][i-1][k])+
                                            b3*(sxx[j][i+3][k]-sxx[j][i-2][k])+
                                            b4*(sxx[j][i+4][k]-sxx[j][i-3][k])+
                                            b5*(sxx[j][i+5][k]-sxx[j][i-4][k]));
                                
                                sxy_y = dy*(b1*(sxy[j][i][k]-sxy[j-1][i][k])+
                                            b2*(sxy[j+1][i][k]-sxy[j-2][i][k])+
                                            b3*(sxy[j+2][i][k]-sxy[j-3][i][k])+
                                            b4*(sxy[j+3][i][k]-sxy[j-4][i][k])+
                                            b5*(sxy[j+4][i][k]-sxy[j-5][i][k]));
                                
                                sxz_z = dz*(b1*(sxz[j][i][k]-sxz[j][i][k-1])+
                                            b2*(sxz[j][i][k+1]-sxz[j][i][k-2])+
                                            b3*(sxz[j][i][k+2]-sxz[j][i][k-3])+
                                            b4*(sxz[j][i][k+3]-sxz[j][i][k-4])+
                                            b5*(sxz[j][i][k+4]-sxz[j][i][k-5]));
                                
                                syy_y = dy*(b1*(syy[j+1][i][k]-syy[j][i][k])+
                                            b2*(syy[j+2][i][k]-syy[j-1][i][k])+
                                            b3*(syy[j+3][i][k]-syy[j-2][i][k])+
                                            b4*(syy[j+4][i][k]-syy[j-3][i][k])+
                                            b5*(syy[j+5][i][k]-syy[j-4][i][k]));
                                
                                sxy_x = dx*(b1*(sxy[j][i][k]-sxy[j][i-1][k])+
                                            b2*(sxy[j][i+1][k]-sxy[j][i-2][k])+
                                            b3*(sxy[j][i+2][k]-sxy[j][i-3][k])+
                                            b4*(sxy[j][i+3][k]-sxy[j][i-4][k])+
                                            b5*(sxy[j][i+4][k]-sxy[j][i-5][k]));
                                
                                syz_z = dz*(b1*(syz[j][i][k]-syz[j][i][k-1])+
                                            b2*(syz[j][i][k+1]-syz[j][i][k-2])+
                                            b3*(syz[j][i][k+2]-syz[j][i][k-3])+
                                            b4*(syz[j][i][k+3]-syz[j][i][k-4])+
                                            b5*(syz[j][i][k+4]-syz[j][i][k-5]));
                                
                                szz_z = dz*(b1*(szz[j][i][k+1]-szz[j][i][k])+
                                            b2*(szz[j][i][k+2]-szz[j][i][k-1])+
                                            b3*(szz[j][i][k+3]-szz[j][i][k-2])+
                                            b4*(szz[j][i][k+4]-szz[j][i][k-3])+
                                            b5*(szz[j][i][k+5]-szz[j][i][k-4]));
                                
                                sxz_x = dx*(b1*(sxz[j][i][k]-sxz[j][i-1][k])+
                                            b2*(sxz[j][i+1][k]-sxz[j][i-2][k])+
                                            b3*(sxz[j][i+2][k]-sxz[j][i-3][k])+
                                            b4*(sxz[j][i+3][k]-sxz[j][i-4][k])+
                                            b5*(sxz[j][i+4][k]-sxz[j][i-5][k]));
                                
                                
                                syz_y = dy*(b1*(syz[j][i][k]-syz[j-1][i][k])+
                                            b2*(syz[j+1][i][k]-syz[j-2][i][k])+
                                            b3*(syz[j+2][i][k]-syz[j-3][i][k])+
                                            b4*(syz[j+3][i][k]-syz[j-4][i][k])+
                                            b5*(syz[j+4][i][k]-syz[j-5][i][k]));
                                
                                *(svx_j_i+k)=sxx_x + sxy_y +sxz_z;
                                *(svy_j_i+k)=syy_y + sxy_x + syz_z;
                                *(svz_j_i+k)=szz_z + sxz_x + syz_y;
                                
                                
                                /* updating components of particle velocities */
                                vx[j][i][k]+= ((c1*(*(svx_j_i+k))+c2*(*(svx_j_i_2+k))+c3*(*(svx_j_i_3+k)))/rip[j][i][k]);
                                vy[j][i][k]+= ((c1*(*(svy_j_i+k))+c2*(*(svy_j_i_2+k))+c3*(*(svy_j_i_3+k)))/rjp[j][i][k]);
                                vz[j][i][k]+= ((c1*(*(svz_j_i+k))+c2*(*(svz_j_i_2+k))+c3*(*(svz_j_i_3+k)))/rkp[j][i][k]);
                                
                            }
                        }
                    }
                    
                    break;
                    
                case 12 :
                    dx=DT/DX;
                    dy=DT/DY;
                    dz=DT/DZ;
                    
                    /* Taylor coefficients */
                    b1=160083.0/131072.0; b2=-12705.0/131072.0; b3=22869.0/1310720.0;
                    b4=-5445.0/1835008.0; b5=847.0/2359296.0; b6=-63.0/2883584;
                    
                    /* Holberg coefficients E=0.1 %*/
                    if(FDCOEFF==2){
                        b1=1.2508; b2=-0.12034; b3=0.032131; b4=-0.010142; b5=0.0029857; b6=-0.00066667;}
                    
                    
                    for (j=ny1;j<=ny2;j++){
                        svx_j=*(svx+j);svy_j=*(svy+j);svz_j=*(svz+j);
                        svx_j_2=*(svx_2+j);svy_j_2=*(svy_2+j);svz_j_2=*(svz_2+j);
                        svx_j_3=*(svx_3+j);svy_j_3=*(svy_3+j);svz_j_3=*(svz_3+j);
                        
                        
                        for (i=nx1;i<=nx2;i++){
                            svx_j_i=*(svx_j+i);svy_j_i=*(svy_j+i);svz_j_i=*(svz_j+i);
                            svx_j_i_2=*(svx_j_2+i);svy_j_i_2=*(svy_j_2+i);svz_j_i_2=*(svz_j_2+i);
                            svx_j_i_3=*(svx_j_3+i);svy_j_i_3=*(svy_j_3+i);svz_j_i_3=*(svz_j_3+i);
                            
                            for (k=nz1;k<=nz2;k++){
                                
                                sxx_x = dx*(b1*(sxx[j][i+1][k]-sxx[j][i][k])+
                                            b2*(sxx[j][i+2][k]-sxx[j][i-1][k])+
                                            b3*(sxx[j][i+3][k]-sxx[j][i-2][k])+
                                            b4*(sxx[j][i+4][k]-sxx[j][i-3][k])+
                                            b5*(sxx[j][i+5][k]-sxx[j][i-4][k])+
                                            b6*(sxx[j][i+6][k]-sxx[j][i-5][k]));
                                
                                sxy_y = dy*(b1*(sxy[j][i][k]-sxy[j-1][i][k])+
                                            b2*(sxy[j+1][i][k]-sxy[j-2][i][k])+
                                            b3*(sxy[j+2][i][k]-sxy[j-3][i][k])+
                                            b4*(sxy[j+3][i][k]-sxy[j-4][i][k])+
                                            b5*(sxy[j+4][i][k]-sxy[j-5][i][k])+
                                            b6*(sxy[j+5][i][k]-sxy[j-6][i][k]));
                                
                                sxz_z = dz*(b1*(sxy[j][i][k]-sxy[j][i][k-1])+
                                            b2*(sxy[j][i][k+1]-sxy[j][i][k-2])+
                                            b3*(sxy[j][i][k+2]-sxy[j][i][k-3])+
                                            b4*(sxy[j][i][k+3]-sxy[j][i][k-4])+
                                            b5*(sxy[j][i][k+4]-sxy[j][i][k-5])+
                                            b6*(sxy[j][i][k+5]-sxy[j][i][k-6]));
                                
                                syy_y = dy*(b1*(syy[j+1][i][k]-syy[j][i][k])+
                                            b2*(syy[j+2][i][k]-syy[j-1][i][k])+
                                            b3*(syy[j+3][i][k]-syy[j-2][i][k])+
                                            b4*(syy[j+4][i][k]-syy[j-3][i][k])+
                                            b5*(syy[j+5][i][k]-syy[j-4][i][k])+
                                            b6*(syy[j+6][i][k]-syy[j-5][i][k]));
                                
                                sxy_x = dx*(b1*(sxy[j][i][k]-sxy[j][i-1][k])+
                                            b2*(sxy[j][i+1][k]-sxy[j][i-2][k])+
                                            b3*(sxy[j][i+2][k]-sxy[j][i-3][k])+
                                            b4*(sxy[j][i+3][k]-sxy[j][i-4][k])+
                                            b5*(sxy[j][i+4][k]-sxy[j][i-5][k])+
                                            b6*(sxy[j][i+5][k]-sxy[j][i-6][k]));
                                
                                syz_z = dz*(b1*(syz[j][i][k]-syz[j][i][k-1])+
                                            b2*(syz[j][i][k+1]-syz[j][i][k-2])+
                                            b3*(syz[j][i][k+2]-syz[j][i][k-3])+
                                            b4*(syz[j][i][k+3]-syz[j][i][k-4])+
                                            b5*(syz[j][i][k+4]-syz[j][i][k-5])+
                                            b6*(syz[j][i][k+5]-syz[j][i][k-6]));
                                
                                szz_z = dz*(b1*(szz[j][i][k+1]-szz[j][i][k])+
                                            b2*(szz[j][i][k+2]-szz[j][i][k-1])+
                                            b3*(szz[j][i][k+3]-szz[j][i][k-2])+
                                            b4*(szz[j][i][k+4]-szz[j][i][k-3])+
                                            b5*(szz[j][i][k+5]-szz[j][i][k-4])+
                                            b6*(szz[j][i][k+6]-szz[j][i][k-5]));
                                
                                sxz_x = dx*(b1*(sxz[j][i][k]-sxz[j][i-1][k])+
                                            b2*(sxz[j][i+1][k]-sxz[j][i-2][k])+
                                            b3*(sxz[j][i+2][k]-sxz[j][i-3][k])+
                                            b4*(sxz[j][i+3][k]-sxz[j][i-4][k])+
                                            b5*(sxz[j][i+4][k]-sxz[j][i-5][k])+
                                            b6*(sxz[j][i+5][k]-sxz[j][i-6][k]));
                                
                                
                                syz_y = dy*(b1*(syz[j][i][k]-syz[j-1][i][k])+
                                            b2*(syz[j+1][i][k]-syz[j-2][i][k])+
                                            b3*(syz[j+2][i][k]-syz[j-3][i][k])+
                                            b4*(syz[j+3][i][k]-syz[j-4][i][k])+
                                            b5*(syz[j+4][i][k]-syz[j-5][i][k])+
                                            b6*(syz[j+5][i][k]-syz[j-6][i][k]));
                                
                                *(svx_j_i+k)=sxx_x + sxy_y +sxz_z;
                                *(svy_j_i+k)=syy_y + sxy_x + syz_z;
                                *(svz_j_i+k)=szz_z + sxz_x + syz_y;
                                
                                
                                /* updating components of particle velocities */
                                vx[j][i][k]+= ((c1*(*(svx_j_i+k))+c2*(*(svx_j_i_2+k))+c3*(*(svx_j_i_3+k)))/rip[j][i][k]);
                                vy[j][i][k]+= ((c1*(*(svy_j_i+k))+c2*(*(svy_j_i_2+k))+c3*(*(svy_j_i_3+k)))/rjp[j][i][k]);
                                vz[j][i][k]+= ((c1*(*(svz_j_i+k))+c2*(*(svz_j_i_2+k))+c3*(*(svz_j_i_3+k)))/rkp[j][i][k]);
                                
                            }
                        }
                    }
                    
                    break;
                    
            }
            break; /* break for FDORDER_TIME=3 */
            
            
        case 4:
            
            c1=13.0/12.0; c2=-5.0/24.0; c3=1.0/6.0; c4=-1.0/24.0;
            
            switch (FDORDER){
                    
                case 2 :
                    
                    dx=DT/DX;
                    dy=DT/DY;
                    dz=DT/DZ;
                    
                    b1=1.0; /* Taylor coefficients*/
                    if(FDCOEFF==2){
                        b1=1.00100; } /* Holberg coefficients E=0.1 %*/
                    
                    for (j=ny1;j<=ny2;j++){
                        svx_j=*(svx+j);svy_j=*(svy+j);svz_j=*(svz+j);
                        svx_j_2=*(svx_2+j);svy_j_2=*(svy_2+j);svz_j_2=*(svz_2+j);
                        svx_j_3=*(svx_3+j);svy_j_3=*(svy_3+j);svz_j_3=*(svz_3+j);
                        svx_j_4=*(svx_4+j);svy_j_4=*(svy_4+j);svz_j_4=*(svz_4+j);
                        
                        for (i=nx1;i<=nx2;i++){
                            svx_j_i=*(svx_j+i);svy_j_i=*(svy_j+i);svz_j_i=*(svz_j+i);
                            svx_j_i_2=*(svx_j_2+i);svy_j_i_2=*(svy_j_2+i);svz_j_i_2=*(svz_j_2+i);
                            svx_j_i_3=*(svx_j_3+i);svy_j_i_3=*(svy_j_3+i);svz_j_i_3=*(svz_j_3+i);
                            svx_j_i_4=*(svx_j_4+i);svy_j_i_4=*(svy_j_4+i);svz_j_i_4=*(svz_j_4+i);
                            
                            for (k=nz1;k<=nz2;k++){
                                
                                sxx_x = dx*b1*(sxx[j][i+1][k]-sxx[j][i][k]);
                                sxy_y = dy*b1*(sxy[j][i][k]-sxy[j-1][i][k]);
                                sxz_z = dz*b1*(sxz[j][i][k]-sxz[j][i][k-1]); /* backward operator */
                                
                                syy_y = dy*b1*(syy[j+1][i][k]-syy[j][i][k]);
                                sxy_x = dx*b1*(sxy[j][i][k]-sxy[j][i-1][k]);
                                syz_z = dz*b1*(syz[j][i][k]-syz[j][i][k-1]);
                                
                                szz_z = dz*b1*(szz[j][i][k+1]-szz[j][i][k]);
                                sxz_x = dx*b1*(sxz[j][i][k]-sxz[j][i-1][k]);
                                syz_y = dy*b1*(syz[j][i][k]-syz[j-1][i][k]);
                                
                                *(svx_j_i+k)=sxx_x + sxy_y +sxz_z;
                                *(svy_j_i+k)=syy_y + sxy_x + syz_z;
                                *(svz_j_i+k)=szz_z + sxz_x + syz_y;
                                
                                
                                /* updating components of particle velocities */
                                vx[j][i][k]+= ((c1*(*(svx_j_i+k))+c2*(*(svx_j_i_2+k))+c3*(*(svx_j_i_3+k))+c4*(*(svx_j_i_4+k)))/rip[j][i][k]);
                                vy[j][i][k]+= ((c1*(*(svy_j_i+k))+c2*(*(svy_j_i_2+k))+c3*(*(svy_j_i_3+k))+c4*(*(svy_j_i_4+k)))/rjp[j][i][k]);
                                vz[j][i][k]+= ((c1*(*(svz_j_i+k))+c2*(*(svz_j_i_2+k))+c3*(*(svz_j_i_3+k))+c4*(*(svz_j_i_4+k)))/rkp[j][i][k]);
                            
                            }
                        }
                    }
                    
                    break;
                    
                case 4 :
                    
                    dx=DT/DX;
                    dy=DT/DY;
                    dz=DT/DZ;
                    
                    b1=9.0/8.0; b2=-1.0/24.0; /* Taylor coefficients*/
                    if(FDCOEFF==2){
                        b1=1.1382; b2=-0.046414;} /* Holberg coefficients E=0.1 %*/
                    for (j=ny1;j<=ny2;j++){
                        svx_j=*(svx+j);svy_j=*(svy+j);svz_j=*(svz+j);
                        svx_j_2=*(svx_2+j);svy_j_2=*(svy_2+j);svz_j_2=*(svz_2+j);
                        svx_j_3=*(svx_3+j);svy_j_3=*(svy_3+j);svz_j_3=*(svz_3+j);
                        svx_j_4=*(svx_4+j);svy_j_4=*(svy_4+j);svz_j_4=*(svz_4+j);
                        
                        for (i=nx1;i<=nx2;i++){
                            svx_j_i=*(svx_j+i);svy_j_i=*(svy_j+i);svz_j_i=*(svz_j+i);
                            svx_j_i_2=*(svx_j_2+i);svy_j_i_2=*(svy_j_2+i);svz_j_i_2=*(svz_j_2+i);
                            svx_j_i_3=*(svx_j_3+i);svy_j_i_3=*(svy_j_3+i);svz_j_i_3=*(svz_j_3+i);
                            svx_j_i_4=*(svx_j_4+i);svy_j_i_4=*(svy_j_4+i);svz_j_i_4=*(svz_j_4+i);
                            
                            for (k=nz1;k<=nz2;k++){
                                sxx_x = dx*(b1*(sxx[j][i+1][k]-sxx[j][i][k])+b2*(sxx[j][i+2][k]-sxx[j][i-1][k]));
                                sxy_y = dy*(b1*(sxy[j][i][k]-sxy[j-1][i][k])+b2*(sxy[j+1][i][k]-sxy[j-2][i][k]));
                                sxz_z = dz*(b1*(sxz[j][i][k]-sxz[j][i][k-1])+b2*(sxz[j][i][k+1]-sxz[j][i][k-2]));
                                
                                syy_y = dy*(b1*(syy[j+1][i][k]-syy[j][i][k])+b2*(syy[j+2][i][k]-syy[j-1][i][k]));
                                sxy_x = dx*(b1*(sxy[j][i][k]-sxy[j][i-1][k])+b2*(sxy[j][i+1][k]-sxy[j][i-2][k]));
                                syz_z = dz*(b1*(syz[j][i][k]-syz[j][i][k-1])+b2*(syz[j][i][k+1]-syz[j][i][k-2]));
                                
                                szz_z = dz*(b1*(szz[j][i][k+1]-szz[j][i][k])+b2*(szz[j][i][k+2]-szz[j][i][k-1]));
                                sxz_x = dx*(b1*(sxz[j][i][k]-sxz[j][i-1][k])+b2*(sxz[j][i+1][k]-sxz[j][i-2][k]));
                                syz_y = dy*(b1*(syz[j][i][k]-syz[j-1][i][k])+b2*(syz[j+1][i][k]-syz[j-2][i][k]));
                                
                                
                                *(svx_j_i+k)=sxx_x + sxy_y +sxz_z;
                                *(svy_j_i+k)=syy_y + sxy_x + syz_z;
                                *(svz_j_i+k)=szz_z + sxz_x + syz_y;
                                
                                
                                /* updating components of particle velocities */
                                vx[j][i][k]+= ((c1*(*(svx_j_i+k))+c2*(*(svx_j_i_2+k))+c3*(*(svx_j_i_3+k))+c4*(*(svx_j_i_4+k)))/rip[j][i][k]);
                                vy[j][i][k]+= ((c1*(*(svy_j_i+k))+c2*(*(svy_j_i_2+k))+c3*(*(svy_j_i_3+k))+c4*(*(svy_j_i_4+k)))/rjp[j][i][k]);
                                vz[j][i][k]+= ((c1*(*(svz_j_i+k))+c2*(*(svz_j_i_2+k))+c3*(*(svz_j_i_3+k))+c4*(*(svz_j_i_4+k)))/rkp[j][i][k]);
                                
                            }
                        }
                    }
                    
                    break;
                    
                case 6 :
                    
                    dx=DT/DX;
                    dy=DT/DY;
                    dz=DT/DZ;
                    
                    b1=75.0/64.0; b2=-25.0/384.0; b3=3.0/640.0; /* Taylor coefficients*/
                    if(FDCOEFF==2){
                        b1=1.1965; b2=-0.078804; b3=0.0081781;}   /* Holberg coefficients E=0.1 %*/
                    for (j=ny1;j<=ny2;j++){
                        svx_j=*(svx+j);svy_j=*(svy+j);svz_j=*(svz+j);
                        svx_j_2=*(svx_2+j);svy_j_2=*(svy_2+j);svz_j_2=*(svz_2+j);
                        svx_j_3=*(svx_3+j);svy_j_3=*(svy_3+j);svz_j_3=*(svz_3+j);
                        svx_j_4=*(svx_4+j);svy_j_4=*(svy_4+j);svz_j_4=*(svz_4+j);
                        
                        for (i=nx1;i<=nx2;i++){
                            svx_j_i=*(svx_j+i);svy_j_i=*(svy_j+i);svz_j_i=*(svz_j+i);
                            svx_j_i_2=*(svx_j_2+i);svy_j_i_2=*(svy_j_2+i);svz_j_i_2=*(svz_j_2+i);
                            svx_j_i_3=*(svx_j_3+i);svy_j_i_3=*(svy_j_3+i);svz_j_i_3=*(svz_j_3+i);
                            svx_j_i_4=*(svx_j_4+i);svy_j_i_4=*(svy_j_4+i);svz_j_i_4=*(svz_j_4+i);
                            
                            for (k=nz1;k<=nz2;k++){
                                sxx_x = dx*(b1*(sxx[j][i+1][k]-sxx[j][i][k])+
                                            b2*(sxx[j][i+2][k]-sxx[j][i-1][k])+
                                            b3*(sxx[j][i+3][k]-sxx[j][i-2][k]));
                                
                                sxy_y = dy*(b1*(sxy[j][i][k]-sxy[j-1][i][k])+
                                            b2*(sxy[j+1][i][k]-sxy[j-2][i][k])+
                                            b3*(sxy[j+2][i][k]-sxy[j-3][i][k]));
                                
                                sxz_z = dz*(b1*(sxz[j][i][k]-sxz[j][i][k-1])+
                                            b2*(sxz[j][i][k+1]-sxz[j][i][k-2])+
                                            b3*(sxz[j][i][k+2]-sxz[j][i][k-3]));
                                
                                syy_y = dy*(b1*(syy[j+1][i][k]-syy[j][i][k])+
                                            b2*(syy[j+2][i][k]-syy[j-1][i][k])+
                                            b3*(syy[j+3][i][k]-syy[j-2][i][k]));
                                
                                sxy_x = dx*(b1*(sxy[j][i][k]-sxy[j][i-1][k])+
                                            b2*(sxy[j][i+1][k]-sxy[j][i-2][k])+
                                            b3*(sxy[j][i+2][k]-sxy[j][i-3][k]));
                                
                                syz_z = dz*(b1*(syz[j][i][k]-syz[j][i][k-1])+
                                            b2*(syz[j][i][k+1]-syz[j][i][k-2])+
                                            b3*(syz[j][i][k+2]-syz[j][i][k-3]));
                                
                                szz_z = dz*(b1*(szz[j][i][k+1]-szz[j][i][k])+
                                            b2*(szz[j][i][k+2]-szz[j][i][k-1])+
                                            b3*(szz[j][i][k+3]-szz[j][i][k-2]));
                                
                                sxz_x = dx*(b1*(sxz[j][i][k]-sxz[j][i-1][k])+
                                            b2*(sxz[j][i+1][k]-sxz[j][i-2][k])+
                                            b3*(sxz[j][i+2][k]-sxz[j][i-3][k]));
                                
                                
                                syz_y = dy*(b1*(syz[j][i][k]-syz[j-1][i][k])+
                                            b2*(syz[j+1][i][k]-syz[j-2][i][k])+
                                            b3*(syz[j+2][i][k]-syz[j-3][i][k]));
                                

                                
                                *(svx_j_i+k)=sxx_x + sxy_y +sxz_z;
                                *(svy_j_i+k)=syy_y + sxy_x + syz_z;
                                *(svz_j_i+k)=szz_z + sxz_x + syz_y;
                                
                                
                                /* updating components of particle velocities */
                                vx[j][i][k]+= ((c1*(*(svx_j_i+k))+c2*(*(svx_j_i_2+k))+c3*(*(svx_j_i_3+k))+c4*(*(svx_j_i_4+k)))/rip[j][i][k]);
                                vy[j][i][k]+= ((c1*(*(svy_j_i+k))+c2*(*(svy_j_i_2+k))+c3*(*(svy_j_i_3+k))+c4*(*(svy_j_i_4+k)))/rjp[j][i][k]);
                                vz[j][i][k]+= ((c1*(*(svz_j_i+k))+c2*(*(svz_j_i_2+k))+c3*(*(svz_j_i_3+k))+c4*(*(svz_j_i_4+k)))/rkp[j][i][k]);
                            }
                        }
                    }
                    
                    break;
                    
                case 8 :
                    
                    dx=DT/DX;
                    dy=DT/DY;
                    dz=DT/DZ;
                    
                    b1=1225.0/1024.0; b2=-245.0/3072.0; b3=49.0/5120.0; b4=-5.0/7168.0; /* Taylor coefficients*/
                    if(FDCOEFF==2){
                        b1=1.2257; b2=-0.099537; b3=0.018063; b4=-0.0026274;} /* Holberg coefficients E=0.1 %*/
                    
                    for (j=ny1;j<=ny2;j++){
                        svx_j=*(svx+j);svy_j=*(svy+j);svz_j=*(svz+j);
                        svx_j_2=*(svx_2+j);svy_j_2=*(svy_2+j);svz_j_2=*(svz_2+j);
                        svx_j_3=*(svx_3+j);svy_j_3=*(svy_3+j);svz_j_3=*(svz_3+j);
                        svx_j_4=*(svx_4+j);svy_j_4=*(svy_4+j);svz_j_4=*(svz_4+j);
                        
                        for (i=nx1;i<=nx2;i++){
                            svx_j_i=*(svx_j+i);svy_j_i=*(svy_j+i);svz_j_i=*(svz_j+i);
                            svx_j_i_2=*(svx_j_2+i);svy_j_i_2=*(svy_j_2+i);svz_j_i_2=*(svz_j_2+i);
                            svx_j_i_3=*(svx_j_3+i);svy_j_i_3=*(svy_j_3+i);svz_j_i_3=*(svz_j_3+i);
                            svx_j_i_4=*(svx_j_4+i);svy_j_i_4=*(svy_j_4+i);svz_j_i_4=*(svz_j_4+i);
                            
                            for (k=nz1;k<=nz2;k++){
                                sxx_x = dx*(b1*(sxx[j][i+1][k]-sxx[j][i][k])+
                                            b2*(sxx[j][i+2][k]-sxx[j][i-1][k])+
                                            b3*(sxx[j][i+3][k]-sxx[j][i-2][k])+
                                            b4*(sxx[j][i+4][k]-sxx[j][i-3][k]));
                                
                                sxy_y = dy*(b1*(sxy[j][i][k]-sxy[j-1][i][k])+
                                            b2*(sxy[j+1][i][k]-sxy[j-2][i][k])+
                                            b3*(sxy[j+2][i][k]-sxy[j-3][i][k])+
                                            b4*(sxy[j+3][i][k]-sxy[j-4][i][k]));
                                
                                sxz_z = dz*(b1*(sxz[j][i][k]-sxz[j][i][k-1])+
                                            b2*(sxz[j][i][k+1]-sxz[j][i][k-2])+
                                            b3*(sxz[j][i][k+2]-sxz[j][i][k-3])+
                                            b4*(sxz[j][i][k+3]-sxz[j][i][k-4]));
                                
                                syy_y = dy*(b1*(syy[j+1][i][k]-syy[j][i][k])+
                                            b2*(syy[j+2][i][k]-syy[j-1][i][k])+
                                            b3*(syy[j+3][i][k]-syy[j-2][i][k])+
                                            b4*(syy[j+4][i][k]-syy[j-3][i][k]));
                                
                                sxy_x = dx*(b1*(sxy[j][i][k]-sxy[j][i-1][k])+
                                            b2*(sxy[j][i+1][k]-sxy[j][i-2][k])+
                                            b3*(sxy[j][i+2][k]-sxy[j][i-3][k])+
                                            b4*(sxy[j][i+3][k]-sxy[j][i-4][k]));
                                
                                syz_z = dz*(b1*(syz[j][i][k]-syz[j][i][k-1])+
                                            b2*(syz[j][i][k+1]-syz[j][i][k-2])+
                                            b3*(syz[j][i][k+2]-syz[j][i][k-3])+
                                            b4*(syz[j][i][k+3]-syz[j][i][k-4]));
                                
                                szz_z = dz*(b1*(szz[j][i][k+1]-szz[j][i][k])+
                                            b2*(szz[j][i][k+2]-szz[j][i][k-1])+
                                            b3*(szz[j][i][k+3]-szz[j][i][k-2])+
                                            b4*(szz[j][i][k+4]-szz[j][i][k-3]));
                                
                                sxz_x = dx*(b1*(sxz[j][i][k]-sxz[j][i-1][k])+
                                            b2*(sxz[j][i+1][k]-sxz[j][i-2][k])+
                                            b3*(sxz[j][i+2][k]-sxz[j][i-3][k])+
                                            b4*(sxz[j][i+3][k]-sxz[j][i-4][k]));
                                
                                
                                syz_y = dy*(b1*(syz[j][i][k]-syz[j-1][i][k])+
                                            b2*(syz[j+1][i][k]-syz[j-2][i][k])+
                                            b3*(syz[j+2][i][k]-syz[j-3][i][k])+
                                            b4*(syz[j+3][i][k]-syz[j-4][i][k]));
                                
                                
                                *(svx_j_i+k)=sxx_x + sxy_y +sxz_z;
                                *(svy_j_i+k)=syy_y + sxy_x + syz_z;
                                *(svz_j_i+k)=szz_z + sxz_x + syz_y;
                                
                                
                                /* updating components of particle velocities */
                                vx[j][i][k]+= ((c1*(*(svx_j_i+k))+c2*(*(svx_j_i_2+k))+c3*(*(svx_j_i_3+k))+c4*(*(svx_j_i_4+k)))/rip[j][i][k]);
                                vy[j][i][k]+= ((c1*(*(svy_j_i+k))+c2*(*(svy_j_i_2+k))+c3*(*(svy_j_i_3+k))+c4*(*(svy_j_i_4+k)))/rjp[j][i][k]);
                                vz[j][i][k]+= ((c1*(*(svz_j_i+k))+c2*(*(svz_j_i_2+k))+c3*(*(svz_j_i_3+k))+c4*(*(svz_j_i_4+k)))/rkp[j][i][k]);
                                
                            }
                        }
                    }
                    
                    break;
                    
                case 10 :
                    
                    dx=DT/DX;
                    dy=DT/DY;
                    dz=DT/DZ;
                    
                    
                    b1=19845.0/16384.0; b2=-735.0/8192.0; b3=567.0/40960.0; b4=-405.0/229376.0; b5=35.0/294912.0; /* Taylor Coefficients*/
                    if(FDCOEFF==2){
                        b1=1.2415; b2=-0.11231; b3=0.026191; b4=-0.0064682; b5=0.001191;} /* Holberg coefficients E=0.1 %*/
                    
                    for (j=ny1;j<=ny2;j++){
                        svx_j=*(svx+j);svy_j=*(svy+j);svz_j=*(svz+j);
                        svx_j_2=*(svx_2+j);svy_j_2=*(svy_2+j);svz_j_2=*(svz_2+j);
                        svx_j_3=*(svx_3+j);svy_j_3=*(svy_3+j);svz_j_3=*(svz_3+j);
                        svx_j_4=*(svx_4+j);svy_j_4=*(svy_4+j);svz_j_4=*(svz_4+j);
                        
                        for (i=nx1;i<=nx2;i++){
                            svx_j_i=*(svx_j+i);svy_j_i=*(svy_j+i);svz_j_i=*(svz_j+i);
                            svx_j_i_2=*(svx_j_2+i);svy_j_i_2=*(svy_j_2+i);svz_j_i_2=*(svz_j_2+i);
                            svx_j_i_3=*(svx_j_3+i);svy_j_i_3=*(svy_j_3+i);svz_j_i_3=*(svz_j_3+i);
                            svx_j_i_4=*(svx_j_4+i);svy_j_i_4=*(svy_j_4+i);svz_j_i_4=*(svz_j_4+i);
                            
                            for (k=nz1;k<=nz2;k++){
                                
                                sxx_x = dx*(b1*(sxx[j][i+1][k]-sxx[j][i][k])+
                                            b2*(sxx[j][i+2][k]-sxx[j][i-1][k])+
                                            b3*(sxx[j][i+3][k]-sxx[j][i-2][k])+
                                            b4*(sxx[j][i+4][k]-sxx[j][i-3][k])+
                                            b5*(sxx[j][i+5][k]-sxx[j][i-4][k]));
                                
                                sxy_y = dy*(b1*(sxy[j][i][k]-sxy[j-1][i][k])+
                                            b2*(sxy[j+1][i][k]-sxy[j-2][i][k])+
                                            b3*(sxy[j+2][i][k]-sxy[j-3][i][k])+
                                            b4*(sxy[j+3][i][k]-sxy[j-4][i][k])+
                                            b5*(sxy[j+4][i][k]-sxy[j-5][i][k]));
                                
                                sxz_z = dz*(b1*(sxz[j][i][k]-sxz[j][i][k-1])+
                                            b2*(sxz[j][i][k+1]-sxz[j][i][k-2])+
                                            b3*(sxz[j][i][k+2]-sxz[j][i][k-3])+
                                            b4*(sxz[j][i][k+3]-sxz[j][i][k-4])+
                                            b5*(sxz[j][i][k+4]-sxz[j][i][k-5]));
                                
                                syy_y = dy*(b1*(syy[j+1][i][k]-syy[j][i][k])+
                                            b2*(syy[j+2][i][k]-syy[j-1][i][k])+
                                            b3*(syy[j+3][i][k]-syy[j-2][i][k])+
                                            b4*(syy[j+4][i][k]-syy[j-3][i][k])+
                                            b5*(syy[j+5][i][k]-syy[j-4][i][k]));
                                
                                sxy_x = dx*(b1*(sxy[j][i][k]-sxy[j][i-1][k])+
                                            b2*(sxy[j][i+1][k]-sxy[j][i-2][k])+
                                            b3*(sxy[j][i+2][k]-sxy[j][i-3][k])+
                                            b4*(sxy[j][i+3][k]-sxy[j][i-4][k])+
                                            b5*(sxy[j][i+4][k]-sxy[j][i-5][k]));
                                
                                syz_z = dz*(b1*(syz[j][i][k]-syz[j][i][k-1])+
                                            b2*(syz[j][i][k+1]-syz[j][i][k-2])+
                                            b3*(syz[j][i][k+2]-syz[j][i][k-3])+
                                            b4*(syz[j][i][k+3]-syz[j][i][k-4])+
                                            b5*(syz[j][i][k+4]-syz[j][i][k-5]));
                                
                                szz_z = dz*(b1*(szz[j][i][k+1]-szz[j][i][k])+
                                            b2*(szz[j][i][k+2]-szz[j][i][k-1])+
                                            b3*(szz[j][i][k+3]-szz[j][i][k-2])+
                                            b4*(szz[j][i][k+4]-szz[j][i][k-3])+
                                            b5*(szz[j][i][k+5]-szz[j][i][k-4]));
                                
                                sxz_x = dx*(b1*(sxz[j][i][k]-sxz[j][i-1][k])+
                                            b2*(sxz[j][i+1][k]-sxz[j][i-2][k])+
                                            b3*(sxz[j][i+2][k]-sxz[j][i-3][k])+
                                            b4*(sxz[j][i+3][k]-sxz[j][i-4][k])+
                                            b5*(sxz[j][i+4][k]-sxz[j][i-5][k]));
                                
                                
                                syz_y = dy*(b1*(syz[j][i][k]-syz[j-1][i][k])+
                                            b2*(syz[j+1][i][k]-syz[j-2][i][k])+
                                            b3*(syz[j+2][i][k]-syz[j-3][i][k])+
                                            b4*(syz[j+3][i][k]-syz[j-4][i][k])+
                                            b5*(syz[j+4][i][k]-syz[j-5][i][k]));
                                
                                
                                *(svx_j_i+k)=sxx_x + sxy_y +sxz_z;
                                *(svy_j_i+k)=syy_y + sxy_x + syz_z;
                                *(svz_j_i+k)=szz_z + sxz_x + syz_y;
                                
                                
                                /* updating components of particle velocities */
                                vx[j][i][k]+= ((c1*(*(svx_j_i+k))+c2*(*(svx_j_i_2+k))+c3*(*(svx_j_i_3+k))+c4*(*(svx_j_i_4+k)))/rip[j][i][k]);
                                vy[j][i][k]+= ((c1*(*(svy_j_i+k))+c2*(*(svy_j_i_2+k))+c3*(*(svy_j_i_3+k))+c4*(*(svy_j_i_4+k)))/rjp[j][i][k]);
                                vz[j][i][k]+= ((c1*(*(svz_j_i+k))+c2*(*(svz_j_i_2+k))+c3*(*(svz_j_i_3+k))+c4*(*(svz_j_i_4+k)))/rkp[j][i][k]);
                                
                            }
                        }
                    }
                    
                    break;
                    
                case 12 :
                    dx=DT/DX;
                    dy=DT/DY;
                    dz=DT/DZ;
                    
                    /* Taylor coefficients */
                    b1=160083.0/131072.0; b2=-12705.0/131072.0; b3=22869.0/1310720.0;
                    b4=-5445.0/1835008.0; b5=847.0/2359296.0; b6=-63.0/2883584;
                    
                    /* Holberg coefficients E=0.1 %*/
                    if(FDCOEFF==2){
                        b1=1.2508; b2=-0.12034; b3=0.032131; b4=-0.010142; b5=0.0029857; b6=-0.00066667;}
                    
                    for (j=ny1;j<=ny2;j++){
                        svx_j=*(svx+j);svy_j=*(svy+j);svz_j=*(svz+j);
                        svx_j_2=*(svx_2+j);svy_j_2=*(svy_2+j);svz_j_2=*(svz_2+j);
                        svx_j_3=*(svx_3+j);svy_j_3=*(svy_3+j);svz_j_3=*(svz_3+j);
                        svx_j_4=*(svx_4+j);svy_j_4=*(svy_4+j);svz_j_4=*(svz_4+j);
                        
                        for (i=nx1;i<=nx2;i++){
                            svx_j_i=*(svx_j+i);svy_j_i=*(svy_j+i);svz_j_i=*(svz_j+i);
                            svx_j_i_2=*(svx_j_2+i);svy_j_i_2=*(svy_j_2+i);svz_j_i_2=*(svz_j_2+i);
                            svx_j_i_3=*(svx_j_3+i);svy_j_i_3=*(svy_j_3+i);svz_j_i_3=*(svz_j_3+i);
                            svx_j_i_4=*(svx_j_4+i);svy_j_i_4=*(svy_j_4+i);svz_j_i_4=*(svz_j_4+i);
                            
                            for (k=nz1;k<=nz2;k++){
                                sxx_x = dx*(b1*(sxx[j][i+1][k]-sxx[j][i][k])+
                                            b2*(sxx[j][i+2][k]-sxx[j][i-1][k])+
                                            b3*(sxx[j][i+3][k]-sxx[j][i-2][k])+
                                            b4*(sxx[j][i+4][k]-sxx[j][i-3][k])+
                                            b5*(sxx[j][i+5][k]-sxx[j][i-4][k])+
                                            b6*(sxx[j][i+6][k]-sxx[j][i-5][k]));
                                
                                sxy_y = dy*(b1*(sxy[j][i][k]-sxy[j-1][i][k])+
                                            b2*(sxy[j+1][i][k]-sxy[j-2][i][k])+
                                            b3*(sxy[j+2][i][k]-sxy[j-3][i][k])+
                                            b4*(sxy[j+3][i][k]-sxy[j-4][i][k])+
                                            b5*(sxy[j+4][i][k]-sxy[j-5][i][k])+
                                            b6*(sxy[j+5][i][k]-sxy[j-6][i][k]));
                                
                                sxz_z = dz*(b1*(sxy[j][i][k]-sxy[j][i][k-1])+
                                            b2*(sxy[j][i][k+1]-sxy[j][i][k-2])+
                                            b3*(sxy[j][i][k+2]-sxy[j][i][k-3])+
                                            b4*(sxy[j][i][k+3]-sxy[j][i][k-4])+
                                            b5*(sxy[j][i][k+4]-sxy[j][i][k-5])+
                                            b6*(sxy[j][i][k+5]-sxy[j][i][k-6]));
                                
                                syy_y = dy*(b1*(syy[j+1][i][k]-syy[j][i][k])+
                                            b2*(syy[j+2][i][k]-syy[j-1][i][k])+
                                            b3*(syy[j+3][i][k]-syy[j-2][i][k])+
                                            b4*(syy[j+4][i][k]-syy[j-3][i][k])+
                                            b5*(syy[j+5][i][k]-syy[j-4][i][k])+
                                            b6*(syy[j+6][i][k]-syy[j-5][i][k]));
                                
                                sxy_x = dx*(b1*(sxy[j][i][k]-sxy[j][i-1][k])+
                                            b2*(sxy[j][i+1][k]-sxy[j][i-2][k])+
                                            b3*(sxy[j][i+2][k]-sxy[j][i-3][k])+
                                            b4*(sxy[j][i+3][k]-sxy[j][i-4][k])+
                                            b5*(sxy[j][i+4][k]-sxy[j][i-5][k])+
                                            b6*(sxy[j][i+5][k]-sxy[j][i-6][k]));
                                
                                syz_z = dz*(b1*(syz[j][i][k]-syz[j][i][k-1])+
                                            b2*(syz[j][i][k+1]-syz[j][i][k-2])+
                                            b3*(syz[j][i][k+2]-syz[j][i][k-3])+
                                            b4*(syz[j][i][k+3]-syz[j][i][k-4])+
                                            b5*(syz[j][i][k+4]-syz[j][i][k-5])+
                                            b6*(syz[j][i][k+5]-syz[j][i][k-6]));
                                
                                szz_z = dz*(b1*(szz[j][i][k+1]-szz[j][i][k])+
                                            b2*(szz[j][i][k+2]-szz[j][i][k-1])+
                                            b3*(szz[j][i][k+3]-szz[j][i][k-2])+
                                            b4*(szz[j][i][k+4]-szz[j][i][k-3])+
                                            b5*(szz[j][i][k+5]-szz[j][i][k-4])+
                                            b6*(szz[j][i][k+6]-szz[j][i][k-5]));
                                
                                sxz_x = dx*(b1*(sxz[j][i][k]-sxz[j][i-1][k])+
                                            b2*(sxz[j][i+1][k]-sxz[j][i-2][k])+
                                            b3*(sxz[j][i+2][k]-sxz[j][i-3][k])+
                                            b4*(sxz[j][i+3][k]-sxz[j][i-4][k])+
                                            b5*(sxz[j][i+4][k]-sxz[j][i-5][k])+
                                            b6*(sxz[j][i+5][k]-sxz[j][i-6][k]));
                                
                                
                                syz_y = dy*(b1*(syz[j][i][k]-syz[j-1][i][k])+
                                            b2*(syz[j+1][i][k]-syz[j-2][i][k])+
                                            b3*(syz[j+2][i][k]-syz[j-3][i][k])+
                                            b4*(syz[j+3][i][k]-syz[j-4][i][k])+
                                            b5*(syz[j+4][i][k]-syz[j-5][i][k])+
                                            b6*(syz[j+5][i][k]-syz[j-6][i][k]));
                                
                                
                                
                                *(svx_j_i+k)=sxx_x + sxy_y +sxz_z;
                                *(svy_j_i+k)=syy_y + sxy_x + syz_z;
                                *(svz_j_i+k)=szz_z + sxz_x + syz_y;
                                
                                
                                /* updating components of particle velocities */
                                vx[j][i][k]+= ((c1*(*(svx_j_i+k))+c2*(*(svx_j_i_2+k))+c3*(*(svx_j_i_3+k))+c4*(*(svx_j_i_4+k)))/rip[j][i][k]);
                                vy[j][i][k]+= ((c1*(*(svy_j_i+k))+c2*(*(svy_j_i_2+k))+c3*(*(svy_j_i_3+k))+c4*(*(svy_j_i_4+k)))/rjp[j][i][k]);
                                vz[j][i][k]+= ((c1*(*(svz_j_i+k))+c2*(*(svz_j_i_2+k))+c3*(*(svz_j_i_3+k))+c4*(*(svz_j_i_4+k)))/rkp[j][i][k]);
                                
                            }
                        }
                    }
                    
                    break;
                    
            }
            break; /* break for FDORDER_TIME=4 */
            
    }
    
    /* Adding body force components to corresponding particle velocities */
    for (l=1;l<=nsrc;l++) {
        i=(int)srcpos_loc[1][l];
        j=(int)srcpos_loc[2][l];
        k=(int)srcpos_loc[3][l];
        //amp=signals[l][nt]; // unscaled force amplitude
        amp=(DT*signals[l][nt])/(DX*DY*DZ);// scaled force amplitude with F= 1N
        
        //fprintf(FP," amp at timestep nt %i = %5.5e with DX=%5.2f DY=%5.2f DZ=%5.2f  DT=%5.8f\n",nt,amp,DX,DY,DZ,DT);
        
        switch (stype[l]){
            case 2 :
                vx[j][i][k] += amp/rip[j][i][k];  /* single force in x */
                break;
            case 3 :
                vy[j][i][k] += amp/rjp[j][i][k];  /* single force in y / vertical direction */
                break;
            case 4 :
                vz[j][i][k] += amp/rkp[j][i][k];  /* single force in z */
                break;
            case 5 :
                alpha_rad=SOURCE_ALPHA*PI/180; /* custom force */
                beta_rad=SOURCE_BETA*PI/180;
                vx[j][i][k]+=cos(alpha_rad)*sin(beta_rad)*amp/rip[j][i][k];
                vy[j][i][k]+=cos(beta_rad)*amp/rjp[j][i][k]; /*vertical component*/
                vz[j][i][k]+=sin(alpha_rad)*sin(beta_rad)*amp/rkp[j][i][k];
                break;
            default:
                break;
        }
    }
    
    
    /* absorbing boundary condition (exponential damping) */
    
    if (ABS_TYPE==2){
        for (j=ny1;j<=ny2;j++){
            for (i=nx1;i<=nx2;i++){
                for (k=nz1;k<=nz2;k++){
                    vx[j][i][k]*=absorb_coeff[j][i][k];
                    vy[j][i][k]*=absorb_coeff[j][i][k];
                    vz[j][i][k]*=absorb_coeff[j][i][k];
                    sxy[j][i][k]*=absorb_coeff[j][i][k];
                    syz[j][i][k]*=absorb_coeff[j][i][k];
                    sxz[j][i][k]*=absorb_coeff[j][i][k];
                    sxx[j][i][k]*=absorb_coeff[j][i][k];
                    syy[j][i][k]*=absorb_coeff[j][i][k];
                    szz[j][i][k]*=absorb_coeff[j][i][k];
                    
                }
            }
        }
    }
    
    
    
    if (LOG)
        if ((MYID==0) && ((nt+(OUTNTIMESTEPINFO-1))%OUTNTIMESTEPINFO)==0) {
            time2=MPI_Wtime();
            time=time2-time1;
            fprintf(FP," Real time for particle velocity update: \t %4.2f s.\n",time);
        }
    return time;
    
}