xref: /petsc/src/mat/utils/gcreate.c (revision d38fa0fb9e96b5504a279959047722757b3fe1d5)
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;
798ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES
80dfbe8321SBarry Smith   PetscErrorCode ierr;
818ba1e511SMatthew Knepley #endif
82273d9f13SBarry Smith 
83273d9f13SBarry Smith   PetscFunctionBegin;
844482741eSBarry Smith   PetscValidPointer(A,6);
8577431f27SBarry 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);
8677431f27SBarry 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);
8797f1f81fSBarry Smith 
888ba1e511SMatthew Knepley   *A = PETSC_NULL;
898ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES
908ba1e511SMatthew Knepley   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
918ba1e511SMatthew Knepley #endif
928ba1e511SMatthew Knepley 
93273d9f13SBarry Smith   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,0,"Mat",comm,MatDestroy,MatView);
94273d9f13SBarry Smith 
95273d9f13SBarry Smith   B->m             = m;
96273d9f13SBarry Smith   B->n             = n;
97273d9f13SBarry Smith   B->M             = M;
98273d9f13SBarry Smith   B->N             = N;
99521d7252SBarry Smith   B->bs            = 1;
100273d9f13SBarry Smith   B->preallocated  = PETSC_FALSE;
10135d8aa7fSBarry Smith   B->bops->publish = MatPublish_Base;
102273d9f13SBarry Smith   *A               = B;
103273d9f13SBarry Smith   PetscFunctionReturn(0);
104273d9f13SBarry Smith }
105273d9f13SBarry Smith 
1064a2ae208SSatish Balay #undef __FUNCT__
1074a2ae208SSatish Balay #define __FUNCT__ "MatSetFromOptions"
108273d9f13SBarry Smith /*@C
109273d9f13SBarry Smith    MatSetFromOptions - Creates a matrix where the type is determined
110273d9f13SBarry Smith    from the options database. Generates a parallel MPI matrix if the
111273d9f13SBarry Smith    communicator has more than one processor.  The default matrix type is
1127e5f4302SBarry Smith    AIJ, using the routines MatCreateSeqAIJ() and MatCreateMPIAIJ() if
1137e5f4302SBarry Smith    you do not select a type in the options database.
114273d9f13SBarry Smith 
115273d9f13SBarry Smith    Collective on Mat
116273d9f13SBarry Smith 
117273d9f13SBarry Smith    Input Parameter:
118273d9f13SBarry Smith .  A - the matrix
119273d9f13SBarry Smith 
120273d9f13SBarry Smith    Options Database Keys:
121273d9f13SBarry Smith +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
122273d9f13SBarry Smith .    -mat_type mpiaij   - AIJ type, uses MatCreateMPIAIJ()
123273d9f13SBarry Smith .    -mat_type seqbdiag - block diagonal type, uses MatCreateSeqBDiag()
124273d9f13SBarry Smith .    -mat_type mpibdiag - block diagonal type, uses MatCreateMPIBDiag()
125273d9f13SBarry Smith .    -mat_type mpirowbs - rowbs type, uses MatCreateMPIRowbs()
126273d9f13SBarry Smith .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
127273d9f13SBarry Smith .    -mat_type mpidense - dense type, uses MatCreateMPIDense()
128273d9f13SBarry Smith .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
129273d9f13SBarry Smith -    -mat_type mpibaij  - block AIJ type, uses MatCreateMPIBAIJ()
130273d9f13SBarry Smith 
131273d9f13SBarry Smith    Even More Options Database Keys:
132273d9f13SBarry Smith    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
133273d9f13SBarry Smith    for additional format-specific options.
134bd9ce289SLois Curfman McInnes 
1351d69843bSLois Curfman McInnes    Level: beginner
1361d69843bSLois Curfman McInnes 
137dc401e71SLois Curfman McInnes .keywords: matrix, create
138e0b365e2SLois Curfman McInnes 
139fafbff53SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
140fafbff53SBarry Smith           MatCreateSeqBDiag(),MatCreateMPIBDiag(),
14139ddd567SLois Curfman McInnes           MatCreateSeqDense(), MatCreateMPIDense(),
142a209d233SLois Curfman McInnes           MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
143a209d233SLois Curfman McInnes           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
144273d9f13SBarry Smith           MatConvert()
1457807a1faSBarry Smith @*/
146dfbe8321SBarry Smith PetscErrorCode MatSetFromOptions(Mat B)
1477807a1faSBarry Smith {
148dfbe8321SBarry Smith   PetscErrorCode ierr;
1492d156eb4SBarry Smith   PetscMPIInt    size;
150273d9f13SBarry Smith   char           mtype[256];
151273d9f13SBarry Smith   PetscTruth     flg;
152dbb450caSBarry Smith 
1533a40ed3dSBarry Smith   PetscFunctionBegin;
154b0a32e0cSBarry Smith   ierr = PetscOptionsGetString(B->prefix,"-mat_type",mtype,256,&flg);CHKERRQ(ierr);
155273d9f13SBarry Smith   if (flg) {
156273d9f13SBarry Smith     ierr = MatSetType(B,mtype);CHKERRQ(ierr);
157273d9f13SBarry Smith   }
158273d9f13SBarry Smith   if (!B->type_name) {
159273d9f13SBarry Smith     ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
1601d9a38edSKris Buschelman     ierr = MatSetType(B,MATAIJ);CHKERRQ(ierr);
161dfa27b74SSatish Balay   }
1623a40ed3dSBarry Smith   PetscFunctionReturn(0);
1637807a1faSBarry Smith }
1647807a1faSBarry Smith 
1654a2ae208SSatish Balay #undef __FUNCT__
1664a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation"
167273d9f13SBarry Smith /*@C
168273d9f13SBarry Smith    MatSetUpPreallocation
169dae03382SLois Curfman McInnes 
170273d9f13SBarry Smith    Collective on Mat
171dae03382SLois Curfman McInnes 
172273d9f13SBarry Smith    Input Parameter:
173273d9f13SBarry Smith .  A - the matrix
174d5d45c9bSBarry Smith 
175273d9f13SBarry Smith    Level: beginner
176d5d45c9bSBarry Smith 
177273d9f13SBarry Smith .keywords: matrix, create
178273d9f13SBarry Smith 
179273d9f13SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
180273d9f13SBarry Smith           MatCreateSeqBDiag(),MatCreateMPIBDiag(),
181273d9f13SBarry Smith           MatCreateSeqDense(), MatCreateMPIDense(),
182273d9f13SBarry Smith           MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
183273d9f13SBarry Smith           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
184273d9f13SBarry Smith           MatConvert()
185273d9f13SBarry Smith @*/
186dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation(Mat B)
187273d9f13SBarry Smith {
188dfbe8321SBarry Smith   PetscErrorCode ierr;
189273d9f13SBarry Smith 
190273d9f13SBarry Smith   PetscFunctionBegin;
191273d9f13SBarry Smith   if (B->ops->setuppreallocation) {
1928ba1e511SMatthew Knepley     PetscLogInfo(B,"MatSetUpPreallocation: Warning not preallocating matrix storage");
193273d9f13SBarry Smith     ierr = (*B->ops->setuppreallocation)(B);CHKERRQ(ierr);
194273d9f13SBarry Smith     B->ops->setuppreallocation = 0;
195273d9f13SBarry Smith     B->preallocated            = PETSC_TRUE;
196273d9f13SBarry Smith   }
197273d9f13SBarry Smith   PetscFunctionReturn(0);
198273d9f13SBarry Smith }
199273d9f13SBarry Smith 
200273d9f13SBarry Smith /*
201273d9f13SBarry Smith         Copies from Cs header to A
202273d9f13SBarry Smith */
2034a2ae208SSatish Balay #undef __FUNCT__
2044a2ae208SSatish Balay #define __FUNCT__ "MatHeaderCopy"
205dfbe8321SBarry Smith PetscErrorCode MatHeaderCopy(Mat A,Mat C)
206273d9f13SBarry Smith {
207dfbe8321SBarry Smith   PetscErrorCode ierr;
208d44834fbSBarry Smith   PetscInt       refct;
209273d9f13SBarry Smith   PetscOps       *Abops;
210273d9f13SBarry Smith   MatOps         Aops;
211273d9f13SBarry Smith   char           *mtype,*mname;
212273d9f13SBarry Smith 
213273d9f13SBarry Smith   PetscFunctionBegin;
214273d9f13SBarry Smith   /* free all the interior data structures from mat */
215273d9f13SBarry Smith   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
216273d9f13SBarry Smith 
2178a124369SBarry Smith   ierr = PetscMapDestroy(A->rmap);CHKERRQ(ierr);
2188a124369SBarry Smith   ierr = PetscMapDestroy(A->cmap);CHKERRQ(ierr);
219273d9f13SBarry Smith 
220273d9f13SBarry Smith   /* save the parts of A we need */
221273d9f13SBarry Smith   Abops = A->bops;
222273d9f13SBarry Smith   Aops  = A->ops;
223273d9f13SBarry Smith   refct = A->refct;
224273d9f13SBarry Smith   mtype = A->type_name;
225273d9f13SBarry Smith   mname = A->name;
226273d9f13SBarry Smith 
227273d9f13SBarry Smith   /* copy C over to A */
228273d9f13SBarry Smith   ierr  = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
229273d9f13SBarry Smith 
230273d9f13SBarry Smith   /* return the parts of A we saved */
231273d9f13SBarry Smith   A->bops      = Abops;
232273d9f13SBarry Smith   A->ops       = Aops;
233273d9f13SBarry Smith   A->qlist     = 0;
234273d9f13SBarry Smith   A->refct     = refct;
235273d9f13SBarry Smith   A->type_name = mtype;
236273d9f13SBarry Smith   A->name      = mname;
237273d9f13SBarry Smith 
238*d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(C);CHKERRQ(ierr);
239273d9f13SBarry Smith   PetscFunctionReturn(0);
240273d9f13SBarry Smith }
241