source: trunk/libf/dyn3dpar/ugeostr.F90 @ 109

Last change on this file since 109 was 66, checked in by emillour, 14 years ago

EM: Mise a niveau par rapport a la version terrestre (LMDZ5V2.0-dev, rev 1487)

  • Mise a jour des scripts (terrestres) 'makegcm' et 'create_make_gcm'
  • Ajout du script 'makelmdz' (version "amelioree", en Bash, de makegcm)
  • Mise a jour des routines dans phylmd (sauf regr_lat_time_climoz_m.F)
  • disvert (dans dyn3d et dyn3dpar): passage au Fortran 90
  • parallel.F90 (dyn3dpar): correction bug
  • etat0_netcdf.F90 (dyn3d et dyn3dpar) : mise a jour mineure
  • ce0l.F90 (dyn3dpar) : correction bug
  • abort_gcm.F (dyn3dpar) : correction bug
  • ugeostr.F90 (dyn3d et dyn3dpar) : passage au Fortran 90
  • fluxstokenc_p.F (dyn3dpar) : correction bug
  • iniacademic.F90 (dyn3d et dyn3dpar) : passage au Fortran 90
  • friction_p.F (dyn3dpar) : correction bug
  • infotrac.F90 (dyn3d et dyn3dpar) : correction bug mineur sur lecture traceurs
  • caladvtrac.F (dyn3d) : modifications cosmetiques
File size: 1.4 KB
RevLine 
[1]1!
[66]2! $Id: ugeostr.F90 1474 2011-01-14 11:04:45Z lguez $
[1]3!
[66]4subroutine ugeostr(phi,ucov)
[1]5
[66]6  ! Calcul du vent covariant geostrophique a partir du champ de
7  ! geopotentiel.
8  ! We actually compute: (1 - cos^8 \phi) u_g
9  ! to have a wind going smoothly to 0 at the equator.
10  ! We assume that the surface pressure is uniform so that model
11  ! levels are pressure levels.
[1]12
[66]13  implicit none
[1]14
[66]15  include "dimensions.h"
16  include "paramet.h"
17  include "comconst.h"
18  include "comgeom2.h"
[1]19
[66]20  real ucov(iip1,jjp1,llm),phi(iip1,jjp1,llm)
21  real um(jjm,llm),fact,u(iip1,jjm,llm)
22  integer i,j,l
[1]23
[66]24  real zlat
[1]25
[66]26  um(:,:)=0 ! initialize um()
[1]27
[66]28  DO j=1,jjm
[1]29
[66]30     if (abs(sin(rlatv(j))).lt.1.e-4) then
31        zlat=1.e-4
32     else
33        zlat=rlatv(j)
34     endif
35     fact=cos(zlat)
36     fact=fact*fact
37     fact=fact*fact
38     fact=fact*fact
39     fact=(1.-fact)/ &
40          (2.*omeg*sin(zlat)*(rlatu(j+1)-rlatu(j)))
41     fact=-fact/rad
42     DO l=1,llm
43        DO i=1,iim
44           u(i,j,l)=fact*(phi(i,j+1,l)-phi(i,j,l))
45           um(j,l)=um(j,l)+u(i,j,l)/REAL(iim)
46        ENDDO
47     ENDDO
48  ENDDO
49  call dump2d(jjm,llm,um,'Vent-u geostrophique')
[1]50
[66]51  !   calcul des champ de vent:
[1]52
[66]53  DO l=1,llm
54     DO i=1,iip1
55        ucov(i,1,l)=0.
56        ucov(i,jjp1,l)=0.
57     end DO
58     DO  j=2,jjm
59        DO  i=1,iim
60           ucov(i,j,l) = 0.5*(u(i,j,l)+u(i,j-1,l))*cu(i,j)
61        end DO
62        ucov(iip1,j,l)=ucov(1,j,l)
63     end DO
64  end DO
[1]65
[66]66  print *, 301
[1]67
[66]68end subroutine ugeostr
Note: See TracBrowser for help on using the repository browser.