1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests Elemental interface.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **args) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat Cdense,B,C,Ct; 9c4762a1bSJed Brown Vec d; 10c4762a1bSJed Brown PetscInt i,j,m = 5,n,nrows,ncols; 11c4762a1bSJed Brown const PetscInt *rows,*cols; 12c4762a1bSJed Brown IS isrows,iscols; 13c4762a1bSJed Brown PetscScalar *v; 14c4762a1bSJed Brown PetscMPIInt rank,size; 15c4762a1bSJed Brown PetscReal Cnorm; 16c4762a1bSJed Brown PetscBool flg,mats_view=PETSC_FALSE; 17c4762a1bSJed Brown 18*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 195f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 205f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 22c4762a1bSJed Brown n = m; 235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view)); 25c4762a1bSJed Brown 265f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(C,m,n,PETSC_DECIDE,PETSC_DECIDE)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(C,MATELEMENTAL)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(C)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(C)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipIS(C,&isrows,&iscols)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(isrows,&nrows)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(isrows,&rows)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(iscols,&ncols)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(iscols,&cols)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nrows*ncols,&v)); 37c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 38c4762a1bSJed Brown PetscRandom rand; 39c4762a1bSJed Brown PetscScalar rval; 405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rand)); 42c4762a1bSJed Brown for (i=0; i<nrows; i++) { 43c4762a1bSJed Brown for (j=0; j<ncols; j++) { 445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rand,&rval)); 45c4762a1bSJed Brown v[i*ncols+j] = rval; 46c4762a1bSJed Brown } 47c4762a1bSJed Brown } 485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rand)); 49c4762a1bSJed Brown #else 50c4762a1bSJed Brown for (i=0; i<nrows; i++) { 51c4762a1bSJed Brown for (j=0; j<ncols; j++) { 52c4762a1bSJed Brown v[i*ncols+j] = (PetscReal)(10000*rank+100*rows[i]+cols[j]); 53c4762a1bSJed Brown } 54c4762a1bSJed Brown } 55c4762a1bSJed Brown #endif 565f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,nrows,rows,ncols,cols,v,INSERT_VALUES)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(isrows,&rows)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(iscols,&cols)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isrows)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iscols)); 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* Test MatView(), MatDuplicate() and out-of-place MatConvert() */ 655f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(C,MAT_COPY_VALUES,&B)); 66c4762a1bSJed Brown if (mats_view) { 675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Duplicated C:\n")); 685f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 69c4762a1bSJed Brown } 705f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(C,MATMPIDENSE,MAT_INITIAL_MATRIX,&Cdense)); 725f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(C,Cdense,5,&flg)); 7328b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Cdense != C. MatConvert() fails"); 74c4762a1bSJed Brown 75c4762a1bSJed Brown /* Test MatNorm() */ 765f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(C,NORM_1,&Cnorm)); 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* Test MatTranspose(), MatZeroEntries() and MatGetDiagonal() */ 795f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(C,MAT_INITIAL_MATRIX,&Ct)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(MatConjugate(Ct)); 81c4762a1bSJed Brown if (mats_view) { 825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"C's Transpose Conjugate:\n")); 835f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(Ct,PETSC_VIEWER_STDOUT_WORLD)); 84c4762a1bSJed Brown } 855f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroEntries(Ct)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&d)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(d,m>n ? n : m,PETSC_DECIDE)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(d)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(C,d)); 90c4762a1bSJed Brown if (mats_view) { 915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Diagonal of C:\n")); 925f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(d,PETSC_VIEWER_STDOUT_WORLD)); 93c4762a1bSJed Brown } 94c4762a1bSJed Brown if (m>n) { 955f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(C,NULL,d)); 96c4762a1bSJed Brown } else { 975f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(C,d,NULL)); 98c4762a1bSJed Brown } 99c4762a1bSJed Brown if (mats_view) { 1005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Diagonal Scaled C:\n")); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 102c4762a1bSJed Brown } 103c4762a1bSJed Brown 104c4762a1bSJed Brown /* Test MatAXPY(), MatAYPX() and in-place MatConvert() */ 1055f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B)); 1065f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(B,m,n,PETSC_DECIDE,PETSC_DECIDE)); 1075f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(B,MATELEMENTAL)); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(B)); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(B)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipIS(B,&isrows,&iscols)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(isrows,&nrows)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(isrows,&rows)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(iscols,&ncols)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(iscols,&cols)); 115c4762a1bSJed Brown for (i=0; i<nrows; i++) { 116c4762a1bSJed Brown for (j=0; j<ncols; j++) { 117c4762a1bSJed Brown v[i*ncols+j] = (PetscReal)(1000*rows[i]+cols[j]); 118c4762a1bSJed Brown } 119c4762a1bSJed Brown } 1205f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,nrows,rows,ncols,cols,v,INSERT_VALUES)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(v)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(isrows,&rows)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(iscols,&cols)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 1255f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 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)); 129c4762a1bSJed Brown 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)); 132c4762a1bSJed Brown } 1335f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isrows)); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iscols)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* 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)); 14228b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MatMatTransposeMult: B != C*B^T"); 143c20d7725SJed Brown 144c4762a1bSJed Brown if (mats_view) { 1455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"C MatMatTransposeMult C:\n")); 1465f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 147c4762a1bSJed Brown } 148c4762a1bSJed Brown 1495f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Cdense)); 1505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(v)); 1515f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 1525f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 1535f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Ct)); 1545f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&d)); 155*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 156*b122ec5aSJacob Faibussowitsch return 0; 157c4762a1bSJed Brown } 158c4762a1bSJed Brown 159c4762a1bSJed Brown /*TEST 160c4762a1bSJed Brown 161c4762a1bSJed Brown test: 162c4762a1bSJed Brown nsize: 2 163c4762a1bSJed Brown args: -m 3 -n 2 164c4762a1bSJed Brown requires: elemental 165c4762a1bSJed Brown 166c4762a1bSJed Brown test: 167c4762a1bSJed Brown suffix: 2 168c4762a1bSJed Brown nsize: 6 169c4762a1bSJed Brown args: -m 2 -n 3 170c4762a1bSJed Brown requires: elemental 171c4762a1bSJed Brown 172c20d7725SJed Brown test: 173c20d7725SJed Brown suffix: 3 174c20d7725SJed Brown nsize: 1 175c20d7725SJed Brown args: -m 2 -n 3 176c20d7725SJed Brown requires: elemental 177c20d7725SJed Brown output_file: output/ex39_1.out 178c20d7725SJed Brown 179c4762a1bSJed Brown TEST*/ 180