1 #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/ 2 #include <petscblaslapack.h> 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "SNESNGMRESUpdateSubspace_Private" 6 PetscErrorCode SNESNGMRESUpdateSubspace_Private(SNES snes,PetscInt ivec,PetscInt l,Vec F,PetscReal fnorm,Vec X) 7 { 8 SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 9 Vec *Fdot = ngmres->Fdot; 10 Vec *Xdot = ngmres->Xdot; 11 PetscScalar *xi = ngmres->xi; 12 PetscInt i; 13 PetscReal nu; 14 PetscErrorCode ierr; 15 16 PetscFunctionBegin; 17 if (ivec > l) SETERRQ2(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot update vector %d with space size %d!",ivec,l); 18 ierr = VecCopy(F,Fdot[ivec]);CHKERRQ(ierr); 19 ierr = VecCopy(X,Xdot[ivec]);CHKERRQ(ierr); 20 21 ngmres->fnorms[ivec] = fnorm; 22 if (l > 0) { 23 ierr = VecMDot(F,l,Fdot,xi);CHKERRQ(ierr); 24 for (i = 0; i < l; i++) { 25 Q(i,ivec) = xi[i]; 26 Q(ivec,i) = xi[i]; 27 } 28 } else { 29 nu = fnorm*fnorm; 30 Q(0,0) = nu; 31 } 32 PetscFunctionReturn(0); 33 } 34 35 #undef __FUNCT__ 36 #define __FUNCT__ "SNESNGMRESFormCombinedSolution_Private" 37 PetscErrorCode SNESNGMRESFormCombinedSolution_Private(SNES snes,PetscInt l,Vec XM,Vec FM,PetscReal fMnorm,Vec X,Vec XA,Vec FA) 38 { 39 SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 40 PetscInt i,j; 41 Vec *Fdot = ngmres->Fdot; 42 Vec *Xdot = ngmres->Xdot; 43 PetscScalar *beta = ngmres->beta; 44 PetscScalar *xi = ngmres->xi; 45 PetscScalar alph_total = 0.; 46 PetscErrorCode ierr; 47 PetscReal nu; 48 Vec Y = snes->work[2]; 49 PetscBool changed_y,changed_w; 50 51 PetscFunctionBegin; 52 nu = fMnorm*fMnorm; 53 54 /* construct the right hand side and xi factors */ 55 ierr = VecMDot(FM,l,Fdot,xi);CHKERRQ(ierr); 56 for (i = 0; i < l; i++) beta[i] = nu - xi[i]; 57 58 /* construct h */ 59 for (j = 0; j < l; j++) { 60 for (i = 0; i < l; i++) { 61 H(i,j) = Q(i,j)-xi[i]-xi[j]+nu; 62 } 63 } 64 if (l == 1) { 65 /* simply set alpha[0] = beta[0] / H[0, 0] */ 66 if (H(0,0) != 0.) beta[0] = beta[0]/H(0,0); 67 else beta[0] = 0.; 68 } else { 69 #if defined(PETSC_MISSING_LAPACK_GELSS) 70 SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_SUP,"NGMRES with LS requires the LAPACK GELSS routine."); 71 #else 72 ierr = PetscBLASIntCast(l,&ngmres->m);CHKERRQ(ierr); 73 ierr = PetscBLASIntCast(l,&ngmres->n);CHKERRQ(ierr); 74 ngmres->info = 0; 75 ngmres->rcond = -1.; 76 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 77 #if defined(PETSC_USE_COMPLEX) 78 LAPACKgelss_(&ngmres->m,&ngmres->n,&ngmres->nrhs,ngmres->h,&ngmres->lda,ngmres->beta,&ngmres->ldb,ngmres->s,&ngmres->rcond,&ngmres->rank,ngmres->work,&ngmres->lwork,ngmres->rwork,&ngmres->info); 79 #else 80 LAPACKgelss_(&ngmres->m,&ngmres->n,&ngmres->nrhs,ngmres->h,&ngmres->lda,ngmres->beta,&ngmres->ldb,ngmres->s,&ngmres->rcond,&ngmres->rank,ngmres->work,&ngmres->lwork,&ngmres->info); 81 #endif 82 ierr = PetscFPTrapPop();CHKERRQ(ierr); 83 if (ngmres->info < 0) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_LIB,"Bad argument to GELSS"); 84 if (ngmres->info > 0) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_LIB,"SVD failed to converge"); 85 #endif 86 } 87 for (i=0; i<l; i++) { 88 if (PetscIsInfOrNanScalar(beta[i])) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_LIB,"SVD generated inconsistent output"); 89 } 90 alph_total = 0.; 91 for (i = 0; i < l; i++) alph_total += beta[i]; 92 93 ierr = VecCopy(XM,XA);CHKERRQ(ierr); 94 ierr = VecScale(XA,1.-alph_total);CHKERRQ(ierr); 95 ierr = VecMAXPY(XA,l,beta,Xdot);CHKERRQ(ierr); 96 /* check the validity of the step */ 97 ierr = VecCopy(XA,Y);CHKERRQ(ierr); 98 ierr = VecAXPY(Y,-1.0,X);CHKERRQ(ierr); 99 ierr = SNESLineSearchPostCheck(snes->linesearch,X,Y,XA,&changed_y,&changed_w);CHKERRQ(ierr); 100 ierr = SNESComputeFunction(snes,XA,FA);CHKERRQ(ierr); 101 PetscFunctionReturn(0); 102 } 103 104 #undef __FUNCT__ 105 #define __FUNCT__ "SNESNGMRESCalculateDifferences_Private" 106 PetscErrorCode SNESNGMRESCalculateDifferences_Private(SNES snes,PetscInt l,Vec X,Vec F,Vec XM,Vec FM,Vec XA,Vec FA,Vec D,PetscReal *dnorm,PetscReal *dminnorm,PetscReal *fAnorm) 107 { 108 PetscErrorCode ierr; 109 SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 110 PetscReal dcurnorm; 111 Vec *Xdot = ngmres->Xdot; 112 PetscInt i; 113 114 PetscFunctionBegin; 115 if (ngmres->singlereduction) { 116 *dminnorm = -1.0; 117 if (fAnorm) { 118 ierr = VecNormBegin(FA,NORM_2,fAnorm);CHKERRQ(ierr); 119 } 120 if (dnorm) { 121 ierr = VecCopy(XA,D);CHKERRQ(ierr); 122 ierr = VecAXPY(D,-1.0,XM);CHKERRQ(ierr); 123 ierr = VecNormBegin(D,NORM_2,dnorm);CHKERRQ(ierr); 124 } 125 if (dminnorm) { 126 for (i=0; i<l; i++) { 127 ierr=VecAXPY(Xdot[i],-1.0,XA);CHKERRQ(ierr); 128 } 129 for (i=0; i<l; i++) { 130 ierr = VecNormBegin(Xdot[i],NORM_2,&ngmres->xnorms[i]);CHKERRQ(ierr); 131 } 132 } 133 if (fAnorm) {ierr = VecNormEnd(FA,NORM_2,fAnorm);CHKERRQ(ierr);} 134 if (dnorm) {ierr = VecNormEnd(D,NORM_2,dnorm);CHKERRQ(ierr);} 135 if (dminnorm) { 136 for (i=0; i<l; i++) { 137 ierr = VecNormEnd(Xdot[i],NORM_2,&ngmres->xnorms[i]);CHKERRQ(ierr); 138 } 139 for (i=0; i<l; i++) { 140 dcurnorm = ngmres->xnorms[i]; 141 if ((dcurnorm < *dminnorm) || (*dminnorm < 0.0)) *dminnorm = dcurnorm; 142 ierr=VecAXPY(Xdot[i],1.0,XA);CHKERRQ(ierr); 143 } 144 } 145 } else { 146 if (dnorm) { 147 ierr=VecCopy(XA,D);CHKERRQ(ierr); 148 ierr=VecAXPY(D,-1.0,XM);CHKERRQ(ierr); 149 ierr=VecNormBegin(D,NORM_2,dnorm);CHKERRQ(ierr); 150 } 151 if (fAnorm) { 152 ierr=VecNormBegin(FA,NORM_2,fAnorm);CHKERRQ(ierr); 153 } 154 if (dnorm) { 155 ierr=VecNormEnd(D,NORM_2,dnorm);CHKERRQ(ierr); 156 } 157 if (fAnorm) { 158 ierr=VecNormEnd(FA,NORM_2,fAnorm);CHKERRQ(ierr); 159 } 160 if (dminnorm) { 161 *dminnorm = -1.0; 162 for (i=0; i<l; i++) { 163 ierr=VecCopy(XA,D);CHKERRQ(ierr); 164 ierr=VecAXPY(D,-1.0,Xdot[i]);CHKERRQ(ierr); 165 ierr=VecNorm(D,NORM_2,&dcurnorm);CHKERRQ(ierr); 166 if ((dcurnorm < *dminnorm) || (*dminnorm < 0.0)) *dminnorm = dcurnorm; 167 } 168 } 169 } 170 PetscFunctionReturn(0); 171 } 172 173 #undef __FUNCT__ 174 #define __FUNCT__ "SNESNGMRESSelect_Private" 175 PetscErrorCode SNESNGMRESSelect_Private(SNES snes,PetscInt k_restart,Vec XM,Vec FM,PetscReal fMnorm,Vec XA,Vec FA,PetscReal fAnorm,PetscReal dnorm,PetscReal fminnorm,PetscReal dminnorm,Vec X,Vec F,Vec Y,PetscReal *fnorm) 176 { 177 SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 178 PetscErrorCode ierr; 179 PetscBool lssucceed,selectA; 180 181 PetscFunctionBegin; 182 if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) { 183 /* X = X + \lambda(XA - X) */ 184 if (ngmres->monitor) { 185 ierr = PetscViewerASCIIPrintf(ngmres->monitor,"||F_A||_2 = %e, ||F_M||_2 = %e\n",fAnorm,fMnorm);CHKERRQ(ierr); 186 } 187 ierr = VecCopy(FM,F);CHKERRQ(ierr); 188 ierr = VecCopy(XM,X);CHKERRQ(ierr); 189 ierr = VecCopy(XA,Y);CHKERRQ(ierr); 190 ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr); 191 *fnorm = fMnorm; 192 ierr = SNESLineSearchApply(ngmres->additive_linesearch,X,F,fnorm,Y);CHKERRQ(ierr); 193 ierr = SNESLineSearchGetSuccess(ngmres->additive_linesearch,&lssucceed);CHKERRQ(ierr); 194 if (!lssucceed) { 195 if (++snes->numFailures >= snes->maxFailures) { 196 snes->reason = SNES_DIVERGED_LINE_SEARCH; 197 PetscFunctionReturn(0); 198 } 199 } 200 if (ngmres->monitor) { 201 ierr = PetscViewerASCIIPrintf(ngmres->monitor,"Additive solution: ||F||_2 = %e\n",*fnorm);CHKERRQ(ierr); 202 } 203 } else if (ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) { 204 selectA = PETSC_TRUE; 205 /* Conditions for choosing the accelerated answer */ 206 /* Criterion A -- the norm of the function isn't increased above the minimum by too much */ 207 if (fAnorm >= ngmres->gammaA*fminnorm) selectA = PETSC_FALSE; 208 209 /* Criterion B -- the choice of x^A isn't too close to some other choice */ 210 if (ngmres->epsilonB*dnorm<dminnorm || PetscSqrtReal(*fnorm)<ngmres->deltaB*PetscSqrtReal(fminnorm)) { 211 } else selectA=PETSC_FALSE; 212 213 if (selectA) { 214 if (ngmres->monitor) { 215 ierr = PetscViewerASCIIPrintf(ngmres->monitor,"picked X_A, ||F_A||_2 = %e, ||F_M||_2 = %e\n",fAnorm,fMnorm);CHKERRQ(ierr); 216 } 217 /* copy it over */ 218 *fnorm = fAnorm; 219 ierr = VecCopy(FA,F);CHKERRQ(ierr); 220 ierr = VecCopy(XA,X);CHKERRQ(ierr); 221 } else { 222 if (ngmres->monitor) { 223 ierr = PetscViewerASCIIPrintf(ngmres->monitor,"picked X_M, ||F_A||_2 = %e, ||F_M||_2 = %e\n",fAnorm,fMnorm);CHKERRQ(ierr); 224 } 225 *fnorm = fMnorm; 226 ierr = VecCopy(XM,Y);CHKERRQ(ierr); 227 ierr = VecAXPY(Y,-1.0,X);CHKERRQ(ierr); 228 ierr = VecCopy(FM,F);CHKERRQ(ierr); 229 ierr = VecCopy(XM,X);CHKERRQ(ierr); 230 } 231 } else { /* none */ 232 *fnorm = fAnorm; 233 ierr = VecCopy(FA,F);CHKERRQ(ierr); 234 ierr = VecCopy(XA,X);CHKERRQ(ierr); 235 } 236 PetscFunctionReturn(0); 237 } 238 239 #undef __FUNCT__ 240 #define __FUNCT__ "SNESNGMRESSelectRestart_Private" 241 PetscErrorCode SNESNGMRESSelectRestart_Private(SNES snes,PetscInt l,PetscReal fAnorm,PetscReal dnorm,PetscReal fminnorm,PetscReal dminnorm,PetscBool *selectRestart) 242 { 243 SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 244 PetscErrorCode ierr; 245 246 PetscFunctionBegin; 247 *selectRestart = PETSC_FALSE; 248 /* difference stagnation restart */ 249 if ((ngmres->epsilonB*dnorm > dminnorm) && (PetscSqrtReal(fAnorm) > ngmres->deltaB*PetscSqrtReal(fminnorm)) && l > 0) { 250 if (ngmres->monitor) { 251 ierr = PetscViewerASCIIPrintf(ngmres->monitor,"difference restart: %e > %e\n",ngmres->epsilonB*dnorm,dminnorm);CHKERRQ(ierr); 252 } 253 *selectRestart = PETSC_TRUE; 254 } 255 /* residual stagnation restart */ 256 if (PetscSqrtReal(fAnorm) > ngmres->gammaC*PetscSqrtReal(fminnorm)) { 257 if (ngmres->monitor) { 258 ierr = PetscViewerASCIIPrintf(ngmres->monitor,"residual restart: %e > %e\n",PetscSqrtReal(fAnorm),ngmres->gammaC*PetscSqrtReal(fminnorm));CHKERRQ(ierr); 259 } 260 *selectRestart = PETSC_TRUE; 261 } 262 PetscFunctionReturn(0); 263 } 264