Actual source code: petscdmmg.h

  1: /*
  2:   Defines the interface functions for the DMMG object.
  3: */
  4: #ifndef __PETSCDMMG_H
 6:  #include petscsnes.h
 7:  #include petscda.h

 10: /*S
 11:      DMMG -  Data structure to easily manage multi-level non-linear solvers on grids managed by DM
 12:           
 13:    Level: intermediate

 15:   Concepts: multigrid, Newton-multigrid

 17: .seealso:  VecPackCreate(), DA, VecPack, DM, DMMGCreate(), DMMGSetKSP(), DMMGSetSNES()
 18: S*/
 19: typedef struct _p_DMMG *DMMG;
 20: struct _p_DMMG {
 21:   DM             dm;                   /* grid information for this level */
 22:   Vec            x,b,r;                /* global vectors used in multigrid preconditioner for this level*/
 23:   Mat            J;                    /* matrix on this level */
 24:   Mat            R;                    /* restriction to next coarser level (not defined on level 0) */
 25:   PetscInt       nlevels;              /* number of levels above this one (total number of levels on level 0)*/
 26:   MPI_Comm       comm;
 27:   PetscErrorCode (*solve)(DMMG*,PetscInt);
 28:   void           *user;
 29:   PetscTruth     galerkin;                  /* for A_c = R*A*R^T */

 31:   /* KSP only */
 32:   KSP            ksp;
 33:   PetscErrorCode (*rhs)(DMMG,Vec);
 34:   PetscTruth     matricesset;               /* User had called DMMGSetKSP() and the matrices have been computed */

 36:   /* SNES only */
 37:   Mat            B;
 38:   Vec            Rscale;                 /* scaling to restriction before computing Jacobian */
 39:   PetscErrorCode (*computejacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
 40:   PetscErrorCode (*computefunction)(SNES,Vec,Vec,void*);

 42:   PetscTruth     updatejacobian;         /* compute new Jacobian when DMMGComputeJacobian_Multigrid() is called */
 43:   PetscInt       updatejacobianperiod;   /* how often, inside a SNES, the Jacobian is recomputed */

 45:   MatFDColoring  fdcoloring;             /* only used with FD coloring for Jacobian */
 46:   SNES           snes;
 47:   PetscErrorCode (*initialguess)(DMMG,Vec);
 48:   Vec            w,work1,work2;         /* global vectors */
 49:   Vec            lwork1;

 51:   /* FAS only */
 52:   NLF            nlf;                   /* FAS smoother object */
 53:   VecScatter     inject;                /* inject from this level to the next coarsest */
 54:   PetscTruth     monitor,monitorall;
 55:   PetscInt       presmooth,postsmooth,coarsesmooth;
 56:   PetscReal      rtol,abstol,rrtol;       /* convergence tolerance */
 57: 
 58: };

 60: EXTERN PetscErrorCode PETSCSNES_DLLEXPORT DMMGCreate(MPI_Comm,PetscInt,void*,DMMG**);
 61: EXTERN PetscErrorCode PETSCSNES_DLLEXPORT DMMGDestroy(DMMG*);
 62: EXTERN PetscErrorCode PETSCSNES_DLLEXPORT DMMGSetUp(DMMG*);
 63: EXTERN PetscErrorCode PETSCSNES_DLLEXPORT DMMGSetKSP(DMMG*,PetscErrorCode (*)(DMMG,Vec),PetscErrorCode (*)(DMMG,Mat));
 64: EXTERN PetscErrorCode PETSCSNES_DLLEXPORT DMMGSetSNES(DMMG*,PetscErrorCode (*)(SNES,Vec,Vec,void*),PetscErrorCode (*)(SNES,Vec,Mat*,Mat*,MatStructure*,void*));
 65: EXTERN PetscErrorCode PETSCSNES_DLLEXPORT DMMGSetInitialGuess(DMMG*,PetscErrorCode (*)(DMMG,Vec));
 66: EXTERN PetscErrorCode PETSCSNES_DLLEXPORT DMMGInitialGuessCurrent(DMMG,Vec);
 67: EXTERN PetscErrorCode PETSCSNES_DLLEXPORT DMMGView(DMMG*,PetscViewer);
 68: EXTERN PetscErrorCode PETSCSNES_DLLEXPORT DMMGSolve(DMMG*);
 69: EXTERN PetscErrorCode PETSCSNES_DLLEXPORT DMMGSetUseMatrixFree(DMMG*);
 70: EXTERN PetscErrorCode PETSCSNES_DLLEXPORT DMMGSetDM(DMMG*,DM);
 71: EXTERN PetscErrorCode PETSCSNES_DLLEXPORT DMMGSetUpLevel(DMMG*,KSP,PetscInt);
 72: EXTERN PetscErrorCode PETSCSNES_DLLEXPORT DMMGSetUseGalerkinCoarse(DMMG*);
 73: EXTERN PetscErrorCode PETSCSNES_DLLEXPORT DMMGSetNullSpace(DMMG*,PetscTruth,PetscInt,PetscErrorCode (*)(DMMG,Vec[]));

 75: EXTERN PetscErrorCode PETSCSNES_DLLEXPORT DMMGSetSNESLocal_Private(DMMG*,DALocalFunction1,DALocalFunction1,DALocalFunction1,DALocalFunction1);
 76: #if defined(PETSC_HAVE_ADIC)
 77: #  define DMMGSetSNESLocal(dmmg,function,jacobian,ad_function,admf_function) \
 78:   DMMGSetSNESLocal_Private(dmmg,(DALocalFunction1)function,(DALocalFunction1)jacobian,(DALocalFunction1)(ad_function),(DALocalFunction1)(admf_function))
 79: #else
 80: #  define DMMGSetSNESLocal(dmmg,function,jacobian,ad_function,admf_function) DMMGSetSNESLocal_Private(dmmg,(DALocalFunction1)function,(DALocalFunction1)jacobian,(DALocalFunction1)0,(DALocalFunction1)0)
 81: #endif

 83: EXTERN PetscErrorCode PETSCSNES_DLLEXPORT DMMGSetSNESLocali_Private(DMMG*,PetscErrorCode(*)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*),PetscErrorCode(*)(DALocalInfo*,MatStencil*,void*,void*,void*),PetscErrorCode(*)(DALocalInfo*,MatStencil*,void*,void*,void*));
 84: #if defined(PETSC_HAVE_ADIC)
 85: #  define DMMGSetSNESLocali(dmmg,function,ad_function,admf_function) DMMGSetSNESLocali_Private(dmmg,(PetscErrorCode(*)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*))function,(PetscErrorCode(*)(DALocalInfo*,MatStencil*,void*,void*,void*))(ad_function),(PetscErrorCode(*)(DALocalInfo*,MatStencil*,void*,void*,void*))(admf_function))
 86: #else
 87: #  define DMMGSetSNESLocali(dmmg,function,ad_function,admf_function) DMMGSetSNESLocali_Private(dmmg,(PetscErrorCode(*)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*))function,0,0)
 88: #endif

 90: /*MC
 91:    DMMGGetRHS - Returns the right hand side vector from a DMMG solve on the finest grid

 93:    Synopsis:
 94:    Vec DMMGGetRHS(DMMG *dmmg)

 96:    Not Collective, but resulting vector is parallel

 98:    Input Parameters:
 99: .   dmmg - DMMG solve context

101:    Level: intermediate

103:    Fortran Usage:
104: .     DMMGGetRHS(DMMG dmmg,Vec b,PetscErrorCode ierr)

106: .seealso: DMMGCreate(), DMMGSetSNES(), DMMGSetKSP(), DMMGSetSNESLocal(), DMMGGetRHS()

108: M*/
109: #define DMMGGetRHS(ctx)              (ctx)[(ctx)[0]->nlevels-1]->b

