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