snek5000.resources#

resources for building#

Contains common snakemake rules and jinja templates.

Snakemake rules#

MAKETOOLS = f"CC={config['CC']} FC={config['FC']} CFLAGS='{config['CFLAGS']}' FFLAGS='{config['FFLAGS']}' ./maketools"


rule tools:
    params:
        maketools=MAKETOOLS,
    shell:
        """
        pushd tools/
        {params.maketools} all
        """


rule _tool:
    input:
        "tools/{tool}",
        "nek5000_make_config.yml",
    output:
        "bin/{tool}",
    params:
        maketools=MAKETOOLS,
    shell:
        """
        pushd tools/
        {params.maketools} {wildcards.tool}
        """


rule tools_clean:
    shell:
        """
        rm -rf 3rd_party/gslib/lib 3rd_party/gslib/gslib
        pushd tools/
        ./maketools clean
        """


# NOTE: May not compile as needed unless proper flags are supplied. The command
# makenek takes care of gslib building using the bash function make_3rd_party
rule build_third_party:
    input:
        "nek5000_make_config.yml",
        nekconfig="bin/nekconfig",
    output:
        "3rd_party/gslib/lib/libgs.a",
        "3rd_party/blasLapack/libblasLapack.a",
    shell:
        """
        export CC="{config[MPICC]}" FC="{config[MPIFC]}" \
            CFLAGS="{config[CFLAGS]}" FFLAGS="{config[FFLAGS]}"
        {input.nekconfig} -build-dep
        """
  • compiler.smk: module for compiling user code

from pprint import pprint

import snek5000
from snek5000.clusters import nproc_available
from snek5000.util import now


NEK_SOURCE_ROOT = snek5000.get_nek_source_root()


# Snakemake configuration
rule show_config:
    run:
        pprint(config)


# compile
rule compile:
    input:
        "SIZE",
        f"{config['CASE']}.usr",
        "makefile_usr.inc",
        "makefile",
        NEK_SOURCE_ROOT + "/3rd_party/gslib/lib/libgs.a",
        NEK_SOURCE_ROOT + "/3rd_party/blasLapack/libblasLapack.a",
    params:
        make="make",
    output:
        exe="nek5000",
    shell:
        """
        {params.make} -j {output.exe} | tee -a build.log
        """


# mpiexec
rule mpiexec:
    input:
        f"{config['CASE']}.re2",
        f"{config['CASE']}.ma2",
        f"{config['CASE']}.par",
        "SESSION.NAME",
        "nek5000",
    output:
        "logs/run_" + now() + ".log",
    resources:
        nproc=nproc_available(),
    params:
        redirect=">",
        end="",
    shell:
        """
        ln -sf {output} {config[CASE]}.log
        echo "Log file:"
        realpath {config[CASE]}.log
        {config[MPIEXEC]} -n {resources.nproc} {config[MPIEXEC_FLAGS]} ./nek5000 {params.redirect} {output} {params.end}
        echo $PWD
        """


# run in background
use rule mpiexec as run with:
    params:
        redirect=">",
        end="&",


# run in foreground
use rule mpiexec as run_fg with:
    params:
        redirect="| tee",
        end="",
  • io.smk: module for cleaning and archiving simulation artefacts

import itertools
from glob import iglob

import snek5000
from snek5000.util.archive import tar_name, archive, clean_simul


NEK_SOURCE_ROOT = snek5000.get_nek_source_root()


# clean compiler output
rule clean:
    shell:
        """
        echo "cleaning Nek5000 ..."
        rm -fv {config[CASE]}.f nek5000
        rm -rf obj
        rm -fv $NEK_SOURCE_ROOT/core/mpif.h
        """


# clean simulation files
rule cleansimul:
    params:
        tarball=tar_name(compress_format=".zst"),
    run:
        clean_simul(config["CASE"], params.tarball)


# clean third party
rule clean3rd:
    params:
        gslib_dir=NEK_SOURCE_ROOT + "/3rd_party/gslib",
    shell:
        """
        nekconfig clean < <(yes)
        rm -rf {params.gslib_dir}/{{gslib-*,*.tar.gz}}
        """


# clean everything
rule cleanall:
    shell:
        """
        rm -rf obj
        snakemake clean clean3rd cleansimul -j
        """


# create an archive with all the results
rule archive:
    params:
        solution=sorted(
            itertools.chain.from_iterable(
                (
                    iglob(f"{prefix}{config['CASE']}0.f*")
                    for prefix in ("", "c2D", "sts")
                )
            )
        ),
        rest=[
            "SESSION.NAME",
            "params_simul.xml",
            "SIZE",
            *iglob(f"rs6{config['CASE']}0.f*"),
            f"{config['CASE']}.re2",
            f"{config['CASE']}.ma2",
            f"{config['CASE']}.par",
            f"{config['CASE']}.usr",
        ],
        tarball=tar_name(),
    run:
        archive(params.tarball, params.solution, remove=True)
        archive(params.tarball, params.rest)
        archive(params.tarball + ".zst", readonly=True)
  • internal.smk: module for generating simulation files, serves as inputs to Snakemake rules above

