1 static char help[] = "Test conversion of ScaLAPACK matrices.\n\n"; 2 3 #include <petscmat.h> 4 5 int main(int argc, char** argv) 6 { 7 Mat A,A_scalapack; 8 PetscInt i,j,M=10,N=5,nloc,mloc,nrows,ncols; 9 PetscErrorCode ierr; 10 PetscMPIInt rank,size; 11 IS isrows,iscols; 12 const PetscInt *rows,*cols; 13 PetscScalar *v; 14 MatType type; 15 PetscBool isDense,isAIJ,flg; 16 17 ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 18 CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 19 CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 20 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 21 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL)); 22 23 /* Create a matrix */ 24 CHKERRQ(MatCreate(PETSC_COMM_WORLD, &A)); 25 mloc = PETSC_DECIDE; 26 CHKERRQ(PetscSplitOwnershipEqual(PETSC_COMM_WORLD,&mloc,&M)); 27 nloc = PETSC_DECIDE; 28 CHKERRQ(PetscSplitOwnershipEqual(PETSC_COMM_WORLD,&nloc,&N)); 29 CHKERRQ(MatSetSizes(A,mloc,nloc,M,N)); 30 CHKERRQ(MatSetType(A,MATDENSE)); 31 CHKERRQ(MatSetFromOptions(A)); 32 CHKERRQ(MatSetUp(A)); 33 34 /* Set local matrix entries */ 35 CHKERRQ(MatGetOwnershipIS(A,&isrows,&iscols)); 36 CHKERRQ(ISGetLocalSize(isrows,&nrows)); 37 CHKERRQ(ISGetIndices(isrows,&rows)); 38 CHKERRQ(ISGetLocalSize(iscols,&ncols)); 39 CHKERRQ(ISGetIndices(iscols,&cols)); 40 CHKERRQ(PetscMalloc1(nrows*ncols,&v)); 41 42 for (i=0; i<nrows; i++) { 43 for (j=0; j<ncols; j++) { 44 if (size == 1) { 45 v[i*ncols+j] = (PetscScalar)(i+j); 46 } else { 47 v[i*ncols+j] = (PetscScalar)rank+j*0.1; 48 } 49 } 50 } 51 CHKERRQ(MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES)); 52 CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 53 CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 54 55 /* Test MatSetValues() by converting A to A_scalapack */ 56 CHKERRQ(MatGetType(A,&type)); 57 if (size == 1) { 58 CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATSEQDENSE,&isDense)); 59 CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isAIJ)); 60 } else { 61 CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATMPIDENSE,&isDense)); 62 CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isAIJ)); 63 } 64 65 if (isDense || isAIJ) { 66 Mat Aexplicit; 67 CHKERRQ(MatConvert(A,MATSCALAPACK,MAT_INITIAL_MATRIX,&A_scalapack)); 68 CHKERRQ(MatComputeOperator(A_scalapack,isAIJ?MATAIJ:MATDENSE,&Aexplicit)); 69 CHKERRQ(MatMultEqual(Aexplicit,A_scalapack,5,&flg)); 70 PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Aexplicit != A_scalapack."); 71 CHKERRQ(MatDestroy(&Aexplicit)); 72 73 /* Test MAT_REUSE_MATRIX which is only supported for inplace conversion */ 74 CHKERRQ(MatConvert(A,MATSCALAPACK,MAT_INPLACE_MATRIX,&A)); 75 CHKERRQ(MatMultEqual(A_scalapack,A,5,&flg)); 76 PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A_scalapack != A."); 77 CHKERRQ(MatDestroy(&A_scalapack)); 78 } 79 80 CHKERRQ(ISRestoreIndices(isrows,&rows)); 81 CHKERRQ(ISRestoreIndices(iscols,&cols)); 82 CHKERRQ(ISDestroy(&isrows)); 83 CHKERRQ(ISDestroy(&iscols)); 84 CHKERRQ(PetscFree(v)); 85 CHKERRQ(MatDestroy(&A)); 86 ierr = PetscFinalize(); 87 return ierr; 88 } 89 90 /*TEST 91 92 build: 93 requires: scalapack 94 95 test: 96 nsize: 6 97 98 test: 99 suffix: 2 100 nsize: 6 101 args: -mat_type aij 102 output_file: output/ex243_1.out 103 104 test: 105 suffix: 3 106 nsize: 6 107 args: -mat_type scalapack 108 output_file: output/ex243_1.out 109 110 TEST*/ 111