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