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(), 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(),  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*);

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

 62:   PetscTruth     getcoloringfrommat;     /* call a graph coloring algorithm on the matrix to get the coloring, instead of getting it from the DM */
 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*));
 86: EXTERN PetscErrorCode  DMMGSetFromOptions(DMMG*);

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

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

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

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


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

130:    Synopsis:
131:    Vec DMMGGetRHS(DMMG *dmmg)

133:    Not Collective, but resulting vector is parallel

135:    Input Parameters:
136: .   dmmg - DMMG solve context

138:    Level: intermediate

140:    Fortran Usage:
141: .     DMMGGetRHS(DMMG dmmg,Vec b,PetscErrorCode ierr)

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

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

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

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

153:    Synopsis:
154:    Vec DMMGGetx(DMMG *dmmg)

156:    Not Collective, but resulting vector is parallel

158:    Input Parameters:
159: .   dmmg - DMMG solve context

161:    Level: intermediate

163:    Fortran Usage:
164: .     DMMGGetx(DMMG dmmg,Vec x,PetscErrorCode ierr)

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

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

171: /*MC
172:    DMMGGetJ - Returns the Jacobian (matrix) for the finest level

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(), DMMGGetB(), DMMGGetRHS()

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

189: /*MC
190:    DMMGGetComm - Returns the MPI_Comm for the finest level

192:    Synopsis:
193:    MPI_Comm DMMGGetJ(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 DMMGGetComm(ctx)           (ctx)[(ctx)[0]->nlevels-1]->comm

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

211:    Synopsis:
212:    Mat DMMGGetJ(DMMG *dmmg)

214:    Not Collective

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

219:    Level: intermediate

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

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

226: /*MC
227:    DMMGGetFine - Returns the DMMG associated with the finest level

229:    Synopsis:
230:    DMMG DMMGGetFine(DMMG *dmmg)

232:    Not Collective

234:    Input Parameter:
235: .   dmmg - DMMG solve context

237:    Level: intermediate

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

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


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

248:    Synopsis:
249:    KSP DMMGGetKSP(DMMG *dmmg)

251:    Not Collective

253:    Input Parameter:
254: .   dmmg - DMMG solve context

256:    Level: intermediate

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

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

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

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

270:    Synopsis:
271:    SNES DMMGGetSNES(DMMG *dmmg)

273:    Not Collective

275:    Input Parameter:
276: .   dmmg - DMMG solve context

278:    Level: intermediate

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

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

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

287: /*MC
288:    DMMGGetDM - Gets the DM object on the finest level

290:    Synopsis:
291:    DM DMMGGetDM(DMMG *dmmg)

293:    Not Collective

295:    Input Parameter:
296: .   dmmg - DMMG solve context

298:    Level: intermediate

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

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

305: /*MC
306:    DMMGGetDA - Gets the DA object on the finest level

308:    Synopsis:
309:    DA DMMGGetDA(DMMG *dmmg)

311:    Not Collective

313:    Input Parameter:
314: .   dmmg - DMMG solve context

316:    Level: intermediate

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

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

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

325: /*MC
326:    DMMGGetDMComposite - Gets the DMComposite object on the finest level

328:    Synopsis:
329:    DMComposite DMMGGetDMComposite(DMMG *dmmg)

331:    Not Collective

333:    Input Parameter:
334: .   dmmg - DMMG solve context

336:    Level: intermediate

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

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

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

345: /*MC
346:    DMMGGetUser - Returns the user context for a particular level

348:    Synopsis:
349:    void* DMMGGetUser(DMMG *dmmg,PetscInt level)

351:    Not Collective

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

357:    Level: intermediate

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

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

364: /*MC
365:    DMMGSetUser - Sets the user context for a particular level

367:    Synopsis:
368:    PetscErrorCode DMMGSetUser(DMMG *dmmg,PetscInt level,void *ctx)

370:    Not Collective

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

377:    Level: intermediate

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

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

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

387: /*MC
388:    DMMGGetLevels - Gets the number of levels in a DMMG object

390:    Synopsis:
391:    PetscInt DMMGGetLevels(DMMG *dmmg)

393:    Not Collective

395:    Input Parameter:
396: .   dmmg - DMMG solve context

398:    Level: intermediate

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

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

406: #endif