1c2fc9fa9SBarry Smith 2c2fc9fa9SBarry Smith #include <../src/snes/impls/vi/ss/vissimpl.h> /*I "petscsnes.h" I*/ 3af0996ceSBarry Smith #include <../include/petsc/private/kspimpl.h> 4af0996ceSBarry Smith #include <../include/petsc/private/matimpl.h> 5af0996ceSBarry Smith #include <../include/petsc/private/dmimpl.h> 6c2fc9fa9SBarry Smith 7c2fc9fa9SBarry Smith 8c2fc9fa9SBarry Smith /* 9c2fc9fa9SBarry Smith SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem. 10c2fc9fa9SBarry Smith 11c2fc9fa9SBarry Smith Input Parameter: 12c2fc9fa9SBarry Smith . phi - the semismooth function 13c2fc9fa9SBarry Smith 14c2fc9fa9SBarry Smith Output Parameter: 15c2fc9fa9SBarry Smith . merit - the merit function 16c2fc9fa9SBarry Smith . phinorm - ||phi|| 17c2fc9fa9SBarry Smith 18c2fc9fa9SBarry Smith Notes: 19c2fc9fa9SBarry Smith The merit function for the mixed complementarity problem is defined as 20c2fc9fa9SBarry Smith merit = 0.5*phi^T*phi 21c2fc9fa9SBarry Smith */ 22c2fc9fa9SBarry Smith #undef __FUNCT__ 23c2fc9fa9SBarry Smith #define __FUNCT__ "SNESVIComputeMeritFunction" 24c2fc9fa9SBarry Smith static PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal *merit,PetscReal *phinorm) 25c2fc9fa9SBarry Smith { 26c2fc9fa9SBarry Smith PetscErrorCode ierr; 27c2fc9fa9SBarry Smith 28c2fc9fa9SBarry Smith PetscFunctionBegin; 29c2fc9fa9SBarry Smith ierr = VecNormBegin(phi,NORM_2,phinorm);CHKERRQ(ierr); 30c2fc9fa9SBarry Smith ierr = VecNormEnd(phi,NORM_2,phinorm);CHKERRQ(ierr); 31c2fc9fa9SBarry Smith 32c2fc9fa9SBarry Smith *merit = 0.5*(*phinorm)*(*phinorm); 33c2fc9fa9SBarry Smith PetscFunctionReturn(0); 34c2fc9fa9SBarry Smith } 35c2fc9fa9SBarry Smith 36c2fc9fa9SBarry Smith PETSC_STATIC_INLINE PetscScalar Phi(PetscScalar a,PetscScalar b) 37c2fc9fa9SBarry Smith { 38c2fc9fa9SBarry Smith return a + b - PetscSqrtScalar(a*a + b*b); 39c2fc9fa9SBarry Smith } 40c2fc9fa9SBarry Smith 41c2fc9fa9SBarry Smith PETSC_STATIC_INLINE PetscScalar DPhi(PetscScalar a,PetscScalar b) 42c2fc9fa9SBarry Smith { 43c2fc9fa9SBarry Smith if ((PetscAbsScalar(a) >= 1.e-6) || (PetscAbsScalar(b) >= 1.e-6)) return 1.0 - a/ PetscSqrtScalar(a*a + b*b); 44c2fc9fa9SBarry Smith else return .5; 45c2fc9fa9SBarry Smith } 46c2fc9fa9SBarry Smith 47c2fc9fa9SBarry Smith /* 48c2fc9fa9SBarry Smith SNESVIComputeFunction - Reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear equations in semismooth form. 49c2fc9fa9SBarry Smith 50c2fc9fa9SBarry Smith Input Parameters: 51c2fc9fa9SBarry Smith . snes - the SNES context 52d5f1b7e6SEd Bueler . X - current iterate 53c2fc9fa9SBarry Smith . functx - user defined function context 54c2fc9fa9SBarry Smith 55c2fc9fa9SBarry Smith Output Parameters: 56c2fc9fa9SBarry Smith . phi - Semismooth function 57c2fc9fa9SBarry Smith 58c2fc9fa9SBarry Smith */ 59c2fc9fa9SBarry Smith #undef __FUNCT__ 60c2fc9fa9SBarry Smith #define __FUNCT__ "SNESVIComputeFunction" 61c2fc9fa9SBarry Smith static PetscErrorCode SNESVIComputeFunction(SNES snes,Vec X,Vec phi,void *functx) 62c2fc9fa9SBarry Smith { 63c2fc9fa9SBarry Smith PetscErrorCode ierr; 64f450aa47SBarry Smith SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS*)snes->data; 65c2fc9fa9SBarry Smith Vec Xl = snes->xl,Xu = snes->xu,F = snes->vec_func; 665edff71fSBarry Smith PetscScalar *phi_arr,*f_arr,*l,*u; 675edff71fSBarry Smith const PetscScalar *x_arr; 68c2fc9fa9SBarry Smith PetscInt i,nlocal; 69c2fc9fa9SBarry Smith 70c2fc9fa9SBarry Smith PetscFunctionBegin; 71c2fc9fa9SBarry Smith ierr = (*vi->computeuserfunction)(snes,X,F,functx);CHKERRQ(ierr); 72c2fc9fa9SBarry Smith ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 735edff71fSBarry Smith ierr = VecGetArrayRead(X,&x_arr);CHKERRQ(ierr); 74c2fc9fa9SBarry Smith ierr = VecGetArray(F,&f_arr);CHKERRQ(ierr); 75c2fc9fa9SBarry Smith ierr = VecGetArray(Xl,&l);CHKERRQ(ierr); 76c2fc9fa9SBarry Smith ierr = VecGetArray(Xu,&u);CHKERRQ(ierr); 77c2fc9fa9SBarry Smith ierr = VecGetArray(phi,&phi_arr);CHKERRQ(ierr); 78c2fc9fa9SBarry Smith 79c2fc9fa9SBarry Smith for (i=0; i < nlocal; i++) { 80e270355aSBarry Smith if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { /* no constraints on variable */ 81c2fc9fa9SBarry Smith phi_arr[i] = f_arr[i]; 82e270355aSBarry Smith } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { /* upper bound on variable only */ 83c2fc9fa9SBarry Smith phi_arr[i] = -Phi(u[i] - x_arr[i],-f_arr[i]); 84e270355aSBarry Smith } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { /* lower bound on variable only */ 85c2fc9fa9SBarry Smith phi_arr[i] = Phi(x_arr[i] - l[i],f_arr[i]); 86c2fc9fa9SBarry Smith } else if (l[i] == u[i]) { 87c2fc9fa9SBarry Smith phi_arr[i] = l[i] - x_arr[i]; 88c2fc9fa9SBarry Smith } else { /* both bounds on variable */ 89c2fc9fa9SBarry Smith phi_arr[i] = Phi(x_arr[i] - l[i],-Phi(u[i] - x_arr[i],-f_arr[i])); 90c2fc9fa9SBarry Smith } 91c2fc9fa9SBarry Smith } 92c2fc9fa9SBarry Smith 935edff71fSBarry Smith ierr = VecRestoreArrayRead(X,&x_arr);CHKERRQ(ierr); 94c2fc9fa9SBarry Smith ierr = VecRestoreArray(F,&f_arr);CHKERRQ(ierr); 95c2fc9fa9SBarry Smith ierr = VecRestoreArray(Xl,&l);CHKERRQ(ierr); 96c2fc9fa9SBarry Smith ierr = VecRestoreArray(Xu,&u);CHKERRQ(ierr); 97c2fc9fa9SBarry Smith ierr = VecRestoreArray(phi,&phi_arr);CHKERRQ(ierr); 98c2fc9fa9SBarry Smith PetscFunctionReturn(0); 99c2fc9fa9SBarry Smith } 100c2fc9fa9SBarry Smith 101c2fc9fa9SBarry Smith /* 102c2fc9fa9SBarry Smith SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the 103c2fc9fa9SBarry Smith the semismooth jacobian. 104c2fc9fa9SBarry Smith */ 105c2fc9fa9SBarry Smith #undef __FUNCT__ 106c2fc9fa9SBarry Smith #define __FUNCT__ "SNESVIComputeBsubdifferentialVectors" 107c2fc9fa9SBarry Smith PetscErrorCode SNESVIComputeBsubdifferentialVectors(SNES snes,Vec X,Vec F,Mat jac,Vec Da,Vec Db) 108c2fc9fa9SBarry Smith { 109c2fc9fa9SBarry Smith PetscErrorCode ierr; 110c2fc9fa9SBarry Smith PetscScalar *l,*u,*x,*f,*da,*db,da1,da2,db1,db2; 111c2fc9fa9SBarry Smith PetscInt i,nlocal; 112c2fc9fa9SBarry Smith 113c2fc9fa9SBarry Smith PetscFunctionBegin; 114c2fc9fa9SBarry Smith ierr = VecGetArray(X,&x);CHKERRQ(ierr); 115c2fc9fa9SBarry Smith ierr = VecGetArray(F,&f);CHKERRQ(ierr); 116c2fc9fa9SBarry Smith ierr = VecGetArray(snes->xl,&l);CHKERRQ(ierr); 117c2fc9fa9SBarry Smith ierr = VecGetArray(snes->xu,&u);CHKERRQ(ierr); 118c2fc9fa9SBarry Smith ierr = VecGetArray(Da,&da);CHKERRQ(ierr); 119c2fc9fa9SBarry Smith ierr = VecGetArray(Db,&db);CHKERRQ(ierr); 120c2fc9fa9SBarry Smith ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 121c2fc9fa9SBarry Smith 122c2fc9fa9SBarry Smith for (i=0; i< nlocal; i++) { 123e270355aSBarry Smith if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { /* no constraints on variable */ 124c2fc9fa9SBarry Smith da[i] = 0; 125c2fc9fa9SBarry Smith db[i] = 1; 126e270355aSBarry Smith } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) { /* upper bound on variable only */ 127c2fc9fa9SBarry Smith da[i] = DPhi(u[i] - x[i], -f[i]); 128c2fc9fa9SBarry Smith db[i] = DPhi(-f[i],u[i] - x[i]); 129e270355aSBarry Smith } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) { /* lower bound on variable only */ 130c2fc9fa9SBarry Smith da[i] = DPhi(x[i] - l[i], f[i]); 131c2fc9fa9SBarry Smith db[i] = DPhi(f[i],x[i] - l[i]); 132c2fc9fa9SBarry Smith } else if (l[i] == u[i]) { /* fixed variable */ 133c2fc9fa9SBarry Smith da[i] = 1; 134c2fc9fa9SBarry Smith db[i] = 0; 135c2fc9fa9SBarry Smith } else { /* upper and lower bounds on variable */ 136c2fc9fa9SBarry Smith da1 = DPhi(x[i] - l[i], -Phi(u[i] - x[i], -f[i])); 137c2fc9fa9SBarry Smith db1 = DPhi(-Phi(u[i] - x[i], -f[i]),x[i] - l[i]); 138c2fc9fa9SBarry Smith da2 = DPhi(u[i] - x[i], -f[i]); 139c2fc9fa9SBarry Smith db2 = DPhi(-f[i],u[i] - x[i]); 140c2fc9fa9SBarry Smith da[i] = da1 + db1*da2; 141c2fc9fa9SBarry Smith db[i] = db1*db2; 142c2fc9fa9SBarry Smith } 143c2fc9fa9SBarry Smith } 144c2fc9fa9SBarry Smith 145c2fc9fa9SBarry Smith ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 146c2fc9fa9SBarry Smith ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 147c2fc9fa9SBarry Smith ierr = VecRestoreArray(snes->xl,&l);CHKERRQ(ierr); 148c2fc9fa9SBarry Smith ierr = VecRestoreArray(snes->xu,&u);CHKERRQ(ierr); 149c2fc9fa9SBarry Smith ierr = VecRestoreArray(Da,&da);CHKERRQ(ierr); 150c2fc9fa9SBarry Smith ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr); 151c2fc9fa9SBarry Smith PetscFunctionReturn(0); 152c2fc9fa9SBarry Smith } 153c2fc9fa9SBarry Smith 154c2fc9fa9SBarry Smith /* 155c2fc9fa9SBarry 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. 156c2fc9fa9SBarry Smith 157c2fc9fa9SBarry Smith Input Parameters: 158c2fc9fa9SBarry Smith . Da - Diagonal shift vector for the semismooth jacobian. 159c2fc9fa9SBarry Smith . Db - Row scaling vector for the semismooth jacobian. 160c2fc9fa9SBarry Smith 161c2fc9fa9SBarry Smith Output Parameters: 162c2fc9fa9SBarry Smith . jac - semismooth jacobian 163c2fc9fa9SBarry Smith . jac_pre - optional preconditioning matrix 164c2fc9fa9SBarry Smith 165c2fc9fa9SBarry Smith Notes: 166c2fc9fa9SBarry Smith The semismooth jacobian matrix is given by 167c2fc9fa9SBarry Smith jac = Da + Db*jacfun 168c2fc9fa9SBarry Smith where Db is the row scaling matrix stored as a vector, 169c2fc9fa9SBarry Smith Da is the diagonal perturbation matrix stored as a vector 170c2fc9fa9SBarry Smith and jacfun is the jacobian of the original nonlinear function. 171c2fc9fa9SBarry Smith */ 172c2fc9fa9SBarry Smith #undef __FUNCT__ 173c2fc9fa9SBarry Smith #define __FUNCT__ "SNESVIComputeJacobian" 174c2fc9fa9SBarry Smith PetscErrorCode SNESVIComputeJacobian(Mat jac, Mat jac_pre,Vec Da, Vec Db) 175c2fc9fa9SBarry Smith { 176c2fc9fa9SBarry Smith PetscErrorCode ierr; 177c2fc9fa9SBarry Smith 178c2fc9fa9SBarry Smith /* Do row scaling and add diagonal perturbation */ 1790298fd71SBarry Smith ierr = MatDiagonalScale(jac,Db,NULL);CHKERRQ(ierr); 180c2fc9fa9SBarry Smith ierr = MatDiagonalSet(jac,Da,ADD_VALUES);CHKERRQ(ierr); 181c2fc9fa9SBarry Smith if (jac != jac_pre) { /* If jac and jac_pre are different */ 1820298fd71SBarry Smith ierr = MatDiagonalScale(jac_pre,Db,NULL);CHKERRQ(ierr); 183c2fc9fa9SBarry Smith ierr = MatDiagonalSet(jac_pre,Da,ADD_VALUES);CHKERRQ(ierr); 184c2fc9fa9SBarry Smith } 185c2fc9fa9SBarry Smith PetscFunctionReturn(0); 186c2fc9fa9SBarry Smith } 187c2fc9fa9SBarry Smith 188c2fc9fa9SBarry Smith /* 189c2fc9fa9SBarry Smith SNESVIComputeMeritFunctionGradient - Computes the gradient of the merit function psi. 190c2fc9fa9SBarry Smith 191c2fc9fa9SBarry Smith Input Parameters: 192c2fc9fa9SBarry Smith phi - semismooth function. 193c2fc9fa9SBarry Smith H - semismooth jacobian 194c2fc9fa9SBarry Smith 195c2fc9fa9SBarry Smith Output Parameters: 196c2fc9fa9SBarry Smith dpsi - merit function gradient 197c2fc9fa9SBarry Smith 198c2fc9fa9SBarry Smith Notes: 199c2fc9fa9SBarry Smith The merit function gradient is computed as follows 200c2fc9fa9SBarry Smith dpsi = H^T*phi 201c2fc9fa9SBarry Smith */ 202c2fc9fa9SBarry Smith #undef __FUNCT__ 203c2fc9fa9SBarry Smith #define __FUNCT__ "SNESVIComputeMeritFunctionGradient" 204c2fc9fa9SBarry Smith PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi) 205c2fc9fa9SBarry Smith { 206c2fc9fa9SBarry Smith PetscErrorCode ierr; 207c2fc9fa9SBarry Smith 208c2fc9fa9SBarry Smith PetscFunctionBegin; 209c2fc9fa9SBarry Smith ierr = MatMultTranspose(H,phi,dpsi);CHKERRQ(ierr); 210c2fc9fa9SBarry Smith PetscFunctionReturn(0); 211c2fc9fa9SBarry Smith } 212c2fc9fa9SBarry Smith 213c2fc9fa9SBarry Smith 214c2fc9fa9SBarry Smith 215c2fc9fa9SBarry Smith /* 216f450aa47SBarry Smith SNESSolve_VINEWTONSSLS - Solves the complementarity problem with a semismooth Newton 217c2fc9fa9SBarry Smith method using a line search. 218c2fc9fa9SBarry Smith 219c2fc9fa9SBarry Smith Input Parameters: 220c2fc9fa9SBarry Smith . snes - the SNES context 221c2fc9fa9SBarry Smith 222c2fc9fa9SBarry Smith Application Interface Routine: SNESSolve() 223c2fc9fa9SBarry Smith 224c2fc9fa9SBarry Smith Notes: 225c2fc9fa9SBarry Smith This implements essentially a semismooth Newton method with a 226d5f1b7e6SEd Bueler line search. The default line search does not do any line search 227d5f1b7e6SEd Bueler but rather takes a full Newton step. 228016d19f9SBarry Smith 22904d7464bSBarry 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 230016d19f9SBarry Smith 231c2fc9fa9SBarry Smith */ 232c2fc9fa9SBarry Smith #undef __FUNCT__ 233f450aa47SBarry Smith #define __FUNCT__ "SNESSolve_VINEWTONSSLS" 234f450aa47SBarry Smith PetscErrorCode SNESSolve_VINEWTONSSLS(SNES snes) 235c2fc9fa9SBarry Smith { 236f450aa47SBarry Smith SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS*)snes->data; 237c2fc9fa9SBarry Smith PetscErrorCode ierr; 238c2fc9fa9SBarry Smith PetscInt maxits,i,lits; 239422a814eSBarry Smith SNESLineSearchReason lssucceed; 240c2fc9fa9SBarry Smith PetscReal gnorm,xnorm=0,ynorm; 2419bd66eb0SPeter Brune Vec Y,X,F; 242c2fc9fa9SBarry Smith KSPConvergedReason kspreason; 2436cab3a1bSJed Brown DM dm; 244942e3340SBarry Smith DMSNES sdm; 245c2fc9fa9SBarry Smith 246c2fc9fa9SBarry Smith PetscFunctionBegin; 2476cab3a1bSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 248942e3340SBarry Smith ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr); 2491aa26658SKarl Rupp 25022c6f798SBarry Smith vi->computeuserfunction = sdm->ops->computefunction; 25122c6f798SBarry Smith sdm->ops->computefunction = SNESVIComputeFunction; 252c2fc9fa9SBarry Smith 253c2fc9fa9SBarry Smith snes->numFailures = 0; 254c2fc9fa9SBarry Smith snes->numLinearSolveFailures = 0; 255c2fc9fa9SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 256c2fc9fa9SBarry Smith 257c2fc9fa9SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 258c2fc9fa9SBarry Smith X = snes->vec_sol; /* solution vector */ 259c2fc9fa9SBarry Smith F = snes->vec_func; /* residual vector */ 260c2fc9fa9SBarry Smith Y = snes->work[0]; /* work vectors */ 261c2fc9fa9SBarry Smith 262e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 263c2fc9fa9SBarry Smith snes->iter = 0; 264c2fc9fa9SBarry Smith snes->norm = 0.0; 265e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 266c2fc9fa9SBarry Smith 267c2fc9fa9SBarry Smith ierr = SNESVIProjectOntoBounds(snes,X);CHKERRQ(ierr); 268c2fc9fa9SBarry Smith ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr); 269c2fc9fa9SBarry Smith if (snes->domainerror) { 270c2fc9fa9SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 27122c6f798SBarry Smith sdm->ops->computefunction = vi->computeuserfunction; 272c2fc9fa9SBarry Smith PetscFunctionReturn(0); 273c2fc9fa9SBarry Smith } 274c2fc9fa9SBarry Smith /* Compute Merit function */ 275c2fc9fa9SBarry Smith ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr); 276c2fc9fa9SBarry Smith 277c2fc9fa9SBarry Smith ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 278c2fc9fa9SBarry Smith ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 279422a814eSBarry Smith SNESCheckFunctionNorm(snes,vi->merit); 280c2fc9fa9SBarry Smith 281e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 282c2fc9fa9SBarry Smith snes->norm = vi->phinorm; 283e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 284a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,vi->phinorm,0);CHKERRQ(ierr); 285c2fc9fa9SBarry Smith ierr = SNESMonitor(snes,0,vi->phinorm);CHKERRQ(ierr); 286c2fc9fa9SBarry Smith 287c2fc9fa9SBarry Smith /* test convergence */ 288c2fc9fa9SBarry Smith ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 289c2fc9fa9SBarry Smith if (snes->reason) { 29022c6f798SBarry Smith sdm->ops->computefunction = vi->computeuserfunction; 291c2fc9fa9SBarry Smith PetscFunctionReturn(0); 292c2fc9fa9SBarry Smith } 293c2fc9fa9SBarry Smith 294c2fc9fa9SBarry Smith for (i=0; i<maxits; i++) { 295c2fc9fa9SBarry Smith 296c2fc9fa9SBarry Smith /* Call general purpose update function */ 297c2fc9fa9SBarry Smith if (snes->ops->update) { 298c2fc9fa9SBarry Smith ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 299c2fc9fa9SBarry Smith } 300c2fc9fa9SBarry Smith 301c2fc9fa9SBarry Smith /* Solve J Y = Phi, where J is the semismooth jacobian */ 302b9600128SPeter Brune 303b9600128SPeter Brune /* Get the jacobian -- note that the function must be the original function for snes_fd and snes_fd_color to work for this*/ 30422c6f798SBarry Smith sdm->ops->computefunction = vi->computeuserfunction; 305d1e9a80fSBarry Smith ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 30622c6f798SBarry Smith sdm->ops->computefunction = SNESVIComputeFunction; 307b9600128SPeter Brune 308c2fc9fa9SBarry Smith /* Get the diagonal shift and row scaling vectors */ 309c2fc9fa9SBarry Smith ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr); 310c2fc9fa9SBarry Smith /* Compute the semismooth jacobian */ 311c2fc9fa9SBarry Smith ierr = SNESVIComputeJacobian(snes->jacobian,snes->jacobian_pre,vi->Da,vi->Db);CHKERRQ(ierr); 312c2fc9fa9SBarry Smith /* Compute the merit function gradient */ 313c2fc9fa9SBarry Smith ierr = SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);CHKERRQ(ierr); 31423ee1639SBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 315d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,vi->phi,Y);CHKERRQ(ierr); 316c2fc9fa9SBarry Smith ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 317c2fc9fa9SBarry Smith 318c2fc9fa9SBarry Smith if (kspreason < 0) { 319c2fc9fa9SBarry Smith if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 320c2fc9fa9SBarry Smith ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 321c2fc9fa9SBarry Smith snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 322c2fc9fa9SBarry Smith break; 323c2fc9fa9SBarry Smith } 324c2fc9fa9SBarry Smith } 325c2fc9fa9SBarry Smith ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 326c2fc9fa9SBarry Smith snes->linear_its += lits; 327c2fc9fa9SBarry Smith ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 328c2fc9fa9SBarry Smith /* 3296b2b7091SBarry Smith if (snes->ops->precheck) { 330c2fc9fa9SBarry Smith PetscBool changed_y = PETSC_FALSE; 3316b2b7091SBarry Smith ierr = (*snes->ops->precheck)(snes,X,Y,snes->precheck,&changed_y);CHKERRQ(ierr); 332c2fc9fa9SBarry Smith } 333c2fc9fa9SBarry Smith 334c2fc9fa9SBarry Smith if (PetscLogPrintInfo) { 335c2fc9fa9SBarry Smith ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 336c2fc9fa9SBarry Smith } 337c2fc9fa9SBarry Smith */ 338c2fc9fa9SBarry Smith /* Compute a (scaled) negative update in the line search routine: 339c2fc9fa9SBarry Smith Y <- X - lambda*Y 340c2fc9fa9SBarry Smith and evaluate G = function(Y) (depends on the line search). 341c2fc9fa9SBarry Smith */ 342c2fc9fa9SBarry Smith ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 343c2fc9fa9SBarry Smith ynorm = 1; gnorm = vi->phinorm; 344f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, vi->phi, &gnorm, Y);CHKERRQ(ierr); 345422a814eSBarry Smith ierr = SNESLineSearchGetReason(snes->linesearch, &lssucceed);CHKERRQ(ierr); 346f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr); 347c2fc9fa9SBarry Smith ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)vi->phinorm,(double)gnorm,(double)ynorm,(int)lssucceed);CHKERRQ(ierr); 348c2fc9fa9SBarry Smith if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 349c2fc9fa9SBarry Smith if (snes->domainerror) { 350c2fc9fa9SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 35122c6f798SBarry Smith sdm->ops->computefunction = vi->computeuserfunction; 352c2fc9fa9SBarry Smith PetscFunctionReturn(0); 353c2fc9fa9SBarry Smith } 354422a814eSBarry Smith if (lssucceed) { 355c2fc9fa9SBarry Smith if (++snes->numFailures >= snes->maxFailures) { 356c2fc9fa9SBarry Smith PetscBool ismin; 357c2fc9fa9SBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 3589bd66eb0SPeter Brune ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,vi->phi,X,gnorm,&ismin);CHKERRQ(ierr); 359c2fc9fa9SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 360c2fc9fa9SBarry Smith break; 361c2fc9fa9SBarry Smith } 362c2fc9fa9SBarry Smith } 363c2fc9fa9SBarry Smith /* Update function and solution vectors */ 364c2fc9fa9SBarry Smith vi->phinorm = gnorm; 365c2fc9fa9SBarry Smith vi->merit = 0.5*vi->phinorm*vi->phinorm; 366c2fc9fa9SBarry Smith /* Monitor convergence */ 367e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 368c2fc9fa9SBarry Smith snes->iter = i+1; 369c2fc9fa9SBarry Smith snes->norm = vi->phinorm; 370e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 371a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr); 372c2fc9fa9SBarry Smith ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 373c2fc9fa9SBarry Smith /* Test for convergence, xnorm = || X || */ 374e2a6519dSDmitry Karpeev if (snes->ops->converged != SNESConvergedSkip) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 375c2fc9fa9SBarry Smith ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 376c2fc9fa9SBarry Smith if (snes->reason) break; 377c2fc9fa9SBarry Smith } 378c2fc9fa9SBarry Smith if (i == maxits) { 379c2fc9fa9SBarry Smith ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 380c2fc9fa9SBarry Smith if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 381c2fc9fa9SBarry Smith } 38222c6f798SBarry Smith sdm->ops->computefunction = vi->computeuserfunction; 383c2fc9fa9SBarry Smith PetscFunctionReturn(0); 384c2fc9fa9SBarry Smith } 385c2fc9fa9SBarry Smith 386c2fc9fa9SBarry Smith /* -------------------------------------------------------------------------- */ 387c2fc9fa9SBarry Smith /* 388f450aa47SBarry Smith SNESSetUp_VINEWTONSSLS - Sets up the internal data structures for the later use 389c2fc9fa9SBarry Smith of the SNES nonlinear solver. 390c2fc9fa9SBarry Smith 391c2fc9fa9SBarry Smith Input Parameter: 392c2fc9fa9SBarry Smith . snes - the SNES context 393c2fc9fa9SBarry Smith 394c2fc9fa9SBarry Smith Application Interface Routine: SNESSetUp() 395c2fc9fa9SBarry Smith 396c2fc9fa9SBarry Smith Notes: 397c2fc9fa9SBarry Smith For basic use of the SNES solvers, the user need not explicitly call 398c2fc9fa9SBarry Smith SNESSetUp(), since these actions will automatically occur during 399c2fc9fa9SBarry Smith the call to SNESSolve(). 400c2fc9fa9SBarry Smith */ 401c2fc9fa9SBarry Smith #undef __FUNCT__ 402f450aa47SBarry Smith #define __FUNCT__ "SNESSetUp_VINEWTONSSLS" 403f450aa47SBarry Smith PetscErrorCode SNESSetUp_VINEWTONSSLS(SNES snes) 404c2fc9fa9SBarry Smith { 405c2fc9fa9SBarry Smith PetscErrorCode ierr; 406f450aa47SBarry Smith SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS*) snes->data; 407c2fc9fa9SBarry Smith 408c2fc9fa9SBarry Smith PetscFunctionBegin; 409c2fc9fa9SBarry Smith ierr = SNESSetUp_VI(snes);CHKERRQ(ierr); 410c2fc9fa9SBarry Smith ierr = VecDuplicate(snes->vec_sol, &vi->dpsi);CHKERRQ(ierr); 411c2fc9fa9SBarry Smith ierr = VecDuplicate(snes->vec_sol, &vi->phi);CHKERRQ(ierr); 412c2fc9fa9SBarry Smith ierr = VecDuplicate(snes->vec_sol, &vi->Da);CHKERRQ(ierr); 413c2fc9fa9SBarry Smith ierr = VecDuplicate(snes->vec_sol, &vi->Db);CHKERRQ(ierr); 414c2fc9fa9SBarry Smith ierr = VecDuplicate(snes->vec_sol, &vi->z);CHKERRQ(ierr); 415c2fc9fa9SBarry Smith ierr = VecDuplicate(snes->vec_sol, &vi->t);CHKERRQ(ierr); 416c2fc9fa9SBarry Smith PetscFunctionReturn(0); 417c2fc9fa9SBarry Smith } 418c2fc9fa9SBarry Smith /* -------------------------------------------------------------------------- */ 419c2fc9fa9SBarry Smith #undef __FUNCT__ 420f450aa47SBarry Smith #define __FUNCT__ "SNESReset_VINEWTONSSLS" 421f450aa47SBarry Smith PetscErrorCode SNESReset_VINEWTONSSLS(SNES snes) 422c2fc9fa9SBarry Smith { 423f450aa47SBarry Smith SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS*) snes->data; 424c2fc9fa9SBarry Smith PetscErrorCode ierr; 425c2fc9fa9SBarry Smith 426c2fc9fa9SBarry Smith PetscFunctionBegin; 427c2fc9fa9SBarry Smith ierr = SNESReset_VI(snes);CHKERRQ(ierr); 428c2fc9fa9SBarry Smith ierr = VecDestroy(&vi->dpsi);CHKERRQ(ierr); 429c2fc9fa9SBarry Smith ierr = VecDestroy(&vi->phi);CHKERRQ(ierr); 430c2fc9fa9SBarry Smith ierr = VecDestroy(&vi->Da);CHKERRQ(ierr); 431c2fc9fa9SBarry Smith ierr = VecDestroy(&vi->Db);CHKERRQ(ierr); 432c2fc9fa9SBarry Smith ierr = VecDestroy(&vi->z);CHKERRQ(ierr); 433c2fc9fa9SBarry Smith ierr = VecDestroy(&vi->t);CHKERRQ(ierr); 434c2fc9fa9SBarry Smith PetscFunctionReturn(0); 435c2fc9fa9SBarry Smith } 436c2fc9fa9SBarry Smith 437c2fc9fa9SBarry Smith /* -------------------------------------------------------------------------- */ 438c2fc9fa9SBarry Smith /* 439f450aa47SBarry Smith SNESSetFromOptions_VINEWTONSSLS - Sets various parameters for the SNESVI method. 440c2fc9fa9SBarry Smith 441c2fc9fa9SBarry Smith Input Parameter: 442c2fc9fa9SBarry Smith . snes - the SNES context 443c2fc9fa9SBarry Smith 444c2fc9fa9SBarry Smith Application Interface Routine: SNESSetFromOptions() 445c2fc9fa9SBarry Smith */ 446c2fc9fa9SBarry Smith #undef __FUNCT__ 447f450aa47SBarry Smith #define __FUNCT__ "SNESSetFromOptions_VINEWTONSSLS" 4484416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_VINEWTONSSLS(PetscOptionItems *PetscOptionsObject,SNES snes) 449c2fc9fa9SBarry Smith { 450c2fc9fa9SBarry Smith PetscErrorCode ierr; 451f1c6b773SPeter Brune SNESLineSearch linesearch; 452c2fc9fa9SBarry Smith 453c2fc9fa9SBarry Smith PetscFunctionBegin; 454e55864a3SBarry Smith ierr = SNESSetFromOptions_VI(PetscOptionsObject,snes);CHKERRQ(ierr); 455e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SNES semismooth method options");CHKERRQ(ierr); 456c2fc9fa9SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 4579bd66eb0SPeter Brune /* set up the default line search */ 4589bd66eb0SPeter Brune if (!snes->linesearch) { 4597601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 4601a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);CHKERRQ(ierr); 461639ea3faSPeter Brune ierr = SNESLineSearchBTSetAlpha(linesearch, 0.0);CHKERRQ(ierr); 4629bd66eb0SPeter Brune } 463c2fc9fa9SBarry Smith PetscFunctionReturn(0); 464c2fc9fa9SBarry Smith } 465c2fc9fa9SBarry Smith 466c2fc9fa9SBarry Smith 467c2fc9fa9SBarry Smith /* -------------------------------------------------------------------------- */ 468c2fc9fa9SBarry Smith /*MC 469f450aa47SBarry Smith SNESVINEWTONSSLS - Semi-smooth solver for variational inequalities based on Newton's method 470c2fc9fa9SBarry Smith 47161589011SJed Brown Options Database: 472d5f1b7e6SEd Bueler + -snes_type <vinewtonssls,vinewtonrsls> a semi-smooth solver, a reduced space active set method 47361589011SJed Brown - -snes_vi_monitor - prints the number of active constraints at each iteration. 47461589011SJed Brown 475c2fc9fa9SBarry Smith Level: beginner 476c2fc9fa9SBarry Smith 477b80f3ac1SShri Abhyankar References: 478*96a0c994SBarry Smith + 1. - T. S. Munson, F. Facchinei, M. C. Ferris, A. Fischer, and C. Kanzow. The semismooth 479b80f3ac1SShri Abhyankar algorithm for large scale complementarity problems. INFORMS Journal on Computing, 13 (2001). 480*96a0c994SBarry Smith - 2. - T. S. Munson, and S. Benson. Flexible Complementarity Solvers for Large Scale 481d5f1b7e6SEd Bueler Applications, Optimization Methods and Software, 21 (2006). 482b80f3ac1SShri Abhyankar 483d5f1b7e6SEd Bueler .seealso: SNESVISetVariableBounds(), SNESVISetComputeVariableBounds(), SNESCreate(), SNES, SNESSetType(), SNESVINEWTONRSLS, SNESNEWTONTR, SNESLineSearchSet(),SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck() 484c2fc9fa9SBarry Smith 485c2fc9fa9SBarry Smith M*/ 486c2fc9fa9SBarry Smith #undef __FUNCT__ 487f450aa47SBarry Smith #define __FUNCT__ "SNESCreate_VINEWTONSSLS" 4888cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_VINEWTONSSLS(SNES snes) 489c2fc9fa9SBarry Smith { 490c2fc9fa9SBarry Smith PetscErrorCode ierr; 491f450aa47SBarry Smith SNES_VINEWTONSSLS *vi; 492c2fc9fa9SBarry Smith 493c2fc9fa9SBarry Smith PetscFunctionBegin; 494f450aa47SBarry Smith snes->ops->reset = SNESReset_VINEWTONSSLS; 495f450aa47SBarry Smith snes->ops->setup = SNESSetUp_VINEWTONSSLS; 496f450aa47SBarry Smith snes->ops->solve = SNESSolve_VINEWTONSSLS; 497c2fc9fa9SBarry Smith snes->ops->destroy = SNESDestroy_VI; 498f450aa47SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_VINEWTONSSLS; 4990298fd71SBarry Smith snes->ops->view = NULL; 500c2fc9fa9SBarry Smith 501c2fc9fa9SBarry Smith snes->usesksp = PETSC_TRUE; 502c2fc9fa9SBarry Smith snes->usespc = PETSC_FALSE; 503c2fc9fa9SBarry Smith 504b00a9115SJed Brown ierr = PetscNewLog(snes,&vi);CHKERRQ(ierr); 505c2fc9fa9SBarry Smith snes->data = (void*)vi; 506c2fc9fa9SBarry Smith 507bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESVISetVariableBounds_C",SNESVISetVariableBounds_VI);CHKERRQ(ierr); 508bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESVISetComputeVariableBounds_C",SNESVISetComputeVariableBounds_VI);CHKERRQ(ierr); 509c2fc9fa9SBarry Smith PetscFunctionReturn(0); 510c2fc9fa9SBarry Smith } 511c2fc9fa9SBarry Smith 512