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*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 215f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 225f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),NULL)); 24c4762a1bSJed Brown 25c4762a1bSJed Brown /* Load matrix A */ 265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATAIJ)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,fd)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&fd)); 31c4762a1bSJed Brown 32c4762a1bSJed Brown /* Create identity matrix B */ 335f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A,&m,&n)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(B,m,m,PETSC_DECIDE,PETSC_DECIDE)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(B,MATAIJ)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(B)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(B)); 39c4762a1bSJed Brown 405f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(B,&rstart,&rend)); 41c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 425f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,1,&i,1,&i,&one,INSERT_VALUES)); 43c4762a1bSJed Brown } 445f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 46c4762a1bSJed Brown 47c4762a1bSJed Brown /* Compute AtA = A^T*B*A, B = identity matrix */ 485f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAP(B,A,MAT_INITIAL_MATRIX,fill,&AtA)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAP(B,A,MAT_REUSE_MATRIX,fill,&AtA)); 50dd400576SPatrick Sanan if (rank == 0) printf("C = A^T*B*A is done...\n"); 515f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 52c4762a1bSJed Brown 53c4762a1bSJed Brown /* Compute C = A^T*A */ 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMult(A,A,MAT_INITIAL_MATRIX,fill,&C)); 55dd400576SPatrick Sanan if (rank == 0) printf("C = A^T*A is done...\n"); 565f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 605f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(C,AtA,20,&equal)); 6128b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"A^T*A != At*A"); 625f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&AtA)); 63c4762a1bSJed Brown 645f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 66*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 67*b122ec5aSJacob 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