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 /* 2401f79c099SShri 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. 2412b4ed282SShri Abhyankar 2422b4ed282SShri Abhyankar Input Parameters: 2432b4ed282SShri Abhyankar . snes - the SNES context 2442b4ed282SShri Abhyankar . X - the current iterate 2452b4ed282SShri Abhyankar . vec_func - nonlinear function evaluated at x 2462b4ed282SShri Abhyankar 2472b4ed282SShri Abhyankar Output Parameters: 2482b4ed282SShri Abhyankar . jac - semismooth jacobian 2492b4ed282SShri Abhyankar . jac_pre - optional preconditioning matrix 2502b4ed282SShri Abhyankar . flag - flag passed on by SNESComputeJacobian. 2512b4ed282SShri Abhyankar . jacctx - user provided jacobian context 2522b4ed282SShri Abhyankar 2532b4ed282SShri Abhyankar Notes: 2542b4ed282SShri Abhyankar The semismooth jacobian matrix is given by 2552b4ed282SShri Abhyankar jac = Da + Db*jacfun 2562b4ed282SShri Abhyankar where Db is the row scaling matrix stored as a vector, 2572b4ed282SShri Abhyankar Da is the diagonal perturbation matrix stored as a vector 2582b4ed282SShri Abhyankar and jacfun is the jacobian of the original nonlinear function. 2592b4ed282SShri Abhyankar */ 2602b4ed282SShri Abhyankar #undef __FUNCT__ 2611f79c099SShri Abhyankar #define __FUNCT__ "SNESVIComputeJacobian" 2621f79c099SShri Abhyankar PetscErrorCode SNESVIComputeJacobian(SNES snes,Vec X,Mat *jac, Mat *jac_pre, MatStructure *flg,void* jacctx) 2632b4ed282SShri Abhyankar { 2642b4ed282SShri Abhyankar PetscErrorCode ierr; 2652b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 2662b4ed282SShri Abhyankar PetscScalar *l,*u,*x,*f,*da,*db,*z,*t,t1,t2,ci,di,ei; 2672b4ed282SShri Abhyankar PetscInt i,nlocal; 2682b4ed282SShri Abhyankar Vec F = snes->vec_func; 2692b4ed282SShri Abhyankar 2702b4ed282SShri Abhyankar PetscFunctionBegin; 2712b4ed282SShri Abhyankar 2722b4ed282SShri Abhyankar ierr = (*vi->computeuserjacobian)(snes,X,jac,jac_pre,flg,jacctx);CHKERRQ(ierr); 2732b4ed282SShri Abhyankar 2742b4ed282SShri Abhyankar ierr = VecGetArray(X,&x);CHKERRQ(ierr); 2752b4ed282SShri Abhyankar ierr = VecGetArray(F,&f);CHKERRQ(ierr); 2762b4ed282SShri Abhyankar ierr = VecGetArray(vi->xl,&l);CHKERRQ(ierr); 2772b4ed282SShri Abhyankar ierr = VecGetArray(vi->xu,&u);CHKERRQ(ierr); 2782b4ed282SShri Abhyankar ierr = VecGetArray(vi->Da,&da);CHKERRQ(ierr); 2792b4ed282SShri Abhyankar ierr = VecGetArray(vi->Db,&db);CHKERRQ(ierr); 2802b4ed282SShri Abhyankar ierr = VecGetArray(vi->z,&z);CHKERRQ(ierr); 2812b4ed282SShri Abhyankar 2822b4ed282SShri Abhyankar ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 2832b4ed282SShri Abhyankar /* Set the elements of the vector z: 2842b4ed282SShri Abhyankar z[i] = 1 if (x[i] - l[i],f[i]) = (0,0) or (u[i] - x[i],f[i]) = (0,0) 2852b4ed282SShri Abhyankar else z[i] = 0 2862b4ed282SShri Abhyankar */ 2872b4ed282SShri Abhyankar for(i=0;i < nlocal;i++) { 2882b4ed282SShri Abhyankar da[i] = db[i] = z[i] = 0; 2892b4ed282SShri Abhyankar if(PetscAbsScalar(f[i]) <= vi->const_tol) { 2902b4ed282SShri Abhyankar if ((l[i] > PETSC_VI_NINF) && (PetscAbsScalar(x[i]-l[i]) <= vi->const_tol)) { 2912b4ed282SShri Abhyankar da[i] = 1; 2922b4ed282SShri Abhyankar z[i] = 1; 2932b4ed282SShri Abhyankar } 2942b4ed282SShri Abhyankar if ((u[i] < PETSC_VI_INF) && (PetscAbsScalar(u[i]-x[i]) <= vi->const_tol)) { 2952b4ed282SShri Abhyankar db[i] = 1; 2962b4ed282SShri Abhyankar z[i] = 1; 2972b4ed282SShri Abhyankar } 2982b4ed282SShri Abhyankar } 2992b4ed282SShri Abhyankar } 3002b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->z,&z);CHKERRQ(ierr); 3012b4ed282SShri Abhyankar ierr = MatMult(*jac,vi->z,vi->t);CHKERRQ(ierr); 3022b4ed282SShri Abhyankar ierr = VecGetArray(vi->t,&t);CHKERRQ(ierr); 3032b4ed282SShri Abhyankar /* Compute the elements of the diagonal perturbation vector Da and row scaling vector Db */ 3042b4ed282SShri Abhyankar for(i=0;i< nlocal;i++) { 3052b4ed282SShri Abhyankar /* Free variables */ 3062b4ed282SShri Abhyankar if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) { 3072b4ed282SShri Abhyankar da[i] = 0; db[i] = -1; 3082b4ed282SShri Abhyankar } 3092b4ed282SShri Abhyankar /* lower bounded variables */ 3102b4ed282SShri Abhyankar else if (u[i] >= PETSC_VI_INF) { 3112b4ed282SShri Abhyankar if (da[i] >= 1) { 3122b4ed282SShri Abhyankar t2 = PetscScalarNorm(1,t[i]); 3132b4ed282SShri Abhyankar da[i] = 1/t2 - 1; 3142b4ed282SShri Abhyankar db[i] = t[i]/t2 - 1; 3152b4ed282SShri Abhyankar } else { 3162b4ed282SShri Abhyankar t1 = x[i] - l[i]; 3172b4ed282SShri Abhyankar t2 = PetscScalarNorm(t1,f[i]); 3182b4ed282SShri Abhyankar da[i] = t1/t2 - 1; 3192b4ed282SShri Abhyankar db[i] = f[i]/t2 - 1; 3202b4ed282SShri Abhyankar } 3212b4ed282SShri Abhyankar } 3222b4ed282SShri Abhyankar /* upper bounded variables */ 3232b4ed282SShri Abhyankar else if (l[i] <= PETSC_VI_NINF) { 3242b4ed282SShri Abhyankar if (db[i] >= 1) { 3252b4ed282SShri Abhyankar t2 = PetscScalarNorm(1,t[i]); 3262b4ed282SShri Abhyankar da[i] = -1/t2 -1; 3272b4ed282SShri Abhyankar db[i] = -t[i]/t2 - 1; 3282b4ed282SShri Abhyankar } 3292b4ed282SShri Abhyankar else { 3302b4ed282SShri Abhyankar t1 = u[i] - x[i]; 3312b4ed282SShri Abhyankar t2 = PetscScalarNorm(t1,f[i]); 3322b4ed282SShri Abhyankar da[i] = t1/t2 - 1; 3332b4ed282SShri Abhyankar db[i] = -f[i]/t2 - 1; 3342b4ed282SShri Abhyankar } 3352b4ed282SShri Abhyankar } 3362b4ed282SShri Abhyankar /* Fixed variables */ 3372b4ed282SShri Abhyankar else if (l[i] == u[i]) { 3382b4ed282SShri Abhyankar da[i] = -1; 3392b4ed282SShri Abhyankar db[i] = 0; 3402b4ed282SShri Abhyankar } 3412b4ed282SShri Abhyankar /* Box constrained variables */ 3422b4ed282SShri Abhyankar else { 3432b4ed282SShri Abhyankar if (db[i] >= 1) { 3442b4ed282SShri Abhyankar t2 = PetscScalarNorm(1,t[i]); 3452b4ed282SShri Abhyankar ci = 1/t2 + 1; 3462b4ed282SShri Abhyankar di = t[i]/t2 + 1; 3472b4ed282SShri Abhyankar } 3482b4ed282SShri Abhyankar else { 3492b4ed282SShri Abhyankar t1 = x[i] - u[i]; 3502b4ed282SShri Abhyankar t2 = PetscScalarNorm(t1,f[i]); 3512b4ed282SShri Abhyankar ci = t1/t2 + 1; 3522b4ed282SShri Abhyankar di = f[i]/t2 + 1; 3532b4ed282SShri Abhyankar } 3542b4ed282SShri Abhyankar 3552b4ed282SShri Abhyankar if (da[i] >= 1) { 3562b4ed282SShri Abhyankar t1 = ci + di*t[i]; 3572b4ed282SShri Abhyankar t2 = PetscScalarNorm(1,t1); 3582b4ed282SShri Abhyankar t1 = t1/t2 - 1; 3592b4ed282SShri Abhyankar t2 = 1/t2 - 1; 3602b4ed282SShri Abhyankar } 3612b4ed282SShri Abhyankar else { 3622b4ed282SShri Abhyankar ierr = ComputeFischerFunction(u[i]-x[i],-f[i],&ei);CHKERRQ(ierr); 3632b4ed282SShri Abhyankar t2 = PetscScalarNorm(x[i]-l[i],ei); 3642b4ed282SShri Abhyankar t1 = ei/t2 - 1; 3652b4ed282SShri Abhyankar t2 = (x[i] - l[i])/t2 - 1; 3662b4ed282SShri Abhyankar } 3672b4ed282SShri Abhyankar 3682b4ed282SShri Abhyankar da[i] = t2 + t1*ci; 3692b4ed282SShri Abhyankar db[i] = t1*di; 3702b4ed282SShri Abhyankar } 3712b4ed282SShri Abhyankar } 3722b4ed282SShri Abhyankar ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 3732b4ed282SShri Abhyankar ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 3742b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->xl,&l);CHKERRQ(ierr); 3752b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->xu,&u);CHKERRQ(ierr); 3762b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->Da,&da);CHKERRQ(ierr); 3772b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->Db,&db);CHKERRQ(ierr); 3782b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->t,&t);CHKERRQ(ierr); 3792b4ed282SShri Abhyankar 3802b4ed282SShri Abhyankar /* Do row scaling and add diagonal perturbation */ 3812b4ed282SShri Abhyankar ierr = MatDiagonalScale(*jac,vi->Db,PETSC_NULL);CHKERRQ(ierr); 3822b4ed282SShri Abhyankar ierr = MatDiagonalSet(*jac,vi->Da,ADD_VALUES);CHKERRQ(ierr); 3832b4ed282SShri Abhyankar if (*jac != *jac_pre) { /* If jac and jac_pre are different */ 3842b4ed282SShri Abhyankar ierr = MatDiagonalScale(*jac_pre,vi->Db,PETSC_NULL); 3852b4ed282SShri Abhyankar ierr = MatDiagonalSet(*jac_pre,vi->Da,ADD_VALUES);CHKERRQ(ierr); 3862b4ed282SShri Abhyankar } 3872b4ed282SShri Abhyankar 3882b4ed282SShri Abhyankar PetscFunctionReturn(0); 3892b4ed282SShri Abhyankar } 3902b4ed282SShri Abhyankar 3912b4ed282SShri Abhyankar /* 3922b4ed282SShri Abhyankar SNESVIComputeMeritFunctionGradient - Computes the gradient of the merit function psi. 3932b4ed282SShri Abhyankar 3942b4ed282SShri Abhyankar Input Parameters: 3952b4ed282SShri Abhyankar . phi - semismooth function. 3962b4ed282SShri Abhyankar . H - semismooth jacobian 3972b4ed282SShri Abhyankar 3982b4ed282SShri Abhyankar Output Parameters: 3992b4ed282SShri Abhyankar . dpsi - merit function gradient 4002b4ed282SShri Abhyankar 4012b4ed282SShri Abhyankar Notes: 4022b4ed282SShri Abhyankar The merit function gradient is computed as follows 4032b4ed282SShri Abhyankar dpsi = H^T*phi 4042b4ed282SShri Abhyankar */ 4052b4ed282SShri Abhyankar #undef __FUNCT__ 4062b4ed282SShri Abhyankar #define __FUNCT__ "SNESVIComputeMeritFunctionGradient" 4072b4ed282SShri Abhyankar PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi) 4082b4ed282SShri Abhyankar { 4092b4ed282SShri Abhyankar PetscErrorCode ierr; 4102b4ed282SShri Abhyankar 4112b4ed282SShri Abhyankar PetscFunctionBegin; 4122b4ed282SShri Abhyankar ierr = MatMultTranspose(H,phi,dpsi); 4132b4ed282SShri Abhyankar 4142b4ed282SShri Abhyankar PetscFunctionReturn(0); 4152b4ed282SShri Abhyankar } 4162b4ed282SShri Abhyankar 4172b4ed282SShri Abhyankar /* 4182b4ed282SShri Abhyankar SNESVISetDescentDirection - Sets the descent direction for the semismooth algorithm 4192b4ed282SShri Abhyankar 4202b4ed282SShri Abhyankar Input Parameters: 4212b4ed282SShri Abhyankar . snes - the SNES context. 4222b4ed282SShri Abhyankar . dpsi - gradient of the merit function. 4232b4ed282SShri Abhyankar 4242b4ed282SShri Abhyankar Output Parameters: 4252b4ed282SShri Abhyankar . flg - PETSC_TRUE if the sufficient descent condition is not satisfied. 4262b4ed282SShri Abhyankar 4272b4ed282SShri Abhyankar Notes: 4282b4ed282SShri Abhyankar The condition for the sufficient descent direction is 4292b4ed282SShri Abhyankar dpsi^T*Y > delta*||Y||^rho 4302b4ed282SShri Abhyankar where rho, delta are as defined in the SNES_VI structure. 4312b4ed282SShri Abhyankar If this condition is satisfied then the existing descent direction i.e. 4322b4ed282SShri 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. 4332b4ed282SShri Abhyankar */ 4342b4ed282SShri Abhyankar #undef __FUNCT__ 4352b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckDescentDirection" 436ace3abfcSBarry Smith PetscErrorCode SNESVICheckDescentDirection(SNES snes,Vec dpsi, Vec Y,PetscBool* flg) 4372b4ed282SShri Abhyankar { 4382b4ed282SShri Abhyankar PetscErrorCode ierr; 4392b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 4402b4ed282SShri Abhyankar PetscScalar dpsidotY; 4412b4ed282SShri Abhyankar PetscReal norm_Y,rhs; 4422b4ed282SShri Abhyankar const PetscReal rho = vi->rho,delta=vi->delta; 4432b4ed282SShri Abhyankar 4442b4ed282SShri Abhyankar PetscFunctionBegin; 4452b4ed282SShri Abhyankar 4462b4ed282SShri Abhyankar *flg = PETSC_FALSE; 4472b4ed282SShri Abhyankar ierr = VecDot(dpsi,Y,&dpsidotY);CHKERRQ(ierr); 4482b4ed282SShri Abhyankar ierr = VecNormBegin(Y,NORM_2,&norm_Y);CHKERRQ(ierr); 4492b4ed282SShri Abhyankar ierr = VecNormEnd(Y,NORM_2,&norm_Y);CHKERRQ(ierr); 4502b4ed282SShri Abhyankar 4512b4ed282SShri Abhyankar rhs = delta*PetscPowScalar(norm_Y,rho); 4522b4ed282SShri Abhyankar 4532b4ed282SShri Abhyankar if (dpsidotY <= rhs) *flg = PETSC_TRUE; 4542b4ed282SShri Abhyankar 4552b4ed282SShri Abhyankar PetscFunctionReturn(0); 4562b4ed282SShri Abhyankar } 4572b4ed282SShri Abhyankar 4582b4ed282SShri Abhyankar /* 4592b4ed282SShri Abhyankar SNESVIAdjustInitialGuess - Readjusts the initial guess to the SNES solver supplied by the user so that the initial guess lies inside the feasible region . 4602b4ed282SShri Abhyankar 4612b4ed282SShri Abhyankar Input Parameters: 4622b4ed282SShri Abhyankar . lb - lower bound. 4632b4ed282SShri Abhyankar . ub - upper bound. 4642b4ed282SShri Abhyankar 4652b4ed282SShri Abhyankar Output Parameters: 4662b4ed282SShri Abhyankar . X - the readjusted initial guess. 4672b4ed282SShri Abhyankar 4682b4ed282SShri Abhyankar Notes: 4692b4ed282SShri Abhyankar The readjusted initial guess X[i] = max(max(min(l[i],X[i]),min(X[i],u[i])),min(l[i],u[i])) 4702b4ed282SShri Abhyankar */ 4712b4ed282SShri Abhyankar #undef __FUNCT__ 4722b4ed282SShri Abhyankar #define __FUNCT__ "SNESVIAdjustInitialGuess" 4732b4ed282SShri Abhyankar PetscErrorCode SNESVIAdjustInitialGuess(Vec X, Vec lb, Vec ub) 4742b4ed282SShri Abhyankar { 4752b4ed282SShri Abhyankar PetscErrorCode ierr; 4762b4ed282SShri Abhyankar PetscInt i,nlocal; 4772b4ed282SShri Abhyankar PetscScalar *x,*l,*u; 4782b4ed282SShri Abhyankar 4792b4ed282SShri Abhyankar PetscFunctionBegin; 4802b4ed282SShri Abhyankar 4812b4ed282SShri Abhyankar ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 4822b4ed282SShri Abhyankar ierr = VecGetArray(X,&x);CHKERRQ(ierr); 4832b4ed282SShri Abhyankar ierr = VecGetArray(lb,&l);CHKERRQ(ierr); 4842b4ed282SShri Abhyankar ierr = VecGetArray(ub,&u);CHKERRQ(ierr); 4852b4ed282SShri Abhyankar 4862b4ed282SShri Abhyankar for(i = 0; i < nlocal; i++) { 4872b4ed282SShri Abhyankar x[i] = PetscMax(PetscMax(PetscMin(x[i],l[i]),PetscMin(x[i],u[i])),PetscMin(l[i],u[i])); 4882b4ed282SShri Abhyankar } 4892b4ed282SShri Abhyankar 4902b4ed282SShri Abhyankar ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 4912b4ed282SShri Abhyankar ierr = VecRestoreArray(lb,&l);CHKERRQ(ierr); 4922b4ed282SShri Abhyankar ierr = VecRestoreArray(ub,&u);CHKERRQ(ierr); 4932b4ed282SShri Abhyankar 4942b4ed282SShri Abhyankar PetscFunctionReturn(0); 4952b4ed282SShri Abhyankar } 4962b4ed282SShri Abhyankar 4972b4ed282SShri Abhyankar /* -------------------------------------------------------------------- 4982b4ed282SShri Abhyankar 4992b4ed282SShri Abhyankar This file implements a semismooth truncated Newton method with a line search, 5002b4ed282SShri Abhyankar for solving a system of nonlinear equations in complementarity form, using the KSP, Vec, 5012b4ed282SShri Abhyankar and Mat interfaces for linear solvers, vectors, and matrices, 5022b4ed282SShri Abhyankar respectively. 5032b4ed282SShri Abhyankar 5042b4ed282SShri Abhyankar The following basic routines are required for each nonlinear solver: 5052b4ed282SShri Abhyankar SNESCreate_XXX() - Creates a nonlinear solver context 5062b4ed282SShri Abhyankar SNESSetFromOptions_XXX() - Sets runtime options 5072b4ed282SShri Abhyankar SNESSolve_XXX() - Solves the nonlinear system 5082b4ed282SShri Abhyankar SNESDestroy_XXX() - Destroys the nonlinear solver context 5092b4ed282SShri Abhyankar The suffix "_XXX" denotes a particular implementation, in this case 5102b4ed282SShri Abhyankar we use _VI (e.g., SNESCreate_VI, SNESSolve_VI) for solving 5112b4ed282SShri Abhyankar systems of nonlinear equations with a line search (LS) method. 5122b4ed282SShri Abhyankar These routines are actually called via the common user interface 5132b4ed282SShri Abhyankar routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 5142b4ed282SShri Abhyankar SNESDestroy(), so the application code interface remains identical 5152b4ed282SShri Abhyankar for all nonlinear solvers. 5162b4ed282SShri Abhyankar 5172b4ed282SShri Abhyankar Another key routine is: 5182b4ed282SShri Abhyankar SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 5192b4ed282SShri Abhyankar by setting data structures and options. The interface routine SNESSetUp() 5202b4ed282SShri Abhyankar is not usually called directly by the user, but instead is called by 5212b4ed282SShri Abhyankar SNESSolve() if necessary. 5222b4ed282SShri Abhyankar 5232b4ed282SShri Abhyankar Additional basic routines are: 5242b4ed282SShri Abhyankar SNESView_XXX() - Prints details of runtime options that 5252b4ed282SShri Abhyankar have actually been used. 5262b4ed282SShri Abhyankar These are called by application codes via the interface routines 5272b4ed282SShri Abhyankar SNESView(). 5282b4ed282SShri Abhyankar 5292b4ed282SShri Abhyankar The various types of solvers (preconditioners, Krylov subspace methods, 5302b4ed282SShri Abhyankar nonlinear solvers, timesteppers) are all organized similarly, so the 5312b4ed282SShri Abhyankar above description applies to these categories also. 5322b4ed282SShri Abhyankar 5332b4ed282SShri Abhyankar -------------------------------------------------------------------- */ 5342b4ed282SShri Abhyankar /* 53569c03620SShri Abhyankar SNESSolveVI_SS - Solves the complementarity problem with a semismooth Newton 5362b4ed282SShri Abhyankar method using a line search. 5372b4ed282SShri Abhyankar 5382b4ed282SShri Abhyankar Input Parameters: 5392b4ed282SShri Abhyankar . snes - the SNES context 5402b4ed282SShri Abhyankar 5412b4ed282SShri Abhyankar Output Parameter: 5422b4ed282SShri Abhyankar . outits - number of iterations until termination 5432b4ed282SShri Abhyankar 5442b4ed282SShri Abhyankar Application Interface Routine: SNESSolve() 5452b4ed282SShri Abhyankar 5462b4ed282SShri Abhyankar Notes: 5472b4ed282SShri Abhyankar This implements essentially a semismooth Newton method with a 5482b4ed282SShri Abhyankar line search. By default a cubic backtracking line search 5492b4ed282SShri Abhyankar is employed, as described in the text "Numerical Methods for 5502b4ed282SShri Abhyankar Unconstrained Optimization and Nonlinear Equations" by Dennis 5512b4ed282SShri Abhyankar and Schnabel. 5522b4ed282SShri Abhyankar */ 5532b4ed282SShri Abhyankar #undef __FUNCT__ 55469c03620SShri Abhyankar #define __FUNCT__ "SNESSolveVI_SS" 55569c03620SShri Abhyankar PetscErrorCode SNESSolveVI_SS(SNES snes) 5562b4ed282SShri Abhyankar { 5572b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 5582b4ed282SShri Abhyankar PetscErrorCode ierr; 5592b4ed282SShri Abhyankar PetscInt maxits,i,lits; 560ace3abfcSBarry Smith PetscBool lssucceed,changedir; 5612b4ed282SShri Abhyankar MatStructure flg = DIFFERENT_NONZERO_PATTERN; 5622b4ed282SShri Abhyankar PetscReal gnorm,xnorm=0,ynorm; 5632b4ed282SShri Abhyankar Vec Y,X,F,G,W; 5642b4ed282SShri Abhyankar KSPConvergedReason kspreason; 5652b4ed282SShri Abhyankar 5662b4ed282SShri Abhyankar PetscFunctionBegin; 5672b4ed282SShri Abhyankar snes->numFailures = 0; 5682b4ed282SShri Abhyankar snes->numLinearSolveFailures = 0; 5692b4ed282SShri Abhyankar snes->reason = SNES_CONVERGED_ITERATING; 5702b4ed282SShri Abhyankar 5712b4ed282SShri Abhyankar maxits = snes->max_its; /* maximum number of iterations */ 5722b4ed282SShri Abhyankar X = snes->vec_sol; /* solution vector */ 5732b4ed282SShri Abhyankar F = snes->vec_func; /* residual vector */ 5742b4ed282SShri Abhyankar Y = snes->work[0]; /* work vectors */ 5752b4ed282SShri Abhyankar G = snes->work[1]; 5762b4ed282SShri Abhyankar W = snes->work[2]; 5772b4ed282SShri Abhyankar 5782b4ed282SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 5792b4ed282SShri Abhyankar snes->iter = 0; 5802b4ed282SShri Abhyankar snes->norm = 0.0; 5812b4ed282SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 5822b4ed282SShri Abhyankar 5832b4ed282SShri Abhyankar ierr = SNESVIAdjustInitialGuess(X,vi->xl,vi->xu);CHKERRQ(ierr); 5842b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr); 5852b4ed282SShri Abhyankar if (snes->domainerror) { 5862b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 5872b4ed282SShri Abhyankar PetscFunctionReturn(0); 5882b4ed282SShri Abhyankar } 5892b4ed282SShri Abhyankar /* Compute Merit function */ 5902b4ed282SShri Abhyankar ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr); 5912b4ed282SShri Abhyankar 5922b4ed282SShri Abhyankar ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 5932b4ed282SShri Abhyankar ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 5942b4ed282SShri Abhyankar if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 5952b4ed282SShri Abhyankar 5962b4ed282SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 5972b4ed282SShri Abhyankar snes->norm = vi->phinorm; 5982b4ed282SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 5992b4ed282SShri Abhyankar SNESLogConvHistory(snes,vi->phinorm,0); 6002b4ed282SShri Abhyankar SNESMonitor(snes,0,vi->phinorm); 6012b4ed282SShri Abhyankar 6022b4ed282SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 6032b4ed282SShri Abhyankar snes->ttol = vi->phinorm*snes->rtol; 6042b4ed282SShri Abhyankar /* test convergence */ 6052b4ed282SShri Abhyankar ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 6062b4ed282SShri Abhyankar if (snes->reason) PetscFunctionReturn(0); 6072b4ed282SShri Abhyankar 6082b4ed282SShri Abhyankar for (i=0; i<maxits; i++) { 6092b4ed282SShri Abhyankar 6102b4ed282SShri Abhyankar /* Call general purpose update function */ 6112b4ed282SShri Abhyankar if (snes->ops->update) { 6122b4ed282SShri Abhyankar ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 6132b4ed282SShri Abhyankar } 6142b4ed282SShri Abhyankar 6152b4ed282SShri Abhyankar /* Solve J Y = Phi, where J is the semismooth jacobian */ 6162b4ed282SShri Abhyankar ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 6172b4ed282SShri Abhyankar 6182b4ed282SShri Abhyankar ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 6192b4ed282SShri Abhyankar ierr = SNES_KSPSolve(snes,snes->ksp,vi->phi,Y);CHKERRQ(ierr); 6202b4ed282SShri Abhyankar ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 6212b4ed282SShri Abhyankar ierr = SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);CHKERRQ(ierr); 6222b4ed282SShri Abhyankar ierr = SNESVICheckDescentDirection(snes,vi->dpsi,Y,&changedir);CHKERRQ(ierr); 6232b4ed282SShri Abhyankar if (kspreason < 0 || changedir) { 6242b4ed282SShri Abhyankar if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 6252b4ed282SShri 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); 6262b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 6272b4ed282SShri Abhyankar break; 6282b4ed282SShri Abhyankar } 6292b4ed282SShri Abhyankar ierr = VecCopy(vi->dpsi,Y);CHKERRQ(ierr); 6302b4ed282SShri Abhyankar } 6312b4ed282SShri Abhyankar ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 6322b4ed282SShri Abhyankar snes->linear_its += lits; 6332b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 6342b4ed282SShri Abhyankar /* 6352b4ed282SShri Abhyankar if (vi->precheckstep) { 636ace3abfcSBarry Smith PetscBool changed_y = PETSC_FALSE; 6372b4ed282SShri Abhyankar ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr); 6382b4ed282SShri Abhyankar } 6392b4ed282SShri Abhyankar 6402b4ed282SShri Abhyankar if (PetscLogPrintInfo){ 6412b4ed282SShri Abhyankar ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 6422b4ed282SShri Abhyankar } 6432b4ed282SShri Abhyankar */ 6442b4ed282SShri Abhyankar /* Compute a (scaled) negative update in the line search routine: 6452b4ed282SShri Abhyankar Y <- X - lambda*Y 6462b4ed282SShri Abhyankar and evaluate G = function(Y) (depends on the line search). 6472b4ed282SShri Abhyankar */ 6482b4ed282SShri Abhyankar ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 6492b4ed282SShri Abhyankar ynorm = 1; gnorm = vi->phinorm; 6502b4ed282SShri Abhyankar ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 6512b4ed282SShri Abhyankar ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 6522b4ed282SShri Abhyankar if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 6532b4ed282SShri Abhyankar if (snes->domainerror) { 6542b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 6552b4ed282SShri Abhyankar PetscFunctionReturn(0); 6562b4ed282SShri Abhyankar } 6572b4ed282SShri Abhyankar if (!lssucceed) { 6582b4ed282SShri Abhyankar if (++snes->numFailures >= snes->maxFailures) { 659ace3abfcSBarry Smith PetscBool ismin; 6602b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_LINE_SEARCH; 6612b4ed282SShri Abhyankar ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 6622b4ed282SShri Abhyankar if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 6632b4ed282SShri Abhyankar break; 6642b4ed282SShri Abhyankar } 6652b4ed282SShri Abhyankar } 6662b4ed282SShri Abhyankar /* Update function and solution vectors */ 6672b4ed282SShri Abhyankar vi->phinorm = gnorm; 6682b4ed282SShri Abhyankar vi->merit = 0.5*vi->phinorm*vi->phinorm; 6692b4ed282SShri Abhyankar ierr = VecCopy(G,vi->phi);CHKERRQ(ierr); 6702b4ed282SShri Abhyankar ierr = VecCopy(W,X);CHKERRQ(ierr); 6712b4ed282SShri Abhyankar /* Monitor convergence */ 6722b4ed282SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 6732b4ed282SShri Abhyankar snes->iter = i+1; 6742b4ed282SShri Abhyankar snes->norm = vi->phinorm; 6752b4ed282SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 6762b4ed282SShri Abhyankar SNESLogConvHistory(snes,snes->norm,lits); 6772b4ed282SShri Abhyankar SNESMonitor(snes,snes->iter,snes->norm); 6782b4ed282SShri Abhyankar /* Test for convergence, xnorm = || X || */ 6792b4ed282SShri Abhyankar if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 6802b4ed282SShri Abhyankar ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 6812b4ed282SShri Abhyankar if (snes->reason) break; 6822b4ed282SShri Abhyankar } 6832b4ed282SShri Abhyankar if (i == maxits) { 6842b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 6852b4ed282SShri Abhyankar if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 6862b4ed282SShri Abhyankar } 6872b4ed282SShri Abhyankar PetscFunctionReturn(0); 6882b4ed282SShri Abhyankar } 68969c03620SShri Abhyankar 69069c03620SShri Abhyankar #undef __FUNCT__ 69169c03620SShri Abhyankar #define __FUNCT__ "SNESVICreateIndexSets_AS" 69269c03620SShri Abhyankar PetscErrorCode SNESVICreateIndexSets_AS(SNES snes,Vec Db,PetscScalar thresh,IS* ISact,IS* ISinact) 69369c03620SShri Abhyankar { 69469c03620SShri Abhyankar PetscErrorCode ierr; 69569c03620SShri Abhyankar PetscInt i,nlocal,ilow,ihigh,nloc_isact=0,nloc_isinact=0; 69669c03620SShri Abhyankar PetscInt *idx_act,*idx_inact,i1=0,i2=0; 69769c03620SShri Abhyankar PetscScalar *db; 69869c03620SShri Abhyankar 69969c03620SShri Abhyankar PetscFunctionBegin; 70069c03620SShri Abhyankar 70169c03620SShri Abhyankar ierr = VecGetLocalSize(Db,&nlocal);CHKERRQ(ierr); 70269c03620SShri Abhyankar ierr = VecGetOwnershipRange(Db,&ilow,&ihigh);CHKERRQ(ierr); 70369c03620SShri Abhyankar ierr = VecGetArray(Db,&db);CHKERRQ(ierr); 70469c03620SShri Abhyankar /* Compute the sizes of the active and inactive sets */ 70569c03620SShri Abhyankar for (i=0; i < nlocal;i++) { 70669c03620SShri Abhyankar if (db[i] <= thresh) nloc_isact++; 70769c03620SShri Abhyankar else nloc_isact++; 70869c03620SShri Abhyankar } 70969c03620SShri Abhyankar 71069c03620SShri Abhyankar ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr); 71169c03620SShri Abhyankar ierr = PetscMalloc(nloc_isinact*sizeof(PetscInt),&idx_inact);CHKERRQ(ierr); 71269c03620SShri Abhyankar 71369c03620SShri Abhyankar /* Creating the indexing arrays */ 71469c03620SShri Abhyankar for(i=ilow; i < ihigh; i++) { 71569c03620SShri Abhyankar if (db[i] <= thresh) idx_act[i1++] = i; 71669c03620SShri Abhyankar else idx_inact[i2++] = i; 71769c03620SShri Abhyankar } 71869c03620SShri Abhyankar 71969c03620SShri Abhyankar /* Create the index sets */ 72069c03620SShri Abhyankar ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isact,idx_act,PETSC_COPY_VALUES,ISact);CHKERRQ(ierr); 72169c03620SShri Abhyankar ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isinact,idx_inact,PETSC_COPY_VALUES,ISinact);CHKERRQ(ierr); 72269c03620SShri Abhyankar 72369c03620SShri Abhyankar ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr); 72469c03620SShri Abhyankar ierr = PetscFree(idx_act);CHKERRQ(ierr); 72569c03620SShri Abhyankar ierr = PetscFree(idx_inact);CHKERRQ(ierr); 72669c03620SShri Abhyankar PetscFunctionReturn(0); 72769c03620SShri Abhyankar } 72869c03620SShri Abhyankar 72969c03620SShri Abhyankar /* Variational Inequality solver using active set method */ 73069c03620SShri Abhyankar #undef __FUNCT__ 73169c03620SShri Abhyankar #define __FUNCT__ "SNESSolveVI_AS" 73269c03620SShri Abhyankar PetscErrorCode SNESSolveVI_AS(SNES snes) 73369c03620SShri Abhyankar { 73469c03620SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 73569c03620SShri Abhyankar PetscErrorCode ierr; 73669c03620SShri Abhyankar PetscInt maxits,i,lits; 73769c03620SShri Abhyankar PetscBool lssucceed,changedir; 73869c03620SShri Abhyankar MatStructure flg = DIFFERENT_NONZERO_PATTERN; 73969c03620SShri Abhyankar PetscReal gnorm,xnorm=0,ynorm; 74069c03620SShri Abhyankar Vec Y,X,F,G,W; 74169c03620SShri Abhyankar KSPConvergedReason kspreason; 74269c03620SShri Abhyankar IS IS_act,IS_inact; /* _act -> active set _inact -> inactive set */ 74369c03620SShri Abhyankar PetscScalar thresh,J_norm1; 74469c03620SShri Abhyankar 74569c03620SShri Abhyankar PetscFunctionBegin; 74669c03620SShri Abhyankar snes->numFailures = 0; 74769c03620SShri Abhyankar snes->numLinearSolveFailures = 0; 74869c03620SShri Abhyankar snes->reason = SNES_CONVERGED_ITERATING; 74969c03620SShri Abhyankar 75069c03620SShri Abhyankar maxits = snes->max_its; /* maximum number of iterations */ 75169c03620SShri Abhyankar X = snes->vec_sol; /* solution vector */ 75269c03620SShri Abhyankar F = snes->vec_func; /* residual vector */ 75369c03620SShri Abhyankar Y = snes->work[0]; /* work vectors */ 75469c03620SShri Abhyankar G = snes->work[1]; 75569c03620SShri Abhyankar W = snes->work[2]; 75669c03620SShri Abhyankar 75769c03620SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 75869c03620SShri Abhyankar snes->iter = 0; 75969c03620SShri Abhyankar snes->norm = 0.0; 76069c03620SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 76169c03620SShri Abhyankar 76269c03620SShri Abhyankar ierr = SNESVIAdjustInitialGuess(X,vi->xl,vi->xu);CHKERRQ(ierr); 76369c03620SShri Abhyankar ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr); 76469c03620SShri Abhyankar if (snes->domainerror) { 76569c03620SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 76669c03620SShri Abhyankar PetscFunctionReturn(0); 76769c03620SShri Abhyankar } 76869c03620SShri Abhyankar /* Compute Merit function */ 76969c03620SShri Abhyankar ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr); 77069c03620SShri Abhyankar 77169c03620SShri Abhyankar ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 77269c03620SShri Abhyankar ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 77369c03620SShri Abhyankar if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 77469c03620SShri Abhyankar 77569c03620SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 77669c03620SShri Abhyankar snes->norm = vi->phinorm; 77769c03620SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 77869c03620SShri Abhyankar SNESLogConvHistory(snes,vi->phinorm,0); 77969c03620SShri Abhyankar SNESMonitor(snes,0,vi->phinorm); 78069c03620SShri Abhyankar 78169c03620SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 78269c03620SShri Abhyankar snes->ttol = vi->phinorm*snes->rtol; 78369c03620SShri Abhyankar /* test convergence */ 78469c03620SShri Abhyankar ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 78569c03620SShri Abhyankar if (snes->reason) PetscFunctionReturn(0); 78669c03620SShri Abhyankar 78769c03620SShri Abhyankar for (i=0; i<maxits; i++) { 78869c03620SShri Abhyankar 78969c03620SShri Abhyankar /* Call general purpose update function */ 79069c03620SShri Abhyankar if (snes->ops->update) { 79169c03620SShri Abhyankar ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 79269c03620SShri Abhyankar } 79369c03620SShri Abhyankar ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 79469c03620SShri Abhyankar /* Compute the threshold value for creating active and inactive sets */ 79569c03620SShri Abhyankar ierr = MatNorm(snes->jacobian,NORM_1,&J_norm1);CHKERRQ(ierr); 79669c03620SShri Abhyankar thresh = PetscMin(vi->merit,1e-2)/(1+J_norm1); 79769c03620SShri Abhyankar /* Create active and inactive index sets */ 79869c03620SShri Abhyankar ierr = SNESVICreateIndexSets_AS(snes,vi->Db,thresh,&IS_act,&IS_inact);CHKERRQ(ierr); 79969c03620SShri Abhyankar SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Active set semismooth algorithm not implemented yet"); 80069c03620SShri Abhyankar ierr = VecView(vi->Db,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 80169c03620SShri Abhyankar ierr = ISView(IS_act,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 80269c03620SShri Abhyankar ierr = ISView(IS_inact,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 80369c03620SShri Abhyankar 804*730c24dcSShri Abhyankar ierr = ISDestroy(IS_act);CHKERRQ(ierr); 805*730c24dcSShri Abhyankar ierr = ISDestroy(IS_inact);CHKERRQ(ierr); 806*730c24dcSShri Abhyankar 80769c03620SShri Abhyankar ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 80869c03620SShri Abhyankar ierr = SNES_KSPSolve(snes,snes->ksp,vi->phi,Y);CHKERRQ(ierr); 80969c03620SShri Abhyankar ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 81069c03620SShri Abhyankar ierr = SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);CHKERRQ(ierr); 81169c03620SShri Abhyankar ierr = SNESVICheckDescentDirection(snes,vi->dpsi,Y,&changedir);CHKERRQ(ierr); 81269c03620SShri Abhyankar if (kspreason < 0 || changedir) { 81369c03620SShri Abhyankar if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 81469c03620SShri 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); 81569c03620SShri Abhyankar snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 81669c03620SShri Abhyankar break; 81769c03620SShri Abhyankar } 81869c03620SShri Abhyankar ierr = VecCopy(vi->dpsi,Y);CHKERRQ(ierr); 81969c03620SShri Abhyankar } 82069c03620SShri Abhyankar ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 82169c03620SShri Abhyankar snes->linear_its += lits; 82269c03620SShri Abhyankar ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 82369c03620SShri Abhyankar /* 82469c03620SShri Abhyankar if (vi->precheckstep) { 82569c03620SShri Abhyankar PetscBool changed_y = PETSC_FALSE; 82669c03620SShri Abhyankar ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr); 82769c03620SShri Abhyankar } 82869c03620SShri Abhyankar 82969c03620SShri Abhyankar if (PetscLogPrintInfo){ 83069c03620SShri Abhyankar ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 83169c03620SShri Abhyankar } 83269c03620SShri Abhyankar */ 83369c03620SShri Abhyankar /* Compute a (scaled) negative update in the line search routine: 83469c03620SShri Abhyankar Y <- X - lambda*Y 83569c03620SShri Abhyankar and evaluate G = function(Y) (depends on the line search). 83669c03620SShri Abhyankar */ 83769c03620SShri Abhyankar ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 83869c03620SShri Abhyankar ynorm = 1; gnorm = vi->phinorm; 83969c03620SShri Abhyankar ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 84069c03620SShri Abhyankar ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 84169c03620SShri Abhyankar if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 84269c03620SShri Abhyankar if (snes->domainerror) { 84369c03620SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 84469c03620SShri Abhyankar PetscFunctionReturn(0); 84569c03620SShri Abhyankar } 84669c03620SShri Abhyankar if (!lssucceed) { 84769c03620SShri Abhyankar if (++snes->numFailures >= snes->maxFailures) { 84869c03620SShri Abhyankar PetscBool ismin; 84969c03620SShri Abhyankar snes->reason = SNES_DIVERGED_LINE_SEARCH; 85069c03620SShri Abhyankar ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 85169c03620SShri Abhyankar if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 85269c03620SShri Abhyankar break; 85369c03620SShri Abhyankar } 85469c03620SShri Abhyankar } 85569c03620SShri Abhyankar /* Update function and solution vectors */ 85669c03620SShri Abhyankar vi->phinorm = gnorm; 85769c03620SShri Abhyankar vi->merit = 0.5*vi->phinorm*vi->phinorm; 85869c03620SShri Abhyankar ierr = VecCopy(G,vi->phi);CHKERRQ(ierr); 85969c03620SShri Abhyankar ierr = VecCopy(W,X);CHKERRQ(ierr); 86069c03620SShri Abhyankar /* Monitor convergence */ 86169c03620SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 86269c03620SShri Abhyankar snes->iter = i+1; 86369c03620SShri Abhyankar snes->norm = vi->phinorm; 86469c03620SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 86569c03620SShri Abhyankar SNESLogConvHistory(snes,snes->norm,lits); 86669c03620SShri Abhyankar SNESMonitor(snes,snes->iter,snes->norm); 86769c03620SShri Abhyankar /* Test for convergence, xnorm = || X || */ 86869c03620SShri Abhyankar if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 86969c03620SShri Abhyankar ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 87069c03620SShri Abhyankar if (snes->reason) break; 87169c03620SShri Abhyankar } 87269c03620SShri Abhyankar if (i == maxits) { 87369c03620SShri Abhyankar ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 87469c03620SShri Abhyankar if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 87569c03620SShri Abhyankar } 87669c03620SShri Abhyankar PetscFunctionReturn(0); 87769c03620SShri Abhyankar } 87869c03620SShri Abhyankar 8792b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 8802b4ed282SShri Abhyankar /* 8812b4ed282SShri Abhyankar SNESSetUp_VI - Sets up the internal data structures for the later use 8822b4ed282SShri Abhyankar of the SNESVI nonlinear solver. 8832b4ed282SShri Abhyankar 8842b4ed282SShri Abhyankar Input Parameter: 8852b4ed282SShri Abhyankar . snes - the SNES context 8862b4ed282SShri Abhyankar . x - the solution vector 8872b4ed282SShri Abhyankar 8882b4ed282SShri Abhyankar Application Interface Routine: SNESSetUp() 8892b4ed282SShri Abhyankar 8902b4ed282SShri Abhyankar Notes: 8912b4ed282SShri Abhyankar For basic use of the SNES solvers, the user need not explicitly call 8922b4ed282SShri Abhyankar SNESSetUp(), since these actions will automatically occur during 8932b4ed282SShri Abhyankar the call to SNESSolve(). 8942b4ed282SShri Abhyankar */ 8952b4ed282SShri Abhyankar #undef __FUNCT__ 8962b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetUp_VI" 8972b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes) 8982b4ed282SShri Abhyankar { 8992b4ed282SShri Abhyankar PetscErrorCode ierr; 9002b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*) snes->data; 9012b4ed282SShri Abhyankar PetscInt i_start[3],i_end[3]; 9022b4ed282SShri Abhyankar 9032b4ed282SShri Abhyankar PetscFunctionBegin; 9042b4ed282SShri Abhyankar if (!snes->vec_sol_update) { 9052b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr); 9062b4ed282SShri Abhyankar ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr); 9072b4ed282SShri Abhyankar } 9082b4ed282SShri Abhyankar if (!snes->work) { 9092b4ed282SShri Abhyankar snes->nwork = 3; 9102b4ed282SShri Abhyankar ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 9112b4ed282SShri Abhyankar ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr); 9122b4ed282SShri Abhyankar } 9132b4ed282SShri Abhyankar 9142b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->phi); CHKERRQ(ierr); 9152b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->dpsi); CHKERRQ(ierr); 9162b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->Da); CHKERRQ(ierr); 9172b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->Db); CHKERRQ(ierr); 9182b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->z);CHKERRQ(ierr); 9192b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->t); CHKERRQ(ierr); 9202b4ed282SShri Abhyankar 9212b4ed282SShri Abhyankar /* If the lower and upper bound on variables are not set, set it to 9222b4ed282SShri Abhyankar -Inf and Inf */ 9232b4ed282SShri Abhyankar if (!vi->xl && !vi->xu) { 9242b4ed282SShri Abhyankar vi->usersetxbounds = PETSC_FALSE; 9252b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->xl); CHKERRQ(ierr); 9262b4ed282SShri Abhyankar ierr = VecSet(vi->xl,PETSC_VI_NINF);CHKERRQ(ierr); 9272b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->xu); CHKERRQ(ierr); 9282b4ed282SShri Abhyankar ierr = VecSet(vi->xu,PETSC_VI_INF);CHKERRQ(ierr); 9292b4ed282SShri Abhyankar } else { 9302b4ed282SShri Abhyankar /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */ 9312b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr); 9322b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(vi->xl,i_start+1,i_end+1);CHKERRQ(ierr); 9332b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(vi->xu,i_start+2,i_end+2);CHKERRQ(ierr); 9342b4ed282SShri 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])) 9352b4ed282SShri 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."); 9362b4ed282SShri Abhyankar } 9372b4ed282SShri Abhyankar 9382b4ed282SShri Abhyankar vi->computeuserfunction = snes->ops->computefunction; 9392b4ed282SShri Abhyankar vi->computeuserjacobian = snes->ops->computejacobian; 9402b4ed282SShri Abhyankar 9411f79c099SShri Abhyankar snes->ops->computefunction = SNESVIComputeFunction; 9421f79c099SShri Abhyankar snes->ops->computejacobian = SNESVIComputeJacobian; 9432b4ed282SShri Abhyankar PetscFunctionReturn(0); 9442b4ed282SShri Abhyankar } 9452b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 9462b4ed282SShri Abhyankar /* 9472b4ed282SShri Abhyankar SNESDestroy_VI - Destroys the private SNES_VI context that was created 9482b4ed282SShri Abhyankar with SNESCreate_VI(). 9492b4ed282SShri Abhyankar 9502b4ed282SShri Abhyankar Input Parameter: 9512b4ed282SShri Abhyankar . snes - the SNES context 9522b4ed282SShri Abhyankar 9532b4ed282SShri Abhyankar Application Interface Routine: SNESDestroy() 9542b4ed282SShri Abhyankar */ 9552b4ed282SShri Abhyankar #undef __FUNCT__ 9562b4ed282SShri Abhyankar #define __FUNCT__ "SNESDestroy_VI" 9572b4ed282SShri Abhyankar PetscErrorCode SNESDestroy_VI(SNES snes) 9582b4ed282SShri Abhyankar { 9592b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*) snes->data; 9602b4ed282SShri Abhyankar PetscErrorCode ierr; 9612b4ed282SShri Abhyankar 9622b4ed282SShri Abhyankar PetscFunctionBegin; 9632b4ed282SShri Abhyankar if (snes->vec_sol_update) { 9642b4ed282SShri Abhyankar ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr); 9652b4ed282SShri Abhyankar snes->vec_sol_update = PETSC_NULL; 9662b4ed282SShri Abhyankar } 9672b4ed282SShri Abhyankar if (snes->nwork) { 9682b4ed282SShri Abhyankar ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 9692b4ed282SShri Abhyankar snes->nwork = 0; 9702b4ed282SShri Abhyankar snes->work = PETSC_NULL; 9712b4ed282SShri Abhyankar } 9722b4ed282SShri Abhyankar 9732b4ed282SShri Abhyankar /* clear vectors */ 9742b4ed282SShri Abhyankar ierr = VecDestroy(vi->phi); CHKERRQ(ierr); 9752b4ed282SShri Abhyankar ierr = VecDestroy(vi->dpsi); CHKERRQ(ierr); 9762b4ed282SShri Abhyankar ierr = VecDestroy(vi->Da); CHKERRQ(ierr); 9772b4ed282SShri Abhyankar ierr = VecDestroy(vi->Db); CHKERRQ(ierr); 9782b4ed282SShri Abhyankar ierr = VecDestroy(vi->z); CHKERRQ(ierr); 9792b4ed282SShri Abhyankar ierr = VecDestroy(vi->t); CHKERRQ(ierr); 9802b4ed282SShri Abhyankar if (!vi->usersetxbounds) { 9812b4ed282SShri Abhyankar ierr = VecDestroy(vi->xl); CHKERRQ(ierr); 9822b4ed282SShri Abhyankar ierr = VecDestroy(vi->xu); CHKERRQ(ierr); 9832b4ed282SShri Abhyankar } 984c92abb8aSShri Abhyankar if (vi->lsmonitor) { 985c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr); 9862b4ed282SShri Abhyankar } 9872b4ed282SShri Abhyankar ierr = PetscFree(snes->data);CHKERRQ(ierr); 9882b4ed282SShri Abhyankar 9892b4ed282SShri Abhyankar /* clear composed functions */ 9902b4ed282SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr); 991c92abb8aSShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr); 9922b4ed282SShri Abhyankar PetscFunctionReturn(0); 9932b4ed282SShri Abhyankar } 9942b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 9952b4ed282SShri Abhyankar #undef __FUNCT__ 9962b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNo_VI" 9972b4ed282SShri Abhyankar 9982b4ed282SShri Abhyankar /* 9992b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchNo routine in snes/impls/ls/ls.c 10002b4ed282SShri Abhyankar 10012b4ed282SShri Abhyankar */ 1002ace3abfcSBarry 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) 10032b4ed282SShri Abhyankar { 10042b4ed282SShri Abhyankar PetscErrorCode ierr; 10052b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1006ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 10072b4ed282SShri Abhyankar 10082b4ed282SShri Abhyankar PetscFunctionBegin; 10092b4ed282SShri Abhyankar *flag = PETSC_TRUE; 10102b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 10112b4ed282SShri Abhyankar ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 10122b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 10132b4ed282SShri Abhyankar if (vi->postcheckstep) { 10142b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 10152b4ed282SShri Abhyankar } 10162b4ed282SShri Abhyankar if (changed_y) { 10172b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 10182b4ed282SShri Abhyankar } 10192b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 10202b4ed282SShri Abhyankar if (!snes->domainerror) { 10212b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 10222b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 10232b4ed282SShri Abhyankar } 1024c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1025c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 1026c92abb8aSShri Abhyankar } 10272b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 10282b4ed282SShri Abhyankar PetscFunctionReturn(0); 10292b4ed282SShri Abhyankar } 10302b4ed282SShri Abhyankar 10312b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 10322b4ed282SShri Abhyankar #undef __FUNCT__ 10332b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNoNorms_VI" 10342b4ed282SShri Abhyankar 10352b4ed282SShri Abhyankar /* 10362b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c 10372b4ed282SShri Abhyankar */ 1038ace3abfcSBarry 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) 10392b4ed282SShri Abhyankar { 10402b4ed282SShri Abhyankar PetscErrorCode ierr; 10412b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1042ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 10432b4ed282SShri Abhyankar 10442b4ed282SShri Abhyankar PetscFunctionBegin; 10452b4ed282SShri Abhyankar *flag = PETSC_TRUE; 10462b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 10472b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 10482b4ed282SShri Abhyankar if (vi->postcheckstep) { 10492b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 10502b4ed282SShri Abhyankar } 10512b4ed282SShri Abhyankar if (changed_y) { 10522b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 10532b4ed282SShri Abhyankar } 10542b4ed282SShri Abhyankar 10552b4ed282SShri Abhyankar /* don't evaluate function the last time through */ 10562b4ed282SShri Abhyankar if (snes->iter < snes->max_its-1) { 10572b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 10582b4ed282SShri Abhyankar } 10592b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 10602b4ed282SShri Abhyankar PetscFunctionReturn(0); 10612b4ed282SShri Abhyankar } 10622b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 10632b4ed282SShri Abhyankar #undef __FUNCT__ 10642b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchCubic_VI" 10652b4ed282SShri Abhyankar /* 10662b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchCubic in snes/impls/ls/ls.c 10672b4ed282SShri Abhyankar */ 1068ace3abfcSBarry 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) 10692b4ed282SShri Abhyankar { 10702b4ed282SShri Abhyankar /* 10712b4ed282SShri Abhyankar Note that for line search purposes we work with with the related 10722b4ed282SShri Abhyankar minimization problem: 10732b4ed282SShri Abhyankar min z(x): R^n -> R, 10742b4ed282SShri Abhyankar where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 10752b4ed282SShri Abhyankar This function z(x) is same as the merit function for the complementarity problem. 10762b4ed282SShri Abhyankar */ 10772b4ed282SShri Abhyankar 10782b4ed282SShri Abhyankar PetscReal initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 10792b4ed282SShri Abhyankar PetscReal minlambda,lambda,lambdatemp; 10802b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 10812b4ed282SShri Abhyankar PetscScalar cinitslope; 10822b4ed282SShri Abhyankar #endif 10832b4ed282SShri Abhyankar PetscErrorCode ierr; 10842b4ed282SShri Abhyankar PetscInt count; 10852b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1086ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 10872b4ed282SShri Abhyankar MPI_Comm comm; 10882b4ed282SShri Abhyankar 10892b4ed282SShri Abhyankar PetscFunctionBegin; 10902b4ed282SShri Abhyankar ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 10912b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 10922b4ed282SShri Abhyankar *flag = PETSC_TRUE; 10932b4ed282SShri Abhyankar 10942b4ed282SShri Abhyankar ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 10952b4ed282SShri Abhyankar if (*ynorm == 0.0) { 1096c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1097c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 10982b4ed282SShri Abhyankar } 10992b4ed282SShri Abhyankar *gnorm = fnorm; 11002b4ed282SShri Abhyankar ierr = VecCopy(x,w);CHKERRQ(ierr); 11012b4ed282SShri Abhyankar ierr = VecCopy(f,g);CHKERRQ(ierr); 11022b4ed282SShri Abhyankar *flag = PETSC_FALSE; 11032b4ed282SShri Abhyankar goto theend1; 11042b4ed282SShri Abhyankar } 11052b4ed282SShri Abhyankar if (*ynorm > vi->maxstep) { /* Step too big, so scale back */ 1106c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1107c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Scaling step by %G old ynorm %G\n",vi->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr); 11082b4ed282SShri Abhyankar } 11092b4ed282SShri Abhyankar ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr); 11102b4ed282SShri Abhyankar *ynorm = vi->maxstep; 11112b4ed282SShri Abhyankar } 11122b4ed282SShri Abhyankar ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 11132b4ed282SShri Abhyankar minlambda = vi->minlambda/rellength; 11142b4ed282SShri Abhyankar /* ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); */ 11152b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 11162b4ed282SShri Abhyankar ierr = VecDot(vi->dpsi,y,&cinitslope);CHKERRQ(ierr); 11172b4ed282SShri Abhyankar initslope = PetscRealPart(cinitslope); 11182b4ed282SShri Abhyankar #else 11192b4ed282SShri Abhyankar ierr = VecDot(vi->dpsi,y,&initslope);CHKERRQ(ierr); 11202b4ed282SShri Abhyankar #endif 11212b4ed282SShri Abhyankar if (initslope > 0.0) initslope = -initslope; 11222b4ed282SShri Abhyankar if (initslope == 0.0) initslope = -1.0; 11232b4ed282SShri Abhyankar 11242b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 11252b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 11262b4ed282SShri Abhyankar ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 11272b4ed282SShri Abhyankar *flag = PETSC_FALSE; 11282b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 11292b4ed282SShri Abhyankar goto theend1; 11302b4ed282SShri Abhyankar } 11312b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 11322b4ed282SShri Abhyankar if (snes->domainerror) { 11332b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 11342b4ed282SShri Abhyankar PetscFunctionReturn(0); 11352b4ed282SShri Abhyankar } 11362b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 11372b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 11382b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 11392b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */ 1140c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1141c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 11422b4ed282SShri Abhyankar } 11432b4ed282SShri Abhyankar goto theend1; 11442b4ed282SShri Abhyankar } 11452b4ed282SShri Abhyankar 11462b4ed282SShri Abhyankar /* Fit points with quadratic */ 11472b4ed282SShri Abhyankar lambda = 1.0; 11482b4ed282SShri Abhyankar lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 11492b4ed282SShri Abhyankar lambdaprev = lambda; 11502b4ed282SShri Abhyankar gnormprev = *gnorm; 11512b4ed282SShri Abhyankar if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 11522b4ed282SShri Abhyankar if (lambdatemp <= .1*lambda) lambda = .1*lambda; 11532b4ed282SShri Abhyankar else lambda = lambdatemp; 11542b4ed282SShri Abhyankar 11552b4ed282SShri Abhyankar ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 11562b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 11572b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 11582b4ed282SShri Abhyankar *flag = PETSC_FALSE; 11592b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 11602b4ed282SShri Abhyankar goto theend1; 11612b4ed282SShri Abhyankar } 11622b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 11632b4ed282SShri Abhyankar if (snes->domainerror) { 11642b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 11652b4ed282SShri Abhyankar PetscFunctionReturn(0); 11662b4ed282SShri Abhyankar } 11672b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 11682b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1169c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1170c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr); 11712b4ed282SShri Abhyankar } 11722b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */ 1173c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1174c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 11752b4ed282SShri Abhyankar } 11762b4ed282SShri Abhyankar goto theend1; 11772b4ed282SShri Abhyankar } 11782b4ed282SShri Abhyankar 11792b4ed282SShri Abhyankar /* Fit points with cubic */ 11802b4ed282SShri Abhyankar count = 1; 11812b4ed282SShri Abhyankar while (PETSC_TRUE) { 11822b4ed282SShri Abhyankar if (lambda <= minlambda) { 1183c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1184c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 1185c92abb8aSShri 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); 11862b4ed282SShri Abhyankar } 11872b4ed282SShri Abhyankar *flag = PETSC_FALSE; 11882b4ed282SShri Abhyankar break; 11892b4ed282SShri Abhyankar } 11902b4ed282SShri Abhyankar t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 11912b4ed282SShri Abhyankar t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 11922b4ed282SShri Abhyankar a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 11932b4ed282SShri Abhyankar b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 11942b4ed282SShri Abhyankar d = b*b - 3*a*initslope; 11952b4ed282SShri Abhyankar if (d < 0.0) d = 0.0; 11962b4ed282SShri Abhyankar if (a == 0.0) { 11972b4ed282SShri Abhyankar lambdatemp = -initslope/(2.0*b); 11982b4ed282SShri Abhyankar } else { 11992b4ed282SShri Abhyankar lambdatemp = (-b + sqrt(d))/(3.0*a); 12002b4ed282SShri Abhyankar } 12012b4ed282SShri Abhyankar lambdaprev = lambda; 12022b4ed282SShri Abhyankar gnormprev = *gnorm; 12032b4ed282SShri Abhyankar if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 12042b4ed282SShri Abhyankar if (lambdatemp <= .1*lambda) lambda = .1*lambda; 12052b4ed282SShri Abhyankar else lambda = lambdatemp; 12062b4ed282SShri Abhyankar ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 12072b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 12082b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 12092b4ed282SShri 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); 12102b4ed282SShri Abhyankar *flag = PETSC_FALSE; 12112b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 12122b4ed282SShri Abhyankar break; 12132b4ed282SShri Abhyankar } 12142b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 12152b4ed282SShri Abhyankar if (snes->domainerror) { 12162b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 12172b4ed282SShri Abhyankar PetscFunctionReturn(0); 12182b4ed282SShri Abhyankar } 12192b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 12202b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 12212b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* is reduction enough? */ 1222c92abb8aSShri Abhyankar if (vi->lsmonitor) { 12232b4ed282SShri Abhyankar ierr = PetscPrintf(comm," Line search: Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 12242b4ed282SShri Abhyankar } 12252b4ed282SShri Abhyankar break; 12262b4ed282SShri Abhyankar } else { 1227c92abb8aSShri Abhyankar if (vi->lsmonitor) { 12282b4ed282SShri Abhyankar ierr = PetscPrintf(comm," Line search: Cubic step no good, shrinking lambda, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 12292b4ed282SShri Abhyankar } 12302b4ed282SShri Abhyankar } 12312b4ed282SShri Abhyankar count++; 12322b4ed282SShri Abhyankar } 12332b4ed282SShri Abhyankar theend1: 12342b4ed282SShri Abhyankar /* Optional user-defined check for line search step validity */ 12352b4ed282SShri Abhyankar if (vi->postcheckstep && *flag) { 12362b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 12372b4ed282SShri Abhyankar if (changed_y) { 12382b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 12392b4ed282SShri Abhyankar } 12402b4ed282SShri Abhyankar if (changed_y || changed_w) { /* recompute the function if the step has changed */ 12412b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 12422b4ed282SShri Abhyankar if (snes->domainerror) { 12432b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 12442b4ed282SShri Abhyankar PetscFunctionReturn(0); 12452b4ed282SShri Abhyankar } 12462b4ed282SShri Abhyankar ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 12472b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 12482b4ed282SShri Abhyankar ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 12492b4ed282SShri Abhyankar ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 12502b4ed282SShri Abhyankar ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 12512b4ed282SShri Abhyankar } 12522b4ed282SShri Abhyankar } 12532b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 12542b4ed282SShri Abhyankar PetscFunctionReturn(0); 12552b4ed282SShri Abhyankar } 12562b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 12572b4ed282SShri Abhyankar #undef __FUNCT__ 12582b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchQuadratic_VI" 12592b4ed282SShri Abhyankar /* 12602b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchQuadratic in snes/impls/ls/ls.c 12612b4ed282SShri Abhyankar */ 1262ace3abfcSBarry 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) 12632b4ed282SShri Abhyankar { 12642b4ed282SShri Abhyankar /* 12652b4ed282SShri Abhyankar Note that for line search purposes we work with with the related 12662b4ed282SShri Abhyankar minimization problem: 12672b4ed282SShri Abhyankar min z(x): R^n -> R, 12682b4ed282SShri Abhyankar where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 12692b4ed282SShri Abhyankar z(x) is the same as the merit function for the complementarity problem 12702b4ed282SShri Abhyankar */ 12712b4ed282SShri Abhyankar PetscReal initslope,minlambda,lambda,lambdatemp,rellength; 12722b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 12732b4ed282SShri Abhyankar PetscScalar cinitslope; 12742b4ed282SShri Abhyankar #endif 12752b4ed282SShri Abhyankar PetscErrorCode ierr; 12762b4ed282SShri Abhyankar PetscInt count; 12772b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1278ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 12792b4ed282SShri Abhyankar 12802b4ed282SShri Abhyankar PetscFunctionBegin; 12812b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 12822b4ed282SShri Abhyankar *flag = PETSC_TRUE; 12832b4ed282SShri Abhyankar 12842b4ed282SShri Abhyankar ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 12852b4ed282SShri Abhyankar if (*ynorm == 0.0) { 1286c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1287c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr); 12882b4ed282SShri Abhyankar } 12892b4ed282SShri Abhyankar *gnorm = fnorm; 12902b4ed282SShri Abhyankar ierr = VecCopy(x,w);CHKERRQ(ierr); 12912b4ed282SShri Abhyankar ierr = VecCopy(f,g);CHKERRQ(ierr); 12922b4ed282SShri Abhyankar *flag = PETSC_FALSE; 12932b4ed282SShri Abhyankar goto theend2; 12942b4ed282SShri Abhyankar } 12952b4ed282SShri Abhyankar if (*ynorm > vi->maxstep) { /* Step too big, so scale back */ 12962b4ed282SShri Abhyankar ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr); 12972b4ed282SShri Abhyankar *ynorm = vi->maxstep; 12982b4ed282SShri Abhyankar } 12992b4ed282SShri Abhyankar ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 13002b4ed282SShri Abhyankar minlambda = vi->minlambda/rellength; 13012b4ed282SShri Abhyankar /* ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); */ 13022b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 13032b4ed282SShri Abhyankar ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 13042b4ed282SShri Abhyankar initslope = PetscRealPart(cinitslope); 13052b4ed282SShri Abhyankar #else 13062b4ed282SShri Abhyankar ierr = VecDot(vi->dpsi,y,&initslope);CHKERRQ(ierr); 13072b4ed282SShri Abhyankar #endif 13082b4ed282SShri Abhyankar if (initslope > 0.0) initslope = -initslope; 13092b4ed282SShri Abhyankar if (initslope == 0.0) initslope = -1.0; 13102b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr); 13112b4ed282SShri Abhyankar 13122b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 13132b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 13142b4ed282SShri Abhyankar ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 13152b4ed282SShri Abhyankar *flag = PETSC_FALSE; 13162b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 13172b4ed282SShri Abhyankar goto theend2; 13182b4ed282SShri Abhyankar } 13192b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 13202b4ed282SShri Abhyankar if (snes->domainerror) { 13212b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 13222b4ed282SShri Abhyankar PetscFunctionReturn(0); 13232b4ed282SShri Abhyankar } 13242b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 13252b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 13262b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */ 1327c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1328c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 13292b4ed282SShri Abhyankar } 13302b4ed282SShri Abhyankar goto theend2; 13312b4ed282SShri Abhyankar } 13322b4ed282SShri Abhyankar 13332b4ed282SShri Abhyankar /* Fit points with quadratic */ 13342b4ed282SShri Abhyankar lambda = 1.0; 13352b4ed282SShri Abhyankar count = 1; 13362b4ed282SShri Abhyankar while (PETSC_TRUE) { 13372b4ed282SShri Abhyankar if (lambda <= minlambda) { /* bad luck; use full step */ 1338c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1339c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr); 1340c92abb8aSShri 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); 13412b4ed282SShri Abhyankar } 13422b4ed282SShri Abhyankar ierr = VecCopy(x,w);CHKERRQ(ierr); 13432b4ed282SShri Abhyankar *flag = PETSC_FALSE; 13442b4ed282SShri Abhyankar break; 13452b4ed282SShri Abhyankar } 13462b4ed282SShri Abhyankar lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 13472b4ed282SShri Abhyankar if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 13482b4ed282SShri Abhyankar if (lambdatemp <= .1*lambda) lambda = .1*lambda; 13492b4ed282SShri Abhyankar else lambda = lambdatemp; 13502b4ed282SShri Abhyankar 13512b4ed282SShri Abhyankar ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 13522b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 13532b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 13542b4ed282SShri 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); 13552b4ed282SShri Abhyankar *flag = PETSC_FALSE; 13562b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 13572b4ed282SShri Abhyankar break; 13582b4ed282SShri Abhyankar } 13592b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 13602b4ed282SShri Abhyankar if (snes->domainerror) { 13612b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 13622b4ed282SShri Abhyankar PetscFunctionReturn(0); 13632b4ed282SShri Abhyankar } 13642b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 13652b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 13662b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */ 1367c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1368c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line Search: Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr); 13692b4ed282SShri Abhyankar } 13702b4ed282SShri Abhyankar break; 13712b4ed282SShri Abhyankar } 13722b4ed282SShri Abhyankar count++; 13732b4ed282SShri Abhyankar } 13742b4ed282SShri Abhyankar theend2: 13752b4ed282SShri Abhyankar /* Optional user-defined check for line search step validity */ 13762b4ed282SShri Abhyankar if (vi->postcheckstep) { 13772b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 13782b4ed282SShri Abhyankar if (changed_y) { 13792b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 13802b4ed282SShri Abhyankar } 13812b4ed282SShri Abhyankar if (changed_y || changed_w) { /* recompute the function if the step has changed */ 13822b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g); 13832b4ed282SShri Abhyankar if (snes->domainerror) { 13842b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 13852b4ed282SShri Abhyankar PetscFunctionReturn(0); 13862b4ed282SShri Abhyankar } 13872b4ed282SShri Abhyankar ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 13882b4ed282SShri Abhyankar ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 13892b4ed282SShri Abhyankar ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 13902b4ed282SShri Abhyankar ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 13912b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 13922b4ed282SShri Abhyankar } 13932b4ed282SShri Abhyankar } 13942b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 13952b4ed282SShri Abhyankar PetscFunctionReturn(0); 13962b4ed282SShri Abhyankar } 13972b4ed282SShri Abhyankar 1398ace3abfcSBarry 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*/ 13992b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 14002b4ed282SShri Abhyankar EXTERN_C_BEGIN 14012b4ed282SShri Abhyankar #undef __FUNCT__ 14022b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchSet_VI" 14032b4ed282SShri Abhyankar PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_VI(SNES snes,FCN2 func,void *lsctx) 14042b4ed282SShri Abhyankar { 14052b4ed282SShri Abhyankar PetscFunctionBegin; 14062b4ed282SShri Abhyankar ((SNES_VI *)(snes->data))->LineSearch = func; 14072b4ed282SShri Abhyankar ((SNES_VI *)(snes->data))->lsP = lsctx; 14082b4ed282SShri Abhyankar PetscFunctionReturn(0); 14092b4ed282SShri Abhyankar } 14102b4ed282SShri Abhyankar EXTERN_C_END 14112b4ed282SShri Abhyankar 14122b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 14132b4ed282SShri Abhyankar EXTERN_C_BEGIN 14142b4ed282SShri Abhyankar #undef __FUNCT__ 14152b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchSetMonitor_VI" 1416ace3abfcSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetMonitor_VI(SNES snes,PetscBool flg) 14172b4ed282SShri Abhyankar { 14182b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 14192b4ed282SShri Abhyankar PetscErrorCode ierr; 14202b4ed282SShri Abhyankar 14212b4ed282SShri Abhyankar PetscFunctionBegin; 1422c92abb8aSShri Abhyankar if (flg && !vi->lsmonitor) { 1423c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",((PetscObject)snes)->tablevel,&vi->lsmonitor);CHKERRQ(ierr); 1424c92abb8aSShri Abhyankar } else if (!flg && vi->lsmonitor) { 1425c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr); 14262b4ed282SShri Abhyankar } 14272b4ed282SShri Abhyankar PetscFunctionReturn(0); 14282b4ed282SShri Abhyankar } 14292b4ed282SShri Abhyankar EXTERN_C_END 14302b4ed282SShri Abhyankar 14312b4ed282SShri Abhyankar /* 14322b4ed282SShri Abhyankar SNESView_VI - Prints info from the SNESVI data structure. 14332b4ed282SShri Abhyankar 14342b4ed282SShri Abhyankar Input Parameters: 14352b4ed282SShri Abhyankar . SNES - the SNES context 14362b4ed282SShri Abhyankar . viewer - visualization context 14372b4ed282SShri Abhyankar 14382b4ed282SShri Abhyankar Application Interface Routine: SNESView() 14392b4ed282SShri Abhyankar */ 14402b4ed282SShri Abhyankar #undef __FUNCT__ 14412b4ed282SShri Abhyankar #define __FUNCT__ "SNESView_VI" 14422b4ed282SShri Abhyankar static PetscErrorCode SNESView_VI(SNES snes,PetscViewer viewer) 14432b4ed282SShri Abhyankar { 14442b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI *)snes->data; 14452b4ed282SShri Abhyankar const char *cstr; 14462b4ed282SShri Abhyankar PetscErrorCode ierr; 1447ace3abfcSBarry Smith PetscBool iascii; 14482b4ed282SShri Abhyankar 14492b4ed282SShri Abhyankar PetscFunctionBegin; 14502b4ed282SShri Abhyankar ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 14512b4ed282SShri Abhyankar if (iascii) { 14522b4ed282SShri Abhyankar if (vi->LineSearch == SNESLineSearchNo_VI) cstr = "SNESLineSearchNo"; 14532b4ed282SShri Abhyankar else if (vi->LineSearch == SNESLineSearchQuadratic_VI) cstr = "SNESLineSearchQuadratic"; 14542b4ed282SShri Abhyankar else if (vi->LineSearch == SNESLineSearchCubic_VI) cstr = "SNESLineSearchCubic"; 14552b4ed282SShri Abhyankar else cstr = "unknown"; 14562b4ed282SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 14572b4ed282SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," alpha=%G, maxstep=%G, minlambda=%G\n",vi->alpha,vi->maxstep,vi->minlambda);CHKERRQ(ierr); 14582b4ed282SShri Abhyankar } else { 14592b4ed282SShri Abhyankar SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ VI",((PetscObject)viewer)->type_name); 14602b4ed282SShri Abhyankar } 14612b4ed282SShri Abhyankar PetscFunctionReturn(0); 14622b4ed282SShri Abhyankar } 14632b4ed282SShri Abhyankar 14642b4ed282SShri Abhyankar /* 14652b4ed282SShri Abhyankar SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. 14662b4ed282SShri Abhyankar 14672b4ed282SShri Abhyankar Input Parameters: 14682b4ed282SShri Abhyankar . snes - the SNES context. 14692b4ed282SShri Abhyankar . xl - lower bound. 14702b4ed282SShri Abhyankar . xu - upper bound. 14712b4ed282SShri Abhyankar 14722b4ed282SShri Abhyankar Notes: 14732b4ed282SShri Abhyankar If this routine is not called then the lower and upper bounds are set to 14742b4ed282SShri Abhyankar -Infinity and Infinity respectively during SNESSetUp. 14752b4ed282SShri Abhyankar */ 14762b4ed282SShri Abhyankar 14772b4ed282SShri Abhyankar #undef __FUNCT__ 14782b4ed282SShri Abhyankar #define __FUNCT__ "SNESVISetVariableBounds" 147969c03620SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu) 14802b4ed282SShri Abhyankar { 14812b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 14822b4ed282SShri Abhyankar 14832b4ed282SShri Abhyankar PetscFunctionBegin; 14842b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 14852b4ed282SShri Abhyankar PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 14862b4ed282SShri Abhyankar PetscValidHeaderSpecific(xu,VEC_CLASSID,3); 14872b4ed282SShri Abhyankar 14882b4ed282SShri Abhyankar /* Check if SNESSetFunction is called */ 14892b4ed282SShri Abhyankar if(!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 14902b4ed282SShri Abhyankar 14912b4ed282SShri Abhyankar /* Check if the vector sizes are compatible for lower and upper bounds */ 14922b4ed282SShri 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); 14932b4ed282SShri 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); 14942b4ed282SShri Abhyankar vi->usersetxbounds = PETSC_TRUE; 14952b4ed282SShri Abhyankar vi->xl = xl; 14962b4ed282SShri Abhyankar vi->xu = xu; 14972b4ed282SShri Abhyankar 14982b4ed282SShri Abhyankar PetscFunctionReturn(0); 14992b4ed282SShri Abhyankar } 15002b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 15012b4ed282SShri Abhyankar /* 15022b4ed282SShri Abhyankar SNESSetFromOptions_VI - Sets various parameters for the SNESVI method. 15032b4ed282SShri Abhyankar 15042b4ed282SShri Abhyankar Input Parameter: 15052b4ed282SShri Abhyankar . snes - the SNES context 15062b4ed282SShri Abhyankar 15072b4ed282SShri Abhyankar Application Interface Routine: SNESSetFromOptions() 15082b4ed282SShri Abhyankar */ 15092b4ed282SShri Abhyankar #undef __FUNCT__ 15102b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetFromOptions_VI" 15112b4ed282SShri Abhyankar static PetscErrorCode SNESSetFromOptions_VI(SNES snes) 15122b4ed282SShri Abhyankar { 15132b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI *)snes->data; 15142b4ed282SShri Abhyankar const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 151569c03620SShri Abhyankar const char *vies[] = {"ss","as"}; 15162b4ed282SShri Abhyankar PetscErrorCode ierr; 15172b4ed282SShri Abhyankar PetscInt indx; 151869c03620SShri Abhyankar PetscBool flg,set,flg2; 15192b4ed282SShri Abhyankar 15202b4ed282SShri Abhyankar PetscFunctionBegin; 15212b4ed282SShri Abhyankar ierr = PetscOptionsHead("SNES semismooth method options");CHKERRQ(ierr); 15222b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_alpha","Function norm must decrease by","None",vi->alpha,&vi->alpha,0);CHKERRQ(ierr); 15232b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_maxstep","Step must be less than","None",vi->maxstep,&vi->maxstep,0);CHKERRQ(ierr); 15242b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_minlambda","Minimum lambda allowed","None",vi->minlambda,&vi->minlambda,0);CHKERRQ(ierr); 15252b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_delta","descent test fraction","None",vi->delta,&vi->delta,0);CHKERRQ(ierr); 15262b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_rho","descent test power","None",vi->rho,&vi->rho,0);CHKERRQ(ierr); 15272b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_const_tol","constraint tolerance","None",vi->const_tol,&vi->const_tol,0);CHKERRQ(ierr); 1528acfcf0e5SJed Brown ierr = PetscOptionsBool("-snes_vi_lsmonitor","Print progress of line searches","SNESLineSearchSetMonitor",vi->lsmonitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 15292b4ed282SShri Abhyankar if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);} 153069c03620SShri Abhyankar ierr = PetscOptionsEList("-snes_vi_type","Semismooth algorithm used","",vies,2,"ss",&indx,&flg2);CHKERRQ(ierr); 153169c03620SShri Abhyankar if (flg2) { 153269c03620SShri Abhyankar switch (indx) { 153369c03620SShri Abhyankar case 0: 153469c03620SShri Abhyankar snes->ops->solve = SNESSolveVI_SS; 153569c03620SShri Abhyankar break; 153669c03620SShri Abhyankar case 1: 153769c03620SShri Abhyankar snes->ops->solve = SNESSolveVI_AS; 153869c03620SShri Abhyankar break; 153969c03620SShri Abhyankar } 154069c03620SShri Abhyankar } 1541c92abb8aSShri Abhyankar ierr = PetscOptionsEList("-snes_vi_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 15422b4ed282SShri Abhyankar if (flg) { 15432b4ed282SShri Abhyankar switch (indx) { 15442b4ed282SShri Abhyankar case 0: 15452b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr); 15462b4ed282SShri Abhyankar break; 15472b4ed282SShri Abhyankar case 1: 15482b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr); 15492b4ed282SShri Abhyankar break; 15502b4ed282SShri Abhyankar case 2: 15512b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr); 15522b4ed282SShri Abhyankar break; 15532b4ed282SShri Abhyankar case 3: 15542b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr); 15552b4ed282SShri Abhyankar break; 15562b4ed282SShri Abhyankar } 15572b4ed282SShri Abhyankar } 15582b4ed282SShri Abhyankar ierr = PetscOptionsTail();CHKERRQ(ierr); 15592b4ed282SShri Abhyankar PetscFunctionReturn(0); 15602b4ed282SShri Abhyankar } 15612b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 15622b4ed282SShri Abhyankar /*MC 15632b4ed282SShri Abhyankar SNESVI - Semismooth newton method based nonlinear solver that uses a line search 15642b4ed282SShri Abhyankar 15652b4ed282SShri Abhyankar Options Database: 15662b4ed282SShri Abhyankar + -snes_vi [cubic,quadratic,basic,basicnonorms] - Selects line search 15672b4ed282SShri Abhyankar . -snes_vi_alpha <alpha> - Sets alpha 15682b4ed282SShri 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) 15692b4ed282SShri Abhyankar . -snes_vi_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 15702b4ed282SShri Abhyankar -snes_vi_delta <delta> - Sets the fraction used in the descent test. 15712b4ed282SShri Abhyankar -snes_vi_rho <rho> - Sets the power used in the descent test. 15722b4ed282SShri Abhyankar For a descent direction to be accepted it has to satisfy the condition dpsi^T*d <= -delta*||d||^rho 15732b4ed282SShri Abhyankar - -snes_vi_monitor - print information about progress of line searches 15742b4ed282SShri Abhyankar 15752b4ed282SShri Abhyankar 15762b4ed282SShri Abhyankar Level: beginner 15772b4ed282SShri Abhyankar 15782b4ed282SShri Abhyankar .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 15792b4ed282SShri Abhyankar SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 15802b4ed282SShri Abhyankar SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams() 15812b4ed282SShri Abhyankar 15822b4ed282SShri Abhyankar M*/ 15832b4ed282SShri Abhyankar EXTERN_C_BEGIN 15842b4ed282SShri Abhyankar #undef __FUNCT__ 15852b4ed282SShri Abhyankar #define __FUNCT__ "SNESCreate_VI" 15862b4ed282SShri Abhyankar PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_VI(SNES snes) 15872b4ed282SShri Abhyankar { 15882b4ed282SShri Abhyankar PetscErrorCode ierr; 15892b4ed282SShri Abhyankar SNES_VI *vi; 15902b4ed282SShri Abhyankar 15912b4ed282SShri Abhyankar PetscFunctionBegin; 15922b4ed282SShri Abhyankar snes->ops->setup = SNESSetUp_VI; 159369c03620SShri Abhyankar snes->ops->solve = SNESSolveVI_SS; 15942b4ed282SShri Abhyankar snes->ops->destroy = SNESDestroy_VI; 15952b4ed282SShri Abhyankar snes->ops->setfromoptions = SNESSetFromOptions_VI; 15962b4ed282SShri Abhyankar snes->ops->view = SNESView_VI; 15972b4ed282SShri Abhyankar snes->ops->converged = SNESDefaultConverged_VI; 15982b4ed282SShri Abhyankar 15992b4ed282SShri Abhyankar ierr = PetscNewLog(snes,SNES_VI,&vi);CHKERRQ(ierr); 16002b4ed282SShri Abhyankar snes->data = (void*)vi; 16012b4ed282SShri Abhyankar vi->alpha = 1.e-4; 16022b4ed282SShri Abhyankar vi->maxstep = 1.e8; 16032b4ed282SShri Abhyankar vi->minlambda = 1.e-12; 16042b4ed282SShri Abhyankar vi->LineSearch = SNESLineSearchCubic_VI; 16052b4ed282SShri Abhyankar vi->lsP = PETSC_NULL; 16062b4ed282SShri Abhyankar vi->postcheckstep = PETSC_NULL; 16072b4ed282SShri Abhyankar vi->postcheck = PETSC_NULL; 16082b4ed282SShri Abhyankar vi->precheckstep = PETSC_NULL; 16092b4ed282SShri Abhyankar vi->precheck = PETSC_NULL; 16102b4ed282SShri Abhyankar vi->rho = 2.1; 16112b4ed282SShri Abhyankar vi->delta = 1e-10; 16122b4ed282SShri Abhyankar vi->const_tol = 2.2204460492503131e-16; 16132b4ed282SShri Abhyankar vi->computessfunction = ComputeFischerFunction; 16142b4ed282SShri Abhyankar 16152b4ed282SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_VI",SNESLineSearchSetMonitor_VI);CHKERRQ(ierr); 16162b4ed282SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_VI",SNESLineSearchSet_VI);CHKERRQ(ierr); 16172b4ed282SShri Abhyankar 16182b4ed282SShri Abhyankar PetscFunctionReturn(0); 16192b4ed282SShri Abhyankar } 16202b4ed282SShri Abhyankar EXTERN_C_END 1621