<div dir="ltr"><div>OK thank you !</div><div><br></div><div>For the PetscDSAddBoundary I am not sure either ... I did this for now :</div><div><br></div><div>#ifdef PETSC_HAVE_FORTRAN_CAPS<br>#define petscdsaddboundary_ PETSCDSADDBOUNDARY<br>#elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) && !defined(FORTRANDOUBLEUNDERSCORE)<br>#define petscdsaddboundary_ petscdsaddboundary<br>#endif<br>PetscFortranCallbackId bocofunc, bocofunc_time;<br>static PetscErrorCode ourbocofunc(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *a_xI, const PetscScalar *a_xG, void *ctx)<br>{<br>    PetscObjectUseFortranCallback((PetscDS)ctx, bocofunc,<br>                                 (PetscReal*, const PetscReal*, const PetscReal*, const PetscScalar*, const PetscScalar*, void*, PetscErrorCode*),<br>                                 (&time, c, n, a_xI, a_xG, ctx, &ierr));<br>}<br>static PetscErrorCode ourbocofunc_time(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *a_xI, const PetscScalar *a_xG, void *ctx)<br>{<br>    PetscObjectUseFortranCallback((PetscDS)ctx, bocofunc_time,<br>                                 (PetscReal*, const PetscReal*, const PetscReal*, const PetscScalar*, const PetscScalar*, void*, PetscErrorCode*),<br>                                 (&time, c, n, a_xI, a_xG, ctx, &ierr));<br>}<br>PETSC_EXTERN void petscdsaddboundary_(PetscDS *prob, DMBoundaryConditionType *type, char *name, char *labelname, PetscInt *field, PetscInt *numcomps, PetscInt *comps,<br>                                      void (*bcFunc)(void),<br>                                      void (*bcFunc_t)(void),<br>                                      PetscInt *numids, const PetscInt *ids, void *ctx, PetscErrorCode *ierr, PETSC_FORTRAN_CHARLEN_T namelen, PETSC_FORTRAN_CHARLEN_T labelnamelen)<br>{<br>    char *newname, *newlabelname;<br>    FIXCHAR(name, namelen, newname);<br>    FIXCHAR(labelname, labelnamelen, newlabelname);<br>    *ierr = PetscObjectSetFortranCallback((PetscObject)*prob, PETSC_FORTRAN_CALLBACK_CLASS, &bocofunc, (PetscVoidFunction)bcFunc, prob);<br>    *ierr = PetscObjectSetFortranCallback((PetscObject)*prob, PETSC_FORTRAN_CALLBACK_CLASS, &bocofunc_time, (PetscVoidFunction)bcFunc_t, prob);<br>    *ierr = PetscDSAddBoundary(*prob, *type, newname, newlabelname, *field, *numcomps, comps, (void (*)(void))ourbocofunc, (void (*)(void))ourbocofunc_time, *numids, ids, ctx);<br>    FREECHAR(name, newname);<br>    FREECHAR(labelname, newlabelname);<br>}<br></div><div><br></div><div><br></div><div>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 ?</div><div><br></div><div>Thanks again,  </div><br clear="all"><div><div dir="ltr" class="gmail_signature" data-smartmail="gmail_signature"><div dir="ltr"><div><div dir="ltr"><div dir="ltr"><div><div><div><div><div>Thibault </div></div></div></div></div></div></div></div></div></div></div><br></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">Le dim. 13 déc. 2020 à 15:52, Matthew Knepley <<a href="mailto:knepley@gmail.com">knepley@gmail.com</a>> a écrit :<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div dir="ltr">On Sun, Dec 13, 2020 at 9:35 AM Thibault Bridel-Bertomeu <<a href="mailto:thibault.bridelbertomeu@gmail.com" target="_blank">thibault.bridelbertomeu@gmail.com</a>> wrote:<br></div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr">Hello Matthew, <div><br></div><div>Thank you for your answer and the pointer to SNESSetFunction(). Is the PETSC_F90_2PTR_PROTO(ptr) necessary for every wrapper of that category ?</div><div><br></div><div>So far, for the PetscDSSetRiemannSolver, I did this :</div><div><br></div><div>#ifdef PETSC_HAVE_FORTRAN_CAPS<br>#define petscdssetriemannsolver_ PETSCDSSETRIEMANNSOLVER<br>#elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) && !defined(FORTRANDOUBLEUNDERSCORE)<br>#define petscdssetriemannsolver_ petscdssetriemannsolver<br>#endif</div><div><br>PetscFortranCallbackId riemannsolver;</div><div><br>static PetscErrorCode ourriemannsolver(PetscInt dim, PetscInt Nf, PetscReal x[], PetscReal n[], PetscScalar uL[], PetscScalar uR[], PetscInt numConstants, PetscScalar constants[], PetscScalar flux[], void *ctx)<br>{<br>    PetscObjectUseFortranCallback((PetscDS)ctx, riemannsolver, (PetscInt*, PetscInt*, PetscReal*, PetscReal*, PetscScalar*, PetscScalar*, PetscInt*, PetscScalar*, PetscScalar*, void*, PetscErrorCode*),<br>                                  (&dim, &Nf, x, n, uL, uR, &numConstants, constants, flux, _ctx, &ierr));<br>}</div><div><br>PETSC_EXTERN void petscdssetriemannsolver_(PetscDS *prob, PetscInt *f,<br>                                           void (*rs)(PetscInt *dim, PetscInt *Nf, PetscReal x[], PetscReal n[], PetscScalar uL[], PetscScalar uR[], PetscInt *numConstants, PetscScalar constants[], PetscScalar flux[], void *ctx, PetscErrorCode *jerr),<br>                                           PetscErrorCode *ierr)<br>{<br>    *ierr = PetscObjectSetFortranCallback((PetscObject)*prob, PETSC_FORTRAN_CALLBACK_CLASS, &riemannsolver, (PetscVoidFunction)rs, prob);<br>    *ierr = PetscDSSetRiemannSolver(*prob, *f, (void*)ourriemannsolver);<br>}<br></div><div><br></div><div>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.</div><div>What do you think ? </div></div></blockquote><div><br></div><div>It looks good for right now. We can help you debug it when the rest of the wrappers are in place.</div><div><br></div><div>  Thanks,</div><div><br></div><div>     Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div><div><div dir="ltr"><div dir="ltr"><div><div dir="ltr"><div dir="ltr"><div><div><div><div><div>Thibault </div></div></div></div></div></div></div></div></div></div></div></div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">Le dim. 13 déc. 2020 à 15:29, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> a écrit :<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div dir="ltr">On Sun, Dec 13, 2020 at 9:17 AM Mark Adams <<a href="mailto:mfadams@lbl.gov" target="_blank">mfadams@lbl.gov</a>> wrote:<br></div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr">I don't think function pointers to PETSc (DS and DM) methods are not going to work in Fortran:<div><br></div><div>  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);<br></div><div><br></div><div>You could write a funcs.c file that you call from your fortran code, like, call setBC1(prob,...,ierr)</div><div><br></div><div>and put PhysicsBoundary_Advect_Inflow and setBC1 in funcs.c, for instance. </div></div></blockquote><div><br></div><div>Wrappers for functions that that function arguments are possible, but more involved. We have to create an internal struct that</div><div>holds Fortran function pointers, and then call them with the right arguments. This has been done for SNESSetFunction() for</div><div>instance, so we would start emulating that.</div><div><br></div><div>  Thanks,</div><div><br></div><div>      Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Sun, Dec 13, 2020 at 5:32 AM Thibault Bridel-Bertomeu <<a href="mailto:thibault.bridelbertomeu@gmail.com" target="_blank">thibault.bridelbertomeu@gmail.com</a>> wrote:<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr">Good morning all, <div><br></div><div>Thank you Barry for your answer.</div><div>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".</div><div>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 ?</div><div><br></div><div>Thank you very much, </div><div><br></div><div><div><div dir="ltr"><div dir="ltr"><div><div dir="ltr"><div dir="ltr"><div><div><div><div><div>Thibault</div></div></div></div></div></div></div></div></div></div></div><br></div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">Le sam. 12 déc. 2020 à 23:28, Barry Smith <<a href="mailto:bsmith@petsc.dev" target="_blank">bsmith@petsc.dev</a>> a écrit :<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><br><div><br><blockquote type="cite"><div>On Dec 12, 2020, at 2:59 PM, Thibault Bridel-Bertomeu <<a href="mailto:thibault.bridelbertomeu@gmail.com" target="_blank">thibault.bridelbertomeu@gmail.com</a>> wrote:</div><br><div><div dir="auto">Dear Jed, dear Barry,</div><div dir="auto"><br></div><div dir="auto">Thank you for the fast answers ! </div><div dir="auto"><br></div><div dir="auto">If I have any success I will make sure to make a pull request to provide some version of ex11 in Fortran. </div><div dir="auto"><br></div><div dir="auto">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. </div><div dir="auto">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 ?</div></div></blockquote><div><br></div>  Yes,</div><div><br><blockquote type="cite"><div><div dir="auto"> However there are a few things I could not find / understand yet. </div><div dir="auto">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 ...</div></div></blockquote><div><br></div>   This secret argument is put automatically by the Fortran compiler.</div><div><br><blockquote type="cite"><div><div dir="auto">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 ? </div></div></blockquote><div><br></div><div>PETSC_EXTERN void dmview_(DM *da,PetscViewer *vin,PetscErrorCode *ierr)</div><div>{</div><div>  PetscViewer v;</div><div>  PetscPatchDefaultViewers_Fortran(vin,v);</div><div>  *ierr = DMView(*da,v);</div><div>}</div><div><br></div><div>dm/interface/ftn-custom/zdmf.c</div><div><br></div><div><br></div><blockquote type="cite"><div><div dir="auto"><br></div><div dir="auto">Thank you very much again,</div><div dir="auto"><br></div><div dir="auto">Thibault Bridel-Bertomeu</div><div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">Le sam. 12 déc. 2020 à 21:48, Barry Smith <<a href="mailto:bsmith@petsc.dev" target="_blank">bsmith@petsc.dev</a>> a écrit :<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><div><br></div>   PETSc Fortran interfaces are a combination of automatically generated and manually generated. <div><br></div><div>   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. </div><div><br></div><div>   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.</div><div><br></div><div>   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.</div></div><div><div><br></div><div>  Barry</div><div><br><div><br><blockquote type="cite"><div>On Dec 12, 2020, at 9:30 AM, Thibault Bridel-Bertomeu <<a href="mailto:thibault.bridelbertomeu@gmail.com" target="_blank">thibault.bridelbertomeu@gmail.com</a>> wrote:</div><br><div><div dir="ltr"><div>Dear all, </div><div><br></div><div>Is there somewhere a version of the TS tutorial ex11.c in Fortran ?</div><div>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 ...</div><div><br></div><div>Thank you very much in advance,</div><br clear="all"><div><div dir="ltr"><div dir="ltr"><div><div dir="ltr"><div dir="ltr"><div><div><div><div><div>Thibault Bridel-Bertomeu<br>—<br></div></div></div></div>Eng, MSc, PhD</div><div>Research Engineer</div><div>CEA/CESTA</div><div>33114 LE BARP</div><div>Tel.: (+33)557046924</div><div>Mob.: (+33)611025322<br></div><div>Mail: <a href="mailto:thibault.bridelbertomeu@gmail.com" target="_blank">thibault.bridelbertomeu@gmail.com</a><br></div></div></div></div></div></div></div></div>
</div></blockquote></div><br></div></div></blockquote></div></div>-- <br><div dir="ltr"><div dir="ltr"><div><div dir="ltr"><div dir="ltr"><div><div><div><div><div>Thibault Bridel-Bertomeu<br>—<br></div></div></div></div>Eng, MSc, PhD</div><div>Research Engineer</div><div>CEA/CESTA</div><div>33114 LE BARP</div><div>Tel.: (+33)557046924</div><div>Mob.: (+33)611025322<br></div><div>Mail: <a href="mailto:thibault.bridelbertomeu@gmail.com" target="_blank">thibault.bridelbertomeu@gmail.com</a><br></div></div></div></div></div></div>
</div></blockquote></div><br></div></blockquote></div>
</blockquote></div>
</blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>-- Norbert Wiener</div><div><br></div><div><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>
</blockquote></div>
</blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>-- Norbert Wiener</div><div><br></div><div><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>
</blockquote></div>