/
module_fr_sfire_atm.F
877 lines (726 loc) · 33 KB
/
module_fr_sfire_atm.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
!WRF:MEDIATION_LAYER:FIRE_MODEL
!*** Jan Mandel August 2007 - February 2010
!*** email: Jan.Mandel@gmail.com
! Routines dealing with the atmosphere
module module_fr_sfire_atm
use module_model_constants, only: cp,xlv,g
use module_fr_sfire_phys, only: fire_wind_height
use module_fr_sfire_util
contains
SUBROUTINE fire_tendency( &
ids,ide, kds,kde, jds,jde, & ! dimensions
ims,ime, kms,kme, jms,jme, &
its,ite, kts,kte, jts,jte, &
grnhfx,grnqfx,canhfx,canqfx, & ! heat fluxes summed up to atm grid
alfg,alfc,z1can, & ! coeffients, properties, geometry
zs,z_at_w,dz8w,mu,rho, &
rthfrten,rqvfrten) ! theta and Qv tendencies
! This routine is atmospheric physics
! it does NOT go into module_fr_sfire_phys because it is not fire physics
! taken from the code by Ned Patton, only order of arguments change to the convention here
! --- this routine takes fire generated heat and moisture fluxes and
! calculates their influence on the theta and water vapor
! --- note that these tendencies are valid at the Arakawa-A location
IMPLICIT NONE
! --- incoming variables
INTEGER , INTENT(in) :: ids,ide, kds,kde, jds,jde, &
ims,ime, kms,kme, jms,jme, &
its,ite, kts,kte, jts,jte
REAL, INTENT(in), DIMENSION( ims:ime,jms:jme ) :: grnhfx,grnqfx ! W/m^2
REAL, INTENT(in), DIMENSION( ims:ime,jms:jme ) :: canhfx,canqfx ! W/m^2
REAL, INTENT(in), DIMENSION( ims:ime,jms:jme ) :: zs ! topography (m abv sealvl)
REAL, INTENT(in), DIMENSION( ims:ime,jms:jme ) :: mu ! dry air mass (Pa)
REAL, INTENT(in), DIMENSION( ims:ime,kms:kme,jms:jme ) :: z_at_w ! m abv sealvl
REAL, INTENT(in), DIMENSION( ims:ime,kms:kme,jms:jme ) :: dz8w ! dz across w-lvl
REAL, INTENT(in), DIMENSION( ims:ime,kms:kme,jms:jme ) :: rho ! density
REAL, INTENT(in) :: alfg ! extinction depth ground fire heat (m)
REAL, INTENT(in) :: alfc ! extinction depth crown fire heat (m)
REAL, INTENT(in) :: z1can ! height of crown fire heat release (m)
! --- outgoing variables
REAL, INTENT(out), DIMENSION( ims:ime,kms:kme,jms:jme ) :: &
rthfrten, & ! theta tendency from fire (in mass units)
rqvfrten ! Qv tendency from fire (in mass units)
! --- local variables
INTEGER :: i,j,k
INTEGER :: i_st,i_en, j_st,j_en, k_st,k_en
REAL :: cp_i
REAL :: rho_i
REAL :: xlv_i
REAL :: z_w
REAL :: fact_g, fact_c
REAL :: alfg_i, alfc_i
REAL, DIMENSION( its:ite,kts:kte,jts:jte ) :: hfx,qfx
!! character(len=128)::msg
do j=jts,jte
do k=kts,min(kte+1,kde)
do i=its,ite
rthfrten(i,k,j)=0.
rqvfrten(i,k,j)=0.
enddo
enddo
enddo
! --- set some local constants
cp_i = 1./cp ! inverse of specific heat
xlv_i = 1./xlv ! inverse of latent heat
alfg_i = 1./alfg
alfc_i = 1./alfc
!!write(msg,'(8e11.3)')cp,cp_i,xlv,xlv_i,alfg,alfc,z1can
!!call message(msg)
call print_2d_stats(its,ite,jts,jte,ims,ime,jms,jme,grnhfx,'fire_tendency:grnhfx')
call print_2d_stats(its,ite,jts,jte,ims,ime,jms,jme,grnqfx,'fire_tendency:grnqfx')
! --- set loop indicies : note that
i_st = MAX(its,ids+1)
i_en = MIN(ite,ide-1)
k_st = kts
k_en = MIN(kte,kde-1)
j_st = MAX(jts,jds+1)
j_en = MIN(jte,jde-1)
! --- distribute fluxes
DO j = j_st,j_en
DO k = k_st,k_en
DO i = i_st,i_en
! --- set z (in meters above ground)
z_w = z_at_w(i,k,j) - zs(i,j) ! should be zero when k=k_st
! --- heat flux
fact_g = cp_i * EXP( - alfg_i * z_w )
IF ( z_w < z1can ) THEN
fact_c = cp_i
ELSE
fact_c = cp_i * EXP( - alfc_i * (z_w - z1can) )
END IF
hfx(i,k,j) = fact_g * grnhfx(i,j) + fact_c * canhfx(i,j)
!! write(msg,2)i,k,j,z_w,grnhfx(i,j),hfx(i,k,j)
!!2 format('hfx:',3i4,6e11.3)
!! call message(msg)
! --- vapor flux
fact_g = xlv_i * EXP( - alfg_i * z_w )
IF (z_w < z1can) THEN
fact_c = xlv_i
ELSE
fact_c = xlv_i * EXP( - alfc_i * (z_w - z1can) )
END IF
qfx(i,k,j) = fact_g * grnqfx(i,j) + fact_c * canqfx(i,j)
!! if(hfx(i,k,j).ne.0. .or. qfx(i,k,j) .ne. 0.)then
!! write(msg,1)i,k,j,hfx(i,k,j),qfx(i,k,j)
!!1 format('tend:',3i6,2e11.3)
!! call message(msg)
! endif
END DO
END DO
END DO
! --- add flux divergence to tendencies
!
! multiply by dry air mass (mu) to eliminate the need to
! call sr. calculate_phy_tend (in dyn_em/module_em.F)
DO j = j_st,j_en
DO k = k_st,k_en-1
DO i = i_st,i_en
rho_i = 1./rho(i,k,j)
rthfrten(i,k,j) = - mu(i,j) * rho_i * (hfx(i,k+1,j)-hfx(i,k,j)) / dz8w(i,k,j)
rqvfrten(i,k,j) = - mu(i,j) * rho_i * (qfx(i,k+1,j)-qfx(i,k,j)) / dz8w(i,k,j)
END DO
END DO
END DO
call print_3d_stats(its,ite,kts,kte,jts,jte,ims,ime,kms,kme,jms,jme,rthfrten,'fire_tendency:rthfrten')
call print_3d_stats(its,ite,kts,kte,jts,jte,ims,ime,kms,kme,jms,jme,rqvfrten,'fire_tendency:rqvfrten')
RETURN
END SUBROUTINE fire_tendency
!
!***
!
subroutine interpolate_atm2fire(id, & ! for debug output, <= 0 no output
ids,ide, kds,kde, jds,jde, & ! atm grid dimensions
ims,ime, kms,kme, jms,jme, &
ips,ipe,jps,jpe, &
its,ite,jts,jte, &
ifds, ifde, jfds, jfde, & ! fire grid dimensions
ifms, ifme, jfms, jfme, &
ifps, ifpe, jfps, jfpe, & ! fire patch bounds
ifts,ifte,jfts,jfte, &
ir,jr, & ! atm/fire grid ratio
u_frame, v_frame, & ! velocity frame correction
u,v, & ! atm grid arrays in
ph,phb, &
z0,zs, &
uah,vah, &
uf,vf) ! fire grid arrays out
implicit none
! Jan Mandel, October 2010
!*** purpose: interpolate winds and height
!*** arguments
integer, intent(in)::id
integer, intent(in):: &
ids,ide, kds,kde, jds,jde, & ! atm domain bounds
ims,ime, kms,kme, jms,jme, & ! atm memory bounds
ips,ipe,jps,jpe, &
its,ite,jts,jte, & ! atm tile bounds
ifds, ifde, jfds, jfde, & ! fire domain bounds
ifms, ifme, jfms, jfme, & ! fire memory bounds
ifps, ifpe, jfps, jfpe, & ! fire patch bounds
ifts,ifte,jfts,jfte, & ! fire tile bounds
ir,jr ! atm/fire grid refinement ratio
real, intent(in):: u_frame, v_frame ! velocity frame correction
real,intent(in),dimension(ims:ime,kms:kme,jms:jme)::&
u,v, & ! atm wind velocity, staggered
ph,phb ! geopotential
real,intent(in),dimension(ims:ime,jms:jme)::&
z0, & ! roughness height
zs ! terrain height
real,intent(out),dimension(ims:ime,jms:jme)::&
uah, & ! atm wind at fire_wind_height, diagnostics
vah !
real,intent(out), dimension(ifms:ifme,jfms:jfme)::&
uf,vf ! wind velocity fire grid nodes
!*** local
character(len=256)::msg
#define TDIMS its-2,ite+2,jts-2,jte+2
real, dimension(its-2:ite+2,jts-2:jte+2):: ua,va ! atm winds, interpolated over height, still staggered grid
real, dimension(its-2:ite+2,kds:kde,jts-2:jte+2):: altw,altub,altvb,hgtu,hgtv ! altitudes
integer:: i,j,k,ifts1,ifte1,jfts1,jfte1,ite1,jte1
integer::itst,itet,jtst,jtet,itsu,iteu,jtsu,jteu,itsv,itev,jtsv,jtev
integer::kdmax,its1,jts1,ips1,jps1
integer::itsou,iteou,jtsou,jteou,itsov,iteov,jtsov,jteov
real:: ground,loght,loglast,logz0,logfwh,ht,zr
real::r_nan
integer::i_nan
equivalence (i_nan,r_nan)
!*** executable
! debug init local arrays
i_nan=2147483647
ua=r_nan
va=r_nan
altw=r_nan
altub=r_nan
hgtu=r_nan
hgtv=r_nan
if(kds.ne.1)then
!$OMP CRITICAL(SFIRE_DRIVER_CRIT)
write(msg,*)'WARNING: bottom index kds=',kds,' should be 1?'
call message(msg)
!$OMP END CRITICAL(SFIRE_DRIVER_CRIT)
endif
! ^ j
! ------------ |
! | | ----> i
! u p |
! | | nodes in cell (i,j)
! ------v----- view from top
!
! v(ide,jde+1)
! -------x------
! | |
! | |
! u(ide,jde) x x x u(ide+1,jde)
! | p(ide,hde) |
! | | p=ph,phb,z0,...
! -------x------
! v(ide,jde)
!
! staggered values set u(ids:ide+1,jds:jde) v(ids:ide,jds:jde+1)
! p=ph+phb set at (ids:ide,jds:jde)
! location of u(i,j) needs p(i-1,j) and p(i,j)
! location of v(i,j) needs p(i,j-1) and p(i,j)
! *** NOTE need HALO in ph, phb
! so we can compute only u(ids+1:ide,jds:jde) v(ids:ide,jds+1,jde)
! unless we extend p at the boundary
! but because we care about the fire way in the inside it does not matter
! if the fire gets close to domain boundary the simulation is over anyway
ite1=snode(ite,ide,1)
jte1=snode(jte,jde,1)
! do this in any case to check for nans
call print_3d_stats(its,ite1,kds,kde,jts,jte,ims,ime,kms,kme,jms,jme,u,'wind U in')
call print_3d_stats(its,ite,kds,kde,jts,jte1,ims,ime,kms,kme,jms,jme,v,'wind V in')
if(fire_print_msg.gt.0)then
!$OMP CRITICAL(SFIRE_DRIVER_CRIT)
write(msg,'(a,f7.2,a)')'interpolate_atm2fire: log-interpolation of wind to',fire_wind_height,' m'
call message(msg)
!$OMP END CRITICAL(SFIRE_DRIVER_CRIT)
endif
! indexing
! for w
itst=ifval(ids.eq.its,its,its-1)
itet=ifval(ide.eq.ite,ite,ite+1)
jtst=ifval(jds.ge.jts,jts,jts-1)
jtet=ifval(jde.eq.jte,jte,jte+1)
! for u
itsu=ifval(ids.eq.its,its+1,its) ! staggered direction
iteu=ifval(ide.eq.ite,ite,ite+1) ! staggered direction
jtsu=ifval(jds.ge.jts,jts,jts-1)
jteu=ifval(jde.eq.jte,jte,jte+1)
! for v
jtsv=ifval(jds.eq.jts,jts+1,jts) ! staggered direction
jtev=ifval(jde.eq.jte,jte,jte+1) ! staggered direction
itsv=ifval(ids.ge.its,its,its-1)
itev=ifval(ide.eq.ite,ite,ite+1)
if(fire_print_msg.ge.1)then
!$OMP CRITICAL(SFIRE_DRIVER_CRIT)
write(msg,7001)'atm input ','tile',its,ite,jts,jte
call message(msg)
write(msg,7001)'altw ','tile',itst,itet,jtst,jtet
call message(msg)
!$OMP END CRITICAL(SFIRE_DRIVER_CRIT)
endif
7001 format(a,' dimensions ',a4,':',i6,' to ',i6,' by ',i6,' to ',i6)
!**********************************************************
!* *
!* find the altitude of the w points *
!* *
!**********************************************************
!! in future, replace by z8w & test if the same
kdmax=kde-1 ! max layer to interpolate from, can be less
do j = jtst,jtet
do k=kds,kdmax+1
do i = itst,itet
altw(i,k,j) = (ph(i,k,j)+phb(i,k,j))/g ! altitude of the bottom w-point
enddo
enddo
enddo
! values at u points
if(fire_print_msg.ge.1)then
!$OMP CRITICAL(SFIRE_DRIVER_CRIT)
write(msg,7001)'u interp at','tile',itsu,iteu,jtsu,jteu
call message(msg)
!$OMP END CRITICAL(SFIRE_DRIVER_CRIT)
endif
!**********************************************************
!* *
!* interpolate the altitude from w points to u points *
!* *
!**********************************************************
do j = jtsu,jteu
do k=kds,kdmax+1
do i = itsu,iteu
altub(i,k,j)= 0.5*(altw(i-1,k,j)+altw(i,k,j)) ! altitude of the bottom point under u-point
enddo
enddo
do k=kds,kdmax
do i = itsu,iteu
hgtu(i,k,j) = 0.5*(altub(i,k,j)+altub(i,k+1,j)) - altub(i,kds,j) ! height of the u-point above the ground
enddo
enddo
enddo
! values at v points
if(fire_print_msg.ge.1)then
!$OMP CRITICAL(SFIRE_DRIVER_CRIT)
write(msg,7001)'v interp at','tile',itsv,itev,jtsv,jtev
call message(msg)
!$OMP END CRITICAL(SFIRE_DRIVER_CRIT)
endif
!**********************************************************
!* *
!* interpolate the altitude from w points to v points *
!* *
!**********************************************************
do j = jtsv,jtev
do k=kds,kdmax+1
do i = itsv,itev
altvb(i,k,j)= 0.5*(altw(i,k,j-1)+altw(i,k,j)) ! altitude of the bottom point under v-point
enddo
enddo
do k=kds,kdmax
do i = itsv,itev
hgtv(i,k,j) = 0.5*(altvb(i,k,j)+altvb(i,k+1,j)) - altvb(i,kds,j) ! height of the v-point above the ground
enddo
enddo
enddo
#ifdef DEBUG_OUT
call write_array_m3(itsu,iteu,kds,kdmax,jtsu,jteu,its-2,ite+2,kds,kde,jts-2,jte+2,altub,'altub',id)
call write_array_m3(itsv,itev,kds,kdmax,jtsv,jtev,its-2,ite+2,kds,kde,jts-2,jte+2,altvb,'altvb',id)
call write_array_m3(itsu,iteu,kds,kdmax,jtsu,jteu,its-2,ite+2,kds,kde,jts-2,jte+2,hgtu,'hgtu',id)
call write_array_m3(itsv,itev,kds,kdmax,jtsv,jtev,its-2,ite+2,kds,kde,jts-2,jte+2,hgtv,'hgtv',id)
#endif
logfwh = log(fire_wind_height)
!**********************************************************
!* *
!* interpolate u vertically to fire_wind_height *
!* *
!**********************************************************
! interpolate u, staggered in X
do j = jtsu,jteu ! compute on domain by 1 smaller, otherwise z0 and ph not available
do i = itsu,iteu ! compute with halo 2
zr = 0.5*(z0(i,j)+z0(i-1,j)) ! interpolated roughness length under this u point
if(fire_wind_height > zr)then !
do k=kds,kdmax
ht = hgtu(i,k,j) ! height of this u point above the ground
if( .not. ht < fire_wind_height) then ! found layer k this point is in
loght = log(ht)
if(k.eq.kds)then ! first layer, log linear interpolation from 0 at zr
logz0 = log(zr)
ua(i,j)= u(i,k,j)*(logfwh-logz0)/(loght-logz0)
else ! log linear interpolation
loglast=log(hgtu(i,k-1,j))
ua(i,j)= u(i,k-1,j) + (u(i,k,j) - u(i,k-1,j)) * ( logfwh - loglast) / (loght - loglast)
endif
goto 10
endif
if(k.eq.kdmax)then ! last layer, still not high enough
ua(i,j)=u(i,k,j)
endif
enddo
10 continue
else ! roughness higher than the fire wind height
ua(i,j)=0.
endif
enddo
enddo
!**********************************************************
!* *
!* interpolate v vertically to fire_wind_height *
!* *
!**********************************************************
! interpolate v, staggered in Y
do j = jtsv,jtev
do i = itsv,itev
zr = 0.5*(z0(i,j-1)+z0(i,j)) ! roughness length under this v point
if(fire_wind_height > zr)then !
do k=kds,kdmax
ht = hgtv(i,k,j) ! height of this u point above the ground
if( .not. ht < fire_wind_height) then ! found layer k this point is in
loght = log(ht)
if(k.eq.kds)then ! first layer, log linear interpolation from 0 at zr
logz0 = log(zr)
va(i,j)= v(i,k,j)*(logfwh-logz0)/(loght-logz0)
else ! log linear interpolation
loglast=log(hgtv(i,k-1,j))
va(i,j)= v(i,k-1,j) + (v(i,k,j) - v(i,k-1,j)) * ( logfwh - loglast) / (loght - loglast)
endif
goto 11
endif
if(k.eq.kdmax)then ! last layer, still not high enough
va(i,j)=v(i,k,j)
endif
enddo
11 continue
else ! roughness higher than the fire wind height
va(i,j)=0.
endif
enddo
enddo
#ifdef DEBUG_OUT
! store the output for diagnostics
do j = jts,jte1
do i = its,ite1
uah(i,j)=ua(i,j)
vah(i,j)=va(i,j)
enddo
enddo
call write_array_m(its,ite1,jts,jte,ims,ime,jms,jme,uah,'uah_n',id) ! no reflection
call write_array_m(its,ite,jts,jte1,ims,ime,jms,jme,vah,'vah_n',id)
#endif
!**********************************************************
!* *
!* interpolate ua,va vertically to the fire mesh *
!* *
!**********************************************************
ips1 = ifval(ips.eq.ids,ips+1,ips)
call continue_at_boundary(1,1,0., & ! x direction
TDIMS, &! memory dims atm grid tile
ids+1,ide,jds,jde, & ! domain dims - where u defined
ips1,ipe,jps,jpe, & ! patch dims
itsu,iteu,jtsu,jteu, & ! tile dims - in nonextended direction one beyond if at patch boundary but not domain
itsou,iteou,jtsou,jteou, & ! out, where set
ua) ! array
jps1 = ifval(jps.eq.jds,jps+1,jps)
call continue_at_boundary(1,1,0., & ! y direction
TDIMS, & ! memory dims atm grid tile
ids,ide,jds+1,jde, & ! domain dims - where v wind defined
ips,ipe,jps1,jpe, & ! patch dims
itsv,itev,jtsv,jtev, & ! tile dims
itsov,iteov,jtsov,jteov, & ! where set
va) ! array
! store the output for diagnostics
do j = jts,jte1
do i = its,ite1
uah(i,j)=ua(i,j)
vah(i,j)=va(i,j)
enddo
enddo
#ifdef DEBUG_OUT
call write_array_m(itsou,iteou,jtsou,jteou,TDIMS,ua,'ua',id)
call write_array_m(itsov,iteov,jtsov,jteov,TDIMS,va,'va',id)
#endif
!$OMP CRITICAL(SFIRE_DRIVER_CRIT)
! don't have all values valid, don't check
write(msg,12)'atm mesh wind U at',fire_wind_height,' m'
call print_2d_stats(itsou,iteou,jtsou,jteou,TDIMS,ua,msg)
write(msg,12)'atm mesh wind V at',fire_wind_height,' m'
call print_2d_stats(itsov,iteov,jtsov,jteov,TDIMS,va,msg)
12 format(a,f6.2,a)
call print_2d_stats(its,ite1,jts,jte,ims,ime,jms,jme,uah,'UAH')
call print_2d_stats(its,ite,jts,jte1,ims,ime,jms,jme,vah,'VAH')
!call write_array_m(its,ite1,jts,jte,ims,ime,jms,jme,uah,'uah',id)
!call write_array_m(its,ite,jts,jte1,ims,ime,jms,jme,vah,'vah',id)
!$OMP END CRITICAL(SFIRE_DRIVER_CRIT)
! ---------------
! | F | F | F | F | Example of atmospheric and fire grid with
! |-------|-------| ir=jr=4.
! | F | F | F | F | Winds are given at the midpoints of the sides of the atmosphere grid,
! ua------z-------| interpolated to midpoints of the cells of the fine fire grid F.
! | F | F | F | F | This is (1,1) cell of atmosphere grid, and [*] is the (1,1) cell of the fire grid.
! |---------------| ua(1,1) <--> uf(0.5,2.5)
! | * | F | F | F | va(1,1) <--> vf(2.5,0.5)
! -------va------ za(1,1) <--> zf(2.5,2.5)
!
! ^ x2
! | --------va(1,2)---------
! | | | | Example of atmospheric and fire grid with
! | | | | ir=jr=1.
! | | za,zf | Winds are given at the midpoints of the sides of the atmosphere grid,
! | ua(1,1)----uf,vf-----ua(2,1) interpolated to midpoints of the cells of the (the same) fire grid
! | | (1,1) | ua(1,1) <--> uf(0.5,1)
! | | | | va(1,1) <--> vf(1,0.5)
! | | | | za(1,1) <--> zf(1,1)
! | --------va(1,1)---------
! |--------------------> x1
!
! Meshes are aligned by the lower left cell of the domain. Then in the above figure
! u = node with the ua component of the wind at (ids,jds), midpoint of side
! v = node with the va component of the wind at (ids,jds), midpoint of side
! * = fire grid node at (ifds,jfds)
! z = node with height, midpoint of cell
!
! ua(ids,jds)=uf(ifds-0.5,jfds+jr*0.5-0.5) = uf(ifds-0.5,jfds+(jr-1)*0.5)
! va(ids,jds)=vf(ifds+ir*0.5-0.5,jfds-0.5) = vf(ifds+(ir-1)*0.5,jfds-0.5)
! za(ids,jds)=zf(ifds+ir*0.5-0.5,jfds+jr*0.5-0.5) = zf(ifds+(ir-1)*0.5,jfds+(jr-1)*0.5)
! ifts1=max(snode(ifts,ifps,-1),ifds) ! go 1 beyond patch boundary but not at domain boundary
! ifte1=min(snode(ifte,ifpe,+1),ifde)
! jfts1=max(snode(jfts,jfps,-1),jfds)
! jfte1=min(snode(jfte,jfpe,+1),jfde)
call interpolate_2d( &
TDIMS, & ! memory dims atm grid tile
itsou,iteou,jtsou,jteou,& ! where set
ifms,ifme,jfms,jfme, & ! array dims fire grid
ifts,ifte,jfts,jfte,& ! dimensions on the fire grid to interpolate to
ir,jr, & ! refinement ratio
real(ids),real(jds),ifds-0.5,jfds+(jr-1)*0.5, & ! line up by lower left corner of domain
ua, & ! in atm grid
uf) ! out fire grid
call interpolate_2d( &
TDIMS, & ! memory dims atm grid tile
itsov,iteov,jtsov,jteov,& ! where set
ifms,ifme,jfms,jfme, & ! array dims fire grid
ifts,ifte,jfts,jfte,& ! dimensions on the fire grid to interpolate to
ir,jr, & ! refinement ratio
real(ids),real(jds),ifds+(ir-1)*0.5,jfds-0.5, & ! line up by lower left corner of domain
va, & ! in atm grid
vf) ! out fire grid
call print_2d_stats_vec(ifts,ifte,jfts,jfte,ifms,ifme,jfms,jfme,uf,vf,'fire wind (m/s)')
! call print_2d_stats_vec(ifts1,ifte1,jfts1,jfte1,ifms,ifme,jfms,jfme,uf,vf,'fire wind extended')
#ifdef DEBUG_OUT
call write_array_m(ifts,ifte,jfts,jfte,ifms,ifme,jfms,jfme,uf,'uf1',id)
call write_array_m(ifts,ifte,jfts,jfte,ifms,ifme,jfms,jfme,vf,'vf1',id)
! call write_array_m(ifts1,ifte1,jfts1,jfte1,ifms,ifme,jfms,jfme,uf,'uf1',id)
! call write_array_m(ifts1,ifte1,jfts1,jfte1,ifms,ifme,jfms,jfme,vf,'vf1',id)
#endif
return
end subroutine interpolate_atm2fire
!
!***
!
subroutine interpolate_wind2fire_height(id, & ! to identify debugging prints and files if needed
ids,ide, kds,kde, jds,jde, & ! atm grid dimensions
ims,ime, kms,kme, jms,jme, &
ips,ipe,jps,jpe, &
its,ite,jts,jte, &
ifds, ifde, jfds, jfde, & ! fire grid dimensions
ifms, ifme, jfms, jfme, &
ifps, ifpe, jfps, jfpe, & ! fire patch bounds
ifts,ifte,jfts,jfte, &
ir,jr, & ! atm/fire grid ratio
u_frame, v_frame, & ! velocity frame correction
u,v,ph,phb, & ! input atmospheric arrays
fz0,fwh, & ! input fire arrays
uf,vf) ! output fire arrays
implicit none
!*** arguments
integer, intent(in):: id, & ! debug identification
ids,ide, kds,kde, jds,jde, & ! atm domain bounds
ims,ime, kms,kme, jms,jme, & ! atm memory bounds
ips,ipe,jps,jpe, &
its,ite,jts,jte, & ! atm tile bounds
ifds, ifde, jfds, jfde, & ! fire domain bounds
ifms, ifme, jfms, jfme, & ! fire memory bounds
ifps, ifpe, jfps, jfpe, & ! fire patch bounds
ifts,ifte,jfts,jfte, & ! fire tile bounds
ir,jr ! atm/fire grid refinement ratio
real, intent(in):: u_frame, v_frame ! velocity frame correction
real,intent(in),dimension(ims:ime,kms:kme,jms:jme)::&
u,v, & ! atm wind velocity, staggered
ph,phb ! geopotential
real,intent(in),dimension(ifms:ifme,jfms:jfme)::&
fz0, & ! roughness height
fwh ! height to read the wind from
real,intent(out),dimension(ifms:ifme,jfms:jfme)::&
uf, & ! atm wind at fire_wind_height, diagnostics
vf !
!*** local
integer:: i,j,k,jcb,jcm,icb,icm,kdmax,kmin,kmax
integer::itst,itet,jtst,jtet
integer::iftst,iftet,jftst,jftet
real:: wjcb,wjcm,wicb,wicm,ht,i_g,loght,zr,ht_last,logwh,wh,loght_last,uk,vk,uk1,vk1,z0,logz0
real, dimension (its-1:ite+1,kds:kde,jts-1:jte+1):: z
character(len=128)::msg
!*** executable
! print *,'interpolate_wind2fire_height start, id=',id
kdmax=kde-1 ! max layer to use
! find the altitude of atm cell centers
! index bounds for cell centers - need to go one beyond if at end of tile but not domain
itst=ifval(ids.eq.its,its,its-1)
itet=ifval(ide.eq.ite,ite,ite+1)
jtst=ifval(jds.ge.jts,jts,jts-1)
jtet=ifval(jde.eq.jte,jte,jte+1)
print *,'its, ite, jts, jte =',its, ite, jts, jte
print *,'itst, itet, jtst, jtet=',itst, itet, jtst, jtet
! get altitudes
i_g = 1./g
do j = jtst,jtet
do k=kds,kdmax+1
do i = itst,itet
z(i,k,j) = (ph(i,k,j)+phb(i,k,j))*i_g ! altitude of the bottom w-point
! print *,'i,k,j=',i,k,j,'ph, phb, z=',ph(i,k,j),phb(i,k,j),z(i,k,j)
enddo
enddo
do k=kds,kdmax
do i = itst,itet
z(i,k,j) = (z(i,k,j)+z(i,k+1,j))*0.5 - z(i,kds,j) ! height of the cell center
enddo
enddo
enddo
! index bounds for fire mesh cell centers
! to prevent setting values from uninitialized memory
iftst=ifval(ifds.eq.ifts,ifts+ir/2,ifts)
iftet=ifval(ifde.eq.ifte,ifte-ir/2,ifte)
jftst=ifval(jfds.ge.jfts,jfts+jr/2,jfts)
jftet=ifval(jfde.eq.jfte,jfte-jr/2,jfte)
print *,'iftst, iftet, jftst, jftet=',iftst, iftet, jftst, jftet
! zero out first, to prevent unitialized values on strips along domain boundaries
! it would be faster but longer code to do just cleanup loop on the strips
do j = jfts,jfte
do i = ifts,ifte
uf(i,j)=0.
vf(i,j)=0.
enddo
enddo
! vertical and horizontal interpolation
kmin=kde ! init stats
kmax=kds
loop_j: do j = jftst,jftet
call coarse(j,jr,-2,jcb,wjcb) ! get interpolation coefficients from the boundary
call coarse(j,jr,ir,jcm,wjcm) ! get interpolation coefficients from the midpoint
loop_i: do i = iftst,iftet
call coarse(i,ir,-2,icb,wicb) ! get interpolation coefficients from the boundary
call coarse(i,ir,ir,icm,wicm) ! get interpolation coefficients from the midpoint
z0 = fz0(i,j) ! roughness length
wh = fwh(i,j) ! wind height
! print *,'i=',i,' j=',j,' icb=',icb,' jcb=',jcb,' z0=',z0,' wh=',wh
if(.not. wh > z0)then
goto 92
endif
ht_last=z0 ! initialize starting height of this layer
loop_k: do k=kds,kdmax ! search for layer k such that ht(k-1)<=wh<ht(k), ht(0)=z0
! interpolate height from atmospheric cell midpoints
ht=interpolate_h(its-1,ite+1,kds,kde,jts-1,jts+1,icm,k,jcm,wicm,wjcm,z)
! print *,'i=',i,' j=',j,'k=',k,' ht=',ht
if(.not. ht < wh) exit loop_k ! found layer k this point is in
ht_last = ht
enddo loop_k
if(k .gt. kdmax) then
goto 91 ! run out of vertical levels, this must be wrong
endif
kmin=min(k,kmin)
kmax=max(k,kmax)
! found layer k, ht_last < wh <= ht
logz0 = log(z0)
logwh= log(wh)
loght_last = log(ht_last)
loght = log(ht)
! interpolate u at level k from staggered coarse grid: boundary in i, midpoints in j
uk=interpolate_h(ims,ime,kms,kme,jms,jme,icb,k,jcm,wicb,wjcm,u) - u_frame
! print *,'i=',i,' j=',j,' uk=',uk
! interpolate v at level k from staggered coarse grid: midpoints in i, boundary in j
vk=interpolate_h(ims,ime,kms,kme,jms,jme,icm,k,jcb,wicm,wjcb,v) - v_frame
! print *,'i=',i,' j=',j,' vk=',vk
if(k.gt.kds)then ! interpolate u,v horizontaly at the previous layer, k-1
! interpolate u at level k-1 from staggered coarse grid: boundary in i, midpoints in j
uk1=interpolate_h(ims,ime,kms,kme,jms,jme,icb,k-1,jcm,wicb,wjcm,u)
! print *,'i=',i,' j=',j,' uk1=',uk1
! interpolate v at level k-1 from staggered coarse grid: midpoints in i, boundary in j
vk1=interpolate_h(ims,ime,kms,kme,jms,jme,icm,k-1,jcb,wicm,wjcb,v)
! print *,'i=',i,' j=',j,' uk1=',uk1
else ! there is no previous layer, wind at roughness is assumed 0
uk1=0.
vk1=0.
endif
! log interpolate the wind to wh
uf(i,j)= uk1 + (uk - uk1) * ( logwh - loght_last) / (loght - loght_last)
vf(i,j)= vk1 + (vk - vk1) * ( logwh - loght_last) / (loght - loght_last)
! print *,'i=',i,' j=',j,' uk=',uk,' vk=',vk,' uk1=',uk1,' vk1=',vk1,' uf=',uf(i,j),' vf=',vf(i,j)
enddo loop_i
enddo loop_j
! print *,'interpolate_wind2fire_height complete, id=',id
!$OMP CRITICAL(SFIRE_ATM_CRIT)
write(msg,*)'wind interpolated from layers',kmin,' to ',kmax
call message(msg)
!$OMP END CRITICAL(SFIRE_ATM_CRIT)
call print_chsum(id,ifms,ifme,1,1,jfms,jfme,ifds,ifde,1,1,jfds,jfde,ifps,ifpe,1,1,jfps,jfpe,0,0,0,uf,'uf')
call print_chsum(id,ifms,ifme,1,1,jfms,jfme,ifds,ifde,1,1,jfds,jfde,ifps,ifpe,1,1,jfps,jfpe,0,0,0,vf,'vf')
return
91 call crash('interpolate_wind2fire_height: fire wind height too large, increase kdmax or atm height')
92 continue
!$OMP CRITICAL(SFIRE_ATM_CRIT)
write(msg,*)'fz0(',i,j,')=',fz0(i,j),'fwh(',i,j,')=',fwh(i,j)
call message(msg)
!$OMP END CRITICAL(SFIRE_ATM_CRIT)
call crash('interpolate_wind2fire_height: fire wind height must be more than the roughness height')
contains
real function interpolate_h(ims,ime,kms,kme,jms,jme,ic,kc,jc,wic,wjc,a)
!*** purpose: bilinear interpolation from a(ic:ic+1,k,jc:jc+1) with weights wicm wjcm
implicit none
!*** arguments
integer, intent(in)::ims,ime,jms,kms,kme,jme,ic,kc,jc
real, intent(in)::wic,wjc,a(ims:ime,kms:kme,jms:jme)
!*** executable
interpolate_h = &
a(ic,kc,jc) *wic *wjc + &
a(ic,kc,jc+1) *wic *(1.-wjc) + &
a(ic+1,kc,jc) *(1.-wic)*wjc + &
a(ic+1,kc,jc+1)*(1.-wic)*(1.-wjc)
end function interpolate_h
subroutine coarse(ix,nr,ia,ic,w)
!*** find coarse mesh index and interpolation weight
!*** arguments
implicit none
integer, intent(in)::ix,nr,ia
integer, intent(out)::ic
real, intent(out)::w
!*** description
! given fine mesh with nr cells for each coarse mesh cell and such that
! coarse mesh node 1 is aligned with the fine mesh at (na+1)/2
! for fine mesh node ix find its coarse mesh coordinate c and return
! ic=floor(c), the index of the nearest coarse mesh node below
! w =1-(c-ic), the interpolation weight from coarse mesh node ic to fine mesh node ix
!
! Intended use:
! fine mesh nodes are always at the middle of their cells
!
! the alignment when the coarse nodes are on the boundary of coarse cells:
! |---1---|---2---|.......|--nr---| fine mesh
! 1-------------------------------2 coarse mesh
! ia = -2 because coarse node 1 is aligned with the fine mesh at -1/2 = (-2 + 1)/2
!
! the alignment when the coarse node is at the midpoint of coarse cell:
! |---1---|---2---|---3---|---4---| fine mesh, here nr=4
! |---------------1---------------| coarse mesh
! ia = nr because coarse node 1 is aligned with the fine mesh at (nr + 1)/2
! here, (4 + 1)/2 = 2.5
!
!*** local
real:: c,a
!*** executable
a = (ia + 1)*0.5 ! the location on the fine grid where coarse node 1 is aligned
c = 1 + (ix - a)/nr ! coarse mesh coordinate of ix
ic= floor(c) ! nearest coarse node to the left
w = (1 + ic) - c ! interpolation weight, 1-(c-ic)
end subroutine coarse
end subroutine interpolate_wind2fire_height
end module module_fr_sfire_atm