source: trunk/LMDZ.MARS/libf/dyn3d/newstart.F @ 222

Last change on this file since 222 was 171, checked in by emillour, 14 years ago

Mars GCM:

  • added modifications (from JYC & FGG) to tracer.h & initracer.F for ions
  • minor improvement to newstart.F (q=x option, check that tracer index provided by user is valid).

EM

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