plugin_src/green/pinel_thickplate.c File Reference

#include <stdio.h>
#include <stdlib.h>
#include <math.h>
#include <gsl/gsl_sf_bessel.h>
#include <gsl/gsl_integration.h>
#include <gsl/gsl_statistics_double.h>
#include "crusde_api.h"

Include dependency graph for pinel_thickplate.c:

Include dependency graph

Defines

#define SIZE_W   100000

Functions

double zeroth_order_bessel_integral (double e, void *params)
double first_order_bessel_integral (double e, void *params)
double f (double x, void *params)
const char * get_name ()
const char * get_version ()
const char * get_authors ()
const PluginCategory get_category ()
const char * get_description ()
void register_parameter ()
 Register parameters this Kernel claims from the input.
int register_output_fields ()
 Register output fields for spatial directions this Green's function calculates.
void request_plugins ()
 Request necessary plugins: crusde_request_green_plugin("elastic halfspace (pinel)").
void init ()
 Initialization of the convolution. Allocation of memory for inputs and outputs.
void run ()
 Performs the fast convolution.
void clear ()
 Clean-up before this plug-in gets unloaded.
int get_value_at (double **result, int x, int y, int t)
 Returns the Green's Function value at Point(x,y) at time 0.

Variables

double g
 gravity [m/s^2]
double E
 Young's Modulus [GPa].
double nu
 Poisson's ratio [-].
double rho_f
 Density of the inviscid fluid [kg/m^3].
double H
 thickness of the elastic plate [m]
double U_vert_const
 Constant part of the vertical displacement equation.
double U_hori_const
 Constant part of the horizontal displacement equation.
double xx = 0.0
double yy = 0.0
double theta = 0.0
double r = 0.0
double local_r = 0.0
double result_j0 = -1.0
double result_j1 = -1.0
double error = 0.0
int quadrant = -1
short sin_sign = -23
short cos_sign = -23
double one_minus_nu_squared = 0.0
double H_squared = 0.0
double H2 = 0.0
double g_rho_f_E_term = 0.0
double g_rho_f_E_term2 = 0.0
double nu4 = 0.0
double sin_h
double cos_h
double A
double B
double D
double global_e
int x_pos
int y_pos
int z_pos
double backpack [3]
 array for return values
double * elastic = NULL
 keeps return values from elastic HS model
int(* elastic_halfspace )(double **, int, int, int)=NULL
gsl_integration_workspace * w
gsl_function j0_integral
gsl_function j1_integral
double lower_integral_bound
double upper_integral_bound
int i = 0

Detailed Description

The Earth is modeled as being made of an elastic layer of arbitrary thickness H lying over an inviscid fluid of density f. The fully relaxed response of the Earth to load changes is modeled. The fluid is considered to be in equilibrium. We calculate surface displacement for the equation of equilibrium (div()=0) in Lagrangian coordinates.

Boundary conditions:<br /> .... The convention is such that a compression corresponds to a positive stress. At the lower surface of of the plate (at z=H) a pressure equal to the hydrostatic pressure in the fluid and no tangential component are applied. Displacements are calculated by integration of Bessel functions. At the surface (z=0) they are expressed by:

\[ G_{z}^{H}(r) = G_{z,elastic}(r) + \frac{(1-\nu^2)g}{E\pi} \int^{\infty}_0{(\frac{B}{D} -1) J_0(\epsilon r)d\epsilon}\\ \]

\[ G_{r}^{H}(r) = G_{r,elastic}(r) + \frac{g}{2\pi E (1-2\nu)} \int^{\infty}_0{\left[(1-\nu^2)(2\epsilon \frac{A}{D} + 4\nu)\right] J_1(\epsilon r)d\epsilon } \]

where $J_0$ and $J_1$ are Bessel functions of the first kind of zeroth- and first-order, respectively. $G_{z}^{H}(r)$ and $G_{r}^{H}(r)$ determine vertical and radial final relaxed displacement at point $r$, respectively. The coefficients $A, B,$ and $D$ are defined by:

\[ A = \epsilon^2 H^2 - \nu (\cosh(2 \epsilon H) - 1) - \frac{2(1-\nu^2)}{E} g \rho_f (\nu \frac{\sinh(2 \epsilon H)}{\epsilon} + H)\\ \]

\[ B = \frac{1}{2} \epsilon \sinh(2 \epsilon H) + \epsilon^2 H^2 + \frac{1-\nu^2}{E} g \rho_f (\cosh(2 \epsilon H) - 1) \]

\[ D = -\epsilon^3 H^2 + \frac{1}{2} \epsilon (\cosh(2 \epsilon H) - 1) + \frac{1-\nu^2}{E} g \rho_f (\sinh(2 \epsilon H) + 2 \epsilon H) \]

