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 /* 240*a79edbefSShri Abhyankar SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the 241*a79edbefSShri Abhyankar the semismooth jacobian. 2422b4ed282SShri Abhyankar */ 2432b4ed282SShri Abhyankar #undef __FUNCT__ 244*a79edbefSShri Abhyankar #define __FUNCT__ "SNESVIComputeBsubdifferentialVectors" 245*a79edbefSShri 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); 258*a79edbefSShri Abhyankar ierr = VecGetArray(Da,&da);CHKERRQ(ierr); 259*a79edbefSShri 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); 281*a79edbefSShri 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); 356*a79edbefSShri Abhyankar ierr = VecRestoreArray(Da,&da);CHKERRQ(ierr); 357*a79edbefSShri Abhyankar ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr); 3582b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->t,&t);CHKERRQ(ierr); 3592b4ed282SShri Abhyankar 360*a79edbefSShri Abhyankar PetscFunctionReturn(0); 361*a79edbefSShri Abhyankar } 362*a79edbefSShri Abhyankar 363*a79edbefSShri Abhyankar /* 364*a79edbefSShri 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. 365*a79edbefSShri Abhyankar 366*a79edbefSShri Abhyankar Input Parameters: 367*a79edbefSShri Abhyankar . Da - Diagonal shift vector for the semismooth jacobian. 368*a79edbefSShri Abhyankar . Db - Row scaling vector for the semismooth jacobian. 369*a79edbefSShri Abhyankar 370*a79edbefSShri Abhyankar Output Parameters: 371*a79edbefSShri Abhyankar . jac - semismooth jacobian 372*a79edbefSShri Abhyankar . jac_pre - optional preconditioning matrix 373*a79edbefSShri Abhyankar 374*a79edbefSShri Abhyankar Notes: 375*a79edbefSShri Abhyankar The semismooth jacobian matrix is given by 376*a79edbefSShri Abhyankar jac = Da + Db*jacfun 377*a79edbefSShri Abhyankar where Db is the row scaling matrix stored as a vector, 378*a79edbefSShri Abhyankar Da is the diagonal perturbation matrix stored as a vector 379*a79edbefSShri Abhyankar and jacfun is the jacobian of the original nonlinear function. 380*a79edbefSShri Abhyankar */ 381*a79edbefSShri Abhyankar #undef __FUNCT__ 382*a79edbefSShri Abhyankar #define __FUNCT__ "SNESVIComputeJacobian" 383*a79edbefSShri Abhyankar PetscErrorCode SNESVIComputeJacobian(Mat jac, Mat jac_pre,Vec Da, Vec Db) 384*a79edbefSShri Abhyankar { 385*a79edbefSShri Abhyankar PetscErrorCode ierr; 386*a79edbefSShri Abhyankar 3872b4ed282SShri Abhyankar /* Do row scaling and add diagonal perturbation */ 388*a79edbefSShri Abhyankar ierr = MatDiagonalScale(jac,Db,PETSC_NULL);CHKERRQ(ierr); 389*a79edbefSShri Abhyankar ierr = MatDiagonalSet(jac,Da,ADD_VALUES);CHKERRQ(ierr); 390*a79edbefSShri Abhyankar if (jac != jac_pre) { /* If jac and jac_pre are different */ 391*a79edbefSShri Abhyankar ierr = MatDiagonalScale(jac_pre,Db,PETSC_NULL); 392*a79edbefSShri 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 4602b4ed282SShri 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 */ 623*a79edbefSShri Abhyankar /* Get the nonlinear function jacobian */ 6242b4ed282SShri Abhyankar ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 625*a79edbefSShri Abhyankar /* Get the diagonal shift and row scaling vectors */ 626*a79edbefSShri Abhyankar ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr); 627*a79edbefSShri Abhyankar /* Compute the semismooth jacobian */ 628*a79edbefSShri 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" 704*a79edbefSShri 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++) { 71869c03620SShri Abhyankar if (db[i] <= thresh) nloc_isact++; 71969c03620SShri Abhyankar else nloc_isact++; 72069c03620SShri Abhyankar } 72169c03620SShri Abhyankar 72269c03620SShri Abhyankar ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr); 72369c03620SShri Abhyankar ierr = PetscMalloc(nloc_isinact*sizeof(PetscInt),&idx_inact);CHKERRQ(ierr); 72469c03620SShri Abhyankar 72569c03620SShri Abhyankar /* Creating the indexing arrays */ 72669c03620SShri Abhyankar for(i=ilow; i < ihigh; i++) { 72769c03620SShri Abhyankar if (db[i] <= thresh) idx_act[i1++] = i; 72869c03620SShri Abhyankar else idx_inact[i2++] = i; 72969c03620SShri Abhyankar } 73069c03620SShri Abhyankar 73169c03620SShri Abhyankar /* Create the index sets */ 73269c03620SShri Abhyankar ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isact,idx_act,PETSC_COPY_VALUES,ISact);CHKERRQ(ierr); 73369c03620SShri Abhyankar ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isinact,idx_inact,PETSC_COPY_VALUES,ISinact);CHKERRQ(ierr); 73469c03620SShri Abhyankar 73569c03620SShri Abhyankar ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr); 73669c03620SShri Abhyankar ierr = PetscFree(idx_act);CHKERRQ(ierr); 73769c03620SShri Abhyankar ierr = PetscFree(idx_inact);CHKERRQ(ierr); 73869c03620SShri Abhyankar PetscFunctionReturn(0); 73969c03620SShri Abhyankar } 74069c03620SShri Abhyankar 74169c03620SShri Abhyankar /* Variational Inequality solver using active set method */ 74269c03620SShri Abhyankar #undef __FUNCT__ 74369c03620SShri Abhyankar #define __FUNCT__ "SNESSolveVI_AS" 74469c03620SShri Abhyankar PetscErrorCode SNESSolveVI_AS(SNES snes) 74569c03620SShri Abhyankar { 74669c03620SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 74769c03620SShri Abhyankar PetscErrorCode ierr; 74869c03620SShri Abhyankar PetscInt maxits,i,lits; 74969c03620SShri Abhyankar PetscBool lssucceed,changedir; 75069c03620SShri Abhyankar MatStructure flg = DIFFERENT_NONZERO_PATTERN; 75169c03620SShri Abhyankar PetscReal gnorm,xnorm=0,ynorm; 75269c03620SShri Abhyankar Vec Y,X,F,G,W; 75369c03620SShri Abhyankar KSPConvergedReason kspreason; 75469c03620SShri Abhyankar 75569c03620SShri Abhyankar PetscFunctionBegin; 75669c03620SShri Abhyankar snes->numFailures = 0; 75769c03620SShri Abhyankar snes->numLinearSolveFailures = 0; 75869c03620SShri Abhyankar snes->reason = SNES_CONVERGED_ITERATING; 75969c03620SShri Abhyankar 76069c03620SShri Abhyankar maxits = snes->max_its; /* maximum number of iterations */ 76169c03620SShri Abhyankar X = snes->vec_sol; /* solution vector */ 76269c03620SShri Abhyankar F = snes->vec_func; /* residual vector */ 76369c03620SShri Abhyankar Y = snes->work[0]; /* work vectors */ 76469c03620SShri Abhyankar G = snes->work[1]; 76569c03620SShri Abhyankar W = snes->work[2]; 76669c03620SShri Abhyankar 76769c03620SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 76869c03620SShri Abhyankar snes->iter = 0; 76969c03620SShri Abhyankar snes->norm = 0.0; 77069c03620SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 77169c03620SShri Abhyankar 77269c03620SShri Abhyankar ierr = SNESVIAdjustInitialGuess(X,vi->xl,vi->xu);CHKERRQ(ierr); 77369c03620SShri Abhyankar ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr); 77469c03620SShri Abhyankar if (snes->domainerror) { 77569c03620SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 77669c03620SShri Abhyankar PetscFunctionReturn(0); 77769c03620SShri Abhyankar } 77869c03620SShri Abhyankar /* Compute Merit function */ 77969c03620SShri Abhyankar ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr); 78069c03620SShri Abhyankar 78169c03620SShri Abhyankar ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 78269c03620SShri Abhyankar ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 78369c03620SShri Abhyankar if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 78469c03620SShri Abhyankar 78569c03620SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 78669c03620SShri Abhyankar snes->norm = vi->phinorm; 78769c03620SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 78869c03620SShri Abhyankar SNESLogConvHistory(snes,vi->phinorm,0); 78969c03620SShri Abhyankar SNESMonitor(snes,0,vi->phinorm); 79069c03620SShri Abhyankar 79169c03620SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 79269c03620SShri Abhyankar snes->ttol = vi->phinorm*snes->rtol; 79369c03620SShri Abhyankar /* test convergence */ 79469c03620SShri Abhyankar ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 79569c03620SShri Abhyankar if (snes->reason) PetscFunctionReturn(0); 79669c03620SShri Abhyankar 79769c03620SShri Abhyankar for (i=0; i<maxits; i++) { 79869c03620SShri Abhyankar 799*a79edbefSShri Abhyankar IS IS_act,IS_inact; /* _act -> active set _inact -> inactive set */ 800*a79edbefSShri Abhyankar PetscReal thresh,J_norm1; 801*a79edbefSShri Abhyankar 80269c03620SShri Abhyankar /* Call general purpose update function */ 80369c03620SShri Abhyankar if (snes->ops->update) { 80469c03620SShri Abhyankar ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 80569c03620SShri Abhyankar } 80669c03620SShri Abhyankar ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 80769c03620SShri Abhyankar /* Compute the threshold value for creating active and inactive sets */ 80869c03620SShri Abhyankar ierr = MatNorm(snes->jacobian,NORM_1,&J_norm1);CHKERRQ(ierr); 80969c03620SShri Abhyankar thresh = PetscMin(vi->merit,1e-2)/(1+J_norm1); 810*a79edbefSShri Abhyankar 811*a79edbefSShri Abhyankar /* Compute B-subdifferential vectors Da and Db */ 812*a79edbefSShri Abhyankar ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr); 81369c03620SShri Abhyankar /* Create active and inactive index sets */ 81469c03620SShri Abhyankar ierr = SNESVICreateIndexSets_AS(snes,vi->Db,thresh,&IS_act,&IS_inact);CHKERRQ(ierr); 81569c03620SShri Abhyankar ierr = VecView(vi->Db,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 81669c03620SShri Abhyankar ierr = ISView(IS_act,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 81769c03620SShri Abhyankar ierr = ISView(IS_inact,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 818*a79edbefSShri Abhyankar SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Active set semismooth algorithm not implemented yet"); 81969c03620SShri Abhyankar 820730c24dcSShri Abhyankar ierr = ISDestroy(IS_act);CHKERRQ(ierr); 821730c24dcSShri Abhyankar ierr = ISDestroy(IS_inact);CHKERRQ(ierr); 822730c24dcSShri Abhyankar 82369c03620SShri Abhyankar ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 82469c03620SShri Abhyankar ierr = SNES_KSPSolve(snes,snes->ksp,vi->phi,Y);CHKERRQ(ierr); 82569c03620SShri Abhyankar ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 82669c03620SShri Abhyankar ierr = SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);CHKERRQ(ierr); 82769c03620SShri Abhyankar ierr = SNESVICheckDescentDirection(snes,vi->dpsi,Y,&changedir);CHKERRQ(ierr); 82869c03620SShri Abhyankar if (kspreason < 0 || changedir) { 82969c03620SShri Abhyankar if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 83069c03620SShri 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); 83169c03620SShri Abhyankar snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 83269c03620SShri Abhyankar break; 83369c03620SShri Abhyankar } 83469c03620SShri Abhyankar ierr = VecCopy(vi->dpsi,Y);CHKERRQ(ierr); 83569c03620SShri Abhyankar } 83669c03620SShri Abhyankar ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 83769c03620SShri Abhyankar snes->linear_its += lits; 83869c03620SShri Abhyankar ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 83969c03620SShri Abhyankar /* 84069c03620SShri Abhyankar if (vi->precheckstep) { 84169c03620SShri Abhyankar PetscBool changed_y = PETSC_FALSE; 84269c03620SShri Abhyankar ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr); 84369c03620SShri Abhyankar } 84469c03620SShri Abhyankar 84569c03620SShri Abhyankar if (PetscLogPrintInfo){ 84669c03620SShri Abhyankar ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 84769c03620SShri Abhyankar } 84869c03620SShri Abhyankar */ 84969c03620SShri Abhyankar /* Compute a (scaled) negative update in the line search routine: 85069c03620SShri Abhyankar Y <- X - lambda*Y 85169c03620SShri Abhyankar and evaluate G = function(Y) (depends on the line search). 85269c03620SShri Abhyankar */ 85369c03620SShri Abhyankar ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 85469c03620SShri Abhyankar ynorm = 1; gnorm = vi->phinorm; 85569c03620SShri Abhyankar ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 85669c03620SShri Abhyankar ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 85769c03620SShri Abhyankar if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 85869c03620SShri Abhyankar if (snes->domainerror) { 85969c03620SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 86069c03620SShri Abhyankar PetscFunctionReturn(0); 86169c03620SShri Abhyankar } 86269c03620SShri Abhyankar if (!lssucceed) { 86369c03620SShri Abhyankar if (++snes->numFailures >= snes->maxFailures) { 86469c03620SShri Abhyankar PetscBool ismin; 86569c03620SShri Abhyankar snes->reason = SNES_DIVERGED_LINE_SEARCH; 86669c03620SShri Abhyankar ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 86769c03620SShri Abhyankar if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 86869c03620SShri Abhyankar break; 86969c03620SShri Abhyankar } 87069c03620SShri Abhyankar } 87169c03620SShri Abhyankar /* Update function and solution vectors */ 87269c03620SShri Abhyankar vi->phinorm = gnorm; 87369c03620SShri Abhyankar vi->merit = 0.5*vi->phinorm*vi->phinorm; 87469c03620SShri Abhyankar ierr = VecCopy(G,vi->phi);CHKERRQ(ierr); 87569c03620SShri Abhyankar ierr = VecCopy(W,X);CHKERRQ(ierr); 87669c03620SShri Abhyankar /* Monitor convergence */ 87769c03620SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 87869c03620SShri Abhyankar snes->iter = i+1; 87969c03620SShri Abhyankar snes->norm = vi->phinorm; 88069c03620SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 88169c03620SShri Abhyankar SNESLogConvHistory(snes,snes->norm,lits); 88269c03620SShri Abhyankar SNESMonitor(snes,snes->iter,snes->norm); 88369c03620SShri Abhyankar /* Test for convergence, xnorm = || X || */ 88469c03620SShri Abhyankar if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 88569c03620SShri Abhyankar ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 88669c03620SShri Abhyankar if (snes->reason) break; 88769c03620SShri Abhyankar } 88869c03620SShri Abhyankar if (i == maxits) { 88969c03620SShri Abhyankar ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 89069c03620SShri Abhyankar if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 89169c03620SShri Abhyankar } 89269c03620SShri Abhyankar PetscFunctionReturn(0); 89369c03620SShri Abhyankar } 89469c03620SShri Abhyankar 8952b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 8962b4ed282SShri Abhyankar /* 8972b4ed282SShri Abhyankar SNESSetUp_VI - Sets up the internal data structures for the later use 8982b4ed282SShri Abhyankar of the SNESVI nonlinear solver. 8992b4ed282SShri Abhyankar 9002b4ed282SShri Abhyankar Input Parameter: 9012b4ed282SShri Abhyankar . snes - the SNES context 9022b4ed282SShri Abhyankar . x - the solution vector 9032b4ed282SShri Abhyankar 9042b4ed282SShri Abhyankar Application Interface Routine: SNESSetUp() 9052b4ed282SShri Abhyankar 9062b4ed282SShri Abhyankar Notes: 9072b4ed282SShri Abhyankar For basic use of the SNES solvers, the user need not explicitly call 9082b4ed282SShri Abhyankar SNESSetUp(), since these actions will automatically occur during 9092b4ed282SShri Abhyankar the call to SNESSolve(). 9102b4ed282SShri Abhyankar */ 9112b4ed282SShri Abhyankar #undef __FUNCT__ 9122b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetUp_VI" 9132b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes) 9142b4ed282SShri Abhyankar { 9152b4ed282SShri Abhyankar PetscErrorCode ierr; 9162b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*) snes->data; 9172b4ed282SShri Abhyankar PetscInt i_start[3],i_end[3]; 9182b4ed282SShri Abhyankar 9192b4ed282SShri Abhyankar PetscFunctionBegin; 9202b4ed282SShri Abhyankar if (!snes->vec_sol_update) { 9212b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr); 9222b4ed282SShri Abhyankar ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr); 9232b4ed282SShri Abhyankar } 9242b4ed282SShri Abhyankar if (!snes->work) { 9252b4ed282SShri Abhyankar snes->nwork = 3; 9262b4ed282SShri Abhyankar ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 9272b4ed282SShri Abhyankar ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr); 9282b4ed282SShri Abhyankar } 9292b4ed282SShri Abhyankar 9302b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->phi); CHKERRQ(ierr); 9312b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->dpsi); CHKERRQ(ierr); 9322b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->Da); CHKERRQ(ierr); 9332b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->Db); CHKERRQ(ierr); 9342b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->z);CHKERRQ(ierr); 9352b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->t); CHKERRQ(ierr); 9362b4ed282SShri Abhyankar 9372b4ed282SShri Abhyankar /* If the lower and upper bound on variables are not set, set it to 9382b4ed282SShri Abhyankar -Inf and Inf */ 9392b4ed282SShri Abhyankar if (!vi->xl && !vi->xu) { 9402b4ed282SShri Abhyankar vi->usersetxbounds = PETSC_FALSE; 9412b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->xl); CHKERRQ(ierr); 9422b4ed282SShri Abhyankar ierr = VecSet(vi->xl,PETSC_VI_NINF);CHKERRQ(ierr); 9432b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->xu); CHKERRQ(ierr); 9442b4ed282SShri Abhyankar ierr = VecSet(vi->xu,PETSC_VI_INF);CHKERRQ(ierr); 9452b4ed282SShri Abhyankar } else { 9462b4ed282SShri Abhyankar /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */ 9472b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr); 9482b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(vi->xl,i_start+1,i_end+1);CHKERRQ(ierr); 9492b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(vi->xu,i_start+2,i_end+2);CHKERRQ(ierr); 9502b4ed282SShri 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])) 9512b4ed282SShri 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."); 9522b4ed282SShri Abhyankar } 9532b4ed282SShri Abhyankar 9542b4ed282SShri Abhyankar vi->computeuserfunction = snes->ops->computefunction; 9551f79c099SShri Abhyankar snes->ops->computefunction = SNESVIComputeFunction; 956*a79edbefSShri Abhyankar 9572b4ed282SShri Abhyankar PetscFunctionReturn(0); 9582b4ed282SShri Abhyankar } 9592b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 9602b4ed282SShri Abhyankar /* 9612b4ed282SShri Abhyankar SNESDestroy_VI - Destroys the private SNES_VI context that was created 9622b4ed282SShri Abhyankar with SNESCreate_VI(). 9632b4ed282SShri Abhyankar 9642b4ed282SShri Abhyankar Input Parameter: 9652b4ed282SShri Abhyankar . snes - the SNES context 9662b4ed282SShri Abhyankar 9672b4ed282SShri Abhyankar Application Interface Routine: SNESDestroy() 9682b4ed282SShri Abhyankar */ 9692b4ed282SShri Abhyankar #undef __FUNCT__ 9702b4ed282SShri Abhyankar #define __FUNCT__ "SNESDestroy_VI" 9712b4ed282SShri Abhyankar PetscErrorCode SNESDestroy_VI(SNES snes) 9722b4ed282SShri Abhyankar { 9732b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*) snes->data; 9742b4ed282SShri Abhyankar PetscErrorCode ierr; 9752b4ed282SShri Abhyankar 9762b4ed282SShri Abhyankar PetscFunctionBegin; 9772b4ed282SShri Abhyankar if (snes->vec_sol_update) { 9782b4ed282SShri Abhyankar ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr); 9792b4ed282SShri Abhyankar snes->vec_sol_update = PETSC_NULL; 9802b4ed282SShri Abhyankar } 9812b4ed282SShri Abhyankar if (snes->nwork) { 9822b4ed282SShri Abhyankar ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 9832b4ed282SShri Abhyankar snes->nwork = 0; 9842b4ed282SShri Abhyankar snes->work = PETSC_NULL; 9852b4ed282SShri Abhyankar } 9862b4ed282SShri Abhyankar 9872b4ed282SShri Abhyankar /* clear vectors */ 9882b4ed282SShri Abhyankar ierr = VecDestroy(vi->phi); CHKERRQ(ierr); 9892b4ed282SShri Abhyankar ierr = VecDestroy(vi->dpsi); CHKERRQ(ierr); 9902b4ed282SShri Abhyankar ierr = VecDestroy(vi->Da); CHKERRQ(ierr); 9912b4ed282SShri Abhyankar ierr = VecDestroy(vi->Db); CHKERRQ(ierr); 9922b4ed282SShri Abhyankar ierr = VecDestroy(vi->z); CHKERRQ(ierr); 9932b4ed282SShri Abhyankar ierr = VecDestroy(vi->t); CHKERRQ(ierr); 9942b4ed282SShri Abhyankar if (!vi->usersetxbounds) { 9952b4ed282SShri Abhyankar ierr = VecDestroy(vi->xl); CHKERRQ(ierr); 9962b4ed282SShri Abhyankar ierr = VecDestroy(vi->xu); CHKERRQ(ierr); 9972b4ed282SShri Abhyankar } 998c92abb8aSShri Abhyankar if (vi->lsmonitor) { 999c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr); 10002b4ed282SShri Abhyankar } 10012b4ed282SShri Abhyankar ierr = PetscFree(snes->data);CHKERRQ(ierr); 10022b4ed282SShri Abhyankar 10032b4ed282SShri Abhyankar /* clear composed functions */ 10042b4ed282SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr); 1005c92abb8aSShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr); 10062b4ed282SShri Abhyankar PetscFunctionReturn(0); 10072b4ed282SShri Abhyankar } 10082b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 10092b4ed282SShri Abhyankar #undef __FUNCT__ 10102b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNo_VI" 10112b4ed282SShri Abhyankar 10122b4ed282SShri Abhyankar /* 10132b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchNo routine in snes/impls/ls/ls.c 10142b4ed282SShri Abhyankar 10152b4ed282SShri Abhyankar */ 1016ace3abfcSBarry 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) 10172b4ed282SShri Abhyankar { 10182b4ed282SShri Abhyankar PetscErrorCode ierr; 10192b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1020ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 10212b4ed282SShri Abhyankar 10222b4ed282SShri Abhyankar PetscFunctionBegin; 10232b4ed282SShri Abhyankar *flag = PETSC_TRUE; 10242b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 10252b4ed282SShri Abhyankar ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 10262b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 10272b4ed282SShri Abhyankar if (vi->postcheckstep) { 10282b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 10292b4ed282SShri Abhyankar } 10302b4ed282SShri Abhyankar if (changed_y) { 10312b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 10322b4ed282SShri Abhyankar } 10332b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 10342b4ed282SShri Abhyankar if (!snes->domainerror) { 10352b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 10362b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 10372b4ed282SShri Abhyankar } 1038c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1039c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 1040c92abb8aSShri Abhyankar } 10412b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 10422b4ed282SShri Abhyankar PetscFunctionReturn(0); 10432b4ed282SShri Abhyankar } 10442b4ed282SShri Abhyankar 10452b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 10462b4ed282SShri Abhyankar #undef __FUNCT__ 10472b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNoNorms_VI" 10482b4ed282SShri Abhyankar 10492b4ed282SShri Abhyankar /* 10502b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c 10512b4ed282SShri Abhyankar */ 1052ace3abfcSBarry 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) 10532b4ed282SShri Abhyankar { 10542b4ed282SShri Abhyankar PetscErrorCode ierr; 10552b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1056ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 10572b4ed282SShri Abhyankar 10582b4ed282SShri Abhyankar PetscFunctionBegin; 10592b4ed282SShri Abhyankar *flag = PETSC_TRUE; 10602b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 10612b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 10622b4ed282SShri Abhyankar if (vi->postcheckstep) { 10632b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 10642b4ed282SShri Abhyankar } 10652b4ed282SShri Abhyankar if (changed_y) { 10662b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 10672b4ed282SShri Abhyankar } 10682b4ed282SShri Abhyankar 10692b4ed282SShri Abhyankar /* don't evaluate function the last time through */ 10702b4ed282SShri Abhyankar if (snes->iter < snes->max_its-1) { 10712b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 10722b4ed282SShri Abhyankar } 10732b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 10742b4ed282SShri Abhyankar PetscFunctionReturn(0); 10752b4ed282SShri Abhyankar } 10762b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 10772b4ed282SShri Abhyankar #undef __FUNCT__ 10782b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchCubic_VI" 10792b4ed282SShri Abhyankar /* 10802b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchCubic in snes/impls/ls/ls.c 10812b4ed282SShri Abhyankar */ 1082ace3abfcSBarry 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) 10832b4ed282SShri Abhyankar { 10842b4ed282SShri Abhyankar /* 10852b4ed282SShri Abhyankar Note that for line search purposes we work with with the related 10862b4ed282SShri Abhyankar minimization problem: 10872b4ed282SShri Abhyankar min z(x): R^n -> R, 10882b4ed282SShri Abhyankar where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 10892b4ed282SShri Abhyankar This function z(x) is same as the merit function for the complementarity problem. 10902b4ed282SShri Abhyankar */ 10912b4ed282SShri Abhyankar 10922b4ed282SShri Abhyankar PetscReal initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 10932b4ed282SShri Abhyankar PetscReal minlambda,lambda,lambdatemp; 10942b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 10952b4ed282SShri Abhyankar PetscScalar cinitslope; 10962b4ed282SShri Abhyankar #endif 10972b4ed282SShri Abhyankar PetscErrorCode ierr; 10982b4ed282SShri Abhyankar PetscInt count; 10992b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1100ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 11012b4ed282SShri Abhyankar MPI_Comm comm; 11022b4ed282SShri Abhyankar 11032b4ed282SShri Abhyankar PetscFunctionBegin; 11042b4ed282SShri Abhyankar ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 11052b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 11062b4ed282SShri Abhyankar *flag = PETSC_TRUE; 11072b4ed282SShri Abhyankar 11082b4ed282SShri Abhyankar ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 11092b4ed282SShri Abhyankar if (*ynorm == 0.0) { 1110c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1111c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 11122b4ed282SShri Abhyankar } 11132b4ed282SShri Abhyankar *gnorm = fnorm; 11142b4ed282SShri Abhyankar ierr = VecCopy(x,w);CHKERRQ(ierr); 11152b4ed282SShri Abhyankar ierr = VecCopy(f,g);CHKERRQ(ierr); 11162b4ed282SShri Abhyankar *flag = PETSC_FALSE; 11172b4ed282SShri Abhyankar goto theend1; 11182b4ed282SShri Abhyankar } 11192b4ed282SShri Abhyankar if (*ynorm > vi->maxstep) { /* Step too big, so scale back */ 1120c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1121c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Scaling step by %G old ynorm %G\n",vi->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr); 11222b4ed282SShri Abhyankar } 11232b4ed282SShri Abhyankar ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr); 11242b4ed282SShri Abhyankar *ynorm = vi->maxstep; 11252b4ed282SShri Abhyankar } 11262b4ed282SShri Abhyankar ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 11272b4ed282SShri Abhyankar minlambda = vi->minlambda/rellength; 11282b4ed282SShri Abhyankar /* ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); */ 11292b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 11302b4ed282SShri Abhyankar ierr = VecDot(vi->dpsi,y,&cinitslope);CHKERRQ(ierr); 11312b4ed282SShri Abhyankar initslope = PetscRealPart(cinitslope); 11322b4ed282SShri Abhyankar #else 11332b4ed282SShri Abhyankar ierr = VecDot(vi->dpsi,y,&initslope);CHKERRQ(ierr); 11342b4ed282SShri Abhyankar #endif 11352b4ed282SShri Abhyankar if (initslope > 0.0) initslope = -initslope; 11362b4ed282SShri Abhyankar if (initslope == 0.0) initslope = -1.0; 11372b4ed282SShri Abhyankar 11382b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 11392b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 11402b4ed282SShri Abhyankar ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 11412b4ed282SShri Abhyankar *flag = PETSC_FALSE; 11422b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 11432b4ed282SShri Abhyankar goto theend1; 11442b4ed282SShri Abhyankar } 11452b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 11462b4ed282SShri Abhyankar if (snes->domainerror) { 11472b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 11482b4ed282SShri Abhyankar PetscFunctionReturn(0); 11492b4ed282SShri Abhyankar } 11502b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 11512b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 11522b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 11532b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */ 1154c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1155c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 11562b4ed282SShri Abhyankar } 11572b4ed282SShri Abhyankar goto theend1; 11582b4ed282SShri Abhyankar } 11592b4ed282SShri Abhyankar 11602b4ed282SShri Abhyankar /* Fit points with quadratic */ 11612b4ed282SShri Abhyankar lambda = 1.0; 11622b4ed282SShri Abhyankar lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 11632b4ed282SShri Abhyankar lambdaprev = lambda; 11642b4ed282SShri Abhyankar gnormprev = *gnorm; 11652b4ed282SShri Abhyankar if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 11662b4ed282SShri Abhyankar if (lambdatemp <= .1*lambda) lambda = .1*lambda; 11672b4ed282SShri Abhyankar else lambda = lambdatemp; 11682b4ed282SShri Abhyankar 11692b4ed282SShri Abhyankar ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 11702b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 11712b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 11722b4ed282SShri Abhyankar *flag = PETSC_FALSE; 11732b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 11742b4ed282SShri Abhyankar goto theend1; 11752b4ed282SShri Abhyankar } 11762b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 11772b4ed282SShri Abhyankar if (snes->domainerror) { 11782b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 11792b4ed282SShri Abhyankar PetscFunctionReturn(0); 11802b4ed282SShri Abhyankar } 11812b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 11822b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1183c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1184c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr); 11852b4ed282SShri Abhyankar } 11862b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */ 1187c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1188c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 11892b4ed282SShri Abhyankar } 11902b4ed282SShri Abhyankar goto theend1; 11912b4ed282SShri Abhyankar } 11922b4ed282SShri Abhyankar 11932b4ed282SShri Abhyankar /* Fit points with cubic */ 11942b4ed282SShri Abhyankar count = 1; 11952b4ed282SShri Abhyankar while (PETSC_TRUE) { 11962b4ed282SShri Abhyankar if (lambda <= minlambda) { 1197c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1198c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 1199c92abb8aSShri 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); 12002b4ed282SShri Abhyankar } 12012b4ed282SShri Abhyankar *flag = PETSC_FALSE; 12022b4ed282SShri Abhyankar break; 12032b4ed282SShri Abhyankar } 12042b4ed282SShri Abhyankar t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 12052b4ed282SShri Abhyankar t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 12062b4ed282SShri Abhyankar a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 12072b4ed282SShri Abhyankar b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 12082b4ed282SShri Abhyankar d = b*b - 3*a*initslope; 12092b4ed282SShri Abhyankar if (d < 0.0) d = 0.0; 12102b4ed282SShri Abhyankar if (a == 0.0) { 12112b4ed282SShri Abhyankar lambdatemp = -initslope/(2.0*b); 12122b4ed282SShri Abhyankar } else { 12132b4ed282SShri Abhyankar lambdatemp = (-b + sqrt(d))/(3.0*a); 12142b4ed282SShri Abhyankar } 12152b4ed282SShri Abhyankar lambdaprev = lambda; 12162b4ed282SShri Abhyankar gnormprev = *gnorm; 12172b4ed282SShri Abhyankar if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 12182b4ed282SShri Abhyankar if (lambdatemp <= .1*lambda) lambda = .1*lambda; 12192b4ed282SShri Abhyankar else lambda = lambdatemp; 12202b4ed282SShri Abhyankar ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 12212b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 12222b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 12232b4ed282SShri 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); 12242b4ed282SShri Abhyankar *flag = PETSC_FALSE; 12252b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 12262b4ed282SShri Abhyankar break; 12272b4ed282SShri Abhyankar } 12282b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 12292b4ed282SShri Abhyankar if (snes->domainerror) { 12302b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 12312b4ed282SShri Abhyankar PetscFunctionReturn(0); 12322b4ed282SShri Abhyankar } 12332b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 12342b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 12352b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* is reduction enough? */ 1236c92abb8aSShri Abhyankar if (vi->lsmonitor) { 12372b4ed282SShri Abhyankar ierr = PetscPrintf(comm," Line search: Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 12382b4ed282SShri Abhyankar } 12392b4ed282SShri Abhyankar break; 12402b4ed282SShri Abhyankar } else { 1241c92abb8aSShri Abhyankar if (vi->lsmonitor) { 12422b4ed282SShri Abhyankar ierr = PetscPrintf(comm," Line search: Cubic step no good, shrinking lambda, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 12432b4ed282SShri Abhyankar } 12442b4ed282SShri Abhyankar } 12452b4ed282SShri Abhyankar count++; 12462b4ed282SShri Abhyankar } 12472b4ed282SShri Abhyankar theend1: 12482b4ed282SShri Abhyankar /* Optional user-defined check for line search step validity */ 12492b4ed282SShri Abhyankar if (vi->postcheckstep && *flag) { 12502b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 12512b4ed282SShri Abhyankar if (changed_y) { 12522b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 12532b4ed282SShri Abhyankar } 12542b4ed282SShri Abhyankar if (changed_y || changed_w) { /* recompute the function if the step has changed */ 12552b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 12562b4ed282SShri Abhyankar if (snes->domainerror) { 12572b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 12582b4ed282SShri Abhyankar PetscFunctionReturn(0); 12592b4ed282SShri Abhyankar } 12602b4ed282SShri Abhyankar ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 12612b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 12622b4ed282SShri Abhyankar ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 12632b4ed282SShri Abhyankar ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 12642b4ed282SShri Abhyankar ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 12652b4ed282SShri Abhyankar } 12662b4ed282SShri Abhyankar } 12672b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 12682b4ed282SShri Abhyankar PetscFunctionReturn(0); 12692b4ed282SShri Abhyankar } 12702b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 12712b4ed282SShri Abhyankar #undef __FUNCT__ 12722b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchQuadratic_VI" 12732b4ed282SShri Abhyankar /* 12742b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchQuadratic in snes/impls/ls/ls.c 12752b4ed282SShri Abhyankar */ 1276ace3abfcSBarry 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) 12772b4ed282SShri Abhyankar { 12782b4ed282SShri Abhyankar /* 12792b4ed282SShri Abhyankar Note that for line search purposes we work with with the related 12802b4ed282SShri Abhyankar minimization problem: 12812b4ed282SShri Abhyankar min z(x): R^n -> R, 12822b4ed282SShri Abhyankar where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 12832b4ed282SShri Abhyankar z(x) is the same as the merit function for the complementarity problem 12842b4ed282SShri Abhyankar */ 12852b4ed282SShri Abhyankar PetscReal initslope,minlambda,lambda,lambdatemp,rellength; 12862b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 12872b4ed282SShri Abhyankar PetscScalar cinitslope; 12882b4ed282SShri Abhyankar #endif 12892b4ed282SShri Abhyankar PetscErrorCode ierr; 12902b4ed282SShri Abhyankar PetscInt count; 12912b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1292ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 12932b4ed282SShri Abhyankar 12942b4ed282SShri Abhyankar PetscFunctionBegin; 12952b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 12962b4ed282SShri Abhyankar *flag = PETSC_TRUE; 12972b4ed282SShri Abhyankar 12982b4ed282SShri Abhyankar ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 12992b4ed282SShri Abhyankar if (*ynorm == 0.0) { 1300c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1301c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr); 13022b4ed282SShri Abhyankar } 13032b4ed282SShri Abhyankar *gnorm = fnorm; 13042b4ed282SShri Abhyankar ierr = VecCopy(x,w);CHKERRQ(ierr); 13052b4ed282SShri Abhyankar ierr = VecCopy(f,g);CHKERRQ(ierr); 13062b4ed282SShri Abhyankar *flag = PETSC_FALSE; 13072b4ed282SShri Abhyankar goto theend2; 13082b4ed282SShri Abhyankar } 13092b4ed282SShri Abhyankar if (*ynorm > vi->maxstep) { /* Step too big, so scale back */ 13102b4ed282SShri Abhyankar ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr); 13112b4ed282SShri Abhyankar *ynorm = vi->maxstep; 13122b4ed282SShri Abhyankar } 13132b4ed282SShri Abhyankar ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 13142b4ed282SShri Abhyankar minlambda = vi->minlambda/rellength; 13152b4ed282SShri Abhyankar /* ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); */ 13162b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 13172b4ed282SShri Abhyankar ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 13182b4ed282SShri Abhyankar initslope = PetscRealPart(cinitslope); 13192b4ed282SShri Abhyankar #else 13202b4ed282SShri Abhyankar ierr = VecDot(vi->dpsi,y,&initslope);CHKERRQ(ierr); 13212b4ed282SShri Abhyankar #endif 13222b4ed282SShri Abhyankar if (initslope > 0.0) initslope = -initslope; 13232b4ed282SShri Abhyankar if (initslope == 0.0) initslope = -1.0; 13242b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr); 13252b4ed282SShri Abhyankar 13262b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 13272b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 13282b4ed282SShri Abhyankar ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 13292b4ed282SShri Abhyankar *flag = PETSC_FALSE; 13302b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 13312b4ed282SShri Abhyankar goto theend2; 13322b4ed282SShri Abhyankar } 13332b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 13342b4ed282SShri Abhyankar if (snes->domainerror) { 13352b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 13362b4ed282SShri Abhyankar PetscFunctionReturn(0); 13372b4ed282SShri Abhyankar } 13382b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 13392b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 13402b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */ 1341c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1342c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 13432b4ed282SShri Abhyankar } 13442b4ed282SShri Abhyankar goto theend2; 13452b4ed282SShri Abhyankar } 13462b4ed282SShri Abhyankar 13472b4ed282SShri Abhyankar /* Fit points with quadratic */ 13482b4ed282SShri Abhyankar lambda = 1.0; 13492b4ed282SShri Abhyankar count = 1; 13502b4ed282SShri Abhyankar while (PETSC_TRUE) { 13512b4ed282SShri Abhyankar if (lambda <= minlambda) { /* bad luck; use full step */ 1352c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1353c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr); 1354c92abb8aSShri 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); 13552b4ed282SShri Abhyankar } 13562b4ed282SShri Abhyankar ierr = VecCopy(x,w);CHKERRQ(ierr); 13572b4ed282SShri Abhyankar *flag = PETSC_FALSE; 13582b4ed282SShri Abhyankar break; 13592b4ed282SShri Abhyankar } 13602b4ed282SShri Abhyankar lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 13612b4ed282SShri Abhyankar if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 13622b4ed282SShri Abhyankar if (lambdatemp <= .1*lambda) lambda = .1*lambda; 13632b4ed282SShri Abhyankar else lambda = lambdatemp; 13642b4ed282SShri Abhyankar 13652b4ed282SShri Abhyankar ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 13662b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 13672b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 13682b4ed282SShri 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); 13692b4ed282SShri Abhyankar *flag = PETSC_FALSE; 13702b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 13712b4ed282SShri Abhyankar break; 13722b4ed282SShri Abhyankar } 13732b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 13742b4ed282SShri Abhyankar if (snes->domainerror) { 13752b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 13762b4ed282SShri Abhyankar PetscFunctionReturn(0); 13772b4ed282SShri Abhyankar } 13782b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 13792b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 13802b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */ 1381c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1382c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line Search: Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr); 13832b4ed282SShri Abhyankar } 13842b4ed282SShri Abhyankar break; 13852b4ed282SShri Abhyankar } 13862b4ed282SShri Abhyankar count++; 13872b4ed282SShri Abhyankar } 13882b4ed282SShri Abhyankar theend2: 13892b4ed282SShri Abhyankar /* Optional user-defined check for line search step validity */ 13902b4ed282SShri Abhyankar if (vi->postcheckstep) { 13912b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 13922b4ed282SShri Abhyankar if (changed_y) { 13932b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 13942b4ed282SShri Abhyankar } 13952b4ed282SShri Abhyankar if (changed_y || changed_w) { /* recompute the function if the step has changed */ 13962b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g); 13972b4ed282SShri Abhyankar if (snes->domainerror) { 13982b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 13992b4ed282SShri Abhyankar PetscFunctionReturn(0); 14002b4ed282SShri Abhyankar } 14012b4ed282SShri Abhyankar ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 14022b4ed282SShri Abhyankar ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 14032b4ed282SShri Abhyankar ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 14042b4ed282SShri Abhyankar ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 14052b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 14062b4ed282SShri Abhyankar } 14072b4ed282SShri Abhyankar } 14082b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 14092b4ed282SShri Abhyankar PetscFunctionReturn(0); 14102b4ed282SShri Abhyankar } 14112b4ed282SShri Abhyankar 1412ace3abfcSBarry 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*/ 14132b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 14142b4ed282SShri Abhyankar EXTERN_C_BEGIN 14152b4ed282SShri Abhyankar #undef __FUNCT__ 14162b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchSet_VI" 14172b4ed282SShri Abhyankar PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_VI(SNES snes,FCN2 func,void *lsctx) 14182b4ed282SShri Abhyankar { 14192b4ed282SShri Abhyankar PetscFunctionBegin; 14202b4ed282SShri Abhyankar ((SNES_VI *)(snes->data))->LineSearch = func; 14212b4ed282SShri Abhyankar ((SNES_VI *)(snes->data))->lsP = lsctx; 14222b4ed282SShri Abhyankar PetscFunctionReturn(0); 14232b4ed282SShri Abhyankar } 14242b4ed282SShri Abhyankar EXTERN_C_END 14252b4ed282SShri Abhyankar 14262b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 14272b4ed282SShri Abhyankar EXTERN_C_BEGIN 14282b4ed282SShri Abhyankar #undef __FUNCT__ 14292b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchSetMonitor_VI" 1430ace3abfcSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetMonitor_VI(SNES snes,PetscBool flg) 14312b4ed282SShri Abhyankar { 14322b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 14332b4ed282SShri Abhyankar PetscErrorCode ierr; 14342b4ed282SShri Abhyankar 14352b4ed282SShri Abhyankar PetscFunctionBegin; 1436c92abb8aSShri Abhyankar if (flg && !vi->lsmonitor) { 1437c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",((PetscObject)snes)->tablevel,&vi->lsmonitor);CHKERRQ(ierr); 1438c92abb8aSShri Abhyankar } else if (!flg && vi->lsmonitor) { 1439c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr); 14402b4ed282SShri Abhyankar } 14412b4ed282SShri Abhyankar PetscFunctionReturn(0); 14422b4ed282SShri Abhyankar } 14432b4ed282SShri Abhyankar EXTERN_C_END 14442b4ed282SShri Abhyankar 14452b4ed282SShri Abhyankar /* 14462b4ed282SShri Abhyankar SNESView_VI - Prints info from the SNESVI data structure. 14472b4ed282SShri Abhyankar 14482b4ed282SShri Abhyankar Input Parameters: 14492b4ed282SShri Abhyankar . SNES - the SNES context 14502b4ed282SShri Abhyankar . viewer - visualization context 14512b4ed282SShri Abhyankar 14522b4ed282SShri Abhyankar Application Interface Routine: SNESView() 14532b4ed282SShri Abhyankar */ 14542b4ed282SShri Abhyankar #undef __FUNCT__ 14552b4ed282SShri Abhyankar #define __FUNCT__ "SNESView_VI" 14562b4ed282SShri Abhyankar static PetscErrorCode SNESView_VI(SNES snes,PetscViewer viewer) 14572b4ed282SShri Abhyankar { 14582b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI *)snes->data; 14592b4ed282SShri Abhyankar const char *cstr; 14602b4ed282SShri Abhyankar PetscErrorCode ierr; 1461ace3abfcSBarry Smith PetscBool iascii; 14622b4ed282SShri Abhyankar 14632b4ed282SShri Abhyankar PetscFunctionBegin; 14642b4ed282SShri Abhyankar ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 14652b4ed282SShri Abhyankar if (iascii) { 14662b4ed282SShri Abhyankar if (vi->LineSearch == SNESLineSearchNo_VI) cstr = "SNESLineSearchNo"; 14672b4ed282SShri Abhyankar else if (vi->LineSearch == SNESLineSearchQuadratic_VI) cstr = "SNESLineSearchQuadratic"; 14682b4ed282SShri Abhyankar else if (vi->LineSearch == SNESLineSearchCubic_VI) cstr = "SNESLineSearchCubic"; 14692b4ed282SShri Abhyankar else cstr = "unknown"; 14702b4ed282SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 14712b4ed282SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," alpha=%G, maxstep=%G, minlambda=%G\n",vi->alpha,vi->maxstep,vi->minlambda);CHKERRQ(ierr); 14722b4ed282SShri Abhyankar } else { 14732b4ed282SShri Abhyankar SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ VI",((PetscObject)viewer)->type_name); 14742b4ed282SShri Abhyankar } 14752b4ed282SShri Abhyankar PetscFunctionReturn(0); 14762b4ed282SShri Abhyankar } 14772b4ed282SShri Abhyankar 14782b4ed282SShri Abhyankar /* 14792b4ed282SShri Abhyankar SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. 14802b4ed282SShri Abhyankar 14812b4ed282SShri Abhyankar Input Parameters: 14822b4ed282SShri Abhyankar . snes - the SNES context. 14832b4ed282SShri Abhyankar . xl - lower bound. 14842b4ed282SShri Abhyankar . xu - upper bound. 14852b4ed282SShri Abhyankar 14862b4ed282SShri Abhyankar Notes: 14872b4ed282SShri Abhyankar If this routine is not called then the lower and upper bounds are set to 14882b4ed282SShri Abhyankar -Infinity and Infinity respectively during SNESSetUp. 14892b4ed282SShri Abhyankar */ 14902b4ed282SShri Abhyankar 14912b4ed282SShri Abhyankar #undef __FUNCT__ 14922b4ed282SShri Abhyankar #define __FUNCT__ "SNESVISetVariableBounds" 149369c03620SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu) 14942b4ed282SShri Abhyankar { 14952b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 14962b4ed282SShri Abhyankar 14972b4ed282SShri Abhyankar PetscFunctionBegin; 14982b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 14992b4ed282SShri Abhyankar PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 15002b4ed282SShri Abhyankar PetscValidHeaderSpecific(xu,VEC_CLASSID,3); 15012b4ed282SShri Abhyankar 15022b4ed282SShri Abhyankar /* Check if SNESSetFunction is called */ 15032b4ed282SShri Abhyankar if(!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 15042b4ed282SShri Abhyankar 15052b4ed282SShri Abhyankar /* Check if the vector sizes are compatible for lower and upper bounds */ 15062b4ed282SShri 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); 15072b4ed282SShri 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); 15082b4ed282SShri Abhyankar vi->usersetxbounds = PETSC_TRUE; 15092b4ed282SShri Abhyankar vi->xl = xl; 15102b4ed282SShri Abhyankar vi->xu = xu; 15112b4ed282SShri Abhyankar 15122b4ed282SShri Abhyankar PetscFunctionReturn(0); 15132b4ed282SShri Abhyankar } 15142b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 15152b4ed282SShri Abhyankar /* 15162b4ed282SShri Abhyankar SNESSetFromOptions_VI - Sets various parameters for the SNESVI method. 15172b4ed282SShri Abhyankar 15182b4ed282SShri Abhyankar Input Parameter: 15192b4ed282SShri Abhyankar . snes - the SNES context 15202b4ed282SShri Abhyankar 15212b4ed282SShri Abhyankar Application Interface Routine: SNESSetFromOptions() 15222b4ed282SShri Abhyankar */ 15232b4ed282SShri Abhyankar #undef __FUNCT__ 15242b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetFromOptions_VI" 15252b4ed282SShri Abhyankar static PetscErrorCode SNESSetFromOptions_VI(SNES snes) 15262b4ed282SShri Abhyankar { 15272b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI *)snes->data; 15282b4ed282SShri Abhyankar const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 152969c03620SShri Abhyankar const char *vies[] = {"ss","as"}; 15302b4ed282SShri Abhyankar PetscErrorCode ierr; 15312b4ed282SShri Abhyankar PetscInt indx; 153269c03620SShri Abhyankar PetscBool flg,set,flg2; 15332b4ed282SShri Abhyankar 15342b4ed282SShri Abhyankar PetscFunctionBegin; 15352b4ed282SShri Abhyankar ierr = PetscOptionsHead("SNES semismooth method options");CHKERRQ(ierr); 15362b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_alpha","Function norm must decrease by","None",vi->alpha,&vi->alpha,0);CHKERRQ(ierr); 15372b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_maxstep","Step must be less than","None",vi->maxstep,&vi->maxstep,0);CHKERRQ(ierr); 15382b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_minlambda","Minimum lambda allowed","None",vi->minlambda,&vi->minlambda,0);CHKERRQ(ierr); 15392b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_delta","descent test fraction","None",vi->delta,&vi->delta,0);CHKERRQ(ierr); 15402b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_rho","descent test power","None",vi->rho,&vi->rho,0);CHKERRQ(ierr); 15412b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_const_tol","constraint tolerance","None",vi->const_tol,&vi->const_tol,0);CHKERRQ(ierr); 1542acfcf0e5SJed Brown ierr = PetscOptionsBool("-snes_vi_lsmonitor","Print progress of line searches","SNESLineSearchSetMonitor",vi->lsmonitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 15432b4ed282SShri Abhyankar if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);} 154469c03620SShri Abhyankar ierr = PetscOptionsEList("-snes_vi_type","Semismooth algorithm used","",vies,2,"ss",&indx,&flg2);CHKERRQ(ierr); 154569c03620SShri Abhyankar if (flg2) { 154669c03620SShri Abhyankar switch (indx) { 154769c03620SShri Abhyankar case 0: 154869c03620SShri Abhyankar snes->ops->solve = SNESSolveVI_SS; 154969c03620SShri Abhyankar break; 155069c03620SShri Abhyankar case 1: 155169c03620SShri Abhyankar snes->ops->solve = SNESSolveVI_AS; 155269c03620SShri Abhyankar break; 155369c03620SShri Abhyankar } 155469c03620SShri Abhyankar } 1555c92abb8aSShri Abhyankar ierr = PetscOptionsEList("-snes_vi_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 15562b4ed282SShri Abhyankar if (flg) { 15572b4ed282SShri Abhyankar switch (indx) { 15582b4ed282SShri Abhyankar case 0: 15592b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr); 15602b4ed282SShri Abhyankar break; 15612b4ed282SShri Abhyankar case 1: 15622b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr); 15632b4ed282SShri Abhyankar break; 15642b4ed282SShri Abhyankar case 2: 15652b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr); 15662b4ed282SShri Abhyankar break; 15672b4ed282SShri Abhyankar case 3: 15682b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr); 15692b4ed282SShri Abhyankar break; 15702b4ed282SShri Abhyankar } 15712b4ed282SShri Abhyankar } 15722b4ed282SShri Abhyankar ierr = PetscOptionsTail();CHKERRQ(ierr); 15732b4ed282SShri Abhyankar PetscFunctionReturn(0); 15742b4ed282SShri Abhyankar } 15752b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 15762b4ed282SShri Abhyankar /*MC 15772b4ed282SShri Abhyankar SNESVI - Semismooth newton method based nonlinear solver that uses a line search 15782b4ed282SShri Abhyankar 15792b4ed282SShri Abhyankar Options Database: 15802b4ed282SShri Abhyankar + -snes_vi [cubic,quadratic,basic,basicnonorms] - Selects line search 15812b4ed282SShri Abhyankar . -snes_vi_alpha <alpha> - Sets alpha 15822b4ed282SShri 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) 15832b4ed282SShri Abhyankar . -snes_vi_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 15842b4ed282SShri Abhyankar -snes_vi_delta <delta> - Sets the fraction used in the descent test. 15852b4ed282SShri Abhyankar -snes_vi_rho <rho> - Sets the power used in the descent test. 15862b4ed282SShri Abhyankar For a descent direction to be accepted it has to satisfy the condition dpsi^T*d <= -delta*||d||^rho 15872b4ed282SShri Abhyankar - -snes_vi_monitor - print information about progress of line searches 15882b4ed282SShri Abhyankar 15892b4ed282SShri Abhyankar 15902b4ed282SShri Abhyankar Level: beginner 15912b4ed282SShri Abhyankar 15922b4ed282SShri Abhyankar .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 15932b4ed282SShri Abhyankar SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 15942b4ed282SShri Abhyankar SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams() 15952b4ed282SShri Abhyankar 15962b4ed282SShri Abhyankar M*/ 15972b4ed282SShri Abhyankar EXTERN_C_BEGIN 15982b4ed282SShri Abhyankar #undef __FUNCT__ 15992b4ed282SShri Abhyankar #define __FUNCT__ "SNESCreate_VI" 16002b4ed282SShri Abhyankar PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_VI(SNES snes) 16012b4ed282SShri Abhyankar { 16022b4ed282SShri Abhyankar PetscErrorCode ierr; 16032b4ed282SShri Abhyankar SNES_VI *vi; 16042b4ed282SShri Abhyankar 16052b4ed282SShri Abhyankar PetscFunctionBegin; 16062b4ed282SShri Abhyankar snes->ops->setup = SNESSetUp_VI; 160769c03620SShri Abhyankar snes->ops->solve = SNESSolveVI_SS; 16082b4ed282SShri Abhyankar snes->ops->destroy = SNESDestroy_VI; 16092b4ed282SShri Abhyankar snes->ops->setfromoptions = SNESSetFromOptions_VI; 16102b4ed282SShri Abhyankar snes->ops->view = SNESView_VI; 16112b4ed282SShri Abhyankar snes->ops->converged = SNESDefaultConverged_VI; 16122b4ed282SShri Abhyankar 16132b4ed282SShri Abhyankar ierr = PetscNewLog(snes,SNES_VI,&vi);CHKERRQ(ierr); 16142b4ed282SShri Abhyankar snes->data = (void*)vi; 16152b4ed282SShri Abhyankar vi->alpha = 1.e-4; 16162b4ed282SShri Abhyankar vi->maxstep = 1.e8; 16172b4ed282SShri Abhyankar vi->minlambda = 1.e-12; 16182b4ed282SShri Abhyankar vi->LineSearch = SNESLineSearchCubic_VI; 16192b4ed282SShri Abhyankar vi->lsP = PETSC_NULL; 16202b4ed282SShri Abhyankar vi->postcheckstep = PETSC_NULL; 16212b4ed282SShri Abhyankar vi->postcheck = PETSC_NULL; 16222b4ed282SShri Abhyankar vi->precheckstep = PETSC_NULL; 16232b4ed282SShri Abhyankar vi->precheck = PETSC_NULL; 16242b4ed282SShri Abhyankar vi->rho = 2.1; 16252b4ed282SShri Abhyankar vi->delta = 1e-10; 16262b4ed282SShri Abhyankar vi->const_tol = 2.2204460492503131e-16; 16272b4ed282SShri Abhyankar vi->computessfunction = ComputeFischerFunction; 16282b4ed282SShri Abhyankar 16292b4ed282SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_VI",SNESLineSearchSetMonitor_VI);CHKERRQ(ierr); 16302b4ed282SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_VI",SNESLineSearchSet_VI);CHKERRQ(ierr); 16312b4ed282SShri Abhyankar 16322b4ed282SShri Abhyankar PetscFunctionReturn(0); 16332b4ed282SShri Abhyankar } 16342b4ed282SShri Abhyankar EXTERN_C_END 1635