xref: /petsc/src/mat/impls/cdiagonal/cdiagonal.c (revision 4d43f095504bf18eecf7d1688d0cb4f6e5ace5bf)
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 
83423f386SBarry Smith 
93423f386SBarry Smith /* ----------------------------------------------------------------------------------------*/
103423f386SBarry Smith static PetscErrorCode MatDestroy_ConstantDiagonal(Mat mat)
113423f386SBarry Smith {
123423f386SBarry Smith   PetscErrorCode       ierr;
133423f386SBarry Smith 
143423f386SBarry Smith   PetscFunctionBegin;
153423f386SBarry Smith   ierr = PetscFree(mat->data);CHKERRQ(ierr);
163423f386SBarry Smith   PetscFunctionReturn(0);
173423f386SBarry Smith }
183423f386SBarry Smith 
193423f386SBarry Smith static PetscErrorCode MatView_ConstantDiagonal(Mat J,PetscViewer viewer)
203423f386SBarry Smith {
213423f386SBarry Smith   PetscErrorCode       ierr;
223423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)J->data;
233423f386SBarry Smith   PetscBool            iascii;
243423f386SBarry Smith 
253423f386SBarry Smith   PetscFunctionBegin;
263423f386SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
273423f386SBarry Smith   if (iascii) {
283423f386SBarry Smith     PetscViewerFormat    format;
293423f386SBarry Smith 
303423f386SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
313423f386SBarry Smith     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) PetscFunctionReturn(0);
323423f386SBarry Smith #if !defined(PETSC_USE_COMPLEX)
333423f386SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Diagonal value: %g\n",ctx->diag);CHKERRQ(ierr);
343423f386SBarry Smith #else
353423f386SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Diagonal value: %g + i %g\n",PetscRealPart(ctx->diag),PetscImaginaryPart(ctx->diag));CHKERRQ(ierr);
363423f386SBarry Smith #endif
373423f386SBarry Smith   }
383423f386SBarry Smith   PetscFunctionReturn(0);
393423f386SBarry Smith }
403423f386SBarry Smith 
413423f386SBarry Smith static PetscErrorCode MatAssemblyEnd_ConstantDiagonal(Mat J,MatAssemblyType mt)
423423f386SBarry Smith {
433423f386SBarry Smith   PetscFunctionBegin;
443423f386SBarry Smith   PetscFunctionReturn(0);
453423f386SBarry Smith }
463423f386SBarry Smith 
473423f386SBarry Smith static PetscErrorCode MatMult_ConstantDiagonal(Mat J,Vec x,Vec y)
483423f386SBarry Smith {
493423f386SBarry Smith   PetscErrorCode       ierr;
503423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)J->data;
513423f386SBarry Smith 
523423f386SBarry Smith   PetscFunctionBegin;
533954eb43SBarry Smith   ierr = VecAXPBY(y,ctx->diag,0.0,x);CHKERRQ(ierr);
543423f386SBarry Smith   PetscFunctionReturn(0);
553423f386SBarry Smith }
563423f386SBarry Smith 
573423f386SBarry Smith PetscErrorCode MatGetDiagonal_ConstantDiagonal(Mat J,Vec x)
583423f386SBarry Smith {
593423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)J->data;
603423f386SBarry Smith   PetscErrorCode       ierr;
613423f386SBarry Smith 
623423f386SBarry Smith   PetscFunctionBegin;
633423f386SBarry Smith   ierr = VecSet(x,ctx->diag);CHKERRQ(ierr);
643423f386SBarry Smith   PetscFunctionReturn(0);
653423f386SBarry Smith }
663423f386SBarry Smith 
673423f386SBarry Smith static PetscErrorCode MatShift_ConstantDiagonal(Mat Y,PetscScalar a)
683423f386SBarry Smith {
693954eb43SBarry Smith   PetscErrorCode       ierr;
703423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)Y->data;
713423f386SBarry Smith 
723423f386SBarry Smith   PetscFunctionBegin;
733954eb43SBarry Smith   if (a != 0.) {ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr);}
743423f386SBarry Smith   ctx->diag += a;
753423f386SBarry Smith   PetscFunctionReturn(0);
763423f386SBarry Smith }
773423f386SBarry Smith 
783423f386SBarry Smith static PetscErrorCode MatScale_ConstantDiagonal(Mat Y,PetscScalar a)
793423f386SBarry Smith {
80*4d43f095SSatish Balay   PetscErrorCode       ierr;
813423f386SBarry Smith   Mat_ConstantDiagonal *ctx  = (Mat_ConstantDiagonal*)Y->data;
823423f386SBarry Smith 
833423f386SBarry Smith   PetscFunctionBegin;
84*4d43f095SSatish Balay   if (a != 1.) {ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr);}
853423f386SBarry Smith   ctx->diag *= a;
863423f386SBarry Smith   PetscFunctionReturn(0);
873423f386SBarry Smith }
883423f386SBarry Smith 
893423f386SBarry Smith static PetscErrorCode MatZeroEntries_ConstantDiagonal(Mat Y)
903423f386SBarry Smith {
91*4d43f095SSatish Balay   PetscErrorCode       ierr;
923423f386SBarry Smith   Mat_ConstantDiagonal *ctx  = (Mat_ConstantDiagonal*)Y->data;
933423f386SBarry Smith 
943423f386SBarry Smith   PetscFunctionBegin;
95*4d43f095SSatish Balay   if (ctx->diag != 0.0) {ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr);}
963423f386SBarry Smith   ctx->diag = 0.0;
973423f386SBarry Smith   PetscFunctionReturn(0);
983423f386SBarry Smith }
993423f386SBarry Smith 
1003423f386SBarry Smith PetscErrorCode MatSOR_ConstantDiagonal(Mat matin,Vec x,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec y)
1013423f386SBarry Smith {
1023423f386SBarry Smith   PetscErrorCode       ierr;
1033423f386SBarry Smith   Mat_ConstantDiagonal *ctx  = (Mat_ConstantDiagonal*)matin->data;
1043423f386SBarry Smith 
1053423f386SBarry Smith   PetscFunctionBegin;
1063423f386SBarry Smith   if (ctx->diag == 0.0) matin->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1073423f386SBarry Smith   else matin->factorerrortype = MAT_FACTOR_NOERROR;
1083954eb43SBarry Smith   ierr = VecAXPBY(y,1.0/ctx->diag,0.0,x);CHKERRQ(ierr);
1093423f386SBarry Smith   PetscFunctionReturn(0);
1103423f386SBarry Smith }
1113423f386SBarry Smith 
1123423f386SBarry Smith PetscErrorCode MatGetInfo_ConstantDiagonal(Mat A,MatInfoType flag,MatInfo *info)
1133423f386SBarry Smith {
1143423f386SBarry Smith   PetscFunctionBegin;
1153423f386SBarry Smith   info->block_size   = 1.0;
1163423f386SBarry Smith   info->nz_allocated = 1.0;
1173423f386SBarry Smith   info->nz_used      = 1.0;
1183423f386SBarry Smith   info->nz_unneeded  = 0.0;
1193423f386SBarry Smith   info->assemblies   = (double)A->num_ass;
1203423f386SBarry Smith   info->mallocs      = 0.0;
1213423f386SBarry Smith   info->memory       = ((PetscObject)A)->mem;
1223423f386SBarry Smith   if (A->factortype) {
1233423f386SBarry Smith     info->fill_ratio_given  = 1.0;
1243423f386SBarry Smith     info->fill_ratio_needed = 1.0;
1253423f386SBarry Smith     info->factor_mallocs    = 0.0;
1263423f386SBarry Smith   } else {
1273423f386SBarry Smith     info->fill_ratio_given  = 0;
1283423f386SBarry Smith     info->fill_ratio_needed = 0;
1293423f386SBarry Smith     info->factor_mallocs    = 0;
1303423f386SBarry Smith   }
1313423f386SBarry Smith   PetscFunctionReturn(0);
1323423f386SBarry Smith }
1333423f386SBarry Smith 
1343423f386SBarry Smith /*@
1353423f386SBarry Smith    MatCreateConstantDiagonal - Creates a matrix with a uniform value along the diagonal
1363423f386SBarry Smith 
1373423f386SBarry Smith    Collective on MPI_Comm
1383423f386SBarry Smith 
1393423f386SBarry Smith    Input Parameters:
1403423f386SBarry Smith +  comm - MPI communicator
1413423f386SBarry Smith .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1423423f386SBarry Smith            This value should be the same as the local size used in creating the
1433423f386SBarry Smith            y vector for the matrix-vector product y = Ax.
1443423f386SBarry Smith .  n - This value should be the same as the local size used in creating the
1453423f386SBarry Smith        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
1463423f386SBarry Smith        calculated if N is given) For square matrices n is almost always m.
1473423f386SBarry Smith .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
1483423f386SBarry Smith .  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
1493423f386SBarry Smith -  diag - the diagonal value
1503423f386SBarry Smith 
1513423f386SBarry Smith    Output Parameter:
1523423f386SBarry Smith .  J - the diagonal matrix
1533423f386SBarry Smith 
1543423f386SBarry Smith    Level: advanced
1553423f386SBarry Smith 
1563423f386SBarry Smith    Notes:
1573423f386SBarry Smith     Only supports square matrices with the same number of local rows and columns
1583423f386SBarry Smith 
1593423f386SBarry Smith .seealso: MatDestroy(), MATCONSTANTDIAGONAL, MatScale(), MatShift(), MatMult(), MatGetDiagonal(), MatGetFactor(), MatSolve()
1603423f386SBarry Smith 
1613423f386SBarry Smith @*/
1623423f386SBarry Smith PetscErrorCode  MatCreateConstantDiagonal(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar diag,Mat *J)
1633423f386SBarry Smith {
1643423f386SBarry Smith   PetscErrorCode ierr;
1653423f386SBarry Smith 
1663423f386SBarry Smith   PetscFunctionBegin;
1673423f386SBarry Smith   ierr = MatCreate(comm,J);CHKERRQ(ierr);
1683423f386SBarry Smith   ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr);
1693423f386SBarry Smith   ierr = MatSetType(*J,MATCONSTANTDIAGONAL);CHKERRQ(ierr);
1703423f386SBarry Smith   ierr = MatShift(*J,diag);CHKERRQ(ierr);
1713423f386SBarry Smith   ierr = MatSetUp(*J);CHKERRQ(ierr);
1723423f386SBarry Smith   PetscFunctionReturn(0);
1733423f386SBarry Smith }
1743423f386SBarry Smith 
1753423f386SBarry Smith PETSC_EXTERN PetscErrorCode  MatCreate_ConstantDiagonal(Mat A)
1763423f386SBarry Smith {
1773423f386SBarry Smith   PetscErrorCode       ierr;
1783423f386SBarry Smith   Mat_ConstantDiagonal *ctx;
1793423f386SBarry Smith 
1803423f386SBarry Smith   PetscFunctionBegin;
1813423f386SBarry Smith   ierr = PetscNew(&ctx);CHKERRQ(ierr);
1823423f386SBarry Smith   ctx->diag = 0.0;
1833423f386SBarry Smith   A->data   = (void*)ctx;
1843423f386SBarry Smith 
1853423f386SBarry Smith   A->assembled        = PETSC_TRUE;
1863423f386SBarry Smith   A->preallocated     = PETSC_TRUE;
1873423f386SBarry Smith   A->ops->mult        = MatMult_ConstantDiagonal;
1883423f386SBarry Smith   A->ops->sor         = MatSOR_ConstantDiagonal;
1893423f386SBarry Smith   A->ops->shift       = MatShift_ConstantDiagonal;
1903423f386SBarry Smith   A->ops->scale       = MatScale_ConstantDiagonal;
1913423f386SBarry Smith   A->ops->getdiagonal = MatGetDiagonal_ConstantDiagonal;
1923423f386SBarry Smith   A->ops->view        = MatView_ConstantDiagonal;
1933423f386SBarry Smith   A->ops->zeroentries = MatZeroEntries_ConstantDiagonal;
1943423f386SBarry Smith   A->ops->assemblyend = MatAssemblyEnd_ConstantDiagonal;
1953423f386SBarry Smith   A->ops->destroy     = MatDestroy_ConstantDiagonal;
1963423f386SBarry Smith   A->ops->getinfo     = MatGetInfo_ConstantDiagonal;
1973423f386SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATCONSTANTDIAGONAL);CHKERRQ(ierr);
1983423f386SBarry Smith   PetscFunctionReturn(0);
1993423f386SBarry Smith }
2003423f386SBarry Smith 
2013423f386SBarry Smith static PetscErrorCode MatFactorNumeric_ConstantDiagonal(Mat fact,Mat A,const MatFactorInfo *info)
2023423f386SBarry Smith {
203*4d43f095SSatish Balay   PetscErrorCode       ierr;
2043423f386SBarry Smith   Mat_ConstantDiagonal *actx = (Mat_ConstantDiagonal*)A->data,*fctx = (Mat_ConstantDiagonal*)fact->data;
2053423f386SBarry Smith 
2063423f386SBarry Smith   PetscFunctionBegin;
2073423f386SBarry Smith   if (actx->diag == 0.0) fact->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
2083423f386SBarry Smith   else fact->factorerrortype = MAT_FACTOR_NOERROR;
2093423f386SBarry Smith   fctx->diag = 1.0/actx->diag;
210*4d43f095SSatish Balay   ierr = PetscObjectStateIncrease((PetscObject)fact);CHKERRQ(ierr);
2113423f386SBarry Smith   fact->ops->solve = MatMult_ConstantDiagonal;
2123423f386SBarry Smith   PetscFunctionReturn(0);
2133423f386SBarry Smith }
2143423f386SBarry Smith 
2153423f386SBarry Smith static PetscErrorCode MatFactorSymbolic_LU_ConstantDiagonal(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
2163423f386SBarry Smith {
2173423f386SBarry Smith   PetscFunctionBegin;
2183423f386SBarry Smith   fact->ops->lufactornumeric = MatFactorNumeric_ConstantDiagonal;
2193423f386SBarry Smith   PetscFunctionReturn(0);
2203423f386SBarry Smith }
2213423f386SBarry Smith 
2223423f386SBarry Smith static PetscErrorCode MatFactorSymbolic_Cholesky_ConstantDiagonal(Mat fact,Mat A,IS isrow,const MatFactorInfo *info)
2233423f386SBarry Smith {
2243423f386SBarry Smith   PetscFunctionBegin;
2253423f386SBarry Smith   fact->ops->choleskyfactornumeric = MatFactorNumeric_ConstantDiagonal;  PetscFunctionReturn(0);
2263423f386SBarry Smith   PetscFunctionReturn(0);
2273423f386SBarry Smith }
2283423f386SBarry Smith 
2293423f386SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_constantdiagonal_petsc(Mat A,MatFactorType ftype,Mat *B)
2303423f386SBarry Smith {
2313423f386SBarry Smith   PetscInt       n = A->rmap->n, N = A->rmap->N;
2323423f386SBarry Smith   PetscErrorCode ierr;
2333423f386SBarry Smith 
2343423f386SBarry Smith   PetscFunctionBegin;
2353423f386SBarry Smith   ierr = MatCreateConstantDiagonal(PetscObjectComm((PetscObject)A),n,n,N,N,0,B);CHKERRQ(ierr);
2363423f386SBarry Smith 
2373423f386SBarry Smith   (*B)->factortype = ftype;
2383423f386SBarry Smith   (*B)->ops->ilufactorsymbolic      = MatFactorSymbolic_LU_ConstantDiagonal;
2393423f386SBarry Smith   (*B)->ops->lufactorsymbolic       = MatFactorSymbolic_LU_ConstantDiagonal;
2403423f386SBarry Smith   (*B)->ops->iccfactorsymbolic      = MatFactorSymbolic_Cholesky_ConstantDiagonal;
2413423f386SBarry Smith   (*B)->ops->choleskyfactorsymbolic = MatFactorSymbolic_Cholesky_ConstantDiagonal;
2423423f386SBarry Smith 
2433423f386SBarry Smith   (*B)->ops->shift       = NULL;
2443423f386SBarry Smith   (*B)->ops->scale       = NULL;
2453423f386SBarry Smith   (*B)->ops->mult        = NULL;
2463423f386SBarry Smith   (*B)->ops->sor         = NULL;
2473423f386SBarry Smith   (*B)->ops->zeroentries = NULL;
2483423f386SBarry Smith 
2493423f386SBarry Smith   ierr = PetscFree((*B)->solvertype);CHKERRQ(ierr);
2503423f386SBarry Smith   ierr = PetscStrallocpy(MATSOLVERPETSC,&(*B)->solvertype);CHKERRQ(ierr);
2513423f386SBarry Smith   PetscFunctionReturn(0);
2523423f386SBarry Smith }
253