File: /nfs/home/0/users/jenkins/mfix.git/model/write_spx1.f

1     !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvC
2     !                                                                      C
3     !  Module name: WRITE_SPX1                                             C
4     !  Purpose: write out the time-dependent restart records (REAL)        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, U_g           C
17     !                        V_g, W_g, U_s, V_s, W_s, ROP_s, T_g, T_s      C
18     !                        IJKMAX2, MMAX                                 C
19     !  Variables modified: None                                            C
20     !                                                                      C
21     !  Local variables:  LC, N, NEXT_REC, NUM_REC                          C
22     !                                                                      C
23     !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^C
24     !
25           SUBROUTINE WRITE_SPX1(L, unit_add)
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 funits
39           USE scalars
40           USE output
41           USE rxns
42           USE cdist
43           USE compar           !//
44           USE mpi_utility      !//
45           USE cutcell
46           use discretelement, only: PRINT_DES_DATA
47           use discretelement, only: DISCRETE_ELEMENT
48           use discretelement, only: PARTICLES, NFACTOR
49     
50     !//       USE tmp_array
51           IMPLICIT NONE
52     !-----------------------------------------------
53     !   G l o b a l   P a r a m e t e r s
54     !-----------------------------------------------
55     !-----------------------------------------------
56     !   D u m m y   A r g u m e n t s
57     !-----------------------------------------------
58     !
59     !             flag whether to write a particular SPx file
60           INTEGER L
61     
62     !              offset for use in post_mfix
63           INTEGER  unit_add
64     !-----------------------------------------------
65     !   L o c a l   P a r a m e t e r s
66     !-----------------------------------------------
67     !-----------------------------------------------
68     !   L o c a l   V a r i a b l e s
69     !-----------------------------------------------
70     !
71     ! local variables
72     !
73     !//
74           double precision, allocatable :: array1(:)     !//
75           double precision, allocatable :: array2(:)     !//
76     
77     !             loop counters
78           INTEGER LC, N
79     !
80     !             Pointer to the next record
81           INTEGER NEXT_REC
82     !
83     !              Number of records written each time step
84           INTEGER  NUM_REC
85     
86           INTEGER  uspx   ! UNIT_SPX + offset from post_mfix
87           CHARACTER(LEN=50), DIMENSION(1) :: LINE   !error message
88           double precision, dimension(:), allocatable :: TMP_VAR
89     
90           allocate(TMP_VAR(DIMENSION_3))
91     
92     !-----------------------------------------------
93           uspx = UNIT_SPX + unit_add
94     
95     !
96           if (myPE .eq.PE_IO) then
97              allocate (array1(ijkmax2))   !//
98              allocate (array2(ijkmax3))   !//
99           else
100              allocate (array1(1))   !//
101              allocate (array2(1))   !//
102           end if
103     !
104     ! ".SP1" FILE         EP_g    [ ROP_g, RO_g  must be calculated ...
105     !                                        not written out ]
106     !
107           SELECT CASE (L)
108           CASE (1)
109              if (myPE.eq.PE_IO.or.bDist_IO) then
110                 READ (uspx + L, REC=3) NEXT_REC, NUM_REC
111                 NUM_REC = NEXT_REC
112                 WRITE (uspx + L, REC=NEXT_REC) REAL(TIME), NSTEP
113                 NEXT_REC = NEXT_REC + 1
114              end if
115              if (bDist_IO) then
116                 IF(RE_INDEXING) THEN
117                    CALL UNSHIFT_DP_ARRAY(EP_g,TMP_VAR)
118                    call OUT_BIN_R(uspx+L,TMP_VAR,size(TMP_VAR),NEXT_REC)
119                 ELSE
120                    call OUT_BIN_R(uspx+L,EP_g,size(EP_g),NEXT_REC)
121                 ENDIF
122     !           call OUT_BIN_R(uspx+L,EP_g,size(EP_g),NEXT_REC)
123              else
124                 call gatherWriteSpx (EP_g,array2, array1, uspx+L, NEXT_REC)   !//
125              end if
126              if (myPE .eq. PE_IO.or.bDist_IO) then
127                 NUM_REC = NEXT_REC - NUM_REC
128                 WRITE (uspx + L, REC=3) NEXT_REC, NUM_REC
129                 if(unit_add == 0) CALL FLUSH_bin (uspx + L)
130              end if
131     
132     
133     ! remove this redundant call here to write_des_data in case of new
134     ! coupled runs that contain at least a particle and involve some initial
135     ! settling of the system.
136     ! the call made in des_time_march is a better call for capturing the
137     ! initial state of such a des continuum coupled system
138              IF(DISCRETE_ELEMENT.AND.PRINT_DES_DATA .AND. &
139                 .NOT.(TRIM(RUN_TYPE)=='NEW' .AND. PARTICLES /=0 .AND. &
140                       NFACTOR >0 .AND. TIME == ZERO)) THEN
141                    CALL WRITE_DES_DATA
142              ENDIF
143     
144     
145     
146     !        call MPI_Barrier(MPI_COMM_WORLD,mpierr)  !//PAR_I/O enforce barrier here
147     !
148     ! ".SP2" FILE         P_g , P_star
149     !
150           CASE (2)
151              if (myPE.eq.PE_IO.or.bDist_IO) then
152                 READ (uspx + L, REC=3) NEXT_REC, NUM_REC
153                 NUM_REC = NEXT_REC
154                 WRITE (uspx + L, REC=NEXT_REC) REAL(TIME), NSTEP
155                 NEXT_REC = NEXT_REC + 1
156              end if
157              if (bDist_IO) then
158                 IF(RE_INDEXING) THEN
159                    CALL UNSHIFT_DP_ARRAY(P_g,TMP_VAR)
160                    call OUT_BIN_R(uspx+L,TMP_VAR,size(TMP_VAR),NEXT_REC)
161                    CALL UNSHIFT_DP_ARRAY(P_star,TMP_VAR)
162                    call OUT_BIN_R(uspx+L,TMP_VAR,size(TMP_VAR),NEXT_REC)
163                 ELSE
164                   call OUT_BIN_R(uspx+L,P_g,size(P_g),NEXT_REC)
165                   call OUT_BIN_R(uspx+L,P_star,size(P_star),NEXT_REC)
166                 ENDIF
167     !           call OUT_BIN_R(uspx+L,P_g,size(P_g),NEXT_REC)
168     !           call OUT_BIN_R(uspx+L,P_star,size(P_star),NEXT_REC)
169              else
170                call gatherWriteSpx (P_g,array2, array1, uspx+L, NEXT_REC)   !//
171                call gatherWriteSpx (P_star,array2, array1, uspx+L, NEXT_REC)   !//
172              end if
173              if (myPE.eq.PE_IO.or.bDist_IO) then
174                 NUM_REC = NEXT_REC - NUM_REC
175                 WRITE (uspx + L, REC=3) NEXT_REC, NUM_REC
176                 if(unit_add == 0) CALL FLUSH_bin (uspx + L)
177              end if
178     !        call MPI_Barrier(MPI_COMM_WORLD,mpierr)  !//PAR_I/O enforce barrier here
179     !
180     ! ".SP3" FILE         U_g , V_g , W_g
181     !
182           CASE (3)
183              if (myPE.eq.PE_IO.or.bDist_IO) then
184                 READ (uspx + L, REC=3) NEXT_REC, NUM_REC
185                 NUM_REC = NEXT_REC
186                 WRITE (uspx + L, REC=NEXT_REC) REAL(TIME), NSTEP
187                 NEXT_REC = NEXT_REC + 1
188              end if
189              if (bDist_IO) then
190                 IF(RE_INDEXING) THEN
191                    CALL UNSHIFT_DP_ARRAY(U_g,TMP_VAR)
192                    call OUT_BIN_R(uspx+L,TMP_VAR,size(TMP_VAR),NEXT_REC)
193                    CALL UNSHIFT_DP_ARRAY(V_g,TMP_VAR)
194                    call OUT_BIN_R(uspx+L,TMP_VAR,size(TMP_VAR),NEXT_REC)
195                    CALL UNSHIFT_DP_ARRAY(W_g,TMP_VAR)
196                    call OUT_BIN_R(uspx+L,TMP_VAR,size(TMP_VAR),NEXT_REC)
197                 ELSE
198                    call OUT_BIN_R(uspx+L,U_g,size(U_g),NEXT_REC)
199                    call OUT_BIN_R(uspx+L,V_g,size(V_g),NEXT_REC)
200                    call OUT_BIN_R(uspx+L,W_g,size(W_g),NEXT_REC)
201                 ENDIF
202     !           call OUT_BIN_R(uspx+L,U_g,size(U_g),NEXT_REC)
203     !           call OUT_BIN_R(uspx+L,V_g,size(V_g),NEXT_REC)
204     !           call OUT_BIN_R(uspx+L,W_g,size(W_g),NEXT_REC)
205              else
206                call gatherWriteSpx (U_g,array2, array1, uspx+L, NEXT_REC)   !//
207                call gatherWriteSpx (V_g,array2, array1, uspx+L, NEXT_REC)   !//
208                call gatherWriteSpx (W_g,array2, array1, uspx+L, NEXT_REC)   !//
209              end if
210              if (myPE.eq.PE_IO.or.bDist_IO) then
211                 NUM_REC = NEXT_REC - NUM_REC
212                 WRITE (uspx + L, REC=3) NEXT_REC, NUM_REC
213                 if(unit_add == 0) CALL FLUSH_bin (uspx + L)
214              end if
215     !        call MPI_Barrier(MPI_COMM_WORLD,mpierr)  !//PAR_I/O enforce barrier here
216     !
217     ! ".SP4" FILE         U_s , V_s , W_s
218     !
219           CASE (4)
220              if (myPE.eq.PE_IO.or.bDist_IO) then
221                 READ (uspx + L, REC=3) NEXT_REC, NUM_REC
222                 NUM_REC = NEXT_REC
223                 WRITE (uspx + L, REC=NEXT_REC) REAL(TIME), NSTEP
224                 NEXT_REC = NEXT_REC + 1
225              end if
226              if (bDist_IO) then
227                 DO LC = 1, MMAX
228                    IF(RE_INDEXING) THEN
229                       CALL UNSHIFT_DP_ARRAY(U_s(:,LC),TMP_VAR)
230                       call OUT_BIN_R(uspx+L,TMP_VAR,size(TMP_VAR),NEXT_REC)
231                       CALL UNSHIFT_DP_ARRAY(V_s(:,LC),TMP_VAR)
232                       call OUT_BIN_R(uspx+L,TMP_VAR,size(TMP_VAR),NEXT_REC)
233                       CALL UNSHIFT_DP_ARRAY(W_s(:,LC),TMP_VAR)
234                       call OUT_BIN_R(uspx+L,TMP_VAR,size(TMP_VAR),NEXT_REC)
235                    ELSE
236                       call OUT_BIN_R(uspx+L,U_s(:,LC),Size(U_s(:,LC)),NEXT_REC)
237                       call OUT_BIN_R(uspx+L,V_s(:,LC),Size(V_s(:,LC)),NEXT_REC)
238                       call OUT_BIN_R(uspx+L,W_s(:,LC),Size(W_s(:,LC)),NEXT_REC)
239                    ENDIF
240                 ENDDO
241     !        DO LC = 1, MMAX
242     !          call OUT_BIN_R(uspx+L,U_s(:,LC),Size(U_s(:,LC)),NEXT_REC)
243     !          call OUT_BIN_R(uspx+L,V_s(:,LC),Size(V_s(:,LC)),NEXT_REC)
244     !          call OUT_BIN_R(uspx+L,W_s(:,LC),Size(W_s(:,LC)),NEXT_REC)
245     !        END DO
246              else
247              DO LC = 1, MMAX
248                 call gatherWriteSpx (U_s(:,LC),array2, array1, uspx+L, NEXT_REC)
249                 call gatherWriteSpx (V_s(:,LC),array2, array1, uspx+L, NEXT_REC)
250                 call gatherWriteSpx (W_s(:,LC),array2, array1, uspx+L, NEXT_REC)
251              END DO
252              end if
253              if (myPE.eq.PE_IO.or.bDist_IO) then
254                 NUM_REC = NEXT_REC - NUM_REC
255                 WRITE (uspx + L, REC=3) NEXT_REC, NUM_REC
256                 if(unit_add == 0) CALL FLUSH_bin (uspx + L)
257              end if
258     !        call MPI_Barrier(MPI_COMM_WORLD,mpierr)  !//PAR_I/O enforce barrier here
259     !
260     ! ".SP5" FILE         ROP_s
261     !
262           CASE (5)
263              if (myPE.eq.PE_IO.or.bDist_IO) then
264                 READ (uspx + L, REC=3) NEXT_REC, NUM_REC
265                 NUM_REC = NEXT_REC
266                 WRITE (uspx + L, REC=NEXT_REC) REAL(TIME), NSTEP
267                 NEXT_REC = NEXT_REC + 1
268              end if
269              if (bDist_IO) then
270                 DO LC = 1, MMAX
271                    IF(RE_INDEXING) THEN
272                       CALL UNSHIFT_DP_ARRAY(ROP_s(:,LC),TMP_VAR)
273                       call OUT_BIN_R(uspx+L,TMP_VAR,size(TMP_VAR),NEXT_REC)
274                    ELSE
275                       call OUT_BIN_R(uspx+L,ROP_s(:,LC),size(ROP_s(:,LC)), NEXT_REC)
276                    ENDIF
277                 ENDDO
278     !          DO LC = 1, MMAX
279     !            call OUT_BIN_R(uspx+L,ROP_s(:,LC),size(ROP_s(:,LC)), NEXT_REC)
280     !         END DO
281              else
282              DO LC = 1, MMAX
283                 call gatherWriteSpx (ROP_s(:,LC),array2, array1, uspx+L, NEXT_REC)
284              END DO
285              end if
286              if (myPE.eq.PE_IO.or.bDist_IO) then
287                 NUM_REC = NEXT_REC - NUM_REC
288                 WRITE (uspx + L, REC=3) NEXT_REC, NUM_REC
289                 if(unit_add == 0) CALL FLUSH_bin (uspx + L)
290              end if
291     !        call MPI_Barrier(MPI_COMM_WORLD,mpierr)  !//PAR_I/O enforce barrier here
292     !
293     ! ".SP6" FILE         T_g  , T_s
294     !
295           CASE (6)
296              if (myPE.eq.PE_IO.or.bDist_IO) then
297                 READ (uspx + L, REC=3) NEXT_REC, NUM_REC
298                 NUM_REC = NEXT_REC
299                 WRITE (uspx + L, REC=NEXT_REC) REAL(TIME), NSTEP
300                 NEXT_REC = NEXT_REC + 1
301              end if
302              if (bDist_IO) then
303                 IF(RE_INDEXING) THEN
304                    CALL UNSHIFT_DP_ARRAY(T_g,TMP_VAR)
305                    call OUT_BIN_R(uspx+L,TMP_VAR,size(TMP_VAR),NEXT_REC)
306                    DO LC = 1, MMAX
307                       CALL UNSHIFT_DP_ARRAY(T_s(:,LC),TMP_VAR)
308                       call OUT_BIN_R(uspx+L,TMP_VAR,size(TMP_VAR), NEXT_REC)
309                    END DO
310                 ELSE
311                    call OUT_BIN_R(uspx+L,T_g,size(T_g), NEXT_REC)
312                    DO LC = 1, MMAX
313                       call OUT_BIN_R(uspx+L,T_s(:,LC),size(T_s(:,LC)), NEXT_REC)
314                    END DO
315                 ENDIF
316     !           call OUT_BIN_R(uspx+L,T_g,size(T_g), NEXT_REC)
317     !          DO LC = 1, MMAX
318     !            call OUT_BIN_R(uspx+L,T_s(:,LC),size(T_s(:,LC)), NEXT_REC)
319     !          END DO
320              else
321              call gatherWriteSpx (T_g,array2, array1, uspx+L, NEXT_REC)   !//
322               DO LC = 1, MMAX
323                 call gatherWriteSpx (T_s(:,LC),array2, array1, uspx+L, NEXT_REC)
324               END DO
325              end if
326              if (myPE.eq.PE_IO.or.bDist_IO) then
327                 NUM_REC = NEXT_REC - NUM_REC
328                 WRITE (uspx + L, REC=3) NEXT_REC, NUM_REC
329                 if(unit_add == 0) CALL FLUSH_bin (uspx + L)
330              end if
331     !
332     ! ".SP7" FILE         X_g, X_s
333     !
334           CASE (7)
335              if (myPE.eq.PE_IO.or.bDist_IO) then
336                 READ (uspx + L, REC=3) NEXT_REC, NUM_REC
337                 NUM_REC = NEXT_REC
338                 WRITE (uspx + L, REC=NEXT_REC) REAL(TIME), NSTEP
339                 NEXT_REC = NEXT_REC + 1
340              end if
341              if (bDist_IO) then
342                 IF(RE_INDEXING) THEN
343                    DO N = 1, NMAX(0)
344                       CALL UNSHIFT_DP_ARRAY(X_G(:,N),TMP_VAR)
345                       call OUT_BIN_R(uspx+L,TMP_VAR,size(TMP_VAR),NEXT_REC)
346                    END DO
347                    DO LC = 1, MMAX
348                       DO N = 1, NMAX(LC)
349                          CALL UNSHIFT_DP_ARRAY(X_s(:,LC,N),TMP_VAR)
350                          call OUT_BIN_R(uspx+L,TMP_VAR,size(TMP_VAR), NEXT_REC)
351                       ENDDO
352                    END DO
353                 ELSE
354                    DO N = 1, NMAX(0)
355                       call OUT_BIN_R(uspx+L,X_G(:,N),size(X_G(:,N)), NEXT_REC)
356                    END DO
357                    DO LC = 1, MMAX
358                       DO N = 1, NMAX(LC)
359                          call OUT_BIN_R(uspx+L,X_s(:,LC,N),size(X_s(:,LC,N)), NEXT_REC)
360                       END DO
361                    END DO
362                 ENDIF
363     
364     !           DO N = 1, NMAX(0)
365     !             call OUT_BIN_R(uspx+L,X_G(:,N),size(X_G(:,N)), NEXT_REC)
366     !           END DO
367     !           DO LC = 1, MMAX
368     !            DO N = 1, NMAX(LC)
369     !               call OUT_BIN_R(uspx+L,X_s(:,LC,N),size(X_s(:,LC,N)), NEXT_REC)
370     !            END DO
371     !           END DO
372     
373              else
374               DO N = 1, NMAX(0)
375                 call gatherWriteSpx (X_G(:,N),array2, array1, uspx+L, NEXT_REC)
376               END DO
377               DO LC = 1, MMAX
378                 DO N = 1, NMAX(LC)
379                    call gatherWriteSpx (X_s(:,LC,N),array2, array1, uspx+L, NEXT_REC)
380                 END DO
381               END DO
382              end if
383              if (myPE.eq.PE_IO.or.bDist_IO) then
384                 NUM_REC = NEXT_REC - NUM_REC
385                 WRITE (uspx + L, REC=3) NEXT_REC, NUM_REC
386                 if(unit_add == 0) CALL FLUSH_bin (uspx + L)
387              end if
388     !        call MPI_Barrier(MPI_COMM_WORLD,mpierr)  !//PAR_I/O enforce barrier here
389     !
390     ! ".SP8" FILE         THETA_m
391     !
392           CASE (8)
393              if (myPE.eq.PE_IO.or.bDist_IO) then
394                 READ (uspx + L, REC=3) NEXT_REC, NUM_REC
395                 NUM_REC = NEXT_REC
396                 WRITE (uspx + L, REC=NEXT_REC) REAL(TIME), NSTEP
397                 NEXT_REC = NEXT_REC + 1
398              end if
399              if (bDist_IO) then
400                 IF(RE_INDEXING) THEN
401                    DO LC = 1, MMAX
402                       CALL UNSHIFT_DP_ARRAY(THETA_m(:,LC),TMP_VAR)
403                       call OUT_BIN_R(uspx+L,TMP_VAR,size(TMP_VAR),NEXT_REC)
404                    ENDDO
405                 ELSE
406                    DO LC = 1, MMAX
407                       call OUT_BIN_R(uspx+L,THETA_m(:,LC),size(THETA_m(:,LC)), NEXT_REC)
408                    END DO
409                 ENDIF
410     !          DO LC = 1, MMAX
411     !            call OUT_BIN_R(uspx+L,THETA_m(:,LC),size(THETA_m(:,LC)), NEXT_REC)
412     !         END DO
413              else
414               DO LC = 1, MMAX
415                 call gatherWriteSpx (THETA_m(:,LC),array2, array1, uspx+L, NEXT_REC)
416               END DO
417              end if
418              if (myPE.eq.PE_IO .or. bDist_IO) then
419                 NUM_REC = NEXT_REC - NUM_REC
420                 WRITE (uspx + L, REC=3) NEXT_REC, NUM_REC
421                 if(unit_add == 0) CALL FLUSH_bin (uspx + L)
422              end if
423     !        call MPI_Barrier(MPI_COMM_WORLD,mpierr)  !//PAR_I/O enforce barrier here
424     !
425     ! ".SP9" FILE         Scalar
426     !
427           CASE (9)
428              if (myPE.eq.PE_IO.or.bDist_IO) then
429                 READ (uspx + L, REC=3) NEXT_REC, NUM_REC
430                 NUM_REC = NEXT_REC
431                 WRITE (uspx + L, REC=NEXT_REC) REAL(TIME), NSTEP
432                 NEXT_REC = NEXT_REC + 1
433              end if
434              if (bDist_IO) then
435                 IF(RE_INDEXING) THEN
436                    DO LC = 1, Nscalar
437                       CALL UNSHIFT_DP_ARRAY(Scalar(:,LC),TMP_VAR)
438                       call OUT_BIN_R(uspx+L,TMP_VAR,size(TMP_VAR),NEXT_REC)
439                    ENDDO
440                 ELSE
441                    DO LC = 1, Nscalar
442                      call OUT_BIN_R(uspx+L,Scalar(:,LC),size(Scalar(:,LC)), NEXT_REC)
443                    END DO
444                 ENDIF
445     !          DO LC = 1, Nscalar
446     !            call OUT_BIN_R(uspx+L,Scalar(:,LC),size(Scalar(:,LC)), NEXT_REC)
447     !          END DO
448              else
449               DO LC = 1, Nscalar
450                 call gatherWriteSpx (Scalar(:,LC),array2, array1, uspx+L, NEXT_REC)
451               END DO
452              end if
453              if (myPE.eq.PE_IO.or.bDist_IO) then
454                 NUM_REC = NEXT_REC - NUM_REC
455                 WRITE (uspx + L, REC=3) NEXT_REC, NUM_REC
456                 if(unit_add == 0) CALL FLUSH_bin (uspx + L)
457              end if
458     !        call MPI_Barrier(MPI_COMM_WORLD,mpierr)  !//PAR_I/O enforce barrier here
459     !
460           CASE (10)  ! Reaction rates
461     
462              if (myPE.eq.PE_IO.or.bDist_IO) then
463                 READ (uspx + L, REC=3) NEXT_REC, NUM_REC
464                 NUM_REC = NEXT_REC
465                 WRITE (uspx + L, REC=NEXT_REC) REAL(TIME), NSTEP
466                 NEXT_REC = NEXT_REC + 1
467              end if
468              if (bDist_IO) then
469                 IF(RE_INDEXING) THEN
470                    DO LC = 1, nRR
471                       CALL UNSHIFT_DP_ARRAY(ReactionRates(:,LC),TMP_VAR)
472                       call OUT_BIN_R(uspx+L,TMP_VAR,size(TMP_VAR),NEXT_REC)
473                    ENDDO
474                 ELSE
475                    DO LC = 1, nRR
476                       call OUT_BIN_R(uspx+L,ReactionRates(:,LC),size(ReactionRates(:,LC)), NEXT_REC)
477                    END DO
478                 ENDIF
479     !          DO LC = 1, nRR
480     !            call OUT_BIN_R(uspx+L,ReactionRates(:,LC),size(ReactionRates(:,LC)), NEXT_REC)
481     !          END DO
482              else
483               DO LC = 1, nRR
484                 call gatherWriteSpx (ReactionRates(:,LC),array2, array1, uspx+L, NEXT_REC)
485               END DO
486              end if
487              if (myPE.eq.PE_IO.or.bDist_IO) then
488                 NUM_REC = NEXT_REC - NUM_REC
489                 WRITE (uspx + L, REC=3) NEXT_REC, NUM_REC
490                 if(unit_add == 0) CALL FLUSH_bin (uspx + L)
491              end if
492     !
493     ! ".SP11" FILE         turbulence
494     !
495           CASE (11)
496              if (myPE.eq.PE_IO.or.bDist_IO) then
497                 READ (uspx + L, REC=3) NEXT_REC, NUM_REC
498                 NUM_REC = NEXT_REC
499                 WRITE (uspx + L, REC=NEXT_REC) REAL(TIME), NSTEP
500                 NEXT_REC = NEXT_REC + 1
501              end if
502              if (K_Epsilon) then
503                 if (bDist_IO) then
504                 IF(RE_INDEXING) THEN
505                    CALL UNSHIFT_DP_ARRAY(K_Turb_G,TMP_VAR)
506                    call OUT_BIN_R(uspx+L,TMP_VAR,size(TMP_VAR),NEXT_REC)
507                    CALL UNSHIFT_DP_ARRAY(E_Turb_G,TMP_VAR)
508                    call OUT_BIN_R(uspx+L,TMP_VAR,size(TMP_VAR),NEXT_REC)
509                 ELSE
510                    call OUT_BIN_R(uspx+L,K_Turb_G,size(K_Turb_G), NEXT_REC)
511                    call OUT_BIN_R(uspx+L,E_Turb_G,size(E_Turb_G), NEXT_REC)
512                 ENDIF
513     !             call OUT_BIN_R(uspx+L,K_Turb_G,size(K_Turb_G), NEXT_REC)
514     !             call OUT_BIN_R(uspx+L,E_Turb_G,size(E_Turb_G), NEXT_REC)
515     
516                 else
517                  call gatherWriteSpx (K_Turb_G,array2, array1, uspx+L, NEXT_REC)
518                  call gatherWriteSpx (E_Turb_G,array2, array1, uspx+L, NEXT_REC)
519                 end if
520               end if
521     
522                if (myPE.eq.PE_IO.or.bDist_IO) then
523                 NUM_REC = NEXT_REC - NUM_REC
524                 WRITE (uspx + L, REC=3) NEXT_REC, NUM_REC
525                 if(unit_add == 0) CALL FLUSH_bin (uspx + L)
526                end if
527     !        call MPI_Barrier(MPI_COMM_WORLD,mpierr)  !//PAR_I/O enforce barrier here
528     !
529     !
530           CASE DEFAULT
531                 LINE(1) = 'Unknown SPx file index'
532                 CALL WRITE_ERROR ('WRITE_SPX1', LINE, 1)
533                 CALL MFIX_EXIT(myPE)
534           END SELECT
535     
536     !//      call unlock_tmp_array
537     !
538           deallocate (array1)
539           deallocate (array2)
540           deallocate (TMP_VAR)
541     !
542           RETURN
543           END SUBROUTINE WRITE_SPX1
544     
545           subroutine gatherWriteSpx(VAR, array2, array1, uspxL, NEXT_REC)
546             USE geometry
547             USE compar           !//
548             USE mpi_utility      !//d pnicol : for gatherWriteSpx
549             USE sendrecv         !//d pnicol : for gatherWriteSpx
550             USE cutcell
551             USE in_binary_512
552             IMPLICIT NONE
553             integer uspxL, NEXT_REC
554             double precision, dimension(ijkmax2) :: array1
555             double precision, dimension(ijkmax3) :: array2
556             double precision, dimension(DIMENSION_3) :: VAR
557             double precision, dimension(:), allocatable :: TMP_VAR
558     
559             allocate(TMP_VAR(DIMENSION_3))
560     
561     !       call MPI_Barrier(MPI_COMM_WORLD,mpierr)  !//PAR_I/O enforce barrier here
562             IF(RE_INDEXING) THEN
563                TMP_VAR = UNDEFINED
564                CALL UNSHIFT_DP_ARRAY(VAR,TMP_VAR)
565                CALL gather (TMP_VAR,array2,root)
566             ELSE
567                CALL gather (VAR,array2,root)
568             ENDIF
569     !        call gather (VAR,array2,root)  !//d pnicol
570     
571     !       call MPI_Barrier(MPI_COMM_WORLD,mpierr)  !//PAR_I/O enforce barrier here
572             if (myPE.eq.PE_IO) then
573                call convert_to_io_dp(array2,array1,ijkmax2)
574                CALL OUT_BIN_R (uspxL, array1, IJKMAX2, NEXT_REC)
575             end if
576     
577             deallocate(TMP_VAR)
578     
579           End subroutine gatherWriteSpx
580     
581     
582     
583           subroutine gatherWriteSpx_netcdf(VAR, arr1, arr2 , arr4d, ncid, varid , &
584                     nx,ny,nz,ijkmax2_use , ijkmax3_use)
585     
586     
587          USE geometry
588          USE compar           !//
589          USE mpi_utility      !//d pnicol : for gatherWriteSpx
590          USE sendrecv         !//d pnicol : for gatherWriteSpx
591          USE MFIX_netcdf
592          USE in_binary_512
593     
594          IMPLICIT NONE
595     
596          integer          :: ncid , varid , nx,ny,nz , ijkmax2_use , ijkmax3_use
597          integer          :: ii , jj , kk , ijk
598     
599          double precision ::  arr1(ijkmax2_use)
600          double precision ::  arr2(ijkmax3_use)
601          double precision ::  arr4d(nx,ny,nz,1)
602          double precision ::  var(dimension_3)
603     
604          call gather(var,arr2,root)
605          if (myPE .eq. PE_IO) then
606             call convert_to_io_dp(arr2,arr1,ijkmax2_use)
607     
608             ijk = 0
609             do kk = 1,nz
610                do jj = 1,ny
611                   do ii = 1,nx
612                      ijk = ijk + 1
613                      arr4d(ii,jj,kk,1) = arr1(ijk)
614                   end do
615                end do
616             end do
617     
618             call MFIX_check_netcdf( MFIX_nf90_put_var(ncid,varid,arr4d) )
619     
620          end if
621     
622     
623         End subroutine gatherWriteSpx_netcdf
624     
625     
626     
627           subroutine gatherWriteSpx_netcdf_int(VAR, arr1, arr2 , arr4d, ncid, &
628                     varid , nx,ny,nz,ijkmax2_use , ijkmax3_use)
629     
630     
631          USE geometry
632          USE compar           !//
633          USE mpi_utility      !//d pnicol : for gatherWriteSpx
634          USE sendrecv         !//d pnicol : for gatherWriteSpx
635          USE MFIX_netcdf
636          USE in_binary_512i
637     
638          IMPLICIT NONE
639     
640          integer          :: ncid , varid , nx,ny,nz , ijkmax2_use , ijkmax3_use
641          integer          :: ii , jj , kk , ijk
642     
643          integer ::  arr1(ijkmax2_use)
644          integer ::  arr2(ijkmax3_use)
645          integer :: arr4d(nx,ny,nz,1)
646          integer ::   var(dimension_3)
647     
648          call gather(var,arr2,root)
649          if (myPE .eq. PE_IO) then
650             call convert_to_io_i(arr2,arr1,ijkmax2_use)
651     
652             ijk = 0
653             do kk = 1,nz
654                do jj = 1,ny
655                   do ii = 1,nx
656                      ijk = ijk + 1
657                      arr4d(ii,jj,kk,1) = arr1(ijk)
658                   end do
659                end do
660             end do
661     
662             call MFIX_check_netcdf( MFIX_nf90_put_var(ncid,varid,arr4d) )
663     
664          end if
665     
666     
667         End subroutine gatherWriteSpx_netcdf_int
668     
669     
670     
671     
672             subroutine copy_d_to_r(darr,rarr,nx,ny,nz)
673             implicit none
674     
675             integer          :: nx , ny , nz
676             double precision :: darr(*)
677             real             :: rarr(nx,ny,*)
678             integer          :: i , j , k , ijk
679     
680     
681             ijk = 0
682     
683             do i = 1,nx
684                do j = 1,ny
685                   do k = 1,nz
686                      ijk = ijk + 1
687                      rarr(i,j,k) = darr(ijk)
688                   end do
689                end do
690             end do
691     
692             return
693             end subroutine copy_d_to_r
694     
695     
696     
697     !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
698     !
699     !                          write_mesh_netcdf
700     !
701     !
702             SUBROUTINE write_mesh_netcdf
703     
704             USE param
705             USE param1
706             USE fldvar
707             USE geometry
708             USE physprop
709             USE run
710     !       USE funits
711             USE scalars
712     !       USE output
713             USE rxns
714             USE cdist
715             USE compar
716             USE mpi_utility
717             USE MFIX_netcdf
718     !       USE tmp_array
719     
720             implicit none
721     
722             integer   :: ncid    , x_dimid , y_dimid , z_dimid  , t_dimid
723             integer   :: varid_x , varid_y , varid_z , L , dimids(4)
724             integer   :: varid_flag , coords_dimid , varid_coords , coords
725     
726             character(LEN=80) :: fname
727     
728             double precision, dimension(:) , allocatable :: xloc
729             double precision, dimension(:) , allocatable :: yloc
730             double precision, dimension(:) , allocatable :: zloc
731     
732             integer, dimension(:) , allocatable :: arr1
733             integer, dimension(:) , allocatable :: arr2
734             integer, dimension(:,:,:,:) , allocatable :: arr4d
735     
736     
737             if (.not. MFIX_usingNETCDF()) return
738     
739             if (.not. bGlobalNetcdf) return  ! no netCDF writes asked for
740     
741             if (.not. bFirst_netcdf_write) return
742     
743             if (myPE.eq.PE_IO) then
744                     allocate ( arr1(ijkmax2))
745                     allocate ( arr2(ijkmax3))
746                     allocate ( arr4d(imax2,jmax2,kmax2,1))
747                     allocate ( xloc(imax2) )
748                     allocate ( yloc(jmax2) )
749                     allocate ( zloc(kmax2) )
750             else
751                     allocate ( arr1(1))
752                     allocate ( arr2(1))
753                     allocate ( arr4d(1,1,1,1))
754                     allocate ( xloc(1) )
755                     allocate ( yloc(1) )
756                     allocate ( zloc(1) )
757             end if
758     
759             if (myPE.eq.PE_IO) then
760                     xloc(1) = -dx(1)
761                     do L = 2,imax2
762                             xloc(L) = xloc(L-1) + dx(L)
763                     end do
764     
765                     yloc(1) = -dy(1)
766                     do L = 2,jmax2
767                             yloc(L) = yloc(L-1) + dy(L)
768                     end do
769     
770                     zloc(1) = -dz(1)
771                     do L = 2,kmax2
772                             zloc(L) = zloc(L-1) + dz(L)
773                     end do
774     
775                     fname = trim(run_name) // "_MESH.nc"
776                     call MFIX_check_netcdf( MFIX_nf90_create(fname, NF90_64BIT_OFFSET, ncid) )
777     
778                     call MFIX_check_netcdf( MFIX_nf90_def_dim(ncid, "x"   , imax2   , x_dimid  ) )
779                     call MFIX_check_netcdf( MFIX_nf90_def_dim(ncid, "y"   , jmax2   , y_dimid  ) )
780                     call MFIX_check_netcdf( MFIX_nf90_def_dim(ncid, "z"   , kmax2   , z_dimid  ) )
781                     call MFIX_check_netcdf( MFIX_nf90_def_dim(ncid, "coordinates"   , 1   , coords_dimid  ) )
782                     call MFIX_check_netcdf( MFIX_nf90_def_dim(ncid, "t"   , 1       , t_dimid  ) )  ! 4
783     
784                     dimids =  (/ x_dimid , y_dimid, z_dimid , t_dimid/)
785     
786     
787                     ! The dimids array is used to pass the IDs of the dimensions of
788                     ! the variables. Note that in fortran arrays are stored in
789                     ! column-major format.
790                     call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, "x"  , NF90_DOUBLE, x_dimid, varid_x  ) )
791                     call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, "y"  , NF90_DOUBLE, y_dimid, varid_y  ) )
792                     call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, "z"  , NF90_DOUBLE, z_dimid, varid_z  ) )
793                     call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, "coordinates"  , NF90_INT, coords_dimid, varid_coords  ) )
794                     call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, "flag"  , NF90_INT, dimids, varid_flag  ) )  ! 9
795     
796     
797                     call MFIX_check_netcdf( MFIX_nf90_enddef(ncid) )
798     
799                     coords = 0
800                     if (COORDINATES .eq. 'CYLINDRICAL') coords = 1
801     
802                     call MFIX_check_netcdf( MFIX_nf90_put_var(ncid,varid_coords,coords) )
803                     call MFIX_check_netcdf( MFIX_nf90_put_var(ncid,varid_x,xloc) )
804                     call MFIX_check_netcdf( MFIX_nf90_put_var(ncid,varid_y,yloc) )
805                     call MFIX_check_netcdf( MFIX_nf90_put_var(ncid,varid_z,zloc) )
806             end if
807     
808     
809             ! needs to be called by all processes
810     
811             call gatherWriteSpx_netcdf_int(flag, arr1, arr2 , arr4d, ncid, &
812                    varid_flag , imax2,jmax2,kmax2,ijkmax2 , ijkmax3)
813     
814     
815     
816             if (myPE.eq.PE_IO) then
817                     call MFIX_check_netcdf( MFIX_nf90_close(ncid) )
818             end if
819     
820             deallocate ( arr1  )
821             deallocate ( arr2  )
822             deallocate ( arr4d )
823             deallocate ( xloc  )
824             deallocate ( yloc  )
825             deallocate ( zloc  )
826     
827             return
828             end subroutine write_mesh_netcdf
829     
830     
831     
832     
833     
834             SUBROUTINE write_netcdf(L, unit_add, the_time)
835     
836             USE param
837             USE param1
838             USE fldvar
839             USE geometry
840             USE physprop
841             USE run
842     !       USE funits
843             USE scalars
844     !       USE output
845             USE rxns
846             USE cdist
847             USE compar
848             USE mpi_utility
849             USE MFIX_netcdf
850     !       USE tmp_array
851     
852     
853             implicit none
854     
855             integer :: L , unit_add , I , n , ii
856     
857             integer   :: ncid , x_dimid , y_dimid , z_dimid
858             integer   :: t_dimid
859             integer   :: dimids(4) , varid_epg , varid_pg
860     
861             integer   :: varid_pstar  , varid_ug , varid_vg , varid_wg
862             integer   :: varid_tg , varid_x , varid_y , varid_z , varid_t
863             integer   :: varid_coords , coords_dimid , coords
864     
865             integer   :: varid_us(20) , varid_vs(20) , varid_ws(20)  !! MMAX
866             integer   :: varid_rops(20)  , varid_ts(20) !! mmax
867             integer   :: varid_thetam(20) !! mmax
868     
869             integer   :: varid_xg(20)  ! nmax(0)
870             integer   :: varid_xs(20,20)  ! mmax , MAX(nmax(1:mmax))
871     
872             integer   :: varid_scalar(20)  ! nscalar
873             integer   :: varid_rr(20)      ! nRR
874     
875             integer   :: varid_kturbg , varid_eturbg
876     
877     
878             character(LEN=80) :: fname, var_name
879             character(LEN=9) :: fname_index
880     
881             double precision, dimension(:) , allocatable :: arr1
882             double precision, dimension(:) , allocatable :: arr2
883     
884             double precision, dimension(:,:,:,:) , allocatable :: arr4d
885     
886     
887             double precision, dimension(:) , allocatable :: xloc
888             double precision, dimension(:) , allocatable :: yloc
889             double precision, dimension(:) , allocatable :: zloc
890     
891             double precision :: the_time
892             logical          :: file_exists
893     
894     
895     ! bWrite_netcdf(1)  : EP_g
896     ! bWrite_netcdf(2)  : P_g
897     ! bWrite_netcdf(3)  : P_star
898     ! bWrite_netcdf(4)  : U_g / V_g / W_g
899     ! bWrite_netcdf(5)  : U_s / V_s / W_s
900     ! bWrite_netcdf(6)  : ROP_s
901     ! bWrite_netcdf(7)  : T_g
902     ! bWrite_netcdf(8)  : T_s
903     ! bWrite_netcdf(9)  : X_g
904     ! bWrite_netcdf(10) : X_s
905     ! bWrite_netcdf(11) : Theta_m
906     ! bWrite_netCDF(12) : Scalar
907     ! bWrite_netCDF(13) : ReactionRates
908     ! bWrite_netCDF(14) : k_turb_g , e_turb_g
909     
910             if (.not. MFIX_usingNETCDF()) return
911             if (.not. bGlobalNetcdf) return
912     
913             call write_mesh_netcdf
914     
915             if (myPE.eq.PE_IO .and. .not.bDist_IO) then
916                allocate (arr1(ijkmax2))
917                allocate (arr2(ijkmax3))
918                allocate (arr4d(imax2,jmax2,kmax2,1))
919                allocate ( xloc(imax2) )
920                allocate ( yloc(jmax2) )
921                allocate ( zloc(kmax2) )
922     
923     
924     
925     
926             else
927                allocate (arr1(1))
928                allocate (arr2(1))
929                allocate (arr4d(1,1,1,1))
930                allocate ( xloc(1) )
931                allocate ( yloc(1) )
932                allocate ( zloc(1) )
933            end if
934     
935     !       CALL MPI_BARRIER(MPI_COMM_WORLD,mpierr)
936             if (myPE .ne. PE_IO) goto 1234
937     
938     
939             xloc(1) = -dx(1)
940             do II = 2,imax2
941                xloc(II) = xloc(II-1) + dx(II)
942             end do
943     
944     
945             yloc(1) = -dy(1)
946             do II = 2,jmax2
947                yloc(II) = yloc(II-1) + dy(II)
948             end do
949     
950             zloc(1) = -dz(1)
951             do II = 2,kmax2
952                zloc(II) = zloc(II-1) + dz(II)
953             end do
954     
955     
956     
957     
958             if (bFirst_netcdf_write .and. MFIX_usingNETCDF()) then
959                bFirst_netcdf_write = .false.
960                fname = trim(run_name) // '_netcdf_index.txt'
961                inquire (file=fname,exist=file_exists)
962     
963                if (file_exists) then
964                    open (unit=11,file=fname,status='old')
965                    read (11,*) netCDF_file_index
966                    close (unit=11)
967                else
968                    netCDF_file_index = 0
969                end if
970             end if
971     
972             fname_index = '_xxxxx.nc'
973             write (fname_index(2:6),'(i5.5)') netCDF_file_index
974             fname = trim(run_name)// fname_index
975             call MFIX_check_netcdf( MFIX_nf90_create(fname, NF90_64BIT_OFFSET, ncid) )
976             netCDF_file_index = netCDF_file_index + 1
977     
978             if (MFIX_usingNETCDF()) then
979                     fname = trim(run_name) // '_netcdf_index.txt'
980                     open (unit=11,file=fname,status='unknown')
981                     write (11,*) netCDF_file_index
982                     close (unit=11)
983             end if
984     
985     
986             call MFIX_check_netcdf( MFIX_nf90_def_dim(ncid, "x"   , imax2   , x_dimid  ) )  ! 1
987             call MFIX_check_netcdf( MFIX_nf90_def_dim(ncid, "y"   , jmax2   , y_dimid  ) )  ! 2
988             call MFIX_check_netcdf( MFIX_nf90_def_dim(ncid, "z"   , kmax2   , z_dimid  ) )  ! 3
989             call MFIX_check_netcdf( MFIX_nf90_def_dim(ncid, "coordinates"   , 1   , coords_dimid  ) )  ! 3
990             call MFIX_check_netcdf( MFIX_nf90_def_dim(ncid, "t"   , 1       , t_dimid  ) )  ! 4
991     
992     
993             call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, "x"  , NF90_DOUBLE, x_dimid, varid_x  ) )  ! 5
994             call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, "y"  , NF90_DOUBLE, y_dimid, varid_y  ) )  ! 6
995             call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, "z"  , NF90_DOUBLE, z_dimid, varid_z  ) )  ! 7
996             call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, "coordinates"  , NF90_INT, coords_dimid, varid_coords  ) )  ! 7
997             call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, "t"  , NF90_DOUBLE, t_dimid, varid_t  ) )  ! 8
998     
999             dimids =  (/ x_dimid , y_dimid, z_dimid , t_dimid/)
1000     
1001             if (bWrite_netcdf(1)) call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, "EP_g"  , NF90_DOUBLE, dimids, varid_epg  ) )  ! 9
1002             if (bWrite_netcdf(2)) call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, "P_g"   , NF90_DOUBLE, dimids, varid_pg  ) )   ! 10
1003     
1004             if (bWrite_netcdf(3)) call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, "P_star", NF90_DOUBLE, dimids, varid_pstar) )
1005             if (bWrite_netcdf(4)) call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, "U_g"   , NF90_DOUBLE, dimids, varid_ug   ) )
1006             if (bWrite_netcdf(4)) call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, "V_g"   , NF90_DOUBLE, dimids, varid_vg   ) )
1007             if (bWrite_netcdf(4)) call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, "W_g"   , NF90_DOUBLE, dimids, varid_wg   ) )
1008             if (bWrite_netcdf(7)) call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, "T_g"   , NF90_DOUBLE, dimids, varid_tg   ) )
1009             do i = 1,1   ! mmax
1010                var_name = 'U_s_xxx'
1011                write (var_name(5:7),'(i3.3)') I
1012                if (bWrite_netcdf(5)) call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, var_name, NF90_DOUBLE, dimids, varid_us(I)) )
1013     
1014                var_name = 'V_s_xxx'
1015                write (var_name(5:7),'(i3.3)') I
1016                if (bWrite_netcdf(5)) call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, var_name, NF90_DOUBLE, dimids, varid_vs(I)) )
1017     
1018                var_name = 'W_s_xxx'
1019                write (var_name(5:7),'(i3.3)') I
1020                if (bWrite_netcdf(5)) call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, var_name, NF90_DOUBLE, dimids, varid_ws(I)) )
1021     
1022                var_name = 'ROP_s_xxx'
1023                write (var_name(7:10),'(i3.3)') I
1024                if (bWrite_netcdf(6)) call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, var_name, NF90_DOUBLE, dimids, varid_rops(I)) )
1025     
1026                var_name = 'T_s_xxx'
1027                write (var_name(5:7),'(i3.3)') I
1028                if (bWrite_netcdf(8)) call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, var_name, NF90_DOUBLE, dimids, varid_ts(I)) )
1029     
1030                var_name = 'Theta_m_xxx'
1031                write (var_name(9:11),'(i3.3)') I
1032                if (bWrite_netcdf(11)) call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, var_name, NF90_DOUBLE, dimids, varid_thetam(I)) )
1033     
1034                DO N = 1, NMAX(i)
1035                   var_name = 'X_s_xxx_xxx'
1036                   write (var_name(5:7) ,'(i3.3)') I
1037                   write (var_name(9:11),'(i3.3)') n
1038                   if (bWrite_netcdf(10)) call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, var_name, NF90_DOUBLE, dimids, varid_xs(I,n)) )
1039                END DO
1040     
1041     
1042             end do
1043     
1044             do i = 1,nmax(0)
1045                var_name = 'X_g_xxx'
1046                write (var_name(5:7),'(i3.3)') I
1047                if (bWrite_netcdf(9)) call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, var_name, NF90_DOUBLE, dimids, varid_xg(I)) )
1048             end do
1049     
1050             do i = 1,nscalar
1051                var_name = 'Scalar_xxx'
1052                write (var_name(8:10),'(i3.3)') I
1053                if (bWrite_netCDF(12)) call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, var_name, NF90_DOUBLE, dimids, varid_scalar(I)) )
1054             end do
1055     
1056             do i = 1,nRR
1057                var_name = 'RRates_xxx'
1058                write (var_name(8:10),'(i3.3)') I
1059                if (bWrite_netCDF(13)) call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, var_name, NF90_DOUBLE, dimids, varid_rr(I)) )
1060             end do
1061     
1062     
1063             if (bWrite_netcdf(14) .and. k_Epsilon) then
1064                call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, 'k_turb_g', NF90_DOUBLE, dimids, varid_kturbg) )
1065                call MFIX_check_netcdf( MFIX_nf90_def_var(ncid, 'e_turb_g', NF90_DOUBLE, dimids, varid_eturbg) )
1066             end if
1067     
1068     
1069     
1070             call MFIX_check_netcdf( MFIX_nf90_enddef(ncid) ) ! 11
1071     
1072      1234   continue
1073                bFirst_netcdf_write = .false.
1074     
1075     !       CALL MPI_BARRIER(MPI_COMM_WORLD,mpierr)
1076     
1077     
1078             if (myPE .eq. PE_IO) then
1079                call MFIX_check_netcdf( MFIX_nf90_put_var(ncid,varid_t,the_time) )   ! 12
1080                call MFIX_check_netcdf( MFIX_nf90_put_var(ncid,varid_x,xloc) )
1081                call MFIX_check_netcdf( MFIX_nf90_put_var(ncid,varid_y,yloc) )
1082                call MFIX_check_netcdf( MFIX_nf90_put_var(ncid,varid_z,zloc) )
1083                coords = 0
1084                if (COORDINATES .eq. 'CYLINDRICAL') coords = 1
1085                call MFIX_check_netcdf( MFIX_nf90_put_var(ncid,varid_coords,coords) )
1086             end if
1087     
1088             if (bWrite_netcdf(1)) then
1089     
1090                 call gatherWriteSpx_netcdf(EP_g, arr1, arr2 , arr4d, ncid, varid_epg , &
1091                     imax2,jmax2,kmax2,ijkmax2 , ijkmax3)
1092     
1093             end if
1094     
1095             if (bWrite_netcdf(2)) then
1096     
1097                 call gatherWriteSpx_netcdf(P_g, arr1, arr2 , arr4d, ncid, varid_pg , &
1098                     imax2,jmax2,kmax2,ijkmax2 , ijkmax3)
1099     
1100             end if
1101     
1102     
1103     
1104     
1105             if (bWrite_netcdf(3)) then
1106     
1107                 call gatherWriteSpx_netcdf(P_star, arr1, arr2 , arr4d, ncid, varid_pstar , &
1108                     imax2,jmax2,kmax2,ijkmax2 , ijkmax3)
1109     
1110             end if
1111     
1112             if (bWrite_netcdf(4)) then
1113     
1114                 call gatherWriteSpx_netcdf(U_g, arr1, arr2 , arr4d, ncid, varid_ug , &
1115                     imax2,jmax2,kmax2,ijkmax2 , ijkmax3)
1116     
1117                 call gatherWriteSpx_netcdf(V_g, arr1, arr2 , arr4d, ncid, varid_vg , &
1118                     imax2,jmax2,kmax2,ijkmax2 , ijkmax3)
1119     
1120                 call gatherWriteSpx_netcdf(W_g, arr1, arr2 , arr4d, ncid, varid_wg , &
1121                     imax2,jmax2,kmax2,ijkmax2 , ijkmax3)
1122     
1123             end if
1124     
1125             if (bWrite_netcdf(7)) then
1126     
1127                 call gatherWriteSpx_netcdf(t_g, arr1, arr2 , arr4d, ncid, varid_tg , &
1128                     imax2,jmax2,kmax2,ijkmax2 , ijkmax3)
1129     
1130             end if
1131     
1132     
1133             do i = 1,1   ! mmax
1134     
1135                if (bWrite_netcdf(5)) then
1136     
1137                   call gatherWriteSpx_netcdf(u_s(:,i) , arr1, arr2 , arr4d, ncid, varid_us(i) , &
1138                        imax2,jmax2,kmax2,ijkmax2 , ijkmax3)
1139     
1140                   call gatherWriteSpx_netcdf(v_s(:,i) , arr1, arr2 , arr4d, ncid, varid_vs(i) , &
1141                        imax2,jmax2,kmax2,ijkmax2 , ijkmax3)
1142     
1143                   call gatherWriteSpx_netcdf(w_s(:,i) , arr1, arr2 , arr4d, ncid, varid_ws(i) , &
1144                        imax2,jmax2,kmax2,ijkmax2 , ijkmax3)
1145     
1146                end if
1147     
1148                if (bWrite_netcdf(6)) then
1149     
1150                   call gatherWriteSpx_netcdf(ROP_s(:,i) , arr1, arr2 , arr4d, ncid, varid_rops(i) , &
1151                        imax2,jmax2,kmax2,ijkmax2 , ijkmax3)
1152     
1153                end if
1154     
1155                if (bWrite_netcdf(8)) then
1156     
1157                   call gatherWriteSpx_netcdf(T_s(:,i) , arr1, arr2 , arr4d, ncid, varid_ts(i) , &
1158                        imax2,jmax2,kmax2,ijkmax2 , ijkmax3)
1159     
1160                end if
1161     
1162                if (bWrite_netcdf(11)) then
1163     
1164                   call gatherWriteSpx_netcdf(Theta_m(:,i) , arr1, arr2 , arr4d, ncid, varid_thetam(i) , &
1165                        imax2,jmax2,kmax2,ijkmax2 , ijkmax3)
1166     
1167                end if
1168     
1169                if (bWrite_netcdf(10)) then
1170                   do N = 1,nmax(i)
1171     
1172                      call gatherWriteSpx_netcdf(X_s(:,i,N) , arr1, arr2 , arr4d, ncid, varid_xs(i,N) , &
1173                           imax2,jmax2,kmax2,ijkmax2 , ijkmax3)
1174     
1175                   end do
1176                end if
1177     
1178     
1179             end do
1180     
1181             if (bWrite_netcdf(9)) then
1182                do i = 1,nmax(0)
1183     
1184                    call gatherWriteSpx_netcdf(X_g(:,i) , arr1, arr2 , arr4d, ncid, varid_xg(i) , &
1185                        imax2,jmax2,kmax2,ijkmax2 , ijkmax3)
1186     
1187                end do
1188             end if
1189     
1190            if (bWrite_netcdf(12)) then
1191               do i = 1,nscalar
1192     
1193                  call gatherWriteSpx_netcdf(Scalar(:,i) , arr1, arr2 , arr4d, ncid, varid_scalar(i) , &
1194                        imax2,jmax2,kmax2,ijkmax2 , ijkmax3)
1195     
1196                end do
1197             end if
1198     
1199     
1200             if (bWrite_netcdf(13)) then
1201                do i = 1,nRR
1202     
1203                  call gatherWriteSpx_netcdf(ReactionRates(:,i) , arr1, arr2 , arr4d, ncid, varid_rr(i) , &
1204                        imax2,jmax2,kmax2,ijkmax2 , ijkmax3)
1205     
1206               end do
1207             end if
1208     
1209             if (bWrite_netcdf(14) .and. k_Epsilon) then
1210     
1211                call gatherWriteSpx_netcdf(k_turb_g , arr1, arr2 , arr4d, ncid, varid_kturbg , &
1212                        imax2,jmax2,kmax2,ijkmax2 , ijkmax3)
1213     
1214                call gatherWriteSpx_netcdf(e_turb_g , arr1, arr2 , arr4d, ncid, varid_eturbg , &
1215                        imax2,jmax2,kmax2,ijkmax2 , ijkmax3)
1216     
1217             end if
1218     
1219     
1220     
1221             ! Close the file. This frees up any internal netCDF resources
1222             ! associated with the file, and flushes any buffers.
1223     !       CALL MPI_BARRIER(MPI_COMM_WORLD,mpierr)
1224             if (myPE .eq. PE_IO) then
1225                call MFIX_check_netcdf( MFIX_nf90_close(ncid) )
1226             end if
1227     
1228     
1229             deallocate (arr1)
1230             deallocate (arr2)
1231             deallocate (arr4d)
1232             deallocate (xloc)
1233             deallocate (yloc)
1234             deallocate (zloc)
1235     
1236             return
1237             end subroutine write_netcdf
1238     
1239