source: trunk/LMDZ.GENERIC/libf/phystd/newtrelax.F90 @ 1351

Last change on this file since 1351 was 1315, checked in by milmd, 10 years ago

LMDZ.GENERIC. OpenMP directives added in generic physic. When running in pure OpenMP or hybrid OpenMP/MPI, may have some bugs with condense_cloud and wstats routines.

File size: 3.6 KB
Line 
1subroutine newtrelax(ngrid,nlayer,mu0,sinlat,popsk,temp,pplay,pplev,dtrad,firstcall)
2       
3  implicit none
4
5!#include "dimensions.h"
6!#include "dimphys.h"
7#include "comcstfi.h"
8#include "callkeys.h"
9#include "netcdf.inc"
10
11!==================================================================
12!     
13!     Purpose
14!     -------
15!     Alternative Newtonian radiative transfer scheme.
16!     
17!     Authors
18!     -------
19!     R. Wordsworth (2010)
20!     
21!==================================================================
22 
23 
24  ! Input
25  integer,intent(in) :: ngrid, nlayer
26  logical,intent(in) :: firstcall
27  real,intent(in) :: mu0(ngrid)            ! cosine of sun incident angle
28  real,intent(in) :: sinlat(ngrid)         ! sine of latitude
29  real,intent(in) :: temp(ngrid,nlayer)    ! temperature at each layer (K)
30  real,intent(in) :: pplay(ngrid,nlayer)   ! pressure at each layer (Pa)
31  real,intent(in) :: pplev(ngrid,nlayer+1) ! pressure at each level (Pa)
32  real,intent(in) :: popsk(ngrid,nlayer)   ! pot. T to T converter
33
34  ! Output
35  real,intent(out) :: dtrad(ngrid,nlayer)
36
37  ! Internal
38  real Trelax_V, Trelax_H
39  real,allocatable,dimension(:,:),save :: Trelax
40!$OMP THREADPRIVATE(Trelax)
41
42  real T_trop ! relaxation temperature at tropopause (K)
43  real T_surf ! relaxation temperature at surface (K)
44  real dT_EP  ! Equator-Pole relaxation temperature difference (K)
45
46  real sig, f_sig, sig_trop
47  integer l,ig
48
49
50  logical tidallocked
51  parameter (tidallocked = .true.)
52
53  ! Setup relaxation temperature 
54  if(firstcall) then
55
56     ALLOCATE(Trelax(ngrid,nlayer))
57
58     print*,'-----------------------------------------------------'
59     print*,'| ATTENTION: You are using a Newtonian cooling scheme'
60     print*,'| for the radiative transfer. This means that ALL'
61     print*,'| other physics subroutines must be switched off.'
62     print*,'-----------------------------------------------------'
63
64     if(tidallocked)then
65        do ig=1,ngrid
66
67           T_surf = 126. + 239.*mu0(ig)
68           T_trop = 140. + 52.*mu0(ig)
69           do l=1,nlayer
70
71              if(mu0(ig).le.0.0)then ! night side
72                 Trelax(ig,l)=0.0
73              else                   ! day side
74                 Trelax(ig,l) = T_surf*popsk(ig,l)
75                 if (Trelax(ig,l).lt.T_trop) Trelax(ig,l) = T_trop
76              endif
77
78           enddo
79        enddo
80
81     else
82
83        T_trop = 200.
84        T_surf = 288.
85        dT_EP  = 70.
86
87        sig_trop=(T_trop/T_surf)**(1./rcp)
88
89        do l=1,nlayer
90           do ig=1,ngrid
91
92              ! vertically varying component
93              Trelax_V = T_surf*popsk(ig,l)
94              if (Trelax_V.lt.T_trop) Trelax_V = T_trop
95             
96              ! horizontally varying component
97              sig = pplay(ig,l)/pplev(ig,1)
98              if(sig.ge.sig_trop)then
99                 f_sig=sin((pi/2)*((sig-sig_trop)/(1-sig_trop)))
100              else
101                 f_sig=0.0
102              endif
103              Trelax_H = -f_sig*dT_EP*(sinlat(ig)**2 - 1./3.)
104             
105              Trelax(ig,l) = Trelax_V + Trelax_H           
106           
107           enddo
108        enddo
109
110     endif
111
112  endif
113
114  ! Calculate radiative forcing
115  do l=1,nlayer
116     do ig=1,ngrid
117        dtrad(ig,l) = -(temp(ig,l) - Trelax(ig,l)) / tau_relax
118        if(temp(ig,l).gt.500.)then ! Trelax(ig,l))then
119           print*,'ig=',ig
120           print*,'l=',l
121           print*,'temp=',temp(ig,l)
122           print*,'Trelax=',Trelax(ig,l)
123        endif
124     enddo
125  enddo
126
127  call writediagfi(ngrid,'Tref','rad forc temp','K',3,Trelax)
128  !call writediagfi(ngrid,'ThetaZ','stellar zenith angle','deg',2,mu0)
129
130end subroutine newtrelax
Note: See TracBrowser for help on using the repository browser.