1c4762a1bSJed Brown static char help[] = "Tests various routines for MATSHELL\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown typedef struct _n_User *User; 6c4762a1bSJed Brown struct _n_User { 7c4762a1bSJed Brown Mat B; 8c4762a1bSJed Brown }; 9c4762a1bSJed Brown 10c4762a1bSJed Brown static PetscErrorCode MatGetDiagonal_User(Mat A,Vec X) 11c4762a1bSJed Brown { 12c4762a1bSJed Brown User user; 13c4762a1bSJed Brown 14c4762a1bSJed Brown PetscFunctionBegin; 15*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(A,&user)); 16*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(user->B,X)); 17c4762a1bSJed Brown PetscFunctionReturn(0); 18c4762a1bSJed Brown } 19c4762a1bSJed Brown 20c4762a1bSJed Brown static PetscErrorCode MatMult_User(Mat A,Vec X,Vec Y) 21c4762a1bSJed Brown { 22c4762a1bSJed Brown User user; 23c4762a1bSJed Brown 24c4762a1bSJed Brown PetscFunctionBegin; 25*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(A,&user)); 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(user->B,X,Y)); 27c4762a1bSJed Brown PetscFunctionReturn(0); 28c4762a1bSJed Brown } 29c4762a1bSJed Brown 30c4762a1bSJed Brown static PetscErrorCode MatMultTranspose_User(Mat A,Vec X,Vec Y) 31c4762a1bSJed Brown { 32c4762a1bSJed Brown User user; 33c4762a1bSJed Brown 34c4762a1bSJed Brown PetscFunctionBegin; 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(A,&user)); 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(user->B,X,Y)); 37c4762a1bSJed Brown PetscFunctionReturn(0); 38c4762a1bSJed Brown } 39c4762a1bSJed Brown 40c4762a1bSJed Brown static PetscErrorCode MatCopy_User(Mat A,Mat X,MatStructure str) 41c4762a1bSJed Brown { 42c4762a1bSJed Brown User user,userX; 43c4762a1bSJed Brown 44c4762a1bSJed Brown PetscFunctionBegin; 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(A,&user)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(X,&userX)); 472c71b3e2SJacob Faibussowitsch PetscCheckFalse(user != userX,PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"This should not happen"); 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)user->B)); 49c4762a1bSJed Brown PetscFunctionReturn(0); 50c4762a1bSJed Brown } 51c4762a1bSJed Brown 52c4762a1bSJed Brown static PetscErrorCode MatDestroy_User(Mat A) 53c4762a1bSJed Brown { 54c4762a1bSJed Brown User user; 55c4762a1bSJed Brown 56c4762a1bSJed Brown PetscFunctionBegin; 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(A,&user)); 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectDereference((PetscObject)user->B)); 59c4762a1bSJed Brown PetscFunctionReturn(0); 60c4762a1bSJed Brown } 61c4762a1bSJed Brown 62c4762a1bSJed Brown int main(int argc,char **args) 63c4762a1bSJed Brown { 64c4762a1bSJed Brown User user; 65c4762a1bSJed Brown Mat A,S; 66c4762a1bSJed Brown PetscScalar *data,diag = 1.3; 67c4762a1bSJed Brown PetscReal tol = PETSC_SMALL; 68c4762a1bSJed Brown PetscInt i,j,m = PETSC_DECIDE,n = PETSC_DECIDE,M = 17,N = 15,s1,s2; 69c4762a1bSJed Brown PetscInt test, ntest = 2; 70c4762a1bSJed Brown PetscMPIInt rank,size; 71c4762a1bSJed Brown PetscBool nc = PETSC_FALSE, cong, flg; 72c4762a1bSJed Brown PetscBool ronl = PETSC_TRUE; 73c4762a1bSJed Brown PetscBool randomize = PETSC_FALSE, submat = PETSC_FALSE; 74c4762a1bSJed Brown PetscBool keep = PETSC_FALSE; 75c4762a1bSJed Brown PetscBool testzerorows = PETSC_TRUE, testdiagscale = PETSC_TRUE, testgetdiag = PETSC_TRUE, testsubmat = PETSC_TRUE; 76c4762a1bSJed Brown PetscBool testshift = PETSC_TRUE, testscale = PETSC_TRUE, testdup = PETSC_TRUE, testreset = PETSC_TRUE; 77c4762a1bSJed Brown PetscBool testaxpy = PETSC_TRUE, testaxpyd = PETSC_TRUE, testaxpyerr = PETSC_FALSE; 78c4762a1bSJed Brown PetscErrorCode ierr; 79c4762a1bSJed Brown 80c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 81*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 82*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL)); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-ml",&m,NULL)); 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nl",&n,NULL)); 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-square_nc",&nc,NULL)); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-rows_only",&ronl,NULL)); 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-randomize",&randomize,NULL)); 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-submat",&submat,NULL)); 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_zerorows",&testzerorows,NULL)); 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_diagscale",&testdiagscale,NULL)); 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_getdiag",&testgetdiag,NULL)); 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_shift",&testshift,NULL)); 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_scale",&testscale,NULL)); 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_dup",&testdup,NULL)); 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_reset",&testreset,NULL)); 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_submat",&testsubmat,NULL)); 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_axpy",&testaxpy,NULL)); 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_axpy_different",&testaxpyd,NULL)); 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_axpy_error",&testaxpyerr,NULL)); 102*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-loop",&ntest,NULL)); 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-tol",&tol,NULL)); 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-diag",&diag,NULL)); 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-keep",&keep,NULL)); 106c4762a1bSJed Brown /* This tests square matrices with different row/col layout */ 107c4762a1bSJed Brown if (nc && size > 1) { 108c4762a1bSJed Brown M = PetscMax(PetscMax(N,M),1); 109c4762a1bSJed Brown N = M; 110c4762a1bSJed Brown m = n = 0; 111c4762a1bSJed Brown if (rank == 0) { m = M-1; n = 1; } 112c4762a1bSJed Brown else if (rank == 1) { m = 1; n = N-1; } 113c4762a1bSJed Brown } 114*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,m,n,M,N,NULL,&A)); 115*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A,&m,&n)); 116*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&M,&N)); 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&s1,NULL)); 118c4762a1bSJed Brown s2 = 1; 119c4762a1bSJed Brown while (s2 < M) s2 *= 10; 120*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArray(A,&data)); 121c4762a1bSJed Brown for (j = 0; j < N; j++) { 122c4762a1bSJed Brown for (i = 0; i < m; i++) { 123c4762a1bSJed Brown data[j*m + i] = s2*j + i + s1 + 1; 124c4762a1bSJed Brown } 125c4762a1bSJed Brown } 126*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 128c4762a1bSJed Brown 129c4762a1bSJed Brown if (submat) { 130c4762a1bSJed Brown Mat A2; 131c4762a1bSJed Brown IS r,c; 132c4762a1bSJed Brown PetscInt rst,ren,cst,cen; 133c4762a1bSJed Brown 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&rst,&ren)); 135*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRangeColumn(A,&cst,&cen)); 136*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PetscObjectComm((PetscObject)A),(ren-rst)/2,rst,1,&r)); 137*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PetscObjectComm((PetscObject)A),(cen-cst)/2,cst,1,&c)); 138*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrix(A,r,c,MAT_INITIAL_MATRIX,&A2)); 139*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&r)); 140*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&c)); 141*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 142c4762a1bSJed Brown A = A2; 143c4762a1bSJed Brown } 144c4762a1bSJed Brown 145*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&M,&N)); 146*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A,&m,&n)); 147*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatHasCongruentLayouts(A,&cong)); 148c4762a1bSJed Brown 149*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATAIJ,MAT_INPLACE_MATRIX,&A)); 150*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,keep)); 151*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)A,"initial")); 152*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(A,NULL,"-view_mat")); 153c4762a1bSJed Brown 154*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&user)); 155*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,m,n,M,N,user,&S)); 156*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(S,MATOP_MULT,(void (*)(void))MatMult_User)); 157*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(S,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMultTranspose_User)); 158c4762a1bSJed Brown if (cong) { 159*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(S,MATOP_GET_DIAGONAL,(void (*)(void))MatGetDiagonal_User)); 160c4762a1bSJed Brown } 161*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(S,MATOP_COPY,(void (*)(void))MatCopy_User)); 162*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(S,MATOP_DESTROY,(void (*)(void))MatDestroy_User)); 163*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&user->B)); 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* Square and rows only scaling */ 166c4762a1bSJed Brown ronl = cong ? ronl : PETSC_TRUE; 167c4762a1bSJed Brown 168c4762a1bSJed Brown for (test = 0; test < ntest; test++) { 169c4762a1bSJed Brown PetscReal err; 170c4762a1bSJed Brown 171*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAddEqual(A,S,10,&flg)); 172c4762a1bSJed Brown if (!flg) { 173*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"[test %" PetscInt_FMT "] Error mult add\n",test)); 174c4762a1bSJed Brown } 175*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAddEqual(A,S,10,&flg)); 176c4762a1bSJed Brown if (!flg) { 177*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"[test %" PetscInt_FMT "] Error mult add (T)\n",test)); 178c4762a1bSJed Brown } 179c4762a1bSJed Brown if (testzerorows) { 180c4762a1bSJed Brown Mat ST,B,C,BT,BTT; 181c4762a1bSJed Brown IS zr; 182c4762a1bSJed Brown Vec x = NULL, b1 = NULL, b2 = NULL; 183c4762a1bSJed Brown PetscInt *idxs = NULL, nr = 0; 184c4762a1bSJed Brown 185c4762a1bSJed Brown if (rank == (test%size)) { 186c4762a1bSJed Brown nr = 1; 187*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nr,&idxs)); 188c4762a1bSJed Brown if (test%2) { 189c4762a1bSJed Brown idxs[0] = (2*M - 1 - test/2)%M; 190c4762a1bSJed Brown } else { 191c4762a1bSJed Brown idxs[0] = (test/2)%M; 192c4762a1bSJed Brown } 193c4762a1bSJed Brown idxs[0] = PetscMax(idxs[0],0); 194c4762a1bSJed Brown } 195*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,nr,idxs,PETSC_OWN_POINTER,&zr)); 196*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)zr,"ZR")); 197*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISViewFromOptions(zr,NULL,"-view_is")); 198*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A,&x,&b1)); 199c4762a1bSJed Brown if (randomize) { 200*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(x,NULL)); 201*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(b1,NULL)); 202c4762a1bSJed Brown } else { 203*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(x,11.4)); 204*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(b1,-14.2)); 205c4762a1bSJed Brown } 206*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(b1,&b2)); 207*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(b1,b2)); 208*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)b1,"A_B1")); 209*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)b2,"A_B2")); 210c4762a1bSJed Brown if (size > 1 && !cong) { /* MATMPIAIJ ZeroRows and ZeroRowsColumns are buggy in this case */ 211*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b1)); 212c4762a1bSJed Brown } 213c4762a1bSJed Brown if (ronl) { 214*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroRowsIS(A,zr,diag,x,b1)); 215*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroRowsIS(S,zr,diag,x,b2)); 216c4762a1bSJed Brown } else { 217*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroRowsColumnsIS(A,zr,diag,x,b1)); 218*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroRowsColumnsIS(S,zr,diag,x,b2)); 219*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&zr)); 220c4762a1bSJed Brown /* Mix zerorows and zerorowscols */ 221c4762a1bSJed Brown nr = 0; 222c4762a1bSJed Brown idxs = NULL; 223dd400576SPatrick Sanan if (rank == 0) { 224c4762a1bSJed Brown nr = 1; 225*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nr,&idxs)); 226c4762a1bSJed Brown if (test%2) { 227c4762a1bSJed Brown idxs[0] = (3*M - 2 - test/2)%M; 228c4762a1bSJed Brown } else { 229c4762a1bSJed Brown idxs[0] = (test/2+1)%M; 230c4762a1bSJed Brown } 231c4762a1bSJed Brown idxs[0] = PetscMax(idxs[0],0); 232c4762a1bSJed Brown } 233*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,nr,idxs,PETSC_OWN_POINTER,&zr)); 234*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)zr,"ZR2")); 235*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISViewFromOptions(zr,NULL,"-view_is")); 236*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroRowsIS(A,zr,diag*2.0+PETSC_SMALL,NULL,NULL)); 237*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroRowsIS(S,zr,diag*2.0+PETSC_SMALL,NULL,NULL)); 238c4762a1bSJed Brown } 239*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&zr)); 240c4762a1bSJed Brown 241c4762a1bSJed Brown if (b1) { 242c4762a1bSJed Brown Vec b; 243c4762a1bSJed Brown 244*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(b1,NULL,"-view_b")); 245*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(b2,NULL,"-view_b")); 246*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(b1,&b)); 247*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(b1,b)); 248*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(b,-1.0,b2)); 249*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(b,NORM_INFINITY,&err)); 250c4762a1bSJed Brown if (err >= tol) { 251*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"[test %" PetscInt_FMT "] Error b %g\n",test,(double)err)); 252c4762a1bSJed Brown } 253*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b)); 254c4762a1bSJed Brown } 255*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b1)); 256*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b2)); 257*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 258*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(S,MATDENSE,MAT_INITIAL_MATRIX,&B)); 259c4762a1bSJed Brown 260*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateTranspose(S,&ST)); 261*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatComputeOperator(ST,MATDENSE,&BT)); 262*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(BT,MAT_INITIAL_MATRIX,&BTT)); 263*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)B,"S")); 264*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)BTT,"STT")); 265*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&C)); 266*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)C,"A")); 267c4762a1bSJed Brown 268*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(C,NULL,"-view_mat")); 269*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(B,NULL,"-view_mat")); 270*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(BTT,NULL,"-view_mat")); 271c4762a1bSJed Brown 272*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(C,-1.0,B,SAME_NONZERO_PATTERN)); 273*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(C,NORM_FROBENIUS,&err)); 274c4762a1bSJed Brown if (err >= tol) { 275*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"[test %" PetscInt_FMT "] Error mat mult after %s %g\n",test,ronl ? "MatZeroRows" : "MatZeroRowsColumns",(double)err)); 276c4762a1bSJed Brown } 277c4762a1bSJed Brown 278*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATDENSE,MAT_REUSE_MATRIX,&C)); 279*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(C,-1.0,BTT,SAME_NONZERO_PATTERN)); 280*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(C,NORM_FROBENIUS,&err)); 281c4762a1bSJed Brown if (err >= tol) { 282*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"[test %" PetscInt_FMT "] Error mat mult transpose after %s %g\n",test,ronl ? "MatZeroRows" : "MatZeroRowsColumns",(double)err)); 283c4762a1bSJed Brown } 284c4762a1bSJed Brown 285*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&ST)); 286*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&BTT)); 287*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&BT)); 288*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 289*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 290c4762a1bSJed Brown } 291c4762a1bSJed Brown if (testdiagscale) { /* MatDiagonalScale() */ 292c4762a1bSJed Brown Vec vr,vl; 293c4762a1bSJed Brown 294*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A,&vr,&vl)); 295c4762a1bSJed Brown if (randomize) { 296*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(vr,NULL)); 297*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(vl,NULL)); 298c4762a1bSJed Brown } else { 299*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(vr,test%2 ? 0.15 : 1.0/0.15)); 300*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(vl,test%2 ? -1.2 : 1.0/-1.2)); 301c4762a1bSJed Brown } 302*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(A,vl,vr)); 303*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(S,vl,vr)); 304*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&vr)); 305*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&vl)); 306c4762a1bSJed Brown } 307c4762a1bSJed Brown 308c4762a1bSJed Brown if (testscale) { /* MatScale() */ 309*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,test%2 ? 1.4 : 1.0/1.4)); 310*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(S,test%2 ? 1.4 : 1.0/1.4)); 311c4762a1bSJed Brown } 312c4762a1bSJed Brown 313c4762a1bSJed Brown if (testshift && cong) { /* MatShift() : MATSHELL shift is broken when row/cols layout are not congruent and left/right scaling have been applied */ 314*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShift(A,test%2 ? -77.5 : 77.5)); 315*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShift(S,test%2 ? -77.5 : 77.5)); 316c4762a1bSJed Brown } 317c4762a1bSJed Brown 318c4762a1bSJed Brown if (testgetdiag && cong) { /* MatGetDiagonal() */ 319c4762a1bSJed Brown Vec dA,dS; 320c4762a1bSJed Brown 321*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A,&dA,NULL)); 322*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(S,&dS,NULL)); 323*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(A,dA)); 324*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(S,dS)); 325*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(dA,-1.0,dS)); 326*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(dA,NORM_INFINITY,&err)); 327c4762a1bSJed Brown if (err >= tol) { 328*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"[test %" PetscInt_FMT "] Error diag %g\n",test,(double)err)); 329c4762a1bSJed Brown } 330*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&dA)); 331*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&dS)); 332c4762a1bSJed Brown } 333c4762a1bSJed Brown 334c4762a1bSJed Brown if (testdup && !test) { 335c4762a1bSJed Brown Mat A2, S2; 336c4762a1bSJed Brown 337*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&A2)); 338*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(S,MAT_COPY_VALUES,&S2)); 339*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 340*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&S)); 341c4762a1bSJed Brown A = A2; 342c4762a1bSJed Brown S = S2; 343c4762a1bSJed Brown } 344c4762a1bSJed Brown 345c4762a1bSJed Brown if (testsubmat) { 346c4762a1bSJed Brown Mat sA,sS,dA,dS,At,St; 347c4762a1bSJed Brown IS r,c; 348c4762a1bSJed Brown PetscInt rst,ren,cst,cen; 349c4762a1bSJed Brown 350*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&rst,&ren)); 351*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRangeColumn(A,&cst,&cen)); 352*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PetscObjectComm((PetscObject)A),(ren-rst)/2,rst,1,&r)); 353*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PetscObjectComm((PetscObject)A),(cen-cst)/2,cst,1,&c)); 354*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrix(A,r,c,MAT_INITIAL_MATRIX,&sA)); 355*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrix(S,r,c,MAT_INITIAL_MATRIX,&sS)); 356*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAddEqual(sA,sS,10,&flg)); 357c4762a1bSJed Brown if (!flg) { 358*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"[test %" PetscInt_FMT "] Error submatrix mult add\n",test)); 359c4762a1bSJed Brown } 360*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAddEqual(sA,sS,10,&flg)); 361c4762a1bSJed Brown if (!flg) { 362*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"[test %" PetscInt_FMT "] Error submatrix mult add (T)\n",test)); 363c4762a1bSJed Brown } 364*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(sA,MATDENSE,MAT_INITIAL_MATRIX,&dA)); 365*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(sS,MATDENSE,MAT_INITIAL_MATRIX,&dS)); 366*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(dA,-1.0,dS,SAME_NONZERO_PATTERN)); 367*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(dA,NORM_FROBENIUS,&err)); 368c4762a1bSJed Brown if (err >= tol) { 369*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"[test %" PetscInt_FMT "] Error mat submatrix %g\n",test,(double)err)); 370c4762a1bSJed Brown } 371*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&sA)); 372*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&sS)); 373*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&dA)); 374*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&dS)); 375*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateTranspose(A,&At)); 376*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateTranspose(S,&St)); 377*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrix(At,c,r,MAT_INITIAL_MATRIX,&sA)); 378*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrix(St,c,r,MAT_INITIAL_MATRIX,&sS)); 379*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAddEqual(sA,sS,10,&flg)); 380c4762a1bSJed Brown if (!flg) { 381*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"[test %" PetscInt_FMT "] Error submatrix (T) mult add\n",test)); 382c4762a1bSJed Brown } 383*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAddEqual(sA,sS,10,&flg)); 384c4762a1bSJed Brown if (!flg) { 385*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"[test %" PetscInt_FMT "] Error submatrix (T) mult add (T)\n",test)); 386c4762a1bSJed Brown } 387*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(sA,MATDENSE,MAT_INITIAL_MATRIX,&dA)); 388*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(sS,MATDENSE,MAT_INITIAL_MATRIX,&dS)); 389*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(dA,-1.0,dS,SAME_NONZERO_PATTERN)); 390*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(dA,NORM_FROBENIUS,&err)); 391c4762a1bSJed Brown if (err >= tol) { 392*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"[test %" PetscInt_FMT "] Error mat submatrix (T) %g\n",test,(double)err)); 393c4762a1bSJed Brown } 394*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&sA)); 395*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&sS)); 396*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&dA)); 397*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&dS)); 398*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&At)); 399*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&St)); 400*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&r)); 401*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&c)); 402c4762a1bSJed Brown } 403c4762a1bSJed Brown 404c4762a1bSJed Brown if (testaxpy) { 405c4762a1bSJed Brown Mat tA,tS,dA,dS; 406c4762a1bSJed Brown MatStructure str[3] = { SAME_NONZERO_PATTERN, SUBSET_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN }; 407c4762a1bSJed Brown 408*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&tA)); 409c4762a1bSJed Brown if (testaxpyd && !(test%2)) { 410*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)tA)); 411c4762a1bSJed Brown tS = tA; 412c4762a1bSJed Brown } else { 413*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)S)); 414c4762a1bSJed Brown tS = S; 415c4762a1bSJed Brown } 416*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(A,0.5,tA,str[test%3])); 417*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(S,0.5,tS,str[test%3])); 418c4762a1bSJed Brown /* this will trigger an error the next MatMult or MatMultTranspose call for S */ 419*5f80ce2aSJacob Faibussowitsch if (testaxpyerr) CHKERRQ(MatScale(tA,0)); 420*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&tA)); 421*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&tS)); 422*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAddEqual(A,S,10,&flg)); 423c4762a1bSJed Brown if (!flg) { 424*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"[test %" PetscInt_FMT "] Error axpy mult add\n",test)); 425c4762a1bSJed Brown } 426*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAddEqual(A,S,10,&flg)); 427c4762a1bSJed Brown if (!flg) { 428*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"[test %" PetscInt_FMT "] Error axpy mult add (T)\n",test)); 429c4762a1bSJed Brown } 430*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&dA)); 431*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(S,MATDENSE,MAT_INITIAL_MATRIX,&dS)); 432*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(dA,-1.0,dS,SAME_NONZERO_PATTERN)); 433*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(dA,NORM_FROBENIUS,&err)); 434c4762a1bSJed Brown if (err >= tol) { 435*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"[test %" PetscInt_FMT "] Error mat submatrix %g\n",test,(double)err)); 436c4762a1bSJed Brown } 437*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&dA)); 438*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&dS)); 439c4762a1bSJed Brown } 440c4762a1bSJed Brown 441c4762a1bSJed Brown if (testreset && (ntest == 1 || test == ntest-2)) { 442c4762a1bSJed Brown /* reset MATSHELL */ 443*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(S,MAT_FINAL_ASSEMBLY)); 444*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(S,MAT_FINAL_ASSEMBLY)); 445c4762a1bSJed Brown /* reset A */ 446*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCopy(user->B,A,DIFFERENT_NONZERO_PATTERN)); 447c4762a1bSJed Brown } 448c4762a1bSJed Brown } 449c4762a1bSJed Brown 450*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 451*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&S)); 452*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(user)); 453c4762a1bSJed Brown ierr = PetscFinalize(); 454c4762a1bSJed Brown return ierr; 455c4762a1bSJed Brown } 456c4762a1bSJed Brown 457c4762a1bSJed Brown /*TEST 458c4762a1bSJed Brown 459c4762a1bSJed Brown testset: 460c4762a1bSJed Brown suffix: rect 461c4762a1bSJed Brown requires: !single 462c4762a1bSJed Brown output_file: output/ex221_1.out 463c4762a1bSJed Brown nsize: {{1 3}} 464c4762a1bSJed Brown args: -loop 3 -keep {{0 1}} -M {{12 19}} -N {{19 12}} -submat {{0 1}} -test_axpy_different {{0 1}} 465c4762a1bSJed Brown 466c4762a1bSJed Brown testset: 467c4762a1bSJed Brown suffix: square 468c4762a1bSJed Brown requires: !single 469c4762a1bSJed Brown output_file: output/ex221_1.out 470c4762a1bSJed Brown nsize: {{1 3}} 471c4762a1bSJed Brown args: -M 21 -N 21 -loop 4 -rows_only {{0 1}} -keep {{0 1}} -submat {{0 1}} -test_axpy_different {{0 1}} 472c4762a1bSJed Brown TEST*/ 473