1c2fc9fa9SBarry Smith 2c2fc9fa9SBarry Smith #include <../src/snes/impls/vi/ss/vissimpl.h> /*I "petscsnes.h" I*/ 3c2fc9fa9SBarry Smith 4*f6dfbefdSBarry Smith /*@ 5c2fc9fa9SBarry Smith SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem. 6c2fc9fa9SBarry Smith 7c2fc9fa9SBarry Smith Input Parameter: 8*f6dfbefdSBarry Smith . phi - the `Vec` holding the evaluation of the semismooth function 9c2fc9fa9SBarry Smith 10*f6dfbefdSBarry Smith Output Parameters: 11*f6dfbefdSBarry Smith + merit - the merit function 1/2 ||phi||^2 12*f6dfbefdSBarry Smith - phinorm - the two-norm of the vector, ||phi|| 13c2fc9fa9SBarry Smith 14*f6dfbefdSBarry Smith .seealso: `SNESVINEWTONSSLS`, `SNESVIComputeFunction()` 15*f6dfbefdSBarry Smith @*/ 16*f6dfbefdSBarry Smith PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal *merit, PetscReal *phinorm) { 17c2fc9fa9SBarry Smith PetscFunctionBegin; 189566063dSJacob Faibussowitsch PetscCall(VecNormBegin(phi, NORM_2, phinorm)); 199566063dSJacob Faibussowitsch PetscCall(VecNormEnd(phi, NORM_2, phinorm)); 20c2fc9fa9SBarry Smith *merit = 0.5 * (*phinorm) * (*phinorm); 21c2fc9fa9SBarry Smith PetscFunctionReturn(0); 22c2fc9fa9SBarry Smith } 23c2fc9fa9SBarry Smith 249371c9d4SSatish Balay static inline PetscScalar Phi(PetscScalar a, PetscScalar b) { 25c2fc9fa9SBarry Smith return a + b - PetscSqrtScalar(a * a + b * b); 26c2fc9fa9SBarry Smith } 27c2fc9fa9SBarry Smith 289371c9d4SSatish Balay static inline PetscScalar DPhi(PetscScalar a, PetscScalar b) { 29c2fc9fa9SBarry Smith if ((PetscAbsScalar(a) >= 1.e-6) || (PetscAbsScalar(b) >= 1.e-6)) return 1.0 - a / PetscSqrtScalar(a * a + b * b); 30c2fc9fa9SBarry Smith else return .5; 31c2fc9fa9SBarry Smith } 32c2fc9fa9SBarry Smith 33*f6dfbefdSBarry Smith /*@ 34*f6dfbefdSBarry Smith SNESVIComputeFunction - Provides the function that reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear 35*f6dfbefdSBarry Smith equations in semismooth form. 36c2fc9fa9SBarry Smith 37c2fc9fa9SBarry Smith Input Parameters: 38*f6dfbefdSBarry Smith + snes - the SNES context 39d5f1b7e6SEd Bueler . X - current iterate 40*f6dfbefdSBarry Smith - functx - user defined function context 41c2fc9fa9SBarry Smith 42*f6dfbefdSBarry Smith Output Parameter: 43*f6dfbefdSBarry Smith . phi - the evaluation of Semismooth function at X 44c2fc9fa9SBarry Smith 45*f6dfbefdSBarry Smith .seealso: `SNESVINEWTONSSLS`, `SNESVIComputeMeritFunction()` 46*f6dfbefdSBarry Smith @*/ 47*f6dfbefdSBarry Smith PetscErrorCode SNESVIComputeFunction(SNES snes, Vec X, Vec phi, void *functx) { 48f450aa47SBarry Smith SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data; 49c2fc9fa9SBarry Smith Vec Xl = snes->xl, Xu = snes->xu, F = snes->vec_func; 505edff71fSBarry Smith PetscScalar *phi_arr, *f_arr, *l, *u; 515edff71fSBarry Smith const PetscScalar *x_arr; 52c2fc9fa9SBarry Smith PetscInt i, nlocal; 53c2fc9fa9SBarry Smith 54c2fc9fa9SBarry Smith PetscFunctionBegin; 559566063dSJacob Faibussowitsch PetscCall((*vi->computeuserfunction)(snes, X, F, functx)); 569566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &nlocal)); 579566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x_arr)); 589566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f_arr)); 599566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xl, &l)); 609566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xu, &u)); 619566063dSJacob Faibussowitsch PetscCall(VecGetArray(phi, &phi_arr)); 62c2fc9fa9SBarry Smith 63c2fc9fa9SBarry Smith for (i = 0; i < nlocal; i++) { 64e270355aSBarry Smith if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { /* no constraints on variable */ 65c2fc9fa9SBarry Smith phi_arr[i] = f_arr[i]; 66e270355aSBarry Smith } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { /* upper bound on variable only */ 67c2fc9fa9SBarry Smith phi_arr[i] = -Phi(u[i] - x_arr[i], -f_arr[i]); 68e270355aSBarry Smith } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { /* lower bound on variable only */ 69c2fc9fa9SBarry Smith phi_arr[i] = Phi(x_arr[i] - l[i], f_arr[i]); 70c2fc9fa9SBarry Smith } else if (l[i] == u[i]) { 71c2fc9fa9SBarry Smith phi_arr[i] = l[i] - x_arr[i]; 72c2fc9fa9SBarry Smith } else { /* both bounds on variable */ 73c2fc9fa9SBarry Smith phi_arr[i] = Phi(x_arr[i] - l[i], -Phi(u[i] - x_arr[i], -f_arr[i])); 74c2fc9fa9SBarry Smith } 75c2fc9fa9SBarry Smith } 76c2fc9fa9SBarry Smith 779566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x_arr)); 789566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f_arr)); 799566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xl, &l)); 809566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xu, &u)); 819566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(phi, &phi_arr)); 82c2fc9fa9SBarry Smith PetscFunctionReturn(0); 83c2fc9fa9SBarry Smith } 84c2fc9fa9SBarry Smith 85c2fc9fa9SBarry Smith /* 86c2fc9fa9SBarry Smith SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the 87c2fc9fa9SBarry Smith the semismooth jacobian. 88c2fc9fa9SBarry Smith */ 899371c9d4SSatish Balay PetscErrorCode SNESVIComputeBsubdifferentialVectors(SNES snes, Vec X, Vec F, Mat jac, Vec Da, Vec Db) { 90c2fc9fa9SBarry Smith PetscScalar *l, *u, *x, *f, *da, *db, da1, da2, db1, db2; 91c2fc9fa9SBarry Smith PetscInt i, nlocal; 92c2fc9fa9SBarry Smith 93c2fc9fa9SBarry Smith PetscFunctionBegin; 949566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 959566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 969566063dSJacob Faibussowitsch PetscCall(VecGetArray(snes->xl, &l)); 979566063dSJacob Faibussowitsch PetscCall(VecGetArray(snes->xu, &u)); 989566063dSJacob Faibussowitsch PetscCall(VecGetArray(Da, &da)); 999566063dSJacob Faibussowitsch PetscCall(VecGetArray(Db, &db)); 1009566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &nlocal)); 101c2fc9fa9SBarry Smith 102c2fc9fa9SBarry Smith for (i = 0; i < nlocal; i++) { 103e270355aSBarry Smith if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { /* no constraints on variable */ 104c2fc9fa9SBarry Smith da[i] = 0; 105c2fc9fa9SBarry Smith db[i] = 1; 106e270355aSBarry Smith } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { /* upper bound on variable only */ 107c2fc9fa9SBarry Smith da[i] = DPhi(u[i] - x[i], -f[i]); 108c2fc9fa9SBarry Smith db[i] = DPhi(-f[i], u[i] - x[i]); 109e270355aSBarry Smith } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { /* lower bound on variable only */ 110c2fc9fa9SBarry Smith da[i] = DPhi(x[i] - l[i], f[i]); 111c2fc9fa9SBarry Smith db[i] = DPhi(f[i], x[i] - l[i]); 112c2fc9fa9SBarry Smith } else if (l[i] == u[i]) { /* fixed variable */ 113c2fc9fa9SBarry Smith da[i] = 1; 114c2fc9fa9SBarry Smith db[i] = 0; 115c2fc9fa9SBarry Smith } else { /* upper and lower bounds on variable */ 116c2fc9fa9SBarry Smith da1 = DPhi(x[i] - l[i], -Phi(u[i] - x[i], -f[i])); 117c2fc9fa9SBarry Smith db1 = DPhi(-Phi(u[i] - x[i], -f[i]), x[i] - l[i]); 118c2fc9fa9SBarry Smith da2 = DPhi(u[i] - x[i], -f[i]); 119c2fc9fa9SBarry Smith db2 = DPhi(-f[i], u[i] - x[i]); 120c2fc9fa9SBarry Smith da[i] = da1 + db1 * da2; 121c2fc9fa9SBarry Smith db[i] = db1 * db2; 122c2fc9fa9SBarry Smith } 123c2fc9fa9SBarry Smith } 124c2fc9fa9SBarry Smith 1259566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 1269566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 1279566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(snes->xl, &l)); 1289566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(snes->xu, &u)); 1299566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Da, &da)); 1309566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Db, &db)); 131c2fc9fa9SBarry Smith PetscFunctionReturn(0); 132c2fc9fa9SBarry Smith } 133c2fc9fa9SBarry Smith 134c2fc9fa9SBarry Smith /* 135c2fc9fa9SBarry Smith 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. 136c2fc9fa9SBarry Smith 137c2fc9fa9SBarry Smith Input Parameters: 138c2fc9fa9SBarry Smith . Da - Diagonal shift vector for the semismooth jacobian. 139c2fc9fa9SBarry Smith . Db - Row scaling vector for the semismooth jacobian. 140c2fc9fa9SBarry Smith 141c2fc9fa9SBarry Smith Output Parameters: 142c2fc9fa9SBarry Smith . jac - semismooth jacobian 143c2fc9fa9SBarry Smith . jac_pre - optional preconditioning matrix 144c2fc9fa9SBarry Smith 145*f6dfbefdSBarry Smith Note: 146c2fc9fa9SBarry Smith The semismooth jacobian matrix is given by 147c2fc9fa9SBarry Smith jac = Da + Db*jacfun 148c2fc9fa9SBarry Smith where Db is the row scaling matrix stored as a vector, 149c2fc9fa9SBarry Smith Da is the diagonal perturbation matrix stored as a vector 150c2fc9fa9SBarry Smith and jacfun is the jacobian of the original nonlinear function. 151c2fc9fa9SBarry Smith */ 1529371c9d4SSatish Balay PetscErrorCode SNESVIComputeJacobian(Mat jac, Mat jac_pre, Vec Da, Vec Db) { 153c2fc9fa9SBarry Smith /* Do row scaling and add diagonal perturbation */ 154362febeeSStefano Zampini PetscFunctionBegin; 1559566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(jac, Db, NULL)); 1569566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(jac, Da, ADD_VALUES)); 157c2fc9fa9SBarry Smith if (jac != jac_pre) { /* If jac and jac_pre are different */ 1589566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(jac_pre, Db, NULL)); 1599566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(jac_pre, Da, ADD_VALUES)); 160c2fc9fa9SBarry Smith } 161c2fc9fa9SBarry Smith PetscFunctionReturn(0); 162c2fc9fa9SBarry Smith } 163c2fc9fa9SBarry Smith 164c2fc9fa9SBarry Smith /* 165c2fc9fa9SBarry Smith SNESVIComputeMeritFunctionGradient - Computes the gradient of the merit function psi. 166c2fc9fa9SBarry Smith 167c2fc9fa9SBarry Smith Input Parameters: 168c2fc9fa9SBarry Smith phi - semismooth function. 169c2fc9fa9SBarry Smith H - semismooth jacobian 170c2fc9fa9SBarry Smith 171*f6dfbefdSBarry Smith Output Parameter: 172c2fc9fa9SBarry Smith dpsi - merit function gradient 173c2fc9fa9SBarry Smith 174*f6dfbefdSBarry Smith Note: 175c2fc9fa9SBarry Smith The merit function gradient is computed as follows 176c2fc9fa9SBarry Smith dpsi = H^T*phi 177c2fc9fa9SBarry Smith */ 1789371c9d4SSatish Balay PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi) { 179c2fc9fa9SBarry Smith PetscFunctionBegin; 1809566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(H, phi, dpsi)); 181c2fc9fa9SBarry Smith PetscFunctionReturn(0); 182c2fc9fa9SBarry Smith } 183c2fc9fa9SBarry Smith 184c2fc9fa9SBarry Smith /* 185f450aa47SBarry Smith SNESSolve_VINEWTONSSLS - Solves the complementarity problem with a semismooth Newton 186c2fc9fa9SBarry Smith method using a line search. 187c2fc9fa9SBarry Smith 188*f6dfbefdSBarry Smith Input Parameter: 189c2fc9fa9SBarry Smith . snes - the SNES context 190c2fc9fa9SBarry Smith 191c2fc9fa9SBarry Smith Application Interface Routine: SNESSolve() 192c2fc9fa9SBarry Smith 193*f6dfbefdSBarry Smith Note: 194c2fc9fa9SBarry Smith This implements essentially a semismooth Newton method with a 195d5f1b7e6SEd Bueler line search. The default line search does not do any line search 196d5f1b7e6SEd Bueler but rather takes a full Newton step. 197016d19f9SBarry Smith 19804d7464bSBarry Smith Developer Note: the code in this file should be slightly modified so that this routine need not exist and the SNESSolve_NEWTONLS() routine is called directly with the appropriate wrapped function and Jacobian evaluations 199016d19f9SBarry Smith 200c2fc9fa9SBarry Smith */ 2019371c9d4SSatish Balay PetscErrorCode SNESSolve_VINEWTONSSLS(SNES snes) { 202f450aa47SBarry Smith SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data; 203c2fc9fa9SBarry Smith PetscInt maxits, i, lits; 204422a814eSBarry Smith SNESLineSearchReason lssucceed; 205c2fc9fa9SBarry Smith PetscReal gnorm, xnorm = 0, ynorm; 2069bd66eb0SPeter Brune Vec Y, X, F; 207c2fc9fa9SBarry Smith KSPConvergedReason kspreason; 2086cab3a1bSJed Brown DM dm; 209942e3340SBarry Smith DMSNES sdm; 210c2fc9fa9SBarry Smith 211c2fc9fa9SBarry Smith PetscFunctionBegin; 2129566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 2139566063dSJacob Faibussowitsch PetscCall(DMGetDMSNES(dm, &sdm)); 2141aa26658SKarl Rupp 21522c6f798SBarry Smith vi->computeuserfunction = sdm->ops->computefunction; 21622c6f798SBarry Smith sdm->ops->computefunction = SNESVIComputeFunction; 217c2fc9fa9SBarry Smith 218c2fc9fa9SBarry Smith snes->numFailures = 0; 219c2fc9fa9SBarry Smith snes->numLinearSolveFailures = 0; 220c2fc9fa9SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 221c2fc9fa9SBarry Smith 222c2fc9fa9SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 223c2fc9fa9SBarry Smith X = snes->vec_sol; /* solution vector */ 224c2fc9fa9SBarry Smith F = snes->vec_func; /* residual vector */ 225c2fc9fa9SBarry Smith Y = snes->work[0]; /* work vectors */ 226c2fc9fa9SBarry Smith 2279566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 228c2fc9fa9SBarry Smith snes->iter = 0; 229c2fc9fa9SBarry Smith snes->norm = 0.0; 2309566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 231c2fc9fa9SBarry Smith 2329566063dSJacob Faibussowitsch PetscCall(SNESVIProjectOntoBounds(snes, X)); 2339566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, vi->phi)); 234c2fc9fa9SBarry Smith if (snes->domainerror) { 235c2fc9fa9SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 23622c6f798SBarry Smith sdm->ops->computefunction = vi->computeuserfunction; 237c2fc9fa9SBarry Smith PetscFunctionReturn(0); 238c2fc9fa9SBarry Smith } 239c2fc9fa9SBarry Smith /* Compute Merit function */ 2409566063dSJacob Faibussowitsch PetscCall(SNESVIComputeMeritFunction(vi->phi, &vi->merit, &vi->phinorm)); 241c2fc9fa9SBarry Smith 2429566063dSJacob Faibussowitsch PetscCall(VecNormBegin(X, NORM_2, &xnorm)); /* xnorm <- ||x|| */ 2439566063dSJacob Faibussowitsch PetscCall(VecNormEnd(X, NORM_2, &xnorm)); 244422a814eSBarry Smith SNESCheckFunctionNorm(snes, vi->merit); 245c2fc9fa9SBarry Smith 2469566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 247c2fc9fa9SBarry Smith snes->norm = vi->phinorm; 2489566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2499566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, vi->phinorm, 0)); 2509566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, vi->phinorm)); 251c2fc9fa9SBarry Smith 252c2fc9fa9SBarry Smith /* test convergence */ 253dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, vi->phinorm, &snes->reason, snes->cnvP); 254c2fc9fa9SBarry Smith if (snes->reason) { 25522c6f798SBarry Smith sdm->ops->computefunction = vi->computeuserfunction; 256c2fc9fa9SBarry Smith PetscFunctionReturn(0); 257c2fc9fa9SBarry Smith } 258c2fc9fa9SBarry Smith 259c2fc9fa9SBarry Smith for (i = 0; i < maxits; i++) { 260c2fc9fa9SBarry Smith /* Call general purpose update function */ 261dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 262c2fc9fa9SBarry Smith 263c2fc9fa9SBarry Smith /* Solve J Y = Phi, where J is the semismooth jacobian */ 264b9600128SPeter Brune 265b9600128SPeter Brune /* Get the jacobian -- note that the function must be the original function for snes_fd and snes_fd_color to work for this*/ 26622c6f798SBarry Smith sdm->ops->computefunction = vi->computeuserfunction; 2679566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 26807b62357SFande Kong SNESCheckJacobianDomainerror(snes); 26922c6f798SBarry Smith sdm->ops->computefunction = SNESVIComputeFunction; 270b9600128SPeter Brune 271c2fc9fa9SBarry Smith /* Get the diagonal shift and row scaling vectors */ 2729566063dSJacob Faibussowitsch PetscCall(SNESVIComputeBsubdifferentialVectors(snes, X, F, snes->jacobian, vi->Da, vi->Db)); 273c2fc9fa9SBarry Smith /* Compute the semismooth jacobian */ 2749566063dSJacob Faibussowitsch PetscCall(SNESVIComputeJacobian(snes->jacobian, snes->jacobian_pre, vi->Da, vi->Db)); 275c2fc9fa9SBarry Smith /* Compute the merit function gradient */ 2769566063dSJacob Faibussowitsch PetscCall(SNESVIComputeMeritFunctionGradient(snes->jacobian, vi->phi, vi->dpsi)); 2779566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 2789566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, vi->phi, Y)); 2799566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(snes->ksp, &kspreason)); 280c2fc9fa9SBarry Smith 281c2fc9fa9SBarry Smith if (kspreason < 0) { 282c2fc9fa9SBarry Smith if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 28363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", number linear solve failures %" PetscInt_FMT " greater than current SNES allowed, stopping solve\n", snes->iter, snes->numLinearSolveFailures)); 284c2fc9fa9SBarry Smith snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 285c2fc9fa9SBarry Smith break; 286c2fc9fa9SBarry Smith } 287c2fc9fa9SBarry Smith } 2889566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 289c2fc9fa9SBarry Smith snes->linear_its += lits; 29063a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits)); 291c2fc9fa9SBarry Smith /* 2926b2b7091SBarry Smith if (snes->ops->precheck) { 293c2fc9fa9SBarry Smith PetscBool changed_y = PETSC_FALSE; 294dbbe0bcdSBarry Smith PetscUseTypeMethod(snes,precheck ,X,Y,snes->precheck,&changed_y); 295c2fc9fa9SBarry Smith } 296c2fc9fa9SBarry Smith 2971baa6e33SBarry Smith if (PetscLogPrintInfo) PetscCall(SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W)); 298c2fc9fa9SBarry Smith */ 299c2fc9fa9SBarry Smith /* Compute a (scaled) negative update in the line search routine: 300c2fc9fa9SBarry Smith Y <- X - lambda*Y 301c2fc9fa9SBarry Smith and evaluate G = function(Y) (depends on the line search). 302c2fc9fa9SBarry Smith */ 3039566063dSJacob Faibussowitsch PetscCall(VecCopy(Y, snes->vec_sol_update)); 3049371c9d4SSatish Balay ynorm = 1; 3059371c9d4SSatish Balay gnorm = vi->phinorm; 3069566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(snes->linesearch, X, vi->phi, &gnorm, Y)); 3079566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed)); 3089566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &gnorm, &ynorm)); 3099566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)vi->phinorm, (double)gnorm, (double)ynorm, (int)lssucceed)); 310c2fc9fa9SBarry Smith if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 311c2fc9fa9SBarry Smith if (snes->domainerror) { 312c2fc9fa9SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 31322c6f798SBarry Smith sdm->ops->computefunction = vi->computeuserfunction; 314c2fc9fa9SBarry Smith PetscFunctionReturn(0); 315c2fc9fa9SBarry Smith } 316422a814eSBarry Smith if (lssucceed) { 317c2fc9fa9SBarry Smith if (++snes->numFailures >= snes->maxFailures) { 318c2fc9fa9SBarry Smith PetscBool ismin; 319c2fc9fa9SBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 3209566063dSJacob Faibussowitsch PetscCall(SNESVICheckLocalMin_Private(snes, snes->jacobian, vi->phi, X, gnorm, &ismin)); 321c2fc9fa9SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 322c2fc9fa9SBarry Smith break; 323c2fc9fa9SBarry Smith } 324c2fc9fa9SBarry Smith } 325c2fc9fa9SBarry Smith /* Update function and solution vectors */ 326c2fc9fa9SBarry Smith vi->phinorm = gnorm; 327c2fc9fa9SBarry Smith vi->merit = 0.5 * vi->phinorm * vi->phinorm; 328c2fc9fa9SBarry Smith /* Monitor convergence */ 3299566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 330c2fc9fa9SBarry Smith snes->iter = i + 1; 331c2fc9fa9SBarry Smith snes->norm = vi->phinorm; 332c1e67a49SFande Kong snes->xnorm = xnorm; 333c1e67a49SFande Kong snes->ynorm = ynorm; 3349566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 3359566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 3369566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 337c2fc9fa9SBarry Smith /* Test for convergence, xnorm = || X || */ 3389566063dSJacob Faibussowitsch if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm)); 339dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, vi->phinorm, &snes->reason, snes->cnvP); 340c2fc9fa9SBarry Smith if (snes->reason) break; 341c2fc9fa9SBarry Smith } 342c2fc9fa9SBarry Smith if (i == maxits) { 34363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits)); 344c2fc9fa9SBarry Smith if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 345c2fc9fa9SBarry Smith } 34622c6f798SBarry Smith sdm->ops->computefunction = vi->computeuserfunction; 347c2fc9fa9SBarry Smith PetscFunctionReturn(0); 348c2fc9fa9SBarry Smith } 349c2fc9fa9SBarry Smith 350c2fc9fa9SBarry Smith /* 351f450aa47SBarry Smith SNESSetUp_VINEWTONSSLS - Sets up the internal data structures for the later use 352c2fc9fa9SBarry Smith of the SNES nonlinear solver. 353c2fc9fa9SBarry Smith 354c2fc9fa9SBarry Smith Input Parameter: 355c2fc9fa9SBarry Smith . snes - the SNES context 356c2fc9fa9SBarry Smith 357c2fc9fa9SBarry Smith Application Interface Routine: SNESSetUp() 358c2fc9fa9SBarry Smith 359*f6dfbefdSBarry Smith Note: 360c2fc9fa9SBarry Smith For basic use of the SNES solvers, the user need not explicitly call 361c2fc9fa9SBarry Smith SNESSetUp(), since these actions will automatically occur during 362c2fc9fa9SBarry Smith the call to SNESSolve(). 363c2fc9fa9SBarry Smith */ 3649371c9d4SSatish Balay PetscErrorCode SNESSetUp_VINEWTONSSLS(SNES snes) { 365f450aa47SBarry Smith SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data; 366c2fc9fa9SBarry Smith 367c2fc9fa9SBarry Smith PetscFunctionBegin; 3689566063dSJacob Faibussowitsch PetscCall(SNESSetUp_VI(snes)); 3699566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol, &vi->dpsi)); 3709566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol, &vi->phi)); 3719566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol, &vi->Da)); 3729566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol, &vi->Db)); 3739566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol, &vi->z)); 3749566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol, &vi->t)); 375c2fc9fa9SBarry Smith PetscFunctionReturn(0); 376c2fc9fa9SBarry Smith } 377*f6dfbefdSBarry Smith 3789371c9d4SSatish Balay PetscErrorCode SNESReset_VINEWTONSSLS(SNES snes) { 379f450aa47SBarry Smith SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data; 380c2fc9fa9SBarry Smith 381c2fc9fa9SBarry Smith PetscFunctionBegin; 3829566063dSJacob Faibussowitsch PetscCall(SNESReset_VI(snes)); 3839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vi->dpsi)); 3849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vi->phi)); 3859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vi->Da)); 3869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vi->Db)); 3879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vi->z)); 3889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vi->t)); 389c2fc9fa9SBarry Smith PetscFunctionReturn(0); 390c2fc9fa9SBarry Smith } 391c2fc9fa9SBarry Smith 392c2fc9fa9SBarry Smith /* 393f450aa47SBarry Smith SNESSetFromOptions_VINEWTONSSLS - Sets various parameters for the SNESVI method. 394c2fc9fa9SBarry Smith 395c2fc9fa9SBarry Smith Input Parameter: 396c2fc9fa9SBarry Smith . snes - the SNES context 397c2fc9fa9SBarry Smith 398c2fc9fa9SBarry Smith Application Interface Routine: SNESSetFromOptions() 399c2fc9fa9SBarry Smith */ 4009371c9d4SSatish Balay static PetscErrorCode SNESSetFromOptions_VINEWTONSSLS(SNES snes, PetscOptionItems *PetscOptionsObject) { 401c2fc9fa9SBarry Smith PetscFunctionBegin; 402dbbe0bcdSBarry Smith PetscCall(SNESSetFromOptions_VI(snes, PetscOptionsObject)); 403d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES semismooth method options"); 404d0609cedSBarry Smith PetscOptionsHeadEnd(); 405c2fc9fa9SBarry Smith PetscFunctionReturn(0); 406c2fc9fa9SBarry Smith } 407c2fc9fa9SBarry Smith 408c2fc9fa9SBarry Smith /*MC 409f450aa47SBarry Smith SNESVINEWTONSSLS - Semi-smooth solver for variational inequalities based on Newton's method 410c2fc9fa9SBarry Smith 411*f6dfbefdSBarry Smith Options Database Keys: 412d5f1b7e6SEd Bueler + -snes_type <vinewtonssls,vinewtonrsls> a semi-smooth solver, a reduced space active set method 41361589011SJed Brown - -snes_vi_monitor - prints the number of active constraints at each iteration. 41461589011SJed Brown 415c2fc9fa9SBarry Smith Level: beginner 416c2fc9fa9SBarry Smith 417b80f3ac1SShri Abhyankar References: 418606c0280SSatish Balay + * - T. S. Munson, F. Facchinei, M. C. Ferris, A. Fischer, and C. Kanzow. The semismooth 419b80f3ac1SShri Abhyankar algorithm for large scale complementarity problems. INFORMS Journal on Computing, 13 (2001). 420606c0280SSatish Balay - * - T. S. Munson, and S. Benson. Flexible Complementarity Solvers for Large Scale 421d5f1b7e6SEd Bueler Applications, Optimization Methods and Software, 21 (2006). 422b80f3ac1SShri Abhyankar 423*f6dfbefdSBarry Smith Notes: 424*f6dfbefdSBarry Smith This family of algorithm is much like an interior point method. 425c2fc9fa9SBarry Smith 426*f6dfbefdSBarry Smith The reduced space active set solvers `SNESVINEWTONRSLS` provide an alternative approach that does not result in extremely ill-conditioned linear systems 427*f6dfbefdSBarry Smith 428*f6dfbefdSBarry Smith .seealso: `SNESVINEWTONRSLS`, `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESVINEWTONRSLS`, `SNESNEWTONTR`, `SNESLineSearchSetType()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()` 429c2fc9fa9SBarry Smith M*/ 4309371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_VINEWTONSSLS(SNES snes) { 431f450aa47SBarry Smith SNES_VINEWTONSSLS *vi; 432d8d34be6SBarry Smith SNESLineSearch linesearch; 433c2fc9fa9SBarry Smith 434c2fc9fa9SBarry Smith PetscFunctionBegin; 435f450aa47SBarry Smith snes->ops->reset = SNESReset_VINEWTONSSLS; 436f450aa47SBarry Smith snes->ops->setup = SNESSetUp_VINEWTONSSLS; 437f450aa47SBarry Smith snes->ops->solve = SNESSolve_VINEWTONSSLS; 438c2fc9fa9SBarry Smith snes->ops->destroy = SNESDestroy_VI; 439f450aa47SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_VINEWTONSSLS; 4400298fd71SBarry Smith snes->ops->view = NULL; 441c2fc9fa9SBarry Smith 442c2fc9fa9SBarry Smith snes->usesksp = PETSC_TRUE; 443efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 444c2fc9fa9SBarry Smith 4459566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 446ec786807SJed Brown if (!((PetscObject)linesearch)->type_name) { 4479566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT)); 4489566063dSJacob Faibussowitsch PetscCall(SNESLineSearchBTSetAlpha(linesearch, 0.0)); 449ec786807SJed Brown } 450d8d34be6SBarry Smith 4514fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 4524fc747eaSLawrence Mitchell 4539566063dSJacob Faibussowitsch PetscCall(PetscNewLog(snes, &vi)); 454c2fc9fa9SBarry Smith snes->data = (void *)vi; 455c2fc9fa9SBarry Smith 4569566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", SNESVISetVariableBounds_VI)); 4579566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", SNESVISetComputeVariableBounds_VI)); 458c2fc9fa9SBarry Smith PetscFunctionReturn(0); 459c2fc9fa9SBarry Smith } 460