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