Ignore:
Timestamp:
Jan 18, 2022, 3:29:00 PM (3 years ago)
Author:
romain.vande
Message:

LMDZ_MARS RV : Open_MP; Reading files in parallel for albedocaps parametrisation
New subroutine specific for the reading (TESicealbedo = .true.)

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/libf/phymars/albedocaps.F90

    r2508 r2609  
    99                     emisice, albedice, watercaptag, albedo_h2o_cap, &
    1010                     emissiv, albedodat
     11USE mod_phys_lmdz_transfert_para, ONLY: bcast
     12USE mod_phys_lmdz_para, ONLY: is_master
    1113implicit none
    1214
     
    2325! local variables:
    2426logical,save :: firstcall=.true.
     27!$OMP THREADPRIVATE(firstcall)
    2528integer :: ig,icap
     29
     30!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
     31
     32! local variables:
     33real,save :: zls_old ! value of zls from a previous call
     34real,save :: pi,radeg ! to convert radians to degrees
     35
     36!$OMP THREADPRIVATE(zls_old,pi,radeg)
     37
     38! TES datasets: (hard coded fixed length/sizes; for now)
     39integer,parameter :: TESlonsize=72
     40! longitudes, in TES files, in degrees, from TESlon(1)=-177.5 to TESlon(72)=177.5
     41real,save :: TESlon(TESlonsize)
     42integer,parameter :: TESlatsize=30
     43! latitudes (north hemisphere file), in degrees, from TESlatn(1)=31,
     44! to TESlatn(30)=89 ; TESlatn(8)=45
     45real,save :: TESlatn(TESlatsize)
     46! latitudes (south hemisphere file), in degrees, from TESlats(1)=-89,
     47! to TESlats(30)=-31 ; TESlats(23)=-45
     48real,save :: TESlats(TESlatsize)
     49integer,parameter :: TESlssize=72
     50! Solar longitude in TES files, TESls(1)=2.5 to TESls(72)=357.5
     51real,save :: TESls(TESlssize)
     52! TES North albedo (=-1 for missing values)
     53real,save :: TESalbn(TESlonsize,TESlatsize,TESlssize)
     54! TES South albedo (=-1 for missing values)
     55real,save :: TESalbs(TESlonsize,TESlatsize,TESlssize)
     56
     57!$OMP THREADPRIVATE(TESlon,TESlatn,TESlats,TESls,TESalbn,TESalbs)
     58
     59
     60!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
    2661
    2762! 1. Initializations
     
    4782  endif
    4883 
     84  call read_TES_icecap_albedo(zls_old,pi,radeg,TESlon,TESlatn,TESlats,TESls,TESalbn,TESalbs)
     85 
    4986  firstcall=.false.
    5087endif ! of if (firstcall)
    51 
     88   
    5289do ig=1,ngrid
    5390  if (latitude(ig).lt.0.) then
     
    5693    icap=1 ! Northern hemisphere
    5794  endif
    58 
     95 
    5996  if (piceco2(ig).gt.0) then
    6097    ! set emissivity of surface to be the ice emissivity
     
    6299    ! set the surface albedo to be the ice albedo
    63100    if (TESicealbedo) then
    64       call TES_icecap_albedo(zls,ig,psolaralb(ig,1),icap)
     101      call TES_icecap_albedo(zls,ig,psolaralb(ig,1),icap,zls_old,pi,radeg,TESlon,TESlatn,TESlats,TESls,TESalbn,TESalbs)
    65102      psolaralb(ig,2)=psolaralb(ig,1)
    66103    else
     
    82119  endif ! of if (piceco2(ig).gt.0)
    83120enddo ! of ig=1,ngrid
     121
    84122end subroutine albedocaps
    85123
    86124!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
    87 subroutine TES_icecap_albedo(zls,ig,alb,icap)
    88 
    89 use geometry_mod, only: latitude, longitude ! in radians
    90 use surfdat_h, only: albedice, TESice_Ncoef, TESice_Scoef
     125
     126subroutine read_TES_icecap_albedo(zls_old,pi,radeg,TESlon,TESlatn,TESlats,TESls,TESalbn,TESalbs)
     127
    91128use datafile_mod, only: datadir
    92129use netcdf, only: nf90_open, NF90_NOWRITE, NF90_NOERR, &
    93130                  nf90_strerror, nf90_inq_varid, nf90_get_var, nf90_close
     131USE mod_phys_lmdz_para, ONLY: is_master
     132USE mod_phys_lmdz_transfert_para, ONLY: bcast
    94133                 
    95134implicit none
    96135
    97136! arguments:
    98 real,intent(in) :: zls ! solar longitude (rad)
    99 integer,intent(in) :: ig ! grid point index
    100 real,intent(out) :: alb ! (interpolated) TES ice albedo at that grid point
    101 integer :: icap ! =1: Northern hemisphere =2: Southern hemisphere
    102137
    103138! local variables:
    104 logical,save :: firstcall=.true.
    105 real,save :: zls_old ! value of zls from a previous call
    106 integer,save :: tinf,tsup ! encompassing time indexes of TES data
    107 real,save :: reltime ! relative position in-between time indexes (in [0;1])
    108 integer :: latinf,latsup ! encompassing latitude indexes of TES data
    109 real :: rellat ! relative position in-between latitude indexes (in[0;1])
    110 integer :: loninf,lonsup ! encompassing longitude indexes of TES data
    111 real :: rellon !relative position in-between longitude indexes (in[0;1])
    112 real,save :: pi,radeg ! to convert radians to degrees
    113 real :: zlsd ! solar longitude, in degrees
    114 real :: latd ! latitude, in degrees
    115 real :: lond ! longitude, in degrees
    116 integer :: i
     139real:: zls_old ! value of zls from a previous call
     140real:: pi,radeg ! to convert radians to degrees
    117141character(len=20),parameter :: modname="TES_icecap_albedo"
    118142
    119143! TES datasets: (hard coded fixed length/sizes; for now)
    120 integer,parameter :: TESlonsize=72
    121 real,parameter :: TESdeltalon=5.0 ! step in longitude in TES files
     144integer,parameter:: TESlonsize=72
    122145! longitudes, in TES files, in degrees, from TESlon(1)=-177.5 to TESlon(72)=177.5
    123 real,save :: TESlon(TESlonsize)
    124 integer,parameter :: TESlatsize=30
    125 real,parameter :: TESdeltalat=2.0 ! step in latitude in TES files
     146real:: TESlon(TESlonsize)
     147integer,parameter:: TESlatsize=30
    126148! latitudes (north hemisphere file), in degrees, from TESlatn(1)=31,
    127149! to TESlatn(30)=89 ; TESlatn(8)=45
    128 real,parameter :: TESlatnmin=45. ! minimum TES latitude (North hemisphere)
    129 real,parameter :: TESlatsmax=-45. ! maximum TES latitude (South hemisphere)
    130 real,save :: TESlatn(TESlatsize)
     150real:: TESlatn(TESlatsize)
    131151! latitudes (south hemisphere file), in degrees, from TESlats(1)=-89,
    132152! to TESlats(30)=-31 ; TESlats(23)=-45
    133 real,save :: TESlats(TESlatsize)
    134 integer,parameter :: TESlssize=72
    135 real,parameter :: TESdeltals=5.0 ! step in solar longitude in TES files
     153real:: TESlats(TESlatsize)
     154integer,parameter:: TESlssize=72
    136155! Solar longitude in TES files, TESls(1)=2.5 to TESls(72)=357.5
    137 real,save :: TESls(TESlssize)
     156real:: TESls(TESlssize)
    138157! TES North albedo (=-1 for missing values)
    139 real,save :: TESalbn(TESlonsize,TESlatsize,TESlssize)
     158real:: TESalbn(TESlonsize,TESlatsize,TESlssize)
    140159! TES South albedo (=-1 for missing values)
    141 real,save :: TESalbs(TESlonsize,TESlatsize,TESlssize)
    142 ! encompassing nodes arranged as follow : 4 3
    143 real :: val(4)                          ! 1 2
     160real:: TESalbs(TESlonsize,TESlatsize,TESlssize)
    144161
    145162!NetCDF variables:
     
    149166
    150167! 0. Preliminary stuff
    151 if (firstcall) then
     168
     169if(is_master) then
     170
    152171! Load TES albedoes for Northern Hemisphere
    153172  ierr=nf90_open(trim(datadir)//"/npsc_albedo.nc",NF90_NOWRITE,nid)
     
    273292
    274293  zls_old=-999 ! dummy initialization
    275 
    276   firstcall=.false.
    277 endif ! of if firstcall
     294   
     295endif !is_master
     296   
     297    call bcast(TESlon)
     298    call bcast(TESlatn)
     299    call bcast(TESlats)
     300    call bcast(TESls)
     301    call bcast(TESalbn)
     302    call bcast(TESalbs)
     303    call bcast(pi)
     304    call bcast(zls_old)
     305    call bcast(radeg)
     306 
     307end subroutine read_TES_icecap_albedo
     308
     309subroutine TES_icecap_albedo(zls,ig,alb,icap,zls_old,pi,radeg,TESlon,TESlatn,TESlats,TESls,TESalbn,TESalbs)
     310
     311use geometry_mod, only: latitude, longitude ! in radians
     312use surfdat_h, only: albedice, TESice_Ncoef, TESice_Scoef
     313use datafile_mod, only: datadir
     314use netcdf, only: nf90_open, NF90_NOWRITE, NF90_NOERR, &
     315                  nf90_strerror, nf90_inq_varid, nf90_get_var, nf90_close
     316USE mod_phys_lmdz_transfert_para, ONLY: bcast
     317                 
     318implicit none
     319
     320! arguments:
     321real,intent(in) :: zls ! solar longitude (rad)
     322integer,intent(in) :: ig ! grid point index
     323real,intent(out) :: alb ! (interpolated) TES ice albedo at that grid point
     324integer :: icap ! =1: Northern hemisphere =2: Southern hemisphere
     325
     326! local variables:
     327integer,save :: tinf,tsup ! encompassing time indexes of TES data
     328real,save :: reltime ! relative position in-between time indexes (in [0;1])
     329integer :: latinf,latsup ! encompassing latitude indexes of TES data
     330real :: rellat ! relative position in-between latitude indexes (in[0;1])
     331integer :: loninf,lonsup ! encompassing longitude indexes of TES data
     332real :: rellon !relative position in-between longitude indexes (in[0;1])
     333real :: zlsd ! solar longitude, in degrees
     334real :: latd ! latitude, in degrees
     335real :: lond ! longitude, in degrees
     336integer :: i
     337
     338!$OMP THREADPRIVATE(tinf,tsup,reltime)
     339
     340! TES datasets: (hard coded fixed length/sizes; for now)
     341real,parameter :: TESdeltalon=5.0 ! step in longitude in TES files
     342real,parameter :: TESdeltalat=2.0 ! step in latitude in TES files
     343real,parameter :: TESlatnmin=45. ! minimum TES latitude (North hemisphere)
     344real,parameter :: TESlatsmax=-45. ! maximum TES latitude (South hemisphere)
     345! to TESlats(30)=-31 ; TESlats(23)=-45
     346real,parameter :: TESdeltals=5.0 ! step in solar longitude in TES files
     347
     348real:: zls_old ! value of zls from a previous call
     349real:: pi,radeg ! to convert radians to degrees
     350
     351! TES datasets: (hard coded fixed length/sizes; for now)
     352integer,parameter:: TESlonsize=72
     353! longitudes, in TES files, in degrees, from TESlon(1)=-177.5 to TESlon(72)=177.5
     354real:: TESlon(TESlonsize)
     355integer,parameter:: TESlatsize=30
     356! latitudes (north hemisphere file), in degrees, from TESlatn(1)=31,
     357! to TESlatn(30)=89 ; TESlatn(8)=45
     358real:: TESlatn(TESlatsize)
     359! latitudes (south hemisphere file), in degrees, from TESlats(1)=-89,
     360! to TESlats(30)=-31 ; TESlats(23)=-45
     361real:: TESlats(TESlatsize)
     362integer,parameter:: TESlssize=72
     363! Solar longitude in TES files, TESls(1)=2.5 to TESls(72)=357.5
     364real:: TESls(TESlssize)
     365! TES North albedo (=-1 for missing values)
     366real:: TESalbn(TESlonsize,TESlatsize,TESlssize)
     367! TES South albedo (=-1 for missing values)
     368real:: TESalbs(TESlonsize,TESlatsize,TESlssize)
     369! encompassing nodes arranged as follow : 4 3
     370real :: val(4)                          ! 1 2
    278371
    279372! 1. Identify encompassing latitudes
     
    281374! Check that latitude is such that there is TES data to use
    282375! (ie: latitude 45 deg and poleward) otherwise use 'default' albedoes
     376
    283377latd=latitude(ig)*radeg ! latitude, in degrees
    284378if (icap.eq.1) then
Note: See TracChangeset for help on using the changeset viewer.