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