Actual source code: mffd.c

petsc-3.7.7 2017-09-25
Report Typos and Errors
  2: #include <petsc/private/matimpl.h>
  3: #include <../src/mat/impls/mffd/mffdimpl.h>   /*I  "petscmat.h"   I*/

  5: PetscFunctionList MatMFFDList              = 0;
  6: PetscBool         MatMFFDRegisterAllCalled = PETSC_FALSE;

  8: PetscClassId  MATMFFD_CLASSID;
  9: PetscLogEvent MATMFFD_Mult;

 11: static PetscBool MatMFFDPackageInitialized = PETSC_FALSE;
 14: /*@C
 15:   MatMFFDFinalizePackage - This function destroys everything in the MatMFFD package. It is
 16:   called from PetscFinalize().

 18:   Level: developer

 20: .keywords: Petsc, destroy, package
 21: .seealso: PetscFinalize(), MatCreateMFFD(), MatCreateSNESMF()
 22: @*/
 23: PetscErrorCode  MatMFFDFinalizePackage(void)
 24: {

 28:   PetscFunctionListDestroy(&MatMFFDList);
 29:   MatMFFDPackageInitialized = PETSC_FALSE;
 30:   MatMFFDRegisterAllCalled  = PETSC_FALSE;
 31:   return(0);
 32: }

 36: /*@C
 37:   MatMFFDInitializePackage - This function initializes everything in the MatMFFD package. It is called
 38:   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to MatCreate_MFFD()
 39:   when using static libraries.

 41:   Level: developer

 43: .keywords: Vec, initialize, package
 44: .seealso: PetscInitialize()
 45: @*/
 46: PetscErrorCode  MatMFFDInitializePackage(void)
 47: {
 48:   char           logList[256];
 49:   char           *className;
 50:   PetscBool      opt;

 54:   if (MatMFFDPackageInitialized) return(0);
 55:   MatMFFDPackageInitialized = PETSC_TRUE;
 56:   /* Register Classes */
 57:   PetscClassIdRegister("MatMFFD",&MATMFFD_CLASSID);
 58:   /* Register Constructors */
 59:   MatMFFDRegisterAll();
 60:   /* Register Events */
 61:   PetscLogEventRegister("MatMult MF",          MATMFFD_CLASSID,&MATMFFD_Mult);

 63:   /* Process info exclusions */
 64:   PetscOptionsGetString(NULL,NULL, "-info_exclude", logList, 256, &opt);
 65:   if (opt) {
 66:     PetscStrstr(logList, "matmffd", &className);
 67:     if (className) {
 68:       PetscInfoDeactivateClass(MATMFFD_CLASSID);
 69:     }
 70:   }
 71:   /* Process summary exclusions */
 72:   PetscOptionsGetString(NULL,NULL, "-log_exclude", logList, 256, &opt);
 73:   if (opt) {
 74:     PetscStrstr(logList, "matmffd", &className);
 75:     if (className) {
 76:       PetscLogEventDeactivateClass(MATMFFD_CLASSID);
 77:     }
 78:   }
 79:   PetscRegisterFinalize(MatMFFDFinalizePackage);
 80:   return(0);
 81: }

 85: /*@C
 86:     MatMFFDSetType - Sets the method that is used to compute the
 87:     differencing parameter for finite differene matrix-free formulations.

 89:     Input Parameters:
 90: +   mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMFFD()
 91:           or MatSetType(mat,MATMFFD);
 92: -   ftype - the type requested, either MATMFFD_WP or MATMFFD_DS

 94:     Level: advanced

 96:     Notes:
 97:     For example, such routines can compute h for use in
 98:     Jacobian-vector products of the form

100:                         F(x+ha) - F(x)
101:           F'(u)a  ~=  ----------------
102:                               h

104: .seealso: MatCreateSNESMF(), MatMFFDRegister(), MatMFFDSetFunction(), MatCreateMFFD()
105: @*/
106: PetscErrorCode  MatMFFDSetType(Mat mat,MatMFFDType ftype)
107: {
108:   PetscErrorCode ierr,(*r)(MatMFFD);
109:   MatMFFD        ctx = (MatMFFD)mat->data;
110:   PetscBool      match;


116:   PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);
117:   if (!match) return(0);

119:   /* already set, so just return */
120:   PetscObjectTypeCompare((PetscObject)ctx,ftype,&match);
121:   if (match) return(0);

123:   /* destroy the old one if it exists */
124:   if (ctx->ops->destroy) {
125:     (*ctx->ops->destroy)(ctx);
126:   }

