1*fe47926eSStefano Zampini #include <petscmat.h> 2*fe47926eSStefano Zampini 3*fe47926eSStefano Zampini static char help[] = "Example of MatMat ops with MatDense in PETSc.\n"; 4*fe47926eSStefano Zampini 5*fe47926eSStefano Zampini int main(int argc, char **args) 6*fe47926eSStefano Zampini { 7*fe47926eSStefano Zampini Mat A, P, PtAP, RARt, ABC, Pt; 8*fe47926eSStefano Zampini PetscInt n = 4, m = 2; // Example dimensions 9*fe47926eSStefano Zampini PetscMPIInt size; 10*fe47926eSStefano Zampini 11*fe47926eSStefano Zampini PetscCall(PetscInitialize(&argc, &args, NULL, help)); 12*fe47926eSStefano Zampini PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 13*fe47926eSStefano Zampini PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This example is for sequential runs only."); 14*fe47926eSStefano Zampini 15*fe47926eSStefano Zampini // Create dense matrix P (n x m) 16*fe47926eSStefano Zampini PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, n, m, NULL, &P)); 17*fe47926eSStefano Zampini PetscScalar P_values[] = {1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0}; 18*fe47926eSStefano Zampini PetscInt rows[] = {0, 1, 2, 3}; 19*fe47926eSStefano Zampini PetscInt cols[] = {0, 1}; 20*fe47926eSStefano Zampini PetscCall(MatSetValues(P, n, rows, m, cols, P_values, INSERT_VALUES)); 21*fe47926eSStefano Zampini PetscCall(MatAssemblyBegin(P, MAT_FINAL_ASSEMBLY)); 22*fe47926eSStefano Zampini PetscCall(MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY)); 23*fe47926eSStefano Zampini PetscCall(MatTranspose(P, MAT_INITIAL_MATRIX, &Pt)); 24*fe47926eSStefano Zampini 25*fe47926eSStefano Zampini // Create dense matrix A (n x n) 26*fe47926eSStefano Zampini PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, n, n, NULL, &A)); 27*fe47926eSStefano Zampini PetscCall(MatSetBlockSize(A, n)); // Set block size for A 28*fe47926eSStefano Zampini PetscScalar A_values[] = {4.0, 1.0, 2.0, 0.0, 1.0, 3.0, 0.0, 1.0, 2.0, 0.0, 5.0, 2.0, 0.0, 1.0, 2.0, 6.0}; 29*fe47926eSStefano Zampini PetscInt indices[] = {0, 1, 2, 3}; 30*fe47926eSStefano Zampini PetscCall(MatSetValues(A, n, indices, n, indices, A_values, INSERT_VALUES)); 31*fe47926eSStefano Zampini PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 32*fe47926eSStefano Zampini PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 33*fe47926eSStefano Zampini 34*fe47926eSStefano Zampini PetscCall(MatPtAP(A, P, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &PtAP)); 35*fe47926eSStefano Zampini PetscCall(MatRARt(A, Pt, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &RARt)); 36*fe47926eSStefano Zampini PetscCall(MatMatMatMult(Pt, A, P, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &ABC)); 37*fe47926eSStefano Zampini 38*fe47926eSStefano Zampini // View matrices 39*fe47926eSStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_SELF, "Matrix P:\n")); 40*fe47926eSStefano Zampini PetscCall(MatView(P, PETSC_VIEWER_STDOUT_SELF)); 41*fe47926eSStefano Zampini 42*fe47926eSStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_SELF, "Matrix A:\n")); 43*fe47926eSStefano Zampini PetscCall(MatView(A, PETSC_VIEWER_STDOUT_SELF)); 44*fe47926eSStefano Zampini 45*fe47926eSStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_SELF, "Matrix PtAP:\n")); 46*fe47926eSStefano Zampini PetscCall(MatView(PtAP, PETSC_VIEWER_STDOUT_SELF)); 47*fe47926eSStefano Zampini 48*fe47926eSStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_SELF, "Matrix RARt:\n")); 49*fe47926eSStefano Zampini PetscCall(MatView(RARt, PETSC_VIEWER_STDOUT_SELF)); 50*fe47926eSStefano Zampini 51*fe47926eSStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_SELF, "Matrix ABC:\n")); 52*fe47926eSStefano Zampini PetscCall(MatView(ABC, PETSC_VIEWER_STDOUT_SELF)); 53*fe47926eSStefano Zampini 54*fe47926eSStefano Zampini // Clean up 55*fe47926eSStefano Zampini PetscCall(MatDestroy(&P)); 56*fe47926eSStefano Zampini PetscCall(MatDestroy(&Pt)); 57*fe47926eSStefano Zampini PetscCall(MatDestroy(&A)); 58*fe47926eSStefano Zampini PetscCall(MatDestroy(&PtAP)); 59*fe47926eSStefano Zampini PetscCall(MatDestroy(&RARt)); 60*fe47926eSStefano Zampini PetscCall(MatDestroy(&ABC)); 61*fe47926eSStefano Zampini 62*fe47926eSStefano Zampini PetscCall(PetscFinalize()); 63*fe47926eSStefano Zampini return 0; 64*fe47926eSStefano Zampini } 65*fe47926eSStefano Zampini 66*fe47926eSStefano Zampini /*TEST 67*fe47926eSStefano Zampini 68*fe47926eSStefano Zampini test: 69*fe47926eSStefano Zampini diff_args: -j 70*fe47926eSStefano Zampini suffix: 1 71*fe47926eSStefano Zampini 72*fe47926eSStefano Zampini TEST*/ 73