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