xref: /petsc/src/mat/utils/gcreate.c (revision 8ab5b3267b42fcbe56c0122fd1434be193b60e5a)
17807a1faSBarry Smith 
2273d9f13SBarry Smith #include "src/mat/matimpl.h"       /*I "petscmat.h"  I*/
3325e03aeSBarry Smith #include "petscsys.h"
47807a1faSBarry Smith 
54a2ae208SSatish Balay #undef __FUNCT__
64a2ae208SSatish Balay #define __FUNCT__ "MatPublish_Base"
76849ba73SBarry Smith static PetscErrorCode MatPublish_Base(PetscObject obj)
835d8aa7fSBarry Smith {
935d8aa7fSBarry Smith   PetscFunctionBegin;
1035d8aa7fSBarry Smith   PetscFunctionReturn(0);
1135d8aa7fSBarry Smith }
1235d8aa7fSBarry Smith 
1335d8aa7fSBarry Smith 
144a2ae208SSatish Balay #undef __FUNCT__
154a2ae208SSatish Balay #define __FUNCT__ "MatCreate"
16325ab940SBarry Smith /*@C
1769dd0797SLois Curfman McInnes    MatCreate - Creates a matrix where the type is determined
187e5f4302SBarry Smith    from either a call to MatSetType() or from the options database
197e5f4302SBarry Smith    with a call to MatSetFromOptions(). The default matrix type is
207e5f4302SBarry Smith    AIJ, using the routines MatCreateSeqAIJ() or MatCreateMPIAIJ()
217e5f4302SBarry Smith    if you do not set a type in the options database. If you never
227e5f4302SBarry Smith    call MatSetType() or MatSetFromOptions() it will generate an
23f8ab6608SSatish Balay    error when you try to use the matrix.
2483e1b59cSLois Curfman McInnes 
25cb13003dSBarry Smith    Collective on MPI_Comm
26cb13003dSBarry Smith 
277807a1faSBarry Smith    Input Parameters:
2882b900a8SBarry Smith +  m - number of local rows (or PETSC_DECIDE)
2982b900a8SBarry Smith .  n - number of local columns (or PETSC_DECIDE)
3082b900a8SBarry Smith .  M - number of global rows (or PETSC_DETERMINE)
3182b900a8SBarry Smith .  N - number of global columns (or PETSC_DETERMINE)
32cb13003dSBarry Smith -  comm - MPI communicator
337807a1faSBarry Smith 
347807a1faSBarry Smith    Output Parameter:
35dc401e71SLois Curfman McInnes .  A - the matrix
36e0b365e2SLois Curfman McInnes 
37273d9f13SBarry Smith    Options Database Keys:
38273d9f13SBarry Smith +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
39273d9f13SBarry Smith .    -mat_type mpiaij   - AIJ type, uses MatCreateMPIAIJ()
40273d9f13SBarry Smith .    -mat_type seqbdiag - block diagonal type, uses MatCreateSeqBDiag()
41273d9f13SBarry Smith .    -mat_type mpibdiag - block diagonal type, uses MatCreateMPIBDiag()
42273d9f13SBarry Smith .    -mat_type mpirowbs - rowbs type, uses MatCreateMPIRowbs()
43273d9f13SBarry Smith .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
44273d9f13SBarry Smith .    -mat_type mpidense - dense type, uses MatCreateMPIDense()
45273d9f13SBarry Smith .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
46273d9f13SBarry Smith -    -mat_type mpibaij  - block AIJ type, uses MatCreateMPIBAIJ()
47e0b365e2SLois Curfman McInnes 
4883e1b59cSLois Curfman McInnes    Even More Options Database Keys:
4983e1b59cSLois Curfman McInnes    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
5083e1b59cSLois Curfman McInnes    for additional format-specific options.
51e0b365e2SLois Curfman McInnes 
52bd9ce289SLois Curfman McInnes    Notes:
53ec6e0d80SSatish Balay    If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the
54ec6e0d80SSatish Balay    user must ensure that they are chosen to be compatible with the
55ec6e0d80SSatish Balay    vectors. To do this, one first considers the matrix-vector product
56ec6e0d80SSatish Balay    'y = A x'. The 'm' that is used in the above routine must match the
57ec6e0d80SSatish Balay    local size used in the vector creation routine VecCreateMPI() for 'y'.
58ec6e0d80SSatish Balay    Likewise, the 'n' used must match that used as the local size in
59ec6e0d80SSatish Balay    VecCreateMPI() for 'x'.
60ec6e0d80SSatish Balay 
61273d9f13SBarry Smith    Level: beginner
62273d9f13SBarry Smith 
63a2fc510eSBarry Smith    User manual sections:
64a2fc510eSBarry Smith +   sec_matcreate
65a2fc510eSBarry Smith -   chapter_matrices
66a2fc510eSBarry Smith 
67273d9f13SBarry Smith .keywords: matrix, create
68273d9f13SBarry Smith 
69273d9f13SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
70273d9f13SBarry Smith           MatCreateSeqBDiag(),MatCreateMPIBDiag(),
71273d9f13SBarry Smith           MatCreateSeqDense(), MatCreateMPIDense(),
72273d9f13SBarry Smith           MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
73273d9f13SBarry Smith           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
74273d9f13SBarry Smith           MatConvert()
75273d9f13SBarry Smith @*/
7638baddfdSBarry Smith PetscErrorCode MatCreate(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A)
77273d9f13SBarry Smith {
78273d9f13SBarry Smith   Mat            B;
79dfbe8321SBarry Smith   PetscErrorCode ierr;
80273d9f13SBarry Smith 
81273d9f13SBarry Smith   PetscFunctionBegin;
824482741eSBarry Smith   PetscValidPointer(A,6);
8377431f27SBarry Smith   if (M > 0 && m > M) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Local column size %D cannot be larger than global column size %D",m,M);
8477431f27SBarry Smith   if (N > 0 && n > N) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Local row size %D cannot be larger than global row size %D",n,N);
8597f1f81fSBarry Smith 
868ba1e511SMatthew Knepley   *A = PETSC_NULL;
878ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES
888ba1e511SMatthew Knepley   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
898ba1e511SMatthew Knepley #endif
908ba1e511SMatthew Knepley 
9152e6d16bSBarry Smith   ierr = PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,0,"Mat",comm,MatDestroy,MatView);CHKERRQ(ierr);
92273d9f13SBarry Smith 
93273d9f13SBarry Smith   B->m             = m;
94273d9f13SBarry Smith   B->n             = n;
95273d9f13SBarry Smith   B->M             = M;
96273d9f13SBarry Smith   B->N             = N;
97521d7252SBarry Smith   B->bs            = 1;
98273d9f13SBarry Smith   B->preallocated  = PETSC_FALSE;
9935d8aa7fSBarry Smith   B->bops->publish = MatPublish_Base;
100273d9f13SBarry Smith   *A               = B;
101273d9f13SBarry Smith   PetscFunctionReturn(0);
102273d9f13SBarry Smith }
103273d9f13SBarry Smith 
1044a2ae208SSatish Balay #undef __FUNCT__
1054a2ae208SSatish Balay #define __FUNCT__ "MatSetFromOptions"
106273d9f13SBarry Smith /*@C
107273d9f13SBarry Smith    MatSetFromOptions - Creates a matrix where the type is determined
108273d9f13SBarry Smith    from the options database. Generates a parallel MPI matrix if the
109273d9f13SBarry Smith    communicator has more than one processor.  The default matrix type is
1107e5f4302SBarry Smith    AIJ, using the routines MatCreateSeqAIJ() and MatCreateMPIAIJ() if
1117e5f4302SBarry Smith    you do not select a type in the options database.
112273d9f13SBarry Smith 
113273d9f13SBarry Smith    Collective on Mat
114273d9f13SBarry Smith 
115273d9f13SBarry Smith    Input Parameter:
116273d9f13SBarry Smith .  A - the matrix
117273d9f13SBarry Smith 
118273d9f13SBarry Smith    Options Database Keys:
119273d9f13SBarry Smith +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
120273d9f13SBarry Smith .    -mat_type mpiaij   - AIJ type, uses MatCreateMPIAIJ()
121273d9f13SBarry Smith .    -mat_type seqbdiag - block diagonal type, uses MatCreateSeqBDiag()
122273d9f13SBarry Smith .    -mat_type mpibdiag - block diagonal type, uses MatCreateMPIBDiag()
123273d9f13SBarry Smith .    -mat_type mpirowbs - rowbs type, uses MatCreateMPIRowbs()
124273d9f13SBarry Smith .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
125273d9f13SBarry Smith .    -mat_type mpidense - dense type, uses MatCreateMPIDense()
126273d9f13SBarry Smith .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
127273d9f13SBarry Smith -    -mat_type mpibaij  - block AIJ type, uses MatCreateMPIBAIJ()
128273d9f13SBarry Smith 
129273d9f13SBarry Smith    Even More Options Database Keys:
130273d9f13SBarry Smith    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
131273d9f13SBarry Smith    for additional format-specific options.
132bd9ce289SLois Curfman McInnes 
1331d69843bSLois Curfman McInnes    Level: beginner
1341d69843bSLois Curfman McInnes 
135dc401e71SLois Curfman McInnes .keywords: matrix, create
136e0b365e2SLois Curfman McInnes 
137fafbff53SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
138fafbff53SBarry Smith           MatCreateSeqBDiag(),MatCreateMPIBDiag(),
13939ddd567SLois Curfman McInnes           MatCreateSeqDense(), MatCreateMPIDense(),
140a209d233SLois Curfman McInnes           MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
141a209d233SLois Curfman McInnes           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
142273d9f13SBarry Smith           MatConvert()
1437807a1faSBarry Smith @*/
144dfbe8321SBarry Smith PetscErrorCode MatSetFromOptions(Mat B)
1457807a1faSBarry Smith {
146dfbe8321SBarry Smith   PetscErrorCode ierr;
147273d9f13SBarry Smith   char           mtype[256];
148273d9f13SBarry Smith   PetscTruth     flg;
149dbb450caSBarry Smith 
1503a40ed3dSBarry Smith   PetscFunctionBegin;
151b0a32e0cSBarry Smith   ierr = PetscOptionsGetString(B->prefix,"-mat_type",mtype,256,&flg);CHKERRQ(ierr);
152273d9f13SBarry Smith   if (flg) {
153273d9f13SBarry Smith     ierr = MatSetType(B,mtype);CHKERRQ(ierr);
154273d9f13SBarry Smith   }
155273d9f13SBarry Smith   if (!B->type_name) {
15630735b05SKris Buschelman     ierr = MatSetType(B,MATAIJ);CHKERRQ(ierr);
157dfa27b74SSatish Balay   }
1583a40ed3dSBarry Smith   PetscFunctionReturn(0);
1597807a1faSBarry Smith }
1607807a1faSBarry Smith 
1614a2ae208SSatish Balay #undef __FUNCT__
1624a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation"
163273d9f13SBarry Smith /*@C
164273d9f13SBarry Smith    MatSetUpPreallocation
165dae03382SLois Curfman McInnes 
166273d9f13SBarry Smith    Collective on Mat
167dae03382SLois Curfman McInnes 
168273d9f13SBarry Smith    Input Parameter:
169273d9f13SBarry Smith .  A - the matrix
170d5d45c9bSBarry Smith 
171273d9f13SBarry Smith    Level: beginner
172d5d45c9bSBarry Smith 
173273d9f13SBarry Smith .keywords: matrix, create
174273d9f13SBarry Smith 
175273d9f13SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
176273d9f13SBarry Smith           MatCreateSeqBDiag(),MatCreateMPIBDiag(),
177273d9f13SBarry Smith           MatCreateSeqDense(), MatCreateMPIDense(),
178273d9f13SBarry Smith           MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
179273d9f13SBarry Smith           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
180273d9f13SBarry Smith           MatConvert()
181273d9f13SBarry Smith @*/
182dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation(Mat B)
183273d9f13SBarry Smith {
184dfbe8321SBarry Smith   PetscErrorCode ierr;
185273d9f13SBarry Smith 
186273d9f13SBarry Smith   PetscFunctionBegin;
187273d9f13SBarry Smith   if (B->ops->setuppreallocation) {
1888ba1e511SMatthew Knepley     PetscLogInfo(B,"MatSetUpPreallocation: Warning not preallocating matrix storage");
189273d9f13SBarry Smith     ierr = (*B->ops->setuppreallocation)(B);CHKERRQ(ierr);
190273d9f13SBarry Smith     B->ops->setuppreallocation = 0;
191273d9f13SBarry Smith     B->preallocated            = PETSC_TRUE;
192273d9f13SBarry Smith   }
193273d9f13SBarry Smith   PetscFunctionReturn(0);
194273d9f13SBarry Smith }
195273d9f13SBarry Smith 
196273d9f13SBarry Smith /*
197273d9f13SBarry Smith         Copies from Cs header to A
198273d9f13SBarry Smith */
1994a2ae208SSatish Balay #undef __FUNCT__
2004a2ae208SSatish Balay #define __FUNCT__ "MatHeaderCopy"
201dfbe8321SBarry Smith PetscErrorCode MatHeaderCopy(Mat A,Mat C)
202273d9f13SBarry Smith {
203dfbe8321SBarry Smith   PetscErrorCode ierr;
204d44834fbSBarry Smith   PetscInt       refct;
205273d9f13SBarry Smith   PetscOps       *Abops;
206273d9f13SBarry Smith   MatOps         Aops;
207273d9f13SBarry Smith   char           *mtype,*mname;
20830735b05SKris Buschelman   void           *spptr;
209273d9f13SBarry Smith 
210273d9f13SBarry Smith   PetscFunctionBegin;
211273d9f13SBarry Smith   /* free all the interior data structures from mat */
212273d9f13SBarry Smith   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
213273d9f13SBarry Smith 
2148a124369SBarry Smith   ierr = PetscMapDestroy(A->rmap);CHKERRQ(ierr);
2158a124369SBarry Smith   ierr = PetscMapDestroy(A->cmap);CHKERRQ(ierr);
216273d9f13SBarry Smith 
217273d9f13SBarry Smith   /* save the parts of A we need */
218273d9f13SBarry Smith   Abops = A->bops;
219273d9f13SBarry Smith   Aops  = A->ops;
220273d9f13SBarry Smith   refct = A->refct;
221273d9f13SBarry Smith   mtype = A->type_name;
222273d9f13SBarry Smith   mname = A->name;
22330735b05SKris Buschelman   spptr = A->spptr;
22430735b05SKris Buschelman 
22530735b05SKris Buschelman   if (C->spptr) {
22630735b05SKris Buschelman     ierr = PetscFree(C->spptr);CHKERRQ(ierr);
22730735b05SKris Buschelman     C->spptr = PETSC_NULL;
22830735b05SKris Buschelman   }
229273d9f13SBarry Smith 
230273d9f13SBarry Smith   /* copy C over to A */
231273d9f13SBarry Smith   ierr  = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
232273d9f13SBarry Smith 
233273d9f13SBarry Smith   /* return the parts of A we saved */
234273d9f13SBarry Smith   A->bops      = Abops;
235273d9f13SBarry Smith   A->ops       = Aops;
236273d9f13SBarry Smith   A->qlist     = 0;
237273d9f13SBarry Smith   A->refct     = refct;
238273d9f13SBarry Smith   A->type_name = mtype;
239273d9f13SBarry Smith   A->name      = mname;
24030735b05SKris Buschelman   A->spptr     = spptr;
241273d9f13SBarry Smith 
242d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(C);CHKERRQ(ierr);
243273d9f13SBarry Smith   PetscFunctionReturn(0);
244273d9f13SBarry Smith }
245*8ab5b326SKris Buschelman /*
246*8ab5b326SKris Buschelman         Replace A's header with that of C
247*8ab5b326SKris Buschelman         This is essentially code moved from MatDestroy
248*8ab5b326SKris Buschelman */
249*8ab5b326SKris Buschelman #undef __FUNCT__
250*8ab5b326SKris Buschelman #define __FUNCT__ "MatHeaderReplace"
251*8ab5b326SKris Buschelman PetscErrorCode MatHeaderReplace(Mat A,Mat C)
252*8ab5b326SKris Buschelman {
253*8ab5b326SKris Buschelman   PetscErrorCode ierr;
254*8ab5b326SKris Buschelman 
255*8ab5b326SKris Buschelman   PetscFunctionBegin;
256*8ab5b326SKris Buschelman   /* free all the interior data structures from mat */
257*8ab5b326SKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
258*8ab5b326SKris Buschelman   ierr = PetscHeaderDestroy_Private(A);CHKERRQ(ierr);
259*8ab5b326SKris Buschelman 
260*8ab5b326SKris Buschelman   /* copy C over to A */
261*8ab5b326SKris Buschelman   if (C) {
262*8ab5b326SKris Buschelman     ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
263*8ab5b326SKris Buschelman 
264*8ab5b326SKris Buschelman     ierr = PetscLogObjectDestroy((PetscObject)C);CHKERRQ(ierr);
265*8ab5b326SKris Buschelman     ierr = PetscFree(C);CHKERRQ(ierr);
266*8ab5b326SKris Buschelman   }
267*8ab5b326SKris Buschelman 
268*8ab5b326SKris Buschelman   PetscFunctionReturn(0);
269*8ab5b326SKris Buschelman }
270