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

Last change on this file since 2613 was 2470, checked in by yjaziri, 4 years ago

Generic GCM:
global1d and szangle for 1D simulation moved from callcorrk to callkeys
to defined a consistent 1D sza in physiq_mod used also in chemistry
+ typo
YJ

  • Property svn:executable set to *
File size: 15.4 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, global1d
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
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  !     added by JVO and YJ to read modern traceur.def
93  character(len=500) :: line ! to store a line of text
94  LOGICAL :: moderntracdef=.false. ! JVO, YJ : modern traceur.def
95
96  ! --------------
97  ! Initialisation
98  ! --------------
99
100  pi=2.E+0*asin(1.E+0)
101
102  cloudfrac(:)  = 0.0
103  totcloudfrac  = 0.0
104
105
106  nlayer=llm
107  nlevel=nlayer+1
108
109  !! this is defined in comsaison_h
110  ALLOCATE(mu0(1))
111  ALLOCATE(fract(1))
112  ALLOCATE(glat(1))
113
114  !  Load parameters from "run.def"
115  ! -------------------------------
116
117  ! check if 'kcm1d.def' file is around (otherwise reading parameters
118  ! from callphys.def via getin() routine won't work.)
119  open(90,file='kcm1d.def',status='old',form='formatted',&
120       iostat=ierr)
121  if (ierr.ne.0) then
122     write(*,*) 'Cannot find required file "kcm1d.def"'
123     write(*,*) '  (which should contain some input parameters'
124     write(*,*) '   along with the following line:'
125     write(*,*) '   INCLUDEDEF=callphys.def'
126     write(*,*) '   )'
127     write(*,*) ' ... might as well stop here ...'
128     stop
129  else
130     close(90)
131  endif
132 
133! now, run.def is needed anyway. so we create a dummy temporary one
134! for ioipsl to work. if a run.def is already here, stop the
135! program and ask the user to do a bit of cleaning
136  open(90,file='run.def',status='old',form='formatted',&
137       iostat=ierr)
138  if (ierr.eq.0) then
139     close(90)
140     write(*,*) 'There is already a run.def file.'
141     write(*,*) 'This is not compatible with 1D runs.'
142     write(*,*) 'Please remove the file and restart the run.'
143     write(*,*) 'Runtime parameters are supposed to be in kcm1d.def'
144     stop
145  else
146     call system('touch run.def')
147     call system("echo 'INCLUDEDEF=callphys.def' >> run.def")
148     call system("echo 'INCLUDEDEF=kcm1d.def' >> run.def")
149  endif
150   
151! check if we are going to run with or without tracers
152  write(*,*) "Run with or without tracer transport ?"
153  tracer=.false. ! default value
154  call getin("tracer",tracer)
155  write(*,*) " tracer = ",tracer
156
157
158
159  ! Tracer initialisation
160  ! ---------------------
161  if (tracer) then
162     ! load tracer names from file 'traceur.def'
163     open(90,file='traceur.def',status='old',form='formatted',&
164          iostat=ierr)
165     if (ierr.eq.0) then
166        write(*,*) "kcm1d: Reading file traceur.def"
167        ! read number of tracers:
168        !! - Modif. by JVO and YJ for modern planetary traceur.def ---------------
169        READ(90,'(A)') line
170        IF (trim(line).ne.'#ModernTrac-v1') THEN ! Test modern traceur.def
171           READ(line,*,iostat=ierr) nq ! Try standard traceur.def
172        ELSE
173           moderntracdef = .true.
174           DO
175             READ(90,'(A)',iostat=ierr) line
176             IF (ierr.eq.0) THEN
177               IF (index(line,'#').ne.1) THEN ! Allows arbitary number of comments lines in the header
178                 READ(line,*,iostat=ierr) nq
179                 EXIT
180               ENDIF
181             ELSE ! If pb, or if reached EOF without having found nbtr
182                write(*,*) "rcm1d: error reading number of tracers"
183                write(*,*) "   (first line of traceur.def) "
184                stop
185             ENDIF
186           ENDDO
187        ENDIF ! if modern or standard traceur.def
188        if (ierr.ne.0) then
189           write(*,*) "kcm1d: error reading number of tracers"
190           write(*,*) "   (first line of traceur.def) "
191           stop
192        endif
193        nqtot=nq
194        ! allocate arrays which depend on number of tracers
195        allocate(nametrac(nq))
196        allocate(q(nlayer,nq))
197        allocate(vmr(nlayer,nq))
198        allocate(qsurf(nq))
199
200        do iq=1,nq
201           ! minimal version, just read in the tracer names, 1 per line
202           read(90,*,iostat=ierr) nametrac(iq)
203           write(*,*) 'tracer here is', nametrac(iq)
204           if (ierr.ne.0) then
205              write(*,*) 'kcm1d: error reading tracer names...'
206              stop
207           endif
208        enddo !of do iq=1,nq
209        close(90)
210     endif
211
212     call initracer(1,nq,nametrac)
213
214  endif
215 
216
217
218
219
220  psurf_n=100000. ! default value for psurf
221  print*,'Dry surface pressure (Pa)?'
222  call getin("psurf",psurf_n)
223  write(*,*) " psurf = ",psurf_n
224
225! OK. now that run.def has been read once -- any variable is in memory.
226! so we can dump the dummy run.def
227  call system("rm -rf run.def")
228
229  tsurf=300.0 ! default value for tsurf
230  print*,'Surface temperature (K)?'
231  call getin("tref",tsurf)
232  write(*,*) " tsurf = ",tsurf
233
234  g=10.0 ! default value for g
235  print*,'Gravity ?'
236  call getin("g",g)
237  write(*,*) " g = ",g
238  glat(1)=g
239
240  periastr = 1.0
241  apoastr  = 1.0
242  print*,'Periastron (AU)?'
243  call getin("periastr",periastr)
244  write(*,*) "periastron = ",periastr
245  dist_star = periastr
246  fract     = 0.5
247  print*,'Apoastron (AU)?'
248  call getin("apoastr",apoastr)
249  write(*,*) "apoastron = ",apoastr
250
251  albedo=0.2 ! default value for albedo
252  print*,'Albedo of bare ground?'
253  call getin("albedo",albedo)
254  write(*,*) " albedo = ",albedo
255  do iw=1,L_NSPECTV
256     albedo_wv(1,iw)=albedo
257  enddo
258
259  emis=1.0 ! default value for emissivity
260  PRINT *,'Emissivity of bare ground ?'
261  call getin("emis",emis)
262  write(*,*) " emis = ",emis
263
264  pceil=100.0 ! Pascals
265  PRINT *,'Ceiling pressure (Pa) ?'
266  call getin("pceil",pceil)
267  write(*,*) " pceil = ", pceil
268
269  check_cpp_match = .false.
270  call getin("check_cpp_match",check_cpp_match)
271  if (check_cpp_match) then
272     print*,"In 1D modeling, check_cpp_match is supposed to be F"
273     print*,"Please correct callphys.def"
274     stop
275  endif
276
277
278!!! GEOGRAPHICAL INITIALIZATIONS
279     !!! AREA. useless in 1D
280  cell_area(1)=1.E+0 !JL+EM to have access to the area in the diagfi.nc files.
281     !!! GEOPOTENTIAL. useless in 1D because control by surface pressure
282  phisfi(1)=0.E+0
283     !!! LATITUDE. can be set.
284  latitude=0 ! default value for latitude
285  PRINT *,'latitude (in degrees) ?'
286  call getin("latitude",latitude)
287  write(*,*) " latitude = ",latitude
288  latitude=latitude*pi/180.E+0
289     !!! LONGITUDE. useless in 1D.
290  longitude=0.E+0
291  longitude=longitude*pi/180.E+0
292
293  rad=6400000.
294  !rad = -99999.
295  !PRINT *,'PLANETARY RADIUS in m ?'
296  !call getin("rad",rad)
297  ! Planetary  radius is needed to compute shadow of the rings
298  !IF (rad.eq.-99999. .and. rings_shadow .eqv. .true.) THEN
299  !   PRINT *,"STOP. I NEED rad IN RCM1D.DEF."
300  !   STOP
301  !ELSE
302  !   PRINT *,"--> rad = ",rad
303  !ENDIF
304
305
306
307  call init_physics_distribution(regular_lonlat,4,1,1,1,nlayer,1)
308  call init_interface_dyn_phys
309  CALL init_regular_lonlat(1,1,longitude,latitude,(/0.,0./),(/0.,0./))
310  call init_geometry(1,longitude,latitude,(/0.,0.,0.,0./),(/0.,0.,0.,0./),cell_area)
311  call init_dimphy(1,nlayer) ! Initialize dimphy module
312
313!!! CALL INIFIS
314!!! - read callphys.def
315!!! - calculate sine and cosine of longitude and latitude
316!!! - calculate mugaz and cp
317!!! NB: some operations are being done dummily in inifis in 1D
318!!! - physical constants: nevermind, things are done allright below
319!!! - physical frequency: nevermind, in inifis this is a simple print
320
321  cpp=1.d-7 !JL because we divide by cpp in inifis, there may be a more elegant solution
322  CALL inifis(1,llm,nq,0,86400.0,1,1.0,latitude,longitude,cell_area,rad,g,r,cpp)
323
324  if(.not.global1d)then
325     print*,'Error, kcm1d must have global1d=.true. in kcm1d.def!'
326     stop
327  end if
328
329  !write(*,*) 'BASE 1'
330  !write(*,*) 1,llm,nq,0,86400.0,1,1.0,latitude,longitude,cell_area,rad,g,r,cpp
331
332  do iq=1,nq
333     do ilay=1,nlayer
334        q(ilay,iq) = 0.
335     enddo
336  enddo
337 
338  do iq=1,nq
339     qsurf(iq) = 0.
340  enddo
341
342  firstcall = .true.
343  lastcall  = .false.
344
345  iter    = 1
346  Tstrat  = 150.0
347  dTstrat = 1000.0
348
349  ! ---------
350  ! Run model
351  ! ---------
352  !do
353  psurf = psurf_n
354
355     !    Create vertical profiles
356  call kcmprof_fn(nlayer,psurf,qsurf(1),tsurf,    &
357       tstrat,play,plev,zlay,temp,q(:,1),muvar(1))
358         
359
360  !write(*,*) 'BASE 2'
361  !write(*,*) nlayer,psurf,qsurf(1),tsurf,    &
362  !           tstrat
363  !write(*,*) play
364  !write(*,*) plev
365  !write(*,*) zlay
366  !write(*,*) temp
367  !write(*,*) q(:,1),muvar(1)
368
369         
370               
371     call iniaerosol()
372     
373
374     
375
376     !    Run radiative transfer
377     call callcorrk(1,nlayer,q,nq,qsurf,                  &
378          albedo_wv,albedo_equivalent,                    &
379          emis,mu0,plev,play,temp,                        &
380          tsurf,fract,dist_star,aerosol,muvar,            &
381          dtlw,dtsw,fluxsurf_lw,fluxsurf_sw,              &
382          fluxsurfabs_sw,fluxtop_lw,                      &
383          fluxabs_sw,fluxtop_dn,OLR_nu,OSR_nu,tau_col,    &
384          cloudfrac,totcloudfrac,.false.,firstcall,lastcall)
385
386     !write(*,*) 'BASE 3'
387     !write(*,*) 1,nlayer,'A',q,'B',nq,'C',qsurf,                      &
388     !     albedo_wv,'D',albedo_equivalent,'E',                    &
389     !     emis,'F',mu0,'G',plev,'H',play,'I',temp,'J',                        &
390     !     tsurf,'K',fract,'L',dist_star,'M',aerosol,'N',muvar,'O',            &
391     !     dtlw,'P',dtsw,'Q',fluxsurf_lw,'R',fluxsurf_sw,'S',              &
392     !     fluxsurfabs_sw,'T',fluxtop_lw,'U',                      &
393     !     fluxabs_sw,'V',fluxtop_dn,'W',OLR_nu,'X',OSR_nu,'Y',tau_col,'Z',    &
394     !     cloudfrac,'A1',totcloudfrac,'A2',.false.,'A3',firstcall,'A4',lastcall
395
396
397
398
399
400
401
402     !    Iterate stratospheric temperature
403     print*,'Tstrat = ',Tstrat
404     dTstrat = Tstrat
405     !Tstrat  = Tsurf*(0.2786*(psurf/100000.)**(-1.123))**0.25
406     ! skin temperature (gray approx.) using analytic pure H2 expression
407     !Tstrat  = (fluxabs_sw/(2*sigma))**0.25 ! skin temperature (gray approx.)
408     Tstrat  = (fluxtop_lw/(2*sigma))**0.25 ! skin temperature (gray approx.)
409     dTstrat = dTstrat-Tstrat
410
411     !if(abs(dTstrat).lt.1.0)then
412     !   print*,'dTstrat = ',dTstrat
413     !   exit
414     !endif
415
416     !iter=iter+1
417     !if(iter.eq.100)then
418     !   print*,'Stratosphere failed to converge after'
419     !   print*,'100 iteration, aborting run.'
420     !   call abort
421     !endif
422
423  !end do
424
425
426  ! Run radiative transfer one last time to get OLR,OSR
427  firstcall=.false.
428  lastcall=.true.
429  call callcorrk(1,nlayer,q,nq,qsurf,                          &
430       albedo_wv,albedo_equivalent,emis,mu0,plev,play,temp,    &
431       tsurf,fract,dist_star,aerosol,muvar,                    &
432       dtlw,dtsw,fluxsurf_lw,fluxsurf_sw,fluxsurfabs_sw,       &
433       fluxtop_lw, fluxabs_sw,fluxtop_dn,OLR_nu,OSR_nu,        &
434       tau_col,                                                &
435       cloudfrac,totcloudfrac,.false.,firstcall,lastcall)
436
437
438     !write(*,*) 'BASE 4'
439     !write(*,*) 1,nlayer,'A',q,'B',nq,'C',qsurf,                      &
440     !     albedo_wv,'D',albedo_equivalent,'E',                    &
441     !     emis,'F',mu0,'G',plev,'H',play,'I',temp,'J',                        &
442     !     tsurf,'K',fract,'L',dist_star,'M',aerosol,'N',muvar,'O',            &
443     !     dtlw,'P',dtsw,'Q',fluxsurf_lw,'R',fluxsurf_sw,'S',              &
444     !     fluxsurfabs_sw,'T',fluxtop_lw,'U',                      &
445     !     fluxabs_sw,'V',fluxtop_dn,'W',OLR_nu,'X',OSR_nu,'Y',tau_col,'Z',    &
446     !     cloudfrac,'A1',totcloudfrac,'A2',.false.,'A3',firstcall,'A4',lastcall
447
448
449
450
451
452
453
454  ! Calculate total atmospheric energy
455  Eatmtot=0.0
456  !  call calcenergy_kcm(nlayer,tsurf,temp,play,plev,qsurf,&
457  !     q(:,1),muvar,Eatmtot)
458
459  ! ------------------------
460  ! Save data to ascii files
461  ! ------------------------
462
463  print*,'Saving profiles...'
464  open(115,file='profpres.out',form='formatted')
465  open(116,file='proftemp.out',form='formatted')
466  open(117,file='profztab.out',form='formatted')
467  open(118,file='profqvar.out',form='formatted')
468  open(119,file='profvvar.out',form='formatted')
469
470  write(115,*) psurf
471  write(116,*) tsurf
472  write(117,*) 0.0
473  write(118,*) qsurf(1)
474  write(119,*) qsurf(1)*(muvar(1)/mH2O)
475  do ilay=1,nlayer
476     vmr(ilay,1) = q(ilay,1)*(muvar(ilay+1)/mH2O)
477     write(115,*) play(ilay)
478     write(116,*) temp(ilay)
479     write(117,*) zlay(ilay)
480     write(118,*) q(ilay,1)
481     write(119,*) vmr(ilay,1)
482  enddo
483  close(115)
484  close(116)
485  close(117)
486  close(118)
487  close(119)
488
489  print*, tsurf,psurf,fluxtop_dn,fluxabs_sw,fluxtop_lw
490
491  print*,'Saving scalars...'
492  open(116,file='surf_vals.out')
493  write(116,*) tsurf,psurf,fluxtop_dn,         &
494       fluxabs_sw,fluxtop_lw
495  close(116)
496  open(111,file='ene_vals.out')
497  write(111,*) tsurf,psurf,Eatmtot,Tstrat
498  close(111)
499
500  print*,'Saving spectra...'
501  open(117,file='OLRnu.out')
502  do iw=1,L_NSPECTI
503     write(117,*) OLR_nu(1,iw)
504  enddo
505  close(117)
506
507  open(127,file='OSRnu.out')
508  do iw=1,L_NSPECTV
509     write(127,*) OSR_nu(1,iw)
510  enddo
511  close(127) 
512
513end program kcm1d
Note: See TracBrowser for help on using the repository browser.