PCG_SH.c 19.8 KB
Newer Older
1
/*-----------------------------------------------------------------------------------------
2
 * Copyright (C) 2016  For the list of authors, see file AUTHORS.
3
 *
Florian Wittkamp's avatar
Florian Wittkamp committed
4
 * This file is part of IFOS.
5
 *
Florian Wittkamp's avatar
Florian Wittkamp committed
6
 * IFOS is free software: you can redistribute it and/or modify
7 8 9
 * it under the terms of the GNU General Public License as published by
 * the Free Software Foundation, version 2.0 of the License only.
 *
Florian Wittkamp's avatar
Florian Wittkamp committed
10
 * IFOS is distributed in the hope that it will be useful,
11 12 13 14 15
 * but WITHOUT ANY WARRANTY; without even the implied warranty of
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 * GNU General Public License for more details.
 *
 * You should have received a copy of the GNU General Public License
Florian Wittkamp's avatar
Florian Wittkamp committed
16
 * along with IFOS. See file COPYING and/or <http://www.gnu.org/licenses/gpl-2.0.html>.
17 18 19 20 21 22 23 24 25 26
 -----------------------------------------------------------------------------------------*/

/*------------------------------------------------------------------------
 * Module for the Preconditioned Conjugate Gradient Method (PCG)
 * for the material parameters vp, vs, rho and lambda, mu, rho respectively
 *
 * ----------------------------------------------------------------------*/

#include "fd.h"

