[1216] | 1 | module phyredem |
---|
[135] | 2 | |
---|
[1216] | 3 | implicit none |
---|
[135] | 4 | |
---|
[1216] | 5 | contains |
---|
[135] | 6 | |
---|
[1216] | 7 | subroutine physdem0(filename,lonfi,latfi,nsoil,ngrid,nlay,nq, & |
---|
| 8 | phystep,day_ini,time,airefi, & |
---|
| 9 | alb,ith,pzmea,pzstd,pzsig,pzgam,pzthe) |
---|
| 10 | ! create physics restart file and write time-independent variables |
---|
[1222] | 11 | use comsoil_h, only: volcapa, mlayer |
---|
[1216] | 12 | use comgeomfi_h, only: area |
---|
[1222] | 13 | use surfdat_h, only: zmea, zstd, zsig, zgam, zthe, & |
---|
[1216] | 14 | albedice, emisice, emissiv, & |
---|
| 15 | iceradius, dtemisice, phisfi |
---|
| 16 | use iostart, only : open_restartphy, close_restartphy, & |
---|
| 17 | put_var, put_field, length |
---|
| 18 | use mod_grid_phy_lmdz, only : klon_glo |
---|
[135] | 19 | |
---|
[1216] | 20 | implicit none |
---|
[135] | 21 | #include "planete.h" |
---|
| 22 | #include "comcstfi.h" |
---|
[1216] | 23 | character(len=*), intent(in) :: filename |
---|
| 24 | real,intent(in) :: lonfi(ngrid) |
---|
| 25 | real,intent(in) :: latfi(ngrid) |
---|
| 26 | integer,intent(in) :: nsoil |
---|
| 27 | integer,intent(in) :: ngrid |
---|
| 28 | integer,intent(in) :: nlay |
---|
| 29 | integer,intent(in) :: nq |
---|
| 30 | real,intent(in) :: phystep |
---|
| 31 | real,intent(in) :: day_ini |
---|
| 32 | real,intent(in) :: time |
---|
| 33 | real,intent(in) :: airefi(ngrid) |
---|
| 34 | real,intent(in) :: alb(ngrid) |
---|
| 35 | real,intent(in) :: ith(ngrid,nsoil) |
---|
| 36 | real,intent(in) :: pzmea(ngrid) |
---|
| 37 | real,intent(in) :: pzstd(ngrid) |
---|
| 38 | real,intent(in) :: pzsig(ngrid) |
---|
| 39 | real,intent(in) :: pzgam(ngrid) |
---|
| 40 | real,intent(in) :: pzthe(ngrid) |
---|
| 41 | |
---|
| 42 | real :: tab_cntrl(length) ! nb "length=100" defined in iostart module |
---|
| 43 | |
---|
| 44 | ! Create physics start file |
---|
| 45 | call open_restartphy(filename) |
---|
[135] | 46 | |
---|
[1216] | 47 | ! tab_cntrl() contains run parameters |
---|
| 48 | tab_cntrl(:)=0 ! initialization |
---|
| 49 | !cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 50 | ! Fill control array tab_cntrl(:) with paramleters for this run |
---|
| 51 | !cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
| 52 | ! Informations on the physics grid |
---|
| 53 | tab_cntrl(1) = float(klon_glo) ! number of nodes on physics grid |
---|
| 54 | tab_cntrl(2) = float(nlay) ! number of atmospheric layers |
---|
| 55 | tab_cntrl(3) = day_ini + int(time) ! final day |
---|
| 56 | tab_cntrl(4) = time -int(time) ! final time of day |
---|
[135] | 57 | |
---|
[1216] | 58 | ! Informations about Mars, used by dynamics and physics |
---|
| 59 | tab_cntrl(5) = rad ! radius of Mars (m) ~3397200 |
---|
| 60 | tab_cntrl(6) = omeg ! rotation rate (rad.s-1) |
---|
| 61 | tab_cntrl(7) = g ! gravity (m.s-2) ~3.72 |
---|
| 62 | tab_cntrl(8) = mugaz ! Molar mass of the atmosphere (g.mol-1) ~43.49 |
---|
| 63 | tab_cntrl(9) = rcp ! = r/cp ~0.256793 (=kappa dans dynamique) |
---|
| 64 | tab_cntrl(10) = daysec ! length of a sol (s) ~88775 |
---|
[135] | 65 | |
---|
[1216] | 66 | tab_cntrl(11) = phystep ! time step in the physics |
---|
| 67 | tab_cntrl(12) = 0. |
---|
| 68 | tab_cntrl(13) = 0. |
---|
[135] | 69 | |
---|
[1216] | 70 | ! Informations about Mars, only for physics |
---|
| 71 | tab_cntrl(14) = year_day ! length of year (sols) ~668.6 |
---|
| 72 | tab_cntrl(15) = periastr ! min. star-planet distance (AU) |
---|
| 73 | tab_cntrl(16) = apoastr ! max. star-planet distance (AU) |
---|
| 74 | tab_cntrl(17) = peri_day ! date of periastron (sols since N. spring) |
---|
| 75 | tab_cntrl(18) = obliquit ! Obliquity of the planet (deg) ~23.98 |
---|
[135] | 76 | |
---|
[1216] | 77 | ! Boundary layer and turbulence |
---|
| 78 | tab_cntrl(19) = z0 ! surface roughness (m) ~0.01 |
---|
| 79 | tab_cntrl(20) = lmixmin ! mixing length ~100 |
---|
| 80 | tab_cntrl(21) = emin_turb ! minimal energy ~1.e-8 |
---|
[135] | 81 | |
---|
[1216] | 82 | ! Optical properties of polar caps and ground emissivity |
---|
| 83 | tab_cntrl(22) = albedice(1) ! Albedo of northern cap ~0.5 |
---|
| 84 | tab_cntrl(23) = albedice(2) ! Albedo of southern cap ~0.5 |
---|
| 85 | tab_cntrl(24) = emisice(1) ! Emissivity of northern cap ~0.95 |
---|
| 86 | tab_cntrl(25) = emisice(2) ! Emissivity of southern cap ~0.95 |
---|
| 87 | tab_cntrl(26) = emissiv ! Emissivity of martian soil ~.95 |
---|
| 88 | tab_cntrl(31) = iceradius(1) ! mean scat radius of CO2 snow (north) |
---|
| 89 | tab_cntrl(32) = iceradius(2) ! mean scat radius of CO2 snow (south) |
---|
| 90 | tab_cntrl(33) = dtemisice(1) ! time scale for snow metamorphism (north) |
---|
| 91 | tab_cntrl(34) = dtemisice(2) ! time scale for snow metamorphism (south) |
---|
[135] | 92 | |
---|
[1216] | 93 | tab_cntrl(28) = 0. |
---|
| 94 | tab_cntrl(29) = 0. |
---|
| 95 | tab_cntrl(30) = 0. |
---|
[588] | 96 | |
---|
[135] | 97 | ! Soil properties: |
---|
[1216] | 98 | tab_cntrl(35) = volcapa ! soil volumetric heat capacity |
---|
[135] | 99 | |
---|
[1216] | 100 | call put_var("controle","Control parameters",tab_cntrl) |
---|
| 101 | |
---|
| 102 | ! Write the mid-layer depths |
---|
| 103 | call put_var("soildepth","Soil mid-layer depth",mlayer) |
---|
| 104 | |
---|
| 105 | ! Write longitudes |
---|
| 106 | call put_field("longitude","Longitudes of physics grid",lonfi) |
---|
| 107 | |
---|
| 108 | ! Write latitudes |
---|
| 109 | call put_field("latitude","Latitudes of physics grid",latfi) |
---|
| 110 | |
---|
| 111 | ! Write mesh areas |
---|
| 112 | call put_field("area","Mesh area",area) |
---|
| 113 | |
---|
| 114 | ! Write surface geopotential |
---|
| 115 | call put_field("phisfi","Geopotential at the surface",phisfi) |
---|
| 116 | |
---|
| 117 | ! Write surface albedo |
---|
[1222] | 118 | call put_field("albedodat","Albedo of bare ground",alb) |
---|
[1216] | 119 | |
---|
| 120 | ! Subgrid topogaphy variables |
---|
| 121 | call put_field("ZMEA","Relief: mean relief",zmea) |
---|
| 122 | call put_field("ZSTD","Relief: standard deviation",zstd) |
---|
| 123 | call put_field("ZSIG","Relief: sigma parameter",zsig) |
---|
| 124 | call put_field("ZGAM","Relief: gamma parameter",zgam) |
---|
| 125 | call put_field("ZTHE","Relief: theta parameter",zthe) |
---|
| 126 | |
---|
| 127 | ! Soil thermal inertia |
---|
[1222] | 128 | call put_field("inertiedat","Soil thermal inertia",ith) |
---|
[1216] | 129 | |
---|
| 130 | ! Close file |
---|
| 131 | call close_restartphy |
---|
| 132 | |
---|
| 133 | end subroutine physdem0 |
---|
[135] | 134 | |
---|
[1216] | 135 | subroutine physdem1(filename,nsoil,ngrid,nlay,nq, & |
---|
| 136 | phystep,time,tsurf,tsoil,emis,q2,qsurf, & |
---|
| 137 | cloudfrac,totcloudfrac,hice) |
---|
| 138 | ! write time-dependent variable to restart file |
---|
| 139 | use iostart, only : open_restartphy, close_restartphy, & |
---|
| 140 | put_var, put_field |
---|
| 141 | use infotrac, only: tname |
---|
| 142 | implicit none |
---|
| 143 | !====================================================================== |
---|
| 144 | !#include "temps.h" |
---|
| 145 | #include "comcstfi.h" |
---|
| 146 | #include "planete.h" |
---|
| 147 | !====================================================================== |
---|
| 148 | character(len=*),intent(in) :: filename |
---|
| 149 | integer,intent(in) :: nsoil |
---|
| 150 | integer,intent(in) :: ngrid |
---|
| 151 | integer,intent(in) :: nlay |
---|
| 152 | integer,intent(in) :: nq |
---|
| 153 | real,intent(in) :: phystep |
---|
| 154 | real,intent(in) :: time |
---|
| 155 | real,intent(in) :: tsurf(ngrid) |
---|
| 156 | real,intent(in) :: tsoil(ngrid,nsoil) |
---|
| 157 | real,intent(in) :: emis(ngrid) |
---|
| 158 | real,intent(in) :: q2(ngrid,nlay+1) |
---|
| 159 | real,intent(in) :: qsurf(ngrid,nq) |
---|
| 160 | real,intent(in) :: cloudfrac(ngrid,nlay) |
---|
| 161 | real,intent(in) :: totcloudfrac(ngrid) |
---|
| 162 | real,intent(in) :: hice(ngrid) |
---|
[135] | 163 | |
---|
[1216] | 164 | integer :: iq |
---|
| 165 | |
---|
| 166 | ! Open file |
---|
| 167 | call open_restartphy(filename) |
---|
[135] | 168 | |
---|
[1216] | 169 | ! First variable to write must be "Time", in order to correctly |
---|
| 170 | ! set time counter in file |
---|
| 171 | !call put_var("Time","Temps de simulation",time) |
---|
| 172 | |
---|
| 173 | ! Surface temperature |
---|
| 174 | call put_field("tsurf","Surface temperature",tsurf) |
---|
| 175 | |
---|
| 176 | ! Soil temperature |
---|
| 177 | call put_field("tsoil","Soil temperature",tsoil) |
---|
| 178 | |
---|
| 179 | ! Emissivity of the surface |
---|
| 180 | call put_field("emis","Surface emissivity",emis) |
---|
| 181 | |
---|
| 182 | ! Planetary Boundary Layer |
---|
| 183 | call put_field("q2","pbl wind variance",q2) |
---|
| 184 | |
---|
| 185 | ! cloud fraction and sea ice (NB: these should be optional... to be improved) |
---|
| 186 | call put_field("cloudfrac","Cloud fraction",cloudfrac) |
---|
| 187 | call put_field("totcloudfrac","Total fraction",totcloudfrac) |
---|
| 188 | call put_field("hice","Height of oceanic ice",hice) |
---|
[135] | 189 | |
---|
[1216] | 190 | ! tracers |
---|
| 191 | if (nq>0) then |
---|
| 192 | do iq=1,nq |
---|
| 193 | call put_field(tname(iq),"tracer on surface",qsurf(:,iq)) |
---|
| 194 | enddo |
---|
| 195 | endif ! of if (nq>0) |
---|
[135] | 196 | |
---|
[1216] | 197 | ! close file |
---|
| 198 | CALL close_restartphy |
---|
| 199 | !$OMP BARRIER |
---|
[135] | 200 | |
---|
[1216] | 201 | end subroutine physdem1 |
---|
[135] | 202 | |
---|
[1216] | 203 | end module phyredem |
---|