DrMZ

DrMZ.advection_diffusion_pde!Method
advection_diffusion_pde!(duhat,uhat,p,t)

RHS for the advection-diffusion equation $u_t = - u_x + ν u_{xx}$ for numerical integration in Fourier space where ν is the viscosity.

source
DrMZ.advection_pde!Method
advection_pde!(duhat,uhat,p,t)

RHS for the advection equation $u_t = - u_x$ for numerical integration in Fourier space.

source
DrMZ.average_errorMethod
average_error(domain,error)

Compute the average error using the trapezoid rule $\frac{1}{T} \int_0^T error(t) dt$.

source
DrMZ.average_ic_errorMethod
average_ic_error(target,prediction)

Compute the two-norm relative error between the prediction and target values for an initial condition.

source
DrMZ.basis_derivativeMethod
basis_derivative(basis,nodes)

Compute the derivatives of each basis function using auto differentiation.

source
DrMZ.basis_evalMethod
basis_eval(basis,nodes)

Evaluate the basis functions for a specified spatial grid or at a specified point. Evaluation for a specified grid returns an $N x M$ matrix with $N$ respresenting the number of nodes and $M$ representing the number of custom basis functions. Evalaution for a specified point returns a $1 x M$ vector.

source
DrMZ.build_basisMethod
build_basis(trunk,L1,L2,M,nodes,weights,L;cutoff=1e-13,tol=1e-12)

Generate a dictionary of orthonormal custom basis functions based on a specified quadrature rule. An orthogonal polynomial expansion is utilized for the evaluation at points away from the specified quadrature nodes.

source
DrMZ.build_basis_expandMethod
build_basis_expand(trunk,L1,L2,M,nodes,weights;cutoff=1e-13,tol=1e-12)

Generate a dictionary of periodic, orthonormal custom basis functions based on the trapezoid rule. B-splines are utilized for the evaluation at points away from the specified quadrature nodes.

To Do: Add time-sampling capabilities...

source
DrMZ.build_dense_modelMethod
function build_dense_model(number_layers,neurons,activations)

Build a feedforward neural network (FFNN) consisting of number_layers of Flux dense layers for the specified number of neurons and activations.

Examples

julia> build_dense_model(2,[(128,128),(128,128)],[relu,relu])
Chain(Dense(128, 128, NNlib.relu), Dense(128, 128, NNlib.relu))
source
DrMZ.burgers_fluxMethod
burgers_flux(u_right,u_left)

Numerical flux for Burgers nonlinearity, $f(u) = \frac{1}{2}u^2$, for use with discontinuous Galerkin method.

source
DrMZ.central_differenceMethod
central_difference(u_j,u_jpos,u_jneg,mu)

Compute the second order central difference for the viscous term of the viscous Burgers equation $u_t = - u u_x + ν u_{xx}$. mu is equal to $ \frac{\nu}{\Delta x^2}$.

source
DrMZ.cheby_diffMethod
cheby_diff(sol,N,L1,L2)

Compute the derivative of using a Chebyshev differentiation matrix on the interval $[L1,L2]$ with N discretization points.

source
DrMZ.cheby_diff_matrixMethod
cheby_diff_matrix(N,L1,L2)

Generate the Chebyshev differentiation matrix for the interval $[L1,L2]$ with N discretization points.

source
DrMZ.cheby_gridMethod
cheby_grid(N,L1,L2)

Generate the grid of Chebyshev points on the interval $[L1,L2]$ with N discretization points.

source
DrMZ.clenshaw_curtisMethod
clenshaw_curtis(N,L1,L2)

Compute the nodes and weights for Clenshaw-Curtis quadrature with N discretization points and on the interval $[L1,L2]$.

source
DrMZ.dlegendre_normMethod
dlegendre_norm(x,L1,L2,n)

Compute the 1st derivative of the n-th shifted Legendre polynomial

source
DrMZ.error_relMethod
error_rel(target,prediction)

Compute the relative error between the prediction and target values.

source
DrMZ.error_seMethod
error_se(target,prediction)

Compute the squared error between the prediction and target values.

source
DrMZ.exp_kernel_periodicMethod
function exp_kernel_periodic(fnc,x_locations;length_scale=0.5)

Covariance kernel for radial basis function (GRF) and IC function for domain.

source
DrMZ.fourier_diffMethod
fourier_diff(sol,N,dL;format="matrix")

