1c4762a1bSJed Brown static char help[] = "Tests MatDenseGetArray() and MatView()/MatLoad() with binary viewers.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown #include <petscviewer.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown static PetscErrorCode CheckValues(Mat A,PetscBool one) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown const PetscScalar *array; 9c4762a1bSJed Brown PetscInt M,N,rstart,rend,lda,i,j; 10c4762a1bSJed Brown 11c4762a1bSJed Brown PetscFunctionBegin; 125f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArrayRead(A,&array)); 135f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetLDA(A,&lda)); 145f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&M,&N)); 155f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend)); 16c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 17c4762a1bSJed Brown for (j=0; j<N; j++) { 18c4762a1bSJed Brown PetscInt ii = i - rstart, jj = j; 19c4762a1bSJed Brown PetscReal v = (PetscReal)(one ? 1 : (1 + i + j*M)); 20c4762a1bSJed Brown PetscReal w = PetscRealPart(array[ii + jj*lda]); 212c71b3e2SJacob Faibussowitsch PetscCheckFalse(PetscAbsReal(v-w) > 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Matrix entry (%" PetscInt_FMT ",%" PetscInt_FMT ") should be %g, got %g",i,j,(double)v,(double)w); 22c4762a1bSJed Brown } 23c4762a1bSJed Brown } 245f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArrayRead(A,&array)); 25c4762a1bSJed Brown PetscFunctionReturn(0); 26c4762a1bSJed Brown } 27c4762a1bSJed Brown 28c4762a1bSJed Brown #define CheckValuesIJ(A) CheckValues(A,PETSC_FALSE) 29c4762a1bSJed Brown #define CheckValuesOne(A) CheckValues(A,PETSC_TRUE) 30c4762a1bSJed Brown 31c4762a1bSJed Brown int main(int argc,char **args) 32c4762a1bSJed Brown { 33c4762a1bSJed Brown Mat A; 34c4762a1bSJed Brown PetscInt i,j,M = 4,N = 3,rstart,rend; 35c4762a1bSJed Brown PetscScalar *array; 36a5225ed3SStefano Zampini char mattype[256]; 37c4762a1bSJed Brown PetscViewer view; 38c4762a1bSJed Brown 39*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,NULL,help)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcpy(mattype,MATMPIDENSE)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-mat_type",mattype,sizeof(mattype),NULL)); 42c4762a1bSJed Brown /* 43c4762a1bSJed Brown Create a parallel dense matrix shared by all processors 44c4762a1bSJed Brown */ 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,M,N,NULL,&A)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 485f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,mattype,MAT_INPLACE_MATRIX,&A)); 49c4762a1bSJed Brown /* 50c4762a1bSJed Brown Set values into the matrix 51c4762a1bSJed Brown */ 52c4762a1bSJed Brown for (i=0; i<M; i++) { 53c4762a1bSJed Brown for (j=0; j<N; j++) { 54c4762a1bSJed Brown PetscScalar v = (PetscReal)(1 + i + j*M); 555f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES)); 56c4762a1bSJed Brown } 57c4762a1bSJed Brown } 585f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,2.0)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,1.0/2.0)); 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* 64c4762a1bSJed Brown Store the binary matrix to a file 65c4762a1bSJed Brown */ 665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_WRITE, &view)); 67c4762a1bSJed Brown for (i=0; i<2; i++) { 685f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,view)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(view,PETSC_VIEWER_NATIVE)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,view)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPopFormat(view)); 72c4762a1bSJed Brown } 735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&view)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* 77c4762a1bSJed Brown Now reload the matrix and check its values 78c4762a1bSJed Brown */ 795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,mattype)); 82c4762a1bSJed Brown for (i=0; i<4; i++) { 835f80ce2aSJacob Faibussowitsch if (i > 0) CHKERRQ(MatZeroEntries(A)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,view)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(CheckValuesIJ(A)); 86c4762a1bSJed Brown } 875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&view)); 88c4762a1bSJed Brown 895f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1((rend-rstart)*N,&array)); 91c4762a1bSJed Brown for (i=0; i<(rend-rstart)*N; i++) array[i] = (PetscReal)1; 925f80ce2aSJacob Faibussowitsch CHKERRQ(MatDensePlaceArray(A,array)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,2.0)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,1.0/2.0)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(CheckValuesOne(A)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_WRITE,&view)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,view)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseResetArray(A)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(array)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(CheckValuesIJ(A)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinarySetSkipHeader(view,PETSC_TRUE)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,view)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinarySetSkipHeader(view,PETSC_FALSE)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&view)); 1055f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 106c4762a1bSJed Brown 1075f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,mattype)); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,view)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(CheckValuesOne(A)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroEntries(A)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinarySetSkipHeader(view,PETSC_TRUE)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,view)); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinarySetSkipHeader(view,PETSC_FALSE)); 1165f80ce2aSJacob Faibussowitsch CHKERRQ(CheckValuesIJ(A)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&view)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 119c4762a1bSJed Brown 120c4762a1bSJed Brown { 121c4762a1bSJed Brown PetscInt m = PETSC_DECIDE, n = PETSC_DECIDE; 1225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSplitOwnership(PETSC_COMM_WORLD,&m,&M)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSplitOwnership(PETSC_COMM_WORLD,&n,&N)); 124c4762a1bSJed Brown /* TODO: MatCreateDense requires data!=NULL at all processes! */ 1255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(m*N+1,&array)); 126c4762a1bSJed Brown 1275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view)); 1285f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,m,n,M,N,array,&A)); 1295f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,view)); 1305f80ce2aSJacob Faibussowitsch CHKERRQ(CheckValuesOne(A)); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinarySetSkipHeader(view,PETSC_TRUE)); 1325f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,view)); 1335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinarySetSkipHeader(view,PETSC_FALSE)); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(CheckValuesIJ(A)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&view)); 137c4762a1bSJed Brown 1385f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,m,n,M,N,array,&A)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(CheckValuesIJ(A)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 141c4762a1bSJed Brown 1425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(array)); 143c4762a1bSJed Brown } 144c4762a1bSJed Brown 145*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 146*b122ec5aSJacob Faibussowitsch return 0; 147c4762a1bSJed Brown } 148c4762a1bSJed Brown 149c4762a1bSJed Brown /*TEST 150c4762a1bSJed Brown 151c4762a1bSJed Brown testset: 152c4762a1bSJed Brown args: -viewer_binary_mpiio 0 153c4762a1bSJed Brown output_file: output/ex16.out 154c4762a1bSJed Brown test: 155c4762a1bSJed Brown suffix: stdio_1 156c4762a1bSJed Brown nsize: 1 157a5225ed3SStefano Zampini args: -mat_type seqdense 158c4762a1bSJed Brown test: 159c4762a1bSJed Brown suffix: stdio_2 160c4762a1bSJed Brown nsize: 2 161c4762a1bSJed Brown test: 162c4762a1bSJed Brown suffix: stdio_3 163c4762a1bSJed Brown nsize: 3 164c4762a1bSJed Brown test: 165c4762a1bSJed Brown suffix: stdio_4 166c4762a1bSJed Brown nsize: 4 167c4762a1bSJed Brown test: 168c4762a1bSJed Brown suffix: stdio_5 169c4762a1bSJed Brown nsize: 5 170a5225ed3SStefano Zampini test: 171a5225ed3SStefano Zampini requires: cuda 172a5225ed3SStefano Zampini args: -mat_type seqdensecuda 173a5225ed3SStefano Zampini suffix: stdio_cuda_1 174a5225ed3SStefano Zampini nsize: 1 175a5225ed3SStefano Zampini test: 176a5225ed3SStefano Zampini requires: cuda 177a5225ed3SStefano Zampini args: -mat_type mpidensecuda 178a5225ed3SStefano Zampini suffix: stdio_cuda_2 179a5225ed3SStefano Zampini nsize: 2 180a5225ed3SStefano Zampini test: 181a5225ed3SStefano Zampini requires: cuda 182a5225ed3SStefano Zampini args: -mat_type mpidensecuda 183a5225ed3SStefano Zampini suffix: stdio_cuda_3 184a5225ed3SStefano Zampini nsize: 3 185a5225ed3SStefano Zampini test: 186a5225ed3SStefano Zampini requires: cuda 187a5225ed3SStefano Zampini args: -mat_type mpidensecuda 188a5225ed3SStefano Zampini suffix: stdio_cuda_4 189a5225ed3SStefano Zampini nsize: 4 190a5225ed3SStefano Zampini test: 191a5225ed3SStefano Zampini requires: cuda 192a5225ed3SStefano Zampini args: -mat_type mpidensecuda 193a5225ed3SStefano Zampini suffix: stdio_cuda_5 194a5225ed3SStefano Zampini nsize: 5 195c4762a1bSJed Brown 196c4762a1bSJed Brown testset: 197c4762a1bSJed Brown requires: mpiio 198c4762a1bSJed Brown args: -viewer_binary_mpiio 1 199c4762a1bSJed Brown output_file: output/ex16.out 200c4762a1bSJed Brown test: 201c4762a1bSJed Brown suffix: mpiio_1 202c4762a1bSJed Brown nsize: 1 203c4762a1bSJed Brown test: 204c4762a1bSJed Brown suffix: mpiio_2 205c4762a1bSJed Brown nsize: 2 206c4762a1bSJed Brown test: 207c4762a1bSJed Brown suffix: mpiio_3 208c4762a1bSJed Brown nsize: 3 209c4762a1bSJed Brown test: 210c4762a1bSJed Brown suffix: mpiio_4 211c4762a1bSJed Brown nsize: 4 212c4762a1bSJed Brown test: 213c4762a1bSJed Brown suffix: mpiio_5 214c4762a1bSJed Brown nsize: 5 215a5225ed3SStefano Zampini test: 216a5225ed3SStefano Zampini requires: cuda 217a5225ed3SStefano Zampini args: -mat_type mpidensecuda 218a5225ed3SStefano Zampini suffix: mpiio_cuda_1 219a5225ed3SStefano Zampini nsize: 1 220a5225ed3SStefano Zampini test: 221a5225ed3SStefano Zampini requires: cuda 222a5225ed3SStefano Zampini args: -mat_type mpidensecuda 223a5225ed3SStefano Zampini suffix: mpiio_cuda_2 224a5225ed3SStefano Zampini nsize: 2 225a5225ed3SStefano Zampini test: 226a5225ed3SStefano Zampini requires: cuda 227a5225ed3SStefano Zampini args: -mat_type mpidensecuda 228a5225ed3SStefano Zampini suffix: mpiio_cuda_3 229a5225ed3SStefano Zampini nsize: 3 230a5225ed3SStefano Zampini test: 231a5225ed3SStefano Zampini requires: cuda 232a5225ed3SStefano Zampini args: -mat_type mpidensecuda 233a5225ed3SStefano Zampini suffix: mpiio_cuda_4 234a5225ed3SStefano Zampini nsize: 4 235a5225ed3SStefano Zampini test: 236a5225ed3SStefano Zampini requires: cuda 237a5225ed3SStefano Zampini args: -mat_type mpidensecuda 238a5225ed3SStefano Zampini suffix: mpiio_cuda_5 239a5225ed3SStefano Zampini nsize: 5 240c4762a1bSJed Brown 241c4762a1bSJed Brown TEST*/ 242