Line data Source code
1 : !--------------------------------------------------------------------------------
2 : ! Species summations for history
3 : !--------------------------------------------------------------------------------
4 : module species_sums_diags
5 :
6 : use shr_kind_mod, only : r8 => shr_kind_r8
7 : use shr_kind_mod, only : CL => SHR_KIND_CL
8 : use cam_history, only : addfld
9 : use cam_history, only : outfld
10 : use ppgrid, only : pver
11 : use spmd_utils, only : masterproc
12 : use cam_abortutils, only : endrun
13 : use cam_logfile, only : iulog
14 : use mo_chem_utls, only : get_spc_ndx
15 : use shr_expr_parser_mod , only : shr_exp_parse, shr_exp_item_t, shr_exp_list_destroy
16 :
17 : implicit none
18 : private
19 : public :: species_sums_init
20 : public :: species_sums_output
21 : public :: species_sums_readnl
22 : public :: species_sums_final
23 :
24 : type(shr_exp_item_t), pointer :: vmr_grps => null()
25 : type(shr_exp_item_t), pointer :: mmr_grps => null()
26 :
27 : integer, parameter :: maxlines = 200
28 : character(len=CL), allocatable :: vmr_sums(:)
29 : character(len=CL), allocatable :: mmr_sums(:)
30 :
31 : contains
32 :
33 : !-------------------------------------------------------------------
34 : !-------------------------------------------------------------------
35 1536 : subroutine species_sums_readnl(nlfile)
36 :
37 : use namelist_utils, only: find_group_name
38 : use units, only: getunit, freeunit
39 : use spmd_utils, only: mpicom, mpi_character, masterprocid
40 :
41 : ! args
42 : character(len=*), intent(in) :: nlfile ! filepath for file containing namelist input
43 :
44 : ! Local variables
45 : integer :: unitn, ierr, n,i
46 :
47 : namelist /species_sums_nl/ vmr_sums, mmr_sums
48 :
49 1536 : allocate( vmr_sums( maxlines ) )
50 308736 : vmr_sums(:) = ' '
51 1536 : allocate( mmr_sums( maxlines ) )
52 308736 : mmr_sums(:) = ' '
53 :
54 : ! Read namelist
55 1536 : if (masterproc) then
56 2 : unitn = getunit()
57 2 : open( unitn, file=trim(nlfile), status='old' )
58 2 : call find_group_name(unitn, 'species_sums_nl', status=ierr)
59 2 : if (ierr == 0) then
60 0 : read(unitn, species_sums_nl, iostat=ierr)
61 0 : if (ierr /= 0) then
62 0 : call endrun('species_sums_readnl:: ERROR reading namelist')
63 : end if
64 : end if
65 2 : close(unitn)
66 2 : call freeunit(unitn)
67 : end if
68 :
69 : ! Broadcast namelist variables
70 1536 : call mpi_bcast(vmr_sums,len(vmr_sums(1))*maxlines, mpi_character, masterprocid, mpicom, ierr)
71 1536 : call mpi_bcast(mmr_sums,len(mmr_sums(1))*maxlines, mpi_character, masterprocid, mpicom, ierr)
72 :
73 1536 : if (masterproc) then
74 :
75 2 : write(iulog,*) ' '
76 2 : write(iulog,*) 'species_sums_readnl -- vmr_sums :'
77 402 : n = count(len_trim(vmr_sums)>0)
78 2 : do i=1,n
79 2 : write(iulog,*) trim(vmr_sums(i))
80 : end do
81 :
82 2 : write(iulog,*) ' '
83 2 : write(iulog,*) 'species_sums_readnl -- mmr_sums :'
84 402 : n = count(len_trim(mmr_sums)>0)
85 2 : do i=1,n
86 2 : write(iulog,*) trim(mmr_sums(i))
87 : end do
88 :
89 : end if
90 :
91 :
92 1536 : end subroutine species_sums_readnl
93 : !--------------------------------------------------------------------------------
94 : !--------------------------------------------------------------------------------
95 1536 : subroutine species_sums_init
96 :
97 : type(shr_exp_item_t), pointer :: grp => null()
98 :
99 : integer :: j
100 :
101 : ! parse the terms of the summations
102 1536 : vmr_grps => shr_exp_parse( vmr_sums )
103 1536 : deallocate( vmr_sums )
104 1536 : mmr_grps => shr_exp_parse( mmr_sums )
105 1536 : deallocate( mmr_sums )
106 :
107 : ! add history fields
108 :
109 1536 : if (masterproc) write(iulog,*) 'species_sums_init -- VMR SUMS:'
110 1536 : grp => vmr_grps
111 1536 : do while(associated(grp))
112 :
113 0 : if (masterproc) then
114 0 : write(iulog,*) ' grp name : ',trim(grp%name)
115 :
116 0 : do j = 1, grp%n_terms
117 0 : write(iulog,'(f12.4,a,a)') grp%coeffs(j),' * ',trim(grp%vars(j))
118 : end do
119 : end if
120 :
121 0 : call addfld( trim(grp%name), (/ 'lev' /),'A', 'mole/mole','summation of species volume mixing ratios')
122 0 : grp => grp%next_item
123 : enddo
124 :
125 1536 : if (masterproc) write(iulog,*) 'species_sums_init -- MMR SUMS:'
126 1536 : grp => mmr_grps
127 1536 : do while(associated(grp))
128 :
129 0 : if (masterproc) then
130 0 : write(iulog,*) ' grp name : ',trim(grp%name)
131 :
132 0 : do j = 1, grp%n_terms
133 0 : write(iulog,'(f12.4,a,a)') grp%coeffs(j),' * ',trim(grp%vars(j))
134 : end do
135 : end if
136 :
137 0 : call addfld( trim(grp%name), (/ 'lev' /),'A', 'kg/kg','summation of species mass mixing ratios')
138 0 : grp => grp%next_item
139 : enddo
140 :
141 1536 : end subroutine species_sums_init
142 :
143 : !--------------------------------------------------------------------------------
144 : !--------------------------------------------------------------------------------
145 1489176 : subroutine species_sums_output( vmr, mmr, ncol, lchnk )
146 :
147 : real(r8), intent(in) :: vmr(:,:,:)
148 : real(r8), intent(in) :: mmr(:,:,:)
149 : integer, intent(in) :: ncol, lchnk
150 :
151 : integer :: j, spc_ndx
152 2978352 : real(r8) :: group_sum(ncol,pver)
153 : character(len=16) :: spc_name
154 : type(shr_exp_item_t), pointer :: grp
155 :
156 : ! output species groups ( or families )
157 1489176 : grp => vmr_grps
158 1489176 : do while(associated(grp))
159 : ! look up the corresponding species index ...
160 0 : group_sum(:,:) = 0._r8
161 0 : do j = 1, grp%n_terms
162 0 : spc_name = grp%vars(j)
163 0 : spc_ndx = get_spc_ndx( spc_name )
164 0 : if ( spc_ndx < 1 ) then
165 0 : call endrun('species_sums_output species name not found : '//trim(spc_name))
166 : endif
167 0 : group_sum(:ncol,:) = group_sum(:ncol,:) + grp%coeffs(j)*vmr(:ncol,:,spc_ndx)
168 : enddo
169 0 : call outfld( trim(grp%name), group_sum(:ncol,:), ncol, lchnk )
170 0 : grp => grp%next_item
171 : end do
172 :
173 1489176 : grp => mmr_grps
174 1489176 : do while(associated(grp))
175 : ! look up the corresponding species index ...
176 0 : group_sum(:,:) = 0._r8
177 0 : do j = 1, grp%n_terms
178 0 : spc_name = grp%vars(j)
179 0 : spc_ndx = get_spc_ndx( spc_name )
180 0 : if ( spc_ndx < 1 ) then
181 0 : call endrun('species_sums_output species name not found : '//trim(spc_name))
182 : endif
183 0 : group_sum(:ncol,:) = group_sum(:ncol,:) + grp%coeffs(j)*mmr(:ncol,:,spc_ndx)
184 : enddo
185 0 : call outfld( trim(grp%name), group_sum(:ncol,:), ncol, lchnk )
186 0 : grp => grp%next_item
187 : end do
188 :
189 1489176 : end subroutine species_sums_output
190 :
191 : !--------------------------------------------------------------------------------
192 : !--------------------------------------------------------------------------------
193 1536 : subroutine species_sums_final
194 :
195 1536 : if (associated(vmr_grps)) call shr_exp_list_destroy(vmr_grps)
196 1536 : if (associated(mmr_grps)) call shr_exp_list_destroy(mmr_grps)
197 :
198 1536 : end subroutine species_sums_final
199 :
200 : end module species_sums_diags
|