1d24d4204SJose E. Roman 2d24d4204SJose E. Roman static char help[] = "Tests ScaLAPACK interface.\n\n"; 3d24d4204SJose E. Roman 4d24d4204SJose E. Roman #include <petscmat.h> 5d24d4204SJose E. Roman 6d24d4204SJose E. Roman int main(int argc,char **args) 7d24d4204SJose E. Roman { 8f7ec113fSDamian Marek Mat Cdense,Caij,B,C,Ct,Asub; 9d24d4204SJose E. Roman Vec d; 10d24d4204SJose E. Roman PetscInt i,j,M = 5,N,mb = 2,nb,nrows,ncols,mloc,nloc; 11d24d4204SJose E. Roman const PetscInt *rows,*cols; 12d24d4204SJose E. Roman IS isrows,iscols; 13d24d4204SJose E. Roman PetscErrorCode ierr; 14d24d4204SJose E. Roman PetscScalar *v; 15f7ec113fSDamian Marek PetscMPIInt rank,color; 16d24d4204SJose E. Roman PetscReal Cnorm; 17d24d4204SJose E. Roman PetscBool flg,mats_view=PETSC_FALSE; 18f7ec113fSDamian Marek MPI_Comm subcomm; 19d24d4204SJose E. Roman 20d24d4204SJose E. Roman ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 215f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 23d24d4204SJose E. Roman N = M; 245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mb",&mb,NULL)); 26d24d4204SJose E. Roman nb = mb; 275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nb",&nb,NULL)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view)); 29d24d4204SJose E. Roman 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(C,MATSCALAPACK)); 32d24d4204SJose E. Roman mloc = PETSC_DECIDE; 335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSplitOwnershipEqual(PETSC_COMM_WORLD,&mloc,&M)); 34d24d4204SJose E. Roman nloc = PETSC_DECIDE; 355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSplitOwnershipEqual(PETSC_COMM_WORLD,&nloc,&N)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(C,mloc,nloc,M,N)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatScaLAPACKSetBlockSizes(C,mb,nb)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(C)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(C)); 405f80ce2aSJacob Faibussowitsch /*CHKERRQ(MatCreateScaLAPACK(PETSC_COMM_WORLD,mb,nb,M,N,0,0,&C)); */ 41d24d4204SJose E. Roman 425f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipIS(C,&isrows,&iscols)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(isrows,&nrows)); 445f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(isrows,&rows)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(iscols,&ncols)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(iscols,&cols)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nrows*ncols,&v)); 48d24d4204SJose E. Roman for (i=0;i<nrows;i++) { 49d24d4204SJose E. Roman for (j=0;j<ncols;j++) v[i*ncols+j] = (PetscReal)(rows[i]+1+(cols[j]+1)*0.01); 50d24d4204SJose E. Roman } 515f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,nrows,rows,ncols,cols,v,INSERT_VALUES)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(v)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(isrows,&rows)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(iscols,&cols)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isrows)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iscols)); 59d24d4204SJose E. Roman 60d24d4204SJose E. Roman /* Test MatView(), MatDuplicate() and out-of-place MatConvert() */ 615f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(C,MAT_COPY_VALUES,&B)); 62d24d4204SJose E. Roman if (mats_view) { 635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Duplicated C:\n")); 645f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 65d24d4204SJose E. Roman } 665f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(C,MATDENSE,MAT_INITIAL_MATRIX,&Cdense)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(C,Cdense,5,&flg)); 69*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Check fails: Cdense != C"); 70d24d4204SJose E. Roman 71d24d4204SJose E. Roman /* Test MatNorm() */ 725f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(C,NORM_1,&Cnorm)); 73d24d4204SJose E. Roman 74d24d4204SJose E. Roman /* Test MatTranspose(), MatZeroEntries() and MatGetDiagonal() */ 755f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(C,MAT_INITIAL_MATRIX,&Ct)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(MatConjugate(Ct)); 77d24d4204SJose E. Roman if (mats_view) { 785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"C's Transpose Conjugate:\n")); 795f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(Ct,PETSC_VIEWER_STDOUT_WORLD)); 80d24d4204SJose E. Roman } 815f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroEntries(Ct)); 825f80ce2aSJacob Faibussowitsch if (M>N) CHKERRQ(MatCreateVecs(C,&d,NULL)); 835f80ce2aSJacob Faibussowitsch else CHKERRQ(MatCreateVecs(C,NULL,&d)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(C,d)); 85d24d4204SJose E. Roman if (mats_view) { 865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Diagonal of C:\n")); 875f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(d,PETSC_VIEWER_STDOUT_WORLD)); 88d24d4204SJose E. Roman } 89d24d4204SJose E. Roman if (M>N) { 905f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(C,NULL,d)); 91d24d4204SJose E. Roman } else { 925f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(C,d,NULL)); 93d24d4204SJose E. Roman } 94d24d4204SJose E. Roman if (mats_view) { 955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Diagonal Scaled C:\n")); 965f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 97d24d4204SJose E. Roman } 98d24d4204SJose E. Roman 99d24d4204SJose E. Roman /* Test MatAXPY(), MatAYPX() and in-place MatConvert() */ 1005f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(B,MATSCALAPACK)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(B,mloc,nloc,PETSC_DECIDE,PETSC_DECIDE)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(MatScaLAPACKSetBlockSizes(B,mb,nb)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(B)); 1055f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(B)); 1065f80ce2aSJacob Faibussowitsch /* CHKERRQ(MatCreateScaLAPACK(PETSC_COMM_WORLD,mb,nb,M,N,0,0,&B)); */ 1075f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipIS(B,&isrows,&iscols)); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(isrows,&nrows)); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(isrows,&rows)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(iscols,&ncols)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(iscols,&cols)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nrows*ncols,&v)); 113d24d4204SJose E. Roman for (i=0;i<nrows;i++) { 114d24d4204SJose E. Roman for (j=0;j<ncols;j++) v[i*ncols+j] = (PetscReal)(1000*rows[i]+cols[j]); 115d24d4204SJose E. Roman } 1165f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,nrows,rows,ncols,cols,v,INSERT_VALUES)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(v)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(isrows,&rows)); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(iscols,&cols)); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 122d24d4204SJose E. Roman if (mats_view) { 1235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"B:\n")); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 125d24d4204SJose E. Roman } 1265f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(B,2.5,C,SAME_NONZERO_PATTERN)); 1275f80ce2aSJacob Faibussowitsch CHKERRQ(MatAYPX(B,3.75,C,SAME_NONZERO_PATTERN)); 1285f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(B,MATDENSE,MAT_INPLACE_MATRIX,&B)); 129d24d4204SJose E. Roman if (mats_view) { 1305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"B after MatAXPY and MatAYPX:\n")); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 132d24d4204SJose E. Roman } 1335f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isrows)); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iscols)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 136d24d4204SJose E. Roman 137d24d4204SJose E. Roman /* Test MatMatTransposeMult(): B = C*C^T */ 1385f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatTransposeMult(C,C,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(C,2.0)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatTransposeMult(C,C,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatTransposeMultEqual(C,C,B,10,&flg)); 142*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Check fails: B != C*C^T"); 143d24d4204SJose E. Roman 144d24d4204SJose E. Roman if (mats_view) { 1455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"C MatMatTransposeMult C:\n")); 1465f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 147d24d4204SJose E. Roman } 148d24d4204SJose E. Roman 149d24d4204SJose E. Roman /* Test MatMult() */ 1505f80ce2aSJacob Faibussowitsch CHKERRQ(MatComputeOperator(C,MATAIJ,&Caij)); 1515f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(C,Caij,5,&flg)); 152*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_NOTSAMETYPE,"C != Caij. MatMultEqual() fails"); 1535f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeEqual(C,Caij,5,&flg)); 154*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_NOTSAMETYPE,"C != Caij. MatMultTransposeEqual() fails"); 155d24d4204SJose E. Roman 156d24d4204SJose E. Roman /* Test MatMultAdd() and MatMultTransposeAddEqual() */ 1575f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAddEqual(C,Caij,5,&flg)); 158*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_NOTSAMETYPE,"C != Caij. MatMultAddEqual() fails"); 1595f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAddEqual(C,Caij,5,&flg)); 160*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_NOTSAMETYPE,"C != Caij. MatMultTransposeAddEqual() fails"); 161d24d4204SJose E. Roman 162d24d4204SJose E. Roman /* Test MatMatMult() */ 1635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-test_matmatmult",&flg)); 164d24d4204SJose E. Roman if (flg) { 165d24d4204SJose E. Roman Mat CC,CCaij; 1665f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(C,C,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&CC)); 1675f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(Caij,Caij,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&CCaij)); 1685f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(CC,CCaij,5,&flg)); 169*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_NOTSAMETYPE,"CC != CCaij. MatMatMult() fails"); 1705f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&CCaij)); 1715f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&CC)); 172d24d4204SJose E. Roman } 173d24d4204SJose E. Roman 174f7ec113fSDamian Marek /* Test MatCreate() on subcomm */ 175f7ec113fSDamian Marek color = rank%2; 1765f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_split(PETSC_COMM_WORLD,color,0,&subcomm)); 177f7ec113fSDamian Marek if (color==0) { 1785f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(subcomm,&Asub)); 1795f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(Asub,MATSCALAPACK)); 180f7ec113fSDamian Marek mloc = PETSC_DECIDE; 1815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSplitOwnershipEqual(subcomm,&mloc,&M)); 182f7ec113fSDamian Marek nloc = PETSC_DECIDE; 1835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSplitOwnershipEqual(subcomm,&nloc,&N)); 1845f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(Asub,mloc,nloc,M,N)); 1855f80ce2aSJacob Faibussowitsch CHKERRQ(MatScaLAPACKSetBlockSizes(Asub,mb,nb)); 1865f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(Asub)); 1875f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(Asub)); 1885f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Asub)); 189f7ec113fSDamian Marek } 190f7ec113fSDamian Marek 1915f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Cdense)); 1925f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Caij)); 1935f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 1945f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 1955f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Ct)); 1965f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&d)); 1975f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_free(&subcomm)); 198d24d4204SJose E. Roman ierr = PetscFinalize(); 199d24d4204SJose E. Roman return ierr; 200d24d4204SJose E. Roman } 201d24d4204SJose E. Roman 202d24d4204SJose E. Roman /*TEST 203d24d4204SJose E. Roman 204d24d4204SJose E. Roman build: 205d24d4204SJose E. Roman requires: scalapack 206d24d4204SJose E. Roman 207d24d4204SJose E. Roman test: 208d24d4204SJose E. Roman nsize: 2 209d24d4204SJose E. Roman args: -mb 5 -nb 5 -M 12 -N 10 210d24d4204SJose E. Roman requires: scalapack 211d24d4204SJose E. Roman 212d24d4204SJose E. Roman test: 213d24d4204SJose E. Roman suffix: 2 214d24d4204SJose E. Roman nsize: 6 215d24d4204SJose E. Roman args: -mb 8 -nb 6 -M 20 -N 50 216d24d4204SJose E. Roman requires: scalapack 217d24d4204SJose E. Roman output_file: output/ex242_1.out 218d24d4204SJose E. Roman 219d24d4204SJose E. Roman test: 220d24d4204SJose E. Roman suffix: 3 221d24d4204SJose E. Roman nsize: 3 222d24d4204SJose E. Roman args: -mb 2 -nb 2 -M 20 -N 20 -test_matmatmult 223d24d4204SJose E. Roman requires: scalapack 224d24d4204SJose E. Roman output_file: output/ex242_1.out 225d24d4204SJose E. Roman 226d24d4204SJose E. Roman TEST*/ 227