where c$ is the variable of integration. If $H \to \infty$, then $\frac{A}{D} \to -\frac{2\nu}{\epsilon}$ and $\frac{A}{D} \to 1$ which gives the solution of the elastic half-space.

The Hankel transform integral is calculated in the bounds $\[0 .. 0.003\]$ using gsl_integration_qag of the GNU Scientific Library (http://www.gnu.org/software/gsl/).


Define Documentation

#define SIZE_W   100000
 


Function Documentation

void clear  ) 
 

Clean-up before this plug-in gets unloaded.

Free the integration workspace!

double f double  x,
void *  params
 

double first_order_bessel_integral double  e,
void *  params
 

const char* get_authors  ) 
 

const PluginCategory get_category  ) 
 

const char* get_description  ) 
 

const char* get_name  ) 
 

int get_value_at double **  result,
int  x,
int  y,
int  t
 

Returns the Green's Function value at Point(x,y) at time 0.

Parameters:
x The x-Coordinate of the requested value.
y The y-Coordinate of the requested value.
t The time-Coordinate of the requested value (set to NO_TIME in non-temporal modelling environment).
Returns:
The value found at Point[t][y][x]

const char* get_version  ) 
 

void init  ) 
 

Initialization of the convolution. Allocation of memory for inputs and outputs.

This function must not be called before register_parameter() unless none of the necessary values depends on parameters provided by the user, which are only set after they have been registered. This function is called some time after register_parameter().

Here U_vert_const and U_hori_const are initialized as follows:

\[ U_{z,const} = \frac{g}{\pi} \frac{(1-\nu ^2)}{E} \]

\[ U_{r, const} = -\frac{g}{2\pi} \frac{(1+\nu)(1-2\nu)}{E} \]

int register_output_fields  ) 
 

Register output fields for spatial directions this Green's function calculates.

This function calls crusde_register_output_field defined in crusde_api.h to register output field in the following order:

  1. x, y, z (horizontal-x, horizontal-y, vertical).

void register_parameter  ) 
 

Register parameters this Kernel claims from the input.

This function calls register_green_param() defined in crusde_api.h to register references to parameters this Green's function will need to operate properly. In case an XML is used to configure the experiment, the reference to this parameter will be identified by the string passed as second argument to register_green_param().

This function registers the references in the following order:

  1. g, XML config identifier is "g" SI-unit: [m/s^2]
  2. E, XML config identifier is "E" SI-unit: [GPa]
  3. nu, XML config identifier is "nu" SI-unit: [-]
  4. rho_f, XML config identifier is "rho_f" SI-unit: [kg/m^3]
  5. H, XML config identifier is "H" SI-unit: [m]

See also:
register_green_param()

register_parameter() (temp_plugin.c.tmp)

void request_plugins  ) 
 

Request necessary plugins: crusde_request_green_plugin("elastic halfspace (pinel)").

empty

void run  ) 
 

Performs the fast convolution.

writes data for actual time step to file

double zeroth_order_bessel_integral double  e,
void *  params
 


Variable Documentation

double A
 

double B
 

double backpack[3]
 

array for return values

double cos_h
 

short cos_sign = -23
 

double D
 

double E
 

Young's Modulus [GPa].

double* elastic = NULL
 

keeps return values from elastic HS model

int(* elastic_halfspace)(double **, int, int, int)=NULL
 

void DOMTreeErrorReporter::error = 0.0
 

double g
 

gravity [m/s^2]

double g_rho_f_E_term = 0.0
 

double g_rho_f_E_term2 = 0.0
 

double global_e
 

double H
 

thickness of the elastic plate [m]

double H2 = 0.0
 

double H_squared = 0.0
 

int i = 0
 

gsl_function j0_integral
 

gsl_function j1_integral
 

double local_r = 0.0
 

double lower_integral_bound
 

double nu
 

Poisson's ratio [-].

double nu4 = 0.0
 

double one_minus_nu_squared = 0.0
 

int quadrant = -1
 

double r = 0.0
 

double result_j0 = -1.0
 

double result_j1 = -1.0
 

double rho_f
 

Density of the inviscid fluid [kg/m^3].

double sin_h
 

short sin_sign = -23
 

double theta = 0.0
 

double U_hori_const
 

Constant part of the horizontal displacement equation.

double U_vert_const
 

Constant part of the vertical displacement equation.

double upper_integral_bound
 

gsl_integration_workspace* w
 

int x_pos
 

double xx = 0.0
 

int y_pos
 

double yy = 0.0
 

int z_pos
 


Generated on Sun Jul 29 08:17:29 2007 for CrusDe by doxygen 1.3.8