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; 11*9566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&B)); 12*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B,6,6,PETSC_DETERMINE,PETSC_DETERMINE)); 13*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 14*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 15*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(B,&ms,&me)); 16c4762a1bSJed Brown for (i=ms; i<me; i++) { 17*9566063dSJacob Faibussowitsch PetscCall(MatSetValue(B,i,i,1.0*i,INSERT_VALUES)); 18c4762a1bSJed Brown } 19*9566063dSJacob Faibussowitsch PetscCall(MatSetValue(B,me-1,me-1,me*me,INSERT_VALUES)); 20*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 21*9566063dSJacob 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; 33*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X[0],&Y)); 34*9566063dSJacob Faibussowitsch PetscCall(VecCopy(X[0],Y)); 35*9566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y,-1.0,X[1])); 36*9566063dSJacob Faibussowitsch PetscCall(VecNorm(Y,NORM_INFINITY,&norm)); 37c4762a1bSJed Brown 38*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-verbose",&verbose,NULL)); 39c4762a1bSJed Brown if (norm < PETSC_SQRT_MACHINE_EPSILON && verbose < 1) { 40*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%30s: norm difference < sqrt(eps_machine)\n",test)); 41c4762a1bSJed Brown } else { 42*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%30s: norm difference %g\n",test,(double)norm)); 43c4762a1bSJed Brown } 44c4762a1bSJed Brown if (verbose > 1) { 45*9566063dSJacob Faibussowitsch PetscCall(VecView(X[0],PETSC_VIEWER_STDOUT_WORLD)); 46*9566063dSJacob Faibussowitsch PetscCall(VecView(X[1],PETSC_VIEWER_STDOUT_WORLD)); 47*9566063dSJacob Faibussowitsch PetscCall(VecView(Y,PETSC_VIEWER_STDOUT_WORLD)); 48c4762a1bSJed Brown } 49*9566063dSJacob 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; 58*9566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(right,2,&rtmp)); 59*9566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(left,2,<mp)); 60*9566063dSJacob Faibussowitsch PetscCall(MatScale(A,PETSC_PI)); 61*9566063dSJacob Faibussowitsch PetscCall(MatScale(B,PETSC_PI)); 62*9566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A,left,right)); 63*9566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(B,left,right)); 64*9566063dSJacob Faibussowitsch PetscCall(MatShift(A,PETSC_PI)); 65*9566063dSJacob Faibussowitsch PetscCall(MatShift(B,PETSC_PI)); 66c4762a1bSJed Brown 67*9566063dSJacob Faibussowitsch PetscCall(MatMult(A,X,ltmp[0])); 68*9566063dSJacob Faibussowitsch PetscCall(MatMult(B,X,ltmp[1])); 69*9566063dSJacob Faibussowitsch PetscCall(Compare2(ltmp,"MatMult")); 70c4762a1bSJed Brown 71*9566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A,Y,rtmp[0])); 72*9566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(B,Y,rtmp[1])); 73*9566063dSJacob Faibussowitsch PetscCall(Compare2(rtmp,"MatMultTranspose")); 74c4762a1bSJed Brown 75*9566063dSJacob Faibussowitsch PetscCall(VecCopy(Y1,ltmp[0])); 76*9566063dSJacob Faibussowitsch PetscCall(VecCopy(Y1,ltmp[1])); 77*9566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A,X,ltmp[0],ltmp[0])); 78*9566063dSJacob Faibussowitsch PetscCall(MatMultAdd(B,X,ltmp[1],ltmp[1])); 79*9566063dSJacob Faibussowitsch PetscCall(Compare2(ltmp,"MatMultAdd v2==v3")); 80c4762a1bSJed Brown 81*9566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A,X,Y1,ltmp[0])); 82*9566063dSJacob Faibussowitsch PetscCall(MatMultAdd(B,X,Y1,ltmp[1])); 83*9566063dSJacob Faibussowitsch PetscCall(Compare2(ltmp,"MatMultAdd v2!=v3")); 84c4762a1bSJed Brown 85*9566063dSJacob Faibussowitsch PetscCall(VecCopy(X1,rtmp[0])); 86*9566063dSJacob Faibussowitsch PetscCall(VecCopy(X1,rtmp[1])); 87*9566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(A,Y,rtmp[0],rtmp[0])); 88*9566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(B,Y,rtmp[1],rtmp[1])); 89*9566063dSJacob Faibussowitsch PetscCall(Compare2(rtmp,"MatMultTransposeAdd v2==v3")); 90c4762a1bSJed Brown 91*9566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(A,Y,X1,rtmp[0])); 92*9566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(B,Y,X1,rtmp[1])); 93*9566063dSJacob Faibussowitsch PetscCall(Compare2(rtmp,"MatMultTransposeAdd v2!=v3")); 94c4762a1bSJed Brown 95*9566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(2,<mp)); 96*9566063dSJacob Faibussowitsch PetscCall(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*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 109*9566063dSJacob Faibussowitsch PetscCall(AssembleMatrix(PETSC_COMM_WORLD,&A)); 110*9566063dSJacob Faibussowitsch PetscCall(AssembleMatrix(PETSC_COMM_WORLD,&B)); 111*9566063dSJacob Faibussowitsch PetscCall(MatSetOperation(B,MATOP_CREATE_SUBMATRIX,NULL)); 112*9566063dSJacob Faibussowitsch PetscCall(MatSetOperation(B,MATOP_CREATE_SUBMATRICES,NULL)); 113*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&ms,NULL)); 114c4762a1bSJed Brown 115c4762a1bSJed Brown idxrow[0] = ms+1; 116c4762a1bSJed Brown idxrow[1] = ms+2; 117c4762a1bSJed Brown idxrow[2] = ms+4; 118*9566063dSJacob Faibussowitsch PetscCall(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; 123*9566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,3,idxcol,PETSC_USE_POINTER,&iscol)); 124c4762a1bSJed Brown 125*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A,isrow,iscol,MAT_INITIAL_MATRIX,&Asub)); 126*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(B,isrow,iscol,MAT_INITIAL_MATRIX,&Bsub)); 127c4762a1bSJed Brown 128*9566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Asub,&right,&left)); 129*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(right,&X)); 130*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(right,&X1)); 131*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(left,&Y)); 132*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(left,&Y1)); 133c4762a1bSJed Brown 134*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-random",&random,NULL)); 135c4762a1bSJed Brown if (random) { 136*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(right,NULL)); 137*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(left,NULL)); 138*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(X,NULL)); 139*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Y,NULL)); 140*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(X1,NULL)); 141*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Y1,NULL)); 142c4762a1bSJed Brown } else { 143*9566063dSJacob Faibussowitsch PetscCall(VecSet(right,1.0)); 144*9566063dSJacob Faibussowitsch PetscCall(VecSet(left,2.0)); 145*9566063dSJacob Faibussowitsch PetscCall(VecSet(X,3.0)); 146*9566063dSJacob Faibussowitsch PetscCall(VecSet(Y,4.0)); 147*9566063dSJacob Faibussowitsch PetscCall(VecSet(X1,3.0)); 148*9566063dSJacob Faibussowitsch PetscCall(VecSet(Y1,4.0)); 149c4762a1bSJed Brown } 150*9566063dSJacob Faibussowitsch PetscCall(CheckMatrices(Asub,Bsub,left,right,X,Y,X1,Y1)); 151*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow)); 152*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol)); 153*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 154*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 155*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Asub)); 156*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Bsub)); 157*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&left)); 158*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&right)); 159*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 160*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 161*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X1)); 162*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y1)); 163*9566063dSJacob 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