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*327415f7SBarry Smith PetscFunctionBeginUser; 199566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 229566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 23c4762a1bSJed Brown n = m; 249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 259566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view)); 26c4762a1bSJed Brown 279566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); 289566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,m,n,PETSC_DECIDE,PETSC_DECIDE)); 299566063dSJacob Faibussowitsch PetscCall(MatSetType(C,MATELEMENTAL)); 309566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 319566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 329566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipIS(C,&isrows,&iscols)); 339566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrows,&nrows)); 349566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrows,&rows)); 359566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscols,&ncols)); 369566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscols,&cols)); 379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows*ncols,&v)); 38c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 39c4762a1bSJed Brown PetscRandom rand; 40c4762a1bSJed Brown PetscScalar rval; 419566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 429566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 43c4762a1bSJed Brown for (i=0; i<nrows; i++) { 44c4762a1bSJed Brown for (j=0; j<ncols; j++) { 459566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&rval)); 46c4762a1bSJed Brown v[i*ncols+j] = rval; 47c4762a1bSJed Brown } 48c4762a1bSJed Brown } 499566063dSJacob Faibussowitsch PetscCall(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 579566063dSJacob Faibussowitsch PetscCall(MatSetValues(C,nrows,rows,ncols,cols,v,INSERT_VALUES)); 589566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrows,&rows)); 599566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscols,&cols)); 609566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 619566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 629566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows)); 639566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscols)); 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* Test MatView(), MatDuplicate() and out-of-place MatConvert() */ 669566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C,MAT_COPY_VALUES,&B)); 67c4762a1bSJed Brown if (mats_view) { 689566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Duplicated C:\n")); 699566063dSJacob Faibussowitsch PetscCall(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 70c4762a1bSJed Brown } 719566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 729566063dSJacob Faibussowitsch PetscCall(MatConvert(C,MATMPIDENSE,MAT_INITIAL_MATRIX,&Cdense)); 739566063dSJacob Faibussowitsch PetscCall(MatMultEqual(C,Cdense,5,&flg)); 7428b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Cdense != C. MatConvert() fails"); 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* Test MatNorm() */ 779566063dSJacob Faibussowitsch PetscCall(MatNorm(C,NORM_1,&Cnorm)); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* Test MatTranspose(), MatZeroEntries() and MatGetDiagonal() */ 809566063dSJacob Faibussowitsch PetscCall(MatTranspose(C,MAT_INITIAL_MATRIX,&Ct)); 819566063dSJacob Faibussowitsch PetscCall(MatConjugate(Ct)); 82c4762a1bSJed Brown if (mats_view) { 839566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"C's Transpose Conjugate:\n")); 849566063dSJacob Faibussowitsch PetscCall(MatView(Ct,PETSC_VIEWER_STDOUT_WORLD)); 85c4762a1bSJed Brown } 869566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Ct)); 879566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&d)); 889566063dSJacob Faibussowitsch PetscCall(VecSetSizes(d,m>n ? n : m,PETSC_DECIDE)); 899566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(d)); 909566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(C,d)); 91c4762a1bSJed Brown if (mats_view) { 929566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Diagonal of C:\n")); 939566063dSJacob Faibussowitsch PetscCall(VecView(d,PETSC_VIEWER_STDOUT_WORLD)); 94c4762a1bSJed Brown } 95c4762a1bSJed Brown if (m>n) { 969566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(C,NULL,d)); 97c4762a1bSJed Brown } else { 989566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(C,d,NULL)); 99c4762a1bSJed Brown } 100c4762a1bSJed Brown if (mats_view) { 1019566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Diagonal Scaled C:\n")); 1029566063dSJacob Faibussowitsch PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 103c4762a1bSJed Brown } 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* Test MatAXPY(), MatAYPX() and in-place MatConvert() */ 1069566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); 1079566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B,m,n,PETSC_DECIDE,PETSC_DECIDE)); 1089566063dSJacob Faibussowitsch PetscCall(MatSetType(B,MATELEMENTAL)); 1099566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 1109566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 1119566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipIS(B,&isrows,&iscols)); 1129566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrows,&nrows)); 1139566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrows,&rows)); 1149566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscols,&ncols)); 1159566063dSJacob Faibussowitsch PetscCall(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 } 1219566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,nrows,rows,ncols,cols,v,INSERT_VALUES)); 1229566063dSJacob Faibussowitsch PetscCall(PetscFree(v)); 1239566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrows,&rows)); 1249566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscols,&cols)); 1259566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 1269566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 1279566063dSJacob Faibussowitsch PetscCall(MatAXPY(B,2.5,C,SAME_NONZERO_PATTERN)); 1289566063dSJacob Faibussowitsch PetscCall(MatAYPX(B,3.75,C,SAME_NONZERO_PATTERN)); 1299566063dSJacob Faibussowitsch PetscCall(MatConvert(B,MATDENSE,MAT_INPLACE_MATRIX,&B)); 130c4762a1bSJed Brown if (mats_view) { 1319566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"B after MatAXPY and MatAYPX:\n")); 1329566063dSJacob Faibussowitsch PetscCall(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 133c4762a1bSJed Brown } 1349566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows)); 1359566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscols)); 1369566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* Test MatMatTransposeMult(): B = C*C^T */ 1399566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMult(C,C,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B)); 1409566063dSJacob Faibussowitsch PetscCall(MatScale(C,2.0)); 1419566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMult(C,C,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B)); 1429566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMultEqual(C,C,B,10,&flg)); 14328b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MatMatTransposeMult: B != C*B^T"); 144c20d7725SJed Brown 145c4762a1bSJed Brown if (mats_view) { 1469566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"C MatMatTransposeMult C:\n")); 1479566063dSJacob Faibussowitsch PetscCall(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 148c4762a1bSJed Brown } 149c4762a1bSJed Brown 1509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Cdense)); 1519566063dSJacob Faibussowitsch PetscCall(PetscFree(v)); 1529566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1539566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1549566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ct)); 1559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&d)); 1569566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 157b122ec5aSJacob Faibussowitsch return 0; 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