128:    PetscFunctionListFind(MatMFFDList,ftype,&r);
129:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown MatMFFD type %s given",ftype);
130:   (*r)(ctx);
131:   PetscObjectChangeTypeName((PetscObject)ctx,ftype);
132:   return(0);
133: }

135: typedef PetscErrorCode (*FCN1)(void*,Vec); /* force argument to next function to not be extern C*/
138: PetscErrorCode  MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func)
139: {
140:   MatMFFD ctx = (MatMFFD)mat->data;

143:   ctx->funcisetbase = func;
144:   return(0);
145: }

147: typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
150: PetscErrorCode  MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci)
151: {
152:   MatMFFD ctx = (MatMFFD)mat->data;

155:   ctx->funci = funci;
156:   return(0);
157: }

161: PetscErrorCode  MatMFFDResetHHistory_MFFD(Mat J)
162: {
163:   MatMFFD ctx = (MatMFFD)J->data;

166:   ctx->ncurrenth = 0;
167:   return(0);
168: }

172: /*@C
173:    MatMFFDRegister - Adds a method to the MatMFFD registry.

175:    Not Collective

177:    Input Parameters:
178: +  name_solver - name of a new user-defined compute-h module
179: -  routine_create - routine to create method context

181:    Level: developer

183:    Notes:
184:    MatMFFDRegister() may be called multiple times to add several user-defined solvers.

186:    Sample usage:
187: .vb
188:    MatMFFDRegister("my_h",MyHCreate);
189: .ve

191:    Then, your solver can be chosen with the procedural interface via
192: $     MatMFFDSetType(mfctx,"my_h")
193:    or at runtime via the option
194: $     -mat_mffd_type my_h

196: .keywords: MatMFFD, register

198: .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
199:  @*/
200: PetscErrorCode  MatMFFDRegister(const char sname[],PetscErrorCode (*function)(MatMFFD))
201: {

205:   PetscFunctionListAdd(&MatMFFDList,sname,function);
206:   return(0);
207: }

209: /* ----------------------------------------------------------------------------------------*/
212: PetscErrorCode MatDestroy_MFFD(Mat mat)
213: {
215:   MatMFFD        ctx = (MatMFFD)mat->data;

218:   VecDestroy(&ctx->w);
219:   VecDestroy(&ctx->drscale);
220:   VecDestroy(&ctx->dlscale);
221:   VecDestroy(&ctx->dshift);
222:   if (ctx->current_f_allocated) {
223:     VecDestroy(&ctx->current_f);
224:   }
225:   if (ctx->ops->destroy) {(*ctx->ops->destroy)(ctx);}
226:   PetscHeaderDestroy(&ctx);
227:   mat->data = 0;

229:   PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C",NULL);
230:   PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C",NULL);
231:   PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C",NULL);
232:   PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunction_C",NULL);
233:   PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctionError_C",NULL);
234:   PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C",NULL);
235:   PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetPeriod_C",NULL);
236:   PetscObjectComposeFunction((PetscObject)mat,"MatMFFDResetHHistory_C",NULL);
237:   return(0);
238: }

242: /*
243:    MatMFFDView_MFFD - Views matrix-free parameters.

245: */
246: PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer)
247: {
249:   MatMFFD        ctx = (MatMFFD)J->data;
250:   PetscBool      iascii, viewbase, viewfunction;
251:   const char     *prefix;

254:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
255:   if (iascii) {
256:     PetscViewerASCIIPrintf(viewer,"Matrix-free approximation:\n");
257:     PetscViewerASCIIPushTab(viewer);
258:     PetscViewerASCIIPrintf(viewer,"err=%g (relative error in function evaluation)\n",(double)ctx->error_rel);
259:     if (!((PetscObject)ctx)->type_name) {
260:       PetscViewerASCIIPrintf(viewer,"The compute h routine has not yet been set\n");
261:     } else {
262:       PetscViewerASCIIPrintf(viewer,"Using %s compute h routine\n",((PetscObject)ctx)->type_name);
263:     }
264:     if (ctx->ops->view) {
265:       (*ctx->ops->view)(ctx,viewer);
266:     }
267:     PetscObjectGetOptionsPrefix((PetscObject)J, &prefix);

269:     PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_base", &viewbase);
270:     if (viewbase) {
271:       PetscViewerASCIIPrintf(viewer, "Base:\n");
272:       VecView(ctx->current_u, viewer);
273:     }
274:     PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_function", &viewfunction);
275:     if (viewfunction) {
276:       PetscViewerASCIIPrintf(viewer, "Function:\n");
277:       VecView(ctx->current_f, viewer);
278:     }
279:     PetscViewerASCIIPopTab(viewer);
280:   }
281:   return(0);
282: }

