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 6*9308c137SBarry Smith #undef __FUNCT__ 7*9308c137SBarry Smith #define __FUNCT__ "SNESMonitorVI" 8*9308c137SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorVI(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 9*9308c137SBarry Smith { 10*9308c137SBarry Smith PetscErrorCode ierr; 11*9308c137SBarry Smith SNES_VI *vi = (SNES_VI*)snes->data; 12*9308c137SBarry Smith PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor) dummy; 13*9308c137SBarry Smith const PetscScalar *x,*xl,*xu,*f; 14*9308c137SBarry Smith Vec ff; 15*9308c137SBarry Smith PetscInt i,n; 16*9308c137SBarry Smith PetscReal rnorm,fnorm; 17*9308c137SBarry Smith 18*9308c137SBarry Smith PetscFunctionBegin; 19*9308c137SBarry Smith if (!dummy) { 20*9308c137SBarry Smith ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",0,&viewer);CHKERRQ(ierr); 21*9308c137SBarry Smith } 22*9308c137SBarry Smith ierr = VecGetLocalSize(snes->vec_sol,&n);CHKERRQ(ierr); 23*9308c137SBarry Smith ierr = VecDuplicate(snes->vec_sol,&ff);CHKERRQ(ierr); 24*9308c137SBarry Smith ierr = SNESComputeFunction(snes,snes->vec_sol,ff);CHKERRQ(ierr); 25*9308c137SBarry Smith ierr = VecGetArrayRead(vi->xl,&xl);CHKERRQ(ierr); 26*9308c137SBarry Smith ierr = VecGetArrayRead(vi->xu,&xu);CHKERRQ(ierr); 27*9308c137SBarry Smith ierr = VecGetArrayRead(snes->vec_sol,&x);CHKERRQ(ierr); 28*9308c137SBarry Smith ierr = VecGetArrayRead(ff,&f);CHKERRQ(ierr); 29*9308c137SBarry Smith 30*9308c137SBarry Smith rnorm = 0.0; 31*9308c137SBarry Smith for (i=0; i<n; i++) { 32*9308c137SBarry Smith if ((x[i] > xl[i] + 1.e-8) && (x[i] < xu[i] - 1.e-8)) rnorm += f[i]*f[i]; 33*9308c137SBarry Smith } 34*9308c137SBarry Smith ierr = VecRestoreArrayRead(ff,&f);CHKERRQ(ierr); 35*9308c137SBarry Smith ierr = VecRestoreArrayRead(vi->xl,&xl);CHKERRQ(ierr); 36*9308c137SBarry Smith ierr = VecRestoreArrayRead(vi->xu,&xu);CHKERRQ(ierr); 37*9308c137SBarry Smith ierr = VecRestoreArrayRead(snes->vec_sol,&x);CHKERRQ(ierr); 38*9308c137SBarry Smith ierr = VecDestroy(ff);CHKERRQ(ierr); 39*9308c137SBarry Smith ierr = MPI_Allreduce(&rnorm,&fnorm,1,MPIU_REAL,MPI_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr); 40*9308c137SBarry Smith fnorm = sqrt(fnorm); 41*9308c137SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(viewer,"%3D SNES VI Function norm %14.12e \n",its,fnorm);CHKERRQ(ierr); 42*9308c137SBarry Smith if (!dummy) { 43*9308c137SBarry Smith ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr); 44*9308c137SBarry Smith } 45*9308c137SBarry Smith PetscFunctionReturn(0); 46*9308c137SBarry Smith } 47*9308c137SBarry Smith 482b4ed282SShri Abhyankar /* 492b4ed282SShri Abhyankar Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function, 502b4ed282SShri 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 512b4ed282SShri Abhyankar 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 522b4ed282SShri Abhyankar for this trick. One assumes that the probability that W is in the null space of J is very, very small. 532b4ed282SShri Abhyankar */ 542b4ed282SShri Abhyankar #undef __FUNCT__ 552b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckLocalMin_Private" 56ace3abfcSBarry Smith PetscErrorCode SNESVICheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin) 572b4ed282SShri Abhyankar { 582b4ed282SShri Abhyankar PetscReal a1; 592b4ed282SShri Abhyankar PetscErrorCode ierr; 60ace3abfcSBarry Smith PetscBool hastranspose; 612b4ed282SShri Abhyankar 622b4ed282SShri Abhyankar PetscFunctionBegin; 632b4ed282SShri Abhyankar *ismin = PETSC_FALSE; 642b4ed282SShri Abhyankar ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 652b4ed282SShri Abhyankar if (hastranspose) { 662b4ed282SShri Abhyankar /* Compute || J^T F|| */ 672b4ed282SShri Abhyankar ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 682b4ed282SShri Abhyankar ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 692b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"|| J^T F|| %G near zero implies found a local minimum\n",a1/fnorm);CHKERRQ(ierr); 702b4ed282SShri Abhyankar if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 712b4ed282SShri Abhyankar } else { 722b4ed282SShri Abhyankar Vec work; 732b4ed282SShri Abhyankar PetscScalar result; 742b4ed282SShri Abhyankar PetscReal wnorm; 752b4ed282SShri Abhyankar 762b4ed282SShri Abhyankar ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr); 772b4ed282SShri Abhyankar ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr); 782b4ed282SShri Abhyankar ierr = VecDuplicate(W,&work);CHKERRQ(ierr); 792b4ed282SShri Abhyankar ierr = MatMult(A,W,work);CHKERRQ(ierr); 802b4ed282SShri Abhyankar ierr = VecDot(F,work,&result);CHKERRQ(ierr); 812b4ed282SShri Abhyankar ierr = VecDestroy(work);CHKERRQ(ierr); 822b4ed282SShri Abhyankar a1 = PetscAbsScalar(result)/(fnorm*wnorm); 832b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %G near zero implies found a local minimum\n",a1);CHKERRQ(ierr); 842b4ed282SShri Abhyankar if (a1 < 1.e-4) *ismin = PETSC_TRUE; 852b4ed282SShri Abhyankar } 862b4ed282SShri Abhyankar PetscFunctionReturn(0); 872b4ed282SShri Abhyankar } 882b4ed282SShri Abhyankar 892b4ed282SShri Abhyankar /* 902b4ed282SShri Abhyankar Checks if J^T(F - J*X) = 0 912b4ed282SShri Abhyankar */ 922b4ed282SShri Abhyankar #undef __FUNCT__ 932b4ed282SShri Abhyankar #define __FUNCT__ "SNESVICheckResidual_Private" 942b4ed282SShri Abhyankar PetscErrorCode SNESVICheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2) 952b4ed282SShri Abhyankar { 962b4ed282SShri Abhyankar PetscReal a1,a2; 972b4ed282SShri Abhyankar PetscErrorCode ierr; 98ace3abfcSBarry Smith PetscBool hastranspose; 992b4ed282SShri Abhyankar 1002b4ed282SShri Abhyankar PetscFunctionBegin; 1012b4ed282SShri Abhyankar ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 1022b4ed282SShri Abhyankar if (hastranspose) { 1032b4ed282SShri Abhyankar ierr = MatMult(A,X,W1);CHKERRQ(ierr); 1042b4ed282SShri Abhyankar ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr); 1052b4ed282SShri Abhyankar 1062b4ed282SShri Abhyankar /* Compute || J^T W|| */ 1072b4ed282SShri Abhyankar ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 1082b4ed282SShri Abhyankar ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 1092b4ed282SShri Abhyankar ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 1102b4ed282SShri Abhyankar if (a1 != 0.0) { 1112b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %G near zero implies inconsistent rhs\n",a2/a1);CHKERRQ(ierr); 1122b4ed282SShri Abhyankar } 1132b4ed282SShri Abhyankar } 1142b4ed282SShri Abhyankar PetscFunctionReturn(0); 1152b4ed282SShri Abhyankar } 1162b4ed282SShri Abhyankar 1172b4ed282SShri Abhyankar /* 1182b4ed282SShri Abhyankar SNESDefaultConverged_VI - Checks the convergence of the semismooth newton algorithm. 1192b4ed282SShri Abhyankar 1202b4ed282SShri Abhyankar Notes: 1212b4ed282SShri Abhyankar The convergence criterion currently implemented is 1222b4ed282SShri Abhyankar merit < abstol 1232b4ed282SShri Abhyankar merit < rtol*merit_initial 1242b4ed282SShri Abhyankar */ 1252b4ed282SShri Abhyankar #undef __FUNCT__ 1262b4ed282SShri Abhyankar #define __FUNCT__ "SNESDefaultConverged_VI" 1272b4ed282SShri Abhyankar PetscErrorCode SNESDefaultConverged_VI(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gradnorm,PetscReal merit,SNESConvergedReason *reason,void *dummy) 1282b4ed282SShri Abhyankar { 1292b4ed282SShri Abhyankar PetscErrorCode ierr; 1302b4ed282SShri Abhyankar 1312b4ed282SShri Abhyankar PetscFunctionBegin; 1322b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1332b4ed282SShri Abhyankar PetscValidPointer(reason,6); 1342b4ed282SShri Abhyankar 1352b4ed282SShri Abhyankar *reason = SNES_CONVERGED_ITERATING; 1362b4ed282SShri Abhyankar 1372b4ed282SShri Abhyankar if (!it) { 1382b4ed282SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 1392b4ed282SShri Abhyankar snes->ttol = merit*snes->rtol; 1402b4ed282SShri Abhyankar } 1412b4ed282SShri Abhyankar if (merit != merit) { 1422b4ed282SShri Abhyankar ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr); 1432b4ed282SShri Abhyankar *reason = SNES_DIVERGED_FNORM_NAN; 1442b4ed282SShri Abhyankar } else if (merit < snes->abstol) { 1452b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"Converged due to merit function %G < %G\n",merit,snes->abstol);CHKERRQ(ierr); 1462b4ed282SShri Abhyankar *reason = SNES_CONVERGED_FNORM_ABS; 1472b4ed282SShri Abhyankar } else if (snes->nfuncs >= snes->max_funcs) { 1482b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr); 1492b4ed282SShri Abhyankar *reason = SNES_DIVERGED_FUNCTION_COUNT; 1502b4ed282SShri Abhyankar } 1512b4ed282SShri Abhyankar 1522b4ed282SShri Abhyankar if (it && !*reason) { 1532b4ed282SShri Abhyankar if (merit < snes->ttol) { 1542b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"Converged due to merit function %G < %G (relative tolerance)\n",merit,snes->ttol);CHKERRQ(ierr); 1552b4ed282SShri Abhyankar *reason = SNES_CONVERGED_FNORM_RELATIVE; 1562b4ed282SShri Abhyankar } 1572b4ed282SShri Abhyankar } 1582b4ed282SShri Abhyankar PetscFunctionReturn(0); 1592b4ed282SShri Abhyankar } 1602b4ed282SShri Abhyankar 1612b4ed282SShri Abhyankar /* 1622b4ed282SShri Abhyankar SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem. 1632b4ed282SShri Abhyankar 1642b4ed282SShri Abhyankar Input Parameter: 1652b4ed282SShri Abhyankar . phi - the semismooth function 1662b4ed282SShri Abhyankar 1672b4ed282SShri Abhyankar Output Parameter: 1682b4ed282SShri Abhyankar . merit - the merit function 1692b4ed282SShri Abhyankar . phinorm - ||phi|| 1702b4ed282SShri Abhyankar 1712b4ed282SShri Abhyankar Notes: 1722b4ed282SShri Abhyankar The merit function for the mixed complementarity problem is defined as 1732b4ed282SShri Abhyankar merit = 0.5*phi^T*phi 1742b4ed282SShri Abhyankar */ 1752b4ed282SShri Abhyankar #undef __FUNCT__ 1762b4ed282SShri Abhyankar #define __FUNCT__ "SNESVIComputeMeritFunction" 1772b4ed282SShri Abhyankar static PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal* merit,PetscReal* phinorm) 1782b4ed282SShri Abhyankar { 1792b4ed282SShri Abhyankar PetscErrorCode ierr; 1802b4ed282SShri Abhyankar 1812b4ed282SShri Abhyankar PetscFunctionBegin; 1822b4ed282SShri Abhyankar ierr = VecNormBegin(phi,NORM_2,phinorm); 1832b4ed282SShri Abhyankar ierr = VecNormEnd(phi,NORM_2,phinorm); 1842b4ed282SShri Abhyankar 1852b4ed282SShri Abhyankar *merit = 0.5*(*phinorm)*(*phinorm); 1862b4ed282SShri Abhyankar PetscFunctionReturn(0); 1872b4ed282SShri Abhyankar } 1882b4ed282SShri Abhyankar 1892b4ed282SShri Abhyankar /* 1902b4ed282SShri Abhyankar ComputeFischerFunction - Computes the semismooth fischer burmeister function for a mixed complementarity equation. 1912b4ed282SShri Abhyankar 1922b4ed282SShri Abhyankar Notes: 1932b4ed282SShri Abhyankar The Fischer-Burmeister function is defined as 1942b4ed282SShri Abhyankar ff(a,b) := sqrt(a*a + b*b) - a - b 1952b4ed282SShri Abhyankar and is used reformulate a complementarity equation as a semismooth equation. 1962b4ed282SShri Abhyankar */ 1972b4ed282SShri Abhyankar 1982b4ed282SShri Abhyankar #undef __FUNCT__ 1992b4ed282SShri Abhyankar #define __FUNCT__ "ComputeFischerFunction" 2002b4ed282SShri Abhyankar static PetscErrorCode ComputeFischerFunction(PetscScalar a, PetscScalar b, PetscScalar* ff) 2012b4ed282SShri Abhyankar { 2022b4ed282SShri Abhyankar PetscFunctionBegin; 2032b4ed282SShri Abhyankar *ff = sqrt(a*a + b*b) - a - b; 2042b4ed282SShri Abhyankar PetscFunctionReturn(0); 2052b4ed282SShri Abhyankar } 2062b4ed282SShri Abhyankar 2072b4ed282SShri Abhyankar /* 2081f79c099SShri Abhyankar SNESVIComputeFunction - Reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear equations in semismooth form. 2092b4ed282SShri Abhyankar 2102b4ed282SShri Abhyankar Input Parameters: 2112b4ed282SShri Abhyankar . snes - the SNES context 2122b4ed282SShri Abhyankar . x - current iterate 2132b4ed282SShri Abhyankar . functx - user defined function context 2142b4ed282SShri Abhyankar 2152b4ed282SShri Abhyankar Output Parameters: 2162b4ed282SShri Abhyankar . phi - Semismooth function 2172b4ed282SShri Abhyankar 2182b4ed282SShri Abhyankar Notes: 2192b4ed282SShri Abhyankar The result of this function is done by cases: 2202b4ed282SShri Abhyankar + l[i] == -infinity, u[i] == infinity -- phi[i] = -f[i] 2212b4ed282SShri Abhyankar . l[i] == -infinity, u[i] finite -- phi[i] = ss(u[i]-x[i], -f[i]) 2222b4ed282SShri Abhyankar . l[i] finite, u[i] == infinity -- phi[i] = ss(x[i]-l[i], f[i]) 2232b4ed282SShri Abhyankar . l[i] finite < u[i] finite -- phi[i] = phi(x[i]-l[i], ss(u[i]-x[i], -f[u])) 2242b4ed282SShri Abhyankar - otherwise l[i] == u[i] -- phi[i] = l[i] - x[i] 2252b4ed282SShri Abhyankar ss is the semismoothing function used to reformulate the nonlinear equation in complementarity 2262b4ed282SShri Abhyankar form to semismooth form 2272b4ed282SShri Abhyankar 2282b4ed282SShri Abhyankar */ 2292b4ed282SShri Abhyankar #undef __FUNCT__ 2301f79c099SShri Abhyankar #define __FUNCT__ "SNESVIComputeFunction" 2311f79c099SShri Abhyankar static PetscErrorCode SNESVIComputeFunction(SNES snes,Vec X,Vec phi,void* functx) 2322b4ed282SShri Abhyankar { 2332b4ed282SShri Abhyankar PetscErrorCode ierr; 2342b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 2352b4ed282SShri Abhyankar Vec Xl = vi->xl,Xu = vi->xu,F = snes->vec_func; 2362b4ed282SShri Abhyankar PetscScalar *phi_arr,*x_arr,*f_arr,*l,*u,t; 2372b4ed282SShri Abhyankar PetscInt i,nlocal; 2382b4ed282SShri Abhyankar 2392b4ed282SShri Abhyankar PetscFunctionBegin; 2402b4ed282SShri Abhyankar ierr = (*vi->computeuserfunction)(snes,X,F,functx);CHKERRQ(ierr); 2412b4ed282SShri Abhyankar ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 2422b4ed282SShri Abhyankar 2432b4ed282SShri Abhyankar ierr = VecGetArray(X,&x_arr);CHKERRQ(ierr); 2442b4ed282SShri Abhyankar ierr = VecGetArray(F,&f_arr);CHKERRQ(ierr); 2452b4ed282SShri Abhyankar ierr = VecGetArray(Xl,&l);CHKERRQ(ierr); 2462b4ed282SShri Abhyankar ierr = VecGetArray(Xu,&u);CHKERRQ(ierr); 2472b4ed282SShri Abhyankar ierr = VecGetArray(phi,&phi_arr);CHKERRQ(ierr); 2482b4ed282SShri Abhyankar 2492b4ed282SShri Abhyankar for (i=0;i < nlocal;i++) { 2502b4ed282SShri Abhyankar if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) { 2512b4ed282SShri Abhyankar phi_arr[i] = -f_arr[i]; 2522b4ed282SShri Abhyankar } 2532b4ed282SShri Abhyankar else if (l[i] <= PETSC_VI_NINF) { 2542b4ed282SShri Abhyankar t = u[i] - x_arr[i]; 2552b4ed282SShri Abhyankar ierr = ComputeFischerFunction(t,-f_arr[i],&phi_arr[i]);CHKERRQ(ierr); 2562b4ed282SShri Abhyankar phi_arr[i] = -phi_arr[i]; 2572b4ed282SShri Abhyankar } 2582b4ed282SShri Abhyankar else if (u[i] >= PETSC_VI_INF) { 2592b4ed282SShri Abhyankar t = x_arr[i] - l[i]; 2602b4ed282SShri Abhyankar ierr = ComputeFischerFunction(t,f_arr[i],&phi_arr[i]);CHKERRQ(ierr); 2612b4ed282SShri Abhyankar } 2622b4ed282SShri Abhyankar else if (l[i] == u[i]) { 2632b4ed282SShri Abhyankar phi_arr[i] = l[i] - x_arr[i]; 2642b4ed282SShri Abhyankar } 2652b4ed282SShri Abhyankar else { 2662b4ed282SShri Abhyankar t = u[i] - x_arr[i]; 2672b4ed282SShri Abhyankar ierr = ComputeFischerFunction(t,-f_arr[i],&phi_arr[i]); 2682b4ed282SShri Abhyankar t = x_arr[i] - l[i]; 2692b4ed282SShri Abhyankar ierr = ComputeFischerFunction(t,phi_arr[i],&phi_arr[i]); 2702b4ed282SShri Abhyankar } 2712b4ed282SShri Abhyankar } 2722b4ed282SShri Abhyankar 2732b4ed282SShri Abhyankar ierr = VecRestoreArray(X,&x_arr);CHKERRQ(ierr); 2742b4ed282SShri Abhyankar ierr = VecRestoreArray(F,&f_arr);CHKERRQ(ierr); 2752b4ed282SShri Abhyankar ierr = VecRestoreArray(Xl,&l);CHKERRQ(ierr); 2762b4ed282SShri Abhyankar ierr = VecRestoreArray(Xu,&u);CHKERRQ(ierr); 2772b4ed282SShri Abhyankar ierr = VecRestoreArray(phi,&phi_arr);CHKERRQ(ierr); 2782b4ed282SShri Abhyankar 2792b4ed282SShri Abhyankar PetscFunctionReturn(0); 2802b4ed282SShri Abhyankar } 2812b4ed282SShri Abhyankar 2822b4ed282SShri Abhyankar /* 283a79edbefSShri Abhyankar SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the 284a79edbefSShri Abhyankar the semismooth jacobian. 2852b4ed282SShri Abhyankar */ 2862b4ed282SShri Abhyankar #undef __FUNCT__ 287a79edbefSShri Abhyankar #define __FUNCT__ "SNESVIComputeBsubdifferentialVectors" 288a79edbefSShri Abhyankar PetscErrorCode SNESVIComputeBsubdifferentialVectors(SNES snes,Vec X,Vec F,Mat jac,Vec Da,Vec Db) 2892b4ed282SShri Abhyankar { 2902b4ed282SShri Abhyankar PetscErrorCode ierr; 2912b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 2922b4ed282SShri Abhyankar PetscScalar *l,*u,*x,*f,*da,*db,*z,*t,t1,t2,ci,di,ei; 2932b4ed282SShri Abhyankar PetscInt i,nlocal; 2942b4ed282SShri Abhyankar 2952b4ed282SShri Abhyankar PetscFunctionBegin; 2962b4ed282SShri Abhyankar 2972b4ed282SShri Abhyankar ierr = VecGetArray(X,&x);CHKERRQ(ierr); 2982b4ed282SShri Abhyankar ierr = VecGetArray(F,&f);CHKERRQ(ierr); 2992b4ed282SShri Abhyankar ierr = VecGetArray(vi->xl,&l);CHKERRQ(ierr); 3002b4ed282SShri Abhyankar ierr = VecGetArray(vi->xu,&u);CHKERRQ(ierr); 301a79edbefSShri Abhyankar ierr = VecGetArray(Da,&da);CHKERRQ(ierr); 302a79edbefSShri Abhyankar ierr = VecGetArray(Db,&db);CHKERRQ(ierr); 3032b4ed282SShri Abhyankar ierr = VecGetArray(vi->z,&z);CHKERRQ(ierr); 3042b4ed282SShri Abhyankar 3052b4ed282SShri Abhyankar ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 3062b4ed282SShri Abhyankar /* Set the elements of the vector z: 3072b4ed282SShri Abhyankar z[i] = 1 if (x[i] - l[i],f[i]) = (0,0) or (u[i] - x[i],f[i]) = (0,0) 3082b4ed282SShri Abhyankar else z[i] = 0 3092b4ed282SShri Abhyankar */ 3102b4ed282SShri Abhyankar for(i=0;i < nlocal;i++) { 3112b4ed282SShri Abhyankar da[i] = db[i] = z[i] = 0; 3122b4ed282SShri Abhyankar if(PetscAbsScalar(f[i]) <= vi->const_tol) { 3132b4ed282SShri Abhyankar if ((l[i] > PETSC_VI_NINF) && (PetscAbsScalar(x[i]-l[i]) <= vi->const_tol)) { 3142b4ed282SShri Abhyankar da[i] = 1; 3152b4ed282SShri Abhyankar z[i] = 1; 3162b4ed282SShri Abhyankar } 3172b4ed282SShri Abhyankar if ((u[i] < PETSC_VI_INF) && (PetscAbsScalar(u[i]-x[i]) <= vi->const_tol)) { 3182b4ed282SShri Abhyankar db[i] = 1; 3192b4ed282SShri Abhyankar z[i] = 1; 3202b4ed282SShri Abhyankar } 3212b4ed282SShri Abhyankar } 3222b4ed282SShri Abhyankar } 3232b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->z,&z);CHKERRQ(ierr); 324a79edbefSShri Abhyankar ierr = MatMult(jac,vi->z,vi->t);CHKERRQ(ierr); 3252b4ed282SShri Abhyankar ierr = VecGetArray(vi->t,&t);CHKERRQ(ierr); 3262b4ed282SShri Abhyankar /* Compute the elements of the diagonal perturbation vector Da and row scaling vector Db */ 3272b4ed282SShri Abhyankar for(i=0;i< nlocal;i++) { 3282b4ed282SShri Abhyankar /* Free variables */ 3292b4ed282SShri Abhyankar if ((l[i] <= PETSC_VI_NINF) && (u[i] >= PETSC_VI_INF)) { 3302b4ed282SShri Abhyankar da[i] = 0; db[i] = -1; 3312b4ed282SShri Abhyankar } 3322b4ed282SShri Abhyankar /* lower bounded variables */ 3332b4ed282SShri Abhyankar else if (u[i] >= PETSC_VI_INF) { 3342b4ed282SShri Abhyankar if (da[i] >= 1) { 3352b4ed282SShri Abhyankar t2 = PetscScalarNorm(1,t[i]); 3362b4ed282SShri Abhyankar da[i] = 1/t2 - 1; 3372b4ed282SShri Abhyankar db[i] = t[i]/t2 - 1; 3382b4ed282SShri Abhyankar } else { 3392b4ed282SShri Abhyankar t1 = x[i] - l[i]; 3402b4ed282SShri Abhyankar t2 = PetscScalarNorm(t1,f[i]); 3412b4ed282SShri Abhyankar da[i] = t1/t2 - 1; 3422b4ed282SShri Abhyankar db[i] = f[i]/t2 - 1; 3432b4ed282SShri Abhyankar } 3442b4ed282SShri Abhyankar } 3452b4ed282SShri Abhyankar /* upper bounded variables */ 3462b4ed282SShri Abhyankar else if (l[i] <= PETSC_VI_NINF) { 3472b4ed282SShri Abhyankar if (db[i] >= 1) { 3482b4ed282SShri Abhyankar t2 = PetscScalarNorm(1,t[i]); 3492b4ed282SShri Abhyankar da[i] = -1/t2 -1; 3502b4ed282SShri Abhyankar db[i] = -t[i]/t2 - 1; 3512b4ed282SShri Abhyankar } 3522b4ed282SShri Abhyankar else { 3532b4ed282SShri Abhyankar t1 = u[i] - x[i]; 3542b4ed282SShri Abhyankar t2 = PetscScalarNorm(t1,f[i]); 3552b4ed282SShri Abhyankar da[i] = t1/t2 - 1; 3562b4ed282SShri Abhyankar db[i] = -f[i]/t2 - 1; 3572b4ed282SShri Abhyankar } 3582b4ed282SShri Abhyankar } 3592b4ed282SShri Abhyankar /* Fixed variables */ 3602b4ed282SShri Abhyankar else if (l[i] == u[i]) { 3612b4ed282SShri Abhyankar da[i] = -1; 3622b4ed282SShri Abhyankar db[i] = 0; 3632b4ed282SShri Abhyankar } 3642b4ed282SShri Abhyankar /* Box constrained variables */ 3652b4ed282SShri Abhyankar else { 3662b4ed282SShri Abhyankar if (db[i] >= 1) { 3672b4ed282SShri Abhyankar t2 = PetscScalarNorm(1,t[i]); 3682b4ed282SShri Abhyankar ci = 1/t2 + 1; 3692b4ed282SShri Abhyankar di = t[i]/t2 + 1; 3702b4ed282SShri Abhyankar } 3712b4ed282SShri Abhyankar else { 3722b4ed282SShri Abhyankar t1 = x[i] - u[i]; 3732b4ed282SShri Abhyankar t2 = PetscScalarNorm(t1,f[i]); 3742b4ed282SShri Abhyankar ci = t1/t2 + 1; 3752b4ed282SShri Abhyankar di = f[i]/t2 + 1; 3762b4ed282SShri Abhyankar } 3772b4ed282SShri Abhyankar 3782b4ed282SShri Abhyankar if (da[i] >= 1) { 3792b4ed282SShri Abhyankar t1 = ci + di*t[i]; 3802b4ed282SShri Abhyankar t2 = PetscScalarNorm(1,t1); 3812b4ed282SShri Abhyankar t1 = t1/t2 - 1; 3822b4ed282SShri Abhyankar t2 = 1/t2 - 1; 3832b4ed282SShri Abhyankar } 3842b4ed282SShri Abhyankar else { 3852b4ed282SShri Abhyankar ierr = ComputeFischerFunction(u[i]-x[i],-f[i],&ei);CHKERRQ(ierr); 3862b4ed282SShri Abhyankar t2 = PetscScalarNorm(x[i]-l[i],ei); 3872b4ed282SShri Abhyankar t1 = ei/t2 - 1; 3882b4ed282SShri Abhyankar t2 = (x[i] - l[i])/t2 - 1; 3892b4ed282SShri Abhyankar } 3902b4ed282SShri Abhyankar 3912b4ed282SShri Abhyankar da[i] = t2 + t1*ci; 3922b4ed282SShri Abhyankar db[i] = t1*di; 3932b4ed282SShri Abhyankar } 3942b4ed282SShri Abhyankar } 3952b4ed282SShri Abhyankar ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 3962b4ed282SShri Abhyankar ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 3972b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->xl,&l);CHKERRQ(ierr); 3982b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->xu,&u);CHKERRQ(ierr); 399a79edbefSShri Abhyankar ierr = VecRestoreArray(Da,&da);CHKERRQ(ierr); 400a79edbefSShri Abhyankar ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr); 4012b4ed282SShri Abhyankar ierr = VecRestoreArray(vi->t,&t);CHKERRQ(ierr); 4022b4ed282SShri Abhyankar 403a79edbefSShri Abhyankar PetscFunctionReturn(0); 404a79edbefSShri Abhyankar } 405a79edbefSShri Abhyankar 406a79edbefSShri Abhyankar /* 407a79edbefSShri 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. 408a79edbefSShri Abhyankar 409a79edbefSShri Abhyankar Input Parameters: 410a79edbefSShri Abhyankar . Da - Diagonal shift vector for the semismooth jacobian. 411a79edbefSShri Abhyankar . Db - Row scaling vector for the semismooth jacobian. 412a79edbefSShri Abhyankar 413a79edbefSShri Abhyankar Output Parameters: 414a79edbefSShri Abhyankar . jac - semismooth jacobian 415a79edbefSShri Abhyankar . jac_pre - optional preconditioning matrix 416a79edbefSShri Abhyankar 417a79edbefSShri Abhyankar Notes: 418a79edbefSShri Abhyankar The semismooth jacobian matrix is given by 419a79edbefSShri Abhyankar jac = Da + Db*jacfun 420a79edbefSShri Abhyankar where Db is the row scaling matrix stored as a vector, 421a79edbefSShri Abhyankar Da is the diagonal perturbation matrix stored as a vector 422a79edbefSShri Abhyankar and jacfun is the jacobian of the original nonlinear function. 423a79edbefSShri Abhyankar */ 424a79edbefSShri Abhyankar #undef __FUNCT__ 425a79edbefSShri Abhyankar #define __FUNCT__ "SNESVIComputeJacobian" 426a79edbefSShri Abhyankar PetscErrorCode SNESVIComputeJacobian(Mat jac, Mat jac_pre,Vec Da, Vec Db) 427a79edbefSShri Abhyankar { 428a79edbefSShri Abhyankar PetscErrorCode ierr; 429a79edbefSShri Abhyankar 4302b4ed282SShri Abhyankar /* Do row scaling and add diagonal perturbation */ 431a79edbefSShri Abhyankar ierr = MatDiagonalScale(jac,Db,PETSC_NULL);CHKERRQ(ierr); 432a79edbefSShri Abhyankar ierr = MatDiagonalSet(jac,Da,ADD_VALUES);CHKERRQ(ierr); 433a79edbefSShri Abhyankar if (jac != jac_pre) { /* If jac and jac_pre are different */ 434a79edbefSShri Abhyankar ierr = MatDiagonalScale(jac_pre,Db,PETSC_NULL); 435a79edbefSShri Abhyankar ierr = MatDiagonalSet(jac_pre,Da,ADD_VALUES);CHKERRQ(ierr); 4362b4ed282SShri Abhyankar } 4372b4ed282SShri Abhyankar PetscFunctionReturn(0); 4382b4ed282SShri Abhyankar } 4392b4ed282SShri Abhyankar 4402b4ed282SShri Abhyankar /* 4412b4ed282SShri Abhyankar SNESVIAdjustInitialGuess - Readjusts the initial guess to the SNES solver supplied by the user so that the initial guess lies inside the feasible region . 4422b4ed282SShri Abhyankar 4432b4ed282SShri Abhyankar Input Parameters: 4442b4ed282SShri Abhyankar . lb - lower bound. 4452b4ed282SShri Abhyankar . ub - upper bound. 4462b4ed282SShri Abhyankar 4472b4ed282SShri Abhyankar Output Parameters: 4482b4ed282SShri Abhyankar . X - the readjusted initial guess. 4492b4ed282SShri Abhyankar 4502b4ed282SShri Abhyankar Notes: 4512b4ed282SShri Abhyankar The readjusted initial guess X[i] = max(max(min(l[i],X[i]),min(X[i],u[i])),min(l[i],u[i])) 4522b4ed282SShri Abhyankar */ 4532b4ed282SShri Abhyankar #undef __FUNCT__ 4542b4ed282SShri Abhyankar #define __FUNCT__ "SNESVIAdjustInitialGuess" 4552b4ed282SShri Abhyankar PetscErrorCode SNESVIAdjustInitialGuess(Vec X, Vec lb, Vec ub) 4562b4ed282SShri Abhyankar { 4572b4ed282SShri Abhyankar PetscErrorCode ierr; 4582b4ed282SShri Abhyankar PetscInt i,nlocal; 4592b4ed282SShri Abhyankar PetscScalar *x,*l,*u; 4602b4ed282SShri Abhyankar 4612b4ed282SShri Abhyankar PetscFunctionBegin; 4622b4ed282SShri Abhyankar 4632b4ed282SShri Abhyankar ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 4642b4ed282SShri Abhyankar ierr = VecGetArray(X,&x);CHKERRQ(ierr); 4652b4ed282SShri Abhyankar ierr = VecGetArray(lb,&l);CHKERRQ(ierr); 4662b4ed282SShri Abhyankar ierr = VecGetArray(ub,&u);CHKERRQ(ierr); 4672b4ed282SShri Abhyankar 4682b4ed282SShri Abhyankar for(i = 0; i < nlocal; i++) { 4692b4ed282SShri Abhyankar x[i] = PetscMax(PetscMax(PetscMin(x[i],l[i]),PetscMin(x[i],u[i])),PetscMin(l[i],u[i])); 4702b4ed282SShri Abhyankar } 4712b4ed282SShri Abhyankar 4722b4ed282SShri Abhyankar ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 4732b4ed282SShri Abhyankar ierr = VecRestoreArray(lb,&l);CHKERRQ(ierr); 4742b4ed282SShri Abhyankar ierr = VecRestoreArray(ub,&u);CHKERRQ(ierr); 4752b4ed282SShri Abhyankar 4762b4ed282SShri Abhyankar PetscFunctionReturn(0); 4772b4ed282SShri Abhyankar } 4782b4ed282SShri Abhyankar 4792b4ed282SShri Abhyankar /* -------------------------------------------------------------------- 4802b4ed282SShri Abhyankar 4812b4ed282SShri Abhyankar This file implements a semismooth truncated Newton method with a line search, 4822b4ed282SShri Abhyankar for solving a system of nonlinear equations in complementarity form, using the KSP, Vec, 4832b4ed282SShri Abhyankar and Mat interfaces for linear solvers, vectors, and matrices, 4842b4ed282SShri Abhyankar respectively. 4852b4ed282SShri Abhyankar 4862b4ed282SShri Abhyankar The following basic routines are required for each nonlinear solver: 4872b4ed282SShri Abhyankar SNESCreate_XXX() - Creates a nonlinear solver context 4882b4ed282SShri Abhyankar SNESSetFromOptions_XXX() - Sets runtime options 4892b4ed282SShri Abhyankar SNESSolve_XXX() - Solves the nonlinear system 4902b4ed282SShri Abhyankar SNESDestroy_XXX() - Destroys the nonlinear solver context 4912b4ed282SShri Abhyankar The suffix "_XXX" denotes a particular implementation, in this case 4922b4ed282SShri Abhyankar we use _VI (e.g., SNESCreate_VI, SNESSolve_VI) for solving 4932b4ed282SShri Abhyankar systems of nonlinear equations with a line search (LS) method. 4942b4ed282SShri Abhyankar These routines are actually called via the common user interface 4952b4ed282SShri Abhyankar routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 4962b4ed282SShri Abhyankar SNESDestroy(), so the application code interface remains identical 4972b4ed282SShri Abhyankar for all nonlinear solvers. 4982b4ed282SShri Abhyankar 4992b4ed282SShri Abhyankar Another key routine is: 5002b4ed282SShri Abhyankar SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 5012b4ed282SShri Abhyankar by setting data structures and options. The interface routine SNESSetUp() 5022b4ed282SShri Abhyankar is not usually called directly by the user, but instead is called by 5032b4ed282SShri Abhyankar SNESSolve() if necessary. 5042b4ed282SShri Abhyankar 5052b4ed282SShri Abhyankar Additional basic routines are: 5062b4ed282SShri Abhyankar SNESView_XXX() - Prints details of runtime options that 5072b4ed282SShri Abhyankar have actually been used. 5082b4ed282SShri Abhyankar These are called by application codes via the interface routines 5092b4ed282SShri Abhyankar SNESView(). 5102b4ed282SShri Abhyankar 5112b4ed282SShri Abhyankar The various types of solvers (preconditioners, Krylov subspace methods, 5122b4ed282SShri Abhyankar nonlinear solvers, timesteppers) are all organized similarly, so the 5132b4ed282SShri Abhyankar above description applies to these categories also. 5142b4ed282SShri Abhyankar 5152b4ed282SShri Abhyankar -------------------------------------------------------------------- */ 5162b4ed282SShri Abhyankar /* 51769c03620SShri Abhyankar SNESSolveVI_SS - Solves the complementarity problem with a semismooth Newton 5182b4ed282SShri Abhyankar method using a line search. 5192b4ed282SShri Abhyankar 5202b4ed282SShri Abhyankar Input Parameters: 5212b4ed282SShri Abhyankar . snes - the SNES context 5222b4ed282SShri Abhyankar 5232b4ed282SShri Abhyankar Output Parameter: 5242b4ed282SShri Abhyankar . outits - number of iterations until termination 5252b4ed282SShri Abhyankar 5262b4ed282SShri Abhyankar Application Interface Routine: SNESSolve() 5272b4ed282SShri Abhyankar 5282b4ed282SShri Abhyankar Notes: 5292b4ed282SShri Abhyankar This implements essentially a semismooth Newton method with a 5302b4ed282SShri Abhyankar line search. By default a cubic backtracking line search 5312b4ed282SShri Abhyankar is employed, as described in the text "Numerical Methods for 5322b4ed282SShri Abhyankar Unconstrained Optimization and Nonlinear Equations" by Dennis 5332b4ed282SShri Abhyankar and Schnabel. 5342b4ed282SShri Abhyankar */ 5352b4ed282SShri Abhyankar #undef __FUNCT__ 53669c03620SShri Abhyankar #define __FUNCT__ "SNESSolveVI_SS" 53769c03620SShri Abhyankar PetscErrorCode SNESSolveVI_SS(SNES snes) 5382b4ed282SShri Abhyankar { 5392b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 5402b4ed282SShri Abhyankar PetscErrorCode ierr; 5412b4ed282SShri Abhyankar PetscInt maxits,i,lits; 5423b336b49SShri Abhyankar PetscBool lssucceed; 5432b4ed282SShri Abhyankar MatStructure flg = DIFFERENT_NONZERO_PATTERN; 5442b4ed282SShri Abhyankar PetscReal gnorm,xnorm=0,ynorm; 5452b4ed282SShri Abhyankar Vec Y,X,F,G,W; 5462b4ed282SShri Abhyankar KSPConvergedReason kspreason; 5472b4ed282SShri Abhyankar 5482b4ed282SShri Abhyankar PetscFunctionBegin; 5492b4ed282SShri Abhyankar snes->numFailures = 0; 5502b4ed282SShri Abhyankar snes->numLinearSolveFailures = 0; 5512b4ed282SShri Abhyankar snes->reason = SNES_CONVERGED_ITERATING; 5522b4ed282SShri Abhyankar 5532b4ed282SShri Abhyankar maxits = snes->max_its; /* maximum number of iterations */ 5542b4ed282SShri Abhyankar X = snes->vec_sol; /* solution vector */ 5552b4ed282SShri Abhyankar F = snes->vec_func; /* residual vector */ 5562b4ed282SShri Abhyankar Y = snes->work[0]; /* work vectors */ 5572b4ed282SShri Abhyankar G = snes->work[1]; 5582b4ed282SShri Abhyankar W = snes->work[2]; 5592b4ed282SShri Abhyankar 5602b4ed282SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 5612b4ed282SShri Abhyankar snes->iter = 0; 5622b4ed282SShri Abhyankar snes->norm = 0.0; 5632b4ed282SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 5642b4ed282SShri Abhyankar 5652b4ed282SShri Abhyankar ierr = SNESVIAdjustInitialGuess(X,vi->xl,vi->xu);CHKERRQ(ierr); 5662b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr); 5672b4ed282SShri Abhyankar if (snes->domainerror) { 5682b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 5692b4ed282SShri Abhyankar PetscFunctionReturn(0); 5702b4ed282SShri Abhyankar } 5712b4ed282SShri Abhyankar /* Compute Merit function */ 5722b4ed282SShri Abhyankar ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr); 5732b4ed282SShri Abhyankar 5742b4ed282SShri Abhyankar ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 5752b4ed282SShri Abhyankar ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 5762b4ed282SShri Abhyankar if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 5772b4ed282SShri Abhyankar 5782b4ed282SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 5792b4ed282SShri Abhyankar snes->norm = vi->phinorm; 5802b4ed282SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 5812b4ed282SShri Abhyankar SNESLogConvHistory(snes,vi->phinorm,0); 5822b4ed282SShri Abhyankar SNESMonitor(snes,0,vi->phinorm); 5832b4ed282SShri Abhyankar 5842b4ed282SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 5852b4ed282SShri Abhyankar snes->ttol = vi->phinorm*snes->rtol; 5862b4ed282SShri Abhyankar /* test convergence */ 5872b4ed282SShri Abhyankar ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 5882b4ed282SShri Abhyankar if (snes->reason) PetscFunctionReturn(0); 5892b4ed282SShri Abhyankar 5902b4ed282SShri Abhyankar for (i=0; i<maxits; i++) { 5912b4ed282SShri Abhyankar 5922b4ed282SShri Abhyankar /* Call general purpose update function */ 5932b4ed282SShri Abhyankar if (snes->ops->update) { 5942b4ed282SShri Abhyankar ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 5952b4ed282SShri Abhyankar } 5962b4ed282SShri Abhyankar 5972b4ed282SShri Abhyankar /* Solve J Y = Phi, where J is the semismooth jacobian */ 598a79edbefSShri Abhyankar /* Get the nonlinear function jacobian */ 5992b4ed282SShri Abhyankar ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 600a79edbefSShri Abhyankar /* Get the diagonal shift and row scaling vectors */ 601a79edbefSShri Abhyankar ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr); 602a79edbefSShri Abhyankar /* Compute the semismooth jacobian */ 603a79edbefSShri Abhyankar ierr = SNESVIComputeJacobian(snes->jacobian,snes->jacobian_pre,vi->Da,vi->Db);CHKERRQ(ierr); 6042b4ed282SShri Abhyankar 6052b4ed282SShri Abhyankar ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 6062b4ed282SShri Abhyankar ierr = SNES_KSPSolve(snes,snes->ksp,vi->phi,Y);CHKERRQ(ierr); 6072b4ed282SShri Abhyankar ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 6083b336b49SShri Abhyankar 6093b336b49SShri Abhyankar if (kspreason < 0) { 6102b4ed282SShri Abhyankar if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 6112b4ed282SShri 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); 6122b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 6132b4ed282SShri Abhyankar break; 6142b4ed282SShri Abhyankar } 6152b4ed282SShri Abhyankar } 6162b4ed282SShri Abhyankar ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 6172b4ed282SShri Abhyankar snes->linear_its += lits; 6182b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 6192b4ed282SShri Abhyankar /* 6202b4ed282SShri Abhyankar if (vi->precheckstep) { 621ace3abfcSBarry Smith PetscBool changed_y = PETSC_FALSE; 6222b4ed282SShri Abhyankar ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr); 6232b4ed282SShri Abhyankar } 6242b4ed282SShri Abhyankar 6252b4ed282SShri Abhyankar if (PetscLogPrintInfo){ 6262b4ed282SShri Abhyankar ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 6272b4ed282SShri Abhyankar } 6282b4ed282SShri Abhyankar */ 6292b4ed282SShri Abhyankar /* Compute a (scaled) negative update in the line search routine: 6302b4ed282SShri Abhyankar Y <- X - lambda*Y 6312b4ed282SShri Abhyankar and evaluate G = function(Y) (depends on the line search). 6322b4ed282SShri Abhyankar */ 6332b4ed282SShri Abhyankar ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 6342b4ed282SShri Abhyankar ynorm = 1; gnorm = vi->phinorm; 6352b4ed282SShri Abhyankar ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 6362b4ed282SShri Abhyankar ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 6372b4ed282SShri Abhyankar if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 6382b4ed282SShri Abhyankar if (snes->domainerror) { 6392b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 6402b4ed282SShri Abhyankar PetscFunctionReturn(0); 6412b4ed282SShri Abhyankar } 6422b4ed282SShri Abhyankar if (!lssucceed) { 6432b4ed282SShri Abhyankar if (++snes->numFailures >= snes->maxFailures) { 644ace3abfcSBarry Smith PetscBool ismin; 6452b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_LINE_SEARCH; 6462b4ed282SShri Abhyankar ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 6472b4ed282SShri Abhyankar if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 6482b4ed282SShri Abhyankar break; 6492b4ed282SShri Abhyankar } 6502b4ed282SShri Abhyankar } 6512b4ed282SShri Abhyankar /* Update function and solution vectors */ 6522b4ed282SShri Abhyankar vi->phinorm = gnorm; 6532b4ed282SShri Abhyankar vi->merit = 0.5*vi->phinorm*vi->phinorm; 6542b4ed282SShri Abhyankar ierr = VecCopy(G,vi->phi);CHKERRQ(ierr); 6552b4ed282SShri Abhyankar ierr = VecCopy(W,X);CHKERRQ(ierr); 6562b4ed282SShri Abhyankar /* Monitor convergence */ 6572b4ed282SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 6582b4ed282SShri Abhyankar snes->iter = i+1; 6592b4ed282SShri Abhyankar snes->norm = vi->phinorm; 6602b4ed282SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 6612b4ed282SShri Abhyankar SNESLogConvHistory(snes,snes->norm,lits); 6622b4ed282SShri Abhyankar SNESMonitor(snes,snes->iter,snes->norm); 6632b4ed282SShri Abhyankar /* Test for convergence, xnorm = || X || */ 6642b4ed282SShri Abhyankar if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 6652b4ed282SShri Abhyankar ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 6662b4ed282SShri Abhyankar if (snes->reason) break; 6672b4ed282SShri Abhyankar } 6682b4ed282SShri Abhyankar if (i == maxits) { 6692b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 6702b4ed282SShri Abhyankar if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 6712b4ed282SShri Abhyankar } 6722b4ed282SShri Abhyankar PetscFunctionReturn(0); 6732b4ed282SShri Abhyankar } 67469c03620SShri Abhyankar 67569c03620SShri Abhyankar #undef __FUNCT__ 67669c03620SShri Abhyankar #define __FUNCT__ "SNESVICreateIndexSets_AS" 677a79edbefSShri Abhyankar PetscErrorCode SNESVICreateIndexSets_AS(SNES snes,Vec Db,PetscReal thresh,IS* ISact,IS* ISinact) 67869c03620SShri Abhyankar { 67969c03620SShri Abhyankar PetscErrorCode ierr; 68069c03620SShri Abhyankar PetscInt i,nlocal,ilow,ihigh,nloc_isact=0,nloc_isinact=0; 68169c03620SShri Abhyankar PetscInt *idx_act,*idx_inact,i1=0,i2=0; 68269c03620SShri Abhyankar PetscScalar *db; 68369c03620SShri Abhyankar 68469c03620SShri Abhyankar PetscFunctionBegin; 68569c03620SShri Abhyankar 68669c03620SShri Abhyankar ierr = VecGetLocalSize(Db,&nlocal);CHKERRQ(ierr); 68769c03620SShri Abhyankar ierr = VecGetOwnershipRange(Db,&ilow,&ihigh);CHKERRQ(ierr); 68869c03620SShri Abhyankar ierr = VecGetArray(Db,&db);CHKERRQ(ierr); 68969c03620SShri Abhyankar /* Compute the sizes of the active and inactive sets */ 69069c03620SShri Abhyankar for (i=0; i < nlocal;i++) { 691fca35906SShri Abhyankar if (PetscAbsScalar(db[i]) <= thresh) nloc_isact++; 692fca35906SShri Abhyankar else nloc_isinact++; 69369c03620SShri Abhyankar } 69469c03620SShri Abhyankar ierr = PetscMalloc(nloc_isact*sizeof(PetscInt),&idx_act);CHKERRQ(ierr); 69569c03620SShri Abhyankar ierr = PetscMalloc(nloc_isinact*sizeof(PetscInt),&idx_inact);CHKERRQ(ierr); 69669c03620SShri Abhyankar 69769c03620SShri Abhyankar /* Creating the indexing arrays */ 698fca35906SShri Abhyankar for(i=0; i < nlocal; i++) { 699fca35906SShri Abhyankar if (PetscAbsScalar(db[i]) <= thresh) idx_act[i1++] = ilow+i; 700fca35906SShri Abhyankar else idx_inact[i2++] = ilow+i; 70169c03620SShri Abhyankar } 70269c03620SShri Abhyankar 70369c03620SShri Abhyankar /* Create the index sets */ 70469c03620SShri Abhyankar ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isact,idx_act,PETSC_COPY_VALUES,ISact);CHKERRQ(ierr); 70569c03620SShri Abhyankar ierr = ISCreateGeneral(((PetscObject)snes)->comm,nloc_isinact,idx_inact,PETSC_COPY_VALUES,ISinact);CHKERRQ(ierr); 70669c03620SShri Abhyankar 70769c03620SShri Abhyankar ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr); 70869c03620SShri Abhyankar ierr = PetscFree(idx_act);CHKERRQ(ierr); 70969c03620SShri Abhyankar ierr = PetscFree(idx_inact);CHKERRQ(ierr); 71069c03620SShri Abhyankar PetscFunctionReturn(0); 71169c03620SShri Abhyankar } 71269c03620SShri Abhyankar 713dbd914b8SShri Abhyankar /* Create active and inactive set vectors. The local size of this vector is set and petsc computes the global size */ 714dbd914b8SShri Abhyankar #undef __FUNCT__ 715dbd914b8SShri Abhyankar #define __FUNCT__ "SNESVICreateVectors_AS" 716dbd914b8SShri Abhyankar PetscErrorCode SNESVICreateVectors_AS(SNES snes,PetscInt n,Vec* newv) 717dbd914b8SShri Abhyankar { 718dbd914b8SShri Abhyankar PetscErrorCode ierr; 719dbd914b8SShri Abhyankar Vec v; 720dbd914b8SShri Abhyankar 721dbd914b8SShri Abhyankar PetscFunctionBegin; 722dbd914b8SShri Abhyankar ierr = VecCreate(((PetscObject)snes)->comm,&v);CHKERRQ(ierr); 723dbd914b8SShri Abhyankar ierr = VecSetSizes(v,n,PETSC_DECIDE);CHKERRQ(ierr); 724dbd914b8SShri Abhyankar ierr = VecSetFromOptions(v);CHKERRQ(ierr); 725dbd914b8SShri Abhyankar *newv = v; 726dbd914b8SShri Abhyankar 727dbd914b8SShri Abhyankar PetscFunctionReturn(0); 728dbd914b8SShri Abhyankar } 729dbd914b8SShri Abhyankar 730dbd914b8SShri Abhyankar 73169c03620SShri Abhyankar /* Variational Inequality solver using active set method */ 73269c03620SShri Abhyankar #undef __FUNCT__ 73369c03620SShri Abhyankar #define __FUNCT__ "SNESSolveVI_AS" 73469c03620SShri Abhyankar PetscErrorCode SNESSolveVI_AS(SNES snes) 73569c03620SShri Abhyankar { 73669c03620SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 73769c03620SShri Abhyankar PetscErrorCode ierr; 738408da460SShri Abhyankar PetscInt maxits,i,lits,Nis_act=0; 7393b336b49SShri Abhyankar PetscBool lssucceed; 74069c03620SShri Abhyankar MatStructure flg = DIFFERENT_NONZERO_PATTERN; 74169c03620SShri Abhyankar PetscReal gnorm,xnorm=0,ynorm; 74269c03620SShri Abhyankar Vec Y,X,F,G,W; 74369c03620SShri Abhyankar KSPConvergedReason kspreason; 74469c03620SShri Abhyankar 74569c03620SShri Abhyankar PetscFunctionBegin; 74669c03620SShri Abhyankar snes->numFailures = 0; 74769c03620SShri Abhyankar snes->numLinearSolveFailures = 0; 74869c03620SShri Abhyankar snes->reason = SNES_CONVERGED_ITERATING; 74969c03620SShri Abhyankar 75069c03620SShri Abhyankar maxits = snes->max_its; /* maximum number of iterations */ 75169c03620SShri Abhyankar X = snes->vec_sol; /* solution vector */ 75269c03620SShri Abhyankar F = snes->vec_func; /* residual vector */ 75369c03620SShri Abhyankar Y = snes->work[0]; /* work vectors */ 75469c03620SShri Abhyankar G = snes->work[1]; 75569c03620SShri Abhyankar W = snes->work[2]; 75669c03620SShri Abhyankar 75769c03620SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 75869c03620SShri Abhyankar snes->iter = 0; 75969c03620SShri Abhyankar snes->norm = 0.0; 76069c03620SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 76169c03620SShri Abhyankar 76269c03620SShri Abhyankar ierr = SNESVIAdjustInitialGuess(X,vi->xl,vi->xu);CHKERRQ(ierr); 76369c03620SShri Abhyankar ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr); 76469c03620SShri Abhyankar if (snes->domainerror) { 76569c03620SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 76669c03620SShri Abhyankar PetscFunctionReturn(0); 76769c03620SShri Abhyankar } 76869c03620SShri Abhyankar /* Compute Merit function */ 76969c03620SShri Abhyankar ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr); 77069c03620SShri Abhyankar 77169c03620SShri Abhyankar ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 77269c03620SShri Abhyankar ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 77369c03620SShri Abhyankar if PetscIsInfOrNanReal(vi->merit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 77469c03620SShri Abhyankar 77569c03620SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 77669c03620SShri Abhyankar snes->norm = vi->phinorm; 77769c03620SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 77869c03620SShri Abhyankar SNESLogConvHistory(snes,vi->phinorm,0); 77969c03620SShri Abhyankar SNESMonitor(snes,0,vi->phinorm); 78069c03620SShri Abhyankar 78169c03620SShri Abhyankar /* set parameter for default relative tolerance convergence test */ 78269c03620SShri Abhyankar snes->ttol = vi->phinorm*snes->rtol; 78369c03620SShri Abhyankar /* test convergence */ 78469c03620SShri Abhyankar ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 78569c03620SShri Abhyankar if (snes->reason) PetscFunctionReturn(0); 78669c03620SShri Abhyankar 78769c03620SShri Abhyankar for (i=0; i<maxits; i++) { 78869c03620SShri Abhyankar 789a79edbefSShri Abhyankar IS IS_act,IS_inact; /* _act -> active set _inact -> inactive set */ 790a79edbefSShri Abhyankar PetscReal thresh,J_norm1; 791dbd914b8SShri Abhyankar VecScatter scat_act,scat_inact; 792408da460SShri Abhyankar PetscInt nis_act,nis_inact,Nis_act_prev; 793dbd914b8SShri Abhyankar Vec Da_act,Da_inact,Db_inact; 794dbd914b8SShri Abhyankar Vec Y_act,Y_inact,phi_act,phi_inact; 795d76c674bSShri Abhyankar Mat jac_inact_inact,jac_inact_act,prejac_inact_inact; 796a79edbefSShri Abhyankar 79769c03620SShri Abhyankar /* Call general purpose update function */ 79869c03620SShri Abhyankar if (snes->ops->update) { 79969c03620SShri Abhyankar ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 80069c03620SShri Abhyankar } 80169c03620SShri Abhyankar ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 80269c03620SShri Abhyankar /* Compute the threshold value for creating active and inactive sets */ 80369c03620SShri Abhyankar ierr = MatNorm(snes->jacobian,NORM_1,&J_norm1);CHKERRQ(ierr); 80469c03620SShri Abhyankar thresh = PetscMin(vi->merit,1e-2)/(1+J_norm1); 805a79edbefSShri Abhyankar 806a79edbefSShri Abhyankar /* Compute B-subdifferential vectors Da and Db */ 807a79edbefSShri Abhyankar ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr); 808dbd914b8SShri Abhyankar 80969c03620SShri Abhyankar /* Create active and inactive index sets */ 81069c03620SShri Abhyankar ierr = SNESVICreateIndexSets_AS(snes,vi->Db,thresh,&IS_act,&IS_inact);CHKERRQ(ierr); 81169c03620SShri Abhyankar 812408da460SShri Abhyankar Nis_act_prev = Nis_act; 813408da460SShri Abhyankar /* Get sizes of active and inactive sets */ 814dbd914b8SShri Abhyankar ierr = ISGetLocalSize(IS_act,&nis_act);CHKERRQ(ierr); 815dbd914b8SShri Abhyankar ierr = ISGetLocalSize(IS_inact,&nis_inact);CHKERRQ(ierr); 816408da460SShri Abhyankar ierr = ISGetSize(IS_act,&Nis_act);CHKERRQ(ierr); 817dbd914b8SShri Abhyankar 818d76c674bSShri Abhyankar ierr = PetscPrintf(PETSC_COMM_WORLD,"Size of active set = %d, size of inactive set = %d\n",nis_act,nis_inact);CHKERRQ(ierr); 819d76c674bSShri Abhyankar 820dbd914b8SShri Abhyankar /* Create active and inactive set vectors */ 821dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_act,&Da_act);CHKERRQ(ierr); 822dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_inact,&Da_inact);CHKERRQ(ierr); 823dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_inact,&Db_inact);CHKERRQ(ierr); 824dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_act,&phi_act);CHKERRQ(ierr); 825dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_inact,&phi_inact);CHKERRQ(ierr); 826dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_act,&Y_act);CHKERRQ(ierr); 827dbd914b8SShri Abhyankar ierr = SNESVICreateVectors_AS(snes,nis_inact,&Y_inact);CHKERRQ(ierr); 828dbd914b8SShri Abhyankar 829dbd914b8SShri Abhyankar /* Create inactive set submatrices */ 830dbd914b8SShri Abhyankar ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_act,MAT_INITIAL_MATRIX,&jac_inact_act);CHKERRQ(ierr); 831dbd914b8SShri Abhyankar ierr = MatGetSubMatrix(snes->jacobian,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&jac_inact_inact);CHKERRQ(ierr); 832dbd914b8SShri Abhyankar 833dbd914b8SShri Abhyankar /* Create scatter contexts */ 834dbd914b8SShri Abhyankar ierr = VecScatterCreate(vi->Da,IS_act,Da_act,PETSC_NULL,&scat_act);CHKERRQ(ierr); 835dbd914b8SShri Abhyankar ierr = VecScatterCreate(vi->Da,IS_inact,Da_inact,PETSC_NULL,&scat_inact);CHKERRQ(ierr); 836dbd914b8SShri Abhyankar 837dbd914b8SShri Abhyankar /* Do a vec scatter to active and inactive set vectors */ 838dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_act,vi->Da,Da_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 839dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_act,vi->Da,Da_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 840dbd914b8SShri Abhyankar 841dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_inact,vi->Da,Da_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 842dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_inact,vi->Da,Da_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 843dbd914b8SShri Abhyankar 844dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_inact,vi->Db,Db_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 845dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_inact,vi->Db,Db_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 846dbd914b8SShri Abhyankar 847dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_act,vi->phi,phi_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 848dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_act,vi->phi,phi_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 849dbd914b8SShri Abhyankar 850dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_inact,vi->phi,phi_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 851dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_inact,vi->phi,phi_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 852dbd914b8SShri Abhyankar 853dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 854dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_act,Y,Y_act,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 855dbd914b8SShri Abhyankar 856dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 857dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_inact,Y,Y_inact,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 858dbd914b8SShri Abhyankar 859dbd914b8SShri Abhyankar /* Active set direction */ 860dbd914b8SShri Abhyankar ierr = VecPointwiseDivide(Y_act,phi_act,Da_act);CHKERRQ(ierr); 861d76c674bSShri Abhyankar /* inactive set jacobian and preconditioner */ 862dbd914b8SShri Abhyankar ierr = VecPointwiseDivide(Da_inact,Da_inact,Db_inact);CHKERRQ(ierr); 863dbd914b8SShri Abhyankar ierr = MatDiagonalSet(jac_inact_inact,Da_inact,ADD_VALUES);CHKERRQ(ierr); 864d76c674bSShri Abhyankar if (snes->jacobian != snes->jacobian_pre) { 865d76c674bSShri Abhyankar ierr = MatGetSubMatrix(snes->jacobian_pre,IS_inact,IS_inact,MAT_INITIAL_MATRIX,&prejac_inact_inact);CHKERRQ(ierr); 866d76c674bSShri Abhyankar ierr = MatDiagonalSet(prejac_inact_inact,Da_inact,ADD_VALUES);CHKERRQ(ierr); 867d76c674bSShri Abhyankar } else prejac_inact_inact = jac_inact_inact; 868d76c674bSShri Abhyankar 869dbd914b8SShri Abhyankar /* right hand side */ 870dbd914b8SShri Abhyankar ierr = VecPointwiseDivide(phi_inact,phi_inact,Db_inact);CHKERRQ(ierr); 871dbd914b8SShri Abhyankar ierr = MatMult(jac_inact_act,Y_act,Db_inact);CHKERRQ(ierr); 872dbd914b8SShri Abhyankar ierr = VecAXPY(phi_inact,-1.0,Db_inact);CHKERRQ(ierr); 873dbd914b8SShri Abhyankar 874408da460SShri Abhyankar if ((i != 0) && (Nis_act != Nis_act_prev)) { 875408da460SShri Abhyankar KSP kspnew,snesksp; 876408da460SShri Abhyankar PC pcnew; 877408da460SShri Abhyankar /* The active and inactive set sizes have changed so need to create a new snes->ksp object */ 878408da460SShri Abhyankar ierr = SNESGetKSP(snes,&snesksp);CHKERRQ(ierr); 879408da460SShri Abhyankar ierr = KSPCreate(((PetscObject)snes)->comm,&kspnew);CHKERRQ(ierr); 880408da460SShri Abhyankar /* Copy over snes->ksp info */ 881408da460SShri Abhyankar kspnew->pc_side = snesksp->pc_side; 882408da460SShri Abhyankar kspnew->rtol = snesksp->rtol; 883408da460SShri Abhyankar kspnew->abstol = snesksp->abstol; 884408da460SShri Abhyankar kspnew->max_it = snesksp->max_it; 885408da460SShri Abhyankar ierr = KSPSetType(kspnew,((PetscObject)snesksp)->type_name);CHKERRQ(ierr); 886408da460SShri Abhyankar ierr = KSPGetPC(kspnew,&pcnew);CHKERRQ(ierr); 887408da460SShri Abhyankar ierr = PCSetType(kspnew->pc,((PetscObject)snesksp->pc)->type_name);CHKERRQ(ierr); 888408da460SShri Abhyankar ierr = KSPDestroy(snesksp);CHKERRQ(ierr); 889408da460SShri Abhyankar snes->ksp = kspnew; 890408da460SShri Abhyankar ierr = PetscLogObjectParent(snes,kspnew);CHKERRQ(ierr); 891408da460SShri Abhyankar ierr = KSPSetFromOptions(kspnew);CHKERRQ(ierr); 892408da460SShri Abhyankar } 893408da460SShri Abhyankar 894d76c674bSShri Abhyankar ierr = KSPSetOperators(snes->ksp,jac_inact_inact,prejac_inact_inact,flg);CHKERRQ(ierr); 895dbd914b8SShri Abhyankar ierr = SNES_KSPSolve(snes,snes->ksp,phi_inact,Y_inact);CHKERRQ(ierr); 896dbd914b8SShri Abhyankar ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 8973b336b49SShri Abhyankar if (kspreason < 0) { 8983b336b49SShri Abhyankar if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 8993b336b49SShri 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); 9003b336b49SShri Abhyankar snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 9013b336b49SShri Abhyankar break; 9023b336b49SShri Abhyankar } 9033b336b49SShri Abhyankar } 904dbd914b8SShri Abhyankar 905dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 906dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_act,Y_act,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 907dbd914b8SShri Abhyankar ierr = VecScatterBegin(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 908dbd914b8SShri Abhyankar ierr = VecScatterEnd(scat_inact,Y_inact,Y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 909dbd914b8SShri Abhyankar 910dbd914b8SShri Abhyankar ierr = VecDestroy(Da_act);CHKERRQ(ierr); 911dbd914b8SShri Abhyankar ierr = VecDestroy(Da_inact);CHKERRQ(ierr); 912dbd914b8SShri Abhyankar ierr = VecDestroy(Db_inact);CHKERRQ(ierr); 913dbd914b8SShri Abhyankar ierr = VecDestroy(phi_act);CHKERRQ(ierr); 914dbd914b8SShri Abhyankar ierr = VecDestroy(phi_inact);CHKERRQ(ierr); 915dbd914b8SShri Abhyankar ierr = VecDestroy(Y_act);CHKERRQ(ierr); 916dbd914b8SShri Abhyankar ierr = VecDestroy(Y_inact);CHKERRQ(ierr); 917dbd914b8SShri Abhyankar ierr = VecScatterDestroy(scat_act);CHKERRQ(ierr); 918dbd914b8SShri Abhyankar ierr = VecScatterDestroy(scat_inact);CHKERRQ(ierr); 919730c24dcSShri Abhyankar ierr = ISDestroy(IS_act);CHKERRQ(ierr); 920730c24dcSShri Abhyankar ierr = ISDestroy(IS_inact);CHKERRQ(ierr); 921bac2f86eSShri Abhyankar ierr = MatDestroy(jac_inact_act);CHKERRQ(ierr); 922bac2f86eSShri Abhyankar ierr = MatDestroy(jac_inact_inact);CHKERRQ(ierr); 923d76c674bSShri Abhyankar if (snes->jacobian != snes->jacobian_pre) { 924d76c674bSShri Abhyankar ierr = MatDestroy(prejac_inact_inact);CHKERRQ(ierr); 925d76c674bSShri Abhyankar } 926730c24dcSShri Abhyankar 92769c03620SShri Abhyankar ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 92869c03620SShri Abhyankar snes->linear_its += lits; 92969c03620SShri Abhyankar ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 93069c03620SShri Abhyankar /* 93169c03620SShri Abhyankar if (vi->precheckstep) { 93269c03620SShri Abhyankar PetscBool changed_y = PETSC_FALSE; 93369c03620SShri Abhyankar ierr = (*vi->precheckstep)(snes,X,Y,vi->precheck,&changed_y);CHKERRQ(ierr); 93469c03620SShri Abhyankar } 93569c03620SShri Abhyankar 93669c03620SShri Abhyankar if (PetscLogPrintInfo){ 93769c03620SShri Abhyankar ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 93869c03620SShri Abhyankar } 93969c03620SShri Abhyankar */ 94069c03620SShri Abhyankar /* Compute a (scaled) negative update in the line search routine: 94169c03620SShri Abhyankar Y <- X - lambda*Y 94269c03620SShri Abhyankar and evaluate G = function(Y) (depends on the line search). 94369c03620SShri Abhyankar */ 94469c03620SShri Abhyankar ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 94569c03620SShri Abhyankar ynorm = 1; gnorm = vi->phinorm; 9463b336b49SShri Abhyankar ierr = SNESVIComputeJacobian(snes->jacobian,snes->jacobian_pre,vi->Da,vi->Db);CHKERRQ(ierr); 94769c03620SShri Abhyankar ierr = (*vi->LineSearch)(snes,vi->lsP,X,vi->phi,G,Y,W,vi->phinorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 94869c03620SShri Abhyankar ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",vi->phinorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 94969c03620SShri Abhyankar if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 95069c03620SShri Abhyankar if (snes->domainerror) { 95169c03620SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 95269c03620SShri Abhyankar PetscFunctionReturn(0); 95369c03620SShri Abhyankar } 95469c03620SShri Abhyankar if (!lssucceed) { 95569c03620SShri Abhyankar if (++snes->numFailures >= snes->maxFailures) { 95669c03620SShri Abhyankar PetscBool ismin; 95769c03620SShri Abhyankar snes->reason = SNES_DIVERGED_LINE_SEARCH; 95869c03620SShri Abhyankar ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 95969c03620SShri Abhyankar if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 96069c03620SShri Abhyankar break; 96169c03620SShri Abhyankar } 96269c03620SShri Abhyankar } 96369c03620SShri Abhyankar /* Update function and solution vectors */ 96469c03620SShri Abhyankar vi->phinorm = gnorm; 96569c03620SShri Abhyankar vi->merit = 0.5*vi->phinorm*vi->phinorm; 96669c03620SShri Abhyankar ierr = VecCopy(G,vi->phi);CHKERRQ(ierr); 96769c03620SShri Abhyankar ierr = VecCopy(W,X);CHKERRQ(ierr); 96869c03620SShri Abhyankar /* Monitor convergence */ 96969c03620SShri Abhyankar ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 97069c03620SShri Abhyankar snes->iter = i+1; 97169c03620SShri Abhyankar snes->norm = vi->phinorm; 97269c03620SShri Abhyankar ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 97369c03620SShri Abhyankar SNESLogConvHistory(snes,snes->norm,lits); 97469c03620SShri Abhyankar SNESMonitor(snes,snes->iter,snes->norm); 97569c03620SShri Abhyankar /* Test for convergence, xnorm = || X || */ 97669c03620SShri Abhyankar if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 97769c03620SShri Abhyankar ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 97869c03620SShri Abhyankar if (snes->reason) break; 97969c03620SShri Abhyankar } 98069c03620SShri Abhyankar if (i == maxits) { 98169c03620SShri Abhyankar ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 98269c03620SShri Abhyankar if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 98369c03620SShri Abhyankar } 98469c03620SShri Abhyankar PetscFunctionReturn(0); 98569c03620SShri Abhyankar } 98669c03620SShri Abhyankar 9872b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 9882b4ed282SShri Abhyankar /* 9892b4ed282SShri Abhyankar SNESSetUp_VI - Sets up the internal data structures for the later use 9902b4ed282SShri Abhyankar of the SNESVI nonlinear solver. 9912b4ed282SShri Abhyankar 9922b4ed282SShri Abhyankar Input Parameter: 9932b4ed282SShri Abhyankar . snes - the SNES context 9942b4ed282SShri Abhyankar . x - the solution vector 9952b4ed282SShri Abhyankar 9962b4ed282SShri Abhyankar Application Interface Routine: SNESSetUp() 9972b4ed282SShri Abhyankar 9982b4ed282SShri Abhyankar Notes: 9992b4ed282SShri Abhyankar For basic use of the SNES solvers, the user need not explicitly call 10002b4ed282SShri Abhyankar SNESSetUp(), since these actions will automatically occur during 10012b4ed282SShri Abhyankar the call to SNESSolve(). 10022b4ed282SShri Abhyankar */ 10032b4ed282SShri Abhyankar #undef __FUNCT__ 10042b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetUp_VI" 10052b4ed282SShri Abhyankar PetscErrorCode SNESSetUp_VI(SNES snes) 10062b4ed282SShri Abhyankar { 10072b4ed282SShri Abhyankar PetscErrorCode ierr; 10082b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*) snes->data; 10092b4ed282SShri Abhyankar PetscInt i_start[3],i_end[3]; 10102b4ed282SShri Abhyankar 10112b4ed282SShri Abhyankar PetscFunctionBegin; 10122b4ed282SShri Abhyankar if (!snes->vec_sol_update) { 10132b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr); 10142b4ed282SShri Abhyankar ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr); 10152b4ed282SShri Abhyankar } 10162b4ed282SShri Abhyankar if (!snes->work) { 10172b4ed282SShri Abhyankar snes->nwork = 3; 10182b4ed282SShri Abhyankar ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 10192b4ed282SShri Abhyankar ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr); 10202b4ed282SShri Abhyankar } 10212b4ed282SShri Abhyankar 10222b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->phi); CHKERRQ(ierr); 10232b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->Da); CHKERRQ(ierr); 10242b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->Db); CHKERRQ(ierr); 10252b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->z);CHKERRQ(ierr); 10262b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->t); CHKERRQ(ierr); 10272b4ed282SShri Abhyankar 10282b4ed282SShri Abhyankar /* If the lower and upper bound on variables are not set, set it to 10292b4ed282SShri Abhyankar -Inf and Inf */ 10302b4ed282SShri Abhyankar if (!vi->xl && !vi->xu) { 10312b4ed282SShri Abhyankar vi->usersetxbounds = PETSC_FALSE; 10322b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->xl); CHKERRQ(ierr); 10332b4ed282SShri Abhyankar ierr = VecSet(vi->xl,PETSC_VI_NINF);CHKERRQ(ierr); 10342b4ed282SShri Abhyankar ierr = VecDuplicate(snes->vec_sol, &vi->xu); CHKERRQ(ierr); 10352b4ed282SShri Abhyankar ierr = VecSet(vi->xu,PETSC_VI_INF);CHKERRQ(ierr); 10362b4ed282SShri Abhyankar } else { 10372b4ed282SShri Abhyankar /* Check if lower bound, upper bound and solution vector distribution across the processors is identical */ 10382b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(snes->vec_sol,i_start,i_end);CHKERRQ(ierr); 10392b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(vi->xl,i_start+1,i_end+1);CHKERRQ(ierr); 10402b4ed282SShri Abhyankar ierr = VecGetOwnershipRange(vi->xu,i_start+2,i_end+2);CHKERRQ(ierr); 10412b4ed282SShri 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])) 10422b4ed282SShri 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."); 10432b4ed282SShri Abhyankar } 10442b4ed282SShri Abhyankar 10452b4ed282SShri Abhyankar vi->computeuserfunction = snes->ops->computefunction; 10461f79c099SShri Abhyankar snes->ops->computefunction = SNESVIComputeFunction; 1047a79edbefSShri Abhyankar 10482b4ed282SShri Abhyankar PetscFunctionReturn(0); 10492b4ed282SShri Abhyankar } 10502b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 10512b4ed282SShri Abhyankar /* 10522b4ed282SShri Abhyankar SNESDestroy_VI - Destroys the private SNES_VI context that was created 10532b4ed282SShri Abhyankar with SNESCreate_VI(). 10542b4ed282SShri Abhyankar 10552b4ed282SShri Abhyankar Input Parameter: 10562b4ed282SShri Abhyankar . snes - the SNES context 10572b4ed282SShri Abhyankar 10582b4ed282SShri Abhyankar Application Interface Routine: SNESDestroy() 10592b4ed282SShri Abhyankar */ 10602b4ed282SShri Abhyankar #undef __FUNCT__ 10612b4ed282SShri Abhyankar #define __FUNCT__ "SNESDestroy_VI" 10622b4ed282SShri Abhyankar PetscErrorCode SNESDestroy_VI(SNES snes) 10632b4ed282SShri Abhyankar { 10642b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*) snes->data; 10652b4ed282SShri Abhyankar PetscErrorCode ierr; 10662b4ed282SShri Abhyankar 10672b4ed282SShri Abhyankar PetscFunctionBegin; 10682b4ed282SShri Abhyankar if (snes->vec_sol_update) { 10692b4ed282SShri Abhyankar ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr); 10702b4ed282SShri Abhyankar snes->vec_sol_update = PETSC_NULL; 10712b4ed282SShri Abhyankar } 10722b4ed282SShri Abhyankar if (snes->nwork) { 10732b4ed282SShri Abhyankar ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 10742b4ed282SShri Abhyankar snes->nwork = 0; 10752b4ed282SShri Abhyankar snes->work = PETSC_NULL; 10762b4ed282SShri Abhyankar } 10772b4ed282SShri Abhyankar 10782b4ed282SShri Abhyankar /* clear vectors */ 10792b4ed282SShri Abhyankar ierr = VecDestroy(vi->phi); CHKERRQ(ierr); 10802b4ed282SShri Abhyankar ierr = VecDestroy(vi->Da); CHKERRQ(ierr); 10812b4ed282SShri Abhyankar ierr = VecDestroy(vi->Db); CHKERRQ(ierr); 10822b4ed282SShri Abhyankar ierr = VecDestroy(vi->z); CHKERRQ(ierr); 10832b4ed282SShri Abhyankar ierr = VecDestroy(vi->t); CHKERRQ(ierr); 10842b4ed282SShri Abhyankar if (!vi->usersetxbounds) { 10852b4ed282SShri Abhyankar ierr = VecDestroy(vi->xl); CHKERRQ(ierr); 10862b4ed282SShri Abhyankar ierr = VecDestroy(vi->xu); CHKERRQ(ierr); 10872b4ed282SShri Abhyankar } 1088c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1089c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr); 10902b4ed282SShri Abhyankar } 10912b4ed282SShri Abhyankar ierr = PetscFree(snes->data);CHKERRQ(ierr); 10922b4ed282SShri Abhyankar 10932b4ed282SShri Abhyankar /* clear composed functions */ 10942b4ed282SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr); 1095c92abb8aSShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr); 10962b4ed282SShri Abhyankar PetscFunctionReturn(0); 10972b4ed282SShri Abhyankar } 10982b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 10992b4ed282SShri Abhyankar #undef __FUNCT__ 11002b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNo_VI" 11012b4ed282SShri Abhyankar 11022b4ed282SShri Abhyankar /* 11032b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchNo routine in snes/impls/ls/ls.c 11042b4ed282SShri Abhyankar 11052b4ed282SShri Abhyankar */ 1106ace3abfcSBarry 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) 11072b4ed282SShri Abhyankar { 11082b4ed282SShri Abhyankar PetscErrorCode ierr; 11092b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1110ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 11112b4ed282SShri Abhyankar 11122b4ed282SShri Abhyankar PetscFunctionBegin; 11132b4ed282SShri Abhyankar *flag = PETSC_TRUE; 11142b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 11152b4ed282SShri Abhyankar ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 11162b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 11172b4ed282SShri Abhyankar if (vi->postcheckstep) { 11182b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 11192b4ed282SShri Abhyankar } 11202b4ed282SShri Abhyankar if (changed_y) { 11212b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 11222b4ed282SShri Abhyankar } 11232b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 11242b4ed282SShri Abhyankar if (!snes->domainerror) { 11252b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 11262b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 11272b4ed282SShri Abhyankar } 1128c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1129c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 1130c92abb8aSShri Abhyankar } 11312b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 11322b4ed282SShri Abhyankar PetscFunctionReturn(0); 11332b4ed282SShri Abhyankar } 11342b4ed282SShri Abhyankar 11352b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 11362b4ed282SShri Abhyankar #undef __FUNCT__ 11372b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchNoNorms_VI" 11382b4ed282SShri Abhyankar 11392b4ed282SShri Abhyankar /* 11402b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchNoNorms in snes/impls/ls/ls.c 11412b4ed282SShri Abhyankar */ 1142ace3abfcSBarry 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) 11432b4ed282SShri Abhyankar { 11442b4ed282SShri Abhyankar PetscErrorCode ierr; 11452b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1146ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 11472b4ed282SShri Abhyankar 11482b4ed282SShri Abhyankar PetscFunctionBegin; 11492b4ed282SShri Abhyankar *flag = PETSC_TRUE; 11502b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 11512b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 11522b4ed282SShri Abhyankar if (vi->postcheckstep) { 11532b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 11542b4ed282SShri Abhyankar } 11552b4ed282SShri Abhyankar if (changed_y) { 11562b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 11572b4ed282SShri Abhyankar } 11582b4ed282SShri Abhyankar 11592b4ed282SShri Abhyankar /* don't evaluate function the last time through */ 11602b4ed282SShri Abhyankar if (snes->iter < snes->max_its-1) { 11612b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 11622b4ed282SShri Abhyankar } 11632b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 11642b4ed282SShri Abhyankar PetscFunctionReturn(0); 11652b4ed282SShri Abhyankar } 11662b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 11672b4ed282SShri Abhyankar #undef __FUNCT__ 11682b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchCubic_VI" 11692b4ed282SShri Abhyankar /* 11702b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchCubic in snes/impls/ls/ls.c 11712b4ed282SShri Abhyankar */ 1172ace3abfcSBarry 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) 11732b4ed282SShri Abhyankar { 11742b4ed282SShri Abhyankar /* 11752b4ed282SShri Abhyankar Note that for line search purposes we work with with the related 11762b4ed282SShri Abhyankar minimization problem: 11772b4ed282SShri Abhyankar min z(x): R^n -> R, 11782b4ed282SShri Abhyankar where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 11792b4ed282SShri Abhyankar This function z(x) is same as the merit function for the complementarity problem. 11802b4ed282SShri Abhyankar */ 11812b4ed282SShri Abhyankar 11822b4ed282SShri Abhyankar PetscReal initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 11832b4ed282SShri Abhyankar PetscReal minlambda,lambda,lambdatemp; 11842b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 11852b4ed282SShri Abhyankar PetscScalar cinitslope; 11862b4ed282SShri Abhyankar #endif 11872b4ed282SShri Abhyankar PetscErrorCode ierr; 11882b4ed282SShri Abhyankar PetscInt count; 11892b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1190ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 11912b4ed282SShri Abhyankar MPI_Comm comm; 11922b4ed282SShri Abhyankar 11932b4ed282SShri Abhyankar PetscFunctionBegin; 11942b4ed282SShri Abhyankar ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 11952b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 11962b4ed282SShri Abhyankar *flag = PETSC_TRUE; 11972b4ed282SShri Abhyankar 11982b4ed282SShri Abhyankar ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 11992b4ed282SShri Abhyankar if (*ynorm == 0.0) { 1200c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1201c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 12022b4ed282SShri Abhyankar } 12032b4ed282SShri Abhyankar *gnorm = fnorm; 12042b4ed282SShri Abhyankar ierr = VecCopy(x,w);CHKERRQ(ierr); 12052b4ed282SShri Abhyankar ierr = VecCopy(f,g);CHKERRQ(ierr); 12062b4ed282SShri Abhyankar *flag = PETSC_FALSE; 12072b4ed282SShri Abhyankar goto theend1; 12082b4ed282SShri Abhyankar } 12092b4ed282SShri Abhyankar if (*ynorm > vi->maxstep) { /* Step too big, so scale back */ 1210c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1211c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Scaling step by %G old ynorm %G\n",vi->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr); 12122b4ed282SShri Abhyankar } 12132b4ed282SShri Abhyankar ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr); 12142b4ed282SShri Abhyankar *ynorm = vi->maxstep; 12152b4ed282SShri Abhyankar } 12162b4ed282SShri Abhyankar ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 12172b4ed282SShri Abhyankar minlambda = vi->minlambda/rellength; 12183b336b49SShri Abhyankar ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 12192b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 12203b336b49SShri Abhyankar ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 12212b4ed282SShri Abhyankar initslope = PetscRealPart(cinitslope); 12222b4ed282SShri Abhyankar #else 12233b336b49SShri Abhyankar ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 12242b4ed282SShri Abhyankar #endif 12252b4ed282SShri Abhyankar if (initslope > 0.0) initslope = -initslope; 12262b4ed282SShri Abhyankar if (initslope == 0.0) initslope = -1.0; 12272b4ed282SShri Abhyankar 12282b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 12292b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 12302b4ed282SShri Abhyankar ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 12312b4ed282SShri Abhyankar *flag = PETSC_FALSE; 12322b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 12332b4ed282SShri Abhyankar goto theend1; 12342b4ed282SShri Abhyankar } 12352b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 12362b4ed282SShri Abhyankar if (snes->domainerror) { 12372b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 12382b4ed282SShri Abhyankar PetscFunctionReturn(0); 12392b4ed282SShri Abhyankar } 12402b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 12412b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 12422b4ed282SShri Abhyankar ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 12432b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */ 1244c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1245c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 12462b4ed282SShri Abhyankar } 12472b4ed282SShri Abhyankar goto theend1; 12482b4ed282SShri Abhyankar } 12492b4ed282SShri Abhyankar 12502b4ed282SShri Abhyankar /* Fit points with quadratic */ 12512b4ed282SShri Abhyankar lambda = 1.0; 12522b4ed282SShri Abhyankar lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 12532b4ed282SShri Abhyankar lambdaprev = lambda; 12542b4ed282SShri Abhyankar gnormprev = *gnorm; 12552b4ed282SShri Abhyankar if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 12562b4ed282SShri Abhyankar if (lambdatemp <= .1*lambda) lambda = .1*lambda; 12572b4ed282SShri Abhyankar else lambda = lambdatemp; 12582b4ed282SShri Abhyankar 12592b4ed282SShri Abhyankar ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 12602b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 12612b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 12622b4ed282SShri Abhyankar *flag = PETSC_FALSE; 12632b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 12642b4ed282SShri Abhyankar goto theend1; 12652b4ed282SShri Abhyankar } 12662b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 12672b4ed282SShri Abhyankar if (snes->domainerror) { 12682b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 12692b4ed282SShri Abhyankar PetscFunctionReturn(0); 12702b4ed282SShri Abhyankar } 12712b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 12722b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1273c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1274c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr); 12752b4ed282SShri Abhyankar } 12762b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */ 1277c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1278c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 12792b4ed282SShri Abhyankar } 12802b4ed282SShri Abhyankar goto theend1; 12812b4ed282SShri Abhyankar } 12822b4ed282SShri Abhyankar 12832b4ed282SShri Abhyankar /* Fit points with cubic */ 12842b4ed282SShri Abhyankar count = 1; 12852b4ed282SShri Abhyankar while (PETSC_TRUE) { 12862b4ed282SShri Abhyankar if (lambda <= minlambda) { 1287c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1288c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 1289c92abb8aSShri 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); 12902b4ed282SShri Abhyankar } 12912b4ed282SShri Abhyankar *flag = PETSC_FALSE; 12922b4ed282SShri Abhyankar break; 12932b4ed282SShri Abhyankar } 12942b4ed282SShri Abhyankar t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 12952b4ed282SShri Abhyankar t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 12962b4ed282SShri Abhyankar a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 12972b4ed282SShri Abhyankar b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 12982b4ed282SShri Abhyankar d = b*b - 3*a*initslope; 12992b4ed282SShri Abhyankar if (d < 0.0) d = 0.0; 13002b4ed282SShri Abhyankar if (a == 0.0) { 13012b4ed282SShri Abhyankar lambdatemp = -initslope/(2.0*b); 13022b4ed282SShri Abhyankar } else { 13032b4ed282SShri Abhyankar lambdatemp = (-b + sqrt(d))/(3.0*a); 13042b4ed282SShri Abhyankar } 13052b4ed282SShri Abhyankar lambdaprev = lambda; 13062b4ed282SShri Abhyankar gnormprev = *gnorm; 13072b4ed282SShri Abhyankar if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 13082b4ed282SShri Abhyankar if (lambdatemp <= .1*lambda) lambda = .1*lambda; 13092b4ed282SShri Abhyankar else lambda = lambdatemp; 13102b4ed282SShri Abhyankar ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 13112b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 13122b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 13132b4ed282SShri 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); 13142b4ed282SShri Abhyankar *flag = PETSC_FALSE; 13152b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 13162b4ed282SShri Abhyankar break; 13172b4ed282SShri Abhyankar } 13182b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 13192b4ed282SShri Abhyankar if (snes->domainerror) { 13202b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 13212b4ed282SShri Abhyankar PetscFunctionReturn(0); 13222b4ed282SShri Abhyankar } 13232b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 13242b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 13252b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* is reduction enough? */ 1326c92abb8aSShri Abhyankar if (vi->lsmonitor) { 13272b4ed282SShri Abhyankar ierr = PetscPrintf(comm," Line search: Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 13282b4ed282SShri Abhyankar } 13292b4ed282SShri Abhyankar break; 13302b4ed282SShri Abhyankar } else { 1331c92abb8aSShri Abhyankar if (vi->lsmonitor) { 13322b4ed282SShri Abhyankar ierr = PetscPrintf(comm," Line search: Cubic step no good, shrinking lambda, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 13332b4ed282SShri Abhyankar } 13342b4ed282SShri Abhyankar } 13352b4ed282SShri Abhyankar count++; 13362b4ed282SShri Abhyankar } 13372b4ed282SShri Abhyankar theend1: 13382b4ed282SShri Abhyankar /* Optional user-defined check for line search step validity */ 13392b4ed282SShri Abhyankar if (vi->postcheckstep && *flag) { 13402b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 13412b4ed282SShri Abhyankar if (changed_y) { 13422b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 13432b4ed282SShri Abhyankar } 13442b4ed282SShri Abhyankar if (changed_y || changed_w) { /* recompute the function if the step has changed */ 13452b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 13462b4ed282SShri Abhyankar if (snes->domainerror) { 13472b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 13482b4ed282SShri Abhyankar PetscFunctionReturn(0); 13492b4ed282SShri Abhyankar } 13502b4ed282SShri Abhyankar ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 13512b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 13522b4ed282SShri Abhyankar ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 13532b4ed282SShri Abhyankar ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 13542b4ed282SShri Abhyankar ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 13552b4ed282SShri Abhyankar } 13562b4ed282SShri Abhyankar } 13572b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 13582b4ed282SShri Abhyankar PetscFunctionReturn(0); 13592b4ed282SShri Abhyankar } 13602b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 13612b4ed282SShri Abhyankar #undef __FUNCT__ 13622b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchQuadratic_VI" 13632b4ed282SShri Abhyankar /* 13642b4ed282SShri Abhyankar This routine is a copy of SNESLineSearchQuadratic in snes/impls/ls/ls.c 13652b4ed282SShri Abhyankar */ 1366ace3abfcSBarry 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) 13672b4ed282SShri Abhyankar { 13682b4ed282SShri Abhyankar /* 13692b4ed282SShri Abhyankar Note that for line search purposes we work with with the related 13702b4ed282SShri Abhyankar minimization problem: 13712b4ed282SShri Abhyankar min z(x): R^n -> R, 13722b4ed282SShri Abhyankar where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 13732b4ed282SShri Abhyankar z(x) is the same as the merit function for the complementarity problem 13742b4ed282SShri Abhyankar */ 13752b4ed282SShri Abhyankar PetscReal initslope,minlambda,lambda,lambdatemp,rellength; 13762b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 13772b4ed282SShri Abhyankar PetscScalar cinitslope; 13782b4ed282SShri Abhyankar #endif 13792b4ed282SShri Abhyankar PetscErrorCode ierr; 13802b4ed282SShri Abhyankar PetscInt count; 13812b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 1382ace3abfcSBarry Smith PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 13832b4ed282SShri Abhyankar 13842b4ed282SShri Abhyankar PetscFunctionBegin; 13852b4ed282SShri Abhyankar ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 13862b4ed282SShri Abhyankar *flag = PETSC_TRUE; 13872b4ed282SShri Abhyankar 13882b4ed282SShri Abhyankar ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 13892b4ed282SShri Abhyankar if (*ynorm == 0.0) { 1390c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1391c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr); 13922b4ed282SShri Abhyankar } 13932b4ed282SShri Abhyankar *gnorm = fnorm; 13942b4ed282SShri Abhyankar ierr = VecCopy(x,w);CHKERRQ(ierr); 13952b4ed282SShri Abhyankar ierr = VecCopy(f,g);CHKERRQ(ierr); 13962b4ed282SShri Abhyankar *flag = PETSC_FALSE; 13972b4ed282SShri Abhyankar goto theend2; 13982b4ed282SShri Abhyankar } 13992b4ed282SShri Abhyankar if (*ynorm > vi->maxstep) { /* Step too big, so scale back */ 14002b4ed282SShri Abhyankar ierr = VecScale(y,vi->maxstep/(*ynorm));CHKERRQ(ierr); 14012b4ed282SShri Abhyankar *ynorm = vi->maxstep; 14022b4ed282SShri Abhyankar } 14032b4ed282SShri Abhyankar ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 14042b4ed282SShri Abhyankar minlambda = vi->minlambda/rellength; 14053b336b49SShri Abhyankar ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 14062b4ed282SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 14072b4ed282SShri Abhyankar ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 14082b4ed282SShri Abhyankar initslope = PetscRealPart(cinitslope); 14092b4ed282SShri Abhyankar #else 14103b336b49SShri Abhyankar ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 14112b4ed282SShri Abhyankar #endif 14122b4ed282SShri Abhyankar if (initslope > 0.0) initslope = -initslope; 14132b4ed282SShri Abhyankar if (initslope == 0.0) initslope = -1.0; 14142b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr); 14152b4ed282SShri Abhyankar 14162b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 14172b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 14182b4ed282SShri Abhyankar ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 14192b4ed282SShri Abhyankar *flag = PETSC_FALSE; 14202b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 14212b4ed282SShri Abhyankar goto theend2; 14222b4ed282SShri Abhyankar } 14232b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 14242b4ed282SShri Abhyankar if (snes->domainerror) { 14252b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 14262b4ed282SShri Abhyankar PetscFunctionReturn(0); 14272b4ed282SShri Abhyankar } 14282b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 14292b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 14302b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + vi->alpha*initslope) { /* Sufficient reduction */ 1431c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1432c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 14332b4ed282SShri Abhyankar } 14342b4ed282SShri Abhyankar goto theend2; 14352b4ed282SShri Abhyankar } 14362b4ed282SShri Abhyankar 14372b4ed282SShri Abhyankar /* Fit points with quadratic */ 14382b4ed282SShri Abhyankar lambda = 1.0; 14392b4ed282SShri Abhyankar count = 1; 14402b4ed282SShri Abhyankar while (PETSC_TRUE) { 14412b4ed282SShri Abhyankar if (lambda <= minlambda) { /* bad luck; use full step */ 1442c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1443c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr); 1444c92abb8aSShri 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); 14452b4ed282SShri Abhyankar } 14462b4ed282SShri Abhyankar ierr = VecCopy(x,w);CHKERRQ(ierr); 14472b4ed282SShri Abhyankar *flag = PETSC_FALSE; 14482b4ed282SShri Abhyankar break; 14492b4ed282SShri Abhyankar } 14502b4ed282SShri Abhyankar lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 14512b4ed282SShri Abhyankar if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 14522b4ed282SShri Abhyankar if (lambdatemp <= .1*lambda) lambda = .1*lambda; 14532b4ed282SShri Abhyankar else lambda = lambdatemp; 14542b4ed282SShri Abhyankar 14552b4ed282SShri Abhyankar ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 14562b4ed282SShri Abhyankar if (snes->nfuncs >= snes->max_funcs) { 14572b4ed282SShri Abhyankar ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 14582b4ed282SShri 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); 14592b4ed282SShri Abhyankar *flag = PETSC_FALSE; 14602b4ed282SShri Abhyankar snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 14612b4ed282SShri Abhyankar break; 14622b4ed282SShri Abhyankar } 14632b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 14642b4ed282SShri Abhyankar if (snes->domainerror) { 14652b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 14662b4ed282SShri Abhyankar PetscFunctionReturn(0); 14672b4ed282SShri Abhyankar } 14682b4ed282SShri Abhyankar ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 14692b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 14702b4ed282SShri Abhyankar if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*vi->alpha*initslope) { /* sufficient reduction */ 1471c92abb8aSShri Abhyankar if (vi->lsmonitor) { 1472c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorPrintf(vi->lsmonitor," Line Search: Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr); 14732b4ed282SShri Abhyankar } 14742b4ed282SShri Abhyankar break; 14752b4ed282SShri Abhyankar } 14762b4ed282SShri Abhyankar count++; 14772b4ed282SShri Abhyankar } 14782b4ed282SShri Abhyankar theend2: 14792b4ed282SShri Abhyankar /* Optional user-defined check for line search step validity */ 14802b4ed282SShri Abhyankar if (vi->postcheckstep) { 14812b4ed282SShri Abhyankar ierr = (*vi->postcheckstep)(snes,x,y,w,vi->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 14822b4ed282SShri Abhyankar if (changed_y) { 14832b4ed282SShri Abhyankar ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 14842b4ed282SShri Abhyankar } 14852b4ed282SShri Abhyankar if (changed_y || changed_w) { /* recompute the function if the step has changed */ 14862b4ed282SShri Abhyankar ierr = SNESComputeFunction(snes,w,g); 14872b4ed282SShri Abhyankar if (snes->domainerror) { 14882b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 14892b4ed282SShri Abhyankar PetscFunctionReturn(0); 14902b4ed282SShri Abhyankar } 14912b4ed282SShri Abhyankar ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 14922b4ed282SShri Abhyankar ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 14932b4ed282SShri Abhyankar ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 14942b4ed282SShri Abhyankar ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 14952b4ed282SShri Abhyankar if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 14962b4ed282SShri Abhyankar } 14972b4ed282SShri Abhyankar } 14982b4ed282SShri Abhyankar ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 14992b4ed282SShri Abhyankar PetscFunctionReturn(0); 15002b4ed282SShri Abhyankar } 15012b4ed282SShri Abhyankar 1502ace3abfcSBarry 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*/ 15032b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 15042b4ed282SShri Abhyankar EXTERN_C_BEGIN 15052b4ed282SShri Abhyankar #undef __FUNCT__ 15062b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchSet_VI" 15072b4ed282SShri Abhyankar PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_VI(SNES snes,FCN2 func,void *lsctx) 15082b4ed282SShri Abhyankar { 15092b4ed282SShri Abhyankar PetscFunctionBegin; 15102b4ed282SShri Abhyankar ((SNES_VI *)(snes->data))->LineSearch = func; 15112b4ed282SShri Abhyankar ((SNES_VI *)(snes->data))->lsP = lsctx; 15122b4ed282SShri Abhyankar PetscFunctionReturn(0); 15132b4ed282SShri Abhyankar } 15142b4ed282SShri Abhyankar EXTERN_C_END 15152b4ed282SShri Abhyankar 15162b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 15172b4ed282SShri Abhyankar EXTERN_C_BEGIN 15182b4ed282SShri Abhyankar #undef __FUNCT__ 15192b4ed282SShri Abhyankar #define __FUNCT__ "SNESLineSearchSetMonitor_VI" 1520ace3abfcSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetMonitor_VI(SNES snes,PetscBool flg) 15212b4ed282SShri Abhyankar { 15222b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 15232b4ed282SShri Abhyankar PetscErrorCode ierr; 15242b4ed282SShri Abhyankar 15252b4ed282SShri Abhyankar PetscFunctionBegin; 1526c92abb8aSShri Abhyankar if (flg && !vi->lsmonitor) { 1527c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",((PetscObject)snes)->tablevel,&vi->lsmonitor);CHKERRQ(ierr); 1528c92abb8aSShri Abhyankar } else if (!flg && vi->lsmonitor) { 1529c92abb8aSShri Abhyankar ierr = PetscViewerASCIIMonitorDestroy(vi->lsmonitor);CHKERRQ(ierr); 15302b4ed282SShri Abhyankar } 15312b4ed282SShri Abhyankar PetscFunctionReturn(0); 15322b4ed282SShri Abhyankar } 15332b4ed282SShri Abhyankar EXTERN_C_END 15342b4ed282SShri Abhyankar 15352b4ed282SShri Abhyankar /* 15362b4ed282SShri Abhyankar SNESView_VI - Prints info from the SNESVI data structure. 15372b4ed282SShri Abhyankar 15382b4ed282SShri Abhyankar Input Parameters: 15392b4ed282SShri Abhyankar . SNES - the SNES context 15402b4ed282SShri Abhyankar . viewer - visualization context 15412b4ed282SShri Abhyankar 15422b4ed282SShri Abhyankar Application Interface Routine: SNESView() 15432b4ed282SShri Abhyankar */ 15442b4ed282SShri Abhyankar #undef __FUNCT__ 15452b4ed282SShri Abhyankar #define __FUNCT__ "SNESView_VI" 15462b4ed282SShri Abhyankar static PetscErrorCode SNESView_VI(SNES snes,PetscViewer viewer) 15472b4ed282SShri Abhyankar { 15482b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI *)snes->data; 154978c4b9d3SShri Abhyankar const char *cstr,*tstr; 15502b4ed282SShri Abhyankar PetscErrorCode ierr; 1551ace3abfcSBarry Smith PetscBool iascii; 15522b4ed282SShri Abhyankar 15532b4ed282SShri Abhyankar PetscFunctionBegin; 15542b4ed282SShri Abhyankar ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 15552b4ed282SShri Abhyankar if (iascii) { 15562b4ed282SShri Abhyankar if (vi->LineSearch == SNESLineSearchNo_VI) cstr = "SNESLineSearchNo"; 15572b4ed282SShri Abhyankar else if (vi->LineSearch == SNESLineSearchQuadratic_VI) cstr = "SNESLineSearchQuadratic"; 15582b4ed282SShri Abhyankar else if (vi->LineSearch == SNESLineSearchCubic_VI) cstr = "SNESLineSearchCubic"; 15592b4ed282SShri Abhyankar else cstr = "unknown"; 156078c4b9d3SShri Abhyankar if (snes->ops->solve == SNESSolveVI_SS) tstr = "Semismooth"; 156178c4b9d3SShri Abhyankar else if (snes->ops->solve == SNESSolveVI_AS) tstr = "Active Set"; 156278c4b9d3SShri Abhyankar else tstr = "unknown"; 156378c4b9d3SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," VI algorithm: %s\n",tstr);CHKERRQ(ierr); 15642b4ed282SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 15652b4ed282SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," alpha=%G, maxstep=%G, minlambda=%G\n",vi->alpha,vi->maxstep,vi->minlambda);CHKERRQ(ierr); 15662b4ed282SShri Abhyankar } else { 15672b4ed282SShri Abhyankar SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ VI",((PetscObject)viewer)->type_name); 15682b4ed282SShri Abhyankar } 15692b4ed282SShri Abhyankar PetscFunctionReturn(0); 15702b4ed282SShri Abhyankar } 15712b4ed282SShri Abhyankar 15722b4ed282SShri Abhyankar /* 15732b4ed282SShri Abhyankar SNESVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu. 15742b4ed282SShri Abhyankar 15752b4ed282SShri Abhyankar Input Parameters: 15762b4ed282SShri Abhyankar . snes - the SNES context. 15772b4ed282SShri Abhyankar . xl - lower bound. 15782b4ed282SShri Abhyankar . xu - upper bound. 15792b4ed282SShri Abhyankar 15802b4ed282SShri Abhyankar Notes: 15812b4ed282SShri Abhyankar If this routine is not called then the lower and upper bounds are set to 158284c105d7SBarry Smith PETSC_VI_INF and PETSC_VI_NINF respectively during SNESSetUp(). 158384c105d7SBarry Smith 15842b4ed282SShri Abhyankar */ 15852b4ed282SShri Abhyankar 15862b4ed282SShri Abhyankar #undef __FUNCT__ 15872b4ed282SShri Abhyankar #define __FUNCT__ "SNESVISetVariableBounds" 158869c03620SShri Abhyankar PetscErrorCode SNESVISetVariableBounds(SNES snes, Vec xl, Vec xu) 15892b4ed282SShri Abhyankar { 15902b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI*)snes->data; 15912b4ed282SShri Abhyankar 15922b4ed282SShri Abhyankar PetscFunctionBegin; 15932b4ed282SShri Abhyankar PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 15942b4ed282SShri Abhyankar PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 15952b4ed282SShri Abhyankar PetscValidHeaderSpecific(xu,VEC_CLASSID,3); 15962b4ed282SShri Abhyankar 15972b4ed282SShri Abhyankar /* Check if SNESSetFunction is called */ 15982b4ed282SShri Abhyankar if(!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 15992b4ed282SShri Abhyankar 16002b4ed282SShri Abhyankar /* Check if the vector sizes are compatible for lower and upper bounds */ 16012b4ed282SShri 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); 16022b4ed282SShri 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); 16032b4ed282SShri Abhyankar vi->usersetxbounds = PETSC_TRUE; 16042b4ed282SShri Abhyankar vi->xl = xl; 16052b4ed282SShri Abhyankar vi->xu = xu; 16062b4ed282SShri Abhyankar 16072b4ed282SShri Abhyankar PetscFunctionReturn(0); 16082b4ed282SShri Abhyankar } 16092b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 16102b4ed282SShri Abhyankar /* 16112b4ed282SShri Abhyankar SNESSetFromOptions_VI - Sets various parameters for the SNESVI method. 16122b4ed282SShri Abhyankar 16132b4ed282SShri Abhyankar Input Parameter: 16142b4ed282SShri Abhyankar . snes - the SNES context 16152b4ed282SShri Abhyankar 16162b4ed282SShri Abhyankar Application Interface Routine: SNESSetFromOptions() 16172b4ed282SShri Abhyankar */ 16182b4ed282SShri Abhyankar #undef __FUNCT__ 16192b4ed282SShri Abhyankar #define __FUNCT__ "SNESSetFromOptions_VI" 16202b4ed282SShri Abhyankar static PetscErrorCode SNESSetFromOptions_VI(SNES snes) 16212b4ed282SShri Abhyankar { 16222b4ed282SShri Abhyankar SNES_VI *vi = (SNES_VI *)snes->data; 16232b4ed282SShri Abhyankar const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 162469c03620SShri Abhyankar const char *vies[] = {"ss","as"}; 16252b4ed282SShri Abhyankar PetscErrorCode ierr; 16262b4ed282SShri Abhyankar PetscInt indx; 162769c03620SShri Abhyankar PetscBool flg,set,flg2; 16282b4ed282SShri Abhyankar 16292b4ed282SShri Abhyankar PetscFunctionBegin; 16302b4ed282SShri Abhyankar ierr = PetscOptionsHead("SNES semismooth method options");CHKERRQ(ierr); 1631*9308c137SBarry Smith ierr = PetscOptionsBool("-snes_vi_monitor","Monitor all non-active variables","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 1632*9308c137SBarry Smith if (flg) { 1633*9308c137SBarry Smith ierr = SNESMonitorSet(snes,SNESMonitorVI,0,0);CHKERRQ(ierr); 1634*9308c137SBarry Smith } 16352b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_alpha","Function norm must decrease by","None",vi->alpha,&vi->alpha,0);CHKERRQ(ierr); 16362b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_maxstep","Step must be less than","None",vi->maxstep,&vi->maxstep,0);CHKERRQ(ierr); 16372b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_minlambda","Minimum lambda allowed","None",vi->minlambda,&vi->minlambda,0);CHKERRQ(ierr); 16382b4ed282SShri Abhyankar ierr = PetscOptionsReal("-snes_vi_const_tol","constraint tolerance","None",vi->const_tol,&vi->const_tol,0);CHKERRQ(ierr); 1639acfcf0e5SJed Brown ierr = PetscOptionsBool("-snes_vi_lsmonitor","Print progress of line searches","SNESLineSearchSetMonitor",vi->lsmonitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 16402b4ed282SShri Abhyankar if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);} 164169c03620SShri Abhyankar ierr = PetscOptionsEList("-snes_vi_type","Semismooth algorithm used","",vies,2,"ss",&indx,&flg2);CHKERRQ(ierr); 164269c03620SShri Abhyankar if (flg2) { 164369c03620SShri Abhyankar switch (indx) { 164469c03620SShri Abhyankar case 0: 164569c03620SShri Abhyankar snes->ops->solve = SNESSolveVI_SS; 164669c03620SShri Abhyankar break; 164769c03620SShri Abhyankar case 1: 164869c03620SShri Abhyankar snes->ops->solve = SNESSolveVI_AS; 164969c03620SShri Abhyankar break; 165069c03620SShri Abhyankar } 165169c03620SShri Abhyankar } 1652c92abb8aSShri Abhyankar ierr = PetscOptionsEList("-snes_vi_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 16532b4ed282SShri Abhyankar if (flg) { 16542b4ed282SShri Abhyankar switch (indx) { 16552b4ed282SShri Abhyankar case 0: 16562b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchNo_VI,PETSC_NULL);CHKERRQ(ierr); 16572b4ed282SShri Abhyankar break; 16582b4ed282SShri Abhyankar case 1: 16592b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_VI,PETSC_NULL);CHKERRQ(ierr); 16602b4ed282SShri Abhyankar break; 16612b4ed282SShri Abhyankar case 2: 16622b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_VI,PETSC_NULL);CHKERRQ(ierr); 16632b4ed282SShri Abhyankar break; 16642b4ed282SShri Abhyankar case 3: 16652b4ed282SShri Abhyankar ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_VI,PETSC_NULL);CHKERRQ(ierr); 16662b4ed282SShri Abhyankar break; 16672b4ed282SShri Abhyankar } 16682b4ed282SShri Abhyankar } 16692b4ed282SShri Abhyankar ierr = PetscOptionsTail();CHKERRQ(ierr); 16702b4ed282SShri Abhyankar PetscFunctionReturn(0); 16712b4ed282SShri Abhyankar } 16722b4ed282SShri Abhyankar /* -------------------------------------------------------------------------- */ 16732b4ed282SShri Abhyankar /*MC 16742b4ed282SShri Abhyankar SNESVI - Semismooth newton method based nonlinear solver that uses a line search 16752b4ed282SShri Abhyankar 16762b4ed282SShri Abhyankar Options Database: 16772b4ed282SShri Abhyankar + -snes_vi [cubic,quadratic,basic,basicnonorms] - Selects line search 16782b4ed282SShri Abhyankar . -snes_vi_alpha <alpha> - Sets alpha 16792b4ed282SShri 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) 16802b4ed282SShri Abhyankar . -snes_vi_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 16812b4ed282SShri Abhyankar - -snes_vi_monitor - print information about progress of line searches 16822b4ed282SShri Abhyankar 16832b4ed282SShri Abhyankar 16842b4ed282SShri Abhyankar Level: beginner 16852b4ed282SShri Abhyankar 16862b4ed282SShri Abhyankar .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 16872b4ed282SShri Abhyankar SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 16882b4ed282SShri Abhyankar SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams() 16892b4ed282SShri Abhyankar 16902b4ed282SShri Abhyankar M*/ 16912b4ed282SShri Abhyankar EXTERN_C_BEGIN 16922b4ed282SShri Abhyankar #undef __FUNCT__ 16932b4ed282SShri Abhyankar #define __FUNCT__ "SNESCreate_VI" 16942b4ed282SShri Abhyankar PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_VI(SNES snes) 16952b4ed282SShri Abhyankar { 16962b4ed282SShri Abhyankar PetscErrorCode ierr; 16972b4ed282SShri Abhyankar SNES_VI *vi; 16982b4ed282SShri Abhyankar 16992b4ed282SShri Abhyankar PetscFunctionBegin; 17002b4ed282SShri Abhyankar snes->ops->setup = SNESSetUp_VI; 170169c03620SShri Abhyankar snes->ops->solve = SNESSolveVI_SS; 17022b4ed282SShri Abhyankar snes->ops->destroy = SNESDestroy_VI; 17032b4ed282SShri Abhyankar snes->ops->setfromoptions = SNESSetFromOptions_VI; 17042b4ed282SShri Abhyankar snes->ops->view = SNESView_VI; 17052b4ed282SShri Abhyankar snes->ops->converged = SNESDefaultConverged_VI; 17062b4ed282SShri Abhyankar 17072b4ed282SShri Abhyankar ierr = PetscNewLog(snes,SNES_VI,&vi);CHKERRQ(ierr); 17082b4ed282SShri Abhyankar snes->data = (void*)vi; 17092b4ed282SShri Abhyankar vi->alpha = 1.e-4; 17102b4ed282SShri Abhyankar vi->maxstep = 1.e8; 17112b4ed282SShri Abhyankar vi->minlambda = 1.e-12; 17122b4ed282SShri Abhyankar vi->LineSearch = SNESLineSearchCubic_VI; 17132b4ed282SShri Abhyankar vi->lsP = PETSC_NULL; 17142b4ed282SShri Abhyankar vi->postcheckstep = PETSC_NULL; 17152b4ed282SShri Abhyankar vi->postcheck = PETSC_NULL; 17162b4ed282SShri Abhyankar vi->precheckstep = PETSC_NULL; 17172b4ed282SShri Abhyankar vi->precheck = PETSC_NULL; 17182b4ed282SShri Abhyankar vi->const_tol = 2.2204460492503131e-16; 17192b4ed282SShri Abhyankar vi->computessfunction = ComputeFischerFunction; 17202b4ed282SShri Abhyankar 17212b4ed282SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_VI",SNESLineSearchSetMonitor_VI);CHKERRQ(ierr); 17222b4ed282SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_VI",SNESLineSearchSet_VI);CHKERRQ(ierr); 17232b4ed282SShri Abhyankar 17242b4ed282SShri Abhyankar PetscFunctionReturn(0); 17252b4ed282SShri Abhyankar } 17262b4ed282SShri Abhyankar EXTERN_C_END 1727