xref: /petsc/src/mat/tests/ex209.c (revision ffa8c5705e8ab2cf85ee1d14dbe507a6e2eb5283)
1 static char help[] = "Test MatTransposeMatMult() \n\n";
2 
3 /* Example:
4   mpiexec -n 8 ./ex209 -f <inputfile> (e.g., inputfile=ceres_solver_iteration_001_A.petsc)
5 */
6 
7 #include <petscmat.h>
8 
9 int main(int argc,char **args)
10 {
11   Mat            A,C,AtA,B;
12   PetscViewer    fd;
13   char           file[PETSC_MAX_PATH_LEN];
14   PetscReal      fill = 4.0;
15   PetscMPIInt    rank,size;
16   PetscBool      equal;
17   PetscInt       i,m,n,rstart,rend;
18   PetscScalar    one=1.0;
19 
20   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
21   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
22   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
23   PetscCall(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),NULL));
24 
25   /* Load matrix A */
26   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
27   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
28   PetscCall(MatSetType(A,MATAIJ));
29   PetscCall(MatLoad(A,fd));
30   PetscCall(PetscViewerDestroy(&fd));
31 
32   /* Create identity matrix B */
33   PetscCall(MatGetLocalSize(A,&m,&n));
34   PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
35   PetscCall(MatSetSizes(B,m,m,PETSC_DECIDE,PETSC_DECIDE));
36   PetscCall(MatSetType(B,MATAIJ));
37   PetscCall(MatSetFromOptions(B));
38   PetscCall(MatSetUp(B));
39 
40   PetscCall(MatGetOwnershipRange(B,&rstart,&rend));
41   for (i=rstart; i<rend; i++) {
42     PetscCall(MatSetValues(B,1,&i,1,&i,&one,INSERT_VALUES));
43   }
44   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
45   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
46 
47   /* Compute AtA = A^T*B*A, B = identity matrix */
48   PetscCall(MatPtAP(B,A,MAT_INITIAL_MATRIX,fill,&AtA));
49   PetscCall(MatPtAP(B,A,MAT_REUSE_MATRIX,fill,&AtA));
50   if (rank == 0) printf("C = A^T*B*A is done...\n");
51   PetscCall(MatDestroy(&B));
52 
53   /* Compute C = A^T*A */
54   PetscCall(MatTransposeMatMult(A,A,MAT_INITIAL_MATRIX,fill,&C));
55   if (rank == 0) printf("C = A^T*A is done...\n");
56   PetscCall(MatTransposeMatMult(A,A,MAT_REUSE_MATRIX,fill,&C));
57   if (rank == 0) printf("REUSE C = A^T*A is done...\n");
58 
59   /* Compare C and AtA */
60   PetscCall(MatMultEqual(C,AtA,20,&equal));
61   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"A^T*A != At*A");
62   PetscCall(MatDestroy(&AtA));
63 
64   PetscCall(MatDestroy(&C));
65   PetscCall(MatDestroy(&A));
66   PetscCall(PetscFinalize());
67   return 0;
68 }
69 
70 /*TEST
71 
72    test:
73       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
74       args: -f ${DATAFILESPATH}/matrices/arco1
75 
76    test:
77       suffix: 2
78       nsize: 4
79       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
80       args: -f ${DATAFILESPATH}/matrices/arco1 -matptap_via nonscalable
81 
82 TEST*/
83