Actual source code: neldermead.c

  1: #include <../src/tao/unconstrained/impls/neldermead/neldermead.h>
  2: #include <petscvec.h>


  5: /*------------------------------------------------------------*/
  6: static PetscErrorCode NelderMeadSort(TAO_NelderMead *nm)
  7: {
  8:   PetscReal *values = nm->f_values;
  9:   PetscInt  *indices = nm->indices;
 10:   PetscInt  dim = nm->N+1;
 11:   PetscInt  i,j,index;
 12:   PetscReal val;

 15:   for (i=1;i<dim;i++) {
 16:     index = indices[i];
 17:     val = values[index];
 18:     for (j=i-1; j>=0 && values[indices[j]] > val; j--) {
 19:       indices[j+1] = indices[j];
 20:     }
 21:     indices[j+1] = index;
 22:   }
 23:   return(0);
 24: }


 27: /*------------------------------------------------------------*/
 28: static PetscErrorCode NelderMeadReplace(TAO_NelderMead *nm, PetscInt index, Vec Xmu, PetscReal f)
 29: {

 33:   /*  Add new vector's fraction of average */
 34:   VecAXPY(nm->Xbar,nm->oneOverN,Xmu);
 35:   VecCopy(Xmu,nm->simplex[index]);
 36:   nm->f_values[index] = f;

 38:   NelderMeadSort(nm);

 40:   /*  Subtract last vector from average */
 41:   VecAXPY(nm->Xbar,-nm->oneOverN,nm->simplex[nm->indices[nm->N]]);
 42:   return(0);
 43: }

 45: /* ---------------------------------------------------------- */
 46: static PetscErrorCode TaoSetUp_NM(Tao tao)
 47: {
 49:   TAO_NelderMead *nm = (TAO_NelderMead *)tao->data;
 50:   PetscInt       n;

 53:   VecGetSize(tao->solution,&n);
 54:   nm->N = n;
 55:   nm->oneOverN = 1.0/n;
 56:   VecDuplicateVecs(tao->solution,nm->N+1,&nm->simplex);
 57:   PetscMalloc1(nm->N+1,&nm->f_values);
 58:   PetscMalloc1(nm->N+1,&nm->indices);
 59:   VecDuplicate(tao->solution,&nm->Xbar);
 60:   VecDuplicate(tao->solution,&nm->Xmur);
 61:   VecDuplicate(tao->solution,&nm->Xmue);
 62:   VecDuplicate(tao->solution,&nm->Xmuc);

 64:   tao->gradient=NULL;
 65:   tao->step=0;
 66:   return(0);
 67: }

 69: /* ---------------------------------------------------------- */
 70: static PetscErrorCode TaoDestroy_NM(Tao tao)
 71: {
 72:   TAO_NelderMead *nm = (TAO_NelderMead*)tao->data;

 76:   if (tao->setupcalled) {
 77:     VecDestroyVecs(nm->N+1,&nm->simplex);
 78:     VecDestroy(&nm->Xmuc);
 79:     VecDestroy(&nm->Xmue);
 80:     VecDestroy(&nm->Xmur);
 81:     VecDestroy(&nm->Xbar);
 82:   }
 83:   PetscFree(nm->indices);
 84:   PetscFree(nm->f_values);
 85:   PetscFree(tao->data);
 86:   return(0);
 87: }

 89: /*------------------------------------------------------------*/
 90: static PetscErrorCode TaoSetFromOptions_NM(PetscOptionItems *PetscOptionsObject,Tao tao)
 91: {
 92:   TAO_NelderMead *nm = (TAO_NelderMead*)tao->data;

 96:   PetscOptionsHead(PetscOptionsObject,"Nelder-Mead options");
 97:   PetscOptionsReal("-tao_nm_lamda","initial step length","",nm->lamda,&nm->lamda,NULL);
 98:   PetscOptionsReal("-tao_nm_mu","mu","",nm->mu_oc,&nm->mu_oc,NULL);
 99:   nm->mu_ic = -nm->mu_oc;
100:   nm->mu_r = nm->mu_oc*2.0;
101:   nm->mu_e = nm->mu_oc*4.0;
102:   PetscOptionsTail();
103:   return(0);
104: }

