xref: /petsc/src/sys/logging/plog.c (revision 0bc1c2e8e3dc4082cdfa268312e338118df56085)
15c6c1daeSBarry Smith /*
25c6c1daeSBarry Smith       PETSc code to log object creation and destruction and PETSc events.
35c6c1daeSBarry Smith 
45c6c1daeSBarry Smith       This provides the public API used by the rest of PETSc and by users.
55c6c1daeSBarry Smith 
65c6c1daeSBarry Smith       These routines use a private API that is not used elsewhere in PETSc and is not
75c6c1daeSBarry Smith       accessible to users. The private API is defined in logimpl.h and the utils directory.
85c6c1daeSBarry Smith 
9b665b14eSToby Isaac       ***
10b665b14eSToby Isaac 
11b665b14eSToby Isaac       This file, and only this file, is for functions that interact with the global logging state
125c6c1daeSBarry Smith */
13af0996ceSBarry Smith #include <petsc/private/logimpl.h> /*I    "petscsys.h"   I*/
1453e0a2f3SToby Isaac #include <petsc/private/loghandlerimpl.h>
155c6c1daeSBarry Smith #include <petsctime.h>
16665c2dedSJed Brown #include <petscviewer.h>
178fe3844cSJunchao Zhang #include <petscdevice.h>
188fe3844cSJunchao Zhang #include <petsc/private/deviceimpl.h>
195c6c1daeSBarry Smith 
20c708d6e3SStefano Zampini #if defined(PETSC_HAVE_THREADSAFETY)
21c708d6e3SStefano Zampini 
22c708d6e3SStefano Zampini PetscInt           petsc_log_gid = -1; /* Global threadId counter */
23c708d6e3SStefano Zampini PETSC_TLS PetscInt petsc_log_tid = -1; /* Local threadId */
24c708d6e3SStefano Zampini 
25c708d6e3SStefano Zampini /* shared variables */
26c708d6e3SStefano Zampini PetscSpinlock PetscLogSpinLock;
27c708d6e3SStefano Zampini 
282611ad71SToby Isaac PetscInt PetscLogGetTid(void)
292611ad71SToby Isaac {
302611ad71SToby Isaac   if (petsc_log_tid < 0) {
312611ad71SToby Isaac     PetscCall(PetscSpinlockLock(&PetscLogSpinLock));
322611ad71SToby Isaac     petsc_log_tid = ++petsc_log_gid;
332611ad71SToby Isaac     PetscCall(PetscSpinlockUnlock(&PetscLogSpinLock));
342611ad71SToby Isaac   }
352611ad71SToby Isaac   return petsc_log_tid;
362611ad71SToby Isaac }
372611ad71SToby Isaac 
38c708d6e3SStefano Zampini #endif
39c708d6e3SStefano Zampini 
405c6c1daeSBarry Smith /* Global counters */
415c6c1daeSBarry Smith PetscLogDouble petsc_BaseTime        = 0.0;
425c6c1daeSBarry Smith PetscLogDouble petsc_TotalFlops      = 0.0; /* The number of flops */
435c6c1daeSBarry Smith PetscLogDouble petsc_send_ct         = 0.0; /* The number of sends */
445c6c1daeSBarry Smith PetscLogDouble petsc_recv_ct         = 0.0; /* The number of receives */
455c6c1daeSBarry Smith PetscLogDouble petsc_send_len        = 0.0; /* The total length of all sent messages */
465c6c1daeSBarry Smith PetscLogDouble petsc_recv_len        = 0.0; /* The total length of all received messages */
475c6c1daeSBarry Smith PetscLogDouble petsc_isend_ct        = 0.0; /* The number of immediate sends */
485c6c1daeSBarry Smith PetscLogDouble petsc_irecv_ct        = 0.0; /* The number of immediate receives */
495c6c1daeSBarry Smith PetscLogDouble petsc_isend_len       = 0.0; /* The total length of all immediate send messages */
505c6c1daeSBarry Smith PetscLogDouble petsc_irecv_len       = 0.0; /* The total length of all immediate receive messages */
515c6c1daeSBarry Smith PetscLogDouble petsc_wait_ct         = 0.0; /* The number of waits */
525c6c1daeSBarry Smith PetscLogDouble petsc_wait_any_ct     = 0.0; /* The number of anywaits */
535c6c1daeSBarry Smith PetscLogDouble petsc_wait_all_ct     = 0.0; /* The number of waitalls */
545c6c1daeSBarry Smith PetscLogDouble petsc_sum_of_waits_ct = 0.0; /* The total number of waits */
555c6c1daeSBarry Smith PetscLogDouble petsc_allreduce_ct    = 0.0; /* The number of reductions */
565c6c1daeSBarry Smith PetscLogDouble petsc_gather_ct       = 0.0; /* The number of gathers and gathervs */
575c6c1daeSBarry Smith PetscLogDouble petsc_scatter_ct      = 0.0; /* The number of scatters and scattervs */
58c708d6e3SStefano Zampini 
59c708d6e3SStefano Zampini /* Thread Local storage */
60c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_TotalFlops_th      = 0.0;
61c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_send_ct_th         = 0.0;
62c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_recv_ct_th         = 0.0;
63c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_send_len_th        = 0.0;
64c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_recv_len_th        = 0.0;
65c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_isend_ct_th        = 0.0;
66c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_irecv_ct_th        = 0.0;
67c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_isend_len_th       = 0.0;
68c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_irecv_len_th       = 0.0;
69c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_wait_ct_th         = 0.0;
70c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_wait_any_ct_th     = 0.0;
71c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_wait_all_ct_th     = 0.0;
72c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_sum_of_waits_ct_th = 0.0;
73c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_allreduce_ct_th    = 0.0;
74c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_gather_ct_th       = 0.0;
75c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_scatter_ct_th      = 0.0;
76c708d6e3SStefano Zampini 
77bec0b493Shannah_mairs PetscLogDouble petsc_ctog_ct        = 0.0; /* The total number of CPU to GPU copies */
78bec0b493Shannah_mairs PetscLogDouble petsc_gtoc_ct        = 0.0; /* The total number of GPU to CPU copies */
79bec0b493Shannah_mairs PetscLogDouble petsc_ctog_sz        = 0.0; /* The total size of CPU to GPU copies */
80bec0b493Shannah_mairs PetscLogDouble petsc_gtoc_sz        = 0.0; /* The total size of GPU to CPU copies */
8145c4b7c1SBarry Smith PetscLogDouble petsc_ctog_ct_scalar = 0.0; /* The total number of CPU to GPU copies */
8245c4b7c1SBarry Smith PetscLogDouble petsc_gtoc_ct_scalar = 0.0; /* The total number of GPU to CPU copies */
8345c4b7c1SBarry Smith PetscLogDouble petsc_ctog_sz_scalar = 0.0; /* The total size of CPU to GPU copies */
8445c4b7c1SBarry Smith PetscLogDouble petsc_gtoc_sz_scalar = 0.0; /* The total size of GPU to CPU copies */
85958c4211Shannah_mairs PetscLogDouble petsc_gflops         = 0.0; /* The flops done on a GPU */
86958c4211Shannah_mairs PetscLogDouble petsc_gtime          = 0.0; /* The time spent on a GPU */
87c708d6e3SStefano Zampini 
88c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_ctog_ct_th        = 0.0;
89c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_gtoc_ct_th        = 0.0;
90c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_ctog_sz_th        = 0.0;
91c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_gtoc_sz_th        = 0.0;
92c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_ctog_ct_scalar_th = 0.0;
93c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_gtoc_ct_scalar_th = 0.0;
94c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_ctog_sz_scalar_th = 0.0;
95c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_gtoc_sz_scalar_th = 0.0;
96c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_gflops_th         = 0.0;
97c708d6e3SStefano Zampini PETSC_TLS PetscLogDouble petsc_gtime_th          = 0.0;
982611ad71SToby Isaac 
992611ad71SToby Isaac PetscBool PetscLogMemory = PETSC_FALSE;
1002611ad71SToby Isaac PetscBool PetscLogSyncOn = PETSC_FALSE;
1012611ad71SToby Isaac 
1022611ad71SToby Isaac PetscBool PetscLogGpuTimeFlag = PETSC_FALSE;
1032611ad71SToby Isaac 
104*0bc1c2e8SToby Isaac PetscInt PetscLogNumViewersCreated   = 0;
105*0bc1c2e8SToby Isaac PetscInt PetscLogNumViewersDestroyed = 0;
106*0bc1c2e8SToby Isaac 
1072611ad71SToby Isaac PetscLogState petsc_log_state = NULL;
1082611ad71SToby Isaac 
10910dd146fSPierre Jolivet #define PETSC_LOG_HANDLER_HOT_BLANK {NULL, NULL, NULL, NULL, NULL, NULL}
1102611ad71SToby Isaac 
1112611ad71SToby Isaac PetscLogHandlerHot PetscLogHandlers[PETSC_LOG_HANDLER_MAX] = {
1122611ad71SToby Isaac   PETSC_LOG_HANDLER_HOT_BLANK,
1132611ad71SToby Isaac   PETSC_LOG_HANDLER_HOT_BLANK,
1142611ad71SToby Isaac   PETSC_LOG_HANDLER_HOT_BLANK,
1152611ad71SToby Isaac   PETSC_LOG_HANDLER_HOT_BLANK,
1162611ad71SToby Isaac };
1172611ad71SToby Isaac 
1182611ad71SToby Isaac #undef PETSC_LOG_HANDLERS_HOT_BLANK
1192611ad71SToby Isaac 
1202611ad71SToby Isaac #if defined(PETSC_USE_LOG)
1212611ad71SToby Isaac   #include <../src/sys/logging/handler/impls/default/logdefault.h>
122c708d6e3SStefano Zampini 
123c708d6e3SStefano Zampini   #if defined(PETSC_HAVE_THREADSAFETY)
124c708d6e3SStefano Zampini PetscErrorCode PetscAddLogDouble(PetscLogDouble *tot, PetscLogDouble *tot_th, PetscLogDouble tmp)
125c708d6e3SStefano Zampini {
126c708d6e3SStefano Zampini   *tot_th += tmp;
1273ba16761SJacob Faibussowitsch   PetscCall(PetscSpinlockLock(&PetscLogSpinLock));
128c708d6e3SStefano Zampini   *tot += tmp;
1293ba16761SJacob Faibussowitsch   PetscCall(PetscSpinlockUnlock(&PetscLogSpinLock));
1303ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
131c708d6e3SStefano Zampini }
132c708d6e3SStefano Zampini 
133c708d6e3SStefano Zampini PetscErrorCode PetscAddLogDoubleCnt(PetscLogDouble *cnt, PetscLogDouble *tot, PetscLogDouble *cnt_th, PetscLogDouble *tot_th, PetscLogDouble tmp)
134c708d6e3SStefano Zampini {
135c708d6e3SStefano Zampini   *cnt_th = *cnt_th + 1;
136c708d6e3SStefano Zampini   *tot_th += tmp;
1373ba16761SJacob Faibussowitsch   PetscCall(PetscSpinlockLock(&PetscLogSpinLock));
13857508eceSPierre Jolivet   *tot += (PetscLogDouble)tmp;
139c708d6e3SStefano Zampini   *cnt += *cnt + 1;
1403ba16761SJacob Faibussowitsch   PetscCall(PetscSpinlockUnlock(&PetscLogSpinLock));
1413ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
142c708d6e3SStefano Zampini }
143c708d6e3SStefano Zampini 
144bec0b493Shannah_mairs   #endif
1455c6c1daeSBarry Smith 
14653e0a2f3SToby Isaac static PetscErrorCode PetscLogTryGetHandler(PetscLogHandlerType type, PetscLogHandler *handler)
14753e0a2f3SToby Isaac {
14853e0a2f3SToby Isaac   PetscFunctionBegin;
14953e0a2f3SToby Isaac   PetscAssertPointer(handler, 2);
15053e0a2f3SToby Isaac   *handler = NULL;
15153e0a2f3SToby Isaac   for (int i = 0; i < PETSC_LOG_HANDLER_MAX; i++) {
15253e0a2f3SToby Isaac     PetscLogHandler h = PetscLogHandlers[i].handler;
15353e0a2f3SToby Isaac     if (h) {
15453e0a2f3SToby Isaac       PetscBool match;
15553e0a2f3SToby Isaac 
15653e0a2f3SToby Isaac       PetscCall(PetscObjectTypeCompare((PetscObject)h, type, &match));
15753e0a2f3SToby Isaac       if (match) {
15853e0a2f3SToby Isaac         *handler = PetscLogHandlers[i].handler;
15953e0a2f3SToby Isaac         PetscFunctionReturn(PETSC_SUCCESS);
16053e0a2f3SToby Isaac       }
16153e0a2f3SToby Isaac     }
16253e0a2f3SToby Isaac   }
16353e0a2f3SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
16453e0a2f3SToby Isaac }
16553e0a2f3SToby Isaac 
16653e0a2f3SToby Isaac /*@
167b665b14eSToby Isaac   PetscLogGetDefaultHandler - Get the default log handler if it is running.
168b665b14eSToby Isaac 
169b665b14eSToby Isaac   Not collective
170b665b14eSToby Isaac 
171b665b14eSToby Isaac   Output Parameter:
172b665b14eSToby Isaac . handler - the default `PetscLogHandler`, or `NULL` if it is not running.
173b665b14eSToby Isaac 
174b665b14eSToby Isaac   Level: developer
175b665b14eSToby Isaac 
176b665b14eSToby Isaac   Notes:
177b665b14eSToby Isaac   The default handler is started with `PetscLogDefaultBegin()`,
178b665b14eSToby Isaac   if the options flags `-log_all` or `-log_view` is given without arguments,
179b665b14eSToby Isaac   or for `-log_view :output:format` if `format` is not `ascii_xml` or `ascii_flamegraph`.
180b665b14eSToby Isaac 
181b665b14eSToby Isaac .seealso: [](ch_profiling)
182b665b14eSToby Isaac @*/
183b665b14eSToby Isaac PetscErrorCode PetscLogGetDefaultHandler(PetscLogHandler *handler)
184b665b14eSToby Isaac {
185b665b14eSToby Isaac   PetscFunctionBegin;
186294de794SToby Isaac   PetscCall(PetscLogTryGetHandler(PETSCLOGHANDLERDEFAULT, handler));
187b665b14eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
188b665b14eSToby Isaac }
189b665b14eSToby Isaac 
190b665b14eSToby Isaac static PetscErrorCode PetscLogGetHandler(PetscLogHandlerType type, PetscLogHandler *handler)
191b665b14eSToby Isaac {
192b665b14eSToby Isaac   PetscFunctionBegin;
193b665b14eSToby Isaac   PetscAssertPointer(handler, 2);
194b665b14eSToby Isaac   PetscCall(PetscLogTryGetHandler(type, handler));
1951943a749SToby Isaac   PetscCheck(*handler != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "A PetscLogHandler of type %s has not been started.", type);
196b665b14eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
197b665b14eSToby Isaac }
198b665b14eSToby Isaac 
199b665b14eSToby Isaac /*@
20053e0a2f3SToby Isaac   PetscLogGetState - Get the `PetscLogState` for PETSc's global logging, used
20153e0a2f3SToby Isaac   by all default log handlers (`PetscLogDefaultBegin()`,
20253e0a2f3SToby Isaac   `PetscLogNestedBegin()`, `PetscLogTraceBegin()`, `PetscLogMPEBegin()`,
20353e0a2f3SToby Isaac   `PetscLogPerfstubsBegin()`).
20453e0a2f3SToby Isaac 
20553e0a2f3SToby Isaac   Collective on `PETSC_COMM_WORLD`
20653e0a2f3SToby Isaac 
20753e0a2f3SToby Isaac   Output Parameter:
208b665b14eSToby Isaac . state - The `PetscLogState` changed by registrations (such as
209b665b14eSToby Isaac           `PetscLogEventRegister()`) and actions (such as `PetscLogEventBegin()` or
210b8004f34SBarry Smith           `PetscLogStagePush()`), or `NULL` if logging is not active
21153e0a2f3SToby Isaac 
21253e0a2f3SToby Isaac   Level: developer
21353e0a2f3SToby Isaac 
21453e0a2f3SToby Isaac .seealso: [](ch_profiling), `PetscLogState`
21553e0a2f3SToby Isaac @*/
21653e0a2f3SToby Isaac PetscErrorCode PetscLogGetState(PetscLogState *state)
21753e0a2f3SToby Isaac {
21853e0a2f3SToby Isaac   PetscFunctionBegin;
21953e0a2f3SToby Isaac   PetscAssertPointer(state, 1);
22053e0a2f3SToby Isaac   *state = petsc_log_state;
22153e0a2f3SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
22253e0a2f3SToby Isaac }
22353e0a2f3SToby Isaac 
22453e0a2f3SToby Isaac static PetscErrorCode PetscLogHandlerCopyToHot(PetscLogHandler h, PetscLogHandlerHot *hot)
22553e0a2f3SToby Isaac {
22653e0a2f3SToby Isaac   PetscFunctionBegin;
22753e0a2f3SToby Isaac   hot->handler       = h;
22853e0a2f3SToby Isaac   hot->eventBegin    = h->ops->eventbegin;
22953e0a2f3SToby Isaac   hot->eventEnd      = h->ops->eventend;
23053e0a2f3SToby Isaac   hot->eventSync     = h->ops->eventsync;
23153e0a2f3SToby Isaac   hot->objectCreate  = h->ops->objectcreate;
23253e0a2f3SToby Isaac   hot->objectDestroy = h->ops->objectdestroy;
23353e0a2f3SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
23453e0a2f3SToby Isaac }
23553e0a2f3SToby Isaac 
23653e0a2f3SToby Isaac /*@
23753e0a2f3SToby Isaac   PetscLogHandlerStart - Connect a log handler to PETSc's global logging stream and state.
23853e0a2f3SToby Isaac 
23953e0a2f3SToby Isaac   Logically collective
24053e0a2f3SToby Isaac 
24153e0a2f3SToby Isaac   Input Parameters:
24253e0a2f3SToby Isaac . h - a `PetscLogHandler`
24353e0a2f3SToby Isaac 
24453e0a2f3SToby Isaac   Level: developer
24553e0a2f3SToby Isaac 
24653e0a2f3SToby Isaac   Notes:
24753e0a2f3SToby Isaac   Users should only need this if they create their own log handlers: handlers that are started
24853e0a2f3SToby Isaac   from the command line (such as `-log_view` and `-log_trace`) or from a function like
24953e0a2f3SToby Isaac   `PetscLogNestedBegin()` will automatically be started.
25053e0a2f3SToby Isaac 
25153e0a2f3SToby Isaac   There is a limit of `PESC_LOG_HANDLER_MAX` handlers that can be active at one time.
25253e0a2f3SToby Isaac 
25353e0a2f3SToby Isaac   To disconnect a handler from the global stream call `PetscLogHandlerStop()`.
25453e0a2f3SToby Isaac 
25553e0a2f3SToby Isaac   When a log handler is started, stages that have already been pushed with `PetscLogStagePush()`,
25653e0a2f3SToby Isaac   will be pushed for the new log handler, but it will not be informed of any events that are
25753e0a2f3SToby Isaac   in progress.  It is recommended to start any user-defined log handlers immediately following
258b8004f34SBarry Smith   `PetscInitialize()`  before any user-defined stages are pushed.
25953e0a2f3SToby Isaac 
260b8004f34SBarry Smith .seealso: [](ch_profiling), `PetscLogHandler`, `PetscLogState`, `PetscLogHandlerStop()`, `PetscInitialize()`
26153e0a2f3SToby Isaac @*/
26253e0a2f3SToby Isaac PetscErrorCode PetscLogHandlerStart(PetscLogHandler h)
26353e0a2f3SToby Isaac {
26453e0a2f3SToby Isaac   PetscFunctionBegin;
26553e0a2f3SToby Isaac   for (PetscInt i = 0; i < PETSC_LOG_HANDLER_MAX; i++) {
26653e0a2f3SToby Isaac     if (PetscLogHandlers[i].handler == h) PetscFunctionReturn(PETSC_SUCCESS);
26753e0a2f3SToby Isaac   }
26853e0a2f3SToby Isaac   for (PetscInt i = 0; i < PETSC_LOG_HANDLER_MAX; i++) {
26953e0a2f3SToby Isaac     if (PetscLogHandlers[i].handler == NULL) {
27053e0a2f3SToby Isaac       PetscCall(PetscObjectReference((PetscObject)h));
27153e0a2f3SToby Isaac       PetscCall(PetscLogHandlerCopyToHot(h, &PetscLogHandlers[i]));
27253e0a2f3SToby Isaac       if (petsc_log_state) {
27353e0a2f3SToby Isaac         PetscLogStage stack_height;
27453e0a2f3SToby Isaac         PetscIntStack orig_stack, temp_stack;
27553e0a2f3SToby Isaac 
27653e0a2f3SToby Isaac         PetscCall(PetscLogHandlerSetState(h, petsc_log_state));
27753e0a2f3SToby Isaac         stack_height = petsc_log_state->stage_stack->top + 1;
27853e0a2f3SToby Isaac         PetscCall(PetscIntStackCreate(&temp_stack));
27953e0a2f3SToby Isaac         orig_stack                     = petsc_log_state->stage_stack;
28053e0a2f3SToby Isaac         petsc_log_state->stage_stack   = temp_stack;
28153e0a2f3SToby Isaac         petsc_log_state->current_stage = -1;
28253e0a2f3SToby Isaac         for (int s = 0; s < stack_height; s++) {
283835f2295SStefano Zampini           PetscLogStage stage = orig_stack->stack[s];
28453e0a2f3SToby Isaac           PetscCall(PetscLogHandlerStagePush(h, stage));
28553e0a2f3SToby Isaac           PetscCall(PetscIntStackPush(temp_stack, stage));
28653e0a2f3SToby Isaac           petsc_log_state->current_stage = stage;
28753e0a2f3SToby Isaac         }
28853e0a2f3SToby Isaac         PetscCall(PetscIntStackDestroy(temp_stack));
28953e0a2f3SToby Isaac         petsc_log_state->stage_stack = orig_stack;
29053e0a2f3SToby Isaac       }
29153e0a2f3SToby Isaac       PetscFunctionReturn(PETSC_SUCCESS);
29253e0a2f3SToby Isaac     }
29353e0a2f3SToby Isaac   }
29453e0a2f3SToby Isaac   SETERRQ(PetscObjectComm((PetscObject)h), PETSC_ERR_ARG_WRONGSTATE, "%d log handlers already started, cannot start another", PETSC_LOG_HANDLER_MAX);
29553e0a2f3SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
29653e0a2f3SToby Isaac }
29753e0a2f3SToby Isaac 
29853e0a2f3SToby Isaac /*@
29953e0a2f3SToby Isaac   PetscLogHandlerStop - Disconnect a log handler from PETSc's global logging stream.
30053e0a2f3SToby Isaac 
30153e0a2f3SToby Isaac   Logically collective
30253e0a2f3SToby Isaac 
30353e0a2f3SToby Isaac   Input Parameters:
30453e0a2f3SToby Isaac . h - a `PetscLogHandler`
30553e0a2f3SToby Isaac 
30653e0a2f3SToby Isaac   Level: developer
30753e0a2f3SToby Isaac 
30853e0a2f3SToby Isaac   Note:
30953e0a2f3SToby Isaac   After `PetscLogHandlerStop()`, the handler can still access the global logging state
31053e0a2f3SToby Isaac   with `PetscLogHandlerGetState()`, so that it can access the registry when post-processing
31153e0a2f3SToby Isaac   (for instance, in `PetscLogHandlerView()`),
31253e0a2f3SToby Isaac 
31353e0a2f3SToby Isaac   When a log handler is stopped, the remaining stages will be popped before it is
31453e0a2f3SToby Isaac   disconnected from the log stream.
31553e0a2f3SToby Isaac 
31653e0a2f3SToby Isaac .seealso: [](ch_profiling), `PetscLogHandler`, `PetscLogState`, `PetscLogHandlerStart()`
31753e0a2f3SToby Isaac @*/
31853e0a2f3SToby Isaac PetscErrorCode PetscLogHandlerStop(PetscLogHandler h)
31953e0a2f3SToby Isaac {
32053e0a2f3SToby Isaac   PetscFunctionBegin;
32153e0a2f3SToby Isaac   for (PetscInt i = 0; i < PETSC_LOG_HANDLER_MAX; i++) {
32253e0a2f3SToby Isaac     if (PetscLogHandlers[i].handler == h) {
32353e0a2f3SToby Isaac       if (petsc_log_state) {
32453e0a2f3SToby Isaac         PetscLogState state;
32553e0a2f3SToby Isaac         PetscLogStage stack_height;
32653e0a2f3SToby Isaac         PetscIntStack orig_stack, temp_stack;
32753e0a2f3SToby Isaac 
32853e0a2f3SToby Isaac         PetscCall(PetscLogHandlerGetState(h, &state));
32953e0a2f3SToby Isaac         PetscCheck(state == petsc_log_state, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONGSTATE, "Called PetscLogHandlerStop() for a PetscLogHander that was not started.");
33053e0a2f3SToby Isaac         stack_height = petsc_log_state->stage_stack->top + 1;
33153e0a2f3SToby Isaac         PetscCall(PetscIntStackCreate(&temp_stack));
33253e0a2f3SToby Isaac         orig_stack                   = petsc_log_state->stage_stack;
33353e0a2f3SToby Isaac         petsc_log_state->stage_stack = temp_stack;
33453e0a2f3SToby Isaac         for (int s = 0; s < stack_height; s++) {
335835f2295SStefano Zampini           PetscLogStage stage = orig_stack->stack[s];
33653e0a2f3SToby Isaac 
33753e0a2f3SToby Isaac           PetscCall(PetscIntStackPush(temp_stack, stage));
33853e0a2f3SToby Isaac         }
33953e0a2f3SToby Isaac         for (int s = 0; s < stack_height; s++) {
34053e0a2f3SToby Isaac           PetscLogStage stage;
34153e0a2f3SToby Isaac           PetscBool     empty;
34253e0a2f3SToby Isaac 
34353e0a2f3SToby Isaac           PetscCall(PetscIntStackPop(temp_stack, &stage));
34453e0a2f3SToby Isaac           PetscCall(PetscIntStackEmpty(temp_stack, &empty));
34553e0a2f3SToby Isaac           if (!empty) {
34653e0a2f3SToby Isaac             PetscCall(PetscIntStackTop(temp_stack, &petsc_log_state->current_stage));
34753e0a2f3SToby Isaac           } else petsc_log_state->current_stage = -1;
34853e0a2f3SToby Isaac           PetscCall(PetscLogHandlerStagePop(h, stage));
34953e0a2f3SToby Isaac         }
35053e0a2f3SToby Isaac         PetscCall(PetscIntStackDestroy(temp_stack));
35153e0a2f3SToby Isaac         petsc_log_state->stage_stack = orig_stack;
35253e0a2f3SToby Isaac         PetscCall(PetscIntStackTop(petsc_log_state->stage_stack, &petsc_log_state->current_stage));
35353e0a2f3SToby Isaac       }
35453e0a2f3SToby Isaac       PetscCall(PetscArrayzero(&PetscLogHandlers[i], 1));
35553e0a2f3SToby Isaac       PetscCall(PetscObjectDereference((PetscObject)h));
35653e0a2f3SToby Isaac     }
35753e0a2f3SToby Isaac   }
35853e0a2f3SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
35953e0a2f3SToby Isaac }
3605c6c1daeSBarry Smith 
361cc4c1da9SBarry Smith /*@
3620b4b7b1cSBarry Smith   PetscLogIsActive - Check if logging (profiling) is currently in progress.
3634dd65854SConnor Ward 
3644dd65854SConnor Ward   Not Collective
3654dd65854SConnor Ward 
3664dd65854SConnor Ward   Output Parameter:
367811af0c4SBarry Smith . isActive - `PETSC_TRUE` if logging is in progress, `PETSC_FALSE` otherwise
3684dd65854SConnor Ward 
3694dd65854SConnor Ward   Level: beginner
3704dd65854SConnor Ward 
371b665b14eSToby Isaac .seealso: [](ch_profiling), `PetscLogDefaultBegin()`
3724dd65854SConnor Ward @*/
373d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogIsActive(PetscBool *isActive)
374d71ae5a4SJacob Faibussowitsch {
3754dd65854SConnor Ward   PetscFunctionBegin;
376b665b14eSToby Isaac   *isActive = PETSC_FALSE;
377b665b14eSToby Isaac   if (petsc_log_state) {
378b665b14eSToby Isaac     for (PetscInt i = 0; i < PETSC_LOG_HANDLER_MAX; i++) {
379b665b14eSToby Isaac       if (PetscLogHandlers[i].handler) {
380b665b14eSToby Isaac         *isActive = PETSC_TRUE;
381b665b14eSToby Isaac         PetscFunctionReturn(PETSC_SUCCESS);
382b665b14eSToby Isaac       }
383b665b14eSToby Isaac     }
384b665b14eSToby Isaac   }
385b665b14eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
386b665b14eSToby Isaac }
387b665b14eSToby Isaac 
388b665b14eSToby Isaac PETSC_UNUSED static PetscErrorCode PetscLogEventBeginIsActive(PetscBool *isActive)
389b665b14eSToby Isaac {
390b665b14eSToby Isaac   PetscFunctionBegin;
391b665b14eSToby Isaac   *isActive = PETSC_FALSE;
392b665b14eSToby Isaac   if (petsc_log_state) {
393b665b14eSToby Isaac     for (PetscInt i = 0; i < PETSC_LOG_HANDLER_MAX; i++) {
394b665b14eSToby Isaac       if (PetscLogHandlers[i].eventBegin) {
395b665b14eSToby Isaac         *isActive = PETSC_TRUE;
396b665b14eSToby Isaac         PetscFunctionReturn(PETSC_SUCCESS);
397b665b14eSToby Isaac       }
398b665b14eSToby Isaac     }
399b665b14eSToby Isaac   }
400b665b14eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
401b665b14eSToby Isaac }
402b665b14eSToby Isaac 
403b665b14eSToby Isaac PETSC_UNUSED static PetscErrorCode PetscLogEventEndIsActive(PetscBool *isActive)
404b665b14eSToby Isaac {
405b665b14eSToby Isaac   PetscFunctionBegin;
406b665b14eSToby Isaac   *isActive = PETSC_FALSE;
407b665b14eSToby Isaac   if (petsc_log_state) {
408b665b14eSToby Isaac     for (PetscInt i = 0; i < PETSC_LOG_HANDLER_MAX; i++) {
409b665b14eSToby Isaac       if (PetscLogHandlers[i].eventEnd) {
410b665b14eSToby Isaac         *isActive = PETSC_TRUE;
411b665b14eSToby Isaac         PetscFunctionReturn(PETSC_SUCCESS);
412b665b14eSToby Isaac       }
413b665b14eSToby Isaac     }
414b665b14eSToby Isaac   }
4153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4164dd65854SConnor Ward }
4174dd65854SConnor Ward 
41861cc7448SToby Isaac PETSC_INTERN PetscErrorCode PetscLogTypeBegin(PetscLogHandlerType type)
41953e0a2f3SToby Isaac {
42053e0a2f3SToby Isaac   PetscLogHandler handler;
42153e0a2f3SToby Isaac 
42253e0a2f3SToby Isaac   PetscFunctionBegin;
42353e0a2f3SToby Isaac   PetscCall(PetscLogTryGetHandler(type, &handler));
42453e0a2f3SToby Isaac   if (handler) PetscFunctionReturn(PETSC_SUCCESS);
42553e0a2f3SToby Isaac   PetscCall(PetscLogHandlerCreate(PETSC_COMM_WORLD, &handler));
42653e0a2f3SToby Isaac   PetscCall(PetscLogHandlerSetType(handler, type));
42753e0a2f3SToby Isaac   PetscCall(PetscLogHandlerStart(handler));
42853e0a2f3SToby Isaac   PetscCall(PetscLogHandlerDestroy(&handler));
42953e0a2f3SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
43053e0a2f3SToby Isaac }
43153e0a2f3SToby Isaac 
432cc4c1da9SBarry Smith /*@
4330b4b7b1cSBarry Smith   PetscLogDefaultBegin - Turns on logging (profiling) of PETSc code using the default log handler (profiler). This logs time, flop
4340b4b7b1cSBarry Smith   rates, and object creation and should not slow programs down too much.
4355c6c1daeSBarry Smith 
4368f14a041SBarry Smith   Logically Collective on `PETSC_COMM_WORLD`
4375c6c1daeSBarry Smith 
438811af0c4SBarry Smith   Options Database Key:
4390b4b7b1cSBarry Smith . -log_view [viewertype:filename:viewerformat] - Prints summary of flop and timing (profiling) information to the
4400b4b7b1cSBarry Smith                                                  screen (for PETSc configured with `--with-log=1` (which is the default)).
4410b4b7b1cSBarry Smith                                                  This option must be provided before `PetscInitialize()`.
4425c6c1daeSBarry Smith 
44310450e9eSJacob Faibussowitsch   Example Usage:
4445c6c1daeSBarry Smith .vb
4455c6c1daeSBarry Smith       PetscInitialize(...);
446bb1d7374SBarry Smith       PetscLogDefaultBegin();
4475c6c1daeSBarry Smith        ... code ...
4485c6c1daeSBarry Smith       PetscLogView(viewer); or PetscLogDump();
4495c6c1daeSBarry Smith       PetscFinalize();
4505c6c1daeSBarry Smith .ve
4515c6c1daeSBarry Smith 
452d1f92df0SBarry Smith   Level: advanced
453d1f92df0SBarry Smith 
4540b4b7b1cSBarry Smith   Notes:
455811af0c4SBarry Smith   `PetscLogView()` or `PetscLogDump()` actually cause the printing of
4565c6c1daeSBarry Smith   the logging information.
4575c6c1daeSBarry Smith 
4580b4b7b1cSBarry Smith   This routine may be called more than once.
4590b4b7b1cSBarry Smith 
4600b4b7b1cSBarry Smith   To provide the `-log_view` option in your source code you must call  PetscCall(PetscOptionsSetValue(NULL, "-log_view", NULL));
4610b4b7b1cSBarry Smith   before you call `PetscInitialize()`
4620b4b7b1cSBarry Smith 
463b665b14eSToby Isaac .seealso: [](ch_profiling), `PetscLogDump()`, `PetscLogView()`, `PetscLogTraceBegin()`
4645c6c1daeSBarry Smith @*/
465d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogDefaultBegin(void)
466d71ae5a4SJacob Faibussowitsch {
4675c6c1daeSBarry Smith   PetscFunctionBegin;
468294de794SToby Isaac   PetscCall(PetscLogTypeBegin(PETSCLOGHANDLERDEFAULT));
4693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4705c6c1daeSBarry Smith }
4715c6c1daeSBarry Smith 
4725c6c1daeSBarry Smith /*@C
473b665b14eSToby Isaac   PetscLogTraceBegin - Begins trace logging.  Every time a PETSc event
4745c6c1daeSBarry Smith   begins or ends, the event name is printed.
4755c6c1daeSBarry Smith 
476cc4c1da9SBarry Smith   Logically Collective on `PETSC_COMM_WORLD`, No Fortran Support
4775c6c1daeSBarry Smith 
4785c6c1daeSBarry Smith   Input Parameter:
4795c6c1daeSBarry Smith . file - The file to print trace in (e.g. stdout)
4805c6c1daeSBarry Smith 
4815c6c1daeSBarry Smith   Options Database Key:
482b665b14eSToby Isaac . -log_trace [filename] - Begins `PetscLogTraceBegin()`
4835c6c1daeSBarry Smith 
484d1f92df0SBarry Smith   Level: intermediate
485d1f92df0SBarry Smith 
4865c6c1daeSBarry Smith   Notes:
487811af0c4SBarry Smith   `PetscLogTraceBegin()` prints the processor number, the execution time (sec),
4885c6c1daeSBarry Smith   then "Event begin:" or "Event end:" followed by the event name.
4895c6c1daeSBarry Smith 
490811af0c4SBarry Smith   `PetscLogTraceBegin()` allows tracing of all PETSc calls, which is useful
4915c6c1daeSBarry Smith   to determine where a program is hanging without running in the
4925c6c1daeSBarry Smith   debugger.  Can be used in conjunction with the -info option.
4935c6c1daeSBarry Smith 
494b665b14eSToby Isaac .seealso: [](ch_profiling), `PetscLogDump()`, `PetscLogView()`, `PetscLogDefaultBegin()`
4955c6c1daeSBarry Smith @*/
496d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogTraceBegin(FILE *file)
497d71ae5a4SJacob Faibussowitsch {
498b665b14eSToby Isaac   PetscLogHandler handler;
4994d86920dSPierre Jolivet 
5005c6c1daeSBarry Smith   PetscFunctionBegin;
501294de794SToby Isaac   PetscCall(PetscLogTryGetHandler(PETSCLOGHANDLERTRACE, &handler));
502b665b14eSToby Isaac   if (handler) PetscFunctionReturn(PETSC_SUCCESS);
503b665b14eSToby Isaac   PetscCall(PetscLogHandlerCreateTrace(PETSC_COMM_WORLD, file, &handler));
504b665b14eSToby Isaac   PetscCall(PetscLogHandlerStart(handler));
505b665b14eSToby Isaac   PetscCall(PetscLogHandlerDestroy(&handler));
506b665b14eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
507b665b14eSToby Isaac }
508a297a907SKarl Rupp 
509b665b14eSToby Isaac PETSC_INTERN PetscErrorCode PetscLogHandlerCreate_Nested(MPI_Comm, PetscLogHandler *);
510b665b14eSToby Isaac 
511cc4c1da9SBarry Smith /*@
512b665b14eSToby Isaac   PetscLogNestedBegin - Turns on nested logging of objects and events. This logs flop
513b665b14eSToby Isaac   rates and object creation and should not slow programs down too much.
514b665b14eSToby Isaac 
515cc4c1da9SBarry Smith   Logically Collective on `PETSC_COMM_WORLD`, No Fortran Support
516b665b14eSToby Isaac 
517b665b14eSToby Isaac   Options Database Keys:
518b665b14eSToby Isaac . -log_view :filename.xml:ascii_xml - Prints an XML summary of flop and timing information to the file
519b665b14eSToby Isaac 
520b665b14eSToby Isaac   Example Usage:
521b665b14eSToby Isaac .vb
522b665b14eSToby Isaac       PetscInitialize(...);
523b665b14eSToby Isaac       PetscLogNestedBegin();
524b665b14eSToby Isaac        ... code ...
525b665b14eSToby Isaac       PetscLogView(viewer);
526b665b14eSToby Isaac       PetscFinalize();
527b665b14eSToby Isaac .ve
528b665b14eSToby Isaac 
529b665b14eSToby Isaac   Level: advanced
530b665b14eSToby Isaac 
531b665b14eSToby Isaac .seealso: `PetscLogDump()`, `PetscLogView()`, `PetscLogTraceBegin()`, `PetscLogDefaultBegin()`
532b665b14eSToby Isaac @*/
533b665b14eSToby Isaac PetscErrorCode PetscLogNestedBegin(void)
534b665b14eSToby Isaac {
535b665b14eSToby Isaac   PetscFunctionBegin;
536294de794SToby Isaac   PetscCall(PetscLogTypeBegin(PETSCLOGHANDLERNESTED));
5373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5385c6c1daeSBarry Smith }
5395c6c1daeSBarry Smith 
54053e0a2f3SToby Isaac /*@C
54153e0a2f3SToby Isaac   PetscLogLegacyCallbacksBegin - Create and start a log handler from callbacks
54253e0a2f3SToby Isaac   matching the now deprecated function pointers `PetscLogPLB`, `PetscLogPLE`,
54353e0a2f3SToby Isaac   `PetscLogPHC`, `PetscLogPHD`.
54453e0a2f3SToby Isaac 
5458f14a041SBarry Smith   Logically Collective on `PETSC_COMM_WORLD`
54653e0a2f3SToby Isaac 
54753e0a2f3SToby Isaac   Input Parameters:
54853e0a2f3SToby Isaac + PetscLogPLB - A callback that will be executed by `PetscLogEventBegin()` (or `NULL`)
54953e0a2f3SToby Isaac . PetscLogPLE - A callback that will be executed by `PetscLogEventEnd()` (or `NULL`)
55053e0a2f3SToby Isaac . PetscLogPHC - A callback that will be executed by `PetscLogObjectCreate()` (or `NULL`)
55153e0a2f3SToby Isaac - PetscLogPHD - A callback that will be executed by `PetscLogObjectCreate()` (or `NULL`)
55253e0a2f3SToby Isaac 
55353e0a2f3SToby Isaac   Calling sequence of `PetscLogPLB`:
55453e0a2f3SToby Isaac + e  - a `PetscLogEvent` that is beginning
55553e0a2f3SToby Isaac . _i - deprecated, unused
55653e0a2f3SToby Isaac . o1 - a `PetscObject` associated with `e` (or `NULL`)
55753e0a2f3SToby Isaac . o2 - a `PetscObject` associated with `e` (or `NULL`)
55853e0a2f3SToby Isaac . o3 - a `PetscObject` associated with `e` (or `NULL`)
55953e0a2f3SToby Isaac - o4 - a `PetscObject` associated with `e` (or `NULL`)
56053e0a2f3SToby Isaac 
56153e0a2f3SToby Isaac   Calling sequence of `PetscLogPLE`:
56253e0a2f3SToby Isaac + e  - a `PetscLogEvent` that is beginning
56353e0a2f3SToby Isaac . _i - deprecated, unused
56453e0a2f3SToby Isaac . o1 - a `PetscObject` associated with `e` (or `NULL`)
56553e0a2f3SToby Isaac . o2 - a `PetscObject` associated with `e` (or `NULL`)
56653e0a2f3SToby Isaac . o3 - a `PetscObject` associated with `e` (or `NULL`)
56753e0a2f3SToby Isaac - o4 - a `PetscObject` associated with `e` (or `NULL`)
56853e0a2f3SToby Isaac 
56953e0a2f3SToby Isaac   Calling sequence of `PetscLogPHC`:
57053e0a2f3SToby Isaac . o - a `PetscObject` that has just been created
57153e0a2f3SToby Isaac 
57253e0a2f3SToby Isaac   Calling sequence of `PetscLogPHD`:
57353e0a2f3SToby Isaac . o - a `PetscObject` that is about to be destroyed
57453e0a2f3SToby Isaac 
57553e0a2f3SToby Isaac   Level: advanced
57653e0a2f3SToby Isaac 
57753e0a2f3SToby Isaac   Notes:
57853e0a2f3SToby Isaac   This is for transitioning from the deprecated function `PetscLogSet()` and should not be used in new code.
57953e0a2f3SToby Isaac 
58053e0a2f3SToby Isaac   This should help migrate external log handlers to use `PetscLogHandler`, but
58153e0a2f3SToby Isaac   callbacks that depend on the deprecated `PetscLogStage` datatype will have to be
58253e0a2f3SToby Isaac   updated.
58353e0a2f3SToby Isaac 
58453e0a2f3SToby Isaac .seealso: [](ch_profiling), `PetscLogHandler`, `PetscLogHandlerStart()`, `PetscLogState`
58553e0a2f3SToby Isaac @*/
58653e0a2f3SToby Isaac PetscErrorCode PetscLogLegacyCallbacksBegin(PetscErrorCode (*PetscLogPLB)(PetscLogEvent e, int _i, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4), PetscErrorCode (*PetscLogPLE)(PetscLogEvent e, int _i, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4), PetscErrorCode (*PetscLogPHC)(PetscObject o), PetscErrorCode (*PetscLogPHD)(PetscObject o))
58753e0a2f3SToby Isaac {
58853e0a2f3SToby Isaac   PetscLogHandler handler;
58953e0a2f3SToby Isaac 
59053e0a2f3SToby Isaac   PetscFunctionBegin;
59153e0a2f3SToby Isaac   PetscCall(PetscLogHandlerCreateLegacy(PETSC_COMM_WORLD, PetscLogPLB, PetscLogPLE, PetscLogPHC, PetscLogPHD, &handler));
59253e0a2f3SToby Isaac   PetscCall(PetscLogHandlerStart(handler));
59353e0a2f3SToby Isaac   PetscCall(PetscLogHandlerDestroy(&handler));
59453e0a2f3SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
59553e0a2f3SToby Isaac }
59653e0a2f3SToby Isaac 
5972611ad71SToby Isaac   #if defined(PETSC_HAVE_MPE)
5982611ad71SToby Isaac     #include <mpe.h>
5992611ad71SToby Isaac static PetscBool PetscBeganMPE = PETSC_FALSE;
6002611ad71SToby Isaac   #endif
6012611ad71SToby Isaac 
6022611ad71SToby Isaac /*@C
6032611ad71SToby Isaac   PetscLogMPEBegin - Turns on MPE logging of events. This creates large log files and slows the
6042611ad71SToby Isaac   program down.
6052611ad71SToby Isaac 
606cc4c1da9SBarry Smith   Collective on `PETSC_COMM_WORLD`, No Fortran Support
6072611ad71SToby Isaac 
6082611ad71SToby Isaac   Options Database Key:
6092611ad71SToby Isaac . -log_mpe - Prints extensive log information
6102611ad71SToby Isaac 
6112611ad71SToby Isaac   Level: advanced
6122611ad71SToby Isaac 
6132611ad71SToby Isaac   Note:
6142611ad71SToby Isaac   A related routine is `PetscLogDefaultBegin()` (with the options key `-log_view`), which is
6152611ad71SToby Isaac   intended for production runs since it logs only flop rates and object creation (and should
6162611ad71SToby Isaac   not significantly slow the programs).
6172611ad71SToby Isaac 
618b665b14eSToby Isaac .seealso: [](ch_profiling), `PetscLogDump()`, `PetscLogDefaultBegin()`, `PetscLogEventActivate()`,
6192611ad71SToby Isaac           `PetscLogEventDeactivate()`
6202611ad71SToby Isaac @*/
6212611ad71SToby Isaac PetscErrorCode PetscLogMPEBegin(void)
6222611ad71SToby Isaac {
6232611ad71SToby Isaac   PetscFunctionBegin;
6242611ad71SToby Isaac   #if defined(PETSC_HAVE_MPE)
6252611ad71SToby Isaac   /* Do MPE initialization */
6262611ad71SToby Isaac   if (!MPE_Initialized_logging()) { /* This function exists in mpich 1.1.2 and higher */
6272611ad71SToby Isaac     PetscCall(PetscInfo(0, "Initializing MPE.\n"));
6282611ad71SToby Isaac     PetscCall(MPE_Init_log());
6292611ad71SToby Isaac 
6302611ad71SToby Isaac     PetscBeganMPE = PETSC_TRUE;
6312611ad71SToby Isaac   } else {
6322611ad71SToby Isaac     PetscCall(PetscInfo(0, "MPE already initialized. Not attempting to reinitialize.\n"));
6332611ad71SToby Isaac   }
634294de794SToby Isaac   PetscCall(PetscLogTypeBegin(PETSCLOGHANDLERMPE));
6352611ad71SToby Isaac   #else
6362611ad71SToby Isaac   SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, "PETSc was configured without MPE support, reconfigure with --with-mpe or --download-mpe");
6372611ad71SToby Isaac   #endif
6382611ad71SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
6392611ad71SToby Isaac }
6402611ad71SToby Isaac 
64153e0a2f3SToby Isaac   #if defined(PETSC_HAVE_TAU_PERFSTUBS)
64253e0a2f3SToby Isaac     #include <../src/sys/perfstubs/timer.h>
64353e0a2f3SToby Isaac   #endif
64453e0a2f3SToby Isaac 
64553e0a2f3SToby Isaac /*@C
64653e0a2f3SToby Isaac   PetscLogPerfstubsBegin - Turns on logging of events using the perfstubs interface.
64753e0a2f3SToby Isaac 
648cc4c1da9SBarry Smith   Collective on `PETSC_COMM_WORLD`, No Fortran Support
64953e0a2f3SToby Isaac 
65053e0a2f3SToby Isaac   Options Database Key:
65153e0a2f3SToby Isaac . -log_perfstubs - use an external log handler through the perfstubs interface
65253e0a2f3SToby Isaac 
65353e0a2f3SToby Isaac   Level: advanced
65453e0a2f3SToby Isaac 
65553e0a2f3SToby Isaac .seealso: [](ch_profiling), `PetscLogDefaultBegin()`, `PetscLogEventActivate()`
65653e0a2f3SToby Isaac @*/
65753e0a2f3SToby Isaac PetscErrorCode PetscLogPerfstubsBegin(void)
65853e0a2f3SToby Isaac {
65953e0a2f3SToby Isaac   PetscFunctionBegin;
66053e0a2f3SToby Isaac   #if defined(PETSC_HAVE_TAU_PERFSTUBS)
661294de794SToby Isaac   PetscCall(PetscLogTypeBegin(PETSCLOGHANDLERPERFSTUBS));
66253e0a2f3SToby Isaac   #else
66353e0a2f3SToby Isaac   SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, "PETSc was configured without perfstubs support, reconfigure with --with-tau-perfstubs");
66453e0a2f3SToby Isaac   #endif
66553e0a2f3SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
66653e0a2f3SToby Isaac }
66753e0a2f3SToby Isaac 
6685c6c1daeSBarry Smith /*@
669b665b14eSToby Isaac   PetscLogActions - Determines whether actions are logged for the default log handler.
6705c6c1daeSBarry Smith 
6715c6c1daeSBarry Smith   Not Collective
6725c6c1daeSBarry Smith 
6735c6c1daeSBarry Smith   Input Parameter:
674811af0c4SBarry Smith . flag - `PETSC_TRUE` if actions are to be logged
675811af0c4SBarry Smith 
676811af0c4SBarry Smith   Options Database Key:
677b665b14eSToby Isaac + -log_exclude_actions - (deprecated) Does nothing
678b665b14eSToby Isaac - -log_include_actions - Turn on action logging
6795c6c1daeSBarry Smith 
6805c6c1daeSBarry Smith   Level: intermediate
6815c6c1daeSBarry Smith 
682811af0c4SBarry Smith   Note:
683811af0c4SBarry Smith   Logging of actions continues to consume more memory as the program
6845c6c1daeSBarry Smith   runs. Long running programs should consider turning this feature off.
685aec76313SJacob Faibussowitsch 
686b665b14eSToby Isaac .seealso: [](ch_profiling), `PetscLogStagePush()`, `PetscLogStagePop()`, `PetscLogGetDefaultHandler()`
6875c6c1daeSBarry Smith @*/
688d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogActions(PetscBool flag)
689d71ae5a4SJacob Faibussowitsch {
6905c6c1daeSBarry Smith   PetscFunctionBegin;
691dff009beSToby Isaac   for (PetscInt i = 0; i < PETSC_LOG_HANDLER_MAX; i++) {
692dff009beSToby Isaac     PetscLogHandler h = PetscLogHandlers[i].handler;
693dff009beSToby Isaac 
694dff009beSToby Isaac     if (h) PetscCall(PetscLogHandlerSetLogActions(h, flag));
695dff009beSToby Isaac   }
6963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6975c6c1daeSBarry Smith }
6985c6c1daeSBarry Smith 
6995c6c1daeSBarry Smith /*@
7005c6c1daeSBarry Smith   PetscLogObjects - Determines whether objects are logged for the graphical viewer.
7015c6c1daeSBarry Smith 
7025c6c1daeSBarry Smith   Not Collective
7035c6c1daeSBarry Smith 
7045c6c1daeSBarry Smith   Input Parameter:
705811af0c4SBarry Smith . flag - `PETSC_TRUE` if objects are to be logged
706811af0c4SBarry Smith 
707811af0c4SBarry Smith   Options Database Key:
708b665b14eSToby Isaac + -log_exclude_objects - (deprecated) Does nothing
709b665b14eSToby Isaac - -log_include_objects - Turns on object logging
7105c6c1daeSBarry Smith 
7115c6c1daeSBarry Smith   Level: intermediate
7125c6c1daeSBarry Smith 
713811af0c4SBarry Smith   Note:
714811af0c4SBarry Smith   Logging of objects continues to consume more memory as the program
7155c6c1daeSBarry Smith   runs. Long running programs should consider turning this feature off.
7165c6c1daeSBarry Smith 
717b665b14eSToby Isaac .seealso: [](ch_profiling), `PetscLogStagePush()`, `PetscLogStagePop()`, `PetscLogGetDefaultHandler()`
7185c6c1daeSBarry Smith @*/
719d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogObjects(PetscBool flag)
720d71ae5a4SJacob Faibussowitsch {
7215c6c1daeSBarry Smith   PetscFunctionBegin;
722dff009beSToby Isaac   for (PetscInt i = 0; i < PETSC_LOG_HANDLER_MAX; i++) {
723dff009beSToby Isaac     PetscLogHandler h = PetscLogHandlers[i].handler;
724dff009beSToby Isaac 
725dff009beSToby Isaac     if (h) PetscCall(PetscLogHandlerSetLogObjects(h, flag));
726dff009beSToby Isaac   }
7273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7285c6c1daeSBarry Smith }
7295c6c1daeSBarry Smith 
7305c6c1daeSBarry Smith /*------------------------------------------------ Stage Functions --------------------------------------------------*/
731cc4c1da9SBarry Smith /*@
73274c0405dSRichard Tran Mills   PetscLogStageRegister - Attaches a character string name to a logging stage.
7335c6c1daeSBarry Smith 
7345c6c1daeSBarry Smith   Not Collective
7355c6c1daeSBarry Smith 
7365c6c1daeSBarry Smith   Input Parameter:
7375c6c1daeSBarry Smith . sname - The name to associate with that stage
7385c6c1daeSBarry Smith 
7395c6c1daeSBarry Smith   Output Parameter:
740b665b14eSToby Isaac . stage - The stage number or -1 if logging is not active (`PetscLogIsActive()`).
7415c6c1daeSBarry Smith 
7425c6c1daeSBarry Smith   Level: intermediate
7435c6c1daeSBarry Smith 
744d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogStagePush()`, `PetscLogStagePop()`
7455c6c1daeSBarry Smith @*/
746d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogStageRegister(const char sname[], PetscLogStage *stage)
747d71ae5a4SJacob Faibussowitsch {
748b665b14eSToby Isaac   PetscLogState state;
7495c6c1daeSBarry Smith 
7505c6c1daeSBarry Smith   PetscFunctionBegin;
751b665b14eSToby Isaac   *stage = -1;
752b665b14eSToby Isaac   PetscCall(PetscLogGetState(&state));
753b665b14eSToby Isaac   if (state) PetscCall(PetscLogStateStageRegister(state, sname, stage));
7543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7555c6c1daeSBarry Smith }
7565c6c1daeSBarry Smith 
757cc4c1da9SBarry Smith /*@
758811af0c4SBarry Smith   PetscLogStagePush - This function pushes a stage on the logging stack. Events started and stopped until `PetscLogStagePop()` will be associated with the stage
7595c6c1daeSBarry Smith 
7605c6c1daeSBarry Smith   Not Collective
7615c6c1daeSBarry Smith 
7625c6c1daeSBarry Smith   Input Parameter:
7635c6c1daeSBarry Smith . stage - The stage on which to log
7645c6c1daeSBarry Smith 
76510450e9eSJacob Faibussowitsch   Example Usage:
766811af0c4SBarry Smith   If the option -log_view is used to run the program containing the
7675c6c1daeSBarry Smith   following code, then 2 sets of summary data will be printed during
7685c6c1daeSBarry Smith   PetscFinalize().
7695c6c1daeSBarry Smith .vb
7705c6c1daeSBarry Smith       PetscInitialize(int *argc,char ***args,0,0);
7715c6c1daeSBarry Smith       [stage 0 of code]
7725c6c1daeSBarry Smith       PetscLogStagePush(1);
7735c6c1daeSBarry Smith       [stage 1 of code]
7745c6c1daeSBarry Smith       PetscLogStagePop();
7755c6c1daeSBarry Smith       PetscBarrier(...);
7765c6c1daeSBarry Smith       [more stage 0 of code]
7775c6c1daeSBarry Smith       PetscFinalize();
7785c6c1daeSBarry Smith .ve
7795c6c1daeSBarry Smith 
780d1f92df0SBarry Smith   Level: intermediate
781d1f92df0SBarry Smith 
782811af0c4SBarry Smith   Note:
783811af0c4SBarry Smith   Use `PetscLogStageRegister()` to register a stage.
7845c6c1daeSBarry Smith 
785d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogStagePop()`, `PetscLogStageRegister()`, `PetscBarrier()`
7865c6c1daeSBarry Smith @*/
787d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogStagePush(PetscLogStage stage)
788d71ae5a4SJacob Faibussowitsch {
789b665b14eSToby Isaac   PetscLogState state;
7905c6c1daeSBarry Smith 
7915c6c1daeSBarry Smith   PetscFunctionBegin;
792b665b14eSToby Isaac   PetscCall(PetscLogGetState(&state));
793b665b14eSToby Isaac   if (!state) PetscFunctionReturn(PETSC_SUCCESS);
794b665b14eSToby Isaac   for (int i = 0; i < PETSC_LOG_HANDLER_MAX; i++) {
795b665b14eSToby Isaac     PetscLogHandler h = PetscLogHandlers[i].handler;
796b665b14eSToby Isaac     if (h) PetscCall(PetscLogHandlerStagePush(h, stage));
797b665b14eSToby Isaac   }
798b665b14eSToby Isaac   PetscCall(PetscLogStateStagePush(state, stage));
7993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8005c6c1daeSBarry Smith }
8015c6c1daeSBarry Smith 
802cc4c1da9SBarry Smith /*@
803811af0c4SBarry Smith   PetscLogStagePop - This function pops a stage from the logging stack that was pushed with `PetscLogStagePush()`
8045c6c1daeSBarry Smith 
8055c6c1daeSBarry Smith   Not Collective
8065c6c1daeSBarry Smith 
80710450e9eSJacob Faibussowitsch   Example Usage:
808811af0c4SBarry Smith   If the option -log_view is used to run the program containing the
8095c6c1daeSBarry Smith   following code, then 2 sets of summary data will be printed during
8105c6c1daeSBarry Smith   PetscFinalize().
8115c6c1daeSBarry Smith .vb
8125c6c1daeSBarry Smith       PetscInitialize(int *argc,char ***args,0,0);
8135c6c1daeSBarry Smith       [stage 0 of code]
8145c6c1daeSBarry Smith       PetscLogStagePush(1);
8155c6c1daeSBarry Smith       [stage 1 of code]
8165c6c1daeSBarry Smith       PetscLogStagePop();
8175c6c1daeSBarry Smith       PetscBarrier(...);
8185c6c1daeSBarry Smith       [more stage 0 of code]
8195c6c1daeSBarry Smith       PetscFinalize();
8205c6c1daeSBarry Smith .ve
8215c6c1daeSBarry Smith 
8225c6c1daeSBarry Smith   Level: intermediate
8235c6c1daeSBarry Smith 
824d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogStagePush()`, `PetscLogStageRegister()`, `PetscBarrier()`
8255c6c1daeSBarry Smith @*/
826d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogStagePop(void)
827d71ae5a4SJacob Faibussowitsch {
828b665b14eSToby Isaac   PetscLogState state;
829b665b14eSToby Isaac   PetscLogStage current_stage;
8305c6c1daeSBarry Smith 
8315c6c1daeSBarry Smith   PetscFunctionBegin;
832b665b14eSToby Isaac   PetscCall(PetscLogGetState(&state));
833b665b14eSToby Isaac   if (!state) PetscFunctionReturn(PETSC_SUCCESS);
834b665b14eSToby Isaac   current_stage = state->current_stage;
835b665b14eSToby Isaac   PetscCall(PetscLogStateStagePop(state));
836b665b14eSToby Isaac   for (int i = 0; i < PETSC_LOG_HANDLER_MAX; i++) {
837b665b14eSToby Isaac     PetscLogHandler h = PetscLogHandlers[i].handler;
838b665b14eSToby Isaac     if (h) PetscCall(PetscLogHandlerStagePop(h, current_stage));
839b665b14eSToby Isaac   }
8403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8415c6c1daeSBarry Smith }
8425c6c1daeSBarry Smith 
8435c6c1daeSBarry Smith /*@
844811af0c4SBarry Smith   PetscLogStageSetActive - Sets if a stage is used for `PetscLogEventBegin()` and `PetscLogEventEnd()`.
8455c6c1daeSBarry Smith 
8465c6c1daeSBarry Smith   Not Collective
8475c6c1daeSBarry Smith 
8485c6c1daeSBarry Smith   Input Parameters:
8495c6c1daeSBarry Smith + stage    - The stage
850811af0c4SBarry Smith - isActive - The activity flag, `PETSC_TRUE` for logging, else `PETSC_FALSE` (defaults to `PETSC_TRUE`)
8515c6c1daeSBarry Smith 
8525c6c1daeSBarry Smith   Level: intermediate
8535c6c1daeSBarry Smith 
854811af0c4SBarry Smith   Note:
855811af0c4SBarry Smith   If this is set to `PETSC_FALSE` the logging acts as if the stage did not exist
856811af0c4SBarry Smith 
857d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogStageRegister()`, `PetscLogStagePush()`, `PetscLogStagePop()`, `PetscLogEventBegin()`, `PetscLogEventEnd()`, `PetscPreLoadBegin()`, `PetscPreLoadEnd()`, `PetscPreLoadStage()`
8585c6c1daeSBarry Smith @*/
859d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogStageSetActive(PetscLogStage stage, PetscBool isActive)
860d71ae5a4SJacob Faibussowitsch {
861b665b14eSToby Isaac   PetscLogState state;
8625c6c1daeSBarry Smith 
8635c6c1daeSBarry Smith   PetscFunctionBegin;
864b665b14eSToby Isaac   PetscCall(PetscLogGetState(&state));
865b665b14eSToby Isaac   if (state) PetscCall(PetscLogStateStageSetActive(state, stage, isActive));
8663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8675c6c1daeSBarry Smith }
8685c6c1daeSBarry Smith 
8695c6c1daeSBarry Smith /*@
870811af0c4SBarry Smith   PetscLogStageGetActive - Checks if a stage is used for `PetscLogEventBegin()` and `PetscLogEventEnd()`.
8715c6c1daeSBarry Smith 
8725c6c1daeSBarry Smith   Not Collective
8735c6c1daeSBarry Smith 
8745c6c1daeSBarry Smith   Input Parameter:
8755c6c1daeSBarry Smith . stage - The stage
8765c6c1daeSBarry Smith 
8775c6c1daeSBarry Smith   Output Parameter:
878811af0c4SBarry Smith . isActive - The activity flag, `PETSC_TRUE` for logging, else `PETSC_FALSE` (defaults to `PETSC_TRUE`)
8795c6c1daeSBarry Smith 
8805c6c1daeSBarry Smith   Level: intermediate
8815c6c1daeSBarry Smith 
882d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogStageRegister()`, `PetscLogStagePush()`, `PetscLogStagePop()`, `PetscLogEventBegin()`, `PetscLogEventEnd()`, `PetscPreLoadBegin()`, `PetscPreLoadEnd()`, `PetscPreLoadStage()`
8835c6c1daeSBarry Smith @*/
884d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogStageGetActive(PetscLogStage stage, PetscBool *isActive)
885d71ae5a4SJacob Faibussowitsch {
886b665b14eSToby Isaac   PetscLogState state;
8875c6c1daeSBarry Smith 
8885c6c1daeSBarry Smith   PetscFunctionBegin;
889b665b14eSToby Isaac   *isActive = PETSC_FALSE;
890b665b14eSToby Isaac   PetscCall(PetscLogGetState(&state));
891b665b14eSToby Isaac   if (state) PetscCall(PetscLogStateStageGetActive(state, stage, isActive));
8923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8935c6c1daeSBarry Smith }
8945c6c1daeSBarry Smith 
8955c6c1daeSBarry Smith /*@
896811af0c4SBarry Smith   PetscLogStageSetVisible - Determines stage visibility in `PetscLogView()`
8975c6c1daeSBarry Smith 
8985c6c1daeSBarry Smith   Not Collective
8995c6c1daeSBarry Smith 
9005c6c1daeSBarry Smith   Input Parameters:
9015c6c1daeSBarry Smith + stage     - The stage
902811af0c4SBarry Smith - isVisible - The visibility flag, `PETSC_TRUE` to print, else `PETSC_FALSE` (defaults to `PETSC_TRUE`)
9035c6c1daeSBarry Smith 
9045c6c1daeSBarry Smith   Level: intermediate
9055c6c1daeSBarry Smith 
906aec76313SJacob Faibussowitsch   Developer Notes:
907b665b14eSToby Isaac   Visibility only affects the default log handler in `PetscLogView()`: stages that are
908b665b14eSToby Isaac   set to invisible are suppressed from output.
909811af0c4SBarry Smith 
910b665b14eSToby Isaac .seealso: [](ch_profiling), `PetscLogStageGetVisible()`, `PetscLogStageRegister()`, `PetscLogStagePush()`, `PetscLogStagePop()`, `PetscLogView()`, `PetscLogGetDefaultHandler()`
9115c6c1daeSBarry Smith @*/
912d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogStageSetVisible(PetscLogStage stage, PetscBool isVisible)
9135c6c1daeSBarry Smith 
914dff009beSToby Isaac {
9155c6c1daeSBarry Smith   PetscFunctionBegin;
916dff009beSToby Isaac   for (PetscInt i = 0; i < PETSC_LOG_HANDLER_MAX; i++) {
917dff009beSToby Isaac     PetscLogHandler h = PetscLogHandlers[i].handler;
918dff009beSToby Isaac 
919dff009beSToby Isaac     if (h) PetscCall(PetscLogHandlerStageSetVisible(h, stage, isVisible));
920dff009beSToby Isaac   }
9213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9225c6c1daeSBarry Smith }
9235c6c1daeSBarry Smith 
9245c6c1daeSBarry Smith /*@
925811af0c4SBarry Smith   PetscLogStageGetVisible - Returns stage visibility in `PetscLogView()`
9265c6c1daeSBarry Smith 
9275c6c1daeSBarry Smith   Not Collective
9285c6c1daeSBarry Smith 
9295c6c1daeSBarry Smith   Input Parameter:
9305c6c1daeSBarry Smith . stage - The stage
9315c6c1daeSBarry Smith 
9325c6c1daeSBarry Smith   Output Parameter:
933811af0c4SBarry Smith . isVisible - The visibility flag, `PETSC_TRUE` to print, else `PETSC_FALSE` (defaults to `PETSC_TRUE`)
9345c6c1daeSBarry Smith 
9355c6c1daeSBarry Smith   Level: intermediate
9365c6c1daeSBarry Smith 
937b665b14eSToby Isaac .seealso: [](ch_profiling), `PetscLogStageSetVisible()`, `PetscLogStageRegister()`, `PetscLogStagePush()`, `PetscLogStagePop()`, `PetscLogView()`, `PetscLogGetDefaultHandler()`
9385c6c1daeSBarry Smith @*/
939d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogStageGetVisible(PetscLogStage stage, PetscBool *isVisible)
940d71ae5a4SJacob Faibussowitsch {
941b665b14eSToby Isaac   PetscLogHandler handler;
9425c6c1daeSBarry Smith 
9435c6c1daeSBarry Smith   PetscFunctionBegin;
944b665b14eSToby Isaac   *isVisible = PETSC_FALSE;
945294de794SToby Isaac   PetscCall(PetscLogTryGetHandler(PETSCLOGHANDLERDEFAULT, &handler));
946dff009beSToby Isaac   if (handler) { PetscCall(PetscLogHandlerStageGetVisible(handler, stage, isVisible)); }
9473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9485c6c1daeSBarry Smith }
9495c6c1daeSBarry Smith 
950cc4c1da9SBarry Smith /*@
9515c6c1daeSBarry Smith   PetscLogStageGetId - Returns the stage id when given the stage name.
9525c6c1daeSBarry Smith 
9535c6c1daeSBarry Smith   Not Collective
9545c6c1daeSBarry Smith 
9555c6c1daeSBarry Smith   Input Parameter:
9565c6c1daeSBarry Smith . name - The stage name
9575c6c1daeSBarry Smith 
9585c6c1daeSBarry Smith   Output Parameter:
9595a4a3fabSBarry Smith . stage - The stage, , or -1 if no stage with that name exists
9605c6c1daeSBarry Smith 
9615c6c1daeSBarry Smith   Level: intermediate
9625c6c1daeSBarry Smith 
963d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogStageRegister()`, `PetscLogStagePush()`, `PetscLogStagePop()`, `PetscPreLoadBegin()`, `PetscPreLoadEnd()`, `PetscPreLoadStage()`
9645c6c1daeSBarry Smith @*/
965d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogStageGetId(const char name[], PetscLogStage *stage)
966d71ae5a4SJacob Faibussowitsch {
967b665b14eSToby Isaac   PetscLogState state;
9685c6c1daeSBarry Smith 
9695c6c1daeSBarry Smith   PetscFunctionBegin;
970b665b14eSToby Isaac   *stage = -1;
971b665b14eSToby Isaac   PetscCall(PetscLogGetState(&state));
972b665b14eSToby Isaac   if (state) PetscCall(PetscLogStateGetStageFromName(state, name, stage));
9733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9745c6c1daeSBarry Smith }
9755c6c1daeSBarry Smith 
976cc4c1da9SBarry Smith /*@
97753e0a2f3SToby Isaac   PetscLogStageGetName - Returns the stage name when given the stage id.
97853e0a2f3SToby Isaac 
97953e0a2f3SToby Isaac   Not Collective
98053e0a2f3SToby Isaac 
98153e0a2f3SToby Isaac   Input Parameter:
98253e0a2f3SToby Isaac . stage - The stage
98353e0a2f3SToby Isaac 
98453e0a2f3SToby Isaac   Output Parameter:
98553e0a2f3SToby Isaac . name - The stage name
98653e0a2f3SToby Isaac 
98753e0a2f3SToby Isaac   Level: intermediate
98853e0a2f3SToby Isaac 
98953e0a2f3SToby Isaac .seealso: [](ch_profiling), `PetscLogStageRegister()`, `PetscLogStagePush()`, `PetscLogStagePop()`, `PetscPreLoadBegin()`, `PetscPreLoadEnd()`, `PetscPreLoadStage()`
99053e0a2f3SToby Isaac @*/
991cc4c1da9SBarry Smith PetscErrorCode PetscLogStageGetName(PetscLogStage stage, const char *name[])
99253e0a2f3SToby Isaac {
99353e0a2f3SToby Isaac   PetscLogStageInfo stage_info;
99453e0a2f3SToby Isaac   PetscLogState     state;
99553e0a2f3SToby Isaac 
99653e0a2f3SToby Isaac   PetscFunctionBegin;
997b665b14eSToby Isaac   *name = NULL;
99853e0a2f3SToby Isaac   PetscCall(PetscLogGetState(&state));
999b665b14eSToby Isaac   if (!state) PetscFunctionReturn(PETSC_SUCCESS);
100053e0a2f3SToby Isaac   PetscCall(PetscLogStateStageGetInfo(state, stage, &stage_info));
100153e0a2f3SToby Isaac   *name = stage_info.name;
100253e0a2f3SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
100353e0a2f3SToby Isaac }
100453e0a2f3SToby Isaac 
10055c6c1daeSBarry Smith /*------------------------------------------------ Event Functions --------------------------------------------------*/
10067a101e5eSJacob Faibussowitsch 
1007cc4c1da9SBarry Smith /*@
1008811af0c4SBarry Smith   PetscLogEventRegister - Registers an event name for logging operations
10095c6c1daeSBarry Smith 
10105c6c1daeSBarry Smith   Not Collective
10115c6c1daeSBarry Smith 
1012d8d19677SJose E. Roman   Input Parameters:
10135c6c1daeSBarry Smith + name    - The name associated with the event
10145c6c1daeSBarry Smith - classid - The classid associated to the class for this event, obtain either with
1015811af0c4SBarry Smith            `PetscClassIdRegister()` or use a predefined one such as `KSP_CLASSID`, `SNES_CLASSID`, the predefined ones
10165c6c1daeSBarry Smith            are only available in C code
10175c6c1daeSBarry Smith 
10185c6c1daeSBarry Smith   Output Parameter:
1019811af0c4SBarry Smith . event - The event id for use with `PetscLogEventBegin()` and `PetscLogEventEnd()`.
10205c6c1daeSBarry Smith 
102110450e9eSJacob Faibussowitsch   Example Usage:
10225c6c1daeSBarry Smith .vb
10235c6c1daeSBarry Smith       PetscLogEvent USER_EVENT;
10245c6c1daeSBarry Smith       PetscClassId classid;
10255c6c1daeSBarry Smith       PetscLogDouble user_event_flops;
10265c6c1daeSBarry Smith       PetscClassIdRegister("class name",&classid);
10275c6c1daeSBarry Smith       PetscLogEventRegister("User event name",classid,&USER_EVENT);
10285c6c1daeSBarry Smith       PetscLogEventBegin(USER_EVENT,0,0,0,0);
10295c6c1daeSBarry Smith          [code segment to monitor]
10305c6c1daeSBarry Smith          PetscLogFlops(user_event_flops);
10315c6c1daeSBarry Smith       PetscLogEventEnd(USER_EVENT,0,0,0,0);
10325c6c1daeSBarry Smith .ve
10335c6c1daeSBarry Smith 
1034d1f92df0SBarry Smith   Level: intermediate
1035d1f92df0SBarry Smith 
10365c6c1daeSBarry Smith   Notes:
10375c6c1daeSBarry Smith   PETSc automatically logs library events if the code has been
1038a2553e36SBarry Smith   configured with --with-log (which is the default) and
1039811af0c4SBarry Smith   -log_view or -log_all is specified.  `PetscLogEventRegister()` is
10405c6c1daeSBarry Smith   intended for logging user events to supplement this PETSc
10415c6c1daeSBarry Smith   information.
10425c6c1daeSBarry Smith 
1043495fc317SBarry Smith   PETSc can gather data for use with the utilities Jumpshot
10445c6c1daeSBarry Smith   (part of the MPICH distribution).  If PETSc has been compiled
10455c6c1daeSBarry Smith   with flag -DPETSC_HAVE_MPE (MPE is an additional utility within
10465c6c1daeSBarry Smith   MPICH), the user can employ another command line option, -log_mpe,
10475c6c1daeSBarry Smith   to create a logfile, "mpe.log", which can be visualized
1048495fc317SBarry Smith   Jumpshot.
10495c6c1daeSBarry Smith 
10505c6c1daeSBarry Smith   The classid is associated with each event so that classes of events
10515c6c1daeSBarry Smith   can be disabled simultaneously, such as all matrix events. The user
1052811af0c4SBarry Smith   can either use an existing classid, such as `MAT_CLASSID`, or create
10535c6c1daeSBarry Smith   their own as shown in the example.
10545c6c1daeSBarry Smith 
1055c5deb1d5SJed Brown   If an existing event with the same name exists, its event handle is
1056c5deb1d5SJed Brown   returned instead of creating a new event.
1057c5deb1d5SJed Brown 
1058d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogStageRegister()`, `PetscLogEventBegin()`, `PetscLogEventEnd()`, `PetscLogFlops()`,
1059db781477SPatrick Sanan           `PetscLogEventActivate()`, `PetscLogEventDeactivate()`, `PetscClassIdRegister()`
10605c6c1daeSBarry Smith @*/
1061d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogEventRegister(const char name[], PetscClassId classid, PetscLogEvent *event)
1062d71ae5a4SJacob Faibussowitsch {
1063b665b14eSToby Isaac   PetscLogState state;
10645c6c1daeSBarry Smith 
10655c6c1daeSBarry Smith   PetscFunctionBegin;
1066b665b14eSToby Isaac   *event = -1;
1067b665b14eSToby Isaac   PetscCall(PetscLogGetState(&state));
1068b665b14eSToby Isaac   if (state) PetscCall(PetscLogStateEventRegister(state, name, classid, event));
10693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10705c6c1daeSBarry Smith }
10715c6c1daeSBarry Smith 
10725c6c1daeSBarry Smith /*@
1073217044c2SLisandro Dalcin   PetscLogEventSetCollective - Indicates that a particular event is collective.
1074217044c2SLisandro Dalcin 
10755aefd447SStefano Zampini   Logically Collective
1076217044c2SLisandro Dalcin 
1077d8d19677SJose E. Roman   Input Parameters:
1078217044c2SLisandro Dalcin + event      - The event id
1079ffbd2f08SBarry Smith - collective - `PetscBool` indicating whether a particular event is collective
1080217044c2SLisandro Dalcin 
1081d1f92df0SBarry Smith   Level: developer
1082d1f92df0SBarry Smith 
1083811af0c4SBarry Smith   Notes:
1084811af0c4SBarry Smith   New events returned from `PetscLogEventRegister()` are collective by default.
1085811af0c4SBarry Smith 
1086ffbd2f08SBarry Smith   Collective events are handled specially if the command line option `-log_sync` is used. In that case the logging saves information about
1087811af0c4SBarry Smith   two parts of the event; the time for all the MPI ranks to synchronize and then the time for the actual computation/communication
10885aefd447SStefano Zampini   to be performed. This option is useful to debug imbalance within the computations or communications.
1089217044c2SLisandro Dalcin 
1090d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogEventBegin()`, `PetscLogEventEnd()`, `PetscLogEventRegister()`
1091217044c2SLisandro Dalcin @*/
1092d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogEventSetCollective(PetscLogEvent event, PetscBool collective)
1093d71ae5a4SJacob Faibussowitsch {
1094b665b14eSToby Isaac   PetscLogState state;
1095217044c2SLisandro Dalcin 
1096217044c2SLisandro Dalcin   PetscFunctionBegin;
1097b665b14eSToby Isaac   PetscCall(PetscLogGetState(&state));
1098b665b14eSToby Isaac   if (state) PetscCall(PetscLogStateEventSetCollective(state, event, collective));
1099b665b14eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1100b665b14eSToby Isaac }
1101b665b14eSToby Isaac 
1102b665b14eSToby Isaac /*
1103b665b14eSToby Isaac   PetscLogClassSetActiveAll - Activate or inactivate logging for all events associated with a PETSc object class in every stage.
1104b665b14eSToby Isaac 
1105b665b14eSToby Isaac   Not Collective
1106b665b14eSToby Isaac 
1107b665b14eSToby Isaac   Input Parameters:
1108b665b14eSToby Isaac + classid - The object class, for example `MAT_CLASSID`, `SNES_CLASSID`, etc.
1109b665b14eSToby Isaac - isActive - if `PETSC_FALSE`, events associated with this class will not be send to log handlers.
1110b665b14eSToby Isaac 
1111b665b14eSToby Isaac   Level: developer
1112b665b14eSToby Isaac 
1113b665b14eSToby Isaac .seealso: [](ch_profiling), `PetscLogEventActivate()`, `PetscLogEventActivateAll()`, `PetscLogStageSetActive()`, `PetscLogEventActivateClass()`
1114b665b14eSToby Isaac */
1115b665b14eSToby Isaac static PetscErrorCode PetscLogClassSetActiveAll(PetscClassId classid, PetscBool isActive)
1116b665b14eSToby Isaac {
1117b665b14eSToby Isaac   PetscLogState state;
1118b665b14eSToby Isaac 
1119b665b14eSToby Isaac   PetscFunctionBegin;
1120b665b14eSToby Isaac   PetscCall(PetscLogGetState(&state));
1121b665b14eSToby Isaac   if (state) PetscCall(PetscLogStateClassSetActiveAll(state, classid, isActive));
11223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1123217044c2SLisandro Dalcin }
1124217044c2SLisandro Dalcin 
1125217044c2SLisandro Dalcin /*@
1126fa2bb9feSLisandro Dalcin   PetscLogEventIncludeClass - Activates event logging for a PETSc object class in every stage.
1127fa2bb9feSLisandro Dalcin 
1128fa2bb9feSLisandro Dalcin   Not Collective
1129fa2bb9feSLisandro Dalcin 
1130fa2bb9feSLisandro Dalcin   Input Parameter:
1131811af0c4SBarry Smith . classid - The object class, for example `MAT_CLASSID`, `SNES_CLASSID`, etc.
1132fa2bb9feSLisandro Dalcin 
1133fa2bb9feSLisandro Dalcin   Level: developer
1134fa2bb9feSLisandro Dalcin 
1135d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogEventActivateClass()`, `PetscLogEventDeactivateClass()`, `PetscLogEventActivate()`, `PetscLogEventDeactivate()`
1136fa2bb9feSLisandro Dalcin @*/
1137d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogEventIncludeClass(PetscClassId classid)
1138d71ae5a4SJacob Faibussowitsch {
1139fa2bb9feSLisandro Dalcin   PetscFunctionBegin;
1140b665b14eSToby Isaac   PetscCall(PetscLogClassSetActiveAll(classid, PETSC_TRUE));
11413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1142fa2bb9feSLisandro Dalcin }
1143fa2bb9feSLisandro Dalcin 
1144fa2bb9feSLisandro Dalcin /*@
1145fa2bb9feSLisandro Dalcin   PetscLogEventExcludeClass - Deactivates event logging for a PETSc object class in every stage.
1146fa2bb9feSLisandro Dalcin 
1147fa2bb9feSLisandro Dalcin   Not Collective
1148fa2bb9feSLisandro Dalcin 
1149fa2bb9feSLisandro Dalcin   Input Parameter:
1150811af0c4SBarry Smith . classid - The object class, for example `MAT_CLASSID`, `SNES_CLASSID`, etc.
1151fa2bb9feSLisandro Dalcin 
1152fa2bb9feSLisandro Dalcin   Level: developer
1153fa2bb9feSLisandro Dalcin 
1154811af0c4SBarry Smith   Note:
1155811af0c4SBarry Smith   If a class is excluded then events associated with that class are not logged.
1156811af0c4SBarry Smith 
1157d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogEventDeactivateClass()`, `PetscLogEventActivateClass()`, `PetscLogEventDeactivate()`, `PetscLogEventActivate()`
1158fa2bb9feSLisandro Dalcin @*/
1159d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogEventExcludeClass(PetscClassId classid)
1160d71ae5a4SJacob Faibussowitsch {
1161b665b14eSToby Isaac   PetscFunctionBegin;
1162b665b14eSToby Isaac   PetscCall(PetscLogClassSetActiveAll(classid, PETSC_FALSE));
1163b665b14eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1164b665b14eSToby Isaac }
1165b665b14eSToby Isaac 
1166b665b14eSToby Isaac /*
1167b665b14eSToby Isaac   PetscLogEventSetActive - Activate or inactivate logging for an event in a given stage
1168b665b14eSToby Isaac 
1169b665b14eSToby Isaac   Not Collective
1170b665b14eSToby Isaac 
1171b665b14eSToby Isaac   Input Parameters:
1172b665b14eSToby Isaac + stage - A registered `PetscLogStage` (or `PETSC_DEFAULT` for the current stage)
1173b665b14eSToby Isaac . event - A `PetscLogEvent`
1174b665b14eSToby Isaac - isActive - If `PETSC_FALSE`, activity from this event (`PetscLogEventBegin()`, `PetscLogEventEnd()`, `PetscLogEventSync()`) will not be sent to log handlers during this stage
1175b665b14eSToby Isaac 
1176b665b14eSToby Isaac   Usage:
1177b665b14eSToby Isaac .vb
1178b665b14eSToby Isaac       PetscLogEventSetActive(VEC_SetValues, PETSC_FALSE);
1179b665b14eSToby Isaac         [code where you do not want to log VecSetValues()]
1180b665b14eSToby Isaac       PetscLogEventSetActive(VEC_SetValues, PETSC_TRUE);
1181b665b14eSToby Isaac         [code where you do want to log VecSetValues()]
1182b665b14eSToby Isaac .ve
1183b665b14eSToby Isaac 
1184b665b14eSToby Isaac   Level: advanced
1185b665b14eSToby Isaac 
1186b665b14eSToby Isaac   Note:
1187b665b14eSToby Isaac   The event may be either a pre-defined PETSc event (found in include/petsclog.h)
1188b665b14eSToby Isaac   or an event number obtained with `PetscLogEventRegister()`.
1189b665b14eSToby Isaac 
1190b665b14eSToby Isaac .seealso: [](ch_profiling), `PetscLogEventDeactivatePush()`, `PetscLogEventDeactivatePop()`
1191b665b14eSToby Isaac */
1192b665b14eSToby Isaac static PetscErrorCode PetscLogEventSetActive(PetscLogStage stage, PetscLogEvent event, PetscBool isActive)
1193b665b14eSToby Isaac {
1194b665b14eSToby Isaac   PetscLogState state;
1195fa2bb9feSLisandro Dalcin 
1196fa2bb9feSLisandro Dalcin   PetscFunctionBegin;
1197b665b14eSToby Isaac   PetscCall(PetscLogGetState(&state));
1198b665b14eSToby Isaac   if (state) PetscCall(PetscLogStateEventSetActive(state, stage, event, isActive));
11993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1200fa2bb9feSLisandro Dalcin }
1201fa2bb9feSLisandro Dalcin 
1202fa2bb9feSLisandro Dalcin /*@
12035c6c1daeSBarry Smith   PetscLogEventActivate - Indicates that a particular event should be logged.
12045c6c1daeSBarry Smith 
12055c6c1daeSBarry Smith   Not Collective
12065c6c1daeSBarry Smith 
12075c6c1daeSBarry Smith   Input Parameter:
12085c6c1daeSBarry Smith . event - The event id
12095c6c1daeSBarry Smith 
121010450e9eSJacob Faibussowitsch   Example Usage:
12115c6c1daeSBarry Smith .vb
12125c6c1daeSBarry Smith       PetscLogEventDeactivate(VEC_SetValues);
12135c6c1daeSBarry Smith         [code where you do not want to log VecSetValues()]
12145c6c1daeSBarry Smith       PetscLogEventActivate(VEC_SetValues);
12155c6c1daeSBarry Smith         [code where you do want to log VecSetValues()]
12165c6c1daeSBarry Smith .ve
12175c6c1daeSBarry Smith 
1218d1f92df0SBarry Smith   Level: advanced
1219d1f92df0SBarry Smith 
12205c6c1daeSBarry Smith   Note:
12215c6c1daeSBarry Smith   The event may be either a pre-defined PETSc event (found in include/petsclog.h)
1222811af0c4SBarry Smith   or an event number obtained with `PetscLogEventRegister()`.
12235c6c1daeSBarry Smith 
1224b665b14eSToby Isaac .seealso: [](ch_profiling), `PetscLogEventDeactivate()`, `PetscLogEventDeactivatePush()`, `PetscLogEventDeactivatePop()`
12255c6c1daeSBarry Smith @*/
1226d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogEventActivate(PetscLogEvent event)
1227d71ae5a4SJacob Faibussowitsch {
12285c6c1daeSBarry Smith   PetscFunctionBegin;
1229b665b14eSToby Isaac   PetscCall(PetscLogEventSetActive(PETSC_DEFAULT, event, PETSC_TRUE));
12303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12315c6c1daeSBarry Smith }
12325c6c1daeSBarry Smith 
12335c6c1daeSBarry Smith /*@
12345c6c1daeSBarry Smith   PetscLogEventDeactivate - Indicates that a particular event should not be logged.
12355c6c1daeSBarry Smith 
12365c6c1daeSBarry Smith   Not Collective
12375c6c1daeSBarry Smith 
12385c6c1daeSBarry Smith   Input Parameter:
12395c6c1daeSBarry Smith . event - The event id
12405c6c1daeSBarry Smith 
124110450e9eSJacob Faibussowitsch   Example Usage:
12425c6c1daeSBarry Smith .vb
12435c6c1daeSBarry Smith       PetscLogEventDeactivate(VEC_SetValues);
12445c6c1daeSBarry Smith         [code where you do not want to log VecSetValues()]
12455c6c1daeSBarry Smith       PetscLogEventActivate(VEC_SetValues);
12465c6c1daeSBarry Smith         [code where you do want to log VecSetValues()]
12475c6c1daeSBarry Smith .ve
12485c6c1daeSBarry Smith 
1249d1f92df0SBarry Smith   Level: advanced
1250d1f92df0SBarry Smith 
12515c6c1daeSBarry Smith   Note:
12525c6c1daeSBarry Smith   The event may be either a pre-defined PETSc event (found in
1253811af0c4SBarry Smith   include/petsclog.h) or an event number obtained with `PetscLogEventRegister()`).
12545c6c1daeSBarry Smith 
1255d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogEventActivate()`, `PetscLogEventDeactivatePush()`, `PetscLogEventDeactivatePop()`
12565c6c1daeSBarry Smith @*/
1257d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogEventDeactivate(PetscLogEvent event)
1258d71ae5a4SJacob Faibussowitsch {
12595c6c1daeSBarry Smith   PetscFunctionBegin;
1260b665b14eSToby Isaac   PetscCall(PetscLogEventSetActive(PETSC_DEFAULT, event, PETSC_FALSE));
12613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12625c6c1daeSBarry Smith }
12635c6c1daeSBarry Smith 
12645c6c1daeSBarry Smith /*@
1265811af0c4SBarry Smith   PetscLogEventDeactivatePush - Indicates that a particular event should not be logged until `PetscLogEventDeactivatePop()` is called
1266c00cb57fSBarry Smith 
1267c00cb57fSBarry Smith   Not Collective
1268c00cb57fSBarry Smith 
1269c00cb57fSBarry Smith   Input Parameter:
1270c00cb57fSBarry Smith . event - The event id
1271c00cb57fSBarry Smith 
127210450e9eSJacob Faibussowitsch   Example Usage:
1273c00cb57fSBarry Smith .vb
1274c00cb57fSBarry Smith       PetscLogEventDeactivatePush(VEC_SetValues);
1275c00cb57fSBarry Smith         [code where you do not want to log VecSetValues()]
1276c00cb57fSBarry Smith       PetscLogEventDeactivatePop(VEC_SetValues);
1277c00cb57fSBarry Smith         [code where you do want to log VecSetValues()]
1278c00cb57fSBarry Smith .ve
1279c00cb57fSBarry Smith 
1280d1f92df0SBarry Smith   Level: advanced
1281d1f92df0SBarry Smith 
1282c00cb57fSBarry Smith   Note:
1283c00cb57fSBarry Smith   The event may be either a pre-defined PETSc event (found in
1284811af0c4SBarry Smith   include/petsclog.h) or an event number obtained with `PetscLogEventRegister()`).
1285c00cb57fSBarry Smith 
1286baca6076SPierre Jolivet   PETSc's default log handler (`PetscLogDefaultBegin()`) respects this function because it can make the output of `PetscLogView()` easier to interpret, but other handlers (such as the nested handler, `PetscLogNestedBegin()`) ignore it because suppressing events is not helpful in their output formats.
1287b665b14eSToby Isaac 
12884b7c4d4dSPierre Jolivet .seealso: [](ch_profiling), `PetscLogEventActivate()`, `PetscLogEventDeactivate()`, `PetscLogEventDeactivatePop()`
1289c00cb57fSBarry Smith @*/
1290d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogEventDeactivatePush(PetscLogEvent event)
1291d71ae5a4SJacob Faibussowitsch {
1292c00cb57fSBarry Smith   PetscFunctionBegin;
1293dff009beSToby Isaac   for (PetscInt i = 0; i < PETSC_LOG_HANDLER_MAX; i++) {
1294dff009beSToby Isaac     PetscLogHandler h = PetscLogHandlers[i].handler;
1295dff009beSToby Isaac 
1296dff009beSToby Isaac     if (h) PetscCall(PetscLogHandlerEventDeactivatePush(h, PETSC_DEFAULT, event));
1297dff009beSToby Isaac   }
12983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1299c00cb57fSBarry Smith }
1300c00cb57fSBarry Smith 
1301c00cb57fSBarry Smith /*@
1302811af0c4SBarry Smith   PetscLogEventDeactivatePop - Indicates that a particular event should again be logged after the logging was turned off with `PetscLogEventDeactivatePush()`
1303c00cb57fSBarry Smith 
1304c00cb57fSBarry Smith   Not Collective
1305c00cb57fSBarry Smith 
1306c00cb57fSBarry Smith   Input Parameter:
1307c00cb57fSBarry Smith . event - The event id
1308c00cb57fSBarry Smith 
130910450e9eSJacob Faibussowitsch   Example Usage:
1310c00cb57fSBarry Smith .vb
1311c00cb57fSBarry Smith       PetscLogEventDeactivatePush(VEC_SetValues);
1312c00cb57fSBarry Smith         [code where you do not want to log VecSetValues()]
1313c00cb57fSBarry Smith       PetscLogEventDeactivatePop(VEC_SetValues);
1314c00cb57fSBarry Smith         [code where you do want to log VecSetValues()]
1315c00cb57fSBarry Smith .ve
1316c00cb57fSBarry Smith 
1317d1f92df0SBarry Smith   Level: advanced
1318d1f92df0SBarry Smith 
1319c00cb57fSBarry Smith   Note:
1320c00cb57fSBarry Smith   The event may be either a pre-defined PETSc event (found in
1321811af0c4SBarry Smith   include/petsclog.h) or an event number obtained with `PetscLogEventRegister()`).
1322c00cb57fSBarry Smith 
1323d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogEventActivate()`, `PetscLogEventDeactivatePush()`
1324c00cb57fSBarry Smith @*/
1325d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogEventDeactivatePop(PetscLogEvent event)
1326d71ae5a4SJacob Faibussowitsch {
1327c00cb57fSBarry Smith   PetscFunctionBegin;
1328dff009beSToby Isaac   for (PetscInt i = 0; i < PETSC_LOG_HANDLER_MAX; i++) {
1329dff009beSToby Isaac     PetscLogHandler h = PetscLogHandlers[i].handler;
1330dff009beSToby Isaac 
1331dff009beSToby Isaac     if (h) PetscCall(PetscLogHandlerEventDeactivatePop(h, PETSC_DEFAULT, event));
1332dff009beSToby Isaac   }
13333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1334c00cb57fSBarry Smith }
1335c00cb57fSBarry Smith 
1336c00cb57fSBarry Smith /*@
1337811af0c4SBarry Smith   PetscLogEventSetActiveAll - Turns on logging of all events
13385c6c1daeSBarry Smith 
13395c6c1daeSBarry Smith   Not Collective
13405c6c1daeSBarry Smith 
13415c6c1daeSBarry Smith   Input Parameters:
13425c6c1daeSBarry Smith + event    - The event id
13435c6c1daeSBarry Smith - isActive - The activity flag determining whether the event is logged
13445c6c1daeSBarry Smith 
13455c6c1daeSBarry Smith   Level: advanced
13465c6c1daeSBarry Smith 
1347b665b14eSToby Isaac .seealso: [](ch_profiling), `PetscLogEventActivate()`, `PetscLogEventDeactivate()`
13485c6c1daeSBarry Smith @*/
1349d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogEventSetActiveAll(PetscLogEvent event, PetscBool isActive)
1350d71ae5a4SJacob Faibussowitsch {
1351b665b14eSToby Isaac   PetscLogState state;
13525c6c1daeSBarry Smith 
13535c6c1daeSBarry Smith   PetscFunctionBegin;
1354b665b14eSToby Isaac   PetscCall(PetscLogGetState(&state));
1355b665b14eSToby Isaac   if (state) PetscCall(PetscLogStateEventSetActiveAll(state, event, isActive));
1356b665b14eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
13575c6c1daeSBarry Smith }
1358b665b14eSToby Isaac 
1359b665b14eSToby Isaac /*
1360b665b14eSToby Isaac   PetscLogClassSetActive - Activates event logging for a PETSc object class for the current stage
1361b665b14eSToby Isaac 
1362b665b14eSToby Isaac   Not Collective
1363b665b14eSToby Isaac 
1364b665b14eSToby Isaac   Input Parameters:
1365b665b14eSToby Isaac + stage - A registered `PetscLogStage` (or `PETSC_DEFAULT` for the current stage)
1366b665b14eSToby Isaac . classid - The event class, for example `MAT_CLASSID`, `SNES_CLASSID`, etc.
1367b665b14eSToby Isaac - isActive - If `PETSC_FALSE`, events associated with this class are not sent to log handlers.
1368b665b14eSToby Isaac 
1369b665b14eSToby Isaac   Level: developer
1370b665b14eSToby Isaac 
1371b665b14eSToby Isaac .seealso: [](ch_profiling), `PetscLogEventIncludeClass()`, `PetscLogEventActivate()`, `PetscLogEventActivateAll()`, `PetscLogStageSetActive()`
1372b665b14eSToby Isaac */
1373b665b14eSToby Isaac static PetscErrorCode PetscLogClassSetActive(PetscLogStage stage, PetscClassId classid, PetscBool isActive)
1374b665b14eSToby Isaac {
1375b665b14eSToby Isaac   PetscLogState state;
1376b665b14eSToby Isaac 
1377b665b14eSToby Isaac   PetscFunctionBegin;
1378b665b14eSToby Isaac   PetscCall(PetscLogGetState(&state));
1379b665b14eSToby Isaac   if (state) PetscCall(PetscLogStateClassSetActive(state, stage, classid, isActive));
13803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13815c6c1daeSBarry Smith }
13825c6c1daeSBarry Smith 
13835c6c1daeSBarry Smith /*@
1384811af0c4SBarry Smith   PetscLogEventActivateClass - Activates event logging for a PETSc object class for the current stage
13855c6c1daeSBarry Smith 
13865c6c1daeSBarry Smith   Not Collective
13875c6c1daeSBarry Smith 
13885c6c1daeSBarry Smith   Input Parameter:
1389811af0c4SBarry Smith . classid - The event class, for example `MAT_CLASSID`, `SNES_CLASSID`, etc.
13905c6c1daeSBarry Smith 
13915c6c1daeSBarry Smith   Level: developer
13925c6c1daeSBarry Smith 
1393d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogEventIncludeClass()`, `PetscLogEventExcludeClass()`, `PetscLogEventDeactivateClass()`, `PetscLogEventActivate()`, `PetscLogEventDeactivate()`
13945c6c1daeSBarry Smith @*/
1395d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogEventActivateClass(PetscClassId classid)
1396d71ae5a4SJacob Faibussowitsch {
13975c6c1daeSBarry Smith   PetscFunctionBegin;
1398b665b14eSToby Isaac   PetscCall(PetscLogClassSetActive(PETSC_DEFAULT, classid, PETSC_TRUE));
13993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14005c6c1daeSBarry Smith }
14015c6c1daeSBarry Smith 
14025c6c1daeSBarry Smith /*@
1403811af0c4SBarry Smith   PetscLogEventDeactivateClass - Deactivates event logging for a PETSc object class for the current stage
14045c6c1daeSBarry Smith 
14055c6c1daeSBarry Smith   Not Collective
14065c6c1daeSBarry Smith 
14075c6c1daeSBarry Smith   Input Parameter:
1408811af0c4SBarry Smith . classid - The event class, for example `MAT_CLASSID`, `SNES_CLASSID`, etc.
14095c6c1daeSBarry Smith 
14105c6c1daeSBarry Smith   Level: developer
14115c6c1daeSBarry Smith 
1412d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogEventIncludeClass()`, `PetscLogEventExcludeClass()`, `PetscLogEventActivateClass()`, `PetscLogEventActivate()`, `PetscLogEventDeactivate()`
14135c6c1daeSBarry Smith @*/
1414d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogEventDeactivateClass(PetscClassId classid)
1415d71ae5a4SJacob Faibussowitsch {
14165c6c1daeSBarry Smith   PetscFunctionBegin;
1417b665b14eSToby Isaac   PetscCall(PetscLogClassSetActive(PETSC_DEFAULT, classid, PETSC_FALSE));
14183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14195c6c1daeSBarry Smith }
14205c6c1daeSBarry Smith 
14215c6c1daeSBarry Smith /*MC
142262872c28SLisandro Dalcin   PetscLogEventSync - Synchronizes the beginning of a user event.
142362872c28SLisandro Dalcin 
142462872c28SLisandro Dalcin   Synopsis:
142562872c28SLisandro Dalcin   #include <petsclog.h>
1426b665b14eSToby Isaac   PetscErrorCode PetscLogEventSync(PetscLogEvent e, MPI_Comm comm)
142762872c28SLisandro Dalcin 
142862872c28SLisandro Dalcin   Collective
142962872c28SLisandro Dalcin 
143062872c28SLisandro Dalcin   Input Parameters:
1431b665b14eSToby Isaac + e    - `PetscLogEvent` obtained from `PetscLogEventRegister()`
143262872c28SLisandro Dalcin - comm - an MPI communicator
143362872c28SLisandro Dalcin 
143410450e9eSJacob Faibussowitsch   Example Usage:
143562872c28SLisandro Dalcin .vb
143662872c28SLisandro Dalcin   PetscLogEvent USER_EVENT;
143710450e9eSJacob Faibussowitsch 
143862872c28SLisandro Dalcin   PetscLogEventRegister("User event", 0, &USER_EVENT);
143962872c28SLisandro Dalcin   PetscLogEventSync(USER_EVENT, PETSC_COMM_WORLD);
144062872c28SLisandro Dalcin   PetscLogEventBegin(USER_EVENT, 0, 0, 0, 0);
144162872c28SLisandro Dalcin   [code segment to monitor]
144262872c28SLisandro Dalcin   PetscLogEventEnd(USER_EVENT, 0, 0, 0 , 0);
144362872c28SLisandro Dalcin .ve
144462872c28SLisandro Dalcin 
1445d1f92df0SBarry Smith   Level: developer
1446d1f92df0SBarry Smith 
1447811af0c4SBarry Smith   Note:
144810450e9eSJacob Faibussowitsch   This routine should be called only if there is not a `PetscObject` available to pass to
144910450e9eSJacob Faibussowitsch   `PetscLogEventBegin()`.
145062872c28SLisandro Dalcin 
1451d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogEventRegister()`, `PetscLogEventBegin()`, `PetscLogEventEnd()`
145262872c28SLisandro Dalcin M*/
145362872c28SLisandro Dalcin 
145462872c28SLisandro Dalcin /*MC
14555c6c1daeSBarry Smith   PetscLogEventBegin - Logs the beginning of a user event.
14565c6c1daeSBarry Smith 
14575c6c1daeSBarry Smith   Synopsis:
1458aaa7dc30SBarry Smith   #include <petsclog.h>
1459b665b14eSToby Isaac   PetscErrorCode PetscLogEventBegin(PetscLogEvent e, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
14605c6c1daeSBarry Smith 
14615c6c1daeSBarry Smith   Not Collective
14625c6c1daeSBarry Smith 
14635c6c1daeSBarry Smith   Input Parameters:
1464b665b14eSToby Isaac + e  - `PetscLogEvent` obtained from `PetscLogEventRegister()`
1465baca6076SPierre Jolivet . o1 - object associated with the event, or `NULL`
1466baca6076SPierre Jolivet . o2 - object associated with the event, or `NULL`
1467baca6076SPierre Jolivet . o3 - object associated with the event, or `NULL`
1468baca6076SPierre Jolivet - o4 - object associated with the event, or `NULL`
14695c6c1daeSBarry Smith 
14705c6c1daeSBarry Smith   Fortran Synopsis:
14715c6c1daeSBarry Smith   void PetscLogEventBegin(int e, PetscErrorCode ierr)
14725c6c1daeSBarry Smith 
147310450e9eSJacob Faibussowitsch   Example Usage:
14745c6c1daeSBarry Smith .vb
14755c6c1daeSBarry Smith   PetscLogEvent USER_EVENT;
147610450e9eSJacob Faibussowitsch 
14775c6c1daeSBarry Smith   PetscLogDouble user_event_flops;
14785c6c1daeSBarry Smith   PetscLogEventRegister("User event",0, &USER_EVENT);
14795c6c1daeSBarry Smith   PetscLogEventBegin(USER_EVENT, 0, 0, 0, 0);
14805c6c1daeSBarry Smith   [code segment to monitor]
14815c6c1daeSBarry Smith   PetscLogFlops(user_event_flops);
14825c6c1daeSBarry Smith   PetscLogEventEnd(USER_EVENT, 0, 0, 0, 0);
14835c6c1daeSBarry Smith .ve
14845c6c1daeSBarry Smith 
1485d1f92df0SBarry Smith   Level: intermediate
1486d1f92df0SBarry Smith 
1487811af0c4SBarry Smith   Developer Note:
148810450e9eSJacob Faibussowitsch   `PetscLogEventBegin()` and `PetscLogEventBegin()` return error codes instead of explicitly
148910450e9eSJacob Faibussowitsch   handling the errors that occur in the macro directly because other packages that use this
149010450e9eSJacob Faibussowitsch   macros have used them in their own functions or methods that do not return error codes and it
149110450e9eSJacob Faibussowitsch   would be disruptive to change the current behavior.
1492d0609cedSBarry Smith 
1493d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogEventRegister()`, `PetscLogEventEnd()`, `PetscLogFlops()`
14945c6c1daeSBarry Smith M*/
14955c6c1daeSBarry Smith 
14965c6c1daeSBarry Smith /*MC
14975c6c1daeSBarry Smith   PetscLogEventEnd - Log the end of a user event.
14985c6c1daeSBarry Smith 
14995c6c1daeSBarry Smith   Synopsis:
1500aaa7dc30SBarry Smith   #include <petsclog.h>
1501b665b14eSToby Isaac   PetscErrorCode PetscLogEventEnd(PetscLogEvent e, PetscObject o1, PetscObject o2, PetscObject o3, PetscObject o4)
15025c6c1daeSBarry Smith 
15035c6c1daeSBarry Smith   Not Collective
15045c6c1daeSBarry Smith 
15055c6c1daeSBarry Smith   Input Parameters:
1506b665b14eSToby Isaac + e  - `PetscLogEvent` obtained from `PetscLogEventRegister()`
1507baca6076SPierre Jolivet . o1 - object associated with the event, or `NULL`
1508baca6076SPierre Jolivet . o2 - object associated with the event, or `NULL`
1509baca6076SPierre Jolivet . o3 - object associated with the event, or `NULL`
1510baca6076SPierre Jolivet - o4 - object associated with the event, or `NULL`
15115c6c1daeSBarry Smith 
15125c6c1daeSBarry Smith   Fortran Synopsis:
15135c6c1daeSBarry Smith   void PetscLogEventEnd(int e, PetscErrorCode ierr)
15145c6c1daeSBarry Smith 
151510450e9eSJacob Faibussowitsch   Example Usage:
15165c6c1daeSBarry Smith .vb
15175c6c1daeSBarry Smith   PetscLogEvent USER_EVENT;
151810450e9eSJacob Faibussowitsch 
15195c6c1daeSBarry Smith   PetscLogDouble user_event_flops;
152010450e9eSJacob Faibussowitsch   PetscLogEventRegister("User event", 0, &USER_EVENT);
15215c6c1daeSBarry Smith   PetscLogEventBegin(USER_EVENT, 0, 0, 0, 0);
15225c6c1daeSBarry Smith   [code segment to monitor]
15235c6c1daeSBarry Smith   PetscLogFlops(user_event_flops);
15245c6c1daeSBarry Smith   PetscLogEventEnd(USER_EVENT, 0, 0, 0, 0);
15255c6c1daeSBarry Smith .ve
15265c6c1daeSBarry Smith 
15275c6c1daeSBarry Smith   Level: intermediate
15285c6c1daeSBarry Smith 
1529d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogEventRegister()`, `PetscLogEventBegin()`, `PetscLogFlops()`
15305c6c1daeSBarry Smith M*/
15315c6c1daeSBarry Smith 
15325c6c1daeSBarry Smith /*@C
15338b08f494SToby Isaac   PetscLogStageGetPerfInfo - Return the performance information about the given stage
15348b08f494SToby Isaac 
1535cc4c1da9SBarry Smith   No Fortran Support
1536cc4c1da9SBarry Smith 
15378b08f494SToby Isaac   Input Parameters:
15388b08f494SToby Isaac . stage - The stage number or `PETSC_DETERMINE` for the current stage
15398b08f494SToby Isaac 
15408b08f494SToby Isaac   Output Parameter:
15418b08f494SToby Isaac . info - This structure is filled with the performance information
15428b08f494SToby Isaac 
15438b08f494SToby Isaac   Level: intermediate
15448b08f494SToby Isaac 
15458b08f494SToby Isaac   Notes:
15468b08f494SToby Isaac   This is a low level routine used by the logging functions in PETSc.
15478b08f494SToby Isaac 
15488b08f494SToby Isaac   A `PETSCLOGHANDLERDEFAULT` must be running for this to work, having been started either with
15495804573cSPierre Jolivet   `PetscLogDefaultBegin()` or from the command line with `-log_view`.  If it was not started,
15500ba5dad9SToby Isaac   all performance statistics in `info` will be zeroed.
15518b08f494SToby Isaac 
15528b08f494SToby Isaac .seealso: [](ch_profiling), `PetscLogEventRegister()`, `PetscLogEventBegin()`, `PetscLogEventEnd()`, `PetscLogGetDefaultHandler()`
15538b08f494SToby Isaac @*/
15548b08f494SToby Isaac PetscErrorCode PetscLogStageGetPerfInfo(PetscLogStage stage, PetscEventPerfInfo *info)
15558b08f494SToby Isaac {
15568b08f494SToby Isaac   PetscLogHandler     handler;
15578b08f494SToby Isaac   PetscEventPerfInfo *event_info;
15588b08f494SToby Isaac 
15598b08f494SToby Isaac   PetscFunctionBegin;
15608b08f494SToby Isaac   PetscAssertPointer(info, 2);
15610ba5dad9SToby Isaac   PetscCall(PetscLogTryGetHandler(PETSCLOGHANDLERDEFAULT, &handler));
15620ba5dad9SToby Isaac   if (handler) {
15638b08f494SToby Isaac     PetscCall(PetscLogHandlerGetStagePerfInfo(handler, stage, &event_info));
15648b08f494SToby Isaac     *info = *event_info;
15650ba5dad9SToby Isaac   } else {
15660ba5dad9SToby Isaac     PetscCall(PetscInfo(NULL, "Default log handler is not running, PetscLogStageGetPerfInfo() returning zeros\n"));
15670ba5dad9SToby Isaac     PetscCall(PetscMemzero(info, sizeof(*info)));
15680ba5dad9SToby Isaac   }
15698b08f494SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
15708b08f494SToby Isaac }
15718b08f494SToby Isaac 
15728b08f494SToby Isaac /*@C
1573b665b14eSToby Isaac   PetscLogEventGetPerfInfo - Return the performance information about the given event in the given stage
1574b665b14eSToby Isaac 
1575cc4c1da9SBarry Smith   No Fortran Support
1576cc4c1da9SBarry Smith 
1577b665b14eSToby Isaac   Input Parameters:
1578b665b14eSToby Isaac + stage - The stage number or `PETSC_DETERMINE` for the current stage
1579b665b14eSToby Isaac - event - The event number
1580b665b14eSToby Isaac 
1581b665b14eSToby Isaac   Output Parameter:
1582b665b14eSToby Isaac . info - This structure is filled with the performance information
1583b665b14eSToby Isaac 
1584b665b14eSToby Isaac   Level: intermediate
1585b665b14eSToby Isaac 
1586b665b14eSToby Isaac   Note:
1587b665b14eSToby Isaac   This is a low level routine used by the logging functions in PETSc
1588b665b14eSToby Isaac 
15890ba5dad9SToby Isaac   A `PETSCLOGHANDLERDEFAULT` must be running for this to work, having been started either with
15905804573cSPierre Jolivet   `PetscLogDefaultBegin()` or from the command line with `-log_view`.  If it was not started,
15910ba5dad9SToby Isaac   all performance statistics in `info` will be zeroed.
15920ba5dad9SToby Isaac 
1593b665b14eSToby Isaac .seealso: [](ch_profiling), `PetscLogEventRegister()`, `PetscLogEventBegin()`, `PetscLogEventEnd()`, `PetscLogGetDefaultHandler()`
1594b665b14eSToby Isaac @*/
1595b665b14eSToby Isaac PetscErrorCode PetscLogEventGetPerfInfo(PetscLogStage stage, PetscLogEvent event, PetscEventPerfInfo *info)
1596b665b14eSToby Isaac {
1597b665b14eSToby Isaac   PetscLogHandler     handler;
1598b665b14eSToby Isaac   PetscEventPerfInfo *event_info;
1599b665b14eSToby Isaac 
1600b665b14eSToby Isaac   PetscFunctionBegin;
1601b665b14eSToby Isaac   PetscAssertPointer(info, 3);
16020ba5dad9SToby Isaac   PetscCall(PetscLogTryGetHandler(PETSCLOGHANDLERDEFAULT, &handler));
16030ba5dad9SToby Isaac   if (handler) {
1604dff009beSToby Isaac     PetscCall(PetscLogHandlerGetEventPerfInfo(handler, stage, event, &event_info));
1605b665b14eSToby Isaac     *info = *event_info;
16060ba5dad9SToby Isaac   } else {
16070ba5dad9SToby Isaac     PetscCall(PetscInfo(NULL, "Default log handler is not running, PetscLogEventGetPerfInfo() returning zeros\n"));
16080ba5dad9SToby Isaac     PetscCall(PetscMemzero(info, sizeof(*info)));
16090ba5dad9SToby Isaac   }
1610b665b14eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1611b665b14eSToby Isaac }
1612b665b14eSToby Isaac 
1613cc4c1da9SBarry Smith /*@
1614b665b14eSToby Isaac   PetscLogEventSetDof - Set the nth number of degrees of freedom of a numerical problem associated with this event
1615b665b14eSToby Isaac 
1616b665b14eSToby Isaac   Not Collective
1617b665b14eSToby Isaac 
1618b665b14eSToby Isaac   Input Parameters:
1619b665b14eSToby Isaac + event - The event id to log
1620b665b14eSToby Isaac . n     - The dof index, in [0, 8)
1621b665b14eSToby Isaac - dof   - The number of dofs
1622b665b14eSToby Isaac 
1623b665b14eSToby Isaac   Options Database Key:
1624b665b14eSToby Isaac . -log_view - Activates log summary
1625b665b14eSToby Isaac 
1626b665b14eSToby Isaac   Level: developer
1627b665b14eSToby Isaac 
1628b665b14eSToby Isaac   Note:
1629b665b14eSToby Isaac   This is to enable logging of convergence
1630b665b14eSToby Isaac 
1631b665b14eSToby Isaac .seealso: `PetscLogEventSetError()`, `PetscLogEventRegister()`, `PetscLogGetDefaultHandler()`
1632b665b14eSToby Isaac @*/
1633b665b14eSToby Isaac PetscErrorCode PetscLogEventSetDof(PetscLogEvent event, PetscInt n, PetscLogDouble dof)
1634b665b14eSToby Isaac {
1635b665b14eSToby Isaac   PetscFunctionBegin;
1636b665b14eSToby Isaac   PetscCheck(!(n < 0) && !(n > 7), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Error index %" PetscInt_FMT " is not in [0, 8)", n);
1637dff009beSToby Isaac   for (PetscInt i = 0; i < PETSC_LOG_HANDLER_MAX; i++) {
1638dff009beSToby Isaac     PetscLogHandler h = PetscLogHandlers[i].handler;
1639dff009beSToby Isaac 
1640dff009beSToby Isaac     if (h) {
1641dff009beSToby Isaac       PetscEventPerfInfo *event_info;
1642dff009beSToby Isaac 
1643dff009beSToby Isaac       PetscCall(PetscLogHandlerGetEventPerfInfo(h, PETSC_DEFAULT, event, &event_info));
1644dff009beSToby Isaac       if (event_info) event_info->dof[n] = dof;
1645dff009beSToby Isaac     }
1646b665b14eSToby Isaac   }
1647b665b14eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1648b665b14eSToby Isaac }
1649b665b14eSToby Isaac 
1650cc4c1da9SBarry Smith /*@
1651b665b14eSToby Isaac   PetscLogEventSetError - Set the nth error associated with a numerical problem associated with this event
1652b665b14eSToby Isaac 
1653b665b14eSToby Isaac   Not Collective
1654b665b14eSToby Isaac 
1655b665b14eSToby Isaac   Input Parameters:
1656b665b14eSToby Isaac + event - The event id to log
1657b665b14eSToby Isaac . n     - The error index, in [0, 8)
1658b665b14eSToby Isaac - error - The error
1659b665b14eSToby Isaac 
1660b665b14eSToby Isaac   Options Database Key:
1661b665b14eSToby Isaac . -log_view - Activates log summary
1662b665b14eSToby Isaac 
1663b665b14eSToby Isaac   Level: developer
1664b665b14eSToby Isaac 
1665b665b14eSToby Isaac   Notes:
1666b665b14eSToby Isaac   This is to enable logging of convergence, and enable users to interpret the errors as they wish. For example,
1667b665b14eSToby Isaac   as different norms, or as errors for different fields
1668b665b14eSToby Isaac 
1669b665b14eSToby Isaac   This is a low level routine used by the logging functions in PETSc
1670b665b14eSToby Isaac 
1671b665b14eSToby Isaac .seealso: `PetscLogEventSetDof()`, `PetscLogEventRegister()`, `PetscLogGetDefaultHandler()`
1672b665b14eSToby Isaac @*/
1673b665b14eSToby Isaac PetscErrorCode PetscLogEventSetError(PetscLogEvent event, PetscInt n, PetscLogDouble error)
1674b665b14eSToby Isaac {
1675b665b14eSToby Isaac   PetscFunctionBegin;
1676b665b14eSToby Isaac   PetscCheck(!(n < 0) && !(n > 7), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Error index %" PetscInt_FMT " is not in [0, 8)", n);
1677dff009beSToby Isaac   for (PetscInt i = 0; i < PETSC_LOG_HANDLER_MAX; i++) {
1678dff009beSToby Isaac     PetscLogHandler h = PetscLogHandlers[i].handler;
1679dff009beSToby Isaac 
1680dff009beSToby Isaac     if (h) {
1681dff009beSToby Isaac       PetscEventPerfInfo *event_info;
1682dff009beSToby Isaac 
1683dff009beSToby Isaac       PetscCall(PetscLogHandlerGetEventPerfInfo(h, PETSC_DEFAULT, event, &event_info));
1684dff009beSToby Isaac       if (event_info) event_info->errors[n] = error;
1685dff009beSToby Isaac     }
1686b665b14eSToby Isaac   }
1687b665b14eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1688b665b14eSToby Isaac }
1689b665b14eSToby Isaac 
1690cc4c1da9SBarry Smith /*@
16915c6c1daeSBarry Smith   PetscLogEventGetId - Returns the event id when given the event name.
16925c6c1daeSBarry Smith 
16935c6c1daeSBarry Smith   Not Collective
16945c6c1daeSBarry Smith 
16955c6c1daeSBarry Smith   Input Parameter:
16965c6c1daeSBarry Smith . name - The event name
16975c6c1daeSBarry Smith 
16985c6c1daeSBarry Smith   Output Parameter:
1699c5deb1d5SJed Brown . event - The event, or -1 if no event with that name exists
17005c6c1daeSBarry Smith 
17015c6c1daeSBarry Smith   Level: intermediate
17025c6c1daeSBarry Smith 
1703d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogEventBegin()`, `PetscLogEventEnd()`, `PetscLogStageGetId()`
17045c6c1daeSBarry Smith @*/
1705d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogEventGetId(const char name[], PetscLogEvent *event)
1706d71ae5a4SJacob Faibussowitsch {
1707b665b14eSToby Isaac   PetscLogState state;
17085c6c1daeSBarry Smith 
17095c6c1daeSBarry Smith   PetscFunctionBegin;
1710b665b14eSToby Isaac   *event = -1;
1711b665b14eSToby Isaac   PetscCall(PetscLogGetState(&state));
1712b665b14eSToby Isaac   if (state) PetscCall(PetscLogStateGetEventFromName(state, name, event));
17133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17145c6c1daeSBarry Smith }
17155c6c1daeSBarry Smith 
1716cc4c1da9SBarry Smith /*@
171753e0a2f3SToby Isaac   PetscLogEventGetName - Returns the event name when given the event id.
171853e0a2f3SToby Isaac 
171953e0a2f3SToby Isaac   Not Collective
172053e0a2f3SToby Isaac 
172153e0a2f3SToby Isaac   Input Parameter:
172253e0a2f3SToby Isaac . event - The event
172353e0a2f3SToby Isaac 
172453e0a2f3SToby Isaac   Output Parameter:
172553e0a2f3SToby Isaac . name - The event name
172653e0a2f3SToby Isaac 
172753e0a2f3SToby Isaac   Level: intermediate
172853e0a2f3SToby Isaac 
172953e0a2f3SToby Isaac .seealso: [](ch_profiling), `PetscLogEventRegister()`, `PetscLogEventBegin()`, `PetscLogEventEnd()`, `PetscPreLoadBegin()`, `PetscPreLoadEnd()`, `PetscPreLoadStage()`
173053e0a2f3SToby Isaac @*/
1731cc4c1da9SBarry Smith PetscErrorCode PetscLogEventGetName(PetscLogEvent event, const char *name[])
173253e0a2f3SToby Isaac {
173353e0a2f3SToby Isaac   PetscLogEventInfo event_info;
173453e0a2f3SToby Isaac   PetscLogState     state;
173553e0a2f3SToby Isaac 
173653e0a2f3SToby Isaac   PetscFunctionBegin;
1737b665b14eSToby Isaac   *name = NULL;
173853e0a2f3SToby Isaac   PetscCall(PetscLogGetState(&state));
1739b665b14eSToby Isaac   if (!state) PetscFunctionReturn(PETSC_SUCCESS);
174053e0a2f3SToby Isaac   PetscCall(PetscLogStateEventGetInfo(state, event, &event_info));
174153e0a2f3SToby Isaac   *name = event_info.name;
174253e0a2f3SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
174353e0a2f3SToby Isaac }
174453e0a2f3SToby Isaac 
174553e0a2f3SToby Isaac /*@
174653e0a2f3SToby Isaac   PetscLogEventsPause - Put event logging into "paused" mode: timers and counters for in-progress events are paused, and any events that happen before logging is resumed with `PetscLogEventsResume()` are logged in the "Main Stage" of execution.
174753e0a2f3SToby Isaac 
174853e0a2f3SToby Isaac   Not collective
174953e0a2f3SToby Isaac 
175053e0a2f3SToby Isaac   Level: advanced
175153e0a2f3SToby Isaac 
175253e0a2f3SToby Isaac   Notes:
175353e0a2f3SToby Isaac   When an external library or runtime has is initialized it can involve lots of setup time that skews the statistics of any unrelated running events: this function is intended to isolate such calls in the default log summary (`PetscLogDefaultBegin()`, `PetscLogView()`).
175453e0a2f3SToby Isaac 
175553e0a2f3SToby Isaac   Other log handlers (such as the nested handler, `PetscLogNestedBegin()`) will ignore this function.
175653e0a2f3SToby Isaac 
1757b665b14eSToby Isaac .seealso: [](ch_profiling), `PetscLogEventDeactivatePush()`, `PetscLogEventDeactivatePop()`, `PetscLogEventsResume()`, `PetscLogGetDefaultHandler()`
175853e0a2f3SToby Isaac @*/
175953e0a2f3SToby Isaac PetscErrorCode PetscLogEventsPause(void)
176053e0a2f3SToby Isaac {
176153e0a2f3SToby Isaac   PetscFunctionBegin;
1762dff009beSToby Isaac   for (PetscInt i = 0; i < PETSC_LOG_HANDLER_MAX; i++) {
1763dff009beSToby Isaac     PetscLogHandler h = PetscLogHandlers[i].handler;
1764dff009beSToby Isaac 
1765dff009beSToby Isaac     if (h) PetscCall(PetscLogHandlerEventsPause(h));
1766dff009beSToby Isaac   }
176753e0a2f3SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
176853e0a2f3SToby Isaac }
176953e0a2f3SToby Isaac 
177053e0a2f3SToby Isaac /*@
177153e0a2f3SToby Isaac   PetscLogEventsResume - Return logging to normal behavior after it was paused with `PetscLogEventsPause()`.
177253e0a2f3SToby Isaac 
177353e0a2f3SToby Isaac   Not collective
177453e0a2f3SToby Isaac 
177553e0a2f3SToby Isaac   Level: advanced
177653e0a2f3SToby Isaac 
1777b665b14eSToby Isaac .seealso: [](ch_profiling), `PetscLogEventDeactivatePush()`, `PetscLogEventDeactivatePop()`, `PetscLogEventsPause()`, `PetscLogGetDefaultHandler()`
177853e0a2f3SToby Isaac @*/
177953e0a2f3SToby Isaac PetscErrorCode PetscLogEventsResume(void)
178053e0a2f3SToby Isaac {
178153e0a2f3SToby Isaac   PetscFunctionBegin;
1782dff009beSToby Isaac   for (PetscInt i = 0; i < PETSC_LOG_HANDLER_MAX; i++) {
1783dff009beSToby Isaac     PetscLogHandler h = PetscLogHandlers[i].handler;
1784dff009beSToby Isaac 
1785dff009beSToby Isaac     if (h) PetscCall(PetscLogHandlerEventsResume(h));
1786dff009beSToby Isaac   }
178753e0a2f3SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
178853e0a2f3SToby Isaac }
178953e0a2f3SToby Isaac 
17901c1ad86eSToby Isaac /*------------------------------------------------ Class Functions --------------------------------------------------*/
17911c1ad86eSToby Isaac 
17921c1ad86eSToby Isaac /*MC
17931c1ad86eSToby Isaac    PetscLogObjectCreate - Log the creation of a `PetscObject`
17941c1ad86eSToby Isaac 
17951c1ad86eSToby Isaac    Synopsis:
17961c1ad86eSToby Isaac    #include <petsclog.h>
17971c1ad86eSToby Isaac    PetscErrorCode PetscLogObjectCreate(PetscObject h)
17981c1ad86eSToby Isaac 
17991c1ad86eSToby Isaac    Not Collective
18001c1ad86eSToby Isaac 
18011c1ad86eSToby Isaac    Input Parameters:
18021c1ad86eSToby Isaac .  h - A `PetscObject`
18031c1ad86eSToby Isaac 
18041c1ad86eSToby Isaac    Level: developer
18051c1ad86eSToby Isaac 
18061c1ad86eSToby Isaac    Developer Note:
18071c1ad86eSToby Isaac      Called internally by PETSc when creating objects: users do not need to call this directly.
1808b665b14eSToby Isaac      Notification of the object creation is sent to each `PetscLogHandler` that is running.
18091c1ad86eSToby Isaac 
1810b665b14eSToby Isaac .seealso: [](ch_profiling), `PetscLogHandler`, `PetscLogObjectDestroy()`
18111c1ad86eSToby Isaac M*/
18121c1ad86eSToby Isaac 
18131c1ad86eSToby Isaac /*MC
18141c1ad86eSToby Isaac    PetscLogObjectDestroy - Logs the destruction of a `PetscObject`
18151c1ad86eSToby Isaac 
18161c1ad86eSToby Isaac    Synopsis:
18171c1ad86eSToby Isaac    #include <petsclog.h>
18181c1ad86eSToby Isaac    PetscErrorCode PetscLogObjectDestroy(PetscObject h)
18191c1ad86eSToby Isaac 
18201c1ad86eSToby Isaac    Not Collective
18211c1ad86eSToby Isaac 
18221c1ad86eSToby Isaac    Input Parameters:
18231c1ad86eSToby Isaac .  h - A `PetscObject`
18241c1ad86eSToby Isaac 
18251c1ad86eSToby Isaac    Level: developer
18261c1ad86eSToby Isaac 
18271c1ad86eSToby Isaac    Developer Note:
18281c1ad86eSToby Isaac      Called internally by PETSc when destroying objects: users do not need to call this directly.
1829b665b14eSToby Isaac      Notification of the object creation is sent to each `PetscLogHandler` that is running.
18301c1ad86eSToby Isaac 
1831b665b14eSToby Isaac .seealso: [](ch_profiling), `PetscLogHandler`, `PetscLogObjectCreate()`
18321c1ad86eSToby Isaac M*/
18331c1ad86eSToby Isaac 
1834cc4c1da9SBarry Smith /*@
183553e0a2f3SToby Isaac   PetscLogClassGetClassId - Returns the `PetscClassId` when given the class name.
183653e0a2f3SToby Isaac 
183753e0a2f3SToby Isaac   Not Collective
183853e0a2f3SToby Isaac 
183953e0a2f3SToby Isaac   Input Parameter:
184053e0a2f3SToby Isaac . name - The class name
184153e0a2f3SToby Isaac 
184253e0a2f3SToby Isaac   Output Parameter:
184353e0a2f3SToby Isaac . classid - The `PetscClassId` id, or -1 if no class with that name exists
184453e0a2f3SToby Isaac 
184553e0a2f3SToby Isaac   Level: intermediate
184653e0a2f3SToby Isaac 
184753e0a2f3SToby Isaac .seealso: [](ch_profiling), `PetscLogEventBegin()`, `PetscLogEventEnd()`, `PetscLogStageGetId()`
184853e0a2f3SToby Isaac @*/
184953e0a2f3SToby Isaac PetscErrorCode PetscLogClassGetClassId(const char name[], PetscClassId *classid)
185053e0a2f3SToby Isaac {
185153e0a2f3SToby Isaac   PetscLogClass     log_class;
185253e0a2f3SToby Isaac   PetscLogClassInfo class_info;
185353e0a2f3SToby Isaac   PetscLogState     state;
185453e0a2f3SToby Isaac 
185553e0a2f3SToby Isaac   PetscFunctionBegin;
1856b665b14eSToby Isaac   *classid = -1;
185753e0a2f3SToby Isaac   PetscCall(PetscLogGetState(&state));
1858b665b14eSToby Isaac   if (!state) PetscFunctionReturn(PETSC_SUCCESS);
185953e0a2f3SToby Isaac   PetscCall(PetscLogStateGetClassFromName(state, name, &log_class));
186053e0a2f3SToby Isaac   if (log_class < 0) {
186153e0a2f3SToby Isaac     *classid = -1;
186253e0a2f3SToby Isaac     PetscFunctionReturn(PETSC_SUCCESS);
186353e0a2f3SToby Isaac   }
186453e0a2f3SToby Isaac   PetscCall(PetscLogStateClassGetInfo(state, log_class, &class_info));
186553e0a2f3SToby Isaac   *classid = class_info.classid;
186653e0a2f3SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
186753e0a2f3SToby Isaac }
186853e0a2f3SToby Isaac 
186953e0a2f3SToby Isaac /*@C
187053e0a2f3SToby Isaac   PetscLogClassIdGetName - Returns a `PetscClassId`'s name.
187153e0a2f3SToby Isaac 
187253e0a2f3SToby Isaac   Not Collective
187353e0a2f3SToby Isaac 
187453e0a2f3SToby Isaac   Input Parameter:
187553e0a2f3SToby Isaac . classid - A `PetscClassId`
187653e0a2f3SToby Isaac 
187753e0a2f3SToby Isaac   Output Parameter:
187853e0a2f3SToby Isaac . name - The class name
187953e0a2f3SToby Isaac 
188053e0a2f3SToby Isaac   Level: intermediate
188153e0a2f3SToby Isaac 
188253e0a2f3SToby Isaac .seealso: [](ch_profiling), `PetscLogClassRegister()`, `PetscLogClassBegin()`, `PetscLogClassEnd()`, `PetscPreLoadBegin()`, `PetscPreLoadEnd()`, `PetscPreLoadClass()`
188353e0a2f3SToby Isaac @*/
188453e0a2f3SToby Isaac PetscErrorCode PetscLogClassIdGetName(PetscClassId classid, const char **name)
188553e0a2f3SToby Isaac {
188653e0a2f3SToby Isaac   PetscLogClass     log_class;
188753e0a2f3SToby Isaac   PetscLogClassInfo class_info;
188853e0a2f3SToby Isaac   PetscLogState     state;
188953e0a2f3SToby Isaac 
189053e0a2f3SToby Isaac   PetscFunctionBegin;
189153e0a2f3SToby Isaac   PetscCall(PetscLogGetState(&state));
189253e0a2f3SToby Isaac   PetscCall(PetscLogStateGetClassFromClassId(state, classid, &log_class));
189353e0a2f3SToby Isaac   PetscCall(PetscLogStateClassGetInfo(state, log_class, &class_info));
189453e0a2f3SToby Isaac   *name = class_info.name;
189553e0a2f3SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
189653e0a2f3SToby Isaac }
189753e0a2f3SToby Isaac 
18985c6c1daeSBarry Smith /*------------------------------------------------ Output Functions -------------------------------------------------*/
1899cc4c1da9SBarry Smith /*@
19005c6c1daeSBarry Smith   PetscLogDump - Dumps logs of objects to a file. This file is intended to
19015c6c1daeSBarry Smith   be read by bin/petscview. This program no longer exists.
19025c6c1daeSBarry Smith 
1903811af0c4SBarry Smith   Collective on `PETSC_COMM_WORLD`
19045c6c1daeSBarry Smith 
19055c6c1daeSBarry Smith   Input Parameter:
1906aec76313SJacob Faibussowitsch . sname - an optional file name
19075c6c1daeSBarry Smith 
190810450e9eSJacob Faibussowitsch   Example Usage:
19095c6c1daeSBarry Smith .vb
19105c6c1daeSBarry Smith   PetscInitialize(...);
1911b665b14eSToby Isaac   PetscLogDefaultBegin();
1912b665b14eSToby Isaac   // ... code ...
19135c6c1daeSBarry Smith   PetscLogDump(filename);
19145c6c1daeSBarry Smith   PetscFinalize();
19155c6c1daeSBarry Smith .ve
19165c6c1daeSBarry Smith 
1917d1f92df0SBarry Smith   Level: advanced
1918d1f92df0SBarry Smith 
1919811af0c4SBarry Smith   Note:
192037fdd005SBarry Smith   The default file name is Log.<rank> where <rank> is the MPI process rank. If no name is specified,
19215c6c1daeSBarry Smith   this file will be used.
19225c6c1daeSBarry Smith 
1923b665b14eSToby Isaac .seealso: [](ch_profiling), `PetscLogDefaultBegin()`, `PetscLogView()`, `PetscLogGetDefaultHandler()`
19245c6c1daeSBarry Smith @*/
1925d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogDump(const char sname[])
1926d71ae5a4SJacob Faibussowitsch {
1927b665b14eSToby Isaac   PetscLogHandler handler;
19285c6c1daeSBarry Smith 
19295c6c1daeSBarry Smith   PetscFunctionBegin;
1930294de794SToby Isaac   PetscCall(PetscLogGetHandler(PETSCLOGHANDLERDEFAULT, &handler));
1931dff009beSToby Isaac   PetscCall(PetscLogHandlerDump(handler, sname));
1932b665b14eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
19335c6c1daeSBarry Smith }
1934b665b14eSToby Isaac 
1935cc4c1da9SBarry Smith /*@
1936b665b14eSToby Isaac   PetscLogMPEDump - Dumps the MPE logging info to file for later use with Jumpshot.
1937b665b14eSToby Isaac 
19388f14a041SBarry Smith   Collective on `PETSC_COMM_WORLD`
1939b665b14eSToby Isaac 
1940b665b14eSToby Isaac   Input Parameter:
1941b665b14eSToby Isaac . sname - filename for the MPE logfile
1942b665b14eSToby Isaac 
1943b665b14eSToby Isaac   Level: advanced
1944b665b14eSToby Isaac 
1945b665b14eSToby Isaac .seealso: [](ch_profiling), `PetscLogDump()`, `PetscLogMPEBegin()`
1946b665b14eSToby Isaac @*/
1947b665b14eSToby Isaac PetscErrorCode PetscLogMPEDump(const char sname[])
1948b665b14eSToby Isaac {
1949b665b14eSToby Isaac   PetscFunctionBegin;
1950b665b14eSToby Isaac   #if defined(PETSC_HAVE_MPE)
1951b665b14eSToby Isaac   if (PetscBeganMPE) {
1952b665b14eSToby Isaac     char name[PETSC_MAX_PATH_LEN];
1953b665b14eSToby Isaac 
1954b665b14eSToby Isaac     PetscCall(PetscInfo(0, "Finalizing MPE.\n"));
1955b665b14eSToby Isaac     if (sname) {
1956b665b14eSToby Isaac       PetscCall(PetscStrncpy(name, sname, sizeof(name)));
19575c6c1daeSBarry Smith     } else {
1958b665b14eSToby Isaac       PetscCall(PetscGetProgramName(name, sizeof(name)));
19595c6c1daeSBarry Smith     }
1960b665b14eSToby Isaac     PetscCall(MPE_Finish_log(name));
19615c6c1daeSBarry Smith   } else {
1962b665b14eSToby Isaac     PetscCall(PetscInfo(0, "Not finalizing MPE (not started by PETSc).\n"));
19635c6c1daeSBarry Smith   }
1964c2a741eeSJunchao Zhang   #else
1965b665b14eSToby Isaac   SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, "PETSc was configured without MPE support, reconfigure with --with-mpe or --download-mpe");
1966c2a741eeSJunchao Zhang   #endif
19673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19685c6c1daeSBarry Smith }
19695c6c1daeSBarry Smith 
1970ffeef943SBarry Smith /*@
19717d6c928cSSatish Balay   PetscLogView - Prints a summary of the logging.
19725c6c1daeSBarry Smith 
19738f14a041SBarry Smith   Collective
19745c6c1daeSBarry Smith 
19755c6c1daeSBarry Smith   Input Parameter:
1976f14045dbSBarry Smith . viewer - an ASCII viewer
19775c6c1daeSBarry Smith 
19785c6c1daeSBarry Smith   Options Database Keys:
1979bb1d7374SBarry Smith + -log_view [:filename]                    - Prints summary of log information
1980bb1d7374SBarry Smith . -log_view :filename.py:ascii_info_detail - Saves logging information from each process as a Python file
1981607d249eSBarry Smith . -log_view :filename.xml:ascii_xml        - Saves a summary of the logging information in a nested format (see below for how to view it)
1982d0a29bd7SConnor Ward . -log_view :filename.txt:ascii_flamegraph - Saves logging information in a format suitable for visualising as a Flame Graph (see below for how to view it)
1983156b51fbSBarry Smith . -log_view_memory                         - Also display memory usage in each event
1984156b51fbSBarry Smith . -log_view_gpu_time                       - Also display time in each event for GPU kernels (Note this may slow the computation)
1985811af0c4SBarry Smith . -log_all                                 - Saves a file Log.rank for each MPI rank with details of each step of the computation
1986bb1d7374SBarry Smith - -log_trace [filename]                    - Displays a trace of what each process is doing
19875c6c1daeSBarry Smith 
1988d1f92df0SBarry Smith   Level: beginner
1989d1f92df0SBarry Smith 
19905c6c1daeSBarry Smith   Notes:
1991da81f932SPierre Jolivet   It is possible to control the logging programmatically but we recommend using the options database approach whenever possible
19925c6c1daeSBarry Smith   By default the summary is printed to stdout.
19935c6c1daeSBarry Smith 
1994bb1d7374SBarry Smith   Before calling this routine you must have called either PetscLogDefaultBegin() or PetscLogNestedBegin()
1995bb1d7374SBarry Smith 
1996bb1d7374SBarry Smith   If PETSc is configured with --with-logging=0 then this functionality is not available
1997bb1d7374SBarry Smith 
1998607d249eSBarry Smith   To view the nested XML format filename.xml first copy  ${PETSC_DIR}/share/petsc/xml/performance_xml2html.xsl to the current
1999607d249eSBarry Smith   directory then open filename.xml with your browser. Specific notes for certain browsers
20001d27aa22SBarry Smith .vb
20011d27aa22SBarry Smith     Firefox and Internet explorer - simply open the file
20021d27aa22SBarry Smith     Google Chrome - you must start up Chrome with the option --allow-file-access-from-files
20031d27aa22SBarry Smith     Safari - see https://ccm.net/faq/36342-safari-how-to-enable-local-file-access
20041d27aa22SBarry Smith .ve
20051d27aa22SBarry Smith   or one can use the package <http://xmlsoft.org/XSLT/xsltproc2.html> to translate the xml file to html and then open it with
2006607d249eSBarry Smith   your browser.
20072add09c0SLisandro Dalcin   Alternatively, use the script ${PETSC_DIR}/lib/petsc/bin/petsc-performance-view to automatically open a new browser
20082add09c0SLisandro Dalcin   window and render the XML log file contents.
2009607d249eSBarry Smith 
2010bb1d7374SBarry Smith   The nested XML format was kindly donated by Koos Huijssen and Christiaan M. Klaij  MARITIME  RESEARCH  INSTITUTE  NETHERLANDS
2011bb1d7374SBarry Smith 
20121d27aa22SBarry Smith   The Flame Graph output can be visualised using either the original Flame Graph script <https://github.com/brendangregg/FlameGraph>
20131d27aa22SBarry Smith   or using speedscope <https://www.speedscope.app>.
2014d0a29bd7SConnor Ward   Old XML profiles may be converted into this format using the script ${PETSC_DIR}/lib/petsc/bin/xml2flamegraph.py.
2015d0a29bd7SConnor Ward 
2016d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogDefaultBegin()`, `PetscLogDump()`
20175c6c1daeSBarry Smith @*/
2018d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogView(PetscViewer viewer)
2019d71ae5a4SJacob Faibussowitsch {
2020f14045dbSBarry Smith   PetscBool         isascii;
2021f14045dbSBarry Smith   PetscViewerFormat format;
2022b665b14eSToby Isaac   int               stage;
2023b665b14eSToby Isaac   PetscLogState     state;
2024b665b14eSToby Isaac   PetscIntStack     temp_stack;
2025b665b14eSToby Isaac   PetscLogHandler   handler;
2026b665b14eSToby Isaac   PetscBool         is_empty;
20275c6c1daeSBarry Smith 
20285c6c1daeSBarry Smith   PetscFunctionBegin;
2029b665b14eSToby Isaac   PetscCall(PetscLogGetState(&state));
203037b78d16SBarry Smith   /* Pop off any stages the user forgot to remove */
2031b665b14eSToby Isaac   PetscCall(PetscIntStackCreate(&temp_stack));
2032b665b14eSToby Isaac   PetscCall(PetscLogStateGetCurrentStage(state, &stage));
203337b78d16SBarry Smith   while (stage >= 0) {
2034b665b14eSToby Isaac     PetscCall(PetscLogStagePop());
2035b665b14eSToby Isaac     PetscCall(PetscIntStackPush(temp_stack, stage));
2036b665b14eSToby Isaac     PetscCall(PetscLogStateGetCurrentStage(state, &stage));
203737b78d16SBarry Smith   }
20389566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
203928b400f6SJacob Faibussowitsch   PetscCheck(isascii, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Currently can only view logging to ASCII");
20409566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
2041b665b14eSToby Isaac   if (format == PETSC_VIEWER_ASCII_XML || format == PETSC_VIEWER_ASCII_FLAMEGRAPH) {
2042294de794SToby Isaac     PetscCall(PetscLogGetHandler(PETSCLOGHANDLERNESTED, &handler));
2043b665b14eSToby Isaac     PetscCall(PetscLogHandlerView(handler, viewer));
2044b665b14eSToby Isaac   } else {
2045294de794SToby Isaac     PetscCall(PetscLogGetHandler(PETSCLOGHANDLERDEFAULT, &handler));
2046b665b14eSToby Isaac     PetscCall(PetscLogHandlerView(handler, viewer));
20475c6c1daeSBarry Smith   }
2048b665b14eSToby Isaac   PetscCall(PetscIntStackEmpty(temp_stack, &is_empty));
2049b665b14eSToby Isaac   while (!is_empty) {
2050b665b14eSToby Isaac     PetscCall(PetscIntStackPop(temp_stack, &stage));
2051b665b14eSToby Isaac     PetscCall(PetscLogStagePush(stage));
2052b665b14eSToby Isaac     PetscCall(PetscIntStackEmpty(temp_stack, &is_empty));
2053b665b14eSToby Isaac   }
2054b665b14eSToby Isaac   PetscCall(PetscIntStackDestroy(temp_stack));
20553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20565c6c1daeSBarry Smith }
20575c6c1daeSBarry Smith 
2058f14045dbSBarry Smith /*@C
2059811af0c4SBarry Smith   PetscLogViewFromOptions - Processes command line options to determine if/how a `PetscLog` is to be viewed.
2060f14045dbSBarry Smith 
2061811af0c4SBarry Smith   Collective on `PETSC_COMM_WORLD`
2062f14045dbSBarry Smith 
2063811af0c4SBarry Smith   Level: developer
2064f14045dbSBarry Smith 
2065d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogView()`
2066f14045dbSBarry Smith @*/
2067d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogViewFromOptions(void)
2068d71ae5a4SJacob Faibussowitsch {
2069ad2e3d55SToby Isaac   PetscInt          n_max = PETSC_LOG_VIEW_FROM_OPTIONS_MAX;
2070ad2e3d55SToby Isaac   PetscViewer       viewers[PETSC_LOG_VIEW_FROM_OPTIONS_MAX];
2071ad2e3d55SToby Isaac   PetscViewerFormat formats[PETSC_LOG_VIEW_FROM_OPTIONS_MAX];
2072f14045dbSBarry Smith   PetscBool         flg;
2073f14045dbSBarry Smith 
2074f14045dbSBarry Smith   PetscFunctionBegin;
2075648c30bcSBarry Smith   PetscCall(PetscOptionsCreateViewers(PETSC_COMM_WORLD, NULL, NULL, "-log_view", &n_max, viewers, formats, &flg));
2076*0bc1c2e8SToby Isaac   /*
2077*0bc1c2e8SToby Isaac      PetscLogHandlerView_Default_Info() wants to be sure that the only objects still around are these viewers, so keep track of how many there are
2078*0bc1c2e8SToby Isaac    */
2079*0bc1c2e8SToby Isaac   PetscLogNumViewersCreated = n_max;
2080ad2e3d55SToby Isaac   for (PetscInt i = 0; i < n_max; i++) {
2081*0bc1c2e8SToby Isaac     PetscInt refct;
2082*0bc1c2e8SToby Isaac 
2083ad2e3d55SToby Isaac     PetscCall(PetscViewerPushFormat(viewers[i], formats[i]));
2084ad2e3d55SToby Isaac     PetscCall(PetscLogView(viewers[i]));
2085ad2e3d55SToby Isaac     PetscCall(PetscViewerPopFormat(viewers[i]));
2086*0bc1c2e8SToby Isaac     PetscCall(PetscObjectGetReference((PetscObject)viewers[i], &refct));
2087648c30bcSBarry Smith     PetscCall(PetscViewerDestroy(&viewers[i]));
2088*0bc1c2e8SToby Isaac     if (refct == 1) PetscLogNumViewersDestroyed++;
2089f14045dbSBarry Smith   }
20903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2091f14045dbSBarry Smith }
2092f14045dbSBarry Smith 
2093b665b14eSToby Isaac PETSC_INTERN PetscErrorCode PetscLogHandlerNestedSetThreshold(PetscLogHandler, PetscLogDouble, PetscLogDouble *);
2094b665b14eSToby Isaac 
2095b665b14eSToby Isaac /*@
2096b665b14eSToby Isaac   PetscLogSetThreshold - Set the threshold time for logging the events; this is a percentage out of 100, so 1. means any event
2097b665b14eSToby Isaac   that takes 1 or more percent of the time.
2098b665b14eSToby Isaac 
20998f14a041SBarry Smith   Logically Collective on `PETSC_COMM_WORLD`
2100b665b14eSToby Isaac 
2101b665b14eSToby Isaac   Input Parameter:
2102b665b14eSToby Isaac . newThresh - the threshold to use
2103b665b14eSToby Isaac 
2104b665b14eSToby Isaac   Output Parameter:
2105b665b14eSToby Isaac . oldThresh - the previously set threshold value
2106b665b14eSToby Isaac 
2107b665b14eSToby Isaac   Options Database Keys:
2108b665b14eSToby Isaac . -log_view :filename.xml:ascii_xml - Prints an XML summary of flop and timing information to the file
2109b665b14eSToby Isaac 
2110b665b14eSToby Isaac   Example Usage:
2111b665b14eSToby Isaac .vb
2112b665b14eSToby Isaac   PetscInitialize(...);
2113b665b14eSToby Isaac   PetscLogNestedBegin();
2114b665b14eSToby Isaac   PetscLogSetThreshold(0.1,&oldthresh);
2115b665b14eSToby Isaac   // ... code ...
2116b665b14eSToby Isaac   PetscLogView(viewer);
2117b665b14eSToby Isaac   PetscFinalize();
2118b665b14eSToby Isaac .ve
2119b665b14eSToby Isaac 
2120b665b14eSToby Isaac   Level: advanced
2121b665b14eSToby Isaac 
2122b665b14eSToby Isaac   Note:
2123b665b14eSToby Isaac   This threshold is only used by the nested log handler
2124b665b14eSToby Isaac 
2125b665b14eSToby Isaac .seealso: `PetscLogDump()`, `PetscLogView()`, `PetscLogTraceBegin()`, `PetscLogDefaultBegin()`,
2126b665b14eSToby Isaac           `PetscLogNestedBegin()`
2127b665b14eSToby Isaac @*/
2128b665b14eSToby Isaac PetscErrorCode PetscLogSetThreshold(PetscLogDouble newThresh, PetscLogDouble *oldThresh)
2129b665b14eSToby Isaac {
2130b665b14eSToby Isaac   PetscLogHandler handler;
2131b665b14eSToby Isaac 
2132b665b14eSToby Isaac   PetscFunctionBegin;
2133294de794SToby Isaac   PetscCall(PetscLogTryGetHandler(PETSCLOGHANDLERNESTED, &handler));
2134b665b14eSToby Isaac   PetscCall(PetscLogHandlerNestedSetThreshold(handler, newThresh, oldThresh));
2135b665b14eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
2136b665b14eSToby Isaac }
2137b665b14eSToby Isaac 
21385c6c1daeSBarry Smith /*----------------------------------------------- Counter Functions -------------------------------------------------*/
2139cc4c1da9SBarry Smith /*@
21405c6c1daeSBarry Smith   PetscGetFlops - Returns the number of flops used on this processor
21415c6c1daeSBarry Smith   since the program began.
21425c6c1daeSBarry Smith 
21435c6c1daeSBarry Smith   Not Collective
21445c6c1daeSBarry Smith 
21455c6c1daeSBarry Smith   Output Parameter:
214610450e9eSJacob Faibussowitsch . flops - number of floating point operations
21475c6c1daeSBarry Smith 
2148d1f92df0SBarry Smith   Level: intermediate
2149d1f92df0SBarry Smith 
21505c6c1daeSBarry Smith   Notes:
21515c6c1daeSBarry Smith   A global counter logs all PETSc flop counts.  The user can use
2152811af0c4SBarry Smith   `PetscLogFlops()` to increment this counter to include flops for the
21535c6c1daeSBarry Smith   application code.
21545c6c1daeSBarry Smith 
21554b7c4d4dSPierre Jolivet   A separate counter `PetscLogGpuFlops()` logs the flops that occur on any GPU associated with this MPI rank
2156811af0c4SBarry Smith 
21574b7c4d4dSPierre Jolivet .seealso: [](ch_profiling), `PetscLogGpuFlops()`, `PetscTime()`, `PetscLogFlops()`
21585c6c1daeSBarry Smith @*/
2159d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscGetFlops(PetscLogDouble *flops)
2160d71ae5a4SJacob Faibussowitsch {
21615c6c1daeSBarry Smith   PetscFunctionBegin;
21625c6c1daeSBarry Smith   *flops = petsc_TotalFlops;
21633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21645c6c1daeSBarry Smith }
21655c6c1daeSBarry Smith 
21661c1ad86eSToby Isaac /*@C
21671c1ad86eSToby Isaac   PetscLogObjectState - Record information about an object with the default log handler
21681c1ad86eSToby Isaac 
21691c1ad86eSToby Isaac   Not Collective
21701c1ad86eSToby Isaac 
21711c1ad86eSToby Isaac   Input Parameters:
21721c1ad86eSToby Isaac + obj    - the `PetscObject`
21731c1ad86eSToby Isaac . format - a printf-style format string
21741c1ad86eSToby Isaac - ...    - printf arguments to format
21751c1ad86eSToby Isaac 
21761c1ad86eSToby Isaac   Level: developer
21771c1ad86eSToby Isaac 
2178b665b14eSToby Isaac .seealso: [](ch_profiling), `PetscLogObjectCreate()`, `PetscLogObjectDestroy()`, `PetscLogGetDefaultHandler()`
21791c1ad86eSToby Isaac @*/
2180d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogObjectState(PetscObject obj, const char format[], ...)
2181d71ae5a4SJacob Faibussowitsch {
21825c6c1daeSBarry Smith   PetscFunctionBegin;
2183dff009beSToby Isaac   for (PetscInt i = 0; i < PETSC_LOG_HANDLER_MAX; i++) {
2184dff009beSToby Isaac     PetscLogHandler h = PetscLogHandlers[i].handler;
2185dff009beSToby Isaac 
2186dff009beSToby Isaac     if (h) {
2187dff009beSToby Isaac       va_list Argp;
21885c6c1daeSBarry Smith       va_start(Argp, format);
2189dff009beSToby Isaac       PetscCall(PetscLogHandlerLogObjectState_Internal(h, obj, format, Argp));
21905c6c1daeSBarry Smith       va_end(Argp);
2191b665b14eSToby Isaac     }
2192dff009beSToby Isaac   }
21933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21945c6c1daeSBarry Smith }
21955c6c1daeSBarry Smith 
21965c6c1daeSBarry Smith /*MC
21975c6c1daeSBarry Smith   PetscLogFlops - Adds floating point operations to the global counter.
21985c6c1daeSBarry Smith 
21995c6c1daeSBarry Smith   Synopsis:
2200aaa7dc30SBarry Smith   #include <petsclog.h>
22015c6c1daeSBarry Smith   PetscErrorCode PetscLogFlops(PetscLogDouble f)
22025c6c1daeSBarry Smith 
22035c6c1daeSBarry Smith   Not Collective
22045c6c1daeSBarry Smith 
22055c6c1daeSBarry Smith   Input Parameter:
22065c6c1daeSBarry Smith . f - flop counter
22075c6c1daeSBarry Smith 
220810450e9eSJacob Faibussowitsch   Example Usage:
22095c6c1daeSBarry Smith .vb
22105c6c1daeSBarry Smith   PetscLogEvent USER_EVENT;
221110450e9eSJacob Faibussowitsch 
22125c6c1daeSBarry Smith   PetscLogEventRegister("User event", 0, &USER_EVENT);
22135c6c1daeSBarry Smith   PetscLogEventBegin(USER_EVENT, 0, 0, 0, 0);
22145c6c1daeSBarry Smith   [code segment to monitor]
22155c6c1daeSBarry Smith   PetscLogFlops(user_flops)
22165c6c1daeSBarry Smith   PetscLogEventEnd(USER_EVENT, 0, 0, 0, 0);
22175c6c1daeSBarry Smith .ve
22185c6c1daeSBarry Smith 
2219d1f92df0SBarry Smith   Level: intermediate
2220d1f92df0SBarry Smith 
2221811af0c4SBarry Smith   Note:
222210450e9eSJacob Faibussowitsch    A global counter logs all PETSc flop counts. The user can use PetscLogFlops() to increment
222310450e9eSJacob Faibussowitsch    this counter to include flops for the application code.
22245c6c1daeSBarry Smith 
22254b7c4d4dSPierre Jolivet .seealso: [](ch_profiling), `PetscLogGpuFlops()`, `PetscLogEventRegister()`, `PetscLogEventBegin()`, `PetscLogEventEnd()`, `PetscGetFlops()`
22265c6c1daeSBarry Smith M*/
22275c6c1daeSBarry Smith 
22285c6c1daeSBarry Smith /*MC
222910450e9eSJacob Faibussowitsch   PetscPreLoadBegin - Begin a segment of code that may be preloaded (run twice) to get accurate
223010450e9eSJacob Faibussowitsch   timings
22315c6c1daeSBarry Smith 
22325c6c1daeSBarry Smith   Synopsis:
2233aaa7dc30SBarry Smith   #include <petsclog.h>
22345c6c1daeSBarry Smith   void PetscPreLoadBegin(PetscBool flag, char *name);
22355c6c1daeSBarry Smith 
22365c6c1daeSBarry Smith   Not Collective
22375c6c1daeSBarry Smith 
2238d8d19677SJose E. Roman   Input Parameters:
223910450e9eSJacob Faibussowitsch + flag - `PETSC_TRUE` to run twice, `PETSC_FALSE` to run once, may be overridden with command
224010450e9eSJacob Faibussowitsch          line option `-preload true|false`
224110450e9eSJacob Faibussowitsch - name - name of first stage (lines of code timed separately with `-log_view`) to be preloaded
22425c6c1daeSBarry Smith 
224310450e9eSJacob Faibussowitsch   Example Usage:
22445c6c1daeSBarry Smith .vb
224510450e9eSJacob Faibussowitsch   PetscPreLoadBegin(PETSC_TRUE, "first stage");
224610450e9eSJacob Faibussowitsch   // lines of code
22475c6c1daeSBarry Smith   PetscPreLoadStage("second stage");
224810450e9eSJacob Faibussowitsch   // lines of code
22495c6c1daeSBarry Smith   PetscPreLoadEnd();
22505c6c1daeSBarry Smith .ve
22515c6c1daeSBarry Smith 
2252d1f92df0SBarry Smith   Level: intermediate
2253d1f92df0SBarry Smith 
2254811af0c4SBarry Smith   Note:
225595452b02SPatrick Sanan   Only works in C/C++, not Fortran
22565c6c1daeSBarry Smith 
225710450e9eSJacob Faibussowitsch   Flags available within the macro\:
225810450e9eSJacob Faibussowitsch + PetscPreLoadingUsed - `PETSC_TRUE` if we are or have done preloading
225910450e9eSJacob Faibussowitsch . PetscPreLoadingOn   - `PETSC_TRUE` if it is CURRENTLY doing preload
226010450e9eSJacob Faibussowitsch . PetscPreLoadIt      - `0` for the first computation (with preloading turned off it is only
226110450e9eSJacob Faibussowitsch                         `0`) `1`  for the second
226210450e9eSJacob Faibussowitsch - PetscPreLoadMax     - number of times it will do the computation, only one when preloading is
226310450e9eSJacob Faibussowitsch                         turned on
226410450e9eSJacob Faibussowitsch 
226510450e9eSJacob Faibussowitsch   The first two variables are available throughout the program, the second two only between the
226610450e9eSJacob Faibussowitsch   `PetscPreLoadBegin()` and `PetscPreLoadEnd()`
22675c6c1daeSBarry Smith 
2268d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogEventRegister()`, `PetscLogEventBegin()`, `PetscLogEventEnd()`, `PetscPreLoadEnd()`, `PetscPreLoadStage()`
22695c6c1daeSBarry Smith M*/
22705c6c1daeSBarry Smith 
22715c6c1daeSBarry Smith /*MC
227210450e9eSJacob Faibussowitsch   PetscPreLoadEnd - End a segment of code that may be preloaded (run twice) to get accurate
227310450e9eSJacob Faibussowitsch   timings
22745c6c1daeSBarry Smith 
22755c6c1daeSBarry Smith   Synopsis:
2276aaa7dc30SBarry Smith   #include <petsclog.h>
22775c6c1daeSBarry Smith   void PetscPreLoadEnd(void);
22785c6c1daeSBarry Smith 
22795c6c1daeSBarry Smith   Not Collective
22805c6c1daeSBarry Smith 
228110450e9eSJacob Faibussowitsch   Example Usage:
22825c6c1daeSBarry Smith .vb
228310450e9eSJacob Faibussowitsch   PetscPreLoadBegin(PETSC_TRUE, "first stage");
228410450e9eSJacob Faibussowitsch   // lines of code
22855c6c1daeSBarry Smith   PetscPreLoadStage("second stage");
228610450e9eSJacob Faibussowitsch   // lines of code
22875c6c1daeSBarry Smith   PetscPreLoadEnd();
22885c6c1daeSBarry Smith .ve
22895c6c1daeSBarry Smith 
2290d1f92df0SBarry Smith   Level: intermediate
2291d1f92df0SBarry Smith 
2292811af0c4SBarry Smith   Note:
2293dd01b7e5SBarry Smith   Only works in C/C++ not Fortran
22945c6c1daeSBarry Smith 
2295d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogEventRegister()`, `PetscLogEventBegin()`, `PetscLogEventEnd()`, `PetscPreLoadBegin()`, `PetscPreLoadStage()`
22965c6c1daeSBarry Smith M*/
22975c6c1daeSBarry Smith 
22985c6c1daeSBarry Smith /*MC
229910450e9eSJacob Faibussowitsch   PetscPreLoadStage - Start a new segment of code to be timed separately to get accurate timings
23005c6c1daeSBarry Smith 
23015c6c1daeSBarry Smith   Synopsis:
2302aaa7dc30SBarry Smith   #include <petsclog.h>
23035c6c1daeSBarry Smith   void PetscPreLoadStage(char *name);
23045c6c1daeSBarry Smith 
23055c6c1daeSBarry Smith   Not Collective
23065c6c1daeSBarry Smith 
230710450e9eSJacob Faibussowitsch   Example Usage:
23085c6c1daeSBarry Smith .vb
230910450e9eSJacob Faibussowitsch   PetscPreLoadBegin(PETSC_TRUE,"first stage");
231010450e9eSJacob Faibussowitsch   // lines of code
23115c6c1daeSBarry Smith   PetscPreLoadStage("second stage");
231210450e9eSJacob Faibussowitsch   // lines of code
23135c6c1daeSBarry Smith   PetscPreLoadEnd();
23145c6c1daeSBarry Smith .ve
23155c6c1daeSBarry Smith 
2316d1f92df0SBarry Smith   Level: intermediate
2317d1f92df0SBarry Smith 
2318811af0c4SBarry Smith   Note:
2319dd01b7e5SBarry Smith   Only works in C/C++ not Fortran
23205c6c1daeSBarry Smith 
2321d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogEventRegister()`, `PetscLogEventBegin()`, `PetscLogEventEnd()`, `PetscPreLoadBegin()`, `PetscPreLoadEnd()`
23225c6c1daeSBarry Smith M*/
23235c6c1daeSBarry Smith 
2324a4af0ceeSJacob Faibussowitsch   #if PetscDefined(HAVE_DEVICE)
2325a4af0ceeSJacob Faibussowitsch     #include <petsc/private/deviceimpl.h>
23269ffd0706SHong Zhang 
2327cc4c1da9SBarry Smith /*@
2328156b51fbSBarry Smith   PetscLogGpuTime - turn on the logging of GPU time for GPU kernels
2329156b51fbSBarry Smith 
2330811af0c4SBarry Smith   Options Database Key:
2331efa05fe8SBarry Smith . -log_view_gpu_time - provide the GPU times for all events in the `-log_view` output
2332156b51fbSBarry Smith 
2333d1f92df0SBarry Smith   Level: advanced
2334d1f92df0SBarry Smith 
2335156b51fbSBarry Smith   Notes:
233610450e9eSJacob Faibussowitsch   Turning on the timing of the GPU kernels can slow down the entire computation and should only
2337efa05fe8SBarry Smith   be used when studying the performance of individual operations on GPU such as vector operations and
233810450e9eSJacob Faibussowitsch   matrix-vector operations.
2339156b51fbSBarry Smith 
2340efa05fe8SBarry Smith   If this option is not used then times for most of the events in the `-log_view` output will be listed as Nan, indicating the times are not available
2341efa05fe8SBarry Smith 
234210450e9eSJacob Faibussowitsch   This routine should only be called once near the beginning of the program. Once it is started
234310450e9eSJacob Faibussowitsch   it cannot be turned off.
2344156b51fbSBarry Smith 
2345d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogView()`, `PetscLogGpuFlops()`, `PetscLogGpuTimeEnd()`, `PetscLogGpuTimeBegin()`
2346156b51fbSBarry Smith @*/
2347d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogGpuTime(void)
2348d71ae5a4SJacob Faibussowitsch {
2349473903fcSJunchao Zhang   PetscFunctionBegin;
2350473903fcSJunchao Zhang   PetscCheck(petsc_gtime == 0.0, PETSC_COMM_SELF, PETSC_ERR_SUP, "GPU logging has already been turned on");
2351156b51fbSBarry Smith   PetscLogGpuTimeFlag = PETSC_TRUE;
2352473903fcSJunchao Zhang   PetscFunctionReturn(PETSC_SUCCESS);
2353156b51fbSBarry Smith }
2354156b51fbSBarry Smith 
2355cc4c1da9SBarry Smith /*@
23569ffd0706SHong Zhang   PetscLogGpuTimeBegin - Start timer for device
23579ffd0706SHong Zhang 
2358d1f92df0SBarry Smith   Level: intermediate
2359d1f92df0SBarry Smith 
23609ffd0706SHong Zhang   Notes:
23617d766218SJunchao Zhang   When GPU is enabled, the timer is run on the GPU, it is a separate logging of time
236210450e9eSJacob Faibussowitsch   devoted to GPU computations (excluding kernel launch times).
2363811af0c4SBarry Smith 
23647d766218SJunchao Zhang   When GPU is not available, the timer is run on the CPU, it is a separate logging of
236510450e9eSJacob Faibussowitsch   time devoted to GPU computations (including kernel launch times).
2366811af0c4SBarry Smith 
236710450e9eSJacob Faibussowitsch   There is no need to call WaitForCUDA() or WaitForHIP() between `PetscLogGpuTimeBegin()` and
236810450e9eSJacob Faibussowitsch   `PetscLogGpuTimeEnd()`
2369811af0c4SBarry Smith 
237010450e9eSJacob Faibussowitsch   This timer should NOT include times for data transfers between the GPU and CPU, nor setup
237110450e9eSJacob Faibussowitsch   actions such as allocating space.
2372811af0c4SBarry Smith 
237310450e9eSJacob Faibussowitsch   The regular logging captures the time for data transfers and any CPU activities during the
237410450e9eSJacob Faibussowitsch   event. It is used to compute the flop rate on the GPU as it is actively engaged in running a
237510450e9eSJacob Faibussowitsch   kernel.
23769ffd0706SHong Zhang 
23779ffd0706SHong Zhang   Developer Notes:
237810450e9eSJacob Faibussowitsch   The GPU event timer captures the execution time of all the kernels launched in the default
237954c05997SPierre Jolivet   stream by the CPU between `PetscLogGpuTimeBegin()` and `PetscLogGpuTimeEnd()`.
2380811af0c4SBarry Smith 
238154c05997SPierre Jolivet   `PetscLogGpuTimeBegin()` and `PetscLogGpuTimeEnd()` insert the begin and end events into the
238210450e9eSJacob Faibussowitsch   default stream (stream 0). The device will record a time stamp for the event when it reaches
238310450e9eSJacob Faibussowitsch   that event in the stream. The function xxxEventSynchronize() is called in
238454c05997SPierre Jolivet   `PetscLogGpuTimeEnd()` to block CPU execution, but not continued GPU execution, until the
238510450e9eSJacob Faibussowitsch   timer event is recorded.
23869ffd0706SHong Zhang 
2387d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogView()`, `PetscLogGpuFlops()`, `PetscLogGpuTimeEnd()`, `PetscLogGpuTime()`
23889ffd0706SHong Zhang @*/
2389d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogGpuTimeBegin(void)
2390d71ae5a4SJacob Faibussowitsch {
2391b665b14eSToby Isaac   PetscBool isActive;
2392b665b14eSToby Isaac 
23939ffd0706SHong Zhang   PetscFunctionBegin;
2394b665b14eSToby Isaac   PetscCall(PetscLogEventBeginIsActive(&isActive));
2395b665b14eSToby Isaac   if (!isActive || !PetscLogGpuTimeFlag) PetscFunctionReturn(PETSC_SUCCESS);
23967d766218SJunchao Zhang     #if defined(PETSC_HAVE_DEVICE) && !defined(PETSC_HAVE_KOKKOS_WITHOUT_GPU)
23977d766218SJunchao Zhang   {
2398a4af0ceeSJacob Faibussowitsch     PetscDeviceContext dctx;
2399a4af0ceeSJacob Faibussowitsch 
24009566063dSJacob Faibussowitsch     PetscCall(PetscDeviceContextGetCurrentContext(&dctx));
24019566063dSJacob Faibussowitsch     PetscCall(PetscDeviceContextBeginTimer_Internal(dctx));
2402a4af0ceeSJacob Faibussowitsch   }
24037d766218SJunchao Zhang     #else
24047d766218SJunchao Zhang   PetscCall(PetscTimeSubtract(&petsc_gtime));
24057d766218SJunchao Zhang     #endif
24063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24079ffd0706SHong Zhang }
24089ffd0706SHong Zhang 
2409cc4c1da9SBarry Smith /*@
24109ffd0706SHong Zhang   PetscLogGpuTimeEnd - Stop timer for device
24119ffd0706SHong Zhang 
24129ffd0706SHong Zhang   Level: intermediate
24139ffd0706SHong Zhang 
2414d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogView()`, `PetscLogGpuFlops()`, `PetscLogGpuTimeBegin()`
24159ffd0706SHong Zhang @*/
2416d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLogGpuTimeEnd(void)
2417d71ae5a4SJacob Faibussowitsch {
2418b665b14eSToby Isaac   PetscBool isActive;
2419b665b14eSToby Isaac 
24209ffd0706SHong Zhang   PetscFunctionBegin;
2421b665b14eSToby Isaac   PetscCall(PetscLogEventEndIsActive(&isActive));
2422b665b14eSToby Isaac   if (!isActive || !PetscLogGpuTimeFlag) PetscFunctionReturn(PETSC_SUCCESS);
24237d766218SJunchao Zhang     #if defined(PETSC_HAVE_DEVICE) && !defined(PETSC_HAVE_KOKKOS_WITHOUT_GPU)
24247d766218SJunchao Zhang   {
2425a4af0ceeSJacob Faibussowitsch     PetscDeviceContext dctx;
2426a4af0ceeSJacob Faibussowitsch     PetscLogDouble     elapsed;
2427a4af0ceeSJacob Faibussowitsch 
24289566063dSJacob Faibussowitsch     PetscCall(PetscDeviceContextGetCurrentContext(&dctx));
24299566063dSJacob Faibussowitsch     PetscCall(PetscDeviceContextEndTimer_Internal(dctx, &elapsed));
2430a4af0ceeSJacob Faibussowitsch     petsc_gtime += (elapsed / 1000.0);
2431a4af0ceeSJacob Faibussowitsch   }
24327d766218SJunchao Zhang     #else
24337d766218SJunchao Zhang   PetscCall(PetscTimeAdd(&petsc_gtime));
24347d766218SJunchao Zhang     #endif
24353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24369ffd0706SHong Zhang }
2437c708d6e3SStefano Zampini 
24389ffd0706SHong Zhang   #endif /* end of PETSC_HAVE_DEVICE */
24399ffd0706SHong Zhang 
2440cb9ef012SToby Isaac #endif /* PETSC_USE_LOG*/
2441cb9ef012SToby Isaac 
2442dd01b7e5SBarry Smith /* -- Utility functions for logging from Fortran -- */
2443b665b14eSToby Isaac 
2444b665b14eSToby Isaac PETSC_EXTERN PetscErrorCode PetscASend(int count, int datatype)
2445b665b14eSToby Isaac {
2446b665b14eSToby Isaac   PetscFunctionBegin;
2447cb9ef012SToby Isaac #if PetscDefined(USE_LOG)
2448b665b14eSToby Isaac   PetscCall(PetscAddLogDouble(&petsc_send_ct, &petsc_send_ct_th, 1));
2449b665b14eSToby Isaac   #if !defined(MPIUNI_H) && !defined(PETSC_HAVE_BROKEN_RECURSIVE_MACRO)
2450b665b14eSToby Isaac   PetscCall(PetscMPITypeSize(count, MPI_Type_f2c((MPI_Fint)datatype), &petsc_send_len, &petsc_send_len_th));
2451b665b14eSToby Isaac   #endif
2452cb9ef012SToby Isaac #endif
2453b665b14eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
2454b665b14eSToby Isaac }
2455b665b14eSToby Isaac 
2456b665b14eSToby Isaac PETSC_EXTERN PetscErrorCode PetscARecv(int count, int datatype)
2457b665b14eSToby Isaac {
2458b665b14eSToby Isaac   PetscFunctionBegin;
2459cb9ef012SToby Isaac #if PetscDefined(USE_LOG)
2460b665b14eSToby Isaac   PetscCall(PetscAddLogDouble(&petsc_recv_ct, &petsc_recv_ct_th, 1));
2461b665b14eSToby Isaac   #if !defined(MPIUNI_H) && !defined(PETSC_HAVE_BROKEN_RECURSIVE_MACRO)
2462b665b14eSToby Isaac   PetscCall(PetscMPITypeSize(count, MPI_Type_f2c((MPI_Fint)datatype), &petsc_recv_len, &petsc_recv_len_th));
2463b665b14eSToby Isaac   #endif
2464cb9ef012SToby Isaac #endif
2465b665b14eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
2466b665b14eSToby Isaac }
2467b665b14eSToby Isaac 
2468b665b14eSToby Isaac PETSC_EXTERN PetscErrorCode PetscAReduce(void)
2469b665b14eSToby Isaac {
2470b665b14eSToby Isaac   PetscFunctionBegin;
2471cb9ef012SToby Isaac   if (PetscDefined(USE_LOG)) PetscCall(PetscAddLogDouble(&petsc_allreduce_ct, &petsc_allreduce_ct_th, 1));
2472b665b14eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
2473b665b14eSToby Isaac }
2474b665b14eSToby Isaac 
24755c6c1daeSBarry Smith PetscClassId PETSC_LARGEST_CLASSID = PETSC_SMALLEST_CLASSID;
24765c6c1daeSBarry Smith PetscClassId PETSC_OBJECT_CLASSID  = 0;
24775c6c1daeSBarry Smith 
24782611ad71SToby Isaac static PetscBool PetscLogInitializeCalled = PETSC_FALSE;
24792611ad71SToby Isaac 
24802611ad71SToby Isaac PETSC_INTERN PetscErrorCode PetscLogInitialize(void)
24812611ad71SToby Isaac {
24822611ad71SToby Isaac   int stage;
24832611ad71SToby Isaac 
24842611ad71SToby Isaac   PetscFunctionBegin;
24852611ad71SToby Isaac   if (PetscLogInitializeCalled) PetscFunctionReturn(PETSC_SUCCESS);
24862611ad71SToby Isaac   PetscLogInitializeCalled = PETSC_TRUE;
24872611ad71SToby Isaac   if (PetscDefined(USE_LOG)) {
24882611ad71SToby Isaac     /* Setup default logging structures */
24892611ad71SToby Isaac     PetscCall(PetscLogStateCreate(&petsc_log_state));
24902611ad71SToby Isaac     for (PetscInt i = 0; i < PETSC_LOG_HANDLER_MAX; i++) {
24912611ad71SToby Isaac       if (PetscLogHandlers[i].handler) PetscCall(PetscLogHandlerSetState(PetscLogHandlers[i].handler, petsc_log_state));
24922611ad71SToby Isaac     }
24932611ad71SToby Isaac     PetscCall(PetscLogStateStageRegister(petsc_log_state, "Main Stage", &stage));
24942611ad71SToby Isaac     PetscCall(PetscSpinlockCreate(&PetscLogSpinLock));
24952611ad71SToby Isaac #if defined(PETSC_HAVE_THREADSAFETY)
24962611ad71SToby Isaac     petsc_log_tid = 0;
24972611ad71SToby Isaac     petsc_log_gid = 0;
24982611ad71SToby Isaac #endif
24992611ad71SToby Isaac 
25002611ad71SToby Isaac     /* All processors sync here for more consistent logging */
25012611ad71SToby Isaac     PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
25022611ad71SToby Isaac     PetscCall(PetscTime(&petsc_BaseTime));
25032611ad71SToby Isaac     PetscCall(PetscLogStagePush(stage));
25042611ad71SToby Isaac   }
25052611ad71SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
25062611ad71SToby Isaac }
25072611ad71SToby Isaac 
25082611ad71SToby Isaac PETSC_INTERN PetscErrorCode PetscLogFinalize(void)
25092611ad71SToby Isaac {
25102611ad71SToby Isaac   PetscFunctionBegin;
25112611ad71SToby Isaac   if (PetscDefined(USE_LOG)) {
2512b665b14eSToby Isaac     /* Resetting phase */
2513b665b14eSToby Isaac     // pop remaining stages
2514b665b14eSToby Isaac     if (petsc_log_state) {
2515b665b14eSToby Isaac       while (petsc_log_state->current_stage >= 0) { PetscCall(PetscLogStagePop()); }
2516b665b14eSToby Isaac     }
25172611ad71SToby Isaac     for (int i = 0; i < PETSC_LOG_HANDLER_MAX; i++) PetscCall(PetscLogHandlerDestroy(&PetscLogHandlers[i].handler));
25182611ad71SToby Isaac     PetscCall(PetscArrayzero(PetscLogHandlers, PETSC_LOG_HANDLER_MAX));
25192611ad71SToby Isaac     PetscCall(PetscLogStateDestroy(&petsc_log_state));
25202611ad71SToby Isaac 
25212611ad71SToby Isaac     petsc_TotalFlops         = 0.0;
25222611ad71SToby Isaac     petsc_BaseTime           = 0.0;
25232611ad71SToby Isaac     petsc_TotalFlops         = 0.0;
25242611ad71SToby Isaac     petsc_send_ct            = 0.0;
25252611ad71SToby Isaac     petsc_recv_ct            = 0.0;
25262611ad71SToby Isaac     petsc_send_len           = 0.0;
25272611ad71SToby Isaac     petsc_recv_len           = 0.0;
25282611ad71SToby Isaac     petsc_isend_ct           = 0.0;
25292611ad71SToby Isaac     petsc_irecv_ct           = 0.0;
25302611ad71SToby Isaac     petsc_isend_len          = 0.0;
25312611ad71SToby Isaac     petsc_irecv_len          = 0.0;
25322611ad71SToby Isaac     petsc_wait_ct            = 0.0;
25332611ad71SToby Isaac     petsc_wait_any_ct        = 0.0;
25342611ad71SToby Isaac     petsc_wait_all_ct        = 0.0;
25352611ad71SToby Isaac     petsc_sum_of_waits_ct    = 0.0;
25362611ad71SToby Isaac     petsc_allreduce_ct       = 0.0;
25372611ad71SToby Isaac     petsc_gather_ct          = 0.0;
25382611ad71SToby Isaac     petsc_scatter_ct         = 0.0;
25392611ad71SToby Isaac     petsc_TotalFlops_th      = 0.0;
25402611ad71SToby Isaac     petsc_send_ct_th         = 0.0;
25412611ad71SToby Isaac     petsc_recv_ct_th         = 0.0;
25422611ad71SToby Isaac     petsc_send_len_th        = 0.0;
25432611ad71SToby Isaac     petsc_recv_len_th        = 0.0;
25442611ad71SToby Isaac     petsc_isend_ct_th        = 0.0;
25452611ad71SToby Isaac     petsc_irecv_ct_th        = 0.0;
25462611ad71SToby Isaac     petsc_isend_len_th       = 0.0;
25472611ad71SToby Isaac     petsc_irecv_len_th       = 0.0;
25482611ad71SToby Isaac     petsc_wait_ct_th         = 0.0;
25492611ad71SToby Isaac     petsc_wait_any_ct_th     = 0.0;
25502611ad71SToby Isaac     petsc_wait_all_ct_th     = 0.0;
25512611ad71SToby Isaac     petsc_sum_of_waits_ct_th = 0.0;
25522611ad71SToby Isaac     petsc_allreduce_ct_th    = 0.0;
25532611ad71SToby Isaac     petsc_gather_ct_th       = 0.0;
25542611ad71SToby Isaac     petsc_scatter_ct_th      = 0.0;
25552611ad71SToby Isaac 
25562611ad71SToby Isaac     petsc_ctog_ct    = 0.0;
25572611ad71SToby Isaac     petsc_gtoc_ct    = 0.0;
25582611ad71SToby Isaac     petsc_ctog_sz    = 0.0;
25592611ad71SToby Isaac     petsc_gtoc_sz    = 0.0;
25602611ad71SToby Isaac     petsc_gflops     = 0.0;
25612611ad71SToby Isaac     petsc_gtime      = 0.0;
25622611ad71SToby Isaac     petsc_ctog_ct_th = 0.0;
25632611ad71SToby Isaac     petsc_gtoc_ct_th = 0.0;
25642611ad71SToby Isaac     petsc_ctog_sz_th = 0.0;
25652611ad71SToby Isaac     petsc_gtoc_sz_th = 0.0;
25662611ad71SToby Isaac     petsc_gflops_th  = 0.0;
25672611ad71SToby Isaac     petsc_gtime_th   = 0.0;
25682611ad71SToby Isaac   }
25692611ad71SToby Isaac   PETSC_LARGEST_CLASSID    = PETSC_SMALLEST_CLASSID;
25702611ad71SToby Isaac   PETSC_OBJECT_CLASSID     = 0;
25712611ad71SToby Isaac   PetscLogInitializeCalled = PETSC_FALSE;
25722611ad71SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
25732611ad71SToby Isaac }
25742611ad71SToby Isaac 
2575cc4c1da9SBarry Smith /*@
25765c6c1daeSBarry Smith   PetscClassIdRegister - Registers a new class name for objects and logging operations in an application code.
25775c6c1daeSBarry Smith 
25785c6c1daeSBarry Smith   Not Collective
25795c6c1daeSBarry Smith 
25805c6c1daeSBarry Smith   Input Parameter:
25815c6c1daeSBarry Smith . name - The class name
25825c6c1daeSBarry Smith 
25835c6c1daeSBarry Smith   Output Parameter:
25845c6c1daeSBarry Smith . oclass - The class id or classid
25855c6c1daeSBarry Smith 
25865c6c1daeSBarry Smith   Level: developer
25875c6c1daeSBarry Smith 
2588d1f92df0SBarry Smith .seealso: [](ch_profiling), `PetscLogEventRegister()`
25895c6c1daeSBarry Smith @*/
2590d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscClassIdRegister(const char name[], PetscClassId *oclass)
2591d71ae5a4SJacob Faibussowitsch {
25925c6c1daeSBarry Smith   PetscFunctionBegin;
25935c6c1daeSBarry Smith   *oclass = ++PETSC_LARGEST_CLASSID;
25945c6c1daeSBarry Smith #if defined(PETSC_USE_LOG)
2595b665b14eSToby Isaac   {
2596b665b14eSToby Isaac     PetscLogState state;
2597b665b14eSToby Isaac     PetscLogClass logclass;
2598b665b14eSToby Isaac 
2599b665b14eSToby Isaac     PetscCall(PetscLogGetState(&state));
2600b665b14eSToby Isaac     if (state) PetscCall(PetscLogStateClassRegister(state, name, *oclass, &logclass));
2601b665b14eSToby Isaac   }
26025c6c1daeSBarry Smith #endif
26033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
26045c6c1daeSBarry Smith }
2605