xref: /petsc/src/mat/tests/ex209.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
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   PetscErrorCode ierr;
15c4762a1bSJed Brown   PetscReal      fill = 4.0;
16c4762a1bSJed Brown   PetscMPIInt    rank,size;
17c4762a1bSJed Brown   PetscBool      equal;
18c4762a1bSJed Brown   PetscInt       i,m,n,rstart,rend;
19c4762a1bSJed Brown   PetscScalar    one=1.0;
20c4762a1bSJed Brown 
21c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
225f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
235f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),NULL));
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   /* Load matrix A */
275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATAIJ));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(A,fd));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   /* Create identity matrix B */
345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetLocalSize(A,&m,&n));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(B,m,m,PETSC_DECIDE,PETSC_DECIDE));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(B,MATAIJ));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(B));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(B));
40c4762a1bSJed Brown 
415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(B,&rstart,&rend));
42c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
435f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(B,1,&i,1,&i,&one,INSERT_VALUES));
44c4762a1bSJed Brown   }
455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
47c4762a1bSJed Brown 
48c4762a1bSJed Brown   /* Compute AtA = A^T*B*A, B = identity matrix */
495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPtAP(B,A,MAT_INITIAL_MATRIX,fill,&AtA));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPtAP(B,A,MAT_REUSE_MATRIX,fill,&AtA));
51dd400576SPatrick Sanan   if (rank == 0) printf("C = A^T*B*A is done...\n");
525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   /* Compute C = A^T*A */
555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTransposeMatMult(A,A,MAT_INITIAL_MATRIX,fill,&C));
56dd400576SPatrick Sanan   if (rank == 0) printf("C = A^T*A is done...\n");
575f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTransposeMatMult(A,A,MAT_REUSE_MATRIX,fill,&C));
58dd400576SPatrick Sanan   if (rank == 0) printf("REUSE C = A^T*A is done...\n");
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /* Compare C and AtA */
615f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultEqual(C,AtA,20,&equal));
62*28b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"A^T*A != At*A");
635f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&AtA));
64c4762a1bSJed Brown 
655f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
67c4762a1bSJed Brown   ierr = PetscFinalize();
68c4762a1bSJed Brown   return ierr;
69c4762a1bSJed Brown }
70c4762a1bSJed Brown 
71c4762a1bSJed Brown /*TEST
72c4762a1bSJed Brown 
73c4762a1bSJed Brown    test:
74dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
75c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/arco1
76c4762a1bSJed Brown 
77c4762a1bSJed Brown    test:
78c4762a1bSJed Brown       suffix: 2
79c4762a1bSJed Brown       nsize: 4
80dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
81c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/arco1 -matptap_via nonscalable
82c4762a1bSJed Brown 
83c4762a1bSJed Brown TEST*/
84