286: /*
287:    MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This
288:    allows the user to indicate the beginning of a new linear solve by calling
289:    MatAssemblyXXX() on the matrix free matrix. This then allows the
290:    MatCreateMFFD_WP() to properly compute ||U|| only the first time
291:    in the linear solver rather than every time.

293:    This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library hence
294:    it must be labeled as PETSC_EXTERN
295: */
296: PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt)
297: {
299:   MatMFFD        j = (MatMFFD)J->data;

302:   MatMFFDResetHHistory(J);
303:   j->vshift = 0.0;
304:   j->vscale = 1.0;
305:   return(0);
306: }

310: /*
311:   MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:

313:         y ~= (F(u + ha) - F(u))/h,
314:   where F = nonlinear function, as set by SNESSetFunction()
315:         u = current iterate
316:         h = difference interval
317: */
318: PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y)
319: {
320:   MatMFFD        ctx = (MatMFFD)mat->data;
321:   PetscScalar    h;
322:   Vec            w,U,F;
324:   PetscBool      zeroa;

327:   if (!ctx->current_u) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MatMFFDSetBase() has not been called, this is often caused by forgetting to call \n\t\tMatAssemblyBegin/End on the first Mat in the SNES compute function");
328:   /* We log matrix-free matrix-vector products separately, so that we can
329:      separate the performance monitoring from the cases that use conventional
330:      storage.  We may eventually modify event logging to associate events
331:      with particular objects, hence alleviating the more general problem. */
332:   PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);

334:   w = ctx->w;
335:   U = ctx->current_u;
336:   F = ctx->current_f;
337:   /*
338:       Compute differencing parameter
339:   */
340:   if (!ctx->ops->compute) {
341:     MatMFFDSetType(mat,MATMFFD_WP);
342:     MatSetFromOptions(mat);
343:   }
344:   (*ctx->ops->compute)(ctx,U,a,&h,&zeroa);
345:   if (zeroa) {
346:     VecSet(y,0.0);
347:     return(0);
348:   }

350:   if (mat->erroriffailure && PetscIsInfOrNanScalar(h)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed Nan differencing parameter h");
351:   if (ctx->checkh) {
352:     (*ctx->checkh)(ctx->checkhctx,U,a,&h);
353:   }

355:   /* keep a record of the current differencing parameter h */
356:   ctx->currenth = h;
357: #if defined(PETSC_USE_COMPLEX)
358:   PetscInfo2(mat,"Current differencing parameter: %g + %g i\n",(double)PetscRealPart(h),(double)PetscImaginaryPart(h));
359: #else
360:   PetscInfo1(mat,"Current differencing parameter: %15.12e\n",h);
361: #endif
362:   if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) {
363:     ctx->historyh[ctx->ncurrenth] = h;
364:   }
365:   ctx->ncurrenth++;

367:   /* w = u + ha */
368:   if (ctx->drscale) {
369:     VecPointwiseMult(ctx->drscale,a,U);
370:     VecAYPX(U,h,w);
371:   } else {
372:     VecWAXPY(w,h,a,U);
373:   }

375:   /* compute func(U) as base for differencing; only needed first time in and not when provided by user */
376:   if (ctx->ncurrenth == 1 && ctx->current_f_allocated) {
377:     (*ctx->func)(ctx->funcctx,U,F);
378:   }
379:   (*ctx->func)(ctx->funcctx,w,y);

381:   VecAXPY(y,-1.0,F);
382:   VecScale(y,1.0/h);

384:   if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
385:     VecAXPBY(y,ctx->vshift,ctx->vscale,a);
386:   }
387:   if (ctx->dlscale) {
388:     VecPointwiseMult(y,ctx->dlscale,y);
389:   }
390:   if (ctx->dshift) {
391:     VecPointwiseMult(ctx->dshift,a,U);
392:     VecAXPY(y,1.0,U);
393:   }

395:   if (mat->nullsp) {MatNullSpaceRemove(mat->nullsp,y);}

397:   PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);
398:   return(0);
399: }

