xref: /petsc/src/mat/utils/gcreate.c (revision be1d678a52e6eff2808b2fa31ae986cdbf03c9fe)
1*be1d678aSKris Buschelman #define PETSCMAT_DLL
27807a1faSBarry Smith 
3273d9f13SBarry Smith #include "src/mat/matimpl.h"       /*I "petscmat.h"  I*/
4325e03aeSBarry Smith #include "petscsys.h"
57807a1faSBarry Smith 
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 }
1335d8aa7fSBarry Smith 
1435d8aa7fSBarry Smith 
154a2ae208SSatish Balay #undef __FUNCT__
164a2ae208SSatish Balay #define __FUNCT__ "MatCreate"
17325ab940SBarry Smith /*@C
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 
287807a1faSBarry Smith    Input Parameters:
2982b900a8SBarry Smith +  m - number of local rows (or PETSC_DECIDE)
3082b900a8SBarry Smith .  n - number of local columns (or PETSC_DECIDE)
3182b900a8SBarry Smith .  M - number of global rows (or PETSC_DETERMINE)
3282b900a8SBarry Smith .  N - number of global columns (or PETSC_DETERMINE)
33cb13003dSBarry Smith -  comm - MPI communicator
347807a1faSBarry Smith 
357807a1faSBarry Smith    Output Parameter:
36dc401e71SLois Curfman McInnes .  A - the matrix
37e0b365e2SLois Curfman McInnes 
38273d9f13SBarry Smith    Options Database Keys:
39273d9f13SBarry Smith +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
40273d9f13SBarry Smith .    -mat_type mpiaij   - AIJ type, uses MatCreateMPIAIJ()
41273d9f13SBarry Smith .    -mat_type seqbdiag - block diagonal type, uses MatCreateSeqBDiag()
42273d9f13SBarry Smith .    -mat_type mpibdiag - block diagonal type, uses MatCreateMPIBDiag()
43273d9f13SBarry Smith .    -mat_type mpirowbs - rowbs type, uses MatCreateMPIRowbs()
44273d9f13SBarry Smith .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
45273d9f13SBarry Smith .    -mat_type mpidense - dense type, uses MatCreateMPIDense()
46273d9f13SBarry Smith .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
47273d9f13SBarry Smith -    -mat_type mpibaij  - block AIJ type, uses MatCreateMPIBAIJ()
48e0b365e2SLois Curfman McInnes 
4983e1b59cSLois Curfman McInnes    Even More Options Database Keys:
5083e1b59cSLois Curfman McInnes    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
5183e1b59cSLois Curfman McInnes    for additional format-specific options.
52e0b365e2SLois Curfman McInnes 
53bd9ce289SLois Curfman McInnes    Notes:
54ec6e0d80SSatish Balay    If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the
55ec6e0d80SSatish Balay    user must ensure that they are chosen to be compatible with the
56ec6e0d80SSatish Balay    vectors. To do this, one first considers the matrix-vector product
57ec6e0d80SSatish Balay    'y = A x'. The 'm' that is used in the above routine must match the
58ec6e0d80SSatish Balay    local size used in the vector creation routine VecCreateMPI() for 'y'.
59ec6e0d80SSatish Balay    Likewise, the 'n' used must match that used as the local size in
60ec6e0d80SSatish Balay    VecCreateMPI() for 'x'.
61ec6e0d80SSatish Balay 
62273d9f13SBarry Smith    Level: beginner
63273d9f13SBarry Smith 
64a2fc510eSBarry Smith    User manual sections:
65a2fc510eSBarry Smith +   sec_matcreate
66a2fc510eSBarry Smith -   chapter_matrices
67a2fc510eSBarry Smith 
68273d9f13SBarry Smith .keywords: matrix, create
69273d9f13SBarry Smith 
70273d9f13SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
71273d9f13SBarry Smith           MatCreateSeqBDiag(),MatCreateMPIBDiag(),
72273d9f13SBarry Smith           MatCreateSeqDense(), MatCreateMPIDense(),
73273d9f13SBarry Smith           MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
74273d9f13SBarry Smith           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
75273d9f13SBarry Smith           MatConvert()
76273d9f13SBarry Smith @*/
77*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A)
78273d9f13SBarry Smith {
79273d9f13SBarry Smith   Mat            B;
80dfbe8321SBarry Smith   PetscErrorCode ierr;
81273d9f13SBarry Smith 
82273d9f13SBarry Smith   PetscFunctionBegin;
834482741eSBarry Smith   PetscValidPointer(A,6);
8477431f27SBarry 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);
8577431f27SBarry 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);
8697f1f81fSBarry Smith 
878ba1e511SMatthew Knepley   *A = PETSC_NULL;
888ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES
898ba1e511SMatthew Knepley   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
908ba1e511SMatthew Knepley #endif
918ba1e511SMatthew Knepley 
9252e6d16bSBarry Smith   ierr = PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,0,"Mat",comm,MatDestroy,MatView);CHKERRQ(ierr);
93273d9f13SBarry Smith 
94273d9f13SBarry Smith   B->m             = m;
95273d9f13SBarry Smith   B->n             = n;
96273d9f13SBarry Smith   B->M             = M;
97273d9f13SBarry Smith   B->N             = N;
98521d7252SBarry Smith   B->bs            = 1;
99273d9f13SBarry Smith   B->preallocated  = PETSC_FALSE;
10035d8aa7fSBarry Smith   B->bops->publish = MatPublish_Base;
101273d9f13SBarry Smith   *A               = B;
102273d9f13SBarry Smith   PetscFunctionReturn(0);
103273d9f13SBarry Smith }
104273d9f13SBarry Smith 
1054a2ae208SSatish Balay #undef __FUNCT__
1064a2ae208SSatish Balay #define __FUNCT__ "MatSetFromOptions"
107273d9f13SBarry Smith /*@C
108273d9f13SBarry Smith    MatSetFromOptions - Creates a matrix where the type is determined
109273d9f13SBarry Smith    from the options database. Generates a parallel MPI matrix if the
110273d9f13SBarry Smith    communicator has more than one processor.  The default matrix type is
1117e5f4302SBarry Smith    AIJ, using the routines MatCreateSeqAIJ() and MatCreateMPIAIJ() if
1127e5f4302SBarry Smith    you do not select a type in the options database.
113273d9f13SBarry Smith 
114273d9f13SBarry Smith    Collective on Mat
115273d9f13SBarry Smith 
116273d9f13SBarry Smith    Input Parameter:
117273d9f13SBarry Smith .  A - the matrix
118273d9f13SBarry Smith 
119273d9f13SBarry Smith    Options Database Keys:
120273d9f13SBarry Smith +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
121273d9f13SBarry Smith .    -mat_type mpiaij   - AIJ type, uses MatCreateMPIAIJ()
122273d9f13SBarry Smith .    -mat_type seqbdiag - block diagonal type, uses MatCreateSeqBDiag()
123273d9f13SBarry Smith .    -mat_type mpibdiag - block diagonal type, uses MatCreateMPIBDiag()
124273d9f13SBarry Smith .    -mat_type mpirowbs - rowbs type, uses MatCreateMPIRowbs()
125273d9f13SBarry Smith .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
126273d9f13SBarry Smith .    -mat_type mpidense - dense type, uses MatCreateMPIDense()
127273d9f13SBarry Smith .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
128273d9f13SBarry Smith -    -mat_type mpibaij  - block AIJ type, uses MatCreateMPIBAIJ()
129273d9f13SBarry Smith 
130273d9f13SBarry Smith    Even More Options Database Keys:
131273d9f13SBarry Smith    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
132273d9f13SBarry Smith    for additional format-specific options.
133bd9ce289SLois Curfman McInnes 
1341d69843bSLois Curfman McInnes    Level: beginner
1351d69843bSLois Curfman McInnes 
136dc401e71SLois Curfman McInnes .keywords: matrix, create
137e0b365e2SLois Curfman McInnes 
138fafbff53SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
139fafbff53SBarry Smith           MatCreateSeqBDiag(),MatCreateMPIBDiag(),
14039ddd567SLois Curfman McInnes           MatCreateSeqDense(), MatCreateMPIDense(),
141a209d233SLois Curfman McInnes           MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
142a209d233SLois Curfman McInnes           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
143273d9f13SBarry Smith           MatConvert()
1447807a1faSBarry Smith @*/
145*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat B)
1467807a1faSBarry Smith {
147dfbe8321SBarry Smith   PetscErrorCode ierr;
148273d9f13SBarry Smith   char           mtype[256];
149273d9f13SBarry Smith   PetscTruth     flg;
150dbb450caSBarry Smith 
1513a40ed3dSBarry Smith   PetscFunctionBegin;
152b0a32e0cSBarry Smith   ierr = PetscOptionsGetString(B->prefix,"-mat_type",mtype,256,&flg);CHKERRQ(ierr);
153273d9f13SBarry Smith   if (flg) {
154273d9f13SBarry Smith     ierr = MatSetType(B,mtype);CHKERRQ(ierr);
155273d9f13SBarry Smith   }
156273d9f13SBarry Smith   if (!B->type_name) {
15730735b05SKris Buschelman     ierr = MatSetType(B,MATAIJ);CHKERRQ(ierr);
158dfa27b74SSatish Balay   }
1593a40ed3dSBarry Smith   PetscFunctionReturn(0);
1607807a1faSBarry Smith }
1617807a1faSBarry Smith 
1624a2ae208SSatish Balay #undef __FUNCT__
1634a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation"
164273d9f13SBarry Smith /*@C
165273d9f13SBarry Smith    MatSetUpPreallocation
166dae03382SLois Curfman McInnes 
167273d9f13SBarry Smith    Collective on Mat
168dae03382SLois Curfman McInnes 
169273d9f13SBarry Smith    Input Parameter:
170273d9f13SBarry Smith .  A - the matrix
171d5d45c9bSBarry Smith 
172273d9f13SBarry Smith    Level: beginner
173d5d45c9bSBarry Smith 
174273d9f13SBarry Smith .keywords: matrix, create
175273d9f13SBarry Smith 
176273d9f13SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
177273d9f13SBarry Smith           MatCreateSeqBDiag(),MatCreateMPIBDiag(),
178273d9f13SBarry Smith           MatCreateSeqDense(), MatCreateMPIDense(),
179273d9f13SBarry Smith           MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
180273d9f13SBarry Smith           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
181273d9f13SBarry Smith           MatConvert()
182273d9f13SBarry Smith @*/
183*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat B)
184273d9f13SBarry Smith {
185dfbe8321SBarry Smith   PetscErrorCode ierr;
186273d9f13SBarry Smith 
187273d9f13SBarry Smith   PetscFunctionBegin;
188273d9f13SBarry Smith   if (B->ops->setuppreallocation) {
18963ba0a88SBarry Smith     ierr = PetscLogInfo((B,"MatSetUpPreallocation: Warning not preallocating matrix storage\n"));CHKERRQ(ierr);
190273d9f13SBarry Smith     ierr = (*B->ops->setuppreallocation)(B);CHKERRQ(ierr);
191273d9f13SBarry Smith     B->ops->setuppreallocation = 0;
192273d9f13SBarry Smith     B->preallocated            = PETSC_TRUE;
193273d9f13SBarry Smith   }
194273d9f13SBarry Smith   PetscFunctionReturn(0);
195273d9f13SBarry Smith }
196273d9f13SBarry Smith 
197273d9f13SBarry Smith /*
198273d9f13SBarry Smith         Copies from Cs header to A
199273d9f13SBarry Smith */
2004a2ae208SSatish Balay #undef __FUNCT__
2014a2ae208SSatish Balay #define __FUNCT__ "MatHeaderCopy"
202dfbe8321SBarry Smith PetscErrorCode MatHeaderCopy(Mat A,Mat C)
203273d9f13SBarry Smith {
204dfbe8321SBarry Smith   PetscErrorCode ierr;
205d44834fbSBarry Smith   PetscInt       refct;
206273d9f13SBarry Smith   PetscOps       *Abops;
207273d9f13SBarry Smith   MatOps         Aops;
208273d9f13SBarry Smith   char           *mtype,*mname;
20930735b05SKris Buschelman   void           *spptr;
210273d9f13SBarry Smith 
211273d9f13SBarry Smith   PetscFunctionBegin;
212273d9f13SBarry Smith   /* free all the interior data structures from mat */
213273d9f13SBarry Smith   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
214273d9f13SBarry Smith 
2158a124369SBarry Smith   ierr = PetscMapDestroy(A->rmap);CHKERRQ(ierr);
2168a124369SBarry Smith   ierr = PetscMapDestroy(A->cmap);CHKERRQ(ierr);
217273d9f13SBarry Smith 
218273d9f13SBarry Smith   /* save the parts of A we need */
219273d9f13SBarry Smith   Abops = A->bops;
220273d9f13SBarry Smith   Aops  = A->ops;
221273d9f13SBarry Smith   refct = A->refct;
222273d9f13SBarry Smith   mtype = A->type_name;
223273d9f13SBarry Smith   mname = A->name;
22430735b05SKris Buschelman   spptr = A->spptr;
22530735b05SKris Buschelman 
22630735b05SKris Buschelman   if (C->spptr) {
22730735b05SKris Buschelman     ierr = PetscFree(C->spptr);CHKERRQ(ierr);
22830735b05SKris Buschelman     C->spptr = PETSC_NULL;
22930735b05SKris Buschelman   }
230273d9f13SBarry Smith 
231273d9f13SBarry Smith   /* copy C over to A */
232273d9f13SBarry Smith   ierr  = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
233273d9f13SBarry Smith 
234273d9f13SBarry Smith   /* return the parts of A we saved */
235273d9f13SBarry Smith   A->bops      = Abops;
236273d9f13SBarry Smith   A->ops       = Aops;
237273d9f13SBarry Smith   A->qlist     = 0;
238273d9f13SBarry Smith   A->refct     = refct;
239273d9f13SBarry Smith   A->type_name = mtype;
240273d9f13SBarry Smith   A->name      = mname;
24130735b05SKris Buschelman   A->spptr     = spptr;
242273d9f13SBarry Smith 
243d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(C);CHKERRQ(ierr);
244273d9f13SBarry Smith   PetscFunctionReturn(0);
245273d9f13SBarry Smith }
2468ab5b326SKris Buschelman /*
2478ab5b326SKris Buschelman         Replace A's header with that of C
2488ab5b326SKris Buschelman         This is essentially code moved from MatDestroy
2498ab5b326SKris Buschelman */
2508ab5b326SKris Buschelman #undef __FUNCT__
2518ab5b326SKris Buschelman #define __FUNCT__ "MatHeaderReplace"
2528ab5b326SKris Buschelman PetscErrorCode MatHeaderReplace(Mat A,Mat C)
2538ab5b326SKris Buschelman {
2548ab5b326SKris Buschelman   PetscErrorCode ierr;
2558ab5b326SKris Buschelman 
2568ab5b326SKris Buschelman   PetscFunctionBegin;
2578ab5b326SKris Buschelman   /* free all the interior data structures from mat */
2588ab5b326SKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
259d2c95936SBarry Smith   ierr = PetscHeaderDestroy_Private((PetscObject)A);CHKERRQ(ierr);
260cb8bee85SHong Zhang   ierr = PetscMapDestroy(A->rmap);CHKERRQ(ierr);
261cb8bee85SHong Zhang   ierr = PetscMapDestroy(A->cmap);CHKERRQ(ierr);
2628ab5b326SKris Buschelman 
2638ab5b326SKris Buschelman   /* copy C over to A */
2648ab5b326SKris Buschelman   if (C) {
2658ab5b326SKris Buschelman     ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
2668ab5b326SKris Buschelman     ierr = PetscLogObjectDestroy((PetscObject)C);CHKERRQ(ierr);
2678ab5b326SKris Buschelman     ierr = PetscFree(C);CHKERRQ(ierr);
2688ab5b326SKris Buschelman   }
2698ab5b326SKris Buschelman   PetscFunctionReturn(0);
2708ab5b326SKris Buschelman }
271