xref: /petsc/src/sys/objects/fcallback.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1af0996ceSBarry Smith #include <petsc/private/petscimpl.h> /*I   "petscsys.h"    I*/
2f6291634SJed Brown 
3f6291634SJed Brown typedef struct _FortranCallbackLink *FortranCallbackLink;
4f6291634SJed Brown struct _FortranCallbackLink {
5f6291634SJed Brown   char                  *type_name;
6f6291634SJed Brown   PetscFortranCallbackId max;
7f6291634SJed Brown   FortranCallbackLink    next;
8f6291634SJed Brown };
9f6291634SJed Brown 
10f6291634SJed Brown typedef struct {
11e0cd13aeSBarry Smith   PetscFortranCallbackId basecount;
12e0cd13aeSBarry Smith   PetscFortranCallbackId maxsubtypecount;
13f6291634SJed Brown   FortranCallbackLink    subtypes;
14f6291634SJed Brown } FortranCallbackBase;
15f6291634SJed Brown 
16f6291634SJed Brown static FortranCallbackBase *_classbase;
17f6291634SJed Brown static PetscClassId         _maxclassid = PETSC_SMALLEST_CLASSID;
18f6291634SJed Brown 
19*9371c9d4SSatish Balay static PetscErrorCode PetscFortranCallbackFinalize(void) {
20f6291634SJed Brown   PetscFunctionBegin;
215f80ce2aSJacob Faibussowitsch   for (PetscInt i = PETSC_SMALLEST_CLASSID; i < _maxclassid; i++) {
22f6291634SJed Brown     FortranCallbackBase *base = &_classbase[i - PETSC_SMALLEST_CLASSID];
23f6291634SJed Brown     FortranCallbackLink  next, link = base->subtypes;
24f6291634SJed Brown     for (; link; link = next) {
25f6291634SJed Brown       next = link->next;
269566063dSJacob Faibussowitsch       PetscCall(PetscFree(link->type_name));
279566063dSJacob Faibussowitsch       PetscCall(PetscFree(link));
28f6291634SJed Brown     }
29f6291634SJed Brown   }
309566063dSJacob Faibussowitsch   PetscCall(PetscFree(_classbase));
31f6291634SJed Brown   _maxclassid = PETSC_SMALLEST_CLASSID;
32f6291634SJed Brown   PetscFunctionReturn(0);
33f6291634SJed Brown }
34f6291634SJed Brown 
35de6d466bSJed Brown /*@C
36f6291634SJed Brown    PetscFortranCallbackRegister - register a type+subtype callback
37f6291634SJed Brown 
38f6291634SJed Brown    Not Collective
39f6291634SJed Brown 
404165533cSJose E. Roman    Input Parameters:
41f6291634SJed Brown +  classid - ID of class on which to register callback
420298fd71SBarry Smith -  subtype - subtype string, or NULL for class ids
43f6291634SJed Brown 
444165533cSJose E. Roman    Output Parameter:
45f6291634SJed Brown .  id - callback id
46f6291634SJed Brown 
47f6291634SJed Brown    Level: developer
48f6291634SJed Brown 
49db781477SPatrick Sanan .seealso: `PetscFortranCallbackGetSizes()`
50f6291634SJed Brown @*/
51*9371c9d4SSatish Balay PetscErrorCode PetscFortranCallbackRegister(PetscClassId classid, const char *subtype, PetscFortranCallbackId *id) {
52f6291634SJed Brown   FortranCallbackBase *base;
53f6291634SJed Brown   FortranCallbackLink  link;
54f6291634SJed Brown 
55f6291634SJed Brown   PetscFunctionBegin;
563ca90d2dSJacob Faibussowitsch   if (subtype) PetscValidCharPointer(subtype, 2);
573ca90d2dSJacob Faibussowitsch   PetscValidPointer(id, 3);
58cc73adaaSBarry Smith   PetscCheck(classid >= PETSC_SMALLEST_CLASSID && classid <= PETSC_LARGEST_CLASSID, PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "ClassId %d corrupt", classid);
59f6291634SJed Brown   *id = 0;
60f6291634SJed Brown   if (classid >= _maxclassid) {
61f6291634SJed Brown     PetscClassId         newmax = PETSC_SMALLEST_CLASSID + 2 * (PETSC_LARGEST_CLASSID - PETSC_SMALLEST_CLASSID);
62f6291634SJed Brown     FortranCallbackBase *newbase;
639566063dSJacob Faibussowitsch     if (!_classbase) PetscCall(PetscRegisterFinalize(PetscFortranCallbackFinalize));
649566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(newmax - PETSC_SMALLEST_CLASSID, &newbase));
659566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(newbase, _classbase, _maxclassid - PETSC_SMALLEST_CLASSID));
669566063dSJacob Faibussowitsch     PetscCall(PetscFree(_classbase));
67a297a907SKarl Rupp 
68f6291634SJed Brown     _classbase  = newbase;
69f6291634SJed Brown     _maxclassid = newmax;
70f6291634SJed Brown   }
71f6291634SJed Brown   base = &_classbase[classid - PETSC_SMALLEST_CLASSID];
72a297a907SKarl Rupp   if (!subtype) *id = PETSC_SMALLEST_FORTRAN_CALLBACK + base->basecount++;
73a297a907SKarl Rupp   else {
74f6291634SJed Brown     for (link = base->subtypes; link; link = link->next) { /* look for either both NULL or matching values (implies both non-NULL) */
75f6291634SJed Brown       PetscBool match;
769566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(subtype, link->type_name, &match));
77f6291634SJed Brown       if (match) { /* base type or matching subtype */
78f6291634SJed Brown         goto found;
79f6291634SJed Brown       }
80f6291634SJed Brown     }
81f6291634SJed Brown     /* Not found. Create node and prepend to class' subtype list */
829566063dSJacob Faibussowitsch     PetscCall(PetscNew(&link));
839566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(subtype, &link->type_name));
84a297a907SKarl Rupp 
85f6291634SJed Brown     link->max      = PETSC_SMALLEST_FORTRAN_CALLBACK;
86f6291634SJed Brown     link->next     = base->subtypes;
87f6291634SJed Brown     base->subtypes = link;
88a297a907SKarl Rupp 
89f6291634SJed Brown   found:
90f6291634SJed Brown     *id = link->max++;
91a297a907SKarl Rupp 
92f6291634SJed Brown     base->maxsubtypecount = PetscMax(base->maxsubtypecount, link->max - PETSC_SMALLEST_FORTRAN_CALLBACK);
93f6291634SJed Brown   }
94f6291634SJed Brown   PetscFunctionReturn(0);
95f6291634SJed Brown }
96f6291634SJed Brown 
97de6d466bSJed Brown /*@C
98f6291634SJed Brown    PetscFortranCallbackGetSizes - get sizes of class and subtype pointer arrays
99f6291634SJed Brown 
100f6291634SJed Brown    Collective
101f6291634SJed Brown 
1024165533cSJose E. Roman    Input Parameter:
103f6291634SJed Brown .  classid - class Id
104f6291634SJed Brown 
1054165533cSJose E. Roman    Output Parameters:
106f6291634SJed Brown +  numbase - number of registered class callbacks
107f6291634SJed Brown -  numsubtype - max number of registered subtype callbacks
108f6291634SJed Brown 
109f6291634SJed Brown    Level: developer
110f6291634SJed Brown 
111db781477SPatrick Sanan .seealso: `PetscFortranCallbackRegister()`
112f6291634SJed Brown @*/
113*9371c9d4SSatish Balay PetscErrorCode PetscFortranCallbackGetSizes(PetscClassId classid, PetscFortranCallbackId *numbase, PetscFortranCallbackId *numsubtype) {
114f6291634SJed Brown   PetscFunctionBegin;
1155f80ce2aSJacob Faibussowitsch   PetscValidPointer(numbase, 2);
1165f80ce2aSJacob Faibussowitsch   PetscValidPointer(numsubtype, 3);
117f6291634SJed Brown   if (classid < _maxclassid) {
118f6291634SJed Brown     FortranCallbackBase *base = &_classbase[classid - PETSC_SMALLEST_CLASSID];
119f6291634SJed Brown     *numbase                  = base->basecount;
120f6291634SJed Brown     *numsubtype               = base->maxsubtypecount;
121f6291634SJed Brown   } else { /* nothing registered */
122f6291634SJed Brown     *numbase    = 0;
123f6291634SJed Brown     *numsubtype = 0;
124f6291634SJed Brown   }
125f6291634SJed Brown   PetscFunctionReturn(0);
126f6291634SJed Brown }
127