Main packageο
Module contentsο
- Module with functions necessary for the creation of faults and the calculation of slip distributions
Author/Github: Alex Villarroel Carrasco / @alexvillarroel University of Concepcion,Chile.
SOME FUNCTIONS ARE COPYRIGHTED BY RODRIGO CIFUENTES LOBOS in his github repository: @RCifuentesLobos : PaleoTsunami
/
- geostochpy.Mw_kauselramirez(ms)[source]ο
calculate the Mw with a ms
- Parameters:
ms ([type]) β [description]
- Returns:
[description]
- Return type:
[type]
- geostochpy.N_max_matriz_covarianza(C)[source]ο
Busca la esquina de la curva l de los valores propios Entrada: C : matriz de covarianza de la falla
- geostochpy.calcula_laplaciano(tamano_sf, F)[source]ο
Entradas tamano_sf: tamano subfalla F: Campo escalar a calcularle el laplaciano
- geostochpy.corr2d_fourier(slip, test, lonslip=[], latslip=[], lontest=[], lattest=[], normalizar=True)[source]ο
Calcula la correlacion cruzada entre 2 arrays 2d utilizando el metodo de la transformada rapida de fourier para obtener indices para filtrar Entradas: slip: array que contiene el slip generado test: arrray de prueba con el cual se busca realzar el filtro (e.g.: anom. de gravedad) lonslip, latslip: longitud y latitud de las subfallas que componen el slip lontest, lattest: longitud y latitud de las subfallas que componen el test
- geostochpy.crea_falla(lats, lons, prof, dip, strike, latini, latfin, area_sf, profundidad, razon_aspecto)[source]ο
Crea una falla a partir de datos del modelo Slab2.0 de Hayes, 2018 Entradas: lats: latitudes de Slab2.0 lons: longitudes de Slab2.0 prof: profundidades de Slab2.0 dip: dip de Slab2.0 strike: strike de Slab2.0 latini: latitud inicial (norte) latfin: latitud final (sur) area_sf: Area deseada para cada subfalla profundidad: profundidad media de la falla, sobre la cual se calculara el punto medio de la distribucion razon_aspecto : Razon de aspecto de la falla (largo/ancho) se calcula el ancho de la falla a partir del largo, para respetar la relacion L/A = const. de Kanamori Anderson, 1975
- geostochpy.dist_sf(lon1, lon2, lat1, lat2)[source]ο
Calcula la distancia en metros entre dos subfallas (sf) i y j utilizando el metodo de karney, 2013 y el elipsoide wgs84 Entradas: lon1: longitud de subfalla i-esima lon2: longitud de subfalla j-esima lat1: latitud de subfalla i-esima lat2: latitud de subfalla j-esima
- geostochpy.dist_sf_alt(lon1, lon2, lat1, lat2)[source]ο
Calcula la distancia en metros entre dos subfallas (sf) i y j resolviendo el problema inverso de la geodesia. Utiliza el elipsoide WGS84. Entradas: lon1: longitud de subfalla i-esima lon2: longitud de subfalla j-esima lat1: latitud de subfalla i-esima lat2: latitud de subfalla j-esima
- geostochpy.dist_sf_vectorized(lon1, lon2, lat1, lat2)[source]ο
Compute the distance between two points on the sphere .
- Parameters:
lon1 ([type]) β [description]
lon2 ([type]) β [description]
lat1 ([type]) β [description]
lat2 ([type]) β [description]
- Returns:
[description]
- Return type:
[type]
- geostochpy.distribucion_slip(C, mu, N)[source]ο
Calcula la distribucion de slip con la expansion de karhunen-loeve log-normal exp(mu)exp(sum(zk*sqrt(lambdak)*vk)) Entradas: C: matriz de covarianza, se le calculan los eigenvalues y vectors mu: matriz con medias N: numero de modos a contar para la sumatoria
- geostochpy.escalar_magnitud_momento(Mw, slip, prof, largo_subfalla, ancho_subfalla, prem=False)[source]ο
Escala la distribucion de slip a la magnitud de momento deseada de entrada. crea una matriz de slip uniforme con magnitud de momento igual a la primitiva para normalizar el slip. Despues escala con matriz uniforme con magnitud de momento deseada Entradas: Mw: magnitud de momento deseada slip: distribucion βprimitivaβ de slip prof: matriz de profundidades de la interfaz lons: longitudes de las subfallas lats: latitudes de las subfallas observacion: los 4 ultimos argumentos son para calcular la magnitud de momento con la funcion βmagnitud_momentoβ definida anteriormente
- geostochpy.estima_rigidez_BilekModel(largo_subfalla, tau, prem=None, profs=None)[source]ο
Estima la rigidez de la falla dado el tiempo de ruptura y el largo de la falla segun la relacion dada en Bilek y Lay, 1999 mu = densidad*largo_falla**2/(0.8**2*tau**2). entradas: largo_falla: largo de la falla en metros tau: tiempo de duracion de la ruptura
- geostochpy.estima_rigidez_prem(profs)[source]ο
Estima la rigidez de cada subfalla en funcion de su profunidad a partir de los datos de velocidad de onda S y densdidad del modelo PREM, dada Vs=sqrt(mu/rho) Entradas: profs: profundidades de las subfallas
- geostochpy.hypocenter(lons, lats, depth, length, width)[source]ο
Generate a random hypocenter from a given grid .
- Parameters:
lons ([type]) β [description]
lats ([type]) β [description]
depth ([type]) β [description]
length ([type]) β [description]
width ([type]) β [description]
- Returns:
[description]
- Return type:
[type]
- geostochpy.interp_slabtofault(lons, lats, nx, ny, slabdep, slabdip, slabstrike, slabrake)[source]ο
- geostochpy.load_PREM()[source]ο
carga los datos de velocidad de onda de corte del modelo prem directorio: directorio donde esta el archivo del modelo slab2
- geostochpy.load_files_slab2(zone='south_america', rake=False)[source]ο
Load data from the hayes 2018 slab2.0 model.
Inputs: zone : str Subduction zone to fech the model. Available zones: -
alaska
: Alaska -calabria
: Calabria -caribbean
: Caribbean -cascadia
: Cascadia -central_america
: Central America -cotabalo
: Cotabalo -halmahera
: Halmahera -hellenic
: Hellenic Arc -himalaya
: Himalaya -hindu_kush
: Hindu Kush -izu_bonin
: Izu-Bonin -kamchatka
: Kamchatka-Kuril Islands-Japan -kermadec
: Kermadec -makran
: Makran -manila_trench
: Manila Trench -muertos_trough
: Muertos Trough -new_guinea
: New Guinea -pamir
: Pamir -philippines
: Philippines -puysegur
: Puysegur -ryukyu
: Ryukyu -scotia_sea
: Scotia Sea -solomon_islands
: Solomon Islands -south_america
: South America -sulawesi
: Sulawesi -sumatra_java
: Sumatra-Java -vanuatu
: Vanuatu If you have and wanna load rake file, invoke function like load_file_slab2(zone,True)
- geostochpy.load_trench(trench_file)[source]ο
upload file with the location of the south american plate trench Input: trench file
- geostochpy.magnitud_momento(slip, prof, largo_subfalla, ancho_subfalla, prem=False)[source]ο
Calcula la magnitud de momento Mw de la distribucion de slip con la ecuacion Mw = 2/3log(Mo)-6.05 donde Mo = mu * area de ruptura * slip Entradas: slip: corresponde a la distribucion de slip prof: matriz con profundidades lons: matriz con longitudes lats: matriz con latitudes
- geostochpy.make_fault_alongtrench(lons_trench, lats_trench, northlat, nx, ny, width, length)[source]ο
Function that creates a fault by calculating the latitudes and longitudes of various subfaults, considering the trench as a boundary. Inputs:
- geostochpy.make_fault_alongtrench_optimized(lons_trench, lats_trench, northlat, nx, ny, width, length)[source]ο
- geostochpy.matriz_covarianza(dip, prof, lons, lats)[source]ο
calcula la matriz de correlacion: Cij = exp(-(dstrike(i,j)/rstrike)-(ddip(i,j)/rdip) donde: dstrike y ddip son estimados de las distancias entre las subfallas i & j ddip = dprof/sin(dip) d es la distancia euclidiana entre dos puntos aproximada segun sus lats y lons rstrike y rdip son los largos de correlacion en cada direccion Entradas: dip: matriz con los dips de cada punto prof: matriz con las profundidades de cada punto lons: matriz con longitud de cada subfalla, necesaria para el calculo de las distancias entre estas lats: matriz con latitud de cada subfalla, necesaria para el calculo de las distancias entre estas Observacion: en esta primera instancia se considera la aproximacion 1 grado approx 111.12 km (generalizar para fallas mas grandes) Observacion 2: aproximacion para calculo de distancia corregida con calculo geodesico de distancia (aumenta el costo computacional bastante)
- geostochpy.matriz_covarianza_optimized(dip, prof, lons, lats, largo_falla, ancho_falla, alpha=0.5)[source]ο
calcula la matriz de correlacion: Cij = exp(-(dstrike(i,j)/rstrike)-(ddip(i,j)/rdip) donde: dstrike y ddip son estimados de las distancias entre las subfallas i & j ddip = dprof/sin(dip) d es la distancia euclidiana entre dos puntos aproximada segun sus lats y lons rstrike y rdip son los largos de correlacion en cada direccion Entradas: dip: matriz con los dips de cada punto prof: matriz con las profundidades de cada punto lons: matriz con longitud de cada subfalla, necesaria para el calculo de las distancias entre estas lats: matriz con latitud de cada subfalla, necesaria para el calculo de las distancias entre estas Observacion: en esta primera instancia se considera la aproximacion 1 grado approx 111.12 km (generalizar para fallas mas grandes) Observacion 2: aproximacion para calculo de distancia corregida con calculo geodesico de distancia (aumenta el costo computacional bastante) Obeservacion 3: el ancho y largo de la falla debe ser en metros
- geostochpy.matriz_medias(media, prof)[source]ο
calcula las medias mu_i para cada subfalla segun mu_i = log(mu*tau_i)-1/2log(alpha**2+1) Entradas: media: valor promedio alrededor del que se desea que se centre el slip prof: matriz con profundidades de cada subfalla, necesaria para el taper
- geostochpy.media_slip(Mw, largo_subfalla, ancho_subfalla, prof)[source]ο
Largo, ancho y profundidad deben estar en metros
- geostochpy.plot_3d(X_array, Y_array, depth, Slip, filename=None)[source]ο
Plot a 3D plot of a 3D Slip grid with a depth of the surface .
- Parameters:
X_array ([type]) β [description]
Y_array ([type]) β [description]
depth ([type]) β [description]
Slip ([type]) β [description]
filename ([type], optional) β [description], defaults to None
- Returns:
[description]
- Return type:
[type]
- geostochpy.plot_deformation(X_grid, Y_grid, lonfosa, latfosa, Deformation, filename, show=False)[source]ο
- geostochpy.plot_deformation_gmt(region, X_grid, Y_grid, dx, dy, lonfosa, latfosa, Deformation, filename)[source]ο
Plot a Slip grid of data on a grid
- Parameters:
region ([type]) β [description]
X_grid ([type]) β [description]
Y_grid ([type]) β [description]
lonfosa ([type]) β [description]
latfosa ([type]) β [description]
Slip ([type]) β [description]
dx ([type]) β [description]
dy ([type]) β [description]
archivo_salida ([type]) β [description]
- geostochpy.plot_slip(X_grid, Y_grid, lonfosa, latfosa, Slip, filename, show=False, cmap='rainbow')[source]ο
- geostochpy.plot_slip_css(region, lons, lats, lonfosa, latfosa, Slip, cmap='rainbow')[source]ο
Plots the CSS projection of a SlipLine
- Parameters:
region ([type]) β [description]
lons ([type]) β [description]
lats ([type]) β [description]
lonfosa ([type]) β [description]
latfosa ([type]) β [description]
Slip ([type]) β [description]
cmap (str, optional) β [description], defaults to βrainbowβ
- geostochpy.plot_slip_gmt(region, X_grid, Y_grid, lonfosa, latfosa, Slip, dx, dy, filename=False)[source]ο
Plot a Slip grid of data on a grid
- Parameters:
region ([type]) β [description]
X_grid ([type]) β [description]
Y_grid ([type]) β [description]
lonfosa ([type]) β [description]
latfosa ([type]) β [description]
Slip ([type]) β [description]
dx ([type]) β [description]
dy ([type]) β [description]
archivo_salida ([type]) β [description]
- geostochpy.plot_slip_gmt_relief(region, X_grid, Y_grid, lonfosa, latfosa, Slip, dx, dy, filename=False)[source]ο
Plot a Slip grid of data on a grid
- Parameters:
region ([type]) β [description]
X_grid ([type]) β [description]
Y_grid ([type]) β [description]
lonfosa ([type]) β [description]
latfosa ([type]) β [description]
Slip ([type]) β [description]
dx ([type]) β [description]
dy ([type]) β [description]
archivo_salida ([type]) β [description]
- geostochpy.taper_LeVeque(d, dmax)[source]ο
realiza el taper propuesto por leveque et al., 2016 t(d)=1-exp(-20|d-dmax|/dmax) donde dmax es la profunidad maxima y de la falla y d es la profundidad de cada subfalla Entrada: d: matriz con profundidades
- geostochpy.taper_except_trench_tukey(depth, dip_taperfunc=<function tukey_window>, strike_taperfunc=<function tukey_window_equal>, alpha_dip=0.5, alpha_strike=0.5)[source]ο
2-D Taper to a Slip. Tukey windows along the strike and along the dip.
- Parameters:
Slip ([type]) β [description]
dip_taperfunc ([type], optional) β [description], defaults to tukey_window
strike_taperfunc ([type], optional) β [description], defaults to windows.tukey
alpha_dip (float, optional) β [description], defaults to 0.5
alpha_strike (float, optional) β [description], defaults to 0.5
- Returns:
[description]
- Return type:
[type]
- geostochpy.taper_slip_fosa(Slip, ventana)[source]ο
Aplica un taper para atenuar el slip en las cercanias de la fosa Entradas: Slip: array con los valores de slip por subfalla ventana: ventana a usarse para aplicar el taper
- geostochpy.tau_ruptura_BilekModel(largo_subfalla, beta=3500, prem=False, profs=None)[source]ο
Para modelo sin PREM, es decir, velocidad de corte constante # Calcula el tiempo de ruptura de la falla dada a partir del largo de falla, considerando la velocidad de ruptura como 0.8 de la velocidad de propagacion de onda S y una velocidad de cizalle de 3.5 km/s utilizando la aproximacion dada en Bilek y Lay, 1999 tau approx L/(0.8beta) Entradas: largo_falla: largo de la falla en metros dado por la funcion que crea la falla beta: velocidad de la onda de cizalle en metros/segundo # Para modelo prem, es decir, velocidad de corte variable Calcula el tiempo de ruptura de la falla dada a partir del largo de falla, considerando la velocidad de ruptura como 0.8 de la velocidad de propagacion de onda S y una velocidad de cizalle dependiente de la profundidad media de la falla y la velocidad vs asociada a esta en el modelo prem, utilizando la aproximacion dada en Bilek y Lay, 1999 tau approx L/(0.8beta) Entradas: largo_falla: largo de la falla en metros dado por la funcion que crea la falla prof_media: profundidad media de la falla prof_prem: profundidades del modelo prem vs: modelo de velocidades del modelo prem asociadas a las profundidades de prof_prem
- geostochpy.tukey_window(N, alpha=0.5)[source]ο
Genera una ventana Tukey.Se ha ajustado los valores de la parte de Hann para que no alcancen completamente cero al comienzo
ParΓ‘metros: - N: Longitud de la ventana. - alpha: ParΓ‘metro de apertura (0 para una ventana rectangular, 1 para una ventana Hann).
Retorna: - Ventana Tukey.
- geostochpy.tukey_window_equal(N, alpha=0.5)[source]ο
Genera una ventana Tukey.Se ha ajustado los valores de la parte de Hann para que no alcancen completamente cero al comienzo
ParΓ‘metros: - N: Longitud de la ventana. - alpha: ParΓ‘metro de apertura (0 para una ventana rectangular, 1 para una ventana Hann).
Retorna: - Ventana Tukey.
Submodulesο
geostochpy.modcomcot moduleο
geostochpy.modfilters moduleο
geostochpy.modokada moduleο
modulo con funciones que calculan la deformacion a partir del slip utilizando las funciones del modulo dtopotools
- geostochpy.modokada.calcula_Mw(lons, lats, razon_aspecto, strike, dip, prof, rake, slip, largo_falla)[source]ο
con la funcion nativa de dtopotools .Mw() se calcula la magnitud de momento de la distribucion de slip Entradas: lons: array con las longitudes de cada subfalla lats: array con las latitudes de cada subfalla razon_aspecto: razon de aspecto de la falla (largo/ancho) strike: array con los strikes de cada subfalla dip: array con los dips de cada subfalla prof: array con la profundidad de cada subfalla en metros rake: array con los rakes de cada subfalla slip: array con el slip de cada subfalla generado con la metodologia de LeVeque, 2016 largo_falla: largo de la falla en metros
- geostochpy.modokada.cargar_deformacion(nombre_archivo, directorio=None, tipo=3)[source]ο
Funcion para cargar los valores del calculo de la deformacion del suelo oceanico obtenidas con dtopo. Por defecto buscara en el cwd, pero se puede especificar un directorio si se desea leer un archivo ubicado en otro directorio. Entradas: nombre_archivo: nombre del archivo tt3 a leer directorio: por defecto se trabajara en el directorio actual, se puede especificar otro si se desea tipo: tipo de archivo dtopo a leer. Por defecto se utilizara en dtopo_type = 3
- geostochpy.modokada.crea_falla_dtopo(lons, lats, razon_aspecto, strike, dip, prof, rake, slip, largo_falla)[source]ο
Crea un objeto Falla generando listas con los atributos de cada subfalla Entradas: lons: array con las longitudes de cada subfalla lats: array con las latitudes de cada subfalla razon_aspecto: razon de aspecto de la falla (largo/ancho) strike: array con los strikes de cada subfalla dip: array con los dips de cada subfalla prof: array con la profundidad de cada subfalla en metros rake: array con los rakes de cada subfalla slip: array con el slip de cada subfalla generado con la metodologia de LeVeque, 2016 largo_falla: largo de la falla en metros
- geostochpy.modokada.guardar_okada(dtopo, nombre_archivo, tipo=3)[source]ο
Funcion para guardar las salidas del calculo de la deformacion del suelo oceanico a partir de la distribucion de slip Entradas: dtopo: deformacion oceanica calculada con el modelo de okada nombre_archivo: nombre del archivo de salida tipo: tipo de data a guardar, por defecto se utiliza dtopo_type = 3
- geostochpy.modokada.leer_okada(nombre_archivo, directorio=None, tipo=3)[source]ο
Funcion para leer las salidas del calculo de la deformacion del suelo oceanico obtenidas con dtopo. Por defecto buscara en el cwd, pero se puede especificar un directorio si se desea leer un archivo ubicado en otro directorio. Entradas: nombre_archivo: nombre del archivo tt3 a leer directorio: por defecto se trabajara en el directorio actual, se puede especificar otro si se desea tipo: tipo de archivo dtopo a leer. Por defecto se utilizara en dtopo_type = 3
- geostochpy.modokada.leer_shape_deformacion(nombre_shape)[source]ο
Funcion para leer shapefiles con datos historicos de deformacion postsismica de terremotos en formato shp de acuerdo a los archivos entregados por cristian s Entradas: nombre de archivo shp. Especificar ruta relativa del archivo shp sin considerar el directorio madre de deformacion
- geostochpy.modokada.leer_xls_observaciones_1960(nombre_xls)[source]ο
Funcion para leer los archivos xls de datos de observaciones del terremoto y tsunami de 1960 otorgados por el dr. patricio winkler el formato de estos archivos es, por columna n sitio , nombre sitio, lat, lon, runup, profundidad de agua, altura de ola, inundacion, cambio de tierra, rango de marea, autor, comentarios Entradas: nombre_xls: nombre del archivo con datos
- geostochpy.modokada.okada_solucion(lons, lats, razon_aspecto, strike, dip, prof, rake, slip, largo_falla, resolucion=0.03333333333333333, tamano_buffer=2.0, verbose=False)[source]ο
Calcula la deformacion provocada por la distribucion de slip generada resolviendo okada utilizando las herramientas de clawpack dtopotools. Se crea un objeto Falla generando listas con los atributos de cada subfalla Entradas: lons: array con las longitudes de cada subfalla lats: array con las latitudes de cada subfalla razon_aspecto: razon de aspecto de la falla (largo/ancho) strike: array con los strikes de cada subfalla dip: array con los dips de cada subfalla prof: array con la profundidad de cada subfalla en metros rake: array con los rakes de cada subfalla slip: array con el slip de cada subfalla generado con la metodologia de LeVeque, 2016 largo_falla: largo de la falla en metros resolucion: resolucion espacial con la que se calculara la deformacion (valor default: 1 minuto) tamano_buffer: tamano (en grados) de espacio βextraβ en los bordes de la falla (valor default: 2 grados) verbose: booleano para mostrar por terminal el estado de la resolucion de okada. True: muestra mensaje False: no muestra mensaje (valor default: False)