Line data Source code
1 : module restart_physics
2 :
3 : use shr_kind_mod, only: r8 => shr_kind_r8
4 : use spmd_utils, only: masterproc
5 : use co2_cycle, only: co2_transport
6 : use constituents, only: pcnst
7 :
8 : use radiation, only: radiation_define_restart, radiation_write_restart, &
9 : radiation_read_restart
10 :
11 : use ioFileMod
12 : use cam_abortutils, only: endrun
13 : use camsrfexch, only: cam_in_t, cam_out_t
14 : use cam_logfile, only: iulog
15 : use pio, only: file_desc_t, io_desc_t, var_desc_t, &
16 : pio_double, pio_int, pio_noerr, &
17 : pio_seterrorhandling, pio_bcast_error, &
18 : pio_inq_varid, &
19 : pio_def_var, pio_def_dim, &
20 : pio_put_var, pio_get_var
21 :
22 : implicit none
23 : private
24 : save
25 : !
26 : ! Public interfaces
27 : !
28 : public :: write_restart_physics ! Write the physics restart info out
29 : public :: read_restart_physics ! Read the physics restart info in
30 : public :: init_restart_physics
31 :
32 : !
33 : ! Private data
34 : !
35 :
36 : type(var_desc_t) :: flwds_desc, &
37 : solld_desc, co2prog_desc, co2diag_desc, sols_desc, soll_desc, &
38 : solsd_desc
39 :
40 : type(var_desc_t) :: bcphidry_desc, bcphodry_desc, ocphidry_desc, ocphodry_desc, &
41 : dstdry1_desc, dstdry2_desc, dstdry3_desc, dstdry4_desc
42 :
43 : type(var_desc_t) :: cflx_desc(pcnst)
44 :
45 : type(var_desc_t) :: wsx_desc
46 : type(var_desc_t) :: wsy_desc
47 : type(var_desc_t) :: shf_desc
48 :
49 : CONTAINS
50 2304 : subroutine init_restart_physics ( File, pbuf2d)
51 :
52 : use physics_buffer, only: pbuf_init_restart, physics_buffer_desc
53 : use ppgrid, only: pver, pverp
54 : use chemistry, only: chem_init_restart
55 : use prescribed_ozone, only: init_prescribed_ozone_restart
56 : use prescribed_ghg, only: init_prescribed_ghg_restart
57 : use prescribed_aero, only: init_prescribed_aero_restart
58 : use prescribed_volcaero, only: init_prescribed_volcaero_restart
59 : use cam_grid_support, only: cam_grid_write_attr, cam_grid_id
60 : use cam_grid_support, only: cam_grid_header_info_t
61 : use cam_pio_utils, only: cam_pio_def_dim
62 : use subcol_utils, only: is_subcol_on
63 : use subcol, only: subcol_init_restart
64 : use carma_intr, only: carma_restart_init
65 :
66 : type(file_desc_t), intent(inout) :: file
67 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
68 :
69 : integer :: grid_id
70 : integer :: hdimcnt, ierr, i
71 : integer :: dimids(4)
72 1536 : integer, allocatable :: hdimids(:)
73 1536 : type(cam_grid_header_info_t) :: info
74 : character(len=4) :: num
75 :
76 1536 : call pio_seterrorhandling(File, PIO_BCAST_ERROR)
77 : ! Probably should have the grid write this out.
78 1536 : grid_id = cam_grid_id('physgrid')
79 1536 : call cam_grid_write_attr(File, grid_id, info)
80 1536 : hdimcnt = info%num_hdims()
81 :
82 4608 : do i = 1, hdimcnt
83 4608 : dimids(i) = info%get_hdimid(i)
84 : end do
85 4608 : allocate(hdimids(hdimcnt))
86 4608 : hdimids(1:hdimcnt) = dimids(1:hdimcnt)
87 :
88 1536 : call pbuf_init_restart(File, pbuf2d)
89 :
90 1536 : call chem_init_restart(File)
91 :
92 1536 : call init_prescribed_ozone_restart(File)
93 1536 : call init_prescribed_ghg_restart(File)
94 1536 : call init_prescribed_aero_restart(File)
95 1536 : call init_prescribed_volcaero_restart(File)
96 :
97 1536 : ierr = pio_def_var(File, 'FLWDS', pio_double, hdimids, flwds_desc)
98 1536 : ierr = pio_def_var(File, 'SOLS', pio_double, hdimids, sols_desc)
99 1536 : ierr = pio_def_var(File, 'SOLL', pio_double, hdimids, soll_desc)
100 1536 : ierr = pio_def_var(File, 'SOLSD', pio_double, hdimids, solsd_desc)
101 1536 : ierr = pio_def_var(File, 'SOLLD', pio_double, hdimids, solld_desc)
102 :
103 1536 : ierr = pio_def_var(File, 'BCPHIDRY', pio_double, hdimids, bcphidry_desc)
104 1536 : ierr = pio_def_var(File, 'BCPHODRY', pio_double, hdimids, bcphodry_desc)
105 1536 : ierr = pio_def_var(File, 'OCPHIDRY', pio_double, hdimids, ocphidry_desc)
106 1536 : ierr = pio_def_var(File, 'OCPHODRY', pio_double, hdimids, ocphodry_desc)
107 1536 : ierr = pio_def_var(File, 'DSTDRY1', pio_double, hdimids, dstdry1_desc)
108 1536 : ierr = pio_def_var(File, 'DSTDRY2', pio_double, hdimids, dstdry2_desc)
109 1536 : ierr = pio_def_var(File, 'DSTDRY3', pio_double, hdimids, dstdry3_desc)
110 1536 : ierr = pio_def_var(File, 'DSTDRY4', pio_double, hdimids, dstdry4_desc)
111 :
112 1536 : if(co2_transport()) then
113 0 : ierr = pio_def_var(File, 'CO2PROG', pio_double, hdimids, co2prog_desc)
114 0 : ierr = pio_def_var(File, 'CO2DIAG', pio_double, hdimids, co2diag_desc)
115 : end if
116 :
117 : ! cam_import variables -- write the constituent surface fluxes as individual 2D arrays
118 : ! rather than as a single variable with a pcnst dimension. Note that the cflx components
119 : ! are only needed for those constituents that are not passed to the coupler. The restart
120 : ! for constituents passed through the coupler are handled by the .rs. restart file. But
121 : ! we don't currently have a mechanism to know whether the constituent is handled by the
122 : ! coupler or not, so we write all of cflx to the CAM restart file.
123 370176 : do i = 1, pcnst
124 368640 : write(num,'(i4.4)') i
125 370176 : ierr = pio_def_var(File, 'CFLX'//num, pio_double, hdimids, cflx_desc(i))
126 : end do
127 :
128 1536 : ierr = pio_def_var(File, 'wsx', pio_double, hdimids, wsx_desc)
129 1536 : ierr = pio_def_var(File, 'wsy', pio_double, hdimids, wsy_desc)
130 1536 : ierr = pio_def_var(File, 'shf', pio_double, hdimids, shf_desc)
131 :
132 1536 : call radiation_define_restart(file)
133 :
134 1536 : if (is_subcol_on()) then
135 0 : call subcol_init_restart(file, hdimids)
136 : end if
137 :
138 1536 : call carma_restart_init(file)
139 :
140 1536 : end subroutine init_restart_physics
141 :
142 3072 : subroutine write_restart_physics (File, cam_in, cam_out, pbuf2d)
143 :
144 : !-----------------------------------------------------------------------
145 1536 : use physics_buffer, only: physics_buffer_desc, pbuf_write_restart
146 : use phys_grid, only: phys_decomp
147 :
148 : use ppgrid, only: begchunk, endchunk, pcols, pverp
149 : use chemistry, only: chem_write_restart
150 : use prescribed_ozone, only: write_prescribed_ozone_restart
151 : use prescribed_ghg, only: write_prescribed_ghg_restart
152 : use prescribed_aero, only: write_prescribed_aero_restart
153 : use prescribed_volcaero, only: write_prescribed_volcaero_restart
154 :
155 : use cam_history_support, only: fillvalue
156 : use spmd_utils, only: iam
157 : use cam_grid_support, only: cam_grid_write_dist_array, cam_grid_id
158 : use cam_grid_support, only: cam_grid_get_decomp, cam_grid_dimensions
159 : use cam_grid_support, only: cam_grid_write_var
160 : use pio, only: pio_write_darray
161 : use subcol_utils, only: is_subcol_on
162 : use subcol, only: subcol_write_restart
163 : use carma_intr, only: carma_restart_write
164 : !
165 : ! Input arguments
166 : !
167 : type(file_desc_t), intent(inout) :: File
168 : type(cam_in_t), intent(in) :: cam_in(begchunk:endchunk)
169 : type(cam_out_t), intent(in) :: cam_out(begchunk:endchunk)
170 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
171 : !
172 : ! Local workspace
173 : !
174 : type(io_desc_t), pointer :: iodesc
175 3072 : real(r8):: tmpfield(pcols, begchunk:endchunk)
176 : integer :: i, m ! loop index
177 : integer :: ncol ! number of vertical columns
178 : integer :: ierr
179 : integer :: physgrid
180 : integer :: dims(3), gdims(3)
181 : integer :: nhdims
182 : !-----------------------------------------------------------------------
183 :
184 : ! Write grid vars
185 1536 : call cam_grid_write_var(File, phys_decomp)
186 :
187 : ! Physics buffer
188 1536 : if (is_subcol_on()) then
189 0 : call subcol_write_restart(File)
190 : end if
191 :
192 1536 : call pbuf_write_restart(File, pbuf2d)
193 :
194 1536 : physgrid = cam_grid_id('physgrid')
195 1536 : call cam_grid_dimensions(physgrid, gdims(1:2), nhdims)
196 :
197 : ! data for chemistry
198 1536 : call chem_write_restart(File)
199 :
200 1536 : call write_prescribed_ozone_restart(File)
201 1536 : call write_prescribed_ghg_restart(File)
202 1536 : call write_prescribed_aero_restart(File)
203 1536 : call write_prescribed_volcaero_restart(File)
204 :
205 : ! cam_in/out variables
206 : ! This is a group of surface variables so can reuse dims
207 1536 : dims(1) = pcols
208 1536 : dims(2) = endchunk - begchunk + 1
209 0 : call cam_grid_get_decomp(physgrid, dims(1:2), gdims(1:nhdims), &
210 1536 : pio_double, iodesc)
211 :
212 9216 : do i = begchunk, endchunk
213 7680 : ncol = cam_out(i)%ncol
214 118272 : tmpfield(:ncol, i) = cam_out(i)%flwds(:ncol)
215 : ! Only have to do this once (cam_in/out vars all same shape)
216 9216 : if (ncol < pcols) then
217 19968 : tmpfield(ncol+1:, i) = fillvalue
218 : end if
219 : end do
220 1536 : call pio_write_darray(File, flwds_desc, iodesc, tmpfield, ierr)
221 :
222 9216 : do i = begchunk, endchunk
223 7680 : ncol = cam_out(i)%ncol
224 119808 : tmpfield(:ncol, i) = cam_out(i)%sols(:ncol)
225 : end do
226 1536 : call pio_write_darray(File, sols_desc, iodesc, tmpfield, ierr)
227 :
228 9216 : do i = begchunk, endchunk
229 7680 : ncol = cam_out(i)%ncol
230 119808 : tmpfield(:ncol, i) = cam_out(i)%soll(:ncol)
231 : end do
232 1536 : call pio_write_darray(File, soll_desc, iodesc, tmpfield, ierr)
233 :
234 9216 : do i = begchunk, endchunk
235 7680 : ncol = cam_out(i)%ncol
236 119808 : tmpfield(:ncol, i) = cam_out(i)%solsd(:ncol)
237 : end do
238 1536 : call pio_write_darray(File, solsd_desc, iodesc, tmpfield, ierr)
239 :
240 9216 : do i = begchunk, endchunk
241 7680 : ncol = cam_out(i)%ncol
242 119808 : tmpfield(:ncol, i) = cam_out(i)%solld(:ncol)
243 : end do
244 1536 : call pio_write_darray(File, solld_desc, iodesc, tmpfield, ierr)
245 :
246 9216 : do i = begchunk, endchunk
247 7680 : ncol = cam_out(i)%ncol
248 119808 : tmpfield(:ncol, i) = cam_out(i)%bcphidry(:ncol)
249 : end do
250 1536 : call pio_write_darray(File, bcphidry_desc, iodesc, tmpfield, ierr)
251 :
252 9216 : do i = begchunk, endchunk
253 7680 : ncol = cam_out(i)%ncol
254 119808 : tmpfield(:ncol, i) = cam_out(i)%bcphodry(:ncol)
255 : end do
256 1536 : call pio_write_darray(File, bcphodry_desc, iodesc, tmpfield, ierr)
257 :
258 9216 : do i = begchunk, endchunk
259 7680 : ncol = cam_out(i)%ncol
260 119808 : tmpfield(:ncol, i) = cam_out(i)%ocphidry(:ncol)
261 : end do
262 1536 : call pio_write_darray(File, ocphidry_desc, iodesc, tmpfield, ierr)
263 :
264 9216 : do i = begchunk, endchunk
265 7680 : ncol = cam_out(i)%ncol
266 119808 : tmpfield(:ncol, i) = cam_out(i)%ocphodry(:ncol)
267 : end do
268 1536 : call pio_write_darray(File, ocphodry_desc, iodesc, tmpfield, ierr)
269 :
270 9216 : do i = begchunk, endchunk
271 7680 : ncol = cam_out(i)%ncol
272 119808 : tmpfield(:ncol, i) = cam_out(i)%dstdry1(:ncol)
273 : end do
274 1536 : call pio_write_darray(File, dstdry1_desc, iodesc, tmpfield, ierr)
275 :
276 9216 : do i = begchunk, endchunk
277 7680 : ncol = cam_out(i)%ncol
278 119808 : tmpfield(:ncol, i) = cam_out(i)%dstdry2(:ncol)
279 : end do
280 1536 : call pio_write_darray(File, dstdry2_desc, iodesc, tmpfield, ierr)
281 :
282 9216 : do i = begchunk, endchunk
283 7680 : ncol = cam_out(i)%ncol
284 119808 : tmpfield(:ncol, i) = cam_out(i)%dstdry3(:ncol)
285 : end do
286 1536 : call pio_write_darray(File, dstdry3_desc, iodesc, tmpfield, ierr)
287 :
288 9216 : do i = begchunk, endchunk
289 7680 : ncol = cam_out(i)%ncol
290 119808 : tmpfield(:ncol, i) = cam_out(i)%dstdry4(:ncol)
291 : end do
292 1536 : call pio_write_darray(File, dstdry4_desc, iodesc, tmpfield, ierr)
293 :
294 1536 : if (co2_transport()) then
295 0 : do i = begchunk, endchunk
296 0 : ncol = cam_out(i)%ncol
297 0 : tmpfield(:ncol, i) = cam_out(i)%co2prog(:ncol)
298 : end do
299 0 : call pio_write_darray(File, co2prog_desc, iodesc, tmpfield, ierr)
300 :
301 0 : do i = begchunk, endchunk
302 0 : ncol = cam_out(i)%ncol
303 0 : tmpfield(:ncol, i) = cam_out(i)%co2diag(:ncol)
304 : end do
305 0 : call pio_write_darray(File, co2diag_desc, iodesc, tmpfield, ierr)
306 : end if
307 :
308 : ! cam_in components
309 370176 : do m = 1, pcnst
310 2211840 : do i = begchunk, endchunk
311 1843200 : ncol = cam_in(i)%ncol
312 28753920 : tmpfield(:ncol, i) = cam_in(i)%cflx(:ncol, m)
313 : end do
314 370176 : call pio_write_darray(File, cflx_desc(m), iodesc, tmpfield, ierr)
315 : end do
316 :
317 9216 : do i = begchunk, endchunk
318 7680 : ncol = cam_in(i)%ncol
319 119808 : tmpfield(:ncol,i) = cam_in(i)%wsx(:ncol)
320 : end do
321 1536 : call pio_write_darray(File, wsx_desc, iodesc, tmpfield, ierr)
322 :
323 9216 : do i = begchunk, endchunk
324 7680 : ncol = cam_in(i)%ncol
325 119808 : tmpfield(:ncol,i) = cam_in(i)%wsy(:ncol)
326 : end do
327 1536 : call pio_write_darray(File, wsy_desc, iodesc, tmpfield, ierr)
328 :
329 9216 : do i = begchunk, endchunk
330 7680 : ncol = cam_in(i)%ncol
331 119808 : tmpfield(:ncol,i) = cam_in(i)%shf(:ncol)
332 : end do
333 1536 : call pio_write_darray(File, shf_desc, iodesc, tmpfield, ierr)
334 :
335 1536 : call radiation_write_restart(file)
336 1536 : call carma_restart_write(file)
337 :
338 1536 : end subroutine write_restart_physics
339 :
340 : !#######################################################################
341 :
342 768 : subroutine read_restart_physics(File, cam_in, cam_out, pbuf2d)
343 :
344 : !-----------------------------------------------------------------------
345 1536 : use physics_buffer, only: physics_buffer_desc, pbuf_read_restart
346 :
347 : use ppgrid, only: begchunk, endchunk, pcols, pver, pverp
348 : use chemistry, only: chem_read_restart
349 : use cam_grid_support, only: cam_grid_read_dist_array, cam_grid_id
350 : use cam_grid_support, only: cam_grid_get_decomp, cam_grid_dimensions
351 : use cam_history_support, only: fillvalue
352 :
353 : use prescribed_ozone, only: read_prescribed_ozone_restart
354 : use prescribed_ghg, only: read_prescribed_ghg_restart
355 : use prescribed_aero, only: read_prescribed_aero_restart
356 : use prescribed_volcaero, only: read_prescribed_volcaero_restart
357 : use subcol_utils, only: is_subcol_on
358 : use subcol, only: subcol_read_restart
359 : use pio, only: pio_read_darray
360 : use carma_intr, only: carma_restart_read
361 : !
362 : ! Arguments
363 : !
364 : type(file_desc_t), intent(inout) :: File
365 : type(cam_in_t), pointer :: cam_in(:)
366 : type(cam_out_t), pointer :: cam_out(:)
367 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
368 : !
369 : ! Local workspace
370 : !
371 768 : real(r8), allocatable :: tmpfield2(:,:)
372 : integer :: i, c, m ! loop index
373 : integer :: ierr ! I/O status
374 : type(io_desc_t), pointer :: iodesc
375 : type(var_desc_t) :: vardesc
376 : integer :: csize, vsize
377 : character(len=4) :: num
378 : integer :: dims(3), gdims(3), nhdims
379 : integer :: err_handling
380 : integer :: physgrid
381 : !-----------------------------------------------------------------------
382 :
383 : ! subcol_read_restart must be called before pbuf_read_restart
384 768 : if (is_subcol_on()) then
385 0 : call subcol_read_restart(File)
386 : end if
387 :
388 768 : call pbuf_read_restart(File, pbuf2d)
389 :
390 768 : csize=endchunk-begchunk+1
391 768 : dims(1) = pcols
392 768 : dims(2) = csize
393 :
394 768 : physgrid = cam_grid_id('physgrid')
395 :
396 768 : call cam_grid_dimensions(physgrid, gdims(1:2))
397 :
398 768 : if (gdims(2) == 1) then
399 : nhdims = 1
400 : else
401 768 : nhdims = 2
402 : end if
403 0 : call cam_grid_get_decomp(physgrid, dims(1:2), gdims(1:nhdims), pio_double, &
404 768 : iodesc)
405 :
406 : ! data for chemistry
407 768 : call chem_read_restart(File)
408 :
409 768 : call read_prescribed_ozone_restart(File)
410 768 : call read_prescribed_ghg_restart(File)
411 768 : call read_prescribed_aero_restart(File)
412 768 : call read_prescribed_volcaero_restart(File)
413 :
414 2304 : allocate(tmpfield2(pcols, begchunk:endchunk))
415 66048 : tmpfield2 = fillvalue
416 :
417 768 : ierr = pio_inq_varid(File, 'FLWDS', vardesc)
418 768 : call pio_read_darray(File, vardesc, iodesc, tmpfield2, ierr)
419 4608 : do c=begchunk,endchunk
420 66048 : do i=1,pcols
421 65280 : cam_out(c)%flwds(i) = tmpfield2(i, c)
422 : end do
423 : end do
424 :
425 768 : ierr = pio_inq_varid(File, 'SOLS', vardesc)
426 768 : call pio_read_darray(File, vardesc, iodesc, tmpfield2, ierr)
427 4608 : do c=begchunk,endchunk
428 66048 : do i=1,pcols
429 65280 : cam_out(c)%sols(i) = tmpfield2(i, c)
430 : end do
431 : end do
432 :
433 768 : ierr = pio_inq_varid(File, 'SOLL', vardesc)
434 768 : call pio_read_darray(File, vardesc, iodesc, tmpfield2, ierr)
435 4608 : do c=begchunk,endchunk
436 66048 : do i=1,pcols
437 65280 : cam_out(c)%soll(i) = tmpfield2(i, c)
438 : end do
439 : end do
440 :
441 768 : ierr = pio_inq_varid(File, 'SOLSD', vardesc)
442 768 : call pio_read_darray(File, vardesc, iodesc, tmpfield2, ierr)
443 4608 : do c=begchunk,endchunk
444 66048 : do i=1,pcols
445 65280 : cam_out(c)%solsd(i) = tmpfield2(i, c)
446 : end do
447 : end do
448 :
449 768 : ierr = pio_inq_varid(File, 'SOLLD', vardesc)
450 768 : call pio_read_darray(File, vardesc, iodesc, tmpfield2, ierr)
451 4608 : do c=begchunk,endchunk
452 66048 : do i=1,pcols
453 65280 : cam_out(c)%solld(i) = tmpfield2(i, c)
454 : end do
455 : end do
456 :
457 768 : ierr = pio_inq_varid(File, 'BCPHIDRY', vardesc)
458 768 : call pio_read_darray(File, vardesc, iodesc, tmpfield2, ierr)
459 4608 : do c=begchunk,endchunk
460 66048 : do i=1,pcols
461 65280 : cam_out(c)%bcphidry(i) = tmpfield2(i, c)
462 : end do
463 : end do
464 :
465 768 : ierr = pio_inq_varid(File, 'BCPHODRY', vardesc)
466 768 : call pio_read_darray(File, vardesc, iodesc, tmpfield2, ierr)
467 4608 : do c=begchunk,endchunk
468 66048 : do i=1,pcols
469 65280 : cam_out(c)%bcphodry(i) = tmpfield2(i, c)
470 : end do
471 : end do
472 :
473 768 : ierr = pio_inq_varid(File, 'OCPHIDRY', vardesc)
474 768 : call pio_read_darray(File, vardesc, iodesc, tmpfield2, ierr)
475 4608 : do c=begchunk,endchunk
476 66048 : do i=1,pcols
477 65280 : cam_out(c)%ocphidry(i) = tmpfield2(i, c)
478 : end do
479 : end do
480 :
481 768 : ierr = pio_inq_varid(File, 'OCPHODRY', vardesc)
482 768 : call pio_read_darray(File, vardesc, iodesc, tmpfield2, ierr)
483 4608 : do c=begchunk,endchunk
484 66048 : do i=1,pcols
485 65280 : cam_out(c)%ocphodry(i) = tmpfield2(i, c)
486 : end do
487 : end do
488 :
489 768 : ierr = pio_inq_varid(File, 'DSTDRY1', vardesc)
490 768 : call pio_read_darray(File, vardesc, iodesc, tmpfield2, ierr)
491 4608 : do c=begchunk,endchunk
492 66048 : do i=1,pcols
493 65280 : cam_out(c)%dstdry1(i) = tmpfield2(i, c)
494 : end do
495 : end do
496 :
497 768 : ierr = pio_inq_varid(File, 'DSTDRY2', vardesc)
498 768 : call pio_read_darray(File, vardesc, iodesc, tmpfield2, ierr)
499 4608 : do c=begchunk,endchunk
500 66048 : do i=1,pcols
501 65280 : cam_out(c)%dstdry2(i) = tmpfield2(i, c)
502 : end do
503 : end do
504 :
505 768 : ierr = pio_inq_varid(File, 'DSTDRY3', vardesc)
506 768 : call pio_read_darray(File, vardesc, iodesc, tmpfield2, ierr)
507 4608 : do c=begchunk,endchunk
508 66048 : do i=1,pcols
509 65280 : cam_out(c)%dstdry3(i) = tmpfield2(i, c)
510 : end do
511 : end do
512 :
513 768 : ierr = pio_inq_varid(File, 'DSTDRY4', vardesc)
514 768 : call pio_read_darray(File, vardesc, iodesc, tmpfield2, ierr)
515 4608 : do c=begchunk,endchunk
516 66048 : do i=1,pcols
517 65280 : cam_out(c)%dstdry4(i) = tmpfield2(i, c)
518 : end do
519 : end do
520 :
521 768 : if (co2_transport()) then
522 0 : ierr = pio_inq_varid(File, 'CO2PROG', vardesc)
523 0 : call pio_read_darray(File, vardesc, iodesc, tmpfield2, ierr)
524 0 : do c=begchunk,endchunk
525 0 : do i=1,pcols
526 0 : cam_out(c)%co2prog(i) = tmpfield2(i, c)
527 : end do
528 : end do
529 :
530 0 : ierr = pio_inq_varid(File, 'CO2DIAG', vardesc)
531 0 : call pio_read_darray(File, vardesc, iodesc, tmpfield2, ierr)
532 0 : do c=begchunk,endchunk
533 0 : do i=1,pcols
534 0 : cam_out(c)%co2diag(i) = tmpfield2(i, c)
535 : end do
536 : end do
537 : end if
538 :
539 : ! Reading the CFLX* components from the restart is optional for
540 : ! backwards compatibility. These fields were not needed for an
541 : ! exact restart until the UNICON scheme was added. More generally,
542 : ! these components are only needed if they are not handled by the
543 : ! coupling layer restart (the ".rs." file), and if the values are
544 : ! used in the tphysbc physics before the tphysac code has a chance
545 : ! to update the values that are coming from boundary datasets.
546 185088 : do m = 1, pcnst
547 :
548 184320 : write(num,'(i4.4)') m
549 :
550 184320 : call pio_seterrorhandling(File, PIO_BCAST_ERROR, err_handling)
551 184320 : ierr = pio_inq_varid(File, 'CFLX'//num, vardesc)
552 184320 : call pio_seterrorhandling(File, err_handling)
553 :
554 185088 : if (ierr == PIO_NOERR) then ! CFLX variable found on restart file
555 184320 : call pio_read_darray(File, vardesc, iodesc, tmpfield2, ierr)
556 1105920 : do c= begchunk, endchunk
557 15851520 : do i = 1, pcols
558 15667200 : cam_in(c)%cflx(i,m) = tmpfield2(i, c)
559 : end do
560 : end do
561 : end if
562 :
563 : end do
564 :
565 768 : call pio_seterrorhandling(File, PIO_BCAST_ERROR, err_handling)
566 768 : ierr = pio_inq_varid(File, 'wsx', vardesc)
567 768 : if (ierr == PIO_NOERR) then ! variable found on restart file
568 768 : call pio_read_darray(File, vardesc, iodesc, tmpfield2, ierr)
569 4608 : do c= begchunk, endchunk
570 66048 : do i = 1, pcols
571 65280 : cam_in(c)%wsx(i) = tmpfield2(i, c)
572 : end do
573 : end do
574 : end if
575 768 : ierr = pio_inq_varid(File, 'wsy', vardesc)
576 768 : if (ierr == PIO_NOERR) then ! variable found on restart file
577 768 : call pio_read_darray(File, vardesc, iodesc, tmpfield2, ierr)
578 4608 : do c= begchunk, endchunk
579 66048 : do i = 1, pcols
580 65280 : cam_in(c)%wsy(i) = tmpfield2(i, c)
581 : end do
582 : end do
583 : end if
584 768 : ierr = pio_inq_varid(File, 'shf', vardesc)
585 768 : if (ierr == PIO_NOERR) then ! variable found on restart file
586 768 : call pio_read_darray(File, vardesc, iodesc, tmpfield2, ierr)
587 4608 : do c= begchunk, endchunk
588 66048 : do i = 1, pcols
589 65280 : cam_in(c)%shf(i) = tmpfield2(i, c)
590 : end do
591 : end do
592 : endif
593 768 : call pio_seterrorhandling(File, err_handling)
594 :
595 768 : deallocate(tmpfield2)
596 :
597 768 : call radiation_read_restart(file)
598 768 : call carma_restart_read(file)
599 :
600 10752 : end subroutine read_restart_physics
601 :
602 : end module restart_physics
|