xref: /petsc/src/dm/tests/ex22.c (revision 92e3ebb5f0dcfe835315f7c1b19ee4ea4682ad38)
1c4762a1bSJed Brown static char help[] = "Tests MatSetValuesBlockedStencil() in 3d.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown #include <petscdm.h>
5c4762a1bSJed Brown #include <petscdmda.h>
6c4762a1bSJed Brown 
7d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
8d71ae5a4SJacob Faibussowitsch {
9c4762a1bSJed Brown   PetscInt        M = 3, N = 4, P = 2, s = 1, w = 2, i, m = PETSC_DECIDE, n = PETSC_DECIDE, p = PETSC_DECIDE;
10c4762a1bSJed Brown   DM              da;
11c4762a1bSJed Brown   Mat             mat;
12c4762a1bSJed Brown   DMDAStencilType stencil_type = DMDA_STENCIL_BOX;
13c4762a1bSJed Brown   PetscBool       flg          = PETSC_FALSE;
14c4762a1bSJed Brown   MatStencil      idx[2], idy[2];
15c4762a1bSJed Brown   PetscScalar    *values;
16c4762a1bSJed Brown 
17327415f7SBarry Smith   PetscFunctionBeginUser;
18c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL));
209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &N, NULL));
219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-P", &P, NULL));
229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-p", &p, NULL));
259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-s", &s, NULL));
269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-w", &w, NULL));
279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-star", &flg, NULL));
28c4762a1bSJed Brown   if (flg) stencil_type = DMDA_STENCIL_STAR;
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   /* Create distributed array and get vectors */
319566063dSJacob Faibussowitsch   PetscCall(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, stencil_type, M, N, P, m, n, p, w, s, 0, 0, 0, &da));
329566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
339566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
349566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(da, &mat));
35*92e3ebb5SRichard Tran Mills   PetscCall(MatSetFromOptions(mat));
36c4762a1bSJed Brown 
379371c9d4SSatish Balay   idx[0].i = 1;
389371c9d4SSatish Balay   idx[0].j = 1;
399371c9d4SSatish Balay   idx[0].k = 0;
409371c9d4SSatish Balay   idx[1].i = 2;
419371c9d4SSatish Balay   idx[1].j = 1;
429371c9d4SSatish Balay   idx[1].k = 0;
439371c9d4SSatish Balay   idy[0].i = 1;
449371c9d4SSatish Balay   idy[0].j = 2;
459371c9d4SSatish Balay   idy[0].k = 0;
469371c9d4SSatish Balay   idy[1].i = 2;
479371c9d4SSatish Balay   idy[1].j = 2;
489371c9d4SSatish Balay   idy[1].k = 0;
499566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * 2 * w * w, &values));
50c4762a1bSJed Brown   for (i = 0; i < 2 * 2 * w * w; i++) values[i] = i;
519566063dSJacob Faibussowitsch   PetscCall(MatSetValuesBlockedStencil(mat, 2, idx, 2, idy, values, INSERT_VALUES));
529566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
539566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
54c4762a1bSJed Brown 
55*92e3ebb5SRichard Tran Mills   flg = PETSC_FALSE;
56*92e3ebb5SRichard Tran Mills   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_matdiagonalscalelocal", &flg, NULL));
57*92e3ebb5SRichard Tran Mills   if (flg) {
58*92e3ebb5SRichard Tran Mills     Vec         vec, scaling;
59*92e3ebb5SRichard Tran Mills     PetscInt    size;
60*92e3ebb5SRichard Tran Mills     PetscMPIInt rank;
61*92e3ebb5SRichard Tran Mills 
62*92e3ebb5SRichard Tran Mills     PetscCall(DMGetLocalVector(da, &vec));
63*92e3ebb5SRichard Tran Mills     /* We need to duplicate vec, since we are not allowed to write to its ghost values, yet we need a vector in which
64*92e3ebb5SRichard Tran Mills        we can write to these locations for testing MatDiagonalScaleLocal(). */
65*92e3ebb5SRichard Tran Mills     PetscCall(VecDuplicate(vec, &scaling));
66*92e3ebb5SRichard Tran Mills     PetscCall(VecGetLocalSize(scaling, &size));
67*92e3ebb5SRichard Tran Mills     PetscCall(PetscFree(values));
68*92e3ebb5SRichard Tran Mills     PetscCall(VecGetArrayWrite(scaling, &values));
69*92e3ebb5SRichard Tran Mills     PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &rank));
70*92e3ebb5SRichard Tran Mills     for (i = 0; i < size; i++) values[i] = (PetscScalar)(rank + 1);
71*92e3ebb5SRichard Tran Mills     PetscCall(VecRestoreArrayWrite(scaling, &values));
72*92e3ebb5SRichard Tran Mills     PetscCall(MatDiagonalScaleLocal(mat, scaling));
73*92e3ebb5SRichard Tran Mills     PetscCall(DMRestoreLocalVector(da, &vec));
74*92e3ebb5SRichard Tran Mills     PetscCall(DMRestoreLocalVector(da, &scaling));
75*92e3ebb5SRichard Tran Mills   }
76*92e3ebb5SRichard Tran Mills 
77c4762a1bSJed Brown   /* Free memory */
789566063dSJacob Faibussowitsch   PetscCall(PetscFree(values));
799566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
809566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
819566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
82b122ec5aSJacob Faibussowitsch   return 0;
83c4762a1bSJed Brown }
84*92e3ebb5SRichard Tran Mills 
85*92e3ebb5SRichard Tran Mills /*TEST
86*92e3ebb5SRichard Tran Mills 
87*92e3ebb5SRichard Tran Mills       test:
88*92e3ebb5SRichard Tran Mills          suffix: baij
89*92e3ebb5SRichard Tran Mills          nsize: 2
90*92e3ebb5SRichard Tran Mills          args: -mat_type baij
91*92e3ebb5SRichard Tran Mills          output_file: output/empty.out
92*92e3ebb5SRichard Tran Mills 
93*92e3ebb5SRichard Tran Mills       test:
94*92e3ebb5SRichard Tran Mills          suffix: aij
95*92e3ebb5SRichard Tran Mills          nsize: 2
96*92e3ebb5SRichard Tran Mills          args: -mat_type aij
97*92e3ebb5SRichard Tran Mills          output_file: output/empty.out
98*92e3ebb5SRichard Tran Mills 
99*92e3ebb5SRichard Tran Mills       test:
100*92e3ebb5SRichard Tran Mills          suffix: baij-diagonalscalelocal
101*92e3ebb5SRichard Tran Mills          nsize: 2
102*92e3ebb5SRichard Tran Mills          args: -mat_type baij -test_matdiagonalscalelocal
103*92e3ebb5SRichard Tran Mills          output_file: output/empty.out
104*92e3ebb5SRichard Tran Mills 
105*92e3ebb5SRichard Tran Mills       test:
106*92e3ebb5SRichard Tran Mills          suffix: aij-diagonalscalelocal
107*92e3ebb5SRichard Tran Mills          nsize: 2
108*92e3ebb5SRichard Tran Mills          args: -mat_type aij -test_matdiagonalscalelocal
109*92e3ebb5SRichard Tran Mills          output_file: output/empty.out
110*92e3ebb5SRichard Tran Mills 
111*92e3ebb5SRichard Tran Mills TEST*/
112