1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Nonlinear DAE benchmark problems.\n"; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown /* 5*c4762a1bSJed Brown Include "petscts.h" so that we can use TS solvers. Note that this 6*c4762a1bSJed Brown file automatically includes: 7*c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 8*c4762a1bSJed Brown petscmat.h - matrices 9*c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 10*c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 11*c4762a1bSJed Brown petscksp.h - linear solvers 12*c4762a1bSJed Brown */ 13*c4762a1bSJed Brown #include <petscts.h> 14*c4762a1bSJed Brown 15*c4762a1bSJed Brown typedef struct _Problem* Problem; 16*c4762a1bSJed Brown struct _Problem { 17*c4762a1bSJed Brown PetscErrorCode (*destroy)(Problem); 18*c4762a1bSJed Brown TSIFunction function; 19*c4762a1bSJed Brown TSIJacobian jacobian; 20*c4762a1bSJed Brown PetscErrorCode (*solution)(PetscReal,Vec,void*); 21*c4762a1bSJed Brown MPI_Comm comm; 22*c4762a1bSJed Brown PetscReal final_time; 23*c4762a1bSJed Brown PetscInt n; 24*c4762a1bSJed Brown PetscBool hasexact; 25*c4762a1bSJed Brown void *data; 26*c4762a1bSJed Brown }; 27*c4762a1bSJed Brown 28*c4762a1bSJed Brown /* 29*c4762a1bSJed Brown Stiff 3-variable system from chemical reactions, due to Robertson (1966), problem ROBER in Hairer&Wanner, ODE 2, 1996 30*c4762a1bSJed Brown */ 31*c4762a1bSJed Brown static PetscErrorCode RoberFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx) 32*c4762a1bSJed Brown { 33*c4762a1bSJed Brown PetscErrorCode ierr; 34*c4762a1bSJed Brown PetscScalar *f; 35*c4762a1bSJed Brown const PetscScalar *x,*xdot; 36*c4762a1bSJed Brown 37*c4762a1bSJed Brown PetscFunctionBeginUser; 38*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 39*c4762a1bSJed Brown ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 41*c4762a1bSJed Brown f[0] = xdot[0] + 0.04*x[0] - 1e4*x[1]*x[2]; 42*c4762a1bSJed Brown f[1] = xdot[1] - 0.04*x[0] + 1e4*x[1]*x[2] + 3e7*PetscSqr(x[1]); 43*c4762a1bSJed Brown f[2] = xdot[2] - 3e7*PetscSqr(x[1]); 44*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 45*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr); 46*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 47*c4762a1bSJed Brown PetscFunctionReturn(0); 48*c4762a1bSJed Brown } 49*c4762a1bSJed Brown 50*c4762a1bSJed Brown static PetscErrorCode RoberJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx) 51*c4762a1bSJed Brown { 52*c4762a1bSJed Brown PetscErrorCode ierr; 53*c4762a1bSJed Brown PetscInt rowcol[] = {0,1,2}; 54*c4762a1bSJed Brown PetscScalar J[3][3]; 55*c4762a1bSJed Brown const PetscScalar *x,*xdot; 56*c4762a1bSJed Brown 57*c4762a1bSJed Brown PetscFunctionBeginUser; 58*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 59*c4762a1bSJed Brown ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr); 60*c4762a1bSJed Brown J[0][0] = a + 0.04; J[0][1] = -1e4*x[2]; J[0][2] = -1e4*x[1]; 61*c4762a1bSJed Brown J[1][0] = -0.04; J[1][1] = a + 1e4*x[2] + 3e7*2*x[1]; J[1][2] = 1e4*x[1]; 62*c4762a1bSJed Brown J[2][0] = 0; J[2][1] = -3e7*2*x[1]; J[2][2] = a; 63*c4762a1bSJed Brown ierr = MatSetValues(B,3,rowcol,3,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr); 66*c4762a1bSJed Brown 67*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 68*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 69*c4762a1bSJed Brown if (A != B) { 70*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 71*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 72*c4762a1bSJed Brown } 73*c4762a1bSJed Brown PetscFunctionReturn(0); 74*c4762a1bSJed Brown } 75*c4762a1bSJed Brown 76*c4762a1bSJed Brown static PetscErrorCode RoberSolution(PetscReal t,Vec X,void *ctx) 77*c4762a1bSJed Brown { 78*c4762a1bSJed Brown PetscErrorCode ierr; 79*c4762a1bSJed Brown PetscScalar *x; 80*c4762a1bSJed Brown 81*c4762a1bSJed Brown PetscFunctionBeginUser; 82*c4762a1bSJed Brown if (t != 0) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"not implemented"); 83*c4762a1bSJed Brown ierr = VecGetArray(X,&x);CHKERRQ(ierr); 84*c4762a1bSJed Brown x[0] = 1; 85*c4762a1bSJed Brown x[1] = 0; 86*c4762a1bSJed Brown x[2] = 0; 87*c4762a1bSJed Brown ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 88*c4762a1bSJed Brown PetscFunctionReturn(0); 89*c4762a1bSJed Brown } 90*c4762a1bSJed Brown 91*c4762a1bSJed Brown static PetscErrorCode RoberCreate(Problem p) 92*c4762a1bSJed Brown { 93*c4762a1bSJed Brown 94*c4762a1bSJed Brown PetscFunctionBeginUser; 95*c4762a1bSJed Brown p->destroy = 0; 96*c4762a1bSJed Brown p->function = &RoberFunction; 97*c4762a1bSJed Brown p->jacobian = &RoberJacobian; 98*c4762a1bSJed Brown p->solution = &RoberSolution; 99*c4762a1bSJed Brown p->final_time = 1e11; 100*c4762a1bSJed Brown p->n = 3; 101*c4762a1bSJed Brown PetscFunctionReturn(0); 102*c4762a1bSJed Brown } 103*c4762a1bSJed Brown 104*c4762a1bSJed Brown /* 105*c4762a1bSJed Brown Stiff scalar valued problem 106*c4762a1bSJed Brown */ 107*c4762a1bSJed Brown 108*c4762a1bSJed Brown typedef struct { 109*c4762a1bSJed Brown PetscReal lambda; 110*c4762a1bSJed Brown } CECtx; 111*c4762a1bSJed Brown 112*c4762a1bSJed Brown static PetscErrorCode CEDestroy(Problem p) 113*c4762a1bSJed Brown { 114*c4762a1bSJed Brown PetscErrorCode ierr; 115*c4762a1bSJed Brown 116*c4762a1bSJed Brown PetscFunctionBeginUser; 117*c4762a1bSJed Brown ierr = PetscFree(p->data);CHKERRQ(ierr); 118*c4762a1bSJed Brown PetscFunctionReturn(0); 119*c4762a1bSJed Brown } 120*c4762a1bSJed Brown 121*c4762a1bSJed Brown static PetscErrorCode CEFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx) 122*c4762a1bSJed Brown { 123*c4762a1bSJed Brown PetscErrorCode ierr; 124*c4762a1bSJed Brown PetscReal l = ((CECtx*)ctx)->lambda; 125*c4762a1bSJed Brown PetscScalar *f; 126*c4762a1bSJed Brown const PetscScalar *x,*xdot; 127*c4762a1bSJed Brown 128*c4762a1bSJed Brown PetscFunctionBeginUser; 129*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 130*c4762a1bSJed Brown ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr); 131*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 132*c4762a1bSJed Brown f[0] = xdot[0] + l*(x[0] - PetscCosReal(t)); 133*c4762a1bSJed Brown #if 0 134*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," f(t=%g,x=%g,xdot=%g) = %g\n",(double)t,(double)x[0],(double)xdot[0],(double)f[0]);CHKERRQ(ierr); 135*c4762a1bSJed Brown #endif 136*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 137*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr); 138*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 139*c4762a1bSJed Brown PetscFunctionReturn(0); 140*c4762a1bSJed Brown } 141*c4762a1bSJed Brown 142*c4762a1bSJed Brown static PetscErrorCode CEJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx) 143*c4762a1bSJed Brown { 144*c4762a1bSJed Brown PetscReal l = ((CECtx*)ctx)->lambda; 145*c4762a1bSJed Brown PetscErrorCode ierr; 146*c4762a1bSJed Brown PetscInt rowcol[] = {0}; 147*c4762a1bSJed Brown PetscScalar J[1][1]; 148*c4762a1bSJed Brown const PetscScalar *x,*xdot; 149*c4762a1bSJed Brown 150*c4762a1bSJed Brown PetscFunctionBeginUser; 151*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 152*c4762a1bSJed Brown ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr); 153*c4762a1bSJed Brown J[0][0] = a + l; 154*c4762a1bSJed Brown ierr = MatSetValues(B,1,rowcol,1,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 155*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 156*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr); 157*c4762a1bSJed Brown 158*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 159*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 160*c4762a1bSJed Brown if (A != B) { 161*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 162*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 163*c4762a1bSJed Brown } 164*c4762a1bSJed Brown PetscFunctionReturn(0); 165*c4762a1bSJed Brown } 166*c4762a1bSJed Brown 167*c4762a1bSJed Brown static PetscErrorCode CESolution(PetscReal t,Vec X,void *ctx) 168*c4762a1bSJed Brown { 169*c4762a1bSJed Brown PetscReal l = ((CECtx*)ctx)->lambda; 170*c4762a1bSJed Brown PetscErrorCode ierr; 171*c4762a1bSJed Brown PetscScalar *x; 172*c4762a1bSJed Brown 173*c4762a1bSJed Brown PetscFunctionBeginUser; 174*c4762a1bSJed Brown ierr = VecGetArray(X,&x);CHKERRQ(ierr); 175*c4762a1bSJed Brown x[0] = l/(l*l+1)*(l*PetscCosReal(t)+PetscSinReal(t)) - l*l/(l*l+1)*PetscExpReal(-l*t); 176*c4762a1bSJed Brown ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 177*c4762a1bSJed Brown PetscFunctionReturn(0); 178*c4762a1bSJed Brown } 179*c4762a1bSJed Brown 180*c4762a1bSJed Brown static PetscErrorCode CECreate(Problem p) 181*c4762a1bSJed Brown { 182*c4762a1bSJed Brown PetscErrorCode ierr; 183*c4762a1bSJed Brown CECtx *ce; 184*c4762a1bSJed Brown 185*c4762a1bSJed Brown PetscFunctionBeginUser; 186*c4762a1bSJed Brown ierr = PetscMalloc(sizeof(CECtx),&ce);CHKERRQ(ierr); 187*c4762a1bSJed Brown p->data = (void*)ce; 188*c4762a1bSJed Brown 189*c4762a1bSJed Brown p->destroy = &CEDestroy; 190*c4762a1bSJed Brown p->function = &CEFunction; 191*c4762a1bSJed Brown p->jacobian = &CEJacobian; 192*c4762a1bSJed Brown p->solution = &CESolution; 193*c4762a1bSJed Brown p->final_time = 10; 194*c4762a1bSJed Brown p->n = 1; 195*c4762a1bSJed Brown p->hasexact = PETSC_TRUE; 196*c4762a1bSJed Brown 197*c4762a1bSJed Brown ce->lambda = 10; 198*c4762a1bSJed Brown ierr = PetscOptionsBegin(p->comm,NULL,"CE options","");CHKERRQ(ierr); 199*c4762a1bSJed Brown { 200*c4762a1bSJed Brown ierr = PetscOptionsReal("-problem_ce_lambda","Parameter controlling stiffness: xdot + lambda*(x - cos(t))","",ce->lambda,&ce->lambda,NULL);CHKERRQ(ierr); 201*c4762a1bSJed Brown } 202*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 203*c4762a1bSJed Brown PetscFunctionReturn(0); 204*c4762a1bSJed Brown } 205*c4762a1bSJed Brown 206*c4762a1bSJed Brown /* 207*c4762a1bSJed Brown * Stiff 3-variable oscillatory system from chemical reactions. problem OREGO in Hairer&Wanner 208*c4762a1bSJed Brown */ 209*c4762a1bSJed Brown static PetscErrorCode OregoFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx) 210*c4762a1bSJed Brown { 211*c4762a1bSJed Brown PetscErrorCode ierr; 212*c4762a1bSJed Brown PetscScalar *f; 213*c4762a1bSJed Brown const PetscScalar *x,*xdot; 214*c4762a1bSJed Brown 215*c4762a1bSJed Brown PetscFunctionBeginUser; 216*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 217*c4762a1bSJed Brown ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr); 218*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 219*c4762a1bSJed Brown f[0] = xdot[0] - 77.27*(x[1] + x[0]*(1. - 8.375e-6*x[0] - x[1])); 220*c4762a1bSJed Brown f[1] = xdot[1] - 1/77.27*(x[2] - (1. + x[0])*x[1]); 221*c4762a1bSJed Brown f[2] = xdot[2] - 0.161*(x[0] - x[2]); 222*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 223*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr); 224*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 225*c4762a1bSJed Brown PetscFunctionReturn(0); 226*c4762a1bSJed Brown } 227*c4762a1bSJed Brown 228*c4762a1bSJed Brown static PetscErrorCode OregoJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx) 229*c4762a1bSJed Brown { 230*c4762a1bSJed Brown PetscErrorCode ierr; 231*c4762a1bSJed Brown PetscInt rowcol[] = {0,1,2}; 232*c4762a1bSJed Brown PetscScalar J[3][3]; 233*c4762a1bSJed Brown const PetscScalar *x,*xdot; 234*c4762a1bSJed Brown 235*c4762a1bSJed Brown PetscFunctionBeginUser; 236*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 237*c4762a1bSJed Brown ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr); 238*c4762a1bSJed Brown J[0][0] = a - 77.27*((1. - 8.375e-6*x[0] - x[1]) - 8.375e-6*x[0]); 239*c4762a1bSJed Brown J[0][1] = -77.27*(1. - x[0]); 240*c4762a1bSJed Brown J[0][2] = 0; 241*c4762a1bSJed Brown J[1][0] = 1./77.27*x[1]; 242*c4762a1bSJed Brown J[1][1] = a + 1./77.27*(1. + x[0]); 243*c4762a1bSJed Brown J[1][2] = -1./77.27; 244*c4762a1bSJed Brown J[2][0] = -0.161; 245*c4762a1bSJed Brown J[2][1] = 0; 246*c4762a1bSJed Brown J[2][2] = a + 0.161; 247*c4762a1bSJed Brown ierr = MatSetValues(B,3,rowcol,3,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 248*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 249*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr); 250*c4762a1bSJed Brown 251*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 252*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 253*c4762a1bSJed Brown if (A != B) { 254*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 255*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 256*c4762a1bSJed Brown } 257*c4762a1bSJed Brown PetscFunctionReturn(0); 258*c4762a1bSJed Brown } 259*c4762a1bSJed Brown 260*c4762a1bSJed Brown static PetscErrorCode OregoSolution(PetscReal t,Vec X,void *ctx) 261*c4762a1bSJed Brown { 262*c4762a1bSJed Brown PetscErrorCode ierr; 263*c4762a1bSJed Brown PetscScalar *x; 264*c4762a1bSJed Brown 265*c4762a1bSJed Brown PetscFunctionBeginUser; 266*c4762a1bSJed Brown if (t != 0) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"not implemented"); 267*c4762a1bSJed Brown ierr = VecGetArray(X,&x);CHKERRQ(ierr); 268*c4762a1bSJed Brown x[0] = 1; 269*c4762a1bSJed Brown x[1] = 2; 270*c4762a1bSJed Brown x[2] = 3; 271*c4762a1bSJed Brown ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 272*c4762a1bSJed Brown PetscFunctionReturn(0); 273*c4762a1bSJed Brown } 274*c4762a1bSJed Brown 275*c4762a1bSJed Brown static PetscErrorCode OregoCreate(Problem p) 276*c4762a1bSJed Brown { 277*c4762a1bSJed Brown 278*c4762a1bSJed Brown PetscFunctionBeginUser; 279*c4762a1bSJed Brown p->destroy = 0; 280*c4762a1bSJed Brown p->function = &OregoFunction; 281*c4762a1bSJed Brown p->jacobian = &OregoJacobian; 282*c4762a1bSJed Brown p->solution = &OregoSolution; 283*c4762a1bSJed Brown p->final_time = 360; 284*c4762a1bSJed Brown p->n = 3; 285*c4762a1bSJed Brown PetscFunctionReturn(0); 286*c4762a1bSJed Brown } 287*c4762a1bSJed Brown 288*c4762a1bSJed Brown 289*c4762a1bSJed Brown /* 290*c4762a1bSJed Brown * User-defined monitor for comparing to exact solutions when possible 291*c4762a1bSJed Brown */ 292*c4762a1bSJed Brown typedef struct { 293*c4762a1bSJed Brown MPI_Comm comm; 294*c4762a1bSJed Brown Problem problem; 295*c4762a1bSJed Brown Vec x; 296*c4762a1bSJed Brown } MonitorCtx; 297*c4762a1bSJed Brown 298*c4762a1bSJed Brown static PetscErrorCode MonitorError(TS ts,PetscInt step,PetscReal t,Vec x,void *ctx) 299*c4762a1bSJed Brown { 300*c4762a1bSJed Brown PetscErrorCode ierr; 301*c4762a1bSJed Brown MonitorCtx *mon = (MonitorCtx*)ctx; 302*c4762a1bSJed Brown PetscReal h,nrm_x,nrm_exact,nrm_diff; 303*c4762a1bSJed Brown 304*c4762a1bSJed Brown PetscFunctionBeginUser; 305*c4762a1bSJed Brown if (!mon->problem->solution) PetscFunctionReturn(0); 306*c4762a1bSJed Brown ierr = (*mon->problem->solution)(t,mon->x,mon->problem->data);CHKERRQ(ierr); 307*c4762a1bSJed Brown ierr = VecNorm(x,NORM_2,&nrm_x);CHKERRQ(ierr); 308*c4762a1bSJed Brown ierr = VecNorm(mon->x,NORM_2,&nrm_exact);CHKERRQ(ierr); 309*c4762a1bSJed Brown ierr = VecAYPX(mon->x,-1,x);CHKERRQ(ierr); 310*c4762a1bSJed Brown ierr = VecNorm(mon->x,NORM_2,&nrm_diff);CHKERRQ(ierr); 311*c4762a1bSJed Brown ierr = TSGetTimeStep(ts,&h);CHKERRQ(ierr); 312*c4762a1bSJed Brown if (step < 0) { 313*c4762a1bSJed Brown ierr = PetscPrintf(mon->comm,"Interpolated final solution ");CHKERRQ(ierr); 314*c4762a1bSJed Brown } 315*c4762a1bSJed Brown ierr = PetscPrintf(mon->comm,"step %4D t=%12.8e h=% 8.2e |x|=%9.2e |x_e|=%9.2e |x-x_e|=%9.2e\n",step,(double)t,(double)h,(double)nrm_x,(double)nrm_exact,(double)nrm_diff);CHKERRQ(ierr); 316*c4762a1bSJed Brown PetscFunctionReturn(0); 317*c4762a1bSJed Brown } 318*c4762a1bSJed Brown 319*c4762a1bSJed Brown 320*c4762a1bSJed Brown int main(int argc,char **argv) 321*c4762a1bSJed Brown { 322*c4762a1bSJed Brown PetscFunctionList plist = NULL; 323*c4762a1bSJed Brown char pname[256]; 324*c4762a1bSJed Brown TS ts; /* nonlinear solver */ 325*c4762a1bSJed Brown Vec x,r; /* solution, residual vectors */ 326*c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 327*c4762a1bSJed Brown Problem problem; 328*c4762a1bSJed Brown PetscBool use_monitor; 329*c4762a1bSJed Brown PetscInt steps,nonlinits,linits,snesfails,rejects; 330*c4762a1bSJed Brown PetscReal ftime; 331*c4762a1bSJed Brown MonitorCtx mon; 332*c4762a1bSJed Brown PetscErrorCode ierr; 333*c4762a1bSJed Brown PetscMPIInt size; 334*c4762a1bSJed Brown 335*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 336*c4762a1bSJed Brown Initialize program 337*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 338*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 339*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 340*c4762a1bSJed Brown if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs"); 341*c4762a1bSJed Brown 342*c4762a1bSJed Brown /* Register the available problems */ 343*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&plist,"rober",&RoberCreate);CHKERRQ(ierr); 344*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&plist,"ce",&CECreate);CHKERRQ(ierr); 345*c4762a1bSJed Brown ierr = PetscFunctionListAdd(&plist,"orego",&OregoCreate);CHKERRQ(ierr); 346*c4762a1bSJed Brown ierr = PetscStrcpy(pname,"ce");CHKERRQ(ierr); 347*c4762a1bSJed Brown 348*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 349*c4762a1bSJed Brown Set runtime options 350*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 351*c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Timestepping benchmark options","");CHKERRQ(ierr); 352*c4762a1bSJed Brown { 353*c4762a1bSJed Brown ierr = PetscOptionsFList("-problem_type","Name of problem to run","",plist,pname,pname,sizeof(pname),NULL);CHKERRQ(ierr); 354*c4762a1bSJed Brown use_monitor = PETSC_FALSE; 355*c4762a1bSJed Brown ierr = PetscOptionsBool("-monitor_error","Display errors relative to exact solutions","",use_monitor,&use_monitor,NULL);CHKERRQ(ierr); 356*c4762a1bSJed Brown } 357*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 358*c4762a1bSJed Brown 359*c4762a1bSJed Brown /* Create the new problem */ 360*c4762a1bSJed Brown ierr = PetscNew(&problem);CHKERRQ(ierr); 361*c4762a1bSJed Brown problem->comm = MPI_COMM_WORLD; 362*c4762a1bSJed Brown { 363*c4762a1bSJed Brown PetscErrorCode (*pcreate)(Problem); 364*c4762a1bSJed Brown 365*c4762a1bSJed Brown ierr = PetscFunctionListFind(plist,pname,&pcreate);CHKERRQ(ierr); 366*c4762a1bSJed Brown if (!pcreate) SETERRQ1(PETSC_COMM_SELF,1,"No problem '%s'",pname); 367*c4762a1bSJed Brown ierr = (*pcreate)(problem);CHKERRQ(ierr); 368*c4762a1bSJed Brown } 369*c4762a1bSJed Brown 370*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 371*c4762a1bSJed Brown Create necessary matrix and vectors 372*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 373*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 374*c4762a1bSJed Brown ierr = MatSetSizes(A,problem->n,problem->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 375*c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 376*c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 377*c4762a1bSJed Brown 378*c4762a1bSJed Brown ierr = MatCreateVecs(A,&x,NULL);CHKERRQ(ierr); 379*c4762a1bSJed Brown ierr = VecDuplicate(x,&r);CHKERRQ(ierr); 380*c4762a1bSJed Brown 381*c4762a1bSJed Brown mon.comm = PETSC_COMM_WORLD; 382*c4762a1bSJed Brown mon.problem = problem; 383*c4762a1bSJed Brown ierr = VecDuplicate(x,&mon.x);CHKERRQ(ierr); 384*c4762a1bSJed Brown 385*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 386*c4762a1bSJed Brown Create timestepping solver context 387*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 388*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 389*c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 390*c4762a1bSJed Brown ierr = TSSetType(ts,TSROSW);CHKERRQ(ierr); /* Rosenbrock-W */ 391*c4762a1bSJed Brown ierr = TSSetIFunction(ts,NULL,problem->function,problem->data);CHKERRQ(ierr); 392*c4762a1bSJed Brown ierr = TSSetIJacobian(ts,A,A,problem->jacobian,problem->data);CHKERRQ(ierr); 393*c4762a1bSJed Brown ierr = TSSetMaxTime(ts,problem->final_time);CHKERRQ(ierr); 394*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 395*c4762a1bSJed Brown ierr = TSSetMaxStepRejections(ts,10);CHKERRQ(ierr); 396*c4762a1bSJed Brown ierr = TSSetMaxSNESFailures(ts,-1);CHKERRQ(ierr); /* unlimited */ 397*c4762a1bSJed Brown if (use_monitor) { 398*c4762a1bSJed Brown ierr = TSMonitorSet(ts,&MonitorError,&mon,NULL);CHKERRQ(ierr); 399*c4762a1bSJed Brown } 400*c4762a1bSJed Brown 401*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 402*c4762a1bSJed Brown Set initial conditions 403*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 404*c4762a1bSJed Brown ierr = (*problem->solution)(0,x,problem->data);CHKERRQ(ierr); 405*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,.001);CHKERRQ(ierr); 406*c4762a1bSJed Brown ierr = TSSetSolution(ts,x);CHKERRQ(ierr); 407*c4762a1bSJed Brown 408*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 409*c4762a1bSJed Brown Set runtime options 410*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 411*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 412*c4762a1bSJed Brown 413*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 414*c4762a1bSJed Brown Solve nonlinear system 415*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 416*c4762a1bSJed Brown ierr = TSSolve(ts,x);CHKERRQ(ierr); 417*c4762a1bSJed Brown ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); 418*c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 419*c4762a1bSJed Brown ierr = TSGetSNESFailures(ts,&snesfails);CHKERRQ(ierr); 420*c4762a1bSJed Brown ierr = TSGetStepRejections(ts,&rejects);CHKERRQ(ierr); 421*c4762a1bSJed Brown ierr = TSGetSNESIterations(ts,&nonlinits);CHKERRQ(ierr); 422*c4762a1bSJed Brown ierr = TSGetKSPIterations(ts,&linits);CHKERRQ(ierr); 423*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"steps %D (%D rejected, %D SNES fails), ftime %g, nonlinits %D, linits %D\n",steps,rejects,snesfails,(double)ftime,nonlinits,linits);CHKERRQ(ierr); 424*c4762a1bSJed Brown 425*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 426*c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 427*c4762a1bSJed Brown are no longer needed. 428*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 429*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 430*c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 431*c4762a1bSJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 432*c4762a1bSJed Brown ierr = VecDestroy(&mon.x);CHKERRQ(ierr); 433*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 434*c4762a1bSJed Brown if (problem->destroy) { 435*c4762a1bSJed Brown ierr = (*problem->destroy)(problem);CHKERRQ(ierr); 436*c4762a1bSJed Brown } 437*c4762a1bSJed Brown ierr = PetscFree(problem);CHKERRQ(ierr); 438*c4762a1bSJed Brown ierr = PetscFunctionListDestroy(&plist);CHKERRQ(ierr); 439*c4762a1bSJed Brown 440*c4762a1bSJed Brown ierr = PetscFinalize(); 441*c4762a1bSJed Brown return ierr; 442*c4762a1bSJed Brown } 443*c4762a1bSJed Brown 444*c4762a1bSJed Brown /*TEST 445*c4762a1bSJed Brown 446*c4762a1bSJed Brown test: 447*c4762a1bSJed Brown requires: !complex 448*c4762a1bSJed Brown args: -monitor_error -ts_atol 1e-2 -ts_rtol 1e-2 -ts_exact_final_time interpolate -ts_type arkimex 449*c4762a1bSJed Brown 450*c4762a1bSJed Brown test: 451*c4762a1bSJed Brown suffix: 2 452*c4762a1bSJed Brown requires: !single !complex 453*c4762a1bSJed Brown args: -ts_atol 1e-2 -ts_rtol 1e-2 -ts_max_time 15 -ts_type arkimex -ts_arkimex_type 2e -problem_type orego -ts_arkimex_initial_guess_extrapolate 0 -ts_adapt_time_step_increase_delay 4 454*c4762a1bSJed Brown 455*c4762a1bSJed Brown test: 456*c4762a1bSJed Brown suffix: 3 457*c4762a1bSJed Brown requires: !single !complex 458*c4762a1bSJed Brown args: -ts_atol 1e-2 -ts_rtol 1e-2 -ts_max_time 15 -ts_type arkimex -ts_arkimex_type 2e -problem_type orego -ts_arkimex_initial_guess_extrapolate 1 459*c4762a1bSJed Brown 460*c4762a1bSJed Brown TEST*/ 461