1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Tests MatGetRowMax(), MatGetRowMin(), MatGetRowMaxAbs()\n"; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown #include <petscmat.h> 5*c4762a1bSJed Brown 6*c4762a1bSJed Brown #define M 5 7*c4762a1bSJed Brown #define N 6 8*c4762a1bSJed Brown 9*c4762a1bSJed Brown int main(int argc,char **args) 10*c4762a1bSJed Brown { 11*c4762a1bSJed Brown Mat A; 12*c4762a1bSJed Brown Vec min,max,maxabs; 13*c4762a1bSJed Brown PetscInt m,n; 14*c4762a1bSJed Brown PetscInt imin[M],imax[M],imaxabs[M],indices[N],row; 15*c4762a1bSJed Brown PetscScalar values[N]; 16*c4762a1bSJed Brown PetscErrorCode ierr; 17*c4762a1bSJed Brown MatType type; 18*c4762a1bSJed Brown PetscMPIInt size; 19*c4762a1bSJed Brown PetscBool doTest=PETSC_TRUE; 20*c4762a1bSJed Brown 21*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 22*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 25*c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 26*c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 27*c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 28*c4762a1bSJed Brown row = 0; 29*c4762a1bSJed Brown 30*c4762a1bSJed Brown indices[0] = 0; indices[1] = 1; indices[2] = 2; indices[3] = 3; indices[4] = 4; indices[5] = 5; 31*c4762a1bSJed Brown values[0] = -1.0; values[1] = 0.0; values[2] = 1.0; values[3] = 3.0; values[4] = 4.0; values[5] = -5.0; 32*c4762a1bSJed Brown 33*c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,6,indices,values,INSERT_VALUES);CHKERRQ(ierr); 34*c4762a1bSJed Brown row = 1; 35*c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES);CHKERRQ(ierr); 36*c4762a1bSJed Brown row = 4; 37*c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,1,indices+4,values+4,INSERT_VALUES);CHKERRQ(ierr); 38*c4762a1bSJed Brown row = 4; 39*c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,2,indices+4,values+4,INSERT_VALUES);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 41*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 42*c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 43*c4762a1bSJed Brown 44*c4762a1bSJed Brown ierr = MatGetLocalSize(A, &m,&n);CHKERRQ(ierr); 45*c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&min);CHKERRQ(ierr); 46*c4762a1bSJed Brown ierr = VecSetSizes(min,m,PETSC_DECIDE);CHKERRQ(ierr); 47*c4762a1bSJed Brown ierr = VecSetFromOptions(min);CHKERRQ(ierr); 48*c4762a1bSJed Brown ierr = VecDuplicate(min,&max);CHKERRQ(ierr); 49*c4762a1bSJed Brown ierr = VecDuplicate(min,&maxabs);CHKERRQ(ierr); 50*c4762a1bSJed Brown 51*c4762a1bSJed Brown /* Test MatGetRowMin, MatGetRowMax and MatGetRowMaxAbs */ 52*c4762a1bSJed Brown if (size == 1) { 53*c4762a1bSJed Brown ierr = MatGetRowMin(A,min,imin);CHKERRQ(ierr); 54*c4762a1bSJed Brown ierr = MatGetRowMax(A,max,imax);CHKERRQ(ierr); 55*c4762a1bSJed Brown ierr = MatGetRowMaxAbs(A,maxabs,imaxabs);CHKERRQ(ierr); 56*c4762a1bSJed Brown 57*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Row Minimums\n");CHKERRQ(ierr); 58*c4762a1bSJed Brown ierr = VecView(min,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 59*c4762a1bSJed Brown ierr = PetscIntView(5,imin,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 60*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Row Maximums\n");CHKERRQ(ierr); 61*c4762a1bSJed Brown ierr = VecView(max,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 62*c4762a1bSJed Brown ierr = PetscIntView(5,imax,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Row Maximum Absolute Values\n");CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = VecView(maxabs,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = PetscIntView(5,imaxabs,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 66*c4762a1bSJed Brown 67*c4762a1bSJed Brown } else { 68*c4762a1bSJed Brown ierr = MatGetType(A,&type);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\nMatrix type: %s\n",type);CHKERRQ(ierr); 70*c4762a1bSJed Brown /* AIJ */ 71*c4762a1bSJed Brown ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&doTest);CHKERRQ(ierr); 72*c4762a1bSJed Brown if (doTest) { 73*c4762a1bSJed Brown ierr = MatGetRowMaxAbs(A,maxabs,NULL);CHKERRQ(ierr); 74*c4762a1bSJed Brown ierr = MatGetRowMaxAbs(A,maxabs,imaxabs);CHKERRQ(ierr); 75*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Row Maximum Absolute Values:\n");CHKERRQ(ierr); 76*c4762a1bSJed Brown ierr = VecView(maxabs,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 77*c4762a1bSJed Brown } 78*c4762a1bSJed Brown /* BAIJ */ 79*c4762a1bSJed Brown ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIBAIJ,&doTest);CHKERRQ(ierr); 80*c4762a1bSJed Brown if (doTest) { 81*c4762a1bSJed Brown ierr = MatGetRowMaxAbs(A,maxabs,NULL);CHKERRQ(ierr); 82*c4762a1bSJed Brown ierr = MatGetRowMaxAbs(A,maxabs,imaxabs);CHKERRQ(ierr); 83*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Row Maximum Absolute Values:\n");CHKERRQ(ierr); 84*c4762a1bSJed Brown ierr = VecView(maxabs,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 85*c4762a1bSJed Brown } 86*c4762a1bSJed Brown } 87*c4762a1bSJed Brown 88*c4762a1bSJed Brown if (size == 1) { 89*c4762a1bSJed Brown ierr = MatConvert(A,MATDENSE,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr); 90*c4762a1bSJed Brown 91*c4762a1bSJed Brown ierr = MatGetRowMin(A,min,imin);CHKERRQ(ierr); 92*c4762a1bSJed Brown ierr = MatGetRowMax(A,max,imax);CHKERRQ(ierr); 93*c4762a1bSJed Brown ierr = MatGetRowMaxAbs(A,maxabs,imaxabs);CHKERRQ(ierr); 94*c4762a1bSJed Brown 95*c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 96*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Row Minimums\n");CHKERRQ(ierr); 97*c4762a1bSJed Brown ierr = VecView(min,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 98*c4762a1bSJed Brown ierr = PetscIntView(5,imin,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 99*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Row Maximums\n");CHKERRQ(ierr); 100*c4762a1bSJed Brown ierr = VecView(max,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 101*c4762a1bSJed Brown ierr = PetscIntView(5,imax,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 102*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Row Maximum Absolute Values\n");CHKERRQ(ierr); 103*c4762a1bSJed Brown ierr = VecView(maxabs,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 104*c4762a1bSJed Brown ierr = PetscIntView(5,imaxabs,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 105*c4762a1bSJed Brown } 106*c4762a1bSJed Brown 107*c4762a1bSJed Brown ierr = VecDestroy(&min);CHKERRQ(ierr); 108*c4762a1bSJed Brown ierr = VecDestroy(&max);CHKERRQ(ierr); 109*c4762a1bSJed Brown ierr = VecDestroy(&maxabs);CHKERRQ(ierr); 110*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 111*c4762a1bSJed Brown ierr = PetscFinalize(); 112*c4762a1bSJed Brown return ierr; 113*c4762a1bSJed Brown } 114*c4762a1bSJed Brown 115*c4762a1bSJed Brown 116*c4762a1bSJed Brown 117*c4762a1bSJed Brown /*TEST 118*c4762a1bSJed Brown 119*c4762a1bSJed Brown test: 120*c4762a1bSJed Brown output_file: output/ex114.out 121*c4762a1bSJed Brown 122*c4762a1bSJed Brown test: 123*c4762a1bSJed Brown suffix: 2 124*c4762a1bSJed Brown nsize: 2 125*c4762a1bSJed Brown 126*c4762a1bSJed Brown test: 127*c4762a1bSJed Brown suffix: 3 128*c4762a1bSJed Brown nsize: 2 129*c4762a1bSJed Brown args: -mat_type baij 130*c4762a1bSJed Brown 131*c4762a1bSJed Brown TEST*/ 132