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 PetscInt iter=0,kspits; 48 TaoConvergedReason reason; 49 TaoLineSearchConvergedReason ls_reason; 50 PetscErrorCode ierr; 51 52 PetscFunctionBegin; 53 /* Assume that Setup has been called! 54 Set the structure for the Jacobian and create a linear solver. */ 55 delta = ssls->delta; 56 rho = ssls->rho; 57 58 ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 59 ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 60 ierr = TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_SSLS_FunctionGradient,tao);CHKERRQ(ierr); 61 ierr = TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao);CHKERRQ(ierr); 62 63 /* Calculate the function value and fischer function value at the 64 current iterate */ 65 ierr = TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,ssls->dpsi);CHKERRQ(ierr); 66 ierr = VecNorm(ssls->dpsi,NORM_2,&ndpsi);CHKERRQ(ierr); 67 68 while (1) { 69 ierr=PetscInfo3(tao, "iter: %D, merit: %g, ndpsi: %g\n",iter, (double)ssls->merit, (double)ndpsi);CHKERRQ(ierr); 70 /* Check the termination criteria */ 71 ierr = TaoMonitor(tao,iter++,ssls->merit,ndpsi,0.0,t,&reason);CHKERRQ(ierr); 72 if (reason!=TAO_CONTINUE_ITERATING) break; 73 74 /* Calculate direction. (Really negative of newton direction. Therefore, 75 rest of the code uses -d.) */ 76 ierr = KSPSetOperators(tao->ksp,tao->jacobian,tao->jacobian_pre);CHKERRQ(ierr); 77 ierr = KSPSolve(tao->ksp,ssls->ff,tao->stepdirection);CHKERRQ(ierr); 78 ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 79 tao->ksp_its+=kspits; 80 ierr = VecNorm(tao->stepdirection,NORM_2,&normd);CHKERRQ(ierr); 81 ierr = VecDot(tao->stepdirection,ssls->dpsi,&innerd);CHKERRQ(ierr); 82 83 /* Make sure that we have a descent direction */ 84 if (innerd <= delta*pow(normd, rho)) { 85 ierr = PetscInfo(tao, "newton direction not descent\n");CHKERRQ(ierr); 86 ierr = VecCopy(ssls->dpsi,tao->stepdirection);CHKERRQ(ierr); 87 ierr = VecDot(tao->stepdirection,ssls->dpsi,&innerd);CHKERRQ(ierr); 88 } 89 90 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 91 innerd = -innerd; 92 93 ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0); 94 ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&psi,ssls->dpsi,tao->stepdirection,&t,&ls_reason);CHKERRQ(ierr); 95 ierr = VecNorm(ssls->dpsi,NORM_2,&ndpsi);CHKERRQ(ierr); 96 } 97 PetscFunctionReturn(0); 98 } 99 100 /* ---------------------------------------------------------- */ 101 /*MC 102 TAOSSILS - semi-smooth infeasible linesearch algorithm for solving 103 complementarity constraints 104 105 Options Database Keys: 106 + -tao_ssls_delta - descent test fraction 107 - -tao_ssls_rho - descent test power 108 109 Level: beginner 110 M*/ 111 EXTERN_C_BEGIN 112 #undef __FUNCT__ 113 #define __FUNCT__ "TaoCreate_SSILS" 114 PetscErrorCode TaoCreate_SSILS(Tao tao) 115 { 116 TAO_SSLS *ssls; 117 PetscErrorCode ierr; 118 const char *armijo_type = TAOLINESEARCHARMIJO; 119 120 PetscFunctionBegin; 121 ierr = PetscNewLog(tao,&ssls);CHKERRQ(ierr); 122 tao->data = (void*)ssls; 123 tao->ops->solve=TaoSolve_SSILS; 124 tao->ops->setup=TaoSetUp_SSILS; 125 tao->ops->view=TaoView_SSLS; 126 tao->ops->setfromoptions=TaoSetFromOptions_SSLS; 127 tao->ops->destroy = TaoDestroy_SSILS; 128 129 ssls->delta = 1e-10; 130 ssls->rho = 2.1; 131 132 ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 133 ierr = TaoLineSearchSetType(tao->linesearch,armijo_type);CHKERRQ(ierr); 134 ierr = TaoLineSearchSetFromOptions(tao->linesearch); 135 /* Note: linesearch objective and objectivegradient routines are set in solve routine */ 136 ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);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 EXTERN_C_END 154 155