Actual source code: gcreate.c

  1: #include <petsc/private/matimpl.h>

  3: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat mat,PetscInt rbs, PetscInt cbs)
  4: {
  5:   if (!mat->preallocated) return 0;
  8:   return 0;
  9: }

 11: PETSC_INTERN PetscErrorCode MatShift_Basic(Mat Y,PetscScalar a)
 12: {
 13:   PetscInt       i,start,end;
 14:   PetscScalar    alpha = a;
 15:   PetscBool      prevoption;

 17:   MatGetOption(Y,MAT_NO_OFF_PROC_ENTRIES,&prevoption);
 18:   MatSetOption(Y,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);
 19:   MatGetOwnershipRange(Y,&start,&end);
 20:   for (i=start; i<end; i++) {
 21:     if (i < Y->cmap->N) {
 22:       MatSetValues(Y,1,&i,1,&i,&alpha,ADD_VALUES);
 23:     }
 24:   }
 25:   MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);
 26:   MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);
 27:   MatSetOption(Y,MAT_NO_OFF_PROC_ENTRIES,prevoption);
 28:   return 0;
 29: }

 31: /*@
 32:    MatCreate - Creates a matrix where the type is determined
 33:    from either a call to MatSetType() or from the options database
 34:    with a call to MatSetFromOptions(). The default matrix type is
 35:    AIJ, using the routines MatCreateSeqAIJ() or MatCreateAIJ()
 36:    if you do not set a type in the options database. If you never
 37:    call MatSetType() or MatSetFromOptions() it will generate an
 38:    error when you try to use the matrix.

 40:    Collective

 42:    Input Parameter:
 43: .  comm - MPI communicator

 45:    Output Parameter:
 46: .  A - the matrix

 48:    Options Database Keys:
 49: +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
 50: .    -mat_type mpiaij   - AIJ type, uses MatCreateAIJ()
 51: .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
 52: .    -mat_type mpidense - dense type, uses MatCreateDense()
 53: .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
 54: -    -mat_type mpibaij  - block AIJ type, uses MatCreateBAIJ()

 56:    Even More Options Database Keys:
 57:    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
 58:    for additional format-specific options.

 60:    Level: beginner

 62: .seealso: MatCreateSeqAIJ(), MatCreateAIJ(),
 63:           MatCreateSeqDense(), MatCreateDense(),
 64:           MatCreateSeqBAIJ(), MatCreateBAIJ(),
 65:           MatCreateSeqSBAIJ(), MatCreateSBAIJ(),
 66:           MatConvert()
 67: @*/
 68: PetscErrorCode  MatCreate(MPI_Comm comm,Mat *A)
 69: {
 70:   Mat            B;


 74:   *A = NULL;
 75:   MatInitializePackage();

 77:   PetscHeaderCreate(B,MAT_CLASSID,"Mat","Matrix","Mat",comm,MatDestroy,MatView);
 78:   PetscLayoutCreate(comm,&B->rmap);
 79:   PetscLayoutCreate(comm,&B->cmap);
 80:   PetscStrallocpy(VECSTANDARD,&B->defaultvectype);

 82:   B->congruentlayouts = PETSC_DECIDE;
 83:   B->preallocated     = PETSC_FALSE;
 84: #if defined(PETSC_HAVE_DEVICE)
 85:   B->boundtocpu       = PETSC_TRUE;
 86: #endif
 87:   *A                  = B;
 88:   return 0;
 89: }

 91: /*@
 92:    MatSetErrorIfFailure - Causes Mat to generate an error, for example a zero pivot, is detected.

 94:    Logically Collective on Mat

 96:    Input Parameters:
 97: +  mat -  matrix obtained from MatCreate()
 98: -  flg - PETSC_TRUE indicates you want the error generated

100:    Level: advanced

102: .seealso: PCSetErrorIfFailure()
103: @*/
104: PetscErrorCode  MatSetErrorIfFailure(Mat mat,PetscBool flg)
105: {
108:   mat->erroriffailure = flg;
109:   return 0;
110: }

