xref: /petsc/src/mat/impls/cdiagonal/cdiagonal.c (revision 3423f3868bef8a3400519521cb1daa60d2113289)
1*3423f386SBarry Smith 
2*3423f386SBarry Smith #include <petsc/private/matimpl.h>  /*I "petscmat.h" I*/
3*3423f386SBarry Smith 
4*3423f386SBarry Smith typedef struct {
5*3423f386SBarry Smith   PetscScalar diag;
6*3423f386SBarry Smith } Mat_ConstantDiagonal;
7*3423f386SBarry Smith 
8*3423f386SBarry Smith 
9*3423f386SBarry Smith /* ----------------------------------------------------------------------------------------*/
10*3423f386SBarry Smith static PetscErrorCode MatDestroy_ConstantDiagonal(Mat mat)
11*3423f386SBarry Smith {
12*3423f386SBarry Smith   PetscErrorCode       ierr;
13*3423f386SBarry Smith 
14*3423f386SBarry Smith   PetscFunctionBegin;
15*3423f386SBarry Smith   ierr = PetscFree(mat->data);CHKERRQ(ierr);
16*3423f386SBarry Smith   PetscFunctionReturn(0);
17*3423f386SBarry Smith }
18*3423f386SBarry Smith 
19*3423f386SBarry Smith static PetscErrorCode MatView_ConstantDiagonal(Mat J,PetscViewer viewer)
20*3423f386SBarry Smith {
21*3423f386SBarry Smith   PetscErrorCode       ierr;
22*3423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)J->data;
23*3423f386SBarry Smith   PetscBool            iascii;
24*3423f386SBarry Smith 
25*3423f386SBarry Smith   PetscFunctionBegin;
26*3423f386SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
27*3423f386SBarry Smith   if (iascii) {
28*3423f386SBarry Smith     PetscViewerFormat    format;
29*3423f386SBarry Smith 
30*3423f386SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
31*3423f386SBarry Smith     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) PetscFunctionReturn(0);
32*3423f386SBarry Smith #if !defined(PETSC_USE_COMPLEX)
33*3423f386SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Diagonal value: %g\n",ctx->diag);CHKERRQ(ierr);
34*3423f386SBarry Smith #else
35*3423f386SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Diagonal value: %g + i %g\n",PetscRealPart(ctx->diag),PetscImaginaryPart(ctx->diag));CHKERRQ(ierr);
36*3423f386SBarry Smith #endif
37*3423f386SBarry Smith   }
38*3423f386SBarry Smith   PetscFunctionReturn(0);
39*3423f386SBarry Smith }
40*3423f386SBarry Smith 
41*3423f386SBarry Smith static PetscErrorCode MatAssemblyEnd_ConstantDiagonal(Mat J,MatAssemblyType mt)
42*3423f386SBarry Smith {
43*3423f386SBarry Smith   PetscFunctionBegin;
44*3423f386SBarry Smith   PetscFunctionReturn(0);
45*3423f386SBarry Smith }
46*3423f386SBarry Smith 
47*3423f386SBarry Smith static PetscErrorCode MatMult_ConstantDiagonal(Mat J,Vec x,Vec y)
48*3423f386SBarry Smith {
49*3423f386SBarry Smith   PetscErrorCode       ierr;
50*3423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)J->data;
51*3423f386SBarry Smith 
52*3423f386SBarry Smith   PetscFunctionBegin;
53*3423f386SBarry Smith   ierr = VecAXPBY(y,ctx->diag,0.0,x);
54*3423f386SBarry Smith   PetscFunctionReturn(0);
55*3423f386SBarry Smith }
56*3423f386SBarry Smith 
57*3423f386SBarry Smith PetscErrorCode MatGetDiagonal_ConstantDiagonal(Mat J,Vec x)
58*3423f386SBarry Smith {
59*3423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)J->data;
60*3423f386SBarry Smith   PetscErrorCode       ierr;
61*3423f386SBarry Smith 
62*3423f386SBarry Smith   PetscFunctionBegin;
63*3423f386SBarry Smith   ierr = VecSet(x,ctx->diag);CHKERRQ(ierr);
64*3423f386SBarry Smith   PetscFunctionReturn(0);
65*3423f386SBarry Smith }
66*3423f386SBarry Smith 
67*3423f386SBarry Smith static PetscErrorCode MatShift_ConstantDiagonal(Mat Y,PetscScalar a)
68*3423f386SBarry Smith {
69*3423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)Y->data;
70*3423f386SBarry Smith 
71*3423f386SBarry Smith   PetscFunctionBegin;
72*3423f386SBarry Smith   if (a != 0.) PetscObjectStateIncrease((PetscObject)Y);
73*3423f386SBarry Smith   ctx->diag += a;
74*3423f386SBarry Smith   PetscFunctionReturn(0);
75*3423f386SBarry Smith }
76*3423f386SBarry Smith 
77*3423f386SBarry Smith static PetscErrorCode MatScale_ConstantDiagonal(Mat Y,PetscScalar a)
78*3423f386SBarry Smith {
79*3423f386SBarry Smith   Mat_ConstantDiagonal *ctx  = (Mat_ConstantDiagonal*)Y->data;
80*3423f386SBarry Smith 
81*3423f386SBarry Smith   PetscFunctionBegin;
82*3423f386SBarry Smith   if (a != 1.) PetscObjectStateIncrease((PetscObject)Y);
83*3423f386SBarry Smith   ctx->diag *= a;
84*3423f386SBarry Smith   PetscFunctionReturn(0);
85*3423f386SBarry Smith }
86*3423f386SBarry Smith 
87*3423f386SBarry Smith static PetscErrorCode MatZeroEntries_ConstantDiagonal(Mat Y)
88*3423f386SBarry Smith {
89*3423f386SBarry Smith   Mat_ConstantDiagonal *ctx  = (Mat_ConstantDiagonal*)Y->data;
90*3423f386SBarry Smith 
91*3423f386SBarry Smith   PetscFunctionBegin;
92*3423f386SBarry Smith   if (ctx->diag != 0.0) PetscObjectStateIncrease((PetscObject)Y);
93*3423f386SBarry Smith   ctx->diag = 0.0;
94*3423f386SBarry Smith   PetscFunctionReturn(0);
95*3423f386SBarry Smith }
96*3423f386SBarry Smith 
97*3423f386SBarry Smith PetscErrorCode MatSOR_ConstantDiagonal(Mat matin,Vec x,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec y)
98*3423f386SBarry Smith {
99*3423f386SBarry Smith   PetscErrorCode       ierr;
100*3423f386SBarry Smith   Mat_ConstantDiagonal *ctx  = (Mat_ConstantDiagonal*)matin->data;
101*3423f386SBarry Smith 
102*3423f386SBarry Smith   PetscFunctionBegin;
103*3423f386SBarry Smith   if (ctx->diag == 0.0) matin->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
104*3423f386SBarry Smith   else matin->factorerrortype = MAT_FACTOR_NOERROR;
105*3423f386SBarry Smith   ierr = VecAXPBY(y,1.0/ctx->diag,0.0,x);
106*3423f386SBarry Smith   PetscFunctionReturn(0);
107*3423f386SBarry Smith }
108*3423f386SBarry Smith 
109*3423f386SBarry Smith PetscErrorCode MatGetInfo_ConstantDiagonal(Mat A,MatInfoType flag,MatInfo *info)
110*3423f386SBarry Smith {
111*3423f386SBarry Smith   PetscFunctionBegin;
112*3423f386SBarry Smith   info->block_size   = 1.0;
113*3423f386SBarry Smith   info->nz_allocated = 1.0;
114*3423f386SBarry Smith   info->nz_used      = 1.0;
115*3423f386SBarry Smith   info->nz_unneeded  = 0.0;
116*3423f386SBarry Smith   info->assemblies   = (double)A->num_ass;
117*3423f386SBarry Smith   info->mallocs      = 0.0;
118*3423f386SBarry Smith   info->memory       = ((PetscObject)A)->mem;
119*3423f386SBarry Smith   if (A->factortype) {
120*3423f386SBarry Smith     info->fill_ratio_given  = 1.0;
121*3423f386SBarry Smith     info->fill_ratio_needed = 1.0;
122*3423f386SBarry Smith     info->factor_mallocs    = 0.0;
123*3423f386SBarry Smith   } else {
124*3423f386SBarry Smith     info->fill_ratio_given  = 0;
125*3423f386SBarry Smith     info->fill_ratio_needed = 0;
126*3423f386SBarry Smith     info->factor_mallocs    = 0;
127*3423f386SBarry Smith   }
128*3423f386SBarry Smith   PetscFunctionReturn(0);
129*3423f386SBarry Smith }
130*3423f386SBarry Smith 
131*3423f386SBarry Smith /*@
132*3423f386SBarry Smith    MatCreateConstantDiagonal - Creates a matrix with a uniform value along the diagonal
133*3423f386SBarry Smith 
134*3423f386SBarry Smith    Collective on MPI_Comm
135*3423f386SBarry Smith 
136*3423f386SBarry Smith    Input Parameters:
137*3423f386SBarry Smith +  comm - MPI communicator
138*3423f386SBarry Smith .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
139*3423f386SBarry Smith            This value should be the same as the local size used in creating the
140*3423f386SBarry Smith            y vector for the matrix-vector product y = Ax.
141*3423f386SBarry Smith .  n - This value should be the same as the local size used in creating the
142*3423f386SBarry Smith        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
143*3423f386SBarry Smith        calculated if N is given) For square matrices n is almost always m.
144*3423f386SBarry Smith .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
145*3423f386SBarry Smith .  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
146*3423f386SBarry Smith -  diag - the diagonal value
147*3423f386SBarry Smith 
148*3423f386SBarry Smith    Output Parameter:
149*3423f386SBarry Smith .  J - the diagonal matrix
150*3423f386SBarry Smith 
151*3423f386SBarry Smith    Level: advanced
152*3423f386SBarry Smith 
153*3423f386SBarry Smith    Notes:
154*3423f386SBarry Smith     Only supports square matrices with the same number of local rows and columns
155*3423f386SBarry Smith 
156*3423f386SBarry Smith .seealso: MatDestroy(), MATCONSTANTDIAGONAL, MatScale(), MatShift(), MatMult(), MatGetDiagonal(), MatGetFactor(), MatSolve()
157*3423f386SBarry Smith 
158*3423f386SBarry Smith @*/
159*3423f386SBarry Smith PetscErrorCode  MatCreateConstantDiagonal(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar diag,Mat *J)
160*3423f386SBarry Smith {
161*3423f386SBarry Smith   PetscErrorCode ierr;
162*3423f386SBarry Smith 
163*3423f386SBarry Smith   PetscFunctionBegin;
164*3423f386SBarry Smith   ierr = MatCreate(comm,J);CHKERRQ(ierr);
165*3423f386SBarry Smith   ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr);
166*3423f386SBarry Smith   ierr = MatSetType(*J,MATCONSTANTDIAGONAL);CHKERRQ(ierr);
167*3423f386SBarry Smith   ierr = MatShift(*J,diag);CHKERRQ(ierr);
168*3423f386SBarry Smith   ierr = MatSetUp(*J);CHKERRQ(ierr);
169*3423f386SBarry Smith   PetscFunctionReturn(0);
170*3423f386SBarry Smith }
171*3423f386SBarry Smith 
172*3423f386SBarry Smith PETSC_EXTERN PetscErrorCode  MatCreate_ConstantDiagonal(Mat A)
173*3423f386SBarry Smith {
174*3423f386SBarry Smith   PetscErrorCode       ierr;
175*3423f386SBarry Smith   Mat_ConstantDiagonal *ctx;
176*3423f386SBarry Smith 
177*3423f386SBarry Smith   PetscFunctionBegin;
178*3423f386SBarry Smith   ierr = PetscNew(&ctx);CHKERRQ(ierr);
179*3423f386SBarry Smith   ctx->diag = 0.0;
180*3423f386SBarry Smith   A->data   = (void*)ctx;
181*3423f386SBarry Smith 
182*3423f386SBarry Smith   A->assembled        = PETSC_TRUE;
183*3423f386SBarry Smith   A->preallocated     = PETSC_TRUE;
184*3423f386SBarry Smith   A->ops->mult        = MatMult_ConstantDiagonal;
185*3423f386SBarry Smith   A->ops->sor         = MatSOR_ConstantDiagonal;
186*3423f386SBarry Smith   A->ops->shift       = MatShift_ConstantDiagonal;
187*3423f386SBarry Smith   A->ops->scale       = MatScale_ConstantDiagonal;
188*3423f386SBarry Smith   A->ops->getdiagonal = MatGetDiagonal_ConstantDiagonal;
189*3423f386SBarry Smith   A->ops->view        = MatView_ConstantDiagonal;
190*3423f386SBarry Smith   A->ops->zeroentries = MatZeroEntries_ConstantDiagonal;
191*3423f386SBarry Smith   A->ops->assemblyend = MatAssemblyEnd_ConstantDiagonal;
192*3423f386SBarry Smith   A->ops->destroy     = MatDestroy_ConstantDiagonal;
193*3423f386SBarry Smith   A->ops->getinfo     = MatGetInfo_ConstantDiagonal;
194*3423f386SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATCONSTANTDIAGONAL);CHKERRQ(ierr);
195*3423f386SBarry Smith   PetscFunctionReturn(0);
196*3423f386SBarry Smith }
197*3423f386SBarry Smith 
198*3423f386SBarry Smith static PetscErrorCode MatFactorNumeric_ConstantDiagonal(Mat fact,Mat A,const MatFactorInfo *info)
199*3423f386SBarry Smith {
200*3423f386SBarry Smith   Mat_ConstantDiagonal *actx = (Mat_ConstantDiagonal*)A->data,*fctx = (Mat_ConstantDiagonal*)fact->data;
201*3423f386SBarry Smith 
202*3423f386SBarry Smith   PetscFunctionBegin;
203*3423f386SBarry Smith   if (actx->diag == 0.0) fact->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
204*3423f386SBarry Smith   else fact->factorerrortype = MAT_FACTOR_NOERROR;
205*3423f386SBarry Smith   fctx->diag = 1.0/actx->diag;
206*3423f386SBarry Smith   PetscObjectStateIncrease((PetscObject)fact);
207*3423f386SBarry Smith   fact->ops->solve = MatMult_ConstantDiagonal;
208*3423f386SBarry Smith   PetscFunctionReturn(0);
209*3423f386SBarry Smith }
210*3423f386SBarry Smith 
211*3423f386SBarry Smith static PetscErrorCode MatFactorSymbolic_LU_ConstantDiagonal(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
212*3423f386SBarry Smith {
213*3423f386SBarry Smith   PetscFunctionBegin;
214*3423f386SBarry Smith   fact->ops->lufactornumeric = MatFactorNumeric_ConstantDiagonal;
215*3423f386SBarry Smith   PetscFunctionReturn(0);
216*3423f386SBarry Smith }
217*3423f386SBarry Smith 
218*3423f386SBarry Smith static PetscErrorCode MatFactorSymbolic_Cholesky_ConstantDiagonal(Mat fact,Mat A,IS isrow,const MatFactorInfo *info)
219*3423f386SBarry Smith {
220*3423f386SBarry Smith   PetscFunctionBegin;
221*3423f386SBarry Smith   fact->ops->choleskyfactornumeric = MatFactorNumeric_ConstantDiagonal;  PetscFunctionReturn(0);
222*3423f386SBarry Smith   PetscFunctionReturn(0);
223*3423f386SBarry Smith }
224*3423f386SBarry Smith 
225*3423f386SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_constantdiagonal_petsc(Mat A,MatFactorType ftype,Mat *B)
226*3423f386SBarry Smith {
227*3423f386SBarry Smith   PetscInt       n = A->rmap->n, N = A->rmap->N;
228*3423f386SBarry Smith   PetscErrorCode ierr;
229*3423f386SBarry Smith 
230*3423f386SBarry Smith   PetscFunctionBegin;
231*3423f386SBarry Smith   ierr = MatCreateConstantDiagonal(PetscObjectComm((PetscObject)A),n,n,N,N,0,B);CHKERRQ(ierr);
232*3423f386SBarry Smith 
233*3423f386SBarry Smith   (*B)->factortype = ftype;
234*3423f386SBarry Smith   (*B)->ops->ilufactorsymbolic      = MatFactorSymbolic_LU_ConstantDiagonal;
235*3423f386SBarry Smith   (*B)->ops->lufactorsymbolic       = MatFactorSymbolic_LU_ConstantDiagonal;
236*3423f386SBarry Smith   (*B)->ops->iccfactorsymbolic      = MatFactorSymbolic_Cholesky_ConstantDiagonal;
237*3423f386SBarry Smith   (*B)->ops->choleskyfactorsymbolic = MatFactorSymbolic_Cholesky_ConstantDiagonal;
238*3423f386SBarry Smith 
239*3423f386SBarry Smith   (*B)->ops->shift       = NULL;
240*3423f386SBarry Smith   (*B)->ops->scale       = NULL;
241*3423f386SBarry Smith   (*B)->ops->mult        = NULL;
242*3423f386SBarry Smith   (*B)->ops->sor         = NULL;
243*3423f386SBarry Smith   (*B)->ops->zeroentries = NULL;
244*3423f386SBarry Smith 
245*3423f386SBarry Smith   ierr = PetscFree((*B)->solvertype);CHKERRQ(ierr);
246*3423f386SBarry Smith   ierr = PetscStrallocpy(MATSOLVERPETSC,&(*B)->solvertype);CHKERRQ(ierr);
247*3423f386SBarry Smith   PetscFunctionReturn(0);
248*3423f386SBarry Smith }
249