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 99371c9d4SSatish Balay int main(int argc, char **args) { 10c4762a1bSJed Brown Mat A, C, AtA, B; 11c4762a1bSJed Brown PetscViewer fd; 12c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; 13c4762a1bSJed Brown PetscReal fill = 4.0; 14c4762a1bSJed Brown PetscMPIInt rank, size; 15c4762a1bSJed Brown PetscBool equal; 16c4762a1bSJed Brown PetscInt i, m, n, rstart, rend; 17c4762a1bSJed Brown PetscScalar one = 1.0; 18c4762a1bSJed Brown 19327415f7SBarry Smith PetscFunctionBeginUser; 209566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 239566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), NULL)); 24c4762a1bSJed Brown 25c4762a1bSJed Brown /* Load matrix A */ 269566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd)); 279566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 289566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATAIJ)); 299566063dSJacob Faibussowitsch PetscCall(MatLoad(A, fd)); 309566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 31c4762a1bSJed Brown 32c4762a1bSJed Brown /* Create identity matrix B */ 339566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, &n)); 349566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &B)); 359566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, m, m, PETSC_DECIDE, PETSC_DECIDE)); 369566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATAIJ)); 379566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 389566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 39c4762a1bSJed Brown 409566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(B, &rstart, &rend)); 41*48a46eb9SPierre Jolivet for (i = rstart; i < rend; i++) PetscCall(MatSetValues(B, 1, &i, 1, &i, &one, INSERT_VALUES)); 429566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 439566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* Compute AtA = A^T*B*A, B = identity matrix */ 469566063dSJacob Faibussowitsch PetscCall(MatPtAP(B, A, MAT_INITIAL_MATRIX, fill, &AtA)); 479566063dSJacob Faibussowitsch PetscCall(MatPtAP(B, A, MAT_REUSE_MATRIX, fill, &AtA)); 48dd400576SPatrick Sanan if (rank == 0) printf("C = A^T*B*A is done...\n"); 499566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* Compute C = A^T*A */ 529566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(A, A, MAT_INITIAL_MATRIX, fill, &C)); 53dd400576SPatrick Sanan if (rank == 0) printf("C = A^T*A is done...\n"); 549566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(A, A, MAT_REUSE_MATRIX, fill, &C)); 55dd400576SPatrick Sanan if (rank == 0) printf("REUSE C = A^T*A is done...\n"); 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* Compare C and AtA */ 589566063dSJacob Faibussowitsch PetscCall(MatMultEqual(C, AtA, 20, &equal)); 5928b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "A^T*A != At*A"); 609566063dSJacob Faibussowitsch PetscCall(MatDestroy(&AtA)); 61c4762a1bSJed Brown 629566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 639566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 649566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 65b122ec5aSJacob Faibussowitsch return 0; 66c4762a1bSJed Brown } 67c4762a1bSJed Brown 68c4762a1bSJed Brown /*TEST 69c4762a1bSJed Brown 70c4762a1bSJed Brown test: 71dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 72c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/arco1 73c4762a1bSJed Brown 74c4762a1bSJed Brown test: 75c4762a1bSJed Brown suffix: 2 76c4762a1bSJed Brown nsize: 4 77dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 78c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/arco1 -matptap_via nonscalable 79c4762a1bSJed Brown 80c4762a1bSJed Brown TEST*/ 81