27
void PCG_SH(float ** taper_coeff, int nsrc, float ** srcpos, int ** recpos, int ntr_glob, int iter, int nfstart_jac, float ** waveconv_u, float C_vs, float ** gradp_u, float ** waveconv_rho, float C_rho, float ** gradp_rho, float Vs_avg, float F_LOW_PASS, int PCG_iter_start){
28 29
    
    extern int NX, NY, IDX, IDY, SPATFILTER, GRAD_FILTER;
30
    extern int FORWARD_ONLY, SWS_TAPER_GRAD_VERT, SWS_TAPER_GRAD_HOR, SWS_TAPER_GRAD_SOURCES, SWS_TAPER_FILE;
31 32 33 34 35 36 37
    extern int POS[3], MYID, ACOUSTIC,WAVETYPE;
    extern char JACOBIAN[STRING_SIZE];
    
    char jac[225], jac2[225];
    int i, j;
    float betaz, betan, gradplastiter, gradclastiter, betar, beta;
    extern FILE *FP;
38
    FILE *FP3, *FP4, *FP6, *FP5 = NULL;
39 40
    int use_conjugate_1=1;
    int use_conjugate_2=1;
41

42
    
43 44 45 46 47 48 49 50
    /* Check if conjugate gradient can be used */
    if( (iter-PCG_iter_start) < 2 ) {
        use_conjugate_2=0;
        if( (iter-PCG_iter_start) < 1 ) {
            use_conjugate_1=0;
        }
    }
    
51 52 53 54
    /* ===================================================================================================================================================== */
    /* ===================================================== GRADIENT Zs ================================================================================== */
    /* ===================================================================================================================================================== */
    
55
    if((FORWARD_ONLY==0)&&(!ACOUSTIC)){
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
        
        /* Preconditioning of the gradient */
        /* ------------------------------- */
        
        /* apply taper on the gradient */
        /* --------------------------- */
        if (SWS_TAPER_GRAD_VERT){   /*vertical gradient taper is applied*/
            taper_grad(waveconv_u,taper_coeff,srcpos,nsrc,recpos,ntr_glob,1);}
        
        if (SWS_TAPER_GRAD_HOR){   /*horizontal gradient taper is applied*/
            taper_grad(waveconv_u,taper_coeff,srcpos,nsrc,recpos,ntr_glob,2);}
        
        if (SWS_TAPER_GRAD_SOURCES){   /*cylindrical taper around sources is applied*/
            taper_grad(waveconv_u,taper_coeff,srcpos,nsrc,recpos,ntr_glob,3);}
        
        if (SWS_TAPER_FILE){   /* read taper from BIN-File*/
            taper_grad(waveconv_u,taper_coeff,srcpos,nsrc,recpos,ntr_glob,5);}
        
        /* save gradient */
        sprintf(jac,"%s_g_u_SH.old.%i.%i",JACOBIAN,POS[1],POS[2]);
        FP3=fopen(jac,"wb");
        
        for (i=1;i<=NX;i=i+IDX){
            for (j=1;j<=NY;j=j+IDY){
                fwrite(&waveconv_u[j][i],sizeof(float),1,FP3);
            }
        }
        
        fclose(FP3);
        
        MPI_Barrier(MPI_COMM_WORLD);
        
        /* merge gradient file */
        sprintf(jac,"%s_g_u_SH.old",JACOBIAN);
        if (MYID==0) mergemod(jac,3);
        
        /* Normalize gradient to maximum value */
        /*norm(waveconv_u,iter,2);*/
        
        /* apply spatial wavelength filter */
        if(SPATFILTER==1){
            if (MYID==0){
                fprintf(FP,"\n Spatial filter is applied to gradient (written by PE %d)\n",MYID);}
            spat_filt(waveconv_u,iter,2);}
        
        /* apply 2D-Gaussian filter*/
102
        if(GRAD_FILTER==1){smooth(waveconv_u,2,1,Vs_avg,F_LOW_PASS);}
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
        
        /* output of the preconditioned gradient */
        for (i=1;i<=NX;i=i+IDX){
            for (j=1;j<=NY;j=j+IDY){
                waveconv_u[j][i] = C_vs * waveconv_u[j][i];
                gradp_u[j][i]=waveconv_u[j][i];
            }
        }
        
        /* save gradient for output as inversion result */
        if(iter==nfstart_jac){
            sprintf(jac,"%s_p_u_SH_it%d.old.%i.%i",JACOBIAN,iter,POS[1],POS[2]);
            FP3=fopen(jac,"wb");
            
            for (i=1;i<=NX;i=i+IDX){
                for (j=1;j<=NY;j=j+IDY){
                    fwrite(&waveconv_u[j][i],sizeof(float),1,FP3);
                }
            }
            
            fclose(FP3);
            
            MPI_Barrier(MPI_COMM_WORLD);
            
            /* merge gradient file */
            sprintf(jac,"%s_p_u_SH_it%d.old",JACOBIAN,iter);
            if (MYID==0) mergemod(jac,3);
            MPI_Barrier(MPI_COMM_WORLD);
            sprintf(jac,"%s_p_u_SH_it%d.old.%i.%i",JACOBIAN,iter,POS[1],POS[2]);
            remove(jac);
        }
        
        /* calculate conjugate gradient direction, if iter > 1 (after Mora 1987) */
        /* --------------------------------------------------------------------- */
        
        if((iter>1)&&(use_conjugate_1)){
            
            sprintf(jac,"%s_p_u_SH.old.%i.%i",JACOBIAN,POS[1],POS[2]);
            FP6=fopen(jac,"rb");
            
            if((iter>2)&&(use_conjugate_2)){
                sprintf(jac2,"%s_c_u_SH.old.%i.%i",JACOBIAN,POS[1],POS[2]);
                FP5=fopen(jac2,"rb");
            }
            
            /* apply scalar product to obtain the coefficient beta */
            betaz = 0.0;
            betan = 0.0;
            for (i=1;i<=NX;i=i+IDX){
                for (j=1;j<=NY;j=j+IDY){
                    
                    fread(&gradplastiter,sizeof(float),1,FP6);

156
                    /*if(gradglastiter==gradg[j][i]) declare_error("TEST1");*/
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
                    /*if (MYID==10)  printf("TEST beta (MYID=%d) bei (j,i)=(%i,%i): gradg(k-1) = %e, gradg(k) = %e\n",MYID,j,i,gradglastiter,gradg[j][i]);*/
                    
                    /*
                     betaz += (1e5*gradp[j][i]) * ( (1e5*gradg[j][i]) - (1e5*gradglastiter) );
                     betan += (1e5*gradplastiter) * (1e5*gradglastiter);
                     */
                    
                    /* Polak and Ribiere */
                    /*betaz += (gradp_u[j][i]) * ( (gradg_u[j][i]) - (gradglastiter) );
                     betan += (gradplastiter) * (gradglastiter);*/
                    
                    /* Polak and Ribiere */
                    betaz += (gradp_u[j][i]) * ( (gradp_u[j][i]) - (gradplastiter) );
                    betan += (gradplastiter) * (gradplastiter);
                    
                    /* Fletcher and Reeves */
                    /*betaz += (gradp[j][i]) * (gradg[j][i]);
                     betan += (gradplastiter) * (gradglastiter);*/
                    
                }
            }
            
           // printf("TEST: vor exchange (MYID=%d): beta = betaz/betan = %e/%e = %e\n",MYID,betaz,betan,betaz/betan);
            
            /*betaz = exchange_L2(betaz,1,1);
             betan = exchange_L2(betan,1,1);*/
            
            betar = 0.0;
            MPI_Allreduce(&betaz,&betar,1,MPI_FLOAT,MPI_SUM,MPI_COMM_WORLD);
            betaz = betar;
            
            betar = 0.0;
            MPI_Allreduce(&betan,&betar,1,MPI_FLOAT,MPI_SUM,MPI_COMM_WORLD);
            betan = betar;
            
            // 			beta = betaz/betan;
            beta = 0.0f;
            if(betan !=0.0f) beta = betaz/betan;
            
            /* direction reset */
            if(beta<0.0){beta = 0.0;}
            
            /*betaVs = beta;*/
           // printf("\n\nTEST: nach exchange (MYID=%d): beta = %e / %e = %e\n",MYID,betaz,betan,beta);
            
            fseek(FP6,0,SEEK_SET);
            
            for (i=1;i<=NX;i=i+IDX){
                for (j=1;j<=NY;j=j+IDY){
                    
                    if(iter==2){
                        fread(&gradplastiter,sizeof(float),1,FP6);
                        waveconv_u[j][i] = gradp_u[j][i] + gradplastiter * beta;
                    }
                    
                    if((iter>2)&&(use_conjugate_2)){
                        fread(&gradclastiter,sizeof(float),1,FP5);
                        waveconv_u[j][i] = gradp_u[j][i] + gradclastiter * beta;
                    }
                    
                
//                    if (iter >= 2){
//                        if (isnan(waveconv_u[j][i]) || isinf(waveconv_u[j][i])){
//                            sum = 0.0;
//                            h = 0;
//                            for (ii=-1;ii<=1;ii++){
//                                for (jj=-1;jj<=1;jj++){
//                                    if (isnan(waveconv_u[j+jj][i+ii]) || isinf(waveconv_u[j+jj][i+ii])) continue;
//                                    sum += waveconv_u[j+jj][i+ii];
//                                    h++;
//                                }
//                            }
//                            if (h>0) waveconv_u[j][i] = sum / h;
//                            else waveconv_u[j][i] = 0.0;
//                        }
//                    }
                    
                }
            }
            
            fclose(FP6);
            
            if((iter>2)&&(use_conjugate_2)){fclose(FP5);}
        }
        
        /* output of the conjugate gradient */
        if((iter>1)&&(use_conjugate_1)){
            sprintf(jac2,"%s_c_u_SH.old.%i.%i",JACOBIAN,POS[1],POS[2]);
            FP5=fopen(jac2,"wb");
            
            for (i=1;i<=NX;i=i+IDX){
                for (j=1;j<=NY;j=j+IDY){
                    fwrite(&waveconv_u[j][i],sizeof(float),1,FP5);
                }
            }
            
            fclose(FP5);
            
            MPI_Barrier(MPI_COMM_WORLD);
            
            /* merge gradient file */
            sprintf(jac2,"%s_c_u_SH.old",JACOBIAN);
            if (MYID==0) mergemod(jac2,3);
        }
        
        sprintf(jac,"%s_p_u_SH.old.%i.%i",JACOBIAN,POS[1],POS[2]);
        FP4=fopen(jac,"wb");
        
        /* output of the preconditioned gradient */
        for (i=1;i<=NX;i=i+IDX){
            for (j=1;j<=NY;j=j+IDY){
                /*fwrite(&waveconv_u[j][i],sizeof(float),1,FP4);*/
                fwrite(&gradp_u[j][i],sizeof(float),1,FP4);
            }
        }
        
        fclose(FP4);
        
        MPI_Barrier(MPI_COMM_WORLD);
        
        /* merge gradient file */
        sprintf(jac,"%s_p_u_SH.old",JACOBIAN);
        if (MYID==0) mergemod(jac,3);
        
    }
    
    /* ===================================================================================================================================================== */
    /* ===================================================== GRADIENT rho ================================================================================== */
    /* ===================================================================================================================================================== */
    
287
    if(FORWARD_ONLY==0){
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
        
        /* Preconditioning of the gradient */
        /* ------------------------------- */
        if (SWS_TAPER_GRAD_VERT){   /*vertical gradient taper is applied*/
            taper_grad(waveconv_rho,taper_coeff,srcpos,nsrc,recpos,ntr_glob,1);}
        
        if (SWS_TAPER_GRAD_HOR){   /*horizontal gradient taper is applied*/
            taper_grad(waveconv_rho,taper_coeff,srcpos,nsrc,recpos,ntr_glob,2);}
        
        if (SWS_TAPER_GRAD_SOURCES){   /*cylindrical taper around sources is applied*/
            taper_grad(waveconv_rho,taper_coeff,srcpos,nsrc,recpos,ntr_glob,3);}
        
        if (SWS_TAPER_FILE){   /* read taper from BIN-File*/
            taper_grad(waveconv_rho,taper_coeff,srcpos,nsrc,recpos,ntr_glob,6);}
        
        /* save gradient */
        sprintf(jac,"%s_g_rho_SH.old.%i.%i",JACOBIAN,POS[1],POS[2]);
        FP3=fopen(jac,"wb");
        
        for (i=1;i<=NX;i=i+IDX){
            for (j=1;j<=NY;j=j+IDY){
                fwrite(&waveconv_rho[j][i],sizeof(float),1,FP3);
            }
        }
        
        fclose(FP3);
        
        MPI_Barrier(MPI_COMM_WORLD);
        
        /* merge gradient file */
        sprintf(jac,"%s_g_rho_SH.old",JACOBIAN);
        if (MYID==0) mergemod(jac,3);
        
        /* Normalize gradient to maximum value */
        /*norm(waveconv_rho,iter,3);*/
        
        /* apply spatial wavelength filter */
        if(SPATFILTER==1){
            if (MYID==0){
                fprintf(FP,"\n Spatial filter is applied to gradient (written by PE %d)\n",MYID);}
            spat_filt(waveconv_rho,iter,3);}
        
        /* apply 2D-Gaussian filter*/
331
        if(GRAD_FILTER==1){smooth(waveconv_rho,3,1,Vs_avg,F_LOW_PASS);}
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
        
        /* output of the preconditioned gradient */
        for (i=1;i<=NX;i=i+IDX){
            for (j=1;j<=NY;j=j+IDY){
                waveconv_rho[j][i] = C_rho * waveconv_rho[j][i];
                gradp_rho[j][i]=waveconv_rho[j][i];
            }
        }
        
        /* save gradient for output as inversion result */
        if(iter==nfstart_jac){
            sprintf(jac,"%s_p_rho_SH_it%d.old.%i.%i",JACOBIAN,iter,POS[1],POS[2]);
            FP3=fopen(jac,"wb");
            
            for (i=1;i<=NX;i=i+IDX){
                for (j=1;j<=NY;j=j+IDY){
                    fwrite(&waveconv_rho[j][i],sizeof(float),1,FP3);
                }
            }
            
            fclose(FP3);
            
            MPI_Barrier(MPI_COMM_WORLD);
            
            /* merge gradient file */ 
            sprintf(jac,"%s_p_rho_SH_it%d.old",JACOBIAN,iter);
            if (MYID==0) mergemod(jac,3);
            MPI_Barrier(MPI_COMM_WORLD);
            sprintf(jac,"%s_p_rho_SH_it%d.old.%i.%i",JACOBIAN,iter,POS[1],POS[2]);
            remove(jac);
        }
        
        /* calculate conjugate gradient direction, if iter > 1 (after Mora 1987) */
        /* --------------------------------------------------------------------- */
        
        if((iter>1)&&(use_conjugate_1)){
            
            sprintf(jac,"%s_p_rho_SH.old.%i.%i",JACOBIAN,POS[1],POS[2]);
            FP6=fopen(jac,"rb");
            
            if((iter>2)&&(use_conjugate_2)){
                sprintf(jac2,"%s_c_rho_SH.old.%i.%i",JACOBIAN,POS[1],POS[2]);
                FP5=fopen(jac2,"rb");
            }
            
            /* apply scalar product to obtain the coefficient beta */
            betaz = 0.0;
            betan = 0.0;
            for (i=1;i<=NX;i=i+IDX){
                for (j=1;j<=NY;j=j+IDY){
                    
                    fread(&gradplastiter,sizeof(float),1,FP6);
                    
385
                    /*if(gradglastiter==gradg[j][i]) declare_error("TEST1");*/
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
                    /*if (MYID==10)  printf("TEST beta (MYID=%d) bei (j,i)=(%i,%i): gradg(k-1) = %e, gradg(k) = %e\n",MYID,j,i,gradglastiter,gradg[j][i]);*/
                    
                    /*
                     betaz += (1e5*gradp[j][i]) * ( (1e5*gradg[j][i]) - (1e5*gradglastiter) );
                     betan += (1e5*gradplastiter) * (1e5*gradglastiter);
                     */
                    
                    /* Polak and Ribiere */
                    /*betaz += (gradp_rho[j][i]) * ( (gradg_rho[j][i]) - (gradglastiter) );
                     betan += (gradplastiter) * (gradglastiter);*/
                    
                    /* Polak and Ribiere */
                    betaz += (gradp_rho[j][i]) * ( (gradp_rho[j][i]) - (gradplastiter) );
                    betan += (gradplastiter) * (gradplastiter);
                    
                    /* Fletcher and Reeves */
                    /*betaz += (gradp[j][i]) * (gradg[j][i]);
                     betan += (gradplastiter) * (gradglastiter);*/
                    
                }
            }
            
            /*printf("TEST: vor exchange (MYID=%d): beta = betaz/betan = %e/%e = %e\n",MYID,betaz,betan,betaz/betan);*/
            
            /*betaz = exchange_L2(betaz,1,1);
             betan = exchange_L2(betan,1,1);*/
            
            betar = 0.0;
            MPI_Allreduce(&betaz,&betar,1,MPI_FLOAT,MPI_SUM,MPI_COMM_WORLD);
            betaz = betar;
            
            betar = 0.0;
            MPI_Allreduce(&betan,&betar,1,MPI_FLOAT,MPI_SUM,MPI_COMM_WORLD);
            betan = betar;
            
            // 			beta = betaz/betan;
            beta = 0.0f;
            if(betan !=0.0f) beta = betaz/betan;
            
            /* direction reset */
            if(beta<0.0){beta = 0.0;}
            
            /*betarho = beta;*/
            //printf("\n\nTEST: nach exchange (MYID=%d): beta = %e / %e = %e\n",MYID,betaz,betan,beta);
            
            fseek(FP6,0,SEEK_SET);
            
            for (i=1;i<=NX;i=i+IDX){
                for (j=1;j<=NY;j=j+IDY){
                    
                    if(iter==2){
                        fread(&gradplastiter,sizeof(float),1,FP6);
                        waveconv_rho[j][i] = gradp_rho[j][i] + gradplastiter * beta;
                    }
                    
                    if((iter>2)&&(use_conjugate_2)){
                        fread(&gradclastiter,sizeof(float),1,FP5);
                        waveconv_rho[j][i] = gradp_rho[j][i] + gradclastiter * beta;
                    }
                    
                    /*if (iter >= 2){
                     if (isnan(waveconv_u[j][i]) || isinf(waveconv_u[j][i])){
                     sum = 0.0;
                     h = 0;
                     for (ii=-1;ii<=1;ii++){
                     for (jj=-1;jj<=1;jj++){
                     if (isnan(waveconv_rho[j+jj][i+ii]) || isinf(waveconv_rho[j+jj][i+ii])) continue;
                     sum += waveconv_rho[j+jj][i+ii];
                     h++;
                     }
                     }
                     if (h>0) waveconv_rho[j][i] = sum / h;
                     else waveconv_rho[j][i] = 0.0;
                     }
                     }*/
                }
            }
            
            fclose(FP6);
            
            if((iter>2)&&(use_conjugate_2)){fclose(FP5);}
        }
        
        /* output of the conjugate gradient */
        if((iter>1)&&(use_conjugate_1)){
            sprintf(jac2,"%s_c_rho_SH.old.%i.%i",JACOBIAN,POS[1],POS[2]);
            FP5=fopen(jac2,"wb");
            
            for (i=1;i<=NX;i=i+IDX){
                for (j=1;j<=NY;j=j+IDY){
                    fwrite(&waveconv_rho[j][i],sizeof(float),1,FP5);
                }
            }
            
            fclose(FP5);
            
            MPI_Barrier(MPI_COMM_WORLD);
            
            /* merge gradient file */ 
            sprintf(jac2,"%s_c_rho_SH.old",JACOBIAN);
            if (MYID==0) mergemod(jac2,3);  
        }
        
        sprintf(jac,"%s_p_rho_SH.old.%i.%i",JACOBIAN,POS[1],POS[2]);
        FP4=fopen(jac,"wb");
        
        /* output of the preconditioned gradient */
        for (i=1;i<=NX;i=i+IDX){
            for (j=1;j<=NY;j=j+IDY){
                /*fwrite(&waveconv_rho[j][i],sizeof(float),1,FP4);*/
                fwrite(&gradp_rho[j][i],sizeof(float),1,FP4);
            }
        }
        
        fclose(FP4);
        
        MPI_Barrier(MPI_COMM_WORLD);
        
        /* merge gradient file */ 
        sprintf(jac,"%s_p_rho_SH.old",JACOBIAN);
        if (MYID==0) mergemod(jac,3);
        
    }
}