111: #define DMMGGetr(ctx)              (ctx)[(ctx)[0]->nlevels-1]->r

113: /*MC
114:    DMMGGetx - Returns the solution vector from a DMMG solve on the finest grid

116:    Synopsis:
117:    Vec DMMGGetx(DMMG *dmmg)

119:    Not Collective, but resulting vector is parallel

121:    Input Parameters:
122: .   dmmg - DMMG solve context

124:    Level: intermediate

126:    Fortran Usage:
127: .     DMMGGetx(DMMG dmmg,Vec x,PetscErrorCode ierr)

129: .seealso: DMMGCreate(), DMMGSetSNES(), DMMGSetKSP(), DMMGSetSNESLocal()

131: M*/
132: #define DMMGGetx(ctx)              (ctx)[(ctx)[0]->nlevels-1]->x

134: /*MC
135:    DMMGGetJ - Returns the Jacobian (matrix) for the finest level

137:    Synopsis:
138:    Mat DMMGGetJ(DMMG *dmmg)

140:    Not Collective

142:    Input Parameter:
143: .   dmmg - DMMG solve context

145:    Level: intermediate

147: .seealso: DMMGCreate(), DMMGSetUser(), DMMGGetB(), DMMGGetRHS()

149: M*/
150: #define DMMGGetJ(ctx)              (ctx)[(ctx)[0]->nlevels-1]->J

