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:      DMMGArray - Fortran only. This is used in the main program when doing DMMGCreate(), DMMGSetDM() etc.
 12:         in the subroutines like FormFunction() one should use DMMG.

 14:         You can use DMMGArrayGetDMMG(DMMGArray,DMMG,ierr) to obtain the DMMG from a DMMG.

 16:    Level: intermediate

 18:   Concepts: multigrid, Newton-multigrid

 20: .seealso:  DMCompositeCreate(), DA, DMComposite, DM, DMMGCreate(), DMMGSetKSP(), DMMGSetSNES(), DMMGSetInitialGuess(),
 21:            DMMGSetNullSpace(), DMMGSetUseGalerkin(), DMMGSetMatType()
 22: S*/

 24: /*S
 25:      DMMG -  Data structure to easily manage multi-level non-linear solvers on grids managed by DM
 26:           
 27:    Level: intermediate

 29:    Fortran Users: see also DMMGArray

 31:   Concepts: multigrid, Newton-multigrid

 33: .seealso:  DMCompositeCreate(), DA, DMComposite, DM, DMMGCreate(), DMMGSetKSP(), DMMGSetSNES(), DMMGSetInitialGuess(),
 34:            DMMGSetNullSpace(), DMMGSetUseGalerkin(), DMMGSetMatType()
 35: S*/
 36: typedef struct _n_DMMG* DMMG;
 37: struct _n_DMMG {
 38:   DM             dm;                   /* grid information for this level */
 39:   Vec            x,b,r;                /* global vectors used in multigrid preconditioner for this level*/
 40:   Mat            J;                    /* matrix on this level */
 41:   Mat            B;
 42:   Mat            R;                    /* restriction to next coarser level (not defined on level 0) */
 43:   PetscInt       nlevels;              /* number of levels above this one (total number of levels on level 0)*/
 44:   MPI_Comm       comm;
 45:   PetscErrorCode (*solve)(DMMG*,PetscInt);
 46:   void           *user;
 47:   PetscTruth     galerkin;                  /* for A_c = R*A*R^T */
 48:   MatType        mtype;                     /* create matrices of this type */
 49:   char           *prefix;

 51:   /* KSP only */
 52:   KSP            ksp;
 53:   PetscErrorCode (*rhs)(DMMG,Vec);

 55:   /* SNES only */
 56:   Vec            Rscale;                 /* scaling to restriction before computing Jacobian */
 57:   PetscErrorCode (*computejacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
 58:   PetscErrorCode (*computefunction)(SNES,Vec,Vec,void*);

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

 63:   ISColoringType isctype;
 64:   MatFDColoring  fdcoloring;             /* only used with FD coloring for Jacobian */
 65:   SNES           snes;
 66:   PetscErrorCode (*initialguess)(DMMG,Vec);
 67:   Vec            w,work1,work2;         /* global vectors */
 68:   Vec            lwork1;

 70:   PetscErrorCode (*lfj)(void);          /* function used when computing Jacobian via FD, usually da->lf */

 72:   /* FAS only */
 73:   NLF            nlf;                   /* FAS smoother object */
 74:   VecScatter     inject;                /* inject from this level to the next coarsest */
 75:   PetscTruth     monitor,monitorall;
 76:   PetscInt       presmooth,postsmooth,coarsesmooth;
 77:   PetscReal      rtol,abstol,rrtol;       /* convergence tolerance */
 78: 
 79: };

 81: EXTERN PetscErrorCode  DMMGCreate(MPI_Comm,PetscInt,void*,DMMG**);
 82: EXTERN PetscErrorCode  DMMGDestroy(DMMG*);
 83: EXTERN PetscErrorCode  DMMGSetUp(DMMG*);
 84: EXTERN PetscErrorCode  DMMGSetKSP(DMMG*,PetscErrorCode (*)(DMMG,Vec),PetscErrorCode (*)(DMMG,Mat,Mat));
 85: EXTERN PetscErrorCode  DMMGSetSNES(DMMG*,PetscErrorCode (*)(SNES,Vec,Vec,void*),PetscErrorCode (*)(SNES,Vec,Mat*,Mat*,MatStructure*,void*));

 87: EXTERN PetscErrorCode  DMMGSetInitialGuessLocal(DMMG*,PetscErrorCode (*)(void));
 88: EXTERN PetscErrorCode  DMMGSetInitialGuess(DMMG*,PetscErrorCode (*)(DMMG,Vec));
 89: EXTERN PetscErrorCode  DMMGInitialGuessCurrent(DMMG,Vec);
 90: EXTERN PetscErrorCode  DMMGView(DMMG*,PetscViewer);
 91: EXTERN PetscErrorCode  DMMGSolve(DMMG*);
 92: EXTERN PetscErrorCode  DMMGSetUseMatrixFree(DMMG*);
 93: EXTERN PetscErrorCode  DMMGSetDM(DMMG*,DM);
 94: EXTERN PetscErrorCode  DMMGSetUpLevel(DMMG*,KSP,PetscInt);
 95: EXTERN PetscErrorCode  DMMGSetUseGalerkinCoarse(DMMG*);
 96: EXTERN PetscErrorCode  DMMGSetNullSpace(DMMG*,PetscTruth,PetscInt,PetscErrorCode (*)(DMMG,Vec[]));
 97: EXTERN PetscErrorCode  DMMGSetMatType(DMMG*,MatType);
 98: EXTERN PetscErrorCode  DMMGSetISColoringType(DMMG*,ISColoringType);
 99: EXTERN PetscErrorCode  DMMGSetPrefix(DMMG*,const char*);

101: EXTERN PetscErrorCode  DMMGSetSNESLocal_Private(DMMG*,DALocalFunction1,DALocalFunction1,DALocalFunction1,DALocalFunction1);
102: #if defined(PETSC_HAVE_ADIC)
103: #  define DMMGSetSNESLocal(dmmg,function,jacobian,ad_function,admf_function) \
104:   DMMGSetSNESLocal_Private(dmmg,(DALocalFunction1)function,(DALocalFunction1)jacobian,(DALocalFunction1)(ad_function),(DALocalFunction1)(admf_function))
105: #else
106: #  define DMMGSetSNESLocal(dmmg,function,jacobian,ad_function,admf_function) DMMGSetSNESLocal_Private(dmmg,(DALocalFunction1)function,(DALocalFunction1)jacobian,(DALocalFunction1)0,(DALocalFunction1)0)
107: #endif

109: EXTERN PetscErrorCode  DMMGSetSNESLocali_Private(DMMG*,PetscErrorCode(*)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*),PetscErrorCode(*)(DALocalInfo*,MatStencil*,void*,void*,void*),PetscErrorCode(*)(DALocalInfo*,MatStencil*,void*,void*,void*));
110: #if defined(PETSC_HAVE_ADIC)
111: #  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))
112: #else
113: #  define DMMGSetSNESLocali(dmmg,function,ad_function,admf_function) DMMGSetSNESLocali_Private(dmmg,(PetscErrorCode(*)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*))function,0,0)
114: #endif

