1f6291634SJed Brown #include <petscsys.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); 37f6291634SJed Brown _maxclassid = PETSC_SMALLEST_CLASSID; 38f6291634SJed Brown PetscFunctionReturn(0); 39f6291634SJed Brown } 40f6291634SJed Brown 41f6291634SJed Brown #undef __FUNCT__ 42f6291634SJed Brown #define __FUNCT__ "PetscFortranCallbackRegister" 43*de6d466bSJed Brown /*@C 44f6291634SJed Brown PetscFortranCallbackRegister - register a type+subtype callback 45f6291634SJed Brown 46f6291634SJed Brown Not Collective 47f6291634SJed Brown 48f6291634SJed Brown Input Arguments: 49f6291634SJed Brown + classid - ID of class on which to register callback 50f6291634SJed Brown - subtype - subtype string, or PETSC_NULL for class ids 51f6291634SJed Brown 52f6291634SJed Brown Output Arguments: 53f6291634SJed Brown . id - callback id 54f6291634SJed Brown 55f6291634SJed Brown Level: developer 56f6291634SJed Brown 57f6291634SJed Brown .seealso: PetscFortranCallbackGetSizes() 58f6291634SJed Brown @*/ 59f6291634SJed Brown PetscErrorCode PetscFortranCallbackRegister(PetscClassId classid,const char *subtype,PetscFortranCallbackId *id) 60f6291634SJed Brown { 61f6291634SJed Brown PetscErrorCode ierr; 62f6291634SJed Brown FortranCallbackBase *base; 63f6291634SJed Brown FortranCallbackLink link; 64f6291634SJed Brown 65f6291634SJed Brown PetscFunctionBegin; 66f6291634SJed Brown *id = 0; 67f6291634SJed Brown if (classid < PETSC_SMALLEST_CLASSID || PETSC_LARGEST_CLASSID <= classid) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"ClassId %D corrupt",classid); 68f6291634SJed Brown if (classid >= _maxclassid) { 69f6291634SJed Brown PetscClassId newmax = PETSC_SMALLEST_CLASSID + 2*(PETSC_LARGEST_CLASSID-PETSC_SMALLEST_CLASSID); 70f6291634SJed Brown FortranCallbackBase *newbase; 71f6291634SJed Brown if (!_classbase) { 72f6291634SJed Brown ierr = PetscRegisterFinalize(PetscFortranCallbackFinalize);CHKERRQ(ierr); 73f6291634SJed Brown } 74f6291634SJed Brown ierr = PetscMalloc((newmax-PETSC_SMALLEST_CLASSID)*sizeof(_classbase[0]),&newbase);CHKERRQ(ierr); 75f6291634SJed Brown ierr = PetscMemzero(newbase,(newmax-PETSC_SMALLEST_CLASSID)*sizeof(_classbase[0]));CHKERRQ(ierr); 76f6291634SJed Brown ierr = PetscMemcpy(newbase,_classbase,(_maxclassid-PETSC_SMALLEST_CLASSID)*sizeof(_classbase[0]));CHKERRQ(ierr); 77f6291634SJed Brown ierr = PetscFree(_classbase);CHKERRQ(ierr); 78f6291634SJed Brown _classbase = newbase; 79f6291634SJed Brown _maxclassid = newmax; 80f6291634SJed Brown } 81f6291634SJed Brown base = &_classbase[classid-PETSC_SMALLEST_CLASSID]; 82f6291634SJed Brown if (!subtype) { 83f6291634SJed Brown *id = PETSC_SMALLEST_FORTRAN_CALLBACK + base->basecount++; 84f6291634SJed Brown } 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 */ 93f6291634SJed Brown ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr); 94f6291634SJed Brown ierr = PetscStrallocpy(subtype,&link->type_name);CHKERRQ(ierr); 95f6291634SJed Brown link->max = PETSC_SMALLEST_FORTRAN_CALLBACK; 96f6291634SJed Brown link->next = base->subtypes; 97f6291634SJed Brown base->subtypes = link; 98f6291634SJed Brown found: 99f6291634SJed Brown *id = link->max++; 100f6291634SJed Brown base->maxsubtypecount = PetscMax(base->maxsubtypecount,link->max-PETSC_SMALLEST_FORTRAN_CALLBACK); 101f6291634SJed Brown } 102f6291634SJed Brown PetscFunctionReturn(0); 103f6291634SJed Brown } 104f6291634SJed Brown 105f6291634SJed Brown #undef __FUNCT__ 106f6291634SJed Brown #define __FUNCT__ "PetscFortranCallbackGetSizes" 107*de6d466bSJed Brown /*@C 108f6291634SJed Brown PetscFortranCallbackGetSizes - get sizes of class and subtype pointer arrays 109f6291634SJed Brown 110f6291634SJed Brown Collective 111f6291634SJed Brown 112f6291634SJed Brown Input Arguments: 113f6291634SJed Brown . classid - class Id 114f6291634SJed Brown 115f6291634SJed Brown Output Arguments: 116f6291634SJed Brown + numbase - number of registered class callbacks 117f6291634SJed Brown - numsubtype - max number of registered subtype callbacks 118f6291634SJed Brown 119f6291634SJed Brown Level: developer 120f6291634SJed Brown 121f6291634SJed Brown .seealso: PetscFortranCallbackRegister() 122f6291634SJed Brown @*/ 123f6291634SJed Brown PetscErrorCode PetscFortranCallbackGetSizes(PetscClassId classid,PetscInt *numbase,PetscInt *numsubtype) 124f6291634SJed Brown { 125f6291634SJed Brown 126f6291634SJed Brown PetscFunctionBegin; 127f6291634SJed Brown if (classid < _maxclassid) { 128f6291634SJed Brown FortranCallbackBase *base = &_classbase[classid-PETSC_SMALLEST_CLASSID]; 129f6291634SJed Brown *numbase = base->basecount; 130f6291634SJed Brown *numsubtype = base->maxsubtypecount; 131f6291634SJed Brown } else { /* nothing registered */ 132f6291634SJed Brown *numbase = 0; 133f6291634SJed Brown *numsubtype = 0; 134f6291634SJed Brown } 135f6291634SJed Brown PetscFunctionReturn(0); 136f6291634SJed Brown } 137