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 */ 5fe8322adSStefano Zampini PETSC_STATIC_INLINE void LagrangeBasisVals(PetscInt n,PetscReal t,const PetscReal T[],PetscScalar L[]) 6fe8322adSStefano Zampini { 7fe8322adSStefano Zampini PetscInt k,j; 8fe8322adSStefano Zampini for (k=0; k<n; k++) 9fe8322adSStefano Zampini for (L[k]=1, j=0; j<n; j++) 10fe8322adSStefano Zampini if (j != k) 11fe8322adSStefano Zampini L[k] *= (t - T[j])/(T[k] - T[j]); 12fe8322adSStefano Zampini } 13fe8322adSStefano Zampini 14fe8322adSStefano Zampini PETSC_STATIC_INLINE void LagrangeBasisDers(PetscInt n,PetscReal t,const PetscReal T[],PetscScalar dL[]) 15fe8322adSStefano Zampini { 16fe8322adSStefano Zampini PetscInt k,j,i; 17fe8322adSStefano Zampini for (k=0; k<n; k++) 18fe8322adSStefano Zampini for (dL[k]=0, j=0; j<n; j++) 19fe8322adSStefano Zampini if (j != k) { 20fe8322adSStefano Zampini PetscReal L = 1/(T[k] - T[j]); 21fe8322adSStefano Zampini for (i=0; i<n; i++) 22fe8322adSStefano Zampini if (i != j && i != k) 23fe8322adSStefano Zampini L *= (t - T[i])/(T[k] - T[i]); 24fe8322adSStefano Zampini dL[k] += L; 25fe8322adSStefano Zampini } 26fe8322adSStefano Zampini } 27fe8322adSStefano Zampini 28fe8322adSStefano Zampini PETSC_STATIC_INLINE PetscInt LagrangeGetId(PetscReal t, PetscInt n, const PetscReal T[], const PetscBool Taken[]) 29fe8322adSStefano Zampini { 30fe8322adSStefano Zampini PetscInt _tid = 0; 31fe8322adSStefano Zampini while (_tid < n && PetscAbsReal(t-T[_tid]) > PETSC_SMALL) _tid++; 32fe8322adSStefano Zampini if (_tid < n && !Taken[_tid]) { 33fe8322adSStefano Zampini return _tid; 34fe8322adSStefano Zampini } else { /* we get back a negative id, where the maximum time is stored, since we use usually reconstruct backward in time */ 35fe8322adSStefano Zampini PetscReal max = PETSC_MIN_REAL; 36fe8322adSStefano Zampini PetscInt maxloc = n; 37fe8322adSStefano Zampini _tid = 0; 38fe8322adSStefano Zampini while (_tid < n) { maxloc = (max < T[_tid] && !Taken[_tid]) ? (max = T[_tid],_tid) : maxloc; _tid++; } 39fe8322adSStefano Zampini return -maxloc-1; 40fe8322adSStefano Zampini } 41fe8322adSStefano Zampini } 42fe8322adSStefano Zampini 43fe8322adSStefano Zampini PetscErrorCode TSTrajectoryReconstruct_Private(TSTrajectory tj,TS ts,PetscReal t,Vec U,Vec Udot) 44fe8322adSStefano Zampini { 45fe8322adSStefano Zampini TSHistory tsh = tj->tsh; 46fe8322adSStefano Zampini const PetscReal *tshhist; 47fe8322adSStefano Zampini const PetscInt *tshhist_id; 48fe8322adSStefano Zampini PetscInt id, cnt, i, tshn; 49fe8322adSStefano Zampini PetscErrorCode ierr; 50fe8322adSStefano Zampini 51fe8322adSStefano Zampini PetscFunctionBegin; 52fe8322adSStefano Zampini ierr = TSHistoryGetLocFromTime(tsh,t,&id);CHKERRQ(ierr); 53fe8322adSStefano Zampini ierr = TSHistoryGetHistory(tsh,&tshn,&tshhist,&tshhist_id,NULL);CHKERRQ(ierr); 54fe8322adSStefano Zampini if (id == -1 || id == -tshn - 1) { 55fe8322adSStefano Zampini PetscReal t0 = tshn ? tshhist[0] : 0.0; 56fe8322adSStefano Zampini PetscReal tf = tshn ? tshhist[tshn-1] : 0.0; 57fe8322adSStefano Zampini SETERRQ4(PetscObjectComm((PetscObject)tj),PETSC_ERR_PLIB,"Requested time %g is outside the history interval [%g, %g] (%d)",(double)t,(double)t0,(double)tf,tshn); 58fe8322adSStefano Zampini } 59fe8322adSStefano Zampini if (tj->monitor) { 60fe8322adSStefano Zampini ierr = PetscViewerASCIIPrintf(tj->monitor,"Reconstructing at time %g, order %D\n",(double)t,tj->lag.order);CHKERRQ(ierr); 61fe8322adSStefano Zampini } 62fe8322adSStefano Zampini if (!tj->lag.T) { 63fe8322adSStefano Zampini PetscInt o = tj->lag.order+1; 64fe8322adSStefano Zampini ierr = PetscMalloc5(o,&tj->lag.L,o,&tj->lag.T,o,&tj->lag.WW,2*o,&tj->lag.TT,o,&tj->lag.TW);CHKERRQ(ierr); 65fe8322adSStefano Zampini for (i = 0; i < o; i++) tj->lag.T[i] = PETSC_MAX_REAL; 66fe8322adSStefano Zampini ierr = VecDuplicateVecs(U ? U : Udot,o,&tj->lag.W);CHKERRQ(ierr); 67fe8322adSStefano Zampini } 68fe8322adSStefano Zampini cnt = 0; 69fe8322adSStefano Zampini ierr = PetscMemzero(tj->lag.TT,2*(tj->lag.order+1)*sizeof(PetscBool));CHKERRQ(ierr); 70fe8322adSStefano Zampini if (id < 0 || Udot) { /* populate snapshots for interpolation */ 71fe8322adSStefano Zampini PetscInt s,nid = id < 0 ? -(id+1) : id; 72fe8322adSStefano Zampini 73fe8322adSStefano Zampini PetscInt up = PetscMin(nid + tj->lag.order/2+1,tshn); 74fe8322adSStefano Zampini PetscInt low = PetscMax(up-tj->lag.order-1,0); 75fe8322adSStefano Zampini up = PetscMin(PetscMax(low + tj->lag.order + 1,up),tshn); 76fe8322adSStefano Zampini if (tj->monitor) { 77fe8322adSStefano Zampini ierr = PetscViewerASCIIPushTab(tj->monitor);CHKERRQ(ierr); 78fe8322adSStefano Zampini } 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) { 86fe8322adSStefano Zampini ierr = PetscViewerASCIIPrintf(tj->monitor,"Reusing snapshot %D, step %D, time %g\n",tid,tshhist_id[s],(double)t);CHKERRQ(ierr); 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]; 98*fe8e0496SStefano Zampini PetscInt tid; 99*fe8e0496SStefano Zampini 100fe8322adSStefano Zampini if (tj->lag.TT[tj->lag.order+1 + s-low]) continue; 101*fe8e0496SStefano Zampini tid = LagrangeGetId(t,tj->lag.order+1,tj->lag.T,tj->lag.TT); 102fe8322adSStefano Zampini if (tid >= 0) SETERRQ(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) { 106fe8322adSStefano Zampini ierr = PetscViewerASCIIPrintf(tj->monitor,"Discarding snapshot %D at time %g\n",tid,(double)tj->lag.T[tid]);CHKERRQ(ierr); 107fe8322adSStefano Zampini } else { 108fe8322adSStefano Zampini ierr = PetscViewerASCIIPrintf(tj->monitor,"New snapshot %D\n",tid);CHKERRQ(ierr); 109fe8322adSStefano Zampini } 110fe8322adSStefano Zampini ierr = PetscViewerASCIIPushTab(tj->monitor);CHKERRQ(ierr); 111fe8322adSStefano Zampini } 112fe8322adSStefano Zampini ierr = TSTrajectoryGetVecs(tj,ts,tshhist_id[s],&t,tj->lag.W[tid],NULL);CHKERRQ(ierr); 113fe8322adSStefano Zampini tj->lag.T[tid] = t; 114fe8322adSStefano Zampini if (tj->monitor) { 115fe8322adSStefano Zampini ierr = PetscViewerASCIIPopTab(tj->monitor);CHKERRQ(ierr); 116fe8322adSStefano Zampini } 117fe8322adSStefano Zampini tj->lag.TT[tid] = PETSC_TRUE; 118fe8322adSStefano Zampini tj->lag.WW[cnt] = tj->lag.W[tid]; 119fe8322adSStefano Zampini tj->lag.TW[cnt] = t; 120fe8322adSStefano Zampini tj->lag.TT[tj->lag.order+1 + s-low] = PETSC_TRUE; 121fe8322adSStefano Zampini cnt++; 122fe8322adSStefano Zampini } 123fe8322adSStefano Zampini if (tj->monitor) { 124fe8322adSStefano Zampini ierr = PetscViewerASCIIPopTab(tj->monitor);CHKERRQ(ierr); 125fe8322adSStefano Zampini } 126fe8322adSStefano Zampini } 127fe8322adSStefano Zampini ierr = PetscMemzero(tj->lag.TT,(tj->lag.order+1)*sizeof(PetscBool));CHKERRQ(ierr); 128fe8322adSStefano Zampini if (id >=0 && U) { /* requested time match */ 129fe8322adSStefano Zampini PetscInt tid = LagrangeGetId(t,tj->lag.order+1,tj->lag.T,tj->lag.TT); 130fe8322adSStefano Zampini if (tj->monitor) { 131fe8322adSStefano Zampini ierr = PetscViewerASCIIPrintf(tj->monitor,"Retrieving solution from exact step\n");CHKERRQ(ierr); 132fe8322adSStefano Zampini ierr = PetscViewerASCIIPushTab(tj->monitor);CHKERRQ(ierr); 133fe8322adSStefano Zampini } 134fe8322adSStefano Zampini if (tid < 0) { 135fe8322adSStefano Zampini tid = -tid-1; 136fe8322adSStefano Zampini if (tj->monitor) { 137fe8322adSStefano Zampini if (tj->lag.T[tid] < PETSC_MAX_REAL) { 138fe8322adSStefano Zampini ierr = PetscViewerASCIIPrintf(tj->monitor,"Discarding snapshot %D at time %g\n",tid,(double)tj->lag.T[tid]);CHKERRQ(ierr); 139fe8322adSStefano Zampini } else { 140fe8322adSStefano Zampini ierr = PetscViewerASCIIPrintf(tj->monitor,"New snapshot %D\n",tid);CHKERRQ(ierr); 141fe8322adSStefano Zampini } 142fe8322adSStefano Zampini ierr = PetscViewerASCIIPushTab(tj->monitor);CHKERRQ(ierr); 143fe8322adSStefano Zampini } 144fe8322adSStefano Zampini ierr = TSTrajectoryGetVecs(tj,ts,tshhist_id[id],&t,tj->lag.W[tid],NULL);CHKERRQ(ierr); 145fe8322adSStefano Zampini if (tj->monitor) { 146fe8322adSStefano Zampini ierr = PetscViewerASCIIPopTab(tj->monitor);CHKERRQ(ierr); 147fe8322adSStefano Zampini } 148fe8322adSStefano Zampini tj->lag.T[tid] = t; 149fe8322adSStefano Zampini } else if (tj->monitor) { 150fe8322adSStefano Zampini ierr = PetscViewerASCIIPrintf(tj->monitor,"Reusing snapshot %D step %D, time %g\n",tid,tshhist_id[id],(double)t);CHKERRQ(ierr); 151fe8322adSStefano Zampini } 152fe8322adSStefano Zampini ierr = VecCopy(tj->lag.W[tid],U);CHKERRQ(ierr); 153fe8322adSStefano Zampini ierr = PetscObjectStateGet((PetscObject)U,&tj->lag.Ucached.state);CHKERRQ(ierr); 154fe8322adSStefano Zampini ierr = PetscObjectGetId((PetscObject)U,&tj->lag.Ucached.id);CHKERRQ(ierr); 155fe8322adSStefano Zampini tj->lag.Ucached.time = t; 156fe8322adSStefano Zampini tj->lag.Ucached.step = tshhist_id[id]; 157fe8322adSStefano Zampini if (tj->monitor) { 158fe8322adSStefano Zampini ierr = PetscViewerASCIIPopTab(tj->monitor);CHKERRQ(ierr); 159fe8322adSStefano Zampini } 160fe8322adSStefano Zampini } 161fe8322adSStefano Zampini if (id < 0 && U) { 162fe8322adSStefano Zampini if (tj->monitor) { 163fe8322adSStefano Zampini ierr = PetscViewerASCIIPrintf(tj->monitor,"Interpolating solution with %D snapshots\n",cnt);CHKERRQ(ierr); 164fe8322adSStefano Zampini } 165fe8322adSStefano Zampini LagrangeBasisVals(cnt,t,tj->lag.TW,tj->lag.L); 166fe8322adSStefano Zampini ierr = VecZeroEntries(U);CHKERRQ(ierr); 167fe8322adSStefano Zampini ierr = VecMAXPY(U,cnt,tj->lag.L,tj->lag.WW);CHKERRQ(ierr); 168fe8322adSStefano Zampini ierr = PetscObjectStateGet((PetscObject)U,&tj->lag.Ucached.state);CHKERRQ(ierr); 169fe8322adSStefano Zampini ierr = PetscObjectGetId((PetscObject)U,&tj->lag.Ucached.id);CHKERRQ(ierr); 170fe8322adSStefano Zampini tj->lag.Ucached.time = t; 171fe8322adSStefano Zampini tj->lag.Ucached.step = PETSC_MIN_INT; 172fe8322adSStefano Zampini } 173fe8322adSStefano Zampini if (Udot) { 174fe8322adSStefano Zampini if (tj->monitor) { 175fe8322adSStefano Zampini ierr = PetscViewerASCIIPrintf(tj->monitor,"Interpolating derivative with %D snapshots\n",cnt);CHKERRQ(ierr); 176fe8322adSStefano Zampini } 177fe8322adSStefano Zampini LagrangeBasisDers(cnt,t,tj->lag.TW,tj->lag.L); 178fe8322adSStefano Zampini ierr = VecZeroEntries(Udot);CHKERRQ(ierr); 179fe8322adSStefano Zampini ierr = VecMAXPY(Udot,cnt,tj->lag.L,tj->lag.WW);CHKERRQ(ierr); 180fe8322adSStefano Zampini ierr = PetscObjectStateGet((PetscObject)Udot,&tj->lag.Udotcached.state);CHKERRQ(ierr); 181fe8322adSStefano Zampini ierr = PetscObjectGetId((PetscObject)Udot,&tj->lag.Udotcached.id);CHKERRQ(ierr); 182fe8322adSStefano Zampini tj->lag.Udotcached.time = t; 183fe8322adSStefano Zampini tj->lag.Udotcached.step = PETSC_MIN_INT; 184fe8322adSStefano Zampini } 185fe8322adSStefano Zampini PetscFunctionReturn(0); 186fe8322adSStefano Zampini } 187