global_cycle 1.14.0
Loading...
Searching...
No Matches
land_increments.F90
Go to the documentation of this file.
1
4
5module land_increments
6
7 private
8
9 public gaussian_to_fv3_interp
10 public add_increment_soil
11 public add_increment_snow
12 public calculate_landinc_mask
13 public apply_land_da_adjustments_soil
14 public apply_land_da_adjustments_snd
15 public lsm_noah, lsm_noahmp
16
17 integer, parameter :: lsm_noah=1
18 integer, parameter :: lsm_noahmp=2
20
21 ! control state for soil analysis:
22
23 real, parameter :: tfreez=273.16
24contains
25
45
46subroutine gaussian_to_fv3_interp(lsoil_incr,rla,rlo, &
47 stcinc,slcinc,soilsnow_tile,lensfc,lsoil,idim,jdim,lsm, myrank)
48
49 use utils
50 use gdswzd_mod
51 use read_write_data, only : idim_gaus, jdim_gaus, &
53 use mpi
54
55 implicit none
56
57 integer, intent(in) :: lsoil_incr, lensfc, lsoil, idim, jdim, myrank, lsm
58
59 integer, intent(in) :: soilsnow_tile(lensfc)
60 real, intent(inout) :: rla(lensfc), rlo(lensfc)
61
62 real, intent(out) :: stcinc(lensfc,lsoil)
63 real, intent(out) :: slcinc(lensfc,lsoil)
64
65 integer :: iopt, nret, kgds_gaus(200)
66 integer :: igaus, jgaus, ij
67 integer :: mask_tile
68 integer :: itile, jtile
69 integer :: j, ierr
70 integer :: igausp1, jgausp1
71 logical :: upd_slc, upd_stc
72 real :: fill
73
74 integer, allocatable :: id1(:,:), id2(:,:), jdc(:,:)
75
76 real :: wsum
77 real :: stc_inc(lsoil)
78 real :: slc_inc(lsoil)
79 real, allocatable :: xpts(:), ypts(:), lats(:), lons(:)
80 real, allocatable :: dum2d(:,:), lats_rad(:), lons_rad(:)
81 real, allocatable :: agrid(:,:,:), s2c(:,:,:)
82
83 integer :: k
84 logical :: gaus_has_soil
85
86 ! this produces the same lat/lon as can be read in from the file
87
88 kgds_gaus = 0
89 kgds_gaus(1) = 4 ! oct 6 - type of grid (gaussian)
90 kgds_gaus(2) = idim_gaus ! oct 7-8 - # pts on latitude circle
91 kgds_gaus(3) = jdim_gaus
92 kgds_gaus(4) = 90000 ! oct 11-13 - lat of origin
93 kgds_gaus(5) = 0 ! oct 14-16 - lon of origin
94 kgds_gaus(6) = 128 ! oct 17 - resolution flag
95 kgds_gaus(7) = -90000 ! oct 18-20 - lat of extreme point
96 kgds_gaus(8) = nint(-360000./float(idim_gaus)) ! oct 21-23 - lon of extreme point
97 kgds_gaus(9) = nint((360.0 / float(idim_gaus))*1000.0)
98 ! oct 24-25 - longitude direction incr.
99 kgds_gaus(10) = jdim_gaus/2 ! oct 26-27 - number of circles pole to equator
100 kgds_gaus(12) = 255 ! oct 29 - reserved
101 kgds_gaus(20) = 255 ! oct 5 - not used, set to 255
102
103
104 if (lsm==lsm_noah) then
105 upd_stc=.true.
106 upd_slc=.false. ! not coded
107 elseif (lsm==lsm_noahmp) then
108 upd_stc=.true.
109 upd_slc=.true.
110 endif
111
112 print*,' Start interpolating first ', lsoil_incr, ' surface layers only'
113
114 !----------------------------------------------------------------------
115 ! call gdswzd to compute the lat/lon of each gsi gaussian grid point.
116 !----------------------------------------------------------------------
117
118 iopt = 0
119 fill = -9999.
120 allocate(xpts(idim_gaus*jdim_gaus))
121 allocate(ypts(idim_gaus*jdim_gaus))
122 allocate(lats(idim_gaus*jdim_gaus))
123 allocate(lons(idim_gaus*jdim_gaus))
124 xpts = fill
125 ypts = fill
126 lats = fill
127 lons = fill
128
129 call gdswzd(kgds_gaus,iopt,(idim_gaus*jdim_gaus),fill,xpts,ypts,lons,lats,nret)
130
131 if (nret /= (idim_gaus*jdim_gaus)) then
132 print*,'fatal error: problem in gdswzd. stop.'
133 call mpi_abort(mpi_comm_world, 12, ierr)
134 endif
135
136 deallocate (xpts, ypts)
137
138 allocate(dum2d(idim_gaus,jdim_gaus))
139 dum2d = reshape(lats, (/idim_gaus,jdim_gaus/) )
140 deallocate(lats)
141
142 allocate(lats_rad(jdim_gaus))
143
144 do j = 1, jdim_gaus
145 lats_rad(j) = dum2d(1,jdim_gaus-j+1) * 3.1415926 / 180.0
146 enddo
147
148 dum2d = reshape(lons, (/idim_gaus,jdim_gaus/) )
149 deallocate(lons)
150 allocate(lons_rad(idim_gaus))
151 lons_rad = dum2d(:,1) * 3.1415926 / 180.0
152
153 deallocate(dum2d)
154
155 allocate(agrid(idim,jdim,2))
156 agrid(:,:,1) = reshape(rlo, (/idim,jdim/) )
157 agrid(:,:,2) = reshape(rla, (/idim,jdim/) )
158 agrid = agrid * 3.1415926 / 180.0
159
160 allocate(id1(idim,jdim))
161 allocate(id2(idim,jdim))
162 allocate(jdc(idim,jdim))
163 allocate(s2c(idim,jdim,4))
164
165 !----------------------------------------------------------------------
166 ! compute bilinear weights for each model point from the nearest
167 ! four gsi/gaussian points. does not account for mask. that
168 ! happens later.
169 !----------------------------------------------------------------------
170
171 call remap_coef( 1, idim, 1, jdim, idim_gaus, jdim_gaus, &
172 lons_rad, lats_rad, id1, id2, jdc, s2c, agrid )
173
174 deallocate(lons_rad, lats_rad, agrid)
175
176 ! initialize matrix
177 stcinc = 0.0
178 slcinc = 0.0
179
180 ij_loop : do ij = 1, lensfc
181
182 mask_tile = soilsnow_tile(ij)
183
184 ! get i,j index on array of (idim,jdim) from known ij
185
186 jtile = (ij-1) / idim + 1
187 itile = mod(ij,idim)
188 if (itile==0) itile = idim
189
190 !----------------------------------------------------------------------
191 ! do update to soil temperature grid cells, using bilinear interp
192 !----------------------------------------------------------------------
193
194 if (mask_tile == 1) then
195 ! these are the four nearest grid cells on the gaus grid
196 igaus = id1(itile,jtile)
197 jgaus = jdc(itile,jtile)
198 igausp1 = id2(itile,jtile)
199 jgausp1 = jdc(itile,jtile)+1
200
201 ! make sure gaus grid has soil nearby
202 gaus_has_soil = .false.
203 if (soilsnow_gaus(igaus,jgaus) == 1 .or. &
204 soilsnow_gaus(igausp1,jgaus) == 1 .or. &
205 soilsnow_gaus(igausp1,jgausp1) == 1 .or. &
206 soilsnow_gaus(igaus,jgausp1) == 1) gaus_has_soil = .true.
207
208 if (.not. gaus_has_soil) then
209 cycle ij_loop
210 endif
211
212 stc_inc = 0.0
213 slc_inc = 0.0
214 wsum = 0.0
215
216 if (soilsnow_gaus(igaus,jgaus) == 1) then
217 do k = 1, lsoil_incr
218 if (upd_stc) &
219 stc_inc(k) = stc_inc(k) + (s2c(itile,jtile,1) * stc_inc_gaus(k,igaus,jgaus))
220 if (upd_slc) &
221 slc_inc(k) = slc_inc(k) + (s2c(itile,jtile,1) * slc_inc_gaus(k,igaus,jgaus))
222 enddo
223 wsum = wsum + s2c(itile,jtile,1)
224 endif
225
226 if (soilsnow_gaus(igausp1,jgaus) == 1) then
227 do k = 1, lsoil_incr
228 if (upd_stc) &
229 stc_inc(k) = stc_inc(k) + (s2c(itile,jtile,2) * stc_inc_gaus(k,igausp1,jgaus))
230 if (upd_slc) &
231 slc_inc(k) = slc_inc(k) + (s2c(itile,jtile,2) * slc_inc_gaus(k,igausp1,jgaus))
232 enddo
233 wsum = wsum + s2c(itile,jtile,2)
234 endif
235
236 if (soilsnow_gaus(igausp1,jgausp1) == 1) then
237 do k = 1, lsoil_incr
238 if (upd_stc) &
239 stc_inc(k) = stc_inc(k) + (s2c(itile,jtile,3) * stc_inc_gaus(k,igausp1,jgausp1))
240 if (upd_slc) &
241 slc_inc(k) = slc_inc(k) + (s2c(itile,jtile,3) * slc_inc_gaus(k,igausp1,jgausp1))
242 enddo
243 wsum = wsum + s2c(itile,jtile,3)
244 endif
245
246 if (soilsnow_gaus(igaus,jgausp1) == 1) then
247 do k = 1, lsoil_incr
248 if (upd_stc) &
249 stc_inc(k) = stc_inc(k) + (s2c(itile,jtile,4) * stc_inc_gaus(k,igaus,jgausp1))
250 if (upd_slc) &
251 slc_inc(k) = slc_inc(k) + (s2c(itile,jtile,4) * slc_inc_gaus(k,igaus,jgausp1))
252 enddo
253 wsum = wsum + s2c(itile,jtile,4)
254 endif
255
256 ! normalize increments
257 do k = 1, lsoil_incr
258 stc_inc(k) = stc_inc(k) / wsum
259 stcinc(ij,k) = stc_inc(k)
260 slc_inc(k) = slc_inc(k) / wsum
261 slcinc(ij,k) = slc_inc(k)
262 enddo
263
264 endif ! if soil/snow point
265
266 enddo ij_loop
267
268 write(*,'(a,i2)') ' Finish soil increments interpolation for rank : ', myrank
269
270 deallocate(id1, id2, jdc, s2c)
271
272end subroutine gaussian_to_fv3_interp
273
294
295subroutine add_increment_soil(lsoil_incr,stcinc,slcinc,stc_state,smc_state,slc_state,stc_updated,&
296 slc_updated,soilsnow_tile,soilsnow_fg_tile,lensfc,lsoil,lsm,myrank)
297
298 use mpi
299
300 implicit none
301
302 integer, intent(in) :: lsoil_incr, lensfc, lsoil, myrank, lsm
303
304 integer, intent(in) :: soilsnow_tile(lensfc), soilsnow_fg_tile(lensfc)
305 real, intent(inout) :: stc_state(lensfc, lsoil)
306 real, intent(inout) :: slc_state(lensfc, lsoil)
307 real, intent(inout) :: smc_state(lensfc, lsoil)
308 integer, intent(out) :: stc_updated(lensfc), slc_updated(lensfc)
309
310
311 integer :: ij
312 integer :: mask_tile, mask_fg_tile
313 logical :: upd_slc, upd_stc
314
315 real :: stcinc(lensfc,lsoil)
316 real :: slcinc(lensfc,lsoil)
317
318 integer :: k, nother, nsnowupd
319 integer :: nstcupd, nslcupd, nfrozen, nfrozen_upd
320 logical :: soil_freeze, soil_ice
321
322 stc_updated=0
323 slc_updated=0
324
325 if (lsm==lsm_noah) then
326 upd_stc=.true.
327 upd_slc=.false. ! not coded
328 elseif (lsm==lsm_noahmp) then
329 upd_stc=.true.
330 upd_slc=.true.
331 endif
332
333 print*
334 print*,'adjust soil using increments on cubed-sphere tiles'
335 print*,'updating soil temps', upd_stc
336 print*,'updating soil moisture', upd_slc
337 print*,'adjusting first ', lsoil_incr, ' surface layers only'
338
339 ! initialize variables for counts statitics to be zeros
340 nother = 0 ! grid cells not land
341 nsnowupd = 0 ! grid cells with snow (temperature not yet updated)
342 nslcupd = 0 ! grid cells that are updated
343 nstcupd = 0 ! grid cells that are updated
344 nfrozen = 0 ! not update as frozen soil
345 nfrozen_upd = 0 ! not update as frozen soil
346
347 ij_loop : do ij = 1, lensfc
348
349 mask_tile = soilsnow_tile(ij)
350 mask_fg_tile = soilsnow_fg_tile(ij)
351
352 !----------------------------------------------------------------------
353 ! mask: 1 - soil, 2 - snow, 0 - land-ice, -1 - not land
354 !----------------------------------------------------------------------
355
356 if (mask_tile <= 0) then ! skip if neither soil nor snow
357 nother = nother + 1
358 cycle ij_loop
359 endif
360
361 !----------------------------------------------------------------------
362 ! if snow is present before or after snow update, skip soil analysis
363 !----------------------------------------------------------------------
364
365 if (mask_fg_tile == 2 .or. mask_tile == 2) then
366 nsnowupd = nsnowupd + 1
367 cycle ij_loop
368 endif
369
370 !----------------------------------------------------------------------
371 ! do update to soil temperature grid cells
372 !----------------------------------------------------------------------
373
374 if (mask_tile == 1) then
375
376 !----------------------------------------------------------------------
377 ! add the interpolated increment to the background
378 !----------------------------------------------------------------------
379
380 soil_freeze=.false.
381 soil_ice=.false.
382 do k = 1, lsoil_incr
383
384 if ( stc_state(ij,k) < tfreez) soil_freeze=.true.
385 if ( smc_state(ij,k) - slc_state(ij,k) > 0.001 ) soil_ice=.true.
386
387 if (upd_stc) then
388 stc_state(ij,k) = stc_state(ij,k) + stcinc(ij,k)
389 if (k==1) then
390 stc_updated(ij) = 1
391 nstcupd = nstcupd + 1
392 endif
393 endif
394
395 if ( (stc_state(ij,k) < tfreez) .and. (.not. soil_freeze) .and. (k==1) )&
396 nfrozen_upd = nfrozen_upd + 1
397
398 ! do not do updates if this layer or any above is frozen
399 if ( (.not. soil_freeze ) .and. (.not. soil_ice ) ) then
400 if (upd_slc) then
401 if (k==1) then
402 nslcupd = nslcupd + 1
403 slc_updated(ij) = 1
404 endif
405 ! apply zero limit here (higher, model-specific limits are
406 ! later)
407 slc_state(ij,k) = max(slc_state(ij,k) + slcinc(ij,k), 0.0)
408 smc_state(ij,k) = max(smc_state(ij,k) + slcinc(ij,k), 0.0)
409 endif
410 else
411 if (k==1) nfrozen = nfrozen+1
412 endif
413
414 enddo
415
416 endif ! if soil/snow point
417
418 enddo ij_loop
419
420 write(*,'(a,i2)') ' statistics of grids number processed for rank : ', myrank
421 write(*,'(a,i8)') ' soil grid total', lensfc
422 write(*,'(a,i8)') ' soil grid cells slc updated = ',nslcupd
423 write(*,'(a,i8)') ' soil grid cells stc updated = ',nstcupd
424 write(*,'(a,i8)') ' soil grid cells not updated, frozen = ',nfrozen
425 write(*,'(a,i8)') ' soil grid cells update, became frozen = ',nfrozen_upd
426 write(*,'(a,i8)') ' (not updated yet) snow grid cells = ', nsnowupd
427 write(*,'(a,i8)') ' grid cells, without soil or snow = ', nother
428
429end subroutine add_increment_soil
430
431
443
444subroutine add_increment_snow(snd_inc,mask,lensfc,snd)
445
446 implicit none
447
448 integer, intent(in) :: lensfc
449 real, intent(in) :: snd_inc(lensfc)
450 integer, intent(in) :: mask(lensfc)
451 real, intent(inout) :: snd(lensfc)
452
453 integer :: i
454
455
456 do i =1, lensfc
457 if (mask(i) > 0) then ! if land
458 snd(i) = max( snd(i) + snd_inc(i), 0.)
459 endif
460 enddo
461
462end subroutine add_increment_snow
463
475subroutine calculate_landinc_mask(swe,vtype,stype,lensfc,veg_type_landice,mask)
476
477 implicit none
478
479 integer, intent(in) :: lensfc, veg_type_landice
480 real, intent(in) :: swe(lensfc)
481 real, intent(in) :: vtype(lensfc),stype(lensfc)
482 integer, intent(out) :: mask(lensfc)
483
484 integer :: i
485
486 mask = -1 ! not land
487
488 ! land (but not land-ice)
489 do i=1,lensfc
490 if (nint(stype(i)) .GT. 0) then
491 if (swe(i) .GT. 0.001) then ! snow covered land
492 mask(i) = 2
493 else ! non-snow covered land
494 mask(i) = 1
495 endif
496 end if ! else should work here too
497 if ( nint(vtype(i)) == veg_type_landice ) then ! land-ice
498 mask(i) = 0
499 endif
500 end do
501
502end subroutine calculate_landinc_mask
503
532
533subroutine apply_land_da_adjustments_soil(lsoil_incr, lsm, isot, ivegsrc,lensfc, &
534 lsoil, rsoiltype, mask, stc_bck, stc_adj, smc_adj, slc_adj, &
535 stc_updated, slc_updated, zsoil)
536
537 use mpi
538 use set_soilveg_snippet_mod, only: set_soilveg_noah,set_soilveg_noahmp
539 use sflx_snippet, only: frh2o
540
541 implicit none
542
543 integer, intent(in) :: lsoil_incr, lsm, lensfc, lsoil, isot, ivegsrc
544 real, intent(in) :: rsoiltype(lensfc) ! soil types, as real
545 integer, intent(in) :: mask(lensfc)
546 real, intent(in) :: stc_bck(lensfc, lsoil)
547 integer, intent(in) :: stc_updated(lensfc), slc_updated(lensfc)
548 real, intent(inout) :: smc_adj(lensfc,lsoil), slc_adj(lensfc,lsoil)
549 real, intent(inout) :: stc_adj(lensfc, lsoil)
550 real(kind=4), intent(in) :: zsoil(lsoil)
551
552
553 logical :: frzn_bck, frzn_anl
554 logical :: soil_freeze, soil_ice
555
556 integer :: i, l, n_freeze, n_thaw, ierr
557 integer :: myrank, soiltype, iret, n_stc, n_slc
558 logical :: upd_slc, upd_stc
559
560 real :: slc_new
561
562 real, parameter :: tfreez=273.16
563 real, dimension(30) :: maxsmc, bb, satpsi
564 real, dimension(4) :: dz ! layer thickness
565
566 real, parameter :: hfus=0.3336e06
567 real, parameter :: grav=9.80616
568 real :: smp
569
570 call mpi_comm_rank(mpi_comm_world, myrank, ierr)
571
572 if (lsm==lsm_noah) then
573 upd_stc=.true.
574 upd_slc=.false.
575 elseif (lsm==lsm_noahmp) then
576 upd_stc=.true.
577 upd_slc=.true.
578 endif
579
580 select case (lsm )
581 case(lsm_noah)
582 ! initialise soil properties
583 call set_soilveg_noah(isot, ivegsrc, maxsmc, bb, satpsi, iret)
584 if (iret < 0) then
585 print *, 'FATAL ERROR: problem in set_soilveg_noah'
586 call mpi_abort(mpi_comm_world, 10, ierr)
587 endif
588
589 print *, 'Adjusting noah model smc after stc DA update'
590
591 n_freeze = 0
592 n_thaw = 0
593
594 do i=1,lensfc
595 if (mask(i) > 0) then ! if soil location
596 do l = 1, lsoil
597 frzn_bck = (stc_bck(i,l) .LT. tfreez )
598 frzn_anl = (stc_adj(i,l) .LT. tfreez )
599
600 if (frzn_bck .eqv. frzn_anl) then
601 cycle
602 elseif (frzn_bck .and. .not. frzn_anl) then
603 n_thaw = n_thaw + 1
604 else
605 n_freeze = n_freeze + 1
606 endif
607
608 ! make adjustment (same routine for both)
609 soiltype = nint(rsoiltype(i))
610 ! bb and maxsmc are in the namelist_soilveg, need soiltype index
611 call frh2o(stc_adj(i,l), smc_adj(i,l),slc_adj(i,l), maxsmc(soiltype), &
612 bb(soiltype), satpsi(soiltype),slc_new)
613
614 slc_adj(i,l) = max( min( slc_new, smc_adj(i,l)), 0.0 )
615 enddo
616 endif
617 enddo
618 print *, 'adjusted: ', n_thaw,' thawed,', n_freeze, ' frozen'
619
620 case (lsm_noahmp)
621
622 if (upd_stc) then
623
624 call set_soilveg_noahmp(isot, ivegsrc, maxsmc, bb, satpsi, iret)
625 if (iret < 0) then
626 print *, 'FATAL ERROR: problem in set_soilveg_noahmp'
627 call mpi_abort(mpi_comm_world, 10, ierr)
628 endif
629
630 n_stc = 0
631 n_slc = 0
632
633 do i=1,lensfc
634 if (stc_updated(i) == 1 ) then ! soil-only location
635 n_stc = n_stc+1
636 soiltype = nint(rsoiltype(i))
637 do l = 1, lsoil_incr
638 !case 1: frz ==> frz, recalculate slc, smc remains
639 !case 2: unfrz ==> frz, recalculate slc, smc remains
640 !both cases are considered in the following if case
641 if (stc_adj(i,l) .LT. tfreez )then
642 !recompute supercool liquid water,smc_anl remain unchanged
643 smp = hfus*(tfreez-stc_adj(i,l))/(grav*stc_adj(i,l)) !(m)
644 slc_new=maxsmc(soiltype)*(smp/satpsi(soiltype))**(-1./bb(soiltype))
645 slc_adj(i,l) = max( min( slc_new, smc_adj(i,l)), 0.0 )
646 endif
647 !case 3: frz ==> unfrz, melt all soil ice (if any)
648 if (stc_adj(i,l) .GT. tfreez )then !do not rely on stc_bck
649 slc_adj(i,l)=smc_adj(i,l)
650 endif
651 enddo
652 endif
653 enddo
654
655 endif
656
657 if (upd_slc) then
658
659 dz(1) = -zsoil(1)
660 do l = 2,lsoil
661 dz(l) = -zsoil(l) + zsoil(l-1)
662 enddo
663 print *, 'Applying soil moisture mins '
664
665 do i=1,lensfc
666 if (slc_updated(i) == 1 ) then
667 n_slc = n_slc+1
668 ! apply SM bounds (later: add upper SMC limit)
669 do l = 1, lsoil_incr
670 ! noah-mp minimum is 1 mm per layer (in SMC)
671 ! no need to maintain frozen amount, would be v. small.
672 slc_adj(i,l) = max( 0.001/dz(l), slc_adj(i,l) )
673 smc_adj(i,l) = max( 0.001/dz(l), smc_adj(i,l) )
674 enddo
675 endif
676 enddo
677 endif
678
679 case default
680 print *, 'FATAL ERROR: unrecognised LSM,', lsm
681 call mpi_abort(mpi_comm_world, 10, ierr)
682 end select
683
684 write(*,'(a,i2)') 'statistics of grids number processed for rank : ', myrank
685 write(*,'(a,i8)') ' soil grid total', lensfc
686 write(*,'(a,i8)') ' soil grid cells with slc update', n_slc
687 write(*,'(a,i8)') ' soil grid cells with stc update', n_stc
688
689end subroutine apply_land_da_adjustments_soil
690
696
705
706subroutine apply_land_da_adjustments_snd(lsm, lensfc, mask, swe_bck, snd_bck, snd_anl, swe_adj)
707
708 use mpi
709 use bulk_snow_module, only: calc_density
710
711 implicit none
712
713 integer, intent(in) :: lsm, lensfc
714 integer, intent(in) :: mask(lensfc)
715 real, intent(in) :: swe_bck(lensfc), snd_bck(lensfc)
716 real, intent(in) :: snd_anl(lensfc)
717 real, intent(inout) :: swe_adj(lensfc)
718
719 integer :: ierr, myrank, i
720
721 real :: density(lensfc)
722
723 call mpi_comm_rank(mpi_comm_world, myrank, ierr)
724
725 if (lsm .NE. lsm_noah) then
726 print *, 'FATAL ERROR: apply_land_da_adjustments not coded for models other than noah', lsm
727 call mpi_abort(mpi_comm_world, 10, ierr)
728 endif
729
730 ! calculate snow density from forecasts
731 call calc_density(lensfc, mask, swe_bck, snd_bck, myrank, density)
732
733 do i =1, lensfc
734 if ( mask(i)>0 ) then
735 swe_adj(i) = snd_anl(i)*density(i)
736 endif
737 enddo
738
739
740end subroutine apply_land_da_adjustments_snd
741
742end module land_increments
This module contains routines that read and write data.
integer, public jdim_gaus
'j' dimension of GSI gaussian grid.
integer, dimension(:,:), allocatable, public soilsnow_gaus
GSI soil / snow mask for land on the gaussian grid.
real, dimension(:,:,:), allocatable, public slc_inc_gaus
GSI soil moisture increments on the gaussian grid.
real, dimension(:,:,:), allocatable, public stc_inc_gaus
GSI soil temperature increments on the gaussian grid.
integer, public idim_gaus
'i' dimension of GSI gaussian grid.
Module containing utility routines.
Definition utils.F90:7
subroutine, public remap_coef(is, ie, js, je, im, jm, lon, lat, id1, id2, jdc, s2c, agrid)
Generate the weights and index of the grids used in the bilinear interpolation.
Definition utils.F90:40