static char help[] = "Solves DAE with integrator only on non-algebraic terms \n"; #include /* \dot{U} = f(U,V) F(U,V) = 0 Same as ex6.c except the user provided functions take input values as a single vector instead of two vectors */ /* f(U,V) = U + V */ PetscErrorCode f(PetscReal t,Vec UV,Vec F) { PetscErrorCode ierr; const PetscScalar *u,*v; PetscScalar *f; PetscInt n,i; PetscFunctionBeginUser; ierr = VecGetLocalSize(UV,&n);CHKERRQ(ierr); n = n/2; ierr = VecGetArrayRead(UV,&u);CHKERRQ(ierr); v = u + n; ierr = VecGetArrayWrite(F,&f);CHKERRQ(ierr); for (i=0; it = t; ierr = VecScatterBegin(ctx->scatterU,U,ctx->UV,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); ierr = VecScatterEnd(ctx->scatterU,U,ctx->UV,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); ierr = SNESSolve(ctx->snes,NULL,ctx->V);CHKERRQ(ierr); ierr = VecScatterBegin(ctx->scatterV,ctx->V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); ierr = VecScatterEnd(ctx->scatterV,ctx->V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); ierr = (*ctx->f)(t,ctx->UV,F);CHKERRQ(ierr); PetscFunctionReturn(0); } /* Defines the nonlinear function that is passed to the nonlinear solver */ PetscErrorCode SNESFunction(SNES snes,Vec V,Vec F,void *actx) { AppCtx *ctx = (AppCtx*)actx; PetscErrorCode ierr; PetscFunctionBeginUser; ierr = VecScatterBegin(ctx->scatterV,V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); ierr = VecScatterEnd(ctx->scatterV,V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); ierr = (*ctx->F)(ctx->t,ctx->UV,F);CHKERRQ(ierr); PetscFunctionReturn(0); } /*TEST test: args: -ts_monitor -ts_view TEST*/