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