Actual source code: daapp.c

  1: #include "taodaapplication.h"     /*I "taodaapplication.h" I*/
  2: #include "taoapp.h"
  3: #include "petsc.h"
  4: #include "daapp_impl.h"

  6: int DAAPP_COOKIE=0;
  7: static char PPetscDAAppXXX[] = "PPetscDAApp";

  9: int TaoAppDAApp(TAO_APPLICATION, DA_APPLICATION *);
 10: int DAAppDestroy(DA_APPLICATION);
 11: int DAAppExtend(TAO_APPLICATION);

 15: /* @C
 16:   DAApplicationCreate - Creates a TaoApplication defined
 17: a PETSc Distributed Array (DA) object.  The routines used to evaluate
 18: the objective function, constraint, and derivative information are
 19: designed specifically for these problems.

 21:    Input Parameters:
 22: +  comm - an MPI communicator
 23: .  tda - an array of Distributed Array objects
 24: -  nnda - the number of Distibuted Array objects

 26:    Output Parameters:
 27: .  newdaapplication - the TaoApplication structure

 29: .seealso TaoDAAppSolve(), TaoApplicationCreate(), TaoAppDestroy();

 31:    Level: beginner

 33: .keywords: Application, DA
 34: @ */
 35: int DAApplicationCreate(MPI_Comm comm, DA* tda, int nnda, TAO_APPLICATION* newdaapplication){
 36:   int info;
 37:   TAO_APPLICATION ttapp;
 39:   info = TaoApplicationCreate(comm,&ttapp); CHKERRQ(info);
 40:   info = TaoAppSetDAApp(ttapp,tda,nnda); CHKERRQ(info);
 41:   info = PetscInfo(tda[0],"Creating a DA_APPLICATION object.\n"); CHKERRQ(info);
 42:   *newdaapplication=ttapp;
 43:   return(0);
 44: }

 48: int DAAppExtend(TAO_APPLICATION daapplication){
 49:   int info;
 50:   MPI_Comm comm;
 51:   DA_APPLICATION daapp;

 54:   if (DAAPP_COOKIE==0){
 55:     info=PetscLogClassRegister(&DAAPP_COOKIE,"TAO DA Application");CHKERRQ(info);
 56:   }
 57:   info=PetscObjectGetComm((PetscObject)daapplication,&comm); CHKERRQ(info);
 58:   info = PetscHeaderCreate(daapp,_p_DA_APPLICATION,int,DAAPP_COOKIE,-1,"DA Application",comm,DAAppDestroy,0); CHKERRQ(info);
 59:   info = PetscLogObjectCreate(daapp); CHKERRQ(info);
 60:   info = PetscLogObjectMemory(daapp, sizeof(struct _p_DA_APPLICATION)); CHKERRQ(info);
 61:   
 62:   daapp->nda=0;
 63:   daapp->ndamax=PETSCDAAPPMAXGRIDS;
 64:   daapp->currentlevel=0;
 65:   daapp->IsComplementarity=PETSC_FALSE;
 66:   daapp->kspflag=SAME_NONZERO_PATTERN;
 67:   daapp->computedafunction=0;
 68:   daapp->computedagradient=0;
 69:   daapp->computedafunctiongradient=0;
 70:   daapp->computedahessian=0;
 71:   daapp->usrdafctx=0; daapp->usrdafgctx=0; daapp->usrdagctx=0; daapp->usrdahctx=0;
 72:   daapp->computedabounds=0;
 73:   daapp->bounddactx=0;
 74:   daapp->nbeforemonitors=0;
 75:   daapp->naftermonitors=0;
 76:   int ii;
 77:   info=TaoAppAddObject(daapplication,PPetscDAAppXXX,(void*)daapp,&ii); CHKERRQ(info);
 78:   info=TaoAppSetDestroyRoutine(daapplication,(int(*)(void*))(TaoAppDestroyDAApp), (void*)daapplication); CHKERRQ(info);
 79:   info=TaoAppSetOptionsRoutine(daapplication,DAAppSetOptions); CHKERRQ(info);

 81:   return(0);
 82: }

 86: /*@
 87:   TaoAppSetDAApp - Extend the functionality of a TaoApplication by 
 88: adding support for applications defined
 89: a PETSc Distributed Array (DA) object.  The routines used to evaluate
 90: the objective function, constraint, and derivative information are
 91: designed specifically for these problems.

 93:    Input Parameters:
 94: +  daappliation - an existing application object
 95: .  tda - an array of Distributed Array objects
 96: -  nnda - the number of Distibuted Array objects

 98: .seealso TaoDAAppSolve(), TaoApplicationCreate();

100:    Level: beginner

102: .keywords: Application, DA
103: @*/
104: int TaoAppSetDAApp(TAO_APPLICATION daapplication, DA* tda, int nnda){
105:   int i,info;
106:   DA_APPLICATION daapp;
107:   PetscScalar zero=0.0;
109:   if (nnda > PETSCDAAPPMAXGRIDS){
110:     SETERRQ1(1,"Number of grids cannot exceed %d.\n",PETSCDAAPPMAXGRIDS); }
111:   info = DAAppExtend(daapplication); CHKERRQ(info);
112:   info = TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
113:   info = PetscInfo(daapplication,"TAO Extending TAO_APPLICATION for DA_APPLICATION object.\n"); CHKERRQ(info);
114:   info = DAAppSetMatType(daapplication,MATAIJ); CHKERRQ(info);
115:   daapp->nda=nnda;
116:   daapp->currentlevel=0;
117:   for (i=0;i<nnda; i++){
119:     info = PetscObjectReference((PetscObject)(tda[i])); CHKERRQ(info);
120:     daapp->grid[i].da=tda[i];
121:     info = DACreateGlobalVector(daapp->grid[i].da,&daapp->grid[i].X); CHKERRQ(info);
122:     info = VecDuplicate(daapp->grid[i].X,&daapp->grid[i].R); CHKERRQ(info);
123:     info = VecDuplicate(daapp->grid[i].X,&daapp->grid[i].RHS); CHKERRQ(info);
124:     daapp->grid[i].XL=0;
125:     daapp->grid[i].XU=0;
126:     daapp->grid[i].H=0;
127:     daapp->grid[i].Interpolate=0;
128:     daapp->grid[i].CScale=0;
129:     daapp->grid[i].mgrid=0;
130:     if (i>0){
131:       info = DAGetInterpolation(daapp->grid[i-1].da,daapp->grid[i].da,&daapp->grid[i].Interpolate,&daapp->grid[i].CScale); CHKERRQ(info);
132:     }
133:   }
134:   info = VecSet(daapp->grid[0].X, zero); CHKERRQ(info);
135:   info = TaoAppSetInitialSolutionVec(daapplication,daapp->grid[0].X); CHKERRQ(info);
136:   info = PetscInfo(daapp,"Create work vectors for DA_APPLICATION object.\n"); CHKERRQ(info);

138:   return(0);
139: }