403: /*
404:   MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix

406:         y ~= (F(u + ha) - F(u))/h,
407:   where F = nonlinear function, as set by SNESSetFunction()
408:         u = current iterate
409:         h = difference interval
410: */
411: PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a)
412: {
413:   MatMFFD        ctx = (MatMFFD)mat->data;
414:   PetscScalar    h,*aa,*ww,v;
415:   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
416:   Vec            w,U;
418:   PetscInt       i,rstart,rend;

421:   if (!ctx->funci) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Requires calling MatMFFDSetFunctioni() first");

423:   w    = ctx->w;
424:   U    = ctx->current_u;
425:   (*ctx->func)(ctx->funcctx,U,a);
426:   (*ctx->funcisetbase)(ctx->funcctx,U);
427:   VecCopy(U,w);

429:   VecGetOwnershipRange(a,&rstart,&rend);
430:   VecGetArray(a,&aa);
431:   for (i=rstart; i<rend; i++) {
432:     VecGetArray(w,&ww);
433:     h    = ww[i-rstart];
434:     if (h == 0.0) h = 1.0;
435:     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
436:     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
437:     h *= epsilon;

439:     ww[i-rstart] += h;
440:     VecRestoreArray(w,&ww);
441:     (*ctx->funci)(ctx->funcctx,i,w,&v);
442:     aa[i-rstart]  = (v - aa[i-rstart])/h;

444:     /* possibly shift and scale result */
445:     if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
446:       aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart];
447:     }

449:     VecGetArray(w,&ww);
450:     ww[i-rstart] -= h;
451:     VecRestoreArray(w,&ww);
452:   }
453:   VecRestoreArray(a,&aa);
454:   return(0);
455: }

459: PetscErrorCode MatDiagonalScale_MFFD(Mat mat,Vec ll,Vec rr)
460: {
461:   MatMFFD        aij = (MatMFFD)mat->data;

465:   if (ll && !aij->dlscale) {
466:     VecDuplicate(ll,&aij->dlscale);
467:   }
468:   if (rr && !aij->drscale) {
469:     VecDuplicate(rr,&aij->drscale);
470:   }
471:   if (ll) {
472:     VecCopy(ll,aij->dlscale);
473:   }
474:   if (rr) {
475:     VecCopy(rr,aij->drscale);
476:   }
477:   return(0);
478: }

482: PetscErrorCode MatDiagonalSet_MFFD(Mat mat,Vec ll,InsertMode mode)
483: {
484:   MatMFFD        aij = (MatMFFD)mat->data;

488:   if (mode == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"No diagonal set with INSERT_VALUES");
489:   if (!aij->dshift) {
490:     VecDuplicate(ll,&aij->dshift);
491:   }
492:   VecAXPY(aij->dshift,1.0,ll);
493:   return(0);
494: }

498: PetscErrorCode MatShift_MFFD(Mat Y,PetscScalar a)
499: {
500:   MatMFFD shell = (MatMFFD)Y->data;

503:   shell->vshift += a;
504:   return(0);
505: }

509: PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a)
510: {
511:   MatMFFD shell = (MatMFFD)Y->data;

514:   shell->vscale *= a;
515:   return(0);
516: }

520: PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F)
521: {
523:   MatMFFD        ctx = (MatMFFD)J->data;

526:   MatMFFDResetHHistory(J);

528:   ctx->current_u = U;
529:   if (F) {
530:     if (ctx->current_f_allocated) {VecDestroy(&ctx->current_f);}
531:     ctx->current_f           = F;
532:     ctx->current_f_allocated = PETSC_FALSE;
533:   } else if (!ctx->current_f_allocated) {
534:     MatCreateVecs(J,NULL,&ctx->current_f);

536:     ctx->current_f_allocated = PETSC_TRUE;
537:   }
538:   if (!ctx->w) {
539:     VecDuplicate(ctx->current_u, &ctx->w);
540:   }
541:   J->assembled = PETSC_TRUE;
542:   return(0);
543: }

545: typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/

549: PetscErrorCode  MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void *ectx)
550: {
551:   MatMFFD ctx = (MatMFFD)J->data;

554:   ctx->checkh    = fun;
555:   ctx->checkhctx = ectx;
556:   return(0);
557: }

561: /*@C
562:    MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all
563:    MatMFFD options in the database.

565:    Collective on Mat

567:    Input Parameter:
568: +  A - the Mat context
569: -  prefix - the prefix to prepend to all option names

571:    Notes:
572:    A hyphen (-) must NOT be given at the beginning of the prefix name.
573:    The first character of all runtime options is AUTOMATICALLY the hyphen.

575:    Level: advanced

577: .keywords: SNES, matrix-free, parameters

579: .seealso: MatSetFromOptions(), MatCreateSNESMF(), MatCreateMFFD()
580: @*/
581: PetscErrorCode  MatMFFDSetOptionsPrefix(Mat mat,const char prefix[])

