1 module rotate_winds_module
6 use misc_definitions_module
9 integer :: orig_selected_projection
14 ! ARW Wind Rotation Code
17 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
20 ! Purpose: Rotate grid-relative winds to Earth-relative winds !
21 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
22 subroutine map_to_met(u, u_mask, v, v_mask, &
25 xlon_u, xlon_v, xlat_u, xlat_v)
30 integer, intent(in) :: us1, us2, ue1, ue2, vs1, vs2, ve1, ve2
31 real, pointer, dimension(:,:) :: u, v, xlon_u, xlon_v, xlat_u, xlat_v
32 type (bitarray), intent(in) :: u_mask, v_mask
34 orig_selected_projection = iget_selected_domain()
35 call select_domain(SOURCE_PROJ)
36 call metmap_xform(u, u_mask, v, v_mask, &
39 xlon_u, xlon_v, xlat_u, xlat_v, 1)
40 call select_domain(orig_selected_projection)
42 end subroutine map_to_met
45 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
48 ! Purpose: Rotate Earth-relative winds to grid-relative winds !
49 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
50 subroutine met_to_map(u, u_mask, v, v_mask, &
53 xlon_u, xlon_v, xlat_u, xlat_v)
58 integer, intent(in) :: us1, us2, ue1, ue2, vs1, vs2, ve1, ve2
59 real, pointer, dimension(:,:) :: u, v, xlon_u, xlon_v, xlat_u, xlat_v
60 type (bitarray), intent(in) :: u_mask, v_mask
62 orig_selected_projection = iget_selected_domain()
64 call metmap_xform(u, u_mask, v, v_mask, &
67 xlon_u, xlon_v, xlat_u, xlat_v, -1)
68 call select_domain(orig_selected_projection)
70 end subroutine met_to_map
73 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
74 ! Name: metmap_xform !
76 ! Purpose: Do the actual work of rotating winds for C grid. !
77 ! If idir= 1, rotate grid-relative winds to Earth-relative winds !
78 ! If idir=-1, rotate Earth-relative winds to grid-relative winds !
80 ! ASSUMPTIONS: 1) MEMORY ORDER IS XY. !
81 ! 2) U ARRAY HAS ONE MORE COLUMN THAN THE V ARRAY, AND V ARRAY !
82 ! HAS ONE MORE ROW THAN U ARRAY. !
83 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
84 subroutine metmap_xform(u, u_mask, v, v_mask, &
85 us1, us2, ue1, ue2, vs1, vs2, ve1, ve2, &
86 xlon_u, xlon_v, xlat_u, xlat_v, idir)
91 integer, intent(in) :: us1, us2, ue1, ue2, vs1, vs2, ve1, ve2, idir
92 real, pointer, dimension(:,:) :: u, v, xlon_u, xlon_v, xlat_u, xlat_v
93 type (bitarray), intent(in) :: u_mask, v_mask
97 real :: u_weight, v_weight
98 real :: u_map, v_map, alpha, diff
99 real, pointer, dimension(:,:) :: u_new, v_new, u_mult, v_mult
100 logical :: do_last_col_u, do_last_row_u, do_last_col_v, do_last_row_v
102 ! If the proj_info structure has not been initialized, we don't have
103 ! information about the projection and standard longitude.
104 if (proj_stack(current_nest_number)%init) then
106 ! Only rotate winds for Lambert conformal, polar stereographic, or Cassini
107 if ((proj_stack(current_nest_number)%code == PROJ_LC) .or. &
108 (proj_stack(current_nest_number)%code == PROJ_PS) .or. &
109 (proj_stack(current_nest_number)%code == PROJ_CASSINI)) then
110 call mprintf((idir == 1),LOGFILE,'Rotating map winds to earth winds.')
111 call mprintf((idir == -1),LOGFILE,'Rotating earth winds to grid winds')
113 allocate(u_mult(us1:ue1,us2:ue2))
114 allocate(v_mult(vs1:ve1,vs2:ve2))
118 if (bitarray_test(u_mask, i-us1+1, j-us2+1)) then
128 if (bitarray_test(v_mask, i-vs1+1, j-vs2+1)) then
136 if (ue1-us1 == ve1-vs1) then
137 do_last_col_u = .false.
138 do_last_col_v = .true.
140 do_last_col_u = .true.
141 do_last_col_v = .false.
144 if (ue2-us2 == ve2-vs2) then
145 do_last_row_u = .true.
146 do_last_row_v = .false.
148 do_last_row_u = .false.
149 do_last_row_v = .true.
152 ! Create arrays to hold rotated winds
153 allocate(u_new(us1:ue1, us2:ue2))
154 allocate(v_new(vs1:ve1, vs2:ve2))
160 diff = idir * (xlon_u(i,j) - proj_stack(current_nest_number)%stdlon)
161 if (diff > 180.) then
163 else if (diff < -180.) then
167 ! Calculate the rotation angle, alpha, in radians
168 if (proj_stack(current_nest_number)%code == PROJ_LC) then
169 alpha = diff * proj_stack(current_nest_number)%cone * rad_per_deg * proj_stack(current_nest_number)%hemi
170 else if (proj_stack(current_nest_number)%code == PROJ_CASSINI) then
172 diff = xlon_u(i,j)-xlon_u(i,j-1)
173 if (diff > 180.) then
175 else if (diff < -180.) then
178 alpha = atan2( (-cos(xlat_u(i,j)*rad_per_deg) * diff*rad_per_deg), &
179 (xlat_u(i,j)-xlat_u(i,j-1))*rad_per_deg &
181 else if (j == us2) then
182 diff = xlon_u(i,j+1)-xlon_u(i,j)
183 if (diff > 180.) then
185 else if (diff < -180.) then
188 alpha = atan2( (-cos(xlat_u(i,j)*rad_per_deg) * diff*rad_per_deg), &
189 (xlat_u(i,j+1)-xlat_u(i,j))*rad_per_deg &
192 diff = xlon_u(i,j+1)-xlon_u(i,j-1)
193 if (diff > 180.) then
195 else if (diff < -180.) then
198 alpha = atan2( (-cos(xlat_u(i,j)*rad_per_deg) * diff*rad_per_deg), &
199 (xlat_u(i,j+1)-xlat_u(i,j-1))*rad_per_deg &
203 alpha = diff * rad_per_deg * proj_stack(current_nest_number)%hemi
208 ! On C grid, take U_ij, and get V value at the same lat/lon
209 ! by averaging the four surrounding V points
210 if (bitarray_test(u_mask, i-us1+1, j-us2+1)) then
213 if (j == ue2 .and. do_last_row_u) then
214 v_weight = v_mult(i,j)
215 v_map = v(i,j)*v_mult(i,j)
217 v_weight = v_mult(i,j) + v_mult(i,j+1)
218 v_map = v(i,j)*v_mult(i,j) + v(i,j+1)*v_mult(i,j+1)
220 else if (i == ue1 .and. do_last_col_u) then
221 if (j == ue2 .and. do_last_row_u) then
222 v_weight = v_mult(i-1,j)
225 v_weight = v_mult(i-1,j) + v_mult(i-1,j+1)
226 v_map = v(i-1,j)*v_mult(i-1,j) + v(i-1,j+1)*v_mult(i-1,j+1)
228 else if (j == ue2 .and. do_last_row_u) then
229 v_weight = v_mult(i-1,j-1) + v_mult(i,j-1)
230 v_map = v(i-1,j-1)*v_mult(i-1,j-1) + v(i,j-1)*v_mult(i,j-1)
232 v_weight = v_mult(i-1,j) + v_mult(i-1,j+1) + v_mult(i,j) + v_mult(i,j+1)
233 v_map = v(i-1,j)*v_mult(i-1,j) + v(i-1,j+1)*v_mult(i-1,j+1) + v(i,j)*v_mult(i,j) + v(i,j+1)*v_mult(i,j+1)
235 if (v_weight > 0.) then
236 u_new(i,j) = cos(alpha)*u_map + sin(alpha)*v_map/v_weight
251 diff = idir * (xlon_v(i,j) - proj_stack(current_nest_number)%stdlon)
252 if (diff > 180.) then
254 else if (diff < -180.) then
258 if (proj_stack(current_nest_number)%code == PROJ_LC) then
259 alpha = diff * proj_stack(current_nest_number)%cone * rad_per_deg * proj_stack(current_nest_number)%hemi
260 else if (proj_stack(current_nest_number)%code == PROJ_CASSINI) then
262 diff = xlon_v(i,j)-xlon_v(i,j-1)
263 if (diff > 180.) then
265 else if (diff < -180.) then
268 alpha = atan2( (-cos(xlat_v(i,j)*rad_per_deg) * diff*rad_per_deg), &
269 (xlat_v(i,j)-xlat_v(i,j-1))*rad_per_deg &
271 else if (j == vs2) then
272 diff = xlon_v(i,j+1)-xlon_v(i,j)
273 if (diff > 180.) then
275 else if (diff < -180.) then
278 alpha = atan2( (-cos(xlat_v(i,j)*rad_per_deg) * diff*rad_per_deg), &
279 (xlat_v(i,j+1)-xlat_v(i,j))*rad_per_deg &
282 diff = xlon_v(i,j+1)-xlon_v(i,j-1)
283 if (diff > 180.) then
285 else if (diff < -180.) then
288 alpha = atan2( (-cos(xlat_v(i,j)*rad_per_deg) * diff*rad_per_deg), &
289 (xlat_v(i,j+1)-xlat_v(i,j-1))*rad_per_deg &
293 alpha = diff * rad_per_deg * proj_stack(current_nest_number)%hemi
298 if (bitarray_test(v_mask, i-vs1+1, j-vs2+1)) then
301 if (i == ve1 .and. do_last_col_v) then
302 u_weight = u_mult(i,j)
303 u_map = u(i,j)*u_mult(i,j)
305 u_weight = u_mult(i,j) + u_mult(i+1,j)
306 u_map = u(i,j)*u_mult(i,j) + u(i+1,j)*u_mult(i+1,j)
308 else if (j == ve2 .and. do_last_row_v) then
309 if (i == ve1 .and. do_last_col_v) then
310 u_weight = u_mult(i,j-1)
311 u_map = u(i,j-1)*u_mult(i,j-1)
313 u_weight = u_mult(i,j-1) + u_mult(i+1,j-1)
314 u_map = u(i,j-1)*u_mult(i,j-1) + u(i+1,j-1)*u_mult(i+1,j-1)
316 else if (i == ve1 .and. do_last_col_v) then
317 u_weight = u_mult(i,j) + u_mult(i,j-1)
318 u_map = u(i,j)*u_mult(i,j) + u(i,j-1)*u_mult(i,j-1)
320 u_weight = u_mult(i,j-1) + u_mult(i,j) + u_mult(i+1,j-1) + u_mult(i+1,j)
321 u_map = u(i,j-1)*u_mult(i,j-1) + u(i,j)*u_mult(i,j) + u(i+1,j-1)*u_mult(i+1,j-1) + u(i+1,j)*u_mult(i+1,j)
323 if (u_weight > 0.) then
324 v_new(i,j) = -sin(alpha)*u_map/u_weight + cos(alpha)*v_map
335 ! Copy rotated winds back into argument arrays
346 call mprintf(.true.,ERROR,'In metmap_xform(), uninitialized proj_info structure.')
349 end subroutine metmap_xform
353 ! NMM Wind Rotation Code
356 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
357 ! Name: map_to_met_nmm !
359 ! Purpose: Rotate grid-relative winds to Earth-relative winds !
360 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
361 subroutine map_to_met_nmm(u, u_mask, v, v_mask, &
362 vs1, vs2, ve1, ve2, &
368 integer, intent(in) :: vs1, vs2, ve1, ve2
369 real, pointer, dimension(:,:) :: u, v, xlat_v, xlon_v
370 type (bitarray), intent(in) :: u_mask, v_mask
372 orig_selected_projection = iget_selected_domain()
373 call select_domain(SOURCE_PROJ)
374 call metmap_xform_nmm(u, u_mask, v, v_mask, &
375 vs1, vs2, ve1, ve2, &
377 call select_domain(orig_selected_projection)
379 end subroutine map_to_met_nmm
382 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
383 ! Name: met_to_map_nmm !
385 ! Purpose: Rotate Earth-relative winds to grid-relative winds !
386 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
387 subroutine met_to_map_nmm(u, u_mask, v, v_mask, &
388 vs1, vs2, ve1, ve2, &
394 integer, intent(in) :: vs1, vs2, ve1, ve2
395 real, pointer, dimension(:,:) :: u, v, xlat_v, xlon_v
396 type (bitarray), intent(in) :: u_mask, v_mask
398 orig_selected_projection = iget_selected_domain()
399 call select_domain(1)
400 call metmap_xform_nmm(u, u_mask, v, v_mask, &
401 vs1, vs2, ve1, ve2, &
403 call select_domain(orig_selected_projection)
405 end subroutine met_to_map_nmm
408 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
409 ! Name: metmap_xform_nmm !
411 ! Purpose: Do the actual work of rotating winds for E grid. !
412 ! If idir= 1, rotate grid-relative winds to Earth-relative winds !
413 ! If idir=-1, rotate Earth-relative winds to grid-relative winds !
415 ! ASSUMPTIONS: 1) MEMORY ORDER IS XY. !
416 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
417 subroutine metmap_xform_nmm(u, u_mask, v, v_mask, &
418 vs1, vs2, ve1, ve2, &
419 xlat_v, xlon_v, idir)
424 integer, intent(in) :: vs1, vs2, ve1, ve2, idir
425 real, pointer, dimension(:,:) :: u, v, xlat_v, xlon_v
426 type (bitarray), intent(in) :: u_mask, v_mask
430 real :: u_map, v_map, diff, alpha
431 real :: phi0, lmbd0, big_denominator, relm, rlat_v,rlon_v, clontemp
432 real :: sin_phi0, cos_phi0, cos_alpha, sin_alpha
433 real, pointer, dimension(:,:) :: u_new, v_new
436 ! If the proj_info structure has not been initialized, we don't have
437 ! information about the projection and standard longitude.
438 if (proj_stack(current_nest_number)%init) then
440 if (proj_stack(current_nest_number)%code == PROJ_ROTLL) then
442 call mprintf((idir == 1),LOGFILE,'Rotating map winds to earth winds.')
443 call mprintf((idir == -1),LOGFILE,'Rotating earth winds to grid winds')
445 ! Create arrays to hold rotated winds
446 allocate(u_new(vs1:ve1, vs2:ve2))
447 allocate(v_new(vs1:ve1, vs2:ve2))
449 phi0 = proj_stack(current_nest_number)%lat1 * rad_per_deg
451 clontemp= proj_stack(current_nest_number)%lon1
453 if (clontemp .lt. 0.) then
454 lmbd0 = (clontemp + 360.) * rad_per_deg
456 lmbd0 = (clontemp) * rad_per_deg
465 ! Calculate the sine and cosine of rotation angle
466 rlat_v = xlat_v(i,j) * rad_per_deg
467 rlon_v = xlon_v(i,j) * rad_per_deg
468 relm = rlon_v - lmbd0
469 big_denominator = cos(asin( &
470 cos_phi0 * sin(rlat_v) - &
471 sin_phi0 * cos(rlat_v) * cos(relm) &
474 sin_alpha = sin_phi0 * sin(relm) / &
477 cos_alpha = (cos_phi0 * cos(rlat_v) + &
478 sin_phi0 * sin(rlat_v) * cos(relm)) / &
482 if (bitarray_test(u_mask, i-vs1+1, j-vs2+1)) then
484 if (bitarray_test(v_mask, i-vs1+1, j-vs2+1)) then
490 u_new(i,j) = cos_alpha*u_map + idir*sin_alpha*v_map
496 if (bitarray_test(v_mask, i-vs1+1, j-vs2+1)) then
498 if (bitarray_test(u_mask, i-vs1+1, j-vs2+1)) then
504 v_new(i,j) = -idir*sin_alpha*u_map + cos_alpha*v_map
512 ! Copy rotated winds back into argument arrays
519 ! Only rotate winds for Lambert conformal, polar stereographic, or Cassini
520 else if ((proj_stack(current_nest_number)%code == PROJ_LC) .or. &
521 (proj_stack(current_nest_number)%code == PROJ_PS) .or. &
522 (proj_stack(current_nest_number)%code == PROJ_CASSINI)) then
524 call mprintf((idir == 1),LOGFILE,'Rotating map winds to earth winds.')
525 call mprintf((idir == -1),LOGFILE,'Rotating earth winds to grid winds')
527 ! Create arrays to hold rotated winds
528 allocate(u_new(vs1:ve1, vs2:ve2))
529 allocate(v_new(vs1:ve1, vs2:ve2))
534 diff = idir * (xlon_v(i,j) - proj_stack(current_nest_number)%stdlon)
535 if (diff > 180.) then
537 else if (diff < -180.) then
541 ! Calculate the rotation angle, alpha, in radians
542 if (proj_stack(current_nest_number)%code == PROJ_LC) then
543 alpha = diff * proj_stack(current_nest_number)%cone * &
544 rad_per_deg * proj_stack(current_nest_number)%hemi
546 alpha = diff * rad_per_deg * proj_stack(current_nest_number)%hemi
550 if (bitarray_test(u_mask, i-vs1+1, j-vs2+1)) then
552 if (bitarray_test(v_mask, i-vs1+1, j-vs2+1)) then
558 u_new(i,j) = cos(alpha)*u_map + idir*sin(alpha)*v_map
564 if (bitarray_test(v_mask, i-vs1+1, j-vs2+1)) then
566 if (bitarray_test(u_mask, i-vs1+1, j-vs2+1)) then
572 v_new(i,j) = -idir*sin(alpha)*u_map + cos(alpha)*v_map
580 ! Copy rotated winds back into argument arrays
590 call mprintf(.true.,ERROR,'In metmap_xform_nmm(), uninitialized proj_info structure.')
593 end subroutine metmap_xform_nmm
595 end module rotate_winds_module