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