粒子の座標ごとに値を持つ計算結果 (複数グループ可)

粒子の座標ごとに値を持つ計算結果を出力する場合、 表 75 に示す関数を 使用します。

ここで示す関数を使うと、複数のグループの粒子を出力することができます。 各グループの出力の最初と最後で、 cg_iric_write_sol_particlegroup_groupbegin_fcg_iric_write_sol_particlegroup_groupend_f を呼び出してください。

出力するプログラムの例は、リスト 107 を参照してください。

注釈

粒子の座標ごとに値を持つ計算結果 で示した関数とは異なり、 ここで示す関数では、一度の関数呼び出しでは粒子一つ分のデータを出力します。

表 75 粒子ごとに値を持つ計算結果の出力に利用する関数 (複数グループ可)

関数

備考

cg_iric_write_sol_particlegroup_groupbegin_f

粒子の計算結果の出力を開始する

cg_iric_write_sol_particlegroup_groupend_f

粒子の計算結果の出力を終了する

cg_iric_write_sol_particlegroup_pos2d_f

粒子の位置を出力する (2次元)

cg_iric_write_sol_particlegroup_pos3d_f

粒子の位置を出力する (3次元)

cg_iric_write_sol_particlegroup_integer_f

整数の粒子ごとに値を持つ計算結果を出力する

cg_iric_write_sol_particlegroup_real_f

倍精度実数の粒子ごとに値を持つ計算結果を出力する

リスト 107 サンプルプログラム (粒子の座標ごとに値を持つ計算結果 (複数グループ可))
 1program SampleProgram
 2  implicit none
 3  include 'cgnslib_f.h'
 4
 5  integer:: fin, ier, isize, jsize
 6  integer:: canceled
 7  integer:: locked
 8  double precision:: time
 9  double precision, dimension(:,:), allocatable::grid_x, grid_y
10  integer:: numparticles = 10
11  double precision, dimension(:), allocatable:: particle_x, particle_y,
12  double precision, dimension(:), allocatable:: velocity_x, velocity_y, temperature
13  integer:: i
14  integer:: status = 1
15  character(len=20):: condFile
16
17  condFile = 'test.cgn'
18
19  ! CGNS ファイルのオープン
20  call cg_open_f(condFile, CG_MODE_MODIFY, fin, ier)
21  if (ier /=0) STOP "*** Open error of CGNS file ***"
22
23  ! 内部変数の初期化
24  call cg_iric_init_f(fin, ier)
25  if (ier /=0) STOP "*** Initialize error of CGNS file ***"
26
27  ! 格子のサイズを調べる
28  call cg_iric_gotogridcoord2d_f(isize, jsize, ier)
29  ! 格子を読み込むためのメモリを確保
30  allocate(grid_x(isize, jsize), grid_y(isize, jsize))
31  ! 計算結果を保持するメモリを確保。
32  allocate(particle_x(numparticles), particle_y(numparticles))
33  allocate(velocity_x(numparticles), velocity_y(numparticles), temperature(numparticles))
34
35  ! 格子を読み込む
36  call cg_iric_getgridcoord2d_f(grid_x, grid_y, ier)
37
38  ! 初期状態の情報を出力
39  time = 0
40  call cg_iric_write_sol_time_f(time, ier)
41
42  call cg_iric_write_sol_particlegroup_groupbegin_f('driftwood', ier)
43  do i = 1, numparticles
44    ! (ここで particle_x, particle_x, velocity_x, velocity_y, temperature に適切な値を設定)
45    call cg_iric_write_sol_particlegroup_pos2d_f(particle_x(i), particle_y(i), ier)
46    call cg_iric_write_sol_particlegroup_real_f('VelocityX', velocity_x(i), ier)
47    call cg_iric_write_sol_particlegroup_real_f('VelocityY', velocity_y(i), ier)
48    call cg_iric_write_sol_particlegroup_real_f('Temperature', temperature(i), ier)
49  end do
50  call cg_iric_write_sol_particlegroup_groupend_f(ier)
51
52  do
53    time = time + 10.0
54
55    ! (ここで計算を実行)
56
57    call iric_check_cancel_f(canceled)
58    if (canceled == 1) exit
59
60    ! 計算結果を出力
61    call iric_write_sol_start_f(condFile, ier)
62    call cg_iric_write_sol_time_f(time, ier)
63    call cg_iric_write_sol_particlegroup_groupbegin_f('driftwood', ier)
64    do i = 1, numparticles
65      ! (ここで particle_x, particle_x, velocity_x, velocity_y, temperature に適切な値を設定)
66      call cg_iric_write_sol_particlegroup_pos2d_f(particle_x(i), particle_y(i), ier)
67      call cg_iric_write_sol_particlegroup_real_f('VelocityX', velocity_x(i), ier)
68      call cg_iric_write_sol_particlegroup_real_f('VelocityY', velocity_y(i), ier)
69      call cg_iric_write_sol_particlegroup_real_f('Temperature', temperature(i), ier)
70    end do
71    call cg_iric_write_sol_particlegroup_groupend_f(ier)
72
73    if (time > 1000) exit
74  end do
75
76  ! CGNS ファイルのクローズ
77  call cg_close_f(fin, ier)
78  stop
79end program SampleProgram