12b4ed282SShri Abhyankar #define PETSCSNES_DLL 22b4ed282SShri Abhyankar 32b4ed282SShri Abhyankar #include "../src/snes/impls/vi/viimpl.h" 42b4ed282SShri Abhyankar 52b4ed282SShri Abhyankar /* 62b4ed282SShri Abhyankar Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function, 72b4ed282SShri 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 82b4ed282SShri Abhyankar 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 92b4ed282SShri Abhyankar for this trick. One assumes that the probability that W is in the null space of J is very, very small. 102b4ed282SShri Abhyankar */ 112b4ed282SShri Abhyankar #undef __FUNCT__ 122b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckLocalMin_Private" 13ace3abfcSBarry Smith PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin) 142b4ed282SShri Abhyankar { 152b4ed282SShri Abhyankar PetscReal a1; 162b4ed282SShri Abhyankar PetscErrorCode ierr; 17ace3abfcSBarry Smith PetscBool hastranspose; 182b4ed282SShri Abhyankar 192b4ed282SShri Abhyankar PetscFunctionBegin; 202b4ed282SShri Abhyankar *ismin = PETSC_FALSE; 212b4ed282SShri Abhyankar ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 222b4ed282SShri Abhyankar if (hastranspose) { 232b4ed282SShri Abhyankar /* Compute || J^T F|| */ 242b4ed282SShri Abhyankar ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 252b4ed282SShri Abhyankar ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 262b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"|| J^T F|| %G near zero implies found a local minimum\n",a1/fnorm);CHKERRQ(ierr); 272b4ed282SShri Abhyankar if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 282b4ed282SShri Abhyankar } else { 292b4ed282SShri Abhyankar Vec work; 302b4ed282SShri Abhyankar PetscScalar result; 312b4ed282SShri Abhyankar PetscReal wnorm; 322b4ed282SShri Abhyankar 332b4ed282SShri Abhyankar ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr); 342b4ed282SShri Abhyankar ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr); 352b4ed282SShri Abhyankar ierr = VecDuplicate(W,&work);CHKERRQ(ierr); 362b4ed282SShri Abhyankar ierr = MatMult(A,W,work);CHKERRQ(ierr); 372b4ed282SShri Abhyankar ierr = VecDot(F,work,&result);CHKERRQ(ierr); 382b4ed282SShri Abhyankar ierr = VecDestroy(work);CHKERRQ(ierr); 392b4ed282SShri Abhyankar a1 = PetscAbsScalar(result)/(fnorm*wnorm); 402b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %G near zero implies found a local minimum\n",a1);CHKERRQ(ierr); 412b4ed282SShri Abhyankar if (a1 < 1.e-4) *ismin = PETSC_TRUE; 422b4ed282SShri Abhyankar } 432b4ed282SShri Abhyankar PetscFunctionReturn(0); 442b4ed282SShri Abhyankar } 452b4ed282SShri Abhyankar 462b4ed282SShri Abhyankar /* 472b4ed282SShri Abhyankar Checks if J^T(F - J*X) = 0 482b4ed282SShri Abhyankar */ 492b4ed282SShri Abhyankar #undef __FUNCT__ 502b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckResidual_Private" 512b4ed282SShri Abhyankar PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2) 522b4ed282SShri Abhyankar { 532b4ed282SShri Abhyankar PetscReal a1,a2; 542b4ed282SShri Abhyankar PetscErrorCode ierr; 55ace3abfcSBarry Smith PetscBool hastranspose; 562b4ed282SShri Abhyankar 572b4ed282SShri Abhyankar PetscFunctionBegin; 582b4ed282SShri Abhyankar ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 592b4ed282SShri Abhyankar if (hastranspose) { 602b4ed282SShri Abhyankar ierr = MatMult(A,X,W1);CHKERRQ(ierr); 612b4ed282SShri Abhyankar ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr); 622b4ed282SShri Abhyankar 632b4ed282SShri Abhyankar /* Compute || J^T W|| */ 642b4ed282SShri Abhyankar ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 652b4ed282SShri Abhyankar ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 662b4ed282SShri Abhyankar ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 672b4ed282SShri Abhyankar if (a1 != 0.0) { 682b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %G near zero implies inconsistent rhs\n",a2/a1);CHKERRQ(ierr); 692b4ed282SShri Abhyankar } 702b4ed282SShri Abhyankar } 712b4ed282SShri Abhyankar PetscFunctionReturn(0); 722b4ed282SShri Abhyankar } 732b4ed282SShri Abhyankar 742b4ed282SShri Abhyankar /* 752b4ed282SShri Abhyankar SNESDefaultConverged_VI - Checks the convergence of the semismooth newton algorithm. 762b4ed282SShri Abhyankar 772b4ed282SShri Abhyankar Notes: 782b4ed282SShri Abhyankar The convergence criterion currently implemented is 792b4ed282SShri Abhyankar merit < abstol 802b4ed282SShri Abhyankar merit < rtol*merit_initial 812b4ed282SShri Abhyankar */ 822b4ed282SShri Abhyankar #undef __FUNCT__ 832b4ed282SShri Abhyankar #define __FUNCT__ "SNESDefaultConverged_VI" 842b4ed282SShri Abhyankar PetscErrorCode SNESDefaultConverged_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal merit,SNESConvergedReason *reason,void *dummy) 852b4ed282SShri Abhyankar { 862b4ed282SShri Abhyankar PetscErrorCode ierr; 872b4ed282SShri Abhyankar 882b4ed282SShri Abhyankar PetscFunctionBegin; 892b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 902b4ed282SShri Abhyankar PetscValidPointer(reason,6); 912b4ed282SShri Abhyankar 922b4ed282SShri Abhyankar *reason = SNES_CONVERGED_ITERATING; 932b4ed282SShri Abhyankar 942b4ed282SShri Abhyankar if (!it) { 952b4ed282SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 962b4ed282SShri Abhyankar snes->ttol = merit*snes->rtol; 972b4ed282SShri Abhyankar } 982b4ed282SShri Abhyankar if (merit != merit) { 992b4ed282SShri Abhyankar ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr); 1002b4ed282SShri Abhyankar *reason = SNES_DIVERGED_FNORM_NAN; 1012b4ed282SShri Abhyankar } else if (merit < snes->abstol) { 1022b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"Converged due to merit function %G < %G\n",merit,snes->abstol);CHKERRQ(ierr); 1032b4ed282SShri Abhyankar *reason = SNES_CONVERGED_FNORM_ABS; 1042b4ed282SShri Abhyankar } else if (snes->nfuncs >= snes->max_funcs) { 1052b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr); 1062b4ed282SShri Abhyankar *reason = SNES_DIVERGED_FUNCTION_COUNT; 1072b4ed282SShri Abhyankar } 1082b4ed282SShri Abhyankar 1092b4ed282SShri Abhyankar if (it && !*reason) { 1102b4ed282SShri Abhyankar if (merit < snes->ttol) { 1112b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"Converged due to merit function %G < %G (relative tolerance)\n",merit,snes->ttol);CHKERRQ(ierr); 1122b4ed282SShri Abhyankar *reason = SNES_CONVERGED_FNORM_RELATIVE; 1132b4ed282SShri Abhyankar } 1142b4ed282SShri Abhyankar } 1152b4ed282SShri Abhyankar PetscFunctionReturn(0); 1162b4ed282SShri Abhyankar } 1172b4ed282SShri Abhyankar 1182b4ed282SShri Abhyankar /* 1192b4ed282SShri Abhyankar SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem. 1202b4ed282SShri Abhyankar 1212b4ed282SShri Abhyankar Input Parameter: 1222b4ed282SShri Abhyankar . phi - the semismooth function 1232b4ed282SShri Abhyankar 1242b4ed282SShri Abhyankar Output Parameter: 1252b4ed282SShri Abhyankar . merit - the merit function 1262b4ed282SShri Abhyankar . phinorm - ||phi|| 1272b4ed282SShri Abhyankar 1282b4ed282SShri Abhyankar Notes: 1292b4ed282SShri Abhyankar The merit function for the mixed complementarity problem is defined as 1302b4ed282SShri Abhyankar merit = 0.5*phi^T*phi 1312b4ed282SShri Abhyankar */ 1322b4ed282SShri Abhyankar #undef __FUNCT__ 1332b4ed282SShri Abhyankar #define __FUNCT__ "SNESVIComputeMeritFunction" 1342b4ed282SShri Abhyankar static PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal* merit,PetscReal* phinorm) 1352b4ed282SShri Abhyankar { 1362b4ed282SShri Abhyankar PetscErrorCode ierr; 1372b4ed282SShri Abhyankar 1382b4ed282SShri Abhyankar PetscFunctionBegin; 1392b4ed282SShri Abhyankar ierr = VecNormBegin(phi,NORM_2,phinorm); 1402b4ed282SShri Abhyankar ierr = VecNormEnd(phi,NORM_2,phinorm); 1412b4ed282SShri Abhyankar 1422b4ed282SShri Abhyankar *merit = 0.5*(*phinorm)*(*phinorm); 1432b4ed282SShri Abhyankar PetscFunctionReturn(0); 1442b4ed282SShri Abhyankar } 1452b4ed282SShri Abhyankar 1462b4ed282SShri Abhyankar /* 1472b4ed282SShri Abhyankar ComputeFischerFunction - Computes the semismooth fischer burmeister function for a mixed complementarity equation. 1482b4ed282SShri Abhyankar 1492b4ed282SShri Abhyankar Notes: 1502b4ed282SShri Abhyankar The Fischer-Burmeister function is defined as 1512b4ed282SShri Abhyankar ff(a,b) := sqrt(a*a + b*b) - a - b 1522b4ed282SShri Abhyankar and is used reformulate a complementarity equation as a semismooth equation. 1532b4ed282SShri Abhyankar */ 1542b4ed282SShri Abhyankar 1552b4ed282SShri Abhyankar #undef __FUNCT__ 1562b4ed282SShri Abhyankar #define __FUNCT__ "ComputeFischerFunction" 1572b4ed282SShri Abhyankar static PetscErrorCode ComputeFischerFunction(PetscScalar a, PetscScalar b, PetscScalar* ff) 1582b4ed282SShri Abhyankar { 1592b4ed282SShri Abhyankar PetscFunctionBegin; 1602b4ed282SShri Abhyankar *ff = sqrt(a*a + b*b) - a - b; 1612b4ed282SShri Abhyankar PetscFunctionReturn(0); 1622b4ed282SShri Abhyankar } 1632b4ed282SShri Abhyankar 1642b4ed282SShri Abhyankar /* 1651f79c099SShri Abhyankar SNESVIComputeFunction - Reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear equations in semismooth form. 1662b4ed282SShri Abhyankar 1672b4ed282SShri Abhyankar Input Parameters: 1682b4ed282SShri Abhyankar . snes - the SNES context 1692b4ed282SShri Abhyankar . x - current iterate 1702b4ed282SShri Abhyankar . functx - user defined function context 1712b4ed282SShri Abhyankar 1722b4ed282SShri Abhyankar Output Parameters: 1732b4ed282SShri Abhyankar . phi - Semismooth function 1742b4ed282SShri Abhyankar 1752b4ed282SShri Abhyankar Notes: 1762b4ed282SShri Abhyankar The result of this function is done by cases: 1772b4ed282SShri Abhyankar + l[i] == -infinity, u[i] == infinity -- phi[i] = -f[i] 1782b4ed282SShri Abhyankar . l[i] == -infinity, u[i] finite -- phi[i] = ss(u[i]-x[i], -f[i]) 1792b4ed282SShri Abhyankar . l[i] finite, u[i] == infinity -- phi[i] = ss(x[i]-l[i], f[i]) 1802b4ed282SShri Abhyankar . l[i] finite < u[i] finite -- phi[i] = phi(x[i]-l[i], ss(u[i]-x[i], -f[u])) 1812b4ed282SShri Abhyankar - otherwise l[i] == u[i] -- phi[i] = l[i] - x[i] 1822b4ed282SShri Abhyankar ss is the semismoothing function used to reformulate the nonlinear equation in complementarity 1832b4ed282SShri Abhyankar form to semismooth form 1842b4ed282SShri Abhyankar 1852b4ed282SShri Abhyankar */ 1862b4ed282SShri Abhyankar #undef __FUNCT__ 1871f79c099SShri Abhyankar #define __FUNCT__ "SNESVIComputeFunction" 1881f79c099SShri Abhyankar static PetscErrorCode SNESVIComputeFunction(SNES snes,Vec X,Vec phi,void* functx) 1892b4ed282SShri Abhyankar { 1902b4ed282SShri Abhyankar PetscErrorCode ierr; 1912b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1922b4ed282SShri Abhyankar Vec Xl = vi->xl,Xu = vi->xu,F = snes->vec_func; 1932b4ed282SShri Abhyankar PetscScalar *phi_arr,*x_arr,*f_arr,*l,*u,t; 1942b4ed282SShri Abhyankar PetscInt i,nlocal; 1952b4ed282SShri Abhyankar 1962b4ed282SShri Abhyankar PetscFunctionBegin; 1972b4ed282SShri Abhyankar ierr = (*vi->computeuserfunction)(snes,X,F,functx);CHKERRQ(ierr); 1982b4ed282SShri Abhyankar ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 1992b4ed282SShri Abhyankar 2002b4ed282SShri Abhyankar ierr = VecGetArray(X,&x_arr);CHKERRQ(ierr); 2012b4ed282SShri Abhyankar ierr = VecGetArray(F,&f_arr);CHKERRQ(ierr); 2022b4ed282SShri Abhyankar ierr = VecGetArray(Xl,&l);CHKERRQ(ierr); 2032b4ed282SShri Abhyankar ierr = VecGetArray(Xu,&u);CHKERRQ(ierr); 2042b4ed282SShri Abhyankar ierr = VecGetArray(phi,&phi_arr);CHKERRQ(ierr); 2052b4ed282SShri Abhyankar 2062b4ed282SShri Abhyankar for (i=0;i < nlocal;i++) { 2072b4ed282SShri Abhyankar if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) { 2082b4ed282SShri Abhyankar phi_arr[i] = -f_arr[i]; 2092b4ed282SShri Abhyankar } 2102b4ed282SShri Abhyankar else if (l[i] <= PETSC_VI_NINF) { 2112b4ed282SShri Abhyankar t = u[i] - x_arr[i]; 2122b4ed282SShri Abhyankar ierr = ComputeFischerFunction(t,-f_arr[i],&phi_arr[i]);CHKERRQ(ierr); 2132b4ed282SShri Abhyankar phi_arr[i] = -phi_arr[i]; 2142b4ed282SShri Abhyankar } 2152b4ed282SShri Abhyankar else if (u[i] >= PETSC_VI_INF) { 2162b4ed282SShri Abhyankar t = x_arr[i] - l[i]; 2172b4ed282SShri Abhyankar ierr = ComputeFischerFunction(t,f_arr[i],&phi_arr[i]);CHKERRQ(ierr); 2182b4ed282SShri Abhyankar } 2192b4ed282SShri Abhyankar else if (l[i] == u[i]) { 2202b4ed282SShri Abhyankar phi_arr[i] = l[i] - x_arr[i]; 2212b4ed282SShri Abhyankar } 2222b4ed282SShri Abhyankar else { 2232b4ed282SShri Abhyankar t = u[i] - x_arr[i]; 2242b4ed282SShri Abhyankar ierr = ComputeFischerFunction(t,-f_arr[i],&phi_arr[i]); 2252b4ed282SShri Abhyankar t = x_arr[i] - l[i]; 2262b4ed282SShri Abhyankar ierr = ComputeFischerFunction(t,phi_arr[i],&phi_arr[i]); 2272b4ed282SShri Abhyankar } 2282b4ed282SShri Abhyankar } 2292b4ed282SShri Abhyankar 2302b4ed282SShri Abhyankar ierr = VecRestoreArray(X,&x_arr);CHKERRQ(ierr); 2312b4ed282SShri Abhyankar ierr = VecRestoreArray(F,&f_arr);CHKERRQ(ierr); 2322b4ed282SShri Abhyankar ierr = VecRestoreArray(Xl,&l);CHKERRQ(ierr); 2332b4ed282SShri Abhyankar ierr = VecRestoreArray(Xu,&u);CHKERRQ(ierr); 2342b4ed282SShri Abhyankar ierr = VecRestoreArray(phi,&phi_arr);CHKERRQ(ierr); 2352b4ed282SShri Abhyankar 2362b4ed282SShri Abhyankar PetscFunctionReturn(0); 2372b4ed282SShri Abhyankar } 2382b4ed282SShri Abhyankar 2392b4ed282SShri Abhyankar /* 240a79edbefSShri Abhyankar SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the 241a79edbefSShri Abhyankar the semismooth jacobian. 2422b4ed282SShri Abhyankar */ 2432b4ed282SShri Abhyankar #undef __FUNCT__ 244a79edbefSShri Abhyankar #define __FUNCT__ "SNESVIComputeBsubdifferentialVectors" 245a79edbefSShri Abhyankar PetscErrorCode SNESVIComputeBsubdifferentialVectors(SNES snes,Vec X,Vec F,Mat jac,Vec Da,Vec Db) 2462b4ed282SShri Abhyankar { 2472b4ed282SShri Abhyankar PetscErrorCode ierr; 2482b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 2492b4ed282SShri Abhyankar PetscScalar *l,*u,*x,*f,*da,*db,*z,*t,t1,t2,ci,di,ei; 2502b4ed282SShri Abhyankar PetscInt i,nlocal; 2512b4ed282SShri Abhyankar 2522b4ed282SShri Abhyankar PetscFunctionBegin; 2532b4ed282SShri Abhyankar 2542b4ed282SShri Abhyankar ierr = VecGetArray(X,&x);CHKERRQ(ierr); 2552b4ed282SShri Abhyankar ierr = VecGetArray(F,&f);CHKERRQ(ierr); 2562b4ed282SShri Abhyankar ierr = VecGetArray(vi->xl,&l);CHKERRQ(ierr); 2572b4ed282SShri Abhyankar ierr = VecGetArray(vi->xu,&u);CHKERRQ(ierr); 258a79edbefSShri Abhyankar ierr = VecGetArray(Da,&da);CHKERRQ(ierr); 259a79edbefSShri Abhyankar ierr = VecGetArray(Db,&db);CHKERRQ(ierr); 2602b4ed282SShri Abhyankar ierr = VecGetArray(vi->z,&z);CHKERRQ(ierr); 2612b4ed282SShri Abhyankar 2622b4ed282SShri Abhyankar ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 2632b4ed282SShri Abhyankar /* Set the elements of the vector z: 2642b4ed282SShri Abhyankar z[i] = 1 if (x[i] - l[i],f[i]) = (0,0) or (u[i] - x[i],f[i]) = (0,0) 2652b4ed282SShri Abhyankar else z[i] = 0 2662b4ed282SShri Abhyankar */ 2672b4ed282SShri Abhyankar for(i=0;i < nlocal;i++) { 2682b4ed282SShri Abhyankar da[i] = db[i] = z[i] = 0; 2692b4ed282SShri Abhyankar if(PetscAbsScalar(f[i]) <= vi->const_tol) { 2702b4ed282SShri Abhyankar if ((l[i] > PETSC_VI_NINF) && (PetscAbsScalar(x[i]-l[i]) <= vi->const_tol)) { 2712b4ed282SShri Abhyankar da[i] = 1; 2722b4ed282SShri Abhyankar z[i] = 1; 2732b4ed282SShri Abhyankar } 2742b4ed282SShri Abhyankar if ((u[i] < PETSC_VI_INF) && (PetscAbsScalar(u[i]-x[i]) <= vi->const_tol)) { 2752b4ed282SShri Abhyankar db[i] = 1; 2762b4ed282SShri Abhyankar z[i] = 1; 2772b4ed282SShri Abhyankar } 2782b4ed282SShri Abhyankar } 2792b4ed282SShri Abhyankar } 2802b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->z,&z);CHKERRQ(ierr); 281a79edbefSShri Abhyankar ierr = MatMult(jac,vi->z,vi->t);CHKERRQ(ierr); 2822b4ed282SShri Abhyankar ierr = VecGetArray(vi->t,&t);CHKERRQ(ierr); 2832b4ed282SShri Abhyankar /* Compute the elements of the diagonal perturbation vector Da and row scaling vector Db */ 2842b4ed282SShri Abhyankar for(i=0;i< nlocal;i++) { 2852b4ed282SShri Abhyankar /* Free variables */ 2862b4ed282SShri Abhyankar if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) { 2872b4ed282SShri Abhyankar da[i] = 0; db[i] = -1; 2882b4ed282SShri Abhyankar } 2892b4ed282SShri Abhyankar /* lower bounded variables */ 2902b4ed282SShri Abhyankar else if (u[i] >= PETSC_VI_INF) { 2912b4ed282SShri Abhyankar if (da[i] >= 1) { 2922b4ed282SShri Abhyankar t2 = PetscScalarNorm(1,t[i]); 2932b4ed282SShri Abhyankar da[i] = 1/t2 - 1; 2942b4ed282SShri Abhyankar db[i] = t[i]/t2 - 1; 2952b4ed282SShri Abhyankar } else { 2962b4ed282SShri Abhyankar t1 = x[i] - l[i]; 2972b4ed282SShri Abhyankar t2 = PetscScalarNorm(t1,f[i]); 2982b4ed282SShri Abhyankar da[i] = t1/t2 - 1; 2992b4ed282SShri Abhyankar db[i] = f[i]/t2 - 1; 3002b4ed282SShri Abhyankar } 3012b4ed282SShri Abhyankar } 3022b4ed282SShri Abhyankar /* upper bounded variables */ 3032b4ed282SShri Abhyankar else if (l[i] <= PETSC_VI_NINF) { 3042b4ed282SShri Abhyankar if (db[i] >= 1) { 3052b4ed282SShri Abhyankar t2 = PetscScalarNorm(1,t[i]); 3062b4ed282SShri Abhyankar da[i] = -1/t2 -1; 3072b4ed282SShri Abhyankar db[i] = -t[i]/t2 - 1; 3082b4ed282SShri Abhyankar } 3092b4ed282SShri Abhyankar else { 3102b4ed282SShri Abhyankar t1 = u[i] - x[i]; 3112b4ed282SShri Abhyankar t2 = PetscScalarNorm(t1,f[i]); 3122b4ed282SShri Abhyankar da[i] = t1/t2 - 1; 3132b4ed282SShri Abhyankar db[i] = -f[i]/t2 - 1; 3142b4ed282SShri Abhyankar } 3152b4ed282SShri Abhyankar } 3162b4ed282SShri Abhyankar /* Fixed variables */ 3172b4ed282SShri Abhyankar else if (l[i] == u[i]) { 3182b4ed282SShri Abhyankar da[i] = -1; 3192b4ed282SShri Abhyankar db[i] = 0; 3202b4ed282SShri Abhyankar } 3212b4ed282SShri Abhyankar /* Box constrained variables */ 3222b4ed282SShri Abhyankar else { 3232b4ed282SShri Abhyankar if (db[i] >= 1) { 3242b4ed282SShri Abhyankar t2 = PetscScalarNorm(1,t[i]); 3252b4ed282SShri Abhyankar ci = 1/t2 + 1; 3262b4ed282SShri Abhyankar di = t[i]/t2 + 1; 3272b4ed282SShri Abhyankar } 3282b4ed282SShri Abhyankar else { 3292b4ed282SShri Abhyankar t1 = x[i] - u[i]; 3302b4ed282SShri Abhyankar t2 = PetscScalarNorm(t1,f[i]); 3312b4ed282SShri Abhyankar ci = t1/t2 + 1; 3322b4ed282SShri Abhyankar di = f[i]/t2 + 1; 3332b4ed282SShri Abhyankar } 3342b4ed282SShri Abhyankar 3352b4ed282SShri Abhyankar if (da[i] >= 1) { 3362b4ed282SShri Abhyankar t1 = ci + di*t[i]; 3372b4ed282SShri Abhyankar t2 = PetscScalarNorm(1,t1); 3382b4ed282SShri Abhyankar t1 = t1/t2 - 1; 3392b4ed282SShri Abhyankar t2 = 1/t2 - 1; 3402b4ed282SShri Abhyankar } 3412b4ed282SShri Abhyankar else { 3422b4ed282SShri Abhyankar ierr = ComputeFischerFunction(u[i]-x[i],-f[i],&ei);CHKERRQ(ierr); 3432b4ed282SShri Abhyankar t2 = PetscScalarNorm(x[i]-l[i],ei); 3442b4ed282SShri Abhyankar t1 = ei/t2 - 1; 3452b4ed282SShri Abhyankar t2 = (x[i] - l[i])/t2 - 1; 3462b4ed282SShri Abhyankar } 3472b4ed282SShri Abhyankar 3482b4ed282SShri Abhyankar da[i] = t2 + t1*ci; 3492b4ed282SShri Abhyankar db[i] = t1*di; 3502b4ed282SShri Abhyankar } 3512b4ed282SShri Abhyankar } 3522b4ed282SShri Abhyankar ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 3532b4ed282SShri Abhyankar ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 3542b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->xl,&l);CHKERRQ(ierr); 3552b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->xu,&u);CHKERRQ(ierr); 356a79edbefSShri Abhyankar ierr = VecRestoreArray(Da,&da);CHKERRQ(ierr); 357a79edbefSShri Abhyankar ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr); 3582b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->t,&t);CHKERRQ(ierr); 3592b4ed282SShri Abhyankar 360a79edbefSShri Abhyankar PetscFunctionReturn(0); 361a79edbefSShri Abhyankar } 362a79edbefSShri Abhyankar 363a79edbefSShri Abhyankar /* 364a79edbefSShri Abhyankar SNESVIComputeJacobian - 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. 365a79edbefSShri Abhyankar 366a79edbefSShri Abhyankar Input Parameters: 367a79edbefSShri Abhyankar . Da - Diagonal shift vector for the semismooth jacobian. 368a79edbefSShri Abhyankar . Db - Row scaling vector for the semismooth jacobian. 369a79edbefSShri Abhyankar 370a79edbefSShri Abhyankar Output Parameters: 371a79edbefSShri Abhyankar . jac - semismooth jacobian 372a79edbefSShri Abhyankar . jac_pre - optional preconditioning matrix 373a79edbefSShri Abhyankar 374a79edbefSShri Abhyankar Notes: 375a79edbefSShri Abhyankar The semismooth jacobian matrix is given by 376a79edbefSShri Abhyankar jac = Da + Db*jacfun 377a79edbefSShri Abhyankar where Db is the row scaling matrix stored as a vector, 378a79edbefSShri Abhyankar Da is the diagonal perturbation matrix stored as a vector 379a79edbefSShri Abhyankar and jacfun is the jacobian of the original nonlinear function. 380a79edbefSShri Abhyankar */ 381a79edbefSShri Abhyankar #undef __FUNCT__ 382a79edbefSShri Abhyankar #define __FUNCT__ "SNESVIComputeJacobian" 383a79edbefSShri Abhyankar PetscErrorCode SNESVIComputeJacobian(Mat jac, Mat jac_pre,Vec Da, Vec Db) 384a79edbefSShri Abhyankar { 385a79edbefSShri Abhyankar PetscErrorCode ierr; 386a79edbefSShri Abhyankar 3872b4ed282SShri Abhyankar /* Do row scaling and add diagonal perturbation */ 388a79edbefSShri Abhyankar ierr = MatDiagonalScale(jac,Db,PETSC_NULL);CHKERRQ(ierr); 389a79edbefSShri Abhyankar ierr = MatDiagonalSet(jac,Da,ADD_VALUES);CHKERRQ(ierr); 390a79edbefSShri Abhyankar if (jac != jac_pre) { /* If jac and jac_pre are different */ 391a79edbefSShri Abhyankar ierr = MatDiagonalScale(jac_pre,Db,PETSC_NULL); 392a79edbefSShri Abhyankar ierr = MatDiagonalSet(jac_pre,Da,ADD_VALUES);CHKERRQ(ierr); 3932b4ed282SShri Abhyankar } 3942b4ed282SShri Abhyankar 3952b4ed282SShri Abhyankar PetscFunctionReturn(0); 3962b4ed282SShri Abhyankar } 3972b4ed282SShri Abhyankar 3982b4ed282SShri Abhyankar /* 3992b4ed282SShri Abhyankar SNESVIComputeMeritFunctionGradient - Computes the gradient of the merit function psi. 4002b4ed282SShri Abhyankar 4012b4ed282SShri Abhyankar Input Parameters: 4022b4ed282SShri Abhyankar . phi - semismooth function. 4032b4ed282SShri Abhyankar . H - semismooth jacobian 4042b4ed282SShri Abhyankar 4052b4ed282SShri Abhyankar Output Parameters: 4062b4ed282SShri Abhyankar . dpsi - merit function gradient 4072b4ed282SShri Abhyankar 4082b4ed282SShri Abhyankar Notes: 4092b4ed282SShri Abhyankar The merit function gradient is computed as follows 4102b4ed282SShri Abhyankar dpsi = H^T*phi 4112b4ed282SShri Abhyankar */ 4122b4ed282SShri Abhyankar #undef __FUNCT__ 4132b4ed282SShri Abhyankar #define __FUNCT__ "SNESVIComputeMeritFunctionGradient" 4142b4ed282SShri Abhyankar PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi) 4152b4ed282SShri Abhyankar { 4162b4ed282SShri Abhyankar PetscErrorCode ierr; 4172b4ed282SShri Abhyankar 4182b4ed282SShri Abhyankar PetscFunctionBegin; 4192b4ed282SShri Abhyankar ierr = MatMultTranspose(H,phi,dpsi); 4202b4ed282SShri Abhyankar 4212b4ed282SShri Abhyankar PetscFunctionReturn(0); 4222b4ed282SShri Abhyankar } 4232b4ed282SShri Abhyankar 4242b4ed282SShri Abhyankar /* 4252b4ed282SShri Abhyankar SNESVISetDescentDirection - Sets the descent direction for the semismooth algorithm 4262b4ed282SShri Abhyankar 4272b4ed282SShri Abhyankar Input Parameters: 4282b4ed282SShri Abhyankar . snes - the SNES context. 4292b4ed282SShri Abhyankar . dpsi - gradient of the merit function. 4302b4ed282SShri Abhyankar 4312b4ed282SShri Abhyankar Output Parameters: 4322b4ed282SShri Abhyankar . flg - PETSC_TRUE if the sufficient descent condition is not satisfied. 4332b4ed282SShri Abhyankar 4342b4ed282SShri Abhyankar Notes: 4352b4ed282SShri Abhyankar The condition for the sufficient descent direction is 4362b4ed282SShri Abhyankar dpsi^T*Y > delta*||Y||^rho 4372b4ed282SShri Abhyankar where rho, delta are as defined in the SNES_VI structure. 4382b4ed282SShri Abhyankar If this condition is satisfied then the existing descent direction i.e. 4392b4ed282SShri 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. 4402b4ed282SShri Abhyankar */ 4412b4ed282SShri Abhyankar #undef __FUNCT__ 4422b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckDescentDirection" 443ace3abfcSBarry Smith PetscErrorCode SNESVICheckDescentDirection(SNES snes,Vec dpsi, Vec Y,PetscBool* flg) 4442b4ed282SShri Abhyankar { 4452b4ed282SShri Abhyankar PetscErrorCode ierr; 4462b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 4472b4ed282SShri Abhyankar PetscScalar dpsidotY; 4482b4ed282SShri Abhyankar PetscReal norm_Y,rhs; 4492b4ed282SShri Abhyankar const PetscReal rho = vi->rho,delta=vi->delta; 4502b4ed282SShri Abhyankar 4512b4ed282SShri Abhyankar PetscFunctionBegin; 4522b4ed282SShri Abhyankar 4532b4ed282SShri Abhyankar *flg = PETSC_FALSE; 4542b4ed282SShri Abhyankar ierr = VecDot(dpsi,Y,&dpsidotY);CHKERRQ(ierr); 4552b4ed282SShri Abhyankar ierr = VecNormBegin(Y,NORM_2,&norm_Y);CHKERRQ(ierr); 4562b4ed282SShri Abhyankar ierr = VecNormEnd(Y,NORM_2,&norm_Y);CHKERRQ(ierr); 4572b4ed282SShri Abhyankar 4582b4ed282SShri Abhyankar rhs = delta*PetscPowScalar(norm_Y,rho); 4592b4ed282SShri Abhyankar 460887080aaSShri Abhyankar if (dpsidotY <= rhs) *flg = PETSC_TRUE; 4612b4ed282SShri Abhyankar 4622b4ed282SShri Abhyankar PetscFunctionReturn(0); 4632b4ed282SShri Abhyankar } 4642b4ed282SShri Abhyankar 4652b4ed282SShri Abhyankar /* 4662b4ed282SShri Abhyankar SNESVIAdjustInitialGuess - Readjusts the initial guess to the SNES solver supplied by the user so that the initial guess lies inside the feasible region . 4672b4ed282SShri Abhyankar 4682b4ed282SShri Abhyankar Input Parameters: 4692b4ed282SShri Abhyankar . lb - lower bound. 4702b4ed282SShri Abhyankar . ub - upper bound. 4712b4ed282SShri Abhyankar 4722b4ed282SShri Abhyankar Output Parameters: 4732b4ed282SShri Abhyankar . X - the readjusted initial guess. 4742b4ed282SShri Abhyankar 4752b4ed282SShri Abhyankar Notes: 4762b4ed282SShri Abhyankar The readjusted initial guess X[i] = max(max(min(l[i],X[i]),min(X[i],u[i])),min(l[i],u[i])) 4772b4ed282SShri Abhyankar */ 4782b4ed282SShri Abhyankar #undef __FUNCT__ 4792b4ed282SShri Abhyankar #define __FUNCT__ "SNESVIAdjustInitialGuess" 4802b4ed282SShri Abhyankar PetscErrorCode SNESVIAdjustInitialGuess(Vec X, Vec lb, Vec ub) 4812b4ed282SShri Abhyankar { 4822b4ed282SShri Abhyankar PetscErrorCode ierr; 4832b4ed282SShri Abhyankar PetscInt i,nlocal; 4842b4ed282SShri Abhyankar PetscScalar *x,*l,*u; 4852b4ed282SShri Abhyankar 4862b4ed282SShri Abhyankar PetscFunctionBegin; 4872b4ed282SShri Abhyankar 4882b4ed282SShri Abhyankar ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 4892b4ed282SShri Abhyankar ierr = VecGetArray(X,&x);CHKERRQ(ierr); 4902b4ed282SShri Abhyankar ierr = VecGetArray(lb,&l);CHKERRQ(ierr); 4912b4ed282SShri Abhyankar ierr = VecGetArray(ub,&u);CHKERRQ(ierr); 4922b4ed282SShri Abhyankar 4932b4ed282SShri Abhyankar for(i = 0; i < nlocal; i++) { 4942b4ed282SShri Abhyankar x[i] = PetscMax(PetscMax(PetscMin(x[i],l[i]),PetscMin(x[i],u[i])),PetscMin(l[i],u[i])); 4952b4ed282SShri Abhyankar } 4962b4ed282SShri Abhyankar 4972b4ed282SShri Abhyankar ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 4982b4ed282SShri Abhyankar ierr = VecRestoreArray(lb,&l);CHKERRQ(ierr); 4992b4ed282SShri Abhyankar ierr = VecRestoreArray(ub,&u);CHKERRQ(ierr); 5002b4ed282SShri Abhyankar 5012b4ed282SShri Abhyankar PetscFunctionReturn(0); 5022b4ed282SShri Abhyankar } 5032b4ed282SShri Abhyankar 5042b4ed282SShri Abhyankar /* -------------------------------------------------------------------- 5052b4ed282SShri Abhyankar 5062b4ed282SShri Abhyankar This file implements a semismooth truncated Newton method with a line search, 5072b4ed282SShri Abhyankar for solving a system of nonlinear equations in complementarity form, using the KSP, Vec, 5082b4ed282SShri Abhyankar and Mat interfaces for linear solvers, vectors, and matrices, 5092b4ed282SShri Abhyankar respectively. 5102b4ed282SShri Abhyankar 5112b4ed282SShri Abhyankar The following basic routines are required for each nonlinear solver: 5122b4ed282SShri Abhyankar SNESCreate_XXX() - Creates a nonlinear solver context 5132b4ed282SShri Abhyankar SNESSetFromOptions_XXX() - Sets runtime options 5142b4ed282SShri Abhyankar SNESSolve_XXX() - Solves the nonlinear system 5152b4ed282SShri Abhyankar SNESDestroy_XXX() - Destroys the nonlinear solver context 5162b4ed282SShri Abhyankar The suffix "_XXX" denotes a particular implementation, in this case 5172b4ed282SShri Abhyankar we use _VI (e.g., SNESCreate_VI, SNESSolve_VI) for solving 5182b4ed282SShri Abhyankar systems of nonlinear equations with a line search (LS) method. 5192b4ed282SShri Abhyankar These routines are actually called via the common user interface 5202b4ed282SShri Abhyankar routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 5212b4ed282SShri Abhyankar SNESDestroy(), so the application code interface remains identical 5222b4ed282SShri Abhyankar for all nonlinear solvers. 5232b4ed282SShri Abhyankar 5242b4ed282SShri Abhyankar Another key routine is: 5252b4ed282SShri Abhyankar SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 5262b4ed282SShri Abhyankar by setting data structures and options. The interface routine SNESSetUp() 5272b4ed282SShri Abhyankar is not usually called directly by the user, but instead is called by 5282b4ed282SShri Abhyankar SNESSolve() if necessary. 5292b4ed282SShri Abhyankar 5302b4ed282SShri Abhyankar Additional basic routines are: 5312b4ed282SShri Abhyankar SNESView_XXX() - Prints details of runtime options that 5322b4ed282SShri Abhyankar have actually been used. 5332b4ed282SShri Abhyankar These are called by application codes via the interface routines 5342b4ed282SShri Abhyankar SNESView(). 5352b4ed282SShri Abhyankar 5362b4ed282SShri Abhyankar The various types of solvers (preconditioners, Krylov subspace methods, 5372b4ed282SShri Abhyankar nonlinear solvers, timesteppers) are all organized similarly, so the 5382b4ed282SShri Abhyankar above description applies to these categories also. 5392b4ed282SShri Abhyankar 5402b4ed282SShri Abhyankar -------------------------------------------------------------------- */ 5412b4ed282SShri Abhyankar /* 54269c03620SShri Abhyankar SNESSolveVI_SS - Solves the complementarity problem with a semismooth Newton 5432b4ed282SShri Abhyankar method using a line search. 5442b4ed282SShri Abhyankar 5452b4ed282SShri Abhyankar Input Parameters: 5462b4ed282SShri Abhyankar . snes - the SNES context 5472b4ed282SShri Abhyankar 5482b4ed282SShri Abhyankar Output Parameter: 5492b4ed282SShri Abhyankar . outits - number of iterations until termination 5502b4ed282SShri Abhyankar 5512b4ed282SShri Abhyankar Application Interface Routine: SNESSolve() 5522b4ed282SShri Abhyankar 5532b4ed282SShri Abhyankar Notes: 5542b4ed282SShri Abhyankar This implements essentially a semismooth Newton method with a 5552b4ed282SShri Abhyankar line search. By default a cubic backtracking line search 5562b4ed282SShri Abhyankar is employed, as described in the text "Numerical Methods for 5572b4ed282SShri Abhyankar Unconstrained Optimization and Nonlinear Equations" by Dennis 5582b4ed282SShri Abhyankar and Schnabel. 5592b4ed282SShri Abhyankar */ 5602b4ed282SShri Abhyankar #undef __FUNCT__ 56169c03620SShri Abhyankar #define __FUNCT__ "SNESSolveVI_SS" 56269c03620SShri Abhyankar PetscErrorCode SNESSolveVI_SS(SNES snes) 5632b4ed282SShri Abhyankar { 5642b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 5652b4ed282SShri Abhyankar PetscErrorCode ierr; 5662b4ed282SShri Abhyankar PetscInt maxits,i,lits; 567ace3abfcSBarry Smith PetscBool lssucceed,changedir; 5682b4ed282SShri Abhyankar MatStructure flg = DIFFERENT_NONZERO_PATTERN; 5692b4ed282SShri Abhyankar PetscReal gnorm,xnorm=0,ynorm; 5702b4ed282SShri Abhyankar Vec Y,X,F,G,W; 5712b4ed282SShri Abhyankar KSPConvergedReason kspreason; 5722b4ed282SShri Abhyankar 5732b4ed282SShri Abhyankar PetscFunctionBegin; 5742b4ed282SShri Abhyankar snes->numFailures = 0; 5752b4ed282SShri Abhyankar snes->numLinearSolveFailures = 0; 5762b4ed282SShri Abhyankar snes->reason = SNES_CONVERGED_ITERATING; 5772b4ed282SShri Abhyankar 5782b4ed282SShri Abhyankar maxits = snes->max_its; /* maximum number of iterations */ 5792b4ed282SShri Abhyankar X = snes->vec_sol; /* solution vector */ 5802b4ed282SShri Abhyankar F = snes->vec_func; /* residual vector */ 5812b4ed282SShri Abhyankar Y = snes->work[0]; /* work vectors */ 5822b4ed282SShri Abhyankar G = snes->work[1]; 5832b4ed282SShri Abhyankar W = snes->work[2]; 5842b4ed282SShri Abhyankar 5852b4ed282SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 5862b4ed282SShri Abhyankar snes->iter = 0; 5872b4ed282SShri Abhyankar snes->norm = 0.0; 5882b4ed282SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 5892b4ed282SShri Abhyankar 5902b4ed282SShri Abhyankar ierr = SNESVIAdjustInitialGuess(X,vi->xl,vi->xu);CHKERRQ(ierr); 5912b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr); 5922b4ed282SShri Abhyankar if (snes->domainerror) { 5932b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 5942b4ed282SShri Abhyankar PetscFunctionReturn(0); 5952b4ed282SShri Abhyankar } 5962b4ed282SShri Abhyankar /* Compute Merit function */ 5972b4ed282SShri Abhyankar ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr); 5982b4ed282SShri Abhyankar 5992b4ed282SShri Abhyankar ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 6002b4ed282SShri Abhyankar ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 6012b4ed282SShri Abhyankar if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 6022b4ed282SShri Abhyankar 6032b4ed282SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 6042b4ed282SShri Abhyankar snes->norm = vi->phinorm; 6052b4ed282SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 6062b4ed282SShri Abhyankar SNESLogConvHistory(snes,vi->phinorm,0); 6072b4ed282SShri Abhyankar SNESMonitor(snes,0,vi->phinorm); 6082b4ed282SShri Abhyankar 6092b4ed282SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 6102b4ed282SShri Abhyankar snes->ttol = vi->phinorm*snes->rtol; 6112b4ed282SShri Abhyankar /* test convergence */ 6122b4ed282SShri Abhyankar ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 6132b4ed282SShri Abhyankar if (snes->reason) PetscFunctionReturn(0); 6142b4ed282SShri Abhyankar 6152b4ed282SShri Abhyankar for (i=0; i<maxits; i++) { 6162b4ed282SShri Abhyankar 6172b4ed282SShri Abhyankar /* Call general purpose update function */ 6182b4ed282SShri Abhyankar if (snes->ops->update) { 6192b4ed282SShri Abhyankar ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 6202b4ed282SShri Abhyankar } 6212b4ed282SShri Abhyankar 6222b4ed282SShri Abhyankar /* Solve J Y = Phi, where J is the semismooth jacobian */ 623a79edbefSShri Abhyankar /* Get the nonlinear function jacobian */ 6242b4ed282SShri Abhyankar ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 625a79edbefSShri Abhyankar /* Get the diagonal shift and row scaling vectors */ 626a79edbefSShri Abhyankar ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr); 627a79edbefSShri Abhyankar /* Compute the semismooth jacobian */ 628a79edbefSShri Abhyankar ierr = SNESVIComputeJacobian(snes->jacobian,snes->jacobian_pre,vi->Da,vi->Db);CHKERRQ(ierr); 6292b4ed282SShri Abhyankar 6302b4ed282SShri Abhyankar ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 6312b4ed282SShri Abhyankar ierr = SNES_KSPSolve(snes,snes->ksp,vi->phi,Y);CHKERRQ(ierr); 6322b4ed282SShri Abhyankar ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 6332b4ed282SShri Abhyankar ierr = SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);CHKERRQ(ierr); 6342b4ed282SShri Abhyankar ierr = SNESVICheckDescentDirection(snes,vi->dpsi,Y,&changedir);CHKERRQ(ierr); 6352b4ed282SShri Abhyankar if (kspreason < 0 || changedir) { 6362b4ed282SShri Abhyankar if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 6372b4ed282SShri 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); 6382b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 6392b4ed282SShri Abhyankar break; 6402b4ed282SShri Abhyankar } 6412b4ed282SShri Abhyankar ierr = VecCopy(vi->dpsi,Y);CHKERRQ(ierr); 6422b4ed282SShri Abhyankar } 6432b4ed282SShri Abhyankar ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 6442b4ed282SShri Abhyankar snes->linear_its += lits; 6452b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 6462b4ed282SShri Abhyankar /* 6472b4ed282SShri Abhyankar if (vi->precheckstep) { 648ace3abfcSBarry Smith PetscBool changed_y = PETSC_FALSE; 6492b4ed282SShri Abhyankar ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr); 6502b4ed282SShri Abhyankar } 6512b4ed282SShri Abhyankar 6522b4ed282SShri Abhyankar if (PetscLogPrintInfo){ 6532b4ed282SShri Abhyankar ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 6542b4ed282SShri Abhyankar } 6552b4ed282SShri Abhyankar */ 6562b4ed282SShri Abhyankar /* Compute a (scaled) negative update in the line search routine: 6572b4ed282SShri Abhyankar Y <- X - lambda*Y 6582b4ed282SShri Abhyankar and evaluate G = function(Y) (depends on the line search). 6592b4ed282SShri Abhyankar */ 6602b4ed282SShri Abhyankar ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 6612b4ed282SShri Abhyankar ynorm = 1; gnorm = vi->phinorm; 6622b4ed282SShri Abhyankar ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 6632b4ed282SShri Abhyankar ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 6642b4ed282SShri Abhyankar if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 6652b4ed282SShri Abhyankar if (snes->domainerror) { 6662b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 6672b4ed282SShri Abhyankar PetscFunctionReturn(0); 6682b4ed282SShri Abhyankar } 6692b4ed282SShri Abhyankar if (!lssucceed) { 6702b4ed282SShri Abhyankar if (++snes->numFailures >= snes->maxFailures) { 671ace3abfcSBarry Smith PetscBool ismin; 6722b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_LINE_SEARCH; 6732b4ed282SShri Abhyankar ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 6742b4ed282SShri Abhyankar if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 6752b4ed282SShri Abhyankar break; 6762b4ed282SShri Abhyankar } 6772b4ed282SShri Abhyankar } 6782b4ed282SShri Abhyankar /* Update function and solution vectors */ 6792b4ed282SShri Abhyankar vi->phinorm = gnorm; 6802b4ed282SShri Abhyankar vi->merit = 0.5*vi->phinorm*vi->phinorm; 6812b4ed282SShri Abhyankar ierr = VecCopy(G,vi->phi);CHKERRQ(ierr); 6822b4ed282SShri Abhyankar ierr = VecCopy(W,X);CHKERRQ(ierr); 6832b4ed282SShri Abhyankar /* Monitor convergence */ 6842b4ed282SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 6852b4ed282SShri Abhyankar snes->iter = i+1; 6862b4ed282SShri Abhyankar snes->norm = vi->phinorm; 6872b4ed282SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 6882b4ed282SShri Abhyankar SNESLogConvHistory(snes,snes->norm,lits); 6892b4ed282SShri Abhyankar SNESMonitor(snes,snes->iter,snes->norm); 6902b4ed282SShri Abhyankar /* Test for convergence, xnorm = || X || */ 6912b4ed282SShri Abhyankar if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 6922b4ed282SShri Abhyankar ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 6932b4ed282SShri Abhyankar if (snes->reason) break; 6942b4ed282SShri Abhyankar } 6952b4ed282SShri Abhyankar if (i == maxits) { 6962b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 6972b4ed282SShri Abhyankar if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 6982b4ed282SShri Abhyankar } 6992b4ed282SShri Abhyankar PetscFunctionReturn(0); 7002b4ed282SShri Abhyankar } 70169c03620SShri Abhyankar 70269c03620SShri Abhyankar #undef __FUNCT__ 70369c03620SShri Abhyankar #define __FUNCT__ "SNESVICreateIndexSets_AS" 704a79edbefSShri Abhyankar PetscErrorCode SNESVICreateIndexSets_AS(SNES snes,Vec Db,PetscReal thresh,IS* ISact,IS* ISinact) 70569c03620SShri Abhyankar { 70669c03620SShri Abhyankar PetscErrorCode ierr; 70769c03620SShri Abhyankar PetscInt i,nlocal,ilow,ihigh,nloc_isact=0,nloc_isinact=0; 70869c03620SShri Abhyankar PetscInt *idx_act,*idx_inact,i1=0,i2=0; 70969c03620SShri Abhyankar PetscScalar *db; 71069c03620SShri Abhyankar 71169c03620SShri Abhyankar PetscFunctionBegin; 71269c03620SShri Abhyankar 71369c03620SShri Abhyankar ierr = VecGetLocalSize(Db,&nlocal);CHKERRQ(ierr); 71469c03620SShri Abhyankar ierr = VecGetOwnershipRange(Db,&ilow,&ihigh);CHKERRQ(ierr); 71569c03620SShri Abhyankar ierr = VecGetArray(Db,&db);CHKERRQ(ierr); 71669c03620SShri Abhyankar /* Compute the sizes of the active and inactive sets */ 71769c03620SShri Abhyankar for (i=0; i < nlocal;i++) { 718fca35906SShri Abhyankar if (PetscAbsScalar(db[i]) <= thresh) nloc_isact++; 719fca35906SShri Abhyankar else nloc_isinact++; 72069c03620SShri Abhyankar } 72169c03620SShri Abhyankar ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr); 72269c03620SShri Abhyankar ierr = PetscMalloc(nloc_isinact*sizeof(PetscInt),&idx_inact);CHKERRQ(ierr); 72369c03620SShri Abhyankar 72469c03620SShri Abhyankar /* Creating the indexing arrays */ 725fca35906SShri Abhyankar for(i=0; i < nlocal; i++) { 726fca35906SShri Abhyankar if (PetscAbsScalar(db[i]) <= thresh) idx_act[i1++] = ilow+i; 727fca35906SShri Abhyankar else idx_inact[i2++] = ilow+i; 72869c03620SShri Abhyankar } 72969c03620SShri Abhyankar 73069c03620SShri Abhyankar /* Create the index sets */ 73169c03620SShri Abhyankar ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isact,idx_act,PETSC_COPY_VALUES,ISact);CHKERRQ(ierr); 73269c03620SShri Abhyankar ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isinact,idx_inact,PETSC_COPY_VALUES,ISinact);CHKERRQ(ierr); 73369c03620SShri Abhyankar 73469c03620SShri Abhyankar ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr); 73569c03620SShri Abhyankar ierr = PetscFree(idx_act);CHKERRQ(ierr); 73669c03620SShri Abhyankar ierr = PetscFree(idx_inact);CHKERRQ(ierr); 73769c03620SShri Abhyankar PetscFunctionReturn(0); 73869c03620SShri Abhyankar } 73969c03620SShri Abhyankar 740dbd914b8SShri Abhyankar /* Create active and inactive set vectors. The local size of this vector is set and petsc computes the global size */ 741dbd914b8SShri Abhyankar #undef __FUNCT__ 742dbd914b8SShri Abhyankar #define __FUNCT__ "SNESVICreateVectors_AS" 743dbd914b8SShri Abhyankar PetscErrorCode SNESVICreateVectors_AS(SNES snes,PetscInt n,Vec* newv) 744dbd914b8SShri Abhyankar { 745dbd914b8SShri Abhyankar PetscErrorCode ierr; 746dbd914b8SShri Abhyankar Vec v; 747dbd914b8SShri Abhyankar 748dbd914b8SShri Abhyankar PetscFunctionBegin; 749dbd914b8SShri Abhyankar ierr = VecCreate(((PetscObject)snes)->comm,&v);CHKERRQ(ierr); 750dbd914b8SShri Abhyankar ierr = VecSetSizes(v,n,PETSC_DECIDE);CHKERRQ(ierr); 751dbd914b8SShri Abhyankar ierr = VecSetFromOptions(v);CHKERRQ(ierr); 752dbd914b8SShri Abhyankar *newv = v; 753dbd914b8SShri Abhyankar 754dbd914b8SShri Abhyankar PetscFunctionReturn(0); 755dbd914b8SShri Abhyankar } 756dbd914b8SShri Abhyankar 757dbd914b8SShri Abhyankar 75869c03620SShri Abhyankar /* Variational Inequality solver using active set method */ 75969c03620SShri Abhyankar #undef __FUNCT__ 76069c03620SShri Abhyankar #define __FUNCT__ "SNESSolveVI_AS" 76169c03620SShri Abhyankar PetscErrorCode SNESSolveVI_AS(SNES snes) 76269c03620SShri Abhyankar { 76369c03620SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 76469c03620SShri Abhyankar PetscErrorCode ierr; 76569c03620SShri Abhyankar PetscInt maxits,i,lits; 76669c03620SShri Abhyankar PetscBool lssucceed,changedir; 76769c03620SShri Abhyankar MatStructure flg = DIFFERENT_NONZERO_PATTERN; 76869c03620SShri Abhyankar PetscReal gnorm,xnorm=0,ynorm; 76969c03620SShri Abhyankar Vec Y,X,F,G,W; 77069c03620SShri Abhyankar KSPConvergedReason kspreason; 77169c03620SShri Abhyankar 77269c03620SShri Abhyankar PetscFunctionBegin; 77369c03620SShri Abhyankar snes->numFailures = 0; 77469c03620SShri Abhyankar snes->numLinearSolveFailures = 0; 77569c03620SShri Abhyankar snes->reason = SNES_CONVERGED_ITERATING; 77669c03620SShri Abhyankar 77769c03620SShri Abhyankar maxits = snes->max_its; /* maximum number of iterations */ 77869c03620SShri Abhyankar X = snes->vec_sol; /* solution vector */ 77969c03620SShri Abhyankar F = snes->vec_func; /* residual vector */ 78069c03620SShri Abhyankar Y = snes->work[0]; /* work vectors */ 78169c03620SShri Abhyankar G = snes->work[1]; 78269c03620SShri Abhyankar W = snes->work[2]; 78369c03620SShri Abhyankar 78469c03620SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 78569c03620SShri Abhyankar snes->iter = 0; 78669c03620SShri Abhyankar snes->norm = 0.0; 78769c03620SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 78869c03620SShri Abhyankar 78969c03620SShri Abhyankar ierr = SNESVIAdjustInitialGuess(X,vi->xl,vi->xu);CHKERRQ(ierr); 79069c03620SShri Abhyankar ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr); 79169c03620SShri Abhyankar if (snes->domainerror) { 79269c03620SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 79369c03620SShri Abhyankar PetscFunctionReturn(0); 79469c03620SShri Abhyankar } 79569c03620SShri Abhyankar /* Compute Merit function */ 79669c03620SShri Abhyankar ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr); 79769c03620SShri Abhyankar 79869c03620SShri Abhyankar ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 79969c03620SShri Abhyankar ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 80069c03620SShri Abhyankar if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 80169c03620SShri Abhyankar 80269c03620SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 80369c03620SShri Abhyankar snes->norm = vi->phinorm; 80469c03620SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 80569c03620SShri Abhyankar SNESLogConvHistory(snes,vi->phinorm,0); 80669c03620SShri Abhyankar SNESMonitor(snes,0,vi->phinorm); 80769c03620SShri Abhyankar 80869c03620SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 80969c03620SShri Abhyankar snes->ttol = vi->phinorm*snes->rtol; 81069c03620SShri Abhyankar /* test convergence */ 81169c03620SShri Abhyankar ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 81269c03620SShri Abhyankar if (snes->reason) PetscFunctionReturn(0); 81369c03620SShri Abhyankar 81469c03620SShri Abhyankar for (i=0; i<maxits; i++) { 81569c03620SShri Abhyankar 816a79edbefSShri Abhyankar IS IS_act,IS_inact; /* _act -> active set _inact -> inactive set */ 817a79edbefSShri Abhyankar PetscReal thresh,J_norm1; 818dbd914b8SShri Abhyankar VecScatter scat_act,scat_inact; 819dbd914b8SShri Abhyankar PetscInt nis_act,nis_inact; 820dbd914b8SShri Abhyankar Vec Da_act,Da_inact,Db_inact; 821dbd914b8SShri Abhyankar Vec Y_act,Y_inact,phi_act,phi_inact; 822dbd914b8SShri Abhyankar Mat jac_inact_inact,jac_inact_act; 823a79edbefSShri Abhyankar 82469c03620SShri Abhyankar /* Call general purpose update function */ 82569c03620SShri Abhyankar if (snes->ops->update) { 82669c03620SShri Abhyankar ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 82769c03620SShri Abhyankar } 82869c03620SShri Abhyankar ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 82969c03620SShri Abhyankar /* Compute the threshold value for creating active and inactive sets */ 83069c03620SShri Abhyankar ierr = MatNorm(snes->jacobian,NORM_1,&J_norm1);CHKERRQ(ierr); 83169c03620SShri Abhyankar thresh = PetscMin(vi->merit,1e-2)/(1+J_norm1); 832a79edbefSShri Abhyankar 833a79edbefSShri Abhyankar /* Compute B-subdifferential vectors Da and Db */ 834a79edbefSShri Abhyankar ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr); 835dbd914b8SShri Abhyankar 83669c03620SShri Abhyankar /* Create active and inactive index sets */ 83769c03620SShri Abhyankar ierr = SNESVICreateIndexSets_AS(snes,vi->Db,thresh,&IS_act,&IS_inact);CHKERRQ(ierr); 83869c03620SShri Abhyankar 839dbd914b8SShri Abhyankar /* Get local sizes of active and inactive sets */ 840dbd914b8SShri Abhyankar ierr = ISGetLocalSize(IS_act,&nis_act);CHKERRQ(ierr); 841dbd914b8SShri Abhyankar ierr = ISGetLocalSize(IS_inact,&nis_inact);CHKERRQ(ierr); 842dbd914b8SShri Abhyankar 843dbd914b8SShri Abhyankar /* Create active and inactive set vectors */ 844dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_act,&Da_act);CHKERRQ(ierr); 845dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_inact,&Da_inact);CHKERRQ(ierr); 846dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_inact,&Db_inact);CHKERRQ(ierr); 847dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_act,&phi_act);CHKERRQ(ierr); 848dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_inact,&phi_inact);CHKERRQ(ierr); 849dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_act,&Y_act);CHKERRQ(ierr); 850dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_inact,&Y_inact);CHKERRQ(ierr); 851dbd914b8SShri Abhyankar 852dbd914b8SShri Abhyankar /* Create inactive set submatrices */ 853dbd914b8SShri Abhyankar ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_act,MAT_INITIAL_MATRIX,&jac_inact_act);CHKERRQ(ierr); 854dbd914b8SShri Abhyankar ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);CHKERRQ(ierr); 855dbd914b8SShri Abhyankar 856dbd914b8SShri Abhyankar /* Create scatter contexts */ 857dbd914b8SShri Abhyankar ierr = VecScatterCreate(vi->Da,IS_act,Da_act,PETSC_NULL,&scat_act);CHKERRQ(ierr); 858dbd914b8SShri Abhyankar ierr = VecScatterCreate(vi->Da,IS_inact,Da_inact,PETSC_NULL,&scat_inact);CHKERRQ(ierr); 859dbd914b8SShri Abhyankar 860dbd914b8SShri Abhyankar /* Do a vec scatter to active and inactive set vectors */ 861dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_act,vi->Da,Da_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 862dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_act,vi->Da,Da_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 863dbd914b8SShri Abhyankar 864dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_inact,vi->Da,Da_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 865dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_inact,vi->Da,Da_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 866dbd914b8SShri Abhyankar 867dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_inact,vi->Db,Db_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 868dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_inact,vi->Db,Db_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 869dbd914b8SShri Abhyankar 870dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_act,vi->phi,phi_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 871dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_act,vi->phi,phi_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 872dbd914b8SShri Abhyankar 873dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_inact,vi->phi,phi_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 874dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_inact,vi->phi,phi_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 875dbd914b8SShri Abhyankar 876dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 877dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 878dbd914b8SShri Abhyankar 879dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 880dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 881dbd914b8SShri Abhyankar 882dbd914b8SShri Abhyankar /* Active set direction */ 883dbd914b8SShri Abhyankar ierr = VecPointwiseDivide(Y_act,phi_act,Da_act);CHKERRQ(ierr); 884dbd914b8SShri Abhyankar /* inactive set jacobian */ 885dbd914b8SShri Abhyankar ierr = VecPointwiseDivide(Da_inact,Da_inact,Db_inact);CHKERRQ(ierr); 886dbd914b8SShri Abhyankar ierr = MatDiagonalSet(jac_inact_inact,Da_inact,ADD_VALUES);CHKERRQ(ierr); 887dbd914b8SShri Abhyankar /* right hand side */ 888dbd914b8SShri Abhyankar ierr = VecPointwiseDivide(phi_inact,phi_inact,Db_inact);CHKERRQ(ierr); 889dbd914b8SShri Abhyankar ierr = MatMult(jac_inact_act,Y_act,Db_inact);CHKERRQ(ierr); 890dbd914b8SShri Abhyankar ierr = VecAXPY(phi_inact,-1.0,Db_inact);CHKERRQ(ierr); 891dbd914b8SShri Abhyankar 892dbd914b8SShri Abhyankar /* USING THE SAME MATRIX AS PRECONDITIONER....NEED TO CHANGE THIS */ 893dbd914b8SShri Abhyankar ierr = KSPSetOperators(snes->ksp,jac_inact_inact,jac_inact_inact,flg);CHKERRQ(ierr); 894dbd914b8SShri Abhyankar ierr = SNES_KSPSolve(snes,snes->ksp,phi_inact,Y_inact);CHKERRQ(ierr); 895dbd914b8SShri Abhyankar ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 896dbd914b8SShri Abhyankar /* Compute the jacobian of the semismooth function which is needed for calculating the merit function 897dbd914b8SShri Abhyankar gradient */ 898dbd914b8SShri Abhyankar ierr = SNESVIComputeJacobian(snes->jacobian,snes->jacobian_pre,vi->Da,vi->Db);CHKERRQ(ierr); 899dbd914b8SShri Abhyankar ierr = SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);CHKERRQ(ierr); 900dbd914b8SShri Abhyankar 901dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 902dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 903dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 904dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 905dbd914b8SShri Abhyankar 906dbd914b8SShri Abhyankar ierr = VecDestroy(Da_act);CHKERRQ(ierr); 907dbd914b8SShri Abhyankar ierr = VecDestroy(Da_inact);CHKERRQ(ierr); 908dbd914b8SShri Abhyankar ierr = VecDestroy(Db_inact);CHKERRQ(ierr); 909dbd914b8SShri Abhyankar ierr = VecDestroy(phi_act);CHKERRQ(ierr); 910dbd914b8SShri Abhyankar ierr = VecDestroy(phi_inact);CHKERRQ(ierr); 911dbd914b8SShri Abhyankar ierr = VecDestroy(Y_act);CHKERRQ(ierr); 912dbd914b8SShri Abhyankar ierr = VecDestroy(Y_inact);CHKERRQ(ierr); 913dbd914b8SShri Abhyankar ierr = VecScatterDestroy(scat_act);CHKERRQ(ierr); 914dbd914b8SShri Abhyankar ierr = VecScatterDestroy(scat_inact);CHKERRQ(ierr); 915730c24dcSShri Abhyankar ierr = ISDestroy(IS_act);CHKERRQ(ierr); 916730c24dcSShri Abhyankar ierr = ISDestroy(IS_inact);CHKERRQ(ierr); 917*bac2f86eSShri Abhyankar ierr = MatDestroy(jac_inact_act);CHKERRQ(ierr); 918*bac2f86eSShri Abhyankar ierr = MatDestroy(jac_inact_inact);CHKERRQ(ierr); 919730c24dcSShri Abhyankar 92069c03620SShri Abhyankar ierr = SNESVICheckDescentDirection(snes,vi->dpsi,Y,&changedir);CHKERRQ(ierr); 92169c03620SShri Abhyankar if (kspreason < 0 || changedir) { 92269c03620SShri Abhyankar if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 92369c03620SShri 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); 92469c03620SShri Abhyankar snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 92569c03620SShri Abhyankar break; 92669c03620SShri Abhyankar } 92769c03620SShri Abhyankar ierr = VecCopy(vi->dpsi,Y);CHKERRQ(ierr); 92869c03620SShri Abhyankar } 92969c03620SShri Abhyankar ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 93069c03620SShri Abhyankar snes->linear_its += lits; 93169c03620SShri Abhyankar ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 93269c03620SShri Abhyankar /* 93369c03620SShri Abhyankar if (vi->precheckstep) { 93469c03620SShri Abhyankar PetscBool changed_y = PETSC_FALSE; 93569c03620SShri Abhyankar ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr); 93669c03620SShri Abhyankar } 93769c03620SShri Abhyankar 93869c03620SShri Abhyankar if (PetscLogPrintInfo){ 93969c03620SShri Abhyankar ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 94069c03620SShri Abhyankar } 94169c03620SShri Abhyankar */ 94269c03620SShri Abhyankar /* Compute a (scaled) negative update in the line search routine: 94369c03620SShri Abhyankar Y <- X - lambda*Y 94469c03620SShri Abhyankar and evaluate G = function(Y) (depends on the line search). 94569c03620SShri Abhyankar */ 94669c03620SShri Abhyankar ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 94769c03620SShri Abhyankar ynorm = 1; gnorm = vi->phinorm; 94869c03620SShri Abhyankar ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 94969c03620SShri Abhyankar ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 95069c03620SShri Abhyankar if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 95169c03620SShri Abhyankar if (snes->domainerror) { 95269c03620SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 95369c03620SShri Abhyankar PetscFunctionReturn(0); 95469c03620SShri Abhyankar } 95569c03620SShri Abhyankar if (!lssucceed) { 95669c03620SShri Abhyankar if (++snes->numFailures >= snes->maxFailures) { 95769c03620SShri Abhyankar PetscBool ismin; 95869c03620SShri Abhyankar snes->reason = SNES_DIVERGED_LINE_SEARCH; 95969c03620SShri Abhyankar ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 96069c03620SShri Abhyankar if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 96169c03620SShri Abhyankar break; 96269c03620SShri Abhyankar } 96369c03620SShri Abhyankar } 96469c03620SShri Abhyankar /* Update function and solution vectors */ 96569c03620SShri Abhyankar vi->phinorm = gnorm; 96669c03620SShri Abhyankar vi->merit = 0.5*vi->phinorm*vi->phinorm; 96769c03620SShri Abhyankar ierr = VecCopy(G,vi->phi);CHKERRQ(ierr); 96869c03620SShri Abhyankar ierr = VecCopy(W,X);CHKERRQ(ierr); 96969c03620SShri Abhyankar /* Monitor convergence */ 97069c03620SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 97169c03620SShri Abhyankar snes->iter = i+1; 97269c03620SShri Abhyankar snes->norm = vi->phinorm; 97369c03620SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 97469c03620SShri Abhyankar SNESLogConvHistory(snes,snes->norm,lits); 97569c03620SShri Abhyankar SNESMonitor(snes,snes->iter,snes->norm); 97669c03620SShri Abhyankar /* Test for convergence, xnorm = || X || */ 97769c03620SShri Abhyankar if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 97869c03620SShri Abhyankar ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 97969c03620SShri Abhyankar if (snes->reason) break; 98069c03620SShri Abhyankar } 98169c03620SShri Abhyankar if (i == maxits) { 98269c03620SShri Abhyankar ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 98369c03620SShri Abhyankar if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 98469c03620SShri Abhyankar } 98569c03620SShri Abhyankar PetscFunctionReturn(0); 98669c03620SShri Abhyankar } 98769c03620SShri Abhyankar 9882b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 9892b4ed282SShri Abhyankar /* 9902b4ed282SShri Abhyankar SNESSetUp_VI - Sets up the internal data structures for the later use 9912b4ed282SShri Abhyankar of the SNESVI nonlinear solver. 9922b4ed282SShri Abhyankar 9932b4ed282SShri Abhyankar Input Parameter: 9942b4ed282SShri Abhyankar . snes - the SNES context 9952b4ed282SShri Abhyankar . x - the solution vector 9962b4ed282SShri Abhyankar 9972b4ed282SShri Abhyankar Application Interface Routine: SNESSetUp() 9982b4ed282SShri Abhyankar 9992b4ed282SShri Abhyankar Notes: 10002b4ed282SShri Abhyankar For basic use of the SNES solvers, the user need not explicitly call 10012b4ed282SShri Abhyankar SNESSetUp(), since these actions will automatically occur during 10022b4ed282SShri Abhyankar the call to SNESSolve(). 10032b4ed282SShri Abhyankar */ 10042b4ed282SShri Abhyankar #undef __FUNCT__ 10052b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetUp_VI" 10062b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes) 10072b4ed282SShri Abhyankar { 10082b4ed282SShri Abhyankar PetscErrorCode ierr; 10092b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*) snes->data; 10102b4ed282SShri Abhyankar PetscInt i_start[3],i_end[3]; 10112b4ed282SShri Abhyankar 10122b4ed282SShri Abhyankar PetscFunctionBegin; 10132b4ed282SShri Abhyankar if (!snes->vec_sol_update) { 10142b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr); 10152b4ed282SShri Abhyankar ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr); 10162b4ed282SShri Abhyankar } 10172b4ed282SShri Abhyankar if (!snes->work) { 10182b4ed282SShri Abhyankar snes->nwork = 3; 10192b4ed282SShri Abhyankar ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 10202b4ed282SShri Abhyankar ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr); 10212b4ed282SShri Abhyankar } 10222b4ed282SShri Abhyankar 10232b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->phi); CHKERRQ(ierr); 10242b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->dpsi); CHKERRQ(ierr); 10252b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->Da); CHKERRQ(ierr); 10262b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->Db); CHKERRQ(ierr); 10272b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->z);CHKERRQ(ierr); 10282b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->t); CHKERRQ(ierr); 10292b4ed282SShri Abhyankar 10302b4ed282SShri Abhyankar /* If the lower and upper bound on variables are not set, set it to 10312b4ed282SShri Abhyankar -Inf and Inf */ 10322b4ed282SShri Abhyankar if (!vi->xl && !vi->xu) { 10332b4ed282SShri Abhyankar vi->usersetxbounds = PETSC_FALSE; 10342b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->xl); CHKERRQ(ierr); 10352b4ed282SShri Abhyankar ierr = VecSet(vi->xl,PETSC_VI_NINF);CHKERRQ(ierr); 10362b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->xu); CHKERRQ(ierr); 10372b4ed282SShri Abhyankar ierr = VecSet(vi->xu,PETSC_VI_INF);CHKERRQ(ierr); 10382b4ed282SShri Abhyankar } else { 10392b4ed282SShri Abhyankar /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */ 10402b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr); 10412b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(vi->xl,i_start+1,i_end+1);CHKERRQ(ierr); 10422b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(vi->xu,i_start+2,i_end+2);CHKERRQ(ierr); 10432b4ed282SShri 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])) 10442b4ed282SShri 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."); 10452b4ed282SShri Abhyankar } 10462b4ed282SShri Abhyankar 10472b4ed282SShri Abhyankar vi->computeuserfunction = snes->ops->computefunction; 10481f79c099SShri Abhyankar snes->ops->computefunction = SNESVIComputeFunction; 1049a79edbefSShri Abhyankar 10502b4ed282SShri Abhyankar PetscFunctionReturn(0); 10512b4ed282SShri Abhyankar } 10522b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 10532b4ed282SShri Abhyankar /* 10542b4ed282SShri Abhyankar SNESDestroy_VI - Destroys the private SNES_VI context that was created 10552b4ed282SShri Abhyankar with SNESCreate_VI(). 10562b4ed282SShri Abhyankar 10572b4ed282SShri Abhyankar Input Parameter: 10582b4ed282SShri Abhyankar . snes - the SNES context 10592b4ed282SShri Abhyankar 10602b4ed282SShri Abhyankar Application Interface Routine: SNESDestroy() 10612b4ed282SShri Abhyankar */ 10622b4ed282SShri Abhyankar #undef __FUNCT__ 10632b4ed282SShri Abhyankar #define __FUNCT__ "SNESDestroy_VI" 10642b4ed282SShri Abhyankar PetscErrorCode SNESDestroy_VI(SNES snes) 10652b4ed282SShri Abhyankar { 10662b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*) snes->data; 10672b4ed282SShri Abhyankar PetscErrorCode ierr; 10682b4ed282SShri Abhyankar 10692b4ed282SShri Abhyankar PetscFunctionBegin; 10702b4ed282SShri Abhyankar if (snes->vec_sol_update) { 10712b4ed282SShri Abhyankar ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr); 10722b4ed282SShri Abhyankar snes->vec_sol_update = PETSC_NULL; 10732b4ed282SShri Abhyankar } 10742b4ed282SShri Abhyankar if (snes->nwork) { 10752b4ed282SShri Abhyankar ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 10762b4ed282SShri Abhyankar snes->nwork = 0; 10772b4ed282SShri Abhyankar snes->work = PETSC_NULL; 10782b4ed282SShri Abhyankar } 10792b4ed282SShri Abhyankar 10802b4ed282SShri Abhyankar /* clear vectors */ 10812b4ed282SShri Abhyankar ierr = VecDestroy(vi->phi); CHKERRQ(ierr); 10822b4ed282SShri Abhyankar ierr = VecDestroy(vi->dpsi); CHKERRQ(ierr); 10832b4ed282SShri Abhyankar ierr = VecDestroy(vi->Da); CHKERRQ(ierr); 10842b4ed282SShri Abhyankar ierr = VecDestroy(vi->Db); CHKERRQ(ierr); 10852b4ed282SShri Abhyankar ierr = VecDestroy(vi->z); CHKERRQ(ierr); 10862b4ed282SShri Abhyankar ierr = VecDestroy(vi->t); CHKERRQ(ierr); 10872b4ed282SShri Abhyankar if (!vi->usersetxbounds) { 10882b4ed282SShri Abhyankar ierr = VecDestroy(vi->xl); CHKERRQ(ierr); 10892b4ed282SShri Abhyankar ierr = VecDestroy(vi->xu); CHKERRQ(ierr); 10902b4ed282SShri Abhyankar } 1091c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1092c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr); 10932b4ed282SShri Abhyankar } 10942b4ed282SShri Abhyankar ierr = PetscFree(snes->data);CHKERRQ(ierr); 10952b4ed282SShri Abhyankar 10962b4ed282SShri Abhyankar /* clear composed functions */ 10972b4ed282SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr); 1098c92abb8aSShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr); 10992b4ed282SShri Abhyankar PetscFunctionReturn(0); 11002b4ed282SShri Abhyankar } 11012b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 11022b4ed282SShri Abhyankar #undef __FUNCT__ 11032b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNo_VI" 11042b4ed282SShri Abhyankar 11052b4ed282SShri Abhyankar /* 11062b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchNo routine in snes/impls/ls/ls.c 11072b4ed282SShri Abhyankar 11082b4ed282SShri Abhyankar */ 1109ace3abfcSBarry Smith 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,PetscBool *flag) 11102b4ed282SShri Abhyankar { 11112b4ed282SShri Abhyankar PetscErrorCode ierr; 11122b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1113ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 11142b4ed282SShri Abhyankar 11152b4ed282SShri Abhyankar PetscFunctionBegin; 11162b4ed282SShri Abhyankar *flag = PETSC_TRUE; 11172b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 11182b4ed282SShri Abhyankar ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 11192b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 11202b4ed282SShri Abhyankar if (vi->postcheckstep) { 11212b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 11222b4ed282SShri Abhyankar } 11232b4ed282SShri Abhyankar if (changed_y) { 11242b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 11252b4ed282SShri Abhyankar } 11262b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 11272b4ed282SShri Abhyankar if (!snes->domainerror) { 11282b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 11292b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 11302b4ed282SShri Abhyankar } 1131c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1132c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 1133c92abb8aSShri Abhyankar } 11342b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 11352b4ed282SShri Abhyankar PetscFunctionReturn(0); 11362b4ed282SShri Abhyankar } 11372b4ed282SShri Abhyankar 11382b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 11392b4ed282SShri Abhyankar #undef __FUNCT__ 11402b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNoNorms_VI" 11412b4ed282SShri Abhyankar 11422b4ed282SShri Abhyankar /* 11432b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c 11442b4ed282SShri Abhyankar */ 1145ace3abfcSBarry Smith 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,PetscBool *flag) 11462b4ed282SShri Abhyankar { 11472b4ed282SShri Abhyankar PetscErrorCode ierr; 11482b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1149ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 11502b4ed282SShri Abhyankar 11512b4ed282SShri Abhyankar PetscFunctionBegin; 11522b4ed282SShri Abhyankar *flag = PETSC_TRUE; 11532b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 11542b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 11552b4ed282SShri Abhyankar if (vi->postcheckstep) { 11562b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 11572b4ed282SShri Abhyankar } 11582b4ed282SShri Abhyankar if (changed_y) { 11592b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 11602b4ed282SShri Abhyankar } 11612b4ed282SShri Abhyankar 11622b4ed282SShri Abhyankar /* don't evaluate function the last time through */ 11632b4ed282SShri Abhyankar if (snes->iter < snes->max_its-1) { 11642b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 11652b4ed282SShri Abhyankar } 11662b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 11672b4ed282SShri Abhyankar PetscFunctionReturn(0); 11682b4ed282SShri Abhyankar } 11692b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 11702b4ed282SShri Abhyankar #undef __FUNCT__ 11712b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchCubic_VI" 11722b4ed282SShri Abhyankar /* 11732b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchCubic in snes/impls/ls/ls.c 11742b4ed282SShri Abhyankar */ 1175ace3abfcSBarry Smith 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,PetscBool *flag) 11762b4ed282SShri Abhyankar { 11772b4ed282SShri Abhyankar /* 11782b4ed282SShri Abhyankar Note that for line search purposes we work with with the related 11792b4ed282SShri Abhyankar minimization problem: 11802b4ed282SShri Abhyankar min z(x): R^n -> R, 11812b4ed282SShri Abhyankar where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 11822b4ed282SShri Abhyankar This function z(x) is same as the merit function for the complementarity problem. 11832b4ed282SShri Abhyankar */ 11842b4ed282SShri Abhyankar 11852b4ed282SShri Abhyankar PetscReal initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 11862b4ed282SShri Abhyankar PetscReal minlambda,lambda,lambdatemp; 11872b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 11882b4ed282SShri Abhyankar PetscScalar cinitslope; 11892b4ed282SShri Abhyankar #endif 11902b4ed282SShri Abhyankar PetscErrorCode ierr; 11912b4ed282SShri Abhyankar PetscInt count; 11922b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1193ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 11942b4ed282SShri Abhyankar MPI_Comm comm; 11952b4ed282SShri Abhyankar 11962b4ed282SShri Abhyankar PetscFunctionBegin; 11972b4ed282SShri Abhyankar ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 11982b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 11992b4ed282SShri Abhyankar *flag = PETSC_TRUE; 12002b4ed282SShri Abhyankar 12012b4ed282SShri Abhyankar ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 12022b4ed282SShri Abhyankar if (*ynorm == 0.0) { 1203c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1204c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 12052b4ed282SShri Abhyankar } 12062b4ed282SShri Abhyankar *gnorm = fnorm; 12072b4ed282SShri Abhyankar ierr = VecCopy(x,w);CHKERRQ(ierr); 12082b4ed282SShri Abhyankar ierr = VecCopy(f,g);CHKERRQ(ierr); 12092b4ed282SShri Abhyankar *flag = PETSC_FALSE; 12102b4ed282SShri Abhyankar goto theend1; 12112b4ed282SShri Abhyankar } 12122b4ed282SShri Abhyankar if (*ynorm > vi->maxstep) { /* Step too big, so scale back */ 1213c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1214c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Scaling step by %G old ynorm %G\n",vi->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr); 12152b4ed282SShri Abhyankar } 12162b4ed282SShri Abhyankar ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr); 12172b4ed282SShri Abhyankar *ynorm = vi->maxstep; 12182b4ed282SShri Abhyankar } 12192b4ed282SShri Abhyankar ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 12202b4ed282SShri Abhyankar minlambda = vi->minlambda/rellength; 12212b4ed282SShri Abhyankar /* ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); */ 12222b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 12232b4ed282SShri Abhyankar ierr = VecDot(vi->dpsi,y,&cinitslope);CHKERRQ(ierr); 12242b4ed282SShri Abhyankar initslope = PetscRealPart(cinitslope); 12252b4ed282SShri Abhyankar #else 12262b4ed282SShri Abhyankar ierr = VecDot(vi->dpsi,y,&initslope);CHKERRQ(ierr); 12272b4ed282SShri Abhyankar #endif 12282b4ed282SShri Abhyankar if (initslope > 0.0) initslope = -initslope; 12292b4ed282SShri Abhyankar if (initslope == 0.0) initslope = -1.0; 12302b4ed282SShri Abhyankar 12312b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 12322b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 12332b4ed282SShri Abhyankar ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 12342b4ed282SShri Abhyankar *flag = PETSC_FALSE; 12352b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 12362b4ed282SShri Abhyankar goto theend1; 12372b4ed282SShri Abhyankar } 12382b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 12392b4ed282SShri Abhyankar if (snes->domainerror) { 12402b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 12412b4ed282SShri Abhyankar PetscFunctionReturn(0); 12422b4ed282SShri Abhyankar } 12432b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 12442b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 12452b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 12462b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */ 1247c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1248c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 12492b4ed282SShri Abhyankar } 12502b4ed282SShri Abhyankar goto theend1; 12512b4ed282SShri Abhyankar } 12522b4ed282SShri Abhyankar 12532b4ed282SShri Abhyankar /* Fit points with quadratic */ 12542b4ed282SShri Abhyankar lambda = 1.0; 12552b4ed282SShri Abhyankar lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 12562b4ed282SShri Abhyankar lambdaprev = lambda; 12572b4ed282SShri Abhyankar gnormprev = *gnorm; 12582b4ed282SShri Abhyankar if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 12592b4ed282SShri Abhyankar if (lambdatemp <= .1*lambda) lambda = .1*lambda; 12602b4ed282SShri Abhyankar else lambda = lambdatemp; 12612b4ed282SShri Abhyankar 12622b4ed282SShri Abhyankar ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 12632b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 12642b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 12652b4ed282SShri Abhyankar *flag = PETSC_FALSE; 12662b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 12672b4ed282SShri Abhyankar goto theend1; 12682b4ed282SShri Abhyankar } 12692b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 12702b4ed282SShri Abhyankar if (snes->domainerror) { 12712b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 12722b4ed282SShri Abhyankar PetscFunctionReturn(0); 12732b4ed282SShri Abhyankar } 12742b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 12752b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1276c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1277c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr); 12782b4ed282SShri Abhyankar } 12792b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */ 1280c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1281c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 12822b4ed282SShri Abhyankar } 12832b4ed282SShri Abhyankar goto theend1; 12842b4ed282SShri Abhyankar } 12852b4ed282SShri Abhyankar 12862b4ed282SShri Abhyankar /* Fit points with cubic */ 12872b4ed282SShri Abhyankar count = 1; 12882b4ed282SShri Abhyankar while (PETSC_TRUE) { 12892b4ed282SShri Abhyankar if (lambda <= minlambda) { 1290c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1291c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 1292c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," 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); 12932b4ed282SShri Abhyankar } 12942b4ed282SShri Abhyankar *flag = PETSC_FALSE; 12952b4ed282SShri Abhyankar break; 12962b4ed282SShri Abhyankar } 12972b4ed282SShri Abhyankar t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 12982b4ed282SShri Abhyankar t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 12992b4ed282SShri Abhyankar a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 13002b4ed282SShri Abhyankar b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 13012b4ed282SShri Abhyankar d = b*b - 3*a*initslope; 13022b4ed282SShri Abhyankar if (d < 0.0) d = 0.0; 13032b4ed282SShri Abhyankar if (a == 0.0) { 13042b4ed282SShri Abhyankar lambdatemp = -initslope/(2.0*b); 13052b4ed282SShri Abhyankar } else { 13062b4ed282SShri Abhyankar lambdatemp = (-b + sqrt(d))/(3.0*a); 13072b4ed282SShri Abhyankar } 13082b4ed282SShri Abhyankar lambdaprev = lambda; 13092b4ed282SShri Abhyankar gnormprev = *gnorm; 13102b4ed282SShri Abhyankar if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 13112b4ed282SShri Abhyankar if (lambdatemp <= .1*lambda) lambda = .1*lambda; 13122b4ed282SShri Abhyankar else lambda = lambdatemp; 13132b4ed282SShri Abhyankar ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 13142b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 13152b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 13162b4ed282SShri 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); 13172b4ed282SShri Abhyankar *flag = PETSC_FALSE; 13182b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 13192b4ed282SShri Abhyankar break; 13202b4ed282SShri Abhyankar } 13212b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 13222b4ed282SShri Abhyankar if (snes->domainerror) { 13232b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 13242b4ed282SShri Abhyankar PetscFunctionReturn(0); 13252b4ed282SShri Abhyankar } 13262b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 13272b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 13282b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* is reduction enough? */ 1329c92abb8aSShri Abhyankar if (vi->lsmonitor) { 13302b4ed282SShri Abhyankar ierr = PetscPrintf(comm," Line search: Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 13312b4ed282SShri Abhyankar } 13322b4ed282SShri Abhyankar break; 13332b4ed282SShri Abhyankar } else { 1334c92abb8aSShri Abhyankar if (vi->lsmonitor) { 13352b4ed282SShri Abhyankar ierr = PetscPrintf(comm," Line search: Cubic step no good, shrinking lambda, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 13362b4ed282SShri Abhyankar } 13372b4ed282SShri Abhyankar } 13382b4ed282SShri Abhyankar count++; 13392b4ed282SShri Abhyankar } 13402b4ed282SShri Abhyankar theend1: 13412b4ed282SShri Abhyankar /* Optional user-defined check for line search step validity */ 13422b4ed282SShri Abhyankar if (vi->postcheckstep && *flag) { 13432b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 13442b4ed282SShri Abhyankar if (changed_y) { 13452b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 13462b4ed282SShri Abhyankar } 13472b4ed282SShri Abhyankar if (changed_y || changed_w) { /* recompute the function if the step has changed */ 13482b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 13492b4ed282SShri Abhyankar if (snes->domainerror) { 13502b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 13512b4ed282SShri Abhyankar PetscFunctionReturn(0); 13522b4ed282SShri Abhyankar } 13532b4ed282SShri Abhyankar ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 13542b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 13552b4ed282SShri Abhyankar ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 13562b4ed282SShri Abhyankar ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 13572b4ed282SShri Abhyankar ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 13582b4ed282SShri Abhyankar } 13592b4ed282SShri Abhyankar } 13602b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 13612b4ed282SShri Abhyankar PetscFunctionReturn(0); 13622b4ed282SShri Abhyankar } 13632b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 13642b4ed282SShri Abhyankar #undef __FUNCT__ 13652b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchQuadratic_VI" 13662b4ed282SShri Abhyankar /* 13672b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchQuadratic in snes/impls/ls/ls.c 13682b4ed282SShri Abhyankar */ 1369ace3abfcSBarry Smith 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,PetscBool *flag) 13702b4ed282SShri Abhyankar { 13712b4ed282SShri Abhyankar /* 13722b4ed282SShri Abhyankar Note that for line search purposes we work with with the related 13732b4ed282SShri Abhyankar minimization problem: 13742b4ed282SShri Abhyankar min z(x): R^n -> R, 13752b4ed282SShri Abhyankar where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 13762b4ed282SShri Abhyankar z(x) is the same as the merit function for the complementarity problem 13772b4ed282SShri Abhyankar */ 13782b4ed282SShri Abhyankar PetscReal initslope,minlambda,lambda,lambdatemp,rellength; 13792b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 13802b4ed282SShri Abhyankar PetscScalar cinitslope; 13812b4ed282SShri Abhyankar #endif 13822b4ed282SShri Abhyankar PetscErrorCode ierr; 13832b4ed282SShri Abhyankar PetscInt count; 13842b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1385ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 13862b4ed282SShri Abhyankar 13872b4ed282SShri Abhyankar PetscFunctionBegin; 13882b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 13892b4ed282SShri Abhyankar *flag = PETSC_TRUE; 13902b4ed282SShri Abhyankar 13912b4ed282SShri Abhyankar ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 13922b4ed282SShri Abhyankar if (*ynorm == 0.0) { 1393c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1394c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr); 13952b4ed282SShri Abhyankar } 13962b4ed282SShri Abhyankar *gnorm = fnorm; 13972b4ed282SShri Abhyankar ierr = VecCopy(x,w);CHKERRQ(ierr); 13982b4ed282SShri Abhyankar ierr = VecCopy(f,g);CHKERRQ(ierr); 13992b4ed282SShri Abhyankar *flag = PETSC_FALSE; 14002b4ed282SShri Abhyankar goto theend2; 14012b4ed282SShri Abhyankar } 14022b4ed282SShri Abhyankar if (*ynorm > vi->maxstep) { /* Step too big, so scale back */ 14032b4ed282SShri Abhyankar ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr); 14042b4ed282SShri Abhyankar *ynorm = vi->maxstep; 14052b4ed282SShri Abhyankar } 14062b4ed282SShri Abhyankar ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 14072b4ed282SShri Abhyankar minlambda = vi->minlambda/rellength; 14082b4ed282SShri Abhyankar /* ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); */ 14092b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 14102b4ed282SShri Abhyankar ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 14112b4ed282SShri Abhyankar initslope = PetscRealPart(cinitslope); 14122b4ed282SShri Abhyankar #else 14132b4ed282SShri Abhyankar ierr = VecDot(vi->dpsi,y,&initslope);CHKERRQ(ierr); 14142b4ed282SShri Abhyankar #endif 14152b4ed282SShri Abhyankar if (initslope > 0.0) initslope = -initslope; 14162b4ed282SShri Abhyankar if (initslope == 0.0) initslope = -1.0; 14172b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr); 14182b4ed282SShri Abhyankar 14192b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 14202b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 14212b4ed282SShri Abhyankar ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 14222b4ed282SShri Abhyankar *flag = PETSC_FALSE; 14232b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 14242b4ed282SShri Abhyankar goto theend2; 14252b4ed282SShri Abhyankar } 14262b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 14272b4ed282SShri Abhyankar if (snes->domainerror) { 14282b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 14292b4ed282SShri Abhyankar PetscFunctionReturn(0); 14302b4ed282SShri Abhyankar } 14312b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 14322b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 14332b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */ 1434c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1435c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 14362b4ed282SShri Abhyankar } 14372b4ed282SShri Abhyankar goto theend2; 14382b4ed282SShri Abhyankar } 14392b4ed282SShri Abhyankar 14402b4ed282SShri Abhyankar /* Fit points with quadratic */ 14412b4ed282SShri Abhyankar lambda = 1.0; 14422b4ed282SShri Abhyankar count = 1; 14432b4ed282SShri Abhyankar while (PETSC_TRUE) { 14442b4ed282SShri Abhyankar if (lambda <= minlambda) { /* bad luck; use full step */ 1445c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1446c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr); 1447c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 14482b4ed282SShri Abhyankar } 14492b4ed282SShri Abhyankar ierr = VecCopy(x,w);CHKERRQ(ierr); 14502b4ed282SShri Abhyankar *flag = PETSC_FALSE; 14512b4ed282SShri Abhyankar break; 14522b4ed282SShri Abhyankar } 14532b4ed282SShri Abhyankar lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 14542b4ed282SShri Abhyankar if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 14552b4ed282SShri Abhyankar if (lambdatemp <= .1*lambda) lambda = .1*lambda; 14562b4ed282SShri Abhyankar else lambda = lambdatemp; 14572b4ed282SShri Abhyankar 14582b4ed282SShri Abhyankar ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 14592b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 14602b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 14612b4ed282SShri 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); 14622b4ed282SShri Abhyankar *flag = PETSC_FALSE; 14632b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 14642b4ed282SShri Abhyankar break; 14652b4ed282SShri Abhyankar } 14662b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 14672b4ed282SShri Abhyankar if (snes->domainerror) { 14682b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 14692b4ed282SShri Abhyankar PetscFunctionReturn(0); 14702b4ed282SShri Abhyankar } 14712b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 14722b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 14732b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */ 1474c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1475c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line Search: Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr); 14762b4ed282SShri Abhyankar } 14772b4ed282SShri Abhyankar break; 14782b4ed282SShri Abhyankar } 14792b4ed282SShri Abhyankar count++; 14802b4ed282SShri Abhyankar } 14812b4ed282SShri Abhyankar theend2: 14822b4ed282SShri Abhyankar /* Optional user-defined check for line search step validity */ 14832b4ed282SShri Abhyankar if (vi->postcheckstep) { 14842b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 14852b4ed282SShri Abhyankar if (changed_y) { 14862b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 14872b4ed282SShri Abhyankar } 14882b4ed282SShri Abhyankar if (changed_y || changed_w) { /* recompute the function if the step has changed */ 14892b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g); 14902b4ed282SShri Abhyankar if (snes->domainerror) { 14912b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 14922b4ed282SShri Abhyankar PetscFunctionReturn(0); 14932b4ed282SShri Abhyankar } 14942b4ed282SShri Abhyankar ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 14952b4ed282SShri Abhyankar ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 14962b4ed282SShri Abhyankar ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 14972b4ed282SShri Abhyankar ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 14982b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 14992b4ed282SShri Abhyankar } 15002b4ed282SShri Abhyankar } 15012b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 15022b4ed282SShri Abhyankar PetscFunctionReturn(0); 15032b4ed282SShri Abhyankar } 15042b4ed282SShri Abhyankar 1505ace3abfcSBarry Smith typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool*); /* force argument to next function to not be extern C*/ 15062b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 15072b4ed282SShri Abhyankar EXTERN_C_BEGIN 15082b4ed282SShri Abhyankar #undef __FUNCT__ 15092b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchSet_VI" 15102b4ed282SShri Abhyankar PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_VI(SNES snes,FCN2 func,void *lsctx) 15112b4ed282SShri Abhyankar { 15122b4ed282SShri Abhyankar PetscFunctionBegin; 15132b4ed282SShri Abhyankar ((SNES_VI *)(snes->data))->LineSearch = func; 15142b4ed282SShri Abhyankar ((SNES_VI *)(snes->data))->lsP = lsctx; 15152b4ed282SShri Abhyankar PetscFunctionReturn(0); 15162b4ed282SShri Abhyankar } 15172b4ed282SShri Abhyankar EXTERN_C_END 15182b4ed282SShri Abhyankar 15192b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 15202b4ed282SShri Abhyankar EXTERN_C_BEGIN 15212b4ed282SShri Abhyankar #undef __FUNCT__ 15222b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchSetMonitor_VI" 1523ace3abfcSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetMonitor_VI(SNES snes,PetscBool flg) 15242b4ed282SShri Abhyankar { 15252b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 15262b4ed282SShri Abhyankar PetscErrorCode ierr; 15272b4ed282SShri Abhyankar 15282b4ed282SShri Abhyankar PetscFunctionBegin; 1529c92abb8aSShri Abhyankar if (flg && !vi->lsmonitor) { 1530c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",((PetscObject)snes)->tablevel,&vi->lsmonitor);CHKERRQ(ierr); 1531c92abb8aSShri Abhyankar } else if (!flg && vi->lsmonitor) { 1532c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr); 15332b4ed282SShri Abhyankar } 15342b4ed282SShri Abhyankar PetscFunctionReturn(0); 15352b4ed282SShri Abhyankar } 15362b4ed282SShri Abhyankar EXTERN_C_END 15372b4ed282SShri Abhyankar 15382b4ed282SShri Abhyankar /* 15392b4ed282SShri Abhyankar SNESView_VI - Prints info from the SNESVI data structure. 15402b4ed282SShri Abhyankar 15412b4ed282SShri Abhyankar Input Parameters: 15422b4ed282SShri Abhyankar . SNES - the SNES context 15432b4ed282SShri Abhyankar . viewer - visualization context 15442b4ed282SShri Abhyankar 15452b4ed282SShri Abhyankar Application Interface Routine: SNESView() 15462b4ed282SShri Abhyankar */ 15472b4ed282SShri Abhyankar #undef __FUNCT__ 15482b4ed282SShri Abhyankar #define __FUNCT__ "SNESView_VI" 15492b4ed282SShri Abhyankar static PetscErrorCode SNESView_VI(SNES snes,PetscViewer viewer) 15502b4ed282SShri Abhyankar { 15512b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI *)snes->data; 15522b4ed282SShri Abhyankar const char *cstr; 15532b4ed282SShri Abhyankar PetscErrorCode ierr; 1554ace3abfcSBarry Smith PetscBool iascii; 15552b4ed282SShri Abhyankar 15562b4ed282SShri Abhyankar PetscFunctionBegin; 15572b4ed282SShri Abhyankar ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 15582b4ed282SShri Abhyankar if (iascii) { 15592b4ed282SShri Abhyankar if (vi->LineSearch == SNESLineSearchNo_VI) cstr = "SNESLineSearchNo"; 15602b4ed282SShri Abhyankar else if (vi->LineSearch == SNESLineSearchQuadratic_VI) cstr = "SNESLineSearchQuadratic"; 15612b4ed282SShri Abhyankar else if (vi->LineSearch == SNESLineSearchCubic_VI) cstr = "SNESLineSearchCubic"; 15622b4ed282SShri Abhyankar else cstr = "unknown"; 15632b4ed282SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 15642b4ed282SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," alpha=%G, maxstep=%G, minlambda=%G\n",vi->alpha,vi->maxstep,vi->minlambda);CHKERRQ(ierr); 15652b4ed282SShri Abhyankar } else { 15662b4ed282SShri Abhyankar SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ VI",((PetscObject)viewer)->type_name); 15672b4ed282SShri Abhyankar } 15682b4ed282SShri Abhyankar PetscFunctionReturn(0); 15692b4ed282SShri Abhyankar } 15702b4ed282SShri Abhyankar 15712b4ed282SShri Abhyankar /* 15722b4ed282SShri Abhyankar SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. 15732b4ed282SShri Abhyankar 15742b4ed282SShri Abhyankar Input Parameters: 15752b4ed282SShri Abhyankar . snes - the SNES context. 15762b4ed282SShri Abhyankar . xl - lower bound. 15772b4ed282SShri Abhyankar . xu - upper bound. 15782b4ed282SShri Abhyankar 15792b4ed282SShri Abhyankar Notes: 15802b4ed282SShri Abhyankar If this routine is not called then the lower and upper bounds are set to 15812b4ed282SShri Abhyankar -Infinity and Infinity respectively during SNESSetUp. 15822b4ed282SShri Abhyankar */ 15832b4ed282SShri Abhyankar 15842b4ed282SShri Abhyankar #undef __FUNCT__ 15852b4ed282SShri Abhyankar #define __FUNCT__ "SNESVISetVariableBounds" 158669c03620SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu) 15872b4ed282SShri Abhyankar { 15882b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 15892b4ed282SShri Abhyankar 15902b4ed282SShri Abhyankar PetscFunctionBegin; 15912b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 15922b4ed282SShri Abhyankar PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 15932b4ed282SShri Abhyankar PetscValidHeaderSpecific(xu,VEC_CLASSID,3); 15942b4ed282SShri Abhyankar 15952b4ed282SShri Abhyankar /* Check if SNESSetFunction is called */ 15962b4ed282SShri Abhyankar if(!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 15972b4ed282SShri Abhyankar 15982b4ed282SShri Abhyankar /* Check if the vector sizes are compatible for lower and upper bounds */ 15992b4ed282SShri 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); 16002b4ed282SShri 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); 16012b4ed282SShri Abhyankar vi->usersetxbounds = PETSC_TRUE; 16022b4ed282SShri Abhyankar vi->xl = xl; 16032b4ed282SShri Abhyankar vi->xu = xu; 16042b4ed282SShri Abhyankar 16052b4ed282SShri Abhyankar PetscFunctionReturn(0); 16062b4ed282SShri Abhyankar } 16072b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 16082b4ed282SShri Abhyankar /* 16092b4ed282SShri Abhyankar SNESSetFromOptions_VI - Sets various parameters for the SNESVI method. 16102b4ed282SShri Abhyankar 16112b4ed282SShri Abhyankar Input Parameter: 16122b4ed282SShri Abhyankar . snes - the SNES context 16132b4ed282SShri Abhyankar 16142b4ed282SShri Abhyankar Application Interface Routine: SNESSetFromOptions() 16152b4ed282SShri Abhyankar */ 16162b4ed282SShri Abhyankar #undef __FUNCT__ 16172b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetFromOptions_VI" 16182b4ed282SShri Abhyankar static PetscErrorCode SNESSetFromOptions_VI(SNES snes) 16192b4ed282SShri Abhyankar { 16202b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI *)snes->data; 16212b4ed282SShri Abhyankar const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 162269c03620SShri Abhyankar const char *vies[] = {"ss","as"}; 16232b4ed282SShri Abhyankar PetscErrorCode ierr; 16242b4ed282SShri Abhyankar PetscInt indx; 162569c03620SShri Abhyankar PetscBool flg,set,flg2; 16262b4ed282SShri Abhyankar 16272b4ed282SShri Abhyankar PetscFunctionBegin; 16282b4ed282SShri Abhyankar ierr = PetscOptionsHead("SNES semismooth method options");CHKERRQ(ierr); 16292b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_alpha","Function norm must decrease by","None",vi->alpha,&vi->alpha,0);CHKERRQ(ierr); 16302b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_maxstep","Step must be less than","None",vi->maxstep,&vi->maxstep,0);CHKERRQ(ierr); 16312b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_minlambda","Minimum lambda allowed","None",vi->minlambda,&vi->minlambda,0);CHKERRQ(ierr); 16322b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_delta","descent test fraction","None",vi->delta,&vi->delta,0);CHKERRQ(ierr); 16332b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_rho","descent test power","None",vi->rho,&vi->rho,0);CHKERRQ(ierr); 16342b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_const_tol","constraint tolerance","None",vi->const_tol,&vi->const_tol,0);CHKERRQ(ierr); 1635acfcf0e5SJed Brown ierr = PetscOptionsBool("-snes_vi_lsmonitor","Print progress of line searches","SNESLineSearchSetMonitor",vi->lsmonitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 16362b4ed282SShri Abhyankar if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);} 163769c03620SShri Abhyankar ierr = PetscOptionsEList("-snes_vi_type","Semismooth algorithm used","",vies,2,"ss",&indx,&flg2);CHKERRQ(ierr); 163869c03620SShri Abhyankar if (flg2) { 163969c03620SShri Abhyankar switch (indx) { 164069c03620SShri Abhyankar case 0: 164169c03620SShri Abhyankar snes->ops->solve = SNESSolveVI_SS; 164269c03620SShri Abhyankar break; 164369c03620SShri Abhyankar case 1: 164469c03620SShri Abhyankar snes->ops->solve = SNESSolveVI_AS; 164569c03620SShri Abhyankar break; 164669c03620SShri Abhyankar } 164769c03620SShri Abhyankar } 1648c92abb8aSShri Abhyankar ierr = PetscOptionsEList("-snes_vi_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 16492b4ed282SShri Abhyankar if (flg) { 16502b4ed282SShri Abhyankar switch (indx) { 16512b4ed282SShri Abhyankar case 0: 16522b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr); 16532b4ed282SShri Abhyankar break; 16542b4ed282SShri Abhyankar case 1: 16552b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr); 16562b4ed282SShri Abhyankar break; 16572b4ed282SShri Abhyankar case 2: 16582b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr); 16592b4ed282SShri Abhyankar break; 16602b4ed282SShri Abhyankar case 3: 16612b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr); 16622b4ed282SShri Abhyankar break; 16632b4ed282SShri Abhyankar } 16642b4ed282SShri Abhyankar } 16652b4ed282SShri Abhyankar ierr = PetscOptionsTail();CHKERRQ(ierr); 16662b4ed282SShri Abhyankar PetscFunctionReturn(0); 16672b4ed282SShri Abhyankar } 16682b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 16692b4ed282SShri Abhyankar /*MC 16702b4ed282SShri Abhyankar SNESVI - Semismooth newton method based nonlinear solver that uses a line search 16712b4ed282SShri Abhyankar 16722b4ed282SShri Abhyankar Options Database: 16732b4ed282SShri Abhyankar + -snes_vi [cubic,quadratic,basic,basicnonorms] - Selects line search 16742b4ed282SShri Abhyankar . -snes_vi_alpha <alpha> - Sets alpha 16752b4ed282SShri 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) 16762b4ed282SShri Abhyankar . -snes_vi_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 16772b4ed282SShri Abhyankar -snes_vi_delta <delta> - Sets the fraction used in the descent test. 16782b4ed282SShri Abhyankar -snes_vi_rho <rho> - Sets the power used in the descent test. 16792b4ed282SShri Abhyankar For a descent direction to be accepted it has to satisfy the condition dpsi^T*d <= -delta*||d||^rho 16802b4ed282SShri Abhyankar - -snes_vi_monitor - print information about progress of line searches 16812b4ed282SShri Abhyankar 16822b4ed282SShri Abhyankar 16832b4ed282SShri Abhyankar Level: beginner 16842b4ed282SShri Abhyankar 16852b4ed282SShri Abhyankar .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 16862b4ed282SShri Abhyankar SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 16872b4ed282SShri Abhyankar SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams() 16882b4ed282SShri Abhyankar 16892b4ed282SShri Abhyankar M*/ 16902b4ed282SShri Abhyankar EXTERN_C_BEGIN 16912b4ed282SShri Abhyankar #undef __FUNCT__ 16922b4ed282SShri Abhyankar #define __FUNCT__ "SNESCreate_VI" 16932b4ed282SShri Abhyankar PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_VI(SNES snes) 16942b4ed282SShri Abhyankar { 16952b4ed282SShri Abhyankar PetscErrorCode ierr; 16962b4ed282SShri Abhyankar SNES_VI *vi; 16972b4ed282SShri Abhyankar 16982b4ed282SShri Abhyankar PetscFunctionBegin; 16992b4ed282SShri Abhyankar snes->ops->setup = SNESSetUp_VI; 170069c03620SShri Abhyankar snes->ops->solve = SNESSolveVI_SS; 17012b4ed282SShri Abhyankar snes->ops->destroy = SNESDestroy_VI; 17022b4ed282SShri Abhyankar snes->ops->setfromoptions = SNESSetFromOptions_VI; 17032b4ed282SShri Abhyankar snes->ops->view = SNESView_VI; 17042b4ed282SShri Abhyankar snes->ops->converged = SNESDefaultConverged_VI; 17052b4ed282SShri Abhyankar 17062b4ed282SShri Abhyankar ierr = PetscNewLog(snes,SNES_VI,&vi);CHKERRQ(ierr); 17072b4ed282SShri Abhyankar snes->data = (void*)vi; 17082b4ed282SShri Abhyankar vi->alpha = 1.e-4; 17092b4ed282SShri Abhyankar vi->maxstep = 1.e8; 17102b4ed282SShri Abhyankar vi->minlambda = 1.e-12; 17112b4ed282SShri Abhyankar vi->LineSearch = SNESLineSearchCubic_VI; 17122b4ed282SShri Abhyankar vi->lsP = PETSC_NULL; 17132b4ed282SShri Abhyankar vi->postcheckstep = PETSC_NULL; 17142b4ed282SShri Abhyankar vi->postcheck = PETSC_NULL; 17152b4ed282SShri Abhyankar vi->precheckstep = PETSC_NULL; 17162b4ed282SShri Abhyankar vi->precheck = PETSC_NULL; 17172b4ed282SShri Abhyankar vi->rho = 2.1; 17182b4ed282SShri Abhyankar vi->delta = 1e-10; 17192b4ed282SShri Abhyankar vi->const_tol = 2.2204460492503131e-16; 17202b4ed282SShri Abhyankar vi->computessfunction = ComputeFischerFunction; 17212b4ed282SShri Abhyankar 17222b4ed282SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_VI",SNESLineSearchSetMonitor_VI);CHKERRQ(ierr); 17232b4ed282SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_VI",SNESLineSearchSet_VI);CHKERRQ(ierr); 17242b4ed282SShri Abhyankar 17252b4ed282SShri Abhyankar PetscFunctionReturn(0); 17262b4ed282SShri Abhyankar } 17272b4ed282SShri Abhyankar EXTERN_C_END 1728