1c4762a1bSJed Brown static char help[] = "Tests the use of MatTranspose_Nest and MatMatMult_Nest_Dense\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown PetscErrorCode TestInitialMatrix(void) 6c4762a1bSJed Brown { 7c4762a1bSJed Brown const PetscInt nr = 2,nc = 3,nk = 10; 8c4762a1bSJed Brown PetscInt n,N,m,M; 9c4762a1bSJed Brown const PetscInt arow[2*3] = { 2,2,2,3,3,3 }; 10c4762a1bSJed Brown const PetscInt acol[2*3] = { 3,2,4,3,2,4 }; 11c4762a1bSJed Brown Mat A,Atranspose,B,C; 12c4762a1bSJed Brown Mat subs[2*3],**block; 13c4762a1bSJed Brown Vec x,y,Ax,ATy; 14c4762a1bSJed Brown PetscInt i,j; 15c20d7725SJed Brown PetscScalar dot1,dot2,zero = 0.0,one = 1.0,*valsB,*valsC; 16c4762a1bSJed Brown PetscReal norm; 17c4762a1bSJed Brown PetscRandom rctx; 18c20d7725SJed Brown PetscBool equal; 19c4762a1bSJed Brown 20c4762a1bSJed Brown PetscFunctionBegin; 21*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rctx)); 22c4762a1bSJed Brown /* Force the random numbers to have imaginary part 0 so printed results are the same for --with-scalar-type=real or --with-scalar-type=complex */ 23*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(rctx,zero,one)); 24*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rctx)); 25c4762a1bSJed Brown for (i=0; i<(nr * nc); i++) { 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqDense(PETSC_COMM_WORLD,arow[i],acol[i],NULL,&subs[i])); 27c4762a1bSJed Brown } 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateNest(PETSC_COMM_WORLD,nr,NULL,nc,NULL,subs,&A)); 29*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A, &x, NULL)); 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A, NULL, &y)); 31*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x, &ATy)); 32*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(y, &Ax)); 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(A,rctx)); 34*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(A, MAT_INITIAL_MATRIX, &Atranspose)); 35c4762a1bSJed Brown 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNestGetSubMats(A, NULL, NULL, &block)); 38c4762a1bSJed Brown for (i=0; i<nr; i++) { 39c4762a1bSJed Brown for (j=0; j<nc; j++) { 40*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(block[i][j], PETSC_VIEWER_STDOUT_WORLD)); 41c4762a1bSJed Brown } 42c4762a1bSJed Brown } 43c4762a1bSJed Brown 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(Atranspose, PETSC_VIEWER_STDOUT_WORLD)); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNestGetSubMats(Atranspose, NULL, NULL, &block)); 46c4762a1bSJed Brown for (i=0; i<nc; i++) { 47c4762a1bSJed Brown for (j=0; j<nr; j++) { 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(block[i][j], PETSC_VIEWER_STDOUT_WORLD)); 49c4762a1bSJed Brown } 50c4762a1bSJed Brown } 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Check <Ax, y> = <x, A^Ty> */ 53c4762a1bSJed Brown for (i=0; i<10; i++) { 54*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(x,rctx)); 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(y,rctx)); 56c4762a1bSJed Brown 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A, x, Ax)); 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDot(Ax, y, &dot1)); 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(Atranspose, y, ATy)); 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDot(ATy, x, &dot2)); 61c4762a1bSJed Brown 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "<Ax, y> = %g\n", (double)PetscRealPart(dot1))); 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "<x, A^Ty> = %g\n",(double)PetscRealPart(dot2))); 64c4762a1bSJed Brown } 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Ax)); 68c4762a1bSJed Brown 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqDense(PETSC_COMM_WORLD,acol[0]+acol[nr]+acol[2*nr],nk,NULL,&B)); 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(B,rctx)); 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C)); 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C)); 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(A,B,C,10,&equal)); 742c71b3e2SJacob Faibussowitsch PetscCheckFalse(!equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in C != A*B"); 75c20d7725SJed Brown 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&M,&N)); 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A,&m,&n)); 78c4762a1bSJed Brown for (i=0; i<nk; i++) { 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetColumn(B,i,&valsB)); 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,N,valsB,&x)); 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A,NULL,&Ax)); 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,x,Ax)); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetColumn(C,i,&valsC)); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,m,M,valsC,&y)); 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,-1.0,Ax)); 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Ax)); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreColumn(C,&valsC)); 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreColumn(B,&valsB)); 91c4762a1bSJed Brown } 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(C,NORM_INFINITY,&norm)); 932c71b3e2SJacob Faibussowitsch PetscCheckFalse(norm > PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in MatMatMult(): %g",(double)norm); 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 96c4762a1bSJed Brown 97c4762a1bSJed Brown for (i=0; i<(nr * nc); i++) { 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&subs[i])); 99c4762a1bSJed Brown } 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Atranspose)); 102*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ATy)); 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rctx)); 104c4762a1bSJed Brown PetscFunctionReturn(0); 105c4762a1bSJed Brown } 106c4762a1bSJed Brown 107c4762a1bSJed Brown PetscErrorCode TestReuseMatrix(void) 108c4762a1bSJed Brown { 109c4762a1bSJed Brown const PetscInt n = 2; 110c4762a1bSJed Brown Mat A; 111c4762a1bSJed Brown Mat subs[2*2],**block; 112c4762a1bSJed Brown PetscInt i,j; 113c4762a1bSJed Brown PetscRandom rctx; 114c4762a1bSJed Brown PetscScalar zero = 0.0, one = 1.0; 115c4762a1bSJed Brown 116c4762a1bSJed Brown PetscFunctionBegin; 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rctx)); 118*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(rctx,zero,one)); 119*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rctx)); 120c4762a1bSJed Brown for (i=0; i<(n * n); i++) { 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqDense(PETSC_COMM_WORLD,n,n,NULL,&subs[i])); 122c4762a1bSJed Brown } 123*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateNest(PETSC_COMM_WORLD,n,NULL,n,NULL,subs,&A)); 124*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(A,rctx)); 125c4762a1bSJed Brown 126*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNestGetSubMats(A, NULL, NULL, &block)); 128c4762a1bSJed Brown for (i=0; i<n; i++) { 129c4762a1bSJed Brown for (j=0; j<n; j++) { 130*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(block[i][j], PETSC_VIEWER_STDOUT_WORLD)); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown } 133*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(A,MAT_INPLACE_MATRIX,&A)); 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 135*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNestGetSubMats(A, NULL, NULL, &block)); 136c4762a1bSJed Brown for (i=0; i<n; i++) { 137c4762a1bSJed Brown for (j=0; j<n; j++) { 138*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(block[i][j], PETSC_VIEWER_STDOUT_WORLD)); 139c4762a1bSJed Brown } 140c4762a1bSJed Brown } 141c4762a1bSJed Brown 142c4762a1bSJed Brown for (i=0; i<(n * n); i++) { 143*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&subs[i])); 144c4762a1bSJed Brown } 145*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 146*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rctx)); 147c4762a1bSJed Brown PetscFunctionReturn(0); 148c4762a1bSJed Brown } 149c4762a1bSJed Brown 150c4762a1bSJed Brown int main(int argc,char **args) 151c4762a1bSJed Brown { 152c4762a1bSJed Brown PetscErrorCode ierr; 153c4762a1bSJed Brown 154c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 155*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestInitialMatrix()); 156*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestReuseMatrix()); 157c4762a1bSJed Brown ierr = PetscFinalize(); 158c4762a1bSJed Brown return ierr; 159c4762a1bSJed Brown } 160c4762a1bSJed Brown 161c4762a1bSJed Brown /*TEST 162c4762a1bSJed Brown 163c4762a1bSJed Brown test: 164c4762a1bSJed Brown args: -malloc_dump 165c4762a1bSJed Brown 166c4762a1bSJed Brown TEST*/ 167