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