1 #include <../src/tao/complementarity/impls/ssls/ssls.h> 2 3 PetscErrorCode TaoSetUp_SSILS(Tao tao) 4 { 5 TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 6 PetscErrorCode ierr; 7 8 PetscFunctionBegin; 9 ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); 10 ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr); 11 ierr = VecDuplicate(tao->solution,&ssls->ff);CHKERRQ(ierr); 12 ierr = VecDuplicate(tao->solution,&ssls->dpsi);CHKERRQ(ierr); 13 ierr = VecDuplicate(tao->solution,&ssls->da);CHKERRQ(ierr); 14 ierr = VecDuplicate(tao->solution,&ssls->db);CHKERRQ(ierr); 15 ierr = VecDuplicate(tao->solution,&ssls->t1);CHKERRQ(ierr); 16 ierr = VecDuplicate(tao->solution,&ssls->t2);CHKERRQ(ierr); 17 PetscFunctionReturn(0); 18 } 19 20 PetscErrorCode TaoDestroy_SSILS(Tao tao) 21 { 22 TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 23 PetscErrorCode ierr; 24 25 PetscFunctionBegin; 26 ierr = VecDestroy(&ssls->ff);CHKERRQ(ierr); 27 ierr = VecDestroy(&ssls->dpsi);CHKERRQ(ierr); 28 ierr = VecDestroy(&ssls->da);CHKERRQ(ierr); 29 ierr = VecDestroy(&ssls->db);CHKERRQ(ierr); 30 ierr = VecDestroy(&ssls->t1);CHKERRQ(ierr); 31 ierr = VecDestroy(&ssls->t2);CHKERRQ(ierr); 32 ierr = PetscFree(tao->data);CHKERRQ(ierr); 33 PetscFunctionReturn(0); 34 } 35 36 static PetscErrorCode TaoSolve_SSILS(Tao tao) 37 { 38 TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 39 PetscReal psi, ndpsi, normd, innerd, t=0; 40 PetscReal delta, rho; 41 TaoConvergedReason reason; 42 TaoLineSearchConvergedReason ls_reason; 43 PetscErrorCode ierr; 44 45 PetscFunctionBegin; 46 /* Assume that Setup has been called! 47 Set the structure for the Jacobian and create a linear solver. */ 48 delta = ssls->delta; 49 rho = ssls->rho; 50 51 ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 52 ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 53 ierr = TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_SSLS_FunctionGradient,tao);CHKERRQ(ierr); 54 ierr = TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao);CHKERRQ(ierr); 55 56 /* Calculate the function value and fischer function value at the 57 current iterate */ 58 ierr = TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,ssls->dpsi);CHKERRQ(ierr); 59 ierr = VecNorm(ssls->dpsi,NORM_2,&ndpsi);CHKERRQ(ierr); 60 61 while (PETSC_TRUE) { 62 ierr=PetscInfo3(tao, "iter: %D, merit: %g, ndpsi: %g\n",tao->niter, (double)ssls->merit, (double)ndpsi);CHKERRQ(ierr); 63 /* Check the termination criteria */ 64 ierr = TaoMonitor(tao,tao->niter,ssls->merit,ndpsi,0.0,t,&reason);CHKERRQ(ierr); 65 if (reason!=TAO_CONTINUE_ITERATING) break; 66 tao->niter++; 67 68 /* Calculate direction. (Really negative of newton direction. Therefore, 69 rest of the code uses -d.) */ 70 ierr = KSPSetOperators(tao->ksp,tao->jacobian,tao->jacobian_pre);CHKERRQ(ierr); 71 ierr = KSPSolve(tao->ksp,ssls->ff,tao->stepdirection);CHKERRQ(ierr); 72 ierr = KSPGetIterationNumber(tao->ksp,&tao->ksp_its);CHKERRQ(ierr); 73 tao->ksp_tot_its+=tao->ksp_its; 74 ierr = VecNorm(tao->stepdirection,NORM_2,&normd);CHKERRQ(ierr); 75 ierr = VecDot(tao->stepdirection,ssls->dpsi,&innerd);CHKERRQ(ierr); 76 77 /* Make sure that we have a descent direction */ 78 if (innerd <= delta*pow(normd, rho)) { 79 ierr = PetscInfo(tao, "newton direction not descent\n");CHKERRQ(ierr); 80 ierr = VecCopy(ssls->dpsi,tao->stepdirection);CHKERRQ(ierr); 81 ierr = VecDot(tao->stepdirection,ssls->dpsi,&innerd);CHKERRQ(ierr); 82 } 83 84 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 85 innerd = -innerd; 86 87 ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr); 88 ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&psi,ssls->dpsi,tao->stepdirection,&t,&ls_reason);CHKERRQ(ierr); 89 ierr = VecNorm(ssls->dpsi,NORM_2,&ndpsi);CHKERRQ(ierr); 90 } 91 PetscFunctionReturn(0); 92 } 93 94 /* ---------------------------------------------------------- */ 95 /*MC 96 TAOSSILS - semi-smooth infeasible linesearch algorithm for solving 97 complementarity constraints 98 99 Options Database Keys: 100 + -tao_ssls_delta - descent test fraction 101 - -tao_ssls_rho - descent test power 102 103 Level: beginner 104 M*/ 105 PETSC_EXTERN PetscErrorCode TaoCreate_SSILS(Tao tao) 106 { 107 TAO_SSLS *ssls; 108 PetscErrorCode ierr; 109 const char *armijo_type = TAOLINESEARCHARMIJO; 110 111 PetscFunctionBegin; 112 ierr = PetscNewLog(tao,&ssls);CHKERRQ(ierr); 113 tao->data = (void*)ssls; 114 tao->ops->solve=TaoSolve_SSILS; 115 tao->ops->setup=TaoSetUp_SSILS; 116 tao->ops->view=TaoView_SSLS; 117 tao->ops->setfromoptions=TaoSetFromOptions_SSLS; 118 tao->ops->destroy = TaoDestroy_SSILS; 119 120 ssls->delta = 1e-10; 121 ssls->rho = 2.1; 122 123 ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 124 ierr = TaoLineSearchSetType(tao->linesearch,armijo_type);CHKERRQ(ierr); 125 ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 126 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 127 /* Note: linesearch objective and objectivegradient routines are set in solve routine */ 128 ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr); 129 ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr); 130 131 /* Override default settings (unless already changed) */ 132 if (!tao->max_it_changed) tao->max_it = 2000; 133 if (!tao->max_funcs_changed) tao->max_funcs = 4000; 134 if (!tao->gttol_changed) tao->gttol = 0; 135 if (!tao->grtol_changed) tao->grtol = 0; 136 #if defined(PETSC_USE_REAL_SINGLE) 137 if (!tao->gatol_changed) tao->gatol = 1.0e-6; 138 if (!tao->fmin_changed) tao->fmin = 1.0e-4; 139 #else 140 if (!tao->gatol_changed) tao->gatol = 1.0e-16; 141 if (!tao->fmin_changed) tao->fmin = 1.0e-8; 142 #endif 143 PetscFunctionReturn(0); 144 } 145 146 147