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