1*c4762a1bSJed Brown static char help[] = "Test MatMatMult() for AIJ and Dense matrices.\n\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscmat.h> 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown int main(int argc,char **argv) 6*c4762a1bSJed Brown { 7*c4762a1bSJed Brown Mat A,B,C,D; 8*c4762a1bSJed Brown PetscInt i,M,N,Istart,Iend,n=7,j,J,Ii,m=8,am,an; 9*c4762a1bSJed Brown PetscScalar v; 10*c4762a1bSJed Brown PetscErrorCode ierr; 11*c4762a1bSJed Brown PetscRandom r; 12*c4762a1bSJed Brown PetscBool equal=PETSC_FALSE; 13*c4762a1bSJed Brown PetscReal fill = 1.0; 14*c4762a1bSJed Brown PetscMPIInt size; 15*c4762a1bSJed Brown 16*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 17*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr); 18*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); 19*c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL);CHKERRQ(ierr); 20*c4762a1bSJed Brown 21*c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&r);CHKERRQ(ierr); 22*c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr); 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown /* Create a aij matrix A */ 25*c4762a1bSJed Brown M = N = m*n; 26*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 27*c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 28*c4762a1bSJed Brown ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr); 29*c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 30*c4762a1bSJed Brown ierr = MatMPIAIJSetPreallocation(A,5,NULL,5,NULL);CHKERRQ(ierr); 31*c4762a1bSJed Brown ierr = MatSeqAIJSetPreallocation(A,5,NULL);CHKERRQ(ierr); 32*c4762a1bSJed Brown 33*c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&Istart,&Iend);CHKERRQ(ierr); 34*c4762a1bSJed Brown am = Iend - Istart; 35*c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 36*c4762a1bSJed Brown v = -1.0; i = Ii/n; j = Ii - i*n; 37*c4762a1bSJed Brown if (i>0) {J = Ii - n; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 38*c4762a1bSJed Brown if (i<m-1) {J = Ii + n; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 39*c4762a1bSJed Brown if (j>0) {J = Ii - 1; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 40*c4762a1bSJed Brown if (j<n-1) {J = Ii + 1; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 41*c4762a1bSJed Brown v = 4.0; ierr = MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr); 42*c4762a1bSJed Brown } 43*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 44*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 45*c4762a1bSJed Brown 46*c4762a1bSJed Brown /* Create a dense matrix B */ 47*c4762a1bSJed Brown ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr); 48*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr); 49*c4762a1bSJed Brown ierr = MatSetSizes(B,an,PETSC_DECIDE,PETSC_DECIDE,M);CHKERRQ(ierr); 50*c4762a1bSJed Brown ierr = MatSetType(B,MATDENSE);CHKERRQ(ierr); 51*c4762a1bSJed Brown ierr = MatSeqDenseSetPreallocation(B,NULL);CHKERRQ(ierr); 52*c4762a1bSJed Brown ierr = MatMPIDenseSetPreallocation(B,NULL);CHKERRQ(ierr); 53*c4762a1bSJed Brown ierr = MatSetFromOptions(B);CHKERRQ(ierr); 54*c4762a1bSJed Brown ierr = MatSetRandom(B,r);CHKERRQ(ierr); 55*c4762a1bSJed Brown ierr = PetscRandomDestroy(&r);CHKERRQ(ierr); 56*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 57*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 58*c4762a1bSJed Brown 59*c4762a1bSJed Brown /* Test C = A*B (aij*dense) */ 60*c4762a1bSJed Brown ierr = MatMatMult(A,B,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr); 61*c4762a1bSJed Brown ierr = MatMatMult(A,B,MAT_REUSE_MATRIX,fill,&C);CHKERRQ(ierr); 62*c4762a1bSJed Brown 63*c4762a1bSJed Brown ierr = MatMatMultSymbolic(A,B,fill,&D);CHKERRQ(ierr); 64*c4762a1bSJed Brown for (i=0; i<2; i++) { 65*c4762a1bSJed Brown ierr = MatMatMultNumeric(A,B,D);CHKERRQ(ierr); 66*c4762a1bSJed Brown } 67*c4762a1bSJed Brown ierr = MatEqual(C,D,&equal);CHKERRQ(ierr); 68*c4762a1bSJed Brown if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"C != D"); 69*c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 70*c4762a1bSJed Brown 71*c4762a1bSJed Brown /* Test D = C*A (dense*aij) */ 72*c4762a1bSJed Brown ierr = MatMatMult(C,A,MAT_INITIAL_MATRIX,fill,&D);CHKERRQ(ierr); 73*c4762a1bSJed Brown ierr = MatMatMult(C,A,MAT_REUSE_MATRIX,fill,&D);CHKERRQ(ierr); 74*c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 75*c4762a1bSJed Brown 76*c4762a1bSJed Brown /* Test D = A*C (aij*dense) */ 77*c4762a1bSJed Brown ierr = MatMatMult(A,C,MAT_INITIAL_MATRIX,fill,&D);CHKERRQ(ierr); 78*c4762a1bSJed Brown ierr = MatMatMult(A,C,MAT_REUSE_MATRIX,fill,&D);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 80*c4762a1bSJed Brown 81*c4762a1bSJed Brown /* Test D = B*C (dense*dense) */ 82*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 83*c4762a1bSJed Brown if (size == 1) { 84*c4762a1bSJed Brown ierr = MatMatMult(B,C,MAT_INITIAL_MATRIX,fill,&D);CHKERRQ(ierr); 85*c4762a1bSJed Brown ierr = MatMatMult(B,C,MAT_REUSE_MATRIX,fill,&D);CHKERRQ(ierr); 86*c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 87*c4762a1bSJed Brown } 88*c4762a1bSJed Brown 89*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 90*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 91*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 92*c4762a1bSJed Brown ierr = PetscFinalize(); 93*c4762a1bSJed Brown return ierr; 94*c4762a1bSJed Brown } 95*c4762a1bSJed Brown 96*c4762a1bSJed Brown 97*c4762a1bSJed Brown /*TEST 98*c4762a1bSJed Brown 99*c4762a1bSJed Brown test: 100*c4762a1bSJed Brown args: -M 10 -N 10 101*c4762a1bSJed Brown output_file: output/ex109.out 102*c4762a1bSJed Brown 103*c4762a1bSJed Brown test: 104*c4762a1bSJed Brown suffix: 2 105*c4762a1bSJed Brown nsize: 3 106*c4762a1bSJed Brown output_file: output/ex109.out 107*c4762a1bSJed Brown 108*c4762a1bSJed Brown TEST*/ 109