152: /*MC
153:    DMMGGetComm - Returns the MPI_Comm for the finest level

155:    Synopsis:
156:    MPI_Comm DMMGGetJ(DMMG *dmmg)

158:    Not Collective

160:    Input Parameter:
161: .   dmmg - DMMG solve context

163:    Level: intermediate

165: .seealso: DMMGCreate(), DMMGSetUser(), DMMGGetJ()

167: M*/
168: #define DMMGGetComm(ctx)           (ctx)[(ctx)[0]->nlevels-1]->comm

170: /*MC
171:    DMMGGetB - Returns the matrix for the finest level used to construct the preconditioner; usually 
172:               the same as the Jacobian

174:    Synopsis:
175:    Mat DMMGGetJ(DMMG *dmmg)

177:    Not Collective

179:    Input Parameter:
180: .   dmmg - DMMG solve context

182:    Level: intermediate

184: .seealso: DMMGCreate(), DMMGSetUser(), DMMGGetJ()

186: M*/
187: #define DMMGGetB(ctx)              (ctx)[(ctx)[0]->nlevels-1]->B

189: /*MC
190:    DMMGGetFine - Returns the DMMG associated with the finest level

192:    Synopsis:
193:    DMMG DMMGGetFine(DMMG *dmmg)

195:    Not Collective

197:    Input Parameter:
198: .   dmmg - DMMG solve context

200:    Level: intermediate

202: .seealso: DMMGCreate(), DMMGSetUser(), DMMGGetJ()

204: M*/
205: #define DMMGGetFine(ctx)           (ctx)[(ctx)[0]->nlevels-1]


208: /*MC
209:    DMMGGetKSP - Gets the KSP object (linear solver object) for the finest level

211:    Synopsis:
212:    KSP DMMGGetKSP(DMMG *dmmg)

214:    Not Collective

216:    Input Parameter:
217: .   dmmg - DMMG solve context

219:    Level: intermediate

221:    Notes: If this is a linear problem (i.e. DMMGSetKSP() was used) then this is the 
222:      master linear solver. If this is a nonlinear problem (i.e. DMMGSetSNES() was used) this
223:      returns the KSP (linear solver) that is associated with the SNES (nonlinear solver)

225: .seealso: DMMGCreate(), DMMGSetUser(), DMMGGetJ(), KSPGetSNES()

227: M*/
228: #define DMMGGetKSP(ctx)            (ctx)[(ctx)[0]->nlevels-1]->ksp

