1c4762a1bSJed Brown static char help[] = "Test MatMatMult() for AIJ and Dense matrices.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown int main(int argc,char **argv) 6c4762a1bSJed Brown { 7c20d7725SJed Brown Mat A,B,C; 8c20d7725SJed Brown PetscInt M=10,N=5; 9c4762a1bSJed Brown PetscRandom r; 10c4762a1bSJed Brown PetscBool equal=PETSC_FALSE; 11c4762a1bSJed Brown PetscReal fill = 1.0; 12c4762a1bSJed Brown PetscInt nza,am,an; 13c4762a1bSJed Brown 14*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL)); 175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL)); 18c4762a1bSJed Brown 195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&r)); 205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(r)); 21c4762a1bSJed Brown 22c4762a1bSJed Brown /* create a aij matrix A */ 235f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 245f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,M)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATAIJ)); 26c4762a1bSJed Brown nza = (PetscInt)(.3*M); /* num of nozeros in each row of A */ 275f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(A,nza,NULL)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJSetPreallocation(A,nza,NULL,nza,NULL)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(A,r)); 30c4762a1bSJed Brown 31c4762a1bSJed Brown /* create a dense matrix B */ 325f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A,&am,&an)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(B,PETSC_DECIDE,am,N,PETSC_DECIDE)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(B,MATDENSE)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqDenseSetPreallocation(B,NULL)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIDenseSetPreallocation(B,NULL)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(B,r)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&r)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* Test MatMatMult() */ 445f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(B,A,MAT_INITIAL_MATRIX,fill,&C)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(B,A,MAT_REUSE_MATRIX,fill,&C)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(B,A,C,10,&equal)); 4728b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"C != B*A"); 48c4762a1bSJed Brown 495f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 505f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 515f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 52*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 53*b122ec5aSJacob Faibussowitsch return 0; 54c4762a1bSJed Brown } 55c4762a1bSJed Brown 56c4762a1bSJed Brown /*TEST 57c4762a1bSJed Brown 58c4762a1bSJed Brown test: 59c4762a1bSJed Brown output_file: output/ex122.out 60c4762a1bSJed Brown 61c4762a1bSJed Brown TEST*/ 62