583: {
584:   MatMFFD        mfctx = mat ? (MatMFFD)mat->data : (MatMFFD)NULL;

590:   PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix);
591:   return(0);
592: }

596: PetscErrorCode  MatSetFromOptions_MFFD(PetscOptionItems *PetscOptionsObject,Mat mat)
597: {
598:   MatMFFD        mfctx = (MatMFFD)mat->data;
600:   PetscBool      flg;
601:   char           ftype[256];

606:   PetscObjectOptionsBegin((PetscObject)mfctx);
607:   PetscOptionsFList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg);
608:   if (flg) {
609:     MatMFFDSetType(mat,ftype);
610:   }

612:   PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);
613:   PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);

615:   flg  = PETSC_FALSE;
616:   PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,NULL);
617:   if (flg) {
618:     MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);
619:   }
620:   if (mfctx->ops->setfromoptions) {
621:     (*mfctx->ops->setfromoptions)(PetscOptionsObject,mfctx);
622:   }
623:   PetscOptionsEnd();
624:   return(0);
625: }

629: PetscErrorCode  MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period)
630: {
631:   MatMFFD ctx = (MatMFFD)mat->data;

635:   ctx->recomputeperiod = period;
636:   return(0);
637: }

641: PetscErrorCode  MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
642: {
643:   MatMFFD ctx = (MatMFFD)mat->data;

646:   ctx->func    = func;
647:   ctx->funcctx = funcctx;
648:   return(0);
649: }

653: PetscErrorCode  MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error)
654: {
655:   MatMFFD ctx = (MatMFFD)mat->data;

659:   if (error != PETSC_DEFAULT) ctx->error_rel = error;
660:   return(0);
661: }

665: static PetscErrorCode MatMissingDiagonal_MFFD(Mat A,PetscBool  *missing,PetscInt *d)
666: {
668:   *missing = PETSC_FALSE;
669:   return(0);
670: }

672: /*MC
673:   MATMFFD - MATMFFD = "mffd" - A matrix free matrix type.

675:   Level: advanced

677: .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction(), MatMFFDSetType(),  
678:           MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction()
679:           MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(),
680:           MatMFFDGetH(),
681: M*/
684: PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A)
685: {
686:   MatMFFD        mfctx;

690:   MatMFFDInitializePackage();

692:   PetscHeaderCreate(mfctx,MATMFFD_CLASSID,"MatMFFD","Matrix-free Finite Differencing","Mat",PetscObjectComm((PetscObject)A),MatDestroy_MFFD,MatView_MFFD);

694:   mfctx->error_rel                = PETSC_SQRT_MACHINE_EPSILON;
695:   mfctx->recomputeperiod          = 1;
696:   mfctx->count                    = 0;
697:   mfctx->currenth                 = 0.0;
698:   mfctx->historyh                 = NULL;
699:   mfctx->ncurrenth                = 0;
700:   mfctx->maxcurrenth              = 0;
701:   ((PetscObject)mfctx)->type_name = 0;

703:   mfctx->vshift = 0.0;
704:   mfctx->vscale = 1.0;

706:   /*
707:      Create the empty data structure to contain compute-h routines.
708:      These will be filled in below from the command line options or
709:      a later call with MatMFFDSetType() or if that is not called
710:      then it will default in the first use of MatMult_MFFD()
711:   */
712:   mfctx->ops->compute        = 0;
713:   mfctx->ops->destroy        = 0;
714:   mfctx->ops->view           = 0;
715:   mfctx->ops->setfromoptions = 0;
716:   mfctx->hctx                = 0;

718:   mfctx->func    = 0;
719:   mfctx->funcctx = 0;
720:   mfctx->w       = NULL;

722:   A->data = mfctx;

724:   A->ops->mult            = MatMult_MFFD;
725:   A->ops->destroy         = MatDestroy_MFFD;
726:   A->ops->view            = MatView_MFFD;
727:   A->ops->assemblyend     = MatAssemblyEnd_MFFD;
728:   A->ops->getdiagonal     = MatGetDiagonal_MFFD;
729:   A->ops->scale           = MatScale_MFFD;
730:   A->ops->shift           = MatShift_MFFD;
731:   A->ops->diagonalscale   = MatDiagonalScale_MFFD;
732:   A->ops->diagonalset     = MatDiagonalSet_MFFD;
733:   A->ops->setfromoptions  = MatSetFromOptions_MFFD;
734:   A->ops->missingdiagonal = MatMissingDiagonal_MFFD;
735:   A->assembled           = PETSC_TRUE;

737:   PetscLayoutSetUp(A->rmap);
738:   PetscLayoutSetUp(A->cmap);

740:   PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetBase_C",MatMFFDSetBase_MFFD);
741:   PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioniBase_C",MatMFFDSetFunctioniBase_MFFD);
742:   PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioni_C",MatMFFDSetFunctioni_MFFD);
743:   PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunction_C",MatMFFDSetFunction_MFFD);
744:   PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetCheckh_C",MatMFFDSetCheckh_MFFD);
745:   PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetPeriod_C",MatMFFDSetPeriod_MFFD);
746:   PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctionError_C",MatMFFDSetFunctionError_MFFD);
747:   PetscObjectComposeFunction((PetscObject)A,"MatMFFDResetHHistory_C",MatMFFDResetHHistory_MFFD);

749:   mfctx->mat = A;

751:   PetscObjectChangeTypeName((PetscObject)A,MATMFFD);
752:   return(0);
753: }

