OpenMP 入門

ナビゲーション:前へ   上へ   次へ

12 nowait指示節

do指示構文、sections指示構文、single指示構文(及び本テキストでは説明しないworkshare指示構文)はその出口で同期が取られます。(暗黙のバリア) このため処理が先に終わったスレッドは他のすべてのスレッドの処理が終わるまで何もせずに待つことになります。
(参考1:暗黙のバリアはこれらの指示構文以外にparallel指示構文、parallel do指示構文、parallel sections指示構文の出口でも適用されます)
(参考2:暗黙のバリアはこれらの各指示構文の出口でのみ適用され、入口では適用されません)

nowait指示節を利用すると、処理が早く終わったスレッドが他のスレッドを待つことなく次の処理へ移行するように指示することができます。

例えば下記のコード例1では、最初のループ(1st Loop部分)の処理がすべてのスレッドで完了したのを待ってから、 2つ目のループ(2nd Loop部分)に処理が移行します。

[コード例1:nowaitを指定しない例]
program exampleNowait
  implicit none
  integer i, a(100), b(50)
!$omp parallel

  !----------------- 1st Loop -----------------
!$omp do
  do i=1, 100
    a(i) = i
  end do
!$omp end do    ! <---ここで同期が取られる
  !--------------- END 1st Loop ---------------

  !----------------- 2nd Loop -----------------
!$omp do
  do i=1, 50
    b(i) = i*2
  end do
!$omp end do
  !--------------- END 2nd Loop ---------------

!$omp end parallel
  print *, sum(a), sum(b)
end program

次のコード例2では、最初のループ(1st Loop部分)の処理において自分が担当する部分の処理を終えたスレッドは、 他のスレッドの完了を待つことなく2つ目のループ処理(2nd Loop部分)に移行します。

[コード例2:nowaitを指定する例]
program exampleNowait
  implicit none
  integer i, a(100), b(50)
!$omp parallel

  !----------------- 1st Loop -----------------
!$omp do
  do i=1, 100
    a(i) = i
  end do
!$omp end do nowait    ! <---同期は取られない
  !--------------- END 1st Loop ---------------

  !----------------- 2nd Loop -----------------
!$omp do
  do i=1, 50
    b(i) = i*2
  end do
!$omp end do
  !--------------- END 2nd Loop ---------------

!$omp end parallel
  print *, sum(a), sum(b)
end program

nowaitは何もしないで無駄に待つ状況を回避する意味で有用ですが、 一方でいとも簡単に意図しない不具合を招いてしまいかねないものです。 nowaitを利用しないコードでまず十分に動作の確認を行った上でnowaitを追加する等、 その利用は慎重に行う必要があります。

12.1 ☆演習課題:nowait指示節

下記に示すプログラムは配列aと配列bを初期化後にそれぞれの合計を計算して表示するプログラムです。 このコードをできるだけ効率よく実行するためにNOWAIT指示節を用いて書き換えて下さい。
[課題のコード]
program kadaiNowait
  use omp_lib
  implicit none
  integer i, a(10), b(20)
  print *, "Before Parallel Region"
!$omp parallel
  print *, "Before 1st Loop :", omp_get_thread_num()
!$omp do
  do i=1,10
    a(i) = i
  end do
!$omp end do
  print *, "Between 1st and 2nd Loops :", omp_get_thread_num()
!$omp do
  do i=1,20
    b(i) = 21-i
  end do
!$omp end do
  print *, "After 2nd Loop :", omp_get_thread_num()
!$omp end parallel
  print *, "After Parallel Region"
  print *, "Sums are :", sum(a), sum(b)
end program
[現状のまま4スレッドで実行した場合の実行例]
 Before Parallel Region
 Before 1st Loop : 1
 Before 1st Loop : 0
 Before 1st Loop : 2
 Before 1st Loop : 3
 Between 1st and 2nd Loops : 3
 Between 1st and 2nd Loops : 2
 Between 1st and 2nd Loops : 0
 Between 1st and 2nd Loops : 1
 After 2nd Loop : 1
 After 2nd Loop : 0
 After 2nd Loop : 2
 After 2nd Loop : 3
 After Parallel Region
 Sums are : 55 210

[書き換え後に4スレッドで実行した場合の実行例]
 Before Parallel Region
 Before 1st Loop : 1
 Before 1st Loop : 2
 Between 1st and 2nd Loops : 2
 After 2nd Loop : 2
 Between 1st and 2nd Loops : 1
 After 2nd Loop : 1
 Before 1st Loop : 3
 Between 1st and 2nd Loops : 3
 After 2nd Loop : 3
 Before 1st Loop : 0
 Between 1st and 2nd Loops : 0
 After 2nd Loop : 0
 After Parallel Region
 Sums are : 55 210
解答例:kadaiNowait.f90


ナビゲーション:前へ   上へ   次へ
関連情報
ご案内
© 日本ニューメリカルアルゴリズムズグループ株式会社 2025
Privacy Policy  /  Trademarks