1900f6b5bSMatthew G. Knepley #include <petscconvest.h> /*I "petscconvest.h" I*/ 2900f6b5bSMatthew G. Knepley #include <petscts.h> 3900f6b5bSMatthew G. Knepley 4900f6b5bSMatthew G. Knepley #include <petsc/private/petscconvestimpl.h> 5900f6b5bSMatthew G. Knepley 6900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstSetTS_Private(PetscConvEst ce, PetscObject solver) 7900f6b5bSMatthew G. Knepley { 8900f6b5bSMatthew G. Knepley PetscClassId id; 9900f6b5bSMatthew G. Knepley PetscErrorCode ierr; 10900f6b5bSMatthew G. Knepley 11900f6b5bSMatthew G. Knepley PetscFunctionBegin; 12900f6b5bSMatthew G. Knepley ierr = PetscObjectGetClassId(ce->solver, &id);CHKERRQ(ierr); 13900f6b5bSMatthew G. Knepley if (id != TS_CLASSID) SETERRQ(PetscObjectComm((PetscObject) ce), PETSC_ERR_ARG_WRONG, "Solver was not a TS"); 14900f6b5bSMatthew G. Knepley ierr = TSGetDM((TS) ce->solver, &ce->idm);CHKERRQ(ierr); 15900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 16900f6b5bSMatthew G. Knepley } 17900f6b5bSMatthew G. Knepley 18900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstInitGuessTS_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u) 19900f6b5bSMatthew G. Knepley { 20900f6b5bSMatthew G. Knepley PetscErrorCode ierr; 21900f6b5bSMatthew G. Knepley 22900f6b5bSMatthew G. Knepley PetscFunctionBegin; 23*2e61be88SMatthew G. Knepley ierr = TSComputeInitialCondition((TS) ce->solver, u);CHKERRQ(ierr); 24900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 25900f6b5bSMatthew G. Knepley } 26900f6b5bSMatthew G. Knepley 27900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstComputeErrorTS_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[]) 28900f6b5bSMatthew G. Knepley { 29900f6b5bSMatthew G. Knepley Vec e; 30900f6b5bSMatthew G. Knepley PetscErrorCode ierr; 31900f6b5bSMatthew G. Knepley 32900f6b5bSMatthew G. Knepley PetscFunctionBegin; 33900f6b5bSMatthew G. Knepley ierr = VecDuplicate(u, &e);CHKERRQ(ierr); 34900f6b5bSMatthew G. Knepley ierr = TSComputeExactError((TS) ce->solver, u, e);CHKERRQ(ierr); 35900f6b5bSMatthew G. Knepley ierr = VecNorm(e, NORM_2, errors);CHKERRQ(ierr); 36900f6b5bSMatthew G. Knepley ierr = VecDestroy(&e);CHKERRQ(ierr); 37900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 38900f6b5bSMatthew G. Knepley } 39900f6b5bSMatthew G. Knepley 40900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstGetConvRateTS_Private(PetscConvEst ce, PetscReal alpha[]) 41900f6b5bSMatthew G. Knepley { 42900f6b5bSMatthew G. Knepley TS ts = (TS) ce->solver; 43900f6b5bSMatthew G. Knepley Vec u; 44900f6b5bSMatthew G. Knepley PetscReal *dt, *x, *y, slope, intercept; 45900f6b5bSMatthew G. Knepley PetscInt Ns, oNs, Nf = ce->Nf, f, Nr = ce->Nr, r; 46900f6b5bSMatthew G. Knepley PetscErrorCode ierr; 47900f6b5bSMatthew G. Knepley 48900f6b5bSMatthew G. Knepley PetscFunctionBegin; 49900f6b5bSMatthew G. Knepley ierr = TSGetSolution(ts, &u);CHKERRQ(ierr); 50900f6b5bSMatthew G. Knepley ierr = PetscMalloc1(Nr+1, &dt);CHKERRQ(ierr); 51900f6b5bSMatthew G. Knepley ierr = TSGetTimeStep(ts, &dt[0]);CHKERRQ(ierr); 52900f6b5bSMatthew G. Knepley ierr = TSGetMaxSteps(ts, &oNs);CHKERRQ(ierr); 53900f6b5bSMatthew G. Knepley Ns = oNs; 54900f6b5bSMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 55900f6b5bSMatthew G. Knepley if (r > 0) { 56*2e61be88SMatthew G. Knepley dt[r] = dt[r-1]/ce->r; 57*2e61be88SMatthew G. Knepley Ns = PetscCeilReal(Ns*ce->r); 58900f6b5bSMatthew G. Knepley } 59900f6b5bSMatthew G. Knepley ierr = TSSetTime(ts, 0.0);CHKERRQ(ierr); 60900f6b5bSMatthew G. Knepley ierr = TSSetStepNumber(ts, 0);CHKERRQ(ierr); 61900f6b5bSMatthew G. Knepley ierr = TSSetTimeStep(ts, dt[r]);CHKERRQ(ierr); 62900f6b5bSMatthew G. Knepley ierr = TSSetMaxSteps(ts, Ns);CHKERRQ(ierr); 63900f6b5bSMatthew G. Knepley ierr = PetscConvEstComputeInitialGuess(ce, r, NULL, u);CHKERRQ(ierr); 64900f6b5bSMatthew G. Knepley ierr = TSSolve(ts, u);CHKERRQ(ierr); 65900f6b5bSMatthew G. Knepley ierr = PetscLogEventBegin(ce->event, ce, 0, 0, 0);CHKERRQ(ierr); 66900f6b5bSMatthew G. Knepley ierr = PetscConvEstComputeError(ce, r, NULL, u, &ce->errors[r*Nf]);CHKERRQ(ierr); 67900f6b5bSMatthew G. Knepley ierr = PetscLogEventEnd(ce->event, ce, 0, 0, 0);CHKERRQ(ierr); 68900f6b5bSMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 69900f6b5bSMatthew G. Knepley ierr = PetscLogEventSetDof(ce->event, f, 1.0/dt[r]);CHKERRQ(ierr); 70900f6b5bSMatthew G. Knepley ierr = PetscLogEventSetError(ce->event, f, ce->errors[r*Nf+f]);CHKERRQ(ierr); 71900f6b5bSMatthew G. Knepley } 72900f6b5bSMatthew G. Knepley } 73900f6b5bSMatthew G. Knepley /* Fit convergence rate */ 74900f6b5bSMatthew G. Knepley if (Nr) { 75900f6b5bSMatthew G. Knepley ierr = PetscMalloc2(Nr+1, &x, Nr+1, &y);CHKERRQ(ierr); 76900f6b5bSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 77900f6b5bSMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 78900f6b5bSMatthew G. Knepley x[r] = PetscLog10Real(dt[r*Nf+f]); 79900f6b5bSMatthew G. Knepley y[r] = PetscLog10Real(ce->errors[r*Nf+f]); 80900f6b5bSMatthew G. Knepley } 81900f6b5bSMatthew G. Knepley ierr = PetscLinearRegression(Nr+1, x, y, &slope, &intercept);CHKERRQ(ierr); 82900f6b5bSMatthew G. Knepley /* Since lg err = s lg dt + b */ 83900f6b5bSMatthew G. Knepley alpha[f] = slope; 84900f6b5bSMatthew G. Knepley } 85900f6b5bSMatthew G. Knepley ierr = PetscFree2(x, y);CHKERRQ(ierr); 86900f6b5bSMatthew G. Knepley } 87900f6b5bSMatthew G. Knepley /* Reset solver */ 88900f6b5bSMatthew G. Knepley ierr = TSSetConvergedReason(ts, TS_CONVERGED_ITERATING);CHKERRQ(ierr); 89900f6b5bSMatthew G. Knepley ierr = TSSetTime(ts, 0.0);CHKERRQ(ierr); 90900f6b5bSMatthew G. Knepley ierr = TSSetStepNumber(ts, 0);CHKERRQ(ierr); 91900f6b5bSMatthew G. Knepley ierr = TSSetTimeStep(ts, dt[0]);CHKERRQ(ierr); 92900f6b5bSMatthew G. Knepley ierr = TSSetMaxSteps(ts, oNs);CHKERRQ(ierr); 93900f6b5bSMatthew G. Knepley ierr = PetscConvEstComputeInitialGuess(ce, 0, NULL, u);CHKERRQ(ierr); 94900f6b5bSMatthew G. Knepley ierr = PetscFree(dt);CHKERRQ(ierr); 95900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 96900f6b5bSMatthew G. Knepley } 97900f6b5bSMatthew G. Knepley 98900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstUseTS(PetscConvEst ce) 99900f6b5bSMatthew G. Knepley { 100900f6b5bSMatthew G. Knepley PetscFunctionBegin; 101900f6b5bSMatthew G. Knepley PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 102900f6b5bSMatthew G. Knepley ce->ops->setsolver = PetscConvEstSetTS_Private; 103900f6b5bSMatthew G. Knepley ce->ops->initguess = PetscConvEstInitGuessTS_Private; 104900f6b5bSMatthew G. Knepley ce->ops->computeerror = PetscConvEstComputeErrorTS_Private; 105900f6b5bSMatthew G. Knepley ce->ops->getconvrate = PetscConvEstGetConvRateTS_Private; 106900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 107900f6b5bSMatthew G. Knepley } 108