Using OpenAD

We split the usage explanation in three parts, first an indroduction to a script that wraps the steps in the sequence of tools into one convenient call, second a more detailed explanation of the individual steps, and third a section on how to compile and link the resulting code.

The openad script

The different components of the OpenAD framework are used to transform the code in a sequence of steps, the "pipeline". Depending on the particular problem there are certain variations to the pipeline that achieve a better performance of the generated code. The most common pipeline setups are encapsulated in a perl script. The script is part of the OpenAD environment that is used to download and build the framework and relies on the same environment, so we cd into the OpenAD download/build environment with
    cd <path_to_OpenAD>
and set up for
The pipeline encapsulating perl script is in
    <path_to_OpenAD>/tools/openad/openad
Invoking it with the -h option displays the script usage of which the following choices are the most important:
Mode:
--bb-rev BasicBlockPreaccumulationReverse [Default]
--bb BasicBlockPreaccumulation
--memopstradeoff MemOpsTradeoffPreaccumulation

Both the --bb and the --memopstradeoff options select forward-only (tangent linear) modes. The -bb version uses the Angel library for determining an elimination sequence and the --memopstradeoff option uses a tradeoff between data locality and operations count heuristics to determine an elimination sequence.
As an example we use the toy problem code assumed to be in a file called head.f and transform this code with the following step:
openad head.f
which produces the adjoint version of head.f emitting the following messages:
Preprocess Fortran: 'head.f' -> 'head.pre.f'
Fortran to WHIRL: 'head.pre.f' -> 'head.pre.B'
WHIRL to XAIF: 'head.pre.B' -> 'head.xaif'
xaifBooster: 'head.xaif' -> 'head.xb.xaif'
XAIF to WHIRL: ('head.pre.B' 'head.xb.xaif') -> 'head.xb.x2w.B'
WHIRL to Fortran: 'head.xb.x2w.B' -> 'head.xb.x2w.w2f.f'
Postprocess Fortran: 'head.xb.x2w.w2f.f' -> 'head.xb.x2w.w2f.pp.f'

OpenAD: 'head.f' --> 'head.xb.x2w.w2f.pp.f'
For larger projects it is more appropriate to customize the sequence by adding the necessary steps to a Makefile. For an example see the Makefile used for the shallow water model case study. More details on the individual steps are given in the following section. Below that we give details on compiling and linking the resulting source file.

 Fortran tool pipeline

Similarly to the previous section we assume the target code to be contained in a file called head.f. In each step we extend the file name by a suitable extension or as  predetermined by the respective tools. This allows an easy construction of makefile targets.
  1. Canonicalizer: There is a Python based canonicalizer which performs a set of simple code transformations to the input code. The top level script is located at
    <path_to_OpenAD>/OpenADFortTk/tools/canonicalize/canon.v1.py
    and takes  the input file as an argument and dumps the output to stdout.
    input (fortran)
    head.f
    output (fortran)
    head.canon.f
    command line
    python canon.v1.py head.f > head.canon.f

  2. Fortran Front-End: mfef90 is the Fortran90 front-end of the Open64 compiler parses the Fortran input and produces a .B file containing Open64's internal representation (known as the "whirl" format) of the input code. The binary (for 32bit linux) is located at
    <path_to_OpenAD>/Open64/osprey1.0/targ_ia32_ia64_linux/crayf90/sgi/mfef90
    and takes a fortran file as an input argument. The -F flag enforces C preprocessing (if needed), the -N132 flag extends the fixed format line length in F77-style codes to 132 characters, the -z flag cleans up certain whirl inconsistencies and obfuscations to simplify the subsequent analysis and transformation
    input (fortran)
    head.canon.f
    output (whirl)
    head.canon.B
    command line
    mfef90 -z -F -N132 head.canon.f

  3. Translation from whirl to XAIF: using the binary (for 32 bit linux) located at
    <path_to_OpenAD>/OpenADFortTk/OpenADFortTk-x86-Linux/bin/whirl2xaif
    and takes the whirl file as an input argument and lets the user specify the output file with the -o switch.
    input (whirl)
    head.canon.B
    output (xaif)
    head.canon.xaif
    command line
    whirl2xaif -o head.canon.xaif head.canon.B

  4. AD-transformation: the xaifBooster component provides several transformation binaries <t> that modify XAIF content. The most important three <t> are located at
    <path_to_OpenAD>/xaifBooster/algorithms/BasicBlockPreaccumulation/test/t (tangent-linear code with basic block preaccumulation)
    <path_to_OpenAD>/xaifBooster/algorithms/MemOpsTradeoffPreaccumulation/test/t (tangent-linear code with data locality/operations tradeoffs)
    <path_to_OpenAD>/xaifBooster/algorithms/BasicBlockPreaccumulationReverse/test/t (adjoint code with basic block preaccumulation)
    and each will display a usage message giving a self explanatory list of options when invoked without any arguments. All of them require an input file, the XAIF schema files and the intrinsics catalogue. The previous step sets up the XAIF input to expect the XAIF schema files (symbolic links are sufficient)  in the working directory. The schema files are all the .xsd files found in
    <path_to_OpenAD>/xaif/schema/
    and the example intrinsics catalogue contained in
    <path_to_OpenAD>/xaif/schema/examples/inlinable_intrinsics.xaif
    works for whirl2xaif. Choose <t> to be one of the above transformations:
    input (xaif, schema and intrinsics catalogue)
    head.canon.xaif, *.xsd, inlinable_intrinsics.xaif
    output (xaif)
    head.canon.xb.xaif
    command line
    <t> -i head.canon.xaif -c <path_to_OpenAD>/xaif/schema/examples/inlinable_intrinsics.xaif -o head.canon.xb.xaif

  5. Translation from XAIF (modified)  to whirl: this step requires the original whirl content, i.e. the head.canon.B file. We use the binary (for 32 bit linux) located at
    <path_to_OpenAD>/OpenADFortTk/OpenADFortTk-x86-Linux/bin/xaif2whirl
    which takes the modified xaif and the original whirl file as input arguments and produces modified whirl as output.
    input (xaif and original whirl)
    head.canon.xb.xaif head.canon.B
    output (whirl)
    head.canon.xb.x2w.B
    command line
    xaif2whirl head.canon.B head.canon.xb.xaif
    The optional switch --structured indicates a well structured control flow graph, i.e. a control flow only consisting of nested loops and branches without alternative return, break, continue, goto. 

  6. Unparsing whirl to Fortran: Open64 provides an unparser as a binary which (for 32 bit linux) is located at
    <path_to_OpenAD>/Open64/osprey1.0/targ_ia32_ia64_linux/whirl2f/whirl2f
    and takes a whirl file as input. The whirl2f binary internally accesses whirl2f_be located in the same directory.
    input (whirl)
    head.canon.xb.x2w.B
    output (fortran)
    head.canon.xb.x2w.w2f.f
    command line
    whirl2f -openad head.canon.xb.x2w.B
    The optional switch -openad slightly modifies the behavior of the unparser in a way suitable for OpenAD.

  7. Fortran postprocessing step: this includes active type member access for all AD transformations and - applicable to reversal schemes in adjoint code generation - subroutine template expansion and inlinable subroutine call expansions. This is done with a perl script contained in
    <path_to_OpenAD>/OpenADFortTk/tools/multiprocess/multi-pp.pl
    and takes a fortran file as input. For tangent-linear code generation the following syntax applies:
    input (fortran)
    head.canon.xb.x2w.w2f.f
    output (fortran)
    head.canon.xb.x2w.w2f.pp.f
    command line
    perl multi-pp.pl -f head.canon.xb.x2w.w2f.f
    The -f switch (read "forward")  indicates postprocessing for tangent-linear code.
    For adjoint code generation the principal reversal scheme on a subroutine level (split, joint, or hybrid) is controlled by subroutine templates. Details on this approach can be found on the templating and inlining page as well as [Naumann2005STf].  The default name for the subroutine template file is ad_template.f and the default name for the inlinable subroutine calls file is ad_inline.f.  Examples can be found on the templating and inlining page.
    input (fortran)
    head.canon.xb.x2w.f,  ad_template.f, ad_inline.f
    output (fortran)
    head.canon.xb.x2w.pp.f
    command line
    perl multi-pp.pl head.canon.xb.x2w.w2f.f

