1c4762a1bSJed Brown static char help[] = "Test MatMatMult() and MatPtAP() for AIJ matrices.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown extern PetscErrorCode testPTAPRectangular(void); 6c4762a1bSJed Brown 7c4762a1bSJed Brown int main(int argc,char **argv) 8c4762a1bSJed Brown { 9c4762a1bSJed Brown Mat A,B,C,D; 10c4762a1bSJed Brown PetscScalar a[] ={1.,1.,0.,0.,1.,1.,0.,0.,1.}; 11c4762a1bSJed Brown PetscInt ij[]={0,1,2}; 12c4762a1bSJed Brown PetscErrorCode ierr; 13c20d7725SJed Brown PetscReal fill=4.0; 14c4762a1bSJed Brown PetscMPIInt size,rank; 15c20d7725SJed Brown PetscBool isequal; 16c20d7725SJed Brown #if defined(PETSC_HAVE_HYPRE) 17c4762a1bSJed Brown PetscBool test_hypre=PETSC_FALSE; 18c20d7725SJed Brown #endif 19c4762a1bSJed Brown 20c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 21c4762a1bSJed Brown #if defined(PETSC_HAVE_HYPRE) 22c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-test_hypre",&test_hypre,NULL);CHKERRQ(ierr); 23c4762a1bSJed Brown #endif 24ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 25ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 26c4762a1bSJed Brown 27c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 28c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,3,3);CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr); 30c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 32c4762a1bSJed Brown ierr = MatSetOption(A,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 33c4762a1bSJed Brown 34dd400576SPatrick Sanan if (rank == 0) { 35c4762a1bSJed Brown ierr = MatSetValues(A,3,ij,3,ij,a,ADD_VALUES);CHKERRQ(ierr); 36c4762a1bSJed Brown } 37c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* Test MatMatMult() */ 41c4762a1bSJed Brown ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); /* B = A^T */ 42c4762a1bSJed Brown ierr = MatMatMult(B,A,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr); /* C = B*A */ 43c4762a1bSJed Brown ierr = MatMatMult(B,A,MAT_REUSE_MATRIX,fill,&C);CHKERRQ(ierr); /* recompute C=B*A */ 44c4762a1bSJed Brown ierr = MatSetOptionsPrefix(C,"C_");CHKERRQ(ierr); 45c20d7725SJed Brown ierr = MatMatMultEqual(B,A,C,10,&isequal);CHKERRQ(ierr); 46c20d7725SJed Brown if (!isequal) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MatMatMult: C != B*A"); 47c4762a1bSJed Brown 48c20d7725SJed Brown ierr = MatMatMult(C,A,MAT_INITIAL_MATRIX,fill,&D);CHKERRQ(ierr); /* D = C*A = (A^T*A)*A */ 49c20d7725SJed Brown ierr = MatMatMult(C,A,MAT_REUSE_MATRIX,fill,&D);CHKERRQ(ierr); 50c20d7725SJed Brown ierr = MatMatMultEqual(C,A,D,10,&isequal);CHKERRQ(ierr); 51c20d7725SJed Brown if (!isequal) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MatMatMult: D != C*A"); 52c4762a1bSJed Brown 53c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 54c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 55c20d7725SJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 56c4762a1bSJed Brown 57c20d7725SJed Brown /* Test MatPtAP */ 58c4762a1bSJed Brown ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); /* B = A */ 59c4762a1bSJed Brown ierr = MatPtAP(A,B,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr); /* C = B^T*A*B */ 60c20d7725SJed Brown ierr = MatPtAPMultEqual(A,B,C,10,&isequal);CHKERRQ(ierr); 61c20d7725SJed Brown if (!isequal) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MatPtAP: C != B^T*A*B"); 62c4762a1bSJed Brown 63c20d7725SJed Brown /* Repeat MatPtAP to test symbolic/numeric separation for reuse of the symbolic product */ 64c20d7725SJed Brown ierr = MatPtAP(A,B,MAT_REUSE_MATRIX,fill,&C);CHKERRQ(ierr); 65c20d7725SJed Brown ierr = MatPtAPMultEqual(A,B,C,10,&isequal);CHKERRQ(ierr); 66c20d7725SJed Brown if (!isequal) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MatPtAP(reuse): C != B^T*A*B"); 67c4762a1bSJed Brown 68c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 69c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 70c4762a1bSJed Brown 71c4762a1bSJed Brown if (size == 1) { 72c4762a1bSJed Brown /* A test contributed by Tobias Neckel <neckel@in.tum.de> */ 73c4762a1bSJed Brown ierr = testPTAPRectangular();CHKERRQ(ierr); 74c4762a1bSJed Brown 75c4762a1bSJed Brown /* test MatMatTransposeMult(): A*B^T */ 76c4762a1bSJed Brown ierr = MatMatTransposeMult(A,A,MAT_INITIAL_MATRIX,fill,&D);CHKERRQ(ierr); /* D = A*A^T */ 77c20d7725SJed Brown ierr = MatScale(A,2.0);CHKERRQ(ierr); 78c20d7725SJed Brown ierr = MatMatTransposeMult(A,A,MAT_REUSE_MATRIX,fill,&D);CHKERRQ(ierr); 79c20d7725SJed Brown ierr = MatMatTransposeMultEqual(A,A,D,10,&isequal);CHKERRQ(ierr); 80c20d7725SJed Brown if (!isequal) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MatMatTranspose: D != A*A^T"); 81c4762a1bSJed Brown } 82c4762a1bSJed Brown 83c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 84c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 85c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 87c4762a1bSJed Brown ierr = PetscFinalize(); 88c4762a1bSJed Brown return ierr; 89c4762a1bSJed Brown } 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* a test contributed by Tobias Neckel <neckel@in.tum.de>, 02 Jul 2008 */ 92c4762a1bSJed Brown PetscErrorCode testPTAPRectangular(void) 93c4762a1bSJed Brown { 94c20d7725SJed Brown const int rows = 3,cols = 5; 95c4762a1bSJed Brown PetscErrorCode ierr; 96c4762a1bSJed Brown int i; 97c20d7725SJed Brown Mat A,P,C; 98c4762a1bSJed Brown 99c4762a1bSJed Brown PetscFunctionBegin; 100c4762a1bSJed Brown /* set up A */ 101c4762a1bSJed Brown ierr = MatCreateSeqAIJ(PETSC_COMM_WORLD, rows, rows,1, NULL, &A);CHKERRQ(ierr); 102c4762a1bSJed Brown for (i=0; i<rows; i++) { 103c4762a1bSJed Brown ierr = MatSetValue(A, i, i, 1.0, INSERT_VALUES);CHKERRQ(ierr); 104c4762a1bSJed Brown } 105c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 106c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* set up P */ 109c4762a1bSJed Brown ierr = MatCreateSeqAIJ(PETSC_COMM_WORLD, rows, cols,5, NULL, &P);CHKERRQ(ierr); 110c4762a1bSJed Brown ierr = MatSetValue(P, 0, 0, 1.0, INSERT_VALUES);CHKERRQ(ierr); 111c4762a1bSJed Brown ierr = MatSetValue(P, 0, 1, 2.0, INSERT_VALUES);CHKERRQ(ierr); 112c4762a1bSJed Brown ierr = MatSetValue(P, 0, 2, 0.0, INSERT_VALUES);CHKERRQ(ierr); 113c4762a1bSJed Brown 114c4762a1bSJed Brown ierr = MatSetValue(P, 0, 3, -1.0, INSERT_VALUES);CHKERRQ(ierr); 115c4762a1bSJed Brown 116c4762a1bSJed Brown ierr = MatSetValue(P, 1, 0, 0.0, INSERT_VALUES);CHKERRQ(ierr); 117c4762a1bSJed Brown ierr = MatSetValue(P, 1, 1, -1.0, INSERT_VALUES);CHKERRQ(ierr); 118c4762a1bSJed Brown ierr = MatSetValue(P, 1, 2, 1.0, INSERT_VALUES);CHKERRQ(ierr); 119c4762a1bSJed Brown 120c4762a1bSJed Brown ierr = MatSetValue(P, 2, 0, 3.0, INSERT_VALUES);CHKERRQ(ierr); 121c4762a1bSJed Brown ierr = MatSetValue(P, 2, 1, 0.0, INSERT_VALUES);CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = MatSetValue(P, 2, 2, -3.0, INSERT_VALUES);CHKERRQ(ierr); 123c4762a1bSJed Brown 124c4762a1bSJed Brown ierr = MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 125c4762a1bSJed Brown ierr = MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 126c4762a1bSJed Brown 127c4762a1bSJed Brown /* compute C */ 128c4762a1bSJed Brown ierr = MatPtAP(A, P, MAT_INITIAL_MATRIX, 1.0, &C);CHKERRQ(ierr); 129c4762a1bSJed Brown 130c4762a1bSJed Brown ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 131c4762a1bSJed Brown ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 132c4762a1bSJed Brown 133c4762a1bSJed Brown /* compare results */ 134c4762a1bSJed Brown /* 135c4762a1bSJed Brown printf("C:\n"); 136c4762a1bSJed Brown ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 137c4762a1bSJed Brown 138c4762a1bSJed Brown blitz::Array<double,2> actualC(cols, cols); 139c4762a1bSJed Brown actualC = 0.0; 140c4762a1bSJed Brown for (int i=0; i<cols; i++) { 141c4762a1bSJed Brown for (int j=0; j<cols; j++) { 142*2f7452b8SBarry Smith ierr = MatGetValues(C, 1, &i, 1, &j, &actualC(i,j));CHKERRQ(ierr); 143c4762a1bSJed Brown } 144c4762a1bSJed Brown } 145c4762a1bSJed Brown blitz::Array<double,2> expectedC(cols, cols); 146c4762a1bSJed Brown expectedC = 0.0; 147c4762a1bSJed Brown 148c4762a1bSJed Brown expectedC(0,0) = 10.0; 149c4762a1bSJed Brown expectedC(0,1) = 2.0; 150c4762a1bSJed Brown expectedC(0,2) = -9.0; 151c4762a1bSJed Brown expectedC(0,3) = -1.0; 152c4762a1bSJed Brown expectedC(1,0) = 2.0; 153c4762a1bSJed Brown expectedC(1,1) = 5.0; 154c4762a1bSJed Brown expectedC(1,2) = -1.0; 155c4762a1bSJed Brown expectedC(1,3) = -2.0; 156c4762a1bSJed Brown expectedC(2,0) = -9.0; 157c4762a1bSJed Brown expectedC(2,1) = -1.0; 158c4762a1bSJed Brown expectedC(2,2) = 10.0; 159c4762a1bSJed Brown expectedC(2,3) = 0.0; 160c4762a1bSJed Brown expectedC(3,0) = -1.0; 161c4762a1bSJed Brown expectedC(3,1) = -2.0; 162c4762a1bSJed Brown expectedC(3,2) = 0.0; 163c4762a1bSJed Brown expectedC(3,3) = 1.0; 164c4762a1bSJed Brown 165c4762a1bSJed Brown int check = areBlitzArrays2NumericallyEqual(actualC,expectedC); 166c4762a1bSJed Brown validateEqualsWithParams3(check, -1 , "testPTAPRectangular()", check, actualC(check), expectedC(check)); 167c4762a1bSJed Brown */ 168c4762a1bSJed Brown 169c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 170c4762a1bSJed Brown ierr = MatDestroy(&P);CHKERRQ(ierr); 171c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 172c4762a1bSJed Brown PetscFunctionReturn(0); 173c4762a1bSJed Brown } 174c4762a1bSJed Brown 175c4762a1bSJed Brown /*TEST 176c4762a1bSJed Brown 177c4762a1bSJed Brown test: 178c4762a1bSJed Brown 179c4762a1bSJed Brown test: 180c4762a1bSJed Brown suffix: 2 181c4762a1bSJed Brown nsize: 2 182c20d7725SJed Brown args: -matmatmult_via nonscalable 183c20d7725SJed Brown output_file: output/ex93_1.out 184c4762a1bSJed Brown 185c4762a1bSJed Brown test: 186c4762a1bSJed Brown suffix: 3 187c4762a1bSJed Brown nsize: 2 188c20d7725SJed Brown output_file: output/ex93_1.out 189c4762a1bSJed Brown 190c4762a1bSJed Brown test: 191c4762a1bSJed Brown suffix: 4 192c4762a1bSJed Brown nsize: 2 193c20d7725SJed Brown args: -matptap_via scalable 194c20d7725SJed Brown output_file: output/ex93_1.out 195c4762a1bSJed Brown 196c4762a1bSJed Brown test: 197c4762a1bSJed Brown suffix: btheap 198c20d7725SJed Brown args: -matmatmult_via btheap -matmattransmult_via color 199c4762a1bSJed Brown output_file: output/ex93_1.out 200c4762a1bSJed Brown 201c4762a1bSJed Brown test: 202c4762a1bSJed Brown suffix: heap 203c20d7725SJed Brown args: -matmatmult_via heap 204c4762a1bSJed Brown output_file: output/ex93_1.out 205c4762a1bSJed Brown 206c4762a1bSJed Brown #HYPRE PtAP is broken for complex numbers 207c4762a1bSJed Brown test: 208c4762a1bSJed Brown suffix: hypre 209c4762a1bSJed Brown nsize: 3 210c4762a1bSJed Brown requires: hypre !complex 211c20d7725SJed Brown args: -matmatmult_via hypre -matptap_via hypre -test_hypre 212c20d7725SJed Brown output_file: output/ex93_hypre.out 213c4762a1bSJed Brown 214c4762a1bSJed Brown test: 215c4762a1bSJed Brown suffix: llcondensed 216c20d7725SJed Brown args: -matmatmult_via llcondensed 217c4762a1bSJed Brown output_file: output/ex93_1.out 218c4762a1bSJed Brown 219c4762a1bSJed Brown test: 220c4762a1bSJed Brown suffix: scalable 221c20d7725SJed Brown args: -matmatmult_via scalable 222c4762a1bSJed Brown output_file: output/ex93_1.out 223c4762a1bSJed Brown 224c4762a1bSJed Brown test: 225c4762a1bSJed Brown suffix: scalable_fast 226c20d7725SJed Brown args: -matmatmult_via scalable_fast 227c4762a1bSJed Brown output_file: output/ex93_1.out 228c4762a1bSJed Brown 229c4762a1bSJed Brown TEST*/ 230