xref: /petsc/src/mat/utils/gcreate.c (revision 91e63d38360eb9bc922f79d792328cc4769c01ac)
1 #include <petsc/private/matimpl.h>       /*I "petscmat.h"  I*/
2 
3 PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat mat,PetscInt rbs, PetscInt cbs)
4 {
5   PetscFunctionBegin;
6   if (!mat->preallocated) PetscFunctionReturn(0);
7   PetscCheckFalse(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);
8   PetscCheckFalse(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);
9   PetscFunctionReturn(0);
10 }
11 
12 PETSC_INTERN PetscErrorCode MatShift_Basic(Mat Y,PetscScalar a)
13 {
14   PetscErrorCode ierr;
15   PetscInt       i,start,end;
16   PetscScalar    alpha = a;
17   PetscBool      prevoption;
18 
19   PetscFunctionBegin;
20   ierr = MatGetOption(Y,MAT_NO_OFF_PROC_ENTRIES,&prevoption);CHKERRQ(ierr);
21   ierr = MatSetOption(Y,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
22   ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr);
23   for (i=start; i<end; i++) {
24     if (i < Y->cmap->N) {
25       ierr = MatSetValues(Y,1,&i,1,&i,&alpha,ADD_VALUES);CHKERRQ(ierr);
26     }
27   }
28   ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
29   ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
30   ierr = MatSetOption(Y,MAT_NO_OFF_PROC_ENTRIES,prevoption);CHKERRQ(ierr);
31   PetscFunctionReturn(0);
32 }
33 
34 /*@
35    MatCreate - Creates a matrix where the type is determined
36    from either a call to MatSetType() or from the options database
37    with a call to MatSetFromOptions(). The default matrix type is
38    AIJ, using the routines MatCreateSeqAIJ() or MatCreateAIJ()
39    if you do not set a type in the options database. If you never
40    call MatSetType() or MatSetFromOptions() it will generate an
41    error when you try to use the matrix.
42 
43    Collective
44 
45    Input Parameter:
46 .  comm - MPI communicator
47 
48    Output Parameter:
49 .  A - the matrix
50 
51    Options Database Keys:
52 +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
53 .    -mat_type mpiaij   - AIJ type, uses MatCreateAIJ()
54 .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
55 .    -mat_type mpidense - dense type, uses MatCreateDense()
56 .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
57 -    -mat_type mpibaij  - block AIJ type, uses MatCreateBAIJ()
58 
59    Even More Options Database Keys:
60    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
61    for additional format-specific options.
62 
63    Level: beginner
64 
65 .seealso: MatCreateSeqAIJ(), MatCreateAIJ(),
66           MatCreateSeqDense(), MatCreateDense(),
67           MatCreateSeqBAIJ(), MatCreateBAIJ(),
68           MatCreateSeqSBAIJ(), MatCreateSBAIJ(),
69           MatConvert()
70 @*/
71 PetscErrorCode  MatCreate(MPI_Comm comm,Mat *A)
72 {
73   Mat            B;
74   PetscErrorCode ierr;
75 
76   PetscFunctionBegin;
77   PetscValidPointer(A,2);
78 
79   *A = NULL;
80   ierr = MatInitializePackage();CHKERRQ(ierr);
81 
82   ierr = PetscHeaderCreate(B,MAT_CLASSID,"Mat","Matrix","Mat",comm,MatDestroy,MatView);CHKERRQ(ierr);
83   ierr = PetscLayoutCreate(comm,&B->rmap);CHKERRQ(ierr);
84   ierr = PetscLayoutCreate(comm,&B->cmap);CHKERRQ(ierr);
85   ierr = PetscStrallocpy(VECSTANDARD,&B->defaultvectype);CHKERRQ(ierr);
86 
87   B->congruentlayouts = PETSC_DECIDE;
88   B->preallocated     = PETSC_FALSE;
89 #if defined(PETSC_HAVE_DEVICE)
90   B->boundtocpu       = PETSC_TRUE;
91 #endif
92   *A                  = B;
93   PetscFunctionReturn(0);
94 }
95 
96 /*@
97    MatSetErrorIfFailure - Causes Mat to generate an error, for example a zero pivot, is detected.
98 
99    Logically Collective on Mat
100 
101    Input Parameters:
102 +  mat -  matrix obtained from MatCreate()
103 -  flg - PETSC_TRUE indicates you want the error generated
104 
105    Level: advanced
106 
107 .seealso: PCSetErrorIfFailure()
108 @*/
109 PetscErrorCode  MatSetErrorIfFailure(Mat mat,PetscBool flg)
110 {
111   PetscFunctionBegin;
112   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
113   PetscValidLogicalCollectiveBool(mat,flg,2);
114   mat->erroriffailure = flg;
115   PetscFunctionReturn(0);
116 }
117 
118 /*@
119   MatSetSizes - Sets the local and global sizes, and checks to determine compatibility
120 
121   Collective on Mat
122 
123   Input Parameters:
124 +  A - the matrix
125 .  m - number of local rows (or PETSC_DECIDE)
126 .  n - number of local columns (or PETSC_DECIDE)
127 .  M - number of global rows (or PETSC_DETERMINE)
128 -  N - number of global columns (or PETSC_DETERMINE)
129 
130    Notes:
131    m (n) and M (N) cannot be both PETSC_DECIDE
132    If one processor calls this with M (N) of PETSC_DECIDE then all processors must, otherwise the program will hang.
133 
134    If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the
135    user must ensure that they are chosen to be compatible with the
136    vectors. To do this, one first considers the matrix-vector product
137    'y = A x'. The 'm' that is used in the above routine must match the
138    local size used in the vector creation routine VecCreateMPI() for 'y'.
139    Likewise, the 'n' used must match that used as the local size in
140    VecCreateMPI() for 'x'.
141 
142    You cannot change the sizes once they have been set.
143 
144    The sizes must be set before MatSetUp() or MatXXXSetPreallocation() is called.
145 
146   Level: beginner
147 
148 .seealso: MatGetSize(), PetscSplitOwnership()
149 @*/
150 PetscErrorCode  MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N)
151 {
152   PetscFunctionBegin;
153   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
154   PetscValidLogicalCollectiveInt(A,M,4);
155   PetscValidLogicalCollectiveInt(A,N,5);
156   PetscCheckFalse(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);
157   PetscCheckFalse(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);
158   PetscCheckFalse((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);
159   PetscCheckFalse((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);
160   A->rmap->n = m;
161   A->cmap->n = n;
162   A->rmap->N = M > -1 ? M : A->rmap->N;
163   A->cmap->N = N > -1 ? N : A->cmap->N;
164   PetscFunctionReturn(0);
165 }
166 
167 /*@
168    MatSetFromOptions - Creates a matrix where the type is determined
169    from the options database. Generates a parallel MPI matrix if the
170    communicator has more than one processor.  The default matrix type is
171    AIJ, using the routines MatCreateSeqAIJ() and MatCreateAIJ() if
172    you do not select a type in the options database.
173 
174    Collective on Mat
175 
176    Input Parameter:
177 .  A - the matrix
178 
179    Options Database Keys:
180 +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
181 .    -mat_type mpiaij   - AIJ type, uses MatCreateAIJ()
182 .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
183 .    -mat_type mpidense - dense type, uses MatCreateDense()
184 .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
185 -    -mat_type mpibaij  - block AIJ type, uses MatCreateBAIJ()
186 
187    Even More Options Database Keys:
188    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
189    for additional format-specific options.
190 
191    Level: beginner
192 
193 .seealso: MatCreateSeqAIJ((), MatCreateAIJ(),
194           MatCreateSeqDense(), MatCreateDense(),
195           MatCreateSeqBAIJ(), MatCreateBAIJ(),
196           MatCreateSeqSBAIJ(), MatCreateSBAIJ(),
197           MatConvert()
198 @*/
199 PetscErrorCode  MatSetFromOptions(Mat B)
200 {
201   PetscErrorCode ierr;
202   const char     *deft = MATAIJ;
203   char           type[256];
204   PetscBool      flg,set;
205   PetscInt       bind_below = 0;
206 
207   PetscFunctionBegin;
208   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
209 
210   ierr = PetscObjectOptionsBegin((PetscObject)B);CHKERRQ(ierr);
211 
212   if (B->rmap->bs < 0) {
213     PetscInt newbs = -1;
214     ierr = PetscOptionsInt("-mat_block_size","Set the blocksize used to store the matrix","MatSetBlockSize",newbs,&newbs,&flg);CHKERRQ(ierr);
215     if (flg) {
216       ierr = PetscLayoutSetBlockSize(B->rmap,newbs);CHKERRQ(ierr);
217       ierr = PetscLayoutSetBlockSize(B->cmap,newbs);CHKERRQ(ierr);
218     }
219   }
220 
221   ierr = PetscOptionsFList("-mat_type","Matrix type","MatSetType",MatList,deft,type,256,&flg);CHKERRQ(ierr);
222   if (flg) {
223     ierr = MatSetType(B,type);CHKERRQ(ierr);
224   } else if (!((PetscObject)B)->type_name) {
225     ierr = MatSetType(B,deft);CHKERRQ(ierr);
226   }
227 
228   ierr = PetscOptionsName("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",&B->checksymmetryonassembly);CHKERRQ(ierr);
229   ierr = PetscOptionsReal("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",B->checksymmetrytol,&B->checksymmetrytol,NULL);CHKERRQ(ierr);
230   ierr = PetscOptionsBool("-mat_null_space_test","Checks if provided null space is correct in MatAssemblyEnd()","MatSetNullSpaceTest",B->checknullspaceonassembly,&B->checknullspaceonassembly,NULL);CHKERRQ(ierr);
231   ierr = PetscOptionsBool("-mat_error_if_failure","Generate an error if an error occurs when factoring the matrix","MatSetErrorIfFailure",B->erroriffailure,&B->erroriffailure,NULL);CHKERRQ(ierr);
232 
233   if (B->ops->setfromoptions) {
234     ierr = (*B->ops->setfromoptions)(PetscOptionsObject,B);CHKERRQ(ierr);
235   }
236 
237   flg  = PETSC_FALSE;
238   ierr = 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);CHKERRQ(ierr);
239   if (set) {ierr = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,flg);CHKERRQ(ierr);}
240   flg  = PETSC_FALSE;
241   ierr = 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);CHKERRQ(ierr);
242   if (set) {ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,flg);CHKERRQ(ierr);}
243   flg  = PETSC_FALSE;
244   ierr = 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);CHKERRQ(ierr);
245   if (set) {ierr = MatSetOption(B,MAT_IGNORE_ZERO_ENTRIES,flg);CHKERRQ(ierr);}
246 
247   flg  = PETSC_FALSE;
248   ierr = PetscOptionsBool("-mat_form_explicit_transpose","Hint to form an explicit transpose for operations like MatMultTranspose","MatSetOption",flg,&flg,&set);CHKERRQ(ierr);
249   if (set) {ierr = MatSetOption(B,MAT_FORM_EXPLICIT_TRANSPOSE,flg);CHKERRQ(ierr);}
250 
251   /* Bind to CPU if below a user-specified size threshold.
252    * This perhaps belongs in the options for the GPU Mat types, but MatBindToCPU() does nothing when called on non-GPU types,
253    * and putting it here makes is more maintainable than duplicating this for all. */
254   ierr = 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);CHKERRQ(ierr);
255   if (flg && B->rmap->n < bind_below) {
256     ierr = MatBindToCPU(B,PETSC_TRUE);CHKERRQ(ierr);
257   }
258 
259   /* process any options handlers added with PetscObjectAddOptionsHandler() */
260   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)B);CHKERRQ(ierr);
261   ierr = PetscOptionsEnd();CHKERRQ(ierr);
262   PetscFunctionReturn(0);
263 }
264 
265 /*@C
266    MatXAIJSetPreallocation - set preallocation for serial and parallel AIJ, BAIJ, and SBAIJ matrices and their unassembled versions.
267 
268    Collective on Mat
269 
270    Input Parameters:
271 +  A - matrix being preallocated
272 .  bs - block size
273 .  dnnz - number of nonzero column blocks per block row of diagonal part of parallel matrix
274 .  onnz - number of nonzero column blocks per block row of off-diagonal part of parallel matrix
275 .  dnnzu - number of nonzero column blocks per block row of upper-triangular part of diagonal part of parallel matrix
276 -  onnzu - number of nonzero column blocks per block row of upper-triangular part of off-diagonal part of parallel matrix
277 
278    Level: beginner
279 
280 .seealso: MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation(),
281           PetscSplitOwnership()
282 @*/
283 PetscErrorCode MatXAIJSetPreallocation(Mat A,PetscInt bs,const PetscInt dnnz[],const PetscInt onnz[],const PetscInt dnnzu[],const PetscInt onnzu[])
284 {
285   PetscErrorCode ierr;
286   PetscInt       cbs;
287   void           (*aij)(void);
288   void           (*is)(void);
289   void           (*hyp)(void) = NULL;
290 
291   PetscFunctionBegin;
292   if (bs != PETSC_DECIDE) { /* don't mess with an already set block size */
293     ierr = MatSetBlockSize(A,bs);CHKERRQ(ierr);
294   }
295   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
296   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
297   ierr = MatGetBlockSizes(A,&bs,&cbs);CHKERRQ(ierr);
298   /* these routines assumes bs == cbs, this should be checked somehow */
299   ierr = MatSeqBAIJSetPreallocation(A,bs,0,dnnz);CHKERRQ(ierr);
300   ierr = MatMPIBAIJSetPreallocation(A,bs,0,dnnz,0,onnz);CHKERRQ(ierr);
301   ierr = MatSeqSBAIJSetPreallocation(A,bs,0,dnnzu);CHKERRQ(ierr);
302   ierr = MatMPISBAIJSetPreallocation(A,bs,0,dnnzu,0,onnzu);CHKERRQ(ierr);
303   /*
304     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
305     good before going on with it.
306   */
307   ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
308   ierr = PetscObjectQueryFunction((PetscObject)A,"MatISSetPreallocation_C",&is);CHKERRQ(ierr);
309 #if defined(PETSC_HAVE_HYPRE)
310   ierr = PetscObjectQueryFunction((PetscObject)A,"MatHYPRESetPreallocation_C",&hyp);CHKERRQ(ierr);
311 #endif
312   if (!aij && !is && !hyp) {
313     ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
314   }
315   if (aij || is || hyp) {
316     if (bs == cbs && bs == 1) {
317       ierr = MatSeqAIJSetPreallocation(A,0,dnnz);CHKERRQ(ierr);
318       ierr = MatMPIAIJSetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr);
319       ierr = MatISSetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr);
320 #if defined(PETSC_HAVE_HYPRE)
321       ierr = MatHYPRESetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr);
322 #endif
323     } else { /* Convert block-row precallocation to scalar-row */
324       PetscInt i,m,*sdnnz,*sonnz;
325       ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
326       ierr = PetscMalloc2((!!dnnz)*m,&sdnnz,(!!onnz)*m,&sonnz);CHKERRQ(ierr);
327       for (i=0; i<m; i++) {
328         if (dnnz) sdnnz[i] = dnnz[i/bs] * cbs;
329         if (onnz) sonnz[i] = onnz[i/bs] * cbs;
330       }
331       ierr = MatSeqAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL);CHKERRQ(ierr);
332       ierr = MatMPIAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr);
333       ierr = MatISSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr);
334 #if defined(PETSC_HAVE_HYPRE)
335       ierr = MatHYPRESetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr);
336 #endif
337       ierr = PetscFree2(sdnnz,sonnz);CHKERRQ(ierr);
338     }
339   }
340   PetscFunctionReturn(0);
341 }
342 
343 /*
344         Merges some information from Cs header to A; the C object is then destroyed
345 
346         This is somewhat different from MatHeaderReplace() it would be nice to merge the code
347 */
348 PetscErrorCode MatHeaderMerge(Mat A,Mat *C)
349 {
350   PetscErrorCode   ierr;
351   PetscInt         refct;
352   PetscOps         Abops;
353   struct _MatOps   Aops;
354   char             *mtype,*mname,*mprefix;
355   Mat_Product      *product;
356   Mat_Redundant    *redundant;
357   PetscObjectState state;
358 
359   PetscFunctionBegin;
360   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
361   PetscValidHeaderSpecific(*C,MAT_CLASSID,2);
362   if (A == *C) PetscFunctionReturn(0);
363   PetscCheckSameComm(A,1,*C,2);
364   /* save the parts of A we need */
365   Abops     = ((PetscObject)A)->bops[0];
366   Aops      = A->ops[0];
367   refct     = ((PetscObject)A)->refct;
368   mtype     = ((PetscObject)A)->type_name;
369   mname     = ((PetscObject)A)->name;
370   state     = ((PetscObject)A)->state;
371   mprefix   = ((PetscObject)A)->prefix;
372   product   = A->product;
373   redundant = A->redundant;
374 
375   /* zero these so the destroy below does not free them */
376   ((PetscObject)A)->type_name = NULL;
377   ((PetscObject)A)->name      = NULL;
378 
379   /* free all the interior data structures from mat */
380   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
381 
382   ierr = PetscFree(A->defaultvectype);CHKERRQ(ierr);
383   ierr = PetscLayoutDestroy(&A->rmap);CHKERRQ(ierr);
384   ierr = PetscLayoutDestroy(&A->cmap);CHKERRQ(ierr);
385   ierr = PetscFunctionListDestroy(&((PetscObject)A)->qlist);CHKERRQ(ierr);
386   ierr = PetscObjectListDestroy(&((PetscObject)A)->olist);CHKERRQ(ierr);
387   ierr = PetscComposedQuantitiesDestroy((PetscObject)A);CHKERRQ(ierr);
388 
389   /* copy C over to A */
390   ierr = PetscMemcpy(A,*C,sizeof(struct _p_Mat));CHKERRQ(ierr);
391 
392   /* return the parts of A we saved */
393   ((PetscObject)A)->bops[0]   = Abops;
394   A->ops[0]                   = Aops;
395   ((PetscObject)A)->refct     = refct;
396   ((PetscObject)A)->type_name = mtype;
397   ((PetscObject)A)->name      = mname;
398   ((PetscObject)A)->prefix    = mprefix;
399   ((PetscObject)A)->state     = state + 1;
400   A->product                  = product;
401   A->redundant                = redundant;
402 
403   /* since these two are copied into A we do not want them destroyed in C */
404   ((PetscObject)*C)->qlist = NULL;
405   ((PetscObject)*C)->olist = NULL;
406 
407   ierr = PetscHeaderDestroy(C);CHKERRQ(ierr);
408   PetscFunctionReturn(0);
409 }
410 /*
411         Replace A's header with that of C; the C object is then destroyed
412 
413         This is essentially code moved from MatDestroy()
414 
415         This is somewhat different from MatHeaderMerge() it would be nice to merge the code
416 
417         Used in DM hence is declared PETSC_EXTERN
418 */
419 PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat A,Mat *C)
420 {
421   PetscErrorCode   ierr;
422   PetscInt         refct;
423   PetscObjectState state;
424   struct _p_Mat    buffer;
425   MatStencilInfo   stencil;
426 
427   PetscFunctionBegin;
428   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
429   PetscValidHeaderSpecific(*C,MAT_CLASSID,2);
430   if (A == *C) PetscFunctionReturn(0);
431   PetscCheckSameComm(A,1,*C,2);
432   PetscCheckFalse(((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);
433 
434   /* swap C and A */
435   refct   = ((PetscObject)A)->refct;
436   state   = ((PetscObject)A)->state;
437   stencil = A->stencil;
438   ierr  = PetscMemcpy(&buffer,A,sizeof(struct _p_Mat));CHKERRQ(ierr);
439   ierr  = PetscMemcpy(A,*C,sizeof(struct _p_Mat));CHKERRQ(ierr);
440   ierr  = PetscMemcpy(*C,&buffer,sizeof(struct _p_Mat));CHKERRQ(ierr);
441   ((PetscObject)A)->refct = refct;
442   ((PetscObject)A)->state = state + 1;
443   A->stencil              = stencil;
444 
445   ((PetscObject)*C)->refct = 1;
446   ierr = MatShellSetOperation(*C,MATOP_DESTROY,(void(*)(void))NULL);CHKERRQ(ierr);
447   ierr = MatDestroy(C);CHKERRQ(ierr);
448   PetscFunctionReturn(0);
449 }
450 
451 /*@
452      MatBindToCPU - marks a matrix to temporarily stay on the CPU and perform computations on the CPU
453 
454    Logically collective on Mat
455 
456    Input Parameters:
457 +   A - the matrix
458 -   flg - bind to the CPU if value of PETSC_TRUE
459 
460    Level: intermediate
461 
462 .seealso: MatBoundToCPU()
463 @*/
464 PetscErrorCode MatBindToCPU(Mat A,PetscBool flg)
465 {
466   PetscFunctionBegin;
467   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
468   PetscValidLogicalCollectiveBool(A,flg,2);
469 #if defined(PETSC_HAVE_DEVICE)
470   if (A->boundtocpu == flg) PetscFunctionReturn(0);
471   A->boundtocpu = flg;
472   if (A->ops->bindtocpu) {
473     PetscErrorCode ierr;
474     ierr = (*A->ops->bindtocpu)(A,flg);CHKERRQ(ierr);
475   }
476 #endif
477   PetscFunctionReturn(0);
478 }
479 
480 /*@
481      MatBoundToCPU - query if a matrix is bound to the CPU
482 
483    Input Parameter:
484 .   A - the matrix
485 
486    Output Parameter:
487 .   flg - the logical flag
488 
489    Level: intermediate
490 
491 .seealso: MatBindToCPU()
492 @*/
493 PetscErrorCode MatBoundToCPU(Mat A,PetscBool *flg)
494 {
495   PetscFunctionBegin;
496   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
497   PetscValidPointer(flg,2);
498 #if defined(PETSC_HAVE_DEVICE)
499   *flg = A->boundtocpu;
500 #else
501   *flg = PETSC_TRUE;
502 #endif
503   PetscFunctionReturn(0);
504 }
505 
506 PetscErrorCode MatSetValuesCOO_Basic(Mat A,const PetscScalar coo_v[],InsertMode imode)
507 {
508   IS             is_coo_i,is_coo_j;
509   const PetscInt *coo_i,*coo_j;
510   PetscInt       n,n_i,n_j;
511   PetscScalar    zero = 0.;
512   PetscErrorCode ierr;
513 
514   PetscFunctionBegin;
515   ierr = PetscObjectQuery((PetscObject)A,"__PETSc_coo_i",(PetscObject*)&is_coo_i);CHKERRQ(ierr);
516   ierr = PetscObjectQuery((PetscObject)A,"__PETSc_coo_j",(PetscObject*)&is_coo_j);CHKERRQ(ierr);
517   PetscCheckFalse(!is_coo_i,PetscObjectComm((PetscObject)A),PETSC_ERR_COR,"Missing coo_i IS");
518   PetscCheckFalse(!is_coo_j,PetscObjectComm((PetscObject)A),PETSC_ERR_COR,"Missing coo_j IS");
519   ierr = ISGetLocalSize(is_coo_i,&n_i);CHKERRQ(ierr);
520   ierr = ISGetLocalSize(is_coo_j,&n_j);CHKERRQ(ierr);
521   PetscCheckFalse(n_i != n_j,PETSC_COMM_SELF,PETSC_ERR_COR,"Wrong local size %" PetscInt_FMT " != %" PetscInt_FMT,n_i,n_j);
522   ierr = ISGetIndices(is_coo_i,&coo_i);CHKERRQ(ierr);
523   ierr = ISGetIndices(is_coo_j,&coo_j);CHKERRQ(ierr);
524   if (imode != ADD_VALUES) {
525     ierr = MatZeroEntries(A);CHKERRQ(ierr);
526   }
527   for (n = 0; n < n_i; n++) {
528     ierr = MatSetValue(A,coo_i[n],coo_j[n],coo_v ? coo_v[n] : zero,ADD_VALUES);CHKERRQ(ierr);
529   }
530   ierr = ISRestoreIndices(is_coo_i,&coo_i);CHKERRQ(ierr);
531   ierr = ISRestoreIndices(is_coo_j,&coo_j);CHKERRQ(ierr);
532   PetscFunctionReturn(0);
533 }
534 
535 PetscErrorCode MatSetPreallocationCOO_Basic(Mat A,PetscCount ncoo,const PetscInt coo_i[],const PetscInt coo_j[])
536 {
537   Mat            preallocator;
538   IS             is_coo_i,is_coo_j;
539   PetscScalar    zero = 0.0;
540   PetscErrorCode ierr;
541 
542   PetscFunctionBegin;
543   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
544   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
545   ierr = MatCreate(PetscObjectComm((PetscObject)A),&preallocator);CHKERRQ(ierr);
546   ierr = MatSetType(preallocator,MATPREALLOCATOR);CHKERRQ(ierr);
547   ierr = MatSetSizes(preallocator,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
548   ierr = MatSetLayouts(preallocator,A->rmap,A->cmap);CHKERRQ(ierr);
549   ierr = MatSetUp(preallocator);CHKERRQ(ierr);
550   for (PetscCount n = 0; n < ncoo; n++) {
551     ierr = MatSetValue(preallocator,coo_i[n],coo_j[n],zero,INSERT_VALUES);CHKERRQ(ierr);
552   }
553   ierr = MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
554   ierr = MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
555   ierr = MatPreallocatorPreallocate(preallocator,PETSC_TRUE,A);CHKERRQ(ierr);
556   ierr = MatDestroy(&preallocator);CHKERRQ(ierr);
557   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);
558   ierr = ISCreateGeneral(PETSC_COMM_SELF,ncoo,coo_i,PETSC_COPY_VALUES,&is_coo_i);CHKERRQ(ierr);
559   ierr = ISCreateGeneral(PETSC_COMM_SELF,ncoo,coo_j,PETSC_COPY_VALUES,&is_coo_j);CHKERRQ(ierr);
560   ierr = PetscObjectCompose((PetscObject)A,"__PETSc_coo_i",(PetscObject)is_coo_i);CHKERRQ(ierr);
561   ierr = PetscObjectCompose((PetscObject)A,"__PETSc_coo_j",(PetscObject)is_coo_j);CHKERRQ(ierr);
562   ierr = ISDestroy(&is_coo_i);CHKERRQ(ierr);
563   ierr = ISDestroy(&is_coo_j);CHKERRQ(ierr);
564   PetscFunctionReturn(0);
565 }
566 
567 /*@
568    MatSetPreallocationCOO - set preallocation for matrices using a coordinate format of the entries with global indices
569 
570    Collective on Mat
571 
572    Input Parameters:
573 +  A - matrix being preallocated
574 .  ncoo - number of entries
575 .  coo_i - row indices
576 -  coo_j - column indices
577 
578    Level: beginner
579 
580    Notes:
581    Entries can be repeated, see MatSetValuesCOO(). Entries with negative row or column indices are allowed
582    but will be ignored. The corresponding entries in MatSetValuesCOO() will be ignored too. Remote entries
583    are allowed and will be properly added or inserted to the matrix, unless the matrix option MAT_IGNORE_OFF_PROC_ENTRIES
584    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.
585 
586    The arrays coo_i and coo_j may be freed immediately after calling this function.
587 
588 .seealso: MatSetValuesCOO(), MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation(), MatSetPreallocationCOOLocal()
589 @*/
590 PetscErrorCode MatSetPreallocationCOO(Mat A,PetscCount ncoo,const PetscInt coo_i[],const PetscInt coo_j[])
591 {
592   PetscErrorCode (*f)(Mat,PetscCount,const PetscInt[],const PetscInt[]) = NULL;
593   PetscErrorCode ierr;
594 
595   PetscFunctionBegin;
596   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
597   PetscValidType(A,1);
598   if (ncoo) PetscValidIntPointer(coo_i,3);
599   if (ncoo) PetscValidIntPointer(coo_j,4);
600   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
601   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
602   ierr = PetscObjectQueryFunction((PetscObject)A,"MatSetPreallocationCOO_C",&f);CHKERRQ(ierr);
603 
604   ierr = PetscLogEventBegin(MAT_PreallCOO,A,0,0,0);CHKERRQ(ierr);
605   if (f) {
606     ierr = (*f)(A,ncoo,coo_i,coo_j);CHKERRQ(ierr);
607   } else { /* allow fallback, very slow */
608     ierr = MatSetPreallocationCOO_Basic(A,ncoo,coo_i,coo_j);CHKERRQ(ierr);
609   }
610   ierr = PetscLogEventEnd(MAT_PreallCOO,A,0,0,0);CHKERRQ(ierr);
611   A->preallocated = PETSC_TRUE;
612   A->nonzerostate++;
613   PetscFunctionReturn(0);
614 }
615 
616 /*@
617    MatSetPreallocationCOOLocal - set preallocation for matrices using a coordinate format of the entries with local indices
618 
619    Collective on Mat
620 
621    Input Parameters:
622 +  A - matrix being preallocated
623 .  ncoo - number of entries
624 .  coo_i - row indices (local numbering; may be modified)
625 -  coo_j - column indices (local numbering; may be modified)
626 
627    Level: beginner
628 
629    Notes:
630    The local indices are translated using the local to global mapping, thus MatSetLocalToGlobalMapping() must have been
631    called prior to this function.
632 
633    The indices coo_i and coo_j may be modified within this function. They might be translated to corresponding global
634    indices, but the caller should not rely on them having any specific value after this function returns. The arrays
635    can be freed or reused immediately after this function returns.
636 
637    Entries can be repeated, see MatSetValuesCOO(). Entries with negative row or column indices are allowed
638    but will be ignored. The corresponding entries in MatSetValuesCOO() will be ignored too. Remote entries
639    are allowed and will be properly added or inserted to the matrix.
640 
641 .seealso: MatSetValuesCOO(), MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation(), MatSetPreallocationCOO()
642 @*/
643 PetscErrorCode MatSetPreallocationCOOLocal(Mat A,PetscCount ncoo,PetscInt coo_i[],PetscInt coo_j[])
644 {
645   PetscErrorCode ierr;
646   PetscErrorCode (*f)(Mat,PetscCount,PetscInt[],PetscInt[]) = NULL;
647 
648   PetscFunctionBegin;
649   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
650   PetscValidType(A,1);
651   if (ncoo) PetscValidIntPointer(coo_i,3);
652   if (ncoo) PetscValidIntPointer(coo_j,4);
653   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);
654   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
655   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
656 
657   ierr = PetscObjectQueryFunction((PetscObject)A,"MatSetPreallocationCOOLocal_C",&f);CHKERRQ(ierr);
658   if (f) {
659     ierr = (*f)(A,ncoo,coo_i,coo_j);CHKERRQ(ierr);
660     A->nonzerostate++;
661   } else {
662     ISLocalToGlobalMapping ltog_row,ltog_col;
663     ierr = MatGetLocalToGlobalMapping(A,&ltog_row,&ltog_col);CHKERRQ(ierr);
664     if (ltog_row) { ierr = ISLocalToGlobalMappingApply(ltog_row,ncoo,coo_i,coo_i);CHKERRQ(ierr); }
665     if (ltog_col) { ierr = ISLocalToGlobalMappingApply(ltog_col,ncoo,coo_j,coo_j);CHKERRQ(ierr); }
666     ierr = MatSetPreallocationCOO(A,ncoo,coo_i,coo_j);CHKERRQ(ierr);
667   }
668   A->preallocated = PETSC_TRUE;
669   PetscFunctionReturn(0);
670 }
671 
672 /*@
673    MatSetValuesCOO - set values at once in a matrix preallocated using MatSetPreallocationCOO()
674 
675    Collective on Mat
676 
677    Input Parameters:
678 +  A - matrix being preallocated
679 .  coo_v - the matrix values (can be NULL)
680 -  imode - the insert mode
681 
682    Level: beginner
683 
684    Notes: The values must follow the order of the indices prescribed with MatSetPreallocationCOO() or MatSetPreallocationCOOLocal().
685           When repeated entries are specified in the COO indices the coo_v values are first properly summed, regardless of the value of imode.
686           The imode flag indicates if coo_v must be added to the current values of the matrix (ADD_VALUES) or overwritten (INSERT_VALUES).
687           MatAssemblyBegin() and MatAssemblyEnd() do not need to be called after this routine. It automatically handles the assembly process.
688 
689 .seealso: MatSetPreallocationCOO(), MatSetPreallocationCOOLocal(), InsertMode, INSERT_VALUES, ADD_VALUES
690 @*/
691 PetscErrorCode MatSetValuesCOO(Mat A, const PetscScalar coo_v[], InsertMode imode)
692 {
693   PetscErrorCode (*f)(Mat,const PetscScalar[],InsertMode) = NULL;
694   PetscErrorCode ierr;
695 
696   PetscFunctionBegin;
697   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
698   PetscValidType(A,1);
699   MatCheckPreallocated(A,1);
700   PetscValidLogicalCollectiveEnum(A,imode,3);
701   ierr = PetscObjectQueryFunction((PetscObject)A,"MatSetValuesCOO_C",&f);CHKERRQ(ierr);
702   ierr = PetscLogEventBegin(MAT_SetVCOO,A,0,0,0);CHKERRQ(ierr);
703   if (f) {
704     ierr = (*f)(A,coo_v,imode);CHKERRQ(ierr);
705   } else { /* allow fallback */
706     ierr = MatSetValuesCOO_Basic(A,coo_v,imode);CHKERRQ(ierr);
707   }
708   ierr = PetscLogEventEnd(MAT_SetVCOO,A,0,0,0);CHKERRQ(ierr);
709   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
710   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
711   PetscFunctionReturn(0);
712 }
713 
714 /*@
715    MatSetBindingPropagates - Sets whether the state of being bound to the CPU for a GPU matrix type propagates to child and some other associated objects
716 
717    Input Parameters:
718 +  A - the matrix
719 -  flg - flag indicating whether the boundtocpu flag should be propagated
720 
721    Level: developer
722 
723    Notes:
724    If the value of flg is set to true, the following will occur:
725 
726    MatCreateSubMatrices() and MatCreateRedundantMatrix() will bind created matrices to CPU if the input matrix is bound to the CPU.
727    MatCreateVecs() will bind created vectors to CPU if the input matrix is bound to the CPU.
728    The bindingpropagates flag itself is also propagated by the above routines.
729 
730    Developer Notes:
731    If the fine-scale DMDA has the -dm_bind_below option set to true, then DMCreateInterpolationScale() calls MatSetBindingPropagates()
732    on the restriction/interpolation operator to set the bindingpropagates flag to true.
733 
734 .seealso: VecSetBindingPropagates(), MatGetBindingPropagates()
735 @*/
736 PetscErrorCode MatSetBindingPropagates(Mat A,PetscBool flg)
737 {
738   PetscFunctionBegin;
739   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
740 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
741   A->bindingpropagates = flg;
742 #endif
743   PetscFunctionReturn(0);
744 }
745 
746 /*@
747    MatGetBindingPropagates - Gets whether the state of being bound to the CPU for a GPU matrix type propagates to child and some other associated objects
748 
749    Input Parameter:
750 .  A - the matrix
751 
752    Output Parameter:
753 .  flg - flag indicating whether the boundtocpu flag will be propagated
754 
755    Level: developer
756 
757 .seealso: MatSetBindingPropagates()
758 @*/
759 PetscErrorCode MatGetBindingPropagates(Mat A,PetscBool *flg)
760 {
761   PetscFunctionBegin;
762   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
763   PetscValidBoolPointer(flg,2);
764 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
765   *flg = A->bindingpropagates;
766 #else
767   *flg = PETSC_FALSE;
768 #endif
769   PetscFunctionReturn(0);
770 }
771