Compiling and Linking

The fortran source resulting from the transformation pipeline relies on a new "active" type whose value component (%v) contains the original variable value and whose derivative component (%d) the corresponding derivative. The active type is defined in file active_module.f90 which can also be found in
<path_to_OpenAD>/xaifBooster/testRoundTrip/
and needs to be compiled prior to any of the transformed fortran sources or the "driver" code which makes use of the active type. The unparser of the Open64 front-end uses refers to specific kind parameters defined in  w2f__types.f90 which can be found in the same directory. Examples for drivers for  tangent-linear code can be found in the case studies.  Using a Fortran90 compiler <f95> of your choice  and assuming the driver code is contained in a file called driver.f, the following steps lead to an executable binary:
<f95> -c w2f__types.f90
<f95> -c active_module.f90
<f95> -c head.canon.xb.x2w.w2f.pp.f  # note that the transformed fortran source is fixed format Fortran90
<f95> w2f__types.o active_module.o head.canon.xb.x2w.w2f.pp.o -o driver driver.f
Examples for drivers for adjoint code can also be found in the case studies. In addition to the modules mentioned above we need another set of modules that work in conjunction with the subroutine templates and inlinable subroutine calls. In
<path_to_OpenAD>/xaifBooster/testRoundTrip/
we have a very simple set1 of these module files:
    OpenAD_checkpoints.f90
    OpenAD_dct.f90
    OpenAD_tape.f90
    OpenAD_rev.f90
which work with the code for inlinable subroutine calls  found  in
    ad_inline.f
and for the subroutine template expansion for joint reversal either
    simple.ad_template_joint.f
or  for split reversal
    simple.ad_template_split.f  
The OpenAD_*.f90 module files need to be compiled prior to compiling head.canon.xb.x2w.w2f.pp.f.To avoid a certain run-time condition the code in ad_inline.f uses a small C routine defined in
    iaddr.c
which needs to be compiled with a C compiler <cc> compatible <f95>. The complete sequence of steps  as follows:
<f95> -c w2f__types.f90
<f95> -c active_module.f90
<f95> -c OpenAD_checkpoints.f90
<f95> -c OpenAD_dct.f90
<f95> -c OpenAD_tape.f90
<f95> -c OpenAD_rev.f90
<cc>  -c iaddr.c
<f95> -c head.canon.xb.x2w.w2f.pp.f  # note that the transformed fortran source is fixed format Fortran90
<f95> w2f__types.o active_module.o OpenAD_checkpoints.o OpenAD_dct.o OpenAD_tape.o OpenAD_rev.o iaddr.c head.canon.xb.x2w.w2f.pp.o -o driver driver.f

1) a more efficient example of the support modules along with other template files can be found in the shallow water model case study.