xref: /petsc/src/mat/utils/gcreate.c (revision 3194b578e315f5cda09358e793b5717ffe550174)
17807a1faSBarry Smith 
2c6db04a5SJed Brown #include <private/matimpl.h>       /*I "petscmat.h"  I*/
37807a1faSBarry Smith 
47adad957SLisandro Dalcin #if 0
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 }
127adad957SLisandro Dalcin #endif
1335d8aa7fSBarry Smith 
144a2ae208SSatish Balay #undef __FUNCT__
154a2ae208SSatish Balay #define __FUNCT__ "MatCreate"
1605869f15SSatish Balay /*@
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 
27f69a0ea3SMatthew Knepley    Input Parameter:
28f69a0ea3SMatthew Knepley .  comm - MPI communicator
297807a1faSBarry Smith 
307807a1faSBarry Smith    Output Parameter:
31dc401e71SLois Curfman McInnes .  A - the matrix
32e0b365e2SLois Curfman McInnes 
33273d9f13SBarry Smith    Options Database Keys:
34273d9f13SBarry Smith +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
35273d9f13SBarry Smith .    -mat_type mpiaij   - AIJ type, uses MatCreateMPIAIJ()
36273d9f13SBarry Smith .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
37273d9f13SBarry Smith .    -mat_type mpidense - dense type, uses MatCreateMPIDense()
38273d9f13SBarry Smith .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
39273d9f13SBarry Smith -    -mat_type mpibaij  - block AIJ type, uses MatCreateMPIBAIJ()
40e0b365e2SLois Curfman McInnes 
4183e1b59cSLois Curfman McInnes    Even More Options Database Keys:
4283e1b59cSLois Curfman McInnes    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
4383e1b59cSLois Curfman McInnes    for additional format-specific options.
44e0b365e2SLois Curfman McInnes 
45bd9ce289SLois Curfman McInnes    Notes:
46ec6e0d80SSatish Balay 
47273d9f13SBarry Smith    Level: beginner
48273d9f13SBarry Smith 
49a2fc510eSBarry Smith    User manual sections:
50a2fc510eSBarry Smith +   sec_matcreate
51a2fc510eSBarry Smith -   chapter_matrices
52a2fc510eSBarry Smith 
53273d9f13SBarry Smith .keywords: matrix, create
54273d9f13SBarry Smith 
55cd05a4c0SHong Zhang .seealso: MatCreateSeqAIJ(), MatCreateMPIAIJ(),
56273d9f13SBarry Smith           MatCreateSeqDense(), MatCreateMPIDense(),
576531c467SSatish Balay           MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
58273d9f13SBarry Smith           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
59273d9f13SBarry Smith           MatConvert()
60273d9f13SBarry Smith @*/
617087cfbeSBarry Smith PetscErrorCode  MatCreate(MPI_Comm comm,Mat *A)
62273d9f13SBarry Smith {
63273d9f13SBarry Smith   Mat            B;
64dfbe8321SBarry Smith   PetscErrorCode ierr;
65273d9f13SBarry Smith 
66273d9f13SBarry Smith   PetscFunctionBegin;
67f69a0ea3SMatthew Knepley   PetscValidPointer(A,2);
6897f1f81fSBarry Smith 
698ba1e511SMatthew Knepley   *A = PETSC_NULL;
708ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES
718ba1e511SMatthew Knepley   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
728ba1e511SMatthew Knepley #endif
738ba1e511SMatthew Knepley 
74*3194b578SJed Brown   ierr = PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_CLASSID,0,"Mat","Matrix","Mat",comm,MatDestroy,MatView);CHKERRQ(ierr);
7526283091SBarry Smith   ierr = PetscLayoutCreate(comm,&B->rmap);CHKERRQ(ierr);
7626283091SBarry Smith   ierr = PetscLayoutCreate(comm,&B->cmap);CHKERRQ(ierr);
77273d9f13SBarry Smith   B->preallocated  = PETSC_FALSE;
78273d9f13SBarry Smith   *A               = B;
79273d9f13SBarry Smith   PetscFunctionReturn(0);
80273d9f13SBarry Smith }
81273d9f13SBarry Smith 
824a2ae208SSatish Balay #undef __FUNCT__
83f69a0ea3SMatthew Knepley #define __FUNCT__ "MatSetSizes"
84f69a0ea3SMatthew Knepley /*@
85f69a0ea3SMatthew Knepley   MatSetSizes - Sets the local and global sizes, and checks to determine compatibility
86f69a0ea3SMatthew Knepley 
87f69a0ea3SMatthew Knepley   Collective on Mat
88f69a0ea3SMatthew Knepley 
89f69a0ea3SMatthew Knepley   Input Parameters:
90f69a0ea3SMatthew Knepley +  A - the matrix
91f69a0ea3SMatthew Knepley .  m - number of local rows (or PETSC_DECIDE)
92f69a0ea3SMatthew Knepley .  n - number of local columns (or PETSC_DECIDE)
93f69a0ea3SMatthew Knepley .  M - number of global rows (or PETSC_DETERMINE)
94f69a0ea3SMatthew Knepley -  N - number of global columns (or PETSC_DETERMINE)
95f69a0ea3SMatthew Knepley 
96f69a0ea3SMatthew Knepley    Notes:
97f69a0ea3SMatthew Knepley    m (n) and M (N) cannot be both PETSC_DECIDE
98f69a0ea3SMatthew Knepley    If one processor calls this with M (N) of PETSC_DECIDE then all processors must, otherwise the program will hang.
99f69a0ea3SMatthew Knepley 
100f69a0ea3SMatthew Knepley    If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the
101f69a0ea3SMatthew Knepley    user must ensure that they are chosen to be compatible with the
102f69a0ea3SMatthew Knepley    vectors. To do this, one first considers the matrix-vector product
103f69a0ea3SMatthew Knepley    'y = A x'. The 'm' that is used in the above routine must match the
104f69a0ea3SMatthew Knepley    local size used in the vector creation routine VecCreateMPI() for 'y'.
105f69a0ea3SMatthew Knepley    Likewise, the 'n' used must match that used as the local size in
106f69a0ea3SMatthew Knepley    VecCreateMPI() for 'x'.
107f69a0ea3SMatthew Knepley 
108f69a0ea3SMatthew Knepley   Level: beginner
109f69a0ea3SMatthew Knepley 
110f69a0ea3SMatthew Knepley .seealso: MatGetSize(), PetscSplitOwnership()
111f69a0ea3SMatthew Knepley @*/
1127087cfbeSBarry Smith PetscErrorCode  MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N)
113f69a0ea3SMatthew Knepley {
114284134d9SBarry Smith   PetscErrorCode ierr;
115284134d9SBarry Smith 
116f69a0ea3SMatthew Knepley   PetscFunctionBegin;
1170700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
118e32f2f54SBarry Smith   if (M > 0 && m > M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local column size %D cannot be larger than global column size %D",m,M);
119e32f2f54SBarry Smith   if (N > 0 && n > N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local row size %D cannot be larger than global row size %D",n,N);
120284134d9SBarry Smith   if (A->ops->setsizes) {
121284134d9SBarry Smith     /* Since this will not be set until the type has been set, this will NOT be called on the initial
122284134d9SBarry Smith        call of MatSetSizes() (which must be called BEFORE MatSetType() */
123284134d9SBarry Smith     ierr = (*A->ops->setsizes)(A,m,n,M,N);CHKERRQ(ierr);
124284134d9SBarry Smith   } else {
125e32f2f54SBarry Smith     if ((A->rmap->n >= 0 || A->rmap->N >= 0) && (A->rmap->n != m || A->rmap->N != M)) SETERRQ4(PETSC_COMM_SELF,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);
126e32f2f54SBarry Smith     if ((A->cmap->n >= 0 || A->cmap->N >= 0) && (A->cmap->n != n || A->cmap->N != N)) SETERRQ4(PETSC_COMM_SELF,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);
127284134d9SBarry Smith   }
128d0f46423SBarry Smith   A->rmap->n = m;
129d0f46423SBarry Smith   A->cmap->n = n;
130d0f46423SBarry Smith   A->rmap->N = M;
131d0f46423SBarry Smith   A->cmap->N = N;
132117016b1SBarry Smith 
133f69a0ea3SMatthew Knepley   PetscFunctionReturn(0);
134f69a0ea3SMatthew Knepley }
135f69a0ea3SMatthew Knepley 
136f69a0ea3SMatthew Knepley #undef __FUNCT__
1374a2ae208SSatish Balay #define __FUNCT__ "MatSetFromOptions"
13805869f15SSatish Balay /*@
139273d9f13SBarry Smith    MatSetFromOptions - Creates a matrix where the type is determined
140273d9f13SBarry Smith    from the options database. Generates a parallel MPI matrix if the
141273d9f13SBarry Smith    communicator has more than one processor.  The default matrix type is
1427e5f4302SBarry Smith    AIJ, using the routines MatCreateSeqAIJ() and MatCreateMPIAIJ() if
1437e5f4302SBarry Smith    you do not select a type in the options database.
144273d9f13SBarry Smith 
145273d9f13SBarry Smith    Collective on Mat
146273d9f13SBarry Smith 
147273d9f13SBarry Smith    Input Parameter:
148273d9f13SBarry Smith .  A - the matrix
149273d9f13SBarry Smith 
150273d9f13SBarry Smith    Options Database Keys:
151273d9f13SBarry Smith +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
152273d9f13SBarry Smith .    -mat_type mpiaij   - AIJ type, uses MatCreateMPIAIJ()
153273d9f13SBarry Smith .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
154273d9f13SBarry Smith .    -mat_type mpidense - dense type, uses MatCreateMPIDense()
155273d9f13SBarry Smith .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
156273d9f13SBarry Smith -    -mat_type mpibaij  - block AIJ type, uses MatCreateMPIBAIJ()
157273d9f13SBarry Smith 
158273d9f13SBarry Smith    Even More Options Database Keys:
159273d9f13SBarry Smith    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
160273d9f13SBarry Smith    for additional format-specific options.
161bd9ce289SLois Curfman McInnes 
1621d69843bSLois Curfman McInnes    Level: beginner
1631d69843bSLois Curfman McInnes 
164dc401e71SLois Curfman McInnes .keywords: matrix, create
165e0b365e2SLois Curfman McInnes 
166fafbff53SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
16739ddd567SLois Curfman McInnes           MatCreateSeqDense(), MatCreateMPIDense(),
1686531c467SSatish Balay           MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
169a209d233SLois Curfman McInnes           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
170273d9f13SBarry Smith           MatConvert()
1717807a1faSBarry Smith @*/
1727087cfbeSBarry Smith PetscErrorCode  MatSetFromOptions(Mat B)
1737807a1faSBarry Smith {
174dfbe8321SBarry Smith   PetscErrorCode ierr;
175f3be49caSLisandro Dalcin   const char     *deft = MATAIJ;
176f3be49caSLisandro Dalcin   char           type[256];
17769df5c0cSJed Brown   PetscBool      flg,set;
178dbb450caSBarry Smith 
1793a40ed3dSBarry Smith   PetscFunctionBegin;
1800700a824SBarry Smith   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
181f3be49caSLisandro Dalcin 
182*3194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)B);CHKERRQ(ierr);
183f3be49caSLisandro Dalcin     ierr = PetscOptionsList("-mat_type","Matrix type","MatSetType",MatList,deft,type,256,&flg);CHKERRQ(ierr);
184273d9f13SBarry Smith     if (flg) {
185f3be49caSLisandro Dalcin       ierr = MatSetType(B,type);CHKERRQ(ierr);
186f3be49caSLisandro Dalcin     } else if (!((PetscObject)B)->type_name) {
187f3be49caSLisandro Dalcin       ierr = MatSetType(B,deft);CHKERRQ(ierr);
188273d9f13SBarry Smith     }
189f3be49caSLisandro Dalcin 
190f3be49caSLisandro Dalcin     if (B->ops->setfromoptions) {
191f3be49caSLisandro Dalcin       ierr = (*B->ops->setfromoptions)(B);CHKERRQ(ierr);
192593c1905SSatish Balay     }
193f3be49caSLisandro Dalcin 
19469df5c0cSJed Brown     flg = PETSC_FALSE;
19569df5c0cSJed Brown     ierr = PetscOptionsBool("-mat_new_nonzero_location_err","Generate an error if new nonzeros are created in the matrix structure (useful to test preallocation)","MatSetOption",flg,&flg,&set);CHKERRQ(ierr);
19669df5c0cSJed Brown     if (set) {ierr = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,flg);CHKERRQ(ierr);}
19769df5c0cSJed Brown     flg = PETSC_FALSE;
19869df5c0cSJed Brown     ierr = PetscOptionsBool("-mat_new_nonzero_allocation_err","Generate an error if new nonzeros are allocated in the matrix structure (useful to test preallocation)","MatSetOption",flg,&flg,&set);CHKERRQ(ierr);
19969df5c0cSJed Brown     if (set) {ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,flg);CHKERRQ(ierr);}
20069df5c0cSJed Brown 
2015d973c19SBarry Smith     /* process any options handlers added with PetscObjectAddOptionsHandler() */
2025d973c19SBarry Smith     ierr = PetscObjectProcessOptionsHandlers((PetscObject)B);CHKERRQ(ierr);
203f3be49caSLisandro Dalcin   ierr = PetscOptionsEnd();CHKERRQ(ierr);
204f3be49caSLisandro Dalcin 
2053a40ed3dSBarry Smith   PetscFunctionReturn(0);
2067807a1faSBarry Smith }
2077807a1faSBarry Smith 
2084a2ae208SSatish Balay #undef __FUNCT__
2094a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation"
210bc08b0f1SBarry Smith /*@
211e0cc20caSBarry Smith    MatSetUpPreallocation - If the user has not set preallocation for this matrix then a default preallocation that is likely to be inefficient is used.
212dae03382SLois Curfman McInnes 
213273d9f13SBarry Smith    Collective on Mat
214dae03382SLois Curfman McInnes 
215273d9f13SBarry Smith    Input Parameter:
216273d9f13SBarry Smith .  A - the matrix
217d5d45c9bSBarry Smith 
218e0cc20caSBarry Smith    Level: advanced
219e0cc20caSBarry Smith 
220e0cc20caSBarry Smith    Notes: See the Performance chapter of the PETSc users manual for how to preallocate matrices
221d5d45c9bSBarry Smith 
222273d9f13SBarry Smith .keywords: matrix, create
223273d9f13SBarry Smith 
224273d9f13SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(),
225273d9f13SBarry Smith           MatCreateSeqDense(), MatCreateMPIDense(),
2266531c467SSatish Balay           MatCreateSeqBAIJ(), MatCreateMPIBAIJ(),
227273d9f13SBarry Smith           MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(),
228273d9f13SBarry Smith           MatConvert()
229273d9f13SBarry Smith @*/
2307087cfbeSBarry Smith PetscErrorCode  MatSetUpPreallocation(Mat B)
231273d9f13SBarry Smith {
232dfbe8321SBarry Smith   PetscErrorCode ierr;
233273d9f13SBarry Smith 
234273d9f13SBarry Smith   PetscFunctionBegin;
23517667f90SBarry Smith   if (!B->preallocated && B->ops->setuppreallocation) {
236ae15b995SBarry Smith     ierr = PetscInfo(B,"Warning not preallocating matrix storage\n");CHKERRQ(ierr);
237273d9f13SBarry Smith     ierr = (*B->ops->setuppreallocation)(B);CHKERRQ(ierr);
238273d9f13SBarry Smith   }
239210f0121SBarry Smith   B->preallocated = PETSC_TRUE;
240273d9f13SBarry Smith   PetscFunctionReturn(0);
241273d9f13SBarry Smith }
242273d9f13SBarry Smith 
243273d9f13SBarry Smith /*
244eb6b5d47SBarry Smith         Merges some information from Cs header to A; the C object is then destroyed
245d0f46423SBarry Smith 
246d0f46423SBarry Smith         This is somewhat different from MatHeaderReplace() it would be nice to merge the code
247273d9f13SBarry Smith */
2484a2ae208SSatish Balay #undef __FUNCT__
249eb6b5d47SBarry Smith #define __FUNCT__ "MatHeaderMerge"
250eb6b5d47SBarry Smith PetscErrorCode MatHeaderMerge(Mat A,Mat C)
251273d9f13SBarry Smith {
252dfbe8321SBarry Smith   PetscErrorCode         ierr;
253d44834fbSBarry Smith   PetscInt               refct;
254273d9f13SBarry Smith   PetscOps               *Abops;
255273d9f13SBarry Smith   MatOps                 Aops;
256273d9f13SBarry Smith   char                   *mtype,*mname;
25730735b05SKris Buschelman   void                   *spptr;
258273d9f13SBarry Smith 
259273d9f13SBarry Smith   PetscFunctionBegin;
260273d9f13SBarry Smith   /* save the parts of A we need */
2617adad957SLisandro Dalcin   Abops = ((PetscObject)A)->bops;
262273d9f13SBarry Smith   Aops  = A->ops;
2637adad957SLisandro Dalcin   refct = ((PetscObject)A)->refct;
2645c9eb25fSBarry Smith   mtype = ((PetscObject)A)->type_name;
2655c9eb25fSBarry Smith   mname = ((PetscObject)A)->name;
26630735b05SKris Buschelman   spptr = A->spptr;
26730735b05SKris Buschelman 
2685c9eb25fSBarry Smith   /* zero these so the destroy below does not free them */
2695c9eb25fSBarry Smith   ((PetscObject)A)->type_name = 0;
2705c9eb25fSBarry Smith   ((PetscObject)A)->name      = 0;
2715c9eb25fSBarry Smith 
2727c99f97cSSatish Balay   /* free all the interior data structures from mat */
2737c99f97cSSatish Balay   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
2747c99f97cSSatish Balay 
27530735b05SKris Buschelman   ierr = PetscFree(C->spptr);CHKERRQ(ierr);
27605b42c5fSBarry Smith 
2776bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->rmap);CHKERRQ(ierr);
2786bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->cmap);CHKERRQ(ierr);
2795c9eb25fSBarry Smith   ierr = PetscFListDestroy(&((PetscObject)A)->qlist);CHKERRQ(ierr);
2806bf464f9SBarry Smith   ierr = PetscOListDestroy(&((PetscObject)A)->olist);CHKERRQ(ierr);
281273d9f13SBarry Smith 
282273d9f13SBarry Smith   /* copy C over to A */
283273d9f13SBarry Smith   ierr  = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
284273d9f13SBarry Smith 
285273d9f13SBarry Smith   /* return the parts of A we saved */
2867adad957SLisandro Dalcin   ((PetscObject)A)->bops      = Abops;
287273d9f13SBarry Smith   A->ops                      = Aops;
2887adad957SLisandro Dalcin   ((PetscObject)A)->refct     = refct;
2897adad957SLisandro Dalcin   ((PetscObject)A)->type_name = mtype;
2907adad957SLisandro Dalcin   ((PetscObject)A)->name      = mname;
29130735b05SKris Buschelman   A->spptr                    = spptr;
292273d9f13SBarry Smith 
2935c9eb25fSBarry Smith   /* since these two are copied into A we do not want them destroyed in C */
2945c9eb25fSBarry Smith   ((PetscObject)C)->qlist = 0;
2955c9eb25fSBarry Smith   ((PetscObject)C)->olist = 0;
2966bf464f9SBarry Smith   ierr = PetscHeaderDestroy(&C);CHKERRQ(ierr);
297273d9f13SBarry Smith   PetscFunctionReturn(0);
298273d9f13SBarry Smith }
2998ab5b326SKris Buschelman /*
300eb6b5d47SBarry Smith         Replace A's header with that of C; the C object is then destroyed
301d0f46423SBarry Smith 
302eb6b5d47SBarry Smith         This is essentially code moved from MatDestroy()
303eb6b5d47SBarry Smith 
304eb6b5d47SBarry Smith         This is somewhat different from MatHeaderMerge() it would be nice to merge the code
3058ab5b326SKris Buschelman */
3068ab5b326SKris Buschelman #undef __FUNCT__
3078ab5b326SKris Buschelman #define __FUNCT__ "MatHeaderReplace"
3088ab5b326SKris Buschelman PetscErrorCode MatHeaderReplace(Mat A,Mat C)
3098ab5b326SKris Buschelman {
3108ab5b326SKris Buschelman   PetscErrorCode ierr;
31127b31e29SJed Brown   PetscInt refct;
3128ab5b326SKris Buschelman 
3138ab5b326SKris Buschelman   PetscFunctionBegin;
31427b31e29SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
31527b31e29SJed Brown   PetscValidHeaderSpecific(C,MAT_CLASSID,2);
3166d7c1e57SBarry Smith   if (A == C) PetscFunctionReturn(0);
317a7e80f87SJed Brown   PetscCheckSameComm(A,1,C,2);
318a7e80f87SJed Brown   if (((PetscObject)C)->refct != 1) SETERRQ1(((PetscObject)C)->comm,PETSC_ERR_ARG_WRONGSTATE,"Object C has refct %D > 1, would leave hanging reference",((PetscObject)C)->refct);
3196d7c1e57SBarry Smith 
3208ab5b326SKris Buschelman   /* free all the interior data structures from mat */
3218ab5b326SKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
322d2c95936SBarry Smith   ierr = PetscHeaderDestroy_Private((PetscObject)A);CHKERRQ(ierr);
323c650bb5fSLisandro Dalcin   ierr = PetscFree(A->ops);CHKERRQ(ierr);
3246bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->rmap);CHKERRQ(ierr);
3256bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->cmap);CHKERRQ(ierr);
326d95db058SSatish Balay   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
3278ab5b326SKris Buschelman 
3288ab5b326SKris Buschelman   /* copy C over to A */
32927b31e29SJed Brown   refct = ((PetscObject)A)->refct;
3308ab5b326SKris Buschelman   ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
33127b31e29SJed Brown   ((PetscObject)A)->refct = refct;
3328ab5b326SKris Buschelman   ierr = PetscLogObjectDestroy((PetscObject)C);CHKERRQ(ierr);
3338ab5b326SKris Buschelman   ierr = PetscFree(C);CHKERRQ(ierr);
3348ab5b326SKris Buschelman   PetscFunctionReturn(0);
3358ab5b326SKris Buschelman }
336