xref: /petsc/src/mat/tutorials/ex2.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown static char help[] = "testing SeqDense matrices with an LDA (leading dimension of the user-allocated arrray) larger than M.\n";
2c4762a1bSJed Brown /*
3c4762a1bSJed Brown  * Example code testing SeqDense matrices with an LDA (leading dimension
4c4762a1bSJed Brown  * of the user-allocated arrray) 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*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,0,help));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&rnd));
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   /*
24c4762a1bSJed Brown    * Create matrix and three vectors: these are all normal
25c4762a1bSJed Brown    */
265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(size*size,&a));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(size*size,&b));
28c4762a1bSJed Brown   for (i=0; i<size; i++) {
29c4762a1bSJed Brown     for (j=0; j<size; j++) {
305f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValue(rnd,&a[i+j*size]));
31c4762a1bSJed Brown       b[i+j*size] = a[i+j*size];
32c4762a1bSJed Brown     }
33c4762a1bSJed Brown   }
345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&A));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,size,size,size,size));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATSEQDENSE));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqDenseSetPreallocation(A,a));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
40c4762a1bSJed Brown 
415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(size,&x));
42c4762a1bSJed Brown   for (i=0; i<size; i++) {
43c4762a1bSJed Brown     x[i] = one;
44c4762a1bSJed Brown   }
455f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size,x,&X));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(X));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(X));
48c4762a1bSJed Brown 
495f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(size,&y));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size,y,&Y));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(Y));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(Y));
53c4762a1bSJed Brown 
545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(size,&z));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size,z,&Z));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(Z));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(Z));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   /*
60c4762a1bSJed Brown    * Now create submatrices and subvectors
61c4762a1bSJed Brown    */
625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&A11));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A11,size1,size1,size1,size1));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A11,MATSEQDENSE));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqDenseSetPreallocation(A11,b));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseSetLDA(A11,size));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A11,MAT_FINAL_ASSEMBLY));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A11,MAT_FINAL_ASSEMBLY));
69c4762a1bSJed Brown 
705f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&A12));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A12,size1,size2,size1,size2));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A12,MATSEQDENSE));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqDenseSetPreallocation(A12,b+size1*size));
745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseSetLDA(A12,size));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A12,MAT_FINAL_ASSEMBLY));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A12,MAT_FINAL_ASSEMBLY));
77c4762a1bSJed Brown 
785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&A21));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A21,size2,size1,size2,size1));
805f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A21,MATSEQDENSE));
815f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqDenseSetPreallocation(A21,b+size1));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseSetLDA(A21,size));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A21,MAT_FINAL_ASSEMBLY));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A21,MAT_FINAL_ASSEMBLY));
85c4762a1bSJed Brown 
865f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&A22));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A22,size2,size2,size2,size2));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A22,MATSEQDENSE));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqDenseSetPreallocation(A22,b+size1*size+size1));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseSetLDA(A22,size));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A22,MAT_FINAL_ASSEMBLY));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A22,MAT_FINAL_ASSEMBLY));
93c4762a1bSJed Brown 
945f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size1,x,&X1));
955f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size2,x+size1,&X2));
965f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size1,z,&Z1));
975f80ce2aSJacob Faibussowitsch   CHKERRQ(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    */
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(A,X,Y));
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(A11,X1,Z1));
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultAdd(A12,X2,Z1,Z1));
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(A22,X2,Z2));
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultAdd(A21,X1,Z2,Z2));
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(Z,-1.0,Y));
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(Z,NORM_2,&nrm));
110c4762a1bSJed Brown   if (nrm > 100.0*PETSC_MACHINE_EPSILON) {
1115f80ce2aSJacob Faibussowitsch     CHKERRQ(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    */
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomGetValue(rnd,&v));
118c4762a1bSJed Brown   i    = 1;
119c4762a1bSJed Brown   j    = size-2;
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
121c4762a1bSJed Brown   j   -= size1;
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A12,1,&i,1,&j,&v,INSERT_VALUES));
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomGetValue(rnd,&v));
124c4762a1bSJed Brown   i    = j = size1+1;
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
126c4762a1bSJed Brown   i     =j = 1;
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A22,1,&i,1,&j,&v,INSERT_VALUES));
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A12,MAT_FINAL_ASSEMBLY));
1315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A12,MAT_FINAL_ASSEMBLY));
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A22,MAT_FINAL_ASSEMBLY));
1335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A22,MAT_FINAL_ASSEMBLY));
134c4762a1bSJed Brown 
1355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(A,X,Y));
1365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(A11,X1,Z1));
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultAdd(A12,X2,Z1,Z1));
1385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(A22,X2,Z2));
1395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultAdd(A21,X1,Z2,Z2));
1405f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(Z,-1.0,Y));
1415f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(Z,NORM_2,&nrm));
142c4762a1bSJed Brown   if (nrm > 100.0*PETSC_MACHINE_EPSILON) {
1435f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Test2; error norm=%g\n",(double)nrm));
144c4762a1bSJed Brown   }
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   /*
147c4762a1bSJed Brown    * Transpose product
148c4762a1bSJed Brown    */
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTranspose(A,X,Y));
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTranspose(A11,X1,Z1));
1515f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTransposeAdd(A21,X2,Z1,Z1));
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTranspose(A22,X2,Z2));
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTransposeAdd(A12,X1,Z2,Z2));
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(Z,-1.0,Y));
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(Z,NORM_2,&nrm));
156c4762a1bSJed Brown   if (nrm > 100.0*PETSC_MACHINE_EPSILON) {
1575f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Test3; error norm=%g\n",(double)nrm));
158c4762a1bSJed Brown   }
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(a));
1605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(b));
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(x));
1625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(y));
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(z));
1645f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
1655f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A11));
1665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A12));
1675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A21));
1685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A22));
169c4762a1bSJed Brown 
1705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&X));
1715f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Y));
1725f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Z));
173c4762a1bSJed Brown 
1745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&X1));
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&X2));
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Z1));
1775f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Z2));
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rnd));
179*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
180*b122ec5aSJacob Faibussowitsch   return 0;
181c4762a1bSJed Brown }
182c4762a1bSJed Brown 
183c4762a1bSJed Brown /*TEST
184c4762a1bSJed Brown 
185c4762a1bSJed Brown    test:
186c4762a1bSJed Brown 
187c4762a1bSJed Brown TEST*/
188