xref: /petsc/src/mat/tutorials/ex2.c (revision 6aad120caa16b1027d343a5f30f73d01448e4dc0)
1*6aad120cSJose 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
4*6aad120cSJose 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 
209566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,0,help));
219566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rnd));
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   /*
24c4762a1bSJed Brown    * Create matrix and three vectors: these are all normal
25c4762a1bSJed Brown    */
269566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size*size,&a));
279566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size*size,&b));
28c4762a1bSJed Brown   for (i=0; i<size; i++) {
29c4762a1bSJed Brown     for (j=0; j<size; j++) {
309566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rnd,&a[i+j*size]));
31c4762a1bSJed Brown       b[i+j*size] = a[i+j*size];
32c4762a1bSJed Brown     }
33c4762a1bSJed Brown   }
349566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&A));
359566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,size,size,size,size));
369566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATSEQDENSE));
379566063dSJacob Faibussowitsch   PetscCall(MatSeqDenseSetPreallocation(A,a));
389566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
399566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
40c4762a1bSJed Brown 
419566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size,&x));
42c4762a1bSJed Brown   for (i=0; i<size; i++) {
43c4762a1bSJed Brown     x[i] = one;
44c4762a1bSJed Brown   }
459566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size,x,&X));
469566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(X));
479566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(X));
48c4762a1bSJed Brown 
499566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size,&y));
509566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size,y,&Y));
519566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(Y));
529566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(Y));
53c4762a1bSJed Brown 
549566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size,&z));
559566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size,z,&Z));
569566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(Z));
579566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(Z));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   /*
60c4762a1bSJed Brown    * Now create submatrices and subvectors
61c4762a1bSJed Brown    */
629566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&A11));
639566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A11,size1,size1,size1,size1));
649566063dSJacob Faibussowitsch   PetscCall(MatSetType(A11,MATSEQDENSE));
659566063dSJacob Faibussowitsch   PetscCall(MatSeqDenseSetPreallocation(A11,b));
669566063dSJacob Faibussowitsch   PetscCall(MatDenseSetLDA(A11,size));
679566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A11,MAT_FINAL_ASSEMBLY));
689566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A11,MAT_FINAL_ASSEMBLY));
69c4762a1bSJed Brown 
709566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&A12));
719566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A12,size1,size2,size1,size2));
729566063dSJacob Faibussowitsch   PetscCall(MatSetType(A12,MATSEQDENSE));
739566063dSJacob Faibussowitsch   PetscCall(MatSeqDenseSetPreallocation(A12,b+size1*size));
749566063dSJacob Faibussowitsch   PetscCall(MatDenseSetLDA(A12,size));
759566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A12,MAT_FINAL_ASSEMBLY));
769566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A12,MAT_FINAL_ASSEMBLY));
77c4762a1bSJed Brown 
789566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&A21));
799566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A21,size2,size1,size2,size1));
809566063dSJacob Faibussowitsch   PetscCall(MatSetType(A21,MATSEQDENSE));
819566063dSJacob Faibussowitsch   PetscCall(MatSeqDenseSetPreallocation(A21,b+size1));
829566063dSJacob Faibussowitsch   PetscCall(MatDenseSetLDA(A21,size));
839566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A21,MAT_FINAL_ASSEMBLY));
849566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A21,MAT_FINAL_ASSEMBLY));
85c4762a1bSJed Brown 
869566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&A22));
879566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A22,size2,size2,size2,size2));
889566063dSJacob Faibussowitsch   PetscCall(MatSetType(A22,MATSEQDENSE));
899566063dSJacob Faibussowitsch   PetscCall(MatSeqDenseSetPreallocation(A22,b+size1*size+size1));
909566063dSJacob Faibussowitsch   PetscCall(MatDenseSetLDA(A22,size));
919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A22,MAT_FINAL_ASSEMBLY));
929566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A22,MAT_FINAL_ASSEMBLY));
93c4762a1bSJed Brown 
949566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size1,x,&X1));
959566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size2,x+size1,&X2));
969566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size1,z,&Z1));
979566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size2,z+size1,&Z2));
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   /*
100c4762a1bSJed Brown    * Now multiple matrix times input in two ways;
101c4762a1bSJed Brown    * compare the results
102c4762a1bSJed Brown    */
1039566063dSJacob Faibussowitsch   PetscCall(MatMult(A,X,Y));
1049566063dSJacob Faibussowitsch   PetscCall(MatMult(A11,X1,Z1));
1059566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(A12,X2,Z1,Z1));
1069566063dSJacob Faibussowitsch   PetscCall(MatMult(A22,X2,Z2));
1079566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(A21,X1,Z2,Z2));
1089566063dSJacob Faibussowitsch   PetscCall(VecAXPY(Z,-1.0,Y));
1099566063dSJacob Faibussowitsch   PetscCall(VecNorm(Z,NORM_2,&nrm));
110c4762a1bSJed Brown   if (nrm > 100.0*PETSC_MACHINE_EPSILON) {
1119566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test1; error norm=%g\n",(double)nrm));
112c4762a1bSJed Brown   }
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   /*
115c4762a1bSJed Brown    * Next test: change both matrices
116c4762a1bSJed Brown    */
1179566063dSJacob Faibussowitsch   PetscCall(PetscRandomGetValue(rnd,&v));
118c4762a1bSJed Brown   i    = 1;
119c4762a1bSJed Brown   j    = size-2;
1209566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
121c4762a1bSJed Brown   j   -= size1;
1229566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A12,1,&i,1,&j,&v,INSERT_VALUES));
1239566063dSJacob Faibussowitsch   PetscCall(PetscRandomGetValue(rnd,&v));
124c4762a1bSJed Brown   i    = j = size1+1;
1259566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
126c4762a1bSJed Brown   i     =j = 1;
1279566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A22,1,&i,1,&j,&v,INSERT_VALUES));
1289566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
1299566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
1309566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A12,MAT_FINAL_ASSEMBLY));
1319566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A12,MAT_FINAL_ASSEMBLY));
1329566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A22,MAT_FINAL_ASSEMBLY));
1339566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A22,MAT_FINAL_ASSEMBLY));
134c4762a1bSJed Brown 
1359566063dSJacob Faibussowitsch   PetscCall(MatMult(A,X,Y));
1369566063dSJacob Faibussowitsch   PetscCall(MatMult(A11,X1,Z1));
1379566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(A12,X2,Z1,Z1));
1389566063dSJacob Faibussowitsch   PetscCall(MatMult(A22,X2,Z2));
1399566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(A21,X1,Z2,Z2));
1409566063dSJacob Faibussowitsch   PetscCall(VecAXPY(Z,-1.0,Y));
1419566063dSJacob Faibussowitsch   PetscCall(VecNorm(Z,NORM_2,&nrm));
142c4762a1bSJed Brown   if (nrm > 100.0*PETSC_MACHINE_EPSILON) {
1439566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test2; error norm=%g\n",(double)nrm));
144c4762a1bSJed Brown   }
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   /*
147c4762a1bSJed Brown    * Transpose product
148c4762a1bSJed Brown    */
1499566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(A,X,Y));
1509566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(A11,X1,Z1));
1519566063dSJacob Faibussowitsch   PetscCall(MatMultTransposeAdd(A21,X2,Z1,Z1));
1529566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(A22,X2,Z2));
1539566063dSJacob Faibussowitsch   PetscCall(MatMultTransposeAdd(A12,X1,Z2,Z2));
1549566063dSJacob Faibussowitsch   PetscCall(VecAXPY(Z,-1.0,Y));
1559566063dSJacob Faibussowitsch   PetscCall(VecNorm(Z,NORM_2,&nrm));
156c4762a1bSJed Brown   if (nrm > 100.0*PETSC_MACHINE_EPSILON) {
1579566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test3; error norm=%g\n",(double)nrm));
158c4762a1bSJed Brown   }
1599566063dSJacob Faibussowitsch   PetscCall(PetscFree(a));
1609566063dSJacob Faibussowitsch   PetscCall(PetscFree(b));
1619566063dSJacob Faibussowitsch   PetscCall(PetscFree(x));
1629566063dSJacob Faibussowitsch   PetscCall(PetscFree(y));
1639566063dSJacob Faibussowitsch   PetscCall(PetscFree(z));
1649566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1659566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A11));
1669566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A12));
1679566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A21));
1689566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A22));
169c4762a1bSJed Brown 
1709566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
1719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Y));
1729566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Z));
173c4762a1bSJed Brown 
1749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X1));
1759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X2));
1769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Z1));
1779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Z2));
1789566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rnd));
1799566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
180b122ec5aSJacob Faibussowitsch   return 0;
181c4762a1bSJed Brown }
182c4762a1bSJed Brown 
183c4762a1bSJed Brown /*TEST
184c4762a1bSJed Brown 
185c4762a1bSJed Brown    test:
186c4762a1bSJed Brown 
187c4762a1bSJed Brown TEST*/
188