Tutorial de 5 minutos#
La manera mas accesible de usar pydap es como cliente de acceso a datos cientificos en servidores remotos de OPeNDAP.
OPeNDAP - la visión original#
La vision original de OPeNDAP (Cornillion, et al 1993) fue el hacer las equivalencias
\( \;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\; \boxed{\text{URL} \approx \text{Dataset Remoto} }\)
\( \;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\;\; \boxed{\text{URL + Expresión de Restricción} \approx \text{Subregion de un Dataset Remoto}} \)
Esa vision original fue la que conllevo el desarrollo del protocolo DAP2. En la actualidad, tanto OPeNDAP como Unidata implementan el protocolo DAP4, el cual es mas moderno y abarca mas typos de information, y cubre todos los elementos esenciales cubiertos por el protocolo DAP2 (para mas informacion vea DAP4 specification).
La aportación de PyDAP:#
The logica interna de PyDAP permite la construccion de expressiones de restriccion (CEs for su siglas en ingles) para cada url, de una manera interactiva, facilitando el accesso a subconjunto de datos remotos atraves de OPeNDAP. Ademas, como pydap es un “backend engine” del paquete de Python Xarray, usuarios pueden escalar su flujo de trabajo con la combinacion Xarray+PyDAP y el uso de parallelismo que Xarray permite. En general, un dominio basico del uso de expressiones de restriction (CEs) es importante para maximizar los protocolos de OPeNDAP.
Objetivos:#
Demonstrar como especificar el protocolo DAP4 al servidor remoto de OPeNDAP.
El uso de
XarrayyPyDAPpara descargar un subcojunto de datos remotos en 2 escenario typicos:a)Un archivo remoto con extensionNcMLque representa un archivo virtua l de aggregacion , yb)dos archivos remotos de format Netcdf.Demonstrat las distintas maneras en que pueden user las Condiciones de Restriccion (
CEs), y como estas se pueden pasar al servidor para que cualquier operacion de extraer subconjuntos sea hecha por el servidor OPeNDAP, de una maneraproximaa los archivos remotos, de una manera eficiente.
Requiremientos#
Archivos expuestos por un servidor OPeNDAP que implemente el protocolo DAP4. Por ejemplo, el servidor: http://test.opendap.org/opendap/.
pydap>=3.5.8
xarray>=2025.0
numpy>=2.0
Note
The gran mayoria de los servidores OPeNDAP de la NASA implementan el protocolo DAP4.
from pydap.client import open_url, consolidate_metadata, create_session
import xarray as xr
import numpy as np
# create a session to inspect downloads. cache_name must have `debug`
session = create_session(use_cache=True, cache_kwargs={"cache_name":'data/debug_case1'})
session.cache.clear()
1) Accesso a subconjuntos de un archivo NcML file#
El archivo que utilizaremos tiene formate NcML representando una aggregacion virtual de muchos archivos dataset, el cual puede ser encontrado en el servidos prueba con nombre: aggExisting.ncml.
Los servidores de OPeNDAP pueden ser configurados para producir estos archivos de aggregation con formato NcML. La ventaja es que el usuario trabaja con un solo URL para toda la informacion posible, mientras que en el escenario cuando los archivos no has sido aggregados, el usuario tiene entonces que trabajar con multiples (a veces miles) de URLs.
ncml_url = "http://test.opendap.org/opendap/data/ncml/agg/aggExisting.ncml"
dap4_ncml_url = ncml_url.replace("http", "dap4")
print("=============================================================\n URL DAP4: \n", dap4_ncml_url, "\n=============================================================")
=============================================================
URL DAP4:
dap4://test.opendap.org/opendap/data/ncml/agg/aggExisting.ncml
=============================================================
Ahora utilizamos Xarray y Pydap para “abrir” el archivo, como si estuviera en tu computadora de trabajo. Para esto ejecte el siguiente bloque de codigo
ds = xr.open_dataset(
dap4_ncml_url,
engine='pydap',
session = session,
chunks={},
)
ds
<xarray.Dataset> Size: 12kB
Dimensions: (time: 59, lat: 3, lon: 4)
Coordinates:
* lat (lat) float32 12B 41.0 40.0 39.0
* lon (lon) float32 16B -109.0 -107.0 -105.0 -103.0
* time (time) datetime64[ns] 472B 2000-01-01 2000-01-02 ... 2000-02-28
Data variables:
P (time, lat, lon) float64 6kB dask.array<chunksize=(59, 3, 4), meta=np.ndarray>
T (time, lat, lon) float64 6kB dask.array<chunksize=(59, 3, 4), meta=np.ndarray>
Attributes:
title: Example DataNote
El archivo sigue estando remoto, pero el programa de Xarray junto con Pydap permiten accesar a la information dentro archivo, como si estuviera accesando un archivo dentro de su computadora.
Como dercargamos un elemento de una variable del archivo remoto?#
Para demostrar como descargar datos, primero inspeccionamos la descripcion de la variable T. Es importante, antes que nada, entender la estructura interna del archivo.
ds['T']
<xarray.DataArray 'T' (time: 59, lat: 3, lon: 4)> Size: 6kB
dask.array<open_dataset-T, shape=(59, 3, 4), dtype=float64, chunksize=(59, 3, 4), chunktype=numpy.ndarray>
Coordinates:
* lat (lat) float32 12B 41.0 40.0 39.0
* lon (lon) float32 16B -109.0 -107.0 -105.0 -103.0
* time (time) datetime64[ns] 472B 2000-01-01 2000-01-02 ... 2000-02-28
Attributes:
long_name: surface temperature
units: degC
Maps: ('/time', '/lat', '/lon')Note
La informacion que describe la variable T implica que toda T esta contenida en un solo chunk. Xarray y OPeNDAP en general transmiten y realizan operaciones dividiento la informacion del archivo en chunks. Y Xarray interpreta a toda variable dentro de un archivo remoto en OPeNDAP, como un solo chunk, incluso cuando el archivo remoto divide la representacion de cada variable en diversos chunks.
# clear the cache to inspect what is being downloaded
session.cache.clear()
ds['T'].isel(time=1, lon=0, lat=0).load()
<xarray.DataArray 'T' ()> Size: 8B
array(100.)
Coordinates:
lat float32 4B 41.0
lon float32 4B -109.0
time datetime64[ns] 8B 2000-01-02
Attributes:
long_name: surface temperature
units: degC
Maps: ('/time', '/lat', '/lon')print("====================================== \n Solicitud enviada al Servidor OPeNDAP \n ", session.cache.urls()[0].split("?")[-1].split("&dap4.checksum")[0].replace("%5B","[").replace("%5D","]").replace("%3A",":").replace("%2F","/"), "\n====================================== ")
======================================
Solicitud enviada al Servidor OPeNDAP
dap4.ce=/T[1:1:1][0:1:0][0:1:0]
======================================
La expression de restriccion (CE) fue contruida por el metodo
.isel de Xarray. Este metodo interno de Xarray fue entonces enviada al servidor OPeNDAP, el cual hizo todo el trabajo por nosotros!
2) Accesando a subconjuntos en 2 archivos remotos que pertenecen al mismo proyecto#
En este escenario, los dos archivos remotos describen informacion contigua del mismo proyecto, bajo la suposicion que estos dos archivos pueden aggregarse a lo largo de una dimension. Por ejemplo, un archivo representa los valores de variables geofisicas en la fecha 10/Sept/2025, y el segundo archivo representa valores de las mismas variables geofisicas en la siguiente fecha disponible.
Utilizaremotes para este ejemplo los siguientes archivos: coads_climatology and coads_climatology2. Estos dos archivos abarcan la misma cobertura espacial, y pueden ser aggregados en tiempo.
Note
Es importante verificar siempre que los conjuntos de datos se puedan agregar. PyDAP y Xarray contienen lógica interna que verifica si dos o más conjuntos de datos se pueden concatenar. Sin embargo, estas comprobaciones de seguridad solo consideran dimensiones y coordenadas.
Un paso importante será el uso de Expresiones de Restricción (CE) para garantizar que solo se concatenen las variables de interés.
Warning
Uno de estos archivos tiene variables adicionales que no están presentes en el otro archivo y que descartaremos mediante el uso de CE.
urls = ["http://test.opendap.org/opendap/data/nc/coads_climatology.nc", "http://test.opendap.org/opendap/data/nc/coads_climatology2.nc"]
dap4_urls = [url.replace("http","dap4") for url in urls]
# Expression de Restriccion (CE)
dap4_CE = "?dap4.ce=" + ";".join(["/SST", "/COADSX", "/COADSY", "/TIME"])
# Final list of OPeNDAP URLs
dap4ce_urls =[url+dap4_CE for url in dap4_urls]
print("====================================================\nURLs de OPeNDAP con protocolo DAP4 \n", dap4ce_urls)
====================================================
URLs de OPeNDAP con protocolo DAP4
['dap4://test.opendap.org/opendap/data/nc/coads_climatology.nc?dap4.ce=/SST;/COADSX;/COADSY;/TIME', 'dap4://test.opendap.org/opendap/data/nc/coads_climatology2.nc?dap4.ce=/SST;/COADSX;/COADSY;/TIME']
Note
Q:¿Por qué usar CEs cuando Xarray tiene un método .drop_variables? Porque Xarray necesita analizar primero todos los metadatos remotos para luego descartar las variables. En algunos archivos, es posible encontrar hasta 1000 variables. Xarray las analizaría todas y luego las descartaría. Con CEs, el servidor envía metadatos restringidos asociados únicamente a las variables deseadas. Asi, Xarray solo procesa las variables de importancia.
Warning
Xarray espera la presencia de dimensiones en los metadatos. Al construir la CE, el usuario debe asegurarse de incluir todas las dimensiones asociadas con las variables de interés. En el ejemplo anterior, COASX, COADSY y TIME son las dimensiones de SST.
Consolidate Metadata acelera el proceso de abrir una serie de archivos.#
consolidate_metadata(dap4ce_urls, session=session, concat_dim="TIME")
datacube has dimensions ['COADSX[0:1:179]', 'COADSY[0:1:89]'] , and concat dim: `['TIME']`
Note
consolidate_metadata(dap4_urls, concat_dim='...', session=session) descarga las dimensiones del archivo remoto y las almacena en formato SQLite, para su reuso. Esto significa que el objecto session permite autentificar y actua como un database manager! ¡Esta práctica puede resultar en una mejora del rendimiento en flujos de trabajo entre 10 y 100 veces más rápidos!
Usamos Xarray como herramienta para abrir, descargar, y almacenar la informacion remota.#
Internamente, Xarray utiliza pydap para comunicar con el servidor de OPeNDAP.
ds = xr.open_mfdataset(
dap4ce_urls,
engine='pydap',
concat_dim='TIME',
session=session,
combine="nested",
parallel=True,
decode_times=False,
)
ds
<xarray.Dataset> Size: 2MB
Dimensions: (TIME: 24, COADSY: 90, COADSX: 180)
Coordinates:
* COADSX (COADSX) float64 1kB 21.0 23.0 25.0 27.0 ... 375.0 377.0 379.0
* COADSY (COADSY) float64 720B -89.0 -87.0 -85.0 -83.0 ... 85.0 87.0 89.0
* TIME (TIME) float64 192B 366.0 1.096e+03 ... 7.671e+03 8.401e+03
Data variables:
SST (TIME, COADSY, COADSX) float32 2MB dask.array<chunksize=(12, 90, 180), meta=np.ndarray>
Attributes:
history: FERRET V4.30 (debug/no GUI) 15-Aug-96
Unlimited_Dimension: TIMEds['SST']
<xarray.DataArray 'SST' (TIME: 24, COADSY: 90, COADSX: 180)> Size: 2MB
dask.array<concatenate, shape=(24, 90, 180), dtype=float32, chunksize=(12, 90, 180), chunktype=numpy.ndarray>
Coordinates:
* COADSX (COADSX) float64 1kB 21.0 23.0 25.0 27.0 ... 375.0 377.0 379.0
* COADSY (COADSY) float64 720B -89.0 -87.0 -85.0 -83.0 ... 85.0 87.0 89.0
* TIME (TIME) float64 192B 366.0 1.096e+03 ... 7.671e+03 8.401e+03
Attributes:
long_name: SEA SURFACE TEMPERATURE
history: From coads_climatology
units: Deg C
Maps: ('/TIME', '/COADSY', '/COADSX')Que pasa si queremos descargar un solo elemento#
session.cache.clear()
%%time
ds['SST'].isel(TIME=0, COADSX=0, COADSY=0).load() # this should download a single point one of the files
CPU times: user 9.94 ms, sys: 5.91 ms, total: 15.8 ms
Wall time: 373 ms
<xarray.DataArray 'SST' ()> Size: 4B
array(-1.e+34, dtype=float32)
Coordinates:
COADSX float64 8B 21.0
COADSY float64 8B -89.0
TIME float64 8B 366.0
Attributes:
long_name: SEA SURFACE TEMPERATURE
history: From coads_climatology
units: Deg C
Maps: ('/TIME', '/COADSY', '/COADSX')print("====================================== \n Solicitud enviada al Servidor OPeNDAP:\n ", session.cache.urls()[0].split("?")[-1].split("&dap4.checksum")[0].replace("%5B","[").replace("%5D","]").replace("%3A",":").replace("%2F","/"), "\n====================================== ")
======================================
Solicitud enviada al Servidor OPeNDAP:
dap4.ce=/SST[0:1:11][0:1:89][0:1:179]
======================================
Toda la variable fue descargada innecessariamente !!#
Lo que queremos, es ver que la solicitud enviada al servidor OPeNDAP contenga la siguiente CE:
dap4.ce=/SST[0:1:0][0:1:0][0:1:0]
xr.open_mfdataset no pasa el argumento de seleccion al servidor, de la misma manera en que xr.open_dataset lo hace. En su lugar, Xarray solicita toda la variable, y ya descargada, Xarray hace la selection localmente, de acuerdo al argumento .isel proporcionado por el usuario.
Como asegurar que la selection is enviada al servidor OPeNDAP?#
La respuesta es proporcional le argumento extra, chunk, cuando abrimos/creamos el dataset con Xarray. Este argumento chunk debe ser igual al tamano de selection que esperamos como resultado final.
Warning
Si el argumento chunk is mas pequeno que el tamano que esperanos de nuestra descarga, Xarray terminara enviando muchas solicitudes de descargas innecesarias al servidor, para luego juntar todos los subconjuntos descargados. Este flujo de trabajo tampoco es idea, pues Xarray termina haciendo trabajo extra. Lo ideal es hacer que el servidor haga todo el trabajo, y Xarray solo proporciona el parallelismo.
A continuacion demostramos lograr que Xarray pase la selection al servidor remoto OPeNDAP, el cual hace la mayoria del trabajo de selection cerca del archivo remoto, y solo envia la informacion requerida.
# consolidate metadata again, since the cached metadata was cleared before
consolidate_metadata(dap4ce_urls, session=session, concat_dim="TIME")
datacube has dimensions ['COADSX[0:1:179]', 'COADSY[0:1:89]'] , and concat dim: `['TIME']`
# For a single element in all dimensions, the expected size of the download is:
expected_sizes = {"TIME":1, "COADSX":1, "COADSY":1}
%%time
ds = xr.open_mfdataset(
dap4ce_urls,
engine='pydap',
concat_dim='TIME',
session=session,
combine="nested",
parallel=True,
decode_times=False,
chunks=expected_sizes, # <---------
)
CPU times: user 205 ms, sys: 59.4 ms, total: 264 ms
Wall time: 380 ms
session.cache.clear()
ds['SST'] # inspect chunks before download
<xarray.DataArray 'SST' (TIME: 24, COADSY: 90, COADSX: 180)> Size: 2MB
dask.array<concatenate, shape=(24, 90, 180), dtype=float32, chunksize=(1, 1, 1), chunktype=numpy.ndarray>
Coordinates:
* COADSX (COADSX) float64 1kB 21.0 23.0 25.0 27.0 ... 375.0 377.0 379.0
* COADSY (COADSY) float64 720B -89.0 -87.0 -85.0 -83.0 ... 85.0 87.0 89.0
* TIME (TIME) float64 192B 366.0 1.096e+03 ... 7.671e+03 8.401e+03
Attributes:
long_name: SEA SURFACE TEMPERATURE
history: From coads_climatology
units: Deg C
Maps: ('/TIME', '/COADSY', '/COADSX')%%time
ds['SST'].isel(TIME=0, COADSX=0, COADSY=0).load() # triggers download of an individual chunk
CPU times: user 4.53 s, sys: 117 ms, total: 4.64 s
Wall time: 5.05 s
<xarray.DataArray 'SST' ()> Size: 4B
array(-1.e+34, dtype=float32)
Coordinates:
COADSX float64 8B 21.0
COADSY float64 8B -89.0
TIME float64 8B 366.0
Attributes:
long_name: SEA SURFACE TEMPERATURE
history: From coads_climatology
units: Deg C
Maps: ('/TIME', '/COADSY', '/COADSX')print("====================================== \n Solicitud enviada al Servidor OPeNDAP:\n ", session.cache.urls()[0].split("?")[-1].split("&dap4.checksum")[0].replace("%5B","[").replace("%5D","]").replace("%3A",":").replace("%2F","/"), "\n====================================== ")
======================================
Solicitud enviada al Servidor OPeNDAP:
dap4.ce=/SST[0:1:0][0:1:0][0:1:0]
======================================
Warning: Be cautious about chunking#
¡Ahora solo descargamos exactamente lo que solicitamos! Sin embargo, en algunos casos, el tiempo de descarga puede ser hasta diez veces más lento que cuando solicitamos más datos. La razón de esta lentitud se puede atribuir a la cantidad de fragmentos (chunks) que generó el gráfico de Dask.
Entonces, estos dos son los escenarios que experimentados al tratar de descargar un subconjunto de datos.
No chunk definido. Xarray descarga toda la variable.Chunk definido. Xarray descarga solo el elemento deseado del archivo remoto. Pero durante este processo,
388800chunks fueron creados!
Idealmente, el Chunk manager solo debería activar la descarga de un único fragmento. Sin embargo, se crearon 388800 para garantizar la transferencia del elemento al servidor. Esto, en ocasiones, puede provocar que el cliente (Xarray) tarde mas de lo normal.
En el escenario anterior, llegamos a extremos. Es mejor encontrar un punto medio para los chunks. Lo demostramos a continuación, pero ahora con subconjuntos a lo largo del tiempo.
consolidate_metadata(dap4ce_urls, session=session, concat_dim="TIME")
datacube has dimensions ['COADSX[0:1:179]', 'COADSY[0:1:89]'] , and concat dim: `['TIME']`
download_sizes = {"COADSY":1} # note that we will subset across all time
%%time
ds = xr.open_mfdataset(
dap4ce_urls,
engine='pydap',
concat_dim='TIME',
session=session,
combine="nested",
parallel=True,
decode_times=False,
chunks=download_sizes,
)
session.cache.clear()
CPU times: user 65.9 ms, sys: 43.9 ms, total: 110 ms
Wall time: 224 ms
ds['SST']
<xarray.DataArray 'SST' (TIME: 24, COADSY: 90, COADSX: 180)> Size: 2MB
dask.array<concatenate, shape=(24, 90, 180), dtype=float32, chunksize=(12, 1, 180), chunktype=numpy.ndarray>
Coordinates:
* COADSX (COADSX) float64 1kB 21.0 23.0 25.0 27.0 ... 375.0 377.0 379.0
* COADSY (COADSY) float64 720B -89.0 -87.0 -85.0 -83.0 ... 85.0 87.0 89.0
* TIME (TIME) float64 192B 366.0 1.096e+03 ... 7.671e+03 8.401e+03
Attributes:
long_name: SEA SURFACE TEMPERATURE
history: From coads_climatology
units: Deg C
Maps: ('/TIME', '/COADSY', '/COADSX')%%time
ds['SST'].isel(COADSX=0, COADSY=0).load()
CPU times: user 18 ms, sys: 12.1 ms, total: 30.1 ms
Wall time: 324 ms
<xarray.DataArray 'SST' (TIME: 24)> Size: 96B
array([-1.e+34, -1.e+34, -1.e+34, -1.e+34, -1.e+34, -1.e+34, -1.e+34,
-1.e+34, -1.e+34, -1.e+34, -1.e+34, -1.e+34, -1.e+34, -1.e+34,
-1.e+34, -1.e+34, -1.e+34, -1.e+34, -1.e+34, -1.e+34, -1.e+34,
-1.e+34, -1.e+34, -1.e+34], dtype=float32)
Coordinates:
COADSX float64 8B 21.0
COADSY float64 8B -89.0
* TIME (TIME) float64 192B 366.0 1.096e+03 ... 7.671e+03 8.401e+03
Attributes:
long_name: SEA SURFACE TEMPERATURE
history: From coads_climatology
units: Deg C
Maps: ('/TIME', '/COADSY', '/COADSX')print("====================================== \n Solicitudes enviadas al Servidor OPeNDAP:\n ", [url.split("?")[-1].split("&dap4.checksum")[0].replace("%5B","[").replace("%5D","]").replace("%3A",":").replace("%2F","/") for url in session.cache.urls()], "\n====================================== ")
======================================
Solicitudes enviadas al Servidor OPeNDAP:
['dap4.ce=/SST[0:1:11][0:1:0][0:1:179]', 'dap4.ce=/SST[0:1:11][0:1:0][0:1:179]']
======================================