SCALE-RM
scale_aetracer_kajino13.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
13 !-------------------------------------------------------------------------------
15  !-----------------------------------------------------------------------------
16  !
17  !++ used modules
18  !
19  use scale_precision
20  use scale_stdio
21  !-----------------------------------------------------------------------------
22  implicit none
23  private
24  !-----------------------------------------------------------------------------
25  !
26  !++ Public procedure
27  !
28  public :: aetracer_kajino13_setup
29 
30  include "inc_aetracer_kajino13.h"
31 
32  !-----------------------------------------------------------------------------
33 contains
34  !-----------------------------------------------------------------------------
35  subroutine aetracer_kajino13_setup
36  use scale_process, only: &
38  implicit none
39 
40  integer, allocatable :: aero_idx(:,:,:,:)
41  integer :: n_kap_max, n_siz_max, ncat_max
42  integer :: NASIZ(3), NAKAP(3)
43  character(len=H_SHORT) :: attribute, catego, aunit
44 
45  namelist / param_tracer_kajino13 / &
46  ae_ctg, &
47  nasiz, &
48  nakap
49 
50  integer :: m, ierr, ik, ic, ia0, is0
51 
52  if( io_l ) write(io_fid_log,*)
53  if( io_l ) write(io_fid_log,*) '+++ READ NUMBER of TRACER for Aerosol'
54 
55  ncat_max = max( ic_mix, ic_sea, ic_dus )
56 
57  nasiz(:) = 64
58  nakap(:) = 1
59 
60  rewind(io_fid_conf)
61  read(io_fid_conf,nml=param_tracer_kajino13,iostat=ierr)
62 
63  if( ierr < 0 ) then !--- missing
64  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
65  elseif( ierr > 0 ) then !--- fatal error
66  write(*,*) 'xxx Not appropriate names in namelist PARAM_TRACER_KAJINO13, Check!'
67  call prc_mpistop
68  end if
69 
70  if( io_l ) write(io_fid_log,nml=param_tracer_kajino13)
71 
72  if( ae_ctg > ncat_max ) then
73  write(*,*) 'xxx AE_CTG should be smaller than', ncat_max+1, 'stop'
74  call prc_mpistop
75  endif
76 
77  allocate( nsiz(ae_ctg) )
78  allocate( nkap(ae_ctg) )
79 
80  nkap(1:ae_ctg) = nakap(1:ae_ctg)
81  nsiz(1:ae_ctg) = nasiz(1:ae_ctg)
82 
83  if( maxval( nkap ) /= 1 .OR. minval( nkap ) /= 1 ) then
84  write(*,*) 'xxx NKAP(:) /= 1 is not supported now, Stop!'
85  call prc_mpistop
86  end if
87 
88  qa_ae = 0
89 ! do ia0 = 1, N_ATR
90  do ic = 1, ae_ctg
91  do ik = 1, nkap(ic)
92  do is0 = 1, nsiz(ic)
93  do ia0 = 1, n_atr
94  qa_ae = qa_ae + 1
95  enddo
96  enddo
97  enddo
98  enddo
99  qa_ae = qa_ae + gas_ctg
100  allocate( aq_ae_name(qa_ae) )
101  allocate( aq_ae_desc(qa_ae) )
102  allocate( aq_ae_unit(qa_ae) )
103 
104  n_siz_max = 0
105  n_kap_max = 0
106  do ic = 1, ae_ctg
107  n_siz_max = max(n_siz_max, nsiz(ic))
108  n_kap_max = max(n_kap_max, nkap(ic))
109  enddo
110 
111  allocate( aero_idx(n_atr,ae_ctg,n_kap_max,n_siz_max) )
112  m = 0
113 ! do ia0 = 1, N_ATR
114  do ic = 1, ae_ctg
115  do ik = 1, nkap(ic)
116  do is0 = 1, nsiz(ic)
117  do ia0 = 1, n_atr
118  m = m+1
119  aero_idx(ia0,ic,ik,is0) = m
120  enddo
121  enddo
122  enddo
123  enddo
124 
125  !-----------------------------------------------------------------------------
126  !
127  !++ calculate each category and aerosol
128  !
129  !-----------------------------------------------------------------------------
130  ic = qa_ae-gas_ctg+ig_h2so4
131  write(aq_ae_unit(ic),'(a)') 'kg/kg'
132  ic = qa_ae-gas_ctg+ig_cgas
133  write(aq_ae_unit(ic),'(a)') 'kg/kg'
134 
135 ! do ia0 = 1, N_ATR
136 ! if( ia0 == 1 ) then
137 ! write(attribute,'(a)') "Number"
138 ! write(aunit,'(a)') "num/kg"
139 ! elseif( ia0 == 2 ) then
140 ! write(attribute,'(a)') "Section"
141 ! write(aunit,'(a)') "m2/kg"
142 ! elseif( ia0 == 3 ) then
143 ! write(attribute,'(a)') "Volume"
144 ! write(aunit,'(a)') "m3/kg"
145 ! elseif( ia0 == 4 ) then
146 ! write(attribute,'(a)') "Mass"
147 ! write(aunit,'(a)') "kg/kg"
148 ! elseif( ia0 == 5 ) then
149 ! write(attribute,'(a)') "kpXmass"
150 ! write(aunit,'(a)') "kg/kg"
151 ! endif
152  do ic = 1, ae_ctg !aerosol category
153  do ik = 1, nkap(ic) !kappa bin
154  do is0 = 1, nsiz(ic)
155  do ia0 = 1, n_atr
156  if( ia0 == 1 ) then
157  write(attribute,'(a)') "Number"
158  write(aunit,'(a)') "num/kg"
159  elseif( ia0 == 2 ) then
160  write(attribute,'(a)') "Section"
161  write(aunit,'(a)') "m2/kg"
162  elseif( ia0 == 3 ) then
163  write(attribute,'(a)') "Volume"
164  write(aunit,'(a)') "m3/kg"
165  elseif( ia0 == 4 ) then
166  write(attribute,'(a)') "Mass"
167  write(aunit,'(a)') "kg/kg"
168  elseif( ia0 == 5 ) then
169  write(attribute,'(a)') "kXm"
170  write(aunit,'(a)') "kg/kg"
171  endif
172  if( ic == ic_mix ) then
173  write(catego,'(a)') "Sulf_"
174  elseif( ic == ic_sea ) then
175  write(catego,'(a)') "Salt_"
176  elseif( ic == ic_dus ) then
177  write(catego,'(a)') "Dust_"
178  endif
179  write(aq_ae_unit(aero_idx(ia0,ic,ik,is0)),'(a)') trim(aunit)
180  write(aq_ae_name(aero_idx(ia0,ic,ik,is0)),'(a,a,i0)') trim(catego), trim(attribute), is0
181  write(aq_ae_desc(aero_idx(ia0,ic,ik,is0)),'(a,a,a,i0)') trim(attribute), ' mixing radio of ', trim(catego), is0
182  enddo
183  enddo
184  enddo
185  enddo
186  ic = qa_ae-gas_ctg+ig_h2so4
187  write(aq_ae_name(ic),'(a)') 'H2SO4_Gas'
188  ic = qa_ae-gas_ctg+ig_cgas
189  write(aq_ae_name(ic),'(a)') 'Condensable_GAS'
190 
191  ic = qa_ae-gas_ctg+ig_h2so4
192  write(aq_ae_desc(ic),'(a)') 'Mixing ratio of H2SO4 Gas'
193  ic = qa_ae-gas_ctg+ig_cgas
194  write(aq_ae_desc(ic),'(a)') 'Mixing ratio of Condensable GAS'
195 
196  deallocate(aero_idx)
197 
198  return
199  end subroutine aetracer_kajino13_setup
200 
201 end module scale_aetracer_kajino13
subroutine, public prc_mpistop
Abort MPI.
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
module TRACER / kajino13
module STDIO
Definition: scale_stdio.F90:12
module PROCESS
subroutine, public aetracer_kajino13_setup
module PRECISION
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56