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