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