112: /*@
113:   MatSetSizes - Sets the local and global sizes, and checks to determine compatibility

115:   Collective on Mat

117:   Input Parameters:
118: +  A - the matrix
119: .  m - number of local rows (or PETSC_DECIDE)
120: .  n - number of local columns (or PETSC_DECIDE)
121: .  M - number of global rows (or PETSC_DETERMINE)
122: -  N - number of global columns (or PETSC_DETERMINE)

124:    Notes:
125:    m (n) and M (N) cannot be both PETSC_DECIDE
126:    If one processor calls this with M (N) of PETSC_DECIDE then all processors must, otherwise the program will hang.

128:    If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the
129:    user must ensure that they are chosen to be compatible with the
130:    vectors. To do this, one first considers the matrix-vector product
131:    'y = A x'. The 'm' that is used in the above routine must match the
132:    local size used in the vector creation routine VecCreateMPI() for 'y'.
133:    Likewise, the 'n' used must match that used as the local size in
134:    VecCreateMPI() for 'x'.

136:    You cannot change the sizes once they have been set.

138:    The sizes must be set before MatSetUp() or MatXXXSetPreallocation() is called.

140:   Level: beginner

142: .seealso: MatGetSize(), PetscSplitOwnership()
143: @*/
144: PetscErrorCode  MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N)
145: {
153:   A->rmap->n = m;
154:   A->cmap->n = n;
155:   A->rmap->N = M > -1 ? M : A->rmap->N;
156:   A->cmap->N = N > -1 ? N : A->cmap->N;
157:   return 0;
158: }

160: /*@
161:    MatSetFromOptions - Creates a matrix where the type is determined
162:    from the options database. Generates a parallel MPI matrix if the
163:    communicator has more than one processor.  The default matrix type is
164:    AIJ, using the routines MatCreateSeqAIJ() and MatCreateAIJ() if
165:    you do not select a type in the options database.

167:    Collective on Mat

169:    Input Parameter:
170: .  A - the matrix

172:    Options Database Keys:
173: +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
174: .    -mat_type mpiaij   - AIJ type, uses MatCreateAIJ()
175: .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
176: .    -mat_type mpidense - dense type, uses MatCreateDense()
177: .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
178: -    -mat_type mpibaij  - block AIJ type, uses MatCreateBAIJ()

180:    Even More Options Database Keys:
181:    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
182:    for additional format-specific options.

184:    Level: beginner

186: .seealso: MatCreateSeqAIJ((), MatCreateAIJ(),
187:           MatCreateSeqDense(), MatCreateDense(),
188:           MatCreateSeqBAIJ(), MatCreateBAIJ(),
189:           MatCreateSeqSBAIJ(), MatCreateSBAIJ(),
190:           MatConvert()
191: @*/
192: PetscErrorCode  MatSetFromOptions(Mat B)
193: {
195:   const char     *deft = MATAIJ;
196:   char           type[256];
197:   PetscBool      flg,set;
198:   PetscInt       bind_below = 0;


202:   PetscObjectOptionsBegin((PetscObject)B);

204:   if (B->rmap->bs < 0) {
205:     PetscInt newbs = -1;
206:     PetscOptionsInt("-mat_block_size","Set the blocksize used to store the matrix","MatSetBlockSize",newbs,&newbs,&flg);
207:     if (flg) {
208:       PetscLayoutSetBlockSize(B->rmap,newbs);
209:       PetscLayoutSetBlockSize(B->cmap,newbs);
210:     }
211:   }

213:   PetscOptionsFList("-mat_type","Matrix type","MatSetType",MatList,deft,type,256,&flg);
214:   if (flg) {
215:     MatSetType(B,type);
216:   } else if (!((PetscObject)B)->type_name) {
217:     MatSetType(B,deft);
218:   }

220:   PetscOptionsName("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",&B->checksymmetryonassembly);
221:   PetscOptionsReal("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",B->checksymmetrytol,&B->checksymmetrytol,NULL);
222:   PetscOptionsBool("-mat_null_space_test","Checks if provided null space is correct in MatAssemblyEnd()","MatSetNullSpaceTest",B->checknullspaceonassembly,&B->checknullspaceonassembly,NULL);
223:   PetscOptionsBool("-mat_error_if_failure","Generate an error if an error occurs when factoring the matrix","MatSetErrorIfFailure",B->erroriffailure,&B->erroriffailure,NULL);

225:   if (B->ops->setfromoptions) {
226:     (*B->ops->setfromoptions)(PetscOptionsObject,B);
227:   }

229:   flg  = PETSC_FALSE;
230:   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);
231:   if (set) MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,flg);
232:   flg  = PETSC_FALSE;
233:   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);
234:   if (set) MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,flg);
235:   flg  = PETSC_FALSE;
236:   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);
237:   if (set) MatSetOption(B,MAT_IGNORE_ZERO_ENTRIES,flg);

