source: trunk/libf/phyvenus/load_psi.F @ 4

Last change on this file since 4 was 3, checked in by slebonnois, 15 years ago

Creation de repertoires:

  • chantiers : pour communiquer sur nos projets de modifs
  • documentation : pour stocker les docs

Ajout de:

  • libf/phytitan : physique de Titan
  • libf/chimtitan: chimie de Titan
  • libf/phyvenus : physique de Venus
File size: 4.2 KB
Line 
1      SUBROUTINE load_psi(
2     S              psurf, ztop, ksive,
3     S              temp, psimap, deltapsimap)
4     
5      IMPLICIT none
6
7#include "dimensions.h"
8#include "dimphy.h"
9#include "raddim.h"
10#include "YOMCST.h"
11#include "comcstVE.h"
12C
13C     ------------------------------------------------------------------
14C
15C     PURPOSE.
16C     --------
17C
18c     This routine loads the longwave matrix of factors Ksi
19c     (interpolated for a given column)
20c     to build the Net Exchange Rates matrix Psi.
21c     Psi(i,j,nu) = Ksi(i,j,nu) * ( B(i,nu)-B(j,nu) )
22c     
23c     The Ksi matrixes have been computed by Vincent Eymet
24C
25c     The NER matrix is then integrated in frequency.
26c
27C     AUTHOR.
28C     -------
29C        Sebastien Lebonnois
30C
31C     MODIFICATIONS.
32C     --------------
33C        version multimatrice (topographie, sommet nuages): 20/12/2006
34C     ------------------------------------------------------------------
35C
36C* ARGUMENTS:
37C
38c inputs
39      real   psurf(klon)           ! Surface pressure
40      real   ztop(klon)            ! Altitude of the top of cloud deck (km)
41      real   ksive(0:kflev+1,0:kflev+1,nnuve,nbmat)  ! ksi matrixes in Vincent's file
42      real   temp(klon,0:kflev+1)  ! Temperature in layer (K)
43c outputs
44      real   psimap(0:kflev+1,0:kflev+1,klon)
45      real   deltapsimap(0:kflev+1,0:kflev+1,klon)
46
47c local variables
48      integer i,j,ig,band
49      integer mat,m,mat0
50      character*100 file
51      real   bplck(0:nlve+1,nnuve)    ! Planck luminances in table layers
52      real   y(0:nlve,nnuve)          ! intermediaire Planck
53      real   zdblay(0:nlve+1,nnuve)   ! gradient en temperature de planck
54      real   ksi
55      real   factp,factz
56
57c -----------------------
58c Main loop on grid point
59c -----------------------
60      do ig=1,klon
61
62c Planck function
63c ---------------
64
65       do band=1,nnuve
66         do j=0,nlve
67c B(T,l) = al/(exp(bl/T)-1)
68           y(j,band) = exp(bl(band)/temp(ig,j))-1.
69           bplck(j,band) = al(band)/(y(j,band))
70           zdblay(j,band) = al(band)*bl(band)*exp(bl(band)/temp(ig,j))/
71     .                 ((temp(ig,j)*temp(ig,j))*(y(j,band)*y(j,band)))
72         enddo
73         bplck(nlve+1,band) = 0.0
74         zdblay(nlve+1,band)= 0.0
75       enddo
76
77c interpolating ksi
78c  and
79c computing psi and deltapsi
80c ---------------------------
81
82c init
83       do j=0,nlve+1
84        do i=0,nlve+1
85         psimap(i,j,ig) = 0.0  ! positif quand nrj de i->j
86         deltapsimap(i,j,ig) = 0.0
87        enddo
88       enddo
89
90c finding the right matrixes
91       mat0 = 0
92       do mat=1,nbmat-nbztopve
93         if (  (psurfve(mat).ge.psurf(ig))
94     .    .and.(psurfve(mat+nbztopve).lt.psurf(ig))
95     .    .and.(ztopve(mat).lt.ztop(ig))
96     .    .and.(ztopve(mat+1).ge.ztop(ig)) ) then
97              mat0  = mat
98c             print*,'ig=',ig,'  mat0=',mat
99              factp = (psurf(ig)            -psurfve(mat))
100     .               /(psurfve(mat+nbztopve)-psurfve(mat))
101              factz = (ztop(ig)     -ztopve(mat))
102     .               /(ztopve(mat+1)-ztopve(mat))
103         endif
104       enddo
105       if (mat0.eq.0) then
106         write(*,*) 'This is subroutine load_psi'
107         print*,'Probleme pour interpolation au point ig=',ig
108         print*,'psurf = ',psurf(ig),' ztop = ',ztop(ig)
109         stop
110       endif
111
112c------------------
113c !!TEST!! Matrice unique fixee: psurf = 90 bar, ztop = 70
114c      mat0 = 24
115c      print*,'MATRICE UNIQUE: ',mat0,' / ps=',psurfve(mat0),
116c    .                                ' / ztop=',ztopve(mat0)
117c------------------
118
119c interpolation of ksi and computation of psi,deltapsi
120       do band=1,nnuve
121        do j=0,nlve+1
122         do i=0,nlve+1
123          ksi = ksive(i,j,band,mat0)*(1-factz)*(1-factp)
124     .         +ksive(i,j,band,mat0+1)*factz  *(1-factp)
125     .         +ksive(i,j,band,mat0+nbztopve)*(1-factz)*factp
126     .         +ksive(i,j,band,mat0+nbztopve+1)*factz  *factp
127         psimap(i,j,ig) = psimap(i,j,ig) +
128     .          ksi*(bplck(i,band)-bplck(j,band))
129         deltapsimap(i,j,ig) = deltapsimap(i,j,ig) +
130     .          ksi*zdblay(i,band)
131         enddo
132        enddo
133       enddo
134
135      enddo !ig
136c -----------------------
137c End loop on grid point
138c -----------------------
139
140c     print*,"LOAD_PSI OK"
141
142      return
143      end
144
Note: See TracBrowser for help on using the repository browser.