143: int TaoAppDAApp(TAO_APPLICATION daapplication, DA_APPLICATION *daapp){
144:   int info;
145:   DA_APPLICATION daapp2;
148:   info=TaoAppQueryForObject(daapplication,PPetscDAAppXXX,(void**)&daapp2); CHKERRQ(info);
149:   if (daapp2){
151:   } else {
152:     SETERRQ(1,"TAO ERROR: Must call TaoAppSetDAApp() first");
153:   }
154:   *daapp=daapp2;
155:   return(0);
156: }

160: int TaoAppDestroyDAApp(TAO_APPLICATION daapplication){
161:   int info;
162:   DA_APPLICATION daapp;
164:   info=TaoAppQueryForObject(daapplication,PPetscDAAppXXX,(void**)&daapp); CHKERRQ(info);
165:   if (daapp){
166:     info=DAAppDestroy(daapp); CHKERRQ(info);
167:     //    info=TaoAppQueryRemoveObject(daapplication,PPetscDAAppXXX); CHKERRQ(info);
168:   }

170:   return(0);
171: }

175: int DAAppDestroy(DA_APPLICATION daapp){
176:   int i,info,nnda=daapp->nda;


181:   if (--daapp->refct > 0) return(0);

183:   info = PetscInfo(daapp,"Destroying work vectors for DA_APPLICATION object.\n"); CHKERRQ(info);
184:   for (i=0;i<nnda; i++){
185:     if (daapp->grid[i].coloring){
186:       info = ISColoringDestroy(daapp->grid[i].coloring); CHKERRQ(info) daapp->grid[i].coloring=0;
187:     }
188:     if (daapp->grid[i].H){
189:       info = MatDestroy(daapp->grid[i].H); CHKERRQ(info);daapp->grid[i].H=0;
190:     }
191:     if (daapp->grid[i].XL && daapp->grid[i].XU){
192:       info = VecDestroy(daapp->grid[i].XL); CHKERRQ(info);daapp->grid[i].XL=0;
193:       info = VecDestroy(daapp->grid[i].XU); CHKERRQ(info);daapp->grid[i].XU=0;
194:     }
195:     if (i>0){
196:       info = MatDestroy(daapp->grid[i].Interpolate); CHKERRQ(info);daapp->grid[i].Interpolate=0;
197:       info = VecDestroy(daapp->grid[i].CScale); CHKERRQ(info);daapp->grid[i].CScale=0;
198:     }
199:     info = VecDestroy(daapp->grid[i].RHS); CHKERRQ(info);daapp->grid[i].RHS=0;
200:     info = VecDestroy(daapp->grid[i].R); CHKERRQ(info);daapp->grid[i].R=0;
201:     info = VecDestroy(daapp->grid[i].X); CHKERRQ(info);daapp->grid[i].X=0;
202:     info = DADestroy(daapp->grid[i].da); CHKERRQ(info);daapp->grid[i].da=0;
203:   }
204:   daapp->nda=0;
205:   daapp->currentlevel=0;
206:   PetscLogObjectDestroy(daapp);
207:   PetscHeaderDestroy(daapp); 
208:   return(0);
209: }

213: /*@
214:   DAAppGetDA - Get the DA on a the specified level.
215:   

217:    Input Parameters:
218: +  daapplication - the TAO DAApplication structure
219: -  n - the level

221:    Output Parameter:
222: .  da - address of the pointer to the DA.


225: .seealso TaoAppSetDAApp();

227:    Level: intermediate

229: .keywords: Application, DA
230: @*/
231: int DAAppGetDA(TAO_APPLICATION daapplication, int n, DA *da){
232:   int info;
233:   DA_APPLICATION daapp;
235:   info=TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
236:   *da=daapp->grid[n].da;
237:   return(0);
238: }

242: /*@
243:   DAAppGetNumberOfDAGrids - Get the number of grids specified on the application
244:   

246:    Input Parameters:
247: .  daapplication - the TAO DAApplication structure

249:    Output Parameter:
250: .  n - number of DA grids specified

252: .seealso TaoAppSetDAApp();

254:    Level: intermediate

256: .keywords: Application, DA
257: @*/
258: int DAAppGetNumberOfDAGrids(TAO_APPLICATION daapplication, int *n){
259:   int info;
260:   DA_APPLICATION daapp;
262:   info=TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
263:   *n=daapp->nda;
264:   return(0);
265: }