757: /*@
758:    MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF()

760:    Collective on Vec

762:    Input Parameters:
763: +  comm - MPI communicator
764: .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
765:            This value should be the same as the local size used in creating the
766:            y vector for the matrix-vector product y = Ax.
767: .  n - This value should be the same as the local size used in creating the
768:        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
769:        calculated if N is given) For square matrices n is almost always m.
770: .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
771: -  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)


774:    Output Parameter:
775: .  J - the matrix-free matrix

777:    Options Database Keys: call MatSetFromOptions() to trigger these
778: +  -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS)
779: -  -mat_mffd_err - square root of estimated relative error in function evaluation
780: -  -mat_mffd_period - how often h is recomputed, defaults to 1, everytime


783:    Level: advanced

785:    Notes:
786:    The matrix-free matrix context merely contains the function pointers
787:    and work space for performing finite difference approximations of
788:    Jacobian-vector products, F'(u)*a,

790:    The default code uses the following approach to compute h

792: .vb
793:      F'(u)*a = [F(u+h*a) - F(u)]/h where
794:      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
795:        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
796:  where
797:      error_rel = square root of relative error in function evaluation
798:      umin = minimum iterate parameter
799: .ve

801:    You can call SNESSetJacobian() with MatMFFDComputeJacobian() if you are using matrix and not a different
802:    preconditioner matrix

804:    The user can set the error_rel via MatMFFDSetFunctionError() and
805:    umin via MatMFFDDSSetUmin(); see Users-Manual: ch_snes for details.

807:    The user should call MatDestroy() when finished with the matrix-free
808:    matrix context.

810:    Options Database Keys:
811: +  -mat_mffd_err <error_rel> - Sets error_rel
812: .  -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only)
813: -  -mat_mffd_check_positivity

815: .keywords: default, matrix-free, create, matrix

817: .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction()
818:           MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(),
819:           MatMFFDGetH(), MatMFFDRegister(), MatMFFDComputeJacobian()

821: @*/
822: PetscErrorCode  MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J)
823: {

827:   MatCreate(comm,J);
828:   MatSetSizes(*J,m,n,M,N);
829:   MatSetType(*J,MATMFFD);
830:   MatSetUp(*J);
831:   return(0);
832: }


837: /*@
838:    MatMFFDGetH - Gets the last value that was used as the differencing
839:    parameter.

841:    Not Collective

843:    Input Parameters:
844: .  mat - the matrix obtained with MatCreateSNESMF()

846:    Output Paramter:
847: .  h - the differencing step size

849:    Level: advanced

851: .keywords: SNES, matrix-free, parameters

853: .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory()
854: @*/
855: PetscErrorCode  MatMFFDGetH(Mat mat,PetscScalar *h)
856: {
857:   MatMFFD        ctx = (MatMFFD)mat->data;
859:   PetscBool      match;

862:   PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);
863:   if (!match) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");

865:   *h = ctx->currenth;
866:   return(0);
867: }

