1c2fc9fa9SBarry Smith 2c2fc9fa9SBarry Smith #include <../src/snes/impls/vi/ss/vissimpl.h> /*I "petscsnes.h" I*/ 3b45d2f2cSJed Brown #include <../include/petsc-private/kspimpl.h> 4b45d2f2cSJed Brown #include <../include/petsc-private/matimpl.h> 5b45d2f2cSJed Brown #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 52c2fc9fa9SBarry Smith . 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; 64c2fc9fa9SBarry Smith SNES_VISS *vi = (SNES_VISS*)snes->data; 65c2fc9fa9SBarry Smith Vec Xl = snes->xl,Xu = snes->xu,F = snes->vec_func; 66c2fc9fa9SBarry Smith PetscScalar *phi_arr,*x_arr,*f_arr,*l,*u; 67c2fc9fa9SBarry Smith PetscInt i,nlocal; 68c2fc9fa9SBarry Smith 69c2fc9fa9SBarry Smith PetscFunctionBegin; 70c2fc9fa9SBarry Smith ierr = (*vi->computeuserfunction)(snes,X,F,functx);CHKERRQ(ierr); 71c2fc9fa9SBarry Smith ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr); 72c2fc9fa9SBarry Smith ierr = VecGetArray(X,&x_arr);CHKERRQ(ierr); 73c2fc9fa9SBarry Smith ierr = VecGetArray(F,&f_arr);CHKERRQ(ierr); 74c2fc9fa9SBarry Smith ierr = VecGetArray(Xl,&l);CHKERRQ(ierr); 75c2fc9fa9SBarry Smith ierr = VecGetArray(Xu,&u);CHKERRQ(ierr); 76c2fc9fa9SBarry Smith ierr = VecGetArray(phi,&phi_arr);CHKERRQ(ierr); 77c2fc9fa9SBarry Smith 78c2fc9fa9SBarry Smith for (i=0;i < nlocal;i++) { 79c2fc9fa9SBarry Smith if ((PetscRealPart(l[i]) <= SNES_VI_NINF) && (PetscRealPart(u[i]) >= SNES_VI_INF)) { /* no constraints on variable */ 80c2fc9fa9SBarry Smith phi_arr[i] = f_arr[i]; 81c2fc9fa9SBarry Smith } else if (PetscRealPart(l[i]) <= SNES_VI_NINF) { /* upper bound on variable only */ 82c2fc9fa9SBarry Smith phi_arr[i] = -Phi(u[i] - x_arr[i],-f_arr[i]); 83c2fc9fa9SBarry Smith } else if (PetscRealPart(u[i]) >= SNES_VI_INF) { /* lower bound on variable only */ 84c2fc9fa9SBarry Smith phi_arr[i] = Phi(x_arr[i] - l[i],f_arr[i]); 85c2fc9fa9SBarry Smith } else if (l[i] == u[i]) { 86c2fc9fa9SBarry Smith phi_arr[i] = l[i] - x_arr[i]; 87c2fc9fa9SBarry Smith } else { /* both bounds on variable */ 88c2fc9fa9SBarry Smith phi_arr[i] = Phi(x_arr[i] - l[i],-Phi(u[i] - x_arr[i],-f_arr[i])); 89c2fc9fa9SBarry Smith } 90c2fc9fa9SBarry Smith } 91c2fc9fa9SBarry Smith 92c2fc9fa9SBarry Smith ierr = VecRestoreArray(X,&x_arr);CHKERRQ(ierr); 93c2fc9fa9SBarry Smith ierr = VecRestoreArray(F,&f_arr);CHKERRQ(ierr); 94c2fc9fa9SBarry Smith ierr = VecRestoreArray(Xl,&l);CHKERRQ(ierr); 95c2fc9fa9SBarry Smith ierr = VecRestoreArray(Xu,&u);CHKERRQ(ierr); 96c2fc9fa9SBarry Smith ierr = VecRestoreArray(phi,&phi_arr);CHKERRQ(ierr); 97c2fc9fa9SBarry Smith PetscFunctionReturn(0); 98c2fc9fa9SBarry Smith } 99c2fc9fa9SBarry Smith 100c2fc9fa9SBarry Smith /* 101c2fc9fa9SBarry Smith SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the 102c2fc9fa9SBarry Smith the semismooth jacobian. 103c2fc9fa9SBarry Smith */ 104c2fc9fa9SBarry Smith #undef __FUNCT__ 105c2fc9fa9SBarry Smith #define __FUNCT__ "SNESVIComputeBsubdifferentialVectors" 106c2fc9fa9SBarry Smith PetscErrorCode SNESVIComputeBsubdifferentialVectors(SNES snes,Vec X,Vec F,Mat jac,Vec Da,Vec Db) 107c2fc9fa9SBarry Smith { 108c2fc9fa9SBarry Smith PetscErrorCode ierr; 109c2fc9fa9SBarry Smith PetscScalar *l,*u,*x,*f,*da,*db,da1,da2,db1,db2; 110c2fc9fa9SBarry Smith PetscInt i,nlocal; 111c2fc9fa9SBarry Smith 112c2fc9fa9SBarry Smith PetscFunctionBegin; 113c2fc9fa9SBarry Smith 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++) { 123c2fc9fa9SBarry Smith if ((PetscRealPart(l[i]) <= SNES_VI_NINF) && (PetscRealPart(u[i]) >= SNES_VI_INF)) {/* no constraints on variable */ 124c2fc9fa9SBarry Smith da[i] = 0; 125c2fc9fa9SBarry Smith db[i] = 1; 126c2fc9fa9SBarry Smith } else if (PetscRealPart(l[i]) <= SNES_VI_NINF) { /* 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]); 129c2fc9fa9SBarry Smith } else if (PetscRealPart(u[i]) >= SNES_VI_INF) { /* 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 */ 179c2fc9fa9SBarry Smith ierr = MatDiagonalScale(jac,Db,PETSC_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 */ 182c2fc9fa9SBarry Smith ierr = MatDiagonalScale(jac_pre,Db,PETSC_NULL); 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 /* 216c2fc9fa9SBarry Smith SNESSolve_VISS - 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 Output Parameter: 223c2fc9fa9SBarry Smith . outits - number of iterations until termination 224c2fc9fa9SBarry Smith 225c2fc9fa9SBarry Smith Application Interface Routine: SNESSolve() 226c2fc9fa9SBarry Smith 227c2fc9fa9SBarry Smith Notes: 228c2fc9fa9SBarry Smith This implements essentially a semismooth Newton method with a 229c2fc9fa9SBarry Smith line search. The default line search does not do any line seach 230c2fc9fa9SBarry Smith but rather takes a full newton step. 231016d19f9SBarry Smith 232016d19f9SBarry Smith Developer Note: the code in this file should be slightly modified so that this routine need not exist and the SNESSolve_LS() routine is called directly with the appropriate wrapped function and Jacobian evaluations 233016d19f9SBarry Smith 234c2fc9fa9SBarry Smith */ 235c2fc9fa9SBarry Smith #undef __FUNCT__ 236c2fc9fa9SBarry Smith #define __FUNCT__ "SNESSolve_VISS" 237c2fc9fa9SBarry Smith PetscErrorCode SNESSolve_VISS(SNES snes) 238c2fc9fa9SBarry Smith { 239c2fc9fa9SBarry Smith SNES_VISS *vi = (SNES_VISS*)snes->data; 240c2fc9fa9SBarry Smith PetscErrorCode ierr; 241c2fc9fa9SBarry Smith PetscInt maxits,i,lits; 242c2fc9fa9SBarry Smith PetscBool lssucceed; 243c2fc9fa9SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 244c2fc9fa9SBarry Smith PetscReal gnorm,xnorm=0,ynorm; 2459bd66eb0SPeter Brune Vec Y,X,F; 246c2fc9fa9SBarry Smith KSPConvergedReason kspreason; 2476cab3a1bSJed Brown DM dm; 248*942e3340SBarry Smith DMSNES sdm; 249c2fc9fa9SBarry Smith 250c2fc9fa9SBarry Smith PetscFunctionBegin; 2516cab3a1bSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 252*942e3340SBarry Smith ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr); 2536cab3a1bSJed Brown vi->computeuserfunction = sdm->computefunction; 2546cab3a1bSJed Brown sdm->computefunction = SNESVIComputeFunction; 255c2fc9fa9SBarry Smith 256c2fc9fa9SBarry Smith snes->numFailures = 0; 257c2fc9fa9SBarry Smith snes->numLinearSolveFailures = 0; 258c2fc9fa9SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 259c2fc9fa9SBarry Smith 260c2fc9fa9SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 261c2fc9fa9SBarry Smith X = snes->vec_sol; /* solution vector */ 262c2fc9fa9SBarry Smith F = snes->vec_func; /* residual vector */ 263c2fc9fa9SBarry Smith Y = snes->work[0]; /* work vectors */ 264c2fc9fa9SBarry Smith 265c2fc9fa9SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 266c2fc9fa9SBarry Smith snes->iter = 0; 267c2fc9fa9SBarry Smith snes->norm = 0.0; 268c2fc9fa9SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 269c2fc9fa9SBarry Smith 270c2fc9fa9SBarry Smith ierr = SNESVIProjectOntoBounds(snes,X);CHKERRQ(ierr); 271c2fc9fa9SBarry Smith ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr); 272c2fc9fa9SBarry Smith if (snes->domainerror) { 273c2fc9fa9SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 2746cab3a1bSJed Brown sdm->computefunction = vi->computeuserfunction; 275c2fc9fa9SBarry Smith PetscFunctionReturn(0); 276c2fc9fa9SBarry Smith } 277c2fc9fa9SBarry Smith /* Compute Merit function */ 278c2fc9fa9SBarry Smith ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr); 279c2fc9fa9SBarry Smith 280c2fc9fa9SBarry Smith ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 281c2fc9fa9SBarry Smith ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 282c2fc9fa9SBarry Smith if (PetscIsInfOrNanReal(vi->merit)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 283c2fc9fa9SBarry Smith 284c2fc9fa9SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 285c2fc9fa9SBarry Smith snes->norm = vi->phinorm; 286c2fc9fa9SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 287c2fc9fa9SBarry Smith SNESLogConvHistory(snes,vi->phinorm,0); 288c2fc9fa9SBarry Smith ierr = SNESMonitor(snes,0,vi->phinorm);CHKERRQ(ierr); 289c2fc9fa9SBarry Smith 290c2fc9fa9SBarry Smith /* set parameter for default relative tolerance convergence test */ 291c2fc9fa9SBarry Smith snes->ttol = vi->phinorm*snes->rtol; 292c2fc9fa9SBarry Smith /* test convergence */ 293c2fc9fa9SBarry Smith ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 294c2fc9fa9SBarry Smith if (snes->reason) { 2956cab3a1bSJed Brown sdm->computefunction = vi->computeuserfunction; 296c2fc9fa9SBarry Smith PetscFunctionReturn(0); 297c2fc9fa9SBarry Smith } 298c2fc9fa9SBarry Smith 299c2fc9fa9SBarry Smith for (i=0; i<maxits; i++) { 300c2fc9fa9SBarry Smith 301c2fc9fa9SBarry Smith /* Call general purpose update function */ 302c2fc9fa9SBarry Smith if (snes->ops->update) { 303c2fc9fa9SBarry Smith ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 304c2fc9fa9SBarry Smith } 305c2fc9fa9SBarry Smith 306c2fc9fa9SBarry Smith /* Solve J Y = Phi, where J is the semismooth jacobian */ 307b9600128SPeter Brune 308b9600128SPeter Brune /* Get the jacobian -- note that the function must be the original function for snes_fd and snes_fd_color to work for this*/ 309b9600128SPeter Brune sdm->computefunction = vi->computeuserfunction; 310c2fc9fa9SBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 311b9600128SPeter Brune sdm->computefunction = SNESVIComputeFunction; 312b9600128SPeter Brune 313c2fc9fa9SBarry Smith /* Get the diagonal shift and row scaling vectors */ 314c2fc9fa9SBarry Smith ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr); 315c2fc9fa9SBarry Smith /* Compute the semismooth jacobian */ 316c2fc9fa9SBarry Smith ierr = SNESVIComputeJacobian(snes->jacobian,snes->jacobian_pre,vi->Da,vi->Db);CHKERRQ(ierr); 317c2fc9fa9SBarry Smith /* Compute the merit function gradient */ 318c2fc9fa9SBarry Smith ierr = SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);CHKERRQ(ierr); 319c2fc9fa9SBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 320c2fc9fa9SBarry Smith ierr = SNES_KSPSolve(snes,snes->ksp,vi->phi,Y);CHKERRQ(ierr); 321c2fc9fa9SBarry Smith ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 322c2fc9fa9SBarry Smith 323c2fc9fa9SBarry Smith if (kspreason < 0) { 324c2fc9fa9SBarry Smith if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 325c2fc9fa9SBarry 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); 326c2fc9fa9SBarry Smith snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 327c2fc9fa9SBarry Smith break; 328c2fc9fa9SBarry Smith } 329c2fc9fa9SBarry Smith } 330c2fc9fa9SBarry Smith ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 331c2fc9fa9SBarry Smith snes->linear_its += lits; 332c2fc9fa9SBarry Smith ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 333c2fc9fa9SBarry Smith /* 334c2fc9fa9SBarry Smith if (snes->ops->precheckstep) { 335c2fc9fa9SBarry Smith PetscBool changed_y = PETSC_FALSE; 336c2fc9fa9SBarry Smith ierr = (*snes->ops->precheckstep)(snes,X,Y,snes->precheck,&changed_y);CHKERRQ(ierr); 337c2fc9fa9SBarry Smith } 338c2fc9fa9SBarry Smith 339c2fc9fa9SBarry Smith if (PetscLogPrintInfo){ 340c2fc9fa9SBarry Smith ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 341c2fc9fa9SBarry Smith } 342c2fc9fa9SBarry Smith */ 343c2fc9fa9SBarry Smith /* Compute a (scaled) negative update in the line search routine: 344c2fc9fa9SBarry Smith Y <- X - lambda*Y 345c2fc9fa9SBarry Smith and evaluate G = function(Y) (depends on the line search). 346c2fc9fa9SBarry Smith */ 347c2fc9fa9SBarry Smith ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 348c2fc9fa9SBarry Smith ynorm = 1; gnorm = vi->phinorm; 3499bd66eb0SPeter Brune /* ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,vi->phi,Y,vi->phinorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); */ 350f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, vi->phi, &gnorm, Y);CHKERRQ(ierr); 351f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr); 352c2fc9fa9SBarry 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); 353c2fc9fa9SBarry Smith if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 354c2fc9fa9SBarry Smith if (snes->domainerror) { 355c2fc9fa9SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 3566cab3a1bSJed Brown sdm->computefunction = vi->computeuserfunction; 357c2fc9fa9SBarry Smith PetscFunctionReturn(0); 358c2fc9fa9SBarry Smith } 359f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr); 360c2fc9fa9SBarry Smith if (!lssucceed) { 361c2fc9fa9SBarry Smith if (++snes->numFailures >= snes->maxFailures) { 362c2fc9fa9SBarry Smith PetscBool ismin; 363c2fc9fa9SBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 3649bd66eb0SPeter Brune ierr = SNESVICheckLocalMin_Private(snes,snes->jacobian,vi->phi,X,gnorm,&ismin);CHKERRQ(ierr); 365c2fc9fa9SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 366c2fc9fa9SBarry Smith break; 367c2fc9fa9SBarry Smith } 368c2fc9fa9SBarry Smith } 369c2fc9fa9SBarry Smith /* Update function and solution vectors */ 370c2fc9fa9SBarry Smith vi->phinorm = gnorm; 371c2fc9fa9SBarry Smith vi->merit = 0.5*vi->phinorm*vi->phinorm; 372c2fc9fa9SBarry Smith /* Monitor convergence */ 373c2fc9fa9SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 374c2fc9fa9SBarry Smith snes->iter = i+1; 375c2fc9fa9SBarry Smith snes->norm = vi->phinorm; 376c2fc9fa9SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 377c2fc9fa9SBarry Smith SNESLogConvHistory(snes,snes->norm,lits); 378c2fc9fa9SBarry Smith ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 379c2fc9fa9SBarry Smith /* Test for convergence, xnorm = || X || */ 380c2fc9fa9SBarry Smith if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 381c2fc9fa9SBarry Smith ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 382c2fc9fa9SBarry Smith if (snes->reason) break; 383c2fc9fa9SBarry Smith } 384c2fc9fa9SBarry Smith if (i == maxits) { 385c2fc9fa9SBarry Smith ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 386c2fc9fa9SBarry Smith if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 387c2fc9fa9SBarry Smith } 3886cab3a1bSJed Brown sdm->computefunction = vi->computeuserfunction; 389c2fc9fa9SBarry Smith PetscFunctionReturn(0); 390c2fc9fa9SBarry Smith } 391c2fc9fa9SBarry Smith 392c2fc9fa9SBarry Smith /* -------------------------------------------------------------------------- */ 393c2fc9fa9SBarry Smith /* 394c2fc9fa9SBarry Smith SNESSetUp_VISS - Sets up the internal data structures for the later use 395c2fc9fa9SBarry Smith of the SNES nonlinear solver. 396c2fc9fa9SBarry Smith 397c2fc9fa9SBarry Smith Input Parameter: 398c2fc9fa9SBarry Smith . snes - the SNES context 399c2fc9fa9SBarry Smith . x - the solution vector 400c2fc9fa9SBarry Smith 401c2fc9fa9SBarry Smith Application Interface Routine: SNESSetUp() 402c2fc9fa9SBarry Smith 403c2fc9fa9SBarry Smith Notes: 404c2fc9fa9SBarry Smith For basic use of the SNES solvers, the user need not explicitly call 405c2fc9fa9SBarry Smith SNESSetUp(), since these actions will automatically occur during 406c2fc9fa9SBarry Smith the call to SNESSolve(). 407c2fc9fa9SBarry Smith */ 408c2fc9fa9SBarry Smith #undef __FUNCT__ 409c2fc9fa9SBarry Smith #define __FUNCT__ "SNESSetUp_VISS" 410c2fc9fa9SBarry Smith PetscErrorCode SNESSetUp_VISS(SNES snes) 411c2fc9fa9SBarry Smith { 412c2fc9fa9SBarry Smith PetscErrorCode ierr; 413c2fc9fa9SBarry Smith SNES_VISS *vi = (SNES_VISS*) snes->data; 414c2fc9fa9SBarry Smith 415c2fc9fa9SBarry Smith PetscFunctionBegin; 416c2fc9fa9SBarry Smith ierr = SNESSetUp_VI(snes);CHKERRQ(ierr); 417c2fc9fa9SBarry Smith ierr = VecDuplicate(snes->vec_sol, &vi->dpsi);CHKERRQ(ierr); 418c2fc9fa9SBarry Smith ierr = VecDuplicate(snes->vec_sol, &vi->phi);CHKERRQ(ierr); 419c2fc9fa9SBarry Smith ierr = VecDuplicate(snes->vec_sol, &vi->Da);CHKERRQ(ierr); 420c2fc9fa9SBarry Smith ierr = VecDuplicate(snes->vec_sol, &vi->Db);CHKERRQ(ierr); 421c2fc9fa9SBarry Smith ierr = VecDuplicate(snes->vec_sol, &vi->z);CHKERRQ(ierr); 422c2fc9fa9SBarry Smith ierr = VecDuplicate(snes->vec_sol, &vi->t);CHKERRQ(ierr); 423c2fc9fa9SBarry Smith PetscFunctionReturn(0); 424c2fc9fa9SBarry Smith } 425c2fc9fa9SBarry Smith /* -------------------------------------------------------------------------- */ 426c2fc9fa9SBarry Smith #undef __FUNCT__ 427c2fc9fa9SBarry Smith #define __FUNCT__ "SNESReset_VISS" 428c2fc9fa9SBarry Smith PetscErrorCode SNESReset_VISS(SNES snes) 429c2fc9fa9SBarry Smith { 430c2fc9fa9SBarry Smith SNES_VISS *vi = (SNES_VISS*) snes->data; 431c2fc9fa9SBarry Smith PetscErrorCode ierr; 432c2fc9fa9SBarry Smith 433c2fc9fa9SBarry Smith PetscFunctionBegin; 434c2fc9fa9SBarry Smith ierr = SNESReset_VI(snes);CHKERRQ(ierr); 435c2fc9fa9SBarry Smith ierr = VecDestroy(&vi->dpsi);CHKERRQ(ierr); 436c2fc9fa9SBarry Smith ierr = VecDestroy(&vi->phi);CHKERRQ(ierr); 437c2fc9fa9SBarry Smith ierr = VecDestroy(&vi->Da);CHKERRQ(ierr); 438c2fc9fa9SBarry Smith ierr = VecDestroy(&vi->Db);CHKERRQ(ierr); 439c2fc9fa9SBarry Smith ierr = VecDestroy(&vi->z);CHKERRQ(ierr); 440c2fc9fa9SBarry Smith ierr = VecDestroy(&vi->t);CHKERRQ(ierr); 441c2fc9fa9SBarry Smith PetscFunctionReturn(0); 442c2fc9fa9SBarry Smith } 443c2fc9fa9SBarry Smith 444c2fc9fa9SBarry Smith /* -------------------------------------------------------------------------- */ 445c2fc9fa9SBarry Smith /* 446c2fc9fa9SBarry Smith SNESSetFromOptions_VISS - Sets various parameters for the SNESVI method. 447c2fc9fa9SBarry Smith 448c2fc9fa9SBarry Smith Input Parameter: 449c2fc9fa9SBarry Smith . snes - the SNES context 450c2fc9fa9SBarry Smith 451c2fc9fa9SBarry Smith Application Interface Routine: SNESSetFromOptions() 452c2fc9fa9SBarry Smith */ 453c2fc9fa9SBarry Smith #undef __FUNCT__ 454c2fc9fa9SBarry Smith #define __FUNCT__ "SNESSetFromOptions_VISS" 455c2fc9fa9SBarry Smith static PetscErrorCode SNESSetFromOptions_VISS(SNES snes) 456c2fc9fa9SBarry Smith { 457c2fc9fa9SBarry Smith PetscErrorCode ierr; 458f1c6b773SPeter Brune SNESLineSearch linesearch; 459c2fc9fa9SBarry Smith 460c2fc9fa9SBarry Smith PetscFunctionBegin; 461c2fc9fa9SBarry Smith ierr = SNESSetFromOptions_VI(snes);CHKERRQ(ierr); 462c2fc9fa9SBarry Smith ierr = PetscOptionsHead("SNES semismooth method options");CHKERRQ(ierr); 463c2fc9fa9SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 4649bd66eb0SPeter Brune /* set up the default line search */ 4659bd66eb0SPeter Brune if (!snes->linesearch) { 466f1c6b773SPeter Brune ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 4671a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);CHKERRQ(ierr); 468639ea3faSPeter Brune ierr = SNESLineSearchBTSetAlpha(linesearch, 0.0);CHKERRQ(ierr); 4699bd66eb0SPeter Brune } 4709bd66eb0SPeter Brune 471c2fc9fa9SBarry Smith PetscFunctionReturn(0); 472c2fc9fa9SBarry Smith } 473c2fc9fa9SBarry Smith 474c2fc9fa9SBarry Smith 475c2fc9fa9SBarry Smith /* -------------------------------------------------------------------------- */ 476c2fc9fa9SBarry Smith /*MC 477c2fc9fa9SBarry Smith SNESVISS - Semi-smooth solver for variational inequalities based on Newton's method 478c2fc9fa9SBarry Smith 47961589011SJed Brown Options Database: 48061589011SJed Brown + -snes_vi_type <ss,rs,rsaug> a semi-smooth solver, a reduced space active set method, and a reduced space active set method that does not eliminate the active constraints from the Jacobian instead augments the Jacobian with additional variables that enforce the constraints 48161589011SJed Brown - -snes_vi_monitor - prints the number of active constraints at each iteration. 48261589011SJed Brown 483c2fc9fa9SBarry Smith Level: beginner 484c2fc9fa9SBarry Smith 485b80f3ac1SShri Abhyankar References: 486b80f3ac1SShri Abhyankar - T. S. Munson, F. Facchinei, M. C. Ferris, A. Fischer, and C. Kanzow. The semismooth 487b80f3ac1SShri Abhyankar algorithm for large scale complementarity problems. INFORMS Journal on Computing, 13 (2001). 488b80f3ac1SShri Abhyankar 48961589011SJed Brown .seealso: SNESVISetVariableBounds(), SNESVISetComputeVariableBounds(), SNESCreate(), SNES, SNESSetType(), SNESVIRS, SNESVISS, SNESTR, SNESLineSearchSet(), 490c2fc9fa9SBarry Smith SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 491c2fc9fa9SBarry Smith SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams() 492c2fc9fa9SBarry Smith 493c2fc9fa9SBarry Smith M*/ 494c2fc9fa9SBarry Smith EXTERN_C_BEGIN 495c2fc9fa9SBarry Smith #undef __FUNCT__ 496c2fc9fa9SBarry Smith #define __FUNCT__ "SNESCreate_VISS" 497c2fc9fa9SBarry Smith PetscErrorCode SNESCreate_VISS(SNES snes) 498c2fc9fa9SBarry Smith { 499c2fc9fa9SBarry Smith PetscErrorCode ierr; 500c2fc9fa9SBarry Smith SNES_VISS *vi; 501c2fc9fa9SBarry Smith 502c2fc9fa9SBarry Smith PetscFunctionBegin; 503c2fc9fa9SBarry Smith snes->ops->reset = SNESReset_VISS; 504c2fc9fa9SBarry Smith snes->ops->setup = SNESSetUp_VISS; 505c2fc9fa9SBarry Smith snes->ops->solve = SNESSolve_VISS; 506c2fc9fa9SBarry Smith snes->ops->destroy = SNESDestroy_VI; 507c2fc9fa9SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_VISS; 508c2fc9fa9SBarry Smith snes->ops->view = PETSC_NULL; 509c2fc9fa9SBarry Smith 510c2fc9fa9SBarry Smith snes->usesksp = PETSC_TRUE; 511c2fc9fa9SBarry Smith snes->usespc = PETSC_FALSE; 512c2fc9fa9SBarry Smith 513c2fc9fa9SBarry Smith ierr = PetscNewLog(snes,SNES_VISS,&vi);CHKERRQ(ierr); 514c2fc9fa9SBarry Smith snes->data = (void*)vi; 515c2fc9fa9SBarry Smith 51661589011SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESVISetVariableBounds_C","SNESVISetVariableBounds_VI",SNESVISetVariableBounds_VI);CHKERRQ(ierr); 51761589011SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESVISetComputeVariableBounds_C","SNESVISetComputeVariableBounds_VI",SNESVISetComputeVariableBounds_VI);CHKERRQ(ierr); 518c2fc9fa9SBarry Smith PetscFunctionReturn(0); 519c2fc9fa9SBarry Smith } 520c2fc9fa9SBarry Smith EXTERN_C_END 521c2fc9fa9SBarry Smith 522