269: /*@
270:   DAAppGetCurrentLevel - Get the number of the current DA grid

272:    Input Parameters:
273: .  daapplication - the TAO DAApplication structure

275:    Output Parameter:
276: .  n - number of DA grids specified

278: .seealso TaoAppSetDAApp();

280:    Level: intermediate

282: .keywords: Application, DA
283: @*/
284: int DAAppGetCurrentLevel(TAO_APPLICATION daapplication, int *n){
285:   int info;
286:   DA_APPLICATION daapp;
288:   info=TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
289:   *n=daapp->currentlevel;
290:   return(0);
291: }


296: static int DAAppPetscHessian(TAO_APPLICATION daapplication , Vec X , Mat *M, Mat *MPre, MatStructure *flag, void*ctx){
297:   int     info, clevel;
298:   Mat     H=*M;
299:   DA_APPLICATION daapp=(DA_APPLICATION)ctx;

304:   PetscStackPush("TAO User DA Hessian");
305:   clevel=daapp->currentlevel;
306:   
307:   //    daapp->FDGrad=daapp->grid[i].RHS;   /* Needed for finite differences gradient vector */
308:   info = PetscInfo1(daapp,"TAO hessian evaluation at DA_APPLICATION object, level %d.\n",clevel); CHKERRQ(info);
309:   info = (*daapp->computedahessian)(daapplication,daapp->grid[clevel].da,X,H,daapp->usrdahctx);CHKERRQ(info);
310:   *flag=daapp->kspflag;

312:   PetscStackPop;

314:   return(0);
315: }


320: /*@
321:   DAAppSetHessianMat - Creates the matrices used to store the Hessian at finer grid levels

323:    Collective on TAO_APPLICATION

325:    Input Parameters:
326: .  daapplication - the DA Application object

328:    Level: advanced


331: .keywords:  DA, hessian

333: .seealso: DAAppSetHessianRoutine(), DAAppSetMatType();

335: @*/
336: int DAAppSetHessianMat(TAO_APPLICATION daapplication){
337:   int i,info,nnda;
338:   DA_APPLICATION daapp;

341:   info=TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
342:   nnda=daapp->nda;
343:   for (i=0;i<nnda; i++){
344:     if (daapp->grid[i].H==0){
345:       info = DAGetMatrix(daapp->grid[i].da,daapp->HessianMatrixType,&daapp->grid[i].H); CHKERRQ(info);
346:       info = MatSetOption(daapp->grid[i].H,MAT_NO_NEW_NONZERO_LOCATIONS); CHKERRQ(info);
347:       info = MatSetOption(daapp->grid[i].H,MAT_STRUCTURALLY_SYMMETRIC); CHKERRQ(info);
348:       info = MatSetOption(daapp->grid[i].H,MAT_SYMMETRIC); CHKERRQ(info);
349:       info = PetscInfo1(daapp->grid[i].da,"TAO create Hessian for DA_APPLICATION object, level %d.\n",i); CHKERRQ(info);
350:       info = DAGetColoring(daapp->grid[i].da,IS_COLORING_GLOBAL,&(daapp->grid[i].coloring));CHKERRQ(info);
351:     }
352:   }
353:   info = TaoAppSetHessianMat(daapplication,daapp->grid[0].H,daapp->grid[0].H);CHKERRQ(info);
354:   
355:   return(0);
356: }
359: /*@C
360:   DAAppSetHessianRoutine - Set a routine that will evaluate the hessian
361:   on the given DA at the given point.

363:    Collective on TAO_APPLICATION

365:    Input Parameters:
366: +  daapplication - the DA Application object
367: .  hess - the function pointer for the hessian evaluation routine
368: -  ctx - the hessian context

370:    Calling sequence of hess:
371: $     hess(TAO_APPLICATION daapplication,DA da, Vec x,Mat H,void *ctx);

373: +  daapplication - the TAO_APPLICATION context
374: .  da - the Distributed Array
375: .  x - input vector
376: .  H - hessian matrix
377: -  ctx - user-defined hessian context set from DAAppSetHessianRoutine()

379:    Level: beginner

381:    Options Database Key:
382: .  -tao_view_hessian - view the hessian after each evaluation using PETSC_VIEWER_STDOUT_WORLD

384: .keywords:  DA, hessian

386: .seealso: DAAppSetObjectiveAndGradientRoutine();
387: @*/
388: int DAAppSetHessianRoutine(TAO_APPLICATION daapplication, int (*hess)(TAO_APPLICATION,DA,Vec,Mat,void*),void *ctx){
389:   int info;
390:   DA_APPLICATION daapp;

393:   info=TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
394:   daapp->computedahessian=hess;
395:   daapp->usrdahctx=ctx;
396:   if (hess){
397:     info = TaoAppSetHessianRoutine(daapplication,DAAppPetscHessian,(void*)daapp);CHKERRQ(info);
398:   } else {
399:     info = TaoAppSetHessianRoutine(daapplication,0,0);CHKERRQ(info);
400:   }
401:   info = DAAppSetHessianMat(daapplication);CHKERRQ(info);
402:   return(0);
403: }

407: static int DAAppPetscFunctionGradient(TAO_APPLICATION daapplication , Vec X ,double *ff, Vec G, void*ctx){
408:   int info;
409:   DA_APPLICATION daapp=(DA_APPLICATION)ctx;

414:   PetscStackPush("TAO User DA Objective and gradient function");
415:   info = (*daapp->computedafunctiongradient)(daapplication,daapp->grid[daapp->currentlevel].da,X,ff,G,daapp->usrdafgctx);
416:   CHKERRQ(info);
417:   info = PetscInfo1(daapplication,"TAO function evaluation at DA_APPLICATION object, level %d.\n",daapp->currentlevel); CHKERRQ(info);
418:   PetscStackPop;

420:   return(0);
421: }


