hexed 0.3.0
 
Loading...
Searching...
No Matches
hexed::Solver Class Reference

The main class that basically runs everything. More...

#include <Solver.hpp>

Public Member Functions

 Solver (int n_dim, int row_size, double root_mesh_size, bool local_time_stepping=false, Transport_model viscosity_model=inviscid, Transport_model thermal_conductivity_model=inviscid, std::shared_ptr< Namespace > space=std::make_shared< Namespace >(), std::shared_ptr< Printer_set > printer=std::make_shared< Printer_set >(), bool implicit=false)
 
setup
Namespacenspace ()
 
Meshmesh ()
 Reference to the namespace which can be used to edit user-defined parameters.
 
void read_mesh (std::string file_name, std::vector< Flow_bc * > extremal_bcs, Surface_geom *=nullptr, Flow_bc *surface_bc=nullptr)
 reads mesh from file
 
void read_state (std::string file_name)
 Reads flow state from file.
 
Storage_params storage_params ()
 
void snap_faces ()
 warps the boundary elements such that the element faces coincide with the boundary at their quadrature points.
 
void calc_jacobian (bool snap_faces=true)
 compute the Jacobian of all elements based on the current position of the vertices and value of any face warping.
 
void initialize (const Spacetime_func &)
 set the flow state
 
bool using_art_visc ()
 returns true if artificial viscosity is currently turned on
 
void set_art_visc_off ()
 turns off artificial viscosity
 
void set_art_visc_constant (double)
 turns on artificial viscosity and initializes coefficient to a uniform value
 
void set_art_visc_row_size (int)
 modify the polynomial order of smoothness-based artificial viscosity
 
void set_fix_admissibility (bool)
 turns on/off the thermodynamic admissibility-preserving scheme
 
void set_uncertainty (const Element_func &func)
 set Element::uncertainty for each element according to func.
 
void set_uncert_surface_rep (int bc_sn)
 Set uncertainty metric based on surface representation quality.
 
void synch_extruded_uncert ()
 set uncertainty of each element to be at least the maximum uncertainty of any elements extruded from it
 
time marching
void update ()
 
void update_implicit ()
 (experimental) performs an implicit time step
 
void compute_residual ()
 
void compute_lts_constraints ()
 Computes the minimum ratio between the local diffusive and convective time steps.
 
bool is_admissible ()
 check whether flowfield is admissible (e.g. density and energy are positive)
 
void update_art_visc_smoothness (double advect_length)
 updates the aritificial viscosity coefficient based on smoothness of the flow variables
 
void update_art_visc_elwise (double width, bool pde_based=false)
 (experimental) sets artificial viscosity based on elementwise smoothness
 
void set_art_visc_admis ()
 set the Laplacian artificial viscosity to a minimal value that will encourage thermodynamic admissibility
 
Iteration_status iteration_status ()
 an object providing all available information about the status of the time marching iteration.
 
void reset_counters ()
 
output

functions that compute some form of output data

std::vector< double > sample (int ref_level, bool is_deformed, int serial_n, int i_qpoint, const Qpoint_func &)
 evaluate arbitrary functions at arbitrary locations
 
std::vector< double > sample (int ref_level, bool is_deformed, int serial_n, const Element_func &)
 
const Stopwatch_treestopwatch_tree ()
 obtain performance data
 
std::vector< double > integral_field (const Qpoint_func &integrand)
 compute an integral over the entire flow field at the current time
 
std::vector< double > integral_surface (const Boundary_func &integrand, int bc_sn)
 compute an integral over all surfaces where a particular boundary condition has been enforced
 
std::vector< std::array< double, 2 > > bounds_field (const Qpoint_func &, int n_sample=20)
 compute the min and max of variables over entire flow field.
 
void visualize_field (std::string format, std::string name, const Qpoint_func &output_variables, int n_sample=10, bool wireframe=false)
 write a visualization file describing the entire flow field (but not identifying surfaces)
 
void visualize_surface (std::string format, std::string name, int bc_sn, const Boundary_func &, int n_sample=10, bool wireframe=false)
 write a visualization file describing all surfaces where a particular boundary condition has been enforced.
 
void visualize_contour (std::string format, std::string name, const Qpoint_func &contour_by, const Qpoint_func &output_variables, int n_sample=10)
 
