GROMACSは、分子動力学シミュレーションのソフトウェアパッケージです。OCTOPUSでは、バッチリクエストによる GROMACS 2021.2 の利用が可能となっております。GROMACSの利用は申請不要です。SQUIDの利用を申し込まれた方は、どなたでも利用可能です。
基本的な利用方法
-
GROMACSの実行は、バッチリクエストによる処理のみ許可しております。
フロントエンドノードに接続し、計算に必要な入力ファイル、ジョブスクリプトを作成後、ジョブ投入して行います。
GROMACSを利用する際のジョブスクリプト例と、ジョブ実行方法について以下で解説いたします。
※GROMACSについてのマニュアルは公式HPをご参照ください。
ジョブスクリプトの作成:汎用CPUノード群向け
-
以下の例は152並列(2ノード使用、ノード当たり76並列)でGROMACSを実行する場合のジョブスクリプト例です。ファイル名に特に指定はありませんが、本項ではgromacs.shとしています。
|
1 2 3 4 5 6 7 8 9 10 11 12 |
#!/bin/bash #PBS -q SQUID #PBS -l cpunum_job=76 #PBS --group=[グループ名] #PBS -b 2 #PBS -T intmpi #PBS -l elapstim_req=01:00:00 #PBS -v OMP_NUM_THREADS=1 module load BaseApp module load gromacs/2021.2 cd $PBS_O_WORKDIR mpirun ${NQSII_MPIOPTS} -np 152 -ppn 76 gmx_mpi_d mdrun (Input file) |
OpenMPを使用するプログラムは、環境変数 OMP_NUM_THREADS に基づいてスレッド数を決定しますが、指定しない場合、過剰なCPUコアが割り当てられてしまいエラーが発生してしまう可能性があります。
そのため、ジョブスクリプト内でスレッド数を明示的に設定することを推奨いたします。
ジョブスクリプトのその他の行についてはこちらをご参照ください。
ジョブスクリプトの作成:GPUノード群向け
以下の例はGPUを使ってGROMACSを実行する場合のジョブスクリプト例です (2ノード使用、計16GPU)。ファイル名に特に指定はありませんが、本項ではgromacs.shとしています。-npmeの値 (長距離相互作用計算を担当するプロセス数) は、全体のプロセス数や、長距離と短距離相互作用計算の演算量のバランスによって最適値が異なるため、実行時間を計測して調整してください。
123456789101112131415161718
#!/bin/bash#PBS -q SQUID#PBS --group=(グループ名)#PBS -l elapstim_req=00:10:00#PBS -b 2#PBS -T openmpi#PBS -l gpunum_job=8#PBS -v NQSV_MPI_MODULE=BaseApp/2025:gromacs/2025.4mpi.GPU#PBS -v OMP_NUM_THREADS=8#PBS -v GMX_GPU_PME_DECOMPOSITION=1 module load BaseAppmodule load gromacs/2025.4mpi.GPU cd $PBS_O_WORKDIR mpirun ${NQSV_MPIOPTS} -np 16 -npernode 8 gmx_mpi mdrun -s poly-ch2.tpr \ -npme 8 -nb gpu -pme gpu -bonded gpu -update gpu
ジョブスクリプトのその他の行についてはこちらをご参照ください。
実行方法
作成したジョブスクリプトを投入します。
% qsub gromacs.sh
実行が成功すると、結果ファイルに計算結果が出力されます。
|
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 |
#!/bin/bash #PBS -q SQUID #PBS --group=(グループ名) #PBS -l elapstim_req=00:10:00 #PBS -b 2 #PBS -T openmpi #PBS -l gpunum_job=8 #PBS -v NQSV_MPI_MODULE=BaseApp/2025:gromacs/2025.4mpi.GPU #PBS -v OMP_NUM_THREADS=8 #PBS -v GMX_GPU_PME_DECOMPOSITION=1 module load BaseApp module load gromacs/2025.4mpi.GPU cd $PBS_O_WORKDIR mpirun ${NQSV_MPIOPTS} -np 16 -npernode 8 gmx_mpi mdrun -s poly-ch2.tpr \ -npme 8 -nb gpu -pme gpu -bonded gpu -update gpu |
ジョブスクリプトのその他の行についてはこちらをご参照ください。
-
作成したジョブスクリプトを投入します。
% qsub gromacs.sh
実行が成功すると、結果ファイルに計算結果が出力されます。

