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

Last change on this file since 3574 was 3574, checked in by jbclement, 32 hours ago

COMMON:
The compilation generates a Fortran subroutine to track compilation and version (SVN or Git) details through the executable file. Put the argument 'version' as an option when executing the code to display these information and create a file "version_diff.txt" containing the diff result.
It can work with every program but it has been implemented only for: 'gcm', 'parallel gcm', 'pem', '1D Mars PCM', 'Mars newstart', 'Mars start2archive' and '1D Generic PCM'.
JBC

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