-
Notifications
You must be signed in to change notification settings - Fork 6
/
Copy pathmodule_optical_averaging.F
5409 lines (5058 loc) · 221 KB
/
module_optical_averaging.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
!************************************************************************
! This computer software was prepared by Battelle Memorial Institute,
! hereinafter the Contractor, under Contract No. DE-AC05-76RL0 1830 with
! the Department of Energy (DOE). NEITHER THE GOVERNMENT NOR THE
! CONTRACTOR MAKES ANY WARRANTY, EXPRESS OR IMPLIED, OR ASSUMES ANY
! LIABILITY FOR THE USE OF THIS SOFTWARE.
!
! Module to Compute Aerosol Optical Properties
! * Author: Jerome D. Fast
! * Originators of parts of code:
! Rahul A. Zaveri, Jim Barnard, Richard C. Easter, William I.
! Gustafson Jr.
! Last update: February 2009
!
! Contact:
! Jerome D. Fast, PhD
! Staff Scientist
! Pacific Northwest National Laboratory
! P.O. Box 999, MSIN K9-30
! Richland, WA, 99352
! Phone: (509) 372-6116
! Email: Jerome.Fast@pnl.gov
!
! Please report any bugs or problems to Jerome Fast, the WRF-chem
! implmentation team leader for PNNL.
!
! Terms of Use:
! 1) Users are requested to consult the primary author prior to
! modifying this module or incorporating it or its submodules in
! another code. This is meant to ensure that the any linkages and/or
! assumptions will not adversely affect the operation of this module.
! 2) The source code in this module is intended for research and
! educational purposes. Users are requested to contact the primary
! author regarding the use of the code for any commercial application.
! 3) Users preparing publications resulting from the usage of this code
! are requested to cite one or more of the references below
! (depending on the application) for proper acknowledgement.
!
! References:
! * Fast, J.D., W.I. Gustafson Jr., R.C. Easter, R.A. Zaveri, J.C.
! Barnard, E.G. Chapman, G.A. Grell, and S.E. Peckham (2005), Evolution
! of ozone, particulates, and aerosol direct radiative forcing in the
! vicinity of Houston using a fully-coupled meteorology-chemistry-
! aerosol model. JGR, 111, doi:10.1029/2005JD006721.
! * Barnard, J.C., J.D. Fast, G. Paredes-Miranda, W.P. Arnott, and
! A. Laskin (2010), Technical note: evaluation of the WRF-Chem
! "aerosol chemical to aerosol optical properties" module using data
! from the MILAGRO campaign, Atmos. Chem. Phys., 10, 7325-7340,
! doi:10.5194/acp-10-7325-2010.
!
! Contact Jerome Fast for updates on the status of manuscripts under
! review.
!
! Additional information:
! * www.pnl.gov/atmospheric/research/wrf-chem
!
! Support:
! Funding for this code development was provided by the U.S. Department
! of Energy under the auspices of Atmospheric Science Program of the
! Office of Biological and Environmental Research the PNNL Laboratory
! Research and Directed Research and Development program.
!************************************************************************
module module_optical_averaging
USE module_data_rrtmgaeropt
implicit none
integer, parameter, private :: lunerr = -1
contains
!----------------------------------------------------------------------------------
! Aerosol optical properties computed using three methods (option_method):
! 1) volume averaging mixing rule: method that assumes internal-mixing of aerosol
! composition that averages the refractive indices for each size bin
! 2) Maxwell-Garnett mixing rule: method that randomly distributes black carbon
! within a particle
! 3) shell-core: method that assumes a "core" composed of black carbon surrounded
! by a "shell" composed of all other compositions
!
! There are two Mie routines included (option_mie):
! 1) subroutine mieaer: Employs a Chebyshev economization (Fast et al. 2006, Ghan
! et al. (2001) so that full Mie computations are called only once and then
! expansion coeffiecients are used for subsequent times to save CPU. This
! method is somewhat less accurate than full Mie calculation.
! 2) subroutine mieaer_sc: Full Mie calculation at each time step that also
! permits computation of shell-core method.
!
! Sectional and modal size distributions are treated similary, but there is
! separate code currrently to handle differences between MOSAIC and MADE/SORGAM.
!
! Methodology for sectional:
! * 3-D arrays for refractive index, wet radius, and aerosol number produced by
! optical_prep_sectional are then passed into mieaer_sectional
! * subroutine mieaer produces vertical profiles of aerosol optical properties for
! 4 wavelengths that are put into 3-D arrays and passed back up to chem_driver.F
! * tauaer*, waer*, gaer* passed to module_ra_gsfcsw.F
! * tauaer*, waer*, gaer*, l2-l7 passed to module_phot_fastj.F
! Methodology for modal:
! * similar to sectional, except divide modal mass into discrete size bins first
! * currently assume same 8 size bins as MOSAIC, but other bins are possible
!
! THIS CODE IS STILL BEING TESTED. USERS ARE ENCOURAGED TO USE ONLY
! AER_OP_OPT=1
!
subroutine optical_averaging(id,curr_secs,dtstep,config_flags, &
nbin_o,haveaer,option_method,option_mie,chem,dz8w,alt, &
relhum,is_gc, &
tauaersw,extaersw,gaersw,waersw,bscoefsw, &
l2aer,l3aer,l4aer,l5aer,l6aer,l7aer, &
tauaerlw,extaerlw, &
ids,ide, jds,jde, kds,kde, &
ims,ime, jms,jme, kms,kme, &
its,ite, jts,jte, kts,kte )
!----------------------------------------------------------------------------------
USE module_configure
USE module_state_description
USE module_model_constants
! USE module_data_mosaic_therm, only: nbin_a, nbin_a_maxd
USE module_peg_util, only: peg_error_fatal, peg_message
IMPLICIT NONE
!
INTEGER, INTENT(IN ) :: id, &
ids,ide, jds,jde, kds,kde, &
ims,ime, jms,jme, kms,kme, &
its,ite, jts,jte, kts,kte
INTEGER, INTENT(IN ) :: nbin_o
REAL(KIND=8), INTENT(IN ) :: curr_secs
REAL, INTENT(IN ) :: dtstep
LOGICAL, INTENT(IN ) :: is_gc
!
! array that holds all advected chemical species
!
REAL, DIMENSION( ims:ime, kms:kme, jms:jme, num_chem ), &
INTENT(IN ) :: chem
!
REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), &
INTENT(IN ) :: relhum,dz8w, alt
integer nspint
parameter ( nspint = 4 ) ! number of spectral interval shortwave bands
REAL, DIMENSION( ims:ime, kms:kme, jms:jme,1:nspint ), &
INTENT(INOUT ) :: &
tauaersw,extaersw,gaersw,waersw,bscoefsw
REAL, DIMENSION( ims:ime, kms:kme, jms:jme, 1:nspint ), &
INTENT(INOUT ) :: &
l2aer, l3aer, l4aer, l5aer, l6aer, l7aer
REAL, DIMENSION( ims:ime, kms:kme, jms:jme,1:nlwbands), &
INTENT(INOUT ) :: &
tauaerlw,extaerlw
!
TYPE(grid_config_rec_type), INTENT(IN ) :: config_flags
LOGICAL, INTENT(IN) :: haveaer
!
! local variables
!
real, dimension( its:ite, kts:kte, jts:jte, 1:nbin_o ) :: &
radius_wet, number_bin, radius_core
real, dimension( 1:nbin_o, kts:kte) :: &
radius_wet_col, number_bin_col, radius_core_col
complex, dimension( its:ite, kts:kte, jts:jte, 1:nbin_o ) :: & !for gocart
refindx0, refindx_core0, refindx_shell0
complex, dimension( its:ite, kts:kte, jts:jte, 1:nbin_o,1:nspint) :: &
swrefindx, swrefindx_core, swrefindx_shell
complex, dimension( 1:nbin_o, kts:kte,1:nspint) :: &
swrefindx_col, swrefindx_core_col, swrefindx_shell_col
complex, dimension( 1:nbin_o, kts:kte) :: &
swrefindx_col1, swrefindx_core_col1, swrefindx_shell_col1
complex, dimension( its:ite, kts:kte, jts:jte, 1:nbin_o,1:nlwbands) :: &
lwrefindx, lwrefindx_core, lwrefindx_shell
complex, dimension( 1:nbin_o, kts:kte,1:nlwbands) :: &
lwrefindx_col, lwrefindx_core_col, lwrefindx_shell_col
real, dimension( kts:kte ) :: dz
!
integer iclm, jclm, k, isize
integer option_method, option_mie
real, dimension( nspint, kts:kte ) :: &
swsizeaer,swextaer,swwaer,swgaer,swtauaer,swbscoef
real, dimension( nspint, kts:kte ) :: &
l2, l3, l4, l5, l6, l7
real, dimension( nlwbands, kts:kte ) :: &
lwtauaer,lwextaer
real refr
integer ns
real fv
complex aa, bb
character*150 msg
!----------------------------------------------------------------------------------
!
! write( msg, '(a, 6i4)' ) &
! 'jdf ', ids, ide, jds, jde, kds, kde
! call peg_message( lunerr, msg )
! write( msg, '(a, 6i4)' ) &
! 'jdf ', ims, ime, jms, jme, kms, kme
! call peg_message( lunerr, msg )
! write( msg, '(a, 6i4)' ) &
! 'jdf ', its, ite, jts, jte, kts, kte
! call peg_message( lunerr, msg )
chem_select: SELECT CASE(config_flags%chem_opt)
!
CASE (233)
call optical_prep_sectional_wrfgc(nbin_o, is_gc, chem, alt,relhum, &
radius_core,radius_wet, number_bin, &
swrefindx,swrefindx_core, swrefindx_shell, &
lwrefindx,lwrefindx_core, lwrefindx_shell, &
ids,ide, jds,jde, kds,kde, &
ims,ime, jms,jme, kms,kme, &
its,ite, jts,jte, kts,kte )
!gocart is now wavelength dependent --SAM 4/25/11
!and for shortwave and longwave - similar to modal --SAM 4/25/11
END SELECT chem_select
do jclm = jts, jte
do iclm = its, ite
do k = kts, kte
dz(k) = dz8w(iclm, k, jclm) ! cell depth (m)
end do
do k = kts, kte
do isize = 1, nbin_o
number_bin_col(isize,k) = number_bin(iclm,k,jclm,isize)
radius_wet_col(isize,k) = radius_wet(iclm,k,jclm,isize)
swrefindx_col(isize,k,:) = swrefindx(iclm,k,jclm,isize,:)
swrefindx_col1(isize,k) = swrefindx(iclm,k,jclm,isize,3) ! at 600 nm
lwrefindx_col(isize,k,:) = lwrefindx(iclm,k,jclm,isize,:)
radius_core_col(isize,k) = radius_core(iclm,k,jclm,isize)
swrefindx_core_col(isize,k,:) = swrefindx_core(iclm,k,jclm,isize,:)
swrefindx_shell_col(isize,k,:) = swrefindx_shell(iclm,k,jclm,isize,:)
swrefindx_core_col1(isize,k) = swrefindx_core(iclm,k,jclm,isize,3)
swrefindx_shell_col1(isize,k) = swrefindx_shell(iclm,k,jclm,isize,3)
lwrefindx_core_col(isize,k,:) = lwrefindx_core(iclm,k,jclm,isize,:)
lwrefindx_shell_col(isize,k,:) = lwrefindx_shell(iclm,k,jclm,isize,:)
! JCB, Feb. 20, 2008: in the case of shell/core and the use of the Mie
! routine, set the refractive index of the shell used in the printout
! equal to the actual refractive index of the shell
if(option_method.eq.3.and.option_mie.eq.2) &
swrefindx_col(isize,k,:) = swrefindx_shell(iclm,k,jclm,isize,:) ! JCB
swrefindx_col1(isize,k) = swrefindx_shell(iclm,k,jclm,isize,3) ! JCB
! JCB, Feb. 20, 2008: set core radius = 0 for very small cores; this
! prevents problems with full-blown Mie calculations that do not deal
! well with very small cores. For very small cores, the amount of
! absorption is negligible, and therefore setting the core radius to zero
! has virtually no effect on calculated optical properties
if(radius_wet_col(isize,k) < 1e-20) then
radius_core_col(isize,k)=0.0
else if(radius_core_col(isize,k)/radius_wet_col(isize,k)**3.le.0.0001) then
radius_core_col(isize,k)=0.0 ! JCB
end if
enddo
enddo
if (option_method .eq. 2) then
do k = kts, kte
do isize = 1, nbin_o
do ns=1,nspint
fv = (radius_core_col(isize,k)/radius_wet_col(isize,k))**3 ! volume fraction
aa=(swrefindx_core_col(isize,k,ns)**2+2.0*swrefindx_shell(iclm,k,jclm,isize,ns)**2)
bb=fv*(swrefindx_core_col(isize,k,ns)**2-swrefindx_shell(iclm,k,jclm,isize,ns)**2)
swrefindx_col(isize,k,ns)= swrefindx_shell(iclm,k,jclm,isize,ns)*sqrt((aa+2.0*bb)/(aa-bb))
if (ns==3) then
swrefindx_col1(isize,k)= swrefindx_shell(iclm,k,jclm,isize,ns)*sqrt((aa+2.0*bb)/(aa-bb))
endif
enddo
enddo
enddo
endif
if (option_method .le. 2) then
do k = kts, kte
do isize = 1, nbin_o
radius_core_col(isize,k) = 0.0
swrefindx_core_col(isize,k,:) = cmplx(0.0,0.0)
swrefindx_core_col1(isize,k) = cmplx(0.0,0.0)
enddo
enddo
endif
!
!!$ if(id.eq.1.and.iclm.eq.84.and.jclm.eq.52) then
!!$ print*,'jdf printout 1'
!!$ do isize = 1, nbin_o
!!$ write(*,888) isize,number_bin_col(isize,1), &
!!$ radius_wet_col(isize,1),radius_core_col(isize,1), &
!!$ real(refindx_col(isize,1)), &
!!$ imag(refindx_col(isize,1)), &
!!$ real(refindx_core_col(isize,1)), &
!!$ imag(refindx_core_col(isize,1)),dz(1)
!!$ enddo
!!$ endif
!!$ if(id.eq.2.and.iclm.eq.59.and.jclm.eq.63) then
!!$ print*,'jdf printout 2'
!!$ do isize = 1, nbin_o
!!$ write(*,888) isize,number_bin_col(isize,1), &
!!$ radius_wet_col(isize,1),radius_core_col(isize,1), &
!!$ real(refindx_col(isize,1)), &
!!$ imag(refindx_col(isize,1)), &
!!$ real(refindx_core_col(isize,1)), &
!!$ imag(refindx_core_col(isize,1)),dz(1)
!!$ enddo
!!$ endif
!!$ 888 format(i3,9e12.5)
!
! Initialize LW vars as not all options compute it
lwtauaer(:,:)=1.e-20
lwextaer(:,:)=1.e-20
if (option_mie .eq. 1) then
call mieaer(id, iclm, jclm, nbin_o, &
number_bin_col, radius_wet_col,swrefindx_col, &
lwrefindx_col, &
dz, curr_secs, kts, kte, &
swsizeaer,swextaer,swwaer,swgaer,swtauaer, &
lwextaer,lwtauaer, &
l2, l3, l4, l5, l6, l7,swbscoef )
endif
if (option_mie .ge. 2 .and. option_method .le. 2) then
call mieaer_sc(id, iclm, jclm, nbin_o, &
number_bin_col, radius_wet_col, swrefindx_col1, &
radius_core_col, swrefindx_core_col1, &
dz, curr_secs, kte, &
swsizeaer, swextaer, swwaer, swgaer, swtauaer, &
l2, l3, l4, l5, l6, l7, swbscoef )
endif
if (option_mie .ge. 2 .and. option_method .eq. 3) then
call mieaer_sc(id, iclm, jclm, nbin_o, &
number_bin_col, radius_wet_col, swrefindx_shell_col1, &
radius_core_col, swrefindx_core_col1, &
dz, curr_secs, kte, &
swsizeaer, swextaer, swwaer, swgaer, swtauaer, &
l2, l3, l4, l5, l6, l7, swbscoef )
endif
!
do k=kts,kte
!jdf
! write( msg, '(a, 5i4)' ) &
! 'jdf sw k', k, kts, kte, iclm, jclm
! call peg_message( lunerr, msg )
!jdf
!shortwave
! tauaersw(iclm,k,jclm,:) = swtauaer(:,k)
! extaersw(iclm,k,jclm,:) = swextaer(:,k)
! gaersw(iclm,k,jclm,:) = swgaer(:,k)
! waersw(iclm,k,jclm,:) = swwaer(:,k)
! bscoefsw(iclm,k,jclm,:) = swbscoef(:,k)
do ns=1,nspint
tauaersw(iclm,k,jclm,ns) = amax1(swtauaer(ns,k),1.e-20)
extaersw(iclm,k,jclm,ns) = amax1(swextaer(ns,k),1.e-20)
gaersw(iclm,k,jclm,ns) = amax1(amin1(swgaer(ns,k),1.0-1.e-8),1.e-20)
waersw(iclm,k,jclm,ns) = amax1(amin1(swwaer(ns,k),1.0-1.e-8),1.e-20)
bscoefsw(iclm,k,jclm,ns) = amax1(swbscoef(ns,k),1.e-20)
enddo
l2aer(iclm,k,jclm,:) = l2(:,k)
l3aer(iclm,k,jclm,:) = l3(:,k)
l4aer(iclm,k,jclm,:) = l4(:,k)
l5aer(iclm,k,jclm,:) = l5(:,k)
l6aer(iclm,k,jclm,:) = l6(:,k)
l7aer(iclm,k,jclm,:) = l7(:,k)
!longwave
do ns=1,nlwbands
tauaerlw(iclm,k,jclm,ns) = amax1(lwtauaer(ns,k),1.e-20)
extaerlw(iclm,k,jclm,ns) = amax1(lwextaer(ns,k),1.e-20)
enddo
enddo
!!$ if(id.eq.1.and.iclm.eq.84.and.jclm.eq.52) then
!!$ write(*,889) sizeaer(1,1),sizeaer(2,1),sizeaer(3,1),sizeaer(4,1)
!!$ write(*,889) extaer(1,1),extaer(2,1),extaer(3,1),extaer(4,1)
!!$ write(*,889) waer(1,1),waer(2,1),waer(3,1),waer(4,1)
!!$ write(*,889) gaer(1,1),gaer(2,1),gaer(3,1),gaer(4,1)
!!$ write(*,889) tauaer(1,1),tauaer(2,1),tauaer(3,1),tauaer(4,1)
!!$ write(*,889) bscoef(1,1),bscoef(2,1),bscoef(3,1),bscoef(4,1)
!!$ write(*,889) l2(1,1),l2(2,1),l2(3,1),l2(4,1)
!!$ write(*,889) l3(1,1),l3(2,1),l3(3,1),l3(4,1)
!!$ write(*,889) l4(1,1),l4(2,1),l4(3,1),l4(4,1)
!!$ write(*,889) l5(1,1),l5(2,1),l5(3,1),l5(4,1)
!!$ write(*,889) l6(1,1),l6(2,1),l6(3,1),l6(4,1)
!!$ write(*,889) l7(1,1),l7(2,1),l7(3,1),l7(4,1)
!!$ endif
!!$ if(id.eq.2.and.iclm.eq.59.and.jclm.eq.63) then
!!$ write(*,889) sizeaer(1,1),sizeaer(2,1),sizeaer(3,1),sizeaer(4,1)
!!$ write(*,889) extaer(1,1),extaer(2,1),extaer(3,1),extaer(4,1)
!!$ write(*,889) waer(1,1),waer(2,1),waer(3,1),waer(4,1)
!!$ write(*,889) gaer(1,1),gaer(2,1),gaer(3,1),gaer(4,1)
!!$ write(*,889) tauaer(1,1),tauaer(2,1),tauaer(3,1),tauaer(4,1)
!!$ write(*,889) bscoef(1,1),bscoef(2,1),bscoef(3,1),bscoef(4,1)
!!$ write(*,889) l2(1,1),l2(2,1),l2(3,1),l2(4,1)
!!$ write(*,889) l3(1,1),l3(2,1),l3(3,1),l3(4,1)
!!$ write(*,889) l4(1,1),l4(2,1),l4(3,1),l4(4,1)
!!$ write(*,889) l5(1,1),l5(2,1),l5(3,1),l5(4,1)
!!$ write(*,889) l6(1,1),l6(2,1),l6(3,1),l6(4,1)
!!$ write(*,889) l7(1,1),l7(2,1),l7(3,1),l7(4,1)
!!$ endif
!!$ 889 format(4e12.5)
enddo
enddo
!
return
!
end subroutine optical_averaging
!
!----------------------------------------------------------------------------------
!----------------------------------------------------------------------------------
! 4/26/11, SAM, added wavelength dependent refractive indexes for shortwave and longwave,
! 4/26/11, similar to what is done for modal CASE
! 9/21/09, SAM a modification of optical_prep_modal subroutine for GOCART aerosol model -
! SAM 7/18/09 - Modal parameters for OC1 (hydrophobic) OC2 (hydrophylic), BC1,BC2,
! and sulfate - just use dginia (meters) and sginia from module_data_sorgam.
! Not using accumulation mode from d'Almedia 1991 Table 7.1 and 7.2 global model
! 10/16/18 - A. Ukhov, bug fix: dust particles having radii in the range 0.1-0.46 microns
! were not accounted in the calculation of the mass redistribution between the GOCART and
! MOZAIC grids.
! 10/24/18 - A. Ukhov, bug fix: mass redistribution between GOCART dust/sea salt and
! MOZAIC bins should be computed using interpolation over the logarithmic axis.
! MOZAIC bins should be computed using interpolation over the logarithmic axis.
! This subroutine computes volume-averaged refractive index and wet radius needed
! by the mie calculations. Aerosol number is also passed into the mie calculations
! in terms of other units.
!
subroutine optical_prep_sectional_wrfgc(nbin_o, is_gc, chem, alt,relhum, &
radius_core,radius_wet, number_bin, &
swrefindx,swrefindx_core, swrefindx_shell, &
lwrefindx,lwrefindx_core, lwrefindx_shell, &
ids,ide, jds,jde, kds,kde, &
ims,ime, jms,jme, kms,kme, &
its,ite, jts,jte, kts,kte )
!
USE module_configure
USE module_model_constants
USE module_data_sorgam, only: dginin, dginia, dginic, sginin, sginia, sginic
! USE module_data_gocart_seas
! USE module_data_gocartchem, only: oc_mfac
! USE module_data_mosaic_asect, only: hygro_msa_aer
USE module_data_gigc_asect
USE module_data_gocart_dust, only: ra_dust, rb_dust
!
INTEGER, INTENT(IN ) :: its,ite, jts,jte, kts,kte, nbin_o
INTEGER, INTENT(IN ) :: ims,ime, jms,jme, kms,kme
INTEGER, INTENT(IN ) :: ids,ide, jds,jde, kds,kde
REAL, DIMENSION( ims:ime, kms:kme, jms:jme, num_chem ), &
INTENT(IN ) :: chem
REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), &
INTENT(IN ) :: alt,relhum
LOGICAL, INTENT( IN ) :: is_gc
REAL, DIMENSION( its:ite, kts:kte, jts:jte, 1:nbin_o), &
INTENT(OUT ) :: &
radius_wet, number_bin, radius_core
COMPLEX, DIMENSION( its:ite, kts:kte, jts:jte,1:nbin_o,nswbands), &
INTENT(OUT ) :: swrefindx, swrefindx_core, swrefindx_shell
COMPLEX, DIMENSION( its:ite, kts:kte, jts:jte,1:nbin_o,nlwbands), &
INTENT(OUT ) :: lwrefindx, lwrefindx_core, lwrefindx_shell
!
! local variables
!
integer i, j, k, l, m, n, ns, isize, itype
complex ref_index_lvcite , ref_index_nh4hso4 , &
ref_index_nh4msa , ref_index_nh4no3 , ref_index_nh4cl , &
ref_index_nano3 , ref_index_na2so4 , &
ref_index_na3hso4, ref_index_nahso4 , ref_index_namsa, &
ref_index_caso4 , ref_index_camsa2 , ref_index_cano3, &
ref_index_cacl2 , ref_index_caco3 , ref_index_h2so4, &
ref_index_hhso4 , ref_index_hno3 , ref_index_hcl, &
ref_index_msa , ref_index_bc, &
ref_index_oin , ref_index_aro1 , ref_index_aro2, &
ref_index_alk1 , ref_index_ole1 , ref_index_api1, &
ref_index_api2 , ref_index_lim1 , ref_index_lim2, &
ri_dum , ri_ave_a
complex, dimension(nswbands) :: & ! now only 5 aerosols have wave-dependent refr
swref_index_oc , swref_index_dust , swref_index_nh4so4, swref_index_nacl,swref_index_h2o
complex, dimension(nlwbands) :: & ! now only 5 aerosols have wave-dependent refr
lwref_index_oc , lwref_index_dust , lwref_index_nh4so4, lwref_index_nacl,lwref_index_h2o
real(KIND=8) mass_so4, mass_nit, mass_nh4, mass_oc, mass_bc, mass_seas, &
mass_dst, mass_soas, mass_water
real(KIND=8) mass_so4i, mass_so4j, mass_niti, mass_nitj, &
mass_nh4i, mass_nh4j, mass_ocpii, mass_ocpij, &
mass_ocpoi, mass_ocpoj, mass_bcpii, mass_bcpij, &
mass_bcpoi, mass_bcpoj, mass_salj, mass_salc, &
mass_soasi, mass_soasj
real(KIND=8) mass_so4_gc, mass_nit_gc, mass_nh4_gc, mass_ocpo_gc, &
mass_ocpi_gc,mass_bcpo_gc, mass_bcpi_gc, mass_dusttmp, &
mass_soas_gc, mass_seas_gc, mass_dst_gc, mass_oc_gc, &
mass_bc_gc
real(KIND=8) vol_h2otmp, vol_wet_aer, volaer_tmp
real(KIND=8) conv2so4, conv2nit, conv2nh4, conv2ocbc, &
conv2dst, conv2seas, conv2soas
real(KIND=8) dlogoc, dhigoc
real(KIND=8) vol_so4, vol_nit, vol_nh4, vol_oc, vol_bc, &
vol_seas, vol_dst, vol_soas, vol_h2o
real(KIND=8) dens_h2o, mole_dryair, conv1a, conv1b
real(KIND=8) mass_dry_a, mass_wet_a, vol_dry_a , vol_wet_a , vol_shell, &
dp_dry_a , dp_wet_a , num_a , dp_bc_a , num_a_lo , &
num_a_hi
real(KIND=8) sixpi, dlo, dhi, xlo, xhi, dx
integer iflag
real, dimension(1:nbin_o) :: xdia_cm
real, dimension(1:nbin_o) :: xnum_secti, xnum_sectj, xnum_sectc, xnum_sect_sna_oc, xnum_sect_bc, xnum_sect_sala, xnum_sect_salc
real, dimension(1:nbin_o) :: xmas_secti, xmas_sectj, xmas_sectc, xmas_sect_sna_oc, xmas_sect_bc, xmas_sect_sala, xmas_sect_salc
real, parameter :: pi = 3.141592653589
real, parameter :: FRAC2Aitken = 0.25 ! Fraction of modal mass in Aitken mode - applied globally to each species
real, parameter :: ndust1 = 4 ! Number of dust bins of GEOS-Chem
real :: dustfrc_gigc4bin(ndust1, nbin_o) ! GEOS-Chem dust size distribution - mass fracs in MOSAIC 4-bins
real :: dgnum_um, dlo_um, dhi_um, duma, dxbin, relh_frc, dgnum_sna_oc, dgnum_bc, sigma_sna_oc, sigma_bc, &
dgnum_sala, dgnum_salc, sigma_sala, sigma_salc
real :: dlo_sectm(nbin_o), dhi_sectm(nbin_o)
!
! real sginin,sginia,sginic from module_data_sorgam.F
!
! Mass from modal distribution is divided into individual sections before
! being passed back into the Mie routine.
! * currently use the same size bins as 8 default MOSAIC size bins
! * dlo_um and dhi_um define the lower and upper bounds of individual sections
! used to compute optical properties
! * sigmas for 3 modes taken from module_sorgam_data.F
! * these parameters are needed by sect02 that is called later
! * sginin=1.7, sginia=2.0, sginic=2.5
!
itype = 1
sixpi = 6.0/pi
dlo_um = 0.0390625 ! low diameter limit (um)
dhi_um = 10.0 ! high diameter limit (um)
duma = 1.0
dlo = dlo_um*1.0E-6 ! low diameter limit (m)
dhi = dhi_um*1.0E-6 ! high diameter limit (m)
dxbin = (log(dhi) - log(dlo))/nbin_o
do n = 1, nbin_o
dlo_sectm(n) = exp( log(dlo) + dxbin * (n-1) )
dhi_sectm(n) = exp( log(dlo) + dxbin * n )
xdia_cm(n) = 0.5 * (dlo_sectm(n) + dhi_sectm(n)) * 1.0E2 ! unit:cm
end do
! Dust bin mass fraction
dustfrc_gigc4bin = 0.
do m = 1, ndust1
dlogoc = ra_dust(m) * 2.E-6 ! low diameter limit (m)
dhigoc = rb_dust(m) * 2.E-6 ! high diameter limit (m)
do n = 1, nbin_o
dustfrc_gigc4bin(m, n) = max(DBLE(0.),min(DBLE(log(dhi_sectm(n))),log(dhigoc))- &
max(log(dlogoc),DBLE(log(dlo_sectm(n)))) )/(log(dhigoc)-log(dlogoc))
end do
end do
!
! Define refractive indicies
! * assume na and cl are the same as nacl
! * assume so4, no3, and nh4 are the same as nh4no3
! * assume ca and co3 are the same as caco3
! * assume msa is just msa
! Further work:
! * to be more precise, need to compute electrolytes to apportion
! so4, no3, nh4, na, cl, msa, ca, co3 among various componds
! as was done previously in module_mosaic_therm.F
!
do ns = 1, nswbands
swref_index_nh4so4(ns) = cmplx(refrsw_sulf(ns),refisw_sulf(ns))
swref_index_oc(ns) = cmplx(refrsw_oc(ns),refisw_oc(ns))
swref_index_dust(ns) = cmplx(refrsw_dust(ns),refisw_dust(ns))
swref_index_nacl(ns) = cmplx(refrsw_seas(ns),refisw_seas(ns))
swref_index_h2o(ns) = cmplx(refrwsw(ns),refiwsw(ns))
enddo
do ns = 1, nlwbands
lwref_index_nh4so4(ns) = cmplx(refrlw_sulf(ns),refilw_sulf(ns))
lwref_index_oc(ns) = cmplx(refrlw_oc(ns),refilw_oc(ns))
lwref_index_dust(ns) = cmplx(refrlw_dust(ns),refilw_dust(ns))
lwref_index_nacl(ns) = cmplx(refrlw_seas(ns),refilw_seas(ns))
lwref_index_h2o(ns) = cmplx(refrwlw(ns),refiwlw(ns))
enddo
! ref_index_nh4so4 = cmplx(1.52,0.)
ref_index_lvcite = cmplx(1.50,0.)
ref_index_nh4hso4= cmplx(1.47,0.)
ref_index_nh4msa = cmplx(1.50,0.) ! assumed
ref_index_nh4no3 = cmplx(1.50,0.)
ref_index_nh4cl = cmplx(1.50,0.)
! ref_index_nacl = cmplx(1.45,0.)
ref_index_nano3 = cmplx(1.50,0.)
ref_index_na2so4 = cmplx(1.50,0.)
ref_index_na3hso4= cmplx(1.50,0.)
ref_index_nahso4 = cmplx(1.50,0.)
ref_index_namsa = cmplx(1.50,0.) ! assumed
ref_index_caso4 = cmplx(1.56,0.006)
ref_index_camsa2 = cmplx(1.56,0.006) ! assumed
ref_index_cano3 = cmplx(1.56,0.006)
ref_index_cacl2 = cmplx(1.52,0.006)
ref_index_caco3 = cmplx(1.68,0.006)
ref_index_h2so4 = cmplx(1.43,0.)
ref_index_hhso4 = cmplx(1.43,0.)
ref_index_hno3 = cmplx(1.50,0.)
ref_index_hcl = cmplx(1.50,0.)
ref_index_msa = cmplx(1.43,0.) ! assumed
! ref_index_oc = cmplx(1.45,0.) ! JCB, Feb. 20, 2008: no complex part?
! JCB, Feb. 20, 2008: set the refractive index of BC equal to the
! midpoint of ranges given in Bond and Bergstrom, Light absorption by
! carboneceous particles: an investigative review 2006, Aerosol Sci.
! and Tech., 40:27-67.
! ref_index_bc = cmplx(1.82,0.74) old value
ref_index_bc = cmplx(1.85,0.71)
ref_index_oin = cmplx(1.55,0.006) ! JCB, Feb. 20, 2008: "other inorganics"
ref_index_aro1 = cmplx(1.45,0.)
ref_index_aro2 = cmplx(1.45,0.)
ref_index_alk1 = cmplx(1.45,0.)
ref_index_ole1 = cmplx(1.45,0.)
ref_index_api1 = cmplx(1.45,0.)
ref_index_api2 = cmplx(1.45,0.)
ref_index_lim1 = cmplx(1.45,0.)
ref_index_lim2 = cmplx(1.45,0.)
! ref_index_h2o = cmplx(1.33,0.)
! JCB, Feb. 20, 2008: the density of BC is updated to reflect values
! published by Bond and Bergstrom, Light absorption by carboneceous
! particles: an investigative review 2006, Aerosol Sci. and Tech., 40:27-67.
swrefindx=0.0
lwrefindx=0.0
radius_wet=0.0
number_bin=0.0
radius_core=0.0
swrefindx_core=0.0
swrefindx_shell=0.0
lwrefindx_core=0.0
lwrefindx_shell=0.0
!
! units:
! * mass - g/cc(air)
! * number - #/cc(air)
! * volume - cc(air)/cc(air)
! * diameter - cm
!
do j = jts, jte
do k = kts, kte
do i = its, ite
mass_so4i = 0.0
mass_so4j = 0.0
mass_niti = 0.0
mass_nitj = 0.0
mass_nh4i = 0.0
mass_nh4j = 0.0
mass_ocpii = 0.0
mass_ocpij = 0.0
mass_ocpoi = 0.0
mass_ocpoj = 0.0
mass_bcpii = 0.0
mass_bcpij = 0.0
mass_bcpoi = 0.0
mass_bcpoj = 0.0
mass_salj = 0.0
mass_salc = 0.0
mass_soasi = 0.0
mass_soasj = 0.0
mass_so4_gc = 0.0
mass_nit_gc = 0.0
mass_nh4_gc = 0.0
mass_oc_gc = 0.0
mass_bc_gc = 0.0
mass_seas_gc = 0.0
mass_dst_gc = 0.0
mass_soas_gc = 0.0
mass_so4 = 0.0
mass_nit = 0.0
mass_nh4 = 0.0
mass_oc = 0.0
mass_bc = 0.0
mass_seas = 0.0
mass_dst = 0.0
mass_water = 0.0
! water densities in g/cm^3
dens_h2o = 1.0
! dry air mole weight g/mol
mole_dryair = 28.97
! convert diagnostic ug/kg aerosols to g/cm^3 (cc => cm^3) air
conv1a = (1.0/alt(i,k,j)) * 1.0E-12
! convert diagnostic #/kg aerosols number concentrations to #/cm^3 (cc => cm^3) air
conv1b = (1.0/alt(i,k,j)) * 1.0E-6
! convert ppmv sulfate to ug/kg
conv2so4 = mw_so4_aer/mole_dryair * 1.0E3
! convert ppmv nitrate to ug/kg
conv2nit = mw_nit_aer/mole_dryair * 1.0E3
! convert ppmv ammonium to ug/kg
conv2nh4 = mw_nh4_aer/mole_dryair * 1.0E3
! convert ppmv oc/bc to ug/kg
conv2ocbc = mw_oc_aer/mole_dryair * 1.0E3
! convert ppmv dust to ug/kg
conv2dst = mw_dst_aer/mole_dryair * 1.0E3
! convert ppmv sea salt ug/kg
conv2seas = mw_seas_aer/mole_dryair * 1.0E3
! convert ppmv simple soa ug/kg
conv2soas = mw_soas_aer/mole_dryair * 1.0E3
! XFENG 11/15/2019 - Convert the units of diagnostic SO4, NIT, NH4, OC, BC, DUST, SEASALT, SOAS masses from ppmv to g/cm^3
if (is_gc) then
mass_so4_gc = chem(i, k, j, p_so4) * conv2so4 * conv1a
mass_nit_gc = chem(i, k, j, p_nit) * conv2nit * conv1a
mass_nh4_gc = chem(i, k, j, p_nh4) * conv2nh4 * conv1a
mass_ocpi_gc = chem(i, k, j, p_ocpi) * conv2ocbc * conv1a
mass_ocpo_gc = chem(i, k, j, p_ocpo) * conv2ocbc * conv1a
mass_bcpi_gc = chem(i, k, j, p_bcpi) * conv2ocbc * conv1a
mass_bcpo_gc = chem(i, k, j, p_bcpo) * conv2ocbc * conv1a
mass_soas_gc = chem(i, k, j, p_soas) * conv2soas * conv1a
else
! Put the part of aerosol mass into Aitken mode
mass_so4i = FRAC2Aitken * chem(i, k, j, p_so4) * conv2so4 * conv1a
mass_niti = FRAC2Aitken * chem(i, k, j, p_nit) * conv2nit * conv1a
mass_nh4i = FRAC2Aitken * chem(i, k, j, p_nh4) * conv2nh4 * conv1a
mass_ocpii = FRAC2Aitken * chem(i, k, j, p_ocpi) * conv2ocbc * conv1a
mass_ocpoi = FRAC2Aitken * chem(i, k, j, p_ocpo) * conv2ocbc * conv1a
mass_bcpii = FRAC2Aitken * chem(i, k, j, p_bcpi) * conv2ocbc * conv1a
mass_bcpoi = FRAC2Aitken * chem(i, k, j, p_bcpo) * conv2ocbc * conv1a
mass_soasi = FRAC2Aitken * chem(i, k, j, p_soas) * conv2soas * conv1a
! Put the rest of aerosol mass into Accumulation mode
mass_so4j = (1.0-FRAC2Aitken) * chem(i, k, j, p_so4) * conv2so4 * conv1a
mass_nitj = (1.0-FRAC2Aitken) * chem(i, k, j, p_nit) * conv2nit * conv1a
mass_nh4j = (1.0-FRAC2Aitken) * chem(i, k, j, p_nh4) * conv2nh4 * conv1a
mass_ocpij = (1.0-FRAC2Aitken) * chem(i, k, j, p_ocpi) * conv2ocbc * conv1a
mass_ocpoj = (1.0-FRAC2Aitken) * chem(i, k, j, p_ocpo) * conv2ocbc * conv1a
mass_bcpij = (1.0-FRAC2Aitken) * chem(i, k, j, p_bcpi) * conv2ocbc * conv1a
mass_bcpoj = (1.0-FRAC2Aitken) * chem(i, k, j, p_bcpo) * conv2ocbc * conv1a
mass_soasj = (1.0-FRAC2Aitken) * chem(i, k, j, p_soas) * conv2soas * conv1a
end if
! sea salt in accumulation and coarse mode
mass_salj = chem(i, k, j, p_sala) * conv2seas * conv1a ! SALA represents the accumulation mode sea salt aerosol
mass_salc = chem(i, k, j, p_salc) * conv2seas * conv1a ! SALC represents the coarse mode sea salt aerosol
itype = 1
do isize = 1, nbin_o
if (is_gc) then
dgnum_sna_oc = 0.07*2 ! unit:um
sigma_sna_oc = 1.6
call sect02_new(dgnum_sna_oc, sigma_sna_oc, duma, nbin_o, dlo_um, dhi_um, &
xnum_sect_sna_oc, xmas_sect_sna_oc)
dgnum_bc = 0.02*2 ! unit:um
sigma_bc = 1.6
call sect02_new(dgnum_bc, sigma_bc, duma, nbin_o, dlo_um, dhi_um, &
xnum_sect_bc, xmas_sect_bc)
! Sulfate, nitrate, ammonium, oc, bc, soas aerosol (g/cm^3)
mass_so4 = mass_so4_gc * xmas_sect_sna_oc(isize)
mass_nit = mass_nit_gc * xmas_sect_sna_oc(isize)
mass_nh4 = mass_nit_gc * xmas_sect_sna_oc(isize)
mass_oc = mass_ocpi_gc * xmas_sect_sna_oc(isize) + mass_ocpo_gc * xmas_sect_sna_oc(isize)
mass_bc = mass_bcpi_gc * xmas_sect_bc(isize) + mass_bcpo_gc * xmas_sect_bc(isize)
mass_soas = mass_soas_gc * xmas_sect_sna_oc(isize)
else
dgnum_um = dginin * 1.E6
call sect02_new(dgnum_um, sginin, duma, nbin_o, dlo_um, dhi_um, &
xnum_secti, xmas_secti)
dgnum_um = dginia * 1.E6
call sect02_new(dgnum_um, sginia, duma, nbin_o, dlo_um, dhi_um, &
xnum_sectj, xmas_sectj)
dgnum_um = dginic * 1.E6
call sect02_new(dgnum_um, sginic, duma, nbin_o, dlo_um, dhi_um, &
xnum_sectc, xmas_sectc)
mass_so4 = mass_so4i * xmas_secti(isize) + mass_so4j * xmas_sectj(isize)
mass_nit = mass_niti * xmas_secti(isize) + mass_nitj * xmas_sectj(isize)
mass_nh4 = mass_nh4i * xmas_secti(isize) + mass_nh4j * xmas_sectj(isize)
mass_oc = mass_ocpii * xmas_secti(isize) + mass_ocpij * xmas_sectj(isize) + &
mass_ocpoi * xmas_secti(isize) + mass_ocpoj * xmas_sectj(isize)
mass_bc = mass_bcpii * xmas_secti(isize) + mass_bcpij * xmas_sectj(isize) + &
mass_bcpoi * xmas_secti(isize) + mass_bcpoj * xmas_sectj(isize)
mass_soas = mass_soasi * xmas_secti(isize) + mass_soasj * xmas_sectj(isize)
end if
! Sea salt aerosol (g/cm^3)
dgnum_sala = 0.09*2 ! unit:um
sigma_sala = 1.5
call sect02_new(dgnum_sala, sigma_sala, duma, nbin_o, dlo_um, dhi_um, &
xnum_sect_sala, xmas_sect_sala)
dgnum_salc = 0.4*2 ! unit:um
sigma_salc = 1.8
call sect02_new(dgnum_salc, sigma_salc, duma, nbin_o, dlo_um, dhi_um, &
xnum_sect_salc, xmas_sect_salc)
mass_seas = mass_salj * xmas_sect_sala(isize) + mass_salc * xmas_sect_salc(isize)
! Dust aerosol (g/cm^3)
n = 0
mass_dusttmp = 0.0
do m = p_dst1, p_dst4
n = n + 1
mass_dusttmp = mass_dusttmp + dustfrc_gigc4bin(n, isize)*chem(i, k, j, m)
end do
mass_dst = mass_dusttmp * conv2dst * conv1a
! Aerosol water
relh_frc = amin1(0.9, relhum(i, k, j)) ! Put in fractional relative humidity, max of 0.9 here
vol_h2otmp = 0.0
vol_h2otmp = mass_so4 * hygro_so4_aer / dens_so4_aer_gc + &
mass_nit * hygro_nit_aer / dens_nit_aer_gc + &
mass_nh4 * hygro_nh4_aer / dens_nh4_aer_gc + &
mass_oc * hygro_oc_aer / dens_oc_aer_gc + &
mass_bc * hygro_bc_aer / dens_bc_aer_gc + &
mass_seas * hygro_seas_aer/ dens_seas_aer_gc + &
mass_soas * hygro_soas_aer/ dens_soas_aer_gc
if (isize .eq. 1) then
vol_h2otmp = vol_h2otmp + mass_dst * hygro_dst_aer / dens_dst1_aer_gc
else
vol_h2otmp = vol_h2otmp + mass_dst * hygro_dst_aer / dens_dst2_aer_gc
end if
vol_h2otmp = relh_frc * vol_h2otmp / (1. - relh_frc)
mass_water = vol_h2otmp * dens_h2o
! Calculate the volume of aerosols
vol_so4 = mass_so4 / dens_so4_aer_gc
vol_nit = mass_nit / dens_nit_aer_gc
vol_nh4 = mass_nh4 / dens_nh4_aer_gc
vol_oc = mass_oc / dens_oc_aer_gc
vol_bc = mass_bc / dens_bc_aer_gc
vol_seas = mass_seas / dens_seas_aer_gc
vol_soas = mass_soas / dens_soas_aer_gc
vol_h2o = mass_water / dens_h2o
if (isize .eq. 1) then
vol_dst = mass_dst / dens_dst1_aer_gc
else
vol_dst = mass_dst / dens_dst2_aer_gc
end if
! 7/23/09 SAM calculate vol_h2o from kappas in Petters and Kreidenweis ACP, 2007, vol. 7, 1961-1971.
! Their kappas are the hygroscopicities used in Abdul-Razzak and Ghan, 2004, JGR, V105, p. 6837-6844.
! These kappas are defined in module_data_sorgam and module_data_mosaic_asect.
! Note that hygroscopicities are at 298K and specific surface tension - further refinement could
! include temperature dependence in Petters and Kreidenweis
! Also, for hygroscopic BC part, assume kappa of OC (how can BC be hydrophylic?)
mass_dry_a = mass_so4 + mass_nit + mass_nh4 + &
mass_oc + mass_bc + mass_seas + &
mass_dst + mass_soas
mass_wet_a = mass_dry_a + mass_water
vol_dry_a = vol_so4 + vol_nit + vol_nh4 + &
vol_oc + vol_bc + vol_seas + &
vol_dst + vol_soas
vol_wet_a = vol_dry_a + vol_h2o
vol_shell = vol_wet_a - vol_bc
! Aerosol number concentrations in each size bin
num_a = sixpi * vol_wet_a / (xdia_cm(isize) * xdia_cm(isize) * xdia_cm(isize))
num_a_lo = sixpi * vol_dry_a / (dlo_sectm(isize) * dlo_sectm(isize) * dlo_sectm(isize) * 1.0E6)
num_a_hi = sixpi * vol_dry_a / (dhi_sectm(isize) * dhi_sectm(isize) * dhi_sectm(isize) * 1.0E6)
if ( vol_dry_a .le. 1.0E-15 .and. num_a .le. 1.0E-10 ) then
if ( num_a .gt. num_a_lo ) then
num_a = num_a_lo
else if (num_a .lt. num_a_hi) then
num_a = num_a_hi
end if
dp_dry_a = dhi_sectm(isize)*1.0E2
dp_wet_a = dhi_sectm(isize)*1.0E2
dp_bc_a = dhi_sectm(isize)*1.0E2
else
if ( num_a .gt. num_a_lo ) then
num_a = num_a_lo
else if ( num_a .lt. num_a_hi ) then
num_a = num_a_hi
end if
dp_dry_a = (sixpi * vol_dry_a / num_a )**0.3333333
dp_wet_a = (sixpi * vol_wet_a / num_a )**0.3333333
dp_bc_a = (sixpi * vol_bc / num_a )**0.3333333
end if
!shortwave
do ns=1,nswbands
ri_dum = (0.0,0.0)
ri_dum = (swref_index_nh4so4(ns) * vol_so4) + &
(ref_index_nh4no3 * vol_nit) + &
(ref_index_nh4no3 * vol_nh4) + &
(swref_index_dust(ns) * vol_dst) + &
(swref_index_oc(ns) * vol_oc) + &
(ref_index_bc * vol_bc) + &
(swref_index_nacl(ns) * vol_seas) + &
(swref_index_oc(ns) * vol_soas) + &
(swref_index_h2o(ns) * vol_h2o)
ri_ave_a = ri_dum/vol_wet_a
ri_dum = (swref_index_nh4so4(ns) * vol_so4) + &
(ref_index_nh4no3 * vol_nit) + &
(ref_index_nh4no3 * vol_nh4) + &
(swref_index_dust(ns) * vol_dst) + &
(swref_index_oc(ns) * vol_oc) + &
(swref_index_nacl(ns) * vol_seas) + &
(swref_index_oc(ns) * vol_soas) + &
(swref_index_h2o(ns) * vol_h2o)
if(dp_wet_a/2.0 .lt. dlo_sectm(isize)*1.0E2/2.0) then
swrefindx(i, k, j, isize, ns) = (1.5,0.0)
radius_wet(i, k, j, isize) = dlo_sectm(isize)*1.0E2/2.0
number_bin(i, k, j, isize) = num_a
radius_core(i, k, j, isize) = 0.0
swrefindx_core(i, k, j, isize, ns) = ref_index_bc
swrefindx_shell(i, k, j, isize, ns) = ref_index_oin
else if (vol_shell .lt. 1.0E-20) then
swrefindx(i, k, j, isize, ns) = (1.5,0)
radius_wet(i, k, j, isize) = dlo_sectm(isize)*1.0E2/2.0
number_bin(i, k, j, isize) = num_a
radius_core(i, k, j, isize) = 0.0
swrefindx_core(i, k, j, isize, ns) = ref_index_bc
swrefindx_shell(i, k, j, isize, ns) = ref_index_oin
else
swrefindx(i, k, j, isize, ns) = ri_ave_a
radius_wet(i, k, j, isize) = dp_wet_a/2.0
number_bin(i, k, j, isize) = num_a
radius_core(i, k, j, isize) = dp_bc_a/2.0
swrefindx_core(i, k, j, isize, ns) = ref_index_bc
swrefindx_shell(i, k, j, isize, ns) = ri_dum/vol_shell
endif
enddo ! ns shortwave
!longwave
do ns=1,nlwbands
ri_dum = (0.0,0.0)
ri_dum = (lwref_index_nh4so4(ns) * vol_so4) + &
(ref_index_nh4no3 * vol_nit) + &
(ref_index_nh4no3 * vol_nh4) + &
(lwref_index_dust(ns) * vol_dst) + &
(lwref_index_oc(ns) * vol_oc) + &
(ref_index_bc * vol_bc ) + &
(lwref_index_nacl(ns) * vol_seas) + &
(lwref_index_oc(ns) * vol_soas ) + &
(lwref_index_h2o(ns) * vol_h2o)
ri_ave_a = ri_dum/vol_wet_a
ri_dum = (lwref_index_nh4so4(ns) * vol_so4) + &
(ref_index_nh4no3 * vol_nit) + &
(ref_index_nh4no3 * vol_nh4) + &
(lwref_index_dust(ns) * vol_dst) + &
(lwref_index_oc(ns) * vol_oc) + &
(lwref_index_nacl(ns) * vol_seas) + &
(lwref_index_oc(ns) * vol_soas ) + &
(lwref_index_h2o(ns) * vol_h2o)
if(dp_wet_a/2.0 .lt. dlo_sectm(isize)*1.0E2/2.0) then
lwrefindx(i, k, j, isize, ns) = (1.5,0.0)
radius_wet(i, k, j, isize) = dlo_sectm(isize)*1.0E2/2.0
number_bin(i, k, j, isize) = num_a
radius_core(i, k, j, isize) = 0.0
lwrefindx_core(i, k, j, isize, ns) = ref_index_bc
lwrefindx_shell(i, k, j, isize, ns) = ref_index_oin
elseif(vol_shell .lt. 1.0E-20) then
lwrefindx(i, k, j, isize, ns) = (1.5,0.0)
radius_wet(i, k, j, isize) = dlo_sectm(isize)*1.0E2/2.0
number_bin(i, k, j, isize) = num_a
radius_core(i, k, j, isize) = 0.0
lwrefindx_core(i, k, j, isize, ns) = ref_index_bc
lwrefindx_shell(i, k, j, isize, ns) = ref_index_oin
else
lwrefindx(i, k, j, isize, ns) = ri_ave_a
radius_wet(i, k, j, isize) = dp_wet_a/2.0
number_bin(i, k, j, isize) = num_a
radius_core(i, k, j, isize) = dp_bc_a/2.0
lwrefindx_core(i, k, j, isize, ns) = ref_index_bc
lwrefindx_shell(i, k, j, isize, ns) = ri_dum/vol_shell
end if
enddo ! ns longwave
enddo !isize
enddo !i
enddo !j
enddo !k
return
end subroutine optical_prep_sectional_wrfgc
!below is the detail calculation for MIE code
!czhao
!
!***********************************************************************
! <1.> subr mieaer