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 */ 59fbee547SJacob Faibussowitsch static inline void LagrangeBasisVals(PetscInt n,PetscReal t,const PetscReal T[],PetscScalar L[]) 6fe8322adSStefano Zampini { 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 159fbee547SJacob Faibussowitsch static inline void LagrangeBasisDers(PetscInt n,PetscReal t,const PetscReal T[],PetscScalar dL[]) 16fe8322adSStefano Zampini { 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 319fbee547SJacob Faibussowitsch static inline PetscInt LagrangeGetId(PetscReal t, PetscInt n, const PetscReal T[], const PetscBool Taken[]) 32fe8322adSStefano Zampini { 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; 41fe8322adSStefano Zampini while (_tid < n) { maxloc = (max < T[_tid] && !Taken[_tid]) ? (max = T[_tid],_tid) : maxloc; _tid++; } 42fe8322adSStefano Zampini return -maxloc-1; 43fe8322adSStefano Zampini } 44fe8322adSStefano Zampini } 45fe8322adSStefano Zampini 46fe8322adSStefano Zampini PetscErrorCode TSTrajectoryReconstruct_Private(TSTrajectory tj,TS ts,PetscReal t,Vec U,Vec Udot) 47fe8322adSStefano Zampini { 48fe8322adSStefano Zampini TSHistory tsh = tj->tsh; 49fe8322adSStefano Zampini const PetscReal *tshhist; 50fe8322adSStefano Zampini const PetscInt *tshhist_id; 51fe8322adSStefano Zampini PetscInt id, cnt, i, tshn; 52fe8322adSStefano Zampini 53fe8322adSStefano Zampini PetscFunctionBegin; 549566063dSJacob Faibussowitsch PetscCall(TSHistoryGetLocFromTime(tsh,t,&id)); 559566063dSJacob Faibussowitsch PetscCall(TSHistoryGetHistory(tsh,&tshn,&tshhist,&tshhist_id,NULL)); 56fe8322adSStefano Zampini if (id == -1 || id == -tshn - 1) { 57fe8322adSStefano Zampini PetscReal t0 = tshn ? tshhist[0] : 0.0; 58fe8322adSStefano Zampini PetscReal tf = tshn ? tshhist[tshn-1] : 0.0; 5963a3b9bcSJacob 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); 60fe8322adSStefano Zampini } 61fe8322adSStefano Zampini if (tj->monitor) { 6263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Reconstructing at time %g, order %" PetscInt_FMT "\n",(double)t,tj->lag.order)); 63fe8322adSStefano Zampini } 64fe8322adSStefano Zampini if (!tj->lag.T) { 65fe8322adSStefano Zampini PetscInt o = tj->lag.order+1; 669566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(o,&tj->lag.L,o,&tj->lag.T,o,&tj->lag.WW,2*o,&tj->lag.TT,o,&tj->lag.TW)); 67fe8322adSStefano Zampini for (i = 0; i < o; i++) tj->lag.T[i] = PETSC_MAX_REAL; 689566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(U ? U : Udot,o,&tj->lag.W)); 69fe8322adSStefano Zampini } 70fe8322adSStefano Zampini cnt = 0; 719566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(tj->lag.TT,2*(tj->lag.order+1))); 72fe8322adSStefano Zampini if (id < 0 || Udot) { /* populate snapshots for interpolation */ 73fe8322adSStefano Zampini PetscInt s,nid = id < 0 ? -(id+1) : id; 74fe8322adSStefano Zampini 75fe8322adSStefano Zampini PetscInt up = PetscMin(nid + tj->lag.order/2+1,tshn); 76fe8322adSStefano Zampini PetscInt low = PetscMax(up-tj->lag.order-1,0); 77fe8322adSStefano Zampini up = PetscMin(PetscMax(low + tj->lag.order + 1,up),tshn); 78*1baa6e33SBarry Smith if (tj->monitor) PetscCall(PetscViewerASCIIPushTab(tj->monitor)); 79fe8322adSStefano Zampini 80fe8322adSStefano Zampini /* first see if we can reuse any */ 81fe8322adSStefano Zampini for (s = up-1; s >= low; s--) { 82fe8322adSStefano Zampini PetscReal t = tshhist[s]; 83fe8322adSStefano Zampini PetscInt tid = LagrangeGetId(t,tj->lag.order+1,tj->lag.T,tj->lag.TT); 84fe8322adSStefano Zampini if (tid < 0) continue; 85fe8322adSStefano Zampini if (tj->monitor) { 8663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Reusing snapshot %" PetscInt_FMT ", step %" PetscInt_FMT ", time %g\n",tid,tshhist_id[s],(double)t)); 87fe8322adSStefano Zampini } 88fe8322adSStefano Zampini tj->lag.TT[tid] = PETSC_TRUE; 89fe8322adSStefano Zampini tj->lag.WW[cnt] = tj->lag.W[tid]; 90fe8322adSStefano Zampini tj->lag.TW[cnt] = t; 91fe8322adSStefano Zampini tj->lag.TT[tj->lag.order+1 + s-low] = PETSC_TRUE; /* tell the next loop to skip it */ 92fe8322adSStefano Zampini cnt++; 93fe8322adSStefano Zampini } 94fe8322adSStefano Zampini 95fe8322adSStefano Zampini /* now load the missing ones */ 96fe8322adSStefano Zampini for (s = up-1; s >= low; s--) { 97fe8322adSStefano Zampini PetscReal t = tshhist[s]; 98fe8e0496SStefano Zampini PetscInt tid; 99fe8e0496SStefano Zampini 100fe8322adSStefano Zampini if (tj->lag.TT[tj->lag.order+1 + s-low]) continue; 101fe8e0496SStefano Zampini tid = LagrangeGetId(t,tj->lag.order+1,tj->lag.T,tj->lag.TT); 1023c633725SBarry Smith PetscCheck(tid < 0,PetscObjectComm((PetscObject)tj),PETSC_ERR_PLIB,"This should not happen"); 103fe8322adSStefano Zampini tid = -tid-1; 104fe8322adSStefano Zampini if (tj->monitor) { 105fe8322adSStefano Zampini if (tj->lag.T[tid] < PETSC_MAX_REAL) { 10663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Discarding snapshot %" PetscInt_FMT " at time %g\n",tid,(double)tj->lag.T[tid])); 107fe8322adSStefano Zampini } else { 10863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(tj->monitor,"New snapshot %" PetscInt_FMT "\n",tid)); 109fe8322adSStefano Zampini } 1109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(tj->monitor)); 111fe8322adSStefano Zampini } 1129566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGetVecs(tj,ts,tshhist_id[s],&t,tj->lag.W[tid],NULL)); 113fe8322adSStefano Zampini tj->lag.T[tid] = t; 114*1baa6e33SBarry Smith if (tj->monitor) PetscCall(PetscViewerASCIIPopTab(tj->monitor)); 115fe8322adSStefano Zampini tj->lag.TT[tid] = PETSC_TRUE; 116fe8322adSStefano Zampini tj->lag.WW[cnt] = tj->lag.W[tid]; 117fe8322adSStefano Zampini tj->lag.TW[cnt] = t; 118fe8322adSStefano Zampini tj->lag.TT[tj->lag.order+1 + s-low] = PETSC_TRUE; 119fe8322adSStefano Zampini cnt++; 120fe8322adSStefano Zampini } 121*1baa6e33SBarry Smith if (tj->monitor) PetscCall(PetscViewerASCIIPopTab(tj->monitor)); 122fe8322adSStefano Zampini } 1239566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(tj->lag.TT,tj->lag.order+1)); 124fe8322adSStefano Zampini if (id >=0 && U) { /* requested time match */ 125fe8322adSStefano Zampini PetscInt tid = LagrangeGetId(t,tj->lag.order+1,tj->lag.T,tj->lag.TT); 126fe8322adSStefano Zampini if (tj->monitor) { 1279566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Retrieving solution from exact step\n")); 1289566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(tj->monitor)); 129fe8322adSStefano Zampini } 130fe8322adSStefano Zampini if (tid < 0) { 131fe8322adSStefano Zampini tid = -tid-1; 132fe8322adSStefano Zampini if (tj->monitor) { 133fe8322adSStefano Zampini if (tj->lag.T[tid] < PETSC_MAX_REAL) { 13463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Discarding snapshot %" PetscInt_FMT " at time %g\n",tid,(double)tj->lag.T[tid])); 135fe8322adSStefano Zampini } else { 13663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(tj->monitor,"New snapshot %" PetscInt_FMT "\n",tid)); 137fe8322adSStefano Zampini } 1389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(tj->monitor)); 139fe8322adSStefano Zampini } 1409566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGetVecs(tj,ts,tshhist_id[id],&t,tj->lag.W[tid],NULL)); 141*1baa6e33SBarry Smith if (tj->monitor) PetscCall(PetscViewerASCIIPopTab(tj->monitor)); 142fe8322adSStefano Zampini tj->lag.T[tid] = t; 143fe8322adSStefano Zampini } else if (tj->monitor) { 14463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Reusing snapshot %" PetscInt_FMT " step %" PetscInt_FMT ", time %g\n",tid,tshhist_id[id],(double)t)); 145fe8322adSStefano Zampini } 1469566063dSJacob Faibussowitsch PetscCall(VecCopy(tj->lag.W[tid],U)); 1479566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)U,&tj->lag.Ucached.state)); 1489566063dSJacob Faibussowitsch PetscCall(PetscObjectGetId((PetscObject)U,&tj->lag.Ucached.id)); 149fe8322adSStefano Zampini tj->lag.Ucached.time = t; 150fe8322adSStefano Zampini tj->lag.Ucached.step = tshhist_id[id]; 151*1baa6e33SBarry Smith if (tj->monitor) PetscCall(PetscViewerASCIIPopTab(tj->monitor)); 152fe8322adSStefano Zampini } 153fe8322adSStefano Zampini if (id < 0 && U) { 154fe8322adSStefano Zampini if (tj->monitor) { 15563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Interpolating solution with %" PetscInt_FMT " snapshots\n",cnt)); 156fe8322adSStefano Zampini } 157fe8322adSStefano Zampini LagrangeBasisVals(cnt,t,tj->lag.TW,tj->lag.L); 1589566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(U)); 1599566063dSJacob Faibussowitsch PetscCall(VecMAXPY(U,cnt,tj->lag.L,tj->lag.WW)); 1609566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)U,&tj->lag.Ucached.state)); 1619566063dSJacob Faibussowitsch PetscCall(PetscObjectGetId((PetscObject)U,&tj->lag.Ucached.id)); 162fe8322adSStefano Zampini tj->lag.Ucached.time = t; 163fe8322adSStefano Zampini tj->lag.Ucached.step = PETSC_MIN_INT; 164fe8322adSStefano Zampini } 165fe8322adSStefano Zampini if (Udot) { 166fe8322adSStefano Zampini if (tj->monitor) { 16763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Interpolating derivative with %" PetscInt_FMT " snapshots\n",cnt)); 168fe8322adSStefano Zampini } 169fe8322adSStefano Zampini LagrangeBasisDers(cnt,t,tj->lag.TW,tj->lag.L); 1709566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Udot)); 1719566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Udot,cnt,tj->lag.L,tj->lag.WW)); 1729566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)Udot,&tj->lag.Udotcached.state)); 1739566063dSJacob Faibussowitsch PetscCall(PetscObjectGetId((PetscObject)Udot,&tj->lag.Udotcached.id)); 174fe8322adSStefano Zampini tj->lag.Udotcached.time = t; 175fe8322adSStefano Zampini tj->lag.Udotcached.step = PETSC_MIN_INT; 176fe8322adSStefano Zampini } 177fe8322adSStefano Zampini PetscFunctionReturn(0); 178fe8322adSStefano Zampini } 179