import snek5000
from snek5000.util.files import create_session


NEK_SOURCE_ROOT = snek5000.get_nek_source_root()


# generate a box mesh
rule generate_box:
    input:
        box=f"{config['CASE']}.box",
        genbox=NEK_SOURCE_ROOT + "/bin/genbox",
    output:
        "box.re2",
    shell:
        "{input.genbox} <<< {input.box}"


# rename mesh file re2
rule move_box:
    input:
        "box.re2",
    output:
        f"{config['CASE']}.re2",
    shell:
        "mv -f {input} {output}"


# generate map / connectivity matrix
rule generate_map:
    input:
        f"{config['CASE']}.re2",
        genmap=NEK_SOURCE_ROOT + "/bin/genmap",
    output:
        f"{config['CASE']}.ma2",
    resources:
        tolerance="0.01",
    shell:
        r'printf "{config[CASE]}\n{resources.tolerance}\n" '
        "  | {input.genmap} "
        "  | (head -n 20; tail -n 10) "


# generate makefile
rule generate_makefile:
    input:
        f"{config['CASE']}.re2",
        f"{config['CASE']}.ma2",
        f"{config['CASE']}.usr",
        cmd=NEK_SOURCE_ROOT + "/bin/nekconfig",
    output:
        "makefile",
        ".state",
    shell:
        """
        set +u  # Do not error on undefined bash variables
        export CC="{config[MPICC]}" FC="{config[MPIFC]}" \
            CFLAGS="{config[CFLAGS]}" FFLAGS="{config[FFLAGS]} {config[includes]}" \
            USR="{config[objects]}"
        {input.cmd} {config[CASE]}
        """


# generate sessionfile
rule generate_session:
    input:
        re2=f"{config['CASE']}.re2",
        ma2=f"{config['CASE']}.ma2",
        par=f"{config['CASE']}.par",
    output:
        "SESSION.NAME",
    run:
        create_session(config["CASE"], input.re2, input.ma2, input.par)

Jinja templates#

