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; 22ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 23ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 24589a23caSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),NULL);CHKERRQ(ierr); 25c4762a1bSJed Brown 26c4762a1bSJed Brown /* Load matrix A */ 27c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); 28c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr); 30c4762a1bSJed Brown ierr = MatLoad(A,fd);CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* Create identity matrix B */ 34c4762a1bSJed Brown ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 35c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr); 36c4762a1bSJed Brown ierr = MatSetSizes(B,m,m,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 37c4762a1bSJed Brown ierr = MatSetType(B,MATAIJ);CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = MatSetFromOptions(B);CHKERRQ(ierr); 39c4762a1bSJed Brown ierr = MatSetUp(B);CHKERRQ(ierr); 40c4762a1bSJed Brown 41c4762a1bSJed Brown ierr = MatGetOwnershipRange(B,&rstart,&rend);CHKERRQ(ierr); 42c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 43c4762a1bSJed Brown ierr = MatSetValues(B,1,&i,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr); 44c4762a1bSJed Brown } 45c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* Compute AtA = A^T*B*A, B = identity matrix */ 49c4762a1bSJed Brown ierr = MatPtAP(B,A,MAT_INITIAL_MATRIX,fill,&AtA);CHKERRQ(ierr); 50c4762a1bSJed Brown ierr = MatPtAP(B,A,MAT_REUSE_MATRIX,fill,&AtA);CHKERRQ(ierr); 51c4762a1bSJed Brown if (!rank) printf("C = A^T*B*A is done...\n"); 52c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* Compute C = A^T*A */ 55c4762a1bSJed Brown ierr = MatTransposeMatMult(A,A,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr); 56c4762a1bSJed Brown if (!rank) printf("C = A^T*A is done...\n"); 57c4762a1bSJed Brown ierr = MatTransposeMatMult(A,A,MAT_REUSE_MATRIX,fill,&C);CHKERRQ(ierr); 58c4762a1bSJed Brown if (!rank) printf("REUSE C = A^T*A is done...\n"); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Compare C and AtA */ 61c4762a1bSJed Brown ierr = MatMultEqual(C,AtA,20,&equal);CHKERRQ(ierr); 62c4762a1bSJed Brown if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"A^T*A != At*A"); 63c4762a1bSJed Brown ierr = MatDestroy(&AtA);CHKERRQ(ierr); 64c4762a1bSJed Brown 65c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 66c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 67c4762a1bSJed Brown ierr = PetscFinalize(); 68c4762a1bSJed Brown return ierr; 69c4762a1bSJed Brown } 70c4762a1bSJed Brown 71c4762a1bSJed Brown /*TEST 72c4762a1bSJed Brown 73c4762a1bSJed Brown test: 74*dfd57a17SPierre 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 80*dfd57a17SPierre 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