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

Last change on this file since 1458 was 1403, checked in by emillour, 10 years ago

All models: Reorganizing the physics/dynamics interface.

  • makelmdz and makelmdz_fcm scripts adapted to handle the new directory settings
  • misc: (replaces what was the "bibio" directory)
  • Should only contain extremely generic (and non physics or dynamics-specific) routines
  • Therefore moved initdynav.F90, initfluxsto.F, inithist.F, writedynav.F90, write_field.F90, writehist.F to "dyn3d_common"
  • dynlonlat_phylonlat: (new interface directory)
  • This directory contains routines relevent to physics/dynamics grid interactions, e.g. routines gr_dyn_fi or gr_fi_dyn and calfis
  • Moreover the dynlonlat_phylonlat contains directories "phy*" corresponding to each physics package "phy*" to be used. These subdirectories should only contain specific interfaces (e.g. iniphysiq) or main programs (e.g. newstart)
  • phy*/dyn1d: this subdirectory contains the 1D model using physics from phy*

EM

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