1a312c225SMatthew G Knepley /* Defines the basic SNES object */ 2*19653cdaSPeter Brune #include <../src/snes/impls/ngmres/snesngmres.h> 3*19653cdaSPeter Brune #include <petscblaslapack.h> 4a312c225SMatthew G Knepley 5087dfb9eSxuemin 6*19653cdaSPeter Brune /*MC 7*19653cdaSPeter Brune SNESNGMRES - The Nonlinear Generalized Minimum Residual (NGMRES) method of Oosterlee and Washio. 8087dfb9eSxuemin 9*19653cdaSPeter Brune Level: beginner 10a312c225SMatthew G Knepley 11*19653cdaSPeter Brune Notes: Supports only left preconditioning 12*19653cdaSPeter Brune 13*19653cdaSPeter Brune "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", C. W. Oosterlee and T. Washio, 14*19653cdaSPeter Brune SIAM Journal on Scientific Computing, 21(5), 2000. 15*19653cdaSPeter Brune 16*19653cdaSPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 17*19653cdaSPeter Brune M*/ 18087dfb9eSxuemin 19a312c225SMatthew G Knepley #undef __FUNCT__ 20a312c225SMatthew G Knepley #define __FUNCT__ "SNESReset_NGMRES" 21a312c225SMatthew G Knepley PetscErrorCode SNESReset_NGMRES(SNES snes) 22a312c225SMatthew G Knepley { 23a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 24a312c225SMatthew G Knepley PetscErrorCode ierr; 25a312c225SMatthew G Knepley 26a312c225SMatthew G Knepley PetscFunctionBegin; 27*19653cdaSPeter Brune ierr = VecDestroyVecs(ngmres->msize, &ngmres->rdot);CHKERRQ(ierr); 28*19653cdaSPeter Brune ierr = VecDestroyVecs(ngmres->msize, &ngmres->xdot);CHKERRQ(ierr); 29a312c225SMatthew G Knepley PetscFunctionReturn(0); 30a312c225SMatthew G Knepley } 31a312c225SMatthew G Knepley 32a312c225SMatthew G Knepley #undef __FUNCT__ 33a312c225SMatthew G Knepley #define __FUNCT__ "SNESDestroy_NGMRES" 34a312c225SMatthew G Knepley PetscErrorCode SNESDestroy_NGMRES(SNES snes) 35a312c225SMatthew G Knepley { 36a312c225SMatthew G Knepley PetscErrorCode ierr; 37a312c225SMatthew G Knepley 38a312c225SMatthew G Knepley PetscFunctionBegin; 39a312c225SMatthew G Knepley ierr = SNESReset_NGMRES(snes);CHKERRQ(ierr); 40a312c225SMatthew G Knepley if (snes->work) {ierr = VecDestroyVecs(snes->nwork, &snes->work);CHKERRQ(ierr);} 41*19653cdaSPeter Brune if (snes->data) { 42*19653cdaSPeter Brune SNES_NGMRES * ngmres = (SNES_NGMRES *)snes->data; 43*19653cdaSPeter Brune ierr = PetscFree4(ngmres->h, ngmres->beta, ngmres->xi, ngmres->q);CHKERRQ(ierr); 44*19653cdaSPeter Brune ierr = PetscFree(ngmres->s);CHKERRQ(ierr); 45*19653cdaSPeter Brune #if PETSC_USE_COMPLEX 46*19653cdaSPeter Brune ierr = PetscFree(ngmres->rwork); 47*19653cdaSPeter Brune #endif 48*19653cdaSPeter Brune ierr = PetscFree(ngmres->work); 49*19653cdaSPeter Brune } 50*19653cdaSPeter Brune ierr = PetscFree(snes->data); 51a312c225SMatthew G Knepley PetscFunctionReturn(0); 52a312c225SMatthew G Knepley } 53a312c225SMatthew G Knepley 54a312c225SMatthew G Knepley #undef __FUNCT__ 55a312c225SMatthew G Knepley #define __FUNCT__ "SNESSetUp_NGMRES" 56a312c225SMatthew G Knepley PetscErrorCode SNESSetUp_NGMRES(SNES snes) 57a312c225SMatthew G Knepley { 58a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 59*19653cdaSPeter Brune PetscInt msize,hsize; 60a312c225SMatthew G Knepley PetscErrorCode ierr; 61a312c225SMatthew G Knepley 62a312c225SMatthew G Knepley PetscFunctionBegin; 63087dfb9eSxuemin msize = ngmres->msize; /* restart size */ 64*19653cdaSPeter Brune hsize = msize * msize; 65087dfb9eSxuemin 66087dfb9eSxuemin 67*19653cdaSPeter Brune //explicit least squares minimization solve 68*19653cdaSPeter Brune ierr = PetscMalloc5(hsize,PetscScalar,&ngmres->h, 69*19653cdaSPeter Brune msize,PetscScalar,&ngmres->beta, 70*19653cdaSPeter Brune msize,PetscScalar,&ngmres->xi, 71*19653cdaSPeter Brune msize,PetscReal,&ngmres->r_norms, 72*19653cdaSPeter Brune hsize,PetscScalar,&ngmres->q);CHKERRQ(ierr); 73*19653cdaSPeter Brune ngmres->nrhs = 1; 74*19653cdaSPeter Brune ngmres->lda = msize; 75*19653cdaSPeter Brune ngmres->ldb = msize; 76*19653cdaSPeter Brune ierr = PetscMalloc(msize*sizeof(PetscScalar),&ngmres->s);CHKERRQ(ierr); 77087dfb9eSxuemin 78*19653cdaSPeter Brune ierr = PetscMemzero(ngmres->h,hsize*sizeof(PetscScalar));CHKERRQ(ierr); 79*19653cdaSPeter Brune ierr = PetscMemzero(ngmres->q,hsize*sizeof(PetscScalar));CHKERRQ(ierr); 80*19653cdaSPeter Brune ierr = PetscMemzero(ngmres->xi,msize*sizeof(PetscScalar));CHKERRQ(ierr); 81*19653cdaSPeter Brune ierr = PetscMemzero(ngmres->beta,msize*sizeof(PetscScalar));CHKERRQ(ierr); 82211b2d39SPeter Brune 83*19653cdaSPeter Brune ngmres->lwork = 12*msize; 84*19653cdaSPeter Brune #if PETSC_USE_COMPLEX 85*19653cdaSPeter Brune ierr = PetscMalloc(sizeof(PetscReal)*ngmres->lwork,&ngmres->rwork); 86*19653cdaSPeter Brune #endif 87*19653cdaSPeter Brune ierr = PetscMalloc(sizeof(PetscScalar)*ngmres->lwork,&ngmres->work); 88211b2d39SPeter Brune 89*19653cdaSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, ngmres->msize, &ngmres->xdot);CHKERRQ(ierr); 90*19653cdaSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, ngmres->msize, &ngmres->rdot);CHKERRQ(ierr); 91211b2d39SPeter Brune ierr = SNESDefaultGetWork(snes, 3);CHKERRQ(ierr); 92a312c225SMatthew G Knepley PetscFunctionReturn(0); 93a312c225SMatthew G Knepley } 94a312c225SMatthew G Knepley 95a312c225SMatthew G Knepley #undef __FUNCT__ 96a312c225SMatthew G Knepley #define __FUNCT__ "SNESSetFromOptions_NGMRES" 97a312c225SMatthew G Knepley PetscErrorCode SNESSetFromOptions_NGMRES(SNES snes) 98a312c225SMatthew G Knepley { 99a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 100a312c225SMatthew G Knepley PetscErrorCode ierr; 101a312c225SMatthew G Knepley 102a312c225SMatthew G Knepley PetscFunctionBegin; 103a312c225SMatthew G Knepley ierr = PetscOptionsHead("SNES NGMRES options");CHKERRQ(ierr); 104*19653cdaSPeter Brune ierr = PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES", ngmres->msize, &ngmres->msize, PETSC_NULL);CHKERRQ(ierr); 105*19653cdaSPeter Brune ierr = PetscOptionsInt("-snes_ngmres_restart", "Maximum iterations before restart.", "SNES", ngmres->k_rmax, &ngmres->k_rmax, PETSC_NULL);CHKERRQ(ierr); 106*19653cdaSPeter Brune ierr = PetscOptionsBool("-snes_ngmres_debug", "Debugging output for NGMRES", "SNES", ngmres->debug, &ngmres->debug, PETSC_NULL);CHKERRQ(ierr); 107a312c225SMatthew G Knepley ierr = PetscOptionsTail();CHKERRQ(ierr); 108a312c225SMatthew G Knepley PetscFunctionReturn(0); 109a312c225SMatthew G Knepley } 110a312c225SMatthew G Knepley 111a312c225SMatthew G Knepley #undef __FUNCT__ 112a312c225SMatthew G Knepley #define __FUNCT__ "SNESView_NGMRES" 113a312c225SMatthew G Knepley PetscErrorCode SNESView_NGMRES(SNES snes, PetscViewer viewer) 114a312c225SMatthew G Knepley { 115a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 116a312c225SMatthew G Knepley PetscBool iascii; 117a312c225SMatthew G Knepley PetscErrorCode ierr; 118a312c225SMatthew G Knepley 119a312c225SMatthew G Knepley PetscFunctionBegin; 120a312c225SMatthew G Knepley ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 121a312c225SMatthew G Knepley if (iascii) { 122a312c225SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer, " Size of space %d\n", ngmres->msize);CHKERRQ(ierr); 123a312c225SMatthew G Knepley } 124a312c225SMatthew G Knepley PetscFunctionReturn(0); 125a312c225SMatthew G Knepley } 126a312c225SMatthew G Knepley 127*19653cdaSPeter Brune 128a312c225SMatthew G Knepley #undef __FUNCT__ 129a312c225SMatthew G Knepley #define __FUNCT__ "SNESSolve_NGMRES" 130211b2d39SPeter Brune 131a312c225SMatthew G Knepley PetscErrorCode SNESSolve_NGMRES(SNES snes) 132a312c225SMatthew G Knepley { 1334a0c5b0cSMatthew G Knepley SNES pc; 134087dfb9eSxuemin SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 135*19653cdaSPeter Brune 136*19653cdaSPeter Brune 137*19653cdaSPeter Brune 138*19653cdaSPeter Brune //present solution, residual, and preconditioned residual 139*19653cdaSPeter Brune Vec x, r, f, b, d; 140*19653cdaSPeter Brune Vec x_A, r_A; 141*19653cdaSPeter Brune 142*19653cdaSPeter Brune //previous iterations to construct the subspace 143*19653cdaSPeter Brune Vec *rdot = ngmres->rdot; 144*19653cdaSPeter Brune Vec *xdot = ngmres->xdot; 145*19653cdaSPeter Brune 146*19653cdaSPeter Brune //coefficients and RHS to the minimization problem 147*19653cdaSPeter Brune PetscScalar *beta = ngmres->beta; 148*19653cdaSPeter Brune PetscScalar *xi = ngmres->xi; 149*19653cdaSPeter Brune PetscReal r_norm, r_A_norm; 150*19653cdaSPeter Brune PetscReal nu; 151*19653cdaSPeter Brune PetscScalar alph_total = 0.; 152*19653cdaSPeter Brune PetscScalar qentry; 153*19653cdaSPeter Brune PetscInt i, j, k, k_restart, l, ivec; 154*19653cdaSPeter Brune 155*19653cdaSPeter Brune //solution selection data 156*19653cdaSPeter Brune PetscBool selectA, selectRestart; 157*19653cdaSPeter Brune PetscReal d_norm, d_min_norm, d_cur_norm; 158*19653cdaSPeter Brune PetscReal r_min_norm; 159*19653cdaSPeter Brune 160a312c225SMatthew G Knepley PetscErrorCode ierr; 161a312c225SMatthew G Knepley 162a312c225SMatthew G Knepley PetscFunctionBegin; 163*19653cdaSPeter Brune //variable initialization 164a312c225SMatthew G Knepley snes->reason = SNES_CONVERGED_ITERATING; 165*19653cdaSPeter Brune x = snes->vec_sol; 166*19653cdaSPeter Brune r = snes->vec_func; 167*19653cdaSPeter Brune f = snes->work[0]; 168*19653cdaSPeter Brune b = snes->vec_rhs; 169*19653cdaSPeter Brune x_A = snes->vec_sol_update; 170*19653cdaSPeter Brune r_A = snes->work[1]; 171*19653cdaSPeter Brune d = snes->work[2]; 172*19653cdaSPeter Brune ierr = VecDuplicate(x,&r);CHKERRQ(ierr); 173a312c225SMatthew G Knepley 1744a0c5b0cSMatthew G Knepley ierr = SNESGetPC(snes, &pc);CHKERRQ(ierr); 175a312c225SMatthew G Knepley ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 176a312c225SMatthew G Knepley snes->iter = 0; 177a312c225SMatthew G Knepley snes->norm = 0.; 178a312c225SMatthew G Knepley ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 179*19653cdaSPeter Brune 180*19653cdaSPeter Brune //initialization -- 181*19653cdaSPeter Brune 182*19653cdaSPeter Brune //r = F(u) 183*19653cdaSPeter Brune ierr = SNESComputeFunction(snes, x, r);CHKERRQ(ierr); /* r = F(x) */ 184a312c225SMatthew G Knepley if (snes->domainerror) { 185a312c225SMatthew G Knepley snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 186a312c225SMatthew G Knepley PetscFunctionReturn(0); 187a312c225SMatthew G Knepley } 188*19653cdaSPeter Brune 189*19653cdaSPeter Brune //nu = (r, r); 190*19653cdaSPeter Brune ierr = VecNorm(r, NORM_2, &r_norm);CHKERRQ(ierr); 191*19653cdaSPeter Brune r_min_norm = r_norm; 192*19653cdaSPeter Brune nu = r_norm*r_norm; 193*19653cdaSPeter Brune if (PetscIsInfOrNanReal(r_norm)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FP, "Infinite or not-a-number generated in norm"); 194*19653cdaSPeter Brune 195*19653cdaSPeter Brune //q_{11} = nu 196*19653cdaSPeter Brune 197*19653cdaSPeter Brune Q(0,0) = nu; 198*19653cdaSPeter Brune ngmres->r_norms[0] = r_norm; 199*19653cdaSPeter Brune //rdot[0] = r 200*19653cdaSPeter Brune ierr = VecCopy(x, xdot[0]);CHKERRQ(ierr); 201*19653cdaSPeter Brune ierr = VecCopy(r, rdot[0]);CHKERRQ(ierr); 202087dfb9eSxuemin 203a312c225SMatthew G Knepley ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 204*19653cdaSPeter Brune snes->norm = r_norm; 205a312c225SMatthew G Knepley ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 206*19653cdaSPeter Brune SNESLogConvHistory(snes, r_norm, 0); 207*19653cdaSPeter Brune ierr = SNESMonitor(snes, 0, r_norm);CHKERRQ(ierr); 208*19653cdaSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,r_norm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 209a312c225SMatthew G Knepley if (snes->reason) PetscFunctionReturn(0); 210a312c225SMatthew G Knepley 211*19653cdaSPeter Brune k_restart = 1; 212*19653cdaSPeter Brune l = 1; 213*19653cdaSPeter Brune for (k=1; k<snes->max_its; k++) { 214087dfb9eSxuemin 215*19653cdaSPeter Brune // 216087dfb9eSxuemin 217*19653cdaSPeter Brune //select which vector of the stored subspace will be updated 218*19653cdaSPeter Brune ivec = k_restart % ngmres->msize; //replace the last used part of the subspace 219*19653cdaSPeter Brune 220*19653cdaSPeter Brune 221*19653cdaSPeter Brune //Computation of x^M 222*19653cdaSPeter Brune ierr = SNESSolve(pc, b, x);CHKERRQ(ierr); 223*19653cdaSPeter Brune //r = F(x) 224*19653cdaSPeter Brune ierr = SNESComputeFunction(snes, x, r);CHKERRQ(ierr); 225*19653cdaSPeter Brune ierr = VecNorm(r, NORM_2, &r_norm);CHKERRQ(ierr); 226*19653cdaSPeter Brune //nu = (r, r) 227*19653cdaSPeter Brune ngmres->r_norms[ivec] = r_norm; 228*19653cdaSPeter Brune nu = r_norm*r_norm; 229*19653cdaSPeter Brune if (r_min_norm > r_norm) r_min_norm = r_norm; //the minimum norm is now of r^M 230*19653cdaSPeter Brune 231*19653cdaSPeter Brune //construct the right hand side and xi factors 232*19653cdaSPeter Brune for (i = 0; i < l; i++) { 233*19653cdaSPeter Brune VecDot(rdot[i], r, &xi[i]); 234*19653cdaSPeter Brune beta[i] = nu - xi[i]; 235*19653cdaSPeter Brune } 236*19653cdaSPeter Brune 237*19653cdaSPeter Brune //construct h 238*19653cdaSPeter Brune for (j = 0; j < l; j++) { 239*19653cdaSPeter Brune for (i = 0; i < l; i++) { 240*19653cdaSPeter Brune H(i, j) = Q(i, j) - xi[i] - xi[j] + nu; 241*19653cdaSPeter Brune } 242*19653cdaSPeter Brune } 243*19653cdaSPeter Brune #ifdef PETSC_MISSING_LAPACK_GELSS 244*19653cdaSPeter Brune SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "NGMRES with LS requires the LAPACK GELSS routine."); 2454a0c5b0cSMatthew G Knepley #else 246*19653cdaSPeter Brune ngmres->m = PetscBLASIntCast(l); 247*19653cdaSPeter Brune ngmres->n = PetscBLASIntCast(l); 248*19653cdaSPeter Brune ngmres->info = PetscBLASIntCast(0); 249*19653cdaSPeter Brune ngmres->rcond = -1.; 250*19653cdaSPeter Brune ngmres->nrhs = 1; 251*19653cdaSPeter Brune #ifdef PETSC_USE_COMPLEX 252*19653cdaSPeter Brune LAPACKgelss_(&ngmres->m, 253*19653cdaSPeter Brune &ngmres->n, 254*19653cdaSPeter Brune &ngmres->nrhs, 255*19653cdaSPeter Brune ngmres->h, 256*19653cdaSPeter Brune &ngmres->lda, 257*19653cdaSPeter Brune ngmres->beta, 258*19653cdaSPeter Brune &ngmres->ldb, 259*19653cdaSPeter Brune ngmres->s, 260*19653cdaSPeter Brune &ngmres->rcond, 261*19653cdaSPeter Brune &ngmres->rank, 262*19653cdaSPeter Brune ngmres->work, 263*19653cdaSPeter Brune &ngmres->lwork, 264*19653cdaSPeter Brune ngmres->rwork, 265*19653cdaSPeter Brune &ngmres->info); 266*19653cdaSPeter Brune #else 267*19653cdaSPeter Brune LAPACKgelss_(&ngmres->m, 268*19653cdaSPeter Brune &ngmres->n, 269*19653cdaSPeter Brune &ngmres->nrhs, 270*19653cdaSPeter Brune ngmres->h, 271*19653cdaSPeter Brune &ngmres->lda, 272*19653cdaSPeter Brune ngmres->beta, 273*19653cdaSPeter Brune &ngmres->ldb, 274*19653cdaSPeter Brune ngmres->s, 275*19653cdaSPeter Brune &ngmres->rcond, 276*19653cdaSPeter Brune &ngmres->rank, 277*19653cdaSPeter Brune ngmres->work, 278*19653cdaSPeter Brune &ngmres->lwork, 279*19653cdaSPeter Brune &ngmres->info); 280*19653cdaSPeter Brune #endif 281*19653cdaSPeter Brune if (ngmres->info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELSS"); 282*19653cdaSPeter Brune if (ngmres->info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SVD failed to converge"); 283a312c225SMatthew G Knepley #endif 284a312c225SMatthew G Knepley 285*19653cdaSPeter Brune alph_total = 0.; 286*19653cdaSPeter Brune for (i = 0; i < l; i++) { 287*19653cdaSPeter Brune alph_total += beta[i]; 288c0bbabecSxuemin } 289*19653cdaSPeter Brune ierr = VecCopy(x, x_A);CHKERRQ(ierr); 290*19653cdaSPeter Brune ierr = VecScale(x_A, 1. - alph_total);CHKERRQ(ierr); 291087dfb9eSxuemin 292*19653cdaSPeter Brune for(i=0;i<l;i++){ 293*19653cdaSPeter Brune ierr= VecAXPY(x_A, beta[i], xdot[i]);CHKERRQ(ierr); 294*19653cdaSPeter Brune } 295*19653cdaSPeter Brune ierr = SNESComputeFunction(snes, x_A, r_A);CHKERRQ(ierr); 296*19653cdaSPeter Brune ierr = VecNorm(r_A, NORM_2, &r_A_norm);CHKERRQ(ierr); 297*19653cdaSPeter Brune 298*19653cdaSPeter Brune selectA = PETSC_TRUE; 299*19653cdaSPeter Brune //Conditions for choosing the accelerated answer -- 300*19653cdaSPeter Brune 301*19653cdaSPeter Brune //Criterion A -- the norm of the function isn't increased above the minimum by too much 302*19653cdaSPeter Brune if (r_A_norm >= ngmres->gammaA*r_min_norm) { 303*19653cdaSPeter Brune selectA = PETSC_FALSE; 304a312c225SMatthew G Knepley } 305087dfb9eSxuemin 306*19653cdaSPeter Brune //Criterion B -- the choice of x^A isn't too close to some other choice 307*19653cdaSPeter Brune ierr=VecCopy(x_A,d);CHKERRQ(ierr); 308*19653cdaSPeter Brune ierr=VecAXPY(d,-1,x);CHKERRQ(ierr); 309*19653cdaSPeter Brune ierr=VecNorm(d,NORM_2,&d_norm);CHKERRQ(ierr); 310*19653cdaSPeter Brune d_min_norm=10000000; 311*19653cdaSPeter Brune for(i=0;i<l;i++) { 312*19653cdaSPeter Brune ierr=VecCopy(x_A,d);CHKERRQ(ierr); 313*19653cdaSPeter Brune ierr=VecAXPY(d,-1,xdot[i]);CHKERRQ(ierr); 314*19653cdaSPeter Brune ierr=VecNorm(d,NORM_2,&d_cur_norm);CHKERRQ(ierr); 315*19653cdaSPeter Brune if(d_cur_norm<d_min_norm) d_min_norm=d_cur_norm; 316*19653cdaSPeter Brune } 317*19653cdaSPeter Brune if (ngmres->epsilonB*d_norm<d_min_norm || sqrt(r_norm)<ngmres->deltaB*sqrt(r_min_norm)) { 318*19653cdaSPeter Brune } else { 319*19653cdaSPeter Brune selectA=PETSC_FALSE; 320*19653cdaSPeter Brune } 321211b2d39SPeter Brune 322087dfb9eSxuemin 323*19653cdaSPeter Brune if (selectA) { 324*19653cdaSPeter Brune if (ngmres->debug) 325*19653cdaSPeter Brune PetscPrintf(PETSC_COMM_WORLD, "picked r_A, ||r_A||_2 = %e, ||r_M||_2 = %e\n", r_A_norm, r_norm); 326*19653cdaSPeter Brune //copy it over 327*19653cdaSPeter Brune r_norm = r_A_norm; 328*19653cdaSPeter Brune nu = r_norm*r_norm; 329*19653cdaSPeter Brune ierr = VecCopy(r_A, r);CHKERRQ(ierr); 330*19653cdaSPeter Brune ierr = VecCopy(x_A, x);CHKERRQ(ierr); 331*19653cdaSPeter Brune } else { 332*19653cdaSPeter Brune if(ngmres->debug) 333*19653cdaSPeter Brune PetscPrintf(PETSC_COMM_WORLD, "picked r_M, ||r_A||_2 = %e, ||r_M||_2 = %e\n", r_A_norm, r_norm); 334*19653cdaSPeter Brune } 335211b2d39SPeter Brune 336*19653cdaSPeter Brune selectRestart = PETSC_FALSE; 337*19653cdaSPeter Brune 338*19653cdaSPeter Brune //maximum iteration criterion 339*19653cdaSPeter Brune if (k_restart > ngmres->k_rmax) { 340*19653cdaSPeter Brune selectRestart = PETSC_TRUE; 341*19653cdaSPeter Brune } 342*19653cdaSPeter Brune 343*19653cdaSPeter Brune //difference stagnation restart 344*19653cdaSPeter Brune if ((ngmres->epsilonB*d_norm > d_min_norm) && 345*19653cdaSPeter Brune (sqrt(r_A_norm) > ngmres->deltaB*sqrt(r_min_norm))) { 346*19653cdaSPeter Brune if (ngmres->debug) 347*19653cdaSPeter Brune PetscPrintf(PETSC_COMM_WORLD, "difference restart: %e > %e\n", ngmres->epsilonB*d_norm, d_min_norm); 348*19653cdaSPeter Brune selectRestart = PETSC_TRUE; 349*19653cdaSPeter Brune } 350*19653cdaSPeter Brune 351*19653cdaSPeter Brune // residual stagnation restart 352*19653cdaSPeter Brune if (sqrt(r_A_norm) > ngmres->gammaC*sqrt(r_min_norm)) { 353*19653cdaSPeter Brune if (ngmres->debug) 354*19653cdaSPeter Brune PetscPrintf(PETSC_COMM_WORLD, "residual restart: %e > %e\n", sqrt(r_A_norm), ngmres->gammaC*sqrt(r_min_norm)); 355*19653cdaSPeter Brune selectRestart = PETSC_TRUE; 356*19653cdaSPeter Brune } 357*19653cdaSPeter Brune 358*19653cdaSPeter Brune if (selectRestart) { 359*19653cdaSPeter Brune if (ngmres->debug) 360*19653cdaSPeter Brune PetscPrintf(PETSC_COMM_WORLD, "Restarted at iteration %d\n", k_restart); 361*19653cdaSPeter Brune k_restart = 1; 362*19653cdaSPeter Brune l = 1; 363*19653cdaSPeter Brune //q_{11} = nu 364*19653cdaSPeter Brune ngmres->r_norms[0] = r_norm; 365*19653cdaSPeter Brune nu = r_norm*r_norm; 366*19653cdaSPeter Brune Q(0,0) = nu; 367*19653cdaSPeter Brune //rdot[0] = r 368*19653cdaSPeter Brune ierr = VecCopy(x, xdot[0]);CHKERRQ(ierr); 369*19653cdaSPeter Brune ierr = VecCopy(r, rdot[0]);CHKERRQ(ierr); 370*19653cdaSPeter Brune } else { 371*19653cdaSPeter Brune //select the current size of the subspace 372*19653cdaSPeter Brune if (l < ngmres->msize) { 373*19653cdaSPeter Brune l++; 374*19653cdaSPeter Brune } 375*19653cdaSPeter Brune k_restart++; 376*19653cdaSPeter Brune //place the current entry in the list of previous entries 377*19653cdaSPeter Brune ierr = VecCopy(r, rdot[ivec]);CHKERRQ(ierr); 378*19653cdaSPeter Brune ierr = VecCopy(x, xdot[ivec]);CHKERRQ(ierr); 379*19653cdaSPeter Brune ngmres->r_norms[ivec] = r_norm; 380*19653cdaSPeter Brune if (r_min_norm > r_norm) r_min_norm = r_norm; //the minimum norm is now of r^A 381*19653cdaSPeter Brune for (i = 0; i < l; i++) { 382*19653cdaSPeter Brune VecDot(r, rdot[i], &qentry); 383*19653cdaSPeter Brune Q(i, ivec) = qentry; 384*19653cdaSPeter Brune Q(ivec, i) = qentry; 385*19653cdaSPeter Brune } 386*19653cdaSPeter Brune } 387*19653cdaSPeter Brune 388*19653cdaSPeter Brune SNESLogConvHistory(snes, r_norm, k); 389*19653cdaSPeter Brune ierr = SNESMonitor(snes, k, r_norm);CHKERRQ(ierr); 390*19653cdaSPeter Brune 391087dfb9eSxuemin snes->iter =k; 392*19653cdaSPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,r_norm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 393087dfb9eSxuemin if (snes->reason) PetscFunctionReturn(0); 394a312c225SMatthew G Knepley } 395a312c225SMatthew G Knepley snes->reason = SNES_DIVERGED_MAX_IT; 396a312c225SMatthew G Knepley PetscFunctionReturn(0); 397a312c225SMatthew G Knepley } 398a312c225SMatthew G Knepley 399a312c225SMatthew G Knepley 400a312c225SMatthew G Knepley 401a312c225SMatthew G Knepley EXTERN_C_BEGIN 402a312c225SMatthew G Knepley #undef __FUNCT__ 403a312c225SMatthew G Knepley #define __FUNCT__ "SNESCreate_NGMRES" 404a312c225SMatthew G Knepley PetscErrorCode SNESCreate_NGMRES(SNES snes) 405a312c225SMatthew G Knepley { 406a312c225SMatthew G Knepley SNES_NGMRES *ngmres; 407a312c225SMatthew G Knepley PetscErrorCode ierr; 408a312c225SMatthew G Knepley 409a312c225SMatthew G Knepley PetscFunctionBegin; 410a312c225SMatthew G Knepley snes->ops->destroy = SNESDestroy_NGMRES; 411a312c225SMatthew G Knepley snes->ops->setup = SNESSetUp_NGMRES; 412a312c225SMatthew G Knepley snes->ops->setfromoptions = SNESSetFromOptions_NGMRES; 413a312c225SMatthew G Knepley snes->ops->view = SNESView_NGMRES; 414a312c225SMatthew G Knepley snes->ops->solve = SNESSolve_NGMRES; 415a312c225SMatthew G Knepley snes->ops->reset = SNESReset_NGMRES; 416a312c225SMatthew G Knepley 417a312c225SMatthew G Knepley ierr = PetscNewLog(snes, SNES_NGMRES, &ngmres);CHKERRQ(ierr); 418a312c225SMatthew G Knepley snes->data = (void*) ngmres; 419*19653cdaSPeter Brune ngmres->msize = 10; 420*19653cdaSPeter Brune ngmres->debug = PETSC_FALSE; 421*19653cdaSPeter Brune 422*19653cdaSPeter Brune ngmres->gammaA = 2.; 423*19653cdaSPeter Brune ngmres->gammaC = 2.; 424*19653cdaSPeter Brune ngmres->deltaB = 0.99; 425*19653cdaSPeter Brune ngmres->epsilonB = 0.01; 426*19653cdaSPeter Brune ngmres->k_rmax = 100; 4274a0c5b0cSMatthew G Knepley 4284a0c5b0cSMatthew G Knepley ierr = SNESGetPC(snes, &snes->pc);CHKERRQ(ierr); 429a312c225SMatthew G Knepley PetscFunctionReturn(0); 430a312c225SMatthew G Knepley } 431a312c225SMatthew G Knepley EXTERN_C_END 432