1 #include <petsctaolinesearch.h> 2 #include <../src/tao/matrix/lmvmmat.h> 3 #include <../src/tao/unconstrained/impls/owlqn/owlqn.h> 4 5 #define OWLQN_BFGS 0 6 #define OWLQN_SCALED_GRADIENT 1 7 #define OWLQN_GRADIENT 2 8 9 static PetscErrorCode ProjDirect_OWLQN(Vec d, Vec g) 10 { 11 PetscErrorCode ierr; 12 const PetscReal *gptr; 13 PetscReal *dptr; 14 PetscInt low,high,low1,high1,i; 15 16 PetscFunctionBegin; 17 ierr=VecGetOwnershipRange(d,&low,&high);CHKERRQ(ierr); 18 ierr=VecGetOwnershipRange(g,&low1,&high1);CHKERRQ(ierr); 19 20 ierr = VecGetArrayRead(g,&gptr);CHKERRQ(ierr); 21 ierr = VecGetArray(d,&dptr);CHKERRQ(ierr); 22 for (i = 0; i < high-low; i++) { 23 if (dptr[i] * gptr[i] <= 0.0 ) { 24 dptr[i] = 0.0; 25 } 26 } 27 ierr = VecRestoreArray(d,&dptr);CHKERRQ(ierr); 28 ierr = VecRestoreArrayRead(g,&gptr);CHKERRQ(ierr); 29 PetscFunctionReturn(0); 30 } 31 32 static PetscErrorCode ComputePseudoGrad_OWLQN(Vec x, Vec gv, PetscReal lambda) 33 { 34 PetscErrorCode ierr; 35 const PetscReal *xptr; 36 PetscReal *gptr; 37 PetscInt low,high,low1,high1,i; 38 39 PetscFunctionBegin; 40 ierr=VecGetOwnershipRange(x,&low,&high);CHKERRQ(ierr); 41 ierr=VecGetOwnershipRange(gv,&low1,&high1);CHKERRQ(ierr); 42 43 ierr = VecGetArrayRead(x,&xptr);CHKERRQ(ierr); 44 ierr = VecGetArray(gv,&gptr);CHKERRQ(ierr); 45 for (i = 0; i < high-low; i++) { 46 if (xptr[i] < 0.0) gptr[i] = gptr[i] - lambda; 47 else if (xptr[i] > 0.0) gptr[i] = gptr[i] + lambda; 48 else if (gptr[i] + lambda < 0.0) gptr[i] = gptr[i] + lambda; 49 else if (gptr[i] - lambda > 0.0) gptr[i] = gptr[i] - lambda; 50 else gptr[i] = 0.0; 51 } 52 ierr = VecRestoreArray(gv,&gptr);CHKERRQ(ierr); 53 ierr = VecRestoreArrayRead(x,&xptr);CHKERRQ(ierr); 54 PetscFunctionReturn(0); 55 } 56 57 static PetscErrorCode TaoSolve_OWLQN(Tao tao) 58 { 59 TAO_OWLQN *lmP = (TAO_OWLQN *)tao->data; 60 PetscReal f, fold, gdx, gnorm; 61 PetscReal step = 1.0; 62 PetscReal delta; 63 PetscErrorCode ierr; 64 PetscInt stepType; 65 PetscInt iter = 0; 66 TaoConvergedReason reason = TAO_CONTINUE_ITERATING; 67 TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 68 69 PetscFunctionBegin; 70 if (tao->XL || tao->XU || tao->ops->computebounds) { 71 ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by owlqn algorithm\n");CHKERRQ(ierr); 72 } 73 74 /* Check convergence criteria */ 75 ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr); 76 77 ierr = VecCopy(tao->gradient, lmP->GV);CHKERRQ(ierr); 78 79 ierr = ComputePseudoGrad_OWLQN(tao->solution,lmP->GV,lmP->lambda);CHKERRQ(ierr); 80 81 ierr = VecNorm(lmP->GV,NORM_2,&gnorm);CHKERRQ(ierr); 82 83 if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 84 85 ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, step, &reason);CHKERRQ(ierr); 86 if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 87 88 /* Set initial scaling for the function */ 89 if (f != 0.0) { 90 delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 91 } else { 92 delta = 2.0 / (gnorm*gnorm); 93 } 94 ierr = MatLMVMSetDelta(lmP->M,delta);CHKERRQ(ierr); 95 96 /* Set counter for gradient/reset steps */ 97 lmP->bfgs = 0; 98 lmP->sgrad = 0; 99 lmP->grad = 0; 100 101 /* Have not converged; continue with Newton method */ 102 while (reason == TAO_CONTINUE_ITERATING) { 103 /* Compute direction */ 104 ierr = MatLMVMUpdate(lmP->M,tao->solution,tao->gradient);CHKERRQ(ierr); 105 ierr = MatLMVMSolve(lmP->M, lmP->GV, lmP->D);CHKERRQ(ierr); 106 107 ierr = ProjDirect_OWLQN(lmP->D,lmP->GV);CHKERRQ(ierr); 108 109 ++lmP->bfgs; 110 111 /* Check for success (descent direction) */ 112 ierr = VecDot(lmP->D, lmP->GV , &gdx);CHKERRQ(ierr); 113 if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) { 114 115 /* Step is not descent or direction produced not a number 116 We can assert bfgsUpdates > 1 in this case because 117 the first solve produces the scaled gradient direction, 118 which is guaranteed to be descent 119 120 Use steepest descent direction (scaled) */ 121 ++lmP->grad; 122 123 if (f != 0.0) { 124 delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 125 } else { 126 delta = 2.0 / (gnorm*gnorm); 127 } 128 ierr = MatLMVMSetDelta(lmP->M, delta);CHKERRQ(ierr); 129 ierr = MatLMVMReset(lmP->M);CHKERRQ(ierr); 130 ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 131 ierr = MatLMVMSolve(lmP->M,lmP->GV, lmP->D);CHKERRQ(ierr); 132 133 ierr = ProjDirect_OWLQN(lmP->D,lmP->GV);CHKERRQ(ierr); 134 135 lmP->bfgs = 1; 136 ++lmP->sgrad; 137 stepType = OWLQN_SCALED_GRADIENT; 138 } else { 139 if (1 == lmP->bfgs) { 140 /* The first BFGS direction is always the scaled gradient */ 141 ++lmP->sgrad; 142 stepType = OWLQN_SCALED_GRADIENT; 143 } else { 144 ++lmP->bfgs; 145 stepType = OWLQN_BFGS; 146 } 147 } 148 149 ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr); 150 151 /* Perform the linesearch */ 152 fold = f; 153 ierr = VecCopy(tao->solution, lmP->Xold);CHKERRQ(ierr); 154 ierr = VecCopy(tao->gradient, lmP->Gold);CHKERRQ(ierr); 155 156 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, lmP->GV, lmP->D, &step,&ls_status);CHKERRQ(ierr); 157 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 158 159 while (((int)ls_status < 0) && (stepType != OWLQN_GRADIENT)) { 160 161 /* Reset factors and use scaled gradient step */ 162 f = fold; 163 ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr); 164 ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr); 165 ierr = VecCopy(tao->gradient, lmP->GV);CHKERRQ(ierr); 166 167 ierr = ComputePseudoGrad_OWLQN(tao->solution,lmP->GV,lmP->lambda);CHKERRQ(ierr); 168 169 switch(stepType) { 170 case OWLQN_BFGS: 171 /* Failed to obtain acceptable iterate with BFGS step 172 Attempt to use the scaled gradient direction */ 173 174 if (f != 0.0) { 175 delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 176 } else { 177 delta = 2.0 / (gnorm*gnorm); 178 } 179 ierr = MatLMVMSetDelta(lmP->M, delta);CHKERRQ(ierr); 180 ierr = MatLMVMReset(lmP->M);CHKERRQ(ierr); 181 ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 182 ierr = MatLMVMSolve(lmP->M, lmP->GV, lmP->D);CHKERRQ(ierr); 183 184 ierr = ProjDirect_OWLQN(lmP->D,lmP->GV);CHKERRQ(ierr); 185 186 lmP->bfgs = 1; 187 ++lmP->sgrad; 188 stepType = OWLQN_SCALED_GRADIENT; 189 break; 190 191 case OWLQN_SCALED_GRADIENT: 192 /* The scaled gradient step did not produce a new iterate; 193 attempt to use the gradient direction. 194 Need to make sure we are not using a different diagonal scaling */ 195 ierr = MatLMVMSetDelta(lmP->M, 1.0);CHKERRQ(ierr); 196 ierr = MatLMVMReset(lmP->M);CHKERRQ(ierr); 197 ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 198 ierr = MatLMVMSolve(lmP->M, lmP->GV, lmP->D);CHKERRQ(ierr); 199 200 ierr = ProjDirect_OWLQN(lmP->D,lmP->GV);CHKERRQ(ierr); 201 202 lmP->bfgs = 1; 203 ++lmP->grad; 204 stepType = OWLQN_GRADIENT; 205 break; 206 } 207 ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr); 208 209 210 /* Perform the linesearch */ 211 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, lmP->GV, lmP->D, &step, &ls_status);CHKERRQ(ierr); 212 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 213 } 214 215 if ((int)ls_status < 0) { 216 /* Failed to find an improving point*/ 217 f = fold; 218 ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr); 219 ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr); 220 ierr = VecCopy(tao->gradient, lmP->GV);CHKERRQ(ierr); 221 step = 0.0; 222 } else { 223 /* a little hack here, because that gv is used to store g */ 224 ierr = VecCopy(lmP->GV, tao->gradient);CHKERRQ(ierr); 225 } 226 227 ierr = ComputePseudoGrad_OWLQN(tao->solution,lmP->GV,lmP->lambda);CHKERRQ(ierr); 228 229 /* Check for termination */ 230 231 ierr = VecNorm(lmP->GV,NORM_2,&gnorm);CHKERRQ(ierr); 232 233 iter++; 234 ierr = TaoMonitor(tao,iter,f,gnorm,0.0,step,&reason);CHKERRQ(ierr); 235 236 if ((int)ls_status < 0) break; 237 } 238 PetscFunctionReturn(0); 239 } 240 241 static PetscErrorCode TaoSetUp_OWLQN(Tao tao) 242 { 243 TAO_OWLQN *lmP = (TAO_OWLQN *)tao->data; 244 PetscInt n,N; 245 PetscErrorCode ierr; 246 247 PetscFunctionBegin; 248 /* Existence of tao->solution checked in TaoSetUp() */ 249 if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); } 250 if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr); } 251 if (!lmP->D) {ierr = VecDuplicate(tao->solution,&lmP->D);CHKERRQ(ierr); } 252 if (!lmP->GV) {ierr = VecDuplicate(tao->solution,&lmP->GV);CHKERRQ(ierr); } 253 if (!lmP->Xold) {ierr = VecDuplicate(tao->solution,&lmP->Xold);CHKERRQ(ierr); } 254 if (!lmP->Gold) {ierr = VecDuplicate(tao->solution,&lmP->Gold);CHKERRQ(ierr); } 255 256 /* Create matrix for the limited memory approximation */ 257 ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 258 ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 259 ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&lmP->M);CHKERRQ(ierr); 260 ierr = MatLMVMAllocateVectors(lmP->M,tao->solution);CHKERRQ(ierr); 261 PetscFunctionReturn(0); 262 } 263 264 /* ---------------------------------------------------------- */ 265 static PetscErrorCode TaoDestroy_OWLQN(Tao tao) 266 { 267 TAO_OWLQN *lmP = (TAO_OWLQN *)tao->data; 268 PetscErrorCode ierr; 269 270 PetscFunctionBegin; 271 if (tao->setupcalled) { 272 ierr = VecDestroy(&lmP->Xold);CHKERRQ(ierr); 273 ierr = VecDestroy(&lmP->Gold);CHKERRQ(ierr); 274 ierr = VecDestroy(&lmP->D);CHKERRQ(ierr); 275 ierr = MatDestroy(&lmP->M);CHKERRQ(ierr); 276 ierr = VecDestroy(&lmP->GV);CHKERRQ(ierr); 277 } 278 ierr = PetscFree(tao->data);CHKERRQ(ierr); 279 PetscFunctionReturn(0); 280 } 281 282 /*------------------------------------------------------------*/ 283 static PetscErrorCode TaoSetFromOptions_OWLQN(PetscOptionItems *PetscOptionsObject,Tao tao) 284 { 285 TAO_OWLQN *lmP = (TAO_OWLQN *)tao->data; 286 PetscErrorCode ierr; 287 288 PetscFunctionBegin; 289 ierr = PetscOptionsHead(PetscOptionsObject,"Orthant-Wise Limited-memory method for Quasi-Newton unconstrained optimization");CHKERRQ(ierr); 290 ierr = PetscOptionsReal("-tao_owlqn_lambda", "regulariser weight","", 100,&lmP->lambda,NULL); CHKERRQ(ierr); 291 ierr = PetscOptionsTail();CHKERRQ(ierr); 292 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 293 PetscFunctionReturn(0); 294 } 295 296 /*------------------------------------------------------------*/ 297 static PetscErrorCode TaoView_OWLQN(Tao tao, PetscViewer viewer) 298 { 299 TAO_OWLQN *lm = (TAO_OWLQN *)tao->data; 300 PetscBool isascii; 301 PetscErrorCode ierr; 302 303 PetscFunctionBegin; 304 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 305 if (isascii) { 306 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 307 ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", lm->bfgs);CHKERRQ(ierr); 308 ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", lm->sgrad);CHKERRQ(ierr); 309 ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lm->grad);CHKERRQ(ierr); 310 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 311 } 312 PetscFunctionReturn(0); 313 } 314 315 /* ---------------------------------------------------------- */ 316 /*MC 317 TAOOWLQN - orthant-wise limited memory quasi-newton algorithm 318 319 . - tao_owlqn_lambda - regulariser weight 320 321 Level: beginner 322 M*/ 323 324 325 PETSC_EXTERN PetscErrorCode TaoCreate_OWLQN(Tao tao) 326 { 327 TAO_OWLQN *lmP; 328 const char *owarmijo_type = TAOLINESEARCHOWARMIJO; 329 PetscErrorCode ierr; 330 331 PetscFunctionBegin; 332 tao->ops->setup = TaoSetUp_OWLQN; 333 tao->ops->solve = TaoSolve_OWLQN; 334 tao->ops->view = TaoView_OWLQN; 335 tao->ops->setfromoptions = TaoSetFromOptions_OWLQN; 336 tao->ops->destroy = TaoDestroy_OWLQN; 337 338 ierr = PetscNewLog(tao,&lmP);CHKERRQ(ierr); 339 lmP->D = 0; 340 lmP->M = 0; 341 lmP->GV = 0; 342 lmP->Xold = 0; 343 lmP->Gold = 0; 344 lmP->lambda = 1.0; 345 346 tao->data = (void*)lmP; 347 /* Override default settings (unless already changed) */ 348 if (!tao->max_it_changed) tao->max_it = 2000; 349 if (!tao->max_funcs_changed) tao->max_funcs = 4000; 350 351 ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 352 ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 353 ierr = TaoLineSearchSetType(tao->linesearch,owarmijo_type);CHKERRQ(ierr); 354 ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 355 ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 356 PetscFunctionReturn(0); 357 } 358