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
Initialize_AppContext(UserCtx * puser)37d71ae5a4SJacob Faibussowitsch PetscErrorCode Initialize_AppContext(UserCtx *puser)
38d71ae5a4SJacob Faibussowitsch {
39c4762a1bSJed Brown UserCtx user;
40c4762a1bSJed Brown
41c4762a1bSJed Brown PetscFunctionBegin;
429566063dSJacob Faibussowitsch PetscCall(PetscNew(&user));
43d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Advection-reaction options", "ex35.cxx");
44c4762a1bSJed Brown {
45c4762a1bSJed Brown user->nvars = 2;
46c4762a1bSJed Brown user->A = 1;
47c4762a1bSJed Brown user->B = 3;
48c4762a1bSJed Brown user->alpha = 0.02;
49c4762a1bSJed Brown user->leftbc.u = 1;
50c4762a1bSJed Brown user->rightbc.u = 1;
51c4762a1bSJed Brown user->leftbc.v = 3;
52c4762a1bSJed Brown user->rightbc.v = 3;
53c4762a1bSJed Brown user->n = 10;
54c4762a1bSJed Brown user->ntsteps = 10000;
55c4762a1bSJed Brown user->io = PETSC_FALSE;
569566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-A", "Reaction rate", "ex35.cxx", user->A, &user->A, NULL));
579566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-B", "Reaction rate", "ex35.cxx", user->B, &user->B, NULL));
589566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-alpha", "Diffusion coefficient", "ex35.cxx", user->alpha, &user->alpha, NULL));
599566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-uleft", "Dirichlet boundary condition", "ex35.cxx", user->leftbc.u, &user->leftbc.u, NULL));
609566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-uright", "Dirichlet boundary condition", "ex35.cxx", user->rightbc.u, &user->rightbc.u, NULL));
619566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-vleft", "Dirichlet boundary condition", "ex35.cxx", user->leftbc.v, &user->leftbc.v, NULL));
629566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-vright", "Dirichlet boundary condition", "ex35.cxx", user->rightbc.v, &user->rightbc.v, NULL));
639566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-n", "Number of 1-D elements", "ex35.cxx", user->n, &user->n, NULL));
649566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ndt", "Number of time steps", "ex35.cxx", user->ntsteps, &user->ntsteps, NULL));
659566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-io", "Write the mesh and solution output to a file.", "ex35.cxx", user->io, &user->io, NULL));
66c4762a1bSJed Brown user->npts = user->n + 1;
67c4762a1bSJed Brown }
68d0609cedSBarry Smith PetscOptionsEnd();
69c4762a1bSJed Brown
70c4762a1bSJed Brown *puser = user;
713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
72c4762a1bSJed Brown }
73c4762a1bSJed Brown
Destroy_AppContext(UserCtx * user)74d71ae5a4SJacob Faibussowitsch PetscErrorCode Destroy_AppContext(UserCtx *user)
75d71ae5a4SJacob Faibussowitsch {
76c4762a1bSJed Brown PetscFunctionBegin;
779566063dSJacob Faibussowitsch PetscCall(PetscFree(*user));
783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
79c4762a1bSJed Brown }
80c4762a1bSJed Brown
81c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS, Vec, void *);
82c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS, PetscReal, Vec, Vec, void *);
83c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS, PetscReal, Vec, Vec, Vec, void *);
84c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
85c4762a1bSJed Brown
86c4762a1bSJed Brown /****************
87c4762a1bSJed Brown * *
88c4762a1bSJed Brown * MAIN *
89c4762a1bSJed Brown * *
90c4762a1bSJed Brown ****************/
main(int argc,char ** argv)91d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
92d71ae5a4SJacob Faibussowitsch {
93c4762a1bSJed Brown TS ts; /* nonlinear solver */
94c4762a1bSJed Brown Vec X; /* solution, residual vectors */
95c4762a1bSJed Brown Mat J; /* Jacobian matrix */
96c4762a1bSJed Brown PetscInt steps, mx;
97c4762a1bSJed Brown PetscReal hx, dt, ftime;
98c4762a1bSJed Brown UserCtx user; /* user-defined work context */
99c4762a1bSJed Brown TSConvergedReason reason;
100c4762a1bSJed Brown DM dm;
101c4762a1bSJed Brown const char *fields[2] = {"U", "V"};
102c4762a1bSJed Brown
103327415f7SBarry Smith PetscFunctionBeginUser;
104c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, nullptr, help));
105c4762a1bSJed Brown
106c4762a1bSJed Brown /* Initialize the user context struct */
1079566063dSJacob Faibussowitsch PetscCall(Initialize_AppContext(&user));
108c4762a1bSJed Brown
109c4762a1bSJed Brown /* Fill in the user defined work context: */
1109566063dSJacob Faibussowitsch PetscCall(DMMoabCreateBoxMesh(PETSC_COMM_WORLD, 1, PETSC_FALSE, NULL, user->n, 1, &dm));
1119566063dSJacob Faibussowitsch PetscCall(DMMoabSetFieldNames(dm, user->nvars, fields));
1129566063dSJacob Faibussowitsch PetscCall(DMMoabSetBlockSize(dm, user->nvars));
1139566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm));
114c4762a1bSJed Brown
115c4762a1bSJed Brown /* SetUp the data structures for DMMOAB */
1169566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm));
117c4762a1bSJed Brown
118c4762a1bSJed Brown /* Create timestepping solver context */
1199566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1209566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, dm));
1219566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSARKIMEX));
1229566063dSJacob Faibussowitsch PetscCall(TSSetEquationType(ts, TS_EQ_DAE_IMPLICIT_INDEX1));
1239566063dSJacob Faibussowitsch PetscCall(DMSetMatType(dm, MATBAIJ));
1249566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &J));
125c4762a1bSJed Brown
1269566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, FormRHSFunction, user));
1279566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, FormIFunction, user));
1289566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, J, J, FormIJacobian, user));
129c4762a1bSJed Brown
130c4762a1bSJed Brown ftime = 10.0;
1319566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, user->ntsteps));
1329566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, ftime));
1339566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
134c4762a1bSJed Brown
135c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136c4762a1bSJed Brown Create the solution vector and set the initial conditions
137c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1389566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &X));
139c4762a1bSJed Brown
1409566063dSJacob Faibussowitsch PetscCall(FormInitialSolution(ts, X, user));
1419566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, X));
1429566063dSJacob Faibussowitsch PetscCall(VecGetSize(X, &mx));
143c4762a1bSJed Brown hx = 1.0 / (PetscReal)(mx / 2 - 1);
144c4762a1bSJed Brown dt = 0.4 * PetscSqr(hx) / user->alpha; /* Diffusive stability limit */
1459566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt));
146c4762a1bSJed Brown
147c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148c4762a1bSJed Brown Set runtime options
149c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1509566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
151c4762a1bSJed Brown
152c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153c4762a1bSJed Brown Solve nonlinear system
154c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1559566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, X));
1569566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime));
1579566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps));
1589566063dSJacob Faibussowitsch PetscCall(TSGetConvergedReason(ts, &reason));
15963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %" PetscInt_FMT " steps\n", TSConvergedReasons[reason], (double)ftime, steps));
160c4762a1bSJed Brown
161c4762a1bSJed Brown if (user->io) {
162c4762a1bSJed Brown /* Print the numerical solution to screen and then dump to file */
1639566063dSJacob Faibussowitsch PetscCall(VecView(X, PETSC_VIEWER_STDOUT_WORLD));
164c4762a1bSJed Brown
165c4762a1bSJed Brown /* Write out the solution along with the mesh */
1669566063dSJacob Faibussowitsch PetscCall(DMMoabSetGlobalFieldVector(dm, X));
167c4762a1bSJed Brown #ifdef MOAB_HAVE_HDF5
1689566063dSJacob Faibussowitsch PetscCall(DMMoabOutput(dm, "ex35.h5m", ""));
169c4762a1bSJed Brown #else
170c4762a1bSJed Brown /* MOAB does not support true parallel writers that aren't HDF5 based
171c4762a1bSJed Brown And so if you are using VTK as the output format in parallel,
172c4762a1bSJed Brown the data could be jumbled due to the order in which the processors
173c4762a1bSJed Brown write out their parts of the mesh and solution tags
174c4762a1bSJed Brown */
1759566063dSJacob Faibussowitsch PetscCall(DMMoabOutput(dm, "ex35.vtk", ""));
176c4762a1bSJed Brown #endif
177c4762a1bSJed Brown }
178c4762a1bSJed Brown
179c4762a1bSJed Brown /* Free work space.
180c4762a1bSJed Brown Free all PETSc related resources: */
1819566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J));
1829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X));
1839566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
1849566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
185c4762a1bSJed Brown
186c4762a1bSJed Brown /* Free all MOAB related resources: */
1879566063dSJacob Faibussowitsch PetscCall(Destroy_AppContext(&user));
1889566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
189b122ec5aSJacob Faibussowitsch return 0;
190c4762a1bSJed Brown }
191c4762a1bSJed Brown
192c4762a1bSJed Brown /*
193c4762a1bSJed Brown IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
194c4762a1bSJed Brown */
FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void * ptr)195d71ae5a4SJacob Faibussowitsch PetscErrorCode FormIJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat J, Mat Jpre, void *ptr)
196d71ae5a4SJacob Faibussowitsch {
197c4762a1bSJed Brown UserCtx user = (UserCtx)ptr;
198c4762a1bSJed Brown PetscInt dof;
199c4762a1bSJed Brown PetscReal hx;
200c4762a1bSJed Brown DM dm;
201c4762a1bSJed Brown const moab::Range *vlocal;
202c4762a1bSJed Brown PetscBool vonboundary;
203c4762a1bSJed Brown
204c4762a1bSJed Brown PetscFunctionBegin;
2059566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm));
206c4762a1bSJed Brown
207c4762a1bSJed Brown /* get the essential MOAB mesh related quantities needed for FEM assembly */
2089566063dSJacob Faibussowitsch PetscCall(DMMoabGetLocalVertices(dm, &vlocal, NULL));
209c4762a1bSJed Brown
210c4762a1bSJed Brown /* compute local element sizes - structured grid */
211c4762a1bSJed Brown hx = 1.0 / user->n;
212c4762a1bSJed Brown
213c4762a1bSJed Brown /* Compute function over the locally owned part of the grid
214c4762a1bSJed Brown Assemble the operator by looping over edges and computing
215c4762a1bSJed Brown contribution for each vertex dof */
216c4762a1bSJed Brown for (moab::Range::iterator iter = vlocal->begin(); iter != vlocal->end(); iter++) {
217c4762a1bSJed Brown const moab::EntityHandle vhandle = *iter;
218c4762a1bSJed Brown
2199566063dSJacob Faibussowitsch PetscCall(DMMoabGetDofsBlocked(dm, 1, &vhandle, &dof));
220c4762a1bSJed Brown
221c4762a1bSJed Brown /* check if vertex is on the boundary */
2229566063dSJacob Faibussowitsch PetscCall(DMMoabIsEntityOnBoundary(dm, vhandle, &vonboundary));
223c4762a1bSJed Brown
224c4762a1bSJed Brown if (vonboundary) {
2259371c9d4SSatish Balay const PetscScalar bcvals[2][2] = {
2269371c9d4SSatish Balay {hx, 0 },
2279371c9d4SSatish Balay {0, hx}
2289371c9d4SSatish Balay };
2299566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(Jpre, 1, &dof, 1, &dof, &bcvals[0][0], INSERT_VALUES));
2309371c9d4SSatish Balay } else {
231c4762a1bSJed Brown const PetscInt row = dof, col[] = {dof - 1, dof, dof + 1};
232c4762a1bSJed Brown const PetscScalar dxxL = -user->alpha / hx, dxx0 = 2. * user->alpha / hx, dxxR = -user->alpha / hx;
2339371c9d4SSatish Balay const PetscScalar vals[2][3][2] = {
2349371c9d4SSatish Balay {{dxxL, 0}, {a * hx + dxx0, 0}, {dxxR, 0}},
2359371c9d4SSatish Balay {{0, dxxL}, {0, a * hx + dxx0}, {0, dxxR}}
2369371c9d4SSatish Balay };
2379566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(Jpre, 1, &row, 3, col, &vals[0][0][0], INSERT_VALUES));
238c4762a1bSJed Brown }
239c4762a1bSJed Brown }
240c4762a1bSJed Brown
2419566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY));
2429566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY));
243c4762a1bSJed Brown if (J != Jpre) {
2449566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
2459566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
246c4762a1bSJed Brown }
2473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
248c4762a1bSJed Brown }
249c4762a1bSJed Brown
FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void * ptr)250d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormRHSFunction(TS ts, PetscReal t, Vec X, Vec F, void *ptr)
251d71ae5a4SJacob Faibussowitsch {
252c4762a1bSJed Brown UserCtx user = (UserCtx)ptr;
253c4762a1bSJed Brown DM dm;
254c4762a1bSJed Brown PetscReal hx;
255c4762a1bSJed Brown const Field *x;
256c4762a1bSJed Brown Field *f;
257c4762a1bSJed Brown PetscInt dof;
258c4762a1bSJed Brown const moab::Range *ownedvtx;
259c4762a1bSJed Brown
260c4762a1bSJed Brown PetscFunctionBegin;
261c4762a1bSJed Brown hx = 1.0 / user->n;
2629566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm));
263c4762a1bSJed Brown
264c4762a1bSJed Brown /* Get pointers to vector data */
2659566063dSJacob Faibussowitsch PetscCall(VecSet(F, 0.0));
266c4762a1bSJed Brown
2679566063dSJacob Faibussowitsch PetscCall(DMMoabVecGetArrayRead(dm, X, &x));
2689566063dSJacob Faibussowitsch PetscCall(DMMoabVecGetArray(dm, F, &f));
269c4762a1bSJed Brown
2709566063dSJacob Faibussowitsch PetscCall(DMMoabGetLocalVertices(dm, &ownedvtx, NULL));
271c4762a1bSJed Brown
272c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */
273c4762a1bSJed Brown for (moab::Range::iterator iter = ownedvtx->begin(); iter != ownedvtx->end(); iter++) {
274c4762a1bSJed Brown const moab::EntityHandle vhandle = *iter;
2759566063dSJacob Faibussowitsch PetscCall(DMMoabGetDofsBlockedLocal(dm, 1, &vhandle, &dof));
276c4762a1bSJed Brown
277c4762a1bSJed Brown PetscScalar u = x[dof].u, v = x[dof].v;
278c4762a1bSJed Brown f[dof].u = hx * (user->A + u * u * v - (user->B + 1) * u);
279c4762a1bSJed Brown f[dof].v = hx * (user->B * u - u * u * v);
280c4762a1bSJed Brown }
281c4762a1bSJed Brown
282c4762a1bSJed Brown /* Restore vectors */
2839566063dSJacob Faibussowitsch PetscCall(DMMoabVecRestoreArrayRead(dm, X, &x));
2849566063dSJacob Faibussowitsch PetscCall(DMMoabVecRestoreArray(dm, F, &f));
2853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
286c4762a1bSJed Brown }
287c4762a1bSJed Brown
FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,PetscCtx ctx)288*2a8381b2SBarry Smith static PetscErrorCode FormIFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, PetscCtx ctx)
289d71ae5a4SJacob Faibussowitsch {
290c4762a1bSJed Brown UserCtx user = (UserCtx)ctx;
291c4762a1bSJed Brown DM dm;
292c4762a1bSJed Brown Field *x, *xdot, *f;
293c4762a1bSJed Brown PetscReal hx;
294c4762a1bSJed Brown Vec Xloc;
295c4762a1bSJed Brown PetscInt i, bcindx;
296c4762a1bSJed Brown PetscBool elem_on_boundary;
297c4762a1bSJed Brown const moab::Range *vlocal;
298c4762a1bSJed Brown
299c4762a1bSJed Brown PetscFunctionBegin;
300c4762a1bSJed Brown hx = 1.0 / user->n;
3019566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm));
302c4762a1bSJed Brown
303c4762a1bSJed Brown /* get the essential MOAB mesh related quantities needed for FEM assembly */
3049566063dSJacob Faibussowitsch PetscCall(DMMoabGetLocalVertices(dm, &vlocal, NULL));
305c4762a1bSJed Brown
306c4762a1bSJed Brown /* reset the residual vector */
3079566063dSJacob Faibussowitsch PetscCall(VecSet(F, 0.0));
308c4762a1bSJed Brown
3099566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &Xloc));
3109566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
3119566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
312c4762a1bSJed Brown
313c4762a1bSJed Brown /* get the local representation of the arrays from Vectors */
3149566063dSJacob Faibussowitsch PetscCall(DMMoabVecGetArrayRead(dm, Xloc, &x));
3159566063dSJacob Faibussowitsch PetscCall(DMMoabVecGetArrayRead(dm, Xdot, &xdot));
3169566063dSJacob Faibussowitsch PetscCall(DMMoabVecGetArray(dm, F, &f));
317c4762a1bSJed Brown
318c4762a1bSJed Brown /* loop over local elements */
319c4762a1bSJed Brown for (moab::Range::iterator iter = vlocal->begin(); iter != vlocal->end(); iter++) {
320c4762a1bSJed Brown const moab::EntityHandle vhandle = *iter;
321c4762a1bSJed Brown
3229566063dSJacob Faibussowitsch PetscCall(DMMoabGetDofsBlockedLocal(dm, 1, &vhandle, &i));
323c4762a1bSJed Brown
324c4762a1bSJed Brown /* check if vertex is on the boundary */
3259566063dSJacob Faibussowitsch PetscCall(DMMoabIsEntityOnBoundary(dm, vhandle, &elem_on_boundary));
326c4762a1bSJed Brown
327c4762a1bSJed Brown if (elem_on_boundary) {
3289566063dSJacob Faibussowitsch PetscCall(DMMoabGetDofsBlocked(dm, 1, &vhandle, &bcindx));
329c4762a1bSJed Brown if (bcindx == 0) { /* Apply left BC */
330c4762a1bSJed Brown f[i].u = hx * (x[i].u - user->leftbc.u);
331c4762a1bSJed Brown f[i].v = hx * (x[i].v - user->leftbc.v);
332c4762a1bSJed Brown } else { /* Apply right BC */
333c4762a1bSJed Brown f[i].u = hx * (x[i].u - user->rightbc.u);
334c4762a1bSJed Brown f[i].v = hx * (x[i].v - user->rightbc.v);
335c4762a1bSJed Brown }
3369371c9d4SSatish Balay } else {
337c4762a1bSJed Brown f[i].u = hx * xdot[i].u - user->alpha * (x[i - 1].u - 2. * x[i].u + x[i + 1].u) / hx;
338c4762a1bSJed Brown f[i].v = hx * xdot[i].v - user->alpha * (x[i - 1].v - 2. * x[i].v + x[i + 1].v) / hx;
339c4762a1bSJed Brown }
340c4762a1bSJed Brown }
341c4762a1bSJed Brown
342c4762a1bSJed Brown /* Restore data */
3439566063dSJacob Faibussowitsch PetscCall(DMMoabVecRestoreArrayRead(dm, Xloc, &x));
3449566063dSJacob Faibussowitsch PetscCall(DMMoabVecRestoreArrayRead(dm, Xdot, &xdot));
3459566063dSJacob Faibussowitsch PetscCall(DMMoabVecRestoreArray(dm, F, &f));
3469566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &Xloc));
3473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
348c4762a1bSJed Brown }
349c4762a1bSJed Brown
FormInitialSolution(TS ts,Vec X,PetscCtx ctx)350*2a8381b2SBarry Smith PetscErrorCode FormInitialSolution(TS ts, Vec X, PetscCtx ctx)
351d71ae5a4SJacob Faibussowitsch {
352c4762a1bSJed Brown UserCtx user = (UserCtx)ctx;
353c4762a1bSJed Brown PetscReal vpos[3];
354c4762a1bSJed Brown DM dm;
355c4762a1bSJed Brown Field *x;
356c4762a1bSJed Brown const moab::Range *vowned;
357c4762a1bSJed Brown PetscInt dof;
358c4762a1bSJed Brown moab::Range::iterator iter;
359c4762a1bSJed Brown
360c4762a1bSJed Brown PetscFunctionBegin;
3619566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm));
362c4762a1bSJed Brown
363c4762a1bSJed Brown /* get the essential MOAB mesh related quantities needed for FEM assembly */
3649566063dSJacob Faibussowitsch PetscCall(DMMoabGetLocalVertices(dm, &vowned, NULL));
365c4762a1bSJed Brown
3669566063dSJacob Faibussowitsch PetscCall(VecSet(X, 0.0));
367c4762a1bSJed Brown
368c4762a1bSJed Brown /* Get pointers to vector data */
3699566063dSJacob Faibussowitsch PetscCall(DMMoabVecGetArray(dm, X, &x));
370c4762a1bSJed Brown
371c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */
372c4762a1bSJed Brown for (moab::Range::iterator iter = vowned->begin(); iter != vowned->end(); iter++) {
373c4762a1bSJed Brown const moab::EntityHandle vhandle = *iter;
3749566063dSJacob Faibussowitsch PetscCall(DMMoabGetDofsBlockedLocal(dm, 1, &vhandle, &dof));
375c4762a1bSJed Brown
376c4762a1bSJed Brown /* compute the mid-point of the element and use a 1-point lumped quadrature */
3779566063dSJacob Faibussowitsch PetscCall(DMMoabGetVertexCoordinates(dm, 1, &vhandle, vpos));
378c4762a1bSJed Brown
379c4762a1bSJed Brown PetscReal xi = vpos[0];
380c4762a1bSJed Brown x[dof].u = user->leftbc.u * (1. - xi) + user->rightbc.u * xi + PetscSinReal(2. * PETSC_PI * xi);
381c4762a1bSJed Brown x[dof].v = user->leftbc.v * (1. - xi) + user->rightbc.v * xi;
382c4762a1bSJed Brown }
383c4762a1bSJed Brown
384c4762a1bSJed Brown /* Restore vectors */
3859566063dSJacob Faibussowitsch PetscCall(DMMoabVecRestoreArray(dm, X, &x));
3863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
387c4762a1bSJed Brown }
388c4762a1bSJed Brown
389c4762a1bSJed Brown /*TEST
390c4762a1bSJed Brown
391c4762a1bSJed Brown build:
392ca4445c7SIlya Fursov requires: moab !complex
393c4762a1bSJed Brown
394c4762a1bSJed Brown test:
395188af4bfSBarry Smith args: -n 20 -ts_type rosw -ts_rosw_type 2p -ts_time_step 5e-2 -ts_adapt_type none
396c4762a1bSJed Brown
397c4762a1bSJed Brown test:
398c4762a1bSJed Brown suffix: 2
399c4762a1bSJed Brown nsize: 2
400188af4bfSBarry Smith args: -n 50 -ts_type glee -ts_adapt_type none -ts_time_step 0.1 -io
401c4762a1bSJed Brown TODO:
402c4762a1bSJed Brown
403c4762a1bSJed Brown TEST*/
404