[petsc-users] TS tutorial ex11 in Fortran

Thibault Bridel-Bertomeu thibault.bridelbertomeu at gmail.com
Sun Dec 13 08:53:51 CST 2020


OK thank you !

For the PetscDSAddBoundary I am not sure either ... I did this for now :

#ifdef PETSC_HAVE_FORTRAN_CAPS
#define petscdsaddboundary_ PETSCDSADDBOUNDARY
#elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) &&
!defined(FORTRANDOUBLEUNDERSCORE)
#define petscdsaddboundary_ petscdsaddboundary
#endif
PetscFortranCallbackId bocofunc, bocofunc_time;
static PetscErrorCode ourbocofunc(PetscReal time, const PetscReal *c, const
PetscReal *n, const PetscScalar *a_xI, const PetscScalar *a_xG, void *ctx)
{
    PetscObjectUseFortranCallback((PetscDS)ctx, bocofunc,
                                 (PetscReal*, const PetscReal*, const
PetscReal*, const PetscScalar*, const PetscScalar*, void*, PetscErrorCode*),
                                 (&time, c, n, a_xI, a_xG, ctx, &ierr));
}
static PetscErrorCode ourbocofunc_time(PetscReal time, const PetscReal *c,
const PetscReal *n, const PetscScalar *a_xI, const PetscScalar *a_xG, void
*ctx)
{
    PetscObjectUseFortranCallback((PetscDS)ctx, bocofunc_time,
                                 (PetscReal*, const PetscReal*, const
PetscReal*, const PetscScalar*, const PetscScalar*, void*, PetscErrorCode*),
                                 (&time, c, n, a_xI, a_xG, ctx, &ierr));
}
PETSC_EXTERN void petscdsaddboundary_(PetscDS *prob,
DMBoundaryConditionType *type, char *name, char *labelname, PetscInt
*field, PetscInt *numcomps, PetscInt *comps,
                                      void (*bcFunc)(void),
                                      void (*bcFunc_t)(void),
                                      PetscInt *numids, const PetscInt
*ids, void *ctx, PetscErrorCode *ierr, PETSC_FORTRAN_CHARLEN_T namelen,
PETSC_FORTRAN_CHARLEN_T labelnamelen)
{
    char *newname, *newlabelname;
    FIXCHAR(name, namelen, newname);
    FIXCHAR(labelname, labelnamelen, newlabelname);
    *ierr = PetscObjectSetFortranCallback((PetscObject)*prob,
PETSC_FORTRAN_CALLBACK_CLASS, &bocofunc, (PetscVoidFunction)bcFunc, prob);
    *ierr = PetscObjectSetFortranCallback((PetscObject)*prob,
PETSC_FORTRAN_CALLBACK_CLASS, &bocofunc_time, (PetscVoidFunction)bcFunc_t,
prob);
    *ierr = PetscDSAddBoundary(*prob, *type, newname, newlabelname, *field,
*numcomps, comps, (void (*)(void))ourbocofunc, (void
(*)(void))ourbocofunc_time, *numids, ids, ctx);
    FREECHAR(name, newname);
    FREECHAR(labelname, newlabelname);
}


But I do not know how to handle the two char* is the argument list : are
there going to be two PETSC_FORTRAN_CHARLEN_T like I wrote ?

Thanks again,

Thibault


Le dim. 13 déc. 2020 à 15:52, Matthew Knepley <knepley at gmail.com> a écrit :

