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