1*900f6b5bSMatthew G. Knepley #include <petscconvest.h> /*I "petscconvest.h" I*/ 2*900f6b5bSMatthew G. Knepley #include <petscts.h> 3*900f6b5bSMatthew G. Knepley 4*900f6b5bSMatthew G. Knepley #include <petsc/private/petscconvestimpl.h> 5*900f6b5bSMatthew G. Knepley 6*900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstSetTS_Private(PetscConvEst ce, PetscObject solver) 7*900f6b5bSMatthew G. Knepley { 8*900f6b5bSMatthew G. Knepley PetscClassId id; 9*900f6b5bSMatthew G. Knepley PetscErrorCode ierr; 10*900f6b5bSMatthew G. Knepley 11*900f6b5bSMatthew G. Knepley PetscFunctionBegin; 12*900f6b5bSMatthew G. Knepley ierr = PetscObjectGetClassId(ce->solver, &id);CHKERRQ(ierr); 13*900f6b5bSMatthew G. Knepley if (id != TS_CLASSID) SETERRQ(PetscObjectComm((PetscObject) ce), PETSC_ERR_ARG_WRONG, "Solver was not a TS"); 14*900f6b5bSMatthew G. Knepley ierr = TSGetDM((TS) ce->solver, &ce->idm);CHKERRQ(ierr); 15*900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 16*900f6b5bSMatthew G. Knepley } 17*900f6b5bSMatthew G. Knepley 18*900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstInitGuessTS_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u) 19*900f6b5bSMatthew G. Knepley { 20*900f6b5bSMatthew G. Knepley PetscErrorCode ierr; 21*900f6b5bSMatthew G. Knepley 22*900f6b5bSMatthew G. Knepley PetscFunctionBegin; 23*900f6b5bSMatthew G. Knepley ierr = TSComputeInitialGuess((TS) ce->solver, u);CHKERRQ(ierr); 24*900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 25*900f6b5bSMatthew G. Knepley } 26*900f6b5bSMatthew G. Knepley 27*900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstComputeErrorTS_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[]) 28*900f6b5bSMatthew G. Knepley { 29*900f6b5bSMatthew G. Knepley Vec e; 30*900f6b5bSMatthew G. Knepley PetscErrorCode ierr; 31*900f6b5bSMatthew G. Knepley 32*900f6b5bSMatthew G. Knepley PetscFunctionBegin; 33*900f6b5bSMatthew G. Knepley ierr = VecDuplicate(u, &e);CHKERRQ(ierr); 34*900f6b5bSMatthew G. Knepley ierr = TSComputeExactError((TS) ce->solver, u, e);CHKERRQ(ierr); 35*900f6b5bSMatthew G. Knepley ierr = VecNorm(e, NORM_2, errors);CHKERRQ(ierr); 36*900f6b5bSMatthew G. Knepley ierr = VecDestroy(&e);CHKERRQ(ierr); 37*900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 38*900f6b5bSMatthew G. Knepley } 39*900f6b5bSMatthew G. Knepley 40*900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstGetConvRateTS_Private(PetscConvEst ce, PetscReal alpha[]) 41*900f6b5bSMatthew G. Knepley { 42*900f6b5bSMatthew G. Knepley TS ts = (TS) ce->solver; 43*900f6b5bSMatthew G. Knepley Vec u; 44*900f6b5bSMatthew G. Knepley PetscReal *dt, *x, *y, slope, intercept; 45*900f6b5bSMatthew G. Knepley PetscInt Ns, oNs, Nf = ce->Nf, f, Nr = ce->Nr, r; 46*900f6b5bSMatthew G. Knepley PetscErrorCode ierr; 47*900f6b5bSMatthew G. Knepley 48*900f6b5bSMatthew G. Knepley PetscFunctionBegin; 49*900f6b5bSMatthew G. Knepley ierr = TSGetSolution(ts, &u);CHKERRQ(ierr); 50*900f6b5bSMatthew G. Knepley ierr = PetscMalloc1(Nr+1, &dt);CHKERRQ(ierr); 51*900f6b5bSMatthew G. Knepley ierr = TSGetTimeStep(ts, &dt[0]);CHKERRQ(ierr); 52*900f6b5bSMatthew G. Knepley ierr = TSGetMaxSteps(ts, &oNs);CHKERRQ(ierr); 53*900f6b5bSMatthew G. Knepley Ns = oNs; 54*900f6b5bSMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 55*900f6b5bSMatthew G. Knepley if (r > 0) { 56*900f6b5bSMatthew G. Knepley dt[r] = dt[r-1]/2.0; 57*900f6b5bSMatthew G. Knepley Ns *= 2; 58*900f6b5bSMatthew G. Knepley } 59*900f6b5bSMatthew G. Knepley ierr = TSSetTime(ts, 0.0);CHKERRQ(ierr); 60*900f6b5bSMatthew G. Knepley ierr = TSSetStepNumber(ts, 0);CHKERRQ(ierr); 61*900f6b5bSMatthew G. Knepley ierr = TSSetTimeStep(ts, dt[r]);CHKERRQ(ierr); 62*900f6b5bSMatthew G. Knepley ierr = TSSetMaxSteps(ts, Ns);CHKERRQ(ierr); 63*900f6b5bSMatthew G. Knepley ierr = PetscConvEstComputeInitialGuess(ce, r, NULL, u);CHKERRQ(ierr); 64*900f6b5bSMatthew G. Knepley ierr = TSSolve(ts, u);CHKERRQ(ierr); 65*900f6b5bSMatthew G. Knepley ierr = PetscLogEventBegin(ce->event, ce, 0, 0, 0);CHKERRQ(ierr); 66*900f6b5bSMatthew G. Knepley ierr = PetscConvEstComputeError(ce, r, NULL, u, &ce->errors[r*Nf]);CHKERRQ(ierr); 67*900f6b5bSMatthew G. Knepley ierr = PetscLogEventEnd(ce->event, ce, 0, 0, 0);CHKERRQ(ierr); 68*900f6b5bSMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 69*900f6b5bSMatthew G. Knepley ierr = PetscLogEventSetDof(ce->event, f, 1.0/dt[r]);CHKERRQ(ierr); 70*900f6b5bSMatthew G. Knepley ierr = PetscLogEventSetError(ce->event, f, ce->errors[r*Nf+f]);CHKERRQ(ierr); 71*900f6b5bSMatthew G. Knepley } 72*900f6b5bSMatthew G. Knepley } 73*900f6b5bSMatthew G. Knepley /* Fit convergence rate */ 74*900f6b5bSMatthew G. Knepley if (Nr) { 75*900f6b5bSMatthew G. Knepley ierr = PetscMalloc2(Nr+1, &x, Nr+1, &y);CHKERRQ(ierr); 76*900f6b5bSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 77*900f6b5bSMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 78*900f6b5bSMatthew G. Knepley x[r] = PetscLog10Real(dt[r*Nf+f]); 79*900f6b5bSMatthew G. Knepley y[r] = PetscLog10Real(ce->errors[r*Nf+f]); 80*900f6b5bSMatthew G. Knepley } 81*900f6b5bSMatthew G. Knepley ierr = PetscLinearRegression(Nr+1, x, y, &slope, &intercept);CHKERRQ(ierr); 82*900f6b5bSMatthew G. Knepley /* Since lg err = s lg dt + b */ 83*900f6b5bSMatthew G. Knepley alpha[f] = slope; 84*900f6b5bSMatthew G. Knepley } 85*900f6b5bSMatthew G. Knepley ierr = PetscFree2(x, y);CHKERRQ(ierr); 86*900f6b5bSMatthew G. Knepley } 87*900f6b5bSMatthew G. Knepley /* Reset solver */ 88*900f6b5bSMatthew G. Knepley ierr = TSSetConvergedReason(ts, TS_CONVERGED_ITERATING);CHKERRQ(ierr); 89*900f6b5bSMatthew G. Knepley ierr = TSSetTime(ts, 0.0);CHKERRQ(ierr); 90*900f6b5bSMatthew G. Knepley ierr = TSSetStepNumber(ts, 0);CHKERRQ(ierr); 91*900f6b5bSMatthew G. Knepley ierr = TSSetTimeStep(ts, dt[0]);CHKERRQ(ierr); 92*900f6b5bSMatthew G. Knepley ierr = TSSetMaxSteps(ts, oNs);CHKERRQ(ierr); 93*900f6b5bSMatthew G. Knepley ierr = PetscConvEstComputeInitialGuess(ce, 0, NULL, u);CHKERRQ(ierr); 94*900f6b5bSMatthew G. Knepley ierr = PetscFree(dt);CHKERRQ(ierr); 95*900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 96*900f6b5bSMatthew G. Knepley } 97*900f6b5bSMatthew G. Knepley 98*900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstUseTS(PetscConvEst ce) 99*900f6b5bSMatthew G. Knepley { 100*900f6b5bSMatthew G. Knepley PetscFunctionBegin; 101*900f6b5bSMatthew G. Knepley PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 102*900f6b5bSMatthew G. Knepley ce->ops->setsolver = PetscConvEstSetTS_Private; 103*900f6b5bSMatthew G. Knepley ce->ops->initguess = PetscConvEstInitGuessTS_Private; 104*900f6b5bSMatthew G. Knepley ce->ops->computeerror = PetscConvEstComputeErrorTS_Private; 105*900f6b5bSMatthew G. Knepley ce->ops->getconvrate = PetscConvEstGetConvRateTS_Private; 106*900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 107*900f6b5bSMatthew G. Knepley } 108