Gromacs和CP2K交叉编译安装
环境:Rocky Linux 8.10, Slurm作业管理系统
安装编译器和相应库¶
根据系统版本,这里我们使用Gromacs 2024.3版本和CP2K 2024.3版本,并从相应官网地址下载对应版本源码压缩包。
由于这2个软件的编译特性,特别是交叉编译,先编译安装好下面软件。
GCC编译器安装¶
集群安装了GCC 14.1.0版本,如果没有或想安装其他版本,可根据编译器调用和安装自行安装。
集群gcc 14.1.0调用: module load gcc14.1.0
FFTW库安装¶
集群安装了单精度和双精度的fftw 3.3.10库,路径在/opt/soft/fftw-3.3.10。
当然也可自行安装,下载http://www.fftw.org/fftw-3.3.10.tar.gz
- 单精度fftw安装
module load gcc14.1.0
mkdir -p /opt/soft/fftw-3.3.10
tar -zxvf fftw-3.3.10.tar.gz
cd fftw-3.3.10
./configure --prefix=/opt/soft/fftw-3.3.10 \
--libdir=/opt/soft/fftw-3.3.10/lib \
--enable-mpi --enable-openmp --disable-shared --enable-static \
--enable-avx --enable-avx2 --enable-avx512 --enable-float
make -j
make install
make clean
./configure --prefix=/opt/soft/fftw-3.3.10 \
--libdir=/opt/soft/fftw-3.3.10/lib \
--enable-mpi --enable-openmp --disable-shared --enable-static \
--enable-avx --enable-avx2 --enable-avx512
make -j
make install
这样,单双精度的库都会安装在/opt/soft/fftw-3.3.10目录下。当然,个人用户是不能安装在/opt目录下的。
假如是userA用户,那就安装在/home/userA/soft/fftw-3.3.10目录下。
安装后的库目录结构(完整):
libfftw3.a (双精度DP)
libfftw3_mpi.a
libfftw3_omp.a
libfftw3f.a (单精度SP)
libfftw3f_mpi.a
libfftw3f_omp.a
OpenMPI安装¶
CP2K默认推荐使用mpich,单独使用mpich版的cp2k没问题,但是和gromacs交叉编译也用mpich,最终编译的二进制版本运行时,会提示内存泄露问题,因而这里我们用openmpi。
在OpenMPI官网下载openmpi-5.0.5.tar.bz2。编译安装方法,见Gromacs编译安装中的“关于编译mpi版本的Gromacs”的openmpi安装。
集群上默认安装了openmpi-5.0.5,安装路径为/opt/mpi/openmpi-5.0.5
调用方法:source /opt/mpi/openmpi-5.0.5/bin/mpivars.sh
CP2K安装¶
编译安装,参看CP2K编译安装。
使用toolchain编译,生成libcp2k,用于Gromacs交叉编译,不需要像以前使用cp2k所有模块编译安装。
module load gcc14.1.0
source /opt/mpi/openmpi-5.0.5/bin/mpivars.sh
#---假如cp2k安装的目录是/home/userA/soft/cp2k-2024.3---#
cd /home/userA/soft/cp2k-2024.3/tools/toolchain
./install_cp2k_toolchain.sh -j 36 --target-cpu=native --with-gcc=system --with-openmpi=system --with-cmake=install --with-fftw=/opt/soft/fftw-3.3.10 --with-quip=no --with-pexsi=no --with-elpa=no --with-sirius=no --with-plumed=install --with-openblas=install --with-scalapack=install
cp /home/userA/soft/cp2k-2024.3/tools/toolchain/install/arch/* /home/userA/soft/cp2k-2024.3/arch/
source /home/userA/soft/cp2k-2024.3/tools/toolchain/install/setup
cd /home/userA/soft/cp2k-2024.3
make -j 36 ARCH=local VERSION=psmp |tee make.log
#--下面这一步是测试例子是否都通过--#
make -j 36 ARCH=local VERSION=psmp test |tee make.psmp.test.log
#--下面这一步是生成libcp2k库--#
make -j 36 ARCH=local VERSION=psmp libcp2k
这个文件很重要,文件里从"-Wl,--enable-new-dtags"开始到最后-lstdc++结束,就是gromacs编译所需的-DCP2K_LINKER_FLAGS。
Gromacs编译¶
编译安装,参看Gromacs编译安装。
module load gcc14.1.0
source /opt/mpi/openmpi-5.0.5/bin/mpivars.sh
source /home/userA/soft/cp2k-2024.3/tools/toolchain/install/setup
tar -jxvf gromacs-2024.3.tar.gz
cd gromacs-2024.3
mkdir build && cd build
#--下面DCMAKE_INSTALL_PREFIX是安装的目录,个人用户不能安装在/opt目录下可以安装在家目录下--#
cmake .. -DCMAKE_PREFIX_PATH=/opt/soft/fftw-3.3.10 -DCMAKE_INSTALL_PREFIX=/opt/soft/gmx_cp2k_2024.3 \
-DBUILD_SHARED_LIBS=OFF -DGMXAPI=OFF -DGMX_INSTALL_NBLIB_API=OFF -DGMX_MPI=ON -DGMX_CP2K=ON \
-DCP2K_DIR=/home/userA/soft/cp2k-2024.3/lib/local/psmp \
-DGMX_BLAS_USER=/home/userA/soft/cp2k-2024.3/tools/toolchain/install/openblas-0.3.27/lib/libopenblas.a \
-DGMX_LAPACK_USER=/home/userA/soft/cp2k-2024.3/tools/toolchain/install/scalapack-2.2.1/lib/libscalapack.a \
-DCP2K_LINKER_FLAGS="-Wl,--enable-new-dtags -L/opt/mpi/openmpi-5.0.5/lib -Wl,-rpath -Wl,/opt/mpi/openmpi-5.0.5/lib -Wl,--enable-new-dtags -L'/home/userA/soft/cp2k-2024.3/tools/toolchain/install/openblas-0.3.27/lib' -Wl,-rpath,'/home/userA/soft/cp2k-2024.3/tools/toolchain/install/openblas-0.3.27/lib' -L'/opt/soft/fftw-3.3.10/lib' -Wl,-rpath,'/opt/soft/fftw-3.3.10/lib' -L'/home/userA/soft/cp2k-2024.3/tools/toolchain/install/libint-v2.6.0-cp2k-lmax-5/lib' -L'/home/userA/soft/cp2k-2024.3/tools/toolchain/install/libxc-6.2.2/lib' -Wl,-rpath,'/home/userA/soft/cp2k-2024.3/tools/toolchain/install/libxc-6.2.2/lib' -L'/home/userA/soft/cp2k-2024.3/tools/toolchain/install/libgrpp-main-20231225/lib' -Wl,-rpath,'/home/userA/soft/cp2k-2024.3/tools/toolchain/install/libgrpp-main-20231225/lib' -L'/home/userA/soft/cp2k-2024.3/tools/toolchain/install/libxsmm-1.17/lib' -Wl,-rpath,'/home/userA/soft/cp2k-2024.3/tools/toolchain/install/libxsmm-1.17/lib' -L'/home/userA/soft/cp2k-2024.3/tools/toolchain/install/scalapack-2.2.1/lib' -Wl,-rpath,'/home/userA/soft/cp2k-2024.3/tools/toolchain/install/scalapack-2.2.1/lib' -L'/home/userA/soft/cp2k-2024.3/tools/toolchain/install/COSMA-2.6.6/lib' -Wl,-rpath,'/home/userA/soft/cp2k-2024.3/tools/toolchain/install/COSMA-2.6.6/lib' -L'/home/userA/soft/cp2k-2024.3/tools/toolchain/install/gsl-2.7/lib' -Wl,-rpath,'/home/userA/soft/cp2k-2024.3/tools/toolchain/install/gsl-2.7/lib' -L'/home/userA/soft/cp2k-2024.3/tools/toolchain/install/plumed-2.9.0/lib' -Wl,-rpath,'/home/userA/soft/cp2k-2024.3/tools/toolchain/install/plumed-2.9.0/lib' -L'/home/userA/soft/cp2k-2024.3/tools/toolchain/install/spglib-2.3.1/lib' -Wl,-rpath,'/home/userA/soft/cp2k-2024.3/tools/toolchain/install/spglib-2.3.1/lib' -L'/home/userA/soft/cp2k-2024.3/tools/toolchain/install/libvori-220621/lib' -Wl,-rpath,'/home/userA/soft/cp2k-2024.3/tools/toolchain/install/libvori-220621/lib' -L/home/userA/soft/cp2k-2024.3/lib/local/psmp/exts/dbcsr -ldbcsr -lsymspg -l:libplumed.a -ldl -lstdc++ -lz -ldl -lgsl -lcosma_prefixed_pxgemm -lcosma -lcosta -lscalapack -lxsmmf -lxsmm -ldl -lpthread -llibgrpp -lxcf03 -lxc -lint2 -lfftw3_mpi -lfftw3 -lfftw3_omp -lmpi -l:libopenblas.a -lvori -lstdc++ -lstdc++ "
make -j 36 |tee make.log
make install
运行及作业提交¶
这里用简单的N-甲基甲酰胺作为例子,使用Gromacs-CP2K接口来做QM/MM模拟。
# 拷贝软件目录examples/nma
cp -r /opt/soft/gmx_cp2k_2024.3/examples/nma ./
cd nma
# 设置刷新软件环境变量
source /opt/soft/gmx_cp2k_2024.3/env.sh
# Make topology for the system using the following command
gmx_mpi pdb2gmx -f nma.pdb # 根据提示,第1步选择1,第2步选择1
# Perform energy minimization first for that molecule using QMMM interface
gmx_mpi grompp -f em.mdp -p topol.top -c conf.gro -o nma-em.tpr
#!/bin/bash
# Job submission script for SLURM:
# Usage: sbatch <this_script>
##############脚本第一部分,声明SLURM作业管理系统的请求资源参数###################
#SBATCH -J nma-em
#SBATCH -o nma-em.o%J
#SBATCH -N 1
#SBATCH --ntasks-per-node=16
#SBATCH -t 2024:03:00
#SBATCH -p public
########################请求资源参数的解释######################################
# -J: 作业名称
# -o: 标准输出stdout
# -e: 标准错误stderr,不写-e,标准输出stdout和标准错误stderr合并到.ojobid文件中
# -N: 请求的节点数目
# -n: 请求的总核数,如果-N后为1,等同于--ntasks-per-node
# -t: 时间限制
# -p: 申请的分区
# -w: 后面跟请求的指定节点名称
##############################################################################
#################脚本第二部分,定义要运行作业软件的环境变量#######################
# go to work dir
cd $SLURM_SUBMIT_DIR
# 设置环境变量
source /opt/soft/gmx_cp2k_2024.3/env.sh
# 设置CP2K的基准数据文件(基组/势函数)目录
export CP2K_DATA_DIR=/opt/soft/gmx_cp2k_2024.3/data
########################脚本第三部分,运行作业的命令#############################
echo "Starting Gromacs-CP2K run at" `date`
# The program we want to execute,根据自己的作业,主要修改的也是下面这一行
time mpirun -np $SLURM_NTASKS gmx_mpi mdrun -s nma-em.tpr
echo "Finished Gromacs-CP2K run at" `date`
sbatch run-em.slurm提交作业。
下一步将在300K NVT系综下执行QM/MM动力学模拟。先准备md-nvt.mdp产生nma-nvt.tpr。
gmx_mpi grompp -f md-nvt.mdp -p topol.top -c conf.gro -o nma-nvt.tpr
执行NVT QM/MM模拟,先准备run-nvt.slurm提交脚本
#!/bin/bash
# Job submission script for SLURM:
# Usage: sbatch <this_script>
##############脚本第一部分,声明SLURM作业管理系统的请求资源参数###################
#SBATCH -J nma-nvt
#SBATCH -o nma-nvt.o%J
#SBATCH -N 1
#SBATCH --ntasks-per-node=16
#SBATCH -t 2024:03:00
#SBATCH -p public
########################请求资源参数的解释######################################
# -J: 作业名称
# -o: 标准输出stdout
# -e: 标准错误stderr,不写-e,标准输出stdout和标准错误stderr合并到.ojobid文件中
# -N: 请求的节点数目
# -n: 请求的总核数,如果-N后为1,等同于--ntasks-per-node
# -t: 时间限制
# -p: 申请的分区
# -w: 后面跟请求的指定节点名称
##############################################################################
#################脚本第二部分,定义要运行作业软件的环境变量#######################
# go to work dir
cd $SLURM_SUBMIT_DIR
# 设置环境变量
source /opt/soft/gmx_cp2k_2024.3/env.sh
# 设置CP2K的基准数据文件(基组/势函数)目录
export CP2K_DATA_DIR=/opt/soft/gmx_cp2k_2024.3/data
########################脚本第三部分,运行作业的命令#############################
echo "Starting Gromacs-CP2K run at" `date`
# The program we want to execute,根据自己的作业,主要修改的也是下面这一行
time mpirun -np $SLURM_NTASKS gmx_mpi mdrun -s nma-nvt.tpr
echo "Finished Gromacs-CP2K run at" `date`
sbatch run-nvt.slurm提交作业。
本文阅读量 次
本站总访问量 次