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