xref: /petsc/src/mat/tutorials/ex2.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
16aad120cSJose E. Roman static char help[] = "testing SeqDense matrices with an LDA (leading dimension of the user-allocated array) larger than M.\n";
2c4762a1bSJed Brown /*
3c4762a1bSJed Brown  * Example code testing SeqDense matrices with an LDA (leading dimension
46aad120cSJose E. Roman  * of the user-allocated array) larger than M.
5c4762a1bSJed Brown  * This example tests the functionality of MatSetValues(), MatMult(),
6c4762a1bSJed Brown  * and MatMultTranspose().
7c4762a1bSJed Brown  */
8c4762a1bSJed Brown 
9c4762a1bSJed Brown #include <petscmat.h>
10c4762a1bSJed Brown 
11c4762a1bSJed Brown int main(int argc,char **argv)
12c4762a1bSJed Brown {
13c4762a1bSJed Brown   Mat            A,A11,A12,A21,A22;
14c4762a1bSJed Brown   Vec            X,X1,X2,Y,Z,Z1,Z2;
15c4762a1bSJed Brown   PetscScalar    *a,*b,*x,*y,*z,v,one=1;
16c4762a1bSJed Brown   PetscReal      nrm;
17c4762a1bSJed Brown   PetscInt       size=8,size1=6,size2=2, i,j;
18c4762a1bSJed Brown   PetscRandom    rnd;
19c4762a1bSJed Brown 
20*327415f7SBarry Smith   PetscFunctionBeginUser;
219566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,0,help));
229566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rnd));
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   /*
25c4762a1bSJed Brown    * Create matrix and three vectors: these are all normal
26c4762a1bSJed Brown    */
279566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size*size,&a));
289566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size*size,&b));
29c4762a1bSJed Brown   for (i=0; i<size; i++) {
30c4762a1bSJed Brown     for (j=0; j<size; j++) {
319566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rnd,&a[i+j*size]));
32c4762a1bSJed Brown       b[i+j*size] = a[i+j*size];
33c4762a1bSJed Brown     }
34c4762a1bSJed Brown   }
359566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&A));
369566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,size,size,size,size));
379566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATSEQDENSE));
389566063dSJacob Faibussowitsch   PetscCall(MatSeqDenseSetPreallocation(A,a));
399566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
409566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
41c4762a1bSJed Brown 
429566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size,&x));
43c4762a1bSJed Brown   for (i=0; i<size; i++) {
44c4762a1bSJed Brown     x[i] = one;
45c4762a1bSJed Brown   }
469566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size,x,&X));
479566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(X));
489566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(X));
49c4762a1bSJed Brown 
509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size,&y));
519566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size,y,&Y));
529566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(Y));
539566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(Y));
54c4762a1bSJed Brown 
559566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size,&z));
569566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size,z,&Z));
579566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(Z));
589566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(Z));
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /*
61c4762a1bSJed Brown    * Now create submatrices and subvectors
62c4762a1bSJed Brown    */
639566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&A11));
649566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A11,size1,size1,size1,size1));
659566063dSJacob Faibussowitsch   PetscCall(MatSetType(A11,MATSEQDENSE));
669566063dSJacob Faibussowitsch   PetscCall(MatSeqDenseSetPreallocation(A11,b));
679566063dSJacob Faibussowitsch   PetscCall(MatDenseSetLDA(A11,size));
689566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A11,MAT_FINAL_ASSEMBLY));
699566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A11,MAT_FINAL_ASSEMBLY));
70c4762a1bSJed Brown 
719566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&A12));
729566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A12,size1,size2,size1,size2));
739566063dSJacob Faibussowitsch   PetscCall(MatSetType(A12,MATSEQDENSE));
749566063dSJacob Faibussowitsch   PetscCall(MatSeqDenseSetPreallocation(A12,b+size1*size));
759566063dSJacob Faibussowitsch   PetscCall(MatDenseSetLDA(A12,size));
769566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A12,MAT_FINAL_ASSEMBLY));
779566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A12,MAT_FINAL_ASSEMBLY));
78c4762a1bSJed Brown 
799566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&A21));
809566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A21,size2,size1,size2,size1));
819566063dSJacob Faibussowitsch   PetscCall(MatSetType(A21,MATSEQDENSE));
829566063dSJacob Faibussowitsch   PetscCall(MatSeqDenseSetPreallocation(A21,b+size1));
839566063dSJacob Faibussowitsch   PetscCall(MatDenseSetLDA(A21,size));
849566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A21,MAT_FINAL_ASSEMBLY));
859566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A21,MAT_FINAL_ASSEMBLY));
86c4762a1bSJed Brown 
879566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&A22));
889566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A22,size2,size2,size2,size2));
899566063dSJacob Faibussowitsch   PetscCall(MatSetType(A22,MATSEQDENSE));
909566063dSJacob Faibussowitsch   PetscCall(MatSeqDenseSetPreallocation(A22,b+size1*size+size1));
919566063dSJacob Faibussowitsch   PetscCall(MatDenseSetLDA(A22,size));
929566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A22,MAT_FINAL_ASSEMBLY));
939566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A22,MAT_FINAL_ASSEMBLY));
94c4762a1bSJed Brown 
959566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size1,x,&X1));
969566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size2,x+size1,&X2));
979566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size1,z,&Z1));
989566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size2,z+size1,&Z2));
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /*
101c4762a1bSJed Brown    * Now multiple matrix times input in two ways;
102c4762a1bSJed Brown    * compare the results
103c4762a1bSJed Brown    */
1049566063dSJacob Faibussowitsch   PetscCall(MatMult(A,X,Y));
1059566063dSJacob Faibussowitsch   PetscCall(MatMult(A11,X1,Z1));
1069566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(A12,X2,Z1,Z1));
1079566063dSJacob Faibussowitsch   PetscCall(MatMult(A22,X2,Z2));
1089566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(A21,X1,Z2,Z2));
1099566063dSJacob Faibussowitsch   PetscCall(VecAXPY(Z,-1.0,Y));
1109566063dSJacob Faibussowitsch   PetscCall(VecNorm(Z,NORM_2,&nrm));
111c4762a1bSJed Brown   if (nrm > 100.0*PETSC_MACHINE_EPSILON) {
1129566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test1; error norm=%g\n",(double)nrm));
113c4762a1bSJed Brown   }
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /*
116c4762a1bSJed Brown    * Next test: change both matrices
117c4762a1bSJed Brown    */
1189566063dSJacob Faibussowitsch   PetscCall(PetscRandomGetValue(rnd,&v));
119c4762a1bSJed Brown   i    = 1;
120c4762a1bSJed Brown   j    = size-2;
1219566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
122c4762a1bSJed Brown   j   -= size1;
1239566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A12,1,&i,1,&j,&v,INSERT_VALUES));
1249566063dSJacob Faibussowitsch   PetscCall(PetscRandomGetValue(rnd,&v));
125c4762a1bSJed Brown   i    = j = size1+1;
1269566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
127c4762a1bSJed Brown   i     =j = 1;
1289566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A22,1,&i,1,&j,&v,INSERT_VALUES));
1299566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
1309566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
1319566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A12,MAT_FINAL_ASSEMBLY));
1329566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A12,MAT_FINAL_ASSEMBLY));
1339566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A22,MAT_FINAL_ASSEMBLY));
1349566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A22,MAT_FINAL_ASSEMBLY));
135c4762a1bSJed Brown 
1369566063dSJacob Faibussowitsch   PetscCall(MatMult(A,X,Y));
1379566063dSJacob Faibussowitsch   PetscCall(MatMult(A11,X1,Z1));
1389566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(A12,X2,Z1,Z1));
1399566063dSJacob Faibussowitsch   PetscCall(MatMult(A22,X2,Z2));
1409566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(A21,X1,Z2,Z2));
1419566063dSJacob Faibussowitsch   PetscCall(VecAXPY(Z,-1.0,Y));
1429566063dSJacob Faibussowitsch   PetscCall(VecNorm(Z,NORM_2,&nrm));
143c4762a1bSJed Brown   if (nrm > 100.0*PETSC_MACHINE_EPSILON) {
1449566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test2; error norm=%g\n",(double)nrm));
145c4762a1bSJed Brown   }
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   /*
148c4762a1bSJed Brown    * Transpose product
149c4762a1bSJed Brown    */
1509566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(A,X,Y));
1519566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(A11,X1,Z1));
1529566063dSJacob Faibussowitsch   PetscCall(MatMultTransposeAdd(A21,X2,Z1,Z1));
1539566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(A22,X2,Z2));
1549566063dSJacob Faibussowitsch   PetscCall(MatMultTransposeAdd(A12,X1,Z2,Z2));
1559566063dSJacob Faibussowitsch   PetscCall(VecAXPY(Z,-1.0,Y));
1569566063dSJacob Faibussowitsch   PetscCall(VecNorm(Z,NORM_2,&nrm));
157c4762a1bSJed Brown   if (nrm > 100.0*PETSC_MACHINE_EPSILON) {
1589566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test3; error norm=%g\n",(double)nrm));
159c4762a1bSJed Brown   }
1609566063dSJacob Faibussowitsch   PetscCall(PetscFree(a));
1619566063dSJacob Faibussowitsch   PetscCall(PetscFree(b));
1629566063dSJacob Faibussowitsch   PetscCall(PetscFree(x));
1639566063dSJacob Faibussowitsch   PetscCall(PetscFree(y));
1649566063dSJacob Faibussowitsch   PetscCall(PetscFree(z));
1659566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1669566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A11));
1679566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A12));
1689566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A21));
1699566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A22));
170c4762a1bSJed Brown 
1719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
1729566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Y));
1739566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Z));
174c4762a1bSJed Brown 
1759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X1));
1769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X2));
1779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Z1));
1789566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Z2));
1799566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rnd));
1809566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
181b122ec5aSJacob Faibussowitsch   return 0;
182c4762a1bSJed Brown }
183c4762a1bSJed Brown 
184c4762a1bSJed Brown /*TEST
185c4762a1bSJed Brown 
186c4762a1bSJed Brown    test:
187c4762a1bSJed Brown 
188c4762a1bSJed Brown TEST*/
189