GEOS–Chem v9–01–02 Online User's Guide

Previous | Next | Printable View (no frames)

Appendix 7: GEOS–Chem Style Guide

In this Appendix, we shall provide our recommended programming style and coding techniques. We ask that you adhere to this when writing new GEOS–Chem source code.

A7.1 Background

(1) The development language of GEOS–Chem is Fortran 90 (F90). This is a recent update of the Fortran programming language. Fortran 90 introduces several new features over its predecessor, Fortran 77 (F77), including allocatable arrays, modules, and new statements. For more information about Fortran 90 you can see the following tutorials:

F90 is totally backwards-compatible with F77. Older code (a.k.a. "legacy code") which adheres to the F77 standard should compile and run under F90 with no problems. Although we have written new GEOS–Chem code using the new language features of F90, many 3rd-party routines included in GEOS–Chem (such as TPCORE, SMVGEAR, FAST–J, and ISORROPIA) were written to the Fortran 77 standard. Using F90 ensures that both the older and newer routines of GEOS–Chem can be smoothly combined together.

In recent years there have also been further updates to Fortran, namely Fortran 95 (F95), and Fortran 2000 (F2K). These new Fortran versions include the new language features which are present in F90, but generally do not include the old, obsolete features found in F77. All modern Fortran compilers are compatible with F77, F90, F95, and F2K, and can automatically invoke the F77 backwards-compatibility functionality in a manner completely transparent to the user.

(2) We write the GEOS–Chem source code using Fortran-90 fixed-format layout. We find that this is the best way to provide backwards compatibility with GEOS–Chem's older "legacy-code" routines.

The main features of the fixed-format layout are:

(3) We usually denote people by their 3-letter email abbreviations (e.g. bmy = Bob Yantosca, mje = Mat Evans, etc). This is needed when we initial and date various model revisions, e.g, (bmy, 5/1/03). Some newer comments use more characters for the comment (e.g. ccarouge = Claire Carouge, etc).

A7.2 Headers, declarations, indentations, and white space

(4) We generate detailed GEOS–Chem reference documentation with ProTeX. ProTeX is a perl script developed by NASA/GMAO that can generate a LaTeX file from GEOS–Chem's subroutine, function, and module comment headers. The LaTeX file can then be compiled to produce both PostScript and PDF output.

To use ProTeX, subroutines, functions and modules should contain a standard header such as:

!------------------------------------------------------------------------------
!            Harvard University Atmospheric Chemistry Modeling Group 
!------------------------------------------------------------------------------
!BOP
!
! !IROUTINE: metero
!
! !DESCRIPTION: Subroutine METERO calculates meteorological constants needed for the
! dry deposition velocity module. (lwh, gmg, djj, 1989, 1994; bmy, 11/21/02)
!
! Arguments as Output:
! ============================================================================
! (1 ) CZ1 (REAL*8) : Midpoint height of first model level [m]
! (2 ) TC0 (REAL*8) : Array for grid box surface temperature [K]
! (3 ) OBK (REAL*8) : Array for the Monin-Obhukov length [m]
!
!\subsection*{Reference}
! (1 ) Wesely, M. L., 1989. 
! (2 ) Jacob, D.J., and S.C. Wofsy, 1990
!
!\\
!\\
! !INTERFACE: 
!
      SUBROUTINE METERO( CZ1, TC0, OBK )
!
! !USES:
! 
      USE DAO_MOD, ONLY : AIRDEN, HFLUX, T, TS, USTAR
      USE PRESSURE_MOD, ONLY : GET_PEDGE 
      
      IMPLICIT NONE

#     include "CMN_SIZE" ! Size parameters
#     include "CMN_GCTM" ! Physical constants
!
! ! OUTPUT PARAMETERS:
!
      REAL*8, INTENT(OUT) :: CZ1(MAXIJ)
      REAL*8, INTENT(OUT) :: TC0(MAXIJ)
      REAL*8, INTENT(OUT) :: OBK(MAXIJ)
