1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Tests MatHYPRE\n"; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown #include <petscmathypre.h> 5*c4762a1bSJed Brown 6*c4762a1bSJed Brown int main(int argc,char **args) 7*c4762a1bSJed Brown { 8*c4762a1bSJed Brown Mat A,B,C,D; 9*c4762a1bSJed Brown Mat pAB,CD,CAB; 10*c4762a1bSJed Brown hypre_ParCSRMatrix *parcsr; 11*c4762a1bSJed Brown PetscReal err; 12*c4762a1bSJed Brown PetscInt i,j,N = 6, M = 6; 13*c4762a1bSJed Brown PetscErrorCode ierr; 14*c4762a1bSJed Brown PetscBool flg,testptap = PETSC_TRUE,testmatmatmult = PETSC_TRUE; 15*c4762a1bSJed Brown PetscReal norm; 16*c4762a1bSJed Brown char file[256]; 17*c4762a1bSJed Brown 18*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 19*c4762a1bSJed Brown ierr = PetscOptionsGetString(NULL,NULL,"-f",file,256,&flg);CHKERRQ(ierr); 20*c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 21*c4762a1bSJed Brown testptap = PETSC_FALSE; 22*c4762a1bSJed Brown testmatmatmult = PETSC_FALSE; 23*c4762a1bSJed Brown ierr = PetscOptionsInsertString(NULL,"-options_left 0");CHKERRQ(ierr); 24*c4762a1bSJed Brown #endif 25*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-ptap",&testptap,NULL);CHKERRQ(ierr); 26*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-matmatmult",&testmatmatmult,NULL);CHKERRQ(ierr); 27*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 28*c4762a1bSJed Brown if (!flg) { /* Create a matrix and test MatSetValues */ 29*c4762a1bSJed Brown PetscMPIInt size; 30*c4762a1bSJed Brown 31*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 32*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr); 33*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);CHKERRQ(ierr); 34*c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 35*c4762a1bSJed Brown ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr); 36*c4762a1bSJed Brown ierr = MatSeqAIJSetPreallocation(A,9,NULL);CHKERRQ(ierr); 37*c4762a1bSJed Brown ierr = MatMPIAIJSetPreallocation(A,9,NULL,9,NULL);CHKERRQ(ierr); 38*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr); 39*c4762a1bSJed Brown ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = MatSetType(B,MATHYPRE);CHKERRQ(ierr); 41*c4762a1bSJed Brown if (M == N) { 42*c4762a1bSJed Brown ierr = MatHYPRESetPreallocation(B,9,NULL,9,NULL);CHKERRQ(ierr); 43*c4762a1bSJed Brown } else { 44*c4762a1bSJed Brown ierr = MatHYPRESetPreallocation(B,6,NULL,6,NULL);CHKERRQ(ierr); 45*c4762a1bSJed Brown } 46*c4762a1bSJed Brown if (M == N) { 47*c4762a1bSJed Brown for (i=0; i<M; i++) { 48*c4762a1bSJed Brown PetscInt cols[] = {0,1,2,3,4,5}; 49*c4762a1bSJed Brown PetscScalar vals[] = {0,1./size,2./size,3./size,4./size,5./size}; 50*c4762a1bSJed Brown for (j=i-2; j<i+1; j++) { 51*c4762a1bSJed Brown if (j >= N) { 52*c4762a1bSJed Brown ierr = MatSetValue(A,i,N-1,(1.*j*N+i)/(3.*N*size),ADD_VALUES);CHKERRQ(ierr); 53*c4762a1bSJed Brown ierr = MatSetValue(B,i,N-1,(1.*j*N+i)/(3.*N*size),ADD_VALUES);CHKERRQ(ierr); 54*c4762a1bSJed Brown } else if (i > j) { 55*c4762a1bSJed Brown ierr = MatSetValue(A,i,PetscMin(j,N-1),(1.*j*N+i)/(2.*N*size),ADD_VALUES);CHKERRQ(ierr); 56*c4762a1bSJed Brown ierr = MatSetValue(B,i,PetscMin(j,N-1),(1.*j*N+i)/(2.*N*size),ADD_VALUES);CHKERRQ(ierr); 57*c4762a1bSJed Brown } else { 58*c4762a1bSJed Brown ierr = MatSetValue(A,i,PetscMin(j,N-1),-1.-(1.*j*N+i)/(4.*N*size),ADD_VALUES);CHKERRQ(ierr); 59*c4762a1bSJed Brown ierr = MatSetValue(B,i,PetscMin(j,N-1),-1.-(1.*j*N+i)/(4.*N*size),ADD_VALUES);CHKERRQ(ierr); 60*c4762a1bSJed Brown } 61*c4762a1bSJed Brown } 62*c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,PetscMin(6,N),cols,vals,ADD_VALUES);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = MatSetValues(B,1,&i,PetscMin(6,N),cols,vals,ADD_VALUES);CHKERRQ(ierr); 64*c4762a1bSJed Brown } 65*c4762a1bSJed Brown } else { 66*c4762a1bSJed Brown PetscInt rows[2]; 67*c4762a1bSJed Brown PetscBool test_offproc = PETSC_FALSE; 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-test_offproc",&test_offproc,NULL); 70*c4762a1bSJed Brown if (test_offproc) { 71*c4762a1bSJed Brown const PetscInt *ranges; 72*c4762a1bSJed Brown PetscMPIInt rank; 73*c4762a1bSJed Brown 74*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 75*c4762a1bSJed Brown ierr = MatGetOwnershipRanges(A,&ranges);CHKERRQ(ierr); 76*c4762a1bSJed Brown rows[0] = ranges[(rank+1)%size]; 77*c4762a1bSJed Brown rows[1] = ranges[(rank+1)%size + 1]; 78*c4762a1bSJed Brown } else { 79*c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&rows[0],&rows[1]);CHKERRQ(ierr); 80*c4762a1bSJed Brown } 81*c4762a1bSJed Brown for (i=rows[0];i<rows[1];i++) { 82*c4762a1bSJed Brown PetscInt cols[] = {0,1,2,3,4,5}; 83*c4762a1bSJed Brown PetscScalar vals[] = {-1,1,-2,2,-3,3}; 84*c4762a1bSJed Brown 85*c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,PetscMin(6,N),cols,vals,INSERT_VALUES);CHKERRQ(ierr); 86*c4762a1bSJed Brown ierr = MatSetValues(B,1,&i,PetscMin(6,N),cols,vals,INSERT_VALUES);CHKERRQ(ierr); 87*c4762a1bSJed Brown } 88*c4762a1bSJed Brown } 89*c4762a1bSJed Brown /* MAT_FLUSH_ASSEMBLY currently not supported */ 90*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 91*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 92*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 93*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 94*c4762a1bSJed Brown 95*c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 96*c4762a1bSJed Brown /* make the matrix imaginary */ 97*c4762a1bSJed Brown ierr = MatScale(A,PETSC_i);CHKERRQ(ierr); 98*c4762a1bSJed Brown ierr = MatScale(B,PETSC_i);CHKERRQ(ierr); 99*c4762a1bSJed Brown #endif 100*c4762a1bSJed Brown 101*c4762a1bSJed Brown /* MatAXPY further exercises MatSetValues_HYPRE */ 102*c4762a1bSJed Brown ierr = MatAXPY(B,-1.,A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 103*c4762a1bSJed Brown ierr = MatConvert(B,MATMPIAIJ,MAT_INITIAL_MATRIX,&C);CHKERRQ(ierr); 104*c4762a1bSJed Brown ierr = MatNorm(C,NORM_INFINITY,&err);CHKERRQ(ierr); 105*c4762a1bSJed Brown if (err > PETSC_SMALL) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error MatSetValues %g",err); 106*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 107*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 108*c4762a1bSJed Brown } else { 109*c4762a1bSJed Brown PetscViewer viewer; 110*c4762a1bSJed Brown 111*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&viewer);CHKERRQ(ierr); 112*c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 113*c4762a1bSJed Brown ierr = MatLoad(A,viewer);CHKERRQ(ierr); 114*c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 115*c4762a1bSJed Brown ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 116*c4762a1bSJed Brown } 117*c4762a1bSJed Brown 118*c4762a1bSJed Brown /* check conversion routines */ 119*c4762a1bSJed Brown ierr = MatConvert(A,MATHYPRE,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 120*c4762a1bSJed Brown ierr = MatConvert(A,MATHYPRE,MAT_REUSE_MATRIX,&B);CHKERRQ(ierr); 121*c4762a1bSJed Brown ierr = MatConvert(B,MATIS,MAT_INITIAL_MATRIX,&D);CHKERRQ(ierr); 122*c4762a1bSJed Brown ierr = MatConvert(B,MATIS,MAT_REUSE_MATRIX,&D);CHKERRQ(ierr); 123*c4762a1bSJed Brown ierr = MatConvert(B,MATAIJ,MAT_INITIAL_MATRIX,&C);CHKERRQ(ierr); 124*c4762a1bSJed Brown ierr = MatConvert(B,MATAIJ,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr); 125*c4762a1bSJed Brown ierr = MatAXPY(C,-1.,A,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 126*c4762a1bSJed Brown ierr = MatNorm(C,NORM_INFINITY,&err);CHKERRQ(ierr); 127*c4762a1bSJed Brown if (err > PETSC_SMALL) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error Mat AIJ %g",err); 128*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 129*c4762a1bSJed Brown ierr = MatConvert(D,MATAIJ,MAT_INITIAL_MATRIX,&C);CHKERRQ(ierr); 130*c4762a1bSJed Brown ierr = MatAXPY(C,-1.,A,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 131*c4762a1bSJed Brown ierr = MatNorm(C,NORM_INFINITY,&err);CHKERRQ(ierr); 132*c4762a1bSJed Brown if (err > PETSC_SMALL) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error Mat IS %g",err); 133*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 134*c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 135*c4762a1bSJed Brown 136*c4762a1bSJed Brown /* check MatCreateFromParCSR */ 137*c4762a1bSJed Brown ierr = MatHYPREGetParCSR(B,&parcsr);CHKERRQ(ierr); 138*c4762a1bSJed Brown ierr = MatCreateFromParCSR(parcsr,MATAIJ,PETSC_COPY_VALUES,&D);CHKERRQ(ierr); 139*c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 140*c4762a1bSJed Brown ierr = MatCreateFromParCSR(parcsr,MATHYPRE,PETSC_USE_POINTER,&C);CHKERRQ(ierr); 141*c4762a1bSJed Brown 142*c4762a1bSJed Brown /* check MatMult operations */ 143*c4762a1bSJed Brown ierr = MatMultEqual(A,B,4,&flg);CHKERRQ(ierr); 144*c4762a1bSJed Brown if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error MatMult B"); 145*c4762a1bSJed Brown ierr = MatMultEqual(A,C,4,&flg);CHKERRQ(ierr); 146*c4762a1bSJed Brown if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error MatMult C"); 147*c4762a1bSJed Brown ierr = MatMultAddEqual(A,B,4,&flg);CHKERRQ(ierr); 148*c4762a1bSJed Brown if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error MatMultAdd B"); 149*c4762a1bSJed Brown ierr = MatMultAddEqual(A,C,4,&flg);CHKERRQ(ierr); 150*c4762a1bSJed Brown if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error MatMultAdd C"); 151*c4762a1bSJed Brown ierr = MatMultTransposeEqual(A,B,4,&flg);CHKERRQ(ierr); 152*c4762a1bSJed Brown if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error MatMultTranspose B"); 153*c4762a1bSJed Brown ierr = MatMultTransposeEqual(A,C,4,&flg);CHKERRQ(ierr); 154*c4762a1bSJed Brown if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error MatMultTranspose C"); 155*c4762a1bSJed Brown ierr = MatMultTransposeAddEqual(A,B,4,&flg);CHKERRQ(ierr); 156*c4762a1bSJed Brown if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error MatMultTransposeAdd B"); 157*c4762a1bSJed Brown ierr = MatMultTransposeAddEqual(A,C,4,&flg);CHKERRQ(ierr); 158*c4762a1bSJed Brown if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error MatMultTransposeAdd C"); 159*c4762a1bSJed Brown 160*c4762a1bSJed Brown /* check PtAP */ 161*c4762a1bSJed Brown if (testptap && M == N) { 162*c4762a1bSJed Brown Mat pP,hP; 163*c4762a1bSJed Brown 164*c4762a1bSJed Brown /* PETSc MatPtAP -> output is a MatAIJ 165*c4762a1bSJed Brown It uses HYPRE functions when -matptap_via hypre is specified at command line */ 166*c4762a1bSJed Brown ierr = MatPtAP(A,A,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&pP);CHKERRQ(ierr); 167*c4762a1bSJed Brown ierr = MatPtAP(A,A,MAT_REUSE_MATRIX,PETSC_DEFAULT,&pP);CHKERRQ(ierr); 168*c4762a1bSJed Brown ierr = MatNorm(pP,NORM_INFINITY,&norm);CHKERRQ(ierr); 169*c4762a1bSJed Brown 170*c4762a1bSJed Brown /* MatPtAP_HYPRE_HYPRE -> output is a MatHYPRE */ 171*c4762a1bSJed Brown ierr = MatPtAP(C,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&hP);CHKERRQ(ierr); 172*c4762a1bSJed Brown ierr = MatPtAP(C,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&hP);CHKERRQ(ierr); 173*c4762a1bSJed Brown ierr = MatAXPY(hP,-1.,pP,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 174*c4762a1bSJed Brown ierr = MatHasOperation(hP,MATOP_NORM,&flg);CHKERRQ(ierr); 175*c4762a1bSJed Brown if (!flg) { /* TODO add MatNorm_HYPRE */ 176*c4762a1bSJed Brown ierr = MatConvert(hP,MATAIJ,MAT_INPLACE_MATRIX,&hP);CHKERRQ(ierr); 177*c4762a1bSJed Brown } 178*c4762a1bSJed Brown ierr = MatNorm(hP,NORM_INFINITY,&err);CHKERRQ(ierr); 179*c4762a1bSJed Brown if (err/norm > PETSC_SMALL) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error MatPtAP %g %g",err,norm); 180*c4762a1bSJed Brown ierr = MatDestroy(&hP);CHKERRQ(ierr); 181*c4762a1bSJed Brown 182*c4762a1bSJed Brown /* MatPtAP_AIJ_HYPRE -> output can be decided at runtime with -matptap_hypre_outtype */ 183*c4762a1bSJed Brown ierr = MatPtAP(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&hP);CHKERRQ(ierr); 184*c4762a1bSJed Brown ierr = MatPtAP(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&hP);CHKERRQ(ierr); 185*c4762a1bSJed Brown ierr = MatAXPY(hP,-1.,pP,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 186*c4762a1bSJed Brown ierr = MatHasOperation(hP,MATOP_NORM,&flg);CHKERRQ(ierr); 187*c4762a1bSJed Brown if (!flg) { /* TODO add MatNorm_HYPRE */ 188*c4762a1bSJed Brown ierr = MatConvert(hP,MATAIJ,MAT_INPLACE_MATRIX,&hP);CHKERRQ(ierr); 189*c4762a1bSJed Brown } 190*c4762a1bSJed Brown ierr = MatNorm(hP,NORM_INFINITY,&err);CHKERRQ(ierr); 191*c4762a1bSJed Brown if (err/norm > PETSC_SMALL) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error MatPtAP mixed %g %g",err,norm); 192*c4762a1bSJed Brown ierr = MatDestroy(&hP);CHKERRQ(ierr); 193*c4762a1bSJed Brown 194*c4762a1bSJed Brown ierr = MatDestroy(&pP);CHKERRQ(ierr); 195*c4762a1bSJed Brown } 196*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 197*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 198*c4762a1bSJed Brown 199*c4762a1bSJed Brown /* check MatMatMult */ 200*c4762a1bSJed Brown if (testmatmatmult) { 201*c4762a1bSJed Brown ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 202*c4762a1bSJed Brown ierr = MatConvert(A,MATHYPRE,MAT_INITIAL_MATRIX,&C);CHKERRQ(ierr); 203*c4762a1bSJed Brown ierr = MatConvert(B,MATHYPRE,MAT_INITIAL_MATRIX,&D);CHKERRQ(ierr); 204*c4762a1bSJed Brown 205*c4762a1bSJed Brown /* PETSc MatMatMult -> output is a MatAIJ 206*c4762a1bSJed Brown It uses HYPRE functions when -matmatmult_via hypre is specified at command line */ 207*c4762a1bSJed Brown ierr = MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&pAB);CHKERRQ(ierr); 208*c4762a1bSJed Brown ierr = MatMatMult(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&pAB);CHKERRQ(ierr); 209*c4762a1bSJed Brown ierr = MatNorm(pAB,NORM_INFINITY,&norm);CHKERRQ(ierr); 210*c4762a1bSJed Brown 211*c4762a1bSJed Brown /* MatMatMult_HYPRE_HYPRE -> output is a MatHYPRE */ 212*c4762a1bSJed Brown ierr = MatMatMult(C,D,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&CD);CHKERRQ(ierr); 213*c4762a1bSJed Brown ierr = MatMatMult(C,D,MAT_REUSE_MATRIX,PETSC_DEFAULT,&CD);CHKERRQ(ierr); 214*c4762a1bSJed Brown ierr = MatAXPY(CD,-1.,pAB,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 215*c4762a1bSJed Brown ierr = MatHasOperation(CD,MATOP_NORM,&flg);CHKERRQ(ierr); 216*c4762a1bSJed Brown if (!flg) { /* TODO add MatNorm_HYPRE */ 217*c4762a1bSJed Brown ierr = MatConvert(CD,MATAIJ,MAT_INPLACE_MATRIX,&CD);CHKERRQ(ierr); 218*c4762a1bSJed Brown } 219*c4762a1bSJed Brown ierr = MatNorm(CD,NORM_INFINITY,&err);CHKERRQ(ierr); 220*c4762a1bSJed Brown if (err/norm > PETSC_SMALL) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error MatMatMult %g %g",err,norm); 221*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 222*c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 223*c4762a1bSJed Brown ierr = MatDestroy(&CD);CHKERRQ(ierr); 224*c4762a1bSJed Brown ierr = MatDestroy(&pAB);CHKERRQ(ierr); 225*c4762a1bSJed Brown 226*c4762a1bSJed Brown /* When configured with HYPRE, MatMatMatMult is available for the triplet transpose(aij)-aij-aij */ 227*c4762a1bSJed Brown ierr = MatCreateTranspose(A,&C);CHKERRQ(ierr); 228*c4762a1bSJed Brown ierr = MatMatMatMult(C,A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&CAB);CHKERRQ(ierr); 229*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 230*c4762a1bSJed Brown ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&C);CHKERRQ(ierr); 231*c4762a1bSJed Brown ierr = MatMatMult(C,A,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&D);CHKERRQ(ierr); 232*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 233*c4762a1bSJed Brown ierr = MatMatMult(D,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr); 234*c4762a1bSJed Brown ierr = MatNorm(C,NORM_INFINITY,&norm);CHKERRQ(ierr); 235*c4762a1bSJed Brown ierr = MatAXPY(C,-1.,CAB,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 236*c4762a1bSJed Brown ierr = MatNorm(C,NORM_INFINITY,&err);CHKERRQ(ierr); 237*c4762a1bSJed Brown if (err/norm > PETSC_SMALL) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error MatMatMatMult %g %g",err,norm); 238*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 239*c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 240*c4762a1bSJed Brown ierr = MatDestroy(&CAB);CHKERRQ(ierr); 241*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 242*c4762a1bSJed Brown } 243*c4762a1bSJed Brown 244*c4762a1bSJed Brown /* Check MatView */ 245*c4762a1bSJed Brown ierr = MatViewFromOptions(A,NULL,"-view_A");CHKERRQ(ierr); 246*c4762a1bSJed Brown ierr = MatConvert(A,MATHYPRE,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 247*c4762a1bSJed Brown ierr = MatViewFromOptions(B,NULL,"-view_B");CHKERRQ(ierr); 248*c4762a1bSJed Brown 249*c4762a1bSJed Brown /* Check MatDuplicate/MatCopy */ 250*c4762a1bSJed Brown for (j=0;j<3;j++) { 251*c4762a1bSJed Brown MatDuplicateOption dop; 252*c4762a1bSJed Brown 253*c4762a1bSJed Brown dop = MAT_COPY_VALUES; 254*c4762a1bSJed Brown if (j==1) dop = MAT_DO_NOT_COPY_VALUES; 255*c4762a1bSJed Brown if (j==2) dop = MAT_SHARE_NONZERO_PATTERN; 256*c4762a1bSJed Brown 257*c4762a1bSJed Brown for (i=0;i<3;i++) { 258*c4762a1bSJed Brown MatStructure str; 259*c4762a1bSJed Brown 260*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Dup/Copy tests: %D %D\n",j,i);CHKERRQ(ierr); 261*c4762a1bSJed Brown 262*c4762a1bSJed Brown str = DIFFERENT_NONZERO_PATTERN; 263*c4762a1bSJed Brown if (i==1) str = SAME_NONZERO_PATTERN; 264*c4762a1bSJed Brown if (i==2) str = SUBSET_NONZERO_PATTERN; 265*c4762a1bSJed Brown 266*c4762a1bSJed Brown ierr = MatDuplicate(A,dop,&C);CHKERRQ(ierr); 267*c4762a1bSJed Brown ierr = MatDuplicate(B,dop,&D);CHKERRQ(ierr); 268*c4762a1bSJed Brown if (dop != MAT_COPY_VALUES) { 269*c4762a1bSJed Brown ierr = MatCopy(A,C,str);CHKERRQ(ierr); 270*c4762a1bSJed Brown ierr = MatCopy(B,D,str);CHKERRQ(ierr); 271*c4762a1bSJed Brown } 272*c4762a1bSJed Brown /* AXPY with AIJ and HYPRE */ 273*c4762a1bSJed Brown ierr = MatAXPY(C,-1.0,D,str);CHKERRQ(ierr); 274*c4762a1bSJed Brown ierr = MatNorm(C,NORM_INFINITY,&err);CHKERRQ(ierr); 275*c4762a1bSJed Brown if (err > PETSC_SMALL) { 276*c4762a1bSJed Brown ierr = MatViewFromOptions(A,NULL,"-view_duplicate_diff");CHKERRQ(ierr); 277*c4762a1bSJed Brown ierr = MatViewFromOptions(B,NULL,"-view_duplicate_diff");CHKERRQ(ierr); 278*c4762a1bSJed Brown ierr = MatViewFromOptions(C,NULL,"-view_duplicate_diff");CHKERRQ(ierr); 279*c4762a1bSJed Brown SETERRQ3(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error test 1 MatDuplicate/MatCopy %g (%D,%D)",err,j,i); 280*c4762a1bSJed Brown } 281*c4762a1bSJed Brown /* AXPY with HYPRE and HYPRE */ 282*c4762a1bSJed Brown ierr = MatAXPY(D,-1.0,B,str);CHKERRQ(ierr); 283*c4762a1bSJed Brown if (err > PETSC_SMALL) { 284*c4762a1bSJed Brown ierr = MatViewFromOptions(A,NULL,"-view_duplicate_diff");CHKERRQ(ierr); 285*c4762a1bSJed Brown ierr = MatViewFromOptions(B,NULL,"-view_duplicate_diff");CHKERRQ(ierr); 286*c4762a1bSJed Brown ierr = MatViewFromOptions(D,NULL,"-view_duplicate_diff");CHKERRQ(ierr); 287*c4762a1bSJed Brown SETERRQ3(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error test 2 MatDuplicate/MatCopy %g (%D,%D)",err,j,i); 288*c4762a1bSJed Brown } 289*c4762a1bSJed Brown /* Copy from HYPRE to AIJ */ 290*c4762a1bSJed Brown ierr = MatCopy(B,C,str);CHKERRQ(ierr); 291*c4762a1bSJed Brown /* Copy from AIJ to HYPRE */ 292*c4762a1bSJed Brown ierr = MatCopy(A,D,str);CHKERRQ(ierr); 293*c4762a1bSJed Brown /* AXPY with HYPRE and AIJ */ 294*c4762a1bSJed Brown ierr = MatAXPY(D,-1.0,C,str);CHKERRQ(ierr); 295*c4762a1bSJed Brown ierr = MatHasOperation(D,MATOP_NORM,&flg);CHKERRQ(ierr); 296*c4762a1bSJed Brown if (!flg) { /* TODO add MatNorm_HYPRE */ 297*c4762a1bSJed Brown ierr = MatConvert(D,MATAIJ,MAT_INPLACE_MATRIX,&D);CHKERRQ(ierr); 298*c4762a1bSJed Brown } 299*c4762a1bSJed Brown ierr = MatNorm(D,NORM_INFINITY,&err);CHKERRQ(ierr); 300*c4762a1bSJed Brown if (err > PETSC_SMALL) { 301*c4762a1bSJed Brown ierr = MatViewFromOptions(A,NULL,"-view_duplicate_diff");CHKERRQ(ierr); 302*c4762a1bSJed Brown ierr = MatViewFromOptions(C,NULL,"-view_duplicate_diff");CHKERRQ(ierr); 303*c4762a1bSJed Brown ierr = MatViewFromOptions(D,NULL,"-view_duplicate_diff");CHKERRQ(ierr); 304*c4762a1bSJed Brown SETERRQ3(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error test 3 MatDuplicate/MatCopy %g (%D,%D)",err,j,i); 305*c4762a1bSJed Brown } 306*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 307*c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 308*c4762a1bSJed Brown } 309*c4762a1bSJed Brown } 310*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 311*c4762a1bSJed Brown 312*c4762a1bSJed Brown ierr = MatHasCongruentLayouts(A,&flg);CHKERRQ(ierr); 313*c4762a1bSJed Brown if (flg) { 314*c4762a1bSJed Brown Vec y,y2; 315*c4762a1bSJed Brown 316*c4762a1bSJed Brown ierr = MatConvert(A,MATHYPRE,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 317*c4762a1bSJed Brown ierr = MatCreateVecs(A,NULL,&y);CHKERRQ(ierr); 318*c4762a1bSJed Brown ierr = MatCreateVecs(B,NULL,&y2);CHKERRQ(ierr); 319*c4762a1bSJed Brown ierr = MatGetDiagonal(A,y);CHKERRQ(ierr); 320*c4762a1bSJed Brown ierr = MatGetDiagonal(B,y2);CHKERRQ(ierr); 321*c4762a1bSJed Brown ierr = VecAXPY(y2,-1.0,y);CHKERRQ(ierr); 322*c4762a1bSJed Brown ierr = VecNorm(y2,NORM_INFINITY,&err);CHKERRQ(ierr); 323*c4762a1bSJed Brown if (err > PETSC_SMALL) { 324*c4762a1bSJed Brown ierr = VecViewFromOptions(y,NULL,"-view_diagonal_diff");CHKERRQ(ierr); 325*c4762a1bSJed Brown ierr = VecViewFromOptions(y2,NULL,"-view_diagonal_diff");CHKERRQ(ierr); 326*c4762a1bSJed Brown SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Error MatGetDiagonal %g",err); 327*c4762a1bSJed Brown } 328*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 329*c4762a1bSJed Brown ierr = VecDestroy(&y);CHKERRQ(ierr); 330*c4762a1bSJed Brown ierr = VecDestroy(&y2);CHKERRQ(ierr); 331*c4762a1bSJed Brown } 332*c4762a1bSJed Brown 333*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 334*c4762a1bSJed Brown 335*c4762a1bSJed Brown ierr = PetscFinalize(); 336*c4762a1bSJed Brown return ierr; 337*c4762a1bSJed Brown } 338*c4762a1bSJed Brown 339*c4762a1bSJed Brown 340*c4762a1bSJed Brown /*TEST 341*c4762a1bSJed Brown 342*c4762a1bSJed Brown build: 343*c4762a1bSJed Brown requires: hypre 344*c4762a1bSJed Brown 345*c4762a1bSJed Brown test: 346*c4762a1bSJed Brown suffix: 1 347*c4762a1bSJed Brown args: -N 11 -M 11 348*c4762a1bSJed Brown output_file: output/ex115_1.out 349*c4762a1bSJed Brown 350*c4762a1bSJed Brown test: 351*c4762a1bSJed Brown suffix: 2 352*c4762a1bSJed Brown nsize: 3 353*c4762a1bSJed Brown args: -N 13 -M 13 -matmatmult_via hypre 354*c4762a1bSJed Brown output_file: output/ex115_1.out 355*c4762a1bSJed Brown 356*c4762a1bSJed Brown test: 357*c4762a1bSJed Brown suffix: 3 358*c4762a1bSJed Brown nsize: 4 359*c4762a1bSJed Brown args: -M 13 -N 7 -matmatmult_via hypre 360*c4762a1bSJed Brown output_file: output/ex115_1.out 361*c4762a1bSJed Brown 362*c4762a1bSJed Brown test: 363*c4762a1bSJed Brown suffix: 4 364*c4762a1bSJed Brown nsize: 2 365*c4762a1bSJed Brown args: -M 12 -N 19 366*c4762a1bSJed Brown output_file: output/ex115_1.out 367*c4762a1bSJed Brown 368*c4762a1bSJed Brown test: 369*c4762a1bSJed Brown suffix: 5 370*c4762a1bSJed Brown nsize: 3 371*c4762a1bSJed Brown args: -M 13 -N 13 -matptap_via hypre -matptap_hypre_outtype hypre 372*c4762a1bSJed Brown output_file: output/ex115_1.out 373*c4762a1bSJed Brown 374*c4762a1bSJed Brown test: 375*c4762a1bSJed Brown suffix: 6 376*c4762a1bSJed Brown nsize: 3 377*c4762a1bSJed Brown args: -M 12 -N 19 -test_offproc 378*c4762a1bSJed Brown output_file: output/ex115_1.out 379*c4762a1bSJed Brown 380*c4762a1bSJed Brown test: 381*c4762a1bSJed Brown suffix: 7 382*c4762a1bSJed Brown nsize: 3 383*c4762a1bSJed Brown args: -M 19 -N 12 -test_offproc -view_B ::ascii_info_detail 384*c4762a1bSJed Brown output_file: output/ex115_7.out 385*c4762a1bSJed Brown 386*c4762a1bSJed Brown test: 387*c4762a1bSJed Brown suffix: 8 388*c4762a1bSJed Brown nsize: 3 389*c4762a1bSJed Brown args: -M 1 -N 12 -test_offproc 390*c4762a1bSJed Brown output_file: output/ex115_1.out 391*c4762a1bSJed Brown 392*c4762a1bSJed Brown test: 393*c4762a1bSJed Brown suffix: 9 394*c4762a1bSJed Brown nsize: 3 395*c4762a1bSJed Brown args: -M 1 -N 2 -test_offproc 396*c4762a1bSJed Brown output_file: output/ex115_1.out 397*c4762a1bSJed Brown 398*c4762a1bSJed Brown TEST*/ 399