Line data Source code
1 : module mo_imp_sol
2 : use shr_kind_mod, only : r8 => shr_kind_r8
3 : use chem_mods, only : clscnt4, gas_pcnst, clsmap, veclen
4 : use cam_logfile, only : iulog
5 : implicit none
6 : private
7 : public :: imp_slv_inti, imp_sol
8 : save
9 : real(r8), parameter :: rel_err = 1.e-3_r8
10 : real(r8), parameter :: high_rel_err = 1.e-4_r8
11 : !-----------------------------------------------------------------------
12 : ! Newton-Raphson iteration limits
13 : !-----------------------------------------------------------------------
14 : integer, parameter :: itermax = 11
15 : integer, parameter :: cut_limit = 5
16 : real(r8), parameter :: sol_min = 1.e-20_r8
17 : real(r8), parameter :: small = 1.e-40_r8
18 : real(r8) :: epsilon(clscnt4)
19 : logical :: factor(itermax)
20 : contains
21 1536 : subroutine imp_slv_inti
22 : !-----------------------------------------------------------------------
23 : ! ... Initialize the implict solver
24 : !-----------------------------------------------------------------------
25 : use mo_chem_utls, only : get_spc_ndx
26 : implicit none
27 : !-----------------------------------------------------------------------
28 : ! ... Local variables
29 : !-----------------------------------------------------------------------
30 : integer :: m, ox_ndx, o3a_ndx
31 : real(r8) :: eps(gas_pcnst)
32 18432 : factor(:) = .true.
33 405504 : eps(:) = rel_err
34 1536 : ox_ndx = get_spc_ndx( 'OX' )
35 1536 : if( ox_ndx < 1 ) then
36 1536 : ox_ndx = get_spc_ndx( 'O3' )
37 : end if
38 1536 : if( ox_ndx > 0 ) then
39 1536 : eps(ox_ndx) = high_rel_err
40 : end if
41 1536 : m = get_spc_ndx( 'NO' )
42 1536 : if( m > 0 ) then
43 1536 : eps(m) = high_rel_err
44 : end if
45 1536 : m = get_spc_ndx( 'NO2' )
46 1536 : if( m > 0 ) then
47 1536 : eps(m) = high_rel_err
48 : end if
49 1536 : m = get_spc_ndx( 'NO3' )
50 1536 : if( m > 0 ) then
51 1536 : eps(m) = high_rel_err
52 : end if
53 1536 : m = get_spc_ndx( 'HNO3' )
54 1536 : if( m > 0 ) then
55 1536 : eps(m) = high_rel_err
56 : end if
57 1536 : m = get_spc_ndx( 'HO2NO2' )
58 1536 : if( m > 0 ) then
59 1536 : eps(m) = high_rel_err
60 : end if
61 1536 : m = get_spc_ndx( 'N2O5' )
62 1536 : if( m > 0 ) then
63 1536 : eps(m) = high_rel_err
64 : end if
65 1536 : m = get_spc_ndx( 'OH' )
66 1536 : if( m > 0 ) then
67 1536 : eps(m) = high_rel_err
68 : end if
69 1536 : m = get_spc_ndx( 'HO2' )
70 1536 : if( m > 0 ) then
71 1536 : eps(m) = high_rel_err
72 : end if
73 1536 : o3a_ndx = get_spc_ndx( 'O3A' )
74 1536 : if( o3a_ndx > 0 ) then
75 0 : eps(m) = high_rel_err
76 : end if
77 1536 : m = get_spc_ndx( 'XNO' )
78 1536 : if( m > 0 ) then
79 0 : eps(m) = high_rel_err
80 : end if
81 1536 : m = get_spc_ndx( 'XNO2' )
82 1536 : if( m > 0 ) then
83 0 : eps(m) = high_rel_err
84 : end if
85 1536 : m = get_spc_ndx( 'XNO3' )
86 1536 : if( m > 0 ) then
87 0 : eps(m) = high_rel_err
88 : end if
89 1536 : m = get_spc_ndx( 'XHNO3' )
90 1536 : if( m > 0 ) then
91 0 : eps(m) = high_rel_err
92 : end if
93 1536 : m = get_spc_ndx( 'XHO2NO2' )
94 1536 : if( m > 0 ) then
95 0 : eps(m) = high_rel_err
96 : end if
97 1536 : m = get_spc_ndx( 'XNO2NO3' )
98 1536 : if( m > 0 ) then
99 0 : eps(m) = high_rel_err
100 : end if
101 1536 : m = get_spc_ndx( 'NO2XNO3' )
102 1536 : if( m > 0 ) then
103 0 : eps(m) = high_rel_err
104 : end if
105 402432 : do m = 1,clscnt4
106 402432 : epsilon(m) = eps(clsmap(m,4))
107 : end do
108 1536 : end subroutine imp_slv_inti
109 72960 : subroutine imp_sol( base_sol, reaction_rates, het_rates, extfrc, delt, &
110 72960 : ncol, nlev, lchnk, prod_out, loss_out )
111 : !-----------------------------------------------------------------------
112 : ! ... imp_sol advances the volumetric mixing ratio
113 : ! forward one time step via the fully implicit euler scheme.
114 : ! this source is meant for vector architectures such as the
115 : ! nec sx6 and cray x1
116 : !-----------------------------------------------------------------------
117 : use chem_mods, only : rxntot, extcnt, nzcnt, permute, cls_rxt_cnt
118 : use mo_tracname, only : solsym
119 : use mo_lin_matrix, only : linmat
120 : use mo_nln_matrix, only : nlnmat
121 : use mo_lu_factor, only : lu_fac
122 : use mo_lu_solve, only : lu_slv
123 : use mo_prod_loss, only : imp_prod_loss
124 : use mo_indprd, only : indprd
125 : use time_manager, only : get_nstep
126 : use perf_mod, only : t_startf, t_stopf
127 : implicit none
128 : !-----------------------------------------------------------------------
129 : ! ... dummy args
130 : !-----------------------------------------------------------------------
131 : integer, intent(in) :: ncol ! columns in chunck
132 : integer, intent(in) :: nlev
133 : integer, intent(in) :: lchnk ! chunk id
134 : real(r8), intent(in) :: delt ! time step (s)
135 : real(r8), intent(in) :: reaction_rates(ncol*nlev,max(1,rxntot)) ! rxt rates (1/cm^3/s)
136 : real(r8), intent(in) :: extfrc(ncol*nlev,max(1,extcnt)) ! external in-situ forcing (1/cm^3/s)
137 : real(r8), intent(in) :: het_rates(ncol*nlev,max(1,gas_pcnst)) ! washout rates (1/s)
138 : real(r8), intent(inout) :: base_sol(ncol*nlev,gas_pcnst) ! species mixing ratios (vmr)
139 : real(r8), intent(out) :: prod_out(ncol*nlev,max(1,clscnt4))
140 : real(r8), intent(out) :: loss_out(ncol*nlev,max(1,clscnt4))
141 : !-----------------------------------------------------------------------
142 : ! ... local variables
143 : !-----------------------------------------------------------------------
144 : integer :: nr_iter
145 : integer :: ofl
146 : integer :: ofu
147 : integer :: avec_len
148 : integer :: bndx ! base index
149 : integer :: cndx ! class index
150 : integer :: pndx ! permuted class index
151 : integer :: i,m
152 : integer :: fail_cnt(veclen)
153 : integer :: cut_cnt(veclen)
154 : integer :: stp_con_cnt(veclen)
155 : integer :: nstep
156 : real(r8) :: interval_done(veclen)
157 : real(r8) :: dt(veclen)
158 : real(r8) :: dti(veclen)
159 : real(r8) :: max_delta(max(1,clscnt4))
160 145920 : real(r8) :: ind_prd(ncol*nlev,max(1,clscnt4))
161 : logical :: convergence
162 : integer :: chnkpnts ! total spatial points in chunk; ncol*ncol
163 145920 : logical :: diags_out(ncol*nlev,max(1,clscnt4))
164 : real(r8) :: sys_jac_blk(veclen,max(1,nzcnt))
165 : real(r8) :: lin_jac_blk(veclen,max(1,nzcnt))
166 : real(r8) :: solution_blk(veclen,max(1,clscnt4))
167 : real(r8) :: forcing_blk(veclen,max(1,clscnt4))
168 : real(r8) :: iter_invariant_blk(veclen,max(1,clscnt4))
169 : real(r8) :: prod_blk(veclen,max(1,clscnt4))
170 : real(r8) :: loss_blk(veclen,max(1,clscnt4))
171 : real(r8) :: ind_prd_blk(veclen,max(1,clscnt4))
172 : real(r8) :: sbase_sol_blk(veclen,gas_pcnst)
173 : real(r8) :: wrk_blk(veclen)
174 : logical :: spc_conv_blk(veclen,max(1,clscnt4))
175 : logical :: cls_conv_blk(veclen)
176 : logical :: time_stp_done_blk(veclen)
177 : real(r8) :: reaction_rates_blk(veclen,max(1,rxntot))
178 : real(r8) :: extfrc_blk(veclen,max(1,extcnt))
179 : real(r8) :: het_rates_blk(veclen,max(1,gas_pcnst))
180 : real(r8) :: base_sol_blk(veclen,gas_pcnst)
181 72960 : chnkpnts = ncol*nlev
182 8793927168 : prod_out = 0._r8
183 8793927168 : loss_out = 0._r8
184 8793927168 : diags_out = .false.
185 : !-----------------------------------------------------------------------
186 : ! ... class independent forcing
187 : !-----------------------------------------------------------------------
188 : if( cls_rxt_cnt(1,4) > 0 .or. extcnt > 0 ) then
189 : call indprd( 4, ind_prd, clscnt4, base_sol, extfrc, &
190 72960 : reaction_rates, chnkpnts )
191 : else
192 : do m = 1,clscnt4
193 : ind_prd(:,m) = 0._r8
194 : end do
195 : end if
196 72960 : nstep = get_nstep()
197 72960 : ofl = 1
198 : chnkpnts_loop : do
199 1050624 : ofu = min( chnkpnts,ofl + veclen - 1 )
200 1050624 : avec_len = (ofu - ofl) + 1
201 19485923328 : reaction_rates_blk(1:avec_len,:) = reaction_rates(ofl:ofu,:)
202 694462464 : extfrc_blk(1:avec_len,:) = extfrc(ofl:ofu,:)
203 9119416320 : het_rates_blk(1:avec_len,:) = het_rates(ofl:ofu,:)
204 9050075136 : ind_prd_blk(1:avec_len,:) = ind_prd(ofl:ofu,:)
205 9119416320 : base_sol_blk(1:avec_len,:) = base_sol(ofl:ofu,:)
206 34670592 : cls_conv_blk(1:avec_len) = .false.
207 34670592 : dt(1:avec_len) = delt
208 34670592 : cut_cnt(1:avec_len) = 0
209 34670592 : fail_cnt(1:avec_len) = 0
210 34670592 : stp_con_cnt(1:avec_len) = 0
211 34670592 : interval_done(1:avec_len) = 0._r8
212 34670592 : time_stp_done_blk(1:avec_len) = .false.
213 : !-----------------------------------------------------------------------
214 : ! ... time step loop
215 : !-----------------------------------------------------------------------
216 : time_step_loop : do
217 35722848 : dti(1:avec_len) = 1._r8 / dt(1:avec_len)
218 : !-----------------------------------------------------------------------
219 : ! ... transfer from base to class array
220 : !-----------------------------------------------------------------------
221 275276064 : do cndx = 1,clscnt4
222 274225392 : bndx = clsmap(cndx,4)
223 274225392 : pndx = permute(cndx,4)
224 9050488608 : do i = 1, avec_len
225 9049437936 : solution_blk(i,pndx) = base_sol_blk(i,bndx)
226 : end do
227 : end do
228 277377408 : do m = 1,gas_pcnst
229 9119832960 : sbase_sol_blk(1:avec_len,m) = base_sol_blk(1:avec_len,m)
230 : end do
231 : !-----------------------------------------------------------------------
232 : ! ... set the iteration invariant part of the function f(y)
233 : !-----------------------------------------------------------------------
234 : if( cls_rxt_cnt(1,4) > 0 .or. extcnt > 0 ) then
235 275276064 : do m = 1,clscnt4
236 9050488608 : do i = 1, avec_len
237 9049437936 : iter_invariant_blk(i,m) = dti(i) * solution_blk(i,m) + ind_prd_blk(i,m)
238 : end do
239 : end do
240 : else
241 : do m = 1,clscnt4
242 : do i = 1, avec_len
243 : iter_invariant_blk(i,m) = dti(i) * solution_blk(i,m)
244 : end do
245 : end do
246 : end if
247 : !-----------------------------------------------------------------------
248 : ! ... the linear component
249 : !-----------------------------------------------------------------------
250 1050672 : if( cls_rxt_cnt(2,4) > 0 ) then
251 1050672 : call t_startf( 'lin_mat' )
252 : call linmat( avec_len, lin_jac_blk, base_sol_blk, &
253 1050672 : reaction_rates_blk, het_rates_blk )
254 1050672 : call t_stopf( 'lin_mat' )
255 : end if
256 : !=======================================================================
257 : ! the newton-raphson iteration for f(y) = 0
258 : !=======================================================================
259 4411728 : iter_loop : do nr_iter = 1,itermax
260 : !-----------------------------------------------------------------------
261 : ! ... the non-linear component
262 : !-----------------------------------------------------------------------
263 4411704 : if( factor(nr_iter) ) then
264 4411704 : call t_startf( 'nln_mat' )
265 : call nlnmat( avec_len, sys_jac_blk, base_sol_blk, &
266 4411704 : reaction_rates_blk, lin_jac_blk, dti )
267 4411704 : call t_stopf( 'nln_mat' )
268 : !-----------------------------------------------------------------------
269 : ! ... factor the "system" matrix
270 : !-----------------------------------------------------------------------
271 4411704 : call t_startf( 'lu_fac' )
272 4411704 : call lu_fac( avec_len, sys_jac_blk )
273 4411704 : call t_stopf( 'lu_fac' )
274 : end if
275 : !-----------------------------------------------------------------------
276 : ! ... form f(y)
277 : !-----------------------------------------------------------------------
278 4411704 : call t_startf( 'prod_loss' )
279 : call imp_prod_loss( avec_len, prod_blk, loss_blk, &
280 4411704 : base_sol_blk, reaction_rates_blk, het_rates_blk )
281 4411704 : call t_stopf( 'prod_loss' )
282 1155866448 : do m = 1,clscnt4
283 38002418256 : do i = 1, avec_len
284 73693103616 : forcing_blk(i,m) = solution_blk(i,m)*dti(i) &
285 >11169*10^7 : - (iter_invariant_blk(i,m) + prod_blk(i,m) - loss_blk(i,m))
286 : end do
287 : end do
288 : !-----------------------------------------------------------------------
289 : ! ... solve for the mixing ratio at t(n+1)
290 : !-----------------------------------------------------------------------
291 4411704 : call t_startf( 'lu_slv' )
292 4411704 : call lu_slv( avec_len, sys_jac_blk, forcing_blk )
293 4411704 : call t_stopf( 'lu_slv' )
294 1155866448 : do m = 1,clscnt4
295 38002418256 : do i = 1, avec_len
296 37998006552 : if( .not. cls_conv_blk(i) )then
297 28737660780 : solution_blk(i,m) = solution_blk(i,m) + forcing_blk(i,m)
298 : else
299 8108891028 : forcing_blk(i,m) = 0._r8
300 : endif
301 : end do
302 : end do
303 : !-----------------------------------------------------------------------
304 : ! ... convergence measures and test
305 : !-----------------------------------------------------------------------
306 5462376 : conv_chk : if( nr_iter > 1 ) then
307 : !-----------------------------------------------------------------------
308 : ! ... check for convergence
309 : !-----------------------------------------------------------------------
310 880590384 : do cndx = 1,clscnt4
311 877229352 : pndx = permute(cndx,4)
312 877229352 : bndx = clsmap(cndx,4)
313 28948568616 : do i = 1, avec_len
314 28948568616 : if ( abs( solution_blk(i,pndx) ) > sol_min ) then
315 21464856886 : wrk_blk(i) = abs( forcing_blk(i,pndx)/solution_blk(i,pndx) )
316 : else
317 6606482378 : wrk_blk(i) = 0._r8
318 : endif
319 : enddo
320 29825797968 : max_delta(cndx) = maxval( wrk_blk(1:avec_len) )
321 28948568616 : do i = 1, avec_len
322 28071339264 : solution_blk(i,pndx) = max( 0._r8,solution_blk(i,pndx) )
323 28071339264 : base_sol_blk(i,bndx) = solution_blk(i,pndx)
324 28948568616 : if ( abs( forcing_blk(i,pndx) ) > small ) then
325 14358533244 : spc_conv_blk(i,cndx) = abs(forcing_blk(i,pndx)) <= epsilon(cndx)*abs(solution_blk(i,pndx))
326 : else
327 13712806020 : spc_conv_blk(i,cndx) = .true.
328 : endif
329 : enddo
330 >11579*10^7 : where( spc_conv_blk(1:avec_len,cndx) .and. .not.diags_out(ofl:ofu,cndx) )
331 : ! capture output production and loss diagnostics at converged ponits
332 877229352 : prod_out(ofl:ofu,cndx) = prod_blk(1:avec_len,cndx) + ind_prd_blk(1:avec_len,cndx)
333 : loss_out(ofl:ofu,cndx) = loss_blk(1:avec_len,cndx)
334 : diags_out(ofl:ofu,cndx) = .true.
335 : endwhere
336 : end do
337 110914056 : do i = 1, avec_len
338 110914056 : if( .not. cls_conv_blk(i) ) then
339 9869994953 : cls_conv_blk(i) = all( spc_conv_blk(i,:) )
340 : end if
341 : end do
342 43230192 : convergence = all( cls_conv_blk(:) )
343 3361032 : if( convergence ) then
344 : exit iter_loop
345 : end if
346 : else conv_chk
347 : !-----------------------------------------------------------------------
348 : ! ... limit iterate
349 : !-----------------------------------------------------------------------
350 275276064 : do m = 1,clscnt4
351 9050488608 : do i = 1, avec_len
352 9049437936 : solution_blk(i,m) = max( 0._r8,solution_blk(i,m) )
353 : end do
354 : end do
355 : !-----------------------------------------------------------------------
356 : ! ... transfer latest solution back to base array
357 : !-----------------------------------------------------------------------
358 275276064 : do cndx = 1,clscnt4
359 274225392 : pndx = permute(cndx,4)
360 274225392 : bndx = clsmap(cndx,4)
361 9050488608 : do i = 1, avec_len
362 9049437936 : base_sol_blk(i,bndx) = solution_blk(i,pndx)
363 : end do
364 : end do
365 : end if conv_chk
366 : end do iter_loop
367 : !-----------------------------------------------------------------------
368 : ! ... check for newton-raphson convergence
369 : !-----------------------------------------------------------------------
370 34672176 : do i = 1,avec_len
371 34672176 : if( .not. cls_conv_blk(i) ) then
372 24 : fail_cnt(i) = fail_cnt(i) + 1
373 : write(iulog,'('' imp_sol: time step '',1p,g15.7,'' failed to converge @ (lchnk,vctrpos,nstep) = '',3i8)') &
374 24 : dt(i),lchnk,ofl+i-1,nstep
375 24 : stp_con_cnt(i) = 0
376 24 : if( cut_cnt(i) < cut_limit ) then
377 24 : cut_cnt(i) = cut_cnt(i) + 1
378 24 : if( cut_cnt(i) < cut_limit ) then
379 24 : dt(i) = .5_r8 * dt(i)
380 : else
381 0 : dt(i) = .1_r8 * dt(i)
382 : end if
383 6336 : base_sol_blk(i,:) = sbase_sol_blk(i,:)
384 : else
385 : write(iulog,'('' imp_sol: step failed to converge @ (lchnk,vctrpos,nstep,dt,time) = '',3i8,1p,2g15.7)') &
386 0 : lchnk,ofl+i-1,nstep,dt(i),interval_done+dt(i)
387 0 : do m = 1,clscnt4
388 0 : if( .not. spc_conv_blk(i,m) ) then
389 0 : write(iulog,'(1x,a16,1x,1pe10.3)') solsym(clsmap(m,4)), max_delta(m)
390 : end if
391 : end do
392 0 : cls_conv_blk(i) = .true.
393 0 : if( .not. time_stp_done_blk(i) ) then
394 0 : interval_done(i) = interval_done(i) + dt(i)
395 0 : time_stp_done_blk(i) = abs( delt - interval_done(i) ) <= .0001_r8
396 : endif
397 : end if
398 33621480 : elseif( .not. time_stp_done_blk(i) ) then
399 33619992 : interval_done(i) = interval_done(i) + dt(i)
400 33619992 : time_stp_done_blk(i) = abs( delt - interval_done(i) ) <= .0001_r8
401 33619992 : stp_con_cnt(i) = stp_con_cnt(i) + 1
402 33619992 : if( .not. time_stp_done_blk(i) ) then
403 24 : if( stp_con_cnt(i) >= 2 ) then
404 0 : dt(i) = 2._r8*dt(i)
405 0 : stp_con_cnt(i) = 0
406 : end if
407 24 : dt(i) = min( dt(i),delt-interval_done(i) )
408 : else
409 8875671552 : base_sol(ofl+i-1,1:gas_pcnst) = base_sol_blk(i,1:gas_pcnst)
410 : endif
411 : endif
412 : end do
413 34672068 : convergence = all( cls_conv_blk(:) )
414 34672176 : do i = 1,avec_len
415 34672176 : if( cls_conv_blk(i) .and. .not. time_stp_done_blk(i) ) then
416 24 : cls_conv_blk(i) = .false.
417 : endif
418 : end do
419 1050672 : if( .not. convergence ) then
420 : cycle time_step_loop
421 : endif
422 : !-----------------------------------------------------------------------
423 : ! ... check for time step done
424 : !-----------------------------------------------------------------------
425 35721924 : if( all( time_stp_done_blk(1:avec_len) ) ) then
426 : exit time_step_loop
427 : end if
428 : end do time_step_loop
429 1050624 : ofl = ofu + 1
430 1050624 : if( ofl > chnkpnts ) then
431 : exit chnkpnts_loop
432 : end if
433 : end do chnkpnts_loop
434 72960 : end subroutine imp_sol
435 : end module mo_imp_sol
|