source: LMDZ.3.3/branches/rel-1-0-patch/libf/dyn3d/interp_horizFF @ 428

Last change on this file since 428 was 2, checked in by lmdz, 25 years ago

Initial revision

  • Property svn:eol-style set to native
  • Property svn:executable set to *
  • Property svn:keywords set to Author Date Id Revision
File size: 4.0 KB
Line 
1      subroutine interp_horiz (varo,varn,imo,jmo,imn,jmn,lm,
2     &  rlonuo,rlatvo,rlonun,rlatvn) 
3
4c===========================================================
5c  Interpolation Horizontales des variables d'une grille LMDZ
6c (des points SCALAIRES au point SCALAIRES)
7c  dans une autre grille LMDZ en conservant la quantite
8c  totale pour les variables intensives (/m2) : ex : Pression au sol
9c
10c Francois Forget (01/1995)
11c===========================================================
12
13      IMPLICIT NONE
14
15c   Declarations:
16c ==============
17c
18c  ARGUMENTS
19c  """""""""
20       
21       integer imo, jmo ! dimensions ancienne grille (input)
22       integer imn,jmn  ! dimensions nouvelle grille (input)
23
24       real rlonuo(imo+1)     !  Latitude et
25       real rlatvo(jmo)       !  longitude des
26       real rlonun(imn+1)     !  bord des
27       real rlatvn(jmn)     !  cases "scalaires" (input)
28
29       integer lm ! dimension verticale (input)
30       real varo (imo+1, jmo+1,lm) ! var dans l'ancienne grille (input)
31       real varn (imn+1,jmn+1,lm) ! var dans la nouvelle grille (output)
32
33c Autres variables
34c """"""""""""""""
35       real airetest(imn+1,jmn+1)
36       integer ii,jj,l
37
38       real airen (imn+1,jmn+1) ! aire dans la nouvelle grille
39c    Info sur les ktotal intersection entre les cases new/old grille
40       integer kllm, k, ktotal
41       parameter (kllm = 100*100*10)
42       integer iik(kllm), jjk(kllm),jk(kllm),ik(kllm)
43       real intersec(kllm)
44       real R
45       real totn, tots
46
47       logical firstcall, firsttest, aire_ok
48       save firsttest
49       data firsttest /.true./
50       data aire_ok /.true./
51
52       
53
54
55
56c initialisation
57c --------------
58c Si c'est le premier appel, on prepare l'interpolation
59c en calculant pour chaque case autour d'un point scalaire de la
60c nouvelle grille, la surface  de intersection avec chaque
61c    case de l'ancienne grille.
62
63
64        call iniinterp_horiz (imo,jmo,imn,jmn ,kllm,
65     &       rlonuo,rlatvo,rlonun,rlatvn,
66     &          ktotal,iik,jjk,jk,ik,intersec,airen)
67
68      do l=1,lm
69       do jj =1 , jmn+1
70        do ii=1, imn+1
71          varn(ii,jj,l) =0.
72        end do
73       end do
74      end do
75       
76c Interpolation horizontale
77c -------------------------
78c boucle sur toute les ktotal intersections entre les cases
79c de l'ancienne et la  nouvelle grille
80c
81     
82      do k=1,ktotal
83        do l=1,lm
84         varn(iik(k),jjk(k),l) = varn(iik(k),jjk(k),l)
85     &        + varo(ik(k), jk(k),l)*intersec(k)/airen(iik(k),jjk(k))
86        end do
87      end do
88
89c Une seule valeur au pole pour les variables ! :
90c -----------------------------------------------
91       do l=1, lm
92         totn =0.
93         tots =0.
94           do ii =1, imn+1
95             totn = totn + varn(ii,1,l)
96             tots = tots + varn (ii,jmn+1,l)
97           end do
98           do ii =1, imn+1
99             varn(ii,1,l) = totn/float(imn+1)
100             varn(ii,jmn+1,l) = tots/float(imn+1)
101           end do
102       end do
103           
104
105c---------------------------------------------------------------
106c  TEST  TEST  TEST  TEST  TEST  TEST  TEST  TEST  TEST  TEST
107      if (.not.(firsttest)) goto 99
108      firsttest = .false.
109      write (*,*) 'INTERP. HORIZ. : TEST SUR LES AIRES:'
110      do jj =1 , jmn+1
111        do ii=1, imn+1
112          airetest(ii,jj) =0.
113        end do
114      end do
115      do k=1,ktotal
116         airetest(iik(k),jjk(k))= airetest(iik(k),jjk(k)) +intersec(k)
117      end do
118      do jj =1 , jmn+1
119       do ii=1, imn+1
120         r = airen(ii,jj)/airetest(ii,jj)
121         if ((r.gt.1.001).or.(r.lt.0.999)) then
122             write (*,*) '********** PROBLEME D'' AIRES !!!',
123     &                   ' DANS L''INTERPOLATION HORIZONTALE'
124             write(*,*)'ii,jj,airen,airetest',
125     &          ii,jj,airen(ii,jj),airetest(ii,jj)
126             aire_ok = .false.
127         end if
128       end do
129      end do
130      if (aire_ok) write(*,*) 'INTERP. HORIZ. : AIRES OK'
131 99   continue
132
133c FIN TEST  FIN TEST  FIN TEST  FIN TEST  FIN TEST  FIN TEST  FIN TEST
134c---------------------------------------------------------------
135
136
137
138
139
140
141
142
143        return
144        end
Note: See TracBrowser for help on using the repository browser.