Compute the derivative using a Fourier differentiation matrix (default) or the spectral derivative for periodic functions for domain length dL and with N discretization points.

source
DrMZ.gauss_legendreMethod
gauss_legendre(N,L1,L2)

Compute the nodes and weights for Gauss-Legendre quadrature with N discretization points and on the interval $[L1,L2]$.

source
DrMZ.generate_basis_solutionMethod
generate_basis_solution(nodes,weights,tspan,initial_condition,basis,params,pde_function;dt=1e-3,rtol=1e-10,atol=1e-14,peak=1.05)

Generate the solution for a given pde_function and initial_condition on a periodic domain using a custom basis function expansion and a RK45 solver.

source
DrMZ.generate_basis_solution_esdirkMethod
generate_basis_solution_esdirk(nodes,weights,tspan,initial_condition,basis,params,pde_function;dt=1e-4,rtol=1e-8,atol=1e-12)

Generate the solution for a given pde_function and initial_condition on a periodic domain using a custom basis function expansion and a ESDIRK4 solver.

source
DrMZ.generate_basis_solution_implicitMethod
generate_basis_solution_implicit(nodes,weights,tspan,initial_condition,basis,params,pde_function;dt=1e-3,rtol=1e-8,atol=1e-12)

Generate the solution for a given pde_function and initial_condition on a periodic domain using a custom basis function expansion and a Crank-Nicolson solver.

source
DrMZ.generate_fourier_solutionMethod
generate_fourier_solution(L1,L2,tspan,N,initial_condition,pde_function;dt=1e-3,nu=0.1,rtol=1e-10,atol=1e-14)

Generate the solution for a given pde_function and initial_condition on a periodic domain using a N mode Fourier expansion and a RK45 solver.

source
DrMZ.generate_fourier_solution_esdirkMethod
generate_fourier_solution_esdirk(L1,L2,tspan,N,initial_condition,pde_function;dt=1e-4,rtol=1e-8,atol=1e-12,nu=0.1)

Generate the solution for a given pde_function and initial_condition on a periodic domain using a N mode Fourier expansion and a ESDIRK4 solver.

source
DrMZ.generate_fourier_solution_implicitMethod
generate_fourier_solution_implicit(L1,L2,tspan,N,initial_condition,pde_function_explicit,pde_function_implicit;dt=1e-4,rtol=1e-8,atol=1e-12,nu=0.1)

Generate the solution for a given pde_function and initial_condition on a periodic domain using a N mode Fourier expansion and a Crank-Nicolson solver.

source
DrMZ.generate_periodic_functionsMethod
generate_periodic_functions(fnc,x_locations,number_functions,length_scale)

Generate a specified number_functions of random periodic vectors using the exp_kernel_periodic function and a multivariate distribution.

source
DrMZ.generate_periodic_train_testMethod
generate_periodic_train_test(t_span,number_sensors,number_train_functions,number_test_functions,number_solution_points,pde_function_handle;L1=0,L2=2*pi,length_scale=0.5,batch=number_solution_points,dt=1e-3,nu_val=0.1,fnc=(x)->sin(x/2)^2)

Generate the training and testing data for a specified pde_function_handle for periodic boundary conditions using a Fourier spectral method. Defaults to IC $f(\sin^2(x/2))$ and $x ∈ [0,1]$.

source
DrMZ.generate_periodic_train_test_esdirkMethod
generate_periodic_train_test_esdirk(t_span,number_sensors,number_train_functions,number_test_functions,number_solution_points,pde_function_handle;L1=0,L2=2*pi,length_scale=0.5,batch=number_solution_points,dt_size=1e-4,nu_val=0.1,domain="periodic",fnc=(x)->sin(x/2)^2,mode_multiplier=4)

Generate the training and testing data for a specified pde_function_handle for periodic boundary conditions using a Fourier spectral method and a ESDIRK ODE solver. Defaults to IC $f(\sin^2(x/2))$ and $x ∈ [0,1]$.

source
DrMZ.generate_periodic_train_test_implicitMethod
generate_periodic_train_test_implicit(t_span,number_sensors,number_train_functions,number_test_functions,number_solution_points,pde_function_handle;L1=0,L2=2*pi,length_scale=0.5,batch=number_solution_points,dt_size=1e-4,nu_val=0.1,domain="periodic",fnc=(x)->sin(x/2)^2,mode_multiplier=4)

