/
AERO_DATA.F
814 lines (680 loc) · 34.9 KB
/
AERO_DATA.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
!------------------------------------------------------------------------!
! The Community Multiscale Air Quality (CMAQ) system software is in !
! continuous development by various groups and is based on information !
! from these groups: Federal Government employees, contractors working !
! within a United States Government contract, and non-Federal sources !
! including research institutions. These groups give the Government !
! permission to use, prepare derivative works of, and distribute copies !
! of their work in the CMAQ system to the public and to permit others !
! to do so. The United States Environmental Protection Agency !
! therefore grants similar permission to use the CMAQ system software, !
! but users are requested to provide copies of derivative works or !
! products designed to operate in the CMAQ system to the United States !
! Government without restrictions as to use by others. Software !
! that is used with the CMAQ system but distributed under the GNU !
! General Public License or the GNU Lesser General Public License is !
! subject to their copyright restrictions. !
!------------------------------------------------------------------------!
C RCS file, release, date & time of last delta, author, state, [and locker]
C $Header: /project/yoj/arc/CCTM/src/aero/aero6/AERO_DATA.F,v 1.10 2012/01/19 13:18:37 yoj Exp $
C:::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
Module aero_data
C Defines aerosol species arrays and the parameters required in aerosol
C processing.
C Contains:
C Subroutine map_aero
C Subroutine map_pmemis
C Subroutine extract_aero
C Subroutine update_aero
C Function findAero
C Revision History:
C First version was coded in April 2010 by Steve Howard with
C Prakash Bhave, Jeff Young, and Sergey Napelenok.
C HS = Heather Simon; GS = Golam Sarwar; SH = Steve Howard; SR = Shawn Roselle;
C JY = Jeff Young; HOTP = Havala Pye; KF = Kathleen Fahey
C HS 01/24/11 Updated to initial AERO6; PMOTHR -> 9 species
C GS 03/02/11 Revised parameters for Mg, K, and Ca
C SH 03/10/11 Inserted functionality of PMEM_DEFN
C -new subroutine map_pmemis
C HS 03/10/11 Made PNCOM a required species
C SR 03/25/11 Replaced I/O API include files with UTILIO_DEFN
C SH 04/04/11 Added sea-salt speciation factors
C GS 04/09/11 Updated sea-salt speciation; replaced ANAK with ASEACAT;
C made MG, K, CA, SEACAT required species;
C JY 04/21/11 Added optional log messages (verbose_aero)
C JY 05/02/11 Added Reshape to aerospc_ssf initialization for pgf90 compiler
C JY 06/08/12 remove full character blank padding put in for GNU Fortran (GCC) 4.1.2
C HOTP 01/16/13 Removed isoprene acid enhanced aerosol and added new isoprene aerosol
C HOTP 01/18/13 Added information for particle phase reaction of isoprene products
C based on Eddingsaas et al. 2010
C KF 09/17/14 Changed emitted modal PM mass fractions according to
C Elleman and Covert (2010)
C HOTP 09/27/14 Added alkane and PAH SOA species. PAH SOA densities follow
C Chan et al.
C JY 08/24/15 Changed visual index factors
C
C References:
C 1. Eddingsaas, N. C., Vandervelde, D. G., Wennberg, P. O., "Kinetics and
C products of the acid-catalyzed ring-opening of atmospherically relevant
C butyl epoxy alcohols," J. Phys. Chem. A., 2010, 114, 8106-8113.
C 2. Chan et al., Secondary organic aerosol formation from photooxidation
C of naphthalene and alkylnaphthalenes: implications for oxidation of
C intermediate volatility organic compounds (IVOCs), Atmos. Chem. Phys.,
C 2009, 9, 3049-3060.
C JY: From Christian Hogrefe...
C Based on Malm and Hand (Atmos. Env. 41, 3407-3427, 2007), the revised
C IMPROVE extinction calculation includes coarse particles, sea salt, and
C a relative humidity correction for sea salt. Also, the factor for "LAC"
C (light absorbing carbon, i.e. AECI and AECJ) should be 10, not 0 since
C both scattering and absorption contribute to total extinction.
C ASEACAT includes all sea-salt cations in coarse mode (Na, Ca, K, and Mg)
C Also note...
C In the Fortran user-derived spcs_type, below, visual_idx is an optimal dry mass
C extinction efficiency [m^2/g], see White, Atmos.Env.,294(10)(1990),pp 2673-1672
C and Malm, et al., JGR,99(D1)(1994),pp 1347-1370
C----------------------------------------------------------------------
Implicit None
C Number of aerosol species and modes
Integer, Parameter :: n_aerospc = 45 ! number of aero species
Integer, Parameter :: n_mode = 3 ! number of modes:
! 1 = Aitken
! 2 = accumulation
! 3 = coarse
C Default minimum concentration
Real, Parameter :: conmin = 1.0E-30 ! [ ug/m^3 ]
Real, Parameter :: cm_set( n_mode ) = (/conmin, conmin, conmin/)
Real, Parameter :: cm_so4( n_mode ) = (/1.0E-12,1.0E-6, conmin/)
Real, Parameter :: cm_cor( n_mode ) = (/conmin, conmin, 1.889544E-05/)
C Emissions splits
Real, Parameter :: es_fin( n_mode ) = (/ 0.100, 0.900, 0.000 /)
Real, Parameter :: es_acc( n_mode ) = (/ 0.000, 1.000, 0.000 /)
Real, Parameter :: es_cor( n_mode ) = (/ 0.000, 0.000, 1.000 /)
Real, Parameter :: es_0 ( n_mode ) = (/ 0.000, 0.000, 0.000 /)
C Sea salt factors
! Real, Parameter :: ss_so4( n_mode ) = (/ 0.0, 0.0776, 0.0776 /)
! Real, Parameter :: ss_cl ( n_mode ) = (/ 0.0, 0.5338, 0.5338 /)
! Real, Parameter :: ss_na ( n_mode ) = (/ 0.0, 0.3086, 0.0 /)
! Real, Parameter :: ss_mg ( n_mode ) = (/ 0.0, 0.0368, 0.0 /)
! Real, Parameter :: ss_k ( n_mode ) = (/ 0.0, 0.0114, 0.0 /)
! Real, Parameter :: ss_ca ( n_mode ) = (/ 0.0, 0.0118, 0.0 /)
! Real, Parameter :: ss_sea( n_mode ) = (/ 0.0, 0.0 , 0.3686 /)
! Real, Parameter :: ss_0 ( n_mode ) = (/ 0.0, 0.0 , 0.0 /)
C Flag to obtain coagulation coefficients
C by analytical approximation (True) or by Gauss-Hermite quadrature (False)
Logical, Parameter :: fastcoag_flag = .True.
C Define Logical values as T and F for the aerospc table
Logical, Parameter, Private :: T = .true.
Logical, Parameter, Private :: F = .false.
Integer, Private, Save :: logdev
Integer, Private, External :: setup_logdev
C-------------------------------------------------------------------------------------------------------
Type spcs_type
Character( 16 ) :: name( n_mode ) ! names of aerosol species for each mode
Real :: min_conc( n_mode ) ! minimum concentration values for each mode
Real :: density ! density [ kg/m^3 ]
Logical :: no_M2Wet ! flag to exclude from 2nd moment during transport
Logical :: nonVol_soa ! non-volatile SOA flag
Logical :: tracer ! tracer flag; does have not mass
Integer :: charge ! electroneutrality charge
Real :: visual_idx ! visual index factor [ m^2/g ]
Character( 16 ) :: optic_surr ! optical surrogate name
Character( 16 ) :: emis ! file emissions names
Real :: emis_split( n_mode ) ! minimum concentration values for each mode
End Type spcs_type
Type( spcs_type ), Parameter :: aerospc( n_aerospc ) = (/
C nonVolSOA
C | Charge
C -----------Name-------------- NoM2Wet |Tracer| Emissions
C Aitken Accum Coarse Min_Concs Density | | | | Visidx OptSurr Emis Splits
C --------- --------- --------- --------- ------- + + + + ------ -------- -------- ------
& spcs_type((/'ASO4I ','ASO4J ','ASO4K '/), cm_so4, 1800.0, F, F, F, -2, 3.0, 'SOLUTE','PSO4 ',es_fin),
& spcs_type((/'ANO3I ','ANO3J ','ANO3K '/), cm_set, 1800.0, F, F, F, -1, 3.0, 'SOLUTE','PNO3 ',es_fin),
& spcs_type((/'ACLI ','ACLJ ','ACLK '/), cm_set, 2200.0, F, F, F, -1, 1.37,'SOLUTE','PCL ',es_fin),
& spcs_type((/'ANH4I ','ANH4J ','ANH4K '/), cm_set, 1800.0, F, F, F, 1, 3.0, 'SOLUTE','PNH4 ',es_fin),
& spcs_type((/'ANAI ','ANAJ ',' '/), cm_set, 2200.0, F, F, F, 1, 1.37,'SOLUTE','PNA ',es_fin),
& spcs_type((/' ','AMGJ ',' '/), cm_set, 2200.0, F, F, F, 2, 1.0, 'DUST ','PMG ',es_acc),
& spcs_type((/' ','AKJ ',' '/), cm_set, 2200.0, F, F, F, 1, 1.0, 'DUST ','PK ',es_acc),
& spcs_type((/' ','ACAJ ',' '/), cm_set, 2200.0, F, F, F, 2, 1.0, 'DUST ','PCA ',es_acc),
& spcs_type((/'APOCI ','APOCJ ',' '/), cm_set, 2000.0, F, F, F, 0, 4.0, 'DUST ','POC ',es_fin),
& spcs_type((/'APNCOMI','APNCOMJ',' '/), cm_set, 2000.0, F, F, F, 0, 4.0, 'DUST ','PNCOM ',es_fin),
& spcs_type((/'AECI ','AECJ ',' '/), cm_set, 2200.0, F, F, F, 0, 10.0, 'SOOT ','PEC ',es_fin),
& spcs_type((/' ','AFEJ ',' '/), cm_set, 2200.0, F, F, F, 0, 1.0, 'DUST ','PFE ',es_acc),
& spcs_type((/' ','AALJ ',' '/), cm_set, 2200.0, F, F, F, 0, 1.0, 'DUST ','PAL ',es_acc),
& spcs_type((/' ','ASIJ ',' '/), cm_set, 2200.0, F, F, F, 0, 1.0, 'DUST ','PSI ',es_acc),
& spcs_type((/' ','ATIJ ',' '/), cm_set, 2200.0, F, F, F, 0, 1.0, 'DUST ','PTI ',es_acc),
& spcs_type((/' ','AMNJ ',' '/), cm_set, 2200.0, F, F, F, 0, 1.0, 'DUST ','PMN ',es_acc),
& spcs_type((/'AH2OI ','AH2OJ ','AH2OK '/), cm_set, 1000.0, T, F, F, 0, 0.0, 'WATER ','PH2O ',es_fin),
& spcs_type((/'AH3OPI ','AH3OPJ ','AH3OPK '/), cm_set, 1000.0, T, F, T, 0, 0.0, 'WATER ',' ',es_0),
& spcs_type((/'AOTHRI ','AOTHRJ ',' '/), cm_set, 2200.0, F, F, F, 0, 1.0, 'DUST ','PMOTHR',es_fin),
& spcs_type((/' ','AALK1J ',' '/), cm_set, 2000.0, T, F, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ','AALK2J ',' '/), cm_set, 2000.0, T, F, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ','AXYL1J ',' '/), cm_set, 2000.0, T, F, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ','AXYL2J ',' '/), cm_set, 2000.0, T, F, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ','AXYL3J ',' '/), cm_set, 2000.0, F, T, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ','ATOL1J ',' '/), cm_set, 2000.0, T, F, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ','ATOL2J ',' '/), cm_set, 2000.0, T, F, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ','ATOL3J ',' '/), cm_set, 2000.0, F, T, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ','ABNZ1J ',' '/), cm_set, 2000.0, T, F, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ','ABNZ2J ',' '/), cm_set, 2000.0, T, F, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ','ABNZ3J ',' '/), cm_set, 2000.0, F, T, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ','ATRP1J ',' '/), cm_set, 2000.0, T, F, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ','ATRP2J ',' '/), cm_set, 2000.0, T, F, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ','AISO1J ',' '/), cm_set, 2000.0, T, F, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ','AISO2J ',' '/), cm_set, 2000.0, T, F, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ','AISO3J ',' '/), cm_set, 2000.0, F, T, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ','ASQTJ ',' '/), cm_set, 2000.0, T, F, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ','APAH1J ',' '/), cm_set, 1480.0, T, F, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ','APAH2J ',' '/), cm_set, 1480.0, T, F, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ','APAH3J ',' '/), cm_set, 1550.0, F, T, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ','AOLGAJ ',' '/), cm_set, 2000.0, F, T, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ','AOLGBJ ',' '/), cm_set, 2000.0, F, T, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ','AORGCJ ',' '/), cm_set, 2000.0, F, T, F, 0, 4.0, 'DUST ',' ',es_0),
& spcs_type((/' ',' ','ASOIL '/), cm_set, 2600.0, F, F, F, 0, 0.6, 'DUST ',' ',es_cor),
& spcs_type((/' ',' ','ACORS '/), cm_cor, 2200.0, F, F, F, 0, 0.6, 'DUST ','PMC ',es_cor),
& spcs_type((/' ',' ','ASEACAT'/), cm_cor, 2200.0, F, F, F, 1, 1.37,'SOLUTE',' ',es_cor)/)
C Sea-Salt Speciation factors
Real, Parameter :: aerospc_ssf( n_mode, n_aerospc ) = Reshape (
C Aitken Accum Coarse
C ------- ------- -------
& (/ 0.0, 0.0776, 0.0776, ! ASO4
& 0.0, 0.0, 0.0, ! ANO3
& 0.0, 0.5538, 0.5538, ! ACL
& 0.0, 0.0, 0.0, ! ANH4
& 0.0, 0.3086, 0.0, ! ANA
& 0.0, 0.0368, 0.0, ! AMG
& 0.0, 0.0114, 0.0, ! AK
& 0.0, 0.0118, 0.0, ! ACA
& 0.0, 0.0, 0.0, ! APOC
& 0.0, 0.0, 0.0, ! APNCOM
& 0.0, 0.0, 0.0, ! AEC
& 0.0, 0.0, 0.0, ! AFE
& 0.0, 0.0, 0.0, ! AAL
& 0.0, 0.0, 0.0, ! ASI
& 0.0, 0.0, 0.0, ! ATI
& 0.0, 0.0, 0.0, ! AMN
& 0.0, 0.0, 0.0, ! AH2O
& 0.0, 0.0, 0.0, ! AH3OP
& 0.0, 0.0, 0.0, ! AOTHR
& 0.0, 0.0, 0.0, ! AALK1
& 0.0, 0.0, 0.0, ! AALK2
& 0.0, 0.0, 0.0, ! AXYL1
& 0.0, 0.0, 0.0, ! AXYL2
& 0.0, 0.0, 0.0, ! AXYL3
& 0.0, 0.0, 0.0, ! ATOL1
& 0.0, 0.0, 0.0, ! ATOL2
& 0.0, 0.0, 0.0, ! ATOL3
& 0.0, 0.0, 0.0, ! ABNZ1
& 0.0, 0.0, 0.0, ! ABNZ2
& 0.0, 0.0, 0.0, ! ABNZ3
& 0.0, 0.0, 0.0, ! ATRP1
& 0.0, 0.0, 0.0, ! ATRP2
& 0.0, 0.0, 0.0, ! AISO1
& 0.0, 0.0, 0.0, ! AISO2
& 0.0, 0.0, 0.0, ! AISO3
& 0.0, 0.0, 0.0, ! ASQT
& 0.0, 0.0, 0.0, ! APAH1
& 0.0, 0.0, 0.0, ! APAH2
& 0.0, 0.0, 0.0, ! APAH3
& 0.0, 0.0, 0.0, ! AOLGA
& 0.0, 0.0, 0.0, ! AOLGB
& 0.0, 0.0, 0.0, ! AORGC
& 0.0, 0.0, 0.0, ! ASOIL
& 0.0, 0.0, 0.0, ! ACORS
& 0.0, 0.0, 0.3686 /), ! ASEACAT
& (/ n_mode,n_aerospc /), Order = (/ 1,2 /) )
C Required species
Character( 16 ), Private, Parameter :: req_so4 = 'ASO4'
Character( 16 ), Private, Parameter :: req_no3 = 'ANO3'
Character( 16 ), Private, Parameter :: req_cl = 'ACL'
Character( 16 ), Private, Parameter :: req_nh4 = 'ANH4'
Character( 16 ), Private, Parameter :: req_na = 'ANA'
Character( 16 ), Private, Parameter :: req_mg = 'AMG'
Character( 16 ), Private, Parameter :: req_k = 'AK'
Character( 16 ), Private, Parameter :: req_ca = 'ACA'
Character( 16 ), Private, Parameter :: req_poc = 'APOC'
Character( 16 ), Private, Parameter :: req_h2o = 'AH2O'
Character( 16 ), Private, Parameter :: req_h3op = 'AH3OP'
Character( 16 ), Private, Parameter :: req_soil = 'ASOIL'
Character( 16 ), Private, Parameter :: req_cors = 'ACORS'
Character( 16 ), Private, Parameter :: req_ncom = 'APNCOM'
Character( 16 ), Private, Parameter :: req_seacat = 'ASEACAT'
C Optional Species
Character( 16 ), Private, Parameter :: req_phgj = 'APHGJ'
C Indices of required species
Integer :: aso4_idx
Integer :: ano3_idx
Integer :: acl_idx
Integer :: anh4_idx
Integer :: ana_idx
Integer :: amg_idx
Integer :: ak_idx
Integer :: aca_idx
Integer :: apoc_idx
Integer :: ah2o_idx
Integer :: ah3op_idx
Integer :: asoil_idx
Integer :: acors_idx
Integer :: apncom_idx
Integer :: aseacat_idx
C Indices of Optional species
Integer :: aphgj_idx
Real :: aerospc_mw( n_aerospc ) ! aero species M.W. (from AE_SPC.EXT) [ g/mol ]
Real :: aerospc_conc( n_aerospc,n_mode ) ! aero species concentration [ ug/m^3 ]
C Common factors
Real( 8 ) :: h2ofac ! converts mass concentrations [ug/m3] to 3rd moment concentrations [m3/m3]
C-------------------------------------------------------------------------------------------------------
Type mode_type
Character( 16 ) :: num_name ! name of aerosol number variable
Character( 16 ) :: srf_name ! name of aerosol surface area variable
Real :: min_numconc ! minimum number concentration
Real :: min_m2conc ! minimum 2nd moment concentration
Real :: min_m3conc ! minimum 3rd moment concentration
End Type mode_type
Type ( mode_type ), Parameter :: aeromode( n_mode ) = (/
C number surface minimum minimum minimum
C name name numconc m2conc m3conc
C ---------- ------- -------- ------- ------
& mode_type('NUMATKN', 'SRFATKN', conmin, conmin, conmin),
& mode_type('NUMACC ', 'SRFACC ', conmin, conmin, conmin),
& mode_type('NUMCOR ', 'SRFCOR ', conmin, conmin, conmin)/)
Real, Parameter :: min_sigma_g = 1.05
Real, Parameter :: max_sigma_g = 2.50
Real, Parameter :: def_sigma_g( n_mode ) = (/ 1.70, 2.0, 2.2 /) ! default background sigma-g for each mode
Real, Parameter :: def_diam( n_mode ) = (/ 1.0E-8, 7.0E-8, 1.0E-6 /) ! default background mean diameter for each mode
Real :: moment0_conc( n_mode ) ! 0th moment concentration
Real :: moment2_conc( n_mode ) ! 2nd moment concentration
Real :: moment3_conc( n_mode ) ! 3rd moment concentration
C Mass concentration (calculated by GETPAR)
Real :: aeromode_mass( n_mode ) ! [ ug/m^3 ]
C Particle density (calculated by GETPAR)
Real :: aeromode_dens( n_mode ) ! [ kg/m^3 ]
C Geometric mean diameter (calculated by GETPAR)
Real :: aeromode_diam( n_mode ) ! [ m ]
C Log of geometric standard deviation (calculated by GETPAR )
Real :: aeromode_sdev( n_mode )
C Minimum number (calculated in map_aero routine)
Real :: aeromode_minNum( n_mode )
C Minimum 2nd moment (calculated in map_aero routine)
Real :: aeromode_minM2( n_mode )
C Mapping for loading from and unloading to CGRID array
Integer :: aerospc_map( n_aerospc,n_mode ) ! indices of aero species to CGRID
Integer :: aeronum_map( n_mode ) ! indices of aero number variable to CGRID
Integer :: aerosrf_map( n_mode ) ! indices of aero surf area variable to CGRID
C Missing aerosol species map
Logical :: aero_missing( n_aerospc,n_mode ) ! indices of aero species to CGRID
C Emissions mapping
Integer :: n_emis_pm ! number of aerospc with emissions
Integer :: pmem_map( n_aerospc ) ! mapping to aerospc array for PM emissions
Character( 16 ) :: pmem_units ! units for PM emissions for all species
C Private variables for loading from and unloading to CGRID array
Logical, Private, Save :: mapped = .False.
Character( 16 ), Private, Save :: pname = 'Aero_Data'
Contains
C-----------------------------------------------------------------------
Subroutine map_aero()
C Defines aerosol mapping from CGRID for species concentration and moments.
C Revision History:
C First version was coded in April 2010 by Steve Howard with
C Prakash Bhave, Jeff Young, and Sergey Napelenok.
C HS 01/24/11 Renamed AORGPA as POC for AERO6
C GS 03/02/11 Find new req`d species for AERO6 (Mg, K, Ca)
C HS 03/10/11 Get index for new required species, PNCOM
C-----------------------------------------------------------------------
Use rxns_data ! chemical mechanism data
Use cgrid_spcs ! CGRID mechanism species
Use aeromet_data
Use utilio_defn
Implicit None
C Local Variables:
Character( 80 ) :: xmsg
Integer m, n, spc
Real so4fac
Real anthfac
logdev = setup_logdev()
if ( mapped ) Return
aerospc_mw = 0.0
aerospc_map = 0
C build mapping to CGRID for each aero spc
Do m = 1, n_mode
Do spc = 1, n_aerospc
! write( logdev,* ) '==1 ', aerospc( spc )%name( m ), ' ', ae_spc( 1 ), ' ', ae_spc( 15 )
! write( logdev,* ) '==2 ', trim( aerospc( spc )%name( m ) ), ' ', ae_spc( 1 ), ' ', ae_spc( 15 )
If ( aerospc( spc )%name( m ) .Ne. ' ' ) Then
aero_missing( spc,m ) = .False.
n = index1( aerospc( spc )%name( m ), n_ae_spc, ae_spc )
If ( n .Eq. 0 ) Then
xmsg = 'Species '// aerospc(spc)%name( m )
& // ' in aerospc name is not in the AE_SPC table.'
Call m3exit( pname, 0, 0, xmsg, xstat3 )
End If
aerospc_map( spc,m ) = ae_strt - 1 + n
If ( aerospc_mw( spc ) .Lt. 0.5 ) Then ! mw=0 means a new species
aerospc_mw( spc ) = ae_molwt( n )
Else If ( aerospc_mw( spc ) .Ne. ae_molwt( n ) ) Then
xmsg = 'molecular weight of ' // Trim( aerospc( spc )%name( m ) )
& // ' is different from that of the same species'
& // ' in the same or another mode.'
Call m3exit( pname, 0, 0, xmsg, xstat3 )
End If
Else
aero_missing( spc,m ) = .True.
End If
End Do
End Do
C Build mapping to CGRID for aero # and surf area variables
aeronum_map = 0
aerosrf_map = 0
Do m = 1, n_mode
n = index1( aeromode( m )%num_name, n_ae_spc, ae_spc )
If ( n .Eq. 0 ) Then
xmsg = 'Species ' // Trim( aeromode( m )%num_name )
& //' in aeronum name is not in AE_SPC'
Call m3exit( pname, 0, 0, xmsg, xstat3 )
Else
aeronum_map( m ) = ae_strt - 1 + n
End If
n = index1( aeromode( m )%srf_name, n_ae_spc, ae_spc )
If ( n .Eq. 0 ) Then
xmsg = 'species ' // Trim( aeromode( m )%srf_name )
& // ' in aerosrf name is not in AE_SPC'
Call m3exit( pname, 0, 0, xmsg, xstat3 )
Else
aerosrf_map( m ) = ae_strt - 1 + n
End If
End Do
C Find indices of required species
aso4_idx = findAero( req_so4, .True. )
ano3_idx = findAero( req_no3, .True. )
acl_idx = findAero( req_cl, .True. )
anh4_idx = findAero( req_nh4, .True. )
ana_idx = findAero( req_na, .True. )
amg_idx = findAero( req_mg, .True. )
ak_idx = findAero( req_k, .True. )
aca_idx = findAero( req_ca, .True. )
apoc_idx = findAero( req_poc, .True. )
ah2o_idx = findAero( req_h2o, .True. )
ah3op_idx = findAero( req_h3op, .True. )
asoil_idx = findAero( req_soil, .True. )
acors_idx = findAero( req_cors, .True. )
apncom_idx = findAero( req_ncom, .True. )
aseacat_idx = findAero( req_seacat, .True. )
aphgj_idx = findAero( req_phgj, .False. )
#ifdef verbose_aero
Write( logdev,'( /5x, a )' ) 'map_aero required species'
Write( logdev,'( 5x, a, i4 )' ) 'aso4_idx: ', aso4_idx
Write( logdev,'( 5x, a, i4 )' ) 'ano3_idx: ', ano3_idx
Write( logdev,'( 5x, a, i4 )' ) 'acl_idx: ', acl_idx
Write( logdev,'( 5x, a, i4 )' ) 'anh4_idx: ', anh4_idx
Write( logdev,'( 5x, a, i4 )' ) 'ana_idx: ', ana_idx
Write( logdev,'( 5x, a, i4 )' ) 'amg_idx: ', amg_idx
Write( logdev,'( 5x, a, i4 )' ) 'ak_idx: ', ak_idx
Write( logdev,'( 5x, a, i4 )' ) 'aca_idx: ', aca_idx
Write( logdev,'( 5x, a, i4 )' ) 'apoc_idx: ', apoc_idx
Write( logdev,'( 5x, a, i4 )' ) 'ah2o_idx: ', ah2o_idx
Write( logdev,'( 5x, a, i4 )' ) 'ah3op_idx: ', ah3op_idx
Write( logdev,'( 5x, a, i4 )' ) 'asoil_idx: ', asoil_idx
Write( logdev,'( 5x, a, i4 )' ) 'acors_idx: ', acors_idx
Write( logdev,'( 5x, a, i4 )' ) 'apncom_idx: ', apncom_idx
Write( logdev,'( 5x, a, i4 )' ) 'aseacat_idx:', aseacat_idx
If ( aphgj_idx .Gt. 0 ) Then
Write( logdev,'( 5x, a, i4 )' ) 'aphgj_idx: ', aphgj_idx
else
Write( logdev,'( 5x, a, i4 )' ) 'PHGJ not found so aphgj_idx: ', aphgj_idx
End If
#endif
C Compute common factors
h2ofac = 1.0D-9 * f6dpi / Real( aerospc( ah2o_idx )%density, 8 )
C compute aeromode_minNum and aeromode_minM2
so4fac = 1.0E-9 * Real( f6dpi, 4 ) / aerospc( aso4_idx )%density
anthfac = 1.0E-9 * Real( f6dpi, 4 ) / aerospc( acors_idx )%density
Do m = 1, n_mode
If ( m .lt. n_mode ) Then
aeromode_minNum( m ) = so4fac * aerospc( aso4_idx )%min_conc( m ) /
& ( def_diam( m )**3 * Exp( 4.5 * Log( def_sigma_g( m ) )**2 ) )
Else
aeromode_minNum( m ) = anthfac * aerospc( acors_idx )%min_conc( m ) /
& ( def_diam( m )**3 * Exp( 4.5 * Log( def_sigma_g( m ) )**2 ) )
End If
aeromode_minM2( m ) = aeromode_minNum( m ) * def_diam( m )**2 *
& Exp( 2.0 * Log( def_sigma_g( m ) )**2 )
End do
mapped = .True.
Write( logdev,'( 5x, a )' ) ' --- Aero Species Mapped ---'
Return
End Subroutine map_aero
C-----------------------------------------------------------------------
Subroutine map_pmemis ( )
C
C Set the emissions units from the header of EMIS_1 file
C Verify that all the units on the file are consistent
Use utilio_defn
Implicit None
Include SUBST_FILES_ID ! file name parameters
C Parameters:
Character( 10 ), Parameter :: blank10 = ' '
C Local Variables:
Character( 16 ), Save :: pname = 'map_pmemis'
Character( 512 ) :: xmsg1
Character( 1024 ) :: xmsg2
Character( 10 ) :: units
Integer :: indx
Integer :: v
Logical :: found, match
Logical, Save :: pm_mapped = .False.
logdev = setup_logdev()
C Create mapping only if first call
If ( pm_mapped ) Return
C Call routine to map aerosol species array
If ( .Not. mapped ) Call map_aero ( )
C Open the gridded emissions file, which contains gas, aerosol, and non-reactive
C species
If ( .Not. open3( emis_1, fsread3, pname ) ) Then
xmsg1 = 'Could not open '// emis_1 // ' file'
Call m3exit( pname, 0, 0, Trim( xmsg1 ), xstat1 )
End If
If ( .Not. desc3( emis_1 ) ) Then
xmsg1 = 'Could not get '// 'EMIS_1' // ' file description'
Call m3exit( pname, 0, 0, Trim( xmsg1 ), xstat2 )
End If
C Search emissions file for emission species names. Verify that their units
C are the same and set pmem_units
n_emis_pm = 0
pmem_units = 'null'
found = .True.
match = .True.
xmsg1 = 'Could not find the following species in emissions file'
xmsg2 = 'PM Units not uniform in EMIS_1 file.'
Do v = 1, n_aerospc
If ( aerospc( v )%emis .Ne. ' ' ) Then
indx = index1( aerospc( v )%emis, nvars3d, vname3d )
If ( indx .Le. 0 ) Then
xmsg1 = Trim( xmsg1 ) // crlf() // blank10
& // Trim( aerospc( v )%emis )
found = .False.
Cycle
End If
n_emis_pm = n_emis_pm + 1
pmem_map( n_emis_pm ) = v
C Change UNITS to upper case
units = units3d( indx )
Call upcase( units )
C Save units on first emissions
If ( pmem_units .Eq. 'null' ) pmem_units = units
C Check that all emissions units match
If ( pmem_units .Ne. units ) Then
xmsg2 = Trim( xmsg2 ) // crlf() // blank10
& // Trim( aerospc( v )%emis )
& // ' [' // Trim( units3d( indx ) ) // ']'
match = .False.
End If
End If
End Do
If ( .Not. found ) Then
Call m3exit( pname, 0, 0, Trim( xmsg1 ), xstat2 )
End If
If ( .Not. match ) Then
Call m3exit( pname, 0, 0, Trim( xmsg2 ), xstat2 )
End If
#ifdef verbose_aero
Write( logdev,'( /5x, a )' ) 'pmem_map to aerospc'
Do v = 1, n_emis_pm
Write( logdev,'( 5x, a, 2i4, 2x, a )' ) 'pmem_map:', v, pmem_map( v ), aerospc( pmem_map( v ) )%emis
End Do
#endif
pm_mapped = .True.
Write( logdev,'( 5x, a )' ) ' --- PM Emis Species Mapped ---'
Return
End Subroutine map_pmemis
C-----------------------------------------------------------------------
Subroutine extract_aero( conc, minchk )
C Extracts aerosol data into the AERO_DATA:aerospc_conc array
C The original idea is that the data for conc comes from CGRID
C Revision History:
C First version was coded in April 2010 by Steve Howard with
C Prakash Bhave, Jeff Young, and Sergey Napelenok.
C-----------------------------------------------------------------------
Implicit None
C Arguments:
Real, Intent( In ) :: conc( : )
Logical, Intent( In ) :: minchk
C Local Variables:
Integer m, n, spc
If ( .Not. mapped ) Then
Call map_aero()
End If
C Copy grid cell concentrations of aero species to aerospc_conc
aerospc_conc = 0.0
If ( minchk ) Then
Do m = 1, n_mode
Do spc = 1, n_aerospc
n = aerospc_map( spc,m )
If ( n .ne. 0 ) Then
aerospc_conc( spc,m ) = Max( conc( n ), aerospc( spc )%min_conc( m ) ) ! [ug/m^3]
End If
End Do
End Do
Else
Do m = 1, n_mode
Do spc = 1, n_aerospc
n = aerospc_map( spc,m )
If ( n .ne. 0 ) Then
aerospc_conc( spc,m ) = conc( n ) ! [ug/m^3]
End If
End Do
End Do
End If
C Copy grid cell concentrations of aero # and surf area
C Convert and assign to moment0_conc and moment2_conc
moment0_conc = 0.0
moment2_conc = 0.0
If ( minchk ) Then
Do m = 1, n_mode
n = aeronum_map( m )
moment0_conc( m ) = Max( conc( n ), aeromode( m )%min_numconc )
n = aerosrf_map( m )
moment2_conc( m ) = Max( conc( n ), aeromode( m )%min_m2conc )
End Do
Else
Do m = 1, n_mode
n = aeronum_map( m )
moment0_conc( m ) = conc( n )
n = aerosrf_map( m )
moment2_conc( m ) = conc( n )
End Do
End If
Return
End Subroutine extract_aero
C-----------------------------------------------------------------------
Subroutine update_aero( conc, minchk )
C Updates conc from the AERO_DATA:aerospc_conc array.
C The original idea is that the data in conc updates CGRID
C Revision History:
C First version was coded in April 2010 by Steve Howard with
C Prakash Bhave, Jeff Young, and Sergey Napelenok.
C-----------------------------------------------------------------------
Use aeromet_data ! fundamental constants, data type definitions, etc.
Use utilio_defn
Implicit None
C Arguments:
Real, Intent( Out ) :: conc( : )
Logical, Intent( In ) :: minchk
C Local variables:
Character( 80 ) :: xmsg
Integer m, n, spc
If ( .Not. mapped ) Then
xmsg = 'CGRID Species has not been mapped'
Call m3exit( pname, 0, 0, xmsg, xstat3 )
End If
C Copy aerospc_conc back to grid cell concentrations
If ( minchk ) Then
Do m = 1, n_mode
Do spc = 1, n_aerospc
n = aerospc_map( spc,m )
If ( n .Ne. 0 ) Then
conc( n ) = Max( aerospc_conc( spc,m ), aerospc( spc )%min_conc( m ) )
End If
End Do
End Do
Else
Do m = 1, n_mode
Do spc = 1, n_aerospc
n = aerospc_map( spc,m )
If ( n .Ne. 0 ) Then
conc( n ) = aerospc_conc( spc,m )
End If
End Do
End Do
End If
C Copy aero number and surface area back to grid cell concentrations
If ( minchk ) Then
Do m = 1, n_mode
n = aeronum_map( m )
conc( n ) = Max( moment0_conc( m ), aeromode( m )%min_numconc )
End Do
Else
Do m = 1, n_mode
n = aeronum_map( m )
conc( n ) = moment0_conc( m )
End Do
End If
Do m = 1, n_mode
n = aerosrf_map( m )
conc( n ) = Real( pi, 4 ) * moment2_conc( m )
End Do
Return
End Subroutine update_aero
C-----------------------------------------------------------------------
! Integer Function findAero( vname, required ) Result ( idx )
Function findAero( vname, required ) Result ( idx )
C Finds the index of 'required' aerosol species in the aerospc list
C Revision History:
C First version was coded in April 2010 by Steve Howard with
C Prakash Bhave, Jeff Young, and Sergey Napelenok.
C-----------------------------------------------------------------------
Use utilio_defn
Implicit None
C Arguments:
Character( 16 ) :: vname
Logical :: required
Integer :: idx
C Local Variables:
Character( 80 ) :: xmsg
Integer spc, n
C Find the substring vname in aerospc( spc )%name( n )
Do n = 1, n_mode
Do spc = 1, n_aerospc
If ( Index( aerospc( spc )%name( n ), Trim( vname ) ) .Gt. 0 ) Then
idx = spc
Return
End If
End Do
End Do
If ( .Not. required ) Then
idx = -1
Return
End If
xmsg = 'Required Species ' // Trim( vname ) // ' Not found in aerospc names array'
Call m3exit( pname, 0, 0, xmsg, xstat3 )
Return
End Function findAero
End Module aero_data