1c4762a1bSJed Brown static const char help[] = "Tests MatCreateSubMatrix with MatSubMatrix versus MatAIJ, square, shifted (copied from ex97)\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; 115f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(comm,&B)); 125f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(B,6,6,PETSC_DETERMINE,PETSC_DETERMINE)); 135f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(B)); 145f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(B)); 155f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(B,&ms,&me)); 16c4762a1bSJed Brown for (i=ms; i<me; i++) { 175f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(B,i,i,1.0*i,INSERT_VALUES)); 18c4762a1bSJed Brown } 195f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(B,me-1,me-1,me*me,INSERT_VALUES)); 205f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 215f80ce2aSJacob 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; 335f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(X[0],&Y)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(X[0],Y)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(VecAYPX(Y,-1.0,X[1])); 365f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(Y,NORM_INFINITY,&norm)); 37c4762a1bSJed Brown 385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-verbose",&verbose,NULL)); 39c4762a1bSJed Brown if (norm < PETSC_SQRT_MACHINE_EPSILON && verbose < 1) { 405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%30s: norm difference < sqrt(eps_machine)\n",test)); 41c4762a1bSJed Brown } else { 425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%30s: norm difference %g\n",test,(double)norm)); 43c4762a1bSJed Brown } 44c4762a1bSJed Brown if (verbose > 1) { 455f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(X[0],PETSC_VIEWER_STDOUT_WORLD)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(X[1],PETSC_VIEWER_STDOUT_WORLD)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(Y,PETSC_VIEWER_STDOUT_WORLD)); 48c4762a1bSJed Brown } 495f80ce2aSJacob 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; 585f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicateVecs(right,2,&rtmp)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicateVecs(left,2,<mp)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,PETSC_PI)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(B,PETSC_PI)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(A,left,right)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(B,left,right)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(MatShift(A,PETSC_PI)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(MatShift(B,PETSC_PI)); 66c4762a1bSJed Brown 675f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,X,ltmp[0])); 685f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(B,X,ltmp[1])); 695f80ce2aSJacob Faibussowitsch CHKERRQ(Compare2(ltmp,"MatMult")); 70c4762a1bSJed Brown 715f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(A,Y,rtmp[0])); 725f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(B,Y,rtmp[1])); 735f80ce2aSJacob Faibussowitsch CHKERRQ(Compare2(rtmp,"MatMultTranspose")); 74c4762a1bSJed Brown 755f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(Y1,ltmp[0])); 765f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(Y1,ltmp[1])); 775f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd(A,X,ltmp[0],ltmp[0])); 785f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd(B,X,ltmp[1],ltmp[1])); 795f80ce2aSJacob Faibussowitsch CHKERRQ(Compare2(ltmp,"MatMultAdd v2==v3")); 80c4762a1bSJed Brown 815f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd(A,X,Y1,ltmp[0])); 825f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd(B,X,Y1,ltmp[1])); 835f80ce2aSJacob Faibussowitsch CHKERRQ(Compare2(ltmp,"MatMultAdd v2!=v3")); 84c4762a1bSJed Brown 855f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(X1,rtmp[0])); 865f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(X1,rtmp[1])); 875f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAdd(A,Y,rtmp[0],rtmp[0])); 885f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAdd(B,Y,rtmp[1],rtmp[1])); 895f80ce2aSJacob Faibussowitsch CHKERRQ(Compare2(rtmp,"MatMultTransposeAdd v2==v3")); 90c4762a1bSJed Brown 915f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAdd(A,Y,X1,rtmp[0])); 925f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAdd(B,Y,X1,rtmp[1])); 935f80ce2aSJacob Faibussowitsch CHKERRQ(Compare2(rtmp,"MatMultTransposeAdd v2!=v3")); 94c4762a1bSJed Brown 955f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroyVecs(2,<mp)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroyVecs(2,&rtmp)); 97c4762a1bSJed Brown PetscFunctionReturn(0); 98c4762a1bSJed Brown } 99c4762a1bSJed Brown 100c4762a1bSJed Brown int main(int argc, char *argv[]) 101c4762a1bSJed Brown { 102c4762a1bSJed Brown Mat A,B,Asub,Bsub; 103c4762a1bSJed Brown PetscInt ms,idxrow[3],idxcol[3]; 104c4762a1bSJed Brown Vec left,right,X,Y,X1,Y1; 105c4762a1bSJed Brown IS isrow,iscol; 106c4762a1bSJed Brown PetscBool random = PETSC_TRUE; 107c4762a1bSJed Brown 108*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,NULL,help)); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(AssembleMatrix(PETSC_COMM_WORLD,&A)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(AssembleMatrix(PETSC_COMM_WORLD,&B)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOperation(B,MATOP_CREATE_SUBMATRIX,NULL)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOperation(B,MATOP_CREATE_SUBMATRICES,NULL)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&ms,NULL)); 114c4762a1bSJed Brown 115c4762a1bSJed Brown idxrow[0] = ms+1; 116c4762a1bSJed Brown idxrow[1] = ms+2; 117c4762a1bSJed Brown idxrow[2] = ms+4; 1185f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,3,idxrow,PETSC_USE_POINTER,&isrow)); 119c4762a1bSJed Brown 120c4762a1bSJed Brown idxcol[0] = ms+1; 121c4762a1bSJed Brown idxcol[1] = ms+2; 122c4762a1bSJed Brown idxcol[2] = ms+4; 1235f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,3,idxcol,PETSC_USE_POINTER,&iscol)); 124c4762a1bSJed Brown 1255f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrix(A,isrow,iscol,MAT_INITIAL_MATRIX,&Asub)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrix(B,isrow,iscol,MAT_INITIAL_MATRIX,&Bsub)); 127c4762a1bSJed Brown 1285f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(Asub,&right,&left)); 1295f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(right,&X)); 1305f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(right,&X1)); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(left,&Y)); 1325f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(left,&Y1)); 133c4762a1bSJed Brown 1345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-random",&random,NULL)); 135c4762a1bSJed Brown if (random) { 1365f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(right,NULL)); 1375f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(left,NULL)); 1385f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(X,NULL)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(Y,NULL)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(X1,NULL)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(Y1,NULL)); 142c4762a1bSJed Brown } else { 1435f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(right,1.0)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(left,2.0)); 1455f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(X,3.0)); 1465f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(Y,4.0)); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(X1,3.0)); 1485f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(Y1,4.0)); 149c4762a1bSJed Brown } 1505f80ce2aSJacob Faibussowitsch CHKERRQ(CheckMatrices(Asub,Bsub,left,right,X,Y,X1,Y1)); 1515f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isrow)); 1525f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iscol)); 1535f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 1545f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 1555f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Asub)); 1565f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Bsub)); 1575f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&left)); 1585f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&right)); 1595f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&X)); 1605f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Y)); 1615f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&X1)); 1625f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Y1)); 163*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 164*b122ec5aSJacob Faibussowitsch return 0; 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 167c4762a1bSJed Brown /*TEST 168c4762a1bSJed Brown 169c4762a1bSJed Brown test: 170c4762a1bSJed Brown nsize: 3 171c4762a1bSJed Brown 172c4762a1bSJed Brown TEST*/ 173