static char help[] = "Tests MatDenseGetArray() and MatView()/MatLoad() with binary viewers.\n\n"; #include #include static PetscErrorCode CheckValues(Mat A,PetscBool one) { const PetscScalar *array; PetscInt M,N,rstart,rend,lda,i,j; PetscFunctionBegin; PetscCall(MatDenseGetArrayRead(A,&array)); PetscCall(MatDenseGetLDA(A,&lda)); PetscCall(MatGetSize(A,&M,&N)); PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); for (i=rstart; i 0) PetscCall(MatZeroEntries(A)); PetscCall(MatLoad(A,view)); PetscCall(CheckValuesIJ(A)); } PetscCall(PetscViewerDestroy(&view)); PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); PetscCall(PetscMalloc1((rend-rstart)*N,&array)); for (i=0; i<(rend-rstart)*N; i++) array[i] = (PetscReal)1; PetscCall(MatDensePlaceArray(A,array)); PetscCall(MatScale(A,2.0)); PetscCall(MatScale(A,1.0/2.0)); PetscCall(CheckValuesOne(A)); PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_WRITE,&view)); PetscCall(MatView(A,view)); PetscCall(MatDenseResetArray(A)); PetscCall(PetscFree(array)); PetscCall(CheckValuesIJ(A)); PetscCall(PetscViewerBinarySetSkipHeader(view,PETSC_TRUE)); PetscCall(MatView(A,view)); PetscCall(PetscViewerBinarySetSkipHeader(view,PETSC_FALSE)); PetscCall(PetscViewerDestroy(&view)); PetscCall(MatDestroy(&A)); PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); PetscCall(MatSetType(A,mattype)); PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view)); PetscCall(MatLoad(A,view)); PetscCall(CheckValuesOne(A)); PetscCall(MatZeroEntries(A)); PetscCall(PetscViewerBinarySetSkipHeader(view,PETSC_TRUE)); PetscCall(MatLoad(A,view)); PetscCall(PetscViewerBinarySetSkipHeader(view,PETSC_FALSE)); PetscCall(CheckValuesIJ(A)); PetscCall(PetscViewerDestroy(&view)); PetscCall(MatDestroy(&A)); { PetscInt m = PETSC_DECIDE, n = PETSC_DECIDE; PetscCall(PetscSplitOwnership(PETSC_COMM_WORLD,&m,&M)); PetscCall(PetscSplitOwnership(PETSC_COMM_WORLD,&n,&N)); /* TODO: MatCreateDense requires data!=NULL at all processes! */ PetscCall(PetscMalloc1(m*N+1,&array)); PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view)); PetscCall(MatCreateDense(PETSC_COMM_WORLD,m,n,M,N,array,&A)); PetscCall(MatLoad(A,view)); PetscCall(CheckValuesOne(A)); PetscCall(PetscViewerBinarySetSkipHeader(view,PETSC_TRUE)); PetscCall(MatLoad(A,view)); PetscCall(PetscViewerBinarySetSkipHeader(view,PETSC_FALSE)); PetscCall(CheckValuesIJ(A)); PetscCall(MatDestroy(&A)); PetscCall(PetscViewerDestroy(&view)); PetscCall(MatCreateDense(PETSC_COMM_WORLD,m,n,M,N,array,&A)); PetscCall(CheckValuesIJ(A)); PetscCall(MatDestroy(&A)); PetscCall(PetscFree(array)); } PetscCall(PetscFinalize()); return 0; } /*TEST testset: args: -viewer_binary_mpiio 0 output_file: output/ex16.out test: suffix: stdio_1 nsize: 1 args: -mat_type seqdense test: suffix: stdio_2 nsize: 2 test: suffix: stdio_3 nsize: 3 test: suffix: stdio_4 nsize: 4 test: suffix: stdio_5 nsize: 5 test: requires: cuda args: -mat_type seqdensecuda suffix: stdio_cuda_1 nsize: 1 test: requires: cuda args: -mat_type mpidensecuda suffix: stdio_cuda_2 nsize: 2 test: requires: cuda args: -mat_type mpidensecuda suffix: stdio_cuda_3 nsize: 3 test: requires: cuda args: -mat_type mpidensecuda suffix: stdio_cuda_4 nsize: 4 test: requires: cuda args: -mat_type mpidensecuda suffix: stdio_cuda_5 nsize: 5 testset: requires: mpiio args: -viewer_binary_mpiio 1 output_file: output/ex16.out test: suffix: mpiio_1 nsize: 1 test: suffix: mpiio_2 nsize: 2 test: suffix: mpiio_3 nsize: 3 test: suffix: mpiio_4 nsize: 4 test: suffix: mpiio_5 nsize: 5 test: requires: cuda args: -mat_type mpidensecuda suffix: mpiio_cuda_1 nsize: 1 test: requires: cuda args: -mat_type mpidensecuda suffix: mpiio_cuda_2 nsize: 2 test: requires: cuda args: -mat_type mpidensecuda suffix: mpiio_cuda_3 nsize: 3 test: requires: cuda args: -mat_type mpidensecuda suffix: mpiio_cuda_4 nsize: 4 test: requires: cuda args: -mat_type mpidensecuda suffix: mpiio_cuda_5 nsize: 5 TEST*/