-
Notifications
You must be signed in to change notification settings - Fork 2
/
Copy pathDERIVE_3D.F
1784 lines (1528 loc) · 51.8 KB
/
DERIVE_3D.F
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
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
#include "AMReX_REAL.H"
#include "AMReX_CONSTANTS.H"
#include "AMReX_ArrayLim.H"
#include "Derived.H"
#ifdef BL_DERIVE_IAMR
#define TRACER_POSITION 5
#else
#define TRACER_POSITION 6
#endif
subroutine FORT_SETGAMMA (gam)
REAL_T gam
#include "xxmeth.fh"
gamma = gam
end
subroutine FORT_SETALPHA (alp)
REAL_T alp
#include "xxmeth.fh"
alpha = alp
end
subroutine FORT_SETQPOS (qpos)
integer qpos
#include "xxmeth.fh"
qposition = qpos
end
c
c derive log density:
c
c Inputs/Outputs
c u => (const) state array
c ulo,uhi => (const) index limits u
c nvar => (const) number of variables in state array
c xlo => (const) physical location of lo end of u
c xhi => (const) physical location of hi end of u
c dat <=> (modify) array holding derived data
c dlo,dhi => (const) index limits of dat array
c ovlo,ovhi => (const) subregion where derivation is done
c
subroutine FORT_DERLGDN (u,uloi1,uloi2,uloi3,uhii1,uhii2,uhii3,
$ nvar,xlo, xhi,dat,dloi1,dloi2,dloi3,dhii1,dhii2,dhii3, ovlo,ovhi)
integer uloi1,uloi2,uloi3
integer uhii1,uhii2,uhii3
integer dloi1,dloi2,dloi3
integer dhii1,dhii2,dhii3
integer ovlo(3), ovhi(3)
integer nvar
REAL_T xlo(3), xhi(3)
REAL_T u(uloi1:uhii1, uloi2:uhii2, uloi3:uhii3, nvar)
REAL_T dat(dloi1:dhii1, dloi2:dhii2, dloi3:dhii3)
integer i,j,k
REAL_T smallr
data smallr /1.0e-6/
do k = ovlo(3),ovhi(3)
do j = ovlo(2),ovhi(2)
do i = ovlo(1),ovhi(1)
dat(i,j,k) = log10(max(u(i,j,k,1),smallr))
end do
end do
end do
end
c
c derive density:
c
c Inputs/Outputs
c u => (const) state array
c ulo,uhi => (const) index limits of u
c nvar => (const) number of variables in state array
c xlo => (const) physical location of lo end of u
c xhi => (const) physical location of hi end of u
c dat <=> (modify) array holding derived data
c dlo,dhi => (const) index limits of dat array
c ovlo,ovhi => (const) subregion where derivation is done
c
subroutine FORT_DERDEN (u,uloi1,uloi2,uloi3,uhii1,uhii2,uhii3,
$ nvar,xlo, xhi,dat,dloi1,dloi2,dloi3,dhii1,dhii2,dhii3, ovlo,ovhi)
integer uloi1,uloi2,uloi3
integer uhii1,uhii2,uhii3
integer dloi1,dloi2,dloi3
integer dhii1,dhii2,dhii3
integer ovlo(3), ovhi(3)
integer nvar
REAL_T xlo(3), xhi(3)
REAL_T u(uloi1:uhii1, uloi2:uhii2, uloi3:uhii3, nvar)
REAL_T dat(dloi1:dhii1, dloi2:dhii2, dloi3:dhii3)
integer i,j,k
do k = ovlo(3),ovhi(3)
do j = ovlo(2),ovhi(2)
do i = ovlo(1),ovhi(1)
dat(i,j,k) = u(i,j,k,1)
end do
end do
end do
end
c
c derive x velocity:
c
c Inputs/Outputs
c u => (const) state array
c ulo,uhi => (const) index limits of u
c nvar => (const) number of variables in state array
c xlo => (const) physical location of lo end of u
c xhi => (const) physical location of hi end of u
c dat <=> (modify) array holding derived data
c dlo,dhi => (const) index limits of dat array
c ovlo,ovhi => (const) subregion where derivation is done
c
subroutine FORT_DERXVEL (u,uloi1,uloi2,uloi3,uhii1,uhii2,uhii3,
$ nvar,xlo, xhi,dat,dloi1,dloi2,dloi3,dhii1,dhii2,dhii3, ovlo,ovhi)
integer uloi1,uloi2,uloi3
integer uhii1,uhii2,uhii3
integer dloi1,dloi2,dloi3
integer dhii1,dhii2,dhii3
integer ovlo(3), ovhi(3)
integer nvar
REAL_T xlo(3), xhi(3)
REAL_T u(uloi1:uhii1, uloi2:uhii2, uloi3:uhii3, nvar)
REAL_T dat(dloi1:dhii1, dloi2:dhii2, dloi3 :dhii3)
integer i,j,k
do k = ovlo(3),ovhi(3)
do j = ovlo(2),ovhi(2)
do i = ovlo(1),ovhi(1)
dat(i,j,k) = u(i,j,k,2)/u(i,j,k,1)
end do
end do
end do
end
c
c derive y velocity:
c
c Inputs/Outputs
c u => (const) state array
c ulo,uhi => (const) index limits of u
c nvar => (const) number of variables in state array
c xlo => (const) physical location of lo end of u
c xhi => (const) physical location of hi end of u
c dat <=> (modify) array holding derived data
c dlo,dhi => (const) index limits of dat array
c ovlo,ovhi => (const) subregion where derivation is done
c
subroutine FORT_DERYVEL (u,uloi1,uloi2,uloi3,uhii1,uhii2,uhii3,
$ nvar,xlo, xhi,dat,dloi1,dloi2,dloi3,dhii1,dhii2,dhii3, ovlo,ovhi)
integer uloi1,uloi2,uloi3
integer uhii1,uhii2,uhii3
integer dloi1,dloi2,dloi3
integer dhii1,dhii2,dhii3
integer ovlo(3), ovhi(3)
integer nvar
REAL_T xlo(3), xhi(3)
REAL_T u(uloi1:uhii1, uloi2:uhii2, uloi3:uhii3, nvar)
REAL_T dat(dloi1:dhii1, dloi2:dhii2, dloi3:dhii3)
integer i,j,k
do k = ovlo(3),ovhi(3)
do j = ovlo(2),ovhi(2)
do i = ovlo(1),ovhi(1)
dat(i,j,k) = u(i,j,k,3)/u(i,j,k,1)
end do
end do
end do
end
c
c derive z velocity:
c
c Inputs/Outputs
c u => (const) state array
c ulo,uhi => (const) index limits of u
c nvar => (const) number of variables in state array
c xlo => (const) physical location of lo end of u
c xhi => (const) physical location of hi end of u
c dat <=> (modify) array holding derived data
c dlo,dhi => (const) index limits of dat array
c ovlo,ovhi => (const) subregion where derivation is done
c
subroutine FORT_DERZVEL (u,uloi1,uloi2,uloi3,uhii1,uhii2,uhii3,
$ nvar,xlo, xhi,dat,dloi1,dloi2,dloi3,dhii1,dhii2,dhii3, ovlo,ovhi)
integer uloi1,uloi2,uloi3
integer uhii1,uhii2,uhii3
integer dloi1,dloi2,dloi3
integer dhii1,dhii2,dhii3
integer ovlo(3), ovhi(3)
integer nvar
REAL_T xlo(3), xhi(3)
REAL_T u(uloi1:uhii1, uloi2:uhii2, uloi3:uhii3, nvar)
REAL_T dat(dloi1:dhii1, dloi2:dhii2, dloi3:dhii3)
integer i,j,k
do k = ovlo(3),ovhi(3)
do j = ovlo(2),ovhi(2)
do i = ovlo(1),ovhi(1)
dat(i,j,k) = u(i,j,k,4)/u(i,j,k,1)
end do
end do
end do
end
c
c derive total energy = (rho*E)/rho:
c
c Inputs/Outputs
c u => (const) state array
c ulo,uhi => (const) index limits of u
c nvar => (const) number of variables in state array
c xlo => (const) physical location of lo end of u
c xhi => (const) physical location of hi end of u
c dat <=> (modify) array holding derived data
c dlo,dhi => (const) index limits of dat array
c ovlo,ovhi => (const) subregion where derivation is done
c
subroutine FORT_DERTENG (u,uloi1,uloi2,uloi3,uhii1,uhii2,uhii3,
$ nvar,xlo,xhi,dat,dloi1,dloi2,dloi3,dhii1,dhii2,dhii3,ovlo,ovhi)
integer uloi1,uloi2,uloi3
integer uhii1,uhii2,uhii3
integer dloi1,dloi2,dloi3
integer dhii1,dhii2,dhii3
integer ovlo(3), ovhi(3)
integer nvar
REAL_T xlo(3), xhi(3)
REAL_T u(uloi1:uhii1, uloi2:uhii2, uloi3:uhii3, nvar)
REAL_T dat(dloi1:dhii1, dloi2:dhii2, dloi3:dhii3)
integer i,j,k
do k = ovlo(3),ovhi(3)
do j = ovlo(2),ovhi(2)
do i = ovlo(1),ovhi(1)
dat(i,j,k) = u(i,j,k,5)/u(i,j,k,1)
end do
end do
end do
end
c
c derive total energy = (rho*(E-(u^2+v^2)/2)/rho:
c
c Inputs/Outputs
c u => (const) state array
c ulo,uhi => (const) index limits of u
c nvar => (const) number of variables in state array
c xlo => (const) physical location of lo end of u
c xhi => (const) physical location of hi end of u
c dat <=> (modify) array holding derived data
c dlo,dhi => (const) index limits of dat array
c ovlo,ovhi => (const) subregion where derivation is done
c
subroutine FORT_DERIENG (u,uloi1,uloi2,uloi3,uhii1,uhii2,uhii3,
$ nvar,xlo, xhi,dat,dloi1,dloi2,dloi3,dhii1,dhii2,dhii3, ovlo,ovhi)
integer uloi1,uloi2,uloi3
integer uhii1,uhii2,uhii3
integer dloi1,dloi2,dloi3
integer dhii1,dhii2,dhii3
integer ovlo(3), ovhi(3)
integer nvar
REAL_T xlo(3), xhi(3)
REAL_T u(uloi1:uhii1, uloi2:uhii2, uloi3:uhii3, nvar)
REAL_T dat(dloi1:dhii1, dloi2:dhii2, dloi3:dhii3)
integer i,j,k
REAL_T vx,vy,vz,eng
do k = ovlo(3),ovhi(3)
do j = ovlo(2),ovhi(2)
do i = ovlo(1),ovhi(1)
vx = u(i,j,k,2)/u(i,j,k,1)
vy = u(i,j,k,3)/u(i,j,k,1)
vz = u(i,j,k,4)/u(i,j,k,1)
eng = u(i,j,k,5)/u(i,j,k,1)
dat(i,j,k) = eng - (vx**2+vy**2+vz**2)*half
end do
end do
end do
end
c
c derive kenetic energy = rho(u**2 + v**2 + w**2)/2
c
c Inputs/Outputs
c u => (const) state array
c ulo,uhi => (const) index limits of u
c nvar => (const) number of variables in state array
c xlo => (const) physical location of lo end of u
c xhi => (const) physical location of hi end of u
c dat <=> (modify) array holding derived data
c dlo,dhi => (const) index limits of dat array
c ovlo,ovhi => (const) subregion where derivation is done
c
subroutine FORT_DERKENG (u,uloi1,uloi2,uloi3,uhii1,uhii2,uhii3,
$ nvar,xlo,xhi,dat,dloi1,dloi2,dloi3,dhii1,dhii2,dhii3,ovlo,ovhi)
integer uloi1,uloi2,uloi3
integer uhii1,uhii2,uhii3
integer dloi1,dloi2,dloi3
integer dhii1,dhii2,dhii3
integer ovlo(3), ovhi(3)
integer nvar
REAL_T xlo(3), xhi(3)
REAL_T u(uloi1:uhii1, uloi2:uhii2, uloi3:uhii3, nvar)
REAL_T dat(dloi1:dhii1, dloi2:dhii2, dloi3:dhii3)
integer i,j,k
REAL_T eng
do k = ovlo(3),ovhi(3)
do j = ovlo(2),ovhi(2)
do i = ovlo(1),ovhi(1)
eng = u(i,j,k,2)**2 + u(i,j,k,3)**2 + u(i,j,k,4)**2
dat(i,j,k) = half*eng/u(i,j,k,1)
end do
end do
end do
end
c
c derive vorticity
c
c Inputs/Outputs
c u => (const) state array
c ulo,uhi => (const) index limits of u
c nvar => (const) number of variables in state array
c xlo => (const) physical location of lo end of u
c xhi => (const) physical location of hi end of u
c dat <=> (modify) array holding derived data
c dlo,dhi => (const) index limits of dat array
c ovlo,ovhi => (const) subregion where derivation is done
c
c Uses one-sided derivatives at edges.
c
subroutine FORT_DERVORT (f,uloi1,uloi2,uloi3,uhii1,uhii2,uhii3,
$ nvar,xlo,xhi,dat,dloi1,dloi2,dloi3,dhii1,dhii2,dhii3,ovlo,ovhi)
integer uloi1,uloi2,uloi3
integer uhii1,uhii2,uhii3
integer dloi1,dloi2,dloi3
integer dhii1,dhii2,dhii3
integer ovlo(3), ovhi(3)
integer nvar
REAL_T xlo(3), xhi(3)
REAL_T f(uloi1:uhii1, uloi2:uhii2, uloi3:uhii3, nvar)
REAL_T dat(dloi1:dhii1, dloi2:dhii2, dloi3:dhii3)
REAL_T uym,uyp,uy,uzm,uzp,uz
REAL_T vxm,vxp,vx,vzm,vzp,vz
REAL_T wxm,wxp,wx,wym,wyp,wy
REAL_T uc,vc,wc
REAL_T uy2m,uy2p,uz2m,uz2p
REAL_T vx2m,vx2p,vz2m,vz2p
REAL_T wx2m,wx2p,wy2m,wy2p
REAL_T vorx,vory,vorz
REAL_T ddx,ddy,ddz
integer i,j,k
ddx = half*float(uhii1-uloi1+1)/(xhi(1)-xlo(1))
ddy = half*float(uhii2-uloi2+1)/(xhi(2)-xlo(2))
ddz = half*float(uhii3-uloi3+1)/(xhi(3)-xlo(3))
do k = ovlo(3),ovhi(3)
do j = ovlo(2),ovhi(2)
do i = ovlo(1),ovhi(1)
uc = f(i,j,k,2)/f(i,j,k,1)
if (j.eq.ovlo(2).and.j+2.le.ovhi(2)) then
uy2p = f(i,j+2,k,2)/f(i,j+2,k,1)
uyp = f(i,j+1,k,2)/f(i,j+1,k,1)
uy = ddy*(-uy2p+4*uyp-3*uc)
else if (j.eq.ovlo(2)) then
uyp = f(i,j+1,k,2)/f(i,j+1,k,1)
uy = 2*ddy*(uyp-uc)
else if (j.eq.ovhi(2).and.j-2.ge.ovlo(2)) then
uy2m = f(i,j-2,k,2)/f(i,j-2,k,1)
uym = f(i,j-1,k,2)/f(i,j-1,k,1)
uy = ddy*(+uy2m-4*uym+3*uc)
else if (j.eq.ovhi(2)) then
uym = f(i,j-1,k,2)/f(i,j-1,k,1)
uy = 2*ddy*(-uym+uc)
else if (j.gt.ovlo(2).and.j.lt.ovhi(2) ) then
uyp = f(i,j+1,k,2)/f(i,j+1,k,1)
uym = f(i,j-1,k,2)/f(i,j-1,k,1)
uy = ddy*(uyp-uym)
else
write(6,*) 'uy'
stop
end if
if (k.eq.ovlo(3).and.k+2.le.ovhi(3)) then
uz2p = f(i,j,k+2,2)/f(i,j,k+2,1)
uzp = f(i,j,k+1,2)/f(i,j,k+1,1)
uz = ddz*(-uz2p+4*uzp-3*uc)
else if (k.eq.ovlo(3)) then
uzp = f(i,j,k+1,2)/f(i,j,k+1,1)
uz = 2*ddz*(uzp-uc)
else if (k.eq.ovhi(3).and.k-2.ge.ovlo(3)) then
uz2m = f(i,j,k-2,2)/f(i,j,k-2,1)
uzm = f(i,j,k-1,2)/f(i,j,k-1,1)
uz = ddz*(+uz2m-4*uzm+3*uc)
else if (k.eq.ovhi(3)) then
uzm = f(i,j,k-1,2)/f(i,j,k-1,1)
uz = 2*ddz*(-uzm+uc)
else if (k.gt.ovlo(3) .and. k.lt.ovhi(3)) then
uzp = f(i,j,k+1,2)/f(i,j,k+1,1)
uzm = f(i,j,k-1,2)/f(i,j,k-1,1)
uz = ddz*(uzp-uzm)
else
write(6,*) 'uz'
stop
end if
vc = f(i,j,k,3)/f(i,j,k,1)
if (i.eq.ovlo(1).and.i+2.le.ovhi(1)) then
vx2p = f(i+2,j,k,3)/f(i+2,j,k,1)
vxp = f(i+1,j,k,3)/f(i+1,j,k,1)
vx = ddx*(-vx2p+4*vxp-3*vc)
else if (i.eq.ovlo(1)) then
vxp = f(i+1,j,k,3)/f(i+1,j,k,1)
vx = 2*ddx*(vxp-vc)
else if (i.eq.ovhi(1).and.i-2.ge.ovlo(1)) then
vx2m = f(i-2,j,k,3)/f(i-2,j,k,1)
vxm = f(i-1,j,k,3)/f(i-1,j,k,1)
vx = ddx*(+vx2m-4*vxm+3*vc)
else if (i.eq.ovhi(1)) then
vxm = f(i-1,j,k,3)/f(i-1,j,k,1)
vx = 2*ddx*(-vxm+vc)
else if (i.gt.ovlo(1).and.i.lt.ovhi(1)) then
vxp = f(i+1,j,k,3)/f(i+1,j,k,1)
vxm = f(i-1,j,k,3)/f(i-1,j,k,1)
vx = ddx*(vxp-vxm)
else
write(6,*) 'vx'
stop
end if
if (k.eq.ovlo(3).and.k+2.le.ovhi(3)) then
vz2p = f(i,j,k+2,3)/f(i,j,k+2,1)
vzp = f(i,j,k+1,3)/f(i,j,k+1,1)
vz = ddz*(-vz2p+4*vzp-3*vc)
else if (k.eq.ovlo(3)) then
vzp = f(i,j,k+1,3)/f(i,j,k+1,1)
vz = 2*ddz*(vzp-vc)
else if (k.eq.ovhi(3).and.k-2.ge.ovlo(3)) then
vz2m = f(i,j,k-2,3)/f(i,j,k-2,1)
vzm = f(i,j,k-1,3)/f(i,j,k-1,1)
vz = ddz*(+vz2m-4*vzm+3*vc)
else if (k.eq.ovhi(3)) then
vzm = f(i,j,k-1,3)/f(i,j,k-1,1)
vz = 2*ddz*(-vzm+vc)
else if (k.gt.ovlo(3).and.k.lt.ovhi(3)) then
vzp = f(i,j,k+1,3)/f(i,j,k+1,1)
vzm = f(i,j,k-1,3)/f(i,j,k-1,1)
vz = ddz*(vzp-vzm)
else
write(6,*) 'vz'
stop
end if
wc = f(i,j,k,4)/f(i,j,k,1)
if (i.eq.ovlo(1).and.i+2.le.ovhi(1)) then
wx2p = f(i+2,j,k,4)/f(i+2,j,k,1)
wxp = f(i+1,j,k,4)/f(i+1,j,k,1)
wx = ddx*(-wx2p+4*wxp-3*wc)
else if (i.eq.ovlo(1)) then
wxp = f(i+1,j,k,4)/f(i+1,j,k,1)
wx = 2*ddx*(wxp-wc)
else if (i.eq.ovhi(1).and.i-2.ge.ovlo(1)) then
wx2m = f(i-2,j,k,4)/f(i-2,j,k,1)
wxm = f(i-1,j,k,4)/f(i-1,j,k,1)
wx = ddx*(+wx2m-4*wxm+3*wc)
else if (i.eq.ovhi(1)) then
wxm = f(i-1,j,k,4)/f(i-1,j,k,1)
wx = 2*ddx*(-wxm+wc)
else if (i.gt.ovlo(1).and.i.lt.ovhi(1)) then
wxp = f(i+1,j,k,4)/f(i+1,j,k,1)
wxm = f(i-1,j,k,4)/f(i-1,j,k,1)
wx = ddx*(wxp-wxm)
else
write(6,*) 'wx'
stop
end if
if (j.eq.ovlo(2).and.j+2.le.ovhi(2)) then
wy2p = f(i,j+2,k,4)/f(i,j+2,k,1)
wyp = f(i,j+1,k,4)/f(i,j+1,k,1)
wy = ddy*(-wy2p+4*wyp-3*wc)
else if (j.eq.ovlo(2)) then
wyp = f(i,j+1,k,4)/f(i,j+1,k,1)
wy = 2*ddy*(wyp-wc)
else if (j.eq.ovhi(2).and.j-2.ge.ovlo(2)) then
wy2m = f(i,j-2,k,4)/f(i,j-2,k,1)
wym = f(i,j-1,k,4)/f(i,j-1,k,1)
wy = ddy*(+wy2m-4*wym+3*wc)
else if (j.eq.ovhi(2)) then
wym = f(i,j-1,k,4)/f(i,j-1,k,1)
wy = 2*ddy*(-wym+wc)
else if (j.gt.ovlo(2).and.j.lt.ovhi(2)) then
wyp = f(i,j+1,k,4)/f(i,j+1,k,1)
wym = f(i,j-1,k,4)/f(i,j-1,k,1)
wy = ddy*(wyp-wym)
else
write(6,*) 'wy'
stop
end if
vorx = wy-vz
vory = uz-wx
vorz = vx-uy
dat(i,j,k) = sqrt(vorx**2+vory**2+vorz**2)
end do
end do
end do
end
c
c derive vorticity components
c
c Inputs/Outputs
c u => (const) state array
c ulo,uhi => (const) index limits of u
c nvar => (const) number of variables in state array
c xlo => (const) physical location of lo end of u
c xhi => (const) physical location of hi end of u
c dat <=> (modify) array holding derived data
c dlo,dhi => (const) index limits of dat array
c ovlo,ovhi => (const) subregion where derivation is done
c
c Uses one-sided derivatives at edges.
c
subroutine FORT_DERVORTCOMPS (f,uloi1,uloi2,uloi3,uhii1,uhii2,uhii3,
$ nvar,xlo,xhi,dat,dloi1,dloi2,dloi3,dhii1,dhii2,dhii3,ovlo,ovhi)
integer uloi1,uloi2,uloi3
integer uhii1,uhii2,uhii3
integer dloi1,dloi2,dloi3
integer dhii1,dhii2,dhii3
integer ovlo(3), ovhi(3)
integer nvar
REAL_T xlo(3), xhi(3)
REAL_T f(uloi1:uhii1, uloi2:uhii2, uloi3:uhii3, nvar)
REAL_T dat(dloi1:dhii1, dloi2:dhii2, dloi3:dhii3, 3)
REAL_T uym,uyp,uy,uzm,uzp,uz
REAL_T vxm,vxp,vx,vzm,vzp,vz
REAL_T wxm,wxp,wx,wym,wyp,wy
REAL_T uc,vc,wc
REAL_T uy2m,uy2p,uz2m,uz2p
REAL_T vx2m,vx2p,vz2m,vz2p
REAL_T wx2m,wx2p,wy2m,wy2p
REAL_T vorx,vory,vorz
REAL_T ddx,ddy,ddz
integer i,j,k
ddx = half*float(uhii1-uloi1+1)/(xhi(1)-xlo(1))
ddy = half*float(uhii2-uloi2+1)/(xhi(2)-xlo(2))
ddz = half*float(uhii3-uloi3+1)/(xhi(3)-xlo(3))
do k = ovlo(3),ovhi(3)
do j = ovlo(2),ovhi(2)
do i = ovlo(1),ovhi(1)
uc = f(i,j,k,2)/f(i,j,k,1)
if (j.eq.ovlo(2).and.j+2.le.ovhi(2)) then
uy2p = f(i,j+2,k,2)/f(i,j+2,k,1)
uyp = f(i,j+1,k,2)/f(i,j+1,k,1)
uy = ddy*(-uy2p+4*uyp-3*uc)
else if (j.eq.ovlo(2)) then
uyp = f(i,j+1,k,2)/f(i,j+1,k,1)
uy = 2*ddy*(uyp-uc)
else if (j.eq.ovhi(2).and.j-2.ge.ovlo(2)) then
uy2m = f(i,j-2,k,2)/f(i,j-2,k,1)
uym = f(i,j-1,k,2)/f(i,j-1,k,1)
uy = ddy*(+uy2m-4*uym+3*uc)
else if (j.eq.ovhi(2)) then
uym = f(i,j-1,k,2)/f(i,j-1,k,1)
uy = 2*ddy*(-uym+uc)
else if (j.gt.ovlo(2).and.j.lt.ovhi(2) ) then
uyp = f(i,j+1,k,2)/f(i,j+1,k,1)
uym = f(i,j-1,k,2)/f(i,j-1,k,1)
uy = ddy*(uyp-uym)
else
write(6,*) 'uy'
stop
end if
if (k.eq.ovlo(3).and.k+2.le.ovhi(3)) then
uz2p = f(i,j,k+2,2)/f(i,j,k+2,1)
uzp = f(i,j,k+1,2)/f(i,j,k+1,1)
uz = ddz*(-uz2p+4*uzp-3*uc)
else if (k.eq.ovlo(3)) then
uzp = f(i,j,k+1,2)/f(i,j,k+1,1)
uz = 2*ddz*(uzp-uc)
else if (k.eq.ovhi(3).and.k-2.ge.ovlo(3)) then
uz2m = f(i,j,k-2,2)/f(i,j,k-2,1)
uzm = f(i,j,k-1,2)/f(i,j,k-1,1)
uz = ddz*(+uz2m-4*uzm+3*uc)
else if (k.eq.ovhi(3)) then
uzm = f(i,j,k-1,2)/f(i,j,k-1,1)
uz = 2*ddz*(-uzm+uc)
else if (k.gt.ovlo(3) .and. k.lt.ovhi(3)) then
uzp = f(i,j,k+1,2)/f(i,j,k+1,1)
uzm = f(i,j,k-1,2)/f(i,j,k-1,1)
uz = ddz*(uzp-uzm)
else
write(6,*) 'uz'
stop
end if
vc = f(i,j,k,3)/f(i,j,k,1)
if (i.eq.ovlo(1).and.i+2.le.ovhi(1)) then
vx2p = f(i+2,j,k,3)/f(i+2,j,k,1)
vxp = f(i+1,j,k,3)/f(i+1,j,k,1)
vx = ddx*(-vx2p+4*vxp-3*vc)
else if (i.eq.ovlo(1)) then
vxp = f(i+1,j,k,3)/f(i+1,j,k,1)
vx = 2*ddx*(vxp-vc)
else if (i.eq.ovhi(1).and.i-2.ge.ovlo(1)) then
vx2m = f(i-2,j,k,3)/f(i-2,j,k,1)
vxm = f(i-1,j,k,3)/f(i-1,j,k,1)
vx = ddx*(+vx2m-4*vxm+3*vc)
else if (i.eq.ovhi(1)) then
vxm = f(i-1,j,k,3)/f(i-1,j,k,1)
vx = 2*ddx*(-vxm+vc)
else if (i.gt.ovlo(1).and.i.lt.ovhi(1)) then
vxp = f(i+1,j,k,3)/f(i+1,j,k,1)
vxm = f(i-1,j,k,3)/f(i-1,j,k,1)
vx = ddx*(vxp-vxm)
else
write(6,*) 'vx'
stop
end if
if (k.eq.ovlo(3).and.k+2.le.ovhi(3)) then
vz2p = f(i,j,k+2,3)/f(i,j,k+2,1)
vzp = f(i,j,k+1,3)/f(i,j,k+1,1)
vz = ddz*(-vz2p+4*vzp-3*vc)
else if (k.eq.ovlo(3)) then
vzp = f(i,j,k+1,3)/f(i,j,k+1,1)
vz = 2*ddz*(vzp-vc)
else if (k.eq.ovhi(3).and.k-2.ge.ovlo(3)) then
vz2m = f(i,j,k-2,3)/f(i,j,k-2,1)
vzm = f(i,j,k-1,3)/f(i,j,k-1,1)
vz = ddz*(+vz2m-4*vzm+3*vc)
else if (k.eq.ovhi(3)) then
vzm = f(i,j,k-1,3)/f(i,j,k-1,1)
vz = 2*ddz*(-vzm+vc)
else if (k.gt.ovlo(3).and.k.lt.ovhi(3)) then
vzp = f(i,j,k+1,3)/f(i,j,k+1,1)
vzm = f(i,j,k-1,3)/f(i,j,k-1,1)
vz = ddz*(vzp-vzm)
else
write(6,*) 'vz'
stop
end if
wc = f(i,j,k,4)/f(i,j,k,1)
if (i.eq.ovlo(1).and.i+2.le.ovhi(1)) then
wx2p = f(i+2,j,k,4)/f(i+2,j,k,1)
wxp = f(i+1,j,k,4)/f(i+1,j,k,1)
wx = ddx*(-wx2p+4*wxp-3*wc)
else if (i.eq.ovlo(1)) then
wxp = f(i+1,j,k,4)/f(i+1,j,k,1)
wx = 2*ddx*(wxp-wc)
else if (i.eq.ovhi(1).and.i-2.ge.ovlo(1)) then
wx2m = f(i-2,j,k,4)/f(i-2,j,k,1)
wxm = f(i-1,j,k,4)/f(i-1,j,k,1)
wx = ddx*(+wx2m-4*wxm+3*wc)
else if (i.eq.ovhi(1)) then
wxm = f(i-1,j,k,4)/f(i-1,j,k,1)
wx = 2*ddx*(-wxm+wc)
else if (i.gt.ovlo(1).and.i.lt.ovhi(1)) then
wxp = f(i+1,j,k,4)/f(i+1,j,k,1)
wxm = f(i-1,j,k,4)/f(i-1,j,k,1)
wx = ddx*(wxp-wxm)
else
write(6,*) 'wx'
stop
end if
if (j.eq.ovlo(2).and.j+2.le.ovhi(2)) then
wy2p = f(i,j+2,k,4)/f(i,j+2,k,1)
wyp = f(i,j+1,k,4)/f(i,j+1,k,1)
wy = ddy*(-wy2p+4*wyp-3*wc)
else if (j.eq.ovlo(2)) then
wyp = f(i,j+1,k,4)/f(i,j+1,k,1)
wy = 2*ddy*(wyp-wc)
else if (j.eq.ovhi(2).and.j-2.ge.ovlo(2)) then
wy2m = f(i,j-2,k,4)/f(i,j-2,k,1)
wym = f(i,j-1,k,4)/f(i,j-1,k,1)
wy = ddy*(+wy2m-4*wym+3*wc)
else if (j.eq.ovhi(2)) then
wym = f(i,j-1,k,4)/f(i,j-1,k,1)
wy = 2*ddy*(-wym+wc)
else if (j.gt.ovlo(2).and.j.lt.ovhi(2)) then
wyp = f(i,j+1,k,4)/f(i,j+1,k,1)
wym = f(i,j-1,k,4)/f(i,j-1,k,1)
wy = ddy*(wyp-wym)
else
write(6,*) 'wy'
stop
end if
vorx = wy-vz
vory = uz-wx
vorz = vx-uy
CCCC dat(i,j,k) = sqrt(vorx**2+vory**2+vorz**2)
CCCC dat(i,j,k,1) = vorx
CCCC dat(i,j,k,2) = vory
CCCC dat(i,j,k,3) = vorz
dat(i,j,k,1) = sqrt(vorx**2)
dat(i,j,k,2) = sqrt(vory**2)
dat(i,j,k,3) = sqrt(vorz**2)
CCC dat(i,j,k,1) = sqrt(vorx**2+vory**2+vorz**2)
CCC dat(i,j,k,2) = sqrt(vorx**2+vory**2+vorz**2)
CCC dat(i,j,k,3) = sqrt(vorx**2+vory**2+vorz**2)
end do
end do
end do
end
c
c derive concentration variable
c
c Inputs/Outputs
c u => (const) state array
c ulo,uhi => (const) index limits of u
c nvar => (const) number of variables in state array
c xlo => (const) physical location of lo end of u
c xhi => (const) physical location of hi end of u
c dat <=> (modify) array holding derived data
c dlo,dhi => (const) index limits of dat array
c ovlo,ovhi => (const) subregion where derivation is done
c
subroutine FORT_DERCON (u,uloi1,uloi2,uloi3,uhii1,uhii2,uhii3,
$ nvar,xlo,xhi,dat,dloi1,dloi2,dloi3,dhii1,dhii2,dhii3,ovlo,ovhi)
integer uloi1,uloi2,uloi3
integer uhii1,uhii2,uhii3
integer dloi1,dloi2,dloi3
integer dhii1,dhii2,dhii3
integer ovlo(3), ovhi(3)
integer nvar
REAL_T xlo(3), xhi(3)
REAL_T u(uloi1:uhii1, uloi2:uhii2, uloi3:uhii3, nvar)
REAL_T dat(dloi1:dhii1, dloi2:dhii2, dloi3:dhii3)
integer i,j,k
do k = ovlo(3),ovhi(3)
do j = ovlo(2),ovhi(2)
do i = ovlo(1),ovhi(1)
c dat(i,j,k) = u(i,j,k,nvar)/u(i,j,k,1)
c JBB
dat(i,j,k) = u(i,j,k,nvar)
end do
end do
end do
end
c
c INTEGRATE (3D version)
c
c integral = sum { field(i,j,k) * dx * dy * dz }
c
c Inputs / Outputs:
c
c field => field of values to integrate
c DIMS(field) => index limits of field
c DIMS(valid) => index limits of valid data in field
c delta => cell size
c integral <= integral
c
subroutine FORT_INTEGRATE
+ (field, DIMS(field), DIMS(valid), delta, integral)
integer DIMDEC(field)
integer DIMDEC(valid)
REAL_T delta(3), integral
REAL_T field(DIMV(field))
integer i, j, k
REAL_T sum
sum = zero
do k = ARG_L3(valid), ARG_H3(valid)
do j = ARG_L2(valid), ARG_H2(valid)
do i = ARG_L1(valid), ARG_H1(valid)
sum = sum + field(i,j,k)
end do
end do
end do
integral = delta(1) * delta(2) * delta(3) * sum
end
c
c JPDF (3D version)
c
c Inputs / Outputs:
c
c field => field of values to bin (1 = horizontal, 2 = vertical)
c visible => flags indicating visibility of cells (0.0 no, 1.0 yes)
c DIMS(field) => index limits of field
c DIMS(valid) => index limits of visible, and of valid data in field
c hBins => number of horizontal bins
c hMax => upper limit for horizontal bins
c hMin => lower limit for horizontal bins
c vBins => number of vertical bins
c vMax => upper limit for vertical bins
c vMin => lower limit for vertical bins
c excessCount <= count of cells outside the bin limits
c jointCount <= count of cells inside the individual, two-way bins
c
subroutine FORT_JPDF
+ (field, visible, DIMS(field), DIMS(valid),
+ hBins, hMax, hMin, vBins, vMax, vMin,
+ excessCount, jointCount)
double precision foobar
integer DIMDEC(field)
integer DIMDEC(valid)
integer hBins, vBins
integer excessCount
integer jointCount
REAL_T hMax, hMin, vMax, vMin
REAL_T field(DIMV(field), 2)
c the visible box is assumed sized to the valid region
REAL_T visible(DIMV(valid))
dimension jointCount (0 : hBins * vBins - 1)
integer bin, hBin, i, j, k, vBin
REAL_T h1, h2, hVal, v1, v2, vVal
h1 = (hBins - half) / (hMax - hMin)
h2 = (- half) / (hMax - hMin)
v1 = (vBins - half) / (vMax - vMin)
v2 = (- half) / (vMax - vMin)
do k = ARG_L3(valid), ARG_H3(valid)
do j = ARG_L2(valid), ARG_H2(valid)
do i = ARG_L1(valid), ARG_H1(valid)
if (visible(i, j, k) .eq. 1) then
hVal = field(i, j, k, 1)
if (hVal .eq. hMin) then
hBin = 0
else if (hVal .eq. hMax) then
hBin = hBins - 1
else
c if either of the Min and Max limits is an extrema, then
c rounding may choose an out of range bin when none is
c intended. hence the need to deal with this case, above.
hBin = nint ((hVal - hMin) * h1 + (hMax - hVal) * h2)
end if
vVal = field(i, j, k, 2)
if (vVal .eq. vMin) then
vBin = 0
else if (vVal .eq. vMax) then
vBin = vBins - 1
else
c if either of the Min and Max limits is an extrema, then
c rounding may choose an out of range bin when none is
c intended. hence the need to deal with this case, above.
vBin = nint ((vVal - vMin) * v1 + (vMax - vVal) * v2)
end if
if (0 .le. hBin .and. hBin .lt. hBins .and.
+ 0 .le. vBin .and. vBin .lt. vBins) then
bin = hBin + hBins * vBin
jointCount(bin) = jointCount(bin) + 1
else
excessCount = excessCount + 1
end if
end if
end do
end do
end do
return
end
c
c LIMITS (3D version)
c
c maximum = max { field(i,j,k) }
c minimum = min { field(i,j,k) }
c
c Inputs / Outputs:
c
c field => field of values to integrate
c visible => flags indicating visibility of cells (0.0 no, 1.0 yes)
c DIMS(field) => index limits of field
c DIMS(valid) => index limits of visible, and of valid data in field
c maximum <= max of starting value and field cells that are visible
c minimum <= min of starting value and field cells that are visible
c
subroutine FORT_LIMITS
+ (field, visible, DIMS(field), DIMS(valid), maximum, minimum)
integer DIMDEC(field)
integer DIMDEC(valid)
integer ifound
REAL_T maximum, minimum
REAL_T field(DIMV(field))
c the visible box is assumed sized to the valid region
REAL_T visible(DIMV(valid))
integer i, j, k
do k = ARG_L3(valid), ARG_H3(valid)
do j = ARG_L2(valid), ARG_H2(valid)
do i = ARG_L1(valid), ARG_H1(valid)
if (visible(i, j, k) .eq. 1) then
maximum = max (maximum, field(i, j, k))
minimum = min (minimum, field(i, j, k))
end if
end do
end do
end do
end
c
c SUMMASS
c
c MASS = sum{ vol(i,j)*rho(i,j) }
c
c Inputs / Outputs:
c
c rho => density field
c rlo,rhi => index limits of rho aray
c lo,hi => index limits of grid interior
c delta => cell size
c sum <= total mass
c tmp => temp column array
c
subroutine FORT_SUMMASS (rho,DIMS(rho),DIMS(grid),delta,sum,tmp)
integer DIMDEC(rho)
integer DIMDEC(grid)
REAL_T sum, delta(3)
REAL_T rho(DIMV(rho))
REAL_T tmp(DIM2(grid))
integer i, j, k
REAL_T dx, dy, dz, vol
dx = delta(1)
dy = delta(2)
dz = delta(3)
vol = dx*dy*dz
do j = ARG_L2(grid), ARG_H2(grid)
tmp(j) = zero
end do
do k = ARG_L3(grid), ARG_H3(grid)
do i = ARG_L1(grid), ARG_H1(grid)
do j = ARG_L2(grid), ARG_H2(grid)
tmp(j) = tmp(j) + vol*rho(i,j,k)
end do
end do
end do
sum = zero
do j = ARG_L2(grid), ARG_H2(grid)
sum = sum + tmp(j)
end do