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

Annotation of /trunk/dyn3d/covnat.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 261 - (hide annotations)
Wed Mar 7 13:33:15 2018 UTC (6 years, 3 months ago) by guez
File size: 1052 byte(s)
Make procedures writehist and writedynav as identical as
possible. Remove output of phis in writedynav: already in histins and
constant. Add output of temp and q in writehist.

Remove unused variable ndm of module dimens_m. Remove unused variables
kftd, mvar, jcfil and jcfllm of module paramet_m.

1 guez 61 module covnat_m
2 guez 3
3 guez 61 IMPLICIT NONE
4 guez 3
5 guez 61 contains
6 guez 3
7 guez 61 SUBROUTINE covnat(klevel, ucov, vcov, unat, vnat)
8 guez 3
9 guez 61 ! From LMDZ4/libf/dyn3d/covnat.F, version 1.1.1.1 2004/05/19 12:53:07
10 guez 3
11 guez 67 USE comgeom, ONLY: cu, cv
12 guez 61 USE paramet_m, ONLY: iip1, iip2, ip1jm, ip1jmp1
13 guez 3
14 guez 61 ! Authors: F. Hourdin, Phu Le Van
15 guez 3
16 guez 261 ! Objet : calcul des composantes naturelles du vent à partir des
17 guez 61 ! composantes covariantes.
18 guez 3
19 guez 61 INTEGER, intent(in):: klevel
20     REAL, intent(in):: ucov(ip1jmp1, klevel), vcov(ip1jm, klevel)
21     REAL, intent(out):: unat(ip1jmp1, klevel), vnat(ip1jm, klevel)
22 guez 45
23 guez 61 ! Local:
24     INTEGER l, ij
25 guez 45
26 guez 61 !------------------------------------------------------------------
27 guez 45
28 guez 61 DO l = 1, klevel
29     DO ij = 1, iip1
30     unat(ij, l) =0.
31     END DO
32    
33     DO ij = iip2, ip1jm
34     unat(ij, l) = ucov(ij, l) / cu(ij)
35     ENDDO
36    
37     DO ij = ip1jm+1, ip1jmp1
38     unat(ij, l) =0.
39     END DO
40    
41     DO ij = 1, ip1jm
42     vnat(ij, l) = vcov(ij, l) / cv(ij)
43     ENDDO
44     ENDDO
45    
46     END SUBROUTINE covnat
47    
48     end module covnat_m

  ViewVC Help
Powered by ViewVC 1.1.21