xref: /petsc/src/mat/utils/gcreate.c (revision 6531c46765a6e92ffcdfc7a2f84e41b5f3cd76f6)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
27807a1faSBarry Smith 
37c4f633dSBarry Smith #include "private/matimpl.h"       /*I "petscmat.h"  I*/
4325e03aeSBarry Smith #include "petscsys.h"
57807a1faSBarry Smith 
67adad957SLisandro Dalcin #if 0
74a2ae208SSatish Balay #undef __FUNCT__
84a2ae208SSatish Balay #define __FUNCT__ "MatPublish_Base"
96849ba73SBarry Smith static PetscErrorCode MatPublish_Base(PetscObject obj)
1035d8aa7fSBarry Smith {
1135d8aa7fSBarry Smith   PetscFunctionBegin;
1235d8aa7fSBarry Smith   PetscFunctionReturn(0);
1335d8aa7fSBarry Smith }
147adad957SLisandro Dalcin #endif
1535d8aa7fSBarry Smith 
164a2ae208SSatish Balay #undef __FUNCT__
174a2ae208SSatish Balay #define __FUNCT__ "MatCreate"
1805869f15SSatish Balay /*@
1969dd0797SLois Curfman McInnes    MatCreate - Creates a matrix where the type is determined
207e5f4302SBarry Smith    from either a call to MatSetType() or from the options database
217e5f4302SBarry Smith    with a call to MatSetFromOptions(). The default matrix type is
227e5f4302SBarry Smith    AIJ, using the routines MatCreateSeqAIJ() or MatCreateMPIAIJ()
237e5f4302SBarry Smith    if you do not set a type in the options database. If you never
247e5f4302SBarry Smith    call MatSetType() or MatSetFromOptions() it will generate an
25f8ab6608SSatish Balay    error when you try to use the matrix.
2683e1b59cSLois Curfman McInnes 
27cb13003dSBarry Smith    Collective on MPI_Comm
28cb13003dSBarry Smith 
29f69a0ea3SMatthew Knepley    Input Parameter:
30f69a0ea3SMatthew Knepley .  comm - MPI communicator
317807a1faSBarry Smith 
327807a1faSBarry Smith    Output Parameter:
33dc401e71SLois Curfman McInnes .  A - the matrix
34e0b365e2SLois Curfman McInnes 
35273d9f13SBarry Smith    Options Database Keys:
36273d9f13SBarry Smith +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
37273d9f13SBarry Smith .    -mat_type mpiaij   - AIJ type, uses MatCreateMPIAIJ()
38273d9f13SBarry Smith .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
39273d9f13SBarry Smith .    -mat_type mpidense - dense type, uses MatCreateMPIDense()
40273d9f13SBarry Smith .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
41273d9f13SBarry Smith -    -mat_type mpibaij  - block AIJ type, uses MatCreateMPIBAIJ()
42e0b365e2SLois Curfman McInnes 
4383e1b59cSLois Curfman McInnes    Even More Options Database Keys:
4483e1b59cSLois Curfman McInnes    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
4583e1b59cSLois Curfman McInnes    for additional format-specific options.
46e0b365e2SLois Curfman McInnes 
47bd9ce289SLois Curfman McInnes    Notes:
48ec6e0d80SSatish Balay 
49273d9f13SBarry Smith    Level: beginner
50273d9f13SBarry Smith 
51a2fc510eSBarry Smith    User manual sections:
52a2fc510eSBarry Smith +   sec_matcreate
53a2fc510eSBarry Smith -   chapter_matrices
54a2fc510eSBarry Smith 
55273d9f13SBarry Smith .keywords: matrix, create
56273d9f13SBarry Smith 
57cd05a4c0SHong Zhang .seealso: MatCreateSeqAIJ(), MatCreateMPIAIJ(),
58273d9f13SBarry Smith           MatCreateSeqDense(), MatCreateMPIDense(),
59*6531c467SSatish Balay           MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
60273d9f13SBarry Smith           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
61273d9f13SBarry Smith           MatConvert()
62273d9f13SBarry Smith @*/
63f69a0ea3SMatthew Knepley PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm comm,Mat *A)
64273d9f13SBarry Smith {
65273d9f13SBarry Smith   Mat            B;
66dfbe8321SBarry Smith   PetscErrorCode ierr;
67273d9f13SBarry Smith 
68273d9f13SBarry Smith   PetscFunctionBegin;
69f69a0ea3SMatthew Knepley   PetscValidPointer(A,2);
7097f1f81fSBarry Smith 
718ba1e511SMatthew Knepley   *A = PETSC_NULL;
728ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES
738ba1e511SMatthew Knepley   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
748ba1e511SMatthew Knepley #endif
758ba1e511SMatthew Knepley 
7652e6d16bSBarry Smith   ierr = PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,0,"Mat",comm,MatDestroy,MatView);CHKERRQ(ierr);
7726283091SBarry Smith   ierr = PetscLayoutCreate(comm,&B->rmap);CHKERRQ(ierr);
7826283091SBarry Smith   ierr = PetscLayoutCreate(comm,&B->cmap);CHKERRQ(ierr);
79273d9f13SBarry Smith   B->preallocated  = PETSC_FALSE;
80273d9f13SBarry Smith   *A               = B;
81273d9f13SBarry Smith   PetscFunctionReturn(0);
82273d9f13SBarry Smith }
83273d9f13SBarry Smith 
844a2ae208SSatish Balay #undef __FUNCT__
85f69a0ea3SMatthew Knepley #define __FUNCT__ "MatSetSizes"
86f69a0ea3SMatthew Knepley /*@
87f69a0ea3SMatthew Knepley   MatSetSizes - Sets the local and global sizes, and checks to determine compatibility
88f69a0ea3SMatthew Knepley 
89f69a0ea3SMatthew Knepley   Collective on Mat
90f69a0ea3SMatthew Knepley 
91f69a0ea3SMatthew Knepley   Input Parameters:
92f69a0ea3SMatthew Knepley +  A - the matrix
93f69a0ea3SMatthew Knepley .  m - number of local rows (or PETSC_DECIDE)
94f69a0ea3SMatthew Knepley .  n - number of local columns (or PETSC_DECIDE)
95f69a0ea3SMatthew Knepley .  M - number of global rows (or PETSC_DETERMINE)
96f69a0ea3SMatthew Knepley -  N - number of global columns (or PETSC_DETERMINE)
97f69a0ea3SMatthew Knepley 
98f69a0ea3SMatthew Knepley    Notes:
99f69a0ea3SMatthew Knepley    m (n) and M (N) cannot be both PETSC_DECIDE
100f69a0ea3SMatthew Knepley    If one processor calls this with M (N) of PETSC_DECIDE then all processors must, otherwise the program will hang.
101f69a0ea3SMatthew Knepley 
102f69a0ea3SMatthew Knepley    If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the
103f69a0ea3SMatthew Knepley    user must ensure that they are chosen to be compatible with the
104f69a0ea3SMatthew Knepley    vectors. To do this, one first considers the matrix-vector product
105f69a0ea3SMatthew Knepley    'y = A x'. The 'm' that is used in the above routine must match the
106f69a0ea3SMatthew Knepley    local size used in the vector creation routine VecCreateMPI() for 'y'.
107f69a0ea3SMatthew Knepley    Likewise, the 'n' used must match that used as the local size in
108f69a0ea3SMatthew Knepley    VecCreateMPI() for 'x'.
109f69a0ea3SMatthew Knepley 
110f69a0ea3SMatthew Knepley   Level: beginner
111f69a0ea3SMatthew Knepley 
112f69a0ea3SMatthew Knepley .seealso: MatGetSize(), PetscSplitOwnership()
113f69a0ea3SMatthew Knepley @*/
114f69a0ea3SMatthew Knepley PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N)
115f69a0ea3SMatthew Knepley {
116284134d9SBarry Smith   PetscErrorCode ierr;
117284134d9SBarry Smith 
118f69a0ea3SMatthew Knepley   PetscFunctionBegin;
119f69a0ea3SMatthew Knepley   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
120f69a0ea3SMatthew 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);
121f69a0ea3SMatthew 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);
122284134d9SBarry Smith   if (A->ops->setsizes) {
123284134d9SBarry Smith     /* Since this will not be set until the type has been set, this will NOT be called on the initial
124284134d9SBarry Smith        call of MatSetSizes() (which must be called BEFORE MatSetType() */
125284134d9SBarry Smith     ierr = (*A->ops->setsizes)(A,m,n,M,N);CHKERRQ(ierr);
126284134d9SBarry Smith   } else {
127d0f46423SBarry 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);
128d0f46423SBarry 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);
129284134d9SBarry Smith   }
130d0f46423SBarry Smith   A->rmap->n = m;
131d0f46423SBarry Smith   A->cmap->n = n;
132d0f46423SBarry Smith   A->rmap->N = M;
133d0f46423SBarry Smith   A->cmap->N = N;
134117016b1SBarry Smith   if (A->ops->create) {
135117016b1SBarry Smith     ierr = (*A->ops->create)(A);CHKERRQ(ierr);
136117016b1SBarry Smith     A->ops->create = 0;
137117016b1SBarry Smith   }
138117016b1SBarry Smith 
139f69a0ea3SMatthew Knepley   PetscFunctionReturn(0);
140f69a0ea3SMatthew Knepley }
141f69a0ea3SMatthew Knepley 
142f69a0ea3SMatthew Knepley #undef __FUNCT__
1434a2ae208SSatish Balay #define __FUNCT__ "MatSetFromOptions"
14405869f15SSatish Balay /*@
145273d9f13SBarry Smith    MatSetFromOptions - Creates a matrix where the type is determined
146273d9f13SBarry Smith    from the options database. Generates a parallel MPI matrix if the
147273d9f13SBarry Smith    communicator has more than one processor.  The default matrix type is
1487e5f4302SBarry Smith    AIJ, using the routines MatCreateSeqAIJ() and MatCreateMPIAIJ() if
1497e5f4302SBarry Smith    you do not select a type in the options database.
150273d9f13SBarry Smith 
151273d9f13SBarry Smith    Collective on Mat
152273d9f13SBarry Smith 
153273d9f13SBarry Smith    Input Parameter:
154273d9f13SBarry Smith .  A - the matrix
155273d9f13SBarry Smith 
156273d9f13SBarry Smith    Options Database Keys:
157273d9f13SBarry Smith +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
158273d9f13SBarry Smith .    -mat_type mpiaij   - AIJ type, uses MatCreateMPIAIJ()
159273d9f13SBarry Smith .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
160273d9f13SBarry Smith .    -mat_type mpidense - dense type, uses MatCreateMPIDense()
161273d9f13SBarry Smith .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
162273d9f13SBarry Smith -    -mat_type mpibaij  - block AIJ type, uses MatCreateMPIBAIJ()
163273d9f13SBarry Smith 
164273d9f13SBarry Smith    Even More Options Database Keys:
165273d9f13SBarry Smith    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
166273d9f13SBarry Smith    for additional format-specific options.
167bd9ce289SLois Curfman McInnes 
1681d69843bSLois Curfman McInnes    Level: beginner
1691d69843bSLois Curfman McInnes 
170dc401e71SLois Curfman McInnes .keywords: matrix, create
171e0b365e2SLois Curfman McInnes 
172fafbff53SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
17339ddd567SLois Curfman McInnes           MatCreateSeqDense(), MatCreateMPIDense(),
174*6531c467SSatish Balay           MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
175a209d233SLois Curfman McInnes           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
176273d9f13SBarry Smith           MatConvert()
1777807a1faSBarry Smith @*/
178be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat B)
1797807a1faSBarry Smith {
180dfbe8321SBarry Smith   PetscErrorCode ierr;
181f3be49caSLisandro Dalcin   const char     *deft = MATAIJ;
182f3be49caSLisandro Dalcin   char           type[256];
183273d9f13SBarry Smith   PetscTruth     flg;
184dbb450caSBarry Smith 
1853a40ed3dSBarry Smith   PetscFunctionBegin;
186f3be49caSLisandro Dalcin   PetscValidHeaderSpecific(B,MAT_COOKIE,1);
187f3be49caSLisandro Dalcin 
188f3be49caSLisandro Dalcin   ierr = PetscOptionsBegin(((PetscObject)B)->comm,((PetscObject)B)->prefix,"Matrix options","Mat");CHKERRQ(ierr);
189f3be49caSLisandro Dalcin     ierr = PetscOptionsList("-mat_type","Matrix type","MatSetType",MatList,deft,type,256,&flg);CHKERRQ(ierr);
190273d9f13SBarry Smith     if (flg) {
191f3be49caSLisandro Dalcin       ierr = MatSetType(B,type);CHKERRQ(ierr);
192f3be49caSLisandro Dalcin     } else if (!((PetscObject)B)->type_name) {
193f3be49caSLisandro Dalcin       ierr = MatSetType(B,deft);CHKERRQ(ierr);
194273d9f13SBarry Smith     }
195f3be49caSLisandro Dalcin 
196f3be49caSLisandro Dalcin     if (B->ops->setfromoptions) {
197f3be49caSLisandro Dalcin       ierr = (*B->ops->setfromoptions)(B);CHKERRQ(ierr);
198593c1905SSatish Balay     }
199f3be49caSLisandro Dalcin 
200f3be49caSLisandro Dalcin   ierr = PetscOptionsEnd();CHKERRQ(ierr);
201f3be49caSLisandro Dalcin 
2023a40ed3dSBarry Smith   PetscFunctionReturn(0);
2037807a1faSBarry Smith }
2047807a1faSBarry Smith 
2054a2ae208SSatish Balay #undef __FUNCT__
2064a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation"
207bc08b0f1SBarry Smith /*@
208273d9f13SBarry Smith    MatSetUpPreallocation
209dae03382SLois Curfman McInnes 
210273d9f13SBarry Smith    Collective on Mat
211dae03382SLois Curfman McInnes 
212273d9f13SBarry Smith    Input Parameter:
213273d9f13SBarry Smith .  A - the matrix
214d5d45c9bSBarry Smith 
215273d9f13SBarry Smith    Level: beginner
216d5d45c9bSBarry Smith 
217273d9f13SBarry Smith .keywords: matrix, create
218273d9f13SBarry Smith 
219273d9f13SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
220273d9f13SBarry Smith           MatCreateSeqDense(), MatCreateMPIDense(),
221*6531c467SSatish Balay           MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
222273d9f13SBarry Smith           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
223273d9f13SBarry Smith           MatConvert()
224273d9f13SBarry Smith @*/
225be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat B)
226273d9f13SBarry Smith {
227dfbe8321SBarry Smith   PetscErrorCode ierr;
228273d9f13SBarry Smith 
229273d9f13SBarry Smith   PetscFunctionBegin;
23017667f90SBarry Smith   if (!B->preallocated && B->ops->setuppreallocation) {
231ae15b995SBarry Smith     ierr = PetscInfo(B,"Warning not preallocating matrix storage\n");CHKERRQ(ierr);
232273d9f13SBarry Smith     ierr = (*B->ops->setuppreallocation)(B);CHKERRQ(ierr);
233273d9f13SBarry Smith   }
234210f0121SBarry Smith   B->preallocated = PETSC_TRUE;
235273d9f13SBarry Smith   PetscFunctionReturn(0);
236273d9f13SBarry Smith }
237273d9f13SBarry Smith 
238273d9f13SBarry Smith /*
239273d9f13SBarry Smith         Copies from Cs header to A
240d0f46423SBarry Smith 
241d0f46423SBarry Smith         This is somewhat different from MatHeaderReplace() it would be nice to merge the code
242273d9f13SBarry Smith */
2434a2ae208SSatish Balay #undef __FUNCT__
2444a2ae208SSatish Balay #define __FUNCT__ "MatHeaderCopy"
245dfbe8321SBarry Smith PetscErrorCode MatHeaderCopy(Mat A,Mat C)
246273d9f13SBarry Smith {
247dfbe8321SBarry Smith   PetscErrorCode ierr;
248d44834fbSBarry Smith   PetscInt       refct;
249273d9f13SBarry Smith   PetscOps       *Abops;
250273d9f13SBarry Smith   MatOps         Aops;
251273d9f13SBarry Smith   char           *mtype,*mname;
25230735b05SKris Buschelman   void           *spptr;
253273d9f13SBarry Smith 
254273d9f13SBarry Smith   PetscFunctionBegin;
255273d9f13SBarry Smith   /* save the parts of A we need */
2567adad957SLisandro Dalcin   Abops = ((PetscObject)A)->bops;
257273d9f13SBarry Smith   Aops  = A->ops;
2587adad957SLisandro Dalcin   refct = ((PetscObject)A)->refct;
2595c9eb25fSBarry Smith   mtype = ((PetscObject)A)->type_name;
2605c9eb25fSBarry Smith   mname = ((PetscObject)A)->name;
26130735b05SKris Buschelman   spptr = A->spptr;
26230735b05SKris Buschelman 
2635c9eb25fSBarry Smith   /* zero these so the destroy below does not free them */
2645c9eb25fSBarry Smith   ((PetscObject)A)->type_name = 0;
2655c9eb25fSBarry Smith   ((PetscObject)A)->name      = 0;
2665c9eb25fSBarry Smith 
2677c99f97cSSatish Balay   /* free all the interior data structures from mat */
2687c99f97cSSatish Balay   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
2697c99f97cSSatish Balay 
27030735b05SKris Buschelman   ierr = PetscFree(C->spptr);CHKERRQ(ierr);
27105b42c5fSBarry Smith 
27226283091SBarry Smith   ierr = PetscLayoutDestroy(A->rmap);CHKERRQ(ierr);
27326283091SBarry Smith   ierr = PetscLayoutDestroy(A->cmap);CHKERRQ(ierr);
2745c9eb25fSBarry Smith   ierr = PetscFListDestroy(&((PetscObject)A)->qlist);CHKERRQ(ierr);
2755c9eb25fSBarry Smith   ierr = PetscOListDestroy(((PetscObject)A)->olist);CHKERRQ(ierr);
276273d9f13SBarry Smith 
277273d9f13SBarry Smith   /* copy C over to A */
278273d9f13SBarry Smith   ierr  = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
279273d9f13SBarry Smith 
280273d9f13SBarry Smith   /* return the parts of A we saved */
2817adad957SLisandro Dalcin   ((PetscObject)A)->bops      = Abops;
282273d9f13SBarry Smith   A->ops                      = Aops;
2837adad957SLisandro Dalcin   ((PetscObject)A)->refct     = refct;
2847adad957SLisandro Dalcin   ((PetscObject)A)->type_name = mtype;
2857adad957SLisandro Dalcin   ((PetscObject)A)->name      = mname;
28630735b05SKris Buschelman   A->spptr                    = spptr;
287273d9f13SBarry Smith 
2885c9eb25fSBarry Smith   /* since these two are copied into A we do not want them destroyed in C */
2895c9eb25fSBarry Smith   ((PetscObject)C)->qlist = 0;
2905c9eb25fSBarry Smith   ((PetscObject)C)->olist = 0;
291d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(C);CHKERRQ(ierr);
292273d9f13SBarry Smith   PetscFunctionReturn(0);
293273d9f13SBarry Smith }
2948ab5b326SKris Buschelman /*
2958ab5b326SKris Buschelman         Replace A's header with that of C
2968ab5b326SKris Buschelman         This is essentially code moved from MatDestroy
297d0f46423SBarry Smith 
298d0f46423SBarry Smith         This is somewhat different from MatHeaderCopy() it would be nice to merge the code
2998ab5b326SKris Buschelman */
3008ab5b326SKris Buschelman #undef __FUNCT__
3018ab5b326SKris Buschelman #define __FUNCT__ "MatHeaderReplace"
3028ab5b326SKris Buschelman PetscErrorCode MatHeaderReplace(Mat A,Mat C)
3038ab5b326SKris Buschelman {
3048ab5b326SKris Buschelman   PetscErrorCode ierr;
3058ab5b326SKris Buschelman 
3068ab5b326SKris Buschelman   PetscFunctionBegin;
3076d7c1e57SBarry Smith   if (A == C) PetscFunctionReturn(0);
3086d7c1e57SBarry Smith 
3098ab5b326SKris Buschelman   /* free all the interior data structures from mat */
3108ab5b326SKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
311d2c95936SBarry Smith   ierr = PetscHeaderDestroy_Private((PetscObject)A);CHKERRQ(ierr);
312c650bb5fSLisandro Dalcin   ierr = PetscFree(A->ops);CHKERRQ(ierr);
31326283091SBarry Smith   ierr = PetscLayoutDestroy(A->rmap);CHKERRQ(ierr);
31426283091SBarry Smith   ierr = PetscLayoutDestroy(A->cmap);CHKERRQ(ierr);
315d95db058SSatish Balay   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
3168ab5b326SKris Buschelman 
3178ab5b326SKris Buschelman   /* copy C over to A */
3188ab5b326SKris Buschelman   if (C) {
3198ab5b326SKris Buschelman     ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
3208ab5b326SKris Buschelman     ierr = PetscLogObjectDestroy((PetscObject)C);CHKERRQ(ierr);
3218ab5b326SKris Buschelman     ierr = PetscFree(C);CHKERRQ(ierr);
3228ab5b326SKris Buschelman   }
3238ab5b326SKris Buschelman   PetscFunctionReturn(0);
3248ab5b326SKris Buschelman }
325