Fortran with coarrays¶
Fortran is perhaps the only programming language with standard support for distributed memory programming: Since the 2008 standard, coarrays are a part of the language standard.
On DelftBlue, coarrays are supported via the GNU fortran compiler and an additional module built on top of MPI.
- Load necessary modules:
- Compile and link your code using these flags (in addition to any others, like
-O3
or-g
): - To run your code interactively (e.g., on the login node for testing):
You will notice that the following warning is printed multiple times:
The warning is likely caused by underlying one-sided MPI routines forcing themselves on the InfiniBand interface library UCX, and nothing to worry about. In the example below we filter the warning to avoid cluttered output.Alternatively, if you intend to stay within the single node only, disabling the UCX functionality altogether might be an option:
Hello World example¶
The following code prints a message from every process ("image") on the first image:
caf_hello.f08
program hello
implicit none
integer :: i
integer :: rank_to_print[*]
! say hello
do i=1,num_images()
if (i==this_image()) then
rank_to_print[1]=this_image()
end if
sync all
if (this_image()==1) then
write(*,*) 'Hello from image',rank_to_print,'out of',num_images()
end if
sync all
end do
sync all
! say goodbye
do i=1,num_images()
if (i==this_image()) then
rank_to_print[1]=this_image()
end if
sync all
if (this_image()==1) then
write(*,*) 'Goodbye from image',rank_to_print,'out of',num_images()
end if
sync all
end do
end program hello
Example submission script that compiles and runs the program:
caf_hello.slurm
#!/bin/bash
#SBATCH --ntasks=5
#SBATCH --cpus-per-task=1
#SBATCH --account=innovation
#SBATCH --partition=compute
#SBATCH --time=00:01:00
#SBATCH --mem-per-cpu=1GB
#SBATCH --output=caf_hello.out
# Load necessary modules:
module load 2024r1
module load openmpi
module load opencoarrays
```
# Compile your code:
gfortran -fcoarray=lib caf_hello.f08 -o caf_hello.x -L${OPENCOARRAYS_ROOT}/lib64 -lcaf_mpi
```
# Run the code on 5 CPU cores:
srun ./caf_hello.x |grep -v UCX