1 #include <../src/mat/impls/shell/shell.h> /*I "petscmat.h" I*/
2 #include <../src/mat/impls/mffd/mffdimpl.h>
3
4 PetscFunctionList MatMFFDList = NULL;
5 PetscBool MatMFFDRegisterAllCalled = PETSC_FALSE;
6
7 PetscClassId MATMFFD_CLASSID;
8 PetscLogEvent MATMFFD_Mult;
9
10 static PetscBool MatMFFDPackageInitialized = PETSC_FALSE;
11
12 /*@C
13 MatMFFDFinalizePackage - This function destroys everything in the MATMFFD` package. It is
14 called from `PetscFinalize()`.
15
16 Level: developer
17
18 .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `PetscFinalize()`, `MatCreateMFFD()`, `MatCreateSNESMF()`
19 @*/
MatMFFDFinalizePackage(void)20 PetscErrorCode MatMFFDFinalizePackage(void)
21 {
22 PetscFunctionBegin;
23 PetscCall(PetscFunctionListDestroy(&MatMFFDList));
24 MatMFFDPackageInitialized = PETSC_FALSE;
25 MatMFFDRegisterAllCalled = PETSC_FALSE;
26 PetscFunctionReturn(PETSC_SUCCESS);
27 }
28
29 /*@C
30 MatMFFDInitializePackage - This function initializes everything in the MATMFFD` package. It is called
31 from `MatInitializePackage()`.
32
33 Level: developer
34
35 .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `PetscInitialize()`
36 @*/
MatMFFDInitializePackage(void)37 PetscErrorCode MatMFFDInitializePackage(void)
38 {
39 char logList[256];
40 PetscBool opt, pkg;
41
42 PetscFunctionBegin;
43 if (MatMFFDPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
44 MatMFFDPackageInitialized = PETSC_TRUE;
45 /* Register Classes */
46 PetscCall(PetscClassIdRegister("MatMFFD", &MATMFFD_CLASSID));
47 /* Register Constructors */
48 PetscCall(MatMFFDRegisterAll());
49 /* Register Events */
50 PetscCall(PetscLogEventRegister("MatMult MF", MATMFFD_CLASSID, &MATMFFD_Mult));
51 /* Process Info */
52 {
53 PetscClassId classids[1];
54
55 classids[0] = MATMFFD_CLASSID;
56 PetscCall(PetscInfoProcessClass("matmffd", 1, classids));
57 }
58 /* Process summary exclusions */
59 PetscCall(PetscOptionsGetString(NULL, NULL, "-log_exclude", logList, sizeof(logList), &opt));
60 if (opt) {
61 PetscCall(PetscStrInList("matmffd", logList, ',', &pkg));
62 if (pkg) PetscCall(PetscLogEventExcludeClass(MATMFFD_CLASSID));
63 }
64 /* Register package finalizer */
65 PetscCall(PetscRegisterFinalize(MatMFFDFinalizePackage));
66 PetscFunctionReturn(PETSC_SUCCESS);
67 }
68
MatMFFDSetType_MFFD(Mat mat,MatMFFDType ftype)69 static PetscErrorCode MatMFFDSetType_MFFD(Mat mat, MatMFFDType ftype)
70 {
71 MatMFFD ctx;
72 PetscBool match;
73 PetscErrorCode (*r)(MatMFFD);
74
75 PetscFunctionBegin;
76 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
77 PetscAssertPointer(ftype, 2);
78 PetscCall(MatShellGetContext(mat, &ctx));
79
80 /* already set, so just return */
81 PetscCall(PetscObjectTypeCompare((PetscObject)ctx, ftype, &match));
82 if (match) PetscFunctionReturn(PETSC_SUCCESS);
83
84 /* destroy the old one if it exists */
85 PetscTryTypeMethod(ctx, destroy);
86
87 PetscCall(PetscFunctionListFind(MatMFFDList, ftype, &r));
88 PetscCheck(r, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown MatMFFD type %s given", ftype);
89 PetscCall((*r)(ctx));
90 PetscCall(PetscObjectChangeTypeName((PetscObject)ctx, ftype));
91 PetscFunctionReturn(PETSC_SUCCESS);
92 }
93
94 /*@
95 MatMFFDSetType - Sets the method that is used to compute the
96 differencing parameter for finite difference matrix-free formulations.
97
98 Input Parameters:
99 + mat - the "matrix-free" matrix created via `MatCreateSNESMF()`, or `MatCreateMFFD()`
100 or `MatSetType`(mat,`MATMFFD`);
101 - ftype - the type requested, either `MATMFFD_WP` or `MATMFFD_DS`
102
103 Level: advanced
104
105 Note:
106 For example, such routines can compute `h` for use in
107 Jacobian-vector products of the form
108 .vb
109
110 F(x+ha) - F(x)
111 F'(u)a ~= ----------------
112 h
113 .ve
114
115 .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MATMFFD_WP`, `MATMFFD_DS`, `MatCreateSNESMF()`, `MatMFFDRegister()`, `MatMFFDSetFunction()`, `MatCreateMFFD()`
116 @*/
MatMFFDSetType(Mat mat,MatMFFDType ftype)117 PetscErrorCode MatMFFDSetType(Mat mat, MatMFFDType ftype)
118 {
119 PetscFunctionBegin;
120 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
121 PetscAssertPointer(ftype, 2);
122 PetscTryMethod(mat, "MatMFFDSetType_C", (Mat, MatMFFDType), (mat, ftype));
123 PetscFunctionReturn(PETSC_SUCCESS);
124 }
125
126 static PetscErrorCode MatGetDiagonal_MFFD(Mat, Vec);
127
128 typedef PetscErrorCode (*FCN1)(void *, Vec); /* force argument to next function to not be extern C*/
MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func)129 static PetscErrorCode MatMFFDSetFunctioniBase_MFFD(Mat mat, FCN1 func)
130 {
131 MatMFFD ctx;
132
133 PetscFunctionBegin;
134 PetscCall(MatShellGetContext(mat, &ctx));
135 ctx->funcisetbase = func;
136 PetscFunctionReturn(PETSC_SUCCESS);
137 }
138
139 typedef PetscErrorCode (*FCN2)(void *, PetscInt, Vec, PetscScalar *); /* force argument to next function to not be extern C*/
MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci)140 static PetscErrorCode MatMFFDSetFunctioni_MFFD(Mat mat, FCN2 funci)
141 {
142 MatMFFD ctx;
143
144 PetscFunctionBegin;
145 PetscCall(MatShellGetContext(mat, &ctx));
146 ctx->funci = funci;
147 PetscCall(MatShellSetOperation(mat, MATOP_GET_DIAGONAL, (PetscErrorCodeFn *)MatGetDiagonal_MFFD));
148 PetscFunctionReturn(PETSC_SUCCESS);
149 }
150
MatMFFDGetH_MFFD(Mat mat,PetscScalar * h)151 static PetscErrorCode MatMFFDGetH_MFFD(Mat mat, PetscScalar *h)
152 {
153 MatMFFD ctx;
154
155 PetscFunctionBegin;
156 PetscCall(MatShellGetContext(mat, &ctx));
157 *h = ctx->currenth;
158 PetscFunctionReturn(PETSC_SUCCESS);
159 }
160
MatMFFDResetHHistory_MFFD(Mat J)161 static PetscErrorCode MatMFFDResetHHistory_MFFD(Mat J)
162 {
163 MatMFFD ctx;
164
165 PetscFunctionBegin;
166 PetscCall(MatShellGetContext(J, &ctx));
167 ctx->ncurrenth = 0;
168 PetscFunctionReturn(PETSC_SUCCESS);
169 }
170
171 /*@C
172 MatMFFDRegister - Adds a method to the `MATMFFD` registry.
173
174 Not Collective, No Fortran Support
175
176 Input Parameters:
177 + sname - name of a new user-defined compute-h module
178 - function - routine to create method context
179
180 Level: developer
181
182 Note:
183 `MatMFFDRegister()` may be called multiple times to add several user-defined solvers.
184
185 Example Usage:
186 .vb
187 MatMFFDRegister("my_h", MyHCreate);
188 .ve
189
190 Then, your solver can be chosen with the procedural interface via `MatMFFDSetType`(mfctx, "my_h")` or at runtime via the option
191 `-mat_mffd_type my_h`
192
193 .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDRegisterAll()`, `MatMFFDRegisterDestroy()`
194 @*/
MatMFFDRegister(const char sname[],PetscErrorCode (* function)(MatMFFD))195 PetscErrorCode MatMFFDRegister(const char sname[], PetscErrorCode (*function)(MatMFFD))
196 {
197 PetscFunctionBegin;
198 PetscCall(MatInitializePackage());
199 PetscCall(PetscFunctionListAdd(&MatMFFDList, sname, function));
200 PetscFunctionReturn(PETSC_SUCCESS);
201 }
202
MatDestroy_MFFD(Mat mat)203 static PetscErrorCode MatDestroy_MFFD(Mat mat)
204 {
205 MatMFFD ctx;
206
207 PetscFunctionBegin;
208 PetscCall(MatShellGetContext(mat, &ctx));
209 PetscCall(VecDestroy(&ctx->w));
210 PetscCall(VecDestroy(&ctx->current_u));
211 if (ctx->current_f_allocated) PetscCall(VecDestroy(&ctx->current_f));
212 PetscTryTypeMethod(ctx, destroy);
213 PetscCall(PetscHeaderDestroy(&ctx));
214
215 PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetBase_C", NULL));
216 PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctioniBase_C", NULL));
217 PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctioni_C", NULL));
218 PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunction_C", NULL));
219 PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetFunctionError_C", NULL));
220 PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetCheckh_C", NULL));
221 PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetPeriod_C", NULL));
222 PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDResetHHistory_C", NULL));
223 PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetHHistory_C", NULL));
224 PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDSetType_C", NULL));
225 PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMFFDGetH_C", NULL));
226 PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatSNESMFSetReuseBase_C", NULL));
227 PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatSNESMFGetReuseBase_C", NULL));
228 PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatShellSetContext_C", NULL));
229 PetscFunctionReturn(PETSC_SUCCESS);
230 }
231
232 /*
233 MatMFFDView_MFFD - Views matrix-free parameters.
234
235 */
MatView_MFFD(Mat J,PetscViewer viewer)236 static PetscErrorCode MatView_MFFD(Mat J, PetscViewer viewer)
237 {
238 MatMFFD ctx;
239 PetscBool isascii, viewbase, viewfunction;
240 const char *prefix;
241
242 PetscFunctionBegin;
243 PetscCall(MatShellGetContext(J, &ctx));
244 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
245 if (isascii) {
246 PetscCall(PetscViewerASCIIPrintf(viewer, "Matrix-free approximation:\n"));
247 PetscCall(PetscViewerASCIIPushTab(viewer));
248 PetscCall(PetscViewerASCIIPrintf(viewer, "err=%g (relative error in function evaluation)\n", (double)ctx->error_rel));
249 if (!((PetscObject)ctx)->type_name) {
250 PetscCall(PetscViewerASCIIPrintf(viewer, "The compute h routine has not yet been set\n"));
251 } else {
252 PetscCall(PetscViewerASCIIPrintf(viewer, "Using %s compute h routine\n", ((PetscObject)ctx)->type_name));
253 }
254 #if defined(PETSC_USE_COMPLEX)
255 if (ctx->usecomplex) PetscCall(PetscViewerASCIIPrintf(viewer, "Using Lyness complex number trick to compute the matrix-vector product\n"));
256 #endif
257 PetscTryTypeMethod(ctx, view, viewer);
258 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)J, &prefix));
259
260 PetscCall(PetscOptionsHasName(((PetscObject)J)->options, prefix, "-mat_mffd_view_base", &viewbase));
261 if (viewbase) {
262 PetscCall(PetscViewerASCIIPrintf(viewer, "Base:\n"));
263 PetscCall(VecView(ctx->current_u, viewer));
264 }
265 PetscCall(PetscOptionsHasName(((PetscObject)J)->options, prefix, "-mat_mffd_view_function", &viewfunction));
266 if (viewfunction) {
267 PetscCall(PetscViewerASCIIPrintf(viewer, "Function:\n"));
268 PetscCall(VecView(ctx->current_f, viewer));
269 }
270 PetscCall(PetscViewerASCIIPopTab(viewer));
271 }
272 PetscFunctionReturn(PETSC_SUCCESS);
273 }
274
275 /*
276 MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This
277 allows the user to indicate the beginning of a new linear solve by calling
278 MatAssemblyXXX() on the matrix-free matrix. This then allows the
279 MatCreateMFFD_WP() to properly compute ||U|| only the first time
280 in the linear solver rather than every time.
281
282 This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library hence
283 it must be labeled as PETSC_EXTERN
284 */
MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt)285 PETSC_SINGLE_LIBRARY_VISIBILITY_INTERNAL PetscErrorCode MatAssemblyEnd_MFFD(Mat J, MatAssemblyType mt)
286 {
287 MatMFFD j;
288
289 PetscFunctionBegin;
290 PetscCall(MatShellGetContext(J, &j));
291 PetscCall(MatMFFDResetHHistory(J));
292 PetscFunctionReturn(PETSC_SUCCESS);
293 }
294
295 /*
296 MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:
297
298 y ~= (F(u + ha) - F(u))/h,
299 where F = nonlinear function, as set by SNESSetFunction()
300 u = current iterate
301 h = difference interval
302 */
MatMult_MFFD(Mat mat,Vec a,Vec y)303 static PetscErrorCode MatMult_MFFD(Mat mat, Vec a, Vec y)
304 {
305 MatMFFD ctx;
306 PetscScalar h;
307 Vec w, U, F;
308 PetscBool zeroa;
309
310 PetscFunctionBegin;
311 PetscCall(MatShellGetContext(mat, &ctx));
312 PetscCheck(ctx->current_u, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "MatMFFDSetBase() has not been called, this is often caused by forgetting to call MatAssemblyBegin/End on the first Mat in the SNES compute function");
313 /* We log matrix-free matrix-vector products separately, so that we can
314 separate the performance monitoring from the cases that use conventional
315 storage. We may eventually modify event logging to associate events
316 with particular objects, hence alleviating the more general problem. */
317 PetscCall(PetscLogEventBegin(MATMFFD_Mult, a, y, 0, 0));
318
319 w = ctx->w;
320 U = ctx->current_u;
321 F = ctx->current_f;
322 /*
323 Compute differencing parameter
324 */
325 if (!((PetscObject)ctx)->type_name) {
326 PetscCall(MatMFFDSetType(mat, MATMFFD_WP));
327 PetscCall(MatSetFromOptions(mat));
328 }
329 PetscUseTypeMethod(ctx, compute, U, a, &h, &zeroa);
330 if (zeroa) {
331 PetscCall(VecSet(y, 0.0));
332 PetscCall(PetscLogEventEnd(MATMFFD_Mult, a, y, 0, 0));
333 PetscFunctionReturn(PETSC_SUCCESS);
334 }
335
336 PetscCheck(!mat->erroriffailure || !PetscIsInfOrNanScalar(h), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Computed NaN differencing parameter h");
337 if (ctx->checkh) PetscCall((*ctx->checkh)(ctx->checkhctx, U, a, &h));
338
339 /* keep a record of the current differencing parameter h */
340 ctx->currenth = h;
341 #if defined(PETSC_USE_COMPLEX)
342 PetscCall(PetscInfo(mat, "Current differencing parameter: %g + %g i\n", (double)PetscRealPart(h), (double)PetscImaginaryPart(h)));
343 #else
344 PetscCall(PetscInfo(mat, "Current differencing parameter: %15.12e\n", (double)PetscRealPart(h)));
345 #endif
346 if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) ctx->historyh[ctx->ncurrenth] = h;
347 ctx->ncurrenth++;
348
349 #if defined(PETSC_USE_COMPLEX)
350 if (ctx->usecomplex) h = PETSC_i * h;
351 #endif
352
353 /* w = u + ha */
354 PetscCall(VecWAXPY(w, h, a, U));
355
356 /* compute func(U) as base for differencing; only needed first time in and not when provided by user */
357 if (ctx->ncurrenth == 1 && ctx->current_f_allocated) PetscCall((*ctx->func)(ctx->funcctx, U, F));
358 PetscCall((*ctx->func)(ctx->funcctx, w, y));
359
360 #if defined(PETSC_USE_COMPLEX)
361 if (ctx->usecomplex) {
362 PetscCall(VecImaginaryPart(y));
363 h = PetscImaginaryPart(h);
364 } else {
365 PetscCall(VecAXPY(y, -1.0, F));
366 }
367 #else
368 PetscCall(VecAXPY(y, -1.0, F));
369 #endif
370 PetscCall(VecScale(y, 1.0 / h));
371 if (mat->nullsp) PetscCall(MatNullSpaceRemove(mat->nullsp, y));
372
373 PetscCall(PetscLogEventEnd(MATMFFD_Mult, a, y, 0, 0));
374 PetscFunctionReturn(PETSC_SUCCESS);
375 }
376
377 /*
378 MatGetDiagonal_MFFD - Gets the diagonal for a matrix-free matrix
379
380 y ~= (F(u + ha) - F(u))/h,
381 where F = nonlinear function, as set by SNESSetFunction()
382 u = current iterate
383 h = difference interval
384 */
MatGetDiagonal_MFFD(Mat mat,Vec a)385 static PetscErrorCode MatGetDiagonal_MFFD(Mat mat, Vec a)
386 {
387 MatMFFD ctx;
388 PetscScalar h, *aa, *ww, v;
389 PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON, umin = 100.0 * PETSC_SQRT_MACHINE_EPSILON;
390 Vec w, U;
391 PetscInt i, rstart, rend;
392
393 PetscFunctionBegin;
394 PetscCall(MatShellGetContext(mat, &ctx));
395 PetscCheck(ctx->func, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Requires calling MatMFFDSetFunction() first");
396 PetscCheck(ctx->funci, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Requires calling MatMFFDSetFunctioni() first");
397 w = ctx->w;
398 U = ctx->current_u;
399 PetscCall((*ctx->func)(ctx->funcctx, U, a));
400 if (ctx->funcisetbase) PetscCall((*ctx->funcisetbase)(ctx->funcctx, U));
401 PetscCall(VecCopy(U, w));
402
403 PetscCall(VecGetOwnershipRange(a, &rstart, &rend));
404 PetscCall(VecGetArray(a, &aa));
405 for (i = rstart; i < rend; i++) {
406 PetscCall(VecGetArray(w, &ww));
407 h = ww[i - rstart];
408 if (h == 0.0) h = 1.0;
409 if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin;
410 else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
411 h *= epsilon;
412
413 ww[i - rstart] += h;
414 PetscCall(VecRestoreArray(w, &ww));
415 PetscCall((*ctx->funci)(ctx->funcctx, i, w, &v));
416 aa[i - rstart] = (v - aa[i - rstart]) / h;
417
418 PetscCall(VecGetArray(w, &ww));
419 ww[i - rstart] -= h;
420 PetscCall(VecRestoreArray(w, &ww));
421 }
422 PetscCall(VecRestoreArray(a, &aa));
423 PetscFunctionReturn(PETSC_SUCCESS);
424 }
425
MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F)426 PETSC_SINGLE_LIBRARY_VISIBILITY_INTERNAL PetscErrorCode MatMFFDSetBase_MFFD(Mat J, Vec U, Vec F)
427 {
428 MatMFFD ctx;
429
430 PetscFunctionBegin;
431 PetscCall(MatShellGetContext(J, &ctx));
432 PetscCall(MatMFFDResetHHistory(J));
433 if (!ctx->current_u) {
434 PetscCall(VecDuplicate(U, &ctx->current_u));
435 PetscCall(VecLockReadPush(ctx->current_u));
436 }
437 PetscCall(VecLockReadPop(ctx->current_u));
438 PetscCall(VecCopy(U, ctx->current_u));
439 PetscCall(VecLockReadPush(ctx->current_u));
440 if (F) {
441 if (ctx->current_f_allocated) PetscCall(VecDestroy(&ctx->current_f));
442 ctx->current_f = F;
443 ctx->current_f_allocated = PETSC_FALSE;
444 } else if (!ctx->current_f_allocated) {
445 PetscCall(MatCreateVecs(J, NULL, &ctx->current_f));
446 ctx->current_f_allocated = PETSC_TRUE;
447 }
448 if (!ctx->w) PetscCall(VecDuplicate(ctx->current_u, &ctx->w));
449 J->assembled = PETSC_TRUE;
450 PetscFunctionReturn(PETSC_SUCCESS);
451 }
452
453 typedef PetscErrorCode (*FCN3)(void *, Vec, Vec, PetscScalar *); /* force argument to next function to not be extern C*/
MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void * ectx)454 static PetscErrorCode MatMFFDSetCheckh_MFFD(Mat J, FCN3 fun, void *ectx)
455 {
456 MatMFFD ctx;
457
458 PetscFunctionBegin;
459 PetscCall(MatShellGetContext(J, &ctx));
460 ctx->checkh = fun;
461 ctx->checkhctx = ectx;
462 PetscFunctionReturn(PETSC_SUCCESS);
463 }
464
465 /*@
466 MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all
467 MATMFFD` options in the database.
468
469 Collective
470
471 Input Parameters:
472 + mat - the `MATMFFD` context
473 - prefix - the prefix to prepend to all option names
474
475 Note:
476 A hyphen (-) must NOT be given at the beginning of the prefix name.
477 The first character of all runtime options is AUTOMATICALLY the hyphen.
478
479 Level: advanced
480
481 .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatSetFromOptions()`, `MatCreateSNESMF()`, `MatCreateMFFD()`
482 @*/
MatMFFDSetOptionsPrefix(Mat mat,const char prefix[])483 PetscErrorCode MatMFFDSetOptionsPrefix(Mat mat, const char prefix[])
484 {
485 MatMFFD mfctx;
486
487 PetscFunctionBegin;
488 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
489 PetscCall(MatShellGetContext(mat, &mfctx));
490 PetscValidHeaderSpecific(mfctx, MATMFFD_CLASSID, 1);
491 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mfctx, prefix));
492 PetscFunctionReturn(PETSC_SUCCESS);
493 }
494
MatSetFromOptions_MFFD(Mat mat,PetscOptionItems PetscOptionsObject)495 static PetscErrorCode MatSetFromOptions_MFFD(Mat mat, PetscOptionItems PetscOptionsObject)
496 {
497 MatMFFD mfctx;
498 PetscBool flg;
499 char ftype[256];
500
501 PetscFunctionBegin;
502 PetscCall(MatShellGetContext(mat, &mfctx));
503 PetscValidHeaderSpecific(mfctx, MATMFFD_CLASSID, 1);
504 PetscObjectOptionsBegin((PetscObject)mfctx);
505 PetscCall(PetscOptionsFList("-mat_mffd_type", "Matrix free type", "MatMFFDSetType", MatMFFDList, ((PetscObject)mfctx)->type_name, ftype, 256, &flg));
506 if (flg) PetscCall(MatMFFDSetType(mat, ftype));
507
508 PetscCall(PetscOptionsReal("-mat_mffd_err", "set sqrt relative error in function", "MatMFFDSetFunctionError", mfctx->error_rel, &mfctx->error_rel, NULL));
509 PetscCall(PetscOptionsInt("-mat_mffd_period", "how often h is recomputed", "MatMFFDSetPeriod", mfctx->recomputeperiod, &mfctx->recomputeperiod, NULL));
510
511 flg = PETSC_FALSE;
512 PetscCall(PetscOptionsBool("-mat_mffd_check_positivity", "Insure that U + h*a is nonnegative", "MatMFFDSetCheckh", flg, &flg, NULL));
513 if (flg) PetscCall(MatMFFDSetCheckh(mat, MatMFFDCheckPositivity, NULL));
514 #if defined(PETSC_USE_COMPLEX)
515 PetscCall(PetscOptionsBool("-mat_mffd_complex", "Use Lyness complex number trick to compute the matrix-vector product", "None", mfctx->usecomplex, &mfctx->usecomplex, NULL));
516 #endif
517 PetscTryTypeMethod(mfctx, setfromoptions, PetscOptionsObject);
518 PetscOptionsEnd();
519 PetscFunctionReturn(PETSC_SUCCESS);
520 }
521
MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period)522 static PetscErrorCode MatMFFDSetPeriod_MFFD(Mat mat, PetscInt period)
523 {
524 MatMFFD ctx;
525
526 PetscFunctionBegin;
527 PetscCall(MatShellGetContext(mat, &ctx));
528 ctx->recomputeperiod = period;
529 PetscFunctionReturn(PETSC_SUCCESS);
530 }
531
MatMFFDSetFunction_MFFD(Mat mat,MatMFFDFn * func,void * funcctx)532 static PetscErrorCode MatMFFDSetFunction_MFFD(Mat mat, MatMFFDFn *func, void *funcctx)
533 {
534 MatMFFD ctx;
535
536 PetscFunctionBegin;
537 PetscCall(MatShellGetContext(mat, &ctx));
538 ctx->func = func;
539 ctx->funcctx = funcctx;
540 PetscFunctionReturn(PETSC_SUCCESS);
541 }
542
MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error)543 static PetscErrorCode MatMFFDSetFunctionError_MFFD(Mat mat, PetscReal error)
544 {
545 PetscFunctionBegin;
546 if (error != (PetscReal)PETSC_DEFAULT) {
547 MatMFFD ctx;
548
549 PetscCall(MatShellGetContext(mat, &ctx));
550 ctx->error_rel = error;
551 }
552 PetscFunctionReturn(PETSC_SUCCESS);
553 }
554
MatMFFDSetHHistory_MFFD(Mat J,PetscScalar history[],PetscInt nhistory)555 static PetscErrorCode MatMFFDSetHHistory_MFFD(Mat J, PetscScalar history[], PetscInt nhistory)
556 {
557 MatMFFD ctx;
558
559 PetscFunctionBegin;
560 PetscCall(MatShellGetContext(J, &ctx));
561 ctx->historyh = history;
562 ctx->maxcurrenth = nhistory;
563 ctx->currenth = 0.;
564 PetscFunctionReturn(PETSC_SUCCESS);
565 }
566
567 /*MC
568 MATMFFD - "mffd" - A matrix-free matrix type.
569
570 Level: advanced
571
572 Developer Notes:
573 This is implemented on top of `MATSHELL` to get support for scaling and shifting without requiring duplicate code
574
575 Users should not MatShell... operations on this class, there is some error checking for that incorrect usage
576
577 .seealso: [](ch_matrices), `Mat`, `MatCreateMFFD()`, `MatCreateSNESMF()`, `MatMFFDSetFunction()`, `MatMFFDSetType()`,
578 `MatMFFDSetFunctionError()`, `MatMFFDDSSetUmin()`, `MatMFFDSetFunction()`
579 `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `MatCreateSNESMF()`,
580 `MatMFFDGetH()`,
581 M*/
MatCreate_MFFD(Mat A)582 PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A)
583 {
584 MatMFFD mfctx;
585
586 PetscFunctionBegin;
587 PetscCall(MatMFFDInitializePackage());
588
589 PetscCall(PetscHeaderCreate(mfctx, MATMFFD_CLASSID, "MatMFFD", "Matrix-free Finite Differencing", "Mat", PetscObjectComm((PetscObject)A), NULL, NULL));
590
591 mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON;
592 mfctx->recomputeperiod = 1;
593 mfctx->count = 0;
594 mfctx->currenth = 0.0;
595 mfctx->historyh = NULL;
596 mfctx->ncurrenth = 0;
597 mfctx->maxcurrenth = 0;
598 ((PetscObject)mfctx)->type_name = NULL;
599
600 /*
601 Create the empty data structure to contain compute-h routines.
602 These will be filled in below from the command line options or
603 a later call with MatMFFDSetType() or if that is not called
604 then it will default in the first use of MatMult_MFFD()
605 */
606 mfctx->ops->compute = NULL;
607 mfctx->ops->destroy = NULL;
608 mfctx->ops->view = NULL;
609 mfctx->ops->setfromoptions = NULL;
610 mfctx->hctx = NULL;
611
612 mfctx->func = NULL;
613 mfctx->funcctx = NULL;
614 mfctx->w = NULL;
615 mfctx->mat = A;
616
617 PetscCall(MatSetType(A, MATSHELL));
618 PetscCall(MatShellSetContext(A, mfctx));
619 PetscCall(MatShellSetOperation(A, MATOP_MULT, (PetscErrorCodeFn *)MatMult_MFFD));
620 PetscCall(MatShellSetOperation(A, MATOP_DESTROY, (PetscErrorCodeFn *)MatDestroy_MFFD));
621 PetscCall(MatShellSetOperation(A, MATOP_VIEW, (PetscErrorCodeFn *)MatView_MFFD));
622 PetscCall(MatShellSetOperation(A, MATOP_ASSEMBLY_END, (PetscErrorCodeFn *)MatAssemblyEnd_MFFD));
623 PetscCall(MatShellSetOperation(A, MATOP_SET_FROM_OPTIONS, (PetscErrorCodeFn *)MatSetFromOptions_MFFD));
624 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContext_C", MatShellSetContext_Immutable));
625 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContextDestroy_C", MatShellSetContextDestroy_Immutable));
626 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetManageScalingShifts_C", MatShellSetManageScalingShifts_Immutable));
627
628 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetBase_C", MatMFFDSetBase_MFFD));
629 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctioniBase_C", MatMFFDSetFunctioniBase_MFFD));
630 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctioni_C", MatMFFDSetFunctioni_MFFD));
631 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunction_C", MatMFFDSetFunction_MFFD));
632 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetCheckh_C", MatMFFDSetCheckh_MFFD));
633 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetPeriod_C", MatMFFDSetPeriod_MFFD));
634 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetFunctionError_C", MatMFFDSetFunctionError_MFFD));
635 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDResetHHistory_C", MatMFFDResetHHistory_MFFD));
636 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetHHistory_C", MatMFFDSetHHistory_MFFD));
637 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDSetType_C", MatMFFDSetType_MFFD));
638 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMFFDGetH_C", MatMFFDGetH_MFFD));
639 PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATMFFD));
640 PetscFunctionReturn(PETSC_SUCCESS);
641 }
642
643 /*@
644 MatCreateMFFD - Creates a matrix-free matrix of type `MATMFFD` that uses finite differences on a provided function to
645 approximately multiply a vector by the matrix (Jacobian) . See also `MatCreateSNESMF()`
646
647 Collective
648
649 Input Parameters:
650 + comm - MPI communicator
651 . m - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given)
652 This value should be the same as the local size used in creating the
653 y vector for the matrix-vector product y = Ax.
654 . n - This value should be the same as the local size used in creating the
655 x vector for the matrix-vector product y = Ax. (or `PETSC_DECIDE` to have
656 calculated if `N` is given) For square matrices `n` is almost always `m`.
657 . M - number of global rows (or `PETSC_DETERMINE` to have calculated if `m` is given)
658 - N - number of global columns (or `PETSC_DETERMINE` to have calculated if `n` is given)
659
660 Output Parameter:
661 . J - the matrix-free matrix
662
663 Options Database Keys:
664 + -mat_mffd_type - wp or ds (see `MATMFFD_WP` or `MATMFFD_DS`)
665 . -mat_mffd_err - square root of estimated relative error in function evaluation
666 . -mat_mffd_period - how often h is recomputed, defaults to 1, every time
667 . -mat_mffd_check_positivity - possibly decrease `h` until U + h*a has only positive values
668 . -mat_mffd_umin <umin> - Sets umin (for default PETSc routine that computes h only)
669 . -mat_mffd_complex - use the Lyness trick with complex numbers to compute the matrix-vector product instead of differencing
670 (requires real valued functions but that PETSc be configured for complex numbers)
671 . -snes_mf - use the finite difference based matrix-free matrix with `SNESSolve()` and no preconditioner
672 - -snes_mf_operator - use the finite difference based matrix-free matrix with `SNESSolve()` but construct a preconditioner
673 using the matrix passed as `pmat` to `SNESSetJacobian()`.
674
675 Level: advanced
676
677 Notes:
678 Use `MatMFFDSetFunction()` to provide the function that will be differenced to compute the matrix-vector product.
679
680 The matrix-free matrix context contains the function pointers
681 and work space for performing finite difference approximations of
682 Jacobian-vector products, F'(u)*a,
683
684 The default code uses the following approach to compute h
685
686 .vb
687 F'(u)*a = [F(u+h*a) - F(u)]/h where
688 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1}
689 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise
690 where
691 error_rel = square root of relative error in function evaluation
692 umin = minimum iterate parameter
693 .ve
694
695 To have `SNES` use the matrix-free finite difference matrix-vector product and not provide a separate matrix
696 from which to compute the preconditioner (the `pmat` argument `SNESSetJacobian()`), then simply call `SNESSetJacobian()`
697 with `NULL` for the matrices and `MatMFFDComputeJacobian()`. Or use the options database option `-snes_mf`
698
699 The user can set `error_rel` via `MatMFFDSetFunctionError()` and `umin` via `MatMFFDDSSetUmin()`.
700
701 Use `MATSHELL` or `MatCreateShell()` to provide your own custom matrix-vector operation.
702
703 .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatDestroy()`, `MatMFFDSetFunctionError()`, `MatMFFDDSSetUmin()`, `MatMFFDSetFunction()`
704 `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `MatCreateSNESMF()`, `MatCreateShell()`, `MATSHELL`,
705 `MatMFFDGetH()`, `MatMFFDRegister()`, `MatMFFDComputeJacobian()`
706 @*/
MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat * J)707 PetscErrorCode MatCreateMFFD(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, Mat *J)
708 {
709 PetscFunctionBegin;
710 PetscCall(MatCreate(comm, J));
711 PetscCall(MatSetSizes(*J, m, n, M, N));
712 PetscCall(MatSetType(*J, MATMFFD));
713 PetscCall(MatSetUp(*J));
714 PetscFunctionReturn(PETSC_SUCCESS);
715 }
716
717 /*@
718 MatMFFDGetH - Gets the last value that was used as the differencing for a `MATMFFD` matrix
719 parameter.
720
721 Not Collective
722
723 Input Parameters:
724 . mat - the `MATMFFD` matrix
725
726 Output Parameter:
727 . h - the differencing step size
728
729 Level: advanced
730
731 .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDSetHHistory()`, `MatCreateMFFD()`, `MatMFFDResetHHistory()`
732 @*/
MatMFFDGetH(Mat mat,PetscScalar * h)733 PetscErrorCode MatMFFDGetH(Mat mat, PetscScalar *h)
734 {
735 PetscFunctionBegin;
736 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
737 PetscAssertPointer(h, 2);
738 PetscUseMethod(mat, "MatMFFDGetH_C", (Mat, PetscScalar *), (mat, h));
739 PetscFunctionReturn(PETSC_SUCCESS);
740 }
741
742 /*@C
743 MatMFFDSetFunction - Sets the function used in applying the matrix-free `MATMFFD` matrix.
744
745 Logically Collective
746
747 Input Parameters:
748 + mat - the matrix-free matrix `MATMFFD` created via `MatCreateSNESMF()` or `MatCreateMFFD()`
749 . func - the function to use
750 - funcctx - optional function context passed to function
751
752 Level: advanced
753
754 Notes:
755 If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix-free
756 matrix inside your compute Jacobian routine
757
758 If this is not set then it will use the function set with `SNESSetFunction()` if `MatCreateSNESMF()` was used.
759
760 .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDFn`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`,
761 `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESSetFunction()`
762 @*/
MatMFFDSetFunction(Mat mat,MatMFFDFn * func,void * funcctx)763 PetscErrorCode MatMFFDSetFunction(Mat mat, MatMFFDFn *func, void *funcctx)
764 {
765 PetscFunctionBegin;
766 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
767 PetscTryMethod(mat, "MatMFFDSetFunction_C", (Mat, MatMFFDFn *, void *), (mat, func, funcctx));
768 PetscFunctionReturn(PETSC_SUCCESS);
769 }
770
771 /*@C
772 MatMFFDSetFunctioni - Sets the function for computing a single component for a `MATMFFD` matrix
773
774 Logically Collective
775
776 Input Parameters:
777 + mat - the matrix-free matrix `MATMFFD`
778 - funci - the function to use
779
780 Level: advanced
781
782 Notes:
783 If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix-free
784 matrix inside your compute Jacobian routine.
785
786 This function is necessary to compute the diagonal of the matrix.
787 `funci` must not contain any MPI call as it is called inside a loop on the local portion of the vector.
788
789 .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDiFn`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`,
790 `SNESSetFunction()`, `MatGetDiagonal()`
791 @*/
MatMFFDSetFunctioni(Mat mat,MatMFFDiFn * funci)792 PetscErrorCode MatMFFDSetFunctioni(Mat mat, MatMFFDiFn *funci)
793 {
794 PetscFunctionBegin;
795 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
796 PetscTryMethod(mat, "MatMFFDSetFunctioni_C", (Mat, MatMFFDiFn *), (mat, funci));
797 PetscFunctionReturn(PETSC_SUCCESS);
798 }
799
800 /*@C
801 MatMFFDSetFunctioniBase - Sets the function to compute the base vector for a single component function evaluation for a `MATMFFD` matrix
802
803 Logically Collective
804
805 Input Parameters:
806 + mat - the `MATMFFD` matrix-free matrix
807 - func - the function to use
808
809 Level: advanced
810
811 Notes:
812 If you use this you MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` on the matrix-free
813 matrix inside your compute Jacobian routine.
814
815 This function is necessary to compute the diagonal of the matrix, used for example with `PCJACOBI`
816
817 .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`,
818 `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESSetFunction()`, `MatGetDiagonal()`
819 @*/
MatMFFDSetFunctioniBase(Mat mat,MatMFFDiBaseFn * func)820 PetscErrorCode MatMFFDSetFunctioniBase(Mat mat, MatMFFDiBaseFn *func)
821 {
822 PetscFunctionBegin;
823 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
824 PetscTryMethod(mat, "MatMFFDSetFunctioniBase_C", (Mat, MatMFFDiBaseFn *), (mat, func));
825 PetscFunctionReturn(PETSC_SUCCESS);
826 }
827
828 /*@
829 MatMFFDSetPeriod - Sets how often the step-size `h` is recomputed for a `MATMFFD` matrix, by default it is every time
830
831 Logically Collective
832
833 Input Parameters:
834 + mat - the `MATMFFD` matrix-free matrix
835 - period - 1 for every time, 2 for every second etc
836
837 Options Database Key:
838 . -mat_mffd_period <period> - Sets how often `h` is recomputed
839
840 Level: advanced
841
842 .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`,
843 `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`
844 @*/
MatMFFDSetPeriod(Mat mat,PetscInt period)845 PetscErrorCode MatMFFDSetPeriod(Mat mat, PetscInt period)
846 {
847 PetscFunctionBegin;
848 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
849 PetscValidLogicalCollectiveInt(mat, period, 2);
850 PetscTryMethod(mat, "MatMFFDSetPeriod_C", (Mat, PetscInt), (mat, period));
851 PetscFunctionReturn(PETSC_SUCCESS);
852 }
853
854 /*@
855 MatMFFDSetFunctionError - Sets the error_rel for the approximation of matrix-vector products using finite differences with the `MATMFFD` matrix
856
857 Logically Collective
858
859 Input Parameters:
860 + mat - the `MATMFFD` matrix-free matrix
861 - error - relative error (should be set to the square root of the relative error in the function evaluations)
862
863 Options Database Key:
864 . -mat_mffd_err <error_rel> - Sets error_rel
865
866 Level: advanced
867
868 Note:
869 The default matrix-free matrix-vector product routine computes
870 .vb
871 F'(u)*a = [F(u+h*a) - F(u)]/h where
872 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1}
873 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else
874 .ve
875
876 .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`,
877 `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`
878 @*/
MatMFFDSetFunctionError(Mat mat,PetscReal error)879 PetscErrorCode MatMFFDSetFunctionError(Mat mat, PetscReal error)
880 {
881 PetscFunctionBegin;
882 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
883 PetscValidLogicalCollectiveReal(mat, error, 2);
884 PetscTryMethod(mat, "MatMFFDSetFunctionError_C", (Mat, PetscReal), (mat, error));
885 PetscFunctionReturn(PETSC_SUCCESS);
886 }
887
888 /*@
889 MatMFFDSetHHistory - Sets an array to collect a history of the
890 differencing values (h) computed for the matrix-free product `MATMFFD` matrix
891
892 Logically Collective
893
894 Input Parameters:
895 + J - the `MATMFFD` matrix-free matrix
896 . history - space to hold the history
897 - nhistory - number of entries in history, if more entries are generated than
898 nhistory, then the later ones are discarded
899
900 Level: advanced
901
902 Note:
903 Use `MatMFFDResetHHistory()` to reset the history counter and collect
904 a new batch of differencing parameters, h.
905
906 .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDGetH()`, `MatCreateSNESMF()`,
907 `MatMFFDResetHHistory()`, `MatMFFDSetFunctionError()`
908 @*/
MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory)909 PetscErrorCode MatMFFDSetHHistory(Mat J, PetscScalar history[], PetscInt nhistory)
910 {
911 PetscFunctionBegin;
912 PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
913 if (history) PetscAssertPointer(history, 2);
914 PetscValidLogicalCollectiveInt(J, nhistory, 3);
915 PetscUseMethod(J, "MatMFFDSetHHistory_C", (Mat, PetscScalar[], PetscInt), (J, history, nhistory));
916 PetscFunctionReturn(PETSC_SUCCESS);
917 }
918
919 /*@
920 MatMFFDResetHHistory - Resets the counter to zero to begin
921 collecting a new set of differencing histories for the `MATMFFD` matrix
922
923 Logically Collective
924
925 Input Parameter:
926 . J - the matrix-free matrix context
927
928 Level: advanced
929
930 Note:
931 Use `MatMFFDSetHHistory()` to create the original history counter.
932
933 .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDGetH()`, `MatCreateSNESMF()`,
934 `MatMFFDSetHHistory()`, `MatMFFDSetFunctionError()`
935 @*/
MatMFFDResetHHistory(Mat J)936 PetscErrorCode MatMFFDResetHHistory(Mat J)
937 {
938 PetscFunctionBegin;
939 PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
940 PetscTryMethod(J, "MatMFFDResetHHistory_C", (Mat), (J));
941 PetscFunctionReturn(PETSC_SUCCESS);
942 }
943
944 /*@
945 MatMFFDSetBase - Sets the vector `U` at which matrix vector products of the
946 Jacobian are computed for the `MATMFFD` matrix
947
948 Logically Collective
949
950 Input Parameters:
951 + J - the `MATMFFD` matrix
952 . U - the vector
953 - F - (optional) vector that contains F(u) if it has been already computed
954
955 Level: advanced
956
957 Notes:
958 This is rarely used directly
959
960 If `F` is provided then it is not recomputed. Otherwise the function is evaluated at the base
961 point during the first `MatMult()` after each call to `MatMFFDSetBase()`.
962
963 .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMult()`
964 @*/
MatMFFDSetBase(Mat J,Vec U,Vec F)965 PetscErrorCode MatMFFDSetBase(Mat J, Vec U, Vec F)
966 {
967 PetscFunctionBegin;
968 PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
969 PetscValidHeaderSpecific(U, VEC_CLASSID, 2);
970 if (F) PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
971 PetscTryMethod(J, "MatMFFDSetBase_C", (Mat, Vec, Vec), (J, U, F));
972 PetscFunctionReturn(PETSC_SUCCESS);
973 }
974
975 /*@C
976 MatMFFDSetCheckh - Sets a function that checks the computed `h` and adjusts
977 it to satisfy some criteria for the `MATMFFD` matrix
978
979 Logically Collective
980
981 Input Parameters:
982 + J - the `MATMFFD` matrix
983 . fun - the function that checks `h`, see `MatMFFDCheckhFn`
984 - ctx - any context needed by the function
985
986 Options Database Keys:
987 . -mat_mffd_check_positivity <bool> - Ensure that $U + h*a $ is non-negative
988
989 Level: advanced
990
991 Notes:
992 For example, `MatMFFDCheckPositivity()` insures that all entries of U + h*a are non-negative
993
994 The function you provide is called after the default `h` has been computed and allows you to
995 modify it.
996
997 .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDCheckhFn`, `MatMFFDCheckPositivity()`
998 @*/
MatMFFDSetCheckh(Mat J,MatMFFDCheckhFn * fun,PetscCtx ctx)999 PetscErrorCode MatMFFDSetCheckh(Mat J, MatMFFDCheckhFn *fun, PetscCtx ctx)
1000 {
1001 PetscFunctionBegin;
1002 PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
1003 PetscTryMethod(J, "MatMFFDSetCheckh_C", (Mat, MatMFFDCheckhFn *, void *), (J, fun, ctx));
1004 PetscFunctionReturn(PETSC_SUCCESS);
1005 }
1006
1007 /*@
1008 MatMFFDCheckPositivity - Checks that all entries in $U + h*a $ are positive or
1009 zero, decreases `h` until this is satisfied for a `MATMFFD` matrix
1010
1011 Logically Collective
1012
1013 Input Parameters:
1014 + dummy - context variable (unused)
1015 . U - base vector that is added to
1016 . a - vector that is added
1017 - h - scaling factor on `a`, may be changed on output
1018
1019 Options Database Keys:
1020 . -mat_mffd_check_positivity <bool> - Ensure that $U + h*a$ is nonnegative
1021
1022 Level: advanced
1023
1024 Note:
1025 This is rarely used directly, rather it is passed as an argument to `MatMFFDSetCheckh()`
1026
1027 .seealso: [](ch_matrices), `Mat`, `MATMFFD`, `MatMFFDSetCheckh()`
1028 @*/
MatMFFDCheckPositivity(void * dummy,Vec U,Vec a,PetscScalar * h)1029 PetscErrorCode MatMFFDCheckPositivity(void *dummy, Vec U, Vec a, PetscScalar *h)
1030 {
1031 PetscReal val, minval;
1032 PetscScalar *u_vec, *a_vec;
1033 PetscInt i, n;
1034 MPI_Comm comm;
1035
1036 PetscFunctionBegin;
1037 PetscValidHeaderSpecific(U, VEC_CLASSID, 2);
1038 PetscValidHeaderSpecific(a, VEC_CLASSID, 3);
1039 PetscAssertPointer(h, 4);
1040 PetscCall(PetscObjectGetComm((PetscObject)U, &comm));
1041 PetscCall(VecGetArray(U, &u_vec));
1042 PetscCall(VecGetArray(a, &a_vec));
1043 PetscCall(VecGetLocalSize(U, &n));
1044 minval = PetscAbsScalar(*h) * PetscRealConstant(1.01);
1045 for (i = 0; i < n; i++) {
1046 if (PetscRealPart(u_vec[i] + *h * a_vec[i]) <= 0.0) {
1047 val = PetscAbsScalar(u_vec[i] / a_vec[i]);
1048 if (val < minval) minval = val;
1049 }
1050 }
1051 PetscCall(VecRestoreArray(U, &u_vec));
1052 PetscCall(VecRestoreArray(a, &a_vec));
1053 PetscCallMPI(MPIU_Allreduce(&minval, &val, 1, MPIU_REAL, MPIU_MIN, comm));
1054 if (val <= PetscAbsScalar(*h)) {
1055 val = 0.99 * val;
1056 PetscCall(PetscInfo(U, "Scaling back h from %g to %g\n", (double)PetscRealPart(*h), (double)val));
1057 if (PetscRealPart(*h) > 0.0) *h = val;
1058 else *h = -val;
1059 }
1060 PetscFunctionReturn(PETSC_SUCCESS);
1061 }
1062