xref: /petsc/src/ts/trajectory/utils/reconstruct.c (revision 3c633725528e547aaaa9b672a746f5d686a276e1)
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   PetscErrorCode  ierr;
53fe8322adSStefano Zampini 
54fe8322adSStefano Zampini   PetscFunctionBegin;
55fe8322adSStefano Zampini   ierr = TSHistoryGetLocFromTime(tsh,t,&id);CHKERRQ(ierr);
56fe8322adSStefano Zampini   ierr = TSHistoryGetHistory(tsh,&tshn,&tshhist,&tshhist_id,NULL);CHKERRQ(ierr);
57fe8322adSStefano Zampini   if (id == -1 || id == -tshn - 1) {
58fe8322adSStefano Zampini     PetscReal t0 = tshn ? tshhist[0]      : 0.0;
59fe8322adSStefano Zampini     PetscReal tf = tshn ? tshhist[tshn-1] : 0.0;
6098921bdaSJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_PLIB,"Requested time %g is outside the history interval [%g, %g] (%d)",(double)t,(double)t0,(double)tf,tshn);
61fe8322adSStefano Zampini   }
62fe8322adSStefano Zampini   if (tj->monitor) {
63fe8322adSStefano Zampini     ierr = PetscViewerASCIIPrintf(tj->monitor,"Reconstructing at time %g, order %D\n",(double)t,tj->lag.order);CHKERRQ(ierr);
64fe8322adSStefano Zampini   }
65fe8322adSStefano Zampini   if (!tj->lag.T) {
66fe8322adSStefano Zampini     PetscInt o = tj->lag.order+1;
67fe8322adSStefano 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);
68fe8322adSStefano Zampini     for (i = 0; i < o; i++) tj->lag.T[i] = PETSC_MAX_REAL;
69fe8322adSStefano Zampini     ierr = VecDuplicateVecs(U ? U : Udot,o,&tj->lag.W);CHKERRQ(ierr);
70fe8322adSStefano Zampini   }
71fe8322adSStefano Zampini   cnt = 0;
72580bdb30SBarry Smith   ierr = PetscArrayzero(tj->lag.TT,2*(tj->lag.order+1));CHKERRQ(ierr);
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);
79fe8322adSStefano Zampini     if (tj->monitor) {
80fe8322adSStefano Zampini       ierr = PetscViewerASCIIPushTab(tj->monitor);CHKERRQ(ierr);
81fe8322adSStefano Zampini     }
82fe8322adSStefano Zampini 
83fe8322adSStefano Zampini     /* first see if we can reuse any */
84fe8322adSStefano Zampini     for (s = up-1; s >= low; s--) {
85fe8322adSStefano Zampini       PetscReal t = tshhist[s];
86fe8322adSStefano Zampini       PetscInt tid = LagrangeGetId(t,tj->lag.order+1,tj->lag.T,tj->lag.TT);
87fe8322adSStefano Zampini       if (tid < 0) continue;
88fe8322adSStefano Zampini       if (tj->monitor) {
89fe8322adSStefano Zampini         ierr = PetscViewerASCIIPrintf(tj->monitor,"Reusing snapshot %D, step %D, time %g\n",tid,tshhist_id[s],(double)t);CHKERRQ(ierr);
90fe8322adSStefano Zampini       }
91fe8322adSStefano Zampini       tj->lag.TT[tid] = PETSC_TRUE;
92fe8322adSStefano Zampini       tj->lag.WW[cnt] = tj->lag.W[tid];
93fe8322adSStefano Zampini       tj->lag.TW[cnt] = t;
94fe8322adSStefano Zampini       tj->lag.TT[tj->lag.order+1 + s-low] = PETSC_TRUE; /* tell the next loop to skip it */
95fe8322adSStefano Zampini       cnt++;
96fe8322adSStefano Zampini     }
97fe8322adSStefano Zampini 
98fe8322adSStefano Zampini     /* now load the missing ones */
99fe8322adSStefano Zampini     for (s = up-1; s >= low; s--) {
100fe8322adSStefano Zampini       PetscReal t = tshhist[s];
101fe8e0496SStefano Zampini       PetscInt tid;
102fe8e0496SStefano Zampini 
103fe8322adSStefano Zampini       if (tj->lag.TT[tj->lag.order+1 + s-low]) continue;
104fe8e0496SStefano Zampini       tid = LagrangeGetId(t,tj->lag.order+1,tj->lag.T,tj->lag.TT);
105*3c633725SBarry Smith       PetscCheck(tid < 0,PetscObjectComm((PetscObject)tj),PETSC_ERR_PLIB,"This should not happen");
106fe8322adSStefano Zampini       tid = -tid-1;
107fe8322adSStefano Zampini       if (tj->monitor) {
108fe8322adSStefano Zampini         if (tj->lag.T[tid] < PETSC_MAX_REAL) {
109fe8322adSStefano Zampini           ierr = PetscViewerASCIIPrintf(tj->monitor,"Discarding snapshot %D at time %g\n",tid,(double)tj->lag.T[tid]);CHKERRQ(ierr);
110fe8322adSStefano Zampini         } else {
111fe8322adSStefano Zampini           ierr = PetscViewerASCIIPrintf(tj->monitor,"New snapshot %D\n",tid);CHKERRQ(ierr);
112fe8322adSStefano Zampini         }
113fe8322adSStefano Zampini         ierr = PetscViewerASCIIPushTab(tj->monitor);CHKERRQ(ierr);
114fe8322adSStefano Zampini       }
115fe8322adSStefano Zampini       ierr = TSTrajectoryGetVecs(tj,ts,tshhist_id[s],&t,tj->lag.W[tid],NULL);CHKERRQ(ierr);
116fe8322adSStefano Zampini       tj->lag.T[tid] = t;
117fe8322adSStefano Zampini       if (tj->monitor) {
118fe8322adSStefano Zampini         ierr = PetscViewerASCIIPopTab(tj->monitor);CHKERRQ(ierr);
119fe8322adSStefano Zampini       }
120fe8322adSStefano Zampini       tj->lag.TT[tid] = PETSC_TRUE;
121fe8322adSStefano Zampini       tj->lag.WW[cnt] = tj->lag.W[tid];
122fe8322adSStefano Zampini       tj->lag.TW[cnt] = t;
123fe8322adSStefano Zampini       tj->lag.TT[tj->lag.order+1 + s-low] = PETSC_TRUE;
124fe8322adSStefano Zampini       cnt++;
125fe8322adSStefano Zampini     }
126fe8322adSStefano Zampini     if (tj->monitor) {
127fe8322adSStefano Zampini       ierr = PetscViewerASCIIPopTab(tj->monitor);CHKERRQ(ierr);
128fe8322adSStefano Zampini     }
129fe8322adSStefano Zampini   }
130580bdb30SBarry Smith   ierr = PetscArrayzero(tj->lag.TT,tj->lag.order+1);CHKERRQ(ierr);
131fe8322adSStefano Zampini   if (id >=0 && U) { /* requested time match */
132fe8322adSStefano Zampini     PetscInt tid = LagrangeGetId(t,tj->lag.order+1,tj->lag.T,tj->lag.TT);
133fe8322adSStefano Zampini     if (tj->monitor) {
134fe8322adSStefano Zampini       ierr = PetscViewerASCIIPrintf(tj->monitor,"Retrieving solution from exact step\n");CHKERRQ(ierr);
135fe8322adSStefano Zampini       ierr = PetscViewerASCIIPushTab(tj->monitor);CHKERRQ(ierr);
136fe8322adSStefano Zampini     }
137fe8322adSStefano Zampini     if (tid < 0) {
138fe8322adSStefano Zampini       tid = -tid-1;
139fe8322adSStefano Zampini       if (tj->monitor) {
140fe8322adSStefano Zampini         if (tj->lag.T[tid] < PETSC_MAX_REAL) {
141fe8322adSStefano Zampini           ierr = PetscViewerASCIIPrintf(tj->monitor,"Discarding snapshot %D at time %g\n",tid,(double)tj->lag.T[tid]);CHKERRQ(ierr);
142fe8322adSStefano Zampini         } else {
143fe8322adSStefano Zampini           ierr = PetscViewerASCIIPrintf(tj->monitor,"New snapshot %D\n",tid);CHKERRQ(ierr);
144fe8322adSStefano Zampini         }
145fe8322adSStefano Zampini         ierr = PetscViewerASCIIPushTab(tj->monitor);CHKERRQ(ierr);
146fe8322adSStefano Zampini       }
147fe8322adSStefano Zampini       ierr = TSTrajectoryGetVecs(tj,ts,tshhist_id[id],&t,tj->lag.W[tid],NULL);CHKERRQ(ierr);
148fe8322adSStefano Zampini       if (tj->monitor) {
149fe8322adSStefano Zampini         ierr = PetscViewerASCIIPopTab(tj->monitor);CHKERRQ(ierr);
150fe8322adSStefano Zampini       }
151fe8322adSStefano Zampini       tj->lag.T[tid] = t;
152fe8322adSStefano Zampini     } else if (tj->monitor) {
153fe8322adSStefano Zampini       ierr = PetscViewerASCIIPrintf(tj->monitor,"Reusing snapshot %D step %D, time %g\n",tid,tshhist_id[id],(double)t);CHKERRQ(ierr);
154fe8322adSStefano Zampini     }
155fe8322adSStefano Zampini     ierr = VecCopy(tj->lag.W[tid],U);CHKERRQ(ierr);
156fe8322adSStefano Zampini     ierr = PetscObjectStateGet((PetscObject)U,&tj->lag.Ucached.state);CHKERRQ(ierr);
157fe8322adSStefano Zampini     ierr = PetscObjectGetId((PetscObject)U,&tj->lag.Ucached.id);CHKERRQ(ierr);
158fe8322adSStefano Zampini     tj->lag.Ucached.time = t;
159fe8322adSStefano Zampini     tj->lag.Ucached.step = tshhist_id[id];
160fe8322adSStefano Zampini     if (tj->monitor) {
161fe8322adSStefano Zampini       ierr = PetscViewerASCIIPopTab(tj->monitor);CHKERRQ(ierr);
162fe8322adSStefano Zampini     }
163fe8322adSStefano Zampini   }
164fe8322adSStefano Zampini   if (id < 0 && U) {
165fe8322adSStefano Zampini     if (tj->monitor) {
166fe8322adSStefano Zampini       ierr = PetscViewerASCIIPrintf(tj->monitor,"Interpolating solution with %D snapshots\n",cnt);CHKERRQ(ierr);
167fe8322adSStefano Zampini     }
168fe8322adSStefano Zampini     LagrangeBasisVals(cnt,t,tj->lag.TW,tj->lag.L);
169fe8322adSStefano Zampini     ierr = VecZeroEntries(U);CHKERRQ(ierr);
170fe8322adSStefano Zampini     ierr = VecMAXPY(U,cnt,tj->lag.L,tj->lag.WW);CHKERRQ(ierr);
171fe8322adSStefano Zampini     ierr = PetscObjectStateGet((PetscObject)U,&tj->lag.Ucached.state);CHKERRQ(ierr);
172fe8322adSStefano Zampini     ierr = PetscObjectGetId((PetscObject)U,&tj->lag.Ucached.id);CHKERRQ(ierr);
173fe8322adSStefano Zampini     tj->lag.Ucached.time = t;
174fe8322adSStefano Zampini     tj->lag.Ucached.step = PETSC_MIN_INT;
175fe8322adSStefano Zampini   }
176fe8322adSStefano Zampini   if (Udot) {
177fe8322adSStefano Zampini     if (tj->monitor) {
178fe8322adSStefano Zampini       ierr = PetscViewerASCIIPrintf(tj->monitor,"Interpolating derivative with %D snapshots\n",cnt);CHKERRQ(ierr);
179fe8322adSStefano Zampini     }
180fe8322adSStefano Zampini     LagrangeBasisDers(cnt,t,tj->lag.TW,tj->lag.L);
181fe8322adSStefano Zampini     ierr = VecZeroEntries(Udot);CHKERRQ(ierr);
182fe8322adSStefano Zampini     ierr = VecMAXPY(Udot,cnt,tj->lag.L,tj->lag.WW);CHKERRQ(ierr);
183fe8322adSStefano Zampini     ierr = PetscObjectStateGet((PetscObject)Udot,&tj->lag.Udotcached.state);CHKERRQ(ierr);
184fe8322adSStefano Zampini     ierr = PetscObjectGetId((PetscObject)Udot,&tj->lag.Udotcached.id);CHKERRQ(ierr);
185fe8322adSStefano Zampini     tj->lag.Udotcached.time = t;
186fe8322adSStefano Zampini     tj->lag.Udotcached.step = PETSC_MIN_INT;
187fe8322adSStefano Zampini   }
188fe8322adSStefano Zampini   PetscFunctionReturn(0);
189fe8322adSStefano Zampini }
190