File: N:\mfix\model\main.f

1     ! -*- f90 -*-
2     !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
3     !                                                                      !
4     !  MODULE: MAIN                                                        !
5     !                                                                      !
6     !  Purpose: Main module for top level mfix subroutines.                !
7     !                                                                      !
8     !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
9     
10           MODULE MAIN
11     
12           use exit, only: mfix_exit
13     
14     !-----------------------------------------------
15     ! Module variables
16     !-----------------------------------------------
17     ! Final value of CPU time.
18           DOUBLE PRECISION :: CPU1
19     ! time used for computations.
20           DOUBLE PRECISION :: CPUTIME_USED, WALLTIME_USED
21     ! DISTIO variable for specifying the mfix version
22           CHARACTER(LEN=512) :: version
23     ! environment variable
24     !$ CHARACTER(LEN=512) :: omp_num_threads
25     !$ INTEGER :: length
26     !$ INTEGER :: status
27     
28     ! Number of iterations
29           INTEGER :: NIT_TOTAL
30     ! used for activating check_data_30
31           INTEGER :: NCHECK, DNCHECK
32     
33     ! Flag to save results and cleanly exit.
34           LOGICAL :: EXIT_SIGNAL = .FALSE.
35     
36           CHARACTER(LEN=80), DIMENSION(100) :: CMD_LINE_ARGS
37           INTEGER :: CMD_LINE_ARGS_COUNT = 0
38     
39           CONTAINS
40     
41           SUBROUTINE INITIALIZE
42     
43     #ifdef MPI
44           USE mpi, only: mpi_comm_world, mpi_barrier  ! ignore-depcomp
45     #endif
46           USE cdist, only: bdoing_postmfix
47           USE cdist, only: bglobalnetcdf, bstart_with_one_res, bdist_io, bwrite_netcdf
48           USE check, only: check_mass_balance
49           USE check_data_cg, only: check_bc_flags, report_best_processor_size
50           USE coeff, only: init_coeff
51           USE compar, only: mpierr, mype, pe_io
52           USE cont, only: do_cont
53           USE cutcell, only: cartesian_grid, re_indexing, set_corner_cells
54           USE discretelement, only: discrete_element
55           USE drag, only: f_gs
56           USE error_manager, only: err_msg, flush_err_msg
57           USE error_manager, only: init_err_msg, finl_err_msg
58           USE fldvar, only: rop_g, rop_s
59           USE funits, only: dmp_log, unit_log, unit_res
60           USE machine, only: start_log, end_log
61           USE machine, only: wall_time, pc_quickwin, machine_cons, get_run_id, start_log, end_log
62           USE mfix_netcdf, only: mfix_usingnetcdf
63           USE output, only: dbgprn_layout
64           USE output_man, only: init_output_vars, output_manager
65           USE parallel_mpi, only: parallel_init, parallel_fin
66           USE param1, only: n_spx, undefined, zero
67           USE pgcor, only: d_e, d_n, d_t, phase_4_p_g, switch_4_p_g
68           USE physprop, only: mmax
69           USE pscor, only: e_e, e_n, e_t, do_p_s, phase_4_p_s, mcp, switch_4_p_s
70           USE qmom_kinetic_equation, only: qmomk
71           USE read_input, only: get_data
72           USE run, only: id_version, ier
73           USE run, only: automatic_restart, call_usr, dem_solids, dt_max, dt_min
74           USE run, only: iter_restart, nstep, pic_solids, run_type, dt, shear, time, v_sh
75           USE time_cpu, only: CPU00, wall0
76           USE time_cpu, only: cpu_io, cpu_nlog, cpu0, cpuos, time_nlog
77           USE vtk, only: write_vtk_files
78     
79           IMPLICIT NONE
80     
81     !$    INTEGER num_threads, threads_specified, omp_id
82     !$    INTEGER omp_get_num_threads
83     !$    INTEGER omp_get_thread_num
84     
85           ! Temporary storage for DT
86           DOUBLE PRECISION :: DT_tmp
87           ! Save TIME in input file for RESTART_2
88           DOUBLE PRECISION :: TIME_SAVE
89     
90           INTEGER :: LL, MM
91     
92     ! DISTIO
93     ! If you change the value below in this subroutine, you must also
94     ! change it in write_res0.f and the value should also be consistent
95     ! with the check in read_res0
96           version = 'RES = 01.6'
97     
98           bDoing_postmfix = .false.
99     
100     ! Invoke MPI initialization routines and get rank info.
101           CALL PARALLEL_INIT
102           CALL GEN_LOG_BASENAME
103     
104     ! we want only PE_IO to write out common error messages
105           DMP_LOG = (myPE == PE_IO)
106     
107     ! set the version.release of the software
108           ID_VERSION = '2016-1'
109     
110     ! set automatic restart flag to false
111     !      AUTOMATIC_RESTART = .FALSE.
112     !      ITER_RESTART      = 1
113     
114     ! specify the number of processors to be used
115     !$        call get_environment_variable("OMP_NUM_THREADS",omp_num_threads,length,status, .true.)
116     !$      if (status.eq.0 .and. length.ne.0) then
117     !$        read(omp_num_threads,*) threads_specified
118     !$      else
119     !$        WRITE(*,'(A,$)') 'Enter the number of threads to be used for SMP: '
120     !$        READ(*,*) threads_specified
121     !$      endif
122     
123     !$      call omp_set_num_threads(threads_specified)
124     
125     ! Find the number of processors used
126     !$omp  parallel
127     !$      num_threads = omp_get_num_threads()
128     !$      omp_id = omp_get_thread_num()
129     !$      if(omp_id.eq.0) Write(*,*)' Number of threads used for SMP = ',  num_threads
130     !$omp  end parallel
131     
132     ! Set machine dependent constants
133           CALL MACHINE_CONS
134     
135     ! Get the date and time. They give the unique run_id in binary output
136     ! files
137           CALL GET_RUN_ID
138     
139     ! AEOLUS: stop trigger mechanism to terminate MFIX normally before batch
140     ! queue terminates. timestep at the beginning of execution
141           CALL CPU_TIME (CPU00)
142           WALL0 = WALL_TIME()
143     
144     ! Read input data, check data, do computations for IC and BC locations
145     ! and flows, and set geometry parameters such as X, X_E, DToDX, etc.
146           CALL GET_DATA
147     
148     ! Write the initial part of the standard output file
149           CALL WRITE_OUT0
150           IF(.NOT.CARTESIAN_GRID)  CALL WRITE_FLAGS
151     
152     ! Write the initial part of the special output file(s)
153           CALL WRITE_USR0
154     
155     !$    CALL START_LOG
156     !$    IF(DMP_LOG)WRITE (UNIT_LOG, *) ' '
157     !$    IF(DMP_LOG)WRITE (UNIT_LOG, *) ' Number of processors used = ', threads_specified
158     !$    IF(DMP_LOG)WRITE (UNIT_LOG, *) ' '
159     !$    CALL END_LOG
160     
161     !  setup for PC quickwin application
162           CALL PC_QUICKWIN
163     
164           CALL INIT_ERR_MSG('MFIX')
165     
166     
167     ! if not netcdf writes asked for ... globally turn off netcdf
168           if(MFIX_usingNETCDF()) then
169              bGlobalNetcdf = .false.
170              do LL = 1,20
171                 if (bWrite_netcdf(LL)) bGlobalNetcdf = .true.
172              enddo
173           endif
174     
175           DT_TMP = DT
176           SELECT CASE (TRIM(RUN_TYPE))
177     
178           CASE ('NEW')
179     ! Write the initial part of the restart files
180              CALL WRITE_RES0
181              DO LL = 1, N_SPX
182                 CALL WRITE_SPX0 (LL, 0)
183              ENDDO
184     
185           CASE ('RESTART_1')
186     ! Read the time-dependent part of the restart file
187              CALL READ_RES1
188              WRITE(ERR_MSG, 1010) TIME, NSTEP
189              CALL FLUSH_ERR_MSG()
190     
191           CASE ('RESTART_2')
192              TIME_SAVE = TIME
193     ! DISTIO
194              if (myPE .ne. PE_IO .and. bDist_IO .and. bStart_with_one_res) then
195                 write (unit_res,rec=1) version
196                 write (unit_res,rec=2) 4
197                 write (unit_res,rec=3) 4
198              endif
199     
200              CALL READ_RES1
201              TIME = TIME_SAVE
202     
203     1010     FORMAT('Message 1010: Read in data from .RES file for TIME = ',&
204                   G12.5,/'Time step number (NSTEP) =',I7)
205     
206              WRITE(ERR_MSG, 1010) TIME, NSTEP
207              CALL FLUSH_ERR_MSG()
208     
209              CALL WRITE_RES0
210     
211     ! Writing the RES1 and SPX1 can only be done here when re-indexing is turned off
212     ! This will be done after the cell re-indexing is done later in this file.
213     ! This allows restarting independently of the re-indexing setting between
214     ! the previous and current run.
215              IF(.NOT.RE_INDEXING) THEN
216                 CALL WRITE_RES1
217                 DO LL = 1, N_SPX
218                    CALL WRITE_SPX0 (LL, 0)
219                    CALL WRITE_SPX1 (LL, 0)
220                 END DO
221                 call write_netcdf(0,0,time)
222              ENDIF
223     
224           CASE DEFAULT
225              CALL START_LOG
226              IF(DMP_LOG)WRITE (UNIT_LOG, *) &
227                   ' MFIX: Do not know how to process'
228              IF(DMP_LOG)WRITE (UNIT_LOG, *) ' RUN_TYPE in data file'
229              CALL END_LOG
230              call mfix_exit(myPE)
231     
232           END SELECT
233     
234     #ifdef MPI
235           call MPI_Barrier(MPI_COMM_WORLD,mpierr)
236     #endif
237     
238           IF (DT_TMP /= UNDEFINED) THEN
239              DT = MAX(DT_MIN,MIN(DT_MAX,DT))
240           ELSE
241              DT = DT_TMP
242           ENDIF
243     
244     ! Set arrays for computing indices. A secondary call is made
245     ! after cut cell-preprocessing to update array indices.
246           IF(CARTESIAN_GRID) THEN
247              CALL SET_INCREMENTS
248              CALL SET_INCREMENTS3
249           ENDIF
250     
251     !      IF(.NOT.RE_INDEXING) CALL WRITE_IJK_VALUES
252     
253     ! Set the flags for wall surfaces impermeable and identify flow
254     ! boundaries using FLAG_E, FLAG_N, and FLAG_T
255           CALL SET_FLAGS1
256     
257     !  Update flags for Cartesian_GRID.
258           IF(CARTESIAN_GRID) CALL CHECK_BC_FLAGS
259     
260     ! Calculate cell volumes and face areas
261           IF(.NOT.CARTESIAN_GRID) CALL SET_GEOMETRY1
262     
263     ! Find corner cells and set their face areas to zero
264           IF(.NOT.CARTESIAN_GRID)  THEN
265              CALL GET_CORNER_CELLS()
266           ELSE
267              IF (SET_CORNER_CELLS)  CALL GET_CORNER_CELLS ()
268           ENDIF
269     
270     ! Set constant physical properties
271           CALL SET_CONSTPROP
272     
273     ! Set initial conditions
274           CALL SET_IC
275     
276     ! Set point sources.
277           CALL SET_PS
278     
279     ! Set boundary conditions
280           CALL ZERO_NORM_VEL
281           CALL SET_BC0
282     
283     ! Cartesian grid implementation
284           IF(CARTESIAN_GRID) CALL CG_SET_BC0
285     
286     ! Set gas mixture molecular weight
287           CALL SET_MW_MIX_G
288     
289     ! Set the pressure field for a fluidized bed
290           IF (RUN_TYPE == 'NEW') CALL SET_FLUIDBED_P
291     
292     ! Initialize densities.
293           IF (RUN_TYPE == 'NEW') CALL SET_RO_G
294           IF (RUN_TYPE == 'NEW') CALL SET_RO_S
295     
296     ! Initialize time dependent boundary conditions
297           CALL SET_BC1
298     
299     ! Check the field variable data and report errors.
300           IF(.NOT.CARTESIAN_GRID)  CALL CHECK_DATA_20
301     
302     !=======================================================================
303     ! JFD: START MODIFICATION FOR RE-INDEXING CELLS
304     !=======================================================================
305           IF(CARTESIAN_GRID.AND.RE_INDEXING) THEN
306     
307              IF(myPE == PE_IO) THEN
308                 WRITE(*,"(72('='))")
309                 WRITE(*,*)' RE-INDEXING CELLS FOR CARTESIAN GRID...'
310              ENDIF
311              CALL RE_INDEX_ARRAYS
312     
313     
314              !IF(myPE == PE_IO)print*,'Calling REPORT_BEST_IJK_SIZE:'
315              !CALL REPORT_BEST_IJK_SIZE
316              CALL REPORT_BEST_PROCESSOR_SIZE
317              !IF(myPE == PE_IO)print*,'Exiting MFIX after REPORT_BEST_IJK_SIZE.'
318     
319              IF(myPE == PE_IO) WRITE(*,"(72('='))")
320     
321     ! In case of a RESTART_2, write the RES1 and SPX1 files here
322     ! This was commented out earlier in this file.
323              IF(RUN_TYPE == 'RESTART_2') THEN
324                 CALL WRITE_RES1
325                 DO LL = 1, N_SPX
326                    CALL WRITE_SPX0 (LL, 0)
327                    CALL WRITE_SPX1 (LL, 0)
328                 END DO
329                 call write_netcdf(0,0,time)
330              ENDIF
331           ENDIF
332     !=======================================================================
333     ! JFD: END MODIFICATION FOR RE-INDEXING CELLS
334     !=======================================================================
335     
336     ! Setup VTK data for regular (no cut cells) grid
337           IF(.NOT.CARTESIAN_GRID.AND.WRITE_VTK_FILES) CALL SETUP_VTK_NO_CUTCELL
338     
339           IF(DISCRETE_ELEMENT) CALL MAKE_ARRAYS_DES
340           IF(QMOMK) CALL QMOMK_MAKE_ARRAYS
341     
342     ! Set the inflow/outflow BCs for DEM solids
343           IF(DEM_SOLIDS) CALL SET_BC_DEM
344     ! Set the inflow/outflow BC for PIC solids
345           IF(PIC_SOLIDS) CALL SET_BC_PIC
346     
347     ! Set the inital properties of each particle.
348           IF(DEM_SOLIDS) CALL SET_IC_DEM
349     
350     ! AEOLUS: debug prints
351           if (DBGPRN_LAYOUT .or. bdist_io) then
352              !write (*,*) myPE , ' E.4 ... version = ' , version(1:33)
353              call debug_write_layout()
354              call write_parallel_info()
355           endif
356     
357     ! Initializations for CPU time calculations in iterate
358           CPUOS = 0.
359           CALL CPU_TIME (CPU1)
360           CPU_NLOG = CPU1
361           TIME_NLOG = TIME - DT
362     
363     ! Get the initial value of CPU time
364           CALL CPU_TIME (CPU0)
365     
366     ! Find the solution of the equations from TIME to TSTOP at
367     ! intervals of DT
368     
369     !-----------------------------------------------
370     
371           NCHECK  = NSTEP
372           DNCHECK = 1
373           CPU_IO  = ZERO
374           NIT_TOTAL = 0
375     
376           CALL INIT_OUTPUT_VARS
377     
378     ! Parse residual strings
379           CALL PARSE_RESID_STRING ()
380     
381     ! Call user-defined subroutine to set constants, check data, etc.
382           IF (CALL_USR) CALL USR0
383     
384           CALL RRATES_INIT()
385     
386     ! Calculate all the coefficients once before entering the time loop
387           CALL INIT_COEFF(IER)
388     
389           DO MM=1, MMAX
390              F_gs(1,MM) = ZERO
391           ENDDO
392     
393     ! Remove undefined values at wall cells for scalars
394           CALL UNDEF_2_0 (ROP_G)
395           DO MM = 1, MMAX
396              CALL UNDEF_2_0 (ROP_S(1,MM))
397           ENDDO
398     
399     ! Initialize d's and e's to zero
400           DO MM = 0, MMAX
401              D_E(1,MM) = ZERO
402              D_N(1,MM) = ZERO
403              D_T(1,MM) = ZERO
404           ENDDO
405           E_E(:) = ZERO
406           E_N(:) = ZERO
407           E_T(:) = ZERO
408     
409     ! calculate shear velocities if periodic shear BCs are used
410           IF(SHEAR) CALL CAL_D(V_sh)
411     
412     ! Initialize check_mass_balance.  This routine is not active by default.
413     ! Specify a reporting interval (hard-wired in the routine) to activate
414     ! the routine.
415           Call check_mass_balance (0)
416     
417     ! sof modification: now it's only needed to do this once before time-loop
418     ! Mark the phase whose continuity will be solved and used to correct
419     ! void/volume fraction in calc_vol_fr (see subroutine for details)
420           CALL MARK_PHASE_4_COR (PHASE_4_P_G, PHASE_4_P_S, DO_CONT, MCP,&
421                DO_P_S, SWITCH_4_P_G, SWITCH_4_P_S)
422     
423           END SUBROUTINE INITIALIZE
424     
425           SUBROUTINE FINALIZE
426     
427           USE cutcell, only: cartesian_grid
428           USE dashboard
429           USE cut_cell_preproc, only: close_cut_cell_files
430           USE error_manager, only: finl_err_msg
431           USE machine, only: wall_time
432           USE parallel_mpi, only: parallel_fin
433           USE run, only: dt, call_usr, dt_min, get_tunit, tunit
434           USE time_cpu
435           IMPLICIT NONE
436     
437     ! Call user-defined subroutine after time-loop.
438           IF (CALL_USR) CALL USR3
439     
440     ! Get the final value of CPU time.  The difference gives the
441     ! CPU time used for the computations.
442           CALL CPU_TIME (CPU1)
443     
444     ! Compute the CPU time and write it out in the .OUT file.
445           CPUTIME_USED = CPU1 - CPU0 - CPU_IO
446           WALLTIME_USED = WALL_TIME() - WALL0
447           CALL WRITE_OUT3 (CPUTIME_USED, WALLTIME_USED, CPU_IO)
448     
449     ! JFD: cartesian grid implementation
450           IF(WRITE_DASHBOARD) THEN
451              IF(DT>=DT_MIN) THEN
452                 RUN_STATUS = 'Complete.'
453              ELSE
454                 RUN_STATUS = 'DT < DT_MIN.  Recovery not possible!'
455              ENDIF
456              CALL GET_TUNIT(CPUTIME_USED,TUNIT)
457              CALL UPDATE_DASHBOARD(0,CPUTIME_USED,TUNIT)
458           ENDIF
459           IF(CARTESIAN_GRID)  CALL CLOSE_CUT_CELL_FILES
460     
461     ! Finalize and terminate MPI
462           call parallel_fin
463     
464           CALL FINL_ERR_MSG
465     
466           END SUBROUTINE FINALIZE
467     
468     !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
469     !                                                                      !
470     !  Subroutine: GEN_LOG_BASENAME                                        !
471     !  Author: Aytekin Gel                                Date: 19-SEP-03  !
472     !                                                                      !
473     !  Purpose: Generate the file base for DMP logs.                       !
474     !                                                                      !
475     !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
476           SUBROUTINE GEN_LOG_BASENAME
477     
478           use compar, only: myPE
479           use compar, only: fbname
480     
481           implicit none
482     
483     ! Variables for generating file basename with processor id
484           INTEGER :: i1, i10, i100, i1000, i10000
485     
486     ! PAR_I/O Generate file basename for LOG files
487           i10000 = int(myPE/10000)
488           i1000  = int((myPE-i10000*10000)/1000)
489           i100   = int((myPE-i10000*10000-i1000*1000)/100)
490           i10    = int((myPE-i10000*10000-i1000*1000-i100*100)/10)
491           i1     = int((myPE-i10000*10000-i1000*1000-i100*100-i10*10)/1)
492     
493           i10000 = i10000 + 48
494           i1000  = i1000  + 48
495           i100   = i100   + 48
496           i10    = i10    + 48
497           i1     = i1     + 48
498     
499           fbname=char(i10000)//char(i1000)//char(i100)//char(i10)//char(i1)
500     
501           RETURN
502           END SUBROUTINE GEN_LOG_BASENAME
503     
504     !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvC
505     !                                                                      C
506     !  Module name: debug_write()                                          C
507     !  Purpose: Write out full geometry index setup information for the
508     !  case
509     !                                                                      C
510     !  Author: Aytekin Gel                                Date: 19-SEP-03  C
511     !  Reviewer:                                          Date:            C
512     !                                                                      C
513     !                                                                      C
514     !  Literature/Document References:                                     C
515     !                                                                      C
516     !  Variables referenced:                                               C
517     !  Variables modified:                                                 C
518     !                                                                      C
519     !  Local variables:                                                    C
520     !                                                                      C
521     !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^C
522     
523        SUBROUTINE debug_write_layout()
524     
525     !-----------------------------------------------
526     ! Modules
527     !-----------------------------------------------
528           USE cdist
529           USE compar
530           USE functions
531           USE funits
532           USE geometry
533           USE indices
534           USE leqsol
535           USE mpi_utility
536           USE parallel
537           USE param
538           USE param1
539           USE run
540           USE sendrecv
541           USE sendrecv3
542           USE time_cpu
543           IMPLICIT NONE
544     !-----------------------------------------------
545     ! Local Variables
546     !-----------------------------------------------
547     ! phase index
548           INTEGER :: M
549     ! indices
550           INTEGER :: i, j, k, ijk, ijk_GL, ijk_PROC, ijk_IO
551     !
552           integer :: indxA, indxA_gl, indxB, indxB_gl, indxC, indxC_gl
553           integer :: indxD, indxD_gl, indxE, indxE_gl, indxF, indxF_gl
554           integer :: indxG, indxG_gl, indxH, indxH_gl
555     !
556           logical :: amgdbg = .TRUE.
557     
558           character(LEN=80) :: fname
559     
560     !DISTIO
561     !      fname = "layout_xxxx.txt"
562     !      write (fname(8:11),'(i4.4)') myPE
563           fname = "layout_xxxxx.txt"
564           write (fname(8:12),'(i5.5)') myPE
565           open (unit=11,file=fname,status='unknown')
566     
567           write (11,*) ' ********************************************'
568           write (11,*) ' ********************************************'
569           write (11,*) ' ********************************************'
570           write (11,*) ' ********************************************'
571           write (11,*) ' '
572           write (11,*) ' '
573           write (11,*) ' myPE =           ' , myPE
574           write (11,*) ' '
575           write (11,*) ' '
576     
577     
578           IF (AMGDBG .OR. bDist_IO) THEN
579              write(11,"('BLK1: Running from istart3,iend3 .AND. jstart3, jend3 .AND. kstart3, kend3')")
580              write(11,"(' (   i ,    j,     k) =>    ijk      ijk_GL     ijk_PROC    ijk_IO')")
581              write(11,"(' ====================      =====     =======    ========    ======')")
582              DO k = kstart3, kend3
583                 DO i = istart3,iend3
584                    DO j = jstart3, jend3
585                       ijk = FUNIJK(i,j,k)
586                       ijk_GL = FUNIJK_GL(i,j,k)
587                       ijk_PROC = FUNIJK_PROC(i,j,k,myPE)
588                       ijk_IO = FUNIJK_IO(i,j,k)
589                       write(11,"(' (',I4,' , ',I4,' , ',I4,') => ',4(I8,' , '))") &
590                            i,j,k,ijk,ijk_GL,ijk_PROC,ijk_IO
591                    ENDDO
592                 ENDDO
593              ENDDO
594     
595              write(11,"(/,/,'BLK2: Print out Bottom, South, West, East, North, Top neighbors')")
596              write(11,"(' (   i ,    j,     k) =>    ijk    ijk_GL    B_of    S_of    W_of    E_of    N_of    T_of')")
597              write(11,"(' ====================      =====   =======  ======  ======  ======  ======  ======  ======')")
598              DO k = kstart3, kend3
599                 DO i = istart3,iend3
600                    DO j = jstart3, jend3
601                       ijk = FUNIJK(i,j,k)
602                       ijk_GL = FUNIJK_GL(i,j,k)
603                       write(11,"(' (',I4,' , ',I4,' , ',I4,') => ',2(I7,' , '),6(I7,2X))") &
604                            i,j,k,ijk,ijk_GL,bottom_of(ijk),south_of(ijk),west_of(ijk),&
605                            east_of(ijk),north_of(ijk),top_of(ijk)
606                    ENDDO
607                 ENDDO
608              ENDDO
609     
610              write(11,"(/,/,'BLK3: Print out km, jm, im, ip, jp, kp neighbors')")
611              write(11,"(' (   i ,    j,     k) =>    ijk    ijk_GL    km_of   jm_of   im_of   ip_of   jp_of   kp_of')")
612              write(11,"(' ====================      =====   =======  ======  ======  ======  ======  ======  ======')")
613              DO k = kstart3, kend3
614                 DO i = istart3,iend3
615                    DO j = jstart3, jend3
616                       ijk = FUNIJK(i,j,k)
617                       ijk_GL = FUNIJK_GL(i,j,k)
618                       write(11,"(' (',I4,' , ',I4,' , ',I4,') => ',2(I7,' , '),6(I7,2X))") &
619                            i,j,k,ijk,ijk_GL,km_of(ijk),jm_of(ijk),im_of(ijk),&
620                            ip_of(ijk),jp_of(ijk),kp_of(ijk)
621                    ENDDO
622                 ENDDO
623              ENDDO
624     
625              write(11,"(/,'BLK4a: Active Fluid Cells:FLUID_AT(ijk)=.T.',/,&
626                   &           ' (   i ,    j,     k) =>    ijk  [   x ,     ,     z]')")
627              write(11,"(' ====================      =====  ====================')")
628              DO ijk = ijkstart3, ijkend3
629                 I = I_OF(IJK)
630                 J = J_OF(IJK)
631                 K = K_OF(IJK)
632     
633                 !         IF (FLOW_AT_E(IJK)) THEN
634                 IF (FLUID_AT(IJK)) THEN
635                    !          write(11,"(' (',I4,' , ',I4,' , ',I4,') => ',I8)") I,J,K,ijk
636                    write(11,"(' (',I4,' , ',I4,' , ',I4,') => ',I8,' [',E12.5,',',E12.5,' ]')") I,J,K,ijk,X(i),Z(k)
637                 ENDIF
638              ENDDO
639     
640              write(11,"(/,'BLK4b: Cells that are (.NOT.WALL_AT(IJK)) = .T.',/,&
641                   &           ' (   i ,    j,     k) =>    ijk  [   x ,     ,     z]')")
642              write(11,"(' ====================      =====  ====================')")
643              DO ijk = ijkstart3, ijkend3
644                 I = I_OF(IJK)
645                 J = J_OF(IJK)
646                 K = K_OF(IJK)
647     
648                 IF (.NOT.WALL_AT(IJK)) THEN
649                    !          write(11,"(' (',I4,' , ',I4,' , ',I4,') => ',I8)") I,J,K,ijk
650                    write(11,"(' (',I4,' , ',I4,' , ',I4,') => ',I8,' [',E12.5,',',E12.5,' ]')") I,J,K,ijk,X(i),Z(k)
651                 ENDIF
652              ENDDO
653     
654              DO k = kstart3, kend3
655                 DO i = istart3,iend3
656                    DO j = jstart3, jend3
657                       ijk = FUNIJK(i,j,k)
658                       ijk_GL = FUNIJK_GL(i,j,k)
659     
660                       if (i == istart2 .AND. j == jstart2) then
661                          indxA = ijk
662                          indxA_gl = ijk_GL
663                       endif
664                       if (i == istart1 .AND. j == jstart1) then
665                          indxE = ijk
666                          indxE_gl = ijk_GL
667                       endif
668                       if (i == istart2 .AND. j == jend2) then
669                          indxB = ijk
670                          indxB_gl = ijk_GL
671                       endif
672                       if (i == istart1 .AND. j == jend1) then
673                          indxF = ijk
674                          indxF_gl = ijk_GL
675                       endif
676                       if (i == iend1 .AND. j == jstart1) then
677                          indxH = ijk
678                          indxH_gl = ijk_GL
679                       endif
680                       if (i == iend2 .AND. j == jstart2) then
681                          indxD = ijk
682                          indxD_gl = ijk_GL
683                       endif
684                       if (i == iend1 .AND. j == jend1) then
685                          indxG = ijk
686                          indxG_gl = ijk_GL
687                       endif
688                       if (i == iend2 .AND. j == jend2) then
689                          indxC = ijk
690                          indxC_gl = ijk_GL
691                       endif
692                    ENDDO
693                 ENDDO
694                 write(11,"('BLK5:')")
695                 write(11,"(57('='))")
696                 write(11,"('k= ',I5,/,57('='))") k
697                 write(11,"('B= ',I5,' (',I7,')',20X,'C= ',I5,' (',I7,')',/)") indxB, indxB_gl, &
698                      indxC, indxC_gl
699                 !        write(UNIT_LOG,"(' \',34X,'/')")
700                 !        write(UNIT_LOG,"(2X,'\',32X,'/')")
701                 write(11,"(3X,'F= ',I5,' (',I7,')',12X,'G= ',I5,' (',I7,')')") indxF, indxF_gl, &
702                      indxG, indxG_gl
703                 write(11,"(4(9X,'|',29X,'|',/))")
704                 write(11,"(3X,'E= ',I5,' (',I7,')',12X,'H= ',I5,' (',I7,')',/)") indxE, indxE_gl, &
705                      indxH, indxH_gl
706                 !        write(UNIT_LOG,"(2X,'/',32X,'\')")
707                 !        write(UNIT_LOG,"('/',34X,'\')")
708                 write(11,"('A= ',I5,' (',I7,')',20X,'D= ',I5,' (',I7,')',/,/)") indxA, indxA_gl, &
709                      indxD, indxD_gl
710     
711                 !        write(UNIT_LOG,"(' (',I4,' , ',I4,' , ',I4,') => ',2(I7,' , '),6(I7,2X))") &
712                 !                                         i,j,k,ijk,ijk_GL,bottom_of(ijk),south_of(ijk),west_of(ijk),&
713                 !                                        east_of(ijk),north_of(ijk),top_of(ijk)
714     
715              ENDDO
716     
717              !      write(UNIT_LOG,"(/,' (   i ,    j,     k) =>    ijk (Active Fluid)')")
718              !      write(UNIT_LOG,"(' ====================      =====')")
719              !       DO ijk = ijkstart3, ijkend3
720              !         I = I_OF(IJK)
721              !         J = J_OF(IJK)
722              !         K = K_OF(IJK)
723     
724              !         IF (FLOW_AT_E(IJK)) THEN
725              !         IF (FLUID_AT(IJK)) THEN
726              !           write(UNIT_LOG,"(' (',I4,' , ',I4,' , ',I4,') => ',I8)") I,J,K,ijk
727              !         ENDIF
728              !      END DO
729     
730     
731           endif   ! end if(amgdbg .or. bdist_io)
732     
733           M = 0
734           !      CALL WRITE_AB_M (A_M, B_M, IJKMAX2, M, IER)
735     
736           IF (AMGDBG .OR. bDist_IO) THEN
737              write(11,"(/,/,'BLK6: ========= ORIGINAL MFIX VARIABLES ===========')")
738              write(11,"('PE ',I5,': imin1  = ',I6,3X,'imax1= ',I6,/,'PE ',I5,': jmin1  = ',I6,3X,'jmax1= ',I6)") &
739                   myPE,imin1,imax1,myPE,jmin1,jmax1
740              write(11,"('PE ',I5,': kmin1  = ',I6,3X,'kmax1= ',I6)") myPE,kmin1,kmax1
741              write(11,"('-----')")
742              write(11,"('PE ',I5,': imin2  = ',I6,3X,'imax2= ',I6,/,'PE ',I5,': jmin2  = ',I6,3X,'jmax2= ',I6)") &
743                   myPE,imin2,imax2,myPE,jmin2,jmax2
744              write(11,"('PE ',I5,': kmin2  = ',I6,3X,'kmax2= ',I6)") myPE,kmin2,kmax2
745              write(11,"('----- Below xxx3 set is DMP extension ------------')")
746              write(11,"('PE ',I5,': imin3  = ',I6,3X,'imax3= ',I6,/,'PE ',I5,': jmin3  = ',I6,3X,'jmax3= ',I6)") &
747                   myPE,imin3,imax3,myPE,jmin3,jmax3
748              write(11,"('PE ',I5,': kmin3  = ',I6,3X,'kmax3= ',I6)") myPE,kmin3,kmax3
749              write(11,"('----- End of Below xxx3 set is DMP extension -----')")
750              !      write(11,"('PE ',I5,': ijkmax2= ',I6)") myPE,ijkmax2
751              write(11,"('PE ',I5,': ijmax2 = ',I6)") myPE,ijmax2
752              write(11,"('PE ',I5,': ijkmin1= ',I6,' ijkmax1= ',I12)") myPE,ijkmin1, ijkmax1
753              write(11,"('PE ',I5,':          ',6X,' ijkmax2= ',I12)") myPE,ijkmax2
754              write(11,"('PE ',I5,':          ',6X,' ijkmax3= ',I12)") myPE,ijkmax3
755              write(11,"('PE ',I5,': ijkmin4= ',I6,' ijkmax4= ',I12)") myPE,ijkmin4, ijkmax4
756     
757     
758              write(11,"(/,/,' ========= DMP EXTENSION VARIABLES ===========')")
759              !      write(UNIT_LOG,"('PE ',I5,': ijksize  = ',I6)") myPE,ijksize
760              write(11,"('PE ',I5,': ijksize3 = ',I6,3X,'ijksize3_all = ',I6)") myPE,ijksize3,ijksize3_all(myPE)
761              write(11,"('PE ',I5,': ijksize4 = ',I6,3X,'ijksize4_all = ',I6)") myPE,ijksize4,ijksize4_all(myPE)
762              write(11,"('PE ',I5,': ijkstart3  = ',I6,3X,'ijkend3  = ',I6)") myPE,ijkstart3, ijkend3
763              write(11,"('PE ',I5,': ijkstart3_all = ',I6,3X,'ijkstart4_all = ',I6)") myPE,ijkstart3_all(myPE),ijkstart4_all(myPE)
764              write(11,"('PE ',I5,': istart_all = ',I6,3X,'iend_all = ',I6,/,'PE ',I5,': jstart_all = ',I6,3X,'jend_all = ',I6)") &
765                   myPE,istart_all(myPE),iend_all(myPE),myPE,jstart_all(myPE),jend_all(myPE)
766              write(11,"('PE ',I5,': kstart_all = ',I6,3X,'kend_all = ',I6,/,'----------------------')") &
767                   myPE,kstart_all(myPE),kend_all(myPE)
768     
769              write(11,"('PE ',I5,': istart1_all= ',I6,3X,'iend1_all= ',I6,/,'PE ',I5,': jstart1_all= ',I6,3X,'jend3_all= ',I6)") &
770                   myPE,istart1_all(myPE),iend1_all(myPE),myPE,jstart1_all(myPE),jend1_all(myPE)
771              write(11,"('PE ',I5,': kstart1_all= ',I6,3X,'kend1_all= ',I6,/,'----------------------')") &
772                   myPE,kstart1_all(myPE),kend1_all(myPE)
773     
774              write(11,"('PE ',I5,': istart2_all= ',I6,3X,'iend2_all= ',I6,/,'PE ',I5,': jstart2_all= ',I6,3X,'jend3_all= ',I6)") &
775                   myPE,istart2_all(myPE),iend2_all(myPE),myPE,jstart2_all(myPE),jend2_all(myPE)
776              write(11,"('PE ',I5,': kstart2_all= ',I6,3X,'kend2_all= ',I6,/,'----------------------')") &
777                   myPE,kstart2_all(myPE),kend2_all(myPE)
778     
779              write(11,"('PE ',I5,': istart3_all= ',I6,3X,'iend3_all= ',I6,/,'PE ',I5,': jstart3_all= ',I6,3X,'jend3_all= ',I6)") &
780                   myPE,istart3_all(myPE),iend3_all(myPE),myPE,jstart3_all(myPE),jend3_all(myPE)
781              write(11,"('PE ',I5,': kstart3_all= ',I6,3X,'kend3_all= ',I6,/,'----------------------')") &
782                   myPE,kstart3_all(myPE),kend3_all(myPE)
783     
784              write(11,"('PE ',I5,': istart1= ',I6,3X,'iend1= ',I6,/,'PE ',I5,': jstart1= ',I6,3X,'jend1= ',I6)") &
785                   myPE,istart1,iend1,myPE,jstart1,jend1
786              write(11,"('PE ',I5,': kstart1= ',I6,3X,'kend1= ',I6,/,'----------------------')") &
787                   myPE,kstart1,kend1
788              write(11,"('PE ',I5,': istart2= ',I6,3X,'iend2= ',I6,/,'PE ',I5,': jstart2= ',I6,3X,'jend2= ',I6)") &
789                   myPE,istart2,iend2,myPE,jstart2,jend2
790              write(11,"('PE ',I5,': kstart2= ',I6,3X,'kend2= ',I6,/,'----------------------')") &
791                   myPE,kstart2,kend2
792              write(11,"('PE ',I5,': istart3= ',I6,3X,'iend3= ',I6,/,'PE ',I5,': jstart3= ',I6,3X,'jend3= ',I6)") &
793                   myPE,istart3,iend3,myPE,jstart3,jend3
794              write(11,"('PE ',I5,': kstart3= ',I6,3X,'kend3= ',I6,/,'----------------------')") &
795                   myPE,kstart3,kend3
796     
797           ENDIF   ! end if(amgdbg .or. bdist_io)
798     
799           close(unit=11)
800     
801     
802           RETURN
803        END SUBROUTINE DEBUG_WRITE_LAYOUT
804     
805     
806     
807        !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvC
808     
809        !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^C
810     
811        SUBROUTINE write_parallel_info()
812     
813           !-----------------------------------------------
814           !   M o d u l e s
815           !-----------------------------------------------
816           USE compar
817           USE functions
818           USE funits
819           USE geometry
820           USE indices
821           USE leqsol
822           USE mpi_utility
823           USE parallel
824           USE param
825           USE param1
826           USE run
827           USE sendrecv
828           USE sendrecv3
829           USE time_cpu
830           IMPLICIT NONE
831           !-----------------------------------------------
832           ! Dummy arguments
833           !-----------------------------------------------
834           ! Local Variables
835           !-----------------------------------------------
836           ! phase index
837           INTEGER :: M
838           ! indices
839           INTEGER :: i, j, k, ijk, ijk_GL, ijk_PROC, ijk_IO
840           !
841           character(LEN=80) :: fname
842           !-----------------------------------------------
843     
844           !DISTIO
845           !      fname = "p_info_xxxx.txt"
846           !      write (fname(8:11),'(i4.4)') myPE
847           fname = "p_info_xxxxx.txt"
848           write (fname(8:12),'(i5.5)') myPE
849           open (unit=11,file=fname,status='unknown')
850     
851           write (11,*) myPe , ' = myPE'
852     
853           write (11,*) myPE , istart3,iend3
854           write (11,*) myPE , jstart3,jend3
855           write (11,*) myPE , kstart3,kend3
856     
857           write(11,"('BLK1: Running from istart3,iend3 .AND. jstart3, jend3 .AND. kstart3, kend3')")
858           write(11,"(' (   i ,    j,     k)       ijk      ijk_GL     ijk_PROC    ijk_IO')")
859           write(11,"(' ====================      =====     =======    ========    ======')")
860           DO k = kstart3, kend3
861              DO i = istart3,iend3
862                 DO j = jstart3, jend3
863                    ijk = FUNIJK(i,j,k)
864                    ijk_GL = FUNIJK_GL(i,j,k)
865                    ijk_PROC = FUNIJK_PROC(i,j,k,myPE)
866                    ijk_IO = FUNIJK_IO(i,j,k)
867                    write(11,"('  ',I4,'   ',I4,'   ',I4,'     ',4(I8,'   '))" ) &
868                         i,j,k,ijk,ijk_GL,ijk_PROC,ijk_IO
869                 ENDDO
870              ENDDO
871           ENDDO
872     
873           M = 0
874           !      CALL WRITE_AB_M (A_M, B_M, IJKMAX2, M, IER)
875     
876           close(unit=11)
877     
878           RETURN
879        END SUBROUTINE write_parallel_info
880     
881        !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
882        !                                                                      !
883        !  SUBROUTINE: DO_MPI_BCAST                                            !
884        !                                                                      !
885        !  Purpose: Used by pymfix for broadcasting commands.                  !
886        !                                                                      !
887        !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
888        function do_mpi_bcast(str)
889     #ifdef MPI
890           use mpi ! ignore-depcomp
891     #endif
892           implicit none
893     
894           ! TODO: return a dynamically allocated string, instead of fixed size of 100,000 bytes
895           character(len=100000),intent(in) :: str
896           character :: aa(100000)
897           character :: do_mpi_bcast(100000)
898           integer :: ii
899           integer :: ierr
900     
901           do ii = 1,len(str)
902              aa(ii) = str(ii:ii)
903           end do
904     
905     #ifdef MPI
906           call mpi_bcast(aa,100000,mpi_character,0,mpi_comm_world,ierr)
907     #endif
908     
909           do ii = 1,100000
910              do_mpi_bcast(ii:ii) = aa(ii)
911           end do
912     
913        end function do_mpi_bcast
914     
915        subroutine do_write_dbg_vtu_and_vtp_files
916           implicit none
917           call write_dbg_vtu_and_vtp_files
918        end subroutine do_write_dbg_vtu_and_vtp_files
919     
920        subroutine do_backupres
921           use output_man, only: backup_res
922           implicit none
923           call backup_res
924        end subroutine do_backupres
925     
926        subroutine do_reinit(filename)
927           use reinit, only: reinitialize
928           implicit none
929           ! filename of uploaded mfix.dat file
930           character(len=*), intent(in) :: filename
931           call reinitialize(filename)
932        end subroutine do_reinit
933     
934        subroutine do_abort
935           use compar, only: mype
936           implicit none
937           call mfix_exit(mype)
938        end subroutine do_abort
939     
940        !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
941        !  Subroutine: ADD_COMMAND_LINE_ARGUMENT                               !
942        !  Author: M.Meredith                                 Date: 03-FEB-16  !
943        !                                                                      !
944        !  Purpose: Save command line arguments in CMD_LINE_ARGS array.        !
945        !           Used by both mfix.f and pymfix.                            !
946        !                                                                      !
947        !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
948        SUBROUTINE ADD_COMMAND_LINE_ARGUMENT(ARG)
949           implicit none
950           CHARACTER(LEN=80), INTENT(IN) :: ARG
951     
952           CMD_LINE_ARGS_COUNT = CMD_LINE_ARGS_COUNT + 1
953     
954           if (CMD_LINE_ARGS_COUNT > 100) THEN
955              print *,"TOO MANY COMMAND LINE ARGUMENTS"
956              stop
957           ENDIF
958     
959           CMD_LINE_ARGS(CMD_LINE_ARGS_COUNT) = arg
960     
961        END SUBROUTINE ADD_COMMAND_LINE_ARGUMENT
962     
963     END MODULE MAIN
964