xref: /petsc/src/dm/tests/ex12.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown /*
3c4762a1bSJed Brown    Simple example to show how PETSc programs can be run from MATLAB.
4c4762a1bSJed Brown   See ex12.m.
5c4762a1bSJed Brown */
6c4762a1bSJed Brown 
7c4762a1bSJed Brown static char help[] = "Solves the one dimensional heat equation.\n\n";
8c4762a1bSJed Brown 
9c4762a1bSJed Brown #include <petscdm.h>
10c4762a1bSJed Brown #include <petscdmda.h>
11c4762a1bSJed Brown 
12*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
13*d71ae5a4SJacob Faibussowitsch {
14c4762a1bSJed Brown   PetscMPIInt  rank, size;
15c4762a1bSJed Brown   PetscInt     M = 14, time_steps = 20, w = 1, s = 1, localsize, j, i, mybase, myend, globalsize;
16c4762a1bSJed Brown   DM           da;
17c4762a1bSJed Brown   Vec          global, local;
18c4762a1bSJed Brown   PetscScalar *globalptr, *localptr;
19c4762a1bSJed Brown   PetscReal    h, k;
20c4762a1bSJed Brown 
21327415f7SBarry Smith   PetscFunctionBeginUser;
229566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL));
249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-time", &time_steps, NULL));
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   /* Set up the array */
279566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, M, w, s, NULL, &da));
289566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
299566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
309566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &global));
319566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
329566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   /* Make copy of local array for doing updates */
359566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(da, &local));
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   /* determine starting point of each processor */
389566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(global, &mybase, &myend));
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   /* Initialize the Array */
419566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(global, &globalsize));
429566063dSJacob Faibussowitsch   PetscCall(VecGetArray(global, &globalptr));
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   for (i = 0; i < globalsize; i++) {
45c4762a1bSJed Brown     j = i + mybase;
46c4762a1bSJed Brown 
47c4762a1bSJed Brown     globalptr[i] = PetscSinReal((PETSC_PI * j * 6) / ((PetscReal)M) + 1.2 * PetscSinReal((PETSC_PI * j * 2) / ((PetscReal)M))) * 4 + 4;
48c4762a1bSJed Brown   }
49c4762a1bSJed Brown 
509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(global, &localptr));
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   /* Assign Parameters */
53c4762a1bSJed Brown   h = 1.0 / M;
54c4762a1bSJed Brown   k = h * h / 2.2;
559566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(local, &localsize));
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   for (j = 0; j < time_steps; j++) {
58c4762a1bSJed Brown     /* Global to Local */
599566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(da, global, INSERT_VALUES, local));
609566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(da, global, INSERT_VALUES, local));
61c4762a1bSJed Brown 
62c4762a1bSJed Brown     /*Extract local array */
639566063dSJacob Faibussowitsch     PetscCall(VecGetArray(local, &localptr));
649566063dSJacob Faibussowitsch     PetscCall(VecGetArray(global, &globalptr));
65c4762a1bSJed Brown 
66c4762a1bSJed Brown     /* Update Locally - Make array of new values */
67c4762a1bSJed Brown     /* Note: I don't do anything for the first and last entry */
68ad540459SPierre Jolivet     for (i = 1; i < localsize - 1; i++) globalptr[i - 1] = localptr[i] + (k / (h * h)) * (localptr[i + 1] - 2.0 * localptr[i] + localptr[i - 1]);
69c4762a1bSJed Brown 
709566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(global, &globalptr));
719566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(local, &localptr));
72c4762a1bSJed Brown 
73c4762a1bSJed Brown     /* View Wave */
74c4762a1bSJed Brown     /* Set Up Display to Show Heat Graph */
75c4762a1bSJed Brown #if defined(PETSC_USE_SOCKET_VIEWER)
769566063dSJacob Faibussowitsch     PetscCall(VecView(global, PETSC_VIEWER_SOCKET_WORLD));
77c4762a1bSJed Brown #endif
78c4762a1bSJed Brown   }
79c4762a1bSJed Brown 
809566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&local));
819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&global));
829566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
839566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
84b122ec5aSJacob Faibussowitsch   return 0;
85c4762a1bSJed Brown }
86