source: LMDZ4/branches/LMDZ4V5.0-dev/libf/dyn3dpar/ugeostr.F @ 2057

Last change on this file since 2057 was 1299, checked in by Laurent Fairhead, 15 years ago

Nettoyage general pour se rapprocher des normes et éviter des erreurs a la
compilation:

  • tous les FLOAT() sont remplacés par des REAL()
  • tous les STOP dans phylmd sont remplacés par des appels à abort_gcm
  • le common control défini dans le fichier control.h est remplacé par le module control_mod pour éviter des messages sur l'alignement des variables dans les déclarations
  • des $Header$ remplacés par des $Id$ pour svn

Quelques remplacements à faire ont pu m'échapper


General cleanup of the code to try and adhere to norms and to prevent some
compilation errors:

  • all FLOAT() instructions have been replaced by REAL() instructions
  • all STOP instructions in phylmd have been replaced by calls to abort_gcm
  • the common block control defined in the control.h file has been replaced by the control_mod to prevent compilation warnings on the alignement of declared variables
  • $Header$ replaced by $Id$ for svn

Some changes which should have been made might have escaped me

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 1.5 KB
Line 
1!
2! $Id: ugeostr.F 1299 2010-01-20 14:27:21Z emillour $
3!
4      subroutine ugeostr(phi,ucov)
5
6
7c  Calcul du vent covariant geostrophique a partir du champs de
8c  geopotentiel en supposant que le vent au sol est nul.
9
10      implicit none
11
12#include "dimensions.h"
13#include "paramet.h"
14#include "comconst.h"
15#include "comgeom2.h"
16
17      real ucov(iip1,jjp1,llm),phi(iip1,jjp1,llm)
18      real um(jjm,llm),fact,u(iip1,jjm,llm)
19      integer i,j,l
20
21      real zlat
22
23      um(:,:)=0 ! initialize um()
24
25      DO j=1,jjm
26
27         if (abs(sin(rlatv(j))).lt.1.e-4) then
28             zlat=1.e-4
29         else
30             zlat=rlatv(j)
31         endif
32         fact=cos(zlat)
33         fact=fact*fact
34         fact=fact*fact
35         fact=fact*fact
36         fact=(1.-fact)/
37     s    (2.*omeg*sin(zlat)*(rlatu(j+1)-rlatu(j)))
38         fact=-fact/rad
39         DO l=1,llm
40            DO i=1,iim
41               u(i,j,l)=fact*(phi(i,j+1,l)-phi(i,j,l))
42               um(j,l)=um(j,l)+u(i,j,l)/REAL(iim)
43            ENDDO
44         ENDDO
45      ENDDO
46      call dump2d(jjm,llm,um,'Vent-u geostrophique')
47
48c
49c-----------------------------------------------------------------------
50c   calcul des champ de vent:
51c   -------------------------
52
53      DO 301 l=1,llm
54         DO 302 i=1,iip1
55            ucov(i,1,l)=0.
56            ucov(i,jjp1,l)=0.
57302      CONTINUE
58         DO 304 j=2,jjm
59            DO 305 i=1,iim
60               ucov(i,j,l) = 0.5*(u(i,j,l)+u(i,j-1,l))*cu(i,j)
61305         CONTINUE
62            ucov(iip1,j,l)=ucov(1,j,l)
63304      CONTINUE
64301   CONTINUE
65
66      print*,301
67
68      return
69      end
Note: See TracBrowser for help on using the repository browser.