Skip to content

Fortran Style Guide

Naming Convention

Ultimately this is a matter of preference. Here is a style guide that we like and that seems to be prevalent in most scientific codes (as well as the Fortran standard library) and you are welcome to follow it.

  • Use lowercase for all Fortran constructs (do, subroutine, module, ...).
  • Follow short mathematical notation for mathematical variables/functions (Ylm, Gamma, gamma, Enl, Rnl, ...).
  • For other names use all lowercase: try to keep names to one or two syllables; if more are required, use underscores to clarify ( meshexp, numstrings, linspace, meshgrid, argsort, spline, spline_interp, stoperr, stop_error).

For example spline interpolation can be shortened to spline_interpolation, spline_interpolate, spline_interp, spline, but not to splineint (“int” could mean integration, integer, etc. — too much ambiguity, even in the clear context of a computational code). This is in contrast to get_argument() where getarg() is perfectly clean and clear.

Indentation

Use 3 spaces indentation.

Floating Point Numbers

Use dp parameter define in params if needed:

    integer, parameter:: dp=kind(0.d0)                   ! double precision
and declare floats as:
    real(dp) :: a, b, c
Always write all floating point constants with the _dp suffix:

1.0_dp, 3.5_dp, 1.34e8_dp

and never any other way (see also the gotcha Floating Point Numbers). Omitting the dot in the literal constant is also incorrect.

