DrMZ
DrMZ.advection_diffusion_pde! — Methodadvection_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.
DrMZ.advection_pde! — Methodadvection_pde!(duhat,uhat,p,t)RHS for the advection equation $u_t = - u_x$ for numerical integration in Fourier space.
DrMZ.average_error — Methodaverage_error(domain,error)Compute the average error using the trapezoid rule $\frac{1}{T} \int_0^T error(t) dt$.
DrMZ.average_ic_error — Methodaverage_ic_error(target,prediction)Compute the two-norm relative error between the prediction and target values for an initial condition.
DrMZ.basis_derivative — Methodbasis_derivative(basis,nodes)Compute the derivatives of each basis function using auto differentiation.
DrMZ.basis_eval — Methodbasis_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.
DrMZ.basis_interpolate — Methodbasis_interpolate(coeffs,x,L1,L2)DrMZ.build_basis — Methodbuild_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.
DrMZ.build_basis_expand — Methodbuild_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...
DrMZ.build_dense_model — Methodfunction 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))DrMZ.burgers_flux — Methodburgers_flux(u_right,u_left)Numerical flux for Burgers nonlinearity, $f(u) = \frac{1}{2}u^2$, for use with discontinuous Galerkin method.
DrMZ.central_difference — Methodcentral_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}$.
DrMZ.cheby_diff — Methodcheby_diff(sol,N,L1,L2)Compute the derivative of using a Chebyshev differentiation matrix on the interval $[L1,L2]$ with N discretization points.
DrMZ.cheby_diff_matrix — Methodcheby_diff_matrix(N,L1,L2)Generate the Chebyshev differentiation matrix for the interval $[L1,L2]$ with N discretization points.
DrMZ.cheby_grid — Methodcheby_grid(N,L1,L2)Generate the grid of Chebyshev points on the interval $[L1,L2]$ with N discretization points.
DrMZ.clenshaw_curtis — Methodclenshaw_curtis(N,L1,L2)Compute the nodes and weights for Clenshaw-Curtis quadrature with N discretization points and on the interval $[L1,L2]$.
DrMZ.dbasis_interpolate — Methoddbasis_interpolate(coeffs,x,L1,L2,m)DrMZ.dlegendre_norm — Methoddlegendre_norm(x,L1,L2,n)Compute the 1st derivative of the n-th shifted Legendre polynomial
DrMZ.dlegendre_norm_basis_build — Methoddlegendre_norm_basis_build(nmax,L1,L2)DrMZ.error_rel — Methoderror_rel(target,prediction)Compute the relative error between the prediction and target values.
DrMZ.error_se — Methoderror_se(target,prediction)Compute the squared error between the prediction and target values.
DrMZ.exp_kernel_periodic — Methodfunction exp_kernel_periodic(fnc,x_locations;length_scale=0.5)Covariance kernel for radial basis function (GRF) and IC function for domain.
DrMZ.expansion_approximation — Methodexpansion_approximation(basis,coefficients,nodes)Compute the expansion approximation...
DrMZ.expansion_coefficients — Methodexpansion_coefficients(basis,fnc,nodes,weights)Compute the expansion coefficients...
DrMZ.feature_expansion_set — Methodfeature_expansion_set(L1,L2,neurons,solution_loc)DrMZ.feature_expansion_set_x — Methodfeature_expansion_set_x(L1,L2,neurons,x_loc)DrMZ.feature_expansion_single — Methodfeature_expansion_single(L1,L2,neurons,x)DrMZ.fft_norm — Methodfft_norm(solution)Compute the FFT normalized by $\frac{1}{N}$.
DrMZ.fl — Methodfl(ulv,urv,ubv)DrMZ.fourier_diff — Methodfourier_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.
DrMZ.gauss_legendre — Methodgauss_legendre(N,L1,L2)Compute the nodes and weights for Gauss-Legendre quadrature with N discretization points and on the interval $[L1,L2]$.
DrMZ.generate_basis_solution — Methodgenerate_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.
DrMZ.generate_basis_solution_esdirk — Methodgenerate_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.
DrMZ.generate_basis_solution_implicit — Methodgenerate_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.
DrMZ.generate_fourier_solution — Methodgenerate_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.
DrMZ.generate_fourier_solution_esdirk — Methodgenerate_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.
DrMZ.generate_fourier_solution_implicit — Methodgenerate_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.
DrMZ.generate_muscl_minmod_solution — Methodgenerate_muscl_minmod_solution(L1,L2,t_end,N,u0,pde_function_handle;dt=1e-4,kappa=-1)DrMZ.generate_periodic_functions — Methodgenerate_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.
DrMZ.generate_periodic_train_test — Methodgenerate_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]$.
DrMZ.generate_periodic_train_test_esdirk — Methodgenerate_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]$.
DrMZ.generate_periodic_train_test_implicit — Methodgenerate_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]$.
DrMZ.generate_periodic_train_test_initial_conditions — Methodgenerate_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]$.
DrMZ.generate_periodic_train_test_muscl — Methodgenerate_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]$.
DrMZ.generate_sinusoidal_functions_2_parameter — Methodgenerate_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]$.
DrMZ.get_1D_energy_custom — Methodget_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.
DrMZ.get_1D_energy_custom_coefficients — Methodget_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.
DrMZ.get_1D_energy_fft — Methodget_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 $.
DrMZ.get_1D_energy_upwind — Methodget_1D_energy_upwind(u_solution_full,u_solution,N)DrMZ.ic_error — Methodic_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}$.
DrMZ.ifft_norm — Methodifft_norm(solution)Compute the IFFT normalized by $N$.
DrMZ.inviscid_burgers_pde! — Methodinviscid_burgers_pde!(duhat,uhat,p,t)RHS for the inviscid Burgers equation $u_t = - u u_x$ for numerical integration in Fourier space.
DrMZ.kdv_explicit_pde! — Methodkdv_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.
DrMZ.kdv_implicit_pde! — Methodkdv_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.
DrMZ.kdv_pde! — Methodkdv_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.
DrMZ.ks_explicit_pde! — Methodks_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.
DrMZ.ks_implicit_pde! — Methodks_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.
DrMZ.ks_pde! — Methodks_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.
DrMZ.legendre_norm — Methodlegendre_norm(x,L1,L2,n)Compute the n-th orthonormal, shifted Legendre polynomial
DrMZ.legendre_norm_basis_build — Methodlegendre_norm_basis_build(nmax,L1,L2)DrMZ.linear_reg — Methodlinear_reg(x,y)DrMZ.load_basis — Methodload_basis(basis,M,pde_function)DrMZ.load_branch — Methodload_branch(n_epoch,pde_function)Load the trained branch neural networks.
DrMZ.load_data — Methodload_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.
DrMZ.load_data_initial_conditions — Methodload_data_initial_conditions(number_train_functions,number_test_functions,pde_function)Load the initial conditions from the train_data and test_data.
DrMZ.load_data_train_test — Methodload_data_train_test(number_train_functions,number_test_functions,pde_function)Load the train_data and test_data.
DrMZ.load_model — Methodload_model(n_epoch,pde_function)Load the trained branch and trunk neural networks.
DrMZ.load_trunk — Methodload_trunk(n_epoch,pde_function)Load the trained trunk neural networks.
DrMZ.loss_all — Methodloss_all(branch,trunk,initial_conditon,solution_location,target_value)Compute the mean squared error (MSE) for a complete dataset.
DrMZ.minmod — Methodminmod(x,y)DrMZ.mode_extractor — Methodmode_extractor(uhat,N)DrMZ.mse_error — Methodmse_error(target,prediction)Compute the mean squared error between the prediction and target values.
DrMZ.muscl_minmod_RHS! — Methodmuscl_minmod_RHS!(du,u,p,t)DrMZ.muscl_minmod_viscous_RHS! — Methodmuscl_minmod_viscous_RHS!(du,u,p,t)DrMZ.norm_infinity_error — Methodnorm_infinity_error(target,prediction)Compute the two-norm relative error between the prediction and target values.
DrMZ.norm_rel_error — Methodnorm_rel_error(target,prediction)Compute the Euclidean two-norm relative error between the prediction and target values.
DrMZ.norm_rel_error_continuous — Methodfunction norm_rel_error_continuous(target,prediction,weights)Compute the continuous two-norm relative error between the prediction and target values.
DrMZ.orthonormal_check — Methodorthonormal_check(basis,weights;tol = 1e-12)DrMZ.periodic_fill_domain — Methodperiodic_fill_domain(x_locations)Output the full domain from periodic domain specified for x_locations.
DrMZ.periodic_fill_solution — Methodperiodic_fill_solution(solution)Output the full $u(t,x)$ solution from periodic solution.
DrMZ.predict — Methodpredict(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.
DrMZ.quadratic_nonlinear — Methodquadratic_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.
DrMZ.quadratic_nonlinear_basis — Methodquadratic_nonlinear_basis(u,nonlinear_triple)DrMZ.quadratic_nonlinear_triple_product_basis — Methodquadratic_nonlinear_triple_product_basis(basis_nodes,Dbasis_nodes,nodes,weights)DrMZ.quadratic_nonlinear_triple_product_basis_galerkin — Methodquadratic_nonlinear_triple_product_basis_galerkin(basis_nodes,Dbasis_nodes,nodes,weights)DrMZ.reduced_initial_condition — Methodreduced_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.
DrMZ.rhs_advection! — Methodrhs_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.
DrMZ.rhs_advection_diffusion! — Methodrhs_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.
DrMZ.rhs_advection_diffusion_dirichlet! — Methodrhs_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.
DrMZ.rhs_advection_diffusion_galerkin! — Methodrhs_advection_diffusion_galerkin!(du,u,p,t)DrMZ.rhs_advection_galerkin! — Methodrhs_advection_galerkin!(du,u,p,t)DrMZ.rhs_explicit_kdv! — Methodrhs_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.
DrMZ.rhs_explicit_ks! — Methodrhs_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.
DrMZ.rhs_implicit_kdv! — Methodrhs_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.
DrMZ.rhs_implicit_ks! — Methodrhs_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.
DrMZ.rhs_inviscid_burgers! — Methodrhs_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.
DrMZ.rhs_inviscid_burgers_galerkin! — Methodrhs_inviscid_burgers_galerkin!(du,u,p,t)DrMZ.rhs_inviscid_burgers_pseudo! — Methodrhs_inviscid_burgers_pseudo!(du,u,p,t)DrMZ.rhs_kdv! — Methodrhs_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.
DrMZ.rhs_kdv_galerkin! — Methodrhs_kdv_galerkin!(du,u,p,t)DrMZ.rhs_kdv_pseudo! — Methodrhs_kdv_pseudo!(du,u,p,t)DrMZ.rhs_ks! — Methodrhs_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.
DrMZ.rhs_ks_galerkin! — Methodrhs_ks_galerkin!(du,u,p,t)DrMZ.rhs_ks_pseudo! — Methodrhs_ks_pseudo!(du,u,p,t)No padding to account for the impact of potential aliasing
DrMZ.rhs_viscous_burgers! — Methodrhs_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.
DrMZ.rhs_viscous_burgers_galerkin! — Methodrhs_viscous_burgers_galerkin!(du,u,p,t)
To Do: Simplify input parameters - typical all Galerkin RHS'sDrMZ.rhs_viscous_burgers_pseudo! — Methodrhs_viscous_burgers_pseudo!(du,u,p,t)DrMZ.save_basis — Methodsave_basis(basis,M,pde_function)DrMZ.save_data — Methodsave_data(train_data,test_data,number_train_functions,number_test_functions,number_solution_points,pde_function)Save the train_data and test_data.
DrMZ.save_data_initial_conditions — Methodsave_data_initial_conditions(number_train_functions,number_test_functions)Saves the initial conditions from the train_ic and test_ic.
DrMZ.save_model — Methodsave_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.
DrMZ.solution_extraction — Methodsolution_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.
DrMZ.solution_interpolation — Methodsolution_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.
DrMZ.solution_spatial_sampling — Methodsolution_spatial_sampling(x_prediction,x_target,solution)Extract the solution values $u(t,x)$ at a reduced number of equally spaced spatial locations.
DrMZ.solution_temporal_sampling — Methodsolution_temporal_sampling(t_prediction,t_target,solution)Extract the solution values $u(t,x)$ at a reduced number of equally spaced temporal locations.
DrMZ.spectral_approximation_fourier — Methodspectral_approximation_fourier(x_locations,k,coefficients)DrMZ.train_model — Methodtrain_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.
DrMZ.trapezoid — Methodtrapezoid(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.
DrMZ.trapz — Methodtrapz(x_range,integrand)Numerical integration using the multi-application trapezoidal rule.
DrMZ.trunk_build — Methodtrunk_build(trunk,M,dtsample)Generate a dictionary of trunk functions evaluated at $t = 0$ or at a range of times spanning (0:dtsample:1)
DrMZ.trunk_ortho_build — Methodtrunk_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.
DrMZ.trunk_ortho_build_expand — Methodtrunk_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.
DrMZ.ub — Methodub(ulv,urv)DrMZ.ulnl — Methodulnl(uj,ujn,ujnn,kappa,omega)DrMZ.ulpl — Methodulpl(ujp,uj,ujn,kappa,omega)DrMZ.urnl — Methodurnl(ujp,uj,ujn,kappa,omega)DrMZ.urpl — Methodurpl(ujpp,ujp,uj,kappa,omega)DrMZ.viscous_burgers_pde! — Methodviscous_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.
DrMZ.advection_diffusion_pde!DrMZ.advection_pde!DrMZ.average_errorDrMZ.average_ic_errorDrMZ.basis_derivativeDrMZ.basis_evalDrMZ.basis_interpolateDrMZ.build_basisDrMZ.build_basis_expandDrMZ.build_dense_modelDrMZ.burgers_fluxDrMZ.central_differenceDrMZ.cheby_diffDrMZ.cheby_diff_matrixDrMZ.cheby_gridDrMZ.clenshaw_curtisDrMZ.dbasis_interpolateDrMZ.dlegendre_normDrMZ.dlegendre_norm_basis_buildDrMZ.error_relDrMZ.error_seDrMZ.exp_kernel_periodicDrMZ.expansion_approximationDrMZ.expansion_coefficientsDrMZ.feature_expansion_setDrMZ.feature_expansion_set_xDrMZ.feature_expansion_singleDrMZ.fft_normDrMZ.flDrMZ.fourier_diffDrMZ.gauss_legendreDrMZ.generate_basis_solutionDrMZ.generate_basis_solution_esdirkDrMZ.generate_basis_solution_implicitDrMZ.generate_fourier_solutionDrMZ.generate_fourier_solution_esdirkDrMZ.generate_fourier_solution_implicitDrMZ.generate_muscl_minmod_solutionDrMZ.generate_periodic_functionsDrMZ.generate_periodic_train_testDrMZ.generate_periodic_train_test_esdirkDrMZ.generate_periodic_train_test_implicitDrMZ.generate_periodic_train_test_initial_conditionsDrMZ.generate_periodic_train_test_musclDrMZ.generate_sinusoidal_functions_2_parameterDrMZ.get_1D_energy_customDrMZ.get_1D_energy_custom_coefficientsDrMZ.get_1D_energy_fftDrMZ.get_1D_energy_upwindDrMZ.ic_errorDrMZ.ifft_normDrMZ.inviscid_burgers_pde!DrMZ.kdv_explicit_pde!DrMZ.kdv_implicit_pde!DrMZ.kdv_pde!DrMZ.ks_explicit_pde!DrMZ.ks_implicit_pde!DrMZ.ks_pde!DrMZ.legendre_normDrMZ.legendre_norm_basis_buildDrMZ.linear_regDrMZ.load_basisDrMZ.load_branchDrMZ.load_dataDrMZ.load_data_initial_conditionsDrMZ.load_data_train_testDrMZ.load_modelDrMZ.load_trunkDrMZ.loss_allDrMZ.minmodDrMZ.mode_extractorDrMZ.mse_errorDrMZ.muscl_minmod_RHS!DrMZ.muscl_minmod_viscous_RHS!DrMZ.norm_infinity_errorDrMZ.norm_rel_errorDrMZ.norm_rel_error_continuousDrMZ.orthonormal_checkDrMZ.periodic_fill_domainDrMZ.periodic_fill_solutionDrMZ.predictDrMZ.quadratic_nonlinearDrMZ.quadratic_nonlinear_basisDrMZ.quadratic_nonlinear_triple_product_basisDrMZ.quadratic_nonlinear_triple_product_basis_galerkinDrMZ.reduced_initial_conditionDrMZ.rhs_advection!DrMZ.rhs_advection_diffusion!DrMZ.rhs_advection_diffusion_dirichlet!DrMZ.rhs_advection_diffusion_galerkin!DrMZ.rhs_advection_galerkin!DrMZ.rhs_explicit_kdv!DrMZ.rhs_explicit_ks!DrMZ.rhs_implicit_kdv!DrMZ.rhs_implicit_ks!DrMZ.rhs_inviscid_burgers!DrMZ.rhs_inviscid_burgers_galerkin!DrMZ.rhs_inviscid_burgers_pseudo!DrMZ.rhs_kdv!DrMZ.rhs_kdv_galerkin!DrMZ.rhs_kdv_pseudo!DrMZ.rhs_ks!DrMZ.rhs_ks_galerkin!DrMZ.rhs_ks_pseudo!DrMZ.rhs_viscous_burgers!DrMZ.rhs_viscous_burgers_galerkin!DrMZ.rhs_viscous_burgers_pseudo!DrMZ.save_basisDrMZ.save_dataDrMZ.save_data_initial_conditionsDrMZ.save_modelDrMZ.solution_extractionDrMZ.solution_interpolationDrMZ.solution_spatial_samplingDrMZ.solution_temporal_samplingDrMZ.spectral_approximation_fourierDrMZ.train_modelDrMZ.trapezoidDrMZ.trapzDrMZ.trunk_buildDrMZ.trunk_ortho_buildDrMZ.trunk_ortho_build_expandDrMZ.ubDrMZ.ulnlDrMZ.ulplDrMZ.urnlDrMZ.urplDrMZ.viscous_burgers_pde!