1 2 #include <petsc/private/taolinesearchimpl.h> 3 #include <../src/tao/linesearch/impls/owarmijo/owarmijo.h> 4 5 #define REPLACE_FIFO 1 6 #define REPLACE_MRU 2 7 8 #define REFERENCE_MAX 1 9 #define REFERENCE_AVE 2 10 #define REFERENCE_MEAN 3 11 12 static PetscErrorCode ProjWork_OWLQN(Vec w,Vec x,Vec gv,PetscReal *gdx) 13 { 14 const PetscReal *xptr,*gptr; 15 PetscReal *wptr; 16 PetscInt low,high,low1,high1,low2,high2,i; 17 18 PetscFunctionBegin; 19 PetscCall(VecGetOwnershipRange(w,&low,&high)); 20 PetscCall(VecGetOwnershipRange(x,&low1,&high1)); 21 PetscCall(VecGetOwnershipRange(gv,&low2,&high2)); 22 23 *gdx=0.0; 24 PetscCall(VecGetArray(w,&wptr)); 25 PetscCall(VecGetArrayRead(x,&xptr)); 26 PetscCall(VecGetArrayRead(gv,&gptr)); 27 28 for (i=0;i<high-low;i++) { 29 if (xptr[i]*wptr[i]<0.0) wptr[i]=0.0; 30 *gdx = *gdx + gptr[i]*(wptr[i]-xptr[i]); 31 } 32 PetscCall(VecRestoreArray(w,&wptr)); 33 PetscCall(VecRestoreArrayRead(x,&xptr)); 34 PetscCall(VecRestoreArrayRead(gv,&gptr)); 35 PetscFunctionReturn(0); 36 } 37 38 static PetscErrorCode TaoLineSearchDestroy_OWArmijo(TaoLineSearch ls) 39 { 40 TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data; 41 42 PetscFunctionBegin; 43 PetscCall(PetscFree(armP->memory)); 44 if (armP->x) { 45 PetscCall(PetscObjectDereference((PetscObject)armP->x)); 46 } 47 PetscCall(VecDestroy(&armP->work)); 48 PetscCall(PetscFree(ls->data)); 49 PetscFunctionReturn(0); 50 } 51 52 static PetscErrorCode TaoLineSearchSetFromOptions_OWArmijo(PetscOptionItems *PetscOptionsObject,TaoLineSearch ls) 53 { 54 TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data; 55 56 PetscFunctionBegin; 57 PetscOptionsHeadBegin(PetscOptionsObject,"OWArmijo linesearch options"); 58 PetscCall(PetscOptionsReal("-tao_ls_OWArmijo_alpha", "initial reference constant", "", armP->alpha, &armP->alpha,NULL)); 59 PetscCall(PetscOptionsReal("-tao_ls_OWArmijo_beta_inf", "decrease constant one", "", armP->beta_inf, &armP->beta_inf,NULL)); 60 PetscCall(PetscOptionsReal("-tao_ls_OWArmijo_beta", "decrease constant", "", armP->beta, &armP->beta,NULL)); 61 PetscCall(PetscOptionsReal("-tao_ls_OWArmijo_sigma", "acceptance constant", "", armP->sigma, &armP->sigma,NULL)); 62 PetscCall(PetscOptionsInt("-tao_ls_OWArmijo_memory_size", "number of historical elements", "", armP->memorySize, &armP->memorySize,NULL)); 63 PetscCall(PetscOptionsInt("-tao_ls_OWArmijo_reference_policy", "policy for updating reference value", "", armP->referencePolicy, &armP->referencePolicy,NULL)); 64 PetscCall(PetscOptionsInt("-tao_ls_OWArmijo_replacement_policy", "policy for updating memory", "", armP->replacementPolicy, &armP->replacementPolicy,NULL)); 65 PetscCall(PetscOptionsBool("-tao_ls_OWArmijo_nondescending","Use nondescending OWArmijo algorithm","",armP->nondescending,&armP->nondescending,NULL)); 66 PetscOptionsHeadEnd(); 67 PetscFunctionReturn(0); 68 } 69 70 static PetscErrorCode TaoLineSearchView_OWArmijo(TaoLineSearch ls, PetscViewer pv) 71 { 72 TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data; 73 PetscBool isascii; 74 75 PetscFunctionBegin; 76 PetscCall(PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &isascii)); 77 if (isascii) { 78 PetscCall(PetscViewerASCIIPrintf(pv," OWArmijo linesearch")); 79 if (armP->nondescending) { 80 PetscCall(PetscViewerASCIIPrintf(pv, " (nondescending)")); 81 } 82 PetscCall(PetscViewerASCIIPrintf(pv,": alpha=%g beta=%g ",(double)armP->alpha,(double)armP->beta)); 83 PetscCall(PetscViewerASCIIPrintf(pv,"sigma=%g ",(double)armP->sigma)); 84 PetscCall(PetscViewerASCIIPrintf(pv,"memsize=%" PetscInt_FMT "\n",armP->memorySize)); 85 } 86 PetscFunctionReturn(0); 87 } 88 89 /* @ TaoApply_OWArmijo - This routine performs a linesearch. It 90 backtracks until the (nonmonotone) OWArmijo conditions are satisfied. 91 92 Input Parameters: 93 + tao - TAO_SOLVER context 94 . X - current iterate (on output X contains new iterate, X + step*S) 95 . S - search direction 96 . f - merit function evaluated at X 97 . G - gradient of merit function evaluated at X 98 . W - work vector 99 - step - initial estimate of step length 100 101 Output parameters: 102 + f - merit function evaluated at new iterate, X + step*S 103 . G - gradient of merit function evaluated at new iterate, X + step*S 104 . X - new iterate 105 - step - final step length 106 107 Info is set to one of: 108 . 0 - the line search succeeds; the sufficient decrease 109 condition and the directional derivative condition hold 110 111 negative number if an input parameter is invalid 112 - -1 - step < 0 113 114 positive number > 1 if the line search otherwise terminates 115 + 1 - Step is at the lower bound, stepmin. 116 @ */ 117 static PetscErrorCode TaoLineSearchApply_OWArmijo(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s) 118 { 119 TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data; 120 PetscInt i, its=0; 121 PetscReal fact, ref, gdx; 122 PetscInt idx; 123 PetscBool g_computed=PETSC_FALSE; /* to prevent extra gradient computation */ 124 Vec g_old; 125 PetscReal owlqn_minstep=0.005; 126 PetscReal partgdx; 127 MPI_Comm comm; 128 129 PetscFunctionBegin; 130 PetscCall(PetscObjectGetComm((PetscObject)ls,&comm)); 131 fact = 0.0; 132 ls->nfeval=0; 133 ls->reason = TAOLINESEARCH_CONTINUE_ITERATING; 134 if (!armP->work) { 135 PetscCall(VecDuplicate(x,&armP->work)); 136 armP->x = x; 137 PetscCall(PetscObjectReference((PetscObject)armP->x)); 138 } else if (x != armP->x) { 139 PetscCall(VecDestroy(&armP->work)); 140 PetscCall(VecDuplicate(x,&armP->work)); 141 PetscCall(PetscObjectDereference((PetscObject)armP->x)); 142 armP->x = x; 143 PetscCall(PetscObjectReference((PetscObject)armP->x)); 144 } 145 146 PetscCall(TaoLineSearchMonitor(ls, 0, *f, 0.0)); 147 148 /* Check linesearch parameters */ 149 if (armP->alpha < 1) { 150 PetscCall(PetscInfo(ls,"OWArmijo line search error: alpha (%g) < 1\n", (double)armP->alpha)); 151 ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 152 } else if ((armP->beta <= 0) || (armP->beta >= 1)) { 153 PetscCall(PetscInfo(ls,"OWArmijo line search error: beta (%g) invalid\n", (double)armP->beta)); 154 ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 155 } else if ((armP->beta_inf <= 0) || (armP->beta_inf >= 1)) { 156 PetscCall(PetscInfo(ls,"OWArmijo line search error: beta_inf (%g) invalid\n", (double)armP->beta_inf)); 157 ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 158 } else if ((armP->sigma <= 0) || (armP->sigma >= 0.5)) { 159 PetscCall(PetscInfo(ls,"OWArmijo line search error: sigma (%g) invalid\n", (double)armP->sigma)); 160 ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 161 } else if (armP->memorySize < 1) { 162 PetscCall(PetscInfo(ls,"OWArmijo line search error: memory_size (%" PetscInt_FMT ") < 1\n", armP->memorySize)); 163 ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 164 } else if ((armP->referencePolicy != REFERENCE_MAX) && (armP->referencePolicy != REFERENCE_AVE) && (armP->referencePolicy != REFERENCE_MEAN)) { 165 PetscCall(PetscInfo(ls,"OWArmijo line search error: reference_policy invalid\n")); 166 ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 167 } else if ((armP->replacementPolicy != REPLACE_FIFO) && (armP->replacementPolicy != REPLACE_MRU)) { 168 PetscCall(PetscInfo(ls,"OWArmijo line search error: replacement_policy invalid\n")); 169 ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 170 } else if (PetscIsInfOrNanReal(*f)) { 171 PetscCall(PetscInfo(ls,"OWArmijo line search error: initial function inf or nan\n")); 172 ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 173 } 174 175 if (ls->reason != TAOLINESEARCH_CONTINUE_ITERATING) PetscFunctionReturn(0); 176 177 /* Check to see of the memory has been allocated. If not, allocate 178 the historical array and populate it with the initial function 179 values. */ 180 if (!armP->memory) { 181 PetscCall(PetscMalloc1(armP->memorySize, &armP->memory)); 182 } 183 184 if (!armP->memorySetup) { 185 for (i = 0; i < armP->memorySize; i++) { 186 armP->memory[i] = armP->alpha*(*f); 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 215 if (armP->nondescending) { 216 fact = armP->sigma; 217 } 218 219 PetscCall(VecDuplicate(g,&g_old)); 220 PetscCall(VecCopy(g,g_old)); 221 222 ls->step = ls->initstep; 223 while (ls->step >= owlqn_minstep && ls->nfeval < ls->max_funcs) { 224 /* Calculate iterate */ 225 ++its; 226 PetscCall(VecWAXPY(armP->work,ls->step,s,x)); 227 228 partgdx=0.0; 229 PetscCall(ProjWork_OWLQN(armP->work,x,g_old,&partgdx)); 230 PetscCall(MPIU_Allreduce(&partgdx,&gdx,1,MPIU_REAL,MPIU_SUM,comm)); 231 232 /* Check the condition of gdx */ 233 if (PetscIsInfOrNanReal(gdx)) { 234 PetscCall(PetscInfo(ls,"Initial Line Search step * g is Inf or Nan (%g)\n",(double)gdx)); 235 ls->reason=TAOLINESEARCH_FAILED_INFORNAN; 236 PetscFunctionReturn(0); 237 } 238 if (gdx >= 0.0) { 239 PetscCall(PetscInfo(ls,"Initial Line Search step is not descent direction (g's=%g)\n",(double)gdx)); 240 ls->reason = TAOLINESEARCH_FAILED_ASCENT; 241 PetscFunctionReturn(0); 242 } 243 244 /* Calculate function at new iterate */ 245 PetscCall(TaoLineSearchComputeObjectiveAndGradient(ls,armP->work,f,g)); 246 g_computed=PETSC_TRUE; 247 248 PetscCall(TaoLineSearchMonitor(ls, its, *f, ls->step)); 249 250 if (ls->step == ls->initstep) { 251 ls->f_fullstep = *f; 252 } 253 254 if (PetscIsInfOrNanReal(*f)) { 255 ls->step *= armP->beta_inf; 256 } else { 257 /* Check descent condition */ 258 if (armP->nondescending && *f <= ref - ls->step*fact*ref) break; 259 if (!armP->nondescending && *f <= ref + armP->sigma * gdx) break; 260 ls->step *= armP->beta; 261 } 262 } 263 PetscCall(VecDestroy(&g_old)); 264 265 /* Check termination */ 266 if (PetscIsInfOrNanReal(*f)) { 267 PetscCall(PetscInfo(ls, "Function is inf or nan.\n")); 268 ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER; 269 } else if (ls->step < owlqn_minstep) { 270 PetscCall(PetscInfo(ls, "Step length is below tolerance.\n")); 271 ls->reason = TAOLINESEARCH_HALTED_RTOL; 272 } else if (ls->nfeval >= ls->max_funcs) { 273 PetscCall(PetscInfo(ls, "Number of line search function evals (%" PetscInt_FMT ") > maximum allowed (%" PetscInt_FMT ")\n",ls->nfeval, ls->max_funcs)); 274 ls->reason = TAOLINESEARCH_HALTED_MAXFCN; 275 } 276 if (ls->reason) PetscFunctionReturn(0); 277 278 /* Successful termination, update memory */ 279 ls->reason = TAOLINESEARCH_SUCCESS; 280 armP->lastReference = ref; 281 if (armP->replacementPolicy == REPLACE_FIFO) { 282 armP->memory[armP->current++] = *f; 283 if (armP->current >= armP->memorySize) { 284 armP->current = 0; 285 } 286 } else { 287 armP->current = idx; 288 armP->memory[idx] = *f; 289 } 290 291 /* Update iterate and compute gradient */ 292 PetscCall(VecCopy(armP->work,x)); 293 if (!g_computed) { 294 PetscCall(TaoLineSearchComputeGradient(ls, x, g)); 295 } 296 PetscCall(PetscInfo(ls, "%" PetscInt_FMT " function evals in line search, step = %10.4f\n",ls->nfeval, (double)ls->step)); 297 PetscFunctionReturn(0); 298 } 299 300 /*MC 301 TAOLINESEARCHOWARMIJO - Special line-search type for the Orthant-Wise Limited Quasi-Newton (TAOOWLQN) algorithm. 302 Should not be used with any other algorithm. 303 304 Level: developer 305 306 .keywords: Tao, linesearch 307 M*/ 308 PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_OWArmijo(TaoLineSearch ls) 309 { 310 TaoLineSearch_OWARMIJO *armP; 311 312 PetscFunctionBegin; 313 PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 314 PetscCall(PetscNewLog(ls,&armP)); 315 316 armP->memory = NULL; 317 armP->alpha = 1.0; 318 armP->beta = 0.25; 319 armP->beta_inf = 0.25; 320 armP->sigma = 1e-4; 321 armP->memorySize = 1; 322 armP->referencePolicy = REFERENCE_MAX; 323 armP->replacementPolicy = REPLACE_MRU; 324 armP->nondescending=PETSC_FALSE; 325 ls->data = (void*)armP; 326 ls->initstep = 0.1; 327 ls->ops->monitor = NULL; 328 ls->ops->setup = NULL; 329 ls->ops->reset = NULL; 330 ls->ops->apply = TaoLineSearchApply_OWArmijo; 331 ls->ops->view = TaoLineSearchView_OWArmijo; 332 ls->ops->destroy = TaoLineSearchDestroy_OWArmijo; 333 ls->ops->setfromoptions = TaoLineSearchSetFromOptions_OWArmijo; 334 PetscFunctionReturn(0); 335 } 336