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

Last change on this file since 2176 was 1714, checked in by mturbet, 7 years ago

kcm1d like a phoenix rises from the ashes, part I

  • Property svn:executable set to *
File size: 14.5 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, glat
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 callkeys_mod, only: check_cpp_match, pceil, tstrat, tracer
11  use inifis_mod, only: inifis
12  use comcstfi_mod
13  use mod_grid_phy_lmdz, only : regular_lonlat
14  use regular_lonlat_mod, only: init_regular_lonlat
15  use physics_distribution_mod, only: init_physics_distribution
16  use regular_lonlat_mod, only: init_regular_lonlat
17  use mod_interface_dyn_phys, only: init_interface_dyn_phys
18  use geometry_mod, only: init_geometry
19  use dimphy, only : init_dimphy
20
21  implicit none
22
23  !==================================================================
24  !     
25  !     Purpose
26  !     -------
27  !     Run the universal model radiative transfer once in a 1D column.
28  !     Useful for climate evolution studies etc.
29  !     
30  !     It can be compiled with a command like (e.g. for 25 layers):
31  !                                  "makegcm -p std -d 25 kcm1d"
32  !     It requires the files "callphys.def", "gases.def"
33  !     "traceur.def", and "run.def" with a line "INCLUDEDEF=callphys.def"
34  !
35  !     Authors
36  !     -------
37  !     - R. Wordsworth
38  !     - updated by M. Turbet (June 2017)
39  !
40  !==================================================================
41
42#include "dimensions.h"
43!#include "dimphys.h"
44
45  ! --------------------------------------------------------------
46  !  Declarations
47  ! --------------------------------------------------------------
48
49  integer nlayer,nlevel,nq
50  integer ilay,ilev,iq,iw,iter
51  real play(llm)     ! Pressure at the middle of the layers [Pa]
52  real zlay(llm)     ! Altitude at middle of the layers [km]
53  real plev(llm+1)   ! Intermediate pressure levels [Pa]
54  real temp(llm)     ! temperature at the middle of the layers [K]
55  real,allocatable :: q(:,:)   ! tracer mixing ratio [kg/kg]
56  real,allocatable :: vmr(:,:) ! tracer mixing ratio [mol/mol]
57  real,allocatable :: qsurf(:)        ! tracer surface budget [kg/kg] ????
58  real psurf,psurf_n,tsurf     
59
60  real emis, albedo, albedo_equivalent
61  real albedo_wv(1,L_NSPECTV)
62
63  real muvar(llm+1)
64
65  real dtsw(llm) ! heating rate (K/s) due to SW
66  real dtlw(llm) ! heating rate (K/s) due to LW
67  real fluxsurf_lw   ! incident LW flux to surf (W/m2)
68  real fluxtop_lw    ! outgoing LW flux to space (W/m2)
69  real fluxsurf_sw   ! incident SW flux to surf (W/m2)
70  real fluxsurfabs_sw ! absorbed SW flux by the surf (W/m2)
71  real fluxabs_sw    ! SW flux absorbed by planet (W/m2)
72  real fluxtop_dn    ! incident top of atmosphere SW flux (W/m2)
73
74  ! not used
75  real cloudfrac(llm)
76  real totcloudfrac
77  real tau_col
78
79  real dTstrat
80  real aerosol(llm,naerkind) ! aerosol tau (kg/kg)
81  real OLR_nu(1,L_NSPECTI)
82  real OSR_nu(1,L_NSPECTV)
83  real Eatmtot
84
85  integer ierr
86  logical firstcall,lastcall,global1d
87
88  character*20,allocatable :: nametrac(:)   ! name of the tracer (no need for adv trac common)
89 
90  real :: latitude(1), longitude(1), cell_area(1), phisfi(1)
91 
92  ! --------------
93  ! Initialisation
94  ! --------------
95
96  pi=2.E+0*asin(1.E+0)
97
98  cloudfrac(:)  = 0.0
99  totcloudfrac  = 0.0
100
101
102  nlayer=llm
103  nlevel=nlayer+1
104
105  !! this is defined in comsaison_h
106  ALLOCATE(mu0(1))
107  ALLOCATE(fract(1))
108  ALLOCATE(glat(1))
109
110  !  Load parameters from "run.def"
111  ! -------------------------------
112
113  ! check if 'kcm1d.def' file is around (otherwise reading parameters
114  ! from callphys.def via getin() routine won't work.)
115  open(90,file='kcm1d.def',status='old',form='formatted',&
116       iostat=ierr)
117  if (ierr.ne.0) then
118     write(*,*) 'Cannot find required file "kcm1d.def"'
119     write(*,*) '  (which should contain some input parameters'
120     write(*,*) '   along with the following line:'
121     write(*,*) '   INCLUDEDEF=callphys.def'
122     write(*,*) '   )'
123     write(*,*) ' ... might as well stop here ...'
124     stop
125  else
126     close(90)
127  endif
128 
129! now, run.def is needed anyway. so we create a dummy temporary one
130! for ioipsl to work. if a run.def is already here, stop the
131! program and ask the user to do a bit of cleaning
132  open(90,file='run.def',status='old',form='formatted',&
133       iostat=ierr)
134  if (ierr.eq.0) then
135     close(90)
136     write(*,*) 'There is already a run.def file.'
137     write(*,*) 'This is not compatible with 1D runs.'
138     write(*,*) 'Please remove the file and restart the run.'
139     write(*,*) 'Runtime parameters are supposed to be in kcm1d.def'
140     stop
141  else
142     call system('touch run.def')
143     call system("echo 'INCLUDEDEF=callphys.def' >> run.def")
144     call system("echo 'INCLUDEDEF=kcm1d.def' >> run.def")
145  endif
146   
147! check if we are going to run with or without tracers
148  write(*,*) "Run with or without tracer transport ?"
149  tracer=.false. ! default value
150  call getin("tracer",tracer)
151  write(*,*) " tracer = ",tracer
152
153
154
155  ! Tracer initialisation
156  ! ---------------------
157  if (tracer) then
158     ! load tracer names from file 'traceur.def'
159     open(90,file='traceur.def',status='old',form='formatted',&
160          iostat=ierr)
161     if (ierr.eq.0) then
162        write(*,*) "kcm1d: Reading file traceur.def"
163        ! read number of tracers:
164        read(90,*,iostat=ierr) nq
165        if (ierr.ne.0) then
166           write(*,*) "kcm1d: error reading number of tracers"
167           write(*,*) "   (first line of traceur.def) "
168           stop
169        endif
170        nqtot=nq
171        ! allocate arrays which depend on number of tracers
172        allocate(nametrac(nq))
173        allocate(q(nlayer,nq))
174        allocate(vmr(nlayer,nq))
175        allocate(qsurf(nq))
176
177        do iq=1,nq
178           ! minimal version, just read in the tracer names, 1 per line
179           read(90,*,iostat=ierr) nametrac(iq)
180           write(*,*) 'tracer here is', nametrac(iq)
181           if (ierr.ne.0) then
182              write(*,*) 'kcm1d: error reading tracer names...'
183              stop
184           endif
185        enddo !of do iq=1,nq
186     endif
187
188     call initracer(1,nq,nametrac)
189
190  endif
191 
192
193
194
195  global1d = .false. ! default value
196  call getin("global1d",global1d)
197  write(*,*) " global1d = ",global1d
198  if(.not.global1d)then
199     print*,'Error, kcm1d must have global1d=.true. in kcm1d.def!'
200     stop
201  end if
202
203  psurf_n=100000. ! default value for psurf
204  print*,'Dry surface pressure (Pa)?'
205  call getin("psurf",psurf_n)
206  write(*,*) " psurf = ",psurf_n
207
208! OK. now that run.def has been read once -- any variable is in memory.
209! so we can dump the dummy run.def
210  call system("rm -rf run.def")
211
212  tsurf=300.0 ! default value for tsurf
213  print*,'Surface temperature (K)?'
214  call getin("tref",tsurf)
215  write(*,*) " tsurf = ",tsurf
216
217  g=10.0 ! default value for g
218  print*,'Gravity ?'
219  call getin("g",g)
220  write(*,*) " g = ",g
221  glat(1)=g
222
223  periastr = 1.0
224  apoastr  = 1.0
225  print*,'Periastron (AU)?'
226  call getin("periastr",periastr)
227  write(*,*) "periastron = ",periastr
228  dist_star = periastr
229  fract     = 0.5
230  print*,'Apoastron (AU)?'
231  call getin("apoastr",apoastr)
232  write(*,*) "apoastron = ",apoastr
233
234  albedo=0.2 ! default value for albedo
235  print*,'Albedo of bare ground?'
236  call getin("albedo",albedo)
237  write(*,*) " albedo = ",albedo
238  do iw=1,L_NSPECTV
239     albedo_wv(1,iw)=albedo
240  enddo
241
242  emis=1.0 ! default value for emissivity
243  PRINT *,'Emissivity of bare ground ?'
244  call getin("emis",emis)
245  write(*,*) " emis = ",emis
246
247  pceil=100.0 ! Pascals
248  PRINT *,'Ceiling pressure (Pa) ?'
249  call getin("pceil",pceil)
250  write(*,*) " pceil = ", pceil
251
252  check_cpp_match = .false.
253  call getin("check_cpp_match",check_cpp_match)
254  if (check_cpp_match) then
255     print*,"In 1D modeling, check_cpp_match is supposed to be F"
256     print*,"Please correct callphys.def"
257     stop
258  endif
259
260
261!!! GEOGRAPHICAL INITIALIZATIONS
262     !!! AREA. useless in 1D
263  cell_area(1)=1.E+0 !JL+EM to have access to the area in the diagfi.nc files.
264     !!! GEOPOTENTIAL. useless in 1D because control by surface pressure
265  phisfi(1)=0.E+0
266     !!! LATITUDE. can be set.
267  latitude=0 ! default value for latitude
268  PRINT *,'latitude (in degrees) ?'
269  call getin("latitude",latitude)
270  write(*,*) " latitude = ",latitude
271  latitude=latitude*pi/180.E+0
272     !!! LONGITUDE. useless in 1D.
273  longitude=0.E+0
274  longitude=longitude*pi/180.E+0
275
276  rad=6400000.
277  !rad = -99999.
278  !PRINT *,'PLANETARY RADIUS in m ?'
279  !call getin("rad",rad)
280  ! Planetary  radius is needed to compute shadow of the rings
281  !IF (rad.eq.-99999. .and. rings_shadow .eqv. .true.) THEN
282  !   PRINT *,"STOP. I NEED rad IN RCM1D.DEF."
283  !   STOP
284  !ELSE
285  !   PRINT *,"--> rad = ",rad
286  !ENDIF
287
288
289
290  call init_physics_distribution(regular_lonlat,4,1,1,1,nlayer,1)
291  call init_interface_dyn_phys
292  CALL init_regular_lonlat(1,1,longitude,latitude,(/0.,0./),(/0.,0./))
293  call init_geometry(1,longitude,latitude,(/0.,0.,0.,0./),(/0.,0.,0.,0./),cell_area)
294  call init_dimphy(1,nlayer) ! Initialize dimphy module
295
296!!! CALL INIFIS
297!!! - read callphys.def
298!!! - calculate sine and cosine of longitude and latitude
299!!! - calculate mugaz and cp
300!!! NB: some operations are being done dummily in inifis in 1D
301!!! - physical constants: nevermind, things are done allright below
302!!! - physical frequency: nevermind, in inifis this is a simple print
303
304  cpp=1.d-7 !JL because we divide by cpp in inifis, there may be a more elegant solution
305  CALL inifis(1,llm,nq,0,86400.0,1,1.0,latitude,longitude,cell_area,rad,g,r,cpp)
306
307  !write(*,*) 'BASE 1'
308  !write(*,*) 1,llm,nq,0,86400.0,1,1.0,latitude,longitude,cell_area,rad,g,r,cpp
309
310  do iq=1,nq
311     do ilay=1,nlayer
312        q(ilay,iq) = 0.
313     enddo
314  enddo
315 
316  do iq=1,nq
317     qsurf(iq) = 0.
318  enddo
319
320  firstcall = .true.
321  lastcall  = .false.
322
323  iter    = 1
324  Tstrat  = 150.0
325  dTstrat = 1000.0
326
327  ! ---------
328  ! Run model
329  ! ---------
330  !do
331  psurf = psurf_n
332
333     !    Create vertical profiles
334  call kcmprof_fn(nlayer,psurf,qsurf(1),tsurf,    &
335       tstrat,play,plev,zlay,temp,q(:,1),muvar(1))
336         
337
338  !write(*,*) 'BASE 2'
339  !write(*,*) nlayer,psurf,qsurf(1),tsurf,    &
340  !           tstrat
341  !write(*,*) play
342  !write(*,*) plev
343  !write(*,*) zlay
344  !write(*,*) temp
345  !write(*,*) q(:,1),muvar(1)
346
347         
348               
349     call iniaerosol()
350     
351
352     
353
354     !    Run radiative transfer
355     call callcorrk(1,nlayer,q,nq,qsurf,                  &
356          albedo_wv,albedo_equivalent,                    &
357          emis,mu0,plev,play,temp,                        &
358          tsurf,fract,dist_star,aerosol,muvar,            &
359          dtlw,dtsw,fluxsurf_lw,fluxsurf_sw,              &
360          fluxsurfabs_sw,fluxtop_lw,                      &
361          fluxabs_sw,fluxtop_dn,OLR_nu,OSR_nu,tau_col,    &
362          cloudfrac,totcloudfrac,.false.,firstcall,lastcall)
363
364     !write(*,*) 'BASE 3'
365     !write(*,*) 1,nlayer,'A',q,'B',nq,'C',qsurf,                      &
366     !     albedo_wv,'D',albedo_equivalent,'E',                    &
367     !     emis,'F',mu0,'G',plev,'H',play,'I',temp,'J',                        &
368     !     tsurf,'K',fract,'L',dist_star,'M',aerosol,'N',muvar,'O',            &
369     !     dtlw,'P',dtsw,'Q',fluxsurf_lw,'R',fluxsurf_sw,'S',              &
370     !     fluxsurfabs_sw,'T',fluxtop_lw,'U',                      &
371     !     fluxabs_sw,'V',fluxtop_dn,'W',OLR_nu,'X',OSR_nu,'Y',tau_col,'Z',    &
372     !     cloudfrac,'A1',totcloudfrac,'A2',.false.,'A3',firstcall,'A4',lastcall
373
374
375
376
377
378
379
380     !    Iterate stratospheric temperature
381     print*,'Tstrat = ',Tstrat
382     dTstrat = Tstrat
383     !Tstrat  = Tsurf*(0.2786*(psurf/100000.)**(-1.123))**0.25
384     ! skin temperature (gray approx.) using analytic pure H2 expression
385     !Tstrat  = (fluxabs_sw/(2*sigma))**0.25 ! skin temperature (gray approx.)
386     Tstrat  = (fluxtop_lw/(2*sigma))**0.25 ! skin temperature (gray approx.)
387     dTstrat = dTstrat-Tstrat
388
389     !if(abs(dTstrat).lt.1.0)then
390     !   print*,'dTstrat = ',dTstrat
391     !   exit
392     !endif
393
394     !iter=iter+1
395     !if(iter.eq.100)then
396     !   print*,'Stratosphere failed to converge after'
397     !   print*,'100 iteration, aborting run.'
398     !   call abort
399     !endif
400
401  !end do
402
403
404  ! Run radiative transfer one last time to get OLR,OSR
405  firstcall=.false.
406  lastcall=.true.
407  call callcorrk(1,nlayer,q,nq,qsurf,                          &
408       albedo_wv,albedo_equivalent,emis,mu0,plev,play,temp,    &
409       tsurf,fract,dist_star,aerosol,muvar,                    &
410       dtlw,dtsw,fluxsurf_lw,fluxsurf_sw,fluxsurfabs_sw,       &
411       fluxtop_lw, fluxabs_sw,fluxtop_dn,OLR_nu,OSR_nu,        &
412       tau_col,                                                &
413       cloudfrac,totcloudfrac,.false.,firstcall,lastcall)
414
415
416     !write(*,*) 'BASE 4'
417     !write(*,*) 1,nlayer,'A',q,'B',nq,'C',qsurf,                      &
418     !     albedo_wv,'D',albedo_equivalent,'E',                    &
419     !     emis,'F',mu0,'G',plev,'H',play,'I',temp,'J',                        &
420     !     tsurf,'K',fract,'L',dist_star,'M',aerosol,'N',muvar,'O',            &
421     !     dtlw,'P',dtsw,'Q',fluxsurf_lw,'R',fluxsurf_sw,'S',              &
422     !     fluxsurfabs_sw,'T',fluxtop_lw,'U',                      &
423     !     fluxabs_sw,'V',fluxtop_dn,'W',OLR_nu,'X',OSR_nu,'Y',tau_col,'Z',    &
424     !     cloudfrac,'A1',totcloudfrac,'A2',.false.,'A3',firstcall,'A4',lastcall
425
426
427
428
429
430
431
432  ! Calculate total atmospheric energy
433  Eatmtot=0.0
434  !  call calcenergy_kcm(nlayer,tsurf,temp,play,plev,qsurf,&
435  !     q(:,1),muvar,Eatmtot)
436
437  ! ------------------------
438  ! Save data to ascii files
439  ! ------------------------
440
441  print*,'Saving profiles...'
442  open(115,file='profpres.out',form='formatted')
443  open(116,file='proftemp.out',form='formatted')
444  open(117,file='profztab.out',form='formatted')
445  open(118,file='profqvar.out',form='formatted')
446  open(119,file='profvvar.out',form='formatted')
447
448  write(115,*) psurf
449  write(116,*) tsurf
450  write(117,*) 0.0
451  write(118,*) qsurf(1)
452  write(119,*) qsurf(1)*(muvar(1)/mH2O)
453  do ilay=1,nlayer
454     vmr(ilay,1) = q(ilay,1)*(muvar(ilay+1)/mH2O)
455     write(115,*) play(ilay)
456     write(116,*) temp(ilay)
457     write(117,*) zlay(ilay)
458     write(118,*) q(ilay,1)
459     write(119,*) vmr(ilay,1)
460  enddo
461  close(115)
462  close(116)
463  close(117)
464  close(118)
465  close(119)
466
467  print*, tsurf,psurf,fluxtop_dn,fluxabs_sw,fluxtop_lw
468
469  print*,'Saving scalars...'
470  open(116,file='surf_vals.out')
471  write(116,*) tsurf,psurf,fluxtop_dn,         &
472       fluxabs_sw,fluxtop_lw
473  close(116)
474  open(111,file='ene_vals.out')
475  write(111,*) tsurf,psurf,Eatmtot,Tstrat
476  close(111)
477
478  print*,'Saving spectra...'
479  open(117,file='OLRnu.out')
480  do iw=1,L_NSPECTI
481     write(117,*) OLR_nu(1,iw)
482  enddo
483  close(117)
484
485  open(127,file='OSRnu.out')
486  do iw=1,L_NSPECTV
487     write(127,*) OSR_nu(1,iw)
488  enddo
489  close(127) 
490
491end program kcm1d
Note: See TracBrowser for help on using the repository browser.