PETSC_EXTERN void PETSC_STDCALL vecpow_(Vec v,PetscScalar *p, int *__ierr ){ *__ierr = VecPow( (Vec)PetscToPointer((v) ),*p); }
extern PetscErrorCode MatLMVMUpdate(Mat M, Vec x, Vec g) { MatLMVMCtx *ctx; PetscReal rhotemp, rhotol; PetscReal y0temp, s0temp; PetscReal yDy, yDs, sDs; PetscReal sigmanew, denom; PetscErrorCode ierr; PetscInt i; PetscBool same; PetscReal yy_sum=0.0, ys_sum=0.0, ss_sum=0.0; PetscFunctionBegin; PetscValidHeaderSpecific(x,VEC_CLASSID,2); PetscValidHeaderSpecific(g,VEC_CLASSID,3); ierr = PetscObjectTypeCompare((PetscObject)M,MATSHELL,&same);CHKERRQ(ierr); if (!same) SETERRQ(PETSC_COMM_SELF,1,"Matrix M is not type MatLMVM"); ierr = MatShellGetContext(M,(void**)&ctx);CHKERRQ(ierr); if (!ctx->allocated) { ierr = MatLMVMAllocateVectors(M, x); CHKERRQ(ierr); } if (0 == ctx->iter) { ierr = MatLMVMReset(M);CHKERRQ(ierr); } else { ierr = VecAYPX(ctx->Gprev,-1.0,g);CHKERRQ(ierr); ierr = VecAYPX(ctx->Xprev,-1.0,x);CHKERRQ(ierr); ierr = VecDot(ctx->Gprev,ctx->Xprev,&rhotemp);CHKERRQ(ierr); ierr = VecDot(ctx->Gprev,ctx->Gprev,&y0temp);CHKERRQ(ierr); rhotol = ctx->eps * y0temp; if (rhotemp > rhotol) { ++ctx->nupdates; ctx->lmnow = PetscMin(ctx->lmnow+1, ctx->lm); ierr=PetscObjectDereference((PetscObject)ctx->S[ctx->lm]);CHKERRQ(ierr); ierr=PetscObjectDereference((PetscObject)ctx->Y[ctx->lm]);CHKERRQ(ierr); for (i = ctx->lm-1; i >= 0; --i) { ctx->S[i+1] = ctx->S[i]; ctx->Y[i+1] = ctx->Y[i]; ctx->rho[i+1] = ctx->rho[i]; } ctx->S[0] = ctx->Xprev; ctx->Y[0] = ctx->Gprev; PetscObjectReference((PetscObject)ctx->S[0]); PetscObjectReference((PetscObject)ctx->Y[0]); ctx->rho[0] = 1.0 / rhotemp; /* Compute the scaling */ switch(ctx->scaleType) { case MatLMVM_Scale_None: break; case MatLMVM_Scale_Scalar: /* Compute s^T s */ ierr = VecDot(ctx->Xprev,ctx->Xprev,&s0temp);CHKERRQ(ierr); /* Scalar is positive; safeguards are not required. */ /* Save information for scalar scaling */ ctx->yy_history[(ctx->nupdates - 1) % ctx->scalar_history] = y0temp; ctx->ys_history[(ctx->nupdates - 1) % ctx->scalar_history] = rhotemp; ctx->ss_history[(ctx->nupdates - 1) % ctx->scalar_history] = s0temp; /* Compute summations for scalar scaling */ yy_sum = 0; /* No safeguard required; y^T y > 0 */ ys_sum = 0; /* No safeguard required; y^T s > 0 */ ss_sum = 0; /* No safeguard required; s^T s > 0 */ for (i = 0; i < PetscMin(ctx->nupdates, ctx->scalar_history); ++i) { yy_sum += ctx->yy_history[i]; ys_sum += ctx->ys_history[i]; ss_sum += ctx->ss_history[i]; } if (0.0 == ctx->s_alpha) { /* Safeguard ys_sum */ if (0.0 == ys_sum) { ys_sum = TAO_ZERO_SAFEGUARD; } sigmanew = ss_sum / ys_sum; } else if (1.0 == ctx->s_alpha) { /* Safeguard yy_sum */ if (0.0 == yy_sum) { yy_sum = TAO_ZERO_SAFEGUARD; } sigmanew = ys_sum / yy_sum; } else { denom = 2*ctx->s_alpha*yy_sum; /* Safeguard denom */ if (0.0 == denom) { denom = TAO_ZERO_SAFEGUARD; } sigmanew = ((2*ctx->s_alpha-1)*ys_sum + PetscSqrtScalar((2*ctx->s_alpha-1)*(2*ctx->s_alpha-1)*ys_sum*ys_sum - 4*(ctx->s_alpha)*(ctx->s_alpha-1)*yy_sum*ss_sum)) / denom; } switch(ctx->limitType) { case MatLMVM_Limit_Average: if (1.0 == ctx->mu) { ctx->sigma = sigmanew; } else if (ctx->mu) { ctx->sigma = ctx->mu * sigmanew + (1.0 - ctx->mu) * ctx->sigma; } break; case MatLMVM_Limit_Relative: if (ctx->mu) { ctx->sigma = TaoMid((1.0 - ctx->mu) * ctx->sigma, sigmanew, (1.0 + ctx->mu) * ctx->sigma); } break; case MatLMVM_Limit_Absolute: if (ctx->nu) { ctx->sigma = TaoMid(ctx->sigma - ctx->nu, sigmanew, ctx->sigma + ctx->nu); } break; default: ctx->sigma = sigmanew; break; } break; case MatLMVM_Scale_Broyden: /* Original version */ /* Combine DFP and BFGS */ /* This code appears to be numerically unstable. We use the */ /* original version because this was used to generate all of */ /* the data and because it may be the least unstable of the */ /* bunch. */ /* P = Q = inv(D); */ ierr = VecCopy(ctx->D,ctx->P);CHKERRQ(ierr); ierr = VecReciprocal(ctx->P);CHKERRQ(ierr); ierr = VecCopy(ctx->P,ctx->Q);CHKERRQ(ierr); /* V = y*y */ ierr = VecPointwiseMult(ctx->V,ctx->Gprev,ctx->Gprev);CHKERRQ(ierr); /* W = inv(D)*s */ ierr = VecPointwiseMult(ctx->W,ctx->Xprev,ctx->P);CHKERRQ(ierr); ierr = VecDot(ctx->W,ctx->Xprev,&sDs);CHKERRQ(ierr); /* Safeguard rhotemp and sDs */ if (0.0 == rhotemp) { rhotemp = TAO_ZERO_SAFEGUARD; } if (0.0 == sDs) { sDs = TAO_ZERO_SAFEGUARD; } if (1.0 != ctx->phi) { /* BFGS portion of the update */ /* U = (inv(D)*s)*(inv(D)*s) */ ierr = VecPointwiseMult(ctx->U,ctx->W,ctx->W);CHKERRQ(ierr); /* Assemble */ ierr = VecAXPY(ctx->P,1.0/rhotemp,ctx->V);CHKERRQ(ierr); ierr = VecAXPY(ctx->P,-1.0/sDs,ctx->U);CHKERRQ(ierr); } if (0.0 != ctx->phi) { /* DFP portion of the update */ /* U = inv(D)*s*y */ ierr = VecPointwiseMult(ctx->U, ctx->W, ctx->Gprev);CHKERRQ(ierr); /* Assemble */ ierr = VecAXPY(ctx->Q,1.0/rhotemp + sDs/(rhotemp*rhotemp), ctx->V);CHKERRQ(ierr); ierr = VecAXPY(ctx->Q,-2.0/rhotemp,ctx->U);CHKERRQ(ierr); } if (0.0 == ctx->phi) { ierr = VecCopy(ctx->P,ctx->U);CHKERRQ(ierr); } else if (1.0 == ctx->phi) { ierr = VecCopy(ctx->Q,ctx->U);CHKERRQ(ierr); } else { /* Broyden update U=(1-phi)*P + phi*Q */ ierr = VecCopy(ctx->Q,ctx->U);CHKERRQ(ierr); ierr = VecAXPBY(ctx->U,1.0-ctx->phi, ctx->phi, ctx->P);CHKERRQ(ierr); } /* Obtain inverse and ensure positive definite */ ierr = VecReciprocal(ctx->U);CHKERRQ(ierr); ierr = VecAbs(ctx->U);CHKERRQ(ierr); switch(ctx->rScaleType) { case MatLMVM_Rescale_None: break; case MatLMVM_Rescale_Scalar: case MatLMVM_Rescale_GL: if (ctx->rScaleType == MatLMVM_Rescale_GL) { /* Gilbert and Lemarachal use the old diagonal */ ierr = VecCopy(ctx->D,ctx->P);CHKERRQ(ierr); } else { /* The default version uses the current diagonal */ ierr = VecCopy(ctx->U,ctx->P);CHKERRQ(ierr); } /* Compute s^T s */ ierr = VecDot(ctx->Xprev,ctx->Xprev,&s0temp);CHKERRQ(ierr); /* Save information for special cases of scalar rescaling */ ctx->yy_rhistory[(ctx->nupdates - 1) % ctx->rescale_history] = y0temp; ctx->ys_rhistory[(ctx->nupdates - 1) % ctx->rescale_history] = rhotemp; ctx->ss_rhistory[(ctx->nupdates - 1) % ctx->rescale_history] = s0temp; if (0.5 == ctx->r_beta) { if (1 == PetscMin(ctx->nupdates, ctx->rescale_history)) { ierr = VecPointwiseMult(ctx->V,ctx->Y[0],ctx->P);CHKERRQ(ierr); ierr = VecDot(ctx->V,ctx->Y[0],&yy_sum);CHKERRQ(ierr); ierr = VecPointwiseDivide(ctx->W,ctx->S[0],ctx->P);CHKERRQ(ierr); ierr = VecDot(ctx->W,ctx->S[0],&ss_sum);CHKERRQ(ierr); ys_sum = ctx->ys_rhistory[0]; } else { ierr = VecCopy(ctx->P,ctx->Q);CHKERRQ(ierr); ierr = VecReciprocal(ctx->Q);CHKERRQ(ierr); /* Compute summations for scalar scaling */ yy_sum = 0; /* No safeguard required */ ys_sum = 0; /* No safeguard required */ ss_sum = 0; /* No safeguard required */ for (i = 0; i < PetscMin(ctx->nupdates, ctx->rescale_history); ++i) { ierr = VecPointwiseMult(ctx->V,ctx->Y[i],ctx->P);CHKERRQ(ierr); ierr = VecDot(ctx->V,ctx->Y[i],&yDy);CHKERRQ(ierr); yy_sum += yDy; ierr = VecPointwiseMult(ctx->W,ctx->S[i],ctx->Q);CHKERRQ(ierr); ierr = VecDot(ctx->W,ctx->S[i],&sDs);CHKERRQ(ierr); ss_sum += sDs; ys_sum += ctx->ys_rhistory[i]; } } } else if (0.0 == ctx->r_beta) { if (1 == PetscMin(ctx->nupdates, ctx->rescale_history)) { /* Compute summations for scalar scaling */ ierr = VecPointwiseDivide(ctx->W,ctx->S[0],ctx->P);CHKERRQ(ierr); ierr = VecDot(ctx->W, ctx->Y[0], &ys_sum);CHKERRQ(ierr); ierr = VecDot(ctx->W, ctx->W, &ss_sum);CHKERRQ(ierr); yy_sum += ctx->yy_rhistory[0]; } else { ierr = VecCopy(ctx->Q, ctx->P);CHKERRQ(ierr); ierr = VecReciprocal(ctx->Q);CHKERRQ(ierr); /* Compute summations for scalar scaling */ yy_sum = 0; /* No safeguard required */ ys_sum = 0; /* No safeguard required */ ss_sum = 0; /* No safeguard required */ for (i = 0; i < PetscMin(ctx->nupdates, ctx->rescale_history); ++i) { ierr = VecPointwiseMult(ctx->W, ctx->S[i], ctx->Q);CHKERRQ(ierr); ierr = VecDot(ctx->W, ctx->Y[i], &yDs);CHKERRQ(ierr); ys_sum += yDs; ierr = VecDot(ctx->W, ctx->W, &sDs);CHKERRQ(ierr); ss_sum += sDs; yy_sum += ctx->yy_rhistory[i]; } } } else if (1.0 == ctx->r_beta) { /* Compute summations for scalar scaling */ yy_sum = 0; /* No safeguard required */ ys_sum = 0; /* No safeguard required */ ss_sum = 0; /* No safeguard required */ for (i = 0; i < PetscMin(ctx->nupdates, ctx->rescale_history); ++i) { ierr = VecPointwiseMult(ctx->V, ctx->Y[i], ctx->P);CHKERRQ(ierr); ierr = VecDot(ctx->V, ctx->S[i], &yDs);CHKERRQ(ierr); ys_sum += yDs; ierr = VecDot(ctx->V, ctx->V, &yDy);CHKERRQ(ierr); yy_sum += yDy; ss_sum += ctx->ss_rhistory[i]; } } else { ierr = VecCopy(ctx->Q, ctx->P);CHKERRQ(ierr); ierr = VecPow(ctx->P, ctx->r_beta);CHKERRQ(ierr); ierr = VecPointwiseDivide(ctx->Q, ctx->P, ctx->Q);CHKERRQ(ierr); /* Compute summations for scalar scaling */ yy_sum = 0; /* No safeguard required */ ys_sum = 0; /* No safeguard required */ ss_sum = 0; /* No safeguard required */ for (i = 0; i < PetscMin(ctx->nupdates, ctx->rescale_history); ++i) { ierr = VecPointwiseMult(ctx->V, ctx->P, ctx->Y[i]);CHKERRQ(ierr); ierr = VecPointwiseMult(ctx->W, ctx->Q, ctx->S[i]);CHKERRQ(ierr); ierr = VecDot(ctx->V, ctx->V, &yDy);CHKERRQ(ierr); ierr = VecDot(ctx->V, ctx->W, &yDs);CHKERRQ(ierr); ierr = VecDot(ctx->W, ctx->W, &sDs);CHKERRQ(ierr); yy_sum += yDy; ys_sum += yDs; ss_sum += sDs; } } if (0.0 == ctx->r_alpha) { /* Safeguard ys_sum */ if (0.0 == ys_sum) { ys_sum = TAO_ZERO_SAFEGUARD; } sigmanew = ss_sum / ys_sum; } else if (1.0 == ctx->r_alpha) { /* Safeguard yy_sum */ if (0.0 == yy_sum) { ys_sum = TAO_ZERO_SAFEGUARD; } sigmanew = ys_sum / yy_sum; } else { denom = 2*ctx->r_alpha*yy_sum; /* Safeguard denom */ if (0.0 == denom) { denom = TAO_ZERO_SAFEGUARD; } sigmanew = ((2*ctx->r_alpha-1)*ys_sum + PetscSqrtScalar((2*ctx->r_alpha-1)*(2*ctx->r_alpha-1)*ys_sum*ys_sum - 4*ctx->r_alpha*(ctx->r_alpha-1)*yy_sum*ss_sum)) / denom; } /* If Q has small values, then Q^(r_beta - 1) */ /* can have very large values. Hence, ys_sum */ /* and ss_sum can be infinity. In this case, */ /* sigmanew can either be not-a-number or infinity. */ if (PetscIsInfOrNanReal(sigmanew)) { /* sigmanew is not-a-number; skip rescaling */ } else if (!sigmanew) { /* sigmanew is zero; this is a bad case; skip rescaling */ } else { /* sigmanew is positive */ ierr = VecScale(ctx->U, sigmanew);CHKERRQ(ierr); } break; } /* Modify for previous information */ switch(ctx->limitType) { case MatLMVM_Limit_Average: if (1.0 == ctx->mu) { ierr = VecCopy(ctx->D, ctx->U);CHKERRQ(ierr); } else if (ctx->mu) { ierr = VecAXPBY(ctx->D,ctx->mu, 1.0-ctx->mu,ctx->U);CHKERRQ(ierr); } break; case MatLMVM_Limit_Relative: if (ctx->mu) { /* P = (1-mu) * D */ ierr = VecAXPBY(ctx->P, 1.0-ctx->mu, 0.0, ctx->D);CHKERRQ(ierr); /* Q = (1+mu) * D */ ierr = VecAXPBY(ctx->Q, 1.0+ctx->mu, 0.0, ctx->D);CHKERRQ(ierr); ierr = VecMedian(ctx->P, ctx->U, ctx->Q, ctx->D);CHKERRQ(ierr); } break; case MatLMVM_Limit_Absolute: if (ctx->nu) { ierr = VecCopy(ctx->P, ctx->D);CHKERRQ(ierr); ierr = VecShift(ctx->P, -ctx->nu);CHKERRQ(ierr); ierr = VecCopy(ctx->D, ctx->Q);CHKERRQ(ierr); ierr = VecShift(ctx->Q, ctx->nu);CHKERRQ(ierr); ierr = VecMedian(ctx->P, ctx->U, ctx->Q, ctx->P);CHKERRQ(ierr); } break; default: ierr = VecCopy(ctx->U, ctx->D);CHKERRQ(ierr); break; } break; } ierr = PetscObjectDereference((PetscObject)ctx->Xprev);CHKERRQ(ierr); ierr = PetscObjectDereference((PetscObject)ctx->Gprev);CHKERRQ(ierr); ctx->Xprev = ctx->S[ctx->lm]; ctx->Gprev = ctx->Y[ctx->lm]; ierr = PetscObjectReference((PetscObject)ctx->S[ctx->lm]);CHKERRQ(ierr); ierr = PetscObjectReference((PetscObject)ctx->Y[ctx->lm]);CHKERRQ(ierr); } else { ++ctx->nrejects; } } ++ctx->iter; ierr = VecCopy(x, ctx->Xprev);CHKERRQ(ierr); ierr = VecCopy(g, ctx->Gprev);CHKERRQ(ierr); PetscFunctionReturn(0); }