1c2fc9fa9SBarry Smith #include <../src/snes/impls/vi/ss/vissimpl.h> /*I "petscsnes.h" I*/ 2c2fc9fa9SBarry Smith 3f6dfbefdSBarry Smith /*@ 4c2fc9fa9SBarry Smith SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem. 5c2fc9fa9SBarry Smith 6c2fc9fa9SBarry Smith Input Parameter: 7f6dfbefdSBarry Smith . phi - the `Vec` holding the evaluation of the semismooth function 8c2fc9fa9SBarry Smith 9f6dfbefdSBarry Smith Output Parameters: 10f6dfbefdSBarry Smith + merit - the merit function 1/2 ||phi||^2 11f6dfbefdSBarry Smith - phinorm - the two-norm of the vector, ||phi|| 12c2fc9fa9SBarry Smith 13fbda9744SBarry Smith Level: developer 14fbda9744SBarry Smith 15420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESVINEWTONSSLS`, `SNESVIComputeFunction()` 16f6dfbefdSBarry Smith @*/ 17d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal *merit, PetscReal *phinorm) 18d71ae5a4SJacob Faibussowitsch { 19c2fc9fa9SBarry Smith PetscFunctionBegin; 209566063dSJacob Faibussowitsch PetscCall(VecNormBegin(phi, NORM_2, phinorm)); 219566063dSJacob Faibussowitsch PetscCall(VecNormEnd(phi, NORM_2, phinorm)); 22c2fc9fa9SBarry Smith *merit = 0.5 * (*phinorm) * (*phinorm); 233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 24c2fc9fa9SBarry Smith } 25c2fc9fa9SBarry Smith 26d71ae5a4SJacob Faibussowitsch static inline PetscScalar Phi(PetscScalar a, PetscScalar b) 27d71ae5a4SJacob Faibussowitsch { 28c2fc9fa9SBarry Smith return a + b - PetscSqrtScalar(a * a + b * b); 29c2fc9fa9SBarry Smith } 30c2fc9fa9SBarry Smith 31d71ae5a4SJacob Faibussowitsch static inline PetscScalar DPhi(PetscScalar a, PetscScalar b) 32d71ae5a4SJacob Faibussowitsch { 33c2fc9fa9SBarry Smith if ((PetscAbsScalar(a) >= 1.e-6) || (PetscAbsScalar(b) >= 1.e-6)) return 1.0 - a / PetscSqrtScalar(a * a + b * b); 34c2fc9fa9SBarry Smith else return .5; 35c2fc9fa9SBarry Smith } 36c2fc9fa9SBarry Smith 37f6dfbefdSBarry Smith /*@ 38f6dfbefdSBarry Smith SNESVIComputeFunction - Provides the function that reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear 39f6dfbefdSBarry Smith equations in semismooth form. 40c2fc9fa9SBarry Smith 41c2fc9fa9SBarry Smith Input Parameters: 42420bcc1bSBarry Smith + snes - the `SNES` context 43d5f1b7e6SEd Bueler . X - current iterate 44f6dfbefdSBarry Smith - functx - user defined function context 45c2fc9fa9SBarry Smith 46f6dfbefdSBarry Smith Output Parameter: 47420bcc1bSBarry Smith . phi - the evaluation of semismooth function at `X` 48c2fc9fa9SBarry Smith 49fbda9744SBarry Smith Level: developer 50fbda9744SBarry Smith 51420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESVINEWTONSSLS`, `SNESVIComputeMeritFunction()` 52f6dfbefdSBarry Smith @*/ 53d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIComputeFunction(SNES snes, Vec X, Vec phi, void *functx) 54d71ae5a4SJacob Faibussowitsch { 55f450aa47SBarry Smith SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data; 56c2fc9fa9SBarry Smith Vec Xl = snes->xl, Xu = snes->xu, F = snes->vec_func; 575edff71fSBarry Smith PetscScalar *phi_arr, *f_arr, *l, *u; 585edff71fSBarry Smith const PetscScalar *x_arr; 59c2fc9fa9SBarry Smith PetscInt i, nlocal; 60c2fc9fa9SBarry Smith 61c2fc9fa9SBarry Smith PetscFunctionBegin; 629566063dSJacob Faibussowitsch PetscCall((*vi->computeuserfunction)(snes, X, F, functx)); 639566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &nlocal)); 649566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x_arr)); 659566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f_arr)); 669566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xl, &l)); 679566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xu, &u)); 689566063dSJacob Faibussowitsch PetscCall(VecGetArray(phi, &phi_arr)); 69c2fc9fa9SBarry Smith 70c2fc9fa9SBarry Smith for (i = 0; i < nlocal; i++) { 71e270355aSBarry Smith if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { /* no constraints on variable */ 72c2fc9fa9SBarry Smith phi_arr[i] = f_arr[i]; 73e270355aSBarry Smith } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { /* upper bound on variable only */ 74c2fc9fa9SBarry Smith phi_arr[i] = -Phi(u[i] - x_arr[i], -f_arr[i]); 75e270355aSBarry Smith } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { /* lower bound on variable only */ 76c2fc9fa9SBarry Smith phi_arr[i] = Phi(x_arr[i] - l[i], f_arr[i]); 77c2fc9fa9SBarry Smith } else if (l[i] == u[i]) { 78c2fc9fa9SBarry Smith phi_arr[i] = l[i] - x_arr[i]; 79c2fc9fa9SBarry Smith } else { /* both bounds on variable */ 80c2fc9fa9SBarry Smith phi_arr[i] = Phi(x_arr[i] - l[i], -Phi(u[i] - x_arr[i], -f_arr[i])); 81c2fc9fa9SBarry Smith } 82c2fc9fa9SBarry Smith } 83c2fc9fa9SBarry Smith 849566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x_arr)); 859566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f_arr)); 869566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xl, &l)); 879566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xu, &u)); 889566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(phi, &phi_arr)); 893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 90c2fc9fa9SBarry Smith } 91c2fc9fa9SBarry Smith 92c2fc9fa9SBarry Smith /* 93c2fc9fa9SBarry Smith SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the 94c2fc9fa9SBarry Smith the semismooth jacobian. 95c2fc9fa9SBarry Smith */ 9666976f2fSJacob Faibussowitsch static PetscErrorCode SNESVIComputeBsubdifferentialVectors(SNES snes, Vec X, Vec F, Mat jac, Vec Da, Vec Db) 97d71ae5a4SJacob Faibussowitsch { 98c2fc9fa9SBarry Smith PetscScalar *l, *u, *x, *f, *da, *db, da1, da2, db1, db2; 99c2fc9fa9SBarry Smith PetscInt i, nlocal; 100c2fc9fa9SBarry Smith 101c2fc9fa9SBarry Smith PetscFunctionBegin; 1029566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 1039566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 1049566063dSJacob Faibussowitsch PetscCall(VecGetArray(snes->xl, &l)); 1059566063dSJacob Faibussowitsch PetscCall(VecGetArray(snes->xu, &u)); 1069566063dSJacob Faibussowitsch PetscCall(VecGetArray(Da, &da)); 1079566063dSJacob Faibussowitsch PetscCall(VecGetArray(Db, &db)); 1089566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &nlocal)); 109c2fc9fa9SBarry Smith 110c2fc9fa9SBarry Smith for (i = 0; i < nlocal; i++) { 111e270355aSBarry Smith if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { /* no constraints on variable */ 112c2fc9fa9SBarry Smith da[i] = 0; 113c2fc9fa9SBarry Smith db[i] = 1; 114e270355aSBarry Smith } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { /* upper bound on variable only */ 115c2fc9fa9SBarry Smith da[i] = DPhi(u[i] - x[i], -f[i]); 116c2fc9fa9SBarry Smith db[i] = DPhi(-f[i], u[i] - x[i]); 117e270355aSBarry Smith } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { /* lower bound on variable only */ 118c2fc9fa9SBarry Smith da[i] = DPhi(x[i] - l[i], f[i]); 119c2fc9fa9SBarry Smith db[i] = DPhi(f[i], x[i] - l[i]); 120c2fc9fa9SBarry Smith } else if (l[i] == u[i]) { /* fixed variable */ 121c2fc9fa9SBarry Smith da[i] = 1; 122c2fc9fa9SBarry Smith db[i] = 0; 123c2fc9fa9SBarry Smith } else { /* upper and lower bounds on variable */ 124c2fc9fa9SBarry Smith da1 = DPhi(x[i] - l[i], -Phi(u[i] - x[i], -f[i])); 125c2fc9fa9SBarry Smith db1 = DPhi(-Phi(u[i] - x[i], -f[i]), x[i] - l[i]); 126c2fc9fa9SBarry Smith da2 = DPhi(u[i] - x[i], -f[i]); 127c2fc9fa9SBarry Smith db2 = DPhi(-f[i], u[i] - x[i]); 128c2fc9fa9SBarry Smith da[i] = da1 + db1 * da2; 129c2fc9fa9SBarry Smith db[i] = db1 * db2; 130c2fc9fa9SBarry Smith } 131c2fc9fa9SBarry Smith } 132c2fc9fa9SBarry Smith 1339566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 1349566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 1359566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(snes->xl, &l)); 1369566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(snes->xu, &u)); 1379566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Da, &da)); 1389566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Db, &db)); 1393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 140c2fc9fa9SBarry Smith } 141c2fc9fa9SBarry Smith 142c2fc9fa9SBarry Smith /* 143c2fc9fa9SBarry 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. 144c2fc9fa9SBarry Smith 145c2fc9fa9SBarry Smith Input Parameters: 1467addb90fSBarry Smith . Da - Diagonal shift vector for the semismooth Jacobian. 1477addb90fSBarry Smith . Db - Row scaling vector for the semismooth Jacobian. 148c2fc9fa9SBarry Smith 149c2fc9fa9SBarry Smith Output Parameters: 1507addb90fSBarry Smith . jac - semismooth Jacobian 1517addb90fSBarry Smith . jac_pre - optional matrix from which to construct the preconditioner 152c2fc9fa9SBarry Smith 153f6dfbefdSBarry Smith Note: 154c2fc9fa9SBarry Smith The semismooth jacobian matrix is given by 1557addb90fSBarry Smith $ jac = Da + Db*jacfun $ 1567addb90fSBarry Smith where `Db` is the row scaling matrix stored as a vector, 1577addb90fSBarry Smith `Da` is the diagonal perturbation matrix stored as a vector 1587addb90fSBarry Smith and `jacfun` is the Jacobian of the original nonlinear function. 159c2fc9fa9SBarry Smith */ 16066976f2fSJacob Faibussowitsch static PetscErrorCode SNESVIComputeJacobian(Mat jac, Mat jac_pre, Vec Da, Vec Db) 161d71ae5a4SJacob Faibussowitsch { 162c2fc9fa9SBarry Smith /* Do row scaling and add diagonal perturbation */ 163362febeeSStefano Zampini PetscFunctionBegin; 1649566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(jac, Db, NULL)); 1659566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(jac, Da, ADD_VALUES)); 166c2fc9fa9SBarry Smith if (jac != jac_pre) { /* If jac and jac_pre are different */ 1679566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(jac_pre, Db, NULL)); 1689566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(jac_pre, Da, ADD_VALUES)); 169c2fc9fa9SBarry Smith } 1703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 171c2fc9fa9SBarry Smith } 172c2fc9fa9SBarry Smith 173ceaaa498SBarry Smith // PetscClangLinter pragma disable: -fdoc-sowing-chars 174c2fc9fa9SBarry Smith /* 175c2fc9fa9SBarry Smith SNESVIComputeMeritFunctionGradient - Computes the gradient of the merit function psi. 176c2fc9fa9SBarry Smith 177c2fc9fa9SBarry Smith Input Parameters: 178c2fc9fa9SBarry Smith phi - semismooth function. 179c2fc9fa9SBarry Smith H - semismooth jacobian 180c2fc9fa9SBarry Smith 181f6dfbefdSBarry Smith Output Parameter: 182c2fc9fa9SBarry Smith dpsi - merit function gradient 183c2fc9fa9SBarry Smith 184f6dfbefdSBarry Smith Note: 185c2fc9fa9SBarry Smith The merit function gradient is computed as follows 186c2fc9fa9SBarry Smith dpsi = H^T*phi 187c2fc9fa9SBarry Smith */ 188ceaaa498SBarry Smith static PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi) 189d71ae5a4SJacob Faibussowitsch { 190c2fc9fa9SBarry Smith PetscFunctionBegin; 1919566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(H, phi, dpsi)); 1923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 193c2fc9fa9SBarry Smith } 194c2fc9fa9SBarry Smith 19566976f2fSJacob Faibussowitsch static PetscErrorCode SNESSolve_VINEWTONSSLS(SNES snes) 196d71ae5a4SJacob Faibussowitsch { 197f450aa47SBarry Smith SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data; 198c2fc9fa9SBarry Smith PetscInt maxits, i, lits; 199422a814eSBarry Smith SNESLineSearchReason lssucceed; 200c2fc9fa9SBarry Smith PetscReal gnorm, xnorm = 0, ynorm; 2019bd66eb0SPeter Brune Vec Y, X, F; 202c2fc9fa9SBarry Smith KSPConvergedReason kspreason; 2036cab3a1bSJed Brown DM dm; 204942e3340SBarry Smith DMSNES sdm; 205c2fc9fa9SBarry Smith 206c2fc9fa9SBarry Smith PetscFunctionBegin; 2079566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 2089566063dSJacob Faibussowitsch PetscCall(DMGetDMSNES(dm, &sdm)); 2091aa26658SKarl Rupp 21022c6f798SBarry Smith vi->computeuserfunction = sdm->ops->computefunction; 21122c6f798SBarry Smith sdm->ops->computefunction = SNESVIComputeFunction; 212c2fc9fa9SBarry Smith 213c2fc9fa9SBarry Smith snes->numFailures = 0; 214c2fc9fa9SBarry Smith snes->numLinearSolveFailures = 0; 215c2fc9fa9SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 216c2fc9fa9SBarry Smith 217c2fc9fa9SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 218c2fc9fa9SBarry Smith X = snes->vec_sol; /* solution vector */ 219c2fc9fa9SBarry Smith F = snes->vec_func; /* residual vector */ 220c2fc9fa9SBarry Smith Y = snes->work[0]; /* work vectors */ 221c2fc9fa9SBarry Smith 2229566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 223c2fc9fa9SBarry Smith snes->iter = 0; 224c2fc9fa9SBarry Smith snes->norm = 0.0; 2259566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 226c2fc9fa9SBarry Smith 2279566063dSJacob Faibussowitsch PetscCall(SNESVIProjectOntoBounds(snes, X)); 2289566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, vi->phi)); 229c2fc9fa9SBarry Smith if (snes->domainerror) { 230c2fc9fa9SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 23122c6f798SBarry Smith sdm->ops->computefunction = vi->computeuserfunction; 2323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 233c2fc9fa9SBarry Smith } 234c2fc9fa9SBarry Smith /* Compute Merit function */ 2359566063dSJacob Faibussowitsch PetscCall(SNESVIComputeMeritFunction(vi->phi, &vi->merit, &vi->phinorm)); 236c2fc9fa9SBarry Smith 2379566063dSJacob Faibussowitsch PetscCall(VecNormBegin(X, NORM_2, &xnorm)); /* xnorm <- ||x|| */ 2389566063dSJacob Faibussowitsch PetscCall(VecNormEnd(X, NORM_2, &xnorm)); 239422a814eSBarry Smith SNESCheckFunctionNorm(snes, vi->merit); 240c2fc9fa9SBarry Smith 2419566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 242c2fc9fa9SBarry Smith snes->norm = vi->phinorm; 2439566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2449566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, vi->phinorm, 0)); 245c2fc9fa9SBarry Smith 246c2fc9fa9SBarry Smith /* test convergence */ 2472d157150SStefano Zampini PetscCall(SNESConverged(snes, 0, 0.0, 0.0, vi->phinorm)); 2482d157150SStefano Zampini PetscCall(SNESMonitor(snes, 0, vi->phinorm)); 249c2fc9fa9SBarry Smith if (snes->reason) { 25022c6f798SBarry Smith sdm->ops->computefunction = vi->computeuserfunction; 2513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 252c2fc9fa9SBarry Smith } 253c2fc9fa9SBarry Smith 254c2fc9fa9SBarry Smith for (i = 0; i < maxits; i++) { 255c2fc9fa9SBarry Smith /* Call general purpose update function */ 256dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 257c2fc9fa9SBarry Smith 258c2fc9fa9SBarry Smith /* Solve J Y = Phi, where J is the semismooth jacobian */ 259b9600128SPeter Brune 260b9600128SPeter Brune /* Get the jacobian -- note that the function must be the original function for snes_fd and snes_fd_color to work for this*/ 26122c6f798SBarry Smith sdm->ops->computefunction = vi->computeuserfunction; 2629566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 26307b62357SFande Kong SNESCheckJacobianDomainerror(snes); 26422c6f798SBarry Smith sdm->ops->computefunction = SNESVIComputeFunction; 265b9600128SPeter Brune 266c2fc9fa9SBarry Smith /* Get the diagonal shift and row scaling vectors */ 2679566063dSJacob Faibussowitsch PetscCall(SNESVIComputeBsubdifferentialVectors(snes, X, F, snes->jacobian, vi->Da, vi->Db)); 268c2fc9fa9SBarry Smith /* Compute the semismooth jacobian */ 2699566063dSJacob Faibussowitsch PetscCall(SNESVIComputeJacobian(snes->jacobian, snes->jacobian_pre, vi->Da, vi->Db)); 270c2fc9fa9SBarry Smith /* Compute the merit function gradient */ 2719566063dSJacob Faibussowitsch PetscCall(SNESVIComputeMeritFunctionGradient(snes->jacobian, vi->phi, vi->dpsi)); 2729566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 2739566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, vi->phi, Y)); 2749566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(snes->ksp, &kspreason)); 275c2fc9fa9SBarry Smith 276c2fc9fa9SBarry Smith if (kspreason < 0) { 277c2fc9fa9SBarry Smith if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 27863a3b9bcSJacob 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)); 279c2fc9fa9SBarry Smith snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 280c2fc9fa9SBarry Smith break; 281c2fc9fa9SBarry Smith } 282c2fc9fa9SBarry Smith } 2839566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 284c2fc9fa9SBarry Smith snes->linear_its += lits; 28563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits)); 286c2fc9fa9SBarry Smith /* 2876b2b7091SBarry Smith if (snes->ops->precheck) { 288c2fc9fa9SBarry Smith PetscBool changed_y = PETSC_FALSE; 289dbbe0bcdSBarry Smith PetscUseTypeMethod(snes,precheck ,X,Y,snes->precheck,&changed_y); 290c2fc9fa9SBarry Smith } 291c2fc9fa9SBarry Smith 2921baa6e33SBarry Smith if (PetscLogPrintInfo) PetscCall(SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W)); 293c2fc9fa9SBarry Smith */ 294c2fc9fa9SBarry Smith /* Compute a (scaled) negative update in the line search routine: 295c2fc9fa9SBarry Smith Y <- X - lambda*Y 296c2fc9fa9SBarry Smith and evaluate G = function(Y) (depends on the line search). 297c2fc9fa9SBarry Smith */ 2989566063dSJacob Faibussowitsch PetscCall(VecCopy(Y, snes->vec_sol_update)); 2999371c9d4SSatish Balay ynorm = 1; 3009371c9d4SSatish Balay gnorm = vi->phinorm; 3019566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(snes->linesearch, X, vi->phi, &gnorm, Y)); 3029566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed)); 3039566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &gnorm, &ynorm)); 3049566063dSJacob 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)); 305c2fc9fa9SBarry Smith if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 306c2fc9fa9SBarry Smith if (snes->domainerror) { 307c2fc9fa9SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 30822c6f798SBarry Smith sdm->ops->computefunction = vi->computeuserfunction; 3093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 310c2fc9fa9SBarry Smith } 311422a814eSBarry Smith if (lssucceed) { 312c2fc9fa9SBarry Smith if (++snes->numFailures >= snes->maxFailures) { 313c2fc9fa9SBarry Smith PetscBool ismin; 314c2fc9fa9SBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 3159566063dSJacob Faibussowitsch PetscCall(SNESVICheckLocalMin_Private(snes, snes->jacobian, vi->phi, X, gnorm, &ismin)); 316c2fc9fa9SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 317c2fc9fa9SBarry Smith break; 318c2fc9fa9SBarry Smith } 319c2fc9fa9SBarry Smith } 320c2fc9fa9SBarry Smith /* Update function and solution vectors */ 321c2fc9fa9SBarry Smith vi->phinorm = gnorm; 322c2fc9fa9SBarry Smith vi->merit = 0.5 * vi->phinorm * vi->phinorm; 323c2fc9fa9SBarry Smith /* Monitor convergence */ 3249566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 325c2fc9fa9SBarry Smith snes->iter = i + 1; 326c2fc9fa9SBarry Smith snes->norm = vi->phinorm; 327c1e67a49SFande Kong snes->xnorm = xnorm; 328c1e67a49SFande Kong snes->ynorm = ynorm; 3299566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 3309566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 331c2fc9fa9SBarry Smith /* Test for convergence, xnorm = || X || */ 3329566063dSJacob Faibussowitsch if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm)); 3332d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, vi->phinorm)); 3342d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 335c2fc9fa9SBarry Smith if (snes->reason) break; 336c2fc9fa9SBarry Smith } 33722c6f798SBarry Smith sdm->ops->computefunction = vi->computeuserfunction; 3383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 339c2fc9fa9SBarry Smith } 340c2fc9fa9SBarry Smith 34166976f2fSJacob Faibussowitsch static PetscErrorCode SNESSetUp_VINEWTONSSLS(SNES snes) 342d71ae5a4SJacob Faibussowitsch { 343f450aa47SBarry Smith SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data; 344c2fc9fa9SBarry Smith 345c2fc9fa9SBarry Smith PetscFunctionBegin; 3469566063dSJacob Faibussowitsch PetscCall(SNESSetUp_VI(snes)); 347*b2ccae6bSStefano Zampini PetscCall(VecDuplicate(snes->work[0], &vi->dpsi)); 348*b2ccae6bSStefano Zampini PetscCall(VecDuplicate(snes->work[0], &vi->phi)); 349*b2ccae6bSStefano Zampini PetscCall(VecDuplicate(snes->work[0], &vi->Da)); 350*b2ccae6bSStefano Zampini PetscCall(VecDuplicate(snes->work[0], &vi->Db)); 351*b2ccae6bSStefano Zampini PetscCall(VecDuplicate(snes->work[0], &vi->z)); 352*b2ccae6bSStefano Zampini PetscCall(VecDuplicate(snes->work[0], &vi->t)); 3533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 354c2fc9fa9SBarry Smith } 355f6dfbefdSBarry Smith 35666976f2fSJacob Faibussowitsch static PetscErrorCode SNESReset_VINEWTONSSLS(SNES snes) 357d71ae5a4SJacob Faibussowitsch { 358f450aa47SBarry Smith SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data; 359c2fc9fa9SBarry Smith 360c2fc9fa9SBarry Smith PetscFunctionBegin; 3619566063dSJacob Faibussowitsch PetscCall(SNESReset_VI(snes)); 3629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vi->dpsi)); 3639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vi->phi)); 3649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vi->Da)); 3659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vi->Db)); 3669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vi->z)); 3679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vi->t)); 3683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 369c2fc9fa9SBarry Smith } 370c2fc9fa9SBarry Smith 371ce78bad3SBarry Smith static PetscErrorCode SNESSetFromOptions_VINEWTONSSLS(SNES snes, PetscOptionItems PetscOptionsObject) 372d71ae5a4SJacob Faibussowitsch { 373c2fc9fa9SBarry Smith PetscFunctionBegin; 374dbbe0bcdSBarry Smith PetscCall(SNESSetFromOptions_VI(snes, PetscOptionsObject)); 375d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES semismooth method options"); 376d0609cedSBarry Smith PetscOptionsHeadEnd(); 3773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 378c2fc9fa9SBarry Smith } 379c2fc9fa9SBarry Smith 380c2fc9fa9SBarry Smith /*MC 381f450aa47SBarry Smith SNESVINEWTONSSLS - Semi-smooth solver for variational inequalities based on Newton's method 382c2fc9fa9SBarry Smith 383f6dfbefdSBarry Smith Options Database Keys: 384d5f1b7e6SEd Bueler + -snes_type <vinewtonssls,vinewtonrsls> a semi-smooth solver, a reduced space active set method 38561589011SJed Brown - -snes_vi_monitor - prints the number of active constraints at each iteration. 38661589011SJed Brown 387c2fc9fa9SBarry Smith Level: beginner 388c2fc9fa9SBarry Smith 389420bcc1bSBarry Smith Notes: 390420bcc1bSBarry Smith This family of algorithms is much like an interior point method. 391420bcc1bSBarry Smith 392420bcc1bSBarry Smith The reduced space active set solvers `SNESVINEWTONRSLS` provide an alternative approach that does not result in extremely ill-conditioned linear systems 393420bcc1bSBarry Smith 3941d27aa22SBarry Smith See {cite}`munson.facchinei.ea:semismooth` and {cite}`benson2006flexible` 395b80f3ac1SShri Abhyankar 396420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESVINEWTONRSLS`, `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESVINEWTONRSLS`, `SNESNEWTONTR`, `SNESLineSearchSetType()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()` 397c2fc9fa9SBarry Smith M*/ 398d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_VINEWTONSSLS(SNES snes) 399d71ae5a4SJacob Faibussowitsch { 400f450aa47SBarry Smith SNES_VINEWTONSSLS *vi; 401d8d34be6SBarry Smith SNESLineSearch linesearch; 402c2fc9fa9SBarry Smith 403c2fc9fa9SBarry Smith PetscFunctionBegin; 404f450aa47SBarry Smith snes->ops->reset = SNESReset_VINEWTONSSLS; 405f450aa47SBarry Smith snes->ops->setup = SNESSetUp_VINEWTONSSLS; 406f450aa47SBarry Smith snes->ops->solve = SNESSolve_VINEWTONSSLS; 407c2fc9fa9SBarry Smith snes->ops->destroy = SNESDestroy_VI; 408f450aa47SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_VINEWTONSSLS; 4090298fd71SBarry Smith snes->ops->view = NULL; 410c2fc9fa9SBarry Smith 411c2fc9fa9SBarry Smith snes->usesksp = PETSC_TRUE; 412efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 413c2fc9fa9SBarry Smith 4149566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 415ec786807SJed Brown if (!((PetscObject)linesearch)->type_name) { 4169566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT)); 4179566063dSJacob Faibussowitsch PetscCall(SNESLineSearchBTSetAlpha(linesearch, 0.0)); 418ec786807SJed Brown } 419d8d34be6SBarry Smith 4204fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 4214fc747eaSLawrence Mitchell 42277e5a1f9SBarry Smith PetscCall(SNESParametersInitialize(snes)); 42377e5a1f9SBarry Smith 4244dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&vi)); 425c2fc9fa9SBarry Smith snes->data = (void *)vi; 426c2fc9fa9SBarry Smith 4279566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", SNESVISetVariableBounds_VI)); 4289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", SNESVISetComputeVariableBounds_VI)); 4293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 430c2fc9fa9SBarry Smith } 431