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