xref: /petsc/src/mat/tests/ex163.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
1 
2 static char help[] = "Tests MatTransposeMatMult() on MatLoad() matrix \n\n";
3 
4 #include <petscmat.h>
5 
6 int main(int argc,char **args)
7 {
8   Mat            A,C,Bdense,Cdense;
9   PetscErrorCode ierr;
10   PetscViewer    fd;              /* viewer */
11   char           file[PETSC_MAX_PATH_LEN]; /* input file name */
12   PetscBool      flg,viewmats=PETSC_FALSE;
13   PetscMPIInt    rank,size;
14   PetscReal      fill=1.0;
15   PetscInt       m,n,i,j,BN=10,rstart,rend,*rows,*cols;
16   PetscScalar    *Barray,*Carray,rval,*array;
17   Vec            x,y;
18   PetscRandom    rand;
19 
20   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
21   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
22 
23   /* Determine file from which we read the matrix A */
24   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg));
25   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option");
26 
27   /* Load matrix A */
28   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
29   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
30   CHKERRQ(MatLoad(A,fd));
31   CHKERRQ(PetscViewerDestroy(&fd));
32 
33   /* Print (for testing only) */
34   CHKERRQ(PetscOptionsHasName(NULL,NULL, "-view_mats", &viewmats));
35   if (viewmats) {
36     if (rank == 0) printf("A_aij:\n");
37     CHKERRQ(MatView(A,0));
38   }
39 
40   /* Test MatTransposeMatMult_aij_aij() */
41   CHKERRQ(MatTransposeMatMult(A,A,MAT_INITIAL_MATRIX,fill,&C));
42   if (viewmats) {
43     if (rank == 0) printf("\nC = A_aij^T * A_aij:\n");
44     CHKERRQ(MatView(C,0));
45   }
46   CHKERRQ(MatDestroy(&C));
47   CHKERRQ(MatGetLocalSize(A,&m,&n));
48 
49   /* create a dense matrix Bdense */
50   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&Bdense));
51   CHKERRQ(MatSetSizes(Bdense,m,PETSC_DECIDE,PETSC_DECIDE,BN));
52   CHKERRQ(MatSetType(Bdense,MATDENSE));
53   CHKERRQ(MatSetFromOptions(Bdense));
54   CHKERRQ(MatSetUp(Bdense));
55   CHKERRQ(MatGetOwnershipRange(Bdense,&rstart,&rend));
56 
57   CHKERRQ(PetscMalloc3(m,&rows,BN,&cols,m*BN,&array));
58   for (i=0; i<m; i++) rows[i] = rstart + i;
59   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rand));
60   CHKERRQ(PetscRandomSetFromOptions(rand));
61   for (j=0; j<BN; j++) {
62     cols[j] = j;
63     for (i=0; i<m; i++) {
64       CHKERRQ(PetscRandomGetValue(rand,&rval));
65       array[m*j+i] = rval;
66     }
67   }
68   CHKERRQ(MatSetValues(Bdense,m,rows,BN,cols,array,INSERT_VALUES));
69   CHKERRQ(MatAssemblyBegin(Bdense,MAT_FINAL_ASSEMBLY));
70   CHKERRQ(MatAssemblyEnd(Bdense,MAT_FINAL_ASSEMBLY));
71   CHKERRQ(PetscRandomDestroy(&rand));
72   CHKERRQ(PetscFree3(rows,cols,array));
73   if (viewmats) {
74     if (rank == 0) printf("\nBdense:\n");
75     CHKERRQ(MatView(Bdense,0));
76   }
77 
78   /* Test MatTransposeMatMult_aij_dense() */
79   CHKERRQ(MatTransposeMatMult(A,Bdense,MAT_INITIAL_MATRIX,fill,&C));
80   CHKERRQ(MatTransposeMatMult(A,Bdense,MAT_REUSE_MATRIX,fill,&C));
81   if (viewmats) {
82     if (rank == 0) printf("\nC=A^T*Bdense:\n");
83     CHKERRQ(MatView(C,0));
84   }
85 
86   /* Check accuracy */
87   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&Cdense));
88   CHKERRQ(MatSetSizes(Cdense,n,PETSC_DECIDE,PETSC_DECIDE,BN));
89   CHKERRQ(MatSetType(Cdense,MATDENSE));
90   CHKERRQ(MatSetFromOptions(Cdense));
91   CHKERRQ(MatSetUp(Cdense));
92   CHKERRQ(MatAssemblyBegin(Cdense,MAT_FINAL_ASSEMBLY));
93   CHKERRQ(MatAssemblyEnd(Cdense,MAT_FINAL_ASSEMBLY));
94 
95   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
96   if (size == 1) {
97     CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,m,NULL,&x));
98     CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,n,NULL,&y));
99   } else {
100     CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,m,PETSC_DECIDE,NULL,&x));
101     CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&y));
102   }
103 
104   /* Cdense[:,j] = A^T * Bdense[:,j] */
105   CHKERRQ(MatDenseGetArray(Bdense,&Barray));
106   CHKERRQ(MatDenseGetArray(Cdense,&Carray));
107   for (j=0; j<BN; j++) {
108     CHKERRQ(VecPlaceArray(x,Barray));
109     CHKERRQ(VecPlaceArray(y,Carray));
110 
111     CHKERRQ(MatMultTranspose(A,x,y));
112 
113     CHKERRQ(VecResetArray(x));
114     CHKERRQ(VecResetArray(y));
115     Barray += m;
116     Carray += n;
117   }
118   CHKERRQ(MatDenseRestoreArray(Bdense,&Barray));
119   CHKERRQ(MatDenseRestoreArray(Cdense,&Carray));
120   if (viewmats) {
121     if (rank == 0) printf("\nCdense:\n");
122     CHKERRQ(MatView(Cdense,0));
123   }
124 
125   CHKERRQ(MatEqual(C,Cdense,&flg));
126   if (!flg) {
127     if (rank == 0) printf(" C != Cdense\n");
128   }
129 
130   /* Free data structures */
131   CHKERRQ(MatDestroy(&A));
132   CHKERRQ(MatDestroy(&C));
133   CHKERRQ(MatDestroy(&Bdense));
134   CHKERRQ(MatDestroy(&Cdense));
135   CHKERRQ(VecDestroy(&x));
136   CHKERRQ(VecDestroy(&y));
137   ierr = PetscFinalize();
138   return ierr;
139 }
140 
141 /*TEST
142 
143    test:
144       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
145       args: -f ${DATAFILESPATH}/matrices/small
146       output_file: output/ex163.out
147 
148    test:
149       suffix: 2
150       nsize: 3
151       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
152       args: -f ${DATAFILESPATH}/matrices/small
153       output_file: output/ex163.out
154 
155 TEST*/
156