source: trunk/LMDZ.PLUTO/libf/phypluto/dyn1d/kcm1d.F90 @ 3586

Last change on this file since 3586 was 3585, checked in by debatzbr, 13 days ago

Connecting microphysics to radiative transfer + miscellaneous cleans

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