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