!
! !REVISION HISTORY: 
!  (1 ) Now reference GET_PEDGE from "pressure_mod.f". Now reference T from 
!   "dao_mod.f". Removed obsolete code & comments, and added new 
!  documentation header. Now force double precision with "D" 
!  exponents. Now compute OBK here as well. Bundled into F90 module
!  "drydep_mod.f" (bmy, 11/20/02)
! 02 Mar 2011 - J. Fisher - Set aerosol dry deposition velocity to 0.03 cm/s
! over snow and ice based on Nilsson & Rannik, 2001
! 02 Mar 2011 - R. Yantosca - Added ProTeX headers !EOP !------------------------------------------------------------------------------ !BOC ! ! ! LOCAL VARIABLES: ! INTEGER :: I, J, IJLOOP REAL*8 :: P1, P2, THIK, NUM, DEN ! ! !DEFINED PARAMETERS: ! REAL*8, PARAMETER :: KAPPA = 0.4d0 ! Von Karman's Constant REAL*8, PARAMETER :: CP = 1000.0d0 ! ! ! EXTERNAL FUNCTIONS: ! REAL*8, EXTERNAL :: XLTMMP !================================================================= ! METERO begins here! !================================================================= ... code goes here ... END SUBROUTINE METERO !EOC

Note that the ProTeX header contains the following:

  1. Indicators (!BOP, !EOP) as to where ProTeX should look for source code header comments
  2. Indicators (!BOC, !EOC) as to where the source code begins and ends. ProTeX will ignore code between !BOC and !EOC unless you tell it otherwise.
  3. A description of the routine, the original date and most recent
  4. Description of when the routine was last modified (with date & name of the programmers)
  5. A list of input & output arguments
  6. A list of references (if applicable)
  7. References to F90 modules (if any) after the USE keyword followed by the IMPLICIT NONE declaration
  8. Any #include declarations for header files. Since we are using the C-preprocessor, the # character must go in column 1.
  9. Modification NOTES. Each time the file is modified the REVISION HISTORY section should be updated.
  10. Declarations for local variables
  11. Declarations for Fortran parameters (aka constants)
  12. Declarations for external functions (if necessary)

Fortran 90 allows you to declare an argument to a subroutine or function with one of the following descriptors:

  1. INTENT(IN) (read-only),
  2. INTENT(OUT) (write-only) or
  3. INTENT(INOUT) (read-and-write)

This helps you to prevent overwriting arguments which should not be overwritten.

We recommend that you separate sections of source code with one or more divider comments:

!=================================================================
! 1st level header
!=================================================================

   !--------------------------------------------------------------
   ! 2nd level header
   !--------------------------------------------------------------

      !-----------------------
      ! 3rd-level header
      !----------------------

You can line up each header with the the code that follows immediately below it. You can mix and match these styles as you wish. You don't always have to extend the header all the way across the screen.

We recommend that you extend major section headers out so that they end at column #72. This will provide a convenient reminder for the last allowable column of executable code in fixed-format Fortran.

(5) Write GEOS–Chem source code in ALL CAPITALS. This prevent confusion between characters (e.g. the lowercase l and the number 1, which in some fonts can look a lot alike). Good editors such as emacs and Xemacs can "colorize" the different words in a program (e.g. statements, variables, and quoted text are rendered in different colors) and thus makes the code more readable.

Because GEOS–Chem contains as significant amount of 3rd-party code, you will see instances of source code routines written in lowercase letters. You can leave this code as-is, because (1) this code comes to us from outside sources, and (2) it can be cumbersome to convert already-written code to uppercase.

In some special instances, you may want to mix uppercase and lowercase letters, particularly if you want to improve readability or if you are using scientific abbreviations in a name (e.g. Rn_Pb_Be_mod.F).

(6) Indent each block of code such that each new block is placed 3 columns deeper than the last block. For example:

         1         2         3         4
1234567890123456789012345678901234567890
      IF ( X == 0 ) THEN
         PRINT*, 'Hello, X is 0!'

         IF ( Y == 0 ) THEN
            PRINT*, 'Hello, X is 0 and Y is also 0!'
         ENDIF
      ENDIF

A good editor such as emacs will indent for you automatically. You can specify the level of indent in your editor setup file (e.g. .emacs).

(7) Omit indentation for each line of a multiple DO loop. This makes the code more readable. Use this syntax:

DO N = 1, NNPAR
DO L = 1, LLPAR
DO J = 1, JJPAR
DO I = 1, IIPAR
   A(I,J,L,N) = A(I,J,L,N) / 2.0d0
ENDDO
ENDDO
ENDDO
ENDDO

instead of letting the DO loops "creep" across the screen:

DO N = 1, NNPAR
   DO L = 1, LLPAR
      DO J = 1, JJPAR
         DO I = 1, IIPAR
            A(I,J,L,N) = A(I,J,L,N) / 2.0d0
         ENDDO
      ENDDO
   ENDDO