239:   flg  = PETSC_FALSE;
240:   PetscOptionsBool("-mat_form_explicit_transpose","Hint to form an explicit transpose for operations like MatMultTranspose","MatSetOption",flg,&flg,&set);
241:   if (set) MatSetOption(B,MAT_FORM_EXPLICIT_TRANSPOSE,flg);

243:   /* Bind to CPU if below a user-specified size threshold.
244:    * This perhaps belongs in the options for the GPU Mat types, but MatBindToCPU() does nothing when called on non-GPU types,
245:    * and putting it here makes is more maintainable than duplicating this for all. */
246:   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);
247:   if (flg && B->rmap->n < bind_below) {
248:     MatBindToCPU(B,PETSC_TRUE);
249:   }

251:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
252:   PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)B);
253:   PetscOptionsEnd();
254:   return 0;
255: }

257: /*@C
258:    MatXAIJSetPreallocation - set preallocation for serial and parallel AIJ, BAIJ, and SBAIJ matrices and their unassembled versions.

260:    Collective on Mat

262:    Input Parameters:
263: +  A - matrix being preallocated
264: .  bs - block size
265: .  dnnz - number of nonzero column blocks per block row of diagonal part of parallel matrix
266: .  onnz - number of nonzero column blocks per block row of off-diagonal part of parallel matrix
267: .  dnnzu - number of nonzero column blocks per block row of upper-triangular part of diagonal part of parallel matrix
268: -  onnzu - number of nonzero column blocks per block row of upper-triangular part of off-diagonal part of parallel matrix

270:    Level: beginner

272: .seealso: MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation(),
273:           PetscSplitOwnership()
274: @*/
275: PetscErrorCode MatXAIJSetPreallocation(Mat A,PetscInt bs,const PetscInt dnnz[],const PetscInt onnz[],const PetscInt dnnzu[],const PetscInt onnzu[])
276: {
277:   PetscInt       cbs;
278:   void           (*aij)(void);
279:   void           (*is)(void);
280:   void           (*hyp)(void) = NULL;

282:   if (bs != PETSC_DECIDE) { /* don't mess with an already set block size */
283:     MatSetBlockSize(A,bs);
284:   }
285:   PetscLayoutSetUp(A->rmap);
286:   PetscLayoutSetUp(A->cmap);
287:   MatGetBlockSizes(A,&bs,&cbs);
288:   /* these routines assumes bs == cbs, this should be checked somehow */
289:   MatSeqBAIJSetPreallocation(A,bs,0,dnnz);
290:   MatMPIBAIJSetPreallocation(A,bs,0,dnnz,0,onnz);
291:   MatSeqSBAIJSetPreallocation(A,bs,0,dnnzu);
292:   MatMPISBAIJSetPreallocation(A,bs,0,dnnzu,0,onnzu);
293:   /*
294:     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
295:     good before going on with it.
296:   */
297:   PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);
298:   PetscObjectQueryFunction((PetscObject)A,"MatISSetPreallocation_C",&is);
299: #if defined(PETSC_HAVE_HYPRE)
300:   PetscObjectQueryFunction((PetscObject)A,"MatHYPRESetPreallocation_C",&hyp);
301: #endif
302:   if (!aij && !is && !hyp) {
303:     PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);
304:   }
305:   if (aij || is || hyp) {
306:     if (bs == cbs && bs == 1) {
307:       MatSeqAIJSetPreallocation(A,0,dnnz);
308:       MatMPIAIJSetPreallocation(A,0,dnnz,0,onnz);
309:       MatISSetPreallocation(A,0,dnnz,0,onnz);
310: #if defined(PETSC_HAVE_HYPRE)
311:       MatHYPRESetPreallocation(A,0,dnnz,0,onnz);
312: #endif
313:     } else { /* Convert block-row precallocation to scalar-row */
314:       PetscInt i,m,*sdnnz,*sonnz;
315:       MatGetLocalSize(A,&m,NULL);
316:       PetscMalloc2((!!dnnz)*m,&sdnnz,(!!onnz)*m,&sonnz);
317:       for (i=0; i<m; i++) {
318:         if (dnnz) sdnnz[i] = dnnz[i/bs] * cbs;
319:         if (onnz) sonnz[i] = onnz[i/bs] * cbs;
320:       }
321:       MatSeqAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL);
322:       MatMPIAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);
323:       MatISSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);
324: #if defined(PETSC_HAVE_HYPRE)
325:       MatHYPRESetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);
326: #endif
327:       PetscFree2(sdnnz,sonnz);
328:     }
329:   }
330:   return 0;
331: }

