[petsc-users] TAO setup with modules in Fortran 90
Randy Michael Churchill
rchurchi at pppl.gov
Tue Oct 10 11:59:54 CDT 2017
Thanks, I have the rosenbrock1f.F90 example compiling and running now with
both 3.7 and 3.8. Transferring this setup knowledge to my codebase, I ran
into a problem, since I setup my module a little differently than the
rosenbrock1f example. In the specification part of the module, I put a Tao
object, e.g using the rosebrock1f as an example I do:
#include <petsc/finclude/petscdef.h>
module commondat
PetscReal :: alpha
PetscInt :: n
Tao :: taotest
end module commondat
The compiler throws errors on this however. :
commondat.F90(22): error #5082: Syntax error, found '::' when expecting one
of: => = . [ % ( :
Tao :: taotest
----------^
commondat.F90(22): error #6274: This statement must not appear in the
specification part of a module.
Tao :: taotest
------^
commondat.F90(22): error #6793: The POINTER attribute is required. [TAO]
Tao :: taotest
------^
What is the right way to use this? Do I have to define a type within the
module for the Tao object to be in? Or is this an issue in 3.7 when using
only the petscdef.h instead of petsc.h?
And thanks for the advice on moving to Petsc 3.8, we work with some Petsc
people, who know our code base, so will discuss with them on working
towards that.
On Fri, Oct 6, 2017 at 5:59 PM, Barry Smith <bsmith at mcs.anl.gov> wrote:
>
> > On Oct 6, 2017, at 11:18 PM, Randy Michael Churchill <rchurchi at pppl.gov>
> wrote:
> >
> > So if I'm limited to petsc 3.7.6 for reasons of eventually using within
> an existing, larger codebase that depends on 3.7.6, is it possible to use
> TAO with a user-defined module in Fortran90 using 3.7.6?
>
> You should really pus this "existing, larger codebase" to transition to
> PETSc 3.8 sooner, rather than later. Especially for developments in Fortran
> it will make life better for everyone. We are always willing to help users,
> once they have read the changes information, with information to make
> transitioning to the latest PETSc version easy. For any parts of the code
> in C, transitioning from 3.7 to 3.8 should be really simple.
>
> Barry
>
> >
> > I had tried the various forms of includes listed in the documentation,
> e.g. see below. I think I now realize this is an issue with the petsc
> installation on Edison, it does not seem to have the petsctao module in the
> library file (confirmed using nm -D on the library file). If I do the same
> include and use statement but with, for example, petscmat, it compiles fine.
> >
> > I built v3.8 from source, and the petsctao module is in the library
> file, and now the make works.
> >
> > commondat.F90
> > module commondat
> > #include <petsc/finclude/petscdef.h>
> > use petsc
> > PetscReal :: alpha
> > PetscInt :: n
> > end module commondat
> >
> > program rosenbrock1f
> > #include <petsc/finclude/petsctaodef.h>
> > use petsctao
> > use commondat
> >
> >
> >
> >
> > On Fri, Oct 6, 2017 at 7:54 AM, Matthew Knepley <knepley at gmail.com>
> wrote:
> > On Fri, Oct 6, 2017 at 7:36 AM, Barry Smith <bsmith at mcs.anl.gov> wrote:
> >
> > Randy,
> >
> > First you absolutely must use version 3.8 or the master development
> copy. We improved and simplified dramatically how Fortran (90) is utilized
> from PETSc.
> >
> > Note that there is only one simple set of include files and modules
> for Fortran; see the newest documentation.
> >
> >
> > http://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/
> Sys/UsingFortran.html
> >
> > Matt
> >
> >
> > Barry
> >
> >
> > > On Oct 5, 2017, at 11:48 PM, Randy Michael Churchill <
> rchurchi at pppl.gov> wrote:
> > >
> > > A simple setup question with TAO: if I were to convert the
> rosenbrock1f.F90 example to use a module instead of common structures, how
> would I setup the include statements? I've tried various combinations
> (using petscXXXdef.h, petscXXX.h, petscXXX.h90, along with use petscXXX),
> but seem to get errors with all.
> > >
> > > file:rosenbrock1f.h:
> > > module commondat
> > > PetscReal :: alpha
> > > PetscInt :: n
> > > end module commondat
> > >
> > > file:rosenbrock1f.90:
> > > program rosenbrock1f
> > > !!include statements??? which and where???!!!
> > > use commondat
> > > ...
> > >
> > > subroutine FormFunctionGradient(tao, X, f, G, dummy, ierr)
> > > use commondat
> > > implicit none
> > > ...
> > >
> > > (https://www.mcs.anl.gov/petsc/petsc-dev/src/tao/
> unconstrained/examples/tutorials/rosenbrock1f.F90.html)
> >
> >
> >
> >
> > --
> > 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/
> >
> >
> >
> > --
> > R. Michael Churchill
>
>
--
R. Michael Churchill
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20171010/d187248b/attachment-0001.html>
More information about the petsc-users
mailing list