ENDDO

(8) Use LOTS of white space to separate code. For example instead of

IF(X.eq.0)CALL DOLOOP(X,Y,Z,A,B,C,1,2,3)

type:

IF ( X == 0 ) CALL DOLOOP( X, Y, Z, A, B, C, 1, 2, 3 )

This makes the code much more readable, which helps to prevent bugs.. After all, if you can't read what you've written, you will never find your mistake! Don't be afraid to continue the statement to the next line if necessary.

In the above example, we also have used the new F90 equality test operator == (see below for more information).

(9) Line up quantities in assignment statements, e.g. instead of:

A=1
THISLOOP=2
C=3
D=THISLOOP*C - A

it is much more more readable to type:

A        = 1
THISLOOP = 2
C        = 3 
D        = ( THISLOOP * C ) - A

Use parentheses to ensure that terms will be grouped correctly. The Fortran compiler will always evaluate expressions in the order the parentheses are placed, from innermost to outermost.

(10) Do not leave white space between array indices:

A = X(I,J,L)

Leave some white space between arguments of functions and WRITE statements:

A = MYFUNCTION( I, J, L )
WRITE( 6, '(a)' ) A

(11) Arguments in subroutine or function calls should be lined up. If there are an even number of arguments then these can be broken up symmetrically:

 CALL MYSUB( THIS_A, THIS_B, THIS_C,
&            THIS_D, THIS_E, THIS_F )

or

 X = MYFUNCTION( THIS_A, THIS_B, THIS_C,
&                THIS_D, THIS_E, THIS_F )

A7.3 Numeric and character data types

(12) GEOS–Chem uses the following data types:

LOGICAL            --> TRUE or FALSE switch
INTEGER            --> 32 byte integer value
REAL*4             --> 32 byte floating-point value
REAL*8             --> 64 byte floating-point value
CHARACTER(LEN=255) --> string w/ 255 characters (the max limit)

where

INTEGER can express numbers from -2e9   to +2e9
REAL*4  can express numbers from -1e38  to +1e38
REAL*8  can express numbers from -1e312 to +1e312

Always use REAL*8 as the default floating-point type when writing new GEOS–Chem code. This will prevent round off and precision truncation errors. However, there are some instances when you can use the REAL*4 data type instead of REAL*8:

Some third-party routines used by GEOS–Chem (such as TPCORE) use the REAL datatype. On Cray machines the type REAL is assumed to be REAL*8. However, most other platforms/compiler combinations assume that REAL is REAL*4. You must use a special compiler switch (-r8 for SGI and Alpha, or -xtypemap=real:64 for Sun/Sparc) when compiling in order to force variables of type REAL to be interpreted as REAL*8. This is already handled for you in the Makefile.

NOTE: See the GEOS–Chem wiki page on floating point math issues for more information.

(13) Always use D (double-precision) exponents when assigning a value to a REAL*8 variable. Always use this syntax:

REAL*8 :: PI
PI = 3.14159265358979323d0

instead of this syntax:

REAL*8 :: PI
PI = 3.14159265359e0

In F90, the E (single-precision) exponent only yields about 7 decimal places of precision, whereas the D exponent yields 15 or 16 decimal places of precision. Using D exponents prevents roundoff errors.

(14) Declare strings with 255 characters, even if you don't know a priori how long the string will be. You can always use the TRIM statement to strip of excess white space at the end. For example:

CHARACTER(LEN=255) :: STR

STR = 'I am the very model of a modern major general...'

WRITE( 6, '(a)' ) TRIM( STR )

A7.4 Converting from number to character (and vice versa)

(15) In many GEOS–Chem routines, you must create a character variable for a file name which contains the current date. This requires you to be able to convert from INTEGER to CHARACTER. You can do this with a Fortran internal write statement. This is similar to a regular WRITE statement, only that the unit number is replaced by a character variable. Think of it as "writing" your number into a character variable instead of to a file.

For example, the following code will make a string containing the date 20010701:

! Declare variables
CHARACTER(LEN=8) :: DATE_STR
INTEGER :: DATE = 20010701

! FORTRAN internal write -- converts number to string
WRITE( DATE_STR, '(i8)' ) DATE

which can then be used in a file name, etc.

Caveat: The Fortran internal write depicted above may not be supported on the PGI Linux compiler. Instead, you can use the ENCODE function to convert from number to character, as follows:

! Declare variables
CHARACTER(LEN=8) :: DATE_STR
INTEGER :: DATE = 20010701

