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