PolyFEM
|
#include <OperatorSplittingSolver.hpp>
Public Member Functions | |
void | initialize_grid (const mesh::Mesh &mesh, const std::vector< basis::ElementBases > &gbases, const std::vector< basis::ElementBases > &bases, const double &density_dx) |
void | initialize_mesh (const mesh::Mesh &mesh, const int shape, const int n_el, const std::vector< mesh::LocalBoundary > &local_boundary) |
void | initialize_hashtable (const mesh::Mesh &mesh) |
OperatorSplittingSolver () | |
void | initialize_solver (const mesh::Mesh &mesh, const int shape_, const int n_el_, const std::vector< mesh::LocalBoundary > &local_boundary, const std::vector< int > &bnd_nodes) |
OperatorSplittingSolver (const mesh::Mesh &mesh, const int shape, const int n_el, const std::vector< mesh::LocalBoundary > &local_boundary, const std::vector< int > &bnd_nodes) | |
OperatorSplittingSolver (const mesh::Mesh &mesh, const int shape, const int n_el, const std::vector< mesh::LocalBoundary > &local_boundary, const std::vector< int > &boundary_nodes_, const std::vector< int > &pressure_boundary_nodes, const std::vector< int > &bnd_nodes, const StiffnessMatrix &mass, const StiffnessMatrix &stiffness_viscosity, const StiffnessMatrix &stiffness_velocity, const StiffnessMatrix &mass_velocity, const double &dt, const double &viscosity_, const json ¶ms) | |
int | handle_boundary_advection (RowVectorNd &pos) |
int | trace_back (const std::vector< basis::ElementBases > &gbases, const std::vector< basis::ElementBases > &bases, const RowVectorNd &pos_1, const RowVectorNd &vel_1, RowVectorNd &pos_2, RowVectorNd &vel_2, Eigen::MatrixXd &local_pos, const Eigen::MatrixXd &sol, const double dt) |
int | interpolator (const std::vector< basis::ElementBases > &gbases, const std::vector< basis::ElementBases > &bases, const RowVectorNd &pos, RowVectorNd &vel, Eigen::MatrixXd &local_pos, const Eigen::MatrixXd &sol) |
void | interpolator (const RowVectorNd &pos, double &val) |
void | advection (const mesh::Mesh &mesh, const std::vector< basis::ElementBases > &gbases, const std::vector< basis::ElementBases > &bases, Eigen::MatrixXd &sol, const double dt, const Eigen::MatrixXd &local_pts, const int order=1, const int RK=1) |
void | advect_density_exact (const std::vector< basis::ElementBases > &gbases, const std::vector< basis::ElementBases > &bases, const std::shared_ptr< assembler::Problem > problem, const double t, const double dt, const int RK=3) |
void | advect_density (const std::vector< basis::ElementBases > &gbases, const std::vector< basis::ElementBases > &bases, const Eigen::MatrixXd &sol, const double dt, const int RK=3) |
void | advection_FLIP (const mesh::Mesh &mesh, const std::vector< basis::ElementBases > &gbases, const std::vector< basis::ElementBases > &bases, Eigen::MatrixXd &sol, const double dt, const Eigen::MatrixXd &local_pts, const int order=1) |
void | advection_PIC (const mesh::Mesh &mesh, const std::vector< basis::ElementBases > &gbases, const std::vector< basis::ElementBases > &bases, Eigen::MatrixXd &sol, const double dt, const Eigen::MatrixXd &local_pts, const int order=1) |
void | solve_diffusion_1st (const StiffnessMatrix &mass, const std::vector< int > &bnd_nodes, Eigen::MatrixXd &sol) |
void | external_force (const mesh::Mesh &mesh, const assembler::Assembler &assembler, const std::vector< basis::ElementBases > &gbases, const std::vector< basis::ElementBases > &bases, const double dt, Eigen::MatrixXd &sol, const Eigen::MatrixXd &local_pts, const std::shared_ptr< assembler::Problem > problem, const double time) |
void | solve_pressure (const StiffnessMatrix &mixed_stiffness, const std::vector< int > &pressure_boundary_nodes, Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure) |
void | projection (const StiffnessMatrix &velocity_mass, const StiffnessMatrix &mixed_stiffness, const std::vector< int > &boundary_nodes_, Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure) |
void | projection (int n_bases, const std::vector< basis::ElementBases > &gbases, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &pressure_bases, const Eigen::MatrixXd &local_pts, Eigen::MatrixXd &pressure, Eigen::MatrixXd &sol) |
void | initialize_density (const std::shared_ptr< assembler::Problem > &problem) |
long | search_cell (const std::vector< basis::ElementBases > &gbases, const RowVectorNd &pos, Eigen::MatrixXd &local_pts) |
bool | outside_quad (const std::vector< RowVectorNd > &vert, const RowVectorNd &pos) |
void | compute_gbase_val (const int elem_idx, const Eigen::MatrixXd &local_pos, Eigen::MatrixXd &pos) |
void | compute_gbase_jacobi (const int elem_idx, const Eigen::MatrixXd &local_pos, Eigen::MatrixXd &jacobi) |
void | calculate_local_pts (const basis::ElementBases &gbase, const int elem_idx, const RowVectorNd &pos, Eigen::MatrixXd &local_pos) |
void | save_density () |
Public Attributes | |
int | dim |
int | n_el |
int | shape |
RowVectorNd | min_domain |
RowVectorNd | max_domain |
Eigen::MatrixXd | V |
Eigen::MatrixXi | T |
std::vector< std::vector< long > > | hash_table |
Eigen::Matrix< long, Eigen::Dynamic, 1, Eigen::ColMajor, 3, 1 > | hash_table_cell_num |
std::vector< Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > > | position_particle |
std::vector< Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > > | velocity_particle |
std::vector< int > | cellI_particle |
Eigen::MatrixXd | new_sol |
Eigen::MatrixXd | new_sol_w |
std::vector< int > | boundary_elem_id |
std::vector< int > | boundary_nodes |
std::unique_ptr< polysolve::linear::Solver > | solver_diffusion |
std::unique_ptr< polysolve::linear::Solver > | solver_projection |
std::unique_ptr< polysolve::linear::Solver > | solver_mass |
StiffnessMatrix | mat_diffusion |
StiffnessMatrix | mat_projection |
Eigen::VectorXd | density |
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > | grid_cell_num |
double | resolution |
Definition at line 21 of file OperatorSplittingSolver.hpp.
|
inline |
Definition at line 35 of file OperatorSplittingSolver.hpp.
polyfem::solver::OperatorSplittingSolver::OperatorSplittingSolver | ( | const mesh::Mesh & | mesh, |
const int | shape, | ||
const int | n_el, | ||
const std::vector< mesh::LocalBoundary > & | local_boundary, | ||
const std::vector< int > & | bnd_nodes | ||
) |
Definition at line 153 of file OperatorSplittingSolver.cpp.
References initialize_solver(), n_el, and shape.
polyfem::solver::OperatorSplittingSolver::OperatorSplittingSolver | ( | const mesh::Mesh & | mesh, |
const int | shape, | ||
const int | n_el, | ||
const std::vector< mesh::LocalBoundary > & | local_boundary, | ||
const std::vector< int > & | boundary_nodes_, | ||
const std::vector< int > & | pressure_boundary_nodes, | ||
const std::vector< int > & | bnd_nodes, | ||
const StiffnessMatrix & | mass, | ||
const StiffnessMatrix & | stiffness_viscosity, | ||
const StiffnessMatrix & | stiffness_velocity, | ||
const StiffnessMatrix & | mass_velocity, | ||
const double & | dt, | ||
const double & | viscosity_, | ||
const json & | params | ||
) |
Definition at line 161 of file OperatorSplittingSolver.cpp.
References initialize_solver(), polyfem::logger(), mat_diffusion, mat_projection, n_el, shape, solver_diffusion, solver_mass, solver_projection, and val.
void polyfem::solver::OperatorSplittingSolver::advect_density | ( | const std::vector< basis::ElementBases > & | gbases, |
const std::vector< basis::ElementBases > & | bases, | ||
const Eigen::MatrixXd & | sol, | ||
const double | dt, | ||
const int | RK = 3 |
||
) |
Definition at line 481 of file OperatorSplittingSolver.cpp.
References density, dim, grid_cell_num, interpolator(), min_domain, and resolution.
void polyfem::solver::OperatorSplittingSolver::advect_density_exact | ( | const std::vector< basis::ElementBases > & | gbases, |
const std::vector< basis::ElementBases > & | bases, | ||
const std::shared_ptr< assembler::Problem > | problem, | ||
const double | t, | ||
const double | dt, | ||
const int | RK = 3 |
||
) |
Definition at line 408 of file OperatorSplittingSolver.cpp.
References density, dim, grid_cell_num, interpolator(), min_domain, and resolution.
void polyfem::solver::OperatorSplittingSolver::advection | ( | const mesh::Mesh & | mesh, |
const std::vector< basis::ElementBases > & | gbases, | ||
const std::vector< basis::ElementBases > & | bases, | ||
Eigen::MatrixXd & | sol, | ||
const double | dt, | ||
const Eigen::MatrixXd & | local_pts, | ||
const int | order = 1 , |
||
const int | RK = 1 |
||
) |
Definition at line 354 of file OperatorSplittingSolver.cpp.
References dim, interpolator(), n_el, and new_sol.
Referenced by polyfem::State::solve_transient_navier_stokes_split().
void polyfem::solver::OperatorSplittingSolver::advection_FLIP | ( | const mesh::Mesh & | mesh, |
const std::vector< basis::ElementBases > & | gbases, | ||
const std::vector< basis::ElementBases > & | bases, | ||
Eigen::MatrixXd & | sol, | ||
const double | dt, | ||
const Eigen::MatrixXd & | local_pts, | ||
const int | order = 1 |
||
) |
Definition at line 554 of file OperatorSplittingSolver.cpp.
References polyfem::assembler::ElementAssemblyValues::basis_values, calculate_local_pts(), cellI_particle, polyfem::assembler::ElementAssemblyValues::compute(), dim, n_el, new_sol, new_sol_w, position_particle, shape, trace_back(), vals, velocity_particle, and x.
Referenced by polyfem::State::solve_transient_navier_stokes_split().
void polyfem::solver::OperatorSplittingSolver::advection_PIC | ( | const mesh::Mesh & | mesh, |
const std::vector< basis::ElementBases > & | gbases, | ||
const std::vector< basis::ElementBases > & | bases, | ||
Eigen::MatrixXd & | sol, | ||
const double | dt, | ||
const Eigen::MatrixXd & | local_pts, | ||
const int | order = 1 |
||
) |
Definition at line 789 of file OperatorSplittingSolver.cpp.
References polyfem::assembler::ElementAssemblyValues::basis_values, polyfem::mesh::Mesh::cell_vertex(), cellI_particle, polyfem::assembler::ElementAssemblyValues::compute(), dim, gvals, n_el, new_sol, new_sol_w, polyfem::mesh::Mesh::point(), position_particle, shape, trace_back(), vals, and velocity_particle.
void polyfem::solver::OperatorSplittingSolver::calculate_local_pts | ( | const basis::ElementBases & | gbase, |
const int | elem_idx, | ||
const RowVectorNd & | pos, | ||
Eigen::MatrixXd & | local_pos | ||
) |
Definition at line 1273 of file OperatorSplittingSolver.cpp.
References dim, polyfem::basis::ElementBases::eval_geom_mapping(), polyfem::basis::ElementBases::eval_geom_mapping_grads(), and shape.
Referenced by advection_FLIP(), interpolator(), and search_cell().
void polyfem::solver::OperatorSplittingSolver::compute_gbase_jacobi | ( | const int | elem_idx, |
const Eigen::MatrixXd & | local_pos, | ||
Eigen::MatrixXd & | jacobi | ||
) |
void polyfem::solver::OperatorSplittingSolver::compute_gbase_val | ( | const int | elem_idx, |
const Eigen::MatrixXd & | local_pos, | ||
Eigen::MatrixXd & | pos | ||
) |
void polyfem::solver::OperatorSplittingSolver::external_force | ( | const mesh::Mesh & | mesh, |
const assembler::Assembler & | assembler, | ||
const std::vector< basis::ElementBases > & | gbases, | ||
const std::vector< basis::ElementBases > & | bases, | ||
const double | dt, | ||
Eigen::MatrixXd & | sol, | ||
const Eigen::MatrixXd & | local_pts, | ||
const std::shared_ptr< assembler::Problem > | problem, | ||
const double | time | ||
) |
Definition at line 947 of file OperatorSplittingSolver.cpp.
References dim, n_el, and val.
Referenced by polyfem::State::solve_transient_navier_stokes_split().
int polyfem::solver::OperatorSplittingSolver::handle_boundary_advection | ( | RowVectorNd & | pos | ) |
Definition at line 233 of file OperatorSplittingSolver.cpp.
References boundary_elem_id, dim, shape, T, and V.
Referenced by interpolator().
void polyfem::solver::OperatorSplittingSolver::initialize_density | ( | const std::shared_ptr< assembler::Problem > & | problem | ) |
Definition at line 1089 of file OperatorSplittingSolver.cpp.
References density, dim, grid_cell_num, min_domain, and resolution.
void polyfem::solver::OperatorSplittingSolver::initialize_grid | ( | const mesh::Mesh & | mesh, |
const std::vector< basis::ElementBases > & | gbases, | ||
const std::vector< basis::ElementBases > & | bases, | ||
const double & | density_dx | ||
) |
Definition at line 14 of file OperatorSplittingSolver.cpp.
References density, dim, grid_cell_num, max_domain, min_domain, and resolution.
void polyfem::solver::OperatorSplittingSolver::initialize_hashtable | ( | const mesh::Mesh & | mesh | ) |
Definition at line 67 of file OperatorSplittingSolver.cpp.
References dim, polyfem::mesh::Mesh::get_edges(), hash_table, hash_table_cell_num, polyfem::logger(), max_domain, min_domain, T, V, x, y, and z.
Referenced by initialize_solver().
void polyfem::solver::OperatorSplittingSolver::initialize_mesh | ( | const mesh::Mesh & | mesh, |
const int | shape, | ||
const int | n_el, | ||
const std::vector< mesh::LocalBoundary > & | local_boundary | ||
) |
Definition at line 32 of file OperatorSplittingSolver.cpp.
References boundary_elem_id, polyfem::mesh::Mesh::bounding_box(), polyfem::mesh::Mesh::cell_vertex(), dim, polyfem::mesh::Mesh::dimension(), max_domain, min_domain, n_el, polyfem::mesh::Mesh::n_vertices(), polyfem::mesh::Mesh::point(), shape, T, and V.
Referenced by initialize_solver().
void polyfem::solver::OperatorSplittingSolver::initialize_solver | ( | const mesh::Mesh & | mesh, |
const int | shape_, | ||
const int | n_el_, | ||
const std::vector< mesh::LocalBoundary > & | local_boundary, | ||
const std::vector< int > & | bnd_nodes | ||
) |
Definition at line 140 of file OperatorSplittingSolver.cpp.
References boundary_nodes, initialize_hashtable(), initialize_mesh(), n_el, and shape.
Referenced by OperatorSplittingSolver(), and OperatorSplittingSolver().
void polyfem::solver::OperatorSplittingSolver::interpolator | ( | const RowVectorNd & | pos, |
double & | val | ||
) |
Definition at line 320 of file OperatorSplittingSolver.cpp.
References density, dim, grid_cell_num, min_domain, resolution, and val.
int polyfem::solver::OperatorSplittingSolver::interpolator | ( | const std::vector< basis::ElementBases > & | gbases, |
const std::vector< basis::ElementBases > & | bases, | ||
const RowVectorNd & | pos, | ||
RowVectorNd & | vel, | ||
Eigen::MatrixXd & | local_pos, | ||
const Eigen::MatrixXd & | sol | ||
) |
Definition at line 285 of file OperatorSplittingSolver.cpp.
References calculate_local_pts(), polyfem::assembler::ElementAssemblyValues::compute(), dim, handle_boundary_advection(), search_cell(), and vals.
Referenced by advect_density(), advect_density_exact(), advection(), and trace_back().
bool polyfem::solver::OperatorSplittingSolver::outside_quad | ( | const std::vector< RowVectorNd > & | vert, |
const RowVectorNd & | pos | ||
) |
Definition at line 1157 of file OperatorSplittingSolver.cpp.
void polyfem::solver::OperatorSplittingSolver::projection | ( | const StiffnessMatrix & | velocity_mass, |
const StiffnessMatrix & | mixed_stiffness, | ||
const std::vector< int > & | boundary_nodes_, | ||
Eigen::MatrixXd & | sol, | ||
const Eigen::MatrixXd & | pressure | ||
) |
Definition at line 1028 of file OperatorSplittingSolver.cpp.
References solver_diffusion, and solver_mass.
Referenced by polyfem::State::solve_transient_navier_stokes_split().
void polyfem::solver::OperatorSplittingSolver::projection | ( | int | n_bases, |
const std::vector< basis::ElementBases > & | gbases, | ||
const std::vector< basis::ElementBases > & | bases, | ||
const std::vector< basis::ElementBases > & | pressure_bases, | ||
const Eigen::MatrixXd & | local_pts, | ||
Eigen::MatrixXd & | pressure, | ||
Eigen::MatrixXd & | sol | ||
) |
Definition at line 1051 of file OperatorSplittingSolver.cpp.
References polyfem::assembler::ElementAssemblyValues::compute(), dim, n_el, and vals.
void polyfem::solver::OperatorSplittingSolver::save_density | ( | ) |
Definition at line 1318 of file OperatorSplittingSolver.cpp.
References density, dim, and grid_cell_num.
long polyfem::solver::OperatorSplittingSolver::search_cell | ( | const std::vector< basis::ElementBases > & | gbases, |
const RowVectorNd & | pos, | ||
Eigen::MatrixXd & | local_pts | ||
) |
Definition at line 1119 of file OperatorSplittingSolver.cpp.
References calculate_local_pts(), dim, hash_table, hash_table_cell_num, max_domain, min_domain, and shape.
Referenced by interpolator().
void polyfem::solver::OperatorSplittingSolver::solve_diffusion_1st | ( | const StiffnessMatrix & | mass, |
const std::vector< int > & | bnd_nodes, | ||
Eigen::MatrixXd & | sol | ||
) |
Definition at line 912 of file OperatorSplittingSolver.cpp.
References dim, mat_diffusion, solver_diffusion, and x.
Referenced by polyfem::State::solve_transient_navier_stokes_split().
void polyfem::solver::OperatorSplittingSolver::solve_pressure | ( | const StiffnessMatrix & | mixed_stiffness, |
const std::vector< int > & | pressure_boundary_nodes, | ||
Eigen::MatrixXd & | sol, | ||
Eigen::MatrixXd & | pressure | ||
) |
Definition at line 988 of file OperatorSplittingSolver.cpp.
References mat_projection, solver_projection, and x.
Referenced by polyfem::State::solve_transient_navier_stokes_split().
int polyfem::solver::OperatorSplittingSolver::trace_back | ( | const std::vector< basis::ElementBases > & | gbases, |
const std::vector< basis::ElementBases > & | bases, | ||
const RowVectorNd & | pos_1, | ||
const RowVectorNd & | vel_1, | ||
RowVectorNd & | pos_2, | ||
RowVectorNd & | vel_2, | ||
Eigen::MatrixXd & | local_pos, | ||
const Eigen::MatrixXd & | sol, | ||
const double | dt | ||
) |
Definition at line 270 of file OperatorSplittingSolver.cpp.
References interpolator().
Referenced by advection_FLIP(), and advection_PIC().
std::vector<int> polyfem::solver::OperatorSplittingSolver::boundary_elem_id |
Definition at line 169 of file OperatorSplittingSolver.hpp.
Referenced by handle_boundary_advection(), and initialize_mesh().
std::vector<int> polyfem::solver::OperatorSplittingSolver::boundary_nodes |
Definition at line 170 of file OperatorSplittingSolver.hpp.
Referenced by initialize_solver().
std::vector<int> polyfem::solver::OperatorSplittingSolver::cellI_particle |
Definition at line 165 of file OperatorSplittingSolver.hpp.
Referenced by advection_FLIP(), and advection_PIC().
Eigen::VectorXd polyfem::solver::OperatorSplittingSolver::density |
Definition at line 179 of file OperatorSplittingSolver.hpp.
Referenced by advect_density(), advect_density_exact(), initialize_density(), initialize_grid(), interpolator(), and save_density().
int polyfem::solver::OperatorSplittingSolver::dim |
Definition at line 150 of file OperatorSplittingSolver.hpp.
Referenced by advect_density(), advect_density_exact(), advection(), advection_FLIP(), advection_PIC(), calculate_local_pts(), compute_gbase_jacobi(), compute_gbase_val(), external_force(), handle_boundary_advection(), initialize_density(), initialize_grid(), initialize_hashtable(), initialize_mesh(), interpolator(), interpolator(), projection(), save_density(), search_cell(), and solve_diffusion_1st().
Eigen::Matrix<double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3> polyfem::solver::OperatorSplittingSolver::grid_cell_num |
Definition at line 182 of file OperatorSplittingSolver.hpp.
Referenced by advect_density(), advect_density_exact(), initialize_density(), initialize_grid(), interpolator(), and save_density().
std::vector<std::vector<long> > polyfem::solver::OperatorSplittingSolver::hash_table |
Definition at line 160 of file OperatorSplittingSolver.hpp.
Referenced by initialize_hashtable(), and search_cell().
Eigen::Matrix<long, Eigen::Dynamic, 1, Eigen::ColMajor, 3, 1> polyfem::solver::OperatorSplittingSolver::hash_table_cell_num |
Definition at line 161 of file OperatorSplittingSolver.hpp.
Referenced by initialize_hashtable(), and search_cell().
StiffnessMatrix polyfem::solver::OperatorSplittingSolver::mat_diffusion |
Definition at line 176 of file OperatorSplittingSolver.hpp.
Referenced by OperatorSplittingSolver(), and solve_diffusion_1st().
StiffnessMatrix polyfem::solver::OperatorSplittingSolver::mat_projection |
Definition at line 177 of file OperatorSplittingSolver.hpp.
Referenced by OperatorSplittingSolver(), and solve_pressure().
RowVectorNd polyfem::solver::OperatorSplittingSolver::max_domain |
Definition at line 155 of file OperatorSplittingSolver.hpp.
Referenced by initialize_grid(), initialize_hashtable(), initialize_mesh(), and search_cell().
RowVectorNd polyfem::solver::OperatorSplittingSolver::min_domain |
Definition at line 154 of file OperatorSplittingSolver.hpp.
Referenced by advect_density(), advect_density_exact(), initialize_density(), initialize_grid(), initialize_hashtable(), initialize_mesh(), interpolator(), and search_cell().
int polyfem::solver::OperatorSplittingSolver::n_el |
Definition at line 151 of file OperatorSplittingSolver.hpp.
Referenced by advection(), advection_FLIP(), advection_PIC(), external_force(), initialize_mesh(), initialize_solver(), OperatorSplittingSolver(), OperatorSplittingSolver(), and projection().
Eigen::MatrixXd polyfem::solver::OperatorSplittingSolver::new_sol |
Definition at line 166 of file OperatorSplittingSolver.hpp.
Referenced by advection(), advection_FLIP(), and advection_PIC().
Eigen::MatrixXd polyfem::solver::OperatorSplittingSolver::new_sol_w |
Definition at line 167 of file OperatorSplittingSolver.hpp.
Referenced by advection_FLIP(), and advection_PIC().
std::vector<Eigen::Matrix<double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3> > polyfem::solver::OperatorSplittingSolver::position_particle |
Definition at line 163 of file OperatorSplittingSolver.hpp.
Referenced by advection_FLIP(), and advection_PIC().
double polyfem::solver::OperatorSplittingSolver::resolution |
Definition at line 183 of file OperatorSplittingSolver.hpp.
Referenced by advect_density(), advect_density_exact(), initialize_density(), initialize_grid(), and interpolator().
int polyfem::solver::OperatorSplittingSolver::shape |
Definition at line 152 of file OperatorSplittingSolver.hpp.
Referenced by advection_FLIP(), advection_PIC(), calculate_local_pts(), compute_gbase_jacobi(), compute_gbase_val(), handle_boundary_advection(), initialize_mesh(), initialize_solver(), OperatorSplittingSolver(), OperatorSplittingSolver(), and search_cell().
std::unique_ptr<polysolve::linear::Solver> polyfem::solver::OperatorSplittingSolver::solver_diffusion |
Definition at line 172 of file OperatorSplittingSolver.hpp.
Referenced by OperatorSplittingSolver(), projection(), and solve_diffusion_1st().
std::unique_ptr<polysolve::linear::Solver> polyfem::solver::OperatorSplittingSolver::solver_mass |
Definition at line 174 of file OperatorSplittingSolver.hpp.
Referenced by OperatorSplittingSolver(), and projection().
std::unique_ptr<polysolve::linear::Solver> polyfem::solver::OperatorSplittingSolver::solver_projection |
Definition at line 173 of file OperatorSplittingSolver.hpp.
Referenced by OperatorSplittingSolver(), and solve_pressure().
Eigen::MatrixXi polyfem::solver::OperatorSplittingSolver::T |
Definition at line 158 of file OperatorSplittingSolver.hpp.
Referenced by compute_gbase_jacobi(), compute_gbase_val(), handle_boundary_advection(), initialize_hashtable(), and initialize_mesh().
Eigen::MatrixXd polyfem::solver::OperatorSplittingSolver::V |
Definition at line 157 of file OperatorSplittingSolver.hpp.
Referenced by compute_gbase_jacobi(), compute_gbase_val(), handle_boundary_advection(), initialize_hashtable(), and initialize_mesh().
std::vector<Eigen::Matrix<double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3> > polyfem::solver::OperatorSplittingSolver::velocity_particle |
Definition at line 164 of file OperatorSplittingSolver.hpp.
Referenced by advection_FLIP(), and advection_PIC().