xref: /petsc/src/mat/utils/gcreate.c (revision f73d5cc4fab9985badaa1a29b152f456bf59ad34)
17807a1faSBarry Smith 
2b45d2f2cSJed Brown #include <petsc-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
2069b1f4b7SBarry Smith    AIJ, using the routines MatCreateSeqAIJ() or MatCreateAIJ()
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()
3569b1f4b7SBarry Smith .    -mat_type mpiaij   - AIJ type, uses MatCreateAIJ()
36273d9f13SBarry Smith .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
3769b1f4b7SBarry Smith .    -mat_type mpidense - dense type, uses MatCreateDense()
38273d9f13SBarry Smith .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
3969b1f4b7SBarry Smith -    -mat_type mpibaij  - block AIJ type, uses MatCreateBAIJ()
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 
5569b1f4b7SBarry Smith .seealso: MatCreateSeqAIJ(), MatCreateAIJ(),
5669b1f4b7SBarry Smith           MatCreateSeqDense(), MatCreateDense(),
5769b1f4b7SBarry Smith           MatCreateSeqBAIJ(), MatCreateBAIJ(),
5869b1f4b7SBarry Smith           MatCreateSeqSBAIJ(), MatCreateSBAIJ(),
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 
743194b578SJed 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 
108*f73d5cc4SBarry Smith    You cannot change the sizes once they have been set.
109*f73d5cc4SBarry Smith 
110*f73d5cc4SBarry Smith    The sizes must be set before MatSetUp() or MatXXXSetPreallocation() is called.
111*f73d5cc4SBarry Smith 
112f69a0ea3SMatthew Knepley   Level: beginner
113f69a0ea3SMatthew Knepley 
114f69a0ea3SMatthew Knepley .seealso: MatGetSize(), PetscSplitOwnership()
115f69a0ea3SMatthew Knepley @*/
1167087cfbeSBarry Smith PetscErrorCode  MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N)
117f69a0ea3SMatthew Knepley {
118f69a0ea3SMatthew Knepley   PetscFunctionBegin;
1190700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
120e32f2f54SBarry 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);
121e32f2f54SBarry 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);
122e32f2f54SBarry 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);
123e32f2f54SBarry 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);
124d0f46423SBarry Smith   A->rmap->n = m;
125d0f46423SBarry Smith   A->cmap->n = n;
126d0f46423SBarry Smith   A->rmap->N = M;
127d0f46423SBarry Smith   A->cmap->N = N;
128f69a0ea3SMatthew Knepley   PetscFunctionReturn(0);
129f69a0ea3SMatthew Knepley }
130f69a0ea3SMatthew Knepley 
131f69a0ea3SMatthew Knepley #undef __FUNCT__
1324a2ae208SSatish Balay #define __FUNCT__ "MatSetFromOptions"
13305869f15SSatish Balay /*@
134273d9f13SBarry Smith    MatSetFromOptions - Creates a matrix where the type is determined
135273d9f13SBarry Smith    from the options database. Generates a parallel MPI matrix if the
136273d9f13SBarry Smith    communicator has more than one processor.  The default matrix type is
13769b1f4b7SBarry Smith    AIJ, using the routines MatCreateSeqAIJ() and MatCreateAIJ() if
1387e5f4302SBarry Smith    you do not select a type in the options database.
139273d9f13SBarry Smith 
140273d9f13SBarry Smith    Collective on Mat
141273d9f13SBarry Smith 
142273d9f13SBarry Smith    Input Parameter:
143273d9f13SBarry Smith .  A - the matrix
144273d9f13SBarry Smith 
145273d9f13SBarry Smith    Options Database Keys:
146273d9f13SBarry Smith +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
14769b1f4b7SBarry Smith .    -mat_type mpiaij   - AIJ type, uses MatCreateAIJ()
148273d9f13SBarry Smith .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
14969b1f4b7SBarry Smith .    -mat_type mpidense - dense type, uses MatCreateDense()
150273d9f13SBarry Smith .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
15169b1f4b7SBarry Smith -    -mat_type mpibaij  - block AIJ type, uses MatCreateBAIJ()
152273d9f13SBarry Smith 
153273d9f13SBarry Smith    Even More Options Database Keys:
154273d9f13SBarry Smith    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
155273d9f13SBarry Smith    for additional format-specific options.
156bd9ce289SLois Curfman McInnes 
1571d69843bSLois Curfman McInnes    Level: beginner
1581d69843bSLois Curfman McInnes 
159dc401e71SLois Curfman McInnes .keywords: matrix, create
160e0b365e2SLois Curfman McInnes 
16169b1f4b7SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateAIJ(),
16269b1f4b7SBarry Smith           MatCreateSeqDense(), MatCreateDense(),
16369b1f4b7SBarry Smith           MatCreateSeqBAIJ(), MatCreateBAIJ(),
16469b1f4b7SBarry Smith           MatCreateSeqSBAIJ(), MatCreateSBAIJ(),
165273d9f13SBarry Smith           MatConvert()
1667807a1faSBarry Smith @*/
1677087cfbeSBarry Smith PetscErrorCode  MatSetFromOptions(Mat B)
1687807a1faSBarry Smith {
169dfbe8321SBarry Smith   PetscErrorCode ierr;
170f3be49caSLisandro Dalcin   const char     *deft = MATAIJ;
171f3be49caSLisandro Dalcin   char           type[256];
17269df5c0cSJed Brown   PetscBool      flg,set;
173dbb450caSBarry Smith 
1743a40ed3dSBarry Smith   PetscFunctionBegin;
1750700a824SBarry Smith   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
176f3be49caSLisandro Dalcin 
1773194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)B);CHKERRQ(ierr);
178f3be49caSLisandro Dalcin     ierr = PetscOptionsList("-mat_type","Matrix type","MatSetType",MatList,deft,type,256,&flg);CHKERRQ(ierr);
179273d9f13SBarry Smith     if (flg) {
180f3be49caSLisandro Dalcin       ierr = MatSetType(B,type);CHKERRQ(ierr);
181f3be49caSLisandro Dalcin     } else if (!((PetscObject)B)->type_name) {
182f3be49caSLisandro Dalcin       ierr = MatSetType(B,deft);CHKERRQ(ierr);
183273d9f13SBarry Smith     }
184f3be49caSLisandro Dalcin 
185f3be49caSLisandro Dalcin     if (B->ops->setfromoptions) {
186f3be49caSLisandro Dalcin       ierr = (*B->ops->setfromoptions)(B);CHKERRQ(ierr);
187593c1905SSatish Balay     }
188f3be49caSLisandro Dalcin 
18969df5c0cSJed Brown     flg = PETSC_FALSE;
19069df5c0cSJed 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);
19169df5c0cSJed Brown     if (set) {ierr = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,flg);CHKERRQ(ierr);}
19269df5c0cSJed Brown     flg = PETSC_FALSE;
19369df5c0cSJed 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);
19469df5c0cSJed Brown     if (set) {ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,flg);CHKERRQ(ierr);}
19569df5c0cSJed Brown 
1965d973c19SBarry Smith     /* process any options handlers added with PetscObjectAddOptionsHandler() */
1975d973c19SBarry Smith     ierr = PetscObjectProcessOptionsHandlers((PetscObject)B);CHKERRQ(ierr);
198f3be49caSLisandro Dalcin   ierr = PetscOptionsEnd();CHKERRQ(ierr);
199f3be49caSLisandro Dalcin 
2003a40ed3dSBarry Smith   PetscFunctionReturn(0);
2017807a1faSBarry Smith }
2027807a1faSBarry Smith 
2034a2ae208SSatish Balay #undef __FUNCT__
20463562e91SJed Brown #define __FUNCT__ "MatXAIJSetPreallocation"
20563562e91SJed Brown /*@
20663562e91SJed Brown    MatXAIJSetPreallocation - set preallocation for serial and parallel AIJ, BAIJ, and SBAIJ matrices
20763562e91SJed Brown 
20863562e91SJed Brown    Collective on Mat
20963562e91SJed Brown 
21063562e91SJed Brown    Input Arguments:
21163562e91SJed Brown +  A - matrix being preallocated
21263562e91SJed Brown .  bs - block size
2133e5f4774SJed Brown .  dnnz - number of nonzero blocks per block row of diagonal part of parallel matrix
2143e5f4774SJed Brown .  onnz - number of nonzero blocks per block row of off-diagonal part of parallel matrix
2153e5f4774SJed Brown .  dnnzu - number of nonzero blocks per block row of upper-triangular part of diagonal part of parallel matrix
2163e5f4774SJed Brown -  onnzu - number of nonzero blocks per block row of upper-triangular part of off-diagonal part of parallel matrix
21763562e91SJed Brown 
21863562e91SJed Brown    Level: beginner
21963562e91SJed Brown 
22063562e91SJed Brown .seealso: MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation()
22163562e91SJed Brown @*/
2223e5f4774SJed Brown PetscErrorCode MatXAIJSetPreallocation(Mat A,PetscInt bs,const PetscInt *dnnz,const PetscInt *onnz,const PetscInt *dnnzu,const PetscInt *onnzu)
22363562e91SJed Brown {
22463562e91SJed Brown   PetscErrorCode ierr;
22563562e91SJed Brown   void (*aij)(void);
22663562e91SJed Brown 
22763562e91SJed Brown   PetscFunctionBegin;
2283e5f4774SJed Brown   ierr = MatSeqBAIJSetPreallocation(A,bs,0,dnnz);CHKERRQ(ierr);
2293e5f4774SJed Brown   ierr = MatMPIBAIJSetPreallocation(A,bs,0,dnnz,0,onnz);CHKERRQ(ierr);
2303e5f4774SJed Brown   ierr = MatSeqSBAIJSetPreallocation(A,bs,0,dnnzu);CHKERRQ(ierr);
2313e5f4774SJed Brown   ierr = MatMPISBAIJSetPreallocation(A,bs,0,dnnzu,0,onnzu);CHKERRQ(ierr);
23263562e91SJed Brown   /*
23363562e91SJed Brown     In general, we have to do extra work to preallocate for scalar (AIJ) matrices so we check whether it will do any
23463562e91SJed Brown     good before going on with it.
23563562e91SJed Brown   */
23663562e91SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
23763562e91SJed Brown   if (!aij) {
23863562e91SJed Brown     ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
23963562e91SJed Brown   }
24063562e91SJed Brown   if (aij) {
24163562e91SJed Brown     if (bs == 1) {
2423e5f4774SJed Brown       ierr = MatSeqAIJSetPreallocation(A,0,dnnz);CHKERRQ(ierr);
2433e5f4774SJed Brown       ierr = MatMPIAIJSetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr);
2443e5f4774SJed Brown     } else {                    /* Convert block-row precallocation to scalar-row */
24563562e91SJed Brown       PetscInt i,m,*sdnnz,*sonnz;
24663562e91SJed Brown       ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr);
24763562e91SJed Brown       ierr = PetscMalloc2((!!dnnz)*m*bs,PetscInt,&sdnnz,(!!onnz)*m*bs,PetscInt,&sonnz);CHKERRQ(ierr);
24863562e91SJed Brown       for (i=0; i<m*bs; i++) {
24963562e91SJed Brown         if (dnnz) sdnnz[i] = dnnz[i/bs] * bs;
25063562e91SJed Brown         if (onnz) sonnz[i] = onnz[i/bs] * bs;
25163562e91SJed Brown       }
2523e5f4774SJed Brown       ierr = MatSeqAIJSetPreallocation(A,0,dnnz?sdnnz:PETSC_NULL);CHKERRQ(ierr);
2533e5f4774SJed Brown       ierr = MatMPIAIJSetPreallocation(A,0,dnnz?sdnnz:PETSC_NULL,0,onnz?sonnz:PETSC_NULL);CHKERRQ(ierr);
25463562e91SJed Brown       ierr = PetscFree2(sdnnz,sonnz);CHKERRQ(ierr);
25563562e91SJed Brown     }
25663562e91SJed Brown   }
25763562e91SJed Brown   ierr = MatSetBlockSize(A,bs);CHKERRQ(ierr);
25863562e91SJed Brown   PetscFunctionReturn(0);
25963562e91SJed Brown }
26063562e91SJed Brown 
261273d9f13SBarry Smith /*
262eb6b5d47SBarry Smith         Merges some information from Cs header to A; the C object is then destroyed
263d0f46423SBarry Smith 
264d0f46423SBarry Smith         This is somewhat different from MatHeaderReplace() it would be nice to merge the code
265273d9f13SBarry Smith */
2664a2ae208SSatish Balay #undef __FUNCT__
267eb6b5d47SBarry Smith #define __FUNCT__ "MatHeaderMerge"
268eb6b5d47SBarry Smith PetscErrorCode MatHeaderMerge(Mat A,Mat C)
269273d9f13SBarry Smith {
270dfbe8321SBarry Smith   PetscErrorCode         ierr;
271d44834fbSBarry Smith   PetscInt               refct;
272273d9f13SBarry Smith   PetscOps               *Abops;
273273d9f13SBarry Smith   MatOps                 Aops;
274273d9f13SBarry Smith   char                   *mtype,*mname;
27530735b05SKris Buschelman   void                   *spptr;
276273d9f13SBarry Smith 
277273d9f13SBarry Smith   PetscFunctionBegin;
278273d9f13SBarry Smith   /* save the parts of A we need */
2797adad957SLisandro Dalcin   Abops = ((PetscObject)A)->bops;
280273d9f13SBarry Smith   Aops  = A->ops;
2817adad957SLisandro Dalcin   refct = ((PetscObject)A)->refct;
2825c9eb25fSBarry Smith   mtype = ((PetscObject)A)->type_name;
2835c9eb25fSBarry Smith   mname = ((PetscObject)A)->name;
28430735b05SKris Buschelman   spptr = A->spptr;
28530735b05SKris Buschelman 
2865c9eb25fSBarry Smith   /* zero these so the destroy below does not free them */
2875c9eb25fSBarry Smith   ((PetscObject)A)->type_name = 0;
2885c9eb25fSBarry Smith   ((PetscObject)A)->name      = 0;
2895c9eb25fSBarry Smith 
2907c99f97cSSatish Balay   /* free all the interior data structures from mat */
2917c99f97cSSatish Balay   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
2927c99f97cSSatish Balay 
29330735b05SKris Buschelman   ierr = PetscFree(C->spptr);CHKERRQ(ierr);
29405b42c5fSBarry Smith 
2956bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->rmap);CHKERRQ(ierr);
2966bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->cmap);CHKERRQ(ierr);
2975c9eb25fSBarry Smith   ierr = PetscFListDestroy(&((PetscObject)A)->qlist);CHKERRQ(ierr);
2986bf464f9SBarry Smith   ierr = PetscOListDestroy(&((PetscObject)A)->olist);CHKERRQ(ierr);
299273d9f13SBarry Smith 
300273d9f13SBarry Smith   /* copy C over to A */
301273d9f13SBarry Smith   ierr  = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
302273d9f13SBarry Smith 
303273d9f13SBarry Smith   /* return the parts of A we saved */
3047adad957SLisandro Dalcin   ((PetscObject)A)->bops      = Abops;
305273d9f13SBarry Smith   A->ops                      = Aops;
3067adad957SLisandro Dalcin   ((PetscObject)A)->refct     = refct;
3077adad957SLisandro Dalcin   ((PetscObject)A)->type_name = mtype;
3087adad957SLisandro Dalcin   ((PetscObject)A)->name      = mname;
30930735b05SKris Buschelman   A->spptr                    = spptr;
310273d9f13SBarry Smith 
3115c9eb25fSBarry Smith   /* since these two are copied into A we do not want them destroyed in C */
3125c9eb25fSBarry Smith   ((PetscObject)C)->qlist = 0;
3135c9eb25fSBarry Smith   ((PetscObject)C)->olist = 0;
3146bf464f9SBarry Smith   ierr = PetscHeaderDestroy(&C);CHKERRQ(ierr);
315273d9f13SBarry Smith   PetscFunctionReturn(0);
316273d9f13SBarry Smith }
3178ab5b326SKris Buschelman /*
318eb6b5d47SBarry Smith         Replace A's header with that of C; the C object is then destroyed
319d0f46423SBarry Smith 
320eb6b5d47SBarry Smith         This is essentially code moved from MatDestroy()
321eb6b5d47SBarry Smith 
322eb6b5d47SBarry Smith         This is somewhat different from MatHeaderMerge() it would be nice to merge the code
3238ab5b326SKris Buschelman */
3248ab5b326SKris Buschelman #undef __FUNCT__
3258ab5b326SKris Buschelman #define __FUNCT__ "MatHeaderReplace"
3268ab5b326SKris Buschelman PetscErrorCode MatHeaderReplace(Mat A,Mat C)
3278ab5b326SKris Buschelman {
3288ab5b326SKris Buschelman   PetscErrorCode ierr;
32927b31e29SJed Brown   PetscInt refct;
3308ab5b326SKris Buschelman 
3318ab5b326SKris Buschelman   PetscFunctionBegin;
33227b31e29SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
33327b31e29SJed Brown   PetscValidHeaderSpecific(C,MAT_CLASSID,2);
3346d7c1e57SBarry Smith   if (A == C) PetscFunctionReturn(0);
335a7e80f87SJed Brown   PetscCheckSameComm(A,1,C,2);
336a7e80f87SJed 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);
3376d7c1e57SBarry Smith 
3388ab5b326SKris Buschelman   /* free all the interior data structures from mat */
3398ab5b326SKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
340d2c95936SBarry Smith   ierr = PetscHeaderDestroy_Private((PetscObject)A);CHKERRQ(ierr);
341c650bb5fSLisandro Dalcin   ierr = PetscFree(A->ops);CHKERRQ(ierr);
3426bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->rmap);CHKERRQ(ierr);
3436bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->cmap);CHKERRQ(ierr);
344d95db058SSatish Balay   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
3458ab5b326SKris Buschelman 
3468ab5b326SKris Buschelman   /* copy C over to A */
34727b31e29SJed Brown   refct = ((PetscObject)A)->refct;
3488ab5b326SKris Buschelman   ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
34927b31e29SJed Brown   ((PetscObject)A)->refct = refct;
3508ab5b326SKris Buschelman   ierr = PetscLogObjectDestroy((PetscObject)C);CHKERRQ(ierr);
3518ab5b326SKris Buschelman   ierr = PetscFree(C);CHKERRQ(ierr);
3528ab5b326SKris Buschelman   PetscFunctionReturn(0);
3538ab5b326SKris Buschelman }
354