871: /*@C
872:    MatMFFDSetFunction - Sets the function used in applying the matrix free.

874:    Logically Collective on Mat

876:    Input Parameters:
877: +  mat - the matrix free matrix created via MatCreateSNESMF() or MatCreateMFFD()
878: .  func - the function to use
879: -  funcctx - optional function context passed to function

881:    Calling Sequence of func:
882: $     func (void *funcctx, Vec x, Vec f)

884: +  funcctx - user provided context
885: .  x - input vector
886: -  f - computed output function

888:    Level: advanced

890:    Notes:
891:     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
892:     matrix inside your compute Jacobian routine

894:     If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used.

896: .keywords: SNES, matrix-free, function

898: .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD,
899:           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
900: @*/
901: PetscErrorCode  MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
902: {

906:   PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx));
907:   return(0);
908: }

912: /*@C
913:    MatMFFDSetFunctioni - Sets the function for a single component

915:    Logically Collective on Mat

917:    Input Parameters:
918: +  mat - the matrix free matrix created via MatCreateSNESMF()
919: -  funci - the function to use

921:    Level: advanced

923:    Notes:
924:     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
925:     matrix inside your compute Jacobian routine


928: .keywords: SNES, matrix-free, function

930: .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()

932: @*/
933: PetscErrorCode  MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*))
934: {

939:   PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));
940:   return(0);
941: }


946: /*@C
947:    MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation

949:    Logically Collective on Mat

951:    Input Parameters:
952: +  mat - the matrix free matrix created via MatCreateSNESMF()
953: -  func - the function to use

955:    Level: advanced

957:    Notes:
958:     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
959:     matrix inside your compute Jacobian routine


962: .keywords: SNES, matrix-free, function

964: .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
965:           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
966: @*/
967: PetscErrorCode  MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec))
968: {

973:   PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));
974:   return(0);
975: }

979: /*@
980:    MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime

982:    Logically Collective on Mat

984:    Input Parameters:
985: +  mat - the matrix free matrix created via MatCreateSNESMF()
986: -  period - 1 for everytime, 2 for every second etc

988:    Options Database Keys:
989: +  -mat_mffd_period <period>

991:    Level: advanced


994: .keywords: SNES, matrix-free, parameters

996: .seealso: MatCreateSNESMF(),MatMFFDGetH(),
997:           MatMFFDSetHHistory(), MatMFFDResetHHistory()
998: @*/
999: PetscErrorCode  MatMFFDSetPeriod(Mat mat,PetscInt period)
1000: {

1004:   PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period));
1005:   return(0);
1006: }

1010: /*@
1011:    MatMFFDSetFunctionError - Sets the error_rel for the approximation of
1012:    matrix-vector products using finite differences.

1014:    Logically Collective on Mat

1016:    Input Parameters:
1017: +  mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF()
1018: -  error_rel - relative error (should be set to the square root of
1019:                the relative error in the function evaluations)

1021:    Options Database Keys:
1022: +  -mat_mffd_err <error_rel> - Sets error_rel

1024:    Level: advanced

1026:    Notes:
1027:    The default matrix-free matrix-vector product routine computes
1028: .vb
1029:      F'(u)*a = [F(u+h*a) - F(u)]/h where
1030:      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
1031:        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
1032: .ve

1034: .keywords: SNES, matrix-free, parameters

1036: .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
1037:           MatMFFDSetHHistory(), MatMFFDResetHHistory()
1038: @*/
1039: PetscErrorCode  MatMFFDSetFunctionError(Mat mat,PetscReal error)
1040: {

1044:   PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error));
1045:   return(0);
1046: }

1050: /*@
1051:    MatMFFDSetHHistory - Sets an array to collect a history of the
1052:    differencing values (h) computed for the matrix-free product.

1054:    Logically Collective on Mat

1056:    Input Parameters:
1057: +  J - the matrix-free matrix context
1058: .  histroy - space to hold the history
1059: -  nhistory - number of entries in history, if more entries are generated than
1060:               nhistory, then the later ones are discarded

1062:    Level: advanced

1064:    Notes:
1065:    Use MatMFFDResetHHistory() to reset the history counter and collect
1066:    a new batch of differencing parameters, h.

1068: .keywords: SNES, matrix-free, h history, differencing history

1070: .seealso: MatMFFDGetH(), MatCreateSNESMF(),
1071:           MatMFFDResetHHistory(), MatMFFDSetFunctionError()

1073: @*/
1074: PetscErrorCode  MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory)
1075: {
1076:   MatMFFD        ctx = (MatMFFD)J->data;
1078:   PetscBool      match;

1081:   PetscObjectTypeCompare((PetscObject)J,MATMFFD,&match);
1082:   if (!match) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
1083:   ctx->historyh    = history;
1084:   ctx->maxcurrenth = nhistory;
1085:   ctx->currenth    = 0.;
1086:   return(0);
1087: }


