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=10,N=5; 9*c4762a1bSJed Brown PetscErrorCode ierr; 10*c4762a1bSJed Brown PetscRandom r; 11*c4762a1bSJed Brown PetscBool equal=PETSC_FALSE; 12*c4762a1bSJed Brown PetscReal fill = 1.0; 13*c4762a1bSJed Brown PetscInt nza,am,an; 14*c4762a1bSJed Brown 15*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 16*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr); 17*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);CHKERRQ(ierr); 18*c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL);CHKERRQ(ierr); 19*c4762a1bSJed Brown 20*c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&r);CHKERRQ(ierr); 21*c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr); 22*c4762a1bSJed Brown 23*c4762a1bSJed Brown /* create a aij matrix A */ 24*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 25*c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,M);CHKERRQ(ierr); 26*c4762a1bSJed Brown ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr); 27*c4762a1bSJed Brown nza = (PetscInt)(.3*M); /* num of nozeros in each row of A */ 28*c4762a1bSJed Brown ierr = MatSeqAIJSetPreallocation(A,nza,NULL);CHKERRQ(ierr); 29*c4762a1bSJed Brown ierr = MatMPIAIJSetPreallocation(A,nza,NULL,nza,NULL);CHKERRQ(ierr); 30*c4762a1bSJed Brown ierr = MatSetRandom(A,r);CHKERRQ(ierr); 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown /* create a dense matrix B */ 33*c4762a1bSJed Brown ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr); 34*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr); 35*c4762a1bSJed Brown ierr = MatSetSizes(B,PETSC_DECIDE,am,N,PETSC_DECIDE);CHKERRQ(ierr); 36*c4762a1bSJed Brown ierr = MatSetType(B,MATDENSE);CHKERRQ(ierr); 37*c4762a1bSJed Brown ierr = MatSeqDenseSetPreallocation(B,NULL);CHKERRQ(ierr); 38*c4762a1bSJed Brown ierr = MatMPIDenseSetPreallocation(B,NULL);CHKERRQ(ierr); 39*c4762a1bSJed Brown ierr = MatSetRandom(B,r);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = PetscRandomDestroy(&r);CHKERRQ(ierr); 41*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 42*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 43*c4762a1bSJed Brown 44*c4762a1bSJed Brown 45*c4762a1bSJed Brown /* Test MatMatMult() */ 46*c4762a1bSJed Brown ierr = MatMatMult(B,A,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr); 47*c4762a1bSJed Brown 48*c4762a1bSJed Brown ierr = MatMatMultSymbolic(B,A,fill,&D);CHKERRQ(ierr); 49*c4762a1bSJed Brown for (i=0; i<2; i++) { 50*c4762a1bSJed Brown ierr = MatMatMultNumeric(B,A,D);CHKERRQ(ierr); 51*c4762a1bSJed Brown } 52*c4762a1bSJed Brown ierr = MatEqual(C,D,&equal);CHKERRQ(ierr); 53*c4762a1bSJed Brown if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"C != D"); 54*c4762a1bSJed Brown 55*c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 56*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 57*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 58*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 59*c4762a1bSJed Brown ierr = PetscFinalize(); 60*c4762a1bSJed Brown return ierr; 61*c4762a1bSJed Brown } 62*c4762a1bSJed Brown 63*c4762a1bSJed Brown 64*c4762a1bSJed Brown /*TEST 65*c4762a1bSJed Brown 66*c4762a1bSJed Brown test: 67*c4762a1bSJed Brown output_file: output/ex122.out 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown TEST*/ 70