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; 435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&user)); 44c4762a1bSJed Brown 451e1ea65dSPierre 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; 585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-A","Reaction rate","ex35.cxx",user->A,&user->A,NULL)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-B","Reaction rate","ex35.cxx",user->B,&user->B,NULL)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-alpha","Diffusion coefficient","ex35.cxx",user->alpha,&user->alpha,NULL)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsScalar("-uleft","Dirichlet boundary condition","ex35.cxx",user->leftbc.u,&user->leftbc.u,NULL)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsScalar("-uright","Dirichlet boundary condition","ex35.cxx",user->rightbc.u,&user->rightbc.u,NULL)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsScalar("-vleft","Dirichlet boundary condition","ex35.cxx",user->leftbc.v,&user->leftbc.v,NULL)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsScalar("-vright","Dirichlet boundary condition","ex35.cxx",user->rightbc.v,&user->rightbc.v,NULL)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-n","Number of 1-D elements","ex35.cxx",user->n,&user->n,NULL)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-ndt","Number of time steps","ex35.cxx",user->ntsteps,&user->ntsteps,NULL)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-io","Write the mesh and solution output to a file.","ex35.cxx",user->io,&user->io,NULL)); 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 PetscFunctionBegin; 795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(*user)); 80c4762a1bSJed Brown PetscFunctionReturn(0); 81c4762a1bSJed Brown } 82c4762a1bSJed Brown 83c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS,Vec,void*); 84c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*); 85c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*); 86c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 87c4762a1bSJed Brown 88c4762a1bSJed Brown /**************** 89c4762a1bSJed Brown * * 90c4762a1bSJed Brown * MAIN * 91c4762a1bSJed Brown * * 92c4762a1bSJed Brown ****************/ 93c4762a1bSJed Brown int main(int argc,char **argv) 94c4762a1bSJed Brown { 95c4762a1bSJed Brown TS ts; /* nonlinear solver */ 96c4762a1bSJed Brown Vec X; /* solution, residual vectors */ 97c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 98c4762a1bSJed Brown PetscInt steps,mx; 99c4762a1bSJed Brown PetscReal hx,dt,ftime; 100c4762a1bSJed Brown UserCtx user; /* user-defined work context */ 101c4762a1bSJed Brown TSConvergedReason reason; 102c4762a1bSJed Brown DM dm; 103c4762a1bSJed Brown const char *fields[2] = {"U","V"}; 104c4762a1bSJed Brown 105*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char *)0,help)); 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* Initialize the user context struct */ 1085f80ce2aSJacob Faibussowitsch CHKERRQ(Initialize_AppContext(&user)); 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* Fill in the user defined work context: */ 1115f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabCreateBoxMesh(PETSC_COMM_WORLD, 1, PETSC_FALSE, NULL, user->n, 1, &dm)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabSetFieldNames(dm, user->nvars, fields)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabSetBlockSize(dm, user->nvars)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(dm)); 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* SetUp the data structures for DMMOAB */ 1175f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(dm)); 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* Create timestepping solver context */ 1205f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetDM(ts, dm)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSARKIMEX)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetEquationType(ts,TS_EQ_DAE_IMPLICIT_INDEX1)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetMatType(dm,MATBAIJ)); 1255f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(dm,&J)); 126c4762a1bSJed Brown 1275f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,NULL,FormRHSFunction,user)); 1285f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIFunction(ts,NULL,FormIFunction,user)); 1295f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIJacobian(ts,J,J,FormIJacobian,user)); 130c4762a1bSJed Brown 131c4762a1bSJed Brown ftime = 10.0; 1325f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxSteps(ts,user->ntsteps)); 1335f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,ftime)); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 135c4762a1bSJed Brown 136c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 137c4762a1bSJed Brown Create the solution vector and set the initial conditions 138c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1395f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dm, &X)); 140c4762a1bSJed Brown 1415f80ce2aSJacob Faibussowitsch CHKERRQ(FormInitialSolution(ts,X,user)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolution(ts,X)); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(X,&mx)); 144c4762a1bSJed Brown hx = 1.0/(PetscReal)(mx/2-1); 145c4762a1bSJed Brown dt = 0.4 * PetscSqr(hx) / user->alpha; /* Diffusive stability limit */ 1465f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,dt)); 147c4762a1bSJed Brown 148c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 149c4762a1bSJed Brown Set runtime options 150c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1515f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 152c4762a1bSJed Brown 153c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 154c4762a1bSJed Brown Solve nonlinear system 155c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1565f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,X)); 1575f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSolveTime(ts,&ftime)); 1585f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetStepNumber(ts,&steps)); 1595f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetConvergedReason(ts,&reason)); 1605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],ftime,steps)); 161c4762a1bSJed Brown 162c4762a1bSJed Brown if (user->io) { 163c4762a1bSJed Brown /* Print the numerical solution to screen and then dump to file */ 1645f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(X,PETSC_VIEWER_STDOUT_WORLD)); 165c4762a1bSJed Brown 166c4762a1bSJed Brown /* Write out the solution along with the mesh */ 1675f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabSetGlobalFieldVector(dm, X)); 168c4762a1bSJed Brown #ifdef MOAB_HAVE_HDF5 1695f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabOutput(dm, "ex35.h5m", "")); 170c4762a1bSJed Brown #else 171c4762a1bSJed Brown /* MOAB does not support true parallel writers that aren't HDF5 based 172c4762a1bSJed Brown And so if you are using VTK as the output format in parallel, 173c4762a1bSJed Brown the data could be jumbled due to the order in which the processors 174c4762a1bSJed Brown write out their parts of the mesh and solution tags 175c4762a1bSJed Brown */ 1765f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabOutput(dm, "ex35.vtk", "")); 177c4762a1bSJed Brown #endif 178c4762a1bSJed Brown } 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* Free work space. 181c4762a1bSJed Brown Free all PETSc related resources: */ 1825f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&J)); 1835f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&X)); 1845f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 1855f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 186c4762a1bSJed Brown 187c4762a1bSJed Brown /* Free all MOAB related resources: */ 1885f80ce2aSJacob Faibussowitsch CHKERRQ(Destroy_AppContext(&user)); 189*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 190*b122ec5aSJacob Faibussowitsch return 0; 191c4762a1bSJed Brown } 192c4762a1bSJed Brown 193c4762a1bSJed Brown /* 194c4762a1bSJed Brown IJacobian - Compute IJacobian = dF/dU + a dF/dUdot 195c4762a1bSJed Brown */ 196c4762a1bSJed Brown PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ptr) 197c4762a1bSJed Brown { 198c4762a1bSJed Brown UserCtx user = (UserCtx)ptr; 199c4762a1bSJed Brown PetscInt dof; 200c4762a1bSJed Brown PetscReal hx; 201c4762a1bSJed Brown DM dm; 202c4762a1bSJed Brown const moab::Range *vlocal; 203c4762a1bSJed Brown PetscBool vonboundary; 204c4762a1bSJed Brown 205c4762a1bSJed Brown PetscFunctionBegin; 2065f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &dm)); 207c4762a1bSJed Brown 208c4762a1bSJed Brown /* get the essential MOAB mesh related quantities needed for FEM assembly */ 2095f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabGetLocalVertices(dm, &vlocal, NULL)); 210c4762a1bSJed Brown 211c4762a1bSJed Brown /* compute local element sizes - structured grid */ 212c4762a1bSJed Brown hx = 1.0/user->n; 213c4762a1bSJed Brown 214c4762a1bSJed Brown /* Compute function over the locally owned part of the grid 215c4762a1bSJed Brown Assemble the operator by looping over edges and computing 216c4762a1bSJed Brown contribution for each vertex dof */ 217c4762a1bSJed Brown for (moab::Range::iterator iter = vlocal->begin(); iter != vlocal->end(); iter++) { 218c4762a1bSJed Brown const moab::EntityHandle vhandle = *iter; 219c4762a1bSJed Brown 2205f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabGetDofsBlocked(dm, 1, &vhandle, &dof)); 221c4762a1bSJed Brown 222c4762a1bSJed Brown /* check if vertex is on the boundary */ 2235f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabIsEntityOnBoundary(dm,vhandle,&vonboundary)); 224c4762a1bSJed Brown 225c4762a1bSJed Brown if (vonboundary) { 226c4762a1bSJed Brown const PetscScalar bcvals[2][2] = {{hx,0},{0,hx}}; 2275f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlocked(Jpre,1,&dof,1,&dof,&bcvals[0][0],INSERT_VALUES)); 228c4762a1bSJed Brown } 229c4762a1bSJed Brown else { 230c4762a1bSJed Brown const PetscInt row = dof,col[] = {dof-1,dof,dof+1}; 231c4762a1bSJed Brown const PetscScalar dxxL = -user->alpha/hx,dxx0 = 2.*user->alpha/hx,dxxR = -user->alpha/hx; 232c4762a1bSJed Brown const PetscScalar vals[2][3][2] = {{{dxxL,0},{a *hx+dxx0,0},{dxxR,0}}, 233c4762a1bSJed Brown {{0,dxxL},{0,a*hx+dxx0},{0,dxxR}}}; 2345f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlocked(Jpre,1,&row,3,col,&vals[0][0][0],INSERT_VALUES)); 235c4762a1bSJed Brown } 236c4762a1bSJed Brown } 237c4762a1bSJed Brown 2385f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY)); 2395f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY)); 240c4762a1bSJed Brown if (J != Jpre) { 2415f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 2425f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 243c4762a1bSJed Brown } 244c4762a1bSJed Brown PetscFunctionReturn(0); 245c4762a1bSJed Brown } 246c4762a1bSJed Brown 247c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr) 248c4762a1bSJed Brown { 249c4762a1bSJed Brown UserCtx user = (UserCtx)ptr; 250c4762a1bSJed Brown DM dm; 251c4762a1bSJed Brown PetscReal hx; 252c4762a1bSJed Brown const Field *x; 253c4762a1bSJed Brown Field *f; 254c4762a1bSJed Brown PetscInt dof; 255c4762a1bSJed Brown const moab::Range *ownedvtx; 256c4762a1bSJed Brown 257c4762a1bSJed Brown PetscFunctionBegin; 258c4762a1bSJed Brown hx = 1.0/user->n; 2595f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&dm)); 260c4762a1bSJed Brown 261c4762a1bSJed Brown /* Get pointers to vector data */ 2625f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(F,0.0)); 263c4762a1bSJed Brown 2645f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabVecGetArrayRead(dm, X, &x)); 2655f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabVecGetArray(dm, F, &f)); 266c4762a1bSJed Brown 2675f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabGetLocalVertices(dm, &ownedvtx, NULL)); 268c4762a1bSJed Brown 269c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 270c4762a1bSJed Brown for (moab::Range::iterator iter = ownedvtx->begin(); iter != ownedvtx->end(); iter++) { 271c4762a1bSJed Brown const moab::EntityHandle vhandle = *iter; 2725f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabGetDofsBlockedLocal(dm, 1, &vhandle, &dof)); 273c4762a1bSJed Brown 274c4762a1bSJed Brown PetscScalar u = x[dof].u,v = x[dof].v; 275c4762a1bSJed Brown f[dof].u = hx*(user->A + u*u*v - (user->B+1)*u); 276c4762a1bSJed Brown f[dof].v = hx*(user->B*u - u*u*v); 277c4762a1bSJed Brown } 278c4762a1bSJed Brown 279c4762a1bSJed Brown /* Restore vectors */ 2805f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabVecRestoreArrayRead(dm, X, &x)); 2815f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabVecRestoreArray(dm, F, &f)); 282c4762a1bSJed Brown PetscFunctionReturn(0); 283c4762a1bSJed Brown } 284c4762a1bSJed Brown 285c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx) 286c4762a1bSJed Brown { 287c4762a1bSJed Brown UserCtx user = (UserCtx)ctx; 288c4762a1bSJed Brown DM dm; 289c4762a1bSJed Brown Field *x,*xdot,*f; 290c4762a1bSJed Brown PetscReal hx; 291c4762a1bSJed Brown Vec Xloc; 292c4762a1bSJed Brown PetscInt i,bcindx; 293c4762a1bSJed Brown PetscBool elem_on_boundary; 294c4762a1bSJed Brown const moab::Range *vlocal; 295c4762a1bSJed Brown 296c4762a1bSJed Brown PetscFunctionBegin; 297c4762a1bSJed Brown hx = 1.0/user->n; 2985f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &dm)); 299c4762a1bSJed Brown 300c4762a1bSJed Brown /* get the essential MOAB mesh related quantities needed for FEM assembly */ 3015f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabGetLocalVertices(dm, &vlocal, NULL)); 302c4762a1bSJed Brown 303c4762a1bSJed Brown /* reset the residual vector */ 3045f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(F,0.0)); 305c4762a1bSJed Brown 3065f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dm,&Xloc)); 3075f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc)); 3085f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc)); 309c4762a1bSJed Brown 310c4762a1bSJed Brown /* get the local representation of the arrays from Vectors */ 3115f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabVecGetArrayRead(dm, Xloc, &x)); 3125f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabVecGetArrayRead(dm, Xdot, &xdot)); 3135f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabVecGetArray(dm, F, &f)); 314c4762a1bSJed Brown 315c4762a1bSJed Brown /* loop over local elements */ 316c4762a1bSJed Brown for (moab::Range::iterator iter = vlocal->begin(); iter != vlocal->end(); iter++) { 317c4762a1bSJed Brown const moab::EntityHandle vhandle = *iter; 318c4762a1bSJed Brown 3195f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabGetDofsBlockedLocal(dm,1,&vhandle,&i)); 320c4762a1bSJed Brown 321c4762a1bSJed Brown /* check if vertex is on the boundary */ 3225f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabIsEntityOnBoundary(dm,vhandle,&elem_on_boundary)); 323c4762a1bSJed Brown 324c4762a1bSJed Brown if (elem_on_boundary) { 3255f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabGetDofsBlocked(dm, 1, &vhandle, &bcindx)); 326c4762a1bSJed Brown if (bcindx == 0) { /* Apply left BC */ 327c4762a1bSJed Brown f[i].u = hx * (x[i].u - user->leftbc.u); 328c4762a1bSJed Brown f[i].v = hx * (x[i].v - user->leftbc.v); 329c4762a1bSJed Brown } else { /* Apply right BC */ 330c4762a1bSJed Brown f[i].u = hx * (x[i].u - user->rightbc.u); 331c4762a1bSJed Brown f[i].v = hx * (x[i].v - user->rightbc.v); 332c4762a1bSJed Brown } 333c4762a1bSJed Brown } 334c4762a1bSJed Brown else { 335c4762a1bSJed Brown f[i].u = hx * xdot[i].u - user->alpha * (x[i-1].u - 2.*x[i].u + x[i+1].u) / hx; 336c4762a1bSJed Brown f[i].v = hx * xdot[i].v - user->alpha * (x[i-1].v - 2.*x[i].v + x[i+1].v) / hx; 337c4762a1bSJed Brown } 338c4762a1bSJed Brown } 339c4762a1bSJed Brown 340c4762a1bSJed Brown /* Restore data */ 3415f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabVecRestoreArrayRead(dm, Xloc, &x)); 3425f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabVecRestoreArrayRead(dm, Xdot, &xdot)); 3435f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabVecRestoreArray(dm, F, &f)); 3445f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dm, &Xloc)); 345c4762a1bSJed Brown PetscFunctionReturn(0); 346c4762a1bSJed Brown } 347c4762a1bSJed Brown 348c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx) 349c4762a1bSJed Brown { 350c4762a1bSJed Brown UserCtx user = (UserCtx)ctx; 351c4762a1bSJed Brown PetscReal vpos[3]; 352c4762a1bSJed Brown DM dm; 353c4762a1bSJed Brown Field *x; 354c4762a1bSJed Brown const moab::Range *vowned; 355c4762a1bSJed Brown PetscInt dof; 356c4762a1bSJed Brown moab::Range::iterator iter; 357c4762a1bSJed Brown 358c4762a1bSJed Brown PetscFunctionBegin; 3595f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &dm)); 360c4762a1bSJed Brown 361c4762a1bSJed Brown /* get the essential MOAB mesh related quantities needed for FEM assembly */ 3625f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabGetLocalVertices(dm, &vowned, NULL)); 363c4762a1bSJed Brown 3645f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(X, 0.0)); 365c4762a1bSJed Brown 366c4762a1bSJed Brown /* Get pointers to vector data */ 3675f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabVecGetArray(dm, X, &x)); 368c4762a1bSJed Brown 369c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 370c4762a1bSJed Brown for (moab::Range::iterator iter = vowned->begin(); iter != vowned->end(); iter++) { 371c4762a1bSJed Brown const moab::EntityHandle vhandle = *iter; 3725f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabGetDofsBlockedLocal(dm, 1, &vhandle, &dof)); 373c4762a1bSJed Brown 374c4762a1bSJed Brown /* compute the mid-point of the element and use a 1-point lumped quadrature */ 3755f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabGetVertexCoordinates(dm,1,&vhandle,vpos)); 376c4762a1bSJed Brown 377c4762a1bSJed Brown PetscReal xi = vpos[0]; 378c4762a1bSJed Brown x[dof].u = user->leftbc.u*(1.-xi) + user->rightbc.u*xi + PetscSinReal(2.*PETSC_PI*xi); 379c4762a1bSJed Brown x[dof].v = user->leftbc.v*(1.-xi) + user->rightbc.v*xi; 380c4762a1bSJed Brown } 381c4762a1bSJed Brown 382c4762a1bSJed Brown /* Restore vectors */ 3835f80ce2aSJacob Faibussowitsch CHKERRQ(DMMoabVecRestoreArray(dm, X, &x)); 384c4762a1bSJed Brown PetscFunctionReturn(0); 385c4762a1bSJed Brown } 386c4762a1bSJed Brown 387c4762a1bSJed Brown /*TEST 388c4762a1bSJed Brown 389c4762a1bSJed Brown build: 390c4762a1bSJed Brown requires: moab 391c4762a1bSJed Brown 392c4762a1bSJed Brown test: 393c4762a1bSJed Brown args: -n 20 -ts_type rosw -ts_rosw_type 2p -ts_dt 5e-2 -ts_adapt_type none 394c4762a1bSJed Brown 395c4762a1bSJed Brown test: 396c4762a1bSJed Brown suffix: 2 397c4762a1bSJed Brown nsize: 2 398c4762a1bSJed Brown args: -n 50 -ts_type glee -ts_adapt_type none -ts_dt 0.1 -io 399c4762a1bSJed Brown TODO: 400c4762a1bSJed Brown 401c4762a1bSJed Brown TEST*/ 402