xref: /petsc/src/mat/impls/cdiagonal/cdiagonal.c (revision e5f3c4be5c9b3f690a675cc36a53ce609624f7fa)
13423f386SBarry Smith 
23423f386SBarry Smith #include <petsc/private/matimpl.h>  /*I "petscmat.h" I*/
33423f386SBarry Smith 
43423f386SBarry Smith typedef struct {
53423f386SBarry Smith   PetscScalar diag;
63423f386SBarry Smith } Mat_ConstantDiagonal;
73423f386SBarry Smith 
8a01b8767SStefano Zampini static PetscErrorCode MatAXPY_ConstantDiagonal(Mat Y, PetscScalar a, Mat X, MatStructure str)
9a01b8767SStefano Zampini {
10a01b8767SStefano Zampini   Mat_ConstantDiagonal *yctx = (Mat_ConstantDiagonal*)Y->data;
11a01b8767SStefano Zampini   Mat_ConstantDiagonal *xctx = (Mat_ConstantDiagonal*)X->data;
12a01b8767SStefano Zampini 
13a01b8767SStefano Zampini   PetscFunctionBegin;
14a01b8767SStefano Zampini   yctx->diag += a*xctx->diag;
15a01b8767SStefano Zampini   PetscFunctionReturn(0);
16a01b8767SStefano Zampini }
17a01b8767SStefano Zampini 
18a8d4b458SStefano Zampini static PetscErrorCode MatGetRow_ConstantDiagonal(Mat A, PetscInt row, PetscInt *ncols, PetscInt *cols[], PetscScalar *vals[])
19a8d4b458SStefano Zampini {
20a8d4b458SStefano Zampini   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)A->data;
21a8d4b458SStefano Zampini   PetscErrorCode       ierr;
223423f386SBarry Smith 
23a8d4b458SStefano Zampini   PetscFunctionBegin;
24a8d4b458SStefano Zampini   if (ncols) *ncols = 1;
25a8d4b458SStefano Zampini   if (cols) {
26a8d4b458SStefano Zampini     ierr = PetscMalloc1(1,cols);CHKERRQ(ierr);
27a8d4b458SStefano Zampini     (*cols)[0] = row;
28a8d4b458SStefano Zampini   }
29a8d4b458SStefano Zampini   if (vals) {
30a8d4b458SStefano Zampini     ierr = PetscMalloc1(1,vals);CHKERRQ(ierr);
31a8d4b458SStefano Zampini     (*vals)[0] = ctx->diag;
32a8d4b458SStefano Zampini   }
33a8d4b458SStefano Zampini   PetscFunctionReturn(0);
34a8d4b458SStefano Zampini }
35a8d4b458SStefano Zampini 
36a8d4b458SStefano Zampini static PetscErrorCode MatRestoreRow_ConstantDiagonal(Mat A, PetscInt row, PetscInt *ncols, PetscInt *cols[], PetscScalar *vals[])
37a8d4b458SStefano Zampini {
38a8d4b458SStefano Zampini   PetscErrorCode ierr;
39a8d4b458SStefano Zampini 
40a8d4b458SStefano Zampini   PetscFunctionBegin;
41a8d4b458SStefano Zampini   if (ncols) *ncols = 0;
42a8d4b458SStefano Zampini   if (cols) {
43a8d4b458SStefano Zampini     ierr = PetscFree(*cols);CHKERRQ(ierr);
44a8d4b458SStefano Zampini   }
45a8d4b458SStefano Zampini   if (vals) {
46a8d4b458SStefano Zampini     ierr = PetscFree(*vals);CHKERRQ(ierr);
47a8d4b458SStefano Zampini   }
48a8d4b458SStefano Zampini   PetscFunctionReturn(0);
49a8d4b458SStefano Zampini }
50a8d4b458SStefano Zampini 
51a8d4b458SStefano Zampini static PetscErrorCode MatMultTranspose_ConstantDiagonal(Mat A, Vec x, Vec y)
52a8d4b458SStefano Zampini {
53a8d4b458SStefano Zampini   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)A->data;
54a8d4b458SStefano Zampini   PetscErrorCode       ierr;
55a8d4b458SStefano Zampini 
56a8d4b458SStefano Zampini   PetscFunctionBegin;
57a8d4b458SStefano Zampini   ierr = VecAXPBY(y,ctx->diag,0.0,x);CHKERRQ(ierr);
58a8d4b458SStefano Zampini   PetscFunctionReturn(0);
59a8d4b458SStefano Zampini }
60a8d4b458SStefano Zampini 
61a8d4b458SStefano Zampini static PetscErrorCode MatMultAdd_ConstantDiagonal(Mat mat,Vec v1,Vec v2,Vec v3)
62a8d4b458SStefano Zampini {
63a8d4b458SStefano Zampini   PetscErrorCode       ierr;
64a8d4b458SStefano Zampini   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)mat->data;
65a8d4b458SStefano Zampini 
66a8d4b458SStefano Zampini   PetscFunctionBegin;
67a8d4b458SStefano Zampini   if (v2 == v3) {
68a8d4b458SStefano Zampini     ierr = VecAXPBY(v3,ctx->diag,1.0,v1);CHKERRQ(ierr);
69a8d4b458SStefano Zampini   } else {
70a8d4b458SStefano Zampini     ierr = VecAXPBYPCZ(v3,ctx->diag,1.0,0.0,v1,v2);CHKERRQ(ierr);
71a8d4b458SStefano Zampini   }
72a8d4b458SStefano Zampini   PetscFunctionReturn(0);
73a8d4b458SStefano Zampini }
74a8d4b458SStefano Zampini 
75a8d4b458SStefano Zampini static PetscErrorCode MatMultTransposeAdd_ConstantDiagonal(Mat mat,Vec v1,Vec v2,Vec v3)
76a8d4b458SStefano Zampini {
77a8d4b458SStefano Zampini   PetscErrorCode       ierr;
78a8d4b458SStefano Zampini   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)mat->data;
79a8d4b458SStefano Zampini 
80a8d4b458SStefano Zampini   PetscFunctionBegin;
81a8d4b458SStefano Zampini   if (v2 == v3) {
82a8d4b458SStefano Zampini     ierr = VecAXPBY(v3,ctx->diag,1.0,v1);CHKERRQ(ierr);
83a8d4b458SStefano Zampini   } else {
84a8d4b458SStefano Zampini     ierr = VecAXPBYPCZ(v3,ctx->diag,1.0,0.0,v1,v2);CHKERRQ(ierr);
85a8d4b458SStefano Zampini   }
86a8d4b458SStefano Zampini   PetscFunctionReturn(0);
87a8d4b458SStefano Zampini }
88a8d4b458SStefano Zampini 
89*e5f3c4beSJose E. Roman static PetscErrorCode MatNorm_ConstantDiagonal(Mat A,NormType type,PetscReal *nrm)
90*e5f3c4beSJose E. Roman {
91*e5f3c4beSJose E. Roman   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)A->data;
92*e5f3c4beSJose E. Roman 
93*e5f3c4beSJose E. Roman   PetscFunctionBegin;
94*e5f3c4beSJose E. Roman   if (type == NORM_FROBENIUS || type == NORM_2 || type == NORM_1 || type == NORM_INFINITY) *nrm = PetscAbsScalar(ctx->diag);
95*e5f3c4beSJose E. Roman   else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Unsupported norm");
96*e5f3c4beSJose E. Roman   PetscFunctionReturn(0);
97*e5f3c4beSJose E. Roman }
98*e5f3c4beSJose E. Roman 
99a8d4b458SStefano Zampini static PetscErrorCode MatDuplicate_ConstantDiagonal(Mat A, MatDuplicateOption op, Mat *B)
100a8d4b458SStefano Zampini {
101a8d4b458SStefano Zampini   PetscErrorCode       ierr;
102a8d4b458SStefano Zampini   Mat_ConstantDiagonal *actx = (Mat_ConstantDiagonal*)A->data;
103a8d4b458SStefano Zampini 
104a8d4b458SStefano Zampini   PetscFunctionBegin;
105a8d4b458SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)A),B);CHKERRQ(ierr);
106a8d4b458SStefano Zampini   ierr = MatSetSizes(*B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
107a8d4b458SStefano Zampini   ierr = MatSetBlockSizesFromMats(*B,A,A);CHKERRQ(ierr);
108a8d4b458SStefano Zampini   ierr = MatSetType(*B,MATCONSTANTDIAGONAL);CHKERRQ(ierr);
109a8d4b458SStefano Zampini   ierr = PetscLayoutReference(A->rmap,&(*B)->rmap);CHKERRQ(ierr);
110a8d4b458SStefano Zampini   ierr = PetscLayoutReference(A->cmap,&(*B)->cmap);CHKERRQ(ierr);
111a8d4b458SStefano Zampini   if (op == MAT_COPY_VALUES) {
112a8d4b458SStefano Zampini     Mat_ConstantDiagonal *bctx = (Mat_ConstantDiagonal*)(*B)->data;
113a8d4b458SStefano Zampini     bctx->diag = actx->diag;
114a8d4b458SStefano Zampini   }
115a8d4b458SStefano Zampini   PetscFunctionReturn(0);
116a8d4b458SStefano Zampini }
117a8d4b458SStefano Zampini 
118a8d4b458SStefano Zampini static PetscErrorCode MatMissingDiagonal_ConstantDiagonal(Mat mat,PetscBool *missing,PetscInt *dd)
119a8d4b458SStefano Zampini {
120a8d4b458SStefano Zampini   PetscFunctionBegin;
121a8d4b458SStefano Zampini   *missing = PETSC_FALSE;
122a8d4b458SStefano Zampini   PetscFunctionReturn(0);
123a8d4b458SStefano Zampini }
124a8d4b458SStefano Zampini 
1253423f386SBarry Smith static PetscErrorCode MatDestroy_ConstantDiagonal(Mat mat)
1263423f386SBarry Smith {
1273423f386SBarry Smith   PetscErrorCode       ierr;
1283423f386SBarry Smith 
1293423f386SBarry Smith   PetscFunctionBegin;
1303423f386SBarry Smith   ierr = PetscFree(mat->data);CHKERRQ(ierr);
1313423f386SBarry Smith   PetscFunctionReturn(0);
1323423f386SBarry Smith }
1333423f386SBarry Smith 
1343423f386SBarry Smith static PetscErrorCode MatView_ConstantDiagonal(Mat J,PetscViewer viewer)
1353423f386SBarry Smith {
1363423f386SBarry Smith   PetscErrorCode       ierr;
1373423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)J->data;
1383423f386SBarry Smith   PetscBool            iascii;
1393423f386SBarry Smith 
1403423f386SBarry Smith   PetscFunctionBegin;
1413423f386SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1423423f386SBarry Smith   if (iascii) {
1433423f386SBarry Smith     PetscViewerFormat    format;
1443423f386SBarry Smith 
1453423f386SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1463423f386SBarry Smith     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) PetscFunctionReturn(0);
1473423f386SBarry Smith #if !defined(PETSC_USE_COMPLEX)
1483423f386SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Diagonal value: %g\n",ctx->diag);CHKERRQ(ierr);
1493423f386SBarry Smith #else
1503423f386SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Diagonal value: %g + i %g\n",PetscRealPart(ctx->diag),PetscImaginaryPart(ctx->diag));CHKERRQ(ierr);
1513423f386SBarry Smith #endif
1523423f386SBarry Smith   }
1533423f386SBarry Smith   PetscFunctionReturn(0);
1543423f386SBarry Smith }
1553423f386SBarry Smith 
1563423f386SBarry Smith static PetscErrorCode MatAssemblyEnd_ConstantDiagonal(Mat J,MatAssemblyType mt)
1573423f386SBarry Smith {
1583423f386SBarry Smith   PetscFunctionBegin;
1593423f386SBarry Smith   PetscFunctionReturn(0);
1603423f386SBarry Smith }
1613423f386SBarry Smith 
1623423f386SBarry Smith static PetscErrorCode MatMult_ConstantDiagonal(Mat J,Vec x,Vec y)
1633423f386SBarry Smith {
1643423f386SBarry Smith   PetscErrorCode       ierr;
1653423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)J->data;
1663423f386SBarry Smith 
1673423f386SBarry Smith   PetscFunctionBegin;
1683954eb43SBarry Smith   ierr = VecAXPBY(y,ctx->diag,0.0,x);CHKERRQ(ierr);
1693423f386SBarry Smith   PetscFunctionReturn(0);
1703423f386SBarry Smith }
1713423f386SBarry Smith 
1723423f386SBarry Smith PetscErrorCode MatGetDiagonal_ConstantDiagonal(Mat J,Vec x)
1733423f386SBarry Smith {
1743423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)J->data;
1753423f386SBarry Smith   PetscErrorCode       ierr;
1763423f386SBarry Smith 
1773423f386SBarry Smith   PetscFunctionBegin;
1783423f386SBarry Smith   ierr = VecSet(x,ctx->diag);CHKERRQ(ierr);
1793423f386SBarry Smith   PetscFunctionReturn(0);
1803423f386SBarry Smith }
1813423f386SBarry Smith 
1823423f386SBarry Smith static PetscErrorCode MatShift_ConstantDiagonal(Mat Y,PetscScalar a)
1833423f386SBarry Smith {
1843423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)Y->data;
1853423f386SBarry Smith 
1863423f386SBarry Smith   PetscFunctionBegin;
1873423f386SBarry Smith   ctx->diag += a;
1883423f386SBarry Smith   PetscFunctionReturn(0);
1893423f386SBarry Smith }
1903423f386SBarry Smith 
1913423f386SBarry Smith static PetscErrorCode MatScale_ConstantDiagonal(Mat Y,PetscScalar a)
1923423f386SBarry Smith {
1933423f386SBarry Smith   Mat_ConstantDiagonal *ctx  = (Mat_ConstantDiagonal*)Y->data;
1943423f386SBarry Smith 
1953423f386SBarry Smith   PetscFunctionBegin;
1963423f386SBarry Smith   ctx->diag *= a;
1973423f386SBarry Smith   PetscFunctionReturn(0);
1983423f386SBarry Smith }
1993423f386SBarry Smith 
2003423f386SBarry Smith static PetscErrorCode MatZeroEntries_ConstantDiagonal(Mat Y)
2013423f386SBarry Smith {
2023423f386SBarry Smith   Mat_ConstantDiagonal *ctx  = (Mat_ConstantDiagonal*)Y->data;
2033423f386SBarry Smith 
2043423f386SBarry Smith   PetscFunctionBegin;
2053423f386SBarry Smith   ctx->diag = 0.0;
2063423f386SBarry Smith   PetscFunctionReturn(0);
2073423f386SBarry Smith }
2083423f386SBarry Smith 
2093423f386SBarry Smith PetscErrorCode MatSOR_ConstantDiagonal(Mat matin,Vec x,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec y)
2103423f386SBarry Smith {
2113423f386SBarry Smith   PetscErrorCode       ierr;
2123423f386SBarry Smith   Mat_ConstantDiagonal *ctx  = (Mat_ConstantDiagonal*)matin->data;
2133423f386SBarry Smith 
2143423f386SBarry Smith   PetscFunctionBegin;
2153423f386SBarry Smith   if (ctx->diag == 0.0) matin->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
2163423f386SBarry Smith   else matin->factorerrortype = MAT_FACTOR_NOERROR;
2173954eb43SBarry Smith   ierr = VecAXPBY(y,1.0/ctx->diag,0.0,x);CHKERRQ(ierr);
2183423f386SBarry Smith   PetscFunctionReturn(0);
2193423f386SBarry Smith }
2203423f386SBarry Smith 
2213423f386SBarry Smith PetscErrorCode MatGetInfo_ConstantDiagonal(Mat A,MatInfoType flag,MatInfo *info)
2223423f386SBarry Smith {
2233423f386SBarry Smith   PetscFunctionBegin;
2243423f386SBarry Smith   info->block_size   = 1.0;
2253423f386SBarry Smith   info->nz_allocated = 1.0;
2263423f386SBarry Smith   info->nz_used      = 1.0;
2273423f386SBarry Smith   info->nz_unneeded  = 0.0;
2283966268fSBarry Smith   info->assemblies   = A->num_ass;
2293423f386SBarry Smith   info->mallocs      = 0.0;
2303423f386SBarry Smith   info->memory       = ((PetscObject)A)->mem;
2313423f386SBarry Smith   if (A->factortype) {
2323423f386SBarry Smith     info->fill_ratio_given  = 1.0;
2333423f386SBarry Smith     info->fill_ratio_needed = 1.0;
2343423f386SBarry Smith     info->factor_mallocs    = 0.0;
2353423f386SBarry Smith   } else {
2363423f386SBarry Smith     info->fill_ratio_given  = 0;
2373423f386SBarry Smith     info->fill_ratio_needed = 0;
2383423f386SBarry Smith     info->factor_mallocs    = 0;
2393423f386SBarry Smith   }
2403423f386SBarry Smith   PetscFunctionReturn(0);
2413423f386SBarry Smith }
2423423f386SBarry Smith 
2433423f386SBarry Smith /*@
2443423f386SBarry Smith    MatCreateConstantDiagonal - Creates a matrix with a uniform value along the diagonal
2453423f386SBarry Smith 
246d083f849SBarry Smith    Collective
2473423f386SBarry Smith 
2483423f386SBarry Smith    Input Parameters:
2493423f386SBarry Smith +  comm - MPI communicator
2503423f386SBarry Smith .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
2513423f386SBarry Smith            This value should be the same as the local size used in creating the
2523423f386SBarry Smith            y vector for the matrix-vector product y = Ax.
2533423f386SBarry Smith .  n - This value should be the same as the local size used in creating the
2543423f386SBarry Smith        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
2553423f386SBarry Smith        calculated if N is given) For square matrices n is almost always m.
2563423f386SBarry Smith .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
2573423f386SBarry Smith .  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
2583423f386SBarry Smith -  diag - the diagonal value
2593423f386SBarry Smith 
2603423f386SBarry Smith    Output Parameter:
2613423f386SBarry Smith .  J - the diagonal matrix
2623423f386SBarry Smith 
2633423f386SBarry Smith    Level: advanced
2643423f386SBarry Smith 
2653423f386SBarry Smith    Notes:
2663423f386SBarry Smith     Only supports square matrices with the same number of local rows and columns
2673423f386SBarry Smith 
2683423f386SBarry Smith .seealso: MatDestroy(), MATCONSTANTDIAGONAL, MatScale(), MatShift(), MatMult(), MatGetDiagonal(), MatGetFactor(), MatSolve()
2693423f386SBarry Smith 
2703423f386SBarry Smith @*/
2713423f386SBarry Smith PetscErrorCode  MatCreateConstantDiagonal(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar diag,Mat *J)
2723423f386SBarry Smith {
2733423f386SBarry Smith   PetscErrorCode ierr;
2743423f386SBarry Smith 
2753423f386SBarry Smith   PetscFunctionBegin;
2763423f386SBarry Smith   ierr = MatCreate(comm,J);CHKERRQ(ierr);
2773423f386SBarry Smith   ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr);
2783423f386SBarry Smith   ierr = MatSetType(*J,MATCONSTANTDIAGONAL);CHKERRQ(ierr);
2793423f386SBarry Smith   ierr = MatShift(*J,diag);CHKERRQ(ierr);
2803423f386SBarry Smith   ierr = MatSetUp(*J);CHKERRQ(ierr);
2813423f386SBarry Smith   PetscFunctionReturn(0);
2823423f386SBarry Smith }
2833423f386SBarry Smith 
2843423f386SBarry Smith PETSC_EXTERN PetscErrorCode  MatCreate_ConstantDiagonal(Mat A)
2853423f386SBarry Smith {
2863423f386SBarry Smith   PetscErrorCode       ierr;
2873423f386SBarry Smith   Mat_ConstantDiagonal *ctx;
2883423f386SBarry Smith 
2893423f386SBarry Smith   PetscFunctionBegin;
2903423f386SBarry Smith   ierr = PetscNew(&ctx);CHKERRQ(ierr);
2913423f386SBarry Smith   ctx->diag = 0.0;
2923423f386SBarry Smith   A->data   = (void*)ctx;
2933423f386SBarry Smith 
2943423f386SBarry Smith   A->assembled    = PETSC_TRUE;
2953423f386SBarry Smith   A->preallocated = PETSC_TRUE;
296a8d4b458SStefano Zampini 
2973423f386SBarry Smith   A->ops->mult             = MatMult_ConstantDiagonal;
298a8d4b458SStefano Zampini   A->ops->multadd          = MatMultAdd_ConstantDiagonal;
299a8d4b458SStefano Zampini   A->ops->multtranspose    = MatMultTranspose_ConstantDiagonal;
300a8d4b458SStefano Zampini   A->ops->multtransposeadd = MatMultTransposeAdd_ConstantDiagonal;
301*e5f3c4beSJose E. Roman   A->ops->norm             = MatNorm_ConstantDiagonal;
302a8d4b458SStefano Zampini   A->ops->duplicate        = MatDuplicate_ConstantDiagonal;
303a8d4b458SStefano Zampini   A->ops->missingdiagonal  = MatMissingDiagonal_ConstantDiagonal;
304a8d4b458SStefano Zampini   A->ops->getrow           = MatGetRow_ConstantDiagonal;
305a8d4b458SStefano Zampini   A->ops->restorerow       = MatRestoreRow_ConstantDiagonal;
3063423f386SBarry Smith   A->ops->sor              = MatSOR_ConstantDiagonal;
3073423f386SBarry Smith   A->ops->shift            = MatShift_ConstantDiagonal;
3083423f386SBarry Smith   A->ops->scale            = MatScale_ConstantDiagonal;
3093423f386SBarry Smith   A->ops->getdiagonal      = MatGetDiagonal_ConstantDiagonal;
3103423f386SBarry Smith   A->ops->view             = MatView_ConstantDiagonal;
3113423f386SBarry Smith   A->ops->zeroentries      = MatZeroEntries_ConstantDiagonal;
3123423f386SBarry Smith   A->ops->assemblyend      = MatAssemblyEnd_ConstantDiagonal;
3133423f386SBarry Smith   A->ops->destroy          = MatDestroy_ConstantDiagonal;
3143423f386SBarry Smith   A->ops->getinfo          = MatGetInfo_ConstantDiagonal;
315a01b8767SStefano Zampini   A->ops->axpy             = MatAXPY_ConstantDiagonal;
316a8d4b458SStefano Zampini 
3173423f386SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATCONSTANTDIAGONAL);CHKERRQ(ierr);
3183423f386SBarry Smith   PetscFunctionReturn(0);
3193423f386SBarry Smith }
3203423f386SBarry Smith 
3213423f386SBarry Smith static PetscErrorCode MatFactorNumeric_ConstantDiagonal(Mat fact,Mat A,const MatFactorInfo *info)
3223423f386SBarry Smith {
3233423f386SBarry Smith   Mat_ConstantDiagonal *actx = (Mat_ConstantDiagonal*)A->data,*fctx = (Mat_ConstantDiagonal*)fact->data;
3243423f386SBarry Smith 
3253423f386SBarry Smith   PetscFunctionBegin;
3263423f386SBarry Smith   if (actx->diag == 0.0) fact->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3273423f386SBarry Smith   else fact->factorerrortype = MAT_FACTOR_NOERROR;
3283423f386SBarry Smith   fctx->diag = 1.0/actx->diag;
3293423f386SBarry Smith   fact->ops->solve = MatMult_ConstantDiagonal;
3303423f386SBarry Smith   PetscFunctionReturn(0);
3313423f386SBarry Smith }
3323423f386SBarry Smith 
3333423f386SBarry Smith static PetscErrorCode MatFactorSymbolic_LU_ConstantDiagonal(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
3343423f386SBarry Smith {
3353423f386SBarry Smith   PetscFunctionBegin;
3363423f386SBarry Smith   fact->ops->lufactornumeric = MatFactorNumeric_ConstantDiagonal;
3373423f386SBarry Smith   PetscFunctionReturn(0);
3383423f386SBarry Smith }
3393423f386SBarry Smith 
3403423f386SBarry Smith static PetscErrorCode MatFactorSymbolic_Cholesky_ConstantDiagonal(Mat fact,Mat A,IS isrow,const MatFactorInfo *info)
3413423f386SBarry Smith {
3423423f386SBarry Smith   PetscFunctionBegin;
343a01b8767SStefano Zampini   fact->ops->choleskyfactornumeric = MatFactorNumeric_ConstantDiagonal;
3443423f386SBarry Smith   PetscFunctionReturn(0);
3453423f386SBarry Smith }
3463423f386SBarry Smith 
3473423f386SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_constantdiagonal_petsc(Mat A,MatFactorType ftype,Mat *B)
3483423f386SBarry Smith {
3493423f386SBarry Smith   PetscInt       n = A->rmap->n, N = A->rmap->N;
3503423f386SBarry Smith   PetscErrorCode ierr;
3513423f386SBarry Smith 
3523423f386SBarry Smith   PetscFunctionBegin;
3533423f386SBarry Smith   ierr = MatCreateConstantDiagonal(PetscObjectComm((PetscObject)A),n,n,N,N,0,B);CHKERRQ(ierr);
3543423f386SBarry Smith 
3553423f386SBarry Smith   (*B)->factortype = ftype;
3563423f386SBarry Smith   (*B)->ops->ilufactorsymbolic      = MatFactorSymbolic_LU_ConstantDiagonal;
3573423f386SBarry Smith   (*B)->ops->lufactorsymbolic       = MatFactorSymbolic_LU_ConstantDiagonal;
3583423f386SBarry Smith   (*B)->ops->iccfactorsymbolic      = MatFactorSymbolic_Cholesky_ConstantDiagonal;
3593423f386SBarry Smith   (*B)->ops->choleskyfactorsymbolic = MatFactorSymbolic_Cholesky_ConstantDiagonal;
3603423f386SBarry Smith 
3613423f386SBarry Smith   (*B)->ops->shift       = NULL;
3623423f386SBarry Smith   (*B)->ops->scale       = NULL;
3633423f386SBarry Smith   (*B)->ops->mult        = NULL;
3643423f386SBarry Smith   (*B)->ops->sor         = NULL;
3653423f386SBarry Smith   (*B)->ops->zeroentries = NULL;
3663423f386SBarry Smith 
3673423f386SBarry Smith   ierr = PetscFree((*B)->solvertype);CHKERRQ(ierr);
3683423f386SBarry Smith   ierr = PetscStrallocpy(MATSOLVERPETSC,&(*B)->solvertype);CHKERRQ(ierr);
3693423f386SBarry Smith   PetscFunctionReturn(0);
3703423f386SBarry Smith }
371