106: /*------------------------------------------------------------*/
107: static PetscErrorCode TaoView_NM(Tao tao,PetscViewer viewer)
108: {
109:   TAO_NelderMead *nm = (TAO_NelderMead*)tao->data;
110:   PetscBool      isascii;

114:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
115:   if (isascii) {
116:     PetscViewerASCIIPushTab(viewer);
117:     PetscViewerASCIIPrintf(viewer,"expansions: %D\n",nm->nexpand);
118:     PetscViewerASCIIPrintf(viewer,"reflections: %D\n",nm->nreflect);
119:     PetscViewerASCIIPrintf(viewer,"inside contractions: %D\n",nm->nincontract);
120:     PetscViewerASCIIPrintf(viewer,"outside contractionss: %D\n",nm->noutcontract);
121:     PetscViewerASCIIPrintf(viewer,"Shrink steps: %D\n",nm->nshrink);
122:     PetscViewerASCIIPopTab(viewer);
123:   }
124:   return(0);
125: }

127: /*------------------------------------------------------------*/
128: static PetscErrorCode TaoSolve_NM(Tao tao)
129: {
130:   PetscErrorCode     ierr;
131:   TAO_NelderMead     *nm = (TAO_NelderMead*)tao->data;
132:   PetscReal          *x;
133:   PetscInt           i;
134:   Vec                Xmur=nm->Xmur, Xmue=nm->Xmue, Xmuc=nm->Xmuc, Xbar=nm->Xbar;
135:   PetscReal          fr,fe,fc;
136:   PetscInt           shrink;
137:   PetscInt           low,high;

140:   nm->nshrink =      0;
141:   nm->nreflect =     0;
142:   nm->nincontract =  0;
143:   nm->noutcontract = 0;
144:   nm->nexpand =      0;

146:   if (tao->XL || tao->XU || tao->ops->computebounds) {
147:     PetscInfo(tao,"WARNING: Variable bounds have been set but will be ignored by NelderMead algorithm\n");
148:   }

150:   VecCopy(tao->solution,nm->simplex[0]);
151:   TaoComputeObjective(tao,nm->simplex[0],&nm->f_values[0]);
152:   nm->indices[0]=0;
153:   for (i=1;i<nm->N+1;i++){
154:     VecCopy(tao->solution,nm->simplex[i]);
155:     VecGetOwnershipRange(nm->simplex[i],&low,&high);
156:     if (i-1 >= low && i-1 < high) {
157:       VecGetArray(nm->simplex[i],&x);
158:       x[i-1-low] += nm->lamda;
159:       VecRestoreArray(nm->simplex[i],&x);
160:     }

162:     TaoComputeObjective(tao,nm->simplex[i],&nm->f_values[i]);
163:     nm->indices[i] = i;
164:   }

166:   /*  Xbar  = (Sum of all simplex vectors - worst vector)/N */
167:   NelderMeadSort(nm);
168:   VecSet(Xbar,0.0);
169:   for (i=0;i<nm->N;i++) {
170:     VecAXPY(Xbar,1.0,nm->simplex[nm->indices[i]]);
171:   }
172:   VecScale(Xbar,nm->oneOverN);
173:   tao->reason = TAO_CONTINUE_ITERATING;
174:   while (1) {
175:     /* Call general purpose update function */
176:     if (tao->ops->update) {
177:       (*tao->ops->update)(tao, tao->niter, tao->user_update);
178:     }
179:     ++tao->niter;
180:     shrink = 0;
181:     VecCopy(nm->simplex[nm->indices[0]],tao->solution);
182:     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);
183:     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);
184:     (*tao->ops->convergencetest)(tao,tao->cnvP);
185:     if (tao->reason != TAO_CONTINUE_ITERATING) break;

187:     /* x(mu) = (1 + mu)Xbar - mu*X_N+1 */
188:     VecAXPBYPCZ(Xmur,1+nm->mu_r,-nm->mu_r,0,Xbar,nm->simplex[nm->indices[nm->N]]);
189:     TaoComputeObjective(tao,Xmur,&fr);

