FV3 Bundle
MieInterface.F90
Go to the documentation of this file.
1 ! $Id: MieInterface.F90,v 1.1.2.1 2015/03/09 15:40:15 drholdaw Exp $
2 
3 
4 #include "MAPL_Generic.h"
5 
7  use esmf
8  use mapl_mod
9  use chem_miemod
10 
11  implicit none
12  private
13 
14  type(chem_mie), save :: mietable
15  logical , save :: mietablecreated=.false.
16 
17  public get_aerooptprop
19  public get_aeroindex
21  public get_aerotaulist
22  public create_aerooptprop
23 
25  module procedure get_aerooptpropallaero3d
26  module procedure get_aerooptpropallaero4d
27  end interface get_aerooptpropallaero
28 
29 contains
30 
31 !=================================================================
32 
33  subroutine get_aerooptprop(idx,ib,qaero,rh,tau,ssa,asy,rc)
34 
35  integer, intent(IN ) :: idx
36  integer, intent(IN ) :: ib
37  real, intent(IN ) :: qaero
38  real, intent(IN ) :: rh
39  real, intent(OUT) :: tau
40  real, intent(OUT) :: ssa
41  real, intent(OUT) :: asy
42  integer, optional, intent(OUT) :: rc
43 
44  integer :: status
45  character(len=ESMF_MAXSTR) :: iam='Get_AeroOptProp'
46 
47  call chem_miequery(mietable,idx,real(ib),qaero,rh, &
48  tau,ssa,asy,rc=status)
49  verify_(status)
50 
51  return_(esmf_success)
52  end subroutine get_aerooptprop
53 
54 !=================================================================
55 
56  subroutine get_aerooptpropallaero3d(na,aerosol,nb,offset,qaero,rh,tau,ssa,asy,rc)
57 
58  integer, intent(in) :: na ! number of aerosols
59  character(len=*), intent(IN ) :: aerosol(:) ! list of aerosols
60  integer, intent(IN ) :: nb
61  integer, intent(IN ) :: offset
62  real, intent(IN ) :: qaero(:,:,:)
63  real, intent(IN ) :: rh(:,:)
64  real, intent(OUT) :: tau(:,:,:)
65  real, intent(OUT) :: ssa(:,:,:)
66  real, intent(OUT) :: asy(:,:,:)
67  integer, optional, intent(OUT) :: rc
68 
69  integer :: status
70  character(len=ESMF_MAXSTR) :: Iam='Get_AeroOptPropAllAero3D'
71 
72  integer :: l, idx
73 
74  real*8 :: tau8(size(tau,1),size(tau,2),size(tau,3))
75  real*8 :: ssa8(size(tau,1),size(tau,2),size(tau,3))
76  real*8 :: asy8(size(tau,1),size(tau,2),size(tau,3))
77 
78  tau8 = 0.0
79  ssa8 = 0.0
80  asy8 = 0.0
81 
82  do l = 1, na
83  idx = chem_miequeryidx( mietable, aerosol(l), status)
84  verify_(status)
85 
86  call chem_miequeryallband3d(mietable,idx,nb,offset,qaero(:,:,l),rh, &
87  tau,ssa,asy,rc=status)
88  verify_(status)
89 
90  tau8 = tau8 + tau
91  ssa8 = ssa8 + (ssa * tau)
92  asy8 = asy8 + asy * (ssa * tau)
93  end do
94 
95  tau = tau8
96  ssa = ssa8
97  asy = asy8
98 
99  return_(esmf_success)
100  end subroutine get_aerooptpropallaero3d
101 
102 !=================================================================
103 
104  subroutine get_aerooptpropallaero4d(na,aerosol,nb,offset,qaero,rh,tau,ssa,asy,rc)
106  integer, intent(IN ) :: na ! number of aerosols
107  character(len=*), intent(IN ) :: aerosol(:) ! list of aerosols
108  integer, intent(IN ) :: nb
109  integer, intent(IN ) :: offset
110  real, intent(IN ) :: qaero(:,:,:,:)
111  real, intent(IN ) :: rh(:,:,:)
112  real, intent(OUT) :: tau(:,:,:,:)
113  real, intent(OUT) :: ssa(:,:,:,:)
114  real, intent(OUT) :: asy(:,:,:,:)
115  integer, optional, intent(OUT) :: rc
116 
117  integer :: status
118  character(len=ESMF_MAXSTR) :: Iam='Get_AeroOptPropAllAero4D'
119 
120  integer :: l, idx
121 
122  real*8 :: tau8(size(tau,1),size(tau,2),size(tau,3),size(tau,4))
123  real*8 :: ssa8(size(tau,1),size(tau,2),size(tau,3),size(tau,4))
124  real*8 :: asy8(size(tau,1),size(tau,2),size(tau,3),size(tau,4))
125 
126  tau8 = 0.0
127  ssa8 = 0.0
128  asy8 = 0.0
129 
130  do l = 1, na
131  idx = chem_miequeryidx( mietable, aerosol(l), status)
132  verify_(status)
133 
134  call chem_miequeryallband4d(mietable,idx,nb,offset,qaero(:,:,:,l),rh, &
135  tau,ssa,asy,rc=status)
136  verify_(status)
137 
138  tau8 = tau8 + tau
139  ssa8 = ssa8 + (ssa * tau)
140  asy8 = asy8 + asy * (ssa * tau)
141  end do
142 
143  tau = tau8
144  ssa = ssa8
145  asy = asy8
146 
147  return_(esmf_success)
148  end subroutine get_aerooptpropallaero4d
149 
150 !=================================================================
151 
152  subroutine get_aerotaulist(idx,ib,nn,qaero,rh,tau,rc)
154  integer, intent(IN ) :: idx
155  integer, intent(IN ) :: ib,nn
156  real, intent(IN ) :: qaero(nn)
157  real, intent(IN ) :: rh(nn)
158  real, intent(OUT) :: tau(nn)
159  integer, optional, intent(OUT) :: rc
160 
161  integer :: status
162  character(len=ESMF_MAXSTR) :: iam='Get_AeroTauList'
163 
164 
165  call chem_miequerytaulist(mietable,idx,real(ib),qaero,rh, &
166  tau,rc=status)
167  verify_(status)
168 
169  return_(esmf_success)
170  end subroutine get_aerotaulist
171 
172 !=================================================================
173 
174  integer function get_aeroindex(aerosol, rc ) result (idx)
175  character (len=*), intent(IN ) :: aerosol
176  integer, optional, intent(OUT) :: rc
177 
178  integer :: status
179  character(len=ESMF_MAXSTR) :: iam='Get_AeroIndex'
180 
181 
182  idx = chem_miequeryidx( mietable, aerosol, status )
183  verify_(status)
184 
185  return_(esmf_success)
186  end function get_aeroindex
187 
188  subroutine get_aeroindextableprops(na,aerosol,rh,routine,&
189  idx,nrh,bext,bsca,gasy,arh,irh,rc)
190 
191  integer, intent(in) :: na ! number of aerosols
192  character(len=*), intent(in) :: aerosol(:) ! list of aerosols
193  real, intent(in) :: rh(:,:) ! relative humidity
194  character(len=5), intent(in) :: routine ! name of routine
195 
196  integer, intent(out) :: idx(:) ! subindex of aerosol table
197  integer, intent(out) :: nrh(:) ! number of RH levels in table
198  real, intent(out) :: bext(:,:,:,:) ! mass extinction efficiency
199  real, intent(out) :: bsca(:,:,:,:) ! mass scattering efficiency
200  real, intent(out) :: gasy(:,:,:,:) ! asymmetry parameter
201  real, intent(out) :: arh(:,:,:) ! mapped slope from hi-res RH map
202  integer, intent(out) :: irh(:,:,:) ! mapped pointer from hi-res RH map
203  integer, optional, intent(out) :: rc
204 
205  integer :: status, i,j,k, ii, jj, isnap
206  integer :: num_radii, num_rh
207  character(len=ESMF_MAXSTR) :: iam='Get_AeroIndexTableProps'
208 
209  ii = size(rh,1)
210  jj = size(rh,2)
211 
212  ! Initialize outgoing matrices
213  idx = 0
214  nrh = 0
215  bext = 0.0
216  bsca = 0.0
217  gasy = 0.0
218  arh = 0.0
219  irh = 0
220 
221  do k = 1, na
222 
223  idx(k) = chem_miequeryidx( mietable, aerosol(k), status)
224  verify_(status)
225 
226  nrh(k) = mietable%vtableUse%nrh
227 
228  num_rh = size(mietable%vtableUse%bext,2)
229  num_radii = size(mietable%vtableUse%bext,3)
230 
231  select case(routine)
232  case('SOLAR') ! Use only the first 8 tables for SOLAR
233  bext(1:8,1:num_rh,1:num_radii,k) = mietable%vtableUse%bext(1:8,:,:)
234  bsca(1:8,1:num_rh,1:num_radii,k) = mietable%vtableUse%bsca(1:8,:,:)
235  gasy(1:8,1:num_rh,1:num_radii,k) = mietable%vtableUse%g( 1:8,:,:)
236  case('IRRAD') ! Use only the last 10 tables for IRRAD
237  bext(1:10,1:num_rh,1:num_radii,k) = mietable%vtableUse%bext(9:18,:,:)
238  bsca(1:10,1:num_rh,1:num_radii,k) = mietable%vtableUse%bsca(9:18,:,:)
239  gasy(1:10,1:num_rh,1:num_radii,k) = mietable%vtableUse%g( 9:18,:,:)
240  case default
241  assert_(.false.) ! We MUST have a case
242  end select
243 
244  do j = 1, jj
245  do i = 1, ii
246  isnap = int(((min(max(rh(i,j),0.0),0.99))+0.001)*1000)
247 
248  arh(i,j,k) = mietable%vtableUse%rha(isnap)
249  irh(i,j,k) = mietable%vtableUse%rhi(isnap)
250  end do
251  end do
252 
253  end do
254 
255  return_(esmf_success)
256  end subroutine get_aeroindextableprops
257 
258  subroutine create_aerooptprop(cf,rc)
260  type(esmf_config), intent(IN ) :: cf
261  integer, optional, intent(OUT) :: rc
262 
263  integer :: status
264  character(len=ESMF_MAXSTR) :: iam='Create_AeroOptProp'
265 
266  mietable = chem_miecreate(cf,status)
267  verify_(status)
268 
269  return_(esmf_success)
270  end subroutine create_aerooptprop
271 
272 
273 end module aerooptproptablemod
subroutine, public get_aerotaulist(idx, ib, nn, qaero, rh, tau, rc)
subroutine, public get_aerooptprop(idx, ib, qaero, rh, tau, ssa, asy, rc)
integer function, public get_aeroindex(aerosol, rc)
subroutine get_aerooptpropallaero3d(na, aerosol, nb, offset, qaero, rh, tau, ssa, asy, rc)
type(chem_mie), save mietable
subroutine, public get_aeroindextableprops(na, aerosol, rh, routine, idx, nrh, bext, bsca, gasy, arh, irh, rc)
logical, save mietablecreated
subroutine get_aerooptpropallaero4d(na, aerosol, nb, offset, qaero, rh, tau, ssa, asy, rc)
subroutine, public create_aerooptprop(cf, rc)
#define max(a, b)
Definition: mosaic_util.h:33
#define min(a, b)
Definition: mosaic_util.h:32