426: /*@C
427:   DAAppSetObjectiveAndGradientRoutine - Set a routine that will evaluate the objective and
428:   gradient functions on the given DA at the given point.

430:    Collective on TAO_APPLICATION

432:    Input Parameters:
433: +  daapplication - the DA Application object
434: .  grad - the function pointer for the gradient evaluation routine
435: -  ctx - the function-gradient context

437:    Calling sequence of funcgrad:
438: $     funcgrad(TAO_APPLICATION daapplication,DA da, Vec x,double *f,Vec g,void *ctx);

440: +  daapplication - the TAO_APPLICATION daapplication context
441: .  da - the Distributed Array
442: .  x - input vector
443: .  f - objective value 
444: .  g - gradient vector 
445: -  ctx - user defined function-gradient context set from DAAppSetObjectiveAndGradientRoutine()

447:    Fortran Note:
448:    If your Fortran compiler does not recognize symbols over 31 characters in length, then
449:    use the identical routine with the shortened name DAAppSetObjectiveAndGradientRou()

451:    Level: beginner

453:    Options Database Key:
454: .  -tao_view_gradient - view the gradient after each evaluation using PETSC_VIEWER_STDOUT_SELF

456: .keywords: DA, Gradient, Objective Function

458: .seealso: DAAppSetObjectiveRoutine(), DAAppSetGradientRoutine();

460: @*/
461: int DAAppSetObjectiveAndGradientRoutine(TAO_APPLICATION daapplication, int (*funcgrad)(TAO_APPLICATION,DA,Vec,double*,Vec, void*),void *ctx){
462:   int info;
463:   DA_APPLICATION daapp;
465:   info=TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
466:   if (funcgrad){ 
467:     info=TaoAppSetObjectiveAndGradientRoutine(daapplication,
468:                                               DAAppPetscFunctionGradient,(void*)daapp);CHKERRQ(info);
469:   } else {
470:     info=TaoAppSetObjectiveAndGradientRoutine(daapplication,0,0);CHKERRQ(info);
471:   }
472:   info=TaoAppSetInitialSolutionVec(daapplication,daapp->grid[0].X);CHKERRQ(info);
473:   daapp->computedafunctiongradient=funcgrad;
474:   daapp->usrdafgctx=ctx;
475:   info = PetscInfo(daapp,"Set objective function pointer for DA_APPLICATION object.\n"); CHKERRQ(info);
476:   return(0);
477: }

481: static int DAAppPetscGradient(TAO_APPLICATION daapplication , Vec X ,Vec G, void*ctx){
482:   int info;
483:   DA_APPLICATION daapp=(DA_APPLICATION)ctx;

488:   PetscStackPush("TAO User DA Gradient Evaluation");
489:   info = (*daapp->computedagradient)(daapplication,daapp->grid[daapp->currentlevel].da,X,G,daapp->usrdafgctx);
490:   CHKERRQ(info);
491:   info = PetscInfo1(daapplication,"TAO gradient evaluation at DA_APPLICATION object, level %d.\n",daapp->currentlevel); CHKERRQ(info);
492:   PetscStackPop;

494:   return(0);
495: }


500: /*@C
501:   DAAppSetGradientRoutine - Set a routine that will evaluate the gradient 
502:   function on the given DA at the given point.

504:    Collective on TAO_APPLICATION

506:    Input Parameters:
507: +  daapplication - the DA Application object
508: .  grad - the function pointer for the gradient evaluation routine
509: -  ctx - the gradient context

511:    Calling sequence of grad:
512: $     grad(TAO_APPLICATION daapplication,DA da, Vec x,Vec g,void *ctx);

514: +  daapplication - the TAO_APPLICATION context
515: .  da - the Distributed Array
516: .  x - input vector
517: .  g - gradient vector 
518: -  ctx - user defined gradient context set from DAAppSetGradientRoutine()

520:    Level: intermediate

522:    Options Database Key:
523: .  -tao_view_gradient - view the gradient after each evaluation using PETSC_VIEWER_STDOUT_SELF

525: .keywords:  DA, gradient

527: .seealso: DAAppSetObjectiveRoutine(), DAAppSetObjectiveAndGradientRoutine();

529: @*/
530: int DAAppSetGradientRoutine(TAO_APPLICATION daapplication, int (*grad)(TAO_APPLICATION,DA,Vec,Vec, void*),void *ctx){
531:   int info;
532:   DA_APPLICATION daapp;
534:   info=TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
535:   if (grad){ 
536:     info=TaoAppSetGradientRoutine(daapplication, DAAppPetscGradient, (void*)daapp);CHKERRQ(info);
537:   } else {
538:     info=TaoAppSetGradientRoutine(daapplication,0,0);CHKERRQ(info);
539:   }
540:   daapp->computedagradient=grad;
541:   daapp->usrdagctx=ctx;
542:   info = PetscInfo(daapp,"Set gradient function pointer for DA_APPLICATION object.\n"); CHKERRQ(info);
543:   return(0);
544: }

548: static int DAAppPetscObjectiveFunction(TAO_APPLICATION daapplication , Vec X ,double* f, void*ctx){
549:   int info;
550:   DA_APPLICATION daapp=(DA_APPLICATION)ctx;

555:   PetscStackPush("TAO User DA Objective function");
556:   info = (*daapp->computedafunction)(daapplication,daapp->grid[daapp->currentlevel].da,X,f,daapp->usrdafgctx);
557:   CHKERRQ(info);
558:   info = PetscInfo1(daapp,"TAO gradient evaluation at DA_APPLICATION object, level %d.\n",daapp->currentlevel); CHKERRQ(info);
559:   PetscStackPop;

561:   return(0);
562: }