116: EXTERN PetscErrorCode  DMMGSetSNESLocalib_Private(DMMG*,PetscErrorCode(*)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*),PetscErrorCode(*)(DALocalInfo*,MatStencil*,void*,void*,void*),PetscErrorCode(*)(DALocalInfo*,MatStencil*,void*,void*,void*));
117: #if defined(PETSC_HAVE_ADIC)
118: #  define DMMGSetSNESLocalib(dmmg,function,ad_function,admf_function) DMMGSetSNESLocalib_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))
119: #else
120: #  define DMMGSetSNESLocalib(dmmg,function,ad_function,admf_function) DMMGSetSNESLocalib_Private(dmmg,(PetscErrorCode(*)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*))function,0,0)
121: #endif


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

128:    Synopsis:
129:    Vec DMMGGetRHS(DMMG *dmmg)

131:    Not Collective, but resulting vector is parallel

133:    Input Parameters:
134: .   dmmg - DMMG solve context

136:    Level: intermediate

138:    Fortran Usage:
139: .     DMMGGetRHS(DMMG dmmg,Vec b,PetscErrorCode ierr)

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

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

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

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

151:    Synopsis:
152:    Vec DMMGGetx(DMMG *dmmg)

154:    Not Collective, but resulting vector is parallel

156:    Input Parameters:
157: .   dmmg - DMMG solve context

159:    Level: intermediate

161:    Fortran Usage:
162: .     DMMGGetx(DMMG dmmg,Vec x,PetscErrorCode ierr)

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

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