Generate the training and testing data for a specified pde_function_handle for periodic boundary conditions using a Fourier spectral method and a Crank-Nicolson solver. Defaults to IC $f(\sin^2(x/2))$ and $x ∈ [0,1]$.

source
DrMZ.generate_periodic_train_test_initial_conditionsMethod
generate_periodic_train_test_initial_conditions(t_span,number_sensors,number_test_functions,number_train_functions,number_solution_points,pde_function_handle;L1=0,L2=2*pi,length_scale=0.5,batch=number_solution_points,dt=1e-3,nu_val=0.1,domain="periodic",fnc=(x)->sin(x/2)^2)

Generate the training and testing data for a specified pde_function_handle for periodic boundary conditions using a Fourier spectral method. Defaults to IC $f(\sin^2(x/2))$ and $x ∈ [0,1]$.

source
DrMZ.generate_periodic_train_test_musclMethod
generate_periodic_train_test_muscl(t_span,number_sensors,number_test_functions,number_train_functions,number_solution_points,pde_function_handle;L1=0,L2=2*pi,length_scale=0.5,batch=number_solution_points,dt_size=1e-4,upwind_solution_points=4096,fnc=(x)->sin(x/2)^2)

Generate the training and testing data for a specified pde_function_handle for periodic boundary conditions using a MUSCL method. Defaults to IC $f(\sin^2(x/2))$ and $x ∈ [0,1]$.

source
DrMZ.generate_sinusoidal_functions_2_parameterMethod
generate_sinusoidal_functions_2_parameter(x_locations,number_functions)

Generate a specified number_functions of random periodic vectors for the distribution $α \sin(x)+β$ for $α ∈ [-1,1]$ and $β ∈ [-1,1]$.

source
DrMZ.get_1D_energy_customMethod
get_1D_energy_custom(basis,u_solution,L1,L2,weights;multiplier=1/(4*pi))

Compute the energy in the custom basis domain: $ \frac{1}{2} \sum \vert \a_k \vert^2 $. Multiplier defaults to $ \frac{1}{4\pi} $ to match the Fourier calculation.

source
DrMZ.get_1D_energy_custom_coefficientsMethod
get_1D_energy_custom_coefficients(u_coefficients;zeta=1/(4*pi))

Compute the energy in the custom basis domain: $ \zeta \sum \vert \a_k \vert^2 $. $ \zeta $ defaults to $ \frac{1}{4\pi} $ to match the Fourier calculation.

source
DrMZ.get_1D_energy_fftMethod
get_1D_energy_fft(u_solution)

Compute the energy in the Fourier domain using the scaling of $ \frac{1}{N} $. Note: this does not include the 2π multiplier found in Parseval's identity for Fourier series and computes $ \frac{1}{2} \sum \vert \hat{u}_k \vert^2 $.

source
DrMZ.ic_errorMethod
ic_error(target,prediction)

Compute the relative error between the prediction and target values for an initial condition. If the target = 0, the target in denomenator is augmented by $\epsilon_{machine}$.

source
DrMZ.inviscid_burgers_pde!Method
inviscid_burgers_pde!(duhat,uhat,p,t)

RHS for the inviscid Burgers equation $u_t = - u u_x$ for numerical integration in Fourier space.

source
DrMZ.kdv_explicit_pde!Method
kdv_explicit_pde!(duhat,uhat,p,t)

Explicit portion of the RHS for the Korteweg-de Vries equation $u_t = - u u_x - \nu^2 u_{xxx}$ for numerical integration numerical integration in Fourier space.

source
DrMZ.kdv_implicit_pde!Method
kdv_implicit_pde!(duhat,uhat,p,t)

Implicit portion of the RHS for the Korteweg-de Vries equation $u_t = - u u_x - \nu^2 u_{xxx}$ for numerical integration numerical integration in Fourier space.

source
DrMZ.kdv_pde!Method
kdv_pde!(duhat,uhat,p,t)

RHS for the Korteweg-de Vries equation $u_t = - u u_x - \nu^2 u_{xxx}$ for numerical integration numerical integration in Fourier space.

source
DrMZ.ks_explicit_pde!Method
ks_explicit_pde!(duhat,uhat,p,t)

Explicit portion of the RHS for the Kuramoto-Sivashinsky equation $u_t = - u u_x -u_{xx} - \nu u_{xxxx}$ for numerical integration numerical integration in Fourier space.

