1 #include "petscvec.h" 2 #include "taosolver.h" 3 #include "tao-private/taolinesearch_impl.h" 4 #include "armijo.h" 5 6 #define REPLACE_FIFO 1 7 #define REPLACE_MRU 2 8 9 #define REFERENCE_MAX 1 10 #define REFERENCE_AVE 2 11 #define REFERENCE_MEAN 3 12 13 #undef __FUNCT__ 14 #define __FUNCT__ "TaoLineSearchDestroy_Armijo" 15 static PetscErrorCode TaoLineSearchDestroy_Armijo(TaoLineSearch ls) 16 { 17 TAOLINESEARCH_ARMIJO_CTX *armP = (TAOLINESEARCH_ARMIJO_CTX *)ls->data; 18 PetscErrorCode ierr; 19 20 PetscFunctionBegin; 21 22 if (armP->memory != PETSC_NULL) { 23 ierr = PetscFree(armP->memory); CHKERRQ(ierr); 24 armP->memory = PETSC_NULL; 25 } 26 if (armP->x) { 27 ierr = PetscObjectDereference((PetscObject)armP->x); CHKERRQ(ierr); 28 } 29 ierr = VecDestroy(&armP->work); CHKERRQ(ierr); 30 ierr = PetscFree(ls->data); CHKERRQ(ierr); 31 ls->data = PETSC_NULL; 32 PetscFunctionReturn(0); 33 } 34 35 #undef __FUNCT__ 36 #define __FUNCT__ "TaoLineSearchReset_Armijo" 37 static PetscErrorCode TaoLineSearchReset_Armijo(TaoLineSearch ls) 38 { 39 TAOLINESEARCH_ARMIJO_CTX *armP = (TAOLINESEARCH_ARMIJO_CTX *)ls->data; 40 PetscErrorCode ierr; 41 42 PetscFunctionBegin; 43 44 if (armP->memory != PETSC_NULL) { 45 ierr = PetscFree(armP->memory); CHKERRQ(ierr); 46 armP->memory = PETSC_NULL; 47 } 48 armP->memorySetup = PETSC_FALSE; 49 PetscFunctionReturn(0); 50 } 51 52 #undef __FUNCT__ 53 #define __FUNCT__ "TaoLineSearchSetFromOptions_Armijo" 54 static PetscErrorCode TaoLineSearchSetFromOptions_Armijo(TaoLineSearch ls) 55 { 56 TAOLINESEARCH_ARMIJO_CTX *armP = (TAOLINESEARCH_ARMIJO_CTX *)ls->data; 57 PetscErrorCode ierr; 58 59 PetscFunctionBegin; 60 ierr = PetscOptionsHead("Armijo linesearch options");CHKERRQ(ierr); 61 ierr = PetscOptionsReal("-tao_ls_armijo_alpha", "initial reference constant", "", armP->alpha, &armP->alpha, 0); CHKERRQ(ierr); 62 ierr = PetscOptionsReal("-tao_ls_armijo_beta_inf", "decrease constant one", "", armP->beta_inf, &armP->beta_inf, 0); CHKERRQ(ierr); 63 ierr = PetscOptionsReal("-tao_ls_armijo_beta", "decrease constant", "", armP->beta, &armP->beta, 0); CHKERRQ(ierr); 64 ierr = PetscOptionsReal("-tao_ls_armijo_sigma", "acceptance constant", "", armP->sigma, &armP->sigma, 0); CHKERRQ(ierr); 65 ierr = PetscOptionsInt("-tao_ls_armijo_memory_size", "number of historical elements", "", armP->memorySize, &armP->memorySize, 0); CHKERRQ(ierr); 66 ierr = PetscOptionsInt("-tao_ls_armijo_reference_policy", "policy for updating reference value", "", armP->referencePolicy, &armP->referencePolicy, 0); CHKERRQ(ierr); 67 ierr = PetscOptionsInt("-tao_ls_armijo_replacement_policy", "policy for updating memory", "", armP->replacementPolicy, &armP->replacementPolicy, 0); CHKERRQ(ierr); 68 ierr = PetscOptionsBool("-tao_ls_armijo_nondescending","Use nondescending armijo algorithm","",armP->nondescending,&armP->nondescending, 0); CHKERRQ(ierr); 69 ierr = PetscOptionsTail();CHKERRQ(ierr); 70 PetscFunctionReturn(0); 71 } 72 73 #undef __FUNCT__ 74 #define __FUNCT__ "TaoLineSearchView_Armijo" 75 static PetscErrorCode TaoLineSearchView_Armijo(TaoLineSearch ls, PetscViewer pv) 76 { 77 TAOLINESEARCH_ARMIJO_CTX *armP = (TAOLINESEARCH_ARMIJO_CTX *)ls->data; 78 PetscBool isascii; 79 PetscErrorCode ierr; 80 81 PetscFunctionBegin; 82 ierr = PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &isascii); CHKERRQ(ierr); 83 if (isascii) { 84 ierr = PetscViewerASCIIPrintf(pv," maxf=%D, ftol=%G, gtol=%G\n",ls->max_funcs, ls->rtol, ls->ftol); CHKERRQ(ierr); 85 ierr=PetscViewerASCIIPrintf(pv," Armijo linesearch",armP->alpha);CHKERRQ(ierr); 86 if (armP->nondescending) { 87 ierr = PetscViewerASCIIPrintf(pv, " (nondescending)"); CHKERRQ(ierr); 88 } 89 if (ls->bounded) { 90 ierr = PetscViewerASCIIPrintf(pv," (projected)"); CHKERRQ(ierr); 91 } 92 ierr=PetscViewerASCIIPrintf(pv,": alpha=%G beta=%G ",armP->alpha,armP->beta);CHKERRQ(ierr); 93 ierr=PetscViewerASCIIPrintf(pv,"sigma=%G ",armP->sigma);CHKERRQ(ierr); 94 ierr=PetscViewerASCIIPrintf(pv,"memsize=%D\n",armP->memorySize);CHKERRQ(ierr); 95 } else { 96 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for Armijo TaoLineSearch",((PetscObject)pv)->type_name); 97 } 98 PetscFunctionReturn(0); 99 } 100 101 #undef __FUNCT__ 102 #define __FUNCT__ "TaoLineSearchApply_Armijo" 103 /* @ TaoApply_Armijo - This routine performs a linesearch. It 104 backtracks until the (nonmonotone) Armijo conditions are satisfied. 105 106 Input Parameters: 107 + tao - TaoSolver context 108 . X - current iterate (on output X contains new iterate, X + step*S) 109 . S - search direction 110 . f - merit function evaluated at X 111 . G - gradient of merit function evaluated at X 112 . W - work vector 113 - step - initial estimate of step length 114 115 Output parameters: 116 + f - merit function evaluated at new iterate, X + step*S 117 . G - gradient of merit function evaluated at new iterate, X + step*S 118 . X - new iterate 119 - step - final step length 120 121 Info is set to one of: 122 . 0 - the line search succeeds; the sufficient decrease 123 condition and the directional derivative condition hold 124 125 negative number if an input parameter is invalid 126 - -1 - step < 0 127 128 positive number > 1 if the line search otherwise terminates 129 + 1 - Step is at the lower bound, stepmin. 130 @ */ 131 static PetscErrorCode TaoLineSearchApply_Armijo(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s) 132 { 133 TAOLINESEARCH_ARMIJO_CTX *armP = (TAOLINESEARCH_ARMIJO_CTX *)ls->data; 134 PetscErrorCode ierr; 135 PetscInt i; 136 PetscReal fact, ref, gdx; 137 PetscInt idx; 138 PetscBool g_computed=PETSC_FALSE; /* to prevent extra gradient computation */ 139 140 PetscFunctionBegin; 141 142 ls->reason = TAOLINESEARCH_CONTINUE_ITERATING; 143 if (!armP->work) { 144 ierr = VecDuplicate(x,&armP->work); CHKERRQ(ierr); 145 armP->x = x; 146 ierr = PetscObjectReference((PetscObject)armP->x); CHKERRQ(ierr); 147 } 148 /* If x has changed, then recreate work */ 149 else if (x != armP->x) { 150 ierr = VecDestroy(&armP->work); CHKERRQ(ierr); 151 ierr = VecDuplicate(x,&armP->work); CHKERRQ(ierr); 152 ierr = PetscObjectDereference((PetscObject)armP->x); CHKERRQ(ierr); 153 armP->x = x; 154 ierr = PetscObjectReference((PetscObject)armP->x); CHKERRQ(ierr); 155 } 156 157 /* Check linesearch parameters */ 158 if (armP->alpha < 1) { 159 ierr = PetscInfo1(ls,"Armijo line search error: alpha (%G) < 1\n", armP->alpha); CHKERRQ(ierr); 160 ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 161 } 162 163 else if ((armP->beta <= 0) || (armP->beta >= 1)) { 164 ierr = PetscInfo1(ls,"Armijo line search error: beta (%G) invalid\n", armP->beta); CHKERRQ(ierr); 165 ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 166 167 } 168 169 else if ((armP->beta_inf <= 0) || (armP->beta_inf >= 1)) { 170 ierr = PetscInfo1(ls,"Armijo line search error: beta_inf (%G) invalid\n", armP->beta_inf); CHKERRQ(ierr); 171 ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 172 } 173 174 else if ((armP->sigma <= 0) || (armP->sigma >= 0.5)) { 175 ierr = PetscInfo1(ls,"Armijo line search error: sigma (%G) invalid\n", armP->sigma); CHKERRQ(ierr); 176 ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 177 } 178 179 else if (armP->memorySize < 1) { 180 ierr = PetscInfo1(ls,"Armijo line search error: memory_size (%D) < 1\n", armP->memorySize); CHKERRQ(ierr); 181 ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 182 } 183 184 else if ((armP->referencePolicy != REFERENCE_MAX) && 185 (armP->referencePolicy != REFERENCE_AVE) && 186 (armP->referencePolicy != REFERENCE_MEAN)) { 187 ierr = PetscInfo(ls,"Armijo line search error: reference_policy invalid\n"); CHKERRQ(ierr); 188 ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 189 190 } 191 192 else if ((armP->replacementPolicy != REPLACE_FIFO) && 193 (armP->replacementPolicy != REPLACE_MRU)) { 194 ierr = PetscInfo(ls,"Armijo line search error: replacement_policy invalid\n"); CHKERRQ(ierr); 195 ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 196 } 197 198 else if (PetscIsInfOrNanReal(*f)) { 199 ierr = PetscInfo(ls,"Armijo line search error: initial function inf or nan\n"); CHKERRQ(ierr); 200 ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 201 } 202 203 if (ls->reason != TAOLINESEARCH_CONTINUE_ITERATING) { 204 PetscFunctionReturn(0); 205 } 206 207 /* Check to see of the memory has been allocated. If not, allocate 208 the historical array and populate it with the initial function 209 values. */ 210 if (armP->memory == PETSC_NULL) { 211 ierr = PetscMalloc(sizeof(PetscReal)*armP->memorySize, &armP->memory ); CHKERRQ(ierr); 212 } 213 214 if (!armP->memorySetup) { 215 for (i = 0; i < armP->memorySize; i++) { 216 armP->memory[i] = armP->alpha*(*f); 217 } 218 219 armP->current = 0; 220 armP->lastReference = armP->memory[0]; 221 armP->memorySetup=PETSC_TRUE; 222 } 223 224 /* Calculate reference value (MAX) */ 225 ref = armP->memory[0]; 226 idx = 0; 227 228 for (i = 1; i < armP->memorySize; i++) { 229 if (armP->memory[i] > ref) { 230 ref = armP->memory[i]; 231 idx = i; 232 } 233 } 234 235 if (armP->referencePolicy == REFERENCE_AVE) { 236 ref = 0; 237 for (i = 0; i < armP->memorySize; i++) { 238 ref += armP->memory[i]; 239 } 240 ref = ref / armP->memorySize; 241 ref = PetscMax(ref, armP->memory[armP->current]); 242 } 243 else if (armP->referencePolicy == REFERENCE_MEAN) { 244 ref = PetscMin(ref, 0.5*(armP->lastReference + armP->memory[armP->current])); 245 } 246 ierr = VecDot(g,s,&gdx); CHKERRQ(ierr); 247 248 if (PetscIsInfOrNanReal(gdx)) { 249 ierr = PetscInfo1(ls,"Initial Line Search step * g is Inf or Nan (%G)\n",gdx); CHKERRQ(ierr); 250 ls->reason=TAOLINESEARCH_FAILED_INFORNAN; 251 PetscFunctionReturn(0); 252 } 253 if (gdx >= 0.0) { 254 ierr = PetscInfo1(ls,"Initial Line Search step is not descent direction (g's=%G)\n",gdx); CHKERRQ(ierr); 255 ls->reason = TAOLINESEARCH_FAILED_ASCENT; 256 PetscFunctionReturn(0); 257 } 258 259 if (armP->nondescending) { 260 fact = armP->sigma; 261 } else { 262 fact = armP->sigma * gdx; 263 } 264 ls->step = ls->initstep; 265 while (ls->step >= ls->stepmin && (ls->nfeval+ls->nfgeval) < ls->max_funcs) { 266 /* Calculate iterate */ 267 ierr = VecCopy(x,armP->work); CHKERRQ(ierr); 268 ierr = VecAXPY(armP->work,ls->step,s); CHKERRQ(ierr); 269 if (ls->bounded) { 270 ierr = VecMedian(ls->lower,armP->work,ls->upper,armP->work); CHKERRQ(ierr); 271 } 272 273 /* Calculate function at new iterate */ 274 if (ls->hasobjective) { 275 ierr = TaoLineSearchComputeObjective(ls,armP->work,f); CHKERRQ(ierr); 276 g_computed=PETSC_FALSE; 277 } else if (ls->usegts) { 278 ierr = TaoLineSearchComputeObjectiveAndGTS(ls,armP->work,f,&gdx); CHKERRQ(ierr); 279 g_computed=PETSC_FALSE; 280 } else { 281 ierr = TaoLineSearchComputeObjectiveAndGradient(ls,armP->work,f,g); CHKERRQ(ierr); 282 g_computed=PETSC_TRUE; 283 } 284 if (ls->step == ls->initstep) { 285 ls->f_fullstep = *f; 286 } 287 288 if (PetscIsInfOrNanReal(*f)) { 289 ls->step *= armP->beta_inf; 290 } 291 else { 292 /* Check descent condition */ 293 if (armP->nondescending && *f <= ref - ls->step*fact*ref) 294 break; 295 if (!armP->nondescending && *f <= ref + ls->step*fact) { 296 break; 297 } 298 299 ls->step *= armP->beta; 300 } 301 } 302 303 /* Check termination */ 304 if (PetscIsInfOrNanReal(*f)) { 305 ierr = PetscInfo(ls, "Function is inf or nan.\n"); CHKERRQ(ierr); 306 ls->reason = TAOLINESEARCH_FAILED_INFORNAN; 307 } else if (ls->step < ls->stepmin) { 308 ierr = PetscInfo(ls, "Step length is below tolerance.\n"); CHKERRQ(ierr); 309 ls->reason = TAOLINESEARCH_HALTED_RTOL; 310 } else if ((ls->nfeval+ls->nfgeval) >= ls->max_funcs) { 311 ierr = PetscInfo2(ls, "Number of line search function evals (%D) > maximum allowed (%D)\n",ls->nfeval+ls->nfgeval, ls->max_funcs); CHKERRQ(ierr); 312 ls->reason = TAOLINESEARCH_HALTED_MAXFCN; 313 } 314 if (ls->reason) { 315 PetscFunctionReturn(0); 316 } 317 318 /* Successful termination, update memory */ 319 armP->lastReference = ref; 320 if (armP->replacementPolicy == REPLACE_FIFO) { 321 armP->memory[armP->current++] = *f; 322 if (armP->current >= armP->memorySize) { 323 armP->current = 0; 324 } 325 } else { 326 armP->current = idx; 327 armP->memory[idx] = *f; 328 } 329 330 /* Update iterate and compute gradient */ 331 ierr = VecCopy(armP->work,x); CHKERRQ(ierr); 332 if (!g_computed) { 333 ierr = TaoLineSearchComputeGradient(ls, x, g); CHKERRQ(ierr); 334 } 335 336 /* Finish computations */ 337 ierr = PetscInfo2(ls, "%D function evals in line search, step = %G\n",ls->nfeval, ls->step); CHKERRQ(ierr); 338 339 PetscFunctionReturn(0); 340 } 341 342 EXTERN_C_BEGIN 343 #undef __FUNCT__ 344 #define __FUNCT__ "TaoLineSearchCreate_Armijo" 345 PetscErrorCode TaoLineSearchCreate_Armijo(TaoLineSearch ls) 346 { 347 TAOLINESEARCH_ARMIJO_CTX *armP; 348 PetscErrorCode ierr; 349 350 PetscFunctionBegin; 351 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 352 ierr = PetscNewLog(ls,&armP);CHKERRQ(ierr); 353 354 armP->memory = PETSC_NULL; 355 armP->alpha = 1.0; 356 armP->beta = 0.5; 357 armP->beta_inf = 0.5; 358 armP->sigma = 1e-4; 359 armP->memorySize = 1; 360 armP->referencePolicy = REFERENCE_MAX; 361 armP->replacementPolicy = REPLACE_MRU; 362 armP->nondescending=PETSC_FALSE; 363 ls->data = (void*)armP; 364 ls->initstep=1.0; 365 ls->ops->setup=0; 366 ls->ops->apply=TaoLineSearchApply_Armijo; 367 ls->ops->view = TaoLineSearchView_Armijo; 368 ls->ops->destroy = TaoLineSearchDestroy_Armijo; 369 ls->ops->reset = TaoLineSearchReset_Armijo; 370 ls->ops->setfromoptions = TaoLineSearchSetFromOptions_Armijo; 371 372 PetscFunctionReturn(0); 373 } 374 EXTERN_C_END 375