333: /*
334:         Merges some information from Cs header to A; the C object is then destroyed

336:         This is somewhat different from MatHeaderReplace() it would be nice to merge the code
337: */
338: PetscErrorCode MatHeaderMerge(Mat A,Mat *C)
339: {
340:   PetscInt         refct;
341:   PetscOps         Abops;
342:   struct _MatOps   Aops;
343:   char             *mtype,*mname,*mprefix;
344:   Mat_Product      *product;
345:   Mat_Redundant    *redundant;
346:   PetscObjectState state;

350:   if (A == *C) return 0;
352:   /* save the parts of A we need */
353:   Abops     = ((PetscObject)A)->bops[0];
354:   Aops      = A->ops[0];
355:   refct     = ((PetscObject)A)->refct;
356:   mtype     = ((PetscObject)A)->type_name;
357:   mname     = ((PetscObject)A)->name;
358:   state     = ((PetscObject)A)->state;
359:   mprefix   = ((PetscObject)A)->prefix;
360:   product   = A->product;
361:   redundant = A->redundant;

363:   /* zero these so the destroy below does not free them */
364:   ((PetscObject)A)->type_name = NULL;
365:   ((PetscObject)A)->name      = NULL;

367:   /* free all the interior data structures from mat */
368:   (*A->ops->destroy)(A);

370:   PetscFree(A->defaultvectype);
371:   PetscLayoutDestroy(&A->rmap);
372:   PetscLayoutDestroy(&A->cmap);
373:   PetscFunctionListDestroy(&((PetscObject)A)->qlist);
374:   PetscObjectListDestroy(&((PetscObject)A)->olist);
375:   PetscComposedQuantitiesDestroy((PetscObject)A);

377:   /* copy C over to A */
378:   PetscMemcpy(A,*C,sizeof(struct _p_Mat));

380:   /* return the parts of A we saved */
381:   ((PetscObject)A)->bops[0]   = Abops;
382:   A->ops[0]                   = Aops;
383:   ((PetscObject)A)->refct     = refct;
384:   ((PetscObject)A)->type_name = mtype;
385:   ((PetscObject)A)->name      = mname;
386:   ((PetscObject)A)->prefix    = mprefix;
387:   ((PetscObject)A)->state     = state + 1;
388:   A->product                  = product;
389:   A->redundant                = redundant;

391:   /* since these two are copied into A we do not want them destroyed in C */
392:   ((PetscObject)*C)->qlist = NULL;
393:   ((PetscObject)*C)->olist = NULL;

395:   PetscHeaderDestroy(C);
396:   return 0;
397: }
398: /*
399:         Replace A's header with that of C; the C object is then destroyed

401:         This is essentially code moved from MatDestroy()

403:         This is somewhat different from MatHeaderMerge() it would be nice to merge the code

405:         Used in DM hence is declared PETSC_EXTERN
406: */
407: PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat A,Mat *C)
408: {
409:   PetscInt         refct;
410:   PetscObjectState state;
411:   struct _p_Mat    buffer;
412:   MatStencilInfo   stencil;

416:   if (A == *C) return 0;

420:   /* swap C and A */
421:   refct   = ((PetscObject)A)->refct;
422:   state   = ((PetscObject)A)->state;
423:   stencil = A->stencil;
424:   PetscMemcpy(&buffer,A,sizeof(struct _p_Mat));
425:   PetscMemcpy(A,*C,sizeof(struct _p_Mat));
426:   PetscMemcpy(*C,&buffer,sizeof(struct _p_Mat));
427:   ((PetscObject)A)->refct = refct;
428:   ((PetscObject)A)->state = state + 1;
429:   A->stencil              = stencil;

431:   ((PetscObject)*C)->refct = 1;
432:   MatShellSetOperation(*C,MATOP_DESTROY,(void(*)(void))NULL);
433:   MatDestroy(C);
434:   return 0;
435: }

