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