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

Last change on this file since 1243 was 1216, checked in by emillour, 11 years ago

Generic model:
Major cleanup, in order to ease the use of LMDZ.GENERIC with (parallel) dynamics
in LMDZ.COMMON: (NB: this will break LMDZ.UNIVERSAL, which should be thrashed
in the near future)

  • Updated makegcm_* scripts (and makdim) and added the "-full" (to enforce full recomputation of the model) option
  • In dyn3d: converted control.h to module control_mod.F90 and converted iniadvtrac.F to module infotrac.F90
  • Added module mod_const_mpi.F90 in dyn3d (not used in serial mode)
  • Rearanged input/outputs routines everywhere to handle serial/MPI cases. physdem.F => phyredem.F90 , phyetat0.F => phyetat0.F90 ; all read/write routines for startfi files are gathered in module iostart.F90
  • added parallelism related routines init_phys_lmdz.F90, comgeomphy.F90, dimphy.F90, iniphysiq.F90, mod_grid_phy_lmdz.F90, mod_phys_lmdz_mpi_data.F90, mod_phys_lmdz_mpi_transfert.F90, mod_phys_lmdz_omp_data.F90, mod_phys_lmdz_omp_transfert.F90, mod_phys_lmdz_para.F90, mod_phys_lmdz_transfert_para.F90 in phymars and mod_const_mpi.F90 in dyn3d (for compliance with parallelism)
  • added created generic routines 'planetwide_maxval' and 'planetwide_minval', in module "planetwide_mod", that enable obtaining the max and min of a field over the whole planet. This should be further imroved with computation of means (possibly area weighed), etc.

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