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

Last change on this file since 4160 was 4144, checked in by debatzbr, 2 weeks ago

Pluto PCM: Add IR opacity diagnostics for haze only.
BBT

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