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