xref: /petsc/src/sys/logging/handler/impls/nested/lognested.h (revision f39ec7874ab38227f85e1a66dd12007450d0487f)
1 #pragma once
2 
3 #include <petsc/private/loghandlerimpl.h>
4 #include <petsc/private/logimpl.h>
5 #include <../src/sys/logging/handler/impls/default/logdefault.h>
6 #include <petsc/private/hashmap.h>
7 
8 typedef int NestedId;
9 
10 static inline PETSC_UNUSED NestedId NestedIdFromEvent(PetscLogEvent e)
11 {
12   return e;
13 }
14 static inline PETSC_UNUSED NestedId NestedIdFromStage(PetscLogStage s)
15 {
16   return -(s + 2);
17 }
18 
19 typedef struct _n_NestedIdPair NestedIdPair;
20 struct _n_NestedIdPair {
21   PetscLogEvent root;
22   NestedId      leaf;
23 };
24 
25 #define NestedIdPairHash(key)     PetscHashCombine(PetscHash_UInt32((PetscHash32_t)((key).root)), PetscHash_UInt32((PetscHash32_t)((key).leaf)))
26 #define NestedIdPairEqual(k1, k2) (((k1).root == (k2).root) && ((k1).leaf == (k2).leaf))
27 
28 PETSC_HASH_MAP(NestedHash, NestedIdPair, PetscLogEvent, NestedIdPairHash, NestedIdPairEqual, -1)
29 
30 typedef struct _n_PetscLogHandler_Nested *PetscLogHandler_Nested;
31 struct _n_PetscLogHandler_Nested {
32   PetscLogState   state;
33   PetscLogHandler handler;
34   PetscNestedHash pair_map;
35   PetscIntStack   stack; // stack of nested ids
36   PetscClassId    nested_stage_id;
37   PetscLogDouble  threshold;
38 };
39 
40 typedef struct {
41   const char *name;
42   PetscInt    id;
43   PetscInt    parent;
44   PetscInt    num_descendants;
45 } PetscNestedEventNode;
46 
47 typedef struct {
48   MPI_Comm              comm;
49   PetscLogGlobalNames   global_events;
50   PetscNestedEventNode *nodes;
51   PetscEventPerfInfo   *perf;
52 } PetscNestedEventTree;
53 
54 typedef enum {
55   PETSC_LOG_NESTED_XML,
56   PETSC_LOG_NESTED_FLAMEGRAPH
57 } PetscLogNestedType;
58 
59 PETSC_INTERN PetscErrorCode PetscLogHandlerView_Nested_XML(PetscLogHandler_Nested, PetscNestedEventTree *, PetscViewer);
60 PETSC_INTERN PetscErrorCode PetscLogHandlerView_Nested_Flamegraph(PetscLogHandler_Nested, PetscNestedEventTree *, PetscViewer);
61