1 subroutine da_wpec_constraint_geoterm (cori,rho, u, v, term_x, term_y)
3 !---------------------------------------------------------------------------
4 ! Purpose: calculates nonlinear geostrophic term in wpec constraint equation.
6 ! method: term is k x rho f u on a single level.
7 !---------------------------------------------------------------------------
10 real, intent(in) :: cori(ims:ime,jms:jme) ! Coriolis factor.
11 real, intent(in) :: rho(ims:ime,jms:jme) ! Density
12 real, intent(in) :: u(:,:) ! u wind comp. (dot pts)
13 real, intent(in) :: v(:,:) ! v wind comp. (dot pts)
14 real, intent(out) :: term_x(:,:) ! x component of term.
15 real, intent(out) :: term_y(:,:) ! y component of term.
16 integer :: is, ie ! 1st dim. end points.
17 integer :: js, je ! 2nd dim. end points.
19 if (trace_use) call da_trace_entry("da_wpec_constraint_geoterm")
21 !---------------------------------------------------------------------------
23 !---------------------------------------------------------------------------
25 ! Computation to check for edge of domain:
26 is = its; ie = ite; js = jts; je = jte
27 if (.not. global .and. its == ids) is = ids+1
28 if (.not. global .and. ite == ide) ie = ide-1
29 if (jts == jds) js = jds+1; if (jte == jde) je = jde-1
31 !---------------------------------------------------------------------------
32 ! [2.0] Calculate term_x = -f rho v~:
33 !---------------------------------------------------------------------------
35 term_x(is:ie, js:je) = term_x(is:ie, js:je)-rho(is:ie, js:je) &
36 * v(is:ie, js:je) * cori(is:ie, js:je)
38 !---------------------------------------------------------------------------
39 ! [3.0] Calculate term_y = f rho u~:
40 !---------------------------------------------------------------------------
42 term_y(is:ie, js:je) = term_y(is:ie, js:je)+rho(is:ie, js:je) &
43 * u(is:ie, js:je) * cori(is:ie, js:je)
45 if (trace_use) call da_trace_exit("da_wpec_constraint_geoterm")
47 end subroutine da_wpec_constraint_geoterm