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