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

Last change on this file was 3929, checked in by debatzbr, 2 months ago

Pluto PCM: Add optical diagnostics for aerosols
BBT

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