! Writes value from DATE into DATE_STR
ENCODE( 8, '(i8)', DATE_STR ) DATE

You must specify the length of the character variable (in this case, 8), the format string, and the character variable when calling ENCODE.

(16) If you want to extract numbers from a character string and store them into an integer variable, we recommend that you do a Fortran internal read. This is the reverse of the FORTRAN internal write described above:

! Declare variables
CHARACTER(LEN=8) :: DATE_STR = '20010101'
INTEGER :: DATE

! FORTRAN internal read -- converts string to number
READ( DATE_STR, '(i8)' ) DATE

The PGI Linux compiler may not support Fortran internal reads. Instead, you can use the DECODE function (which is the reverse of the ENCODE function described above). DECODE is called as follows:

! Declare variables
CHARACTER(LEN=8) :: DATE_STR = '20010101'
INTEGER :: DATE

! Reads string value from DATE_STR into DATE
DECODE( 8, '(i8)', DATE_STR ) DATE

As with ENCODE, you must specify the length of the character variable (in this case, 8), the format string, and the character variable in the call to DECODE.

(17) In practice, when writing GEOS–Chem code, you may have to use both the internal read / internal write and ENCODE / DECODE methods simultaneously, separated by the appropriate C–preprocessor compiler flags.

! Declare variables
CHARACTER(LEN=8) :: DATE_STR
INTEGER :: DATE = 20010701
 
# include "define.h" ! C-preprocessor flags
 
#if defined( LINUX_PGI )
   ! Write numeric value from DATE into DATE_STR for PGI Linux
   ENCODE( 8, '(i8)', DATE_STR ) DATE
#else
   ! FORTRAN Internal Write for other platforms
   WRITE( DATE_STR, '(i8)' ) DATE
#endif

A7.5 New language features of Fortran 90

(18) Always use the new F90-style operators instead of the older F77-style operators:

== instead of .EQ.
/= instead of .NE.
>  instead of .GT.
>= instead of .GE.
<  instead of .LT.
>= instead of .LE.

These are easier to read and take up only one or two spaces instead of four spaces.

The operators .AND. .OR. .NOT. are still the same in F90 as in F77.

(19) Always use the F90 comment character ! instead of the F77 comment character C. You can place the ! comment in any column of code, which produces much more legible comments.

(20) Use the official F90 continuation character, &, in column 6. This will make it easier to merge free-format and fixed-format code if necessary.

(21) Use the F90 array assignment functionality to assign a value to all elements of an array without using array indices. For example:

INTEGER :: ARRAY(10,10)
ARRAY(:,:) = 0

You can also leave off the default array mask (:,:), so

ARRAY = 0

is also acceptable.

(22) Always initialize PARAMETER constants on the same line where they are declared. Use this syntax:

REAL*8, PARAMETER :: PI = 3.14159265358979323d0

instead of this obsolete F77 syntax:

REAL*8 PI
PARAMETER( PI = 3.14159265358979323 )

(23) Always use F90 direct assignment statements instead of the obsolete F77 DATA statements. Use this syntax:

REAL*8 :: A(2) = (/ 1.0d0, 2.0d0 /)

instead of this obsolete F77 syntax:

REAL*8 A(2)
DATA A / 1.0, 2.0 /

(24) Always include the SAVE attribute on the same line where the variable is declared. Use this syntax:

LOGICAL, SAVE :: FIRSTTIME = .TRUE.

instead of this obsolete F77 syntax:

LOGICAL FIRSTTIME 
SAVE FIRSTTIME
DATA FIRSTTIME / .TRUE. /
Note: SAVEd variables within a subroutine or function keep their values from one call to the next. This allows you to do some initialization only on the first call to a subroutine, for example:
LOGICAL, SAVE :: FIRSTTIME = .TRUE. 

IF ( FIRSTTIME ) THEN
   PRINT*, 'This is the first time this routine is called!'
   FIRSTTIME = .FALSE.
ENDIF

Since the value of FIRSTTIME is saved between calls, the sentence above will only be printed out on the first call to the subroutine.

(25) Always include the EXTERNAL attribute on the same line where a function's type is declared. Use this syntax:

INTEGER, EXTERNAL :: MYFUNC

instead of this obsolete F77 syntax.

INTEGER MYFUNC
EXTERNAL MYFUNC

(26) Use the F90 SELECT CASE statement to pick from a list of options wherever expedient.

