12#define DefaultBCDirichlet \
15 BCType::DIRICHLET, { 0, 0, 0 } \
27 std::string(
"q_vector size mismatch. ") +
28 std::to_string(q_vector.size()) +
" vs " +
29 std::to_string(num_local_dofs));
31 const std::string fname =
"lbs::acceleration::DiffusionMIPSolver::"
33 if (
A_ ==
nullptr or
rhs_ ==
nullptr or
ksp_ ==
nullptr)
34 throw std::logic_error(fname +
": Some or all PETSc elements are null. "
35 "Check that Initialize has been called.");
38 <<
" Starting assembly";
45 const size_t num_faces = cell.faces_.size();
47 const size_t num_nodes = cell_mapping.
NumNodes();
48 const auto cc_nodes = cell_mapping.GetNodeLocations();
51 const auto& cell_K_matrix = unit_cell_matrices.K_matrix;
52 const auto& cell_M_matrix = unit_cell_matrices.M_matrix;
53 const auto& cell_Vi = unit_cell_matrices.Vi_vectors;
58 typedef std::pair<bool, double> DirichFlagVal;
59 std::vector<DirichFlagVal> node_is_dirichlet(num_nodes, {
false, 0.0});
60 for (
size_t f = 0; f < num_faces; ++f)
62 const auto& face = cell.faces_[f];
63 if (not face.has_neighbor_)
66 if (
bcs_.count(face.neighbor_id_) > 0) bc =
bcs_.at(face.neighbor_id_);
70 const size_t num_face_nodes = cell_mapping.NumFaceNodes(f);
71 for (
size_t fi = 0; fi < num_face_nodes; ++fi)
72 node_is_dirichlet[cell_mapping.MapFaceNode(f, fi)] = {
true,
77 for (
size_t g = 0; g < num_groups; ++g)
80 const double Dg = xs.Dg[g];
81 const double sigr_g = xs.sigR[g];
83 std::vector<double> qg(num_nodes, 0.0);
84 for (
size_t j = 0; j < num_nodes; j++)
88 for (
size_t i = 0; i < num_nodes; i++)
90 if (node_is_dirichlet[i].first)
continue;
92 double entry_rhs_i = 0.0;
93 for (
size_t j = 0; j < num_nodes; j++)
97 const double entry_aij =
98 Dg * cell_K_matrix[i][j] + sigr_g * cell_M_matrix[i][j];
100 if (not node_is_dirichlet[j].first)
101 MatSetValue(
A_, imap, jmap, entry_aij, ADD_VALUES);
104 const double bcvalue = node_is_dirichlet[j].second;
105 VecSetValue(
rhs_, imap, -entry_aij * bcvalue, ADD_VALUES);
108 entry_rhs_i += qg[j] * cell_M_matrix[i][j];
111 VecSetValue(
rhs_, imap, entry_rhs_i, ADD_VALUES);
115 for (
size_t f = 0; f < num_faces; ++f)
117 const auto& face = cell.faces_[f];
118 const size_t num_face_nodes = cell_mapping.NumFaceNodes(f);
120 const auto& face_M = unit_cell_matrices.face_M_matrices[f];
121 const auto& face_Si = unit_cell_matrices.face_Si_vectors[f];
123 if (not face.has_neighbor_)
126 if (
bcs_.count(face.neighbor_id_) > 0)
127 bc =
bcs_.at(face.neighbor_id_);
131 const double bc_value = bc.values[0];
133 for (
size_t fi = 0; fi < num_face_nodes; ++fi)
135 const int i = cell_mapping.MapFaceNode(f, fi);
140 MatSetValue(
A_, imap, imap, 1.0, ADD_VALUES);
141 VecSetValue(
rhs_, imap, bc_value, ADD_VALUES);
147 const double aval = bc.values[0];
148 const double bval = bc.values[1];
149 const double fval = bc.values[2];
151 if (std::fabs(bval) < 1.0e-12)
continue;
153 for (
size_t fi = 0; fi < num_face_nodes; fi++)
155 const int i = cell_mapping.MapFaceNode(f, fi);
158 if (std::fabs(aval) >= 1.0e-12)
160 for (
size_t fj = 0; fj < num_face_nodes; fj++)
162 const int j = cell_mapping.MapFaceNode(f, fj);
165 const double aij = (aval / bval) * face_M[i][j];
167 MatSetValue(
A_, ir, jr, aij, ADD_VALUES);
171 if (std::fabs(fval) >= 1.0e-12)
173 const double rhs_val = (fval / bval) * face_Si[i];
175 VecSetValue(
rhs_, ir, rhs_val, ADD_VALUES);
184 MatAssemblyBegin(
A_, MAT_FINAL_ASSEMBLY);
185 MatAssemblyEnd(
A_, MAT_FINAL_ASSEMBLY);
186 VecAssemblyBegin(
rhs_);
187 VecAssemblyEnd(
rhs_);
192 MatGetInfo(
A_, MAT_GLOBAL_SUM, &info);
194 Chi::log.
Log() <<
"Number of mallocs used = " << info.mallocs
195 <<
"\nNumber of non-zeros allocated = " << info.nz_allocated
196 <<
"\nNumber of non-zeros used = " << info.nz_used
197 <<
"\nNumber of unneeded non-zeros = " << info.nz_unneeded;
202 PetscBool symmetry = PETSC_FALSE;
203 MatIsSymmetric(
A_, 1.0e-6, &symmetry);
204 if (symmetry == PETSC_FALSE)
205 throw std::logic_error(fname +
":Symmetry check failed");
212 <<
" Assembly completed";
#define ChiInvalidArgumentIf(condition, message)
static chi::Timer program_timer
LogStream Log(LOG_LVL level=LOG_0)
std::string GetTimeString() const
size_t GetNumLocalAndGhostDOFs(const UnknownManager &unknown_manager) const
const CellMapping & GetCellMapping(const chi_mesh::Cell &cell) const
virtual int64_t MapDOF(const chi_mesh::Cell &cell, unsigned int node, const UnknownManager &unknown_manager, unsigned int unknown_id, unsigned int component) const =0
virtual int64_t MapDOFLocal(const chi_mesh::Cell &cell, unsigned int node, const UnknownManager &unknown_manager, unsigned int unknown_id, unsigned int component) const =0
std::vector< Unknown > unknowns_
LocalCellHandler local_cells
void AssembleAand_b(const std::vector< double > &q_vector) override
struct lbs::acceleration::DiffusionSolver::Options options
const chi_math::SpatialDiscretization & sdm_
const chi_math::UnknownManager uk_man_
const chi_mesh::MeshContinuum & grid_
const std::map< uint64_t, BoundaryCondition > bcs_
const MatID2XSMap mat_id_2_xs_map_
const std::vector< UnitCellMatrices > & unit_cell_matrices_
#define DefaultBCDirichlet
bool verbose
Verbosity flag.
bool perform_symmetry_check
For debugging only (very expensive)