xref: /petsc/src/mat/tutorials/ex2.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "testing SeqDense matrices with an LDA (leading dimension of the user-allocated arrray) larger than M.\n";
2*c4762a1bSJed Brown /*
3*c4762a1bSJed Brown  * Example code testing SeqDense matrices with an LDA (leading dimension
4*c4762a1bSJed Brown  * of the user-allocated arrray) larger than M.
5*c4762a1bSJed Brown  * This example tests the functionality of MatSetValues(), MatMult(),
6*c4762a1bSJed Brown  * and MatMultTranspose().
7*c4762a1bSJed Brown  */
8*c4762a1bSJed Brown 
9*c4762a1bSJed Brown #include <petscmat.h>
10*c4762a1bSJed Brown 
11*c4762a1bSJed Brown int main(int argc,char **argv)
12*c4762a1bSJed Brown {
13*c4762a1bSJed Brown   Mat            A,A11,A12,A21,A22;
14*c4762a1bSJed Brown   Vec            X,X1,X2,Y,Z,Z1,Z2;
15*c4762a1bSJed Brown   PetscScalar    *a,*b,*x,*y,*z,v,one=1;
16*c4762a1bSJed Brown   PetscReal      nrm;
17*c4762a1bSJed Brown   PetscErrorCode ierr;
18*c4762a1bSJed Brown   PetscInt       size=8,size1=6,size2=2, i,j;
19*c4762a1bSJed Brown   PetscRandom    rnd;
20*c4762a1bSJed Brown 
21*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr;
22*c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_SELF,&rnd);CHKERRQ(ierr);
23*c4762a1bSJed Brown 
24*c4762a1bSJed Brown   /*
25*c4762a1bSJed Brown    * Create matrix and three vectors: these are all normal
26*c4762a1bSJed Brown    */
27*c4762a1bSJed Brown   ierr = PetscMalloc1(size*size,&a);CHKERRQ(ierr);
28*c4762a1bSJed Brown   ierr = PetscMalloc1(size*size,&b);CHKERRQ(ierr);
29*c4762a1bSJed Brown   for (i=0; i<size; i++) {
30*c4762a1bSJed Brown     for (j=0; j<size; j++) {
31*c4762a1bSJed Brown       ierr = PetscRandomGetValue(rnd,&a[i+j*size]);CHKERRQ(ierr);
32*c4762a1bSJed Brown       b[i+j*size] = a[i+j*size];
33*c4762a1bSJed Brown     }
34*c4762a1bSJed Brown   }
35*c4762a1bSJed Brown   ierr = MatCreate(MPI_COMM_SELF,&A);CHKERRQ(ierr);
36*c4762a1bSJed Brown   ierr = MatSetSizes(A,size,size,size,size);CHKERRQ(ierr);
37*c4762a1bSJed Brown   ierr = MatSetType(A,MATSEQDENSE);CHKERRQ(ierr);
38*c4762a1bSJed Brown   ierr = MatSeqDenseSetPreallocation(A,a);CHKERRQ(ierr);
39*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
40*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
41*c4762a1bSJed Brown 
42*c4762a1bSJed Brown   ierr = PetscMalloc1(size,&x);CHKERRQ(ierr);
43*c4762a1bSJed Brown   for (i=0; i<size; i++) {
44*c4762a1bSJed Brown     x[i] = one;
45*c4762a1bSJed Brown   }
46*c4762a1bSJed Brown   ierr = VecCreateSeqWithArray(MPI_COMM_SELF,1,size,x,&X);CHKERRQ(ierr);
47*c4762a1bSJed Brown   ierr = VecAssemblyBegin(X);CHKERRQ(ierr);
48*c4762a1bSJed Brown   ierr = VecAssemblyEnd(X);CHKERRQ(ierr);
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown   ierr = PetscMalloc1(size,&y);CHKERRQ(ierr);
51*c4762a1bSJed Brown   ierr = VecCreateSeqWithArray(MPI_COMM_SELF,1,size,y,&Y);CHKERRQ(ierr);
52*c4762a1bSJed Brown   ierr = VecAssemblyBegin(Y);CHKERRQ(ierr);
53*c4762a1bSJed Brown   ierr = VecAssemblyEnd(Y);CHKERRQ(ierr);
54*c4762a1bSJed Brown 
55*c4762a1bSJed Brown   ierr = PetscMalloc1(size,&z);CHKERRQ(ierr);
56*c4762a1bSJed Brown   ierr = VecCreateSeqWithArray(MPI_COMM_SELF,1,size,z,&Z);CHKERRQ(ierr);
57*c4762a1bSJed Brown   ierr = VecAssemblyBegin(Z);CHKERRQ(ierr);
58*c4762a1bSJed Brown   ierr = VecAssemblyEnd(Z);CHKERRQ(ierr);
59*c4762a1bSJed Brown 
60*c4762a1bSJed Brown   /*
61*c4762a1bSJed Brown    * Now create submatrices and subvectors
62*c4762a1bSJed Brown    */
63*c4762a1bSJed Brown   ierr = MatCreate(MPI_COMM_SELF,&A11);CHKERRQ(ierr);
64*c4762a1bSJed Brown   ierr = MatSetSizes(A11,size1,size1,size1,size1);CHKERRQ(ierr);
65*c4762a1bSJed Brown   ierr = MatSetType(A11,MATSEQDENSE);CHKERRQ(ierr);
66*c4762a1bSJed Brown   ierr = MatSeqDenseSetPreallocation(A11,b);CHKERRQ(ierr);
67*c4762a1bSJed Brown   ierr = MatSeqDenseSetLDA(A11,size);CHKERRQ(ierr);
68*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A11,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
69*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A11,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
70*c4762a1bSJed Brown 
71*c4762a1bSJed Brown   ierr = MatCreate(MPI_COMM_SELF,&A12);CHKERRQ(ierr);
72*c4762a1bSJed Brown   ierr = MatSetSizes(A12,size1,size2,size1,size2);CHKERRQ(ierr);
73*c4762a1bSJed Brown   ierr = MatSetType(A12,MATSEQDENSE);CHKERRQ(ierr);
74*c4762a1bSJed Brown   ierr = MatSeqDenseSetPreallocation(A12,b+size1*size);CHKERRQ(ierr);
75*c4762a1bSJed Brown   ierr = MatSeqDenseSetLDA(A12,size);CHKERRQ(ierr);
76*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A12,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
77*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A12,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
78*c4762a1bSJed Brown 
79*c4762a1bSJed Brown   ierr = MatCreate(MPI_COMM_SELF,&A21);CHKERRQ(ierr);
80*c4762a1bSJed Brown   ierr = MatSetSizes(A21,size2,size1,size2,size1);CHKERRQ(ierr);
81*c4762a1bSJed Brown   ierr = MatSetType(A21,MATSEQDENSE);CHKERRQ(ierr);
82*c4762a1bSJed Brown   ierr = MatSeqDenseSetPreallocation(A21,b+size1);CHKERRQ(ierr);
83*c4762a1bSJed Brown   ierr = MatSeqDenseSetLDA(A21,size);CHKERRQ(ierr);
84*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A21,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
85*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A21,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
86*c4762a1bSJed Brown 
87*c4762a1bSJed Brown   ierr = MatCreate(MPI_COMM_SELF,&A22);CHKERRQ(ierr);
88*c4762a1bSJed Brown   ierr = MatSetSizes(A22,size2,size2,size2,size2);CHKERRQ(ierr);
89*c4762a1bSJed Brown   ierr = MatSetType(A22,MATSEQDENSE);CHKERRQ(ierr);
90*c4762a1bSJed Brown   ierr = MatSeqDenseSetPreallocation(A22,b+size1*size+size1);CHKERRQ(ierr);
91*c4762a1bSJed Brown   ierr = MatSeqDenseSetLDA(A22,size);CHKERRQ(ierr);
92*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A22,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
93*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A22,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
94*c4762a1bSJed Brown 
95*c4762a1bSJed Brown   ierr = VecCreateSeqWithArray(MPI_COMM_SELF,1,size1,x,&X1);CHKERRQ(ierr);
96*c4762a1bSJed Brown   ierr = VecCreateSeqWithArray(MPI_COMM_SELF,1,size2,x+size1,&X2);CHKERRQ(ierr);
97*c4762a1bSJed Brown   ierr = VecCreateSeqWithArray(MPI_COMM_SELF,1,size1,z,&Z1);CHKERRQ(ierr);
98*c4762a1bSJed Brown   ierr = VecCreateSeqWithArray(MPI_COMM_SELF,1,size2,z+size1,&Z2);CHKERRQ(ierr);
99*c4762a1bSJed Brown 
100*c4762a1bSJed Brown   /*
101*c4762a1bSJed Brown    * Now multiple matrix times input in two ways;
102*c4762a1bSJed Brown    * compare the results
103*c4762a1bSJed Brown    */
104*c4762a1bSJed Brown   ierr = MatMult(A,X,Y);CHKERRQ(ierr);
105*c4762a1bSJed Brown   ierr = MatMult(A11,X1,Z1);CHKERRQ(ierr);
106*c4762a1bSJed Brown   ierr = MatMultAdd(A12,X2,Z1,Z1);CHKERRQ(ierr);
107*c4762a1bSJed Brown   ierr = MatMult(A22,X2,Z2);CHKERRQ(ierr);
108*c4762a1bSJed Brown   ierr = MatMultAdd(A21,X1,Z2,Z2);CHKERRQ(ierr);
109*c4762a1bSJed Brown   ierr = VecAXPY(Z,-1.0,Y);CHKERRQ(ierr);
110*c4762a1bSJed Brown   ierr = VecNorm(Z,NORM_2,&nrm);CHKERRQ(ierr);
111*c4762a1bSJed Brown   if (nrm > 100.0*PETSC_MACHINE_EPSILON) {
112*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Test1; error norm=%g\n",(double)nrm);CHKERRQ(ierr);
113*c4762a1bSJed Brown   }
114*c4762a1bSJed Brown 
115*c4762a1bSJed Brown   /*
116*c4762a1bSJed Brown    * Next test: change both matrices
117*c4762a1bSJed Brown    */
118*c4762a1bSJed Brown   ierr = PetscRandomGetValue(rnd,&v);CHKERRQ(ierr);
119*c4762a1bSJed Brown   i    = 1;
120*c4762a1bSJed Brown   j    = size-2;
121*c4762a1bSJed Brown   ierr = MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr);
122*c4762a1bSJed Brown   j   -= size1;
123*c4762a1bSJed Brown   ierr = MatSetValues(A12,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr);
124*c4762a1bSJed Brown   ierr = PetscRandomGetValue(rnd,&v);CHKERRQ(ierr);
125*c4762a1bSJed Brown   i    = j = size1+1;
126*c4762a1bSJed Brown   ierr = MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr);
127*c4762a1bSJed Brown   i     =j = 1;
128*c4762a1bSJed Brown   ierr = MatSetValues(A22,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr);
129*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
130*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
131*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A12,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
132*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A12,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
133*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A22,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
134*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A22,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
135*c4762a1bSJed Brown 
136*c4762a1bSJed Brown   ierr = MatMult(A,X,Y);CHKERRQ(ierr);
137*c4762a1bSJed Brown   ierr = MatMult(A11,X1,Z1);CHKERRQ(ierr);
138*c4762a1bSJed Brown   ierr = MatMultAdd(A12,X2,Z1,Z1);CHKERRQ(ierr);
139*c4762a1bSJed Brown   ierr = MatMult(A22,X2,Z2);CHKERRQ(ierr);
140*c4762a1bSJed Brown   ierr = MatMultAdd(A21,X1,Z2,Z2);CHKERRQ(ierr);
141*c4762a1bSJed Brown   ierr = VecAXPY(Z,-1.0,Y);CHKERRQ(ierr);
142*c4762a1bSJed Brown   ierr = VecNorm(Z,NORM_2,&nrm);CHKERRQ(ierr);
143*c4762a1bSJed Brown   if (nrm > 100.0*PETSC_MACHINE_EPSILON) {
144*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Test2; error norm=%g\n",(double)nrm);CHKERRQ(ierr);
145*c4762a1bSJed Brown   }
146*c4762a1bSJed Brown 
147*c4762a1bSJed Brown   /*
148*c4762a1bSJed Brown    * Transpose product
149*c4762a1bSJed Brown    */
150*c4762a1bSJed Brown   ierr = MatMultTranspose(A,X,Y);CHKERRQ(ierr);
151*c4762a1bSJed Brown   ierr = MatMultTranspose(A11,X1,Z1);CHKERRQ(ierr);
152*c4762a1bSJed Brown   ierr = MatMultTransposeAdd(A21,X2,Z1,Z1);CHKERRQ(ierr);
153*c4762a1bSJed Brown   ierr = MatMultTranspose(A22,X2,Z2);CHKERRQ(ierr);
154*c4762a1bSJed Brown   ierr = MatMultTransposeAdd(A12,X1,Z2,Z2);CHKERRQ(ierr);
155*c4762a1bSJed Brown   ierr = VecAXPY(Z,-1.0,Y);CHKERRQ(ierr);
156*c4762a1bSJed Brown   ierr = VecNorm(Z,NORM_2,&nrm);CHKERRQ(ierr);
157*c4762a1bSJed Brown   if (nrm > 100.0*PETSC_MACHINE_EPSILON) {
158*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Test3; error norm=%g\n",(double)nrm);CHKERRQ(ierr);
159*c4762a1bSJed Brown   }
160*c4762a1bSJed Brown   ierr = PetscFree(a);CHKERRQ(ierr);
161*c4762a1bSJed Brown   ierr = PetscFree(b);CHKERRQ(ierr);
162*c4762a1bSJed Brown   ierr = PetscFree(x);CHKERRQ(ierr);
163*c4762a1bSJed Brown   ierr = PetscFree(y);CHKERRQ(ierr);
164*c4762a1bSJed Brown   ierr = PetscFree(z);CHKERRQ(ierr);
165*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
166*c4762a1bSJed Brown   ierr = MatDestroy(&A11);CHKERRQ(ierr);
167*c4762a1bSJed Brown   ierr = MatDestroy(&A12);CHKERRQ(ierr);
168*c4762a1bSJed Brown   ierr = MatDestroy(&A21);CHKERRQ(ierr);
169*c4762a1bSJed Brown   ierr = MatDestroy(&A22);CHKERRQ(ierr);
170*c4762a1bSJed Brown 
171*c4762a1bSJed Brown   ierr = VecDestroy(&X);CHKERRQ(ierr);
172*c4762a1bSJed Brown   ierr = VecDestroy(&Y);CHKERRQ(ierr);
173*c4762a1bSJed Brown   ierr = VecDestroy(&Z);CHKERRQ(ierr);
174*c4762a1bSJed Brown 
175*c4762a1bSJed Brown   ierr = VecDestroy(&X1);CHKERRQ(ierr);
176*c4762a1bSJed Brown   ierr = VecDestroy(&X2);CHKERRQ(ierr);
177*c4762a1bSJed Brown   ierr = VecDestroy(&Z1);CHKERRQ(ierr);
178*c4762a1bSJed Brown   ierr = VecDestroy(&Z2);CHKERRQ(ierr);
179*c4762a1bSJed Brown   ierr = PetscRandomDestroy(&rnd);CHKERRQ(ierr);
180*c4762a1bSJed Brown   ierr = PetscFinalize();
181*c4762a1bSJed Brown   return ierr;
182*c4762a1bSJed Brown }
183*c4762a1bSJed Brown 
184*c4762a1bSJed Brown 
185*c4762a1bSJed Brown /*TEST
186*c4762a1bSJed Brown 
187*c4762a1bSJed Brown    test:
188*c4762a1bSJed Brown 
189*c4762a1bSJed Brown TEST*/
190