xref: /petsc/src/mat/tests/ex97.c (revision 9566063d113dddea24716c546802770db7481bc0)
1c4762a1bSJed Brown static const char help[] = "Tests MatCreateSubMatrix with MatSubMatrix versus MatAIJ, non-square\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown static PetscErrorCode AssembleMatrix(MPI_Comm comm,Mat *A)
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   Mat            B;
8c4762a1bSJed Brown   PetscInt       i,ms,me;
9c4762a1bSJed Brown 
10c4762a1bSJed Brown   PetscFunctionBegin;
11*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm,&B));
12*9566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B,5,6,PETSC_DETERMINE,PETSC_DETERMINE));
13*9566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(B));
14*9566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
15*9566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(B,&ms,&me));
16c4762a1bSJed Brown   for (i=ms; i<me; i++) {
17*9566063dSJacob Faibussowitsch     PetscCall(MatSetValue(B,i,i,1.0*i,INSERT_VALUES));
18c4762a1bSJed Brown   }
19*9566063dSJacob Faibussowitsch   PetscCall(MatSetValue(B,me-1,me,me*me,INSERT_VALUES));
20*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
21*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
22c4762a1bSJed Brown   *A   = B;
23c4762a1bSJed Brown   PetscFunctionReturn(0);
24c4762a1bSJed Brown }
25c4762a1bSJed Brown 
26c4762a1bSJed Brown static PetscErrorCode Compare2(Vec *X,const char *test)
27c4762a1bSJed Brown {
28c4762a1bSJed Brown   PetscReal      norm;
29c4762a1bSJed Brown   Vec            Y;
30c4762a1bSJed Brown   PetscInt       verbose = 0;
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   PetscFunctionBegin;
33*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X[0],&Y));
34*9566063dSJacob Faibussowitsch   PetscCall(VecCopy(X[0],Y));
35*9566063dSJacob Faibussowitsch   PetscCall(VecAYPX(Y,-1.0,X[1]));
36*9566063dSJacob Faibussowitsch   PetscCall(VecNorm(Y,NORM_INFINITY,&norm));
37c4762a1bSJed Brown 
38*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-verbose",&verbose,NULL));
39c4762a1bSJed Brown   if (norm < PETSC_SQRT_MACHINE_EPSILON && verbose < 1) {
40*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%30s: norm difference < sqrt(eps_machine)\n",test));
41c4762a1bSJed Brown   } else {
42*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%30s: norm difference %g\n",test,(double)norm));
43c4762a1bSJed Brown   }
44c4762a1bSJed Brown   if (verbose > 1) {
45*9566063dSJacob Faibussowitsch     PetscCall(VecView(X[0],PETSC_VIEWER_STDOUT_WORLD));
46*9566063dSJacob Faibussowitsch     PetscCall(VecView(X[1],PETSC_VIEWER_STDOUT_WORLD));
47*9566063dSJacob Faibussowitsch     PetscCall(VecView(Y,PETSC_VIEWER_STDOUT_WORLD));
48c4762a1bSJed Brown   }
49*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Y));
50c4762a1bSJed Brown   PetscFunctionReturn(0);
51c4762a1bSJed Brown }
52c4762a1bSJed Brown 
53c4762a1bSJed Brown static PetscErrorCode CheckMatrices(Mat A,Mat B,Vec left,Vec right,Vec X,Vec Y,Vec X1,Vec Y1)
54c4762a1bSJed Brown {
55c4762a1bSJed Brown   Vec            *ltmp,*rtmp;
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   PetscFunctionBegin;
58*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(right,2,&rtmp));
59*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(left,2,&ltmp));
60*9566063dSJacob Faibussowitsch   PetscCall(MatScale(A,PETSC_PI));
61*9566063dSJacob Faibussowitsch   PetscCall(MatScale(B,PETSC_PI));
62*9566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(A,left,right));
63*9566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(B,left,right));
64c4762a1bSJed Brown 
65*9566063dSJacob Faibussowitsch   PetscCall(MatMult(A,X,ltmp[0]));
66*9566063dSJacob Faibussowitsch   PetscCall(MatMult(B,X,ltmp[1]));
67*9566063dSJacob Faibussowitsch   PetscCall(Compare2(ltmp,"MatMult"));
68c4762a1bSJed Brown 
69*9566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(A,Y,rtmp[0]));
70*9566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(B,Y,rtmp[1]));
71*9566063dSJacob Faibussowitsch   PetscCall(Compare2(rtmp,"MatMultTranspose"));
72c4762a1bSJed Brown 
73*9566063dSJacob Faibussowitsch   PetscCall(VecCopy(Y1,ltmp[0]));
74*9566063dSJacob Faibussowitsch   PetscCall(VecCopy(Y1,ltmp[1]));
75*9566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(A,X,ltmp[0],ltmp[0]));
76*9566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(B,X,ltmp[1],ltmp[1]));
77*9566063dSJacob Faibussowitsch   PetscCall(Compare2(ltmp,"MatMultAdd v2==v3"));
78c4762a1bSJed Brown 
79*9566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(A,X,Y1,ltmp[0]));
80*9566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(B,X,Y1,ltmp[1]));
81*9566063dSJacob Faibussowitsch   PetscCall(Compare2(ltmp,"MatMultAdd v2!=v3"));
82c4762a1bSJed Brown 
83*9566063dSJacob Faibussowitsch   PetscCall(VecCopy(X1,rtmp[0]));
84*9566063dSJacob Faibussowitsch   PetscCall(VecCopy(X1,rtmp[1]));
85*9566063dSJacob Faibussowitsch   PetscCall(MatMultTransposeAdd(A,Y,rtmp[0],rtmp[0]));
86*9566063dSJacob Faibussowitsch   PetscCall(MatMultTransposeAdd(B,Y,rtmp[1],rtmp[1]));
87*9566063dSJacob Faibussowitsch   PetscCall(Compare2(rtmp,"MatMultTransposeAdd v2==v3"));
88c4762a1bSJed Brown 
89*9566063dSJacob Faibussowitsch   PetscCall(MatMultTransposeAdd(A,Y,X1,rtmp[0]));
90*9566063dSJacob Faibussowitsch   PetscCall(MatMultTransposeAdd(B,Y,X1,rtmp[1]));
91*9566063dSJacob Faibussowitsch   PetscCall(Compare2(rtmp,"MatMultTransposeAdd v2!=v3"));
92c4762a1bSJed Brown 
93*9566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(2,&ltmp));
94*9566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(2,&rtmp));
95c4762a1bSJed Brown   PetscFunctionReturn(0);
96c4762a1bSJed Brown }
97c4762a1bSJed Brown 
98c4762a1bSJed Brown int main(int argc, char *argv[])
99c4762a1bSJed Brown {
100c4762a1bSJed Brown   Mat            A,B,Asub,Bsub;
101c4762a1bSJed Brown   PetscInt       ms,idxrow[3],idxcol[4];
102c4762a1bSJed Brown   Vec            left,right,X,Y,X1,Y1;
103c4762a1bSJed Brown   IS             isrow,iscol;
104c4762a1bSJed Brown   PetscBool      random = PETSC_TRUE;
105c4762a1bSJed Brown 
106*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
107*9566063dSJacob Faibussowitsch   PetscCall(AssembleMatrix(PETSC_COMM_WORLD,&A));
108*9566063dSJacob Faibussowitsch   PetscCall(AssembleMatrix(PETSC_COMM_WORLD,&B));
109*9566063dSJacob Faibussowitsch   PetscCall(MatSetOperation(B,MATOP_CREATE_SUBMATRIX,NULL));
110*9566063dSJacob Faibussowitsch   PetscCall(MatSetOperation(B,MATOP_CREATE_SUBMATRICES,NULL));
111*9566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A,&ms,NULL));
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   idxrow[0] = ms+1;
114c4762a1bSJed Brown   idxrow[1] = ms+2;
115c4762a1bSJed Brown   idxrow[2] = ms+4;
116*9566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,3,idxrow,PETSC_USE_POINTER,&isrow));
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   idxcol[0] = ms+1;
119c4762a1bSJed Brown   idxcol[1] = ms+2;
120c4762a1bSJed Brown   idxcol[2] = ms+4;
121c4762a1bSJed Brown   idxcol[3] = ms+5;
122*9566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,4,idxcol,PETSC_USE_POINTER,&iscol));
123c4762a1bSJed Brown 
124*9566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(A,isrow,iscol,MAT_INITIAL_MATRIX,&Asub));
125*9566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(B,isrow,iscol,MAT_INITIAL_MATRIX,&Bsub));
126c4762a1bSJed Brown 
127*9566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Asub,&right,&left));
128*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(right,&X));
129*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(right,&X1));
130*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(left,&Y));
131*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(left,&Y1));
132c4762a1bSJed Brown 
133*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-random",&random,NULL));
134c4762a1bSJed Brown   if (random) {
135*9566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(right,NULL));
136*9566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(left,NULL));
137*9566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(X,NULL));
138*9566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(Y,NULL));
139*9566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(X1,NULL));
140*9566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(Y1,NULL));
141c4762a1bSJed Brown   } else {
142*9566063dSJacob Faibussowitsch     PetscCall(VecSet(right,1.0));
143*9566063dSJacob Faibussowitsch     PetscCall(VecSet(left,2.0));
144*9566063dSJacob Faibussowitsch     PetscCall(VecSet(X,3.0));
145*9566063dSJacob Faibussowitsch     PetscCall(VecSet(Y,4.0));
146*9566063dSJacob Faibussowitsch     PetscCall(VecSet(X1,3.0));
147*9566063dSJacob Faibussowitsch     PetscCall(VecSet(Y1,4.0));
148c4762a1bSJed Brown   }
149*9566063dSJacob Faibussowitsch   PetscCall(CheckMatrices(Asub,Bsub,left,right,X,Y,X1,Y1));
150*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrow));
151*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscol));
152*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
153*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
154*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Asub));
155*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Bsub));
156*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&left));
157*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&right));
158*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
159*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Y));
160*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X1));
161*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Y1));
162*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
163b122ec5aSJacob Faibussowitsch   return 0;
164c4762a1bSJed Brown }
165c4762a1bSJed Brown 
166c4762a1bSJed Brown /*TEST
167c4762a1bSJed Brown 
168c4762a1bSJed Brown    test:
169c4762a1bSJed Brown       nsize: 3
170c4762a1bSJed Brown 
171c4762a1bSJed Brown TEST*/
172