Value defined at particles (groups supported)¶
When you output value defined at particles, please use the functions in Table 83.
Using the functions explained here, you can output multiple groups of particles. To output data please call cg_iric_write_sol_particlegroup_groupbegin and cg_iric_write_sol_particlegroup_groupend, before and after outputting data.
List 174 shows an example of the process to output value defined at particles.
Note
In functions explained here, data for one particle is output with each function call.
Subroutine |
Remarks |
---|---|
cg_iric_write_sol_particlegroup_groupbegin |
Start outputting calculation result defined at particles. |
cg_iric_write_sol_particlegroup_groupend |
Finish outputting calculation result defined at particles. |
cg_iric_write_sol_particlegroup_pos2d |
Outputs particle positions (two-dimensions) |
cg_iric_write_sol_particlegroup_pos3d |
Outputs particle positions (three-dimensions) |
cg_iric_write_sol_particlegroup_integer |
Outputs integer-type calculation results, giving a value for each particle. |
cg_iric_write_sol_particlegroup_real |
Outputs double-precision real-type calculation results, giving a value for each particle. |
1program SampleProgram
2 use iric
3 implicit none
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 ! Open CGNS file
20 call cg_iric_open(condFile, IRIC_MODE_MODIFY, fin, ier)
21 if (ier /=0) STOP "*** Open error of CGNS file ***"
22
23 ! Check the grid size
24 call cg_iric_read_grid2d_str_size(fin, isize, jsize, ier)
25 ! Allocate memory for loading the grid
26 allocate(grid_x(isize, jsize), grid_y(isize, jsize))
27 ! Allocate memory for calculation result.
28 allocate(particle_x(numparticles), particle_y(numparticles))
29 allocate(velocity_x(numparticles), velocity_y(numparticles), temperature(numparticles))
30
31 ! Read the grid into memory
32 call cg_iric_read_grid2d_coords(fin, grid_x, grid_y, ier)
33
34 ! Output the initial state information.
35 time = 0
36 call cg_iric_write_sol_start(fin, ier)
37 call cg_iric_write_sol_time(fin, time, ier)
38
39 call cg_iric_write_sol_particlegroup_groupbegin(fin, 'driftwood', ier)
40 do i = 1, numparticles
41 ! (Input values to particle_x, particle_x, velocity_x, velocity_y, temperature here)
42 call cg_iric_write_sol_particlegroup_pos2d(fin, particle_x(i), particle_y(i), ier)
43 call cg_iric_write_sol_particlegroup_real(fin, 'VelocityX', velocity_x(i), ier)
44 call cg_iric_write_sol_particlegroup_real(fin, 'VelocityY', velocity_y(i), ier)
45 call cg_iric_write_sol_particlegroup_real(fin, 'Temperature', temperature(i), ier)
46 end do
47 call cg_iric_write_sol_particlegroup_groupend(fin, ier)
48 call cg_iric_write_sol_end(fin, ier)
49
50 do
51 time = time + 10.0
52
53 ! (Perform calculation here)
54
55 call iric_check_cancel(canceled)
56 if (canceled == 1) exit
57
58 ! Output calculation results
59 call cg_iric_write_sol_start(fin, ier)
60 call cg_iric_write_sol_time(fin, time, ier)
61 call cg_iric_write_sol_particlegroup_groupbegin(fin, 'driftwood', ier)
62 do i = 1, numparticles
63 ! (Input values to particle_x, particle_x, velocity_x, velocity_y, temperature here)
64 call cg_iric_write_sol_particlegroup_pos2d(fin, particle_x(i), particle_y(i), ier)
65 call cg_iric_write_sol_particlegroup_real(fin, 'VelocityX', velocity_x(i), ier)
66 call cg_iric_write_sol_particlegroup_real(fin, 'VelocityY', velocity_y(i), ier)
67 call cg_iric_write_sol_particlegroup_real(fin, 'Temperature', temperature(i), ier)
68 end do
69 call cg_iric_write_sol_particlegroup_groupend(fin, ier)
70 call cg_iric_write_sol_end(fin, ier)
71
72 if (time > 1000) exit
73 end do
74
75 ! Close CGNS file
76 call cg_iric_close(fin, ier)
77 stop
78end program SampleProgram