xref: /petsc/src/mat/utils/gcreate.c (revision bce096a405b23f07feb1783633a56a7a90eceb89)
1 #define PETSCMAT_DLL
2 
3 #include "src/mat/matimpl.h"       /*I "petscmat.h"  I*/
4 #include "petscsys.h"
5 
6 #undef __FUNCT__
7 #define __FUNCT__ "MatPublish_Base"
8 static PetscErrorCode MatPublish_Base(PetscObject obj)
9 {
10   PetscFunctionBegin;
11   PetscFunctionReturn(0);
12 }
13 
14 
15 #undef __FUNCT__
16 #define __FUNCT__ "MatCreate"
17 /*@C
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.
25 
26    Collective on MPI_Comm
27 
28    Input Parameter:
29 .  comm - MPI communicator
30 
31    Output Parameter:
32 .  A - the matrix
33 
34    Options Database Keys:
35 +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
36 .    -mat_type mpiaij   - AIJ type, uses MatCreateMPIAIJ()
37 .    -mat_type seqbdiag - block diagonal type, uses MatCreateSeqBDiag()
38 .    -mat_type mpibdiag - block diagonal type, uses MatCreateMPIBDiag()
39 .    -mat_type mpirowbs - rowbs type, uses MatCreateMPIRowbs()
40 .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
41 .    -mat_type mpidense - dense type, uses MatCreateMPIDense()
42 .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
43 -    -mat_type mpibaij  - block AIJ type, uses MatCreateMPIBAIJ()
44 
45    Even More Options Database Keys:
46    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
47    for additional format-specific options.
48 
49    Notes:
50 
51    Level: beginner
52 
53    User manual sections:
54 +   sec_matcreate
55 -   chapter_matrices
56 
57 .keywords: matrix, create
58 
59 .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
60           MatCreateSeqBDiag(),MatCreateMPIBDiag(),
61           MatCreateSeqDense(), MatCreateMPIDense(),
62           MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
63           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
64           MatConvert()
65 @*/
66 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm comm,Mat *A)
67 {
68   Mat            B;
69   PetscErrorCode ierr;
70 
71   PetscFunctionBegin;
72   PetscValidPointer(A,2);
73 
74   *A = PETSC_NULL;
75 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
76   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
77 #endif
78 
79   ierr = PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,0,"Mat",comm,MatDestroy,MatView);CHKERRQ(ierr);
80   B->m             = -1;
81   B->M             = -1;
82   B->n             = -1;
83   B->N             = -1;
84   B->bs            = 1;
85   B->preallocated  = PETSC_FALSE;
86   B->bops->publish = MatPublish_Base;
87   *A               = B;
88   PetscFunctionReturn(0);
89 }
90 
91 #undef __FUNCT__
92 #define __FUNCT__ "MatSetSizes"
93 /*@
94   MatSetSizes - Sets the local and global sizes, and checks to determine compatibility
95 
96   Collective on Mat
97 
98   Input Parameters:
99 +  A - the matrix
100 .  m - number of local rows (or PETSC_DECIDE)
101 .  n - number of local columns (or PETSC_DECIDE)
102 .  M - number of global rows (or PETSC_DETERMINE)
103 -  N - number of global columns (or PETSC_DETERMINE)
104 
105    Notes:
106    m (n) and M (N) cannot be both PETSC_DECIDE
107    If one processor calls this with M (N) of PETSC_DECIDE then all processors must, otherwise the program will hang.
108 
109    If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the
110    user must ensure that they are chosen to be compatible with the
111    vectors. To do this, one first considers the matrix-vector product
112    'y = A x'. The 'm' that is used in the above routine must match the
113    local size used in the vector creation routine VecCreateMPI() for 'y'.
114    Likewise, the 'n' used must match that used as the local size in
115    VecCreateMPI() for 'x'.
116 
117   Level: beginner
118 
119 .seealso: MatGetSize(), PetscSplitOwnership()
120 @*/
121 PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N)
122 {
123   PetscFunctionBegin;
124   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
125   if (M > 0 && m > M) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Local column size %D cannot be larger than global column size %D",m,M);
126   if (N > 0 && n > N) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Local row size %D cannot be larger than global row size %D",n,N);
127   if ((A->m >= 0 || A->M >= 0) && (A->m != m || A->M != 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->m,A->M);
128   if ((A->n >= 0 || A->N >= 0) && (A->n != n || A->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->n,A->N);
129   A->m = m;
130   A->n = n;
131   A->M = M;
132   A->N = N;
133   PetscFunctionReturn(0);
134 }
135 
136 #undef __FUNCT__
137 #define __FUNCT__ "MatSetFromOptions"
138 /*@C
139    MatSetFromOptions - Creates a matrix where the type is determined
140    from the options database. Generates a parallel MPI matrix if the
141    communicator has more than one processor.  The default matrix type is
142    AIJ, using the routines MatCreateSeqAIJ() and MatCreateMPIAIJ() if
143    you do not select a type in the options database.
144 
145    Collective on Mat
146 
147    Input Parameter:
148 .  A - the matrix
149 
150    Options Database Keys:
151 +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
152 .    -mat_type mpiaij   - AIJ type, uses MatCreateMPIAIJ()
153 .    -mat_type seqbdiag - block diagonal type, uses MatCreateSeqBDiag()
154 .    -mat_type mpibdiag - block diagonal type, uses MatCreateMPIBDiag()
155 .    -mat_type mpirowbs - rowbs type, uses MatCreateMPIRowbs()
156 .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
157 .    -mat_type mpidense - dense type, uses MatCreateMPIDense()
158 .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
159 -    -mat_type mpibaij  - block AIJ type, uses MatCreateMPIBAIJ()
160 
161    Even More Options Database Keys:
162    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
163    for additional format-specific options.
164 
165    Level: beginner
166 
167 .keywords: matrix, create
168 
169 .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
170           MatCreateSeqBDiag(),MatCreateMPIBDiag(),
171           MatCreateSeqDense(), MatCreateMPIDense(),
172           MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
173           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
174           MatConvert()
175 @*/
176 PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat B)
177 {
178   PetscErrorCode ierr;
179   char           mtype[256];
180   PetscTruth     flg;
181 
182   PetscFunctionBegin;
183   ierr = PetscOptionsGetString(B->prefix,"-mat_type",mtype,256,&flg);CHKERRQ(ierr);
184   if (flg) {
185     ierr = MatSetType(B,mtype);CHKERRQ(ierr);
186   }
187   if (!B->type_name) {
188     ierr = MatSetType(B,MATAIJ);CHKERRQ(ierr);
189   }
190   PetscFunctionReturn(0);
191 }
192 
193 #undef __FUNCT__
194 #define __FUNCT__ "MatSetUpPreallocation"
195 /*@C
196    MatSetUpPreallocation
197 
198    Collective on Mat
199 
200    Input Parameter:
201 .  A - the matrix
202 
203    Level: beginner
204 
205 .keywords: matrix, create
206 
207 .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
208           MatCreateSeqBDiag(),MatCreateMPIBDiag(),
209           MatCreateSeqDense(), MatCreateMPIDense(),
210           MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
211           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
212           MatConvert()
213 @*/
214 PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat B)
215 {
216   PetscErrorCode ierr;
217 
218   PetscFunctionBegin;
219   if (B->ops->setuppreallocation) {
220     ierr = PetscLogInfo((B,"MatSetUpPreallocation: Warning not preallocating matrix storage\n"));CHKERRQ(ierr);
221     ierr = (*B->ops->setuppreallocation)(B);CHKERRQ(ierr);
222     B->ops->setuppreallocation = 0;
223     B->preallocated            = PETSC_TRUE;
224   }
225   PetscFunctionReturn(0);
226 }
227 
228 /*
229         Copies from Cs header to A
230 */
231 #undef __FUNCT__
232 #define __FUNCT__ "MatHeaderCopy"
233 PetscErrorCode MatHeaderCopy(Mat A,Mat C)
234 {
235   PetscErrorCode ierr;
236   PetscInt       refct;
237   PetscOps       *Abops;
238   MatOps         Aops;
239   char           *mtype,*mname;
240   void           *spptr;
241 
242   PetscFunctionBegin;
243   /* free all the interior data structures from mat */
244   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
245 
246   ierr = PetscMapDestroy(A->rmap);CHKERRQ(ierr);
247   ierr = PetscMapDestroy(A->cmap);CHKERRQ(ierr);
248 
249   /* save the parts of A we need */
250   Abops = A->bops;
251   Aops  = A->ops;
252   refct = A->refct;
253   mtype = A->type_name;
254   mname = A->name;
255   spptr = A->spptr;
256 
257   if (C->spptr) {
258     ierr = PetscFree(C->spptr);CHKERRQ(ierr);
259     C->spptr = PETSC_NULL;
260   }
261 
262   /* copy C over to A */
263   ierr  = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
264 
265   /* return the parts of A we saved */
266   A->bops      = Abops;
267   A->ops       = Aops;
268   A->qlist     = 0;
269   A->refct     = refct;
270   A->type_name = mtype;
271   A->name      = mname;
272   A->spptr     = spptr;
273 
274   ierr = PetscHeaderDestroy(C);CHKERRQ(ierr);
275   PetscFunctionReturn(0);
276 }
277 /*
278         Replace A's header with that of C
279         This is essentially code moved from MatDestroy
280 */
281 #undef __FUNCT__
282 #define __FUNCT__ "MatHeaderReplace"
283 PetscErrorCode MatHeaderReplace(Mat A,Mat C)
284 {
285   PetscErrorCode ierr;
286 
287   PetscFunctionBegin;
288   /* free all the interior data structures from mat */
289   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
290   ierr = PetscHeaderDestroy_Private((PetscObject)A);CHKERRQ(ierr);
291   ierr = PetscMapDestroy(A->rmap);CHKERRQ(ierr);
292   ierr = PetscMapDestroy(A->cmap);CHKERRQ(ierr);
293 
294   /* copy C over to A */
295   if (C) {
296     ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
297     ierr = PetscLogObjectDestroy((PetscObject)C);CHKERRQ(ierr);
298     ierr = PetscFree(C);CHKERRQ(ierr);
299   }
300   PetscFunctionReturn(0);
301 }
302