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 PetscErrorCode ierr; 14c4762a1bSJed Brown PetscScalar *v; 15c4762a1bSJed Brown PetscMPIInt rank,size; 16c4762a1bSJed Brown PetscReal Cnorm; 17c4762a1bSJed Brown PetscBool flg,mats_view=PETSC_FALSE; 18c4762a1bSJed Brown 19c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 205f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 215f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 23c4762a1bSJed Brown n = m; 245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view)); 26c4762a1bSJed Brown 275f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(C,m,n,PETSC_DECIDE,PETSC_DECIDE)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(C,MATELEMENTAL)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(C)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(C)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipIS(C,&isrows,&iscols)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(isrows,&nrows)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(isrows,&rows)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(iscols,&ncols)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(iscols,&cols)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nrows*ncols,&v)); 38c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 39c4762a1bSJed Brown PetscRandom rand; 40c4762a1bSJed Brown PetscScalar rval; 415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rand)); 43c4762a1bSJed Brown for (i=0; i<nrows; i++) { 44c4762a1bSJed Brown for (j=0; j<ncols; j++) { 455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rand,&rval)); 46c4762a1bSJed Brown v[i*ncols+j] = rval; 47c4762a1bSJed Brown } 48c4762a1bSJed Brown } 495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rand)); 50c4762a1bSJed Brown #else 51c4762a1bSJed Brown for (i=0; i<nrows; i++) { 52c4762a1bSJed Brown for (j=0; j<ncols; j++) { 53c4762a1bSJed Brown v[i*ncols+j] = (PetscReal)(10000*rank+100*rows[i]+cols[j]); 54c4762a1bSJed Brown } 55c4762a1bSJed Brown } 56c4762a1bSJed Brown #endif 575f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,nrows,rows,ncols,cols,v,INSERT_VALUES)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(isrows,&rows)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(iscols,&cols)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isrows)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iscols)); 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* Test MatView(), MatDuplicate() and out-of-place MatConvert() */ 665f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(C,MAT_COPY_VALUES,&B)); 67c4762a1bSJed Brown if (mats_view) { 685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Duplicated C:\n")); 695f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 70c4762a1bSJed Brown } 715f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 725f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(C,MATMPIDENSE,MAT_INITIAL_MATRIX,&Cdense)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(C,Cdense,5,&flg)); 74*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Cdense != C. MatConvert() fails"); 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* Test MatNorm() */ 775f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(C,NORM_1,&Cnorm)); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* Test MatTranspose(), MatZeroEntries() and MatGetDiagonal() */ 805f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(C,MAT_INITIAL_MATRIX,&Ct)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(MatConjugate(Ct)); 82c4762a1bSJed Brown if (mats_view) { 835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"C's Transpose Conjugate:\n")); 845f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(Ct,PETSC_VIEWER_STDOUT_WORLD)); 85c4762a1bSJed Brown } 865f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroEntries(Ct)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&d)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(d,m>n ? n : m,PETSC_DECIDE)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(d)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(C,d)); 91c4762a1bSJed Brown if (mats_view) { 925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Diagonal of C:\n")); 935f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(d,PETSC_VIEWER_STDOUT_WORLD)); 94c4762a1bSJed Brown } 95c4762a1bSJed Brown if (m>n) { 965f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(C,NULL,d)); 97c4762a1bSJed Brown } else { 985f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(C,d,NULL)); 99c4762a1bSJed Brown } 100c4762a1bSJed Brown if (mats_view) { 1015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Diagonal Scaled C:\n")); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 103c4762a1bSJed Brown } 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* Test MatAXPY(), MatAYPX() and in-place MatConvert() */ 1065f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B)); 1075f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(B,m,n,PETSC_DECIDE,PETSC_DECIDE)); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(B,MATELEMENTAL)); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(B)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(B)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipIS(B,&isrows,&iscols)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(isrows,&nrows)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(isrows,&rows)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(iscols,&ncols)); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(iscols,&cols)); 116c4762a1bSJed Brown for (i=0; i<nrows; i++) { 117c4762a1bSJed Brown for (j=0; j<ncols; j++) { 118c4762a1bSJed Brown v[i*ncols+j] = (PetscReal)(1000*rows[i]+cols[j]); 119c4762a1bSJed Brown } 120c4762a1bSJed Brown } 1215f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,nrows,rows,ncols,cols,v,INSERT_VALUES)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(v)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(isrows,&rows)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(iscols,&cols)); 1255f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 1275f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(B,2.5,C,SAME_NONZERO_PATTERN)); 1285f80ce2aSJacob Faibussowitsch CHKERRQ(MatAYPX(B,3.75,C,SAME_NONZERO_PATTERN)); 1295f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(B,MATDENSE,MAT_INPLACE_MATRIX,&B)); 130c4762a1bSJed Brown if (mats_view) { 1315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"B after MatAXPY and MatAYPX:\n")); 1325f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 133c4762a1bSJed Brown } 1345f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isrows)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iscols)); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* Test MatMatTransposeMult(): B = C*C^T */ 1395f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatTransposeMult(C,C,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(C,2.0)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatTransposeMult(C,C,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatTransposeMultEqual(C,C,B,10,&flg)); 143*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MatMatTransposeMult: B != C*B^T"); 144c20d7725SJed Brown 145c4762a1bSJed Brown if (mats_view) { 1465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"C MatMatTransposeMult C:\n")); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 148c4762a1bSJed Brown } 149c4762a1bSJed Brown 1505f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Cdense)); 1515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(v)); 1525f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 1535f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 1545f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Ct)); 1555f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&d)); 156c4762a1bSJed Brown ierr = PetscFinalize(); 157c4762a1bSJed Brown return ierr; 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown /*TEST 161c4762a1bSJed Brown 162c4762a1bSJed Brown test: 163c4762a1bSJed Brown nsize: 2 164c4762a1bSJed Brown args: -m 3 -n 2 165c4762a1bSJed Brown requires: elemental 166c4762a1bSJed Brown 167c4762a1bSJed Brown test: 168c4762a1bSJed Brown suffix: 2 169c4762a1bSJed Brown nsize: 6 170c4762a1bSJed Brown args: -m 2 -n 3 171c4762a1bSJed Brown requires: elemental 172c4762a1bSJed Brown 173c20d7725SJed Brown test: 174c20d7725SJed Brown suffix: 3 175c20d7725SJed Brown nsize: 1 176c20d7725SJed Brown args: -m 2 -n 3 177c20d7725SJed Brown requires: elemental 178c20d7725SJed Brown output_file: output/ex39_1.out 179c20d7725SJed Brown 180c4762a1bSJed Brown TEST*/ 181