source: trunk/LMDZ.GENERIC/libf/phystd/newstart.F @ 995

Last change on this file since 995 was 993, checked in by emillour, 12 years ago

Generic GCM:

  • Some more cleanup in dynamics:
    • Moved "start2archive" (and auxilliary routines) to phystd
    • removed unused (obsolete) testharm.F , para_netcdf.h , readhead_NC.F , angtot.h from dyn3d
    • removed obsolete addit.F (and change corresponding lines in gcm)
    • remove unused "description.h" (and many places where it was "included")

EM

File size: 50.9 KB
Line 
1C======================================================================
2      PROGRAM newstart
3c=======================================================================
4c
5c
6c   Auteur:   Christophe Hourdin/Francois Forget/Yann Wanherdrick
7c   ------
8c             Derniere modif : 12/03
9c
10c
11c   Objet:  Create or modify the initial state for the LMD Mars GCM
12c   -----           (fichiers NetCDF start et startfi)
13c
14c
15c=======================================================================
16
17      USE tracer_h
18      USE comsoil_h
19      USE surfdat_h
20      USE comgeomfi_h
21      use datafile_mod, only: datadir
22! to use  'getin'
23      USE ioipsl_getincom, only: getin
24      implicit none
25
26#include "dimensions.h"
27#include "dimphys.h"
28#include "planete.h"
29#include "paramet.h"
30#include "comconst.h"
31#include "comvert.h"
32#include "comgeom2.h"
33#include "control.h"
34#include "logic.h"
35#include "ener.h"
36#include "temps.h"
37#include "comdissnew.h"
38#include "serre.h"
39#include "netcdf.inc"
40#include "advtrac.h"
41c=======================================================================
42c   Declarations
43c=======================================================================
44
45c Variables dimension du fichier "start_archive"
46c------------------------------------
47      CHARACTER relief*3
48
49
50c Variables pour les lectures NetCDF des fichiers "start_archive"
51c--------------------------------------------------
52      INTEGER nid_dyn, nid_fi,nid,nvarid
53      INTEGER length
54      parameter (length = 100)
55      INTEGER tab0
56      INTEGER   NB_ETATMAX
57      parameter (NB_ETATMAX = 100)
58
59      REAL  date
60      REAL p_rad,p_omeg,p_g,p_cpp,p_mugaz,p_daysec
61
62c Variable histoire
63c------------------
64      REAL vcov(iip1,jjm,llm),ucov(iip1,jjp1,llm) ! vents covariants
65      REAL phis(iip1,jjp1)
66      REAL q(iip1,jjp1,llm,nqmx)               ! champs advectes
67
68c autre variables dynamique nouvelle grille
69c------------------------------------------
70      REAL pks(iip1,jjp1)
71      REAL w(iip1,jjp1,llm+1)
72      REAL pbaru(ip1jmp1,llm),pbarv(ip1jm,llm)
73!      REAL dv(ip1jm,llm),du(ip1jmp1,llm)
74!      REAL dh(ip1jmp1,llm),dp(ip1jmp1)
75      REAL phi(iip1,jjp1,llm)
76
77      integer klatdat,klongdat
78      PARAMETER (klatdat=180,klongdat=360)
79
80c Physique sur grille scalaire
81c----------------------------
82      real zmeaS(iip1,jjp1),zstdS(iip1,jjp1)
83      real zsigS(iip1,jjp1),zgamS(iip1,jjp1),ztheS(iip1,jjp1)
84
85c variable physique
86c------------------
87      REAL tsurf(ngridmx)       ! surface temperature
88      REAL tsoil(ngridmx,nsoilmx) ! soil temperature
89!      REAL co2ice(ngridmx)     ! CO2 ice layer
90      REAL emis(ngridmx)        ! surface emissivity
91      real emisread             ! added by RW
92      REAL qsurf(ngridmx,nqmx)
93      REAL q2(ngridmx,nlayermx+1)
94!      REAL rnaturfi(ngridmx)
95      real alb(iip1,jjp1),albfi(ngridmx) ! albedos
96      real ith(iip1,jjp1,nsoilmx),ithfi(ngridmx,nsoilmx) ! thermal inertia (3D)
97      real surfith(iip1,jjp1),surfithfi(ngridmx) ! surface thermal inertia (2D)
98      REAL latfi(ngridmx),lonfi(ngridmx),airefi(ngridmx)
99
100      INTEGER i,j,l,isoil,ig,idum
101      real mugaz ! molar mass of the atmosphere
102
103      integer ierr
104
105c Variables on the new grid along scalar points
106c------------------------------------------------------
107!      REAL p(iip1,jjp1)
108      REAL t(iip1,jjp1,llm)
109      REAL tset(iip1,jjp1,llm)
110      real phisold_newgrid(iip1,jjp1)
111      REAL :: teta(iip1, jjp1, llm)
112      REAL :: pk(iip1,jjp1,llm)
113      REAL :: pkf(iip1,jjp1,llm)
114      REAL :: ps(iip1, jjp1)
115      REAL :: masse(iip1,jjp1,llm)
116      REAL :: xpn,xps,xppn(iim),xpps(iim)
117      REAL :: p3d(iip1, jjp1, llm+1)
118      REAL :: beta(iip1,jjp1,llm)
119!      REAL dteta(ip1jmp1,llm)
120
121c Variable de l'ancienne grille
122c------------------------------
123      real time
124      real tab_cntrl(100)
125      real tab_cntrl_bis(100)
126
127c variables diverses
128c-------------------
129      real choix_1,pp
130      character*80      fichnom
131      character*250  filestring
132      integer Lmodif,iq,thermo
133      character modif*20
134      real z_reel(iip1,jjp1)
135      real tsud,albsud,alb_bb,ith_bb,Tiso,Tabove
136      real ptoto,pcap,patm,airetot,ptotn,patmn,psea
137!      real ssum
138      character*1 yes
139      logical :: flagtset=.false. ,  flagps0=.false.
140      real val, val2, val3 ! to store temporary variables
141      real :: iceith=2000 ! thermal inertia of subterranean ice
142      integer iref,jref
143
144      INTEGER :: itau
145     
146      INTEGER :: nq,numvanle
147      character(len=20) :: txt ! to store some text
148      character(len=50) :: surfacefile ! "surface.nc" (or equivalent file)
149      character(len=150) :: longline
150      integer :: count
151      real :: profile(llm+1) ! to store an atmospheric profile + surface value
152
153!     added by RW for test
154      real pmean, phi0
155
156!     added by BC for equilibrium temperature startup
157      real teque
158
159!     added by BC for cloud fraction setup
160      REAL hice(ngridmx),cloudfrac(ngridmx,nlayermx)
161      REAL totalfrac(ngridmx)
162
163!     added by RW for nuketharsis
164      real fact1
165      real fact2
166
167
168c sortie visu pour les champs dynamiques
169c---------------------------------------
170!      INTEGER :: visuid
171!      real :: time_step,t_ops,t_wrt
172!      CHARACTER*80 :: visu_file
173
174      cpp    = 0.
175      preff  = 0.
176      pa     = 0. ! to ensure disaster rather than minor error if we don`t
177                  ! make deliberate choice of these values elsewhere.
178
179c=======================================================================
180c   Choice of the start file(s) to use
181c=======================================================================
182
183      write(*,*) 'From which kind of files do you want to create new',
184     .  'start and startfi files'
185      write(*,*) '    0 - from a file start_archive'
186      write(*,*) '    1 - from files start and startfi'
187 
188c-----------------------------------------------------------------------
189c   Open file(s) to modify (start or start_archive)
190c-----------------------------------------------------------------------
191
192      DO
193         read(*,*,iostat=ierr) choix_1
194         if ((choix_1 /= 0).OR.(choix_1 /=1)) EXIT
195      ENDDO
196
197c     Open start_archive
198c     ~~~~~~~~~~~~~~~~~~~~~~~~~~
199      if (choix_1.eq.0) then
200
201        write(*,*) 'Creating start files from:'
202        write(*,*) './start_archive.nc'
203        write(*,*)
204        fichnom = 'start_archive.nc'
205        ierr = NF_OPEN (fichnom, NF_NOWRITE,nid)
206        IF (ierr.NE.NF_NOERR) THEN
207          write(6,*)' Problem opening file:',fichnom
208          write(6,*)' ierr = ', ierr
209          CALL ABORT
210        ENDIF
211        tab0 = 50
212        Lmodif = 1
213
214c     OR open start and startfi files
215c     ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
216      else
217        write(*,*) 'Creating start files from:'
218        write(*,*) './start.nc and ./startfi.nc'
219        write(*,*)
220        fichnom = 'start.nc'
221        ierr = NF_OPEN (fichnom, NF_NOWRITE,nid_dyn)
222        IF (ierr.NE.NF_NOERR) THEN
223          write(6,*)' Problem opening file:',fichnom
224          write(6,*)' ierr = ', ierr
225          CALL ABORT
226        ENDIF
227 
228        fichnom = 'startfi.nc'
229        ierr = NF_OPEN (fichnom, NF_NOWRITE,nid_fi)
230        IF (ierr.NE.NF_NOERR) THEN
231          write(6,*)' Problem opening file:',fichnom
232          write(6,*)' ierr = ', ierr
233          CALL ABORT
234        ENDIF
235
236        tab0 = 0
237        Lmodif = 0
238
239      endif
240
241
242c=======================================================================
243c  INITIALISATIONS DIVERSES
244c=======================================================================
245! Load tracer names:
246      call iniadvtrac(nq,numvanle)
247      ! tnom(:) now contains tracer names
248! JL12 we will need the tracer names to read start in dyneta0
249
250! Initialize global tracer indexes (stored in tracer.h)
251! ... this has to be done before phyetat0
252      call initracer(ngridmx,nqmx,tnom)
253
254! Take care of arrays in common modules
255      ! ALLOCATE ARRAYS in surfdat_h (if not already done, e.g. when using start_archive)
256      IF (.not. ALLOCATED(albedodat)) ALLOCATE(albedodat(ngridmx))
257      IF (.not. ALLOCATED(phisfi)) ALLOCATE(phisfi(ngridmx))
258      IF (.not. ALLOCATED(zmea)) ALLOCATE(zmea(ngridmx))
259      IF (.not. ALLOCATED(zstd)) ALLOCATE(zstd(ngridmx))
260      IF (.not. ALLOCATED(zsig)) ALLOCATE(zsig(ngridmx))
261      IF (.not. ALLOCATED(zgam)) ALLOCATE(zgam(ngridmx))
262      IF (.not. ALLOCATED(zthe)) ALLOCATE(zthe(ngridmx))
263      ! ALLOCATE ARRAYS in comsoil_h (if not already done)
264      IF (.not.ALLOCATED(layer))
265     . ALLOCATE(layer(nsoilmx))
266      IF (.not.ALLOCATED(mlayer))
267     . ALLOCATE(mlayer(0:nsoilmx-1))
268      IF (.not.ALLOCATED(inertiedat))
269     . ALLOCATE(inertiedat(ngridmx,nsoilmx))
270      ! ALLOCATE ARRAYS IN comgeomfi_h (done in inifis usually)
271      IF (.not. ALLOCATED(lati)) ALLOCATE(lati(ngridmx))
272      IF (.not. ALLOCATED(long)) ALLOCATE(long(ngridmx))
273      IF (.not. ALLOCATED(area)) ALLOCATE(area(ngridmx))
274      cursor = 1 ! added by AS in dimphys. 1 for sequential runs.
275
276c-----------------------------------------------------------------------
277c Lecture du tableau des parametres du run (pour la dynamique)
278c-----------------------------------------------------------------------
279
280      if (choix_1.eq.0) then
281
282        write(*,*) 'reading tab_cntrl START_ARCHIVE'
283c
284        ierr = NF_INQ_VARID (nid, "controle", nvarid)
285#ifdef NC_DOUBLE
286        ierr = NF_GET_VAR_DOUBLE(nid, nvarid, tab_cntrl)
287#else
288        ierr = NF_GET_VAR_REAL(nid, nvarid, tab_cntrl)
289#endif
290c
291      else if (choix_1.eq.1) then
292
293        write(*,*) 'reading tab_cntrl START'
294c
295        ierr = NF_INQ_VARID (nid_dyn, "controle", nvarid)
296#ifdef NC_DOUBLE
297        ierr = NF_GET_VAR_DOUBLE(nid_dyn, nvarid, tab_cntrl)
298#else
299        ierr = NF_GET_VAR_REAL(nid_dyn, nvarid, tab_cntrl)
300#endif
301c
302        write(*,*) 'reading tab_cntrl STARTFI'
303c
304        ierr = NF_INQ_VARID (nid_fi, "controle", nvarid)
305#ifdef NC_DOUBLE
306        ierr = NF_GET_VAR_DOUBLE(nid_fi, nvarid, tab_cntrl_bis)
307#else
308        ierr = NF_GET_VAR_REAL(nid_fi, nvarid, tab_cntrl_bis)
309#endif
310c
311        do i=1,50
312          tab_cntrl(i+50)=tab_cntrl_bis(i)
313        enddo
314        write(*,*) 'printing tab_cntrl', tab_cntrl
315        do i=1,100
316          write(*,*) i,tab_cntrl(i)
317        enddo
318       
319        ! Lmodif set to 0 to disable modifications possibility in phyeta0                           
320        write(*,*) 'Reading file START'
321        fichnom = 'start.nc'
322        CALL dynetat0(fichnom,nqmx,vcov,ucov,teta,q,masse,
323     .       ps,phis,time)
324
325        write(*,*) 'Reading file STARTFI'
326        fichnom = 'startfi.nc'
327        CALL phyetat0 (ngridmx,fichnom,tab0,Lmodif,nsoilmx,nqmx,
328     .        day_ini,time,
329     .        tsurf,tsoil,emis,q2,qsurf,   !) ! temporary modif by RDW
330     .        cloudfrac,totalfrac,hice)
331
332        ! copy albedo and soil thermal inertia on (local) physics grid
333        do i=1,ngridmx
334          albfi(i) = albedodat(i)
335          do j=1,nsoilmx
336           ithfi(i,j) = inertiedat(i,j)
337          enddo
338        ! build a surfithfi(:) using 1st layer of ithfi(:), which might
339        ! be needed later on if reinitializing soil thermal inertia
340          surfithfi(i)=ithfi(i,1)
341        enddo
342        ! also copy albedo and soil thermal inertia on (local) dynamics grid
343        ! so that options below can manipulate either (but must then ensure
344        ! to correctly recast things on physics grid)
345        call gr_fi_dyn(1,ngridmx,iip1,jjp1,albfi,alb)
346        call gr_fi_dyn(nsoilmx,ngridmx,iip1,jjp1,ithfi,ith)
347        call gr_fi_dyn(1,ngridmx,iip1,jjp1,surfithfi,surfith)
348     
349      endif
350c-----------------------------------------------------------------------
351c               Initialisation des constantes dynamique
352c-----------------------------------------------------------------------
353
354      kappa = tab_cntrl(9)
355      etot0 = tab_cntrl(12)
356      ptot0 = tab_cntrl(13)
357      ztot0 = tab_cntrl(14)
358      stot0 = tab_cntrl(15)
359      ang0 = tab_cntrl(16)
360      write(*,*) "Newstart: kappa,etot0,ptot0,ztot0,stot0,ang0"
361      write(*,*) kappa,etot0,ptot0,ztot0,stot0,ang0
362
363      ! for vertical coordinate
364      preff=tab_cntrl(18) ! reference surface pressure
365      pa=tab_cntrl(17)  ! reference pressure at which coord is purely pressure
366      !NB: in start_archive files tab_cntrl(17)=tab_cntrl(18)=0
367      write(*,*) "Newstart: preff=",preff," pa=",pa
368      yes=' '
369      do while ((yes.ne.'y').and.(yes.ne.'n'))
370        write(*,*) "Change the values of preff and pa? (y/n)"
371        read(*,fmt='(a)') yes
372      end do
373      if (yes.eq.'y') then
374        write(*,*)"New value of reference surface pressure preff?"
375        write(*,*)"   (for Mars, typically preff=610)"
376        read(*,*) preff
377        write(*,*)"New value of reference pressure pa for purely"
378        write(*,*)"pressure levels (for hybrid coordinates)?"
379        write(*,*)"   (for Mars, typically pa=20)"
380        read(*,*) pa
381      endif
382c-----------------------------------------------------------------------
383c   Lecture du tab_cntrl et initialisation des constantes physiques
384c  - pour start:  Lmodif = 0 => pas de modifications possibles
385c                  (modif dans le tabfi de readfi + loin)
386c  - pour start_archive:  Lmodif = 1 => modifications possibles
387c-----------------------------------------------------------------------
388      if (choix_1.eq.0) then
389         call tabfi (ngridmx,nid,Lmodif,tab0,day_ini,lllm,p_rad,
390     .            p_omeg,p_g,p_cpp,p_mugaz,p_daysec,time)
391      else if (choix_1.eq.1) then
392         Lmodif=1 ! Lmodif set to 1 to allow modifications in phyeta0                           
393         call tabfi (ngridmx,nid_fi,Lmodif,tab0,day_ini,lllm,p_rad,
394     .            p_omeg,p_g,p_cpp,p_mugaz,p_daysec,time)
395      endif
396
397      rad = p_rad
398      omeg = p_omeg
399      g = p_g
400      cpp = p_cpp
401      mugaz = p_mugaz
402      daysec = p_daysec
403
404      kappa = 8.314*1000./(p_mugaz*p_cpp) ! added by RDW
405
406c=======================================================================
407c  INITIALISATIONS DIVERSES
408c=======================================================================
409! Load parameters from run.def file
410      CALL defrun_new( 99, .TRUE. )
411      CALL iniconst
412      CALL inigeom
413      idum=-1
414      idum=0
415
416c Initialisation coordonnees /aires
417c -------------------------------
418! Note: rlatu(:) and rlonv(:) are commons defined in "comgeom.h"
419!       rlatu() and rlonv() are given in radians
420      latfi(1)=rlatu(1)
421      lonfi(1)=0.
422      DO j=2,jjm
423         DO i=1,iim
424            latfi((j-2)*iim+1+i)=rlatu(j)
425            lonfi((j-2)*iim+1+i)=rlonv(i)
426         ENDDO
427      ENDDO
428      latfi(ngridmx)=rlatu(jjp1)
429      lonfi(ngridmx)=0.
430       
431      ! build airefi(), mesh area on physics grid
432      CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,aire,airefi)
433      ! Poles are single points on physics grid
434      airefi(1)=airefi(1)*iim
435      airefi(ngridmx)=airefi(ngridmx)*iim
436
437c=======================================================================
438c   lecture topographie, albedo, inertie thermique, relief sous-maille
439c=======================================================================
440
441      if (choix_1.eq.0) then  ! for start_archive files,
442                              ! where an external "surface.nc" file is needed
443
444c do while((relief(1:3).ne.'mol').AND.(relief(1:3).ne.'pla'))
445c       write(*,*)
446c       write(*,*) 'choix du relief (mola,pla)'
447c       write(*,*) '(Topographie MGS MOLA, plat)'
448c       read(*,fmt='(a3)') relief
449        relief="mola"
450c     enddo
451
452!    First get the correct value of datadir, if not already done:
453        ! default 'datadir' is set in "datafile_mod"
454        call getin("datadir",datadir)
455        write(*,*) 'Available surface data files are:'
456        filestring='ls '//trim(datadir)//' | grep .nc'
457        call system(filestring)
458
459        write(*,*) ''
460        write(*,*) 'Please choose the relevant file',
461     &  ' (e.g. "surface_mars.nc")'
462        write(*,*) ' or "none" to not use any (and set planetary'
463        write(*,*) '  albedo and surface thermal inertia)'
464        read(*,fmt='(a50)') surfacefile
465
466        if (surfacefile.ne."none") then
467
468          CALL datareadnc(relief,surfacefile,phis,alb,surfith,
469     &          zmeaS,zstdS,zsigS,zgamS,ztheS)
470        else
471        ! specific case when not using a "surface.nc" file
472          phis(:,:)=0
473          zmeaS(:,:)=0
474          zstdS(:,:)=0
475          zsigS(:,:)=0
476          zgamS(:,:)=0
477          ztheS(:,:)=0
478         
479          write(*,*) "Enter value of albedo of the bare ground:"
480          read(*,*) alb(1,1)
481          alb(:,:)=alb(1,1)
482         
483          write(*,*) "Enter value of thermal inertia of soil:"
484          read(*,*) surfith(1,1)
485          surfith(:,:)=surfith(1,1)
486       
487        endif ! of if (surfacefile.ne."none")
488
489        CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,phis,phisfi)
490        CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,surfith,surfithfi)
491        CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,alb,albfi)
492        CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,zmeaS,zmea)
493        CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,zstdS,zstd)
494        CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,zsigS,zsig)
495        CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,zgamS,zgam)
496        CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,ztheS,zthe)
497
498      endif ! of if (choix_1.eq.0)
499
500
501c=======================================================================
502c  Lecture des fichiers (start ou start_archive)
503c=======================================================================
504
505      if (choix_1.eq.0) then
506
507        write(*,*) 'Reading file START_ARCHIVE'
508        CALL lect_start_archive(date,tsurf,tsoil,emis,q2,
509     .   t,ucov,vcov,ps,teta,phisold_newgrid,q,qsurf,
510     &   surfith,nid)
511        write(*,*) "OK, read start_archive file"
512        ! copy soil thermal inertia
513        ithfi(:,:)=inertiedat(:,:)
514       
515        ierr= NF_CLOSE(nid)
516
517      else if (choix_1.eq.1) then
518         !do nothing, start and startfi have already been read
519      else
520        CALL exit(1)
521      endif
522
523      dtvr   = daysec/FLOAT(day_step)
524      dtphys   = dtvr * FLOAT(iphysiq)
525
526c=======================================================================
527c
528c=======================================================================
529
530      do ! infinite loop on list of changes
531
532      write(*,*)
533      write(*,*)
534      write(*,*) 'List of possible changes :'
535      write(*,*) '~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~'
536      write(*,*)
537      write(*,*) 'flat : no topography ("aquaplanet")'
538      write(*,*) 'nuketharsis : no Tharsis bulge'
539      write(*,*) 'bilball : uniform albedo and thermal inertia'
540      write(*,*) 'coldspole : cold subsurface and high albedo at S.pole'
541      write(*,*) 'qname : change tracer name'
542      write(*,*) 'q=0 : ALL tracer =zero'
543      write(*,*) 'q=x : give a specific uniform value to one tracer'
544      write(*,*) 'q=profile    : specify a profile for a tracer'
545!      write(*,*) 'ini_q : tracers initialisation for chemistry, water an
546!     $d ice   '
547!      write(*,*) 'ini_q-H2O : tracers initialisation for chemistry and
548!     $ice '
549!      write(*,*) 'ini_q-iceH2O : tracers initialisation for chemistry on
550!     $ly '
551      write(*,*) 'noglacier : Remove tropical H2O ice if |lat|<45'
552      write(*,*) 'watercapn : H20 ice on permanent N polar cap '
553      write(*,*) 'watercaps : H20 ice on permanent S polar cap '
554      write(*,*) 'noacglac  : H2O ice across Noachis Terra'
555      write(*,*) 'oborealis : H2O ice across Vastitas Borealis'
556      write(*,*) 'iceball   : Thick ice layer all over surface'
557      write(*,*) 'wetstart  : start with a wet atmosphere'
558      write(*,*) 'isotherm  : Isothermal Temperatures, wind set to zero'
559      write(*,*) 'radequi   : Earth-like radiative equilibrium temperature
560     $ profile (lat-alt) and winds set to zero'
561      write(*,*) 'coldstart : Start X K above the CO2 frost point and
562     $set wind to zero (assumes 100% CO2)'
563      write(*,*) 'co2ice=0 : remove CO2 polar cap'
564      write(*,*) 'ptot : change total pressure'
565      write(*,*) 'emis : change surface emissivity'
566      write(*,*) 'therm_ini_s : Set soil thermal inertia to reference sur
567     &face values'
568
569        write(*,*)
570        write(*,*) 'Change to perform ?'
571        write(*,*) '   (enter keyword or return to end)'
572        write(*,*)
573
574        read(*,fmt='(a20)') modif
575        if (modif(1:1) .eq. ' ') exit ! exit loop on changes
576
577        write(*,*)
578        write(*,*) trim(modif) , ' : '
579
580c       'flat : no topography ("aquaplanet")'
581c       -------------------------------------
582        if (trim(modif) .eq. 'flat') then
583c         set topo to zero
584          z_reel(1:iip1,1:jjp1)=0
585          phis(1:iip1,1:jjp1)=g*z_reel(1:iip1,1:jjp1)
586          CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,phis,phisfi)
587          write(*,*) 'topography set to zero.'
588          write(*,*) 'WARNING : the subgrid topography parameters',
589     &    ' were not set to zero ! => set calllott to F'                   
590
591c        Choice of surface pressure
592         yes=' '
593         do while ((yes.ne.'y').and.(yes.ne.'n'))
594            write(*,*) 'Do you wish to choose homogeneous surface',
595     &                 'pressure (y) or let newstart interpolate ',
596     &                 ' the previous field  (n)?'
597             read(*,fmt='(a)') yes
598         end do
599         if (yes.eq.'y') then
600           flagps0=.true.
601           write(*,*) 'New value for ps (Pa) ?'
602 201       read(*,*,iostat=ierr) patm
603            if(ierr.ne.0) goto 201
604             write(*,*)
605             write(*,*) ' new ps everywhere (Pa) = ', patm
606             write(*,*)
607             do j=1,jjp1
608               do i=1,iip1
609                 ps(i,j)=patm
610               enddo
611             enddo
612         end if
613
614c       'nuketharsis : no tharsis bulge for Early Mars'
615c       ---------------------------------------------
616        else if (trim(modif) .eq. 'nuketharsis') then
617
618           DO j=1,jjp1       
619              DO i=1,iim
620                 ig=1+(j-2)*iim +i
621                 if(j.eq.1) ig=1
622                 if(j.eq.jjp1) ig=ngridmx
623
624                 fact1=(((rlonv(i)*180./pi)+100)**2 +
625     &                (rlatu(j)*180./pi)**2)/65**2
626                 fact2=exp( -fact1**2.5 )
627
628                 phis(i,j) = phis(i,j) - (phis(i,j)+4000.*g)*fact2
629
630!                 if(phis(i,j).gt.2500.*g)then
631!                    if(rlatu(j)*180./pi.gt.-80.)then ! avoid chopping south polar cap
632!                       phis(i,j)=2500.*g
633!                    endif
634!                 endif
635
636              ENDDO
637           ENDDO
638          CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,phis,phisfi)
639
640
641c       bilball : uniform albedo, thermal inertia
642c       -----------------------------------------
643        else if (trim(modif) .eq. 'bilball') then
644          write(*,*) 'constante albedo and iner.therm:'
645          write(*,*) 'New value for albedo (ex: 0.25) ?'
646 101      read(*,*,iostat=ierr) alb_bb
647          if(ierr.ne.0) goto 101
648          write(*,*)
649          write(*,*) ' uniform albedo (new value):',alb_bb
650          write(*,*)
651
652          write(*,*) 'New value for thermal inertia (eg: 247) ?'
653 102      read(*,*,iostat=ierr) ith_bb
654          if(ierr.ne.0) goto 102
655          write(*,*) 'uniform thermal inertia (new value):',ith_bb
656          DO j=1,jjp1
657             DO i=1,iip1
658                alb(i,j) = alb_bb       ! albedo
659                do isoil=1,nsoilmx
660                  ith(i,j,isoil) = ith_bb       ! thermal inertia
661                enddo
662             END DO
663          END DO
664!          CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,ith,ithfi)
665          CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi)
666          CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,alb,albfi)
667
668c       coldspole : sous-sol de la calotte sud toujours froid
669c       -----------------------------------------------------
670        else if (trim(modif) .eq. 'coldspole') then
671          write(*,*)'new value for the subsurface temperature',
672     &   ' beneath the permanent southern polar cap ? (eg: 141 K)'
673 103      read(*,*,iostat=ierr) tsud
674          if(ierr.ne.0) goto 103
675          write(*,*)
676          write(*,*) ' new value of the subsurface temperature:',tsud
677c         nouvelle temperature sous la calotte permanente
678          do l=2,nsoilmx
679               tsoil(ngridmx,l) =  tsud
680          end do
681
682
683          write(*,*)'new value for the albedo',
684     &   'of the permanent southern polar cap ? (eg: 0.75)'
685 104      read(*,*,iostat=ierr) albsud
686          if(ierr.ne.0) goto 104
687          write(*,*)
688
689c         ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
690c         Option 1:  only the albedo of the pole is modified :   
691          albfi(ngridmx)=albsud
692          write(*,*) 'ig=',ngridmx,'   albedo perennial cap ',
693     &    albfi(ngridmx)
694
695c         ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
696c          Option 2 A haute resolution : coordonnee de la vrai calotte ~   
697c           DO j=1,jjp1
698c             DO i=1,iip1
699c                ig=1+(j-2)*iim +i
700c                if(j.eq.1) ig=1
701c                if(j.eq.jjp1) ig=ngridmx
702c                if ((rlatu(j)*180./pi.lt.-84.).and.
703c     &            (rlatu(j)*180./pi.gt.-91.).and.
704c     &            (rlonv(i)*180./pi.gt.-91.).and.
705c     &            (rlonv(i)*180./pi.lt.0.))         then
706cc    albedo de la calotte permanente fixe a albsud
707c                   alb(i,j)=albsud
708c                   write(*,*) 'lat=',rlatu(j)*180./pi,
709c     &                      ' lon=',rlonv(i)*180./pi
710cc     fin de la condition sur les limites de la calotte permanente
711c                end if
712c             ENDDO
713c          ENDDO
714c      ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
715
716c         CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,alb,albfi)
717
718
719c       ptot : Modification of the total pressure: ice + current atmosphere
720c       -------------------------------------------------------------------
721        else if (trim(modif).eq.'ptot') then
722
723          ! check if we have a co2_ice surface tracer:
724          if (igcm_co2_ice.eq.0) then
725            write(*,*) " No surface CO2 ice !"
726            write(*,*) " only atmospheric pressure will be considered!"
727          endif
728c         calcul de la pression totale glace + atm actuelle
729          patm=0.
730          airetot=0.
731          pcap=0.
732          DO j=1,jjp1
733             DO i=1,iim
734                ig=1+(j-2)*iim +i
735                if(j.eq.1) ig=1
736                if(j.eq.jjp1) ig=ngridmx
737                patm = patm + ps(i,j)*aire(i,j)
738                airetot= airetot + aire(i,j)
739                if (igcm_co2_ice.ne.0) then
740                  !pcap = pcap + aire(i,j)*co2ice(ig)*g
741                  pcap = pcap + aire(i,j)*qsurf(ig,igcm_co2_ice)*g
742                endif
743             ENDDO
744          ENDDO
745          ptoto = pcap + patm
746
747          print*,'Current total pressure at surface (co2 ice + atm) ',
748     &       ptoto/airetot
749
750          print*,'new value?'
751          read(*,*) ptotn
752          ptotn=ptotn*airetot
753          patmn=ptotn-pcap
754          print*,'ptoto,patm,ptotn,patmn'
755          print*,ptoto,patm,ptotn,patmn
756          print*,'Mult. factor for pressure (atm only)', patmn/patm
757          do j=1,jjp1
758             do i=1,iip1
759                ps(i,j)=ps(i,j)*patmn/patm
760             enddo
761          enddo
762
763
764
765c        Correction pour la conservation des traceurs
766         yes=' '
767         do while ((yes.ne.'y').and.(yes.ne.'n'))
768            write(*,*) 'Do you wish to conserve tracer total mass (y)',
769     &              ' or tracer mixing ratio (n) ?'
770             read(*,fmt='(a)') yes
771         end do
772
773         if (yes.eq.'y') then
774           write(*,*) 'OK : conservation of tracer total mass'
775           DO iq =1, nqmx
776             DO l=1,llm
777               DO j=1,jjp1
778                  DO i=1,iip1
779                    q(i,j,l,iq)=q(i,j,l,iq)*patm/patmn
780                  ENDDO
781               ENDDO
782             ENDDO
783           ENDDO
784          else
785            write(*,*) 'OK : conservation of tracer mixing ratio'
786          end if
787
788c        Correction pour la pression au niveau de la mer
789         yes=' '
790         do while ((yes.ne.'y').and.(yes.ne.'n'))
791            write(*,*) 'Do you wish fix pressure at sea level (y)',
792     &              ' or not (n) ?'
793             read(*,fmt='(a)') yes
794         end do
795
796         if (yes.eq.'y') then
797           write(*,*) 'Value?'
798                read(*,*,iostat=ierr) psea
799             DO i=1,iip1
800               DO j=1,jjp1
801                    ps(i,j)=psea
802
803                  ENDDO
804               ENDDO
805                write(*,*) 'psea=',psea
806          else
807            write(*,*) 'no'
808          end if
809
810
811c       emis : change surface emissivity (added by RW)
812c       ----------------------------------------------
813        else if (trim(modif).eq.'emis') then
814
815           print*,'new value?'
816           read(*,*) emisread
817
818           do i=1,ngridmx
819              emis(i)=emisread
820           enddo
821
822c       qname : change tracer name
823c       --------------------------
824        else if (trim(modif).eq.'qname') then
825         yes='y'
826         do while (yes.eq.'y')
827          write(*,*) 'Which tracer name do you want to change ?'
828          do iq=1,nqmx
829            write(*,'(i3,a3,a20)')iq,' : ',trim(tnom(iq))
830          enddo
831          write(*,'(a35,i3)')
832     &            '(enter tracer number; between 1 and ',nqmx
833          write(*,*)' or any other value to quit this option)'
834          read(*,*) iq
835          if ((iq.ge.1).and.(iq.le.nqmx)) then
836            write(*,*)'Change tracer name ',trim(tnom(iq)),' to ?'
837            read(*,*) txt
838            tnom(iq)=txt
839            write(*,*)'Do you want to change another tracer name (y/n)?'
840            read(*,'(a)') yes
841          else
842! inapropiate value of iq; quit this option
843            yes='n'
844          endif ! of if ((iq.ge.1).and.(iq.le.nqmx))
845         enddo ! of do while (yes.ne.'y')
846
847c       q=0 : set tracers to zero
848c       -------------------------
849        else if (trim(modif).eq.'q=0') then
850c          mise a 0 des q (traceurs)
851          write(*,*) 'Tracers set to 0 (1.E-30 in fact)'
852           DO iq =1, nqmx
853             DO l=1,llm
854               DO j=1,jjp1
855                  DO i=1,iip1
856                    q(i,j,l,iq)=1.e-30
857                  ENDDO
858               ENDDO
859             ENDDO
860           ENDDO
861
862c          set surface tracers to zero
863           DO iq =1, nqmx
864             DO ig=1,ngridmx
865                 qsurf(ig,iq)=0.
866             ENDDO
867           ENDDO
868
869c       q=x : initialise tracer manually
870c       --------------------------------
871        else if (trim(modif).eq.'q=x') then
872             write(*,*) 'Which tracer do you want to modify ?'
873             do iq=1,nqmx
874               write(*,*)iq,' : ',trim(tnom(iq))
875             enddo
876             write(*,*) '(choose between 1 and ',nqmx,')'
877             read(*,*) iq
878             write(*,*)'mixing ratio of tracer ',trim(tnom(iq)),
879     &                 ' ? (kg/kg)'
880             read(*,*) val
881             DO l=1,llm
882               DO j=1,jjp1
883                  DO i=1,iip1
884                    q(i,j,l,iq)=val
885                  ENDDO
886               ENDDO
887             ENDDO
888             write(*,*) 'SURFACE value of tracer ',trim(tnom(iq)),
889     &                   ' ? (kg/m2)'
890             read(*,*) val
891             DO ig=1,ngridmx
892                 qsurf(ig,iq)=val
893             ENDDO
894
895c       q=profile : initialize tracer with a given profile
896c       --------------------------------------------------
897        else if (trim(modif) .eq. 'q=profile') then
898             write(*,*) 'Tracer profile will be sought in ASCII file'
899             write(*,*) "'profile_tracer' where 'tracer' is tracer name"
900             write(*,*) "(one value per line in file; starting with"
901             write(*,*) "surface value, the 1st atmospheric layer"
902             write(*,*) "followed by 2nd, etc. up to top of atmosphere)"
903             write(*,*) 'Which tracer do you want to set?'
904             do iq=1,nqmx
905               write(*,*)iq,' : ',trim(tnom(iq))
906             enddo
907             write(*,*) '(choose between 1 and ',nqmx,')'
908             read(*,*) iq
909             if ((iq.lt.1).or.(iq.gt.nqmx)) then
910               ! wrong value for iq, go back to menu
911               write(*,*) "wrong input value:",iq
912               cycle
913             endif
914             ! look for input file 'profile_tracer'
915             txt="profile_"//trim(tnom(iq))
916             open(41,file=trim(txt),status='old',form='formatted',
917     &            iostat=ierr)
918             if (ierr.eq.0) then
919               ! OK, found file 'profile_...', load the profile
920               do l=1,llm+1
921                 read(41,*,iostat=ierr) profile(l)
922                 if (ierr.ne.0) then ! something went wrong
923                   exit ! quit loop
924                 endif
925               enddo
926               if (ierr.eq.0) then
927                 ! initialize tracer values
928                 qsurf(:,iq)=profile(1)
929                 do l=1,llm
930                   q(:,:,l,iq)=profile(l+1)
931                 enddo
932                 write(*,*)'OK, tracer ',trim(tnom(iq)),' initialized ',
933     &                     'using values from file ',trim(txt)
934               else
935                 write(*,*)'problem reading file ',trim(txt),' !'
936                 write(*,*)'No modifications to tracer ',trim(tnom(iq))
937               endif
938             else
939               write(*,*)'Could not find file ',trim(txt),' !'
940               write(*,*)'No modifications to tracer ',trim(tnom(iq))
941             endif
942             
943
944c      wetstart : wet atmosphere with a north to south gradient
945c      --------------------------------------------------------
946       else if (trim(modif) .eq. 'wetstart') then
947        ! check that there is indeed a water vapor tracer
948        if (igcm_h2o_vap.eq.0) then
949          write(*,*) "No water vapour tracer! Can't use this option"
950          stop
951        endif
952          DO l=1,llm
953            DO j=1,jjp1
954              DO i=1,iip1
955                q(i,j,l,igcm_h2o_vap)=150.e-6 * (rlatu(j)+pi/2.) / pi
956              ENDDO
957            ENDDO
958          ENDDO
959
960         write(*,*) 'Water mass mixing ratio at north pole='
961     *               ,q(1,1,1,igcm_h2o_vap)
962         write(*,*) '---------------------------south pole='
963     *               ,q(1,jjp1,1,igcm_h2o_vap)
964
965c      noglacier : remove tropical water ice (to initialize high res sim)
966c      --------------------------------------------------
967        else if (trim(modif) .eq. 'noglacier') then
968           if (igcm_h2o_ice.eq.0) then
969             write(*,*) "No water ice tracer! Can't use this option"
970             stop
971           endif
972           do ig=1,ngridmx
973             j=(ig-2)/iim +2
974              if(ig.eq.1) j=1
975              write(*,*) 'OK: remove surface ice for |lat|<45'
976              if (abs(rlatu(j)*180./pi).lt.45.) then
977                   qsurf(ig,igcm_h2o_ice)=0.
978              end if
979           end do
980
981
982c      watercapn : H20 ice on permanent northern cap
983c      --------------------------------------------------
984        else if (trim(modif) .eq. 'watercapn') then
985           if (igcm_h2o_ice.eq.0) then
986             write(*,*) "No water ice tracer! Can't use this option"
987             stop
988           endif
989
990           DO j=1,jjp1       
991              DO i=1,iim
992                 ig=1+(j-2)*iim +i
993                 if(j.eq.1) ig=1
994                 if(j.eq.jjp1) ig=ngridmx
995
996                 if (rlatu(j)*180./pi.gt.80.) then
997                    qsurf(ig,igcm_h2o_ice)=3.4e3
998                    !do isoil=1,nsoilmx
999                    !   ith(i,j,isoil) = 18000. ! thermal inertia
1000                    !enddo
1001                   write(*,*)'     ==> Ice mesh North boundary (deg)= ',
1002     &                   rlatv(j-1)*180./pi
1003                 end if
1004              ENDDO
1005           ENDDO
1006           CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi)
1007
1008c$$$           do ig=1,ngridmx
1009c$$$             j=(ig-2)/iim +2
1010c$$$              if(ig.eq.1) j=1
1011c$$$              if (rlatu(j)*180./pi.gt.80.) then
1012c$$$
1013c$$$                   qsurf(ig,igcm_h2o_ice)=1.e5
1014c$$$                   qsurf(ig,igcm_h2o_vap)=0.0!1.e5
1015c$$$
1016c$$$                   write(*,*) 'ig=',ig,'    H2O ice mass (kg/m2)= ',
1017c$$$     &              qsurf(ig,igcm_h2o_ice)
1018c$$$
1019c$$$                   write(*,*)'     ==> Ice mesh South boundary (deg)= ',
1020c$$$     &              rlatv(j)*180./pi
1021c$$$                end if
1022c$$$           enddo
1023
1024c      watercaps : H20 ice on permanent southern cap
1025c      -------------------------------------------------
1026        else if (trim(modif) .eq. 'watercaps') then
1027           if (igcm_h2o_ice.eq.0) then
1028              write(*,*) "No water ice tracer! Can't use this option"
1029              stop
1030           endif
1031
1032           DO j=1,jjp1       
1033              DO i=1,iim
1034                 ig=1+(j-2)*iim +i
1035                 if(j.eq.1) ig=1
1036                 if(j.eq.jjp1) ig=ngridmx
1037
1038                 if (rlatu(j)*180./pi.lt.-80.) then
1039                    qsurf(ig,igcm_h2o_ice)=3.4e3
1040                    !do isoil=1,nsoilmx
1041                    !   ith(i,j,isoil) = 18000. ! thermal inertia
1042                    !enddo
1043                   write(*,*)'     ==> Ice mesh South boundary (deg)= ',
1044     &                   rlatv(j-1)*180./pi
1045                 end if
1046              ENDDO
1047           ENDDO
1048           CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi)
1049
1050c$$$           do ig=1,ngridmx
1051c$$$               j=(ig-2)/iim +2
1052c$$$               if(ig.eq.1) j=1
1053c$$$               if (rlatu(j)*180./pi.lt.-80.) then
1054c$$$                  qsurf(ig,igcm_h2o_ice)=1.e5
1055c$$$                  qsurf(ig,igcm_h2o_vap)=0.0 !1.e5
1056c$$$
1057c$$$                  write(*,*) 'ig=',ig,'   H2O ice mass (kg/m2)= ',
1058c$$$     &                 qsurf(ig,igcm_h2o_ice)
1059c$$$                  write(*,*)'     ==> Ice mesh North boundary (deg)= ',
1060c$$$     &                 rlatv(j-1)*180./pi
1061c$$$               end if
1062c$$$           enddo
1063
1064
1065c       noacglac : H2O ice across highest terrain
1066c       --------------------------------------------
1067        else if (trim(modif) .eq. 'noacglac') then
1068           if (igcm_h2o_ice.eq.0) then
1069             write(*,*) "No water ice tracer! Can't use this option"
1070             stop
1071           endif
1072          DO j=1,jjp1       
1073             DO i=1,iim
1074                ig=1+(j-2)*iim +i
1075                if(j.eq.1) ig=1
1076                if(j.eq.jjp1) ig=ngridmx
1077
1078                if(phis(i,j).gt.1000.*g)then
1079                    alb(i,j) = 0.5 ! snow value
1080                    do isoil=1,nsoilmx
1081                       ith(i,j,isoil) = 18000. ! thermal inertia
1082                       ! this leads to rnat set to 'ocean' in physiq.F90
1083                       ! actually no, because it is soil not surface
1084                    enddo
1085                endif
1086             ENDDO
1087          ENDDO
1088          CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi)
1089          CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,alb,albfi)
1090          CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,phis,phisfi)
1091
1092
1093
1094c       oborealis : H2O oceans across Vastitas Borealis
1095c       -----------------------------------------------
1096        else if (trim(modif) .eq. 'oborealis') then
1097           if (igcm_h2o_ice.eq.0) then
1098             write(*,*) "No water ice tracer! Can't use this option"
1099             stop
1100           endif
1101          DO j=1,jjp1       
1102             DO i=1,iim
1103                ig=1+(j-2)*iim +i
1104                if(j.eq.1) ig=1
1105                if(j.eq.jjp1) ig=ngridmx
1106
1107                if(phis(i,j).lt.-4000.*g)then
1108!                if( (phis(i,j).lt.-4000.*g)
1109!     &               .and. (rlatu(j)*180./pi.lt.0.) )then ! south hemisphere only
1110!                    phis(i,j)=-8000.0*g ! proper ocean
1111                    phis(i,j)=-4000.0*g ! scanty ocean
1112
1113                    alb(i,j) = 0.07 ! oceanic value
1114                    do isoil=1,nsoilmx
1115                       ith(i,j,isoil) = 18000. ! thermal inertia
1116                       ! this leads to rnat set to 'ocean' in physiq.F90
1117                       ! actually no, because it is soil not surface
1118                    enddo
1119                endif
1120             ENDDO
1121          ENDDO
1122          CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi)
1123          CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,alb,albfi)
1124          CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,phis,phisfi)
1125
1126c       iborealis : H2O ice in Northern plains
1127c       --------------------------------------
1128        else if (trim(modif) .eq. 'iborealis') then
1129           if (igcm_h2o_ice.eq.0) then
1130             write(*,*) "No water ice tracer! Can't use this option"
1131             stop
1132           endif
1133          DO j=1,jjp1       
1134             DO i=1,iim
1135                ig=1+(j-2)*iim +i
1136                if(j.eq.1) ig=1
1137                if(j.eq.jjp1) ig=ngridmx
1138
1139                if(phis(i,j).lt.-4000.*g)then
1140                   !qsurf(ig,igcm_h2o_ice)=1.e3
1141                   qsurf(ig,igcm_h2o_ice)=241.4 ! to make total 33 kg m^-2
1142                endif
1143             ENDDO
1144          ENDDO
1145          CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi)
1146          CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,alb,albfi)
1147          CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,phis,phisfi)
1148
1149
1150c       oceanball : H2O liquid everywhere
1151c       ----------------------------
1152        else if (trim(modif) .eq. 'oceanball') then
1153           if (igcm_h2o_ice.eq.0) then
1154             write(*,*) "No water ice tracer! Can't use this option"
1155             stop
1156           endif
1157          DO j=1,jjp1       
1158             DO i=1,iim
1159                ig=1+(j-2)*iim +i
1160                if(j.eq.1) ig=1
1161                if(j.eq.jjp1) ig=ngridmx
1162
1163                qsurf(ig,igcm_h2o_vap)=0.0    ! for ocean, this is infinite source
1164                qsurf(ig,igcm_h2o_ice)=0.0
1165                alb(i,j) = 0.07 ! ocean value
1166
1167                do isoil=1,nsoilmx
1168                   ith(i,j,isoil) = 18000. ! thermal inertia
1169                   !ith(i,j,isoil) = 50. ! extremely low for test
1170                   ! this leads to rnat set to 'ocean' in physiq.F90
1171                enddo
1172
1173             ENDDO
1174          ENDDO
1175          CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi)
1176          CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,alb,albfi)
1177          CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,phis,phisfi)
1178
1179c       iceball : H2O ice everywhere
1180c       ----------------------------
1181        else if (trim(modif) .eq. 'iceball') then
1182           if (igcm_h2o_ice.eq.0) then
1183             write(*,*) "No water ice tracer! Can't use this option"
1184             stop
1185           endif
1186          DO j=1,jjp1       
1187             DO i=1,iim
1188                ig=1+(j-2)*iim +i
1189                if(j.eq.1) ig=1
1190                if(j.eq.jjp1) ig=ngridmx
1191
1192                qsurf(ig,igcm_h2o_vap)=-50.    ! for ocean, this is infinite source
1193                qsurf(ig,igcm_h2o_ice)=50.     ! == to 0.05 m of oceanic ice
1194                alb(i,j) = 0.6 ! ice albedo value
1195
1196                do isoil=1,nsoilmx
1197                   !ith(i,j,isoil) = 18000. ! thermal inertia
1198                   ! this leads to rnat set to 'ocean' in physiq.F90
1199                enddo
1200
1201             ENDDO
1202          ENDDO
1203          CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi)
1204          CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,alb,albfi)
1205
1206c       isotherm : Isothermal temperatures and no winds
1207c       -----------------------------------------------
1208        else if (trim(modif) .eq. 'isotherm') then
1209
1210          write(*,*)'Isothermal temperature of the atmosphere,
1211     &           surface and subsurface'
1212          write(*,*) 'Value of this temperature ? :'
1213 203      read(*,*,iostat=ierr) Tiso
1214          if(ierr.ne.0) goto 203
1215
1216          tsurf(1:ngridmx)=Tiso
1217         
1218          tsoil(1:ngridmx,1:nsoilmx)=Tiso
1219         
1220          Tset(1:iip1,1:jjp1,1:llm)=Tiso
1221          flagtset=.true.
1222         
1223          ucov(1:iip1,1:jjp1,1:llm)=0
1224          vcov(1:iip1,1:jjm,1:llm)=0
1225          q2(1:ngridmx,1:nlayermx+1)=0
1226
1227c       radequi : Radiative equilibrium profile of temperatures and no winds
1228c       --------------------------------------------------------------------
1229        else if (trim(modif) .eq. 'radequi') then
1230
1231          write(*,*)'radiative equilibrium temperature profile'       
1232
1233          do ig=1, ngridmx
1234             teque= 335.0-60.0*sin(latfi(ig))*sin(latfi(ig))-
1235     &            10.0*cos(latfi(ig))*cos(latfi(ig))
1236             tsurf(ig) = MAX(220.0,teque)
1237          end do
1238          do l=2,nsoilmx
1239             do ig=1, ngridmx
1240                tsoil(ig,l) = tsurf(ig)
1241             end do
1242          end do
1243          DO j=1,jjp1
1244             DO i=1,iim
1245                Do l=1,llm
1246                   teque=335.-60.0*sin(rlatu(j))*sin(rlatu(j))-
1247     &                  10.0*cos(rlatu(j))*cos(rlatu(j))
1248                   Tset(i,j,l)=MAX(220.0,teque)
1249                end do
1250             end do
1251          end do
1252          flagtset=.true.
1253          ucov(1:iip1,1:jjp1,1:llm)=0
1254          vcov(1:iip1,1:jjm,1:llm)=0
1255          q2(1:ngridmx,1:nlayermx+1)=0
1256
1257c       coldstart : T set 1K above CO2 frost point and no winds
1258c       ------------------------------------------------
1259        else if (trim(modif) .eq. 'coldstart') then
1260
1261          write(*,*)'set temperature of the atmosphere,'
1262     &,'surface and subsurface how many degrees above CO2 frost point?'
1263 204      read(*,*,iostat=ierr) Tabove
1264          if(ierr.ne.0) goto 204
1265
1266            DO j=1,jjp1
1267             DO i=1,iim
1268                ig=1+(j-2)*iim +i
1269                if(j.eq.1) ig=1
1270                if(j.eq.jjp1) ig=ngridmx
1271            tsurf(ig) = (-3167.8)/(log(.01*ps(i,j))-23.23)+Tabove
1272             END DO
1273            END DO
1274          do l=1,nsoilmx
1275            do ig=1, ngridmx
1276              tsoil(ig,l) = tsurf(ig)
1277            end do
1278          end do
1279          DO j=1,jjp1
1280           DO i=1,iim
1281            Do l=1,llm
1282               pp = aps(l) +bps(l)*ps(i,j)
1283               Tset(i,j,l)=(-3167.8)/(log(.01*pp)-23.23)+Tabove
1284            end do
1285           end do
1286          end do
1287
1288          flagtset=.true.
1289          ucov(1:iip1,1:jjp1,1:llm)=0
1290          vcov(1:iip1,1:jjm,1:llm)=0
1291          q2(1:ngridmx,1:nlayermx+1)=0
1292
1293
1294c       co2ice=0 : remove CO2 polar ice caps'
1295c       ------------------------------------------------
1296        else if (trim(modif) .eq. 'co2ice=0') then
1297         ! check that there is indeed a co2_ice tracer ...
1298          if (igcm_co2_ice.ne.0) then
1299           do ig=1,ngridmx
1300              !co2ice(ig)=0
1301              qsurf(ig,igcm_co2_ice)=0
1302              emis(ig)=emis(ngridmx/2)
1303           end do
1304          else
1305            write(*,*) "Can't remove CO2 ice!! (no co2_ice tracer)"
1306          endif
1307       
1308!       therm_ini_s: (re)-set soil thermal inertia to reference surface values
1309!       ----------------------------------------------------------------------
1310
1311        else if (trim(modif).eq.'therm_ini_s') then
1312!          write(*,*)"surfithfi(1):",surfithfi(1)
1313          do isoil=1,nsoilmx
1314            inertiedat(1:ngridmx,isoil)=surfithfi(1:ngridmx)
1315          enddo
1316          write(*,*)'OK: Soil thermal inertia has been reset to referenc
1317     &e surface values'
1318!         write(*,*)"inertiedat(1,1):",inertiedat(1,1)
1319          ithfi(:,:)=inertiedat(:,:)
1320         ! recast ithfi() onto ith()
1321         call gr_fi_dyn(nsoilmx,ngridmx,iip1,jjp1,ithfi,ith)
1322! Check:
1323!         do i=1,iip1
1324!           do j=1,jjp1
1325!             do isoil=1,nsoilmx
1326!               write(77,*) i,j,isoil,"  ",ith(i,j,isoil)
1327!             enddo
1328!           enddo
1329!        enddo
1330
1331
1332        else
1333          write(*,*) '       Unknown (misspelled?) option!!!'
1334        end if ! of if (trim(modif) .eq. '...') elseif ...
1335       
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346       enddo ! of do ! infinite loop on liste of changes
1347
1348 999  continue
1349
1350 
1351c=======================================================================
1352c   Initialisation for cloud fraction and oceanic ice (added by BC 2010)
1353c=======================================================================
1354      DO ig=1, ngridmx
1355         totalfrac(ig)=0.5
1356         DO l=1,llm
1357            cloudfrac(ig,l)=0.5
1358         ENDDO
1359! Ehouarn, march 2012: also add some initialisation for hice
1360         hice(ig)=0.0
1361      ENDDO
1362     
1363c========================================================
1364
1365!      DO ig=1,ngridmx
1366!         IF(tsurf(ig) .LT. 273.16-1.8) THEN
1367!            hice(ig)=(273.16-1.8-tsurf(ig))/(273.16-1.8-240)*1.
1368!            hice(ig)=min(hice(ig),1.0)
1369!         ENDIF
1370!      ENDDO
1371
1372c=======================================================================
1373c   Correct pressure on the new grid (menu 0)
1374c=======================================================================
1375
1376
1377      if ((choix_1.eq.0).and.(.not.flagps0)) then
1378        r = 1000.*8.31/mugaz
1379
1380        do j=1,jjp1
1381          do i=1,iip1
1382             ps(i,j) = ps(i,j) *
1383     .            exp((phisold_newgrid(i,j)-phis(i,j)) /
1384     .                                  (t(i,j,1) * r))
1385          end do
1386        end do
1387
1388c periodicite de ps en longitude
1389        do j=1,jjp1
1390          ps(1,j) = ps(iip1,j)
1391        end do
1392      end if
1393
1394         
1395c=======================================================================
1396c=======================================================================
1397
1398c=======================================================================
1399c    Initialisation de la physique / ecriture de newstartfi :
1400c=======================================================================
1401
1402
1403      CALL inifilr
1404      CALL pression(ip1jmp1, ap, bp, ps, p3d)
1405
1406c-----------------------------------------------------------------------
1407c   Initialisation  pks:
1408c-----------------------------------------------------------------------
1409
1410      CALL exner_hyb(ip1jmp1, ps, p3d, beta, pks, pk, pkf)
1411! Calcul de la temperature potentielle teta
1412
1413      if (flagtset) then
1414          DO l=1,llm
1415             DO j=1,jjp1
1416                DO i=1,iim
1417                   teta(i,j,l) = Tset(i,j,l) * cpp/pk(i,j,l)
1418                ENDDO
1419                teta (iip1,j,l)= teta (1,j,l)
1420             ENDDO
1421          ENDDO
1422      else if (choix_1.eq.0) then
1423         DO l=1,llm
1424            DO j=1,jjp1
1425               DO i=1,iim
1426                  teta(i,j,l) = t(i,j,l) * cpp/pk(i,j,l)
1427               ENDDO
1428               teta (iip1,j,l)= teta (1,j,l)
1429            ENDDO
1430         ENDDO
1431      endif
1432
1433C Calcul intermediaire
1434c
1435      if (choix_1.eq.0) then
1436         CALL massdair( p3d, masse  )
1437c
1438         print *,' ALPHAX ',alphax
1439
1440         DO  l = 1, llm
1441           DO  i    = 1, iim
1442             xppn(i) = aire( i, 1   ) * masse(  i     ,  1   , l )
1443             xpps(i) = aire( i,jjp1 ) * masse(  i     , jjp1 , l )
1444           ENDDO
1445             xpn      = SUM(xppn)/apoln
1446             xps      = SUM(xpps)/apols
1447           DO i   = 1, iip1
1448             masse(   i   ,   1     ,  l )   = xpn
1449             masse(   i   ,   jjp1  ,  l )   = xps
1450           ENDDO
1451         ENDDO
1452      endif
1453      phis(iip1,:) = phis(1,:)
1454
1455      CALL inidissip ( lstardis, nitergdiv, nitergrot, niterh,
1456     *                tetagdiv, tetagrot , tetatemp  )
1457      itau=0
1458      if (choix_1.eq.0) then
1459         day_ini=int(date)
1460      endif
1461c
1462      CALL geopot  ( ip1jmp1, teta  , pk , pks,  phis  , phi   )
1463
1464      CALL caldyn0( itau,ucov,vcov,teta,ps,masse,pk,phis ,
1465     *                phi,w, pbaru,pbarv,day_ini+time )
1466
1467         
1468      CALL dynredem0("restart.nc",day_ini,anneeref,phis,nqmx)
1469      CALL dynredem1("restart.nc",0.0,vcov,ucov,teta,q,nqmx,masse,ps)
1470C
1471C Ecriture etat initial physique
1472C
1473
1474!      do ig=1,ngridmx
1475!         print*,'surface ice in newstart=',qsurf(ig,igcm_h2o_ice)
1476!         print*,'surface liq in newstart=',qsurf(ig,igcm_h2o_vap)
1477!      enddo
1478!      stop
1479
1480
1481      call physdem1(ngridmx,"restartfi.nc",lonfi,latfi,nsoilmx,nqmx,
1482     .              dtphys,float(day_ini),
1483     .              time,tsurf,tsoil,emis,q2,qsurf,
1484     .              airefi,albfi,ithfi,zmea,zstd,zsig,zgam,zthe,
1485     .              cloudfrac,totalfrac,hice,tnom)
1486
1487c=======================================================================
1488c       Formats
1489c=======================================================================
1490
1491   1  FORMAT(//10x,'la valeur de im =',i4,2x,'lue sur le fichier de dema
1492     *rrage est differente de la valeur parametree iim =',i4//)
1493   2  FORMAT(//10x,'la valeur de jm =',i4,2x,'lue sur le fichier de dema
1494     *rrage est differente de la valeur parametree jjm =',i4//)
1495   3  FORMAT(//10x,'la valeur de lllm =',i4,2x,'lue sur le fichier demar
1496     *rage est differente de la valeur parametree llm =',i4//)
1497
1498      end
1499
Note: See TracBrowser for help on using the repository browser.