xref: /petsc/src/mat/utils/gcreate.c (revision 73107ff1babc248481eafc1deee82c9703b7a0a0)
17807a1faSBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/matimpl.h>       /*I "petscmat.h"  I*/
37807a1faSBarry Smith 
44a2ae208SSatish Balay #undef __FUNCT__
57d68702bSBarry Smith #define __FUNCT__ "MatShift_Basic"
67d68702bSBarry Smith PetscErrorCode MatShift_Basic(Mat Y,PetscScalar a)
77d68702bSBarry Smith {
87d68702bSBarry Smith   PetscErrorCode ierr;
97d68702bSBarry Smith   PetscInt       i,start,end;
107d68702bSBarry Smith   PetscScalar    alpha = a;
117d68702bSBarry Smith   PetscBool      prevoption;
127d68702bSBarry Smith 
137d68702bSBarry Smith   PetscFunctionBegin;
147d68702bSBarry Smith   ierr = MatGetOption(Y,MAT_NO_OFF_PROC_ENTRIES,&prevoption);CHKERRQ(ierr);
157d68702bSBarry Smith   ierr = MatSetOption(Y,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
167d68702bSBarry Smith   ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr);
177d68702bSBarry Smith   for (i=start; i<end; i++) {
187d68702bSBarry Smith     ierr = MatSetValues(Y,1,&i,1,&i,&alpha,ADD_VALUES);CHKERRQ(ierr);
197d68702bSBarry Smith   }
207d68702bSBarry Smith   ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
217d68702bSBarry Smith   ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
227d68702bSBarry Smith   ierr = MatSetOption(Y,MAT_NO_OFF_PROC_ENTRIES,prevoption);CHKERRQ(ierr);
237d68702bSBarry Smith   PetscFunctionReturn(0);
247d68702bSBarry Smith }
257d68702bSBarry Smith 
267d68702bSBarry Smith #undef __FUNCT__
274a2ae208SSatish Balay #define __FUNCT__ "MatCreate"
2805869f15SSatish Balay /*@
2969dd0797SLois Curfman McInnes    MatCreate - Creates a matrix where the type is determined
307e5f4302SBarry Smith    from either a call to MatSetType() or from the options database
317e5f4302SBarry Smith    with a call to MatSetFromOptions(). The default matrix type is
3269b1f4b7SBarry Smith    AIJ, using the routines MatCreateSeqAIJ() or MatCreateAIJ()
337e5f4302SBarry Smith    if you do not set a type in the options database. If you never
347e5f4302SBarry Smith    call MatSetType() or MatSetFromOptions() it will generate an
35f8ab6608SSatish Balay    error when you try to use the matrix.
3683e1b59cSLois Curfman McInnes 
37cb13003dSBarry Smith    Collective on MPI_Comm
38cb13003dSBarry Smith 
39f69a0ea3SMatthew Knepley    Input Parameter:
40f69a0ea3SMatthew Knepley .  comm - MPI communicator
417807a1faSBarry Smith 
427807a1faSBarry Smith    Output Parameter:
43dc401e71SLois Curfman McInnes .  A - the matrix
44e0b365e2SLois Curfman McInnes 
45273d9f13SBarry Smith    Options Database Keys:
46273d9f13SBarry Smith +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
4769b1f4b7SBarry Smith .    -mat_type mpiaij   - AIJ type, uses MatCreateAIJ()
48273d9f13SBarry Smith .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
4969b1f4b7SBarry Smith .    -mat_type mpidense - dense type, uses MatCreateDense()
50273d9f13SBarry Smith .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
5169b1f4b7SBarry Smith -    -mat_type mpibaij  - block AIJ type, uses MatCreateBAIJ()
52e0b365e2SLois Curfman McInnes 
5383e1b59cSLois Curfman McInnes    Even More Options Database Keys:
5483e1b59cSLois Curfman McInnes    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
5583e1b59cSLois Curfman McInnes    for additional format-specific options.
56e0b365e2SLois Curfman McInnes 
57bd9ce289SLois Curfman McInnes    Notes:
58ec6e0d80SSatish Balay 
59273d9f13SBarry Smith    Level: beginner
60273d9f13SBarry Smith 
61a2fc510eSBarry Smith    User manual sections:
62a2fc510eSBarry Smith +   sec_matcreate
63a2fc510eSBarry Smith -   chapter_matrices
64a2fc510eSBarry Smith 
65273d9f13SBarry Smith .keywords: matrix, create
66273d9f13SBarry Smith 
6769b1f4b7SBarry Smith .seealso: MatCreateSeqAIJ(), MatCreateAIJ(),
6869b1f4b7SBarry Smith           MatCreateSeqDense(), MatCreateDense(),
6969b1f4b7SBarry Smith           MatCreateSeqBAIJ(), MatCreateBAIJ(),
7069b1f4b7SBarry Smith           MatCreateSeqSBAIJ(), MatCreateSBAIJ(),
71273d9f13SBarry Smith           MatConvert()
72273d9f13SBarry Smith @*/
737087cfbeSBarry Smith PetscErrorCode  MatCreate(MPI_Comm comm,Mat *A)
74273d9f13SBarry Smith {
75273d9f13SBarry Smith   Mat            B;
76dfbe8321SBarry Smith   PetscErrorCode ierr;
77273d9f13SBarry Smith 
78273d9f13SBarry Smith   PetscFunctionBegin;
79f69a0ea3SMatthew Knepley   PetscValidPointer(A,2);
8097f1f81fSBarry Smith 
810298fd71SBarry Smith   *A = NULL;
82607a6623SBarry Smith   ierr = MatInitializePackage();CHKERRQ(ierr);
838ba1e511SMatthew Knepley 
84*73107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(B,MAT_CLASSID,"Mat","Matrix","Mat",comm,MatDestroy,MatView);CHKERRQ(ierr);
8526283091SBarry Smith   ierr = PetscLayoutCreate(comm,&B->rmap);CHKERRQ(ierr);
8626283091SBarry Smith   ierr = PetscLayoutCreate(comm,&B->cmap);CHKERRQ(ierr);
8726fbe8dcSKarl Rupp 
88273d9f13SBarry Smith   B->preallocated = PETSC_FALSE;
89273d9f13SBarry Smith   *A              = B;
90273d9f13SBarry Smith   PetscFunctionReturn(0);
91273d9f13SBarry Smith }
92273d9f13SBarry Smith 
934a2ae208SSatish Balay #undef __FUNCT__
94f69a0ea3SMatthew Knepley #define __FUNCT__ "MatSetSizes"
95f69a0ea3SMatthew Knepley /*@
96f69a0ea3SMatthew Knepley   MatSetSizes - Sets the local and global sizes, and checks to determine compatibility
97f69a0ea3SMatthew Knepley 
98f69a0ea3SMatthew Knepley   Collective on Mat
99f69a0ea3SMatthew Knepley 
100f69a0ea3SMatthew Knepley   Input Parameters:
101f69a0ea3SMatthew Knepley +  A - the matrix
102f69a0ea3SMatthew Knepley .  m - number of local rows (or PETSC_DECIDE)
103f69a0ea3SMatthew Knepley .  n - number of local columns (or PETSC_DECIDE)
104f69a0ea3SMatthew Knepley .  M - number of global rows (or PETSC_DETERMINE)
105f69a0ea3SMatthew Knepley -  N - number of global columns (or PETSC_DETERMINE)
106f69a0ea3SMatthew Knepley 
107f69a0ea3SMatthew Knepley    Notes:
108f69a0ea3SMatthew Knepley    m (n) and M (N) cannot be both PETSC_DECIDE
109f69a0ea3SMatthew Knepley    If one processor calls this with M (N) of PETSC_DECIDE then all processors must, otherwise the program will hang.
110f69a0ea3SMatthew Knepley 
111f69a0ea3SMatthew Knepley    If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the
112f69a0ea3SMatthew Knepley    user must ensure that they are chosen to be compatible with the
113f69a0ea3SMatthew Knepley    vectors. To do this, one first considers the matrix-vector product
114f69a0ea3SMatthew Knepley    'y = A x'. The 'm' that is used in the above routine must match the
115f69a0ea3SMatthew Knepley    local size used in the vector creation routine VecCreateMPI() for 'y'.
116f69a0ea3SMatthew Knepley    Likewise, the 'n' used must match that used as the local size in
117f69a0ea3SMatthew Knepley    VecCreateMPI() for 'x'.
118f69a0ea3SMatthew Knepley 
119f73d5cc4SBarry Smith    You cannot change the sizes once they have been set.
120f73d5cc4SBarry Smith 
121f73d5cc4SBarry Smith    The sizes must be set before MatSetUp() or MatXXXSetPreallocation() is called.
122f73d5cc4SBarry Smith 
123f69a0ea3SMatthew Knepley   Level: beginner
124f69a0ea3SMatthew Knepley 
125f69a0ea3SMatthew Knepley .seealso: MatGetSize(), PetscSplitOwnership()
126f69a0ea3SMatthew Knepley @*/
1277087cfbeSBarry Smith PetscErrorCode  MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N)
128f69a0ea3SMatthew Knepley {
129f69a0ea3SMatthew Knepley   PetscFunctionBegin;
1300700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
131b28e6c9fSJed Brown   if (M > 0) PetscValidLogicalCollectiveInt(A,M,4);
132b28e6c9fSJed Brown   if (N > 0) PetscValidLogicalCollectiveInt(A,N,5);
133e32f2f54SBarry 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);
134e32f2f54SBarry 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);
135461878b2SBarry 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);
136461878b2SBarry 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);
137d0f46423SBarry Smith   A->rmap->n = m;
138d0f46423SBarry Smith   A->cmap->n = n;
139d0f46423SBarry Smith   A->rmap->N = M;
140d0f46423SBarry Smith   A->cmap->N = N;
141f69a0ea3SMatthew Knepley   PetscFunctionReturn(0);
142f69a0ea3SMatthew Knepley }
143f69a0ea3SMatthew Knepley 
144f69a0ea3SMatthew Knepley #undef __FUNCT__
1454a2ae208SSatish Balay #define __FUNCT__ "MatSetFromOptions"
14605869f15SSatish Balay /*@
147273d9f13SBarry Smith    MatSetFromOptions - Creates a matrix where the type is determined
148273d9f13SBarry Smith    from the options database. Generates a parallel MPI matrix if the
149273d9f13SBarry Smith    communicator has more than one processor.  The default matrix type is
15069b1f4b7SBarry Smith    AIJ, using the routines MatCreateSeqAIJ() and MatCreateAIJ() if
1517e5f4302SBarry Smith    you do not select a type in the options database.
152273d9f13SBarry Smith 
153273d9f13SBarry Smith    Collective on Mat
154273d9f13SBarry Smith 
155273d9f13SBarry Smith    Input Parameter:
156273d9f13SBarry Smith .  A - the matrix
157273d9f13SBarry Smith 
158273d9f13SBarry Smith    Options Database Keys:
159273d9f13SBarry Smith +    -mat_type seqaij   - AIJ type, uses MatCreateSeqAIJ()
16069b1f4b7SBarry Smith .    -mat_type mpiaij   - AIJ type, uses MatCreateAIJ()
161273d9f13SBarry Smith .    -mat_type seqdense - dense type, uses MatCreateSeqDense()
16269b1f4b7SBarry Smith .    -mat_type mpidense - dense type, uses MatCreateDense()
163273d9f13SBarry Smith .    -mat_type seqbaij  - block AIJ type, uses MatCreateSeqBAIJ()
16469b1f4b7SBarry Smith -    -mat_type mpibaij  - block AIJ type, uses MatCreateBAIJ()
165273d9f13SBarry Smith 
166273d9f13SBarry Smith    Even More Options Database Keys:
167273d9f13SBarry Smith    See the manpages for particular formats (e.g., MatCreateSeqAIJ())
168273d9f13SBarry Smith    for additional format-specific options.
169bd9ce289SLois Curfman McInnes 
1701d69843bSLois Curfman McInnes    Level: beginner
1711d69843bSLois Curfman McInnes 
172dc401e71SLois Curfman McInnes .keywords: matrix, create
173e0b365e2SLois Curfman McInnes 
17469b1f4b7SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateAIJ(),
17569b1f4b7SBarry Smith           MatCreateSeqDense(), MatCreateDense(),
17669b1f4b7SBarry Smith           MatCreateSeqBAIJ(), MatCreateBAIJ(),
17769b1f4b7SBarry Smith           MatCreateSeqSBAIJ(), MatCreateSBAIJ(),
178273d9f13SBarry Smith           MatConvert()
1797807a1faSBarry Smith @*/
1807087cfbeSBarry Smith PetscErrorCode  MatSetFromOptions(Mat B)
1817807a1faSBarry Smith {
182dfbe8321SBarry Smith   PetscErrorCode ierr;
183f3be49caSLisandro Dalcin   const char     *deft = MATAIJ;
184f3be49caSLisandro Dalcin   char           type[256];
18569df5c0cSJed Brown   PetscBool      flg,set;
186dbb450caSBarry Smith 
1873a40ed3dSBarry Smith   PetscFunctionBegin;
1880700a824SBarry Smith   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
189f3be49caSLisandro Dalcin 
1903194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)B);CHKERRQ(ierr);
191535b19f3SBarry Smith 
192535b19f3SBarry Smith   if (B->rmap->bs < 0) {
193535b19f3SBarry Smith     PetscInt newbs = -1;
194535b19f3SBarry Smith     ierr = PetscOptionsInt("-mat_block_size","Set the blocksize used to store the matrix","MatSetBlockSize",newbs,&newbs,&flg);CHKERRQ(ierr);
195535b19f3SBarry Smith     if (flg) {
196535b19f3SBarry Smith       ierr = PetscLayoutSetBlockSize(B->rmap,newbs);CHKERRQ(ierr);
197535b19f3SBarry Smith       ierr = PetscLayoutSetBlockSize(B->cmap,newbs);CHKERRQ(ierr);
198535b19f3SBarry Smith     }
199535b19f3SBarry Smith   }
200535b19f3SBarry Smith 
201a264d7a6SBarry Smith   ierr = PetscOptionsFList("-mat_type","Matrix type","MatSetType",MatList,deft,type,256,&flg);CHKERRQ(ierr);
202273d9f13SBarry Smith   if (flg) {
203f3be49caSLisandro Dalcin     ierr = MatSetType(B,type);CHKERRQ(ierr);
204f3be49caSLisandro Dalcin   } else if (!((PetscObject)B)->type_name) {
205f3be49caSLisandro Dalcin     ierr = MatSetType(B,deft);CHKERRQ(ierr);
206273d9f13SBarry Smith   }
207f3be49caSLisandro Dalcin 
208840d65ccSBarry Smith   ierr = PetscOptionsName("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",&B->checksymmetryonassembly);CHKERRQ(ierr);
20994ae4db5SBarry Smith   ierr = PetscOptionsReal("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",B->checksymmetrytol,&B->checksymmetrytol,NULL);CHKERRQ(ierr);
21094ae4db5SBarry Smith   ierr = PetscOptionsBool("-mat_null_space_test","Checks if provided null space is correct in MatAssemblyEnd()","MatSetNullSpaceTest",B->checknullspaceonassembly,&B->checknullspaceonassembly,NULL);CHKERRQ(ierr);
211840d65ccSBarry Smith 
212f3be49caSLisandro Dalcin   if (B->ops->setfromoptions) {
213e55864a3SBarry Smith     ierr = (*B->ops->setfromoptions)(PetscOptionsObject,B);CHKERRQ(ierr);
214593c1905SSatish Balay   }
215f3be49caSLisandro Dalcin 
21669df5c0cSJed Brown   flg  = PETSC_FALSE;
21769df5c0cSJed 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);
21869df5c0cSJed Brown   if (set) {ierr = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,flg);CHKERRQ(ierr);}
21969df5c0cSJed Brown   flg  = PETSC_FALSE;
22069df5c0cSJed 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);
22169df5c0cSJed Brown   if (set) {ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,flg);CHKERRQ(ierr);}
22269df5c0cSJed Brown 
2235d973c19SBarry Smith   /* process any options handlers added with PetscObjectAddOptionsHandler() */
2245d973c19SBarry Smith   ierr = PetscObjectProcessOptionsHandlers((PetscObject)B);CHKERRQ(ierr);
225f3be49caSLisandro Dalcin   ierr = PetscOptionsEnd();CHKERRQ(ierr);
2263a40ed3dSBarry Smith   PetscFunctionReturn(0);
2277807a1faSBarry Smith }
2287807a1faSBarry Smith 
2294a2ae208SSatish Balay #undef __FUNCT__
23063562e91SJed Brown #define __FUNCT__ "MatXAIJSetPreallocation"
23163562e91SJed Brown /*@
23263562e91SJed Brown    MatXAIJSetPreallocation - set preallocation for serial and parallel AIJ, BAIJ, and SBAIJ matrices
23363562e91SJed Brown 
23463562e91SJed Brown    Collective on Mat
23563562e91SJed Brown 
23663562e91SJed Brown    Input Arguments:
23763562e91SJed Brown +  A - matrix being preallocated
23863562e91SJed Brown .  bs - block size
2393e5f4774SJed Brown .  dnnz - number of nonzero blocks per block row of diagonal part of parallel matrix
2403e5f4774SJed Brown .  onnz - number of nonzero blocks per block row of off-diagonal part of parallel matrix
2413e5f4774SJed Brown .  dnnzu - number of nonzero blocks per block row of upper-triangular part of diagonal part of parallel matrix
2423e5f4774SJed Brown -  onnzu - number of nonzero blocks per block row of upper-triangular part of off-diagonal part of parallel matrix
24363562e91SJed Brown 
24463562e91SJed Brown    Level: beginner
24563562e91SJed Brown 
246ab978733SBarry Smith .seealso: MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation(),
247ab978733SBarry Smith           PetscSplitOwnership()
24863562e91SJed Brown @*/
2494cce697cSJed Brown PetscErrorCode MatXAIJSetPreallocation(Mat A,PetscInt bs,const PetscInt dnnz[],const PetscInt onnz[],const PetscInt dnnzu[],const PetscInt onnzu[])
25063562e91SJed Brown {
25163562e91SJed Brown   PetscErrorCode ierr;
25263562e91SJed Brown   void           (*aij)(void);
25363562e91SJed Brown 
25463562e91SJed Brown   PetscFunctionBegin;
255dec54756SJed Brown   ierr = MatSetBlockSize(A,bs);CHKERRQ(ierr);
25669bbac97SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
257dec54756SJed Brown   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
258dec54756SJed Brown   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
2593e5f4774SJed Brown   ierr = MatSeqBAIJSetPreallocation(A,bs,0,dnnz);CHKERRQ(ierr);
2603e5f4774SJed Brown   ierr = MatMPIBAIJSetPreallocation(A,bs,0,dnnz,0,onnz);CHKERRQ(ierr);
2613e5f4774SJed Brown   ierr = MatSeqSBAIJSetPreallocation(A,bs,0,dnnzu);CHKERRQ(ierr);
2623e5f4774SJed Brown   ierr = MatMPISBAIJSetPreallocation(A,bs,0,dnnzu,0,onnzu);CHKERRQ(ierr);
26363562e91SJed Brown   /*
26463562e91SJed Brown     In general, we have to do extra work to preallocate for scalar (AIJ) matrices so we check whether it will do any
26563562e91SJed Brown     good before going on with it.
26663562e91SJed Brown   */
26763562e91SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
26863562e91SJed Brown   if (!aij) {
26963562e91SJed Brown     ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
27063562e91SJed Brown   }
27163562e91SJed Brown   if (aij) {
27263562e91SJed Brown     if (bs == 1) {
2733e5f4774SJed Brown       ierr = MatSeqAIJSetPreallocation(A,0,dnnz);CHKERRQ(ierr);
2743e5f4774SJed Brown       ierr = MatMPIAIJSetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr);
2753e5f4774SJed Brown     } else {                    /* Convert block-row precallocation to scalar-row */
27663562e91SJed Brown       PetscInt i,m,*sdnnz,*sonnz;
2770298fd71SBarry Smith       ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
278dcca6d9dSJed Brown       ierr = PetscMalloc2((!!dnnz)*m,&sdnnz,(!!onnz)*m,&sonnz);CHKERRQ(ierr);
279dec54756SJed Brown       for (i=0; i<m; i++) {
28063562e91SJed Brown         if (dnnz) sdnnz[i] = dnnz[i/bs] * bs;
28163562e91SJed Brown         if (onnz) sonnz[i] = onnz[i/bs] * bs;
28263562e91SJed Brown       }
2830298fd71SBarry Smith       ierr = MatSeqAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL);CHKERRQ(ierr);
2840298fd71SBarry Smith       ierr = MatMPIAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr);
28563562e91SJed Brown       ierr = PetscFree2(sdnnz,sonnz);CHKERRQ(ierr);
28663562e91SJed Brown     }
28763562e91SJed Brown   }
28863562e91SJed Brown   PetscFunctionReturn(0);
28963562e91SJed Brown }
29063562e91SJed Brown 
291273d9f13SBarry Smith /*
292eb6b5d47SBarry Smith         Merges some information from Cs header to A; the C object is then destroyed
293d0f46423SBarry Smith 
294d0f46423SBarry Smith         This is somewhat different from MatHeaderReplace() it would be nice to merge the code
295273d9f13SBarry Smith */
2964a2ae208SSatish Balay #undef __FUNCT__
297eb6b5d47SBarry Smith #define __FUNCT__ "MatHeaderMerge"
298eb6b5d47SBarry Smith PetscErrorCode MatHeaderMerge(Mat A,Mat C)
299273d9f13SBarry Smith {
300dfbe8321SBarry Smith   PetscErrorCode ierr;
301d44834fbSBarry Smith   PetscInt       refct;
302*73107ff1SLisandro Dalcin   PetscOps       Abops;
303*73107ff1SLisandro Dalcin   struct _MatOps Aops;
304273d9f13SBarry Smith   char           *mtype,*mname;
30530735b05SKris Buschelman   void           *spptr;
306273d9f13SBarry Smith 
307273d9f13SBarry Smith   PetscFunctionBegin;
308273d9f13SBarry Smith   /* save the parts of A we need */
309*73107ff1SLisandro Dalcin   Abops = ((PetscObject)A)->bops[0];
310*73107ff1SLisandro Dalcin   Aops  = A->ops[0];
3117adad957SLisandro Dalcin   refct = ((PetscObject)A)->refct;
3125c9eb25fSBarry Smith   mtype = ((PetscObject)A)->type_name;
3135c9eb25fSBarry Smith   mname = ((PetscObject)A)->name;
31430735b05SKris Buschelman   spptr = A->spptr;
31530735b05SKris Buschelman 
3165c9eb25fSBarry Smith   /* zero these so the destroy below does not free them */
3175c9eb25fSBarry Smith   ((PetscObject)A)->type_name = 0;
3185c9eb25fSBarry Smith   ((PetscObject)A)->name      = 0;
3195c9eb25fSBarry Smith 
3207c99f97cSSatish Balay   /* free all the interior data structures from mat */
3217c99f97cSSatish Balay   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
3227c99f97cSSatish Balay 
32330735b05SKris Buschelman   ierr = PetscFree(C->spptr);CHKERRQ(ierr);
32405b42c5fSBarry Smith 
3256bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->rmap);CHKERRQ(ierr);
3266bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->cmap);CHKERRQ(ierr);
327140e18c1SBarry Smith   ierr = PetscFunctionListDestroy(&((PetscObject)A)->qlist);CHKERRQ(ierr);
328140e18c1SBarry Smith   ierr = PetscObjectListDestroy(&((PetscObject)A)->olist);CHKERRQ(ierr);
329273d9f13SBarry Smith 
330273d9f13SBarry Smith   /* copy C over to A */
331273d9f13SBarry Smith   ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
332273d9f13SBarry Smith 
333273d9f13SBarry Smith   /* return the parts of A we saved */
334*73107ff1SLisandro Dalcin   ((PetscObject)A)->bops[0]   = Abops;
335*73107ff1SLisandro Dalcin   A->ops[0]                   = Aops;
3367adad957SLisandro Dalcin   ((PetscObject)A)->refct     = refct;
3377adad957SLisandro Dalcin   ((PetscObject)A)->type_name = mtype;
3387adad957SLisandro Dalcin   ((PetscObject)A)->name      = mname;
33930735b05SKris Buschelman   A->spptr                    = spptr;
340273d9f13SBarry Smith 
3415c9eb25fSBarry Smith   /* since these two are copied into A we do not want them destroyed in C */
3425c9eb25fSBarry Smith   ((PetscObject)C)->qlist = 0;
3435c9eb25fSBarry Smith   ((PetscObject)C)->olist = 0;
34426fbe8dcSKarl Rupp 
3456bf464f9SBarry Smith   ierr = PetscHeaderDestroy(&C);CHKERRQ(ierr);
346273d9f13SBarry Smith   PetscFunctionReturn(0);
347273d9f13SBarry Smith }
3488ab5b326SKris Buschelman /*
349eb6b5d47SBarry Smith         Replace A's header with that of C; the C object is then destroyed
350d0f46423SBarry Smith 
351eb6b5d47SBarry Smith         This is essentially code moved from MatDestroy()
352eb6b5d47SBarry Smith 
353eb6b5d47SBarry Smith         This is somewhat different from MatHeaderMerge() it would be nice to merge the code
354b30237c6SBarry Smith 
355b30237c6SBarry Smith         Used in DM hence is declared PETSC_EXTERN
3568ab5b326SKris Buschelman */
3578ab5b326SKris Buschelman #undef __FUNCT__
3588ab5b326SKris Buschelman #define __FUNCT__ "MatHeaderReplace"
359b30237c6SBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat A,Mat C)
3608ab5b326SKris Buschelman {
3618ab5b326SKris Buschelman   PetscErrorCode   ierr;
36227b31e29SJed Brown   PetscInt         refct;
363fefd9316SJose E. Roman   PetscObjectState state;
3648ab5b326SKris Buschelman 
3658ab5b326SKris Buschelman   PetscFunctionBegin;
36627b31e29SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
36727b31e29SJed Brown   PetscValidHeaderSpecific(C,MAT_CLASSID,2);
3686d7c1e57SBarry Smith   if (A == C) PetscFunctionReturn(0);
369a7e80f87SJed Brown   PetscCheckSameComm(A,1,C,2);
370ce94432eSBarry 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);
3716d7c1e57SBarry Smith 
3728ab5b326SKris Buschelman   /* free all the interior data structures from mat */
3738ab5b326SKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
374d2c95936SBarry Smith   ierr = PetscHeaderDestroy_Private((PetscObject)A);CHKERRQ(ierr);
3756bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->rmap);CHKERRQ(ierr);
3766bf464f9SBarry Smith   ierr = PetscLayoutDestroy(&A->cmap);CHKERRQ(ierr);
377d95db058SSatish Balay   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
3788ab5b326SKris Buschelman 
3798ab5b326SKris Buschelman   /* copy C over to A */
38027b31e29SJed Brown   refct = ((PetscObject)A)->refct;
381fefd9316SJose E. Roman   state = ((PetscObject)A)->state;
3828ab5b326SKris Buschelman   ierr  = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
38326fbe8dcSKarl Rupp 
38427b31e29SJed Brown   ((PetscObject)A)->refct = refct;
385fefd9316SJose E. Roman   ((PetscObject)A)->state = state + 1;
38626fbe8dcSKarl Rupp 
3878ab5b326SKris Buschelman   ierr = PetscFree(C);CHKERRQ(ierr);
3888ab5b326SKris Buschelman   PetscFunctionReturn(0);
3898ab5b326SKris Buschelman }
390