xref: /petsc/src/tao/unconstrained/impls/neldermead/neldermead.c (revision 03047865b8d8757cf1cf9cda45785c1537b01dc1)
1aaa7dc30SBarry Smith #include <../src/tao/unconstrained/impls/neldermead/neldermead.h>
2aaa7dc30SBarry Smith #include <petscvec.h>
3a7e14dcfSSatish Balay 
4d71ae5a4SJacob Faibussowitsch static PetscErrorCode NelderMeadSort(TAO_NelderMead *nm)
5d71ae5a4SJacob Faibussowitsch {
640244768SBarry Smith   PetscReal *values  = nm->f_values;
740244768SBarry Smith   PetscInt  *indices = nm->indices;
840244768SBarry Smith   PetscInt   dim     = nm->N + 1;
940244768SBarry Smith   PetscInt   i, j, index;
1040244768SBarry Smith   PetscReal  val;
1140244768SBarry Smith 
1240244768SBarry Smith   PetscFunctionBegin;
1340244768SBarry Smith   for (i = 1; i < dim; i++) {
1440244768SBarry Smith     index = indices[i];
1540244768SBarry Smith     val   = values[index];
16ad540459SPierre Jolivet     for (j = i - 1; j >= 0 && values[indices[j]] > val; j--) indices[j + 1] = indices[j];
1740244768SBarry Smith     indices[j + 1] = index;
1840244768SBarry Smith   }
193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2040244768SBarry Smith }
2140244768SBarry Smith 
22d71ae5a4SJacob Faibussowitsch static PetscErrorCode NelderMeadReplace(TAO_NelderMead *nm, PetscInt index, Vec Xmu, PetscReal f)
23d71ae5a4SJacob Faibussowitsch {
2440244768SBarry Smith   PetscFunctionBegin;
2540244768SBarry Smith   /*  Add new vector's fraction of average */
269566063dSJacob Faibussowitsch   PetscCall(VecAXPY(nm->Xbar, nm->oneOverN, Xmu));
279566063dSJacob Faibussowitsch   PetscCall(VecCopy(Xmu, nm->simplex[index]));
2840244768SBarry Smith   nm->f_values[index] = f;
2940244768SBarry Smith 
309566063dSJacob Faibussowitsch   PetscCall(NelderMeadSort(nm));
3140244768SBarry Smith 
3240244768SBarry Smith   /*  Subtract last vector from average */
339566063dSJacob Faibussowitsch   PetscCall(VecAXPY(nm->Xbar, -nm->oneOverN, nm->simplex[nm->indices[nm->N]]));
343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3540244768SBarry Smith }
3640244768SBarry Smith 
37d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_NM(Tao tao)
38d71ae5a4SJacob Faibussowitsch {
39a7e14dcfSSatish Balay   TAO_NelderMead *nm = (TAO_NelderMead *)tao->data;
4053506e15SBarry Smith   PetscInt        n;
41a7e14dcfSSatish Balay 
42a7e14dcfSSatish Balay   PetscFunctionBegin;
439566063dSJacob Faibussowitsch   PetscCall(VecGetSize(tao->solution, &n));
4453506e15SBarry Smith   nm->N        = n;
4553506e15SBarry Smith   nm->oneOverN = 1.0 / n;
469566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(tao->solution, nm->N + 1, &nm->simplex));
479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nm->N + 1, &nm->f_values));
489566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nm->N + 1, &nm->indices));
499566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &nm->Xbar));
509566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &nm->Xmur));
519566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &nm->Xmue));
529566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &nm->Xmuc));
53a7e14dcfSSatish Balay 
5483c8fe1dSLisandro Dalcin   tao->gradient = NULL;
55a7e14dcfSSatish Balay   tao->step     = 0;
563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
57a7e14dcfSSatish Balay }
58a7e14dcfSSatish Balay 
59d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_NM(Tao tao)
60d71ae5a4SJacob Faibussowitsch {
61a7e14dcfSSatish Balay   TAO_NelderMead *nm = (TAO_NelderMead *)tao->data;
6253506e15SBarry Smith 
63a7e14dcfSSatish Balay   PetscFunctionBegin;
64a7e14dcfSSatish Balay   if (tao->setupcalled) {
659566063dSJacob Faibussowitsch     PetscCall(VecDestroyVecs(nm->N + 1, &nm->simplex));
669566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&nm->Xmuc));
679566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&nm->Xmue));
689566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&nm->Xmur));
699566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&nm->Xbar));
70a7e14dcfSSatish Balay   }
719566063dSJacob Faibussowitsch   PetscCall(PetscFree(nm->indices));
729566063dSJacob Faibussowitsch   PetscCall(PetscFree(nm->f_values));
739566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
75a7e14dcfSSatish Balay }
76a7e14dcfSSatish Balay 
77ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_NM(Tao tao, PetscOptionItems PetscOptionsObject)
78d71ae5a4SJacob Faibussowitsch {
79a7e14dcfSSatish Balay   TAO_NelderMead *nm = (TAO_NelderMead *)tao->data;
80a7e14dcfSSatish Balay 
81a7e14dcfSSatish Balay   PetscFunctionBegin;
82d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Nelder-Mead options");
83a8d3b578SPierre Jolivet   PetscCall(PetscOptionsDeprecated("-tao_nm_lamda", "-tao_nm_lambda", "3.18.4", NULL));
84a8d3b578SPierre Jolivet   PetscCall(PetscOptionsReal("-tao_nm_lambda", "initial step length", "", nm->lambda, &nm->lambda, NULL));
859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nm_mu", "mu", "", nm->mu_oc, &nm->mu_oc, NULL));
86a7e14dcfSSatish Balay   nm->mu_ic = -nm->mu_oc;
87a7e14dcfSSatish Balay   nm->mu_r  = nm->mu_oc * 2.0;
88a7e14dcfSSatish Balay   nm->mu_e  = nm->mu_oc * 4.0;
89d0609cedSBarry Smith   PetscOptionsHeadEnd();
903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
91a7e14dcfSSatish Balay }
92a7e14dcfSSatish Balay 
93d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_NM(Tao tao, PetscViewer viewer)
94d71ae5a4SJacob Faibussowitsch {
95a7e14dcfSSatish Balay   TAO_NelderMead *nm = (TAO_NelderMead *)tao->data;
96a7e14dcfSSatish Balay   PetscBool       isascii;
97a7e14dcfSSatish Balay 
98a7e14dcfSSatish Balay   PetscFunctionBegin;
999566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
100a7e14dcfSSatish Balay   if (isascii) {
1019566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
10263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "expansions: %" PetscInt_FMT "\n", nm->nexpand));
10363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "reflections: %" PetscInt_FMT "\n", nm->nreflect));
10463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "inside contractions: %" PetscInt_FMT "\n", nm->nincontract));
10563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "outside contractionss: %" PetscInt_FMT "\n", nm->noutcontract));
10663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Shrink steps: %" PetscInt_FMT "\n", nm->nshrink));
1079566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
108a7e14dcfSSatish Balay   }
1093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
110a7e14dcfSSatish Balay }
111a7e14dcfSSatish Balay 
112d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_NM(Tao tao)
113d71ae5a4SJacob Faibussowitsch {
114a7e14dcfSSatish Balay   TAO_NelderMead *nm = (TAO_NelderMead *)tao->data;
115a7e14dcfSSatish Balay   PetscReal      *x;
1168931d482SJason Sarich   PetscInt        i;
117a7e14dcfSSatish Balay   Vec             Xmur = nm->Xmur, Xmue = nm->Xmue, Xmuc = nm->Xmuc, Xbar = nm->Xbar;
118a7e14dcfSSatish Balay   PetscReal       fr, fe, fc;
119a7e14dcfSSatish Balay   PetscInt        shrink;
120a7e14dcfSSatish Balay   PetscInt        low, high;
121a7e14dcfSSatish Balay 
122a7e14dcfSSatish Balay   PetscFunctionBegin;
123a7e14dcfSSatish Balay   nm->nshrink      = 0;
124a7e14dcfSSatish Balay   nm->nreflect     = 0;
125a7e14dcfSSatish Balay   nm->nincontract  = 0;
126a7e14dcfSSatish Balay   nm->noutcontract = 0;
127a7e14dcfSSatish Balay   nm->nexpand      = 0;
128a7e14dcfSSatish Balay 
12948a46eb9SPierre Jolivet   if (tao->XL || tao->XU || tao->ops->computebounds) PetscCall(PetscInfo(tao, "WARNING: Variable bounds have been set but will be ignored by NelderMead algorithm\n"));
130a7e14dcfSSatish Balay 
1319566063dSJacob Faibussowitsch   PetscCall(VecCopy(tao->solution, nm->simplex[0]));
1329566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjective(tao, nm->simplex[0], &nm->f_values[0]));
133a7e14dcfSSatish Balay   nm->indices[0] = 0;
134a7e14dcfSSatish Balay   for (i = 1; i < nm->N + 1; i++) {
1359566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->solution, nm->simplex[i]));
1369566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(nm->simplex[i], &low, &high));
137a7e14dcfSSatish Balay     if (i - 1 >= low && i - 1 < high) {
1389566063dSJacob Faibussowitsch       PetscCall(VecGetArray(nm->simplex[i], &x));
139a8d3b578SPierre Jolivet       x[i - 1 - low] += nm->lambda;
1409566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(nm->simplex[i], &x));
141a7e14dcfSSatish Balay     }
142a7e14dcfSSatish Balay 
1439566063dSJacob Faibussowitsch     PetscCall(TaoComputeObjective(tao, nm->simplex[i], &nm->f_values[i]));
144a7e14dcfSSatish Balay     nm->indices[i] = i;
145a7e14dcfSSatish Balay   }
146a7e14dcfSSatish Balay 
147a7e14dcfSSatish Balay   /*  Xbar  = (Sum of all simplex vectors - worst vector)/N */
1489566063dSJacob Faibussowitsch   PetscCall(NelderMeadSort(nm));
1499566063dSJacob Faibussowitsch   PetscCall(VecSet(Xbar, 0.0));
15048a46eb9SPierre Jolivet   for (i = 0; i < nm->N; i++) PetscCall(VecAXPY(Xbar, 1.0, nm->simplex[nm->indices[i]]));
1519566063dSJacob Faibussowitsch   PetscCall(VecScale(Xbar, nm->oneOverN));
1523ecd9318SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
153a7e14dcfSSatish Balay   while (1) {
154e1e80dc8SAlp Dener     /* Call general purpose update function */
155dbbe0bcdSBarry Smith     PetscTryTypeMethod(tao, update, tao->niter, tao->user_update);
156d643b47cSToby Isaac     ++tao->niter;
157a7e14dcfSSatish Balay     shrink = 0;
1589566063dSJacob Faibussowitsch     PetscCall(VecCopy(nm->simplex[nm->indices[0]], tao->solution));
1599566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, nm->f_values[nm->indices[0]], nm->f_values[nm->indices[nm->N]] - nm->f_values[nm->indices[0]], 0.0, tao->ksp_its));
1609566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, nm->f_values[nm->indices[0]], nm->f_values[nm->indices[nm->N]] - nm->f_values[nm->indices[0]], 0.0, 1.0));
161dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
1623ecd9318SAlp Dener     if (tao->reason != TAO_CONTINUE_ITERATING) break;
163a7e14dcfSSatish Balay 
164a7e14dcfSSatish Balay     /* x(mu) = (1 + mu)Xbar - mu*X_N+1 */
1659566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(Xmur, 1 + nm->mu_r, -nm->mu_r, 0, Xbar, nm->simplex[nm->indices[nm->N]]));
1669566063dSJacob Faibussowitsch     PetscCall(TaoComputeObjective(tao, Xmur, &fr));
167a7e14dcfSSatish Balay 
168a7e14dcfSSatish Balay     if (nm->f_values[nm->indices[0]] <= fr && fr < nm->f_values[nm->indices[nm->N - 1]]) {
169a7e14dcfSSatish Balay       /*  reflect */
170a7e14dcfSSatish Balay       nm->nreflect++;
1719566063dSJacob Faibussowitsch       PetscCall(PetscInfo(0, "Reflect\n"));
1729566063dSJacob Faibussowitsch       PetscCall(NelderMeadReplace(nm, nm->indices[nm->N], Xmur, fr));
17353506e15SBarry Smith     } else if (fr < nm->f_values[nm->indices[0]]) {
174a7e14dcfSSatish Balay       /*  expand */
175a7e14dcfSSatish Balay       nm->nexpand++;
1769566063dSJacob Faibussowitsch       PetscCall(PetscInfo(0, "Expand\n"));
1779566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(Xmue, 1 + nm->mu_e, -nm->mu_e, 0, Xbar, nm->simplex[nm->indices[nm->N]]));
1789566063dSJacob Faibussowitsch       PetscCall(TaoComputeObjective(tao, Xmue, &fe));
179a7e14dcfSSatish Balay       if (fe < fr) {
1809566063dSJacob Faibussowitsch         PetscCall(NelderMeadReplace(nm, nm->indices[nm->N], Xmue, fe));
181a7e14dcfSSatish Balay       } else {
1829566063dSJacob Faibussowitsch         PetscCall(NelderMeadReplace(nm, nm->indices[nm->N], Xmur, fr));
183a7e14dcfSSatish Balay       }
184a7e14dcfSSatish Balay     } else if (nm->f_values[nm->indices[nm->N - 1]] <= fr && fr < nm->f_values[nm->indices[nm->N]]) {
185a7e14dcfSSatish Balay       /* outside contraction */
186a7e14dcfSSatish Balay       nm->noutcontract++;
1879566063dSJacob Faibussowitsch       PetscCall(PetscInfo(0, "Outside Contraction\n"));
1889566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(Xmuc, 1 + nm->mu_oc, -nm->mu_oc, 0, Xbar, nm->simplex[nm->indices[nm->N]]));
189a7e14dcfSSatish Balay 
1909566063dSJacob Faibussowitsch       PetscCall(TaoComputeObjective(tao, Xmuc, &fc));
191*ac530a7eSPierre Jolivet       if (fc <= fr) PetscCall(NelderMeadReplace(nm, nm->indices[nm->N], Xmuc, fc));
192*ac530a7eSPierre Jolivet       else shrink = 1;
193a7e14dcfSSatish Balay     } else {
194a7e14dcfSSatish Balay       /* inside contraction */
195a7e14dcfSSatish Balay       nm->nincontract++;
1969566063dSJacob Faibussowitsch       PetscCall(PetscInfo(0, "Inside Contraction\n"));
1979566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(Xmuc, 1 + nm->mu_ic, -nm->mu_ic, 0, Xbar, nm->simplex[nm->indices[nm->N]]));
1989566063dSJacob Faibussowitsch       PetscCall(TaoComputeObjective(tao, Xmuc, &fc));
199*ac530a7eSPierre Jolivet       if (fc < nm->f_values[nm->indices[nm->N]]) PetscCall(NelderMeadReplace(nm, nm->indices[nm->N], Xmuc, fc));
200*ac530a7eSPierre Jolivet       else shrink = 1;
201a7e14dcfSSatish Balay     }
202a7e14dcfSSatish Balay 
203a7e14dcfSSatish Balay     if (shrink) {
204a7e14dcfSSatish Balay       nm->nshrink++;
2059566063dSJacob Faibussowitsch       PetscCall(PetscInfo(0, "Shrink\n"));
206a7e14dcfSSatish Balay 
207a7e14dcfSSatish Balay       for (i = 1; i < nm->N + 1; i++) {
2089566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(nm->simplex[nm->indices[i]], 1.5, -0.5, nm->simplex[nm->indices[0]]));
2099566063dSJacob Faibussowitsch         PetscCall(TaoComputeObjective(tao, nm->simplex[nm->indices[i]], &nm->f_values[nm->indices[i]]));
210a7e14dcfSSatish Balay       }
2119566063dSJacob Faibussowitsch       PetscCall(VecAXPBY(Xbar, 1.5 * nm->oneOverN, -0.5, nm->simplex[nm->indices[0]]));
212a7e14dcfSSatish Balay 
213a7e14dcfSSatish Balay       /*  Add last vector's fraction of average */
2149566063dSJacob Faibussowitsch       PetscCall(VecAXPY(Xbar, nm->oneOverN, nm->simplex[nm->indices[nm->N]]));
2159566063dSJacob Faibussowitsch       PetscCall(NelderMeadSort(nm));
216a7e14dcfSSatish Balay       /*  Subtract new last vector from average */
2179566063dSJacob Faibussowitsch       PetscCall(VecAXPY(Xbar, -nm->oneOverN, nm->simplex[nm->indices[nm->N]]));
218a7e14dcfSSatish Balay     }
219a7e14dcfSSatish Balay   }
2203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
221a7e14dcfSSatish Balay }
222a7e14dcfSSatish Balay 
2231eb8069cSJason Sarich /*MC
2241eb8069cSJason Sarich  TAONM - Nelder-Mead solver for derivative free, unconstrained minimization
2251eb8069cSJason Sarich 
2261eb8069cSJason Sarich  Options Database Keys:
227a8d3b578SPierre Jolivet + -tao_nm_lambda - initial step length
228a2b725a8SWilliam Gropp - -tao_nm_mu - expansion/contraction factor
2291eb8069cSJason Sarich 
2301eb8069cSJason Sarich  Level: beginner
2311eb8069cSJason Sarich M*/
2321eb8069cSJason Sarich 
233d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_NM(Tao tao)
234d71ae5a4SJacob Faibussowitsch {
235a7e14dcfSSatish Balay   TAO_NelderMead *nm;
236a7e14dcfSSatish Balay 
237a7e14dcfSSatish Balay   PetscFunctionBegin;
2384dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&nm));
239a7e14dcfSSatish Balay   tao->data = (void *)nm;
240a7e14dcfSSatish Balay 
241a7e14dcfSSatish Balay   tao->ops->setup          = TaoSetUp_NM;
242a7e14dcfSSatish Balay   tao->ops->solve          = TaoSolve_NM;
243a7e14dcfSSatish Balay   tao->ops->view           = TaoView_NM;
244a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_NM;
245a7e14dcfSSatish Balay   tao->ops->destroy        = TaoDestroy_NM;
246a7e14dcfSSatish Balay 
2476552cf8aSJason Sarich   /* Override default settings (unless already changed) */
248606f75f6SBarry Smith   PetscCall(TaoParametersInitialize(tao));
249606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, max_it, 2000);
250606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, max_funcs, 4000);
251a7e14dcfSSatish Balay 
25283c8fe1dSLisandro Dalcin   nm->simplex = NULL;
253a8d3b578SPierre Jolivet   nm->lambda  = 1;
254a7e14dcfSSatish Balay 
255a7e14dcfSSatish Balay   nm->mu_ic = -0.5;
256a7e14dcfSSatish Balay   nm->mu_oc = 0.5;
257a7e14dcfSSatish Balay   nm->mu_r  = 1.0;
258a7e14dcfSSatish Balay   nm->mu_e  = 2.0;
2593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
260a7e14dcfSSatish Balay }
261