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