void vis_cart_surf (std::string format, std::string name, int bc_sn, const Boundary_func &func=Uncertainty())
 visualize the Cartesian surface which theoretically exists after element deletion but before any vertex snapping
 
void vis_lts_constraints (std::string format, std::string name, int n_sample=10)
 visualize the local time step constraints imposed by convection and diffusion, respectively
 
void write_state (std::string file_name)
 Writes flow state to file.
 
Array< double > skews ()
 get a list of the Equiangle_skewness for each element
 

Detailed Description

The main class that basically runs everything.

If you want to run a simulation with hexed through the C++ API, your workflow should be roughly the following:

  1. construct a Solver object
  2. interact with the mesh() object to build the mesh topology and/or snap vertices
  3. call calc_jacobian() to initialize internal parameters based on the mesh
  4. call initialize() to initialize the flow state
  5. call update() repeatedly to progress the simulation
  6. call the some of the functions in output section to get the data you want from the simulation

Constructor & Destructor Documentation

◆ Solver()

hexed::Solver::Solver ( int n_dim,
int row_size,
double root_mesh_size,
bool local_time_stepping = false,
Transport_model viscosity_model = inviscid,
Transport_model thermal_conductivity_model = inviscid,
std::shared_ptr< Namespace > space = std::make_shared<Namespace>(),
std::shared_ptr< Printer_set > printer = std::make_shared<Printer_set>(),
bool implicit = false )
Parameters
n_dimnumber of dimensions
row_sizerow size of the basis (see Terminology)
root_mesh_sizesets the value of Mesh::root_mesh_size()
local_time_steppingwhether to use local or global time stepping
viscosity_modeldetermines whether the flow has viscosity (natural, not artificial) and if so, how it depends on temperature
thermal_conductivity_modeldetermines whether the flow has thermal conductivity and if so, how it depends on temperature
spaceNamespace containing any user-defined parameters affecting the behavior of the solver. If no namespace is provided, a new blank namespace is creqated. Any optional parameters which are not found in the namespace shall be created with their default values.
printerwhat to do with any information the solver wants to print for the user to see
implicitif true, allocate storage for solving with an implicit method (experimental feature—not ready for production use)

If viscosity_model and thermal_conductivity_model are both inviscid and you don't turn on artificial viscosity, you will be solving the pure inviscid flow equations. Otherwise, you will be solving the viscous flow equations using the LDG scheme, potentially with some of the diffusion coefficients (artificial viscosity, natural viscosity, thermal conductivity) set to zero.

Member Function Documentation

◆ bounds_field()

std::vector< std::array< double, 2 > > hexed::Solver::bounds_field ( const Qpoint_func & func,
int n_sample = 20 )

compute the min and max of variables over entire flow field.

Layout: {{var0_min, var0_max}, {var1_min, var1_max}, ...}. Bounds are approximated by uniformly sampling a block n_sample-on-a-side in each element.

◆ calc_jacobian()

void hexed::Solver::calc_jacobian ( bool snap_faces = true)

compute the Jacobian of all elements based on the current position of the vertices and value of any face warping.

Mesh topology must be valid (no duplicate or missing connections) before calling this function.

Parameters
snap_facesif true, this function will go ahead and perform face snapping for you

◆ compute_lts_constraints()

void hexed::Solver::compute_lts_constraints ( )

Computes the minimum ratio between the local diffusive and convective time steps.

Assumes no Chebyshev acceleration. Result is written to min_lts_dc_ratio in the HIL namespace.

◆ iteration_status()

Iteration_status hexed::Solver::iteration_status ( )

an object providing all available information about the status of the time marching iteration.

The Iteration_status::start_time member will refer to when the Solver object was created (specifically at the start of the Solver::Solver body).

◆ mesh()

Mesh & hexed::Solver::mesh ( )

Reference to the namespace which can be used to edit user-defined parameters.

fetch the Mesh.

An object the user can use to build the mesh. Note that whenever elements are added, the flow state, and Jacobian are uninitialized, the time step scale is uniformly 1, and the mesh quality may be poor. The functions below must be used to complete the setup before any flow calculation can begin.

◆ read_mesh()

void hexed::Solver::read_mesh ( std::string file_name,
std::vector< Flow_bc * > extremal_bcs,
Surface_geom * geom = nullptr,
Flow_bc * surface_bc = nullptr )

