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 6d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckValues(Mat A, PetscBool one) 7d71ae5a4SJacob Faibussowitsch { 8c4762a1bSJed Brown const PetscScalar *array; 9c4762a1bSJed Brown PetscInt M, N, rstart, rend, lda, i, j; 10c4762a1bSJed Brown 11c4762a1bSJed Brown PetscFunctionBegin; 129566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(A, &array)); 139566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(A, &lda)); 149566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 159566063dSJacob Faibussowitsch PetscCall(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]); 2108401ef6SPierre Jolivet PetscCheck(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 } 249566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(A, &array)); 253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 26c4762a1bSJed Brown } 27c4762a1bSJed Brown 28c4762a1bSJed Brown #define CheckValuesIJ(A) CheckValues(A, PETSC_FALSE) 29c4762a1bSJed Brown #define CheckValuesOne(A) CheckValues(A, PETSC_TRUE) 30c4762a1bSJed Brown 31d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 32d71ae5a4SJacob Faibussowitsch { 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 39327415f7SBarry Smith PetscFunctionBeginUser; 409566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, NULL, help)); 41c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(mattype, MATMPIDENSE, sizeof(mattype))); 429566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-mat_type", mattype, sizeof(mattype), NULL)); 43c4762a1bSJed Brown /* 44c4762a1bSJed Brown Create a parallel dense matrix shared by all processors 45c4762a1bSJed Brown */ 469566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, M, N, NULL, &A)); 479566063dSJacob Faibussowitsch PetscCall(MatConvert(A, mattype, MAT_INPLACE_MATRIX, &A)); 48c4762a1bSJed Brown /* 49c4762a1bSJed Brown Set values into the matrix 50c4762a1bSJed Brown */ 51c4762a1bSJed Brown for (i = 0; i < M; i++) { 52c4762a1bSJed Brown for (j = 0; j < N; j++) { 53c4762a1bSJed Brown PetscScalar v = (PetscReal)(1 + i + j * M); 549566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES)); 55c4762a1bSJed Brown } 56c4762a1bSJed Brown } 579566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 589566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 599566063dSJacob Faibussowitsch PetscCall(MatScale(A, 2.0)); 609566063dSJacob Faibussowitsch PetscCall(MatScale(A, 1.0 / 2.0)); 61c4762a1bSJed Brown 62c4762a1bSJed Brown /* 63c4762a1bSJed Brown Store the binary matrix to a file 64c4762a1bSJed Brown */ 659566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_WRITE, &view)); 66c4762a1bSJed Brown for (i = 0; i < 2; i++) { 679566063dSJacob Faibussowitsch PetscCall(MatView(A, view)); 689566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(view, PETSC_VIEWER_NATIVE)); 699566063dSJacob Faibussowitsch PetscCall(MatView(A, view)); 709566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(view)); 71c4762a1bSJed Brown } 729566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 74c4762a1bSJed Brown 75c4762a1bSJed Brown /* 76c4762a1bSJed Brown Now reload the matrix and check its values 77c4762a1bSJed Brown */ 789566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_READ, &view)); 799566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 809566063dSJacob Faibussowitsch PetscCall(MatSetType(A, mattype)); 81c4762a1bSJed Brown for (i = 0; i < 4; i++) { 829566063dSJacob Faibussowitsch if (i > 0) PetscCall(MatZeroEntries(A)); 839566063dSJacob Faibussowitsch PetscCall(MatLoad(A, view)); 849566063dSJacob Faibussowitsch PetscCall(CheckValuesIJ(A)); 85c4762a1bSJed Brown } 869566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 87c4762a1bSJed Brown 889566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((rend - rstart) * N, &array)); 90c4762a1bSJed Brown for (i = 0; i < (rend - rstart) * N; i++) array[i] = (PetscReal)1; 919566063dSJacob Faibussowitsch PetscCall(MatDensePlaceArray(A, array)); 929566063dSJacob Faibussowitsch PetscCall(MatScale(A, 2.0)); 939566063dSJacob Faibussowitsch PetscCall(MatScale(A, 1.0 / 2.0)); 949566063dSJacob Faibussowitsch PetscCall(CheckValuesOne(A)); 959566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_WRITE, &view)); 969566063dSJacob Faibussowitsch PetscCall(MatView(A, view)); 979566063dSJacob Faibussowitsch PetscCall(MatDenseResetArray(A)); 989566063dSJacob Faibussowitsch PetscCall(PetscFree(array)); 999566063dSJacob Faibussowitsch PetscCall(CheckValuesIJ(A)); 1009566063dSJacob Faibussowitsch PetscCall(PetscViewerBinarySetSkipHeader(view, PETSC_TRUE)); 1019566063dSJacob Faibussowitsch PetscCall(MatView(A, view)); 1029566063dSJacob Faibussowitsch PetscCall(PetscViewerBinarySetSkipHeader(view, PETSC_FALSE)); 1039566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 1049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 105c4762a1bSJed Brown 1069566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 1079566063dSJacob Faibussowitsch PetscCall(MatSetType(A, mattype)); 1089566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_READ, &view)); 1099566063dSJacob Faibussowitsch PetscCall(MatLoad(A, view)); 1109566063dSJacob Faibussowitsch PetscCall(CheckValuesOne(A)); 1119566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(A)); 1129566063dSJacob Faibussowitsch PetscCall(PetscViewerBinarySetSkipHeader(view, PETSC_TRUE)); 1139566063dSJacob Faibussowitsch PetscCall(MatLoad(A, view)); 1149566063dSJacob Faibussowitsch PetscCall(PetscViewerBinarySetSkipHeader(view, PETSC_FALSE)); 1159566063dSJacob Faibussowitsch PetscCall(CheckValuesIJ(A)); 1169566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 1179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 118c4762a1bSJed Brown 119c4762a1bSJed Brown { 120c4762a1bSJed Brown PetscInt m = PETSC_DECIDE, n = PETSC_DECIDE; 1219566063dSJacob Faibussowitsch PetscCall(PetscSplitOwnership(PETSC_COMM_WORLD, &m, &M)); 1229566063dSJacob Faibussowitsch PetscCall(PetscSplitOwnership(PETSC_COMM_WORLD, &n, &N)); 123c4762a1bSJed Brown /* TODO: MatCreateDense requires data!=NULL at all processes! */ 1249566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m * N + 1, &array)); 125c4762a1bSJed Brown 1269566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_READ, &view)); 1279566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, m, n, M, N, array, &A)); 1289566063dSJacob Faibussowitsch PetscCall(MatLoad(A, view)); 1299566063dSJacob Faibussowitsch PetscCall(CheckValuesOne(A)); 1309566063dSJacob Faibussowitsch PetscCall(PetscViewerBinarySetSkipHeader(view, PETSC_TRUE)); 1319566063dSJacob Faibussowitsch PetscCall(MatLoad(A, view)); 1329566063dSJacob Faibussowitsch PetscCall(PetscViewerBinarySetSkipHeader(view, PETSC_FALSE)); 1339566063dSJacob Faibussowitsch PetscCall(CheckValuesIJ(A)); 1349566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1359566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 136c4762a1bSJed Brown 1379566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, m, n, M, N, array, &A)); 1389566063dSJacob Faibussowitsch PetscCall(CheckValuesIJ(A)); 1399566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 140c4762a1bSJed Brown 1419566063dSJacob Faibussowitsch PetscCall(PetscFree(array)); 142c4762a1bSJed Brown } 143c4762a1bSJed Brown 1449566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 145b122ec5aSJacob Faibussowitsch return 0; 146c4762a1bSJed Brown } 147c4762a1bSJed Brown 148c4762a1bSJed Brown /*TEST 149c4762a1bSJed Brown 150c4762a1bSJed Brown testset: 151c4762a1bSJed Brown args: -viewer_binary_mpiio 0 152*3886731fSPierre Jolivet output_file: output/empty.out 153c4762a1bSJed Brown test: 154c4762a1bSJed Brown suffix: stdio_1 155c4762a1bSJed Brown nsize: 1 156a5225ed3SStefano Zampini args: -mat_type seqdense 157c4762a1bSJed Brown test: 158c4762a1bSJed Brown suffix: stdio_2 159c4762a1bSJed Brown nsize: 2 160c4762a1bSJed Brown test: 161c4762a1bSJed Brown suffix: stdio_3 162c4762a1bSJed Brown nsize: 3 163c4762a1bSJed Brown test: 164c4762a1bSJed Brown suffix: stdio_4 165c4762a1bSJed Brown nsize: 4 166c4762a1bSJed Brown test: 167c4762a1bSJed Brown suffix: stdio_5 168c4762a1bSJed Brown nsize: 5 169a5225ed3SStefano Zampini test: 170a5225ed3SStefano Zampini requires: cuda 171a5225ed3SStefano Zampini args: -mat_type seqdensecuda 172a5225ed3SStefano Zampini suffix: stdio_cuda_1 173a5225ed3SStefano Zampini nsize: 1 174a5225ed3SStefano Zampini test: 175a5225ed3SStefano Zampini requires: cuda 176a5225ed3SStefano Zampini args: -mat_type mpidensecuda 177a5225ed3SStefano Zampini suffix: stdio_cuda_2 178a5225ed3SStefano Zampini nsize: 2 179a5225ed3SStefano Zampini test: 180a5225ed3SStefano Zampini requires: cuda 181a5225ed3SStefano Zampini args: -mat_type mpidensecuda 182a5225ed3SStefano Zampini suffix: stdio_cuda_3 183a5225ed3SStefano Zampini nsize: 3 184a5225ed3SStefano Zampini test: 185a5225ed3SStefano Zampini requires: cuda 186a5225ed3SStefano Zampini args: -mat_type mpidensecuda 187a5225ed3SStefano Zampini suffix: stdio_cuda_4 188a5225ed3SStefano Zampini nsize: 4 189a5225ed3SStefano Zampini test: 190a5225ed3SStefano Zampini requires: cuda 191a5225ed3SStefano Zampini args: -mat_type mpidensecuda 192a5225ed3SStefano Zampini suffix: stdio_cuda_5 193a5225ed3SStefano Zampini nsize: 5 194c4762a1bSJed Brown 195c4762a1bSJed Brown testset: 196c4762a1bSJed Brown requires: mpiio 197c4762a1bSJed Brown args: -viewer_binary_mpiio 1 198*3886731fSPierre Jolivet output_file: output/empty.out 199c4762a1bSJed Brown test: 200c4762a1bSJed Brown suffix: mpiio_1 201c4762a1bSJed Brown nsize: 1 202c4762a1bSJed Brown test: 203c4762a1bSJed Brown suffix: mpiio_2 204c4762a1bSJed Brown nsize: 2 205c4762a1bSJed Brown test: 206c4762a1bSJed Brown suffix: mpiio_3 207c4762a1bSJed Brown nsize: 3 208c4762a1bSJed Brown test: 209c4762a1bSJed Brown suffix: mpiio_4 210c4762a1bSJed Brown nsize: 4 211c4762a1bSJed Brown test: 212c4762a1bSJed Brown suffix: mpiio_5 213c4762a1bSJed Brown nsize: 5 214a5225ed3SStefano Zampini test: 215a5225ed3SStefano Zampini requires: cuda 216a5225ed3SStefano Zampini args: -mat_type mpidensecuda 217a5225ed3SStefano Zampini suffix: mpiio_cuda_1 218a5225ed3SStefano Zampini nsize: 1 219a5225ed3SStefano Zampini test: 220a5225ed3SStefano Zampini requires: cuda 221a5225ed3SStefano Zampini args: -mat_type mpidensecuda 222a5225ed3SStefano Zampini suffix: mpiio_cuda_2 223a5225ed3SStefano Zampini nsize: 2 224a5225ed3SStefano Zampini test: 225a5225ed3SStefano Zampini requires: cuda 226a5225ed3SStefano Zampini args: -mat_type mpidensecuda 227a5225ed3SStefano Zampini suffix: mpiio_cuda_3 228a5225ed3SStefano Zampini nsize: 3 229a5225ed3SStefano Zampini test: 230a5225ed3SStefano Zampini requires: cuda 231a5225ed3SStefano Zampini args: -mat_type mpidensecuda 232a5225ed3SStefano Zampini suffix: mpiio_cuda_4 233a5225ed3SStefano Zampini nsize: 4 234a5225ed3SStefano Zampini test: 235a5225ed3SStefano Zampini requires: cuda 236a5225ed3SStefano Zampini args: -mat_type mpidensecuda 237a5225ed3SStefano Zampini suffix: mpiio_cuda_5 238a5225ed3SStefano Zampini nsize: 5 239c4762a1bSJed Brown 240c4762a1bSJed Brown TEST*/ 241