1af0996ceSBarry Smith #include <petsc/private/taolinesearchimpl.h> 2aaa7dc30SBarry Smith #include <../src/tao/linesearch/impls/owarmijo/owarmijo.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 ProjWork_OWLQN(Vec w, Vec x, Vec gv, PetscReal *gdx) 12d71ae5a4SJacob Faibussowitsch { 135e081366SBarry Smith const PetscReal *xptr, *gptr; 145e081366SBarry Smith PetscReal *wptr; 15a7e14dcfSSatish Balay PetscInt low, high, low1, high1, low2, high2, i; 16a7e14dcfSSatish Balay 17a7e14dcfSSatish Balay PetscFunctionBegin; 189566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(w, &low, &high)); 199566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x, &low1, &high1)); 209566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(gv, &low2, &high2)); 21a7e14dcfSSatish Balay 22a7e14dcfSSatish Balay *gdx = 0.0; 239566063dSJacob Faibussowitsch PetscCall(VecGetArray(w, &wptr)); 249566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xptr)); 259566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gv, &gptr)); 26a7e14dcfSSatish Balay 2753506e15SBarry Smith for (i = 0; i < high - low; i++) { 285e081366SBarry Smith if (xptr[i] * wptr[i] < 0.0) wptr[i] = 0.0; 29a7e14dcfSSatish Balay *gdx = *gdx + gptr[i] * (wptr[i] - xptr[i]); 30a7e14dcfSSatish Balay } 319566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(w, &wptr)); 329566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xptr)); 339566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gv, &gptr)); 343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 35a7e14dcfSSatish Balay } 36a7e14dcfSSatish Balay 37d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchDestroy_OWArmijo(TaoLineSearch ls) 38d71ae5a4SJacob Faibussowitsch { 398caf6e8cSBarry Smith TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data; 40a7e14dcfSSatish Balay 41a7e14dcfSSatish Balay PetscFunctionBegin; 429566063dSJacob Faibussowitsch PetscCall(PetscFree(armP->memory)); 431baa6e33SBarry Smith if (armP->x) PetscCall(PetscObjectDereference((PetscObject)armP->x)); 449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&armP->work)); 459566063dSJacob Faibussowitsch PetscCall(PetscFree(ls->data)); 463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 47a7e14dcfSSatish Balay } 48a7e14dcfSSatish Balay 49ce78bad3SBarry Smith static PetscErrorCode TaoLineSearchSetFromOptions_OWArmijo(TaoLineSearch ls, PetscOptionItems PetscOptionsObject) 50d71ae5a4SJacob Faibussowitsch { 518caf6e8cSBarry Smith TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data; 52a7e14dcfSSatish Balay 53a7e14dcfSSatish Balay PetscFunctionBegin; 54d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "OWArmijo linesearch options"); 559566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_OWArmijo_alpha", "initial reference constant", "", armP->alpha, &armP->alpha, NULL)); 569566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_OWArmijo_beta_inf", "decrease constant one", "", armP->beta_inf, &armP->beta_inf, NULL)); 579566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_OWArmijo_beta", "decrease constant", "", armP->beta, &armP->beta, NULL)); 589566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ls_OWArmijo_sigma", "acceptance constant", "", armP->sigma, &armP->sigma, NULL)); 599566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-tao_ls_OWArmijo_memory_size", "number of historical elements", "", armP->memorySize, &armP->memorySize, NULL)); 609566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-tao_ls_OWArmijo_reference_policy", "policy for updating reference value", "", armP->referencePolicy, &armP->referencePolicy, NULL)); 619566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-tao_ls_OWArmijo_replacement_policy", "policy for updating memory", "", armP->replacementPolicy, &armP->replacementPolicy, NULL)); 629566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_ls_OWArmijo_nondescending", "Use nondescending OWArmijo algorithm", "", armP->nondescending, &armP->nondescending, NULL)); 63d0609cedSBarry Smith PetscOptionsHeadEnd(); 643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 65a7e14dcfSSatish Balay } 66a7e14dcfSSatish Balay 67d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchView_OWArmijo(TaoLineSearch ls, PetscViewer pv) 68d71ae5a4SJacob Faibussowitsch { 698caf6e8cSBarry Smith TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data; 70a7e14dcfSSatish Balay PetscBool isascii; 71a7e14dcfSSatish Balay 72a7e14dcfSSatish Balay PetscFunctionBegin; 739566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &isascii)); 74a7e14dcfSSatish Balay if (isascii) { 7563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, " OWArmijo linesearch")); 7648a46eb9SPierre Jolivet if (armP->nondescending) PetscCall(PetscViewerASCIIPrintf(pv, " (nondescending)")); 779566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, ": alpha=%g beta=%g ", (double)armP->alpha, (double)armP->beta)); 789566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "sigma=%g ", (double)armP->sigma)); 7963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "memsize=%" PetscInt_FMT "\n", armP->memorySize)); 80a7e14dcfSSatish Balay } 813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 82a7e14dcfSSatish Balay } 83a7e14dcfSSatish Balay 84a7e14dcfSSatish Balay /* @ TaoApply_OWArmijo - This routine performs a linesearch. It 85a7e14dcfSSatish Balay backtracks until the (nonmonotone) OWArmijo conditions are satisfied. 86a7e14dcfSSatish Balay 87a7e14dcfSSatish Balay Input Parameters: 88a7e14dcfSSatish Balay + tao - TAO_SOLVER context 89a7e14dcfSSatish Balay . X - current iterate (on output X contains new iterate, X + step*S) 90a7e14dcfSSatish Balay . S - search direction 91a7e14dcfSSatish Balay . f - merit function evaluated at X 92a7e14dcfSSatish Balay . G - gradient of merit function evaluated at X 93a7e14dcfSSatish Balay . W - work vector 94a7e14dcfSSatish Balay - step - initial estimate of step length 95a7e14dcfSSatish Balay 96a7e14dcfSSatish Balay Output parameters: 97a7e14dcfSSatish Balay + f - merit function evaluated at new iterate, X + step*S 98a7e14dcfSSatish Balay . G - gradient of merit function evaluated at new iterate, X + step*S 99a7e14dcfSSatish Balay . X - new iterate 100a7e14dcfSSatish Balay - step - final step length 101a7e14dcfSSatish Balay 102a7e14dcfSSatish Balay Info is set to one of: 103a7e14dcfSSatish Balay . 0 - the line search succeeds; the sufficient decrease 104a7e14dcfSSatish Balay condition and the directional derivative condition hold 105a7e14dcfSSatish Balay 106a7e14dcfSSatish Balay negative number if an input parameter is invalid 107a7e14dcfSSatish Balay - -1 - step < 0 108a7e14dcfSSatish Balay 109a7e14dcfSSatish Balay positive number > 1 if the line search otherwise terminates 110a7e14dcfSSatish Balay + 1 - Step is at the lower bound, stepmin. 111a7e14dcfSSatish Balay @ */ 112d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchApply_OWArmijo(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s) 113d71ae5a4SJacob Faibussowitsch { 1148caf6e8cSBarry Smith TaoLineSearch_OWARMIJO *armP = (TaoLineSearch_OWARMIJO *)ls->data; 1152a0dac07SAlp Dener PetscInt i, its = 0; 116a7e14dcfSSatish Balay PetscReal fact, ref, gdx; 117a7e14dcfSSatish Balay PetscInt idx; 118a7e14dcfSSatish Balay PetscBool g_computed = PETSC_FALSE; /* to prevent extra gradient computation */ 119a7e14dcfSSatish Balay Vec g_old; 120a7e14dcfSSatish Balay PetscReal owlqn_minstep = 0.005; 121a7e14dcfSSatish Balay PetscReal partgdx; 122ba4b436cSBarry Smith MPI_Comm comm; 123a7e14dcfSSatish Balay 124a7e14dcfSSatish Balay PetscFunctionBegin; 1259566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ls, &comm)); 126a7e14dcfSSatish Balay fact = 0.0; 127a7e14dcfSSatish Balay ls->nfeval = 0; 128a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_CONTINUE_ITERATING; 129a7e14dcfSSatish Balay if (!armP->work) { 1309566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &armP->work)); 131a7e14dcfSSatish Balay armP->x = x; 1329566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)armP->x)); 13353506e15SBarry Smith } else if (x != armP->x) { 1349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&armP->work)); 1359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &armP->work)); 1369566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)armP->x)); 137a7e14dcfSSatish Balay armP->x = x; 1389566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)armP->x)); 139a7e14dcfSSatish Balay } 140a7e14dcfSSatish Balay 1419566063dSJacob Faibussowitsch PetscCall(TaoLineSearchMonitor(ls, 0, *f, 0.0)); 1422a0dac07SAlp Dener 143a7e14dcfSSatish Balay /* Check linesearch parameters */ 144a7e14dcfSSatish Balay if (armP->alpha < 1) { 1459566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "OWArmijo line search error: alpha (%g) < 1\n", (double)armP->alpha)); 146a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER; 14753506e15SBarry Smith } else if ((armP->beta <= 0) || (armP->beta >= 1)) { 1489566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "OWArmijo line search error: beta (%g) invalid\n", (double)armP->beta)); 149a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER; 15053506e15SBarry Smith } else if ((armP->beta_inf <= 0) || (armP->beta_inf >= 1)) { 1519566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "OWArmijo line search error: beta_inf (%g) invalid\n", (double)armP->beta_inf)); 152a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER; 15353506e15SBarry Smith } else if ((armP->sigma <= 0) || (armP->sigma >= 0.5)) { 1549566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "OWArmijo line search error: sigma (%g) invalid\n", (double)armP->sigma)); 155a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER; 15653506e15SBarry Smith } else if (armP->memorySize < 1) { 15763a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ls, "OWArmijo line search error: memory_size (%" PetscInt_FMT ") < 1\n", armP->memorySize)); 158a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER; 15953506e15SBarry Smith } else if ((armP->referencePolicy != REFERENCE_MAX) && (armP->referencePolicy != REFERENCE_AVE) && (armP->referencePolicy != REFERENCE_MEAN)) { 1609566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "OWArmijo line search error: reference_policy invalid\n")); 161a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER; 16253506e15SBarry Smith } else if ((armP->replacementPolicy != REPLACE_FIFO) && (armP->replacementPolicy != REPLACE_MRU)) { 1639566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "OWArmijo line search error: replacement_policy invalid\n")); 164a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER; 16553506e15SBarry Smith } else if (PetscIsInfOrNanReal(*f)) { 1669566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "OWArmijo line search error: initial function inf or nan\n")); 167a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER; 168a7e14dcfSSatish Balay } 169a7e14dcfSSatish Balay 1703ba16761SJacob Faibussowitsch if (ls->reason != TAOLINESEARCH_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS); 171a7e14dcfSSatish Balay 172a7e14dcfSSatish Balay /* Check to see of the memory has been allocated. If not, allocate 173a7e14dcfSSatish Balay the historical array and populate it with the initial function 174a7e14dcfSSatish Balay values. */ 17548a46eb9SPierre Jolivet if (!armP->memory) PetscCall(PetscMalloc1(armP->memorySize, &armP->memory)); 176a7e14dcfSSatish Balay 177a7e14dcfSSatish Balay if (!armP->memorySetup) { 178ad540459SPierre Jolivet for (i = 0; i < armP->memorySize; i++) armP->memory[i] = armP->alpha * (*f); 179a7e14dcfSSatish Balay armP->current = 0; 180a7e14dcfSSatish Balay armP->lastReference = armP->memory[0]; 181a7e14dcfSSatish Balay armP->memorySetup = PETSC_TRUE; 182a7e14dcfSSatish Balay } 183a7e14dcfSSatish Balay 184a7e14dcfSSatish Balay /* Calculate reference value (MAX) */ 185a7e14dcfSSatish Balay ref = armP->memory[0]; 186a7e14dcfSSatish Balay idx = 0; 187a7e14dcfSSatish Balay 188a7e14dcfSSatish Balay for (i = 1; i < armP->memorySize; i++) { 189a7e14dcfSSatish Balay if (armP->memory[i] > ref) { 190a7e14dcfSSatish Balay ref = armP->memory[i]; 191a7e14dcfSSatish Balay idx = i; 192a7e14dcfSSatish Balay } 193a7e14dcfSSatish Balay } 194a7e14dcfSSatish Balay 195a7e14dcfSSatish Balay if (armP->referencePolicy == REFERENCE_AVE) { 196a7e14dcfSSatish Balay ref = 0; 197ad540459SPierre Jolivet for (i = 0; i < armP->memorySize; i++) ref += armP->memory[i]; 198a7e14dcfSSatish Balay ref = ref / armP->memorySize; 199a7e14dcfSSatish Balay ref = PetscMax(ref, armP->memory[armP->current]); 20053506e15SBarry Smith } else if (armP->referencePolicy == REFERENCE_MEAN) { 201a7e14dcfSSatish Balay ref = PetscMin(ref, 0.5 * (armP->lastReference + armP->memory[armP->current])); 202a7e14dcfSSatish Balay } 203a7e14dcfSSatish Balay 204ad540459SPierre Jolivet if (armP->nondescending) fact = armP->sigma; 205a7e14dcfSSatish Balay 2069566063dSJacob Faibussowitsch PetscCall(VecDuplicate(g, &g_old)); 2079566063dSJacob Faibussowitsch PetscCall(VecCopy(g, g_old)); 208a7e14dcfSSatish Balay 209a7e14dcfSSatish Balay ls->step = ls->initstep; 210a7e14dcfSSatish Balay while (ls->step >= owlqn_minstep && ls->nfeval < ls->max_funcs) { 211a7e14dcfSSatish Balay /* Calculate iterate */ 2122a0dac07SAlp Dener ++its; 213ef46b1a6SStefano Zampini PetscCall(VecWAXPY(armP->work, ls->step, s, x)); 214a7e14dcfSSatish Balay 215a7e14dcfSSatish Balay partgdx = 0.0; 2169566063dSJacob Faibussowitsch PetscCall(ProjWork_OWLQN(armP->work, x, g_old, &partgdx)); 217462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&partgdx, &gdx, 1, MPIU_REAL, MPIU_SUM, comm)); 218a7e14dcfSSatish Balay 219a7e14dcfSSatish Balay /* Check the condition of gdx */ 220a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(gdx)) { 221*76c63389SBarry Smith PetscCall(PetscInfo(ls, "Initial Line Search step * g is infinity or NaN (%g)\n", (double)gdx)); 222a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_INFORNAN; 2233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 224a7e14dcfSSatish Balay } 225a7e14dcfSSatish Balay if (gdx >= 0.0) { 2269566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Initial Line Search step is not descent direction (g's=%g)\n", (double)gdx)); 227a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_ASCENT; 2283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 229a7e14dcfSSatish Balay } 230a7e14dcfSSatish Balay 231a7e14dcfSSatish Balay /* Calculate function at new iterate */ 232e6994092SAlex Lindsay PetscCall(VecLockReadPush(x)); 2339566063dSJacob Faibussowitsch PetscCall(TaoLineSearchComputeObjectiveAndGradient(ls, armP->work, f, g)); 234e6994092SAlex Lindsay PetscCall(VecLockReadPop(x)); 235a7e14dcfSSatish Balay g_computed = PETSC_TRUE; 236a7e14dcfSSatish Balay 2379566063dSJacob Faibussowitsch PetscCall(TaoLineSearchMonitor(ls, its, *f, ls->step)); 2382a0dac07SAlp Dener 239ad540459SPierre Jolivet if (ls->step == ls->initstep) ls->f_fullstep = *f; 240a7e14dcfSSatish Balay 241a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(*f)) { 242a7e14dcfSSatish Balay ls->step *= armP->beta_inf; 24353506e15SBarry Smith } else { 244a7e14dcfSSatish Balay /* Check descent condition */ 24553506e15SBarry Smith if (armP->nondescending && *f <= ref - ls->step * fact * ref) break; 24653506e15SBarry Smith if (!armP->nondescending && *f <= ref + armP->sigma * gdx) break; 247a7e14dcfSSatish Balay ls->step *= armP->beta; 248a7e14dcfSSatish Balay } 249a7e14dcfSSatish Balay } 2509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&g_old)); 251a7e14dcfSSatish Balay 252a7e14dcfSSatish Balay /* Check termination */ 253a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(*f)) { 2549566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Function is inf or nan.\n")); 255a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_BADPARAMETER; 256a7e14dcfSSatish Balay } else if (ls->step < owlqn_minstep) { 2579566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Step length is below tolerance.\n")); 258a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_RTOL; 259a7e14dcfSSatish Balay } else if (ls->nfeval >= ls->max_funcs) { 26063a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Number of line search function evals (%" PetscInt_FMT ") > maximum allowed (%" PetscInt_FMT ")\n", ls->nfeval, ls->max_funcs)); 261a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_MAXFCN; 262a7e14dcfSSatish Balay } 2633ba16761SJacob Faibussowitsch if (ls->reason) PetscFunctionReturn(PETSC_SUCCESS); 264a7e14dcfSSatish Balay 265a7e14dcfSSatish Balay /* Successful termination, update memory */ 266a31dc708SJason Sarich ls->reason = TAOLINESEARCH_SUCCESS; 267a7e14dcfSSatish Balay armP->lastReference = ref; 268a7e14dcfSSatish Balay if (armP->replacementPolicy == REPLACE_FIFO) { 269a7e14dcfSSatish Balay armP->memory[armP->current++] = *f; 270ad540459SPierre Jolivet if (armP->current >= armP->memorySize) armP->current = 0; 271a7e14dcfSSatish Balay } else { 272a7e14dcfSSatish Balay armP->current = idx; 273a7e14dcfSSatish Balay armP->memory[idx] = *f; 274a7e14dcfSSatish Balay } 275a7e14dcfSSatish Balay 276a7e14dcfSSatish Balay /* Update iterate and compute gradient */ 2779566063dSJacob Faibussowitsch PetscCall(VecCopy(armP->work, x)); 27848a46eb9SPierre Jolivet if (!g_computed) PetscCall(TaoLineSearchComputeGradient(ls, x, g)); 27963a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ls, "%" PetscInt_FMT " function evals in line search, step = %10.4f\n", ls->nfeval, (double)ls->step)); 2803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 281a7e14dcfSSatish Balay } 282a7e14dcfSSatish Balay 28390b6438dSAlp Dener /*MC 2841d27aa22SBarry Smith TAOLINESEARCHOWARMIJO - Special line-search type for the Orthant-Wise Limited Quasi-Newton (`TAOOWLQN`) algorithm. 28590b6438dSAlp Dener Should not be used with any other algorithm. 28690b6438dSAlp Dener 28790b6438dSAlp Dener Level: developer 28890b6438dSAlp Dener 289f8d70eaaSPierre Jolivet .seealso: `TaoLineSearch`, `TAOOWLQN`, `Tao` 29090b6438dSAlp Dener M*/ 291d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_OWArmijo(TaoLineSearch ls) 292d71ae5a4SJacob Faibussowitsch { 2938caf6e8cSBarry Smith TaoLineSearch_OWARMIJO *armP; 294a7e14dcfSSatish Balay 295a7e14dcfSSatish Balay PetscFunctionBegin; 296a7e14dcfSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1); 2974dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&armP)); 298a7e14dcfSSatish Balay 2996c23d075SBarry Smith armP->memory = NULL; 300a7e14dcfSSatish Balay armP->alpha = 1.0; 301a7e14dcfSSatish Balay armP->beta = 0.25; 302a7e14dcfSSatish Balay armP->beta_inf = 0.25; 303a7e14dcfSSatish Balay armP->sigma = 1e-4; 304a7e14dcfSSatish Balay armP->memorySize = 1; 305a7e14dcfSSatish Balay armP->referencePolicy = REFERENCE_MAX; 306a7e14dcfSSatish Balay armP->replacementPolicy = REPLACE_MRU; 307a7e14dcfSSatish Balay armP->nondescending = PETSC_FALSE; 308a7e14dcfSSatish Balay ls->data = (void *)armP; 309a7e14dcfSSatish Balay ls->initstep = 0.1; 3104664e3ffSStefano Zampini ls->ops->monitor = NULL; 31183c8fe1dSLisandro Dalcin ls->ops->setup = NULL; 31283c8fe1dSLisandro Dalcin ls->ops->reset = NULL; 313a7e14dcfSSatish Balay ls->ops->apply = TaoLineSearchApply_OWArmijo; 314a7e14dcfSSatish Balay ls->ops->view = TaoLineSearchView_OWArmijo; 315a7e14dcfSSatish Balay ls->ops->destroy = TaoLineSearchDestroy_OWArmijo; 316a7e14dcfSSatish Balay ls->ops->setfromoptions = TaoLineSearchSetFromOptions_OWArmijo; 3173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 318a7e14dcfSSatish Balay } 319