Actual source code: gcreate.c

petsc-master 2019-03-19
Report Typos and Errors

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

  4: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat mat,PetscInt rbs, PetscInt cbs)
  5: {
  7:   if (!mat->preallocated) return(0);
  8:   if (mat->rmap->bs > 0 && mat->rmap->bs != rbs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot change row block size %D to %D\n",mat->rmap->bs,rbs);
  9:   if (mat->cmap->bs > 0 && mat->cmap->bs != cbs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot change column block size %D to %D\n",mat->cmap->bs,cbs);
 10:   return(0);
 11: }

 13: PETSC_INTERN PetscErrorCode MatShift_Basic(Mat Y,PetscScalar a)
 14: {
 16:   PetscInt       i,start,end;
 17:   PetscScalar    alpha = a;
 18:   PetscBool      prevoption;

 21:   MatGetOption(Y,MAT_NO_OFF_PROC_ENTRIES,&prevoption);
 22:   MatSetOption(Y,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);
 23:   MatGetOwnershipRange(Y,&start,&end);
 24:   for (i=start; i<end; i++) {
 25:     MatSetValues(Y,1,&i,1,&i,&alpha,ADD_VALUES);
 26:   }
 27:   MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);
 28:   MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);
 29:   MatSetOption(Y,MAT_NO_OFF_PROC_ENTRIES,prevoption);
 30:   return(0);
 31: }

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

 42:    Collective on MPI_Comm

 44:    Input Parameter:
 45: .  comm - MPI communicator

 47:    Output Parameter:
 48: .  A - the matrix

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

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

 62:    Notes:

 64:    Level: beginner

 66:    User manual sections:
 67: +   Section 3.1 Creating and Assembling Matrices
 68: -   Chapter 3 Matrices

 70: .keywords: matrix, create

 72: .seealso: MatCreateSeqAIJ(), MatCreateAIJ(),
 73:           MatCreateSeqDense(), MatCreateDense(),
 74:           MatCreateSeqBAIJ(), MatCreateBAIJ(),
 75:           MatCreateSeqSBAIJ(), MatCreateSBAIJ(),
 76:           MatConvert()
 77: @*/
 78: PetscErrorCode  MatCreate(MPI_Comm comm,Mat *A)
 79: {
 80:   Mat            B;


 86:   *A = NULL;
 87:   MatInitializePackage();

 89:   PetscHeaderCreate(B,MAT_CLASSID,"Mat","Matrix","Mat",comm,MatDestroy,MatView);
 90:   PetscLayoutCreate(comm,&B->rmap);
 91:   PetscLayoutCreate(comm,&B->cmap);
 92:   PetscStrallocpy(VECSTANDARD,&B->defaultvectype);

 94:   B->congruentlayouts = PETSC_DECIDE;
 95:   B->preallocated     = PETSC_FALSE;
 96:   *A                  = B;
 97:   return(0);
 98: }

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

103:    Logically Collective on Mat

105:    Input Parameters:
106: +  mat -  matrix obtained from MatCreate()
107: -  flg - PETSC_TRUE indicates you want the error generated

109:    Level: advanced

111: .keywords: Mat, set, initial guess, nonzero

113: .seealso: PCSetErrorIfFailure()
114: @*/
115: PetscErrorCode  MatSetErrorIfFailure(Mat mat,PetscBool flg)
116: {
120:   mat->erroriffailure = flg;
121:   return(0);
122: }

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

127:   Collective on Mat

129:   Input Parameters:
130: +  A - the matrix
131: .  m - number of local rows (or PETSC_DECIDE)
132: .  n - number of local columns (or PETSC_DECIDE)
133: .  M - number of global rows (or PETSC_DETERMINE)
134: -  N - number of global columns (or PETSC_DETERMINE)

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

140:    If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the
141:    user must ensure that they are chosen to be compatible with the
142:    vectors. To do this, one first considers the matrix-vector product
143:    'y = A x'. The 'm' that is used in the above routine must match the
144:    local size used in the vector creation routine VecCreateMPI() for 'y'.
145:    Likewise, the 'n' used must match that used as the local size in
146:    VecCreateMPI() for 'x'.

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

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

152:   Level: beginner