> On Sun, Dec 13, 2020 at 9:35 AM Thibault Bridel-Bertomeu <
> thibault.bridelbertomeu at gmail.com> wrote:
>
>> Hello Matthew,
>>
>> Thank you for your answer and the pointer to SNESSetFunction(). Is the
>> PETSC_F90_2PTR_PROTO(ptr) necessary for every wrapper of that category ?
>>
>> So far, for the PetscDSSetRiemannSolver, I did this :
>>
>> #ifdef PETSC_HAVE_FORTRAN_CAPS
>> #define petscdssetriemannsolver_ PETSCDSSETRIEMANNSOLVER
>> #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) &&
>> !defined(FORTRANDOUBLEUNDERSCORE)
>> #define petscdssetriemannsolver_ petscdssetriemannsolver
>> #endif
>>
>> PetscFortranCallbackId riemannsolver;
>>
>> static PetscErrorCode ourriemannsolver(PetscInt dim, PetscInt Nf,
>> PetscReal x[], PetscReal n[], PetscScalar uL[], PetscScalar uR[], PetscInt
>> numConstants, PetscScalar constants[], PetscScalar flux[], void *ctx)
>> {
>>     PetscObjectUseFortranCallback((PetscDS)ctx, riemannsolver,
>> (PetscInt*, PetscInt*, PetscReal*, PetscReal*, PetscScalar*, PetscScalar*,
>> PetscInt*, PetscScalar*, PetscScalar*, void*, PetscErrorCode*),
>>                                   (&dim, &Nf, x, n, uL, uR,
>> &numConstants, constants, flux, _ctx, &ierr));
>> }
>>
>> PETSC_EXTERN void petscdssetriemannsolver_(PetscDS *prob, PetscInt *f,
>>                                            void (*rs)(PetscInt *dim,
>> PetscInt *Nf, PetscReal x[], PetscReal n[], PetscScalar uL[], PetscScalar
>> uR[], PetscInt *numConstants, PetscScalar constants[], PetscScalar flux[],
>> void *ctx, PetscErrorCode *jerr),
>>                                            PetscErrorCode *ierr)
>> {
>>     *ierr = PetscObjectSetFortranCallback((PetscObject)*prob,
>> PETSC_FORTRAN_CALLBACK_CLASS, &riemannsolver, (PetscVoidFunction)rs, prob);
>>     *ierr = PetscDSSetRiemannSolver(*prob, *f, (void*)ourriemannsolver);
>> }
>>
>> It compiles, but I have no idea whether it works or not, and I won't have
>> until the program is entirely built (with the TS and all) and it actually
>> has to call a Riemann solver added with that routine.
>> What do you think ?
>>
>
> It looks good for right now. We can help you debug it when the rest of the
> wrappers are in place.
>
>   Thanks,
>
>      Matt
>
>
>> Thibault
>>
>> Le dim. 13 déc. 2020 à 15:29, Matthew Knepley <knepley at gmail.com> a
>> écrit :
>>
>>> On Sun, Dec 13, 2020 at 9:17 AM Mark Adams <mfadams at lbl.gov> wrote:
>>>
>>>> I don't think function pointers to PETSc (DS and DM) methods are not
>>>> going to work in Fortran:
>>>>
>>>>   ierr = PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "inflow",
>>>>  "Face Sets", 0, 0, NULL, (void (*)(void)) PhysicsBoundary_Advect_Inflow,
>>>> NULL,  ALEN(inflowids),  inflowids,  phys);CHKERRQ(ierr);
>>>>
>>>> You could write a funcs.c file that you call from your fortran code,
>>>> like, call setBC1(prob,...,ierr)
>>>>
>>>> and put PhysicsBoundary_Advect_Inflow and setBC1 in funcs.c, for
>>>> instance.
>>>>
>>>
>>> Wrappers for functions that that function arguments are possible, but
>>> more involved. We have to create an internal struct that
>>> holds Fortran function pointers, and then call them with the right
>>> arguments. This has been done for SNESSetFunction() for
>>> instance, so we would start emulating that.
>>>
>>>   Thanks,
>>>
>>>       Matt
>>>
>>>
>>>>
>>>> On Sun, Dec 13, 2020 at 5:32 AM Thibault Bridel-Bertomeu <
>>>> thibault.bridelbertomeu at gmail.com> wrote:
>>>>
>>>>> Good morning all,
>>>>>
>>>>> Thank you Barry for your answer.
>>>>> I started adding some interfaces (PetscFVSetComponentName, PetscFVView
>>>>> & PetscFVSetType so far) and I have to say I think it is working quite
>>>>> well, but the prototypes of those functions are still quite "simple".
>>>>> I am stuck at how to implement the wrappers for
>>>>> PetscDSSetRiemannSolver and PetscDSSetContext though, especially on how to
>>>>> pass a function as an argument to PetscDSSetRiemannSolver ... Are there any
>>>>> similar functions that may already have their wrappers ?
>>>>>
>>>>> Thank you very much,
>>>>>
>>>>> Thibault
>>>>>
>>>>>
>>>>> Le sam. 12 déc. 2020 à 23:28, Barry Smith <bsmith at petsc.dev> a écrit :
>>>>>
>>>>>>
>>>>>>
>>>>>> On Dec 12, 2020, at 2:59 PM, Thibault Bridel-Bertomeu <
>>>>>> thibault.bridelbertomeu at gmail.com> wrote:
>>>>>>
>>>>>> Dear Jed, dear Barry,
>>>>>>
>>>>>> Thank you for the fast answers !
>>>>>>
>>>>>> If I have any success I will make sure to make a pull request to
>>>>>> provide some version of ex11 in Fortran.
>>>>>>
>>>>>> Regarding the stubs, I admit I started looking in that direction to
>>>>>> add the missing wrappers but I am not sure I fully understand the process
>>>>>> yet.
>>>>>> For each C function, I gotta provide a Fortran interface in a .h90
>>>>>> file as well as a C function that has a Fortran-like prototype and calls
>>>>>> the C function - right ?
>>>>>>
>>>>>>
>>>>>>   Yes,
>>>>>>
>>>>>> However there are a few things I could not find / understand yet.
>>>>>> For instance, it appears that for C functions that have character
>>>>>> string arguments take an extra argument in their
>>>>>> Fortran-like-prototype-wrapper, namely the length of the string. Is that
>>>>>> passed automatically ? I couldn’t find where it could come from ...
>>>>>>
>>>>>>
>>>>>>    This secret argument is put automatically by the Fortran compiler.
>>>>>>
>>>>>> Another thing is for functions like PetscFVView. I guess the wrapping
>>>>>> is less straightforward because I tried a quick something and it
>>>>>> segfault’ed. I couldnt find the wrapper for DMView although there is such a
>>>>>> routine in Fortran too. Could you please detail how to wrap such functions
>>>>>> ?
>>>>>>
>>>>>>
>>>>>> PETSC_EXTERN void dmview_(DM *da,PetscViewer *vin,PetscErrorCode
>>>>>> *ierr)
>>>>>> {
>>>>>>   PetscViewer v;
>>>>>>   PetscPatchDefaultViewers_Fortran(vin,v);
>>>>>>   *ierr = DMView(*da,v);
>>>>>> }
>>>>>>
>>>>>> dm/interface/ftn-custom/zdmf.c
>>>>>>
>>>>>>
>>>>>>
>>>>>> Thank you very much again,
>>>>>>
>>>>>> Thibault Bridel-Bertomeu
>>>>>>
>>>>>> Le sam. 12 déc. 2020 à 21:48, Barry Smith <bsmith at petsc.dev> a
>>>>>> écrit :
>>>>>>
>>>>>>>
>>>>>>>    PETSc Fortran interfaces are a combination of automatically
>>>>>>> generated and manually generated.
>>>>>>>
>>>>>>>    For any C PETSc function if the manual page begins with /*@  it
>>>>>>> generates the Fortran interface automatically (make allfortranstubs).  If
>>>>>>> it begins /*@C then either the Fortran interface is done manually or is
>>>>>>> missing.
>>>>>>>
>>>>>>>    C functions that have character string arguments or function
>>>>>>> arguments (or a few other special cases) need to be manually provided.  The
>>>>>>> automatically generated stubs go in the directory ftn-auto while manually
>>>>>>> generated ones go in the directory fin-custom.
>>>>>>>
>>>>>>>    Perhaps you could first generate a list of "missing" Fortran
>>>>>>> stubs and then for each stub determine why it is missing and if it can be
>>>>>>> provided. Some are likely easy to provide but a few (involving function
>>>>>>> arguments) will be more involved. Once you have all the stubs available
>>>>>>> translating ex11.c becomes straightforward.
>>>>>>>
>>>>>>>   Barry
>>>>>>>
>>>>>>>
>>>>>>> On Dec 12, 2020, at 9:30 AM, Thibault Bridel-Bertomeu <
>>>>>>> thibault.bridelbertomeu at gmail.com> wrote:
>>>>>>>
>>>>>>> Dear all,
>>>>>>>
>>>>>>> Is there somewhere a version of the TS tutorial ex11.c in Fortran ?
>>>>>>> I am looking into building in F90 (let's say that it is an
>>>>>>> unavoidable constraint) an unstructured 3D solver of the Euler equations
>>>>>>> using the "new" features of PETSc - mostly DMPlex & PetscFV - but I think
>>>>>>> there are some interfaces missing and I find it hard to find workarounds in
>>>>>>> Fortran. I would be grateful if anyone could please give me some pointers
>>>>>>> ...
>>>>>>>
>>>>>>> Thank you very much in advance,
>>>>>>>
>>>>>>> Thibault Bridel-Bertomeu
>>>>>>>>>>>>>> Eng, MSc, PhD
>>>>>>> Research Engineer
>>>>>>> CEA/CESTA
>>>>>>> 33114 LE BARP
>>>>>>> Tel.: (+33)557046924
>>>>>>> Mob.: (+33)611025322
>>>>>>> Mail: thibault.bridelbertomeu at gmail.com
>>>>>>>
>>>>>>>
>>>>>>> --
>>>>>> Thibault Bridel-Bertomeu
>>>>>>>>>>>> Eng, MSc, PhD
>>>>>> Research Engineer
>>>>>> CEA/CESTA
>>>>>> 33114 LE BARP
>>>>>> Tel.: (+33)557046924
>>>>>> Mob.: (+33)611025322
>>>>>> Mail: thibault.bridelbertomeu at gmail.com
>>>>>>
>>>>>>
>>>>>>
>>>
>>> --
>>> What most experimenters take for granted before they begin their
>>> experiments is infinitely more interesting than any results to which their
>>> experiments lead.
>>> -- Norbert Wiener
>>>
>>> https://www.cse.buffalo.edu/~knepley/
>>> <http://www.cse.buffalo.edu/~knepley/>
>>>
>>
>
> --
> What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which their
> experiments lead.
> -- Norbert Wiener
>
> https://www.cse.buffalo.edu/~knepley/
> <http://www.cse.buffalo.edu/~knepley/>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20201213/4c73b906/attachment.html>


More information about the petsc-users mailing list