1d24d4204SJose E. Roman static char help[] = "Test conversion of ScaLAPACK matrices.\n\n"; 2d24d4204SJose E. Roman 3d24d4204SJose E. Roman #include <petscmat.h> 4d24d4204SJose E. Roman 5d24d4204SJose E. Roman int main(int argc, char** argv) 6d24d4204SJose E. Roman { 7d24d4204SJose E. Roman Mat A,A_scalapack; 8d24d4204SJose E. Roman PetscInt i,j,M=10,N=5,nloc,mloc,nrows,ncols; 9d24d4204SJose E. Roman PetscErrorCode ierr; 10d24d4204SJose E. Roman PetscMPIInt rank,size; 11d24d4204SJose E. Roman IS isrows,iscols; 12d24d4204SJose E. Roman const PetscInt *rows,*cols; 13d24d4204SJose E. Roman PetscScalar *v; 14d24d4204SJose E. Roman MatType type; 15d24d4204SJose E. Roman PetscBool isDense,isAIJ,flg; 16d24d4204SJose E. Roman 17d24d4204SJose E. Roman ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 185f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 195f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL)); 22d24d4204SJose E. Roman 23d24d4204SJose E. Roman /* Create a matrix */ 245f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD, &A)); 25d24d4204SJose E. Roman mloc = PETSC_DECIDE; 265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSplitOwnershipEqual(PETSC_COMM_WORLD,&mloc,&M)); 27d24d4204SJose E. Roman nloc = PETSC_DECIDE; 285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSplitOwnershipEqual(PETSC_COMM_WORLD,&nloc,&N)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,mloc,nloc,M,N)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATDENSE)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 33d24d4204SJose E. Roman 34d24d4204SJose E. Roman /* Set local matrix entries */ 355f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipIS(A,&isrows,&iscols)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(isrows,&nrows)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(isrows,&rows)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(iscols,&ncols)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(iscols,&cols)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nrows*ncols,&v)); 41d24d4204SJose E. Roman 42d24d4204SJose E. Roman for (i=0; i<nrows; i++) { 43d24d4204SJose E. Roman for (j=0; j<ncols; j++) { 44d24d4204SJose E. Roman if (size == 1) { 45d24d4204SJose E. Roman v[i*ncols+j] = (PetscScalar)(i+j); 46d24d4204SJose E. Roman } else { 47d24d4204SJose E. Roman v[i*ncols+j] = (PetscScalar)rank+j*0.1; 48d24d4204SJose E. Roman } 49d24d4204SJose E. Roman } 50d24d4204SJose E. Roman } 515f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 54d24d4204SJose E. Roman 55d24d4204SJose E. Roman /* Test MatSetValues() by converting A to A_scalapack */ 565f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetType(A,&type)); 57d24d4204SJose E. Roman if (size == 1) { 585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATSEQDENSE,&isDense)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isAIJ)); 60d24d4204SJose E. Roman } else { 615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATMPIDENSE,&isDense)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isAIJ)); 63d24d4204SJose E. Roman } 64d24d4204SJose E. Roman 65d24d4204SJose E. Roman if (isDense || isAIJ) { 66d24d4204SJose E. Roman Mat Aexplicit; 675f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATSCALAPACK,MAT_INITIAL_MATRIX,&A_scalapack)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(MatComputeOperator(A_scalapack,isAIJ?MATAIJ:MATDENSE,&Aexplicit)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(Aexplicit,A_scalapack,5,&flg)); 70*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Aexplicit != A_scalapack."); 715f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Aexplicit)); 72d24d4204SJose E. Roman 73d24d4204SJose E. Roman /* Test MAT_REUSE_MATRIX which is only supported for inplace conversion */ 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATSCALAPACK,MAT_INPLACE_MATRIX,&A)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(A_scalapack,A,5,&flg)); 76*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A_scalapack != A."); 775f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A_scalapack)); 78d24d4204SJose E. Roman } 79d24d4204SJose E. Roman 805f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(isrows,&rows)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(iscols,&cols)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isrows)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iscols)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(v)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 86d24d4204SJose E. Roman ierr = PetscFinalize(); 87d24d4204SJose E. Roman return ierr; 88d24d4204SJose E. Roman } 89d24d4204SJose E. Roman 90d24d4204SJose E. Roman /*TEST 91d24d4204SJose E. Roman 92d24d4204SJose E. Roman build: 93d24d4204SJose E. Roman requires: scalapack 94d24d4204SJose E. Roman 95d24d4204SJose E. Roman test: 96d24d4204SJose E. Roman nsize: 6 97d24d4204SJose E. Roman 98d24d4204SJose E. Roman test: 99d24d4204SJose E. Roman suffix: 2 100d24d4204SJose E. Roman nsize: 6 101d24d4204SJose E. Roman args: -mat_type aij 102d24d4204SJose E. Roman output_file: output/ex243_1.out 103d24d4204SJose E. Roman 104d24d4204SJose E. Roman test: 105d24d4204SJose E. Roman suffix: 3 106d24d4204SJose E. Roman nsize: 6 107d24d4204SJose E. Roman args: -mat_type scalapack 108d24d4204SJose E. Roman output_file: output/ex243_1.out 109d24d4204SJose E. Roman 110d24d4204SJose E. Roman TEST*/ 111