In F77 you would have had to write:

IF ( X .EQ. 1 ) THEN
   CALL MYSUB( X1 )
ELSE IF ( X .EQ. 2 ) THEN 
   CALL MYSUB( X2 )
ELSE
   CALL MYSUB( X0 )
ENDIF

but in F90, you may write this as:

SELECT CASE ( X )
   CASE( 1 )
      CALL MYSUB( X1 )
   CASE( 2 ) 
      CALL MYSUB( X2 )
   CASE DEFAULT
      CALL MYSUB( X0 )
END SELECT

SELECT CASE also works with character constants:

SELECT CASE ( TRIM( TRACERNAME ) )
   CASE( 'NOx' )
      CALL MYSUB( 1 )
   CASE( 'Ox' )
      CALL MYSUB( 2 )
   CASE DEFAULT
      CALL MYSUB( 3 )
END SELECT

A note of warning: You cannot specify CASE( X ) where X is a variable, or else the compiler will balk. In that case you must resort to the IF - ELSE IF block structure as shown above.

(27) Always use the F90 intrinsic functions LOG( X ) and LOG10( X ) to take the natural and common logarithms of X. Here X may be declared either as REAL*4 or REAL*8.

The F77 intrinsic functions ALOG( X ) and ALOG10( X ) are not supported on some compilers, such as on the Alpha platform. This can lead to a compile time error.

(28) You may use F90's keyword argument capability to make subroutine calls more clear. For example, the subroutine:

 SUBROUTINE READ_A1( NYMD,     NHMS, 
&                    ALBEDO,   CLDTOT,   EFLUX,    EVAP,    
&                    ... etc ...                         )

can be called like this:

      CALL READ_A1( NYMD     = NYMD, 
     &              NHMS     = NHMS,
     &              ALBEDO   = ALBD,
     &              ... etc ...      )

In each pair separated by an equals sign such as:

ALBEDO   = ALBD,

The name on the left of the equals sign is the name of the argument as defined in the subroutine (i.e. ALBEDO). This is often called the "dummy argument" because it is just a name for the memory that is getting passed to it from outside of the subroutine.