169: /*MC
170:    DMMGGetJ - Returns the Jacobian (matrix) for the finest level

172:    Synopsis:
173:    Mat DMMGGetJ(DMMG *dmmg)

175:    Not Collective

177:    Input Parameter:
178: .   dmmg - DMMG solve context

180:    Level: intermediate

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

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

187: /*MC
188:    DMMGGetComm - Returns the MPI_Comm for the finest level

190:    Synopsis:
191:    MPI_Comm DMMGGetJ(DMMG *dmmg)

193:    Not Collective

195:    Input Parameter:
196: .   dmmg - DMMG solve context

198:    Level: intermediate

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

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

205: /*MC
206:    DMMGGetB - Returns the matrix for the finest level used to construct the preconditioner; usually 
207:               the same as the Jacobian

209:    Synopsis:
210:    Mat DMMGGetJ(DMMG *dmmg)

212:    Not Collective

214:    Input Parameter:
215: .   dmmg - DMMG solve context

217:    Level: intermediate

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

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

224: /*MC
225:    DMMGGetFine - Returns the DMMG associated with the finest level

227:    Synopsis:
228:    DMMG DMMGGetFine(DMMG *dmmg)

230:    Not Collective

232:    Input Parameter:
233: .   dmmg - DMMG solve context

235:    Level: intermediate

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

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


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

246:    Synopsis:
247:    KSP DMMGGetKSP(DMMG *dmmg)

249:    Not Collective

251:    Input Parameter:
252: .   dmmg - DMMG solve context

254:    Level: intermediate

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

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

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

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

268:    Synopsis:
269:    SNES DMMGGetSNES(DMMG *dmmg)

271:    Not Collective

273:    Input Parameter:
274: .   dmmg - DMMG solve context

276:    Level: intermediate

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

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

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

285: /*MC
286:    DMMGGetDM - Gets the DM object on the finest level

288:    Synopsis:
289:    DM DMMGGetDM(DMMG *dmmg)

291:    Not Collective

293:    Input Parameter:
294: .   dmmg - DMMG solve context

296:    Level: intermediate

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

300: M*/
301: #define DMMGGetDM(ctx)             ((ctx)[(ctx)[0]->nlevels-1]->dm)

303: /*MC
304:    DMMGGetDA - Gets the DA object on the finest level

306:    Synopsis:
307:    DA DMMGGetDA(DMMG *dmmg)

309:    Not Collective

311:    Input Parameter:
312: .   dmmg - DMMG solve context

314:    Level: intermediate

316:    Notes: Use only if the DMMG was created with a DA, not a DMComposite

318: .seealso: DMMGCreate(), DMMGSetUser(), DMMGGetJ(), KSPGetKSP(), DMMGGetDMComposite()

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

323: /*MC
324:    DMMGGetDMComposite - Gets the DMComposite object on the finest level

326:    Synopsis:
327:    DMComposite DMMGGetDMComposite(DMMG *dmmg)

329:    Not Collective

331:    Input Parameter:
332: .   dmmg - DMMG solve context

334:    Level: intermediate

336:    Notes: Use only if the DMMG was created with a DA, not a DMComposite

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

340: M*/
341: #define DMMGGetDMComposite(ctx)        (DMComposite)((ctx)[(ctx)[0]->nlevels-1]->dm)

343: /*MC
344:    DMMGGetUser - Returns the user context for a particular level

346:    Synopsis:
347:    void* DMMGGetUser(DMMG *dmmg,PetscInt level)

349:    Not Collective

351:    Input Parameters:
352: +   dmmg - DMMG solve context
353: -   level - the number of the level you want the context for

355:    Level: intermediate

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

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

362: /*MC
363:    DMMGSetUser - Sets the user context for a particular level

365:    Synopsis:
366:    PetscErrorCode DMMGSetUser(DMMG *dmmg,PetscInt level,void *ctx)

368:    Not Collective

370:    Input Parameters:
371: +   dmmg - DMMG solve context
372: .   level - the number of the level you want the context for
373: -   ctx - the context

375:    Level: intermediate

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

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

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

385: /*MC
386:    DMMGGetLevels - Gets the number of levels in a DMMG object

388:    Synopsis:
389:    PetscInt DMMGGetLevels(DMMG *dmmg)

391:    Not Collective

393:    Input Parameter:
394: .   dmmg - DMMG solve context

396:    Level: intermediate

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

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

403: /*MC
404:    DMMGGetDMMG - Returns the DMMG struct for the finest level

406:    Synopsis:
407:    DMMG DMMGGetDMMG(DMMG *dmmg)

409:    Not Collective

411:    Input Parameter:
412: .   dmmg - DMMG solve context

414:    Level: intermediate

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

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

421: #define PCDMMG      "pcdmmg"
422: EXTERN PetscErrorCode PCDMMGSetDMMG(PC,DMMG*);

425: #endif