C====================================================================== PROGRAM newstart c======================================================================= c c c Auteur: Christophe Hourdin/Francois Forget/Yann Wanherdrick c ------ c Derniere modif : 12/03 c c c Objet: Create or modify the initial state for the LMD Mars GCM c ----- (fichiers NetCDF start et startfi) c c c======================================================================= implicit none #include "dimensions.h" #include "dimphys.h" #include "surfdat.h" #include "comsoil.h" #include "planete.h" #include "paramet.h" #include "comconst.h" #include "comvert.h" #include "comgeom2.h" #include "control.h" #include "logic.h" #include "description.h" #include "ener.h" #include "temps.h" #include "lmdstd.h" #include "comdissnew.h" #include "clesph0.h" #include "serre.h" #include "netcdf.inc" #include "advtrac.h" #include "tracer.h" c======================================================================= c Declarations c======================================================================= c Variables dimension du fichier "start_archive" c------------------------------------ CHARACTER relief*3 c Variables pour les lectures NetCDF des fichiers "start_archive" c-------------------------------------------------- INTEGER nid_dyn, nid_fi,nid,nvarid INTEGER length parameter (length = 100) INTEGER tab0 INTEGER NB_ETATMAX parameter (NB_ETATMAX = 100) REAL date REAL p_rad,p_omeg,p_g,p_cpp,p_mugaz,p_daysec c Variable histoire c------------------ REAL vcov(iip1,jjm,llm),ucov(iip1,jjp1,llm) ! vents covariants REAL phis(iip1,jjp1) REAL q(iip1,jjp1,llm,nqmx) ! champs advectes c autre variables dynamique nouvelle grille c------------------------------------------ REAL pks(iip1,jjp1) REAL w(iip1,jjp1,llm+1) REAL pbaru(ip1jmp1,llm),pbarv(ip1jm,llm) ! REAL dv(ip1jm,llm),du(ip1jmp1,llm) ! REAL dh(ip1jmp1,llm),dp(ip1jmp1) REAL phi(iip1,jjp1,llm) integer klatdat,klongdat PARAMETER (klatdat=180,klongdat=360) c Physique sur grille scalaire c---------------------------- real zmeaS(iip1,jjp1),zstdS(iip1,jjp1) real zsigS(iip1,jjp1),zgamS(iip1,jjp1),ztheS(iip1,jjp1) c variable physique c------------------ REAL tsurf(ngridmx) ! surface temperature REAL tsoil(ngridmx,nsoilmx) ! soil temperature ! REAL co2ice(ngridmx) ! CO2 ice layer REAL emis(ngridmx) ! surface emissivity real emisread ! added by RW REAL qsurf(ngridmx,nqmx) REAL q2(ngridmx,nlayermx+1) ! REAL rnaturfi(ngridmx) real alb(iip1,jjp1),albfi(ngridmx) ! albedos real ith(iip1,jjp1,nsoilmx),ithfi(ngridmx,nsoilmx) ! thermal inertia (3D) real surfith(iip1,jjp1),surfithfi(ngridmx) ! surface thermal inertia (2D) REAL latfi(ngridmx),lonfi(ngridmx),airefi(ngridmx) INTEGER i,j,l,isoil,ig,idum real mugaz ! molar mass of the atmosphere integer ierr c Variables on the new grid along scalar points c------------------------------------------------------ ! REAL p(iip1,jjp1) REAL t(iip1,jjp1,llm) REAL tset(iip1,jjp1,llm) real phisold_newgrid(iip1,jjp1) REAL :: teta(iip1, jjp1, llm) REAL :: pk(iip1,jjp1,llm) REAL :: pkf(iip1,jjp1,llm) REAL :: ps(iip1, jjp1) REAL :: masse(iip1,jjp1,llm) REAL :: xpn,xps,xppn(iim),xpps(iim) REAL :: p3d(iip1, jjp1, llm+1) REAL :: beta(iip1,jjp1,llm) ! REAL dteta(ip1jmp1,llm) c Variable de l'ancienne grille c------------------------------ real time real tab_cntrl(100) real tab_cntrl_bis(100) c variables diverses c------------------- real choix_1,pp character*80 fichnom integer Lmodif,iq,thermo character modif*20 real z_reel(iip1,jjp1) real tsud,albsud,alb_bb,ith_bb,Tiso,Tabove real ptoto,pcap,patm,airetot,ptotn,patmn,psea ! real ssum character*1 yes logical :: flagtset=.false. , flagps0=.false. real val, val2, val3 ! to store temporary variables real :: iceith=2000 ! thermal inertia of subterranean ice integer iref,jref INTEGER :: itau INTEGER :: nq,numvanle character(len=20) :: txt ! to store some text integer :: count ! MONS data: real :: MONS_Hdn(iip1,jjp1) ! Hdn: %WEH=Mass fraction of H2O real :: MONS_d21(iip1,jjp1) ! ice table "depth" (in kg/m2) ! coefficient to apply to convert d21 to 'true' depth (m) real :: MONS_coeff real :: MONS_coeffS ! coeff for southern hemisphere real :: MONS_coeffN ! coeff for northern hemisphere ! real,parameter :: icedepthmin=1.e-3 ! Ice begins at most at that depth ! added by RW for test real pmean, phi0 ! added by BC for equilibrium temperature startup real teque ! added by BC for cloud fraction setup REAL hice(ngridmx),cloudfrac(ngridmx,nlayermx) REAL totalfrac(ngridmx) ! added by RW for nuketharsis real fact1 real fact2 c sortie visu pour les champs dynamiques c--------------------------------------- ! INTEGER :: visuid ! real :: time_step,t_ops,t_wrt ! CHARACTER*80 :: visu_file ! cpp = 744.499 ! for Mars, instead of 1004.70885 (Earth) ! preff = 610. ! for Mars, instead of 101325. (Earth) ! pa = 20 ! for Mars, instead of 500 (Earth) cpp = 0. preff = 0. pa = 0. ! to ensure disaster rather than minor error if we don`t ! make deliberate choice of these values elsewhere. c======================================================================= c Choice of the start file(s) to use c======================================================================= write(*,*) 'From which kind of files do you want to create new', . 'start and startfi files' write(*,*) ' 0 - from a file start_archive' write(*,*) ' 1 - from files start and startfi' c----------------------------------------------------------------------- c Open file(s) to modify (start or start_archive) c----------------------------------------------------------------------- DO read(*,*,iostat=ierr) choix_1 if ((choix_1 /= 0).OR.(choix_1 /=1)) EXIT ENDDO c Open start_archive c ~~~~~~~~~~~~~~~~~~~~~~~~~~ if (choix_1.eq.0) then write(*,*) 'Creating start files from:' write(*,*) './start_archive.nc' write(*,*) fichnom = 'start_archive.nc' ierr = NF_OPEN (fichnom, NF_NOWRITE,nid) IF (ierr.NE.NF_NOERR) THEN write(6,*)' Problem opening file:',fichnom write(6,*)' ierr = ', ierr CALL ABORT ENDIF tab0 = 50 Lmodif = 1 c OR open start and startfi files c ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ else write(*,*) 'Creating start files from:' write(*,*) './start.nc and ./startfi.nc' write(*,*) fichnom = 'start.nc' ierr = NF_OPEN (fichnom, NF_NOWRITE,nid_dyn) IF (ierr.NE.NF_NOERR) THEN write(6,*)' Problem opening file:',fichnom write(6,*)' ierr = ', ierr CALL ABORT ENDIF fichnom = 'startfi.nc' ierr = NF_OPEN (fichnom, NF_NOWRITE,nid_fi) IF (ierr.NE.NF_NOERR) THEN write(6,*)' Problem opening file:',fichnom write(6,*)' ierr = ', ierr CALL ABORT ENDIF tab0 = 0 Lmodif = 0 endif c----------------------------------------------------------------------- c Lecture du tableau des parametres du run (pour la dynamique) c----------------------------------------------------------------------- if (choix_1.eq.0) then write(*,*) 'reading tab_cntrl START_ARCHIVE' c ierr = NF_INQ_VARID (nid, "controle", nvarid) #ifdef NC_DOUBLE ierr = NF_GET_VAR_DOUBLE(nid, nvarid, tab_cntrl) #else ierr = NF_GET_VAR_REAL(nid, nvarid, tab_cntrl) #endif c else if (choix_1.eq.1) then write(*,*) 'reading tab_cntrl START' c ierr = NF_INQ_VARID (nid_dyn, "controle", nvarid) #ifdef NC_DOUBLE ierr = NF_GET_VAR_DOUBLE(nid_dyn, nvarid, tab_cntrl) #else ierr = NF_GET_VAR_REAL(nid_dyn, nvarid, tab_cntrl) #endif c write(*,*) 'reading tab_cntrl STARTFI' c ierr = NF_INQ_VARID (nid_fi, "controle", nvarid) #ifdef NC_DOUBLE ierr = NF_GET_VAR_DOUBLE(nid_fi, nvarid, tab_cntrl_bis) #else ierr = NF_GET_VAR_REAL(nid_fi, nvarid, tab_cntrl_bis) #endif c do i=1,50 tab_cntrl(i+50)=tab_cntrl_bis(i) enddo write(*,*) 'printing tab_cntrl', tab_cntrl do i=1,100 write(*,*) i,tab_cntrl(i) enddo endif c----------------------------------------------------------------------- c Initialisation des constantes dynamique c----------------------------------------------------------------------- kappa = tab_cntrl(9) etot0 = tab_cntrl(12) ptot0 = tab_cntrl(13) ztot0 = tab_cntrl(14) stot0 = tab_cntrl(15) ang0 = tab_cntrl(16) write(*,*) "Newstart: kappa,etot0,ptot0,ztot0,stot0,ang0" write(*,*) kappa,etot0,ptot0,ztot0,stot0,ang0 ! for vertical coordinate preff=tab_cntrl(18) ! reference surface pressure pa=tab_cntrl(17) ! reference pressure at which coord is purely pressure !NB: in start_archive files tab_cntrl(17)=tab_cntrl(18)=0 if (preff.eq.0) then preff=610 pa=20 endif write(*,*) "Newstart: preff=",preff," pa=",pa yes=' ' do while ((yes.ne.'y').and.(yes.ne.'n')) write(*,*) "Change the values of preff and pa? (y/n)" read(*,fmt='(a)') yes end do if (yes.eq.'y') then write(*,*)"New value of reference surface pressure preff?" write(*,*)" (for Mars, typically preff=610)" read(*,*) preff write(*,*)"New value of reference pressure pa for purely" write(*,*)"pressure levels (for hybrid coordinates)?" write(*,*)" (for Mars, typically pa=20)" read(*,*) pa endif c----------------------------------------------------------------------- c Lecture du tab_cntrl et initialisation des constantes physiques c - pour start: Lmodif = 0 => pas de modifications possibles c (modif dans le tabfi de readfi + loin) c - pour start_archive: Lmodif = 1 => modifications possibles c----------------------------------------------------------------------- if (choix_1.eq.0) then call tabfi (nid,Lmodif,tab0,day_ini,lllm,p_rad, . p_omeg,p_g,p_cpp,p_mugaz,p_daysec,time) else if (choix_1.eq.1) then call tabfi (nid_fi,Lmodif,tab0,day_ini,lllm,p_rad, . p_omeg,p_g,p_cpp,p_mugaz,p_daysec,time) endif rad = p_rad omeg = p_omeg g = p_g cpp = p_cpp mugaz = p_mugaz daysec = p_daysec kappa = 8.314*1000./(p_mugaz*p_cpp) ! added by RDW c======================================================================= c INITIALISATIONS DIVERSES c======================================================================= ! Load tracer names: call iniadvtrac(nq,numvanle) ! tnom(:) now contains tracer names ! Initialize global tracer indexes (stored in tracer.h) call initracer() ! Load parameters from run.def file CALL defrun_new( 99, .TRUE. ) CALL iniconst CALL inigeom idum=-1 idum=0 c Initialisation coordonnees /aires c ------------------------------- ! Note: rlatu(:) and rlonv(:) are commons defined in "comgeom.h" ! rlatu() and rlonv() are given in radians latfi(1)=rlatu(1) lonfi(1)=0. DO j=2,jjm DO i=1,iim latfi((j-2)*iim+1+i)=rlatu(j) lonfi((j-2)*iim+1+i)=rlonv(i) ENDDO ENDDO latfi(ngridmx)=rlatu(jjp1) lonfi(ngridmx)=0. CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,aire,airefi) c======================================================================= c lecture topographie, albedo, inertie thermique, relief sous-maille c======================================================================= if (choix_1.ne.1) then ! pour ne pas avoir besoin du fichier ! surface.dat dans le cas des start c do while((relief(1:3).ne.'mol').AND.(relief(1:3).ne.'pla')) c write(*,*) c write(*,*) 'choix du relief (mola,pla)' c write(*,*) '(Topographie MGS MOLA, plat)' c read(*,fmt='(a3)') relief relief="mola" c enddo CALL datareadnc(relief,phis,alb,surfith,zmeaS,zstdS,zsigS,zgamS, . ztheS) CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,phis,phisfi) CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,surfith,surfithfi) CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,alb,albfi) CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,zmeaS,zmea) CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,zstdS,zstd) CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,zsigS,zsig) CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,zgamS,zgam) CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,ztheS,zthe) endif ! of if (choix_1.ne.1) c======================================================================= c Lecture des fichiers (start ou start_archive) c======================================================================= if (choix_1.eq.0) then write(*,*) 'Reading file START_ARCHIVE' CALL lect_start_archive(date,tsurf,tsoil,emis,q2, . t,ucov,vcov,ps,teta,phisold_newgrid,q,qsurf, & surfith,nid) write(*,*) "OK, read start_archive file" ! copy soil thermal inertia ithfi(:,:)=inertiedat(:,:) ierr= NF_CLOSE(nid) else if (choix_1.eq.1) then ! c'est l'appel a tabfi de phyeta0 qui ! permet de changer les valeurs du ! tab_cntrl Lmodif=1 tab0=0 Lmodif=1 ! Lmodif set to 1 to allow modifications in phyeta0 write(*,*) 'Reading file START' fichnom = 'start.nc' CALL dynetat0(fichnom,nqmx,vcov,ucov,teta,q,masse, . ps,phis,time) write(*,*) 'Reading file STARTFI' fichnom = 'startfi.nc' CALL phyetat0 (fichnom,tab0,Lmodif,nsoilmx,nqmx, . day_ini,time, . tsurf,tsoil,emis,q2,qsurf, !) ! temporary modif by RDW . cloudfrac,totalfrac,hice) ! copy albedo and soil thermal inertia do i=1,ngridmx albfi(i) = albedodat(i) do j=1,nsoilmx ithfi(i,j) = inertiedat(i,j) enddo ! build a surfithfi(:) using 1st layer of ithfi(:), which might ! be neede later on if reinitializing soil thermal inertia surfithfi(i)=ithfi(i,1) enddo else CALL exit(1) endif dtvr = daysec/FLOAT(day_step) dtphys = dtvr * FLOAT(iphysiq) c======================================================================= c c======================================================================= do ! infinite loop on list of changes write(*,*) write(*,*) write(*,*) 'List of possible changes :' write(*,*) '~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~' write(*,*) write(*,*) 'flat : no topography ("aquaplanet")' write(*,*) 'nuketharsis : no Tharsis bulge' write(*,*) 'bilball : uniform albedo and thermal inertia' write(*,*) 'coldspole : cold subsurface and high albedo at S.pole' write(*,*) 'qname : change tracer name' write(*,*) 'q=0 : ALL tracer =zero' write(*,*) 'q=x : give a specific uniform value to one tracer' write(*,*) 'ini_q : tracers initialisation for chemistry, water an $d ice ' write(*,*) 'ini_q-H2O : tracers initialisation for chemistry and $ice ' write(*,*) 'ini_q-iceH2O : tracers initialisation for chemistry on $ly ' write(*,*) 'noglacier : Remove tropical H2O ice if |lat|<45' write(*,*) 'watercapn : H20 ice on permanent N polar cap ' write(*,*) 'watercaps : H20 ice on permanent S polar cap ' write(*,*) 'noacglac : H2O ice across Noachis Terra' write(*,*) 'oborealis : H2O ice across Vastitas Borealis' write(*,*) 'iceball : Thick ice layer all over surface' write(*,*) 'wetstart : start with a wet atmosphere' write(*,*) 'isotherm : Isothermal Temperatures, wind set to zero' write(*,*) 'radequi : Earth-like radiative equilibrium temperature $ profile (lat-alt) and winds set to zero' write(*,*) 'coldstart : Start X K above the CO2 frost point and $set wind to zero (assumes 100% CO2)' write(*,*) 'co2ice=0 : remove CO2 polar cap' write(*,*) 'ptot : change total pressure' write(*,*) 'emis : change surface emissivity' write(*,*) 'therm_ini_s : Set soil thermal inertia to reference sur &face values' ! write(*,*) 'subsoilice_n : Put deep underground ice layer in northe ! &rn hemisphere' ! write(*,*) 'subsoilice_s : Put deep underground ice layer in southe ! &rn hemisphere' ! write(*,*) 'mons_ice : Put underground ice layer according to MONS- ! &derived data' write(*,*) write(*,*) 'Change to perform ?' write(*,*) ' (enter keyword or return to end)' write(*,*) read(*,fmt='(a20)') modif if (modif(1:1) .eq. ' ') exit ! exit loop on changes write(*,*) write(*,*) trim(modif) , ' : ' c 'flat : no topography ("aquaplanet")' c ------------------------------------- if (modif(1:len_trim(modif)) .eq. 'flat') then c set topo to zero CALL initial0(ip1jmp1,z_reel) CALL multscal(ip1jmp1,z_reel,g,phis) CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,phis,phisfi) write(*,*) 'topography set to zero.' write(*,*) 'WARNING : the subgrid topography parameters', & ' were not set to zero ! => set calllott to F' c Choice of surface pressure yes=' ' do while ((yes.ne.'y').and.(yes.ne.'n')) write(*,*) 'Do you wish to choose homogeneous surface', & 'pressure (y) or let newstart interpolate ', & ' the previous field (n)?' read(*,fmt='(a)') yes end do if (yes.eq.'y') then flagps0=.true. write(*,*) 'New value for ps (Pa) ?' 201 read(*,*,iostat=ierr) patm if(ierr.ne.0) goto 201 write(*,*) write(*,*) ' new ps everywhere (Pa) = ', patm write(*,*) do j=1,jjp1 do i=1,iip1 ps(i,j)=patm enddo enddo end if c 'nuketharsis : no tharsis bulge for Early Mars' c --------------------------------------------- else if (modif(1:len_trim(modif)) .eq. 'nuketharsis') then DO j=1,jjp1 DO i=1,iim ig=1+(j-2)*iim +i if(j.eq.1) ig=1 if(j.eq.jjp1) ig=ngridmx fact1=(((rlonv(i)*180./pi)+100)**2 + & (rlatu(j)*180./pi)**2)/65**2 fact2=exp( -fact1**2.5 ) phis(i,j) = phis(i,j) - (phis(i,j)+4000.*g)*fact2 ! if(phis(i,j).gt.2500.*g)then ! if(rlatu(j)*180./pi.gt.-80.)then ! avoid chopping south polar cap ! phis(i,j)=2500.*g ! endif ! endif ENDDO ENDDO CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,phis,phisfi) c bilball : uniform albedo, thermal inertia c ----------------------------------------- else if (modif(1:len_trim(modif)) .eq. 'bilball') then write(*,*) 'constante albedo and iner.therm:' write(*,*) 'New value for albedo (ex: 0.25) ?' 101 read(*,*,iostat=ierr) alb_bb if(ierr.ne.0) goto 101 write(*,*) write(*,*) ' uniform albedo (new value):',alb_bb write(*,*) write(*,*) 'New value for thermal inertia (eg: 247) ?' 102 read(*,*,iostat=ierr) ith_bb if(ierr.ne.0) goto 102 write(*,*) 'uniform thermal inertia (new value):',ith_bb DO j=1,jjp1 DO i=1,iip1 alb(i,j) = alb_bb ! albedo do isoil=1,nsoilmx ith(i,j,isoil) = ith_bb ! thermal inertia enddo END DO END DO ! CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,ith,ithfi) CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi) CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,alb,albfi) c coldspole : sous-sol de la calotte sud toujours froid c ----------------------------------------------------- else if (modif(1:len_trim(modif)) .eq. 'coldspole') then write(*,*)'new value for the subsurface temperature', & ' beneath the permanent southern polar cap ? (eg: 141 K)' 103 read(*,*,iostat=ierr) tsud if(ierr.ne.0) goto 103 write(*,*) write(*,*) ' new value of the subsurface temperature:',tsud c nouvelle temperature sous la calotte permanente do l=2,nsoilmx tsoil(ngridmx,l) = tsud end do write(*,*)'new value for the albedo', & 'of the permanent southern polar cap ? (eg: 0.75)' 104 read(*,*,iostat=ierr) albsud if(ierr.ne.0) goto 104 write(*,*) c ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ c Option 1: only the albedo of the pole is modified : albfi(ngridmx)=albsud write(*,*) 'ig=',ngridmx,' albedo perennial cap ', & albfi(ngridmx) c ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ c Option 2 A haute resolution : coordonnee de la vrai calotte ~ c DO j=1,jjp1 c DO i=1,iip1 c ig=1+(j-2)*iim +i c if(j.eq.1) ig=1 c if(j.eq.jjp1) ig=ngridmx c if ((rlatu(j)*180./pi.lt.-84.).and. c & (rlatu(j)*180./pi.gt.-91.).and. c & (rlonv(i)*180./pi.gt.-91.).and. c & (rlonv(i)*180./pi.lt.0.)) then cc albedo de la calotte permanente fixe a albsud c alb(i,j)=albsud c write(*,*) 'lat=',rlatu(j)*180./pi, c & ' lon=',rlonv(i)*180./pi cc fin de la condition sur les limites de la calotte permanente c end if c ENDDO c ENDDO c ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ c CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,alb,albfi) c ptot : Modification of the total pressure: ice + current atmosphere c ------------------------------------------------------------------- else if (modif(1:len_trim(modif)).eq.'ptot') then ! check if we have a co2_ice surface tracer: if (igcm_co2_ice.eq.0) then write(*,*) " No surface CO2 ice !" write(*,*) " only atmospheric pressure will be considered!" endif c calcul de la pression totale glace + atm actuelle patm=0. airetot=0. pcap=0. DO j=1,jjp1 DO i=1,iim ig=1+(j-2)*iim +i if(j.eq.1) ig=1 if(j.eq.jjp1) ig=ngridmx patm = patm + ps(i,j)*aire(i,j) airetot= airetot + aire(i,j) if (igcm_co2_ice.ne.0) then !pcap = pcap + aire(i,j)*co2ice(ig)*g pcap = pcap + aire(i,j)*qsurf(ig,igcm_co2_ice)*g endif ENDDO ENDDO ptoto = pcap + patm print*,'Current total pressure at surface (co2 ice + atm) ', & ptoto/airetot print*,'new value?' read(*,*) ptotn ptotn=ptotn*airetot patmn=ptotn-pcap print*,'ptoto,patm,ptotn,patmn' print*,ptoto,patm,ptotn,patmn print*,'Mult. factor for pressure (atm only)', patmn/patm do j=1,jjp1 do i=1,iip1 ps(i,j)=ps(i,j)*patmn/patm enddo enddo c Correction pour la conservation des traceurs yes=' ' do while ((yes.ne.'y').and.(yes.ne.'n')) write(*,*) 'Do you wish to conserve tracer total mass (y)', & ' or tracer mixing ratio (n) ?' read(*,fmt='(a)') yes end do if (yes.eq.'y') then write(*,*) 'OK : conservation of tracer total mass' DO iq =1, nqmx DO l=1,llm DO j=1,jjp1 DO i=1,iip1 q(i,j,l,iq)=q(i,j,l,iq)*patm/patmn ENDDO ENDDO ENDDO ENDDO else write(*,*) 'OK : conservation of tracer mixing ratio' end if c Correction pour la pression au niveau de la mer yes=' ' do while ((yes.ne.'y').and.(yes.ne.'n')) write(*,*) 'Do you wish fix pressure at sea level (y)', & ' or not (n) ?' read(*,fmt='(a)') yes end do if (yes.eq.'y') then write(*,*) 'Value?' read(*,*,iostat=ierr) psea DO i=1,iip1 DO j=1,jjp1 ps(i,j)=psea ENDDO ENDDO write(*,*) 'psea=',psea else write(*,*) 'no' end if c emis : change surface emissivity (added by RW) c ---------------------------------------------- else if (trim(modif).eq.'emis') then print*,'new value?' read(*,*) emisread do i=1,ngridmx emis(i)=emisread enddo c qname : change tracer name c -------------------------- else if (trim(modif).eq.'qname') then yes='y' do while (yes.eq.'y') write(*,*) 'Which tracer name do you want to change ?' do iq=1,nqmx write(*,'(i3,a3,a20)')iq,' : ',trim(tnom(iq)) enddo write(*,'(a35,i3)') & '(enter tracer number; between 1 and ',nqmx write(*,*)' or any other value to quit this option)' read(*,*) iq if ((iq.ge.1).and.(iq.le.nqmx)) then write(*,*)'Change tracer name ',trim(tnom(iq)),' to ?' read(*,*) txt tnom(iq)=txt write(*,*)'Do you want to change another tracer name (y/n)?' read(*,'(a)') yes else ! inapropiate value of iq; quit this option yes='n' endif ! of if ((iq.ge.1).and.(iq.le.nqmx)) enddo ! of do while (yes.ne.'y') c q=0 : set tracers to zero c ------------------------- else if (modif(1:len_trim(modif)).eq.'q=0') then c mise a 0 des q (traceurs) write(*,*) 'Tracers set to 0 (1.E-30 in fact)' DO iq =1, nqmx DO l=1,llm DO j=1,jjp1 DO i=1,iip1 q(i,j,l,iq)=1.e-30 ENDDO ENDDO ENDDO ENDDO c set surface tracers to zero DO iq =1, nqmx DO ig=1,ngridmx qsurf(ig,iq)=0. ENDDO ENDDO c q=x : initialise tracer manually c -------------------------------- else if (modif(1:len_trim(modif)).eq.'q=x') then write(*,*) 'Which tracer do you want to modify ?' do iq=1,nqmx write(*,*)iq,' : ',trim(tnom(iq)) enddo write(*,*) '(choose between 1 and ',nqmx,')' read(*,*) iq write(*,*)'mixing ratio of tracer ',trim(tnom(iq)), & ' ? (kg/kg)' read(*,*) val DO l=1,llm DO j=1,jjp1 DO i=1,iip1 q(i,j,l,iq)=val ENDDO ENDDO ENDDO write(*,*) 'SURFACE value of tracer ',trim(tnom(iq)), & ' ? (kg/m2)' read(*,*) val DO ig=1,ngridmx qsurf(ig,iq)=val ENDDO c ini_q : Initialize tracers for chemistry c ----------------------------------------------- else if (modif(1:len_trim(modif)) .eq. 'ini_q') then c For more than 32 layers, possible to initiate thermosphere only thermo=0 yes=' ' if (llm.gt.32) then do while ((yes.ne.'y').and.(yes.ne.'n')) write(*,*)'', & 'initialisation for thermosphere only? (y/n)' read(*,fmt='(a)') yes if (yes.eq.'y') then thermo=1 else thermo=0 endif enddo endif c call inichim_newstart(q,ps,sig,nqmx,latfi,lonfi,airefi, c $ thermo,qsurf) write(*,*) 'Chemical species initialized' if (thermo.eq.0) then c mise a 0 des qsurf (traceurs a la surface) DO iq =1, nqmx DO ig=1,ngridmx qsurf(ig,iq)=0. ENDDO ENDDO endif c ini_q-H2O : as above except for the water vapour tracer c ------------------------------------------------------ else if (modif(1:len_trim(modif)) .eq. 'ini_q-H2O') then ! for more than 32 layers, possible to initiate thermosphere only thermo=0 yes=' ' if(llm.gt.32) then do while ((yes.ne.'y').and.(yes.ne.'n')) write(*,*)'', & 'initialisation for thermosphere only? (y/n)' read(*,fmt='(a)') yes if (yes.eq.'y') then thermo=1 else thermo=0 endif enddo endif c call inichim_newstart(q,ps,sig,nqmx-1,latfi,lonfi,airefi, c $ thermo,qsurf) c write(*,*) 'Initialized chem. species exept last (H2O)' if (thermo.eq.0) then c set surface tracers to zero, except water ice DO iq =1, nqmx if (iq.ne.igcm_h2o_ice) then DO ig=1,ngridmx qsurf(ig,iq)=0. ENDDO endif ENDDO endif c ini_q-iceH2O : as above exept for ice et H2O c ----------------------------------------------- else if (modif(1:len_trim(modif)) .eq. 'ini_q-iceH2O') then c For more than 32 layers, possible to initiate thermosphere only thermo=0 yes=' ' if(llm.gt.32) then do while ((yes.ne.'y').and.(yes.ne.'n')) write(*,*)'', & 'initialisation for thermosphere only? (y/n)' read(*,fmt='(a)') yes if (yes.eq.'y') then thermo=1 else thermo=0 endif enddo endif c call inichim_newstart(q,ps,sig,nqmx-2,latfi,lonfi,airefi, c $ thermo,qsurf) c write(*,*) 'Initialized chem. species exept ice and H2O' if (thermo.eq.0) then c set surface tracers to zero, except water ice DO iq =1, nqmx if (iq.ne.igcm_h2o_ice) then DO ig=1,ngridmx qsurf(ig,iq)=0. ENDDO endif ENDDO endif c wetstart : wet atmosphere with a north to south gradient c -------------------------------------------------------- else if (modif(1:len_trim(modif)) .eq. 'wetstart') then ! check that there is indeed a water vapor tracer if (igcm_h2o_vap.eq.0) then write(*,*) "No water vapour tracer! Can't use this option" stop endif DO l=1,llm DO j=1,jjp1 DO i=1,iip1 q(i,j,l,igcm_h2o_vap)=150.e-6 * (rlatu(j)+pi/2.) / pi ENDDO ENDDO ENDDO write(*,*) 'Water mass mixing ratio at north pole=' * ,q(1,1,1,igcm_h2o_vap) write(*,*) '---------------------------south pole=' * ,q(1,jjp1,1,igcm_h2o_vap) c noglacier : remove tropical water ice (to initialize high res sim) c -------------------------------------------------- else if (modif(1:len_trim(modif)) .eq. 'noglacier') then if (igcm_h2o_ice.eq.0) then write(*,*) "No water ice tracer! Can't use this option" stop endif do ig=1,ngridmx j=(ig-2)/iim +2 if(ig.eq.1) j=1 write(*,*) 'OK: remove surface ice for |lat|<45' if (abs(rlatu(j)*180./pi).lt.45.) then qsurf(ig,igcm_h2o_ice)=0. end if end do c watercapn : H20 ice on permanent northern cap c -------------------------------------------------- else if (modif(1:len_trim(modif)) .eq. 'watercapn') then if (igcm_h2o_ice.eq.0) then write(*,*) "No water ice tracer! Can't use this option" stop endif DO j=1,jjp1 DO i=1,iim ig=1+(j-2)*iim +i if(j.eq.1) ig=1 if(j.eq.jjp1) ig=ngridmx if (rlatu(j)*180./pi.gt.80.) then qsurf(ig,igcm_h2o_ice)=3.4e3 !do isoil=1,nsoilmx ! ith(i,j,isoil) = 18000. ! thermal inertia !enddo write(*,*)' ==> Ice mesh North boundary (deg)= ', & rlatv(j-1)*180./pi end if ENDDO ENDDO CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi) c$$$ do ig=1,ngridmx c$$$ j=(ig-2)/iim +2 c$$$ if(ig.eq.1) j=1 c$$$ if (rlatu(j)*180./pi.gt.80.) then c$$$ c$$$ qsurf(ig,igcm_h2o_ice)=1.e5 c$$$ qsurf(ig,igcm_h2o_vap)=0.0!1.e5 c$$$ c$$$ write(*,*) 'ig=',ig,' H2O ice mass (kg/m2)= ', c$$$ & qsurf(ig,igcm_h2o_ice) c$$$ c$$$ write(*,*)' ==> Ice mesh South boundary (deg)= ', c$$$ & rlatv(j)*180./pi c$$$ end if c$$$ enddo c watercaps : H20 ice on permanent southern cap c ------------------------------------------------- else if (modif(1:len_trim(modif)) .eq. 'watercaps') then if (igcm_h2o_ice.eq.0) then write(*,*) "No water ice tracer! Can't use this option" stop endif DO j=1,jjp1 DO i=1,iim ig=1+(j-2)*iim +i if(j.eq.1) ig=1 if(j.eq.jjp1) ig=ngridmx if (rlatu(j)*180./pi.lt.-80.) then qsurf(ig,igcm_h2o_ice)=3.4e3 !do isoil=1,nsoilmx ! ith(i,j,isoil) = 18000. ! thermal inertia !enddo write(*,*)' ==> Ice mesh South boundary (deg)= ', & rlatv(j-1)*180./pi end if ENDDO ENDDO CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi) c$$$ do ig=1,ngridmx c$$$ j=(ig-2)/iim +2 c$$$ if(ig.eq.1) j=1 c$$$ if (rlatu(j)*180./pi.lt.-80.) then c$$$ qsurf(ig,igcm_h2o_ice)=1.e5 c$$$ qsurf(ig,igcm_h2o_vap)=0.0 !1.e5 c$$$ c$$$ write(*,*) 'ig=',ig,' H2O ice mass (kg/m2)= ', c$$$ & qsurf(ig,igcm_h2o_ice) c$$$ write(*,*)' ==> Ice mesh North boundary (deg)= ', c$$$ & rlatv(j-1)*180./pi c$$$ end if c$$$ enddo c noacglac : H2O ice across highest terrain c -------------------------------------------- else if (modif(1:len_trim(modif)) .eq. 'noacglac') then if (igcm_h2o_ice.eq.0) then write(*,*) "No water ice tracer! Can't use this option" stop endif DO j=1,jjp1 DO i=1,iim ig=1+(j-2)*iim +i if(j.eq.1) ig=1 if(j.eq.jjp1) ig=ngridmx if(phis(i,j).gt.1000.*g)then alb(i,j) = 0.5 ! snow value do isoil=1,nsoilmx ith(i,j,isoil) = 18000. ! thermal inertia ! this leads to rnat set to 'ocean' in physiq.F90 ! actually no, because it is soil not surface enddo endif ENDDO ENDDO CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi) CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,alb,albfi) CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,phis,phisfi) c oborealis : H2O oceans across Vastitas Borealis c ----------------------------------------------- else if (modif(1:len_trim(modif)) .eq. 'oborealis') then if (igcm_h2o_ice.eq.0) then write(*,*) "No water ice tracer! Can't use this option" stop endif DO j=1,jjp1 DO i=1,iim ig=1+(j-2)*iim +i if(j.eq.1) ig=1 if(j.eq.jjp1) ig=ngridmx if(phis(i,j).lt.-4000.*g)then ! if( (phis(i,j).lt.-4000.*g) ! & .and. (rlatu(j)*180./pi.lt.0.) )then ! south hemisphere only ! phis(i,j)=-8000.0*g ! proper ocean phis(i,j)=-4000.0*g ! scanty ocean alb(i,j) = 0.07 ! oceanic value do isoil=1,nsoilmx ith(i,j,isoil) = 18000. ! thermal inertia ! this leads to rnat set to 'ocean' in physiq.F90 ! actually no, because it is soil not surface enddo endif ENDDO ENDDO CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi) CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,alb,albfi) CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,phis,phisfi) c iborealis : H2O ice in Northern plains c -------------------------------------- else if (modif(1:len_trim(modif)) .eq. 'iborealis') then if (igcm_h2o_ice.eq.0) then write(*,*) "No water ice tracer! Can't use this option" stop endif DO j=1,jjp1 DO i=1,iim ig=1+(j-2)*iim +i if(j.eq.1) ig=1 if(j.eq.jjp1) ig=ngridmx if(phis(i,j).lt.-4000.*g)then !qsurf(ig,igcm_h2o_ice)=1.e3 qsurf(ig,igcm_h2o_ice)=241.4 ! to make total 33 kg m^-2 endif ENDDO ENDDO CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi) CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,alb,albfi) CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,phis,phisfi) c oceanball : H2O liquid everywhere c ---------------------------- else if (modif(1:len_trim(modif)) .eq. 'oceanball') then if (igcm_h2o_ice.eq.0) then write(*,*) "No water ice tracer! Can't use this option" stop endif DO j=1,jjp1 DO i=1,iim ig=1+(j-2)*iim +i if(j.eq.1) ig=1 if(j.eq.jjp1) ig=ngridmx qsurf(ig,igcm_h2o_vap)=0.0 ! for ocean, this is infinite source qsurf(ig,igcm_h2o_ice)=0.0 alb(i,j) = 0.07 ! ocean value do isoil=1,nsoilmx ith(i,j,isoil) = 18000. ! thermal inertia !ith(i,j,isoil) = 50. ! extremely low for test ! this leads to rnat set to 'ocean' in physiq.F90 enddo ENDDO ENDDO CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi) CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,alb,albfi) CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,phis,phisfi) c iceball : H2O ice everywhere c ---------------------------- else if (modif(1:len_trim(modif)) .eq. 'iceball') then if (igcm_h2o_ice.eq.0) then write(*,*) "No water ice tracer! Can't use this option" stop endif DO j=1,jjp1 DO i=1,iim ig=1+(j-2)*iim +i if(j.eq.1) ig=1 if(j.eq.jjp1) ig=ngridmx qsurf(ig,igcm_h2o_vap)=-50. ! for ocean, this is infinite source qsurf(ig,igcm_h2o_ice)=50. ! == to 0.05 m of oceanic ice alb(i,j) = 0.6 ! ice albedo value do isoil=1,nsoilmx !ith(i,j,isoil) = 18000. ! thermal inertia ! this leads to rnat set to 'ocean' in physiq.F90 enddo ENDDO ENDDO CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi) CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,alb,albfi) c isotherm : Isothermal temperatures and no winds c ----------------------------------------------- else if (modif(1:len_trim(modif)) .eq. 'isotherm') then write(*,*)'Isothermal temperature of the atmosphere, & surface and subsurface' write(*,*) 'Value of this temperature ? :' 203 read(*,*,iostat=ierr) Tiso if(ierr.ne.0) goto 203 do ig=1, ngridmx tsurf(ig) = Tiso end do do l=2,nsoilmx do ig=1, ngridmx tsoil(ig,l) = Tiso end do end do DO j=1,jjp1 DO i=1,iim Do l=1,llm Tset(i,j,l)=Tiso end do end do end do flagtset=.true. call initial0(llm*ip1jmp1,ucov) call initial0(llm*ip1jm,vcov) call initial0(ngridmx*(llm+1),q2) c radequi : Radiative equilibrium profile of temperatures and no winds c -------------------------------------------------------------------- else if (modif(1:len_trim(modif)) .eq. 'radequi') then write(*,*)'radiative equilibrium temperature profile' do ig=1, ngridmx teque= 335.0-60.0*sin(latfi(ig))*sin(latfi(ig))- & 10.0*cos(latfi(ig))*cos(latfi(ig)) tsurf(ig) = MAX(220.0,teque) end do do l=2,nsoilmx do ig=1, ngridmx tsoil(ig,l) = tsurf(ig) end do end do DO j=1,jjp1 DO i=1,iim Do l=1,llm teque=335.-60.0*sin(rlatu(j))*sin(rlatu(j))- & 10.0*cos(rlatu(j))*cos(rlatu(j)) Tset(i,j,l)=MAX(220.0,teque) end do end do end do flagtset=.true. call initial0(llm*ip1jmp1,ucov) call initial0(llm*ip1jm,vcov) call initial0(ngridmx*(llm+1),q2) c coldstart : T set 1K above CO2 frost point and no winds c ------------------------------------------------ else if (modif(1:len_trim(modif)) .eq. 'coldstart') then write(*,*)'set temperature of the atmosphere,' &,'surface and subsurface how many degrees above CO2 frost point?' 204 read(*,*,iostat=ierr) Tabove if(ierr.ne.0) goto 204 DO j=1,jjp1 DO i=1,iim ig=1+(j-2)*iim +i if(j.eq.1) ig=1 if(j.eq.jjp1) ig=ngridmx tsurf(ig) = (-3167.8)/(log(.01*ps(i,j))-23.23)+Tabove END DO END DO do l=1,nsoilmx do ig=1, ngridmx tsoil(ig,l) = tsurf(ig) end do end do DO j=1,jjp1 DO i=1,iim Do l=1,llm pp = aps(l) +bps(l)*ps(i,j) Tset(i,j,l)=(-3167.8)/(log(.01*pp)-23.23)+Tabove end do end do end do flagtset=.true. call initial0(llm*ip1jmp1,ucov) call initial0(llm*ip1jm,vcov) call initial0(ngridmx*(llm+1),q2) c co2ice=0 : remove CO2 polar ice caps' c ------------------------------------------------ else if (modif(1:len_trim(modif)) .eq. 'co2ice=0') then ! check that there is indeed a co2_ice tracer ... if (igcm_co2_ice.ne.0) then do ig=1,ngridmx !co2ice(ig)=0 qsurf(ig,igcm_co2_ice)=0 emis(ig)=emis(ngridmx/2) end do else write(*,*) "Can't remove CO2 ice!! (no co2_ice tracer)" endif ! therm_ini_s: (re)-set soil thermal inertia to reference surface values ! ---------------------------------------------------------------------- else if (modif(1:len_trim(modif)).eq.'therm_ini_s') then ! write(*,*)"surfithfi(1):",surfithfi(1) do isoil=1,nsoilmx inertiedat(1:ngridmx,isoil)=surfithfi(1:ngridmx) enddo write(*,*)'OK: Soil thermal inertia has been reset to referenc &e surface values' ! write(*,*)"inertiedat(1,1):",inertiedat(1,1) ithfi(:,:)=inertiedat(:,:) ! recast ithfi() onto ith() call gr_fi_dyn(nsoilmx,ngridmx,iip1,jjp1,ithfi,ith) ! Check: ! do i=1,iip1 ! do j=1,jjp1 ! do isoil=1,nsoilmx ! write(77,*) i,j,isoil," ",ith(i,j,isoil) ! enddo ! enddo ! enddo c$$$! subsoilice_n: Put deep ice layer in northern hemisphere soil c$$$! ------------------------------------------------------------ c$$$ c$$$ else if (modif(1:len_trim(modif)).eq.'subsoilice_n') then c$$$ c$$$ write(*,*)'From which latitude (in deg.), up to the north pole, c$$$ &should we put subterranean ice?' c$$$ ierr=1 c$$$ do while (ierr.ne.0) c$$$ read(*,*,iostat=ierr) val c$$$ if (ierr.eq.0) then ! got a value c$$$ ! do a sanity check c$$$ if((val.lt.0.).or.(val.gt.90)) then c$$$ write(*,*)'Latitude should be between 0 and 90 deg. !!!' c$$$ ierr=1 c$$$ else ! find corresponding jref (nearest latitude) c$$$ ! note: rlatu(:) contains decreasing values of latitude c$$$ ! starting from PI/2 to -PI/2 c$$$ do j=1,jjp1 c$$$ if ((rlatu(j)*180./pi.ge.val).and. c$$$ & (rlatu(j+1)*180./pi.le.val)) then c$$$ ! find which grid point is nearest to val: c$$$ if (abs(rlatu(j)*180./pi-val).le. c$$$ & abs((rlatu(j+1)*180./pi-val))) then c$$$ jref=j c$$$ else c$$$ jref=j+1 c$$$ endif c$$$ c$$$ write(*,*)'Will use nearest grid latitude which is:', c$$$ & rlatu(jref)*180./pi c$$$ endif c$$$ enddo ! of do j=1,jjp1 c$$$ endif ! of if((val.lt.0.).or.(val.gt.90)) c$$$ endif !of if (ierr.eq.0) c$$$ enddo ! of do while c$$$ c$$$ ! Build layers() (as in soil_settings.F) c$$$ val2=sqrt(mlayer(0)*mlayer(1)) c$$$ val3=mlayer(1)/mlayer(0) c$$$ do isoil=1,nsoilmx c$$$ layer(isoil)=val2*(val3**(isoil-1)) c$$$ enddo c$$$ c$$$ write(*,*)'At which depth (in m.) does the ice layer begin?' c$$$ write(*,*)'(currently, the deepest soil layer extends down to:' c$$$ & ,layer(nsoilmx),')' c$$$ ierr=1 c$$$ do while (ierr.ne.0) c$$$ read(*,*,iostat=ierr) val2 c$$$! write(*,*)'val2:',val2,'ierr=',ierr c$$$ if (ierr.eq.0) then ! got a value, but do a sanity check c$$$ if(val2.gt.layer(nsoilmx)) then c$$$ write(*,*)'Depth should be less than ',layer(nsoilmx) c$$$ ierr=1 c$$$ endif c$$$ if(val2.lt.layer(1)) then c$$$ write(*,*)'Depth should be more than ',layer(1) c$$$ ierr=1 c$$$ endif c$$$ endif c$$$ enddo ! of do while c$$$ c$$$ ! find the reference index iref the depth corresponds to c$$$! if (val2.lt.layer(1)) then c$$$! iref=1 c$$$! else c$$$ do isoil=1,nsoilmx-1 c$$$ if((val2.gt.layer(isoil)).and.(val2.lt.layer(isoil+1))) c$$$ & then c$$$ iref=isoil c$$$ exit c$$$ endif c$$$ enddo c$$$! endif c$$$ c$$$! write(*,*)'iref:',iref,' jref:',jref c$$$! write(*,*)'layer',layer c$$$! write(*,*)'mlayer',mlayer c$$$ c$$$ ! thermal inertia of the ice: c$$$ ierr=1 c$$$ do while (ierr.ne.0) c$$$ write(*,*)'What is the value of subterranean ice thermal inert c$$$ &ia? (e.g.: 2000)' c$$$ read(*,*,iostat=ierr)iceith c$$$ enddo ! of do while c$$$ c$$$ ! recast ithfi() onto ith() c$$$ call gr_fi_dyn(nsoilmx,ngridmx,iip1,jjp1,ithfi,ith) c$$$ c$$$ do j=1,jref c$$$! write(*,*)'j:',j,'rlatu(j)*180./pi:',rlatu(j)*180./pi c$$$ do i=1,iip1 ! loop on longitudes c$$$ ! Build "equivalent" thermal inertia for the mixed layer c$$$ ith(i,j,iref+1)=sqrt((layer(iref+1)-layer(iref))/ c$$$ & (((val2-layer(iref))/(ith(i,j,iref)**2))+ c$$$ & ((layer(iref+1)-val2)/(iceith)**2))) c$$$ ! Set thermal inertia of lower layers c$$$ do isoil=iref+2,nsoilmx c$$$ ith(i,j,isoil)=iceith ! ice c$$$ enddo c$$$ enddo ! of do i=1,iip1 c$$$ enddo ! of do j=1,jjp1 c$$$ c$$$ c$$$ CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi) c$$$ c$$$! do i=1,nsoilmx c$$$! write(*,*)'i:',i,'ithfi(1,i):',ithfi(1,i) c$$$! enddo c$$$! subsoilice_s: Put deep ice layer in southern hemisphere soil c$$$! ------------------------------------------------------------ c$$$ c$$$ else if (modif(1:len_trim(modif)).eq.'subsoilice_s') then c$$$ c$$$ write(*,*)'From which latitude (in deg.), down to the south pol c$$$ &e, should we put subterranean ice?' c$$$ ierr=1 c$$$ do while (ierr.ne.0) c$$$ read(*,*,iostat=ierr) val c$$$ if (ierr.eq.0) then ! got a value c$$$ ! do a sanity check c$$$ if((val.gt.0.).or.(val.lt.-90)) then c$$$ write(*,*)'Latitude should be between 0 and -90 deg. !!!' c$$$ ierr=1 c$$$ else ! find corresponding jref (nearest latitude) c$$$ ! note: rlatu(:) contains decreasing values of latitude c$$$ ! starting from PI/2 to -PI/2 c$$$ do j=1,jjp1 c$$$ if ((rlatu(j)*180./pi.ge.val).and. c$$$ & (rlatu(j+1)*180./pi.le.val)) then c$$$ ! find which grid point is nearest to val: c$$$ if (abs(rlatu(j)*180./pi-val).le. c$$$ & abs((rlatu(j+1)*180./pi-val))) then c$$$ jref=j c$$$ else c$$$ jref=j+1 c$$$ endif c$$$ c$$$ write(*,*)'Will use nearest grid latitude which is:', c$$$ & rlatu(jref)*180./pi c$$$ endif c$$$ enddo ! of do j=1,jjp1 c$$$ endif ! of if((val.lt.0.).or.(val.gt.90)) c$$$ endif !of if (ierr.eq.0) c$$$ enddo ! of do while c$$$ c$$$ ! Build layers() (as in soil_settings.F) c$$$ val2=sqrt(mlayer(0)*mlayer(1)) c$$$ val3=mlayer(1)/mlayer(0) c$$$ do isoil=1,nsoilmx c$$$ layer(isoil)=val2*(val3**(isoil-1)) c$$$ enddo c$$$ c$$$ write(*,*)'At which depth (in m.) does the ice layer begin?' c$$$ write(*,*)'(currently, the deepest soil layer extends down to:' c$$$ & ,layer(nsoilmx),')' c$$$ ierr=1 c$$$ do while (ierr.ne.0) c$$$ read(*,*,iostat=ierr) val2 c$$$! write(*,*)'val2:',val2,'ierr=',ierr c$$$ if (ierr.eq.0) then ! got a value, but do a sanity check c$$$ if(val2.gt.layer(nsoilmx)) then c$$$ write(*,*)'Depth should be less than ',layer(nsoilmx) c$$$ ierr=1 c$$$ endif c$$$ if(val2.lt.layer(1)) then c$$$ write(*,*)'Depth should be more than ',layer(1) c$$$ ierr=1 c$$$ endif c$$$ endif c$$$ enddo ! of do while c$$$ c$$$ ! find the reference index iref the depth corresponds to c$$$ do isoil=1,nsoilmx-1 c$$$ if((val2.gt.layer(isoil)).and.(val2.lt.layer(isoil+1))) c$$$ & then c$$$ iref=isoil c$$$ exit c$$$ endif c$$$ enddo c$$$ c$$$! write(*,*)'iref:',iref,' jref:',jref c$$$ c$$$ ! thermal inertia of the ice: c$$$ ierr=1 c$$$ do while (ierr.ne.0) c$$$ write(*,*)'What is the value of subterranean ice thermal inert c$$$ &ia? (e.g.: 2000)' c$$$ read(*,*,iostat=ierr)iceith c$$$ enddo ! of do while c$$$ c$$$ ! recast ithfi() onto ith() c$$$ call gr_fi_dyn(nsoilmx,ngridmx,iip1,jjp1,ithfi,ith) c$$$ c$$$ do j=jref,jjp1 c$$$! write(*,*)'j:',j,'rlatu(j)*180./pi:',rlatu(j)*180./pi c$$$ do i=1,iip1 ! loop on longitudes c$$$ ! Build "equivalent" thermal inertia for the mixed layer c$$$ ith(i,j,iref+1)=sqrt((layer(iref+1)-layer(iref))/ c$$$ & (((val2-layer(iref))/(ith(i,j,iref)**2))+ c$$$ & ((layer(iref+1)-val2)/(iceith)**2))) c$$$ ! Set thermal inertia of lower layers c$$$ do isoil=iref+2,nsoilmx c$$$ ith(i,j,isoil)=iceith ! ice c$$$ enddo c$$$ enddo ! of do i=1,iip1 c$$$ enddo ! of do j=jref,jjp1 c$$$ c$$$ c$$$ CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi) c$$$c 'mons_ice' : use MONS data to build subsurface ice table c$$$c -------------------------------------------------------- c$$$ else if (modif(1:len_trim(modif)).eq.'mons_ice') then c$$$ c$$$ ! 1. Load MONS data c$$$ call load_MONS_data(MONS_Hdn,MONS_d21) c$$$ c$$$ ! 2. Get parameters from user c$$$ ierr=1 c$$$ do while (ierr.ne.0) c$$$ write(*,*) "Coefficient to apply to MONS 'depth' in Northern", c$$$ & " Hemisphere?" c$$$ write(*,*) " (should be somewhere between 3.2e-4 and 1.3e-3)" c$$$ read(*,*,iostat=ierr) MONS_coeffN c$$$ enddo c$$$ ierr=1 c$$$ do while (ierr.ne.0) c$$$ write(*,*) "Coefficient to apply to MONS 'depth' in Southern", c$$$ & " Hemisphere?" c$$$ write(*,*) " (should be somewhere between 3.2e-4 and 1.3e-3)" c$$$ read(*,*,iostat=ierr) MONS_coeffS c$$$ enddo c$$$ ierr=1 c$$$ do while (ierr.ne.0) c$$$ write(*,*) "Value of subterranean ice thermal inertia?" c$$$ write(*,*) " (e.g.: 2000, or perhaps 2290)" c$$$ read(*,*,iostat=ierr) iceith c$$$ enddo c$$$ c$$$ ! 3. Build subterranean thermal inertia c$$$ c$$$ ! initialise subsurface inertia with reference surface values c$$$ do isoil=1,nsoilmx c$$$ ithfi(1:ngridmx,isoil)=surfithfi(1:ngridmx) c$$$ enddo c$$$ ! recast ithfi() onto ith() c$$$ call gr_fi_dyn(nsoilmx,ngridmx,iip1,jjp1,ithfi,ith) c$$$ c$$$ do i=1,iip1 ! loop on longitudes c$$$ do j=1,jjp1 ! loop on latitudes c$$$ ! set MONS_coeff c$$$ if (rlatu(j).ge.0) then ! northern hemisphere c$$$ ! N.B: rlatu(:) contains decreasing values of latitude c$$$ ! starting from PI/2 to -PI/2 c$$$ MONS_coeff=MONS_coeffN c$$$ else ! southern hemisphere c$$$ MONS_coeff=MONS_coeffS c$$$ endif c$$$ ! check if we should put subterranean ice c$$$ if (MONS_Hdn(i,j).ge.14.0) then ! no ice if Hdn<14% c$$$ ! compute depth at which ice lies: c$$$ val=MONS_d21(i,j)*MONS_coeff c$$$ ! compute val2= the diurnal skin depth of surface inertia c$$$ ! assuming a volumetric heat cap. of C=1.e6 J.m-3.K-1 c$$$ val2=ith(i,j,1)*1.e-6*sqrt(88775./3.14159) c$$$ if (val.lt.val2) then c$$$ ! ice must be below the (surface inertia) diurnal skin depth c$$$ val=val2 c$$$ endif c$$$ if (val.lt.layer(nsoilmx)) then ! subterranean ice c$$$ ! find the reference index iref that depth corresponds to c$$$ iref=0 c$$$ do isoil=1,nsoilmx-1 c$$$ if ((val.ge.layer(isoil)).and.(val.lt.layer(isoil+1))) c$$$ & then c$$$ iref=isoil c$$$ exit c$$$ endif c$$$ enddo c$$$ ! Build "equivalent" thermal inertia for the mixed layer c$$$ ith(i,j,iref+1)=sqrt((layer(iref+1)-layer(iref))/ c$$$ & (((val-layer(iref))/(ith(i,j,iref+1)**2))+ c$$$ & ((layer(iref+1)-val)/(iceith)**2))) c$$$ ! Set thermal inertia of lower layers c$$$ do isoil=iref+2,nsoilmx c$$$ ith(i,j,isoil)=iceith c$$$ enddo c$$$ endif ! of if (val.lt.layer(nsoilmx)) c$$$ endif ! of if (MONS_Hdn(i,j).lt.14.0) c$$$ enddo ! do j=1,jjp1 c$$$ enddo ! do i=1,iip1 c$$$ c$$$! Check: c$$$! do i=1,iip1 c$$$! do j=1,jjp1 c$$$! do isoil=1,nsoilmx c$$$! write(77,*) i,j,isoil," ",ith(i,j,isoil) c$$$! enddo c$$$! enddo c$$$! enddo c$$$ c$$$ ! recast ith() into ithfi() c$$$ CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi) c$$$ c$$$ else c$$$ write(*,*) ' Unknown (misspelled?) option!!!' end if ! of if (modif(1:len_trim(modif)) .eq. '...') elseif ... enddo ! of do ! infinite loop on liste of changes 999 continue c======================================================================= c Initialisation for cloud fraction and oceanic ice (added by BC 2010) c======================================================================= DO ig=1, ngridmx totalfrac(ig)=0.5 DO l=1,llm cloudfrac(ig,l)=0.5 ENDDO ENDDO c======================================================== ! DO ig=1,ngridmx ! IF(tsurf(ig) .LT. 273.16-1.8) THEN ! hice(ig)=(273.16-1.8-tsurf(ig))/(273.16-1.8-240)*1. ! hice(ig)=min(hice(ig),1.0) ! ENDIF ! ENDDO c======================================================================= c Correct pressure on the new grid (menu 0) c======================================================================= if ((choix_1.eq.0).and.(.not.flagps0)) then r = 1000.*8.31/mugaz phi0=0.0 do j=1,jjp1 do i=1,iip1 phi0 = phi0+phis(i,j)*aire(i,j) end do end do phi0=phi0/airetot do j=1,jjp1 do i=1,iip1 ps(i,j) = ps(i,j) * ! . exp((phisold_newgrid(i,j)-phis(i,j)) / . exp((phi0-phis(i,j)) / . (t(i,j,1) * r)) end do end do ! we must renormalise pressures AGAIN, because exp(-phi) is nonlinear ptot=0.0 do j=1,jjp1 do i=1,iip1 ptot=ptot+ps(i,j)*aire(i,j) enddo enddo do j=1,jjp1 do i=1,iip1 ps(i,j)=ps(i,j)*ptotn/ptot enddo enddo ! periodicity of surface ps in longitude do j=1,jjp1 ps(1,j) = ps(iip1,j) end do end if c======================================================================= c======================================================================= c======================================================================= c Initialisation de la physique / ecriture de newstartfi : c======================================================================= CALL inifilr CALL pression(ip1jmp1, ap, bp, ps, p3d) c----------------------------------------------------------------------- c Initialisation pks: c----------------------------------------------------------------------- CALL exner_hyb(ip1jmp1, ps, p3d, beta, pks, pk, pkf) ! Calcul de la temperature potentielle teta if (flagtset) then DO l=1,llm DO j=1,jjp1 DO i=1,iim teta(i,j,l) = Tset(i,j,l) * cpp/pk(i,j,l) ENDDO teta (iip1,j,l)= teta (1,j,l) ENDDO ENDDO else if (choix_1.eq.0) then DO l=1,llm DO j=1,jjp1 DO i=1,iim teta(i,j,l) = t(i,j,l) * cpp/pk(i,j,l) ENDDO teta (iip1,j,l)= teta (1,j,l) ENDDO ENDDO endif C Calcul intermediaire c if (choix_1.eq.0) then CALL massdair( p3d, masse ) c print *,' ALPHAX ',alphax DO l = 1, llm DO i = 1, iim xppn(i) = aire( i, 1 ) * masse( i , 1 , l ) xpps(i) = aire( i,jjp1 ) * masse( i , jjp1 , l ) ENDDO xpn = SUM(xppn)/apoln xps = SUM(xpps)/apols DO i = 1, iip1 masse( i , 1 , l ) = xpn masse( i , jjp1 , l ) = xps ENDDO ENDDO endif phis(iip1,:) = phis(1,:) CALL inidissip ( lstardis, nitergdiv, nitergrot, niterh, * tetagdiv, tetagrot , tetatemp ) itau=0 if (choix_1.eq.0) then day_ini=int(date) endif c CALL geopot ( ip1jmp1, teta , pk , pks, phis , phi ) CALL caldyn0( itau,ucov,vcov,teta,ps,masse,pk,phis , * phi,w, pbaru,pbarv,day_ini+time ) c CALL caldyn c $ ( itau,ucov,vcov,teta,ps,masse,pk,pkf,phis , c $ phi,conser,du,dv,dteta,dp,w, pbaru,pbarv, day_ini ) CALL dynredem0("restart.nc",day_ini,anneeref,phis,nqmx) CALL dynredem1("restart.nc",0.0,vcov,ucov,teta,q,nqmx,masse,ps) C C Ecriture etat initial physique C ! do ig=1,ngridmx ! print*,'surface ice in newstart=',qsurf(ig,igcm_h2o_ice) ! print*,'surface liq in newstart=',qsurf(ig,igcm_h2o_vap) ! enddo ! stop call physdem1("restartfi.nc",lonfi,latfi,nsoilmx,nqmx, . dtphys,float(day_ini), . time,tsurf,tsoil,emis,q2,qsurf, . airefi,albfi,ithfi,zmea,zstd,zsig,zgam,zthe, . cloudfrac,totalfrac,hice) c======================================================================= c Formats c======================================================================= 1 FORMAT(//10x,'la valeur de im =',i4,2x,'lue sur le fichier de dema *rrage est differente de la valeur parametree iim =',i4//) 2 FORMAT(//10x,'la valeur de jm =',i4,2x,'lue sur le fichier de dema *rrage est differente de la valeur parametree jjm =',i4//) 3 FORMAT(//10x,'la valeur de lllm =',i4,2x,'lue sur le fichier demar *rage est differente de la valeur parametree llm =',i4//) end !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! subroutine load_MONS_data(MONS_Hdn,MONS_d21) use datafile_mod, only: datadir implicit none ! routine to load Benedicte Diez MONS dataset, fill in date in southern ! polar region, and interpolate the result onto the GCM grid #include"dimensions.h" #include"paramet.h" #include"comgeom.h" ! arguments: real,intent(out) :: MONS_Hdn(iip1,jjp1) ! Hdn: %WEH=Mass fraction of H2O real,intent(out) :: MONS_d21(iip1,jjp1) ! ice table "depth" (in kg/m2) ! N.B MONS datasets should be of dimension (iip1,jjp1) ! local variables: character(len=88) :: filename="results_MONS_lat_lon_H_depth.txt" character(len=88) :: txt ! to store some text integer :: ierr,i,j integer,parameter :: nblon=180 ! number of longitudes of MONS datasets integer,parameter :: nblat=90 ! number of latitudes of MONS datasets real :: pi real :: longitudes(nblon) ! MONS dataset longitudes real :: latitudes(nblat) ! MONS dataset latitudes ! MONS dataset: mass fraction of H2O where H is assumed to be in H2O real :: Hdn(nblon,nblat) real :: d21(nblon,nblat)! MONS dataset "depth" (g/cm2) ! Extended MONS dataset (for interp_horiz) real :: Hdnx(nblon+1,nblat) real :: d21x(nblon+1,nblat) real :: lon_bound(nblon+1) ! longitude boundaries real :: lat_bound(nblat-1) ! latitude boundaries ! 1. Initializations: write(*,*) "Loading MONS data" ! Open MONS datafile: open(42,file=trim(datadir)//"/"//trim(filename), & status="old",iostat=ierr) if (ierr/=0) then write(*,*) "Error in load_MONS_data:" write(*,*) "Failed opening file ", & trim(datadir)//"/"//trim(filename) write(*,*)'Check that your path to datagcm:',trim(datadir) write(*,*)' is correct. You can change it in callphys.def with:' write(*,*)' datadir = /absolute/path/to/datagcm' write(*,*)'If necessary ',trim(filename), & ' (and other datafiles)' write(*,*)' can be obtained online at:' write(*,*)' http://www.lmd.jussieu.fr/~forget/datagcm/datafile' CALL ABORT else ! skip first line of file (dummy read) read(42,*) txt endif pi=2.*asin(1.) !2. Load MONS data (on MONS grid) do j=1,nblat do i=1,nblon ! swap latitude index so latitudes go from north pole to south pole: read(42,*) latitudes(nblat-j+1),longitudes(i), & Hdn(i,nblat-j+1),d21(i,nblat-j+1) ! multiply d21 by 10 to convert from g/cm2 to kg/m2 d21(i,nblat-j+1)=d21(i,nblat-j+1)*10.0 enddo enddo close(42) ! there is unfortunately no d21 data for latitudes -77 to -90 ! so we build some by linear interpolation between values at -75 ! and assuming d21=0 at the pole do j=84,90 ! latitudes(84)=-77 ; latitudes(83)=-75 do i=1,nblon d21(i,j)=d21(i,83)*((latitudes(j)+90)/15.0) enddo enddo ! 3. Build extended MONS dataset & boundaries (for interp_horiz) ! longitude boundaries (in radians): do i=1,nblon ! NB: MONS data is every 2 degrees in longitude lon_bound(i)=(longitudes(i)+1.0)*pi/180.0 enddo ! extra 'modulo' value lon_bound(nblon+1)=lon_bound(1)+2.0*pi ! latitude boundaries (in radians): do j=1,nblat-1 ! NB: Mons data is every 2 degrees in latitude lat_bound(j)=(latitudes(j)-1.0)*pi/180.0 enddo ! MONS datasets: do j=1,nblat Hdnx(1:nblon,j)=Hdn(1:nblon,j) Hdnx(nblon+1,j)=Hdnx(1,j) d21x(1:nblon,j)=d21(1:nblon,j) d21x(nblon+1,j)=d21x(1,j) enddo ! Interpolate onto GCM grid call interp_horiz(Hdnx,MONS_Hdn,nblon,nblat-1,iim,jjm,1, & lon_bound,lat_bound,rlonu,rlatv) call interp_horiz(d21x,MONS_d21,nblon,nblat-1,iim,jjm,1, & lon_bound,lat_bound,rlonu,rlatv) end subroutine