1af0996ceSBarry Smith #include <petsc/private/taolinesearchimpl.h> 2aaa7dc30SBarry Smith #include <../src/tao/linesearch/impls/armijo/armijo.h> 3a7e14dcfSSatish Balay 4a7e14dcfSSatish Balay #define REPLACE_FIFO 1 5a7e14dcfSSatish Balay #define REPLACE_MRU 2 6a7e14dcfSSatish Balay 7a7e14dcfSSatish Balay #define REFERENCE_MAX 1 8a7e14dcfSSatish Balay #define REFERENCE_AVE 2 9a7e14dcfSSatish Balay #define REFERENCE_MEAN 3 10a7e14dcfSSatish Balay 11d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchDestroy_Armijo(TaoLineSearch ls) 12d71ae5a4SJacob Faibussowitsch { 138caf6e8cSBarry Smith TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data; 14a7e14dcfSSatish Balay 15a7e14dcfSSatish Balay PetscFunctionBegin; 169566063dSJacob Faibussowitsch PetscCall(PetscFree(armP->memory)); 179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&armP->x)); 189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&armP->work)); 199566063dSJacob Faibussowitsch PetscCall(PetscFree(ls->data)); 203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 21a7e14dcfSSatish Balay } 22a7e14dcfSSatish Balay 23d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchReset_Armijo(TaoLineSearch ls) 24d71ae5a4SJacob Faibussowitsch { 258caf6e8cSBarry Smith TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data; 26a7e14dcfSSatish Balay 27a7e14dcfSSatish Balay PetscFunctionBegin; 289566063dSJacob Faibussowitsch PetscCall(PetscFree(armP->memory)); 29a7e14dcfSSatish Balay armP->memorySetup = PETSC_FALSE; 303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 31a7e14dcfSSatish Balay } 32a7e14dcfSSatish Balay 33ce78bad3SBarry Smith static PetscErrorCode TaoLineSearchSetFromOptions_Armijo(TaoLineSearch ls, PetscOptionItems PetscOptionsObject) 34d71ae5a4SJacob Faibussowitsch { 358caf6e8cSBarry Smith TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data; 36a7e14dcfSSatish Balay 37a7e14dcfSSatish Balay PetscFunctionBegin; 38d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Armijo linesearch options"); 399566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_armijo_alpha", "initial reference constant", "", armP->alpha, &armP->alpha, NULL)); 409566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_armijo_beta_inf", "decrease constant one", "", armP->beta_inf, &armP->beta_inf, NULL)); 419566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_armijo_beta", "decrease constant", "", armP->beta, &armP->beta, NULL)); 429566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_armijo_sigma", "acceptance constant", "", armP->sigma, &armP->sigma, NULL)); 439566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-tao_ls_armijo_memory_size", "number of historical elements", "", armP->memorySize, &armP->memorySize, NULL)); 449566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-tao_ls_armijo_reference_policy", "policy for updating reference value", "", armP->referencePolicy, &armP->referencePolicy, NULL)); 459566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-tao_ls_armijo_replacement_policy", "policy for updating memory", "", armP->replacementPolicy, &armP->replacementPolicy, NULL)); 469566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_ls_armijo_nondescending", "Use nondescending armijo algorithm", "", armP->nondescending, &armP->nondescending, NULL)); 47d0609cedSBarry Smith PetscOptionsHeadEnd(); 483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49a7e14dcfSSatish Balay } 50a7e14dcfSSatish Balay 51d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchView_Armijo(TaoLineSearch ls, PetscViewer pv) 52d71ae5a4SJacob Faibussowitsch { 538caf6e8cSBarry Smith TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data; 54a7e14dcfSSatish Balay PetscBool isascii; 55a7e14dcfSSatish Balay 56a7e14dcfSSatish Balay PetscFunctionBegin; 579566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &isascii)); 58a7e14dcfSSatish Balay if (isascii) { 5963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, " Armijo linesearch")); 6048a46eb9SPierre Jolivet if (armP->nondescending) PetscCall(PetscViewerASCIIPrintf(pv, " (nondescending)")); 6148a46eb9SPierre Jolivet if (ls->bounded) PetscCall(PetscViewerASCIIPrintf(pv, " (projected)")); 629566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, ": alpha=%g beta=%g ", (double)armP->alpha, (double)armP->beta)); 639566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "sigma=%g ", (double)armP->sigma)); 6463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "memsize=%" PetscInt_FMT "\n", armP->memorySize)); 65a7e14dcfSSatish Balay } 663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 67a7e14dcfSSatish Balay } 68a7e14dcfSSatish Balay 69a7e14dcfSSatish Balay /* @ TaoApply_Armijo - This routine performs a linesearch. It 70a7e14dcfSSatish Balay backtracks until the (nonmonotone) Armijo conditions are satisfied. 71a7e14dcfSSatish Balay 72a7e14dcfSSatish Balay Input Parameters: 73441846f8SBarry Smith + tao - Tao context 74a7e14dcfSSatish Balay . X - current iterate (on output X contains new iterate, X + step*S) 75a7e14dcfSSatish Balay . S - search direction 76a7e14dcfSSatish Balay . f - merit function evaluated at X 77a7e14dcfSSatish Balay . G - gradient of merit function evaluated at X 78a7e14dcfSSatish Balay . W - work vector 79a7e14dcfSSatish Balay - step - initial estimate of step length 80a7e14dcfSSatish Balay 81a7e14dcfSSatish Balay Output parameters: 82a7e14dcfSSatish Balay + f - merit function evaluated at new iterate, X + step*S 83a7e14dcfSSatish Balay . G - gradient of merit function evaluated at new iterate, X + step*S 84a7e14dcfSSatish Balay . X - new iterate 85a7e14dcfSSatish Balay - step - final step length 86a7e14dcfSSatish Balay 87a7e14dcfSSatish Balay @ */ 88d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchApply_Armijo(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s) 89d71ae5a4SJacob Faibussowitsch { 908caf6e8cSBarry Smith TaoLineSearch_ARMIJO *armP = (TaoLineSearch_ARMIJO *)ls->data; 912a0dac07SAlp Dener PetscInt i, its = 0; 92a7e14dcfSSatish Balay PetscReal fact, ref, gdx; 93a7e14dcfSSatish Balay PetscInt idx; 94a7e14dcfSSatish Balay PetscBool g_computed = PETSC_FALSE; /* to prevent extra gradient computation */ 95a7e14dcfSSatish Balay 96a7e14dcfSSatish Balay PetscFunctionBegin; 979566063dSJacob Faibussowitsch PetscCall(TaoLineSearchMonitor(ls, 0, *f, 0.0)); 982a0dac07SAlp Dener 99a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_CONTINUE_ITERATING; 100a7e14dcfSSatish Balay if (!armP->work) { 1019566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &armP->work)); 102a7e14dcfSSatish Balay armP->x = x; 1039566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)armP->x)); 10447a47007SBarry Smith } else if (x != armP->x) { 105a7e14dcfSSatish Balay /* If x has changed, then recreate work */ 1069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&armP->work)); 1079566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &armP->work)); 1089566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)armP->x)); 109a7e14dcfSSatish Balay armP->x = x; 1109566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)armP->x)); 111a7e14dcfSSatish Balay } 112a7e14dcfSSatish Balay 113a7e14dcfSSatish Balay /* Check linesearch parameters */ 114a7e14dcfSSatish Balay if (armP->alpha < 1) { 1159566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Armijo line search error: alpha (%g) < 1\n", (double)armP->alpha)); 116a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER; 11747a47007SBarry Smith } else if ((armP->beta <= 0) || (armP->beta >= 1)) { 1189566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Armijo line search error: beta (%g) invalid\n", (double)armP->beta)); 119a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER; 12047a47007SBarry Smith } else if ((armP->beta_inf <= 0) || (armP->beta_inf >= 1)) { 1219566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Armijo line search error: beta_inf (%g) invalid\n", (double)armP->beta_inf)); 122a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER; 12347a47007SBarry Smith } else if ((armP->sigma <= 0) || (armP->sigma >= 0.5)) { 1249566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Armijo line search error: sigma (%g) invalid\n", (double)armP->sigma)); 125a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER; 12647a47007SBarry Smith } else if (armP->memorySize < 1) { 12763a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Armijo line search error: memory_size (%" PetscInt_FMT ") < 1\n", armP->memorySize)); 128a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER; 12947a47007SBarry Smith } else if ((armP->referencePolicy != REFERENCE_MAX) && (armP->referencePolicy != REFERENCE_AVE) && (armP->referencePolicy != REFERENCE_MEAN)) { 1309566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Armijo line search error: reference_policy invalid\n")); 131a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER; 13247a47007SBarry Smith } else if ((armP->replacementPolicy != REPLACE_FIFO) && (armP->replacementPolicy != REPLACE_MRU)) { 1339566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Armijo line search error: replacement_policy invalid\n")); 134a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER; 13547a47007SBarry Smith } else if (PetscIsInfOrNanReal(*f)) { 1369566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Armijo line search error: initial function inf or nan\n")); 137a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER; 138a7e14dcfSSatish Balay } 139a7e14dcfSSatish Balay 1403ba16761SJacob Faibussowitsch if (ls->reason != TAOLINESEARCH_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS); 141a7e14dcfSSatish Balay 142a7e14dcfSSatish Balay /* Check to see of the memory has been allocated. If not, allocate 143a7e14dcfSSatish Balay the historical array and populate it with the initial function 144a7e14dcfSSatish Balay values. */ 14548a46eb9SPierre Jolivet if (!armP->memory) PetscCall(PetscMalloc1(armP->memorySize, &armP->memory)); 146a7e14dcfSSatish Balay 147a7e14dcfSSatish Balay if (!armP->memorySetup) { 148ad540459SPierre Jolivet for (i = 0; i < armP->memorySize; i++) armP->memory[i] = armP->alpha * (*f); 149a7e14dcfSSatish Balay 150a7e14dcfSSatish Balay armP->current = 0; 151a7e14dcfSSatish Balay armP->lastReference = armP->memory[0]; 152a7e14dcfSSatish Balay armP->memorySetup = PETSC_TRUE; 153a7e14dcfSSatish Balay } 154a7e14dcfSSatish Balay 155a7e14dcfSSatish Balay /* Calculate reference value (MAX) */ 1562d631c23SStefano Zampini ref = *f; 1572d631c23SStefano Zampini idx = armP->current; 1582d631c23SStefano Zampini for (i = 0; i < armP->memorySize; i++) { 159a7e14dcfSSatish Balay if (armP->memory[i] > ref) { 160a7e14dcfSSatish Balay ref = armP->memory[i]; 161a7e14dcfSSatish Balay idx = i; 162a7e14dcfSSatish Balay } 163a7e14dcfSSatish Balay } 164a7e14dcfSSatish Balay 165a7e14dcfSSatish Balay if (armP->referencePolicy == REFERENCE_AVE) { 166a7e14dcfSSatish Balay ref = 0; 167ad540459SPierre Jolivet for (i = 0; i < armP->memorySize; i++) ref += armP->memory[i]; 168a7e14dcfSSatish Balay ref = ref / armP->memorySize; 169a7e14dcfSSatish Balay ref = PetscMax(ref, armP->memory[armP->current]); 17047a47007SBarry Smith } else if (armP->referencePolicy == REFERENCE_MEAN) { 171a7e14dcfSSatish Balay ref = PetscMin(ref, 0.5 * (armP->lastReference + armP->memory[armP->current])); 172a7e14dcfSSatish Balay } 1739566063dSJacob Faibussowitsch PetscCall(VecDot(g, s, &gdx)); 174a7e14dcfSSatish Balay 175a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(gdx)) { 176*76c63389SBarry Smith PetscCall(PetscInfo(ls, "Initial Line Search step * g is infinity or NaN (%g)\n", (double)gdx)); 177a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_INFORNAN; 1783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 179a7e14dcfSSatish Balay } 180a7e14dcfSSatish Balay if (gdx >= 0.0) { 1819566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Initial Line Search step is not descent direction (g's=%g)\n", (double)gdx)); 182a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_ASCENT; 1833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 184a7e14dcfSSatish Balay } 185a7e14dcfSSatish Balay 186a7e14dcfSSatish Balay if (armP->nondescending) { 187a7e14dcfSSatish Balay fact = armP->sigma; 188a7e14dcfSSatish Balay } else { 189a7e14dcfSSatish Balay fact = armP->sigma * gdx; 190a7e14dcfSSatish Balay } 191a7e14dcfSSatish Balay ls->step = ls->initstep; 192a7e14dcfSSatish Balay while (ls->step >= ls->stepmin && (ls->nfeval + ls->nfgeval) < ls->max_funcs) { 193a7e14dcfSSatish Balay /* Calculate iterate */ 1942a0dac07SAlp Dener ++its; 195ef46b1a6SStefano Zampini PetscCall(VecWAXPY(armP->work, ls->step, s, x)); 1961baa6e33SBarry Smith if (ls->bounded) PetscCall(VecMedian(ls->lower, armP->work, ls->upper, armP->work)); 197a7e14dcfSSatish Balay 198a7e14dcfSSatish Balay /* Calculate function at new iterate */ 199e6994092SAlex Lindsay PetscCall(VecLockReadPush(x)); 200a7e14dcfSSatish Balay if (ls->hasobjective) { 2019566063dSJacob Faibussowitsch PetscCall(TaoLineSearchComputeObjective(ls, armP->work, f)); 202a7e14dcfSSatish Balay g_computed = PETSC_FALSE; 203a7e14dcfSSatish Balay } else if (ls->usegts) { 2049566063dSJacob Faibussowitsch PetscCall(TaoLineSearchComputeObjectiveAndGTS(ls, armP->work, f, &gdx)); 205a7e14dcfSSatish Balay g_computed = PETSC_FALSE; 206a7e14dcfSSatish Balay } else { 2079566063dSJacob Faibussowitsch PetscCall(TaoLineSearchComputeObjectiveAndGradient(ls, armP->work, f, g)); 208a7e14dcfSSatish Balay g_computed = PETSC_TRUE; 209a7e14dcfSSatish Balay } 210e6994092SAlex Lindsay PetscCall(VecLockReadPop(x)); 211ad540459SPierre Jolivet if (ls->step == ls->initstep) ls->f_fullstep = *f; 212a7e14dcfSSatish Balay 2139566063dSJacob Faibussowitsch PetscCall(TaoLineSearchMonitor(ls, its, *f, ls->step)); 2142a0dac07SAlp Dener 215a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(*f)) { 216a7e14dcfSSatish Balay ls->step *= armP->beta_inf; 21747a47007SBarry Smith } else { 218a7e14dcfSSatish Balay /* Check descent condition */ 2199371c9d4SSatish Balay if (armP->nondescending && *f <= ref - ls->step * fact * ref) break; 220ad540459SPierre Jolivet if (!armP->nondescending && *f <= ref + ls->step * fact) break; 221a7e14dcfSSatish Balay 222a7e14dcfSSatish Balay ls->step *= armP->beta; 223a7e14dcfSSatish Balay } 224a7e14dcfSSatish Balay } 225a7e14dcfSSatish Balay 226a7e14dcfSSatish Balay /* Check termination */ 227a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(*f)) { 2289566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Function is inf or nan.\n")); 229a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_INFORNAN; 230a7e14dcfSSatish Balay } else if (ls->step < ls->stepmin) { 2319566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Step length is below tolerance.\n")); 232a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_RTOL; 233a7e14dcfSSatish Balay } else if ((ls->nfeval + ls->nfgeval) >= ls->max_funcs) { 23463a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Number of line search function evals (%" PetscInt_FMT ") > maximum allowed (%" PetscInt_FMT ")\n", ls->nfeval + ls->nfgeval, ls->max_funcs)); 235a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_MAXFCN; 236a7e14dcfSSatish Balay } 2373ba16761SJacob Faibussowitsch if (ls->reason) PetscFunctionReturn(PETSC_SUCCESS); 238a7e14dcfSSatish Balay 239a7e14dcfSSatish Balay /* Successful termination, update memory */ 240a31dc708SJason Sarich ls->reason = TAOLINESEARCH_SUCCESS; 241a7e14dcfSSatish Balay armP->lastReference = ref; 242a7e14dcfSSatish Balay if (armP->replacementPolicy == REPLACE_FIFO) { 243a7e14dcfSSatish Balay armP->memory[armP->current++] = *f; 244ad540459SPierre Jolivet if (armP->current >= armP->memorySize) armP->current = 0; 245a7e14dcfSSatish Balay } else { 246a7e14dcfSSatish Balay armP->current = idx; 247a7e14dcfSSatish Balay armP->memory[idx] = *f; 248a7e14dcfSSatish Balay } 249a7e14dcfSSatish Balay 250a7e14dcfSSatish Balay /* Update iterate and compute gradient */ 2519566063dSJacob Faibussowitsch PetscCall(VecCopy(armP->work, x)); 25248a46eb9SPierre Jolivet if (!g_computed) PetscCall(TaoLineSearchComputeGradient(ls, x, g)); 25363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ls, "%" PetscInt_FMT " function evals in line search, step = %g\n", ls->nfeval, (double)ls->step)); 2543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 255a7e14dcfSSatish Balay } 256a7e14dcfSSatish Balay 25790b6438dSAlp Dener /*MC 25890b6438dSAlp Dener TAOLINESEARCHARMIJO - Backtracking line-search that satisfies only the (nonmonotone) Armijo condition 25990b6438dSAlp Dener (i.e., sufficient decrease). 26090b6438dSAlp Dener 26190b6438dSAlp Dener Armijo line-search type can be selected with "-tao_ls_type armijo". 26290b6438dSAlp Dener 26390b6438dSAlp Dener Level: developer 26490b6438dSAlp Dener 265db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetType()`, `TaoLineSearchApply()` 26690b6438dSAlp Dener M*/ 267d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_Armijo(TaoLineSearch ls) 268d71ae5a4SJacob Faibussowitsch { 2698caf6e8cSBarry Smith TaoLineSearch_ARMIJO *armP; 270a7e14dcfSSatish Balay 271a7e14dcfSSatish Balay PetscFunctionBegin; 272a7e14dcfSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 2734dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&armP)); 274a7e14dcfSSatish Balay 2756c23d075SBarry Smith armP->memory = NULL; 276a7e14dcfSSatish Balay armP->alpha = 1.0; 277a7e14dcfSSatish Balay armP->beta = 0.5; 278a7e14dcfSSatish Balay armP->beta_inf = 0.5; 279a7e14dcfSSatish Balay armP->sigma = 1e-4; 280a7e14dcfSSatish Balay armP->memorySize = 1; 281a7e14dcfSSatish Balay armP->referencePolicy = REFERENCE_MAX; 282a7e14dcfSSatish Balay armP->replacementPolicy = REPLACE_MRU; 283a7e14dcfSSatish Balay armP->nondescending = PETSC_FALSE; 284a7e14dcfSSatish Balay ls->data = (void *)armP; 285a7e14dcfSSatish Balay ls->initstep = 1.0; 28683c8fe1dSLisandro Dalcin ls->ops->setup = NULL; 2874664e3ffSStefano Zampini ls->ops->monitor = NULL; 288a7e14dcfSSatish Balay ls->ops->apply = TaoLineSearchApply_Armijo; 289a7e14dcfSSatish Balay ls->ops->view = TaoLineSearchView_Armijo; 290a7e14dcfSSatish Balay ls->ops->destroy = TaoLineSearchDestroy_Armijo; 291a7e14dcfSSatish Balay ls->ops->reset = TaoLineSearchReset_Armijo; 292a7e14dcfSSatish Balay ls->ops->setfromoptions = TaoLineSearchSetFromOptions_Armijo; 2933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 294a7e14dcfSSatish Balay } 295