MFIX  2016-1
write_res1.f
Go to the documentation of this file.
1 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvC
2 ! C
3 ! Module name: WRITE_RES1 C
4 ! Purpose: write out the time-dependent restart records C
5 ! C
6 ! Author: P. Nicoletti Date: 13-DEC-91 C
7 ! Reviewer: P. Nicoletti, W. Rogers, M. Syamlal Date: 24-JAN-92 C
8 ! C
9 ! Revision Number: C
10 ! Purpose: C
11 ! Author: Date: dd-mmm-yy C
12 ! Reviewer: Date: dd-mmm-yy C
13 ! C
14 ! Literature/Document References: C
15 ! C
16 ! Variables referenced: TIME, NSTEP, EP_g, P_g, P_star, RO_g, ROP_g C
17 ! T_g, T_s, U_g, V_g, W_g, ROP_s, U_s C
18 ! V_s, W_s, IJKMAX2 C
19 ! Variables modified: None C
20 ! C
21 ! Local variables: LC, N, NEXT_REC C
22 ! C
23 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^C
24 !
25  SUBROUTINE write_res1
26 !...Translated by Pacific-Sierra Research VAST-90 2.06G5 12:17:31 12/09/98
27 !...Switches: -xf
28 !
29 !-----------------------------------------------
30 ! M o d u l e s
31 !-----------------------------------------------
32  USE cdist
33  USE compar !//
34  USE energy
35  USE fldvar
36  USE funits
37  USE geometry
38  USE machine, only: flush_res
39  USE mpi_utility !//d pnicol : for gather
40  USE output
41  USE param
42  USE param1
43  USE physprop
44  USE run
45  USE rxns
46  USE scalars
47  USE sendrecv !//d pnicol : for gather
48 !//d pnicol ... not needed USE tmp_array
49  IMPLICIT NONE
50 !-----------------------------------------------
51 ! G l o b a l P a r a m e t e r s
52 !-----------------------------------------------
53 !-----------------------------------------------
54 ! L o c a l P a r a m e t e r s
55 !-----------------------------------------------
56 !-----------------------------------------------
57 ! L o c a l V a r i a b l e s
58 !-----------------------------------------------
59 !
60 !
61 !
62 !
63 !//d pnicol : allocate arrays for gather/convert_to_io_dp
64  double precision, allocatable :: array1(:)
65  double precision, allocatable :: array2(:)
66 
67 
68 ! loop counter
69  INTEGER :: LC, NN
70 !
71 ! pointer to first time-dependent record in restart file
72  INTEGER :: NEXT_REC
73 !-----------------------------------------------
74 !
75 
76 !//d pnicol
77 ! if (myPE.eq.PE_IO .and. .not.distio) then
78  allocate (array1(ijkmax2))
79  allocate (array2(ijkmax3))
80 ! else
81 ! allocate (array1(1))
82 ! allocate (array2(1))
83 ! end if
84 
85 
86  if (mype.eq.pe_io .or. bdist_io) then
87  READ (unit_res, rec=3) next_rec
88  WRITE (unit_res, rec=next_rec) time, dt, nstep
89  next_rec = next_rec + 1
90  end if
91 !
92 !\\SP Local Send Receive - need to be moved to source later!!
93 
94  if (.not. bdist_io) then
95 
96  call send_recv(ep_g,2)
97  call send_recv(p_g,2)
98  call send_recv(p_star,2)
99  call send_recv(ro_g,2)
100  call send_recv(rop_g,2)
101  call send_recv(x_g,2)
102  call send_recv(t_g,2)
103  call send_recv(u_g,2)
104  call send_recv(v_g,2)
105  call send_recv(w_g,2)
106  call send_recv(rop_s,2)
107  call send_recv(t_s,2)
108  call send_recv(u_s,2)
109  call send_recv(v_s,2)
110  call send_recv(w_s,2)
111  call send_recv(theta_m,2)
112  call send_recv(x_s,2)
113  if(nscalar > 0)call send_recv(scalar,2)
114  if(k_epsilon) THEN
115  call send_recv(k_turb_g,2)
116  call send_recv(e_turb_g,2)
117  endif
118  call send_recv(gama_rg,2)
119  call send_recv(t_rg,2)
120  call send_recv(gama_rs,2)
121  call send_recv(t_rs,2)
122  if(nrr > 0)call send_recv(reactionrates,2)
123 
124  end if
125 
126 
127  call gatherwriteres (ep_g,array2, array1, next_rec) !//d pnicol
128 !
129  call gatherwriteres (p_g,array2, array1, next_rec) !//d pnicol
130 !
131  call gatherwriteres (p_star,array2, array1, next_rec) !//d pnicol
132 !
133  call gatherwriteres (ro_g,array2, array1, next_rec) !//d pnicol
134 !
135  call gatherwriteres (rop_g,array2, array1, next_rec) !//d pnicol
136 !
137  call gatherwriteres (t_g,array2, array1, next_rec) !//d pnicol
138 !
139  DO nn = 1, nmax(0)
140  call gatherwriteres (x_g(:,nn),array2, array1, next_rec) !//d pnicol
141  END DO
142 !
143  call gatherwriteres (u_g,array2, array1, next_rec) !//d pnicol
144 !
145  call gatherwriteres (v_g,array2, array1, next_rec) !//d pnicol
146 !
147  call gatherwriteres (w_g,array2, array1, next_rec) !//d pnicol
148 !
149  DO lc = 1, mmax
150 !
151  call gatherwriteres (rop_s(:,lc),array2, array1, next_rec) !//d pnicol
152 
153  IF(any(solve_ros)) &
154  call gatherwriteres (ro_s(:,lc),array2, array1, next_rec) !//d pnicol
155 !
156  call gatherwriteres (t_s(:,lc),array2, array1, next_rec) !//d pnicol
157 !
158  call gatherwriteres (u_s(:,lc),array2, array1, next_rec) !//d pnicol
159 !
160  call gatherwriteres (v_s(:,lc),array2, array1, next_rec) !//d pnicol
161 !
162  call gatherwriteres (w_s(:,lc),array2, array1, next_rec) !//d pnicol
163 !
164  call gatherwriteres (theta_m(:,lc),array2, array1, next_rec) !//d pnicol
165 !
166  DO nn = 1, nmax(lc)
167  call gatherwriteres (x_s(:,lc,nn),array2, array1, next_rec) !//d pnicol
168  END DO
169  END DO
170 !
171 ! Version 1.3
172 
173  DO lc = 1, nscalar
174  call gatherwriteres (scalar(:,lc),array2, array1, next_rec) !//d pnicol
175  END DO
176 !
177 ! Version 1.4 -- write radiation variables in write_res1
178  call gatherwriteres (gama_rg,array2, array1, next_rec) !//d pnicol
179 
180  call gatherwriteres (t_rg,array2, array1, next_rec) !//d pnicol
181 
182  DO lc = 1, mmax
183  call gatherwriteres (gama_rs(1,lc),array2, array1, next_rec) !//d pnicol
184 
185  call gatherwriteres (t_rs(1,lc),array2, array1, next_rec) !//d pnicol
186  ENDDO
187 
188 !
189 ! Version 1.5
190  DO lc = 1, nrr
191  call gatherwriteres (reactionrates(:,lc),array2, array1, next_rec) !//d pnicol
192  END DO
193 !
194 ! Version 1.6
195 
196  if (k_epsilon) then
197  call gatherwriteres (k_turb_g,array2, array1, next_rec) !//d pnicol
198  call gatherwriteres (e_turb_g,array2, array1, next_rec)
199  endif
200 !---------------------------------------------------------------------
201 
202  if ( (mype.eq.pe_io .and. .not.bdist_io) .or. bdist_io) then
203  CALL flush_res (unit_res)
204  end if
205 
206 ! call MPI_Barrier(MPI_COMM_WORLD,mpierr) !//PAR_I/O enforce barrier here
207 !
208  deallocate (array1) !//d pnicol
209  deallocate (array2) !//d pnicol
210 ! call MPI_Barrier(MPI_COMM_WORLD,mpierr) !//PAR_I/O enforce barrier here
211 
212  call write_res1_netcdf
213 !
214  RETURN
215  END SUBROUTINE write_res1
216 
217  subroutine gatherwriteres(VAR, array2, array1, NEXT_REC)
219  USE geometry
220  USE funits
221  USE cdist
222  USE compar !//
223  USE mpi_utility !//d pnicol : for gather
224  USE sendrecv !//d pnicol : for gather
225 
226  USE cutcell
227  USE in_binary_512
228  USE param, only: dimension_3
229 
230  IMPLICIT NONE
231 
232  double precision, dimension(ijkmax2) :: array1
233  double precision, dimension(ijkmax3) :: array2
234  double precision, dimension(DIMENSION_3) :: VAR,TMP_VAR
235 
236  INTEGER :: NEXT_REC
237 
238 ! call MPI_Barrier(MPI_COMM_WORLD,mpierr) !//PAR_I/O enforce barrier here
239  if (.not.bdist_io) then
240 
241 
242  IF(re_indexing) THEN
243  CALL unshift_dp_array(var,tmp_var)
244  CALL gather (tmp_var,array2,root)
245  ELSE
246  CALL gather (var,array2,root)
247  ENDIF
248 ! call gather (VAR,array2,root) !//d pnicol
249 
250 ! call MPI_Barrier(MPI_COMM_WORLD,mpierr) !//PAR_I/O enforce barrier here
251  if (mype.eq.pe_io) then
252  call convert_to_io_dp(array2,array1,ijkmax2)
253  CALL out_bin_512 (unit_res, array1, ijkmax2, next_rec)
254  end if
255 
256  else
257 
258  IF(re_indexing) THEN
259  CALL unshift_dp_array(var,tmp_var)
260  CALL out_bin_512 (unit_res, tmp_var, size(tmp_var), next_rec)
261  ELSE
262  CALL out_bin_512 (unit_res, var, size(var), next_rec)
263  ENDIF
264 ! CALL OUT_BIN_512 (UNIT_RES, var, size(var), NEXT_REC)
265 
266  end if
267 
268  End subroutine gatherwriteres
269 
270 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
271 ! !
272 ! write_res1_netcdf !
273 ! !
274  subroutine write_res1_netcdf
276  USE param
277  USE param1
278  USE fldvar
279  USE geometry
280  USE physprop
281  USE run
282  USE scalars
283  USE rxns
284  USE cdist
285  USE compar
286  USE mpi_utility
287  USE mfix_netcdf
288  USE energy
289  USE in_binary_512
290 ! USE tmp_array
291 
292 
293  implicit none
294 
295  integer :: I , nn
296 
297  integer :: ncid , xyz_dimid
298  integer :: varid_time , t_dimid
299  integer :: dimids(1) , varid_epg , varid_pg,dims_time(1)
300  integer :: varid_pstar , varid_ug , varid_vg , varid_wg
301  integer :: varid_tg , varid_ropg
302  integer :: varid_rog , varid_gamaRG , varid_TRG
303  integer :: varid_gamaRS(20) , varid_TRS(20)
304 
305  integer :: varid_us(20) , varid_vs(20) , varid_ws(20) !! MMAX
306  integer :: varid_rops(20) , varid_ts(20) !! mmax
307  integer :: varid_thetam(20) !! mmax
308 
309  integer :: varid_xg(20) ! nmax(0)
310  integer :: varid_xs(20,20) ! mmax , MAX(nmax(1:mmax))
311 
312  integer :: varid_scalar(20) ! nscalar
313  integer :: varid_rr(20) ! nRR
314 
315  integer :: varid_kturbg , varid_eturbg
316 
317 
318  character(LEN=80) :: fname, var_name
319 
320  double precision, allocatable :: arr1(:)
321  double precision, allocatable :: arr2(:)
322 
323 ! bWrite_netcdf(1) : EP_g
324 ! bWrite_netcdf(2) : P_g
325 ! bWrite_netcdf(3) : P_star
326 ! bWrite_netcdf(4) : U_g / V_g / W_g
327 ! bWrite_netcdf(5) : U_s / V_s / W_s
328 ! bWrite_netcdf(6) : ROP_s
329 ! bWrite_netcdf(7) : T_g
330 ! bWrite_netcdf(8) : T_s
331 ! bWrite_netcdf(9) : X_g
332 ! bWrite_netcdf(10) : X_s
333 ! bWrite_netcdf(11) : Theta_m
334 ! bWrite_netCDF(12) : Scalar
335 ! bWrite_netCDF(13) : ReactionRates
336 ! bWrite_netCDF(14) : k_turb_g , e_turb_g
337 
338  if (.not. mfix_usingnetcdf()) return
339  if (.not. bglobalnetcdf) return
340 
341  return ! for now ... do not write restart files using netCDF
342 
343  if (mype .eq. pe_io) then
344  allocate (arr1(ijkmax2))
345  allocate (arr2(ijkmax3))
346  else
347  allocate (arr1(1))
348  allocate (arr2(1))
349  end if
350 
351 ! call mpi$barrier(MPI_COMM_WORLD,mpierr)
352  if (mype .ne. pe_io) goto 1234
353 
354  fname = trim(run_name) // "_RES1.nc"
355  call mfix_check_netcdf( mfix_nf90_create(fname, nf90_clobber, ncid) )
356 
357  call mfix_check_netcdf( mfix_nf90_def_dim(ncid, "xyz", ijkmax2, xyz_dimid) )
358  call mfix_check_netcdf( mfix_nf90_def_dim(ncid, "t" , 1 , t_dimid) )
359 
360  ! The dimids array is used to pass the IDs of the dimensions of
361  ! the variables. Note that in fortran arrays are stored in
362  ! column-major format.
363  dims_time(1) = t_dimid
364  dimids = (/ xyz_dimid /)
365 
366  call mfix_check_netcdf( mfix_nf90_def_var(ncid ,"time" , nf90_double ,dims_time , varid_time ) )
367  call mfix_check_netcdf( mfix_nf90_def_var(ncid, "EP_g" , nf90_double, dimids , varid_epg ) )
368  call mfix_check_netcdf( mfix_nf90_def_var(ncid, "P_g" , nf90_double, dimids , varid_pg ) )
369  call mfix_check_netcdf( mfix_nf90_def_var(ncid, "P_star", nf90_double, dimids , varid_pstar) )
370  call mfix_check_netcdf( mfix_nf90_def_var(ncid, "U_g" , nf90_double, dimids , varid_ug ) )
371  call mfix_check_netcdf( mfix_nf90_def_var(ncid, "V_g" , nf90_double, dimids , varid_vg ) )
372  call mfix_check_netcdf( mfix_nf90_def_var(ncid, "W_g" , nf90_double, dimids , varid_wg ) )
373  call mfix_check_netcdf( mfix_nf90_def_var(ncid, "T_g" , nf90_double, dimids , varid_tg ) )
374  call mfix_check_netcdf( mfix_nf90_def_var(ncid, "ROP_g" , nf90_double, dimids , varid_ropg ) )
375  call mfix_check_netcdf( mfix_nf90_def_var(ncid, "RO_g" , nf90_double, dimids , varid_rog ) )
376  call mfix_check_netcdf( mfix_nf90_def_var(ncid, "gamaRG" , nf90_double, dimids , varid_gamarg ) )
377  call mfix_check_netcdf( mfix_nf90_def_var(ncid, "TRG" , nf90_double, dimids , varid_trg ) )
378 
379  do i = 1,1 ! mmax
380  var_name = 'U_s_xxx'
381  write (var_name(5:7),'(i3.3)') i
382  call mfix_check_netcdf( mfix_nf90_def_var(ncid, var_name, nf90_double, dimids, varid_us(i)) )
383 
384  var_name = 'V_s_xxx'
385  write (var_name(5:7),'(i3.3)') i
386  call mfix_check_netcdf( mfix_nf90_def_var(ncid, var_name, nf90_double, dimids, varid_vs(i)) )
387 
388  var_name = 'W_s_xxx'
389  write (var_name(5:7),'(i3.3)') i
390  call mfix_check_netcdf( mfix_nf90_def_var(ncid, var_name, nf90_double, dimids, varid_ws(i)) )
391 
392  var_name = 'ROP_s_xxx'
393  write (var_name(7:10),'(i3.3)') i
394  call mfix_check_netcdf( mfix_nf90_def_var(ncid, var_name, nf90_double, dimids, varid_rops(i)) )
395 
396  var_name = 'T_s_xxx'
397  write (var_name(5:7),'(i3.3)') i
398  call mfix_check_netcdf( mfix_nf90_def_var(ncid, var_name, nf90_double, dimids, varid_ts(i)) )
399 
400  var_name = 'Theta_m_xxx'
401  write (var_name(9:11),'(i3.3)') i
402  call mfix_check_netcdf( mfix_nf90_def_var(ncid, var_name, nf90_double, dimids, varid_thetam(i)) )
403 
404  var_name = 'gamaRS_xxx'
405  write (var_name(8:10),'(i3.3)') i
406  call mfix_check_netcdf( mfix_nf90_def_var(ncid, var_name, nf90_double, dimids, varid_gamars(i)) )
407 
408  var_name = 'TRS_xxx'
409  write (var_name(5:7),'(i3.3)') i
410  call mfix_check_netcdf( mfix_nf90_def_var(ncid, var_name, nf90_double, dimids, varid_trs(i)) )
411 
412  DO nn = 1, nmax(i)
413  var_name = 'X_s_xxx_xxx'
414  write (var_name(5:7) ,'(i3.3)') i
415  write (var_name(9:11),'(i3.3)') nn
416  call mfix_check_netcdf( mfix_nf90_def_var(ncid, var_name, nf90_double, dimids, varid_xs(i,nn)) )
417  END DO
418 
419 
420  end do
421 
422  do i = 1,nmax(0)
423  var_name = 'X_g_xxx'
424  write (var_name(5:7),'(i3.3)') i
425  call mfix_check_netcdf( mfix_nf90_def_var(ncid, var_name, nf90_double, dimids, varid_xg(i)) )
426  end do
427 
428  do i = 1,nscalar
429  var_name = 'Scalar_xxx'
430  write (var_name(8:10),'(i3.3)') i
431  call mfix_check_netcdf( mfix_nf90_def_var(ncid, var_name, nf90_double, dimids, varid_scalar(i)) )
432  end do
433 
434  do i = 1,nrr
435  var_name = 'RRates_xxx'
436  write (var_name(8:10),'(i3.3)') i
437  call mfix_check_netcdf( mfix_nf90_def_var(ncid, var_name, nf90_double, dimids, varid_rr(i)) )
438  end do
439 
440 
441  if (k_epsilon) then
442  call mfix_check_netcdf( mfix_nf90_def_var(ncid, 'k_turb_g', nf90_double, dimids, varid_kturbg) )
443  call mfix_check_netcdf( mfix_nf90_def_var(ncid, 'e_turb_g', nf90_double, dimids, varid_eturbg) )
444  end if
445 
446 
447  call mfix_check_netcdf( mfix_nf90_enddef(ncid) )
448 
449  1234 continue
450 ! call mpi$barrier(MPI_COMM_WORLD,mpierr)
451 
452  if (mype .eq. pe_io) then
453 ! call check_netcdf( nf90_put_var(ncid,varid_time,the_time) ) ! ???
454  call mfix_check_netcdf( mfix_nf90_put_var(ncid,varid_time,time) ) ! ???
455  end if
456 
457  call gather(ep_g,arr2,root)
458 
459  if (mype .eq. pe_io) then
460  call convert_to_io_dp(arr2,arr1,ijkmax2)
461  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_epg, arr1) )
462  end if
463 
464  call gather(p_g,arr2,root)
465  if (mype .eq. pe_io) then
466  call convert_to_io_dp(arr2,arr1,ijkmax2)
467  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_pg, arr1) )
468  end if
469 
470  call gather(p_star,arr2,root)
471  if (mype .eq. pe_io) then
472  call convert_to_io_dp(arr2,arr1,ijkmax2)
473  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_pstar, arr1) )
474  end if
475 
476  call gather(u_g,arr2,root)
477  if (mype .eq. pe_io) then
478  call convert_to_io_dp(arr2,arr1,ijkmax2)
479  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_ug, arr1) )
480  end if
481 
482  call gather(v_g,arr2,root)
483  if (mype .eq. pe_io) then
484  call convert_to_io_dp(arr2,arr1,ijkmax2)
485  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_vg, arr1) )
486  end if
487 
488  call gather(w_g,arr2,root)
489  if (mype .eq. pe_io) then
490  call convert_to_io_dp(arr2,arr1,ijkmax2)
491  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_wg, arr1) )
492  end if
493 
494  call gather(t_g,arr2,root)
495  if (mype .eq. pe_io) then
496  call convert_to_io_dp(arr2,arr1,ijkmax2)
497  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_tg, arr1) )
498  end if
499 
500  call gather(gama_rg,arr2,root)
501  if (mype .eq. pe_io) then
502  call convert_to_io_dp(arr2,arr1,ijkmax2)
503  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_gamarg, arr1) )
504  end if
505 
506  call gather(ro_g,arr2,root)
507  if (mype .eq. pe_io) then
508  call convert_to_io_dp(arr2,arr1,ijkmax2)
509  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_rog, arr1) )
510  end if
511 
512  call gather(rop_g,arr2,root)
513  if (mype .eq. pe_io) then
514  call convert_to_io_dp(arr2,arr1,ijkmax2)
515  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_ropg, arr1) )
516  end if
517 
518  call gather(t_rg,arr2,root)
519  if (mype .eq. pe_io) then
520  call convert_to_io_dp(arr2,arr1,ijkmax2)
521  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_trg, arr1) )
522  end if
523 
524  do i = 1,1 ! mmax
525 
526  call gather(u_s(:,i),arr2,root)
527  if (mype .eq. pe_io) then
528  call convert_to_io_dp(arr2,arr1,ijkmax2)
529  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_us(i), arr1) )
530  end if
531 
532  call gather(v_s(:,i),arr2,root)
533  if (mype .eq. pe_io) then
534  call convert_to_io_dp(arr2,arr1,ijkmax2)
535  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_vs(i), arr1) )
536  end if
537 
538 
539  call gather(w_s(:,i),arr2,root)
540  if (mype .eq. pe_io) then
541  call convert_to_io_dp(arr2,arr1,ijkmax2)
542  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_ws(i), arr1) )
543  end if
544 
545  call gather(rop_s(:,i),arr2,root)
546  if (mype .eq. pe_io) then
547  call convert_to_io_dp(arr2,arr1,ijkmax2)
548  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_rops(i), arr1 ) )
549  end if
550 
551  call gather(t_s(:,i),arr2,root)
552  if (mype .eq. pe_io) then
553  call convert_to_io_dp(arr2,arr1,ijkmax2)
554  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_ts(i), arr1 ) )
555  end if
556 
557  call gather(theta_m(:,i),arr2,root)
558  if (mype .eq. pe_io) then
559  call convert_to_io_dp(arr2,arr1,ijkmax2)
560  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_thetam(i), arr1 ) )
561  end if
562 
563  call gather(gama_rs(:,i),arr2,root)
564  if (mype .eq. pe_io) then
565  call convert_to_io_dp(arr2,arr1,ijkmax2)
566  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_gamars(i), arr1 ) )
567  end if
568 
569  call gather(t_rs(:,i),arr2,root)
570  if (mype .eq. pe_io) then
571  call convert_to_io_dp(arr2,arr1,ijkmax2)
572  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_trs(i) , arr1 ) )
573  end if
574 
575  do nn = 1,nmax(i)
576  call gather(x_s(:,i,nn),arr2,root)
577  if (mype .eq. pe_io) then
578  call convert_to_io_dp(arr2,arr1,ijkmax2)
579  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_xs(i,nn), arr1 ) )
580  end if
581  end do
582 
583  end do
584 
585  do i = 1,nmax(0)
586  call gather(x_g(:,i),arr2,root)
587  if (mype .eq. pe_io) then
588  call convert_to_io_dp(arr2,arr1,ijkmax2)
589  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_xg(i), arr1 ) )
590  end if
591  end do
592 
593  do i = 1,nscalar
594  call gather(scalar(:,i),arr2,root)
595  if (mype .eq. pe_io) then
596  call convert_to_io_dp(arr2,arr1,ijkmax2)
597  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_scalar(i), arr1 ) )
598  end if
599  end do
600 
601  do i = 1,nrr
602  call gather(reactionrates(:,i),arr2,root)
603  if (mype .eq. pe_io) then
604  call convert_to_io_dp(arr2,arr1,ijkmax2)
605  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_rr(i), arr1 ) )
606  end if
607  end do
608 
609  if (k_epsilon) then
610  call gather(k_turb_g,arr2,root)
611  if (mype .eq. pe_io) then
612  call convert_to_io_dp(arr2,arr1,ijkmax2)
613  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_kturbg, arr1) )
614  end if
615 
616  call gather(e_turb_g,arr2,root)
617  if (mype .eq. pe_io) then
618  call convert_to_io_dp(arr2,arr1,ijkmax2)
619  call mfix_check_netcdf( mfix_nf90_put_var(ncid, varid_eturbg, arr1) )
620  end if
621  end if
622 
623  2222 continue
624 
625  ! Close the file. This frees up any internal netCDF resources
626  ! associated with the file, and flushes any buffers.
627  if (mype .eq. pe_io) then
628  call mfix_check_netcdf( mfix_nf90_close(ncid) )
629  end if
630 
631  deallocate (arr1)
632  deallocate (arr2)
633 
634  return
635 
636  end subroutine write_res1_netcdf
637 
638 
639 
double precision, dimension(:,:), allocatable v_s
Definition: fldvar_mod.f:105
logical re_indexing
Definition: cutcell_mod.f:16
integer function mfix_nf90_def_dim(ncid, name, len, dimid)
double precision, dimension(:), allocatable ep_g
Definition: fldvar_mod.f:15
integer function mfix_nf90_create(path, cmode, ncid, initialsize, chunksize)
logical bdist_io
Definition: cdist_mod.f:4
subroutine mfix_check_netcdf(status)
double precision, dimension(:), allocatable k_turb_g
Definition: fldvar_mod.f:161
subroutine gatherwriteres(VAR, array2, array1, NEXT_REC)
Definition: write_res1.f:218
integer dimension_3
Definition: param_mod.f:11
integer ijkmax2
Definition: geometry_mod.f:80
double precision, dimension(:,:), allocatable w_s
Definition: fldvar_mod.f:117
Definition: rxns_mod.f:1
character(len=60) run_name
Definition: run_mod.f:24
double precision, dimension(:), allocatable t_g
Definition: fldvar_mod.f:63
integer nf90_clobber
integer function mfix_nf90_close(ncid)
logical, dimension(dim_m) solve_ros
Definition: run_mod.f:250
double precision, dimension(:,:), allocatable scalar
Definition: fldvar_mod.f:155
double precision dt
Definition: run_mod.f:51
subroutine flush_res(iunit)
Definition: machine_mod.f:235
double precision, dimension(:), allocatable p_g
Definition: fldvar_mod.f:26
double precision, dimension(:,:), allocatable u_s
Definition: fldvar_mod.f:93
integer ijkmax3
Definition: geometry_mod.f:82
integer pe_io
Definition: compar_mod.f:30
integer mmax
Definition: physprop_mod.f:19
double precision, dimension(:,:,:), allocatable x_s
Definition: fldvar_mod.f:78
subroutine write_res1
Definition: write_res1.f:26
double precision, dimension(:,:), allocatable t_s
Definition: fldvar_mod.f:66
double precision, dimension(:,:), allocatable x_g
Definition: fldvar_mod.f:75
integer root
Definition: compar_mod.f:41
double precision, dimension(:,:), allocatable t_rs
Definition: energy_mod.f:24
double precision, dimension(:), allocatable t_rg
Definition: energy_mod.f:21
Definition: cdist_mod.f:2
subroutine write_res1_netcdf
Definition: write_res1.f:275
double precision, dimension(:,:), allocatable theta_m
Definition: fldvar_mod.f:149
subroutine out_bin_512(IUNIT, ARRAY, N, NEXT_REC)
Definition: out_bin_512.f:24
double precision, dimension(:), allocatable v_g
Definition: fldvar_mod.f:99
double precision, dimension(:), allocatable w_g
Definition: fldvar_mod.f:111
Definition: run_mod.f:13
integer nrr
Definition: rxns_mod.f:10
Definition: param_mod.f:2
double precision, dimension(:,:), allocatable ro_s
Definition: fldvar_mod.f:45
integer, dimension(0:dim_m) nmax
Definition: physprop_mod.f:119
logical k_epsilon
Definition: run_mod.f:97
integer mype
Definition: compar_mod.f:24
double precision, dimension(:,:), allocatable reactionrates
Definition: rxns_mod.f:7
double precision, dimension(:), allocatable p_star
Definition: fldvar_mod.f:142
double precision, dimension(:), allocatable gama_rg
Definition: energy_mod.f:15
integer nstep
Definition: run_mod.f:60
logical function mfix_usingnetcdf()
double precision, dimension(:), allocatable u_g
Definition: fldvar_mod.f:87
integer nscalar
Definition: scalars_mod.f:7
integer, parameter unit_res
Definition: funits_mod.f:27
integer function mfix_nf90_enddef(ncid, h_minfree, v_align, v_minfree, r_align)
double precision, dimension(:,:), allocatable rop_s
Definition: fldvar_mod.f:51
logical bglobalnetcdf
Definition: cdist_mod.f:14
double precision, dimension(:,:), allocatable gama_rs
Definition: energy_mod.f:18
double precision time
Definition: run_mod.f:45
subroutine unshift_dp_array(ARRAY_1, ARRAY_2)
double precision, dimension(:), allocatable e_turb_g
Definition: fldvar_mod.f:162
subroutine convert_to_io_dp(arr_internal, arr_io, nn)
double precision, dimension(:), allocatable ro_g
Definition: fldvar_mod.f:32
double precision, dimension(:), allocatable rop_g
Definition: fldvar_mod.f:38
integer nf90_double