This example  job computes the value of pi using a Monte Carlo method implemented  in Fortran 

The Fortran program has not been parallelised and the job  is designed to run in the ns serial queue on the Atos HPC.

Examples are provided both for using the Intel compiler and the GNU compiler.

#!/bin/bash
#
# **************************** LICENSE START ***********************************
#
# Copyright 2021 ECMWF. This software is distributed under the terms
# of the Apache License version 2.0. In applying this license, ECMWF does not
# waive the privileges and immunities granted to it by virtue of its status as
# an Intergovernmental Organization or submit itself to any jurisdiction.
#
# ***************************** LICENSE END ************************************
#
# pi-serial-intel    USER SERVICES  NOVEMBER 2021 - ECMWF 
#
#       This shell-script:
#
#               - Creates a Fortran program to compute the value
#                 of pi using a Monte Carlo method
#               - Compiles the program using the Intel compiler
#               - Runs the program in serial in the nf queue
#
#       This shell script produces the standard output file 
#
#         pi-serial-intel.<JOB-ID>.out
#
#       in the workdir directory (the user's $SCRATCH).
#
#-------------------------------
# setting options for SLURM
#-------------------------------
# Options that are specified within the script file should precede the
# first executable shell command in the file. 
# All options, which are case sensitive, are preceded by #SBATCH. 
# These lines become active only when the script is submitted
# using the "sbatch" command. 
# All job output is written to the workdir directory, by default. 

#SBATCH --qos=nf

        # Specifies that your job will run in the queue (Quality Of
        # Service - QoS) 'nf'. This is the default QoS and is suitable
	# for serial and small parallel jobs.

#SBATCH --job-name=pi-serial-intel

        # Assigns the specified name to the request

#SBATCH --output=%x.%j.out

        # Specifies the name and location of STDOUT where %x is the job
        # name and  %j is the job-id. The file will be # written in the
        # workdir directory if it is a relative path. If not given, the 
        # default is slurm-%j.out in the workdir.

#SBATCH --error=%x.%j.out

        # Specifies the name and location of STDERR where %x is the job
        # name and  %j is the job-id. The file will be # written in the
        # workdir directory if it is a relative path. If not given, the 
        # default is slurm-%j.err in the workdir.

#SBATCH --chdir=/scratch/...

        # Sets the working directory of the batch script before it is
        # executed.

#SBATCH --mail-type=FAIL
        
        # Specifies that an email should be sent in case the job fails.
        # Other options include BEGIN, END, REQUEUE and ALL (any state
        # change).

#SBATCH --time=00:05:00

        # Specifies that your job my run up to HH:MM:SS of wall clock
        # time. The job will be killed if it exceeds this limit. If
        # time is not defined, the default limit for the queue (qos)
        # will be used.

set -e

cd $SCRATCHDIR

# cat is used to read the Fortran program from the input stream
# and to write it to the file 'pi.f90'.
# cat reads line by line until it reaches the line EOF.

cat << EOF > pi.f90
!
! compute PI using a Monte Carlo simulation
!
program pi
  implicit none
  real*8             :: qrandom, x, y
  integer            :: iseed, inside, i
  integer, parameter :: npoints=100000000

  inside = 0
  iseed = 0
  do i = 1, npoints
    x = qrandom(iseed)
    y = qrandom(iseed)
    if ( x * x + y * y < 1.0 ) inside = inside + 1
  enddo
  print *,'pi is ',4.0d0 * dble(inside) / npoints
  stop
end program pi
!
! 'standard' random number generators are not thread safe.
! i.e. they can't be used within an OpenMP program.
!
! Here's a simple one that is.
!
! Returns a pseudo-random double precision value in the
! range [0,1) (i.e. 0 to not quite 1).
!
! I'm not sure how great the choices of a, c and m are although it has
! a period of about 260,000,000 which should be sufficient for our purposes.
! (note that a, c and m are prime)
!
function qrandom(iseed)
  implicit none
  real*8               :: qrandom
  integer              :: iseed
  integer*8, parameter :: a=874737473
  integer*8, parameter :: c=328592863
  integer*8, parameter :: m=2134884953

  if ( iseed < 0 ) then
    print *,'qrandom:  iseed is negative!!!!'
    stop
  endif

  iseed = mod( a * iseed + c, m )
  qrandom = dble(iseed) / dble(m)
  return
end function qrandom
EOF

module load prgenv/intel

ifort -O2 -march=core-avx2 pi.f90 -o pi_serial

./pi_serial

#-------------------------------
# tidy up and terminate
#-------------------------------

exit 0          # terminate the program, returning 0 (default) as return code
                # to the system

