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

Last change on this file since 1149 was 787, checked in by aslmd, 12 years ago

LMDZ.GENERIC. (Sorry for long text but this is a quite major commit)

Paving the path for parallel computations. And moving towards a more flexible code.

Automatic allocation is used within all routines in phystd. No further mention to ngridmx and nqmx.

  1. ngridmx and nqmx are still used in LMDZ.GENERIC in the dyn3d part
  2. if the LMDZ4/LMDZ5 dynamical core is used, there is no more fixed dimensions ngridmx and nqmx --> a fully flexible parallel implementation is now possible (e.g. no need to recompile when changing numbers of processors)

The important stuff :

  • Compilation checked with ifort. OK with and without debug mode. No errors. Checked for: gcm, newstart, rcm1d, kcm1d
  • RUN GCM: Running an Earth test case. Comparison with previous revision --> debug mode : perfect match. bit by bit (diff command). checked with plots --> O1 mode : close match (checked with plots) --> O2 mode : sometimes up to 0.5 K departure.... BUT in this new version O2 and O1 are quite close while in previous version O1 and O2 differed by about, well, typically 0.5 K (pictures available on request)
  • RUN NEWSTART : perfect match (bit-by-bit) in either debug or normal mode.
  • RUN RCM1D : perfect match in normal mode.
  • RUN KCM1D : not tested (I don't know what is the use of kcm1d)

List of main changes :

  • Additional arguments to some subroutines (ngrid and nq)
  • F77 include strategy is obsolete and replaced by F90 module strategy In this new strategy arrays are allocatable and allocated once at first use This has to be done for all common featuring arrays defined with ngridmx or nqmx

surfdat.h >> surfdat_h.F90
tracer.h >> tracer_h.F90
comsaison.h >> comsaison_h.F90
comgeomfi.h >> comgeomfi_h.F90
comsoil.h >> comsoil_h.F90
comdiurn.h >> comdiurn_h.F90
fisice.h >> DELETED. was not used. probably a fossil.
watercap.h >> DELETED. variable put in surfdat_h.F90

  • F77 'save' strategy is obsolete and replaced by F90 'allocatable save' strategy (see previous point and e.g. new version of physiq.F90)
  • Suppressing any mention to advtrac.h which is a common in the dynamics and needs nqmx This was easily solved by adding an argument with tracer names, coming from the dynamics This is probably not a definitive solution, ... but this allows for generic physics to work easily with either LMDZ.GENERIC or LMDZ dynamical cores
  • Removing consistency tests between nq and nqmx ; and ngrid and ngridmx. No use now!
  • Adaptation of rcm1d, kcm1d, newstart given above-mentioned changes

A note on phyetat0 and soil_setting:

  • Now written so that a slice of horizontal size 'ngrid' starting at grid point 'cursor' is read in startfi.nc 'cursor' is defined in dimphys.h and initialized by inifis (or in newstart) this is useful for parallel computations. default behavior is the usual one : sequential runs, cursor is 1, size ngrid is the whole global domain

A note on an additional change :

  • nueffrad is now an argument to callcorrk as is the case for reffrad both are saved in physiq this is for consistency and lisibility (previously nueffrad was saved in callcorrk) ... but there is a call to a function which modifies nueffrad in physiq ... previously this was not modifying nueffrad (although it was quite cumbersome to detect this) ... to be conservative I kept this behaviour and highlighted it with an array nueffrad_dummy ... I added a comment because someone might want to change this
  • Property svn:executable set to *
File size: 10.3 KB
Line 
1program kcm1d
2
3  use radinc_h,      only: NAERKIND
4  use radcommon_h,   only: L_NSPECTI, L_NSPECTV, sigma
5  use watercommon_h, only: mH2O
6  use ioipsl_getincom
7  use comsaison_h
8
9  implicit none
10
11  !==================================================================
12  !     
13  !     Purpose
14  !     -------
15  !     Run the universal model radiative transfer once in a 1D column.
16  !     Useful for climate evolution studies etc.
17  !     
18  !     It can be compiled with a command like (e.g. for 25 layers):
19  !                                  "makegcm -p std -d 25 kcm1d"
20  !     It requires the files "callphys.def", "gases.def"
21  !     "traceur.def", and "run.def" with a line "INCLUDEDEF=callphys.def"
22  !
23  !     Authors
24  !     -------
25  !     R. Wordsworth
26  !
27  !==================================================================
28
29#include "dimensions.h"
30#include "dimphys.h"
31#include "callkeys.h"
32#include "comcstfi.h"
33#include "planete.h"
34#include "control.h"
35
36  ! --------------------------------------------------------------
37  !  Declarations
38  ! --------------------------------------------------------------
39
40  integer nlayer,nlevel,nq
41  integer ilay,ilev,iq,iw,iter
42  real play(nlayermx)     ! Pressure at the middle of the layers [Pa]
43  real zlay(nlayermx)     ! Altitude at middle of the layers [km]
44  real plev(nlayermx+1)   ! Intermediate pressure levels [Pa]
45  real temp(nlayermx)     ! temperature at the middle of the layers [K]
46  real q(nlayermx,nqmx)   ! tracer mixing ratio [kg/kg]
47  real vmr(nlayermx,nqmx) ! tracer mixing ratio [mol/mol]
48  real qsurf(nqmx)        ! tracer surface budget [kg/kg] ????
49  real psurf,psurf_n,tsurf     
50
51  real emis, albedo
52
53  real muvar(nlayermx+1)
54
55  real dtsw(nlayermx) ! heating rate (K/s) due to SW
56  real dtlw(nlayermx) ! heating rate (K/s) due to LW
57  real fluxsurf_lw   ! incident LW flux to surf (W/m2)
58  real fluxtop_lw    ! outgoing LW flux to space (W/m2)
59  real fluxsurf_sw   ! incident SW flux to surf (W/m2)
60  real fluxabs_sw    ! SW flux absorbed by planet (W/m2)
61  real fluxtop_dn    ! incident top of atmosphere SW flux (W/m2)
62
63  ! not used
64  real reffrad(nlayermx,naerkind)
65  real nueffrad(nlayermx,naerkind)
66  real cloudfrac(nlayermx)
67  real totcloudfrac
68  real tau_col
69
70  real dTstrat
71  real aerosol(nlayermx,naerkind) ! aerosol tau (kg/kg)
72  real OLR_nu(1,L_NSPECTI)
73  real OSR_nu(1,L_NSPECTV)
74  real Eatmtot
75
76  integer ierr
77  logical firstcall,lastcall,global1d
78
79  character*20 nametrac(nqmx)   ! name of the tracer (no need for adv trac common)
80
81  ! --------------
82  ! Initialisation
83  ! --------------
84
85  pi=2.E+0*asin(1.E+0)
86
87  reffrad(:,:)  = 0.0
88  nueffrad(:,:) = 0.0
89  cloudfrac(:)  = 0.0
90  totcloudfrac  = 0.0
91
92
93  nlayer=nlayermx
94  nlevel=nlayer+1
95
96  !! this is defined in comsaison_h
97  ALLOCATE(mu0(1))
98  ALLOCATE(fract(1))
99
100
101
102  !  Load parameters from "run.def"
103  ! -------------------------------
104
105  ! check if 'kcm1d.def' file is around (otherwise reading parameters
106  ! from callphys.def via getin() routine won't work.)
107  open(90,file='kcm1d.def',status='old',form='formatted',&
108       iostat=ierr)
109  if (ierr.ne.0) then
110     write(*,*) 'Cannot find required file "kcm1d.def"'
111     write(*,*) '  (which should contain some input parameters'
112     write(*,*) '   along with the following line:'
113     write(*,*) '   INCLUDEDEF=callphys.def'
114     write(*,*) '   )'
115     write(*,*) ' ... might as well stop here ...'
116     stop
117  else
118     close(90)
119  endif
120
121! now, run.def is needed anyway. so we create a dummy temporary one
122! for ioipsl to work. if a run.def is already here, stop the
123! program and ask the user to do a bit of cleaning
124  open(90,file='run.def',status='old',form='formatted',&
125       iostat=ierr)
126  if (ierr.eq.0) then
127     close(90)
128     write(*,*) 'There is already a run.def file.'
129     write(*,*) 'This is not compatible with 1D runs.'
130     write(*,*) 'Please remove the file and restart the run.'
131     write(*,*) 'Runtime parameters are supposed to be in kcm1d.def'
132     stop
133  else
134     call system('touch run.def')
135     call system("echo 'INCLUDEDEF=callphys.def' >> run.def")
136     call system("echo 'INCLUDEDEF=kcm1d.def' >> run.def")
137  endif
138
139  global1d = .false. ! default value
140  call getin("global1d",global1d)
141  if(.not.global1d)then
142     print*,'Error, kcm1d must have global1d=.true. in kcm1d.def!'
143     stop
144  end if
145
146  psurf_n=100000. ! default value for psurf
147  print*,'Dry surface pressure (Pa)?'
148  call getin("psurf",psurf_n)
149  write(*,*) " psurf = ",psurf_n
150
151! OK. now that run.def has been read once -- any variable is in memory.
152! so we can dump the dummy run.def
153  call system("rm -rf run.def")
154
155  tsurf=300.0 ! default value for tsurf
156  print*,'Surface temperature (K)?'
157  call getin("tref",tsurf)
158  write(*,*) " tsurf = ",tsurf
159
160  g=10.0 ! default value for g
161  print*,'Gravity ?'
162  call getin("g",g)
163  write(*,*) " g = ",g
164
165  periastr = 1.0
166  apoastr  = 1.0
167  print*,'Periastron (AU)?'
168  call getin("periastr",periastr)
169  write(*,*) "periastron = ",periastr
170  dist_star = periastr
171  fract     = 0.5
172  print*,'Apoastron (AU)?'
173  call getin("apoastr",apoastr)
174  write(*,*) "apoastron = ",apoastr
175
176  albedo=0.2 ! default value for albedo
177  print*,'Albedo of bare ground?'
178  call getin("albedo",albedo)
179  write(*,*) " albedo = ",albedo
180
181  emis=1.0 ! default value for emissivity
182  PRINT *,'Emissivity of bare ground ?'
183  call getin("emis",emis)
184  write(*,*) " emis = ",emis
185
186  pceil=100.0 ! Pascals
187  PRINT *,'Ceiling pressure (Pa) ?'
188  call getin("pceil",pceil)
189  write(*,*) " pceil = ", pceil
190
191  mugaz=0.
192  cpp=0.
193
194  check_cpp_match = .false.
195  call getin("check_cpp_match",check_cpp_match)
196  if (check_cpp_match) then
197     print*,"In 1D modeling, check_cpp_match is supposed to be F"
198     print*,"Please correct callphys.def"
199     stop
200  endif
201
202  call su_gases
203  call calc_cpp_mugaz
204
205  call inifis(1,llm,0,86400.0,1.0,0.0,0.0,1.0,rad,g,r,cpp)
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        read(90,*,iostat=ierr) nq
217        if (ierr.ne.0) then
218           write(*,*) "kcm1d: error reading number of tracers"
219           write(*,*) "   (first line of traceur.def) "
220           stop
221        else
222           ! check that the number of tracers is indeed nqmx
223           if (nq.ne.nqmx) then
224              write(*,*) "kcm1d: error, wrong number of tracers:"
225              write(*,*) "nq=",nq," whereas nqmx=",nqmx
226              stop
227           endif
228        endif
229
230        do iq=1,nq
231           ! minimal version, just read in the tracer names, 1 per line
232           read(90,*,iostat=ierr) nametrac(iq)
233           if (ierr.ne.0) then
234              write(*,*) 'kcm1d: error reading tracer names...'
235              stop
236           endif
237        enddo !of do iq=1,nq
238     endif
239
240     call initracer(1,nq,nametrac)
241
242  endif
243
244
245  do iq=1,nq
246     do ilay=1,nlayer
247        q(ilay,iq) = 0.
248     enddo
249  enddo
250
251  do iq=1,nq
252     qsurf(iq) = 0.
253  enddo
254
255  firstcall = .true.
256  lastcall  = .false.
257
258  iter    = 1
259  Tstrat  = 150.0
260  dTstrat = 1000.0
261
262  ! ---------
263  ! Run model
264  ! ---------
265  !do
266     psurf = psurf_n
267
268     !    Create vertical profiles
269     call kcmprof_fn(psurf,qsurf(1),tsurf,    &
270          tstrat,play,plev,zlay,temp,q(:,1),muvar(1))
271
272     !    Run radiative transfer
273     call callcorrk(1,nlayer,q,nq,qsurf,      &
274          albedo,emis,mu0,plev,play,temp,                    &
275          tsurf,fract,dist_star,aerosol,muvar,         &
276          dtlw,dtsw,fluxsurf_lw,fluxsurf_sw,fluxtop_lw,    &
277          fluxabs_sw,fluxtop_dn,OLR_nu,OSR_nu,reffrad,nueffrad,tau_col,  &
278          cloudfrac,totcloudfrac,.false.,firstcall,lastcall)
279
280     !    Iterate stratospheric temperature
281     print*,'Tstrat = ',Tstrat
282     dTstrat = Tstrat
283     !Tstrat  = Tsurf*(0.2786*(psurf/100000.)**(-1.123))**0.25
284     ! skin temperature (gray approx.) using analytic pure H2 expression
285     !Tstrat  = (fluxabs_sw/(2*sigma))**0.25 ! skin temperature (gray approx.)
286     Tstrat  = (fluxtop_lw/(2*sigma))**0.25 ! skin temperature (gray approx.)
287     dTstrat = dTstrat-Tstrat
288
289     !if(abs(dTstrat).lt.1.0)then
290     !   print*,'dTstrat = ',dTstrat
291     !   exit
292     !endif
293
294     !iter=iter+1
295     !if(iter.eq.100)then
296     !   print*,'Stratosphere failed to converge after'
297     !   print*,'100 iteration, aborting run.'
298     !   call abort
299     !endif
300
301  !end do
302
303  ! Run radiative transfer one last time to get OLR,OSR
304  firstcall=.false.
305  lastcall=.true.
306  call callcorrk(1,nlayer,q,nq,qsurf,      &
307       albedo,emis,mu0,plev,play,temp,                    &
308       tsurf,fract,dist_star,aerosol,muvar,         &
309       dtlw,dtsw,fluxsurf_lw,fluxsurf_sw,fluxtop_lw,    &
310       fluxabs_sw,fluxtop_dn,OLR_nu,OSR_nu,             &
311       reffrad,nueffrad,tau_col,  &
312       cloudfrac,totcloudfrac,.false.,firstcall,lastcall)
313
314
315  ! Calculate total atmospheric energy
316  Eatmtot=0.0
317  !  call calcenergy_kcm(tsurf,temp,play,plev,qsurf,&
318  !     q(:,1),muvar,Eatmtot)
319
320  ! ------------------------
321  ! Save data to ascii files
322  ! ------------------------
323
324  print*,'Saving profiles...'
325  open(115,file='profpres.out',form='formatted')
326  open(116,file='proftemp.out',form='formatted')
327  open(117,file='profztab.out',form='formatted')
328  open(118,file='profqvar.out',form='formatted')
329  open(119,file='profvvar.out',form='formatted')
330
331  write(115,*) psurf
332  write(116,*) tsurf
333  write(117,*) 0.0
334  write(118,*) qsurf(1)
335  write(119,*) qsurf(1)*(muvar(1)/mH2O)
336  do ilay=1,nlayer
337     vmr(ilay,1) = q(ilay,1)*(muvar(ilay+1)/mH2O)
338     write(115,*) play(ilay)
339     write(116,*) temp(ilay)
340     write(117,*) zlay(ilay)
341     write(118,*) q(ilay,1)
342     write(119,*) vmr(ilay,1)
343  enddo
344  close(115)
345  close(116)
346  close(117)
347  close(118)
348  close(119)
349
350  print*, tsurf,psurf,fluxtop_dn,fluxabs_sw,fluxtop_lw
351
352  print*,'Saving scalars...'
353  open(116,file='surf_vals.out')
354  write(116,*) tsurf,psurf,fluxtop_dn,         &
355       fluxabs_sw,fluxtop_lw
356  close(116)
357  open(111,file='ene_vals.out')
358  write(111,*) tsurf,psurf,Eatmtot,Tstrat
359  close(111)
360
361  print*,'Saving spectra...'
362  open(117,file='OLRnu.out')
363  do iw=1,L_NSPECTI
364     write(117,*) OLR_nu(1,iw)
365  enddo
366  close(117)
367
368  open(127,file='OSRnu.out')
369  do iw=1,L_NSPECTV
370     write(127,*) OSR_nu(1,iw)
371  enddo
372  close(127) 
373
374end program kcm1d
Note: See TracBrowser for help on using the repository browser.