xref: /petsc/src/mat/utils/gcreate.c (revision dfe00d7b09a80b04879b1d4346c83210468b990e)
1af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h"  I*/
27807a1faSBarry Smith 
3cc1836a6SJunchao Zhang #include <../src/mat/impls/aij/seq/aij.h>
4cc1836a6SJunchao Zhang #include <../src/mat/impls/aij/mpi/mpiaij.h>
5cc1836a6SJunchao Zhang 
MatSetBlockSizes_Default(Mat mat,PetscInt rbs,PetscInt cbs)69de2952eSStefano Zampini PetscErrorCode MatSetBlockSizes_Default(Mat mat, PetscInt rbs, PetscInt cbs)
7d71ae5a4SJacob Faibussowitsch {
846533700Sstefano_zampini   PetscFunctionBegin;
93ba16761SJacob Faibussowitsch   if (!mat->preallocated) PetscFunctionReturn(PETSC_SUCCESS);
10aed4548fSBarry Smith   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);
11aed4548fSBarry Smith   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);
123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1346533700Sstefano_zampini }
1446533700Sstefano_zampini 
MatShift_Basic(Mat Y,PetscScalar a)159de2952eSStefano Zampini PetscErrorCode MatShift_Basic(Mat Y, PetscScalar a)
16d71ae5a4SJacob Faibussowitsch {
17cc1836a6SJunchao Zhang   PetscInt    i, start, end, oldValA = 0, oldValB = 0;
187d68702bSBarry Smith   PetscScalar alpha = a;
197d68702bSBarry Smith   PetscBool   prevoption;
20cc1836a6SJunchao Zhang   PetscBool   isSeqAIJDerived, isMPIAIJDerived; // all classes sharing SEQAIJHEADER or MPIAIJHEADER
21cc1836a6SJunchao Zhang   Mat         A = NULL, B = NULL;
227d68702bSBarry Smith 
237d68702bSBarry Smith   PetscFunctionBegin;
249566063dSJacob Faibussowitsch   PetscCall(MatGetOption(Y, MAT_NO_OFF_PROC_ENTRIES, &prevoption));
259566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Y, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
26cc1836a6SJunchao Zhang   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &isSeqAIJDerived, MATSEQAIJ, MATSEQBAIJ, MATSEQSBAIJ, ""));
27cc1836a6SJunchao Zhang   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &isMPIAIJDerived, MATMPIAIJ, MATMPIBAIJ, MATMPISBAIJ, ""));
28cc1836a6SJunchao Zhang 
29cc1836a6SJunchao Zhang   if (isSeqAIJDerived) A = Y;
30cc1836a6SJunchao Zhang   else if (isMPIAIJDerived) {
31cc1836a6SJunchao Zhang     Mat_MPIAIJ *mpiaij = (Mat_MPIAIJ *)Y->data;
32cc1836a6SJunchao Zhang     A                  = mpiaij->A;
33cc1836a6SJunchao Zhang     B                  = mpiaij->B;
34cc1836a6SJunchao Zhang   }
35cc1836a6SJunchao Zhang 
36cc1836a6SJunchao Zhang   if (A) {
37f4f49eeaSPierre Jolivet     oldValA                        = ((Mat_SeqAIJ *)A->data)->nonew;
38f4f49eeaSPierre Jolivet     ((Mat_SeqAIJ *)A->data)->nonew = 0; // so that new nonzero locations are allowed
39cc1836a6SJunchao Zhang   }
40cc1836a6SJunchao Zhang   if (B) {
41f4f49eeaSPierre Jolivet     oldValB                        = ((Mat_SeqAIJ *)B->data)->nonew;
42f4f49eeaSPierre Jolivet     ((Mat_SeqAIJ *)B->data)->nonew = 0;
43cc1836a6SJunchao Zhang   }
44cc1836a6SJunchao Zhang 
459566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Y, &start, &end));
467d68702bSBarry Smith   for (i = start; i < end; i++) {
4748a46eb9SPierre Jolivet     if (i < Y->cmap->N) PetscCall(MatSetValues(Y, 1, &i, 1, &i, &alpha, ADD_VALUES));
48ab6153dcSStefano Zampini   }
499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Y, MAT_FINAL_ASSEMBLY));
509566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Y, MAT_FINAL_ASSEMBLY));
519566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Y, MAT_NO_OFF_PROC_ENTRIES, prevoption));
52f4f49eeaSPierre Jolivet   if (A) ((Mat_SeqAIJ *)A->data)->nonew = oldValA;
53f4f49eeaSPierre Jolivet   if (B) ((Mat_SeqAIJ *)B->data)->nonew = oldValB;
543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
557d68702bSBarry Smith }
567d68702bSBarry Smith 
5705869f15SSatish Balay /*@
5869dd0797SLois Curfman McInnes   MatCreate - Creates a matrix where the type is determined
5911a5261eSBarry Smith   from either a call to `MatSetType()` or from the options database
602920cce0SJacob Faibussowitsch   with a call to `MatSetFromOptions()`.
6183e1b59cSLois Curfman McInnes 
62d083f849SBarry Smith   Collective
63cb13003dSBarry Smith 
64f69a0ea3SMatthew Knepley   Input Parameter:
65f69a0ea3SMatthew Knepley . comm - MPI communicator
667807a1faSBarry Smith 
677807a1faSBarry Smith   Output Parameter:
68dc401e71SLois Curfman McInnes . A - the matrix
69e0b365e2SLois Curfman McInnes 
70273d9f13SBarry Smith   Options Database Keys:
7111a5261eSBarry Smith + -mat_type seqaij   - `MATSEQAIJ` type, uses `MatCreateSeqAIJ()`
7211a5261eSBarry Smith . -mat_type mpiaij   - `MATMPIAIJ` type, uses `MatCreateAIJ()`
7311a5261eSBarry Smith . -mat_type seqdense - `MATSEQDENSE`, uses `MatCreateSeqDense()`
7411a5261eSBarry Smith . -mat_type mpidense - `MATMPIDENSE` type, uses `MatCreateDense()`
7511a5261eSBarry Smith . -mat_type seqbaij  - `MATSEQBAIJ` type, uses `MatCreateSeqBAIJ()`
7611a5261eSBarry Smith - -mat_type mpibaij  - `MATMPIBAIJ` type, uses `MatCreateBAIJ()`
77e0b365e2SLois Curfman McInnes 
782ef1f0ffSBarry Smith    See the manpages for particular formats (e.g., `MATSEQAIJ`)
7983e1b59cSLois Curfman McInnes    for additional format-specific options.
80e0b365e2SLois Curfman McInnes 
81273d9f13SBarry Smith   Level: beginner
82273d9f13SBarry Smith 
832920cce0SJacob Faibussowitsch   Notes:
842920cce0SJacob Faibussowitsch   The default matrix type is `MATAIJ`, using the routines `MatCreateSeqAIJ()` or
852920cce0SJacob Faibussowitsch   `MatCreateAIJ()` if you do not set a type in the options database. If you never call
862920cce0SJacob Faibussowitsch   `MatSetType()` or `MatSetFromOptions()` it will generate an error when you try to use the
872920cce0SJacob Faibussowitsch   matrix.
882920cce0SJacob Faibussowitsch 
891cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateSeqAIJ()`, `MatCreateAIJ()`,
90db781477SPatrick Sanan           `MatCreateSeqDense()`, `MatCreateDense()`,
91db781477SPatrick Sanan           `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`,
92db781477SPatrick Sanan           `MatCreateSeqSBAIJ()`, `MatCreateSBAIJ()`,
93db781477SPatrick Sanan           `MatConvert()`
94273d9f13SBarry Smith @*/
MatCreate(MPI_Comm comm,Mat * A)95d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreate(MPI_Comm comm, Mat *A)
96d71ae5a4SJacob Faibussowitsch {
97273d9f13SBarry Smith   Mat B;
98273d9f13SBarry Smith 
99273d9f13SBarry Smith   PetscFunctionBegin;
1004f572ea9SToby Isaac   PetscAssertPointer(A, 2);
1019566063dSJacob Faibussowitsch   PetscCall(MatInitializePackage());
1028ba1e511SMatthew Knepley 
1039566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(B, MAT_CLASSID, "Mat", "Matrix", "Mat", comm, MatDestroy, MatView));
1049566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &B->rmap));
1059566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &B->cmap));
1069566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(VECSTANDARD, &B->defaultvectype));
1073faff063SStefano Zampini   PetscCall(PetscStrallocpy(PETSCRANDER48, &B->defaultrandtype));
10826fbe8dcSKarl Rupp 
109b94d7dedSBarry Smith   B->symmetric                   = PETSC_BOOL3_UNKNOWN;
110b94d7dedSBarry Smith   B->hermitian                   = PETSC_BOOL3_UNKNOWN;
111b94d7dedSBarry Smith   B->structurally_symmetric      = PETSC_BOOL3_UNKNOWN;
112b94d7dedSBarry Smith   B->spd                         = PETSC_BOOL3_UNKNOWN;
113b94d7dedSBarry Smith   B->symmetry_eternal            = PETSC_FALSE;
114b94d7dedSBarry Smith   B->structural_symmetry_eternal = PETSC_FALSE;
115b94d7dedSBarry Smith 
11694342113SStefano Zampini   B->congruentlayouts = PETSC_DECIDE;
117273d9f13SBarry Smith   B->preallocated     = PETSC_FALSE;
1186f3d89d0SStefano Zampini #if defined(PETSC_HAVE_DEVICE)
1196f3d89d0SStefano Zampini   B->boundtocpu = PETSC_TRUE;
1206f3d89d0SStefano Zampini #endif
121273d9f13SBarry Smith   *A = B;
1223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
123273d9f13SBarry Smith }
124273d9f13SBarry Smith 
125cc4c1da9SBarry Smith /*@
12677433607SBarry Smith   MatCreateFromOptions - Creates a matrix whose type is set from the options database
12777433607SBarry Smith 
12877433607SBarry Smith   Collective
12977433607SBarry Smith 
13077433607SBarry Smith   Input Parameters:
13177433607SBarry Smith + comm   - MPI communicator
13277433607SBarry Smith . prefix - [optional] prefix for the options database
13377433607SBarry Smith . bs     - the blocksize (commonly 1)
13477433607SBarry Smith . m      - the local number of rows (or `PETSC_DECIDE`)
13577433607SBarry Smith . n      - the local number of columns (or `PETSC_DECIDE` or `PETSC_DETERMINE`)
13677433607SBarry Smith . M      - the global number of rows (or `PETSC_DETERMINE`)
13777433607SBarry Smith - N      - the global number of columns (or `PETSC_DETERMINE`)
13877433607SBarry Smith 
13977433607SBarry Smith   Output Parameter:
14077433607SBarry Smith . A - the matrix
14177433607SBarry Smith 
14277433607SBarry Smith   Options Database Key:
143362f43d5SPierre Jolivet . -mat_type - see `MatType`, for example `aij`, `aijcusparse`, `baij`, `sbaij`, `dense`, defaults to `aij`
14477433607SBarry Smith 
14577433607SBarry Smith   Level: beginner
14677433607SBarry Smith 
14777433607SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateSeqAIJ()`, `MatCreateAIJ()`,
14877433607SBarry Smith           `MatCreateSeqDense()`, `MatCreateDense()`,
14977433607SBarry Smith           `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`,
15077433607SBarry Smith           `MatCreateSeqSBAIJ()`, `MatCreateSBAIJ()`,
15177433607SBarry Smith           `MatConvert()`, `MatCreate()`
15277433607SBarry Smith @*/
MatCreateFromOptions(MPI_Comm comm,const char * prefix,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat * A)15377433607SBarry Smith PetscErrorCode MatCreateFromOptions(MPI_Comm comm, const char *prefix, PetscInt bs, PetscInt m, PetscInt n, PetscInt M, PetscInt N, Mat *A)
15477433607SBarry Smith {
15577433607SBarry Smith   PetscFunctionBegin;
15677433607SBarry Smith   PetscAssertPointer(A, 8);
15777433607SBarry Smith   PetscCall(MatCreate(comm, A));
15877433607SBarry Smith   if (prefix) PetscCall(MatSetOptionsPrefix(*A, prefix));
15977433607SBarry Smith   PetscCall(MatSetBlockSize(*A, bs));
16077433607SBarry Smith   PetscCall(MatSetSizes(*A, m, n, M, N));
16177433607SBarry Smith   PetscCall(MatSetFromOptions(*A));
16277433607SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
16377433607SBarry Smith }
16477433607SBarry Smith 
165422a814eSBarry Smith /*@
16611a5261eSBarry Smith   MatSetErrorIfFailure - Causes `Mat` to generate an immediate error, for example a zero pivot, is detected.
167422a814eSBarry Smith 
168c3339decSBarry Smith   Logically Collective
169422a814eSBarry Smith 
170422a814eSBarry Smith   Input Parameters:
17111a5261eSBarry Smith + mat - matrix obtained from `MatCreate()`
17211a5261eSBarry Smith - flg - `PETSC_TRUE` indicates you want the error generated
173422a814eSBarry Smith 
174422a814eSBarry Smith   Level: advanced
175422a814eSBarry Smith 
17611a5261eSBarry Smith   Note:
17711a5261eSBarry Smith   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
17811a5261eSBarry Smith   or result in a `KSPConvergedReason` indicating the method did not converge.
17911a5261eSBarry Smith 
1801cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `PCSetErrorIfFailure()`, `KSPConvergedReason`, `SNESConvergedReason`
181422a814eSBarry Smith @*/
MatSetErrorIfFailure(Mat mat,PetscBool flg)182d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetErrorIfFailure(Mat mat, PetscBool flg)
183d71ae5a4SJacob Faibussowitsch {
184422a814eSBarry Smith   PetscFunctionBegin;
185422a814eSBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
186422a814eSBarry Smith   PetscValidLogicalCollectiveBool(mat, flg, 2);
18784d44b13SHong Zhang   mat->erroriffailure = flg;
1883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
189422a814eSBarry Smith }
190422a814eSBarry Smith 
191f69a0ea3SMatthew Knepley /*@
192f69a0ea3SMatthew Knepley   MatSetSizes - Sets the local and global sizes, and checks to determine compatibility
193f69a0ea3SMatthew Knepley 
194c3339decSBarry Smith   Collective
195f69a0ea3SMatthew Knepley 
196f69a0ea3SMatthew Knepley   Input Parameters:
197f69a0ea3SMatthew Knepley + A - the matrix
19811a5261eSBarry Smith . m - number of local rows (or `PETSC_DECIDE`)
19911a5261eSBarry Smith . n - number of local columns (or `PETSC_DECIDE`)
20011a5261eSBarry Smith . M - number of global rows (or `PETSC_DETERMINE`)
20111a5261eSBarry Smith - N - number of global columns (or `PETSC_DETERMINE`)
202f69a0ea3SMatthew Knepley 
2032fe279fdSBarry Smith   Level: beginner
2042fe279fdSBarry Smith 
205f69a0ea3SMatthew Knepley   Notes:
2062ef1f0ffSBarry Smith   `m` (`n`) and `M` (`N`) cannot be both `PETSC_DECIDE`
2072ef1f0ffSBarry Smith   If one processor calls this with `M` (`N`) of `PETSC_DECIDE` then all processors must, otherwise the program will hang.
208f69a0ea3SMatthew Knepley 
20911a5261eSBarry Smith   If `PETSC_DECIDE` is not used for the arguments 'm' and 'n', then the
210f69a0ea3SMatthew Knepley   user must ensure that they are chosen to be compatible with the
211f69a0ea3SMatthew Knepley   vectors. To do this, one first considers the matrix-vector product
2122ef1f0ffSBarry Smith   'y = A x'. The `m` that is used in the above routine must match the
213d8bfcecaSStefano Zampini   local size of 'y'. Likewise, the `n` used must match the local size of 'x'.
214f69a0ea3SMatthew Knepley 
215727bdf9bSBarry Smith   If `m` and `n` are not `PETSC_DECIDE`, then the values determine the `PetscLayout` of the matrix and the ranges returned by
216727bdf9bSBarry Smith   `MatGetOwnershipRange()`,  `MatGetOwnershipRanges()`, `MatGetOwnershipRangeColumn()`, and `MatGetOwnershipRangesColumn()`.
217727bdf9bSBarry Smith 
218f73d5cc4SBarry Smith   You cannot change the sizes once they have been set.
219f73d5cc4SBarry Smith 
22011a5261eSBarry Smith   The sizes must be set before `MatSetUp()` or MatXXXSetPreallocation() is called.
221f73d5cc4SBarry Smith 
222727bdf9bSBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetSize()`, `PetscSplitOwnership()`, `MatGetOwnershipRange()`, `MatGetOwnershipRanges()`,
223727bdf9bSBarry Smith           `MatGetOwnershipRangeColumn()`, `MatGetOwnershipRangesColumn()`, `PetscLayout`, `VecSetSizes()`
224f69a0ea3SMatthew Knepley @*/
MatSetSizes(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N)225d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N)
226d71ae5a4SJacob Faibussowitsch {
227f69a0ea3SMatthew Knepley   PetscFunctionBegin;
2280700a824SBarry Smith   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
229a69c7061SStefano Zampini   PetscValidLogicalCollectiveInt(A, M, 4);
230a69c7061SStefano Zampini   PetscValidLogicalCollectiveInt(A, N, 5);
231aed4548fSBarry Smith   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);
232aed4548fSBarry Smith   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);
2339371c9d4SSatish Balay   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,
2349371c9d4SSatish Balay              A->rmap->n, A->rmap->N);
2359371c9d4SSatish Balay   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,
2369371c9d4SSatish Balay              A->cmap->n, A->cmap->N);
237d0f46423SBarry Smith   A->rmap->n = m;
238d0f46423SBarry Smith   A->cmap->n = n;
23959cb773eSBarry Smith   A->rmap->N = M > -1 ? M : A->rmap->N;
24059cb773eSBarry Smith   A->cmap->N = N > -1 ? N : A->cmap->N;
2413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
242f69a0ea3SMatthew Knepley }
243f69a0ea3SMatthew Knepley 
24405869f15SSatish Balay /*@
245273d9f13SBarry Smith   MatSetFromOptions - Creates a matrix where the type is determined
2462920cce0SJacob Faibussowitsch   from the options database.
247273d9f13SBarry Smith 
248c3339decSBarry Smith   Collective
249273d9f13SBarry Smith 
250273d9f13SBarry Smith   Input Parameter:
251fe59aa6dSJacob Faibussowitsch . B - the matrix
252273d9f13SBarry Smith 
253273d9f13SBarry Smith   Options Database Keys:
25411a5261eSBarry Smith + -mat_type seqaij   - `MATSEQAIJ` type, uses `MatCreateSeqAIJ()`
25511a5261eSBarry Smith . -mat_type mpiaij   - `MATMPIAIJ` type, uses `MatCreateAIJ()`
25611a5261eSBarry Smith . -mat_type seqdense - `MATSEQDENSE` type, uses `MatCreateSeqDense()`
25711a5261eSBarry Smith . -mat_type mpidense - `MATMPIDENSE`, uses `MatCreateDense()`
25811a5261eSBarry Smith . -mat_type seqbaij  - `MATSEQBAIJ`, uses `MatCreateSeqBAIJ()`
25911a5261eSBarry Smith - -mat_type mpibaij  - `MATMPIBAIJ`, uses `MatCreateBAIJ()`
260273d9f13SBarry Smith 
2612ef1f0ffSBarry Smith    See the manpages for particular formats (e.g., `MATSEQAIJ`)
262273d9f13SBarry Smith    for additional format-specific options.
263bd9ce289SLois Curfman McInnes 
2641d69843bSLois Curfman McInnes   Level: beginner
2651d69843bSLois Curfman McInnes 
2662920cce0SJacob Faibussowitsch   Notes:
2672920cce0SJacob Faibussowitsch   Generates a parallel MPI matrix if the communicator has more than one processor.  The default
2682920cce0SJacob Faibussowitsch   matrix type is `MATAIJ`, using the routines `MatCreateSeqAIJ()` and `MatCreateAIJ()` if you
2692920cce0SJacob Faibussowitsch   do not select a type in the options database.
2702920cce0SJacob Faibussowitsch 
2711cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateSeqAIJ()`, `MatCreateAIJ()`,
272db781477SPatrick Sanan           `MatCreateSeqDense()`, `MatCreateDense()`,
273db781477SPatrick Sanan           `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`,
274db781477SPatrick Sanan           `MatCreateSeqSBAIJ()`, `MatCreateSBAIJ()`,
275db781477SPatrick Sanan           `MatConvert()`
2767807a1faSBarry Smith @*/
MatSetFromOptions(Mat B)277d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetFromOptions(Mat B)
278d71ae5a4SJacob Faibussowitsch {
279f3be49caSLisandro Dalcin   const char *deft = MATAIJ;
280f3be49caSLisandro Dalcin   char        type[256];
28169df5c0cSJed Brown   PetscBool   flg, set;
28258b7e2c1SStefano Zampini   PetscInt    bind_below = 0, newbs = -1;
283dbb450caSBarry Smith 
2843a40ed3dSBarry Smith   PetscFunctionBegin;
2850700a824SBarry Smith   PetscValidHeaderSpecific(B, MAT_CLASSID, 1);
286f3be49caSLisandro Dalcin 
287d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)B);
288535b19f3SBarry Smith 
2899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_block_size", "Set the blocksize used to store the matrix", "MatSetBlockSize", newbs, &newbs, &flg));
290535b19f3SBarry Smith   if (flg) {
2919566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetBlockSize(B->rmap, newbs));
2929566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetBlockSize(B->cmap, newbs));
293535b19f3SBarry Smith   }
294535b19f3SBarry Smith 
295*ff19ca73SPierre Jolivet   PetscCall(PetscOptionsFList("-mat_type", "Matrix type", "MatSetType", MatList, deft, type, PETSC_STATIC_ARRAY_LENGTH(type), &flg));
296*ff19ca73SPierre Jolivet   if (flg) PetscCall(MatSetType(B, type));
297*ff19ca73SPierre Jolivet   else if (!((PetscObject)B)->type_name) PetscCall(MatSetType(B, deft));
298*ff19ca73SPierre Jolivet 
299*ff19ca73SPierre Jolivet   if (newbs > 0) PetscTryTypeMethod(B, setblocksizes, newbs, newbs);
300f3be49caSLisandro Dalcin 
3019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsName("-mat_is_symmetric", "Checks if mat is symmetric on MatAssemblyEnd()", "MatIsSymmetric", &B->checksymmetryonassembly));
3029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_is_symmetric", "Checks if mat is symmetric on MatAssemblyEnd()", "MatIsSymmetric", B->checksymmetrytol, &B->checksymmetrytol, NULL));
3039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_null_space_test", "Checks if provided null space is correct in MatAssemblyEnd()", "MatSetNullSpaceTest", B->checknullspaceonassembly, &B->checknullspaceonassembly, NULL));
3049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_error_if_failure", "Generate an error if an error occurs when factoring the matrix", "MatSetErrorIfFailure", B->erroriffailure, &B->erroriffailure, NULL));
305840d65ccSBarry Smith 
306dbbe0bcdSBarry Smith   PetscTryTypeMethod(B, setfromoptions, PetscOptionsObject);
307f3be49caSLisandro Dalcin 
30869df5c0cSJed Brown   flg = PETSC_FALSE;
3090b4b7b1cSBarry Smith   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));
3109566063dSJacob Faibussowitsch   if (set) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, flg));
31169df5c0cSJed Brown   flg = PETSC_FALSE;
3120b4b7b1cSBarry Smith   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));
3139566063dSJacob Faibussowitsch   if (set) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, flg));
314478db826SMatthew G. Knepley   flg = PETSC_FALSE;
3159566063dSJacob Faibussowitsch   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));
3169566063dSJacob Faibussowitsch   if (set) PetscCall(MatSetOption(B, MAT_IGNORE_ZERO_ENTRIES, flg));
31769df5c0cSJed Brown 
3181a2c6b5cSJunchao Zhang   flg = PETSC_FALSE;
3199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_form_explicit_transpose", "Hint to form an explicit transpose for operations like MatMultTranspose", "MatSetOption", flg, &flg, &set));
3209566063dSJacob Faibussowitsch   if (set) PetscCall(MatSetOption(B, MAT_FORM_EXPLICIT_TRANSPOSE, flg));
3211a2c6b5cSJunchao Zhang 
32216e04d98SRichard Tran Mills   /* Bind to CPU if below a user-specified size threshold.
32316e04d98SRichard Tran Mills    * This perhaps belongs in the options for the GPU Mat types, but MatBindToCPU() does nothing when called on non-GPU types,
32416e04d98SRichard Tran Mills    * and putting it here makes is more maintainable than duplicating this for all. */
3259566063dSJacob Faibussowitsch   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));
32648a46eb9SPierre Jolivet   if (flg && B->rmap->n < bind_below) PetscCall(MatBindToCPU(B, PETSC_TRUE));
32716e04d98SRichard Tran Mills 
3285d973c19SBarry Smith   /* process any options handlers added with PetscObjectAddOptionsHandler() */
329dbbe0bcdSBarry Smith   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)B, PetscOptionsObject));
330d0609cedSBarry Smith   PetscOptionsEnd();
3313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3327807a1faSBarry Smith }
3337807a1faSBarry Smith 
3345d83a8b1SBarry Smith /*@
33511a5261eSBarry Smith   MatXAIJSetPreallocation - set preallocation for serial and parallel `MATAIJ`, `MATBAIJ`, and `MATSBAIJ` matrices and their unassembled versions.
33663562e91SJed Brown 
337c3339decSBarry Smith   Collective
33863562e91SJed Brown 
3394165533cSJose E. Roman   Input Parameters:
34063562e91SJed Brown + A     - matrix being preallocated
34163562e91SJed Brown . bs    - block size
34241319c1dSStefano Zampini . dnnz  - number of nonzero column blocks per block row of diagonal part of parallel matrix
34341319c1dSStefano Zampini . onnz  - number of nonzero column blocks per block row of off-diagonal part of parallel matrix
34441319c1dSStefano Zampini . dnnzu - number of nonzero column blocks per block row of upper-triangular part of diagonal part of parallel matrix
34541319c1dSStefano Zampini - onnzu - number of nonzero column blocks per block row of upper-triangular part of off-diagonal part of parallel matrix
34663562e91SJed Brown 
34763562e91SJed Brown   Level: beginner
34863562e91SJed Brown 
3491cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, `MatMPIBAIJSetPreallocation()`,
3502fe279fdSBarry Smith           `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`,
351db781477SPatrick Sanan           `PetscSplitOwnership()`
35263562e91SJed Brown @*/
MatXAIJSetPreallocation(Mat A,PetscInt bs,const PetscInt dnnz[],const PetscInt onnz[],const PetscInt dnnzu[],const PetscInt onnzu[])353d71ae5a4SJacob Faibussowitsch PetscErrorCode MatXAIJSetPreallocation(Mat A, PetscInt bs, const PetscInt dnnz[], const PetscInt onnz[], const PetscInt dnnzu[], const PetscInt onnzu[])
354d71ae5a4SJacob Faibussowitsch {
35541319c1dSStefano Zampini   PetscInt  cbs;
3560cd8b6e2SStefano Zampini   PetscBool aij, is, hyp;
35763562e91SJed Brown 
35863562e91SJed Brown   PetscFunctionBegin;
35941319c1dSStefano Zampini   if (bs != PETSC_DECIDE) { /* don't mess with an already set block size */
3609566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(A, bs));
36141319c1dSStefano Zampini   }
3629566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->rmap));
3639566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->cmap));
3649566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSizes(A, &bs, &cbs));
36541319c1dSStefano Zampini   /* these routines assumes bs == cbs, this should be checked somehow */
3669566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(A, bs, 0, dnnz));
3679566063dSJacob Faibussowitsch   PetscCall(MatMPIBAIJSetPreallocation(A, bs, 0, dnnz, 0, onnz));
3689566063dSJacob Faibussowitsch   PetscCall(MatSeqSBAIJSetPreallocation(A, bs, 0, dnnzu));
3699566063dSJacob Faibussowitsch   PetscCall(MatMPISBAIJSetPreallocation(A, bs, 0, dnnzu, 0, onnzu));
37063562e91SJed Brown   /*
371e8bd9bafSStefano Zampini     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
37263562e91SJed Brown     good before going on with it.
37363562e91SJed Brown   */
3740cd8b6e2SStefano Zampini   PetscCall(PetscObjectHasFunction((PetscObject)A, "MatMPIAIJSetPreallocation_C", &aij));
3750cd8b6e2SStefano Zampini   PetscCall(PetscObjectHasFunction((PetscObject)A, "MatISSetPreallocation_C", &is));
3760cd8b6e2SStefano Zampini   PetscCall(PetscObjectHasFunction((PetscObject)A, "MatHYPRESetPreallocation_C", &hyp));
3770cd8b6e2SStefano Zampini   if (!aij && !is && !hyp) PetscCall(PetscObjectHasFunction((PetscObject)A, "MatSeqAIJSetPreallocation_C", &aij));
378990279feSStefano Zampini   if (aij || is || hyp) {
37941319c1dSStefano Zampini     if (bs == cbs && bs == 1) {
3809566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJSetPreallocation(A, 0, dnnz));
3819566063dSJacob Faibussowitsch       PetscCall(MatMPIAIJSetPreallocation(A, 0, dnnz, 0, onnz));
3829566063dSJacob Faibussowitsch       PetscCall(MatISSetPreallocation(A, 0, dnnz, 0, onnz));
383990279feSStefano Zampini #if defined(PETSC_HAVE_HYPRE)
3849566063dSJacob Faibussowitsch       PetscCall(MatHYPRESetPreallocation(A, 0, dnnz, 0, onnz));
385990279feSStefano Zampini #endif
3863e5f4774SJed Brown     } else { /* Convert block-row precallocation to scalar-row */
38763562e91SJed Brown       PetscInt i, m, *sdnnz, *sonnz;
3889566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(A, &m, NULL));
3899566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2((!!dnnz) * m, &sdnnz, (!!onnz) * m, &sonnz));
390dec54756SJed Brown       for (i = 0; i < m; i++) {
39141319c1dSStefano Zampini         if (dnnz) sdnnz[i] = dnnz[i / bs] * cbs;
39241319c1dSStefano Zampini         if (onnz) sonnz[i] = onnz[i / bs] * cbs;
39363562e91SJed Brown       }
3949566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJSetPreallocation(A, 0, dnnz ? sdnnz : NULL));
3959566063dSJacob Faibussowitsch       PetscCall(MatMPIAIJSetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL));
3969566063dSJacob Faibussowitsch       PetscCall(MatISSetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL));
397990279feSStefano Zampini #if defined(PETSC_HAVE_HYPRE)
3989566063dSJacob Faibussowitsch       PetscCall(MatHYPRESetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL));
399990279feSStefano Zampini #endif
4009566063dSJacob Faibussowitsch       PetscCall(PetscFree2(sdnnz, sonnz));
40163562e91SJed Brown     }
40263562e91SJed Brown   }
4033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
40463562e91SJed Brown }
40563562e91SJed Brown 
40602218ef1SBarry Smith /*@C
40702218ef1SBarry Smith   MatHeaderMerge - Merges some information from the header of `C` to `A`; the `C` object is then destroyed
408d0f46423SBarry Smith 
40902218ef1SBarry Smith   Collective, No Fortran Support
41002218ef1SBarry Smith 
41102218ef1SBarry Smith   Input Parameters:
41202218ef1SBarry Smith + A - a `Mat` being merged into
41302218ef1SBarry Smith - C - the `Mat` providing the merge information
41402218ef1SBarry Smith 
41502218ef1SBarry Smith   Level: developer
41602218ef1SBarry Smith 
417f2e50d50SStefano Zampini   Notes:
418f2e50d50SStefano Zampini   `A` and `C` must be of the same type.
419f2e50d50SStefano Zampini   The object list and query function list in `A` are retained, as well as the object name, and prefix.
420f2e50d50SStefano Zampini   The object state of `A` is increased by 1.
421f2e50d50SStefano Zampini 
42202218ef1SBarry Smith   Developer Note:
42302218ef1SBarry Smith   This is somewhat different from `MatHeaderReplace()`, it would be nice to merge the code
42402218ef1SBarry Smith 
42502218ef1SBarry Smith .seealso: `Mat`, `MatHeaderReplace()`
42602218ef1SBarry Smith  @*/
MatHeaderMerge(Mat A,Mat * C)427d71ae5a4SJacob Faibussowitsch PetscErrorCode MatHeaderMerge(Mat A, Mat *C)
428d71ae5a4SJacob Faibussowitsch {
429d44834fbSBarry Smith   PetscInt          refct;
43073107ff1SLisandro Dalcin   PetscOps          Abops;
43173107ff1SLisandro Dalcin   struct _MatOps    Aops;
4324768301cSVaclav Hapla   char             *mtype, *mname, *mprefix;
4334222ddf1SHong Zhang   Mat_Product      *product;
43433e6eea4SJose E. Roman   Mat_Redundant    *redundant;
435d4a972cbSStefano Zampini   PetscObjectState  state;
436f2e50d50SStefano Zampini   PetscObjectList   olist;
437f2e50d50SStefano Zampini   PetscFunctionList qlist;
438273d9f13SBarry Smith 
439273d9f13SBarry Smith   PetscFunctionBegin;
4401dc04de0SStefano Zampini   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
4411dc04de0SStefano Zampini   PetscValidHeaderSpecific(*C, MAT_CLASSID, 2);
4423ba16761SJacob Faibussowitsch   if (A == *C) PetscFunctionReturn(PETSC_SUCCESS);
443f2e50d50SStefano Zampini   PetscCheckSameTypeAndComm(A, 1, *C, 2);
444273d9f13SBarry Smith   /* save the parts of A we need */
44573107ff1SLisandro Dalcin   Abops     = ((PetscObject)A)->bops[0];
44673107ff1SLisandro Dalcin   Aops      = A->ops[0];
4477adad957SLisandro Dalcin   refct     = ((PetscObject)A)->refct;
4485c9eb25fSBarry Smith   mtype     = ((PetscObject)A)->type_name;
4495c9eb25fSBarry Smith   mname     = ((PetscObject)A)->name;
450d4a972cbSStefano Zampini   state     = ((PetscObject)A)->state;
4514768301cSVaclav Hapla   mprefix   = ((PetscObject)A)->prefix;
4524222ddf1SHong Zhang   product   = A->product;
45333e6eea4SJose E. Roman   redundant = A->redundant;
454f2e50d50SStefano Zampini   qlist     = ((PetscObject)A)->qlist;
455f2e50d50SStefano Zampini   olist     = ((PetscObject)A)->olist;
45630735b05SKris Buschelman 
4575c9eb25fSBarry Smith   /* zero these so the destroy below does not free them */
458f4259b30SLisandro Dalcin   ((PetscObject)A)->type_name = NULL;
459f4259b30SLisandro Dalcin   ((PetscObject)A)->name      = NULL;
460f2e50d50SStefano Zampini   ((PetscObject)A)->qlist     = NULL;
461f2e50d50SStefano Zampini   ((PetscObject)A)->olist     = NULL;
4625c9eb25fSBarry Smith 
463dbbe0bcdSBarry Smith   /*
464dbbe0bcdSBarry Smith      free all the interior data structures from mat
465dbbe0bcdSBarry Smith      cannot use PetscUseTypeMethod(A,destroy); because compiler
466dbbe0bcdSBarry Smith      thinks it may print NULL type_name and name
467dbbe0bcdSBarry Smith   */
468dbbe0bcdSBarry Smith   PetscTryTypeMethod(A, destroy);
4697c99f97cSSatish Balay 
4709566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->defaultvectype));
4713faff063SStefano Zampini   PetscCall(PetscFree(A->defaultrandtype));
4729566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&A->rmap));
4739566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&A->cmap));
4749566063dSJacob Faibussowitsch   PetscCall(PetscComposedQuantitiesDestroy((PetscObject)A));
475273d9f13SBarry Smith 
476273d9f13SBarry Smith   /* copy C over to A */
47726cc229bSBarry Smith   PetscCall(PetscFree(A->factorprefix));
4789566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(A, *C, sizeof(struct _p_Mat)));
479273d9f13SBarry Smith 
480273d9f13SBarry Smith   /* return the parts of A we saved */
48173107ff1SLisandro Dalcin   ((PetscObject)A)->bops[0]   = Abops;
48273107ff1SLisandro Dalcin   A->ops[0]                   = Aops;
4837adad957SLisandro Dalcin   ((PetscObject)A)->refct     = refct;
4847adad957SLisandro Dalcin   ((PetscObject)A)->type_name = mtype;
4857adad957SLisandro Dalcin   ((PetscObject)A)->name      = mname;
4864768301cSVaclav Hapla   ((PetscObject)A)->prefix    = mprefix;
487d4a972cbSStefano Zampini   ((PetscObject)A)->state     = state + 1;
4884222ddf1SHong Zhang   A->product                  = product;
48933e6eea4SJose E. Roman   A->redundant                = redundant;
490273d9f13SBarry Smith 
491f2e50d50SStefano Zampini   /* Append the saved lists */
492f2e50d50SStefano Zampini   PetscCall(PetscFunctionListDuplicate(qlist, &((PetscObject)A)->qlist));
493f2e50d50SStefano Zampini   PetscCall(PetscObjectListDuplicate(olist, &((PetscObject)A)->olist));
494f2e50d50SStefano Zampini   PetscCall(PetscFunctionListDestroy(&qlist));
495f2e50d50SStefano Zampini   PetscCall(PetscObjectListDestroy(&olist));
496f2e50d50SStefano Zampini 
4975c9eb25fSBarry Smith   /* since these two are copied into A we do not want them destroyed in C */
498f4259b30SLisandro Dalcin   ((PetscObject)*C)->qlist = NULL;
499f4259b30SLisandro Dalcin   ((PetscObject)*C)->olist = NULL;
5009566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(C));
5013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
502273d9f13SBarry Smith }
503d0f46423SBarry Smith 
50402218ef1SBarry Smith /*@
50502218ef1SBarry Smith   MatHeaderReplace - Replaces the internal data of matrix `A` by the internal data of matrix `C` while deleting the outer wrapper of `C`
506eb6b5d47SBarry Smith 
50702218ef1SBarry Smith   Input Parameters:
50802218ef1SBarry Smith + A - a `Mat` whose internal data is to be replaced
50902218ef1SBarry Smith - C - the `Mat` providing new internal data for `A`
510b30237c6SBarry Smith 
51102218ef1SBarry Smith   Level: advanced
51202218ef1SBarry Smith 
51302218ef1SBarry Smith   Example Usage\:
51402218ef1SBarry Smith .vb
51502218ef1SBarry Smith   Mat C;
51602218ef1SBarry Smith   MatCreateSeqAIJWithArrays(..., &C);
51702218ef1SBarry Smith   MatHeaderReplace(A, &C);
51802218ef1SBarry Smith   // C has been destroyed and A contains the matrix entries of C
51902218ef1SBarry Smith .ve
52002218ef1SBarry Smith 
52102218ef1SBarry Smith   Note:
5220b4b7b1cSBarry Smith   This can be used inside a function provided to `SNESSetJacobian()`, `TSSetRHSJacobian()`, or `TSSetIJacobian()` in cases where the user code
5230b4b7b1cSBarry Smith   computes an entirely new sparse matrix  (generally with a different matrix nonzero structure/pattern) for each Newton update.
5240b4b7b1cSBarry Smith   It is usually better to reuse the matrix nonzero structure of `A` instead of constructing an entirely new one.
52502218ef1SBarry Smith 
52602218ef1SBarry Smith   Developer Note:
52702218ef1SBarry Smith   This is somewhat different from `MatHeaderMerge()` it would be nice to merge the code
52802218ef1SBarry Smith 
52902218ef1SBarry Smith .seealso: `Mat`, `MatHeaderMerge()`
53002218ef1SBarry Smith  @*/
MatHeaderReplace(Mat A,Mat * C)53102218ef1SBarry Smith PetscErrorCode MatHeaderReplace(Mat A, Mat *C)
532d71ae5a4SJacob Faibussowitsch {
53327b31e29SJed Brown   PetscInt         refct;
534fefd9316SJose E. Roman   PetscObjectState state;
53528be2f97SBarry Smith   struct _p_Mat    buffer;
53681fa06acSBarry Smith   MatStencilInfo   stencil;
5378ab5b326SKris Buschelman 
5388ab5b326SKris Buschelman   PetscFunctionBegin;
53927b31e29SJed Brown   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
54028be2f97SBarry Smith   PetscValidHeaderSpecific(*C, MAT_CLASSID, 2);
5413ba16761SJacob Faibussowitsch   if (A == *C) PetscFunctionReturn(PETSC_SUCCESS);
54228be2f97SBarry Smith   PetscCheckSameComm(A, 1, *C, 2);
543aed4548fSBarry Smith   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);
5446d7c1e57SBarry Smith 
54528be2f97SBarry Smith   /* swap C and A */
54627b31e29SJed Brown   refct   = ((PetscObject)A)->refct;
547fefd9316SJose E. Roman   state   = ((PetscObject)A)->state;
54881fa06acSBarry Smith   stencil = A->stencil;
5499566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(&buffer, A, sizeof(struct _p_Mat)));
5509566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(A, *C, sizeof(struct _p_Mat)));
5519566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(*C, &buffer, sizeof(struct _p_Mat)));
55227b31e29SJed Brown   ((PetscObject)A)->refct = refct;
553fefd9316SJose E. Roman   ((PetscObject)A)->state = state + 1;
55481fa06acSBarry Smith   A->stencil              = stencil;
55526fbe8dcSKarl Rupp 
556c32d4117SBarry Smith   ((PetscObject)*C)->refct = 1;
5579566063dSJacob Faibussowitsch   PetscCall(MatDestroy(C));
5583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5598ab5b326SKris Buschelman }
560e7e92044SBarry Smith 
561e7e92044SBarry Smith /*@
562b470e4b4SRichard Tran Mills   MatBindToCPU - marks a matrix to temporarily stay on the CPU and perform computations on the CPU
563e7e92044SBarry Smith 
5642ef1f0ffSBarry Smith   Logically Collective
5652216c58aSStefano Zampini 
566e7e92044SBarry Smith   Input Parameters:
567e7e92044SBarry Smith + A   - the matrix
56811a5261eSBarry Smith - flg - bind to the CPU if value of `PETSC_TRUE`
569e7e92044SBarry Smith 
57090ea27d8SSatish Balay   Level: intermediate
5712216c58aSStefano Zampini 
5721cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatBoundToCPU()`
573e7e92044SBarry Smith @*/
MatBindToCPU(Mat A,PetscBool flg)574d71ae5a4SJacob Faibussowitsch PetscErrorCode MatBindToCPU(Mat A, PetscBool flg)
575d71ae5a4SJacob Faibussowitsch {
5767d871021SStefano Zampini   PetscFunctionBegin;
5772ffa8ee7SStefano Zampini   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
5782ffa8ee7SStefano Zampini   PetscValidLogicalCollectiveBool(A, flg, 2);
5792216c58aSStefano Zampini #if defined(PETSC_HAVE_DEVICE)
5803ba16761SJacob Faibussowitsch   if (A->boundtocpu == flg) PetscFunctionReturn(PETSC_SUCCESS);
581b470e4b4SRichard Tran Mills   A->boundtocpu = flg;
582dbbe0bcdSBarry Smith   PetscTryTypeMethod(A, bindtocpu, flg);
5832216c58aSStefano Zampini #endif
5843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5852216c58aSStefano Zampini }
5862216c58aSStefano Zampini 
5872216c58aSStefano Zampini /*@
5882216c58aSStefano Zampini   MatBoundToCPU - query if a matrix is bound to the CPU
5892216c58aSStefano Zampini 
5902216c58aSStefano Zampini   Input Parameter:
5912216c58aSStefano Zampini . A - the matrix
5922216c58aSStefano Zampini 
5932216c58aSStefano Zampini   Output Parameter:
5942216c58aSStefano Zampini . flg - the logical flag
5952216c58aSStefano Zampini 
5962216c58aSStefano Zampini   Level: intermediate
5972216c58aSStefano Zampini 
5981cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatBindToCPU()`
5992216c58aSStefano Zampini @*/
MatBoundToCPU(Mat A,PetscBool * flg)600d71ae5a4SJacob Faibussowitsch PetscErrorCode MatBoundToCPU(Mat A, PetscBool *flg)
601d71ae5a4SJacob Faibussowitsch {
6022ffa8ee7SStefano Zampini   PetscFunctionBegin;
6032ffa8ee7SStefano Zampini   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
6044f572ea9SToby Isaac   PetscAssertPointer(flg, 2);
6052216c58aSStefano Zampini #if defined(PETSC_HAVE_DEVICE)
6062216c58aSStefano Zampini   *flg = A->boundtocpu;
6072216c58aSStefano Zampini #else
6082216c58aSStefano Zampini   *flg = PETSC_TRUE;
6097d871021SStefano Zampini #endif
6103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
611e7e92044SBarry Smith }
6127e8381f9SStefano Zampini 
MatSetValuesCOO_Basic(Mat A,const PetscScalar coo_v[],InsertMode imode)613d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetValuesCOO_Basic(Mat A, const PetscScalar coo_v[], InsertMode imode)
614d71ae5a4SJacob Faibussowitsch {
6157e8381f9SStefano Zampini   IS              is_coo_i, is_coo_j;
6167e8381f9SStefano Zampini   const PetscInt *coo_i, *coo_j;
6177e8381f9SStefano Zampini   PetscInt        n, n_i, n_j;
6187e8381f9SStefano Zampini   PetscScalar     zero = 0.;
6197e8381f9SStefano Zampini 
6207e8381f9SStefano Zampini   PetscFunctionBegin;
6219566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_coo_i", (PetscObject *)&is_coo_i));
6229566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_coo_j", (PetscObject *)&is_coo_j));
62328b400f6SJacob Faibussowitsch   PetscCheck(is_coo_i, PetscObjectComm((PetscObject)A), PETSC_ERR_COR, "Missing coo_i IS");
62428b400f6SJacob Faibussowitsch   PetscCheck(is_coo_j, PetscObjectComm((PetscObject)A), PETSC_ERR_COR, "Missing coo_j IS");
6259566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is_coo_i, &n_i));
6269566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is_coo_j, &n_j));
62708401ef6SPierre Jolivet   PetscCheck(n_i == n_j, PETSC_COMM_SELF, PETSC_ERR_COR, "Wrong local size %" PetscInt_FMT " != %" PetscInt_FMT, n_i, n_j);
6289566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is_coo_i, &coo_i));
6299566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is_coo_j, &coo_j));
63048a46eb9SPierre Jolivet   if (imode != ADD_VALUES) PetscCall(MatZeroEntries(A));
63148a46eb9SPierre Jolivet   for (n = 0; n < n_i; n++) PetscCall(MatSetValue(A, coo_i[n], coo_j[n], coo_v ? coo_v[n] : zero, ADD_VALUES));
6329566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is_coo_i, &coo_i));
6339566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is_coo_j, &coo_j));
6343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6357e8381f9SStefano Zampini }
6367e8381f9SStefano Zampini 
MatSetPreallocationCOO_Basic(Mat A,PetscCount ncoo,PetscInt coo_i[],PetscInt coo_j[])6378063e3c8SPierre Jolivet PetscErrorCode MatSetPreallocationCOO_Basic(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[])
638d71ae5a4SJacob Faibussowitsch {
6397e8381f9SStefano Zampini   Mat         preallocator;
6407e8381f9SStefano Zampini   IS          is_coo_i, is_coo_j;
641835f2295SStefano Zampini   PetscInt    ncoo_i;
6427e8381f9SStefano Zampini   PetscScalar zero = 0.0;
6437e8381f9SStefano Zampini 
6447e8381f9SStefano Zampini   PetscFunctionBegin;
645835f2295SStefano Zampini   PetscCall(PetscIntCast(ncoo, &ncoo_i));
6469566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->rmap));
6479566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->cmap));
6489566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &preallocator));
6499566063dSJacob Faibussowitsch   PetscCall(MatSetType(preallocator, MATPREALLOCATOR));
6509566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(preallocator, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
6519566063dSJacob Faibussowitsch   PetscCall(MatSetLayouts(preallocator, A->rmap, A->cmap));
6529566063dSJacob Faibussowitsch   PetscCall(MatSetUp(preallocator));
65348a46eb9SPierre Jolivet   for (PetscCount n = 0; n < ncoo; n++) PetscCall(MatSetValue(preallocator, coo_i[n], coo_j[n], zero, INSERT_VALUES));
6549566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY));
6559566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY));
6569566063dSJacob Faibussowitsch   PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_TRUE, A));
6579566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&preallocator));
658835f2295SStefano Zampini   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ncoo_i, coo_i, PETSC_COPY_VALUES, &is_coo_i));
659835f2295SStefano Zampini   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ncoo_i, coo_j, PETSC_COPY_VALUES, &is_coo_j));
6609566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc_coo_i", (PetscObject)is_coo_i));
6619566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc_coo_j", (PetscObject)is_coo_j));
6629566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_coo_i));
6639566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_coo_j));
6643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6657e8381f9SStefano Zampini }
6667e8381f9SStefano Zampini 
66756856777SBarry Smith /*@C
668c3dd2894SJed Brown   MatSetPreallocationCOO - set preallocation for matrices using a coordinate format of the entries with global indices
6697e8381f9SStefano Zampini 
670c3339decSBarry Smith   Collective
6717e8381f9SStefano Zampini 
6724165533cSJose E. Roman   Input Parameters:
6737e8381f9SStefano Zampini + A     - matrix being preallocated
67442550becSJunchao Zhang . ncoo  - number of entries
6757e8381f9SStefano Zampini . coo_i - row indices
6767e8381f9SStefano Zampini - coo_j - column indices
6777e8381f9SStefano Zampini 
6787e8381f9SStefano Zampini   Level: beginner
6797e8381f9SStefano Zampini 
680394ed5ebSJunchao Zhang   Notes:
681f13dfd9eSBarry Smith   The indices within `coo_i` and `coo_j` may be modified within this function. The caller should not rely on them
682e8729f6fSJunchao Zhang   having any specific value after this function returns. The arrays can be freed or reused immediately
683e8729f6fSJunchao Zhang   after this function returns.
684e8729f6fSJunchao Zhang 
68511a5261eSBarry Smith   Entries can be repeated, see `MatSetValuesCOO()`. Entries with negative row or column indices are allowed
68611a5261eSBarry Smith   but will be ignored. The corresponding entries in `MatSetValuesCOO()` will be ignored too. Remote entries
687d7547e51SJunchao Zhang   are allowed and will be properly added or inserted to the matrix, unless the matrix option `MAT_IGNORE_OFF_PROC_ENTRIES`
68811a5261eSBarry Smith   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.
6897e8381f9SStefano Zampini 
690d7547e51SJunchao Zhang   If you just want to create a sequential AIJ matrix (`MATSEQAIJ`), and your matrix entries in COO format are unique, you can also use
691d7547e51SJunchao Zhang   `MatCreateSeqAIJFromTriple()`. But that is not recommended for iterative applications.
692d7547e51SJunchao Zhang 
6931cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`,
6942ef1f0ffSBarry Smith           `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOOLocal()`,
6952ef1f0ffSBarry Smith           `DMSetMatrixPreallocateSkip()`, `MatCreateSeqAIJFromTriple()`
6967e8381f9SStefano Zampini @*/
MatSetPreallocationCOO(Mat A,PetscCount ncoo,PetscInt coo_i[],PetscInt coo_j[])697d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetPreallocationCOO(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[])
698d71ae5a4SJacob Faibussowitsch {
6998063e3c8SPierre Jolivet   PetscErrorCode (*f)(Mat, PetscCount, PetscInt[], PetscInt[]) = NULL;
7007e8381f9SStefano Zampini 
7017e8381f9SStefano Zampini   PetscFunctionBegin;
7027e8381f9SStefano Zampini   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
7037e8381f9SStefano Zampini   PetscValidType(A, 1);
7044f572ea9SToby Isaac   if (ncoo) PetscAssertPointer(coo_i, 3);
7054f572ea9SToby Isaac   if (ncoo) PetscAssertPointer(coo_j, 4);
7069566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->rmap));
7079566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->cmap));
7089566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetPreallocationCOO_C", &f));
709cbc6b225SStefano Zampini 
7109566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_PreallCOO, A, 0, 0, 0));
7117e8381f9SStefano Zampini   if (f) {
7129566063dSJacob Faibussowitsch     PetscCall((*f)(A, ncoo, coo_i, coo_j));
7137e8381f9SStefano Zampini   } else { /* allow fallback, very slow */
7149566063dSJacob Faibussowitsch     PetscCall(MatSetPreallocationCOO_Basic(A, ncoo, coo_i, coo_j));
7157e8381f9SStefano Zampini   }
7169566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_PreallCOO, A, 0, 0, 0));
7176834774dSStefano Zampini   A->preallocated = PETSC_TRUE;
718cbc6b225SStefano Zampini   A->nonzerostate++;
7193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7207e8381f9SStefano Zampini }
7217e8381f9SStefano Zampini 
72256856777SBarry Smith /*@C
723c3dd2894SJed Brown   MatSetPreallocationCOOLocal - set preallocation for matrices using a coordinate format of the entries with local indices
724c3dd2894SJed Brown 
725c3339decSBarry Smith   Collective
726c3dd2894SJed Brown 
727c3dd2894SJed Brown   Input Parameters:
728c3dd2894SJed Brown + A     - matrix being preallocated
729c3dd2894SJed Brown . ncoo  - number of entries
730c3dd2894SJed Brown . coo_i - row indices (local numbering; may be modified)
731c3dd2894SJed Brown - coo_j - column indices (local numbering; may be modified)
732c3dd2894SJed Brown 
733c3dd2894SJed Brown   Level: beginner
734c3dd2894SJed Brown 
735c3dd2894SJed Brown   Notes:
73611a5261eSBarry Smith   The local indices are translated using the local to global mapping, thus `MatSetLocalToGlobalMapping()` must have been
73711a5261eSBarry Smith   called prior to this function. For matrices created with `DMCreateMatrix()` the local to global mapping is often already provided.
738c3dd2894SJed Brown 
7392ef1f0ffSBarry Smith   The indices `coo_i` and `coo_j` may be modified within this function. They might be translated to corresponding global
740735d7f90SBarry Smith   indices, but the caller should not rely on them having any specific value after this function returns. The arrays
741735d7f90SBarry Smith   can be freed or reused immediately after this function returns.
742c3dd2894SJed Brown 
74311a5261eSBarry Smith   Entries can be repeated, see `MatSetValuesCOO()`. Entries with negative row or column indices are allowed
74411a5261eSBarry Smith   but will be ignored. The corresponding entries in `MatSetValuesCOO()` will be ignored too. Remote entries
745394ed5ebSJunchao Zhang   are allowed and will be properly added or inserted to the matrix.
746c3dd2894SJed Brown 
7471cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`,
7482ef1f0ffSBarry Smith           `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOO()`,
7492ef1f0ffSBarry Smith           `DMSetMatrixPreallocateSkip()`
750c3dd2894SJed Brown @*/
MatSetPreallocationCOOLocal(Mat A,PetscCount ncoo,PetscInt coo_i[],PetscInt coo_j[])751d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetPreallocationCOOLocal(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[])
752d71ae5a4SJacob Faibussowitsch {
7536834774dSStefano Zampini   PetscErrorCode (*f)(Mat, PetscCount, PetscInt[], PetscInt[]) = NULL;
754c3dd2894SJed Brown 
755c3dd2894SJed Brown   PetscFunctionBegin;
756c3dd2894SJed Brown   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
757c3dd2894SJed Brown   PetscValidType(A, 1);
7584f572ea9SToby Isaac   if (ncoo) PetscAssertPointer(coo_i, 3);
7594f572ea9SToby Isaac   if (ncoo) PetscAssertPointer(coo_j, 4);
7609566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->rmap));
7619566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->cmap));
762cbc6b225SStefano Zampini 
7639566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetPreallocationCOOLocal_C", &f));
7646834774dSStefano Zampini   if (f) {
7659566063dSJacob Faibussowitsch     PetscCall((*f)(A, ncoo, coo_i, coo_j));
766cbc6b225SStefano Zampini     A->nonzerostate++;
7676834774dSStefano Zampini   } else {
768835f2295SStefano Zampini     PetscInt               ncoo_i;
769cbc6b225SStefano Zampini     ISLocalToGlobalMapping ltog_row, ltog_col;
7706497c311SBarry Smith 
7719566063dSJacob Faibussowitsch     PetscCall(MatGetLocalToGlobalMapping(A, &ltog_row, &ltog_col));
772835f2295SStefano Zampini     if (ltog_row) {
773835f2295SStefano Zampini       PetscCall(PetscIntCast(ncoo, &ncoo_i));
774835f2295SStefano Zampini       PetscCall(ISLocalToGlobalMappingApply(ltog_row, ncoo_i, coo_i, coo_i));
775835f2295SStefano Zampini     }
776835f2295SStefano Zampini     if (ltog_col) {
777835f2295SStefano Zampini       PetscCall(PetscIntCast(ncoo, &ncoo_i));
778835f2295SStefano Zampini       PetscCall(ISLocalToGlobalMappingApply(ltog_col, ncoo_i, coo_j, coo_j));
779835f2295SStefano Zampini     }
7809566063dSJacob Faibussowitsch     PetscCall(MatSetPreallocationCOO(A, ncoo, coo_i, coo_j));
7816834774dSStefano Zampini   }
7826834774dSStefano Zampini   A->preallocated = PETSC_TRUE;
7833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
784c3dd2894SJed Brown }
785c3dd2894SJed Brown 
786c3dd2894SJed Brown /*@
78711a5261eSBarry Smith   MatSetValuesCOO - set values at once in a matrix preallocated using `MatSetPreallocationCOO()`
7887e8381f9SStefano Zampini 
789c3339decSBarry Smith   Collective
7907e8381f9SStefano Zampini 
7914165533cSJose E. Roman   Input Parameters:
7927e8381f9SStefano Zampini + A     - matrix being preallocated
7932ef1f0ffSBarry Smith . coo_v - the matrix values (can be `NULL`)
7947e8381f9SStefano Zampini - imode - the insert mode
7957e8381f9SStefano Zampini 
7967e8381f9SStefano Zampini   Level: beginner
7977e8381f9SStefano Zampini 
79811a5261eSBarry Smith   Notes:
79911a5261eSBarry Smith   The values must follow the order of the indices prescribed with `MatSetPreallocationCOO()` or `MatSetPreallocationCOOLocal()`.
80011a5261eSBarry Smith 
8012ef1f0ffSBarry Smith   When repeated entries are specified in the COO indices the `coo_v` values are first properly summed, regardless of the value of imode.
80211a5261eSBarry Smith   The imode flag indicates if coo_v must be added to the current values of the matrix (`ADD_VALUES`) or overwritten (`INSERT_VALUES`).
80311a5261eSBarry Smith 
80411a5261eSBarry Smith   `MatAssemblyBegin()` and `MatAssemblyEnd()` do not need to be called after this routine. It automatically handles the assembly process.
8057e8381f9SStefano Zampini 
8061cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatSetPreallocationCOO()`, `MatSetPreallocationCOOLocal()`, `InsertMode`, `INSERT_VALUES`, `ADD_VALUES`
8077e8381f9SStefano Zampini @*/
MatSetValuesCOO(Mat A,const PetscScalar coo_v[],InsertMode imode)808d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetValuesCOO(Mat A, const PetscScalar coo_v[], InsertMode imode)
809d71ae5a4SJacob Faibussowitsch {
8107e8381f9SStefano Zampini   PetscErrorCode (*f)(Mat, const PetscScalar[], InsertMode) = NULL;
81135cef55dSJunchao Zhang   PetscBool oldFlg;
8127e8381f9SStefano Zampini 
8137e8381f9SStefano Zampini   PetscFunctionBegin;
8147e8381f9SStefano Zampini   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
8157e8381f9SStefano Zampini   PetscValidType(A, 1);
8167e8381f9SStefano Zampini   MatCheckPreallocated(A, 1);
817bfcc3627SStefano Zampini   PetscValidLogicalCollectiveEnum(A, imode, 3);
8189566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetValuesCOO_C", &f));
8199566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_SetVCOO, A, 0, 0, 0));
8207e8381f9SStefano Zampini   if (f) {
82135cef55dSJunchao Zhang     PetscCall((*f)(A, coo_v, imode)); // all known COO implementations do not use MatStash. They do their own off-proc communication
82235cef55dSJunchao Zhang     PetscCall(MatGetOption(A, MAT_NO_OFF_PROC_ENTRIES, &oldFlg));
82335cef55dSJunchao Zhang     PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); // set A->nooffprocentries to avoid costly MatStash scatter in MatAssembly
82435cef55dSJunchao Zhang   } else {
82535cef55dSJunchao Zhang     PetscCall(MatSetValuesCOO_Basic(A, coo_v, imode)); // fall back to MatSetValues, which might use MatStash
8267e8381f9SStefano Zampini   }
8279566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
8289566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
82935cef55dSJunchao Zhang   if (f) PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, oldFlg));
83035cef55dSJunchao Zhang   PetscCall(PetscLogEventEnd(MAT_SetVCOO, A, 0, 0, 0));
8313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8327e8381f9SStefano Zampini }
83365a9ecf2SRichard Tran Mills 
83465a9ecf2SRichard Tran Mills /*@
83565a9ecf2SRichard Tran Mills   MatSetBindingPropagates - Sets whether the state of being bound to the CPU for a GPU matrix type propagates to child and some other associated objects
83665a9ecf2SRichard Tran Mills 
83765a9ecf2SRichard Tran Mills   Input Parameters:
83865a9ecf2SRichard Tran Mills + A   - the matrix
83965a9ecf2SRichard Tran Mills - flg - flag indicating whether the boundtocpu flag should be propagated
84065a9ecf2SRichard Tran Mills 
84165a9ecf2SRichard Tran Mills   Level: developer
84265a9ecf2SRichard Tran Mills 
84365a9ecf2SRichard Tran Mills   Notes:
8442fe279fdSBarry Smith   If the value of flg is set to true, the following will occur
8452fe279fdSBarry Smith +   `MatCreateSubMatrices()` and `MatCreateRedundantMatrix()` - bind created matrices to CPU if the input matrix is bound to the CPU.
8462fe279fdSBarry Smith -   `MatCreateVecs()` - bind created vectors to CPU if the input matrix is bound to the CPU.
84765a9ecf2SRichard Tran Mills 
84865a9ecf2SRichard Tran Mills   The bindingpropagates flag itself is also propagated by the above routines.
84965a9ecf2SRichard Tran Mills 
850fe59aa6dSJacob Faibussowitsch   Developer Notes:
851aa624791SPierre Jolivet   If the fine-scale `DMDA` has the `-dm_bind_below` option set to true, then `DMCreateInterpolationScale()` calls `MatSetBindingPropagates()`
85265a9ecf2SRichard Tran Mills   on the restriction/interpolation operator to set the bindingpropagates flag to true.
85365a9ecf2SRichard Tran Mills 
8541cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `VecSetBindingPropagates()`, `MatGetBindingPropagates()`
85565a9ecf2SRichard Tran Mills @*/
MatSetBindingPropagates(Mat A,PetscBool flg)856d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetBindingPropagates(Mat A, PetscBool flg)
857d71ae5a4SJacob Faibussowitsch {
85865a9ecf2SRichard Tran Mills   PetscFunctionBegin;
85965a9ecf2SRichard Tran Mills   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
860d5e393b6SSuyash Tandon #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
86165a9ecf2SRichard Tran Mills   A->bindingpropagates = flg;
86265a9ecf2SRichard Tran Mills #endif
8633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
86465a9ecf2SRichard Tran Mills }
865e9c74fd6SRichard Tran Mills 
866e9c74fd6SRichard Tran Mills /*@
867e9c74fd6SRichard Tran Mills   MatGetBindingPropagates - Gets whether the state of being bound to the CPU for a GPU matrix type propagates to child and some other associated objects
868e9c74fd6SRichard Tran Mills 
869e9c74fd6SRichard Tran Mills   Input Parameter:
870e9c74fd6SRichard Tran Mills . A - the matrix
871e9c74fd6SRichard Tran Mills 
872e9c74fd6SRichard Tran Mills   Output Parameter:
873e9c74fd6SRichard Tran Mills . flg - flag indicating whether the boundtocpu flag will be propagated
874e9c74fd6SRichard Tran Mills 
875e9c74fd6SRichard Tran Mills   Level: developer
876e9c74fd6SRichard Tran Mills 
8771cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatSetBindingPropagates()`
878e9c74fd6SRichard Tran Mills @*/
MatGetBindingPropagates(Mat A,PetscBool * flg)879d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetBindingPropagates(Mat A, PetscBool *flg)
880d71ae5a4SJacob Faibussowitsch {
881e9c74fd6SRichard Tran Mills   PetscFunctionBegin;
882e9c74fd6SRichard Tran Mills   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
8834f572ea9SToby Isaac   PetscAssertPointer(flg, 2);
884d5e393b6SSuyash Tandon #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
885e9c74fd6SRichard Tran Mills   *flg = A->bindingpropagates;
886e9c74fd6SRichard Tran Mills #else
887e9c74fd6SRichard Tran Mills   *flg = PETSC_FALSE;
888e9c74fd6SRichard Tran Mills #endif
8893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
890e9c74fd6SRichard Tran Mills }
891