-
Notifications
You must be signed in to change notification settings - Fork 148
/
mo_photo.F90
1411 lines (1261 loc) · 56.4 KB
/
mo_photo.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
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
module mo_photo
!----------------------------------------------------------------------
! ... photolysis interp table and related arrays
!----------------------------------------------------------------------
use shr_kind_mod, only : r8 => shr_kind_r8
use ppgrid, only : pcols, pver, begchunk, endchunk
use cam_abortutils, only : endrun
use mo_constants, only : r2d,d2r
use ref_pres, only : num_pr_lev, ptop_ref
use pio
use cam_pio_utils, only : cam_pio_openfile
use spmd_utils, only : masterproc
use cam_logfile, only : iulog
use solar_parms_data, only : f107=>solar_parms_f107, f107a=>solar_parms_f107a
implicit none
private
public :: photo_inti, table_photo
public :: set_ub_col
public :: setcol
public :: photo_timestep_init
public :: photo_register
save
real(r8), parameter :: kg2g = 1.e3_r8
integer, parameter :: pverm = pver - 1
integer :: jno_ndx
integer :: jonitr_ndx
integer :: jho2no2_ndx
integer :: jo2_a_ndx, jo2_b_ndx
integer :: ox_ndx, o3_ndx, o3_inv_ndx, o3rad_ndx
integer :: n2_ndx, no_ndx, o2_ndx, o_ndx
integer, allocatable :: lng_indexer(:)
integer, allocatable :: sht_indexer(:)
integer, allocatable :: euv_indexer(:)
integer :: ki
integer :: last
integer :: next
integer :: n_exo_levs
real(r8) :: delp
real(r8) :: dels
real(r8), allocatable :: days(:)
real(r8), allocatable :: levs(:)
real(r8), allocatable :: o2_exo_coldens(:,:,:,:)
real(r8), allocatable :: o3_exo_coldens(:,:,:,:)
logical :: o_is_inv
logical :: o2_is_inv
logical :: n2_is_inv
logical :: o3_is_inv
logical :: no_is_inv
logical :: has_o2_col
logical :: has_o3_col
logical :: has_fixed_press
real(r8) :: max_zen_angle ! degrees
integer :: jo1d_ndx, jo3p_ndx, jno2_ndx, jn2o5_ndx
integer :: jhno3_ndx, jno3_ndx, jpan_ndx, jmpan_ndx
integer :: jo1da_ndx, jo3pa_ndx, jno2a_ndx, jn2o5a_ndx, jn2o5b_ndx
integer :: jhno3a_ndx, jno3a_ndx, jpana_ndx, jmpana_ndx, jho2no2a_ndx
integer :: jonitra_ndx
integer :: jppi_ndx, jepn1_ndx, jepn2_ndx, jepn3_ndx, jepn4_ndx, jepn6_ndx
integer :: jepn7_ndx, jpni1_ndx, jpni2_ndx, jpni3_ndx, jpni4_ndx, jpni5_ndx
logical :: do_jeuv = .false.
logical :: do_jshort = .false.
integer :: ion_rates_idx = -1
contains
!----------------------------------------------------------------------
!----------------------------------------------------------------------
subroutine photo_register
use mo_jeuv, only : nIonRates
use physics_buffer,only : pbuf_add_field, dtype_r8
! add photo-ionization rates to phys buffer for waccmx ionosphere module
call pbuf_add_field('IonRates' , 'physpkg', dtype_r8, (/pcols,pver,nIonRates/), ion_rates_idx) ! Ionization rates for O+,O2+,N+,N2+,NO+
endsubroutine photo_register
!----------------------------------------------------------------------
!----------------------------------------------------------------------
subroutine photo_inti( xs_coef_file, xs_short_file, xs_long_file, rsf_file, &
photon_file, electron_file, &
exo_coldens_file, maxzen )
!----------------------------------------------------------------------
! ... initialize photolysis module
!----------------------------------------------------------------------
use interpolate_data, only: lininterp_init, lininterp, lininterp_finish, interp_type
use chem_mods, only : phtcnt
use chem_mods, only : ncol_abs => nabscol
use chem_mods, only : rxt_tag_lst, pht_alias_lst, pht_alias_mult
use time_manager, only : get_calday
use ioFileMod, only : getfil
use mo_chem_utls, only : get_spc_ndx, get_rxt_ndx, get_inv_ndx
use mo_jlong, only : jlong_init
use mo_jshort, only : jshort_init
use mo_jeuv, only : jeuv_init, neuv
use phys_grid, only : get_ncols_p, get_rlat_all_p
use solar_irrad_data,only : has_spectrum
use photo_bkgrnd, only : photo_bkgrnd_init
use cam_history, only : addfld
implicit none
!----------------------------------------------------------------------
! ... dummy arguments
!----------------------------------------------------------------------
character(len=*), intent(in) :: xs_long_file, rsf_file
character(len=*), intent(in) :: exo_coldens_file
real(r8), intent(in) :: maxzen
! waccm
character(len=*), intent(in) :: xs_coef_file
character(len=*), intent(in) :: xs_short_file
character(len=*), intent(in) :: photon_file
character(len=*), intent(in) :: electron_file
!----------------------------------------------------------------------
! ... local variables
!----------------------------------------------------------------------
real(r8), parameter :: hPa2Pa = 100._r8
integer :: k, n
type(file_desc_t) :: ncid
type(var_desc_t) :: vid
type(interp_type) :: lat_wgts
integer :: dimid
integer :: nlat
integer :: ntimes
integer :: astat
integer :: ndx
integer :: spc_ndx
integer :: ierr
integer :: c, ncols
integer, allocatable :: dates(:)
real(r8) :: pinterp
real(r8), allocatable :: lats(:)
real(r8), allocatable :: coldens(:,:,:)
character(len=256) :: locfn
character(len=256) :: filespec
real(r8), parameter :: trop_thrshld = 1._r8 ! Pa
real(r8) :: to_lats(pcols)
if( phtcnt < 1 ) then
return
end if
! maximum solar zenith angle for which photo-chemical rates are computed
max_zen_angle = maxzen
if (.not. max_zen_angle>0._r8) then
write(iulog,*) 'photo_inti: max_zen_angle = ',max_zen_angle
call endrun ('photo_inti: max_zen_angle must be greater then zero')
end if
! jeuv_1,,, jeuv_25 --> need euv calculations --> must be waccm
! how to determine if shrt calc is needed ?? -- use top level pressure => waccm = true ? false
if ( .not. has_spectrum ) then
write(iulog,*) 'photo_inti: solar_irrad_data file needs to contain irradiance spectrum'
call endrun('photo_inti: ERROR -- solar irradiance spectrum is missing')
endif
!----------------------------------------------------------------------
! ... allocate indexers
!----------------------------------------------------------------------
allocate( lng_indexer(phtcnt),stat=astat )
if( astat /= 0 ) then
write(iulog,*) 'photo_inti: lng_indexer allocation error = ',astat
call endrun
end if
lng_indexer(:) = 0
allocate( sht_indexer(phtcnt),stat=astat )
if( astat /= 0 ) then
write(iulog,*) 'photo_inti: Failed to allocate sht_indexer; error = ',astat
call endrun
end if
sht_indexer(:) = 0
allocate( euv_indexer(neuv),stat=astat )
if( astat /= 0 ) then
write(iulog,*) 'photo_inti: Failed to allocate euv_indexer; error = ',astat
call endrun
end if
euv_indexer(:) = 0
jno_ndx = get_rxt_ndx( 'jno' )
jo2_a_ndx = get_rxt_ndx( 'jo2_a' )
jo2_b_ndx = get_rxt_ndx( 'jo2_b' )
jo1da_ndx = get_rxt_ndx( 'jo1da' )
jo3pa_ndx = get_rxt_ndx( 'jo3pa' )
jno2a_ndx = get_rxt_ndx( 'jno2a' )
jn2o5a_ndx = get_rxt_ndx( 'jn2o5a' )
jn2o5b_ndx = get_rxt_ndx( 'jn2o5b' )
jhno3a_ndx = get_rxt_ndx( 'jhno3a' )
jno3a_ndx = get_rxt_ndx( 'jno3a' )
jpana_ndx = get_rxt_ndx( 'jpana' )
jmpana_ndx = get_rxt_ndx( 'jmpana' )
jho2no2a_ndx = get_rxt_ndx( 'jho2no2a' )
jonitra_ndx = get_rxt_ndx( 'jonitra' )
jo1d_ndx = get_rxt_ndx( 'jo1d' )
jo3p_ndx = get_rxt_ndx( 'jo3p' )
jno2_ndx = get_rxt_ndx( 'jno2' )
jn2o5_ndx = get_rxt_ndx( 'jn2o5' )
jn2o5_ndx = get_rxt_ndx( 'jn2o5' )
jhno3_ndx = get_rxt_ndx( 'jhno3' )
jno3_ndx = get_rxt_ndx( 'jno3' )
jpan_ndx = get_rxt_ndx( 'jpan' )
jmpan_ndx = get_rxt_ndx( 'jmpan' )
jho2no2_ndx = get_rxt_ndx( 'jho2no2' )
jonitr_ndx = get_rxt_ndx( 'jonitr' )
jppi_ndx = get_rxt_ndx( 'jppi' )
jepn1_ndx = get_rxt_ndx( 'jepn1' )
jepn2_ndx = get_rxt_ndx( 'jepn2' )
jepn3_ndx = get_rxt_ndx( 'jepn3' )
jepn4_ndx = get_rxt_ndx( 'jepn4' )
jepn6_ndx = get_rxt_ndx( 'jepn6' )
jepn7_ndx = get_rxt_ndx( 'jepn7' )
jpni1_ndx = get_rxt_ndx( 'jpni1' )
jpni2_ndx = get_rxt_ndx( 'jpni2' )
jpni3_ndx = get_rxt_ndx( 'jpni3' )
jpni4_ndx = get_rxt_ndx( 'jpni4' )
! added to v02
jpni5_ndx = get_rxt_ndx( 'jpni5' )
ox_ndx = get_spc_ndx( 'OX' )
if( ox_ndx < 1 ) then
ox_ndx = get_spc_ndx( 'O3' )
end if
o3_ndx = get_spc_ndx( 'O3' )
o3rad_ndx = get_spc_ndx( 'O3RAD' )
o3_inv_ndx = get_inv_ndx( 'O3' )
n2_ndx = get_inv_ndx( 'N2' )
n2_is_inv = n2_ndx > 0
if( .not. n2_is_inv ) then
n2_ndx = get_spc_ndx( 'N2' )
end if
o2_ndx = get_inv_ndx( 'O2' )
o2_is_inv = o2_ndx > 0
if( .not. o2_is_inv ) then
o2_ndx = get_spc_ndx( 'O2' )
end if
no_ndx = get_spc_ndx( 'NO' )
no_is_inv = no_ndx < 1
if( no_is_inv ) then
no_ndx = get_inv_ndx( 'NO' )
end if
o3_is_inv = o3_ndx < 1
o_ndx = get_spc_ndx( 'O' )
o_is_inv = o_ndx < 1
if( o_is_inv ) then
o_ndx = get_inv_ndx( 'O' )
end if
do_jshort = o_ndx>0 .and. o2_ndx>0 .and. (o3_ndx>0.or.o3_inv_ndx>0) .and. n2_ndx>0 .and. no_ndx>0
call jeuv_init( photon_file, electron_file, euv_indexer )
do_jeuv = any(euv_indexer(:)>0)
!----------------------------------------------------------------------
! ... call module initializers
!----------------------------------------------------------------------
call jlong_init( xs_long_file, rsf_file, lng_indexer )
if (do_jeuv) then
call photo_bkgrnd_init()
call addfld('Qbkgndtot', (/ 'lev' /), 'A','cm-3 sec-1', 'background ionization rate ' )
call addfld('Qbkgnd_o1', (/ 'lev' /), 'A','cm-3 sec-1', 'background ionization rate ' )
call addfld('Qbkgnd_o2', (/ 'lev' /), 'A','cm-3 sec-1', 'background ionization rate ' )
call addfld('Qbkgnd_n2', (/ 'lev' /), 'A','cm-3 sec-1', 'background ionization rate ' )
call addfld('Qbkgnd_n1', (/ 'lev' /), 'A','cm-3 sec-1', 'background ionization rate ' )
call addfld('Qbkgnd_no', (/ 'lev' /), 'A','cm-3 sec-1', 'background ionization rate ' )
endif
if (do_jshort) then
call jshort_init( xs_coef_file, xs_short_file, sht_indexer )
endif
jho2no2_ndx = get_rxt_ndx( 'jho2no2_b' )
!----------------------------------------------------------------------
! ... check that each photorate is in short or long datasets
!----------------------------------------------------------------------
if( any( ( abs(sht_indexer(:)) + abs(lng_indexer(:)) ) == 0 ) ) then
write(iulog,*) ' '
write(iulog,*) 'photo_inti: the following photorate(s) are not in'
write(iulog,*) ' either the short or long datasets'
write(iulog,*) ' '
do ndx = 1,phtcnt
if( abs(sht_indexer(ndx)) + abs(lng_indexer(ndx)) == 0 ) then
write(iulog,*) ' ',trim( rxt_tag_lst(ndx) )
end if
end do
call endrun
end if
!----------------------------------------------------------------------
! ... output any aliased photorates
!----------------------------------------------------------------------
if( masterproc ) then
if( any( pht_alias_lst(:,1) /= ' ' ) ) then
write(iulog,*) ' '
write(iulog,*) 'photo_inti: the following short photorate(s) are aliased'
write(iulog,*) ' '
do ndx = 1,phtcnt
if( pht_alias_lst(ndx,1) /= ' ' ) then
if( pht_alias_mult(ndx,1) == 1._r8 ) then
write(iulog,*) ' ',trim(rxt_tag_lst(ndx)),' -> ',trim(pht_alias_lst(ndx,1))
else
write(iulog,*) ' ',trim(rxt_tag_lst(ndx)),' -> ',pht_alias_mult(ndx,1),'*',trim(pht_alias_lst(ndx,1))
end if
end if
end do
end if
if( any( pht_alias_lst(:,2) /= ' ' ) ) then
write(iulog,*) ' '
write(iulog,*) 'photo_inti: the following long photorate(s) are aliased'
write(iulog,*) ' '
do ndx = 1,phtcnt
if( pht_alias_lst(ndx,2) /= ' ' ) then
if( pht_alias_mult(ndx,2) == 1._r8 ) then
write(iulog,*) ' ',trim(rxt_tag_lst(ndx)),' -> ',trim(pht_alias_lst(ndx,2))
else
write(iulog,*) ' ',trim(rxt_tag_lst(ndx)),' -> ',pht_alias_mult(ndx,2),'*',trim(pht_alias_lst(ndx,2))
end if
end if
end do
end if
write(iulog,*) ' '
write(iulog,*) '*********************************************'
write(iulog,*) 'photo_inti: euv_indexer'
write(iulog,'(10i6)') euv_indexer(:)
write(iulog,*) 'photo_inti: sht_indexer'
write(iulog,'(10i6)') sht_indexer(:)
write(iulog,*) 'photo_inti: lng_indexer'
write(iulog,'(10i6)') lng_indexer(:)
write(iulog,*) '*********************************************'
write(iulog,*) ' '
endif
!----------------------------------------------------------------------
! ... check for o2, o3 absorber columns
!----------------------------------------------------------------------
if( ncol_abs > 0 ) then
spc_ndx = ox_ndx
if( spc_ndx < 1 ) then
spc_ndx = o3_ndx
end if
if( spc_ndx > 0 ) then
has_o3_col = .true.
else
has_o3_col = .false.
end if
if( ncol_abs > 1 ) then
if( o2_ndx > 1 ) then
has_o2_col = .true.
else
has_o2_col = .false.
end if
else
has_o2_col = .false.
end if
else
has_o2_col = .false.
has_o3_col = .false.
end if
if ( len_trim(exo_coldens_file) == 0 ) then
has_o2_col = .false.
has_o3_col = .false.
endif
has_abs_columns : if( has_o2_col .or. has_o3_col ) then
!-----------------------------------------------------------------------
! ... open exo coldens file
!-----------------------------------------------------------------------
filespec = trim( exo_coldens_file )
call getfil( filespec, locfn, 0 )
call cam_pio_openfile( ncid, trim(locfn), PIO_NOWRITE )
!-----------------------------------------------------------------------
! ... get grid dimensions from file
!-----------------------------------------------------------------------
! ... timing
!-----------------------------------------------------------------------
ierr = pio_inq_dimid( ncid, 'month', dimid )
ierr = pio_inq_dimlen( ncid, dimid, ntimes )
if( ntimes /= 12 ) then
call endrun('photo_inti: exo coldens is not annual period')
end if
allocate( dates(ntimes),days(ntimes),stat=astat )
if( astat /= 0 ) then
write(iulog,*) 'photo_inti: dates,days allocation error = ',astat
call endrun
end if
dates(:) = (/ 116, 214, 316, 415, 516, 615, &
716, 816, 915, 1016, 1115, 1216 /)
!-----------------------------------------------------------------------
! ... initialize the monthly day of year times
!-----------------------------------------------------------------------
do n = 1,ntimes
days(n) = get_calday( dates(n), 0 )
end do
deallocate( dates )
!-----------------------------------------------------------------------
! ... latitudes
!-----------------------------------------------------------------------
ierr = pio_inq_dimid( ncid, 'lat', dimid )
ierr = pio_inq_dimlen( ncid, dimid, nlat )
allocate( lats(nlat), stat=astat )
if( astat /= 0 ) then
write(iulog,*) 'photo_inti: lats allocation error = ',astat
call endrun
end if
ierr = pio_inq_varid( ncid, 'lat', vid )
ierr = pio_get_var( ncid, vid, lats )
lats(:nlat) = lats(:nlat) * d2r
!-----------------------------------------------------------------------
! ... levels
!-----------------------------------------------------------------------
ierr = pio_inq_dimid( ncid, 'lev', dimid )
ierr = pio_inq_dimlen( ncid, dimid, n_exo_levs )
allocate( levs(n_exo_levs), stat=astat )
if( astat /= 0 ) then
write(iulog,*) 'photo_inti: levs allocation error = ',astat
call endrun
end if
ierr = pio_inq_varid( ncid, 'lev', vid )
ierr = pio_get_var( ncid, vid, levs )
levs(:n_exo_levs) = levs(:n_exo_levs) * hPa2Pa
!-----------------------------------------------------------------------
! ... set up regridding
!-----------------------------------------------------------------------
allocate( coldens(nlat,n_exo_levs,ntimes),stat=astat )
if( astat /= 0 ) then
write(iulog,*) 'photo_inti: coldens allocation error = ',astat
call endrun
end if
if( has_o2_col ) then
allocate( o2_exo_coldens(n_exo_levs,pcols,begchunk:endchunk,ntimes),stat=astat )
if( astat /= 0 ) then
write(iulog,*) 'photo_inti: o2_exo_coldens allocation error = ',astat
call endrun
end if
ierr = pio_inq_varid( ncid, 'O2_column_density', vid )
ierr = pio_get_var( ncid, vid,coldens )
do c=begchunk,endchunk
ncols = get_ncols_p(c)
call get_rlat_all_p(c, pcols, to_lats)
call lininterp_init(lats, nlat, to_lats, ncols, 1, lat_wgts)
do n=1,ntimes
do k = 1,n_exo_levs
call lininterp(coldens(:,k,n), nlat, o2_exo_coldens(k,:,c,n), ncols, lat_wgts)
end do
end do
call lininterp_finish(lat_wgts)
enddo
end if
if( has_o3_col ) then
allocate( o3_exo_coldens(n_exo_levs,pcols,begchunk:endchunk,ntimes),stat=astat )
if( astat /= 0 ) then
write(iulog,*) 'photo_inti: o3_exo_coldens allocation error = ',astat
call endrun
end if
ierr = pio_inq_varid( ncid, 'O3_column_density', vid )
ierr = pio_get_var( ncid, vid,coldens )
do c=begchunk,endchunk
ncols = get_ncols_p(c)
call get_rlat_all_p(c, pcols, to_lats)
call lininterp_init(lats, nlat, to_lats, ncols, 1, lat_wgts)
do n=1,ntimes
do k = 1,n_exo_levs
call lininterp(coldens(:,k,n), nlat, o3_exo_coldens(k,:,c,n), ncols, lat_wgts)
end do
end do
call lininterp_finish(lat_wgts)
enddo
end if
call pio_closefile (ncid)
deallocate( coldens,stat=astat )
if( astat /= 0 ) then
write(iulog,*) 'photo_inti: failed to deallocate coldens; error = ',astat
call endrun
end if
has_fixed_press = (num_pr_lev .ne. 0)
!-----------------------------------------------------------------------
! ... setup the pressure interpolation
!-----------------------------------------------------------------------
if( has_fixed_press ) then
pinterp = ptop_ref
if( pinterp <= levs(1) ) then
ki = 1
delp = 0._r8
else
do ki = 2,n_exo_levs
if( pinterp <= levs(ki) ) then
delp = log( pinterp/levs(ki-1) )/log( levs(ki)/levs(ki-1) )
exit
end if
end do
end if
#ifdef DEBUG
if (masterproc) then
write(iulog,*) '-----------------------------------'
write(iulog,*) 'photo_inti: diagnostics'
write(iulog,*) 'ki, delp = ',ki,delp
if (ki>1) then
write(iulog,*) 'pinterp,levs(ki-1:ki) = ',pinterp,levs(ki-1:ki)
else
write(iulog,*) 'pinterp,levs(ki) = ',pinterp,levs(ki)
end if
write(iulog,*) '-----------------------------------'
endif
#endif
end if
end if has_abs_columns
end subroutine photo_inti
subroutine table_photo( photos, pmid, pdel, temper, zmid, zint, &
col_dens, zen_angle, srf_alb, lwc, clouds, &
esfact, vmr, invariants, ncol, lchnk, pbuf )
!-----------------------------------------------------------------
! ... table photorates for wavelengths > 200nm
!-----------------------------------------------------------------
use chem_mods, only : ncol_abs => nabscol, phtcnt, gas_pcnst, nfs
use chem_mods, only : pht_alias_mult, indexm
use mo_jshort, only : nsht => nj, jshort
use mo_jlong, only : nlng => numj, jlong
use mo_jeuv, only : neuv, jeuv, nIonRates
use physics_buffer, only : physics_buffer_desc, pbuf_get_field
use photo_bkgrnd, only : photo_bkgrnd_calc
use cam_history, only : outfld
use infnan, only : nan, assignment(=)
implicit none
!-----------------------------------------------------------------
! ... dummy arguments
!-----------------------------------------------------------------
integer, intent(in) :: lchnk
integer, intent(in) :: ncol
real(r8), intent(in) :: esfact ! earth sun distance factor
real(r8), intent(in) :: vmr(ncol,pver,max(1,gas_pcnst)) ! vmr
real(r8), intent(in) :: col_dens(ncol,pver,ncol_abs) ! column densities (molecules/cm^2)
real(r8), intent(in) :: zen_angle(ncol) ! solar zenith angle (radians)
real(r8), intent(in) :: srf_alb(pcols) ! surface albedo
real(r8), intent(in) :: lwc(ncol,pver) ! liquid water content (kg/kg)
real(r8), intent(in) :: clouds(ncol,pver) ! cloud fraction
real(r8), intent(in) :: pmid(pcols,pver) ! midpoint pressure (Pa)
real(r8), intent(in) :: pdel(pcols,pver) ! pressure delta about midpoint (Pa)
real(r8), intent(in) :: temper(pcols,pver) ! midpoint temperature (K)
real(r8), intent(in) :: zmid(ncol,pver) ! midpoint height (km)
real(r8), intent(in) :: zint(ncol,pver) ! interface height (km)
real(r8), intent(in) :: invariants(ncol,pver,max(1,nfs)) ! invariant densities (molecules/cm^3)
real(r8), intent(inout) :: photos(ncol,pver,phtcnt) ! photodissociation rates (1/s)
type(physics_buffer_desc),pointer :: pbuf(:)
!-----------------------------------------------------------------
! ... local variables
!-----------------------------------------------------------------
real(r8), parameter :: Pa2mb = 1.e-2_r8 ! pascals to mb
integer :: i, k, m ! indicies
integer :: astat
real(r8) :: sza
real(r8) :: alias_factor
real(r8) :: fac1(pver) ! work space for j(no) calc
real(r8) :: fac2(pver) ! work space for j(no) calc
real(r8) :: colo3(pver) ! vertical o3 column density
real(r8) :: parg(pver) ! vertical pressure array (hPa)
real(r8) :: cld_line(pver) ! vertical cloud array
real(r8) :: lwc_line(pver) ! vertical lwc array
real(r8) :: eff_alb(pver) ! effective albedo from cloud modifications
real(r8) :: cld_mult(pver) ! clould multiplier
real(r8), allocatable :: lng_prates(:,:) ! photorates matrix (1/s)
real(r8), allocatable :: sht_prates(:,:) ! photorates matrix (1/s)
real(r8), allocatable :: euv_prates(:,:) ! photorates matrix (1/s)
real(r8), allocatable :: zarg(:)
real(r8), allocatable :: tline(:) ! vertical temperature array
real(r8), allocatable :: o_den(:) ! o density (molecules/cm^3)
real(r8), allocatable :: o2_den(:) ! o2 density (molecules/cm^3)
real(r8), allocatable :: o3_den(:) ! o3 density (molecules/cm^3)
real(r8), allocatable :: no_den(:) ! no density (molecules/cm^3)
real(r8), allocatable :: n2_den(:) ! n2 density (molecules/cm^3)
real(r8), allocatable :: jno_sht(:) ! no short photorate
real(r8), allocatable :: jo2_sht(:,:) ! o2 short photorate
real(r8), pointer :: ionRates(:,:,:) ! Pointer to ionization rates for O+,O2+,N+,N2+,NO+ in pbuf (s-1 from modules mo_jeuv and mo_jshort)
integer :: n_jshrt_levs, p1, p2
real(r8) :: ideltaZkm
real(r8) :: qbktot(ncol,pver)
real(r8) :: qbko1(ncol,pver)
real(r8) :: qbko2(ncol,pver)
real(r8) :: qbkn2(ncol,pver)
real(r8) :: qbkn1(ncol,pver)
real(r8) :: qbkno(ncol,pver)
qbktot(:,:) = nan
qbko1(:,:) = nan
qbko2(:,:) = nan
qbkn2(:,:) = nan
qbkn1(:,:) = nan
qbkno(:,:) = nan
if( phtcnt < 1 ) then
return
end if
if ((.not.do_jshort) .or. (ptop_ref < 10._r8)) then
n_jshrt_levs = pver
p1 = 1
p2 = pver
else
n_jshrt_levs = pver+1
p1 = 2
p2 = pver+1
endif
allocate( zarg(n_jshrt_levs) )
allocate( tline(n_jshrt_levs) )
if (do_jshort) then
allocate( o_den(n_jshrt_levs) )
allocate( o2_den(n_jshrt_levs) )
allocate( o3_den(n_jshrt_levs) )
allocate( no_den(n_jshrt_levs) )
allocate( n2_den(n_jshrt_levs) )
allocate( jno_sht(n_jshrt_levs) )
allocate( jo2_sht(n_jshrt_levs,2) )
endif
!-----------------------------------------------------------------
! ... allocate short, long temp arrays
!-----------------------------------------------------------------
if ( do_jeuv ) then
if (neuv>0) then
allocate( euv_prates(pver,neuv),stat=astat )
if( astat /= 0 ) then
write(iulog,*) 'photo: Failed to allocate euv_prates; error = ',astat
call endrun
end if
endif
endif
if (nsht>0) then
allocate( sht_prates(n_jshrt_levs,nsht),stat=astat )
if( astat /= 0 ) then
write(iulog,*) 'photo: Failed to allocate sht_prates; error = ',astat
call endrun
end if
endif
if (nlng>0) then
allocate( lng_prates(nlng,pver),stat=astat )
if( astat /= 0 ) then
write(iulog,*) 'photo: Failed to allocate lng_prates; error = ',astat
call endrun
end if
endif
!-----------------------------------------------------------------
! ... zero all photorates
!-----------------------------------------------------------------
do m = 1,max(1,phtcnt)
do k = 1,pver
photos(:,k,m) = 0._r8
end do
end do
!------------------------------------------------------------------------------------------------------------
! Point to production rates array in physics buffer where rates will be stored for ionosphere module
! access. Also, initialize rates to zero before column loop since only daylight values are filled
!------------------------------------------------------------------------------------------------------------
if (ion_rates_idx>0) then
call pbuf_get_field(pbuf, ion_rates_idx, ionRates)
ionRates(:,:,:) = 0._r8
endif
col_loop : do i = 1,ncol
if (do_jshort) then
if ( o_is_inv ) then
o_den(p1:p2) = invariants(i,:pver,o_ndx)
else
o_den(p1:p2) = vmr(i,:pver,o_ndx) * invariants(i,:pver,indexm)
endif
if ( o2_is_inv ) then
o2_den(p1:p2) = invariants(i,:pver,o2_ndx)
else
o2_den(p1:p2) = vmr(i,:pver,o2_ndx) * invariants(i,:pver,indexm)
endif
if ( o3_is_inv ) then
o3_den(p1:p2) = invariants(i,:pver,o3_inv_ndx)
else
o3_den(p1:p2) = vmr(i,:,o3_ndx) * invariants(i,:pver,indexm)
endif
if ( n2_is_inv ) then
n2_den(p1:p2) = invariants(i,:,n2_ndx)
else
n2_den(p1:p2) = vmr(i,:pver,n2_ndx) * invariants(i,:pver,indexm)
endif
if ( no_is_inv ) then
no_den(p1:p2) = invariants(i,:pver,no_ndx)
else
no_den(p1:p2) = vmr(i,:pver,no_ndx) * invariants(i,:pver,indexm)
endif
endif
sza = zen_angle(i)*r2d
daylight : if( sza >= 0._r8 .and. sza < max_zen_angle ) then
parg(:) = Pa2mb*pmid(i,:)
colo3(:) = col_dens(i,:,1)
fac1(:) = pdel(i,:)
lwc_line(:) = lwc(i,:)
cld_line(:) = clouds(i,:)
tline(p1:p2) = temper(i,:pver)
zarg(p1:p2) = zmid(i,:pver)
if ( ptop_ref > 10._r8 ) then
if (jppi_ndx > 0 ) photos(i,:,jppi_ndx) = photos(i,:,jppi_ndx) + esfact * 0.42_r8
if (jepn1_ndx > 0 ) photos(i,:,jepn1_ndx) = photos(i,:,jepn1_ndx) + esfact * 1.4_r8
if (jepn2_ndx > 0 ) photos(i,:,jepn2_ndx) = photos(i,:,jepn2_ndx) + esfact * 3.8e-1_r8
if (jepn3_ndx > 0 ) photos(i,:,jepn3_ndx) = photos(i,:,jepn3_ndx) + esfact * 4.7e-2_r8
if (jepn4_ndx > 0 ) photos(i,:,jepn4_ndx) = photos(i,:,jepn4_ndx) + esfact * 1.1_r8
if (jepn6_ndx > 0 ) photos(i,:,jepn6_ndx) = photos(i,:,jepn6_ndx) + esfact * 8.0e-4_r8
if (jepn7_ndx > 0 ) photos(i,:,jepn7_ndx) = photos(i,:,jepn7_ndx) + esfact * 5.2e-2_r8
if (jpni1_ndx > 0 ) photos(i,:,jpni1_ndx) = photos(i,:,jpni1_ndx) + esfact * 0.47_r8
if (jpni2_ndx > 0 ) photos(i,:,jpni2_ndx) = photos(i,:,jpni2_ndx) + esfact * 0.24_r8
if (jpni3_ndx > 0 ) photos(i,:,jpni3_ndx) = photos(i,:,jpni3_ndx) + esfact * 0.15_r8
if (jpni4_ndx > 0 ) photos(i,:,jpni4_ndx) = photos(i,:,jpni4_ndx) + esfact * 6.2e-3_r8
! added to v02
if (jpni5_ndx > 0 ) photos(i,:,jpni5_ndx) = photos(i,:,jpni5_ndx) + esfact * 1.0_r8
endif
if (do_jshort) then
if ( ptop_ref > 10._r8 ) then
!-----------------------------------------------------------------
! Only for lower lid versions of CAM (i.e., not for WACCM)
! Column O3 and O2 above the top of the model
! DEK 20110224
!-----------------------------------------------------------------
ideltaZkm = 1._r8/(zint(i,1) - zint(i,2))
!-----------------------------------------------------------------
!... set density (units: molecules cm-3)
!... used for jshort
!....... Assuming a scale height of 7km for ozone
!....... Assuming a scale height of 7km for O2
!-----------------------------------------------------------------
o3_den(1) = o3_den(2)*7.0_r8 * ideltaZkm
o2_den(1) = o2_den(2)*7.0_r8 * ideltaZkm
no_den(1) = no_den(2)*0.9_r8
n2_den(1) = n2_den(2)*0.9_r8
tline(1) = tline(2) + 5.0_r8
zarg(1) = zarg(2) + (zint(i,1) - zint(i,2))
endif
!-----------------------------------------------------------------
! ... short wave length component
!-----------------------------------------------------------------
call jshort( n_jshrt_levs, sza, n2_den, o2_den, o3_den, &
no_den, tline, zarg, jo2_sht, jno_sht, sht_prates )
do m = 1,phtcnt
if( sht_indexer(m) > 0 ) then
alias_factor = pht_alias_mult(m,1)
if( alias_factor == 1._r8 ) then
photos(i,pver:1:-1,m) = sht_prates(1:pver,sht_indexer(m))
else
photos(i,pver:1:-1,m) = alias_factor * sht_prates(1:pver,sht_indexer(m))
end if
end if
end do
if( jno_ndx > 0 ) photos(i,pver:1:-1,jno_ndx) = jno_sht(1:pver)
if( jo2_a_ndx > 0 ) photos(i,pver:1:-1,jo2_a_ndx) = jo2_sht(1:pver,2)
if( jo2_b_ndx > 0 ) photos(i,pver:1:-1,jo2_b_ndx) = jo2_sht(1:pver,1)
endif
if ( do_jeuv ) then
!-----------------------------------------------------------------
! ... euv photorates do not include cloud effects ??
!-----------------------------------------------------------------
call jeuv( pver, sza, o_den, o2_den, n2_den, zarg, euv_prates )
do m = 1,neuv
if( euv_indexer(m) > 0 ) then
photos(i,:,euv_indexer(m)) = esfact * euv_prates(:,m)
endif
enddo
endif
!-----------------------------------------------------------------
! ... compute eff_alb and cld_mult -- needs to be before jlong
!-----------------------------------------------------------------
call cloud_mod( zen_angle(i), cld_line, lwc_line, fac1, srf_alb(i), &
eff_alb, cld_mult )
cld_mult(:) = esfact * cld_mult(:)
!-----------------------------------------------------------------
! ... long wave length component
!-----------------------------------------------------------------
call jlong( pver, sza, eff_alb, parg, tline, colo3, lng_prates )
do m = 1,phtcnt
if( lng_indexer(m) > 0 ) then
alias_factor = pht_alias_mult(m,2)
if( alias_factor == 1._r8 ) then
photos(i,:,m) = (photos(i,:,m) + lng_prates(lng_indexer(m),:))*cld_mult(:)
else
photos(i,:,m) = (photos(i,:,m) + alias_factor * lng_prates(lng_indexer(m),:))*cld_mult(:)
end if
end if
end do
!-----------------------------------------------------------------
! ... calculate j(no) from formula
!-----------------------------------------------------------------
if( (jno_ndx > 0) .and. (.not.do_jshort)) then
if( has_o2_col .and. has_o3_col ) then
fac1(:) = 1.e-8_r8 * (abs(col_dens(i,:,2)/cos(zen_angle(i))))**.38_r8
fac2(:) = 5.e-19_r8 * abs(col_dens(i,:,1)/cos(zen_angle(i)))
photos(i,:,jno_ndx) = photos(i,:,jno_ndx) + 4.5e-6_r8 * exp( -(fac1(:) + fac2(:)) )
end if
end if
!-----------------------------------------------------------------
! ... add near IR correction to ho2no2
!-----------------------------------------------------------------
if( jho2no2_ndx > 0 ) then
photos(i,:,jho2no2_ndx) = photos(i,:,jho2no2_ndx) + 1.e-5_r8*cld_mult(:)
endif
! Save photo-ionization rates to physics buffer accessed in ionosphere module for WACCMX
if (ion_rates_idx>0) then
ionRates(i,1:pver,1:nIonRates) = esfact * euv_prates(1:pver,1:nIonRates)
endif
end if daylight
if (do_jeuv) then
!-----------------------------------------------------------------
! include background ionization ...
! outside daylight block so this is applied in all columns
!-----------------------------------------------------------------
call photo_bkgrnd_calc( f107, o_den, o2_den, n2_den, no_den, zint(i,:),&
photos(i,:,:), qbko1_out=qbko1(i,:), qbko2_out=qbko2(i,:), &
qbkn2_out=qbkn2(i,:), qbkn1_out=qbkn1(i,:), qbkno_out=qbkno(i,:) )
endif
end do col_loop
if ( do_jeuv ) then
qbktot(:ncol,:) = qbko1(:ncol,:)+qbko2(:ncol,:)+qbkn2(:ncol,:)+qbkn1(:ncol,:)+qbkno(:ncol,:)
call outfld( 'Qbkgndtot', qbktot(:ncol,:),ncol, lchnk )
call outfld( 'Qbkgnd_o1', qbko1(:ncol,:), ncol, lchnk )
call outfld( 'Qbkgnd_o2', qbko2(:ncol,:), ncol, lchnk )
call outfld( 'Qbkgnd_n2', qbkn2(:ncol,:), ncol, lchnk )
call outfld( 'Qbkgnd_n1', qbkn1(:ncol,:), ncol, lchnk )
call outfld( 'Qbkgnd_no', qbkno(:ncol,:), ncol, lchnk )
endif
if ( allocated(lng_prates) ) deallocate( lng_prates )
if ( allocated(sht_prates) ) deallocate( sht_prates )
if ( allocated(euv_prates) ) deallocate( euv_prates )
if ( allocated(zarg) ) deallocate( zarg )
if ( allocated(tline) ) deallocate( tline )
if ( allocated(o_den) ) deallocate( o_den )
if ( allocated(o2_den) ) deallocate( o2_den )
if ( allocated(o3_den) ) deallocate( o3_den )
if ( allocated(no_den) ) deallocate( no_den )
if ( allocated(n2_den) ) deallocate( n2_den )
if ( allocated(jno_sht) ) deallocate( jno_sht )
if ( allocated(jo2_sht) ) deallocate( jo2_sht )
call set_xnox_photo( photos, ncol )
end subroutine table_photo
subroutine cloud_mod( zen_angle, clouds, lwc, delp, srf_alb, &
eff_alb, cld_mult )
!-----------------------------------------------------------------------
! ... cloud alteration factors for photorates and albedo
!-----------------------------------------------------------------------
implicit none
!-----------------------------------------------------------------------
! ... dummy arguments
!-----------------------------------------------------------------------
real(r8), intent(in) :: zen_angle ! zenith angle
real(r8), intent(in) :: srf_alb ! surface albedo
real(r8), intent(in) :: clouds(pver) ! cloud fraction
real(r8), intent(in) :: lwc(pver) ! liquid water content (mass mr)
real(r8), intent(in) :: delp(pver) ! del press about midpoint in pascals
real(r8), intent(out) :: eff_alb(pver) ! effective albedo
real(r8), intent(out) :: cld_mult(pver) ! photolysis mult factor
!-----------------------------------------------------------------------
! ... local variables
!-----------------------------------------------------------------------
integer :: k
real(r8) :: coschi
real(r8) :: del_lwp(pver)
real(r8) :: del_tau(pver)
real(r8) :: above_tau(pver)
real(r8) :: below_tau(pver)
real(r8) :: above_cld(pver)
real(r8) :: below_cld(pver)
real(r8) :: above_tra(pver)
real(r8) :: below_tra(pver)
real(r8) :: fac1(pver)
real(r8) :: fac2(pver)
real(r8), parameter :: rgrav = 1._r8/9.80616_r8
!---------------------------------------------------------
! ... modify lwc for cloud fraction and form
! liquid water path for each layer
!---------------------------------------------------------
where( clouds(:) /= 0._r8 )
del_lwp(:) = rgrav * lwc(:) * delp(:) * 1.e3_r8 / clouds(:)
elsewhere
del_lwp(:) = 0._r8
endwhere
!---------------------------------------------------------
! ... form tau for each model layer
!---------------------------------------------------------
where( clouds(:) /= 0._r8 )
del_tau(:) = del_lwp(:) *.155_r8 * clouds(:)**1.5_r8
elsewhere
del_tau(:) = 0._r8
end where
!---------------------------------------------------------
! ... form integrated tau from top down
!---------------------------------------------------------
above_tau(1) = 0._r8
do k = 1,pverm
above_tau(k+1) = del_tau(k) + above_tau(k)
end do
!---------------------------------------------------------
! ... form integrated tau from bottom up
!---------------------------------------------------------
below_tau(pver) = 0._r8
do k = pverm,1,-1
below_tau(k) = del_tau(k+1) + below_tau(k+1)
end do
!---------------------------------------------------------
! ... form vertically averaged cloud cover above and below
!---------------------------------------------------------
above_cld(1) = 0._r8
do k = 1,pverm
above_cld(k+1) = clouds(k) * del_tau(k) + above_cld(k)
end do
do k = 2,pver
if( above_tau(k) /= 0._r8 ) then
above_cld(k) = above_cld(k) / above_tau(k)
else
above_cld(k) = above_cld(k-1)
end if
end do
below_cld(pver) = 0._r8
do k = pverm,1,-1
below_cld(k) = clouds(k+1) * del_tau(k+1) + below_cld(k+1)
end do
do k = pverm,1,-1
if( below_tau(k) /= 0._r8 ) then
below_cld(k) = below_cld(k) / below_tau(k)