source
DrMZ.ks_implicit_pde!Method
ks_implicit_pde!(duhat,uhat,p,t)

Implicit portion of the RHS for the Kuramoto-Sivashinsky equation $u_t = - u u_x -u_{xx} - \nu u_{xxxx}$ for numerical integration numerical integration in Fourier space.

source
DrMZ.ks_pde!Method
ks_pde!(duhat,uhat,p,t)

RHS for the Kuramoto-Sivashinsky equation $u_t = - u u_x -u_{xx} - \nu u_{xxxx}$ for numerical integration numerical integration in Fourier space.

source
DrMZ.load_dataMethod
load_data(n_epoch,number_train_functions,number_test_functions,pde_function)

Load the trained branch and trunk neural networks along with the train_data and test_data.

source
DrMZ.load_data_initial_conditionsMethod
load_data_initial_conditions(number_train_functions,number_test_functions,pde_function)

Load the initial conditions from the train_data and test_data.

source
DrMZ.load_data_train_testMethod
load_data_train_test(number_train_functions,number_test_functions,pde_function)

Load the train_data and test_data.

source
DrMZ.load_modelMethod
load_model(n_epoch,pde_function)

Load the trained branch and trunk neural networks.

source
DrMZ.loss_allMethod
loss_all(branch,trunk,initial_conditon,solution_location,target_value)

Compute the mean squared error (MSE) for a complete dataset.

source
DrMZ.mse_errorMethod
mse_error(target,prediction)

Compute the mean squared error between the prediction and target values.

source
DrMZ.norm_infinity_errorMethod
norm_infinity_error(target,prediction)

Compute the two-norm relative error between the prediction and target values.

source
DrMZ.norm_rel_errorMethod
norm_rel_error(target,prediction)

Compute the Euclidean two-norm relative error between the prediction and target values.

source
DrMZ.norm_rel_error_continuousMethod
function norm_rel_error_continuous(target,prediction,weights)

Compute the continuous two-norm relative error between the prediction and target values.

source
DrMZ.predictMethod
predict(branch,trunk,initial_condition,x_locations,t_values)

Predict solution $u(t,x)$ at specified output locations using trained operator neural network branch and trunk.

source
DrMZ.quadratic_nonlinearMethod
quadratic_nonlinear!(uhat,N,dL,alpha)

Compute the convolution sum $\frac{ik}{2}\sum_{p+q=k} u_p u_q$ resulting from the quadratic nonlinearity of Burgers equation $u u_x$ in Fourier space. Convolution sum is padded with the 3/2 rule for dealiasing.

source
DrMZ.reduced_initial_conditionMethod
reduced_initial_condition(L1,L2,x_reduced,x_locations,initial_condition)

Extract the $x$ locations and intial condition values $u(x)$ at a reduced number of equally spaced spatial locations.

source
DrMZ.rhs_advection!Method
rhs_advection!(du,u,p,t)

RHS for the advection equation $u_t = - u_x$ for numerical integration in custom basis space using discontinuous Galerkin method.

source
DrMZ.rhs_advection_diffusion!Method
rhs_advection_diffusion!(du,u,p,t)

RHS for the advection-diffusion equation $u_t = - u_x + ν u_{xx}$ for numerical integration in custom basis space using discontinuous Galerkin method and where ν is the viscosity.

source
DrMZ.rhs_advection_diffusion_dirichlet!Method
rhs_advection_diffusion_dirichlet!(du,u,p,t)

RHS for the advection-diffusion equation $u_t = - u_x + ν u_{xx}$ with Dirichlet boundary conditions for numerical integration in custom basis space using discontinuous Galerkin method and where ν is the viscosity.

source
DrMZ.rhs_explicit_kdv!Method
rhs_explicit_kdv!(du,u,p,t)

Explicit portion of the RHS for the Korteweg-de Vries equation $u_t = - u u_x - \nu^2 u_{xxx}$ for numerical integration in custom basis space using discontinuous Galerkin method.

source
DrMZ.rhs_explicit_ks!Method
rhs_explicit_ks!(du,u,p,t)

Explicit portion of the RHS for the Kuramoto-Sivashinsky equation $u_t = - u u_x -u_{xx} - \nu u_{xxxx}$ for numerical integration in custom basis space using discontinuous Galerkin method.

source
DrMZ.rhs_implicit_kdv!Method
rhs_implicit_kdv!(du,u,p,t)

