1c2fc9fa9SBarry Smith 2c2fc9fa9SBarry Smith #include <../src/snes/impls/vi/ss/vissimpl.h> /*I "petscsnes.h" I*/ 3c2fc9fa9SBarry Smith 4f6dfbefdSBarry Smith /*@ 5c2fc9fa9SBarry Smith SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem. 6c2fc9fa9SBarry Smith 7c2fc9fa9SBarry Smith Input Parameter: 8f6dfbefdSBarry Smith . phi - the `Vec` holding the evaluation of the semismooth function 9c2fc9fa9SBarry Smith 10f6dfbefdSBarry Smith Output Parameters: 11f6dfbefdSBarry Smith + merit - the merit function 1/2 ||phi||^2 12f6dfbefdSBarry Smith - phinorm - the two-norm of the vector, ||phi|| 13c2fc9fa9SBarry Smith 14fbda9744SBarry Smith Level: developer 15fbda9744SBarry Smith 16f6dfbefdSBarry Smith .seealso: `SNESVINEWTONSSLS`, `SNESVIComputeFunction()` 17f6dfbefdSBarry Smith @*/ 18d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal *merit, PetscReal *phinorm) 19d71ae5a4SJacob Faibussowitsch { 20c2fc9fa9SBarry Smith PetscFunctionBegin; 219566063dSJacob Faibussowitsch PetscCall(VecNormBegin(phi, NORM_2, phinorm)); 229566063dSJacob Faibussowitsch PetscCall(VecNormEnd(phi, NORM_2, phinorm)); 23c2fc9fa9SBarry Smith *merit = 0.5 * (*phinorm) * (*phinorm); 24*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 25c2fc9fa9SBarry Smith } 26c2fc9fa9SBarry Smith 27d71ae5a4SJacob Faibussowitsch static inline PetscScalar Phi(PetscScalar a, PetscScalar b) 28d71ae5a4SJacob Faibussowitsch { 29c2fc9fa9SBarry Smith return a + b - PetscSqrtScalar(a * a + b * b); 30c2fc9fa9SBarry Smith } 31c2fc9fa9SBarry Smith 32d71ae5a4SJacob Faibussowitsch static inline PetscScalar DPhi(PetscScalar a, PetscScalar b) 33d71ae5a4SJacob Faibussowitsch { 34c2fc9fa9SBarry Smith if ((PetscAbsScalar(a) >= 1.e-6) || (PetscAbsScalar(b) >= 1.e-6)) return 1.0 - a / PetscSqrtScalar(a * a + b * b); 35c2fc9fa9SBarry Smith else return .5; 36c2fc9fa9SBarry Smith } 37c2fc9fa9SBarry Smith 38f6dfbefdSBarry Smith /*@ 39f6dfbefdSBarry Smith SNESVIComputeFunction - Provides the function that reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear 40f6dfbefdSBarry Smith equations in semismooth form. 41c2fc9fa9SBarry Smith 42c2fc9fa9SBarry Smith Input Parameters: 43f6dfbefdSBarry Smith + snes - the SNES context 44d5f1b7e6SEd Bueler . X - current iterate 45f6dfbefdSBarry Smith - functx - user defined function context 46c2fc9fa9SBarry Smith 47f6dfbefdSBarry Smith Output Parameter: 48f6dfbefdSBarry Smith . phi - the evaluation of Semismooth function at X 49c2fc9fa9SBarry Smith 50fbda9744SBarry Smith Level: developer 51fbda9744SBarry Smith 52f6dfbefdSBarry Smith .seealso: `SNESVINEWTONSSLS`, `SNESVIComputeMeritFunction()` 53f6dfbefdSBarry Smith @*/ 54d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIComputeFunction(SNES snes, Vec X, Vec phi, void *functx) 55d71ae5a4SJacob Faibussowitsch { 56f450aa47SBarry Smith SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data; 57c2fc9fa9SBarry Smith Vec Xl = snes->xl, Xu = snes->xu, F = snes->vec_func; 585edff71fSBarry Smith PetscScalar *phi_arr, *f_arr, *l, *u; 595edff71fSBarry Smith const PetscScalar *x_arr; 60c2fc9fa9SBarry Smith PetscInt i, nlocal; 61c2fc9fa9SBarry Smith 62c2fc9fa9SBarry Smith PetscFunctionBegin; 639566063dSJacob Faibussowitsch PetscCall((*vi->computeuserfunction)(snes, X, F, functx)); 649566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &nlocal)); 659566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x_arr)); 669566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f_arr)); 679566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xl, &l)); 689566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xu, &u)); 699566063dSJacob Faibussowitsch PetscCall(VecGetArray(phi, &phi_arr)); 70c2fc9fa9SBarry Smith 71c2fc9fa9SBarry Smith for (i = 0; i < nlocal; i++) { 72e270355aSBarry Smith if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { /* no constraints on variable */ 73c2fc9fa9SBarry Smith phi_arr[i] = f_arr[i]; 74e270355aSBarry Smith } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { /* upper bound on variable only */ 75c2fc9fa9SBarry Smith phi_arr[i] = -Phi(u[i] - x_arr[i], -f_arr[i]); 76e270355aSBarry Smith } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { /* lower bound on variable only */ 77c2fc9fa9SBarry Smith phi_arr[i] = Phi(x_arr[i] - l[i], f_arr[i]); 78c2fc9fa9SBarry Smith } else if (l[i] == u[i]) { 79c2fc9fa9SBarry Smith phi_arr[i] = l[i] - x_arr[i]; 80c2fc9fa9SBarry Smith } else { /* both bounds on variable */ 81c2fc9fa9SBarry Smith phi_arr[i] = Phi(x_arr[i] - l[i], -Phi(u[i] - x_arr[i], -f_arr[i])); 82c2fc9fa9SBarry Smith } 83c2fc9fa9SBarry Smith } 84c2fc9fa9SBarry Smith 859566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x_arr)); 869566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f_arr)); 879566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xl, &l)); 889566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xu, &u)); 899566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(phi, &phi_arr)); 90*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 91c2fc9fa9SBarry Smith } 92c2fc9fa9SBarry Smith 93c2fc9fa9SBarry Smith /* 94c2fc9fa9SBarry Smith SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the 95c2fc9fa9SBarry Smith the semismooth jacobian. 96c2fc9fa9SBarry Smith */ 97d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIComputeBsubdifferentialVectors(SNES snes, Vec X, Vec F, Mat jac, Vec Da, Vec Db) 98d71ae5a4SJacob Faibussowitsch { 99c2fc9fa9SBarry Smith PetscScalar *l, *u, *x, *f, *da, *db, da1, da2, db1, db2; 100c2fc9fa9SBarry Smith PetscInt i, nlocal; 101c2fc9fa9SBarry Smith 102c2fc9fa9SBarry Smith PetscFunctionBegin; 1039566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 1049566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 1059566063dSJacob Faibussowitsch PetscCall(VecGetArray(snes->xl, &l)); 1069566063dSJacob Faibussowitsch PetscCall(VecGetArray(snes->xu, &u)); 1079566063dSJacob Faibussowitsch PetscCall(VecGetArray(Da, &da)); 1089566063dSJacob Faibussowitsch PetscCall(VecGetArray(Db, &db)); 1099566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &nlocal)); 110c2fc9fa9SBarry Smith 111c2fc9fa9SBarry Smith for (i = 0; i < nlocal; i++) { 112e270355aSBarry Smith if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { /* no constraints on variable */ 113c2fc9fa9SBarry Smith da[i] = 0; 114c2fc9fa9SBarry Smith db[i] = 1; 115e270355aSBarry Smith } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { /* upper bound on variable only */ 116c2fc9fa9SBarry Smith da[i] = DPhi(u[i] - x[i], -f[i]); 117c2fc9fa9SBarry Smith db[i] = DPhi(-f[i], u[i] - x[i]); 118e270355aSBarry Smith } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { /* lower bound on variable only */ 119c2fc9fa9SBarry Smith da[i] = DPhi(x[i] - l[i], f[i]); 120c2fc9fa9SBarry Smith db[i] = DPhi(f[i], x[i] - l[i]); 121c2fc9fa9SBarry Smith } else if (l[i] == u[i]) { /* fixed variable */ 122c2fc9fa9SBarry Smith da[i] = 1; 123c2fc9fa9SBarry Smith db[i] = 0; 124c2fc9fa9SBarry Smith } else { /* upper and lower bounds on variable */ 125c2fc9fa9SBarry Smith da1 = DPhi(x[i] - l[i], -Phi(u[i] - x[i], -f[i])); 126c2fc9fa9SBarry Smith db1 = DPhi(-Phi(u[i] - x[i], -f[i]), x[i] - l[i]); 127c2fc9fa9SBarry Smith da2 = DPhi(u[i] - x[i], -f[i]); 128c2fc9fa9SBarry Smith db2 = DPhi(-f[i], u[i] - x[i]); 129c2fc9fa9SBarry Smith da[i] = da1 + db1 * da2; 130c2fc9fa9SBarry Smith db[i] = db1 * db2; 131c2fc9fa9SBarry Smith } 132c2fc9fa9SBarry Smith } 133c2fc9fa9SBarry Smith 1349566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 1359566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 1369566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(snes->xl, &l)); 1379566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(snes->xu, &u)); 1389566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Da, &da)); 1399566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Db, &db)); 140*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 141c2fc9fa9SBarry Smith } 142c2fc9fa9SBarry Smith 143c2fc9fa9SBarry Smith /* 144c2fc9fa9SBarry 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. 145c2fc9fa9SBarry Smith 146c2fc9fa9SBarry Smith Input Parameters: 147c2fc9fa9SBarry Smith . Da - Diagonal shift vector for the semismooth jacobian. 148c2fc9fa9SBarry Smith . Db - Row scaling vector for the semismooth jacobian. 149c2fc9fa9SBarry Smith 150c2fc9fa9SBarry Smith Output Parameters: 151c2fc9fa9SBarry Smith . jac - semismooth jacobian 152c2fc9fa9SBarry Smith . jac_pre - optional preconditioning matrix 153c2fc9fa9SBarry Smith 154f6dfbefdSBarry Smith Note: 155c2fc9fa9SBarry Smith The semismooth jacobian matrix is given by 156c2fc9fa9SBarry Smith jac = Da + Db*jacfun 157c2fc9fa9SBarry Smith where Db is the row scaling matrix stored as a vector, 158c2fc9fa9SBarry Smith Da is the diagonal perturbation matrix stored as a vector 159c2fc9fa9SBarry Smith and jacfun is the jacobian of the original nonlinear function. 160c2fc9fa9SBarry Smith */ 161d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIComputeJacobian(Mat jac, Mat jac_pre, Vec Da, Vec Db) 162d71ae5a4SJacob Faibussowitsch { 163c2fc9fa9SBarry Smith /* Do row scaling and add diagonal perturbation */ 164362febeeSStefano Zampini PetscFunctionBegin; 1659566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(jac, Db, NULL)); 1669566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(jac, Da, ADD_VALUES)); 167c2fc9fa9SBarry Smith if (jac != jac_pre) { /* If jac and jac_pre are different */ 1689566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(jac_pre, Db, NULL)); 1699566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(jac_pre, Da, ADD_VALUES)); 170c2fc9fa9SBarry Smith } 171*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 172c2fc9fa9SBarry Smith } 173c2fc9fa9SBarry Smith 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 */ 188d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi) 189d71ae5a4SJacob Faibussowitsch { 190c2fc9fa9SBarry Smith PetscFunctionBegin; 1919566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(H, phi, dpsi)); 192*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 193c2fc9fa9SBarry Smith } 194c2fc9fa9SBarry Smith 195c2fc9fa9SBarry Smith /* 196f450aa47SBarry Smith SNESSolve_VINEWTONSSLS - Solves the complementarity problem with a semismooth Newton 197c2fc9fa9SBarry Smith method using a line search. 198c2fc9fa9SBarry Smith 199f6dfbefdSBarry Smith Input Parameter: 200c2fc9fa9SBarry Smith . snes - the SNES context 201c2fc9fa9SBarry Smith 202c2fc9fa9SBarry Smith Application Interface Routine: SNESSolve() 203c2fc9fa9SBarry Smith 204f6dfbefdSBarry Smith Note: 205c2fc9fa9SBarry Smith This implements essentially a semismooth Newton method with a 206d5f1b7e6SEd Bueler line search. The default line search does not do any line search 207d5f1b7e6SEd Bueler but rather takes a full Newton step. 208016d19f9SBarry Smith 20904d7464bSBarry 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 210016d19f9SBarry Smith 211c2fc9fa9SBarry Smith */ 212d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSolve_VINEWTONSSLS(SNES snes) 213d71ae5a4SJacob Faibussowitsch { 214f450aa47SBarry Smith SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data; 215c2fc9fa9SBarry Smith PetscInt maxits, i, lits; 216422a814eSBarry Smith SNESLineSearchReason lssucceed; 217c2fc9fa9SBarry Smith PetscReal gnorm, xnorm = 0, ynorm; 2189bd66eb0SPeter Brune Vec Y, X, F; 219c2fc9fa9SBarry Smith KSPConvergedReason kspreason; 2206cab3a1bSJed Brown DM dm; 221942e3340SBarry Smith DMSNES sdm; 222c2fc9fa9SBarry Smith 223c2fc9fa9SBarry Smith PetscFunctionBegin; 2249566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 2259566063dSJacob Faibussowitsch PetscCall(DMGetDMSNES(dm, &sdm)); 2261aa26658SKarl Rupp 22722c6f798SBarry Smith vi->computeuserfunction = sdm->ops->computefunction; 22822c6f798SBarry Smith sdm->ops->computefunction = SNESVIComputeFunction; 229c2fc9fa9SBarry Smith 230c2fc9fa9SBarry Smith snes->numFailures = 0; 231c2fc9fa9SBarry Smith snes->numLinearSolveFailures = 0; 232c2fc9fa9SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 233c2fc9fa9SBarry Smith 234c2fc9fa9SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 235c2fc9fa9SBarry Smith X = snes->vec_sol; /* solution vector */ 236c2fc9fa9SBarry Smith F = snes->vec_func; /* residual vector */ 237c2fc9fa9SBarry Smith Y = snes->work[0]; /* work vectors */ 238c2fc9fa9SBarry Smith 2399566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 240c2fc9fa9SBarry Smith snes->iter = 0; 241c2fc9fa9SBarry Smith snes->norm = 0.0; 2429566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 243c2fc9fa9SBarry Smith 2449566063dSJacob Faibussowitsch PetscCall(SNESVIProjectOntoBounds(snes, X)); 2459566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, vi->phi)); 246c2fc9fa9SBarry Smith if (snes->domainerror) { 247c2fc9fa9SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 24822c6f798SBarry Smith sdm->ops->computefunction = vi->computeuserfunction; 249*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 250c2fc9fa9SBarry Smith } 251c2fc9fa9SBarry Smith /* Compute Merit function */ 2529566063dSJacob Faibussowitsch PetscCall(SNESVIComputeMeritFunction(vi->phi, &vi->merit, &vi->phinorm)); 253c2fc9fa9SBarry Smith 2549566063dSJacob Faibussowitsch PetscCall(VecNormBegin(X, NORM_2, &xnorm)); /* xnorm <- ||x|| */ 2559566063dSJacob Faibussowitsch PetscCall(VecNormEnd(X, NORM_2, &xnorm)); 256422a814eSBarry Smith SNESCheckFunctionNorm(snes, vi->merit); 257c2fc9fa9SBarry Smith 2589566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 259c2fc9fa9SBarry Smith snes->norm = vi->phinorm; 2609566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2619566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, vi->phinorm, 0)); 2629566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, vi->phinorm)); 263c2fc9fa9SBarry Smith 264c2fc9fa9SBarry Smith /* test convergence */ 265dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, vi->phinorm, &snes->reason, snes->cnvP); 266c2fc9fa9SBarry Smith if (snes->reason) { 26722c6f798SBarry Smith sdm->ops->computefunction = vi->computeuserfunction; 268*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 269c2fc9fa9SBarry Smith } 270c2fc9fa9SBarry Smith 271c2fc9fa9SBarry Smith for (i = 0; i < maxits; i++) { 272c2fc9fa9SBarry Smith /* Call general purpose update function */ 273dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 274c2fc9fa9SBarry Smith 275c2fc9fa9SBarry Smith /* Solve J Y = Phi, where J is the semismooth jacobian */ 276b9600128SPeter Brune 277b9600128SPeter Brune /* Get the jacobian -- note that the function must be the original function for snes_fd and snes_fd_color to work for this*/ 27822c6f798SBarry Smith sdm->ops->computefunction = vi->computeuserfunction; 2799566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 28007b62357SFande Kong SNESCheckJacobianDomainerror(snes); 28122c6f798SBarry Smith sdm->ops->computefunction = SNESVIComputeFunction; 282b9600128SPeter Brune 283c2fc9fa9SBarry Smith /* Get the diagonal shift and row scaling vectors */ 2849566063dSJacob Faibussowitsch PetscCall(SNESVIComputeBsubdifferentialVectors(snes, X, F, snes->jacobian, vi->Da, vi->Db)); 285c2fc9fa9SBarry Smith /* Compute the semismooth jacobian */ 2869566063dSJacob Faibussowitsch PetscCall(SNESVIComputeJacobian(snes->jacobian, snes->jacobian_pre, vi->Da, vi->Db)); 287c2fc9fa9SBarry Smith /* Compute the merit function gradient */ 2889566063dSJacob Faibussowitsch PetscCall(SNESVIComputeMeritFunctionGradient(snes->jacobian, vi->phi, vi->dpsi)); 2899566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 2909566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, vi->phi, Y)); 2919566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(snes->ksp, &kspreason)); 292c2fc9fa9SBarry Smith 293c2fc9fa9SBarry Smith if (kspreason < 0) { 294c2fc9fa9SBarry Smith if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 29563a3b9bcSJacob 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)); 296c2fc9fa9SBarry Smith snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 297c2fc9fa9SBarry Smith break; 298c2fc9fa9SBarry Smith } 299c2fc9fa9SBarry Smith } 3009566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 301c2fc9fa9SBarry Smith snes->linear_its += lits; 30263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits)); 303c2fc9fa9SBarry Smith /* 3046b2b7091SBarry Smith if (snes->ops->precheck) { 305c2fc9fa9SBarry Smith PetscBool changed_y = PETSC_FALSE; 306dbbe0bcdSBarry Smith PetscUseTypeMethod(snes,precheck ,X,Y,snes->precheck,&changed_y); 307c2fc9fa9SBarry Smith } 308c2fc9fa9SBarry Smith 3091baa6e33SBarry Smith if (PetscLogPrintInfo) PetscCall(SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W)); 310c2fc9fa9SBarry Smith */ 311c2fc9fa9SBarry Smith /* Compute a (scaled) negative update in the line search routine: 312c2fc9fa9SBarry Smith Y <- X - lambda*Y 313c2fc9fa9SBarry Smith and evaluate G = function(Y) (depends on the line search). 314c2fc9fa9SBarry Smith */ 3159566063dSJacob Faibussowitsch PetscCall(VecCopy(Y, snes->vec_sol_update)); 3169371c9d4SSatish Balay ynorm = 1; 3179371c9d4SSatish Balay gnorm = vi->phinorm; 3189566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(snes->linesearch, X, vi->phi, &gnorm, Y)); 3199566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed)); 3209566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &gnorm, &ynorm)); 3219566063dSJacob 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)); 322c2fc9fa9SBarry Smith if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 323c2fc9fa9SBarry Smith if (snes->domainerror) { 324c2fc9fa9SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 32522c6f798SBarry Smith sdm->ops->computefunction = vi->computeuserfunction; 326*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 327c2fc9fa9SBarry Smith } 328422a814eSBarry Smith if (lssucceed) { 329c2fc9fa9SBarry Smith if (++snes->numFailures >= snes->maxFailures) { 330c2fc9fa9SBarry Smith PetscBool ismin; 331c2fc9fa9SBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 3329566063dSJacob Faibussowitsch PetscCall(SNESVICheckLocalMin_Private(snes, snes->jacobian, vi->phi, X, gnorm, &ismin)); 333c2fc9fa9SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 334c2fc9fa9SBarry Smith break; 335c2fc9fa9SBarry Smith } 336c2fc9fa9SBarry Smith } 337c2fc9fa9SBarry Smith /* Update function and solution vectors */ 338c2fc9fa9SBarry Smith vi->phinorm = gnorm; 339c2fc9fa9SBarry Smith vi->merit = 0.5 * vi->phinorm * vi->phinorm; 340c2fc9fa9SBarry Smith /* Monitor convergence */ 3419566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 342c2fc9fa9SBarry Smith snes->iter = i + 1; 343c2fc9fa9SBarry Smith snes->norm = vi->phinorm; 344c1e67a49SFande Kong snes->xnorm = xnorm; 345c1e67a49SFande Kong snes->ynorm = ynorm; 3469566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 3479566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 3489566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 349c2fc9fa9SBarry Smith /* Test for convergence, xnorm = || X || */ 3509566063dSJacob Faibussowitsch if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm)); 351dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, vi->phinorm, &snes->reason, snes->cnvP); 352c2fc9fa9SBarry Smith if (snes->reason) break; 353c2fc9fa9SBarry Smith } 354c2fc9fa9SBarry Smith if (i == maxits) { 35563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits)); 356c2fc9fa9SBarry Smith if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 357c2fc9fa9SBarry Smith } 35822c6f798SBarry Smith sdm->ops->computefunction = vi->computeuserfunction; 359*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 360c2fc9fa9SBarry Smith } 361c2fc9fa9SBarry Smith 362c2fc9fa9SBarry Smith /* 363f450aa47SBarry Smith SNESSetUp_VINEWTONSSLS - Sets up the internal data structures for the later use 364c2fc9fa9SBarry Smith of the SNES nonlinear solver. 365c2fc9fa9SBarry Smith 366c2fc9fa9SBarry Smith Input Parameter: 367c2fc9fa9SBarry Smith . snes - the SNES context 368c2fc9fa9SBarry Smith 369c2fc9fa9SBarry Smith Application Interface Routine: SNESSetUp() 370c2fc9fa9SBarry Smith 371f6dfbefdSBarry Smith Note: 372c2fc9fa9SBarry Smith For basic use of the SNES solvers, the user need not explicitly call 373c2fc9fa9SBarry Smith SNESSetUp(), since these actions will automatically occur during 374c2fc9fa9SBarry Smith the call to SNESSolve(). 375c2fc9fa9SBarry Smith */ 376d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetUp_VINEWTONSSLS(SNES snes) 377d71ae5a4SJacob Faibussowitsch { 378f450aa47SBarry Smith SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data; 379c2fc9fa9SBarry Smith 380c2fc9fa9SBarry Smith PetscFunctionBegin; 3819566063dSJacob Faibussowitsch PetscCall(SNESSetUp_VI(snes)); 3829566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol, &vi->dpsi)); 3839566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol, &vi->phi)); 3849566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol, &vi->Da)); 3859566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol, &vi->Db)); 3869566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol, &vi->z)); 3879566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol, &vi->t)); 388*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 389c2fc9fa9SBarry Smith } 390f6dfbefdSBarry Smith 391d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESReset_VINEWTONSSLS(SNES snes) 392d71ae5a4SJacob Faibussowitsch { 393f450aa47SBarry Smith SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS *)snes->data; 394c2fc9fa9SBarry Smith 395c2fc9fa9SBarry Smith PetscFunctionBegin; 3969566063dSJacob Faibussowitsch PetscCall(SNESReset_VI(snes)); 3979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vi->dpsi)); 3989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vi->phi)); 3999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vi->Da)); 4009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vi->Db)); 4019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vi->z)); 4029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vi->t)); 403*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 404c2fc9fa9SBarry Smith } 405c2fc9fa9SBarry Smith 406c2fc9fa9SBarry Smith /* 407f450aa47SBarry Smith SNESSetFromOptions_VINEWTONSSLS - Sets various parameters for the SNESVI method. 408c2fc9fa9SBarry Smith 409c2fc9fa9SBarry Smith Input Parameter: 410c2fc9fa9SBarry Smith . snes - the SNES context 411c2fc9fa9SBarry Smith 412c2fc9fa9SBarry Smith Application Interface Routine: SNESSetFromOptions() 413c2fc9fa9SBarry Smith */ 414d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_VINEWTONSSLS(SNES snes, PetscOptionItems *PetscOptionsObject) 415d71ae5a4SJacob Faibussowitsch { 416c2fc9fa9SBarry Smith PetscFunctionBegin; 417dbbe0bcdSBarry Smith PetscCall(SNESSetFromOptions_VI(snes, PetscOptionsObject)); 418d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES semismooth method options"); 419d0609cedSBarry Smith PetscOptionsHeadEnd(); 420*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 421c2fc9fa9SBarry Smith } 422c2fc9fa9SBarry Smith 423c2fc9fa9SBarry Smith /*MC 424f450aa47SBarry Smith SNESVINEWTONSSLS - Semi-smooth solver for variational inequalities based on Newton's method 425c2fc9fa9SBarry Smith 426f6dfbefdSBarry Smith Options Database Keys: 427d5f1b7e6SEd Bueler + -snes_type <vinewtonssls,vinewtonrsls> a semi-smooth solver, a reduced space active set method 42861589011SJed Brown - -snes_vi_monitor - prints the number of active constraints at each iteration. 42961589011SJed Brown 430c2fc9fa9SBarry Smith Level: beginner 431c2fc9fa9SBarry Smith 432b80f3ac1SShri Abhyankar References: 433606c0280SSatish Balay + * - T. S. Munson, F. Facchinei, M. C. Ferris, A. Fischer, and C. Kanzow. The semismooth 434b80f3ac1SShri Abhyankar algorithm for large scale complementarity problems. INFORMS Journal on Computing, 13 (2001). 435606c0280SSatish Balay - * - T. S. Munson, and S. Benson. Flexible Complementarity Solvers for Large Scale 436d5f1b7e6SEd Bueler Applications, Optimization Methods and Software, 21 (2006). 437b80f3ac1SShri Abhyankar 438f6dfbefdSBarry Smith Notes: 439f6dfbefdSBarry Smith This family of algorithm is much like an interior point method. 440c2fc9fa9SBarry Smith 441f6dfbefdSBarry Smith The reduced space active set solvers `SNESVINEWTONRSLS` provide an alternative approach that does not result in extremely ill-conditioned linear systems 442f6dfbefdSBarry Smith 443f6dfbefdSBarry Smith .seealso: `SNESVINEWTONRSLS`, `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESVINEWTONRSLS`, `SNESNEWTONTR`, `SNESLineSearchSetType()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()` 444c2fc9fa9SBarry Smith M*/ 445d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_VINEWTONSSLS(SNES snes) 446d71ae5a4SJacob Faibussowitsch { 447f450aa47SBarry Smith SNES_VINEWTONSSLS *vi; 448d8d34be6SBarry Smith SNESLineSearch linesearch; 449c2fc9fa9SBarry Smith 450c2fc9fa9SBarry Smith PetscFunctionBegin; 451f450aa47SBarry Smith snes->ops->reset = SNESReset_VINEWTONSSLS; 452f450aa47SBarry Smith snes->ops->setup = SNESSetUp_VINEWTONSSLS; 453f450aa47SBarry Smith snes->ops->solve = SNESSolve_VINEWTONSSLS; 454c2fc9fa9SBarry Smith snes->ops->destroy = SNESDestroy_VI; 455f450aa47SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_VINEWTONSSLS; 4560298fd71SBarry Smith snes->ops->view = NULL; 457c2fc9fa9SBarry Smith 458c2fc9fa9SBarry Smith snes->usesksp = PETSC_TRUE; 459efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 460c2fc9fa9SBarry Smith 4619566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 462ec786807SJed Brown if (!((PetscObject)linesearch)->type_name) { 4639566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT)); 4649566063dSJacob Faibussowitsch PetscCall(SNESLineSearchBTSetAlpha(linesearch, 0.0)); 465ec786807SJed Brown } 466d8d34be6SBarry Smith 4674fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 4684fc747eaSLawrence Mitchell 4694dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&vi)); 470c2fc9fa9SBarry Smith snes->data = (void *)vi; 471c2fc9fa9SBarry Smith 4729566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetVariableBounds_C", SNESVISetVariableBounds_VI)); 4739566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESVISetComputeVariableBounds_C", SNESVISetComputeVariableBounds_VI)); 474*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 475c2fc9fa9SBarry Smith } 476