xref: /petsc/src/mat/tests/ex202.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
1c4762a1bSJed Brown static char help[] = "Tests the use of MatTranspose_Nest and MatMatMult_Nest_Dense\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown PetscErrorCode TestInitialMatrix(void)
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   const PetscInt  nr = 2,nc = 3,nk = 10;
8c4762a1bSJed Brown   PetscInt        n,N,m,M;
9c4762a1bSJed Brown   const PetscInt  arow[2*3] = { 2,2,2,3,3,3 };
10c4762a1bSJed Brown   const PetscInt  acol[2*3] = { 3,2,4,3,2,4 };
11c4762a1bSJed Brown   Mat             A,Atranspose,B,C;
12c4762a1bSJed Brown   Mat             subs[2*3],**block;
13c4762a1bSJed Brown   Vec             x,y,Ax,ATy;
14c4762a1bSJed Brown   PetscInt        i,j;
15c20d7725SJed Brown   PetscScalar     dot1,dot2,zero = 0.0,one = 1.0,*valsB,*valsC;
16c4762a1bSJed Brown   PetscReal       norm;
17c4762a1bSJed Brown   PetscRandom     rctx;
18c4762a1bSJed Brown   PetscErrorCode  ierr;
19c20d7725SJed Brown   PetscBool       equal;
20c4762a1bSJed Brown 
21c4762a1bSJed Brown   PetscFunctionBegin;
22c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rctx);CHKERRQ(ierr);
23c4762a1bSJed Brown   /* Force the random numbers to have imaginary part 0 so printed results are the same for --with-scalar-type=real or --with-scalar-type=complex */
24c4762a1bSJed Brown   ierr = PetscRandomSetInterval(rctx,zero,one);CHKERRQ(ierr);
25c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
26c4762a1bSJed Brown   for (i=0; i<(nr * nc); i++) {
27c4762a1bSJed Brown     ierr = MatCreateSeqDense(PETSC_COMM_WORLD,arow[i],acol[i],NULL,&subs[i]);CHKERRQ(ierr);
28c4762a1bSJed Brown   }
29c4762a1bSJed Brown   ierr = MatCreateNest(PETSC_COMM_WORLD,nr,NULL,nc,NULL,subs,&A);CHKERRQ(ierr);
30c4762a1bSJed Brown   ierr = MatCreateVecs(A, &x, NULL);CHKERRQ(ierr);
31c4762a1bSJed Brown   ierr = MatCreateVecs(A, NULL, &y);CHKERRQ(ierr);
32c4762a1bSJed Brown   ierr = VecDuplicate(x, &ATy);CHKERRQ(ierr);
33c4762a1bSJed Brown   ierr = VecDuplicate(y, &Ax);CHKERRQ(ierr);
34c4762a1bSJed Brown   ierr = MatSetRandom(A,rctx);CHKERRQ(ierr);
35c4762a1bSJed Brown   ierr = MatTranspose(A, MAT_INITIAL_MATRIX, &Atranspose);CHKERRQ(ierr);
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   ierr = MatView(A, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
38c4762a1bSJed Brown   ierr = MatNestGetSubMats(A, NULL, NULL, &block);CHKERRQ(ierr);
39c4762a1bSJed Brown   for (i=0; i<nr; i++) {
40c4762a1bSJed Brown     for (j=0; j<nc; j++) {
41c4762a1bSJed Brown       ierr = MatView(block[i][j], PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
42c4762a1bSJed Brown     }
43c4762a1bSJed Brown   }
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   ierr = MatView(Atranspose, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
46c4762a1bSJed Brown   ierr = MatNestGetSubMats(Atranspose, NULL, NULL, &block);CHKERRQ(ierr);
47c4762a1bSJed Brown   for (i=0; i<nc; i++) {
48c4762a1bSJed Brown     for (j=0; j<nr; j++) {
49c4762a1bSJed Brown       ierr = MatView(block[i][j], PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
50c4762a1bSJed Brown     }
51c4762a1bSJed Brown   }
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /* Check <Ax, y> = <x, A^Ty> */
54c4762a1bSJed Brown   for (i=0; i<10; i++) {
55c4762a1bSJed Brown     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
56c4762a1bSJed Brown     ierr = VecSetRandom(y,rctx);CHKERRQ(ierr);
57c4762a1bSJed Brown 
58c4762a1bSJed Brown     ierr = MatMult(A, x, Ax);CHKERRQ(ierr);
59c4762a1bSJed Brown     ierr = VecDot(Ax, y, &dot1);CHKERRQ(ierr);
60c4762a1bSJed Brown     ierr = MatMult(Atranspose, y, ATy);CHKERRQ(ierr);
61c4762a1bSJed Brown     ierr = VecDot(ATy, x, &dot2);CHKERRQ(ierr);
62c4762a1bSJed Brown 
63c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD, "<Ax, y> = %g\n", (double)PetscRealPart(dot1));CHKERRQ(ierr);
64c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD, "<x, A^Ty> = %g\n",(double)PetscRealPart(dot2));CHKERRQ(ierr);
65c4762a1bSJed Brown   }
66c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
67c4762a1bSJed Brown   ierr = VecDestroy(&y);CHKERRQ(ierr);
68c4762a1bSJed Brown   ierr = VecDestroy(&Ax);CHKERRQ(ierr);
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   ierr = MatCreateSeqDense(PETSC_COMM_WORLD,acol[0]+acol[nr]+acol[2*nr],nk,NULL,&B);CHKERRQ(ierr);
71c4762a1bSJed Brown   ierr = MatSetRandom(B,rctx);CHKERRQ(ierr);
72c4762a1bSJed Brown   ierr = MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr);
73c4762a1bSJed Brown   ierr = MatMatMult(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr);
74c20d7725SJed Brown   ierr = MatMatMultEqual(A,B,C,10,&equal);CHKERRQ(ierr);
75*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in C != A*B");
76c20d7725SJed Brown 
77c4762a1bSJed Brown   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
78c4762a1bSJed Brown   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
79c4762a1bSJed Brown   for (i=0; i<nk; i++) {
80c4762a1bSJed Brown     ierr = MatDenseGetColumn(B,i,&valsB);CHKERRQ(ierr);
81c4762a1bSJed Brown     ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,N,valsB,&x);CHKERRQ(ierr);
82c4762a1bSJed Brown     ierr = MatCreateVecs(A,NULL,&Ax);CHKERRQ(ierr);
83c4762a1bSJed Brown     ierr = MatMult(A,x,Ax);CHKERRQ(ierr);
84c4762a1bSJed Brown     ierr = VecDestroy(&x);CHKERRQ(ierr);
85c4762a1bSJed Brown     ierr = MatDenseGetColumn(C,i,&valsC);CHKERRQ(ierr);
86c4762a1bSJed Brown     ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,m,M,valsC,&y);CHKERRQ(ierr);
87c4762a1bSJed Brown     ierr = VecAXPY(y,-1.0,Ax);CHKERRQ(ierr);
88c4762a1bSJed Brown     ierr = VecDestroy(&Ax);CHKERRQ(ierr);
89c4762a1bSJed Brown     ierr = VecDestroy(&y);CHKERRQ(ierr);
90c4762a1bSJed Brown     ierr = MatDenseRestoreColumn(C,&valsC);CHKERRQ(ierr);
91c4762a1bSJed Brown     ierr = MatDenseRestoreColumn(B,&valsB);CHKERRQ(ierr);
92c4762a1bSJed Brown   }
93c4762a1bSJed Brown   ierr = MatNorm(C,NORM_INFINITY,&norm);CHKERRQ(ierr);
94*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(norm > PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in MatMatMult(): %g",(double)norm);
95c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
96c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   for (i=0; i<(nr * nc); i++) {
99c4762a1bSJed Brown     ierr = MatDestroy(&subs[i]);CHKERRQ(ierr);
100c4762a1bSJed Brown   }
101c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
102c4762a1bSJed Brown   ierr = MatDestroy(&Atranspose);CHKERRQ(ierr);
103c4762a1bSJed Brown   ierr = VecDestroy(&ATy);CHKERRQ(ierr);
104c4762a1bSJed Brown   ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
105c4762a1bSJed Brown   PetscFunctionReturn(0);
106c4762a1bSJed Brown }
107c4762a1bSJed Brown 
108c4762a1bSJed Brown PetscErrorCode TestReuseMatrix(void)
109c4762a1bSJed Brown {
110c4762a1bSJed Brown   const PetscInt  n = 2;
111c4762a1bSJed Brown   Mat             A;
112c4762a1bSJed Brown   Mat             subs[2*2],**block;
113c4762a1bSJed Brown   PetscInt        i,j;
114c4762a1bSJed Brown   PetscRandom     rctx;
115c4762a1bSJed Brown   PetscErrorCode  ierr;
116c4762a1bSJed Brown   PetscScalar     zero = 0.0, one = 1.0;
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   PetscFunctionBegin;
119c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rctx);CHKERRQ(ierr);
120c4762a1bSJed Brown   ierr = PetscRandomSetInterval(rctx,zero,one);CHKERRQ(ierr);
121c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
122c4762a1bSJed Brown   for (i=0; i<(n * n); i++) {
123c4762a1bSJed Brown     ierr = MatCreateSeqDense(PETSC_COMM_WORLD,n,n,NULL,&subs[i]);CHKERRQ(ierr);
124c4762a1bSJed Brown   }
125c4762a1bSJed Brown   ierr = MatCreateNest(PETSC_COMM_WORLD,n,NULL,n,NULL,subs,&A);CHKERRQ(ierr);
126c4762a1bSJed Brown   ierr = MatSetRandom(A,rctx);CHKERRQ(ierr);
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   ierr = MatView(A, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
129c4762a1bSJed Brown   ierr = MatNestGetSubMats(A, NULL, NULL, &block);CHKERRQ(ierr);
130c4762a1bSJed Brown   for (i=0; i<n; i++) {
131c4762a1bSJed Brown     for (j=0; j<n; j++) {
132c4762a1bSJed Brown       ierr = MatView(block[i][j], PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
133c4762a1bSJed Brown     }
134c4762a1bSJed Brown   }
135c4762a1bSJed Brown   ierr = MatTranspose(A,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr);
136c4762a1bSJed Brown   ierr = MatView(A, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
137c4762a1bSJed Brown   ierr = MatNestGetSubMats(A, NULL, NULL, &block);CHKERRQ(ierr);
138c4762a1bSJed Brown   for (i=0; i<n; i++) {
139c4762a1bSJed Brown     for (j=0; j<n; j++) {
140c4762a1bSJed Brown       ierr = MatView(block[i][j], PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
141c4762a1bSJed Brown     }
142c4762a1bSJed Brown   }
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   for (i=0; i<(n * n); i++) {
145c4762a1bSJed Brown     ierr = MatDestroy(&subs[i]);CHKERRQ(ierr);
146c4762a1bSJed Brown   }
147c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
148c4762a1bSJed Brown   ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
149c4762a1bSJed Brown   PetscFunctionReturn(0);
150c4762a1bSJed Brown }
151c4762a1bSJed Brown 
152c4762a1bSJed Brown int main(int argc,char **args)
153c4762a1bSJed Brown {
154c4762a1bSJed Brown   PetscErrorCode      ierr;
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
157c4762a1bSJed Brown   ierr = TestInitialMatrix();CHKERRQ(ierr);
158c4762a1bSJed Brown   ierr = TestReuseMatrix();CHKERRQ(ierr);
159c4762a1bSJed Brown   ierr = PetscFinalize();
160c4762a1bSJed Brown   return ierr;
161c4762a1bSJed Brown }
162c4762a1bSJed Brown 
163c4762a1bSJed Brown /*TEST
164c4762a1bSJed Brown 
165c4762a1bSJed Brown    test:
166c4762a1bSJed Brown       args: -malloc_dump
167c4762a1bSJed Brown 
168c4762a1bSJed Brown TEST*/
169