Implicit portion of the RHS for the Korteweg-de Vries equation $u_t = - u u_x - \nu^2 u_{xxx}$ for numerical integration in custom basis space using discontinuous Galerkin method.

source
DrMZ.rhs_implicit_ks!Method
rhs_implicit_ks!(du,u,p,t)

Implicit portion of the RHS for the Kuramoto-Sivashinsky equation $u_t = - u u_x -u_{xx} - \nu u_{xxxx}$ for numerical integration in custom basis space using discontinuous Galerkin method.

source
DrMZ.rhs_inviscid_burgers!Method
rhs_inviscid_burgers!(du,u,p,t)

RHS for the inviscid Burgers equation $u_t = - u u_x$ for numerical integration in custom basis space using discontinuous Galerkin method.

source
DrMZ.rhs_kdv!Method
rhs_kdv!(du,u,p,t)

RHS for the Korteweg-de Vries equation $u_t = - u u_x - \nu^2 u_{xxx}$ for numerical integration in custom basis space using discontinuous Galerkin method.

source
DrMZ.rhs_ks!Method
rhs_ks!(du,u,p,t)

RHS for the Kuramoto-Sivashinsky equation $u_t = - u u_x -u_{xx} - \nu u_{xxxx}$ for numerical integration in custom basis space using discontinuous Galerkin method.

source
DrMZ.rhs_viscous_burgers!Method
rhs_viscous_burgers!(du,u,p,t)

RHS for the viscous Burgers equation $u_t = - u u_x + ν u_{xx}$ for numerical integration in custom basis space using Discontinuous Galerking method and where ν is the viscosity.

source
DrMZ.save_dataMethod
save_data(train_data,test_data,number_train_functions,number_test_functions,number_solution_points,pde_function)

Save the train_data and test_data.

source
DrMZ.save_modelMethod
save_model(branch,trunk,n_epoch,loss_all_train,loss_all_test,pde_function)

Save the trained branch and trunk neural networks and the training and testing loss history.

source
DrMZ.solution_extractionMethod
solution_extraction(x_locations,t_values,solution,initial_condition,number_solution_points)

Extract the specified number_solution_points randomly from the $u(t,x)$ solution space.

source
DrMZ.solution_interpolationMethod
solution_interpolation(t_span_original,x_locations_original,t_span_interpolate,x_locations_interpolate,solution)

Compute the $u(t,x)$ solution at intermediate $(t,x)$ locations using linear interpolation.

source
DrMZ.solution_spatial_samplingMethod
solution_spatial_sampling(x_prediction,x_target,solution)

Extract the solution values $u(t,x)$ at a reduced number of equally spaced spatial locations.

source
DrMZ.solution_temporal_samplingMethod
solution_temporal_sampling(t_prediction,t_target,solution)

Extract the solution values $u(t,x)$ at a reduced number of equally spaced temporal locations.

source
DrMZ.train_modelMethod
train_model(branch,trunk,n_epoch,train_data;learning_rate=0.00001,save_at=2500,starting_epoch=0)

Train the operator neural network using the mean squared error (MSE) and Adam optimization for n_epochs epochs.

source
DrMZ.trapezoidMethod
trapezoid(N,L1,L2)

Compute the nodes and weights for trapezoid rule with N discretization points and on the interval $[0,L2]$. TO DO: Revise for [L1,L2] support.

source
DrMZ.trapzMethod
trapz(x_range,integrand)

Numerical integration using the multi-application trapezoidal rule.

source
DrMZ.trunk_buildMethod
trunk_build(trunk,M,dtsample)

Generate a dictionary of trunk functions evaluated at $t = 0$ or at a range of times spanning (0:dtsample:1)

source
DrMZ.trunk_ortho_buildMethod
trunk_ortho_build(utilde,L1,L2,nodes,weights,L)

Generate a dictionary of orthonormal custom basis functions which use an orthogonal polynomial expansion for evaluation away from the specified quadrature nodes.

source
DrMZ.trunk_ortho_build_expandMethod
trunk_ortho_build_expand(utilde,nodes;p=17)

Generate a dictionary of periodic, orthonormal custom basis functions which use B-splines for evaluation away from the specified quadrature nodes.

source
DrMZ.viscous_burgers_pde!Method
viscous_burgers_pde!(duhat,uhat,p,t)

RHS for the viscous Burgers equation $u_t = - u u_x + ν u_{xx}$ for numerical integration in Fourier space where ν is the viscosity.

source