xref: /petsc/src/mat/utils/gcreate.c (revision 0700a8246d308f50502909ba325e6169d3ee27eb)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
27807a1faSBarry Smith 
37c4f633dSBarry Smith #include "private/matimpl.h"       /*I "petscmat.h"  I*/
47807a1faSBarry Smith 
57adad957SLisandro Dalcin #if 0
64a2ae208SSatish Balay #undef __FUNCT__
74a2ae208SSatish Balay #define __FUNCT__ "MatPublish_Base"
86849ba73SBarry Smith static PetscErrorCode MatPublish_Base(PetscObject obj)
935d8aa7fSBarry Smith {
1035d8aa7fSBarry Smith   PetscFunctionBegin;
1135d8aa7fSBarry Smith   PetscFunctionReturn(0);
1235d8aa7fSBarry Smith }
137adad957SLisandro Dalcin #endif
1435d8aa7fSBarry Smith 
154a2ae208SSatish Balay #undef __FUNCT__
164a2ae208SSatish Balay #define __FUNCT__ "MatCreate"
1705869f15SSatish Balay /*@
1869dd0797SLois Curfman McInnes    MatCreate - Creates a matrix where the type is determined
197e5f4302SBarry Smith    from either a call to MatSetType() or from the options database
207e5f4302SBarry Smith    with a call to MatSetFromOptions(). The default matrix type is
217e5f4302SBarry Smith    AIJ, using the routines MatCreateSeqAIJ() or MatCreateMPIAIJ()
227e5f4302SBarry Smith    if you do not set a type in the options database. If you never
237e5f4302SBarry Smith    call MatSetType() or MatSetFromOptions() it will generate an
24f8ab6608SSatish Balay    error when you try to use the matrix.
2583e1b59cSLois Curfman McInnes 
26cb13003dSBarry Smith    Collective on MPI_Comm
27cb13003dSBarry Smith 
28f69a0ea3SMatthew Knepley    Input Parameter:
29f69a0ea3SMatthew Knepley .  comm - MPI communicator
307807a1faSBarry Smith 
317807a1faSBarry Smith    Output Parameter:
32dc401e71SLois Curfman McInnes .  A - the matrix
33e0b365e2SLois Curfman McInnes 
34273d9f13SBarry Smith    Options Database Keys:
35273d9f13SBarry Smith +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
36273d9f13SBarry Smith .    -mat_type mpiaij   - AIJ type, uses MatCreateMPIAIJ()
37273d9f13SBarry Smith .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
38273d9f13SBarry Smith .    -mat_type mpidense - dense type, uses MatCreateMPIDense()
39273d9f13SBarry Smith .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
40273d9f13SBarry Smith -    -mat_type mpibaij  - block AIJ type, uses MatCreateMPIBAIJ()
41e0b365e2SLois Curfman McInnes 
4283e1b59cSLois Curfman McInnes    Even More Options Database Keys:
4383e1b59cSLois Curfman McInnes    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
4483e1b59cSLois Curfman McInnes    for additional format-specific options.
45e0b365e2SLois Curfman McInnes 
46bd9ce289SLois Curfman McInnes    Notes:
47ec6e0d80SSatish Balay 
48273d9f13SBarry Smith    Level: beginner
49273d9f13SBarry Smith 
50a2fc510eSBarry Smith    User manual sections:
51a2fc510eSBarry Smith +   sec_matcreate
52a2fc510eSBarry Smith -   chapter_matrices
53a2fc510eSBarry Smith 
54273d9f13SBarry Smith .keywords: matrix, create
55273d9f13SBarry Smith 
56cd05a4c0SHong Zhang .seealso: MatCreateSeqAIJ(), MatCreateMPIAIJ(),
57273d9f13SBarry Smith           MatCreateSeqDense(), MatCreateMPIDense(),
586531c467SSatish Balay           MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
59273d9f13SBarry Smith           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
60273d9f13SBarry Smith           MatConvert()
61273d9f13SBarry Smith @*/
62f69a0ea3SMatthew Knepley PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm comm,Mat *A)
63273d9f13SBarry Smith {
64273d9f13SBarry Smith   Mat            B;
65dfbe8321SBarry Smith   PetscErrorCode ierr;
66273d9f13SBarry Smith 
67273d9f13SBarry Smith   PetscFunctionBegin;
68f69a0ea3SMatthew Knepley   PetscValidPointer(A,2);
6997f1f81fSBarry Smith 
708ba1e511SMatthew Knepley   *A = PETSC_NULL;
718ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES
728ba1e511SMatthew Knepley   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
738ba1e511SMatthew Knepley #endif
748ba1e511SMatthew Knepley 
75*0700a824SBarry Smith   ierr = PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_CLASSID,0,"Mat",comm,MatDestroy,MatView);CHKERRQ(ierr);
7626283091SBarry Smith   ierr = PetscLayoutCreate(comm,&B->rmap);CHKERRQ(ierr);
7726283091SBarry Smith   ierr = PetscLayoutCreate(comm,&B->cmap);CHKERRQ(ierr);
78273d9f13SBarry Smith   B->preallocated  = PETSC_FALSE;
79273d9f13SBarry Smith   *A               = B;
80273d9f13SBarry Smith   PetscFunctionReturn(0);
81273d9f13SBarry Smith }
82273d9f13SBarry Smith 
834a2ae208SSatish Balay #undef __FUNCT__
84f69a0ea3SMatthew Knepley #define __FUNCT__ "MatSetSizes"
85f69a0ea3SMatthew Knepley /*@
86f69a0ea3SMatthew Knepley   MatSetSizes - Sets the local and global sizes, and checks to determine compatibility
87f69a0ea3SMatthew Knepley 
88f69a0ea3SMatthew Knepley   Collective on Mat
89f69a0ea3SMatthew Knepley 
90f69a0ea3SMatthew Knepley   Input Parameters:
91f69a0ea3SMatthew Knepley +  A - the matrix
92f69a0ea3SMatthew Knepley .  m - number of local rows (or PETSC_DECIDE)
93f69a0ea3SMatthew Knepley .  n - number of local columns (or PETSC_DECIDE)
94f69a0ea3SMatthew Knepley .  M - number of global rows (or PETSC_DETERMINE)
95f69a0ea3SMatthew Knepley -  N - number of global columns (or PETSC_DETERMINE)
96f69a0ea3SMatthew Knepley 
97f69a0ea3SMatthew Knepley    Notes:
98f69a0ea3SMatthew Knepley    m (n) and M (N) cannot be both PETSC_DECIDE
99f69a0ea3SMatthew Knepley    If one processor calls this with M (N) of PETSC_DECIDE then all processors must, otherwise the program will hang.
100f69a0ea3SMatthew Knepley 
101f69a0ea3SMatthew Knepley    If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the
102f69a0ea3SMatthew Knepley    user must ensure that they are chosen to be compatible with the
103f69a0ea3SMatthew Knepley    vectors. To do this, one first considers the matrix-vector product
104f69a0ea3SMatthew Knepley    'y = A x'. The 'm' that is used in the above routine must match the
105f69a0ea3SMatthew Knepley    local size used in the vector creation routine VecCreateMPI() for 'y'.
106f69a0ea3SMatthew Knepley    Likewise, the 'n' used must match that used as the local size in
107f69a0ea3SMatthew Knepley    VecCreateMPI() for 'x'.
108f69a0ea3SMatthew Knepley 
109f69a0ea3SMatthew Knepley   Level: beginner
110f69a0ea3SMatthew Knepley 
111f69a0ea3SMatthew Knepley .seealso: MatGetSize(), PetscSplitOwnership()
112f69a0ea3SMatthew Knepley @*/
113f69a0ea3SMatthew Knepley PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N)
114f69a0ea3SMatthew Knepley {
115284134d9SBarry Smith   PetscErrorCode ierr;
116284134d9SBarry Smith 
117f69a0ea3SMatthew Knepley   PetscFunctionBegin;
118*0700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
119f69a0ea3SMatthew Knepley   if (M > 0 && m > M) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Local column size %D cannot be larger than global column size %D",m,M);
120f69a0ea3SMatthew Knepley   if (N > 0 && n > N) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Local row size %D cannot be larger than global row size %D",n,N);
121284134d9SBarry Smith   if (A->ops->setsizes) {
122284134d9SBarry Smith     /* Since this will not be set until the type has been set, this will NOT be called on the initial
123284134d9SBarry Smith        call of MatSetSizes() (which must be called BEFORE MatSetType() */
124284134d9SBarry Smith     ierr = (*A->ops->setsizes)(A,m,n,M,N);CHKERRQ(ierr);
125284134d9SBarry Smith   } else {
126d0f46423SBarry Smith     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);
127d0f46423SBarry Smith     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);
128284134d9SBarry Smith   }
129d0f46423SBarry Smith   A->rmap->n = m;
130d0f46423SBarry Smith   A->cmap->n = n;
131d0f46423SBarry Smith   A->rmap->N = M;
132d0f46423SBarry Smith   A->cmap->N = N;
133117016b1SBarry Smith   if (A->ops->create) {
134117016b1SBarry Smith     ierr = (*A->ops->create)(A);CHKERRQ(ierr);
135117016b1SBarry Smith     A->ops->create = 0;
136117016b1SBarry Smith   }
137117016b1SBarry Smith 
138f69a0ea3SMatthew Knepley   PetscFunctionReturn(0);
139f69a0ea3SMatthew Knepley }
140f69a0ea3SMatthew Knepley 
141f69a0ea3SMatthew Knepley #undef __FUNCT__
1424a2ae208SSatish Balay #define __FUNCT__ "MatSetFromOptions"
14305869f15SSatish Balay /*@
144273d9f13SBarry Smith    MatSetFromOptions - Creates a matrix where the type is determined
145273d9f13SBarry Smith    from the options database. Generates a parallel MPI matrix if the
146273d9f13SBarry Smith    communicator has more than one processor.  The default matrix type is
1477e5f4302SBarry Smith    AIJ, using the routines MatCreateSeqAIJ() and MatCreateMPIAIJ() if
1487e5f4302SBarry Smith    you do not select a type in the options database.
149273d9f13SBarry Smith 
150273d9f13SBarry Smith    Collective on Mat
151273d9f13SBarry Smith 
152273d9f13SBarry Smith    Input Parameter:
153273d9f13SBarry Smith .  A - the matrix
154273d9f13SBarry Smith 
155273d9f13SBarry Smith    Options Database Keys:
156273d9f13SBarry Smith +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
157273d9f13SBarry Smith .    -mat_type mpiaij   - AIJ type, uses MatCreateMPIAIJ()
158273d9f13SBarry Smith .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
159273d9f13SBarry Smith .    -mat_type mpidense - dense type, uses MatCreateMPIDense()
160273d9f13SBarry Smith .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
161273d9f13SBarry Smith -    -mat_type mpibaij  - block AIJ type, uses MatCreateMPIBAIJ()
162273d9f13SBarry Smith 
163273d9f13SBarry Smith    Even More Options Database Keys:
164273d9f13SBarry Smith    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
165273d9f13SBarry Smith    for additional format-specific options.
166bd9ce289SLois Curfman McInnes 
1671d69843bSLois Curfman McInnes    Level: beginner
1681d69843bSLois Curfman McInnes 
169dc401e71SLois Curfman McInnes .keywords: matrix, create
170e0b365e2SLois Curfman McInnes 
171fafbff53SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
17239ddd567SLois Curfman McInnes           MatCreateSeqDense(), MatCreateMPIDense(),
1736531c467SSatish Balay           MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
174a209d233SLois Curfman McInnes           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
175273d9f13SBarry Smith           MatConvert()
1767807a1faSBarry Smith @*/
177be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat B)
1787807a1faSBarry Smith {
179dfbe8321SBarry Smith   PetscErrorCode ierr;
180f3be49caSLisandro Dalcin   const char     *deft = MATAIJ;
181f3be49caSLisandro Dalcin   char           type[256];
182273d9f13SBarry Smith   PetscTruth     flg;
183dbb450caSBarry Smith 
1843a40ed3dSBarry Smith   PetscFunctionBegin;
185*0700a824SBarry Smith   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
186f3be49caSLisandro Dalcin 
187f3be49caSLisandro Dalcin   ierr = PetscOptionsBegin(((PetscObject)B)->comm,((PetscObject)B)->prefix,"Matrix options","Mat");CHKERRQ(ierr);
188f3be49caSLisandro Dalcin     ierr = PetscOptionsList("-mat_type","Matrix type","MatSetType",MatList,deft,type,256,&flg);CHKERRQ(ierr);
189273d9f13SBarry Smith     if (flg) {
190f3be49caSLisandro Dalcin       ierr = MatSetType(B,type);CHKERRQ(ierr);
191f3be49caSLisandro Dalcin     } else if (!((PetscObject)B)->type_name) {
192f3be49caSLisandro Dalcin       ierr = MatSetType(B,deft);CHKERRQ(ierr);
193273d9f13SBarry Smith     }
194f3be49caSLisandro Dalcin 
195f3be49caSLisandro Dalcin     if (B->ops->setfromoptions) {
196f3be49caSLisandro Dalcin       ierr = (*B->ops->setfromoptions)(B);CHKERRQ(ierr);
197593c1905SSatish Balay     }
198f3be49caSLisandro Dalcin 
199f3be49caSLisandro Dalcin   ierr = PetscOptionsEnd();CHKERRQ(ierr);
200f3be49caSLisandro Dalcin 
2013a40ed3dSBarry Smith   PetscFunctionReturn(0);
2027807a1faSBarry Smith }
2037807a1faSBarry Smith 
2044a2ae208SSatish Balay #undef __FUNCT__
2054a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation"
206bc08b0f1SBarry Smith /*@
207273d9f13SBarry Smith    MatSetUpPreallocation
208dae03382SLois Curfman McInnes 
209273d9f13SBarry Smith    Collective on Mat
210dae03382SLois Curfman McInnes 
211273d9f13SBarry Smith    Input Parameter:
212273d9f13SBarry Smith .  A - the matrix
213d5d45c9bSBarry Smith 
214273d9f13SBarry Smith    Level: beginner
215d5d45c9bSBarry Smith 
216273d9f13SBarry Smith .keywords: matrix, create
217273d9f13SBarry Smith 
218273d9f13SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
219273d9f13SBarry Smith           MatCreateSeqDense(), MatCreateMPIDense(),
2206531c467SSatish Balay           MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
221273d9f13SBarry Smith           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
222273d9f13SBarry Smith           MatConvert()
223273d9f13SBarry Smith @*/
224be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat B)
225273d9f13SBarry Smith {
226dfbe8321SBarry Smith   PetscErrorCode ierr;
227273d9f13SBarry Smith 
228273d9f13SBarry Smith   PetscFunctionBegin;
22917667f90SBarry Smith   if (!B->preallocated && B->ops->setuppreallocation) {
230ae15b995SBarry Smith     ierr = PetscInfo(B,"Warning not preallocating matrix storage\n");CHKERRQ(ierr);
231273d9f13SBarry Smith     ierr = (*B->ops->setuppreallocation)(B);CHKERRQ(ierr);
232273d9f13SBarry Smith   }
233210f0121SBarry Smith   B->preallocated = PETSC_TRUE;
234273d9f13SBarry Smith   PetscFunctionReturn(0);
235273d9f13SBarry Smith }
236273d9f13SBarry Smith 
237273d9f13SBarry Smith /*
238273d9f13SBarry Smith         Copies from Cs header to A
239d0f46423SBarry Smith 
240d0f46423SBarry Smith         This is somewhat different from MatHeaderReplace() it would be nice to merge the code
241273d9f13SBarry Smith */
2424a2ae208SSatish Balay #undef __FUNCT__
2434a2ae208SSatish Balay #define __FUNCT__ "MatHeaderCopy"
244dfbe8321SBarry Smith PetscErrorCode MatHeaderCopy(Mat A,Mat C)
245273d9f13SBarry Smith {
246dfbe8321SBarry Smith   PetscErrorCode ierr;
247d44834fbSBarry Smith   PetscInt       refct;
248273d9f13SBarry Smith   PetscOps       *Abops;
249273d9f13SBarry Smith   MatOps         Aops;
250273d9f13SBarry Smith   char           *mtype,*mname;
25130735b05SKris Buschelman   void           *spptr;
252273d9f13SBarry Smith 
253273d9f13SBarry Smith   PetscFunctionBegin;
254273d9f13SBarry Smith   /* save the parts of A we need */
2557adad957SLisandro Dalcin   Abops = ((PetscObject)A)->bops;
256273d9f13SBarry Smith   Aops  = A->ops;
2577adad957SLisandro Dalcin   refct = ((PetscObject)A)->refct;
2585c9eb25fSBarry Smith   mtype = ((PetscObject)A)->type_name;
2595c9eb25fSBarry Smith   mname = ((PetscObject)A)->name;
26030735b05SKris Buschelman   spptr = A->spptr;
26130735b05SKris Buschelman 
2625c9eb25fSBarry Smith   /* zero these so the destroy below does not free them */
2635c9eb25fSBarry Smith   ((PetscObject)A)->type_name = 0;
2645c9eb25fSBarry Smith   ((PetscObject)A)->name      = 0;
2655c9eb25fSBarry Smith 
2667c99f97cSSatish Balay   /* free all the interior data structures from mat */
2677c99f97cSSatish Balay   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
2687c99f97cSSatish Balay 
26930735b05SKris Buschelman   ierr = PetscFree(C->spptr);CHKERRQ(ierr);
27005b42c5fSBarry Smith 
27126283091SBarry Smith   ierr = PetscLayoutDestroy(A->rmap);CHKERRQ(ierr);
27226283091SBarry Smith   ierr = PetscLayoutDestroy(A->cmap);CHKERRQ(ierr);
2735c9eb25fSBarry Smith   ierr = PetscFListDestroy(&((PetscObject)A)->qlist);CHKERRQ(ierr);
2745c9eb25fSBarry Smith   ierr = PetscOListDestroy(((PetscObject)A)->olist);CHKERRQ(ierr);
275273d9f13SBarry Smith 
276273d9f13SBarry Smith   /* copy C over to A */
277273d9f13SBarry Smith   ierr  = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
278273d9f13SBarry Smith 
279273d9f13SBarry Smith   /* return the parts of A we saved */
2807adad957SLisandro Dalcin   ((PetscObject)A)->bops      = Abops;
281273d9f13SBarry Smith   A->ops                      = Aops;
2827adad957SLisandro Dalcin   ((PetscObject)A)->refct     = refct;
2837adad957SLisandro Dalcin   ((PetscObject)A)->type_name = mtype;
2847adad957SLisandro Dalcin   ((PetscObject)A)->name      = mname;
28530735b05SKris Buschelman   A->spptr                    = spptr;
286273d9f13SBarry Smith 
2875c9eb25fSBarry Smith   /* since these two are copied into A we do not want them destroyed in C */
2885c9eb25fSBarry Smith   ((PetscObject)C)->qlist = 0;
2895c9eb25fSBarry Smith   ((PetscObject)C)->olist = 0;
290d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(C);CHKERRQ(ierr);
291273d9f13SBarry Smith   PetscFunctionReturn(0);
292273d9f13SBarry Smith }
2938ab5b326SKris Buschelman /*
2948ab5b326SKris Buschelman         Replace A's header with that of C
2958ab5b326SKris Buschelman         This is essentially code moved from MatDestroy
296d0f46423SBarry Smith 
297d0f46423SBarry Smith         This is somewhat different from MatHeaderCopy() it would be nice to merge the code
2988ab5b326SKris Buschelman */
2998ab5b326SKris Buschelman #undef __FUNCT__
3008ab5b326SKris Buschelman #define __FUNCT__ "MatHeaderReplace"
3018ab5b326SKris Buschelman PetscErrorCode MatHeaderReplace(Mat A,Mat C)
3028ab5b326SKris Buschelman {
3038ab5b326SKris Buschelman   PetscErrorCode ierr;
3048ab5b326SKris Buschelman 
3058ab5b326SKris Buschelman   PetscFunctionBegin;
3066d7c1e57SBarry Smith   if (A == C) PetscFunctionReturn(0);
3076d7c1e57SBarry Smith 
3088ab5b326SKris Buschelman   /* free all the interior data structures from mat */
3098ab5b326SKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
310d2c95936SBarry Smith   ierr = PetscHeaderDestroy_Private((PetscObject)A);CHKERRQ(ierr);
311c650bb5fSLisandro Dalcin   ierr = PetscFree(A->ops);CHKERRQ(ierr);
31226283091SBarry Smith   ierr = PetscLayoutDestroy(A->rmap);CHKERRQ(ierr);
31326283091SBarry Smith   ierr = PetscLayoutDestroy(A->cmap);CHKERRQ(ierr);
314d95db058SSatish Balay   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
3158ab5b326SKris Buschelman 
3168ab5b326SKris Buschelman   /* copy C over to A */
3178ab5b326SKris Buschelman   if (C) {
3188ab5b326SKris Buschelman     ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
3198ab5b326SKris Buschelman     ierr = PetscLogObjectDestroy((PetscObject)C);CHKERRQ(ierr);
3208ab5b326SKris Buschelman     ierr = PetscFree(C);CHKERRQ(ierr);
3218ab5b326SKris Buschelman   }
3228ab5b326SKris Buschelman   PetscFunctionReturn(0);
3238ab5b326SKris Buschelman }
324