xref: /petsc/src/mat/utils/gcreate.c (revision 69bbac976351bc191e2527ca4ac8ff0093de2015)
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;
60607a6623SBarry Smith   ierr = MatInitializePackage();CHKERRQ(ierr);
618ba1e511SMatthew Knepley 
6267c2884eSBarry Smith   ierr = PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_CLASSID,"Mat","Matrix","Mat",comm,MatDestroy,MatView);CHKERRQ(ierr);
6326283091SBarry Smith   ierr = PetscLayoutCreate(comm,&B->rmap);CHKERRQ(ierr);
6426283091SBarry Smith   ierr = PetscLayoutCreate(comm,&B->cmap);CHKERRQ(ierr);
6526fbe8dcSKarl Rupp 
66273d9f13SBarry Smith   B->preallocated = PETSC_FALSE;
67273d9f13SBarry Smith   *A              = B;
68273d9f13SBarry Smith   PetscFunctionReturn(0);
69273d9f13SBarry Smith }
70273d9f13SBarry Smith 
714a2ae208SSatish Balay #undef __FUNCT__
72f69a0ea3SMatthew Knepley #define __FUNCT__ "MatSetSizes"
73f69a0ea3SMatthew Knepley /*@
74f69a0ea3SMatthew Knepley   MatSetSizes - Sets the local and global sizes, and checks to determine compatibility
75f69a0ea3SMatthew Knepley 
76f69a0ea3SMatthew Knepley   Collective on Mat
77f69a0ea3SMatthew Knepley 
78f69a0ea3SMatthew Knepley   Input Parameters:
79f69a0ea3SMatthew Knepley +  A - the matrix
80f69a0ea3SMatthew Knepley .  m - number of local rows (or PETSC_DECIDE)
81f69a0ea3SMatthew Knepley .  n - number of local columns (or PETSC_DECIDE)
82f69a0ea3SMatthew Knepley .  M - number of global rows (or PETSC_DETERMINE)
83f69a0ea3SMatthew Knepley -  N - number of global columns (or PETSC_DETERMINE)
84f69a0ea3SMatthew Knepley 
85f69a0ea3SMatthew Knepley    Notes:
86f69a0ea3SMatthew Knepley    m (n) and M (N) cannot be both PETSC_DECIDE
87f69a0ea3SMatthew Knepley    If one processor calls this with M (N) of PETSC_DECIDE then all processors must, otherwise the program will hang.
88f69a0ea3SMatthew Knepley 
89f69a0ea3SMatthew Knepley    If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the
90f69a0ea3SMatthew Knepley    user must ensure that they are chosen to be compatible with the
91f69a0ea3SMatthew Knepley    vectors. To do this, one first considers the matrix-vector product
92f69a0ea3SMatthew Knepley    'y = A x'. The 'm' that is used in the above routine must match the
93f69a0ea3SMatthew Knepley    local size used in the vector creation routine VecCreateMPI() for 'y'.
94f69a0ea3SMatthew Knepley    Likewise, the 'n' used must match that used as the local size in
95f69a0ea3SMatthew Knepley    VecCreateMPI() for 'x'.
96f69a0ea3SMatthew Knepley 
97f73d5cc4SBarry Smith    You cannot change the sizes once they have been set.
98f73d5cc4SBarry Smith 
99f73d5cc4SBarry Smith    The sizes must be set before MatSetUp() or MatXXXSetPreallocation() is called.
100f73d5cc4SBarry Smith 
101f69a0ea3SMatthew Knepley   Level: beginner
102f69a0ea3SMatthew Knepley 
103f69a0ea3SMatthew Knepley .seealso: MatGetSize(), PetscSplitOwnership()
104f69a0ea3SMatthew Knepley @*/
1057087cfbeSBarry Smith PetscErrorCode  MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N)
106f69a0ea3SMatthew Knepley {
107f69a0ea3SMatthew Knepley   PetscFunctionBegin;
1080700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
109b28e6c9fSJed Brown   if (M > 0) PetscValidLogicalCollectiveInt(A,M,4);
110b28e6c9fSJed Brown   if (N > 0) PetscValidLogicalCollectiveInt(A,N,5);
111e32f2f54SBarry 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);
112e32f2f54SBarry 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);
113e32f2f54SBarry 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);
114e32f2f54SBarry 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);
115d0f46423SBarry Smith   A->rmap->n = m;
116d0f46423SBarry Smith   A->cmap->n = n;
117d0f46423SBarry Smith   A->rmap->N = M;
118d0f46423SBarry Smith   A->cmap->N = N;
119f69a0ea3SMatthew Knepley   PetscFunctionReturn(0);
120f69a0ea3SMatthew Knepley }
121f69a0ea3SMatthew Knepley 
122f69a0ea3SMatthew Knepley #undef __FUNCT__
1234a2ae208SSatish Balay #define __FUNCT__ "MatSetFromOptions"
12405869f15SSatish Balay /*@
125273d9f13SBarry Smith    MatSetFromOptions - Creates a matrix where the type is determined
126273d9f13SBarry Smith    from the options database. Generates a parallel MPI matrix if the
127273d9f13SBarry Smith    communicator has more than one processor.  The default matrix type is
12869b1f4b7SBarry Smith    AIJ, using the routines MatCreateSeqAIJ() and MatCreateAIJ() if
1297e5f4302SBarry Smith    you do not select a type in the options database.
130273d9f13SBarry Smith 
131273d9f13SBarry Smith    Collective on Mat
132273d9f13SBarry Smith 
133273d9f13SBarry Smith    Input Parameter:
134273d9f13SBarry Smith .  A - the matrix
135273d9f13SBarry Smith 
136273d9f13SBarry Smith    Options Database Keys:
137273d9f13SBarry Smith +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
13869b1f4b7SBarry Smith .    -mat_type mpiaij   - AIJ type, uses MatCreateAIJ()
139273d9f13SBarry Smith .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
14069b1f4b7SBarry Smith .    -mat_type mpidense - dense type, uses MatCreateDense()
141273d9f13SBarry Smith .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
14269b1f4b7SBarry Smith -    -mat_type mpibaij  - block AIJ type, uses MatCreateBAIJ()
143273d9f13SBarry Smith 
144273d9f13SBarry Smith    Even More Options Database Keys:
145273d9f13SBarry Smith    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
146273d9f13SBarry Smith    for additional format-specific options.
147bd9ce289SLois Curfman McInnes 
1481d69843bSLois Curfman McInnes    Level: beginner
1491d69843bSLois Curfman McInnes 
150dc401e71SLois Curfman McInnes .keywords: matrix, create
151e0b365e2SLois Curfman McInnes 
15269b1f4b7SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateAIJ(),
15369b1f4b7SBarry Smith           MatCreateSeqDense(), MatCreateDense(),
15469b1f4b7SBarry Smith           MatCreateSeqBAIJ(), MatCreateBAIJ(),
15569b1f4b7SBarry Smith           MatCreateSeqSBAIJ(), MatCreateSBAIJ(),
156273d9f13SBarry Smith           MatConvert()
1577807a1faSBarry Smith @*/
1587087cfbeSBarry Smith PetscErrorCode  MatSetFromOptions(Mat B)
1597807a1faSBarry Smith {
160dfbe8321SBarry Smith   PetscErrorCode ierr;
161f3be49caSLisandro Dalcin   const char     *deft = MATAIJ;
162f3be49caSLisandro Dalcin   char           type[256];
16369df5c0cSJed Brown   PetscBool      flg,set;
164dbb450caSBarry Smith 
1653a40ed3dSBarry Smith   PetscFunctionBegin;
1660700a824SBarry Smith   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
167f3be49caSLisandro Dalcin 
1683194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)B);CHKERRQ(ierr);
169535b19f3SBarry Smith 
170535b19f3SBarry Smith   if (B->rmap->bs < 0) {
171535b19f3SBarry Smith     PetscInt newbs = -1;
172535b19f3SBarry Smith     ierr = PetscOptionsInt("-mat_block_size","Set the blocksize used to store the matrix","MatSetBlockSize",newbs,&newbs,&flg);CHKERRQ(ierr);
173535b19f3SBarry Smith     if (flg) {
174535b19f3SBarry Smith       ierr = PetscLayoutSetBlockSize(B->rmap,newbs);CHKERRQ(ierr);
175535b19f3SBarry Smith       ierr = PetscLayoutSetBlockSize(B->cmap,newbs);CHKERRQ(ierr);
176535b19f3SBarry Smith     }
177535b19f3SBarry Smith   }
178535b19f3SBarry Smith 
179a264d7a6SBarry Smith   ierr = PetscOptionsFList("-mat_type","Matrix type","MatSetType",MatList,deft,type,256,&flg);CHKERRQ(ierr);
180273d9f13SBarry Smith   if (flg) {
181f3be49caSLisandro Dalcin     ierr = MatSetType(B,type);CHKERRQ(ierr);
182f3be49caSLisandro Dalcin   } else if (!((PetscObject)B)->type_name) {
183f3be49caSLisandro Dalcin     ierr = MatSetType(B,deft);CHKERRQ(ierr);
184273d9f13SBarry Smith   }
185f3be49caSLisandro Dalcin 
186840d65ccSBarry Smith   ierr = PetscOptionsName("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",&B->checksymmetryonassembly);CHKERRQ(ierr);
1870298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",0.0,&B->checksymmetrytol,NULL);CHKERRQ(ierr);
1880298fd71SBarry Smith   ierr = PetscOptionsBool("-mat_null_space_test","Checks if provided null space is correct in MatAssemblyEnd()","MatSetNullSpaceTest",PETSC_FALSE,&B->checknullspaceonassembly,NULL);CHKERRQ(ierr);
189840d65ccSBarry Smith 
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);
2043a40ed3dSBarry Smith   PetscFunctionReturn(0);
2057807a1faSBarry Smith }
2067807a1faSBarry Smith 
2074a2ae208SSatish Balay #undef __FUNCT__
20863562e91SJed Brown #define __FUNCT__ "MatXAIJSetPreallocation"
20963562e91SJed Brown /*@
21063562e91SJed Brown    MatXAIJSetPreallocation - set preallocation for serial and parallel AIJ, BAIJ, and SBAIJ matrices
21163562e91SJed Brown 
21263562e91SJed Brown    Collective on Mat
21363562e91SJed Brown 
21463562e91SJed Brown    Input Arguments:
21563562e91SJed Brown +  A - matrix being preallocated
21663562e91SJed Brown .  bs - block size
2173e5f4774SJed Brown .  dnnz - number of nonzero blocks per block row of diagonal part of parallel matrix
2183e5f4774SJed Brown .  onnz - number of nonzero blocks per block row of off-diagonal part of parallel matrix
2193e5f4774SJed Brown .  dnnzu - number of nonzero blocks per block row of upper-triangular part of diagonal part of parallel matrix
2203e5f4774SJed Brown -  onnzu - number of nonzero blocks per block row of upper-triangular part of off-diagonal part of parallel matrix
22163562e91SJed Brown 
22263562e91SJed Brown    Level: beginner
22363562e91SJed Brown 
224ab978733SBarry Smith .seealso: MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation(),
225ab978733SBarry Smith           PetscSplitOwnership()
22663562e91SJed Brown @*/
2274cce697cSJed Brown PetscErrorCode MatXAIJSetPreallocation(Mat A,PetscInt bs,const PetscInt dnnz[],const PetscInt onnz[],const PetscInt dnnzu[],const PetscInt onnzu[])
22863562e91SJed Brown {
22963562e91SJed Brown   PetscErrorCode ierr;
23063562e91SJed Brown   void           (*aij)(void);
23163562e91SJed Brown 
23263562e91SJed Brown   PetscFunctionBegin;
233dec54756SJed Brown   ierr = MatSetBlockSize(A,bs);CHKERRQ(ierr);
234*69bbac97SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
235dec54756SJed Brown   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
236dec54756SJed Brown   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
2373e5f4774SJed Brown   ierr = MatSeqBAIJSetPreallocation(A,bs,0,dnnz);CHKERRQ(ierr);
2383e5f4774SJed Brown   ierr = MatMPIBAIJSetPreallocation(A,bs,0,dnnz,0,onnz);CHKERRQ(ierr);
2393e5f4774SJed Brown   ierr = MatSeqSBAIJSetPreallocation(A,bs,0,dnnzu);CHKERRQ(ierr);
2403e5f4774SJed Brown   ierr = MatMPISBAIJSetPreallocation(A,bs,0,dnnzu,0,onnzu);CHKERRQ(ierr);
24163562e91SJed Brown   /*
24263562e91SJed Brown     In general, we have to do extra work to preallocate for scalar (AIJ) matrices so we check whether it will do any
24363562e91SJed Brown     good before going on with it.
24463562e91SJed Brown   */
24563562e91SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
24663562e91SJed Brown   if (!aij) {
24763562e91SJed Brown     ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
24863562e91SJed Brown   }
24963562e91SJed Brown   if (aij) {
25063562e91SJed Brown     if (bs == 1) {
2513e5f4774SJed Brown       ierr = MatSeqAIJSetPreallocation(A,0,dnnz);CHKERRQ(ierr);
2523e5f4774SJed Brown       ierr = MatMPIAIJSetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr);
2533e5f4774SJed Brown     } else {                    /* Convert block-row precallocation to scalar-row */
25463562e91SJed Brown       PetscInt i,m,*sdnnz,*sonnz;
2550298fd71SBarry Smith       ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
256dcca6d9dSJed Brown       ierr = PetscMalloc2((!!dnnz)*m,&sdnnz,(!!onnz)*m,&sonnz);CHKERRQ(ierr);
257dec54756SJed Brown       for (i=0; i<m; i++) {
25863562e91SJed Brown         if (dnnz) sdnnz[i] = dnnz[i/bs] * bs;
25963562e91SJed Brown         if (onnz) sonnz[i] = onnz[i/bs] * bs;
26063562e91SJed Brown       }
2610298fd71SBarry Smith       ierr = MatSeqAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL);CHKERRQ(ierr);
2620298fd71SBarry Smith       ierr = MatMPIAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr);
26363562e91SJed Brown       ierr = PetscFree2(sdnnz,sonnz);CHKERRQ(ierr);
26463562e91SJed Brown     }
26563562e91SJed Brown   }
26663562e91SJed Brown   PetscFunctionReturn(0);
26763562e91SJed Brown }
26863562e91SJed Brown 
269273d9f13SBarry Smith /*
270eb6b5d47SBarry Smith         Merges some information from Cs header to A; the C object is then destroyed
271d0f46423SBarry Smith 
272d0f46423SBarry Smith         This is somewhat different from MatHeaderReplace() it would be nice to merge the code
273273d9f13SBarry Smith */
2744a2ae208SSatish Balay #undef __FUNCT__
275eb6b5d47SBarry Smith #define __FUNCT__ "MatHeaderMerge"
276eb6b5d47SBarry Smith PetscErrorCode MatHeaderMerge(Mat A,Mat C)
277273d9f13SBarry Smith {
278dfbe8321SBarry Smith   PetscErrorCode ierr;
279d44834fbSBarry Smith   PetscInt       refct;
280273d9f13SBarry Smith   PetscOps       *Abops;
281273d9f13SBarry Smith   MatOps         Aops;
282273d9f13SBarry Smith   char           *mtype,*mname;
28330735b05SKris Buschelman   void           *spptr;
284273d9f13SBarry Smith 
285273d9f13SBarry Smith   PetscFunctionBegin;
286273d9f13SBarry Smith   /* save the parts of A we need */
2877adad957SLisandro Dalcin   Abops = ((PetscObject)A)->bops;
288273d9f13SBarry Smith   Aops  = A->ops;
2897adad957SLisandro Dalcin   refct = ((PetscObject)A)->refct;
2905c9eb25fSBarry Smith   mtype = ((PetscObject)A)->type_name;
2915c9eb25fSBarry Smith   mname = ((PetscObject)A)->name;
29230735b05SKris Buschelman   spptr = A->spptr;
29330735b05SKris Buschelman 
2945c9eb25fSBarry Smith   /* zero these so the destroy below does not free them */
2955c9eb25fSBarry Smith   ((PetscObject)A)->type_name = 0;
2965c9eb25fSBarry Smith   ((PetscObject)A)->name      = 0;
2975c9eb25fSBarry Smith 
2987c99f97cSSatish Balay   /* free all the interior data structures from mat */
2997c99f97cSSatish Balay   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
3007c99f97cSSatish Balay 
30130735b05SKris Buschelman   ierr = PetscFree(C->spptr);CHKERRQ(ierr);
30205b42c5fSBarry Smith 
3036bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->rmap);CHKERRQ(ierr);
3046bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->cmap);CHKERRQ(ierr);
305140e18c1SBarry Smith   ierr = PetscFunctionListDestroy(&((PetscObject)A)->qlist);CHKERRQ(ierr);
306140e18c1SBarry Smith   ierr = PetscObjectListDestroy(&((PetscObject)A)->olist);CHKERRQ(ierr);
307273d9f13SBarry Smith 
308273d9f13SBarry Smith   /* copy C over to A */
309273d9f13SBarry Smith   ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
310273d9f13SBarry Smith 
311273d9f13SBarry Smith   /* return the parts of A we saved */
3127adad957SLisandro Dalcin   ((PetscObject)A)->bops      = Abops;
313273d9f13SBarry Smith   A->ops                      = Aops;
3147adad957SLisandro Dalcin   ((PetscObject)A)->refct     = refct;
3157adad957SLisandro Dalcin   ((PetscObject)A)->type_name = mtype;
3167adad957SLisandro Dalcin   ((PetscObject)A)->name      = mname;
31730735b05SKris Buschelman   A->spptr                    = spptr;
318273d9f13SBarry Smith 
3195c9eb25fSBarry Smith   /* since these two are copied into A we do not want them destroyed in C */
3205c9eb25fSBarry Smith   ((PetscObject)C)->qlist = 0;
3215c9eb25fSBarry Smith   ((PetscObject)C)->olist = 0;
32226fbe8dcSKarl Rupp 
3236bf464f9SBarry Smith   ierr = PetscHeaderDestroy(&C);CHKERRQ(ierr);
324273d9f13SBarry Smith   PetscFunctionReturn(0);
325273d9f13SBarry Smith }
3268ab5b326SKris Buschelman /*
327eb6b5d47SBarry Smith         Replace A's header with that of C; the C object is then destroyed
328d0f46423SBarry Smith 
329eb6b5d47SBarry Smith         This is essentially code moved from MatDestroy()
330eb6b5d47SBarry Smith 
331eb6b5d47SBarry Smith         This is somewhat different from MatHeaderMerge() it would be nice to merge the code
332b30237c6SBarry Smith 
333b30237c6SBarry Smith         Used in DM hence is declared PETSC_EXTERN
3348ab5b326SKris Buschelman */
3358ab5b326SKris Buschelman #undef __FUNCT__
3368ab5b326SKris Buschelman #define __FUNCT__ "MatHeaderReplace"
337b30237c6SBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat A,Mat C)
3388ab5b326SKris Buschelman {
3398ab5b326SKris Buschelman   PetscErrorCode ierr;
34027b31e29SJed Brown   PetscInt       refct;
3418ab5b326SKris Buschelman 
3428ab5b326SKris Buschelman   PetscFunctionBegin;
34327b31e29SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
34427b31e29SJed Brown   PetscValidHeaderSpecific(C,MAT_CLASSID,2);
3456d7c1e57SBarry Smith   if (A == C) PetscFunctionReturn(0);
346a7e80f87SJed Brown   PetscCheckSameComm(A,1,C,2);
347ce94432eSBarry 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);
3486d7c1e57SBarry Smith 
3498ab5b326SKris Buschelman   /* free all the interior data structures from mat */
3508ab5b326SKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
351d2c95936SBarry Smith   ierr = PetscHeaderDestroy_Private((PetscObject)A);CHKERRQ(ierr);
352c650bb5fSLisandro Dalcin   ierr = PetscFree(A->ops);CHKERRQ(ierr);
3536bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->rmap);CHKERRQ(ierr);
3546bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->cmap);CHKERRQ(ierr);
355d95db058SSatish Balay   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
3568ab5b326SKris Buschelman 
3578ab5b326SKris Buschelman   /* copy C over to A */
35827b31e29SJed Brown   refct = ((PetscObject)A)->refct;
3598ab5b326SKris Buschelman   ierr  = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
36026fbe8dcSKarl Rupp 
36127b31e29SJed Brown   ((PetscObject)A)->refct = refct;
36226fbe8dcSKarl Rupp 
3638ab5b326SKris Buschelman   ierr = PetscFree(C);CHKERRQ(ierr);
3648ab5b326SKris Buschelman   PetscFunctionReturn(0);
3658ab5b326SKris Buschelman }
366