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