xref: /petsc/src/dm/tests/noflux_check.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31) !
118f812f4SBarry Smith static char help[] = "Check to see of DM_BOUNDARY_MIRROR works in 3D for DMDA with star stencil\n";
218f812f4SBarry Smith 
318f812f4SBarry Smith #include "petscdmda.h"
418f812f4SBarry Smith 
518f812f4SBarry Smith /* Contributed by Gourav Kumbhojkar */
618f812f4SBarry Smith 
globalKMat_3d(Mat K,DMDALocalInfo info)718f812f4SBarry Smith PetscErrorCode globalKMat_3d(Mat K, DMDALocalInfo info)
818f812f4SBarry Smith {
918f812f4SBarry Smith   MatStencil  row, col[7];
1018f812f4SBarry Smith   PetscScalar vals[7];
1118f812f4SBarry Smith   PetscInt    ncols;
1218f812f4SBarry Smith 
1318f812f4SBarry Smith   PetscFunctionBeginUser;
1418f812f4SBarry Smith   for (PetscInt i = info.xs; i < info.xs + info.xm; i++) {
1518f812f4SBarry Smith     for (PetscInt j = info.ys; j < info.ys + info.ym; j++) {
1618f812f4SBarry Smith       for (PetscInt k = info.zs; k < info.zs + info.zm; k++) {
1718f812f4SBarry Smith         ncols = 0;
1818f812f4SBarry Smith         row.i = i;
1918f812f4SBarry Smith         row.j = j;
2018f812f4SBarry Smith         row.k = k;
2118f812f4SBarry Smith 
2218f812f4SBarry Smith         col[0].i      = i;
2318f812f4SBarry Smith         col[0].j      = j;
2418f812f4SBarry Smith         col[0].k      = k;
2518f812f4SBarry Smith         vals[ncols++] = -6.; //ncols=1
2618f812f4SBarry Smith 
2718f812f4SBarry Smith         col[ncols].i  = i - 1;
2818f812f4SBarry Smith         col[ncols].j  = j;
2918f812f4SBarry Smith         col[ncols].k  = k;
3018f812f4SBarry Smith         vals[ncols++] = 1.; //ncols=2
3118f812f4SBarry Smith 
3218f812f4SBarry Smith         col[ncols].i  = i + 1;
3318f812f4SBarry Smith         col[ncols].j  = j;
3418f812f4SBarry Smith         col[ncols].k  = k;
3518f812f4SBarry Smith         vals[ncols++] = 1.; //ncols=3
3618f812f4SBarry Smith 
3718f812f4SBarry Smith         col[ncols].i  = i;
3818f812f4SBarry Smith         col[ncols].j  = j - 1;
3918f812f4SBarry Smith         col[ncols].k  = k;
4018f812f4SBarry Smith         vals[ncols++] = 1.; //ncols=4
4118f812f4SBarry Smith 
4218f812f4SBarry Smith         col[ncols].i  = i;
4318f812f4SBarry Smith         col[ncols].j  = j + 1;
4418f812f4SBarry Smith         col[ncols].k  = k;
4518f812f4SBarry Smith         vals[ncols++] = 1.; //ncols=5
4618f812f4SBarry Smith 
4718f812f4SBarry Smith         col[ncols].i  = i;
4818f812f4SBarry Smith         col[ncols].j  = j;
4918f812f4SBarry Smith         col[ncols].k  = k + 1;
5018f812f4SBarry Smith         vals[ncols++] = 1.; //ncols=6
5118f812f4SBarry Smith 
5218f812f4SBarry Smith         col[ncols].i  = i;
5318f812f4SBarry Smith         col[ncols].j  = j;
5418f812f4SBarry Smith         col[ncols].k  = k - 1;
5518f812f4SBarry Smith         vals[ncols++] = 1.; //ncols=7
5618f812f4SBarry Smith 
5718f812f4SBarry Smith         PetscCall(MatSetValuesStencil(K, 1, &row, ncols, col, vals, ADD_VALUES));
5818f812f4SBarry Smith       }
5918f812f4SBarry Smith     }
6018f812f4SBarry Smith   }
6118f812f4SBarry Smith   PetscCall(MatAssemblyBegin(K, MAT_FINAL_ASSEMBLY));
6218f812f4SBarry Smith   PetscCall(MatAssemblyEnd(K, MAT_FINAL_ASSEMBLY));
6318f812f4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
6418f812f4SBarry Smith }
6518f812f4SBarry Smith 
globalKMat_2d(Mat K,DMDALocalInfo info)6618f812f4SBarry Smith PetscErrorCode globalKMat_2d(Mat K, DMDALocalInfo info)
6718f812f4SBarry Smith {
6818f812f4SBarry Smith   MatStencil  row, col[5];
6918f812f4SBarry Smith   PetscScalar vals[5];
7018f812f4SBarry Smith   PetscInt    ncols;
7118f812f4SBarry Smith 
7218f812f4SBarry Smith   PetscFunctionBeginUser;
7318f812f4SBarry Smith   for (PetscInt i = info.xs; i < info.xs + info.xm; i++) {
7418f812f4SBarry Smith     for (PetscInt j = info.ys; j < info.ys + info.ym; j++) {
7518f812f4SBarry Smith       ncols = 0;
7618f812f4SBarry Smith       row.i = i;
7718f812f4SBarry Smith       row.j = j;
7818f812f4SBarry Smith 
7918f812f4SBarry Smith       col[0].i      = i;
8018f812f4SBarry Smith       col[0].j      = j;
8118f812f4SBarry Smith       vals[ncols++] = -4.; //ncols=1
8218f812f4SBarry Smith 
8318f812f4SBarry Smith       col[ncols].i  = i - 1;
8418f812f4SBarry Smith       col[ncols].j  = j;
8518f812f4SBarry Smith       vals[ncols++] = 1.; //ncols=2
8618f812f4SBarry Smith 
8718f812f4SBarry Smith       col[ncols].i  = i;
8818f812f4SBarry Smith       col[ncols].j  = j - 1;
8918f812f4SBarry Smith       vals[ncols++] = 1.; //ncols=3
9018f812f4SBarry Smith 
9118f812f4SBarry Smith       col[ncols].i  = i + 1;
9218f812f4SBarry Smith       col[ncols].j  = j;
9318f812f4SBarry Smith       vals[ncols++] = 1.; //ncols=4
9418f812f4SBarry Smith 
9518f812f4SBarry Smith       col[ncols].i  = i;
9618f812f4SBarry Smith       col[ncols].j  = j + 1;
9718f812f4SBarry Smith       vals[ncols++] = 1.; //ncols=5
9818f812f4SBarry Smith 
9918f812f4SBarry Smith       PetscCall(MatSetValuesStencil(K, 1, &row, ncols, col, vals, ADD_VALUES));
10018f812f4SBarry Smith     }
10118f812f4SBarry Smith   }
10218f812f4SBarry Smith   PetscCall(MatAssemblyBegin(K, MAT_FINAL_ASSEMBLY));
10318f812f4SBarry Smith   PetscCall(MatAssemblyEnd(K, MAT_FINAL_ASSEMBLY));
10418f812f4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
10518f812f4SBarry Smith }
10618f812f4SBarry Smith 
main(int argc,char ** argv)10718f812f4SBarry Smith int main(int argc, char **argv)
10818f812f4SBarry Smith {
10918f812f4SBarry Smith   DM                     da3d, da2d;
11018f812f4SBarry Smith   DMDALocalInfo          info3d, info2d;
11118f812f4SBarry Smith   Mat                    K3d, K2d;
11218f812f4SBarry Smith   PetscInt               ne, num_pts;
11318f812f4SBarry Smith   ISLocalToGlobalMapping ltgm3d, ltgm2d;
11418f812f4SBarry Smith   Vec                    row2d, row3d;
11518f812f4SBarry Smith   PetscReal              norm2d, norm3d;
11618f812f4SBarry Smith 
117*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
11818f812f4SBarry Smith   ne      = 8;
11918f812f4SBarry Smith   num_pts = ne + 1;
12018f812f4SBarry Smith 
12118f812f4SBarry Smith   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_MIRROR, DM_BOUNDARY_MIRROR, DMDA_STENCIL_STAR, num_pts, num_pts + 1, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da2d));
12218f812f4SBarry Smith   PetscCall(DMSetUp(da2d));
12318f812f4SBarry Smith   PetscCall(DMDAGetLocalInfo(da2d, &info2d));
12418f812f4SBarry Smith   PetscCall(DMCreateMatrix(da2d, &K2d));
12518f812f4SBarry Smith   PetscCall(DMGetLocalToGlobalMapping(da2d, &ltgm2d));
12618f812f4SBarry Smith   PetscCall(ISLocalToGlobalMappingView(ltgm2d, PETSC_VIEWER_STDOUT_WORLD));
12718f812f4SBarry Smith   //PetscFinalize();
12818f812f4SBarry Smith   PetscCall(globalKMat_2d(K2d, info2d));
12918f812f4SBarry Smith   PetscCall(MatView(K2d, PETSC_VIEWER_STDOUT_WORLD));
13018f812f4SBarry Smith   PetscCall(MatCreateVecs(K2d, &row2d, NULL));
13118f812f4SBarry Smith 
13218f812f4SBarry Smith   PetscCall(MatGetRowSum(K2d, row2d));
13318f812f4SBarry Smith   PetscCall(VecNorm(row2d, NORM_2, &norm2d));
13418f812f4SBarry Smith 
13518f812f4SBarry Smith   PetscCheck(norm2d == 0, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "2D atrix row sum should be zero");
13618f812f4SBarry Smith   PetscCall(VecDestroy(&row2d));
13718f812f4SBarry Smith   PetscCall(MatDestroy(&K2d));
13818f812f4SBarry Smith   PetscCall(DMDestroy(&da2d));
13918f812f4SBarry Smith 
14018f812f4SBarry Smith   PetscCall(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_MIRROR, DM_BOUNDARY_MIRROR, DM_BOUNDARY_MIRROR, DMDA_STENCIL_STAR, num_pts, num_pts + 1, num_pts + 2, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, NULL, &da3d));
14118f812f4SBarry Smith   PetscCall(DMSetUp(da3d));
14218f812f4SBarry Smith   PetscCall(DMCreateMatrix(da3d, &K3d));
14318f812f4SBarry Smith   PetscCall(DMDAGetLocalInfo(da3d, &info3d));
14418f812f4SBarry Smith   PetscCall(DMGetLocalToGlobalMapping(da3d, &ltgm3d));
14518f812f4SBarry Smith   PetscCall(ISLocalToGlobalMappingView(ltgm3d, PETSC_VIEWER_STDOUT_WORLD));
14618f812f4SBarry Smith   PetscCall(globalKMat_3d(K3d, info3d));
14718f812f4SBarry Smith   PetscCall(MatView(K3d, PETSC_VIEWER_STDOUT_WORLD));
14818f812f4SBarry Smith   PetscCall(MatCreateVecs(K3d, &row3d, NULL));
14918f812f4SBarry Smith   PetscCall(MatGetRowSum(K3d, row3d));
15018f812f4SBarry Smith   PetscCall(VecNorm(row3d, NORM_2, &norm3d));
15118f812f4SBarry Smith   PetscCheck(norm3d == 0, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "3D atrix row sum should be zero");
15218f812f4SBarry Smith   PetscCall(VecDestroy(&row3d));
15318f812f4SBarry Smith 
15418f812f4SBarry Smith   PetscCall(DMDestroy(&da3d));
15518f812f4SBarry Smith   PetscCall(MatDestroy(&K3d));
15618f812f4SBarry Smith   return PetscFinalize();
15718f812f4SBarry Smith }
15818f812f4SBarry Smith 
15918f812f4SBarry Smith /*TEST
16018f812f4SBarry Smith 
16118f812f4SBarry Smith    test:
16218f812f4SBarry Smith 
16318f812f4SBarry Smith    test:
16418f812f4SBarry Smith       suffix: 2
16518f812f4SBarry Smith       nsize: 2
16618f812f4SBarry Smith 
16718f812f4SBarry Smith    test:
16818f812f4SBarry Smith       suffix: 4
16918f812f4SBarry Smith       nsize: 4
17018f812f4SBarry Smith 
17118f812f4SBarry Smith    test:
17218f812f4SBarry Smith       suffix: 8
17318f812f4SBarry Smith       nsize: 8
17418f812f4SBarry Smith 
17518f812f4SBarry Smith TEST*/
176