# There is one output file produced by this job: 
#
#         pi-serial-intel.$jobid.out
#
# in the working directory specified by the "--chdir" SBATCH option
#
# End of example job 'pi-serial-intel'

#!/bin/bash
#
# **************************** LICENSE START ***********************************
#
# Copyright 2021 ECMWF. This software is distributed under the terms
# of the Apache License version 2.0. In applying this license, ECMWF does not
# waive the privileges and immunities granted to it by virtue of its status as
# an Intergovernmental Organization or submit itself to any jurisdiction.
#
# ***************************** LICENSE END ************************************
#
# pi-serial-gnu      USER SERVICES  NOVEMBER 2021 - ECMWF 
#
#       This shell-script:
#
#               - Creates a Fortran program to compute the value
#                 of pi using a Monte Carlo method
#               - Compiles the program using the GNU compiler
#               - Runs the program in serial in the nf queue
#
#       This shell script produces the standard output file 
#
#         pi-serial-gnu.<JOB-ID>.out
#
#       in the workdir directory (the user's $SCRATCH).
#
#-------------------------------
# setting options for SLURM
#-------------------------------
# Options that are specified within the script file should precede the
# first executable shell command in the file. 
# All options, which are case sensitive, are preceded by #SBATCH. 
# These lines become active only when the script is submitted
# using the "sbatch" command. 
# All job output is written to the workdir directory, by default. 

#SBATCH --qos=nf

        # Specifies that your job will run in the queue (Quality Of
        # Service - QoS) 'nf'. This is the default QoS and is suitable
	# for serial and small parallel jobs.

#SBATCH --job-name=pi-serial-gnu

        # Assigns the specified name to the request

#SBATCH --output=%x.%j.out

        # Specifies the name and location of STDOUT where %x is the job
        # name and  %j is the job-id. The file will be # written in the
        # workdir directory if it is a relative path. If not given, the 
        # default is slurm-%j.out in the workdir.

#SBATCH --error=%x.%j.out

        # Specifies the name and location of STDERR where %x is the job
        # name and  %j is the job-id. The file will be # written in the
        # workdir directory if it is a relative path. If not given, the 
        # default is slurm-%j.err in the workdir.

#SBATCH --chdir=/scratch/...

        # Sets the working directory of the batch script before it is
        # executed.

#SBATCH --mail-type=FAIL
        
        # Specifies that an email should be sent in case the job fails.
        # Other options include BEGIN, END, REQUEUE and ALL (any state
        # change).

#SBATCH --time=00:05:00

        # Specifies that your job my run up to HH:MM:SS of wall clock
        # time. The job will be killed if it exceeds this limit. If
        # time is not defined, the default limit for the queue (qos)
        # will be used.

set -e

cd $SCRATCHDIR

# cat is used to read the Fortran program from the input stream
# and to write it to the file 'pi.f90'.
# cat reads line by line until it reaches the line EOF.

cat << EOF > pi.f90
!
! compute PI using a Monte Carlo simulation
!
program pi
  implicit none
  real*8             :: qrandom, x, y
  integer            :: iseed, inside, i
  integer, parameter :: npoints=100000000

  inside = 0
  iseed = 0
  do i = 1, npoints
    x = qrandom(iseed)
    y = qrandom(iseed)
    if ( x * x + y * y < 1.0 ) inside = inside + 1
  enddo
  print *,'pi is ',4.0d0 * dble(inside) / npoints
  stop
end program pi
!
! 'standard' random number generators are not thread safe.
! i.e. they can't be used within an OpenMP program.
!
! Here's a simple one that is.
!
! Returns a pseudo-random double precision value in the
! range [0,1) (i.e. 0 to not quite 1).
!
! I'm not sure how great the choices of a, c and m are although it has
! a period of about 260,000,000 which should be sufficient for our purposes.
! (note that a, c and m are prime)
!
function qrandom(iseed)
  implicit none
  real*8               :: qrandom
  integer              :: iseed
  integer*8, parameter :: a=874737473
  integer*8, parameter :: c=328592863
  integer*8, parameter :: m=2134884953

  if ( iseed < 0 ) then
    print *,'qrandom:  iseed is negative!!!!'
    stop
  endif

  iseed = mod( a * iseed + c, m )
  qrandom = dble(iseed) / dble(m)
  return
end function qrandom
EOF

module load pg

gfortran -O2 -march=core-avx2 pi.f90 -o pi_serial

./pi_serial

#-------------------------------
# tidy up and terminate
#-------------------------------

exit 0          # terminate the program, returning 0 (default) as return code 
                # to the system 

# There is one output file produced by this job: 
#
#         pi-serial-gnu.$jobid.out
#
# in the working directory specified by the "--chdir" SBATCH option
#
# End of example job 'pi-serial-gnu'