Line data Source code
1 : module stepon
2 :
3 : use shr_kind_mod, only: r8 => shr_kind_r8
4 : use spmd_utils, only: iam, mpicom
5 : use ppgrid, only: begchunk, endchunk
6 :
7 : use physics_types, only: physics_state, physics_tend
8 : use dyn_comp, only: dyn_import_t, dyn_export_t
9 :
10 : use perf_mod, only: t_startf, t_stopf, t_barrierf
11 : use cam_abortutils, only: endrun
12 :
13 : use parallel_mod, only: par
14 : use dimensions_mod, only: nelemd
15 :
16 : use aerosol_properties_mod, only: aerosol_properties
17 : use aerosol_state_mod, only: aerosol_state
18 : use microp_aero, only: aerosol_state_object, aerosol_properties_object
19 :
20 : implicit none
21 : private
22 : save
23 :
24 : public stepon_init
25 : public stepon_run1
26 : public stepon_run2
27 : public stepon_run3
28 : public stepon_final
29 :
30 : class(aerosol_properties), pointer :: aero_props_obj => null()
31 : logical :: aerosols_transported = .false.
32 :
33 : !=========================================================================================
34 : contains
35 : !=========================================================================================
36 :
37 17664 : subroutine stepon_init(dyn_in, dyn_out )
38 :
39 : use cam_history, only: addfld, add_default, horiz_only
40 : use constituents, only: pcnst, cnst_name, cnst_longname
41 : use dimensions_mod, only: fv_nphys, cnst_name_gll, cnst_longname_gll, qsize
42 :
43 : ! arguments
44 : type (dyn_import_t), intent(inout) :: dyn_in ! Dynamics import container
45 : type (dyn_export_t), intent(inout) :: dyn_out ! Dynamics export container
46 :
47 : ! local variables
48 : integer :: m, m_cnst
49 : !----------------------------------------------------------------------------
50 : ! These fields on dynamics grid are output before the call to d_p_coupling.
51 10752 : do m_cnst = 1, qsize
52 0 : call addfld(trim(cnst_name_gll(m_cnst))//'_gll', (/ 'lev' /), 'I', 'kg/kg', &
53 18432 : trim(cnst_longname_gll(m_cnst)), gridname='GLL')
54 0 : call addfld(trim(cnst_name_gll(m_cnst))//'dp_gll', (/ 'lev' /), 'I', 'kg/kg', &
55 19968 : trim(cnst_longname_gll(m_cnst))//'*dp', gridname='GLL')
56 : end do
57 3072 : call addfld('U_gll' ,(/ 'lev' /), 'I', 'm/s ','U wind on gll grid',gridname='GLL')
58 3072 : call addfld('V_gll' ,(/ 'lev' /), 'I', 'm/s ','V wind on gll grid',gridname='GLL')
59 3072 : call addfld('T_gll' ,(/ 'lev' /), 'I', 'K ' ,'T on gll grid' ,gridname='GLL')
60 3072 : call addfld('dp_ref_gll' ,(/ 'lev' /), 'I', ' ' ,'dp dry / dp_ref on gll grid' ,gridname='GLL')
61 1536 : call addfld('PSDRY_gll' ,horiz_only , 'I', 'Pa ' ,'psdry on gll grid' ,gridname='GLL')
62 1536 : call addfld('PS_gll' ,horiz_only , 'I', 'Pa ' ,'ps on gll grid' ,gridname='GLL')
63 1536 : call addfld('PHIS_gll' ,horiz_only , 'I', 'Pa ' ,'PHIS on gll grid' ,gridname='GLL')
64 :
65 : ! Fields for initial condition files
66 3072 : call addfld('U&IC', (/ 'lev' /), 'I', 'm/s', 'Zonal wind', gridname='GLL' )
67 3072 : call addfld('V&IC', (/ 'lev' /), 'I', 'm/s', 'Meridional wind',gridname='GLL' )
68 : ! Don't need to register U&IC V&IC as vector components since we don't interpolate IC files
69 1536 : call add_default('U&IC',0, 'I')
70 1536 : call add_default('V&IC',0, 'I')
71 :
72 1536 : call addfld('PS&IC', horiz_only, 'I', 'Pa', 'Surface pressure', gridname='GLL')
73 3072 : call addfld('T&IC', (/ 'lev' /), 'I', 'K', 'Temperature', gridname='GLL')
74 1536 : call add_default('PS&IC', 0, 'I')
75 1536 : call add_default('T&IC', 0, 'I')
76 :
77 64512 : do m_cnst = 1,pcnst
78 62976 : call addfld(trim(cnst_name(m_cnst))//'&IC', (/ 'lev' /), 'I', 'kg/kg', &
79 188928 : trim(cnst_longname(m_cnst)), gridname='GLL')
80 64512 : call add_default(trim(cnst_name(m_cnst))//'&IC', 0, 'I')
81 : end do
82 :
83 : ! get aerosol properties
84 1536 : aero_props_obj => aerosol_properties_object()
85 :
86 1536 : if (associated(aero_props_obj)) then
87 : ! determine if there are transported aerosol contistuents
88 1536 : aerosols_transported = aero_props_obj%number_transported()>0
89 : end if
90 :
91 1536 : end subroutine stepon_init
92 :
93 : !=========================================================================================
94 :
95 32256 : subroutine stepon_run1( dtime_out, phys_state, phys_tend, &
96 : pbuf2d, dyn_in, dyn_out )
97 :
98 1536 : use time_manager, only: get_step_size
99 : use dp_coupling, only: d_p_coupling
100 : use physics_buffer, only: physics_buffer_desc
101 :
102 : use se_dyn_time_mod,only: tstep ! dynamics timestep
103 :
104 : real(r8), intent(out) :: dtime_out ! Time-step
105 : type(physics_state), intent(inout) :: phys_state(begchunk:endchunk)
106 : type(physics_tend), intent(inout) :: phys_tend(begchunk:endchunk)
107 : type (dyn_import_t), intent(inout) :: dyn_in ! Dynamics import container
108 : type (dyn_export_t), intent(inout) :: dyn_out ! Dynamics export container
109 : type (physics_buffer_desc), pointer :: pbuf2d(:,:)
110 : !----------------------------------------------------------------------------
111 :
112 : integer :: c
113 : class(aerosol_state), pointer :: aero_state_obj
114 16128 : nullify(aero_state_obj)
115 :
116 16128 : dtime_out = get_step_size()
117 :
118 16128 : if (iam < par%nprocs) then
119 16128 : if (tstep <= 0) call endrun('stepon_run1: bad tstep')
120 16128 : if (dtime_out <= 0) call endrun('stepon_run1: bad dtime')
121 :
122 : ! write diagnostic fields on gll grid and initial file
123 16128 : call diag_dynvar_ic(dyn_out%elem, dyn_out%fvm)
124 : end if
125 :
126 16128 : call t_barrierf('sync_d_p_coupling', mpicom)
127 16128 : call t_startf('d_p_coupling')
128 : ! Move data into phys_state structure.
129 16128 : call d_p_coupling(phys_state, phys_tend, pbuf2d, dyn_out )
130 16128 : call t_stopf('d_p_coupling')
131 :
132 : !----------------------------------------------------------
133 : ! update aerosol state object from CAM physics state constituents
134 : !----------------------------------------------------------
135 16128 : if (aerosols_transported) then
136 :
137 0 : do c = begchunk,endchunk
138 0 : aero_state_obj => aerosol_state_object(c)
139 : ! pass number mass or number mixing ratios of aerosol constituents
140 : ! to aerosol state object
141 0 : call aero_state_obj%set_transported(phys_state(c)%q)
142 : end do
143 :
144 : end if
145 :
146 16128 : end subroutine stepon_run1
147 :
148 : !=========================================================================================
149 :
150 29184 : subroutine stepon_run2(phys_state, phys_tend, dyn_in, dyn_out)
151 :
152 16128 : use dp_coupling, only: p_d_coupling
153 : use dyn_grid, only: TimeLevel
154 :
155 : use se_dyn_time_mod, only: TimeLevel_Qdp
156 : use control_mod, only: qsplit
157 : use prim_advance_mod, only: tot_energy_dyn
158 :
159 :
160 : ! arguments
161 : type(physics_state), intent(inout) :: phys_state(begchunk:endchunk)
162 : type(physics_tend), intent(inout) :: phys_tend(begchunk:endchunk)
163 : type (dyn_import_t), intent(inout) :: dyn_in ! Dynamics import container
164 : type (dyn_export_t), intent(inout) :: dyn_out ! Dynamics export container
165 :
166 : ! local variables
167 : integer :: tl_f, tl_fQdp
168 :
169 : integer :: c
170 : class(aerosol_state), pointer :: aero_state_obj
171 :
172 : !----------------------------------------------------------------------------
173 :
174 14592 : tl_f = TimeLevel%n0 ! timelevel which was adjusted by physics
175 14592 : call TimeLevel_Qdp(TimeLevel, qsplit, tl_fQdp)
176 :
177 : !----------------------------------------------------------
178 : ! update physics state with aerosol constituents
179 : !----------------------------------------------------------
180 14592 : nullify(aero_state_obj)
181 :
182 14592 : if (aerosols_transported) then
183 0 : do c = begchunk,endchunk
184 0 : aero_state_obj => aerosol_state_object(c)
185 : ! get mass or number mixing ratios of aerosol constituents
186 0 : call aero_state_obj%get_transported(phys_state(c)%q)
187 : end do
188 : end if
189 :
190 14592 : call t_barrierf('sync_p_d_coupling', mpicom)
191 14592 : call t_startf('p_d_coupling')
192 : ! copy from phys structures -> dynamics structures
193 14592 : call p_d_coupling(phys_state, phys_tend, dyn_in, tl_f, tl_fQdp)
194 14592 : call t_stopf('p_d_coupling')
195 :
196 14592 : if (iam < par%nprocs) then
197 14592 : call tot_energy_dyn(dyn_in%elem,dyn_in%fvm, 1, nelemd, tl_f, tl_fQdp,'dED')
198 : end if
199 :
200 14592 : end subroutine stepon_run2
201 :
202 : !=========================================================================================
203 :
204 14592 : subroutine stepon_run3(dtime, cam_out, phys_state, dyn_in, dyn_out)
205 :
206 14592 : use camsrfexch, only: cam_out_t
207 : use dyn_comp, only: dyn_run
208 : use advect_tend, only: compute_adv_tends_xyz
209 : use dyn_grid, only: TimeLevel
210 : use se_dyn_time_mod,only: TimeLevel_Qdp
211 : use control_mod, only: qsplit
212 : ! arguments
213 : real(r8), intent(in) :: dtime ! Time-step
214 : type(cam_out_t), intent(inout) :: cam_out(:) ! Output from CAM to surface
215 : type(physics_state), intent(inout) :: phys_state(begchunk:endchunk)
216 : type (dyn_import_t), intent(inout) :: dyn_in ! Dynamics import container
217 : type (dyn_export_t), intent(inout) :: dyn_out ! Dynamics export container
218 :
219 : integer :: tl_f, tl_fQdp
220 : !--------------------------------------------------------------------------------------
221 :
222 14592 : call t_startf('comp_adv_tends1')
223 14592 : tl_f = TimeLevel%n0
224 14592 : call TimeLevel_Qdp(TimeLevel, qsplit, tl_fQdp)
225 14592 : call compute_adv_tends_xyz(dyn_in%elem,dyn_in%fvm,1,nelemd,tl_fQdp,tl_f)
226 14592 : call t_stopf('comp_adv_tends1')
227 :
228 14592 : call t_barrierf('sync_dyn_run', mpicom)
229 14592 : call t_startf('dyn_run')
230 14592 : call dyn_run(dyn_out)
231 14592 : call t_stopf('dyn_run')
232 :
233 14592 : call t_startf('comp_adv_tends2')
234 14592 : tl_f = TimeLevel%n0
235 14592 : call TimeLevel_Qdp(TimeLevel, qsplit, tl_fQdp)
236 14592 : call compute_adv_tends_xyz(dyn_in%elem,dyn_in%fvm,1,nelemd,tl_fQdp,tl_f)
237 14592 : call t_stopf('comp_adv_tends2')
238 :
239 14592 : end subroutine stepon_run3
240 :
241 : !=========================================================================================
242 :
243 1536 : subroutine stepon_final(dyn_in, dyn_out)
244 :
245 : type (dyn_import_t), intent(inout) :: dyn_in ! Dynamics import container
246 : type (dyn_export_t), intent(inout) :: dyn_out ! Dynamics export container
247 :
248 14592 : end subroutine stepon_final
249 :
250 : !=========================================================================================
251 :
252 16128 : subroutine diag_dynvar_ic(elem, fvm)
253 : use constituents, only: cnst_type
254 : use cam_history, only: write_inithist, outfld, hist_fld_active, fieldname_len
255 : use dyn_grid, only: TimeLevel
256 :
257 : use se_dyn_time_mod, only: TimeLevel_Qdp ! dynamics typestep
258 : use control_mod, only: qsplit
259 : use hybrid_mod, only: config_thread_region, get_loop_ranges
260 : use hybrid_mod, only: hybrid_t
261 : use dimensions_mod, only: np, npsq, nc, nhc, fv_nphys, qsize, ntrac, nlev
262 : use dimensions_mod, only: cnst_name_gll
263 : use constituents, only: cnst_name
264 : use element_mod, only: element_t
265 : use fvm_control_volume_mod, only: fvm_struct
266 : use fvm_mapping, only: fvm2dyn
267 : use cam_thermo, only: get_sum_species, get_dp_ref, get_ps
268 : use air_composition, only: thermodynamic_active_species_idx
269 : use air_composition, only: thermodynamic_active_species_idx_dycore
270 : use hycoef, only: hyai, hybi, ps0
271 : ! arguments
272 : type(element_t) , intent(in) :: elem(1:nelemd)
273 : type(fvm_struct), intent(inout) :: fvm(:)
274 :
275 : ! local variables
276 : integer :: ie, i, j, k, m, m_cnst, nq
277 : integer :: tl_f, tl_qdp
278 : character(len=fieldname_len) :: tfname
279 :
280 : type(hybrid_t) :: hybrid
281 : integer :: nets, nete
282 16128 : real(r8), allocatable :: ftmp(:,:,:)
283 16128 : real(r8), allocatable :: fld_fvm(:,:,:,:,:), fld_gll(:,:,:,:,:)
284 16128 : real(r8), allocatable :: fld_2d(:,:)
285 : logical :: llimiter(1)
286 : real(r8) :: qtmp(np,np,nlev), dp_ref(np,np,nlev), ps_ref(np,np)
287 16128 : real(r8), allocatable :: factor_array(:,:,:)
288 : integer :: astat
289 : character(len=*), parameter :: prefix = 'diag_dynvar_ic: '
290 : !----------------------------------------------------------------------------
291 :
292 16128 : tl_f = timelevel%n0
293 16128 : call TimeLevel_Qdp(TimeLevel, qsplit, tl_Qdp)
294 :
295 16128 : allocate(ftmp(npsq,nlev,2))
296 :
297 : ! Output tracer fields for analysis of advection schemes
298 112896 : do m_cnst = 1, qsize
299 96768 : tfname = trim(cnst_name_gll(m_cnst))//'_gll'
300 112896 : if (hist_fld_active(tfname)) then
301 0 : do ie = 1, nelemd
302 0 : qtmp(:,:,:) = elem(ie)%state%Qdp(:,:,:,m_cnst,tl_qdp)/&
303 0 : elem(ie)%state%dp3d(:,:,:,tl_f)
304 0 : do j = 1, np
305 0 : do i = 1, np
306 0 : ftmp(i+(j-1)*np,:,1) = elem(ie)%state%Qdp(i,j,:,m_cnst,tl_qdp)/&
307 0 : elem(ie)%state%dp3d(i,j,:,tl_f)
308 : end do
309 : end do
310 0 : call outfld(tfname, ftmp(:,:,1), npsq, ie)
311 : end do
312 : end if
313 : end do
314 :
315 112896 : do m_cnst = 1, qsize
316 96768 : tfname = trim(cnst_name_gll(m_cnst))//'dp_gll'
317 112896 : if (hist_fld_active(tfname)) then
318 0 : do ie = 1, nelemd
319 0 : do j = 1, np
320 0 : do i = 1, np
321 0 : ftmp(i+(j-1)*np,:,1) = elem(ie)%state%Qdp(i,j,:,m_cnst,tl_qdp)
322 : end do
323 : end do
324 0 : call outfld(tfname, ftmp(:,:,1), npsq, ie)
325 : end do
326 : end if
327 : end do
328 :
329 16128 : if (hist_fld_active('U_gll') .or. hist_fld_active('V_gll')) then
330 0 : do ie = 1, nelemd
331 0 : do j = 1, np
332 0 : do i = 1, np
333 0 : ftmp(i+(j-1)*np,:,1) = elem(ie)%state%v(i,j,1,:,tl_f)
334 0 : ftmp(i+(j-1)*np,:,2) = elem(ie)%state%v(i,j,2,:,tl_f)
335 : end do
336 : end do
337 0 : call outfld('U_gll', ftmp(:,:,1), npsq, ie)
338 0 : call outfld('V_gll', ftmp(:,:,2), npsq, ie)
339 : end do
340 : end if
341 :
342 16128 : if (hist_fld_active('T_gll')) then
343 0 : do ie = 1, nelemd
344 0 : do j = 1, np
345 0 : do i = 1, np
346 0 : ftmp(i+(j-1)*np,:,1) = elem(ie)%state%T(i,j,:,tl_f)
347 : end do
348 : end do
349 0 : call outfld('T_gll', ftmp(:,:,1), npsq, ie)
350 : end do
351 : end if
352 :
353 16128 : if (hist_fld_active('dp_ref_gll')) then
354 0 : do ie = 1, nelemd
355 0 : call get_dp_ref(hyai, hybi, ps0, elem(ie)%state%phis(:,:), dp_ref(:,:,:), ps_ref(:,:))
356 0 : do j = 1, np
357 0 : do i = 1, np
358 0 : ftmp(i+(j-1)*np,:,1) = elem(ie)%state%dp3d(i,j,:,tl_f)/dp_ref(i,j,:)
359 : end do
360 : end do
361 0 : call outfld('dp_ref_gll', ftmp(:,:,1), npsq, ie)
362 : end do
363 : end if
364 :
365 16128 : if (hist_fld_active('PSDRY_gll')) then
366 0 : do ie = 1, nelemd
367 0 : do j = 1, np
368 0 : do i = 1, np
369 0 : ftmp(i+(j-1)*np,1,1) = elem(ie)%state%psdry(i,j)
370 : end do
371 : end do
372 0 : call outfld('PSDRY_gll', ftmp(:,1,1), npsq, ie)
373 : end do
374 : end if
375 :
376 16128 : if (hist_fld_active('PS_gll')) then
377 0 : allocate(fld_2d(np,np))
378 0 : do ie = 1, nelemd
379 0 : call get_ps(elem(ie)%state%Qdp(:,:,:,:,tl_Qdp), thermodynamic_active_species_idx_dycore,&
380 0 : elem(ie)%state%dp3d(:,:,:,tl_f),fld_2d,hyai(1)*ps0)
381 0 : do j = 1, np
382 0 : do i = 1, np
383 0 : ftmp(i+(j-1)*np,1,1) = fld_2d(i,j)
384 : end do
385 : end do
386 0 : call outfld('PS_gll', ftmp(:,1,1), npsq, ie)
387 : end do
388 0 : deallocate(fld_2d)
389 : end if
390 :
391 16128 : if (hist_fld_active('PHIS_gll')) then
392 0 : do ie = 1, nelemd
393 0 : call outfld('PHIS_gll', RESHAPE(elem(ie)%state%phis, (/np*np/)), np*np, ie)
394 : end do
395 : end if
396 :
397 16128 : if (write_inithist()) then
398 0 : allocate(fld_2d(np,np))
399 0 : do ie = 1, nelemd
400 0 : call get_ps(elem(ie)%state%Qdp(:,:,:,:,tl_Qdp), thermodynamic_active_species_idx_dycore,&
401 0 : elem(ie)%state%dp3d(:,:,:,tl_f),fld_2d,hyai(1)*ps0)
402 0 : do j = 1, np
403 0 : do i = 1, np
404 0 : ftmp(i+(j-1)*np,1,1) = fld_2d(i,j)
405 : end do
406 : end do
407 0 : call outfld('PS&IC', ftmp(:,1,1), npsq, ie)
408 : end do
409 0 : deallocate(fld_2d)
410 : endif
411 :
412 16128 : deallocate(ftmp)
413 :
414 16128 : if (write_inithist()) then
415 :
416 0 : if (fv_nphys < 1) then
417 0 : allocate(factor_array(np,np,nlev),stat=astat)
418 0 : if (astat /= 0) call endrun(prefix//"Allocate factor_array failed")
419 : endif
420 :
421 0 : do ie = 1, nelemd
422 0 : call outfld('T&IC', RESHAPE(elem(ie)%state%T(:,:,:,tl_f), (/npsq,nlev/)), npsq, ie)
423 0 : call outfld('U&IC', RESHAPE(elem(ie)%state%v(:,:,1,:,tl_f), (/npsq,nlev/)), npsq, ie)
424 0 : call outfld('V&IC', RESHAPE(elem(ie)%state%v(:,:,2,:,tl_f), (/npsq,nlev/)), npsq, ie)
425 :
426 0 : if (fv_nphys < 1) then
427 0 : call get_sum_species(elem(ie)%state%Qdp(:,:,:,:,tl_qdp), &
428 0 : thermodynamic_active_species_idx_dycore, factor_array,dp_dry=elem(ie)%state%dp3d(:,:,:,tl_f))
429 0 : factor_array(:,:,:) = 1.0_r8/factor_array(:,:,:)
430 0 : do m_cnst = 1, qsize
431 0 : if (cnst_type(m_cnst) == 'wet') then
432 : call outfld(trim(cnst_name(m_cnst))//'&IC', &
433 0 : RESHAPE(factor_array(:,:,:)*elem(ie)%state%Qdp(:,:,:,m_cnst,tl_qdp)/&
434 0 : elem(ie)%state%dp3d(:,:,:,tl_f), (/npsq,nlev/)), npsq, ie)
435 : else
436 : call outfld(trim(cnst_name(m_cnst))//'&IC', &
437 0 : RESHAPE(elem(ie)%state%Qdp(:,:,:,m_cnst,tl_qdp)/&
438 0 : elem(ie)%state%dp3d(:,:,:,tl_f), (/npsq,nlev/)), npsq, ie)
439 : end if
440 : end do
441 : end if
442 : end do
443 :
444 0 : if (fv_nphys > 0) then
445 : !JMD $OMP PARALLEL NUM_THREADS(horz_num_threads), DEFAULT(SHARED), PRIVATE(hybrid,nets,nete,n)
446 : !JMD hybrid = config_thread_region(par,'horizontal')
447 0 : hybrid = config_thread_region(par,'serial')
448 0 : call get_loop_ranges(hybrid, ibeg=nets, iend=nete)
449 :
450 0 : allocate(fld_fvm(1-nhc:nc+nhc,1-nhc:nc+nhc,nlev,1,nets:nete),stat=astat)
451 0 : if (astat /= 0) call endrun(prefix//"Allocate fld_fvm failed")
452 0 : allocate(fld_gll(np,np,nlev,1,nets:nete),stat=astat)
453 0 : if (astat /= 0) call endrun(prefix//"Allocate fld_gll failed")
454 0 : allocate(factor_array(nc,nc,nlev),stat=astat)
455 0 : if (astat /= 0) call endrun(prefix//"Allocate factor_array failed")
456 :
457 0 : llimiter = .true.
458 :
459 0 : do m_cnst = 1, ntrac
460 0 : do ie = nets, nete
461 :
462 0 : call get_sum_species(fvm(ie)%c(1:nc,1:nc,:,:),thermodynamic_active_species_idx,factor_array)
463 0 : factor_array(:,:,:) = 1.0_r8/factor_array(:,:,:)
464 :
465 0 : if (cnst_type(m_cnst) == 'wet') then
466 0 : fld_fvm(1:nc,1:nc,:,1,ie) = fvm(ie)%c(1:nc,1:nc,:,m_cnst)*factor_array(:,:,:)
467 : else
468 0 : fld_fvm(1:nc,1:nc,:,1,ie) = fvm(ie)%c(1:nc,1:nc,:,m_cnst)
469 : end if
470 : end do
471 :
472 0 : call fvm2dyn(fld_fvm, fld_gll, hybrid, nets, nete, nlev, fvm(nets:nete), llimiter)
473 :
474 0 : do ie = nets, nete
475 0 : call outfld(trim(cnst_name(m_cnst))//'&IC', &
476 0 : RESHAPE(fld_gll(:,:,:,:,ie), (/npsq,nlev/)), npsq, ie)
477 : end do
478 : end do
479 :
480 0 : deallocate(fld_fvm)
481 0 : deallocate(fld_gll)
482 : end if
483 :
484 0 : deallocate(factor_array)
485 :
486 : end if ! if (write_inithist)
487 :
488 16128 : end subroutine diag_dynvar_ic
489 :
490 : !=========================================================================================
491 :
492 : end module stepon
|