The name on the right of the equals sign is the name of the variable that we are passing down to the subroutine. Here we are using the ALBD variable (from GEOS–Chem module dao_mod.F and passing that down to READ_A1 as the ALBEDO argument.

A7.6 Math optimizations

(29) Logarithms and exponentiation are the most computationally expensive mathematical functions. Therefore, if you wish to code a polynomial expression, do not use the traditional formulation, i.e.:

 Y = A0 + A1*X + A2*X**2 + A3*X**3 + A4*X**4 + A5*X**5 

but instead write the expression as:

 Y = A0 + X* ( A1 + X* ( A2 + X* ( A3 + X* ( A4 + X* ( A5 )))))

This modified polynomial expression has the benefit of eliminating the costly exponentiations. A Nth order polynomial will now only require N multiplications. This will execute much faster than the traditional formulation.

(30) Although it is tempting to replace EXP( –x ) with the first-order approximation 1 – x, this can result in negative values for certain values of x. This can cause negative tracer concentrations.

A7.7 DO loops

(31) Always use the CYCLE statement to skip to the next iteration in a DO-loop. Use this syntax:

DO I = 1, 1000

   ! Skip to next iteration
   IF ( X(I) == 0 ) CYCLE
 
   ! Print
   PRINT*, X(I)
ENDDO

instead of this obsolete F77 syntax:

DO 100 I = 1, 1000
   IF ( X(I) .EQ. 0 ) GOTO 100
   PRINT*, X(I)
100 CONTINUE

Notice that the F90 ENDDO statement replaces the labeled CONTINUE statement. This results in cleaner code.

(32) Always use the EXIT statement to completely exit from a DO loop. Use this syntax:

DO I = 1, 1000
 
   ! Break out of this loop if X==0!
   IF ( X(I) == 0 ) EXIT
 
   ! Print info if we have not exited
   WRITE( 6, '(''Iteration: '', i5, f13.6)' ) I, X(I)
ENDDO 

PRINT*, 'outside loop'

instead of this obsolete F77 syntax:

C Do-loop
DO 100 I = 1, 1000
   IF ( X(I) .EQ. 0 ) GOTO 200
   PRINT*, 'inside loop'
100 CONTINUE

C Continue outside loop
200 CONTINUE
PRINT*, 'outside loop'

(33) Structure your DO loops so that they access memory in the most optimal manner. Fortran is a column-major language, which means that arrays are stored in memory by columns first, then by rows. If you have declared an array such as:

INTEGER :: I, J, L, N
REAL*8  :: ARRAY(IIPAR,JJPAR,LLPAR,NNPAR)

then for optimal efficiency, the leftmost dimension (I) needs to vary the fastest, and needs to be accessed by the innermost DO-loop. Then the next leftmost dimension (J) should be accessed by the next innermost DO-loop, and so on.

Therefore, the proper way to loop over this array is:

DO N = 1, NNPAR
DO L = 1, LLPAR
DO J = 1, JJPAR
DO I = 1, IIPAR
   ARRAY(I,J,L,N) = ARRAY(I,J,L,N) * 2.0d0
ENDDO
ENDDO
ENDDO
ENDDO 

Note that the I index is varying most often, since it is the innermost DO-loop, then J, L, and N. This is opposite to how a car's speedometer varies (i.e. rightmost index varies the fastest).

If you loop through an array in this fashion, with leftmost indicies varying fastest, then the code minimizes the number of times it has to load subsections of the array into cache memory. In this optimal manner of execution, all of the array elements sitting in the cache memory are read in the proper order before the next array subsection needs to be loaded into the cache. However, if you step through array elements in the wrong order, the number of cache loads is proportionally increased. Since it takes a finite amount of time to reload array elements into cache memory, the more times you have to access the cache, the longer it will take the code to execute. This can slow down the code dramatically.

(34) Use "infinite" DO loops in situations where you need to iterate for an unknown number of iterations before some exit criterion is reached. This is most often used to read data from a file whose length you do not know a priori.

Here is an example of an infinite DO loop used to read data from disk:

INTEGER :: I, IOS, IUNIT
REAL*8 :: A(10000) 

! Open file
OPEN( IUNIT, FILE='myfile.txt', IOSTAT=IOS )

! Quit if we can't open the file
IF ( IOS /= 0 ) THEN
   PRINT*, 'Error opening the file!'
   STOP
ENDIF

! Infinite DO loop for reading data from "myfile.txt"
DO

   ! Read I and A(I) from the file
   READ( IUNIT, '(i3,f11.3)', IOSTAT=IOS ) I, A(I)
 
   ! IOS < 0 is end-of-file, so we exit the loop
   IF ( IOS < 0 ) EXIT
	
   ! If IOS < 0, then this is an I/O error,
   ! so we print an error msg and quit
   IF ( IOS > 0 ) THEN
      PRINT*, 'I/O error reading from file!
      STOP
   ENDIF
ENDDO
 
! Close the file after we exit the loop
CLOSE( IUNIT )

In this example, the DO loop wants to execute forever, but it is stopped by an external criterion—the end of the file. The OPEN, READ, and CLOSE functions can return the I/O status to a variable if you specify via the IOSTAT keyword. If the I/O status variable (in this case, named IOS) is negative, then this is a normal end-of-file condition, and so we can just exit the loop and close the file. However, if IOS is a nonzero, positive number, then this means that we have encountered a real error condition.

(35) Reserve the following variable names for DO-loop and array indices:

Therefore, you can assume that any DO-loop that uses I, J, and L is looping over grid box longitude, latitude, and altitude dimensions, and that any DO-loop which uses N is looping over tracer, species, or some other quantity. These special indices should NOT be used to refer other quantities. In this way, if you should get an error message such as:

 Error at grid box (I,J) = (23,34)!

you will immediately understand (I,J) are the longitude and latitude indices of the box where the error happened.

NOTE: In some 3rd-party routines, you will see K used instead of L to denote the altitude index. We recommend leaving these as-is, so as to avoid having to rewrite entire sections of 3rd-party code. However, you should use L for the altitude index in any new GEOS–Chem code that you write.

A7.8 Fortran-90 modules

(36) Place all of your new GEOS–Chem code into Fortran 90 modules, if possible. These allow you to save both variables and routines within a single program unit.

A sample module is provided below.

!------------------------------------------------------------------------------
!          Harvard University Atmospheric Chemistry Modeling Group            !
!------------------------------------------------------------------------------
!BOP
!
! !MODULE: MY\_MODULE
!
! !DESCRIPTION: \subsection*{Overview}
! Module MY\_MODULE contains variables and routines to do something.
! (bmy, 5/1/03)
!\\
!\\
!
! !INTERFACE
       MODULE MY_MODULE
!
! !USES:
! 
      IMPLICIT NONE

 ! Make everything private...
      PRIVATE
 !
 ! !PUBLIC MEMBER FUNCTIONS: 
 !
      PUBLIC :: DRIVER
      PUBLIC :: INIT_MY_MODULE
      PUBLIC :: CLEANUP_MY_MODULE
 !
 ! !REVISION HISTORY:
 !  (1 ) Bug fix in MY_SUBROUTINE (bmy, 5/1/03)
 !EOP
 !------------------------------------------------------------------------------
 !
 ! !PRIVATE DATA MEMBERS:
 !
      ! Argument to the SIN function
      REAL*8, ALLOCATABLE :: B(:)

      CONTAINS

!------------------------------------------------------------------------------
!          Harvard University Atmospheric Chemistry Modeling Group            !
!------------------------------------------------------------------------------
!BOP
!
! !IROUTINE: DRIVER
!
! !DESCRIPTION: \subsection*{Overview}
! Subroutine DRIVER is the driver routine for MY_MODULE. (bmy, 5/1/03)
!\\
!\\
!
! !INTERFACE:
      SUBROUTINE DRIVER
!
! !REVISION HISTORY:
!
!EOP
!------------------------------------------------------------------------------
!BOC
!
! !LOCAL VARIABLES:
!
      LOGICAL, SAVE :: FIRST = .TRUE.
      INTEGER, :: I
      REAL*8 :: X, Y1, Y2
      REAL*8, PARAMETER :: PI180 = 180d0 / 3.14159265358979323d0
 
      !=================================================================
      ! DRIVER begins here!
      !=================================================================
 
      ! Initialize 
      IF ( FIRST ) THEN
         CALL INIT_MY_MODULE
      ENDIF
 
      !=================================================================
      ! Loop over degrees
      !=================================================================
      DO I = 1, 360
 
         ! Convert to radians
         B = DBLE( I ) / PI180
 
         ! Call subroutine
         CALL MY_SUBROUTINE( B(I), Y1 )
 
         ! Call function
         Y2 = MY_FUNCTION( B(I) )
 
         ! Write output
         WRITE( 6, '(i3, 1x, 3(f13.6,1x))' ) I, B(I), Y1, Y2
      ENDDO 
 
      ! Return to calling program
      END SUBROUTINE DRIVER
!EOC
!-----------------------------------------------------------------------------
!          Harvard University Atmospheric Chemistry Modeling Group            !
!------------------------------------------------------------------------------
!BOP
!
! !IROUTINE: MY\_SUBROUTINE
!
! !DESCRIPTION: \subsection*{Overview}
! Function MY_SUBROUTINE returns the sine of a number. (bmy, 5/1/03)
!\\
!\\
!
! !INTERFACE:

      SUBROUTINE MY_SUBROUTINE( X, Y )
!
! !INPUT PARAMETERS: 
!
      ! Argument for the sine function
      REAL*8, INTENT(IN) :: X
!
! !RETURN VALUE:
!     
      REAL*8, INTENT(OUT) :: Y 
!
! !REVISION HISTORY:
! (1 ) Corrected typographic error (bmy, 5/1/03)
!EOP
!------------------------------------------------------------------------------
!BOC
 
      !=================================================================
      ! MY_SUBROUTINE begins here!
      !=================================================================
      Y = SIN( X )
 
      ! Return to calling program
      END SUBROUTINE MY_SUBROUTINE
!EOC
!-----------------------------------------------------------------------------
!          Harvard University Atmospheric Chemistry Modeling Group            !
!------------------------------------------------------------------------------
!BOP
!
! !IROUTINE: MY\_FUNCTION
!
! !DESCRIPTION: \subsection*{Overview}
! Function MY_FUNCTION returns the sine of a number. (bmy, 5/1/03)
!\\
!\\
!
! !INTERFACE:

      FUNCTION MY_FUNCTION( X ) RESULT( Y )
!
! !INPUT PARAMETERS: 
!
      ! Argument for the sine function
      REAL*8, INTENT(IN) :: X
!
! !RETURN VALUE:
!     
      REAL*8, INTENT(OUT) :: Y 
!
! !REVISION HISTORY:
!
!EOP
!------------------------------------------------------------------------------
!BOC
 
      !=================================================================
      ! MY_FUNCTION begins here!
      !=================================================================
      Y = SIN( X )
 
      ! Return to calling program
      END FUNCTION MY_FUNCTION
!EOC
!------------------------------------------------------------------------------
!!          Harvard University Atmospheric Chemistry Modeling Group            !
!------------------------------------------------------------------------------
!BOP
!
! !IROUTINE: INIT\_MY\_MODULE
!
! !DESCRIPTION: \subsection*{Overview}
! Subroutine INIT\_MY\_MODULE allocates and zeroes the B array. (bmy, 5/1/03)
!\\
!\\
!
! !INTERFACE:
      SUBROUTINE INIT_MY_MODULE
!
! !USES
!
      USE ERROR_MOD, ONLY : ALLOC_ERR
!
! !REVISION HISTORY:
!
!EOP
!------------------------------------------------------------------------------
!BOC
!
! !LOCAL VARIABLES
      INTEGER :: AS
 
      !=================================================================
      ! INIT_MY_MODULE begins here!
      !=================================================================
 
      ! Allocate B, check status
      ALLOCATE( B(360), STAT=AS )
 
      ! Error check
      IF ( AS /= 0 ) CALL ALLOC_ERR( 'B' )
 
      ! Zero B
      B(:) = 0d0
 
      ! Return to calling program
      END SUBROUTINE INIT_MY_MODULE
!EOC
!-----------------------------------------------------------------------------
!!          Harvard University Atmospheric Chemistry Modeling Group            !
!------------------------------------------------------------------------------
!BOP
!
! !IROUTINE: CLEANUP\_MY\_MODULE
!
! !DESCRIPTION: \subsection*{Overview}
! Subroutine CLEANUP_MY_MODULE deallocates all module arrays. (bmy, 5/1/03)
!\\
!\\
!
! !INTERFACE:
      SUBROUTINE CLEANUP_MY_MODULE 
!
! !REVISION HISTORY:
!
!EOP
!------------------------------------------------------------------------------
!BOC

      !=================================================================
      ! CLEANUP_MY_MODULE begins here!
      !=================================================================
      IF ( ALLOCATED( B ) ) DEALLOCATE( B )

      ! Return to calling program
      END SUBROUTINE CLEANUP_MY_MODULE
!EOC
 END MODULE MY_MODULE

Fortran 90 modules written for GEOS–Chem contain the following features:

  1. Module header: This is similar to the subroutine header. Contains a list of all module variables, routines, and modification notes, in ProTeX style.

  2. PRIVATE declarations: You can "hide" certain routines or variables within a F90 module from other routines or modules. Ideally your module will be like a "black box" with only a few routines or variables that are publicly accessible. By convention, in GEOS–Chem, we declare everything PRIVATE first then indicate the public routines.

  3. IMPLICIT NONE: If you declare IMPLICIT NONE once at the top of the module, then you don't have to declare it in the individual subroutines; the declaration "carries through" below.

  4. Module variables and arrays: Any variables declared above the CONTAINS statement are as if they were stored in F77 common blocks; that is, they are preserved between calls to the various module routines. Also, you should declare module arrays as ALLOCATABLE whenever possible. This allows you to only set aside memory for that array if a certain routine is called. This results in more efficient memory management.

  5. CONTAINS statement: All module variables and PRIVATE declarations go above the CONTAINS statement. Module routines and functions go below the CONTAINS statement.

  6. Subroutines and functions: Your module will contain various subroutines and functions depending on the particular needs at hand.

  7. An INIT routine: Your module should have a subroutine named INIT_modulename, which initializes and allocates all module arrays.

  8. A CLEANUP routine: Your module should have a subroutine named CLEANUP_modulename, which deallocates module arrays. (You only need this if you have allocatable arrays).

F90 modules are very similar to classes in Java and C++; they allow you to group data and the routines which work on the data in a single package.

Theoretically, each GEOS–Chem routine and/or variable should belong to a module. In practice, this is not true, as we do have common blocks and separate subroutine and function files from some of the older routines which were written by 3rd party sources. However, you should write all new GEOS–Chem code in module form.

Also, it is OK for modules to reference other modules. If Module B references Module A, all that you have to do is to make sure that Module B declaration in the Makefile refers to Module A. The GEOS–Chem Makefiles have been prepared for you by the GEOS–Chem Support Team, so in most instances, you will not have to worry about this.

NOTE: With the development of the column code, we strongly recommend that new modules, subroutines and functions (even in existing modules) do not refer to different modules for gridded variables (dependent on longitude, latitude, altitude) and functions. These should be passed through the argument list instead.

Previous | Next | Printable View (no frames)