xref: /petsc/src/ts/tutorials/ex35.cxx (revision 1e1ea65d8de51fde77ce8a787efbef25e407badc)
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 
37c4762a1bSJed Brown PetscErrorCode Initialize_AppContext(UserCtx *puser)
38c4762a1bSJed Brown {
39c4762a1bSJed Brown   UserCtx           user;
40c4762a1bSJed Brown   PetscErrorCode    ierr;
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   PetscFunctionBegin;
43c4762a1bSJed Brown   ierr = PetscNew(&user);CHKERRQ(ierr);
44c4762a1bSJed Brown 
45*1e1ea65dSPierre Jolivet   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Advection-reaction options","ex35.cxx");CHKERRQ(ierr);
46c4762a1bSJed Brown   {
47c4762a1bSJed Brown     user->nvars  = 2;
48c4762a1bSJed Brown     user->A      = 1;
49c4762a1bSJed Brown     user->B      = 3;
50c4762a1bSJed Brown     user->alpha  = 0.02;
51c4762a1bSJed Brown     user->leftbc.u  = 1;
52c4762a1bSJed Brown     user->rightbc.u = 1;
53c4762a1bSJed Brown     user->leftbc.v  = 3;
54c4762a1bSJed Brown     user->rightbc.v = 3;
55c4762a1bSJed Brown     user->n      = 10;
56c4762a1bSJed Brown     user->ntsteps = 10000;
57c4762a1bSJed Brown     user->io = PETSC_FALSE;
58c4762a1bSJed Brown     ierr = PetscOptionsReal("-A","Reaction rate","ex35.cxx",user->A,&user->A,NULL);CHKERRQ(ierr);
59c4762a1bSJed Brown     ierr = PetscOptionsReal("-B","Reaction rate","ex35.cxx",user->B,&user->B,NULL);CHKERRQ(ierr);
60c4762a1bSJed Brown     ierr = PetscOptionsReal("-alpha","Diffusion coefficient","ex35.cxx",user->alpha,&user->alpha,NULL);CHKERRQ(ierr);
61c4762a1bSJed Brown     ierr = PetscOptionsScalar("-uleft","Dirichlet boundary condition","ex35.cxx",user->leftbc.u,&user->leftbc.u,NULL);CHKERRQ(ierr);
62c4762a1bSJed Brown     ierr = PetscOptionsScalar("-uright","Dirichlet boundary condition","ex35.cxx",user->rightbc.u,&user->rightbc.u,NULL);CHKERRQ(ierr);
63c4762a1bSJed Brown     ierr = PetscOptionsScalar("-vleft","Dirichlet boundary condition","ex35.cxx",user->leftbc.v,&user->leftbc.v,NULL);CHKERRQ(ierr);
64c4762a1bSJed Brown     ierr = PetscOptionsScalar("-vright","Dirichlet boundary condition","ex35.cxx",user->rightbc.v,&user->rightbc.v,NULL);CHKERRQ(ierr);
65c4762a1bSJed Brown     ierr = PetscOptionsInt("-n","Number of 1-D elements","ex35.cxx",user->n,&user->n,NULL);CHKERRQ(ierr);
66c4762a1bSJed Brown     ierr = PetscOptionsInt("-ndt","Number of time steps","ex35.cxx",user->ntsteps,&user->ntsteps,NULL);CHKERRQ(ierr);
67c4762a1bSJed Brown     ierr = PetscOptionsBool("-io","Write the mesh and solution output to a file.","ex35.cxx",user->io,&user->io,NULL);CHKERRQ(ierr);
68c4762a1bSJed Brown     user->npts   = user->n+1;
69c4762a1bSJed Brown   }
70c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   *puser = user;
73c4762a1bSJed Brown   PetscFunctionReturn(0);
74c4762a1bSJed Brown }
75c4762a1bSJed Brown 
76c4762a1bSJed Brown PetscErrorCode Destroy_AppContext(UserCtx *user)
77c4762a1bSJed Brown {
78c4762a1bSJed Brown   PetscErrorCode ierr;
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   PetscFunctionBegin;
81c4762a1bSJed Brown   ierr = PetscFree(*user);CHKERRQ(ierr);
82c4762a1bSJed Brown   PetscFunctionReturn(0);
83c4762a1bSJed Brown }
84c4762a1bSJed Brown 
85c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS,Vec,void*);
86c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*);
87c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
88c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
89c4762a1bSJed Brown 
90c4762a1bSJed Brown /****************
91c4762a1bSJed Brown  *              *
92c4762a1bSJed Brown  *     MAIN     *
93c4762a1bSJed Brown  *              *
94c4762a1bSJed Brown  ****************/
95c4762a1bSJed Brown int main(int argc,char **argv)
96c4762a1bSJed Brown {
97c4762a1bSJed Brown   TS                ts;         /* nonlinear solver */
98c4762a1bSJed Brown   Vec               X;          /* solution, residual vectors */
99c4762a1bSJed Brown   Mat               J;          /* Jacobian matrix */
100c4762a1bSJed Brown   PetscInt          steps,mx;
101c4762a1bSJed Brown   PetscErrorCode    ierr;
102c4762a1bSJed Brown   PetscReal         hx,dt,ftime;
103c4762a1bSJed Brown   UserCtx           user;       /* user-defined work context */
104c4762a1bSJed Brown   TSConvergedReason reason;
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   DM                dm;
107c4762a1bSJed Brown   const char        *fields[2] = {"U","V"};
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char *)0,help);if (ierr) return ierr;
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   /* Initialize the user context struct */
112c4762a1bSJed Brown   ierr = Initialize_AppContext(&user);CHKERRQ(ierr);
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   /* Fill in the user defined work context: */
115c4762a1bSJed Brown   ierr = DMMoabCreateBoxMesh(PETSC_COMM_WORLD, 1, PETSC_FALSE, NULL, user->n, 1, &dm);CHKERRQ(ierr);
116c4762a1bSJed Brown   ierr = DMMoabSetFieldNames(dm, user->nvars, fields);CHKERRQ(ierr);
117c4762a1bSJed Brown   ierr = DMMoabSetBlockSize(dm, user->nvars);CHKERRQ(ierr);
118c4762a1bSJed Brown   ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   /* SetUp the data structures for DMMOAB */
121c4762a1bSJed Brown   ierr = DMSetUp(dm);CHKERRQ(ierr);
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   /*  Create timestepping solver context */
124c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
125c4762a1bSJed Brown   ierr = TSSetDM(ts, dm);CHKERRQ(ierr);
126c4762a1bSJed Brown   ierr = TSSetType(ts,TSARKIMEX);CHKERRQ(ierr);
127c4762a1bSJed Brown   ierr = TSSetEquationType(ts,TS_EQ_DAE_IMPLICIT_INDEX1);CHKERRQ(ierr);
128c4762a1bSJed Brown   ierr = DMSetMatType(dm,MATBAIJ);CHKERRQ(ierr);
129c4762a1bSJed Brown   ierr = DMCreateMatrix(dm,&J);CHKERRQ(ierr);
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts,NULL,FormRHSFunction,user);CHKERRQ(ierr);
132c4762a1bSJed Brown   ierr = TSSetIFunction(ts,NULL,FormIFunction,user);CHKERRQ(ierr);
133c4762a1bSJed Brown   ierr = TSSetIJacobian(ts,J,J,FormIJacobian,user);CHKERRQ(ierr);
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   ftime = 10.0;
136c4762a1bSJed Brown   ierr = TSSetMaxSteps(ts,user->ntsteps);CHKERRQ(ierr);
137c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr);
138c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141c4762a1bSJed Brown      Create the solution vector and set the initial conditions
142c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
143c4762a1bSJed Brown   ierr = DMCreateGlobalVector(dm, &X);CHKERRQ(ierr);
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   ierr = FormInitialSolution(ts,X,user);CHKERRQ(ierr);
146c4762a1bSJed Brown   ierr = TSSetSolution(ts,X);CHKERRQ(ierr);
147c4762a1bSJed Brown   ierr = VecGetSize(X,&mx);CHKERRQ(ierr);
148c4762a1bSJed Brown   hx = 1.0/(PetscReal)(mx/2-1);
149c4762a1bSJed Brown   dt = 0.4 * PetscSqr(hx) / user->alpha; /* Diffusive stability limit */
150c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153c4762a1bSJed Brown      Set runtime options
154c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158c4762a1bSJed Brown      Solve nonlinear system
159c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
160c4762a1bSJed Brown   ierr = TSSolve(ts,X);CHKERRQ(ierr);
161c4762a1bSJed Brown   ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
162c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
163c4762a1bSJed Brown   ierr = TSGetConvergedReason(ts,&reason);CHKERRQ(ierr);
164c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],ftime,steps);CHKERRQ(ierr);
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   if (user->io) {
167c4762a1bSJed Brown     /* Print the numerical solution to screen and then dump to file */
168c4762a1bSJed Brown     ierr = VecView(X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
169c4762a1bSJed Brown 
170c4762a1bSJed Brown     /* Write out the solution along with the mesh */
171c4762a1bSJed Brown     ierr = DMMoabSetGlobalFieldVector(dm, X);CHKERRQ(ierr);
172c4762a1bSJed Brown #ifdef MOAB_HAVE_HDF5
173c4762a1bSJed Brown     ierr = DMMoabOutput(dm, "ex35.h5m", "");CHKERRQ(ierr);
174c4762a1bSJed Brown #else
175c4762a1bSJed Brown     /* MOAB does not support true parallel writers that aren't HDF5 based
176c4762a1bSJed Brown        And so if you are using VTK as the output format in parallel,
177c4762a1bSJed Brown        the data could be jumbled due to the order in which the processors
178c4762a1bSJed Brown        write out their parts of the mesh and solution tags
179c4762a1bSJed Brown     */
180c4762a1bSJed Brown     ierr = DMMoabOutput(dm, "ex35.vtk", "");CHKERRQ(ierr);
181c4762a1bSJed Brown #endif
182c4762a1bSJed Brown   }
183c4762a1bSJed Brown 
184c4762a1bSJed Brown   /* Free work space.
185c4762a1bSJed Brown      Free all PETSc related resources: */
186c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
187c4762a1bSJed Brown   ierr = VecDestroy(&X);CHKERRQ(ierr);
188c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
189c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
190c4762a1bSJed Brown 
191c4762a1bSJed Brown   /* Free all MOAB related resources: */
192c4762a1bSJed Brown   ierr = Destroy_AppContext(&user);CHKERRQ(ierr);
193c4762a1bSJed Brown   ierr = PetscFinalize();
194c4762a1bSJed Brown   return ierr;
195c4762a1bSJed Brown }
196c4762a1bSJed Brown 
197c4762a1bSJed Brown /*
198c4762a1bSJed Brown   IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
199c4762a1bSJed Brown */
200c4762a1bSJed Brown PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ptr)
201c4762a1bSJed Brown {
202c4762a1bSJed Brown   UserCtx             user = (UserCtx)ptr;
203c4762a1bSJed Brown   PetscErrorCode      ierr;
204c4762a1bSJed Brown   PetscInt            dof;
205c4762a1bSJed Brown   PetscReal           hx;
206c4762a1bSJed Brown   DM                  dm;
207c4762a1bSJed Brown   const moab::Range   *vlocal;
208c4762a1bSJed Brown   PetscBool           vonboundary;
209c4762a1bSJed Brown 
210c4762a1bSJed Brown   PetscFunctionBegin;
211c4762a1bSJed Brown   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
212c4762a1bSJed Brown 
213c4762a1bSJed Brown   /* get the essential MOAB mesh related quantities needed for FEM assembly */
214c4762a1bSJed Brown   ierr = DMMoabGetLocalVertices(dm, &vlocal, NULL);CHKERRQ(ierr);
215c4762a1bSJed Brown 
216c4762a1bSJed Brown   /* compute local element sizes - structured grid */
217c4762a1bSJed Brown   hx = 1.0/user->n;
218c4762a1bSJed Brown 
219c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid
220c4762a1bSJed Brown      Assemble the operator by looping over edges and computing
221c4762a1bSJed Brown      contribution for each vertex dof                         */
222c4762a1bSJed Brown   for (moab::Range::iterator iter = vlocal->begin(); iter != vlocal->end(); iter++) {
223c4762a1bSJed Brown     const moab::EntityHandle vhandle = *iter;
224c4762a1bSJed Brown 
225c4762a1bSJed Brown     ierr = DMMoabGetDofsBlocked(dm, 1, &vhandle, &dof);CHKERRQ(ierr);
226c4762a1bSJed Brown 
227c4762a1bSJed Brown     /* check if vertex is on the boundary */
228c4762a1bSJed Brown     ierr = DMMoabIsEntityOnBoundary(dm,vhandle,&vonboundary);CHKERRQ(ierr);
229c4762a1bSJed Brown 
230c4762a1bSJed Brown     if (vonboundary) {
231c4762a1bSJed Brown       const PetscScalar bcvals[2][2] = {{hx,0},{0,hx}};
232c4762a1bSJed Brown       ierr = MatSetValuesBlocked(Jpre,1,&dof,1,&dof,&bcvals[0][0],INSERT_VALUES);CHKERRQ(ierr);
233c4762a1bSJed Brown     }
234c4762a1bSJed Brown     else {
235c4762a1bSJed Brown       const PetscInt    row           = dof,col[] = {dof-1,dof,dof+1};
236c4762a1bSJed Brown       const PetscScalar dxxL          = -user->alpha/hx,dxx0 = 2.*user->alpha/hx,dxxR = -user->alpha/hx;
237c4762a1bSJed Brown       const PetscScalar vals[2][3][2] = {{{dxxL,0},{a *hx+dxx0,0},{dxxR,0}},
238c4762a1bSJed Brown                                          {{0,dxxL},{0,a*hx+dxx0},{0,dxxR}}};
239c4762a1bSJed Brown       ierr = MatSetValuesBlocked(Jpre,1,&row,3,col,&vals[0][0][0],INSERT_VALUES);CHKERRQ(ierr);
240c4762a1bSJed Brown     }
241c4762a1bSJed Brown   }
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   ierr = MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
244c4762a1bSJed Brown   ierr = MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
245c4762a1bSJed Brown   if (J != Jpre) {
246c4762a1bSJed Brown     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
247c4762a1bSJed Brown     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
248c4762a1bSJed Brown   }
249c4762a1bSJed Brown   PetscFunctionReturn(0);
250c4762a1bSJed Brown }
251c4762a1bSJed Brown 
252c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
253c4762a1bSJed Brown {
254c4762a1bSJed Brown   UserCtx           user = (UserCtx)ptr;
255c4762a1bSJed Brown   DM                dm;
256c4762a1bSJed Brown   PetscReal         hx;
257c4762a1bSJed Brown   const Field       *x;
258c4762a1bSJed Brown   Field             *f;
259c4762a1bSJed Brown   PetscInt          dof;
260c4762a1bSJed Brown   const moab::Range *ownedvtx;
261c4762a1bSJed Brown   PetscErrorCode    ierr;
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   PetscFunctionBegin;
264c4762a1bSJed Brown   hx = 1.0/user->n;
265c4762a1bSJed Brown   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
266c4762a1bSJed Brown 
267c4762a1bSJed Brown   /* Get pointers to vector data */
268c4762a1bSJed Brown   ierr = VecSet(F,0.0);CHKERRQ(ierr);
269c4762a1bSJed Brown 
270c4762a1bSJed Brown   ierr = DMMoabVecGetArrayRead(dm, X, &x);CHKERRQ(ierr);
271c4762a1bSJed Brown   ierr = DMMoabVecGetArray(dm, F, &f);CHKERRQ(ierr);
272c4762a1bSJed Brown 
273c4762a1bSJed Brown   ierr = DMMoabGetLocalVertices(dm, &ownedvtx, NULL);CHKERRQ(ierr);
274c4762a1bSJed Brown 
275c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
276c4762a1bSJed Brown   for (moab::Range::iterator iter = ownedvtx->begin(); iter != ownedvtx->end(); iter++) {
277c4762a1bSJed Brown     const moab::EntityHandle vhandle = *iter;
278c4762a1bSJed Brown     ierr = DMMoabGetDofsBlockedLocal(dm, 1, &vhandle, &dof);CHKERRQ(ierr);
279c4762a1bSJed Brown 
280c4762a1bSJed Brown     PetscScalar u = x[dof].u,v = x[dof].v;
281c4762a1bSJed Brown     f[dof].u = hx*(user->A + u*u*v - (user->B+1)*u);
282c4762a1bSJed Brown     f[dof].v = hx*(user->B*u - u*u*v);
283c4762a1bSJed Brown   }
284c4762a1bSJed Brown 
285c4762a1bSJed Brown   /* Restore vectors */
286c4762a1bSJed Brown   ierr = DMMoabVecRestoreArrayRead(dm, X, &x);CHKERRQ(ierr);
287c4762a1bSJed Brown   ierr = DMMoabVecRestoreArray(dm, F, &f);CHKERRQ(ierr);
288c4762a1bSJed Brown   PetscFunctionReturn(0);
289c4762a1bSJed Brown }
290c4762a1bSJed Brown 
291c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
292c4762a1bSJed Brown {
293c4762a1bSJed Brown   UserCtx         user = (UserCtx)ctx;
294c4762a1bSJed Brown   DM              dm;
295c4762a1bSJed Brown   Field           *x,*xdot,*f;
296c4762a1bSJed Brown   PetscReal       hx;
297c4762a1bSJed Brown   Vec             Xloc;
298c4762a1bSJed Brown   PetscErrorCode  ierr;
299c4762a1bSJed Brown   PetscInt        i,bcindx;
300c4762a1bSJed Brown   PetscBool       elem_on_boundary;
301c4762a1bSJed Brown   const moab::Range   *vlocal;
302c4762a1bSJed Brown 
303c4762a1bSJed Brown   PetscFunctionBegin;
304c4762a1bSJed Brown   hx = 1.0/user->n;
305c4762a1bSJed Brown   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
306c4762a1bSJed Brown 
307c4762a1bSJed Brown   /* get the essential MOAB mesh related quantities needed for FEM assembly */
308c4762a1bSJed Brown   ierr = DMMoabGetLocalVertices(dm, &vlocal, NULL);CHKERRQ(ierr);
309c4762a1bSJed Brown 
310c4762a1bSJed Brown   /* reset the residual vector */
311c4762a1bSJed Brown   ierr = VecSet(F,0.0);CHKERRQ(ierr);
312c4762a1bSJed Brown 
313c4762a1bSJed Brown   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
314c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
315c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
316c4762a1bSJed Brown 
317c4762a1bSJed Brown   /* get the local representation of the arrays from Vectors */
318c4762a1bSJed Brown   ierr = DMMoabVecGetArrayRead(dm, Xloc, &x);CHKERRQ(ierr);
319c4762a1bSJed Brown   ierr = DMMoabVecGetArrayRead(dm, Xdot, &xdot);CHKERRQ(ierr);
320c4762a1bSJed Brown   ierr = DMMoabVecGetArray(dm, F, &f);CHKERRQ(ierr);
321c4762a1bSJed Brown 
322c4762a1bSJed Brown   /* loop over local elements */
323c4762a1bSJed Brown   for (moab::Range::iterator iter = vlocal->begin(); iter != vlocal->end(); iter++) {
324c4762a1bSJed Brown     const moab::EntityHandle vhandle = *iter;
325c4762a1bSJed Brown 
326c4762a1bSJed Brown     ierr = DMMoabGetDofsBlockedLocal(dm,1,&vhandle,&i);CHKERRQ(ierr);
327c4762a1bSJed Brown 
328c4762a1bSJed Brown     /* check if vertex is on the boundary */
329c4762a1bSJed Brown     ierr = DMMoabIsEntityOnBoundary(dm,vhandle,&elem_on_boundary);CHKERRQ(ierr);
330c4762a1bSJed Brown 
331c4762a1bSJed Brown     if (elem_on_boundary) {
332c4762a1bSJed Brown       ierr = DMMoabGetDofsBlocked(dm, 1, &vhandle, &bcindx);CHKERRQ(ierr);
333c4762a1bSJed Brown       if (bcindx == 0) {  /* Apply left BC */
334c4762a1bSJed Brown         f[i].u = hx * (x[i].u - user->leftbc.u);
335c4762a1bSJed Brown         f[i].v = hx * (x[i].v - user->leftbc.v);
336c4762a1bSJed Brown       } else {       /* Apply right BC */
337c4762a1bSJed Brown         f[i].u = hx * (x[i].u - user->rightbc.u);
338c4762a1bSJed Brown         f[i].v = hx * (x[i].v - user->rightbc.v);
339c4762a1bSJed Brown       }
340c4762a1bSJed Brown     }
341c4762a1bSJed Brown     else {
342c4762a1bSJed Brown       f[i].u = hx * xdot[i].u - user->alpha * (x[i-1].u - 2.*x[i].u + x[i+1].u) / hx;
343c4762a1bSJed Brown       f[i].v = hx * xdot[i].v - user->alpha * (x[i-1].v - 2.*x[i].v + x[i+1].v) / hx;
344c4762a1bSJed Brown     }
345c4762a1bSJed Brown   }
346c4762a1bSJed Brown 
347c4762a1bSJed Brown   /* Restore data */
348c4762a1bSJed Brown   ierr = DMMoabVecRestoreArrayRead(dm, Xloc, &x);CHKERRQ(ierr);
349c4762a1bSJed Brown   ierr = DMMoabVecRestoreArrayRead(dm, Xdot, &xdot);CHKERRQ(ierr);
350c4762a1bSJed Brown   ierr = DMMoabVecRestoreArray(dm, F, &f);CHKERRQ(ierr);
351c4762a1bSJed Brown   ierr = DMRestoreLocalVector(dm, &Xloc);CHKERRQ(ierr);
352c4762a1bSJed Brown   PetscFunctionReturn(0);
353c4762a1bSJed Brown }
354c4762a1bSJed Brown 
355c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
356c4762a1bSJed Brown {
357c4762a1bSJed Brown   UserCtx           user = (UserCtx)ctx;
358c4762a1bSJed Brown   PetscReal         vpos[3];
359c4762a1bSJed Brown   DM                dm;
360c4762a1bSJed Brown   Field             *x;
361c4762a1bSJed Brown   PetscErrorCode    ierr;
362c4762a1bSJed Brown   const moab::Range *vowned;
363c4762a1bSJed Brown   PetscInt          dof;
364c4762a1bSJed Brown   moab::Range::iterator iter;
365c4762a1bSJed Brown 
366c4762a1bSJed Brown   PetscFunctionBegin;
367c4762a1bSJed Brown   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
368c4762a1bSJed Brown 
369c4762a1bSJed Brown   /* get the essential MOAB mesh related quantities needed for FEM assembly */
370c4762a1bSJed Brown   ierr = DMMoabGetLocalVertices(dm, &vowned, NULL);CHKERRQ(ierr);
371c4762a1bSJed Brown 
372c4762a1bSJed Brown   ierr = VecSet(X, 0.0);CHKERRQ(ierr);
373c4762a1bSJed Brown 
374c4762a1bSJed Brown   /* Get pointers to vector data */
375c4762a1bSJed Brown   ierr = DMMoabVecGetArray(dm, X, &x);CHKERRQ(ierr);
376c4762a1bSJed Brown 
377c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
378c4762a1bSJed Brown   for (moab::Range::iterator iter = vowned->begin(); iter != vowned->end(); iter++) {
379c4762a1bSJed Brown     const moab::EntityHandle vhandle = *iter;
380c4762a1bSJed Brown     ierr = DMMoabGetDofsBlockedLocal(dm, 1, &vhandle, &dof);CHKERRQ(ierr);
381c4762a1bSJed Brown 
382c4762a1bSJed Brown     /* compute the mid-point of the element and use a 1-point lumped quadrature */
383c4762a1bSJed Brown     ierr = DMMoabGetVertexCoordinates(dm,1,&vhandle,vpos);CHKERRQ(ierr);
384c4762a1bSJed Brown 
385c4762a1bSJed Brown     PetscReal xi = vpos[0];
386c4762a1bSJed Brown     x[dof].u = user->leftbc.u*(1.-xi) + user->rightbc.u*xi + PetscSinReal(2.*PETSC_PI*xi);
387c4762a1bSJed Brown     x[dof].v = user->leftbc.v*(1.-xi) + user->rightbc.v*xi;
388c4762a1bSJed Brown   }
389c4762a1bSJed Brown 
390c4762a1bSJed Brown   /* Restore vectors */
391c4762a1bSJed Brown   ierr = DMMoabVecRestoreArray(dm, X, &x);CHKERRQ(ierr);
392c4762a1bSJed Brown   PetscFunctionReturn(0);
393c4762a1bSJed Brown }
394c4762a1bSJed Brown 
395c4762a1bSJed Brown /*TEST
396c4762a1bSJed Brown 
397c4762a1bSJed Brown     build:
398c4762a1bSJed Brown       requires: moab
399c4762a1bSJed Brown 
400c4762a1bSJed Brown     test:
401c4762a1bSJed Brown       args: -n 20 -ts_type rosw -ts_rosw_type 2p -ts_dt 5e-2 -ts_adapt_type none
402c4762a1bSJed Brown 
403c4762a1bSJed Brown     test:
404c4762a1bSJed Brown       suffix: 2
405c4762a1bSJed Brown       nsize: 2
406c4762a1bSJed Brown       args: -n 50 -ts_type glee -ts_adapt_type none -ts_dt 0.1 -io
407c4762a1bSJed Brown       TODO:
408c4762a1bSJed Brown 
409c4762a1bSJed Brown TEST*/
410