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; 119566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&B)); 129566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B,5,6,PETSC_DETERMINE,PETSC_DETERMINE)); 139566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 149566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 159566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(B,&ms,&me)); 16c4762a1bSJed Brown for (i=ms; i<me; i++) { 179566063dSJacob Faibussowitsch PetscCall(MatSetValue(B,i,i,1.0*i,INSERT_VALUES)); 18c4762a1bSJed Brown } 199566063dSJacob Faibussowitsch PetscCall(MatSetValue(B,me-1,me,me*me,INSERT_VALUES)); 209566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 219566063dSJacob Faibussowitsch PetscCall(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; 339566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X[0],&Y)); 349566063dSJacob Faibussowitsch PetscCall(VecCopy(X[0],Y)); 359566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y,-1.0,X[1])); 369566063dSJacob Faibussowitsch PetscCall(VecNorm(Y,NORM_INFINITY,&norm)); 37c4762a1bSJed Brown 389566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-verbose",&verbose,NULL)); 39c4762a1bSJed Brown if (norm < PETSC_SQRT_MACHINE_EPSILON && verbose < 1) { 409566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%30s: norm difference < sqrt(eps_machine)\n",test)); 41c4762a1bSJed Brown } else { 429566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%30s: norm difference %g\n",test,(double)norm)); 43c4762a1bSJed Brown } 44c4762a1bSJed Brown if (verbose > 1) { 459566063dSJacob Faibussowitsch PetscCall(VecView(X[0],PETSC_VIEWER_STDOUT_WORLD)); 469566063dSJacob Faibussowitsch PetscCall(VecView(X[1],PETSC_VIEWER_STDOUT_WORLD)); 479566063dSJacob Faibussowitsch PetscCall(VecView(Y,PETSC_VIEWER_STDOUT_WORLD)); 48c4762a1bSJed Brown } 499566063dSJacob Faibussowitsch PetscCall(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; 589566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(right,2,&rtmp)); 599566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(left,2,<mp)); 609566063dSJacob Faibussowitsch PetscCall(MatScale(A,PETSC_PI)); 619566063dSJacob Faibussowitsch PetscCall(MatScale(B,PETSC_PI)); 629566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A,left,right)); 639566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(B,left,right)); 64c4762a1bSJed Brown 659566063dSJacob Faibussowitsch PetscCall(MatMult(A,X,ltmp[0])); 669566063dSJacob Faibussowitsch PetscCall(MatMult(B,X,ltmp[1])); 679566063dSJacob Faibussowitsch PetscCall(Compare2(ltmp,"MatMult")); 68c4762a1bSJed Brown 699566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A,Y,rtmp[0])); 709566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(B,Y,rtmp[1])); 719566063dSJacob Faibussowitsch PetscCall(Compare2(rtmp,"MatMultTranspose")); 72c4762a1bSJed Brown 739566063dSJacob Faibussowitsch PetscCall(VecCopy(Y1,ltmp[0])); 749566063dSJacob Faibussowitsch PetscCall(VecCopy(Y1,ltmp[1])); 759566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A,X,ltmp[0],ltmp[0])); 769566063dSJacob Faibussowitsch PetscCall(MatMultAdd(B,X,ltmp[1],ltmp[1])); 779566063dSJacob Faibussowitsch PetscCall(Compare2(ltmp,"MatMultAdd v2==v3")); 78c4762a1bSJed Brown 799566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A,X,Y1,ltmp[0])); 809566063dSJacob Faibussowitsch PetscCall(MatMultAdd(B,X,Y1,ltmp[1])); 819566063dSJacob Faibussowitsch PetscCall(Compare2(ltmp,"MatMultAdd v2!=v3")); 82c4762a1bSJed Brown 839566063dSJacob Faibussowitsch PetscCall(VecCopy(X1,rtmp[0])); 849566063dSJacob Faibussowitsch PetscCall(VecCopy(X1,rtmp[1])); 859566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(A,Y,rtmp[0],rtmp[0])); 869566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(B,Y,rtmp[1],rtmp[1])); 879566063dSJacob Faibussowitsch PetscCall(Compare2(rtmp,"MatMultTransposeAdd v2==v3")); 88c4762a1bSJed Brown 899566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(A,Y,X1,rtmp[0])); 909566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(B,Y,X1,rtmp[1])); 919566063dSJacob Faibussowitsch PetscCall(Compare2(rtmp,"MatMultTransposeAdd v2!=v3")); 92c4762a1bSJed Brown 939566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(2,<mp)); 949566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(2,&rtmp)); 95c4762a1bSJed Brown PetscFunctionReturn(0); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown 98c4762a1bSJed Brown int main(int argc, char *argv[]) 99c4762a1bSJed Brown { 100c4762a1bSJed Brown Mat A,B,Asub,Bsub; 101c4762a1bSJed Brown PetscInt ms,idxrow[3],idxcol[4]; 102c4762a1bSJed Brown Vec left,right,X,Y,X1,Y1; 103c4762a1bSJed Brown IS isrow,iscol; 104c4762a1bSJed Brown PetscBool random = PETSC_TRUE; 105c4762a1bSJed Brown 106*327415f7SBarry Smith PetscFunctionBeginUser; 1079566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 1089566063dSJacob Faibussowitsch PetscCall(AssembleMatrix(PETSC_COMM_WORLD,&A)); 1099566063dSJacob Faibussowitsch PetscCall(AssembleMatrix(PETSC_COMM_WORLD,&B)); 1109566063dSJacob Faibussowitsch PetscCall(MatSetOperation(B,MATOP_CREATE_SUBMATRIX,NULL)); 1119566063dSJacob Faibussowitsch PetscCall(MatSetOperation(B,MATOP_CREATE_SUBMATRICES,NULL)); 1129566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&ms,NULL)); 113c4762a1bSJed Brown 114c4762a1bSJed Brown idxrow[0] = ms+1; 115c4762a1bSJed Brown idxrow[1] = ms+2; 116c4762a1bSJed Brown idxrow[2] = ms+4; 1179566063dSJacob Faibussowitsch PetscCall(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; 1239566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,4,idxcol,PETSC_USE_POINTER,&iscol)); 124c4762a1bSJed Brown 1259566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A,isrow,iscol,MAT_INITIAL_MATRIX,&Asub)); 1269566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(B,isrow,iscol,MAT_INITIAL_MATRIX,&Bsub)); 127c4762a1bSJed Brown 1289566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Asub,&right,&left)); 1299566063dSJacob Faibussowitsch PetscCall(VecDuplicate(right,&X)); 1309566063dSJacob Faibussowitsch PetscCall(VecDuplicate(right,&X1)); 1319566063dSJacob Faibussowitsch PetscCall(VecDuplicate(left,&Y)); 1329566063dSJacob Faibussowitsch PetscCall(VecDuplicate(left,&Y1)); 133c4762a1bSJed Brown 1349566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-random",&random,NULL)); 135c4762a1bSJed Brown if (random) { 1369566063dSJacob Faibussowitsch PetscCall(VecSetRandom(right,NULL)); 1379566063dSJacob Faibussowitsch PetscCall(VecSetRandom(left,NULL)); 1389566063dSJacob Faibussowitsch PetscCall(VecSetRandom(X,NULL)); 1399566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Y,NULL)); 1409566063dSJacob Faibussowitsch PetscCall(VecSetRandom(X1,NULL)); 1419566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Y1,NULL)); 142c4762a1bSJed Brown } else { 1439566063dSJacob Faibussowitsch PetscCall(VecSet(right,1.0)); 1449566063dSJacob Faibussowitsch PetscCall(VecSet(left,2.0)); 1459566063dSJacob Faibussowitsch PetscCall(VecSet(X,3.0)); 1469566063dSJacob Faibussowitsch PetscCall(VecSet(Y,4.0)); 1479566063dSJacob Faibussowitsch PetscCall(VecSet(X1,3.0)); 1489566063dSJacob Faibussowitsch PetscCall(VecSet(Y1,4.0)); 149c4762a1bSJed Brown } 1509566063dSJacob Faibussowitsch PetscCall(CheckMatrices(Asub,Bsub,left,right,X,Y,X1,Y1)); 1519566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow)); 1529566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol)); 1539566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1549566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1559566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Asub)); 1569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Bsub)); 1579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&left)); 1589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&right)); 1599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 1609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 1619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X1)); 1629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y1)); 1639566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 164b122ec5aSJacob 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