1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Solves biharmonic equation in 1d.\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /* 5c4762a1bSJed Brown Solves the equation biharmonic equation in split form 6c4762a1bSJed Brown 7c4762a1bSJed Brown w = -kappa \Delta u 8c4762a1bSJed Brown u_t = \Delta w 9c4762a1bSJed Brown -1 <= u <= 1 10c4762a1bSJed Brown Periodic boundary conditions 11c4762a1bSJed Brown 12c4762a1bSJed Brown Evolve the biharmonic heat equation with bounds: (same as biharmonic) 13c4762a1bSJed Brown --------------- 14c4762a1bSJed Brown ./biharmonic2 -ts_monitor -snes_monitor -ts_monitor_draw_solution -pc_type lu -draw_pause .1 -snes_converged_reason -ts_type beuler -da_refine 5 -draw_fields 1 -ts_dt 9.53674e-9 15c4762a1bSJed Brown 16c4762a1bSJed Brown w = -kappa \Delta u + u^3 - u 17c4762a1bSJed Brown u_t = \Delta w 18c4762a1bSJed Brown -1 <= u <= 1 19c4762a1bSJed Brown Periodic boundary conditions 20c4762a1bSJed Brown 21c4762a1bSJed Brown Evolve the Cahn-Hillard equations: (this fails after a few timesteps 12/17/2017) 22c4762a1bSJed Brown --------------- 23c4762a1bSJed Brown ./biharmonic2 -ts_monitor -snes_monitor -ts_monitor_draw_solution -pc_type lu -draw_pause .1 -snes_converged_reason -ts_type beuler -da_refine 6 -draw_fields 1 -kappa .00001 -ts_dt 5.96046e-06 -cahn-hillard 24c4762a1bSJed Brown 25c4762a1bSJed Brown */ 26c4762a1bSJed Brown #include <petscdm.h> 27c4762a1bSJed Brown #include <petscdmda.h> 28c4762a1bSJed Brown #include <petscts.h> 29c4762a1bSJed Brown #include <petscdraw.h> 30c4762a1bSJed Brown 31c4762a1bSJed Brown /* 32c4762a1bSJed Brown User-defined routines 33c4762a1bSJed Brown */ 34c4762a1bSJed Brown extern PetscErrorCode FormFunction(TS,PetscReal,Vec,Vec,Vec,void*),FormInitialSolution(DM,Vec,PetscReal); 35c4762a1bSJed Brown typedef struct {PetscBool cahnhillard;PetscReal kappa;PetscInt energy;PetscReal tol;PetscReal theta;PetscReal theta_c;} UserCtx; 36c4762a1bSJed Brown 37c4762a1bSJed Brown int main(int argc,char **argv) 38c4762a1bSJed Brown { 39c4762a1bSJed Brown TS ts; /* nonlinear solver */ 40c4762a1bSJed Brown Vec x,r; /* solution, residual vectors */ 41c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 42c4762a1bSJed Brown PetscInt steps,Mx; 43c4762a1bSJed Brown DM da; 44c4762a1bSJed Brown MatFDColoring matfdcoloring; 45c4762a1bSJed Brown ISColoring iscoloring; 46c4762a1bSJed Brown PetscReal dt; 47c4762a1bSJed Brown PetscReal vbounds[] = {-100000,100000,-1.1,1.1}; 48c4762a1bSJed Brown SNES snes; 49c4762a1bSJed Brown UserCtx ctx; 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 52c4762a1bSJed Brown Initialize program 53c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 54*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 55c4762a1bSJed Brown ctx.kappa = 1.0; 56*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-kappa",&ctx.kappa,NULL)); 57c4762a1bSJed Brown ctx.cahnhillard = PETSC_FALSE; 58c4762a1bSJed Brown 59*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-cahn-hillard",&ctx.cahnhillard,NULL)); 60*9566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),2,vbounds)); 61*9566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),600,600)); 62c4762a1bSJed Brown ctx.energy = 1; 63*9566063dSJacob Faibussowitsch /*PetscCall(PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL));*/ 64*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL)); 65c4762a1bSJed Brown ctx.tol = 1.0e-8; 66*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-tol",&ctx.tol,NULL)); 67c4762a1bSJed Brown ctx.theta = .001; 68c4762a1bSJed Brown ctx.theta_c = 1.0; 69*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-theta",&ctx.theta,NULL)); 70*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-theta_c",&ctx.theta_c,NULL)); 71c4762a1bSJed Brown 72c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 73c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 74c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 75*9566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10,2,2,NULL,&da)); 76*9566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 77*9566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 78*9566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da,0,"Biharmonic heat equation: w = -kappa*u_xx")); 79*9566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da,1,"Biharmonic heat equation: u")); 80*9566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0,&Mx,0,0,0,0,0,0,0,0,0,0,0)); 81c4762a1bSJed Brown dt = 1.0/(10.*ctx.kappa*Mx*Mx*Mx*Mx); 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 84c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining 85c4762a1bSJed Brown vectors that are the same types 86c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 87*9566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da,&x)); 88*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&r)); 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 91c4762a1bSJed Brown Create timestepping solver context 92c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 93*9566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 94*9566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts,da)); 95*9566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts,TS_NONLINEAR)); 96*9566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts,NULL,FormFunction,&ctx)); 97*9566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,.02)); 98*9566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_INTERPOLATE)); 99c4762a1bSJed Brown 100c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 101c4762a1bSJed Brown Create matrix data structure; set Jacobian evaluation routine 102c4762a1bSJed Brown 103c4762a1bSJed Brown < Set Jacobian matrix data structure and default Jacobian evaluation 104c4762a1bSJed Brown routine. User can override with: 105c4762a1bSJed Brown -snes_mf : matrix-free Newton-Krylov method with no preconditioning 106c4762a1bSJed Brown (unless user explicitly sets preconditioner) 107c4762a1bSJed Brown -snes_mf_operator : form preconditioning matrix as set by the user, 108c4762a1bSJed Brown but use matrix-free approx for Jacobian-vector 109c4762a1bSJed Brown products within Newton-Krylov method 110c4762a1bSJed Brown 111c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 112*9566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts,&snes)); 113*9566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(da,IS_COLORING_GLOBAL,&iscoloring)); 114*9566063dSJacob Faibussowitsch PetscCall(DMSetMatType(da,MATAIJ)); 115*9566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(da,&J)); 116*9566063dSJacob Faibussowitsch PetscCall(MatFDColoringCreate(J,iscoloring,&matfdcoloring)); 117*9566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts)); 118*9566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFromOptions(matfdcoloring)); 119*9566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetUp(J,iscoloring,matfdcoloring)); 120*9566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&iscoloring)); 121*9566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring)); 122c4762a1bSJed Brown 123c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 124c4762a1bSJed Brown Customize nonlinear solver 125c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 126*9566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSBEULER)); 127c4762a1bSJed Brown 128c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 129c4762a1bSJed Brown Set initial conditions 130c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 131*9566063dSJacob Faibussowitsch PetscCall(FormInitialSolution(da,x,ctx.kappa)); 132*9566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,dt)); 133*9566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts,x)); 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 136c4762a1bSJed Brown Set runtime options 137c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 138*9566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 139c4762a1bSJed Brown 140c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 141c4762a1bSJed Brown Solve nonlinear system 142c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 143*9566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,x)); 144*9566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&steps)); 145c4762a1bSJed Brown 146c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 147c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 148c4762a1bSJed Brown are no longer needed. 149c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 150*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 151*9566063dSJacob Faibussowitsch PetscCall(MatFDColoringDestroy(&matfdcoloring)); 152*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 153*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 154*9566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 155*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 156c4762a1bSJed Brown 157*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 158b122ec5aSJacob Faibussowitsch return 0; 159c4762a1bSJed Brown } 160c4762a1bSJed Brown 161c4762a1bSJed Brown typedef struct {PetscScalar w,u;} Field; 162c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 163c4762a1bSJed Brown /* 164c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x). 165c4762a1bSJed Brown 166c4762a1bSJed Brown Input Parameters: 167c4762a1bSJed Brown . ts - the TS context 168c4762a1bSJed Brown . X - input vector 169c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetFunction() 170c4762a1bSJed Brown 171c4762a1bSJed Brown Output Parameter: 172c4762a1bSJed Brown . F - function vector 173c4762a1bSJed Brown */ 174c4762a1bSJed Brown PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec Xdot,Vec F,void *ptr) 175c4762a1bSJed Brown { 176c4762a1bSJed Brown DM da; 177c4762a1bSJed Brown PetscInt i,Mx,xs,xm; 178c4762a1bSJed Brown PetscReal hx,sx; 179c4762a1bSJed Brown Field *x,*xdot,*f; 180c4762a1bSJed Brown Vec localX,localXdot; 181c4762a1bSJed Brown UserCtx *ctx = (UserCtx*)ptr; 182c4762a1bSJed Brown 183c4762a1bSJed Brown PetscFunctionBegin; 184*9566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 185*9566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&localX)); 186*9566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&localXdot)); 187*9566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE)); 188c4762a1bSJed Brown 189c4762a1bSJed Brown hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx); 190c4762a1bSJed Brown 191c4762a1bSJed Brown /* 192c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 193c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 194c4762a1bSJed Brown By placing code between these two statements, computations can be 195c4762a1bSJed Brown done while messages are in transition. 196c4762a1bSJed Brown */ 197*9566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX)); 198*9566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX)); 199*9566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,Xdot,INSERT_VALUES,localXdot)); 200*9566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,Xdot,INSERT_VALUES,localXdot)); 201c4762a1bSJed Brown 202c4762a1bSJed Brown /* 203c4762a1bSJed Brown Get pointers to vector data 204c4762a1bSJed Brown */ 205*9566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da,localX,&x)); 206*9566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da,localXdot,&xdot)); 207*9566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,F,&f)); 208c4762a1bSJed Brown 209c4762a1bSJed Brown /* 210c4762a1bSJed Brown Get local grid boundaries 211c4762a1bSJed Brown */ 212*9566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 213c4762a1bSJed Brown 214c4762a1bSJed Brown /* 215c4762a1bSJed Brown Compute function over the locally owned part of the grid 216c4762a1bSJed Brown */ 217c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 218c4762a1bSJed Brown f[i].w = x[i].w + ctx->kappa*(x[i-1].u + x[i+1].u - 2.0*x[i].u)*sx; 219c4762a1bSJed Brown if (ctx->cahnhillard) { 220c4762a1bSJed Brown switch (ctx->energy) { 221c4762a1bSJed Brown case 1: /* double well */ 222c4762a1bSJed Brown f[i].w += -x[i].u*x[i].u*x[i].u + x[i].u; 223c4762a1bSJed Brown break; 224c4762a1bSJed Brown case 2: /* double obstacle */ 225c4762a1bSJed Brown f[i].w += x[i].u; 226c4762a1bSJed Brown break; 227c4762a1bSJed Brown case 3: /* logarithmic */ 228c4762a1bSJed Brown if (PetscRealPart(x[i].u) < -1.0 + 2.0*ctx->tol) f[i].w += .5*ctx->theta*(-PetscLogReal(ctx->tol) + PetscLogScalar((1.0-x[i].u)/2.0)) + ctx->theta_c*x[i].u; 229c4762a1bSJed Brown else if (PetscRealPart(x[i].u) > 1.0 - 2.0*ctx->tol) f[i].w += .5*ctx->theta*(-PetscLogScalar((1.0+x[i].u)/2.0) + PetscLogReal(ctx->tol)) + ctx->theta_c*x[i].u; 230c4762a1bSJed Brown else f[i].w += .5*ctx->theta*(-PetscLogScalar((1.0+x[i].u)/2.0) + PetscLogScalar((1.0-x[i].u)/2.0)) + ctx->theta_c*x[i].u; 231c4762a1bSJed Brown break; 232c4762a1bSJed Brown } 233c4762a1bSJed Brown } 234c4762a1bSJed Brown f[i].u = xdot[i].u - (x[i-1].w + x[i+1].w - 2.0*x[i].w)*sx; 235c4762a1bSJed Brown } 236c4762a1bSJed Brown 237c4762a1bSJed Brown /* 238c4762a1bSJed Brown Restore vectors 239c4762a1bSJed Brown */ 240*9566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da,localXdot,&xdot)); 241*9566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da,localX,&x)); 242*9566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,F,&f)); 243*9566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&localX)); 244*9566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&localXdot)); 245c4762a1bSJed Brown PetscFunctionReturn(0); 246c4762a1bSJed Brown } 247c4762a1bSJed Brown 248c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 249c4762a1bSJed Brown PetscErrorCode FormInitialSolution(DM da,Vec X,PetscReal kappa) 250c4762a1bSJed Brown { 251c4762a1bSJed Brown PetscInt i,xs,xm,Mx,xgs,xgm; 252c4762a1bSJed Brown Field *x; 253c4762a1bSJed Brown PetscReal hx,xx,r,sx; 254c4762a1bSJed Brown Vec Xg; 255c4762a1bSJed Brown 256c4762a1bSJed Brown PetscFunctionBegin; 257*9566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE)); 258c4762a1bSJed Brown 259c4762a1bSJed Brown hx = 1.0/(PetscReal)Mx; 260c4762a1bSJed Brown sx = 1.0/(hx*hx); 261c4762a1bSJed Brown 262c4762a1bSJed Brown /* 263c4762a1bSJed Brown Get pointers to vector data 264c4762a1bSJed Brown */ 265*9566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(da,&Xg)); 266*9566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,Xg,&x)); 267c4762a1bSJed Brown 268c4762a1bSJed Brown /* 269c4762a1bSJed Brown Get local grid boundaries 270c4762a1bSJed Brown */ 271*9566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 272*9566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da,&xgs,NULL,NULL,&xgm,NULL,NULL)); 273c4762a1bSJed Brown 274c4762a1bSJed Brown /* 275c4762a1bSJed Brown Compute u function over the locally owned part of the grid including ghost points 276c4762a1bSJed Brown */ 277c4762a1bSJed Brown for (i=xgs; i<xgs+xgm; i++) { 278c4762a1bSJed Brown xx = i*hx; 279c4762a1bSJed Brown r = PetscSqrtReal((xx-.5)*(xx-.5)); 280c4762a1bSJed Brown if (r < .125) x[i].u = 1.0; 281c4762a1bSJed Brown else x[i].u = -.50; 282c4762a1bSJed Brown /* fill in x[i].w so that valgrind doesn't detect use of uninitialized memory */ 283c4762a1bSJed Brown x[i].w = 0; 284c4762a1bSJed Brown } 285c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) x[i].w = -kappa*(x[i-1].u + x[i+1].u - 2.0*x[i].u)*sx; 286c4762a1bSJed Brown 287c4762a1bSJed Brown /* 288c4762a1bSJed Brown Restore vectors 289c4762a1bSJed Brown */ 290*9566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Xg,&x)); 291c4762a1bSJed Brown 292c4762a1bSJed Brown /* Grab only the global part of the vector */ 293*9566063dSJacob Faibussowitsch PetscCall(VecSet(X,0)); 294*9566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(da,Xg,ADD_VALUES,X)); 295*9566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(da,Xg,ADD_VALUES,X)); 296*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Xg)); 297c4762a1bSJed Brown PetscFunctionReturn(0); 298c4762a1bSJed Brown } 299c4762a1bSJed Brown 300c4762a1bSJed Brown /*TEST 301c4762a1bSJed Brown 302c4762a1bSJed Brown build: 303c4762a1bSJed Brown requires: !complex !single 304c4762a1bSJed Brown 305c4762a1bSJed Brown test: 306c4762a1bSJed Brown args: -ts_monitor -snes_monitor -pc_type lu -snes_converged_reason -ts_type beuler -da_refine 5 -ts_dt 9.53674e-9 -ts_max_steps 50 307c4762a1bSJed Brown requires: x 308c4762a1bSJed Brown 309c4762a1bSJed Brown TEST*/ 310