source: LMDZ6/branches/Amaury_dev/libf/dyn3d/iniinterp_horiz.F90

Last change on this file was 5159, checked in by abarral, 7 weeks ago

Put dimensions.h and paramet.h into modules

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 4.3 KB
Line 
1
2! $Header$
3
4SUBROUTINE iniinterp_horiz(imo, jmo, imn, jmn, kllm, &
5        rlonuo, rlatvo, rlonun, rlatvn, &
6        ktotal, iik, jjk, jk, ik, intersec, airen)
7
8  IMPLICIT NONE
9
10
11
12  ! ---------------------------------------------------------
13  ! Prepare l' interpolation des variables d'une grille LMDZ
14  !  dans une autre grille LMDZ en conservant la quantite
15  !  totale pour les variables intensives (/m2) : ex : Pression au sol
16
17  !   (Pour chaque case autour d'un point scalaire de la nouvelle
18  !    grille, on calcule la surface (en m2)en intersection avec chaque
19  !    case de l'ancienne grille , pour la future interpolation)
20
21  ! on calcule aussi l' aire dans la nouvelle grille
22
23
24  !   Auteur:  F.Forget 01/1995
25  !   -------
26
27  ! ---------------------------------------------------------
28  !   Declarations:
29  ! ==============
30
31  !  ARGUMENTS
32  !  """""""""
33  ! INPUT
34  INTEGER :: imo, jmo ! dimensions ancienne grille
35  INTEGER :: imn, jmn  ! dimensions nouvelle grille
36  INTEGER :: kllm ! taille du tableau des intersections
37  REAL :: rlonuo(imo + 1)     !  Latitude et
38  REAL :: rlatvo(jmo)       !  longitude des
39  REAL :: rlonun(imn + 1)     !  bord des
40  REAL :: rlatvn(jmn)     !  cases "scalaires" (input)
41
42  ! OUTPUT
43  INTEGER :: ktotal ! nombre totale d'intersections reperees
44  INTEGER :: iik(kllm), jjk(kllm), jk(kllm), ik(kllm)
45  REAL :: intersec(kllm)  ! surface des intersections (m2)
46  REAL :: airen (imn + 1, jmn + 1) ! aire dans la nouvelle grille
47
48
49
50
51  ! Autres variables
52  ! """"""""""""""""
53  INTEGER :: i, j, ii, jj
54  REAL :: a(imo + 1), b(imo + 1), c(jmo + 1), d(jmo + 1)
55  REAL :: an(imn + 1), bn(imn + 1), cn(jmn + 1), dn(jmn + 1)
56  REAL :: aa, bb, cc, dd
57  REAL :: pi
58
59  pi = 2. * ASIN(1.)
60
61
62
63  ! On repere les frontieres des cases :
64  ! ===================================
65
66  ! Attention, on ruse avec des latitudes = 90 deg au pole.
67
68
69  !  ANcienne grile
70  !  """"""""""""""
71  a(1) = - rlonuo(imo + 1)
72  b(1) = rlonuo(1)
73  DO i = 2, imo + 1
74    a(i) = rlonuo(i - 1)
75    b(i) = rlonuo(i)
76  END DO
77
78  d(1) = pi / 2
79  DO j = 1, jmo
80    c(j) = rlatvo(j)
81    d(j + 1) = rlatvo(j)
82  END DO
83  c(jmo + 1) = -pi / 2
84
85
86  !  Nouvelle grille
87  !  """""""""""""""
88  an(1) = - rlonun(imn + 1)
89  bn(1) = rlonun(1)
90  DO i = 2, imn + 1
91    an(i) = rlonun(i - 1)
92    bn(i) = rlonun(i)
93  END DO
94
95  dn(1) = pi / 2
96  DO j = 1, jmn
97    cn(j) = rlatvn(j)
98    dn(j + 1) = rlatvn(j)
99  END DO
100  cn(jmn + 1) = -pi / 2
101
102  ! Calcul de la surface des cases scalaires de la nouvelle grille
103  ! ==============================================================
104  DO ii = 1, imn + 1
105    DO jj = 1, jmn + 1
106      airen(ii, jj) = (bn(ii) - an(ii)) * (sin(dn(jj)) - sin(cn(jj)))
107    END DO
108  END DO
109
110  ! Calcul de la surface des intersections
111  ! ======================================
112
113  ! boucle sur la nouvelle grille
114  ! """"""""""""""""""""""""""""
115  ktotal = 0
116  DO jj = 1, jmn + 1
117    DO j = 1, jmo + 1
118      IF((cn(jj)<d(j)).AND.(dn(jj)>c(j)))THEN
119        DO ii = 1, imn + 1
120          DO i = 1, imo + 1
121            IF (((an(ii)<b(i)).AND.(bn(ii)>a(i))) &
122                    .OR. ((an(ii)<b(i) - 2 * pi).AND.(bn(ii)>a(i) - 2 * pi) &
123                            .AND.(b(i) - 2 * pi<-pi)) &
124                    .OR. ((an(ii)<b(i) + 2 * pi).AND.(bn(ii)>a(i) + 2 * pi) &
125                            .AND.(a(i) + 2 * pi>pi)) &
126                    )THEN
127              ktotal = ktotal + 1
128              iik(ktotal) = ii
129              jjk(ktotal) = jj
130              ik(ktotal) = i
131              jk(ktotal) = j
132
133              dd = min(d(j), dn(jj))
134              cc = cn(jj)
135              IF (cc< c(j))cc = c(j)
136              IF((an(ii)<b(i) - 2 * pi).AND. &
137                      (bn(ii)>a(i) - 2 * pi)) THEN
138                bb = min(b(i) - 2 * pi, bn(ii))
139                aa = an(ii)
140                IF (aa<a(i) - 2 * pi) aa = a(i) - 2 * pi
141              ELSE IF((an(ii)<b(i) + 2 * pi).AND. &
142                      (bn(ii)>a(i) + 2 * pi)) THEN
143                bb = min(b(i) + 2 * pi, bn(ii))
144                aa = an(ii)
145                IF (aa<a(i) + 2 * pi) aa = a(i) + 2 * pi
146              else
147                bb = min(b(i), bn(ii))
148                aa = an(ii)
149                IF (aa<a(i)) aa = a(i)
150              end if
151              intersec(ktotal) = (bb - aa) * (sin(dd) - sin(cc))
152            end if
153          END DO
154        END DO
155      end if
156    END DO
157  END DO
158END SUBROUTINE  iniinterp_horiz
Note: See TracBrowser for help on using the repository browser.