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