{# A template for generating box files.

A box file is a description of a simple box mesh. It describes the spatial
dimension of the mesh, the number of fields, the domain extents and boundary
conditions. The box file is rendered using the write_box method of the Output
class.
-#}
{% set indent = 42 -%}
{{ dim.ljust(indent) }} spatial dimension (dim<0 will create box.re2)
{{ nb_fields.ljust(indent) }} number of fields
#========================================================
#
{% for comment in comments.splitlines() -%}
#    {{ comment }}
{% endfor -%}
#
#========================================================
#
Box
{% for comment, grid in grid_info.items() -%}
{{ grid.ljust(indent) }}{{comment}}
{% endfor -%}
{# A template for generating user makefiles.

By using a user makefile, extra Fortran files can be included in the main
makefile. This file takes care of generating object files (*.o) from a list of
fortran files. The user makefile is rendered using the write_makefile_usr
method of the Output class.
#}
# vim: set ft=make
{% for comment in comments.splitlines() -%}
# {{ comment }}
{% endfor -%}
{% macro make(sources) -%}
{# A cosmetic macro to indent comments nicely #}
{%- set target = sources[0].replace(".f", ".o").split("/")[-1] -%}
{%- set rule_string = "$(OBJDIR)/{}:{}\n\t$(F77) -c $(FL2) -I./ $< -o $@".format(
		target.ljust(20),
		" ".join(sources),
	)
-%}
{{ rule_string }}
{% endmacro -%}

{% for sources in list_of_sources -%}
{{ make(sources) }}
{% endfor -%}
{# A template for generating SIZE files.

A SIZE file contains a collection of Fortran parameters (fortran_var) and their
values (python_var). The values are initialized as parameters in Python. The
template is rendered by the write_size method of the Output class.
#}
{% macro parameter(fortran_var, python_var) -%}
{# A cosmetic macro to indent comments nicely #}
{%- set parameter_string = "parameter ({}={})".format(fortran_var, python_var) -%}
{{ parameter_string.ljust(32) -}}
{% endmacro -%}
c
c     Include file to dimension static arrays
c     and to set some hardwired run-time parameters
c
{% for comment in comments.splitlines() -%}
c     {{ comment }}
{% endfor -%}
c
c
      integer ldim,lx1,lxd,lx2,lx1m,lelg,lelt,lpmin,lpmax,ldimt
      integer lpelt,lbelt,toteq,lcvelt
      integer lelx,lely,lelz,mxprev,lgmres,lorder,lhis
      integer maxobj,lpert,nsessmax,lxo
      integer lfdm,ldimt_proj

      ! BASIC
      {{ parameter("ldim", params.oper.dim) }}  ! domain dimension (2 or 3)
      {{ parameter("lx1", params.oper.elem.order) }}  ! GLL points per element along each direction
      {{ parameter("lxd", order_dealiasing) }}  ! GL  points for over-integration (dealiasing)
      {{ parameter("lx2", order_pressure) }}  ! GLL points for pressure (lx1 or lx1-2)

      {{ parameter("lelg", max_n_seq) }}  ! max total number of elements
      {{ parameter("lpmin", params.oper.nproc_min) }}  ! min MPI ranks
      {{ parameter("lpmax", params.oper.nproc_max) }}  ! max MPI ranks
      {{ parameter("ldimt", (1, params.oper.scalars) | max) }}  ! max auxiliary fields (temperature + scalars)

      ! OPTIONAL
      {{ parameter("ldimt_proj", params.oper.max.scalars_proj) }}  ! max auxiliary fields residual projection
      {{ parameter("lhis", params.oper.max.hist) }}  ! max history/monitoring points
      {{ parameter("maxobj", params.oper.max.obj) }}  ! max number of objects
      {{ parameter("lpert", params.oper.max.perturb) }}  ! max number of perturbations
      {{ parameter("toteq", params.oper.max.scalars_cons) }}  ! max number of conserved scalars in CMT
      {{ parameter("nsessmax", params.oper.max.sessions) }}  ! max sessions to NEKNEK
      {{ parameter("lxo", order_out) }}  ! max GLL points on output (lxo>=lx1)
      {{ parameter("mxprev", params.oper.max.dim_proj) | trim }}
      {{ parameter("lgmres", params.oper.max.dim_krylov) }}  ! max dim of projection & Krylov space
      {{ parameter("lorder", max_order_time) }}  ! max order in time
      {{ parameter("lx1m", order_mesh_solver) }}  ! GLL points mesh solver
      {{ parameter("lfdm", params.oper.misc.fast_diag | int) }}  ! set to 1 for fast diagonalization method
      {{ parameter("lelx", max_nx) }}  ! global tensor mesh dimensions
      {{ parameter("lely", max_ny) | trim }}
{%- if params.oper.dim == 3 %}
      {{ parameter("lelz", max_nz) | trim }}
{%- endif %}
      {{ parameter("lelt", max_n_loc) }}  ! max number of local elements per MPI rank
      {{ parameter("lbelt", order_mhd) }}  ! set to lelt for mhd
      {{ parameter("lpelt", order_linear) }}  ! set to lelt for linear stability
      {{ parameter("lcvelt", order_cvode) }}  ! set to lelt for cvode
{% block user_size %}{% endblock %}
      ! INTERNALS
      include 'SIZE.inc'
c vim: set ft=fortran
#!/bin/bash

function error_quit {
    echo -e "$@"
    echo
    echo -e 'Usage:'
    echo -e './compile.sh --clean'
    echo -e '   To clean build directory. Makenek will ask for cleaning 3rd-party libraries.'
    echo
    echo -e './compile.sh --all'
    echo -e '   To compile the code.'
    exit 1
}

#parameters
export CASE="{{ CASE }}"
export NEK_SOURCE_ROOT=${NEK_SOURCE_ROOT:"../../Nek5000"}
export FC="{{ MPIFC }}"
export CC="{{ MPICC }}"
export CFLAGS="{{ CFLAGS }}"
export FFLAGS="{{ FFLAGS }} {{ INC }}"
export PPLIST=""
export USR="{{ USR }}"
export USR_LFLAGS=""

# arguments
args=("$@")
argsnr=$#

# check arguments
# parameters number check
if [ $[argsnr] -ne 1 ]; then
    error_quit 'Wrong arguments number!'
fi

for il in "$@"
do
case $il in
    --clean)
        ${NEK_SOURCE_ROOT}/bin/makenek clean
        shift
        ;;
    --all)
        ${NEK_SOURCE_ROOT}/bin/makenek ${CASE}
        shift
        ;;
    *)
        error_quit 'Wrong argument.'
        ;;
esac
done

Functions

get_base_template(name)

Get a template from snek5000.resources.

get_base_templates()

Get templates (box, makefile_usr and size) from snek5000.resources.

Classes

class snek5000.resources.BaseTemplates[source]#
property env#
get_base_templates()[source]#

Get templates (box, makefile_usr and size) from snek5000.resources.

get_base_template(name)[source]#

Get a template from snek5000.resources.

snek5000.resources.get_base_templates()#

Get templates (box, makefile_usr and size) from snek5000.resources.

snek5000.resources.get_base_template(name)#

Get a template from snek5000.resources.