437: /*@
438:      MatBindToCPU - marks a matrix to temporarily stay on the CPU and perform computations on the CPU

440:    Logically collective on Mat

442:    Input Parameters:
443: +   A - the matrix
444: -   flg - bind to the CPU if value of PETSC_TRUE

446:    Level: intermediate

448: .seealso: MatBoundToCPU()
449: @*/
450: PetscErrorCode MatBindToCPU(Mat A,PetscBool flg)
451: {
454: #if defined(PETSC_HAVE_DEVICE)
455:   if (A->boundtocpu == flg) return 0;
456:   A->boundtocpu = flg;
457:   if (A->ops->bindtocpu) {
458:     (*A->ops->bindtocpu)(A,flg);
459:   }
460: #endif
461:   return 0;
462: }

464: /*@
465:      MatBoundToCPU - query if a matrix is bound to the CPU

467:    Input Parameter:
468: .   A - the matrix

470:    Output Parameter:
471: .   flg - the logical flag

473:    Level: intermediate

475: .seealso: MatBindToCPU()
476: @*/
477: PetscErrorCode MatBoundToCPU(Mat A,PetscBool *flg)
478: {
481: #if defined(PETSC_HAVE_DEVICE)
482:   *flg = A->boundtocpu;
483: #else
484:   *flg = PETSC_TRUE;
485: #endif
486:   return 0;
487: }

489: PetscErrorCode MatSetValuesCOO_Basic(Mat A,const PetscScalar coo_v[],InsertMode imode)
490: {
491:   IS             is_coo_i,is_coo_j;
492:   const PetscInt *coo_i,*coo_j;
493:   PetscInt       n,n_i,n_j;
494:   PetscScalar    zero = 0.;

496:   PetscObjectQuery((PetscObject)A,"__PETSc_coo_i",(PetscObject*)&is_coo_i);
497:   PetscObjectQuery((PetscObject)A,"__PETSc_coo_j",(PetscObject*)&is_coo_j);
500:   ISGetLocalSize(is_coo_i,&n_i);
501:   ISGetLocalSize(is_coo_j,&n_j);
503:   ISGetIndices(is_coo_i,&coo_i);
504:   ISGetIndices(is_coo_j,&coo_j);
505:   if (imode != ADD_VALUES) {
506:     MatZeroEntries(A);
507:   }
508:   for (n = 0; n < n_i; n++) {
509:     MatSetValue(A,coo_i[n],coo_j[n],coo_v ? coo_v[n] : zero,ADD_VALUES);
510:   }
511:   ISRestoreIndices(is_coo_i,&coo_i);
512:   ISRestoreIndices(is_coo_j,&coo_j);
513:   return 0;
514: }

516: PetscErrorCode MatSetPreallocationCOO_Basic(Mat A,PetscCount ncoo,const PetscInt coo_i[],const PetscInt coo_j[])
517: {
518:   Mat            preallocator;
519:   IS             is_coo_i,is_coo_j;
520:   PetscScalar    zero = 0.0;

522:   PetscLayoutSetUp(A->rmap);
523:   PetscLayoutSetUp(A->cmap);
524:   MatCreate(PetscObjectComm((PetscObject)A),&preallocator);
525:   MatSetType(preallocator,MATPREALLOCATOR);
526:   MatSetSizes(preallocator,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
527:   MatSetLayouts(preallocator,A->rmap,A->cmap);
528:   MatSetUp(preallocator);
529:   for (PetscCount n = 0; n < ncoo; n++) {
530:     MatSetValue(preallocator,coo_i[n],coo_j[n],zero,INSERT_VALUES);
531:   }
532:   MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY);
533:   MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY);
534:   MatPreallocatorPreallocate(preallocator,PETSC_TRUE,A);
535:   MatDestroy(&preallocator);
537:   ISCreateGeneral(PETSC_COMM_SELF,ncoo,coo_i,PETSC_COPY_VALUES,&is_coo_i);
538:   ISCreateGeneral(PETSC_COMM_SELF,ncoo,coo_j,PETSC_COPY_VALUES,&is_coo_j);
539:   PetscObjectCompose((PetscObject)A,"__PETSc_coo_i",(PetscObject)is_coo_i);
540:   PetscObjectCompose((PetscObject)A,"__PETSc_coo_j",(PetscObject)is_coo_j);
541:   ISDestroy(&is_coo_i);
542:   ISDestroy(&is_coo_j);
543:   return 0;
544: }

