xref: /petsc/src/mat/utils/gcreate.c (revision 52e6d16ba989af9362d0fcebb12e73dd7c9666ed) !
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,PetscInt m,PetscInt n,PetscInt M,PetscInt 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   if (M > 0 && m > M) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Local column size %D cannot be larger than global column size %D",m,M);
86   if (N > 0 && n > N) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Local row size %D cannot be larger than global row size %D",n,N);
87 
88   *A = PETSC_NULL;
89 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
90   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
91 #endif
92 
93   ierr = PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,0,"Mat",comm,MatDestroy,MatView);CHKERRQ(ierr);
94 
95   B->m             = m;
96   B->n             = n;
97   B->M             = M;
98   B->N             = N;
99   B->bs            = 1;
100   B->preallocated  = PETSC_FALSE;
101   B->bops->publish = MatPublish_Base;
102   *A               = B;
103   PetscFunctionReturn(0);
104 }
105 
106 #undef __FUNCT__
107 #define __FUNCT__ "MatSetFromOptions"
108 /*@C
109    MatSetFromOptions - Creates a matrix where the type is determined
110    from the options database. Generates a parallel MPI matrix if the
111    communicator has more than one processor.  The default matrix type is
112    AIJ, using the routines MatCreateSeqAIJ() and MatCreateMPIAIJ() if
113    you do not select a type in the options database.
114 
115    Collective on Mat
116 
117    Input Parameter:
118 .  A - the matrix
119 
120    Options Database Keys:
121 +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
122 .    -mat_type mpiaij   - AIJ type, uses MatCreateMPIAIJ()
123 .    -mat_type seqbdiag - block diagonal type, uses MatCreateSeqBDiag()
124 .    -mat_type mpibdiag - block diagonal type, uses MatCreateMPIBDiag()
125 .    -mat_type mpirowbs - rowbs type, uses MatCreateMPIRowbs()
126 .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
127 .    -mat_type mpidense - dense type, uses MatCreateMPIDense()
128 .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
129 -    -mat_type mpibaij  - block AIJ type, uses MatCreateMPIBAIJ()
130 
131    Even More Options Database Keys:
132    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
133    for additional format-specific options.
134 
135    Level: beginner
136 
137 .keywords: matrix, create
138 
139 .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
140           MatCreateSeqBDiag(),MatCreateMPIBDiag(),
141           MatCreateSeqDense(), MatCreateMPIDense(),
142           MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
143           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
144           MatConvert()
145 @*/
146 PetscErrorCode MatSetFromOptions(Mat B)
147 {
148   PetscErrorCode ierr;
149   char           mtype[256];
150   PetscTruth     flg;
151 
152   PetscFunctionBegin;
153   ierr = PetscOptionsGetString(B->prefix,"-mat_type",mtype,256,&flg);CHKERRQ(ierr);
154   if (flg) {
155     ierr = MatSetType(B,mtype);CHKERRQ(ierr);
156   }
157   if (!B->type_name) {
158     ierr = MatSetType(B,MATAIJ);CHKERRQ(ierr);
159   }
160   PetscFunctionReturn(0);
161 }
162 
163 #undef __FUNCT__
164 #define __FUNCT__ "MatSetUpPreallocation"
165 /*@C
166    MatSetUpPreallocation
167 
168    Collective on Mat
169 
170    Input Parameter:
171 .  A - the matrix
172 
173    Level: beginner
174 
175 .keywords: matrix, create
176 
177 .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
178           MatCreateSeqBDiag(),MatCreateMPIBDiag(),
179           MatCreateSeqDense(), MatCreateMPIDense(),
180           MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
181           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
182           MatConvert()
183 @*/
184 PetscErrorCode MatSetUpPreallocation(Mat B)
185 {
186   PetscErrorCode ierr;
187 
188   PetscFunctionBegin;
189   if (B->ops->setuppreallocation) {
190     PetscLogInfo(B,"MatSetUpPreallocation: Warning not preallocating matrix storage");
191     ierr = (*B->ops->setuppreallocation)(B);CHKERRQ(ierr);
192     B->ops->setuppreallocation = 0;
193     B->preallocated            = PETSC_TRUE;
194   }
195   PetscFunctionReturn(0);
196 }
197 
198 /*
199         Copies from Cs header to A
200 */
201 #undef __FUNCT__
202 #define __FUNCT__ "MatHeaderCopy"
203 PetscErrorCode MatHeaderCopy(Mat A,Mat C)
204 {
205   PetscErrorCode ierr;
206   PetscInt       refct;
207   PetscOps       *Abops;
208   MatOps         Aops;
209   char           *mtype,*mname;
210   void           *spptr;
211 
212   PetscFunctionBegin;
213   /* free all the interior data structures from mat */
214   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
215 
216   ierr = PetscMapDestroy(A->rmap);CHKERRQ(ierr);
217   ierr = PetscMapDestroy(A->cmap);CHKERRQ(ierr);
218 
219   /* save the parts of A we need */
220   Abops = A->bops;
221   Aops  = A->ops;
222   refct = A->refct;
223   mtype = A->type_name;
224   mname = A->name;
225   spptr = A->spptr;
226 
227   if (C->spptr) {
228     ierr = PetscFree(C->spptr);CHKERRQ(ierr);
229     C->spptr = PETSC_NULL;
230   }
231 
232   /* copy C over to A */
233   ierr  = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
234 
235   /* return the parts of A we saved */
236   A->bops      = Abops;
237   A->ops       = Aops;
238   A->qlist     = 0;
239   A->refct     = refct;
240   A->type_name = mtype;
241   A->name      = mname;
242   A->spptr     = spptr;
243 
244   ierr = PetscHeaderDestroy(C);CHKERRQ(ierr);
245   PetscFunctionReturn(0);
246 }
247