xref: /petsc/src/ts/interface/tshistory.c (revision e69f218b27fb14433472a1cd01a89601f1c93073)
189818f9dSStefano Zampini #include <petsc/private/tshistoryimpl.h>
289818f9dSStefano Zampini 
389818f9dSStefano Zampini struct _n_TSHistory {
489818f9dSStefano Zampini   MPI_Comm   comm;    /* used for runtime collective checks */
589818f9dSStefano Zampini   PetscReal *hist;    /* time history */
689818f9dSStefano Zampini   PetscInt  *hist_id; /* stores the stepid in time history */
76497c311SBarry Smith   PetscCount n;       /* current number of steps registered */
889818f9dSStefano Zampini   PetscBool  sorted;  /* if the history is sorted in ascending order */
96497c311SBarry Smith   PetscCount c;       /* current capacity of history */
106497c311SBarry Smith   PetscCount s;       /* reallocation size */
1189818f9dSStefano Zampini };
1289818f9dSStefano Zampini 
TSHistoryGetNumSteps(TSHistory tsh,PetscInt * n)13d71ae5a4SJacob Faibussowitsch PetscErrorCode TSHistoryGetNumSteps(TSHistory tsh, PetscInt *n)
14d71ae5a4SJacob Faibussowitsch {
1589818f9dSStefano Zampini   PetscFunctionBegin;
164f572ea9SToby Isaac   PetscAssertPointer(n, 2);
176497c311SBarry Smith   PetscCall(PetscIntCast(tsh->n, n));
183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1989818f9dSStefano Zampini }
2089818f9dSStefano Zampini 
TSHistoryUpdate(TSHistory tsh,PetscInt id,PetscReal time)21d71ae5a4SJacob Faibussowitsch PetscErrorCode TSHistoryUpdate(TSHistory tsh, PetscInt id, PetscReal time)
22d71ae5a4SJacob Faibussowitsch {
2389818f9dSStefano Zampini   PetscFunctionBegin;
2489818f9dSStefano Zampini   if (tsh->n == tsh->c) { /* reallocation */
2589818f9dSStefano Zampini     tsh->c += tsh->s;
269566063dSJacob Faibussowitsch     PetscCall(PetscRealloc(tsh->c * sizeof(*tsh->hist), &tsh->hist));
279566063dSJacob Faibussowitsch     PetscCall(PetscRealloc(tsh->c * sizeof(*tsh->hist_id), &tsh->hist_id));
2889818f9dSStefano Zampini   }
29d070fe2eSStefano Zampini   tsh->sorted = (PetscBool)(tsh->sorted && (tsh->n ? (PetscBool)(time >= tsh->hist[tsh->n - 1]) : PETSC_TRUE));
3089818f9dSStefano Zampini #if defined(PETSC_USE_DEBUG)
3189818f9dSStefano Zampini   if (tsh->n) { /* id should be unique */
3289818f9dSStefano Zampini     PetscInt loc, *ids;
3389818f9dSStefano Zampini 
349566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(tsh->n, &ids));
359566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(ids, tsh->hist_id, tsh->n));
369566063dSJacob Faibussowitsch     PetscCall(PetscSortInt(tsh->n, ids));
379566063dSJacob Faibussowitsch     PetscCall(PetscFindInt(id, tsh->n, ids, &loc));
389566063dSJacob Faibussowitsch     PetscCall(PetscFree(ids));
393c633725SBarry Smith     PetscCheck(loc < 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "History id should be unique");
4089818f9dSStefano Zampini   }
4189818f9dSStefano Zampini #endif
4289818f9dSStefano Zampini   tsh->hist[tsh->n]    = time;
4389818f9dSStefano Zampini   tsh->hist_id[tsh->n] = id;
4489818f9dSStefano Zampini   tsh->n += 1;
453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4689818f9dSStefano Zampini }
4789818f9dSStefano Zampini 
TSHistoryGetTime(TSHistory tsh,PetscBool backward,PetscInt step,PetscReal * t)48d71ae5a4SJacob Faibussowitsch PetscErrorCode TSHistoryGetTime(TSHistory tsh, PetscBool backward, PetscInt step, PetscReal *t)
49d71ae5a4SJacob Faibussowitsch {
508a32c442SStefano Zampini   PetscFunctionBegin;
513ba16761SJacob Faibussowitsch   if (!t) PetscFunctionReturn(PETSC_SUCCESS);
524f572ea9SToby Isaac   PetscAssertPointer(t, 4);
538a32c442SStefano Zampini   if (!tsh->sorted) {
549566063dSJacob Faibussowitsch     PetscCall(PetscSortRealWithArrayInt(tsh->n, tsh->hist, tsh->hist_id));
558a32c442SStefano Zampini     tsh->sorted = PETSC_TRUE;
568a32c442SStefano Zampini   }
57*835f2295SStefano Zampini   PetscCheck(step >= 0 && step < tsh->n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Given time step %" PetscInt_FMT " does not match any in history [0,%" PetscCount_FMT "]", step, tsh->n);
588a32c442SStefano Zampini   if (!backward) *t = tsh->hist[step];
598a32c442SStefano Zampini   else *t = tsh->hist[tsh->n - step - 1];
603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
618a32c442SStefano Zampini }
628a32c442SStefano Zampini 
TSHistoryGetTimeStep(TSHistory tsh,PetscBool backward,PetscInt step,PetscReal * dt)63d71ae5a4SJacob Faibussowitsch PetscErrorCode TSHistoryGetTimeStep(TSHistory tsh, PetscBool backward, PetscInt step, PetscReal *dt)
64d71ae5a4SJacob Faibussowitsch {
6589818f9dSStefano Zampini   PetscFunctionBegin;
663ba16761SJacob Faibussowitsch   if (!dt) PetscFunctionReturn(PETSC_SUCCESS);
674f572ea9SToby Isaac   PetscAssertPointer(dt, 4);
6889818f9dSStefano Zampini   if (!tsh->sorted) {
699566063dSJacob Faibussowitsch     PetscCall(PetscSortRealWithArrayInt(tsh->n, tsh->hist, tsh->hist_id));
7089818f9dSStefano Zampini     tsh->sorted = PETSC_TRUE;
7189818f9dSStefano Zampini   }
72*835f2295SStefano Zampini   PetscCheck(step >= 0 && step <= tsh->n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Given time step %" PetscInt_FMT " does not match any in history [0,%" PetscCount_FMT "]", step, tsh->n);
73*835f2295SStefano Zampini   if (!backward) *dt = tsh->hist[PetscMin(step + 1, tsh->n - 1)] - tsh->hist[PetscMin(step, tsh->n - 1)];
74*835f2295SStefano Zampini   else *dt = tsh->hist[PetscMax(tsh->n - step - 1, 0)] - tsh->hist[PetscMax(tsh->n - step - 2, 0)];
753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7689818f9dSStefano Zampini }
7789818f9dSStefano Zampini 
TSHistoryGetLocFromTime(TSHistory tsh,PetscReal time,PetscInt * loc)78d71ae5a4SJacob Faibussowitsch PetscErrorCode TSHistoryGetLocFromTime(TSHistory tsh, PetscReal time, PetscInt *loc)
79d71ae5a4SJacob Faibussowitsch {
8089818f9dSStefano Zampini   PetscFunctionBegin;
814f572ea9SToby Isaac   PetscAssertPointer(loc, 3);
8289818f9dSStefano Zampini   if (!tsh->sorted) {
839566063dSJacob Faibussowitsch     PetscCall(PetscSortRealWithArrayInt(tsh->n, tsh->hist, tsh->hist_id));
8489818f9dSStefano Zampini     tsh->sorted = PETSC_TRUE;
8589818f9dSStefano Zampini   }
869566063dSJacob Faibussowitsch   PetscCall(PetscFindReal(time, tsh->n, tsh->hist, PETSC_SMALL, loc));
873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8889818f9dSStefano Zampini }
8989818f9dSStefano Zampini 
TSHistorySetHistory(TSHistory tsh,PetscInt n,PetscReal hist[],PetscInt hist_id[],PetscBool sorted)90d71ae5a4SJacob Faibussowitsch PetscErrorCode TSHistorySetHistory(TSHistory tsh, PetscInt n, PetscReal hist[], PetscInt hist_id[], PetscBool sorted)
91d71ae5a4SJacob Faibussowitsch {
9289818f9dSStefano Zampini   PetscFunctionBegin;
9389818f9dSStefano Zampini   PetscValidLogicalCollectiveIntComm(tsh->comm, n, 2);
9454c59aa7SJacob Faibussowitsch   PetscCheck(n >= 0, tsh->comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot request a negative size for history storage");
954f572ea9SToby Isaac   if (n) PetscAssertPointer(hist, 3);
969566063dSJacob Faibussowitsch   PetscCall(PetscFree(tsh->hist));
979566063dSJacob Faibussowitsch   PetscCall(PetscFree(tsh->hist_id));
98115dd874SBarry Smith   tsh->n = (size_t)n;
99115dd874SBarry Smith   tsh->c = (size_t)n;
1009566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tsh->n, &tsh->hist));
1019566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tsh->n, &tsh->hist_id));
102*835f2295SStefano Zampini   for (PetscInt i = 0; i < n; i++) {
10389818f9dSStefano Zampini     tsh->hist[i]    = hist[i];
10489818f9dSStefano Zampini     tsh->hist_id[i] = hist_id ? hist_id[i] : i;
10589818f9dSStefano Zampini   }
1069566063dSJacob Faibussowitsch   if (!sorted) PetscCall(PetscSortRealWithArrayInt(tsh->n, tsh->hist, tsh->hist_id));
10789818f9dSStefano Zampini   tsh->sorted = PETSC_TRUE;
1083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10989818f9dSStefano Zampini }
11089818f9dSStefano Zampini 
TSHistoryGetHistory(TSHistory tsh,PetscInt * n,const PetscReal * hist[],const PetscInt * hist_id[],PetscBool * sorted)111d71ae5a4SJacob Faibussowitsch PetscErrorCode TSHistoryGetHistory(TSHistory tsh, PetscInt *n, const PetscReal *hist[], const PetscInt *hist_id[], PetscBool *sorted)
112d71ae5a4SJacob Faibussowitsch {
11389818f9dSStefano Zampini   PetscFunctionBegin;
1146497c311SBarry Smith   if (n) PetscCall(PetscIntCast(tsh->n, n));
11589818f9dSStefano Zampini   if (hist) *hist = tsh->hist;
11689818f9dSStefano Zampini   if (hist_id) *hist_id = tsh->hist_id;
11789818f9dSStefano Zampini   if (sorted) *sorted = tsh->sorted;
1183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11989818f9dSStefano Zampini }
12089818f9dSStefano Zampini 
TSHistoryDestroy(TSHistory * tsh)121d71ae5a4SJacob Faibussowitsch PetscErrorCode TSHistoryDestroy(TSHistory *tsh)
122d71ae5a4SJacob Faibussowitsch {
12389818f9dSStefano Zampini   PetscFunctionBegin;
1243ba16761SJacob Faibussowitsch   if (!*tsh) PetscFunctionReturn(PETSC_SUCCESS);
1259566063dSJacob Faibussowitsch   PetscCall(PetscFree((*tsh)->hist));
1269566063dSJacob Faibussowitsch   PetscCall(PetscFree((*tsh)->hist_id));
12757508eceSPierre Jolivet   PetscCall(PetscCommDestroy(&(*tsh)->comm));
128f4f49eeaSPierre Jolivet   PetscCall(PetscFree(*tsh));
12989818f9dSStefano Zampini   *tsh = NULL;
1303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13189818f9dSStefano Zampini }
13289818f9dSStefano Zampini 
TSHistoryCreate(MPI_Comm comm,TSHistory * hst)133d71ae5a4SJacob Faibussowitsch PetscErrorCode TSHistoryCreate(MPI_Comm comm, TSHistory *hst)
134d71ae5a4SJacob Faibussowitsch {
13589818f9dSStefano Zampini   TSHistory tsh;
13689818f9dSStefano Zampini 
13789818f9dSStefano Zampini   PetscFunctionBegin;
1384f572ea9SToby Isaac   PetscAssertPointer(hst, 2);
1399566063dSJacob Faibussowitsch   PetscCall(PetscNew(&tsh));
1409566063dSJacob Faibussowitsch   PetscCall(PetscCommDuplicate(comm, &tsh->comm, NULL));
14189818f9dSStefano Zampini 
14289818f9dSStefano Zampini   tsh->c      = 1024; /* capacity */
14389818f9dSStefano Zampini   tsh->s      = 1024; /* reallocation size */
14489818f9dSStefano Zampini   tsh->sorted = PETSC_TRUE;
14589818f9dSStefano Zampini 
1469566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tsh->c, &tsh->hist));
1479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tsh->c, &tsh->hist_id));
14889818f9dSStefano Zampini   *hst = tsh;
1493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15089818f9dSStefano Zampini }
151