154: .seealso: MatGetSize(), PetscSplitOwnership()
155: @*/
156: PetscErrorCode  MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N)
157: {
162:   if (M > 0 && m > M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local column size %D cannot be larger than global column size %D",m,M);
163:   if (N > 0 && n > N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local row size %D cannot be larger than global row size %D",n,N);
164:   if ((A->rmap->n >= 0 && A->rmap->N >= 0) && (A->rmap->n != m || (M > 0 && A->rmap->N != M))) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot change/reset row sizes to %D local %D global after previously setting them to %D local %D global",m,M,A->rmap->n,A->rmap->N);
165:   if ((A->cmap->n >= 0 && A->cmap->N >= 0) && (A->cmap->n != n || (N > 0 && A->cmap->N != N))) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot change/reset column sizes to %D local %D global after previously setting them to %D local %D global",n,N,A->cmap->n,A->cmap->N);
166:   A->rmap->n = m;
167:   A->cmap->n = n;
168:   A->rmap->N = M > -1 ? M : A->rmap->N;
169:   A->cmap->N = N > -1 ? N : A->cmap->N;
170:   return(0);
171: }

173: /*@
174:    MatSetFromOptions - Creates a matrix where the type is determined
175:    from the options database. Generates a parallel MPI matrix if the
176:    communicator has more than one processor.  The default matrix type is
177:    AIJ, using the routines MatCreateSeqAIJ() and MatCreateAIJ() if
178:    you do not select a type in the options database.

180:    Collective on Mat

182:    Input Parameter:
183: .  A - the matrix

185:    Options Database Keys:
186: +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
187: .    -mat_type mpiaij   - AIJ type, uses MatCreateAIJ()
188: .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
189: .    -mat_type mpidense - dense type, uses MatCreateDense()
190: .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
191: -    -mat_type mpibaij  - block AIJ type, uses MatCreateBAIJ()

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

197:    Level: beginner

199: .keywords: matrix, create

201: .seealso: MatCreateSeqAIJ((), MatCreateAIJ(),
202:           MatCreateSeqDense(), MatCreateDense(),
203:           MatCreateSeqBAIJ(), MatCreateBAIJ(),
204:           MatCreateSeqSBAIJ(), MatCreateSBAIJ(),
205:           MatConvert()
206: @*/
207: PetscErrorCode  MatSetFromOptions(Mat B)
208: {
210:   const char     *deft = MATAIJ;
211:   char           type[256];
212:   PetscBool      flg,set;


217:   PetscObjectOptionsBegin((PetscObject)B);

219:   if (B->rmap->bs < 0) {
220:     PetscInt newbs = -1;
221:     PetscOptionsInt("-mat_block_size","Set the blocksize used to store the matrix","MatSetBlockSize",newbs,&newbs,&flg);
222:     if (flg) {
223:       PetscLayoutSetBlockSize(B->rmap,newbs);
224:       PetscLayoutSetBlockSize(B->cmap,newbs);
225:     }
226:   }

228:   PetscOptionsFList("-mat_type","Matrix type","MatSetType",MatList,deft,type,256,&flg);
229:   if (flg) {
230:     MatSetType(B,type);
231:   } else if (!((PetscObject)B)->type_name) {
232:     MatSetType(B,deft);
233:   }

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

240:   if (B->ops->setfromoptions) {
241:     (*B->ops->setfromoptions)(PetscOptionsObject,B);
242:   }

244:   flg  = PETSC_FALSE;
245:   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);
246:   if (set) {MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,flg);}
247:   flg  = PETSC_FALSE;
248:   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);
249:   if (set) {MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,flg);}

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 Arguments:
263: +  A - matrix being preallocated
264: .  bs - block size
265: .  dnnz - number of nonzero blocks per block row of diagonal part of parallel matrix
266: .  onnz - number of nonzero blocks per block row of off-diagonal part of parallel matrix
267: .  dnnzu - number of nonzero blocks per block row of upper-triangular part of diagonal part of parallel matrix
268: -  onnzu - number of nonzero 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: {
278:   void           (*aij)(void);
279:   void           (*is)(void);
280:   void           (*hyp)(void) = NULL;

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

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

334:         This is somewhat different from MatHeaderReplace() it would be nice to merge the code
335: */
336: PetscErrorCode MatHeaderMerge(Mat A,Mat *C)
337: {
339:   PetscInt       refct;
340:   PetscOps       Abops;
341:   struct _MatOps Aops;
342:   char           *mtype,*mname;

345:   /* save the parts of A we need */
346:   Abops = ((PetscObject)A)->bops[0];
347:   Aops  = A->ops[0];
348:   refct = ((PetscObject)A)->refct;
349:   mtype = ((PetscObject)A)->type_name;
350:   mname = ((PetscObject)A)->name;

352:   /* zero these so the destroy below does not free them */
353:   ((PetscObject)A)->type_name = 0;
354:   ((PetscObject)A)->name      = 0;

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

359:   PetscFree(A->defaultvectype);
360:   PetscLayoutDestroy(&A->rmap);
361:   PetscLayoutDestroy(&A->cmap);
362:   PetscFunctionListDestroy(&((PetscObject)A)->qlist);
363:   PetscObjectListDestroy(&((PetscObject)A)->olist);

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

368:   /* return the parts of A we saved */
369:   ((PetscObject)A)->bops[0]   = Abops;
370:   A->ops[0]                   = Aops;
371:   ((PetscObject)A)->refct     = refct;
372:   ((PetscObject)A)->type_name = mtype;
373:   ((PetscObject)A)->name      = mname;

375:   /* since these two are copied into A we do not want them destroyed in C */
376:   ((PetscObject)*C)->qlist = 0;
377:   ((PetscObject)*C)->olist = 0;

379:   PetscHeaderDestroy(C);
380:   return(0);
381: }
382: /*
383:         Replace A's header with that of C; the C object is then destroyed

385:         This is essentially code moved from MatDestroy()

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

389:         Used in DM hence is declared PETSC_EXTERN
390: */
391: PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat A,Mat *C)
392: {
393:   PetscErrorCode   ierr;
394:   PetscInt         refct;
395:   PetscObjectState state;
396:   struct _p_Mat    buffer;

401:   if (A == *C) return(0);
403:   if (((PetscObject)*C)->refct != 1) SETERRQ1(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"Object C has refct %D > 1, would leave hanging reference",((PetscObject)*C)->refct);

405:   /* swap C and A */
406:   refct = ((PetscObject)A)->refct;
407:   state = ((PetscObject)A)->state;
408:   PetscMemcpy(&buffer,A,sizeof(struct _p_Mat));
409:   PetscMemcpy(A,*C,sizeof(struct _p_Mat));
410:   PetscMemcpy(*C,&buffer,sizeof(struct _p_Mat));
411:   ((PetscObject)A)->refct = refct;
412:   ((PetscObject)A)->state = state + 1;

414:   ((PetscObject)*C)->refct = 1;
415:   MatDestroy(C);
416:   return(0);
417: }