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

Last change on this file since 1351 was 1308, checked in by emillour, 10 years ago

Generic GCM:
Some cleanup to simplify dynamics/physics interactions by getting rid
of dimphys.h (i.e. the nlayermx parameter) and minimizing use of
dimension.h in the physics.
EM

  • Property svn:executable set to *
File size: 10.3 KB
Line 
1program kcm1d
2
3  use infotrac, only: nqtot
4  use radinc_h,      only: NAERKIND
5  use radcommon_h,   only: L_NSPECTI, L_NSPECTV, sigma
6  use watercommon_h, only: mH2O
7  use ioipsl_getincom, only: getin
8  use comsaison_h, only: mu0, fract, dist_star
9  use planete_mod
10!  use control_mod
11
12  implicit none
13
14  !==================================================================
15  !     
16  !     Purpose
17  !     -------
18  !     Run the universal model radiative transfer once in a 1D column.
19  !     Useful for climate evolution studies etc.
20  !     
21  !     It can be compiled with a command like (e.g. for 25 layers):
22  !                                  "makegcm -p std -d 25 kcm1d"
23  !     It requires the files "callphys.def", "gases.def"
24  !     "traceur.def", and "run.def" with a line "INCLUDEDEF=callphys.def"
25  !
26  !     Authors
27  !     -------
28  !     R. Wordsworth
29  !
30  !==================================================================
31
32#include "dimensions.h"
33!#include "dimphys.h"
34#include "callkeys.h"
35#include "comcstfi.h"
36!#include "planete.h"
37!#include "control.h"
38
39  ! --------------------------------------------------------------
40  !  Declarations
41  ! --------------------------------------------------------------
42
43  integer nlayer,nlevel,nq
44  integer ilay,ilev,iq,iw,iter
45  real play(llm)     ! Pressure at the middle of the layers [Pa]
46  real zlay(llm)     ! Altitude at middle of the layers [km]
47  real plev(llm+1)   ! Intermediate pressure levels [Pa]
48  real temp(llm)     ! temperature at the middle of the layers [K]
49  real,allocatable :: q(:,:)   ! tracer mixing ratio [kg/kg]
50  real,allocatable :: vmr(:,:) ! tracer mixing ratio [mol/mol]
51  real,allocatable :: qsurf(:)        ! tracer surface budget [kg/kg] ????
52  real psurf,psurf_n,tsurf     
53
54  real emis, albedo
55
56  real muvar(llm+1)
57
58  real dtsw(llm) ! heating rate (K/s) due to SW
59  real dtlw(llm) ! heating rate (K/s) due to LW
60  real fluxsurf_lw   ! incident LW flux to surf (W/m2)
61  real fluxtop_lw    ! outgoing LW flux to space (W/m2)
62  real fluxsurf_sw   ! incident SW flux to surf (W/m2)
63  real fluxabs_sw    ! SW flux absorbed by planet (W/m2)
64  real fluxtop_dn    ! incident top of atmosphere SW flux (W/m2)
65
66  ! not used
67  real reffrad(llm,naerkind)
68  real nueffrad(llm,naerkind)
69  real cloudfrac(llm)
70  real totcloudfrac
71  real tau_col
72
73  real dTstrat
74  real aerosol(llm,naerkind) ! aerosol tau (kg/kg)
75  real OLR_nu(1,L_NSPECTI)
76  real OSR_nu(1,L_NSPECTV)
77  real Eatmtot
78
79  integer ierr
80  logical firstcall,lastcall,global1d
81
82  character*20,allocatable :: nametrac(:)   ! name of the tracer (no need for adv trac common)
83
84  ! --------------
85  ! Initialisation
86  ! --------------
87
88  pi=2.E+0*asin(1.E+0)
89
90  reffrad(:,:)  = 0.0
91  nueffrad(:,:) = 0.0
92  cloudfrac(:)  = 0.0
93  totcloudfrac  = 0.0
94
95
96  nlayer=llm
97  nlevel=nlayer+1
98
99  !! this is defined in comsaison_h
100  ALLOCATE(mu0(1))
101  ALLOCATE(fract(1))
102
103
104
105  !  Load parameters from "run.def"
106  ! -------------------------------
107
108  ! check if 'kcm1d.def' file is around (otherwise reading parameters
109  ! from callphys.def via getin() routine won't work.)
110  open(90,file='kcm1d.def',status='old',form='formatted',&
111       iostat=ierr)
112  if (ierr.ne.0) then
113     write(*,*) 'Cannot find required file "kcm1d.def"'
114     write(*,*) '  (which should contain some input parameters'
115     write(*,*) '   along with the following line:'
116     write(*,*) '   INCLUDEDEF=callphys.def'
117     write(*,*) '   )'
118     write(*,*) ' ... might as well stop here ...'
119     stop
120  else
121     close(90)
122  endif
123
124! now, run.def is needed anyway. so we create a dummy temporary one
125! for ioipsl to work. if a run.def is already here, stop the
126! program and ask the user to do a bit of cleaning
127  open(90,file='run.def',status='old',form='formatted',&
128       iostat=ierr)
129  if (ierr.eq.0) then
130     close(90)
131     write(*,*) 'There is already a run.def file.'
132     write(*,*) 'This is not compatible with 1D runs.'
133     write(*,*) 'Please remove the file and restart the run.'
134     write(*,*) 'Runtime parameters are supposed to be in kcm1d.def'
135     stop
136  else
137     call system('touch run.def')
138     call system("echo 'INCLUDEDEF=callphys.def' >> run.def")
139     call system("echo 'INCLUDEDEF=kcm1d.def' >> run.def")
140  endif
141
142  global1d = .false. ! default value
143  call getin("global1d",global1d)
144  if(.not.global1d)then
145     print*,'Error, kcm1d must have global1d=.true. in kcm1d.def!'
146     stop
147  end if
148
149  psurf_n=100000. ! default value for psurf
150  print*,'Dry surface pressure (Pa)?'
151  call getin("psurf",psurf_n)
152  write(*,*) " psurf = ",psurf_n
153
154! OK. now that run.def has been read once -- any variable is in memory.
155! so we can dump the dummy run.def
156  call system("rm -rf run.def")
157
158  tsurf=300.0 ! default value for tsurf
159  print*,'Surface temperature (K)?'
160  call getin("tref",tsurf)
161  write(*,*) " tsurf = ",tsurf
162
163  g=10.0 ! default value for g
164  print*,'Gravity ?'
165  call getin("g",g)
166  write(*,*) " g = ",g
167
168  periastr = 1.0
169  apoastr  = 1.0
170  print*,'Periastron (AU)?'
171  call getin("periastr",periastr)
172  write(*,*) "periastron = ",periastr
173  dist_star = periastr
174  fract     = 0.5
175  print*,'Apoastron (AU)?'
176  call getin("apoastr",apoastr)
177  write(*,*) "apoastron = ",apoastr
178
179  albedo=0.2 ! default value for albedo
180  print*,'Albedo of bare ground?'
181  call getin("albedo",albedo)
182  write(*,*) " albedo = ",albedo
183
184  emis=1.0 ! default value for emissivity
185  PRINT *,'Emissivity of bare ground ?'
186  call getin("emis",emis)
187  write(*,*) " emis = ",emis
188
189  pceil=100.0 ! Pascals
190  PRINT *,'Ceiling pressure (Pa) ?'
191  call getin("pceil",pceil)
192  write(*,*) " pceil = ", pceil
193
194  mugaz=0.
195  cpp=0.
196
197  check_cpp_match = .false.
198  call getin("check_cpp_match",check_cpp_match)
199  if (check_cpp_match) then
200     print*,"In 1D modeling, check_cpp_match is supposed to be F"
201     print*,"Please correct callphys.def"
202     stop
203  endif
204
205  call su_gases
206  call calc_cpp_mugaz
207
208  call inifis(1,llm,0,86400.0,1.0,0.0,0.0,1.0,rad,g,r,cpp)
209
210  ! Tracer initialisation
211  ! ---------------------
212  if (tracer) then
213     ! load tracer names from file 'traceur.def'
214     open(90,file='traceur.def',status='old',form='formatted',&
215          iostat=ierr)
216     if (ierr.eq.0) then
217        write(*,*) "kcm1d: Reading file traceur.def"
218        ! read number of tracers:
219        read(90,*,iostat=ierr) nq
220        if (ierr.ne.0) then
221           write(*,*) "kcm1d: error reading number of tracers"
222           write(*,*) "   (first line of traceur.def) "
223           stop
224        endif
225        nqtot=nq
226        ! allocate arrays which depend on number of tracers
227        allocate(nametrac(nq))
228        allocate(q(nlayer,nq))
229        allocate(vmr(nlayer,nq))
230        allocate(qsurf(nq))
231
232        do iq=1,nq
233           ! minimal version, just read in the tracer names, 1 per line
234           read(90,*,iostat=ierr) nametrac(iq)
235           if (ierr.ne.0) then
236              write(*,*) 'kcm1d: error reading tracer names...'
237              stop
238           endif
239        enddo !of do iq=1,nq
240     endif
241
242     call initracer(1,nq,nametrac)
243
244  endif
245
246
247  do iq=1,nq
248     do ilay=1,nlayer
249        q(ilay,iq) = 0.
250     enddo
251  enddo
252
253  do iq=1,nq
254     qsurf(iq) = 0.
255  enddo
256
257  firstcall = .true.
258  lastcall  = .false.
259
260  iter    = 1
261  Tstrat  = 150.0
262  dTstrat = 1000.0
263
264  ! ---------
265  ! Run model
266  ! ---------
267  !do
268     psurf = psurf_n
269
270     !    Create vertical profiles
271     call kcmprof_fn(nlayer,psurf,qsurf(1),tsurf,    &
272          tstrat,play,plev,zlay,temp,q(:,1),muvar(1))
273
274     !    Run radiative transfer
275     call callcorrk(1,nlayer,q,nq,qsurf,      &
276          albedo,emis,mu0,plev,play,temp,                    &
277          tsurf,fract,dist_star,aerosol,muvar,         &
278          dtlw,dtsw,fluxsurf_lw,fluxsurf_sw,fluxtop_lw,    &
279          fluxabs_sw,fluxtop_dn,OLR_nu,OSR_nu,reffrad,nueffrad,tau_col,  &
280          cloudfrac,totcloudfrac,.false.,firstcall,lastcall)
281
282     !    Iterate stratospheric temperature
283     print*,'Tstrat = ',Tstrat
284     dTstrat = Tstrat
285     !Tstrat  = Tsurf*(0.2786*(psurf/100000.)**(-1.123))**0.25
286     ! skin temperature (gray approx.) using analytic pure H2 expression
287     !Tstrat  = (fluxabs_sw/(2*sigma))**0.25 ! skin temperature (gray approx.)
288     Tstrat  = (fluxtop_lw/(2*sigma))**0.25 ! skin temperature (gray approx.)
289     dTstrat = dTstrat-Tstrat
290
291     !if(abs(dTstrat).lt.1.0)then
292     !   print*,'dTstrat = ',dTstrat
293     !   exit
294     !endif
295
296     !iter=iter+1
297     !if(iter.eq.100)then
298     !   print*,'Stratosphere failed to converge after'
299     !   print*,'100 iteration, aborting run.'
300     !   call abort
301     !endif
302
303  !end do
304
305  ! Run radiative transfer one last time to get OLR,OSR
306  firstcall=.false.
307  lastcall=.true.
308  call callcorrk(1,nlayer,q,nq,qsurf,      &
309       albedo,emis,mu0,plev,play,temp,                    &
310       tsurf,fract,dist_star,aerosol,muvar,         &
311       dtlw,dtsw,fluxsurf_lw,fluxsurf_sw,fluxtop_lw,    &
312       fluxabs_sw,fluxtop_dn,OLR_nu,OSR_nu,             &
313       reffrad,nueffrad,tau_col,  &
314       cloudfrac,totcloudfrac,.false.,firstcall,lastcall)
315
316
317  ! Calculate total atmospheric energy
318  Eatmtot=0.0
319  !  call calcenergy_kcm(nlayer,tsurf,temp,play,plev,qsurf,&
320  !     q(:,1),muvar,Eatmtot)
321
322  ! ------------------------
323  ! Save data to ascii files
324  ! ------------------------
325
326  print*,'Saving profiles...'
327  open(115,file='profpres.out',form='formatted')
328  open(116,file='proftemp.out',form='formatted')
329  open(117,file='profztab.out',form='formatted')
330  open(118,file='profqvar.out',form='formatted')
331  open(119,file='profvvar.out',form='formatted')
332
333  write(115,*) psurf
334  write(116,*) tsurf
335  write(117,*) 0.0
336  write(118,*) qsurf(1)
337  write(119,*) qsurf(1)*(muvar(1)/mH2O)
338  do ilay=1,nlayer
339     vmr(ilay,1) = q(ilay,1)*(muvar(ilay+1)/mH2O)
340     write(115,*) play(ilay)
341     write(116,*) temp(ilay)
342     write(117,*) zlay(ilay)
343     write(118,*) q(ilay,1)
344     write(119,*) vmr(ilay,1)
345  enddo
346  close(115)
347  close(116)
348  close(117)
349  close(118)
350  close(119)
351
352  print*, tsurf,psurf,fluxtop_dn,fluxabs_sw,fluxtop_lw
353
354  print*,'Saving scalars...'
355  open(116,file='surf_vals.out')
356  write(116,*) tsurf,psurf,fluxtop_dn,         &
357       fluxabs_sw,fluxtop_lw
358  close(116)
359  open(111,file='ene_vals.out')
360  write(111,*) tsurf,psurf,Eatmtot,Tstrat
361  close(111)
362
363  print*,'Saving spectra...'
364  open(117,file='OLRnu.out')
365  do iw=1,L_NSPECTI
366     write(117,*) OLR_nu(1,iw)
367  enddo
368  close(117)
369
370  open(127,file='OSRnu.out')
371  do iw=1,L_NSPECTV
372     write(127,*) OSR_nu(1,iw)
373  enddo
374  close(127) 
375
376end program kcm1d
Note: See TracBrowser for help on using the repository browser.