1f273b952SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/taolinesearchimpl.h> 3aaa7dc30SBarry Smith #include <../src/tao/linesearch/impls/owarmijo/owarmijo.h> 4a7e14dcfSSatish Balay 5a7e14dcfSSatish Balay #define REPLACE_FIFO 1 6a7e14dcfSSatish Balay #define REPLACE_MRU 2 7a7e14dcfSSatish Balay 8a7e14dcfSSatish Balay #define REFERENCE_MAX 1 9a7e14dcfSSatish Balay #define REFERENCE_AVE 2 10a7e14dcfSSatish Balay #define REFERENCE_MEAN 3 11a7e14dcfSSatish Balay 12a7e14dcfSSatish Balay static PetscErrorCode ProjWork_OWLQN(Vec w,Vec x,Vec gv,PetscReal *gdx) 13a7e14dcfSSatish Balay { 145e081366SBarry Smith const PetscReal *xptr,*gptr; 155e081366SBarry Smith PetscReal *wptr; 16a7e14dcfSSatish Balay PetscInt low,high,low1,high1,low2,high2,i; 17a7e14dcfSSatish Balay 18a7e14dcfSSatish Balay PetscFunctionBegin; 199566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(w,&low,&high)); 209566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x,&low1,&high1)); 219566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(gv,&low2,&high2)); 22a7e14dcfSSatish Balay 23a7e14dcfSSatish Balay *gdx=0.0; 249566063dSJacob Faibussowitsch PetscCall(VecGetArray(w,&wptr)); 259566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xptr)); 269566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gv,&gptr)); 27a7e14dcfSSatish Balay 2853506e15SBarry Smith for (i=0;i<high-low;i++) { 295e081366SBarry Smith if (xptr[i]*wptr[i]<0.0) wptr[i]=0.0; 30a7e14dcfSSatish Balay *gdx = *gdx + gptr[i]*(wptr[i]-xptr[i]); 31a7e14dcfSSatish Balay } 329566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(w,&wptr)); 339566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xptr)); 349566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gv,&gptr)); 35a7e14dcfSSatish Balay PetscFunctionReturn(0); 36a7e14dcfSSatish Balay } 37a7e14dcfSSatish Balay 38a7e14dcfSSatish Balay static PetscErrorCode TaoLineSearchDestroy_OWArmijo(TaoLineSearch ls) 39a7e14dcfSSatish Balay { 408caf6e8cSBarry Smith TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data; 41a7e14dcfSSatish Balay 42a7e14dcfSSatish Balay PetscFunctionBegin; 439566063dSJacob Faibussowitsch PetscCall(PetscFree(armP->memory)); 44a7e14dcfSSatish Balay if (armP->x) { 459566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)armP->x)); 46a7e14dcfSSatish Balay } 479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&armP->work)); 489566063dSJacob Faibussowitsch PetscCall(PetscFree(ls->data)); 49a7e14dcfSSatish Balay PetscFunctionReturn(0); 50a7e14dcfSSatish Balay } 51a7e14dcfSSatish Balay 524416b707SBarry Smith static PetscErrorCode TaoLineSearchSetFromOptions_OWArmijo(PetscOptionItems *PetscOptionsObject,TaoLineSearch ls) 53a7e14dcfSSatish Balay { 548caf6e8cSBarry Smith TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data; 55a7e14dcfSSatish Balay 56a7e14dcfSSatish Balay PetscFunctionBegin; 57d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"OWArmijo linesearch options"); 589566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_OWArmijo_alpha", "initial reference constant", "", armP->alpha, &armP->alpha,NULL)); 599566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_OWArmijo_beta_inf", "decrease constant one", "", armP->beta_inf, &armP->beta_inf,NULL)); 609566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_OWArmijo_beta", "decrease constant", "", armP->beta, &armP->beta,NULL)); 619566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_OWArmijo_sigma", "acceptance constant", "", armP->sigma, &armP->sigma,NULL)); 629566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-tao_ls_OWArmijo_memory_size", "number of historical elements", "", armP->memorySize, &armP->memorySize,NULL)); 639566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-tao_ls_OWArmijo_reference_policy", "policy for updating reference value", "", armP->referencePolicy, &armP->referencePolicy,NULL)); 649566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-tao_ls_OWArmijo_replacement_policy", "policy for updating memory", "", armP->replacementPolicy, &armP->replacementPolicy,NULL)); 659566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_ls_OWArmijo_nondescending","Use nondescending OWArmijo algorithm","",armP->nondescending,&armP->nondescending,NULL)); 66d0609cedSBarry Smith PetscOptionsHeadEnd(); 67a7e14dcfSSatish Balay PetscFunctionReturn(0); 68a7e14dcfSSatish Balay } 69a7e14dcfSSatish Balay 70a7e14dcfSSatish Balay static PetscErrorCode TaoLineSearchView_OWArmijo(TaoLineSearch ls, PetscViewer pv) 71a7e14dcfSSatish Balay { 728caf6e8cSBarry Smith TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data; 73a7e14dcfSSatish Balay PetscBool isascii; 74a7e14dcfSSatish Balay 75a7e14dcfSSatish Balay PetscFunctionBegin; 769566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &isascii)); 77a7e14dcfSSatish Balay if (isascii) { 789566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv," OWArmijo linesearch",armP->alpha)); 79a7e14dcfSSatish Balay if (armP->nondescending) { 809566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, " (nondescending)")); 81a7e14dcfSSatish Balay } 829566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,": alpha=%g beta=%g ",(double)armP->alpha,(double)armP->beta)); 839566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"sigma=%g ",(double)armP->sigma)); 849566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv,"memsize=%D\n",armP->memorySize)); 85a7e14dcfSSatish Balay } 86a7e14dcfSSatish Balay PetscFunctionReturn(0); 87a7e14dcfSSatish Balay } 88a7e14dcfSSatish Balay 89a7e14dcfSSatish Balay /* @ TaoApply_OWArmijo - This routine performs a linesearch. It 90a7e14dcfSSatish Balay backtracks until the (nonmonotone) OWArmijo conditions are satisfied. 91a7e14dcfSSatish Balay 92a7e14dcfSSatish Balay Input Parameters: 93a7e14dcfSSatish Balay + tao - TAO_SOLVER context 94a7e14dcfSSatish Balay . X - current iterate (on output X contains new iterate, X + step*S) 95a7e14dcfSSatish Balay . S - search direction 96a7e14dcfSSatish Balay . f - merit function evaluated at X 97a7e14dcfSSatish Balay . G - gradient of merit function evaluated at X 98a7e14dcfSSatish Balay . W - work vector 99a7e14dcfSSatish Balay - step - initial estimate of step length 100a7e14dcfSSatish Balay 101a7e14dcfSSatish Balay Output parameters: 102a7e14dcfSSatish Balay + f - merit function evaluated at new iterate, X + step*S 103a7e14dcfSSatish Balay . G - gradient of merit function evaluated at new iterate, X + step*S 104a7e14dcfSSatish Balay . X - new iterate 105a7e14dcfSSatish Balay - step - final step length 106a7e14dcfSSatish Balay 107a7e14dcfSSatish Balay Info is set to one of: 108a7e14dcfSSatish Balay . 0 - the line search succeeds; the sufficient decrease 109a7e14dcfSSatish Balay condition and the directional derivative condition hold 110a7e14dcfSSatish Balay 111a7e14dcfSSatish Balay negative number if an input parameter is invalid 112a7e14dcfSSatish Balay - -1 - step < 0 113a7e14dcfSSatish Balay 114a7e14dcfSSatish Balay positive number > 1 if the line search otherwise terminates 115a7e14dcfSSatish Balay + 1 - Step is at the lower bound, stepmin. 116a7e14dcfSSatish Balay @ */ 117a7e14dcfSSatish Balay static PetscErrorCode TaoLineSearchApply_OWArmijo(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s) 118a7e14dcfSSatish Balay { 1198caf6e8cSBarry Smith TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data; 1202a0dac07SAlp Dener PetscInt i, its=0; 121a7e14dcfSSatish Balay PetscReal fact, ref, gdx; 122a7e14dcfSSatish Balay PetscInt idx; 123a7e14dcfSSatish Balay PetscBool g_computed=PETSC_FALSE; /* to prevent extra gradient computation */ 124a7e14dcfSSatish Balay Vec g_old; 125a7e14dcfSSatish Balay PetscReal owlqn_minstep=0.005; 126a7e14dcfSSatish Balay PetscReal partgdx; 127ba4b436cSBarry Smith MPI_Comm comm; 128a7e14dcfSSatish Balay 129a7e14dcfSSatish Balay PetscFunctionBegin; 1309566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ls,&comm)); 131a7e14dcfSSatish Balay fact = 0.0; 132a7e14dcfSSatish Balay ls->nfeval=0; 133a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_CONTINUE_ITERATING; 134a7e14dcfSSatish Balay if (!armP->work) { 1359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&armP->work)); 136a7e14dcfSSatish Balay armP->x = x; 1379566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)armP->x)); 13853506e15SBarry Smith } else if (x != armP->x) { 1399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&armP->work)); 1409566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&armP->work)); 1419566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)armP->x)); 142a7e14dcfSSatish Balay armP->x = x; 1439566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)armP->x)); 144a7e14dcfSSatish Balay } 145a7e14dcfSSatish Balay 1469566063dSJacob Faibussowitsch PetscCall(TaoLineSearchMonitor(ls, 0, *f, 0.0)); 1472a0dac07SAlp Dener 148a7e14dcfSSatish Balay /* Check linesearch parameters */ 149a7e14dcfSSatish Balay if (armP->alpha < 1) { 1509566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls,"OWArmijo line search error: alpha (%g) < 1\n", (double)armP->alpha)); 151a7e14dcfSSatish Balay ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 15253506e15SBarry Smith } else if ((armP->beta <= 0) || (armP->beta >= 1)) { 1539566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls,"OWArmijo line search error: beta (%g) invalid\n", (double)armP->beta)); 154a7e14dcfSSatish Balay ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 15553506e15SBarry Smith } else if ((armP->beta_inf <= 0) || (armP->beta_inf >= 1)) { 1569566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls,"OWArmijo line search error: beta_inf (%g) invalid\n", (double)armP->beta_inf)); 157a7e14dcfSSatish Balay ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 15853506e15SBarry Smith } else if ((armP->sigma <= 0) || (armP->sigma >= 0.5)) { 1599566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls,"OWArmijo line search error: sigma (%g) invalid\n", (double)armP->sigma)); 160a7e14dcfSSatish Balay ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 16153506e15SBarry Smith } else if (armP->memorySize < 1) { 1629566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls,"OWArmijo line search error: memory_size (%D) < 1\n", armP->memorySize)); 163a7e14dcfSSatish Balay ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 16453506e15SBarry Smith } else if ((armP->referencePolicy != REFERENCE_MAX) && (armP->referencePolicy != REFERENCE_AVE) && (armP->referencePolicy != REFERENCE_MEAN)) { 1659566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls,"OWArmijo line search error: reference_policy invalid\n")); 166a7e14dcfSSatish Balay ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 16753506e15SBarry Smith } else if ((armP->replacementPolicy != REPLACE_FIFO) && (armP->replacementPolicy != REPLACE_MRU)) { 1689566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls,"OWArmijo line search error: replacement_policy invalid\n")); 169a7e14dcfSSatish Balay ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 17053506e15SBarry Smith } else if (PetscIsInfOrNanReal(*f)) { 1719566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls,"OWArmijo line search error: initial function inf or nan\n")); 172a7e14dcfSSatish Balay ls->reason=TAOLINESEARCH_FAILED_BADPARAMETER; 173a7e14dcfSSatish Balay } 174a7e14dcfSSatish Balay 17553506e15SBarry Smith if (ls->reason != TAOLINESEARCH_CONTINUE_ITERATING) PetscFunctionReturn(0); 176a7e14dcfSSatish Balay 177a7e14dcfSSatish Balay /* Check to see of the memory has been allocated. If not, allocate 178a7e14dcfSSatish Balay the historical array and populate it with the initial function 179a7e14dcfSSatish Balay values. */ 1806c23d075SBarry Smith if (!armP->memory) { 1819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(armP->memorySize, &armP->memory)); 182a7e14dcfSSatish Balay } 183a7e14dcfSSatish Balay 184a7e14dcfSSatish Balay if (!armP->memorySetup) { 185a7e14dcfSSatish Balay for (i = 0; i < armP->memorySize; i++) { 186a7e14dcfSSatish Balay armP->memory[i] = armP->alpha*(*f); 187a7e14dcfSSatish Balay } 188a7e14dcfSSatish Balay armP->current = 0; 189a7e14dcfSSatish Balay armP->lastReference = armP->memory[0]; 190a7e14dcfSSatish Balay armP->memorySetup=PETSC_TRUE; 191a7e14dcfSSatish Balay } 192a7e14dcfSSatish Balay 193a7e14dcfSSatish Balay /* Calculate reference value (MAX) */ 194a7e14dcfSSatish Balay ref = armP->memory[0]; 195a7e14dcfSSatish Balay idx = 0; 196a7e14dcfSSatish Balay 197a7e14dcfSSatish Balay for (i = 1; i < armP->memorySize; i++) { 198a7e14dcfSSatish Balay if (armP->memory[i] > ref) { 199a7e14dcfSSatish Balay ref = armP->memory[i]; 200a7e14dcfSSatish Balay idx = i; 201a7e14dcfSSatish Balay } 202a7e14dcfSSatish Balay } 203a7e14dcfSSatish Balay 204a7e14dcfSSatish Balay if (armP->referencePolicy == REFERENCE_AVE) { 205a7e14dcfSSatish Balay ref = 0; 206a7e14dcfSSatish Balay for (i = 0; i < armP->memorySize; i++) { 207a7e14dcfSSatish Balay ref += armP->memory[i]; 208a7e14dcfSSatish Balay } 209a7e14dcfSSatish Balay ref = ref / armP->memorySize; 210a7e14dcfSSatish Balay ref = PetscMax(ref, armP->memory[armP->current]); 21153506e15SBarry Smith } else if (armP->referencePolicy == REFERENCE_MEAN) { 212a7e14dcfSSatish Balay ref = PetscMin(ref, 0.5*(armP->lastReference + armP->memory[armP->current])); 213a7e14dcfSSatish Balay } 214a7e14dcfSSatish Balay 215a7e14dcfSSatish Balay if (armP->nondescending) { 216a7e14dcfSSatish Balay fact = armP->sigma; 217a7e14dcfSSatish Balay } 218a7e14dcfSSatish Balay 2199566063dSJacob Faibussowitsch PetscCall(VecDuplicate(g,&g_old)); 2209566063dSJacob Faibussowitsch PetscCall(VecCopy(g,g_old)); 221a7e14dcfSSatish Balay 222a7e14dcfSSatish Balay ls->step = ls->initstep; 223a7e14dcfSSatish Balay while (ls->step >= owlqn_minstep && ls->nfeval < ls->max_funcs) { 224a7e14dcfSSatish Balay /* Calculate iterate */ 2252a0dac07SAlp Dener ++its; 226*ef46b1a6SStefano Zampini PetscCall(VecWAXPY(armP->work,ls->step,s,x)); 227a7e14dcfSSatish Balay 228a7e14dcfSSatish Balay partgdx=0.0; 2299566063dSJacob Faibussowitsch PetscCall(ProjWork_OWLQN(armP->work,x,g_old,&partgdx)); 2301c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&partgdx,&gdx,1,MPIU_REAL,MPIU_SUM,comm)); 231a7e14dcfSSatish Balay 232a7e14dcfSSatish Balay /* Check the condition of gdx */ 233a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(gdx)) { 2349566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls,"Initial Line Search step * g is Inf or Nan (%g)\n",(double)gdx)); 235a7e14dcfSSatish Balay ls->reason=TAOLINESEARCH_FAILED_INFORNAN; 236a7e14dcfSSatish Balay PetscFunctionReturn(0); 237a7e14dcfSSatish Balay } 238a7e14dcfSSatish Balay if (gdx >= 0.0) { 2399566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls,"Initial Line Search step is not descent direction (g's=%g)\n",(double)gdx)); 240a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_ASCENT; 241a7e14dcfSSatish Balay PetscFunctionReturn(0); 242a7e14dcfSSatish Balay } 243a7e14dcfSSatish Balay 244a7e14dcfSSatish Balay /* Calculate function at new iterate */ 2459566063dSJacob Faibussowitsch PetscCall(TaoLineSearchComputeObjectiveAndGradient(ls,armP->work,f,g)); 246a7e14dcfSSatish Balay g_computed=PETSC_TRUE; 247a7e14dcfSSatish Balay 2489566063dSJacob Faibussowitsch PetscCall(TaoLineSearchMonitor(ls, its, *f, ls->step)); 2492a0dac07SAlp Dener 250a7e14dcfSSatish Balay if (ls->step == ls->initstep) { 251a7e14dcfSSatish Balay ls->f_fullstep = *f; 252a7e14dcfSSatish Balay } 253a7e14dcfSSatish Balay 254a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(*f)) { 255a7e14dcfSSatish Balay ls->step *= armP->beta_inf; 25653506e15SBarry Smith } else { 257a7e14dcfSSatish Balay /* Check descent condition */ 25853506e15SBarry Smith if (armP->nondescending && *f <= ref - ls->step*fact*ref) break; 25953506e15SBarry Smith if (!armP->nondescending && *f <= ref + armP->sigma * gdx) break; 260a7e14dcfSSatish Balay ls->step *= armP->beta; 261a7e14dcfSSatish Balay } 262a7e14dcfSSatish Balay } 2639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&g_old)); 264a7e14dcfSSatish Balay 265a7e14dcfSSatish Balay /* Check termination */ 266a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(*f)) { 2679566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Function is inf or nan.\n")); 268a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER; 269a7e14dcfSSatish Balay } else if (ls->step < owlqn_minstep) { 2709566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Step length is below tolerance.\n")); 271a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_RTOL; 272a7e14dcfSSatish Balay } else if (ls->nfeval >= ls->max_funcs) { 2739566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Number of line search function evals (%D) > maximum allowed (%D)\n",ls->nfeval, ls->max_funcs)); 274a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_MAXFCN; 275a7e14dcfSSatish Balay } 27653506e15SBarry Smith if (ls->reason) PetscFunctionReturn(0); 277a7e14dcfSSatish Balay 278a7e14dcfSSatish Balay /* Successful termination, update memory */ 279a31dc708SJason Sarich ls->reason = TAOLINESEARCH_SUCCESS; 280a7e14dcfSSatish Balay armP->lastReference = ref; 281a7e14dcfSSatish Balay if (armP->replacementPolicy == REPLACE_FIFO) { 282a7e14dcfSSatish Balay armP->memory[armP->current++] = *f; 283a7e14dcfSSatish Balay if (armP->current >= armP->memorySize) { 284a7e14dcfSSatish Balay armP->current = 0; 285a7e14dcfSSatish Balay } 286a7e14dcfSSatish Balay } else { 287a7e14dcfSSatish Balay armP->current = idx; 288a7e14dcfSSatish Balay armP->memory[idx] = *f; 289a7e14dcfSSatish Balay } 290a7e14dcfSSatish Balay 291a7e14dcfSSatish Balay /* Update iterate and compute gradient */ 2929566063dSJacob Faibussowitsch PetscCall(VecCopy(armP->work,x)); 293a7e14dcfSSatish Balay if (!g_computed) { 2949566063dSJacob Faibussowitsch PetscCall(TaoLineSearchComputeGradient(ls, x, g)); 295a7e14dcfSSatish Balay } 2969566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "%D function evals in line search, step = %10.4f\n",ls->nfeval, (double)ls->step)); 297a7e14dcfSSatish Balay PetscFunctionReturn(0); 298a7e14dcfSSatish Balay } 299a7e14dcfSSatish Balay 30090b6438dSAlp Dener /*MC 30190b6438dSAlp Dener TAOLINESEARCHOWARMIJO - Special line-search type for the Orthant-Wise Limited Quasi-Newton (TAOOWLQN) algorithm. 30290b6438dSAlp Dener Should not be used with any other algorithm. 30390b6438dSAlp Dener 30490b6438dSAlp Dener Level: developer 30590b6438dSAlp Dener 30690b6438dSAlp Dener .keywords: Tao, linesearch 30790b6438dSAlp Dener M*/ 308728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_OWArmijo(TaoLineSearch ls) 309a7e14dcfSSatish Balay { 3108caf6e8cSBarry Smith TaoLineSearch_OWARMIJO *armP; 311a7e14dcfSSatish Balay 312a7e14dcfSSatish Balay PetscFunctionBegin; 313a7e14dcfSSatish Balay PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1); 3149566063dSJacob Faibussowitsch PetscCall(PetscNewLog(ls,&armP)); 315a7e14dcfSSatish Balay 3166c23d075SBarry Smith armP->memory = NULL; 317a7e14dcfSSatish Balay armP->alpha = 1.0; 318a7e14dcfSSatish Balay armP->beta = 0.25; 319a7e14dcfSSatish Balay armP->beta_inf = 0.25; 320a7e14dcfSSatish Balay armP->sigma = 1e-4; 321a7e14dcfSSatish Balay armP->memorySize = 1; 322a7e14dcfSSatish Balay armP->referencePolicy = REFERENCE_MAX; 323a7e14dcfSSatish Balay armP->replacementPolicy = REPLACE_MRU; 324a7e14dcfSSatish Balay armP->nondescending=PETSC_FALSE; 325a7e14dcfSSatish Balay ls->data = (void*)armP; 326a7e14dcfSSatish Balay ls->initstep = 0.1; 3274664e3ffSStefano Zampini ls->ops->monitor = NULL; 32883c8fe1dSLisandro Dalcin ls->ops->setup = NULL; 32983c8fe1dSLisandro Dalcin ls->ops->reset = NULL; 330a7e14dcfSSatish Balay ls->ops->apply = TaoLineSearchApply_OWArmijo; 331a7e14dcfSSatish Balay ls->ops->view = TaoLineSearchView_OWArmijo; 332a7e14dcfSSatish Balay ls->ops->destroy = TaoLineSearchDestroy_OWArmijo; 333a7e14dcfSSatish Balay ls->ops->setfromoptions = TaoLineSearchSetFromOptions_OWArmijo; 334a7e14dcfSSatish Balay PetscFunctionReturn(0); 335a7e14dcfSSatish Balay } 336