230: /*MC
231:    DMMGGetSNES - Gets the SNES object (nonlinear solver) for the finest level

233:    Synopsis:
234:    SNES DMMGGetSNES(DMMG *dmmg)

236:    Not Collective

238:    Input Parameter:
239: .   dmmg - DMMG solve context

241:    Level: intermediate

243:    Notes: If this is a linear problem (i.e. DMMGSetKSP() was used) then this returns PETSC_NULL

245: .seealso: DMMGCreate(), DMMGSetUser(), DMMGGetJ(), KSPGetKSP()

247: M*/
248: #define DMMGGetSNES(ctx)           (ctx)[(ctx)[0]->nlevels-1]->snes

250: /*MC
251:    DMMGGetDA - Gets the DA object on the finest level

253:    Synopsis:
254:    DA DMMGGetDA(DMMG *dmmg)

256:    Not Collective

258:    Input Parameter:
259: .   dmmg - DMMG solve context

261:    Level: intermediate

263:    Notes: Use only if the DMMG was created with a DA, not a VecPack

265: .seealso: DMMGCreate(), DMMGSetUser(), DMMGGetJ(), KSPGetKSP(), DMMGGetVecPack()

267: M*/
268: #define DMMGGetDA(ctx)             (DA)((ctx)[(ctx)[0]->nlevels-1]->dm)

270: /*MC
271:    DMMGGetVecPack - Gets the VecPack object on the finest level

273:    Synopsis:
274:    VecPack DMMGGetVecPack(DMMG *dmmg)

276:    Not Collective

278:    Input Parameter:
279: .   dmmg - DMMG solve context

281:    Level: intermediate

283:    Notes: Use only if the DMMG was created with a DA, not a VecPack

285: .seealso: DMMGCreate(), DMMGSetUser(), DMMGGetJ(), KSPGetKSP(), DMMGGetDA()

287: M*/
288: #define DMMGGetVecPack(ctx)        (VecPack)((ctx)[(ctx)[0]->nlevels-1]->dm)

290: /*MC
291:    DMMGGetUser - Returns the user context for a particular level

293:    Synopsis:
294:    void* DMMGGetUser(DMMG *dmmg,PetscInt level)

296:    Not Collective

298:    Input Parameters:
299: +   dmmg - DMMG solve context
300: -   level - the number of the level you want the context for

302:    Level: intermediate

304: .seealso: DMMGCreate(), DMMGSetUser()

306: M*/
307: #define DMMGGetUser(ctx,level)     ((ctx)[level]->user)

309: /*MC
310:    DMMGSetUser - Sets the user context for a particular level

312:    Synopsis:
313:    PetscErrorCode DMMGSetUser(DMMG *dmmg,PetscInt level,void *ctx)

315:    Not Collective

317:    Input Parameters:
318: +   dmmg - DMMG solve context
319: .   level - the number of the level you want the context for
320: -   ctx - the context

322:    Level: intermediate

324:    Note: if the context is the same for each level just pass it in with 
325:          DMMGCreate() and don't call this macro

327: .seealso: DMMGCreate(), DMMGGetUser()

329: M*/
330: #define DMMGSetUser(ctx,level,usr) ((ctx)[level]->user = usr,0)

332: /*MC
333:    DMMGGetLevels - Gets the number of levels in a DMMG object

335:    Synopsis:
336:    PetscInt DMMGGetLevels(DMMG *dmmg)

338:    Not Collective

340:    Input Parameter:
341: .   dmmg - DMMG solve context

343:    Level: intermediate

345: .seealso: DMMGCreate(), DMMGGetUser()

347: M*/
348: #define DMMGGetLevels(ctx)         (ctx)[0]->nlevels

350: /*MC
351:    DMMGGetDMMG - Returns the DMMG struct for the finest level

353:    Synopsis:
354:    DMMG DMMGGetDMMG(DMMG *dmmg)

356:    Not Collective

358:    Input Parameter:
359: .   dmmg - DMMG solve context

361:    Level: intermediate

363: .seealso: DMMGCreate(), DMMGSetUser(), DMMGGetB()

365: M*/
366: #define DMMGGetDMMG(ctx)              (ctx)[(ctx)[0]->nlevels-1]

369: #endif