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

Last change on this file since 3995 was 3893, checked in by gmilcareck, 8 months ago

Remove all "call abort" and "stop" statement in LMDZ.GENERIC and replacing them by call abort_physic().

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