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