1*c4762a1bSJed Brown static char help[] = "Tests MatComputeOperator() and MatComputeOperatorTranspose()\n\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscmat.h> 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown int main(int argc,char **argv) 6*c4762a1bSJed Brown { 7*c4762a1bSJed Brown Mat A,Ae,Aet; 8*c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN]; 9*c4762a1bSJed Brown char expltype[128],*etype = NULL; 10*c4762a1bSJed Brown PetscInt bs = 1; 11*c4762a1bSJed Brown PetscBool flg, check = PETSC_TRUE; 12*c4762a1bSJed Brown PetscErrorCode ierr; 13*c4762a1bSJed Brown 14*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*) 0,help);if (ierr) return ierr; 15*c4762a1bSJed Brown 16*c4762a1bSJed Brown ierr = PetscOptionsGetString(NULL,NULL,"-expl_type",expltype,sizeof(expltype),&flg);CHKERRQ(ierr); 17*c4762a1bSJed Brown if (flg) { 18*c4762a1bSJed Brown ierr = PetscStrallocpy(expltype,&etype);CHKERRQ(ierr); 19*c4762a1bSJed Brown } 20*c4762a1bSJed Brown ierr = PetscOptionsGetString(NULL,NULL,"-f",filename,sizeof(filename),&flg);CHKERRQ(ierr); 21*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);CHKERRQ(ierr); 22*c4762a1bSJed Brown if (!flg) { 23*c4762a1bSJed Brown PetscInt M = 13,N = 6; 24*c4762a1bSJed Brown 25*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr); 26*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);CHKERRQ(ierr); 27*c4762a1bSJed Brown ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,M,N,NULL,&A);CHKERRQ(ierr); 28*c4762a1bSJed Brown ierr = MatSetBlockSize(A,bs);CHKERRQ(ierr); 29*c4762a1bSJed Brown ierr = MatSetRandom(A,NULL);CHKERRQ(ierr); 30*c4762a1bSJed Brown } else { 31*c4762a1bSJed Brown PetscViewer viewer; 32*c4762a1bSJed Brown 33*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr); 34*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 35*c4762a1bSJed Brown ierr = MatSetBlockSize(A,bs);CHKERRQ(ierr); 36*c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 37*c4762a1bSJed Brown ierr = MatLoad(A,viewer);CHKERRQ(ierr); 38*c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 39*c4762a1bSJed Brown } 40*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)A,"Matrix");CHKERRQ(ierr); 41*c4762a1bSJed Brown ierr = MatViewFromOptions(A,NULL,"-view_expl");CHKERRQ(ierr); 42*c4762a1bSJed Brown 43*c4762a1bSJed Brown ierr = MatComputeOperator(A,etype,&Ae);CHKERRQ(ierr); 44*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)Ae,"Explicit matrix");CHKERRQ(ierr); 45*c4762a1bSJed Brown ierr = MatViewFromOptions(Ae,NULL,"-view_expl");CHKERRQ(ierr); 46*c4762a1bSJed Brown 47*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-check",&check,NULL);CHKERRQ(ierr); 48*c4762a1bSJed Brown if (check) { 49*c4762a1bSJed Brown Mat A2; 50*c4762a1bSJed Brown PetscReal err,tol = PETSC_SMALL; 51*c4762a1bSJed Brown 52*c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-tol",&tol,NULL);CHKERRQ(ierr); 53*c4762a1bSJed Brown ierr = MatConvert(A,etype,MAT_INITIAL_MATRIX,&A2);CHKERRQ(ierr); 54*c4762a1bSJed Brown ierr = MatAXPY(A2,-1.0,Ae,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 55*c4762a1bSJed Brown ierr = MatNorm(A2,NORM_FROBENIUS,&err);CHKERRQ(ierr); 56*c4762a1bSJed Brown if (err > tol) { 57*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Error %g > %g (type %s)\n",(double)err,(double)tol,etype);CHKERRQ(ierr); 58*c4762a1bSJed Brown } 59*c4762a1bSJed Brown ierr = MatDestroy(&A2);CHKERRQ(ierr); 60*c4762a1bSJed Brown } 61*c4762a1bSJed Brown 62*c4762a1bSJed Brown ierr = MatComputeOperatorTranspose(A,etype,&Aet);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)Aet,"Explicit matrix transpose");CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = MatViewFromOptions(Aet,NULL,"-view_expl");CHKERRQ(ierr); 65*c4762a1bSJed Brown 66*c4762a1bSJed Brown ierr = PetscFree(etype);CHKERRQ(ierr); 67*c4762a1bSJed Brown ierr = MatDestroy(&Ae);CHKERRQ(ierr); 68*c4762a1bSJed Brown ierr = MatDestroy(&Aet);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = PetscFinalize(); 71*c4762a1bSJed Brown return ierr; 72*c4762a1bSJed Brown } 73*c4762a1bSJed Brown 74*c4762a1bSJed Brown 75*c4762a1bSJed Brown 76*c4762a1bSJed Brown /*TEST 77*c4762a1bSJed Brown 78*c4762a1bSJed Brown test: 79*c4762a1bSJed Brown output_file: output/ex222_null.out 80*c4762a1bSJed Brown 81*c4762a1bSJed Brown testset: 82*c4762a1bSJed Brown suffix: matexpl_rect 83*c4762a1bSJed Brown output_file: output/ex222_null.out 84*c4762a1bSJed Brown nsize: {{1 3}} 85*c4762a1bSJed Brown args: -expl_type {{dense aij baij}} 86*c4762a1bSJed Brown 87*c4762a1bSJed Brown testset: 88*c4762a1bSJed Brown suffix: matexpl_square 89*c4762a1bSJed Brown output_file: output/ex222_null.out 90*c4762a1bSJed Brown nsize: {{1 3}} 91*c4762a1bSJed Brown args: -bs {{1 2 3}} -M 36 -N 36 -expl_type {{dense aij baij sbaij}} 92*c4762a1bSJed Brown 93*c4762a1bSJed Brown TEST*/ 94