Actual source code: mffd.c
2: #include <private/matimpl.h>
3: #include <../src/mat/impls/mffd/mffdimpl.h> /*I "petscmat.h" I*/
5: PetscFList 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()
22: @*/
23: PetscErrorCode MatMFFDFinalizePackage(void)
24: {
26: MatMFFDPackageInitialized = PETSC_FALSE;
27: MatMFFDRegisterAllCalled = PETSC_FALSE;
28: MatMFFDList = PETSC_NULL;
29: return(0);
30: }
34: /*@C
35: MatMFFDInitializePackage - This function initializes everything in the MatMFFD package. It is called
36: from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to MatCreate_MFFD()
37: when using static libraries.
39: Input Parameter:
40: . path - The dynamic library path, or PETSC_NULL
42: Level: developer
44: .keywords: Vec, initialize, package
45: .seealso: PetscInitialize()
46: @*/
47: PetscErrorCode MatMFFDInitializePackage(const char path[])
48: {
49: char logList[256];
50: char *className;
51: PetscBool opt;
52: PetscErrorCode ierr;
55: if (MatMFFDPackageInitialized) return(0);
56: MatMFFDPackageInitialized = PETSC_TRUE;
57: /* Register Classes */
58: PetscClassIdRegister("MatMFFD",&MATMFFD_CLASSID);
59: /* Register Constructors */
60: MatMFFDRegisterAll(path);
61: /* Register Events */
62: PetscLogEventRegister("MatMult MF", MATMFFD_CLASSID,&MATMFFD_Mult);
64: /* Process info exclusions */
65: PetscOptionsGetString(PETSC_NULL, "-info_exclude", logList, 256, &opt);
66: if (opt) {
67: PetscStrstr(logList, "matmffd", &className);
68: if (className) {
69: PetscInfoDeactivateClass(MATMFFD_CLASSID);
70: }
71: }
72: /* Process summary exclusions */
73: PetscOptionsGetString(PETSC_NULL, "-log_summary_exclude", logList, 256, &opt);
74: if (opt) {
75: PetscStrstr(logList, "matmffd", &className);
76: if (className) {
77: PetscLogEventDeactivateClass(MATMFFD_CLASSID);
78: }
79: }
80: PetscRegisterFinalize(MatMFFDFinalizePackage);
81: return(0);
82: }
86: /*@C
87: MatMFFDSetType - Sets the method that is used to compute the
88: differencing parameter for finite differene matrix-free formulations.
90: Input Parameters:
91: + mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMFFD()
92: or MatSetType(mat,MATMFFD);
93: - ftype - the type requested, either MATMFFD_WP or MATMFFD_DS
95: Level: advanced
97: Notes:
98: For example, such routines can compute h for use in
99: Jacobian-vector products of the form
101: F(x+ha) - F(x)
102: F'(u)a ~= ----------------
103: h
105: .seealso: MatCreateSNESMF(), MatMFFDRegisterDynamic(), MatMFFDSetFunction()
106: @*/
107: PetscErrorCode MatMFFDSetType(Mat mat,const MatMFFDType ftype)
108: {
109: PetscErrorCode ierr,(*r)(MatMFFD);
110: MatMFFD ctx = (MatMFFD)mat->data;
111: PetscBool match;
112:
117: PetscTypeCompare((PetscObject)mat,MATMFFD,&match);
118: if (!match) return(0);
120: /* already set, so just return */
121: PetscTypeCompare((PetscObject)ctx,ftype,&match);
122: if (match) return(0);
124: /* destroy the old one if it exists */
125: if (ctx->ops->destroy) {
126: (*ctx->ops->destroy)(ctx);
127: }
129: PetscFListFind(MatMFFDList,((PetscObject)ctx)->comm,ftype,PETSC_TRUE,(void (**)(void)) &r);
130: if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown MatMFFD type %s given",ftype);
131: (*r)(ctx);
132: PetscObjectChangeTypeName((PetscObject)ctx,ftype);
133: return(0);
134: }
140: PetscErrorCode MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func)
141: {
142: MatMFFD ctx = (MatMFFD)mat->data;
145: ctx->funcisetbase = func;
146: return(0);
147: }
154: PetscErrorCode MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci)
155: {
156: MatMFFD ctx = (MatMFFD)mat->data;
159: ctx->funci = funci;
160: return(0);
161: }
167: PetscErrorCode MatMFFDResetHHistory_MFFD(Mat J)
168: {
169: MatMFFD ctx = (MatMFFD)J->data;
172: ctx->ncurrenth = 0;
173: return(0);
174: }
179: PetscErrorCode MatMFFDRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(MatMFFD))
180: {
182: char fullname[PETSC_MAX_PATH_LEN];
185: PetscFListConcat(path,name,fullname);
186: PetscFListAdd(&MatMFFDList,sname,fullname,(void (*)(void))function);
187: return(0);
188: }
193: /*@C
194: MatMFFDRegisterDestroy - Frees the list of MatMFFD methods that were
195: registered by MatMFFDRegisterDynamic).
197: Not Collective
199: Level: developer
201: .keywords: MatMFFD, register, destroy
203: .seealso: MatMFFDRegisterDynamic), MatMFFDRegisterAll()
204: @*/
205: PetscErrorCode MatMFFDRegisterDestroy(void)
206: {
210: PetscFListDestroy(&MatMFFDList);
211: MatMFFDRegisterAllCalled = PETSC_FALSE;
212: return(0);
213: }
218: PetscErrorCode MatMFFDAddNullSpace_MFFD(Mat J,MatNullSpace nullsp)
219: {
221: MatMFFD ctx = (MatMFFD)J->data;
224: PetscObjectReference((PetscObject)nullsp);
225: if (ctx->sp) { MatNullSpaceDestroy(&ctx->sp); }
226: ctx->sp = nullsp;
227: return(0);
228: }
231: /* ----------------------------------------------------------------------------------------*/
234: PetscErrorCode MatDestroy_MFFD(Mat mat)
235: {
237: MatMFFD ctx = (MatMFFD)mat->data;
240: VecDestroy(&ctx->w);
241: VecDestroy(&ctx->drscale);
242: VecDestroy(&ctx->dlscale);
243: VecDestroy(&ctx->dshift);
244: if (ctx->current_f_allocated) {
245: VecDestroy(&ctx->current_f);
246: }
247: if (ctx->ops->destroy) {(*ctx->ops->destroy)(ctx);}
248: MatNullSpaceDestroy(&ctx->sp);
249: PetscHeaderDestroy(&ctx);
250: mat->data = 0;
252: PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C","",PETSC_NULL);
253: PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C","",PETSC_NULL);
254: PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C","",PETSC_NULL);
255: PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunction_C","",PETSC_NULL);
256: PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctionError_C","",PETSC_NULL);
257: PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C","",PETSC_NULL);
258: PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetPeriod_C","",PETSC_NULL);
259: PetscObjectComposeFunction((PetscObject)mat,"MatMFFDResetHHistory_C","",PETSC_NULL);
260: PetscObjectComposeFunction((PetscObject)mat,"MatMFFDAddNullSpace_C","",PETSC_NULL);
262: return(0);
263: }
267: /*
268: MatMFFDView_MFFD - Views matrix-free parameters.
270: */
271: PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer)
272: {
274: MatMFFD ctx = (MatMFFD)J->data;
275: PetscBool iascii;
278: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
279: if (iascii) {
280: PetscViewerASCIIPrintf(viewer," matrix-free approximation:\n");
281: PetscViewerASCIIPrintf(viewer," err=%G (relative error in function evaluation)\n",ctx->error_rel);
282: if (!((PetscObject)ctx)->type_name) {
283: PetscViewerASCIIPrintf(viewer," The compute h routine has not yet been set\n");
284: } else {
285: PetscViewerASCIIPrintf(viewer," Using %s compute h routine\n",((PetscObject)ctx)->type_name);
286: }
287: if (ctx->ops->view) {
288: (*ctx->ops->view)(ctx,viewer);
289: }
290: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for matrix-free matrix",((PetscObject)viewer)->type_name);
291: return(0);
292: }
296: /*
297: MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This
298: allows the user to indicate the beginning of a new linear solve by calling
299: MatAssemblyXXX() on the matrix free matrix. This then allows the
300: MatCreateMFFD_WP() to properly compute ||U|| only the first time
301: in the linear solver rather than every time.
302: */
303: PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt)
304: {
306: MatMFFD j = (MatMFFD)J->data;
309: MatMFFDResetHHistory(J);
310: j->vshift = 0.0;
311: j->vscale = 1.0;
312: return(0);
313: }
317: /*
318: MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:
320: y ~= (F(u + ha) - F(u))/h,
321: where F = nonlinear function, as set by SNESSetFunction()
322: u = current iterate
323: h = difference interval
324: */
325: PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y)
326: {
327: MatMFFD ctx = (MatMFFD)mat->data;
328: PetscScalar h;
329: Vec w,U,F;
331: PetscBool zeroa;
334: if (!ctx->current_u) SETERRQ(((PetscObject)mat)->comm,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");
335: /* We log matrix-free matrix-vector products separately, so that we can
336: separate the performance monitoring from the cases that use conventional
337: storage. We may eventually modify event logging to associate events
338: with particular objects, hence alleviating the more general problem. */
339: PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);
341: w = ctx->w;
342: U = ctx->current_u;
343: F = ctx->current_f;
344: /*
345: Compute differencing parameter
346: */
347: if (!ctx->ops->compute) {
348: MatMFFDSetType(mat,MATMFFD_WP);
349: MatSetFromOptions(mat);
350: }
351: (*ctx->ops->compute)(ctx,U,a,&h,&zeroa);
352: if (zeroa) {
353: VecSet(y,0.0);
354: return(0);
355: }
357: if (PetscIsInfOrNanScalar(h)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed Nan differencing parameter h");
358: if (ctx->checkh) {
359: (*ctx->checkh)(ctx->checkhctx,U,a,&h);
360: }
362: /* keep a record of the current differencing parameter h */
363: ctx->currenth = h;
364: #if defined(PETSC_USE_COMPLEX)
365: PetscInfo2(mat,"Current differencing parameter: %G + %G i\n",PetscRealPart(h),PetscImaginaryPart(h));
366: #else
367: PetscInfo1(mat,"Current differencing parameter: %15.12e\n",h);
368: #endif
369: if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) {
370: ctx->historyh[ctx->ncurrenth] = h;
371: }
372: ctx->ncurrenth++;
374: /* w = u + ha */
375: if (ctx->drscale) {
376: VecPointwiseMult(ctx->drscale,a,U);
377: VecAYPX(U,h,w);
378: } else {
379: VecWAXPY(w,h,a,U);
380: }
382: /* compute func(U) as base for differencing; only needed first time in and not when provided by user */
383: if (ctx->ncurrenth == 1 && ctx->current_f_allocated) {
384: (*ctx->func)(ctx->funcctx,U,F);
385: }
386: (*ctx->func)(ctx->funcctx,w,y);
388: VecAXPY(y,-1.0,F);
389: VecScale(y,1.0/h);
391: if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
392: VecAXPBY(y,ctx->vshift,ctx->vscale,a);
393: }
394: if (ctx->dlscale) {
395: VecPointwiseMult(y,ctx->dlscale,y);
396: }
397: if (ctx->dshift) {
398: VecPointwiseMult(ctx->dshift,a,U);
399: VecAXPY(y,1.0,U);
400: }
402: if (ctx->sp) {MatNullSpaceRemove(ctx->sp,y,PETSC_NULL);}
404: PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);
405: return(0);
406: }
410: /*
411: MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix
413: y ~= (F(u + ha) - F(u))/h,
414: where F = nonlinear function, as set by SNESSetFunction()
415: u = current iterate
416: h = difference interval
417: */
418: PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a)
419: {
420: MatMFFD ctx = (MatMFFD)mat->data;
421: PetscScalar h,*aa,*ww,v;
422: PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
423: Vec w,U;
425: PetscInt i,rstart,rend;
428: if (!ctx->funci) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Requires calling MatMFFDSetFunctioni() first");
430: w = ctx->w;
431: U = ctx->current_u;
432: (*ctx->func)(ctx->funcctx,U,a);
433: (*ctx->funcisetbase)(ctx->funcctx,U);
434: VecCopy(U,w);
436: VecGetOwnershipRange(a,&rstart,&rend);
437: VecGetArray(a,&aa);
438: for (i=rstart; i<rend; i++) {
439: VecGetArray(w,&ww);
440: h = ww[i-rstart];
441: if (h == 0.0) h = 1.0;
442: #if !defined(PETSC_USE_COMPLEX)
443: if (h < umin && h >= 0.0) h = umin;
444: else if (h < 0.0 && h > -umin) h = -umin;
445: #else
446: if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin;
447: else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
448: #endif
449: h *= epsilon;
450:
451: ww[i-rstart] += h;
452: VecRestoreArray(w,&ww);
453: (*ctx->funci)(ctx->funcctx,i,w,&v);
454: aa[i-rstart] = (v - aa[i-rstart])/h;
456: /* possibly shift and scale result */
457: if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
458: aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart];
459: }
461: VecGetArray(w,&ww);
462: ww[i-rstart] -= h;
463: VecRestoreArray(w,&ww);
464: }
465: VecRestoreArray(a,&aa);
466: return(0);
467: }
471: PetscErrorCode MatDiagonalScale_MFFD(Mat mat,Vec ll,Vec rr)
472: {
473: MatMFFD aij = (MatMFFD)mat->data;
477: if (ll && !aij->dlscale) {
478: VecDuplicate(ll,&aij->dlscale);
479: }
480: if (rr && !aij->drscale) {
481: VecDuplicate(rr,&aij->drscale);
482: }
483: if (ll) {
484: VecCopy(ll,aij->dlscale);
485: }
486: if (rr) {
487: VecCopy(rr,aij->drscale);
488: }
489: return(0);
490: }
494: PetscErrorCode MatDiagonalSet_MFFD(Mat mat,Vec ll,InsertMode mode)
495: {
496: MatMFFD aij = (MatMFFD)mat->data;
500: if (mode == INSERT_VALUES) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"No diagonal set with INSERT_VALUES");
501: if (!aij->dshift) {
502: VecDuplicate(ll,&aij->dshift);
503: }
504: VecAXPY(aij->dshift,1.0,ll);
505: return(0);
506: }
510: PetscErrorCode MatShift_MFFD(Mat Y,PetscScalar a)
511: {
512: MatMFFD shell = (MatMFFD)Y->data;
514: shell->vshift += a;
515: return(0);
516: }
520: PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a)
521: {
522: MatMFFD shell = (MatMFFD)Y->data;
524: shell->vscale *= a;
525: return(0);
526: }
531: PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F)
532: {
534: MatMFFD ctx = (MatMFFD)J->data;
537: MatMFFDResetHHistory(J);
538: ctx->current_u = U;
539: if (F) {
540: if (ctx->current_f_allocated) {VecDestroy(&ctx->current_f);}
541: ctx->current_f = F;
542: ctx->current_f_allocated = PETSC_FALSE;
543: } else if (!ctx->current_f_allocated) {
544: VecDuplicate(ctx->current_u, &ctx->current_f);
545: ctx->current_f_allocated = PETSC_TRUE;
546: }
547: if (!ctx->w) {
548: VecDuplicate(ctx->current_u, &ctx->w);
549: }
550: J->assembled = PETSC_TRUE;
551: return(0);
552: }
558: PetscErrorCode MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void*ectx)
559: {
560: MatMFFD ctx = (MatMFFD)J->data;
563: ctx->checkh = fun;
564: ctx->checkhctx = ectx;
565: return(0);
566: }
571: /*@C
572: MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all
573: MatMFFD options in the database.
575: Collective on Mat
577: Input Parameter:
578: + A - the Mat context
579: - prefix - the prefix to prepend to all option names
581: Notes:
582: A hyphen (-) must NOT be given at the beginning of the prefix name.
583: The first character of all runtime options is AUTOMATICALLY the hyphen.
585: Level: advanced
587: .keywords: SNES, matrix-free, parameters
589: .seealso: MatSetFromOptions(), MatCreateSNESMF()
590: @*/
591: PetscErrorCode MatMFFDSetOptionsPrefix(Mat mat,const char prefix[])
593: {
594: MatMFFD mfctx = mat ? (MatMFFD)mat->data : PETSC_NULL;
599: PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix);
600: return(0);
601: }
605: PetscErrorCode MatSetFromOptions_MFFD(Mat mat)
606: {
607: MatMFFD mfctx = (MatMFFD)mat->data;
609: PetscBool flg;
610: char ftype[256];
615: PetscObjectOptionsBegin((PetscObject)mfctx);
616: PetscOptionsList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg);
617: if (flg) {
618: MatMFFDSetType(mat,ftype);
619: }
621: PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);
622: PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);
624: flg = PETSC_FALSE;
625: PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,PETSC_NULL);
626: if (flg) {
627: MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);
628: }
629: if (mfctx->ops->setfromoptions) {
630: (*mfctx->ops->setfromoptions)(mfctx);
631: }
632: PetscOptionsEnd();
633: return(0);
634: }
639: PetscErrorCode MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period)
640: {
641: MatMFFD ctx = (MatMFFD)mat->data;
645: ctx->recomputeperiod = period;
646: return(0);
647: }
653: PetscErrorCode MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
654: {
655: MatMFFD ctx = (MatMFFD)mat->data;
658: ctx->func = func;
659: ctx->funcctx = funcctx;
660: return(0);
661: }
667: PetscErrorCode MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error)
668: {
669: MatMFFD ctx = (MatMFFD)mat->data;
673: if (error != PETSC_DEFAULT) ctx->error_rel = error;
674: return(0);
675: }
678: /*MC
679: MATMFFD - MATMFFD = "mffd" - A matrix free matrix type.
681: Level: advanced
683: .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction()
684: M*/
688: PetscErrorCode MatCreate_MFFD(Mat A)
689: {
690: MatMFFD mfctx;
691: PetscErrorCode ierr;
694: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
695: MatMFFDInitializePackage(PETSC_NULL);
696: #endif
698: PetscHeaderCreate(mfctx,_p_MatMFFD,struct _MFOps,MATMFFD_CLASSID,0,"MatMFFD","Matrix-free Finite Differencing","Mat",((PetscObject)A)->comm,MatDestroy_MFFD,MatView_MFFD);
699: mfctx->sp = 0;
700: mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON;
701: mfctx->recomputeperiod = 1;
702: mfctx->count = 0;
703: mfctx->currenth = 0.0;
704: mfctx->historyh = PETSC_NULL;
705: mfctx->ncurrenth = 0;
706: mfctx->maxcurrenth = 0;
707: ((PetscObject)mfctx)->type_name = 0;
709: mfctx->vshift = 0.0;
710: mfctx->vscale = 1.0;
712: /*
713: Create the empty data structure to contain compute-h routines.
714: These will be filled in below from the command line options or
715: a later call with MatMFFDSetType() or if that is not called
716: then it will default in the first use of MatMult_MFFD()
717: */
718: mfctx->ops->compute = 0;
719: mfctx->ops->destroy = 0;
720: mfctx->ops->view = 0;
721: mfctx->ops->setfromoptions = 0;
722: mfctx->hctx = 0;
724: mfctx->func = 0;
725: mfctx->funcctx = 0;
726: mfctx->w = PETSC_NULL;
728: A->data = mfctx;
730: A->ops->mult = MatMult_MFFD;
731: A->ops->destroy = MatDestroy_MFFD;
732: A->ops->view = MatView_MFFD;
733: A->ops->assemblyend = MatAssemblyEnd_MFFD;
734: A->ops->getdiagonal = MatGetDiagonal_MFFD;
735: A->ops->scale = MatScale_MFFD;
736: A->ops->shift = MatShift_MFFD;
737: A->ops->diagonalscale = MatDiagonalScale_MFFD;
738: A->ops->diagonalset = MatDiagonalSet_MFFD;
739: A->ops->setfromoptions = MatSetFromOptions_MFFD;
740: A->assembled = PETSC_TRUE;
742: PetscLayoutSetBlockSize(A->rmap,1);
743: PetscLayoutSetBlockSize(A->cmap,1);
744: PetscLayoutSetUp(A->rmap);
745: PetscLayoutSetUp(A->cmap);
747: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetBase_C","MatMFFDSetBase_MFFD",MatMFFDSetBase_MFFD);
748: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunctioniBase_C","MatMFFDSetFunctioniBase_MFFD",MatMFFDSetFunctioniBase_MFFD);
749: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunctioni_C","MatMFFDSetFunctioni_MFFD",MatMFFDSetFunctioni_MFFD);
750: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunction_C","MatMFFDSetFunction_MFFD",MatMFFDSetFunction_MFFD);
751: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetCheckh_C","MatMFFDSetCheckh_MFFD",MatMFFDSetCheckh_MFFD);
752: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetPeriod_C","MatMFFDSetPeriod_MFFD",MatMFFDSetPeriod_MFFD);
753: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunctionError_C","MatMFFDSetFunctionError_MFFD",MatMFFDSetFunctionError_MFFD);
754: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDResetHHistory_C","MatMFFDResetHHistory_MFFD",MatMFFDResetHHistory_MFFD);
755: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDAddNullSpace_C","MatMFFDAddNullSpace_MFFD",MatMFFDAddNullSpace_MFFD);
756: mfctx->mat = A;
757: PetscObjectChangeTypeName((PetscObject)A,MATMFFD);
758: return(0);
759: }
764: /*@
765: MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF()
767: Collective on Vec
769: Input Parameters:
770: + comm - MPI communicator
771: . m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
772: This value should be the same as the local size used in creating the
773: y vector for the matrix-vector product y = Ax.
774: . n - This value should be the same as the local size used in creating the
775: x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
776: calculated if N is given) For square matrices n is almost always m.
777: . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
778: - N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
781: Output Parameter:
782: . J - the matrix-free matrix
784: Options Database Keys: call MatSetFromOptions() to trigger these
785: + -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS)
786: - -mat_mffd_err - square root of estimated relative error in function evaluation
787: - -mat_mffd_period - how often h is recomputed, defaults to 1, everytime
790: Level: advanced
792: Notes:
793: The matrix-free matrix context merely contains the function pointers
794: and work space for performing finite difference approximations of
795: Jacobian-vector products, F'(u)*a,
797: The default code uses the following approach to compute h
799: .vb
800: F'(u)*a = [F(u+h*a) - F(u)]/h where
801: h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1}
802: = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise
803: where
804: error_rel = square root of relative error in function evaluation
805: umin = minimum iterate parameter
806: .ve
808: The user can set the error_rel via MatMFFDSetFunctionError() and
809: umin via MatMFFDDSSetUmin(); see the <A href="../../docs/manual.pdf#nameddest=ch_snes">SNES chapter of the users manual</A> for details.
811: The user should call MatDestroy() when finished with the matrix-free
812: matrix context.
814: Options Database Keys:
815: + -mat_mffd_err <error_rel> - Sets error_rel
816: . -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only)
817: - -mat_mffd_check_positivity
819: .keywords: default, matrix-free, create, matrix
821: .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction()
822: MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(),
823: MatMFFDGetH(), MatMFFDRegisterDynamic), MatMFFDComputeJacobian()
824:
825: @*/
826: PetscErrorCode MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J)
827: {
831: MatCreate(comm,J);
832: MatSetSizes(*J,m,n,M,N);
833: MatSetType(*J,MATMFFD);
834: return(0);
835: }
840: /*@
841: MatMFFDGetH - Gets the last value that was used as the differencing
842: parameter.
844: Not Collective
846: Input Parameters:
847: . mat - the matrix obtained with MatCreateSNESMF()
849: Output Paramter:
850: . h - the differencing step size
852: Level: advanced
854: .keywords: SNES, matrix-free, parameters
856: .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory()
857: @*/
858: PetscErrorCode MatMFFDGetH(Mat mat,PetscScalar *h)
859: {
860: MatMFFD ctx = (MatMFFD)mat->data;
862: PetscBool match;
865: PetscTypeCompare((PetscObject)mat,MATMFFD,&match);
866: if (!match) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
868: *h = ctx->currenth;
869: return(0);
870: }
874: /*@C
875: MatMFFDSetFunction - Sets the function used in applying the matrix free.
877: Logically Collective on Mat
879: Input Parameters:
880: + mat - the matrix free matrix created via MatCreateSNESMF()
881: . func - the function to use
882: - funcctx - optional function context passed to function
884: Level: advanced
886: Notes:
887: If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
888: matrix inside your compute Jacobian routine
890: If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used.
892: .keywords: SNES, matrix-free, function
894: .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD,
895: MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
896: @*/
897: PetscErrorCode MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
898: {
902: PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx));
903: return(0);
904: }
908: /*@C
909: MatMFFDSetFunctioni - Sets the function for a single component
911: Logically Collective on Mat
913: Input Parameters:
914: + mat - the matrix free matrix created via MatCreateSNESMF()
915: - funci - the function to use
917: Level: advanced
919: Notes:
920: If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
921: matrix inside your compute Jacobian routine
924: .keywords: SNES, matrix-free, function
926: .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
928: @*/
929: PetscErrorCode MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*))
930: {
935: PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));
936: return(0);
937: }
942: /*@C
943: MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation
945: Logically Collective on Mat
947: Input Parameters:
948: + mat - the matrix free matrix created via MatCreateSNESMF()
949: - func - the function to use
951: Level: advanced
953: Notes:
954: If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
955: matrix inside your compute Jacobian routine
958: .keywords: SNES, matrix-free, function
960: .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
961: MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
962: @*/
963: PetscErrorCode MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec))
964: {
969: PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));
970: return(0);
971: }
975: /*@
976: MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime
978: Logically Collective on Mat
980: Input Parameters:
981: + mat - the matrix free matrix created via MatCreateSNESMF()
982: - period - 1 for everytime, 2 for every second etc
984: Options Database Keys:
985: + -mat_mffd_period <period>
987: Level: advanced
990: .keywords: SNES, matrix-free, parameters
992: .seealso: MatCreateSNESMF(),MatMFFDGetH(),
993: MatMFFDSetHHistory(), MatMFFDResetHHistory()
994: @*/
995: PetscErrorCode MatMFFDSetPeriod(Mat mat,PetscInt period)
996: {
1000: PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period));
1001: return(0);
1002: }
1006: /*@
1007: MatMFFDSetFunctionError - Sets the error_rel for the approximation of
1008: matrix-vector products using finite differences.
1010: Logically Collective on Mat
1012: Input Parameters:
1013: + mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF()
1014: - error_rel - relative error (should be set to the square root of
1015: the relative error in the function evaluations)
1017: Options Database Keys:
1018: + -mat_mffd_err <error_rel> - Sets error_rel
1020: Level: advanced
1022: Notes:
1023: The default matrix-free matrix-vector product routine computes
1024: .vb
1025: F'(u)*a = [F(u+h*a) - F(u)]/h where
1026: h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1}
1027: = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else
1028: .ve
1030: .keywords: SNES, matrix-free, parameters
1032: .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
1033: MatMFFDSetHHistory(), MatMFFDResetHHistory()
1034: @*/
1035: PetscErrorCode MatMFFDSetFunctionError(Mat mat,PetscReal error)
1036: {
1040: PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error));
1041: return(0);
1042: }
1046: /*@
1047: MatMFFDAddNullSpace - Provides a null space that an operator is
1048: supposed to have. Since roundoff will create a small component in
1049: the null space, if you know the null space you may have it
1050: automatically removed.
1052: Logically Collective on Mat
1054: Input Parameters:
1055: + J - the matrix-free matrix context
1056: - nullsp - object created with MatNullSpaceCreate()
1058: Level: advanced
1060: .keywords: SNES, matrix-free, null space
1062: .seealso: MatNullSpaceCreate(), MatMFFDGetH(), MatCreateSNESMF(), MatCreateMFFD(), MATMFFD
1063: MatMFFDSetHHistory(), MatMFFDResetHHistory()
1064: @*/
1065: PetscErrorCode MatMFFDAddNullSpace(Mat J,MatNullSpace nullsp)
1066: {
1070: PetscTryMethod(J,"MatMFFDAddNullSpace_C",(Mat,MatNullSpace),(J,nullsp));
1071: return(0);
1072: }
1076: /*@
1077: MatMFFDSetHHistory - Sets an array to collect a history of the
1078: differencing values (h) computed for the matrix-free product.
1080: Logically Collective on Mat
1082: Input Parameters:
1083: + J - the matrix-free matrix context
1084: . histroy - space to hold the history
1085: - nhistory - number of entries in history, if more entries are generated than
1086: nhistory, then the later ones are discarded
1088: Level: advanced
1090: Notes:
1091: Use MatMFFDResetHHistory() to reset the history counter and collect
1092: a new batch of differencing parameters, h.
1094: .keywords: SNES, matrix-free, h history, differencing history
1096: .seealso: MatMFFDGetH(), MatCreateSNESMF(),
1097: MatMFFDResetHHistory(), MatMFFDSetFunctionError()
1099: @*/
1100: PetscErrorCode MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory)
1101: {
1102: MatMFFD ctx = (MatMFFD)J->data;
1104: PetscBool match;
1107: PetscTypeCompare((PetscObject)J,MATMFFD,&match);
1108: if (!match) SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
1109: ctx->historyh = history;
1110: ctx->maxcurrenth = nhistory;
1111: ctx->currenth = 0.;
1112: return(0);
1113: }
1118: /*@
1119: MatMFFDResetHHistory - Resets the counter to zero to begin
1120: collecting a new set of differencing histories.
1122: Logically Collective on Mat
1124: Input Parameters:
1125: . J - the matrix-free matrix context
1127: Level: advanced
1129: Notes:
1130: Use MatMFFDSetHHistory() to create the original history counter.
1132: .keywords: SNES, matrix-free, h history, differencing history
1134: .seealso: MatMFFDGetH(), MatCreateSNESMF(),
1135: MatMFFDSetHHistory(), MatMFFDSetFunctionError()
1137: @*/
1138: PetscErrorCode MatMFFDResetHHistory(Mat J)
1139: {
1143: PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J));
1144: return(0);
1145: }
1150: /*@
1151: MatMFFDSetBase - Sets the vector U at which matrix vector products of the
1152: Jacobian are computed
1154: Logically Collective on Mat
1156: Input Parameters:
1157: + J - the MatMFFD matrix
1158: . U - the vector
1159: - F - (optional) vector that contains F(u) if it has been already computed
1161: Notes: This is rarely used directly
1163: If F is provided then it is not recomputed. Otherwise the function is evaluated at the base
1164: point during the first MatMult() after each call to MatMFFDSetBase().
1166: Level: advanced
1168: @*/
1169: PetscErrorCode MatMFFDSetBase(Mat J,Vec U,Vec F)
1170: {
1177: PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));
1178: return(0);
1179: }
1183: /*@C
1184: MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts
1185: it to satisfy some criteria
1187: Logically Collective on Mat
1189: Input Parameters:
1190: + J - the MatMFFD matrix
1191: . fun - the function that checks h
1192: - ctx - any context needed by the function
1194: Options Database Keys:
1195: . -mat_mffd_check_positivity
1197: Level: advanced
1199: Notes: For example, MatMFFDSetCheckPositivity() insures that all entries
1200: of U + h*a are non-negative
1202: .seealso: MatMFFDSetCheckPositivity()
1203: @*/
1204: PetscErrorCode MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void* ctx)
1205: {
1210: PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));
1211: return(0);
1212: }
1216: /*@
1217: MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or
1218: zero, decreases h until this is satisfied.
1220: Logically Collective on Vec
1222: Input Parameters:
1223: + U - base vector that is added to
1224: . a - vector that is added
1225: . h - scaling factor on a
1226: - dummy - context variable (unused)
1228: Options Database Keys:
1229: . -mat_mffd_check_positivity
1231: Level: advanced
1233: Notes: This is rarely used directly, rather it is passed as an argument to
1234: MatMFFDSetCheckh()
1236: .seealso: MatMFFDSetCheckh()
1237: @*/
1238: PetscErrorCode MatMFFDCheckPositivity(void* dummy,Vec U,Vec a,PetscScalar *h)
1239: {
1240: PetscReal val, minval;
1241: PetscScalar *u_vec, *a_vec;
1243: PetscInt i,n;
1244: MPI_Comm comm;
1247: PetscObjectGetComm((PetscObject)U,&comm);
1248: VecGetArray(U,&u_vec);
1249: VecGetArray(a,&a_vec);
1250: VecGetLocalSize(U,&n);
1251: minval = PetscAbsScalar(*h*1.01);
1252: for(i=0;i<n;i++) {
1253: if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) {
1254: val = PetscAbsScalar(u_vec[i]/a_vec[i]);
1255: if (val < minval) minval = val;
1256: }
1257: }
1258: VecRestoreArray(U,&u_vec);
1259: VecRestoreArray(a,&a_vec);
1260: MPI_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);
1261: if (val <= PetscAbsScalar(*h)) {
1262: PetscInfo2(U,"Scaling back h from %G to %G\n",PetscRealPart(*h),.99*val);
1263: if (PetscRealPart(*h) > 0.0) *h = 0.99*val;
1264: else *h = -0.99*val;
1265: }
1266: return(0);
1267: }