MFIX  2016-1
read_res1.f
Go to the documentation of this file.
1 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvC
2 ! C
3 ! Module name: READ_RES1 C
4 ! Purpose: read in the time-dependent restart records C
5 ! C
6 ! Author: P. Nicoletti Date: 03-JAN-92 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: IJKMAX2, MMAX, DT C
17 ! Variables modified: TIME, NSTEP, EP_g, P_g, P_star, RO_g C
18 ! ROP_g, T_g, T_s, U_g, V_g, W_g, ROP_s C
19 ! U_s, V_s, W_s C
20 ! C
21 ! Local variables: TIME_READ, LC, NEXT_REC C
22 ! C
23 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^C
24 !
25  SUBROUTINE read_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 param
33  USE param1
34  USE fldvar
35  USE geometry
36  USE physprop
37  USE run
38  USE rxns
39  USE scalars
40  USE funits
41  USE energy
42  USE compar
43  USE cdist
44  USE mpi_utility
45  USE sendrecv
46  USE in_binary_512
47 
48  IMPLICIT NONE
49 !-----------------------------------------------
50 ! G l o b a l P a r a m e t e r s
51 !-----------------------------------------------
52 !-----------------------------------------------
53 ! L o c a l P a r a m e t e r s
54 !-----------------------------------------------
55 !-----------------------------------------------
56 ! L o c a l V a r i a b l e s
57 !-----------------------------------------------
58 !
59 ! loop counter
60  INTEGER LC
61 !
62 ! Local species index
63  INTEGER NN
64 !
65 ! pointer to the next record
66  INTEGER NEXT_REC
67 !
68 ! file version id
69  CHARACTER(LEN=512) :: VERSION
70 !
71 ! version number
72  REAL VERSION_NUMBER
73 !
74 ! Dummy array
75  DOUBLE PRECISION DT_SAVE
76 
77 !//PAR_I/O declare global scratch arrays
78  double precision, allocatable :: array1(:)
79  double precision, allocatable :: array2(:)
80 !-----------------------------------------------
81 !
82  if (mype .eq. pe_io .or. .not.bstart_with_one_res) then
83  allocate (array1(ijkmax2))
84  allocate (array2(ijkmax3))
85  else
86  allocate (array1(1))
87  allocate (array2(1))
88  end if
89 
90 ! call MPI_barrier(MPI_COMM_WORLD,mpierr)
91 !
92 ! Use DT from data file if DT_FAC is set to 1.0
93  IF (dt_fac == one) dt_save = dt
94 !
95 
96 !//PAR_I/O only PE_IO reads the restart file
97  if (mype == pe_io .or. (bdist_io .and. .not.bstart_with_one_res)) then
98  READ (unit_res, rec=1) version
99  READ (version(6:512), *) version_number
100 
101  READ (unit_res, rec=3) next_rec
102  IF (version_number >= 1.12) THEN
103  READ (unit_res, rec=next_rec) time, dt, nstep
104  ELSE
105  READ (unit_res, rec=next_rec) time, nstep
106  ENDIF
107  next_rec = next_rec + 1
108  end if
109 
110 
111 
112  if (.not.bdist_io .or. bstart_with_one_res) then
113 ! call MPI_barrier(MPI_COMM_WORLD,mpierr)
114 
115  call bcast(version, pe_io) !//PAR_I/O BCAST0c
116  call bcast(version_number, pe_io) !//PAR_I/O BCAST0r
117  call bcast(time, pe_io) !//PAR_I/O BCAST0d
118  call bcast(nstep, pe_io) !//PAR_I/O BCAST0i
119  if (version_number >= 1.12) call bcast(dt, pe_io) !//PAR_I/O BCAST0d
120  end if
121 ! call MPI_barrier(MPI_COMM_WORLD,mpierr)
122 
123 !AE TIME 091501 Store the timestep counter level at the begin of RESTART run
124  nsteprst = nstep
125 
126 ! for now ... do not do RES1 file in netCDF format
127 ! call read_res1_netcdf
128 ! goto 999
129 !
130 
131 
132 !
133  call readscatterres(ep_g,array2, array1, 0, next_rec)
134 
135  call readscatterres(p_g,array2, array1, 0, next_rec)
136 
137  call readscatterres(p_star,array2, array1, 1, next_rec)
138 
139  call readscatterres(ro_g,array2, array1, 0, next_rec)
140 
141  call readscatterres(rop_g,array2, array1, 0, next_rec)
142 
143  call readscatterres(t_g,array2, array1, 1, next_rec)
144 !
145 
146  IF (version_number < 1.15) THEN
147  call readscatterres (t_s(:,1),array2, array1, 1, next_rec)
148 
149  IF (mmax >= 2) THEN
150  call readscatterres (t_s(:,2),array2, array1, 1, next_rec)
151  ELSE
152  if (mype == pe_io) &
153  CALL in_bin_512 (unit_res, array1, ijkmax2, next_rec)
154  ENDIF
155  ENDIF
156 
157  IF (version_number >= 1.05) THEN
158  DO nn = 1, nmax(0)
159  call readscatterres (x_g(:,nn),array2, array1, 1, next_rec)
160  END DO
161  ENDIF
162 
163  call readscatterres(u_g, array2, array1, 0, next_rec)
164  call readscatterres(v_g, array2, array1, 0, next_rec)
165  call readscatterres(w_g, array2, array1, 0, next_rec)
166 !
167  DO lc = 1, mmax
168  call readscatterres(rop_s(:,lc), array2, array1, 0, next_rec)
169 
170  IF(any(solve_ros)) &
171  CALL readscatterres(ro_s(:,lc), array2, array1, 0, next_rec)
172 
173  IF (version_number >= 1.15) THEN
174  call readscatterres(t_s(:,lc), array2, array1, 1, next_rec)
175  END IF
176  call readscatterres(u_s(:,lc), array2, array1, 0, next_rec)
177  call readscatterres(v_s(:,lc), array2, array1, 0, next_rec)
178  call readscatterres(w_s(:,lc), array2, array1, 0, next_rec)
179 
180  IF (version_number >= 1.2) then
181  call readscatterres(theta_m(:,lc), array2, array1, 1, next_rec)
182  end if
183  IF (version_number >= 1.05) THEN
184  DO nn = 1, nmax(lc)
185  call readscatterres(x_s(:,lc,nn), array2, array1, 1, next_rec)
186  END DO
187  ENDIF
188  END DO
189 
190  IF (version_number >= 1.3) THEN
191  DO nn = 1, nscalar
192  call readscatterres(scalar(:,nn), array2, array1, 1, next_rec)
193  END DO
194  ENDIF
195 
196  IF (version_number >= 1.4) THEN
197  call readscatterres(gama_rg, array2, array1, 1, next_rec)
198 
199  call readscatterres(t_rg, array2, array1, 0, next_rec)
200 
201  DO lc = 1, mmax
202  call readscatterres(gama_rs(1,lc), array2, array1, 1, next_rec)
203 
204  call readscatterres(t_rs(1,lc), array2, array1, 0, next_rec)
205 
206  ENDDO
207  ELSE
208  gama_rg(:) = zero
209  t_rg(:) = zero
210  gama_rs(:,:) = zero
211  t_rs(:,:) = zero
212  ENDIF
213 
214 
215  IF (version_number >= 1.5) THEN
216  DO nn = 1, nrr
217  call readscatterres(reactionrates(:,nn), array2, array1, 1, next_rec)
218  END DO
219  ENDIF
220 
221  IF (version_number >= 1.6 .AND. k_epsilon) THEN
222  call readscatterres(k_turb_g, array2, array1, 1, next_rec)
223  call readscatterres(e_turb_g, array2, array1, 1, next_rec)
224  ENDIF
225 !------------------------------------------------------------------------
226 !
227 
228 ! call MPI_barrier(MPI_COMM_WORLD,mpierr)
229  deallocate( array1 )
230  deallocate( array2 )
231 ! call MPI_barrier(MPI_COMM_WORLD,mpierr)
232 
233  if (.not.bdist_io .or. bstart_with_one_res) then
234  call send_recv(rop_g)
235  call send_recv(ro_g)
236  call send_recv(rop_s)
237  end if
238 
239 
240  IF (dt_fac == one) dt = dt_save
241 !
242 
243 ! We may no longer need PATCH_AFTER_RESTART
244 ! CALL PATCH_AFTER_RESTART
245 
246  RETURN
247  END SUBROUTINE read_res1
248 
249  subroutine readscatterres(VAR, array2, array1, init, NEXT_REC)
250  use param1, only: zero, undefined
251  use param, only: dimension_3
252  USE geometry
253  USE funits
254  USE compar
255  USE cdist
256  USE mpi_utility
257  USE sendrecv
258  USE in_binary_512
259  IMPLICIT NONE
260  double precision, dimension(ijkmax2) :: array1
261  double precision, dimension(ijkmax3) :: array2
262  double precision, dimension(DIMENSION_3) :: VAR
263  INTEGER :: init ! define VAR initialization, 0: undefin, 1: zero
264  INTEGER :: NEXT_REC
265 
266 !// Reset global scratch arrays
267  if( init==0 ) then
268  array1(:) = undefined
269  array2(:) = undefined
270  else
271  array1(:) = zero
272  array2(:) = zero
273  endif
274 
275  if (.not.bdist_io .or. bstart_with_one_res) then
276  if (mype == pe_io) then
277  CALL in_bin_512 (unit_res, array1, ijkmax2, next_rec)
278  CALL convert_from_io_dp(array1, array2, ijkmax2)
279  end if
280 ! call MPI_barrier(MPI_COMM_WORLD,mpierr)
281  call scatter(var, array2, pe_io)
282 ! call MPI_barrier(MPI_COMM_WORLD,mpierr)
283  else
284  CALL in_bin_512 (unit_res, var, size(var) , next_rec)
285  end if
286 
287  End subroutine readscatterres
288 
289 
290  subroutine readscatterres_netcdf(VAR, array2, array1, ncid , varid)
291  USE param, only: dimension_3
292  USE geometry
293  USE compar
294  USE cdist
295  USE mpi_utility
296  USE sendrecv
297  USE mfix_netcdf
298  USE in_binary_512
299 
300  IMPLICIT NONE
301 
302  double precision, dimension(ijkmax2) :: array1
303  double precision, dimension(ijkmax3) :: array2
304  double precision, dimension(DIMENSION_3) :: VAR
305 
306  integer :: ncid , varid
307 
308 
309 
310  if (mype .eq. pe_io) then
311  call mfix_check_netcdf( mfix_nf90_get_var(ncid , varid , array1 ) )
312  call convert_from_io_dp(array1,array2,ijkmax2)
313  end if
314 
315 ! call MPI_barrier(MPI_COMM_WORLD,mpierr)
316  call scatter(var,array2,pe_io)
317 ! call MPI_barrier(MPI_COMM_WORLD,mpierr)
318 
319  End subroutine readscatterres_netcdf
320 
321 
322 
323 
324 
325 !// Comments on the modifications for DMP version implementation
326 !// 001 Include header file and common declarations for parallelization
327 !// 020 New local variables for parallelization: array1, array2
328 !// 400 Added sendrecv module and send_recv calls for COMMunication
329 !// 400 Added mpi_utility module and other global reduction (bcast) calls
330 
331 
332 
333 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
334 ! !
335 ! read_res1_netcdf !
336 ! !
337  subroutine read_res1_netcdf
339  USE param
340  USE param1
341  USE fldvar
342  USE geometry
343  USE physprop
344  USE run
345 ! USE funits
346  USE scalars
347 ! USE output
348  USE rxns
349  USE cdist
350  USE compar
351  USE mpi_utility
352  USE mfix_netcdf
353 ! USE tmp_array
354  USE energy
355 
356 
357  implicit none
358 
359  integer :: I , nn
360 
361  integer :: ncid
362  integer :: varid_time
363  integer :: varid_epg , varid_pg
364  integer :: varid_pstar , varid_ug , varid_vg , varid_wg
365  integer :: varid_tg
366  integer :: varid_rog , varid_gamaRG , varid_TRG
367  integer :: varid_gamaRS(20) , varid_TRS(20) , varid_ropg
368 
369  integer :: varid_us(20) , varid_vs(20) , varid_ws(20) !! MMAX
370  integer :: varid_rops(20) , varid_ts(20) !! mmax
371  integer :: varid_thetam(20) !! mmax
372 
373  integer :: varid_xg(20) ! nmax(0)
374  integer :: varid_xs(20,20) ! mmax , MAX(nmax(1:mmax))
375 
376  integer :: varid_scalar(20) ! nscalar
377  integer :: varid_rr(20) ! nRR
378 
379  integer :: varid_kturbg , varid_eturbg
380 
381 
382  character(LEN=80) :: fname, var_name
383 
384  integer nDim , nVars , nAttr , unID , formatNUM
385 
386  integer xyz_id , xyz_dim
387 
388  character(LEN=80) :: varname
389  integer vartype,nvdims,vdims(10),nvatts,rcode
390 
391  double precision, allocatable :: array1(:)
392  double precision, allocatable :: array2(:)
393 
394 
395 ! integer :: MFIX_nf90_create
396 ! integer :: MFIX_nf90_def_dim
397 ! integer :: MFIX_nf90_close
398 ! integer :: MFIX_nf90_open
399 ! integer :: MFIX_nf90_inquire
400 ! integer :: MFIX_nf90_inq_dimid
401 ! integer :: MFIX_nf90_inquire_dimension
402 ! integer :: MFIX_nf90_inq_varid
403 
404 
405 
406 
407 
408 
409 
410 
411 !
412 !
413 
414 
415 ! bWrite_netcdf(1) : EP_g
416 ! bWrite_netcdf(2) : P_g
417 ! bWrite_netcdf(3) : P_star
418 ! bWrite_netcdf(4) : U_g / V_g / W_g
419 ! bWrite_netcdf(5) : U_s / V_s / W_s
420 ! bWrite_netcdf(6) : ROP_s
421 ! bWrite_netcdf(7) : T_g
422 ! bWrite_netcdf(8) : T_s
423 ! bWrite_netcdf(9) : X_g
424 ! bWrite_netcdf(10) : X_s
425 ! bWrite_netcdf(11) : Theta_m
426 ! bWrite_netCDF(12) : Scalar
427 ! bWrite_netCDF(13) : ReactionRates
428 ! bWrite_netCDF(14) : k_turb_g , e_turb_g
429 
430  ! return
431 
432  if (.not. mfix_usingnetcdf()) return
433 
434  if (mype .eq. pe_io) then
435  allocate (array1(ijkmax2))
436  allocate (array2(ijkmax3))
437  else
438  allocate (array1(1))
439  allocate (array2(1))
440  end if
441 
442 
443 
444 ! call MPI_barrier(MPI_COMM_WORLD,mpierr)
445  if (mype .eq. pe_io) then
446  fname = trim(run_name) // "_RES1.nc"
447 
448  call mfix_check_netcdf( mfix_nf90_open(fname, nf90_nowrite, ncid) )
449 
450  call mfix_check_netcdf( mfix_nf90_inquire(ncid, ndim , nvars , nattr , unid , formatnum) )
451 
452 ! write (*,*) ' nDim = ' , nDim
453 ! write (*,*) ' nVars = ' , nVars
454 ! write (*,*) ' nAttr = ' , nAttr
455 ! write (*,*) ' unID = ' , unID
456 ! write (*,*) ' formatNum = ' , formatNum
457 
458  call mfix_check_netcdf( mfix_nf90_inq_dimid(ncid,"xyz",xyz_id) )
459  call mfix_check_netcdf( mfix_nf90_inquire_dimension(ncid,xyz_id,len=xyz_dim) )
460 
461 
462  do i = 1,nvars
463  call mfix_ncvinq(ncid,i,varname,vartype,nvdims,vdims,nvatts,rcode)
464  end do
465 
466  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid ,"time" , varid_time ) )
467  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, "EP_g" , varid_epg ) )
468  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, "P_g" , varid_pg ) )
469  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, "P_star", varid_pstar ) )
470  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, "U_g" , varid_ug ) )
471  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, "V_g" , varid_vg ) )
472  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, "W_g" , varid_wg ) )
473  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, "T_g" , varid_tg ) )
474  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, "RO_g" , varid_rog ) )
475  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, "ROP_g" , varid_ropg ) )
476  call mfix_check_netcdf( mfix_nf90_get_var(ncid,varid_time,time) )
477 
478  end if
479 
480 
481 ! if (myPE .eq. PE_IO) then
482 ! call check_netcdf( nf90_get_var(ncid , varid_epg , array1 ) )
483 ! call convert_from_io_dp(array1,array2,ijkmax2)
484 ! end if
485 ! call scatter(EP_g,array2,PE_IO)
486 
487  call readscatterres_netcdf(ep_g , array2, array1, ncid , varid_epg)
488  call readscatterres_netcdf(p_g , array2, array1, ncid , varid_pg)
489  call readscatterres_netcdf(p_star , array2, array1, ncid , varid_pstar)
490  call readscatterres_netcdf(ro_g , array2, array1, ncid , varid_rog)
491  call readscatterres_netcdf(rop_g , array2, array1, ncid , varid_ropg)
492  call readscatterres_netcdf(u_g , array2, array1, ncid , varid_ug)
493  call readscatterres_netcdf(v_g , array2, array1, ncid , varid_vg)
494  call readscatterres_netcdf(w_g , array2, array1, ncid , varid_wg)
495  call readscatterres_netcdf(t_g , array2, array1, ncid , varid_tg)
496 
497  do i = 1,1 ! mmax
498 
499  if (mype .eq. pe_io) then
500  var_name = 'U_s_xxx'
501  write (var_name(5:7),'(i3.3)') i
502  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, var_name, varid_us(i)) )
503  end if
504  call readscatterres_netcdf(u_s(:,i) , array2, array1, ncid , varid_us(i))
505 
506  if (mype .eq. pe_io) then
507  var_name = 'V_s_xxx'
508  write (var_name(5:7),'(i3.3)') i
509  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, var_name, varid_vs(i)) )
510  end if
511  call readscatterres_netcdf(v_s(:,i) , array2, array1, ncid , varid_vs(i))
512 
513  if (mype .eq. pe_io) then
514  var_name = 'W_s_xxx'
515  write (var_name(5:7),'(i3.3)') i
516  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, var_name, varid_ws(i)) )
517  end if
518  call readscatterres_netcdf(w_s(:,i) , array2, array1, ncid , varid_ws(i))
519 
520  if (mype .eq. pe_io) then
521  var_name = 'ROP_s_xxx'
522  write (var_name(7:10),'(i3.3)') i
523  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, var_name, varid_rops(i)) )
524  end if
525  call readscatterres_netcdf(rop_s(:,i) , array2, array1, ncid , varid_rops(i))
526 
527  if (mype .eq. pe_io) then
528  var_name = 'T_s_xxx'
529  write (var_name(5:7),'(i3.3)') i
530  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, var_name, varid_ts(i)) )
531  end if
532  call readscatterres_netcdf(t_s(:,i) , array2, array1, ncid , varid_ts(i))
533 
534  if (mype .eq. pe_io) then
535  var_name = 'Theta_m_xxx'
536  write (var_name(9:11),'(i3.3)') i
537  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, var_name, varid_thetam(i)) )
538  end if
539  call readscatterres_netcdf(theta_m(:,i) , array2, array1, ncid , varid_thetam(i))
540 
541  if (mype .eq. pe_io) then
542  var_name = 'gamaRS_xxx'
543  write (var_name(8:10),'(i3.3)') i
544  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, var_name, varid_gamars(i)) )
545  end if
546  call readscatterres_netcdf(gama_rs(:,i) , array2, array1, ncid , varid_gamars(i))
547 
548  if (mype .eq. pe_io) then
549  var_name = 'TRS_xxx'
550  write (var_name(5:7),'(i3.3)') i
551  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, var_name, varid_trs(i)) )
552  end if
553  call readscatterres_netcdf(t_rs(:,i) , array2, array1, ncid , varid_trs(i))
554 
555  DO nn = 1, nmax(i)
556  if (mype .eq. pe_io) then
557  var_name = 'X_s_xxx_xxx'
558  write (var_name(5:7) ,'(i3.3)') i
559  write (var_name(9:11),'(i3.3)') nn
560  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, var_name, varid_xs(i,nn)) )
561  end if
562  call readscatterres_netcdf(x_s(:,i,nn) , array2, array1, ncid , varid_xs(i,nn))
563  END DO
564 
565  end do
566 
567  do i = 1,nmax(0)
568  if (mype .eq. pe_io) then
569  var_name = 'X_g_xxx'
570  write (var_name(5:7),'(i3.3)') i
571  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, var_name, varid_xg(i)) )
572  end if
573  call readscatterres_netcdf(x_g(:,i) , array2, array1, ncid , varid_xg(i))
574  end do
575 
576  do i = 1,nscalar
577  if (mype .eq. pe_io) then
578  var_name = 'Scalar_xxx'
579  write (var_name(8:10),'(i3.3)') i
580  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, var_name, varid_scalar(i)) )
581  end if
582  call readscatterres_netcdf(scalar(:,i) , array2, array1, ncid , varid_scalar(i))
583  end do
584 
585  do i = 1,nrr
586  if (mype .eq. pe_io) then
587  var_name = 'RRates_xxx'
588  write (var_name(8:10),'(i3.3)') i
589  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, var_name, varid_rr(i)) )
590  end if
591  call readscatterres_netcdf(reactionrates(:,i) , array2, array1, ncid , varid_rr(i))
592  end do
593 
594  if (k_epsilon) then
595  if (mype .eq. pe_io) then
596  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, 'k_turb_g', varid_kturbg) )
597  end if
598  call readscatterres_netcdf(k_turb_g , array2, array1, ncid , varid_kturbg)
599 
600  if (mype .eq. pe_io) then
601  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, 'e_turb_g', varid_eturbg) )
602  end if
603  call readscatterres_netcdf(e_turb_g , array2, array1, ncid , varid_eturbg)
604  end if
605 
606  if (mype .eq. pe_io) then
607  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, 'gamaRG', varid_gamarg) )
608  end if
609  call readscatterres_netcdf(gama_rg , array2, array1, ncid , varid_gamarg)
610 
611  if (mype .eq. pe_io) then
612  call mfix_check_netcdf( mfix_nf90_inq_varid(ncid, 'TRG', varid_trg) )
613  end if
614  call readscatterres_netcdf(t_rg , array2, array1, ncid , varid_trg)
615 
616  ! Close the file. This frees up any internal netCDF resources
617  ! associated with the file, and flushes any buffers.
618 ! call MPI_barrier(MPI_COMM_WORLD,mpierr)
619  if (mype .eq. pe_io) then
620  call mfix_check_netcdf( mfix_nf90_close(ncid) )
621  end if
622 ! call MPI_barrier(MPI_COMM_WORLD,mpierr)
623 
624  ! stop
625 
626 
627  return
628 
629  end subroutine read_res1_netcdf
630 
631 
632 
633 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvC
634 ! C
635 ! Module name: PATCH_AFTER_RESTART C
636 ! Purpose: Patch new fluid cells after a restart C
637 ! This could occur when restarting with a different C
638 ! grid partition when the RESTART file was generated C
639 ! prior to the Dec. 4th 2014 bug fix. C
640 ! C
641 ! Author: Jeff Dietiker Date: 14-APR-15 C
642 ! C
643 ! Revision Number: C
644 ! Purpose: C
645 ! Author: Date: dd-mmm-yy C
646 ! Reviewer: Date: dd-mmm-yy C
647 ! C
648 ! Literature/Document References: C
649 ! C
650 ! C
651 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^C
652 
653  SUBROUTINE patch_after_restart
654 !...Translated by Pacific-Sierra Research VAST-90 2.06G5 12:17:31 12/09/98
655 !...Switches: -xf
656 !
657 !-----------------------------------------------
658 ! M o d u l e s
659 !-----------------------------------------------
660  USE param
661  USE param1
662  USE fldvar
663  USE geometry
664  USE physprop
665  USE run
666  USE rxns
667  USE scalars
668  USE funits
669  USE energy
670  USE compar
671  USE cdist
672  USE mpi_utility
673  USE sendrecv
674  USE cutcell
675  use functions
676 
677  IMPLICIT NONE
678 !-----------------------------------------------
679 ! G l o b a l P a r a m e t e r s
680 !-----------------------------------------------
681 !-----------------------------------------------
682 ! L o c a l P a r a m e t e r s
683 !-----------------------------------------------
684 !-----------------------------------------------
685 ! L o c a l V a r i a b l e s
686 !-----------------------------------------------
687 !
688 !
689 !-----------------------------------------------
690 ! Local variables
691 !-----------------------------------------------
692 ! indices
693  INTEGER :: I,J,K, IJK, IJKNB
694  INTEGER :: M,NN
695  INTEGER :: NB
696  INTEGER, DIMENSION(6) :: NBCELL
697  LOGICAL :: NB_FOUND
698 
699 
700 !-----------------------------------------------
701 
702  DO ijk = ijkstart3, ijkend3
703 
704  IF (fluid_at(ijk).AND.ep_g(ijk)==undefined) THEN
705 
706 ! Detects new fluid cells that used to be blocked cells with undefined
707 ! values. When a fluid cell has undefined void fraction, this means all
708 ! variables need to be patched. Typically, this will be a fairly small
709 ! cut cell that was flagged as BLOCKED cell with a different partition.
710 ! If a valid fluid cell is found next to this undefined cell, all field
711 ! variables will be copied over. If no valid fluid cell is found, the
712 ! code will continue and will likely stop during the check_data_30
713 ! (zero species mass fractions will yield a zero specific heat).
714  i = i_of(ijk)
715  j = j_of(ijk)
716  k = k_of(ijk)
717 
718  nbcell(1) = im_of(ijk)
719  nbcell(2) = jm_of(ijk)
720  nbcell(3) = km_of(ijk)
721  nbcell(4) = ip_of(ijk)
722  nbcell(5) = jp_of(ijk)
723  nbcell(6) = kp_of(ijk)
724 
725  nb_found = .false.
726 
727  DO nb = 1,6
728 
729  ijknb = nbcell(nb)
730 
731  IF(fluid_at(ijknb).AND.ep_g(ijknb)/=undefined) THEN
732  nb_found = .true.
733  WRITE (*, 1010) mype, i,j,k
734 
735  ep_g(ijk) = ep_g(ijknb)
736  p_g(ijk) = p_g(ijknb)
737  p_star(ijk) = p_star(ijknb)
738  ro_g(ijk) = ro_g(ijknb)
739  rop_g(ijk) = rop_g(ijknb)
740  t_g(ijk) = t_g(ijknb)
741 
742  t_s(ijk,1:mmax) = t_s(ijknb,1:mmax)
743 
744  u_s(ijk,1:mmax) = u_s(ijknb,1:mmax)
745  v_s(ijk,1:mmax) = v_s(ijknb,1:mmax)
746  w_s(ijk,1:mmax) = w_s(ijknb,1:mmax)
747 
748  x_g(ijk,1:nmax(0)) = x_g(ijknb,1:nmax(0))
749 
750  u_g(ijk) = u_g(ijknb)
751  v_g(ijk) = v_g(ijknb)
752  w_g(ijk) = w_g(ijknb)
753 
754  rop_s(ijk,1:mmax) = rop_s(ijknb,1:mmax)
755 
756  IF(any(solve_ros)) ro_s(ijk,1:mmax) = ro_s(ijknb,1:mmax)
757 
758 
759  theta_m(ijk,1:mmax) = theta_m(ijknb,1:mmax)
760 
761  DO m = 1,mmax
762  DO nn = 1, nmax(m)
763  x_s(ijk,m,nn)= x_s(ijknb,m,nn)
764  ENDDO
765  ENDDO
766 
767 
768  DO nn = 1, nscalar
769  scalar(ijk,nn) = scalar(ijknb,nn)
770  END DO
771 
772  gama_rg(ijk) = gama_rg(ijknb)
773  t_rg(ijk) = t_rg(ijknb)
774 
775  gama_rs(ijk,1:mmax) = gama_rs(ijknb,1:mmax)
776  t_rs(ijk,1:mmax) = t_rs(ijknb,1:mmax)
777 
778 
779  DO nn = 1, nrr
780  reactionrates(ijk,nn) = reactionrates(ijknb,nn)
781  END DO
782 
783  IF (k_epsilon) THEN
784  k_turb_g(ijk) = k_turb_g(ijknb)
785  e_turb_g(ijk) = e_turb_g(ijknb)
786  ENDIF
787 
788  EXIT ! Exit as soon as first valid neighbor cell is found
789  ENDIF ! NB is a fluid cell
790 
791  ENDDO ! NB Loop
792 
793  IF(.NOT.nb_found) WRITE (*, 1020) mype, i,j,k ! NO FLUID CELL AMONG NEIGBHORS
794 
795  ENDIF ! New fuid cell
796 
797  ENDDO ! IJK loop
798 
799 1010 FORMAT(1x,'PATCHING NEW FLUID CELL UPON RESTART: MyPE,I,J,K =' ,i6,i6,i6,i6)
800 1020 FORMAT(1x,'UNABLE TO PATCH NEW FLUID CELL UPON RESTART: MyPE,I,J,K =' ,i6,i6,i6,i6)
801  END SUBROUTINE patch_after_restart
double precision, dimension(:,:), allocatable v_s
Definition: fldvar_mod.f:105
subroutine patch_after_restart
Definition: read_res1.f:654
integer ijkend3
Definition: compar_mod.f:80
double precision, dimension(:), allocatable ep_g
Definition: fldvar_mod.f:15
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
double precision dt_fac
Definition: run_mod.f:226
double precision, parameter one
Definition: param1_mod.f:29
integer dimension_3
Definition: param_mod.f:11
subroutine readscatterres_netcdf(VAR, array2, array1, ncid, varid
Definition: read_res1.f:291
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
subroutine mfix_ncvinq(ncid, varid, varnam, vartyp, nvdims, vdims, nvatts, rcode)
double precision, dimension(:), allocatable t_g
Definition: fldvar_mod.f:63
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
double precision, dimension(:), allocatable p_g
Definition: fldvar_mod.f:26
subroutine in_bin_512(IUNIT, ARRAY, N, NEXT_REC)
subroutine read_res1
Definition: read_res1.f:26
double precision, parameter undefined
Definition: param1_mod.f:18
double precision, dimension(:,:), allocatable u_s
Definition: fldvar_mod.f:93
integer nf90_nowrite
logical bstart_with_one_res
Definition: cdist_mod.f:5
integer ijkmax3
Definition: geometry_mod.f:82
integer pe_io
Definition: compar_mod.f:30
subroutine readscatterres(VAR, array2, array1, init, NEXT_REC)
Definition: read_res1.f:250
subroutine convert_from_io_dp(arr_io, arr_internal, nn)
integer mmax
Definition: physprop_mod.f:19
double precision, dimension(:,:,:), allocatable x_s
Definition: fldvar_mod.f:78
integer function mfix_nf90_inquire_dimension(ncid, dimid, name, len)
double precision, dimension(:,:), allocatable t_s
Definition: fldvar_mod.f:66
double precision, dimension(:,:), allocatable x_g
Definition: fldvar_mod.f:75
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
double precision, dimension(:,:), allocatable theta_m
Definition: fldvar_mod.f:149
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
integer nsteprst
Definition: run_mod.f:64
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
integer ijkstart3
Definition: compar_mod.f:80
logical function mfix_usingnetcdf()
double precision, dimension(:), allocatable u_g
Definition: fldvar_mod.f:87
integer nscalar
Definition: scalars_mod.f:7
integer function mfix_nf90_inq_dimid(ncid, name, dimid)
integer, parameter unit_res
Definition: funits_mod.f:27
subroutine read_res1_netcdf
Definition: read_res1.f:338
double precision, dimension(:,:), allocatable rop_s
Definition: fldvar_mod.f:51
double precision, dimension(:,:), allocatable gama_rs
Definition: energy_mod.f:18
double precision time
Definition: run_mod.f:45
double precision, dimension(:), allocatable e_turb_g
Definition: fldvar_mod.f:162
double precision, dimension(:), allocatable ro_g
Definition: fldvar_mod.f:32
integer function mfix_nf90_inquire(ncid, nDimensions, nVariables, nAttributes, unlimitedDimId, formatNum)
double precision, dimension(:), allocatable rop_g
Definition: fldvar_mod.f:38
double precision, dimension(:), allocatable x
Definition: geometry_mod.f:129
integer function mfix_nf90_inq_varid(ncid, name, varid)
double precision, parameter zero
Definition: param1_mod.f:27
integer function mfix_nf90_open(path, mode, ncid, chunksize)