xref: /petsc/src/sys/objects/fcallback.c (revision af0996ce37bc06907c37d8d91773840993d61e62)
1*af0996ceSBarry 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 {
11f6291634SJed Brown   PetscInt            basecount;
12f6291634SJed Brown   PetscInt            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 
19f6291634SJed Brown #undef __FUNCT__
20f6291634SJed Brown #define __FUNCT__ "PetscFortranCallbackFinalize"
21f6291634SJed Brown static PetscErrorCode PetscFortranCallbackFinalize(void)
22f6291634SJed Brown {
23f6291634SJed Brown   PetscErrorCode ierr;
24f6291634SJed Brown   PetscClassId   i;
25f6291634SJed Brown 
26f6291634SJed Brown   PetscFunctionBegin;
27f6291634SJed Brown   for (i=PETSC_SMALLEST_CLASSID; i<_maxclassid; i++) {
28f6291634SJed Brown     FortranCallbackBase *base = &_classbase[i-PETSC_SMALLEST_CLASSID];
29f6291634SJed Brown     FortranCallbackLink next,link = base->subtypes;
30f6291634SJed Brown     for (; link; link=next) {
31f6291634SJed Brown       next = link->next;
32f6291634SJed Brown       ierr = PetscFree(link->type_name);CHKERRQ(ierr);
33f6291634SJed Brown       ierr = PetscFree(link);CHKERRQ(ierr);
34f6291634SJed Brown     }
35f6291634SJed Brown   }
36f6291634SJed Brown   ierr = PetscFree(_classbase);CHKERRQ(ierr);
37a297a907SKarl Rupp 
38f6291634SJed Brown   _maxclassid = PETSC_SMALLEST_CLASSID;
39f6291634SJed Brown   PetscFunctionReturn(0);
40f6291634SJed Brown }
41f6291634SJed Brown 
42f6291634SJed Brown #undef __FUNCT__
43f6291634SJed Brown #define __FUNCT__ "PetscFortranCallbackRegister"
44de6d466bSJed Brown /*@C
45f6291634SJed Brown    PetscFortranCallbackRegister - register a type+subtype callback
46f6291634SJed Brown 
47f6291634SJed Brown    Not Collective
48f6291634SJed Brown 
49f6291634SJed Brown    Input Arguments:
50f6291634SJed Brown +  classid - ID of class on which to register callback
510298fd71SBarry Smith -  subtype - subtype string, or NULL for class ids
52f6291634SJed Brown 
53f6291634SJed Brown    Output Arguments:
54f6291634SJed Brown .  id - callback id
55f6291634SJed Brown 
56f6291634SJed Brown    Level: developer
57f6291634SJed Brown 
58f6291634SJed Brown .seealso: PetscFortranCallbackGetSizes()
59f6291634SJed Brown @*/
60f6291634SJed Brown PetscErrorCode PetscFortranCallbackRegister(PetscClassId classid,const char *subtype,PetscFortranCallbackId *id)
61f6291634SJed Brown {
62f6291634SJed Brown   PetscErrorCode      ierr;
63f6291634SJed Brown   FortranCallbackBase *base;
64f6291634SJed Brown   FortranCallbackLink link;
65f6291634SJed Brown 
66f6291634SJed Brown   PetscFunctionBegin;
67f6291634SJed Brown   *id = 0;
68f6291634SJed Brown   if (classid < PETSC_SMALLEST_CLASSID || PETSC_LARGEST_CLASSID <= classid) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"ClassId %D corrupt",classid);
69f6291634SJed Brown   if (classid >= _maxclassid) {
70f6291634SJed Brown     PetscClassId        newmax = PETSC_SMALLEST_CLASSID + 2*(PETSC_LARGEST_CLASSID-PETSC_SMALLEST_CLASSID);
71f6291634SJed Brown     FortranCallbackBase *newbase;
72f6291634SJed Brown     if (!_classbase) {
73f6291634SJed Brown       ierr = PetscRegisterFinalize(PetscFortranCallbackFinalize);CHKERRQ(ierr);
74f6291634SJed Brown     }
75854ce69bSBarry Smith     ierr = PetscCalloc1(newmax-PETSC_SMALLEST_CLASSID,&newbase);CHKERRQ(ierr);
76f6291634SJed Brown     ierr = PetscMemcpy(newbase,_classbase,(_maxclassid-PETSC_SMALLEST_CLASSID)*sizeof(_classbase[0]));CHKERRQ(ierr);
77f6291634SJed Brown     ierr = PetscFree(_classbase);CHKERRQ(ierr);
78a297a907SKarl Rupp 
79f6291634SJed Brown     _classbase = newbase;
80f6291634SJed Brown     _maxclassid = newmax;
81f6291634SJed Brown   }
82f6291634SJed Brown   base = &_classbase[classid-PETSC_SMALLEST_CLASSID];
83a297a907SKarl Rupp   if (!subtype) *id = PETSC_SMALLEST_FORTRAN_CALLBACK + base->basecount++;
84a297a907SKarl Rupp   else {
85f6291634SJed Brown     for (link=base->subtypes; link; link=link->next) { /* look for either both NULL or matching values (implies both non-NULL) */
86f6291634SJed Brown       PetscBool match;
87f6291634SJed Brown       ierr = PetscStrcmp(subtype,link->type_name,&match);CHKERRQ(ierr);
88f6291634SJed Brown       if (match) { /* base type or matching subtype */
89f6291634SJed Brown         goto found;
90f6291634SJed Brown       }
91f6291634SJed Brown     }
92f6291634SJed Brown     /* Not found. Create node and prepend to class' subtype list */
93854ce69bSBarry Smith     ierr = PetscNew(&link);CHKERRQ(ierr);
94f6291634SJed Brown     ierr = PetscStrallocpy(subtype,&link->type_name);CHKERRQ(ierr);
95a297a907SKarl Rupp 
96f6291634SJed Brown     link->max      = PETSC_SMALLEST_FORTRAN_CALLBACK;
97f6291634SJed Brown     link->next     = base->subtypes;
98f6291634SJed Brown     base->subtypes = link;
99a297a907SKarl Rupp 
100f6291634SJed Brown found:
101f6291634SJed Brown     *id = link->max++;
102a297a907SKarl Rupp 
103f6291634SJed Brown     base->maxsubtypecount = PetscMax(base->maxsubtypecount,link->max-PETSC_SMALLEST_FORTRAN_CALLBACK);
104f6291634SJed Brown   }
105f6291634SJed Brown   PetscFunctionReturn(0);
106f6291634SJed Brown }
107f6291634SJed Brown 
108f6291634SJed Brown #undef __FUNCT__
109f6291634SJed Brown #define __FUNCT__ "PetscFortranCallbackGetSizes"
110de6d466bSJed Brown /*@C
111f6291634SJed Brown    PetscFortranCallbackGetSizes - get sizes of class and subtype pointer arrays
112f6291634SJed Brown 
113f6291634SJed Brown    Collective
114f6291634SJed Brown 
115f6291634SJed Brown    Input Arguments:
116f6291634SJed Brown .  classid - class Id
117f6291634SJed Brown 
118f6291634SJed Brown    Output Arguments:
119f6291634SJed Brown +  numbase - number of registered class callbacks
120f6291634SJed Brown -  numsubtype - max number of registered subtype callbacks
121f6291634SJed Brown 
122f6291634SJed Brown    Level: developer
123f6291634SJed Brown 
124f6291634SJed Brown .seealso: PetscFortranCallbackRegister()
125f6291634SJed Brown @*/
126f6291634SJed Brown PetscErrorCode PetscFortranCallbackGetSizes(PetscClassId classid,PetscInt *numbase,PetscInt *numsubtype)
127f6291634SJed Brown {
128f6291634SJed Brown 
129f6291634SJed Brown   PetscFunctionBegin;
130f6291634SJed Brown   if (classid < _maxclassid) {
131f6291634SJed Brown     FortranCallbackBase *base = &_classbase[classid-PETSC_SMALLEST_CLASSID];
132f6291634SJed Brown     *numbase    = base->basecount;
133f6291634SJed Brown     *numsubtype = base->maxsubtypecount;
134f6291634SJed Brown   } else {                      /* nothing registered */
135f6291634SJed Brown     *numbase    = 0;
136f6291634SJed Brown     *numsubtype = 0;
137f6291634SJed Brown   }
138f6291634SJed Brown   PetscFunctionReturn(0);
139f6291634SJed Brown }
140