xref: /petsc/src/mat/utils/gcreate.c (revision b94d7ded0a05f1bbd5e48daa6f92b28259c75b44)
1af0996ceSBarry Smith #include <petsc/private/matimpl.h>       /*I "petscmat.h"  I*/
27807a1faSBarry Smith 
346533700Sstefano_zampini PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat mat,PetscInt rbs, PetscInt cbs)
446533700Sstefano_zampini {
546533700Sstefano_zampini   PetscFunctionBegin;
65c577a9aSstefano_zampini   if (!mat->preallocated) PetscFunctionReturn(0);
7aed4548fSBarry 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);
8aed4548fSBarry 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);
946533700Sstefano_zampini   PetscFunctionReturn(0);
1046533700Sstefano_zampini }
1146533700Sstefano_zampini 
12db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatShift_Basic(Mat Y,PetscScalar a)
137d68702bSBarry Smith {
147d68702bSBarry Smith   PetscInt       i,start,end;
157d68702bSBarry Smith   PetscScalar    alpha = a;
167d68702bSBarry Smith   PetscBool      prevoption;
177d68702bSBarry Smith 
187d68702bSBarry Smith   PetscFunctionBegin;
199566063dSJacob Faibussowitsch   PetscCall(MatGetOption(Y,MAT_NO_OFF_PROC_ENTRIES,&prevoption));
209566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Y,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE));
219566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Y,&start,&end));
227d68702bSBarry Smith   for (i=start; i<end; i++) {
23ab6153dcSStefano Zampini     if (i < Y->cmap->N) {
249566063dSJacob Faibussowitsch       PetscCall(MatSetValues(Y,1,&i,1,&i,&alpha,ADD_VALUES));
257d68702bSBarry Smith     }
26ab6153dcSStefano Zampini   }
279566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY));
289566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY));
299566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Y,MAT_NO_OFF_PROC_ENTRIES,prevoption));
307d68702bSBarry Smith   PetscFunctionReturn(0);
317d68702bSBarry Smith }
327d68702bSBarry Smith 
3305869f15SSatish Balay /*@
3469dd0797SLois Curfman McInnes    MatCreate - Creates a matrix where the type is determined
357e5f4302SBarry Smith    from either a call to MatSetType() or from the options database
367e5f4302SBarry Smith    with a call to MatSetFromOptions(). The default matrix type is
3769b1f4b7SBarry Smith    AIJ, using the routines MatCreateSeqAIJ() or MatCreateAIJ()
387e5f4302SBarry Smith    if you do not set a type in the options database. If you never
397e5f4302SBarry Smith    call MatSetType() or MatSetFromOptions() it will generate an
40f8ab6608SSatish Balay    error when you try to use the matrix.
4183e1b59cSLois Curfman McInnes 
42d083f849SBarry Smith    Collective
43cb13003dSBarry Smith 
44f69a0ea3SMatthew Knepley    Input Parameter:
45f69a0ea3SMatthew Knepley .  comm - MPI communicator
467807a1faSBarry Smith 
477807a1faSBarry Smith    Output Parameter:
48dc401e71SLois Curfman McInnes .  A - the matrix
49e0b365e2SLois Curfman McInnes 
50273d9f13SBarry Smith    Options Database Keys:
51273d9f13SBarry Smith +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
5269b1f4b7SBarry Smith .    -mat_type mpiaij   - AIJ type, uses MatCreateAIJ()
53273d9f13SBarry Smith .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
5469b1f4b7SBarry Smith .    -mat_type mpidense - dense type, uses MatCreateDense()
55273d9f13SBarry Smith .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
5669b1f4b7SBarry Smith -    -mat_type mpibaij  - block AIJ type, uses MatCreateBAIJ()
57e0b365e2SLois Curfman McInnes 
5883e1b59cSLois Curfman McInnes    Even More Options Database Keys:
5983e1b59cSLois Curfman McInnes    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
6083e1b59cSLois Curfman McInnes    for additional format-specific options.
61e0b365e2SLois Curfman McInnes 
62273d9f13SBarry Smith    Level: beginner
63273d9f13SBarry Smith 
64db781477SPatrick Sanan .seealso: `MatCreateSeqAIJ()`, `MatCreateAIJ()`,
65db781477SPatrick Sanan           `MatCreateSeqDense()`, `MatCreateDense()`,
66db781477SPatrick Sanan           `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`,
67db781477SPatrick Sanan           `MatCreateSeqSBAIJ()`, `MatCreateSBAIJ()`,
68db781477SPatrick Sanan           `MatConvert()`
69273d9f13SBarry Smith @*/
707087cfbeSBarry Smith PetscErrorCode  MatCreate(MPI_Comm comm,Mat *A)
71273d9f13SBarry Smith {
72273d9f13SBarry Smith   Mat            B;
73273d9f13SBarry Smith 
74273d9f13SBarry Smith   PetscFunctionBegin;
75f69a0ea3SMatthew Knepley   PetscValidPointer(A,2);
7697f1f81fSBarry Smith 
770298fd71SBarry Smith   *A = NULL;
789566063dSJacob Faibussowitsch   PetscCall(MatInitializePackage());
798ba1e511SMatthew Knepley 
809566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(B,MAT_CLASSID,"Mat","Matrix","Mat",comm,MatDestroy,MatView));
819566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm,&B->rmap));
829566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm,&B->cmap));
839566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(VECSTANDARD,&B->defaultvectype));
8426fbe8dcSKarl Rupp 
85*b94d7dedSBarry Smith   B->symmetric                   = PETSC_BOOL3_UNKNOWN;
86*b94d7dedSBarry Smith   B->hermitian                   = PETSC_BOOL3_UNKNOWN;
87*b94d7dedSBarry Smith   B->structurally_symmetric      = PETSC_BOOL3_UNKNOWN;
88*b94d7dedSBarry Smith   B->spd                         = PETSC_BOOL3_UNKNOWN;
89*b94d7dedSBarry Smith   B->symmetry_eternal            = PETSC_FALSE;
90*b94d7dedSBarry Smith   B->structural_symmetry_eternal = PETSC_FALSE;
91*b94d7dedSBarry Smith 
9294342113SStefano Zampini   B->congruentlayouts = PETSC_DECIDE;
93273d9f13SBarry Smith   B->preallocated     = PETSC_FALSE;
946f3d89d0SStefano Zampini #if defined(PETSC_HAVE_DEVICE)
956f3d89d0SStefano Zampini   B->boundtocpu       = PETSC_TRUE;
966f3d89d0SStefano Zampini #endif
97273d9f13SBarry Smith   *A                  = B;
98273d9f13SBarry Smith   PetscFunctionReturn(0);
99273d9f13SBarry Smith }
100273d9f13SBarry Smith 
101422a814eSBarry Smith /*@
10284d44b13SHong Zhang    MatSetErrorIfFailure - Causes Mat to generate an error, for example a zero pivot, is detected.
103422a814eSBarry Smith 
104422a814eSBarry Smith    Logically Collective on Mat
105422a814eSBarry Smith 
106422a814eSBarry Smith    Input Parameters:
107422a814eSBarry Smith +  mat -  matrix obtained from MatCreate()
108422a814eSBarry Smith -  flg - PETSC_TRUE indicates you want the error generated
109422a814eSBarry Smith 
110422a814eSBarry Smith    Level: advanced
111422a814eSBarry Smith 
112db781477SPatrick Sanan .seealso: `PCSetErrorIfFailure()`
113422a814eSBarry Smith @*/
11484d44b13SHong Zhang PetscErrorCode  MatSetErrorIfFailure(Mat mat,PetscBool flg)
115422a814eSBarry Smith {
116422a814eSBarry Smith   PetscFunctionBegin;
117422a814eSBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
118422a814eSBarry Smith   PetscValidLogicalCollectiveBool(mat,flg,2);
11984d44b13SHong Zhang   mat->erroriffailure = flg;
120422a814eSBarry Smith   PetscFunctionReturn(0);
121422a814eSBarry Smith }
122422a814eSBarry Smith 
123f69a0ea3SMatthew Knepley /*@
124f69a0ea3SMatthew Knepley   MatSetSizes - Sets the local and global sizes, and checks to determine compatibility
125f69a0ea3SMatthew Knepley 
126f69a0ea3SMatthew Knepley   Collective on Mat
127f69a0ea3SMatthew Knepley 
128f69a0ea3SMatthew Knepley   Input Parameters:
129f69a0ea3SMatthew Knepley +  A - the matrix
130f69a0ea3SMatthew Knepley .  m - number of local rows (or PETSC_DECIDE)
131f69a0ea3SMatthew Knepley .  n - number of local columns (or PETSC_DECIDE)
132f69a0ea3SMatthew Knepley .  M - number of global rows (or PETSC_DETERMINE)
133f69a0ea3SMatthew Knepley -  N - number of global columns (or PETSC_DETERMINE)
134f69a0ea3SMatthew Knepley 
135f69a0ea3SMatthew Knepley    Notes:
136f69a0ea3SMatthew Knepley    m (n) and M (N) cannot be both PETSC_DECIDE
137f69a0ea3SMatthew Knepley    If one processor calls this with M (N) of PETSC_DECIDE then all processors must, otherwise the program will hang.
138f69a0ea3SMatthew Knepley 
139f69a0ea3SMatthew Knepley    If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the
140f69a0ea3SMatthew Knepley    user must ensure that they are chosen to be compatible with the
141f69a0ea3SMatthew Knepley    vectors. To do this, one first considers the matrix-vector product
142f69a0ea3SMatthew Knepley    'y = A x'. The 'm' that is used in the above routine must match the
143f69a0ea3SMatthew Knepley    local size used in the vector creation routine VecCreateMPI() for 'y'.
144f69a0ea3SMatthew Knepley    Likewise, the 'n' used must match that used as the local size in
145f69a0ea3SMatthew Knepley    VecCreateMPI() for 'x'.
146f69a0ea3SMatthew Knepley 
147f73d5cc4SBarry Smith    You cannot change the sizes once they have been set.
148f73d5cc4SBarry Smith 
149f73d5cc4SBarry Smith    The sizes must be set before MatSetUp() or MatXXXSetPreallocation() is called.
150f73d5cc4SBarry Smith 
151f69a0ea3SMatthew Knepley   Level: beginner
152f69a0ea3SMatthew Knepley 
153db781477SPatrick Sanan .seealso: `MatGetSize()`, `PetscSplitOwnership()`
154f69a0ea3SMatthew Knepley @*/
1557087cfbeSBarry Smith PetscErrorCode  MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N)
156f69a0ea3SMatthew Knepley {
157f69a0ea3SMatthew Knepley   PetscFunctionBegin;
1580700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
159a69c7061SStefano Zampini   PetscValidLogicalCollectiveInt(A,M,4);
160a69c7061SStefano Zampini   PetscValidLogicalCollectiveInt(A,N,5);
161aed4548fSBarry 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);
162aed4548fSBarry 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);
163aed4548fSBarry Smith   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,A->rmap->n,A->rmap->N);
164aed4548fSBarry Smith   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,A->cmap->n,A->cmap->N);
165d0f46423SBarry Smith   A->rmap->n = m;
166d0f46423SBarry Smith   A->cmap->n = n;
16759cb773eSBarry Smith   A->rmap->N = M > -1 ? M : A->rmap->N;
16859cb773eSBarry Smith   A->cmap->N = N > -1 ? N : A->cmap->N;
169f69a0ea3SMatthew Knepley   PetscFunctionReturn(0);
170f69a0ea3SMatthew Knepley }
171f69a0ea3SMatthew Knepley 
17205869f15SSatish Balay /*@
173273d9f13SBarry Smith    MatSetFromOptions - Creates a matrix where the type is determined
174273d9f13SBarry Smith    from the options database. Generates a parallel MPI matrix if the
175273d9f13SBarry Smith    communicator has more than one processor.  The default matrix type is
17669b1f4b7SBarry Smith    AIJ, using the routines MatCreateSeqAIJ() and MatCreateAIJ() if
1777e5f4302SBarry Smith    you do not select a type in the options database.
178273d9f13SBarry Smith 
179273d9f13SBarry Smith    Collective on Mat
180273d9f13SBarry Smith 
181273d9f13SBarry Smith    Input Parameter:
182273d9f13SBarry Smith .  A - the matrix
183273d9f13SBarry Smith 
184273d9f13SBarry Smith    Options Database Keys:
185273d9f13SBarry Smith +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
18669b1f4b7SBarry Smith .    -mat_type mpiaij   - AIJ type, uses MatCreateAIJ()
187273d9f13SBarry Smith .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
18869b1f4b7SBarry Smith .    -mat_type mpidense - dense type, uses MatCreateDense()
189273d9f13SBarry Smith .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
19069b1f4b7SBarry Smith -    -mat_type mpibaij  - block AIJ type, uses MatCreateBAIJ()
191273d9f13SBarry Smith 
192273d9f13SBarry Smith    Even More Options Database Keys:
193273d9f13SBarry Smith    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
194273d9f13SBarry Smith    for additional format-specific options.
195bd9ce289SLois Curfman McInnes 
1961d69843bSLois Curfman McInnes    Level: beginner
1971d69843bSLois Curfman McInnes 
198db781477SPatrick Sanan .seealso: `MatCreateSeqAIJ(()`, `MatCreateAIJ()`,
199db781477SPatrick Sanan           `MatCreateSeqDense()`, `MatCreateDense()`,
200db781477SPatrick Sanan           `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`,
201db781477SPatrick Sanan           `MatCreateSeqSBAIJ()`, `MatCreateSBAIJ()`,
202db781477SPatrick Sanan           `MatConvert()`
2037807a1faSBarry Smith @*/
2047087cfbeSBarry Smith PetscErrorCode  MatSetFromOptions(Mat B)
2057807a1faSBarry Smith {
206f3be49caSLisandro Dalcin   const char     *deft = MATAIJ;
207f3be49caSLisandro Dalcin   char           type[256];
20869df5c0cSJed Brown   PetscBool      flg,set;
20916e04d98SRichard Tran Mills   PetscInt       bind_below = 0;
210dbb450caSBarry Smith 
2113a40ed3dSBarry Smith   PetscFunctionBegin;
2120700a824SBarry Smith   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
213f3be49caSLisandro Dalcin 
214d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)B);
215535b19f3SBarry Smith 
216535b19f3SBarry Smith   if (B->rmap->bs < 0) {
217535b19f3SBarry Smith     PetscInt newbs = -1;
2189566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-mat_block_size","Set the blocksize used to store the matrix","MatSetBlockSize",newbs,&newbs,&flg));
219535b19f3SBarry Smith     if (flg) {
2209566063dSJacob Faibussowitsch       PetscCall(PetscLayoutSetBlockSize(B->rmap,newbs));
2219566063dSJacob Faibussowitsch       PetscCall(PetscLayoutSetBlockSize(B->cmap,newbs));
222535b19f3SBarry Smith     }
223535b19f3SBarry Smith   }
224535b19f3SBarry Smith 
2259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-mat_type","Matrix type","MatSetType",MatList,deft,type,256,&flg));
226273d9f13SBarry Smith   if (flg) {
2279566063dSJacob Faibussowitsch     PetscCall(MatSetType(B,type));
228f3be49caSLisandro Dalcin   } else if (!((PetscObject)B)->type_name) {
2299566063dSJacob Faibussowitsch     PetscCall(MatSetType(B,deft));
230273d9f13SBarry Smith   }
231f3be49caSLisandro Dalcin 
2329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsName("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",&B->checksymmetryonassembly));
2339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",B->checksymmetrytol,&B->checksymmetrytol,NULL));
2349566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_null_space_test","Checks if provided null space is correct in MatAssemblyEnd()","MatSetNullSpaceTest",B->checknullspaceonassembly,&B->checknullspaceonassembly,NULL));
2359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_error_if_failure","Generate an error if an error occurs when factoring the matrix","MatSetErrorIfFailure",B->erroriffailure,&B->erroriffailure,NULL));
236840d65ccSBarry Smith 
2371baa6e33SBarry Smith   if (B->ops->setfromoptions) PetscCall((*B->ops->setfromoptions)(PetscOptionsObject,B));
238f3be49caSLisandro Dalcin 
23969df5c0cSJed Brown   flg  = PETSC_FALSE;
2409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_new_nonzero_location_err","Generate an error if new nonzeros are created in the matrix structure (useful to test preallocation)","MatSetOption",flg,&flg,&set));
2419566063dSJacob Faibussowitsch   if (set) PetscCall(MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,flg));
24269df5c0cSJed Brown   flg  = PETSC_FALSE;
2439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_new_nonzero_allocation_err","Generate an error if new nonzeros are allocated in the matrix structure (useful to test preallocation)","MatSetOption",flg,&flg,&set));
2449566063dSJacob Faibussowitsch   if (set) PetscCall(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,flg));
245478db826SMatthew G. Knepley   flg  = PETSC_FALSE;
2469566063dSJacob 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));
2479566063dSJacob Faibussowitsch   if (set) PetscCall(MatSetOption(B,MAT_IGNORE_ZERO_ENTRIES,flg));
24869df5c0cSJed Brown 
2491a2c6b5cSJunchao Zhang   flg  = PETSC_FALSE;
2509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_form_explicit_transpose","Hint to form an explicit transpose for operations like MatMultTranspose","MatSetOption",flg,&flg,&set));
2519566063dSJacob Faibussowitsch   if (set) PetscCall(MatSetOption(B,MAT_FORM_EXPLICIT_TRANSPOSE,flg));
2521a2c6b5cSJunchao Zhang 
25316e04d98SRichard Tran Mills   /* Bind to CPU if below a user-specified size threshold.
25416e04d98SRichard Tran Mills    * This perhaps belongs in the options for the GPU Mat types, but MatBindToCPU() does nothing when called on non-GPU types,
25516e04d98SRichard Tran Mills    * and putting it here makes is more maintainable than duplicating this for all. */
2569566063dSJacob 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));
25716e04d98SRichard Tran Mills   if (flg && B->rmap->n < bind_below) {
2589566063dSJacob Faibussowitsch     PetscCall(MatBindToCPU(B,PETSC_TRUE));
25916e04d98SRichard Tran Mills   }
26016e04d98SRichard Tran Mills 
2615d973c19SBarry Smith   /* process any options handlers added with PetscObjectAddOptionsHandler() */
2629566063dSJacob Faibussowitsch   PetscCall(PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)B));
263d0609cedSBarry Smith   PetscOptionsEnd();
2643a40ed3dSBarry Smith   PetscFunctionReturn(0);
2657807a1faSBarry Smith }
2667807a1faSBarry Smith 
267987010e7SBarry Smith /*@C
268e8bd9bafSStefano Zampini    MatXAIJSetPreallocation - set preallocation for serial and parallel AIJ, BAIJ, and SBAIJ matrices and their unassembled versions.
26963562e91SJed Brown 
27063562e91SJed Brown    Collective on Mat
27163562e91SJed Brown 
2724165533cSJose E. Roman    Input Parameters:
27363562e91SJed Brown +  A - matrix being preallocated
27463562e91SJed Brown .  bs - block size
27541319c1dSStefano Zampini .  dnnz - number of nonzero column blocks per block row of diagonal part of parallel matrix
27641319c1dSStefano Zampini .  onnz - number of nonzero column blocks per block row of off-diagonal part of parallel matrix
27741319c1dSStefano Zampini .  dnnzu - number of nonzero column blocks per block row of upper-triangular part of diagonal part of parallel matrix
27841319c1dSStefano Zampini -  onnzu - number of nonzero column blocks per block row of upper-triangular part of off-diagonal part of parallel matrix
27963562e91SJed Brown 
28063562e91SJed Brown    Level: beginner
28163562e91SJed Brown 
282db781477SPatrick Sanan .seealso: `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`,
283db781477SPatrick Sanan           `PetscSplitOwnership()`
28463562e91SJed Brown @*/
2854cce697cSJed Brown PetscErrorCode MatXAIJSetPreallocation(Mat A,PetscInt bs,const PetscInt dnnz[],const PetscInt onnz[],const PetscInt dnnzu[],const PetscInt onnzu[])
28663562e91SJed Brown {
28741319c1dSStefano Zampini   PetscInt       cbs;
28863562e91SJed Brown   void           (*aij)(void);
289e8bd9bafSStefano Zampini   void           (*is)(void);
290990279feSStefano Zampini   void           (*hyp)(void) = NULL;
29163562e91SJed Brown 
29263562e91SJed Brown   PetscFunctionBegin;
29341319c1dSStefano Zampini   if (bs != PETSC_DECIDE) { /* don't mess with an already set block size */
2949566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(A,bs));
29541319c1dSStefano Zampini   }
2969566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->rmap));
2979566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->cmap));
2989566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSizes(A,&bs,&cbs));
29941319c1dSStefano Zampini   /* these routines assumes bs == cbs, this should be checked somehow */
3009566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(A,bs,0,dnnz));
3019566063dSJacob Faibussowitsch   PetscCall(MatMPIBAIJSetPreallocation(A,bs,0,dnnz,0,onnz));
3029566063dSJacob Faibussowitsch   PetscCall(MatSeqSBAIJSetPreallocation(A,bs,0,dnnzu));
3039566063dSJacob Faibussowitsch   PetscCall(MatMPISBAIJSetPreallocation(A,bs,0,dnnzu,0,onnzu));
30463562e91SJed Brown   /*
305e8bd9bafSStefano 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
30663562e91SJed Brown     good before going on with it.
30763562e91SJed Brown   */
3089566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij));
3099566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)A,"MatISSetPreallocation_C",&is));
310990279feSStefano Zampini #if defined(PETSC_HAVE_HYPRE)
3119566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)A,"MatHYPRESetPreallocation_C",&hyp));
312990279feSStefano Zampini #endif
313990279feSStefano Zampini   if (!aij && !is && !hyp) {
3149566063dSJacob Faibussowitsch     PetscCall(PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij));
31563562e91SJed Brown   }
316990279feSStefano Zampini   if (aij || is || hyp) {
31741319c1dSStefano Zampini     if (bs == cbs && bs == 1) {
3189566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJSetPreallocation(A,0,dnnz));
3199566063dSJacob Faibussowitsch       PetscCall(MatMPIAIJSetPreallocation(A,0,dnnz,0,onnz));
3209566063dSJacob Faibussowitsch       PetscCall(MatISSetPreallocation(A,0,dnnz,0,onnz));
321990279feSStefano Zampini #if defined(PETSC_HAVE_HYPRE)
3229566063dSJacob Faibussowitsch       PetscCall(MatHYPRESetPreallocation(A,0,dnnz,0,onnz));
323990279feSStefano Zampini #endif
3243e5f4774SJed Brown     } else { /* Convert block-row precallocation to scalar-row */
32563562e91SJed Brown       PetscInt i,m,*sdnnz,*sonnz;
3269566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(A,&m,NULL));
3279566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2((!!dnnz)*m,&sdnnz,(!!onnz)*m,&sonnz));
328dec54756SJed Brown       for (i=0; i<m; i++) {
32941319c1dSStefano Zampini         if (dnnz) sdnnz[i] = dnnz[i/bs] * cbs;
33041319c1dSStefano Zampini         if (onnz) sonnz[i] = onnz[i/bs] * cbs;
33163562e91SJed Brown       }
3329566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL));
3339566063dSJacob Faibussowitsch       PetscCall(MatMPIAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL));
3349566063dSJacob Faibussowitsch       PetscCall(MatISSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL));
335990279feSStefano Zampini #if defined(PETSC_HAVE_HYPRE)
3369566063dSJacob Faibussowitsch       PetscCall(MatHYPRESetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL));
337990279feSStefano Zampini #endif
3389566063dSJacob Faibussowitsch       PetscCall(PetscFree2(sdnnz,sonnz));
33963562e91SJed Brown     }
34063562e91SJed Brown   }
34163562e91SJed Brown   PetscFunctionReturn(0);
34263562e91SJed Brown }
34363562e91SJed Brown 
344273d9f13SBarry Smith /*
345eb6b5d47SBarry Smith         Merges some information from Cs header to A; the C object is then destroyed
346d0f46423SBarry Smith 
347d0f46423SBarry Smith         This is somewhat different from MatHeaderReplace() it would be nice to merge the code
348273d9f13SBarry Smith */
34928be2f97SBarry Smith PetscErrorCode MatHeaderMerge(Mat A,Mat *C)
350273d9f13SBarry Smith {
351d44834fbSBarry Smith   PetscInt         refct;
35273107ff1SLisandro Dalcin   PetscOps         Abops;
35373107ff1SLisandro Dalcin   struct _MatOps   Aops;
3544768301cSVaclav Hapla   char             *mtype,*mname,*mprefix;
3554222ddf1SHong Zhang   Mat_Product      *product;
35633e6eea4SJose E. Roman   Mat_Redundant    *redundant;
357d4a972cbSStefano Zampini   PetscObjectState state;
358273d9f13SBarry Smith 
359273d9f13SBarry Smith   PetscFunctionBegin;
3601dc04de0SStefano Zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
3611dc04de0SStefano Zampini   PetscValidHeaderSpecific(*C,MAT_CLASSID,2);
3621dc04de0SStefano Zampini   if (A == *C) PetscFunctionReturn(0);
3631dc04de0SStefano Zampini   PetscCheckSameComm(A,1,*C,2);
364273d9f13SBarry Smith   /* save the parts of A we need */
36573107ff1SLisandro Dalcin   Abops     = ((PetscObject)A)->bops[0];
36673107ff1SLisandro Dalcin   Aops      = A->ops[0];
3677adad957SLisandro Dalcin   refct     = ((PetscObject)A)->refct;
3685c9eb25fSBarry Smith   mtype     = ((PetscObject)A)->type_name;
3695c9eb25fSBarry Smith   mname     = ((PetscObject)A)->name;
370d4a972cbSStefano Zampini   state     = ((PetscObject)A)->state;
3714768301cSVaclav Hapla   mprefix   = ((PetscObject)A)->prefix;
3724222ddf1SHong Zhang   product   = A->product;
37333e6eea4SJose E. Roman   redundant = A->redundant;
37430735b05SKris Buschelman 
3755c9eb25fSBarry Smith   /* zero these so the destroy below does not free them */
376f4259b30SLisandro Dalcin   ((PetscObject)A)->type_name = NULL;
377f4259b30SLisandro Dalcin   ((PetscObject)A)->name      = NULL;
3785c9eb25fSBarry Smith 
3797c99f97cSSatish Balay   /* free all the interior data structures from mat */
3809566063dSJacob Faibussowitsch   PetscCall((*A->ops->destroy)(A));
3817c99f97cSSatish Balay 
3829566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->defaultvectype));
3839566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&A->rmap));
3849566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&A->cmap));
3859566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&((PetscObject)A)->qlist));
3869566063dSJacob Faibussowitsch   PetscCall(PetscObjectListDestroy(&((PetscObject)A)->olist));
3879566063dSJacob Faibussowitsch   PetscCall(PetscComposedQuantitiesDestroy((PetscObject)A));
388273d9f13SBarry Smith 
389273d9f13SBarry Smith   /* copy C over to A */
39026cc229bSBarry Smith   PetscCall(PetscFree(A->factorprefix));
3919566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(A,*C,sizeof(struct _p_Mat)));
392273d9f13SBarry Smith 
393273d9f13SBarry Smith   /* return the parts of A we saved */
39473107ff1SLisandro Dalcin   ((PetscObject)A)->bops[0]   = Abops;
39573107ff1SLisandro Dalcin   A->ops[0]                   = Aops;
3967adad957SLisandro Dalcin   ((PetscObject)A)->refct     = refct;
3977adad957SLisandro Dalcin   ((PetscObject)A)->type_name = mtype;
3987adad957SLisandro Dalcin   ((PetscObject)A)->name      = mname;
3994768301cSVaclav Hapla   ((PetscObject)A)->prefix    = mprefix;
400d4a972cbSStefano Zampini   ((PetscObject)A)->state     = state + 1;
4014222ddf1SHong Zhang   A->product                  = product;
40233e6eea4SJose E. Roman   A->redundant                = redundant;
403273d9f13SBarry Smith 
4045c9eb25fSBarry Smith   /* since these two are copied into A we do not want them destroyed in C */
405f4259b30SLisandro Dalcin   ((PetscObject)*C)->qlist = NULL;
406f4259b30SLisandro Dalcin   ((PetscObject)*C)->olist = NULL;
40726fbe8dcSKarl Rupp 
4089566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(C));
409273d9f13SBarry Smith   PetscFunctionReturn(0);
410273d9f13SBarry Smith }
4118ab5b326SKris Buschelman /*
412eb6b5d47SBarry Smith         Replace A's header with that of C; the C object is then destroyed
413d0f46423SBarry Smith 
414eb6b5d47SBarry Smith         This is essentially code moved from MatDestroy()
415eb6b5d47SBarry Smith 
416eb6b5d47SBarry Smith         This is somewhat different from MatHeaderMerge() it would be nice to merge the code
417b30237c6SBarry Smith 
418b30237c6SBarry Smith         Used in DM hence is declared PETSC_EXTERN
4198ab5b326SKris Buschelman */
42028be2f97SBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat A,Mat *C)
4218ab5b326SKris Buschelman {
42227b31e29SJed Brown   PetscInt         refct;
423fefd9316SJose E. Roman   PetscObjectState state;
42428be2f97SBarry Smith   struct _p_Mat    buffer;
42581fa06acSBarry Smith   MatStencilInfo   stencil;
4268ab5b326SKris Buschelman 
4278ab5b326SKris Buschelman   PetscFunctionBegin;
42827b31e29SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
42928be2f97SBarry Smith   PetscValidHeaderSpecific(*C,MAT_CLASSID,2);
43028be2f97SBarry Smith   if (A == *C) PetscFunctionReturn(0);
43128be2f97SBarry Smith   PetscCheckSameComm(A,1,*C,2);
432aed4548fSBarry 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);
4336d7c1e57SBarry Smith 
43428be2f97SBarry Smith   /* swap C and A */
43527b31e29SJed Brown   refct   = ((PetscObject)A)->refct;
436fefd9316SJose E. Roman   state   = ((PetscObject)A)->state;
43781fa06acSBarry Smith   stencil = A->stencil;
4389566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(&buffer,A,sizeof(struct _p_Mat)));
4399566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(A,*C,sizeof(struct _p_Mat)));
4409566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(*C,&buffer,sizeof(struct _p_Mat)));
44127b31e29SJed Brown   ((PetscObject)A)->refct = refct;
442fefd9316SJose E. Roman   ((PetscObject)A)->state = state + 1;
44381fa06acSBarry Smith   A->stencil              = stencil;
44426fbe8dcSKarl Rupp 
445c32d4117SBarry Smith   ((PetscObject)*C)->refct = 1;
4469566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*C,MATOP_DESTROY,(void(*)(void))NULL));
4479566063dSJacob Faibussowitsch   PetscCall(MatDestroy(C));
4488ab5b326SKris Buschelman   PetscFunctionReturn(0);
4498ab5b326SKris Buschelman }
450e7e92044SBarry Smith 
451e7e92044SBarry Smith /*@
452b470e4b4SRichard Tran Mills      MatBindToCPU - marks a matrix to temporarily stay on the CPU and perform computations on the CPU
453e7e92044SBarry Smith 
4542216c58aSStefano Zampini    Logically collective on Mat
4552216c58aSStefano Zampini 
456e7e92044SBarry Smith    Input Parameters:
457e7e92044SBarry Smith +   A - the matrix
458b470e4b4SRichard Tran Mills -   flg - bind to the CPU if value of PETSC_TRUE
459e7e92044SBarry Smith 
46090ea27d8SSatish Balay    Level: intermediate
4612216c58aSStefano Zampini 
462db781477SPatrick Sanan .seealso: `MatBoundToCPU()`
463e7e92044SBarry Smith @*/
464b470e4b4SRichard Tran Mills PetscErrorCode MatBindToCPU(Mat A,PetscBool flg)
465e7e92044SBarry Smith {
4667d871021SStefano Zampini   PetscFunctionBegin;
4672ffa8ee7SStefano Zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
4682ffa8ee7SStefano Zampini   PetscValidLogicalCollectiveBool(A,flg,2);
4692216c58aSStefano Zampini #if defined(PETSC_HAVE_DEVICE)
470b470e4b4SRichard Tran Mills   if (A->boundtocpu == flg) PetscFunctionReturn(0);
471b470e4b4SRichard Tran Mills   A->boundtocpu = flg;
4721baa6e33SBarry Smith   if (A->ops->bindtocpu) PetscCall((*A->ops->bindtocpu)(A,flg));
4732216c58aSStefano Zampini #endif
474e7e92044SBarry Smith   PetscFunctionReturn(0);
4752216c58aSStefano Zampini }
4762216c58aSStefano Zampini 
4772216c58aSStefano Zampini /*@
4782216c58aSStefano Zampini      MatBoundToCPU - query if a matrix is bound to the CPU
4792216c58aSStefano Zampini 
4802216c58aSStefano Zampini    Input Parameter:
4812216c58aSStefano Zampini .   A - the matrix
4822216c58aSStefano Zampini 
4832216c58aSStefano Zampini    Output Parameter:
4842216c58aSStefano Zampini .   flg - the logical flag
4852216c58aSStefano Zampini 
4862216c58aSStefano Zampini    Level: intermediate
4872216c58aSStefano Zampini 
488db781477SPatrick Sanan .seealso: `MatBindToCPU()`
4892216c58aSStefano Zampini @*/
4902216c58aSStefano Zampini PetscErrorCode MatBoundToCPU(Mat A,PetscBool *flg)
4912216c58aSStefano Zampini {
4922ffa8ee7SStefano Zampini   PetscFunctionBegin;
4932ffa8ee7SStefano Zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
494dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(flg,2);
4952216c58aSStefano Zampini #if defined(PETSC_HAVE_DEVICE)
4962216c58aSStefano Zampini   *flg = A->boundtocpu;
4972216c58aSStefano Zampini #else
4982216c58aSStefano Zampini   *flg = PETSC_TRUE;
4997d871021SStefano Zampini #endif
5002216c58aSStefano Zampini   PetscFunctionReturn(0);
501e7e92044SBarry Smith }
5027e8381f9SStefano Zampini 
5037e8381f9SStefano Zampini PetscErrorCode MatSetValuesCOO_Basic(Mat A,const PetscScalar coo_v[],InsertMode imode)
5047e8381f9SStefano Zampini {
5057e8381f9SStefano Zampini   IS             is_coo_i,is_coo_j;
5067e8381f9SStefano Zampini   const PetscInt *coo_i,*coo_j;
5077e8381f9SStefano Zampini   PetscInt       n,n_i,n_j;
5087e8381f9SStefano Zampini   PetscScalar    zero = 0.;
5097e8381f9SStefano Zampini 
5107e8381f9SStefano Zampini   PetscFunctionBegin;
5119566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)A,"__PETSc_coo_i",(PetscObject*)&is_coo_i));
5129566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)A,"__PETSc_coo_j",(PetscObject*)&is_coo_j));
51328b400f6SJacob Faibussowitsch   PetscCheck(is_coo_i,PetscObjectComm((PetscObject)A),PETSC_ERR_COR,"Missing coo_i IS");
51428b400f6SJacob Faibussowitsch   PetscCheck(is_coo_j,PetscObjectComm((PetscObject)A),PETSC_ERR_COR,"Missing coo_j IS");
5159566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is_coo_i,&n_i));
5169566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is_coo_j,&n_j));
51708401ef6SPierre Jolivet   PetscCheck(n_i == n_j,PETSC_COMM_SELF,PETSC_ERR_COR,"Wrong local size %" PetscInt_FMT " != %" PetscInt_FMT,n_i,n_j);
5189566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is_coo_i,&coo_i));
5199566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is_coo_j,&coo_j));
520e61fc153SStefano Zampini   if (imode != ADD_VALUES) {
5219566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(A));
522e61fc153SStefano Zampini   }
5237e8381f9SStefano Zampini   for (n = 0; n < n_i; n++) {
5249566063dSJacob Faibussowitsch     PetscCall(MatSetValue(A,coo_i[n],coo_j[n],coo_v ? coo_v[n] : zero,ADD_VALUES));
5257e8381f9SStefano Zampini   }
5269566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is_coo_i,&coo_i));
5279566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is_coo_j,&coo_j));
5287e8381f9SStefano Zampini   PetscFunctionReturn(0);
5297e8381f9SStefano Zampini }
5307e8381f9SStefano Zampini 
53182a78a4eSJed Brown PetscErrorCode MatSetPreallocationCOO_Basic(Mat A,PetscCount ncoo,const PetscInt coo_i[],const PetscInt coo_j[])
5327e8381f9SStefano Zampini {
5337e8381f9SStefano Zampini   Mat            preallocator;
5347e8381f9SStefano Zampini   IS             is_coo_i,is_coo_j;
5357e8381f9SStefano Zampini   PetscScalar    zero = 0.0;
5367e8381f9SStefano Zampini 
5377e8381f9SStefano Zampini   PetscFunctionBegin;
5389566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->rmap));
5399566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->cmap));
5409566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&preallocator));
5419566063dSJacob Faibussowitsch   PetscCall(MatSetType(preallocator,MATPREALLOCATOR));
5429566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(preallocator,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N));
5439566063dSJacob Faibussowitsch   PetscCall(MatSetLayouts(preallocator,A->rmap,A->cmap));
5449566063dSJacob Faibussowitsch   PetscCall(MatSetUp(preallocator));
54582a78a4eSJed Brown   for (PetscCount n = 0; n < ncoo; n++) {
5469566063dSJacob Faibussowitsch     PetscCall(MatSetValue(preallocator,coo_i[n],coo_j[n],zero,INSERT_VALUES));
5477e8381f9SStefano Zampini   }
5489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY));
5499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY));
5509566063dSJacob Faibussowitsch   PetscCall(MatPreallocatorPreallocate(preallocator,PETSC_TRUE,A));
5519566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&preallocator));
5522c71b3e2SJacob Faibussowitsch   PetscCheck(ncoo <= PETSC_MAX_INT,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"ncoo %" PetscCount_FMT " overflowed PetscInt; configure --with-64-bit-indices or request support",ncoo);
5539566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,ncoo,coo_i,PETSC_COPY_VALUES,&is_coo_i));
5549566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,ncoo,coo_j,PETSC_COPY_VALUES,&is_coo_j));
5559566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)A,"__PETSc_coo_i",(PetscObject)is_coo_i));
5569566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)A,"__PETSc_coo_j",(PetscObject)is_coo_j));
5579566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_coo_i));
5589566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_coo_j));
5597e8381f9SStefano Zampini   PetscFunctionReturn(0);
5607e8381f9SStefano Zampini }
5617e8381f9SStefano Zampini 
56256856777SBarry Smith /*@C
563c3dd2894SJed Brown    MatSetPreallocationCOO - set preallocation for matrices using a coordinate format of the entries with global indices
5647e8381f9SStefano Zampini 
5657e8381f9SStefano Zampini    Collective on Mat
5667e8381f9SStefano Zampini 
5674165533cSJose E. Roman    Input Parameters:
5687e8381f9SStefano Zampini +  A - matrix being preallocated
56942550becSJunchao Zhang .  ncoo - number of entries
5707e8381f9SStefano Zampini .  coo_i - row indices
5717e8381f9SStefano Zampini -  coo_j - column indices
5727e8381f9SStefano Zampini 
5737e8381f9SStefano Zampini    Level: beginner
5747e8381f9SStefano Zampini 
575394ed5ebSJunchao Zhang    Notes:
576394ed5ebSJunchao Zhang    Entries can be repeated, see MatSetValuesCOO(). Entries with negative row or column indices are allowed
577394ed5ebSJunchao Zhang    but will be ignored. The corresponding entries in MatSetValuesCOO() will be ignored too. Remote entries
578394ed5ebSJunchao Zhang    are allowed and will be properly added or inserted to the matrix, unless the matrix option MAT_IGNORE_OFF_PROC_ENTRIES
579394ed5ebSJunchao Zhang    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.
5807e8381f9SStefano Zampini 
581735d7f90SBarry Smith    The arrays coo_i and coo_j may be freed immediately after calling this function.
582735d7f90SBarry Smith 
583db781477SPatrick Sanan .seealso: `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOOLocal()`, `DMSetMatrixPreallocateSkip()`
5847e8381f9SStefano Zampini @*/
58582a78a4eSJed Brown PetscErrorCode MatSetPreallocationCOO(Mat A,PetscCount ncoo,const PetscInt coo_i[],const PetscInt coo_j[])
5867e8381f9SStefano Zampini {
58782a78a4eSJed Brown   PetscErrorCode (*f)(Mat,PetscCount,const PetscInt[],const PetscInt[]) = NULL;
5887e8381f9SStefano Zampini 
5897e8381f9SStefano Zampini   PetscFunctionBegin;
5907e8381f9SStefano Zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
5917e8381f9SStefano Zampini   PetscValidType(A,1);
5927e8381f9SStefano Zampini   if (ncoo) PetscValidIntPointer(coo_i,3);
5937e8381f9SStefano Zampini   if (ncoo) PetscValidIntPointer(coo_j,4);
5949566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->rmap));
5959566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->cmap));
5969566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)A,"MatSetPreallocationCOO_C",&f));
597cbc6b225SStefano Zampini 
5989566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_PreallCOO,A,0,0,0));
5997e8381f9SStefano Zampini   if (f) {
6009566063dSJacob Faibussowitsch     PetscCall((*f)(A,ncoo,coo_i,coo_j));
6017e8381f9SStefano Zampini   } else { /* allow fallback, very slow */
6029566063dSJacob Faibussowitsch     PetscCall(MatSetPreallocationCOO_Basic(A,ncoo,coo_i,coo_j));
6037e8381f9SStefano Zampini   }
6049566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_PreallCOO,A,0,0,0));
6056834774dSStefano Zampini   A->preallocated = PETSC_TRUE;
606cbc6b225SStefano Zampini   A->nonzerostate++;
6077e8381f9SStefano Zampini   PetscFunctionReturn(0);
6087e8381f9SStefano Zampini }
6097e8381f9SStefano Zampini 
61056856777SBarry Smith /*@C
611c3dd2894SJed Brown    MatSetPreallocationCOOLocal - set preallocation for matrices using a coordinate format of the entries with local indices
612c3dd2894SJed Brown 
613c3dd2894SJed Brown    Collective on Mat
614c3dd2894SJed Brown 
615c3dd2894SJed Brown    Input Parameters:
616c3dd2894SJed Brown +  A - matrix being preallocated
617c3dd2894SJed Brown .  ncoo - number of entries
618c3dd2894SJed Brown .  coo_i - row indices (local numbering; may be modified)
619c3dd2894SJed Brown -  coo_j - column indices (local numbering; may be modified)
620c3dd2894SJed Brown 
621c3dd2894SJed Brown    Level: beginner
622c3dd2894SJed Brown 
623c3dd2894SJed Brown    Notes:
624c3dd2894SJed Brown    The local indices are translated using the local to global mapping, thus MatSetLocalToGlobalMapping() must have been
625c3dd2894SJed Brown    called prior to this function.
626c3dd2894SJed Brown 
627c3dd2894SJed Brown    The indices coo_i and coo_j may be modified within this function. They might be translated to corresponding global
628735d7f90SBarry Smith    indices, but the caller should not rely on them having any specific value after this function returns. The arrays
629735d7f90SBarry Smith    can be freed or reused immediately after this function returns.
630c3dd2894SJed Brown 
631394ed5ebSJunchao Zhang    Entries can be repeated, see MatSetValuesCOO(). Entries with negative row or column indices are allowed
632394ed5ebSJunchao Zhang    but will be ignored. The corresponding entries in MatSetValuesCOO() will be ignored too. Remote entries
633394ed5ebSJunchao Zhang    are allowed and will be properly added or inserted to the matrix.
634c3dd2894SJed Brown 
635db781477SPatrick Sanan .seealso: `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOO()`, `DMSetMatrixPreallocateSkip()`
636c3dd2894SJed Brown @*/
63782a78a4eSJed Brown PetscErrorCode MatSetPreallocationCOOLocal(Mat A,PetscCount ncoo,PetscInt coo_i[],PetscInt coo_j[])
638c3dd2894SJed Brown {
6396834774dSStefano Zampini   PetscErrorCode (*f)(Mat,PetscCount,PetscInt[],PetscInt[]) = NULL;
640c3dd2894SJed Brown 
641c3dd2894SJed Brown   PetscFunctionBegin;
642c3dd2894SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
643c3dd2894SJed Brown   PetscValidType(A,1);
644c3dd2894SJed Brown   if (ncoo) PetscValidIntPointer(coo_i,3);
645c3dd2894SJed Brown   if (ncoo) PetscValidIntPointer(coo_j,4);
6466834774dSStefano Zampini   PetscCheck(ncoo <= PETSC_MAX_INT,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"ncoo %" PetscCount_FMT " overflowed PetscInt; configure --with-64-bit-indices or request support",ncoo);
6479566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->rmap));
6489566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->cmap));
649cbc6b225SStefano Zampini 
6509566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)A,"MatSetPreallocationCOOLocal_C",&f));
6516834774dSStefano Zampini   if (f) {
6529566063dSJacob Faibussowitsch     PetscCall((*f)(A,ncoo,coo_i,coo_j));
653cbc6b225SStefano Zampini     A->nonzerostate++;
6546834774dSStefano Zampini   } else {
655cbc6b225SStefano Zampini     ISLocalToGlobalMapping ltog_row,ltog_col;
6569566063dSJacob Faibussowitsch     PetscCall(MatGetLocalToGlobalMapping(A,&ltog_row,&ltog_col));
6579566063dSJacob Faibussowitsch     if (ltog_row) PetscCall(ISLocalToGlobalMappingApply(ltog_row,ncoo,coo_i,coo_i));
6589566063dSJacob Faibussowitsch     if (ltog_col) PetscCall(ISLocalToGlobalMappingApply(ltog_col,ncoo,coo_j,coo_j));
6599566063dSJacob Faibussowitsch     PetscCall(MatSetPreallocationCOO(A,ncoo,coo_i,coo_j));
6606834774dSStefano Zampini   }
6616834774dSStefano Zampini   A->preallocated = PETSC_TRUE;
662c3dd2894SJed Brown   PetscFunctionReturn(0);
663c3dd2894SJed Brown }
664c3dd2894SJed Brown 
665c3dd2894SJed Brown /*@
666bfcc3627SStefano Zampini    MatSetValuesCOO - set values at once in a matrix preallocated using MatSetPreallocationCOO()
6677e8381f9SStefano Zampini 
6687e8381f9SStefano Zampini    Collective on Mat
6697e8381f9SStefano Zampini 
6704165533cSJose E. Roman    Input Parameters:
6717e8381f9SStefano Zampini +  A - matrix being preallocated
672bfcc3627SStefano Zampini .  coo_v - the matrix values (can be NULL)
6737e8381f9SStefano Zampini -  imode - the insert mode
6747e8381f9SStefano Zampini 
6757e8381f9SStefano Zampini    Level: beginner
6767e8381f9SStefano Zampini 
677c3dd2894SJed Brown    Notes: The values must follow the order of the indices prescribed with MatSetPreallocationCOO() or MatSetPreallocationCOOLocal().
678735d7f90SBarry Smith           When repeated entries are specified in the COO indices the coo_v values are first properly summed, regardless of the value of imode.
679bfcc3627SStefano Zampini           The imode flag indicates if coo_v must be added to the current values of the matrix (ADD_VALUES) or overwritten (INSERT_VALUES).
680735d7f90SBarry Smith           MatAssemblyBegin() and MatAssemblyEnd() do not need to be called after this routine. It automatically handles the assembly process.
6817e8381f9SStefano Zampini 
682db781477SPatrick Sanan .seealso: `MatSetPreallocationCOO()`, `MatSetPreallocationCOOLocal()`, `InsertMode`, `INSERT_VALUES`, `ADD_VALUES`
6837e8381f9SStefano Zampini @*/
6847e8381f9SStefano Zampini PetscErrorCode MatSetValuesCOO(Mat A, const PetscScalar coo_v[], InsertMode imode)
6857e8381f9SStefano Zampini {
6867e8381f9SStefano Zampini   PetscErrorCode (*f)(Mat,const PetscScalar[],InsertMode) = NULL;
6877e8381f9SStefano Zampini 
6887e8381f9SStefano Zampini   PetscFunctionBegin;
6897e8381f9SStefano Zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
6907e8381f9SStefano Zampini   PetscValidType(A,1);
6917e8381f9SStefano Zampini   MatCheckPreallocated(A,1);
692bfcc3627SStefano Zampini   PetscValidLogicalCollectiveEnum(A,imode,3);
6939566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)A,"MatSetValuesCOO_C",&f));
6949566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_SetVCOO,A,0,0,0));
6957e8381f9SStefano Zampini   if (f) {
6969566063dSJacob Faibussowitsch     PetscCall((*f)(A,coo_v,imode));
6977e8381f9SStefano Zampini   } else { /* allow fallback */
6989566063dSJacob Faibussowitsch     PetscCall(MatSetValuesCOO_Basic(A,coo_v,imode));
6997e8381f9SStefano Zampini   }
7009566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_SetVCOO,A,0,0,0));
7019566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
7029566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
7037e8381f9SStefano Zampini   PetscFunctionReturn(0);
7047e8381f9SStefano Zampini }
70565a9ecf2SRichard Tran Mills 
70665a9ecf2SRichard Tran Mills /*@
70765a9ecf2SRichard 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
70865a9ecf2SRichard Tran Mills 
70965a9ecf2SRichard Tran Mills    Input Parameters:
71065a9ecf2SRichard Tran Mills +  A - the matrix
71165a9ecf2SRichard Tran Mills -  flg - flag indicating whether the boundtocpu flag should be propagated
71265a9ecf2SRichard Tran Mills 
71365a9ecf2SRichard Tran Mills    Level: developer
71465a9ecf2SRichard Tran Mills 
71565a9ecf2SRichard Tran Mills    Notes:
71665a9ecf2SRichard Tran Mills    If the value of flg is set to true, the following will occur:
71765a9ecf2SRichard Tran Mills 
71865a9ecf2SRichard Tran Mills    MatCreateSubMatrices() and MatCreateRedundantMatrix() will bind created matrices to CPU if the input matrix is bound to the CPU.
71965a9ecf2SRichard Tran Mills    MatCreateVecs() will bind created vectors to CPU if the input matrix is bound to the CPU.
72065a9ecf2SRichard Tran Mills    The bindingpropagates flag itself is also propagated by the above routines.
72165a9ecf2SRichard Tran Mills 
72265a9ecf2SRichard Tran Mills    Developer Notes:
72365a9ecf2SRichard Tran Mills    If the fine-scale DMDA has the -dm_bind_below option set to true, then DMCreateInterpolationScale() calls MatSetBindingPropagates()
72465a9ecf2SRichard Tran Mills    on the restriction/interpolation operator to set the bindingpropagates flag to true.
72565a9ecf2SRichard Tran Mills 
726db781477SPatrick Sanan .seealso: `VecSetBindingPropagates()`, `MatGetBindingPropagates()`
72765a9ecf2SRichard Tran Mills @*/
72865a9ecf2SRichard Tran Mills PetscErrorCode MatSetBindingPropagates(Mat A,PetscBool flg)
72965a9ecf2SRichard Tran Mills {
73065a9ecf2SRichard Tran Mills   PetscFunctionBegin;
73165a9ecf2SRichard Tran Mills   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
73265a9ecf2SRichard Tran Mills #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
73365a9ecf2SRichard Tran Mills   A->bindingpropagates = flg;
73465a9ecf2SRichard Tran Mills #endif
73565a9ecf2SRichard Tran Mills   PetscFunctionReturn(0);
73665a9ecf2SRichard Tran Mills }
737e9c74fd6SRichard Tran Mills 
738e9c74fd6SRichard Tran Mills /*@
739e9c74fd6SRichard 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
740e9c74fd6SRichard Tran Mills 
741e9c74fd6SRichard Tran Mills    Input Parameter:
742e9c74fd6SRichard Tran Mills .  A - the matrix
743e9c74fd6SRichard Tran Mills 
744e9c74fd6SRichard Tran Mills    Output Parameter:
745e9c74fd6SRichard Tran Mills .  flg - flag indicating whether the boundtocpu flag will be propagated
746e9c74fd6SRichard Tran Mills 
747e9c74fd6SRichard Tran Mills    Level: developer
748e9c74fd6SRichard Tran Mills 
749db781477SPatrick Sanan .seealso: `MatSetBindingPropagates()`
750e9c74fd6SRichard Tran Mills @*/
751e9c74fd6SRichard Tran Mills PetscErrorCode MatGetBindingPropagates(Mat A,PetscBool *flg)
752e9c74fd6SRichard Tran Mills {
753e9c74fd6SRichard Tran Mills   PetscFunctionBegin;
754e9c74fd6SRichard Tran Mills   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
755e9c74fd6SRichard Tran Mills   PetscValidBoolPointer(flg,2);
756e9c74fd6SRichard Tran Mills #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
757e9c74fd6SRichard Tran Mills   *flg = A->bindingpropagates;
758e9c74fd6SRichard Tran Mills #else
759e9c74fd6SRichard Tran Mills   *flg = PETSC_FALSE;
760e9c74fd6SRichard Tran Mills #endif
761e9c74fd6SRichard Tran Mills   PetscFunctionReturn(0);
762e9c74fd6SRichard Tran Mills }
763