126cec326SBarry Smith static char help[] = "Test of setting values in a matrix without preallocation\n\n"; 226cec326SBarry Smith 326cec326SBarry Smith #include <petscmat.h> 426cec326SBarry Smith 526cec326SBarry Smith PetscErrorCode ex1_nonsquare_bs1(void) 626cec326SBarry Smith { 726cec326SBarry Smith Mat A; 826cec326SBarry Smith PetscInt M, N, m, n, bs = 1; 926cec326SBarry Smith char type[16]; 1026cec326SBarry Smith PetscBool flg; 1126cec326SBarry Smith 1226cec326SBarry Smith /* 1326cec326SBarry Smith Create the matrix 1426cec326SBarry Smith */ 1526cec326SBarry Smith PetscFunctionBegin; 1626cec326SBarry Smith M = 10; 1726cec326SBarry Smith N = 12; 1826cec326SBarry Smith PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 1926cec326SBarry Smith PetscCall(PetscOptionsGetString(NULL, NULL, "-type", type, sizeof(type), &flg)); 2026cec326SBarry Smith if (flg) PetscCall(MatSetType(A, type)); 2126cec326SBarry Smith else PetscCall(MatSetType(A, MATAIJ)); 2226cec326SBarry Smith PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, M, N)); 2326cec326SBarry Smith PetscCall(PetscOptionsGetInt(NULL, NULL, "-bs", &bs, NULL)); 2426cec326SBarry Smith PetscCall(MatSetBlockSize(A, bs)); 2526cec326SBarry Smith PetscCall(MatSetFromOptions(A)); 2626cec326SBarry Smith 2726cec326SBarry Smith /* 2826cec326SBarry Smith Get the sizes of the matrix 2926cec326SBarry Smith */ 3026cec326SBarry Smith PetscCall(MatGetLocalSize(A, &m, &n)); 3126cec326SBarry Smith 3226cec326SBarry Smith /* 3326cec326SBarry Smith Insert non-zero pattern (e.g. perform a sweep over the grid). 3426cec326SBarry Smith You can use MatSetValues(), MatSetValuesBlocked() or MatSetValue(). 3526cec326SBarry Smith */ 3626cec326SBarry Smith { 3726cec326SBarry Smith PetscInt ii, jj; 3826cec326SBarry Smith PetscScalar vv = 22.0; 3926cec326SBarry Smith 4026cec326SBarry Smith ii = 3; 4126cec326SBarry Smith jj = 3; 4226cec326SBarry Smith PetscCall(MatSetValue(A, ii, jj, vv, INSERT_VALUES)); 4326cec326SBarry Smith 4426cec326SBarry Smith ii = 7; 4526cec326SBarry Smith jj = 4; 4626cec326SBarry Smith PetscCall(MatSetValue(A, ii, jj, vv, INSERT_VALUES)); 4726cec326SBarry Smith PetscCall(MatSetValue(A, jj, ii, vv, INSERT_VALUES)); 4826cec326SBarry Smith 4926cec326SBarry Smith ii = 9; 5026cec326SBarry Smith jj = 7; 5126cec326SBarry Smith PetscCall(MatSetValue(A, ii, jj, vv, INSERT_VALUES)); 5226cec326SBarry Smith PetscCall(MatSetValue(A, jj, ii, vv, INSERT_VALUES)); 5326cec326SBarry Smith } 5426cec326SBarry Smith PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 5526cec326SBarry Smith PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 5626cec326SBarry Smith 5726cec326SBarry Smith PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_COMMON)); 5826cec326SBarry Smith PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 5926cec326SBarry Smith 6026cec326SBarry Smith /* 6126cec326SBarry Smith Insert same location non-zero values into A. 6226cec326SBarry Smith */ 6326cec326SBarry Smith { 6426cec326SBarry Smith PetscInt ii, jj; 6526cec326SBarry Smith PetscScalar vv; 6626cec326SBarry Smith 6726cec326SBarry Smith ii = 3; 6826cec326SBarry Smith jj = 3; 6926cec326SBarry Smith vv = 0.3; 7026cec326SBarry Smith PetscCall(MatSetValue(A, ii, jj, vv, INSERT_VALUES)); 7126cec326SBarry Smith 7226cec326SBarry Smith ii = 7; 7326cec326SBarry Smith jj = 4; 7426cec326SBarry Smith vv = 3.3; 7526cec326SBarry Smith PetscCall(MatSetValue(A, ii, jj, vv, INSERT_VALUES)); 7626cec326SBarry Smith PetscCall(MatSetValue(A, jj, ii, vv, INSERT_VALUES)); 7726cec326SBarry Smith 7826cec326SBarry Smith ii = 9; 7926cec326SBarry Smith jj = 7; 8026cec326SBarry Smith vv = 4.3; 8126cec326SBarry Smith PetscCall(MatSetValue(A, ii, jj, vv, INSERT_VALUES)); 8226cec326SBarry Smith PetscCall(MatSetValue(A, jj, ii, vv, INSERT_VALUES)); 8326cec326SBarry Smith } 8426cec326SBarry Smith PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 8526cec326SBarry Smith PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 8626cec326SBarry Smith 8726cec326SBarry Smith PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 8826cec326SBarry Smith 8926cec326SBarry Smith PetscCall(MatDestroy(&A)); 9026cec326SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 9126cec326SBarry Smith } 9226cec326SBarry Smith 9326cec326SBarry Smith int main(int argc, char **args) 9426cec326SBarry Smith { 9526cec326SBarry Smith PetscFunctionBeginUser; 96*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 9726cec326SBarry Smith PetscCall(ex1_nonsquare_bs1()); 9826cec326SBarry Smith PetscCall(PetscFinalize()); 9926cec326SBarry Smith return 0; 10026cec326SBarry Smith } 10126cec326SBarry Smith 10226cec326SBarry Smith /*TEST 10326cec326SBarry Smith 10426cec326SBarry Smith testset: 10526cec326SBarry Smith args: -bs {{1 2}} -type {{aij baij sbaij}} 10626cec326SBarry Smith filter: grep -v "type:" 10726cec326SBarry Smith test: 10826cec326SBarry Smith test: 10926cec326SBarry Smith suffix: 2 11026cec326SBarry Smith nsize: 2 11126cec326SBarry Smith 11226cec326SBarry Smith TEST*/ 113