1092: /*@
1093:    MatMFFDResetHHistory - Resets the counter to zero to begin
1094:    collecting a new set of differencing histories.

1096:    Logically Collective on Mat

1098:    Input Parameters:
1099: .  J - the matrix-free matrix context

1101:    Level: advanced

1103:    Notes:
1104:    Use MatMFFDSetHHistory() to create the original history counter.

1106: .keywords: SNES, matrix-free, h history, differencing history

1108: .seealso: MatMFFDGetH(), MatCreateSNESMF(),
1109:           MatMFFDSetHHistory(), MatMFFDSetFunctionError()

1111: @*/
1112: PetscErrorCode  MatMFFDResetHHistory(Mat J)
1113: {

1117:   PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J));
1118:   return(0);
1119: }


1124: /*@
1125:     MatMFFDSetBase - Sets the vector U at which matrix vector products of the
1126:         Jacobian are computed

1128:     Logically Collective on Mat

1130:     Input Parameters:
1131: +   J - the MatMFFD matrix
1132: .   U - the vector
1133: -   F - (optional) vector that contains F(u) if it has been already computed

1135:     Notes: This is rarely used directly

1137:     If F is provided then it is not recomputed. Otherwise the function is evaluated at the base
1138:     point during the first MatMult() after each call to MatMFFDSetBase().

1140:     Level: advanced

1142: @*/
1143: PetscErrorCode  MatMFFDSetBase(Mat J,Vec U,Vec F)
1144: {

1151:   PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));
1152:   return(0);
1153: }

1157: /*@C
1158:     MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts
1159:         it to satisfy some criteria

1161:     Logically Collective on Mat

1163:     Input Parameters:
1164: +   J - the MatMFFD matrix
1165: .   fun - the function that checks h
1166: -   ctx - any context needed by the function

1168:     Options Database Keys:
1169: .   -mat_mffd_check_positivity

1171:     Level: advanced

1173:     Notes: For example, MatMFFDCheckPositivity() insures that all entries
1174:        of U + h*a are non-negative

1176:      The function you provide is called after the default h has been computed and allows you to
1177:      modify it.

1179: .seealso:  MatMFFDCheckPositivity()
1180: @*/
1181: PetscErrorCode  MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void *ctx)
1182: {

1187:   PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));
1188:   return(0);
1189: }

1193: /*@
1194:     MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or
1195:         zero, decreases h until this is satisfied.

1197:     Logically Collective on Vec

1199:     Input Parameters:
1200: +   U - base vector that is added to
1201: .   a - vector that is added
1202: .   h - scaling factor on a
1203: -   dummy - context variable (unused)

1205:     Options Database Keys:
1206: .   -mat_mffd_check_positivity

1208:     Level: advanced

1210:     Notes: This is rarely used directly, rather it is passed as an argument to
1211:            MatMFFDSetCheckh()

1213: .seealso:  MatMFFDSetCheckh()
1214: @*/
1215: PetscErrorCode  MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h)
1216: {
1217:   PetscReal      val, minval;
1218:   PetscScalar    *u_vec, *a_vec;
1220:   PetscInt       i,n;
1221:   MPI_Comm       comm;

1224:   PetscObjectGetComm((PetscObject)U,&comm);
1225:   VecGetArray(U,&u_vec);
1226:   VecGetArray(a,&a_vec);
1227:   VecGetLocalSize(U,&n);
1228:   minval = PetscAbsScalar(*h*1.01);
1229:   for (i=0; i<n; i++) {
1230:     if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) {
1231:       val = PetscAbsScalar(u_vec[i]/a_vec[i]);
1232:       if (val < minval) minval = val;
1233:     }
1234:   }
1235:   VecRestoreArray(U,&u_vec);
1236:   VecRestoreArray(a,&a_vec);
1237:   MPIU_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);
1238:   if (val <= PetscAbsScalar(*h)) {
1239:     PetscInfo2(U,"Scaling back h from %g to %g\n",(double)PetscRealPart(*h),(double)(.99*val));
1240:     if (PetscRealPart(*h) > 0.0) *h =  0.99*val;
1241:     else                         *h = -0.99*val;
1242:   }
1243:   return(0);
1244: }