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 /* 3992b4ed282SShri Abhyankar SNESVIAdjustInitialGuess - Readjusts the initial guess to the SNES solver supplied by the user so that the initial guess lies inside the feasible region . 4002b4ed282SShri Abhyankar 4012b4ed282SShri Abhyankar Input Parameters: 4022b4ed282SShri Abhyankar . lb - lower bound. 4032b4ed282SShri Abhyankar . ub - upper bound. 4042b4ed282SShri Abhyankar 4052b4ed282SShri Abhyankar Output Parameters: 4062b4ed282SShri Abhyankar . X - the readjusted initial guess. 4072b4ed282SShri Abhyankar 4082b4ed282SShri Abhyankar Notes: 4092b4ed282SShri Abhyankar The readjusted initial guess X[i] = max(max(min(l[i],X[i]),min(X[i],u[i])),min(l[i],u[i])) 4102b4ed282SShri Abhyankar */ 4112b4ed282SShri Abhyankar #undef __FUNCT__ 4122b4ed282SShri Abhyankar #define __FUNCT__ "SNESVIAdjustInitialGuess" 4132b4ed282SShri Abhyankar PetscErrorCode SNESVIAdjustInitialGuess(Vec X, Vec lb, Vec ub) 4142b4ed282SShri Abhyankar { 4152b4ed282SShri Abhyankar PetscErrorCode ierr; 4162b4ed282SShri Abhyankar PetscInt i,nlocal; 4172b4ed282SShri Abhyankar PetscScalar *x,*l,*u; 4182b4ed282SShri Abhyankar 4192b4ed282SShri Abhyankar PetscFunctionBegin; 4202b4ed282SShri Abhyankar 4212b4ed282SShri Abhyankar ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 4222b4ed282SShri Abhyankar ierr = VecGetArray(X,&x);CHKERRQ(ierr); 4232b4ed282SShri Abhyankar ierr = VecGetArray(lb,&l);CHKERRQ(ierr); 4242b4ed282SShri Abhyankar ierr = VecGetArray(ub,&u);CHKERRQ(ierr); 4252b4ed282SShri Abhyankar 4262b4ed282SShri Abhyankar for(i = 0; i < nlocal; i++) { 4272b4ed282SShri Abhyankar x[i] = PetscMax(PetscMax(PetscMin(x[i],l[i]),PetscMin(x[i],u[i])),PetscMin(l[i],u[i])); 4282b4ed282SShri Abhyankar } 4292b4ed282SShri Abhyankar 4302b4ed282SShri Abhyankar ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 4312b4ed282SShri Abhyankar ierr = VecRestoreArray(lb,&l);CHKERRQ(ierr); 4322b4ed282SShri Abhyankar ierr = VecRestoreArray(ub,&u);CHKERRQ(ierr); 4332b4ed282SShri Abhyankar 4342b4ed282SShri Abhyankar PetscFunctionReturn(0); 4352b4ed282SShri Abhyankar } 4362b4ed282SShri Abhyankar 4372b4ed282SShri Abhyankar /* -------------------------------------------------------------------- 4382b4ed282SShri Abhyankar 4392b4ed282SShri Abhyankar This file implements a semismooth truncated Newton method with a line search, 4402b4ed282SShri Abhyankar for solving a system of nonlinear equations in complementarity form, using the KSP, Vec, 4412b4ed282SShri Abhyankar and Mat interfaces for linear solvers, vectors, and matrices, 4422b4ed282SShri Abhyankar respectively. 4432b4ed282SShri Abhyankar 4442b4ed282SShri Abhyankar The following basic routines are required for each nonlinear solver: 4452b4ed282SShri Abhyankar SNESCreate_XXX() - Creates a nonlinear solver context 4462b4ed282SShri Abhyankar SNESSetFromOptions_XXX() - Sets runtime options 4472b4ed282SShri Abhyankar SNESSolve_XXX() - Solves the nonlinear system 4482b4ed282SShri Abhyankar SNESDestroy_XXX() - Destroys the nonlinear solver context 4492b4ed282SShri Abhyankar The suffix "_XXX" denotes a particular implementation, in this case 4502b4ed282SShri Abhyankar we use _VI (e.g., SNESCreate_VI, SNESSolve_VI) for solving 4512b4ed282SShri Abhyankar systems of nonlinear equations with a line search (LS) method. 4522b4ed282SShri Abhyankar These routines are actually called via the common user interface 4532b4ed282SShri Abhyankar routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 4542b4ed282SShri Abhyankar SNESDestroy(), so the application code interface remains identical 4552b4ed282SShri Abhyankar for all nonlinear solvers. 4562b4ed282SShri Abhyankar 4572b4ed282SShri Abhyankar Another key routine is: 4582b4ed282SShri Abhyankar SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 4592b4ed282SShri Abhyankar by setting data structures and options. The interface routine SNESSetUp() 4602b4ed282SShri Abhyankar is not usually called directly by the user, but instead is called by 4612b4ed282SShri Abhyankar SNESSolve() if necessary. 4622b4ed282SShri Abhyankar 4632b4ed282SShri Abhyankar Additional basic routines are: 4642b4ed282SShri Abhyankar SNESView_XXX() - Prints details of runtime options that 4652b4ed282SShri Abhyankar have actually been used. 4662b4ed282SShri Abhyankar These are called by application codes via the interface routines 4672b4ed282SShri Abhyankar SNESView(). 4682b4ed282SShri Abhyankar 4692b4ed282SShri Abhyankar The various types of solvers (preconditioners, Krylov subspace methods, 4702b4ed282SShri Abhyankar nonlinear solvers, timesteppers) are all organized similarly, so the 4712b4ed282SShri Abhyankar above description applies to these categories also. 4722b4ed282SShri Abhyankar 4732b4ed282SShri Abhyankar -------------------------------------------------------------------- */ 4742b4ed282SShri Abhyankar /* 47569c03620SShri Abhyankar SNESSolveVI_SS - Solves the complementarity problem with a semismooth Newton 4762b4ed282SShri Abhyankar method using a line search. 4772b4ed282SShri Abhyankar 4782b4ed282SShri Abhyankar Input Parameters: 4792b4ed282SShri Abhyankar . snes - the SNES context 4802b4ed282SShri Abhyankar 4812b4ed282SShri Abhyankar Output Parameter: 4822b4ed282SShri Abhyankar . outits - number of iterations until termination 4832b4ed282SShri Abhyankar 4842b4ed282SShri Abhyankar Application Interface Routine: SNESSolve() 4852b4ed282SShri Abhyankar 4862b4ed282SShri Abhyankar Notes: 4872b4ed282SShri Abhyankar This implements essentially a semismooth Newton method with a 4882b4ed282SShri Abhyankar line search. By default a cubic backtracking line search 4892b4ed282SShri Abhyankar is employed, as described in the text "Numerical Methods for 4902b4ed282SShri Abhyankar Unconstrained Optimization and Nonlinear Equations" by Dennis 4912b4ed282SShri Abhyankar and Schnabel. 4922b4ed282SShri Abhyankar */ 4932b4ed282SShri Abhyankar #undef __FUNCT__ 49469c03620SShri Abhyankar #define __FUNCT__ "SNESSolveVI_SS" 49569c03620SShri Abhyankar PetscErrorCode SNESSolveVI_SS(SNES snes) 4962b4ed282SShri Abhyankar { 4972b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 4982b4ed282SShri Abhyankar PetscErrorCode ierr; 4992b4ed282SShri Abhyankar PetscInt maxits,i,lits; 5003b336b49SShri Abhyankar PetscBool lssucceed; 5012b4ed282SShri Abhyankar MatStructure flg = DIFFERENT_NONZERO_PATTERN; 5022b4ed282SShri Abhyankar PetscReal gnorm,xnorm=0,ynorm; 5032b4ed282SShri Abhyankar Vec Y,X,F,G,W; 5042b4ed282SShri Abhyankar KSPConvergedReason kspreason; 5052b4ed282SShri Abhyankar 5062b4ed282SShri Abhyankar PetscFunctionBegin; 5072b4ed282SShri Abhyankar snes->numFailures = 0; 5082b4ed282SShri Abhyankar snes->numLinearSolveFailures = 0; 5092b4ed282SShri Abhyankar snes->reason = SNES_CONVERGED_ITERATING; 5102b4ed282SShri Abhyankar 5112b4ed282SShri Abhyankar maxits = snes->max_its; /* maximum number of iterations */ 5122b4ed282SShri Abhyankar X = snes->vec_sol; /* solution vector */ 5132b4ed282SShri Abhyankar F = snes->vec_func; /* residual vector */ 5142b4ed282SShri Abhyankar Y = snes->work[0]; /* work vectors */ 5152b4ed282SShri Abhyankar G = snes->work[1]; 5162b4ed282SShri Abhyankar W = snes->work[2]; 5172b4ed282SShri Abhyankar 5182b4ed282SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 5192b4ed282SShri Abhyankar snes->iter = 0; 5202b4ed282SShri Abhyankar snes->norm = 0.0; 5212b4ed282SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 5222b4ed282SShri Abhyankar 5232b4ed282SShri Abhyankar ierr = SNESVIAdjustInitialGuess(X,vi->xl,vi->xu);CHKERRQ(ierr); 5242b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr); 5252b4ed282SShri Abhyankar if (snes->domainerror) { 5262b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 5272b4ed282SShri Abhyankar PetscFunctionReturn(0); 5282b4ed282SShri Abhyankar } 5292b4ed282SShri Abhyankar /* Compute Merit function */ 5302b4ed282SShri Abhyankar ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr); 5312b4ed282SShri Abhyankar 5322b4ed282SShri Abhyankar ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 5332b4ed282SShri Abhyankar ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 5342b4ed282SShri Abhyankar if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 5352b4ed282SShri Abhyankar 5362b4ed282SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 5372b4ed282SShri Abhyankar snes->norm = vi->phinorm; 5382b4ed282SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 5392b4ed282SShri Abhyankar SNESLogConvHistory(snes,vi->phinorm,0); 5402b4ed282SShri Abhyankar SNESMonitor(snes,0,vi->phinorm); 5412b4ed282SShri Abhyankar 5422b4ed282SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 5432b4ed282SShri Abhyankar snes->ttol = vi->phinorm*snes->rtol; 5442b4ed282SShri Abhyankar /* test convergence */ 5452b4ed282SShri Abhyankar ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 5462b4ed282SShri Abhyankar if (snes->reason) PetscFunctionReturn(0); 5472b4ed282SShri Abhyankar 5482b4ed282SShri Abhyankar for (i=0; i<maxits; i++) { 5492b4ed282SShri Abhyankar 5502b4ed282SShri Abhyankar /* Call general purpose update function */ 5512b4ed282SShri Abhyankar if (snes->ops->update) { 5522b4ed282SShri Abhyankar ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 5532b4ed282SShri Abhyankar } 5542b4ed282SShri Abhyankar 5552b4ed282SShri Abhyankar /* Solve J Y = Phi, where J is the semismooth jacobian */ 556a79edbefSShri Abhyankar /* Get the nonlinear function jacobian */ 5572b4ed282SShri Abhyankar ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 558a79edbefSShri Abhyankar /* Get the diagonal shift and row scaling vectors */ 559a79edbefSShri Abhyankar ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr); 560a79edbefSShri Abhyankar /* Compute the semismooth jacobian */ 561a79edbefSShri Abhyankar ierr = SNESVIComputeJacobian(snes->jacobian,snes->jacobian_pre,vi->Da,vi->Db);CHKERRQ(ierr); 5622b4ed282SShri Abhyankar 5632b4ed282SShri Abhyankar ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 5642b4ed282SShri Abhyankar ierr = SNES_KSPSolve(snes,snes->ksp,vi->phi,Y);CHKERRQ(ierr); 5652b4ed282SShri Abhyankar ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 5663b336b49SShri Abhyankar 5673b336b49SShri Abhyankar if (kspreason < 0) { 5682b4ed282SShri Abhyankar if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 5692b4ed282SShri 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); 5702b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 5712b4ed282SShri Abhyankar break; 5722b4ed282SShri Abhyankar } 5732b4ed282SShri Abhyankar } 5742b4ed282SShri Abhyankar ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 5752b4ed282SShri Abhyankar snes->linear_its += lits; 5762b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 5772b4ed282SShri Abhyankar /* 5782b4ed282SShri Abhyankar if (vi->precheckstep) { 579ace3abfcSBarry Smith PetscBool changed_y = PETSC_FALSE; 5802b4ed282SShri Abhyankar ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr); 5812b4ed282SShri Abhyankar } 5822b4ed282SShri Abhyankar 5832b4ed282SShri Abhyankar if (PetscLogPrintInfo){ 5842b4ed282SShri Abhyankar ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 5852b4ed282SShri Abhyankar } 5862b4ed282SShri Abhyankar */ 5872b4ed282SShri Abhyankar /* Compute a (scaled) negative update in the line search routine: 5882b4ed282SShri Abhyankar Y <- X - lambda*Y 5892b4ed282SShri Abhyankar and evaluate G = function(Y) (depends on the line search). 5902b4ed282SShri Abhyankar */ 5912b4ed282SShri Abhyankar ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 5922b4ed282SShri Abhyankar ynorm = 1; gnorm = vi->phinorm; 5932b4ed282SShri Abhyankar ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 5942b4ed282SShri Abhyankar ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 5952b4ed282SShri Abhyankar if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 5962b4ed282SShri Abhyankar if (snes->domainerror) { 5972b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 5982b4ed282SShri Abhyankar PetscFunctionReturn(0); 5992b4ed282SShri Abhyankar } 6002b4ed282SShri Abhyankar if (!lssucceed) { 6012b4ed282SShri Abhyankar if (++snes->numFailures >= snes->maxFailures) { 602ace3abfcSBarry Smith PetscBool ismin; 6032b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_LINE_SEARCH; 6042b4ed282SShri Abhyankar ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 6052b4ed282SShri Abhyankar if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 6062b4ed282SShri Abhyankar break; 6072b4ed282SShri Abhyankar } 6082b4ed282SShri Abhyankar } 6092b4ed282SShri Abhyankar /* Update function and solution vectors */ 6102b4ed282SShri Abhyankar vi->phinorm = gnorm; 6112b4ed282SShri Abhyankar vi->merit = 0.5*vi->phinorm*vi->phinorm; 6122b4ed282SShri Abhyankar ierr = VecCopy(G,vi->phi);CHKERRQ(ierr); 6132b4ed282SShri Abhyankar ierr = VecCopy(W,X);CHKERRQ(ierr); 6142b4ed282SShri Abhyankar /* Monitor convergence */ 6152b4ed282SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 6162b4ed282SShri Abhyankar snes->iter = i+1; 6172b4ed282SShri Abhyankar snes->norm = vi->phinorm; 6182b4ed282SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 6192b4ed282SShri Abhyankar SNESLogConvHistory(snes,snes->norm,lits); 6202b4ed282SShri Abhyankar SNESMonitor(snes,snes->iter,snes->norm); 6212b4ed282SShri Abhyankar /* Test for convergence, xnorm = || X || */ 6222b4ed282SShri Abhyankar if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 6232b4ed282SShri Abhyankar ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 6242b4ed282SShri Abhyankar if (snes->reason) break; 6252b4ed282SShri Abhyankar } 6262b4ed282SShri Abhyankar if (i == maxits) { 6272b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 6282b4ed282SShri Abhyankar if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 6292b4ed282SShri Abhyankar } 6302b4ed282SShri Abhyankar PetscFunctionReturn(0); 6312b4ed282SShri Abhyankar } 63269c03620SShri Abhyankar 63369c03620SShri Abhyankar #undef __FUNCT__ 63469c03620SShri Abhyankar #define __FUNCT__ "SNESVICreateIndexSets_AS" 635a79edbefSShri Abhyankar PetscErrorCode SNESVICreateIndexSets_AS(SNES snes,Vec Db,PetscReal thresh,IS* ISact,IS* ISinact) 63669c03620SShri Abhyankar { 63769c03620SShri Abhyankar PetscErrorCode ierr; 63869c03620SShri Abhyankar PetscInt i,nlocal,ilow,ihigh,nloc_isact=0,nloc_isinact=0; 63969c03620SShri Abhyankar PetscInt *idx_act,*idx_inact,i1=0,i2=0; 64069c03620SShri Abhyankar PetscScalar *db; 64169c03620SShri Abhyankar 64269c03620SShri Abhyankar PetscFunctionBegin; 64369c03620SShri Abhyankar 64469c03620SShri Abhyankar ierr = VecGetLocalSize(Db,&nlocal);CHKERRQ(ierr); 64569c03620SShri Abhyankar ierr = VecGetOwnershipRange(Db,&ilow,&ihigh);CHKERRQ(ierr); 64669c03620SShri Abhyankar ierr = VecGetArray(Db,&db);CHKERRQ(ierr); 64769c03620SShri Abhyankar /* Compute the sizes of the active and inactive sets */ 64869c03620SShri Abhyankar for (i=0; i < nlocal;i++) { 649fca35906SShri Abhyankar if (PetscAbsScalar(db[i]) <= thresh) nloc_isact++; 650fca35906SShri Abhyankar else nloc_isinact++; 65169c03620SShri Abhyankar } 65269c03620SShri Abhyankar ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr); 65369c03620SShri Abhyankar ierr = PetscMalloc(nloc_isinact*sizeof(PetscInt),&idx_inact);CHKERRQ(ierr); 65469c03620SShri Abhyankar 65569c03620SShri Abhyankar /* Creating the indexing arrays */ 656fca35906SShri Abhyankar for(i=0; i < nlocal; i++) { 657fca35906SShri Abhyankar if (PetscAbsScalar(db[i]) <= thresh) idx_act[i1++] = ilow+i; 658fca35906SShri Abhyankar else idx_inact[i2++] = ilow+i; 65969c03620SShri Abhyankar } 66069c03620SShri Abhyankar 66169c03620SShri Abhyankar /* Create the index sets */ 66269c03620SShri Abhyankar ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isact,idx_act,PETSC_COPY_VALUES,ISact);CHKERRQ(ierr); 66369c03620SShri Abhyankar ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isinact,idx_inact,PETSC_COPY_VALUES,ISinact);CHKERRQ(ierr); 66469c03620SShri Abhyankar 66569c03620SShri Abhyankar ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr); 66669c03620SShri Abhyankar ierr = PetscFree(idx_act);CHKERRQ(ierr); 66769c03620SShri Abhyankar ierr = PetscFree(idx_inact);CHKERRQ(ierr); 66869c03620SShri Abhyankar PetscFunctionReturn(0); 66969c03620SShri Abhyankar } 67069c03620SShri Abhyankar 671*d950fb63SShri Abhyankar #undef __FUNCT__ 672*d950fb63SShri Abhyankar #define __FUNCT__ "SNESVICreateIndexSets_RS" 673*d950fb63SShri Abhyankar PetscErrorCode SNESVICreateIndexSets_RS(SNES snes,Vec X,Vec Xl, Vec Xu,IS* ISact,IS* ISinact) 674*d950fb63SShri Abhyankar { 675*d950fb63SShri Abhyankar PetscErrorCode ierr; 676*d950fb63SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 677*d950fb63SShri Abhyankar PetscInt i,nlocal,ilow,ihigh,nloc_isact=0,nloc_isinact=0; 678*d950fb63SShri Abhyankar PetscInt *idx_act,*idx_inact,i1=0,i2=0; 679*d950fb63SShri Abhyankar PetscScalar *x,*l,*u,*f; 680*d950fb63SShri Abhyankar Vec F = snes->vec_func; 681*d950fb63SShri Abhyankar 682*d950fb63SShri Abhyankar PetscFunctionBegin; 683*d950fb63SShri Abhyankar 684*d950fb63SShri Abhyankar ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 685*d950fb63SShri Abhyankar ierr = VecGetOwnershipRange(X,&ilow,&ihigh);CHKERRQ(ierr); 686*d950fb63SShri Abhyankar ierr = VecGetArray(X,&x);CHKERRQ(ierr); 687*d950fb63SShri Abhyankar ierr = VecGetArray(Xl,&l);CHKERRQ(ierr); 688*d950fb63SShri Abhyankar ierr = VecGetArray(Xu,&u);CHKERRQ(ierr); 689*d950fb63SShri Abhyankar ierr = VecGetArray(F,&f);CHKERRQ(ierr); 690*d950fb63SShri Abhyankar /* Compute the sizes of the active and inactive sets */ 691*d950fb63SShri Abhyankar for (i=0; i < nlocal;i++) { 692*d950fb63SShri Abhyankar if (PetscAbsScalar(x[i] - l[i]) <= vi->const_tol && f[i] > vi->const_tol) nloc_isact++; 693*d950fb63SShri Abhyankar else if (PetscAbsScalar(u[i] - x[i]) <= vi->const_tol && f[i] < vi->const_tol) nloc_isact++; 694*d950fb63SShri Abhyankar else nloc_isinact++; 695*d950fb63SShri Abhyankar } 696*d950fb63SShri Abhyankar ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr); 697*d950fb63SShri Abhyankar ierr = PetscMalloc(nloc_isinact*sizeof(PetscInt),&idx_inact);CHKERRQ(ierr); 698*d950fb63SShri Abhyankar 699*d950fb63SShri Abhyankar /* Creating the indexing arrays */ 700*d950fb63SShri Abhyankar for(i=0; i < nlocal; i++) { 701*d950fb63SShri Abhyankar if (PetscAbsScalar(x[i] - l[i]) <= vi->const_tol && f[i] > vi->const_tol) idx_act[i1++] = ilow+i; 702*d950fb63SShri Abhyankar else if (PetscAbsScalar(u[i] - x[i]) <= vi->const_tol && f[i] < vi->const_tol) idx_act[i1++] = ilow+i; 703*d950fb63SShri Abhyankar else idx_inact[i2++] = ilow+i; 704*d950fb63SShri Abhyankar } 705*d950fb63SShri Abhyankar 706*d950fb63SShri Abhyankar /* Create the index sets */ 707*d950fb63SShri Abhyankar ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isact,idx_act,PETSC_COPY_VALUES,ISact);CHKERRQ(ierr); 708*d950fb63SShri Abhyankar ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isinact,idx_inact,PETSC_COPY_VALUES,ISinact);CHKERRQ(ierr); 709*d950fb63SShri Abhyankar 710*d950fb63SShri Abhyankar ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 711*d950fb63SShri Abhyankar ierr = VecRestoreArray(Xl,&l);CHKERRQ(ierr); 712*d950fb63SShri Abhyankar ierr = VecRestoreArray(Xu,&u);CHKERRQ(ierr); 713*d950fb63SShri Abhyankar ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 714*d950fb63SShri Abhyankar ierr = PetscFree(idx_act);CHKERRQ(ierr); 715*d950fb63SShri Abhyankar ierr = PetscFree(idx_inact);CHKERRQ(ierr); 716*d950fb63SShri Abhyankar PetscFunctionReturn(0); 717*d950fb63SShri Abhyankar } 718*d950fb63SShri Abhyankar 719dbd914b8SShri Abhyankar /* Create active and inactive set vectors. The local size of this vector is set and petsc computes the global size */ 720dbd914b8SShri Abhyankar #undef __FUNCT__ 721dbd914b8SShri Abhyankar #define __FUNCT__ "SNESVICreateVectors_AS" 722dbd914b8SShri Abhyankar PetscErrorCode SNESVICreateVectors_AS(SNES snes,PetscInt n,Vec* newv) 723dbd914b8SShri Abhyankar { 724dbd914b8SShri Abhyankar PetscErrorCode ierr; 725dbd914b8SShri Abhyankar Vec v; 726dbd914b8SShri Abhyankar 727dbd914b8SShri Abhyankar PetscFunctionBegin; 728dbd914b8SShri Abhyankar ierr = VecCreate(((PetscObject)snes)->comm,&v);CHKERRQ(ierr); 729dbd914b8SShri Abhyankar ierr = VecSetSizes(v,n,PETSC_DECIDE);CHKERRQ(ierr); 730dbd914b8SShri Abhyankar ierr = VecSetFromOptions(v);CHKERRQ(ierr); 731dbd914b8SShri Abhyankar *newv = v; 732dbd914b8SShri Abhyankar 733dbd914b8SShri Abhyankar PetscFunctionReturn(0); 734dbd914b8SShri Abhyankar } 735dbd914b8SShri Abhyankar 736dbd914b8SShri Abhyankar 73769c03620SShri Abhyankar /* Variational Inequality solver using active set method */ 73869c03620SShri Abhyankar #undef __FUNCT__ 73969c03620SShri Abhyankar #define __FUNCT__ "SNESSolveVI_AS" 74069c03620SShri Abhyankar PetscErrorCode SNESSolveVI_AS(SNES snes) 74169c03620SShri Abhyankar { 74269c03620SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 74369c03620SShri Abhyankar PetscErrorCode ierr; 744408da460SShri Abhyankar PetscInt maxits,i,lits,Nis_act=0; 7453b336b49SShri Abhyankar PetscBool lssucceed; 74669c03620SShri Abhyankar MatStructure flg = DIFFERENT_NONZERO_PATTERN; 74769c03620SShri Abhyankar PetscReal gnorm,xnorm=0,ynorm; 74869c03620SShri Abhyankar Vec Y,X,F,G,W; 74969c03620SShri Abhyankar KSPConvergedReason kspreason; 75069c03620SShri Abhyankar 75169c03620SShri Abhyankar PetscFunctionBegin; 75269c03620SShri Abhyankar snes->numFailures = 0; 75369c03620SShri Abhyankar snes->numLinearSolveFailures = 0; 75469c03620SShri Abhyankar snes->reason = SNES_CONVERGED_ITERATING; 75569c03620SShri Abhyankar 75669c03620SShri Abhyankar maxits = snes->max_its; /* maximum number of iterations */ 75769c03620SShri Abhyankar X = snes->vec_sol; /* solution vector */ 75869c03620SShri Abhyankar F = snes->vec_func; /* residual vector */ 75969c03620SShri Abhyankar Y = snes->work[0]; /* work vectors */ 76069c03620SShri Abhyankar G = snes->work[1]; 76169c03620SShri Abhyankar W = snes->work[2]; 76269c03620SShri Abhyankar 76369c03620SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 76469c03620SShri Abhyankar snes->iter = 0; 76569c03620SShri Abhyankar snes->norm = 0.0; 76669c03620SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 76769c03620SShri Abhyankar 76869c03620SShri Abhyankar ierr = SNESVIAdjustInitialGuess(X,vi->xl,vi->xu);CHKERRQ(ierr); 76969c03620SShri Abhyankar ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr); 77069c03620SShri Abhyankar if (snes->domainerror) { 77169c03620SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 77269c03620SShri Abhyankar PetscFunctionReturn(0); 77369c03620SShri Abhyankar } 77469c03620SShri Abhyankar /* Compute Merit function */ 77569c03620SShri Abhyankar ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr); 77669c03620SShri Abhyankar 77769c03620SShri Abhyankar ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 77869c03620SShri Abhyankar ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 77969c03620SShri Abhyankar if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 78069c03620SShri Abhyankar 78169c03620SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 78269c03620SShri Abhyankar snes->norm = vi->phinorm; 78369c03620SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 78469c03620SShri Abhyankar SNESLogConvHistory(snes,vi->phinorm,0); 78569c03620SShri Abhyankar SNESMonitor(snes,0,vi->phinorm); 78669c03620SShri Abhyankar 78769c03620SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 78869c03620SShri Abhyankar snes->ttol = vi->phinorm*snes->rtol; 78969c03620SShri Abhyankar /* test convergence */ 79069c03620SShri Abhyankar ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 79169c03620SShri Abhyankar if (snes->reason) PetscFunctionReturn(0); 79269c03620SShri Abhyankar 79369c03620SShri Abhyankar for (i=0; i<maxits; i++) { 79469c03620SShri Abhyankar 795a79edbefSShri Abhyankar IS IS_act,IS_inact; /* _act -> active set _inact -> inactive set */ 796a79edbefSShri Abhyankar PetscReal thresh,J_norm1; 797dbd914b8SShri Abhyankar VecScatter scat_act,scat_inact; 798408da460SShri Abhyankar PetscInt nis_act,nis_inact,Nis_act_prev; 799dbd914b8SShri Abhyankar Vec Da_act,Da_inact,Db_inact; 800dbd914b8SShri Abhyankar Vec Y_act,Y_inact,phi_act,phi_inact; 801d76c674bSShri Abhyankar Mat jac_inact_inact,jac_inact_act,prejac_inact_inact; 802a79edbefSShri Abhyankar 80369c03620SShri Abhyankar /* Call general purpose update function */ 80469c03620SShri Abhyankar if (snes->ops->update) { 80569c03620SShri Abhyankar ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 80669c03620SShri Abhyankar } 80769c03620SShri Abhyankar ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 80869c03620SShri Abhyankar /* Compute the threshold value for creating active and inactive sets */ 80969c03620SShri Abhyankar ierr = MatNorm(snes->jacobian,NORM_1,&J_norm1);CHKERRQ(ierr); 81069c03620SShri Abhyankar thresh = PetscMin(vi->merit,1e-2)/(1+J_norm1); 811a79edbefSShri Abhyankar 812a79edbefSShri Abhyankar /* Compute B-subdifferential vectors Da and Db */ 813a79edbefSShri Abhyankar ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr); 814dbd914b8SShri Abhyankar 81569c03620SShri Abhyankar /* Create active and inactive index sets */ 81669c03620SShri Abhyankar ierr = SNESVICreateIndexSets_AS(snes,vi->Db,thresh,&IS_act,&IS_inact);CHKERRQ(ierr); 81769c03620SShri Abhyankar 818408da460SShri Abhyankar Nis_act_prev = Nis_act; 819408da460SShri Abhyankar /* Get sizes of active and inactive sets */ 820dbd914b8SShri Abhyankar ierr = ISGetLocalSize(IS_act,&nis_act);CHKERRQ(ierr); 821dbd914b8SShri Abhyankar ierr = ISGetLocalSize(IS_inact,&nis_inact);CHKERRQ(ierr); 822408da460SShri Abhyankar ierr = ISGetSize(IS_act,&Nis_act);CHKERRQ(ierr); 823dbd914b8SShri Abhyankar 824d76c674bSShri Abhyankar ierr = PetscPrintf(PETSC_COMM_WORLD,"Size of active set = %d, size of inactive set = %d\n",nis_act,nis_inact);CHKERRQ(ierr); 825d76c674bSShri Abhyankar 826dbd914b8SShri Abhyankar /* Create active and inactive set vectors */ 827dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_act,&Da_act);CHKERRQ(ierr); 828dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_inact,&Da_inact);CHKERRQ(ierr); 829dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_inact,&Db_inact);CHKERRQ(ierr); 830dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_act,&phi_act);CHKERRQ(ierr); 831dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_inact,&phi_inact);CHKERRQ(ierr); 832dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_act,&Y_act);CHKERRQ(ierr); 833dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_inact,&Y_inact);CHKERRQ(ierr); 834dbd914b8SShri Abhyankar 835dbd914b8SShri Abhyankar /* Create inactive set submatrices */ 836dbd914b8SShri Abhyankar ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_act,MAT_INITIAL_MATRIX,&jac_inact_act);CHKERRQ(ierr); 837dbd914b8SShri Abhyankar ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);CHKERRQ(ierr); 838dbd914b8SShri Abhyankar 839dbd914b8SShri Abhyankar /* Create scatter contexts */ 840dbd914b8SShri Abhyankar ierr = VecScatterCreate(vi->Da,IS_act,Da_act,PETSC_NULL,&scat_act);CHKERRQ(ierr); 841dbd914b8SShri Abhyankar ierr = VecScatterCreate(vi->Da,IS_inact,Da_inact,PETSC_NULL,&scat_inact);CHKERRQ(ierr); 842dbd914b8SShri Abhyankar 843dbd914b8SShri Abhyankar /* Do a vec scatter to active and inactive set vectors */ 844dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_act,vi->Da,Da_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 845dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_act,vi->Da,Da_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 846dbd914b8SShri Abhyankar 847dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_inact,vi->Da,Da_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 848dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_inact,vi->Da,Da_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 849dbd914b8SShri Abhyankar 850dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_inact,vi->Db,Db_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 851dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_inact,vi->Db,Db_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 852dbd914b8SShri Abhyankar 853dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_act,vi->phi,phi_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 854dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_act,vi->phi,phi_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 855dbd914b8SShri Abhyankar 856dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_inact,vi->phi,phi_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 857dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_inact,vi->phi,phi_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 858dbd914b8SShri Abhyankar 859dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 860dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 861dbd914b8SShri Abhyankar 862dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 863dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 864dbd914b8SShri Abhyankar 865dbd914b8SShri Abhyankar /* Active set direction */ 866dbd914b8SShri Abhyankar ierr = VecPointwiseDivide(Y_act,phi_act,Da_act);CHKERRQ(ierr); 867d76c674bSShri Abhyankar /* inactive set jacobian and preconditioner */ 868dbd914b8SShri Abhyankar ierr = VecPointwiseDivide(Da_inact,Da_inact,Db_inact);CHKERRQ(ierr); 869dbd914b8SShri Abhyankar ierr = MatDiagonalSet(jac_inact_inact,Da_inact,ADD_VALUES);CHKERRQ(ierr); 870d76c674bSShri Abhyankar if (snes->jacobian != snes->jacobian_pre) { 871d76c674bSShri Abhyankar ierr = MatGetSubMatrix(snes->jacobian_pre,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&prejac_inact_inact);CHKERRQ(ierr); 872d76c674bSShri Abhyankar ierr = MatDiagonalSet(prejac_inact_inact,Da_inact,ADD_VALUES);CHKERRQ(ierr); 873d76c674bSShri Abhyankar } else prejac_inact_inact = jac_inact_inact; 874d76c674bSShri Abhyankar 875dbd914b8SShri Abhyankar /* right hand side */ 876dbd914b8SShri Abhyankar ierr = VecPointwiseDivide(phi_inact,phi_inact,Db_inact);CHKERRQ(ierr); 877dbd914b8SShri Abhyankar ierr = MatMult(jac_inact_act,Y_act,Db_inact);CHKERRQ(ierr); 878dbd914b8SShri Abhyankar ierr = VecAXPY(phi_inact,-1.0,Db_inact);CHKERRQ(ierr); 879dbd914b8SShri Abhyankar 880408da460SShri Abhyankar if ((i != 0) && (Nis_act != Nis_act_prev)) { 881408da460SShri Abhyankar KSP kspnew,snesksp; 882408da460SShri Abhyankar PC pcnew; 883408da460SShri Abhyankar /* The active and inactive set sizes have changed so need to create a new snes->ksp object */ 884408da460SShri Abhyankar ierr = SNESGetKSP(snes,&snesksp);CHKERRQ(ierr); 885408da460SShri Abhyankar ierr = KSPCreate(((PetscObject)snes)->comm,&kspnew);CHKERRQ(ierr); 886408da460SShri Abhyankar /* Copy over snes->ksp info */ 887408da460SShri Abhyankar kspnew->pc_side = snesksp->pc_side; 888408da460SShri Abhyankar kspnew->rtol = snesksp->rtol; 889408da460SShri Abhyankar kspnew->abstol = snesksp->abstol; 890408da460SShri Abhyankar kspnew->max_it = snesksp->max_it; 891408da460SShri Abhyankar ierr = KSPSetType(kspnew,((PetscObject)snesksp)->type_name);CHKERRQ(ierr); 892408da460SShri Abhyankar ierr = KSPGetPC(kspnew,&pcnew);CHKERRQ(ierr); 893408da460SShri Abhyankar ierr = PCSetType(kspnew->pc,((PetscObject)snesksp->pc)->type_name);CHKERRQ(ierr); 894408da460SShri Abhyankar ierr = KSPDestroy(snesksp);CHKERRQ(ierr); 895408da460SShri Abhyankar snes->ksp = kspnew; 896408da460SShri Abhyankar ierr = PetscLogObjectParent(snes,kspnew);CHKERRQ(ierr); 897408da460SShri Abhyankar ierr = KSPSetFromOptions(kspnew);CHKERRQ(ierr); 898408da460SShri Abhyankar } 899408da460SShri Abhyankar 900d76c674bSShri Abhyankar ierr = KSPSetOperators(snes->ksp,jac_inact_inact,prejac_inact_inact,flg);CHKERRQ(ierr); 901dbd914b8SShri Abhyankar ierr = SNES_KSPSolve(snes,snes->ksp,phi_inact,Y_inact);CHKERRQ(ierr); 902dbd914b8SShri Abhyankar ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 9033b336b49SShri Abhyankar if (kspreason < 0) { 9043b336b49SShri Abhyankar if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 9053b336b49SShri 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); 9063b336b49SShri Abhyankar snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 9073b336b49SShri Abhyankar break; 9083b336b49SShri Abhyankar } 9093b336b49SShri Abhyankar } 910dbd914b8SShri Abhyankar 911dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 912dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 913dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 914dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 915dbd914b8SShri Abhyankar 916dbd914b8SShri Abhyankar ierr = VecDestroy(Da_act);CHKERRQ(ierr); 917dbd914b8SShri Abhyankar ierr = VecDestroy(Da_inact);CHKERRQ(ierr); 918dbd914b8SShri Abhyankar ierr = VecDestroy(Db_inact);CHKERRQ(ierr); 919dbd914b8SShri Abhyankar ierr = VecDestroy(phi_act);CHKERRQ(ierr); 920dbd914b8SShri Abhyankar ierr = VecDestroy(phi_inact);CHKERRQ(ierr); 921dbd914b8SShri Abhyankar ierr = VecDestroy(Y_act);CHKERRQ(ierr); 922dbd914b8SShri Abhyankar ierr = VecDestroy(Y_inact);CHKERRQ(ierr); 923dbd914b8SShri Abhyankar ierr = VecScatterDestroy(scat_act);CHKERRQ(ierr); 924dbd914b8SShri Abhyankar ierr = VecScatterDestroy(scat_inact);CHKERRQ(ierr); 925730c24dcSShri Abhyankar ierr = ISDestroy(IS_act);CHKERRQ(ierr); 926730c24dcSShri Abhyankar ierr = ISDestroy(IS_inact);CHKERRQ(ierr); 927bac2f86eSShri Abhyankar ierr = MatDestroy(jac_inact_act);CHKERRQ(ierr); 928bac2f86eSShri Abhyankar ierr = MatDestroy(jac_inact_inact);CHKERRQ(ierr); 929d76c674bSShri Abhyankar if (snes->jacobian != snes->jacobian_pre) { 930d76c674bSShri Abhyankar ierr = MatDestroy(prejac_inact_inact);CHKERRQ(ierr); 931d76c674bSShri Abhyankar } 932730c24dcSShri Abhyankar 93369c03620SShri Abhyankar ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 93469c03620SShri Abhyankar snes->linear_its += lits; 93569c03620SShri Abhyankar ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 93669c03620SShri Abhyankar /* 93769c03620SShri Abhyankar if (vi->precheckstep) { 93869c03620SShri Abhyankar PetscBool changed_y = PETSC_FALSE; 93969c03620SShri Abhyankar ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr); 94069c03620SShri Abhyankar } 94169c03620SShri Abhyankar 94269c03620SShri Abhyankar if (PetscLogPrintInfo){ 94369c03620SShri Abhyankar ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 94469c03620SShri Abhyankar } 94569c03620SShri Abhyankar */ 94669c03620SShri Abhyankar /* Compute a (scaled) negative update in the line search routine: 94769c03620SShri Abhyankar Y <- X - lambda*Y 94869c03620SShri Abhyankar and evaluate G = function(Y) (depends on the line search). 94969c03620SShri Abhyankar */ 95069c03620SShri Abhyankar ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 95169c03620SShri Abhyankar ynorm = 1; gnorm = vi->phinorm; 9523b336b49SShri Abhyankar ierr = SNESVIComputeJacobian(snes->jacobian,snes->jacobian_pre,vi->Da,vi->Db);CHKERRQ(ierr); 95369c03620SShri Abhyankar ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 95469c03620SShri Abhyankar ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 95569c03620SShri Abhyankar if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 95669c03620SShri Abhyankar if (snes->domainerror) { 95769c03620SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 95869c03620SShri Abhyankar PetscFunctionReturn(0); 95969c03620SShri Abhyankar } 96069c03620SShri Abhyankar if (!lssucceed) { 96169c03620SShri Abhyankar if (++snes->numFailures >= snes->maxFailures) { 96269c03620SShri Abhyankar PetscBool ismin; 96369c03620SShri Abhyankar snes->reason = SNES_DIVERGED_LINE_SEARCH; 96469c03620SShri Abhyankar ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 96569c03620SShri Abhyankar if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 96669c03620SShri Abhyankar break; 96769c03620SShri Abhyankar } 96869c03620SShri Abhyankar } 96969c03620SShri Abhyankar /* Update function and solution vectors */ 97069c03620SShri Abhyankar vi->phinorm = gnorm; 97169c03620SShri Abhyankar vi->merit = 0.5*vi->phinorm*vi->phinorm; 97269c03620SShri Abhyankar ierr = VecCopy(G,vi->phi);CHKERRQ(ierr); 97369c03620SShri Abhyankar ierr = VecCopy(W,X);CHKERRQ(ierr); 97469c03620SShri Abhyankar /* Monitor convergence */ 97569c03620SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 97669c03620SShri Abhyankar snes->iter = i+1; 97769c03620SShri Abhyankar snes->norm = vi->phinorm; 97869c03620SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 97969c03620SShri Abhyankar SNESLogConvHistory(snes,snes->norm,lits); 98069c03620SShri Abhyankar SNESMonitor(snes,snes->iter,snes->norm); 98169c03620SShri Abhyankar /* Test for convergence, xnorm = || X || */ 98269c03620SShri Abhyankar if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 98369c03620SShri Abhyankar ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 98469c03620SShri Abhyankar if (snes->reason) break; 98569c03620SShri Abhyankar } 98669c03620SShri Abhyankar if (i == maxits) { 98769c03620SShri Abhyankar ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 98869c03620SShri Abhyankar if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 98969c03620SShri Abhyankar } 99069c03620SShri Abhyankar PetscFunctionReturn(0); 99169c03620SShri Abhyankar } 99269c03620SShri Abhyankar 993*d950fb63SShri Abhyankar /* Variational Inequality solver using active set method */ 994*d950fb63SShri Abhyankar #undef __FUNCT__ 995*d950fb63SShri Abhyankar #define __FUNCT__ "SNESSolveVI_RS" 996*d950fb63SShri Abhyankar PetscErrorCode SNESSolveVI_RS(SNES snes) 997*d950fb63SShri Abhyankar { 998*d950fb63SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 999*d950fb63SShri Abhyankar PetscErrorCode ierr; 1000*d950fb63SShri Abhyankar PetscInt maxits,i,lits,Nis_act=0; 1001*d950fb63SShri Abhyankar PetscBool lssucceed; 1002*d950fb63SShri Abhyankar MatStructure flg = DIFFERENT_NONZERO_PATTERN; 1003*d950fb63SShri Abhyankar PetscReal gnorm,xnorm=0,ynorm; 1004*d950fb63SShri Abhyankar Vec Y,X,F,G,W; 1005*d950fb63SShri Abhyankar KSPConvergedReason kspreason; 1006*d950fb63SShri Abhyankar 1007*d950fb63SShri Abhyankar PetscFunctionBegin; 1008*d950fb63SShri Abhyankar snes->numFailures = 0; 1009*d950fb63SShri Abhyankar snes->numLinearSolveFailures = 0; 1010*d950fb63SShri Abhyankar snes->reason = SNES_CONVERGED_ITERATING; 1011*d950fb63SShri Abhyankar 1012*d950fb63SShri Abhyankar maxits = snes->max_its; /* maximum number of iterations */ 1013*d950fb63SShri Abhyankar X = snes->vec_sol; /* solution vector */ 1014*d950fb63SShri Abhyankar F = snes->vec_func; /* residual vector */ 1015*d950fb63SShri Abhyankar Y = snes->work[0]; /* work vectors */ 1016*d950fb63SShri Abhyankar G = snes->work[1]; 1017*d950fb63SShri Abhyankar W = snes->work[2]; 1018*d950fb63SShri Abhyankar 1019*d950fb63SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1020*d950fb63SShri Abhyankar snes->iter = 0; 1021*d950fb63SShri Abhyankar snes->norm = 0.0; 1022*d950fb63SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1023*d950fb63SShri Abhyankar 1024*d950fb63SShri Abhyankar ierr = SNESVIAdjustInitialGuess(X,vi->xl,vi->xu);CHKERRQ(ierr); 1025*d950fb63SShri Abhyankar ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr); 1026*d950fb63SShri Abhyankar if (snes->domainerror) { 1027*d950fb63SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1028*d950fb63SShri Abhyankar PetscFunctionReturn(0); 1029*d950fb63SShri Abhyankar } 1030*d950fb63SShri Abhyankar /* Compute Merit function */ 1031*d950fb63SShri Abhyankar ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr); 1032*d950fb63SShri Abhyankar 1033*d950fb63SShri Abhyankar ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 1034*d950fb63SShri Abhyankar ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 1035*d950fb63SShri Abhyankar if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1036*d950fb63SShri Abhyankar 1037*d950fb63SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1038*d950fb63SShri Abhyankar snes->norm = vi->phinorm; 1039*d950fb63SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1040*d950fb63SShri Abhyankar SNESLogConvHistory(snes,vi->phinorm,0); 1041*d950fb63SShri Abhyankar SNESMonitor(snes,0,vi->phinorm); 1042*d950fb63SShri Abhyankar 1043*d950fb63SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 1044*d950fb63SShri Abhyankar snes->ttol = vi->phinorm*snes->rtol; 1045*d950fb63SShri Abhyankar /* test convergence */ 1046*d950fb63SShri Abhyankar ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1047*d950fb63SShri Abhyankar if (snes->reason) PetscFunctionReturn(0); 1048*d950fb63SShri Abhyankar 1049*d950fb63SShri Abhyankar for (i=0; i<maxits; i++) { 1050*d950fb63SShri Abhyankar 1051*d950fb63SShri Abhyankar IS IS_act,IS_inact; /* _act -> active set _inact -> inactive set */ 1052*d950fb63SShri Abhyankar VecScatter scat_act,scat_inact; 1053*d950fb63SShri Abhyankar PetscInt nis_act,nis_inact,Nis_act_prev; 1054*d950fb63SShri Abhyankar Vec Y_act,Y_inact,phi_inact; 1055*d950fb63SShri Abhyankar Mat jac_inact_inact,prejac_inact_inact; 1056*d950fb63SShri Abhyankar 1057*d950fb63SShri Abhyankar /* Call general purpose update function */ 1058*d950fb63SShri Abhyankar if (snes->ops->update) { 1059*d950fb63SShri Abhyankar ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1060*d950fb63SShri Abhyankar } 1061*d950fb63SShri Abhyankar ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 1062*d950fb63SShri Abhyankar ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr); 1063*d950fb63SShri Abhyankar ierr = SNESVIComputeJacobian(snes->jacobian,snes->jacobian_pre,vi->Da,vi->Db);CHKERRQ(ierr); 1064*d950fb63SShri Abhyankar /* Create active and inactive index sets */ 1065*d950fb63SShri Abhyankar ierr = SNESVICreateIndexSets_RS(snes,X,vi->xl,vi->xu,&IS_act,&IS_inact);CHKERRQ(ierr); 1066*d950fb63SShri Abhyankar 1067*d950fb63SShri Abhyankar Nis_act_prev = Nis_act; 1068*d950fb63SShri Abhyankar /* Get sizes of active and inactive sets */ 1069*d950fb63SShri Abhyankar ierr = ISGetLocalSize(IS_act,&nis_act);CHKERRQ(ierr); 1070*d950fb63SShri Abhyankar ierr = ISGetLocalSize(IS_inact,&nis_inact);CHKERRQ(ierr); 1071*d950fb63SShri Abhyankar ierr = ISGetSize(IS_act,&Nis_act);CHKERRQ(ierr); 1072*d950fb63SShri Abhyankar 1073*d950fb63SShri Abhyankar ierr = PetscPrintf(PETSC_COMM_WORLD,"Size of active set = %d, size of inactive set = %d\n",nis_act,nis_inact);CHKERRQ(ierr); 1074*d950fb63SShri Abhyankar 1075*d950fb63SShri Abhyankar /* Create active and inactive set vectors */ 1076*d950fb63SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_inact,&phi_inact);CHKERRQ(ierr); 1077*d950fb63SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_act,&Y_act);CHKERRQ(ierr); 1078*d950fb63SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_inact,&Y_inact);CHKERRQ(ierr); 1079*d950fb63SShri Abhyankar 1080*d950fb63SShri Abhyankar /* Create inactive set submatrices */ 1081*d950fb63SShri Abhyankar ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);CHKERRQ(ierr); 1082*d950fb63SShri Abhyankar 1083*d950fb63SShri Abhyankar /* Create scatter contexts */ 1084*d950fb63SShri Abhyankar ierr = VecScatterCreate(Y,IS_act,Y_act,PETSC_NULL,&scat_act);CHKERRQ(ierr); 1085*d950fb63SShri Abhyankar ierr = VecScatterCreate(Y,IS_inact,Y_inact,PETSC_NULL,&scat_inact);CHKERRQ(ierr); 1086*d950fb63SShri Abhyankar 1087*d950fb63SShri Abhyankar /* Do a vec scatter to active and inactive set vectors */ 1088*d950fb63SShri Abhyankar ierr = VecScatterBegin(scat_inact,vi->phi,phi_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1089*d950fb63SShri Abhyankar ierr = VecScatterEnd(scat_inact,vi->phi,phi_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1090*d950fb63SShri Abhyankar 1091*d950fb63SShri Abhyankar ierr = VecScatterBegin(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1092*d950fb63SShri Abhyankar ierr = VecScatterEnd(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1093*d950fb63SShri Abhyankar 1094*d950fb63SShri Abhyankar ierr = VecScatterBegin(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1095*d950fb63SShri Abhyankar ierr = VecScatterEnd(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1096*d950fb63SShri Abhyankar 1097*d950fb63SShri Abhyankar /* Active set direction = 0*/ 1098*d950fb63SShri Abhyankar ierr = VecSet(Y_act,0);CHKERRQ(ierr); 1099*d950fb63SShri Abhyankar if (snes->jacobian != snes->jacobian_pre) { 1100*d950fb63SShri Abhyankar ierr = MatGetSubMatrix(snes->jacobian_pre,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&prejac_inact_inact);CHKERRQ(ierr); 1101*d950fb63SShri Abhyankar } else prejac_inact_inact = jac_inact_inact; 1102*d950fb63SShri Abhyankar 1103*d950fb63SShri Abhyankar if ((i != 0) && (Nis_act != Nis_act_prev)) { 1104*d950fb63SShri Abhyankar KSP kspnew,snesksp; 1105*d950fb63SShri Abhyankar PC pcnew; 1106*d950fb63SShri Abhyankar /* The active and inactive set sizes have changed so need to create a new snes->ksp object */ 1107*d950fb63SShri Abhyankar ierr = SNESGetKSP(snes,&snesksp);CHKERRQ(ierr); 1108*d950fb63SShri Abhyankar ierr = KSPCreate(((PetscObject)snes)->comm,&kspnew);CHKERRQ(ierr); 1109*d950fb63SShri Abhyankar /* Copy over snes->ksp info */ 1110*d950fb63SShri Abhyankar kspnew->pc_side = snesksp->pc_side; 1111*d950fb63SShri Abhyankar kspnew->rtol = snesksp->rtol; 1112*d950fb63SShri Abhyankar kspnew->abstol = snesksp->abstol; 1113*d950fb63SShri Abhyankar kspnew->max_it = snesksp->max_it; 1114*d950fb63SShri Abhyankar ierr = KSPSetType(kspnew,((PetscObject)snesksp)->type_name);CHKERRQ(ierr); 1115*d950fb63SShri Abhyankar ierr = KSPGetPC(kspnew,&pcnew);CHKERRQ(ierr); 1116*d950fb63SShri Abhyankar ierr = PCSetType(kspnew->pc,((PetscObject)snesksp->pc)->type_name);CHKERRQ(ierr); 1117*d950fb63SShri Abhyankar ierr = KSPDestroy(snesksp);CHKERRQ(ierr); 1118*d950fb63SShri Abhyankar snes->ksp = kspnew; 1119*d950fb63SShri Abhyankar ierr = PetscLogObjectParent(snes,kspnew);CHKERRQ(ierr); 1120*d950fb63SShri Abhyankar ierr = KSPSetFromOptions(kspnew);CHKERRQ(ierr); 1121*d950fb63SShri Abhyankar } 1122*d950fb63SShri Abhyankar 1123*d950fb63SShri Abhyankar ierr = KSPSetOperators(snes->ksp,jac_inact_inact,prejac_inact_inact,flg);CHKERRQ(ierr); 1124*d950fb63SShri Abhyankar ierr = SNES_KSPSolve(snes,snes->ksp,phi_inact,Y_inact);CHKERRQ(ierr); 1125*d950fb63SShri Abhyankar ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 1126*d950fb63SShri Abhyankar if (kspreason < 0) { 1127*d950fb63SShri Abhyankar if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 1128*d950fb63SShri 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); 1129*d950fb63SShri Abhyankar snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 1130*d950fb63SShri Abhyankar break; 1131*d950fb63SShri Abhyankar } 1132*d950fb63SShri Abhyankar } 1133*d950fb63SShri Abhyankar 1134*d950fb63SShri Abhyankar ierr = VecScatterBegin(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1135*d950fb63SShri Abhyankar ierr = VecScatterEnd(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1136*d950fb63SShri Abhyankar ierr = VecScatterBegin(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1137*d950fb63SShri Abhyankar ierr = VecScatterEnd(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1138*d950fb63SShri Abhyankar 1139*d950fb63SShri Abhyankar ierr = VecDestroy(phi_inact);CHKERRQ(ierr); 1140*d950fb63SShri Abhyankar ierr = VecDestroy(Y_act);CHKERRQ(ierr); 1141*d950fb63SShri Abhyankar ierr = VecDestroy(Y_inact);CHKERRQ(ierr); 1142*d950fb63SShri Abhyankar ierr = VecScatterDestroy(scat_act);CHKERRQ(ierr); 1143*d950fb63SShri Abhyankar ierr = VecScatterDestroy(scat_inact);CHKERRQ(ierr); 1144*d950fb63SShri Abhyankar ierr = ISDestroy(IS_act);CHKERRQ(ierr); 1145*d950fb63SShri Abhyankar ierr = ISDestroy(IS_inact);CHKERRQ(ierr); 1146*d950fb63SShri Abhyankar ierr = MatDestroy(jac_inact_inact);CHKERRQ(ierr); 1147*d950fb63SShri Abhyankar if (snes->jacobian != snes->jacobian_pre) { 1148*d950fb63SShri Abhyankar ierr = MatDestroy(prejac_inact_inact);CHKERRQ(ierr); 1149*d950fb63SShri Abhyankar } 1150*d950fb63SShri Abhyankar 1151*d950fb63SShri Abhyankar ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 1152*d950fb63SShri Abhyankar snes->linear_its += lits; 1153*d950fb63SShri Abhyankar ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 1154*d950fb63SShri Abhyankar /* 1155*d950fb63SShri Abhyankar if (vi->precheckstep) { 1156*d950fb63SShri Abhyankar PetscBool changed_y = PETSC_FALSE; 1157*d950fb63SShri Abhyankar ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr); 1158*d950fb63SShri Abhyankar } 1159*d950fb63SShri Abhyankar 1160*d950fb63SShri Abhyankar if (PetscLogPrintInfo){ 1161*d950fb63SShri Abhyankar ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 1162*d950fb63SShri Abhyankar } 1163*d950fb63SShri Abhyankar */ 1164*d950fb63SShri Abhyankar /* Compute a (scaled) negative update in the line search routine: 1165*d950fb63SShri Abhyankar Y <- X - lambda*Y 1166*d950fb63SShri Abhyankar and evaluate G = function(Y) (depends on the line search). 1167*d950fb63SShri Abhyankar */ 1168*d950fb63SShri Abhyankar ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 1169*d950fb63SShri Abhyankar ynorm = 1; gnorm = vi->phinorm; 1170*d950fb63SShri Abhyankar ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 1171*d950fb63SShri Abhyankar ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 1172*d950fb63SShri Abhyankar if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 1173*d950fb63SShri Abhyankar if (snes->domainerror) { 1174*d950fb63SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1175*d950fb63SShri Abhyankar PetscFunctionReturn(0); 1176*d950fb63SShri Abhyankar } 1177*d950fb63SShri Abhyankar if (!lssucceed) { 1178*d950fb63SShri Abhyankar if (++snes->numFailures >= snes->maxFailures) { 1179*d950fb63SShri Abhyankar PetscBool ismin; 1180*d950fb63SShri Abhyankar snes->reason = SNES_DIVERGED_LINE_SEARCH; 1181*d950fb63SShri Abhyankar ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 1182*d950fb63SShri Abhyankar if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 1183*d950fb63SShri Abhyankar break; 1184*d950fb63SShri Abhyankar } 1185*d950fb63SShri Abhyankar } 1186*d950fb63SShri Abhyankar /* Update function and solution vectors */ 1187*d950fb63SShri Abhyankar vi->phinorm = gnorm; 1188*d950fb63SShri Abhyankar vi->merit = 0.5*vi->phinorm*vi->phinorm; 1189*d950fb63SShri Abhyankar ierr = VecCopy(G,vi->phi);CHKERRQ(ierr); 1190*d950fb63SShri Abhyankar ierr = VecCopy(W,X);CHKERRQ(ierr); 1191*d950fb63SShri Abhyankar /* Monitor convergence */ 1192*d950fb63SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1193*d950fb63SShri Abhyankar snes->iter = i+1; 1194*d950fb63SShri Abhyankar snes->norm = vi->phinorm; 1195*d950fb63SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1196*d950fb63SShri Abhyankar SNESLogConvHistory(snes,snes->norm,lits); 1197*d950fb63SShri Abhyankar SNESMonitor(snes,snes->iter,snes->norm); 1198*d950fb63SShri Abhyankar /* Test for convergence, xnorm = || X || */ 1199*d950fb63SShri Abhyankar if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 1200*d950fb63SShri Abhyankar ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1201*d950fb63SShri Abhyankar if (snes->reason) break; 1202*d950fb63SShri Abhyankar } 1203*d950fb63SShri Abhyankar if (i == maxits) { 1204*d950fb63SShri Abhyankar ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 1205*d950fb63SShri Abhyankar if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1206*d950fb63SShri Abhyankar } 1207*d950fb63SShri Abhyankar PetscFunctionReturn(0); 1208*d950fb63SShri Abhyankar } 1209*d950fb63SShri Abhyankar 12102b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 12112b4ed282SShri Abhyankar /* 12122b4ed282SShri Abhyankar SNESSetUp_VI - Sets up the internal data structures for the later use 12132b4ed282SShri Abhyankar of the SNESVI nonlinear solver. 12142b4ed282SShri Abhyankar 12152b4ed282SShri Abhyankar Input Parameter: 12162b4ed282SShri Abhyankar . snes - the SNES context 12172b4ed282SShri Abhyankar . x - the solution vector 12182b4ed282SShri Abhyankar 12192b4ed282SShri Abhyankar Application Interface Routine: SNESSetUp() 12202b4ed282SShri Abhyankar 12212b4ed282SShri Abhyankar Notes: 12222b4ed282SShri Abhyankar For basic use of the SNES solvers, the user need not explicitly call 12232b4ed282SShri Abhyankar SNESSetUp(), since these actions will automatically occur during 12242b4ed282SShri Abhyankar the call to SNESSolve(). 12252b4ed282SShri Abhyankar */ 12262b4ed282SShri Abhyankar #undef __FUNCT__ 12272b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetUp_VI" 12282b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes) 12292b4ed282SShri Abhyankar { 12302b4ed282SShri Abhyankar PetscErrorCode ierr; 12312b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*) snes->data; 12322b4ed282SShri Abhyankar PetscInt i_start[3],i_end[3]; 12332b4ed282SShri Abhyankar 12342b4ed282SShri Abhyankar PetscFunctionBegin; 12352b4ed282SShri Abhyankar if (!snes->vec_sol_update) { 12362b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr); 12372b4ed282SShri Abhyankar ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr); 12382b4ed282SShri Abhyankar } 12392b4ed282SShri Abhyankar if (!snes->work) { 12402b4ed282SShri Abhyankar snes->nwork = 3; 12412b4ed282SShri Abhyankar ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 12422b4ed282SShri Abhyankar ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr); 12432b4ed282SShri Abhyankar } 12442b4ed282SShri Abhyankar 12452b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->phi); CHKERRQ(ierr); 12462b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->Da); CHKERRQ(ierr); 12472b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->Db); CHKERRQ(ierr); 12482b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->z);CHKERRQ(ierr); 12492b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->t); CHKERRQ(ierr); 12502b4ed282SShri Abhyankar 12512b4ed282SShri Abhyankar /* If the lower and upper bound on variables are not set, set it to 12522b4ed282SShri Abhyankar -Inf and Inf */ 12532b4ed282SShri Abhyankar if (!vi->xl && !vi->xu) { 12542b4ed282SShri Abhyankar vi->usersetxbounds = PETSC_FALSE; 12552b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->xl); CHKERRQ(ierr); 12562b4ed282SShri Abhyankar ierr = VecSet(vi->xl,PETSC_VI_NINF);CHKERRQ(ierr); 12572b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->xu); CHKERRQ(ierr); 12582b4ed282SShri Abhyankar ierr = VecSet(vi->xu,PETSC_VI_INF);CHKERRQ(ierr); 12592b4ed282SShri Abhyankar } else { 12602b4ed282SShri Abhyankar /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */ 12612b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr); 12622b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(vi->xl,i_start+1,i_end+1);CHKERRQ(ierr); 12632b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(vi->xu,i_start+2,i_end+2);CHKERRQ(ierr); 12642b4ed282SShri 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])) 12652b4ed282SShri 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."); 12662b4ed282SShri Abhyankar } 12672b4ed282SShri Abhyankar 12682b4ed282SShri Abhyankar vi->computeuserfunction = snes->ops->computefunction; 12691f79c099SShri Abhyankar snes->ops->computefunction = SNESVIComputeFunction; 1270a79edbefSShri Abhyankar 12712b4ed282SShri Abhyankar PetscFunctionReturn(0); 12722b4ed282SShri Abhyankar } 12732b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 12742b4ed282SShri Abhyankar /* 12752b4ed282SShri Abhyankar SNESDestroy_VI - Destroys the private SNES_VI context that was created 12762b4ed282SShri Abhyankar with SNESCreate_VI(). 12772b4ed282SShri Abhyankar 12782b4ed282SShri Abhyankar Input Parameter: 12792b4ed282SShri Abhyankar . snes - the SNES context 12802b4ed282SShri Abhyankar 12812b4ed282SShri Abhyankar Application Interface Routine: SNESDestroy() 12822b4ed282SShri Abhyankar */ 12832b4ed282SShri Abhyankar #undef __FUNCT__ 12842b4ed282SShri Abhyankar #define __FUNCT__ "SNESDestroy_VI" 12852b4ed282SShri Abhyankar PetscErrorCode SNESDestroy_VI(SNES snes) 12862b4ed282SShri Abhyankar { 12872b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*) snes->data; 12882b4ed282SShri Abhyankar PetscErrorCode ierr; 12892b4ed282SShri Abhyankar 12902b4ed282SShri Abhyankar PetscFunctionBegin; 12912b4ed282SShri Abhyankar if (snes->vec_sol_update) { 12922b4ed282SShri Abhyankar ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr); 12932b4ed282SShri Abhyankar snes->vec_sol_update = PETSC_NULL; 12942b4ed282SShri Abhyankar } 12952b4ed282SShri Abhyankar if (snes->nwork) { 12962b4ed282SShri Abhyankar ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 12972b4ed282SShri Abhyankar snes->nwork = 0; 12982b4ed282SShri Abhyankar snes->work = PETSC_NULL; 12992b4ed282SShri Abhyankar } 13002b4ed282SShri Abhyankar 13012b4ed282SShri Abhyankar /* clear vectors */ 13022b4ed282SShri Abhyankar ierr = VecDestroy(vi->phi); CHKERRQ(ierr); 13032b4ed282SShri Abhyankar ierr = VecDestroy(vi->Da); CHKERRQ(ierr); 13042b4ed282SShri Abhyankar ierr = VecDestroy(vi->Db); CHKERRQ(ierr); 13052b4ed282SShri Abhyankar ierr = VecDestroy(vi->z); CHKERRQ(ierr); 13062b4ed282SShri Abhyankar ierr = VecDestroy(vi->t); CHKERRQ(ierr); 13072b4ed282SShri Abhyankar if (!vi->usersetxbounds) { 13082b4ed282SShri Abhyankar ierr = VecDestroy(vi->xl); CHKERRQ(ierr); 13092b4ed282SShri Abhyankar ierr = VecDestroy(vi->xu); CHKERRQ(ierr); 13102b4ed282SShri Abhyankar } 1311c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1312c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr); 13132b4ed282SShri Abhyankar } 13142b4ed282SShri Abhyankar ierr = PetscFree(snes->data);CHKERRQ(ierr); 13152b4ed282SShri Abhyankar 13162b4ed282SShri Abhyankar /* clear composed functions */ 13172b4ed282SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr); 1318c92abb8aSShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr); 13192b4ed282SShri Abhyankar PetscFunctionReturn(0); 13202b4ed282SShri Abhyankar } 13212b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 13222b4ed282SShri Abhyankar #undef __FUNCT__ 13232b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNo_VI" 13242b4ed282SShri Abhyankar 13252b4ed282SShri Abhyankar /* 13262b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchNo routine in snes/impls/ls/ls.c 13272b4ed282SShri Abhyankar 13282b4ed282SShri Abhyankar */ 1329ace3abfcSBarry 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) 13302b4ed282SShri Abhyankar { 13312b4ed282SShri Abhyankar PetscErrorCode ierr; 13322b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1333ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 13342b4ed282SShri Abhyankar 13352b4ed282SShri Abhyankar PetscFunctionBegin; 13362b4ed282SShri Abhyankar *flag = PETSC_TRUE; 13372b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 13382b4ed282SShri Abhyankar ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 13392b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 13402b4ed282SShri Abhyankar if (vi->postcheckstep) { 13412b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 13422b4ed282SShri Abhyankar } 13432b4ed282SShri Abhyankar if (changed_y) { 13442b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 13452b4ed282SShri Abhyankar } 13462b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 13472b4ed282SShri Abhyankar if (!snes->domainerror) { 13482b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 13492b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 13502b4ed282SShri Abhyankar } 1351c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1352c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 1353c92abb8aSShri Abhyankar } 13542b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 13552b4ed282SShri Abhyankar PetscFunctionReturn(0); 13562b4ed282SShri Abhyankar } 13572b4ed282SShri Abhyankar 13582b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 13592b4ed282SShri Abhyankar #undef __FUNCT__ 13602b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNoNorms_VI" 13612b4ed282SShri Abhyankar 13622b4ed282SShri Abhyankar /* 13632b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c 13642b4ed282SShri Abhyankar */ 1365ace3abfcSBarry 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) 13662b4ed282SShri Abhyankar { 13672b4ed282SShri Abhyankar PetscErrorCode ierr; 13682b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1369ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 13702b4ed282SShri Abhyankar 13712b4ed282SShri Abhyankar PetscFunctionBegin; 13722b4ed282SShri Abhyankar *flag = PETSC_TRUE; 13732b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 13742b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 13752b4ed282SShri Abhyankar if (vi->postcheckstep) { 13762b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 13772b4ed282SShri Abhyankar } 13782b4ed282SShri Abhyankar if (changed_y) { 13792b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 13802b4ed282SShri Abhyankar } 13812b4ed282SShri Abhyankar 13822b4ed282SShri Abhyankar /* don't evaluate function the last time through */ 13832b4ed282SShri Abhyankar if (snes->iter < snes->max_its-1) { 13842b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 13852b4ed282SShri Abhyankar } 13862b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 13872b4ed282SShri Abhyankar PetscFunctionReturn(0); 13882b4ed282SShri Abhyankar } 13892b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 13902b4ed282SShri Abhyankar #undef __FUNCT__ 13912b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchCubic_VI" 13922b4ed282SShri Abhyankar /* 13932b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchCubic in snes/impls/ls/ls.c 13942b4ed282SShri Abhyankar */ 1395ace3abfcSBarry 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) 13962b4ed282SShri Abhyankar { 13972b4ed282SShri Abhyankar /* 13982b4ed282SShri Abhyankar Note that for line search purposes we work with with the related 13992b4ed282SShri Abhyankar minimization problem: 14002b4ed282SShri Abhyankar min z(x): R^n -> R, 14012b4ed282SShri Abhyankar where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 14022b4ed282SShri Abhyankar This function z(x) is same as the merit function for the complementarity problem. 14032b4ed282SShri Abhyankar */ 14042b4ed282SShri Abhyankar 14052b4ed282SShri Abhyankar PetscReal initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 14062b4ed282SShri Abhyankar PetscReal minlambda,lambda,lambdatemp; 14072b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 14082b4ed282SShri Abhyankar PetscScalar cinitslope; 14092b4ed282SShri Abhyankar #endif 14102b4ed282SShri Abhyankar PetscErrorCode ierr; 14112b4ed282SShri Abhyankar PetscInt count; 14122b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1413ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 14142b4ed282SShri Abhyankar MPI_Comm comm; 14152b4ed282SShri Abhyankar 14162b4ed282SShri Abhyankar PetscFunctionBegin; 14172b4ed282SShri Abhyankar ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 14182b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 14192b4ed282SShri Abhyankar *flag = PETSC_TRUE; 14202b4ed282SShri Abhyankar 14212b4ed282SShri Abhyankar ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 14222b4ed282SShri Abhyankar if (*ynorm == 0.0) { 1423c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1424c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 14252b4ed282SShri Abhyankar } 14262b4ed282SShri Abhyankar *gnorm = fnorm; 14272b4ed282SShri Abhyankar ierr = VecCopy(x,w);CHKERRQ(ierr); 14282b4ed282SShri Abhyankar ierr = VecCopy(f,g);CHKERRQ(ierr); 14292b4ed282SShri Abhyankar *flag = PETSC_FALSE; 14302b4ed282SShri Abhyankar goto theend1; 14312b4ed282SShri Abhyankar } 14322b4ed282SShri Abhyankar if (*ynorm > vi->maxstep) { /* Step too big, so scale back */ 1433c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1434c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Scaling step by %G old ynorm %G\n",vi->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr); 14352b4ed282SShri Abhyankar } 14362b4ed282SShri Abhyankar ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr); 14372b4ed282SShri Abhyankar *ynorm = vi->maxstep; 14382b4ed282SShri Abhyankar } 14392b4ed282SShri Abhyankar ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 14402b4ed282SShri Abhyankar minlambda = vi->minlambda/rellength; 14413b336b49SShri Abhyankar ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 14422b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 14433b336b49SShri Abhyankar ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 14442b4ed282SShri Abhyankar initslope = PetscRealPart(cinitslope); 14452b4ed282SShri Abhyankar #else 14463b336b49SShri Abhyankar ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 14472b4ed282SShri Abhyankar #endif 14482b4ed282SShri Abhyankar if (initslope > 0.0) initslope = -initslope; 14492b4ed282SShri Abhyankar if (initslope == 0.0) initslope = -1.0; 14502b4ed282SShri Abhyankar 14512b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 14522b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 14532b4ed282SShri Abhyankar ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 14542b4ed282SShri Abhyankar *flag = PETSC_FALSE; 14552b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 14562b4ed282SShri Abhyankar goto theend1; 14572b4ed282SShri Abhyankar } 14582b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 14592b4ed282SShri Abhyankar if (snes->domainerror) { 14602b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 14612b4ed282SShri Abhyankar PetscFunctionReturn(0); 14622b4ed282SShri Abhyankar } 14632b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 14642b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 14652b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 14662b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */ 1467c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1468c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 14692b4ed282SShri Abhyankar } 14702b4ed282SShri Abhyankar goto theend1; 14712b4ed282SShri Abhyankar } 14722b4ed282SShri Abhyankar 14732b4ed282SShri Abhyankar /* Fit points with quadratic */ 14742b4ed282SShri Abhyankar lambda = 1.0; 14752b4ed282SShri Abhyankar lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 14762b4ed282SShri Abhyankar lambdaprev = lambda; 14772b4ed282SShri Abhyankar gnormprev = *gnorm; 14782b4ed282SShri Abhyankar if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 14792b4ed282SShri Abhyankar if (lambdatemp <= .1*lambda) lambda = .1*lambda; 14802b4ed282SShri Abhyankar else lambda = lambdatemp; 14812b4ed282SShri Abhyankar 14822b4ed282SShri Abhyankar ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 14832b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 14842b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 14852b4ed282SShri Abhyankar *flag = PETSC_FALSE; 14862b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 14872b4ed282SShri Abhyankar goto theend1; 14882b4ed282SShri Abhyankar } 14892b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 14902b4ed282SShri Abhyankar if (snes->domainerror) { 14912b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 14922b4ed282SShri Abhyankar PetscFunctionReturn(0); 14932b4ed282SShri Abhyankar } 14942b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 14952b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1496c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1497c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr); 14982b4ed282SShri Abhyankar } 14992b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */ 1500c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1501c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 15022b4ed282SShri Abhyankar } 15032b4ed282SShri Abhyankar goto theend1; 15042b4ed282SShri Abhyankar } 15052b4ed282SShri Abhyankar 15062b4ed282SShri Abhyankar /* Fit points with cubic */ 15072b4ed282SShri Abhyankar count = 1; 15082b4ed282SShri Abhyankar while (PETSC_TRUE) { 15092b4ed282SShri Abhyankar if (lambda <= minlambda) { 1510c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1511c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 1512c92abb8aSShri 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); 15132b4ed282SShri Abhyankar } 15142b4ed282SShri Abhyankar *flag = PETSC_FALSE; 15152b4ed282SShri Abhyankar break; 15162b4ed282SShri Abhyankar } 15172b4ed282SShri Abhyankar t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 15182b4ed282SShri Abhyankar t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 15192b4ed282SShri Abhyankar a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 15202b4ed282SShri Abhyankar b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 15212b4ed282SShri Abhyankar d = b*b - 3*a*initslope; 15222b4ed282SShri Abhyankar if (d < 0.0) d = 0.0; 15232b4ed282SShri Abhyankar if (a == 0.0) { 15242b4ed282SShri Abhyankar lambdatemp = -initslope/(2.0*b); 15252b4ed282SShri Abhyankar } else { 15262b4ed282SShri Abhyankar lambdatemp = (-b + sqrt(d))/(3.0*a); 15272b4ed282SShri Abhyankar } 15282b4ed282SShri Abhyankar lambdaprev = lambda; 15292b4ed282SShri Abhyankar gnormprev = *gnorm; 15302b4ed282SShri Abhyankar if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 15312b4ed282SShri Abhyankar if (lambdatemp <= .1*lambda) lambda = .1*lambda; 15322b4ed282SShri Abhyankar else lambda = lambdatemp; 15332b4ed282SShri Abhyankar ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 15342b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 15352b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 15362b4ed282SShri 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); 15372b4ed282SShri Abhyankar *flag = PETSC_FALSE; 15382b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 15392b4ed282SShri Abhyankar break; 15402b4ed282SShri Abhyankar } 15412b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 15422b4ed282SShri Abhyankar if (snes->domainerror) { 15432b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 15442b4ed282SShri Abhyankar PetscFunctionReturn(0); 15452b4ed282SShri Abhyankar } 15462b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 15472b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 15482b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* is reduction enough? */ 1549c92abb8aSShri Abhyankar if (vi->lsmonitor) { 15502b4ed282SShri Abhyankar ierr = PetscPrintf(comm," Line search: Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 15512b4ed282SShri Abhyankar } 15522b4ed282SShri Abhyankar break; 15532b4ed282SShri Abhyankar } else { 1554c92abb8aSShri Abhyankar if (vi->lsmonitor) { 15552b4ed282SShri Abhyankar ierr = PetscPrintf(comm," Line search: Cubic step no good, shrinking lambda, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 15562b4ed282SShri Abhyankar } 15572b4ed282SShri Abhyankar } 15582b4ed282SShri Abhyankar count++; 15592b4ed282SShri Abhyankar } 15602b4ed282SShri Abhyankar theend1: 15612b4ed282SShri Abhyankar /* Optional user-defined check for line search step validity */ 15622b4ed282SShri Abhyankar if (vi->postcheckstep && *flag) { 15632b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 15642b4ed282SShri Abhyankar if (changed_y) { 15652b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 15662b4ed282SShri Abhyankar } 15672b4ed282SShri Abhyankar if (changed_y || changed_w) { /* recompute the function if the step has changed */ 15682b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 15692b4ed282SShri Abhyankar if (snes->domainerror) { 15702b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 15712b4ed282SShri Abhyankar PetscFunctionReturn(0); 15722b4ed282SShri Abhyankar } 15732b4ed282SShri Abhyankar ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 15742b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 15752b4ed282SShri Abhyankar ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 15762b4ed282SShri Abhyankar ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 15772b4ed282SShri Abhyankar ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 15782b4ed282SShri Abhyankar } 15792b4ed282SShri Abhyankar } 15802b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 15812b4ed282SShri Abhyankar PetscFunctionReturn(0); 15822b4ed282SShri Abhyankar } 15832b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 15842b4ed282SShri Abhyankar #undef __FUNCT__ 15852b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchQuadratic_VI" 15862b4ed282SShri Abhyankar /* 15872b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchQuadratic in snes/impls/ls/ls.c 15882b4ed282SShri Abhyankar */ 1589ace3abfcSBarry 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) 15902b4ed282SShri Abhyankar { 15912b4ed282SShri Abhyankar /* 15922b4ed282SShri Abhyankar Note that for line search purposes we work with with the related 15932b4ed282SShri Abhyankar minimization problem: 15942b4ed282SShri Abhyankar min z(x): R^n -> R, 15952b4ed282SShri Abhyankar where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 15962b4ed282SShri Abhyankar z(x) is the same as the merit function for the complementarity problem 15972b4ed282SShri Abhyankar */ 15982b4ed282SShri Abhyankar PetscReal initslope,minlambda,lambda,lambdatemp,rellength; 15992b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 16002b4ed282SShri Abhyankar PetscScalar cinitslope; 16012b4ed282SShri Abhyankar #endif 16022b4ed282SShri Abhyankar PetscErrorCode ierr; 16032b4ed282SShri Abhyankar PetscInt count; 16042b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1605ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 16062b4ed282SShri Abhyankar 16072b4ed282SShri Abhyankar PetscFunctionBegin; 16082b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 16092b4ed282SShri Abhyankar *flag = PETSC_TRUE; 16102b4ed282SShri Abhyankar 16112b4ed282SShri Abhyankar ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 16122b4ed282SShri Abhyankar if (*ynorm == 0.0) { 1613c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1614c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr); 16152b4ed282SShri Abhyankar } 16162b4ed282SShri Abhyankar *gnorm = fnorm; 16172b4ed282SShri Abhyankar ierr = VecCopy(x,w);CHKERRQ(ierr); 16182b4ed282SShri Abhyankar ierr = VecCopy(f,g);CHKERRQ(ierr); 16192b4ed282SShri Abhyankar *flag = PETSC_FALSE; 16202b4ed282SShri Abhyankar goto theend2; 16212b4ed282SShri Abhyankar } 16222b4ed282SShri Abhyankar if (*ynorm > vi->maxstep) { /* Step too big, so scale back */ 16232b4ed282SShri Abhyankar ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr); 16242b4ed282SShri Abhyankar *ynorm = vi->maxstep; 16252b4ed282SShri Abhyankar } 16262b4ed282SShri Abhyankar ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 16272b4ed282SShri Abhyankar minlambda = vi->minlambda/rellength; 16283b336b49SShri Abhyankar ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 16292b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 16302b4ed282SShri Abhyankar ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 16312b4ed282SShri Abhyankar initslope = PetscRealPart(cinitslope); 16322b4ed282SShri Abhyankar #else 16333b336b49SShri Abhyankar ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 16342b4ed282SShri Abhyankar #endif 16352b4ed282SShri Abhyankar if (initslope > 0.0) initslope = -initslope; 16362b4ed282SShri Abhyankar if (initslope == 0.0) initslope = -1.0; 16372b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr); 16382b4ed282SShri Abhyankar 16392b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 16402b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 16412b4ed282SShri Abhyankar ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 16422b4ed282SShri Abhyankar *flag = PETSC_FALSE; 16432b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 16442b4ed282SShri Abhyankar goto theend2; 16452b4ed282SShri Abhyankar } 16462b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 16472b4ed282SShri Abhyankar if (snes->domainerror) { 16482b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 16492b4ed282SShri Abhyankar PetscFunctionReturn(0); 16502b4ed282SShri Abhyankar } 16512b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 16522b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 16532b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */ 1654c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1655c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 16562b4ed282SShri Abhyankar } 16572b4ed282SShri Abhyankar goto theend2; 16582b4ed282SShri Abhyankar } 16592b4ed282SShri Abhyankar 16602b4ed282SShri Abhyankar /* Fit points with quadratic */ 16612b4ed282SShri Abhyankar lambda = 1.0; 16622b4ed282SShri Abhyankar count = 1; 16632b4ed282SShri Abhyankar while (PETSC_TRUE) { 16642b4ed282SShri Abhyankar if (lambda <= minlambda) { /* bad luck; use full step */ 1665c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1666c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr); 1667c92abb8aSShri 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); 16682b4ed282SShri Abhyankar } 16692b4ed282SShri Abhyankar ierr = VecCopy(x,w);CHKERRQ(ierr); 16702b4ed282SShri Abhyankar *flag = PETSC_FALSE; 16712b4ed282SShri Abhyankar break; 16722b4ed282SShri Abhyankar } 16732b4ed282SShri Abhyankar lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 16742b4ed282SShri Abhyankar if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 16752b4ed282SShri Abhyankar if (lambdatemp <= .1*lambda) lambda = .1*lambda; 16762b4ed282SShri Abhyankar else lambda = lambdatemp; 16772b4ed282SShri Abhyankar 16782b4ed282SShri Abhyankar ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 16792b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 16802b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 16812b4ed282SShri 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); 16822b4ed282SShri Abhyankar *flag = PETSC_FALSE; 16832b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 16842b4ed282SShri Abhyankar break; 16852b4ed282SShri Abhyankar } 16862b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 16872b4ed282SShri Abhyankar if (snes->domainerror) { 16882b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 16892b4ed282SShri Abhyankar PetscFunctionReturn(0); 16902b4ed282SShri Abhyankar } 16912b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 16922b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 16932b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */ 1694c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1695c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line Search: Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr); 16962b4ed282SShri Abhyankar } 16972b4ed282SShri Abhyankar break; 16982b4ed282SShri Abhyankar } 16992b4ed282SShri Abhyankar count++; 17002b4ed282SShri Abhyankar } 17012b4ed282SShri Abhyankar theend2: 17022b4ed282SShri Abhyankar /* Optional user-defined check for line search step validity */ 17032b4ed282SShri Abhyankar if (vi->postcheckstep) { 17042b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 17052b4ed282SShri Abhyankar if (changed_y) { 17062b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 17072b4ed282SShri Abhyankar } 17082b4ed282SShri Abhyankar if (changed_y || changed_w) { /* recompute the function if the step has changed */ 17092b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g); 17102b4ed282SShri Abhyankar if (snes->domainerror) { 17112b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 17122b4ed282SShri Abhyankar PetscFunctionReturn(0); 17132b4ed282SShri Abhyankar } 17142b4ed282SShri Abhyankar ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 17152b4ed282SShri Abhyankar ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 17162b4ed282SShri Abhyankar ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 17172b4ed282SShri Abhyankar ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 17182b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 17192b4ed282SShri Abhyankar } 17202b4ed282SShri Abhyankar } 17212b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 17222b4ed282SShri Abhyankar PetscFunctionReturn(0); 17232b4ed282SShri Abhyankar } 17242b4ed282SShri Abhyankar 1725ace3abfcSBarry 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*/ 17262b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 17272b4ed282SShri Abhyankar EXTERN_C_BEGIN 17282b4ed282SShri Abhyankar #undef __FUNCT__ 17292b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchSet_VI" 17302b4ed282SShri Abhyankar PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_VI(SNES snes,FCN2 func,void *lsctx) 17312b4ed282SShri Abhyankar { 17322b4ed282SShri Abhyankar PetscFunctionBegin; 17332b4ed282SShri Abhyankar ((SNES_VI *)(snes->data))->LineSearch = func; 17342b4ed282SShri Abhyankar ((SNES_VI *)(snes->data))->lsP = lsctx; 17352b4ed282SShri Abhyankar PetscFunctionReturn(0); 17362b4ed282SShri Abhyankar } 17372b4ed282SShri Abhyankar EXTERN_C_END 17382b4ed282SShri Abhyankar 17392b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 17402b4ed282SShri Abhyankar EXTERN_C_BEGIN 17412b4ed282SShri Abhyankar #undef __FUNCT__ 17422b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchSetMonitor_VI" 1743ace3abfcSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetMonitor_VI(SNES snes,PetscBool flg) 17442b4ed282SShri Abhyankar { 17452b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 17462b4ed282SShri Abhyankar PetscErrorCode ierr; 17472b4ed282SShri Abhyankar 17482b4ed282SShri Abhyankar PetscFunctionBegin; 1749c92abb8aSShri Abhyankar if (flg && !vi->lsmonitor) { 1750c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",((PetscObject)snes)->tablevel,&vi->lsmonitor);CHKERRQ(ierr); 1751c92abb8aSShri Abhyankar } else if (!flg && vi->lsmonitor) { 1752c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr); 17532b4ed282SShri Abhyankar } 17542b4ed282SShri Abhyankar PetscFunctionReturn(0); 17552b4ed282SShri Abhyankar } 17562b4ed282SShri Abhyankar EXTERN_C_END 17572b4ed282SShri Abhyankar 17582b4ed282SShri Abhyankar /* 17592b4ed282SShri Abhyankar SNESView_VI - Prints info from the SNESVI data structure. 17602b4ed282SShri Abhyankar 17612b4ed282SShri Abhyankar Input Parameters: 17622b4ed282SShri Abhyankar . SNES - the SNES context 17632b4ed282SShri Abhyankar . viewer - visualization context 17642b4ed282SShri Abhyankar 17652b4ed282SShri Abhyankar Application Interface Routine: SNESView() 17662b4ed282SShri Abhyankar */ 17672b4ed282SShri Abhyankar #undef __FUNCT__ 17682b4ed282SShri Abhyankar #define __FUNCT__ "SNESView_VI" 17692b4ed282SShri Abhyankar static PetscErrorCode SNESView_VI(SNES snes,PetscViewer viewer) 17702b4ed282SShri Abhyankar { 17712b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI *)snes->data; 177278c4b9d3SShri Abhyankar const char *cstr,*tstr; 17732b4ed282SShri Abhyankar PetscErrorCode ierr; 1774ace3abfcSBarry Smith PetscBool iascii; 17752b4ed282SShri Abhyankar 17762b4ed282SShri Abhyankar PetscFunctionBegin; 17772b4ed282SShri Abhyankar ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 17782b4ed282SShri Abhyankar if (iascii) { 17792b4ed282SShri Abhyankar if (vi->LineSearch == SNESLineSearchNo_VI) cstr = "SNESLineSearchNo"; 17802b4ed282SShri Abhyankar else if (vi->LineSearch == SNESLineSearchQuadratic_VI) cstr = "SNESLineSearchQuadratic"; 17812b4ed282SShri Abhyankar else if (vi->LineSearch == SNESLineSearchCubic_VI) cstr = "SNESLineSearchCubic"; 17822b4ed282SShri Abhyankar else cstr = "unknown"; 178378c4b9d3SShri Abhyankar if (snes->ops->solve == SNESSolveVI_SS) tstr = "Semismooth"; 178478c4b9d3SShri Abhyankar else if (snes->ops->solve == SNESSolveVI_AS) tstr = "Active Set"; 178578c4b9d3SShri Abhyankar else tstr = "unknown"; 178678c4b9d3SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," VI algorithm: %s\n",tstr);CHKERRQ(ierr); 17872b4ed282SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 17882b4ed282SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," alpha=%G, maxstep=%G, minlambda=%G\n",vi->alpha,vi->maxstep,vi->minlambda);CHKERRQ(ierr); 17892b4ed282SShri Abhyankar } else { 17902b4ed282SShri Abhyankar SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ VI",((PetscObject)viewer)->type_name); 17912b4ed282SShri Abhyankar } 17922b4ed282SShri Abhyankar PetscFunctionReturn(0); 17932b4ed282SShri Abhyankar } 17942b4ed282SShri Abhyankar 17952b4ed282SShri Abhyankar /* 17962b4ed282SShri Abhyankar SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. 17972b4ed282SShri Abhyankar 17982b4ed282SShri Abhyankar Input Parameters: 17992b4ed282SShri Abhyankar . snes - the SNES context. 18002b4ed282SShri Abhyankar . xl - lower bound. 18012b4ed282SShri Abhyankar . xu - upper bound. 18022b4ed282SShri Abhyankar 18032b4ed282SShri Abhyankar Notes: 18042b4ed282SShri Abhyankar If this routine is not called then the lower and upper bounds are set to 180584c105d7SBarry Smith PETSC_VI_INF and PETSC_VI_NINF respectively during SNESSetUp(). 180684c105d7SBarry Smith 18072b4ed282SShri Abhyankar */ 18082b4ed282SShri Abhyankar 18092b4ed282SShri Abhyankar #undef __FUNCT__ 18102b4ed282SShri Abhyankar #define __FUNCT__ "SNESVISetVariableBounds" 181169c03620SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu) 18122b4ed282SShri Abhyankar { 18132b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 18142b4ed282SShri Abhyankar 18152b4ed282SShri Abhyankar PetscFunctionBegin; 18162b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 18172b4ed282SShri Abhyankar PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 18182b4ed282SShri Abhyankar PetscValidHeaderSpecific(xu,VEC_CLASSID,3); 18192b4ed282SShri Abhyankar 18202b4ed282SShri Abhyankar /* Check if SNESSetFunction is called */ 18212b4ed282SShri Abhyankar if(!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 18222b4ed282SShri Abhyankar 18232b4ed282SShri Abhyankar /* Check if the vector sizes are compatible for lower and upper bounds */ 18242b4ed282SShri 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); 18252b4ed282SShri 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); 18262b4ed282SShri Abhyankar vi->usersetxbounds = PETSC_TRUE; 18272b4ed282SShri Abhyankar vi->xl = xl; 18282b4ed282SShri Abhyankar vi->xu = xu; 18292b4ed282SShri Abhyankar 18302b4ed282SShri Abhyankar PetscFunctionReturn(0); 18312b4ed282SShri Abhyankar } 18322b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 18332b4ed282SShri Abhyankar /* 18342b4ed282SShri Abhyankar SNESSetFromOptions_VI - Sets various parameters for the SNESVI method. 18352b4ed282SShri Abhyankar 18362b4ed282SShri Abhyankar Input Parameter: 18372b4ed282SShri Abhyankar . snes - the SNES context 18382b4ed282SShri Abhyankar 18392b4ed282SShri Abhyankar Application Interface Routine: SNESSetFromOptions() 18402b4ed282SShri Abhyankar */ 18412b4ed282SShri Abhyankar #undef __FUNCT__ 18422b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetFromOptions_VI" 18432b4ed282SShri Abhyankar static PetscErrorCode SNESSetFromOptions_VI(SNES snes) 18442b4ed282SShri Abhyankar { 18452b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI *)snes->data; 18462b4ed282SShri Abhyankar const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 1847*d950fb63SShri Abhyankar const char *vies[] = {"ss","as","rs"}; 18482b4ed282SShri Abhyankar PetscErrorCode ierr; 18492b4ed282SShri Abhyankar PetscInt indx; 185069c03620SShri Abhyankar PetscBool flg,set,flg2; 18512b4ed282SShri Abhyankar 18522b4ed282SShri Abhyankar PetscFunctionBegin; 18532b4ed282SShri Abhyankar ierr = PetscOptionsHead("SNES semismooth method options");CHKERRQ(ierr); 18542b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_alpha","Function norm must decrease by","None",vi->alpha,&vi->alpha,0);CHKERRQ(ierr); 18552b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_maxstep","Step must be less than","None",vi->maxstep,&vi->maxstep,0);CHKERRQ(ierr); 18562b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_minlambda","Minimum lambda allowed","None",vi->minlambda,&vi->minlambda,0);CHKERRQ(ierr); 18572b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_const_tol","constraint tolerance","None",vi->const_tol,&vi->const_tol,0);CHKERRQ(ierr); 1858acfcf0e5SJed Brown ierr = PetscOptionsBool("-snes_vi_lsmonitor","Print progress of line searches","SNESLineSearchSetMonitor",vi->lsmonitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 18592b4ed282SShri Abhyankar if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);} 1860*d950fb63SShri Abhyankar ierr = PetscOptionsEList("-snes_vi_type","Semismooth algorithm used","",vies,3,"ss",&indx,&flg2);CHKERRQ(ierr); 186169c03620SShri Abhyankar if (flg2) { 186269c03620SShri Abhyankar switch (indx) { 186369c03620SShri Abhyankar case 0: 186469c03620SShri Abhyankar snes->ops->solve = SNESSolveVI_SS; 186569c03620SShri Abhyankar break; 186669c03620SShri Abhyankar case 1: 186769c03620SShri Abhyankar snes->ops->solve = SNESSolveVI_AS; 186869c03620SShri Abhyankar break; 1869*d950fb63SShri Abhyankar case 2: 1870*d950fb63SShri Abhyankar snes->ops->solve = SNESSolveVI_RS; 1871*d950fb63SShri Abhyankar break; 187269c03620SShri Abhyankar } 187369c03620SShri Abhyankar } 1874c92abb8aSShri Abhyankar ierr = PetscOptionsEList("-snes_vi_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 18752b4ed282SShri Abhyankar if (flg) { 18762b4ed282SShri Abhyankar switch (indx) { 18772b4ed282SShri Abhyankar case 0: 18782b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr); 18792b4ed282SShri Abhyankar break; 18802b4ed282SShri Abhyankar case 1: 18812b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr); 18822b4ed282SShri Abhyankar break; 18832b4ed282SShri Abhyankar case 2: 18842b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr); 18852b4ed282SShri Abhyankar break; 18862b4ed282SShri Abhyankar case 3: 18872b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr); 18882b4ed282SShri Abhyankar break; 18892b4ed282SShri Abhyankar } 18902b4ed282SShri Abhyankar } 18912b4ed282SShri Abhyankar ierr = PetscOptionsTail();CHKERRQ(ierr); 18922b4ed282SShri Abhyankar PetscFunctionReturn(0); 18932b4ed282SShri Abhyankar } 18942b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 18952b4ed282SShri Abhyankar /*MC 18962b4ed282SShri Abhyankar SNESVI - Semismooth newton method based nonlinear solver that uses a line search 18972b4ed282SShri Abhyankar 18982b4ed282SShri Abhyankar Options Database: 18992b4ed282SShri Abhyankar + -snes_vi [cubic,quadratic,basic,basicnonorms] - Selects line search 19002b4ed282SShri Abhyankar . -snes_vi_alpha <alpha> - Sets alpha 19012b4ed282SShri 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) 19022b4ed282SShri Abhyankar . -snes_vi_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 19032b4ed282SShri Abhyankar - -snes_vi_monitor - print information about progress of line searches 19042b4ed282SShri Abhyankar 19052b4ed282SShri Abhyankar 19062b4ed282SShri Abhyankar Level: beginner 19072b4ed282SShri Abhyankar 19082b4ed282SShri Abhyankar .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 19092b4ed282SShri Abhyankar SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 19102b4ed282SShri Abhyankar SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams() 19112b4ed282SShri Abhyankar 19122b4ed282SShri Abhyankar M*/ 19132b4ed282SShri Abhyankar EXTERN_C_BEGIN 19142b4ed282SShri Abhyankar #undef __FUNCT__ 19152b4ed282SShri Abhyankar #define __FUNCT__ "SNESCreate_VI" 19162b4ed282SShri Abhyankar PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_VI(SNES snes) 19172b4ed282SShri Abhyankar { 19182b4ed282SShri Abhyankar PetscErrorCode ierr; 19192b4ed282SShri Abhyankar SNES_VI *vi; 19202b4ed282SShri Abhyankar 19212b4ed282SShri Abhyankar PetscFunctionBegin; 19222b4ed282SShri Abhyankar snes->ops->setup = SNESSetUp_VI; 192369c03620SShri Abhyankar snes->ops->solve = SNESSolveVI_SS; 19242b4ed282SShri Abhyankar snes->ops->destroy = SNESDestroy_VI; 19252b4ed282SShri Abhyankar snes->ops->setfromoptions = SNESSetFromOptions_VI; 19262b4ed282SShri Abhyankar snes->ops->view = SNESView_VI; 19272b4ed282SShri Abhyankar snes->ops->converged = SNESDefaultConverged_VI; 19282b4ed282SShri Abhyankar 19292b4ed282SShri Abhyankar ierr = PetscNewLog(snes,SNES_VI,&vi);CHKERRQ(ierr); 19302b4ed282SShri Abhyankar snes->data = (void*)vi; 19312b4ed282SShri Abhyankar vi->alpha = 1.e-4; 19322b4ed282SShri Abhyankar vi->maxstep = 1.e8; 19332b4ed282SShri Abhyankar vi->minlambda = 1.e-12; 19342b4ed282SShri Abhyankar vi->LineSearch = SNESLineSearchCubic_VI; 19352b4ed282SShri Abhyankar vi->lsP = PETSC_NULL; 19362b4ed282SShri Abhyankar vi->postcheckstep = PETSC_NULL; 19372b4ed282SShri Abhyankar vi->postcheck = PETSC_NULL; 19382b4ed282SShri Abhyankar vi->precheckstep = PETSC_NULL; 19392b4ed282SShri Abhyankar vi->precheck = PETSC_NULL; 19402b4ed282SShri Abhyankar vi->const_tol = 2.2204460492503131e-16; 19412b4ed282SShri Abhyankar vi->computessfunction = ComputeFischerFunction; 19422b4ed282SShri Abhyankar 19432b4ed282SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_VI",SNESLineSearchSetMonitor_VI);CHKERRQ(ierr); 19442b4ed282SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_VI",SNESLineSearchSet_VI);CHKERRQ(ierr); 19452b4ed282SShri Abhyankar 19462b4ed282SShri Abhyankar PetscFunctionReturn(0); 19472b4ed282SShri Abhyankar } 19482b4ed282SShri Abhyankar EXTERN_C_END 1949