void VectorMatrixAssembler::preAssemble( const std::size_t mesh_item_id, LocalAssemblerInterface& local_assembler, const NumLib::LocalToGlobalIndexMap& dof_table, const double t, const GlobalVector& x) { auto const indices = NumLib::getIndices(mesh_item_id, dof_table); auto const local_x = x.get(indices); local_assembler.preAssemble(t, local_x); }
Eigen::Vector3d HTProcess::getFlux(std::size_t element_id, MathLib::Point3d const& p, double const t, GlobalVector const& x) const { // fetch local_x from primary variable std::vector<GlobalIndexType> indices_cache; auto const r_c_indices = NumLib::getRowColumnIndices( element_id, *_local_to_global_index_map, indices_cache); std::vector<double> local_x(x.get(r_c_indices.rows)); return _local_assemblers[element_id]->getFlux(p, t, local_x); }
// TODO that essentially duplicates code which is also present in ProcessOutput. double getNodalValue(GlobalVector const& x, MeshLib::Mesh const& mesh, NumLib::LocalToGlobalIndexMap const& dof_table, std::size_t const node_id, std::size_t const global_component_id) { MeshLib::Location const l{mesh.getID(), MeshLib::MeshItemType::Node, node_id}; auto const index = dof_table.getLocalIndex( l, global_component_id, x.getRangeBegin(), x.getRangeEnd()); return x.get(index); }
void getVectorValues( GlobalVector const& x, NumLib::LocalToGlobalIndexMap::RowColumnIndices const& r_c_indices, std::vector<double>& local_x) { auto const& indices = r_c_indices.rows; local_x.clear(); local_x.reserve(indices.size()); for (auto i : indices) { local_x.emplace_back(x.get(i)); } }
NumLib::IterationResult TESProcess::postIterationConcreteProcess( GlobalVector const& x) { bool check_passed = true; if (!Trafo::constrained) { // bounds checking only has to happen if the vapour mass fraction is // non-logarithmic. std::vector<GlobalIndexType> indices_cache; std::vector<double> local_x_cache; std::vector<double> local_x_prev_ts_cache; auto check_variable_bounds = [&](std::size_t id, TESLocalAssemblerInterface& loc_asm) { auto const r_c_indices = NumLib::getRowColumnIndices( id, *this->_local_to_global_index_map, indices_cache); local_x_cache = x.get(r_c_indices.rows); local_x_prev_ts_cache = _x_previous_timestep->get(r_c_indices.rows); if (!loc_asm.checkBounds(local_x_cache, local_x_prev_ts_cache)) check_passed = false; }; GlobalExecutor::executeDereferenced(check_variable_bounds, _local_assemblers); } if (!check_passed) return NumLib::IterationResult::REPEAT_ITERATION; // TODO remove DBUG("ts %lu iteration %lu (in current ts: %lu) try %u accepted", _assembly_params.timestep, _assembly_params.total_iteration, _assembly_params.iteration_in_current_timestep, _assembly_params.number_of_try_of_iteration); _assembly_params.number_of_try_of_iteration = 0; return NumLib::IterationResult::SUCCESS; }
void VectorMatrixAssembler::assemble( const std::size_t mesh_item_id, LocalAssemblerInterface& local_assembler, std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>> const& dof_tables, const double t, const GlobalVector& x, GlobalMatrix& M, GlobalMatrix& K, GlobalVector& b, CoupledSolutionsForStaggeredScheme const* const cpl_xs) { std::vector<std::vector<GlobalIndexType>> indices_of_processes; indices_of_processes.reserve(dof_tables.size()); for (std::size_t i = 0; i < dof_tables.size(); i++) { indices_of_processes.emplace_back( NumLib::getIndices(mesh_item_id, dof_tables[i].get())); } auto const& indices = (cpl_xs == nullptr) ? indices_of_processes[0] : indices_of_processes[cpl_xs->process_id]; _local_M_data.clear(); _local_K_data.clear(); _local_b_data.clear(); if (cpl_xs == nullptr) { auto const local_x = x.get(indices); local_assembler.assemble(t, local_x, _local_M_data, _local_K_data, _local_b_data); } else { auto local_coupled_xs0 = getPreviousLocalSolutions(*cpl_xs, indices_of_processes); auto local_coupled_xs = getCurrentLocalSolutions(*cpl_xs, indices_of_processes); ProcessLib::LocalCoupledSolutions local_coupled_solutions( cpl_xs->dt, cpl_xs->process_id, std::move(local_coupled_xs0), std::move(local_coupled_xs)); local_assembler.assembleForStaggeredScheme(t, _local_M_data, _local_K_data, _local_b_data, local_coupled_solutions); } auto const num_r_c = indices.size(); auto const r_c_indices = NumLib::LocalToGlobalIndexMap::RowColumnIndices(indices, indices); if (!_local_M_data.empty()) { auto const local_M = MathLib::toMatrix(_local_M_data, num_r_c, num_r_c); M.add(r_c_indices, local_M); } if (!_local_K_data.empty()) { auto const local_K = MathLib::toMatrix(_local_K_data, num_r_c, num_r_c); K.add(r_c_indices, local_K); } if (!_local_b_data.empty()) { assert(_local_b_data.size() == num_r_c); b.add(indices, _local_b_data); } }
void VectorMatrixAssembler::assembleWithJacobian( std::size_t const mesh_item_id, LocalAssemblerInterface& local_assembler, std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>> const& dof_tables, const double t, GlobalVector const& x, GlobalVector const& xdot, const double dxdot_dx, const double dx_dx, GlobalMatrix& M, GlobalMatrix& K, GlobalVector& b, GlobalMatrix& Jac, CoupledSolutionsForStaggeredScheme const* const cpl_xs) { std::vector<std::vector<GlobalIndexType>> indices_of_processes; indices_of_processes.reserve(dof_tables.size()); for (std::size_t i = 0; i < dof_tables.size(); i++) { indices_of_processes.emplace_back( NumLib::getIndices(mesh_item_id, dof_tables[i].get())); } auto const& indices = (cpl_xs == nullptr) ? indices_of_processes[0] : indices_of_processes[cpl_xs->process_id]; auto const local_xdot = xdot.get(indices); _local_M_data.clear(); _local_K_data.clear(); _local_b_data.clear(); _local_Jac_data.clear(); if (cpl_xs == nullptr) { auto const local_x = x.get(indices); _jacobian_assembler->assembleWithJacobian( local_assembler, t, local_x, local_xdot, dxdot_dx, dx_dx, _local_M_data, _local_K_data, _local_b_data, _local_Jac_data); } else { auto local_coupled_xs0 = getPreviousLocalSolutions(*cpl_xs, indices_of_processes); auto local_coupled_xs = getCurrentLocalSolutions(*cpl_xs, indices_of_processes); ProcessLib::LocalCoupledSolutions local_coupled_solutions( cpl_xs->dt, cpl_xs->process_id, std::move(local_coupled_xs0), std::move(local_coupled_xs)); _jacobian_assembler->assembleWithJacobianForStaggeredScheme( local_assembler, t, local_xdot, dxdot_dx, dx_dx, _local_M_data, _local_K_data, _local_b_data, _local_Jac_data, local_coupled_solutions); } auto const num_r_c = indices.size(); auto const r_c_indices = NumLib::LocalToGlobalIndexMap::RowColumnIndices(indices, indices); if (!_local_M_data.empty()) { auto const local_M = MathLib::toMatrix(_local_M_data, num_r_c, num_r_c); M.add(r_c_indices, local_M); } if (!_local_K_data.empty()) { auto const local_K = MathLib::toMatrix(_local_K_data, num_r_c, num_r_c); K.add(r_c_indices, local_K); } if (!_local_b_data.empty()) { assert(_local_b_data.size() == num_r_c); b.add(indices, _local_b_data); } if (!_local_Jac_data.empty()) { auto const local_Jac = MathLib::toMatrix(_local_Jac_data, num_r_c, num_r_c); Jac.add(r_c_indices, local_Jac); } else { OGS_FATAL( "No Jacobian has been assembled! This might be due to programming " "errors in the local assembler of the current process."); } }