xref: /petsc/src/mat/utils/gcreate.c (revision 4bf112e79d3f87ddad6dbf7b300a6f889c1e7d34)
1 
2 #include "src/mat/matimpl.h"       /*I "petscmat.h"  I*/
3 #include "petscsys.h"
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "MatPublish_Base"
7 static PetscErrorCode MatPublish_Base(PetscObject obj)
8 {
9   PetscFunctionBegin;
10   PetscFunctionReturn(0);
11 }
12 
13 
14 #undef __FUNCT__
15 #define __FUNCT__ "MatCreate"
16 /*@C
17    MatCreate - Creates a matrix where the type is determined
18    from either a call to MatSetType() or from the options database
19    with a call to MatSetFromOptions(). The default matrix type is
20    AIJ, using the routines MatCreateSeqAIJ() or MatCreateMPIAIJ()
21    if you do not set a type in the options database. If you never
22    call MatSetType() or MatSetFromOptions() it will generate an
23    error when you try to use the matrix.
24 
25    Collective on MPI_Comm
26 
27    Input Parameters:
28 +  m - number of local rows (or PETSC_DECIDE)
29 .  n - number of local columns (or PETSC_DECIDE)
30 .  M - number of global rows (or PETSC_DETERMINE)
31 .  N - number of global columns (or PETSC_DETERMINE)
32 -  comm - MPI communicator
33 
34    Output Parameter:
35 .  A - the matrix
36 
37    Options Database Keys:
38 +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
39 .    -mat_type mpiaij   - AIJ type, uses MatCreateMPIAIJ()
40 .    -mat_type seqbdiag - block diagonal type, uses MatCreateSeqBDiag()
41 .    -mat_type mpibdiag - block diagonal type, uses MatCreateMPIBDiag()
42 .    -mat_type mpirowbs - rowbs type, uses MatCreateMPIRowbs()
43 .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
44 .    -mat_type mpidense - dense type, uses MatCreateMPIDense()
45 .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
46 -    -mat_type mpibaij  - block AIJ type, uses MatCreateMPIBAIJ()
47 
48    Even More Options Database Keys:
49    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
50    for additional format-specific options.
51 
52    Notes:
53    If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the
54    user must ensure that they are chosen to be compatible with the
55    vectors. To do this, one first considers the matrix-vector product
56    'y = A x'. The 'm' that is used in the above routine must match the
57    local size used in the vector creation routine VecCreateMPI() for 'y'.
58    Likewise, the 'n' used must match that used as the local size in
59    VecCreateMPI() for 'x'.
60 
61    Level: beginner
62 
63    User manual sections:
64 +   sec_matcreate
65 -   chapter_matrices
66 
67 .keywords: matrix, create
68 
69 .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
70           MatCreateSeqBDiag(),MatCreateMPIBDiag(),
71           MatCreateSeqDense(), MatCreateMPIDense(),
72           MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
73           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
74           MatConvert()
75 @*/
76 PetscErrorCode MatCreate(MPI_Comm comm,int m,int n,int M,int N,Mat *A)
77 {
78   Mat B;
79 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
80   PetscErrorCode ierr;
81 #endif
82 
83   PetscFunctionBegin;
84   PetscValidPointer(A,6);
85   *A = PETSC_NULL;
86 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
87   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
88 #endif
89 
90   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,0,"Mat",comm,MatDestroy,MatView);
91   PetscLogObjectCreate(B);
92 
93   B->m = m;
94   B->n = n;
95   B->M = M;
96   B->N = N;
97 
98   B->preallocated  = PETSC_FALSE;
99   B->bops->publish = MatPublish_Base;
100   *A = B;
101   PetscFunctionReturn(0);
102 }
103 
104 #undef __FUNCT__
105 #define __FUNCT__ "MatSetFromOptions"
106 /*@C
107    MatSetFromOptions - Creates a matrix where the type is determined
108    from the options database. Generates a parallel MPI matrix if the
109    communicator has more than one processor.  The default matrix type is
110    AIJ, using the routines MatCreateSeqAIJ() and MatCreateMPIAIJ() if
111    you do not select a type in the options database.
112 
113    Collective on Mat
114 
115    Input Parameter:
116 .  A - the matrix
117 
118    Options Database Keys:
119 +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
120 .    -mat_type mpiaij   - AIJ type, uses MatCreateMPIAIJ()
121 .    -mat_type seqbdiag - block diagonal type, uses MatCreateSeqBDiag()
122 .    -mat_type mpibdiag - block diagonal type, uses MatCreateMPIBDiag()
123 .    -mat_type mpirowbs - rowbs type, uses MatCreateMPIRowbs()
124 .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
125 .    -mat_type mpidense - dense type, uses MatCreateMPIDense()
126 .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
127 -    -mat_type mpibaij  - block AIJ type, uses MatCreateMPIBAIJ()
128 
129    Even More Options Database Keys:
130    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
131    for additional format-specific options.
132 
133    Level: beginner
134 
135 .keywords: matrix, create
136 
137 .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
138           MatCreateSeqBDiag(),MatCreateMPIBDiag(),
139           MatCreateSeqDense(), MatCreateMPIDense(),
140           MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
141           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
142           MatConvert()
143 @*/
144 PetscErrorCode MatSetFromOptions(Mat B)
145 {
146   PetscErrorCode ierr;
147   int size;
148   char       mtype[256];
149   PetscTruth flg;
150 
151   PetscFunctionBegin;
152   ierr = PetscOptionsGetString(B->prefix,"-mat_type",mtype,256,&flg);CHKERRQ(ierr);
153   if (flg) {
154     ierr = MatSetType(B,mtype);CHKERRQ(ierr);
155   }
156   if (!B->type_name) {
157     ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
158     ierr = MatSetType(B,MATAIJ);CHKERRQ(ierr);
159   }
160 #if defined(__cplusplus) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && defined(PETSC_HAVE_CXX_NAMESPACE)
161   ierr = MatESISetFromOptions(B);CHKERRQ(ierr);
162 #endif
163   PetscFunctionReturn(0);
164 }
165 
166 #undef __FUNCT__
167 #define __FUNCT__ "MatSetUpPreallocation"
168 /*@C
169    MatSetUpPreallocation
170 
171    Collective on Mat
172 
173    Input Parameter:
174 .  A - the matrix
175 
176    Level: beginner
177 
178 .keywords: matrix, create
179 
180 .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
181           MatCreateSeqBDiag(),MatCreateMPIBDiag(),
182           MatCreateSeqDense(), MatCreateMPIDense(),
183           MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
184           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
185           MatConvert()
186 @*/
187 PetscErrorCode MatSetUpPreallocation(Mat B)
188 {
189   PetscErrorCode ierr;
190 
191   PetscFunctionBegin;
192   if (B->ops->setuppreallocation) {
193     PetscLogInfo(B,"MatSetUpPreallocation: Warning not preallocating matrix storage");
194     ierr = (*B->ops->setuppreallocation)(B);CHKERRQ(ierr);
195     B->ops->setuppreallocation = 0;
196     B->preallocated            = PETSC_TRUE;
197   }
198   PetscFunctionReturn(0);
199 }
200 
201 /*
202         Copies from Cs header to A
203 */
204 #undef __FUNCT__
205 #define __FUNCT__ "MatHeaderCopy"
206 PetscErrorCode MatHeaderCopy(Mat A,Mat C)
207 {
208   PetscErrorCode ierr;
209   int refct;
210   PetscOps    *Abops;
211   MatOps      Aops;
212   char        *mtype,*mname;
213 
214   PetscFunctionBegin;
215   /* free all the interior data structures from mat */
216   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
217 
218   ierr = PetscMapDestroy(A->rmap);CHKERRQ(ierr);
219   ierr = PetscMapDestroy(A->cmap);CHKERRQ(ierr);
220 
221   /* save the parts of A we need */
222   Abops = A->bops;
223   Aops  = A->ops;
224   refct = A->refct;
225   mtype = A->type_name;
226   mname = A->name;
227 
228   /* copy C over to A */
229   ierr  = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
230 
231   /* return the parts of A we saved */
232   A->bops      = Abops;
233   A->ops       = Aops;
234   A->qlist     = 0;
235   A->refct     = refct;
236   A->type_name = mtype;
237   A->name      = mname;
238 
239   PetscLogObjectDestroy(C);
240   PetscHeaderDestroy(C);
241   PetscFunctionReturn(0);
242 }
243