xref: /petsc/src/dm/tests/ex3.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
1c4762a1bSJed Brown static char help[] = "Solves the 1-dimensional wave equation.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdm.h>
4c4762a1bSJed Brown #include <petscdmda.h>
5c4762a1bSJed Brown #include <petscdraw.h>
6c4762a1bSJed Brown 
main(int argc,char ** argv)7d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
8d71ae5a4SJacob Faibussowitsch {
9c4762a1bSJed Brown   PetscMPIInt  rank, size;
106497c311SBarry Smith   PetscInt     M = 60, time_steps = 100, localsize, j, i, mybase, myend, *localnodes = NULL;
116497c311SBarry Smith   int          width, xbase;
12c4762a1bSJed Brown   DM           da;
13c4762a1bSJed Brown   PetscViewer  viewer, viewer_private;
14c4762a1bSJed Brown   PetscDraw    draw;
15c4762a1bSJed Brown   Vec          local, global;
16c4762a1bSJed Brown   PetscScalar *localptr, *globalptr;
17c4762a1bSJed Brown   PetscReal    a, h, k;
18c4762a1bSJed Brown   PetscBool    flg = PETSC_FALSE;
19c4762a1bSJed Brown 
20327415f7SBarry Smith   PetscFunctionBeginUser;
21c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
24c4762a1bSJed Brown 
259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL));
269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-time", &time_steps, NULL));
27c4762a1bSJed Brown   /*
28c4762a1bSJed Brown       Test putting two nodes on each processor, exact last processor gets the rest
29c4762a1bSJed Brown   */
309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-distribute", &flg, NULL));
31c4762a1bSJed Brown   if (flg) {
329566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &localnodes));
33c4762a1bSJed Brown     for (i = 0; i < size - 1; i++) localnodes[i] = 2;
34c4762a1bSJed Brown     localnodes[size - 1] = M - 2 * (size - 1);
35c4762a1bSJed Brown   }
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   /* Set up the array */
389566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, M, 1, 1, localnodes, &da));
399566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
409566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
419566063dSJacob Faibussowitsch   PetscCall(PetscFree(localnodes));
429566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &global));
439566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(da, &local));
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /* Set up display to show combined wave graph */
469566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD, 0, "Entire Solution", 20, 480, 800, 200, &viewer));
479566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
489566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetDoubleBuffer(draw));
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   /* determine starting point of each processor */
519566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(global, &mybase, &myend));
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /* set up display to show my portion of the wave */
54c4762a1bSJed Brown   xbase = (int)((mybase) * ((800.0 - 4.0 * size) / M) + 4.0 * rank);
55c4762a1bSJed Brown   width = (int)((myend - mybase) * 800. / M);
569566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF, 0, "Local Portion of Solution", xbase, 200, width, 200, &viewer_private));
579566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(viewer_private, 0, &draw));
589566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetDoubleBuffer(draw));
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /* Initialize the array */
619566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(local, &localsize));
629566063dSJacob Faibussowitsch   PetscCall(VecGetArray(global, &globalptr));
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   for (i = 1; i < localsize - 1; i++) {
65c4762a1bSJed Brown     j                = (i - 1) + mybase;
66c4762a1bSJed Brown     globalptr[i - 1] = PetscSinReal((PETSC_PI * j * 6) / ((PetscReal)M) + 1.2 * PetscSinReal((PETSC_PI * j * 2) / ((PetscReal)M))) * 2;
67c4762a1bSJed Brown   }
68c4762a1bSJed Brown 
699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(global, &globalptr));
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   /* Assign Parameters */
72c4762a1bSJed Brown   a = 1.0;
73c4762a1bSJed Brown   h = 1.0 / M;
74c4762a1bSJed Brown   k = h;
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   for (j = 0; j < time_steps; j++) {
77c4762a1bSJed Brown     /* Global to Local */
789566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(da, global, INSERT_VALUES, local));
799566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(da, global, INSERT_VALUES, local));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown     /*Extract local array */
829566063dSJacob Faibussowitsch     PetscCall(VecGetArray(local, &localptr));
839566063dSJacob Faibussowitsch     PetscCall(VecGetArray(global, &globalptr));
84c4762a1bSJed Brown 
85c4762a1bSJed Brown     /* Update Locally - Make array of new values */
86c4762a1bSJed Brown     /* Note: I don't do anything for the first and last entry */
87ad540459SPierre Jolivet     for (i = 1; i < localsize - 1; i++) globalptr[i - 1] = .5 * (localptr[i + 1] + localptr[i - 1]) - (k / (2.0 * a * h)) * (localptr[i + 1] - localptr[i - 1]);
889566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(global, &globalptr));
899566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(local, &localptr));
90c4762a1bSJed Brown 
91c4762a1bSJed Brown     /* View my part of Wave */
929566063dSJacob Faibussowitsch     PetscCall(VecView(global, viewer_private));
93c4762a1bSJed Brown 
94c4762a1bSJed Brown     /* View global Wave */
959566063dSJacob Faibussowitsch     PetscCall(VecView(global, viewer));
96c4762a1bSJed Brown   }
97c4762a1bSJed Brown 
989566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
999566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
1009566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer_private));
1019566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&local));
1029566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&global));
103c4762a1bSJed Brown 
1049566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
105b122ec5aSJacob Faibussowitsch   return 0;
106c4762a1bSJed Brown }
107c4762a1bSJed Brown 
108c4762a1bSJed Brown /*TEST
109c4762a1bSJed Brown 
110c4762a1bSJed Brown     test:
111c4762a1bSJed Brown       nsize: 3
112c4762a1bSJed Brown       args: -time 50 -nox
113c4762a1bSJed Brown       requires: x
114*3886731fSPierre Jolivet       output_file: output/empty.out
115c4762a1bSJed Brown 
116c4762a1bSJed Brown TEST*/
117