[P]arallel [Hi]gh-order [Li]brary for [P]DEs
Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
|
Files for the baseline physics. More...
Classes | |
class | AdaptiveSampling |
POD adaptive sampling. More... | |
class | AdaptiveSamplingBase |
class | Adjoint |
Adjoint class. More... | |
class | AnisotropicMeshAdaptation |
class | ArtificialDissipationBase |
Class to add artificial dissipation with an option to add one of the 3 dissipation types: 1. Laplacian, 2. Physical and 3. Enthalpy Laplacian. More... | |
class | ArtificialDissipationFactory |
Creates artificial dissipation pointer. More... | |
class | BaseParameterization |
Abstract class for design parameterization. Objective function and the constraints take this class's pointer as an input to parameterize design variables w.r.t. volume nodes. More... | |
class | BoundPreservingLimiter |
Base Class for implementation of bound preserving limiters. More... | |
class | BoundPreservingLimiterFactory |
This class creates a new BoundPreservingLimiter object based on input parameters. More... | |
class | BoundPreservingLimiterState |
Base Class for bound preserving limiters templated on state. More... | |
class | ConstraintFromObjective_SimOpt |
Creates a constraint from an objective function and a offset value. More... | |
class | DGBase |
DGBase is independent of the number of state variables. More... | |
class | DGBaseState |
Abstract class templated on the number of state variables. More... | |
class | DGFactory |
This class creates a new DGBase object. More... | |
class | DGStrong |
DGStrong class templated on the number of state variables. More... | |
class | DGWeak |
DGWeak class templated on the number of state variables. More... | |
class | DualWeightedResidualError |
DualWeightedResidualError class. More... | |
class | EnthalpyConservingArtificialDissipation |
Adds enthalpy laplacian artificial dissipation (from G.E. Barter and D.L. Darmofal, 2009). More... | |
class | ExactSolutionFactory |
Exact solution function factory. More... | |
class | ExactSolutionFunction |
Exact solution function used for a particular flow setup/case. More... | |
class | ExactSolutionFunction_1DSine |
Exact Solution Function: 1D Sine Function; used for temporal convergence. More... | |
class | ExactSolutionFunction_IsentropicVortex |
Exact Solution Function: Isentropic vortex. More... | |
class | ExactSolutionFunction_Zero |
Exact Solution Function: Zero Function; used as a placeholder when there is no exact solution. More... | |
class | FEValuesShapeHessian |
Class to evaluate hessians of shape functions w.r.t. physical quadrature points. More... | |
class | FlowConstraints |
Use DGBase as a Simulation Constraint ROL::Constraint_SimOpt. More... | |
class | FreeFormDeformation |
Free form deformation class from Sederberg 1986. More... | |
class | FreeFormDeformationParameterization |
FFD design parameterization. Holds an object of FreeFormDeformation and uses it to update the mesh when the control variables are updated. More... | |
class | FreeStreamInitialConditions |
Function used to evaluate farfield conservative solution. More... | |
class | Functional |
Functional base class. More... | |
class | FunctionalErrorNormLpBoundary |
Lp boundary error norm functional class. More... | |
class | FunctionalErrorNormLpVolume |
Lp volume error norm functional class. More... | |
class | FunctionalFactory |
Factory class to construct default functional types. More... | |
class | FunctionalNormLpBoundary |
Lp boundary norm functional class. More... | |
class | FunctionalNormLpVolume |
Lp volume norm functional class. More... | |
class | FunctionalWeightedIntegralBoundary |
Weighted boundary norm functional class. More... | |
class | FunctionalWeightedIntegralVolume |
Weighted volume norm functional class. More... | |
class | GridPostprocessor |
Postprocessor used to output the grid. More... | |
class | HaltonSampling |
Halton sampling. More... | |
class | HighOrderGrid |
class | HyperreducedAdaptiveSampling |
Hyperreduced adaptive sampling. More... | |
class | HyperreducedSamplingErrorUpdated |
Hyperreduced Adaptive Sampling with the updated error indicator. More... | |
class | IdentityParameterization |
Identity design parameterization. Control variables are all volume nodes. More... | |
class | InitialConditionFactory |
Initial condition function factory. More... | |
class | InitialConditionFunction |
Initial condition function used to initialize a particular flow setup/case. More... | |
class | InitialConditionFunction_1DSine |
Initial Condition Function: 1D Sine Function; used for temporal convergence. More... | |
class | InitialConditionFunction_Advection |
Initial Condition Function: 1D Burgers Inviscid. More... | |
class | InitialConditionFunction_AdvectionEnergy |
Initial Condition Function: Advection Energy. More... | |
class | InitialConditionFunction_BurgersInviscid |
Initial Condition Function: 1D Burgers Inviscid. More... | |
class | InitialConditionFunction_BurgersInviscidEnergy |
Initial Condition Function: 1D Burgers Inviscid Energy. More... | |
class | InitialConditionFunction_BurgersRewienski |
Initial Condition Function: 1D Burgers Rewienski. More... | |
class | InitialConditionFunction_BurgersViscous |
Initial Condition Function: 1D Burgers Viscous. More... | |
class | InitialConditionFunction_ConvDiff |
Initial Condition Function: Convection Diffusion Orders of Accuracy. More... | |
class | InitialConditionFunction_ConvDiffEnergy |
Initial Condition Function: Convection Diffusion Energy. More... | |
class | InitialConditionFunction_EulerBase |
Initial Condition Function: Euler Equations (primitive values) More... | |
class | InitialConditionFunction_IsentropicVortex |
Initial Condition Function: Isentropic vortex. More... | |
class | InitialConditionFunction_KHI |
Kelvin-Helmholtz Instability, parametrized by Atwood number. More... | |
class | InitialConditionFunction_LeblancShockTube |
Initial Condition Function: 1D Leblanc Shock Tube. More... | |
class | InitialConditionFunction_LowDensity2D |
Initial Condition Function: 2D Low Density Euler. More... | |
class | InitialConditionFunction_ShuOsherProblem |
Initial Condition Function: 1D Shu Osher Problem. More... | |
class | InitialConditionFunction_SodShockTube |
Initial Condition Function: 1D Sod Shock Tube. More... | |
class | InitialConditionFunction_TaylorGreenVortex |
Initial Condition Function: Taylor Green Vortex (uniform density) More... | |
class | InitialConditionFunction_TaylorGreenVortex_Isothermal |
Initial Condition Function: Taylor Green Vortex (isothermal density) More... | |
class | InitialConditionFunction_Zero |
Initial condition 0. More... | |
class | InnerVolParameterization |
Design parameterization w.r.t. inner volume nodes (i.e. volume nodes excluding those on the boundary). More... | |
class | LaplacianArtificialDissipation |
Adds Laplacian artificial dissipation (from Persson and Peraire, 2008). More... | |
class | LiftDragFunctional |
class | LocalSolution |
Class to store local solution coefficients and provide evaluation functions. More... | |
class | ManufacturedSolutionAdd |
Sum of sine waves manufactured solution. More... | |
class | ManufacturedSolutionAtan |
Hump manufactured solution based on arctangent functions. More... | |
class | ManufacturedSolutionBoundaryLayer |
Scalar boundary layer manufactured solution. More... | |
class | ManufacturedSolutionCosine |
Product of cosine waves manufactured solution. More... | |
class | ManufacturedSolutionEvenPoly |
Sum of even order polynomial functions manufactured solution. More... | |
class | ManufacturedSolutionExample |
class | ManufacturedSolutionExp |
Sum of exponential functions manufactured solution. More... | |
class | ManufacturedSolutionFactory |
Manufactured solution function factory. More... | |
class | ManufacturedSolutionFunction |
Manufactured solution used for grid studies to check convergence orders. More... | |
class | ManufacturedSolutionNavah_MS1 |
Navah and Nadarajah free flows manufactured solution: MS1. More... | |
class | ManufacturedSolutionNavah_MS2 |
Navah and Nadarajah free flows manufactured solution: MS2. More... | |
class | ManufacturedSolutionNavah_MS3 |
Navah and Nadarajah free flows manufactured solution: MS3. More... | |
class | ManufacturedSolutionNavah_MS4 |
Navah and Nadarajah free flows manufactured solution: MS4. More... | |
class | ManufacturedSolutionNavah_MS5 |
Navah and Nadarajah free flows manufactured solution: MS5. More... | |
class | ManufacturedSolutionNavahBase |
class | ManufacturedSolutionPoly |
Sum of polynomial manufactured solution. More... | |
class | ManufacturedSolutionQuadratic |
Quadratic function manufactured solution. More... | |
class | ManufacturedSolutionSine |
Product of sine waves manufactured solution. More... | |
class | ManufacturedSolutionSShock |
S-Shock manufactured solution. More... | |
class | ManufacturedSolutionZero |
Product of zero waves manufactured solution. More... | |
class | MaximumPrincipleLimiter |
Class for implementation of Maximum-Principle-Satisfying limiter derived from BoundPreservingLimiterState class. More... | |
class | MeshAdaptation |
class | MeshErrorEstimateBase |
Abstract class to estimate error for mesh adaptation. More... | |
class | MeshErrorFactory |
Returns pointer to appropriate mesh error class depending on the input parameters. More... | |
class | MeshTypeHelper |
class | MeshTypeHelper< dealii::parallel::distributed::Triangulation< PHILIP_DIM > > |
class | MeshTypeHelper< dealii::parallel::shared::Triangulation< PHILIP_DIM > > |
class | MeshTypeHelper< dealii::Triangulation< PHILIP_DIM > > |
class | MetricToMeshGenerator |
Class to convert metric field to mesh using BAMG. More... | |
class | NNLSSolver |
class | OutletPressureIntegral |
class | PhysicalArtificialDissipation |
Adds Physical artificial dissipation (from Persson and Peraire, 2008). More... | |
class | PositivityPreservingLimiter |
Class for implementation of two forms of the Positivity-Preserving limiter derived from BoundPreservingLimiterState class. More... | |
class | PreconditionILU_RAS |
class | ResidualErrorEstimate |
Class to compute residual based error. More... | |
class | ROLObjectiveSimOpt |
Interface between the ROL::Objective_SimOpt PHiLiP::Functional. More... | |
class | SetInitialCondition |
Class for setting/applying the initial condition. More... | |
class | SolutionIntegral |
Functional to take the integral of the solution. More... | |
class | TargetBoundaryFunctional |
class | TargetFunctional |
TargetFunctional base class. More... | |
class | TargetWallPressure |
class | TVBLimiter |
Class for implementation of a TVD/TVB limiter derived from BoundPreservingLimiterState class. More... | |
Typedefs | |
using | FadType = Sacado::Fad::DFad< double > |
Sacado AD type for first derivatives. | |
using | FadFadType = Sacado::Fad::DFad< FadType > |
Sacado AD type that allows 2nd derivatives. | |
using | codi_FadType = codi::RealForwardGen< double, codi::Direction< double, dimForwardAD > > |
Tapeless forward mode. | |
using | codi_JacobianComputationType = codi::RealReverseIndexVec< dimReverseAD > |
Reverse mode type for Jacobian computation using TapeHelper. | |
using | codi_HessianComputationType = codi::RealReversePrimalIndexGen< codi::RealForwardVec< dimForwardAD >, codi::Direction< codi::RealForwardVec< dimForwardAD >, dimReverseAD > > |
Nested reverse-forward mode type for Jacobian and Hessian computation using TapeHelper. | |
using | RadType = codi_JacobianComputationType |
CoDiPaco reverse-AD type for first derivatives. | |
using | RadFadType = codi_HessianComputationType |
Nested reverse-forward mode type for Jacobian and Hessian computation using TapeHelper. | |
using | dealii_Vector = dealii::LinearAlgebra::distributed::Vector< double > |
using | AdaptVector = dealii::Rol::VectorAdaptor< dealii_Vector > |
using | ROL_Vector = ROL::Vector< double > |
using | DealiiVector = dealii::LinearAlgebra::distributed::Vector< double > |
Functions | |
template<int dim, typename real > | |
std::vector< real > | project_function (const std::vector< real > &function_coeff, const dealii::FESystem< dim, dim > &fe_input, const dealii::FESystem< dim, dim > &fe_output, const dealii::QGauss< dim > &projection_quadrature) |
Get the coefficients of a function projected onto a set of basis (to be replaced with operators->projection_operator). | |
template<int dim, int nstate, typename real2 > | |
void | compute_br2_correction (const dealii::FESystem< dim, dim > &fe_soln, const LocalSolution< real2, dim, dim > &metric_solution, const std::vector< State< real2, nstate >> &lifting_op_R_rhs, std::vector< State< real2, nstate >> &soln_grad_correction) |
template<int dim, int nstate, typename real2 > | |
void | correct_the_gradient (const std::vector< State< real2, nstate >> &soln_grad_corr, const dealii::FESystem< dim, dim > &fe_soln, const std::vector< DirectionalState< real2, dim, nstate >> &soln_jump, const dealii::FullMatrix< double > &interpolation_operator, const std::array< dealii::FullMatrix< real2 >, dim > &gradient_operator, std::vector< DirectionalState< real2, dim, nstate >> &soln_grad) |
template<int dim> | |
dealii::Quadrature< dim > | project_face_quadrature (const dealii::Quadrature< dim - 1 > &face_quadrature_lower_dim, const std::pair< unsigned int, int > face_subface_pair, const typename dealii::QProjector< dim >::DataSetDescriptor face_data_set) |
template<int dim, typename real1 , typename real2 > | |
dealii::Tensor< 1, dim, real1 > | vmult (const dealii::Tensor< 2, dim, real1 > A, const dealii::Tensor< 1, dim, real2 > x) |
template<int dim, typename real1 > | |
real1 | norm (const dealii::Tensor< 1, dim, real1 > x) |
Returns norm of dealii::Tensor<1,dim,real> More... | |
std::pair< unsigned int, double > | solve_linear3 (const dealii::TrilinosWrappers::SparseMatrix &system_matrix, dealii::LinearAlgebra::distributed::Vector< double > &right_hand_side, dealii::LinearAlgebra::distributed::Vector< double > &solution, const Parameters::LinearSolverParam ¶m) |
std::pair< unsigned int, double > | solve_linear (const dealii::TrilinosWrappers::SparseMatrix &system_matrix, dealii::LinearAlgebra::distributed::Vector< double > &right_hand_side, dealii::LinearAlgebra::distributed::Vector< double > &solution, const Parameters::LinearSolverParam ¶m) |
std::pair< unsigned int, double > | solve_linear_2 (const dealii::TrilinosWrappers::SparseMatrix &system_matrix, const dealii::LinearAlgebra::distributed::Vector< double > &right_hand_side, dealii::LinearAlgebra::distributed::Vector< double > &solution, const Parameters::LinearSolverParam ¶m) |
template<int spacedim> | |
void | assign_1d_boundary_ids (const std::map< unsigned int, dealii::types::boundary_id > &boundary_ids, dealii::Triangulation< 1, spacedim > &triangulation) |
template<int dim> | |
void | rotate_indices (std::vector< unsigned int > &numbers, const unsigned int n_indices_per_direction, const char direction, const bool mesh_reader_verbose_output) |
template<int dim, int spacedim> | |
void | assign_1d_boundary_ids (const std::map< unsigned int, dealii::types::boundary_id > &, dealii::Triangulation< dim, spacedim > &) |
void | open_file_toRead (const std::string filepath, std::ifstream &file_in) |
void | read_gmsh_entities (std::ifstream &infile, std::array< std::map< int, int >, 4 > &tag_maps) |
template<int spacedim> | |
void | read_gmsh_nodes (std::ifstream &infile, std::vector< dealii::Point< spacedim >> &vertices, std::map< int, int > &vertex_indices, const bool mesh_reader_verbose_output) |
unsigned int | gmsh_cell_type_to_order (unsigned int cell_type) |
template<int dim> | |
unsigned int | find_grid_order (std::ifstream &infile, const bool mesh_reader_verbose_output) |
unsigned int | ijk_to_num (const unsigned int i, const unsigned int j, const unsigned int k, const unsigned int n_per_line) |
unsigned int | ij_to_num (const unsigned int i, const unsigned int j, const unsigned int n_per_line) |
std::vector< unsigned int > | face_node_finder (const unsigned int degree, const bool mesh_reader_verbose_output) |
template<int dim> | |
std::vector< unsigned int > | gmsh_hierarchic_to_lexicographic (const unsigned int degree, const bool mesh_reader_verbose_output) |
unsigned int | dealii_node_number (const unsigned int i, const unsigned int j, const unsigned int k) |
void | fe_q_node_number (const unsigned int index, unsigned int &i, unsigned int &j, unsigned int &k) |
template<int dim, int spacedim> | |
bool | get_new_rotated_indices (const dealii::CellAccessor< dim, spacedim > &cell, const std::vector< dealii::Point< spacedim >> &all_vertices, const std::vector< unsigned int > &deal_h2l, const std::vector< unsigned int > &rotate_z90degree, std::vector< unsigned int > &high_order_vertices_id) |
template<int dim, int spacedim> | |
bool | get_new_rotated_indices_3D (const dealii::CellAccessor< dim, spacedim > &cell, const std::vector< dealii::Point< spacedim >> &all_vertices, const std::vector< unsigned int > &deal_h2l, const std::vector< unsigned int > &rotate_x90degree_3D, const std::vector< unsigned int > &rotate_y90degree_3D, const std::vector< unsigned int > &rotate_z90degree_3D, std::vector< unsigned int > &high_order_vertices_id_rotated, const bool) |
template<int dim, int spacedim> | |
std::shared_ptr< HighOrderGrid< dim, double > > | read_gmsh (std::string filename, const bool periodic_x, const bool periodic_y, const bool periodic_z, const int x_periodic_1, const int x_periodic_2, const int y_periodic_1, const int y_periodic_2, const int z_periodic_1, const int z_periodic_2, const bool mesh_reader_verbose_output, const bool do_renumber_dofs, int requested_grid_order, const bool use_mesh_smoothing) |
template<int dim, int spacedim> | |
std::shared_ptr< HighOrderGrid< dim, double > > | read_gmsh (std::string filename, const bool do_renumber_dofs, int requested_grid_order=0, const bool use_mesh_smoothing=true) |
Reads Gmsh grid from file at a given requested_grid_order and use_mesh_smoothing input. | |
template std::shared_ptr< HighOrderGrid< PHILIP_DIM, double > > | read_gmsh< PHILIP_DIM, PHILIP_DIM > (std::string filename, const bool periodic_x, const bool periodic_y, const bool periodic_z, const int x_periodic_1, const int x_periodic_2, const int y_periodic_1, const int y_periodic_2, const int z_periodic_1, const int z_periodic_2, const bool mesh_reader_verbose_output, const bool do_renumber_dofs, int requested_grid_order, const bool use_mesh_smoothing) |
template std::shared_ptr< HighOrderGrid< PHILIP_DIM, double > > | read_gmsh< PHILIP_DIM, PHILIP_DIM > (std::string filename, const bool do_renumber_dofs, int requested_grid_order, const bool use_mesh_smoothing) |
template<typename real > | |
std::vector< real > | matrix_stdvector_mult (const dealii::FullMatrix< double > &matrix, const std::vector< real > &vector_in) |
template<typename T > | |
std::vector< T > | flatten (const std::vector< std::vector< T >> &v) |
const dealii::LinearAlgebra::distributed::Vector< double > & | ROL_vector_to_dealii_vector_reference (const ROL::Vector< double > &x) |
Access the read-write deali.II Vector stored within the ROL::Vector. More... | |
dealii::LinearAlgebra::distributed::Vector< double > & | ROL_vector_to_dealii_vector_reference (ROL::Vector< double > &x) |
Access the read-only deali.II Vector stored within the ROL::Vector. More... | |
std::string | get_padded_mpi_rank_string (const int mpi_rank_input) |
bool | isfinite (double value) |
< Provide isfinite for double. More... | |
bool | isfinite (Sacado::Fad::DFad< double > value) |
Provide isfinite for FadFadType. | |
bool | isfinite (Sacado::Fad::DFad< Sacado::Fad::DFad< double >> value) |
Provide isfinite for RadFadType. | |
bool | isfinite (Sacado::Rad::ADvar< Sacado::Fad::DFad< double >> value) |
Variables | |
static constexpr int | dimForwardAD = 1 |
Size of the forward vector mode for CoDiPack. | |
static constexpr int | dimReverseAD = 1 |
Size of the reverse vector mode for CoDiPack. | |
static constexpr int | dimForwardAD = 1 |
Size of the forward vector mode for CoDiPack. | |
static constexpr int | dimReverseAD = 1 |
Size of the reverse vector mode for CoDiPack. | |
Files for the baseline physics.
void PHiLiP::assign_1d_boundary_ids | ( | const std::map< unsigned int, dealii::types::boundary_id > & | boundary_ids, |
dealii::Triangulation< 1, spacedim > & | triangulation | ||
) |
In 1d, boundary indicators are associated with vertices, but this is not currently passed through the SubcellData structure. This function sets boundary indicators on vertices after the triangulation has already been created.
TODO: Fix this properly via SubcellData
Definition at line 35 of file gmsh_reader.cpp.
std::vector<unsigned int> PHiLiP::face_node_finder | ( | const unsigned int | degree, |
const bool | mesh_reader_verbose_output | ||
) |
Solely used in the 3D case. Helps finding the face_nodes during recursive call.
degree | (Degree of the inner cube) -> Degree is -2 in this function call because we lose two points when we deal with the inner faces/cube. -> 1D example: i.e. (1-2-3-4-5), 4th order line with 5 points, if we switch to inner line, we obtain (2-3-4), hence, we obtain a 2nd order line with 3 points. |
For Debug output
Definition at line 582 of file gmsh_reader.cpp.
unsigned int PHiLiP::find_grid_order | ( | std::ifstream & | infile, |
const bool | mesh_reader_verbose_output | ||
) |
Finds the grid order
Definition at line 502 of file gmsh_reader.cpp.
bool PHiLiP::get_new_rotated_indices | ( | const dealii::CellAccessor< dim, spacedim > & | cell, |
const std::vector< dealii::Point< spacedim >> & | all_vertices, | ||
const std::vector< unsigned int > & | deal_h2l, | ||
const std::vector< unsigned int > & | rotate_z90degree, | ||
std::vector< unsigned int > & | high_order_vertices_id | ||
) |
Function to get rotated indices in 2D
Definition at line 1063 of file gmsh_reader.cpp.
bool PHiLiP::get_new_rotated_indices_3D | ( | const dealii::CellAccessor< dim, spacedim > & | cell, |
const std::vector< dealii::Point< spacedim >> & | all_vertices, | ||
const std::vector< unsigned int > & | deal_h2l, | ||
const std::vector< unsigned int > & | rotate_x90degree_3D, | ||
const std::vector< unsigned int > & | rotate_y90degree_3D, | ||
const std::vector< unsigned int > & | rotate_z90degree_3D, | ||
std::vector< unsigned int > & | high_order_vertices_id_rotated, | ||
const bool | |||
) |
Function to get rotated indices in 3D
TECHNICALLY, THE MATCHING VECTOR SHOULD BE ALL 0 IF THEY ALL MATCH
Definition at line 1115 of file gmsh_reader.cpp.
std::vector<unsigned int> PHiLiP::gmsh_hierarchic_to_lexicographic | ( | const unsigned int | degree, |
const bool | mesh_reader_verbose_output | ||
) |
For Debug output
Now, we have an inside hex for hex of order 3 and more (2 has a single point) Idea now is to use recursion and apply the same logic to the inner block
Once this is out, we need to do some node processing, since the nodes are not at the correct locations Use the global index information to track it, i.e., get the transformed indices, and allocate them to the global index. This would make sense since the global index are always true for both GMSH and DEAL.II
Definition at line 675 of file gmsh_reader.cpp.
bool PHiLiP::isfinite | ( | double | value | ) |
< Provide isfinite for double.
Provide isfinite for FadType
Definition at line 15 of file manufactured_solution.cpp.
real1 PHiLiP::norm | ( | const dealii::Tensor< 1, dim, real1 > | x | ) |
Returns norm of dealii::Tensor<1,dim,real>
Had to rewrite this instead of x.norm() because norm() doesn't allow the use of codi variables.
Definition at line 41 of file target_functional.cpp.
std::shared_ptr< HighOrderGrid< dim, double > > PHiLiP::read_gmsh | ( | std::string | filename, |
const bool | periodic_x, | ||
const bool | periodic_y, | ||
const bool | periodic_z, | ||
const int | x_periodic_1, | ||
const int | x_periodic_2, | ||
const int | y_periodic_1, | ||
const int | y_periodic_2, | ||
const int | z_periodic_1, | ||
const int | z_periodic_2, | ||
const bool | mesh_reader_verbose_output, | ||
const bool | do_renumber_dofs, | ||
int | requested_grid_order = 0 , |
||
const bool | use_mesh_smoothing = true |
||
) |
Reads Gmsh grid from input file, supports periodic boundaries. Can request to convert the input grid's order to the requested_grid_order, which will simply interpolate the high-order nodes. Dealii's mesh smoothing can be set to none while using goal oriented mesh adaptation.
When dimEntity == dim, this means we found a Face (2D) or Cell (3D)
For 3D ROTATIONS (INITIATING ONCE HERE, SO WE DON'T RECOMPUTE FOR EACH ROTATION/CELL)
Go through all cells and perform rotations to match gmsh with deal.ii
Convert the equidistant points from Gmsh into the GLL points used by FE_Q in deal.II.
Convert the mesh by interpolating from one order to another.
Definition at line 1247 of file gmsh_reader.cpp.
const dealii::LinearAlgebra::distributed::Vector< double > & PHiLiP::ROL_vector_to_dealii_vector_reference | ( | const ROL::Vector< double > & | x | ) |
Access the read-write deali.II Vector stored within the ROL::Vector.
Note that the ROL::Vector<double> x
should actually be a dealii::Rol::VectorAdaptor, which is derived from ROL::Vector<double>. x
is dynamically downcasted into the VectorAdaptor, which in turn returns a reference to the stored dealii::LinearAlgebra::distributed::Vector<double>.
Definition at line 6 of file rol_to_dealii_vector.cpp.
dealii::LinearAlgebra::distributed::Vector< double > & PHiLiP::ROL_vector_to_dealii_vector_reference | ( | ROL::Vector< double > & | x | ) |
Access the read-only deali.II Vector stored within the ROL::Vector.
Note that the ROL::Vector<double> x
should actually be a dealii::Rol::VectorAdaptor, which is derived from ROL::Vector<double>. x
is dynamically downcasted into the VectorAdaptor, which in turn returns a reference to the stored dealii::LinearAlgebra::distributed::Vector<double>.
Definition at line 15 of file rol_to_dealii_vector.cpp.
void PHiLiP::rotate_indices | ( | std::vector< unsigned int > & | numbers, |
const unsigned int | n_indices_per_direction, | ||
const char | direction, | ||
const bool | mesh_reader_verbose_output | ||
) |
Can have 4 possibilities: 1) Rotate about Z-Axis in clockwise direction 2) Rotate about X-Axis in clockwise direction 3) Rotate about Y-Axis in clockwise direction 4) Flip the entire cube, i.e.,
Rotate Z-Axis - - - -> - - - - - . . . 7-----------8 - -- - - - ^ . 3----------4 - . . - - - . . - - 6 ....... . ..... Rotate Y-Axis - - - . . - - - . 1----------2 . .- - - -> . . . . . . - - - - - . .
Rotate X-Axis
/flip to below, similar to
Rotate Z-Axis - - - -> - - - - - . . . 8-----------7 - -- - - - ^ . 4----------3 - . . - - - . . - - 5 ....... . ..... Rotate Y-Axis - - - . . - - - . 2----------1 . .- - - -> . . . . . . - - - - - . .
Rotate X-Axis
Definition at line 57 of file gmsh_reader.cpp.
std::pair< unsigned int, double > PHiLiP::solve_linear | ( | const dealii::TrilinosWrappers::SparseMatrix & | system_matrix, |
dealii::LinearAlgebra::distributed::Vector< double > & | right_hand_side, | ||
dealii::LinearAlgebra::distributed::Vector< double > & | solution, | ||
const Parameters::LinearSolverParam & | param | ||
) |
Still need to make a LinearSolver class for our problems Note that right hand side should be const however, the Trilinos wrapper gives and error when trying to map it. This is probably because the Trilinos function does not take right_hand_side as a const
Definition at line 167 of file linear_solver.cpp.