567: /*@C
568:   DAAppSetObjectiveRoutine - Set a routine that will evaluate the objective 
569:   function on the given DA at the given point.

571:    Collective on TAO_APPLICATION

573:    Input Parameters:
574: +  daapplication - the DA Application object
575: .  func - the function pointer for the objecive evaluation routine
576: -  ctx - the monitor context

578:    Calling sequence of func:
579: $     func(TAO_APPLICATION daapplication,DA da,Vec x,double *f,void *ctx);

581: +  daapplication - the TAO_APPLICATION context
582: .  da - the Distributed Array
583: .  x - input vector
584: .  f - application sets equal to the function value 
585: -  ctx - user-defined function context set from DAAppSetObjectiveRoutine()

587:    Level: beginner

589: .keywords:  DA, objective

591: .seealso: DAAppSetObjectiveAndGradientRoutine();

593: @*/
594: int DAAppSetObjectiveRoutine(TAO_APPLICATION daapplication, int (*func)(TAO_APPLICATION,DA,Vec,double*, void*),void *ctx){
595:   int info;
596:   DA_APPLICATION daapp;
598:   info=TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
599:   if (func){
600:     info=TaoAppSetObjectiveRoutine(daapplication, DAAppPetscObjectiveFunction, (void*)daapp);CHKERRQ(info);
601:   } else {
602:     info=TaoAppSetObjectiveRoutine(daapplication,0,0);CHKERRQ(info);
603:   }
604:   info=TaoAppSetInitialSolutionVec(daapplication, daapp->grid[0].X);CHKERRQ(info);
605:   daapp->computedafunction=func;
606:   daapp->usrdafctx=ctx;
607:   info = PetscInfo(daapp,"Set objective function pointer for DA_APPLICATION object.\n"); CHKERRQ(info);
608:   return(0);
609: }

613: int DAApplicationEvaluateVariableBounds(TAO_APPLICATION daapplication, Vec XL, Vec XU, void*ctx){
614:   DA_APPLICATION daapp=(DA_APPLICATION)ctx;
615:   int level=daapp->currentlevel, info;

619:   info = PetscInfo1(daapp->grid[level].da,"Compute variable bounds in level %d of DA_APPLICATION object.\n",level); CHKERRQ(info);
620:   info = PetscObjectReference((PetscObject)XL);
621:   info = PetscObjectReference((PetscObject)XU);
622:   daapp->grid[level].XL=XL;
623:   daapp->grid[level].XU=XU;
624:   if (daapp->computedabounds){
625:     info = (*daapp->computedabounds)(daapplication,daapp->grid[level].da,XL,XU,daapp->bounddactx); CHKERRQ(info);
626:   }
627:   return(0);
628: }


633: /*@C
634:   DAAppSetVariableBoundsRoutine - Set a routine that will evaluate the bounds
635:   on the variables.

637:    Collective on TAO_APPLICATION

639:    Input Parameters:
640: +  daapplication - the DA Application object
641: .  bound - the function pointer for the bound evaluation routine
642: -  ctx - the monitor context

644:    Calling sequence of bound:
645: $     bound(TAO_APPLICATION daapplication, DA da,Vec xl, Vec xu, void *ctx);

647: +  daapplication - the TAO_APPLICATION context
648: .  da - the Distributed Array
649: .  xl - vector of upper bounds
650: .  xu - vector of lower bounds
651: -  ctx - user-defined monitor context set from DAAppSetVariableBoundsRoutine()

653:    Level: beginner

655: .keywords:  DA, bounds

657: .seealso: TaoDAAppSolve();

659: @*/
660: int DAAppSetVariableBoundsRoutine(TAO_APPLICATION daapplication, int (*bound)(TAO_APPLICATION,DA,Vec,Vec, void*),void *ctx){
661:   int info;
662:   DA_APPLICATION daapp;
664:   info = TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
665:   info=TaoAppSetVariableBoundsRoutine(daapplication,DAApplicationEvaluateVariableBounds,(void*)daapp);CHKERRQ(info);
666:   daapp->computedabounds=bound;
667:   daapp->bounddactx=ctx;
668:   info = PetscInfo(daapp,"Set variable bounds in DA_APPLICATION object.\n"); CHKERRQ(info);
669:   return(0);
670: }



676: /*@C
677:   DAAppSetConstraintRoutine - Set a routine that will evaluate the constraint
678:   function on the given DA at the given point.

680:    Collective on TAO_APPLICATION

682:    Input Parameters:
683: +  daapplication - the TAO Application object
684: .  grad - the function pointer for the gradient evaluation routine
685: -  ctx - the gradient context

687:    Calling sequence of grad:
688: $     f(TAO_APPLICATION daapplication,DA da, Vec x,Vec r,void *ctx);

690: +  daapplication - the DA_APPLICATION context
691: .  da - the Distributed Array
692: .  x - input vector
693: .  r - constraint vector 
694: -  ctx - user defined gradient context set from DAAppSetGradientRoutine()

696:    Level: intermediate

698:    Options Database Key:
699: .  -tao_view - view the gradient after each evaluation using PETSC_VIEWER_STDOUT_SELF

701: .keywords:  DA, gradient

703: .seealso: DAAppSetJacobianRoutine();
704: @*/
705: int DAAppSetConstraintRoutine(TAO_APPLICATION daapplication, int (*f)(TAO_APPLICATION,DA,Vec,Vec, void*),void *ctx){
706:   int info;
707:   DA_APPLICATION daapp;
709:   info=TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
710:   if (f){ 
711:     info=TaoAppSetConstraintRoutine(daapplication, DAAppPetscGradient, (void*)daapp);CHKERRQ(info);
712:   } else {
713:     info=TaoAppSetConstraintRoutine(daapplication,0,0);CHKERRQ(info);
714:   }
715:   daapp->IsComplementarity=PETSC_TRUE;
716:   info=DAAppSetGradientRoutine(daapplication,f,ctx); CHKERRQ(info);
717:   daapp->computedagradient=f;
718:   daapp->usrdagctx=ctx;
719:   daapp->usrdafgctx=ctx;

721:   info = PetscInfo(daapp,"Set constraint function pointer for DA_APPLICATION object.\n"); CHKERRQ(info);
722:   return(0);
723: }


