xref: /petsc/src/ts/interface/tshistory.c (revision f4f49eeac7efa77fffa46b7ff95a3ed169f659ed)
189818f9dSStefano Zampini #include <petsc/private/tshistoryimpl.h>
289818f9dSStefano Zampini 
389818f9dSStefano Zampini /* These macros can be moved to petscimpl.h eventually */
489818f9dSStefano Zampini #if defined(PETSC_USE_DEBUG)
589818f9dSStefano Zampini 
689818f9dSStefano Zampini   #define PetscValidLogicalCollectiveIntComm(a, b, c) \
789818f9dSStefano Zampini     do { \
889818f9dSStefano Zampini       PetscInt b1[2], b2[2]; \
99371c9d4SSatish Balay       b1[0] = -b; \
109371c9d4SSatish Balay       b1[1] = b; \
111c2dc1cbSBarry Smith       PetscCall(MPIU_Allreduce(b1, b2, 2, MPIU_INT, MPI_MAX, a)); \
123c633725SBarry Smith       PetscCheck(-b2[0] == b2[1], a, PETSC_ERR_ARG_WRONG, "Int value must be same on all processes, argument # %d", c); \
1389818f9dSStefano Zampini     } while (0)
1489818f9dSStefano Zampini 
1589818f9dSStefano Zampini   #define PetscValidLogicalCollectiveBoolComm(a, b, c) \
1689818f9dSStefano Zampini     do { \
1789818f9dSStefano Zampini       PetscMPIInt b1[2], b2[2]; \
189371c9d4SSatish Balay       b1[0] = -(PetscMPIInt)b; \
199371c9d4SSatish Balay       b1[1] = (PetscMPIInt)b; \
201c2dc1cbSBarry Smith       PetscCall(MPIU_Allreduce(b1, b2, 2, MPI_INT, MPI_MAX, a)); \
213c633725SBarry Smith       PetscCheck(-b2[0] == b2[1], a, PETSC_ERR_ARG_WRONG, "Bool value must be same on all processes, argument # %d", c); \
2289818f9dSStefano Zampini     } while (0)
2389818f9dSStefano Zampini 
2489818f9dSStefano Zampini   #define PetscValidLogicalCollectiveRealComm(a, b, c) \
2589818f9dSStefano Zampini     do { \
2689818f9dSStefano Zampini       PetscReal b1[3], b2[3]; \
279371c9d4SSatish Balay       if (PetscIsNanReal(b)) { \
289371c9d4SSatish Balay         b1[2] = 1; \
299371c9d4SSatish Balay       } else { \
309371c9d4SSatish Balay         b1[2] = 0; \
319371c9d4SSatish Balay       }; \
329371c9d4SSatish Balay       b1[0] = -b; \
339371c9d4SSatish Balay       b1[1] = b; \
34712fec58SPierre Jolivet       PetscCall(MPIU_Allreduce(b1, b2, 3, MPIU_REAL, MPIU_MAX, a)); \
353c633725SBarry Smith       PetscCheck((b2[2] == 1) || PetscEqualReal(-b2[0], b2[1]), a, PETSC_ERR_ARG_WRONG, "Real value must be same on all processes, argument # %d", c); \
3689818f9dSStefano Zampini     } while (0)
3789818f9dSStefano Zampini 
3889818f9dSStefano Zampini #else
3989818f9dSStefano Zampini 
409371c9d4SSatish Balay   #define PetscValidLogicalCollectiveRealComm(a, b, c) \
419371c9d4SSatish Balay     do { \
429371c9d4SSatish Balay     } while (0)
439371c9d4SSatish Balay   #define PetscValidLogicalCollectiveIntComm(a, b, c) \
449371c9d4SSatish Balay     do { \
459371c9d4SSatish Balay     } while (0)
469371c9d4SSatish Balay   #define PetscValidLogicalCollectiveBoolComm(a, b, c) \
479371c9d4SSatish Balay     do { \
489371c9d4SSatish Balay     } while (0)
4989818f9dSStefano Zampini 
5089818f9dSStefano Zampini #endif
5189818f9dSStefano Zampini 
5289818f9dSStefano Zampini struct _n_TSHistory {
5389818f9dSStefano Zampini   MPI_Comm   comm;    /* used for runtime collective checks */
5489818f9dSStefano Zampini   PetscReal *hist;    /* time history */
5589818f9dSStefano Zampini   PetscInt  *hist_id; /* stores the stepid in time history */
56115dd874SBarry Smith   size_t     n;       /* current number of steps registered */
5789818f9dSStefano Zampini   PetscBool  sorted;  /* if the history is sorted in ascending order */
58115dd874SBarry Smith   size_t     c;       /* current capacity of history */
59115dd874SBarry Smith   size_t     s;       /* reallocation size */
6089818f9dSStefano Zampini };
6189818f9dSStefano Zampini 
62d71ae5a4SJacob Faibussowitsch PetscErrorCode TSHistoryGetNumSteps(TSHistory tsh, PetscInt *n)
63d71ae5a4SJacob Faibussowitsch {
6489818f9dSStefano Zampini   PetscFunctionBegin;
654f572ea9SToby Isaac   PetscAssertPointer(n, 2);
6689818f9dSStefano Zampini   *n = tsh->n;
673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6889818f9dSStefano Zampini }
6989818f9dSStefano Zampini 
70d71ae5a4SJacob Faibussowitsch PetscErrorCode TSHistoryUpdate(TSHistory tsh, PetscInt id, PetscReal time)
71d71ae5a4SJacob Faibussowitsch {
7289818f9dSStefano Zampini   PetscFunctionBegin;
7389818f9dSStefano Zampini   if (tsh->n == tsh->c) { /* reallocation */
7489818f9dSStefano Zampini     tsh->c += tsh->s;
759566063dSJacob Faibussowitsch     PetscCall(PetscRealloc(tsh->c * sizeof(*tsh->hist), &tsh->hist));
769566063dSJacob Faibussowitsch     PetscCall(PetscRealloc(tsh->c * sizeof(*tsh->hist_id), &tsh->hist_id));
7789818f9dSStefano Zampini   }
7889818f9dSStefano Zampini   tsh->sorted = (PetscBool)(tsh->sorted && (tsh->n ? time >= tsh->hist[tsh->n - 1] : PETSC_TRUE));
7989818f9dSStefano Zampini #if defined(PETSC_USE_DEBUG)
8089818f9dSStefano Zampini   if (tsh->n) { /* id should be unique */
8189818f9dSStefano Zampini     PetscInt loc, *ids;
8289818f9dSStefano Zampini 
839566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(tsh->n, &ids));
849566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(ids, tsh->hist_id, tsh->n));
859566063dSJacob Faibussowitsch     PetscCall(PetscSortInt(tsh->n, ids));
869566063dSJacob Faibussowitsch     PetscCall(PetscFindInt(id, tsh->n, ids, &loc));
879566063dSJacob Faibussowitsch     PetscCall(PetscFree(ids));
883c633725SBarry Smith     PetscCheck(loc < 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "History id should be unique");
8989818f9dSStefano Zampini   }
9089818f9dSStefano Zampini #endif
9189818f9dSStefano Zampini   tsh->hist[tsh->n]    = time;
9289818f9dSStefano Zampini   tsh->hist_id[tsh->n] = id;
9389818f9dSStefano Zampini   tsh->n += 1;
943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9589818f9dSStefano Zampini }
9689818f9dSStefano Zampini 
97d71ae5a4SJacob Faibussowitsch PetscErrorCode TSHistoryGetTime(TSHistory tsh, PetscBool backward, PetscInt step, PetscReal *t)
98d71ae5a4SJacob Faibussowitsch {
998a32c442SStefano Zampini   PetscFunctionBegin;
1003ba16761SJacob Faibussowitsch   if (!t) PetscFunctionReturn(PETSC_SUCCESS);
1014f572ea9SToby Isaac   PetscAssertPointer(t, 4);
1028a32c442SStefano Zampini   if (!tsh->sorted) {
1039566063dSJacob Faibussowitsch     PetscCall(PetscSortRealWithArrayInt(tsh->n, tsh->hist, tsh->hist_id));
1048a32c442SStefano Zampini     tsh->sorted = PETSC_TRUE;
1058a32c442SStefano Zampini   }
10663a3b9bcSJacob Faibussowitsch   PetscCheck(step >= 0 && step < (PetscInt)tsh->n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Given time step %" PetscInt_FMT " does not match any in history [0,%" PetscInt_FMT "]", step, (PetscInt)tsh->n);
1078a32c442SStefano Zampini   if (!backward) *t = tsh->hist[step];
1088a32c442SStefano Zampini   else *t = tsh->hist[tsh->n - step - 1];
1093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1108a32c442SStefano Zampini }
1118a32c442SStefano Zampini 
112d71ae5a4SJacob Faibussowitsch PetscErrorCode TSHistoryGetTimeStep(TSHistory tsh, PetscBool backward, PetscInt step, PetscReal *dt)
113d71ae5a4SJacob Faibussowitsch {
11489818f9dSStefano Zampini   PetscFunctionBegin;
1153ba16761SJacob Faibussowitsch   if (!dt) PetscFunctionReturn(PETSC_SUCCESS);
1164f572ea9SToby Isaac   PetscAssertPointer(dt, 4);
11789818f9dSStefano Zampini   if (!tsh->sorted) {
1189566063dSJacob Faibussowitsch     PetscCall(PetscSortRealWithArrayInt(tsh->n, tsh->hist, tsh->hist_id));
11989818f9dSStefano Zampini     tsh->sorted = PETSC_TRUE;
12089818f9dSStefano Zampini   }
12163a3b9bcSJacob Faibussowitsch   PetscCheck(step >= 0 && step <= (PetscInt)tsh->n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Given time step %" PetscInt_FMT " does not match any in history [0,%" PetscInt_FMT "]", step, (PetscInt)tsh->n);
122115dd874SBarry Smith   if (!backward) *dt = tsh->hist[PetscMin(step + 1, (PetscInt)tsh->n - 1)] - tsh->hist[PetscMin(step, (PetscInt)tsh->n - 1)];
123115dd874SBarry Smith   else *dt = tsh->hist[PetscMax((PetscInt)tsh->n - step - 1, 0)] - tsh->hist[PetscMax((PetscInt)tsh->n - step - 2, 0)];
1243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12589818f9dSStefano Zampini }
12689818f9dSStefano Zampini 
127d71ae5a4SJacob Faibussowitsch PetscErrorCode TSHistoryGetLocFromTime(TSHistory tsh, PetscReal time, PetscInt *loc)
128d71ae5a4SJacob Faibussowitsch {
12989818f9dSStefano Zampini   PetscFunctionBegin;
1304f572ea9SToby Isaac   PetscAssertPointer(loc, 3);
13189818f9dSStefano Zampini   if (!tsh->sorted) {
1329566063dSJacob Faibussowitsch     PetscCall(PetscSortRealWithArrayInt(tsh->n, tsh->hist, tsh->hist_id));
13389818f9dSStefano Zampini     tsh->sorted = PETSC_TRUE;
13489818f9dSStefano Zampini   }
1359566063dSJacob Faibussowitsch   PetscCall(PetscFindReal(time, tsh->n, tsh->hist, PETSC_SMALL, loc));
1363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13789818f9dSStefano Zampini }
13889818f9dSStefano Zampini 
139d71ae5a4SJacob Faibussowitsch PetscErrorCode TSHistorySetHistory(TSHistory tsh, PetscInt n, PetscReal hist[], PetscInt hist_id[], PetscBool sorted)
140d71ae5a4SJacob Faibussowitsch {
14189818f9dSStefano Zampini   PetscFunctionBegin;
14289818f9dSStefano Zampini   PetscValidLogicalCollectiveIntComm(tsh->comm, n, 2);
14354c59aa7SJacob Faibussowitsch   PetscCheck(n >= 0, tsh->comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot request a negative size for history storage");
1444f572ea9SToby Isaac   if (n) PetscAssertPointer(hist, 3);
1459566063dSJacob Faibussowitsch   PetscCall(PetscFree(tsh->hist));
1469566063dSJacob Faibussowitsch   PetscCall(PetscFree(tsh->hist_id));
147115dd874SBarry Smith   tsh->n = (size_t)n;
148115dd874SBarry Smith   tsh->c = (size_t)n;
1499566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tsh->n, &tsh->hist));
1509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tsh->n, &tsh->hist_id));
1515f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0; i < (PetscInt)tsh->n; i++) {
15289818f9dSStefano Zampini     tsh->hist[i]    = hist[i];
15389818f9dSStefano Zampini     tsh->hist_id[i] = hist_id ? hist_id[i] : i;
15489818f9dSStefano Zampini   }
1559566063dSJacob Faibussowitsch   if (!sorted) PetscCall(PetscSortRealWithArrayInt(tsh->n, tsh->hist, tsh->hist_id));
15689818f9dSStefano Zampini   tsh->sorted = PETSC_TRUE;
1573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15889818f9dSStefano Zampini }
15989818f9dSStefano Zampini 
160d71ae5a4SJacob Faibussowitsch PetscErrorCode TSHistoryGetHistory(TSHistory tsh, PetscInt *n, const PetscReal *hist[], const PetscInt *hist_id[], PetscBool *sorted)
161d71ae5a4SJacob Faibussowitsch {
16289818f9dSStefano Zampini   PetscFunctionBegin;
16389818f9dSStefano Zampini   if (n) *n = tsh->n;
16489818f9dSStefano Zampini   if (hist) *hist = tsh->hist;
16589818f9dSStefano Zampini   if (hist_id) *hist_id = tsh->hist_id;
16689818f9dSStefano Zampini   if (sorted) *sorted = tsh->sorted;
1673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16889818f9dSStefano Zampini }
16989818f9dSStefano Zampini 
170d71ae5a4SJacob Faibussowitsch PetscErrorCode TSHistoryDestroy(TSHistory *tsh)
171d71ae5a4SJacob Faibussowitsch {
17289818f9dSStefano Zampini   PetscFunctionBegin;
1733ba16761SJacob Faibussowitsch   if (!*tsh) PetscFunctionReturn(PETSC_SUCCESS);
1749566063dSJacob Faibussowitsch   PetscCall(PetscFree((*tsh)->hist));
1759566063dSJacob Faibussowitsch   PetscCall(PetscFree((*tsh)->hist_id));
1769566063dSJacob Faibussowitsch   PetscCall(PetscCommDestroy(&((*tsh)->comm)));
177*f4f49eeaSPierre Jolivet   PetscCall(PetscFree(*tsh));
17889818f9dSStefano Zampini   *tsh = NULL;
1793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18089818f9dSStefano Zampini }
18189818f9dSStefano Zampini 
182d71ae5a4SJacob Faibussowitsch PetscErrorCode TSHistoryCreate(MPI_Comm comm, TSHistory *hst)
183d71ae5a4SJacob Faibussowitsch {
18489818f9dSStefano Zampini   TSHistory tsh;
18589818f9dSStefano Zampini 
18689818f9dSStefano Zampini   PetscFunctionBegin;
1874f572ea9SToby Isaac   PetscAssertPointer(hst, 2);
18889818f9dSStefano Zampini   *hst = NULL;
1899566063dSJacob Faibussowitsch   PetscCall(PetscNew(&tsh));
1909566063dSJacob Faibussowitsch   PetscCall(PetscCommDuplicate(comm, &tsh->comm, NULL));
19189818f9dSStefano Zampini 
19289818f9dSStefano Zampini   tsh->c      = 1024; /* capacity */
19389818f9dSStefano Zampini   tsh->s      = 1024; /* reallocation size */
19489818f9dSStefano Zampini   tsh->sorted = PETSC_TRUE;
19589818f9dSStefano Zampini 
1969566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tsh->c, &tsh->hist));
1979566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tsh->c, &tsh->hist_id));
19889818f9dSStefano Zampini   *hst = tsh;
1993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20089818f9dSStefano Zampini }
201