xref: /petsc/src/mat/utils/gcreate.c (revision b0a32e0c6855ee6a6cd3495fa7da12ea9885bc5d)
1 /*$Id: gcreate.c,v 1.125 2000/11/28 17:29:54 bsmith Exp bsmith $*/
2 
3 #include "petscsys.h"
4 #include "src/mat/matimpl.h"       /*I "petscmat.h"  I*/
5 
6 #undef __FUNC__
7 #define __FUNC__ "MatPublish_Base"
8 static int MatPublish_Base(PetscObject obj)
9 {
10 #if defined(PETSC_HAVE_AMS)
11   Mat mat = (Mat)obj;
12   int ierr;
13 #endif
14 
15   PetscFunctionBegin;
16 #if defined(PETSC_HAVE_AMS)
17   /* if it is already published then return */
18   if (mat->amem >=0) PetscFunctionReturn(0);
19 
20   ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr);
21   ierr = AMS_Memory_add_field((AMS_Memory)mat->amem,"globalsizes",&mat->M,2,AMS_INT,AMS_READ,
22                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
23   ierr = AMS_Memory_add_field((AMS_Memory)mat->amem,"localsizes",&mat->m,2,AMS_INT,AMS_READ,
24                                 AMS_DISTRIBUTED,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
25   ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr);
26 #endif
27 
28   PetscFunctionReturn(0);
29 }
30 
31 
32 #undef __FUNC__
33 #define __FUNC__ "MatCreate"
34 /*@C
35    MatCreate - Creates a matrix where the type is determined
36    from the options database. Generates a parallel MPI matrix if the
37    communicator has more than one processor.  The default matrix type is
38    AIJ, using the routines MatCreateSeqAIJ() and MatCreateMPIAIJ().
39 
40    Collective on MPI_Comm
41 
42    Input Parameters:
43 +  m - number of local rows (or PETSC_DECIDE)
44 .  n - number of local columns (or PETSC_DECIDE)
45 .  M - number of global rows (or PETSC_DETERMINE)
46 .  N - number of global columns (or PETSC_DETERMINE)
47 -  comm - MPI communicator
48 
49    Output Parameter:
50 .  A - the matrix
51 
52    Options Database Keys:
53 +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
54 .    -mat_type mpiaij   - AIJ type, uses MatCreateMPIAIJ()
55 .    -mat_type seqbdiag - block diagonal type, uses MatCreateSeqBDiag()
56 .    -mat_type mpibdiag - block diagonal type, uses MatCreateMPIBDiag()
57 .    -mat_type mpirowbs - rowbs type, uses MatCreateMPIRowbs()
58 .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
59 .    -mat_type mpidense - dense type, uses MatCreateMPIDense()
60 .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
61 -    -mat_type mpibaij  - block AIJ type, uses MatCreateMPIBAIJ()
62 
63    Even More Options Database Keys:
64    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
65    for additional format-specific options.
66 
67    Notes:
68    If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the
69    user must ensure that they are chosen to be compatible with the
70    vectors. To do this, one first considers the matrix-vector product
71    'y = A x'. The 'm' that is used in the above routine must match the
72    local size used in the vector creation routine VecCreateMPI() for 'y'.
73    Likewise, the 'n' used must match that used as the local size in
74    VecCreateMPI() for 'x'.
75 
76    Level: beginner
77 
78 .keywords: matrix, create
79 
80 .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
81           MatCreateSeqBDiag(),MatCreateMPIBDiag(),
82           MatCreateSeqDense(), MatCreateMPIDense(),
83           MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
84           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
85           MatConvert()
86 @*/
87 int MatCreate(MPI_Comm comm,int m,int n,int M,int N,Mat *A)
88 {
89   Mat B;
90 
91   PetscFunctionBegin;
92   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,0,"Mat",comm,MatDestroy,MatView);
93   PetscLogObjectCreate(B);
94 
95   B->m = m;
96   B->n = n;
97   B->M = M;
98   B->N = N;
99 
100   B->preallocated  = PETSC_FALSE;
101   B->bops->publish = MatPublish_Base;
102   *A = B;
103   PetscFunctionReturn(0);
104 }
105 
106 #undef __FUNC__
107 #define __FUNC__ "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 MatSetFromOptionsSeqAIJ() and MatSetFromOptionsMPIAIJ().
113 
114    Collective on Mat
115 
116    Input Parameter:
117 .  A - the matrix
118 
119    Options Database Keys:
120 +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
121 .    -mat_type mpiaij   - AIJ type, uses MatCreateMPIAIJ()
122 .    -mat_type seqbdiag - block diagonal type, uses MatCreateSeqBDiag()
123 .    -mat_type mpibdiag - block diagonal type, uses MatCreateMPIBDiag()
124 .    -mat_type mpirowbs - rowbs type, uses MatCreateMPIRowbs()
125 .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
126 .    -mat_type mpidense - dense type, uses MatCreateMPIDense()
127 .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
128 -    -mat_type mpibaij  - block AIJ type, uses MatCreateMPIBAIJ()
129 
130    Even More Options Database Keys:
131    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
132    for additional format-specific options.
133 
134    Level: beginner
135 
136 .keywords: matrix, create
137 
138 .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
139           MatCreateSeqBDiag(),MatCreateMPIBDiag(),
140           MatCreateSeqDense(), MatCreateMPIDense(),
141           MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
142           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
143           MatConvert()
144 @*/
145 int MatSetFromOptions(Mat B)
146 {
147   int        ierr,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     if (size == 1) {
159       ierr = MatSetType(B,MATSEQAIJ);CHKERRQ(ierr);
160     } else {
161       ierr = MatSetType(B,MATMPIAIJ);CHKERRQ(ierr);
162     }
163   }
164   PetscFunctionReturn(0);
165 }
166 
167 #undef __FUNC__
168 #define __FUNC__ "MatSetUpPreallocation"
169 /*@C
170    MatSetUpPreallocation
171 
172    Collective on Mat
173 
174    Input Parameter:
175 .  A - the matrix
176 
177    Level: beginner
178 
179 .keywords: matrix, create
180 
181 .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
182           MatCreateSeqBDiag(),MatCreateMPIBDiag(),
183           MatCreateSeqDense(), MatCreateMPIDense(),
184           MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
185           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
186           MatConvert()
187 @*/
188 int MatSetUpPreallocation(Mat B)
189 {
190   int        ierr;
191 
192   PetscFunctionBegin;
193   if (B->ops->setuppreallocation) {
194     PetscLogInfo(B,"MatSetTpPreallocation: Warning not preallocating matrix storage");
195     ierr = (*B->ops->setuppreallocation)(B);CHKERRQ(ierr);
196     B->ops->setuppreallocation = 0;
197     B->preallocated            = PETSC_TRUE;
198   }
199   PetscFunctionReturn(0);
200 }
201 
202 /*
203         Copies from Cs header to A
204 */
205 #undef __FUNC__
206 #define __FUNC__ "MatHeaderCopy"
207 int MatHeaderCopy(Mat A,Mat C)
208 {
209   int         ierr,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 = MapDestroy(A->rmap);CHKERRQ(ierr);
219   ierr = MapDestroy(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   PetscHeaderDestroy(C);
240   PetscFunctionReturn(0);
241 }
242