xref: /petsc/src/mat/utils/gcreate.c (revision 4cce697cf3e7f389152c2af9531432dd57fb072f)
17807a1faSBarry Smith 
2b45d2f2cSJed Brown #include <petsc-private/matimpl.h>       /*I "petscmat.h"  I*/
37807a1faSBarry Smith 
44a2ae208SSatish Balay #undef __FUNCT__
54a2ae208SSatish Balay #define __FUNCT__ "MatCreate"
605869f15SSatish Balay /*@
769dd0797SLois Curfman McInnes    MatCreate - Creates a matrix where the type is determined
87e5f4302SBarry Smith    from either a call to MatSetType() or from the options database
97e5f4302SBarry Smith    with a call to MatSetFromOptions(). The default matrix type is
1069b1f4b7SBarry Smith    AIJ, using the routines MatCreateSeqAIJ() or MatCreateAIJ()
117e5f4302SBarry Smith    if you do not set a type in the options database. If you never
127e5f4302SBarry Smith    call MatSetType() or MatSetFromOptions() it will generate an
13f8ab6608SSatish Balay    error when you try to use the matrix.
1483e1b59cSLois Curfman McInnes 
15cb13003dSBarry Smith    Collective on MPI_Comm
16cb13003dSBarry Smith 
17f69a0ea3SMatthew Knepley    Input Parameter:
18f69a0ea3SMatthew Knepley .  comm - MPI communicator
197807a1faSBarry Smith 
207807a1faSBarry Smith    Output Parameter:
21dc401e71SLois Curfman McInnes .  A - the matrix
22e0b365e2SLois Curfman McInnes 
23273d9f13SBarry Smith    Options Database Keys:
24273d9f13SBarry Smith +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
2569b1f4b7SBarry Smith .    -mat_type mpiaij   - AIJ type, uses MatCreateAIJ()
26273d9f13SBarry Smith .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
2769b1f4b7SBarry Smith .    -mat_type mpidense - dense type, uses MatCreateDense()
28273d9f13SBarry Smith .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
2969b1f4b7SBarry Smith -    -mat_type mpibaij  - block AIJ type, uses MatCreateBAIJ()
30e0b365e2SLois Curfman McInnes 
3183e1b59cSLois Curfman McInnes    Even More Options Database Keys:
3283e1b59cSLois Curfman McInnes    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
3383e1b59cSLois Curfman McInnes    for additional format-specific options.
34e0b365e2SLois Curfman McInnes 
35bd9ce289SLois Curfman McInnes    Notes:
36ec6e0d80SSatish Balay 
37273d9f13SBarry Smith    Level: beginner
38273d9f13SBarry Smith 
39a2fc510eSBarry Smith    User manual sections:
40a2fc510eSBarry Smith +   sec_matcreate
41a2fc510eSBarry Smith -   chapter_matrices
42a2fc510eSBarry Smith 
43273d9f13SBarry Smith .keywords: matrix, create
44273d9f13SBarry Smith 
4569b1f4b7SBarry Smith .seealso: MatCreateSeqAIJ(), MatCreateAIJ(),
4669b1f4b7SBarry Smith           MatCreateSeqDense(), MatCreateDense(),
4769b1f4b7SBarry Smith           MatCreateSeqBAIJ(), MatCreateBAIJ(),
4869b1f4b7SBarry Smith           MatCreateSeqSBAIJ(), MatCreateSBAIJ(),
49273d9f13SBarry Smith           MatConvert()
50273d9f13SBarry Smith @*/
517087cfbeSBarry Smith PetscErrorCode  MatCreate(MPI_Comm comm,Mat *A)
52273d9f13SBarry Smith {
53273d9f13SBarry Smith   Mat            B;
54dfbe8321SBarry Smith   PetscErrorCode ierr;
55273d9f13SBarry Smith 
56273d9f13SBarry Smith   PetscFunctionBegin;
57f69a0ea3SMatthew Knepley   PetscValidPointer(A,2);
5897f1f81fSBarry Smith 
590298fd71SBarry Smith   *A = NULL;
60519f805aSKarl Rupp #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
61607a6623SBarry Smith   ierr = MatInitializePackage();CHKERRQ(ierr);
628ba1e511SMatthew Knepley #endif
638ba1e511SMatthew Knepley 
6467c2884eSBarry Smith   ierr = PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_CLASSID,"Mat","Matrix","Mat",comm,MatDestroy,MatView);CHKERRQ(ierr);
6526283091SBarry Smith   ierr = PetscLayoutCreate(comm,&B->rmap);CHKERRQ(ierr);
6626283091SBarry Smith   ierr = PetscLayoutCreate(comm,&B->cmap);CHKERRQ(ierr);
6726fbe8dcSKarl Rupp 
68273d9f13SBarry Smith   B->preallocated = PETSC_FALSE;
69273d9f13SBarry Smith   *A              = B;
70273d9f13SBarry Smith   PetscFunctionReturn(0);
71273d9f13SBarry Smith }
72273d9f13SBarry Smith 
734a2ae208SSatish Balay #undef __FUNCT__
74f69a0ea3SMatthew Knepley #define __FUNCT__ "MatSetSizes"
75f69a0ea3SMatthew Knepley /*@
76f69a0ea3SMatthew Knepley   MatSetSizes - Sets the local and global sizes, and checks to determine compatibility
77f69a0ea3SMatthew Knepley 
78f69a0ea3SMatthew Knepley   Collective on Mat
79f69a0ea3SMatthew Knepley 
80f69a0ea3SMatthew Knepley   Input Parameters:
81f69a0ea3SMatthew Knepley +  A - the matrix
82f69a0ea3SMatthew Knepley .  m - number of local rows (or PETSC_DECIDE)
83f69a0ea3SMatthew Knepley .  n - number of local columns (or PETSC_DECIDE)
84f69a0ea3SMatthew Knepley .  M - number of global rows (or PETSC_DETERMINE)
85f69a0ea3SMatthew Knepley -  N - number of global columns (or PETSC_DETERMINE)
86f69a0ea3SMatthew Knepley 
87f69a0ea3SMatthew Knepley    Notes:
88f69a0ea3SMatthew Knepley    m (n) and M (N) cannot be both PETSC_DECIDE
89f69a0ea3SMatthew Knepley    If one processor calls this with M (N) of PETSC_DECIDE then all processors must, otherwise the program will hang.
90f69a0ea3SMatthew Knepley 
91f69a0ea3SMatthew Knepley    If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the
92f69a0ea3SMatthew Knepley    user must ensure that they are chosen to be compatible with the
93f69a0ea3SMatthew Knepley    vectors. To do this, one first considers the matrix-vector product
94f69a0ea3SMatthew Knepley    'y = A x'. The 'm' that is used in the above routine must match the
95f69a0ea3SMatthew Knepley    local size used in the vector creation routine VecCreateMPI() for 'y'.
96f69a0ea3SMatthew Knepley    Likewise, the 'n' used must match that used as the local size in
97f69a0ea3SMatthew Knepley    VecCreateMPI() for 'x'.
98f69a0ea3SMatthew Knepley 
99f73d5cc4SBarry Smith    You cannot change the sizes once they have been set.
100f73d5cc4SBarry Smith 
101f73d5cc4SBarry Smith    The sizes must be set before MatSetUp() or MatXXXSetPreallocation() is called.
102f73d5cc4SBarry Smith 
103f69a0ea3SMatthew Knepley   Level: beginner
104f69a0ea3SMatthew Knepley 
105f69a0ea3SMatthew Knepley .seealso: MatGetSize(), PetscSplitOwnership()
106f69a0ea3SMatthew Knepley @*/
1077087cfbeSBarry Smith PetscErrorCode  MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N)
108f69a0ea3SMatthew Knepley {
109f69a0ea3SMatthew Knepley   PetscFunctionBegin;
1100700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
111b28e6c9fSJed Brown   if (M > 0) PetscValidLogicalCollectiveInt(A,M,4);
112b28e6c9fSJed Brown   if (N > 0) PetscValidLogicalCollectiveInt(A,N,5);
113e32f2f54SBarry 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);
114e32f2f54SBarry 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);
115e32f2f54SBarry 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);
116e32f2f54SBarry 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);
117d0f46423SBarry Smith   A->rmap->n = m;
118d0f46423SBarry Smith   A->cmap->n = n;
119d0f46423SBarry Smith   A->rmap->N = M;
120d0f46423SBarry Smith   A->cmap->N = N;
121f69a0ea3SMatthew Knepley   PetscFunctionReturn(0);
122f69a0ea3SMatthew Knepley }
123f69a0ea3SMatthew Knepley 
124f69a0ea3SMatthew Knepley #undef __FUNCT__
1254a2ae208SSatish Balay #define __FUNCT__ "MatSetFromOptions"
12605869f15SSatish Balay /*@
127273d9f13SBarry Smith    MatSetFromOptions - Creates a matrix where the type is determined
128273d9f13SBarry Smith    from the options database. Generates a parallel MPI matrix if the
129273d9f13SBarry Smith    communicator has more than one processor.  The default matrix type is
13069b1f4b7SBarry Smith    AIJ, using the routines MatCreateSeqAIJ() and MatCreateAIJ() if
1317e5f4302SBarry Smith    you do not select a type in the options database.
132273d9f13SBarry Smith 
133273d9f13SBarry Smith    Collective on Mat
134273d9f13SBarry Smith 
135273d9f13SBarry Smith    Input Parameter:
136273d9f13SBarry Smith .  A - the matrix
137273d9f13SBarry Smith 
138273d9f13SBarry Smith    Options Database Keys:
139273d9f13SBarry Smith +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
14069b1f4b7SBarry Smith .    -mat_type mpiaij   - AIJ type, uses MatCreateAIJ()
141273d9f13SBarry Smith .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
14269b1f4b7SBarry Smith .    -mat_type mpidense - dense type, uses MatCreateDense()
143273d9f13SBarry Smith .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
14469b1f4b7SBarry Smith -    -mat_type mpibaij  - block AIJ type, uses MatCreateBAIJ()
145273d9f13SBarry Smith 
146273d9f13SBarry Smith    Even More Options Database Keys:
147273d9f13SBarry Smith    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
148273d9f13SBarry Smith    for additional format-specific options.
149bd9ce289SLois Curfman McInnes 
1501d69843bSLois Curfman McInnes    Level: beginner
1511d69843bSLois Curfman McInnes 
152dc401e71SLois Curfman McInnes .keywords: matrix, create
153e0b365e2SLois Curfman McInnes 
15469b1f4b7SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateAIJ(),
15569b1f4b7SBarry Smith           MatCreateSeqDense(), MatCreateDense(),
15669b1f4b7SBarry Smith           MatCreateSeqBAIJ(), MatCreateBAIJ(),
15769b1f4b7SBarry Smith           MatCreateSeqSBAIJ(), MatCreateSBAIJ(),
158273d9f13SBarry Smith           MatConvert()
1597807a1faSBarry Smith @*/
1607087cfbeSBarry Smith PetscErrorCode  MatSetFromOptions(Mat B)
1617807a1faSBarry Smith {
162dfbe8321SBarry Smith   PetscErrorCode ierr;
163f3be49caSLisandro Dalcin   const char     *deft = MATAIJ;
164f3be49caSLisandro Dalcin   char           type[256];
16569df5c0cSJed Brown   PetscBool      flg,set;
166dbb450caSBarry Smith 
1673a40ed3dSBarry Smith   PetscFunctionBegin;
1680700a824SBarry Smith   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
169f3be49caSLisandro Dalcin 
1703194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)B);CHKERRQ(ierr);
171535b19f3SBarry Smith 
172535b19f3SBarry Smith   if (B->rmap->bs < 0) {
173535b19f3SBarry Smith     PetscInt newbs = -1;
174535b19f3SBarry Smith     ierr = PetscOptionsInt("-mat_block_size","Set the blocksize used to store the matrix","MatSetBlockSize",newbs,&newbs,&flg);CHKERRQ(ierr);
175535b19f3SBarry Smith     if (flg) {
176535b19f3SBarry Smith       ierr = PetscLayoutSetBlockSize(B->rmap,newbs);CHKERRQ(ierr);
177535b19f3SBarry Smith       ierr = PetscLayoutSetBlockSize(B->cmap,newbs);CHKERRQ(ierr);
178535b19f3SBarry Smith     }
179535b19f3SBarry Smith   }
180535b19f3SBarry Smith 
181f3be49caSLisandro Dalcin   ierr = PetscOptionsList("-mat_type","Matrix type","MatSetType",MatList,deft,type,256,&flg);CHKERRQ(ierr);
182273d9f13SBarry Smith   if (flg) {
183f3be49caSLisandro Dalcin     ierr = MatSetType(B,type);CHKERRQ(ierr);
184f3be49caSLisandro Dalcin   } else if (!((PetscObject)B)->type_name) {
185f3be49caSLisandro Dalcin     ierr = MatSetType(B,deft);CHKERRQ(ierr);
186273d9f13SBarry Smith   }
187f3be49caSLisandro Dalcin 
188840d65ccSBarry Smith   ierr = PetscViewerDestroy(&B->viewonassembly);CHKERRQ(ierr);
1890298fd71SBarry Smith   ierr = PetscOptionsViewer("-mat_view","Display mat with the viewer on MatAssemblyEnd()","MatView",&B->viewonassembly,&B->viewformatonassembly,NULL);CHKERRQ(ierr);
190840d65ccSBarry Smith   ierr = PetscOptionsName("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",&B->checksymmetryonassembly);CHKERRQ(ierr);
1910298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",0.0,&B->checksymmetrytol,NULL);CHKERRQ(ierr);
1920298fd71SBarry Smith   ierr = PetscOptionsBool("-mat_null_space_test","Checks if provided null space is correct in MatAssemblyEnd()","MatSetNullSpaceTest",PETSC_FALSE,&B->checknullspaceonassembly,NULL);CHKERRQ(ierr);
193840d65ccSBarry Smith 
194f3be49caSLisandro Dalcin   if (B->ops->setfromoptions) {
195f3be49caSLisandro Dalcin     ierr = (*B->ops->setfromoptions)(B);CHKERRQ(ierr);
196593c1905SSatish Balay   }
197f3be49caSLisandro Dalcin 
19869df5c0cSJed Brown   flg  = PETSC_FALSE;
19969df5c0cSJed 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);
20069df5c0cSJed Brown   if (set) {ierr = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,flg);CHKERRQ(ierr);}
20169df5c0cSJed Brown   flg  = PETSC_FALSE;
20269df5c0cSJed 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);
20369df5c0cSJed Brown   if (set) {ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,flg);CHKERRQ(ierr);}
20469df5c0cSJed Brown 
2055d973c19SBarry Smith   /* process any options handlers added with PetscObjectAddOptionsHandler() */
2065d973c19SBarry Smith   ierr = PetscObjectProcessOptionsHandlers((PetscObject)B);CHKERRQ(ierr);
207f3be49caSLisandro Dalcin   ierr = PetscOptionsEnd();CHKERRQ(ierr);
2083a40ed3dSBarry Smith   PetscFunctionReturn(0);
2097807a1faSBarry Smith }
2107807a1faSBarry Smith 
2114a2ae208SSatish Balay #undef __FUNCT__
21263562e91SJed Brown #define __FUNCT__ "MatXAIJSetPreallocation"
21363562e91SJed Brown /*@
21463562e91SJed Brown    MatXAIJSetPreallocation - set preallocation for serial and parallel AIJ, BAIJ, and SBAIJ matrices
21563562e91SJed Brown 
21663562e91SJed Brown    Collective on Mat
21763562e91SJed Brown 
21863562e91SJed Brown    Input Arguments:
21963562e91SJed Brown +  A - matrix being preallocated
22063562e91SJed Brown .  bs - block size
2213e5f4774SJed Brown .  dnnz - number of nonzero blocks per block row of diagonal part of parallel matrix
2223e5f4774SJed Brown .  onnz - number of nonzero blocks per block row of off-diagonal part of parallel matrix
2233e5f4774SJed Brown .  dnnzu - number of nonzero blocks per block row of upper-triangular part of diagonal part of parallel matrix
2243e5f4774SJed Brown -  onnzu - number of nonzero blocks per block row of upper-triangular part of off-diagonal part of parallel matrix
22563562e91SJed Brown 
22663562e91SJed Brown    Level: beginner
22763562e91SJed Brown 
228ab978733SBarry Smith .seealso: MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation(),
229ab978733SBarry Smith           PetscSplitOwnership()
23063562e91SJed Brown @*/
231*4cce697cSJed Brown PetscErrorCode MatXAIJSetPreallocation(Mat A,PetscInt bs,const PetscInt dnnz[],const PetscInt onnz[],const PetscInt dnnzu[],const PetscInt onnzu[])
23263562e91SJed Brown {
23363562e91SJed Brown   PetscErrorCode ierr;
23463562e91SJed Brown   void           (*aij)(void);
23563562e91SJed Brown 
23663562e91SJed Brown   PetscFunctionBegin;
237dec54756SJed Brown   ierr = MatSetBlockSize(A,bs);CHKERRQ(ierr);
238dec54756SJed Brown   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
239dec54756SJed Brown   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
2403e5f4774SJed Brown   ierr = MatSeqBAIJSetPreallocation(A,bs,0,dnnz);CHKERRQ(ierr);
2413e5f4774SJed Brown   ierr = MatMPIBAIJSetPreallocation(A,bs,0,dnnz,0,onnz);CHKERRQ(ierr);
2423e5f4774SJed Brown   ierr = MatSeqSBAIJSetPreallocation(A,bs,0,dnnzu);CHKERRQ(ierr);
2433e5f4774SJed Brown   ierr = MatMPISBAIJSetPreallocation(A,bs,0,dnnzu,0,onnzu);CHKERRQ(ierr);
24463562e91SJed Brown   /*
24563562e91SJed Brown     In general, we have to do extra work to preallocate for scalar (AIJ) matrices so we check whether it will do any
24663562e91SJed Brown     good before going on with it.
24763562e91SJed Brown   */
24863562e91SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
24963562e91SJed Brown   if (!aij) {
25063562e91SJed Brown     ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
25163562e91SJed Brown   }
25263562e91SJed Brown   if (aij) {
25363562e91SJed Brown     if (bs == 1) {
2543e5f4774SJed Brown       ierr = MatSeqAIJSetPreallocation(A,0,dnnz);CHKERRQ(ierr);
2553e5f4774SJed Brown       ierr = MatMPIAIJSetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr);
2563e5f4774SJed Brown     } else {                    /* Convert block-row precallocation to scalar-row */
25763562e91SJed Brown       PetscInt i,m,*sdnnz,*sonnz;
2580298fd71SBarry Smith       ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
259dec54756SJed Brown       ierr = PetscMalloc2((!!dnnz)*m,PetscInt,&sdnnz,(!!onnz)*m,PetscInt,&sonnz);CHKERRQ(ierr);
260dec54756SJed Brown       for (i=0; i<m; i++) {
26163562e91SJed Brown         if (dnnz) sdnnz[i] = dnnz[i/bs] * bs;
26263562e91SJed Brown         if (onnz) sonnz[i] = onnz[i/bs] * bs;
26363562e91SJed Brown       }
2640298fd71SBarry Smith       ierr = MatSeqAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL);CHKERRQ(ierr);
2650298fd71SBarry Smith       ierr = MatMPIAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr);
26663562e91SJed Brown       ierr = PetscFree2(sdnnz,sonnz);CHKERRQ(ierr);
26763562e91SJed Brown     }
26863562e91SJed Brown   }
26963562e91SJed Brown   PetscFunctionReturn(0);
27063562e91SJed Brown }
27163562e91SJed Brown 
272273d9f13SBarry Smith /*
273eb6b5d47SBarry Smith         Merges some information from Cs header to A; the C object is then destroyed
274d0f46423SBarry Smith 
275d0f46423SBarry Smith         This is somewhat different from MatHeaderReplace() it would be nice to merge the code
276273d9f13SBarry Smith */
2774a2ae208SSatish Balay #undef __FUNCT__
278eb6b5d47SBarry Smith #define __FUNCT__ "MatHeaderMerge"
279eb6b5d47SBarry Smith PetscErrorCode MatHeaderMerge(Mat A,Mat C)
280273d9f13SBarry Smith {
281dfbe8321SBarry Smith   PetscErrorCode ierr;
282d44834fbSBarry Smith   PetscInt       refct;
283273d9f13SBarry Smith   PetscOps       *Abops;
284273d9f13SBarry Smith   MatOps         Aops;
285273d9f13SBarry Smith   char           *mtype,*mname;
28630735b05SKris Buschelman   void           *spptr;
287273d9f13SBarry Smith 
288273d9f13SBarry Smith   PetscFunctionBegin;
289273d9f13SBarry Smith   /* save the parts of A we need */
2907adad957SLisandro Dalcin   Abops = ((PetscObject)A)->bops;
291273d9f13SBarry Smith   Aops  = A->ops;
2927adad957SLisandro Dalcin   refct = ((PetscObject)A)->refct;
2935c9eb25fSBarry Smith   mtype = ((PetscObject)A)->type_name;
2945c9eb25fSBarry Smith   mname = ((PetscObject)A)->name;
29530735b05SKris Buschelman   spptr = A->spptr;
29630735b05SKris Buschelman 
2975c9eb25fSBarry Smith   /* zero these so the destroy below does not free them */
2985c9eb25fSBarry Smith   ((PetscObject)A)->type_name = 0;
2995c9eb25fSBarry Smith   ((PetscObject)A)->name      = 0;
3005c9eb25fSBarry Smith 
3017c99f97cSSatish Balay   /* free all the interior data structures from mat */
3027c99f97cSSatish Balay   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
3037c99f97cSSatish Balay 
30430735b05SKris Buschelman   ierr = PetscFree(C->spptr);CHKERRQ(ierr);
30505b42c5fSBarry Smith 
3066bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->rmap);CHKERRQ(ierr);
3076bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->cmap);CHKERRQ(ierr);
308140e18c1SBarry Smith   ierr = PetscFunctionListDestroy(&((PetscObject)A)->qlist);CHKERRQ(ierr);
309140e18c1SBarry Smith   ierr = PetscObjectListDestroy(&((PetscObject)A)->olist);CHKERRQ(ierr);
310273d9f13SBarry Smith 
311273d9f13SBarry Smith   /* copy C over to A */
312273d9f13SBarry Smith   ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
313273d9f13SBarry Smith 
314273d9f13SBarry Smith   /* return the parts of A we saved */
3157adad957SLisandro Dalcin   ((PetscObject)A)->bops      = Abops;
316273d9f13SBarry Smith   A->ops                      = Aops;
3177adad957SLisandro Dalcin   ((PetscObject)A)->refct     = refct;
3187adad957SLisandro Dalcin   ((PetscObject)A)->type_name = mtype;
3197adad957SLisandro Dalcin   ((PetscObject)A)->name      = mname;
32030735b05SKris Buschelman   A->spptr                    = spptr;
321273d9f13SBarry Smith 
3225c9eb25fSBarry Smith   /* since these two are copied into A we do not want them destroyed in C */
3235c9eb25fSBarry Smith   ((PetscObject)C)->qlist = 0;
3245c9eb25fSBarry Smith   ((PetscObject)C)->olist = 0;
32526fbe8dcSKarl Rupp 
3266bf464f9SBarry Smith   ierr = PetscHeaderDestroy(&C);CHKERRQ(ierr);
327273d9f13SBarry Smith   PetscFunctionReturn(0);
328273d9f13SBarry Smith }
3298ab5b326SKris Buschelman /*
330eb6b5d47SBarry Smith         Replace A's header with that of C; the C object is then destroyed
331d0f46423SBarry Smith 
332eb6b5d47SBarry Smith         This is essentially code moved from MatDestroy()
333eb6b5d47SBarry Smith 
334eb6b5d47SBarry Smith         This is somewhat different from MatHeaderMerge() it would be nice to merge the code
335b30237c6SBarry Smith 
336b30237c6SBarry Smith         Used in DM hence is declared PETSC_EXTERN
3378ab5b326SKris Buschelman */
3388ab5b326SKris Buschelman #undef __FUNCT__
3398ab5b326SKris Buschelman #define __FUNCT__ "MatHeaderReplace"
340b30237c6SBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat A,Mat C)
3418ab5b326SKris Buschelman {
3428ab5b326SKris Buschelman   PetscErrorCode ierr;
34327b31e29SJed Brown   PetscInt       refct;
3448ab5b326SKris Buschelman 
3458ab5b326SKris Buschelman   PetscFunctionBegin;
34627b31e29SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
34727b31e29SJed Brown   PetscValidHeaderSpecific(C,MAT_CLASSID,2);
3486d7c1e57SBarry Smith   if (A == C) PetscFunctionReturn(0);
349a7e80f87SJed Brown   PetscCheckSameComm(A,1,C,2);
350ce94432eSBarry Smith   if (((PetscObject)C)->refct != 1) SETERRQ1(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"Object C has refct %D > 1, would leave hanging reference",((PetscObject)C)->refct);
3516d7c1e57SBarry Smith 
3528ab5b326SKris Buschelman   /* free all the interior data structures from mat */
3538ab5b326SKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
354d2c95936SBarry Smith   ierr = PetscHeaderDestroy_Private((PetscObject)A);CHKERRQ(ierr);
355c650bb5fSLisandro Dalcin   ierr = PetscFree(A->ops);CHKERRQ(ierr);
3566bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->rmap);CHKERRQ(ierr);
3576bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->cmap);CHKERRQ(ierr);
358d95db058SSatish Balay   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
3598ab5b326SKris Buschelman 
3608ab5b326SKris Buschelman   /* copy C over to A */
36127b31e29SJed Brown   refct = ((PetscObject)A)->refct;
3628ab5b326SKris Buschelman   ierr  = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
36326fbe8dcSKarl Rupp 
36427b31e29SJed Brown   ((PetscObject)A)->refct = refct;
36526fbe8dcSKarl Rupp 
3668ab5b326SKris Buschelman   ierr = PetscLogObjectDestroy((PetscObject)C);CHKERRQ(ierr);
3678ab5b326SKris Buschelman   ierr = PetscFree(C);CHKERRQ(ierr);
3688ab5b326SKris Buschelman   PetscFunctionReturn(0);
3698ab5b326SKris Buschelman }
370