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