Actual source code: gcreate.c
1: #define PETSCMAT_DLL
3: #include private/matimpl.h
5: #if 0
8: static PetscErrorCode MatPublish_Base(PetscObject obj)
9: {
11: return(0);
12: }
13: #endif
17: /*@
18: MatCreate - Creates a matrix where the type is determined
19: from either a call to MatSetType() or from the options database
20: with a call to MatSetFromOptions(). The default matrix type is
21: AIJ, using the routines MatCreateSeqAIJ() or MatCreateMPIAIJ()
22: if you do not set a type in the options database. If you never
23: call MatSetType() or MatSetFromOptions() it will generate an
24: error when you try to use the matrix.
26: Collective on MPI_Comm
28: Input Parameter:
29: . comm - MPI communicator
30:
31: Output Parameter:
32: . A - the matrix
34: Options Database Keys:
35: + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ()
36: . -mat_type mpiaij - AIJ type, uses MatCreateMPIAIJ()
37: . -mat_type seqdense - dense type, uses MatCreateSeqDense()
38: . -mat_type mpidense - dense type, uses MatCreateMPIDense()
39: . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ()
40: - -mat_type mpibaij - block AIJ type, uses MatCreateMPIBAIJ()
42: Even More Options Database Keys:
43: See the manpages for particular formats (e.g., MatCreateSeqAIJ())
44: for additional format-specific options.
46: Notes:
48: Level: beginner
50: User manual sections:
51: + Section 3.1 Creating and Assembling Matrices
52: - Chapter 3 Matrices
54: .keywords: matrix, create
56: .seealso: MatCreateSeqAIJ(), MatCreateMPIAIJ(),
57: MatCreateSeqDense(), MatCreateMPIDense(),
58: MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
59: MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
60: MatConvert()
61: @*/
62: PetscErrorCode MatCreate(MPI_Comm comm,Mat *A)
63: {
64: Mat B;
70: *A = PETSC_NULL;
71: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
72: MatInitializePackage(PETSC_NULL);
73: #endif
75: PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,0,"Mat",comm,MatDestroy,MatView);
76: PetscLayoutCreate(comm,&B->rmap);
77: PetscLayoutCreate(comm,&B->cmap);
78: B->preallocated = PETSC_FALSE;
79: *A = B;
80: return(0);
81: }
85: /*@
86: MatSetSizes - Sets the local and global sizes, and checks to determine compatibility
88: Collective on Mat
90: Input Parameters:
91: + A - the matrix
92: . m - number of local rows (or PETSC_DECIDE)
93: . n - number of local columns (or PETSC_DECIDE)
94: . M - number of global rows (or PETSC_DETERMINE)
95: - N - number of global columns (or PETSC_DETERMINE)
97: Notes:
98: m (n) and M (N) cannot be both PETSC_DECIDE
99: If one processor calls this with M (N) of PETSC_DECIDE then all processors must, otherwise the program will hang.
101: If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the
102: user must ensure that they are chosen to be compatible with the
103: vectors. To do this, one first considers the matrix-vector product
104: 'y = A x'. The 'm' that is used in the above routine must match the
105: local size used in the vector creation routine VecCreateMPI() for 'y'.
106: Likewise, the 'n' used must match that used as the local size in
107: VecCreateMPI() for 'x'.
109: Level: beginner
111: .seealso: MatGetSize(), PetscSplitOwnership()
112: @*/
113: PetscErrorCode MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N)
114: {
119: if (M > 0 && m > M) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Local column size %D cannot be larger than global column size %D",m,M);
120: if (N > 0 && n > N) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Local row size %D cannot be larger than global row size %D",n,N);
121: if (A->ops->setsizes) {
122: /* Since this will not be set until the type has been set, this will NOT be called on the initial
123: call of MatSetSizes() (which must be called BEFORE MatSetType() */
124: (*A->ops->setsizes)(A,m,n,M,N);
125: } else {
126: if ((A->rmap->n >= 0 || A->rmap->N >= 0) && (A->rmap->n != m || A->rmap->N != M)) SETERRQ4(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);
127: if ((A->cmap->n >= 0 || A->cmap->N >= 0) && (A->cmap->n != n || A->cmap->N != N)) SETERRQ4(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);
128: }
129: A->rmap->n = m;
130: A->cmap->n = n;
131: A->rmap->N = M;
132: A->cmap->N = N;
133: if (A->ops->create) {
134: (*A->ops->create)(A);
135: A->ops->create = 0;
136: }
138: return(0);
139: }
143: /*@
144: MatSetFromOptions - Creates a matrix where the type is determined
145: from the options database. Generates a parallel MPI matrix if the
146: communicator has more than one processor. The default matrix type is
147: AIJ, using the routines MatCreateSeqAIJ() and MatCreateMPIAIJ() if
148: you do not select a type in the options database.
150: Collective on Mat
152: Input Parameter:
153: . A - the matrix
155: Options Database Keys:
156: + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ()
157: . -mat_type mpiaij - AIJ type, uses MatCreateMPIAIJ()
158: . -mat_type seqdense - dense type, uses MatCreateSeqDense()
159: . -mat_type mpidense - dense type, uses MatCreateMPIDense()
160: . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ()
161: - -mat_type mpibaij - block AIJ type, uses MatCreateMPIBAIJ()
163: Even More Options Database Keys:
164: See the manpages for particular formats (e.g., MatCreateSeqAIJ())
165: for additional format-specific options.
167: Level: beginner
169: .keywords: matrix, create
171: .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
172: MatCreateSeqDense(), MatCreateMPIDense(),
173: MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
174: MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
175: MatConvert()
176: @*/
177: PetscErrorCode MatSetFromOptions(Mat B)
178: {
180: const char *deft = MATAIJ;
181: char type[256];
182: PetscTruth flg;
187: PetscOptionsBegin(((PetscObject)B)->comm,((PetscObject)B)->prefix,"Matrix options","Mat");
188: PetscOptionsList("-mat_type","Matrix type","MatSetType",MatList,deft,type,256,&flg);
189: if (flg) {
190: MatSetType(B,type);
191: } else if (!((PetscObject)B)->type_name) {
192: MatSetType(B,deft);
193: }
195: if (B->ops->setfromoptions) {
196: (*B->ops->setfromoptions)(B);
197: }
199: PetscOptionsEnd();
201: return(0);
202: }
206: /*@
207: MatSetUpPreallocation
209: Collective on Mat
211: Input Parameter:
212: . A - the matrix
214: Level: beginner
216: .keywords: matrix, create
218: .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
219: MatCreateSeqDense(), MatCreateMPIDense(),
220: MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
221: MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
222: MatConvert()
223: @*/
224: PetscErrorCode MatSetUpPreallocation(Mat B)
225: {
229: if (!B->preallocated && B->ops->setuppreallocation) {
230: PetscInfo(B,"Warning not preallocating matrix storage\n");
231: (*B->ops->setuppreallocation)(B);
232: }
233: B->preallocated = PETSC_TRUE;
234: return(0);
235: }
237: /*
238: Copies from Cs header to A
240: This is somewhat different from MatHeaderReplace() it would be nice to merge the code
241: */
244: PetscErrorCode MatHeaderCopy(Mat A,Mat C)
245: {
247: PetscInt refct;
248: PetscOps *Abops;
249: MatOps Aops;
250: char *mtype,*mname;
251: void *spptr;
254: /* save the parts of A we need */
255: Abops = ((PetscObject)A)->bops;
256: Aops = A->ops;
257: refct = ((PetscObject)A)->refct;
258: mtype = ((PetscObject)A)->type_name;
259: mname = ((PetscObject)A)->name;
260: spptr = A->spptr;
262: /* zero these so the destroy below does not free them */
263: ((PetscObject)A)->type_name = 0;
264: ((PetscObject)A)->name = 0;
266: /* free all the interior data structures from mat */
267: (*A->ops->destroy)(A);
269: PetscFree(C->spptr);
271: PetscLayoutDestroy(A->rmap);
272: PetscLayoutDestroy(A->cmap);
273: PetscFListDestroy(&((PetscObject)A)->qlist);
274: PetscOListDestroy(((PetscObject)A)->olist);
276: /* copy C over to A */
277: PetscMemcpy(A,C,sizeof(struct _p_Mat));
279: /* return the parts of A we saved */
280: ((PetscObject)A)->bops = Abops;
281: A->ops = Aops;
282: ((PetscObject)A)->refct = refct;
283: ((PetscObject)A)->type_name = mtype;
284: ((PetscObject)A)->name = mname;
285: A->spptr = spptr;
287: /* since these two are copied into A we do not want them destroyed in C */
288: ((PetscObject)C)->qlist = 0;
289: ((PetscObject)C)->olist = 0;
290: PetscHeaderDestroy(C);
291: return(0);
292: }
293: /*
294: Replace A's header with that of C
295: This is essentially code moved from MatDestroy
297: This is somewhat different from MatHeaderCopy() it would be nice to merge the code
298: */
301: PetscErrorCode MatHeaderReplace(Mat A,Mat C)
302: {
306: if (A == C) return(0);
308: /* free all the interior data structures from mat */
309: (*A->ops->destroy)(A);
310: PetscHeaderDestroy_Private((PetscObject)A);
311: PetscFree(A->ops);
312: PetscLayoutDestroy(A->rmap);
313: PetscLayoutDestroy(A->cmap);
314: PetscFree(A->spptr);
315:
316: /* copy C over to A */
317: if (C) {
318: PetscMemcpy(A,C,sizeof(struct _p_Mat));
319: PetscLogObjectDestroy((PetscObject)C);
320: PetscFree(C);
321: }
322: return(0);
323: }