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 199371c9d4SSatish 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 36*811af0c4SBarry Smith PetscFortranCallbackRegister - register a type+subtype callback. This is used by the PETSc Fortran interface to allow the use of user Fortran functions as arguments 37*811af0c4SBarry Smith to PETSc functions that take function pointers 38f6291634SJed Brown 39f6291634SJed Brown Not Collective 40f6291634SJed Brown 414165533cSJose E. Roman Input Parameters: 42f6291634SJed Brown + classid - ID of class on which to register callback 430298fd71SBarry Smith - subtype - subtype string, or NULL for class ids 44f6291634SJed Brown 454165533cSJose E. Roman Output Parameter: 46f6291634SJed Brown . id - callback id 47f6291634SJed Brown 48f6291634SJed Brown Level: developer 49f6291634SJed Brown 50db781477SPatrick Sanan .seealso: `PetscFortranCallbackGetSizes()` 51f6291634SJed Brown @*/ 529371c9d4SSatish Balay PetscErrorCode PetscFortranCallbackRegister(PetscClassId classid, const char *subtype, PetscFortranCallbackId *id) { 53f6291634SJed Brown FortranCallbackBase *base; 54f6291634SJed Brown FortranCallbackLink link; 55f6291634SJed Brown 56f6291634SJed Brown PetscFunctionBegin; 573ca90d2dSJacob Faibussowitsch if (subtype) PetscValidCharPointer(subtype, 2); 583ca90d2dSJacob Faibussowitsch PetscValidPointer(id, 3); 59cc73adaaSBarry Smith PetscCheck(classid >= PETSC_SMALLEST_CLASSID && classid <= PETSC_LARGEST_CLASSID, PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "ClassId %d corrupt", classid); 60f6291634SJed Brown *id = 0; 61f6291634SJed Brown if (classid >= _maxclassid) { 62f6291634SJed Brown PetscClassId newmax = PETSC_SMALLEST_CLASSID + 2 * (PETSC_LARGEST_CLASSID - PETSC_SMALLEST_CLASSID); 63f6291634SJed Brown FortranCallbackBase *newbase; 649566063dSJacob Faibussowitsch if (!_classbase) PetscCall(PetscRegisterFinalize(PetscFortranCallbackFinalize)); 659566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(newmax - PETSC_SMALLEST_CLASSID, &newbase)); 669566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(newbase, _classbase, _maxclassid - PETSC_SMALLEST_CLASSID)); 679566063dSJacob Faibussowitsch PetscCall(PetscFree(_classbase)); 68a297a907SKarl Rupp 69f6291634SJed Brown _classbase = newbase; 70f6291634SJed Brown _maxclassid = newmax; 71f6291634SJed Brown } 72f6291634SJed Brown base = &_classbase[classid - PETSC_SMALLEST_CLASSID]; 73a297a907SKarl Rupp if (!subtype) *id = PETSC_SMALLEST_FORTRAN_CALLBACK + base->basecount++; 74a297a907SKarl Rupp else { 75f6291634SJed Brown for (link = base->subtypes; link; link = link->next) { /* look for either both NULL or matching values (implies both non-NULL) */ 76f6291634SJed Brown PetscBool match; 779566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(subtype, link->type_name, &match)); 78f6291634SJed Brown if (match) { /* base type or matching subtype */ 79f6291634SJed Brown goto found; 80f6291634SJed Brown } 81f6291634SJed Brown } 82f6291634SJed Brown /* Not found. Create node and prepend to class' subtype list */ 839566063dSJacob Faibussowitsch PetscCall(PetscNew(&link)); 849566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(subtype, &link->type_name)); 85a297a907SKarl Rupp 86f6291634SJed Brown link->max = PETSC_SMALLEST_FORTRAN_CALLBACK; 87f6291634SJed Brown link->next = base->subtypes; 88f6291634SJed Brown base->subtypes = link; 89a297a907SKarl Rupp 90f6291634SJed Brown found: 91f6291634SJed Brown *id = link->max++; 92a297a907SKarl Rupp 93f6291634SJed Brown base->maxsubtypecount = PetscMax(base->maxsubtypecount, link->max - PETSC_SMALLEST_FORTRAN_CALLBACK); 94f6291634SJed Brown } 95f6291634SJed Brown PetscFunctionReturn(0); 96f6291634SJed Brown } 97f6291634SJed Brown 98de6d466bSJed Brown /*@C 99f6291634SJed Brown PetscFortranCallbackGetSizes - get sizes of class and subtype pointer arrays 100f6291634SJed Brown 101f6291634SJed Brown Collective 102f6291634SJed Brown 1034165533cSJose E. Roman Input Parameter: 104f6291634SJed Brown . classid - class Id 105f6291634SJed Brown 1064165533cSJose E. Roman Output Parameters: 107f6291634SJed Brown + numbase - number of registered class callbacks 108f6291634SJed Brown - numsubtype - max number of registered subtype callbacks 109f6291634SJed Brown 110f6291634SJed Brown Level: developer 111f6291634SJed Brown 112db781477SPatrick Sanan .seealso: `PetscFortranCallbackRegister()` 113f6291634SJed Brown @*/ 1149371c9d4SSatish Balay PetscErrorCode PetscFortranCallbackGetSizes(PetscClassId classid, PetscFortranCallbackId *numbase, PetscFortranCallbackId *numsubtype) { 115f6291634SJed Brown PetscFunctionBegin; 1165f80ce2aSJacob Faibussowitsch PetscValidPointer(numbase, 2); 1175f80ce2aSJacob Faibussowitsch PetscValidPointer(numsubtype, 3); 118f6291634SJed Brown if (classid < _maxclassid) { 119f6291634SJed Brown FortranCallbackBase *base = &_classbase[classid - PETSC_SMALLEST_CLASSID]; 120f6291634SJed Brown *numbase = base->basecount; 121f6291634SJed Brown *numsubtype = base->maxsubtypecount; 122f6291634SJed Brown } else { /* nothing registered */ 123f6291634SJed Brown *numbase = 0; 124f6291634SJed Brown *numsubtype = 0; 125f6291634SJed Brown } 126f6291634SJed Brown PetscFunctionReturn(0); 127f6291634SJed Brown } 128