To print floating point double precision numbers without losing precision, use the (es23.16) format (see http://stackoverflow.com/questions/6118231/why-do-i-need-17-significant-digits-and-not-16-to-represent-a-double/).

Modules and Programs

Only use modules and programs. Always setup a module in the following way:

    module ode1d
        use types, only: dp, pi
        use utils, only: stop_error
    implicit none
        private
        public integrate, normalize, parsefunction, get_val, rk4step, eulerstep, &
        rk4step2, get_midpoints, rk4_integrate, rk4_integrate_inward, &
        protected local_array

    contains

    subroutine get_val(...)
     ...
    end subroutine get_val
    ...
    end module ode1d
The implicit none statement works for the whole module (so you don’t need to worry about it). By keeping the private empty, all your subroutines/data types will be private to the module by default. Then you export things by putting it into the public clause.

Setup programs in the following way:

    program uranium
        use fmesh, only: mesh_exp
        use utils, only: stop_error, dp
        use dft, only: atom
        implicit none

    integer, parameter  :: Z = 92
    real(dp), parameter :: r_min = 8e-9_dp, r_max = 50.0_dp, a = 1e7_dp
    ...
    print *, "I am running"
    end program

Notice the “explicit imports” in the use statements. You can also use “implicit imports” like:

use fmesh

This should be avoided (“explicit is better than implicit”) in most cases.

Arrays

When passing arrays in and out of a subroutine/function, use the following pattern for 1D arrays (it is called assumed-shape):

subroutine f(r)
    real(dp), intent(out) :: r(:)
    integer :: n, i
    n = size(r)
    do i = 1, n
        r(i) = 1.0_dp / i**2
    enddo
end subroutine

2D arrays:

    subroutine g(A)
        real(dp), intent(in) :: A(:, :)
        ...
    end subroutine

and call it like this:

    real(dp) :: r(5)
    call f(r)

No array copying is done above. It has the following advantages:

* the shape and size of the array is passed in automatically
* the shape is checked at compile time, the size optionally at runtime
* allows to use strides and all kinds of array arithmetic without actually copying any data.

In order for the array to start with different index than 1, do:

    subroutine print_eigenvalues(kappa_min, lam)
        integer, intent(in) :: kappa_min
        real(dp), intent(in) :: lam(kappa_min:)

        integer :: kappa
        do kappa = kappa_min, ubound(lam, 1)
            print *, kappa, lam(kappa)
        end do
    end subroutine

Multidimensional Arrays

Always access slices as V(:, 1), V(:, 2), or V(:, :, 1), e.g. the colons should be on the left. That way the stride is contiguous and it will be fast. So when you need some slice in your algorithm, always setup the array in a way, so that you call it as above. If you put the colon on the right, it will be slow.

Example:

dydx = matmul(C(:, :, i), y) ! fast
dydx = matmul(C(i, :, :), y) ! slow

In other words, the “fortran storage order” is: smallest/fastest changing/innermost-loop index first, largest/slowest/outermost-loop index last (“Inner-most are left-most.”). So the elements of a 3D array A(N1,N2,N3) are stored, and thus most efficiently accessed, as:

do k = zs,ze
    do j = ys, ye
        do i = xs, xe
            A(i, j, k)
        end do
    end do
end do

Associated array of vectors would then be most efficiently accessed as:

do k = zs,ze
    do j = ys, ye
        A(:, j, k)
    end do
end do

Source Formatting

In general, follow existing style (e.g., from radiation module). To accomodate ease of parsing just described, we will insist on a uniform format.

The structure of a Bifrost routine includes: * the subroutine or function statement, followed by statements declaring the type and dimensions of the arguments; * an implicit none declaration; * a summary of the purpose of the routine; * descriptions of each of the Arguments in the order of the argument list; * (optionally) Further Details e.g., reference to jurnal paper etc.

Example:

 subroutine rhunt(n,xx,x,jlo)

  ! Purpose
  ! -------
  ! Given an array xx(1:N), and given a value x, returns a value jlo
  ! such that x is between xx(jlo) and xx(jlo+1). xx must be
  ! monotonicly increasing, jlo=0 or jlo=N is returned to indicate
  ! that x is out of range.
  ! the value of jlo should be used as initial guess, otherwise HUNT
  ! is not faster than LOCATE.
  ! Arguments
  ! ---------
  ! n, intent(in), scalar, integer      - size of the vector
  ! x, intent(in), scalar, real(4)      - input value 
  ! xx, intent(in), 1D array, real(4)   - input array of size n
  ! jlo, intent(inout), scalar, integer - initial guess on input, return jlo
  !                                     - such that x is between xx(jlo) and xx(jlo+1)
  ! Notes
  ! -----
  ! May 28 2006 : Copied from Numerical Recipes.

  implicit none

  integer                   , intent(in)    :: n
  real(kind=4)              , intent(in)    :: x
  real(kind=4), dimension(:), intent(in)    :: xx
  integer                   , intent(inout) :: jlo

  integer :: inc,jhi,jm
  logical :: ascnd

  ascnd=(xx(n).ge.xx(1))  ! True if ascending order of table, false otherwise.
  if (jlo <= 0 .or. jlo > n) then
  ! ... 
  do !Hunt is done, so begin the final bisection phase:
    if (jhi-jlo <= 1) then
      if (x== xx(n)) jlo=n-1
      if (x== xx(1)) jlo=1
      exit
    else
      jm=(jhi+jlo)/2
      if (x>= xx(jm) .eqv. ascnd) then
        jlo=jm
      else
        jhi=jm
      end if
    end if
  end do

 end subroutine rhunt

More general style suggestions:

  • When naming your variables and program units, always keep in mind that Fortran is a case-insensitive language (e.g. EditOrExit is the same as EditorExit.)

  • Use only characters in the Fortran character set. In particular, accent characters and tabs are not allowed in code, although they are usually OK in comments. If your editor inserts tabs automatically, you should configure it to switch off the functionality when you are editing Fortran source files.

  • Be generous with comments. State the reason for doing something, instead of repeating the Fortran logic in words.

  • Use the new and clearer syntax for LOGICAL comparisons, i.e.:

    == instead of .EQ.
    /= instead of .NE.
    > instead of .GT.
    < instead of .LT.
    >= instead of .GE.
    <= instead of .LE.
  • Where appropriate, simplify your LOGICAL assignments, for example:
IF (my_var == some_value) THEN
  something      = .TRUE.
  something_else = .FALSE.
ELSE
  something      = .FALSE.
  something_else = .TRUE.
END IF

! ... Better approach

something      = (my_var == some_value)
something_else = (my_var /= some_value)
  • Positive logic is usually easier to understand. When using an IF-ELSE-END IF construct you should use positive logic in the IF test, provided that the positive and the negative blocks are about the same length. It may be more appropriate to use negative logic if the negative block is significantly longer than the positive block. For example:
if (my_var /= some_value) then
  call do_this()
else
  call do_that()
end if
Better approach
if (my_var == some_value) then
  call do_that()
else
  call do_this()
end if
  • To improve readability, you should always use the optional space to separate the following Fortran keywords:
else if       end do           end forall   end function
end if        end interface    end module   end program
end select    end subroutine   end type     end where
select case
  • If you have a large or complex code block embedding other code blocks, you may consider naming some or all of them to improve readability.

  • If you have a large or complex interface block or if you have one or more sub-program units in the contains section, you can improve readability by using the full version of the end statement (i.e. end subroutine <name> or end function <name> instead of just end) at the end of each sub-program unit. For readability in general, the full version of the end statement is recommended over the simple end.

Example:

    over_species: do mf_ispecies = 1, mf_nspecies

      ! Due to symmetry we can loop over upper triangular matrix formed from level idices
      ! Eg for 4 levels - x marks pairs we have include in our computation:
      !
      !    1 2 3 4
      !  1 - x x x
      !  2 - - x x
      !  3 - - - x
      !  4 - - - -

      if (qrecion .and. mf_epf) then ! Calculate species specific constants  for Energy terms
        m_particle = mf_species(mf_ispecies)%awgt * AMU ! ~ in [g]
        m_en_ratio = M_ELECTRON / m_particle ! $m_e / m_n$ ~ [dimensionless]
      end if

      over_neutrals: do ilvl = 1, mf_species(mf_ispecies)%nlevel

        gidx_npos =  mf_species(mf_ispecies)%level_to_gidx(ilvl)   ! ilvl -> neutral level
        gidx_enpos = mf_species(mf_ispecies)%level_to_e_gidx(ilvl) ! ilvl -> energy equation

        over_ions: do jlvl = ilvl + 1, mf_species(mf_ispecies)%nlevel

          gidx_ipos  = mf_species(mf_ispecies)%level_to_gidx(jlvl)           ! jlvl -> ionized level
          gidx_eipos = mf_species(mf_ispecies)%level_to_e_gidx(jlvl)         ! jlvl -> energy equation

          has_rec_model = mf_species(mf_ispecies)%has_transitions(jlvl,ilvl) ! if positive model is ready to deal with rec j->i rates
          has_ion_model = mf_species(mf_ispecies)%has_transitions(ilvl,jlvl) ! if positive model is ready to deal with ion i->j rates

          ! ... long code here 

        end do over_ions
      end do over_neutrals
    end do over_species
  • Where possible, consider using CYCLE, EXIT or a WHERE-construct to simplify complicated DO-loops.

  • When writing a real literal with an integer value, put a 0 after the decimal point (i.e. 1.0 as opposed to 1.) to improve readability.

  • Where reasonable and sensible to do so, you should try to match the names of dummy and actual arguments to a SUBROUTINE / FUNCTION.

  • In an array assignment, it is recommended that you use array notations to improve readability. E.g.:

integer :: array1(10, 20), array2(10, 20)
integer :: scalar

! ...

array1 = 1
array2 = array1 * scalar

better approach

integer :: array1(10, 20), array2(10, 20)
integer :: scalar

! ...

array1(:, :) = 1
array2(:, :) = array1(:, :) * scalar
  • where appropriate, use parentheses to improve readability. e.g.:
a = (b * i) + (c / n)
! is easier to read than
a = b * i + c / n
  • to improve readability add brackets to function/subroutines colls without arguments. e.g.:
call update()
! is easier to read then
call update