/[lmdze]/trunk/dyn3d/principal_cshift.f
ViewVC logotype

Contents of /trunk/dyn3d/principal_cshift.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 124 - (show annotations)
Thu Feb 5 15:19:37 2015 UTC (9 years, 3 months ago) by guez
File size: 956 byte(s)
Moved some processing from fxhyp_loop_ik to fxhyp. Now fxhyp_loop_ik
does not necessarily give longitudes near [-pi, pi]. In fxhyp, we look
in rlonm025 whether we need to move the array toward [-pi, pi]. If so,
we apply the same move to all grids: rlonm025, rlonv, rlonp025, rlonu
and the corresponding derivatives. The move itself is done by the new
procedure principal_cshift. This revision makes the logic
clearer. (For example, we do not have a saved variable is2 in
fxhyp_loop_ik any longer and we remove a test on ik in fxhyp_loop_ik.)

Fixed a bad error message in fxhyp_loop_ik: talked about rlonu when
xvrai is not always rlonu.

1 module principal_cshift_m
2
3 implicit none
4
5 contains
6
7 subroutine principal_cshift(is2, xlon, xprimm)
8
9 USE dimens_m, ONLY: iim
10 use nr_util, only: twopi
11 use serre, only: clon
12
13 integer, intent(in):: is2
14 real, intent(inout):: xlon(:), xprimm(:) ! (iim + 1)
15
16 !-----------------------------------------------------
17
18 if (is2 /= 0) then
19 IF (clon <= 0.) THEN
20 IF (is2 /= 1) THEN
21 xlon(:is2 - 1) = xlon(:is2 - 1) + twopi
22 xlon(:iim) = cshift(xlon(:iim), shift = is2 - 1)
23 xprimm(:iim) = cshift(xprimm(:iim), shift = is2 - 1)
24 END IF
25 else
26 xlon(is2 + 1:iim) = xlon(is2 + 1:iim) - twopi
27 xlon(:iim) = cshift(xlon(:iim), shift = is2)
28 xprimm(:iim) = cshift(xprimm(:iim), shift = is2)
29 end IF
30 end if
31
32 xlon(iim + 1) = xlon(1) + twopi
33 xprimm(iim + 1) = xprimm(1)
34
35 end subroutine principal_cshift
36
37 end module principal_cshift_m

  ViewVC Help
Powered by ViewVC 1.1.21