546: /*@
547:    MatSetPreallocationCOO - set preallocation for matrices using a coordinate format of the entries with global indices

549:    Collective on Mat

551:    Input Parameters:
552: +  A - matrix being preallocated
553: .  ncoo - number of entries
554: .  coo_i - row indices
555: -  coo_j - column indices

557:    Level: beginner

559:    Notes:
560:    Entries can be repeated, see MatSetValuesCOO(). Entries with negative row or column indices are allowed
561:    but will be ignored. The corresponding entries in MatSetValuesCOO() will be ignored too. Remote entries
562:    are allowed and will be properly added or inserted to the matrix, unless the matrix option MAT_IGNORE_OFF_PROC_ENTRIES
563:    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.

565:    The arrays coo_i and coo_j may be freed immediately after calling this function.

567: .seealso: MatSetValuesCOO(), MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation(), MatSetPreallocationCOOLocal(), DMSetMatrixPreallocateSkip()
568: @*/
569: PetscErrorCode MatSetPreallocationCOO(Mat A,PetscCount ncoo,const PetscInt coo_i[],const PetscInt coo_j[])
570: {
571:   PetscErrorCode (*f)(Mat,PetscCount,const PetscInt[],const PetscInt[]) = NULL;

577:   PetscLayoutSetUp(A->rmap);
578:   PetscLayoutSetUp(A->cmap);
579:   PetscObjectQueryFunction((PetscObject)A,"MatSetPreallocationCOO_C",&f);

581:   PetscLogEventBegin(MAT_PreallCOO,A,0,0,0);
582:   if (f) {
583:     (*f)(A,ncoo,coo_i,coo_j);
584:   } else { /* allow fallback, very slow */
585:     MatSetPreallocationCOO_Basic(A,ncoo,coo_i,coo_j);
586:   }
587:   PetscLogEventEnd(MAT_PreallCOO,A,0,0,0);
588:   A->preallocated = PETSC_TRUE;
589:   A->nonzerostate++;
590:   return 0;
591: }

593: /*@
594:    MatSetPreallocationCOOLocal - set preallocation for matrices using a coordinate format of the entries with local indices

596:    Collective on Mat

598:    Input Parameters:
599: +  A - matrix being preallocated
600: .  ncoo - number of entries
601: .  coo_i - row indices (local numbering; may be modified)
602: -  coo_j - column indices (local numbering; may be modified)

604:    Level: beginner

606:    Notes:
607:    The local indices are translated using the local to global mapping, thus MatSetLocalToGlobalMapping() must have been
608:    called prior to this function.

610:    The indices coo_i and coo_j may be modified within this function. They might be translated to corresponding global
611:    indices, but the caller should not rely on them having any specific value after this function returns. The arrays
612:    can be freed or reused immediately after this function returns.

614:    Entries can be repeated, see MatSetValuesCOO(). Entries with negative row or column indices are allowed
615:    but will be ignored. The corresponding entries in MatSetValuesCOO() will be ignored too. Remote entries
616:    are allowed and will be properly added or inserted to the matrix.

618: .seealso: MatSetValuesCOO(), MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation(), MatSetPreallocationCOO(), DMSetMatrixPreallocateSkip()
619: @*/
620: PetscErrorCode MatSetPreallocationCOOLocal(Mat A,PetscCount ncoo,PetscInt coo_i[],PetscInt coo_j[])
621: {
622:   PetscErrorCode (*f)(Mat,PetscCount,PetscInt[],PetscInt[]) = NULL;

629:   PetscLayoutSetUp(A->rmap);
630:   PetscLayoutSetUp(A->cmap);

632:   PetscObjectQueryFunction((PetscObject)A,"MatSetPreallocationCOOLocal_C",&f);
633:   if (f) {
634:     (*f)(A,ncoo,coo_i,coo_j);
635:     A->nonzerostate++;
636:   } else {
637:     ISLocalToGlobalMapping ltog_row,ltog_col;
638:     MatGetLocalToGlobalMapping(A,&ltog_row,&ltog_col);
639:     if (ltog_row) ISLocalToGlobalMappingApply(ltog_row,ncoo,coo_i,coo_i);
640:     if (ltog_col) ISLocalToGlobalMappingApply(ltog_col,ncoo,coo_j,coo_j);
641:     MatSetPreallocationCOO(A,ncoo,coo_i,coo_j);
642:   }
643:   A->preallocated = PETSC_TRUE;
644:   return 0;
645: }

