1fe8322adSStefano Zampini #include <petsc/private/tshistoryimpl.h> 2fe8322adSStefano Zampini #include <petscts.h> 3fe8322adSStefano Zampini 4fe8322adSStefano Zampini /* these two functions have been stolen from bdf.c */ 5d71ae5a4SJacob Faibussowitsch static inline void LagrangeBasisVals(PetscInt n, PetscReal t, const PetscReal T[], PetscScalar L[]) 6d71ae5a4SJacob Faibussowitsch { 7fe8322adSStefano Zampini PetscInt k, j; 82da392ccSBarry Smith for (k = 0; k < n; k++) { 92da392ccSBarry Smith for (L[k] = 1, j = 0; j < n; j++) { 102da392ccSBarry Smith if (j != k) L[k] *= (t - T[j]) / (T[k] - T[j]); 112da392ccSBarry Smith } 122da392ccSBarry Smith } 13fe8322adSStefano Zampini } 14fe8322adSStefano Zampini 15d71ae5a4SJacob Faibussowitsch static inline void LagrangeBasisDers(PetscInt n, PetscReal t, const PetscReal T[], PetscScalar dL[]) 16d71ae5a4SJacob Faibussowitsch { 17fe8322adSStefano Zampini PetscInt k, j, i; 182da392ccSBarry Smith for (k = 0; k < n; k++) { 192da392ccSBarry Smith for (dL[k] = 0, j = 0; j < n; j++) { 20fe8322adSStefano Zampini if (j != k) { 21fe8322adSStefano Zampini PetscReal L = 1 / (T[k] - T[j]); 222da392ccSBarry Smith for (i = 0; i < n; i++) { 232da392ccSBarry Smith if (i != j && i != k) L *= (t - T[i]) / (T[k] - T[i]); 242da392ccSBarry Smith } 25fe8322adSStefano Zampini dL[k] += L; 26fe8322adSStefano Zampini } 27fe8322adSStefano Zampini } 282da392ccSBarry Smith } 292da392ccSBarry Smith } 30fe8322adSStefano Zampini 31d71ae5a4SJacob Faibussowitsch static inline PetscInt LagrangeGetId(PetscReal t, PetscInt n, const PetscReal T[], const PetscBool Taken[]) 32d71ae5a4SJacob Faibussowitsch { 33fe8322adSStefano Zampini PetscInt _tid = 0; 34fe8322adSStefano Zampini while (_tid < n && PetscAbsReal(t - T[_tid]) > PETSC_SMALL) _tid++; 35fe8322adSStefano Zampini if (_tid < n && !Taken[_tid]) { 36fe8322adSStefano Zampini return _tid; 37fe8322adSStefano Zampini } else { /* we get back a negative id, where the maximum time is stored, since we use usually reconstruct backward in time */ 38fe8322adSStefano Zampini PetscReal max = PETSC_MIN_REAL; 39fe8322adSStefano Zampini PetscInt maxloc = n; 40fe8322adSStefano Zampini _tid = 0; 419371c9d4SSatish Balay while (_tid < n) { 429371c9d4SSatish Balay maxloc = (max < T[_tid] && !Taken[_tid]) ? (max = T[_tid], _tid) : maxloc; 439371c9d4SSatish Balay _tid++; 449371c9d4SSatish Balay } 45fe8322adSStefano Zampini return -maxloc - 1; 46fe8322adSStefano Zampini } 47fe8322adSStefano Zampini } 48fe8322adSStefano Zampini 49d71ae5a4SJacob Faibussowitsch PetscErrorCode TSTrajectoryReconstruct_Private(TSTrajectory tj, TS ts, PetscReal t, Vec U, Vec Udot) 50d71ae5a4SJacob Faibussowitsch { 51fe8322adSStefano Zampini TSHistory tsh = tj->tsh; 52fe8322adSStefano Zampini const PetscReal *tshhist; 53fe8322adSStefano Zampini const PetscInt *tshhist_id; 54fe8322adSStefano Zampini PetscInt id, cnt, i, tshn; 55fe8322adSStefano Zampini 56fe8322adSStefano Zampini PetscFunctionBegin; 579566063dSJacob Faibussowitsch PetscCall(TSHistoryGetLocFromTime(tsh, t, &id)); 589566063dSJacob Faibussowitsch PetscCall(TSHistoryGetHistory(tsh, &tshn, &tshhist, &tshhist_id, NULL)); 59fe8322adSStefano Zampini if (id == -1 || id == -tshn - 1) { 60fe8322adSStefano Zampini PetscReal t0 = tshn ? tshhist[0] : 0.0; 61fe8322adSStefano Zampini PetscReal tf = tshn ? tshhist[tshn - 1] : 0.0; 6263a3b9bcSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)tj), PETSC_ERR_PLIB, "Requested time %g is outside the history interval [%g, %g] (%" PetscInt_FMT ")", (double)t, (double)t0, (double)tf, tshn); 63fe8322adSStefano Zampini } 6448a46eb9SPierre Jolivet if (tj->monitor) PetscCall(PetscViewerASCIIPrintf(tj->monitor, "Reconstructing at time %g, order %" PetscInt_FMT "\n", (double)t, tj->lag.order)); 65fe8322adSStefano Zampini if (!tj->lag.T) { 66fe8322adSStefano Zampini PetscInt o = tj->lag.order + 1; 679566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(o, &tj->lag.L, o, &tj->lag.T, o, &tj->lag.WW, 2 * o, &tj->lag.TT, o, &tj->lag.TW)); 68fe8322adSStefano Zampini for (i = 0; i < o; i++) tj->lag.T[i] = PETSC_MAX_REAL; 699566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(U ? U : Udot, o, &tj->lag.W)); 70fe8322adSStefano Zampini } 71fe8322adSStefano Zampini cnt = 0; 729566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(tj->lag.TT, 2 * (tj->lag.order + 1))); 73fe8322adSStefano Zampini if (id < 0 || Udot) { /* populate snapshots for interpolation */ 74fe8322adSStefano Zampini PetscInt s, nid = id < 0 ? -(id + 1) : id; 75fe8322adSStefano Zampini 76fe8322adSStefano Zampini PetscInt up = PetscMin(nid + tj->lag.order / 2 + 1, tshn); 77fe8322adSStefano Zampini PetscInt low = PetscMax(up - tj->lag.order - 1, 0); 78fe8322adSStefano Zampini up = PetscMin(PetscMax(low + tj->lag.order + 1, up), tshn); 791baa6e33SBarry Smith if (tj->monitor) PetscCall(PetscViewerASCIIPushTab(tj->monitor)); 80fe8322adSStefano Zampini 81fe8322adSStefano Zampini /* first see if we can reuse any */ 82fe8322adSStefano Zampini for (s = up - 1; s >= low; s--) { 83fe8322adSStefano Zampini PetscReal t = tshhist[s]; 84fe8322adSStefano Zampini PetscInt tid = LagrangeGetId(t, tj->lag.order + 1, tj->lag.T, tj->lag.TT); 85fe8322adSStefano Zampini if (tid < 0) continue; 8648a46eb9SPierre Jolivet if (tj->monitor) PetscCall(PetscViewerASCIIPrintf(tj->monitor, "Reusing snapshot %" PetscInt_FMT ", step %" PetscInt_FMT ", time %g\n", tid, tshhist_id[s], (double)t)); 87fe8322adSStefano Zampini tj->lag.TT[tid] = PETSC_TRUE; 88fe8322adSStefano Zampini tj->lag.WW[cnt] = tj->lag.W[tid]; 89fe8322adSStefano Zampini tj->lag.TW[cnt] = t; 90fe8322adSStefano Zampini tj->lag.TT[tj->lag.order + 1 + s - low] = PETSC_TRUE; /* tell the next loop to skip it */ 91fe8322adSStefano Zampini cnt++; 92fe8322adSStefano Zampini } 93fe8322adSStefano Zampini 94fe8322adSStefano Zampini /* now load the missing ones */ 95fe8322adSStefano Zampini for (s = up - 1; s >= low; s--) { 96fe8322adSStefano Zampini PetscReal t = tshhist[s]; 97fe8e0496SStefano Zampini PetscInt tid; 98fe8e0496SStefano Zampini 99fe8322adSStefano Zampini if (tj->lag.TT[tj->lag.order + 1 + s - low]) continue; 100fe8e0496SStefano Zampini tid = LagrangeGetId(t, tj->lag.order + 1, tj->lag.T, tj->lag.TT); 1013c633725SBarry Smith PetscCheck(tid < 0, PetscObjectComm((PetscObject)tj), PETSC_ERR_PLIB, "This should not happen"); 102fe8322adSStefano Zampini tid = -tid - 1; 103fe8322adSStefano Zampini if (tj->monitor) { 104fe8322adSStefano Zampini if (tj->lag.T[tid] < PETSC_MAX_REAL) { 10563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(tj->monitor, "Discarding snapshot %" PetscInt_FMT " at time %g\n", tid, (double)tj->lag.T[tid])); 106fe8322adSStefano Zampini } else { 10763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(tj->monitor, "New snapshot %" PetscInt_FMT "\n", tid)); 108fe8322adSStefano Zampini } 1099566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(tj->monitor)); 110fe8322adSStefano Zampini } 1119566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGetVecs(tj, ts, tshhist_id[s], &t, tj->lag.W[tid], NULL)); 112fe8322adSStefano Zampini tj->lag.T[tid] = t; 1131baa6e33SBarry Smith if (tj->monitor) PetscCall(PetscViewerASCIIPopTab(tj->monitor)); 114fe8322adSStefano Zampini tj->lag.TT[tid] = PETSC_TRUE; 115fe8322adSStefano Zampini tj->lag.WW[cnt] = tj->lag.W[tid]; 116fe8322adSStefano Zampini tj->lag.TW[cnt] = t; 117fe8322adSStefano Zampini tj->lag.TT[tj->lag.order + 1 + s - low] = PETSC_TRUE; 118fe8322adSStefano Zampini cnt++; 119fe8322adSStefano Zampini } 1201baa6e33SBarry Smith if (tj->monitor) PetscCall(PetscViewerASCIIPopTab(tj->monitor)); 121fe8322adSStefano Zampini } 1229566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(tj->lag.TT, tj->lag.order + 1)); 123fe8322adSStefano Zampini if (id >= 0 && U) { /* requested time match */ 124fe8322adSStefano Zampini PetscInt tid = LagrangeGetId(t, tj->lag.order + 1, tj->lag.T, tj->lag.TT); 125fe8322adSStefano Zampini if (tj->monitor) { 1269566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(tj->monitor, "Retrieving solution from exact step\n")); 1279566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(tj->monitor)); 128fe8322adSStefano Zampini } 129fe8322adSStefano Zampini if (tid < 0) { 130fe8322adSStefano Zampini tid = -tid - 1; 131fe8322adSStefano Zampini if (tj->monitor) { 132fe8322adSStefano Zampini if (tj->lag.T[tid] < PETSC_MAX_REAL) { 13363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(tj->monitor, "Discarding snapshot %" PetscInt_FMT " at time %g\n", tid, (double)tj->lag.T[tid])); 134fe8322adSStefano Zampini } else { 13563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(tj->monitor, "New snapshot %" PetscInt_FMT "\n", tid)); 136fe8322adSStefano Zampini } 1379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(tj->monitor)); 138fe8322adSStefano Zampini } 1399566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGetVecs(tj, ts, tshhist_id[id], &t, tj->lag.W[tid], NULL)); 1401baa6e33SBarry Smith if (tj->monitor) PetscCall(PetscViewerASCIIPopTab(tj->monitor)); 141fe8322adSStefano Zampini tj->lag.T[tid] = t; 142fe8322adSStefano Zampini } else if (tj->monitor) { 14363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(tj->monitor, "Reusing snapshot %" PetscInt_FMT " step %" PetscInt_FMT ", time %g\n", tid, tshhist_id[id], (double)t)); 144fe8322adSStefano Zampini } 1459566063dSJacob Faibussowitsch PetscCall(VecCopy(tj->lag.W[tid], U)); 1469566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)U, &tj->lag.Ucached.state)); 1479566063dSJacob Faibussowitsch PetscCall(PetscObjectGetId((PetscObject)U, &tj->lag.Ucached.id)); 148fe8322adSStefano Zampini tj->lag.Ucached.time = t; 149fe8322adSStefano Zampini tj->lag.Ucached.step = tshhist_id[id]; 1501baa6e33SBarry Smith if (tj->monitor) PetscCall(PetscViewerASCIIPopTab(tj->monitor)); 151fe8322adSStefano Zampini } 152fe8322adSStefano Zampini if (id < 0 && U) { 15348a46eb9SPierre Jolivet if (tj->monitor) PetscCall(PetscViewerASCIIPrintf(tj->monitor, "Interpolating solution with %" PetscInt_FMT " snapshots\n", cnt)); 154fe8322adSStefano Zampini LagrangeBasisVals(cnt, t, tj->lag.TW, tj->lag.L); 1559566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(U)); 1569566063dSJacob Faibussowitsch PetscCall(VecMAXPY(U, cnt, tj->lag.L, tj->lag.WW)); 1579566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)U, &tj->lag.Ucached.state)); 1589566063dSJacob Faibussowitsch PetscCall(PetscObjectGetId((PetscObject)U, &tj->lag.Ucached.id)); 159fe8322adSStefano Zampini tj->lag.Ucached.time = t; 160*1690c2aeSBarry Smith tj->lag.Ucached.step = PETSC_INT_MIN; 161fe8322adSStefano Zampini } 162fe8322adSStefano Zampini if (Udot) { 16348a46eb9SPierre Jolivet if (tj->monitor) PetscCall(PetscViewerASCIIPrintf(tj->monitor, "Interpolating derivative with %" PetscInt_FMT " snapshots\n", cnt)); 164fe8322adSStefano Zampini LagrangeBasisDers(cnt, t, tj->lag.TW, tj->lag.L); 1659566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Udot)); 1669566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Udot, cnt, tj->lag.L, tj->lag.WW)); 1679566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)Udot, &tj->lag.Udotcached.state)); 1689566063dSJacob Faibussowitsch PetscCall(PetscObjectGetId((PetscObject)Udot, &tj->lag.Udotcached.id)); 169fe8322adSStefano Zampini tj->lag.Udotcached.time = t; 170*1690c2aeSBarry Smith tj->lag.Udotcached.step = PETSC_INT_MIN; 171fe8322adSStefano Zampini } 1723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 173fe8322adSStefano Zampini } 174