Interpolation tools¶
Netcdf_add_massif_info script¶
Script interpolation/Netcdf_add_massif_info.py
help:
usage: Netcdf_add_massif_info.py [-h] [-s PATH_SHAPE] [-p] [-v]
[-o PATH_OUTPUT]
path_forcing
Script to add massif info into the netcdf forcing 2D with a variable ZS.
Dimension of ZS = (latitude, longitude). Example of use: python3
Netcdf_add_massif_info /home/reveilletm/alpha4.nc -p python3
Netcdf_add_massif_info.py input_4_test_add_massif_info.nc -p -v -o
output_4_test_add_massif_info.nc
positional arguments:
path_forcing Path to forcing file (mandatory)
options:
-h, --help show this help message and exit
-s PATH_SHAPE, --shape PATH_SHAPE
Path to shapefile file
-p, --plot Make a plot to check
-v, --verbose Print index in loop
-o PATH_OUTPUT, --output PATH_OUTPUT
Path to output file
shapefile2NetCDF tool¶
This script needs a shapefile (in Lambert 93 coordinates) of points. It is useful for simulations on several differents points.
If you want a simulation on a 2D domain delimited by a shapefile: see shapefile2NetCDF_2D.py instead.
General documentation of shapefile2NetCDF script¶
This script creates a NetCDF with the same points than in shapefile. The NetCDF file is necessary to launch Surfex-Crocus simulations. It allows to interpolate FORCING values from SAFRAN reanalysis to your shapefile points. It also allows to take into account the solar mask from mountains around.
Precisely, this script writes skylines in snowtools/DATA/METADATA.xml file. These skylines are necessary for the –addmask option in SURFEX Crocus simulations. This script also allows to keep the skyline views if you want to check them.
WHAT IS PROJECT NUMBER ?¶
Because we add some points in METADATA.xml file, we want to avoid getting several points with same station code During the script conception, we tought that:
each project gets less than 10000 points
there will be less than 100 projects with the necessity to keep the points (and station code) in METADATA.xml
If you don’t push METADATA.xml, you don’t need to think too much about that
The aim of Project Number is to avoid to get twice the same station code (8 numbers as OMM codes). The idea is:
numbers between 10000001 and 10009999 are for project 0
numbers between 10010001 and 10019999 are for project 1
etc
RECORD FOR PROJECT NUMBER¶
ORCHAMP = project 0 = orchamp geometry in vortex/conf/geometries.ini (164 points)
projet 1 used without references ? (15 points)
TOP_CBNA = project 2 = orchamp geometry in vortex/conf/geometries.ini (169 points)
ORCHAMP_MAJ_2022 = project 3 = orchamp geometry in vortex/conf/geometries.ini (26 points)
Reminder: if your project doesn’t need to stay a longtime, don’t push METADATA.xml
EXAMPLES OF USE (script launch)¶
python3 shapefile2NetCDF.py path_shapefile station_name_in_shapefile station_id_in_shapefile project_number
[--name_alt alti_name_in_shapefile] [-o path_name_of_NetCDF_output] [--MNT_alt path_of_MNT_altitude]
[--confirm_overwrite] [--list_skyline all or 1 5 6 if you want skyline for your id_station number 1, 5 and 6]
python3 shapefile2NetCDF.py /home/fructusm/Téléchargements/Plots2020/plots codeplot idplot 0 --name_alt alti
python3 shapefile2NetCDF.py /home/fructusm/Téléchargements/Plots2020/plots codeplot idplot 0 --name_alt alti
--confirm_overwrite (if you have already work on this project)
python3 shapefile2NetCDF.py /home/fructusm/Téléchargements/Plots2020/plots codeplot idplot 0 --name_alt alti
--list_skyline 1 34 47 (in folder output, you'll have skylines for stations 1, 34 and 47)
python3 shapefile2NetCDF.py /home/fructusm/Téléchargements/Plots2021/PlotsMaJ2021 codeplot idplot 0
--confirm_overwrite
TOP_CBNA:
python3 shapefile2NetCDF.py /home/fructusm/Bureau/Shapefile_simu/cn_maille_points/cn_maille_points cd50m ORIG_FID 1
MAJ Orchamp_2022
python3 shapefile2NetCDF.py /home/fructusm/Travail_Orchamp/Orchamp_nouveaux_sites_2022/Plots_new_2022 codeplot
idplot 3
EXAMPLES OF USE OF THE NETCDF FILE IN LOCAL SIMULATION¶
Get the interpolation of the FORCING file, interpol_FORCING:
s2m research -f path_FORCING -b begin_date -e end_date -r path_netcdf.nc -o output_name -g --extractforcing
Obtain the PRO file for the shapefile points:
s2m research -f path_interpol_FORCING -b begin_date -e end_date -o output_name_pour_PRO -g --addmask
Example:
s2m research -f /rd/cenfic3/cenmod/era40/vortex/s2m/alp_flat/reanalysis/meteo/FORCING_2017080106_2018080106.nc
-b 20170801 -e 20180801 -r /home/fructusm/git/snowtools_git/interpolation/NetCDF_from_shapefile.nc
-o output_test_s2m -g --extractforcing
s2m research -f /home/fructusm/OUTPUT_et_PRO/output_test_s2m/meteo/FORCING_2017080106_2018080106.nc
-b 20170801 -e 20180801 -o output_test_s2m_Safran -g --addmask
EXAMPLES OF USE OF THE NETCDF FILE IN HPC SIMULATION¶
s2m research -b 19580801 -e 20200801 -m s2m -f reanalysis2020.2@lafaysse
-r alp_flat:orchamp:/home/cnrm_other/cen/mrns/fructusm/NetCDF_from_shapefile.nc
-o TEST1 -n OPTIONS_V8.1_NEW_OUTPUTS_NC_reanalysis.nam -g --addmask -a 400
Options:
-m model
-f forcing files -> -f reanalysis in order to get the forcing from reanalysis
-r région: add geometry in vortex/conf/geometries.ini
-n namelist (get the same options than reanalysis)
-g if you don’t have a prep -> a spinup has to be made
–addmask in order to use the masks created at the same time than the NetCDF file
-a 400 In order to limit Snow Water Equivalent to 400kg/m3 at the end of the year (1rst of august)
–ntasks 8 if you have only 8 points -> otherwise crash
HPC SIMULATION: CLIMATE SIMULATIONS¶
Forcing files from SMHI_RCA4_MOHC_HadGEM2_ES_RCP85 (for example) are in Raphaelle Samacoits folders
s2m research -b 20100801 -e 20990801 -m adamont -f CLMcom_CCLM4_8_17_CNRM_CERFACS_CNRM_CM5_RCP45@samacoitsr
-r alp_flat:orchamp:/home/cnrm_other/cen/mrns/fructusm/NetCDF_from_shapefile.nc -o TEST_Raphaelle
-n OPTIONS_V8.1_NEW_OUTPUTS_NC_reanalysis.nam -x 20200801 --addmask -a 400
HPC SIMULATION: SXCEN TRANSFERT¶
Transfer files from Belenos to sxcen: Use get_reanalysis which is in snowtools/scripts/extract/vortex
On SXCEN:
python3 get_reanalysis.py --geometry=orchamp --snow --xpid=SPINUP@fructusm
python3 get_reanalysis.py --geometry=orchamp --snow --xpid=TEST_Raphaelle@fructusm
--byear=2010 --eyear=2099
DO NOT FORGET¶
SPINUP
SEND METADATA.xml ON BELENOS IN ORDER TO USE THE MASKS
CHANGE THE NAMELIST: in CSELECT, you can’t have station AND massif_num AT THE SAME TIME
REMOVE IN NAMELIST AVALANCHES DIAG IN CSELECT (in that case, massif_num is necessary)
In Surfex-Crocus code (normally OK now but in order to remember): /SURFEX/src/OFFLINE init_output_oln.F90 (line 130) and ol_write_coord.F90
- interpolation.shapefile2NetCDF.check_id_station_in_Metadata(all_lists)[source]¶
Check that id from stations (creation via make_dict_list) are not in METADATA.xml. Skip this verification if confirm_overwrite is set to True.
- Parameters:
all_lists – Dictionary of lists in order to get id_station list
- Returns:
nothing (in a pythonic way) Screen printings if id stations in METADATA.xml. Othewise: nothing happens
- interpolation.shapefile2NetCDF.create_NetCDF(all_lists, output_name)[source]¶
Creation of 1D NetCDF for simulation reanalysis-projection from dictionary of lists and output name. Dictionary of lists is coming from the shapefile and is built with make_dict_list.
- Parameters:
all_lists – Dictionary of lists in order to fill NetCDF variables
output_name (str) – output name for the NetCDF file
- Returns:
nothing (in a pythonic way). Write a netcdf file.
- interpolation.shapefile2NetCDF.create_skyline(all_lists, path_MNT_alt, path_shapefile, list_skyline)[source]¶
Add in METADATA.xml the points from shapefile with skylines masks (ie the angles of view for the azimuts).
- Parameters:
all_lists – Dictionary of lists
path_MNT_alt (str) – MNT path for altitude values
path_shapefile – shapefile path
path_shapefile – str
list_skyline – List of identification points of shapefile from which we want the skyline.
- Returns:
nothing (in a pythonic way). Plot the skyline and add lines in METADATA.xml
- interpolation.shapefile2NetCDF.main(args=None)[source]¶
Main program: parse argument then launch the creation of NetCDF and skylines
- interpolation.shapefile2NetCDF.make_dict_list(path_shapefile, id_station, nom_station, nom_alt, nom_asp, nom_slop, path_MNT_alt, path_MNT_asp, path_MNT_slop, add_for_METADATA)[source]¶
Create a dictionary of lists from datas in shapefile and data in MNT
- Parameters:
path_shapefile (str) – shapefile path
id_station (str) – name of parameter in shapefile for point number (each point has a unique number)
nom_station (str) – name of parameter in shapefile for point name (each point has a unique name)
nom_alt (str) – if present in shapefile, name of altitude parameter. If not present, value is None.
nom_asp (str) – if present in shapefile, name of aspect parameter. If not present, value is None.
nom_slop (str) – if present in shapefile, name of slope parameter. If not present, value is None.
path_MNT_alt (str) – MNT path for altitude
path_MNT_asp (str) – MNT path for aspect
path_MNT_slop (str) – MNT path for slopes
add_for_METADATA (int) – integer (8 figures) for station code in METADATA.xml file
- Returns:
dictionary of lists:{ ‘lat’: list_latitude, ‘lon’:list_longitude, ‘alt’: list_altitude, ‘asp’: list_aspect, ‘m’: list_massif, ‘slop’: list_slope, ‘id’: list_id_station, ‘nom’: list_nom_station }
- interpolation.shapefile2NetCDF.parseArguments(args)[source]¶
Parsing the arguments when you call the main program. :param args: The list of arguments when you call the main program (typically sys.argv[1:] )
This module is also an executable script:
Script interpolation/shapefile2NetCDF.py
help:
usage: shapefile2NetCDF.py [-h] [--name_alt NAME_ALT] [--name_asp NAME_ASP]
[--name_slop NAME_SLOP] [-o OUTPUT]
[--MNT_alt MNT_ALT] [--MNT_asp MNT_ASP]
[--MNT_slop MNT_SLOP]
[--list_skyline [LIST_SKYLINE ...]]
[--confirm_overwrite]
path_shape name_station id_station project_number
positional arguments:
path_shape Path to shapefile
name_station Shapefile Field Name containing a unique reference for
points
id_station Shapefile Field Number containing a unique reference
for points
project_number Project Number in order to get unique reference in
METADATA
options:
-h, --help show this help message and exit
--name_alt NAME_ALT Shapefile Field Name containing altitude, if it exists
--name_asp NAME_ASP Shapefile Field Name containing aspect, if it exists
--name_slop NAME_SLOP
Shapefile Field Name containing slope, if it exists
-o OUTPUT, --output OUTPUT
Name for NetCDF file to save
--MNT_alt MNT_ALT Path for MNT altitude
--MNT_asp MNT_ASP Path for MNT altitude
--MNT_slop MNT_SLOP Path for MNT altitude
--list_skyline [LIST_SKYLINE ...]
The skyline plot you want
--confirm_overwrite Confirm you want to overwrite
shapefile2NetCDF_2D tool¶
This script needs a shapefile (Lambert 93 or (lat-lon) coordinates). This shapefile can be in .shp or .kml format. This shapefile is delimitating a 2D domain on which you want to make a SURFEX-Crocus simulation.
If you want to do a simulation for several distincts differents points and not a domain: see shapefile2NetCDF.py instead.
General documentation of shapefile2NetCDF_2D script¶
This script creates a NetCDF file using the shapefile. The NetCDF file is necessary to launch Surfex-Crocus 2D simulations. The NetCDF will define a rectangular 2D domain which contains your shapefile. If your 2D domain is belonging to several “SAFRAN massif”, by default: If more than 50% of points are in one massif, all the points are artificially associated with this massif After the execution you have some infos for the NAM_IGN part of the namelist.
Options :
* -o Path to shapefile (with extension .kml or .shp)
* -rlon Longitude Resolution for 2D grid (250 or 30)
* -rlat Latitude Resolution for 2D grid (250 or 30)
* --MNT_alt Path for MNT altitude
* -m Massif number if you want to choose the massif that will be applied to your zone
(Massif number can be found in snowtools/DATA/METADATA.xml)
EXAMPLES OF USE (script launch)¶
python3 shapefile2NetCDF_2D.py SP_Abries_Ristolas.kml -o NetCDF2D_Ristolas.nc
EXAMPLES OF USE OF NETCDF FILE IN HPC SIMULATION¶
First, push on Belenos with scp:
scp NetCDF2D_* fructusm@belenos:/home/fructusm/SIMU_TOP2D/.
Then s2m command:
s2m research --walltime='2:00:00' -b 20190801 -e 20200801 -m s2m -f reanalysis2020.2@lafaysse
-r alp_flat:Aussois250:/home/cnrm_other/cen/mrns/fructusm/SIMU_TOP2D/NetCDF2D_Aussois.nc
-n /home/cnrm_other/cen/mrns/fructusm/SIMU_TOP2D/Aussois.nam --geotype grid --ntasks=52
-g -o Aussois_2D
Warning
DO NOT FORGET THE SPINUP
Options:
* -m model
* -f forcing files -> -f reanalysis in order to get the forcing from reanalysis
* -r région: add geometry in vortex/conf/geometries.ini
* -n namelist (get the same options than reanalysis)
* -g if you don't have a prep -> a spinup has to be made
* -a 400 In order to limit Snow Water Equivalent to 400kg/m3 at the end of the year (1rst of august)
* --geotype grid Type of simulation geometry: grids must be defined in namelist while namelists are automatically
updated for unstructured geometries according to the forcing file. (default:
unstructured) (from s2m research -h)
* --ntasks 50 if 50 < min(column, lines) of netcdf 2D domain
File transfert from Belenos to sxcen: use get_reanalysis which is in snowtools/scripts/extract/vortex
On SXCEN:
python3 get_reanalysis.py --geometry=orchamp --snow --xpid=TEST_Raphaelle@fructusm
-byear=2010 --eyear=2099
- interpolation.shapefile2NetCDF_2D.clean_the_mess(path_shape_shp, clean_all)[source]¶
Remove all the temporary files created with the script
- Parameters:
path_shape_shp (str) – le chemin vers le fichier de shapefile
- Returns:
In a pythonic point of view, nothing. Some files are deleted.
- interpolation.shapefile2NetCDF_2D.conversion_to_L93_if_lat_lon(list_of_four)[source]¶
Idea: we have a list bbox [lon_min, lat_min, lon_max, lat_max]. We convert it in Lambert93 By default, we raise a little in case of Criteria: if in bbox, this is not lower than 200, we are already in L93 Because coordinates are in WGS84 (ie lat-lon -> everything below 200), or in L93
REMINDER: L93 = EPSG 2154 WGS84 = EPSG 4326
- Parameters:
list_of_four (list) – list of 4 floats ([lon_min, lat_min, lon_max, lat_max])
- Returns:
a list of 4 values: the same in L93 coordinates
- interpolation.shapefile2NetCDF_2D.convert_kml2shp(path_shape_kml)[source]¶
kml files (Keyhole Markup Language, google format) are transform in .shp (format GIS Geographic Information System from ESRI) in order to avoid imports from “non standard” librairies in CEN actual configuration of PC.
- Parameters:
path_shape_kml (str) – path to kml shapefile
- Returns:
a string: the path to the .shp file after conversion from the .kml file.
- interpolation.shapefile2NetCDF_2D.create_MNT_tif(dict_info, path_MNT_given=None)[source]¶
Use the dictionary with all the information (coordinates of the lower-left and upper-right points, number of pixel in x and y direction, resolution in x and y direction) and and a MNT (with the good resolution) to create a .tif file with altitude values on the grid define by the information.
- Parameters:
dict_info (dict) – dictionary with all the information for the grid
path_MNT_given (str) – optional, a path for the MNT (without a path, the standard MNT in cenfic3 is used)
- Returns:
In a pythonic point of view, nothing. A .tif file is created.
- interpolation.shapefile2NetCDF_2D.create_dict_all_infos(list_of_four, resol_x=250, resol_y=250)[source]¶
Create a dictionnary with all the important information (coordinates of the lower-left and upper-right points, number of pixel in x and y direction, resolution in x and y direction).
- Parameters:
list_of_four (list) – list of four coordinates in L93 format
resol_x (int) – resolution in x direction
resol_y (int) – resolution in y direction
- Returns:
a dictionary with all the information
- interpolation.shapefile2NetCDF_2D.create_netcdf(massif_number, file_out='NetCDF2D_from_shapefile.nc')[source]¶
Creation of the 2D_NetCDF with 2 values: ‘ZS’ for the altitude and ‘Massif number’ for the massif number.
- Parameters:
massif_number – an integer or a numpy array depending of the percentage of points in the most common massif
file_out (str) – the name of the output for the netcdf file
- Returns:
In a pythonic point of view, nothing. A .nc file is created.
- interpolation.shapefile2NetCDF_2D.find_most_common_massif(dict_info)[source]¶
Use the grid info coming from the dictionary and the massif info coming from snowtools to detect what is the most common massif on the grid. If there are more than 50% of points in the most common massif, all the points are artificially associated with the most common massif. Otherwise, all the points in the grid keep the massif they belong.
- Parameters:
dict_info (dict) – dictionary with all the information for the grid
- Returns:
an integer (the massif number of the most common massif) or a numpy array (if there is not common massif).
- interpolation.shapefile2NetCDF_2D.get_bounds(path_shape_shp)[source]¶
sf.bbox return [6.61446151100006, 45.189366518, 6.8154706149867, 45.2979826160001] ie [lon_min, lat_min, lon_max, lat_max]
- Parameters:
path_shape_shp (str) – path to shapefile
- Returns:
a list of 4 values [lon_min, lat_min, lon_max, lat_max]
- interpolation.shapefile2NetCDF_2D.main(args=None)[source]¶
Main program: parse argument then launch the creation of NetCDF
- interpolation.shapefile2NetCDF_2D.parseArguments(args)[source]¶
Parsing the arguments when you call the main program. :param args: The list of arguments when you call the main program (typically sys.argv[1:] )
- interpolation.shapefile2NetCDF_2D.print_for_namelist(dict_info, output_name)[source]¶
Just to have a print on the screen in order to have all the relevant information for the namelist used later in the Surfex simulations.
- Parameters:
dict_info (dict) – dictionary with all the information for the grid
output_name (str) – the name of the output for the netcdf file
- Returns:
In a pythonic point of view, nothing. Some information are printed on the screen.
This module is also an executable script:
Script interpolation/shapefile2NetCDF_2D.py
help:
usage: shapefile2NetCDF_2D.py [-h] [-o OUTPUT] [-rlon RESOL_LON]
[-rlat RESOL_LAT] [--MNT_alt MNT_ALT]
[-m MASSIF_NUMBER]
path_shape
positional arguments:
path_shape Path to shapefile (with extension .kml or .shp)
options:
-h, --help show this help message and exit
-o OUTPUT, --output OUTPUT
Name for NetCDF file to save
-rlon RESOL_LON, --resol_lon RESOL_LON
Longitude Resolution for 2D grid (250 or 30)
-rlat RESOL_LAT, --resol_lat RESOL_LAT
Latitude Resolution for 2D grid (250 or 30)
--MNT_alt MNT_ALT Path for MNT altitude
-m MASSIF_NUMBER, --massif_number MASSIF_NUMBER
Massif number if you want to choose the massif that
will be applied to your zone
The interpolation tool¶
The interpolation tool is written in FORTRAN.
Makefiles for mageia7 PCs and for MF supercomputers are provided. Needs MPI and parallel netcdf libraries.
- program interpolate_safran¶
This program interpolates SAFRAN meteorological data.
If a namelist called interpolate_safran.nam is present, the configuration information from the namelist is used. Otherwise the default configuration is used. An example namelist can be found in sonwtools/DATA.
- Default configuration:
A file named “GRID.nc” is read, containing the output grid (2D case) or station (1D case) information.
A file named “input.nc” is supposed to contain the data in massif geometry to be interpolated on some grid or station.
A file named “ouput.nc” is written containing the interpolated version of the data.
- Possible configurations via the namelist:
Possibility to specify custom filenames (mandatory for multiinput and multioutput options).
Iterating over several input files
Providing a grid definition file for each of them (set
lmultiinput
=.TRUE. andlmultioutput
=.TRUE.) and producing a separate output file for each of themCombining the data from them on a single grid in a single output file (set
lmultiinput
=.TRUE. andlmultioutput
=.FALSE.)
Treat time steps sepatately in order to save memory in the case of very large fields by setting
ltimechunk
= .TRUE.customize the NetCDF chunksize of the spatial output dimensions in order to optimise write performance for large fields. (
lspatialchunk
= .TRUE. andnlonchunksize
,nlatchunksize
)Select variables to be interpolated. (set
lselectvar
= .TRUE., assign a list of variables tohvar_list
and give the number of selected variables tonnumber_of_variables
)
Input variables can be of any type in input.nc. In output variables, the type of variable of floating points is always single precision, integers types are not modified.
Module SUBS
Quick access
- Variables:
comm
,dim_name_out
,hfilenameg
,ierr
,ipstart
,ixstart
,iystart
,npoint_proc
,nproc
,ntime
,nx
,nx_proc
,ny
,ny_proc
,proc_id
,xundef
- Routines:
abort_interpolate()
,check()
,evaluate_aspect()
,evaluate_slope()
,explicit_slope()
,explicit_slope_lat_lon()
,fuse_fld_3d()
,indices2d()
,interpolzs2d()
,interpolzs2ddimbefore()
,interpolzs2dnotime()
,latlon_ign()
,npoint_proc_distributor()
,read_nc()
,read_nml()
,read_output_grid()
,xy_ign()
,xy_proc_distributor()
Needed modules
netcdf
modn_interpol_safran
: Module with namelist declarations.mpi
Variables
- subs/comm [integer]¶
MPI parameter
- subs/dim_name_out (*) [character(len=20),allocatable]¶
name of output dimensions
- subs/hfilenameg [character(len=100)]¶
Filename of the grid definition file
- subs/ierr [integer]¶
MPI error code
- subs/ipstart [integer]¶
start index in the spatial dimension (points) for the given processor
- subs/ixstart [integer]¶
start index in the x-dimension for the given processor
- subs/iystart [integer]¶
start index in the y-dimension for the given processor
- subs/npoint_proc [integer]¶
number of indices in the spatial dimension to be treated by the processor
- subs/nproc [integer]¶
Number of MPI processors
- subs/ntime [integer]¶
length of time dimension
- subs/nx [integer]¶
length of first spatial dimension
- subs/nx_proc [integer]¶
number of indices in the x-dimension to be treated by the processor
- subs/ny [integer]¶
length of second spatial dimension
- subs/ny_proc [integer]¶
number of indices in the y-dimension to be treated by the processor
- subs/proc_id [integer]¶
MPI rank
- subs/xundef [real,parameter=-9999999.0]¶
Fill_Value for netcdf
Subroutines and functions
- subroutine subs/abort_interpolate(ytext)¶
Crash the program
- Parameters:
ytext [character(len=*),in]
- Called from:
- subroutine subs/xy_proc_distributor(ixydim, iproc_id, inproc, iixstart, inx_proc, iiystart, iny_proc)¶
Routine used for parallelisation with MPI.
Defines indices of the two spatial dimensions to be treated by each of the processors.
- Parameters:
ixydim (2) [integer ,in] :: Dimension lengths of the spatial dimensions
iproc_id [integer ,in]
inproc [integer ,in]
iixstart [integer ,out] :: start index in the x-dimension for the given processor
inx_proc [integer ,out] :: number of indices in the x-dimension to be treated by the processor
iiystart [integer ,out] :: start index in the y-dimension for the given processor
iny_proc [integer ,out] :: number of indices in the y-dimension to be treated by the processor
- Called from:
- subroutine subs/npoint_proc_distributor(ipdim, iproc_id, inproc, iipstart, inp_proc)¶
Routine used for parallelisation with MPI.
Defines indices of the unique spatial dimension in the 1D case to be treated by each of the processors.
- Parameters:
ipdim (2) [integer ,in] :: Dimension lengths of the spatial dimensions
iproc_id [integer ,in]
inproc [integer ,in]
iipstart [integer ,out] :: start index in the spatial dimension (points) for the given processor
inp_proc [integer ,out] :: number of indices in the spatial dimension to be treated by the processor
- Called from:
- subroutine subs/read_output_grid(file_id_geo, pzsout, kmassifout, paspectout, irank, gt, ilendim, platout, plonout, pyout, pxout, pslopeout, kstationout, kmassifgather)¶
Read the netcdf file containing the output grid definition (default: GRID.nc)
- Parameters:
file_id_geo [integer ,in] :: Netcdf id of the grid file
pzsout (*,*) [real ,out,allocatable] :: Elevation of output grid
kmassifout (*,*) [integer ,out,allocatable] :: Massif number of the output points
paspectout (*,*) [real ,out,allocatable] :: aspect of output grid
irank [integer ,out] :: rank of output grid (1 for vector or 2 for matrix)
gt [character(len=2),out] :: GRID TYPE
ilendim (2) [integer ,out] :: lenghts of dimensions
platout (*) [real ,out,allocatable] :: output latitudes of 2D grid or 1D stations
plonout (*) [real ,out,allocatable] :: output longitudes of 2D grid or 1D Stations
pyout (*) [real ,out,allocatable] :: output y coordinates in the case of 2D grid with a XY grid type
pxout (*) [real ,out,allocatable] :: output x coordinates in the case of 2D grid with a XY grid type
pslopeout (*,*) [real ,out,allocatable] :: slope of output grid
kstationout (*) [integer ,out,allocatable] :: Station number of the output points (used in the 1D case only)
kmassifgather (*,*) [integer ,out,allocatable] :: Massif number of the output points gathered from all procs
- Called from:
- Call to:
check()
,abort_interpolate()
,xy_proc_distributor()
,npoint_proc_distributor()
,explicit_slope()
,xy_ign()
,explicit_slope_lat_lon()
- subroutine subs/read_nc(idficnc, kzsin, kmassifin, kaspectin, pslopein)¶
Read NetCDF file with forcing data to be interpolated.
- Parameters:
idficnc [integer ,in] :: NetCDF id of the forcing file
kzsin (*) [integer ,out,allocatable] :: Elevation of the input data points
kmassifin (*) [integer ,out,allocatable] :: Massif of the input data points
kaspectin (*) [integer ,out,allocatable] :: Aspect of the input data points
pslopein (*) [real ,out,allocatable] :: Slope of the input data points
- Called from:
- Call to:
- subroutine subs/indices2d(pzsout, kmassifout, paspectout, pslopeout, kzsin, kmassifin, kaspectin, pslopein, kindicesbas, kindiceshaut)¶
Determines the indexes of input collection of points to interpolate between in order to obtain the output grid.
Selects the hight level above and below the hight of the output grid point, at the corresponding massif number and aspect
- Parameters:
pzsout (*,*) [real ,in] :: Elevation of output collection of points
kmassifout (*,*) [integer ,in] :: Massif of output collection of points
paspectout (*,*) [real ,in] :: Aspect of output collection of points
pslopeout (*,*) [real ,in] :: Slope of output collection of points
kzsin (*) [integer ,in] :: Elevation of input collection of points
kmassifin (*) [integer ,in] :: Massif of input collection of points
kaspectin (*) [integer ,in] :: Aspect of input collection of points
pslopein (*) [real ,in] :: slope of input collection of points
kindicesbas (*,*) [integer ,out,allocatable] :: Indices of the input collection of points corresponding to the points below the target point
kindiceshaut (*,*) [integer ,out,allocatable] :: Indices of the input collection of points corresponding to the points above the target point
- Called from:
- Call to:
- subroutine subs/evaluate_aspect(kaspectin, paspectout, laspect)¶
Evaluate if the aspect of the current input point correspond to the aspect of the current output point
- Parameters:
kaspectin [integer ,in] :: Aspect of input points
paspectout [real ,in] :: Aspect of output points
laspect [logical ,out] :: Aspect correspondence indicator
- Called from:
- subroutine subs/evaluate_slope(pslopein, pslopeout, lslope)¶
Evaluate if the slope of the current input point correspond to the slope of the current output point
- Parameters:
pslopein [real ,in] :: Aspect of input points
pslopeout [real ,in] :: Aspect of output points
lslope [logical ,out] :: Aspect correspondence indicator
- Called from:
- subroutine subs/interpolzs2d(pvarout, pvarin, kindbas, kindhaut, pzsgrid, kzssafran)¶
Linear interpolation between SAFRAN elevation levels over 2D grids
- Parameters:
pvarout (*,*,*,*) [real ,out] :: Interpolated variable on 2D grid. The first 2 dimensions are the spatial dimensions, the last dimension is the time dimension.
pvarin (*,*,*) [real ,in] :: Variable on SAFRAN elevation levels
kindbas (*,*) [integer ,in] :: Index of pvarin containing the values at the SAFRAN elevation level below the grid point elevation at the corresponding massif.
kindhaut (*,*) [integer ,in] :: Index of pvarin containing the values at the SAFRAN elevation level above the grid point elevation at the corresponding massif.
pzsgrid (*,*) [real ,in] :: Gridpoint elevation data
kzssafran (*) [integer ,in] :: SAFRAN elevation data
- Called from:
- subroutine subs/interpolzs2ddimbefore(pvarout, pvarin, kindbas, kindhaut, pzsgrid, kzssafran)¶
Linear interpolation between SAFRAN elevation levels over 2D grids
- Parameters:
pvarout (*,*,*,*) [real ,out] :: Interpolated variable on 2D grid. The second and third dimension are the spatial dimensions, the last dimension is the time dimension.
pvarin (*,*,*) [real ,in] :: Variable on SAFRAN elevation levels
kindbas (*,*) [integer ,in] :: Index of pvarin containing the values at the SAFRAN elevation level below the grid point elevation at the corresponding massif.
kindhaut (*,*) [integer ,in] :: Index of pvarin containing the values at the SAFRAN elevation level above the grid point elevation at the corresponding massif.
pzsgrid (*,*) [real ,in] :: Gridpoint elevation data
kzssafran (*) [integer ,in] :: SAFRAN elevation data
- Called from:
- subroutine subs/interpolzs2dnotime(pvarout, pvarin, kindbas, kindhaut, pzsgrid, kzssafran)¶
Linear interpolation between SAFRAN elevation levels over 2D grids
- Parameters:
pvarout (*,*,*) [real ,out] :: Interpolated variable on 2D grid. The first two dimensions are the spatial dimensions. There is no time dimension.
pvarin (*,*) [real ,in] :: Variable on SAFRAN elevation levels without time dimension.
kindbas (*,*) [integer ,in] :: Index of pvarin containing the values at the SAFRAN elevation level below the grid point elevation at the corresponding massif.
kindhaut (*,*) [integer ,in] :: Index of pvarin containing the values at the SAFRAN elevation level above the grid point elevation at the corresponding massif.
pzsgrid (*,*) [real ,in] :: Gridpoint elevation data
kzssafran (*) [integer ,in] :: SAFRAN elevation data
- Called from:
- subroutine subs/fuse_fld_3d(pvarout_old, pvarout, xundef)¶
Merge two 3D fields.
- Actions :
Put values of PVAROUT into PVAROUT_OLD where PVAROUT_OLD == XUNDEF
- Parameters:
pvarout_old (*,*,*) [real ,inout] :: Variable on 2D grid with existing values to be completed.
pvarout (*,*,*) [real ,in] :: Variable on 2D grid with new values to be added.
xundef [real ,in] :: Missing value indicator.
- subroutine subs/latlon_ign(px, py[, plat[, plon]])¶
Routine to compute geographical coordinates from lambert conformal coordinates
- Parameters:
px (*) [real ,in] :: given conformal x-coordinates of the processed points (meters).
py (*) [real ,in] :: given conformal y-coordinates of the processed points (meters).
- Options:
plat (*,*) [real ,out,] :: returned geographic latitudes of the processed point (degrees).
plon (*,*) [real ,out,] :: returned geographic longitudes of the processed point (degrees).
- Called from:
- subroutine subs/xy_ign(plat, plon, px, py)¶
Routine to compute Lambert coordinates from logitude and latitude
- Parameters:
plat (*) [real ,in] :: given geographic latitudes of the processed points (degrees).
plon (*) [real ,in] :: given geographic longitudes of the processed points (degrees).
px (*,*) [real ,out] :: returned lambert conformal x-coordinates of the processed points (meters).
py (*,*) [real ,out] :: returned lambert conformal y-coordinates of the processed points (meters).
- Called from:
- subroutine subs/explicit_slope(px, py, pzs, psso_slope, psso_dir)¶
calculate resolved slope and aspect on a Lambert conformal grid.
- Parameters:
px (*) [real ,in] :: given Lambert conformal x-coordinates of the processed points (meters).
py (*) [real ,in] :: given Lambert conformal y-coordinates of the processed points (meters).
pzs (*,*) [real ,in] :: resolved model orography
psso_slope (*,*) [real ,out] :: resolved slope tangent
psso_dir (*,*) [real ,out] :: resolved aspect
- Called from:
- subroutine subs/explicit_slope_lat_lon(px, py, pzs, psso_slope, psso_dir)¶
calculate resolved slope and aspect on a lat/lon grid.
- Parameters:
px (*,*) [real ,in] :: given longitude coordinates of the processed points (degrees).
py (*,*) [real ,in] :: given latitude coordinates of the processed points (degrees).
pzs (*,*) [real ,in] :: resolved model orography
psso_slope (*,*) [real ,out] :: resolved slope tangent
psso_dir (*,*) [real ,out] :: resolved aspect
- Called from:
- subroutine subs/check(status, line)¶
Handles error messages from the netcdf API
- Parameters:
status [integer ,in] :: error number
line [character(len=*),in] :: description line
- Called from:
- Call to:
- subroutine subs/read_nml(knamunit)¶
reads namelist interpolate_safran.nam
- Parameters:
knamunit [integer ,in] :: unit number
- Called from:
- Call to:
Module MODN_INTERPOL_SAFRAN
Description
Module with namelist declarations.
- modn_interpol_safran/NAM_SWITCHES_INT [namelist]¶
- modn_interpol_safran/NAM_MULTIIN_SETTING [namelist]¶
- modn_interpol_safran/NAM_MULTIOUT_SETTING [namelist]¶
- modn_interpol_safran/NAM_SELECT_VARS_SETTING [namelist]¶
- modn_interpol_safran/NAM_OTHER_STUFF [namelist]¶
Needed modules
modd_interpol_safran
: Declarations of namelist variables
Module MODD_INTERPOL_SAFRAN
Description
Declarations of namelist variables
Quick access
- Variables:
hfilein
,hfileout
,hfilesin
,hfilesout
,hgridin
,hgridsin
,hvar_list
,lmultiinput
,lmultioutput
,lselectvar
,lspatialchunk
,ltimechunk
,nlatchunksize
,nlonchunksize
,nnumber_input_files
,nnumber_of_variables
,nnumber_output_grids
Variables
- modd_interpol_safran/hfilein [character(len=100)]¶
Filename of the input file. Used if LMULTIINPUT=.FALSE.
- modd_interpol_safran/hfileout [character(len=100)]¶
Filename of the output file. Used if LMULTIOUTPUT=.FALSE.
- modd_interpol_safran/hfilesin (*) [character(len=100),allocatable]¶
Filenames of the input files. Used if LMULTIINPUT=.TRUE.. The number of filenames given must be equal to NNUMBER_INPUT_FILES.
- modd_interpol_safran/hfilesout (*) [character(len=100),allocatable]¶
Filenames of the output files. Used if LMULTIOUTPUT=.TRUE.. The number of filenames given must be equal to NNUMBER_OUTPUT_GRIDS.
- modd_interpol_safran/hgridin [character(len=100)]¶
Filename of the grid definition file. Used if LMULTIOUTPUT=.FALSE.
- modd_interpol_safran/hgridsin (*) [character(len=100),allocatable]¶
Filenames of the grid definition files. Used if LMULTIOUTPUT=.TRUE.. The number of filenames given must be equal to NNUMBER_OUTPUT_GRIDS.
- modd_interpol_safran/hvar_list (*) [character(len=20),allocatable]¶
List of variable names of variables to interpolate. Used if LSELECTVAR=.TRUE..
- modd_interpol_safran/lmultiinput [logical,optional/default=.false.]¶
Set to .TRUE. in order to handle multiple input files
- modd_interpol_safran/lmultioutput [logical,optional/default=.false.]¶
Set to .TRUE. in order to handle multiple output files
- modd_interpol_safran/lselectvar [logical,optional/default=.false.]¶
Set to .TRUE. to choose variables to interpolate from the input file
- modd_interpol_safran/lspatialchunk [logical,optional/default=.false.]¶
Set to .TRUE. to manualy fix chuncksizes for the first and second spatial dimension when writing the output netcdf file(s).
- modd_interpol_safran/ltimechunk [logical,optional/default=.false.]¶
Set to .TRUE. to write one time step at a time to safe memory when working with larger files
- modd_interpol_safran/nlatchunksize [integer,optional/default=115]¶
Chunksize for the second spatial dimension for writing the output netcdf file(s).
- modd_interpol_safran/nlonchunksize [integer,optional/default=115]¶
Chunksize for the first spatial dimension for writing the output netcdf file(s).
- modd_interpol_safran/nnumber_input_files [integer,optional/default=1]¶
Number of input files. Used if LMULTIINPUT=.TRUE..
- modd_interpol_safran/nnumber_of_variables [integer,optional/default=1]¶
Number of variables to interpolate. Used if LSELECTVAR=.TRUE..
- modd_interpol_safran/nnumber_output_grids [integer,optional/default=1]¶
Number of output grid files. Used if LMULTIOUTPUT=.TRUE..