728: /*@C
729:   DAAppSetJacobianRoutine - Set a routine that will evaluate the Jacobian
730:   on the given DA at the given point.

732:    Collective on TAO_APPLICATION

734:    Input Parameters:
735: +  daapplication - the DA Application object
736: .  jac - the function pointer for the Jacobian evaluation routine
737: -  ctx - the jacobian context

739:    Calling sequence of hess:
740: $     jac(TAO_APPLICATION daapplication, DA da, Vec x,Mat J,void *ctx);

742: +  daapplication - the TAO_APPLICATION context
743: .  da - the Distributed Array
744: .  x - input vector
745: .  J - Jacobian matrix
746: -  ctx - user-defined hessian context set from DAAppSetHessianRoutine()

748:    Level: intermediate

750:    Options Database Key:
751: .  -tao_view_jacobian - view the jacobian after each evaluation using PETSC_VIEWER_STDOUT_WORLD

753: .keywords:  DA, hessian

755: .seealso: DAAppSetConstraintRoutine();

757: @*/
758: int DAAppSetJacobianRoutine(TAO_APPLICATION daapplication, int (*jac)(TAO_APPLICATION,DA,Vec,Mat,void*),void *ctx){
759:   int i,info,nnda;
760:   DA_APPLICATION daapp;
762:   info = TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
763:   if (!jac){ SETERRQ(1,"No DA Hessian routine provided\n");}
764:   daapp->computedahessian=jac;
765:   daapp->usrdahctx=ctx;
766:   nnda=daapp->nda;;
767:   for (i=0;i<nnda; i++){
768:     if (daapp->grid[i].H==0){
769:       info = DAGetMatrix(daapp->grid[i].da,daapp->HessianMatrixType,&daapp->grid[i].H); CHKERRQ(info);
770:       info = MatSetOption(daapp->grid[i].H,MAT_NO_NEW_NONZERO_LOCATIONS); CHKERRQ(info);
771:       info = MatSetOption(daapp->grid[i].H,MAT_STRUCTURALLY_SYMMETRIC); CHKERRQ(info);
772:       info = PetscInfo1(daapp->grid[i].da,"TAO create Jacobian for DA_APPLICATION object, level %d.\n",i); CHKERRQ(info);
773:     }
774:   }
775:   
776:   info = TaoAppSetJacobianRoutine(daapplication,DAAppPetscHessian,(void*)daapp);CHKERRQ(info);
777:   info = TaoAppSetJacobianMat(daapplication,daapp->grid[0].H,daapp->grid[0].H);CHKERRQ(info);
778:   
779:   return(0);
780: }

784: /*@C
785:   DAAppSetBeforeMonitor - Set a routine that will be called before the
786:   optimization on each grid.

788:    Collective on TAO_APPLICATION

790:    Input Parameters:
791: +  daapplication - the DA Application object
792: .  beforemonitor - a monitor routine called before solving the application on each DA 
793: -  ctx - the monitor context

795:    Calling sequence of monitor:
796: $     beforemonitor(TAO_APPLICATION daapplication, DA da,int level, void *ctx);

798: +  daapplication - this TAO_APPLICATION context
799: .  da - the Distributed Array
800: .  level - the grid level that will be solved next (level 0 is coarsest)
801: -  ctx - user-defined function context set from DAAppSetBeforeMonitor()

803: Note:
804:    These monitors are different from the monitors that can be called after
805:    each iteration of the optimization algorithm.

807: Note:
808:    The beforemonitor and aftermonitor are good for setting up and destroying the application data.

810:    Level: intermediate

812: .keywords:  DA, monitor

814: .seealso: DAAppSetAfterMonitor(), TaoSetMonitor(), TaoAppSetMonitor();

816: @*/
817: int DAAppSetBeforeMonitor(TAO_APPLICATION daapplication, int (*beforemonitor)(TAO_APPLICATION,DA,int, void*), void *ctx){
818:   int n,info;
819:   DA_APPLICATION daapp;
821:   info=TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
822:   if (beforemonitor){
823:     n=daapp->nbeforemonitors;
824:     if (n>=MAX_DAAP_MONITORS){
825:       SETERRQ(1,"TAO ERROR: Too many beforemonitors set in DAAPP");
826:     }
827:     daapp->beforemonitor[n]=beforemonitor;
828:     daapp->beforemonitorctx[n]=ctx;
829:     daapp->nbeforemonitors++;
830:     info = PetscInfo(daapplication,"TAO: Set user beforemonitor for DA_APPLICATION object.\n"); CHKERRQ(info);
831:   }
832:   return(0);
833: }

837: /*@C
838:   DAAppSetAfterMonitor - Set a routine that will be called after the
839:   optimization on each grid.

841:    Collective on TAO_APPLICATION

843:    Input Parameters:
844: +  daapplication - the DA Application object
845: .  aftermonitor - a monitor routine called after solving the application on each DA 
846: -  ctx - the monitor context

848:    Calling sequence of monitor:
849: $     aftermonitor(TAO_APPLICATION daapplication, DA da, int level, void *ctx);

851: +  daapplication - this TAO_APPLICATION context
852: .  da - the Distributed Array
853: .  level - the grid level that will be solved next (level 0 is coarsest)
854: -  ctx - user-defined function context set from DAAppSetAfterMonitor()

856: Note:
857:    These monitors are different from the monitors that can be called after
858:    each iteration of the optimization algorithm.

860: Note:
861:    The beforemonitor and aftermonitor are good for setting up and destroying the application data.

863:    Level: intermediate

865: .keywords:  DA, monitor

867: .seealso: DAAppSetBeforeMonitor(), TaoSetMonitor(), TaoAppSetMonitor();

869: @*/
870: int DAAppSetAfterMonitor(TAO_APPLICATION daapplication, int (*aftermonitor)(TAO_APPLICATION,DA,int, void*), void *ctx){
871:   int info,n;
872:   DA_APPLICATION daapp;
874:   info=TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
875:   if (aftermonitor){
876:     n=daapp->naftermonitors;
877:     if (n>=MAX_DAAP_MONITORS){
878:       SETERRQ(1,"TAO ERROR: Too many aftermonitors set in DAAPP");
879:     }
880:     daapp->aftermonitor[n]=aftermonitor;
881:     daapp->aftermonitorctx[n]=ctx;
882:     daapp->naftermonitors++;
883:     info = PetscInfo(daapp,"TAO: Set user aftermonitor for DA_APPLICATION object.\n"); CHKERRQ(info);
884:   }
885:   return(0);
886: }



892: /*@
893:   DAAppGetSolution - Sets a pointer to an existing vector that contains
894:   the solution.

896:    Collective on TAO_APPLICATION

898:    Input Parameters:
899: +  daapplication - the DA Application object
900: -  level - the DA on which the current solution is desired (Level 0 is coarsest)

902:    Output Parameters:
903: .  X - the current solution

905:    Level: beginner

907: .seealso: DAAppSetBeforeMonitor(), DAAppSetAfterMonitor();

909: .keywords: Solution, DA
910: @*/
911: int DAAppGetSolution(TAO_APPLICATION daapplication, int level, Vec *X){
912:   int info;
913:   DA_APPLICATION daapp;
915:   info=TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
916:   if (level<0 || level>=daapp->nda){
917:     SETERRQ(1,"No Solution available.  This grid does not exist");
918:   } else if (X && daapp && daapp->grid && level<daapp->nda){
919:     *X=daapp->grid[level].X;
920:   }
921:   return(0);
922: }


927: /*@C
928:   DAAppGetHessianMat - Sets a pointer to an existing matrix that contains
929:   the Hessian at the specified level.

931:    Collective on TAO_APPLICATION

933:    Input Parameters:
934: +  daapplication - the DA Application object
935: -  level - the DA on which the current Hessian is desired (Level 0 is coarsest)

937:    Output Parameters:
938: .  H - the current Hessian

940:    Level: intermediate

942: .seealso: DAAppSetHessianRoutine;

944: .keywords: Hessian, DA
945: @*/
946: int DAAppGetHessianMat(TAO_APPLICATION daapplication, int level, Mat *H){
947:   int info;
948:   DA_APPLICATION daapp;
950:   info=TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
951:   if (level<0 || level>=daapp->nda){
952:     SETERRQ(1,"No Hessian available.  This grid does not exist");
953:   } else if (H && daapp && daapp->grid && level<daapp->nda){
954:     *H=daapp->grid[level].H;
955:   }
956:   return(0);
957: }

961: /*@
962:   DAAppGetInterpolationMatrix - Sets pointers to existing vectors
963:   containg upper and lower bounds on the variables.

965:    Collective on TAO_APPLICATION

967:    Input Parameters:
968: +  daapplication - the DA Application object
969: -  level - the DA on which the current solution is desired (Level 0 is coarsest)

971:    Output Parameters:
972: +  Interpolate - the interpolating matrix
973: -  CScale - the scaling matrix for restriction Matrix

975:    Level: intermediate

977: .seealso: DAAppSetBeforeMonitor(), DAAppSetAfterMonitor();

979: .keywords: Interpolate, DA
980: @*/
981: int DAAppGetInterpolationMatrix(TAO_APPLICATION daapplication, int level, Mat *Interpolate, Vec *CScale){
982:   int info;
983:   DA_APPLICATION daapp;
985:   info=TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
986:   if ( !daapp || level<1 || level>=daapp->nda){
987:     SETERRQ1(1,"No Interpolator available for level %d.  This grid does not exist",level);
988:   } 
989:   if (Interpolate){
990:     *Interpolate=daapp->grid[level].Interpolate;
991:   }
992:   if (CScale){
993:     *CScale=daapp->grid[level].CScale;
994:   }
995:   return(0);
996: }


1001: /*@
1002:   DAAppGetVariableBounds - Sets pointers to existing vectors
1003:   containg upper and lower bounds on the variables.

1005:    Collective on TAO_APPLICATION

1007:    Input Parameters:
1008: +  daapplication - the DA Application object
1009: -  level - the DA on which the current solution is desired (Level 0 is coarsest)

1011:    Output Parameters:
1012: +  XL - the current solution
1013: -  XU - the current solution

1015:    Level: intermediate

1017: .seealso: DAAppSetBeforeMonitor(), DAAppSetAfterMonitor();

1019: .keywords: Solution, DA
1020: @*/
1021: int DAAppGetVariableBounds(TAO_APPLICATION daapplication, int level, Vec *XL, Vec *XU){
1022:   int info;
1023:   DA_APPLICATION daapp;
1025:   info=TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
1026:   if ( !daapp || level<0 || level>=daapp->nda){
1027:     SETERRQ(1,"No Solution available.  This grid does not exist.");
1028:   } 
1029:   if (XL){
1030:     *XL=daapp->grid[level].XL;
1031:   }
1032:   if (XU){
1033:     *XU=daapp->grid[level].XU;
1034:   }
1035:   return(0);
1036: }



1042: /*@
1043:   DAAppSetInitialSolution - Sets the initial solution for
1044:   the grid application.

1046:    Collective on TAO_APPLICATION

1048:    Input Parameters:
1049: +  daapplication - the DA Application object
1050: -  X0 - the initial solution vector

1052:    Level: beginner

1054: .seealso: DAAppGetSolution();

1056: .keywords: Initial Solution, DA
1057: @*/
1058: int DAAppSetInitialSolution(TAO_APPLICATION daapplication, Vec X0){
1059:   int info;
1060:   DA_APPLICATION daapp;
1063:   info=TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
1064:   info=VecCopy(X0,daapp->grid[0].X); CHKERRQ(info);
1065:   return(0);
1066: }



1072: /*
1073:   DAAppSetMatStructure - Set the MatStructure flag to be used in KSPSetOperators:

1075:   Collective on TAO_APPLICATION
1076:   
1077:   Input Parameters:
1078: +  daapplication - the DA Application object
1079: -  flag - indicates the KSP flag

1081:    Level: intermediate

1083:    Note:
1084:    The default value is SAME_NONZERO_PATTERN

1086: .seealso: KSPSetOperators(), TaoAppGetKSP();

1088: .keywords: Linear Solver, Multigrid, DA, KSP

1090: @*/
1091: int DAAppSetMatStructure(TAO_APPLICATION daapplication, MatStructure pflag){
1092:   int info;
1093:   DA_APPLICATION daapp;

1096:   info=TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
1097:   daapp->kspflag=pflag;
1098:   info = PetscInfo(daapp,"Set MatStructure flag in TAO solver.\n"); CHKERRQ(info);

1100:   return(0);
1101: }

1105: /*@
1106:   DAAppSetMatType - set the matrix type to be used for the Hessian matrix

1108:   Collective on TAO_APPLICATION
1109:   
1110:   Input Parameters:
1111: +  daapplication - the DA Application object
1112: -  mattype - the matrix type

1114:    Level: intermediate

1116:    Options Database Key:
1117: .  -tao_da_mattype - select matrix type

1119:    Note:
1120:    The default value is MATMPIAIJ.

1122:    Note:
1123:    Call before DAAppSetHessianRoutine()

1125: .seealso: DAAppSetMatStructure(),  DAAppSetHessianRoutine(), DAGetMatrix();

1127: .keywords: DA, Hessian

1129: @*/
1130: int DAAppSetMatType(TAO_APPLICATION daapplication, MatType mattype){
1131:   int info;
1132:   DA_APPLICATION daapp;

1135:   info=TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
1136:   info=PetscStrcpy(daapp->HessianMatrixType,mattype); CHKERRQ(info);

1138:   return(0);
1139: }


1144: /*@
1145:   DAAppSetOptions - Sets various parameters to be used in this application
1146:   and the TAO solver.

1148:    Collective on TAO_APPLICATION

1150:    Input Parameters:
1151: .  daapplication - the DA Application object

1153:    Level: intermediate

1155: Note:  This routine will be called by the TaoAppSetFromOptions().

1157: .keywords:  DA, options

1159: .seealso: TaoDAAppSolve();

1161: @*/
1162: int DAAppSetOptions(TAO_APPLICATION daapplication){
1163:   int info;
1164:   //  char *prefix;
1165:   PetscTruth flg1=PETSC_FALSE,flg2=PETSC_FALSE;
1166:   DA_APPLICATION daapp;
1167:   //  MPI_Comm comm;

1171:   info = PetscInfo(daapplication,"TaoDAAppSetOptions(): Reading command line for options\n"); CHKERRQ(info);
1172:   info = TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
1173:   /*
1174:   info = PetscObjectGetOptionsPrefix((PetscObject)daapplication,&prefix); CHKERRQ(info);
1175:   info = PetscObjectGetComm((PetscObject)daapplication,&comm); CHKERRQ(info);
1176:   info = PetscOptionsBegin(comm,prefix,"TAO PETSC DA APPLICATIONS","solver");CHKERRQ(info);
1177:   */
1178:   info = DAAppSetMultiGridOptions(daapplication);CHKERRQ(info);


1181:   flg1=PETSC_FALSE,flg2=PETSC_FALSE;
1182:   info = PetscOptionsTruth("-tao_da_samepreconditioner","Use same preconditioner in KSP","DAAppSetMatStructure",
1183:                              PETSC_FALSE,&flg2,&flg1);CHKERRQ(info);
1184:   if (flg1 && flg2) {
1185:     info = DAAppSetMatStructure(daapplication, SAME_PRECONDITIONER);CHKERRQ(info); 
1186:   } else if (flg1){
1187:     info = DAAppSetMatStructure(daapplication, SAME_NONZERO_PATTERN);CHKERRQ(info); 
1188:   }

1190:   info = PetscOptionsString("-tao_da_mattype","the hessian should have matrix type","DAAppSetMatType",
1191:                             daapp->HessianMatrixType,daapp->HessianMatrixType,18,&flg1);CHKERRQ(info);

1193:   /*
1194:   info = PetscOptionsEnd();CHKERRQ(info);
1195:   */
1196:   return(0);
1197: }