MFC
Exascale flow solver
Loading...
Searching...
No Matches
m_data_output.fpp.f90
Go to the documentation of this file.
1# 1 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2!>
3!! @file
4!! @brief Contains module m_data_output
5
6# 1 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp" 1
7# 1 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp" 1
8# 1 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp" 1
9# 2 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
10# 3 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
11# 4 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
12# 5 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
13# 6 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
14
15# 8 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
16# 9 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
17# 10 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
18
19# 17 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
20
21# 46 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
22
23# 58 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
24
25# 68 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
26
27# 98 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
28
29# 110 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
30
31# 120 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
32! New line at end of file is required for FYPP
33# 2 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp" 2
34# 1 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp" 1
35# 1 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp" 1
36# 2 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
37# 3 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
38# 4 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
39# 5 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
40# 6 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
41
42# 8 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
43# 9 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
44# 10 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
45
46# 17 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
47
48# 46 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
49
50# 58 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
51
52# 68 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
53
54# 98 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
55
56# 110 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
57
58# 120 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
59! New line at end of file is required for FYPP
60# 2 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp" 2
61
62# 4 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
63# 5 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
64# 6 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
65# 7 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
66# 8 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
67
68# 20 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
69
70# 43 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
71
72# 48 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
73
74# 53 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
75
76# 58 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
77
78# 63 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
79
80# 68 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
81
82# 76 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
83
84# 81 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
85
86# 86 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
87
88# 91 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
89
90# 96 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
91
92# 101 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
93
94# 106 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
95
96# 111 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
97
98# 116 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
99
100# 121 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
101
102# 151 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
103
104# 192 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
105
106# 206 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
107
108# 231 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
109
110# 242 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
111
112# 244 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
113# 255 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
114
115# 284 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
116
117# 294 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
118
119# 304 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
120
121# 313 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
122
123# 330 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
124
125# 340 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
126
127# 347 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
128
129# 353 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
130
131# 359 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
132
133# 365 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
134
135# 371 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
136
137# 377 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
138! New line at end of file is required for FYPP
139# 3 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp" 2
140# 1 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp" 1
141# 1 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp" 1
142# 2 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
143# 3 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
144# 4 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
145# 5 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
146# 6 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
147
148# 8 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
149# 9 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
150# 10 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
151
152# 17 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
153
154# 46 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
155
156# 58 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
157
158# 68 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
159
160# 98 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
161
162# 110 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
163
164# 120 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
165! New line at end of file is required for FYPP
166# 2 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp" 2
167
168# 7 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
169
170# 17 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
171
172# 22 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
173
174# 27 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
175
176# 32 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
177
178# 37 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
179
180# 42 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
181
182# 47 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
183
184# 52 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
185
186# 57 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
187
188# 62 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
189
190# 73 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
191
192# 78 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
193
194# 83 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
195
196# 88 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
197
198# 103 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
199
200# 131 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
201
202# 160 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
203
204# 175 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
205
206# 193 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
207
208# 215 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
209
210# 244 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
211
212# 259 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
213
214# 269 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
215
216# 278 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
217
218# 294 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
219
220# 304 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
221
222# 311 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
223! New line at end of file is required for FYPP
224# 4 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp" 2
225
226! GPU parallel region (scalar reductions, maxval/minval)
227# 23 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
228
229! GPU parallel loop over threads (most common GPU macro)
230# 43 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
231
232! Required closing for GPU_PARALLEL_LOOP
233# 55 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
234
235! Mark routine for device compilation
236# 112 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
237
238! Declare device-resident data
239# 130 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
240
241! Inner loop within a GPU parallel region
242# 145 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
243
244! Scoped GPU data region
245# 164 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
246
247! Host code with device pointers (for MPI with GPU buffers)
248# 193 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
249
250! Allocate device memory (unscoped)
251# 207 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
252
253! Free device memory
254# 219 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
255
256! Atomic operation on device
257# 231 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
258
259! End atomic capture block
260# 242 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
261
262! Copy data between host and device
263# 254 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
264
265! Synchronization barrier
266# 266 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
267
268! Import GPU library module (openacc or omp_lib)
269# 275 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
270
271! Emit code only for AMD compiler
272# 282 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
273
274! Emit code for non-Cray compilers
275# 289 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
276
277! Emit code only for Cray compiler
278# 296 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
279
280! Emit code for non-NVIDIA compilers
281# 303 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
282
283# 305 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
284# 306 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
285! New line at end of file is required for FYPP
286# 2 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp" 2
287
288# 14 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
289
290! Caution: This macro requires the use of a binding script to set CUDA_VISIBLE_DEVICES, such that we have one GPU device per MPI
291! rank. That's because for both cudaMemAdvise (preferred location) and cudaMemPrefetchAsync we use location = device_id = 0. For an
292! example see misc/nvidia_uvm/bind.sh. NVIDIA unified memory page placement hint
293# 57 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
294
295! Allocate and create GPU device memory
296# 77 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
297
298! Free GPU device memory and deallocate
299# 85 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
300
301! Cray-specific GPU pointer setup for vector fields
302# 109 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
303
304! Cray-specific GPU pointer setup for scalar fields
305# 125 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
306
307! Cray-specific GPU pointer setup for acoustic source spatials
308# 150 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
309
310# 156 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
311
312# 163 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
313! New line at end of file is required for FYPP
314# 6 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp" 2
315# 1 "/home/runner/work/MFC/MFC/src/common/include/case.fpp" 1
316! This file exists so that Fypp can be run without generating case.fpp files for
317! each target. This is useful when generating documentation, for example. This
318! should also let MFC be built with CMake directly, without invoking mfc.sh.
319
320! For pre-process.
321# 9 "/home/runner/work/MFC/MFC/src/common/include/case.fpp"
322
323! For moving immersed boundaries in simulation
324# 14 "/home/runner/work/MFC/MFC/src/common/include/case.fpp"
325# 7 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp" 2
326
327!> @brief Writes solution data, run-time stability diagnostics (ICFL, VCFL, CCFL, Rc), and probe/center-of-mass files
329
332 use m_mpi_proxy
335 use m_helper
337 use m_sim_helpers
339 use m_ibm
341
342 implicit none
343
344 private
350
351 integer :: ib_state_unit = -1 !< I/O unit for IB state binary file
352 real(wp), allocatable, dimension(:,:,:) :: icfl_sf !< ICFL stability criterion
353 real(wp), allocatable, dimension(:,:,:) :: vcfl_sf !< VCFL stability criterion
354 real(wp), allocatable, dimension(:,:,:) :: ccfl_sf !< CCFL stability criterion
355 real(wp), allocatable, dimension(:,:,:) :: rc_sf !< Rc stability criterion
356 real(wp), public, allocatable, dimension(:,:) :: c_mass
357
358# 38 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
359#if defined(MFC_OpenACC)
360# 38 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
361!$acc declare create(icfl_sf, vcfl_sf, ccfl_sf, Rc_sf, c_mass)
362# 38 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
363#elif defined(MFC_OpenMP)
364# 38 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
365!$omp declare target (icfl_sf, vcfl_sf, ccfl_sf, Rc_sf, c_mass)
366# 38 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
367#endif
368
369 real(wp) :: icfl_max_loc, icfl_max_glb !< ICFL stability extrema on local and global grids
370 real(wp) :: vcfl_max_loc, vcfl_max_glb !< VCFL stability extrema on local and global grids
371 real(wp) :: ccfl_max_loc, ccfl_max_glb !< CCFL stability extrema on local and global grids
372 real(wp) :: rc_min_loc, rc_min_glb !< Rc stability extrema on local and global grids
373
374# 44 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
375#if defined(MFC_OpenACC)
376# 44 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
377!$acc declare create(icfl_max_loc, icfl_max_glb, vcfl_max_loc, vcfl_max_glb)
378# 44 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
379#elif defined(MFC_OpenMP)
380# 44 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
381!$omp declare target (icfl_max_loc, icfl_max_glb, vcfl_max_loc, vcfl_max_glb)
382# 44 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
383#endif
384
385# 45 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
386#if defined(MFC_OpenACC)
387# 45 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
388!$acc declare create(ccfl_max_loc, ccfl_max_glb, Rc_min_loc, Rc_min_glb)
389# 45 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
390#elif defined(MFC_OpenMP)
391# 45 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
392!$omp declare target (ccfl_max_loc, ccfl_max_glb, Rc_min_loc, Rc_min_glb)
393# 45 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
394#endif
395
396 !> @name ICFL, VCFL, CCFL and Rc stability criteria extrema over all the time-steps
397 !> @{
398 real(wp) :: icfl_max !< ICFL criterion maximum
399 real(wp) :: vcfl_max !< VCFL criterion maximum
400 real(wp) :: ccfl_max !< CCFL criterion maximum
401 real(wp) :: rc_min !< Rc criterion maximum
402 !> @}
403
404 type(scalar_field), allocatable, dimension(:) :: q_cons_temp_ds
405
406contains
407
408 !> Write data files. Dispatch subroutine that replaces procedure pointer.
409 impure subroutine s_write_data_files(q_cons_vf, q_T_sf, q_prim_vf, t_step, bc_type, beta)
410
411 type(scalar_field), dimension(sys_size), intent(inout) :: q_cons_vf
412 type(scalar_field), intent(inout) :: q_t_sf
413 type(scalar_field), dimension(sys_size), intent(inout) :: q_prim_vf
414 integer, intent(in) :: t_step
415 type(scalar_field), intent(inout), optional :: beta
416 type(integer_field), dimension(1:num_dims,-1:1), intent(in) :: bc_type
417
418 if (.not. parallel_io) then
419 call s_write_serial_data_files(q_cons_vf, q_t_sf, q_prim_vf, t_step, bc_type, beta)
420 else
421 call s_write_parallel_data_files(q_cons_vf, t_step, bc_type, beta)
422 end if
423
424 end subroutine s_write_data_files
425
426 !> Open the run-time information file and write the stability criteria table header
428
429 character(LEN=name_len), parameter :: file_name = 'run_time.inf' !< Name of the run-time information file
430 character(LEN=path_len + name_len) :: file_path !< Relative path to a file in the case directory
431 character(LEN=8) :: file_date !< Creation date of the run-time information file
432
433 file_path = trim(case_dir) // '/' // trim(file_name)
434
435 open (3, file=trim(file_path), form='formatted', status='replace')
436
437 write (3, '(A)') 'Description: Stability information at ' // 'each time-step of the simulation. This'
438 write (3, '(13X,A)') 'data is composed of the inviscid ' // 'Courant-Friedrichs-Lewy (ICFL)'
439 write (3, '(13X,A)') 'number, the viscous CFL (VCFL) number, ' // 'the capillary CFL (CCFL)'
440 write (3, '(13X,A)') 'number and the cell Reynolds (Rc) ' // 'number. Please note that only'
441 write (3, '(13X,A)') 'those stability conditions pertinent ' // 'to the physics included in'
442 write (3, '(13X,A)') 'the current computation are displayed.'
443
444 call date_and_time(date=file_date)
445
446 write (3, '(A)') 'Date: ' // file_date(5:6) // '/' // file_date(7:8) // '/' // file_date(3:4)
447
448 write (3, '(A)') ''; write (3, '(A)') ''
449
450 write (3, '(13X,A9,13X,A10,13X,A10,13X,A10)', advance="no") trim('Time-step'), trim('dt'), trim('Time'), trim('ICFL Max')
451
452 if (viscous) then
453 write (3, '(13X,A10,13X,A16)', advance="no") trim('VCFL Max'), trim('Rc Min')
454 end if
455
456 write (3, *) ! new line
457
459
460 !> Open center-of-mass data files for writing
461 impure subroutine s_open_com_files()
462
463 character(len=path_len + 3*name_len) :: file_path !< Relative path to the CoM file in the case directory
464 integer :: i !< Generic loop iterator
465
466 do i = 1, num_fluids
467 write (file_path, '(A,I0,A)') '/fluid', i, '_com.dat'
468 file_path = trim(case_dir) // trim(file_path)
469 open (i + 120, file=trim(file_path), form='formatted', position='append', status='unknown')
470 if (n == 0) then
471 write (i + 120, '(A)') ' Non-Dimensional Time ' // ' Total Mass ' // ' x-loc ' // ' Total Volume '
472 else if (p == 0) then
473 write (i + 120, &
474 & '(A)') ' Non-Dimensional Time ' // ' Total Mass ' // ' x-loc ' // ' y-loc ' &
475 & // ' Total Volume '
476 else
477 write (i + 120, &
478 & '(A)') ' Non-Dimensional Time ' // ' Total Mass ' // ' x-loc ' // ' y-loc ' // ' z-loc ' &
479 & // ' Total Volume '
480 end if
481 end do
482
483 end subroutine s_open_com_files
484
485 !> Open flow probe data files for writing
486 impure subroutine s_open_probe_files
487
488 character(LEN=path_len + 3*name_len) :: file_path !< Relative path to the probe data file in the case directory
489 integer :: i !< Generic loop iterator
490 logical :: file_exist
491
492 do i = 1, num_probes
493 write (file_path, '(A,I0,A)') '/D/probe', i, '_prim.dat'
494 file_path = trim(case_dir) // trim(file_path)
495
496 inquire (file=trim(file_path), exist=file_exist)
497
498 if (file_exist) then
499 open (i + 30, file=trim(file_path), form='formatted', status='old', position='append')
500 else
501 open (i + 30, file=trim(file_path), form='formatted', status='unknown')
502 end if
503 end do
504
505 if (integral_wrt) then
506 do i = 1, num_integrals
507 write (file_path, '(A,I0,A)') '/D/integral', i, '_prim.dat'
508 file_path = trim(case_dir) // trim(file_path)
509
510 open (i + 70, file=trim(file_path), form='formatted', position='append', status='unknown')
511 end do
512 end if
513
514 end subroutine s_open_probe_files
515
516 !> Open the immersed boundary state file for binary output
517 impure subroutine s_open_ib_state_file
518
519 character(len=path_len + 2*name_len) :: file_loc
520 integer :: ios
521
522 write (file_loc, '(A)') 'ib_state.dat'
523 file_loc = trim(case_dir) // '/D/' // trim(file_loc)
524 open (newunit=ib_state_unit, file=trim(file_loc), form='unformatted', access='stream', status='replace', iostat=ios)
525 if (ios /= 0) call s_mpi_abort('Cannot open IB state output file: ' // trim(file_loc))
526
527 end subroutine s_open_ib_state_file
528
529 !> Write stability criteria extrema to the run-time information file at the given time step
530 impure subroutine s_write_run_time_information(q_prim_vf, t_step)
531
532 type(scalar_field), dimension(sys_size), intent(in) :: q_prim_vf
533 integer, intent(in) :: t_step
534 real(wp) :: rho !< Cell-avg. density
535
536# 191 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
537 real(wp), dimension(num_fluids) :: alpha !< Cell-avg. volume fraction
538 real(wp), dimension(num_vels) :: vel !< Cell-avg. velocity
539# 194 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
540 real(wp) :: vel_sum !< Cell-avg. velocity sum
541 real(wp) :: pres !< Cell-avg. pressure
542 real(wp) :: gamma !< Cell-avg. sp. heat ratio
543 real(wp) :: pi_inf !< Cell-avg. liquid stiffness function
544 real(wp) :: qv !< Cell-avg. internal energy reference value
545 real(wp) :: c !< Cell-avg. sound speed
546 real(wp) :: h !< Cell-avg. enthalpy
547 real(wp), dimension(2) :: re !< Cell-avg. Reynolds numbers
548 integer :: j, k, l
549
550 ! Computing Stability Criteria at Current Time-step
551
552
553# 206 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
554
555# 206 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
556#if defined(MFC_OpenACC)
557# 206 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
558!$acc parallel loop collapse(3) gang vector default(present) private(j, k, l, vel, alpha, Re, rho, vel_sum, pres, gamma, pi_inf, c, H, qv)
559# 206 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
560#elif defined(MFC_OpenMP)
561# 206 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
562
563# 206 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
564
565# 206 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
566
567# 206 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
568!$omp target teams loop defaultmap(firstprivate:scalar) bind(teams,parallel) collapse(3) defaultmap(tofrom:aggregate) defaultmap(tofrom:allocatable) defaultmap(tofrom:pointer) private(j, k, l, vel, alpha, Re, rho, vel_sum, pres, gamma, pi_inf, c, H, qv)
569# 206 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
570#endif
571 do l = 0, p
572 do k = 0, n
573 do j = 0, m
574 call s_compute_enthalpy(q_prim_vf, pres, rho, gamma, pi_inf, re, h, alpha, vel, vel_sum, qv, j, k, l)
575
576 call s_compute_speed_of_sound(pres, rho, gamma, pi_inf, h, alpha, vel_sum, 0._wp, c, qv)
577
578 if (viscous) then
579 call s_compute_stability_from_dt(vel, c, rho, re, j, k, l, icfl_sf, vcfl_sf, rc_sf)
580 else
581 call s_compute_stability_from_dt(vel, c, rho, re, j, k, l, icfl_sf)
582 end if
583 end do
584 end do
585 end do
586
587# 222 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
588#if defined(MFC_OpenACC)
589# 222 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
590!$acc end parallel loop
591# 222 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
592#elif defined(MFC_OpenMP)
593# 222 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
594
595# 222 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
596!$omp end target teams loop
597# 222 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
598#endif
599
600#ifdef _CRAYFTN
601
602# 225 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
603#if defined(MFC_OpenACC)
604# 225 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
605!$acc update host(icfl_sf)
606# 225 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
607#elif defined(MFC_OpenMP)
608# 225 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
609!$omp target update from(icfl_sf)
610# 225 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
611#endif
612
613 if (viscous) then
614
615# 228 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
616#if defined(MFC_OpenACC)
617# 228 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
618!$acc update host(vcfl_sf, Rc_sf)
619# 228 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
620#elif defined(MFC_OpenMP)
621# 228 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
622!$omp target update from(vcfl_sf, Rc_sf)
623# 228 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
624#endif
625 end if
626
627 icfl_max_loc = maxval(icfl_sf)
628
629 if (viscous) then
630 vcfl_max_loc = maxval(vcfl_sf)
631 rc_min_loc = minval(rc_sf)
632 end if
633#else
634
635# 238 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
636
637# 238 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
638#if defined(MFC_OpenACC)
639# 238 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
640!$acc parallel default(present) copyin(icfl_sf) copyout(icfl_max_loc)
641# 238 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
642 icfl_max_loc = maxval(icfl_sf)
643# 238 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
644!$acc end parallel
645# 238 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
646#elif defined(MFC_OpenMP)
647# 238 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
648
649# 238 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
650
651# 238 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
652!$omp target teams defaultmap(tofrom:aggregate) defaultmap(tofrom:allocatable) defaultmap(tofrom:pointer) map(to:icfl_sf) map(from:icfl_max_loc)
653# 238 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
654 icfl_max_loc = maxval(icfl_sf)
655# 238 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
656!$omp end target teams
657# 238 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
658#else
659# 238 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
660 icfl_max_loc = maxval(icfl_sf)
661# 238 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
662#endif
663# 241 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
664 if (viscous .or. dummy) then
665
666# 242 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
667
668# 242 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
669#if defined(MFC_OpenACC)
670# 242 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
671!$acc parallel default(present) copyin(vcfl_sf, Rc_sf) copyout(vcfl_max_loc, Rc_min_loc)
672# 242 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
673 vcfl_max_loc = maxval(vcfl_sf)
674# 242 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
675 rc_min_loc = minval(rc_sf)
676# 242 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
677!$acc end parallel
678# 242 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
679#elif defined(MFC_OpenMP)
680# 242 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
681
682# 242 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
683
684# 242 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
685!$omp target teams defaultmap(tofrom:aggregate) defaultmap(tofrom:allocatable) defaultmap(tofrom:pointer) map(to:vcfl_sf, Rc_sf) map(from:vcfl_max_loc, Rc_min_loc)
686# 242 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
687 vcfl_max_loc = maxval(vcfl_sf)
688# 242 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
689 rc_min_loc = minval(rc_sf)
690# 242 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
691!$omp end target teams
692# 242 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
693#else
694# 242 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
695 vcfl_max_loc = maxval(vcfl_sf)
696# 242 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
697 rc_min_loc = minval(rc_sf)
698# 242 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
699#endif
700# 246 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
701 end if
702#endif
703
704 if (num_procs > 1) then
705 call s_mpi_reduce_stability_criteria_extrema(icfl_max_loc, vcfl_max_loc, rc_min_loc, icfl_max_glb, vcfl_max_glb, &
706 & rc_min_glb)
707 else
711 end if
712
714
715 if (viscous) then
718 end if
719
720 if (proc_rank == 0) then
721 write (3, '(13X,I9,13X,F10.6,13X,F10.6,13X,F10.6)', advance="no") t_step, dt, mytime, icfl_max_glb
722
723 if (viscous) then
724 write (3, '(13X,F10.6,13X,ES16.6)', advance="no") vcfl_max_glb, rc_min_glb
725 end if
726
727 write (3, *) ! new line
728
730 call s_mpi_abort('ICFL is NaN. Exiting.')
731 else if (icfl_max_glb > 1._wp) then
732 print *, 'icfl', icfl_max_glb
733 call s_mpi_abort('ICFL is greater than 1.0. Exiting.')
734 end if
735
736 if (viscous) then
738 call s_mpi_abort('VCFL is NaN. Exiting.')
739 else if (vcfl_max_glb > 1._wp) then
740 print *, 'vcfl', vcfl_max_glb
741 call s_mpi_abort('VCFL is greater than 1.0. Exiting.')
742 end if
743 end if
744 end if
745
746 call s_mpi_barrier()
747
748 end subroutine s_write_run_time_information
749
750 !> Write grid and conservative variable data files in serial format
751 impure subroutine s_write_serial_data_files(q_cons_vf, q_T_sf, q_prim_vf, t_step, bc_type, beta)
752
753 type(scalar_field), dimension(sys_size), intent(inout) :: q_cons_vf
754 type(scalar_field), intent(inout) :: q_t_sf
755 type(scalar_field), dimension(sys_size), intent(inout) :: q_prim_vf
756 integer, intent(in) :: t_step
757 type(scalar_field), intent(inout), optional :: beta
758 type(integer_field), dimension(1:num_dims,-1:1), intent(in) :: bc_type
759 character(LEN=path_len + 2*name_len) :: t_step_dir !< Relative path to the current time-step directory
760 character(LEN=path_len + 3*name_len) :: file_path !< Relative path to the grid and conservative variables data files
761 logical :: file_exist !< Logical used to check existence of current time-step directory
762 character(LEN=15) :: fmt
763 integer :: i, j, k, l, r
764 real(wp) :: gamma, lit_gamma, pi_inf, qv !< Temporary EOS params
765
766 write (t_step_dir, '(A,I0,A,I0)') trim(case_dir) // '/p_all'
767 write (t_step_dir, '(a,i0,a,i0)') trim(case_dir) // '/p_all/p', proc_rank, '/', t_step
768
769 file_path = trim(t_step_dir) // '/.'
770 call my_inquire(file_path, file_exist)
771 if (file_exist) call s_delete_directory(trim(t_step_dir))
772 call s_create_directory(trim(t_step_dir))
773
774 file_path = trim(t_step_dir) // '/x_cb.dat'
775
776 open (2, file=trim(file_path), form='unformatted', status='new')
777 write (2) x_cb(-1:m); close (2)
778
779 if (n > 0) then
780 file_path = trim(t_step_dir) // '/y_cb.dat'
781
782 open (2, file=trim(file_path), form='unformatted', status='new')
783 write (2) y_cb(-1:n); close (2)
784
785 if (p > 0) then
786 file_path = trim(t_step_dir) // '/z_cb.dat'
787
788 open (2, file=trim(file_path), form='unformatted', status='new')
789 write (2) z_cb(-1:p); close (2)
790 end if
791 end if
792
793 do i = 1, sys_size
794 write (file_path, '(A,I0,A)') trim(t_step_dir) // '/q_cons_vf', i, '.dat'
795
796 open (2, file=trim(file_path), form='unformatted', status='new')
797
798 write (2) q_cons_vf(i)%sf(0:m,0:n,0:p); close (2)
799 end do
800
801 ! Lagrangian beta (void fraction) written as q_cons_vf(sys_size+1) to match the parallel I/O path and allow post_process to
802 ! read it.
803 if (bubbles_lagrange) then
804 write (file_path, '(A,I0,A)') trim(t_step_dir) // '/q_cons_vf', sys_size + 1, '.dat'
805
806 open (2, file=trim(file_path), form='unformatted', status='new')
807
808 write (2) beta%sf(0:m,0:n,0:p); close (2)
809 end if
810
811 if (qbmm .and. .not. polytropic) then
812 do i = 1, nb
813 do r = 1, nnode
814 write (file_path, '(A,I0,A)') trim(t_step_dir) // '/pb', sys_size + (i - 1)*nnode + r, '.dat'
815
816 open (2, file=trim(file_path), form='unformatted', status='new')
817
818 write (2) pb_ts(1)%sf(0:m,0:n,0:p,r, i); close (2)
819 end do
820 end do
821
822 do i = 1, nb
823 do r = 1, nnode
824 write (file_path, '(A,I0,A)') trim(t_step_dir) // '/mv', sys_size + (i - 1)*nnode + r, '.dat'
825
826 open (2, file=trim(file_path), form='unformatted', status='new')
827
828 write (2) mv_ts(1)%sf(0:m,0:n,0:p,r, i); close (2)
829 end do
830 end do
831 end if
832
833 ! Writing the IB markers
834 if (ib) then
835 call s_write_serial_ib_data(t_step)
836 end if
837
838 gamma = gammas(1)
839 lit_gamma = gs_min(1)
840 pi_inf = pi_infs(1)
841 qv = qvs(1)
842
843 if (precision == 1) then
844 fmt = "(2F30.3)"
845 else
846 fmt = "(2F40.14)"
847 end if
848
849 write (t_step_dir, '(A,I0,A,I0)') trim(case_dir) // '/D'
850 file_path = trim(t_step_dir) // '/.'
851
852 inquire (file=trim(file_path), exist=file_exist)
853
854 if (.not. file_exist) call s_create_directory(trim(t_step_dir))
855
856 if ((prim_vars_wrt .or. (n == 0 .and. p == 0)) .and. (.not. igr)) then
858 do i = 1, sys_size
859
860# 404 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
861#if defined(MFC_OpenACC)
862# 404 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
863!$acc update host(q_prim_vf(i)%sf(:, :, :))
864# 404 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
865#elif defined(MFC_OpenMP)
866# 404 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
867!$omp target update from(q_prim_vf(i)%sf(:, :, :))
868# 404 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
869#endif
870 end do
871 ! q_prim_vf(bubxb) stores the value of nb needed in riemann solvers, so replace with true primitive value (=1._wp)
872 if (qbmm) then
873 q_prim_vf(bubxb)%sf = 1._wp
874 end if
875 end if
876
877 if (n == 0 .and. p == 0) then
878 if (model_eqns == 2 .and. (.not. igr)) then
879 do i = 1, sys_size
880 write (file_path, '(A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/prim.', i, '.', proc_rank, '.', t_step, '.dat'
881
882 open (2, file=trim(file_path))
883 do j = 0, m
884 ! todo: revisit change here
885 if (((i >= adv_idx%beg) .and. (i <= adv_idx%end))) then
886 write (2, fmt) x_cb(j), q_cons_vf(i)%sf(j, 0, 0)
887 else
888 write (2, fmt) x_cb(j), q_prim_vf(i)%sf(j, 0, 0)
889 end if
890 end do
891 close (2)
892 end do
893 end if
894
895 do i = 1, sys_size
896 write (file_path, '(A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/cons.', i, '.', proc_rank, '.', t_step, '.dat'
897
898 open (2, file=trim(file_path))
899 do j = 0, m
900 write (2, fmt) x_cb(j), q_cons_vf(i)%sf(j, 0, 0)
901 end do
902 close (2)
903 end do
904
905 if (qbmm .and. .not. polytropic) then
906 do i = 1, nb
907 do r = 1, nnode
908 write (file_path, '(A,I0,A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/pres.', i, '.', r, '.', proc_rank, &
909 & '.', t_step, '.dat'
910
911 open (2, file=trim(file_path))
912 do j = 0, m
913 write (2, fmt) x_cb(j), pb_ts(1)%sf(j, 0, 0, r, i)
914 end do
915 close (2)
916 end do
917 end do
918 do i = 1, nb
919 do r = 1, nnode
920 write (file_path, '(A,I0,A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/mv.', i, '.', r, '.', proc_rank, &
921 & '.', t_step, '.dat'
922
923 open (2, file=trim(file_path))
924 do j = 0, m
925 write (2, fmt) x_cb(j), mv_ts(1)%sf(j, 0, 0, r, i)
926 end do
927 close (2)
928 end do
929 end do
930 end if
931 end if
932
933 if (precision == 1) then
934 fmt = "(3F30.7)"
935 else
936 fmt = "(3F40.14)"
937 end if
938
939 if ((n > 0) .and. (p == 0)) then
940 do i = 1, sys_size
941 write (file_path, '(A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/cons.', i, '.', proc_rank, '.', t_step, '.dat'
942 open (2, file=trim(file_path))
943 do j = 0, m
944 do k = 0, n
945 write (2, fmt) x_cb(j), y_cb(k), q_cons_vf(i)%sf(j, k, 0)
946 end do
947 write (2, *)
948 end do
949 close (2)
950 end do
951
952 if (present(beta)) then
953 write (file_path, '(A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/beta.', i, '.', proc_rank, '.', t_step, '.dat'
954 open (2, file=trim(file_path))
955 do j = 0, m
956 do k = 0, n
957 write (2, fmt) x_cb(j), y_cb(k), beta%sf(j, k, 0)
958 end do
959 write (2, *)
960 end do
961 close (2)
962 end if
963
964 if (qbmm .and. .not. polytropic) then
965 do i = 1, nb
966 do r = 1, nnode
967 write (file_path, '(A,I0,A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/pres.', i, '.', r, '.', proc_rank, &
968 & '.', t_step, '.dat'
969
970 open (2, file=trim(file_path))
971 do j = 0, m
972 do k = 0, n
973 write (2, fmt) x_cb(j), y_cb(k), pb_ts(1)%sf(j, k, 0, r, i)
974 end do
975 end do
976 close (2)
977 end do
978 end do
979 do i = 1, nb
980 do r = 1, nnode
981 write (file_path, '(A,I0,A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/mv.', i, '.', r, '.', proc_rank, &
982 & '.', t_step, '.dat'
983
984 open (2, file=trim(file_path))
985 do j = 0, m
986 do k = 0, n
987 write (2, fmt) x_cb(j), y_cb(k), mv_ts(1)%sf(j, k, 0, r, i)
988 end do
989 end do
990 close (2)
991 end do
992 end do
993 end if
994
995 if (prim_vars_wrt .and. (.not. igr)) then
996 do i = 1, sys_size
997 write (file_path, '(A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/prim.', i, '.', proc_rank, '.', t_step, '.dat'
998
999 open (2, file=trim(file_path))
1000
1001 do j = 0, m
1002 do k = 0, n
1003 if (((i >= cont_idx%beg) .and. (i <= cont_idx%end)) .or. ((i >= adv_idx%beg) .and. (i <= adv_idx%end)) &
1004 & ) then
1005 write (2, fmt) x_cb(j), y_cb(k), q_cons_vf(i)%sf(j, k, 0)
1006 else
1007 write (2, fmt) x_cb(j), y_cb(k), q_prim_vf(i)%sf(j, k, 0)
1008 end if
1009 end do
1010 write (2, *)
1011 end do
1012 close (2)
1013 end do
1014 end if
1015 end if
1016
1017 if (precision == 1) then
1018 fmt = "(4F30.7)"
1019 else
1020 fmt = "(4F40.14)"
1021 end if
1022
1023 if (p > 0) then
1024 do i = 1, sys_size
1025 write (file_path, '(A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/cons.', i, '.', proc_rank, '.', t_step, '.dat'
1026 open (2, file=trim(file_path))
1027 do j = 0, m
1028 do k = 0, n
1029 do l = 0, p
1030 write (2, fmt) x_cb(j), y_cb(k), z_cb(l), q_cons_vf(i)%sf(j, k, l)
1031 end do
1032 write (2, *)
1033 end do
1034 write (2, *)
1035 end do
1036 close (2)
1037 end do
1038
1039 if (present(beta)) then
1040 write (file_path, '(A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/beta.', i, '.', proc_rank, '.', t_step, '.dat'
1041 open (2, file=trim(file_path))
1042 do j = 0, m
1043 do k = 0, n
1044 do l = 0, p
1045 write (2, fmt) x_cb(j), y_cb(k), z_cb(l), beta%sf(j, k, l)
1046 end do
1047 write (2, *)
1048 end do
1049 write (2, *)
1050 end do
1051 close (2)
1052 end if
1053
1054 if (qbmm .and. .not. polytropic) then
1055 do i = 1, nb
1056 do r = 1, nnode
1057 write (file_path, '(A,I0,A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/pres.', i, '.', r, '.', proc_rank, &
1058 & '.', t_step, '.dat'
1059
1060 open (2, file=trim(file_path))
1061 do j = 0, m
1062 do k = 0, n
1063 do l = 0, p
1064 write (2, fmt) x_cb(j), y_cb(k), z_cb(l), pb_ts(1)%sf(j, k, l, r, i)
1065 end do
1066 end do
1067 end do
1068 close (2)
1069 end do
1070 end do
1071 do i = 1, nb
1072 do r = 1, nnode
1073 write (file_path, '(A,I0,A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/mv.', i, '.', r, '.', proc_rank, &
1074 & '.', t_step, '.dat'
1075
1076 open (2, file=trim(file_path))
1077 do j = 0, m
1078 do k = 0, n
1079 do l = 0, p
1080 write (2, fmt) x_cb(j), y_cb(k), z_cb(l), mv_ts(1)%sf(j, k, l, r, i)
1081 end do
1082 end do
1083 end do
1084 close (2)
1085 end do
1086 end do
1087 end if
1088
1089 if (prim_vars_wrt .and. (.not. igr)) then
1090 do i = 1, sys_size
1091 write (file_path, '(A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/prim.', i, '.', proc_rank, '.', t_step, '.dat'
1092
1093 open (2, file=trim(file_path))
1094
1095 do j = 0, m
1096 do k = 0, n
1097 do l = 0, p
1098 if (((i >= cont_idx%beg) .and. (i <= cont_idx%end)) .or. ((i >= adv_idx%beg) &
1099 & .and. (i <= adv_idx%end)) .or. ((i >= chemxb) .and. (i <= chemxe))) then
1100 write (2, fmt) x_cb(j), y_cb(k), z_cb(l), q_cons_vf(i)%sf(j, k, l)
1101 else
1102 write (2, fmt) x_cb(j), y_cb(k), z_cb(l), q_prim_vf(i)%sf(j, k, l)
1103 end if
1104 end do
1105 write (2, *)
1106 end do
1107 write (2, *)
1108 end do
1109 close (2)
1110 end do
1111 end if
1112 end if
1113
1114 end subroutine s_write_serial_data_files
1115
1116 !> Write grid and conservative variable data files in parallel via MPI I/O
1117 impure subroutine s_write_parallel_data_files(q_cons_vf, t_step, bc_type, beta)
1118
1119 type(scalar_field), dimension(sys_size), intent(inout) :: q_cons_vf
1120 integer, intent(in) :: t_step
1121 type(scalar_field), intent(inout), optional :: beta
1122 type(integer_field), dimension(1:num_dims,-1:1), intent(in) :: bc_type
1123
1124#ifdef MFC_MPI
1125 integer :: ifile, ierr, data_size
1126 integer, dimension(MPI_STATUS_SIZE) :: status
1127 integer(kind=MPI_OFFSET_kind) :: disp
1128 integer(kind=MPI_OFFSET_kind) :: m_mok, n_mok, p_mok
1129 integer(kind=MPI_OFFSET_kind) :: wp_mok, var_mok, str_mok
1130 integer(kind=MPI_OFFSET_kind) :: nvars_mok
1131 integer(kind=MPI_OFFSET_kind) :: mok
1132 character(LEN=path_len + 2*name_len) :: file_loc
1133 logical :: file_exist, dir_check
1134 character(len=10) :: t_step_string
1135 integer :: i !< Generic loop iterator
1136 integer :: alt_sys !< Altered system size for the lagrangian subgrid bubble model
1137 ! Down sampling variables
1138 integer :: m_ds, n_ds, p_ds
1139 integer :: m_glb_ds, n_glb_ds, p_glb_ds
1140 integer :: m_glb_save, n_glb_save, p_glb_save !< Global save size
1141
1142 if (down_sample) then
1143 call s_downsample_data(q_cons_vf, q_cons_temp_ds, m_ds, n_ds, p_ds, m_glb_ds, n_glb_ds, p_glb_ds)
1144 end if
1145
1146 if (present(beta)) then
1147 alt_sys = sys_size + 1
1148 else
1149 alt_sys = sys_size
1150 end if
1151
1152 if (file_per_process) then
1153 call s_int_to_str(t_step, t_step_string)
1154
1155 if (down_sample) then
1156 call s_initialize_mpi_data_ds(q_cons_temp_ds)
1157 else
1158 if (ib) then
1159 call s_initialize_mpi_data(q_cons_vf, ib_markers)
1160 else
1161 call s_initialize_mpi_data(q_cons_vf)
1162 end if
1163 end if
1164
1165 if (proc_rank == 0) then
1166 file_loc = trim(case_dir) // '/restart_data/lustre_' // trim(t_step_string)
1167 call my_inquire(file_loc, dir_check)
1168 if (dir_check .neqv. .true.) then
1169 call s_create_directory(trim(file_loc))
1170 end if
1171 call s_create_directory(trim(file_loc))
1172 end if
1173 call s_mpi_barrier()
1175
1176 call s_initialize_mpi_data(q_cons_vf)
1177
1178 write (file_loc, '(I0,A,i7.7,A)') t_step, '_', proc_rank, '.dat'
1179 file_loc = trim(case_dir) // '/restart_data/lustre_' // trim(t_step_string) // trim(mpiiofs) // trim(file_loc)
1180 inquire (file=trim(file_loc), exist=file_exist)
1181 if (file_exist .and. proc_rank == 0) then
1182 call mpi_file_delete(file_loc, mpi_info_int, ierr)
1183 end if
1184 call mpi_file_open(mpi_comm_self, file_loc, ior(mpi_mode_wronly, mpi_mode_create), mpi_info_int, ifile, ierr)
1185
1186 if (down_sample) then
1187 data_size = (m_ds + 3)*(n_ds + 3)*(p_ds + 3)
1188 m_glb_save = m_glb_ds + 1
1189 n_glb_save = n_glb_ds + 1
1190 p_glb_save = p_glb_ds + 1
1191 else
1192 data_size = (m + 1)*(n + 1)*(p + 1)
1193 m_glb_save = m_glb + 1
1194 n_glb_save = n_glb + 1
1195 p_glb_save = p_glb + 1
1196 end if
1197
1198 m_mok = int(m_glb_save + 1, mpi_offset_kind)
1199 n_mok = int(n_glb_save + 1, mpi_offset_kind)
1200 p_mok = int(p_glb_save + 1, mpi_offset_kind)
1201 wp_mok = int(storage_size(0._stp)/8, mpi_offset_kind)
1202 mok = int(1._wp, mpi_offset_kind)
1203 str_mok = int(name_len, mpi_offset_kind)
1204 nvars_mok = int(sys_size, mpi_offset_kind)
1205
1206 if (bubbles_euler) then
1207 do i = 1, sys_size
1208 var_mok = int(i, mpi_offset_kind)
1209
1210 call mpi_file_write_all(ifile, mpi_io_data%var(i)%sf, data_size*mpi_io_type, mpi_io_p, status, ierr)
1211 end do
1212 if (qbmm .and. .not. polytropic) then
1213 do i = sys_size + 1, sys_size + 2*nb*nnode
1214 var_mok = int(i, mpi_offset_kind)
1215
1216 call mpi_file_write_all(ifile, mpi_io_data%var(i)%sf, data_size*mpi_io_type, mpi_io_p, status, ierr)
1217 end do
1218 end if
1219 else
1220 if (down_sample) then
1221 do i = 1, sys_size ! TODO: check if sys_size is correct
1222 var_mok = int(i, mpi_offset_kind)
1223
1224 call mpi_file_write_all(ifile, q_cons_temp_ds(i)%sf, data_size*mpi_io_type, mpi_io_p, status, ierr)
1225 end do
1226 else
1227 do i = 1, sys_size ! TODO: check if sys_size is correct
1228 var_mok = int(i, mpi_offset_kind)
1229
1230 call mpi_file_write_all(ifile, mpi_io_data%var(i)%sf, data_size*mpi_io_type, mpi_io_p, status, ierr)
1231 end do
1232 end if
1233 end if
1234
1235 call mpi_file_close(ifile, ierr)
1236 else
1237 if (ib) then
1238 call s_initialize_mpi_data(q_cons_vf, ib_markers)
1239 else if (present(beta)) then
1240 call s_initialize_mpi_data(q_cons_vf, beta=beta)
1241 else
1242 call s_initialize_mpi_data(q_cons_vf)
1243 end if
1244
1245 write (file_loc, '(I0,A)') t_step, '.dat'
1246 file_loc = trim(case_dir) // '/restart_data' // trim(mpiiofs) // trim(file_loc)
1247 inquire (file=trim(file_loc), exist=file_exist)
1248 if (file_exist .and. proc_rank == 0) then
1249 call mpi_file_delete(file_loc, mpi_info_int, ierr)
1250 end if
1251 call mpi_file_open(mpi_comm_world, file_loc, ior(mpi_mode_wronly, mpi_mode_create), mpi_info_int, ifile, ierr)
1252
1253 data_size = (m + 1)*(n + 1)*(p + 1)
1254
1255 m_mok = int(m_glb + 1, mpi_offset_kind)
1256 n_mok = int(n_glb + 1, mpi_offset_kind)
1257 p_mok = int(p_glb + 1, mpi_offset_kind)
1258 wp_mok = int(storage_size(0._stp)/8, mpi_offset_kind)
1259 mok = int(1._wp, mpi_offset_kind)
1260 str_mok = int(name_len, mpi_offset_kind)
1261 nvars_mok = int(alt_sys, mpi_offset_kind)
1262
1263 if (bubbles_euler) then
1264 do i = 1, sys_size
1265 var_mok = int(i, mpi_offset_kind)
1266
1267 disp = m_mok*max(mok, n_mok)*max(mok, p_mok)*wp_mok*(var_mok - 1)
1268
1269 call mpi_file_set_view(ifile, disp, mpi_p, mpi_io_data%view(i), 'native', mpi_info_int, ierr)
1270 call mpi_file_write_all(ifile, mpi_io_data%var(i)%sf, data_size*mpi_io_type, mpi_io_p, status, ierr)
1271 end do
1272 if (qbmm .and. .not. polytropic) then
1273 do i = sys_size + 1, sys_size + 2*nb*nnode
1274 var_mok = int(i, mpi_offset_kind)
1275
1276 disp = m_mok*max(mok, n_mok)*max(mok, p_mok)*wp_mok*(var_mok - 1)
1277
1278 call mpi_file_set_view(ifile, disp, mpi_p, mpi_io_data%view(i), 'native', mpi_info_int, ierr)
1279 call mpi_file_write_all(ifile, mpi_io_data%var(i)%sf, data_size*mpi_io_type, mpi_io_p, status, ierr)
1280 end do
1281 end if
1282 else
1283 do i = 1, sys_size ! TODO: check if sys_size is correct
1284 var_mok = int(i, mpi_offset_kind)
1285
1286 disp = m_mok*max(mok, n_mok)*max(mok, p_mok)*wp_mok*(var_mok - 1)
1287
1288 call mpi_file_set_view(ifile, disp, mpi_p, mpi_io_data%view(i), 'native', mpi_info_int, ierr)
1289 call mpi_file_write_all(ifile, mpi_io_data%var(i)%sf, data_size*mpi_io_type, mpi_io_p, status, ierr)
1290 end do
1291 end if
1292
1293 if (present(beta)) then
1294 var_mok = int(sys_size + 1, mpi_offset_kind)
1295
1296 disp = m_mok*max(mok, n_mok)*max(mok, p_mok)*wp_mok*(var_mok - 1)
1297
1298 call mpi_file_set_view(ifile, disp, mpi_p, mpi_io_data%view(sys_size + 1), 'native', mpi_info_int, ierr)
1299 call mpi_file_write_all(ifile, mpi_io_data%var(sys_size + 1)%sf, data_size*mpi_io_type, mpi_io_p, status, ierr)
1300 end if
1301
1302 call mpi_file_close(ifile, ierr)
1303
1304 if (ib) then
1305 call s_write_parallel_ib_data(t_step)
1306 end if
1307 end if
1308#endif
1309
1310 end subroutine s_write_parallel_data_files
1311
1312 !> Write immersed boundary marker data to a serial (per-processor) unformatted file
1313 subroutine s_write_serial_ib_data(time_step)
1314
1315 integer, intent(in) :: time_step
1316 character(LEN=path_len + 2*name_len) :: file_path
1317 character(LEN=path_len + 2*name_len) :: t_step_dir
1318
1319 write (t_step_dir, '(A,I0,A,I0)') trim(case_dir) // '/p_all'
1320 write (t_step_dir, '(a,i0,a,i0)') trim(case_dir) // '/p_all/p', proc_rank, '/', time_step
1321 write (file_path, '(A,I0,A)') trim(t_step_dir) // '/ib_data.dat'
1322
1323 open (2, file=trim(file_path), form='unformatted', status='new')
1324
1325
1326# 860 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1327#if defined(MFC_OpenACC)
1328# 860 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1329!$acc update host(ib_markers%sf)
1330# 860 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1331#elif defined(MFC_OpenMP)
1332# 860 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1333!$omp target update from(ib_markers%sf)
1334# 860 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1335#endif
1336 write (2) ib_markers%sf(0:m,0:n,0:p); close (2)
1337
1338 end subroutine s_write_serial_ib_data
1339
1340 !> Write immersed boundary marker data in parallel using MPI I/O
1341 subroutine s_write_parallel_ib_data(time_step)
1342
1343 integer, intent(in) :: time_step
1344
1345#ifdef MFC_MPI
1346 character(LEN=path_len + 2*name_len) :: file_loc
1347 integer(kind=MPI_OFFSET_kind) :: disp
1348 integer(kind=MPI_OFFSET_kind) :: m_MOK, n_MOK, p_MOK
1349 integer(kind=MPI_OFFSET_kind) :: WP_MOK, var_MOK, MOK
1350 integer :: ifile, ierr, data_size
1351 integer, dimension(MPI_STATUS_SIZE) :: status
1352
1353
1354# 878 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1355#if defined(MFC_OpenACC)
1356# 878 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1357!$acc update host(ib_markers%sf)
1358# 878 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1359#elif defined(MFC_OpenMP)
1360# 878 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1361!$omp target update from(ib_markers%sf)
1362# 878 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1363#endif
1364
1365 data_size = (m + 1)*(n + 1)*(p + 1)
1366 m_mok = int(m_glb + 1, mpi_offset_kind)
1367 n_mok = int(n_glb + 1, mpi_offset_kind)
1368 p_mok = int(p_glb + 1, mpi_offset_kind)
1369 wp_mok = int(storage_size(0._stp)/8, mpi_offset_kind)
1370 mok = int(1._wp, mpi_offset_kind)
1371
1372 write (file_loc, '(A)') 'ib.dat'
1373 file_loc = trim(case_dir) // '/restart_data' // trim(mpiiofs) // trim(file_loc)
1374 call mpi_file_open(mpi_comm_world, file_loc, ior(mpi_mode_wronly, mpi_mode_create), mpi_info_int, ifile, ierr)
1375
1376 var_mok = int(sys_size + 1, mpi_offset_kind)
1377 disp = m_mok*max(mok, n_mok)*max(mok, p_mok)*wp_mok*(var_mok - 1 + int(time_step/t_step_save))
1378 if (time_step == 0) disp = 0
1379
1380 call mpi_file_set_view(ifile, disp, mpi_integer, mpi_io_ib_data%view, 'native', mpi_info_int, ierr)
1381 call mpi_file_write_all(ifile, mpi_io_ib_data%var%sf, data_size, mpi_integer, status, ierr)
1382 call mpi_file_close(ifile, ierr)
1383#endif
1384
1385 end subroutine s_write_parallel_ib_data
1386
1387 !> Dispatch immersed boundary data output to the serial or parallel writer
1388 subroutine s_write_ib_data_file(time_step)
1389
1390 integer, intent(in) :: time_step
1391
1392 if (parallel_io) then
1393 call s_write_parallel_ib_data(time_step)
1394 else
1395 call s_write_serial_ib_data(time_step)
1396 end if
1397
1398 end subroutine s_write_ib_data_file
1399
1400 !> Write IB state records to D/ib_state.dat (rank 0 only)
1401 impure subroutine s_write_ib_state_file()
1402
1403 integer :: i
1404
1405 do i = 1, num_ibs
1406 write (ib_state_unit) mytime, i, patch_ib(i)%force, patch_ib(i)%torque, patch_ib(i)%vel, patch_ib(i)%angular_vel, &
1407 & patch_ib(i)%angles, patch_ib(i)%x_centroid, patch_ib(i)%y_centroid, patch_ib(i)%z_centroid
1408 end do
1409
1410 end subroutine s_write_ib_state_file
1411
1412 !> Write center-of-mass data at the current time step
1413 impure subroutine s_write_com_files(t_step, c_mass_in)
1414
1415 integer, intent(in) :: t_step
1416 real(wp), dimension(num_fluids, 5), intent(in) :: c_mass_in
1417 integer :: i !< Generic loop iterator
1418 real(wp) :: nondim_time !< Non-dimensional time
1419
1420 if (t_step_old /= dflt_int) then
1421 nondim_time = real(t_step + t_step_old, wp)*dt
1422 else
1423 nondim_time = real(t_step, wp)*dt
1424 end if
1425
1426 if (proc_rank == 0) then
1427 if (n == 0) then
1428 do i = 1, num_fluids
1429 write (i + 120, '(6X,4F24.12)') nondim_time, c_mass_in(i, 1), c_mass_in(i, 2), c_mass_in(i, 5)
1430 end do
1431 else if (p == 0) then
1432 do i = 1, num_fluids
1433 write (i + 120, '(6X,5F24.12)') nondim_time, c_mass_in(i, 1), c_mass_in(i, 2), c_mass_in(i, 3), c_mass_in(i, 5)
1434 end do
1435 else
1436 do i = 1, num_fluids
1437 write (i + 120, '(6X,6F24.12)') nondim_time, c_mass_in(i, 1), c_mass_in(i, 2), c_mass_in(i, 3), c_mass_in(i, &
1438 & 4), c_mass_in(i, 5)
1439 end do
1440 end if
1441 end if
1442
1443 end subroutine s_write_com_files
1444
1445 !> Write flow probe data at the current time step
1446 impure subroutine s_write_probe_files(t_step, q_cons_vf, accel_mag)
1447
1448 integer, intent(in) :: t_step
1449 type(scalar_field), dimension(sys_size), intent(in) :: q_cons_vf
1450 real(wp), dimension(0:m,0:n,0:p), intent(in) :: accel_mag
1451 real(wp), dimension(-1:m) :: distx
1452 real(wp), dimension(-1:n) :: disty
1453 real(wp), dimension(-1:p) :: distz
1454
1455 ! The cell-averaged partial densities, density, velocity, pressure, volume fractions, specific heat ratio function, liquid
1456 ! stiffness function, and sound speed.
1457 real(wp) :: lit_gamma, nbub
1458 real(wp) :: rho
1459 real(wp), dimension(num_vels) :: vel
1460 real(wp) :: pres
1461 real(wp) :: ptilde
1462 real(wp) :: ptot
1463 real(wp) :: alf
1464 real(wp) :: alfgr
1465 real(wp), dimension(num_fluids) :: alpha
1466 real(wp) :: gamma
1467 real(wp) :: pi_inf
1468 real(wp) :: qv
1469 real(wp) :: c
1470 real(wp) :: m00, m10, m01, m20, m11, m02
1471 real(wp) :: varr, varv
1472 real(wp), dimension(Nb) :: nr, r, nrdot, rdot
1473 real(wp) :: nr3
1474 real(wp) :: accel
1475 real(wp) :: int_pres
1476 real(wp) :: max_pres
1477 real(wp), dimension(2) :: re
1478 real(wp), dimension(6) :: tau_e
1479 real(wp) :: g_local
1480 real(wp) :: dyn_p, t
1481 real(wp) :: damage_state
1482 integer :: i, j, k, l, s, d !< Generic loop iterator
1483 real(wp) :: nondim_time !< Non-dimensional time
1484 real(wp) :: tmp !< Temporary variable to store quantity for mpi_allreduce
1485 integer :: npts !< Number of included integral points
1486 real(wp) :: rad, thickness !< For integral quantities
1487 logical :: trigger !< For integral quantities
1488 real(wp) :: rhoyks(1:num_species)
1489
1490 t = dflt_t_guess
1491
1492 if (time_stepper == 23) then
1493 nondim_time = mytime
1494 else
1495 if (t_step_old /= dflt_int) then
1496 nondim_time = real(t_step + t_step_old, wp)*dt
1497 else
1498 nondim_time = real(t_step, wp)*dt
1499 end if
1500 end if
1501
1502 do i = 1, num_probes
1503 rho = 0._wp
1504 do s = 1, num_vels
1505 vel(s) = 0._wp
1506 end do
1507 pres = 0._wp
1508 gamma = 0._wp
1509 pi_inf = 0._wp
1510 qv = 0._wp
1511 c = 0._wp
1512 accel = 0._wp
1513 nr = 0._wp; r = 0._wp
1514 nrdot = 0._wp; rdot = 0._wp
1515 nbub = 0._wp
1516 m00 = 0._wp
1517 m10 = 0._wp
1518 m01 = 0._wp
1519 m20 = 0._wp
1520 m11 = 0._wp
1521 m02 = 0._wp
1522 varr = 0._wp; varv = 0._wp
1523 alf = 0._wp
1524 do s = 1, (num_dims*(num_dims + 1))/2
1525 tau_e(s) = 0._wp
1526 end do
1527 damage_state = 0._wp
1528
1529 if (n == 0) then
1530 if ((probe(i)%x >= x_cb(-1)) .and. (probe(i)%x <= x_cb(m))) then
1531 do s = -1, m
1532 distx(s) = x_cb(s) - probe(i)%x
1533 if (distx(s) < 0._wp) distx(s) = 1000._wp
1534 end do
1535 j = minloc(distx, 1)
1536 if (j == 1) j = 2 ! Pick first point if probe is at edge
1537 k = 0
1538 l = 0
1539
1540 if (chemistry) then
1541 do d = 1, num_species
1542 rhoyks(d) = q_cons_vf(chemxb + d - 1)%sf(j - 2, k, l)
1543 end do
1544 end if
1545
1546 ! Computing/Sharing necessary state variables
1547 if (elasticity) then
1548 call s_convert_to_mixture_variables(q_cons_vf, j - 2, k, l, rho, gamma, pi_inf, qv, re, g_local, &
1549 & fluid_pp(:)%G)
1550 else
1551 call s_convert_to_mixture_variables(q_cons_vf, j - 2, k, l, rho, gamma, pi_inf, qv)
1552 end if
1553 do s = 1, num_vels
1554 vel(s) = q_cons_vf(cont_idx%end + s)%sf(j - 2, k, l)/rho
1555 end do
1556
1557 dyn_p = 0.5_wp*rho*dot_product(vel, vel)
1558
1559 if (elasticity) then
1560 if (cont_damage) then
1561 damage_state = q_cons_vf(damage_idx)%sf(j - 2, k, l)
1562 g_local = g_local*max((1._wp - damage_state), 0._wp)
1563 end if
1564
1565 call s_compute_pressure(q_cons_vf(1)%sf(j - 2, k, l), q_cons_vf(alf_idx)%sf(j - 2, k, l), dyn_p, pi_inf, &
1566 & gamma, rho, qv, rhoyks(:), pres, t, q_cons_vf(stress_idx%beg)%sf(j - 2, k, l), &
1567 & q_cons_vf(mom_idx%beg)%sf(j - 2, k, l), g_local)
1568 else
1569 call s_compute_pressure(q_cons_vf(e_idx)%sf(j - 2, k, l), q_cons_vf(alf_idx)%sf(j - 2, k, l), dyn_p, &
1570 & pi_inf, gamma, rho, qv, rhoyks, pres, t)
1571 end if
1572
1573 if (model_eqns == 4) then
1574 lit_gamma = gammas(1)
1575 else if (elasticity) then
1576 tau_e(1) = q_cons_vf(stress_idx%end)%sf(j - 2, k, l)/rho
1577 end if
1578
1579 if (bubbles_euler) then
1580 alf = q_cons_vf(alf_idx)%sf(j - 2, k, l)
1581 if (num_fluids == 3) then
1582 alfgr = q_cons_vf(alf_idx - 1)%sf(j - 2, k, l)
1583 end if
1584 do s = 1, nb
1585 nr(s) = q_cons_vf(bub_idx%rs(s))%sf(j - 2, k, l)
1586 nrdot(s) = q_cons_vf(bub_idx%vs(s))%sf(j - 2, k, l)
1587 end do
1588
1589 if (adv_n) then
1590 nbub = q_cons_vf(n_idx)%sf(j - 2, k, l)
1591 else
1592 nr3 = 0._wp
1593 do s = 1, nb
1594 nr3 = nr3 + weight(s)*(nr(s)**3._wp)
1595 end do
1596
1597 nbub = sqrt((4._wp*pi/3._wp)*nr3/alf)
1598 end if
1599#ifdef DEBUG
1600 print *, 'In probe, nbub: ', nbub
1601#endif
1602 if (qbmm) then
1603 m00 = q_cons_vf(bub_idx%moms(1, 1))%sf(j - 2, k, l)/nbub
1604 m10 = q_cons_vf(bub_idx%moms(1, 2))%sf(j - 2, k, l)/nbub
1605 m01 = q_cons_vf(bub_idx%moms(1, 3))%sf(j - 2, k, l)/nbub
1606 m20 = q_cons_vf(bub_idx%moms(1, 4))%sf(j - 2, k, l)/nbub
1607 m11 = q_cons_vf(bub_idx%moms(1, 5))%sf(j - 2, k, l)/nbub
1608 m02 = q_cons_vf(bub_idx%moms(1, 6))%sf(j - 2, k, l)/nbub
1609
1610 m10 = m10/m00
1611 m01 = m01/m00
1612 m20 = m20/m00
1613 m11 = m11/m00
1614 m02 = m02/m00
1615
1616 varr = m20 - m10**2._wp
1617 varv = m02 - m01**2._wp
1618 end if
1619 r(:) = nr(:)/nbub
1620 rdot(:) = nrdot(:)/nbub
1621
1622 ptilde = ptil(j - 2, k, l)
1623 ptot = pres - ptilde
1624 end if
1625
1626 ! Compute mixture sound Speed
1627 call s_compute_speed_of_sound(pres, rho, gamma, pi_inf, ((gamma + 1._wp)*pres + pi_inf)/rho, alpha, 0._wp, &
1628 & 0._wp, c, qv)
1629
1630 accel = accel_mag(j - 2, k, l)
1631 end if
1632 else if (p == 0) then
1633 if (chemistry) then
1634 do d = 1, num_species
1635 rhoyks(d) = q_cons_vf(chemxb + d - 1)%sf(j - 2, k - 2, l)
1636 end do
1637 end if
1638
1639 if ((probe(i)%x >= x_cb(-1)) .and. (probe(i)%x <= x_cb(m))) then
1640 if ((probe(i)%y >= y_cb(-1)) .and. (probe(i)%y <= y_cb(n))) then
1641 do s = -1, m
1642 distx(s) = x_cb(s) - probe(i)%x
1643 if (distx(s) < 0._wp) distx(s) = 1000._wp
1644 end do
1645 do s = -1, n
1646 disty(s) = y_cb(s) - probe(i)%y
1647 if (disty(s) < 0._wp) disty(s) = 1000._wp
1648 end do
1649 j = minloc(distx, 1)
1650 k = minloc(disty, 1)
1651 if (j == 1) j = 2 ! Pick first point if probe is at edge
1652 if (k == 1) k = 2 ! Pick first point if probe is at edge
1653 l = 0
1654
1655 ! Computing/Sharing necessary state variables
1656 call s_convert_to_mixture_variables(q_cons_vf, j - 2, k - 2, l, rho, gamma, pi_inf, qv, re, g_local, &
1657 & fluid_pp(:)%G)
1658 do s = 1, num_vels
1659 vel(s) = q_cons_vf(cont_idx%end + s)%sf(j - 2, k - 2, l)/rho
1660 end do
1661
1662 dyn_p = 0.5_wp*rho*dot_product(vel, vel)
1663
1664 if (elasticity) then
1665 if (cont_damage) then
1666 damage_state = q_cons_vf(damage_idx)%sf(j - 2, k - 2, l)
1667 g_local = g_local*max((1._wp - damage_state), 0._wp)
1668 end if
1669
1670 call s_compute_pressure(q_cons_vf(1)%sf(j - 2, k - 2, l), q_cons_vf(alf_idx)%sf(j - 2, k - 2, l), &
1671 & dyn_p, pi_inf, gamma, rho, qv, rhoyks, pres, t, &
1672 & q_cons_vf(stress_idx%beg)%sf(j - 2, k - 2, l), &
1673 & q_cons_vf(mom_idx%beg)%sf(j - 2, k - 2, l), g_local)
1674 else
1675 call s_compute_pressure(q_cons_vf(e_idx)%sf(j - 2, k - 2, l), q_cons_vf(alf_idx)%sf(j - 2, k - 2, l), &
1676 & dyn_p, pi_inf, gamma, rho, qv, rhoyks, pres, t)
1677 end if
1678
1679 if (model_eqns == 4) then
1680 lit_gamma = gs_min(1)
1681 else if (elasticity) then
1682 do s = 1, 3
1683 tau_e(s) = q_cons_vf(s)%sf(j - 2, k - 2, l)/rho
1684 end do
1685 end if
1686
1687 if (bubbles_euler) then
1688 alf = q_cons_vf(alf_idx)%sf(j - 2, k - 2, l)
1689 do s = 1, nb
1690 nr(s) = q_cons_vf(bub_idx%rs(s))%sf(j - 2, k - 2, l)
1691 nrdot(s) = q_cons_vf(bub_idx%vs(s))%sf(j - 2, k - 2, l)
1692 end do
1693
1694 if (adv_n) then
1695 nbub = q_cons_vf(n_idx)%sf(j - 2, k - 2, l)
1696 else
1697 nr3 = 0._wp
1698 do s = 1, nb
1699 nr3 = nr3 + weight(s)*(nr(s)**3._wp)
1700 end do
1701
1702 nbub = sqrt((4._wp*pi/3._wp)*nr3/alf)
1703 end if
1704
1705 r(:) = nr(:)/nbub
1706 rdot(:) = nrdot(:)/nbub
1707 end if
1708 ! Compute mixture sound speed
1709 call s_compute_speed_of_sound(pres, rho, gamma, pi_inf, ((gamma + 1._wp)*pres + pi_inf)/rho, alpha, &
1710 & 0._wp, 0._wp, c, qv)
1711 end if
1712 end if
1713 else
1714 if ((probe(i)%x >= x_cb(-1)) .and. (probe(i)%x <= x_cb(m))) then
1715 if ((probe(i)%y >= y_cb(-1)) .and. (probe(i)%y <= y_cb(n))) then
1716 if ((probe(i)%z >= z_cb(-1)) .and. (probe(i)%z <= z_cb(p))) then
1717 do s = -1, m
1718 distx(s) = x_cb(s) - probe(i)%x
1719 if (distx(s) < 0._wp) distx(s) = 1000._wp
1720 end do
1721 do s = -1, n
1722 disty(s) = y_cb(s) - probe(i)%y
1723 if (disty(s) < 0._wp) disty(s) = 1000._wp
1724 end do
1725 do s = -1, p
1726 distz(s) = z_cb(s) - probe(i)%z
1727 if (distz(s) < 0._wp) distz(s) = 1000._wp
1728 end do
1729 j = minloc(distx, 1)
1730 k = minloc(disty, 1)
1731 l = minloc(distz, 1)
1732 if (j == 1) j = 2 ! Pick first point if probe is at edge
1733 if (k == 1) k = 2 ! Pick first point if probe is at edge
1734 if (l == 1) l = 2 ! Pick first point if probe is at edge
1735
1736 ! Computing/Sharing necessary state variables
1737 call s_convert_to_mixture_variables(q_cons_vf, j - 2, k - 2, l - 2, rho, gamma, pi_inf, qv, re, &
1738 & g_local, fluid_pp(:)%G)
1739 do s = 1, num_vels
1740 vel(s) = q_cons_vf(cont_idx%end + s)%sf(j - 2, k - 2, l - 2)/rho
1741 end do
1742
1743 dyn_p = 0.5_wp*rho*dot_product(vel, vel)
1744
1745 if (chemistry) then
1746 do d = 1, num_species
1747 rhoyks(d) = q_cons_vf(chemxb + d - 1)%sf(j - 2, k - 2, l - 2)
1748 end do
1749 end if
1750
1751 if (elasticity) then
1752 if (cont_damage) then
1753 damage_state = q_cons_vf(damage_idx)%sf(j - 2, k - 2, l - 2)
1754 g_local = g_local*max((1._wp - damage_state), 0._wp)
1755 end if
1756
1757 call s_compute_pressure(q_cons_vf(1)%sf(j - 2, k - 2, l - 2), q_cons_vf(alf_idx)%sf(j - 2, k - 2, &
1758 & l - 2), dyn_p, pi_inf, gamma, rho, qv, rhoyks, pres, t, &
1759 & q_cons_vf(stress_idx%beg)%sf(j - 2, k - 2, l - 2), &
1760 & q_cons_vf(mom_idx%beg)%sf(j - 2, k - 2, l - 2), g_local)
1761 else
1762 call s_compute_pressure(q_cons_vf(e_idx)%sf(j - 2, k - 2, l - 2), q_cons_vf(alf_idx)%sf(j - 2, &
1763 & k - 2, l - 2), dyn_p, pi_inf, gamma, rho, qv, rhoyks, pres, t)
1764 end if
1765
1766 ! Compute mixture sound speed
1767 call s_compute_speed_of_sound(pres, rho, gamma, pi_inf, ((gamma + 1._wp)*pres + pi_inf)/rho, alpha, &
1768 & 0._wp, 0._wp, c, qv)
1769
1770 accel = accel_mag(j - 2, k - 2, l - 2)
1771 end if
1772 end if
1773 end if
1774 end if
1775 if (num_procs > 1) then
1776# 1292 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1777 tmp = rho
1778 call s_mpi_allreduce_sum(tmp, rho)
1779# 1292 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1780 tmp = pres
1781 call s_mpi_allreduce_sum(tmp, pres)
1782# 1292 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1783 tmp = gamma
1784 call s_mpi_allreduce_sum(tmp, gamma)
1785# 1292 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1786 tmp = pi_inf
1787 call s_mpi_allreduce_sum(tmp, pi_inf)
1788# 1292 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1789 tmp = qv
1790 call s_mpi_allreduce_sum(tmp, qv)
1791# 1292 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1792 tmp = c
1793 call s_mpi_allreduce_sum(tmp, c)
1794# 1292 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1795 tmp = accel
1796 call s_mpi_allreduce_sum(tmp, accel)
1797# 1295 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1798
1799 do s = 1, num_vels
1800 tmp = vel(s)
1801 call s_mpi_allreduce_sum(tmp, vel(s))
1802 end do
1803
1804 if (bubbles_euler) then
1805# 1303 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1806 tmp = alf
1807 call s_mpi_allreduce_sum(tmp, alf)
1808# 1303 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1809 tmp = alfgr
1810 call s_mpi_allreduce_sum(tmp, alfgr)
1811# 1303 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1812 tmp = nbub
1813 call s_mpi_allreduce_sum(tmp, nbub)
1814# 1303 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1815 tmp = nr(1)
1816 call s_mpi_allreduce_sum(tmp, nr(1))
1817# 1303 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1818 tmp = nrdot(1)
1819 call s_mpi_allreduce_sum(tmp, nrdot(1))
1820# 1303 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1821 tmp = m00
1822 call s_mpi_allreduce_sum(tmp, m00)
1823# 1303 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1824 tmp = r(1)
1825 call s_mpi_allreduce_sum(tmp, r(1))
1826# 1303 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1827 tmp = rdot(1)
1828 call s_mpi_allreduce_sum(tmp, rdot(1))
1829# 1303 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1830 tmp = ptilde
1831 call s_mpi_allreduce_sum(tmp, ptilde)
1832# 1303 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1833 tmp = ptot
1834 call s_mpi_allreduce_sum(tmp, ptot)
1835# 1306 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1836
1837 if (qbmm) then
1838# 1309 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1839 tmp = varr
1840 call s_mpi_allreduce_sum(tmp, varr)
1841# 1309 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1842 tmp = varv
1843 call s_mpi_allreduce_sum(tmp, varv)
1844# 1309 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1845 tmp = m10
1846 call s_mpi_allreduce_sum(tmp, m10)
1847# 1309 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1848 tmp = m01
1849 call s_mpi_allreduce_sum(tmp, m01)
1850# 1309 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1851 tmp = m20
1852 call s_mpi_allreduce_sum(tmp, m20)
1853# 1309 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1854 tmp = m02
1855 call s_mpi_allreduce_sum(tmp, m02)
1856# 1312 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1857 end if
1858 end if
1859
1860 if (elasticity) then
1861 do s = 1, (num_dims*(num_dims + 1))/2
1862 tmp = tau_e(s)
1863 call s_mpi_allreduce_sum(tmp, tau_e(s))
1864 end do
1865 end if
1866
1867 if (cont_damage) then
1868 tmp = damage_state
1869 call s_mpi_allreduce_sum(tmp, damage_state)
1870 end if
1871 end if
1872 if (proc_rank == 0) then
1873 if (n == 0) then
1874 if (bubbles_euler .and. (num_fluids <= 2)) then
1875 if (qbmm) then
1876 write (i + 30, '(6x,f12.6,14f28.16)') nondim_time, rho, vel(1), pres, alf, r(1), rdot(1), nr(1), &
1877 & nrdot(1), varr, varv, m10, m01, m20, m02
1878 else
1879 write (i + 30, '(6x,f12.6,8f24.8)') nondim_time, rho, vel(1), pres, alf, r(1), rdot(1), nr(1), nrdot(1)
1880 ! ptilde, & ptot
1881 end if
1882 else if (bubbles_euler .and. (num_fluids == 3)) then
1883 write (i + 30, &
1884 & '(6x,f12.6,f24.8,f24.8,f24.8,f24.8,f24.8,' // 'f24.8,f24.8,f24.8,f24.8,f24.8, f24.8)') &
1885 & nondim_time, rho, vel(1), pres, alf, alfgr, nr(1), nrdot(1), r(1), rdot(1), ptilde, ptot
1886 else if (bubbles_euler .and. num_fluids == 4) then
1887 write (i + 30, &
1888 & '(6x,f12.6,f24.8,f24.8,f24.8,f24.8,' // 'f24.8,f24.8,f24.8,f24.8,f24.8,f24.8,f24.8,f24.8,f24.8)') &
1889 & nondim_time, q_cons_vf(1)%sf(j - 2, 0, 0), q_cons_vf(2)%sf(j - 2, 0, 0), q_cons_vf(3)%sf(j - 2, &
1890 & 0, 0), q_cons_vf(4)%sf(j - 2, 0, 0), q_cons_vf(5)%sf(j - 2, 0, 0), q_cons_vf(6)%sf(j - 2, 0, 0), &
1891 & q_cons_vf(7)%sf(j - 2, 0, 0), q_cons_vf(8)%sf(j - 2, 0, 0), q_cons_vf(9)%sf(j - 2, 0, 0), &
1892 & q_cons_vf(10)%sf(j - 2, 0, 0), nbub, r(1), rdot(1)
1893 else
1894 write (i + 30, '(6X,F12.6,F24.8,F24.8,F24.8)') nondim_time, rho, vel(1), pres
1895 end if
1896 else if (p == 0) then
1897 if (bubbles_euler) then
1898# 1354 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1899 write (i + 30, '(6X,10F24.8)') nondim_time, rho, vel(1), vel(2), pres, alf, nr(1), nrdot(1), r(1), &
1900 & rdot(1)
1901# 1357 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1902 else if (elasticity) then
1903# 1359 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1904 write (i + 30, '(6X,F12.6,F24.8,F24.8,F24.8,F24.8,' // 'F24.8,F24.8,F24.8)') nondim_time, rho, &
1905 & vel(1), vel(2), pres, tau_e(1), tau_e(2), tau_e(3)
1906# 1362 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1907 else
1908 write (i + 30, '(6X,F12.6,F24.8,F24.8,F24.8)') nondim_time, rho, vel(1), pres
1909 print *, 'time =', nondim_time, 'rho =', rho, 'pres =', pres
1910 end if
1911 else
1912# 1368 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1913 write (i + 30, &
1914 & '(6X,F12.6,F24.8,F24.8,F24.8,F24.8,' // 'F24.8,F24.8,F24.8,F24.8,F24.8,' // 'F24.8)') &
1915 & nondim_time, rho, vel(1), vel(2), vel(3), pres, gamma, pi_inf, qv, c, accel
1916# 1372 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
1917 end if
1918 end if
1919 end do
1920
1921 if (integral_wrt .and. bubbles_euler) then
1922 if (n == 0) then
1923 do i = 1, num_integrals
1924 int_pres = 0._wp
1925 max_pres = 0._wp
1926 k = 0; l = 0
1927 npts = 0
1928 do j = 1, m
1929 pres = 0._wp
1930 do s = 1, num_vels
1931 vel(s) = 0._wp
1932 end do
1933 rho = 0._wp
1934 pres = 0._wp
1935 gamma = 0._wp
1936 pi_inf = 0._wp
1937 qv = 0._wp
1938
1939 if ((integral(i)%xmin <= x_cb(j)) .and. (integral(i)%xmax >= x_cb(j))) then
1940 npts = npts + 1
1941 call s_convert_to_mixture_variables(q_cons_vf, j, k, l, rho, gamma, pi_inf, qv, re)
1942 do s = 1, num_vels
1943 vel(s) = q_cons_vf(cont_idx%end + s)%sf(j, k, l)/rho
1944 end do
1945
1946 pres = ((q_cons_vf(e_idx)%sf(j, k, l) - 0.5_wp*(q_cons_vf(mom_idx%beg)%sf(j, k, &
1947 & l)**2._wp)/rho)/(1._wp - q_cons_vf(alf_idx)%sf(j, k, l)) - pi_inf - qv)/gamma
1948 int_pres = int_pres + (pres - 1._wp)**2._wp
1949 end if
1950 end do
1951 int_pres = sqrt(int_pres/(1._wp*npts))
1952
1953 if (num_procs > 1) then
1954 tmp = int_pres
1955 call s_mpi_allreduce_sum(tmp, int_pres)
1956 end if
1957
1958 if (proc_rank == 0) then
1959 if (bubbles_euler .and. (num_fluids <= 2)) then
1960 write (i + 70, '(6x,f12.6,f24.8)') nondim_time, int_pres
1961 end if
1962 end if
1963 end do
1964 else if (p == 0) then
1965 if (num_integrals /= 3) then
1966 call s_mpi_abort('Incorrect number of integrals')
1967 end if
1968
1969 rad = integral(1)%xmax
1970 thickness = integral(1)%xmin
1971
1972 do i = 1, num_integrals
1973 int_pres = 0._wp
1974 max_pres = 0._wp
1975 l = 0
1976 npts = 0
1977 do j = 1, m
1978 do k = 1, n
1979 trigger = .false.
1980 if (i == 1) then
1981 ! inner portion
1982 if (sqrt(x_cb(j)**2._wp + y_cb(k)**2._wp) < (rad - 0.5_wp*thickness)) trigger = .true.
1983 else if (i == 2) then
1984 ! net region
1985 if (sqrt(x_cb(j)**2._wp + y_cb(k)**2._wp) > (rad - 0.5_wp*thickness) .and. sqrt(x_cb(j)**2._wp &
1986 & + y_cb(k)**2._wp) < (rad + 0.5_wp*thickness)) trigger = .true.
1987 else if (i == 3) then
1988 ! everything else
1989 if (sqrt(x_cb(j)**2._wp + y_cb(k)**2._wp) > (rad + 0.5_wp*thickness)) trigger = .true.
1990 end if
1991
1992 pres = 0._wp
1993 do s = 1, num_vels
1994 vel(s) = 0._wp
1995 end do
1996 rho = 0._wp
1997 pres = 0._wp
1998 gamma = 0._wp
1999 pi_inf = 0._wp
2000 qv = 0._wp
2001
2002 if (trigger) then
2003 npts = npts + 1
2004 call s_convert_to_mixture_variables(q_cons_vf, j, k, l, rho, gamma, pi_inf, qv, re)
2005 do s = 1, num_vels
2006 vel(s) = q_cons_vf(cont_idx%end + s)%sf(j, k, l)/rho
2007 end do
2008
2009 pres = ((q_cons_vf(e_idx)%sf(j, k, l) - 0.5_wp*(q_cons_vf(mom_idx%beg)%sf(j, k, &
2010 & l)**2._wp)/rho)/(1._wp - q_cons_vf(alf_idx)%sf(j, k, l)) - pi_inf - qv)/gamma
2011 int_pres = int_pres + abs(pres - 1._wp)
2012 max_pres = max(max_pres, abs(pres - 1._wp))
2013 end if
2014 end do
2015 end do
2016
2017 if (npts > 0) then
2018 int_pres = int_pres/(1._wp*npts)
2019 else
2020 int_pres = 0._wp
2021 end if
2022
2023 if (num_procs > 1) then
2024 tmp = int_pres
2025 call s_mpi_allreduce_sum(tmp, int_pres)
2026
2027 tmp = max_pres
2028 call s_mpi_allreduce_max(tmp, max_pres)
2029 end if
2030
2031 if (proc_rank == 0) then
2032 if (bubbles_euler .and. (num_fluids <= 2)) then
2033 write (i + 70, '(6x,f12.6,f24.8,f24.8)') nondim_time, int_pres, max_pres
2034 end if
2035 end if
2036 end do
2037 end if
2038 end if
2039
2040 end subroutine s_write_probe_files
2041
2042 !> Write footer with stability criteria extrema and run-time to the information file, then close it
2044
2045 real(wp) :: run_time !< Run-time of the simulation
2046
2047 write (3, '(A)') ' '
2048 write (3, '(A)') ''
2049
2050 write (3, '(A,F9.6)') 'ICFL Max: ', icfl_max
2051 if (viscous) write (3, '(A,F9.6)') 'VCFL Max: ', vcfl_max
2052 if (viscous) write (3, '(A,F10.6)') 'Rc Min: ', rc_min
2053
2054 call cpu_time(run_time)
2055
2056 write (3, '(A)') ''
2057 write (3, '(A,I0,A)') 'Run-time: ', int(anint(run_time)), 's'
2058 write (3, '(A)') ' '
2059 close (3)
2060
2062
2063 !> Closes communication files
2064 impure subroutine s_close_com_files()
2065
2066 integer :: i !< Generic loop iterator
2067
2068 do i = 1, num_fluids
2069 close (i + 120)
2070 end do
2071
2072 end subroutine s_close_com_files
2073
2074 !> Closes probe files
2075 impure subroutine s_close_probe_files
2076
2077 integer :: i !< Generic loop iterator
2078
2079 do i = 1, num_probes
2080 close (i + 30)
2081 end do
2082
2083 end subroutine s_close_probe_files
2084
2085 !> Close the immersed boundary state file
2086 impure subroutine s_close_ib_state_file
2087
2088 close (ib_state_unit)
2089
2090 end subroutine s_close_ib_state_file
2091
2092 !> Initialize the data output module
2094
2095 integer :: i, m_ds, n_ds, p_ds
2096
2097 if (run_time_info) then
2098#ifdef MFC_DEBUG
2099# 1553 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2100 block
2101# 1553 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2102 use iso_fortran_env, only: output_unit
2103# 1553 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2104
2105# 1553 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2106 print *, 'm_data_output.fpp:1553: ', '@:ALLOCATE(icfl_sf(0:m, 0:n, 0:p))'
2107# 1553 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2108
2109# 1553 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2110 call flush (output_unit)
2111# 1553 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2112 end block
2113# 1553 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2114#endif
2115# 1553 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2116 allocate (icfl_sf(0:m, 0:n, 0:p))
2117# 1553 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2118
2119# 1553 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2120
2121# 1553 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2122#if defined(MFC_OpenACC)
2123# 1553 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2124!$acc enter data create(icfl_sf)
2125# 1553 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2126#elif defined(MFC_OpenMP)
2127# 1553 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2128!$omp target enter data map(always,alloc:icfl_sf)
2129# 1553 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2130#endif
2131 icfl_max = 0._wp
2132
2133 if (viscous) then
2134#ifdef MFC_DEBUG
2135# 1557 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2136 block
2137# 1557 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2138 use iso_fortran_env, only: output_unit
2139# 1557 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2140
2141# 1557 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2142 print *, 'm_data_output.fpp:1557: ', '@:ALLOCATE(vcfl_sf(0:m, 0:n, 0:p))'
2143# 1557 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2144
2145# 1557 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2146 call flush (output_unit)
2147# 1557 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2148 end block
2149# 1557 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2150#endif
2151# 1557 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2152 allocate (vcfl_sf(0:m, 0:n, 0:p))
2153# 1557 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2154
2155# 1557 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2156
2157# 1557 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2158#if defined(MFC_OpenACC)
2159# 1557 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2160!$acc enter data create(vcfl_sf)
2161# 1557 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2162#elif defined(MFC_OpenMP)
2163# 1557 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2164!$omp target enter data map(always,alloc:vcfl_sf)
2165# 1557 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2166#endif
2167#ifdef MFC_DEBUG
2168# 1558 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2169 block
2170# 1558 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2171 use iso_fortran_env, only: output_unit
2172# 1558 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2173
2174# 1558 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2175 print *, 'm_data_output.fpp:1558: ', '@:ALLOCATE(Rc_sf (0:m, 0:n, 0:p))'
2176# 1558 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2177
2178# 1558 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2179 call flush (output_unit)
2180# 1558 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2181 end block
2182# 1558 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2183#endif
2184# 1558 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2185 allocate (rc_sf(0:m, 0:n, 0:p))
2186# 1558 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2187
2188# 1558 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2189
2190# 1558 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2191#if defined(MFC_OpenACC)
2192# 1558 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2193!$acc enter data create(Rc_sf)
2194# 1558 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2195#elif defined(MFC_OpenMP)
2196# 1558 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2197!$omp target enter data map(always,alloc:Rc_sf)
2198# 1558 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2199#endif
2200
2201 vcfl_max = 0._wp
2202 rc_min = 1.e3_wp
2203 end if
2204 end if
2205
2206 if (probe_wrt) then
2207#ifdef MFC_DEBUG
2208# 1566 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2209 block
2210# 1566 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2211 use iso_fortran_env, only: output_unit
2212# 1566 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2213
2214# 1566 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2215 print *, 'm_data_output.fpp:1566: ', '@:ALLOCATE(c_mass(num_fluids,5))'
2216# 1566 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2217
2218# 1566 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2219 call flush (output_unit)
2220# 1566 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2221 end block
2222# 1566 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2223#endif
2224# 1566 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2225 allocate (c_mass(num_fluids,5))
2226# 1566 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2227
2228# 1566 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2229
2230# 1566 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2231#if defined(MFC_OpenACC)
2232# 1566 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2233!$acc enter data create(c_mass)
2234# 1566 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2235#elif defined(MFC_OpenMP)
2236# 1566 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2237!$omp target enter data map(always,alloc:c_mass)
2238# 1566 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2239#endif
2240 end if
2241
2242 if (down_sample) then
2243 m_ds = int((m + 1)/3) - 1
2244 n_ds = int((n + 1)/3) - 1
2245 p_ds = int((p + 1)/3) - 1
2246
2247 allocate (q_cons_temp_ds(1:sys_size))
2248 do i = 1, sys_size
2249 allocate (q_cons_temp_ds(i)%sf(-1:m_ds + 1,-1:n_ds + 1,-1:p_ds + 1))
2250 end do
2251 end if
2252
2253 end subroutine s_initialize_data_output_module
2254
2255 !> Module deallocation and/or disassociation procedures
2257
2258 integer :: i
2259
2260 if (probe_wrt) then
2261#ifdef MFC_DEBUG
2262# 1588 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2263 block
2264# 1588 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2265 use iso_fortran_env, only: output_unit
2266# 1588 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2267
2268# 1588 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2269 print *, 'm_data_output.fpp:1588: ', '@:DEALLOCATE(c_mass)'
2270# 1588 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2271
2272# 1588 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2273 call flush (output_unit)
2274# 1588 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2275 end block
2276# 1588 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2277#endif
2278# 1588 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2279
2280# 1588 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2281#if defined(MFC_OpenACC)
2282# 1588 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2283!$acc exit data delete(c_mass)
2284# 1588 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2285#elif defined(MFC_OpenMP)
2286# 1588 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2287!$omp target exit data map(release:c_mass)
2288# 1588 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2289#endif
2290# 1588 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2291 deallocate (c_mass)
2292 end if
2293
2294 if (run_time_info) then
2295#ifdef MFC_DEBUG
2296# 1592 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2297 block
2298# 1592 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2299 use iso_fortran_env, only: output_unit
2300# 1592 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2301
2302# 1592 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2303 print *, 'm_data_output.fpp:1592: ', '@:DEALLOCATE(icfl_sf)'
2304# 1592 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2305
2306# 1592 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2307 call flush (output_unit)
2308# 1592 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2309 end block
2310# 1592 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2311#endif
2312# 1592 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2313
2314# 1592 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2315#if defined(MFC_OpenACC)
2316# 1592 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2317!$acc exit data delete(icfl_sf)
2318# 1592 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2319#elif defined(MFC_OpenMP)
2320# 1592 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2321!$omp target exit data map(release:icfl_sf)
2322# 1592 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2323#endif
2324# 1592 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2325 deallocate (icfl_sf)
2326 if (viscous) then
2327#ifdef MFC_DEBUG
2328# 1594 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2329 block
2330# 1594 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2331 use iso_fortran_env, only: output_unit
2332# 1594 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2333
2334# 1594 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2335 print *, 'm_data_output.fpp:1594: ', '@:DEALLOCATE(vcfl_sf, Rc_sf)'
2336# 1594 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2337
2338# 1594 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2339 call flush (output_unit)
2340# 1594 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2341 end block
2342# 1594 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2343#endif
2344# 1594 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2345
2346# 1594 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2347#if defined(MFC_OpenACC)
2348# 1594 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2349!$acc exit data delete(vcfl_sf, Rc_sf)
2350# 1594 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2351#elif defined(MFC_OpenMP)
2352# 1594 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2353!$omp target exit data map(release:vcfl_sf, Rc_sf)
2354# 1594 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2355#endif
2356# 1594 "/home/runner/work/MFC/MFC/src/simulation/m_data_output.fpp"
2357 deallocate (vcfl_sf, rc_sf)
2358 end if
2359 end if
2360
2361 if (down_sample) then
2362 do i = 1, sys_size
2363 deallocate (q_cons_temp_ds(i)%sf)
2364 end do
2365 deallocate (q_cons_temp_ds)
2366 end if
2367
2368 end subroutine s_finalize_data_output_module
2369
2370end module m_data_output
type(scalar_field), dimension(sys_size), intent(inout) q_cons_vf
integer, intent(in) k
integer, intent(in) j
integer, intent(in) l
Noncharacteristic and processor boundary condition application for ghost cells and buffer regions.
Platform-specific file and directory operations: create, delete, inquire, getcwd, and basename.
impure subroutine s_delete_directory(dir_name)
Recursively delete a directory using a platform-specific system command.
impure subroutine my_inquire(fileloc, dircheck)
Inquires on the existence of a directory.
impure subroutine s_create_directory(dir_name)
Create a directory and all its parents if it does not exist.
Writes solution data, run-time stability diagnostics (ICFL, VCFL, CCFL, Rc), and probe/center-of-mass...
real(wp), dimension(:,:), allocatable, public c_mass
real(wp), dimension(:,:,:), allocatable ccfl_sf
CCFL stability criterion.
impure subroutine, public s_open_probe_files
Open flow probe data files for writing.
impure subroutine, public s_write_ib_state_file()
Write IB state records to D/ib_state.dat (rank 0 only).
real(wp) rc_min
Rc criterion maximum.
real(wp) vcfl_max
VCFL criterion maximum.
real(wp), dimension(:,:,:), allocatable vcfl_sf
VCFL stability criterion.
impure subroutine, public s_write_probe_files(t_step, q_cons_vf, accel_mag)
Write flow probe data at the current time step.
real(wp) rc_min_glb
Rc stability extrema on local and global grids.
impure subroutine, public s_write_com_files(t_step, c_mass_in)
Write center-of-mass data at the current time step.
impure subroutine, public s_write_parallel_data_files(q_cons_vf, t_step, bc_type, beta)
Write grid and conservative variable data files in parallel via MPI I/O.
impure subroutine, public s_initialize_data_output_module
Initialize the data output module.
impure subroutine, public s_finalize_data_output_module
Module deallocation and/or disassociation procedures.
subroutine s_write_serial_ib_data(time_step)
Write immersed boundary marker data to a serial (per-processor) unformatted file.
impure subroutine, public s_close_run_time_information_file
Write footer with stability criteria extrema and run-time to the information file,...
real(wp) vcfl_max_glb
VCFL stability extrema on local and global grids.
integer ib_state_unit
I/O unit for IB state binary file.
impure subroutine, public s_write_data_files(q_cons_vf, q_t_sf, q_prim_vf, t_step, bc_type, beta)
Write data files. Dispatch subroutine that replaces procedure pointer.
impure subroutine, public s_write_serial_data_files(q_cons_vf, q_t_sf, q_prim_vf, t_step, bc_type, beta)
Write grid and conservative variable data files in serial format.
real(wp) ccfl_max_glb
CCFL stability extrema on local and global grids.
real(wp), dimension(:,:,:), allocatable rc_sf
Rc stability criterion.
impure subroutine, public s_open_ib_state_file
Open the immersed boundary state file for binary output.
impure subroutine, public s_close_probe_files
Closes probe files.
impure subroutine, public s_close_com_files()
Closes communication files.
real(wp) icfl_max_glb
ICFL stability extrema on local and global grids.
impure subroutine, public s_write_run_time_information(q_prim_vf, t_step)
Write stability criteria extrema to the run-time information file at the given time step.
impure subroutine, public s_open_run_time_information_file
Open the run-time information file and write the stability criteria table header.
real(wp) icfl_max
ICFL criterion maximum.
type(scalar_field), dimension(:), allocatable q_cons_temp_ds
subroutine s_write_parallel_ib_data(time_step)
Write immersed boundary marker data in parallel using MPI I/O.
impure subroutine, public s_close_ib_state_file
Close the immersed boundary state file.
real(wp), dimension(:,:,:), allocatable icfl_sf
ICFL stability criterion.
impure subroutine, public s_open_com_files()
Open center-of-mass data files for writing.
subroutine, public s_write_ib_data_file(time_step)
Dispatch immersed boundary data output to the serial or parallel writer.
real(wp) ccfl_max
CCFL criterion maximum.
Rank-staggered file access delays to prevent I/O contention on parallel file systems.
impure subroutine, public delayfileaccess(processrank)
Introduce a rank-dependent busy-wait delay to stagger parallel file access and reduce I/O contention.
Shared derived types for field data, patch geometry, bubble dynamics, and MPI I/O structures.
Global parameters for the computational domain, fluid properties, and simulation algorithm configurat...
real(wp) mytime
Current simulation time.
logical bubbles_euler
Bubbles euler on/off.
logical cont_damage
Continuum damage modeling.
logical igr
Use information geometric regularization.
type(int_bounds_info), dimension(1:3) idwint
real(wp), dimension(:), allocatable, target z_cb
logical, parameter chemistry
Chemistry modeling.
type(ib_patch_parameters), dimension(num_patches_max) patch_ib
Immersed boundary patch parameters.
type(int_bounds_info) mom_idx
Indexes of first & last momentum eqns.
integer num_fluids
number of fluids in the simulation
type(int_bounds_info) stress_idx
Indexes of first and last shear stress eqns.
integer proc_rank
Rank of the local processor.
integer n_idx
Index of number density.
type(mpi_io_ib_var), public mpi_io_ib_data
logical dummy
AMDFlang workaround for case-optimization + GPU-kernel bug.
character(len=name_len) mpiiofs
type(vec3_dt), dimension(num_probes_max) probe
integer sys_size
Number of unknowns in system of eqns.
real(wp), dimension(:), allocatable weight
Simpson quadrature weights.
integer t_step_old
Existing IC/grid folder.
type(int_bounds_info) cont_idx
Indexes of first & last continuity eqns.
logical viscous
Viscous effects.
logical run_time_info
Run-time output flag.
integer model_eqns
Multicomponent flow model.
integer precision
Precision of output files.
type(physical_parameters), dimension(num_fluids_max) fluid_pp
Stiffened gas EOS parameters and Reynolds numbers per fluid.
integer num_dims
Number of spatial dimensions.
type(pres_field), dimension(:), allocatable pb_ts
type(pres_field), dimension(:), allocatable mv_ts
logical bubbles_lagrange
Lagrangian subgrid bubble model switch.
integer num_vels
Number of velocity components (different from num_dims for mhd).
logical polytropic
Polytropic switch.
logical qbmm
Quadrature moment method.
type(bub_bounds_info) bub_idx
Indexes of first & last bubble variable eqns.
integer damage_idx
Index of damage state variable (D) for continuum damage model.
real(wp), dimension(:), allocatable qvs
real(wp), dimension(:), allocatable pi_infs
logical adv_n
Solve the number density equation and compute alpha from number density.
integer num_procs
Number of processors.
character(len=path_len) case_dir
Case folder location.
type(int_bounds_info) adv_idx
Indexes of first & last advection eqns.
logical parallel_io
Format of the data files.
type(integral_parameters), dimension(num_probes_max) integral
real(wp), dimension(:), allocatable, target y_cb
real(wp), dimension(:,:,:), allocatable ptil
Pressure modification.
integer e_idx
Index of energy equation.
logical down_sample
down sample the output files
logical file_per_process
shared file or not when using parallel io
logical elasticity
elasticity modeling, true for hyper or hypo
integer nb
Number of eq. bubble sizes.
type(mpi_io_var), public mpi_io_data
real(wp) dt
Size of the time-step.
integer time_stepper
Time-stepper algorithm.
real(wp), dimension(:), allocatable gammas
real(wp), dimension(:), allocatable gs_min
integer alf_idx
Index of void fraction.
real(wp), dimension(:), allocatable, target x_cb
Basic floating-point utilities: approximate equality, default detection, and coordinate bounds.
logical elemental function, public f_approx_equal(a, b, tol_input)
Check if two floating point numbers of wp are within tolerance.
Utility routines for bubble model setup, coordinate transforms, array sampling, and special functions...
subroutine, public s_downsample_data(q_cons_vf, q_cons_temp, m_ds, n_ds, p_ds, m_glb_ds, n_glb_ds, p_glb_ds)
Downsample conservative variable fields by a factor of 3 in each direction using volume averaging.
elemental subroutine, public s_int_to_str(i, res)
Convert an integer to its trimmed string representation.
Ghost-node immersed boundary method: locates ghost/image points, computes interpolation coefficients,...
type(integer_field), public ib_markers
MPI halo exchange, domain decomposition, and buffer packing/unpacking for the simulation solver.
Simulation helper routines for enthalpy computation, CFL calculation, and stability checks.
subroutine, public s_compute_enthalpy(q_prim_vf, pres, rho, gamma, pi_inf, re, h, alpha, vel, vel_sum, qv, j, k, l)
Computes enthalpy.
subroutine, public s_compute_stability_from_dt(vel, c, rho, re_l, j, k, l, icfl_sf, vcfl_sf, rc_sf)
Computes stability criterion for a specified dt.
Conservative-to-primitive variable conversion, mixture property evaluation, and pressure computation.
subroutine s_compute_speed_of_sound(pres, rho, gamma, pi_inf, h, adv, vel_sum, c_c, c, qv)
Compute the speed of sound from thermodynamic state variables, supporting multiple equation-of-state ...
subroutine, public s_compute_pressure(energy, alf, dyn_p, pi_inf, gamma, rho, qv, rhoyks, pres, t, stress, mom, g, pres_mag)
Compute the pressure from the appropriate equation of state.
subroutine, public s_convert_conservative_to_primitive_variables(qk_cons_vf, q_t_sf, qk_prim_vf, ibounds)
Convert conserved variables (rho*alpha, rho*u, E, alpha) to primitives (rho, u, p,...
subroutine, public s_convert_to_mixture_variables(q_vf, i, j, k, rho, gamma, pi_inf, qv, re_k, g_k, g)
Dispatch to the s_convert_mixture_to_mixture_variables and s_convert_species_to_mixture_variables sub...
Derived type annexing an integer scalar field (SF).
Derived type annexing a scalar field (SF).