647: /*@
648:    MatSetValuesCOO - set values at once in a matrix preallocated using MatSetPreallocationCOO()

650:    Collective on Mat

652:    Input Parameters:
653: +  A - matrix being preallocated
654: .  coo_v - the matrix values (can be NULL)
655: -  imode - the insert mode

657:    Level: beginner

659:    Notes: The values must follow the order of the indices prescribed with MatSetPreallocationCOO() or MatSetPreallocationCOOLocal().
660:           When repeated entries are specified in the COO indices the coo_v values are first properly summed, regardless of the value of imode.
661:           The imode flag indicates if coo_v must be added to the current values of the matrix (ADD_VALUES) or overwritten (INSERT_VALUES).
662:           MatAssemblyBegin() and MatAssemblyEnd() do not need to be called after this routine. It automatically handles the assembly process.

664: .seealso: MatSetPreallocationCOO(), MatSetPreallocationCOOLocal(), InsertMode, INSERT_VALUES, ADD_VALUES
665: @*/
666: PetscErrorCode MatSetValuesCOO(Mat A, const PetscScalar coo_v[], InsertMode imode)
667: {
668:   PetscErrorCode (*f)(Mat,const PetscScalar[],InsertMode) = NULL;

672:   MatCheckPreallocated(A,1);
674:   PetscObjectQueryFunction((PetscObject)A,"MatSetValuesCOO_C",&f);
675:   PetscLogEventBegin(MAT_SetVCOO,A,0,0,0);
676:   if (f) {
677:     (*f)(A,coo_v,imode);
678:   } else { /* allow fallback */
679:     MatSetValuesCOO_Basic(A,coo_v,imode);
680:   }
681:   PetscLogEventEnd(MAT_SetVCOO,A,0,0,0);
682:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
683:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
684:   return 0;
685: }

687: /*@
688:    MatSetBindingPropagates - Sets whether the state of being bound to the CPU for a GPU matrix type propagates to child and some other associated objects

690:    Input Parameters:
691: +  A - the matrix
692: -  flg - flag indicating whether the boundtocpu flag should be propagated

694:    Level: developer

696:    Notes:
697:    If the value of flg is set to true, the following will occur:

699:    MatCreateSubMatrices() and MatCreateRedundantMatrix() will bind created matrices to CPU if the input matrix is bound to the CPU.
700:    MatCreateVecs() will bind created vectors to CPU if the input matrix is bound to the CPU.
701:    The bindingpropagates flag itself is also propagated by the above routines.

703:    Developer Notes:
704:    If the fine-scale DMDA has the -dm_bind_below option set to true, then DMCreateInterpolationScale() calls MatSetBindingPropagates()
705:    on the restriction/interpolation operator to set the bindingpropagates flag to true.

707: .seealso: VecSetBindingPropagates(), MatGetBindingPropagates()
708: @*/
709: PetscErrorCode MatSetBindingPropagates(Mat A,PetscBool flg)
710: {
712: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
713:   A->bindingpropagates = flg;
714: #endif
715:   return 0;
716: }

718: /*@
719:    MatGetBindingPropagates - Gets whether the state of being bound to the CPU for a GPU matrix type propagates to child and some other associated objects

721:    Input Parameter:
722: .  A - the matrix

724:    Output Parameter:
725: .  flg - flag indicating whether the boundtocpu flag will be propagated

727:    Level: developer

729: .seealso: MatSetBindingPropagates()
730: @*/
731: PetscErrorCode MatGetBindingPropagates(Mat A,PetscBool *flg)
732: {
735: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
736:   *flg = A->bindingpropagates;
737: #else
738:   *flg = PETSC_FALSE;
739: #endif
740:   return 0;
741: }