xref: /petsc/src/mat/tests/ex265.c (revision f13dfd9ea68e0ddeee984e65c377a1819eab8a8a) !
1 static char help[] = "Tests inserting new block into SBAIJ and BAIJ matrix \n ";
2 
3 #include <petscdmda.h>
4 
5 int main(int argc, char **argv)
6 {
7   DM          dm;
8   Mat         A;
9   PetscInt    idm = 0, idn = 8;
10   PetscScalar v[] = {1, 2, 3, 4};
11 
12   PetscFunctionBeginUser;
13   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
14   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, 4, 4, PETSC_DECIDE, PETSC_DECIDE, 2, 1, NULL, NULL, &dm));
15   PetscCall(DMSetFromOptions(dm));
16   PetscCall(DMSetUp(dm));
17   PetscCall(DMCreateMatrix(dm, &A));
18   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_FALSE));
19   PetscCall(MatSetValuesBlocked(A, 1, &idm, 1, &idn, v, INSERT_VALUES));
20   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
21   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
22   PetscCall(MatDestroy(&A));
23   PetscCall(DMDestroy(&dm));
24   PetscCall(PetscFinalize());
25   return 0;
26 }
27 
28 /*TEST
29 
30    test:
31       args: -dm_mat_type {{aij baij sbaij}separate output} -mat_view
32 
33    test:
34      suffix: 2
35      nsize: 2
36      args: -dm_mat_type {{aij baij sbaij}separate output} -mat_view
37 
38    test:
39      suffix: 3
40      nsize: 3
41      args: -dm_mat_type {{aij baij sbaij}separate output} -mat_view
42 
43 TEST*/
44