xref: /petsc/src/tao/matrix/adamat.c (revision 4fc747eaadbeca11629f314a99edccbc2ed7b3d3)
1 #include <petscmat.h>              /*I  "mat.h"  I*/
2 
3 typedef struct{
4   Mat      A;
5   Vec      D1;
6   Vec      D2;
7   Vec      W;
8   Vec      W2;
9   Vec      ADADiag;
10   PetscInt GotDiag;
11 } _p_TaoMatADACtx;
12 typedef  _p_TaoMatADACtx* TaoMatADACtx;
13 
14 extern PetscErrorCode MatCreateADA(Mat,Vec, Vec, Mat*);
15 
16 #undef __FUNCT__
17 #define __FUNCT__ "MatMult_ADA"
18 static PetscErrorCode MatMult_ADA(Mat mat,Vec a,Vec y)
19 {
20   TaoMatADACtx   ctx;
21   PetscReal      one = 1.0;
22   PetscErrorCode ierr;
23 
24   PetscFunctionBegin;
25   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
26   ierr = MatMult(ctx->A,a,ctx->W);CHKERRQ(ierr);
27   if (ctx->D1){
28     ierr = VecPointwiseMult(ctx->W,ctx->D1,ctx->W);CHKERRQ(ierr);
29   }
30   ierr = MatMultTranspose(ctx->A,ctx->W,y);CHKERRQ(ierr);
31   if (ctx->D2){
32     ierr = VecPointwiseMult(ctx->W2, ctx->D2, a);CHKERRQ(ierr);
33     ierr = VecAXPY(y, one, ctx->W2);CHKERRQ(ierr);
34   }
35   PetscFunctionReturn(0);
36 }
37 
38 #undef __FUNCT__
39 #define __FUNCT__ "MatMultTranspose_ADA"
40 static PetscErrorCode MatMultTranspose_ADA(Mat mat,Vec a,Vec y)
41 {
42   PetscErrorCode ierr;
43 
44   PetscFunctionBegin;
45   ierr = MatMult_ADA(mat,a,y);CHKERRQ(ierr);
46   PetscFunctionReturn(0);
47 }
48 
49 #undef __FUNCT__
50 #define __FUNCT__ "MatDiagonalSet_ADA"
51 static PetscErrorCode MatDiagonalSet_ADA(Vec D, Mat M)
52 {
53   TaoMatADACtx   ctx;
54   PetscReal      zero=0.0,one = 1.0;
55   PetscErrorCode ierr;
56 
57   PetscFunctionBegin;
58   ierr = MatShellGetContext(M,(void **)&ctx);CHKERRQ(ierr);
59   if (!ctx->D2){
60     ierr = VecDuplicate(D,&ctx->D2);CHKERRQ(ierr);
61     ierr = VecSet(ctx->D2, zero);CHKERRQ(ierr);
62   }
63   ierr = VecAXPY(ctx->D2, one, D);CHKERRQ(ierr);
64   PetscFunctionReturn(0);
65 }
66 
67 #undef __FUNCT__
68 #define __FUNCT__ "MatDestroy_ADA"
69 static PetscErrorCode MatDestroy_ADA(Mat mat)
70 {
71   PetscErrorCode ierr;
72   TaoMatADACtx   ctx;
73 
74   PetscFunctionBegin;
75   ierr=MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
76   ierr=VecDestroy(&ctx->W);CHKERRQ(ierr);
77   ierr=VecDestroy(&ctx->W2);CHKERRQ(ierr);
78   ierr=VecDestroy(&ctx->ADADiag);CHKERRQ(ierr);
79   ierr=MatDestroy(&ctx->A);CHKERRQ(ierr);
80   ierr=VecDestroy(&ctx->D1);CHKERRQ(ierr);
81   ierr=VecDestroy(&ctx->D2);CHKERRQ(ierr);
82   ierr = PetscFree(ctx);CHKERRQ(ierr);
83   PetscFunctionReturn(0);
84 }
85 
86 #undef __FUNCT__
87 #define __FUNCT__ "MatView_ADA"
88 static PetscErrorCode MatView_ADA(Mat mat,PetscViewer viewer)
89 {
90   PetscFunctionBegin;
91   PetscFunctionReturn(0);
92 }
93 
94 #undef __FUNCT__
95 #define __FUNCT__ "MatShift_ADA"
96 static PetscErrorCode MatShift_ADA(Mat Y, PetscReal a)
97 {
98   PetscErrorCode ierr;
99   TaoMatADACtx   ctx;
100 
101   PetscFunctionBegin;
102   ierr = MatShellGetContext(Y,(void **)&ctx);CHKERRQ(ierr);
103   ierr = VecShift(ctx->D2,a);CHKERRQ(ierr);
104   PetscFunctionReturn(0);
105 }
106 
107 #undef __FUNCT__
108 #define __FUNCT__ "MatDuplicate_ADA"
109 static PetscErrorCode MatDuplicate_ADA(Mat mat,MatDuplicateOption op,Mat *M)
110 {
111   PetscErrorCode    ierr;
112   TaoMatADACtx      ctx;
113   Mat               A2;
114   Vec               D1b=NULL,D2b;
115 
116   PetscFunctionBegin;
117   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
118   ierr = MatDuplicate(ctx->A,op,&A2);CHKERRQ(ierr);
119   if (ctx->D1){
120     ierr = VecDuplicate(ctx->D1,&D1b);CHKERRQ(ierr);
121     ierr = VecCopy(ctx->D1,D1b);CHKERRQ(ierr);
122   }
123   ierr = VecDuplicate(ctx->D2,&D2b);CHKERRQ(ierr);
124   ierr = VecCopy(ctx->D2,D2b);CHKERRQ(ierr);
125   ierr = MatCreateADA(A2,D1b,D2b,M);CHKERRQ(ierr);
126   if (ctx->D1){
127     ierr = PetscObjectDereference((PetscObject)D1b);CHKERRQ(ierr);
128   }
129   ierr = PetscObjectDereference((PetscObject)D2b);CHKERRQ(ierr);
130   ierr = PetscObjectDereference((PetscObject)A2);CHKERRQ(ierr);
131   PetscFunctionReturn(0);
132 }
133 
134 #undef __FUNCT__
135 #define __FUNCT__ "MatEqual_ADA"
136 static PetscErrorCode MatEqual_ADA(Mat A,Mat B,PetscBool *flg)
137 {
138   PetscErrorCode ierr;
139   TaoMatADACtx   ctx1,ctx2;
140 
141   PetscFunctionBegin;
142   ierr = MatShellGetContext(A,(void **)&ctx1);CHKERRQ(ierr);
143   ierr = MatShellGetContext(B,(void **)&ctx2);CHKERRQ(ierr);
144   ierr = VecEqual(ctx1->D2,ctx2->D2,flg);CHKERRQ(ierr);
145   if (*flg==PETSC_TRUE){
146     ierr = VecEqual(ctx1->D1,ctx2->D1,flg);CHKERRQ(ierr);
147   }
148   if (*flg==PETSC_TRUE){
149     ierr = MatEqual(ctx1->A,ctx2->A,flg);CHKERRQ(ierr);
150   }
151   PetscFunctionReturn(0);
152 }
153 
154 #undef __FUNCT__
155 #define __FUNCT__ "MatScale_ADA"
156 static PetscErrorCode MatScale_ADA(Mat mat, PetscReal a)
157 {
158   PetscErrorCode ierr;
159   TaoMatADACtx   ctx;
160 
161   PetscFunctionBegin;
162   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
163   ierr = VecScale(ctx->D1,a);CHKERRQ(ierr);
164   if (ctx->D2){
165     ierr = VecScale(ctx->D2,a);CHKERRQ(ierr);
166   }
167   PetscFunctionReturn(0);
168 }
169 
170 #undef __FUNCT__
171 #define __FUNCT__ "MatTranspose_ADA"
172 static PetscErrorCode MatTranspose_ADA(Mat mat,Mat *B)
173 {
174   PetscErrorCode ierr;
175   TaoMatADACtx   ctx;
176 
177   PetscFunctionBegin;
178   if (*B){
179     ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
180     ierr = MatDuplicate(mat,MAT_COPY_VALUES,B);CHKERRQ(ierr);
181   }
182   PetscFunctionReturn(0);
183 }
184 
185 #undef __FUNCT__
186 #define __FUNCT__ "MatADAComputeDiagonal"
187 PetscErrorCode MatADAComputeDiagonal(Mat mat)
188 {
189   PetscErrorCode ierr;
190   PetscInt       i,m,n,low,high;
191   PetscScalar    *dtemp,*dptr;
192   TaoMatADACtx   ctx;
193 
194   PetscFunctionBegin;
195   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
196   ierr = MatGetOwnershipRange(mat, &low, &high);CHKERRQ(ierr);
197   ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr);
198 
199   ierr = PetscMalloc1(n,&dtemp );CHKERRQ(ierr);
200 
201   for (i=0; i<n; i++){
202     ierr = MatGetColumnVector(ctx->A, ctx->W, i);CHKERRQ(ierr);
203     ierr = VecPointwiseMult(ctx->W,ctx->W,ctx->W);CHKERRQ(ierr);
204     ierr = VecDotBegin(ctx->D1, ctx->W,dtemp+i);CHKERRQ(ierr);
205   }
206   for (i=0; i<n; i++){
207     ierr = VecDotEnd(ctx->D1, ctx->W,dtemp+i);CHKERRQ(ierr);
208   }
209 
210   ierr = VecGetArray(ctx->ADADiag,&dptr);CHKERRQ(ierr);
211   for (i=low; i<high; i++){
212     dptr[i-low]= dtemp[i];
213   }
214   ierr = VecRestoreArray(ctx->ADADiag,&dptr);CHKERRQ(ierr);
215   ierr = PetscFree(dtemp);CHKERRQ(ierr);
216   PetscFunctionReturn(0);
217 }
218 
219 #undef __FUNCT__
220 #define __FUNCT__ "MatGetDiagonal_ADA"
221 static PetscErrorCode MatGetDiagonal_ADA(Mat mat,Vec v)
222 {
223   PetscErrorCode  ierr;
224   PetscReal       one=1.0;
225   TaoMatADACtx    ctx;
226 
227   PetscFunctionBegin;
228   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
229   ierr = MatADAComputeDiagonal(mat);CHKERRQ(ierr);
230   ierr=VecCopy(ctx->ADADiag,v);CHKERRQ(ierr);
231   if (ctx->D2){
232     ierr=VecAXPY(v, one, ctx->D2);CHKERRQ(ierr);
233   }
234   PetscFunctionReturn(0);
235 }
236 
237 #undef __FUNCT__
238 #define __FUNCT__ "MatGetSubMatrix_ADA"
239 static PetscErrorCode MatGetSubMatrix_ADA(Mat mat,IS isrow,IS iscol,MatReuse cll, Mat *newmat)
240 {
241   PetscErrorCode    ierr;
242   PetscInt          low,high;
243   IS                ISrow;
244   Vec               D1,D2;
245   Mat               Atemp;
246   TaoMatADACtx      ctx;
247   PetscBool         isequal;
248 
249   PetscFunctionBegin;
250   ierr = ISEqual(isrow,iscol,&isequal);CHKERRQ(ierr);
251   if (!isequal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for idential column and row indices");
252   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
253 
254   ierr = MatGetOwnershipRange(ctx->A,&low,&high);CHKERRQ(ierr);
255   ierr = ISCreateStride(PetscObjectComm((PetscObject)mat),high-low,low,1,&ISrow);CHKERRQ(ierr);
256   ierr = MatGetSubMatrix(ctx->A,ISrow,iscol,cll,&Atemp);CHKERRQ(ierr);
257   ierr = ISDestroy(&ISrow);CHKERRQ(ierr);
258 
259   if (ctx->D1){
260     ierr=VecDuplicate(ctx->D1,&D1);CHKERRQ(ierr);
261     ierr=VecCopy(ctx->D1,D1);CHKERRQ(ierr);
262   } else {
263     D1 = NULL;
264   }
265 
266   if (ctx->D2){
267     Vec D2sub;
268 
269     ierr=VecGetSubVector(ctx->D2,isrow,&D2sub);CHKERRQ(ierr);
270     ierr=VecDuplicate(D2sub,&D2);CHKERRQ(ierr);
271     ierr=VecCopy(D2sub,D2);CHKERRQ(ierr);
272     ierr=VecRestoreSubVector(ctx->D2,isrow,&D2sub);CHKERRQ(ierr);
273   } else {
274     D2 = NULL;
275   }
276 
277   ierr = MatCreateADA(Atemp,D1,D2,newmat);CHKERRQ(ierr);
278   ierr = MatShellGetContext(*newmat,(void **)&ctx);CHKERRQ(ierr);
279   ierr = PetscObjectDereference((PetscObject)Atemp);CHKERRQ(ierr);
280   if (ctx->D1){
281     ierr = PetscObjectDereference((PetscObject)D1);CHKERRQ(ierr);
282   }
283   if (ctx->D2){
284     ierr = PetscObjectDereference((PetscObject)D2);CHKERRQ(ierr);
285   }
286   PetscFunctionReturn(0);
287 }
288 
289 #undef __FUNCT__
290 #define __FUNCT__ "MatGetSubMatrices_ADA"
291 static PetscErrorCode MatGetSubMatrices_ADA(Mat A,PetscInt n, IS *irow,IS *icol,MatReuse scall,Mat **B)
292 {
293   PetscErrorCode ierr;
294   PetscInt       i;
295 
296   PetscFunctionBegin;
297   if (scall == MAT_INITIAL_MATRIX) {
298     ierr = PetscMalloc1(n+1,B );CHKERRQ(ierr);
299   }
300   for ( i=0; i<n; i++ ) {
301     ierr = MatGetSubMatrix_ADA(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
302   }
303   PetscFunctionReturn(0);
304 }
305 
306 #undef __FUNCT__
307 #define __FUNCT__ "MatGetColumnVector_ADA"
308 static PetscErrorCode MatGetColumnVector_ADA(Mat mat,Vec Y, PetscInt col)
309 {
310   PetscErrorCode ierr;
311   PetscInt       low,high;
312   PetscScalar    zero=0.0,one=1.0;
313 
314   PetscFunctionBegin;
315   ierr=VecSet(Y, zero);CHKERRQ(ierr);
316   ierr=VecGetOwnershipRange(Y,&low,&high);CHKERRQ(ierr);
317   if (col>=low && col<high){
318     ierr=VecSetValue(Y,col,one,INSERT_VALUES);CHKERRQ(ierr);
319   }
320   ierr=VecAssemblyBegin(Y);CHKERRQ(ierr);
321   ierr=VecAssemblyEnd(Y);CHKERRQ(ierr);
322   ierr=MatMult_ADA(mat,Y,Y);CHKERRQ(ierr);
323   PetscFunctionReturn(0);
324 }
325 
326 #undef __FUNCT__
327 #define __FUNCT__ "MatConvert_ADA"
328 PETSC_INTERN PetscErrorCode MatConvert_ADA(Mat mat,MatType newtype,Mat *NewMat)
329 {
330   PetscErrorCode ierr;
331   PetscMPIInt    size;
332   PetscBool      sametype, issame, isdense, isseqdense;
333   TaoMatADACtx   ctx;
334 
335   PetscFunctionBegin;
336   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
337   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr);
338 
339   ierr = PetscObjectTypeCompare((PetscObject)mat,newtype,&sametype);CHKERRQ(ierr);
340   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSAME,&issame);CHKERRQ(ierr);
341   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIDENSE,&isdense);CHKERRQ(ierr);
342   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQDENSE,&isseqdense);CHKERRQ(ierr);
343 
344   if (sametype || issame) {
345     ierr=MatDuplicate(mat,MAT_COPY_VALUES,NewMat);CHKERRQ(ierr);
346   } else if (isdense) {
347     PetscInt          i,j,low,high,m,n,M,N;
348     const PetscScalar *dptr;
349     Vec               X;
350 
351     ierr = VecDuplicate(ctx->D2,&X);CHKERRQ(ierr);
352     ierr=MatGetSize(mat,&M,&N);CHKERRQ(ierr);
353     ierr=MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
354     ierr = MatCreateDense(PetscObjectComm((PetscObject)mat),m,m,N,N,NULL,NewMat);CHKERRQ(ierr);
355     ierr = MatGetOwnershipRange(*NewMat,&low,&high);CHKERRQ(ierr);
356     for (i=0;i<M;i++){
357       ierr = MatGetColumnVector_ADA(mat,X,i);CHKERRQ(ierr);
358       ierr = VecGetArrayRead(X,&dptr);CHKERRQ(ierr);
359       for (j=0; j<high-low; j++){
360         ierr = MatSetValue(*NewMat,low+j,i,dptr[j],INSERT_VALUES);CHKERRQ(ierr);
361       }
362       ierr = VecRestoreArrayRead(X,&dptr);CHKERRQ(ierr);
363     }
364     ierr = MatAssemblyBegin(*NewMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
365     ierr = MatAssemblyEnd(*NewMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
366     ierr = VecDestroy(&X);CHKERRQ(ierr);
367   } else if (isseqdense && size==1){
368     PetscInt          i,j,low,high,m,n,M,N;
369     const PetscScalar *dptr;
370     Vec               X;
371 
372     ierr = VecDuplicate(ctx->D2,&X);CHKERRQ(ierr);
373     ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
374     ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
375     ierr = MatCreateSeqDense(PetscObjectComm((PetscObject)mat),N,N,NULL,NewMat);CHKERRQ(ierr);
376     ierr = MatGetOwnershipRange(*NewMat,&low,&high);CHKERRQ(ierr);
377     for (i=0;i<M;i++){
378       ierr = MatGetColumnVector_ADA(mat,X,i);CHKERRQ(ierr);
379       ierr = VecGetArrayRead(X,&dptr);CHKERRQ(ierr);
380       for (j=0; j<high-low; j++){
381         ierr = MatSetValue(*NewMat,low+j,i,dptr[j],INSERT_VALUES);CHKERRQ(ierr);
382       }
383       ierr = VecRestoreArrayRead(X,&dptr);CHKERRQ(ierr);
384     }
385     ierr = MatAssemblyBegin(*NewMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
386     ierr = MatAssemblyEnd(*NewMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
387     ierr = VecDestroy(&X);CHKERRQ(ierr);
388   } else SETERRQ(PETSC_COMM_SELF,1,"No support to convert objects to that type");
389   PetscFunctionReturn(0);
390 }
391 
392 #undef __FUNCT__
393 #define __FUNCT__ "MatNorm_ADA"
394 static PetscErrorCode MatNorm_ADA(Mat mat,NormType type,PetscReal *norm)
395 {
396   PetscErrorCode ierr;
397   TaoMatADACtx   ctx;
398 
399   PetscFunctionBegin;
400   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
401   if (type == NORM_FROBENIUS) {
402     *norm = 1.0;
403   } else if (type == NORM_1 || type == NORM_INFINITY) {
404     *norm = 1.0;
405   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm");
406   PetscFunctionReturn(0);
407 }
408 
409 #undef __FUNCT__
410 #define __FUNCT__ "MatCreateADA"
411 /*@C
412    MatCreateADA - Creates a matrix M=A^T D1 A + D2 where D1, D2 are diagonal
413 
414    Collective on matrix
415 
416    Input Parameters:
417 +  mat - matrix of arbitrary type
418 .  d1 - A vector with diagonal elements of D1
419 -  d2 - A vector
420 
421    Output Parameters:
422 .  J - New matrix whose operations are defined in terms of mat, D1, and D2.
423 
424    Notes:
425    The user provides the input data and is responsible for destroying
426    this data after matrix J has been destroyed.
427    The operation MatMult(A,D2,D1) must be well defined.
428    Before calling the operation MatGetDiagonal(), the function
429    MatADAComputeDiagonal() must be called.  The matrices A and D1 must
430    be the same during calls to MatADAComputeDiagonal() and
431    MatGetDiagonal().
432 
433    Level: developer
434 
435 .seealso: MatCreate()
436 @*/
437 PetscErrorCode MatCreateADA(Mat mat,Vec d1, Vec d2, Mat *J)
438 {
439   MPI_Comm       comm=PetscObjectComm((PetscObject)mat);
440   TaoMatADACtx   ctx;
441   PetscErrorCode ierr;
442   PetscInt       nloc,n;
443 
444   PetscFunctionBegin;
445   ierr = PetscNew(&ctx);CHKERRQ(ierr);
446   ctx->A=mat;
447   ctx->D1=d1;
448   ctx->D2=d2;
449   if (d1){
450     ierr = VecDuplicate(d1,&ctx->W);CHKERRQ(ierr);
451     ierr = PetscObjectReference((PetscObject)d1);CHKERRQ(ierr);
452   } else {
453     ctx->W = NULL;
454   }
455   if (d2){
456     ierr = VecDuplicate(d2,&ctx->W2);CHKERRQ(ierr);
457     ierr = VecDuplicate(d2,&ctx->ADADiag);CHKERRQ(ierr);
458     ierr =  PetscObjectReference((PetscObject)d2);CHKERRQ(ierr);
459   } else {
460     ctx->W2      = NULL;
461     ctx->ADADiag = NULL;
462   }
463 
464   ctx->GotDiag=0;
465   ierr =  PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
466 
467   ierr=VecGetLocalSize(d2,&nloc);CHKERRQ(ierr);
468   ierr=VecGetSize(d2,&n);CHKERRQ(ierr);
469 
470   ierr = MatCreateShell(comm,nloc,nloc,n,n,ctx,J);CHKERRQ(ierr);
471 
472   ierr = MatShellSetOperation(*J,MATOP_MULT,(void(*)(void))MatMult_ADA);CHKERRQ(ierr);
473   ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void(*)(void))MatDestroy_ADA);CHKERRQ(ierr);
474   ierr = MatShellSetOperation(*J,MATOP_VIEW,(void(*)(void))MatView_ADA);CHKERRQ(ierr);
475   ierr = MatShellSetOperation(*J,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_ADA);CHKERRQ(ierr);
476   ierr = MatShellSetOperation(*J,MATOP_DIAGONAL_SET,(void(*)(void))MatDiagonalSet_ADA);CHKERRQ(ierr);
477   ierr = MatShellSetOperation(*J,MATOP_SHIFT,(void(*)(void))MatShift_ADA);CHKERRQ(ierr);
478   ierr = MatShellSetOperation(*J,MATOP_EQUAL,(void(*)(void))MatEqual_ADA);CHKERRQ(ierr);
479   ierr = MatShellSetOperation(*J,MATOP_SCALE,(void(*)(void))MatScale_ADA);CHKERRQ(ierr);
480   ierr = MatShellSetOperation(*J,MATOP_TRANSPOSE,(void(*)(void))MatTranspose_ADA);CHKERRQ(ierr);
481   ierr = MatShellSetOperation(*J,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_ADA);CHKERRQ(ierr);
482   ierr = MatShellSetOperation(*J,MATOP_GET_SUBMATRICES,(void(*)(void))MatGetSubMatrices_ADA);CHKERRQ(ierr);
483   ierr = MatShellSetOperation(*J,MATOP_NORM,(void(*)(void))MatNorm_ADA);CHKERRQ(ierr);
484   ierr = MatShellSetOperation(*J,MATOP_DUPLICATE,(void(*)(void))MatDuplicate_ADA);CHKERRQ(ierr);
485   ierr = MatShellSetOperation(*J,MATOP_GET_SUBMATRIX,(void(*)(void))MatGetSubMatrix_ADA);CHKERRQ(ierr);
486 
487   ierr = PetscLogObjectParent((PetscObject)(*J),(PetscObject)ctx->W);CHKERRQ(ierr);
488   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)(*J));CHKERRQ(ierr);
489 
490   ierr = MatSetOption(*J,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
491   PetscFunctionReturn(0);
492 }
493