1*c4762a1bSJed Brown static const char help[] = "Tests MatCreateSubMatrix with MatSubMatrix versus MatAIJ, non-square\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscmat.h> 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown static PetscErrorCode AssembleMatrix(MPI_Comm comm,Mat *A) 6*c4762a1bSJed Brown { 7*c4762a1bSJed Brown PetscErrorCode ierr; 8*c4762a1bSJed Brown Mat B; 9*c4762a1bSJed Brown PetscInt i,ms,me; 10*c4762a1bSJed Brown 11*c4762a1bSJed Brown PetscFunctionBegin; 12*c4762a1bSJed Brown ierr = MatCreate(comm,&B);CHKERRQ(ierr); 13*c4762a1bSJed Brown ierr = MatSetSizes(B,5,6,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 14*c4762a1bSJed Brown ierr = MatSetFromOptions(B);CHKERRQ(ierr); 15*c4762a1bSJed Brown ierr = MatSetUp(B);CHKERRQ(ierr); 16*c4762a1bSJed Brown ierr = MatGetOwnershipRange(B,&ms,&me);CHKERRQ(ierr); 17*c4762a1bSJed Brown for (i=ms; i<me; i++) { 18*c4762a1bSJed Brown ierr = MatSetValue(B,i,i,1.0*i,INSERT_VALUES);CHKERRQ(ierr); 19*c4762a1bSJed Brown } 20*c4762a1bSJed Brown ierr = MatSetValue(B,me-1,me,me*me,INSERT_VALUES);CHKERRQ(ierr); 21*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 22*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 23*c4762a1bSJed Brown *A = B; 24*c4762a1bSJed Brown PetscFunctionReturn(0); 25*c4762a1bSJed Brown } 26*c4762a1bSJed Brown 27*c4762a1bSJed Brown static PetscErrorCode Compare2(Vec *X,const char *test) 28*c4762a1bSJed Brown { 29*c4762a1bSJed Brown PetscErrorCode ierr; 30*c4762a1bSJed Brown PetscReal norm; 31*c4762a1bSJed Brown Vec Y; 32*c4762a1bSJed Brown PetscInt verbose = 0; 33*c4762a1bSJed Brown 34*c4762a1bSJed Brown PetscFunctionBegin; 35*c4762a1bSJed Brown ierr = VecDuplicate(X[0],&Y);CHKERRQ(ierr); 36*c4762a1bSJed Brown ierr = VecCopy(X[0],Y);CHKERRQ(ierr); 37*c4762a1bSJed Brown ierr = VecAYPX(Y,-1.0,X[1]);CHKERRQ(ierr); 38*c4762a1bSJed Brown ierr = VecNorm(Y,NORM_INFINITY,&norm);CHKERRQ(ierr); 39*c4762a1bSJed Brown 40*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-verbose",&verbose,NULL);CHKERRQ(ierr); 41*c4762a1bSJed Brown if (norm < PETSC_SQRT_MACHINE_EPSILON && verbose < 1) { 42*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%30s: norm difference < sqrt(eps_machine)\n",test);CHKERRQ(ierr); 43*c4762a1bSJed Brown } else { 44*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%30s: norm difference %g\n",test,(double)norm);CHKERRQ(ierr); 45*c4762a1bSJed Brown } 46*c4762a1bSJed Brown if (verbose > 1) { 47*c4762a1bSJed Brown ierr = VecView(X[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 48*c4762a1bSJed Brown ierr = VecView(X[1],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 49*c4762a1bSJed Brown ierr = VecView(Y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 50*c4762a1bSJed Brown } 51*c4762a1bSJed Brown ierr = VecDestroy(&Y);CHKERRQ(ierr); 52*c4762a1bSJed Brown PetscFunctionReturn(0); 53*c4762a1bSJed Brown } 54*c4762a1bSJed Brown 55*c4762a1bSJed Brown static PetscErrorCode CheckMatrices(Mat A,Mat B,Vec left,Vec right,Vec X,Vec Y,Vec X1,Vec Y1) 56*c4762a1bSJed Brown { 57*c4762a1bSJed Brown PetscErrorCode ierr; 58*c4762a1bSJed Brown Vec *ltmp,*rtmp; 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown PetscFunctionBegin; 61*c4762a1bSJed Brown ierr = VecDuplicateVecs(right,2,&rtmp);CHKERRQ(ierr); 62*c4762a1bSJed Brown ierr = VecDuplicateVecs(left,2,<mp);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = MatScale(A,PETSC_PI);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = MatScale(B,PETSC_PI);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = MatDiagonalScale(A,left,right);CHKERRQ(ierr); 66*c4762a1bSJed Brown ierr = MatDiagonalScale(B,left,right);CHKERRQ(ierr); 67*c4762a1bSJed Brown 68*c4762a1bSJed Brown ierr = MatMult(A,X,ltmp[0]);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = MatMult(B,X,ltmp[1]);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = Compare2(ltmp,"MatMult");CHKERRQ(ierr); 71*c4762a1bSJed Brown 72*c4762a1bSJed Brown ierr = MatMultTranspose(A,Y,rtmp[0]);CHKERRQ(ierr); 73*c4762a1bSJed Brown ierr = MatMultTranspose(B,Y,rtmp[1]);CHKERRQ(ierr); 74*c4762a1bSJed Brown ierr = Compare2(rtmp,"MatMultTranspose");CHKERRQ(ierr); 75*c4762a1bSJed Brown 76*c4762a1bSJed Brown ierr = VecCopy(Y1,ltmp[0]);CHKERRQ(ierr); 77*c4762a1bSJed Brown ierr = VecCopy(Y1,ltmp[1]);CHKERRQ(ierr); 78*c4762a1bSJed Brown ierr = MatMultAdd(A,X,ltmp[0],ltmp[0]);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = MatMultAdd(B,X,ltmp[1],ltmp[1]);CHKERRQ(ierr); 80*c4762a1bSJed Brown ierr = Compare2(ltmp,"MatMultAdd v2==v3");CHKERRQ(ierr); 81*c4762a1bSJed Brown 82*c4762a1bSJed Brown ierr = MatMultAdd(A,X,Y1,ltmp[0]);CHKERRQ(ierr); 83*c4762a1bSJed Brown ierr = MatMultAdd(B,X,Y1,ltmp[1]);CHKERRQ(ierr); 84*c4762a1bSJed Brown ierr = Compare2(ltmp,"MatMultAdd v2!=v3");CHKERRQ(ierr); 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown ierr = VecCopy(X1,rtmp[0]);CHKERRQ(ierr); 87*c4762a1bSJed Brown ierr = VecCopy(X1,rtmp[1]);CHKERRQ(ierr); 88*c4762a1bSJed Brown ierr = MatMultTransposeAdd(A,Y,rtmp[0],rtmp[0]);CHKERRQ(ierr); 89*c4762a1bSJed Brown ierr = MatMultTransposeAdd(B,Y,rtmp[1],rtmp[1]);CHKERRQ(ierr); 90*c4762a1bSJed Brown ierr = Compare2(rtmp,"MatMultTransposeAdd v2==v3");CHKERRQ(ierr); 91*c4762a1bSJed Brown 92*c4762a1bSJed Brown ierr = MatMultTransposeAdd(A,Y,X1,rtmp[0]);CHKERRQ(ierr); 93*c4762a1bSJed Brown ierr = MatMultTransposeAdd(B,Y,X1,rtmp[1]);CHKERRQ(ierr); 94*c4762a1bSJed Brown ierr = Compare2(rtmp,"MatMultTransposeAdd v2!=v3");CHKERRQ(ierr); 95*c4762a1bSJed Brown 96*c4762a1bSJed Brown ierr = VecDestroyVecs(2,<mp);CHKERRQ(ierr); 97*c4762a1bSJed Brown ierr = VecDestroyVecs(2,&rtmp);CHKERRQ(ierr); 98*c4762a1bSJed Brown PetscFunctionReturn(0); 99*c4762a1bSJed Brown } 100*c4762a1bSJed Brown 101*c4762a1bSJed Brown int main(int argc, char *argv[]) 102*c4762a1bSJed Brown { 103*c4762a1bSJed Brown PetscErrorCode ierr; 104*c4762a1bSJed Brown Mat A,B,Asub,Bsub; 105*c4762a1bSJed Brown PetscInt ms,idxrow[3],idxcol[4]; 106*c4762a1bSJed Brown Vec left,right,X,Y,X1,Y1; 107*c4762a1bSJed Brown IS isrow,iscol; 108*c4762a1bSJed Brown PetscBool random = PETSC_TRUE; 109*c4762a1bSJed Brown 110*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 111*c4762a1bSJed Brown ierr = AssembleMatrix(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 112*c4762a1bSJed Brown ierr = AssembleMatrix(PETSC_COMM_WORLD,&B);CHKERRQ(ierr); 113*c4762a1bSJed Brown ierr = MatSetOperation(B,MATOP_CREATE_SUBMATRIX,NULL);CHKERRQ(ierr); 114*c4762a1bSJed Brown ierr = MatSetOperation(B,MATOP_CREATE_SUBMATRICES,NULL);CHKERRQ(ierr); 115*c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&ms,NULL);CHKERRQ(ierr); 116*c4762a1bSJed Brown 117*c4762a1bSJed Brown idxrow[0] = ms+1; 118*c4762a1bSJed Brown idxrow[1] = ms+2; 119*c4762a1bSJed Brown idxrow[2] = ms+4; 120*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_WORLD,3,idxrow,PETSC_USE_POINTER,&isrow);CHKERRQ(ierr); 121*c4762a1bSJed Brown 122*c4762a1bSJed Brown idxcol[0] = ms+1; 123*c4762a1bSJed Brown idxcol[1] = ms+2; 124*c4762a1bSJed Brown idxcol[2] = ms+4; 125*c4762a1bSJed Brown idxcol[3] = ms+5; 126*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_WORLD,4,idxcol,PETSC_USE_POINTER,&iscol);CHKERRQ(ierr); 127*c4762a1bSJed Brown 128*c4762a1bSJed Brown ierr = MatCreateSubMatrix(A,isrow,iscol,MAT_INITIAL_MATRIX,&Asub);CHKERRQ(ierr); 129*c4762a1bSJed Brown ierr = MatCreateSubMatrix(B,isrow,iscol,MAT_INITIAL_MATRIX,&Bsub);CHKERRQ(ierr); 130*c4762a1bSJed Brown 131*c4762a1bSJed Brown ierr = MatCreateVecs(Asub,&right,&left);CHKERRQ(ierr); 132*c4762a1bSJed Brown ierr = VecDuplicate(right,&X);CHKERRQ(ierr); 133*c4762a1bSJed Brown ierr = VecDuplicate(right,&X1);CHKERRQ(ierr); 134*c4762a1bSJed Brown ierr = VecDuplicate(left,&Y);CHKERRQ(ierr); 135*c4762a1bSJed Brown ierr = VecDuplicate(left,&Y1);CHKERRQ(ierr); 136*c4762a1bSJed Brown 137*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-random",&random,NULL);CHKERRQ(ierr); 138*c4762a1bSJed Brown if (random) { 139*c4762a1bSJed Brown ierr = VecSetRandom(right,NULL);CHKERRQ(ierr); 140*c4762a1bSJed Brown ierr = VecSetRandom(left,NULL);CHKERRQ(ierr); 141*c4762a1bSJed Brown ierr = VecSetRandom(X,NULL);CHKERRQ(ierr); 142*c4762a1bSJed Brown ierr = VecSetRandom(Y,NULL);CHKERRQ(ierr); 143*c4762a1bSJed Brown ierr = VecSetRandom(X1,NULL);CHKERRQ(ierr); 144*c4762a1bSJed Brown ierr = VecSetRandom(Y1,NULL);CHKERRQ(ierr); 145*c4762a1bSJed Brown } else { 146*c4762a1bSJed Brown ierr = VecSet(right,1.0);CHKERRQ(ierr); 147*c4762a1bSJed Brown ierr = VecSet(left,2.0);CHKERRQ(ierr); 148*c4762a1bSJed Brown ierr = VecSet(X,3.0);CHKERRQ(ierr); 149*c4762a1bSJed Brown ierr = VecSet(Y,4.0);CHKERRQ(ierr); 150*c4762a1bSJed Brown ierr = VecSet(X1,3.0);CHKERRQ(ierr); 151*c4762a1bSJed Brown ierr = VecSet(Y1,4.0);CHKERRQ(ierr); 152*c4762a1bSJed Brown } 153*c4762a1bSJed Brown ierr = CheckMatrices(Asub,Bsub,left,right,X,Y,X1,Y1);CHKERRQ(ierr); 154*c4762a1bSJed Brown ierr = ISDestroy(&isrow);CHKERRQ(ierr); 155*c4762a1bSJed Brown ierr = ISDestroy(&iscol);CHKERRQ(ierr); 156*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 157*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 158*c4762a1bSJed Brown ierr = MatDestroy(&Asub);CHKERRQ(ierr); 159*c4762a1bSJed Brown ierr = MatDestroy(&Bsub);CHKERRQ(ierr); 160*c4762a1bSJed Brown ierr = VecDestroy(&left);CHKERRQ(ierr); 161*c4762a1bSJed Brown ierr = VecDestroy(&right);CHKERRQ(ierr); 162*c4762a1bSJed Brown ierr = VecDestroy(&X);CHKERRQ(ierr); 163*c4762a1bSJed Brown ierr = VecDestroy(&Y);CHKERRQ(ierr); 164*c4762a1bSJed Brown ierr = VecDestroy(&X1);CHKERRQ(ierr); 165*c4762a1bSJed Brown ierr = VecDestroy(&Y1);CHKERRQ(ierr); 166*c4762a1bSJed Brown ierr = PetscFinalize(); 167*c4762a1bSJed Brown return ierr; 168*c4762a1bSJed Brown } 169*c4762a1bSJed Brown 170*c4762a1bSJed Brown 171*c4762a1bSJed Brown 172*c4762a1bSJed Brown /*TEST 173*c4762a1bSJed Brown 174*c4762a1bSJed Brown test: 175*c4762a1bSJed Brown nsize: 3 176*c4762a1bSJed Brown 177*c4762a1bSJed Brown TEST*/ 178