xref: /petsc/src/ts/tutorials/ex35.cxx (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown static const char help[] = "Time-dependent Brusselator reaction-diffusion PDE in 1d. Demonstrates IMEX methods and uses MOAB.\n";
2c4762a1bSJed Brown /*
3c4762a1bSJed Brown    u_t - alpha u_xx = A + u^2 v - (B+1) u
4c4762a1bSJed Brown    v_t - alpha v_xx = B u - u^2 v
5c4762a1bSJed Brown    0 < x < 1;
6c4762a1bSJed Brown    A = 1, B = 3, alpha = 1/50
7c4762a1bSJed Brown 
8c4762a1bSJed Brown    Initial conditions:
9c4762a1bSJed Brown    u(x,0) = 1 + sin(2 pi x)
10c4762a1bSJed Brown    v(x,0) = 3
11c4762a1bSJed Brown 
12c4762a1bSJed Brown    Boundary conditions:
13c4762a1bSJed Brown    u(0,t) = u(1,t) = 1
14c4762a1bSJed Brown    v(0,t) = v(1,t) = 3
15c4762a1bSJed Brown */
16c4762a1bSJed Brown 
17c4762a1bSJed Brown // PETSc includes:
18c4762a1bSJed Brown #include <petscts.h>
19c4762a1bSJed Brown #include <petscdmmoab.h>
20c4762a1bSJed Brown 
21c4762a1bSJed Brown typedef struct {
22c4762a1bSJed Brown   PetscScalar u, v;
23c4762a1bSJed Brown } Field;
24c4762a1bSJed Brown 
25c4762a1bSJed Brown struct pUserCtx {
26c4762a1bSJed Brown   PetscReal A, B;    /* Reaction coefficients */
27c4762a1bSJed Brown   PetscReal alpha;   /* Diffusion coefficient */
28c4762a1bSJed Brown   Field     leftbc;  /* Dirichlet boundary conditions at left boundary */
29c4762a1bSJed Brown   Field     rightbc; /* Dirichlet boundary conditions at right boundary */
30c4762a1bSJed Brown   PetscInt  n, npts; /* Number of mesh points */
31c4762a1bSJed Brown   PetscInt  ntsteps; /* Number of time steps */
32c4762a1bSJed Brown   PetscInt  nvars;   /* Number of variables in the equation system */
33c4762a1bSJed Brown   PetscBool io;
34c4762a1bSJed Brown };
35c4762a1bSJed Brown typedef pUserCtx *UserCtx;
36c4762a1bSJed Brown 
37*9371c9d4SSatish Balay PetscErrorCode Initialize_AppContext(UserCtx *puser) {
38c4762a1bSJed Brown   UserCtx user;
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   PetscFunctionBegin;
419566063dSJacob Faibussowitsch   PetscCall(PetscNew(&user));
42d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Advection-reaction options", "ex35.cxx");
43c4762a1bSJed Brown   {
44c4762a1bSJed Brown     user->nvars     = 2;
45c4762a1bSJed Brown     user->A         = 1;
46c4762a1bSJed Brown     user->B         = 3;
47c4762a1bSJed Brown     user->alpha     = 0.02;
48c4762a1bSJed Brown     user->leftbc.u  = 1;
49c4762a1bSJed Brown     user->rightbc.u = 1;
50c4762a1bSJed Brown     user->leftbc.v  = 3;
51c4762a1bSJed Brown     user->rightbc.v = 3;
52c4762a1bSJed Brown     user->n         = 10;
53c4762a1bSJed Brown     user->ntsteps   = 10000;
54c4762a1bSJed Brown     user->io        = PETSC_FALSE;
559566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-A", "Reaction rate", "ex35.cxx", user->A, &user->A, NULL));
569566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-B", "Reaction rate", "ex35.cxx", user->B, &user->B, NULL));
579566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-alpha", "Diffusion coefficient", "ex35.cxx", user->alpha, &user->alpha, NULL));
589566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-uleft", "Dirichlet boundary condition", "ex35.cxx", user->leftbc.u, &user->leftbc.u, NULL));
599566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-uright", "Dirichlet boundary condition", "ex35.cxx", user->rightbc.u, &user->rightbc.u, NULL));
609566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-vleft", "Dirichlet boundary condition", "ex35.cxx", user->leftbc.v, &user->leftbc.v, NULL));
619566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-vright", "Dirichlet boundary condition", "ex35.cxx", user->rightbc.v, &user->rightbc.v, NULL));
629566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-n", "Number of 1-D elements", "ex35.cxx", user->n, &user->n, NULL));
639566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-ndt", "Number of time steps", "ex35.cxx", user->ntsteps, &user->ntsteps, NULL));
649566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-io", "Write the mesh and solution output to a file.", "ex35.cxx", user->io, &user->io, NULL));
65c4762a1bSJed Brown     user->npts = user->n + 1;
66c4762a1bSJed Brown   }
67d0609cedSBarry Smith   PetscOptionsEnd();
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   *puser = user;
70c4762a1bSJed Brown   PetscFunctionReturn(0);
71c4762a1bSJed Brown }
72c4762a1bSJed Brown 
73*9371c9d4SSatish Balay PetscErrorCode Destroy_AppContext(UserCtx *user) {
74c4762a1bSJed Brown   PetscFunctionBegin;
759566063dSJacob Faibussowitsch   PetscCall(PetscFree(*user));
76c4762a1bSJed Brown   PetscFunctionReturn(0);
77c4762a1bSJed Brown }
78c4762a1bSJed Brown 
79c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS, Vec, void *);
80c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS, PetscReal, Vec, Vec, void *);
81c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS, PetscReal, Vec, Vec, Vec, void *);
82c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
83c4762a1bSJed Brown 
84c4762a1bSJed Brown /****************
85c4762a1bSJed Brown  *              *
86c4762a1bSJed Brown  *     MAIN     *
87c4762a1bSJed Brown  *              *
88c4762a1bSJed Brown  ****************/
89*9371c9d4SSatish Balay int main(int argc, char **argv) {
90c4762a1bSJed Brown   TS                ts; /* nonlinear solver */
91c4762a1bSJed Brown   Vec               X;  /* solution, residual vectors */
92c4762a1bSJed Brown   Mat               J;  /* Jacobian matrix */
93c4762a1bSJed Brown   PetscInt          steps, mx;
94c4762a1bSJed Brown   PetscReal         hx, dt, ftime;
95c4762a1bSJed Brown   UserCtx           user; /* user-defined work context */
96c4762a1bSJed Brown   TSConvergedReason reason;
97c4762a1bSJed Brown   DM                dm;
98c4762a1bSJed Brown   const char       *fields[2] = {"U", "V"};
99c4762a1bSJed Brown 
100327415f7SBarry Smith   PetscFunctionBeginUser;
1019566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   /* Initialize the user context struct */
1049566063dSJacob Faibussowitsch   PetscCall(Initialize_AppContext(&user));
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /* Fill in the user defined work context: */
1079566063dSJacob Faibussowitsch   PetscCall(DMMoabCreateBoxMesh(PETSC_COMM_WORLD, 1, PETSC_FALSE, NULL, user->n, 1, &dm));
1089566063dSJacob Faibussowitsch   PetscCall(DMMoabSetFieldNames(dm, user->nvars, fields));
1099566063dSJacob Faibussowitsch   PetscCall(DMMoabSetBlockSize(dm, user->nvars));
1109566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   /* SetUp the data structures for DMMOAB */
1139566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dm));
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /*  Create timestepping solver context */
1169566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1179566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, dm));
1189566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSARKIMEX));
1199566063dSJacob Faibussowitsch   PetscCall(TSSetEquationType(ts, TS_EQ_DAE_IMPLICIT_INDEX1));
1209566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(dm, MATBAIJ));
1219566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(dm, &J));
122c4762a1bSJed Brown 
1239566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, FormRHSFunction, user));
1249566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, FormIFunction, user));
1259566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, J, J, FormIJacobian, user));
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   ftime = 10.0;
1289566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, user->ntsteps));
1299566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, ftime));
1309566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133c4762a1bSJed Brown      Create the solution vector and set the initial conditions
134c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1359566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &X));
136c4762a1bSJed Brown 
1379566063dSJacob Faibussowitsch   PetscCall(FormInitialSolution(ts, X, user));
1389566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, X));
1399566063dSJacob Faibussowitsch   PetscCall(VecGetSize(X, &mx));
140c4762a1bSJed Brown   hx = 1.0 / (PetscReal)(mx / 2 - 1);
141c4762a1bSJed Brown   dt = 0.4 * PetscSqr(hx) / user->alpha; /* Diffusive stability limit */
1429566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145c4762a1bSJed Brown      Set runtime options
146c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1479566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150c4762a1bSJed Brown      Solve nonlinear system
151c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1529566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, X));
1539566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
1549566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
1559566063dSJacob Faibussowitsch   PetscCall(TSGetConvergedReason(ts, &reason));
15663a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %" PetscInt_FMT " steps\n", TSConvergedReasons[reason], (double)ftime, steps));
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   if (user->io) {
159c4762a1bSJed Brown     /* Print the numerical solution to screen and then dump to file */
1609566063dSJacob Faibussowitsch     PetscCall(VecView(X, PETSC_VIEWER_STDOUT_WORLD));
161c4762a1bSJed Brown 
162c4762a1bSJed Brown     /* Write out the solution along with the mesh */
1639566063dSJacob Faibussowitsch     PetscCall(DMMoabSetGlobalFieldVector(dm, X));
164c4762a1bSJed Brown #ifdef MOAB_HAVE_HDF5
1659566063dSJacob Faibussowitsch     PetscCall(DMMoabOutput(dm, "ex35.h5m", ""));
166c4762a1bSJed Brown #else
167c4762a1bSJed Brown     /* MOAB does not support true parallel writers that aren't HDF5 based
168c4762a1bSJed Brown        And so if you are using VTK as the output format in parallel,
169c4762a1bSJed Brown        the data could be jumbled due to the order in which the processors
170c4762a1bSJed Brown        write out their parts of the mesh and solution tags
171c4762a1bSJed Brown     */
1729566063dSJacob Faibussowitsch     PetscCall(DMMoabOutput(dm, "ex35.vtk", ""));
173c4762a1bSJed Brown #endif
174c4762a1bSJed Brown   }
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   /* Free work space.
177c4762a1bSJed Brown      Free all PETSc related resources: */
1789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
1799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
1809566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1819566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
182c4762a1bSJed Brown 
183c4762a1bSJed Brown   /* Free all MOAB related resources: */
1849566063dSJacob Faibussowitsch   PetscCall(Destroy_AppContext(&user));
1859566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
186b122ec5aSJacob Faibussowitsch   return 0;
187c4762a1bSJed Brown }
188c4762a1bSJed Brown 
189c4762a1bSJed Brown /*
190c4762a1bSJed Brown   IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
191c4762a1bSJed Brown */
192*9371c9d4SSatish Balay PetscErrorCode FormIJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat J, Mat Jpre, void *ptr) {
193c4762a1bSJed Brown   UserCtx            user = (UserCtx)ptr;
194c4762a1bSJed Brown   PetscInt           dof;
195c4762a1bSJed Brown   PetscReal          hx;
196c4762a1bSJed Brown   DM                 dm;
197c4762a1bSJed Brown   const moab::Range *vlocal;
198c4762a1bSJed Brown   PetscBool          vonboundary;
199c4762a1bSJed Brown 
200c4762a1bSJed Brown   PetscFunctionBegin;
2019566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   /* get the essential MOAB mesh related quantities needed for FEM assembly */
2049566063dSJacob Faibussowitsch   PetscCall(DMMoabGetLocalVertices(dm, &vlocal, NULL));
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   /* compute local element sizes - structured grid */
207c4762a1bSJed Brown   hx = 1.0 / user->n;
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid
210c4762a1bSJed Brown      Assemble the operator by looping over edges and computing
211c4762a1bSJed Brown      contribution for each vertex dof                         */
212c4762a1bSJed Brown   for (moab::Range::iterator iter = vlocal->begin(); iter != vlocal->end(); iter++) {
213c4762a1bSJed Brown     const moab::EntityHandle vhandle = *iter;
214c4762a1bSJed Brown 
2159566063dSJacob Faibussowitsch     PetscCall(DMMoabGetDofsBlocked(dm, 1, &vhandle, &dof));
216c4762a1bSJed Brown 
217c4762a1bSJed Brown     /* check if vertex is on the boundary */
2189566063dSJacob Faibussowitsch     PetscCall(DMMoabIsEntityOnBoundary(dm, vhandle, &vonboundary));
219c4762a1bSJed Brown 
220c4762a1bSJed Brown     if (vonboundary) {
221*9371c9d4SSatish Balay       const PetscScalar bcvals[2][2] = {
222*9371c9d4SSatish Balay         {hx, 0 },
223*9371c9d4SSatish Balay         {0,  hx}
224*9371c9d4SSatish Balay       };
2259566063dSJacob Faibussowitsch       PetscCall(MatSetValuesBlocked(Jpre, 1, &dof, 1, &dof, &bcvals[0][0], INSERT_VALUES));
226*9371c9d4SSatish Balay     } else {
227c4762a1bSJed Brown       const PetscInt    row = dof, col[] = {dof - 1, dof, dof + 1};
228c4762a1bSJed Brown       const PetscScalar dxxL = -user->alpha / hx, dxx0 = 2. * user->alpha / hx, dxxR = -user->alpha / hx;
229*9371c9d4SSatish Balay       const PetscScalar vals[2][3][2] = {
230*9371c9d4SSatish Balay         {{dxxL, 0}, {a * hx + dxx0, 0}, {dxxR, 0}},
231*9371c9d4SSatish Balay         {{0, dxxL}, {0, a * hx + dxx0}, {0, dxxR}}
232*9371c9d4SSatish Balay       };
2339566063dSJacob Faibussowitsch       PetscCall(MatSetValuesBlocked(Jpre, 1, &row, 3, col, &vals[0][0][0], INSERT_VALUES));
234c4762a1bSJed Brown     }
235c4762a1bSJed Brown   }
236c4762a1bSJed Brown 
2379566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY));
2389566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY));
239c4762a1bSJed Brown   if (J != Jpre) {
2409566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
2419566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
242c4762a1bSJed Brown   }
243c4762a1bSJed Brown   PetscFunctionReturn(0);
244c4762a1bSJed Brown }
245c4762a1bSJed Brown 
246*9371c9d4SSatish Balay static PetscErrorCode FormRHSFunction(TS ts, PetscReal t, Vec X, Vec F, void *ptr) {
247c4762a1bSJed Brown   UserCtx            user = (UserCtx)ptr;
248c4762a1bSJed Brown   DM                 dm;
249c4762a1bSJed Brown   PetscReal          hx;
250c4762a1bSJed Brown   const Field       *x;
251c4762a1bSJed Brown   Field             *f;
252c4762a1bSJed Brown   PetscInt           dof;
253c4762a1bSJed Brown   const moab::Range *ownedvtx;
254c4762a1bSJed Brown 
255c4762a1bSJed Brown   PetscFunctionBegin;
256c4762a1bSJed Brown   hx = 1.0 / user->n;
2579566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   /* Get pointers to vector data */
2609566063dSJacob Faibussowitsch   PetscCall(VecSet(F, 0.0));
261c4762a1bSJed Brown 
2629566063dSJacob Faibussowitsch   PetscCall(DMMoabVecGetArrayRead(dm, X, &x));
2639566063dSJacob Faibussowitsch   PetscCall(DMMoabVecGetArray(dm, F, &f));
264c4762a1bSJed Brown 
2659566063dSJacob Faibussowitsch   PetscCall(DMMoabGetLocalVertices(dm, &ownedvtx, NULL));
266c4762a1bSJed Brown 
267c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
268c4762a1bSJed Brown   for (moab::Range::iterator iter = ownedvtx->begin(); iter != ownedvtx->end(); iter++) {
269c4762a1bSJed Brown     const moab::EntityHandle vhandle = *iter;
2709566063dSJacob Faibussowitsch     PetscCall(DMMoabGetDofsBlockedLocal(dm, 1, &vhandle, &dof));
271c4762a1bSJed Brown 
272c4762a1bSJed Brown     PetscScalar u = x[dof].u, v = x[dof].v;
273c4762a1bSJed Brown     f[dof].u = hx * (user->A + u * u * v - (user->B + 1) * u);
274c4762a1bSJed Brown     f[dof].v = hx * (user->B * u - u * u * v);
275c4762a1bSJed Brown   }
276c4762a1bSJed Brown 
277c4762a1bSJed Brown   /* Restore vectors */
2789566063dSJacob Faibussowitsch   PetscCall(DMMoabVecRestoreArrayRead(dm, X, &x));
2799566063dSJacob Faibussowitsch   PetscCall(DMMoabVecRestoreArray(dm, F, &f));
280c4762a1bSJed Brown   PetscFunctionReturn(0);
281c4762a1bSJed Brown }
282c4762a1bSJed Brown 
283*9371c9d4SSatish Balay static PetscErrorCode FormIFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) {
284c4762a1bSJed Brown   UserCtx            user = (UserCtx)ctx;
285c4762a1bSJed Brown   DM                 dm;
286c4762a1bSJed Brown   Field             *x, *xdot, *f;
287c4762a1bSJed Brown   PetscReal          hx;
288c4762a1bSJed Brown   Vec                Xloc;
289c4762a1bSJed Brown   PetscInt           i, bcindx;
290c4762a1bSJed Brown   PetscBool          elem_on_boundary;
291c4762a1bSJed Brown   const moab::Range *vlocal;
292c4762a1bSJed Brown 
293c4762a1bSJed Brown   PetscFunctionBegin;
294c4762a1bSJed Brown   hx = 1.0 / user->n;
2959566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
296c4762a1bSJed Brown 
297c4762a1bSJed Brown   /* get the essential MOAB mesh related quantities needed for FEM assembly */
2989566063dSJacob Faibussowitsch   PetscCall(DMMoabGetLocalVertices(dm, &vlocal, NULL));
299c4762a1bSJed Brown 
300c4762a1bSJed Brown   /* reset the residual vector */
3019566063dSJacob Faibussowitsch   PetscCall(VecSet(F, 0.0));
302c4762a1bSJed Brown 
3039566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xloc));
3049566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
3059566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
306c4762a1bSJed Brown 
307c4762a1bSJed Brown   /* get the local representation of the arrays from Vectors */
3089566063dSJacob Faibussowitsch   PetscCall(DMMoabVecGetArrayRead(dm, Xloc, &x));
3099566063dSJacob Faibussowitsch   PetscCall(DMMoabVecGetArrayRead(dm, Xdot, &xdot));
3109566063dSJacob Faibussowitsch   PetscCall(DMMoabVecGetArray(dm, F, &f));
311c4762a1bSJed Brown 
312c4762a1bSJed Brown   /* loop over local elements */
313c4762a1bSJed Brown   for (moab::Range::iterator iter = vlocal->begin(); iter != vlocal->end(); iter++) {
314c4762a1bSJed Brown     const moab::EntityHandle vhandle = *iter;
315c4762a1bSJed Brown 
3169566063dSJacob Faibussowitsch     PetscCall(DMMoabGetDofsBlockedLocal(dm, 1, &vhandle, &i));
317c4762a1bSJed Brown 
318c4762a1bSJed Brown     /* check if vertex is on the boundary */
3199566063dSJacob Faibussowitsch     PetscCall(DMMoabIsEntityOnBoundary(dm, vhandle, &elem_on_boundary));
320c4762a1bSJed Brown 
321c4762a1bSJed Brown     if (elem_on_boundary) {
3229566063dSJacob Faibussowitsch       PetscCall(DMMoabGetDofsBlocked(dm, 1, &vhandle, &bcindx));
323c4762a1bSJed Brown       if (bcindx == 0) { /* Apply left BC */
324c4762a1bSJed Brown         f[i].u = hx * (x[i].u - user->leftbc.u);
325c4762a1bSJed Brown         f[i].v = hx * (x[i].v - user->leftbc.v);
326c4762a1bSJed Brown       } else { /* Apply right BC */
327c4762a1bSJed Brown         f[i].u = hx * (x[i].u - user->rightbc.u);
328c4762a1bSJed Brown         f[i].v = hx * (x[i].v - user->rightbc.v);
329c4762a1bSJed Brown       }
330*9371c9d4SSatish Balay     } else {
331c4762a1bSJed Brown       f[i].u = hx * xdot[i].u - user->alpha * (x[i - 1].u - 2. * x[i].u + x[i + 1].u) / hx;
332c4762a1bSJed Brown       f[i].v = hx * xdot[i].v - user->alpha * (x[i - 1].v - 2. * x[i].v + x[i + 1].v) / hx;
333c4762a1bSJed Brown     }
334c4762a1bSJed Brown   }
335c4762a1bSJed Brown 
336c4762a1bSJed Brown   /* Restore data */
3379566063dSJacob Faibussowitsch   PetscCall(DMMoabVecRestoreArrayRead(dm, Xloc, &x));
3389566063dSJacob Faibussowitsch   PetscCall(DMMoabVecRestoreArrayRead(dm, Xdot, &xdot));
3399566063dSJacob Faibussowitsch   PetscCall(DMMoabVecRestoreArray(dm, F, &f));
3409566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xloc));
341c4762a1bSJed Brown   PetscFunctionReturn(0);
342c4762a1bSJed Brown }
343c4762a1bSJed Brown 
344*9371c9d4SSatish Balay PetscErrorCode FormInitialSolution(TS ts, Vec X, void *ctx) {
345c4762a1bSJed Brown   UserCtx               user = (UserCtx)ctx;
346c4762a1bSJed Brown   PetscReal             vpos[3];
347c4762a1bSJed Brown   DM                    dm;
348c4762a1bSJed Brown   Field                *x;
349c4762a1bSJed Brown   const moab::Range    *vowned;
350c4762a1bSJed Brown   PetscInt              dof;
351c4762a1bSJed Brown   moab::Range::iterator iter;
352c4762a1bSJed Brown 
353c4762a1bSJed Brown   PetscFunctionBegin;
3549566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
355c4762a1bSJed Brown 
356c4762a1bSJed Brown   /* get the essential MOAB mesh related quantities needed for FEM assembly */
3579566063dSJacob Faibussowitsch   PetscCall(DMMoabGetLocalVertices(dm, &vowned, NULL));
358c4762a1bSJed Brown 
3599566063dSJacob Faibussowitsch   PetscCall(VecSet(X, 0.0));
360c4762a1bSJed Brown 
361c4762a1bSJed Brown   /* Get pointers to vector data */
3629566063dSJacob Faibussowitsch   PetscCall(DMMoabVecGetArray(dm, X, &x));
363c4762a1bSJed Brown 
364c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
365c4762a1bSJed Brown   for (moab::Range::iterator iter = vowned->begin(); iter != vowned->end(); iter++) {
366c4762a1bSJed Brown     const moab::EntityHandle vhandle = *iter;
3679566063dSJacob Faibussowitsch     PetscCall(DMMoabGetDofsBlockedLocal(dm, 1, &vhandle, &dof));
368c4762a1bSJed Brown 
369c4762a1bSJed Brown     /* compute the mid-point of the element and use a 1-point lumped quadrature */
3709566063dSJacob Faibussowitsch     PetscCall(DMMoabGetVertexCoordinates(dm, 1, &vhandle, vpos));
371c4762a1bSJed Brown 
372c4762a1bSJed Brown     PetscReal xi = vpos[0];
373c4762a1bSJed Brown     x[dof].u     = user->leftbc.u * (1. - xi) + user->rightbc.u * xi + PetscSinReal(2. * PETSC_PI * xi);
374c4762a1bSJed Brown     x[dof].v     = user->leftbc.v * (1. - xi) + user->rightbc.v * xi;
375c4762a1bSJed Brown   }
376c4762a1bSJed Brown 
377c4762a1bSJed Brown   /* Restore vectors */
3789566063dSJacob Faibussowitsch   PetscCall(DMMoabVecRestoreArray(dm, X, &x));
379c4762a1bSJed Brown   PetscFunctionReturn(0);
380c4762a1bSJed Brown }
381c4762a1bSJed Brown 
382c4762a1bSJed Brown /*TEST
383c4762a1bSJed Brown 
384c4762a1bSJed Brown     build:
385c4762a1bSJed Brown       requires: moab
386c4762a1bSJed Brown 
387c4762a1bSJed Brown     test:
388c4762a1bSJed Brown       args: -n 20 -ts_type rosw -ts_rosw_type 2p -ts_dt 5e-2 -ts_adapt_type none
389c4762a1bSJed Brown 
390c4762a1bSJed Brown     test:
391c4762a1bSJed Brown       suffix: 2
392c4762a1bSJed Brown       nsize: 2
393c4762a1bSJed Brown       args: -n 50 -ts_type glee -ts_adapt_type none -ts_dt 0.1 -io
394c4762a1bSJed Brown       TODO:
395c4762a1bSJed Brown 
396c4762a1bSJed Brown TEST*/
397