Actual source code: gcreate.c

petsc-3.8.4 2018-03-24
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);

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

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

102:    Logically Collective on Mat

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

108:    Level: advanced

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

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

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

126:   Collective on Mat

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

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

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

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

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

151:   Level: beginner

153: .seealso: MatGetSize(), PetscSplitOwnership()
154: @*/
155: PetscErrorCode  MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N)
156: {
161:   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);
162:   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);
163:   if ((A->rmap->n >= 0 && A->rmap->N >= 0) && (A->rmap->n != m || 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);
164:   if ((A->cmap->n >= 0 && A->cmap->N >= 0) && (A->cmap->n != n || 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);
165:   A->rmap->n = m;
166:   A->cmap->n = n;
167:   A->rmap->N = M;
168:   A->cmap->N = N;
169:   return(0);
170: }

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

179:    Collective on Mat

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

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

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

196:    Level: beginner

198: .keywords: matrix, create

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


216:   PetscObjectOptionsBegin((PetscObject)B);

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

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

234:   PetscOptionsName("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",&B->checksymmetryonassembly);
235:   PetscOptionsReal("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",B->checksymmetrytol,&B->checksymmetrytol,NULL);
236:   PetscOptionsBool("-mat_null_space_test","Checks if provided null space is correct in MatAssemblyEnd()","MatSetNullSpaceTest",B->checknullspaceonassembly,&B->checknullspaceonassembly,NULL);

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

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

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

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

258:    Collective on Mat

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

268:    Level: beginner

270: .seealso: MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation(),
271:           PetscSplitOwnership()
272: @*/
273: PetscErrorCode MatXAIJSetPreallocation(Mat A,PetscInt bs,const PetscInt dnnz[],const PetscInt onnz[],const PetscInt dnnzu[],const PetscInt onnzu[])
274: {
276:   void           (*aij)(void);
277:   void           (*is)(void);

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

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

322:         This is somewhat different from MatHeaderReplace() it would be nice to merge the code
323: */
324: PetscErrorCode MatHeaderMerge(Mat A,Mat *C)
325: {
327:   PetscInt       refct;
328:   PetscOps       Abops;
329:   struct _MatOps Aops;
330:   char           *mtype,*mname;

333:   /* save the parts of A we need */
334:   Abops = ((PetscObject)A)->bops[0];
335:   Aops  = A->ops[0];
336:   refct = ((PetscObject)A)->refct;
337:   mtype = ((PetscObject)A)->type_name;
338:   mname = ((PetscObject)A)->name;

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

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

347:   PetscLayoutDestroy(&A->rmap);
348:   PetscLayoutDestroy(&A->cmap);
349:   PetscFunctionListDestroy(&((PetscObject)A)->qlist);
350:   PetscObjectListDestroy(&((PetscObject)A)->olist);

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

355:   /* return the parts of A we saved */
356:   ((PetscObject)A)->bops[0]   = Abops;
357:   A->ops[0]                   = Aops;
358:   ((PetscObject)A)->refct     = refct;
359:   ((PetscObject)A)->type_name = mtype;
360:   ((PetscObject)A)->name      = mname;

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

366:   PetscHeaderDestroy(C);
367:   return(0);
368: }
369: /*
370:         Replace A's header with that of C; the C object is then destroyed

372:         This is essentially code moved from MatDestroy()

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

376:         Used in DM hence is declared PETSC_EXTERN
377: */
378: PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat A,Mat *C)
379: {
380:   PetscErrorCode   ierr;
381:   PetscInt         refct;
382:   PetscObjectState state;
383:   struct _p_Mat    buffer;

388:   if (A == *C) return(0);
390:   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);

392:   /* swap C and A */
393:   refct = ((PetscObject)A)->refct;
394:   state = ((PetscObject)A)->state;
395:   PetscMemcpy(&buffer,A,sizeof(struct _p_Mat));
396:   PetscMemcpy(A,*C,sizeof(struct _p_Mat));
397:   PetscMemcpy(*C,&buffer,sizeof(struct _p_Mat));
398:   ((PetscObject)A)->refct = refct;
399:   ((PetscObject)A)->state = state + 1;

401:   ((PetscObject)*C)->refct = 1;
402:   MatDestroy(C);
403:   return(0);
404: }