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*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 19*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 20*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 21*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 22c4762a1bSJed Brown n = m; 23*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 24*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view)); 25c4762a1bSJed Brown 26*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); 27*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,m,n,PETSC_DECIDE,PETSC_DECIDE)); 28*9566063dSJacob Faibussowitsch PetscCall(MatSetType(C,MATELEMENTAL)); 29*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 30*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 31*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipIS(C,&isrows,&iscols)); 32*9566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrows,&nrows)); 33*9566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrows,&rows)); 34*9566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscols,&ncols)); 35*9566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscols,&cols)); 36*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows*ncols,&v)); 37c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 38c4762a1bSJed Brown PetscRandom rand; 39c4762a1bSJed Brown PetscScalar rval; 40*9566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 41*9566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 42c4762a1bSJed Brown for (i=0; i<nrows; i++) { 43c4762a1bSJed Brown for (j=0; j<ncols; j++) { 44*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&rval)); 45c4762a1bSJed Brown v[i*ncols+j] = rval; 46c4762a1bSJed Brown } 47c4762a1bSJed Brown } 48*9566063dSJacob Faibussowitsch PetscCall(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 56*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(C,nrows,rows,ncols,cols,v,INSERT_VALUES)); 57*9566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrows,&rows)); 58*9566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscols,&cols)); 59*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 60*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 61*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows)); 62*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscols)); 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* Test MatView(), MatDuplicate() and out-of-place MatConvert() */ 65*9566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C,MAT_COPY_VALUES,&B)); 66c4762a1bSJed Brown if (mats_view) { 67*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Duplicated C:\n")); 68*9566063dSJacob Faibussowitsch PetscCall(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 69c4762a1bSJed Brown } 70*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 71*9566063dSJacob Faibussowitsch PetscCall(MatConvert(C,MATMPIDENSE,MAT_INITIAL_MATRIX,&Cdense)); 72*9566063dSJacob Faibussowitsch PetscCall(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() */ 76*9566063dSJacob Faibussowitsch PetscCall(MatNorm(C,NORM_1,&Cnorm)); 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* Test MatTranspose(), MatZeroEntries() and MatGetDiagonal() */ 79*9566063dSJacob Faibussowitsch PetscCall(MatTranspose(C,MAT_INITIAL_MATRIX,&Ct)); 80*9566063dSJacob Faibussowitsch PetscCall(MatConjugate(Ct)); 81c4762a1bSJed Brown if (mats_view) { 82*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"C's Transpose Conjugate:\n")); 83*9566063dSJacob Faibussowitsch PetscCall(MatView(Ct,PETSC_VIEWER_STDOUT_WORLD)); 84c4762a1bSJed Brown } 85*9566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Ct)); 86*9566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&d)); 87*9566063dSJacob Faibussowitsch PetscCall(VecSetSizes(d,m>n ? n : m,PETSC_DECIDE)); 88*9566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(d)); 89*9566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(C,d)); 90c4762a1bSJed Brown if (mats_view) { 91*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Diagonal of C:\n")); 92*9566063dSJacob Faibussowitsch PetscCall(VecView(d,PETSC_VIEWER_STDOUT_WORLD)); 93c4762a1bSJed Brown } 94c4762a1bSJed Brown if (m>n) { 95*9566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(C,NULL,d)); 96c4762a1bSJed Brown } else { 97*9566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(C,d,NULL)); 98c4762a1bSJed Brown } 99c4762a1bSJed Brown if (mats_view) { 100*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Diagonal Scaled C:\n")); 101*9566063dSJacob Faibussowitsch PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 102c4762a1bSJed Brown } 103c4762a1bSJed Brown 104c4762a1bSJed Brown /* Test MatAXPY(), MatAYPX() and in-place MatConvert() */ 105*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); 106*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B,m,n,PETSC_DECIDE,PETSC_DECIDE)); 107*9566063dSJacob Faibussowitsch PetscCall(MatSetType(B,MATELEMENTAL)); 108*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 109*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 110*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipIS(B,&isrows,&iscols)); 111*9566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrows,&nrows)); 112*9566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrows,&rows)); 113*9566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscols,&ncols)); 114*9566063dSJacob Faibussowitsch PetscCall(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 } 120*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,nrows,rows,ncols,cols,v,INSERT_VALUES)); 121*9566063dSJacob Faibussowitsch PetscCall(PetscFree(v)); 122*9566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrows,&rows)); 123*9566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscols,&cols)); 124*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 125*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 126*9566063dSJacob Faibussowitsch PetscCall(MatAXPY(B,2.5,C,SAME_NONZERO_PATTERN)); 127*9566063dSJacob Faibussowitsch PetscCall(MatAYPX(B,3.75,C,SAME_NONZERO_PATTERN)); 128*9566063dSJacob Faibussowitsch PetscCall(MatConvert(B,MATDENSE,MAT_INPLACE_MATRIX,&B)); 129c4762a1bSJed Brown if (mats_view) { 130*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"B after MatAXPY and MatAYPX:\n")); 131*9566063dSJacob Faibussowitsch PetscCall(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 132c4762a1bSJed Brown } 133*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows)); 134*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscols)); 135*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* Test MatMatTransposeMult(): B = C*C^T */ 138*9566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMult(C,C,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B)); 139*9566063dSJacob Faibussowitsch PetscCall(MatScale(C,2.0)); 140*9566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMult(C,C,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B)); 141*9566063dSJacob Faibussowitsch PetscCall(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) { 145*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"C MatMatTransposeMult C:\n")); 146*9566063dSJacob Faibussowitsch PetscCall(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 147c4762a1bSJed Brown } 148c4762a1bSJed Brown 149*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Cdense)); 150*9566063dSJacob Faibussowitsch PetscCall(PetscFree(v)); 151*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 152*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 153*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ct)); 154*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&d)); 155*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 156b122ec5aSJacob 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