forked from xcompact3d/Incompact3d
-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathcase.f90
582 lines (389 loc) · 17.6 KB
/
case.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
!Copyright (c) 2012-2022, Xcompact3d
!This file is part of Xcompact3d (xcompact3d.com)
!SPDX-License-Identifier: BSD 3-Clause
module case
use param
use decomp_2d_constants
use variables
use case_generic
use tgv
use cyl
use hill
use channel
use mixlayer
use gravitycur
use tbl
use abl
use uniform
use sandbox
use cavity
use pipe
use var, only : nzmsize
implicit none
logical :: case_visu_init = .false.
private ! All functions/subroutines private by default
public :: init, boundary_conditions, &
momentum_forcing, scalar_forcing, set_fluid_properties, &
test_flow, preprocessing, postprocessing, visu_case, visu_case_init
contains
!##################################################################
subroutine init (rho1, ux1, uy1, uz1, ep1, phi1, drho1, dux1, duy1, duz1, dphi1, &
pp3, px1, py1, pz1)
real(mytype),dimension(xsize(1),xsize(2),xsize(3)) :: ux1,uy1,uz1,ep1
real(mytype),dimension(xsize(1),xsize(2),xsize(3),nrhotime) :: rho1
real(mytype),dimension(xsize(1),xsize(2),xsize(3),numscalar) :: phi1
real(mytype),dimension(xsize(1),xsize(2),xsize(3),ntime) :: dux1,duy1,duz1,drho1
real(mytype),dimension(xsize(1),xsize(2),xsize(3),ntime,numscalar) :: dphi1
real(mytype),dimension(ph1%zst(1):ph1%zen(1), ph1%zst(2):ph1%zen(2), nzmsize, npress) :: pp3
real(mytype),dimension(xsize(1),xsize(2),xsize(3)) :: px1, py1, pz1
INTEGER :: it, is
!! Zero out the pressure field
pp3(:,:,:,1) = zero
px1(:,:,:) = zero
py1(:,:,:) = zero
pz1(:,:,:) = zero
!! Default density and pressure0 to one
pressure0 = one
rho1(:,:,:,:) = one
if (itype.eq.itype_generic) then
call init_generic (ux1, uy1, uz1, ep1, phi1)
elseif (itype.eq.itype_gravitycur) then
call init_gravitycur(rho1, ux1, uy1, uz1, ep1, phi1)
elseif (itype.eq.itype_tgv) then
call init_tgv (ux1, uy1, uz1, ep1, phi1)
elseif (itype.eq.itype_channel) then
call init_channel (ux1, uy1, uz1, ep1, phi1)
elseif (itype.eq.itype_hill) then
call init_hill (ux1,uy1,uz1,ep1,phi1)
elseif (itype.eq.itype_cyl) then
call init_cyl (ux1, uy1, uz1, phi1)
elseif (itype.eq.itype_mixlayer) then
call init_mixlayer(rho1, ux1, uy1, uz1)
elseif (itype.eq.itype_tbl) then
call init_tbl (ux1, uy1, uz1, ep1, phi1)
elseif (itype.eq.itype_abl) then
call init_abl (ux1, uy1, uz1, ep1, phi1)
elseif (itype.eq.itype_uniform) then
call init_uniform (ux1, uy1, uz1, ep1, phi1)
elseif (itype.EQ.itype_sandbox) THEN
call init_sandbox (ux1, uy1, uz1, ep1, phi1, 0)
elseif (itype.eq.itype_cavity) then
call init_cavity(ux1, uy1, uz1, ep1, phi1)
elseif (itype.eq.itype_pipe) then
call init_pipe(ux1, uy1, uz1, ep1, phi1)
else
if (nrank.eq.0) then
print *, "ERROR: Unknown itype: ", itype
STOP
endif
endif
!! Setup old arrays
do it = 1, ntime
drho1(:,:,:,it) = rho1(:,:,:,1)
dux1(:,:,:,it)=ux1(:,:,:)
duy1(:,:,:,it)=uy1(:,:,:)
duz1(:,:,:,it)=uz1(:,:,:)
enddo
do it = 2, nrhotime
rho1(:,:,:,it) = rho1(:,:,:,1)
enddo
do is = 1, numscalar
do it = 1, ntime
dphi1(:,:,:,it,is) = phi1(:,:,:,is)
enddo
enddo
end subroutine init
!##################################################################
!##################################################################
subroutine boundary_conditions (rho,ux,uy,uz,phi,ep)
real(mytype),dimension(xsize(1),xsize(2),xsize(3)) :: ux,uy,uz,ep
real(mytype),dimension(xsize(1),xsize(2),xsize(3),numscalar) :: phi
real(mytype),dimension(xsize(1),xsize(2),xsize(3),nrhotime) :: rho
if (itype == itype_generic) then
call boundary_conditions_generic (ux,uy,uz,phi,ep)
elseif (itype.eq.itype_gravitycur) then
call boundary_conditions_gravitycur(rho, phi)
elseif (itype.eq.itype_tgv) then
call boundary_conditions_tgv (ux, uy, uz, phi)
elseif (itype.eq.itype_channel) then
call boundary_conditions_channel (ux, uy, uz, phi)
elseif (itype.eq.itype_hill) then
call boundary_conditions_hill (ux,uy,uz,phi,ep)
elseif (itype.eq.itype_cyl) then
call boundary_conditions_cyl (ux, uy, uz, phi)
elseif (itype.eq.itype_tbl) then
call boundary_conditions_tbl (ux, uy, uz, phi)
elseif (itype.eq.itype_abl) then
call boundary_conditions_abl (ux, uy, uz, phi)
elseif (itype.eq.itype_uniform) then
call boundary_conditions_uniform (ux, uy, uz, phi)
elseif (itype.EQ.itype_sandbox) THEN
call boundary_conditions_sandbox (ux, uy, uz, phi)
elseif (itype.eq.itype_cavity) then
call boundary_conditions_cavity(ux, uy, uz, phi)
elseif (itype.eq.itype_pipe) then
call boundary_conditions_pipe (ux, uy, uz, phi)
endif
end subroutine boundary_conditions
!##################################################################
!##################################################################
subroutine preprocessing(rho1, ux1, uy1, uz1, pp3, phi1, ep1)
use decomp_2d, only : xsize, ph1
use visu, only : write_snapshot
use stats, only : overall_statistic
use var, only : nzmsize
use var, only : itime
use var, only : numscalar, nrhotime, npress
real(mytype),dimension(xsize(1),xsize(2),xsize(3)), intent(in) :: ux1, uy1, uz1
real(mytype),dimension(xsize(1),xsize(2),xsize(3),numscalar), intent(in) :: phi1
real(mytype),dimension(xsize(1),xsize(2),xsize(3),nrhotime), intent(in) :: rho1
real(mytype),dimension(xsize(1),xsize(2),xsize(3)), intent(in) :: ep1
real(mytype), dimension(ph1%zst(1):ph1%zen(1), ph1%zst(2):ph1%zen(2), nzmsize, npress), intent(in) :: pp3
!call write_snapshot(rho1, ux1, uy1, uz1, pp3, phi1, ep1, itime)
!call postprocess_case(rho1, ux1, uy1, uz1, pp3, phi1, ep1)
!call overall_statistic(ux1, uy1, uz1, phi1, pp3, ep1)
end subroutine preprocessing
!##################################################################
!##################################################################
subroutine postprocessing(rho1, ux1, uy1, uz1, pp3, phi1, ep1)
use decomp_2d, only : xsize, ph1
use var, only : nzmsize, numscalar, nrhotime, npress, abl_T
real(mytype),dimension(xsize(1),xsize(2),xsize(3)), intent(in) :: ux1, uy1, uz1
real(mytype),dimension(xsize(1),xsize(2),xsize(3),numscalar), intent(in) :: phi1
real(mytype),dimension(xsize(1),xsize(2),xsize(3),nrhotime), intent(in) :: rho1
real(mytype),dimension(xsize(1),xsize(2),xsize(3)), intent(in) :: ep1
real(mytype),dimension(ph1%zst(1):ph1%zen(1), ph1%zst(2):ph1%zen(2), nzmsize, npress), intent(in) :: pp3
integer :: j
! Recover temperature when decomposed (pressure to be recovered externally)
if (itype.eq.itype_abl.and.ibuoyancy.eq.1) then
do j=1,xsize(2)
abl_T(:,j,:,1) = phi1(:,j,:,1) + Tstat(j,1)
enddo
call run_postprocessing(rho1, ux1, uy1, uz1, pp3, abl_T, ep1)
else
call run_postprocessing(rho1, ux1, uy1, uz1, pp3, phi1, ep1)
endif
end subroutine postprocessing
!##################################################################
!##################################################################
subroutine run_postprocessing(rho1, ux1, uy1, uz1, pp3, phi1, ep1)
use decomp_2d, only : xsize, ph1
use visu, only : write_snapshot, end_snapshot
use stats, only : overall_statistic
use var, only : nzmsize
use var, only : itime
use var, only : numscalar, nrhotime, npress
use turbine, only : turbine_output
use probes, only : write_probes
real(mytype),dimension(xsize(1),xsize(2),xsize(3)), intent(in) :: ux1, uy1, uz1
real(mytype),dimension(xsize(1),xsize(2),xsize(3),numscalar), intent(in) :: phi1
real(mytype),dimension(xsize(1),xsize(2),xsize(3),nrhotime), intent(in) :: rho1
real(mytype),dimension(xsize(1),xsize(2),xsize(3)), intent(in) :: ep1
real(mytype),dimension(ph1%zst(1):ph1%zen(1), ph1%zst(2):ph1%zen(2), nzmsize, npress), intent(in) :: pp3
integer :: num
if ((ivisu.ne.0).and.(mod(itime, ioutput).eq.0)) then
call write_snapshot(rho1, ux1, uy1, uz1, pp3, phi1, ep1, itime, num)
! XXX: Ultimate goal for ADIOS2 is to pass do all postproc online - do we need this?
! Currently, needs some way to "register" variables for IO
call visu_case(rho1, ux1, uy1, uz1, pp3, phi1, ep1, num)
call end_snapshot(itime, num)
end if
call postprocess_case(rho1, ux1, uy1, uz1, pp3, phi1, ep1)
call overall_statistic(ux1, uy1, uz1, phi1, pp3, ep1)
if (iturbine.ne.0) then
call turbine_output()
endif
call write_probes(ux1, uy1, uz1, pp3, phi1)
end subroutine run_postprocessing
!##################################################################
!##################################################################
subroutine postprocess_case(rho,ux,uy,uz,pp,phi,ep)
use forces
use var, only : nzmsize
use param, only : npress
real(mytype),dimension(xsize(1),xsize(2),xsize(3)) :: ux,uy,uz
real(mytype),dimension(xsize(1),xsize(2),xsize(3),numscalar) :: phi
real(mytype),dimension(xsize(1),xsize(2),xsize(3),nrhotime) :: rho
real(mytype),dimension(xsize(1),xsize(2),xsize(3)) :: ep
real(mytype), dimension(ph1%zst(1):ph1%zen(1), ph1%zst(2):ph1%zen(2), nzmsize, npress), intent(in) :: pp
if (itype.eq.itype_generic) then
call postprocess_generic (ux, uy, uz, phi, ep)
elseif (itype.eq.itype_gravitycur) then
call postprocess_gravitycur(rho, ux, uy, uz, phi, ep)
elseif (itype.eq.itype_tgv) then
call postprocess_tgv (ux, uy, uz, phi, ep)
elseif (itype.eq.itype_channel) then
call postprocess_channel (ux, uy, uz, pp, phi, ep)
elseif (itype.eq.itype_hill) then
call postprocess_hill(ux, uy, uz, pp, phi, ep)
elseif (itype.eq.itype_cyl) then
call postprocess_cyl (ux, uy, uz, ep)
elseif (itype.eq.itype_tbl) then
call postprocess_tbl (ux, uy, uz, ep)
elseif (itype.eq.itype_abl) then
call postprocess_abl (ux, uy, uz, ep)
elseif (itype.eq.itype_uniform) then
call postprocess_uniform (ux, uy, uz, ep)
elseif (itype.EQ.itype_sandbox) THEN
call postprocess_sandbox (ux, uy, uz, phi, ep)
elseif (itype.eq.itype_cavity) then
call postprocess_cavity(ux, uy, uz, phi)
elseif (itype.eq.itype_pipe) then
call postprocess_pipe(ux, uy, uz, pp, phi, ep)
endif
if (iforces.eq.1) then
call force(ux,uy,uz,ep)
call restart_forces(1)
endif
end subroutine postprocess_case
!##################################################################
!!
!! SUBROUTINE: visu_case_init
!! AUTHOR: PB
!! DESCRIPTION: Initialise case-specific visualization
!!
!##################################################################
subroutine visu_case_init
implicit none
if (itype .eq. itype_tgv) then
call visu_tgv_init(case_visu_init)
else if (itype .eq. itype_channel) then
call visu_channel_init(case_visu_init)
else if (itype .eq. itype_hill) then
call visu_hill_init(case_visu_init)
else if (itype .eq. itype_cyl) then
call visu_cyl_init(case_visu_init)
else if (itype .eq. itype_tbl) then
call visu_tbl_init(case_visu_init)
else if (itype .eq. itype_gravitycur) then
call visu_gravitycur_init(case_visu_init)
else if (itype .eq. itype_uniform) then
call visu_uniform_init(case_visu_init)
end if
end subroutine visu_case_init
!##################################################################
!!
!! SUBROUTINE: visu_case
!! AUTHOR: CF
!! DESCRIPTION: Call case-specific visualization
!!
!##################################################################
subroutine visu_case(rho1,ux1,uy1,uz1,pp3,phi1,ep1,num)
use var, only : nzmsize
use param, only : npress
real(mytype), intent(in), dimension(xsize(1),xsize(2),xsize(3),nrhotime) :: rho1
real(mytype), intent(in), dimension(xsize(1),xsize(2),xsize(3)) :: ux1,uy1,uz1
real(mytype), intent(in), dimension(ph1%zst(1):ph1%zen(1), ph1%zst(2):ph1%zen(2), nzmsize, npress) :: pp3
real(mytype), intent(in), dimension(xsize(1),xsize(2),xsize(3),numscalar) :: phi1
real(mytype), intent(in), dimension(xsize(1),xsize(2),xsize(3)) :: ep1
integer, intent(in) :: num
logical :: called_visu = .false.
if (itype.eq.itype_generic) then
call visu_generic(ux1, uy1, uz1, pp3, phi1, ep1, num)
called_visu = .true.
elseif (itype.eq.itype_tgv) then
call visu_tgv(ux1, uy1, uz1, pp3, phi1, ep1, num)
called_visu = .true.
elseif (itype.eq.itype_channel) then
call visu_channel(ux1, uy1, uz1, pp3, phi1, ep1, num)
called_visu = .true.
elseif (itype.eq.itype_hill) then
call visu_hill(ux1, uy1, uz1, pp3, phi1, ep1, num)
called_visu = .true.
elseif (itype.eq.itype_cyl) then
call visu_cyl(ux1, uy1, uz1, pp3, phi1, ep1, num)
called_visu = .true.
elseif (itype.eq.itype_tbl) then
call visu_tbl(ux1, uy1, uz1, pp3, phi1, ep1, num)
called_visu = .true.
elseif (itype.eq.itype_uniform) then
call visu_uniform(ux1, uy1, uz1, pp3, phi1, ep1, num)
called_visu = .true.
endif
if (called_visu .and. (.not. case_visu_init)) then
print *, "ERROR: tried to run case-specific visu without initialisation!"
print *, " See the TGV case initialisation for example."
STOP
endif
end subroutine visu_case
!##################################################################
!##################################################################
!!
!! SUBROUTINE: momentum_forcing
!! AUTHOR: Paul Bartholomew
!! DESCRIPTION: Calls case-specific forcing functions for the
!! momentum equations.
!!
!##################################################################
subroutine momentum_forcing(dux1, duy1, duz1, rho1, ux1, uy1, uz1, phi1)
use mhd, only: mhd_active,momentum_forcing_mhd
implicit none
real(mytype), intent(in), dimension(xsize(1), xsize(2), xsize(3)) :: ux1, uy1, uz1
real(mytype), intent(in), dimension(xsize(1), xsize(2), xsize(3), nrhotime) :: rho1
real(mytype), intent(in), dimension(xsize(1), xsize(2), xsize(3), numscalar) :: phi1
real(mytype), dimension(xsize(1), xsize(2), xsize(3), ntime) :: dux1, duy1, duz1
if (itype.eq.itype_channel) then
call momentum_forcing_channel(dux1, duy1, duz1, ux1, uy1, uz1)
elseif (itype.eq.itype_abl) then
call momentum_forcing_abl(dux1, duy1, duz1, ux1, uy1, uz1, phi1)
endif
if(mhd_active) then
call momentum_forcing_mhd(dux1(:,:,:,1),duy1(:,:,:,1),duz1(:,:,:,1),ux1,uy1,uz1)
endif
end subroutine momentum_forcing
!##################################################################
!##################################################################
!!
!! SUBROUTINE: scalar_forcing
!! AUTHOR: Kay Schäfer
!! DESCRIPTION: Calls case-specific forcing functions for the
!! scalar transport equations.
!!
!##################################################################
subroutine scalar_forcing(dphi1, rho1, ux1, uy1, uz1, phi1)
implicit none
real(mytype), intent(in), dimension(xsize(1), xsize(2), xsize(3)) :: ux1, uy1, uz1, phi1
real(mytype), intent(in), dimension(xsize(1), xsize(2), xsize(3), nrhotime) :: rho1
real(mytype), dimension(xsize(1),xsize(2),xsize(3),ntime) :: dphi1
if (itype.eq.itype_abl) then
call scalar_forcing_abl(uy1, dphi1, phi1)
endif
end subroutine scalar_forcing
!##################################################################
!##################################################################
subroutine set_fluid_properties(rho1, mu1)
implicit none
real(mytype), dimension(xsize(1), xsize(2), xsize(3)), intent(in) :: rho1
real(mytype), dimension(xsize(1), xsize(2), xsize(3)) :: mu1
if (itype.eq.itype_gravitycur) then
if (ilmn) then
call set_fluid_properties_gravitycur(rho1, mu1)
end if
endif
endsubroutine set_fluid_properties
!##################################################################
!##################################################################
subroutine test_flow(rho1,ux1,uy1,uz1,phi1,ep1,drho1,divu3)
use decomp_2d
use param
use navier, only : divergence
use var, only : numscalar, dv3
use tools, only : test_speed_min_max, compute_cfl, test_scalar_min_max
implicit none
real(mytype), dimension(xsize(1), xsize(2), xsize(3)), intent(in) :: ux1, uy1, uz1, ep1
real(mytype), dimension(xsize(1), xsize(2), xsize(3), nrhotime), intent(in) :: rho1
real(mytype), dimension(xsize(1), xsize(2), xsize(3), numscalar), intent(in) :: phi1
real(mytype), dimension(xsize(1), xsize(2), xsize(3), ntime), intent(in) :: drho1
real(mytype), dimension(zsize(1), zsize(2), zsize(3)), intent(in) :: divu3
if ((mod(itime,ilist)==0 .or. itime == ifirst .or. itime == ilast)) then
call divergence(dv3,rho1,ux1,uy1,uz1,ep1,drho1,divu3,2)
call test_speed_min_max(ux1,uy1,uz1)
call compute_cfl(ux1,uy1,uz1)
if (iscalar==1) call test_scalar_min_max(phi1)
endif
end subroutine test_flow
!##################################################################
!##################################################################
end module case
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!! case.f90 ends here