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