reads mesh from file

Wipes old mesh and flow state. File must be in the native (HDF5-based) mesh format, which you can generate from a previous simulation with Mesh::write. .mesh.h5 will be automatically appended to file_name. New mesh must match the Storage_params of the current one, but the root mesh size will be replaced with that of the new mesh. You still have to initialize (even if you already did before reading the new mesh), but you don't have to calc_jacobian unless you further modify the mesh.

◆ read_state()

void hexed::Solver::read_state ( std::string file_name)

Reads flow state from file.

Essentially a substitute for initialize. .state.h5 will be appended to file name

◆ reset_counters()

void hexed::Solver::reset_counters ( )

reset any variables in iteration_status() that count something since the last call to reset_counters() (e.g. number of artificial viscosity iterations)

◆ sample()

std::vector< double > hexed::Solver::sample ( int ref_level,
bool is_deformed,
int serial_n,
const Element_func & func )

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ set_art_visc_row_size()

void hexed::Solver::set_art_visc_row_size ( int row_size)

modify the polynomial order of smoothness-based artificial viscosity

must be <= row size of discretization (which is the default)

◆ set_fix_admissibility()

void hexed::Solver::set_fix_admissibility ( bool value)

turns on/off the thermodynamic admissibility-preserving scheme

increases robustness at some computational overhead

◆ set_uncert_surface_rep()

void hexed::Solver::set_uncert_surface_rep ( int bc_sn)

Set uncertainty metric based on surface representation quality.

For all deformed elements contacting the boundary specified by bc_sn, computes the unit surface normals at the faces and compares to neighboring elements. Element::uncertainty is set to the total difference between unit normals with all neighbors, where in 3D the total on each edge is computed by Gaussian quadrature in reference space.

Note
Connections between deformed elements and Cartesian elements are ignored. Fully Cartesian connections trivially contribute zero to this metric of uncertainty.

◆ set_uncertainty()

void hexed::Solver::set_uncertainty ( const Element_func & func)

set Element::uncertainty for each element according to func.

Uncertainty metric can be evaluated via sample(ref_level, is_deformed, serial_n, Uncertainty()). This function does some additional work to enforce some conditions on the uncertainty of neighboring elements. Thus, use this function rather than just sample(ref_level, is_deformed, serial_n, func) directly.

◆ update()

void hexed::Solver::update ( )

March the simulation forward by a time step equal to time_step or max_safety times the estimated maximum stable time step, whichever is smaller. Also, the safety factor is not the same as the CFL number (it is scaled by the max allowable CFL for the chosen DG scheme which is often O(1e-2)).

◆ update_art_visc_elwise()

void hexed::Solver::update_art_visc_elwise ( double width,
bool pde_based = false )

(experimental) sets artificial viscosity based on elementwise smoothness

Based on the work of Persson et al. on artificial viscosity with elementwise smoothness indicators. Implemented primarily for evaluating the difference between grid-independent artificial viscosity and conventional methods.

Warning
Not recommended for use in practice. Use Solver::update_art_visc_smoothness, which was developed to supplant this type of approach.

◆ update_implicit()

void hexed::Solver::update_implicit ( )

(experimental) performs an implicit time step

Warning
Experimental! Interesting for reasearch, not effective in practice (yet, anyway).

◆ vis_lts_constraints()

void hexed::Solver::vis_lts_constraints ( std::string format,
std::string name,
int n_sample = 10 )

visualize the local time step constraints imposed by convection and diffusion, respectively

Warning
This function overwrites the reference state, which will invalidate any residual evaluation until update is called again.

◆ visualize_field()

void hexed::Solver::visualize_field ( std::string format,
std::string name,
const Qpoint_func & output_variables,
int n_sample = 10,
bool wireframe = false )

write a visualization file describing the entire flow field (but not identifying surfaces)

Parameters
formatWhich format to write the visualization file in. Accepted values are "xdmf" and "tecplot"
namename of file to write (not including extension)
output_variableswhat variables to write
n_sampleeach element will contain an n_sample by n_sample array of uniformly-spaced sample points
wireframeif true, visualize the mesh edges as a wireframe instead of the filled surface/solid

◆ write_state()

void hexed::Solver::write_state ( std::string file_name)

Writes flow state to file.

.state.h5 will be appended to file_name.


The documentation for this class was generated from the following files: