Actual source code: daapp_solve.c

  1: #include "taodaapplication.h"   /*I  "taodaapplication.h"  I*/
  2: //#include "taoapp.h"
  3: #include "tao.h"
  4: #include "petsc.h"
  5: #include "src/petsctao/application/petscapp/tao_app_impl.h"
  6: #include "daapp_impl.h"

  8: static int Tao_DA_Solve=0;
 10: int TaoAppDAApp(TAO_APPLICATION, DA_APPLICATION *);

 14: /*@
 15:   TaoDAAppSolve - Solve the PETSC DA application.

 17:    Input Parameters:
 18: .  daapplication - the TAO DAApplication structure

 20:    Level: beginner

 22: .keywords: Application, DA, Solve
 23: @*/
 24: int TaoDAAppSolve(TAO_APPLICATION daapplication, TAO_SOLVER tao){
 25:   int i,j,mx,my,mz,iter,info;
 26:   double fval,residual;
 27:   TaoTerminateReason reason;
 28:   DA_APPLICATION daapp;


 33:   info=TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);

 35:   if (Tao_DA_Solve==0){
 36:     info=PetscLogEventRegister(&Tao_DA_Solve,"TaoSolveDAApp",DAAPP_COOKIE); CHKERRQ(info);
 37:   }

 39:   info = PetscLogEventBegin(Tao_DA_Solve,tao,daapp,0,0);CHKERRQ(info);
 40:   for (i=0;i<daapp->nda; i++){
 41:     
 42:     
 43:     info=TaoAppResetCounters(daapplication);CHKERRQ(info);
 44:     info = DAGetInfo(daapp->grid[i].da,PETSC_NULL,&mx,&my,&mz,PETSC_NULL,
 45:                      PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,
 46:                      PETSC_NULL,PETSC_NULL);CHKERRQ(info);

 48:     info = PetscInfo2(daapp->grid[i].da,"Level %d of %d in TAO_DA_APPLICATION object.\n",i,daapp->nda); CHKERRQ(info);
 49:     info = PetscInfo3(daapp->grid[i].da,"Global dimensions of DA:  mx=%d, my=%d, mz=%d .\n",mx,my,mz); CHKERRQ(info);

 51:     daapp->currentlevel=i;
 53:     info=TaoAppSetColoring(daapplication, daapp->grid[i].coloring);CHKERRQ(info);

 55:     if (i>0){
 56:       info = TaoSetDown(tao); CHKERRQ(info);
 57:       info = MatMult(daapp->grid[i].Interpolate,daapp->grid[i-1].X,daapp->grid[i].X); CHKERRQ(info);
 58:     } 

 60:     info = TaoAppSetInitialSolutionVec(daapplication,daapp->grid[i].X);CHKERRQ(info);

 62:     if (daapp->grid[i].H){
 63:       if (!daapplication->ksp) {
 64:         MPI_Comm comm;
 65:         info = PetscObjectGetComm((PetscObject)daapp->grid[i].H,&comm); CHKERRQ(info);
 66:         info = KSPCreate(comm,&daapp->grid[i].ksp); CHKERRQ(info);
 67:         info = KSPSetFromOptions(daapp->grid[i].ksp); CHKERRQ(info);
 68:         daapplication->ksp = daapp->grid[i].ksp;
 69:         /*
 70:         if (tao->ksp) {
 71:           tao->ksp->SetKSP(daapp->grid[i].ksp);
 72:         }
 73:         info = TaoWrapKSP(daapp->grid[i].ksp,&tao->ksp); CHKERRQ(info);
 74:         */
 75:       }
 76:       if (daapp->IsComplementarity==PETSC_FALSE){
 77:         info = TaoAppSetHessianMat(daapplication,daapp->grid[i].H,daapp->grid[i].H);CHKERRQ(info);
 78:       } else {
 79:         info = TaoAppSetJacobianMat(daapplication,daapp->grid[i].H,daapp->grid[i].H);CHKERRQ(info);
 80:       }
 81:     }

 83:     info = TaoSetupApplicationSolver(daapplication, tao); CHKERRQ(info);
 84:     for (j=0;j<daapp->nbeforemonitors;j++){
 85:       info = PetscInfo(daapp->grid[i].da,"TAO: Call before user monitor for DA_APPLICATION object.\n"); CHKERRQ(info);
 86:       info = (*daapp->beforemonitor[j])(daapplication,daapp->grid[i].da,i,
 87:                                         daapp->beforemonitorctx[j]); CHKERRQ(info);
 88:     }

 90:     //    info = TaoSetUp(tao); CHKERRQ(info);

 92:     info = PetscInfo1(daapp->grid[i].da,"TAO: Begin solving level %d of DA_APPLICATION object.\n",i); CHKERRQ(info);
 93:     info = TaoSolveApplication(daapplication,tao);CHKERRQ(info);

 95:     info = TaoGetSolutionStatus(tao,&iter,&fval,&residual,0,0,&reason); CHKERRQ(info);
 96:     if (reason <= 0 ){
 97:       info = PetscInfo1(daapp->grid[i].da,"FAILURE TO FIND SOLUTION at level %d of DA_APPLICATION.\n",i+1); CHKERRQ(info);    
 98:       info = PetscInfo1(daapp->grid[i].da,"  TAO Reason for termination: %d.\n",reason); CHKERRQ(info);
 99:     } else {
100:       info = PetscInfo1(daapp->grid[i].da,"Found solution to DA_APPLICATION at level: %d.\n",i+1);CHKERRQ(info);
101:       info = PetscInfo3(daapp->grid[i].da,"Iterations: %d, Objective Value: %10.8e, Residual: %4.2e.\n",
102:                            iter,fval,residual);CHKERRQ(info);
103:     }

105:     for (j=0;j<daapp->naftermonitors;j++){
106:       info = PetscInfo(daapp->grid[i].da,"TAO: Call after user monitor for DA_APPLICATION object.\n"); CHKERRQ(info);
107:       info = (*daapp->aftermonitor[j])(daapplication,daapp->grid[i].da,i,
108:                                        daapp->aftermonitorctx[j]); CHKERRQ(info);
109:     }

111:     if (daapplication->ksp) {
112:       info = KSPDestroy(daapplication->ksp); CHKERRQ(info);
113:       daapplication->ksp=0;
114:     }
115:     if (i < daapp->nda-1){
116:       info = TaoSetDown(tao); CHKERRQ(info);
117:     } 
118:     
119:   }

121:   info = PetscLogEventEnd(Tao_DA_Solve,tao,daapp,0,0);CHKERRQ(info);
122:   return(0);
123: }

125: