1*2b4ed282SShri Abhyankar #define PETSCSNES_DLL 2*2b4ed282SShri Abhyankar 3*2b4ed282SShri Abhyankar #include "../src/snes/impls/vi/viimpl.h" 4*2b4ed282SShri Abhyankar 5*2b4ed282SShri Abhyankar /* 6*2b4ed282SShri Abhyankar Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function, 7*2b4ed282SShri Abhyankar || F(u) ||_2 but not a zero, F(u) = 0. In the case when one cannot compute J^T F we use the fact that 8*2b4ed282SShri Abhyankar 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 9*2b4ed282SShri Abhyankar for this trick. One assumes that the probability that W is in the null space of J is very, very small. 10*2b4ed282SShri Abhyankar */ 11*2b4ed282SShri Abhyankar #undef __FUNCT__ 12*2b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckLocalMin_Private" 13*2b4ed282SShri Abhyankar PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin) 14*2b4ed282SShri Abhyankar { 15*2b4ed282SShri Abhyankar PetscReal a1; 16*2b4ed282SShri Abhyankar PetscErrorCode ierr; 17*2b4ed282SShri Abhyankar PetscTruth hastranspose; 18*2b4ed282SShri Abhyankar 19*2b4ed282SShri Abhyankar PetscFunctionBegin; 20*2b4ed282SShri Abhyankar *ismin = PETSC_FALSE; 21*2b4ed282SShri Abhyankar ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 22*2b4ed282SShri Abhyankar if (hastranspose) { 23*2b4ed282SShri Abhyankar /* Compute || J^T F|| */ 24*2b4ed282SShri Abhyankar ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 25*2b4ed282SShri Abhyankar ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 26*2b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"|| J^T F|| %G near zero implies found a local minimum\n",a1/fnorm);CHKERRQ(ierr); 27*2b4ed282SShri Abhyankar if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 28*2b4ed282SShri Abhyankar } else { 29*2b4ed282SShri Abhyankar Vec work; 30*2b4ed282SShri Abhyankar PetscScalar result; 31*2b4ed282SShri Abhyankar PetscReal wnorm; 32*2b4ed282SShri Abhyankar 33*2b4ed282SShri Abhyankar ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr); 34*2b4ed282SShri Abhyankar ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr); 35*2b4ed282SShri Abhyankar ierr = VecDuplicate(W,&work);CHKERRQ(ierr); 36*2b4ed282SShri Abhyankar ierr = MatMult(A,W,work);CHKERRQ(ierr); 37*2b4ed282SShri Abhyankar ierr = VecDot(F,work,&result);CHKERRQ(ierr); 38*2b4ed282SShri Abhyankar ierr = VecDestroy(work);CHKERRQ(ierr); 39*2b4ed282SShri Abhyankar a1 = PetscAbsScalar(result)/(fnorm*wnorm); 40*2b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %G near zero implies found a local minimum\n",a1);CHKERRQ(ierr); 41*2b4ed282SShri Abhyankar if (a1 < 1.e-4) *ismin = PETSC_TRUE; 42*2b4ed282SShri Abhyankar } 43*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 44*2b4ed282SShri Abhyankar } 45*2b4ed282SShri Abhyankar 46*2b4ed282SShri Abhyankar /* 47*2b4ed282SShri Abhyankar Checks if J^T(F - J*X) = 0 48*2b4ed282SShri Abhyankar */ 49*2b4ed282SShri Abhyankar #undef __FUNCT__ 50*2b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckResidual_Private" 51*2b4ed282SShri Abhyankar PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2) 52*2b4ed282SShri Abhyankar { 53*2b4ed282SShri Abhyankar PetscReal a1,a2; 54*2b4ed282SShri Abhyankar PetscErrorCode ierr; 55*2b4ed282SShri Abhyankar PetscTruth hastranspose; 56*2b4ed282SShri Abhyankar 57*2b4ed282SShri Abhyankar PetscFunctionBegin; 58*2b4ed282SShri Abhyankar ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 59*2b4ed282SShri Abhyankar if (hastranspose) { 60*2b4ed282SShri Abhyankar ierr = MatMult(A,X,W1);CHKERRQ(ierr); 61*2b4ed282SShri Abhyankar ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr); 62*2b4ed282SShri Abhyankar 63*2b4ed282SShri Abhyankar /* Compute || J^T W|| */ 64*2b4ed282SShri Abhyankar ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 65*2b4ed282SShri Abhyankar ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 66*2b4ed282SShri Abhyankar ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 67*2b4ed282SShri Abhyankar if (a1 != 0.0) { 68*2b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %G near zero implies inconsistent rhs\n",a2/a1);CHKERRQ(ierr); 69*2b4ed282SShri Abhyankar } 70*2b4ed282SShri Abhyankar } 71*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 72*2b4ed282SShri Abhyankar } 73*2b4ed282SShri Abhyankar 74*2b4ed282SShri Abhyankar /* 75*2b4ed282SShri Abhyankar SNESDefaultConverged_VI - Checks the convergence of the semismooth newton algorithm. 76*2b4ed282SShri Abhyankar 77*2b4ed282SShri Abhyankar Notes: 78*2b4ed282SShri Abhyankar The convergence criterion currently implemented is 79*2b4ed282SShri Abhyankar merit < abstol 80*2b4ed282SShri Abhyankar merit < rtol*merit_initial 81*2b4ed282SShri Abhyankar */ 82*2b4ed282SShri Abhyankar #undef __FUNCT__ 83*2b4ed282SShri Abhyankar #define __FUNCT__ "SNESDefaultConverged_VI" 84*2b4ed282SShri Abhyankar PetscErrorCode SNESDefaultConverged_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal merit,SNESConvergedReason *reason,void *dummy) 85*2b4ed282SShri Abhyankar { 86*2b4ed282SShri Abhyankar PetscErrorCode ierr; 87*2b4ed282SShri Abhyankar 88*2b4ed282SShri Abhyankar PetscFunctionBegin; 89*2b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 90*2b4ed282SShri Abhyankar PetscValidPointer(reason,6); 91*2b4ed282SShri Abhyankar 92*2b4ed282SShri Abhyankar *reason = SNES_CONVERGED_ITERATING; 93*2b4ed282SShri Abhyankar 94*2b4ed282SShri Abhyankar if (!it) { 95*2b4ed282SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 96*2b4ed282SShri Abhyankar snes->ttol = merit*snes->rtol; 97*2b4ed282SShri Abhyankar } 98*2b4ed282SShri Abhyankar if (merit != merit) { 99*2b4ed282SShri Abhyankar ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr); 100*2b4ed282SShri Abhyankar *reason = SNES_DIVERGED_FNORM_NAN; 101*2b4ed282SShri Abhyankar } else if (merit < snes->abstol) { 102*2b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"Converged due to merit function %G < %G\n",merit,snes->abstol);CHKERRQ(ierr); 103*2b4ed282SShri Abhyankar *reason = SNES_CONVERGED_FNORM_ABS; 104*2b4ed282SShri Abhyankar } else if (snes->nfuncs >= snes->max_funcs) { 105*2b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr); 106*2b4ed282SShri Abhyankar *reason = SNES_DIVERGED_FUNCTION_COUNT; 107*2b4ed282SShri Abhyankar } 108*2b4ed282SShri Abhyankar 109*2b4ed282SShri Abhyankar if (it && !*reason) { 110*2b4ed282SShri Abhyankar if (merit < snes->ttol) { 111*2b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"Converged due to merit function %G < %G (relative tolerance)\n",merit,snes->ttol);CHKERRQ(ierr); 112*2b4ed282SShri Abhyankar *reason = SNES_CONVERGED_FNORM_RELATIVE; 113*2b4ed282SShri Abhyankar } 114*2b4ed282SShri Abhyankar } 115*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 116*2b4ed282SShri Abhyankar } 117*2b4ed282SShri Abhyankar 118*2b4ed282SShri Abhyankar /* 119*2b4ed282SShri Abhyankar SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem. 120*2b4ed282SShri Abhyankar 121*2b4ed282SShri Abhyankar Input Parameter: 122*2b4ed282SShri Abhyankar . phi - the semismooth function 123*2b4ed282SShri Abhyankar 124*2b4ed282SShri Abhyankar Output Parameter: 125*2b4ed282SShri Abhyankar . merit - the merit function 126*2b4ed282SShri Abhyankar . phinorm - ||phi|| 127*2b4ed282SShri Abhyankar 128*2b4ed282SShri Abhyankar Notes: 129*2b4ed282SShri Abhyankar The merit function for the mixed complementarity problem is defined as 130*2b4ed282SShri Abhyankar merit = 0.5*phi^T*phi 131*2b4ed282SShri Abhyankar */ 132*2b4ed282SShri Abhyankar #undef __FUNCT__ 133*2b4ed282SShri Abhyankar #define __FUNCT__ "SNESVIComputeMeritFunction" 134*2b4ed282SShri Abhyankar static PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal* merit,PetscReal* phinorm) 135*2b4ed282SShri Abhyankar { 136*2b4ed282SShri Abhyankar PetscErrorCode ierr; 137*2b4ed282SShri Abhyankar 138*2b4ed282SShri Abhyankar PetscFunctionBegin; 139*2b4ed282SShri Abhyankar ierr = VecNormBegin(phi,NORM_2,phinorm); 140*2b4ed282SShri Abhyankar ierr = VecNormEnd(phi,NORM_2,phinorm); 141*2b4ed282SShri Abhyankar 142*2b4ed282SShri Abhyankar *merit = 0.5*(*phinorm)*(*phinorm); 143*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 144*2b4ed282SShri Abhyankar } 145*2b4ed282SShri Abhyankar 146*2b4ed282SShri Abhyankar /* 147*2b4ed282SShri Abhyankar ComputeFischerFunction - Computes the semismooth fischer burmeister function for a mixed complementarity equation. 148*2b4ed282SShri Abhyankar 149*2b4ed282SShri Abhyankar Notes: 150*2b4ed282SShri Abhyankar The Fischer-Burmeister function is defined as 151*2b4ed282SShri Abhyankar ff(a,b) := sqrt(a*a + b*b) - a - b 152*2b4ed282SShri Abhyankar and is used reformulate a complementarity equation as a semismooth equation. 153*2b4ed282SShri Abhyankar */ 154*2b4ed282SShri Abhyankar 155*2b4ed282SShri Abhyankar #undef __FUNCT__ 156*2b4ed282SShri Abhyankar #define __FUNCT__ "ComputeFischerFunction" 157*2b4ed282SShri Abhyankar static PetscErrorCode ComputeFischerFunction(PetscScalar a, PetscScalar b, PetscScalar* ff) 158*2b4ed282SShri Abhyankar { 159*2b4ed282SShri Abhyankar PetscFunctionBegin; 160*2b4ed282SShri Abhyankar *ff = sqrt(a*a + b*b) - a - b; 161*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 162*2b4ed282SShri Abhyankar } 163*2b4ed282SShri Abhyankar 164*2b4ed282SShri Abhyankar /* 165*2b4ed282SShri Abhyankar SNESVIComputeSSFunction - Reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear equations in semismooth form. 166*2b4ed282SShri Abhyankar 167*2b4ed282SShri Abhyankar Input Parameters: 168*2b4ed282SShri Abhyankar . snes - the SNES context 169*2b4ed282SShri Abhyankar . x - current iterate 170*2b4ed282SShri Abhyankar . functx - user defined function context 171*2b4ed282SShri Abhyankar 172*2b4ed282SShri Abhyankar Output Parameters: 173*2b4ed282SShri Abhyankar . phi - Semismooth function 174*2b4ed282SShri Abhyankar 175*2b4ed282SShri Abhyankar Notes: 176*2b4ed282SShri Abhyankar The result of this function is done by cases: 177*2b4ed282SShri Abhyankar + l[i] == -infinity, u[i] == infinity -- phi[i] = -f[i] 178*2b4ed282SShri Abhyankar . l[i] == -infinity, u[i] finite -- phi[i] = ss(u[i]-x[i], -f[i]) 179*2b4ed282SShri Abhyankar . l[i] finite, u[i] == infinity -- phi[i] = ss(x[i]-l[i], f[i]) 180*2b4ed282SShri Abhyankar . l[i] finite < u[i] finite -- phi[i] = phi(x[i]-l[i], ss(u[i]-x[i], -f[u])) 181*2b4ed282SShri Abhyankar - otherwise l[i] == u[i] -- phi[i] = l[i] - x[i] 182*2b4ed282SShri Abhyankar ss is the semismoothing function used to reformulate the nonlinear equation in complementarity 183*2b4ed282SShri Abhyankar form to semismooth form 184*2b4ed282SShri Abhyankar 185*2b4ed282SShri Abhyankar */ 186*2b4ed282SShri Abhyankar #undef __FUNCT__ 187*2b4ed282SShri Abhyankar #define __FUNCT__ "SNESVIComputeSSFunction" 188*2b4ed282SShri Abhyankar static PetscErrorCode SNESVIComputeSSFunction(SNES snes,Vec X,Vec phi,void* functx) 189*2b4ed282SShri Abhyankar { 190*2b4ed282SShri Abhyankar PetscErrorCode ierr; 191*2b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 192*2b4ed282SShri Abhyankar Vec Xl = vi->xl,Xu = vi->xu,F = snes->vec_func; 193*2b4ed282SShri Abhyankar PetscScalar *phi_arr,*x_arr,*f_arr,*l,*u,t; 194*2b4ed282SShri Abhyankar PetscInt i,nlocal; 195*2b4ed282SShri Abhyankar 196*2b4ed282SShri Abhyankar PetscFunctionBegin; 197*2b4ed282SShri Abhyankar ierr = (*vi->computeuserfunction)(snes,X,F,functx);CHKERRQ(ierr); 198*2b4ed282SShri Abhyankar ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 199*2b4ed282SShri Abhyankar 200*2b4ed282SShri Abhyankar ierr = VecGetArray(X,&x_arr);CHKERRQ(ierr); 201*2b4ed282SShri Abhyankar ierr = VecGetArray(F,&f_arr);CHKERRQ(ierr); 202*2b4ed282SShri Abhyankar ierr = VecGetArray(Xl,&l);CHKERRQ(ierr); 203*2b4ed282SShri Abhyankar ierr = VecGetArray(Xu,&u);CHKERRQ(ierr); 204*2b4ed282SShri Abhyankar ierr = VecGetArray(phi,&phi_arr);CHKERRQ(ierr); 205*2b4ed282SShri Abhyankar 206*2b4ed282SShri Abhyankar for (i=0;i < nlocal;i++) { 207*2b4ed282SShri Abhyankar if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) { 208*2b4ed282SShri Abhyankar phi_arr[i] = -f_arr[i]; 209*2b4ed282SShri Abhyankar } 210*2b4ed282SShri Abhyankar else if (l[i] <= PETSC_VI_NINF) { 211*2b4ed282SShri Abhyankar t = u[i] - x_arr[i]; 212*2b4ed282SShri Abhyankar ierr = ComputeFischerFunction(t,-f_arr[i],&phi_arr[i]);CHKERRQ(ierr); 213*2b4ed282SShri Abhyankar phi_arr[i] = -phi_arr[i]; 214*2b4ed282SShri Abhyankar } 215*2b4ed282SShri Abhyankar else if (u[i] >= PETSC_VI_INF) { 216*2b4ed282SShri Abhyankar t = x_arr[i] - l[i]; 217*2b4ed282SShri Abhyankar ierr = ComputeFischerFunction(t,f_arr[i],&phi_arr[i]);CHKERRQ(ierr); 218*2b4ed282SShri Abhyankar } 219*2b4ed282SShri Abhyankar else if (l[i] == u[i]) { 220*2b4ed282SShri Abhyankar phi_arr[i] = l[i] - x_arr[i]; 221*2b4ed282SShri Abhyankar } 222*2b4ed282SShri Abhyankar else { 223*2b4ed282SShri Abhyankar t = u[i] - x_arr[i]; 224*2b4ed282SShri Abhyankar ierr = ComputeFischerFunction(t,-f_arr[i],&phi_arr[i]); 225*2b4ed282SShri Abhyankar t = x_arr[i] - l[i]; 226*2b4ed282SShri Abhyankar ierr = ComputeFischerFunction(t,phi_arr[i],&phi_arr[i]); 227*2b4ed282SShri Abhyankar } 228*2b4ed282SShri Abhyankar } 229*2b4ed282SShri Abhyankar 230*2b4ed282SShri Abhyankar ierr = VecRestoreArray(X,&x_arr);CHKERRQ(ierr); 231*2b4ed282SShri Abhyankar ierr = VecRestoreArray(F,&f_arr);CHKERRQ(ierr); 232*2b4ed282SShri Abhyankar ierr = VecRestoreArray(Xl,&l);CHKERRQ(ierr); 233*2b4ed282SShri Abhyankar ierr = VecRestoreArray(Xu,&u);CHKERRQ(ierr); 234*2b4ed282SShri Abhyankar ierr = VecRestoreArray(phi,&phi_arr);CHKERRQ(ierr); 235*2b4ed282SShri Abhyankar 236*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 237*2b4ed282SShri Abhyankar } 238*2b4ed282SShri Abhyankar 239*2b4ed282SShri Abhyankar /* 240*2b4ed282SShri Abhyankar SNESVIComputeSSJacobian - Computes the jacobian of the semismooth function.The Jacobian for the semismooth function is an element of the B-subdifferential of the Fischer-Burmeister function for complementarity problems. 241*2b4ed282SShri Abhyankar 242*2b4ed282SShri Abhyankar Input Parameters: 243*2b4ed282SShri Abhyankar . snes - the SNES context 244*2b4ed282SShri Abhyankar . X - the current iterate 245*2b4ed282SShri Abhyankar . vec_func - nonlinear function evaluated at x 246*2b4ed282SShri Abhyankar 247*2b4ed282SShri Abhyankar Output Parameters: 248*2b4ed282SShri Abhyankar . jac - semismooth jacobian 249*2b4ed282SShri Abhyankar . jac_pre - optional preconditioning matrix 250*2b4ed282SShri Abhyankar . flag - flag passed on by SNESComputeJacobian. 251*2b4ed282SShri Abhyankar . jacctx - user provided jacobian context 252*2b4ed282SShri Abhyankar 253*2b4ed282SShri Abhyankar Notes: 254*2b4ed282SShri Abhyankar The semismooth jacobian matrix is given by 255*2b4ed282SShri Abhyankar jac = Da + Db*jacfun 256*2b4ed282SShri Abhyankar where Db is the row scaling matrix stored as a vector, 257*2b4ed282SShri Abhyankar Da is the diagonal perturbation matrix stored as a vector 258*2b4ed282SShri Abhyankar and jacfun is the jacobian of the original nonlinear function. 259*2b4ed282SShri Abhyankar */ 260*2b4ed282SShri Abhyankar #undef __FUNCT__ 261*2b4ed282SShri Abhyankar #define __FUNCT__ "SNESVIComputeSSJacobian" 262*2b4ed282SShri Abhyankar PetscErrorCode SNESVIComputeSSJacobian(SNES snes,Vec X,Mat *jac, Mat *jac_pre, MatStructure *flg,void* jacctx) 263*2b4ed282SShri Abhyankar { 264*2b4ed282SShri Abhyankar PetscErrorCode ierr; 265*2b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 266*2b4ed282SShri Abhyankar PetscScalar *l,*u,*x,*f,*da,*db,*z,*t,t1,t2,ci,di,ei; 267*2b4ed282SShri Abhyankar PetscInt i,nlocal; 268*2b4ed282SShri Abhyankar Vec F = snes->vec_func; 269*2b4ed282SShri Abhyankar 270*2b4ed282SShri Abhyankar PetscFunctionBegin; 271*2b4ed282SShri Abhyankar 272*2b4ed282SShri Abhyankar ierr = (*vi->computeuserjacobian)(snes,X,jac,jac_pre,flg,jacctx);CHKERRQ(ierr); 273*2b4ed282SShri Abhyankar 274*2b4ed282SShri Abhyankar ierr = VecGetArray(X,&x);CHKERRQ(ierr); 275*2b4ed282SShri Abhyankar ierr = VecGetArray(F,&f);CHKERRQ(ierr); 276*2b4ed282SShri Abhyankar ierr = VecGetArray(vi->xl,&l);CHKERRQ(ierr); 277*2b4ed282SShri Abhyankar ierr = VecGetArray(vi->xu,&u);CHKERRQ(ierr); 278*2b4ed282SShri Abhyankar ierr = VecGetArray(vi->Da,&da);CHKERRQ(ierr); 279*2b4ed282SShri Abhyankar ierr = VecGetArray(vi->Db,&db);CHKERRQ(ierr); 280*2b4ed282SShri Abhyankar ierr = VecGetArray(vi->z,&z);CHKERRQ(ierr); 281*2b4ed282SShri Abhyankar 282*2b4ed282SShri Abhyankar ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 283*2b4ed282SShri Abhyankar /* Set the elements of the vector z: 284*2b4ed282SShri Abhyankar z[i] = 1 if (x[i] - l[i],f[i]) = (0,0) or (u[i] - x[i],f[i]) = (0,0) 285*2b4ed282SShri Abhyankar else z[i] = 0 286*2b4ed282SShri Abhyankar */ 287*2b4ed282SShri Abhyankar for(i=0;i < nlocal;i++) { 288*2b4ed282SShri Abhyankar da[i] = db[i] = z[i] = 0; 289*2b4ed282SShri Abhyankar if(PetscAbsScalar(f[i]) <= vi->const_tol) { 290*2b4ed282SShri Abhyankar if ((l[i] > PETSC_VI_NINF) && (PetscAbsScalar(x[i]-l[i]) <= vi->const_tol)) { 291*2b4ed282SShri Abhyankar da[i] = 1; 292*2b4ed282SShri Abhyankar z[i] = 1; 293*2b4ed282SShri Abhyankar } 294*2b4ed282SShri Abhyankar if ((u[i] < PETSC_VI_INF) && (PetscAbsScalar(u[i]-x[i]) <= vi->const_tol)) { 295*2b4ed282SShri Abhyankar db[i] = 1; 296*2b4ed282SShri Abhyankar z[i] = 1; 297*2b4ed282SShri Abhyankar } 298*2b4ed282SShri Abhyankar } 299*2b4ed282SShri Abhyankar } 300*2b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->z,&z);CHKERRQ(ierr); 301*2b4ed282SShri Abhyankar ierr = MatMult(*jac,vi->z,vi->t);CHKERRQ(ierr); 302*2b4ed282SShri Abhyankar ierr = VecGetArray(vi->t,&t);CHKERRQ(ierr); 303*2b4ed282SShri Abhyankar /* Compute the elements of the diagonal perturbation vector Da and row scaling vector Db */ 304*2b4ed282SShri Abhyankar for(i=0;i< nlocal;i++) { 305*2b4ed282SShri Abhyankar /* Free variables */ 306*2b4ed282SShri Abhyankar if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) { 307*2b4ed282SShri Abhyankar da[i] = 0; db[i] = -1; 308*2b4ed282SShri Abhyankar } 309*2b4ed282SShri Abhyankar /* lower bounded variables */ 310*2b4ed282SShri Abhyankar else if (u[i] >= PETSC_VI_INF) { 311*2b4ed282SShri Abhyankar if (da[i] >= 1) { 312*2b4ed282SShri Abhyankar t2 = PetscScalarNorm(1,t[i]); 313*2b4ed282SShri Abhyankar da[i] = 1/t2 - 1; 314*2b4ed282SShri Abhyankar db[i] = t[i]/t2 - 1; 315*2b4ed282SShri Abhyankar } else { 316*2b4ed282SShri Abhyankar t1 = x[i] - l[i]; 317*2b4ed282SShri Abhyankar t2 = PetscScalarNorm(t1,f[i]); 318*2b4ed282SShri Abhyankar da[i] = t1/t2 - 1; 319*2b4ed282SShri Abhyankar db[i] = f[i]/t2 - 1; 320*2b4ed282SShri Abhyankar } 321*2b4ed282SShri Abhyankar } 322*2b4ed282SShri Abhyankar /* upper bounded variables */ 323*2b4ed282SShri Abhyankar else if (l[i] <= PETSC_VI_NINF) { 324*2b4ed282SShri Abhyankar if (db[i] >= 1) { 325*2b4ed282SShri Abhyankar t2 = PetscScalarNorm(1,t[i]); 326*2b4ed282SShri Abhyankar da[i] = -1/t2 -1; 327*2b4ed282SShri Abhyankar db[i] = -t[i]/t2 - 1; 328*2b4ed282SShri Abhyankar } 329*2b4ed282SShri Abhyankar else { 330*2b4ed282SShri Abhyankar t1 = u[i] - x[i]; 331*2b4ed282SShri Abhyankar t2 = PetscScalarNorm(t1,f[i]); 332*2b4ed282SShri Abhyankar da[i] = t1/t2 - 1; 333*2b4ed282SShri Abhyankar db[i] = -f[i]/t2 - 1; 334*2b4ed282SShri Abhyankar } 335*2b4ed282SShri Abhyankar } 336*2b4ed282SShri Abhyankar /* Fixed variables */ 337*2b4ed282SShri Abhyankar else if (l[i] == u[i]) { 338*2b4ed282SShri Abhyankar da[i] = -1; 339*2b4ed282SShri Abhyankar db[i] = 0; 340*2b4ed282SShri Abhyankar } 341*2b4ed282SShri Abhyankar /* Box constrained variables */ 342*2b4ed282SShri Abhyankar else { 343*2b4ed282SShri Abhyankar if (db[i] >= 1) { 344*2b4ed282SShri Abhyankar t2 = PetscScalarNorm(1,t[i]); 345*2b4ed282SShri Abhyankar ci = 1/t2 + 1; 346*2b4ed282SShri Abhyankar di = t[i]/t2 + 1; 347*2b4ed282SShri Abhyankar } 348*2b4ed282SShri Abhyankar else { 349*2b4ed282SShri Abhyankar t1 = x[i] - u[i]; 350*2b4ed282SShri Abhyankar t2 = PetscScalarNorm(t1,f[i]); 351*2b4ed282SShri Abhyankar ci = t1/t2 + 1; 352*2b4ed282SShri Abhyankar di = f[i]/t2 + 1; 353*2b4ed282SShri Abhyankar } 354*2b4ed282SShri Abhyankar 355*2b4ed282SShri Abhyankar if (da[i] >= 1) { 356*2b4ed282SShri Abhyankar t1 = ci + di*t[i]; 357*2b4ed282SShri Abhyankar t2 = PetscScalarNorm(1,t1); 358*2b4ed282SShri Abhyankar t1 = t1/t2 - 1; 359*2b4ed282SShri Abhyankar t2 = 1/t2 - 1; 360*2b4ed282SShri Abhyankar } 361*2b4ed282SShri Abhyankar else { 362*2b4ed282SShri Abhyankar ierr = ComputeFischerFunction(u[i]-x[i],-f[i],&ei);CHKERRQ(ierr); 363*2b4ed282SShri Abhyankar t2 = PetscScalarNorm(x[i]-l[i],ei); 364*2b4ed282SShri Abhyankar t1 = ei/t2 - 1; 365*2b4ed282SShri Abhyankar t2 = (x[i] - l[i])/t2 - 1; 366*2b4ed282SShri Abhyankar } 367*2b4ed282SShri Abhyankar 368*2b4ed282SShri Abhyankar da[i] = t2 + t1*ci; 369*2b4ed282SShri Abhyankar db[i] = t1*di; 370*2b4ed282SShri Abhyankar } 371*2b4ed282SShri Abhyankar } 372*2b4ed282SShri Abhyankar ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 373*2b4ed282SShri Abhyankar ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 374*2b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->xl,&l);CHKERRQ(ierr); 375*2b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->xu,&u);CHKERRQ(ierr); 376*2b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->Da,&da);CHKERRQ(ierr); 377*2b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->Db,&db);CHKERRQ(ierr); 378*2b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->t,&t);CHKERRQ(ierr); 379*2b4ed282SShri Abhyankar 380*2b4ed282SShri Abhyankar /* Do row scaling and add diagonal perturbation */ 381*2b4ed282SShri Abhyankar ierr = MatDiagonalScale(*jac,vi->Db,PETSC_NULL);CHKERRQ(ierr); 382*2b4ed282SShri Abhyankar ierr = MatDiagonalSet(*jac,vi->Da,ADD_VALUES);CHKERRQ(ierr); 383*2b4ed282SShri Abhyankar if (*jac != *jac_pre) { /* If jac and jac_pre are different */ 384*2b4ed282SShri Abhyankar ierr = MatDiagonalScale(*jac_pre,vi->Db,PETSC_NULL); 385*2b4ed282SShri Abhyankar ierr = MatDiagonalSet(*jac_pre,vi->Da,ADD_VALUES);CHKERRQ(ierr); 386*2b4ed282SShri Abhyankar } 387*2b4ed282SShri Abhyankar 388*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 389*2b4ed282SShri Abhyankar } 390*2b4ed282SShri Abhyankar 391*2b4ed282SShri Abhyankar /* 392*2b4ed282SShri Abhyankar SNESVIComputeMeritFunctionGradient - Computes the gradient of the merit function psi. 393*2b4ed282SShri Abhyankar 394*2b4ed282SShri Abhyankar Input Parameters: 395*2b4ed282SShri Abhyankar . phi - semismooth function. 396*2b4ed282SShri Abhyankar . H - semismooth jacobian 397*2b4ed282SShri Abhyankar 398*2b4ed282SShri Abhyankar Output Parameters: 399*2b4ed282SShri Abhyankar . dpsi - merit function gradient 400*2b4ed282SShri Abhyankar 401*2b4ed282SShri Abhyankar Notes: 402*2b4ed282SShri Abhyankar The merit function gradient is computed as follows 403*2b4ed282SShri Abhyankar dpsi = H^T*phi 404*2b4ed282SShri Abhyankar */ 405*2b4ed282SShri Abhyankar #undef __FUNCT__ 406*2b4ed282SShri Abhyankar #define __FUNCT__ "SNESVIComputeMeritFunctionGradient" 407*2b4ed282SShri Abhyankar PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi) 408*2b4ed282SShri Abhyankar { 409*2b4ed282SShri Abhyankar PetscErrorCode ierr; 410*2b4ed282SShri Abhyankar 411*2b4ed282SShri Abhyankar PetscFunctionBegin; 412*2b4ed282SShri Abhyankar ierr = MatMultTranspose(H,phi,dpsi); 413*2b4ed282SShri Abhyankar 414*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 415*2b4ed282SShri Abhyankar } 416*2b4ed282SShri Abhyankar 417*2b4ed282SShri Abhyankar /* 418*2b4ed282SShri Abhyankar SNESVISetDescentDirection - Sets the descent direction for the semismooth algorithm 419*2b4ed282SShri Abhyankar 420*2b4ed282SShri Abhyankar Input Parameters: 421*2b4ed282SShri Abhyankar . snes - the SNES context. 422*2b4ed282SShri Abhyankar . dpsi - gradient of the merit function. 423*2b4ed282SShri Abhyankar 424*2b4ed282SShri Abhyankar Output Parameters: 425*2b4ed282SShri Abhyankar . flg - PETSC_TRUE if the sufficient descent condition is not satisfied. 426*2b4ed282SShri Abhyankar 427*2b4ed282SShri Abhyankar Notes: 428*2b4ed282SShri Abhyankar The condition for the sufficient descent direction is 429*2b4ed282SShri Abhyankar dpsi^T*Y > delta*||Y||^rho 430*2b4ed282SShri Abhyankar where rho, delta are as defined in the SNES_VI structure. 431*2b4ed282SShri Abhyankar If this condition is satisfied then the existing descent direction i.e. 432*2b4ed282SShri Abhyankar the direction given by the linear solve should be used otherwise it should be set to the negative of the merit function gradient i.e -dpsi. 433*2b4ed282SShri Abhyankar */ 434*2b4ed282SShri Abhyankar #undef __FUNCT__ 435*2b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckDescentDirection" 436*2b4ed282SShri Abhyankar PetscErrorCode SNESVICheckDescentDirection(SNES snes,Vec dpsi, Vec Y,PetscTruth* flg) 437*2b4ed282SShri Abhyankar { 438*2b4ed282SShri Abhyankar PetscErrorCode ierr; 439*2b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 440*2b4ed282SShri Abhyankar PetscScalar dpsidotY; 441*2b4ed282SShri Abhyankar PetscReal norm_Y,rhs; 442*2b4ed282SShri Abhyankar const PetscReal rho = vi->rho,delta=vi->delta; 443*2b4ed282SShri Abhyankar 444*2b4ed282SShri Abhyankar PetscFunctionBegin; 445*2b4ed282SShri Abhyankar 446*2b4ed282SShri Abhyankar *flg = PETSC_FALSE; 447*2b4ed282SShri Abhyankar ierr = VecDot(dpsi,Y,&dpsidotY);CHKERRQ(ierr); 448*2b4ed282SShri Abhyankar ierr = VecNormBegin(Y,NORM_2,&norm_Y);CHKERRQ(ierr); 449*2b4ed282SShri Abhyankar ierr = VecNormEnd(Y,NORM_2,&norm_Y);CHKERRQ(ierr); 450*2b4ed282SShri Abhyankar 451*2b4ed282SShri Abhyankar rhs = delta*PetscPowScalar(norm_Y,rho); 452*2b4ed282SShri Abhyankar 453*2b4ed282SShri Abhyankar if (dpsidotY <= rhs) *flg = PETSC_TRUE; 454*2b4ed282SShri Abhyankar 455*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 456*2b4ed282SShri Abhyankar } 457*2b4ed282SShri Abhyankar 458*2b4ed282SShri Abhyankar /* 459*2b4ed282SShri Abhyankar SNESVIAdjustInitialGuess - Readjusts the initial guess to the SNES solver supplied by the user so that the initial guess lies inside the feasible region . 460*2b4ed282SShri Abhyankar 461*2b4ed282SShri Abhyankar Input Parameters: 462*2b4ed282SShri Abhyankar . lb - lower bound. 463*2b4ed282SShri Abhyankar . ub - upper bound. 464*2b4ed282SShri Abhyankar 465*2b4ed282SShri Abhyankar Output Parameters: 466*2b4ed282SShri Abhyankar . X - the readjusted initial guess. 467*2b4ed282SShri Abhyankar 468*2b4ed282SShri Abhyankar Notes: 469*2b4ed282SShri Abhyankar The readjusted initial guess X[i] = max(max(min(l[i],X[i]),min(X[i],u[i])),min(l[i],u[i])) 470*2b4ed282SShri Abhyankar */ 471*2b4ed282SShri Abhyankar #undef __FUNCT__ 472*2b4ed282SShri Abhyankar #define __FUNCT__ "SNESVIAdjustInitialGuess" 473*2b4ed282SShri Abhyankar PetscErrorCode SNESVIAdjustInitialGuess(Vec X, Vec lb, Vec ub) 474*2b4ed282SShri Abhyankar { 475*2b4ed282SShri Abhyankar PetscErrorCode ierr; 476*2b4ed282SShri Abhyankar PetscInt i,nlocal; 477*2b4ed282SShri Abhyankar PetscScalar *x,*l,*u; 478*2b4ed282SShri Abhyankar 479*2b4ed282SShri Abhyankar PetscFunctionBegin; 480*2b4ed282SShri Abhyankar 481*2b4ed282SShri Abhyankar ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 482*2b4ed282SShri Abhyankar ierr = VecGetArray(X,&x);CHKERRQ(ierr); 483*2b4ed282SShri Abhyankar ierr = VecGetArray(lb,&l);CHKERRQ(ierr); 484*2b4ed282SShri Abhyankar ierr = VecGetArray(ub,&u);CHKERRQ(ierr); 485*2b4ed282SShri Abhyankar 486*2b4ed282SShri Abhyankar for(i = 0; i < nlocal; i++) { 487*2b4ed282SShri Abhyankar x[i] = PetscMax(PetscMax(PetscMin(x[i],l[i]),PetscMin(x[i],u[i])),PetscMin(l[i],u[i])); 488*2b4ed282SShri Abhyankar } 489*2b4ed282SShri Abhyankar 490*2b4ed282SShri Abhyankar ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 491*2b4ed282SShri Abhyankar ierr = VecRestoreArray(lb,&l);CHKERRQ(ierr); 492*2b4ed282SShri Abhyankar ierr = VecRestoreArray(ub,&u);CHKERRQ(ierr); 493*2b4ed282SShri Abhyankar 494*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 495*2b4ed282SShri Abhyankar } 496*2b4ed282SShri Abhyankar 497*2b4ed282SShri Abhyankar /* -------------------------------------------------------------------- 498*2b4ed282SShri Abhyankar 499*2b4ed282SShri Abhyankar This file implements a semismooth truncated Newton method with a line search, 500*2b4ed282SShri Abhyankar for solving a system of nonlinear equations in complementarity form, using the KSP, Vec, 501*2b4ed282SShri Abhyankar and Mat interfaces for linear solvers, vectors, and matrices, 502*2b4ed282SShri Abhyankar respectively. 503*2b4ed282SShri Abhyankar 504*2b4ed282SShri Abhyankar The following basic routines are required for each nonlinear solver: 505*2b4ed282SShri Abhyankar SNESCreate_XXX() - Creates a nonlinear solver context 506*2b4ed282SShri Abhyankar SNESSetFromOptions_XXX() - Sets runtime options 507*2b4ed282SShri Abhyankar SNESSolve_XXX() - Solves the nonlinear system 508*2b4ed282SShri Abhyankar SNESDestroy_XXX() - Destroys the nonlinear solver context 509*2b4ed282SShri Abhyankar The suffix "_XXX" denotes a particular implementation, in this case 510*2b4ed282SShri Abhyankar we use _VI (e.g., SNESCreate_VI, SNESSolve_VI) for solving 511*2b4ed282SShri Abhyankar systems of nonlinear equations with a line search (LS) method. 512*2b4ed282SShri Abhyankar These routines are actually called via the common user interface 513*2b4ed282SShri Abhyankar routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 514*2b4ed282SShri Abhyankar SNESDestroy(), so the application code interface remains identical 515*2b4ed282SShri Abhyankar for all nonlinear solvers. 516*2b4ed282SShri Abhyankar 517*2b4ed282SShri Abhyankar Another key routine is: 518*2b4ed282SShri Abhyankar SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 519*2b4ed282SShri Abhyankar by setting data structures and options. The interface routine SNESSetUp() 520*2b4ed282SShri Abhyankar is not usually called directly by the user, but instead is called by 521*2b4ed282SShri Abhyankar SNESSolve() if necessary. 522*2b4ed282SShri Abhyankar 523*2b4ed282SShri Abhyankar Additional basic routines are: 524*2b4ed282SShri Abhyankar SNESView_XXX() - Prints details of runtime options that 525*2b4ed282SShri Abhyankar have actually been used. 526*2b4ed282SShri Abhyankar These are called by application codes via the interface routines 527*2b4ed282SShri Abhyankar SNESView(). 528*2b4ed282SShri Abhyankar 529*2b4ed282SShri Abhyankar The various types of solvers (preconditioners, Krylov subspace methods, 530*2b4ed282SShri Abhyankar nonlinear solvers, timesteppers) are all organized similarly, so the 531*2b4ed282SShri Abhyankar above description applies to these categories also. 532*2b4ed282SShri Abhyankar 533*2b4ed282SShri Abhyankar -------------------------------------------------------------------- */ 534*2b4ed282SShri Abhyankar /* 535*2b4ed282SShri Abhyankar SNESSolve_VI - Solves the complementarity problem with a semismooth Newton 536*2b4ed282SShri Abhyankar method using a line search. 537*2b4ed282SShri Abhyankar 538*2b4ed282SShri Abhyankar Input Parameters: 539*2b4ed282SShri Abhyankar . snes - the SNES context 540*2b4ed282SShri Abhyankar 541*2b4ed282SShri Abhyankar Output Parameter: 542*2b4ed282SShri Abhyankar . outits - number of iterations until termination 543*2b4ed282SShri Abhyankar 544*2b4ed282SShri Abhyankar Application Interface Routine: SNESSolve() 545*2b4ed282SShri Abhyankar 546*2b4ed282SShri Abhyankar Notes: 547*2b4ed282SShri Abhyankar This implements essentially a semismooth Newton method with a 548*2b4ed282SShri Abhyankar line search. By default a cubic backtracking line search 549*2b4ed282SShri Abhyankar is employed, as described in the text "Numerical Methods for 550*2b4ed282SShri Abhyankar Unconstrained Optimization and Nonlinear Equations" by Dennis 551*2b4ed282SShri Abhyankar and Schnabel. 552*2b4ed282SShri Abhyankar */ 553*2b4ed282SShri Abhyankar #undef __FUNCT__ 554*2b4ed282SShri Abhyankar #define __FUNCT__ "SNESSolve_VI" 555*2b4ed282SShri Abhyankar PetscErrorCode SNESSolve_VI(SNES snes) 556*2b4ed282SShri Abhyankar { 557*2b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 558*2b4ed282SShri Abhyankar PetscErrorCode ierr; 559*2b4ed282SShri Abhyankar PetscInt maxits,i,lits; 560*2b4ed282SShri Abhyankar PetscTruth lssucceed,changedir; 561*2b4ed282SShri Abhyankar MatStructure flg = DIFFERENT_NONZERO_PATTERN; 562*2b4ed282SShri Abhyankar PetscReal gnorm,xnorm=0,ynorm; 563*2b4ed282SShri Abhyankar Vec Y,X,F,G,W; 564*2b4ed282SShri Abhyankar KSPConvergedReason kspreason; 565*2b4ed282SShri Abhyankar 566*2b4ed282SShri Abhyankar PetscFunctionBegin; 567*2b4ed282SShri Abhyankar snes->numFailures = 0; 568*2b4ed282SShri Abhyankar snes->numLinearSolveFailures = 0; 569*2b4ed282SShri Abhyankar snes->reason = SNES_CONVERGED_ITERATING; 570*2b4ed282SShri Abhyankar 571*2b4ed282SShri Abhyankar maxits = snes->max_its; /* maximum number of iterations */ 572*2b4ed282SShri Abhyankar X = snes->vec_sol; /* solution vector */ 573*2b4ed282SShri Abhyankar F = snes->vec_func; /* residual vector */ 574*2b4ed282SShri Abhyankar Y = snes->work[0]; /* work vectors */ 575*2b4ed282SShri Abhyankar G = snes->work[1]; 576*2b4ed282SShri Abhyankar W = snes->work[2]; 577*2b4ed282SShri Abhyankar 578*2b4ed282SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 579*2b4ed282SShri Abhyankar snes->iter = 0; 580*2b4ed282SShri Abhyankar snes->norm = 0.0; 581*2b4ed282SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 582*2b4ed282SShri Abhyankar 583*2b4ed282SShri Abhyankar ierr = SNESVIAdjustInitialGuess(X,vi->xl,vi->xu);CHKERRQ(ierr); 584*2b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr); 585*2b4ed282SShri Abhyankar if (snes->domainerror) { 586*2b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 587*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 588*2b4ed282SShri Abhyankar } 589*2b4ed282SShri Abhyankar /* Compute Merit function */ 590*2b4ed282SShri Abhyankar ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr); 591*2b4ed282SShri Abhyankar 592*2b4ed282SShri Abhyankar ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 593*2b4ed282SShri Abhyankar ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 594*2b4ed282SShri Abhyankar if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 595*2b4ed282SShri Abhyankar 596*2b4ed282SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 597*2b4ed282SShri Abhyankar snes->norm = vi->phinorm; 598*2b4ed282SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 599*2b4ed282SShri Abhyankar SNESLogConvHistory(snes,vi->phinorm,0); 600*2b4ed282SShri Abhyankar SNESMonitor(snes,0,vi->phinorm); 601*2b4ed282SShri Abhyankar 602*2b4ed282SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 603*2b4ed282SShri Abhyankar snes->ttol = vi->phinorm*snes->rtol; 604*2b4ed282SShri Abhyankar /* test convergence */ 605*2b4ed282SShri Abhyankar ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 606*2b4ed282SShri Abhyankar if (snes->reason) PetscFunctionReturn(0); 607*2b4ed282SShri Abhyankar 608*2b4ed282SShri Abhyankar for (i=0; i<maxits; i++) { 609*2b4ed282SShri Abhyankar 610*2b4ed282SShri Abhyankar /* Call general purpose update function */ 611*2b4ed282SShri Abhyankar if (snes->ops->update) { 612*2b4ed282SShri Abhyankar ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 613*2b4ed282SShri Abhyankar } 614*2b4ed282SShri Abhyankar 615*2b4ed282SShri Abhyankar /* Solve J Y = Phi, where J is the semismooth jacobian */ 616*2b4ed282SShri Abhyankar ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 617*2b4ed282SShri Abhyankar 618*2b4ed282SShri Abhyankar ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 619*2b4ed282SShri Abhyankar ierr = SNES_KSPSolve(snes,snes->ksp,vi->phi,Y);CHKERRQ(ierr); 620*2b4ed282SShri Abhyankar ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 621*2b4ed282SShri Abhyankar ierr = SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);CHKERRQ(ierr); 622*2b4ed282SShri Abhyankar ierr = SNESVICheckDescentDirection(snes,vi->dpsi,Y,&changedir);CHKERRQ(ierr); 623*2b4ed282SShri Abhyankar if (kspreason < 0 || changedir) { 624*2b4ed282SShri Abhyankar if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 625*2b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 626*2b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 627*2b4ed282SShri Abhyankar break; 628*2b4ed282SShri Abhyankar } 629*2b4ed282SShri Abhyankar ierr = VecCopy(vi->dpsi,Y);CHKERRQ(ierr); 630*2b4ed282SShri Abhyankar } 631*2b4ed282SShri Abhyankar ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 632*2b4ed282SShri Abhyankar snes->linear_its += lits; 633*2b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 634*2b4ed282SShri Abhyankar /* 635*2b4ed282SShri Abhyankar if (vi->precheckstep) { 636*2b4ed282SShri Abhyankar PetscTruth changed_y = PETSC_FALSE; 637*2b4ed282SShri Abhyankar ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr); 638*2b4ed282SShri Abhyankar } 639*2b4ed282SShri Abhyankar 640*2b4ed282SShri Abhyankar if (PetscLogPrintInfo){ 641*2b4ed282SShri Abhyankar ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 642*2b4ed282SShri Abhyankar } 643*2b4ed282SShri Abhyankar */ 644*2b4ed282SShri Abhyankar /* Compute a (scaled) negative update in the line search routine: 645*2b4ed282SShri Abhyankar Y <- X - lambda*Y 646*2b4ed282SShri Abhyankar and evaluate G = function(Y) (depends on the line search). 647*2b4ed282SShri Abhyankar */ 648*2b4ed282SShri Abhyankar ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 649*2b4ed282SShri Abhyankar ynorm = 1; gnorm = vi->phinorm; 650*2b4ed282SShri Abhyankar ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 651*2b4ed282SShri Abhyankar ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 652*2b4ed282SShri Abhyankar if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 653*2b4ed282SShri Abhyankar if (snes->domainerror) { 654*2b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 655*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 656*2b4ed282SShri Abhyankar } 657*2b4ed282SShri Abhyankar if (!lssucceed) { 658*2b4ed282SShri Abhyankar if (++snes->numFailures >= snes->maxFailures) { 659*2b4ed282SShri Abhyankar PetscTruth ismin; 660*2b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_LINE_SEARCH; 661*2b4ed282SShri Abhyankar ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 662*2b4ed282SShri Abhyankar if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 663*2b4ed282SShri Abhyankar break; 664*2b4ed282SShri Abhyankar } 665*2b4ed282SShri Abhyankar } 666*2b4ed282SShri Abhyankar /* Update function and solution vectors */ 667*2b4ed282SShri Abhyankar vi->phinorm = gnorm; 668*2b4ed282SShri Abhyankar vi->merit = 0.5*vi->phinorm*vi->phinorm; 669*2b4ed282SShri Abhyankar ierr = VecCopy(G,vi->phi);CHKERRQ(ierr); 670*2b4ed282SShri Abhyankar ierr = VecCopy(W,X);CHKERRQ(ierr); 671*2b4ed282SShri Abhyankar /* Monitor convergence */ 672*2b4ed282SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 673*2b4ed282SShri Abhyankar snes->iter = i+1; 674*2b4ed282SShri Abhyankar snes->norm = vi->phinorm; 675*2b4ed282SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 676*2b4ed282SShri Abhyankar SNESLogConvHistory(snes,snes->norm,lits); 677*2b4ed282SShri Abhyankar SNESMonitor(snes,snes->iter,snes->norm); 678*2b4ed282SShri Abhyankar /* Test for convergence, xnorm = || X || */ 679*2b4ed282SShri Abhyankar if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 680*2b4ed282SShri Abhyankar ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 681*2b4ed282SShri Abhyankar if (snes->reason) break; 682*2b4ed282SShri Abhyankar } 683*2b4ed282SShri Abhyankar if (i == maxits) { 684*2b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 685*2b4ed282SShri Abhyankar if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 686*2b4ed282SShri Abhyankar } 687*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 688*2b4ed282SShri Abhyankar } 689*2b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 690*2b4ed282SShri Abhyankar /* 691*2b4ed282SShri Abhyankar SNESSetUp_VI - Sets up the internal data structures for the later use 692*2b4ed282SShri Abhyankar of the SNESVI nonlinear solver. 693*2b4ed282SShri Abhyankar 694*2b4ed282SShri Abhyankar Input Parameter: 695*2b4ed282SShri Abhyankar . snes - the SNES context 696*2b4ed282SShri Abhyankar . x - the solution vector 697*2b4ed282SShri Abhyankar 698*2b4ed282SShri Abhyankar Application Interface Routine: SNESSetUp() 699*2b4ed282SShri Abhyankar 700*2b4ed282SShri Abhyankar Notes: 701*2b4ed282SShri Abhyankar For basic use of the SNES solvers, the user need not explicitly call 702*2b4ed282SShri Abhyankar SNESSetUp(), since these actions will automatically occur during 703*2b4ed282SShri Abhyankar the call to SNESSolve(). 704*2b4ed282SShri Abhyankar */ 705*2b4ed282SShri Abhyankar #undef __FUNCT__ 706*2b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetUp_VI" 707*2b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes) 708*2b4ed282SShri Abhyankar { 709*2b4ed282SShri Abhyankar PetscErrorCode ierr; 710*2b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*) snes->data; 711*2b4ed282SShri Abhyankar PetscInt i_start[3],i_end[3]; 712*2b4ed282SShri Abhyankar 713*2b4ed282SShri Abhyankar PetscFunctionBegin; 714*2b4ed282SShri Abhyankar if (!snes->vec_sol_update) { 715*2b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr); 716*2b4ed282SShri Abhyankar ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr); 717*2b4ed282SShri Abhyankar } 718*2b4ed282SShri Abhyankar if (!snes->work) { 719*2b4ed282SShri Abhyankar snes->nwork = 3; 720*2b4ed282SShri Abhyankar ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 721*2b4ed282SShri Abhyankar ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr); 722*2b4ed282SShri Abhyankar } 723*2b4ed282SShri Abhyankar 724*2b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->phi); CHKERRQ(ierr); 725*2b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->dpsi); CHKERRQ(ierr); 726*2b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->Da); CHKERRQ(ierr); 727*2b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->Db); CHKERRQ(ierr); 728*2b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->z);CHKERRQ(ierr); 729*2b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->t); CHKERRQ(ierr); 730*2b4ed282SShri Abhyankar 731*2b4ed282SShri Abhyankar /* If the lower and upper bound on variables are not set, set it to 732*2b4ed282SShri Abhyankar -Inf and Inf */ 733*2b4ed282SShri Abhyankar if (!vi->xl && !vi->xu) { 734*2b4ed282SShri Abhyankar vi->usersetxbounds = PETSC_FALSE; 735*2b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->xl); CHKERRQ(ierr); 736*2b4ed282SShri Abhyankar ierr = VecSet(vi->xl,PETSC_VI_NINF);CHKERRQ(ierr); 737*2b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->xu); CHKERRQ(ierr); 738*2b4ed282SShri Abhyankar ierr = VecSet(vi->xu,PETSC_VI_INF);CHKERRQ(ierr); 739*2b4ed282SShri Abhyankar } else { 740*2b4ed282SShri Abhyankar /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */ 741*2b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr); 742*2b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(vi->xl,i_start+1,i_end+1);CHKERRQ(ierr); 743*2b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(vi->xu,i_start+2,i_end+2);CHKERRQ(ierr); 744*2b4ed282SShri Abhyankar if ((i_start[0] != i_start[1]) || (i_start[0] != i_start[2]) || (i_end[0] != i_end[1]) || (i_end[0] != i_end[2])) 745*2b4ed282SShri Abhyankar SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Distribution of lower bound, upper bound and the solution vector should be identical across all the processors."); 746*2b4ed282SShri Abhyankar } 747*2b4ed282SShri Abhyankar 748*2b4ed282SShri Abhyankar vi->computeuserfunction = snes->ops->computefunction; 749*2b4ed282SShri Abhyankar vi->computeuserjacobian = snes->ops->computejacobian; 750*2b4ed282SShri Abhyankar 751*2b4ed282SShri Abhyankar snes->ops->computefunction = SNESVIComputeSSFunction; 752*2b4ed282SShri Abhyankar snes->ops->computejacobian = SNESVIComputeSSJacobian; 753*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 754*2b4ed282SShri Abhyankar } 755*2b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 756*2b4ed282SShri Abhyankar /* 757*2b4ed282SShri Abhyankar SNESDestroy_VI - Destroys the private SNES_VI context that was created 758*2b4ed282SShri Abhyankar with SNESCreate_VI(). 759*2b4ed282SShri Abhyankar 760*2b4ed282SShri Abhyankar Input Parameter: 761*2b4ed282SShri Abhyankar . snes - the SNES context 762*2b4ed282SShri Abhyankar 763*2b4ed282SShri Abhyankar Application Interface Routine: SNESDestroy() 764*2b4ed282SShri Abhyankar */ 765*2b4ed282SShri Abhyankar #undef __FUNCT__ 766*2b4ed282SShri Abhyankar #define __FUNCT__ "SNESDestroy_VI" 767*2b4ed282SShri Abhyankar PetscErrorCode SNESDestroy_VI(SNES snes) 768*2b4ed282SShri Abhyankar { 769*2b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*) snes->data; 770*2b4ed282SShri Abhyankar PetscErrorCode ierr; 771*2b4ed282SShri Abhyankar 772*2b4ed282SShri Abhyankar PetscFunctionBegin; 773*2b4ed282SShri Abhyankar if (snes->vec_sol_update) { 774*2b4ed282SShri Abhyankar ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr); 775*2b4ed282SShri Abhyankar snes->vec_sol_update = PETSC_NULL; 776*2b4ed282SShri Abhyankar } 777*2b4ed282SShri Abhyankar if (snes->nwork) { 778*2b4ed282SShri Abhyankar ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 779*2b4ed282SShri Abhyankar snes->nwork = 0; 780*2b4ed282SShri Abhyankar snes->work = PETSC_NULL; 781*2b4ed282SShri Abhyankar } 782*2b4ed282SShri Abhyankar 783*2b4ed282SShri Abhyankar /* clear vectors */ 784*2b4ed282SShri Abhyankar ierr = VecDestroy(vi->phi); CHKERRQ(ierr); 785*2b4ed282SShri Abhyankar ierr = VecDestroy(vi->dpsi); CHKERRQ(ierr); 786*2b4ed282SShri Abhyankar ierr = VecDestroy(vi->Da); CHKERRQ(ierr); 787*2b4ed282SShri Abhyankar ierr = VecDestroy(vi->Db); CHKERRQ(ierr); 788*2b4ed282SShri Abhyankar ierr = VecDestroy(vi->z); CHKERRQ(ierr); 789*2b4ed282SShri Abhyankar ierr = VecDestroy(vi->t); CHKERRQ(ierr); 790*2b4ed282SShri Abhyankar if (!vi->usersetxbounds) { 791*2b4ed282SShri Abhyankar ierr = VecDestroy(vi->xl); CHKERRQ(ierr); 792*2b4ed282SShri Abhyankar ierr = VecDestroy(vi->xu); CHKERRQ(ierr); 793*2b4ed282SShri Abhyankar } 794*2b4ed282SShri Abhyankar if (vi->monitor) { 795*2b4ed282SShri Abhyankar ierr = PetscViewerASCIIMonitorDestroy(vi->monitor);CHKERRQ(ierr); 796*2b4ed282SShri Abhyankar } 797*2b4ed282SShri Abhyankar ierr = PetscFree(snes->data);CHKERRQ(ierr); 798*2b4ed282SShri Abhyankar 799*2b4ed282SShri Abhyankar /* clear composed functions */ 800*2b4ed282SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr); 801*2b4ed282SShri Abhyankar 802*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 803*2b4ed282SShri Abhyankar } 804*2b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 805*2b4ed282SShri Abhyankar #undef __FUNCT__ 806*2b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNo_VI" 807*2b4ed282SShri Abhyankar 808*2b4ed282SShri Abhyankar /* 809*2b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchNo routine in snes/impls/ls/ls.c 810*2b4ed282SShri Abhyankar 811*2b4ed282SShri Abhyankar */ 812*2b4ed282SShri Abhyankar PetscErrorCode SNESLineSearchNo_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 813*2b4ed282SShri Abhyankar { 814*2b4ed282SShri Abhyankar PetscErrorCode ierr; 815*2b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 816*2b4ed282SShri Abhyankar PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 817*2b4ed282SShri Abhyankar 818*2b4ed282SShri Abhyankar PetscFunctionBegin; 819*2b4ed282SShri Abhyankar *flag = PETSC_TRUE; 820*2b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 821*2b4ed282SShri Abhyankar ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 822*2b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 823*2b4ed282SShri Abhyankar if (vi->postcheckstep) { 824*2b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 825*2b4ed282SShri Abhyankar } 826*2b4ed282SShri Abhyankar if (changed_y) { 827*2b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 828*2b4ed282SShri Abhyankar } 829*2b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 830*2b4ed282SShri Abhyankar if (!snes->domainerror) { 831*2b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 832*2b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 833*2b4ed282SShri Abhyankar } 834*2b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 835*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 836*2b4ed282SShri Abhyankar } 837*2b4ed282SShri Abhyankar 838*2b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 839*2b4ed282SShri Abhyankar #undef __FUNCT__ 840*2b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNoNorms_VI" 841*2b4ed282SShri Abhyankar 842*2b4ed282SShri Abhyankar /* 843*2b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c 844*2b4ed282SShri Abhyankar */ 845*2b4ed282SShri Abhyankar PetscErrorCode SNESLineSearchNoNorms_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 846*2b4ed282SShri Abhyankar { 847*2b4ed282SShri Abhyankar PetscErrorCode ierr; 848*2b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 849*2b4ed282SShri Abhyankar PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 850*2b4ed282SShri Abhyankar 851*2b4ed282SShri Abhyankar PetscFunctionBegin; 852*2b4ed282SShri Abhyankar *flag = PETSC_TRUE; 853*2b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 854*2b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 855*2b4ed282SShri Abhyankar if (vi->postcheckstep) { 856*2b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 857*2b4ed282SShri Abhyankar } 858*2b4ed282SShri Abhyankar if (changed_y) { 859*2b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 860*2b4ed282SShri Abhyankar } 861*2b4ed282SShri Abhyankar 862*2b4ed282SShri Abhyankar /* don't evaluate function the last time through */ 863*2b4ed282SShri Abhyankar if (snes->iter < snes->max_its-1) { 864*2b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 865*2b4ed282SShri Abhyankar } 866*2b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 867*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 868*2b4ed282SShri Abhyankar } 869*2b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 870*2b4ed282SShri Abhyankar #undef __FUNCT__ 871*2b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchCubic_VI" 872*2b4ed282SShri Abhyankar /* 873*2b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchCubic in snes/impls/ls/ls.c 874*2b4ed282SShri Abhyankar */ 875*2b4ed282SShri Abhyankar PetscErrorCode SNESLineSearchCubic_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 876*2b4ed282SShri Abhyankar { 877*2b4ed282SShri Abhyankar /* 878*2b4ed282SShri Abhyankar Note that for line search purposes we work with with the related 879*2b4ed282SShri Abhyankar minimization problem: 880*2b4ed282SShri Abhyankar min z(x): R^n -> R, 881*2b4ed282SShri Abhyankar where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 882*2b4ed282SShri Abhyankar This function z(x) is same as the merit function for the complementarity problem. 883*2b4ed282SShri Abhyankar */ 884*2b4ed282SShri Abhyankar 885*2b4ed282SShri Abhyankar PetscReal initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 886*2b4ed282SShri Abhyankar PetscReal minlambda,lambda,lambdatemp; 887*2b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 888*2b4ed282SShri Abhyankar PetscScalar cinitslope; 889*2b4ed282SShri Abhyankar #endif 890*2b4ed282SShri Abhyankar PetscErrorCode ierr; 891*2b4ed282SShri Abhyankar PetscInt count; 892*2b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 893*2b4ed282SShri Abhyankar PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 894*2b4ed282SShri Abhyankar MPI_Comm comm; 895*2b4ed282SShri Abhyankar 896*2b4ed282SShri Abhyankar PetscFunctionBegin; 897*2b4ed282SShri Abhyankar ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 898*2b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 899*2b4ed282SShri Abhyankar *flag = PETSC_TRUE; 900*2b4ed282SShri Abhyankar 901*2b4ed282SShri Abhyankar ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 902*2b4ed282SShri Abhyankar if (*ynorm == 0.0) { 903*2b4ed282SShri Abhyankar if (vi->monitor) { 904*2b4ed282SShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->monitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 905*2b4ed282SShri Abhyankar } 906*2b4ed282SShri Abhyankar *gnorm = fnorm; 907*2b4ed282SShri Abhyankar ierr = VecCopy(x,w);CHKERRQ(ierr); 908*2b4ed282SShri Abhyankar ierr = VecCopy(f,g);CHKERRQ(ierr); 909*2b4ed282SShri Abhyankar *flag = PETSC_FALSE; 910*2b4ed282SShri Abhyankar goto theend1; 911*2b4ed282SShri Abhyankar } 912*2b4ed282SShri Abhyankar if (*ynorm > vi->maxstep) { /* Step too big, so scale back */ 913*2b4ed282SShri Abhyankar if (vi->monitor) { 914*2b4ed282SShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->monitor," Line search: Scaling step by %G old ynorm %G\n",vi->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr); 915*2b4ed282SShri Abhyankar } 916*2b4ed282SShri Abhyankar ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr); 917*2b4ed282SShri Abhyankar *ynorm = vi->maxstep; 918*2b4ed282SShri Abhyankar } 919*2b4ed282SShri Abhyankar ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 920*2b4ed282SShri Abhyankar minlambda = vi->minlambda/rellength; 921*2b4ed282SShri Abhyankar /* ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); */ 922*2b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 923*2b4ed282SShri Abhyankar ierr = VecDot(vi->dpsi,y,&cinitslope);CHKERRQ(ierr); 924*2b4ed282SShri Abhyankar initslope = PetscRealPart(cinitslope); 925*2b4ed282SShri Abhyankar #else 926*2b4ed282SShri Abhyankar ierr = VecDot(vi->dpsi,y,&initslope);CHKERRQ(ierr); 927*2b4ed282SShri Abhyankar #endif 928*2b4ed282SShri Abhyankar if (initslope > 0.0) initslope = -initslope; 929*2b4ed282SShri Abhyankar if (initslope == 0.0) initslope = -1.0; 930*2b4ed282SShri Abhyankar 931*2b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 932*2b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 933*2b4ed282SShri Abhyankar ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 934*2b4ed282SShri Abhyankar *flag = PETSC_FALSE; 935*2b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 936*2b4ed282SShri Abhyankar goto theend1; 937*2b4ed282SShri Abhyankar } 938*2b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 939*2b4ed282SShri Abhyankar if (snes->domainerror) { 940*2b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 941*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 942*2b4ed282SShri Abhyankar } 943*2b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 944*2b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 945*2b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 946*2b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */ 947*2b4ed282SShri Abhyankar if (vi->monitor) { 948*2b4ed282SShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->monitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 949*2b4ed282SShri Abhyankar } 950*2b4ed282SShri Abhyankar goto theend1; 951*2b4ed282SShri Abhyankar } 952*2b4ed282SShri Abhyankar 953*2b4ed282SShri Abhyankar /* Fit points with quadratic */ 954*2b4ed282SShri Abhyankar lambda = 1.0; 955*2b4ed282SShri Abhyankar lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 956*2b4ed282SShri Abhyankar lambdaprev = lambda; 957*2b4ed282SShri Abhyankar gnormprev = *gnorm; 958*2b4ed282SShri Abhyankar if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 959*2b4ed282SShri Abhyankar if (lambdatemp <= .1*lambda) lambda = .1*lambda; 960*2b4ed282SShri Abhyankar else lambda = lambdatemp; 961*2b4ed282SShri Abhyankar 962*2b4ed282SShri Abhyankar ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 963*2b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 964*2b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 965*2b4ed282SShri Abhyankar *flag = PETSC_FALSE; 966*2b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 967*2b4ed282SShri Abhyankar goto theend1; 968*2b4ed282SShri Abhyankar } 969*2b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 970*2b4ed282SShri Abhyankar if (snes->domainerror) { 971*2b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 972*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 973*2b4ed282SShri Abhyankar } 974*2b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 975*2b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 976*2b4ed282SShri Abhyankar if (vi->monitor) { 977*2b4ed282SShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->monitor," Line search: gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr); 978*2b4ed282SShri Abhyankar } 979*2b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */ 980*2b4ed282SShri Abhyankar if (vi->monitor) { 981*2b4ed282SShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->monitor," Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 982*2b4ed282SShri Abhyankar } 983*2b4ed282SShri Abhyankar goto theend1; 984*2b4ed282SShri Abhyankar } 985*2b4ed282SShri Abhyankar 986*2b4ed282SShri Abhyankar /* Fit points with cubic */ 987*2b4ed282SShri Abhyankar count = 1; 988*2b4ed282SShri Abhyankar while (PETSC_TRUE) { 989*2b4ed282SShri Abhyankar if (lambda <= minlambda) { 990*2b4ed282SShri Abhyankar if (vi->monitor) { 991*2b4ed282SShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->monitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 992*2b4ed282SShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->monitor," Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,minlambda,lambda,initslope);CHKERRQ(ierr); 993*2b4ed282SShri Abhyankar } 994*2b4ed282SShri Abhyankar *flag = PETSC_FALSE; 995*2b4ed282SShri Abhyankar break; 996*2b4ed282SShri Abhyankar } 997*2b4ed282SShri Abhyankar t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 998*2b4ed282SShri Abhyankar t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 999*2b4ed282SShri Abhyankar a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 1000*2b4ed282SShri Abhyankar b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 1001*2b4ed282SShri Abhyankar d = b*b - 3*a*initslope; 1002*2b4ed282SShri Abhyankar if (d < 0.0) d = 0.0; 1003*2b4ed282SShri Abhyankar if (a == 0.0) { 1004*2b4ed282SShri Abhyankar lambdatemp = -initslope/(2.0*b); 1005*2b4ed282SShri Abhyankar } else { 1006*2b4ed282SShri Abhyankar lambdatemp = (-b + sqrt(d))/(3.0*a); 1007*2b4ed282SShri Abhyankar } 1008*2b4ed282SShri Abhyankar lambdaprev = lambda; 1009*2b4ed282SShri Abhyankar gnormprev = *gnorm; 1010*2b4ed282SShri Abhyankar if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 1011*2b4ed282SShri Abhyankar if (lambdatemp <= .1*lambda) lambda = .1*lambda; 1012*2b4ed282SShri Abhyankar else lambda = lambdatemp; 1013*2b4ed282SShri Abhyankar ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 1014*2b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 1015*2b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 1016*2b4ed282SShri Abhyankar ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 1017*2b4ed282SShri Abhyankar *flag = PETSC_FALSE; 1018*2b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 1019*2b4ed282SShri Abhyankar break; 1020*2b4ed282SShri Abhyankar } 1021*2b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1022*2b4ed282SShri Abhyankar if (snes->domainerror) { 1023*2b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1024*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 1025*2b4ed282SShri Abhyankar } 1026*2b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 1027*2b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1028*2b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* is reduction enough? */ 1029*2b4ed282SShri Abhyankar if (vi->monitor) { 1030*2b4ed282SShri Abhyankar ierr = PetscPrintf(comm," Line search: Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 1031*2b4ed282SShri Abhyankar } 1032*2b4ed282SShri Abhyankar break; 1033*2b4ed282SShri Abhyankar } else { 1034*2b4ed282SShri Abhyankar if (vi->monitor) { 1035*2b4ed282SShri Abhyankar ierr = PetscPrintf(comm," Line search: Cubic step no good, shrinking lambda, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 1036*2b4ed282SShri Abhyankar } 1037*2b4ed282SShri Abhyankar } 1038*2b4ed282SShri Abhyankar count++; 1039*2b4ed282SShri Abhyankar } 1040*2b4ed282SShri Abhyankar theend1: 1041*2b4ed282SShri Abhyankar /* Optional user-defined check for line search step validity */ 1042*2b4ed282SShri Abhyankar if (vi->postcheckstep && *flag) { 1043*2b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 1044*2b4ed282SShri Abhyankar if (changed_y) { 1045*2b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 1046*2b4ed282SShri Abhyankar } 1047*2b4ed282SShri Abhyankar if (changed_y || changed_w) { /* recompute the function if the step has changed */ 1048*2b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1049*2b4ed282SShri Abhyankar if (snes->domainerror) { 1050*2b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1051*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 1052*2b4ed282SShri Abhyankar } 1053*2b4ed282SShri Abhyankar ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 1054*2b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1055*2b4ed282SShri Abhyankar ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 1056*2b4ed282SShri Abhyankar ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 1057*2b4ed282SShri Abhyankar ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 1058*2b4ed282SShri Abhyankar } 1059*2b4ed282SShri Abhyankar } 1060*2b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1061*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 1062*2b4ed282SShri Abhyankar } 1063*2b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 1064*2b4ed282SShri Abhyankar #undef __FUNCT__ 1065*2b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchQuadratic_VI" 1066*2b4ed282SShri Abhyankar /* 1067*2b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchQuadratic in snes/impls/ls/ls.c 1068*2b4ed282SShri Abhyankar */ 1069*2b4ed282SShri Abhyankar PetscErrorCode SNESLineSearchQuadratic_VI(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 1070*2b4ed282SShri Abhyankar { 1071*2b4ed282SShri Abhyankar /* 1072*2b4ed282SShri Abhyankar Note that for line search purposes we work with with the related 1073*2b4ed282SShri Abhyankar minimization problem: 1074*2b4ed282SShri Abhyankar min z(x): R^n -> R, 1075*2b4ed282SShri Abhyankar where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 1076*2b4ed282SShri Abhyankar z(x) is the same as the merit function for the complementarity problem 1077*2b4ed282SShri Abhyankar */ 1078*2b4ed282SShri Abhyankar PetscReal initslope,minlambda,lambda,lambdatemp,rellength; 1079*2b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 1080*2b4ed282SShri Abhyankar PetscScalar cinitslope; 1081*2b4ed282SShri Abhyankar #endif 1082*2b4ed282SShri Abhyankar PetscErrorCode ierr; 1083*2b4ed282SShri Abhyankar PetscInt count; 1084*2b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1085*2b4ed282SShri Abhyankar PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 1086*2b4ed282SShri Abhyankar 1087*2b4ed282SShri Abhyankar PetscFunctionBegin; 1088*2b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1089*2b4ed282SShri Abhyankar *flag = PETSC_TRUE; 1090*2b4ed282SShri Abhyankar 1091*2b4ed282SShri Abhyankar ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 1092*2b4ed282SShri Abhyankar if (*ynorm == 0.0) { 1093*2b4ed282SShri Abhyankar if (vi->monitor) { 1094*2b4ed282SShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->monitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr); 1095*2b4ed282SShri Abhyankar } 1096*2b4ed282SShri Abhyankar *gnorm = fnorm; 1097*2b4ed282SShri Abhyankar ierr = VecCopy(x,w);CHKERRQ(ierr); 1098*2b4ed282SShri Abhyankar ierr = VecCopy(f,g);CHKERRQ(ierr); 1099*2b4ed282SShri Abhyankar *flag = PETSC_FALSE; 1100*2b4ed282SShri Abhyankar goto theend2; 1101*2b4ed282SShri Abhyankar } 1102*2b4ed282SShri Abhyankar if (*ynorm > vi->maxstep) { /* Step too big, so scale back */ 1103*2b4ed282SShri Abhyankar ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr); 1104*2b4ed282SShri Abhyankar *ynorm = vi->maxstep; 1105*2b4ed282SShri Abhyankar } 1106*2b4ed282SShri Abhyankar ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 1107*2b4ed282SShri Abhyankar minlambda = vi->minlambda/rellength; 1108*2b4ed282SShri Abhyankar /* ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); */ 1109*2b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 1110*2b4ed282SShri Abhyankar ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 1111*2b4ed282SShri Abhyankar initslope = PetscRealPart(cinitslope); 1112*2b4ed282SShri Abhyankar #else 1113*2b4ed282SShri Abhyankar ierr = VecDot(vi->dpsi,y,&initslope);CHKERRQ(ierr); 1114*2b4ed282SShri Abhyankar #endif 1115*2b4ed282SShri Abhyankar if (initslope > 0.0) initslope = -initslope; 1116*2b4ed282SShri Abhyankar if (initslope == 0.0) initslope = -1.0; 1117*2b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr); 1118*2b4ed282SShri Abhyankar 1119*2b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 1120*2b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 1121*2b4ed282SShri Abhyankar ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 1122*2b4ed282SShri Abhyankar *flag = PETSC_FALSE; 1123*2b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 1124*2b4ed282SShri Abhyankar goto theend2; 1125*2b4ed282SShri Abhyankar } 1126*2b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1127*2b4ed282SShri Abhyankar if (snes->domainerror) { 1128*2b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1129*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 1130*2b4ed282SShri Abhyankar } 1131*2b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 1132*2b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1133*2b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */ 1134*2b4ed282SShri Abhyankar if (vi->monitor) { 1135*2b4ed282SShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->monitor," Line Search: Using full step\n");CHKERRQ(ierr); 1136*2b4ed282SShri Abhyankar } 1137*2b4ed282SShri Abhyankar goto theend2; 1138*2b4ed282SShri Abhyankar } 1139*2b4ed282SShri Abhyankar 1140*2b4ed282SShri Abhyankar /* Fit points with quadratic */ 1141*2b4ed282SShri Abhyankar lambda = 1.0; 1142*2b4ed282SShri Abhyankar count = 1; 1143*2b4ed282SShri Abhyankar while (PETSC_TRUE) { 1144*2b4ed282SShri Abhyankar if (lambda <= minlambda) { /* bad luck; use full step */ 1145*2b4ed282SShri Abhyankar if (vi->monitor) { 1146*2b4ed282SShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->monitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr); 1147*2b4ed282SShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->monitor,"Line search: fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 1148*2b4ed282SShri Abhyankar } 1149*2b4ed282SShri Abhyankar ierr = VecCopy(x,w);CHKERRQ(ierr); 1150*2b4ed282SShri Abhyankar *flag = PETSC_FALSE; 1151*2b4ed282SShri Abhyankar break; 1152*2b4ed282SShri Abhyankar } 1153*2b4ed282SShri Abhyankar lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 1154*2b4ed282SShri Abhyankar if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 1155*2b4ed282SShri Abhyankar if (lambdatemp <= .1*lambda) lambda = .1*lambda; 1156*2b4ed282SShri Abhyankar else lambda = lambdatemp; 1157*2b4ed282SShri Abhyankar 1158*2b4ed282SShri Abhyankar ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 1159*2b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 1160*2b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 1161*2b4ed282SShri Abhyankar ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 1162*2b4ed282SShri Abhyankar *flag = PETSC_FALSE; 1163*2b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 1164*2b4ed282SShri Abhyankar break; 1165*2b4ed282SShri Abhyankar } 1166*2b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 1167*2b4ed282SShri Abhyankar if (snes->domainerror) { 1168*2b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1169*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 1170*2b4ed282SShri Abhyankar } 1171*2b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 1172*2b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1173*2b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */ 1174*2b4ed282SShri Abhyankar if (vi->monitor) { 1175*2b4ed282SShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->monitor," Line Search: Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr); 1176*2b4ed282SShri Abhyankar } 1177*2b4ed282SShri Abhyankar break; 1178*2b4ed282SShri Abhyankar } 1179*2b4ed282SShri Abhyankar count++; 1180*2b4ed282SShri Abhyankar } 1181*2b4ed282SShri Abhyankar theend2: 1182*2b4ed282SShri Abhyankar /* Optional user-defined check for line search step validity */ 1183*2b4ed282SShri Abhyankar if (vi->postcheckstep) { 1184*2b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 1185*2b4ed282SShri Abhyankar if (changed_y) { 1186*2b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 1187*2b4ed282SShri Abhyankar } 1188*2b4ed282SShri Abhyankar if (changed_y || changed_w) { /* recompute the function if the step has changed */ 1189*2b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g); 1190*2b4ed282SShri Abhyankar if (snes->domainerror) { 1191*2b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1192*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 1193*2b4ed282SShri Abhyankar } 1194*2b4ed282SShri Abhyankar ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 1195*2b4ed282SShri Abhyankar ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 1196*2b4ed282SShri Abhyankar ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 1197*2b4ed282SShri Abhyankar ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 1198*2b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1199*2b4ed282SShri Abhyankar } 1200*2b4ed282SShri Abhyankar } 1201*2b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 1202*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 1203*2b4ed282SShri Abhyankar } 1204*2b4ed282SShri Abhyankar 1205*2b4ed282SShri Abhyankar typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscTruth*); /* force argument to next function to not be extern C*/ 1206*2b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 1207*2b4ed282SShri Abhyankar EXTERN_C_BEGIN 1208*2b4ed282SShri Abhyankar #undef __FUNCT__ 1209*2b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchSet_VI" 1210*2b4ed282SShri Abhyankar PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_VI(SNES snes,FCN2 func,void *lsctx) 1211*2b4ed282SShri Abhyankar { 1212*2b4ed282SShri Abhyankar PetscFunctionBegin; 1213*2b4ed282SShri Abhyankar ((SNES_VI *)(snes->data))->LineSearch = func; 1214*2b4ed282SShri Abhyankar ((SNES_VI *)(snes->data))->lsP = lsctx; 1215*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 1216*2b4ed282SShri Abhyankar } 1217*2b4ed282SShri Abhyankar EXTERN_C_END 1218*2b4ed282SShri Abhyankar 1219*2b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 1220*2b4ed282SShri Abhyankar EXTERN_C_BEGIN 1221*2b4ed282SShri Abhyankar #undef __FUNCT__ 1222*2b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchSetMonitor_VI" 1223*2b4ed282SShri Abhyankar PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetMonitor_VI(SNES snes,PetscTruth flg) 1224*2b4ed282SShri Abhyankar { 1225*2b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1226*2b4ed282SShri Abhyankar PetscErrorCode ierr; 1227*2b4ed282SShri Abhyankar 1228*2b4ed282SShri Abhyankar PetscFunctionBegin; 1229*2b4ed282SShri Abhyankar if (flg && !vi->monitor) { 1230*2b4ed282SShri Abhyankar ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",((PetscObject)snes)->tablevel,&vi->monitor);CHKERRQ(ierr); 1231*2b4ed282SShri Abhyankar } else if (!flg && vi->monitor) { 1232*2b4ed282SShri Abhyankar ierr = PetscViewerASCIIMonitorDestroy(vi->monitor);CHKERRQ(ierr); 1233*2b4ed282SShri Abhyankar } 1234*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 1235*2b4ed282SShri Abhyankar } 1236*2b4ed282SShri Abhyankar EXTERN_C_END 1237*2b4ed282SShri Abhyankar 1238*2b4ed282SShri Abhyankar /* 1239*2b4ed282SShri Abhyankar SNESView_VI - Prints info from the SNESVI data structure. 1240*2b4ed282SShri Abhyankar 1241*2b4ed282SShri Abhyankar Input Parameters: 1242*2b4ed282SShri Abhyankar . SNES - the SNES context 1243*2b4ed282SShri Abhyankar . viewer - visualization context 1244*2b4ed282SShri Abhyankar 1245*2b4ed282SShri Abhyankar Application Interface Routine: SNESView() 1246*2b4ed282SShri Abhyankar */ 1247*2b4ed282SShri Abhyankar #undef __FUNCT__ 1248*2b4ed282SShri Abhyankar #define __FUNCT__ "SNESView_VI" 1249*2b4ed282SShri Abhyankar static PetscErrorCode SNESView_VI(SNES snes,PetscViewer viewer) 1250*2b4ed282SShri Abhyankar { 1251*2b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI *)snes->data; 1252*2b4ed282SShri Abhyankar const char *cstr; 1253*2b4ed282SShri Abhyankar PetscErrorCode ierr; 1254*2b4ed282SShri Abhyankar PetscTruth iascii; 1255*2b4ed282SShri Abhyankar 1256*2b4ed282SShri Abhyankar PetscFunctionBegin; 1257*2b4ed282SShri Abhyankar ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1258*2b4ed282SShri Abhyankar if (iascii) { 1259*2b4ed282SShri Abhyankar if (vi->LineSearch == SNESLineSearchNo_VI) cstr = "SNESLineSearchNo"; 1260*2b4ed282SShri Abhyankar else if (vi->LineSearch == SNESLineSearchQuadratic_VI) cstr = "SNESLineSearchQuadratic"; 1261*2b4ed282SShri Abhyankar else if (vi->LineSearch == SNESLineSearchCubic_VI) cstr = "SNESLineSearchCubic"; 1262*2b4ed282SShri Abhyankar else cstr = "unknown"; 1263*2b4ed282SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 1264*2b4ed282SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," alpha=%G, maxstep=%G, minlambda=%G\n",vi->alpha,vi->maxstep,vi->minlambda);CHKERRQ(ierr); 1265*2b4ed282SShri Abhyankar } else { 1266*2b4ed282SShri Abhyankar SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ VI",((PetscObject)viewer)->type_name); 1267*2b4ed282SShri Abhyankar } 1268*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 1269*2b4ed282SShri Abhyankar } 1270*2b4ed282SShri Abhyankar 1271*2b4ed282SShri Abhyankar /* 1272*2b4ed282SShri Abhyankar SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. 1273*2b4ed282SShri Abhyankar 1274*2b4ed282SShri Abhyankar Input Parameters: 1275*2b4ed282SShri Abhyankar . snes - the SNES context. 1276*2b4ed282SShri Abhyankar . xl - lower bound. 1277*2b4ed282SShri Abhyankar . xu - upper bound. 1278*2b4ed282SShri Abhyankar 1279*2b4ed282SShri Abhyankar Notes: 1280*2b4ed282SShri Abhyankar If this routine is not called then the lower and upper bounds are set to 1281*2b4ed282SShri Abhyankar -Infinity and Infinity respectively during SNESSetUp. 1282*2b4ed282SShri Abhyankar */ 1283*2b4ed282SShri Abhyankar 1284*2b4ed282SShri Abhyankar #undef __FUNCT__ 1285*2b4ed282SShri Abhyankar #define __FUNCT__ "SNESVISetVariableBounds" 1286*2b4ed282SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu) 1287*2b4ed282SShri Abhyankar { 1288*2b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1289*2b4ed282SShri Abhyankar 1290*2b4ed282SShri Abhyankar PetscFunctionBegin; 1291*2b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1292*2b4ed282SShri Abhyankar PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 1293*2b4ed282SShri Abhyankar PetscValidHeaderSpecific(xu,VEC_CLASSID,3); 1294*2b4ed282SShri Abhyankar 1295*2b4ed282SShri Abhyankar /* Check if SNESSetFunction is called */ 1296*2b4ed282SShri Abhyankar if(!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1297*2b4ed282SShri Abhyankar 1298*2b4ed282SShri Abhyankar /* Check if the vector sizes are compatible for lower and upper bounds */ 1299*2b4ed282SShri Abhyankar if (xl->map->N != snes->vec_func->map->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths lower bound = %D solution vector = %D",xl->map->N,snes->vec_func->map->N); 1300*2b4ed282SShri Abhyankar if (xu->map->N != snes->vec_func->map->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector lengths: upper bound = %D solution vector = %D",xu->map->N,snes->vec_func->map->N); 1301*2b4ed282SShri Abhyankar vi->usersetxbounds = PETSC_TRUE; 1302*2b4ed282SShri Abhyankar vi->xl = xl; 1303*2b4ed282SShri Abhyankar vi->xu = xu; 1304*2b4ed282SShri Abhyankar 1305*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 1306*2b4ed282SShri Abhyankar } 1307*2b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 1308*2b4ed282SShri Abhyankar /* 1309*2b4ed282SShri Abhyankar SNESSetFromOptions_VI - Sets various parameters for the SNESVI method. 1310*2b4ed282SShri Abhyankar 1311*2b4ed282SShri Abhyankar Input Parameter: 1312*2b4ed282SShri Abhyankar . snes - the SNES context 1313*2b4ed282SShri Abhyankar 1314*2b4ed282SShri Abhyankar Application Interface Routine: SNESSetFromOptions() 1315*2b4ed282SShri Abhyankar */ 1316*2b4ed282SShri Abhyankar #undef __FUNCT__ 1317*2b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetFromOptions_VI" 1318*2b4ed282SShri Abhyankar static PetscErrorCode SNESSetFromOptions_VI(SNES snes) 1319*2b4ed282SShri Abhyankar { 1320*2b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI *)snes->data; 1321*2b4ed282SShri Abhyankar const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 1322*2b4ed282SShri Abhyankar PetscErrorCode ierr; 1323*2b4ed282SShri Abhyankar PetscInt indx; 1324*2b4ed282SShri Abhyankar PetscTruth flg,set; 1325*2b4ed282SShri Abhyankar 1326*2b4ed282SShri Abhyankar PetscFunctionBegin; 1327*2b4ed282SShri Abhyankar ierr = PetscOptionsHead("SNES semismooth method options");CHKERRQ(ierr); 1328*2b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_alpha","Function norm must decrease by","None",vi->alpha,&vi->alpha,0);CHKERRQ(ierr); 1329*2b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_maxstep","Step must be less than","None",vi->maxstep,&vi->maxstep,0);CHKERRQ(ierr); 1330*2b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_minlambda","Minimum lambda allowed","None",vi->minlambda,&vi->minlambda,0);CHKERRQ(ierr); 1331*2b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_delta","descent test fraction","None",vi->delta,&vi->delta,0);CHKERRQ(ierr); 1332*2b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_rho","descent test power","None",vi->rho,&vi->rho,0);CHKERRQ(ierr); 1333*2b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_const_tol","constraint tolerance","None",vi->const_tol,&vi->const_tol,0);CHKERRQ(ierr); 1334*2b4ed282SShri Abhyankar ierr = PetscOptionsTruth("-snes_vi_monitor","Print progress of line searches","SNESLineSearchSetMonitor",vi->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 1335*2b4ed282SShri Abhyankar if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);} 1336*2b4ed282SShri Abhyankar 1337*2b4ed282SShri Abhyankar ierr = PetscOptionsEList("-snes_vi","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 1338*2b4ed282SShri Abhyankar if (flg) { 1339*2b4ed282SShri Abhyankar switch (indx) { 1340*2b4ed282SShri Abhyankar case 0: 1341*2b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr); 1342*2b4ed282SShri Abhyankar break; 1343*2b4ed282SShri Abhyankar case 1: 1344*2b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr); 1345*2b4ed282SShri Abhyankar break; 1346*2b4ed282SShri Abhyankar case 2: 1347*2b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr); 1348*2b4ed282SShri Abhyankar break; 1349*2b4ed282SShri Abhyankar case 3: 1350*2b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr); 1351*2b4ed282SShri Abhyankar break; 1352*2b4ed282SShri Abhyankar } 1353*2b4ed282SShri Abhyankar } 1354*2b4ed282SShri Abhyankar ierr = PetscOptionsTail();CHKERRQ(ierr); 1355*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 1356*2b4ed282SShri Abhyankar } 1357*2b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 1358*2b4ed282SShri Abhyankar /*MC 1359*2b4ed282SShri Abhyankar SNESVI - Semismooth newton method based nonlinear solver that uses a line search 1360*2b4ed282SShri Abhyankar 1361*2b4ed282SShri Abhyankar Options Database: 1362*2b4ed282SShri Abhyankar + -snes_vi [cubic,quadratic,basic,basicnonorms] - Selects line search 1363*2b4ed282SShri Abhyankar . -snes_vi_alpha <alpha> - Sets alpha 1364*2b4ed282SShri Abhyankar . -snes_vi_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y) 1365*2b4ed282SShri Abhyankar . -snes_vi_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 1366*2b4ed282SShri Abhyankar -snes_vi_delta <delta> - Sets the fraction used in the descent test. 1367*2b4ed282SShri Abhyankar -snes_vi_rho <rho> - Sets the power used in the descent test. 1368*2b4ed282SShri Abhyankar For a descent direction to be accepted it has to satisfy the condition dpsi^T*d <= -delta*||d||^rho 1369*2b4ed282SShri Abhyankar - -snes_vi_monitor - print information about progress of line searches 1370*2b4ed282SShri Abhyankar 1371*2b4ed282SShri Abhyankar 1372*2b4ed282SShri Abhyankar Level: beginner 1373*2b4ed282SShri Abhyankar 1374*2b4ed282SShri Abhyankar .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 1375*2b4ed282SShri Abhyankar SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 1376*2b4ed282SShri Abhyankar SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams() 1377*2b4ed282SShri Abhyankar 1378*2b4ed282SShri Abhyankar M*/ 1379*2b4ed282SShri Abhyankar EXTERN_C_BEGIN 1380*2b4ed282SShri Abhyankar #undef __FUNCT__ 1381*2b4ed282SShri Abhyankar #define __FUNCT__ "SNESCreate_VI" 1382*2b4ed282SShri Abhyankar PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_VI(SNES snes) 1383*2b4ed282SShri Abhyankar { 1384*2b4ed282SShri Abhyankar PetscErrorCode ierr; 1385*2b4ed282SShri Abhyankar SNES_VI *vi; 1386*2b4ed282SShri Abhyankar 1387*2b4ed282SShri Abhyankar PetscFunctionBegin; 1388*2b4ed282SShri Abhyankar snes->ops->setup = SNESSetUp_VI; 1389*2b4ed282SShri Abhyankar snes->ops->solve = SNESSolve_VI; 1390*2b4ed282SShri Abhyankar snes->ops->destroy = SNESDestroy_VI; 1391*2b4ed282SShri Abhyankar snes->ops->setfromoptions = SNESSetFromOptions_VI; 1392*2b4ed282SShri Abhyankar snes->ops->view = SNESView_VI; 1393*2b4ed282SShri Abhyankar snes->ops->converged = SNESDefaultConverged_VI; 1394*2b4ed282SShri Abhyankar 1395*2b4ed282SShri Abhyankar ierr = PetscNewLog(snes,SNES_VI,&vi);CHKERRQ(ierr); 1396*2b4ed282SShri Abhyankar snes->data = (void*)vi; 1397*2b4ed282SShri Abhyankar vi->alpha = 1.e-4; 1398*2b4ed282SShri Abhyankar vi->maxstep = 1.e8; 1399*2b4ed282SShri Abhyankar vi->minlambda = 1.e-12; 1400*2b4ed282SShri Abhyankar vi->LineSearch = SNESLineSearchCubic_VI; 1401*2b4ed282SShri Abhyankar vi->lsP = PETSC_NULL; 1402*2b4ed282SShri Abhyankar vi->postcheckstep = PETSC_NULL; 1403*2b4ed282SShri Abhyankar vi->postcheck = PETSC_NULL; 1404*2b4ed282SShri Abhyankar vi->precheckstep = PETSC_NULL; 1405*2b4ed282SShri Abhyankar vi->precheck = PETSC_NULL; 1406*2b4ed282SShri Abhyankar vi->rho = 2.1; 1407*2b4ed282SShri Abhyankar vi->delta = 1e-10; 1408*2b4ed282SShri Abhyankar vi->const_tol = 2.2204460492503131e-16; 1409*2b4ed282SShri Abhyankar vi->computessfunction = ComputeFischerFunction; 1410*2b4ed282SShri Abhyankar 1411*2b4ed282SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_VI",SNESLineSearchSetMonitor_VI);CHKERRQ(ierr); 1412*2b4ed282SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_VI",SNESLineSearchSet_VI);CHKERRQ(ierr); 1413*2b4ed282SShri Abhyankar 1414*2b4ed282SShri Abhyankar PetscFunctionReturn(0); 1415*2b4ed282SShri Abhyankar } 1416*2b4ed282SShri Abhyankar EXTERN_C_END 1417