source: trunk/LMDZ.MARS/libf/dynphy_lonlat/phymars/newstart.F @ 2879

Last change on this file since 2879 was 2832, checked in by emillour, 3 years ago

Mars GCM:

  • follow-up to r2826-2829: adapt newstart and lect_start_archive to handle "old" start_archive files where surface co2 ice is stored as "co2ice".

EM

File size: 67.6 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 ioipsl_getincom, only: getin
18      use mod_phys_lmdz_para, only: is_parallel, is_sequential,
19     &                              is_mpi_root, is_omp_root,
20     &                              is_master
21      use infotrac, only: infotrac_init, nqtot, tname
22      use tracer_mod, only: noms, mmol,
23     &                      igcm_dust_number, igcm_dust_mass,
24     &                      igcm_ccn_number, igcm_ccn_mass,
25     &                      igcm_h2o_vap, igcm_h2o_ice, igcm_co2,
26     &                      igcm_hdo_vap, igcm_hdo_ice,
27     &                      igcm_n2, igcm_ar, igcm_o2, igcm_co,
28     &                      igcm_o, igcm_h2
29      use surfdat_h, only: phisfi, z0, zmea, zstd, zsig, zgam, zthe,
30     &                     albedodat, z0_default, qsurf, tsurf,
31     &                     emis, hmons, summit, base, watercap
32      use comsoil_h, only: inertiedat, layer, mlayer, nsoilmx, tsoil
33      use control_mod, only: day_step, iphysiq, anneeref, planet_type
34      use geometry_mod, only: longitude,latitude,cell_area
35      use phyetat0_mod, only: phyetat0
36      use phyredem, only: physdem0, physdem1
37      use iostart, only: open_startphy
38      use dimradmars_mod, only: albedo
39      use dust_param_mod, only: tauscaling
40      use turb_mod, only: q2, wstar
41      use filtreg_mod, only: inifilr
42      USE mod_const_mpi, ONLY: COMM_LMDZ
43      USE comvert_mod, ONLY: ap,bp,pa,preff
44      USE comconst_mod, ONLY: lllm,daysec,dtphys,dtvr,
45     .                  cpp,kappa,rad,omeg,g,r,pi
46      USE serre_mod, ONLY: alphax
47      USE temps_mod, ONLY: day_ini,hour_ini
48      USE ener_mod, ONLY: etot0,ptot0,ztot0,stot0,ang0
49      USE iniphysiq_mod, ONLY: iniphysiq
50      USE exner_hyb_m, ONLY: exner_hyb
51      USE inichim_newstart_mod, ONLY: inichim_newstart
52
53      implicit none
54
55      include "dimensions.h"
56      integer, parameter :: ngridmx = (2+(jjm-1)*iim - 1/jjm)
57      include "paramet.h"
58      include "comgeom2.h"
59      include "comdissnew.h"
60      include "clesph0.h"
61      include "netcdf.inc"
62c=======================================================================
63c   Declarations
64c=======================================================================
65
66c Variables dimension du fichier "start_archive"
67c------------------------------------
68      CHARACTER relief*3
69
70c et autres:
71c----------
72
73c Variables pour les lectures NetCDF des fichiers "start_archive"
74c--------------------------------------------------
75      INTEGER nid_dyn, nid_fi,nid,nvarid
76      INTEGER tab0
77
78      REAL  date
79      REAL p_rad,p_omeg,p_g,p_mugaz,p_daysec
80
81c Variable histoire
82c------------------
83      REAL vcov(iip1,jjm,llm),ucov(iip1,jjp1,llm) ! vents covariants
84      REAL phis(iip1,jjp1)
85      REAL,ALLOCATABLE :: q(:,:,:,:)               ! champs advectes
86
87c autre variables dynamique nouvelle grille
88c------------------------------------------
89      REAL pks(iip1,jjp1)
90      REAL w(iip1,jjp1,llm+1)
91      REAL pbaru(ip1jmp1,llm),pbarv(ip1jm,llm)
92!      REAL dv(ip1jm,llm),du(ip1jmp1,llm)
93!      REAL dh(ip1jmp1,llm),dp(ip1jmp1)
94      REAL phi(iip1,jjp1,llm)
95
96      integer klatdat,klongdat
97      PARAMETER (klatdat=180,klongdat=360)
98
99c Physique sur grille scalaire
100c----------------------------
101      real zmeaS(iip1,jjp1),zstdS(iip1,jjp1)
102      real zsigS(iip1,jjp1),zgamS(iip1,jjp1),ztheS(iip1,jjp1)
103      real hmonsS(iip1,jjp1)
104      real summitS(iip1,jjp1)
105      real baseS(iip1,jjp1)
106      real zavgS(iip1,jjp1)
107      real z0S(iip1,jjp1)
108
109c variable physique
110c------------------
111      REAL tauscadyn(iip1,jjp1) ! dust conversion factor on the dynamics grid
112      real alb(iip1,jjp1),albfi(ngridmx) ! albedos
113      real ith(iip1,jjp1,nsoilmx),ithfi(ngridmx,nsoilmx) ! thermal inertia (3D)
114      real surfith(iip1,jjp1),surfithfi(ngridmx) ! surface thermal inertia (2D)
115!      REAL latfi(ngridmx),lonfi(ngridmx),airefi(ngridmx)
116
117      INTEGER i,j,l,isoil,ig,idum
118      real mugaz ! molar mass of the atmosphere
119
120      integer ierr  !, nbetat
121
122c Variables on the new grid along scalar points
123c------------------------------------------------------
124!      REAL p(iip1,jjp1)
125      REAL t(iip1,jjp1,llm)
126      real phisold_newgrid(iip1,jjp1)
127      REAL :: teta(iip1, jjp1, llm)
128      REAL :: pk(iip1,jjp1,llm)
129      REAL :: pkf(iip1,jjp1,llm)
130      REAL :: ps(iip1, jjp1)
131      REAL :: masse(iip1,jjp1,llm)
132      REAL :: xpn,xps,xppn(iim),xpps(iim)
133      REAL :: p3d(iip1, jjp1, llm+1)
134!      REAL dteta(ip1jmp1,llm)
135
136c Variable de l'ancienne grille
137c------------------------------
138      real time
139      real tab_cntrl(100)
140      real tab_cntrl_bis(100)
141
142c variables diverses
143c-------------------
144      real choix_1 ! ==0 : read start_archive file ; ==1: read start files
145      character*80      fichnom
146      integer Lmodif,iq
147      integer flagthermo, flagh2o
148      character modif*20
149      real tsud,albsud,alb_bb,ith_bb,Tiso
150      real ptoto,pcap,patm,airetot,ptotn,patmn
151!      real ssum
152      character*1 yes
153      logical :: flagiso=.false. ,  flagps0=.false.
154      real val, val2, val3 ! to store temporary variables
155      real :: iceith=2000 ! thermal inertia of subterranean ice
156      real :: iceithN,iceithS ! values of thermal inertias in N & S hemispheres
157      integer iref,jref
158
159      INTEGER :: itau
160      real DoverH !D/H ratio
161     
162      INTEGER :: numvanle
163      character(len=50) :: txt ! to store some text
164      integer :: count
165      real :: profile(llm+1) ! to store an atmospheric profile + surface value
166
167! MONS data:
168      real :: MONS_Hdn(iip1,jjp1) ! Hdn: %WEH=Mass fraction of H2O
169      real :: MONS_d21(iip1,jjp1) ! ice table "depth" (in kg/m2)
170      ! coefficient to apply to convert d21 to 'true' depth (m)
171      real :: MONS_coeff
172      real :: MONS_coeffS ! coeff for southern hemisphere
173      real :: MONS_coeffN ! coeff for northern hemisphere
174!      real,parameter :: icedepthmin=1.e-3 ! Ice begins at most at that depth
175! Reference position for composition
176      real :: latref,lonref,dlatmin,dlonmin
177! Variable used to change composition
178      real :: Svmr,Smmr,Smmr_old,Smmr_new,n,Sn
179      real :: Mair_old,Mair_new,vmr_old,vmr_new
180      real,allocatable :: coefvmr(:)  ! Correction coefficient when changing composition
181      real :: maxq
182      integer :: iloc(1), iqmax
183! sub-grid cloud fraction
184      real :: totcloudfrac(ngridmx)
185
186c sortie visu pour les champs dynamiques
187c---------------------------------------
188!      INTEGER :: visuid
189!      real :: time_step,t_ops,t_wrt
190!      CHARACTER*80 :: visu_file
191
192      cpp    = 744.499 ! for Mars, instead of 1004.70885 (Earth)
193      preff  = 610.    ! for Mars, instead of 101325. (Earth)
194      pa= 20           ! for Mars, instead of 500 (Earth)
195      planet_type="mars"
196
197! initialize "serial/parallel" related stuff:
198! (required because we call tabfi() below, before calling iniphysiq)
199      is_sequential=.true.
200      is_parallel=.false.
201      is_mpi_root=.true.
202      is_omp_root=.true.
203      is_master=.true.
204     
205! Load tracer number and names:
206      call infotrac_init
207! allocate arrays
208      allocate(q(iip1,jjp1,llm,nqtot))
209      allocate(coefvmr(nqtot))
210
211c=======================================================================
212c   Choice of the start file(s) to use
213c=======================================================================
214
215      write(*,*) 'From which kind of files do you want to create new',
216     .  'start and startfi files'
217      write(*,*) '    0 - from a file start_archive'
218      write(*,*) '    1 - from files start and startfi'
219 
220c-----------------------------------------------------------------------
221c   Open file(s) to modify (start or start_archive)
222c-----------------------------------------------------------------------
223
224      DO
225         read(*,*,iostat=ierr) choix_1
226         if ((choix_1 /= 0).OR.(choix_1 /=1)) EXIT
227      ENDDO
228
229c     Open start_archive
230c     ~~~~~~~~~~~~~~~~~~~~~~~~~~
231      if (choix_1.eq.0) then
232
233        write(*,*) 'Creating start files from:'
234        write(*,*) './start_archive.nc'
235        write(*,*)
236        fichnom = 'start_archive.nc'
237        ierr = NF_OPEN (fichnom, NF_NOWRITE,nid)
238        IF (ierr.NE.NF_NOERR) THEN
239          write(6,*)' Problem opening file:',fichnom
240          write(6,*)' ierr = ', ierr
241          CALL ABORT
242        ENDIF
243        tab0 = 50
244        Lmodif = 1
245
246c     OR open start and startfi files
247c     ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
248      else
249        write(*,*) 'Creating start files from:'
250        write(*,*) './start.nc and ./startfi.nc'
251        write(*,*)
252        fichnom = 'start.nc'
253        ierr = NF_OPEN (fichnom, NF_NOWRITE,nid_dyn)
254        IF (ierr.NE.NF_NOERR) THEN
255          write(6,*)' Problem opening file:',fichnom
256          write(6,*)' ierr = ', ierr
257          CALL ABORT
258        ENDIF
259 
260        fichnom = 'startfi.nc'
261        ierr = NF_OPEN (fichnom, NF_NOWRITE,nid_fi)
262        IF (ierr.NE.NF_NOERR) THEN
263          write(6,*)' Problem opening file:',fichnom
264          write(6,*)' ierr = ', ierr
265          CALL ABORT
266        ENDIF
267
268        tab0 = 0
269        Lmodif = 0
270
271      endif
272
273c-----------------------------------------------------------------------
274c Lecture du tableau des parametres du run (pour la dynamique)
275c-----------------------------------------------------------------------
276
277      if (choix_1.eq.0) then
278
279        write(*,*) 'reading tab_cntrl START_ARCHIVE'
280c
281        ierr = NF_INQ_VARID (nid, "controle", nvarid)
282#ifdef NC_DOUBLE
283        ierr = NF_GET_VAR_DOUBLE(nid, nvarid, tab_cntrl)
284#else
285        ierr = NF_GET_VAR_REAL(nid, nvarid, tab_cntrl)
286#endif
287c
288      else if (choix_1.eq.1) then
289
290        write(*,*) 'reading tab_cntrl START'
291c
292        ierr = NF_INQ_VARID (nid_dyn, "controle", nvarid)
293#ifdef NC_DOUBLE
294        ierr = NF_GET_VAR_DOUBLE(nid_dyn, nvarid, tab_cntrl)
295#else
296        ierr = NF_GET_VAR_REAL(nid_dyn, nvarid, tab_cntrl)
297#endif
298c
299        write(*,*) 'reading tab_cntrl STARTFI'
300c
301        ierr = NF_INQ_VARID (nid_fi, "controle", nvarid)
302#ifdef NC_DOUBLE
303        ierr = NF_GET_VAR_DOUBLE(nid_fi, nvarid, tab_cntrl_bis)
304#else
305        ierr = NF_GET_VAR_REAL(nid_fi, nvarid, tab_cntrl_bis)
306#endif
307c
308        do i=1,50
309          tab_cntrl(i+50)=tab_cntrl_bis(i)
310        enddo
311      write(*,*) 'printing tab_cntrl', tab_cntrl
312      do i=1,100
313        write(*,*) i,tab_cntrl(i)
314      enddo
315     
316      endif
317c-----------------------------------------------------------------------
318c               Initialisation des constantes dynamique
319c-----------------------------------------------------------------------
320
321      kappa = tab_cntrl(9)
322      etot0 = tab_cntrl(12)
323      ptot0 = tab_cntrl(13)
324      ztot0 = tab_cntrl(14)
325      stot0 = tab_cntrl(15)
326      ang0 = tab_cntrl(16)
327      write(*,*) "Newstart: kappa,etot0,ptot0,ztot0,stot0,ang0"
328      write(*,*) kappa,etot0,ptot0,ztot0,stot0,ang0
329
330c-----------------------------------------------------------------------
331c   Lecture du tab_cntrl et initialisation des constantes physiques
332c  - pour start:  Lmodif = 0 => pas de modifications possibles
333c                  (modif dans le tabfi de readfi + loin)
334c  - pour start_archive:  Lmodif = 1 => modifications possibles
335c-----------------------------------------------------------------------
336      if (choix_1.eq.0) then
337         ! tabfi requires that input file be first opened by open_startphy(fichnom)
338         fichnom = 'start_archive.nc'
339         call open_startphy(fichnom)
340         call tabfi (nid,Lmodif,tab0,day_ini,lllm,p_rad,
341     .            p_omeg,p_g,p_mugaz,p_daysec,time)
342      else if (choix_1.eq.1) then
343         fichnom = 'startfi.nc'
344         call open_startphy(fichnom)
345         call tabfi (nid_fi,Lmodif,tab0,day_ini,lllm,p_rad,
346     .            p_omeg,p_g,p_mugaz,p_daysec,time)
347      endif
348
349      rad = p_rad
350      omeg = p_omeg
351      g = p_g
352      mugaz = p_mugaz
353      daysec = p_daysec
354
355
356c=======================================================================
357c  INITIALISATIONS DIVERSES
358c=======================================================================
359
360      day_step=180 !?! Note: day_step is a common in "control.h"
361      CALL defrun_new( 99, .TRUE. )
362      dtvr    = daysec/REAL(day_step)
363      CALL iniconst
364      CALL inigeom
365      idum=-1
366      idum=0
367
368! Initialize the physics
369         CALL iniphysiq(iim,jjm,llm,
370     &                  (jjm-1)*iim+2,comm_lmdz,
371     &                  daysec,day_ini,dtphys,
372     &                  rlatu,rlatv,rlonu,rlonv,
373     &                  aire,cu,cv,rad,g,r,cpp,
374     &                  1)
375
376c=======================================================================
377c   lecture topographie, albedo, inertie thermique, relief sous-maille
378c=======================================================================
379
380      if (choix_1.ne.1) then  ! pour ne pas avoir besoin du fichier
381                              ! surface.dat dans le cas des start
382
383c do while((relief(1:3).ne.'mol').AND.(relief(1:3).ne.'pla'))
384c       write(*,*)
385c       write(*,*) 'choix du relief (mola,pla)'
386c       write(*,*) '(Topographie MGS MOLA, plat)'
387c       read(*,fmt='(a3)') relief
388        relief="mola"
389c     enddo
390
391      CALL datareadnc(relief,phis,alb,surfith,z0S,
392     &          zmeaS,zstdS,zsigS,zgamS,ztheS,
393     &          hmonsS,summitS,baseS,zavgS)
394
395      CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,phis,phisfi)
396!      CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,ith,ithfi)
397      CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,surfith,surfithfi)
398      CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,alb,albfi)
399      CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,z0S,z0)
400      CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,zmeaS,zmea)
401      CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,zstdS,zstd)
402      CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,zsigS,zsig)
403      CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,zgamS,zgam)
404      CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,ztheS,zthe)
405      CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,hmonsS,hmons)
406      CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,summitS,summit)
407      CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,baseS,base)
408!      CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,zavgS,zavg)
409
410      endif ! of if (choix_1.ne.1)
411
412
413c=======================================================================
414c  Lecture des fichiers (start ou start_archive)
415c=======================================================================
416
417      if (choix_1.eq.0) then
418
419        write(*,*) 'Reading file START_ARCHIVE'
420        CALL lect_start_archive(ngridmx,llm,nqtot,
421     &   date,tsurf,tsoil,albedo,emis,q2,
422     &   t,ucov,vcov,ps,teta,phisold_newgrid,q,qsurf,
423     &   tauscaling,totcloudfrac,surfith,nid,watercap)
424        write(*,*) "OK, read start_archive file"
425        ! copy soil thermal inertia
426        ithfi(:,:)=inertiedat(:,:)
427       
428        ierr= NF_CLOSE(nid)
429
430      else if (choix_1.eq.1) then !  c'est l'appel a tabfi de phyeta0 qui
431                                  !  permet de changer les valeurs du
432                                  !  tab_cntrl Lmodif=1
433        tab0=0
434        Lmodif=1 ! Lmodif set to 1 to allow modifications in phyeta0                           
435        write(*,*) 'Reading file START'
436        fichnom = 'start.nc'
437        CALL dynetat0(fichnom,vcov,ucov,teta,q,masse,
438     &       ps,phis,time)
439
440        write(*,*) 'Reading file STARTFI'
441        fichnom = 'startfi.nc'
442        CALL phyetat0 (fichnom,tab0,Lmodif,nsoilmx,ngridmx,llm,nqtot,
443     &        day_ini,time,tsurf,tsoil,albedo,emis,
444     &        q2,qsurf,tauscaling,totcloudfrac,
445     &        wstar,watercap)
446       
447        ! copy albedo and soil thermal inertia
448        do i=1,ngridmx
449          albfi(i) = albedodat(i)
450          do j=1,nsoilmx
451           ithfi(i,j) = inertiedat(i,j)
452          enddo
453        ! build a surfithfi(:) using 1st layer of ithfi(:), which might
454        ! be neede later on if reinitializing soil thermal inertia
455          surfithfi(i)=ithfi(i,1)
456        enddo
457
458      else
459        CALL exit(1)
460      endif
461
462      dtvr   = daysec/REAL(day_step)
463      dtphys   = dtvr * REAL(iphysiq)
464
465c=======================================================================
466c
467c=======================================================================
468! If tracer names follow 'old' convention (q01, q02, ...) then
469! rename them
470      count=0
471      do iq=1,nqtot
472        txt=" "
473        write(txt,'(a1,i2.2)') 'q',iq
474        if (txt.eq.tname(iq)) then
475          count=count+1
476        endif
477      enddo ! of do iq=1,nqtot
478     
479      ! initialize tracer names noms(:) and indexes (igcm_co2, igcm_h2o_vap, ...)
480      call initracer(ngridmx,nqtot,qsurf)
481     
482      if (count.eq.nqtot) then
483        write(*,*) 'Newstart: updating tracer names'
484        ! copy noms(:) to tname(:) to have matching tracer names in physics
485        ! and dynamics
486        tname(1:nqtot)=noms(1:nqtot)
487      endif
488
489c=======================================================================
490c
491c=======================================================================
492
493      do ! infinite loop on list of changes
494
495      write(*,*)
496      write(*,*)
497      write(*,*) 'List of possible changes :'
498      write(*,*) '~~~~~~~~~~~~~~~~~~~~~~~~~~'
499      write(*,*)
500      write(*,*) 'flat         : no topography ("aquaplanet")'
501      write(*,*) 'bilball      : uniform albedo and thermal inertia'
502      write(*,*) 'z0           : set a uniform surface roughness length'
503      write(*,*) 'coldspole    : cold subsurface and high albedo at
504     $ S.Pole'
505      write(*,*) 'qname        : change tracer name'
506      write(*,*) 'q=0          : ALL tracer =zero'
507      write(*,*) 'q=factor     : change tracer value by a multiplicative
508     & factor'
509      write(*,*) 'q=x          : give a specific uniform value to one
510     $ tracer'
511      write(*,*) 'q=profile    : specify a profile for a tracer'
512      write(*,*) 'freedust     : rescale dust to a true value'
513      write(*,*) 'ini_q        : tracers initialization for chemistry
514     $ and water vapour'
515      write(*,*) 'ini_q-h2o    : tracers initialization for chemistry
516     $ only'
517      write(*,*) 'composition  : change atm main composition: CO2,N2,Ar,
518     $ O2,CO'
519      write(*,*) 'inihdo       : initialize HDO'
520      write(*,*) 'ini_h2osurf  : reinitialize surface water ice '
521      write(*,*) 'noglacier    : Remove tropical H2O ice if |lat|<45'
522      write(*,*) 'watercapn    : H20 ice on permanent N polar cap '
523      write(*,*) 'watercaps    : H20 ice on permanent S polar cap '
524      write(*,*) 'wetstart     : start with a wet atmosphere'
525      write(*,*) 'isotherm     : Isothermal Temperatures, wind set to
526     $ zero'
527      write(*,*) 'co2ice=0     : remove CO2 polar cap i.e.
528     $ qsurf(co2)=0 '
529      write(*,*) 'ptot         : change total pressure'
530      write(*,*) 'therm_ini_s  : set soil thermal inertia to reference
531     $ surface values'
532      write(*,*) 'subsoilice_n : put deep underground ice layer in
533     $ northern hemisphere'
534      write(*,*) 'subsoilice_s : put deep underground ice layer in
535     $ southern hemisphere'
536      write(*,*) 'mons_ice     : put underground ice layer according
537     $ to MONS derived data'
538
539        write(*,*)
540        write(*,*) 'Change to perform ?'
541        write(*,*) '   (enter keyword or return to end)'
542        write(*,*)
543
544        read(*,fmt='(a20)') modif
545        if (modif(1:1) .eq. ' ') exit ! exit loop on changes
546
547        write(*,*)
548        write(*,*) trim(modif) , ' : '
549
550c       'flat : no topography ("aquaplanet")'
551c       -------------------------------------
552        if (trim(modif) .eq. 'flat') then
553c         set topo to zero
554          phis(:,:)=0
555          CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,phis,phisfi)
556          write(*,*) 'topography set to zero.'
557          write(*,*) 'WARNING : the subgrid topography parameters',
558     &    ' were not set to zero ! => set calllott to F'                   
559
560c        Choice for surface pressure
561         yes=' '
562         do while ((yes.ne.'y').and.(yes.ne.'n'))
563            write(*,*) 'Do you wish to choose homogeneous surface',
564     &                 'pressure (y) or let newstart interpolate ',
565     &                 ' the previous field  (n)?'
566             read(*,fmt='(a)') yes
567         end do
568         if (yes.eq.'y') then
569           flagps0=.true.
570           write(*,*) 'New value for ps (Pa) ?'
571 201       read(*,*,iostat=ierr) patm
572            if(ierr.ne.0) goto 201
573             write(*,*)
574             write(*,*) ' new ps everywhere (Pa) = ', patm
575             write(*,*)
576             do j=1,jjp1
577               do i=1,iip1
578                 ps(i,j)=patm
579               enddo
580             enddo
581         end if
582
583c       bilball : albedo, inertie thermique uniforme
584c       --------------------------------------------
585        else if (trim(modif) .eq. 'bilball') then
586          write(*,*) 'constante albedo and iner.therm:'
587          write(*,*) 'New value for albedo (ex: 0.25) ?'
588 101      read(*,*,iostat=ierr) alb_bb
589          if(ierr.ne.0) goto 101
590          write(*,*)
591          write(*,*) ' uniform albedo (new value):',alb_bb
592          write(*,*)
593
594          write(*,*) 'New value for thermal inertia (eg: 247) ?'
595 102      read(*,*,iostat=ierr) ith_bb
596          if(ierr.ne.0) goto 102
597          write(*,*) 'uniform thermal inertia (new value):',ith_bb
598          DO j=1,jjp1
599             DO i=1,iip1
600                alb(i,j) = alb_bb       ! albedo
601                do isoil=1,nsoilmx
602                  ith(i,j,isoil) = ith_bb       ! thermal inertia
603                enddo
604             END DO
605          END DO
606!          CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,ith,ithfi)
607          CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi)
608          CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,alb,albfi)
609       
610         ! also reset surface roughness length to default value
611         write(*,*) 'surface roughness length set to:',z0_default,' m'
612         z0(:)=z0_default
613
614!       z0 : set surface roughness length to a constant value
615!       -----------------------------------------------------
616        else if (trim(modif) .eq. 'z0') then
617          write(*,*) 'set a uniform surface roughness length'
618          write(*,*) ' value for z0_default (ex: ',z0_default,')?'
619          ierr=1
620          do while (ierr.ne.0)
621            read(*,*,iostat=ierr) z0_default
622          enddo
623          z0(:)=z0_default
624
625c       coldspole : sous-sol de la calotte sud toujours froid
626c       -----------------------------------------------------
627        else if (trim(modif) .eq. 'coldspole') then
628          write(*,*)'new value for the subsurface temperature',
629     &   ' beneath the permanent southern polar cap ? (eg: 141 K)'
630 103      read(*,*,iostat=ierr) tsud
631          if(ierr.ne.0) goto 103
632          write(*,*)
633          write(*,*) ' new value of the subsurface temperature:',tsud
634c         nouvelle temperature sous la calotte permanente
635          do l=2,nsoilmx
636               tsoil(ngridmx,l) =  tsud
637          end do
638
639
640          write(*,*)'new value for the albedo',
641     &   'of the permanent southern polar cap ? (eg: 0.75)'
642 104      read(*,*,iostat=ierr) albsud
643          if(ierr.ne.0) goto 104
644          write(*,*)
645
646c         ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
647c         Option 1:  only the albedo of the pole is modified :   
648          albfi(ngridmx)=albsud
649          write(*,*) 'ig=',ngridmx,'   albedo perennial cap ',
650     &    albfi(ngridmx)
651
652c         ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
653c          Option 2 A haute resolution : coordonnee de la vrai calotte ~   
654c           DO j=1,jjp1
655c             DO i=1,iip1
656c                ig=1+(j-2)*iim +i
657c                if(j.eq.1) ig=1
658c                if(j.eq.jjp1) ig=ngridmx
659c                if ((rlatu(j)*180./pi.lt.-84.).and.
660c     &            (rlatu(j)*180./pi.gt.-91.).and.
661c     &            (rlonv(i)*180./pi.gt.-91.).and.
662c     &            (rlonv(i)*180./pi.lt.0.))         then
663cc    albedo de la calotte permanente fixe a albsud
664c                   alb(i,j)=albsud
665c                   write(*,*) 'lat=',rlatu(j)*180./pi,
666c     &                      ' lon=',rlonv(i)*180./pi
667cc     fin de la condition sur les limites de la calotte permanente
668c                end if
669c             ENDDO
670c          ENDDO
671c      ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
672
673c         CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,alb,albfi)
674
675
676c       ptot : Modification of the total pressure: ice + current atmosphere
677c       -------------------------------------------------------------------
678        else if (trim(modif) .eq. 'ptot') then
679
680c         calcul de la pression totale glace + atm actuelle
681          patm=0.
682          airetot=0.
683          pcap=0.
684          DO j=1,jjp1
685             DO i=1,iim
686                ig=1+(j-2)*iim +i
687                if(j.eq.1) ig=1
688                if(j.eq.jjp1) ig=ngridmx
689                patm = patm + ps(i,j)*aire(i,j)
690                airetot= airetot + aire(i,j)
691                pcap = pcap + aire(i,j)*qsurf(ig,igcm_co2)*g
692             ENDDO
693          ENDDO
694          ptoto = pcap + patm
695
696          print*,'Current total pressure at surface (co2 ice + atm) ',
697     &       ptoto/airetot
698
699          print*,'new value?'
700          read(*,*) ptotn
701          ptotn=ptotn*airetot
702          patmn=ptotn-pcap
703          print*,'ptoto,patm,ptotn,patmn'
704          print*,ptoto,patm,ptotn,patmn
705          print*,'Mult. factor for pressure (atm only)', patmn/patm
706          do j=1,jjp1
707             do i=1,iip1
708                ps(i,j)=ps(i,j)*patmn/patm
709             enddo
710          enddo
711
712c        Correction pour la conservation des traceurs
713         yes=' '
714         do while ((yes.ne.'y').and.(yes.ne.'n'))
715            write(*,*) 'Do you wish to conserve tracer total mass (y)',
716     &              ' or tracer mixing ratio (n) ?'
717             read(*,fmt='(a)') yes
718         end do
719
720         if (yes.eq.'y') then
721           write(*,*) 'OK : conservation of tracer total mass'
722           DO iq =1, nqtot
723             DO l=1,llm
724               DO j=1,jjp1
725                  DO i=1,iip1
726                    q(i,j,l,iq)=q(i,j,l,iq)*patm/patmn
727                  ENDDO
728               ENDDO
729             ENDDO
730           ENDDO
731          else
732            write(*,*) 'OK : conservation of tracer mixing ratio'
733          end if
734
735c       qname : change tracer name
736c       --------------------------
737        else if (trim(modif).eq.'qname') then
738         yes='y'
739         do while (yes.eq.'y')
740          write(*,*) 'Which tracer name do you want to change ?'
741          do iq=1,nqtot
742            write(*,'(i3,a3,a20)')iq,' : ',trim(tname(iq))
743          enddo
744          write(*,'(a35,i3)')
745     &            '(enter tracer number; between 1 and ',nqtot
746          write(*,*)' or any other value to quit this option)'
747          read(*,*) iq
748          if ((iq.ge.1).and.(iq.le.nqtot)) then
749            write(*,*)'Change tracer name ',trim(tname(iq)),' to ?'
750            read(*,*) txt
751            tname(iq)=txt
752            write(*,*)'Do you want to change another tracer name (y/n)?'
753            read(*,'(a)') yes
754          else
755! inapropiate value of iq; quit this option
756            yes='n'
757          endif ! of if ((iq.ge.1).and.(iq.le.nqtot))
758         enddo ! of do while (yes.ne.'y')
759
760c       q=0 : set tracers to zero
761c       -------------------------
762        else if (trim(modif) .eq. 'q=0') then
763c          mise a 0 des q (traceurs)
764          write(*,*) 'Tracers set to 0 (1.E-30 in fact)'
765           DO iq =1, nqtot
766             DO l=1,llm
767               DO j=1,jjp1
768                  DO i=1,iip1
769                    q(i,j,l,iq)=1.e-30
770                  ENDDO
771               ENDDO
772             ENDDO
773           ENDDO
774
775c          set surface tracers to zero
776           DO iq =1, nqtot
777             DO ig=1,ngridmx
778                 qsurf(ig,iq)=0.
779             ENDDO
780           ENDDO
781
782c       q=factor : change value of tracer by a multiplicative factor
783c       ------------------------------------------------------------
784        else if (trim(modif) .eq. 'q=factor') then
785             write(*,*) 'Which tracer do you want to modify ?'
786             do iq=1,nqtot
787               write(*,*)iq,' : ',trim(tname(iq))
788             enddo
789             write(*,*) '(choose between 1 and ',nqtot,')'
790             read(*,*) iq
791             if ((iq.lt.1).or.(iq.gt.nqtot)) then
792               ! wrong value for iq, go back to menu
793               write(*,*) "wrong input value:",iq
794               cycle
795             endif
796             write(*,*)"factor to multiply current mixing ratio by?"
797             read(*,*) val
798             
799             q(1:iip1,1:jjp1,1:llm,iq)=q(1:iip1,1:jjp1,1:llm,iq)*val
800             qsurf(1:ngridmx,iq)=qsurf(1:ngridmx,iq)*val
801
802c       q=x : initialise tracer manually
803c       --------------------------------
804        else if (trim(modif) .eq. 'q=x') then
805             write(*,*) 'Which tracer do you want to modify ?'
806             do iq=1,nqtot
807               write(*,*)iq,' : ',trim(tname(iq))
808             enddo
809             write(*,*) '(choose between 1 and ',nqtot,')'
810             read(*,*) iq
811             if ((iq.lt.1).or.(iq.gt.nqtot)) then
812               ! wrong value for iq, go back to menu
813               write(*,*) "wrong input value:",iq
814               cycle
815             endif
816             write(*,*)'mixing ratio of tracer ',trim(tname(iq)),
817     &                 ' ? (kg/kg)'
818             read(*,*) val
819             DO l=1,llm
820               DO j=1,jjp1
821                  DO i=1,iip1
822                    q(i,j,l,iq)=val
823                  ENDDO
824               ENDDO
825             ENDDO
826             write(*,*) 'SURFACE value of tracer ',trim(tname(iq)),
827     &                   ' ? (kg/m2)'
828             read(*,*) val
829             DO ig=1,ngridmx
830                 qsurf(ig,iq)=val
831             ENDDO
832
833c       q=profile : initialize tracer with a given profile
834c       --------------------------------------------------
835        else if (trim(modif) .eq. 'q=profile') then
836             write(*,*) 'Tracer profile will be sought in ASCII file'
837             write(*,*) "'profile_tracer' where 'tracer' is tracer name"
838             write(*,*) "(one value per line in file; starting with"
839             write(*,*) "surface value, the 1st atmospheric layer"
840             write(*,*) "followed by 2nd, etc. up to top of atmosphere)"
841             write(*,*) 'Which tracer do you want to set?'
842             do iq=1,nqtot
843               write(*,*)iq,' : ',trim(tname(iq))
844             enddo
845             write(*,*) '(choose between 1 and ',nqtot,')'
846             read(*,*) iq
847             if ((iq.lt.1).or.(iq.gt.nqtot)) then
848               ! wrong value for iq, go back to menu
849               write(*,*) "wrong input value:",iq
850               cycle
851             endif
852             ! look for input file 'profile_tracer'
853             txt="profile_"//trim(tname(iq))
854             open(41,file=trim(txt),status='old',form='formatted',
855     &            iostat=ierr)
856             if (ierr.eq.0) then
857               ! OK, found file 'profile_...', load the profile
858               do l=1,llm+1
859                 read(41,*,iostat=ierr) profile(l)
860                 if (ierr.ne.0) then ! something went wrong
861                   exit ! quit loop
862                 endif
863               enddo
864               if (ierr.eq.0) then
865                 ! initialize tracer values
866                 qsurf(:,iq)=profile(1)
867                 do l=1,llm
868                   q(:,:,l,iq)=profile(l+1)
869                 enddo
870                 write(*,*)'OK, tracer ',trim(tname(iq)),
871     &               ' initialized ','using values from file ',trim(txt)
872               else
873                 write(*,*)'problem reading file ',trim(txt),' !'
874                 write(*,*)'No modifications to tracer ',trim(tname(iq))
875               endif
876             else
877               write(*,*)'Could not find file ',trim(txt),' !'
878               write(*,*)'No modifications to tracer ',trim(tname(iq))
879             endif
880             
881c       convert dust from virtual to true values
882c       --------------------------------------------------
883        else if (trim(modif) .eq. 'freedust') then
884         if (minval(tauscaling) .lt. 0) then
885           write(*,*) 'WARNING conversion factor negative'
886           write(*,*) 'This is probably because it was not present
887     &in the file'
888           write(*,*) 'A constant conversion is used instead.'
889           tauscaling(:) = 1.e-3
890         endif
891         CALL gr_fi_dyn(1,ngridmx,iip1,jjp1,tauscaling,tauscadyn)
892          do l=1,llm
893            do j=1,jjp1
894              do i=1,iip1
895                if (igcm_dust_number .ne. 0)
896     &            q(i,j,l,igcm_dust_number) =
897     &            q(i,j,l,igcm_dust_number) * tauscadyn(i,j)
898                if (igcm_dust_mass .ne. 0)
899     &            q(i,j,l,igcm_dust_mass) =
900     &            q(i,j,l,igcm_dust_mass) * tauscadyn(i,j)
901                if (igcm_ccn_number .ne. 0)
902     &            q(i,j,l,igcm_ccn_number) =
903     &            q(i,j,l,igcm_ccn_number) * tauscadyn(i,j)
904                if (igcm_ccn_mass .ne. 0)
905     &            q(i,j,l,igcm_ccn_mass) =
906     &            q(i,j,l,igcm_ccn_mass) * tauscadyn(i,j)
907              end do
908            end do
909          end do
910
911          tauscaling(:) = 1.
912
913         ! We want to have the very same value at lon -180 and lon 180
914          do l = 1,llm
915             do j = 1,jjp1
916                do iq = 1,nqtot
917                   q(iip1,j,l,iq) = q(1,j,l,iq)
918                end do
919             end do
920          end do
921
922          write(*,*) 'done rescaling to true vale'
923
924c       ini_q : Initialize tracers for chemistry
925c       -----------------------------------------------
926        else if (trim(modif) .eq. 'ini_q') then
927          flagh2o    = 1
928          flagthermo = 0
929          yes=' '
930c         For more than 32 layers, possible to initiate thermosphere only     
931          if (llm.gt.32) then
932            do while ((yes.ne.'y').and.(yes.ne.'n'))
933            write(*,*)'',
934     &     'initialisation for thermosphere only? (y/n)'
935            read(*,fmt='(a)') yes
936            if (yes.eq.'y') then
937            flagthermo=1
938            else
939            flagthermo=0
940            endif
941            enddo 
942          endif
943         
944          call inichim_newstart(ngridmx, nqtot, q, qsurf, ps,
945     &                          flagh2o, flagthermo)
946
947         ! We want to have the very same value at lon -180 and lon 180
948          do l = 1,llm
949             do j = 1,jjp1
950                do iq = 1,nqtot
951                   q(iip1,j,l,iq) = q(1,j,l,iq)
952                end do
953             end do
954          end do
955
956          write(*,*) 'inichim_newstart: chemical species and
957     $ water vapour initialised'
958
959c       ini_q-h2o : as above except for the water vapour tracer
960c       ------------------------------------------------------
961        else if (trim(modif) .eq. 'ini_q-h2o') then
962          flagh2o    = 0
963          flagthermo = 0
964          yes=' '
965          ! for more than 32 layers, possible to initiate thermosphere only     
966          if(llm.gt.32) then
967            do while ((yes.ne.'y').and.(yes.ne.'n'))
968            write(*,*)'',
969     &      'initialisation for thermosphere only? (y/n)'
970            read(*,fmt='(a)') yes
971            if (yes.eq.'y') then
972            flagthermo=1
973            else
974            flagthermo=0
975            endif
976            enddo
977          endif
978
979          call inichim_newstart(ngridmx, nqtot, q, qsurf, ps,
980     &                          flagh2o, flagthermo)
981
982         ! We want to have the very same value at lon -180 and lon 180
983          do l = 1,llm
984             do j = 1,jjp1
985                do iq = 1,nqtot
986                   q(iip1,j,l,iq) = q(1,j,l,iq)
987                end do
988             end do
989          end do
990
991          write(*,*) 'inichim_newstart: chemical species initialised
992     $ (except water vapour)'
993
994c      inihdo : initialize HDO with user D/H value
995c      --------------------------------------------------------
996       else if (trim(modif) .eq. 'inihdo') then
997        ! check that there is indeed a water vapor tracer
998        if (igcm_h2o_vap.eq.0) then
999          write(*,*) "No water vapour tracer! Can't use this option"
1000          stop
1001        endif
1002
1003         write(*,*)'Input D/H ratio (in SMOW)'
1004         write(*,*)'If value is <0 then HDO=H2O'
1005 303     read(*,*, iostat=ierr) DoverH
1006         if(ierr.ne.0) goto 303
1007
1008        DoverH = DoverH * 2 * 155.76e-6
1009        if (DoverH.lt.0.0) then
1010           DoverH = 1.
1011        endif
1012        !D/H (SMOW) = 155.76e-6 so HDO/H2O is twice that
1013
1014           do ig=1,ngridmx
1015           qsurf(ig,igcm_h2o_ice)=max(0.,qsurf(ig,igcm_h2o_ice))
1016           end do
1017
1018        ! Update the hdo tracers
1019         q(1:iip1,1:jjp1,1:llm,igcm_hdo_vap)
1020     &      =q(1:iip1,1:jjp1,1:llm,igcm_h2o_vap)* DoverH
1021         q(1:iip1,1:jjp1,1:llm,igcm_hdo_ice)
1022     &      =q(1:iip1,1:jjp1,1:llm,igcm_h2o_ice)* DoverH
1023
1024         qsurf(1:ngridmx,igcm_hdo_ice)
1025     &      =qsurf(1:ngridmx,igcm_h2o_ice)*DoverH
1026
1027
1028
1029c      composition : change main composition: CO2,N2,Ar,O2,CO (FF 03/2014)
1030c      --------------------------------------------------------
1031       else if (trim(modif) .eq. 'composition') then
1032          write(*,*) "Lat (degN)  lon (degE) of the reference site ?"
1033          write(*,*) "e.g. MSL : -4.5  137.  "
1034 301      read(*,*,iostat=ierr) latref, lonref
1035          if(ierr.ne.0) goto 301
1036
1037
1038        !  Select GCM point close to reference site
1039          dlonmin =90.
1040          DO i=1,iip1-1
1041             if (abs(rlonv(i)*180./pi -lonref).lt.dlonmin)then
1042                iref=i
1043                dlonmin=abs(rlonv(i)*180./pi -lonref)
1044             end if   
1045          ENDDO
1046          dlatmin =45.
1047          DO j=1,jjp1
1048             if (abs(rlatu(j)*180./pi -latref).lt.dlatmin)then
1049                jref=j
1050                dlatmin=abs(rlatu(j)*180./pi -latref)
1051             end if   
1052          ENDDO
1053          write(*,*) "In GCM : lat= " ,  rlatu(jref)*180./pi
1054          write(*,*) "In GCM : lon= " ,  rlonv(iref)*180./pi
1055          write(*,*)
1056
1057        ! Compute air molar mass at reference site
1058          Smmr=0.
1059          Sn = 0.
1060          write(*,*) 'igcm_co2 = ', igcm_co2
1061          write(*,*) 'igcm_n2 = ', igcm_n2
1062          write(*,*) 'igcm_ar = ', igcm_ar
1063          write(*,*) 'igcm_o2 = ', igcm_o2
1064          write(*,*) 'igcm_co = ', igcm_co
1065          write(*,*) noms
1066          do iq=1,nqtot
1067             if ((iq.eq.igcm_co2).or.(iq.eq.igcm_n2)
1068     &      .or. (iq.eq.igcm_ar).or.(iq.eq.igcm_o2)
1069     &      .or. (iq.eq.igcm_co)) then
1070                 Smmr=Smmr+q(iref,jref,1,iq)
1071                 Sn=Sn+q(iref,jref,1,iq)/mmol(iq)
1072             end if
1073          end do
1074        ! Special case : implicit non-co2 gases ! JN 11/2019
1075          if ((igcm_n2.eq.0) .or. (igcm_ar.eq.0)) then
1076           write(*,*) "Warning : non-co2 gases are implicit :  "
1077           write(*,*) "At reference site :  "
1078       !    write(*,*) "q= ", q(iref, jref, 1,igcm_co2)
1079           write(*,*) "Sum of mass mix. ratio (ie MMR(co2))=",Smmr
1080           Mair_old = 44.0*Smmr  + 33.87226017157708*(1-Smmr)
1081     
1082       !  33.87226017157708 is the
1083       !   molar mass of non-co2 atmosphere measured by MSL at Ls ~184
1084       
1085          else
1086            ! Assume co2/n2/ar/o2/co are available
1087            Mair_old=(q(iref,jref,1,igcm_co2)*mmol(igcm_co2)
1088     &               +q(iref,jref,1,igcm_n2)*mmol(igcm_n2)
1089     &               +q(iref,jref,1,igcm_ar)*mmol(igcm_ar)
1090     &               +q(iref,jref,1,igcm_o2)*mmol(igcm_o2)
1091     &               +q(iref,jref,1,igcm_co)*mmol(igcm_co))/Smmr
1092          end if
1093
1094          write(*,*)
1095     &      "Air molar mass (g/mol) at reference site= ",Mair_old
1096
1097        ! Ask for new volume mixing ratio at reference site
1098          Svmr =0.
1099          Sn =0.
1100          coefvmr(igcm_co2)=1.
1101
1102          do iq=1,nqtot
1103           coefvmr(iq) = 1.
1104           if ((iq.eq.igcm_n2).or.(iq.eq.igcm_ar)
1105     &     .or. (iq.eq.igcm_o2).or.(iq.eq.igcm_co)) then
1106
1107             vmr_old=q(iref,jref,1,iq)*Mair_old/mmol(iq) 
1108             write(*,*) "Previous vmr("//trim(tname(iq))//")= ", vmr_old
1109
1110              if (iq.eq.igcm_n2) then
1111                write(*,*) "New vmr(n2)? (MSL: 2.8e-02 at Ls~180,",
1112     &           " Trainer et al. 2019)"
1113              endif
1114              if (iq.eq.igcm_ar) then
1115                write(*,*) "New vmr(ar)? (MSL: 2.1e-02 at Ls~180)"
1116              endif
1117              if (iq.eq.igcm_o2) then
1118                write(*,*) "New vmr(o2)? (MSL: 1.7e-03 at Ls~180)"
1119              endif
1120              if (iq.eq.igcm_co) then
1121                write(*,*) "New vmr(co)? (ACS: 1e-03 at Ls~180)"
1122              endif
1123 302          read(*,*,iostat=ierr) vmr_new
1124              if(ierr.ne.0) goto 302
1125              write(*,*) "New vmr("//trim(tname(iq))//")= ",vmr_new
1126              write(*,*)
1127              coefvmr(iq) = vmr_new/vmr_old
1128              Svmr=Svmr+vmr_new
1129              Sn=Sn+vmr_new*mmol(iq)
1130           end if
1131          enddo ! of do iq=1,nqtot
1132
1133        ! Special case : implicit non-co2 gases JN 11/2019
1134          if ((igcm_n2.eq.0) .or. (igcm_ar.eq.0)) then
1135            write(*,*) "Warning : non-co2 gases are implicit"
1136            vmr_old=q(iref,jref,1,igcm_co2)*Mair_old/mmol(igcm_co2) 
1137            write(*,*) "Previous vmr(co2)=", vmr_old
1138            write(*,*) "New vmr(co2) ? (MSL: 0.947 at Ls~180)",
1139     &                  " Trainer et al. 2019)"
1140 666          read(*,*,iostat=ierr) vmr_new
1141              if(ierr.ne.0) goto 666
1142              coefvmr(igcm_co2) = vmr_new/vmr_old
1143              Svmr=Svmr+vmr_new
1144              Sn=vmr_new*mmol(igcm_co2) + (1-vmr_new)
1145     &         *33.87226017157708 ! Molar mass of non-co2 atm from MSL
1146          end if
1147      !  Estimation of new Air molar mass at reference site (assuming vmr_co2 = 1-Svmr)
1148          Mair_new = Sn + (1-Svmr)*mmol(igcm_co2)
1149      !  Estimation of new Air molar mass when non-co2 gases are implicit
1150          if ((igcm_n2.eq.0) .or. (igcm_ar.eq.0)) then
1151              Mair_new=vmr_new*mmol(igcm_co2) + (1-vmr_new)
1152     &         *33.87226017157708 ! Molar mass of non-co2 atm from MSL
1153           write(*,*)
1154     &     "We consider non-co2 gases vmr measured from Curiosity"
1155          end if
1156          write(*,*)
1157     &     "NEW Air molar mass (g/mol) at reference site= ",Mair_new
1158
1159        ! Compute mass mixing ratio changes 
1160          do iq=1,nqtot 
1161            if ((iq.eq.igcm_co2).or.(iq.eq.igcm_n2).or.(iq.eq.igcm_ar)
1162     &          .or. (iq.eq.igcm_o2).or.(iq.eq.igcm_co)) then
1163             write(*,*) "Everywhere mmr("//trim(tname(iq))//
1164     &        ") is multiplied by ",coefvmr(iq)*Mair_old/Mair_new
1165            end if
1166          end do
1167
1168        ! Recompute mass mixing ratios everywhere, and adjust mmr of most abundant species
1169        ! to keep sum of mmr constant.
1170          do l=1,llm
1171           do j=1,jjp1
1172            do i=1,iip1
1173              Smmr_old = 0.
1174              Smmr_new = 0.
1175              do iq=1,nqtot 
1176                if ((iq.eq.igcm_co2).or.(iq.eq.igcm_n2)
1177     &          .or.(iq.eq.igcm_ar)
1178     &          .or. (iq.eq.igcm_o2).or.(iq.eq.igcm_co)
1179     &          .or. (iq.eq.igcm_o) .or. (iq.eq. igcm_h2) ) then
1180                   Smmr_old = Smmr_old + q(i,j,l,iq) ! sum of old mmr
1181                   q(i,j,l,iq)=q(i,j,l,iq)*coefvmr(iq)*Mair_old/Mair_new
1182                   Smmr_new = Smmr_new + q(i,j,l,iq) ! sum of new mmr
1183                end if
1184              enddo
1185              !iloc = maxloc(q(i,j,l,:))
1186              iqmax=0 ; maxq=0
1187              do iq=1,nqtot
1188                if ((iq.eq.igcm_co2).or.(iq.eq.igcm_n2)
1189     &          .or.(iq.eq.igcm_ar)
1190     &          .or. (iq.eq.igcm_o2).or.(iq.eq.igcm_co)
1191     &          .or. (iq.eq.igcm_o) .or. (iq.eq. igcm_h2) ) then
1192                  if (q(i,j,l,iq).gt.maxq) then
1193                    maxq=q(i,j,l,iq)
1194                    iqmax=iq
1195                  endif
1196                endif
1197              enddo
1198              !iqmax = iloc(1)
1199              q(i,j,l,iqmax) = q(i,j,l,iqmax) + Smmr_old - Smmr_new
1200            enddo
1201           enddo
1202          enddo
1203
1204          write(*,*)
1205     &   "The most abundant species is modified everywhere to keep "//
1206     &   "sum of mmr constant"
1207          write(*,*) 'At reference site vmr(CO2)=',
1208     &        q(iref,jref,1,igcm_co2)*Mair_new/mmol(igcm_co2)
1209          write(*,*) "Compared to MSL observation: vmr(CO2)= 0.947 "//
1210     &   "at Ls=180"
1211
1212          Sn = q(iref,jref,1,igcm_co2)*Mair_new/mmol(igcm_co2)
1213     &       + q(iref,jref,1,igcm_n2)*Mair_new/mmol(igcm_n2)
1214     &       + q(iref,jref,1,igcm_ar)*Mair_new/mmol(igcm_ar)
1215     &       + q(iref,jref,1,igcm_o2)*Mair_new/mmol(igcm_o2)
1216     &       + q(iref,jref,1,igcm_co)*Mair_new/mmol(igcm_co)
1217
1218          write(*,*) 'Sum of volume mixing ratios = ', Sn
1219
1220c      wetstart : wet atmosphere with a north to south gradient
1221c      --------------------------------------------------------
1222       else if (trim(modif) .eq. 'wetstart') then
1223        ! check that there is indeed a water vapor tracer
1224        if (igcm_h2o_vap.eq.0) then
1225          write(*,*) "No water vapour tracer! Can't use this option"
1226          stop
1227        endif
1228          DO l=1,llm
1229            DO j=1,jjp1
1230              DO i=1,iip1-1
1231                q(i,j,l,igcm_h2o_vap)=150.e-6 * (rlatu(j)+pi/2.) / pi
1232              ENDDO
1233              ! We want to have the very same value at lon -180 and lon 180
1234              q(iip1,j,l,igcm_h2o_vap) = q(1,j,l,igcm_h2o_vap)
1235            ENDDO
1236          ENDDO
1237
1238         write(*,*) 'Water mass mixing ratio at north pole='
1239     *               ,q(1,1,1,igcm_h2o_vap)
1240         write(*,*) '---------------------------south pole='
1241     *               ,q(1,jjp1,1,igcm_h2o_vap)
1242
1243c      ini_h2osurf : reinitialize surface water ice
1244c      --------------------------------------------------
1245        else if (trim(modif) .eq. 'ini_h2osurf') then
1246          write(*,*)'max surface ice left?(e.g. 0.2 kg/m2=200microns)'
1247 207      read(*,*,iostat=ierr) val
1248          if(ierr.ne.0) goto 207
1249          write(*,*)'also set negative values of surf ice to 0'
1250           do ig=1,ngridmx
1251              qsurf(ig,igcm_h2o_ice)=min(val,qsurf(ig,igcm_h2o_ice))
1252              qsurf(ig,igcm_h2o_ice)=max(0.,qsurf(ig,igcm_h2o_ice))
1253           end do
1254
1255c      noglacier : remove tropical water ice (to initialize high res sim)
1256c      --------------------------------------------------
1257        else if (trim(modif) .eq. 'noglacier') then
1258           do ig=1,ngridmx
1259             j=(ig-2)/iim +2
1260              if(ig.eq.1) j=1
1261              write(*,*) 'OK: remove surface ice for |lat|<45'
1262              if (abs(rlatu(j)*180./pi).lt.45.) then
1263                   qsurf(ig,igcm_h2o_ice)=0.
1264              end if
1265           end do
1266
1267
1268c      watercapn : H20 ice on permanent northern cap
1269c      --------------------------------------------------
1270        else if (trim(modif) .eq. 'watercapn') then
1271           do ig=1,ngridmx
1272             j=(ig-2)/iim +2
1273              if(ig.eq.1) j=1
1274              if (rlatu(j)*180./pi.gt.80.) then
1275                   qsurf(ig,igcm_h2o_ice)=1.e5
1276                   write(*,*) 'ig=',ig,'    H2O ice mass (kg/m2)= ',
1277     &              qsurf(ig,igcm_h2o_ice)
1278                   write(*,*)'     ==> Ice mesh South boundary (deg)= ',
1279     &              rlatv(j)*180./pi
1280                end if
1281           enddo
1282
1283c      watercaps : H20 ice on permanent southern cap
1284c      -------------------------------------------------
1285        else if (trim(modif) .eq. 'watercaps') then
1286           do ig=1,ngridmx
1287               j=(ig-2)/iim +2
1288               if(ig.eq.1) j=1
1289               if (rlatu(j)*180./pi.lt.-80.) then
1290                   qsurf(ig,igcm_h2o_ice)=1.e5
1291                   write(*,*) 'ig=',ig,'   H2O ice mass (kg/m2)= ',
1292     &              qsurf(ig,igcm_h2o_ice)
1293                   write(*,*)'     ==> Ice mesh North boundary (deg)= ',
1294     &              rlatv(j-1)*180./pi
1295               end if
1296           enddo
1297
1298c       isotherm : Isothermal temperatures and no winds
1299c       ------------------------------------------------
1300        else if (trim(modif) .eq. 'isotherm') then
1301
1302          write(*,*)'Isothermal temperature of the atmosphere,
1303     &           surface and subsurface'
1304          write(*,*) 'Value of this temperature ? :'
1305 203      read(*,*,iostat=ierr) Tiso
1306          if(ierr.ne.0) goto 203
1307
1308          do ig=1, ngridmx
1309            tsurf(ig) = Tiso
1310          end do
1311          do l=2,nsoilmx
1312            do ig=1, ngridmx
1313              tsoil(ig,l) = Tiso
1314            end do
1315          end do
1316          flagiso=.true.
1317          call initial0(llm*ip1jmp1,ucov)
1318          call initial0(llm*ip1jm,vcov)
1319          call initial0(ngridmx*(llm+1),q2)
1320
1321c       co2ice=0 : remove CO2 polar ice caps'
1322c       ------------------------------------------------
1323        else if (trim(modif) .eq. 'co2ice=0') then
1324           do ig=1,ngridmx
1325              qsurf(ig,igcm_co2)=0
1326              emis(ig)=emis(ngridmx/2)
1327           end do
1328       
1329!       therm_ini_s: (re)-set soil thermal inertia to reference surface values
1330!       ----------------------------------------------------------------------
1331
1332        else if (trim(modif).eq.'therm_ini_s') then
1333!          write(*,*)"surfithfi(1):",surfithfi(1)
1334          do isoil=1,nsoilmx
1335            inertiedat(1:ngridmx,isoil)=surfithfi(1:ngridmx)
1336          enddo
1337          write(*,*)'OK: Soil thermal inertia has been reset to referenc
1338     &e surface values'
1339!         write(*,*)"inertiedat(1,1):",inertiedat(1,1)
1340          ithfi(:,:)=inertiedat(:,:)
1341         ! recast ithfi() onto ith()
1342         call gr_fi_dyn(nsoilmx,ngridmx,iip1,jjp1,ithfi,ith)
1343! Check:
1344!         do i=1,iip1
1345!           do j=1,jjp1
1346!             do isoil=1,nsoilmx
1347!               write(77,*) i,j,isoil,"  ",ith(i,j,isoil)
1348!             enddo
1349!           enddo
1350!        enddo
1351
1352!       subsoilice_n: Put deep ice layer in northern hemisphere soil
1353!       ------------------------------------------------------------
1354
1355        else if (trim(modif).eq.'subsoilice_n') then
1356
1357         write(*,*)'From which latitude (in deg.), up to the north pole,
1358     &should we put subterranean ice?'
1359         ierr=1
1360         do while (ierr.ne.0)
1361          read(*,*,iostat=ierr) val
1362          if (ierr.eq.0) then ! got a value
1363            ! do a sanity check
1364            if((val.lt.0.).or.(val.gt.90)) then
1365              write(*,*)'Latitude should be between 0 and 90 deg. !!!'
1366              ierr=1
1367            else ! find corresponding jref (nearest latitude)
1368              ! note: rlatu(:) contains decreasing values of latitude
1369              !       starting from PI/2 to -PI/2
1370              do j=1,jjp1
1371                if ((rlatu(j)*180./pi.ge.val).and.
1372     &              (rlatu(j+1)*180./pi.le.val)) then
1373                  ! find which grid point is nearest to val:
1374                  if (abs(rlatu(j)*180./pi-val).le.
1375     &                abs((rlatu(j+1)*180./pi-val))) then
1376                   jref=j
1377                  else
1378                   jref=j+1
1379                  endif
1380                 
1381                 write(*,*)'Will use nearest grid latitude which is:',
1382     &                     rlatu(jref)*180./pi
1383                endif
1384              enddo ! of do j=1,jjp1
1385            endif ! of if((val.lt.0.).or.(val.gt.90))
1386          endif !of if (ierr.eq.0)
1387         enddo ! of do while
1388
1389         ! Build layers() (as in soil_settings.F)
1390         val2=sqrt(mlayer(0)*mlayer(1))
1391         val3=mlayer(1)/mlayer(0)
1392         do isoil=1,nsoilmx
1393           layer(isoil)=val2*(val3**(isoil-1))
1394         enddo
1395
1396         write(*,*)'At which depth (in m.) does the ice layer begin?'
1397         write(*,*)'(currently, the deepest soil layer extends down to:'
1398     &              ,layer(nsoilmx),')'
1399         ierr=1
1400         do while (ierr.ne.0)
1401          read(*,*,iostat=ierr) val2
1402!         write(*,*)'val2:',val2,'ierr=',ierr
1403          if (ierr.eq.0) then ! got a value, but do a sanity check
1404            if(val2.gt.layer(nsoilmx)) then
1405              write(*,*)'Depth should be less than ',layer(nsoilmx)
1406              ierr=1
1407            endif
1408            if(val2.lt.layer(1)) then
1409              write(*,*)'Depth should be more than ',layer(1)
1410              ierr=1
1411            endif
1412          endif
1413         enddo ! of do while
1414         
1415         ! find the reference index iref the depth corresponds to
1416!        if (val2.lt.layer(1)) then
1417!         iref=1
1418!        else
1419          do isoil=1,nsoilmx-1
1420           if((val2.gt.layer(isoil)).and.(val2.lt.layer(isoil+1)))
1421     &       then
1422             iref=isoil
1423             exit
1424           endif
1425          enddo
1426!        endif
1427         
1428!        write(*,*)'iref:',iref,'  jref:',jref
1429!        write(*,*)'layer',layer
1430!        write(*,*)'mlayer',mlayer
1431         
1432         ! thermal inertia of the ice:
1433         ierr=1
1434         do while (ierr.ne.0)
1435          write(*,*)'What is the value of subterranean ice thermal inert
1436     &ia? (e.g.: 2000)'
1437          read(*,*,iostat=ierr)iceith
1438         enddo ! of do while
1439         
1440         ! recast ithfi() onto ith()
1441         call gr_fi_dyn(nsoilmx,ngridmx,iip1,jjp1,ithfi,ith)
1442         
1443         do j=1,jref
1444!           write(*,*)'j:',j,'rlatu(j)*180./pi:',rlatu(j)*180./pi
1445            do i=1,iip1 ! loop on longitudes
1446             ! Build "equivalent" thermal inertia for the mixed layer
1447             ith(i,j,iref+1)=sqrt((layer(iref+1)-layer(iref))/
1448     &                     (((val2-layer(iref))/(ith(i,j,iref)**2))+
1449     &                      ((layer(iref+1)-val2)/(iceith)**2)))
1450             ! Set thermal inertia of lower layers
1451             do isoil=iref+2,nsoilmx
1452              ith(i,j,isoil)=iceith ! ice
1453             enddo
1454            enddo ! of do i=1,iip1
1455         enddo ! of do j=1,jjp1
1456         
1457
1458         CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi)
1459
1460!         do i=1,nsoilmx
1461!         write(*,*)'i:',i,'ithfi(1,i):',ithfi(1,i)
1462!        enddo
1463
1464       
1465!       subsoilice_s: Put deep ice layer in southern hemisphere soil
1466!       ------------------------------------------------------------
1467
1468        else if (trim(modif).eq.'subsoilice_s') then
1469
1470         write(*,*)'From which latitude (in deg.), down to the south pol
1471     &e, should we put subterranean ice?'
1472         ierr=1
1473         do while (ierr.ne.0)
1474          read(*,*,iostat=ierr) val
1475          if (ierr.eq.0) then ! got a value
1476            ! do a sanity check
1477            if((val.gt.0.).or.(val.lt.-90)) then
1478              write(*,*)'Latitude should be between 0 and -90 deg. !!!'
1479              ierr=1
1480            else ! find corresponding jref (nearest latitude)
1481              ! note: rlatu(:) contains decreasing values of latitude
1482              !       starting from PI/2 to -PI/2
1483              do j=1,jjp1
1484                if ((rlatu(j)*180./pi.ge.val).and.
1485     &              (rlatu(j+1)*180./pi.le.val)) then
1486                  ! find which grid point is nearest to val:
1487                  if (abs(rlatu(j)*180./pi-val).le.
1488     &                abs((rlatu(j+1)*180./pi-val))) then
1489                   jref=j
1490                  else
1491                   jref=j+1
1492                  endif
1493                 
1494                 write(*,*)'Will use nearest grid latitude which is:',
1495     &                     rlatu(jref)*180./pi
1496                endif
1497              enddo ! of do j=1,jjp1
1498            endif ! of if((val.lt.0.).or.(val.gt.90))
1499          endif !of if (ierr.eq.0)
1500         enddo ! of do while
1501
1502         ! Build layers() (as in soil_settings.F)
1503         val2=sqrt(mlayer(0)*mlayer(1))
1504         val3=mlayer(1)/mlayer(0)
1505         do isoil=1,nsoilmx
1506           layer(isoil)=val2*(val3**(isoil-1))
1507         enddo
1508
1509         write(*,*)'At which depth (in m.) does the ice layer begin?'
1510         write(*,*)'(currently, the deepest soil layer extends down to:'
1511     &              ,layer(nsoilmx),')'
1512         ierr=1
1513         do while (ierr.ne.0)
1514          read(*,*,iostat=ierr) val2
1515!         write(*,*)'val2:',val2,'ierr=',ierr
1516          if (ierr.eq.0) then ! got a value, but do a sanity check
1517            if(val2.gt.layer(nsoilmx)) then
1518              write(*,*)'Depth should be less than ',layer(nsoilmx)
1519              ierr=1
1520            endif
1521            if(val2.lt.layer(1)) then
1522              write(*,*)'Depth should be more than ',layer(1)
1523              ierr=1
1524            endif
1525          endif
1526         enddo ! of do while
1527         
1528         ! find the reference index iref the depth corresponds to
1529          do isoil=1,nsoilmx-1
1530           if((val2.gt.layer(isoil)).and.(val2.lt.layer(isoil+1)))
1531     &       then
1532             iref=isoil
1533             exit
1534           endif
1535          enddo
1536         
1537!        write(*,*)'iref:',iref,'  jref:',jref
1538         
1539         ! thermal inertia of the ice:
1540         ierr=1
1541         do while (ierr.ne.0)
1542          write(*,*)'What is the value of subterranean ice thermal inert
1543     &ia? (e.g.: 2000)'
1544          read(*,*,iostat=ierr)iceith
1545         enddo ! of do while
1546         
1547         ! recast ithfi() onto ith()
1548         call gr_fi_dyn(nsoilmx,ngridmx,iip1,jjp1,ithfi,ith)
1549         
1550         do j=jref,jjp1
1551!           write(*,*)'j:',j,'rlatu(j)*180./pi:',rlatu(j)*180./pi
1552            do i=1,iip1 ! loop on longitudes
1553             ! Build "equivalent" thermal inertia for the mixed layer
1554             ith(i,j,iref+1)=sqrt((layer(iref+1)-layer(iref))/
1555     &                     (((val2-layer(iref))/(ith(i,j,iref)**2))+
1556     &                      ((layer(iref+1)-val2)/(iceith)**2)))
1557             ! Set thermal inertia of lower layers
1558             do isoil=iref+2,nsoilmx
1559              ith(i,j,isoil)=iceith ! ice
1560             enddo
1561            enddo ! of do i=1,iip1
1562         enddo ! of do j=jref,jjp1
1563         
1564
1565         CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi)
1566
1567c       'mons_ice' : use MONS data to build subsurface ice table
1568c       --------------------------------------------------------
1569        else if (trim(modif).eq.'mons_ice') then
1570       
1571       ! 1. Load MONS data
1572        call load_MONS_data(MONS_Hdn,MONS_d21)
1573       
1574        ! 2. Get parameters from user
1575        ierr=1
1576        do while (ierr.ne.0)
1577          write(*,*) "Coefficient to apply to MONS 'depth' in Northern",
1578     &               " Hemisphere?"
1579          write(*,*) " (should be somewhere between 3.2e-4 and 1.3e-3)"
1580          read(*,*,iostat=ierr) MONS_coeffN
1581        enddo
1582        ierr=1
1583        do while (ierr.ne.0)
1584          write(*,*) "Coefficient to apply to MONS 'depth' in Southern",
1585     &               " Hemisphere?"
1586          write(*,*) " (should be somewhere between 3.2e-4 and 1.3e-3)"
1587          read(*,*,iostat=ierr) MONS_coeffS
1588        enddo
1589        ierr=1
1590        do while (ierr.ne.0)
1591          write(*,*) "Value of subterranean ice thermal inertia ",
1592     &               " in Northern hemisphere?"
1593          write(*,*) " (e.g.: 2000, or perhaps 2290)"
1594!          read(*,*,iostat=ierr) iceith
1595          read(*,*,iostat=ierr) iceithN
1596        enddo
1597        ierr=1
1598        do while (ierr.ne.0)
1599          write(*,*) "Value of subterranean ice thermal inertia ",
1600     &               " in Southern hemisphere?"
1601          write(*,*) " (e.g.: 2000, or perhaps 2290)"
1602!          read(*,*,iostat=ierr) iceith
1603          read(*,*,iostat=ierr) iceithS
1604        enddo
1605       
1606        ! 3. Build subterranean thermal inertia
1607       
1608        ! initialise subsurface inertia with reference surface values
1609        do isoil=1,nsoilmx
1610          ithfi(1:ngridmx,isoil)=surfithfi(1:ngridmx)
1611        enddo
1612        ! recast ithfi() onto ith()
1613        call gr_fi_dyn(nsoilmx,ngridmx,iip1,jjp1,ithfi,ith)
1614       
1615        do i=1,iip1 ! loop on longitudes
1616          do j=1,jjp1 ! loop on latitudes
1617            ! set MONS_coeff
1618            if (rlatu(j).ge.0) then ! northern hemisphere
1619              ! N.B: rlatu(:) contains decreasing values of latitude
1620              !       starting from PI/2 to -PI/2
1621              MONS_coeff=MONS_coeffN
1622              iceith=iceithN
1623            else ! southern hemisphere
1624              MONS_coeff=MONS_coeffS
1625              iceith=iceithS
1626            endif
1627            ! check if we should put subterranean ice
1628            if (MONS_Hdn(i,j).ge.14.0) then ! no ice if Hdn<14%
1629              ! compute depth at which ice lies:
1630              val=MONS_d21(i,j)*MONS_coeff
1631              ! compute val2= the diurnal skin depth of surface inertia
1632              ! assuming a volumetric heat cap. of C=1.e6 J.m-3.K-1
1633              val2=ith(i,j,1)*1.e-6*sqrt(88775./3.14159)
1634              if (val.lt.val2) then
1635                ! ice must be below the (surface inertia) diurnal skin depth
1636                val=val2
1637              endif
1638              if (val.lt.layer(nsoilmx)) then ! subterranean ice
1639                ! find the reference index iref that depth corresponds to
1640                iref=0
1641                do isoil=1,nsoilmx-1
1642                 if ((val.ge.layer(isoil)).and.(val.lt.layer(isoil+1)))
1643     &             then
1644                   iref=isoil
1645                   exit
1646                 endif
1647                enddo
1648                ! Build "equivalent" thermal inertia for the mixed layer
1649                ith(i,j,iref+1)=sqrt((layer(iref+1)-layer(iref))/
1650     &                     (((val-layer(iref))/(ith(i,j,iref+1)**2))+
1651     &                      ((layer(iref+1)-val)/(iceith)**2)))
1652                ! Set thermal inertia of lower layers
1653                do isoil=iref+2,nsoilmx
1654                  ith(i,j,isoil)=iceith
1655                enddo
1656              endif ! of if (val.lt.layer(nsoilmx))
1657            endif ! of if (MONS_Hdn(i,j).lt.14.0)
1658          enddo ! do j=1,jjp1
1659        enddo ! do i=1,iip1
1660       
1661! Check:
1662!         do i=1,iip1
1663!           do j=1,jjp1
1664!             do isoil=1,nsoilmx
1665!               write(77,*) i,j,isoil,"  ",ith(i,j,isoil)
1666!             enddo
1667!           enddo
1668!        enddo
1669
1670        ! recast ith() into ithfi()
1671        CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi)
1672       
1673        else
1674          write(*,*) '       Unknown (misspelled?) option!!!'
1675        end if ! of if (trim(modif) .eq. '...') elseif ...
1676       
1677       enddo ! of do ! infinite loop on liste of changes
1678
1679 999  continue
1680
1681 
1682c=======================================================================
1683c   Correct pressure on the new grid (menu 0)
1684c=======================================================================
1685
1686      if (choix_1.eq.0) then
1687        r = 1000.*8.31/mugaz
1688
1689        do j=1,jjp1
1690          do i=1,iip1
1691             ps(i,j) = ps(i,j) *
1692     .            exp((phisold_newgrid(i,j)-phis(i,j)) /
1693     .                                  (t(i,j,1) * r))
1694          end do
1695        end do
1696 
1697c periodicity of surface ps in longitude
1698        do j=1,jjp1
1699          ps(1,j) = ps(iip1,j)
1700        end do
1701      end if
1702
1703c=======================================================================
1704c=======================================================================
1705
1706c=======================================================================
1707c    Initialisation de la physique / ecriture de newstartfi :
1708c=======================================================================
1709
1710
1711      CALL inifilr
1712      CALL pression(ip1jmp1, ap, bp, ps, p3d)
1713
1714c-----------------------------------------------------------------------
1715c   Initialisation  pks:
1716c-----------------------------------------------------------------------
1717
1718      CALL exner_hyb(ip1jmp1, ps, p3d, pks, pk, pkf)
1719! Calcul de la temperature potentielle teta
1720
1721      if (flagiso) then
1722          DO l=1,llm
1723             DO j=1,jjp1
1724                DO i=1,iim
1725                   teta(i,j,l) = Tiso * cpp/pk(i,j,l)
1726                ENDDO
1727                teta (iip1,j,l)= teta (1,j,l)
1728             ENDDO
1729          ENDDO
1730      else if (choix_1.eq.0) then
1731         DO l=1,llm
1732            DO j=1,jjp1
1733               DO i=1,iim
1734                  teta(i,j,l) = t(i,j,l) * cpp/pk(i,j,l)
1735               ENDDO
1736               teta (iip1,j,l)= teta (1,j,l)
1737            ENDDO
1738         ENDDO
1739      endif
1740
1741C Calcul intermediaire
1742c
1743      if (choix_1.eq.0) then
1744         CALL massdair( p3d, masse  )
1745c
1746         print *,' ALPHAX ',alphax
1747
1748         DO  l = 1, llm
1749           DO  i    = 1, iim
1750             xppn(i) = aire( i, 1   ) * masse(  i     ,  1   , l )
1751             xpps(i) = aire( i,jjp1 ) * masse(  i     , jjp1 , l )
1752           ENDDO
1753             xpn      = SUM(xppn)/apoln
1754             xps      = SUM(xpps)/apols
1755           DO i   = 1, iip1
1756             masse(   i   ,   1     ,  l )   = xpn
1757             masse(   i   ,   jjp1  ,  l )   = xps
1758           ENDDO
1759         ENDDO
1760      endif
1761      phis(iip1,:) = phis(1,:)
1762
1763c      CALL inidissip ( lstardis, nitergdiv, nitergrot, niterh,
1764c     *                tetagdiv, tetagrot , tetatemp  )
1765      itau=0
1766      if (choix_1.eq.0) then
1767         day_ini=int(date)
1768         hour_ini=date-int(date)
1769      endif
1770c
1771      CALL geopot  ( ip1jmp1, teta  , pk , pks,  phis  , phi   )
1772
1773      CALL caldyn0( itau,ucov,vcov,teta,ps,masse,pk,phis ,
1774     *                phi,w, pbaru,pbarv,day_ini+time )
1775c     CALL caldyn
1776c    $  ( itau,ucov,vcov,teta,ps,masse,pk,pkf,phis ,
1777c    $    phi,conser,du,dv,dteta,dp,w, pbaru,pbarv, day_ini )
1778
1779      CALL dynredem0("restart.nc",day_ini,phis)
1780      CALL dynredem1("restart.nc",hour_ini,vcov,ucov,teta,q,
1781     .               masse,ps)
1782C
1783C Ecriture etat initial physique
1784C
1785      call physdem0("restartfi.nc",longitude,latitude,
1786     &              nsoilmx,ngridmx,llm,
1787     &              nqtot,dtphys,real(day_ini),0.0,cell_area,
1788     &              albfi,ithfi)
1789      call physdem1("restartfi.nc",nsoilmx,ngridmx,llm,nqtot,
1790     &              dtphys,hour_ini,
1791     &              tsurf,tsoil,albedo,emis,q2,qsurf,tauscaling,
1792     &              totcloudfrac,wstar,watercap)
1793
1794c=======================================================================
1795c       Formats
1796c=======================================================================
1797
1798   1  FORMAT(//10x,'la valeur de im =',i4,2x,'lue sur le fichier de dema
1799     *rrage est differente de la valeur parametree iim =',i4//)
1800   2  FORMAT(//10x,'la valeur de jm =',i4,2x,'lue sur le fichier de dema
1801     *rrage est differente de la valeur parametree jjm =',i4//)
1802   3  FORMAT(//10x,'la valeur de lllm =',i4,2x,'lue sur le fichier demar
1803     *rage est differente de la valeur parametree llm =',i4//)
1804
1805      write(*,*) "newstart: All is well that ends well."
1806
1807      end
1808
1809!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1810      subroutine load_MONS_data(MONS_Hdn,MONS_d21)
1811     
1812      use datafile_mod, only:datadir
1813     
1814      implicit none
1815      ! routine to load Benedicte Diez MONS dataset, fill in date in southern
1816      ! polar region, and interpolate the result onto the GCM grid
1817      include"dimensions.h"
1818      include"paramet.h"
1819      include"comgeom.h"
1820      ! arguments:
1821      real,intent(out) :: MONS_Hdn(iip1,jjp1) ! Hdn: %WEH=Mass fraction of H2O
1822      real,intent(out) :: MONS_d21(iip1,jjp1) ! ice table "depth" (in kg/m2)
1823      ! N.B MONS datasets should be of dimension (iip1,jjp1)
1824      ! local variables:
1825      character(len=88) :: filename="results_MONS_lat_lon_H_depth.txt"
1826      character(len=88) :: txt ! to store some text
1827      integer :: ierr,i,j
1828      integer,parameter :: nblon=180 ! number of longitudes of MONS datasets
1829      integer,parameter :: nblat=90 ! number of latitudes of MONS datasets
1830      real :: pi
1831      real :: longitudes(nblon) ! MONS dataset longitudes
1832      real :: latitudes(nblat)  ! MONS dataset latitudes
1833      ! MONS dataset: mass fraction of H2O where H is assumed to be in H2O
1834      real :: Hdn(nblon,nblat)
1835      real :: d21(nblon,nblat)! MONS dataset "depth" (g/cm2)
1836
1837      ! Extended MONS dataset (for interp_horiz)
1838      real :: Hdnx(nblon+1,nblat)
1839      real :: d21x(nblon+1,nblat)
1840      real :: lon_bound(nblon+1) ! longitude boundaries
1841      real :: lat_bound(nblat-1) ! latitude boundaries
1842
1843      ! 1. Initializations:
1844
1845      write(*,*) "Loading MONS data"
1846
1847      ! Open MONS datafile:
1848      open(42,file=trim(datadir)//"/"//trim(filename),
1849     &     status="old",iostat=ierr)
1850      if (ierr/=0) then
1851        write(*,*) "Error in load_MONS_data:"
1852        write(*,*) "Failed opening file ",
1853     &             trim(datadir)//"/"//trim(filename)
1854        write(*,*)'1) You can change this directory address in ',
1855     &            'callfis.def with'
1856        write(*,*)'   datadir=/path/to/datafiles'
1857        write(*,*)'2) If necessary ',trim(filename),
1858     &                 ' (and other datafiles)'
1859        write(*,*)'   can be obtained online at:'
1860        write(*,*)'http://www.lmd.jussieu.fr/~lmdz/planets/mars/datadir'
1861        CALL ABORT
1862      else ! skip first line of file (dummy read)
1863         read(42,*) txt
1864      endif
1865
1866      pi=2.*asin(1.)
1867     
1868      !2. Load MONS data (on MONS grid)
1869      do j=1,nblat
1870        do i=1,nblon
1871        ! swap latitude index so latitudes go from north pole to south pole:
1872          read(42,*) latitudes(nblat-j+1),longitudes(i),
1873     &               Hdn(i,nblat-j+1),d21(i,nblat-j+1)
1874        ! multiply d21 by 10 to convert from g/cm2 to kg/m2
1875          d21(i,nblat-j+1)=d21(i,nblat-j+1)*10.0
1876        enddo
1877      enddo
1878      close(42)
1879     
1880      ! there is unfortunately no d21 data for latitudes -77 to -90
1881      ! so we build some by linear interpolation between values at -75
1882      ! and assuming d21=0 at the pole
1883      do j=84,90 ! latitudes(84)=-77 ; latitudes(83)=-75
1884        do i=1,nblon
1885          d21(i,j)=d21(i,83)*((latitudes(j)+90)/15.0)
1886        enddo
1887      enddo
1888
1889      ! 3. Build extended MONS dataset & boundaries (for interp_horiz)
1890      ! longitude boundaries (in radians):
1891      do i=1,nblon
1892        ! NB: MONS data is every 2 degrees in longitude
1893        lon_bound(i)=(longitudes(i)+1.0)*pi/180.0
1894      enddo
1895      ! extra 'modulo' value
1896      lon_bound(nblon+1)=lon_bound(1)+2.0*pi
1897     
1898      ! latitude boundaries (in radians):
1899      do j=1,nblat-1
1900        ! NB: Mons data is every 2 degrees in latitude
1901        lat_bound(j)=(latitudes(j)-1.0)*pi/180.0
1902      enddo
1903     
1904      ! MONS datasets:
1905      do j=1,nblat
1906        Hdnx(1:nblon,j)=Hdn(1:nblon,j)
1907        Hdnx(nblon+1,j)=Hdnx(1,j)
1908        d21x(1:nblon,j)=d21(1:nblon,j)
1909        d21x(nblon+1,j)=d21x(1,j)
1910      enddo
1911     
1912      ! Interpolate onto GCM grid
1913      call interp_horiz(Hdnx,MONS_Hdn,nblon,nblat-1,iim,jjm,1,
1914     &                  lon_bound,lat_bound,rlonu,rlatv)
1915      call interp_horiz(d21x,MONS_d21,nblon,nblat-1,iim,jjm,1,
1916     &                  lon_bound,lat_bound,rlonu,rlatv)
1917     
1918      end subroutine
Note: See TracBrowser for help on using the repository browser.