191:     if (nm->f_values[nm->indices[0]] <= fr && fr < nm->f_values[nm->indices[nm->N-1]]) {
192:       /*  reflect */
193:       nm->nreflect++;
194:       PetscInfo(0,"Reflect\n");
195:       NelderMeadReplace(nm,nm->indices[nm->N],Xmur,fr);
196:     } else if (fr < nm->f_values[nm->indices[0]]) {
197:       /*  expand */
198:       nm->nexpand++;
199:       PetscInfo(0,"Expand\n");
200:       VecAXPBYPCZ(Xmue,1+nm->mu_e,-nm->mu_e,0,Xbar,nm->simplex[nm->indices[nm->N]]);
201:       TaoComputeObjective(tao,Xmue,&fe);
202:       if (fe < fr) {
203:         NelderMeadReplace(nm,nm->indices[nm->N],Xmue,fe);
204:       } else {
205:         NelderMeadReplace(nm,nm->indices[nm->N],Xmur,fr);
206:       }
207:     } else if (nm->f_values[nm->indices[nm->N-1]] <= fr && fr < nm->f_values[nm->indices[nm->N]]) {
208:       /* outside contraction */
209:       nm->noutcontract++;
210:       PetscInfo(0,"Outside Contraction\n");
211:       VecAXPBYPCZ(Xmuc,1+nm->mu_oc,-nm->mu_oc,0,Xbar,nm->simplex[nm->indices[nm->N]]);

213:       TaoComputeObjective(tao,Xmuc,&fc);
214:       if (fc <= fr) {
215:         NelderMeadReplace(nm,nm->indices[nm->N],Xmuc,fc);
216:       } else shrink=1;
217:     } else {
218:       /* inside contraction */
219:       nm->nincontract++;
220:       PetscInfo(0,"Inside Contraction\n");
221:       VecAXPBYPCZ(Xmuc,1+nm->mu_ic,-nm->mu_ic,0,Xbar,nm->simplex[nm->indices[nm->N]]);
222:       TaoComputeObjective(tao,Xmuc,&fc);
223:       if (fc < nm->f_values[nm->indices[nm->N]]) {
224:         NelderMeadReplace(nm,nm->indices[nm->N],Xmuc,fc);
225:       } else shrink = 1;
226:     }

228:     if (shrink) {
229:       nm->nshrink++;
230:       PetscInfo(0,"Shrink\n");

232:       for (i=1;i<nm->N+1;i++) {
233:         VecAXPBY(nm->simplex[nm->indices[i]],1.5,-0.5,nm->simplex[nm->indices[0]]);
234:         TaoComputeObjective(tao,nm->simplex[nm->indices[i]], &nm->f_values[nm->indices[i]]);
235:       }
236:       VecAXPBY(Xbar,1.5*nm->oneOverN,-0.5,nm->simplex[nm->indices[0]]);

238:       /*  Add last vector's fraction of average */
239:       VecAXPY(Xbar,nm->oneOverN,nm->simplex[nm->indices[nm->N]]);
240:       NelderMeadSort(nm);
241:       /*  Subtract new last vector from average */
242:       VecAXPY(Xbar,-nm->oneOverN,nm->simplex[nm->indices[nm->N]]);
243:     }
244:   }
245:   return(0);
246: }

248: /* ---------------------------------------------------------- */
249: /*MC
250:  TAONM - Nelder-Mead solver for derivative free, unconstrained minimization

252:  Options Database Keys:
253: + -tao_nm_lamda - initial step length
254: - -tao_nm_mu - expansion/contraction factor

256:  Level: beginner
257: M*/

259: PETSC_EXTERN PetscErrorCode TaoCreate_NM(Tao tao)
260: {
261:   TAO_NelderMead *nm;

265:   PetscNewLog(tao,&nm);
266:   tao->data = (void*)nm;

268:   tao->ops->setup = TaoSetUp_NM;
269:   tao->ops->solve = TaoSolve_NM;
270:   tao->ops->view = TaoView_NM;
271:   tao->ops->setfromoptions = TaoSetFromOptions_NM;
272:   tao->ops->destroy = TaoDestroy_NM;

274:   /* Override default settings (unless already changed) */
275:   if (!tao->max_it_changed) tao->max_it = 2000;
276:   if (!tao->max_funcs_changed) tao->max_funcs = 4000;

278:   nm->simplex = NULL;
279:   nm->lamda = 1;

281:   nm->mu_ic = -0.5;
282:   nm->mu_oc = 0.5;
283:   nm->mu_r = 1.0;
284:   nm->mu_e = 2.0;

286:   return(0);
287: }