xref: /petsc/src/mat/utils/gcreate.c (revision 26fbe8dc1a3c99fb8dddfa572c8c6b3b4ce3ca53)
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;
70519f805aSKarl Rupp #if !defined(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);
77*26fbe8dcSKarl Rupp 
78273d9f13SBarry Smith   B->preallocated = PETSC_FALSE;
79273d9f13SBarry Smith   *A              = B;
80273d9f13SBarry Smith   PetscFunctionReturn(0);
81273d9f13SBarry Smith }
82273d9f13SBarry Smith 
834a2ae208SSatish Balay #undef __FUNCT__
84f69a0ea3SMatthew Knepley #define __FUNCT__ "MatSetSizes"
85f69a0ea3SMatthew Knepley /*@
86f69a0ea3SMatthew Knepley   MatSetSizes - Sets the local and global sizes, and checks to determine compatibility
87f69a0ea3SMatthew Knepley 
88f69a0ea3SMatthew Knepley   Collective on Mat
89f69a0ea3SMatthew Knepley 
90f69a0ea3SMatthew Knepley   Input Parameters:
91f69a0ea3SMatthew Knepley +  A - the matrix
92f69a0ea3SMatthew Knepley .  m - number of local rows (or PETSC_DECIDE)
93f69a0ea3SMatthew Knepley .  n - number of local columns (or PETSC_DECIDE)
94f69a0ea3SMatthew Knepley .  M - number of global rows (or PETSC_DETERMINE)
95f69a0ea3SMatthew Knepley -  N - number of global columns (or PETSC_DETERMINE)
96f69a0ea3SMatthew Knepley 
97f69a0ea3SMatthew Knepley    Notes:
98f69a0ea3SMatthew Knepley    m (n) and M (N) cannot be both PETSC_DECIDE
99f69a0ea3SMatthew Knepley    If one processor calls this with M (N) of PETSC_DECIDE then all processors must, otherwise the program will hang.
100f69a0ea3SMatthew Knepley 
101f69a0ea3SMatthew Knepley    If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the
102f69a0ea3SMatthew Knepley    user must ensure that they are chosen to be compatible with the
103f69a0ea3SMatthew Knepley    vectors. To do this, one first considers the matrix-vector product
104f69a0ea3SMatthew Knepley    'y = A x'. The 'm' that is used in the above routine must match the
105f69a0ea3SMatthew Knepley    local size used in the vector creation routine VecCreateMPI() for 'y'.
106f69a0ea3SMatthew Knepley    Likewise, the 'n' used must match that used as the local size in
107f69a0ea3SMatthew Knepley    VecCreateMPI() for 'x'.
108f69a0ea3SMatthew Knepley 
109f73d5cc4SBarry Smith    You cannot change the sizes once they have been set.
110f73d5cc4SBarry Smith 
111f73d5cc4SBarry Smith    The sizes must be set before MatSetUp() or MatXXXSetPreallocation() is called.
112f73d5cc4SBarry Smith 
113f69a0ea3SMatthew Knepley   Level: beginner
114f69a0ea3SMatthew Knepley 
115f69a0ea3SMatthew Knepley .seealso: MatGetSize(), PetscSplitOwnership()
116f69a0ea3SMatthew Knepley @*/
1177087cfbeSBarry Smith PetscErrorCode  MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N)
118f69a0ea3SMatthew Knepley {
119f69a0ea3SMatthew Knepley   PetscFunctionBegin;
1200700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
121e32f2f54SBarry 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);
122e32f2f54SBarry 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);
123e32f2f54SBarry 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);
124e32f2f54SBarry 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);
125d0f46423SBarry Smith   A->rmap->n = m;
126d0f46423SBarry Smith   A->cmap->n = n;
127d0f46423SBarry Smith   A->rmap->N = M;
128d0f46423SBarry Smith   A->cmap->N = N;
129f69a0ea3SMatthew Knepley   PetscFunctionReturn(0);
130f69a0ea3SMatthew Knepley }
131f69a0ea3SMatthew Knepley 
132f69a0ea3SMatthew Knepley #undef __FUNCT__
1334a2ae208SSatish Balay #define __FUNCT__ "MatSetFromOptions"
13405869f15SSatish Balay /*@
135273d9f13SBarry Smith    MatSetFromOptions - Creates a matrix where the type is determined
136273d9f13SBarry Smith    from the options database. Generates a parallel MPI matrix if the
137273d9f13SBarry Smith    communicator has more than one processor.  The default matrix type is
13869b1f4b7SBarry Smith    AIJ, using the routines MatCreateSeqAIJ() and MatCreateAIJ() if
1397e5f4302SBarry Smith    you do not select a type in the options database.
140273d9f13SBarry Smith 
141273d9f13SBarry Smith    Collective on Mat
142273d9f13SBarry Smith 
143273d9f13SBarry Smith    Input Parameter:
144273d9f13SBarry Smith .  A - the matrix
145273d9f13SBarry Smith 
146273d9f13SBarry Smith    Options Database Keys:
147273d9f13SBarry Smith +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
14869b1f4b7SBarry Smith .    -mat_type mpiaij   - AIJ type, uses MatCreateAIJ()
149273d9f13SBarry Smith .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
15069b1f4b7SBarry Smith .    -mat_type mpidense - dense type, uses MatCreateDense()
151273d9f13SBarry Smith .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
15269b1f4b7SBarry Smith -    -mat_type mpibaij  - block AIJ type, uses MatCreateBAIJ()
153273d9f13SBarry Smith 
154273d9f13SBarry Smith    Even More Options Database Keys:
155273d9f13SBarry Smith    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
156273d9f13SBarry Smith    for additional format-specific options.
157bd9ce289SLois Curfman McInnes 
1581d69843bSLois Curfman McInnes    Level: beginner
1591d69843bSLois Curfman McInnes 
160dc401e71SLois Curfman McInnes .keywords: matrix, create
161e0b365e2SLois Curfman McInnes 
16269b1f4b7SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateAIJ(),
16369b1f4b7SBarry Smith           MatCreateSeqDense(), MatCreateDense(),
16469b1f4b7SBarry Smith           MatCreateSeqBAIJ(), MatCreateBAIJ(),
16569b1f4b7SBarry Smith           MatCreateSeqSBAIJ(), MatCreateSBAIJ(),
166273d9f13SBarry Smith           MatConvert()
1677807a1faSBarry Smith @*/
1687087cfbeSBarry Smith PetscErrorCode  MatSetFromOptions(Mat B)
1697807a1faSBarry Smith {
170dfbe8321SBarry Smith   PetscErrorCode ierr;
171f3be49caSLisandro Dalcin   const char     *deft = MATAIJ;
172f3be49caSLisandro Dalcin   char           type[256];
17369df5c0cSJed Brown   PetscBool      flg,set;
174dbb450caSBarry Smith 
1753a40ed3dSBarry Smith   PetscFunctionBegin;
1760700a824SBarry Smith   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
177f3be49caSLisandro Dalcin 
1783194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)B);CHKERRQ(ierr);
179535b19f3SBarry Smith 
180535b19f3SBarry Smith   if (B->rmap->bs < 0) {
181535b19f3SBarry Smith     PetscInt newbs = -1;
182535b19f3SBarry Smith     ierr = PetscOptionsInt("-mat_block_size","Set the blocksize used to store the matrix","MatSetBlockSize",newbs,&newbs,&flg);CHKERRQ(ierr);
183535b19f3SBarry Smith     if (flg) {
184535b19f3SBarry Smith       ierr = PetscLayoutSetBlockSize(B->rmap,newbs);CHKERRQ(ierr);
185535b19f3SBarry Smith       ierr = PetscLayoutSetBlockSize(B->cmap,newbs);CHKERRQ(ierr);
186535b19f3SBarry Smith     }
187535b19f3SBarry Smith   }
188535b19f3SBarry Smith 
189f3be49caSLisandro Dalcin   ierr = PetscOptionsList("-mat_type","Matrix type","MatSetType",MatList,deft,type,256,&flg);CHKERRQ(ierr);
190273d9f13SBarry Smith   if (flg) {
191f3be49caSLisandro Dalcin     ierr = MatSetType(B,type);CHKERRQ(ierr);
192f3be49caSLisandro Dalcin   } else if (!((PetscObject)B)->type_name) {
193f3be49caSLisandro Dalcin     ierr = MatSetType(B,deft);CHKERRQ(ierr);
194273d9f13SBarry Smith   }
195f3be49caSLisandro Dalcin 
196840d65ccSBarry Smith   ierr = PetscViewerDestroy(&B->viewonassembly);CHKERRQ(ierr);
1973b52ad7dSJed Brown   ierr = PetscOptionsViewer("-mat_view","Display mat with the viewer on MatAssemblyEnd()","MatView",&B->viewonassembly,&B->viewformatonassembly,PETSC_NULL);CHKERRQ(ierr);
198840d65ccSBarry Smith   ierr = PetscOptionsName("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",&B->checksymmetryonassembly);CHKERRQ(ierr);
199840d65ccSBarry Smith   ierr = PetscOptionsReal("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",0.0,&B->checksymmetrytol,PETSC_NULL);CHKERRQ(ierr);
200840d65ccSBarry Smith   ierr = PetscOptionsBool("-mat_null_space_test","Checks if provided null space is correct in MatAssemblyEnd()","MatSetNullSpaceTest",PETSC_FALSE,&B->checknullspaceonassembly,PETSC_NULL);CHKERRQ(ierr);
201840d65ccSBarry Smith 
202f3be49caSLisandro Dalcin   if (B->ops->setfromoptions) {
203f3be49caSLisandro Dalcin     ierr = (*B->ops->setfromoptions)(B);CHKERRQ(ierr);
204593c1905SSatish Balay   }
205f3be49caSLisandro Dalcin 
20669df5c0cSJed Brown   flg  = PETSC_FALSE;
20769df5c0cSJed 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);
20869df5c0cSJed Brown   if (set) {ierr = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,flg);CHKERRQ(ierr);}
20969df5c0cSJed Brown   flg  = PETSC_FALSE;
21069df5c0cSJed 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);
21169df5c0cSJed Brown   if (set) {ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,flg);CHKERRQ(ierr);}
21269df5c0cSJed Brown 
2135d973c19SBarry Smith   /* process any options handlers added with PetscObjectAddOptionsHandler() */
2145d973c19SBarry Smith   ierr = PetscObjectProcessOptionsHandlers((PetscObject)B);CHKERRQ(ierr);
215f3be49caSLisandro Dalcin   ierr = PetscOptionsEnd();CHKERRQ(ierr);
2163a40ed3dSBarry Smith   PetscFunctionReturn(0);
2177807a1faSBarry Smith }
2187807a1faSBarry Smith 
2194a2ae208SSatish Balay #undef __FUNCT__
22063562e91SJed Brown #define __FUNCT__ "MatXAIJSetPreallocation"
22163562e91SJed Brown /*@
22263562e91SJed Brown    MatXAIJSetPreallocation - set preallocation for serial and parallel AIJ, BAIJ, and SBAIJ matrices
22363562e91SJed Brown 
22463562e91SJed Brown    Collective on Mat
22563562e91SJed Brown 
22663562e91SJed Brown    Input Arguments:
22763562e91SJed Brown +  A - matrix being preallocated
22863562e91SJed Brown .  bs - block size
2293e5f4774SJed Brown .  dnnz - number of nonzero blocks per block row of diagonal part of parallel matrix
2303e5f4774SJed Brown .  onnz - number of nonzero blocks per block row of off-diagonal part of parallel matrix
2313e5f4774SJed Brown .  dnnzu - number of nonzero blocks per block row of upper-triangular part of diagonal part of parallel matrix
2323e5f4774SJed Brown -  onnzu - number of nonzero blocks per block row of upper-triangular part of off-diagonal part of parallel matrix
23363562e91SJed Brown 
23463562e91SJed Brown    Level: beginner
23563562e91SJed Brown 
236ab978733SBarry Smith .seealso: MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation(),
237ab978733SBarry Smith           PetscSplitOwnership()
23863562e91SJed Brown @*/
2393e5f4774SJed Brown PetscErrorCode MatXAIJSetPreallocation(Mat A,PetscInt bs,const PetscInt *dnnz,const PetscInt *onnz,const PetscInt *dnnzu,const PetscInt *onnzu)
24063562e91SJed Brown {
24163562e91SJed Brown   PetscErrorCode ierr;
24263562e91SJed Brown   void           (*aij)(void);
24363562e91SJed Brown 
24463562e91SJed Brown   PetscFunctionBegin;
245dec54756SJed Brown   ierr = MatSetBlockSize(A,bs);CHKERRQ(ierr);
246dec54756SJed Brown   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
247dec54756SJed Brown   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
2483e5f4774SJed Brown   ierr = MatSeqBAIJSetPreallocation(A,bs,0,dnnz);CHKERRQ(ierr);
2493e5f4774SJed Brown   ierr = MatMPIBAIJSetPreallocation(A,bs,0,dnnz,0,onnz);CHKERRQ(ierr);
2503e5f4774SJed Brown   ierr = MatSeqSBAIJSetPreallocation(A,bs,0,dnnzu);CHKERRQ(ierr);
2513e5f4774SJed Brown   ierr = MatMPISBAIJSetPreallocation(A,bs,0,dnnzu,0,onnzu);CHKERRQ(ierr);
25263562e91SJed Brown   /*
25363562e91SJed Brown     In general, we have to do extra work to preallocate for scalar (AIJ) matrices so we check whether it will do any
25463562e91SJed Brown     good before going on with it.
25563562e91SJed Brown   */
25663562e91SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
25763562e91SJed Brown   if (!aij) {
25863562e91SJed Brown     ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
25963562e91SJed Brown   }
26063562e91SJed Brown   if (aij) {
26163562e91SJed Brown     if (bs == 1) {
2623e5f4774SJed Brown       ierr = MatSeqAIJSetPreallocation(A,0,dnnz);CHKERRQ(ierr);
2633e5f4774SJed Brown       ierr = MatMPIAIJSetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr);
2643e5f4774SJed Brown     } else {                    /* Convert block-row precallocation to scalar-row */
26563562e91SJed Brown       PetscInt i,m,*sdnnz,*sonnz;
26663562e91SJed Brown       ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr);
267dec54756SJed Brown       ierr = PetscMalloc2((!!dnnz)*m,PetscInt,&sdnnz,(!!onnz)*m,PetscInt,&sonnz);CHKERRQ(ierr);
268dec54756SJed Brown       for (i=0; i<m; i++) {
26963562e91SJed Brown         if (dnnz) sdnnz[i] = dnnz[i/bs] * bs;
27063562e91SJed Brown         if (onnz) sonnz[i] = onnz[i/bs] * bs;
27163562e91SJed Brown       }
2723e5f4774SJed Brown       ierr = MatSeqAIJSetPreallocation(A,0,dnnz ? sdnnz : PETSC_NULL);CHKERRQ(ierr);
2733e5f4774SJed Brown       ierr = MatMPIAIJSetPreallocation(A,0,dnnz ? sdnnz : PETSC_NULL,0,onnz ? sonnz : PETSC_NULL);CHKERRQ(ierr);
27463562e91SJed Brown       ierr = PetscFree2(sdnnz,sonnz);CHKERRQ(ierr);
27563562e91SJed Brown     }
27663562e91SJed Brown   }
27763562e91SJed Brown   PetscFunctionReturn(0);
27863562e91SJed Brown }
27963562e91SJed Brown 
280273d9f13SBarry Smith /*
281eb6b5d47SBarry Smith         Merges some information from Cs header to A; the C object is then destroyed
282d0f46423SBarry Smith 
283d0f46423SBarry Smith         This is somewhat different from MatHeaderReplace() it would be nice to merge the code
284273d9f13SBarry Smith */
2854a2ae208SSatish Balay #undef __FUNCT__
286eb6b5d47SBarry Smith #define __FUNCT__ "MatHeaderMerge"
287eb6b5d47SBarry Smith PetscErrorCode MatHeaderMerge(Mat A,Mat C)
288273d9f13SBarry Smith {
289dfbe8321SBarry Smith   PetscErrorCode ierr;
290d44834fbSBarry Smith   PetscInt       refct;
291273d9f13SBarry Smith   PetscOps       *Abops;
292273d9f13SBarry Smith   MatOps         Aops;
293273d9f13SBarry Smith   char           *mtype,*mname;
29430735b05SKris Buschelman   void           *spptr;
295273d9f13SBarry Smith 
296273d9f13SBarry Smith   PetscFunctionBegin;
297273d9f13SBarry Smith   /* save the parts of A we need */
2987adad957SLisandro Dalcin   Abops = ((PetscObject)A)->bops;
299273d9f13SBarry Smith   Aops  = A->ops;
3007adad957SLisandro Dalcin   refct = ((PetscObject)A)->refct;
3015c9eb25fSBarry Smith   mtype = ((PetscObject)A)->type_name;
3025c9eb25fSBarry Smith   mname = ((PetscObject)A)->name;
30330735b05SKris Buschelman   spptr = A->spptr;
30430735b05SKris Buschelman 
3055c9eb25fSBarry Smith   /* zero these so the destroy below does not free them */
3065c9eb25fSBarry Smith   ((PetscObject)A)->type_name = 0;
3075c9eb25fSBarry Smith   ((PetscObject)A)->name      = 0;
3085c9eb25fSBarry Smith 
3097c99f97cSSatish Balay   /* free all the interior data structures from mat */
3107c99f97cSSatish Balay   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
3117c99f97cSSatish Balay 
31230735b05SKris Buschelman   ierr = PetscFree(C->spptr);CHKERRQ(ierr);
31305b42c5fSBarry Smith 
3146bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->rmap);CHKERRQ(ierr);
3156bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->cmap);CHKERRQ(ierr);
316140e18c1SBarry Smith   ierr = PetscFunctionListDestroy(&((PetscObject)A)->qlist);CHKERRQ(ierr);
317140e18c1SBarry Smith   ierr = PetscObjectListDestroy(&((PetscObject)A)->olist);CHKERRQ(ierr);
318273d9f13SBarry Smith 
319273d9f13SBarry Smith   /* copy C over to A */
320273d9f13SBarry Smith   ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
321273d9f13SBarry Smith 
322273d9f13SBarry Smith   /* return the parts of A we saved */
3237adad957SLisandro Dalcin   ((PetscObject)A)->bops      = Abops;
324273d9f13SBarry Smith   A->ops                      = Aops;
3257adad957SLisandro Dalcin   ((PetscObject)A)->refct     = refct;
3267adad957SLisandro Dalcin   ((PetscObject)A)->type_name = mtype;
3277adad957SLisandro Dalcin   ((PetscObject)A)->name      = mname;
32830735b05SKris Buschelman   A->spptr                    = spptr;
329273d9f13SBarry Smith 
3305c9eb25fSBarry Smith   /* since these two are copied into A we do not want them destroyed in C */
3315c9eb25fSBarry Smith   ((PetscObject)C)->qlist = 0;
3325c9eb25fSBarry Smith   ((PetscObject)C)->olist = 0;
333*26fbe8dcSKarl Rupp 
3346bf464f9SBarry Smith   ierr = PetscHeaderDestroy(&C);CHKERRQ(ierr);
335273d9f13SBarry Smith   PetscFunctionReturn(0);
336273d9f13SBarry Smith }
3378ab5b326SKris Buschelman /*
338eb6b5d47SBarry Smith         Replace A's header with that of C; the C object is then destroyed
339d0f46423SBarry Smith 
340eb6b5d47SBarry Smith         This is essentially code moved from MatDestroy()
341eb6b5d47SBarry Smith 
342eb6b5d47SBarry Smith         This is somewhat different from MatHeaderMerge() it would be nice to merge the code
3438ab5b326SKris Buschelman */
3448ab5b326SKris Buschelman #undef __FUNCT__
3458ab5b326SKris Buschelman #define __FUNCT__ "MatHeaderReplace"
3468ab5b326SKris Buschelman PetscErrorCode MatHeaderReplace(Mat A,Mat C)
3478ab5b326SKris Buschelman {
3488ab5b326SKris Buschelman   PetscErrorCode ierr;
34927b31e29SJed Brown   PetscInt refct;
3508ab5b326SKris Buschelman 
3518ab5b326SKris Buschelman   PetscFunctionBegin;
35227b31e29SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
35327b31e29SJed Brown   PetscValidHeaderSpecific(C,MAT_CLASSID,2);
3546d7c1e57SBarry Smith   if (A == C) PetscFunctionReturn(0);
355a7e80f87SJed Brown   PetscCheckSameComm(A,1,C,2);
356a7e80f87SJed 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);
3576d7c1e57SBarry Smith 
3588ab5b326SKris Buschelman   /* free all the interior data structures from mat */
3598ab5b326SKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
360d2c95936SBarry Smith   ierr = PetscHeaderDestroy_Private((PetscObject)A);CHKERRQ(ierr);
361c650bb5fSLisandro Dalcin   ierr = PetscFree(A->ops);CHKERRQ(ierr);
3626bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->rmap);CHKERRQ(ierr);
3636bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->cmap);CHKERRQ(ierr);
364d95db058SSatish Balay   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
3658ab5b326SKris Buschelman 
3668ab5b326SKris Buschelman   /* copy C over to A */
36727b31e29SJed Brown   refct = ((PetscObject)A)->refct;
3688ab5b326SKris Buschelman   ierr  = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
369*26fbe8dcSKarl Rupp 
37027b31e29SJed Brown   ((PetscObject)A)->refct = refct;
371*26fbe8dcSKarl Rupp 
3728ab5b326SKris Buschelman   ierr = PetscLogObjectDestroy((PetscObject)C);CHKERRQ(ierr);
3738ab5b326SKris Buschelman   ierr = PetscFree(C);CHKERRQ(ierr);
3748ab5b326SKris Buschelman   PetscFunctionReturn(0);
3758ab5b326SKris Buschelman }
376