Actual source code: gcreate.c
1: #include <petsc/private/matimpl.h>
3: #include <../src/mat/impls/aij/seq/aij.h>
4: #include <../src/mat/impls/aij/mpi/mpiaij.h>
6: PetscErrorCode MatSetBlockSizes_Default(Mat mat, PetscInt rbs, PetscInt cbs)
7: {
8: PetscFunctionBegin;
9: if (!mat->preallocated) PetscFunctionReturn(PETSC_SUCCESS);
10: PetscCheck(mat->rmap->bs <= 0 || mat->rmap->bs == rbs, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot change row block size %" PetscInt_FMT " to %" PetscInt_FMT, mat->rmap->bs, rbs);
11: PetscCheck(mat->cmap->bs <= 0 || mat->cmap->bs == cbs, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot change column block size %" PetscInt_FMT " to %" PetscInt_FMT, mat->cmap->bs, cbs);
12: PetscFunctionReturn(PETSC_SUCCESS);
13: }
15: PetscErrorCode MatShift_Basic(Mat Y, PetscScalar a)
16: {
17: PetscInt i, start, end, oldValA = 0, oldValB = 0;
18: PetscScalar alpha = a;
19: PetscBool prevoption;
20: PetscBool isSeqAIJDerived, isMPIAIJDerived; // all classes sharing SEQAIJHEADER or MPIAIJHEADER
21: Mat A = NULL, B = NULL;
23: PetscFunctionBegin;
24: PetscCall(MatGetOption(Y, MAT_NO_OFF_PROC_ENTRIES, &prevoption));
25: PetscCall(MatSetOption(Y, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
26: PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &isSeqAIJDerived, MATSEQAIJ, MATSEQBAIJ, MATSEQSBAIJ, ""));
27: PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &isMPIAIJDerived, MATMPIAIJ, MATMPIBAIJ, MATMPISBAIJ, ""));
29: if (isSeqAIJDerived) A = Y;
30: else if (isMPIAIJDerived) {
31: Mat_MPIAIJ *mpiaij = (Mat_MPIAIJ *)Y->data;
32: A = mpiaij->A;
33: B = mpiaij->B;
34: }
36: if (A) {
37: oldValA = ((Mat_SeqAIJ *)A->data)->nonew;
38: ((Mat_SeqAIJ *)A->data)->nonew = 0; // so that new nonzero locations are allowed
39: }
40: if (B) {
41: oldValB = ((Mat_SeqAIJ *)B->data)->nonew;
42: ((Mat_SeqAIJ *)B->data)->nonew = 0;
43: }
45: PetscCall(MatGetOwnershipRange(Y, &start, &end));
46: for (i = start; i < end; i++) {
47: if (i < Y->cmap->N) PetscCall(MatSetValues(Y, 1, &i, 1, &i, &alpha, ADD_VALUES));
48: }
49: PetscCall(MatAssemblyBegin(Y, MAT_FINAL_ASSEMBLY));
50: PetscCall(MatAssemblyEnd(Y, MAT_FINAL_ASSEMBLY));
51: PetscCall(MatSetOption(Y, MAT_NO_OFF_PROC_ENTRIES, prevoption));
52: if (A) ((Mat_SeqAIJ *)A->data)->nonew = oldValA;
53: if (B) ((Mat_SeqAIJ *)B->data)->nonew = oldValB;
54: PetscFunctionReturn(PETSC_SUCCESS);
55: }
57: /*@
58: MatCreate - Creates a matrix where the type is determined
59: from either a call to `MatSetType()` or from the options database
60: with a call to `MatSetFromOptions()`.
62: Collective
64: Input Parameter:
65: . comm - MPI communicator
67: Output Parameter:
68: . A - the matrix
70: Options Database Keys:
71: + -mat_type seqaij - `MATSEQAIJ` type, uses `MatCreateSeqAIJ()`
72: . -mat_type mpiaij - `MATMPIAIJ` type, uses `MatCreateAIJ()`
73: . -mat_type seqdense - `MATSEQDENSE`, uses `MatCreateSeqDense()`
74: . -mat_type mpidense - `MATMPIDENSE` type, uses `MatCreateDense()`
75: . -mat_type seqbaij - `MATSEQBAIJ` type, uses `MatCreateSeqBAIJ()`
76: - -mat_type mpibaij - `MATMPIBAIJ` type, uses `MatCreateBAIJ()`
78: See the manpages for particular formats (e.g., `MATSEQAIJ`)
79: for additional format-specific options.
81: Level: beginner
83: Notes:
84: The default matrix type is `MATAIJ`, using the routines `MatCreateSeqAIJ()` or
85: `MatCreateAIJ()` if you do not set a type in the options database. If you never call
86: `MatSetType()` or `MatSetFromOptions()` it will generate an error when you try to use the
87: matrix.
89: .seealso: [](ch_matrices), `Mat`, `MatCreateSeqAIJ()`, `MatCreateAIJ()`,
90: `MatCreateSeqDense()`, `MatCreateDense()`,
91: `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`,
92: `MatCreateSeqSBAIJ()`, `MatCreateSBAIJ()`,
93: `MatConvert()`
94: @*/
95: PetscErrorCode MatCreate(MPI_Comm comm, Mat *A)
96: {
97: Mat B;
99: PetscFunctionBegin;
100: PetscAssertPointer(A, 2);
101: PetscCall(MatInitializePackage());
103: PetscCall(PetscHeaderCreate(B, MAT_CLASSID, "Mat", "Matrix", "Mat", comm, MatDestroy, MatView));
104: PetscCall(PetscLayoutCreate(comm, &B->rmap));
105: PetscCall(PetscLayoutCreate(comm, &B->cmap));
106: PetscCall(PetscStrallocpy(VECSTANDARD, &B->defaultvectype));
107: PetscCall(PetscStrallocpy(PETSCRANDER48, &B->defaultrandtype));
109: B->symmetric = PETSC_BOOL3_UNKNOWN;
110: B->hermitian = PETSC_BOOL3_UNKNOWN;
111: B->structurally_symmetric = PETSC_BOOL3_UNKNOWN;
112: B->spd = PETSC_BOOL3_UNKNOWN;
113: B->symmetry_eternal = PETSC_FALSE;
114: B->structural_symmetry_eternal = PETSC_FALSE;
116: B->congruentlayouts = PETSC_DECIDE;
117: B->preallocated = PETSC_FALSE;
118: #if defined(PETSC_HAVE_DEVICE)
119: B->boundtocpu = PETSC_TRUE;
120: #endif
121: *A = B;
122: PetscFunctionReturn(PETSC_SUCCESS);
123: }
125: /*@
126: MatCreateFromOptions - Creates a matrix whose type is set from the options database
128: Collective
130: Input Parameters:
131: + comm - MPI communicator
132: . prefix - [optional] prefix for the options database
133: . bs - the blocksize (commonly 1)
134: . m - the local number of rows (or `PETSC_DECIDE`)
135: . n - the local number of columns (or `PETSC_DECIDE` or `PETSC_DETERMINE`)
136: . M - the global number of rows (or `PETSC_DETERMINE`)
137: - N - the global number of columns (or `PETSC_DETERMINE`)
139: Output Parameter:
140: . A - the matrix
142: Options Database Key:
143: . -mat_type - see `MatType`, for example `aij`, `aijcusparse`, `baij`, `sbaij`, `dense`, defaults to `aij`
145: Level: beginner
147: .seealso: [](ch_matrices), `Mat`, `MatCreateSeqAIJ()`, `MatCreateAIJ()`,
148: `MatCreateSeqDense()`, `MatCreateDense()`,
149: `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`,
150: `MatCreateSeqSBAIJ()`, `MatCreateSBAIJ()`,
151: `MatConvert()`, `MatCreate()`
152: @*/
153: PetscErrorCode MatCreateFromOptions(MPI_Comm comm, const char *prefix, PetscInt bs, PetscInt m, PetscInt n, PetscInt M, PetscInt N, Mat *A)
154: {
155: PetscFunctionBegin;
156: PetscAssertPointer(A, 8);
157: PetscCall(MatCreate(comm, A));
158: if (prefix) PetscCall(MatSetOptionsPrefix(*A, prefix));
159: PetscCall(MatSetBlockSize(*A, bs));
160: PetscCall(MatSetSizes(*A, m, n, M, N));
161: PetscCall(MatSetFromOptions(*A));
162: PetscFunctionReturn(PETSC_SUCCESS);
163: }
165: /*@
166: MatSetErrorIfFailure - Causes `Mat` to generate an immediate error, for example a zero pivot, is detected.
168: Logically Collective
170: Input Parameters:
171: + mat - matrix obtained from `MatCreate()`
172: - flg - `PETSC_TRUE` indicates you want the error generated
174: Level: advanced
176: Note:
177: If this flag is not set then the matrix operation will note the error and continue. The error may cause a later `PC` or `KSP` error
178: or result in a `KSPConvergedReason` indicating the method did not converge.
180: .seealso: [](ch_matrices), `Mat`, `PCSetErrorIfFailure()`, `KSPConvergedReason`, `SNESConvergedReason`
181: @*/
182: PetscErrorCode MatSetErrorIfFailure(Mat mat, PetscBool flg)
183: {
184: PetscFunctionBegin;
187: mat->erroriffailure = flg;
188: PetscFunctionReturn(PETSC_SUCCESS);
189: }
191: /*@
192: MatSetSizes - Sets the local and global sizes, and checks to determine compatibility
194: Collective
196: Input Parameters:
197: + A - the matrix
198: . m - number of local rows (or `PETSC_DECIDE`)
199: . n - number of local columns (or `PETSC_DECIDE`)
200: . M - number of global rows (or `PETSC_DETERMINE`)
201: - N - number of global columns (or `PETSC_DETERMINE`)
203: Level: beginner
205: Notes:
206: `m` (`n`) and `M` (`N`) cannot be both `PETSC_DECIDE`
207: If one processor calls this with `M` (`N`) of `PETSC_DECIDE` then all processors must, otherwise the program will hang.
209: If `PETSC_DECIDE` is not used for the arguments 'm' and 'n', then the
210: user must ensure that they are chosen to be compatible with the
211: vectors. To do this, one first considers the matrix-vector product
212: 'y = A x'. The `m` that is used in the above routine must match the
213: local size of 'y'. Likewise, the `n` used must match the local size of 'x'.
215: If `m` and `n` are not `PETSC_DECIDE`, then the values determine the `PetscLayout` of the matrix and the ranges returned by
216: `MatGetOwnershipRange()`, `MatGetOwnershipRanges()`, `MatGetOwnershipRangeColumn()`, and `MatGetOwnershipRangesColumn()`.
218: You cannot change the sizes once they have been set.
220: The sizes must be set before `MatSetUp()` or MatXXXSetPreallocation() is called.
222: .seealso: [](ch_matrices), `Mat`, `MatGetSize()`, `PetscSplitOwnership()`, `MatGetOwnershipRange()`, `MatGetOwnershipRanges()`,
223: `MatGetOwnershipRangeColumn()`, `MatGetOwnershipRangesColumn()`, `PetscLayout`, `VecSetSizes()`
224: @*/
225: PetscErrorCode MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N)
226: {
227: PetscFunctionBegin;
231: PetscCheck(M <= 0 || m <= M, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Local row size %" PetscInt_FMT " cannot be larger than global row size %" PetscInt_FMT, m, M);
232: PetscCheck(N <= 0 || n <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Local column size %" PetscInt_FMT " cannot be larger than global column size %" PetscInt_FMT, n, N);
233: PetscCheck((A->rmap->n < 0 || A->rmap->N < 0) || (A->rmap->n == m && (M <= 0 || A->rmap->N == M)), PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot change/reset row sizes to %" PetscInt_FMT " local %" PetscInt_FMT " global after previously setting them to %" PetscInt_FMT " local %" PetscInt_FMT " global", m, M,
234: A->rmap->n, A->rmap->N);
235: PetscCheck((A->cmap->n < 0 || A->cmap->N < 0) || (A->cmap->n == n && (N <= 0 || A->cmap->N == N)), PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot change/reset column sizes to %" PetscInt_FMT " local %" PetscInt_FMT " global after previously setting them to %" PetscInt_FMT " local %" PetscInt_FMT " global", n, N,
236: A->cmap->n, A->cmap->N);
237: A->rmap->n = m;
238: A->cmap->n = n;
239: A->rmap->N = M > -1 ? M : A->rmap->N;
240: A->cmap->N = N > -1 ? N : A->cmap->N;
241: PetscFunctionReturn(PETSC_SUCCESS);
242: }
244: /*@
245: MatSetFromOptions - Creates a matrix where the type is determined
246: from the options database.
248: Collective
250: Input Parameter:
251: . B - the matrix
253: Options Database Keys:
254: + -mat_type seqaij - `MATSEQAIJ` type, uses `MatCreateSeqAIJ()`
255: . -mat_type mpiaij - `MATMPIAIJ` type, uses `MatCreateAIJ()`
256: . -mat_type seqdense - `MATSEQDENSE` type, uses `MatCreateSeqDense()`
257: . -mat_type mpidense - `MATMPIDENSE`, uses `MatCreateDense()`
258: . -mat_type seqbaij - `MATSEQBAIJ`, uses `MatCreateSeqBAIJ()`
259: . -mat_type mpibaij - `MATMPIBAIJ`, uses `MatCreateBAIJ()`
260: - -mat_vec_type - the `VecType` used by `MatCreateVecs()`, see `MatSetVecType()`
262: See the manpages for particular formats (e.g., `MATSEQAIJ`)
263: for additional format-specific options.
265: Level: beginner
267: Notes:
268: Generates a parallel MPI matrix if the communicator has more than one processor. The default
269: matrix type is `MATAIJ`, using the routines `MatCreateSeqAIJ()` and `MatCreateAIJ()` if you
270: do not select a type in the options database.
272: .seealso: [](ch_matrices), `Mat`, `MatCreateSeqAIJ()`, `MatCreateAIJ()`,
273: `MatCreateSeqDense()`, `MatCreateDense()`,
274: `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`,
275: `MatCreateSeqSBAIJ()`, `MatCreateSBAIJ()`,
276: `MatConvert()`
277: @*/
278: PetscErrorCode MatSetFromOptions(Mat B)
279: {
280: const char *deft = MATAIJ;
281: char type[256];
282: PetscBool flg, set;
283: PetscInt bind_below = 0, newbs = -1;
285: PetscFunctionBegin;
288: PetscObjectOptionsBegin((PetscObject)B);
290: PetscCall(PetscOptionsInt("-mat_block_size", "Set the blocksize used to store the matrix", "MatSetBlockSize", newbs, &newbs, &flg));
291: if (flg) {
292: PetscCall(PetscLayoutSetBlockSize(B->rmap, newbs));
293: PetscCall(PetscLayoutSetBlockSize(B->cmap, newbs));
294: }
296: PetscCall(PetscOptionsFList("-mat_type", "Matrix type", "MatSetType", MatList, deft, type, PETSC_STATIC_ARRAY_LENGTH(type), &flg));
297: if (flg) PetscCall(MatSetType(B, type));
298: else if (!((PetscObject)B)->type_name) PetscCall(MatSetType(B, deft));
299: PetscCall(PetscOptionsFList("-mat_vec_type", "Vec type used by MatCreateVecs()", "MatSetVecType", VecList, B->defaultvectype, type, PETSC_STATIC_ARRAY_LENGTH(type), &flg));
300: if (flg) PetscCall(MatSetVecType(B, type));
302: if (newbs > 0) PetscTryTypeMethod(B, setblocksizes, newbs, newbs);
304: PetscCall(PetscOptionsName("-mat_is_symmetric", "Checks if mat is symmetric on MatAssemblyEnd()", "MatIsSymmetric", &B->checksymmetryonassembly));
305: PetscCall(PetscOptionsReal("-mat_is_symmetric", "Checks if mat is symmetric on MatAssemblyEnd()", "MatIsSymmetric", B->checksymmetrytol, &B->checksymmetrytol, NULL));
306: PetscCall(PetscOptionsBool("-mat_null_space_test", "Checks if provided null space is correct in MatAssemblyEnd()", "MatSetNullSpaceTest", B->checknullspaceonassembly, &B->checknullspaceonassembly, NULL));
307: PetscCall(PetscOptionsBool("-mat_error_if_failure", "Generate an error if an error occurs when factoring the matrix", "MatSetErrorIfFailure", B->erroriffailure, &B->erroriffailure, NULL));
309: PetscTryTypeMethod(B, setfromoptions, PetscOptionsObject);
311: flg = PETSC_FALSE;
312: PetscCall(PetscOptionsBool("-mat_new_nonzero_location_err", "Generate an error if new nonzeros are created in the matrix nonzero structure (useful to test preallocation)", "MatSetOption", flg, &flg, &set));
313: if (set) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, flg));
314: flg = PETSC_FALSE;
315: PetscCall(PetscOptionsBool("-mat_new_nonzero_allocation_err", "Generate an error if new nonzeros are allocated in the matrix nonzero structure (useful to test preallocation)", "MatSetOption", flg, &flg, &set));
316: if (set) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, flg));
317: flg = PETSC_FALSE;
318: PetscCall(PetscOptionsBool("-mat_ignore_zero_entries", "For AIJ/IS matrices this will stop zero values from creating a zero location in the matrix", "MatSetOption", flg, &flg, &set));
319: if (set) PetscCall(MatSetOption(B, MAT_IGNORE_ZERO_ENTRIES, flg));
321: flg = PETSC_FALSE;
322: PetscCall(PetscOptionsBool("-mat_form_explicit_transpose", "Hint to form an explicit transpose for operations like MatMultTranspose", "MatSetOption", flg, &flg, &set));
323: if (set) PetscCall(MatSetOption(B, MAT_FORM_EXPLICIT_TRANSPOSE, flg));
325: /* Bind to CPU if below a user-specified size threshold.
326: * This perhaps belongs in the options for the GPU Mat types, but MatBindToCPU() does nothing when called on non-GPU types,
327: * and putting it here makes is more maintainable than duplicating this for all. */
328: PetscCall(PetscOptionsInt("-mat_bind_below", "Set the size threshold (in local rows) below which the Mat is bound to the CPU", "MatBindToCPU", bind_below, &bind_below, &flg));
329: if (flg && B->rmap->n < bind_below) PetscCall(MatBindToCPU(B, PETSC_TRUE));
331: /* process any options handlers added with PetscObjectAddOptionsHandler() */
332: PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)B, PetscOptionsObject));
333: PetscOptionsEnd();
334: PetscFunctionReturn(PETSC_SUCCESS);
335: }
337: /*@
338: MatXAIJSetPreallocation - set preallocation for serial and parallel `MATAIJ`, `MATBAIJ`, and `MATSBAIJ` matrices and their unassembled versions.
340: Collective
342: Input Parameters:
343: + A - matrix being preallocated
344: . bs - block size
345: . dnnz - number of nonzero column blocks per block row of diagonal part of parallel matrix
346: . onnz - number of nonzero column blocks per block row of off-diagonal part of parallel matrix
347: . dnnzu - number of nonzero column blocks per block row of upper-triangular part of diagonal part of parallel matrix
348: - onnzu - number of nonzero column blocks per block row of upper-triangular part of off-diagonal part of parallel matrix
350: Level: beginner
352: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, `MatMPIBAIJSetPreallocation()`,
353: `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`,
354: `PetscSplitOwnership()`
355: @*/
356: PetscErrorCode MatXAIJSetPreallocation(Mat A, PetscInt bs, const PetscInt dnnz[], const PetscInt onnz[], const PetscInt dnnzu[], const PetscInt onnzu[])
357: {
358: PetscInt cbs;
359: PetscBool aij, is, hyp;
361: PetscFunctionBegin;
362: if (bs != PETSC_DECIDE) { /* don't mess with an already set block size */
363: PetscCall(MatSetBlockSize(A, bs));
364: }
365: PetscCall(PetscLayoutSetUp(A->rmap));
366: PetscCall(PetscLayoutSetUp(A->cmap));
367: PetscCall(MatGetBlockSizes(A, &bs, &cbs));
368: /* these routines assumes bs == cbs, this should be checked somehow */
369: PetscCall(MatSeqBAIJSetPreallocation(A, bs, 0, dnnz));
370: PetscCall(MatMPIBAIJSetPreallocation(A, bs, 0, dnnz, 0, onnz));
371: PetscCall(MatSeqSBAIJSetPreallocation(A, bs, 0, dnnzu));
372: PetscCall(MatMPISBAIJSetPreallocation(A, bs, 0, dnnzu, 0, onnzu));
373: /*
374: In general, we have to do extra work to preallocate for scalar (AIJ) or unassembled (IS) matrices so we check whether it will do any
375: good before going on with it.
376: */
377: PetscCall(PetscObjectHasFunction((PetscObject)A, "MatMPIAIJSetPreallocation_C", &aij));
378: PetscCall(PetscObjectHasFunction((PetscObject)A, "MatISSetPreallocation_C", &is));
379: PetscCall(PetscObjectHasFunction((PetscObject)A, "MatHYPRESetPreallocation_C", &hyp));
380: if (!aij && !is && !hyp) PetscCall(PetscObjectHasFunction((PetscObject)A, "MatSeqAIJSetPreallocation_C", &aij));
381: if (aij || is || hyp) {
382: if (bs == cbs && bs == 1) {
383: PetscCall(MatSeqAIJSetPreallocation(A, 0, dnnz));
384: PetscCall(MatMPIAIJSetPreallocation(A, 0, dnnz, 0, onnz));
385: PetscCall(MatISSetPreallocation(A, 0, dnnz, 0, onnz));
386: #if defined(PETSC_HAVE_HYPRE)
387: PetscCall(MatHYPRESetPreallocation(A, 0, dnnz, 0, onnz));
388: #endif
389: } else { /* Convert block-row precallocation to scalar-row */
390: PetscInt i, m, *sdnnz, *sonnz;
391: PetscCall(MatGetLocalSize(A, &m, NULL));
392: PetscCall(PetscMalloc2((!!dnnz) * m, &sdnnz, (!!onnz) * m, &sonnz));
393: for (i = 0; i < m; i++) {
394: if (dnnz) sdnnz[i] = dnnz[i / bs] * cbs;
395: if (onnz) sonnz[i] = onnz[i / bs] * cbs;
396: }
397: PetscCall(MatSeqAIJSetPreallocation(A, 0, dnnz ? sdnnz : NULL));
398: PetscCall(MatMPIAIJSetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL));
399: PetscCall(MatISSetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL));
400: #if defined(PETSC_HAVE_HYPRE)
401: PetscCall(MatHYPRESetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL));
402: #endif
403: PetscCall(PetscFree2(sdnnz, sonnz));
404: }
405: }
406: PetscFunctionReturn(PETSC_SUCCESS);
407: }
409: /*@C
410: MatHeaderMerge - Merges some information from the header of `C` to `A`; the `C` object is then destroyed
412: Collective, No Fortran Support
414: Input Parameters:
415: + A - a `Mat` being merged into
416: - C - the `Mat` providing the merge information
418: Level: developer
420: Notes:
421: `A` and `C` must be of the same type.
422: The object list and query function list in `A` are retained, as well as the object name, and prefix.
423: The object state of `A` is increased by 1.
425: Developer Note:
426: This is somewhat different from `MatHeaderReplace()`, it would be nice to merge the code
428: .seealso: `Mat`, `MatHeaderReplace()`
429: @*/
430: PetscErrorCode MatHeaderMerge(Mat A, Mat *C)
431: {
432: PetscInt refct;
433: PetscOps Abops;
434: struct _MatOps Aops;
435: char *mtype, *mname, *mprefix;
436: Mat_Product *product;
437: Mat_Redundant *redundant;
438: PetscObjectState state;
439: PetscObjectList olist;
440: PetscFunctionList qlist;
442: PetscFunctionBegin;
445: if (A == *C) PetscFunctionReturn(PETSC_SUCCESS);
446: PetscCheckSameTypeAndComm(A, 1, *C, 2);
447: /* save the parts of A we need */
448: Abops = ((PetscObject)A)->bops[0];
449: Aops = A->ops[0];
450: refct = ((PetscObject)A)->refct;
451: mtype = ((PetscObject)A)->type_name;
452: mname = ((PetscObject)A)->name;
453: state = ((PetscObject)A)->state;
454: mprefix = ((PetscObject)A)->prefix;
455: product = A->product;
456: redundant = A->redundant;
457: qlist = ((PetscObject)A)->qlist;
458: olist = ((PetscObject)A)->olist;
460: /* zero these so the destroy below does not free them */
461: ((PetscObject)A)->type_name = NULL;
462: ((PetscObject)A)->name = NULL;
463: ((PetscObject)A)->qlist = NULL;
464: ((PetscObject)A)->olist = NULL;
466: /*
467: free all the interior data structures from mat
468: cannot use PetscUseTypeMethod(A,destroy); because compiler
469: thinks it may print NULL type_name and name
470: */
471: PetscTryTypeMethod(A, destroy);
473: PetscCall(PetscFree(A->defaultvectype));
474: PetscCall(PetscFree(A->defaultrandtype));
475: PetscCall(PetscLayoutDestroy(&A->rmap));
476: PetscCall(PetscLayoutDestroy(&A->cmap));
477: PetscCall(PetscComposedQuantitiesDestroy((PetscObject)A));
479: /* copy C over to A */
480: PetscCall(PetscFree(A->factorprefix));
481: PetscCall(PetscMemcpy(A, *C, sizeof(struct _p_Mat)));
483: /* return the parts of A we saved */
484: ((PetscObject)A)->bops[0] = Abops;
485: A->ops[0] = Aops;
486: ((PetscObject)A)->refct = refct;
487: ((PetscObject)A)->type_name = mtype;
488: ((PetscObject)A)->name = mname;
489: ((PetscObject)A)->prefix = mprefix;
490: ((PetscObject)A)->state = state + 1;
491: A->product = product;
492: A->redundant = redundant;
494: /* Append the saved lists */
495: PetscCall(PetscFunctionListDuplicate(qlist, &((PetscObject)A)->qlist));
496: PetscCall(PetscObjectListDuplicate(olist, &((PetscObject)A)->olist));
497: PetscCall(PetscFunctionListDestroy(&qlist));
498: PetscCall(PetscObjectListDestroy(&olist));
500: /* since these two are copied into A we do not want them destroyed in C */
501: ((PetscObject)*C)->qlist = NULL;
502: ((PetscObject)*C)->olist = NULL;
503: PetscCall(PetscHeaderDestroy(C));
504: PetscFunctionReturn(PETSC_SUCCESS);
505: }
507: /*@
508: MatHeaderReplace - Replaces the internal data of matrix `A` by the internal data of matrix `C` while deleting the outer wrapper of `C`
510: Input Parameters:
511: + A - a `Mat` whose internal data is to be replaced
512: - C - the `Mat` providing new internal data for `A`
514: Level: advanced
516: Example Usage\:
517: .vb
518: Mat C;
519: MatCreateSeqAIJWithArrays(..., &C);
520: MatHeaderReplace(A, &C);
521: // C has been destroyed and A contains the matrix entries of C
522: .ve
524: Note:
525: This can be used inside a function provided to `SNESSetJacobian()`, `TSSetRHSJacobian()`, or `TSSetIJacobian()` in cases where the user code
526: computes an entirely new sparse matrix (generally with a different matrix nonzero structure/pattern) for each Newton update.
527: It is usually better to reuse the matrix nonzero structure of `A` instead of constructing an entirely new one.
529: Developer Note:
530: This is somewhat different from `MatHeaderMerge()` it would be nice to merge the code
532: .seealso: `Mat`, `MatHeaderMerge()`
533: @*/
534: PetscErrorCode MatHeaderReplace(Mat A, Mat *C)
535: {
536: PetscInt refct;
537: PetscObjectState state;
538: struct _p_Mat buffer;
539: MatStencilInfo stencil;
541: PetscFunctionBegin;
544: if (A == *C) PetscFunctionReturn(PETSC_SUCCESS);
545: PetscCheckSameComm(A, 1, *C, 2);
546: PetscCheck(((PetscObject)*C)->refct == 1, PetscObjectComm((PetscObject)C), PETSC_ERR_ARG_WRONGSTATE, "Object C has refct %" PetscInt_FMT " > 1, would leave hanging reference", ((PetscObject)*C)->refct);
548: /* swap C and A */
549: refct = ((PetscObject)A)->refct;
550: state = ((PetscObject)A)->state;
551: stencil = A->stencil;
552: PetscCall(PetscMemcpy(&buffer, A, sizeof(struct _p_Mat)));
553: PetscCall(PetscMemcpy(A, *C, sizeof(struct _p_Mat)));
554: PetscCall(PetscMemcpy(*C, &buffer, sizeof(struct _p_Mat)));
555: ((PetscObject)A)->refct = refct;
556: ((PetscObject)A)->state = state + 1;
557: A->stencil = stencil;
559: ((PetscObject)*C)->refct = 1;
560: PetscCall(MatDestroy(C));
561: PetscFunctionReturn(PETSC_SUCCESS);
562: }
564: /*@
565: MatBindToCPU - marks a matrix to temporarily stay on the CPU and perform computations on the CPU
567: Logically Collective
569: Input Parameters:
570: + A - the matrix
571: - flg - bind to the CPU if value of `PETSC_TRUE`
573: Level: intermediate
575: Note:
576: `MATAIJKOKKOS` has yet to implement CPU binding. If Kokkos is configured without GPU support,
577: we deem a `MATAIJKOKKOS` matrix as bound to the CPU. Different from `MATAIJVIENNACL` with a CPU
578: backend, `MATAIJKOKKOS` always use its own operation implementations (in constrast to CPU-bound
579: `MATAIJVIENNACL`, which uses `MATAIJ`'s CPU operations).
581: .seealso: [](ch_matrices), `Mat`, `MatBoundToCPU()`
582: @*/
583: PetscErrorCode MatBindToCPU(Mat A, PetscBool flg)
584: {
585: PetscFunctionBegin;
588: #if defined(PETSC_HAVE_DEVICE)
589: if (A->boundtocpu == flg) PetscFunctionReturn(PETSC_SUCCESS);
590: A->boundtocpu = flg;
591: PetscTryTypeMethod(A, bindtocpu, flg);
592: #endif
593: PetscFunctionReturn(PETSC_SUCCESS);
594: }
596: /*@
597: MatBoundToCPU - query if a matrix is bound to the CPU
599: Input Parameter:
600: . A - the matrix
602: Output Parameter:
603: . flg - the logical flag
605: Level: intermediate
607: .seealso: [](ch_matrices), `Mat`, `MatBindToCPU()`
608: @*/
609: PetscErrorCode MatBoundToCPU(Mat A, PetscBool *flg)
610: {
611: PetscFunctionBegin;
613: PetscAssertPointer(flg, 2);
614: #if defined(PETSC_HAVE_DEVICE)
615: *flg = A->boundtocpu;
616: #else
617: *flg = PETSC_TRUE;
618: #endif
619: PetscFunctionReturn(PETSC_SUCCESS);
620: }
622: PetscErrorCode MatSetValuesCOO_Basic(Mat A, const PetscScalar coo_v[], InsertMode imode)
623: {
624: IS is_coo_i, is_coo_j;
625: const PetscInt *coo_i, *coo_j;
626: PetscInt n, n_i, n_j;
627: PetscScalar zero = 0.;
629: PetscFunctionBegin;
630: PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_coo_i", (PetscObject *)&is_coo_i));
631: PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_coo_j", (PetscObject *)&is_coo_j));
632: PetscCheck(is_coo_i, PetscObjectComm((PetscObject)A), PETSC_ERR_COR, "Missing coo_i IS");
633: PetscCheck(is_coo_j, PetscObjectComm((PetscObject)A), PETSC_ERR_COR, "Missing coo_j IS");
634: PetscCall(ISGetLocalSize(is_coo_i, &n_i));
635: PetscCall(ISGetLocalSize(is_coo_j, &n_j));
636: PetscCheck(n_i == n_j, PETSC_COMM_SELF, PETSC_ERR_COR, "Wrong local size %" PetscInt_FMT " != %" PetscInt_FMT, n_i, n_j);
637: PetscCall(ISGetIndices(is_coo_i, &coo_i));
638: PetscCall(ISGetIndices(is_coo_j, &coo_j));
639: if (imode != ADD_VALUES) PetscCall(MatZeroEntries(A));
640: for (n = 0; n < n_i; n++) PetscCall(MatSetValue(A, coo_i[n], coo_j[n], coo_v ? coo_v[n] : zero, ADD_VALUES));
641: PetscCall(ISRestoreIndices(is_coo_i, &coo_i));
642: PetscCall(ISRestoreIndices(is_coo_j, &coo_j));
643: PetscFunctionReturn(PETSC_SUCCESS);
644: }
646: PetscErrorCode MatSetPreallocationCOO_Basic(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[])
647: {
648: Mat preallocator;
649: IS is_coo_i, is_coo_j;
650: PetscInt ncoo_i;
651: PetscScalar zero = 0.0;
653: PetscFunctionBegin;
654: PetscCall(PetscIntCast(ncoo, &ncoo_i));
655: PetscCall(PetscLayoutSetUp(A->rmap));
656: PetscCall(PetscLayoutSetUp(A->cmap));
657: PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &preallocator));
658: PetscCall(MatSetType(preallocator, MATPREALLOCATOR));
659: PetscCall(MatSetSizes(preallocator, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
660: PetscCall(MatSetLayouts(preallocator, A->rmap, A->cmap));
661: PetscCall(MatSetUp(preallocator));
662: for (PetscCount n = 0; n < ncoo; n++) PetscCall(MatSetValue(preallocator, coo_i[n], coo_j[n], zero, INSERT_VALUES));
663: PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY));
664: PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY));
665: PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_TRUE, A));
666: PetscCall(MatDestroy(&preallocator));
667: PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ncoo_i, coo_i, PETSC_COPY_VALUES, &is_coo_i));
668: PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ncoo_i, coo_j, PETSC_COPY_VALUES, &is_coo_j));
669: PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc_coo_i", (PetscObject)is_coo_i));
670: PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc_coo_j", (PetscObject)is_coo_j));
671: PetscCall(ISDestroy(&is_coo_i));
672: PetscCall(ISDestroy(&is_coo_j));
673: PetscFunctionReturn(PETSC_SUCCESS);
674: }
676: /*@C
677: MatSetPreallocationCOO - set preallocation for matrices using a coordinate format of the entries with global indices
679: Collective
681: Input Parameters:
682: + A - matrix being preallocated
683: . ncoo - number of entries
684: . coo_i - row indices
685: - coo_j - column indices
687: Level: beginner
689: Notes:
690: The indices within `coo_i` and `coo_j` may be modified within this function. The caller should not rely on them
691: having any specific value after this function returns. The arrays can be freed or reused immediately
692: after this function returns.
694: Entries can be repeated, see `MatSetValuesCOO()`. Entries with negative row or column indices are allowed
695: but will be ignored. The corresponding entries in `MatSetValuesCOO()` will be ignored too. Remote entries
696: are allowed and will be properly added or inserted to the matrix, unless the matrix option `MAT_IGNORE_OFF_PROC_ENTRIES`
697: is set, in which case remote entries are ignored, or `MAT_NO_OFF_PROC_ENTRIES` is set, in which case an error will be generated.
699: If you just want to create a sequential AIJ matrix (`MATSEQAIJ`), and your matrix entries in COO format are unique, you can also use
700: `MatCreateSeqAIJFromTriple()`. But that is not recommended for iterative applications.
702: .seealso: [](ch_matrices), `Mat`, `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`,
703: `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOOLocal()`,
704: `DMSetMatrixPreallocateSkip()`, `MatCreateSeqAIJFromTriple()`
705: @*/
706: PetscErrorCode MatSetPreallocationCOO(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[])
707: {
708: PetscErrorCode (*f)(Mat, PetscCount, PetscInt[], PetscInt[]) = NULL;
710: PetscFunctionBegin;
713: if (ncoo) PetscAssertPointer(coo_i, 3);
714: if (ncoo) PetscAssertPointer(coo_j, 4);
715: PetscCall(PetscLayoutSetUp(A->rmap));
716: PetscCall(PetscLayoutSetUp(A->cmap));
717: PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetPreallocationCOO_C", &f));
719: PetscCall(PetscLogEventBegin(MAT_PreallCOO, A, 0, 0, 0));
720: if (f) PetscCall((*f)(A, ncoo, coo_i, coo_j));
721: else PetscCall(MatSetPreallocationCOO_Basic(A, ncoo, coo_i, coo_j)); /* allow fallback, very slow */
722: PetscCall(PetscLogEventEnd(MAT_PreallCOO, A, 0, 0, 0));
723: A->preallocated = PETSC_TRUE;
724: A->nonzerostate++;
725: PetscFunctionReturn(PETSC_SUCCESS);
726: }
728: /*@C
729: MatSetPreallocationCOOLocal - set preallocation for matrices using a coordinate format of the entries with local indices
731: Collective
733: Input Parameters:
734: + A - matrix being preallocated
735: . ncoo - number of entries
736: . coo_i - row indices (local numbering; may be modified)
737: - coo_j - column indices (local numbering; may be modified)
739: Level: beginner
741: Notes:
742: The local indices are translated using the local to global mapping, thus `MatSetLocalToGlobalMapping()` must have been
743: called prior to this function. For matrices created with `DMCreateMatrix()` the local to global mapping is often already provided.
745: The indices `coo_i` and `coo_j` may be modified within this function. They might be translated to corresponding global
746: indices, but the caller should not rely on them having any specific value after this function returns. The arrays
747: can be freed or reused immediately after this function returns.
749: Entries can be repeated, see `MatSetValuesCOO()`. Entries with negative row or column indices are allowed
750: but will be ignored. The corresponding entries in `MatSetValuesCOO()` will be ignored too. Remote entries
751: are allowed and will be properly added or inserted to the matrix.
753: .seealso: [](ch_matrices), `Mat`, `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`,
754: `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOO()`,
755: `DMSetMatrixPreallocateSkip()`
756: @*/
757: PetscErrorCode MatSetPreallocationCOOLocal(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[])
758: {
759: PetscErrorCode (*f)(Mat, PetscCount, PetscInt[], PetscInt[]) = NULL;
761: PetscFunctionBegin;
764: if (ncoo) PetscAssertPointer(coo_i, 3);
765: if (ncoo) PetscAssertPointer(coo_j, 4);
766: PetscCall(PetscLayoutSetUp(A->rmap));
767: PetscCall(PetscLayoutSetUp(A->cmap));
769: PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetPreallocationCOOLocal_C", &f));
770: if (f) {
771: PetscCall((*f)(A, ncoo, coo_i, coo_j));
772: A->nonzerostate++;
773: } else {
774: PetscInt ncoo_i;
775: ISLocalToGlobalMapping ltog_row, ltog_col;
777: PetscCall(MatGetLocalToGlobalMapping(A, <og_row, <og_col));
778: if (ltog_row) {
779: PetscCall(PetscIntCast(ncoo, &ncoo_i));
780: PetscCall(ISLocalToGlobalMappingApply(ltog_row, ncoo_i, coo_i, coo_i));
781: }
782: if (ltog_col) {
783: PetscCall(PetscIntCast(ncoo, &ncoo_i));
784: PetscCall(ISLocalToGlobalMappingApply(ltog_col, ncoo_i, coo_j, coo_j));
785: }
786: PetscCall(MatSetPreallocationCOO(A, ncoo, coo_i, coo_j));
787: }
788: A->preallocated = PETSC_TRUE;
789: PetscFunctionReturn(PETSC_SUCCESS);
790: }
792: /*@
793: MatSetValuesCOO - set values at once in a matrix preallocated using `MatSetPreallocationCOO()`
795: Collective
797: Input Parameters:
798: + A - matrix being preallocated
799: . coo_v - the matrix values (can be `NULL`)
800: - imode - the insert mode
802: Level: beginner
804: Notes:
805: The values must follow the order of the indices prescribed with `MatSetPreallocationCOO()` or `MatSetPreallocationCOOLocal()`.
807: When repeated entries are specified in the COO indices the `coo_v` values are first properly summed, regardless of the value of imode.
808: The imode flag indicates if coo_v must be added to the current values of the matrix (`ADD_VALUES`) or overwritten (`INSERT_VALUES`).
810: `MatAssemblyBegin()` and `MatAssemblyEnd()` do not need to be called after this routine. It automatically handles the assembly process.
812: .seealso: [](ch_matrices), `Mat`, `MatSetPreallocationCOO()`, `MatSetPreallocationCOOLocal()`, `InsertMode`, `INSERT_VALUES`, `ADD_VALUES`
813: @*/
814: PetscErrorCode MatSetValuesCOO(Mat A, const PetscScalar coo_v[], InsertMode imode)
815: {
816: PetscErrorCode (*f)(Mat, const PetscScalar[], InsertMode) = NULL;
817: PetscBool oldFlg;
819: PetscFunctionBegin;
822: MatCheckPreallocated(A, 1);
824: PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetValuesCOO_C", &f));
825: PetscCall(PetscLogEventBegin(MAT_SetVCOO, A, 0, 0, 0));
826: if (f) {
827: PetscCall((*f)(A, coo_v, imode)); // all known COO implementations do not use MatStash. They do their own off-proc communication
828: PetscCall(MatGetOption(A, MAT_NO_OFF_PROC_ENTRIES, &oldFlg));
829: PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); // set A->nooffprocentries to avoid costly MatStash scatter in MatAssembly
830: } else {
831: PetscCall(MatSetValuesCOO_Basic(A, coo_v, imode)); // fall back to MatSetValues, which might use MatStash
832: }
833: PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
834: PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
835: if (f) PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, oldFlg));
836: PetscCall(PetscLogEventEnd(MAT_SetVCOO, A, 0, 0, 0));
837: PetscFunctionReturn(PETSC_SUCCESS);
838: }
840: /*@
841: MatSetBindingPropagates - Sets whether the state of being bound to the CPU for a GPU matrix type propagates to child and some other associated objects
843: Input Parameters:
844: + A - the matrix
845: - flg - flag indicating whether the boundtocpu flag should be propagated
847: Level: developer
849: Notes:
850: If the value of flg is set to true, the following will occur
851: + `MatCreateSubMatrices()` and `MatCreateRedundantMatrix()` - bind created matrices to CPU if the input matrix is bound to the CPU.
852: - `MatCreateVecs()` - bind created vectors to CPU if the input matrix is bound to the CPU.
854: The bindingpropagates flag itself is also propagated by the above routines.
856: Developer Notes:
857: If the fine-scale `DMDA` has the `-dm_bind_below` option set to true, then `DMCreateInterpolationScale()` calls `MatSetBindingPropagates()`
858: on the restriction/interpolation operator to set the bindingpropagates flag to true.
860: .seealso: [](ch_matrices), `Mat`, `VecSetBindingPropagates()`, `MatGetBindingPropagates()`
861: @*/
862: PetscErrorCode MatSetBindingPropagates(Mat A, PetscBool flg)
863: {
864: PetscFunctionBegin;
866: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
867: A->bindingpropagates = flg;
868: #endif
869: PetscFunctionReturn(PETSC_SUCCESS);
870: }
872: /*@
873: MatGetBindingPropagates - Gets whether the state of being bound to the CPU for a GPU matrix type propagates to child and some other associated objects
875: Input Parameter:
876: . A - the matrix
878: Output Parameter:
879: . flg - flag indicating whether the boundtocpu flag will be propagated
881: Level: developer
883: .seealso: [](ch_matrices), `Mat`, `MatSetBindingPropagates()`
884: @*/
885: PetscErrorCode MatGetBindingPropagates(Mat A, PetscBool *flg)
886: {
887: PetscFunctionBegin;
889: PetscAssertPointer(flg, 2);
890: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
891: *flg = A->bindingpropagates;
892: #else
893: *flg = PETSC_FALSE;
894: #endif
895: PetscFunctionReturn(PETSC_SUCCESS);
896: }