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 41c4762a1bSJed Brown PetscFunctionBegin; 429566063dSJacob Faibussowitsch PetscCall(PetscNew(&user)); 43*d0609cedSBarry 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 } 68*d0609cedSBarry Smith PetscOptionsEnd(); 69c4762a1bSJed Brown 70c4762a1bSJed Brown *puser = user; 71c4762a1bSJed Brown PetscFunctionReturn(0); 72c4762a1bSJed Brown } 73c4762a1bSJed Brown 74c4762a1bSJed Brown PetscErrorCode Destroy_AppContext(UserCtx *user) 75c4762a1bSJed Brown { 76c4762a1bSJed Brown PetscFunctionBegin; 779566063dSJacob Faibussowitsch PetscCall(PetscFree(*user)); 78c4762a1bSJed Brown PetscFunctionReturn(0); 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 ****************/ 91c4762a1bSJed Brown int main(int argc,char **argv) 92c4762a1bSJed Brown { 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 1039566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char *)0,help)); 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* Initialize the user context struct */ 1069566063dSJacob Faibussowitsch PetscCall(Initialize_AppContext(&user)); 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* Fill in the user defined work context: */ 1099566063dSJacob Faibussowitsch PetscCall(DMMoabCreateBoxMesh(PETSC_COMM_WORLD, 1, PETSC_FALSE, NULL, user->n, 1, &dm)); 1109566063dSJacob Faibussowitsch PetscCall(DMMoabSetFieldNames(dm, user->nvars, fields)); 1119566063dSJacob Faibussowitsch PetscCall(DMMoabSetBlockSize(dm, user->nvars)); 1129566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* SetUp the data structures for DMMOAB */ 1159566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* Create timestepping solver context */ 1189566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 1199566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, dm)); 1209566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSARKIMEX)); 1219566063dSJacob Faibussowitsch PetscCall(TSSetEquationType(ts,TS_EQ_DAE_IMPLICIT_INDEX1)); 1229566063dSJacob Faibussowitsch PetscCall(DMSetMatType(dm,MATBAIJ)); 1239566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm,&J)); 124c4762a1bSJed Brown 1259566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,NULL,FormRHSFunction,user)); 1269566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts,NULL,FormIFunction,user)); 1279566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts,J,J,FormIJacobian,user)); 128c4762a1bSJed Brown 129c4762a1bSJed Brown ftime = 10.0; 1309566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts,user->ntsteps)); 1319566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,ftime)); 1329566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 135c4762a1bSJed Brown Create the solution vector and set the initial conditions 136c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1379566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &X)); 138c4762a1bSJed Brown 1399566063dSJacob Faibussowitsch PetscCall(FormInitialSolution(ts,X,user)); 1409566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts,X)); 1419566063dSJacob Faibussowitsch PetscCall(VecGetSize(X,&mx)); 142c4762a1bSJed Brown hx = 1.0/(PetscReal)(mx/2-1); 143c4762a1bSJed Brown dt = 0.4 * PetscSqr(hx) / user->alpha; /* Diffusive stability limit */ 1449566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,dt)); 145c4762a1bSJed Brown 146c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 147c4762a1bSJed Brown Set runtime options 148c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1499566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 150c4762a1bSJed Brown 151c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 152c4762a1bSJed Brown Solve nonlinear system 153c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1549566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,X)); 1559566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts,&ftime)); 1569566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&steps)); 1579566063dSJacob Faibussowitsch PetscCall(TSGetConvergedReason(ts,&reason)); 1589566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],ftime,steps)); 159c4762a1bSJed Brown 160c4762a1bSJed Brown if (user->io) { 161c4762a1bSJed Brown /* Print the numerical solution to screen and then dump to file */ 1629566063dSJacob Faibussowitsch PetscCall(VecView(X,PETSC_VIEWER_STDOUT_WORLD)); 163c4762a1bSJed Brown 164c4762a1bSJed Brown /* Write out the solution along with the mesh */ 1659566063dSJacob Faibussowitsch PetscCall(DMMoabSetGlobalFieldVector(dm, X)); 166c4762a1bSJed Brown #ifdef MOAB_HAVE_HDF5 1679566063dSJacob Faibussowitsch PetscCall(DMMoabOutput(dm, "ex35.h5m", "")); 168c4762a1bSJed Brown #else 169c4762a1bSJed Brown /* MOAB does not support true parallel writers that aren't HDF5 based 170c4762a1bSJed Brown And so if you are using VTK as the output format in parallel, 171c4762a1bSJed Brown the data could be jumbled due to the order in which the processors 172c4762a1bSJed Brown write out their parts of the mesh and solution tags 173c4762a1bSJed Brown */ 1749566063dSJacob Faibussowitsch PetscCall(DMMoabOutput(dm, "ex35.vtk", "")); 175c4762a1bSJed Brown #endif 176c4762a1bSJed Brown } 177c4762a1bSJed Brown 178c4762a1bSJed Brown /* Free work space. 179c4762a1bSJed Brown Free all PETSc related resources: */ 1809566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 1819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 1829566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1839566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 184c4762a1bSJed Brown 185c4762a1bSJed Brown /* Free all MOAB related resources: */ 1869566063dSJacob Faibussowitsch PetscCall(Destroy_AppContext(&user)); 1879566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 188b122ec5aSJacob Faibussowitsch return 0; 189c4762a1bSJed Brown } 190c4762a1bSJed Brown 191c4762a1bSJed Brown /* 192c4762a1bSJed Brown IJacobian - Compute IJacobian = dF/dU + a dF/dUdot 193c4762a1bSJed Brown */ 194c4762a1bSJed Brown PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ptr) 195c4762a1bSJed Brown { 196c4762a1bSJed Brown UserCtx user = (UserCtx)ptr; 197c4762a1bSJed Brown PetscInt dof; 198c4762a1bSJed Brown PetscReal hx; 199c4762a1bSJed Brown DM dm; 200c4762a1bSJed Brown const moab::Range *vlocal; 201c4762a1bSJed Brown PetscBool vonboundary; 202c4762a1bSJed Brown 203c4762a1bSJed Brown PetscFunctionBegin; 2049566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 205c4762a1bSJed Brown 206c4762a1bSJed Brown /* get the essential MOAB mesh related quantities needed for FEM assembly */ 2079566063dSJacob Faibussowitsch PetscCall(DMMoabGetLocalVertices(dm, &vlocal, NULL)); 208c4762a1bSJed Brown 209c4762a1bSJed Brown /* compute local element sizes - structured grid */ 210c4762a1bSJed Brown hx = 1.0/user->n; 211c4762a1bSJed Brown 212c4762a1bSJed Brown /* Compute function over the locally owned part of the grid 213c4762a1bSJed Brown Assemble the operator by looping over edges and computing 214c4762a1bSJed Brown contribution for each vertex dof */ 215c4762a1bSJed Brown for (moab::Range::iterator iter = vlocal->begin(); iter != vlocal->end(); iter++) { 216c4762a1bSJed Brown const moab::EntityHandle vhandle = *iter; 217c4762a1bSJed Brown 2189566063dSJacob Faibussowitsch PetscCall(DMMoabGetDofsBlocked(dm, 1, &vhandle, &dof)); 219c4762a1bSJed Brown 220c4762a1bSJed Brown /* check if vertex is on the boundary */ 2219566063dSJacob Faibussowitsch PetscCall(DMMoabIsEntityOnBoundary(dm,vhandle,&vonboundary)); 222c4762a1bSJed Brown 223c4762a1bSJed Brown if (vonboundary) { 224c4762a1bSJed Brown const PetscScalar bcvals[2][2] = {{hx,0},{0,hx}}; 2259566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(Jpre,1,&dof,1,&dof,&bcvals[0][0],INSERT_VALUES)); 226c4762a1bSJed Brown } 227c4762a1bSJed Brown else { 228c4762a1bSJed Brown const PetscInt row = dof,col[] = {dof-1,dof,dof+1}; 229c4762a1bSJed Brown const PetscScalar dxxL = -user->alpha/hx,dxx0 = 2.*user->alpha/hx,dxxR = -user->alpha/hx; 230c4762a1bSJed Brown const PetscScalar vals[2][3][2] = {{{dxxL,0},{a *hx+dxx0,0},{dxxR,0}}, 231c4762a1bSJed Brown {{0,dxxL},{0,a*hx+dxx0},{0,dxxR}}}; 2329566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(Jpre,1,&row,3,col,&vals[0][0][0],INSERT_VALUES)); 233c4762a1bSJed Brown } 234c4762a1bSJed Brown } 235c4762a1bSJed Brown 2369566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY)); 2379566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY)); 238c4762a1bSJed Brown if (J != Jpre) { 2399566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 2409566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 241c4762a1bSJed Brown } 242c4762a1bSJed Brown PetscFunctionReturn(0); 243c4762a1bSJed Brown } 244c4762a1bSJed Brown 245c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr) 246c4762a1bSJed Brown { 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 283c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx) 284c4762a1bSJed Brown { 285c4762a1bSJed Brown UserCtx user = (UserCtx)ctx; 286c4762a1bSJed Brown DM dm; 287c4762a1bSJed Brown Field *x,*xdot,*f; 288c4762a1bSJed Brown PetscReal hx; 289c4762a1bSJed Brown Vec Xloc; 290c4762a1bSJed Brown PetscInt i,bcindx; 291c4762a1bSJed Brown PetscBool elem_on_boundary; 292c4762a1bSJed Brown const moab::Range *vlocal; 293c4762a1bSJed Brown 294c4762a1bSJed Brown PetscFunctionBegin; 295c4762a1bSJed Brown hx = 1.0/user->n; 2969566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 297c4762a1bSJed Brown 298c4762a1bSJed Brown /* get the essential MOAB mesh related quantities needed for FEM assembly */ 2999566063dSJacob Faibussowitsch PetscCall(DMMoabGetLocalVertices(dm, &vlocal, NULL)); 300c4762a1bSJed Brown 301c4762a1bSJed Brown /* reset the residual vector */ 3029566063dSJacob Faibussowitsch PetscCall(VecSet(F,0.0)); 303c4762a1bSJed Brown 3049566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm,&Xloc)); 3059566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc)); 3069566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc)); 307c4762a1bSJed Brown 308c4762a1bSJed Brown /* get the local representation of the arrays from Vectors */ 3099566063dSJacob Faibussowitsch PetscCall(DMMoabVecGetArrayRead(dm, Xloc, &x)); 3109566063dSJacob Faibussowitsch PetscCall(DMMoabVecGetArrayRead(dm, Xdot, &xdot)); 3119566063dSJacob Faibussowitsch PetscCall(DMMoabVecGetArray(dm, F, &f)); 312c4762a1bSJed Brown 313c4762a1bSJed Brown /* loop over local elements */ 314c4762a1bSJed Brown for (moab::Range::iterator iter = vlocal->begin(); iter != vlocal->end(); iter++) { 315c4762a1bSJed Brown const moab::EntityHandle vhandle = *iter; 316c4762a1bSJed Brown 3179566063dSJacob Faibussowitsch PetscCall(DMMoabGetDofsBlockedLocal(dm,1,&vhandle,&i)); 318c4762a1bSJed Brown 319c4762a1bSJed Brown /* check if vertex is on the boundary */ 3209566063dSJacob Faibussowitsch PetscCall(DMMoabIsEntityOnBoundary(dm,vhandle,&elem_on_boundary)); 321c4762a1bSJed Brown 322c4762a1bSJed Brown if (elem_on_boundary) { 3239566063dSJacob Faibussowitsch PetscCall(DMMoabGetDofsBlocked(dm, 1, &vhandle, &bcindx)); 324c4762a1bSJed Brown if (bcindx == 0) { /* Apply left BC */ 325c4762a1bSJed Brown f[i].u = hx * (x[i].u - user->leftbc.u); 326c4762a1bSJed Brown f[i].v = hx * (x[i].v - user->leftbc.v); 327c4762a1bSJed Brown } else { /* Apply right BC */ 328c4762a1bSJed Brown f[i].u = hx * (x[i].u - user->rightbc.u); 329c4762a1bSJed Brown f[i].v = hx * (x[i].v - user->rightbc.v); 330c4762a1bSJed Brown } 331c4762a1bSJed Brown } 332c4762a1bSJed Brown else { 333c4762a1bSJed Brown f[i].u = hx * xdot[i].u - user->alpha * (x[i-1].u - 2.*x[i].u + x[i+1].u) / hx; 334c4762a1bSJed Brown f[i].v = hx * xdot[i].v - user->alpha * (x[i-1].v - 2.*x[i].v + x[i+1].v) / hx; 335c4762a1bSJed Brown } 336c4762a1bSJed Brown } 337c4762a1bSJed Brown 338c4762a1bSJed Brown /* Restore data */ 3399566063dSJacob Faibussowitsch PetscCall(DMMoabVecRestoreArrayRead(dm, Xloc, &x)); 3409566063dSJacob Faibussowitsch PetscCall(DMMoabVecRestoreArrayRead(dm, Xdot, &xdot)); 3419566063dSJacob Faibussowitsch PetscCall(DMMoabVecRestoreArray(dm, F, &f)); 3429566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &Xloc)); 343c4762a1bSJed Brown PetscFunctionReturn(0); 344c4762a1bSJed Brown } 345c4762a1bSJed Brown 346c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx) 347c4762a1bSJed Brown { 348c4762a1bSJed Brown UserCtx user = (UserCtx)ctx; 349c4762a1bSJed Brown PetscReal vpos[3]; 350c4762a1bSJed Brown DM dm; 351c4762a1bSJed Brown Field *x; 352c4762a1bSJed Brown const moab::Range *vowned; 353c4762a1bSJed Brown PetscInt dof; 354c4762a1bSJed Brown moab::Range::iterator iter; 355c4762a1bSJed Brown 356c4762a1bSJed Brown PetscFunctionBegin; 3579566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 358c4762a1bSJed Brown 359c4762a1bSJed Brown /* get the essential MOAB mesh related quantities needed for FEM assembly */ 3609566063dSJacob Faibussowitsch PetscCall(DMMoabGetLocalVertices(dm, &vowned, NULL)); 361c4762a1bSJed Brown 3629566063dSJacob Faibussowitsch PetscCall(VecSet(X, 0.0)); 363c4762a1bSJed Brown 364c4762a1bSJed Brown /* Get pointers to vector data */ 3659566063dSJacob Faibussowitsch PetscCall(DMMoabVecGetArray(dm, X, &x)); 366c4762a1bSJed Brown 367c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 368c4762a1bSJed Brown for (moab::Range::iterator iter = vowned->begin(); iter != vowned->end(); iter++) { 369c4762a1bSJed Brown const moab::EntityHandle vhandle = *iter; 3709566063dSJacob Faibussowitsch PetscCall(DMMoabGetDofsBlockedLocal(dm, 1, &vhandle, &dof)); 371c4762a1bSJed Brown 372c4762a1bSJed Brown /* compute the mid-point of the element and use a 1-point lumped quadrature */ 3739566063dSJacob Faibussowitsch PetscCall(DMMoabGetVertexCoordinates(dm,1,&vhandle,vpos)); 374c4762a1bSJed Brown 375c4762a1bSJed Brown PetscReal xi = vpos[0]; 376c4762a1bSJed Brown x[dof].u = user->leftbc.u*(1.-xi) + user->rightbc.u*xi + PetscSinReal(2.*PETSC_PI*xi); 377c4762a1bSJed Brown x[dof].v = user->leftbc.v*(1.-xi) + user->rightbc.v*xi; 378c4762a1bSJed Brown } 379c4762a1bSJed Brown 380c4762a1bSJed Brown /* Restore vectors */ 3819566063dSJacob Faibussowitsch PetscCall(DMMoabVecRestoreArray(dm, X, &x)); 382c4762a1bSJed Brown PetscFunctionReturn(0); 383c4762a1bSJed Brown } 384c4762a1bSJed Brown 385c4762a1bSJed Brown /*TEST 386c4762a1bSJed Brown 387c4762a1bSJed Brown build: 388c4762a1bSJed Brown requires: moab 389c4762a1bSJed Brown 390c4762a1bSJed Brown test: 391c4762a1bSJed Brown args: -n 20 -ts_type rosw -ts_rosw_type 2p -ts_dt 5e-2 -ts_adapt_type none 392c4762a1bSJed Brown 393c4762a1bSJed Brown test: 394c4762a1bSJed Brown suffix: 2 395c4762a1bSJed Brown nsize: 2 396c4762a1bSJed Brown args: -n 50 -ts_type glee -ts_adapt_type none -ts_dt 0.1 -io 397c4762a1bSJed Brown TODO: 398c4762a1bSJed Brown 399c4762a1bSJed Brown TEST*/ 400