12b4ed282SShri Abhyankar #define PETSCSNES_DLL 22b4ed282SShri Abhyankar 32b4ed282SShri Abhyankar #include "../src/snes/impls/vi/viimpl.h" 4408da460SShri Abhyankar #include "../include/private/kspimpl.h" 52b4ed282SShri Abhyankar 62b4ed282SShri Abhyankar /* 72b4ed282SShri Abhyankar Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function, 82b4ed282SShri 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 92b4ed282SShri Abhyankar 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 102b4ed282SShri Abhyankar for this trick. One assumes that the probability that W is in the null space of J is very, very small. 112b4ed282SShri Abhyankar */ 122b4ed282SShri Abhyankar #undef __FUNCT__ 132b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckLocalMin_Private" 14ace3abfcSBarry Smith PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin) 152b4ed282SShri Abhyankar { 162b4ed282SShri Abhyankar PetscReal a1; 172b4ed282SShri Abhyankar PetscErrorCode ierr; 18ace3abfcSBarry Smith PetscBool hastranspose; 192b4ed282SShri Abhyankar 202b4ed282SShri Abhyankar PetscFunctionBegin; 212b4ed282SShri Abhyankar *ismin = PETSC_FALSE; 222b4ed282SShri Abhyankar ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 232b4ed282SShri Abhyankar if (hastranspose) { 242b4ed282SShri Abhyankar /* Compute || J^T F|| */ 252b4ed282SShri Abhyankar ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 262b4ed282SShri Abhyankar ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 272b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"|| J^T F|| %G near zero implies found a local minimum\n",a1/fnorm);CHKERRQ(ierr); 282b4ed282SShri Abhyankar if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 292b4ed282SShri Abhyankar } else { 302b4ed282SShri Abhyankar Vec work; 312b4ed282SShri Abhyankar PetscScalar result; 322b4ed282SShri Abhyankar PetscReal wnorm; 332b4ed282SShri Abhyankar 342b4ed282SShri Abhyankar ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr); 352b4ed282SShri Abhyankar ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr); 362b4ed282SShri Abhyankar ierr = VecDuplicate(W,&work);CHKERRQ(ierr); 372b4ed282SShri Abhyankar ierr = MatMult(A,W,work);CHKERRQ(ierr); 382b4ed282SShri Abhyankar ierr = VecDot(F,work,&result);CHKERRQ(ierr); 392b4ed282SShri Abhyankar ierr = VecDestroy(work);CHKERRQ(ierr); 402b4ed282SShri Abhyankar a1 = PetscAbsScalar(result)/(fnorm*wnorm); 412b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %G near zero implies found a local minimum\n",a1);CHKERRQ(ierr); 422b4ed282SShri Abhyankar if (a1 < 1.e-4) *ismin = PETSC_TRUE; 432b4ed282SShri Abhyankar } 442b4ed282SShri Abhyankar PetscFunctionReturn(0); 452b4ed282SShri Abhyankar } 462b4ed282SShri Abhyankar 472b4ed282SShri Abhyankar /* 482b4ed282SShri Abhyankar Checks if J^T(F - J*X) = 0 492b4ed282SShri Abhyankar */ 502b4ed282SShri Abhyankar #undef __FUNCT__ 512b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckResidual_Private" 522b4ed282SShri Abhyankar PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2) 532b4ed282SShri Abhyankar { 542b4ed282SShri Abhyankar PetscReal a1,a2; 552b4ed282SShri Abhyankar PetscErrorCode ierr; 56ace3abfcSBarry Smith PetscBool hastranspose; 572b4ed282SShri Abhyankar 582b4ed282SShri Abhyankar PetscFunctionBegin; 592b4ed282SShri Abhyankar ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 602b4ed282SShri Abhyankar if (hastranspose) { 612b4ed282SShri Abhyankar ierr = MatMult(A,X,W1);CHKERRQ(ierr); 622b4ed282SShri Abhyankar ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr); 632b4ed282SShri Abhyankar 642b4ed282SShri Abhyankar /* Compute || J^T W|| */ 652b4ed282SShri Abhyankar ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 662b4ed282SShri Abhyankar ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 672b4ed282SShri Abhyankar ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 682b4ed282SShri Abhyankar if (a1 != 0.0) { 692b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %G near zero implies inconsistent rhs\n",a2/a1);CHKERRQ(ierr); 702b4ed282SShri Abhyankar } 712b4ed282SShri Abhyankar } 722b4ed282SShri Abhyankar PetscFunctionReturn(0); 732b4ed282SShri Abhyankar } 742b4ed282SShri Abhyankar 752b4ed282SShri Abhyankar /* 762b4ed282SShri Abhyankar SNESDefaultConverged_VI - Checks the convergence of the semismooth newton algorithm. 772b4ed282SShri Abhyankar 782b4ed282SShri Abhyankar Notes: 792b4ed282SShri Abhyankar The convergence criterion currently implemented is 802b4ed282SShri Abhyankar merit < abstol 812b4ed282SShri Abhyankar merit < rtol*merit_initial 822b4ed282SShri Abhyankar */ 832b4ed282SShri Abhyankar #undef __FUNCT__ 842b4ed282SShri Abhyankar #define __FUNCT__ "SNESDefaultConverged_VI" 852b4ed282SShri Abhyankar PetscErrorCode SNESDefaultConverged_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal merit,SNESConvergedReason *reason,void *dummy) 862b4ed282SShri Abhyankar { 872b4ed282SShri Abhyankar PetscErrorCode ierr; 882b4ed282SShri Abhyankar 892b4ed282SShri Abhyankar PetscFunctionBegin; 902b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 912b4ed282SShri Abhyankar PetscValidPointer(reason,6); 922b4ed282SShri Abhyankar 932b4ed282SShri Abhyankar *reason = SNES_CONVERGED_ITERATING; 942b4ed282SShri Abhyankar 952b4ed282SShri Abhyankar if (!it) { 962b4ed282SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 972b4ed282SShri Abhyankar snes->ttol = merit*snes->rtol; 982b4ed282SShri Abhyankar } 992b4ed282SShri Abhyankar if (merit != merit) { 1002b4ed282SShri Abhyankar ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr); 1012b4ed282SShri Abhyankar *reason = SNES_DIVERGED_FNORM_NAN; 1022b4ed282SShri Abhyankar } else if (merit < snes->abstol) { 1032b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"Converged due to merit function %G < %G\n",merit,snes->abstol);CHKERRQ(ierr); 1042b4ed282SShri Abhyankar *reason = SNES_CONVERGED_FNORM_ABS; 1052b4ed282SShri Abhyankar } else if (snes->nfuncs >= snes->max_funcs) { 1062b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr); 1072b4ed282SShri Abhyankar *reason = SNES_DIVERGED_FUNCTION_COUNT; 1082b4ed282SShri Abhyankar } 1092b4ed282SShri Abhyankar 1102b4ed282SShri Abhyankar if (it && !*reason) { 1112b4ed282SShri Abhyankar if (merit < snes->ttol) { 1122b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"Converged due to merit function %G < %G (relative tolerance)\n",merit,snes->ttol);CHKERRQ(ierr); 1132b4ed282SShri Abhyankar *reason = SNES_CONVERGED_FNORM_RELATIVE; 1142b4ed282SShri Abhyankar } 1152b4ed282SShri Abhyankar } 1162b4ed282SShri Abhyankar PetscFunctionReturn(0); 1172b4ed282SShri Abhyankar } 1182b4ed282SShri Abhyankar 1192b4ed282SShri Abhyankar /* 1202b4ed282SShri Abhyankar SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem. 1212b4ed282SShri Abhyankar 1222b4ed282SShri Abhyankar Input Parameter: 1232b4ed282SShri Abhyankar . phi - the semismooth function 1242b4ed282SShri Abhyankar 1252b4ed282SShri Abhyankar Output Parameter: 1262b4ed282SShri Abhyankar . merit - the merit function 1272b4ed282SShri Abhyankar . phinorm - ||phi|| 1282b4ed282SShri Abhyankar 1292b4ed282SShri Abhyankar Notes: 1302b4ed282SShri Abhyankar The merit function for the mixed complementarity problem is defined as 1312b4ed282SShri Abhyankar merit = 0.5*phi^T*phi 1322b4ed282SShri Abhyankar */ 1332b4ed282SShri Abhyankar #undef __FUNCT__ 1342b4ed282SShri Abhyankar #define __FUNCT__ "SNESVIComputeMeritFunction" 1352b4ed282SShri Abhyankar static PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal* merit,PetscReal* phinorm) 1362b4ed282SShri Abhyankar { 1372b4ed282SShri Abhyankar PetscErrorCode ierr; 1382b4ed282SShri Abhyankar 1392b4ed282SShri Abhyankar PetscFunctionBegin; 1402b4ed282SShri Abhyankar ierr = VecNormBegin(phi,NORM_2,phinorm); 1412b4ed282SShri Abhyankar ierr = VecNormEnd(phi,NORM_2,phinorm); 1422b4ed282SShri Abhyankar 1432b4ed282SShri Abhyankar *merit = 0.5*(*phinorm)*(*phinorm); 1442b4ed282SShri Abhyankar PetscFunctionReturn(0); 1452b4ed282SShri Abhyankar } 1462b4ed282SShri Abhyankar 1472b4ed282SShri Abhyankar /* 1482b4ed282SShri Abhyankar ComputeFischerFunction - Computes the semismooth fischer burmeister function for a mixed complementarity equation. 1492b4ed282SShri Abhyankar 1502b4ed282SShri Abhyankar Notes: 1512b4ed282SShri Abhyankar The Fischer-Burmeister function is defined as 1522b4ed282SShri Abhyankar ff(a,b) := sqrt(a*a + b*b) - a - b 1532b4ed282SShri Abhyankar and is used reformulate a complementarity equation as a semismooth equation. 1542b4ed282SShri Abhyankar */ 1552b4ed282SShri Abhyankar 1562b4ed282SShri Abhyankar #undef __FUNCT__ 1572b4ed282SShri Abhyankar #define __FUNCT__ "ComputeFischerFunction" 1582b4ed282SShri Abhyankar static PetscErrorCode ComputeFischerFunction(PetscScalar a, PetscScalar b, PetscScalar* ff) 1592b4ed282SShri Abhyankar { 1602b4ed282SShri Abhyankar PetscFunctionBegin; 1612b4ed282SShri Abhyankar *ff = sqrt(a*a + b*b) - a - b; 1622b4ed282SShri Abhyankar PetscFunctionReturn(0); 1632b4ed282SShri Abhyankar } 1642b4ed282SShri Abhyankar 1652b4ed282SShri Abhyankar /* 1661f79c099SShri Abhyankar SNESVIComputeFunction - Reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear equations in semismooth form. 1672b4ed282SShri Abhyankar 1682b4ed282SShri Abhyankar Input Parameters: 1692b4ed282SShri Abhyankar . snes - the SNES context 1702b4ed282SShri Abhyankar . x - current iterate 1712b4ed282SShri Abhyankar . functx - user defined function context 1722b4ed282SShri Abhyankar 1732b4ed282SShri Abhyankar Output Parameters: 1742b4ed282SShri Abhyankar . phi - Semismooth function 1752b4ed282SShri Abhyankar 1762b4ed282SShri Abhyankar Notes: 1772b4ed282SShri Abhyankar The result of this function is done by cases: 1782b4ed282SShri Abhyankar + l[i] == -infinity, u[i] == infinity -- phi[i] = -f[i] 1792b4ed282SShri Abhyankar . l[i] == -infinity, u[i] finite -- phi[i] = ss(u[i]-x[i], -f[i]) 1802b4ed282SShri Abhyankar . l[i] finite, u[i] == infinity -- phi[i] = ss(x[i]-l[i], f[i]) 1812b4ed282SShri Abhyankar . l[i] finite < u[i] finite -- phi[i] = phi(x[i]-l[i], ss(u[i]-x[i], -f[u])) 1822b4ed282SShri Abhyankar - otherwise l[i] == u[i] -- phi[i] = l[i] - x[i] 1832b4ed282SShri Abhyankar ss is the semismoothing function used to reformulate the nonlinear equation in complementarity 1842b4ed282SShri Abhyankar form to semismooth form 1852b4ed282SShri Abhyankar 1862b4ed282SShri Abhyankar */ 1872b4ed282SShri Abhyankar #undef __FUNCT__ 1881f79c099SShri Abhyankar #define __FUNCT__ "SNESVIComputeFunction" 1891f79c099SShri Abhyankar static PetscErrorCode SNESVIComputeFunction(SNES snes,Vec X,Vec phi,void* functx) 1902b4ed282SShri Abhyankar { 1912b4ed282SShri Abhyankar PetscErrorCode ierr; 1922b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1932b4ed282SShri Abhyankar Vec Xl = vi->xl,Xu = vi->xu,F = snes->vec_func; 1942b4ed282SShri Abhyankar PetscScalar *phi_arr,*x_arr,*f_arr,*l,*u,t; 1952b4ed282SShri Abhyankar PetscInt i,nlocal; 1962b4ed282SShri Abhyankar 1972b4ed282SShri Abhyankar PetscFunctionBegin; 1982b4ed282SShri Abhyankar ierr = (*vi->computeuserfunction)(snes,X,F,functx);CHKERRQ(ierr); 1992b4ed282SShri Abhyankar ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 2002b4ed282SShri Abhyankar 2012b4ed282SShri Abhyankar ierr = VecGetArray(X,&x_arr);CHKERRQ(ierr); 2022b4ed282SShri Abhyankar ierr = VecGetArray(F,&f_arr);CHKERRQ(ierr); 2032b4ed282SShri Abhyankar ierr = VecGetArray(Xl,&l);CHKERRQ(ierr); 2042b4ed282SShri Abhyankar ierr = VecGetArray(Xu,&u);CHKERRQ(ierr); 2052b4ed282SShri Abhyankar ierr = VecGetArray(phi,&phi_arr);CHKERRQ(ierr); 2062b4ed282SShri Abhyankar 2072b4ed282SShri Abhyankar for (i=0;i < nlocal;i++) { 2082b4ed282SShri Abhyankar if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) { 2092b4ed282SShri Abhyankar phi_arr[i] = -f_arr[i]; 2102b4ed282SShri Abhyankar } 2112b4ed282SShri Abhyankar else if (l[i] <= PETSC_VI_NINF) { 2122b4ed282SShri Abhyankar t = u[i] - x_arr[i]; 2132b4ed282SShri Abhyankar ierr = ComputeFischerFunction(t,-f_arr[i],&phi_arr[i]);CHKERRQ(ierr); 2142b4ed282SShri Abhyankar phi_arr[i] = -phi_arr[i]; 2152b4ed282SShri Abhyankar } 2162b4ed282SShri Abhyankar else if (u[i] >= PETSC_VI_INF) { 2172b4ed282SShri Abhyankar t = x_arr[i] - l[i]; 2182b4ed282SShri Abhyankar ierr = ComputeFischerFunction(t,f_arr[i],&phi_arr[i]);CHKERRQ(ierr); 2192b4ed282SShri Abhyankar } 2202b4ed282SShri Abhyankar else if (l[i] == u[i]) { 2212b4ed282SShri Abhyankar phi_arr[i] = l[i] - x_arr[i]; 2222b4ed282SShri Abhyankar } 2232b4ed282SShri Abhyankar else { 2242b4ed282SShri Abhyankar t = u[i] - x_arr[i]; 2252b4ed282SShri Abhyankar ierr = ComputeFischerFunction(t,-f_arr[i],&phi_arr[i]); 2262b4ed282SShri Abhyankar t = x_arr[i] - l[i]; 2272b4ed282SShri Abhyankar ierr = ComputeFischerFunction(t,phi_arr[i],&phi_arr[i]); 2282b4ed282SShri Abhyankar } 2292b4ed282SShri Abhyankar } 2302b4ed282SShri Abhyankar 2312b4ed282SShri Abhyankar ierr = VecRestoreArray(X,&x_arr);CHKERRQ(ierr); 2322b4ed282SShri Abhyankar ierr = VecRestoreArray(F,&f_arr);CHKERRQ(ierr); 2332b4ed282SShri Abhyankar ierr = VecRestoreArray(Xl,&l);CHKERRQ(ierr); 2342b4ed282SShri Abhyankar ierr = VecRestoreArray(Xu,&u);CHKERRQ(ierr); 2352b4ed282SShri Abhyankar ierr = VecRestoreArray(phi,&phi_arr);CHKERRQ(ierr); 2362b4ed282SShri Abhyankar 2372b4ed282SShri Abhyankar PetscFunctionReturn(0); 2382b4ed282SShri Abhyankar } 2392b4ed282SShri Abhyankar 2402b4ed282SShri Abhyankar /* 241a79edbefSShri Abhyankar SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the 242a79edbefSShri Abhyankar the semismooth jacobian. 2432b4ed282SShri Abhyankar */ 2442b4ed282SShri Abhyankar #undef __FUNCT__ 245a79edbefSShri Abhyankar #define __FUNCT__ "SNESVIComputeBsubdifferentialVectors" 246a79edbefSShri Abhyankar PetscErrorCode SNESVIComputeBsubdifferentialVectors(SNES snes,Vec X,Vec F,Mat jac,Vec Da,Vec Db) 2472b4ed282SShri Abhyankar { 2482b4ed282SShri Abhyankar PetscErrorCode ierr; 2492b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 2502b4ed282SShri Abhyankar PetscScalar *l,*u,*x,*f,*da,*db,*z,*t,t1,t2,ci,di,ei; 2512b4ed282SShri Abhyankar PetscInt i,nlocal; 2522b4ed282SShri Abhyankar 2532b4ed282SShri Abhyankar PetscFunctionBegin; 2542b4ed282SShri Abhyankar 2552b4ed282SShri Abhyankar ierr = VecGetArray(X,&x);CHKERRQ(ierr); 2562b4ed282SShri Abhyankar ierr = VecGetArray(F,&f);CHKERRQ(ierr); 2572b4ed282SShri Abhyankar ierr = VecGetArray(vi->xl,&l);CHKERRQ(ierr); 2582b4ed282SShri Abhyankar ierr = VecGetArray(vi->xu,&u);CHKERRQ(ierr); 259a79edbefSShri Abhyankar ierr = VecGetArray(Da,&da);CHKERRQ(ierr); 260a79edbefSShri Abhyankar ierr = VecGetArray(Db,&db);CHKERRQ(ierr); 2612b4ed282SShri Abhyankar ierr = VecGetArray(vi->z,&z);CHKERRQ(ierr); 2622b4ed282SShri Abhyankar 2632b4ed282SShri Abhyankar ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 2642b4ed282SShri Abhyankar /* Set the elements of the vector z: 2652b4ed282SShri Abhyankar z[i] = 1 if (x[i] - l[i],f[i]) = (0,0) or (u[i] - x[i],f[i]) = (0,0) 2662b4ed282SShri Abhyankar else z[i] = 0 2672b4ed282SShri Abhyankar */ 2682b4ed282SShri Abhyankar for(i=0;i < nlocal;i++) { 2692b4ed282SShri Abhyankar da[i] = db[i] = z[i] = 0; 2702b4ed282SShri Abhyankar if(PetscAbsScalar(f[i]) <= vi->const_tol) { 2712b4ed282SShri Abhyankar if ((l[i] > PETSC_VI_NINF) && (PetscAbsScalar(x[i]-l[i]) <= vi->const_tol)) { 2722b4ed282SShri Abhyankar da[i] = 1; 2732b4ed282SShri Abhyankar z[i] = 1; 2742b4ed282SShri Abhyankar } 2752b4ed282SShri Abhyankar if ((u[i] < PETSC_VI_INF) && (PetscAbsScalar(u[i]-x[i]) <= vi->const_tol)) { 2762b4ed282SShri Abhyankar db[i] = 1; 2772b4ed282SShri Abhyankar z[i] = 1; 2782b4ed282SShri Abhyankar } 2792b4ed282SShri Abhyankar } 2802b4ed282SShri Abhyankar } 2812b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->z,&z);CHKERRQ(ierr); 282a79edbefSShri Abhyankar ierr = MatMult(jac,vi->z,vi->t);CHKERRQ(ierr); 2832b4ed282SShri Abhyankar ierr = VecGetArray(vi->t,&t);CHKERRQ(ierr); 2842b4ed282SShri Abhyankar /* Compute the elements of the diagonal perturbation vector Da and row scaling vector Db */ 2852b4ed282SShri Abhyankar for(i=0;i< nlocal;i++) { 2862b4ed282SShri Abhyankar /* Free variables */ 2872b4ed282SShri Abhyankar if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) { 2882b4ed282SShri Abhyankar da[i] = 0; db[i] = -1; 2892b4ed282SShri Abhyankar } 2902b4ed282SShri Abhyankar /* lower bounded variables */ 2912b4ed282SShri Abhyankar else if (u[i] >= PETSC_VI_INF) { 2922b4ed282SShri Abhyankar if (da[i] >= 1) { 2932b4ed282SShri Abhyankar t2 = PetscScalarNorm(1,t[i]); 2942b4ed282SShri Abhyankar da[i] = 1/t2 - 1; 2952b4ed282SShri Abhyankar db[i] = t[i]/t2 - 1; 2962b4ed282SShri Abhyankar } else { 2972b4ed282SShri Abhyankar t1 = x[i] - l[i]; 2982b4ed282SShri Abhyankar t2 = PetscScalarNorm(t1,f[i]); 2992b4ed282SShri Abhyankar da[i] = t1/t2 - 1; 3002b4ed282SShri Abhyankar db[i] = f[i]/t2 - 1; 3012b4ed282SShri Abhyankar } 3022b4ed282SShri Abhyankar } 3032b4ed282SShri Abhyankar /* upper bounded variables */ 3042b4ed282SShri Abhyankar else if (l[i] <= PETSC_VI_NINF) { 3052b4ed282SShri Abhyankar if (db[i] >= 1) { 3062b4ed282SShri Abhyankar t2 = PetscScalarNorm(1,t[i]); 3072b4ed282SShri Abhyankar da[i] = -1/t2 -1; 3082b4ed282SShri Abhyankar db[i] = -t[i]/t2 - 1; 3092b4ed282SShri Abhyankar } 3102b4ed282SShri Abhyankar else { 3112b4ed282SShri Abhyankar t1 = u[i] - x[i]; 3122b4ed282SShri Abhyankar t2 = PetscScalarNorm(t1,f[i]); 3132b4ed282SShri Abhyankar da[i] = t1/t2 - 1; 3142b4ed282SShri Abhyankar db[i] = -f[i]/t2 - 1; 3152b4ed282SShri Abhyankar } 3162b4ed282SShri Abhyankar } 3172b4ed282SShri Abhyankar /* Fixed variables */ 3182b4ed282SShri Abhyankar else if (l[i] == u[i]) { 3192b4ed282SShri Abhyankar da[i] = -1; 3202b4ed282SShri Abhyankar db[i] = 0; 3212b4ed282SShri Abhyankar } 3222b4ed282SShri Abhyankar /* Box constrained variables */ 3232b4ed282SShri Abhyankar else { 3242b4ed282SShri Abhyankar if (db[i] >= 1) { 3252b4ed282SShri Abhyankar t2 = PetscScalarNorm(1,t[i]); 3262b4ed282SShri Abhyankar ci = 1/t2 + 1; 3272b4ed282SShri Abhyankar di = t[i]/t2 + 1; 3282b4ed282SShri Abhyankar } 3292b4ed282SShri Abhyankar else { 3302b4ed282SShri Abhyankar t1 = x[i] - u[i]; 3312b4ed282SShri Abhyankar t2 = PetscScalarNorm(t1,f[i]); 3322b4ed282SShri Abhyankar ci = t1/t2 + 1; 3332b4ed282SShri Abhyankar di = f[i]/t2 + 1; 3342b4ed282SShri Abhyankar } 3352b4ed282SShri Abhyankar 3362b4ed282SShri Abhyankar if (da[i] >= 1) { 3372b4ed282SShri Abhyankar t1 = ci + di*t[i]; 3382b4ed282SShri Abhyankar t2 = PetscScalarNorm(1,t1); 3392b4ed282SShri Abhyankar t1 = t1/t2 - 1; 3402b4ed282SShri Abhyankar t2 = 1/t2 - 1; 3412b4ed282SShri Abhyankar } 3422b4ed282SShri Abhyankar else { 3432b4ed282SShri Abhyankar ierr = ComputeFischerFunction(u[i]-x[i],-f[i],&ei);CHKERRQ(ierr); 3442b4ed282SShri Abhyankar t2 = PetscScalarNorm(x[i]-l[i],ei); 3452b4ed282SShri Abhyankar t1 = ei/t2 - 1; 3462b4ed282SShri Abhyankar t2 = (x[i] - l[i])/t2 - 1; 3472b4ed282SShri Abhyankar } 3482b4ed282SShri Abhyankar 3492b4ed282SShri Abhyankar da[i] = t2 + t1*ci; 3502b4ed282SShri Abhyankar db[i] = t1*di; 3512b4ed282SShri Abhyankar } 3522b4ed282SShri Abhyankar } 3532b4ed282SShri Abhyankar ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 3542b4ed282SShri Abhyankar ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 3552b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->xl,&l);CHKERRQ(ierr); 3562b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->xu,&u);CHKERRQ(ierr); 357a79edbefSShri Abhyankar ierr = VecRestoreArray(Da,&da);CHKERRQ(ierr); 358a79edbefSShri Abhyankar ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr); 3592b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->t,&t);CHKERRQ(ierr); 3602b4ed282SShri Abhyankar 361a79edbefSShri Abhyankar PetscFunctionReturn(0); 362a79edbefSShri Abhyankar } 363a79edbefSShri Abhyankar 364a79edbefSShri Abhyankar /* 365a79edbefSShri 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. 366a79edbefSShri Abhyankar 367a79edbefSShri Abhyankar Input Parameters: 368a79edbefSShri Abhyankar . Da - Diagonal shift vector for the semismooth jacobian. 369a79edbefSShri Abhyankar . Db - Row scaling vector for the semismooth jacobian. 370a79edbefSShri Abhyankar 371a79edbefSShri Abhyankar Output Parameters: 372a79edbefSShri Abhyankar . jac - semismooth jacobian 373a79edbefSShri Abhyankar . jac_pre - optional preconditioning matrix 374a79edbefSShri Abhyankar 375a79edbefSShri Abhyankar Notes: 376a79edbefSShri Abhyankar The semismooth jacobian matrix is given by 377a79edbefSShri Abhyankar jac = Da + Db*jacfun 378a79edbefSShri Abhyankar where Db is the row scaling matrix stored as a vector, 379a79edbefSShri Abhyankar Da is the diagonal perturbation matrix stored as a vector 380a79edbefSShri Abhyankar and jacfun is the jacobian of the original nonlinear function. 381a79edbefSShri Abhyankar */ 382a79edbefSShri Abhyankar #undef __FUNCT__ 383a79edbefSShri Abhyankar #define __FUNCT__ "SNESVIComputeJacobian" 384a79edbefSShri Abhyankar PetscErrorCode SNESVIComputeJacobian(Mat jac, Mat jac_pre,Vec Da, Vec Db) 385a79edbefSShri Abhyankar { 386a79edbefSShri Abhyankar PetscErrorCode ierr; 387a79edbefSShri Abhyankar 3882b4ed282SShri Abhyankar /* Do row scaling and add diagonal perturbation */ 389a79edbefSShri Abhyankar ierr = MatDiagonalScale(jac,Db,PETSC_NULL);CHKERRQ(ierr); 390a79edbefSShri Abhyankar ierr = MatDiagonalSet(jac,Da,ADD_VALUES);CHKERRQ(ierr); 391a79edbefSShri Abhyankar if (jac != jac_pre) { /* If jac and jac_pre are different */ 392a79edbefSShri Abhyankar ierr = MatDiagonalScale(jac_pre,Db,PETSC_NULL); 393a79edbefSShri Abhyankar ierr = MatDiagonalSet(jac_pre,Da,ADD_VALUES);CHKERRQ(ierr); 3942b4ed282SShri Abhyankar } 3952b4ed282SShri Abhyankar PetscFunctionReturn(0); 3962b4ed282SShri Abhyankar } 3972b4ed282SShri Abhyankar 3982b4ed282SShri Abhyankar /* 399ee657d29SShri Abhyankar SNESVIComputeMeritFunctionGradient - Computes the gradient of the merit function psi. 400ee657d29SShri Abhyankar 401ee657d29SShri Abhyankar Input Parameters: 402ee657d29SShri Abhyankar phi - semismooth function. 403ee657d29SShri Abhyankar H - semismooth jacobian 404ee657d29SShri Abhyankar 405ee657d29SShri Abhyankar Output Parameters: 406ee657d29SShri Abhyankar dpsi - merit function gradient 407ee657d29SShri Abhyankar 408ee657d29SShri Abhyankar Notes: 409ee657d29SShri Abhyankar The merit function gradient is computed as follows 410ee657d29SShri Abhyankar dpsi = H^T*phi 411ee657d29SShri Abhyankar */ 412ee657d29SShri Abhyankar #undef __FUNCT__ 413ee657d29SShri Abhyankar #define __FUNCT__ "SNESVIComputeMeritFunctionGradient" 414ee657d29SShri Abhyankar PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi) 415ee657d29SShri Abhyankar { 416ee657d29SShri Abhyankar PetscErrorCode ierr; 417ee657d29SShri Abhyankar 418ee657d29SShri Abhyankar PetscFunctionBegin; 419ee657d29SShri Abhyankar ierr = MatMultTranspose(H,phi,dpsi); 420ee657d29SShri Abhyankar PetscFunctionReturn(0); 421ee657d29SShri Abhyankar } 422ee657d29SShri Abhyankar 423ee657d29SShri Abhyankar /* 4242b4ed282SShri Abhyankar SNESVIAdjustInitialGuess - Readjusts the initial guess to the SNES solver supplied by the user so that the initial guess lies inside the feasible region . 4252b4ed282SShri Abhyankar 4262b4ed282SShri Abhyankar Input Parameters: 4272b4ed282SShri Abhyankar . lb - lower bound. 4282b4ed282SShri Abhyankar . ub - upper bound. 4292b4ed282SShri Abhyankar 4302b4ed282SShri Abhyankar Output Parameters: 4312b4ed282SShri Abhyankar . X - the readjusted initial guess. 4322b4ed282SShri Abhyankar 4332b4ed282SShri Abhyankar Notes: 4342b4ed282SShri Abhyankar The readjusted initial guess X[i] = max(max(min(l[i],X[i]),min(X[i],u[i])),min(l[i],u[i])) 4352b4ed282SShri Abhyankar */ 4362b4ed282SShri Abhyankar #undef __FUNCT__ 4372b4ed282SShri Abhyankar #define __FUNCT__ "SNESVIAdjustInitialGuess" 4382b4ed282SShri Abhyankar PetscErrorCode SNESVIAdjustInitialGuess(Vec X, Vec lb, Vec ub) 4392b4ed282SShri Abhyankar { 4402b4ed282SShri Abhyankar PetscErrorCode ierr; 4412b4ed282SShri Abhyankar PetscInt i,nlocal; 4422b4ed282SShri Abhyankar PetscScalar *x,*l,*u; 4432b4ed282SShri Abhyankar 4442b4ed282SShri Abhyankar PetscFunctionBegin; 4452b4ed282SShri Abhyankar 4462b4ed282SShri Abhyankar ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 4472b4ed282SShri Abhyankar ierr = VecGetArray(X,&x);CHKERRQ(ierr); 4482b4ed282SShri Abhyankar ierr = VecGetArray(lb,&l);CHKERRQ(ierr); 4492b4ed282SShri Abhyankar ierr = VecGetArray(ub,&u);CHKERRQ(ierr); 4502b4ed282SShri Abhyankar 4512b4ed282SShri Abhyankar for(i = 0; i < nlocal; i++) { 4522b4ed282SShri Abhyankar x[i] = PetscMax(PetscMax(PetscMin(x[i],l[i]),PetscMin(x[i],u[i])),PetscMin(l[i],u[i])); 4532b4ed282SShri Abhyankar } 4542b4ed282SShri Abhyankar 4552b4ed282SShri Abhyankar ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 4562b4ed282SShri Abhyankar ierr = VecRestoreArray(lb,&l);CHKERRQ(ierr); 4572b4ed282SShri Abhyankar ierr = VecRestoreArray(ub,&u);CHKERRQ(ierr); 4582b4ed282SShri Abhyankar 4592b4ed282SShri Abhyankar PetscFunctionReturn(0); 4602b4ed282SShri Abhyankar } 4612b4ed282SShri Abhyankar 4622b4ed282SShri Abhyankar /* -------------------------------------------------------------------- 4632b4ed282SShri Abhyankar 4642b4ed282SShri Abhyankar This file implements a semismooth truncated Newton method with a line search, 4652b4ed282SShri Abhyankar for solving a system of nonlinear equations in complementarity form, using the KSP, Vec, 4662b4ed282SShri Abhyankar and Mat interfaces for linear solvers, vectors, and matrices, 4672b4ed282SShri Abhyankar respectively. 4682b4ed282SShri Abhyankar 4692b4ed282SShri Abhyankar The following basic routines are required for each nonlinear solver: 4702b4ed282SShri Abhyankar SNESCreate_XXX() - Creates a nonlinear solver context 4712b4ed282SShri Abhyankar SNESSetFromOptions_XXX() - Sets runtime options 4722b4ed282SShri Abhyankar SNESSolve_XXX() - Solves the nonlinear system 4732b4ed282SShri Abhyankar SNESDestroy_XXX() - Destroys the nonlinear solver context 4742b4ed282SShri Abhyankar The suffix "_XXX" denotes a particular implementation, in this case 4752b4ed282SShri Abhyankar we use _VI (e.g., SNESCreate_VI, SNESSolve_VI) for solving 4762b4ed282SShri Abhyankar systems of nonlinear equations with a line search (LS) method. 4772b4ed282SShri Abhyankar These routines are actually called via the common user interface 4782b4ed282SShri Abhyankar routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 4792b4ed282SShri Abhyankar SNESDestroy(), so the application code interface remains identical 4802b4ed282SShri Abhyankar for all nonlinear solvers. 4812b4ed282SShri Abhyankar 4822b4ed282SShri Abhyankar Another key routine is: 4832b4ed282SShri Abhyankar SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 4842b4ed282SShri Abhyankar by setting data structures and options. The interface routine SNESSetUp() 4852b4ed282SShri Abhyankar is not usually called directly by the user, but instead is called by 4862b4ed282SShri Abhyankar SNESSolve() if necessary. 4872b4ed282SShri Abhyankar 4882b4ed282SShri Abhyankar Additional basic routines are: 4892b4ed282SShri Abhyankar SNESView_XXX() - Prints details of runtime options that 4902b4ed282SShri Abhyankar have actually been used. 4912b4ed282SShri Abhyankar These are called by application codes via the interface routines 4922b4ed282SShri Abhyankar SNESView(). 4932b4ed282SShri Abhyankar 4942b4ed282SShri Abhyankar The various types of solvers (preconditioners, Krylov subspace methods, 4952b4ed282SShri Abhyankar nonlinear solvers, timesteppers) are all organized similarly, so the 4962b4ed282SShri Abhyankar above description applies to these categories also. 4972b4ed282SShri Abhyankar 4982b4ed282SShri Abhyankar -------------------------------------------------------------------- */ 4992b4ed282SShri Abhyankar /* 50069c03620SShri Abhyankar SNESSolveVI_SS - Solves the complementarity problem with a semismooth Newton 5012b4ed282SShri Abhyankar method using a line search. 5022b4ed282SShri Abhyankar 5032b4ed282SShri Abhyankar Input Parameters: 5042b4ed282SShri Abhyankar . snes - the SNES context 5052b4ed282SShri Abhyankar 5062b4ed282SShri Abhyankar Output Parameter: 5072b4ed282SShri Abhyankar . outits - number of iterations until termination 5082b4ed282SShri Abhyankar 5092b4ed282SShri Abhyankar Application Interface Routine: SNESSolve() 5102b4ed282SShri Abhyankar 5112b4ed282SShri Abhyankar Notes: 5122b4ed282SShri Abhyankar This implements essentially a semismooth Newton method with a 5132b4ed282SShri Abhyankar line search. By default a cubic backtracking line search 5142b4ed282SShri Abhyankar is employed, as described in the text "Numerical Methods for 5152b4ed282SShri Abhyankar Unconstrained Optimization and Nonlinear Equations" by Dennis 5162b4ed282SShri Abhyankar and Schnabel. 5172b4ed282SShri Abhyankar */ 5182b4ed282SShri Abhyankar #undef __FUNCT__ 51969c03620SShri Abhyankar #define __FUNCT__ "SNESSolveVI_SS" 52069c03620SShri Abhyankar PetscErrorCode SNESSolveVI_SS(SNES snes) 5212b4ed282SShri Abhyankar { 5222b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 5232b4ed282SShri Abhyankar PetscErrorCode ierr; 5242b4ed282SShri Abhyankar PetscInt maxits,i,lits; 5253b336b49SShri Abhyankar PetscBool lssucceed; 5262b4ed282SShri Abhyankar MatStructure flg = DIFFERENT_NONZERO_PATTERN; 5272b4ed282SShri Abhyankar PetscReal gnorm,xnorm=0,ynorm; 5282b4ed282SShri Abhyankar Vec Y,X,F,G,W; 5292b4ed282SShri Abhyankar KSPConvergedReason kspreason; 5302b4ed282SShri Abhyankar 5312b4ed282SShri Abhyankar PetscFunctionBegin; 5322b4ed282SShri Abhyankar snes->numFailures = 0; 5332b4ed282SShri Abhyankar snes->numLinearSolveFailures = 0; 5342b4ed282SShri Abhyankar snes->reason = SNES_CONVERGED_ITERATING; 5352b4ed282SShri Abhyankar 5362b4ed282SShri Abhyankar maxits = snes->max_its; /* maximum number of iterations */ 5372b4ed282SShri Abhyankar X = snes->vec_sol; /* solution vector */ 5382b4ed282SShri Abhyankar F = snes->vec_func; /* residual vector */ 5392b4ed282SShri Abhyankar Y = snes->work[0]; /* work vectors */ 5402b4ed282SShri Abhyankar G = snes->work[1]; 5412b4ed282SShri Abhyankar W = snes->work[2]; 5422b4ed282SShri Abhyankar 5432b4ed282SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 5442b4ed282SShri Abhyankar snes->iter = 0; 5452b4ed282SShri Abhyankar snes->norm = 0.0; 5462b4ed282SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 5472b4ed282SShri Abhyankar 5482b4ed282SShri Abhyankar ierr = SNESVIAdjustInitialGuess(X,vi->xl,vi->xu);CHKERRQ(ierr); 5492b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr); 5502b4ed282SShri Abhyankar if (snes->domainerror) { 5512b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 5522b4ed282SShri Abhyankar PetscFunctionReturn(0); 5532b4ed282SShri Abhyankar } 5542b4ed282SShri Abhyankar /* Compute Merit function */ 5552b4ed282SShri Abhyankar ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr); 5562b4ed282SShri Abhyankar 5572b4ed282SShri Abhyankar ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 5582b4ed282SShri Abhyankar ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 5592b4ed282SShri Abhyankar if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 5602b4ed282SShri Abhyankar 5612b4ed282SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 5622b4ed282SShri Abhyankar snes->norm = vi->phinorm; 5632b4ed282SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 5642b4ed282SShri Abhyankar SNESLogConvHistory(snes,vi->phinorm,0); 5652b4ed282SShri Abhyankar SNESMonitor(snes,0,vi->phinorm); 5662b4ed282SShri Abhyankar 5672b4ed282SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 5682b4ed282SShri Abhyankar snes->ttol = vi->phinorm*snes->rtol; 5692b4ed282SShri Abhyankar /* test convergence */ 5702b4ed282SShri Abhyankar ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 5712b4ed282SShri Abhyankar if (snes->reason) PetscFunctionReturn(0); 5722b4ed282SShri Abhyankar 5732b4ed282SShri Abhyankar for (i=0; i<maxits; i++) { 5742b4ed282SShri Abhyankar 5752b4ed282SShri Abhyankar /* Call general purpose update function */ 5762b4ed282SShri Abhyankar if (snes->ops->update) { 5772b4ed282SShri Abhyankar ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 5782b4ed282SShri Abhyankar } 5792b4ed282SShri Abhyankar 5802b4ed282SShri Abhyankar /* Solve J Y = Phi, where J is the semismooth jacobian */ 581a79edbefSShri Abhyankar /* Get the nonlinear function jacobian */ 5822b4ed282SShri Abhyankar ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 583a79edbefSShri Abhyankar /* Get the diagonal shift and row scaling vectors */ 584a79edbefSShri Abhyankar ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr); 585a79edbefSShri Abhyankar /* Compute the semismooth jacobian */ 586a79edbefSShri Abhyankar ierr = SNESVIComputeJacobian(snes->jacobian,snes->jacobian_pre,vi->Da,vi->Db);CHKERRQ(ierr); 587ee657d29SShri Abhyankar /* Compute the merit function gradient */ 588ee657d29SShri Abhyankar ierr = SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);CHKERRQ(ierr); 5892b4ed282SShri Abhyankar ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 5902b4ed282SShri Abhyankar ierr = SNES_KSPSolve(snes,snes->ksp,vi->phi,Y);CHKERRQ(ierr); 5912b4ed282SShri Abhyankar ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 5923b336b49SShri Abhyankar 5933b336b49SShri Abhyankar if (kspreason < 0) { 5942b4ed282SShri Abhyankar if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 5952b4ed282SShri 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); 5962b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 5972b4ed282SShri Abhyankar break; 5982b4ed282SShri Abhyankar } 5992b4ed282SShri Abhyankar } 6002b4ed282SShri Abhyankar ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 6012b4ed282SShri Abhyankar snes->linear_its += lits; 6022b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 6032b4ed282SShri Abhyankar /* 6042b4ed282SShri Abhyankar if (vi->precheckstep) { 605ace3abfcSBarry Smith PetscBool changed_y = PETSC_FALSE; 6062b4ed282SShri Abhyankar ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr); 6072b4ed282SShri Abhyankar } 6082b4ed282SShri Abhyankar 6092b4ed282SShri Abhyankar if (PetscLogPrintInfo){ 6102b4ed282SShri Abhyankar ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 6112b4ed282SShri Abhyankar } 6122b4ed282SShri Abhyankar */ 6132b4ed282SShri Abhyankar /* Compute a (scaled) negative update in the line search routine: 6142b4ed282SShri Abhyankar Y <- X - lambda*Y 6152b4ed282SShri Abhyankar and evaluate G = function(Y) (depends on the line search). 6162b4ed282SShri Abhyankar */ 6172b4ed282SShri Abhyankar ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 6182b4ed282SShri Abhyankar ynorm = 1; gnorm = vi->phinorm; 6192b4ed282SShri Abhyankar ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 6202b4ed282SShri Abhyankar ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 6212b4ed282SShri Abhyankar if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 6222b4ed282SShri Abhyankar if (snes->domainerror) { 6232b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 6242b4ed282SShri Abhyankar PetscFunctionReturn(0); 6252b4ed282SShri Abhyankar } 6262b4ed282SShri Abhyankar if (!lssucceed) { 6272b4ed282SShri Abhyankar if (++snes->numFailures >= snes->maxFailures) { 628ace3abfcSBarry Smith PetscBool ismin; 6292b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_LINE_SEARCH; 6302b4ed282SShri Abhyankar ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 6312b4ed282SShri Abhyankar if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 6322b4ed282SShri Abhyankar break; 6332b4ed282SShri Abhyankar } 6342b4ed282SShri Abhyankar } 6352b4ed282SShri Abhyankar /* Update function and solution vectors */ 6362b4ed282SShri Abhyankar vi->phinorm = gnorm; 6372b4ed282SShri Abhyankar vi->merit = 0.5*vi->phinorm*vi->phinorm; 6382b4ed282SShri Abhyankar ierr = VecCopy(G,vi->phi);CHKERRQ(ierr); 6392b4ed282SShri Abhyankar ierr = VecCopy(W,X);CHKERRQ(ierr); 6402b4ed282SShri Abhyankar /* Monitor convergence */ 6412b4ed282SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 6422b4ed282SShri Abhyankar snes->iter = i+1; 6432b4ed282SShri Abhyankar snes->norm = vi->phinorm; 6442b4ed282SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 6452b4ed282SShri Abhyankar SNESLogConvHistory(snes,snes->norm,lits); 6462b4ed282SShri Abhyankar SNESMonitor(snes,snes->iter,snes->norm); 6472b4ed282SShri Abhyankar /* Test for convergence, xnorm = || X || */ 6482b4ed282SShri Abhyankar if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 6492b4ed282SShri Abhyankar ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 6502b4ed282SShri Abhyankar if (snes->reason) break; 6512b4ed282SShri Abhyankar } 6522b4ed282SShri Abhyankar if (i == maxits) { 6532b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 6542b4ed282SShri Abhyankar if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 6552b4ed282SShri Abhyankar } 6562b4ed282SShri Abhyankar PetscFunctionReturn(0); 6572b4ed282SShri Abhyankar } 65869c03620SShri Abhyankar 65969c03620SShri Abhyankar #undef __FUNCT__ 66069c03620SShri Abhyankar #define __FUNCT__ "SNESVICreateIndexSets_AS" 661a79edbefSShri Abhyankar PetscErrorCode SNESVICreateIndexSets_AS(SNES snes,Vec Db,PetscReal thresh,IS* ISact,IS* ISinact) 66269c03620SShri Abhyankar { 66369c03620SShri Abhyankar PetscErrorCode ierr; 66469c03620SShri Abhyankar PetscInt i,nlocal,ilow,ihigh,nloc_isact=0,nloc_isinact=0; 66569c03620SShri Abhyankar PetscInt *idx_act,*idx_inact,i1=0,i2=0; 66669c03620SShri Abhyankar PetscScalar *db; 66769c03620SShri Abhyankar 66869c03620SShri Abhyankar PetscFunctionBegin; 66969c03620SShri Abhyankar 67069c03620SShri Abhyankar ierr = VecGetLocalSize(Db,&nlocal);CHKERRQ(ierr); 67169c03620SShri Abhyankar ierr = VecGetOwnershipRange(Db,&ilow,&ihigh);CHKERRQ(ierr); 67269c03620SShri Abhyankar ierr = VecGetArray(Db,&db);CHKERRQ(ierr); 67369c03620SShri Abhyankar /* Compute the sizes of the active and inactive sets */ 67469c03620SShri Abhyankar for (i=0; i < nlocal;i++) { 675fca35906SShri Abhyankar if (PetscAbsScalar(db[i]) <= thresh) nloc_isact++; 676fca35906SShri Abhyankar else nloc_isinact++; 67769c03620SShri Abhyankar } 67869c03620SShri Abhyankar ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr); 67969c03620SShri Abhyankar ierr = PetscMalloc(nloc_isinact*sizeof(PetscInt),&idx_inact);CHKERRQ(ierr); 68069c03620SShri Abhyankar 68169c03620SShri Abhyankar /* Creating the indexing arrays */ 682fca35906SShri Abhyankar for(i=0; i < nlocal; i++) { 683fca35906SShri Abhyankar if (PetscAbsScalar(db[i]) <= thresh) idx_act[i1++] = ilow+i; 684fca35906SShri Abhyankar else idx_inact[i2++] = ilow+i; 68569c03620SShri Abhyankar } 68669c03620SShri Abhyankar 68769c03620SShri Abhyankar /* Create the index sets */ 68869c03620SShri Abhyankar ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isact,idx_act,PETSC_COPY_VALUES,ISact);CHKERRQ(ierr); 68969c03620SShri Abhyankar ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isinact,idx_inact,PETSC_COPY_VALUES,ISinact);CHKERRQ(ierr); 69069c03620SShri Abhyankar 69169c03620SShri Abhyankar ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr); 69269c03620SShri Abhyankar ierr = PetscFree(idx_act);CHKERRQ(ierr); 69369c03620SShri Abhyankar ierr = PetscFree(idx_inact);CHKERRQ(ierr); 69469c03620SShri Abhyankar PetscFunctionReturn(0); 69569c03620SShri Abhyankar } 69669c03620SShri Abhyankar 697d950fb63SShri Abhyankar #undef __FUNCT__ 698d950fb63SShri Abhyankar #define __FUNCT__ "SNESVICreateIndexSets_RS" 699d950fb63SShri Abhyankar PetscErrorCode SNESVICreateIndexSets_RS(SNES snes,Vec X,Vec Xl, Vec Xu,IS* ISact,IS* ISinact) 700d950fb63SShri Abhyankar { 701d950fb63SShri Abhyankar PetscErrorCode ierr; 702d950fb63SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 703d950fb63SShri Abhyankar PetscInt i,nlocal,ilow,ihigh,nloc_isact=0,nloc_isinact=0; 704d950fb63SShri Abhyankar PetscInt *idx_act,*idx_inact,i1=0,i2=0; 705d950fb63SShri Abhyankar PetscScalar *x,*l,*u,*f; 706d950fb63SShri Abhyankar Vec F = snes->vec_func; 707d950fb63SShri Abhyankar 708d950fb63SShri Abhyankar PetscFunctionBegin; 709d950fb63SShri Abhyankar 710d950fb63SShri Abhyankar ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 711d950fb63SShri Abhyankar ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr); 712d950fb63SShri Abhyankar ierr = VecGetArray(X,&x);CHKERRQ(ierr); 713d950fb63SShri Abhyankar ierr = VecGetArray(Xl,&l);CHKERRQ(ierr); 714d950fb63SShri Abhyankar ierr = VecGetArray(Xu,&u);CHKERRQ(ierr); 715d950fb63SShri Abhyankar ierr = VecGetArray(F,&f);CHKERRQ(ierr); 716d950fb63SShri Abhyankar /* Compute the sizes of the active and inactive sets */ 717d950fb63SShri Abhyankar for (i=0; i < nlocal;i++) { 718d950fb63SShri Abhyankar if (PetscAbsScalar(x[i] - l[i]) <= vi->const_tol && f[i] > vi->const_tol) nloc_isact++; 719d950fb63SShri Abhyankar else if (PetscAbsScalar(u[i] - x[i]) <= vi->const_tol && f[i] < vi->const_tol) nloc_isact++; 720d950fb63SShri Abhyankar else nloc_isinact++; 721d950fb63SShri Abhyankar } 722d950fb63SShri Abhyankar ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr); 723d950fb63SShri Abhyankar ierr = PetscMalloc(nloc_isinact*sizeof(PetscInt),&idx_inact);CHKERRQ(ierr); 724d950fb63SShri Abhyankar 725d950fb63SShri Abhyankar /* Creating the indexing arrays */ 726d950fb63SShri Abhyankar for(i=0; i < nlocal; i++) { 727d950fb63SShri Abhyankar if (PetscAbsScalar(x[i] - l[i]) <= vi->const_tol && f[i] > vi->const_tol) idx_act[i1++] = ilow+i; 728d950fb63SShri Abhyankar else if (PetscAbsScalar(u[i] - x[i]) <= vi->const_tol && f[i] < vi->const_tol) idx_act[i1++] = ilow+i; 729d950fb63SShri Abhyankar else idx_inact[i2++] = ilow+i; 730d950fb63SShri Abhyankar } 731d950fb63SShri Abhyankar 732d950fb63SShri Abhyankar /* Create the index sets */ 733d950fb63SShri Abhyankar ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isact,idx_act,PETSC_COPY_VALUES,ISact);CHKERRQ(ierr); 734d950fb63SShri Abhyankar ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isinact,idx_inact,PETSC_COPY_VALUES,ISinact);CHKERRQ(ierr); 735d950fb63SShri Abhyankar 736d950fb63SShri Abhyankar ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 737d950fb63SShri Abhyankar ierr = VecRestoreArray(Xl,&l);CHKERRQ(ierr); 738d950fb63SShri Abhyankar ierr = VecRestoreArray(Xu,&u);CHKERRQ(ierr); 739d950fb63SShri Abhyankar ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 740d950fb63SShri Abhyankar ierr = PetscFree(idx_act);CHKERRQ(ierr); 741d950fb63SShri Abhyankar ierr = PetscFree(idx_inact);CHKERRQ(ierr); 742d950fb63SShri Abhyankar PetscFunctionReturn(0); 743d950fb63SShri Abhyankar } 744d950fb63SShri Abhyankar 745dbd914b8SShri Abhyankar /* Create active and inactive set vectors. The local size of this vector is set and petsc computes the global size */ 746dbd914b8SShri Abhyankar #undef __FUNCT__ 747dbd914b8SShri Abhyankar #define __FUNCT__ "SNESVICreateVectors_AS" 748dbd914b8SShri Abhyankar PetscErrorCode SNESVICreateVectors_AS(SNES snes,PetscInt n,Vec* newv) 749dbd914b8SShri Abhyankar { 750dbd914b8SShri Abhyankar PetscErrorCode ierr; 751dbd914b8SShri Abhyankar Vec v; 752dbd914b8SShri Abhyankar 753dbd914b8SShri Abhyankar PetscFunctionBegin; 754dbd914b8SShri Abhyankar ierr = VecCreate(((PetscObject)snes)->comm,&v);CHKERRQ(ierr); 755dbd914b8SShri Abhyankar ierr = VecSetSizes(v,n,PETSC_DECIDE);CHKERRQ(ierr); 756dbd914b8SShri Abhyankar ierr = VecSetFromOptions(v);CHKERRQ(ierr); 757dbd914b8SShri Abhyankar *newv = v; 758dbd914b8SShri Abhyankar 759dbd914b8SShri Abhyankar PetscFunctionReturn(0); 760dbd914b8SShri Abhyankar } 761dbd914b8SShri Abhyankar 762732bb160SShri Abhyankar /* Resets the snes PC and KSP when the active set sizes change */ 763732bb160SShri Abhyankar #undef __FUNCT__ 764732bb160SShri Abhyankar #define __FUNCT__ "SNESVIResetPCandKSP" 765732bb160SShri Abhyankar PetscErrorCode SNESVIResetPCandKSP(SNES snes,Mat Amat,Mat Pmat) 766732bb160SShri Abhyankar { 767732bb160SShri Abhyankar PetscErrorCode ierr; 768732bb160SShri Abhyankar KSP kspnew,snesksp; 769732bb160SShri Abhyankar PC pcnew; 770732bb160SShri Abhyankar const MatSolverPackage stype; 771dbd914b8SShri Abhyankar 772732bb160SShri Abhyankar PetscFunctionBegin; 773732bb160SShri Abhyankar /* The active and inactive set sizes have changed so need to create a new snes->ksp object */ 774732bb160SShri Abhyankar ierr = SNESGetKSP(snes,&snesksp);CHKERRQ(ierr); 775732bb160SShri Abhyankar ierr = KSPCreate(((PetscObject)snes)->comm,&kspnew);CHKERRQ(ierr); 776732bb160SShri Abhyankar /* Copy over snes->ksp info */ 777732bb160SShri Abhyankar kspnew->pc_side = snesksp->pc_side; 778732bb160SShri Abhyankar kspnew->rtol = snesksp->rtol; 779732bb160SShri Abhyankar kspnew->abstol = snesksp->abstol; 780732bb160SShri Abhyankar kspnew->max_it = snesksp->max_it; 781732bb160SShri Abhyankar ierr = KSPSetType(kspnew,((PetscObject)snesksp)->type_name);CHKERRQ(ierr); 782732bb160SShri Abhyankar ierr = KSPGetPC(kspnew,&pcnew);CHKERRQ(ierr); 783732bb160SShri Abhyankar ierr = PCSetType(kspnew->pc,((PetscObject)snesksp->pc)->type_name);CHKERRQ(ierr); 784732bb160SShri Abhyankar ierr = PCSetOperators(kspnew->pc,Amat,Pmat,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 785732bb160SShri Abhyankar ierr = PCFactorGetMatSolverPackage(snesksp->pc,&stype);CHKERRQ(ierr); 786732bb160SShri Abhyankar ierr = PCFactorSetMatSolverPackage(kspnew->pc,stype);CHKERRQ(ierr); 787732bb160SShri Abhyankar ierr = KSPDestroy(snesksp);CHKERRQ(ierr); 788732bb160SShri Abhyankar snes->ksp = kspnew; 789732bb160SShri Abhyankar ierr = PetscLogObjectParent(snes,kspnew);CHKERRQ(ierr); 790732bb160SShri Abhyankar ierr = KSPSetFromOptions(kspnew);CHKERRQ(ierr); 791732bb160SShri Abhyankar PetscFunctionReturn(0); 792732bb160SShri Abhyankar } 79369c03620SShri Abhyankar /* Variational Inequality solver using active set method */ 79469c03620SShri Abhyankar #undef __FUNCT__ 79569c03620SShri Abhyankar #define __FUNCT__ "SNESSolveVI_AS" 79669c03620SShri Abhyankar PetscErrorCode SNESSolveVI_AS(SNES snes) 79769c03620SShri Abhyankar { 79869c03620SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 79969c03620SShri Abhyankar PetscErrorCode ierr; 800408da460SShri Abhyankar PetscInt maxits,i,lits,Nis_act=0; 8013b336b49SShri Abhyankar PetscBool lssucceed; 80269c03620SShri Abhyankar MatStructure flg = DIFFERENT_NONZERO_PATTERN; 80369c03620SShri Abhyankar PetscReal gnorm,xnorm=0,ynorm; 80469c03620SShri Abhyankar Vec Y,X,F,G,W; 80569c03620SShri Abhyankar KSPConvergedReason kspreason; 80669c03620SShri Abhyankar 80769c03620SShri Abhyankar PetscFunctionBegin; 80869c03620SShri Abhyankar snes->numFailures = 0; 80969c03620SShri Abhyankar snes->numLinearSolveFailures = 0; 81069c03620SShri Abhyankar snes->reason = SNES_CONVERGED_ITERATING; 81169c03620SShri Abhyankar 81269c03620SShri Abhyankar maxits = snes->max_its; /* maximum number of iterations */ 81369c03620SShri Abhyankar X = snes->vec_sol; /* solution vector */ 81469c03620SShri Abhyankar F = snes->vec_func; /* residual vector */ 81569c03620SShri Abhyankar Y = snes->work[0]; /* work vectors */ 81669c03620SShri Abhyankar G = snes->work[1]; 81769c03620SShri Abhyankar W = snes->work[2]; 81869c03620SShri Abhyankar 81969c03620SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 82069c03620SShri Abhyankar snes->iter = 0; 82169c03620SShri Abhyankar snes->norm = 0.0; 82269c03620SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 82369c03620SShri Abhyankar 82469c03620SShri Abhyankar ierr = SNESVIAdjustInitialGuess(X,vi->xl,vi->xu);CHKERRQ(ierr); 82569c03620SShri Abhyankar ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr); 82669c03620SShri Abhyankar if (snes->domainerror) { 82769c03620SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 82869c03620SShri Abhyankar PetscFunctionReturn(0); 82969c03620SShri Abhyankar } 83069c03620SShri Abhyankar /* Compute Merit function */ 83169c03620SShri Abhyankar ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr); 83269c03620SShri Abhyankar 83369c03620SShri Abhyankar ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 83469c03620SShri Abhyankar ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 83569c03620SShri Abhyankar if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 83669c03620SShri Abhyankar 83769c03620SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 83869c03620SShri Abhyankar snes->norm = vi->phinorm; 83969c03620SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 84069c03620SShri Abhyankar SNESLogConvHistory(snes,vi->phinorm,0); 84169c03620SShri Abhyankar SNESMonitor(snes,0,vi->phinorm); 84269c03620SShri Abhyankar 84369c03620SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 84469c03620SShri Abhyankar snes->ttol = vi->phinorm*snes->rtol; 84569c03620SShri Abhyankar /* test convergence */ 84669c03620SShri Abhyankar ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 84769c03620SShri Abhyankar if (snes->reason) PetscFunctionReturn(0); 84869c03620SShri Abhyankar 84969c03620SShri Abhyankar for (i=0; i<maxits; i++) { 85069c03620SShri Abhyankar 851a79edbefSShri Abhyankar IS IS_act,IS_inact; /* _act -> active set _inact -> inactive set */ 852a79edbefSShri Abhyankar PetscReal thresh,J_norm1; 853dbd914b8SShri Abhyankar VecScatter scat_act,scat_inact; 854408da460SShri Abhyankar PetscInt nis_act,nis_inact,Nis_act_prev; 855dbd914b8SShri Abhyankar Vec Da_act,Da_inact,Db_inact; 856dbd914b8SShri Abhyankar Vec Y_act,Y_inact,phi_act,phi_inact; 857d76c674bSShri Abhyankar Mat jac_inact_inact,jac_inact_act,prejac_inact_inact; 858a79edbefSShri Abhyankar 85969c03620SShri Abhyankar /* Call general purpose update function */ 86069c03620SShri Abhyankar if (snes->ops->update) { 86169c03620SShri Abhyankar ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 86269c03620SShri Abhyankar } 86369c03620SShri Abhyankar ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 86469c03620SShri Abhyankar /* Compute the threshold value for creating active and inactive sets */ 86569c03620SShri Abhyankar ierr = MatNorm(snes->jacobian,NORM_1,&J_norm1);CHKERRQ(ierr); 86669c03620SShri Abhyankar thresh = PetscMin(vi->merit,1e-2)/(1+J_norm1); 867a79edbefSShri Abhyankar 868a79edbefSShri Abhyankar /* Compute B-subdifferential vectors Da and Db */ 869a79edbefSShri Abhyankar ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr); 870dbd914b8SShri Abhyankar 87169c03620SShri Abhyankar /* Create active and inactive index sets */ 87269c03620SShri Abhyankar ierr = SNESVICreateIndexSets_AS(snes,vi->Db,thresh,&IS_act,&IS_inact);CHKERRQ(ierr); 87369c03620SShri Abhyankar 874408da460SShri Abhyankar Nis_act_prev = Nis_act; 875408da460SShri Abhyankar /* Get sizes of active and inactive sets */ 876dbd914b8SShri Abhyankar ierr = ISGetLocalSize(IS_act,&nis_act);CHKERRQ(ierr); 877dbd914b8SShri Abhyankar ierr = ISGetLocalSize(IS_inact,&nis_inact);CHKERRQ(ierr); 878408da460SShri Abhyankar ierr = ISGetSize(IS_act,&Nis_act);CHKERRQ(ierr); 879dbd914b8SShri Abhyankar 880d76c674bSShri Abhyankar ierr = PetscPrintf(PETSC_COMM_WORLD,"Size of active set = %d, size of inactive set = %d\n",nis_act,nis_inact);CHKERRQ(ierr); 881d76c674bSShri Abhyankar 882dbd914b8SShri Abhyankar /* Create active and inactive set vectors */ 883dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_act,&Da_act);CHKERRQ(ierr); 884dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_inact,&Da_inact);CHKERRQ(ierr); 885dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_inact,&Db_inact);CHKERRQ(ierr); 886dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_act,&phi_act);CHKERRQ(ierr); 887dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_inact,&phi_inact);CHKERRQ(ierr); 888dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_act,&Y_act);CHKERRQ(ierr); 889dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_inact,&Y_inact);CHKERRQ(ierr); 890dbd914b8SShri Abhyankar 891dbd914b8SShri Abhyankar /* Create inactive set submatrices */ 892dbd914b8SShri Abhyankar ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_act,MAT_INITIAL_MATRIX,&jac_inact_act);CHKERRQ(ierr); 893dbd914b8SShri Abhyankar ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);CHKERRQ(ierr); 894dbd914b8SShri Abhyankar 895dbd914b8SShri Abhyankar /* Create scatter contexts */ 896dbd914b8SShri Abhyankar ierr = VecScatterCreate(vi->Da,IS_act,Da_act,PETSC_NULL,&scat_act);CHKERRQ(ierr); 897dbd914b8SShri Abhyankar ierr = VecScatterCreate(vi->Da,IS_inact,Da_inact,PETSC_NULL,&scat_inact);CHKERRQ(ierr); 898dbd914b8SShri Abhyankar 899dbd914b8SShri Abhyankar /* Do a vec scatter to active and inactive set vectors */ 900dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_act,vi->Da,Da_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 901dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_act,vi->Da,Da_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 902dbd914b8SShri Abhyankar 903dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_inact,vi->Da,Da_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 904dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_inact,vi->Da,Da_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 905dbd914b8SShri Abhyankar 906dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_inact,vi->Db,Db_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 907dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_inact,vi->Db,Db_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 908dbd914b8SShri Abhyankar 909dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_act,vi->phi,phi_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 910dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_act,vi->phi,phi_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 911dbd914b8SShri Abhyankar 912dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_inact,vi->phi,phi_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 913dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_inact,vi->phi,phi_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 914dbd914b8SShri Abhyankar 915dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 916dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 917dbd914b8SShri Abhyankar 918dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 919dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 920dbd914b8SShri Abhyankar 921dbd914b8SShri Abhyankar /* Active set direction */ 922dbd914b8SShri Abhyankar ierr = VecPointwiseDivide(Y_act,phi_act,Da_act);CHKERRQ(ierr); 923d76c674bSShri Abhyankar /* inactive set jacobian and preconditioner */ 924dbd914b8SShri Abhyankar ierr = VecPointwiseDivide(Da_inact,Da_inact,Db_inact);CHKERRQ(ierr); 925dbd914b8SShri Abhyankar ierr = MatDiagonalSet(jac_inact_inact,Da_inact,ADD_VALUES);CHKERRQ(ierr); 926d76c674bSShri Abhyankar if (snes->jacobian != snes->jacobian_pre) { 927d76c674bSShri Abhyankar ierr = MatGetSubMatrix(snes->jacobian_pre,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&prejac_inact_inact);CHKERRQ(ierr); 928d76c674bSShri Abhyankar ierr = MatDiagonalSet(prejac_inact_inact,Da_inact,ADD_VALUES);CHKERRQ(ierr); 929d76c674bSShri Abhyankar } else prejac_inact_inact = jac_inact_inact; 930d76c674bSShri Abhyankar 931dbd914b8SShri Abhyankar /* right hand side */ 932dbd914b8SShri Abhyankar ierr = VecPointwiseDivide(phi_inact,phi_inact,Db_inact);CHKERRQ(ierr); 933dbd914b8SShri Abhyankar ierr = MatMult(jac_inact_act,Y_act,Db_inact);CHKERRQ(ierr); 934dbd914b8SShri Abhyankar ierr = VecAXPY(phi_inact,-1.0,Db_inact);CHKERRQ(ierr); 935dbd914b8SShri Abhyankar 936408da460SShri Abhyankar if ((i != 0) && (Nis_act != Nis_act_prev)) { 937732bb160SShri Abhyankar ierr = SNESVIResetPCandKSP(snes,jac_inact_inact,prejac_inact_inact);CHKERRQ(ierr); 938408da460SShri Abhyankar } 939408da460SShri Abhyankar 940d76c674bSShri Abhyankar ierr = KSPSetOperators(snes->ksp,jac_inact_inact,prejac_inact_inact,flg);CHKERRQ(ierr); 941dbd914b8SShri Abhyankar ierr = SNES_KSPSolve(snes,snes->ksp,phi_inact,Y_inact);CHKERRQ(ierr); 942dbd914b8SShri Abhyankar ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 9433b336b49SShri Abhyankar if (kspreason < 0) { 9443b336b49SShri Abhyankar if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 9453b336b49SShri 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); 9463b336b49SShri Abhyankar snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 9473b336b49SShri Abhyankar break; 9483b336b49SShri Abhyankar } 9493b336b49SShri Abhyankar } 950dbd914b8SShri Abhyankar 951dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 952dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 953dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 954dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 955dbd914b8SShri Abhyankar 956dbd914b8SShri Abhyankar ierr = VecDestroy(Da_act);CHKERRQ(ierr); 957dbd914b8SShri Abhyankar ierr = VecDestroy(Da_inact);CHKERRQ(ierr); 958dbd914b8SShri Abhyankar ierr = VecDestroy(Db_inact);CHKERRQ(ierr); 959dbd914b8SShri Abhyankar ierr = VecDestroy(phi_act);CHKERRQ(ierr); 960dbd914b8SShri Abhyankar ierr = VecDestroy(phi_inact);CHKERRQ(ierr); 961dbd914b8SShri Abhyankar ierr = VecDestroy(Y_act);CHKERRQ(ierr); 962dbd914b8SShri Abhyankar ierr = VecDestroy(Y_inact);CHKERRQ(ierr); 963dbd914b8SShri Abhyankar ierr = VecScatterDestroy(scat_act);CHKERRQ(ierr); 964dbd914b8SShri Abhyankar ierr = VecScatterDestroy(scat_inact);CHKERRQ(ierr); 965730c24dcSShri Abhyankar ierr = ISDestroy(IS_act);CHKERRQ(ierr); 966730c24dcSShri Abhyankar ierr = ISDestroy(IS_inact);CHKERRQ(ierr); 967bac2f86eSShri Abhyankar ierr = MatDestroy(jac_inact_act);CHKERRQ(ierr); 968bac2f86eSShri Abhyankar ierr = MatDestroy(jac_inact_inact);CHKERRQ(ierr); 969d76c674bSShri Abhyankar if (snes->jacobian != snes->jacobian_pre) { 970d76c674bSShri Abhyankar ierr = MatDestroy(prejac_inact_inact);CHKERRQ(ierr); 971d76c674bSShri Abhyankar } 972730c24dcSShri Abhyankar 97369c03620SShri Abhyankar ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 97469c03620SShri Abhyankar snes->linear_its += lits; 97569c03620SShri Abhyankar ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 97669c03620SShri Abhyankar /* 97769c03620SShri Abhyankar if (vi->precheckstep) { 97869c03620SShri Abhyankar PetscBool changed_y = PETSC_FALSE; 97969c03620SShri Abhyankar ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr); 98069c03620SShri Abhyankar } 98169c03620SShri Abhyankar 98269c03620SShri Abhyankar if (PetscLogPrintInfo){ 98369c03620SShri Abhyankar ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 98469c03620SShri Abhyankar } 98569c03620SShri Abhyankar */ 98669c03620SShri Abhyankar /* Compute a (scaled) negative update in the line search routine: 98769c03620SShri Abhyankar Y <- X - lambda*Y 98869c03620SShri Abhyankar and evaluate G = function(Y) (depends on the line search). 98969c03620SShri Abhyankar */ 99069c03620SShri Abhyankar ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 99169c03620SShri Abhyankar ynorm = 1; gnorm = vi->phinorm; 9923b336b49SShri Abhyankar ierr = SNESVIComputeJacobian(snes->jacobian,snes->jacobian_pre,vi->Da,vi->Db);CHKERRQ(ierr); 993ee657d29SShri Abhyankar ierr = SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);CHKERRQ(ierr); 99469c03620SShri Abhyankar ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 99569c03620SShri Abhyankar ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 99669c03620SShri Abhyankar if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 99769c03620SShri Abhyankar if (snes->domainerror) { 99869c03620SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 99969c03620SShri Abhyankar PetscFunctionReturn(0); 100069c03620SShri Abhyankar } 100169c03620SShri Abhyankar if (!lssucceed) { 100269c03620SShri Abhyankar if (++snes->numFailures >= snes->maxFailures) { 100369c03620SShri Abhyankar PetscBool ismin; 100469c03620SShri Abhyankar snes->reason = SNES_DIVERGED_LINE_SEARCH; 100569c03620SShri Abhyankar ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 100669c03620SShri Abhyankar if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 100769c03620SShri Abhyankar break; 100869c03620SShri Abhyankar } 100969c03620SShri Abhyankar } 101069c03620SShri Abhyankar /* Update function and solution vectors */ 101169c03620SShri Abhyankar vi->phinorm = gnorm; 101269c03620SShri Abhyankar vi->merit = 0.5*vi->phinorm*vi->phinorm; 101369c03620SShri Abhyankar ierr = VecCopy(G,vi->phi);CHKERRQ(ierr); 101469c03620SShri Abhyankar ierr = VecCopy(W,X);CHKERRQ(ierr); 101569c03620SShri Abhyankar /* Monitor convergence */ 101669c03620SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 101769c03620SShri Abhyankar snes->iter = i+1; 101869c03620SShri Abhyankar snes->norm = vi->phinorm; 101969c03620SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 102069c03620SShri Abhyankar SNESLogConvHistory(snes,snes->norm,lits); 102169c03620SShri Abhyankar SNESMonitor(snes,snes->iter,snes->norm); 102269c03620SShri Abhyankar /* Test for convergence, xnorm = || X || */ 102369c03620SShri Abhyankar if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 102469c03620SShri Abhyankar ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 102569c03620SShri Abhyankar if (snes->reason) break; 102669c03620SShri Abhyankar } 102769c03620SShri Abhyankar if (i == maxits) { 102869c03620SShri Abhyankar ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 102969c03620SShri Abhyankar if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 103069c03620SShri Abhyankar } 103169c03620SShri Abhyankar PetscFunctionReturn(0); 103269c03620SShri Abhyankar } 103369c03620SShri Abhyankar 1034d950fb63SShri Abhyankar /* Variational Inequality solver using active set method */ 1035d950fb63SShri Abhyankar #undef __FUNCT__ 1036d950fb63SShri Abhyankar #define __FUNCT__ "SNESSolveVI_RS" 1037d950fb63SShri Abhyankar PetscErrorCode SNESSolveVI_RS(SNES snes) 1038d950fb63SShri Abhyankar { 1039d950fb63SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1040d950fb63SShri Abhyankar PetscErrorCode ierr; 1041d950fb63SShri Abhyankar PetscInt maxits,i,lits,Nis_act=0; 1042d950fb63SShri Abhyankar PetscBool lssucceed; 1043d950fb63SShri Abhyankar MatStructure flg = DIFFERENT_NONZERO_PATTERN; 1044d950fb63SShri Abhyankar PetscReal gnorm,xnorm=0,ynorm; 1045d950fb63SShri Abhyankar Vec Y,X,F,G,W; 1046d950fb63SShri Abhyankar KSPConvergedReason kspreason; 1047d950fb63SShri Abhyankar 1048d950fb63SShri Abhyankar PetscFunctionBegin; 1049d950fb63SShri Abhyankar snes->numFailures = 0; 1050d950fb63SShri Abhyankar snes->numLinearSolveFailures = 0; 1051d950fb63SShri Abhyankar snes->reason = SNES_CONVERGED_ITERATING; 1052d950fb63SShri Abhyankar 1053d950fb63SShri Abhyankar maxits = snes->max_its; /* maximum number of iterations */ 1054d950fb63SShri Abhyankar X = snes->vec_sol; /* solution vector */ 1055d950fb63SShri Abhyankar F = snes->vec_func; /* residual vector */ 1056d950fb63SShri Abhyankar Y = snes->work[0]; /* work vectors */ 1057d950fb63SShri Abhyankar G = snes->work[1]; 1058d950fb63SShri Abhyankar W = snes->work[2]; 1059d950fb63SShri Abhyankar 1060d950fb63SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1061d950fb63SShri Abhyankar snes->iter = 0; 1062d950fb63SShri Abhyankar snes->norm = 0.0; 1063d950fb63SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1064d950fb63SShri Abhyankar 1065d950fb63SShri Abhyankar ierr = SNESVIAdjustInitialGuess(X,vi->xl,vi->xu);CHKERRQ(ierr); 1066d950fb63SShri Abhyankar ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr); 1067d950fb63SShri Abhyankar if (snes->domainerror) { 1068d950fb63SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1069d950fb63SShri Abhyankar PetscFunctionReturn(0); 1070d950fb63SShri Abhyankar } 1071d950fb63SShri Abhyankar /* Compute Merit function */ 1072d950fb63SShri Abhyankar ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr); 1073d950fb63SShri Abhyankar 1074d950fb63SShri Abhyankar ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 1075d950fb63SShri Abhyankar ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 1076d950fb63SShri Abhyankar if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1077d950fb63SShri Abhyankar 1078d950fb63SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1079d950fb63SShri Abhyankar snes->norm = vi->phinorm; 1080d950fb63SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1081d950fb63SShri Abhyankar SNESLogConvHistory(snes,vi->phinorm,0); 1082d950fb63SShri Abhyankar SNESMonitor(snes,0,vi->phinorm); 1083d950fb63SShri Abhyankar 1084d950fb63SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 1085d950fb63SShri Abhyankar snes->ttol = vi->phinorm*snes->rtol; 1086d950fb63SShri Abhyankar /* test convergence */ 1087d950fb63SShri Abhyankar ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1088d950fb63SShri Abhyankar if (snes->reason) PetscFunctionReturn(0); 1089d950fb63SShri Abhyankar 1090d950fb63SShri Abhyankar for (i=0; i<maxits; i++) { 1091d950fb63SShri Abhyankar 1092d950fb63SShri Abhyankar IS IS_act,IS_inact; /* _act -> active set _inact -> inactive set */ 1093d950fb63SShri Abhyankar VecScatter scat_act,scat_inact; 1094d950fb63SShri Abhyankar PetscInt nis_act,nis_inact,Nis_act_prev; 1095d950fb63SShri Abhyankar Vec Y_act,Y_inact,phi_inact; 1096d950fb63SShri Abhyankar Mat jac_inact_inact,prejac_inact_inact; 1097d950fb63SShri Abhyankar 1098d950fb63SShri Abhyankar /* Call general purpose update function */ 1099d950fb63SShri Abhyankar if (snes->ops->update) { 1100d950fb63SShri Abhyankar ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1101d950fb63SShri Abhyankar } 1102d950fb63SShri Abhyankar ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 1103d950fb63SShri Abhyankar ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr); 1104d950fb63SShri Abhyankar ierr = SNESVIComputeJacobian(snes->jacobian,snes->jacobian_pre,vi->Da,vi->Db);CHKERRQ(ierr); 1105ee657d29SShri Abhyankar ierr = SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);CHKERRQ(ierr); 1106d950fb63SShri Abhyankar /* Create active and inactive index sets */ 1107d950fb63SShri Abhyankar ierr = SNESVICreateIndexSets_RS(snes,X,vi->xl,vi->xu,&IS_act,&IS_inact);CHKERRQ(ierr); 1108d950fb63SShri Abhyankar 1109d950fb63SShri Abhyankar Nis_act_prev = Nis_act; 1110d950fb63SShri Abhyankar /* Get sizes of active and inactive sets */ 1111d950fb63SShri Abhyankar ierr = ISGetLocalSize(IS_act,&nis_act);CHKERRQ(ierr); 1112d950fb63SShri Abhyankar ierr = ISGetLocalSize(IS_inact,&nis_inact);CHKERRQ(ierr); 1113d950fb63SShri Abhyankar ierr = ISGetSize(IS_act,&Nis_act);CHKERRQ(ierr); 1114d950fb63SShri Abhyankar 1115d950fb63SShri Abhyankar ierr = PetscPrintf(PETSC_COMM_WORLD,"Size of active set = %d, size of inactive set = %d\n",nis_act,nis_inact);CHKERRQ(ierr); 1116d950fb63SShri Abhyankar 1117d950fb63SShri Abhyankar /* Create active and inactive set vectors */ 1118d950fb63SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_inact,&phi_inact);CHKERRQ(ierr); 1119d950fb63SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_act,&Y_act);CHKERRQ(ierr); 1120d950fb63SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_inact,&Y_inact);CHKERRQ(ierr); 1121d950fb63SShri Abhyankar 1122d950fb63SShri Abhyankar /* Create inactive set submatrices */ 1123d950fb63SShri Abhyankar ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);CHKERRQ(ierr); 1124d950fb63SShri Abhyankar 1125d950fb63SShri Abhyankar /* Create scatter contexts */ 1126d950fb63SShri Abhyankar ierr = VecScatterCreate(Y,IS_act,Y_act,PETSC_NULL,&scat_act);CHKERRQ(ierr); 1127d950fb63SShri Abhyankar ierr = VecScatterCreate(Y,IS_inact,Y_inact,PETSC_NULL,&scat_inact);CHKERRQ(ierr); 1128d950fb63SShri Abhyankar 1129d950fb63SShri Abhyankar /* Do a vec scatter to active and inactive set vectors */ 1130d950fb63SShri Abhyankar ierr = VecScatterBegin(scat_inact,vi->phi,phi_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1131d950fb63SShri Abhyankar ierr = VecScatterEnd(scat_inact,vi->phi,phi_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1132d950fb63SShri Abhyankar 1133d950fb63SShri Abhyankar ierr = VecScatterBegin(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1134d950fb63SShri Abhyankar ierr = VecScatterEnd(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1135d950fb63SShri Abhyankar 1136d950fb63SShri Abhyankar ierr = VecScatterBegin(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1137d950fb63SShri Abhyankar ierr = VecScatterEnd(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1138d950fb63SShri Abhyankar 1139d950fb63SShri Abhyankar /* Active set direction = 0*/ 1140d950fb63SShri Abhyankar ierr = VecSet(Y_act,0);CHKERRQ(ierr); 1141d950fb63SShri Abhyankar if (snes->jacobian != snes->jacobian_pre) { 1142d950fb63SShri Abhyankar ierr = MatGetSubMatrix(snes->jacobian_pre,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&prejac_inact_inact);CHKERRQ(ierr); 1143d950fb63SShri Abhyankar } else prejac_inact_inact = jac_inact_inact; 1144d950fb63SShri Abhyankar 1145d950fb63SShri Abhyankar if ((i != 0) && (Nis_act != Nis_act_prev)) { 1146732bb160SShri Abhyankar ierr = SNESVIResetPCandKSP(snes,jac_inact_inact,prejac_inact_inact);CHKERRQ(ierr); 1147d950fb63SShri Abhyankar } 1148d950fb63SShri Abhyankar 1149d950fb63SShri Abhyankar ierr = KSPSetOperators(snes->ksp,jac_inact_inact,prejac_inact_inact,flg);CHKERRQ(ierr); 1150d950fb63SShri Abhyankar ierr = SNES_KSPSolve(snes,snes->ksp,phi_inact,Y_inact);CHKERRQ(ierr); 1151d950fb63SShri Abhyankar ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 1152d950fb63SShri Abhyankar if (kspreason < 0) { 1153d950fb63SShri Abhyankar if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 1154d950fb63SShri 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); 1155d950fb63SShri Abhyankar snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 1156d950fb63SShri Abhyankar break; 1157d950fb63SShri Abhyankar } 1158d950fb63SShri Abhyankar } 1159d950fb63SShri Abhyankar 1160d950fb63SShri Abhyankar ierr = VecScatterBegin(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1161d950fb63SShri Abhyankar ierr = VecScatterEnd(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1162d950fb63SShri Abhyankar ierr = VecScatterBegin(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1163d950fb63SShri Abhyankar ierr = VecScatterEnd(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1164d950fb63SShri Abhyankar 1165d950fb63SShri Abhyankar ierr = VecDestroy(phi_inact);CHKERRQ(ierr); 1166d950fb63SShri Abhyankar ierr = VecDestroy(Y_act);CHKERRQ(ierr); 1167d950fb63SShri Abhyankar ierr = VecDestroy(Y_inact);CHKERRQ(ierr); 1168d950fb63SShri Abhyankar ierr = VecScatterDestroy(scat_act);CHKERRQ(ierr); 1169d950fb63SShri Abhyankar ierr = VecScatterDestroy(scat_inact);CHKERRQ(ierr); 1170d950fb63SShri Abhyankar ierr = ISDestroy(IS_act);CHKERRQ(ierr); 1171d950fb63SShri Abhyankar ierr = ISDestroy(IS_inact);CHKERRQ(ierr); 1172d950fb63SShri Abhyankar ierr = MatDestroy(jac_inact_inact);CHKERRQ(ierr); 1173d950fb63SShri Abhyankar if (snes->jacobian != snes->jacobian_pre) { 1174d950fb63SShri Abhyankar ierr = MatDestroy(prejac_inact_inact);CHKERRQ(ierr); 1175d950fb63SShri Abhyankar } 1176d950fb63SShri Abhyankar 1177d950fb63SShri Abhyankar ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 1178d950fb63SShri Abhyankar snes->linear_its += lits; 1179d950fb63SShri Abhyankar ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 1180d950fb63SShri Abhyankar /* 1181d950fb63SShri Abhyankar if (vi->precheckstep) { 1182d950fb63SShri Abhyankar PetscBool changed_y = PETSC_FALSE; 1183d950fb63SShri Abhyankar ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr); 1184d950fb63SShri Abhyankar } 1185d950fb63SShri Abhyankar 1186d950fb63SShri Abhyankar if (PetscLogPrintInfo){ 1187d950fb63SShri Abhyankar ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 1188d950fb63SShri Abhyankar } 1189d950fb63SShri Abhyankar */ 1190d950fb63SShri Abhyankar /* Compute a (scaled) negative update in the line search routine: 1191d950fb63SShri Abhyankar Y <- X - lambda*Y 1192d950fb63SShri Abhyankar and evaluate G = function(Y) (depends on the line search). 1193d950fb63SShri Abhyankar */ 1194d950fb63SShri Abhyankar ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 1195d950fb63SShri Abhyankar ynorm = 1; gnorm = vi->phinorm; 1196d950fb63SShri Abhyankar ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 1197d950fb63SShri Abhyankar ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 1198d950fb63SShri Abhyankar if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 1199d950fb63SShri Abhyankar if (snes->domainerror) { 1200d950fb63SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1201d950fb63SShri Abhyankar PetscFunctionReturn(0); 1202d950fb63SShri Abhyankar } 1203d950fb63SShri Abhyankar if (!lssucceed) { 1204d950fb63SShri Abhyankar if (++snes->numFailures >= snes->maxFailures) { 1205d950fb63SShri Abhyankar PetscBool ismin; 1206d950fb63SShri Abhyankar snes->reason = SNES_DIVERGED_LINE_SEARCH; 1207d950fb63SShri Abhyankar ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 1208d950fb63SShri Abhyankar if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 1209d950fb63SShri Abhyankar break; 1210d950fb63SShri Abhyankar } 1211d950fb63SShri Abhyankar } 1212d950fb63SShri Abhyankar /* Update function and solution vectors */ 1213d950fb63SShri Abhyankar vi->phinorm = gnorm; 1214d950fb63SShri Abhyankar vi->merit = 0.5*vi->phinorm*vi->phinorm; 1215d950fb63SShri Abhyankar ierr = VecCopy(G,vi->phi);CHKERRQ(ierr); 1216d950fb63SShri Abhyankar ierr = VecCopy(W,X);CHKERRQ(ierr); 1217d950fb63SShri Abhyankar /* Monitor convergence */ 1218d950fb63SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1219d950fb63SShri Abhyankar snes->iter = i+1; 1220d950fb63SShri Abhyankar snes->norm = vi->phinorm; 1221d950fb63SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1222d950fb63SShri Abhyankar SNESLogConvHistory(snes,snes->norm,lits); 1223d950fb63SShri Abhyankar SNESMonitor(snes,snes->iter,snes->norm); 1224d950fb63SShri Abhyankar /* Test for convergence, xnorm = || X || */ 1225d950fb63SShri Abhyankar if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 1226d950fb63SShri Abhyankar ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1227d950fb63SShri Abhyankar if (snes->reason) break; 1228d950fb63SShri Abhyankar } 1229d950fb63SShri Abhyankar if (i == maxits) { 1230d950fb63SShri Abhyankar ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 1231d950fb63SShri Abhyankar if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1232d950fb63SShri Abhyankar } 1233d950fb63SShri Abhyankar PetscFunctionReturn(0); 1234d950fb63SShri Abhyankar } 1235d950fb63SShri Abhyankar 12362b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 12372b4ed282SShri Abhyankar /* 12382b4ed282SShri Abhyankar SNESSetUp_VI - Sets up the internal data structures for the later use 12392b4ed282SShri Abhyankar of the SNESVI nonlinear solver. 12402b4ed282SShri Abhyankar 12412b4ed282SShri Abhyankar Input Parameter: 12422b4ed282SShri Abhyankar . snes - the SNES context 12432b4ed282SShri Abhyankar . x - the solution vector 12442b4ed282SShri Abhyankar 12452b4ed282SShri Abhyankar Application Interface Routine: SNESSetUp() 12462b4ed282SShri Abhyankar 12472b4ed282SShri Abhyankar Notes: 12482b4ed282SShri Abhyankar For basic use of the SNES solvers, the user need not explicitly call 12492b4ed282SShri Abhyankar SNESSetUp(), since these actions will automatically occur during 12502b4ed282SShri Abhyankar the call to SNESSolve(). 12512b4ed282SShri Abhyankar */ 12522b4ed282SShri Abhyankar #undef __FUNCT__ 12532b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetUp_VI" 12542b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes) 12552b4ed282SShri Abhyankar { 12562b4ed282SShri Abhyankar PetscErrorCode ierr; 12572b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*) snes->data; 12582b4ed282SShri Abhyankar PetscInt i_start[3],i_end[3]; 12592b4ed282SShri Abhyankar 12602b4ed282SShri Abhyankar PetscFunctionBegin; 12612b4ed282SShri Abhyankar if (!snes->vec_sol_update) { 12622b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr); 12632b4ed282SShri Abhyankar ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr); 12642b4ed282SShri Abhyankar } 12652b4ed282SShri Abhyankar if (!snes->work) { 12662b4ed282SShri Abhyankar snes->nwork = 3; 12672b4ed282SShri Abhyankar ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 12682b4ed282SShri Abhyankar ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr); 12692b4ed282SShri Abhyankar } 1270ee657d29SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->dpsi);CHKERRQ(ierr); 12712b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->phi); CHKERRQ(ierr); 12722b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->Da); CHKERRQ(ierr); 12732b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->Db); CHKERRQ(ierr); 12742b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->z);CHKERRQ(ierr); 12752b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->t); CHKERRQ(ierr); 12762b4ed282SShri Abhyankar 12772b4ed282SShri Abhyankar /* If the lower and upper bound on variables are not set, set it to 12782b4ed282SShri Abhyankar -Inf and Inf */ 12792b4ed282SShri Abhyankar if (!vi->xl && !vi->xu) { 12802b4ed282SShri Abhyankar vi->usersetxbounds = PETSC_FALSE; 12812b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->xl); CHKERRQ(ierr); 12822b4ed282SShri Abhyankar ierr = VecSet(vi->xl,PETSC_VI_NINF);CHKERRQ(ierr); 12832b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->xu); CHKERRQ(ierr); 12842b4ed282SShri Abhyankar ierr = VecSet(vi->xu,PETSC_VI_INF);CHKERRQ(ierr); 12852b4ed282SShri Abhyankar } else { 12862b4ed282SShri Abhyankar /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */ 12872b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr); 12882b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(vi->xl,i_start+1,i_end+1);CHKERRQ(ierr); 12892b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(vi->xu,i_start+2,i_end+2);CHKERRQ(ierr); 12902b4ed282SShri 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])) 12912b4ed282SShri 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."); 12922b4ed282SShri Abhyankar } 12932b4ed282SShri Abhyankar 12942b4ed282SShri Abhyankar vi->computeuserfunction = snes->ops->computefunction; 12951f79c099SShri Abhyankar snes->ops->computefunction = SNESVIComputeFunction; 1296a79edbefSShri Abhyankar 12972b4ed282SShri Abhyankar PetscFunctionReturn(0); 12982b4ed282SShri Abhyankar } 12992b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 13002b4ed282SShri Abhyankar /* 13012b4ed282SShri Abhyankar SNESDestroy_VI - Destroys the private SNES_VI context that was created 13022b4ed282SShri Abhyankar with SNESCreate_VI(). 13032b4ed282SShri Abhyankar 13042b4ed282SShri Abhyankar Input Parameter: 13052b4ed282SShri Abhyankar . snes - the SNES context 13062b4ed282SShri Abhyankar 13072b4ed282SShri Abhyankar Application Interface Routine: SNESDestroy() 13082b4ed282SShri Abhyankar */ 13092b4ed282SShri Abhyankar #undef __FUNCT__ 13102b4ed282SShri Abhyankar #define __FUNCT__ "SNESDestroy_VI" 13112b4ed282SShri Abhyankar PetscErrorCode SNESDestroy_VI(SNES snes) 13122b4ed282SShri Abhyankar { 13132b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*) snes->data; 13142b4ed282SShri Abhyankar PetscErrorCode ierr; 13152b4ed282SShri Abhyankar 13162b4ed282SShri Abhyankar PetscFunctionBegin; 13172b4ed282SShri Abhyankar if (snes->vec_sol_update) { 13182b4ed282SShri Abhyankar ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr); 13192b4ed282SShri Abhyankar snes->vec_sol_update = PETSC_NULL; 13202b4ed282SShri Abhyankar } 13212b4ed282SShri Abhyankar if (snes->nwork) { 13222b4ed282SShri Abhyankar ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 13232b4ed282SShri Abhyankar snes->nwork = 0; 13242b4ed282SShri Abhyankar snes->work = PETSC_NULL; 13252b4ed282SShri Abhyankar } 13262b4ed282SShri Abhyankar 13272b4ed282SShri Abhyankar /* clear vectors */ 1328ee657d29SShri Abhyankar ierr = VecDestroy(vi->dpsi);CHKERRQ(ierr); 13292b4ed282SShri Abhyankar ierr = VecDestroy(vi->phi); CHKERRQ(ierr); 13302b4ed282SShri Abhyankar ierr = VecDestroy(vi->Da); CHKERRQ(ierr); 13312b4ed282SShri Abhyankar ierr = VecDestroy(vi->Db); CHKERRQ(ierr); 13322b4ed282SShri Abhyankar ierr = VecDestroy(vi->z); CHKERRQ(ierr); 13332b4ed282SShri Abhyankar ierr = VecDestroy(vi->t); CHKERRQ(ierr); 13342b4ed282SShri Abhyankar if (!vi->usersetxbounds) { 13352b4ed282SShri Abhyankar ierr = VecDestroy(vi->xl); CHKERRQ(ierr); 13362b4ed282SShri Abhyankar ierr = VecDestroy(vi->xu); CHKERRQ(ierr); 13372b4ed282SShri Abhyankar } 1338c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1339c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr); 13402b4ed282SShri Abhyankar } 13412b4ed282SShri Abhyankar ierr = PetscFree(snes->data);CHKERRQ(ierr); 13422b4ed282SShri Abhyankar 13432b4ed282SShri Abhyankar /* clear composed functions */ 13442b4ed282SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr); 1345c92abb8aSShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr); 13462b4ed282SShri Abhyankar PetscFunctionReturn(0); 13472b4ed282SShri Abhyankar } 13482b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 13492b4ed282SShri Abhyankar #undef __FUNCT__ 13502b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNo_VI" 13512b4ed282SShri Abhyankar 13522b4ed282SShri Abhyankar /* 13532b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchNo routine in snes/impls/ls/ls.c 13542b4ed282SShri Abhyankar 13552b4ed282SShri Abhyankar */ 1356ace3abfcSBarry 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) 13572b4ed282SShri Abhyankar { 13582b4ed282SShri Abhyankar PetscErrorCode ierr; 13592b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1360ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 13612b4ed282SShri Abhyankar 13622b4ed282SShri Abhyankar PetscFunctionBegin; 13632b4ed282SShri Abhyankar *flag = PETSC_TRUE; 13642b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 13652b4ed282SShri Abhyankar ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 13662b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 13672b4ed282SShri Abhyankar if (vi->postcheckstep) { 13682b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 13692b4ed282SShri Abhyankar } 13702b4ed282SShri Abhyankar if (changed_y) { 13712b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 13722b4ed282SShri Abhyankar } 13732b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 13742b4ed282SShri Abhyankar if (!snes->domainerror) { 13752b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 13762b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 13772b4ed282SShri Abhyankar } 1378c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1379c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 1380c92abb8aSShri Abhyankar } 13812b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 13822b4ed282SShri Abhyankar PetscFunctionReturn(0); 13832b4ed282SShri Abhyankar } 13842b4ed282SShri Abhyankar 13852b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 13862b4ed282SShri Abhyankar #undef __FUNCT__ 13872b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNoNorms_VI" 13882b4ed282SShri Abhyankar 13892b4ed282SShri Abhyankar /* 13902b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c 13912b4ed282SShri Abhyankar */ 1392ace3abfcSBarry 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) 13932b4ed282SShri Abhyankar { 13942b4ed282SShri Abhyankar PetscErrorCode ierr; 13952b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1396ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 13972b4ed282SShri Abhyankar 13982b4ed282SShri Abhyankar PetscFunctionBegin; 13992b4ed282SShri Abhyankar *flag = PETSC_TRUE; 14002b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 14012b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 14022b4ed282SShri Abhyankar if (vi->postcheckstep) { 14032b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 14042b4ed282SShri Abhyankar } 14052b4ed282SShri Abhyankar if (changed_y) { 14062b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 14072b4ed282SShri Abhyankar } 14082b4ed282SShri Abhyankar 14092b4ed282SShri Abhyankar /* don't evaluate function the last time through */ 14102b4ed282SShri Abhyankar if (snes->iter < snes->max_its-1) { 14112b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 14122b4ed282SShri Abhyankar } 14132b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 14142b4ed282SShri Abhyankar PetscFunctionReturn(0); 14152b4ed282SShri Abhyankar } 14162b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 14172b4ed282SShri Abhyankar #undef __FUNCT__ 14182b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchCubic_VI" 14192b4ed282SShri Abhyankar /* 14202b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchCubic in snes/impls/ls/ls.c 14212b4ed282SShri Abhyankar */ 1422ace3abfcSBarry 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) 14232b4ed282SShri Abhyankar { 14242b4ed282SShri Abhyankar /* 14252b4ed282SShri Abhyankar Note that for line search purposes we work with with the related 14262b4ed282SShri Abhyankar minimization problem: 14272b4ed282SShri Abhyankar min z(x): R^n -> R, 14282b4ed282SShri Abhyankar where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 14292b4ed282SShri Abhyankar This function z(x) is same as the merit function for the complementarity problem. 14302b4ed282SShri Abhyankar */ 14312b4ed282SShri Abhyankar 14322b4ed282SShri Abhyankar PetscReal initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 14332b4ed282SShri Abhyankar PetscReal minlambda,lambda,lambdatemp; 14342b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 14352b4ed282SShri Abhyankar PetscScalar cinitslope; 14362b4ed282SShri Abhyankar #endif 14372b4ed282SShri Abhyankar PetscErrorCode ierr; 14382b4ed282SShri Abhyankar PetscInt count; 14392b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1440ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 14412b4ed282SShri Abhyankar MPI_Comm comm; 14422b4ed282SShri Abhyankar 14432b4ed282SShri Abhyankar PetscFunctionBegin; 14442b4ed282SShri Abhyankar ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 14452b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 14462b4ed282SShri Abhyankar *flag = PETSC_TRUE; 14472b4ed282SShri Abhyankar 14482b4ed282SShri Abhyankar ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 14492b4ed282SShri Abhyankar if (*ynorm == 0.0) { 1450c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1451c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 14522b4ed282SShri Abhyankar } 14532b4ed282SShri Abhyankar *gnorm = fnorm; 14542b4ed282SShri Abhyankar ierr = VecCopy(x,w);CHKERRQ(ierr); 14552b4ed282SShri Abhyankar ierr = VecCopy(f,g);CHKERRQ(ierr); 14562b4ed282SShri Abhyankar *flag = PETSC_FALSE; 14572b4ed282SShri Abhyankar goto theend1; 14582b4ed282SShri Abhyankar } 14592b4ed282SShri Abhyankar if (*ynorm > vi->maxstep) { /* Step too big, so scale back */ 1460c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1461c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Scaling step by %G old ynorm %G\n",vi->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr); 14622b4ed282SShri Abhyankar } 14632b4ed282SShri Abhyankar ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr); 14642b4ed282SShri Abhyankar *ynorm = vi->maxstep; 14652b4ed282SShri Abhyankar } 14662b4ed282SShri Abhyankar ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 14672b4ed282SShri Abhyankar minlambda = vi->minlambda/rellength; 14682b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 1469*89d9e806SShri Abhyankar ierr = VecDot(vi->dpsi,y,&cinitslope);CHKERRQ(ierr); 14702b4ed282SShri Abhyankar initslope = PetscRealPart(cinitslope); 14712b4ed282SShri Abhyankar #else 1472*89d9e806SShri Abhyankar ierr = VecDot(vi->dpsi,y,&initslope);CHKERRQ(ierr); 14732b4ed282SShri Abhyankar #endif 14742b4ed282SShri Abhyankar if (initslope > 0.0) initslope = -initslope; 14752b4ed282SShri Abhyankar if (initslope == 0.0) initslope = -1.0; 14762b4ed282SShri Abhyankar 14772b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 14782b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 14792b4ed282SShri Abhyankar ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 14802b4ed282SShri Abhyankar *flag = PETSC_FALSE; 14812b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 14822b4ed282SShri Abhyankar goto theend1; 14832b4ed282SShri Abhyankar } 14842b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 14852b4ed282SShri Abhyankar if (snes->domainerror) { 14862b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 14872b4ed282SShri Abhyankar PetscFunctionReturn(0); 14882b4ed282SShri Abhyankar } 14892b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 14902b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 14912b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 14922b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */ 1493c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1494c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 14952b4ed282SShri Abhyankar } 14962b4ed282SShri Abhyankar goto theend1; 14972b4ed282SShri Abhyankar } 14982b4ed282SShri Abhyankar 14992b4ed282SShri Abhyankar /* Fit points with quadratic */ 15002b4ed282SShri Abhyankar lambda = 1.0; 15012b4ed282SShri Abhyankar lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 15022b4ed282SShri Abhyankar lambdaprev = lambda; 15032b4ed282SShri Abhyankar gnormprev = *gnorm; 15042b4ed282SShri Abhyankar if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 15052b4ed282SShri Abhyankar if (lambdatemp <= .1*lambda) lambda = .1*lambda; 15062b4ed282SShri Abhyankar else lambda = lambdatemp; 15072b4ed282SShri Abhyankar 15082b4ed282SShri Abhyankar ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 15092b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 15102b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 15112b4ed282SShri Abhyankar *flag = PETSC_FALSE; 15122b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 15132b4ed282SShri Abhyankar goto theend1; 15142b4ed282SShri Abhyankar } 15152b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 15162b4ed282SShri Abhyankar if (snes->domainerror) { 15172b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 15182b4ed282SShri Abhyankar PetscFunctionReturn(0); 15192b4ed282SShri Abhyankar } 15202b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 15212b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1522c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1523c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr); 15242b4ed282SShri Abhyankar } 15252b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */ 1526c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1527c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 15282b4ed282SShri Abhyankar } 15292b4ed282SShri Abhyankar goto theend1; 15302b4ed282SShri Abhyankar } 15312b4ed282SShri Abhyankar 15322b4ed282SShri Abhyankar /* Fit points with cubic */ 15332b4ed282SShri Abhyankar count = 1; 15342b4ed282SShri Abhyankar while (PETSC_TRUE) { 15352b4ed282SShri Abhyankar if (lambda <= minlambda) { 1536c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1537c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 1538c92abb8aSShri 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); 15392b4ed282SShri Abhyankar } 15402b4ed282SShri Abhyankar *flag = PETSC_FALSE; 15412b4ed282SShri Abhyankar break; 15422b4ed282SShri Abhyankar } 15432b4ed282SShri Abhyankar t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 15442b4ed282SShri Abhyankar t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 15452b4ed282SShri Abhyankar a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 15462b4ed282SShri Abhyankar b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 15472b4ed282SShri Abhyankar d = b*b - 3*a*initslope; 15482b4ed282SShri Abhyankar if (d < 0.0) d = 0.0; 15492b4ed282SShri Abhyankar if (a == 0.0) { 15502b4ed282SShri Abhyankar lambdatemp = -initslope/(2.0*b); 15512b4ed282SShri Abhyankar } else { 15522b4ed282SShri Abhyankar lambdatemp = (-b + sqrt(d))/(3.0*a); 15532b4ed282SShri Abhyankar } 15542b4ed282SShri Abhyankar lambdaprev = lambda; 15552b4ed282SShri Abhyankar gnormprev = *gnorm; 15562b4ed282SShri Abhyankar if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 15572b4ed282SShri Abhyankar if (lambdatemp <= .1*lambda) lambda = .1*lambda; 15582b4ed282SShri Abhyankar else lambda = lambdatemp; 15592b4ed282SShri Abhyankar ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 15602b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 15612b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 15622b4ed282SShri 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); 15632b4ed282SShri Abhyankar *flag = PETSC_FALSE; 15642b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 15652b4ed282SShri Abhyankar break; 15662b4ed282SShri Abhyankar } 15672b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 15682b4ed282SShri Abhyankar if (snes->domainerror) { 15692b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 15702b4ed282SShri Abhyankar PetscFunctionReturn(0); 15712b4ed282SShri Abhyankar } 15722b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 15732b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 15742b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* is reduction enough? */ 1575c92abb8aSShri Abhyankar if (vi->lsmonitor) { 15762b4ed282SShri Abhyankar ierr = PetscPrintf(comm," Line search: Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 15772b4ed282SShri Abhyankar } 15782b4ed282SShri Abhyankar break; 15792b4ed282SShri Abhyankar } else { 1580c92abb8aSShri Abhyankar if (vi->lsmonitor) { 15812b4ed282SShri Abhyankar ierr = PetscPrintf(comm," Line search: Cubic step no good, shrinking lambda, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 15822b4ed282SShri Abhyankar } 15832b4ed282SShri Abhyankar } 15842b4ed282SShri Abhyankar count++; 15852b4ed282SShri Abhyankar } 15862b4ed282SShri Abhyankar theend1: 15872b4ed282SShri Abhyankar /* Optional user-defined check for line search step validity */ 15882b4ed282SShri Abhyankar if (vi->postcheckstep && *flag) { 15892b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 15902b4ed282SShri Abhyankar if (changed_y) { 15912b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 15922b4ed282SShri Abhyankar } 15932b4ed282SShri Abhyankar if (changed_y || changed_w) { /* recompute the function if the step has changed */ 15942b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 15952b4ed282SShri Abhyankar if (snes->domainerror) { 15962b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 15972b4ed282SShri Abhyankar PetscFunctionReturn(0); 15982b4ed282SShri Abhyankar } 15992b4ed282SShri Abhyankar ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 16002b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 16012b4ed282SShri Abhyankar ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 16022b4ed282SShri Abhyankar ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 16032b4ed282SShri Abhyankar ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 16042b4ed282SShri Abhyankar } 16052b4ed282SShri Abhyankar } 16062b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 16072b4ed282SShri Abhyankar PetscFunctionReturn(0); 16082b4ed282SShri Abhyankar } 16092b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 16102b4ed282SShri Abhyankar #undef __FUNCT__ 16112b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchQuadratic_VI" 16122b4ed282SShri Abhyankar /* 16132b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchQuadratic in snes/impls/ls/ls.c 16142b4ed282SShri Abhyankar */ 1615ace3abfcSBarry 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) 16162b4ed282SShri Abhyankar { 16172b4ed282SShri Abhyankar /* 16182b4ed282SShri Abhyankar Note that for line search purposes we work with with the related 16192b4ed282SShri Abhyankar minimization problem: 16202b4ed282SShri Abhyankar min z(x): R^n -> R, 16212b4ed282SShri Abhyankar where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 16222b4ed282SShri Abhyankar z(x) is the same as the merit function for the complementarity problem 16232b4ed282SShri Abhyankar */ 16242b4ed282SShri Abhyankar PetscReal initslope,minlambda,lambda,lambdatemp,rellength; 16252b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 16262b4ed282SShri Abhyankar PetscScalar cinitslope; 16272b4ed282SShri Abhyankar #endif 16282b4ed282SShri Abhyankar PetscErrorCode ierr; 16292b4ed282SShri Abhyankar PetscInt count; 16302b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1631ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 16322b4ed282SShri Abhyankar 16332b4ed282SShri Abhyankar PetscFunctionBegin; 16342b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 16352b4ed282SShri Abhyankar *flag = PETSC_TRUE; 16362b4ed282SShri Abhyankar 16372b4ed282SShri Abhyankar ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 16382b4ed282SShri Abhyankar if (*ynorm == 0.0) { 1639c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1640c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr); 16412b4ed282SShri Abhyankar } 16422b4ed282SShri Abhyankar *gnorm = fnorm; 16432b4ed282SShri Abhyankar ierr = VecCopy(x,w);CHKERRQ(ierr); 16442b4ed282SShri Abhyankar ierr = VecCopy(f,g);CHKERRQ(ierr); 16452b4ed282SShri Abhyankar *flag = PETSC_FALSE; 16462b4ed282SShri Abhyankar goto theend2; 16472b4ed282SShri Abhyankar } 16482b4ed282SShri Abhyankar if (*ynorm > vi->maxstep) { /* Step too big, so scale back */ 16492b4ed282SShri Abhyankar ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr); 16502b4ed282SShri Abhyankar *ynorm = vi->maxstep; 16512b4ed282SShri Abhyankar } 16522b4ed282SShri Abhyankar ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 16532b4ed282SShri Abhyankar minlambda = vi->minlambda/rellength; 16542b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 1655*89d9e806SShri Abhyankar ierr = VecDot(vi->dpsi,y,&cinitslope);CHKERRQ(ierr); 16562b4ed282SShri Abhyankar initslope = PetscRealPart(cinitslope); 16572b4ed282SShri Abhyankar #else 1658*89d9e806SShri Abhyankar ierr = VecDot(vi->dpsi,y,&initslope);CHKERRQ(ierr); 16592b4ed282SShri Abhyankar #endif 16602b4ed282SShri Abhyankar if (initslope > 0.0) initslope = -initslope; 16612b4ed282SShri Abhyankar if (initslope == 0.0) initslope = -1.0; 16622b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr); 16632b4ed282SShri Abhyankar 16642b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 16652b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 16662b4ed282SShri Abhyankar ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 16672b4ed282SShri Abhyankar *flag = PETSC_FALSE; 16682b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 16692b4ed282SShri Abhyankar goto theend2; 16702b4ed282SShri Abhyankar } 16712b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 16722b4ed282SShri Abhyankar if (snes->domainerror) { 16732b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 16742b4ed282SShri Abhyankar PetscFunctionReturn(0); 16752b4ed282SShri Abhyankar } 16762b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 16772b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 16782b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */ 1679c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1680c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 16812b4ed282SShri Abhyankar } 16822b4ed282SShri Abhyankar goto theend2; 16832b4ed282SShri Abhyankar } 16842b4ed282SShri Abhyankar 16852b4ed282SShri Abhyankar /* Fit points with quadratic */ 16862b4ed282SShri Abhyankar lambda = 1.0; 16872b4ed282SShri Abhyankar count = 1; 16882b4ed282SShri Abhyankar while (PETSC_TRUE) { 16892b4ed282SShri Abhyankar if (lambda <= minlambda) { /* bad luck; use full step */ 1690c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1691c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr); 1692c92abb8aSShri 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); 16932b4ed282SShri Abhyankar } 16942b4ed282SShri Abhyankar ierr = VecCopy(x,w);CHKERRQ(ierr); 16952b4ed282SShri Abhyankar *flag = PETSC_FALSE; 16962b4ed282SShri Abhyankar break; 16972b4ed282SShri Abhyankar } 16982b4ed282SShri Abhyankar lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 16992b4ed282SShri Abhyankar if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 17002b4ed282SShri Abhyankar if (lambdatemp <= .1*lambda) lambda = .1*lambda; 17012b4ed282SShri Abhyankar else lambda = lambdatemp; 17022b4ed282SShri Abhyankar 17032b4ed282SShri Abhyankar ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 17042b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 17052b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 17062b4ed282SShri 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); 17072b4ed282SShri Abhyankar *flag = PETSC_FALSE; 17082b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 17092b4ed282SShri Abhyankar break; 17102b4ed282SShri Abhyankar } 17112b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 17122b4ed282SShri Abhyankar if (snes->domainerror) { 17132b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 17142b4ed282SShri Abhyankar PetscFunctionReturn(0); 17152b4ed282SShri Abhyankar } 17162b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 17172b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 17182b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */ 1719c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1720c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line Search: Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr); 17212b4ed282SShri Abhyankar } 17222b4ed282SShri Abhyankar break; 17232b4ed282SShri Abhyankar } 17242b4ed282SShri Abhyankar count++; 17252b4ed282SShri Abhyankar } 17262b4ed282SShri Abhyankar theend2: 17272b4ed282SShri Abhyankar /* Optional user-defined check for line search step validity */ 17282b4ed282SShri Abhyankar if (vi->postcheckstep) { 17292b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 17302b4ed282SShri Abhyankar if (changed_y) { 17312b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 17322b4ed282SShri Abhyankar } 17332b4ed282SShri Abhyankar if (changed_y || changed_w) { /* recompute the function if the step has changed */ 17342b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g); 17352b4ed282SShri Abhyankar if (snes->domainerror) { 17362b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 17372b4ed282SShri Abhyankar PetscFunctionReturn(0); 17382b4ed282SShri Abhyankar } 17392b4ed282SShri Abhyankar ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 17402b4ed282SShri Abhyankar ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 17412b4ed282SShri Abhyankar ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 17422b4ed282SShri Abhyankar ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 17432b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 17442b4ed282SShri Abhyankar } 17452b4ed282SShri Abhyankar } 17462b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 17472b4ed282SShri Abhyankar PetscFunctionReturn(0); 17482b4ed282SShri Abhyankar } 17492b4ed282SShri Abhyankar 1750ace3abfcSBarry 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*/ 17512b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 17522b4ed282SShri Abhyankar EXTERN_C_BEGIN 17532b4ed282SShri Abhyankar #undef __FUNCT__ 17542b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchSet_VI" 17552b4ed282SShri Abhyankar PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_VI(SNES snes,FCN2 func,void *lsctx) 17562b4ed282SShri Abhyankar { 17572b4ed282SShri Abhyankar PetscFunctionBegin; 17582b4ed282SShri Abhyankar ((SNES_VI *)(snes->data))->LineSearch = func; 17592b4ed282SShri Abhyankar ((SNES_VI *)(snes->data))->lsP = lsctx; 17602b4ed282SShri Abhyankar PetscFunctionReturn(0); 17612b4ed282SShri Abhyankar } 17622b4ed282SShri Abhyankar EXTERN_C_END 17632b4ed282SShri Abhyankar 17642b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 17652b4ed282SShri Abhyankar EXTERN_C_BEGIN 17662b4ed282SShri Abhyankar #undef __FUNCT__ 17672b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchSetMonitor_VI" 1768ace3abfcSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetMonitor_VI(SNES snes,PetscBool flg) 17692b4ed282SShri Abhyankar { 17702b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 17712b4ed282SShri Abhyankar PetscErrorCode ierr; 17722b4ed282SShri Abhyankar 17732b4ed282SShri Abhyankar PetscFunctionBegin; 1774c92abb8aSShri Abhyankar if (flg && !vi->lsmonitor) { 1775c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",((PetscObject)snes)->tablevel,&vi->lsmonitor);CHKERRQ(ierr); 1776c92abb8aSShri Abhyankar } else if (!flg && vi->lsmonitor) { 1777c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr); 17782b4ed282SShri Abhyankar } 17792b4ed282SShri Abhyankar PetscFunctionReturn(0); 17802b4ed282SShri Abhyankar } 17812b4ed282SShri Abhyankar EXTERN_C_END 17822b4ed282SShri Abhyankar 17832b4ed282SShri Abhyankar /* 17842b4ed282SShri Abhyankar SNESView_VI - Prints info from the SNESVI data structure. 17852b4ed282SShri Abhyankar 17862b4ed282SShri Abhyankar Input Parameters: 17872b4ed282SShri Abhyankar . SNES - the SNES context 17882b4ed282SShri Abhyankar . viewer - visualization context 17892b4ed282SShri Abhyankar 17902b4ed282SShri Abhyankar Application Interface Routine: SNESView() 17912b4ed282SShri Abhyankar */ 17922b4ed282SShri Abhyankar #undef __FUNCT__ 17932b4ed282SShri Abhyankar #define __FUNCT__ "SNESView_VI" 17942b4ed282SShri Abhyankar static PetscErrorCode SNESView_VI(SNES snes,PetscViewer viewer) 17952b4ed282SShri Abhyankar { 17962b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI *)snes->data; 179778c4b9d3SShri Abhyankar const char *cstr,*tstr; 17982b4ed282SShri Abhyankar PetscErrorCode ierr; 1799ace3abfcSBarry Smith PetscBool iascii; 18002b4ed282SShri Abhyankar 18012b4ed282SShri Abhyankar PetscFunctionBegin; 18022b4ed282SShri Abhyankar ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 18032b4ed282SShri Abhyankar if (iascii) { 18042b4ed282SShri Abhyankar if (vi->LineSearch == SNESLineSearchNo_VI) cstr = "SNESLineSearchNo"; 18052b4ed282SShri Abhyankar else if (vi->LineSearch == SNESLineSearchQuadratic_VI) cstr = "SNESLineSearchQuadratic"; 18062b4ed282SShri Abhyankar else if (vi->LineSearch == SNESLineSearchCubic_VI) cstr = "SNESLineSearchCubic"; 18072b4ed282SShri Abhyankar else cstr = "unknown"; 180878c4b9d3SShri Abhyankar if (snes->ops->solve == SNESSolveVI_SS) tstr = "Semismooth"; 180978c4b9d3SShri Abhyankar else if (snes->ops->solve == SNESSolveVI_AS) tstr = "Active Set"; 18100a7c7af0SShri Abhyankar else if (snes->ops->solve == SNESSolveVI_RS) tstr = "Reduced Space"; 181178c4b9d3SShri Abhyankar else tstr = "unknown"; 181278c4b9d3SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," VI algorithm: %s\n",tstr);CHKERRQ(ierr); 18132b4ed282SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 18142b4ed282SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," alpha=%G, maxstep=%G, minlambda=%G\n",vi->alpha,vi->maxstep,vi->minlambda);CHKERRQ(ierr); 18152b4ed282SShri Abhyankar } else { 18162b4ed282SShri Abhyankar SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ VI",((PetscObject)viewer)->type_name); 18172b4ed282SShri Abhyankar } 18182b4ed282SShri Abhyankar PetscFunctionReturn(0); 18192b4ed282SShri Abhyankar } 18202b4ed282SShri Abhyankar 18212b4ed282SShri Abhyankar /* 18222b4ed282SShri Abhyankar SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. 18232b4ed282SShri Abhyankar 18242b4ed282SShri Abhyankar Input Parameters: 18252b4ed282SShri Abhyankar . snes - the SNES context. 18262b4ed282SShri Abhyankar . xl - lower bound. 18272b4ed282SShri Abhyankar . xu - upper bound. 18282b4ed282SShri Abhyankar 18292b4ed282SShri Abhyankar Notes: 18302b4ed282SShri Abhyankar If this routine is not called then the lower and upper bounds are set to 183184c105d7SBarry Smith PETSC_VI_INF and PETSC_VI_NINF respectively during SNESSetUp(). 183284c105d7SBarry Smith 18332b4ed282SShri Abhyankar */ 18342b4ed282SShri Abhyankar 18352b4ed282SShri Abhyankar #undef __FUNCT__ 18362b4ed282SShri Abhyankar #define __FUNCT__ "SNESVISetVariableBounds" 183769c03620SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu) 18382b4ed282SShri Abhyankar { 18392b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 18402b4ed282SShri Abhyankar 18412b4ed282SShri Abhyankar PetscFunctionBegin; 18422b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 18432b4ed282SShri Abhyankar PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 18442b4ed282SShri Abhyankar PetscValidHeaderSpecific(xu,VEC_CLASSID,3); 18452b4ed282SShri Abhyankar 18462b4ed282SShri Abhyankar /* Check if SNESSetFunction is called */ 18472b4ed282SShri Abhyankar if(!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 18482b4ed282SShri Abhyankar 18492b4ed282SShri Abhyankar /* Check if the vector sizes are compatible for lower and upper bounds */ 18502b4ed282SShri 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); 18512b4ed282SShri 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); 18522b4ed282SShri Abhyankar vi->usersetxbounds = PETSC_TRUE; 18532b4ed282SShri Abhyankar vi->xl = xl; 18542b4ed282SShri Abhyankar vi->xu = xu; 18552b4ed282SShri Abhyankar 18562b4ed282SShri Abhyankar PetscFunctionReturn(0); 18572b4ed282SShri Abhyankar } 18582b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 18592b4ed282SShri Abhyankar /* 18602b4ed282SShri Abhyankar SNESSetFromOptions_VI - Sets various parameters for the SNESVI method. 18612b4ed282SShri Abhyankar 18622b4ed282SShri Abhyankar Input Parameter: 18632b4ed282SShri Abhyankar . snes - the SNES context 18642b4ed282SShri Abhyankar 18652b4ed282SShri Abhyankar Application Interface Routine: SNESSetFromOptions() 18662b4ed282SShri Abhyankar */ 18672b4ed282SShri Abhyankar #undef __FUNCT__ 18682b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetFromOptions_VI" 18692b4ed282SShri Abhyankar static PetscErrorCode SNESSetFromOptions_VI(SNES snes) 18702b4ed282SShri Abhyankar { 18712b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI *)snes->data; 18722b4ed282SShri Abhyankar const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 1873d950fb63SShri Abhyankar const char *vies[] = {"ss","as","rs"}; 18742b4ed282SShri Abhyankar PetscErrorCode ierr; 18752b4ed282SShri Abhyankar PetscInt indx; 187669c03620SShri Abhyankar PetscBool flg,set,flg2; 18772b4ed282SShri Abhyankar 18782b4ed282SShri Abhyankar PetscFunctionBegin; 18792b4ed282SShri Abhyankar ierr = PetscOptionsHead("SNES semismooth method options");CHKERRQ(ierr); 18802b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_alpha","Function norm must decrease by","None",vi->alpha,&vi->alpha,0);CHKERRQ(ierr); 18812b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_maxstep","Step must be less than","None",vi->maxstep,&vi->maxstep,0);CHKERRQ(ierr); 18822b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_minlambda","Minimum lambda allowed","None",vi->minlambda,&vi->minlambda,0);CHKERRQ(ierr); 18832b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_const_tol","constraint tolerance","None",vi->const_tol,&vi->const_tol,0);CHKERRQ(ierr); 1884acfcf0e5SJed Brown ierr = PetscOptionsBool("-snes_vi_lsmonitor","Print progress of line searches","SNESLineSearchSetMonitor",vi->lsmonitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 18852b4ed282SShri Abhyankar if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);} 1886d950fb63SShri Abhyankar ierr = PetscOptionsEList("-snes_vi_type","Semismooth algorithm used","",vies,3,"ss",&indx,&flg2);CHKERRQ(ierr); 188769c03620SShri Abhyankar if (flg2) { 188869c03620SShri Abhyankar switch (indx) { 188969c03620SShri Abhyankar case 0: 189069c03620SShri Abhyankar snes->ops->solve = SNESSolveVI_SS; 189169c03620SShri Abhyankar break; 189269c03620SShri Abhyankar case 1: 189369c03620SShri Abhyankar snes->ops->solve = SNESSolveVI_AS; 189469c03620SShri Abhyankar break; 1895d950fb63SShri Abhyankar case 2: 1896d950fb63SShri Abhyankar snes->ops->solve = SNESSolveVI_RS; 1897d950fb63SShri Abhyankar break; 189869c03620SShri Abhyankar } 189969c03620SShri Abhyankar } 1900c92abb8aSShri Abhyankar ierr = PetscOptionsEList("-snes_vi_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 19012b4ed282SShri Abhyankar if (flg) { 19022b4ed282SShri Abhyankar switch (indx) { 19032b4ed282SShri Abhyankar case 0: 19042b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr); 19052b4ed282SShri Abhyankar break; 19062b4ed282SShri Abhyankar case 1: 19072b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr); 19082b4ed282SShri Abhyankar break; 19092b4ed282SShri Abhyankar case 2: 19102b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr); 19112b4ed282SShri Abhyankar break; 19122b4ed282SShri Abhyankar case 3: 19132b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr); 19142b4ed282SShri Abhyankar break; 19152b4ed282SShri Abhyankar } 19162b4ed282SShri Abhyankar } 19172b4ed282SShri Abhyankar ierr = PetscOptionsTail();CHKERRQ(ierr); 19182b4ed282SShri Abhyankar PetscFunctionReturn(0); 19192b4ed282SShri Abhyankar } 19202b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 19212b4ed282SShri Abhyankar /*MC 19222b4ed282SShri Abhyankar SNESVI - Semismooth newton method based nonlinear solver that uses a line search 19232b4ed282SShri Abhyankar 19242b4ed282SShri Abhyankar Options Database: 19252b4ed282SShri Abhyankar + -snes_vi [cubic,quadratic,basic,basicnonorms] - Selects line search 19262b4ed282SShri Abhyankar . -snes_vi_alpha <alpha> - Sets alpha 19272b4ed282SShri 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) 19282b4ed282SShri Abhyankar . -snes_vi_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 19292b4ed282SShri Abhyankar - -snes_vi_monitor - print information about progress of line searches 19302b4ed282SShri Abhyankar 19312b4ed282SShri Abhyankar 19322b4ed282SShri Abhyankar Level: beginner 19332b4ed282SShri Abhyankar 19342b4ed282SShri Abhyankar .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 19352b4ed282SShri Abhyankar SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 19362b4ed282SShri Abhyankar SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams() 19372b4ed282SShri Abhyankar 19382b4ed282SShri Abhyankar M*/ 19392b4ed282SShri Abhyankar EXTERN_C_BEGIN 19402b4ed282SShri Abhyankar #undef __FUNCT__ 19412b4ed282SShri Abhyankar #define __FUNCT__ "SNESCreate_VI" 19422b4ed282SShri Abhyankar PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_VI(SNES snes) 19432b4ed282SShri Abhyankar { 19442b4ed282SShri Abhyankar PetscErrorCode ierr; 19452b4ed282SShri Abhyankar SNES_VI *vi; 19462b4ed282SShri Abhyankar 19472b4ed282SShri Abhyankar PetscFunctionBegin; 19482b4ed282SShri Abhyankar snes->ops->setup = SNESSetUp_VI; 194969c03620SShri Abhyankar snes->ops->solve = SNESSolveVI_SS; 19502b4ed282SShri Abhyankar snes->ops->destroy = SNESDestroy_VI; 19512b4ed282SShri Abhyankar snes->ops->setfromoptions = SNESSetFromOptions_VI; 19522b4ed282SShri Abhyankar snes->ops->view = SNESView_VI; 19532b4ed282SShri Abhyankar snes->ops->converged = SNESDefaultConverged_VI; 19542b4ed282SShri Abhyankar 19552b4ed282SShri Abhyankar ierr = PetscNewLog(snes,SNES_VI,&vi);CHKERRQ(ierr); 19562b4ed282SShri Abhyankar snes->data = (void*)vi; 19572b4ed282SShri Abhyankar vi->alpha = 1.e-4; 19582b4ed282SShri Abhyankar vi->maxstep = 1.e8; 19592b4ed282SShri Abhyankar vi->minlambda = 1.e-12; 19602b4ed282SShri Abhyankar vi->LineSearch = SNESLineSearchCubic_VI; 19612b4ed282SShri Abhyankar vi->lsP = PETSC_NULL; 19622b4ed282SShri Abhyankar vi->postcheckstep = PETSC_NULL; 19632b4ed282SShri Abhyankar vi->postcheck = PETSC_NULL; 19642b4ed282SShri Abhyankar vi->precheckstep = PETSC_NULL; 19652b4ed282SShri Abhyankar vi->precheck = PETSC_NULL; 19662b4ed282SShri Abhyankar vi->const_tol = 2.2204460492503131e-16; 19672b4ed282SShri Abhyankar vi->computessfunction = ComputeFischerFunction; 19682b4ed282SShri Abhyankar 19692b4ed282SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_VI",SNESLineSearchSetMonitor_VI);CHKERRQ(ierr); 19702b4ed282SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_VI",SNESLineSearchSet_VI);CHKERRQ(ierr); 19712b4ed282SShri Abhyankar 19722b4ed282SShri Abhyankar PetscFunctionReturn(0); 19732b4ed282SShri Abhyankar } 19742b4ed282SShri Abhyankar EXTERN_C_END 1975