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