Chi-Tech
diffusion_PWLC_02c_assembleAand_b.cc
Go to the documentation of this file.
1#include "diffusion_PWLC.h"
2
5
7
8#include "chi_runtime.h"
9#include "chi_log.h"
10#include "utils/chi_timer.h"
11
12#define DefaultBCDirichlet \
13 BoundaryCondition \
14 { \
15 BCType::DIRICHLET, { 0, 0, 0 } \
16 }
17
18namespace lbs::acceleration
19{
20// ###################################################################
21/**Assembles both the matrix and the RHS using unit cell-matrices. These are
22 * the routines used in the production versions.*/
23void DiffusionPWLCSolver::AssembleAand_b(const std::vector<double>& q_vector)
24{
25 const size_t num_local_dofs = sdm_.GetNumLocalAndGhostDOFs(uk_man_);
26 ChiInvalidArgumentIf(q_vector.size() != num_local_dofs,
27 std::string("q_vector size mismatch. ") +
28 std::to_string(q_vector.size()) + " vs " +
29 std::to_string(num_local_dofs));
30
31 const std::string fname = "lbs::acceleration::DiffusionMIPSolver::"
32 "AssembleAand_b";
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.");
36 if (options.verbose)
38 << " Starting assembly";
39
40 const size_t num_groups = uk_man_.unknowns_.front().num_components_;
41
42 VecSet(rhs_, 0.0);
43 for (const auto& cell : grid_.local_cells)
44 {
45 const size_t num_faces = cell.faces_.size();
46 const auto& cell_mapping = sdm_.GetCellMapping(cell);
47 const size_t num_nodes = cell_mapping.NumNodes();
48 const auto cc_nodes = cell_mapping.GetNodeLocations();
49 const auto& unit_cell_matrices = unit_cell_matrices_[cell.local_id_];
50
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;
54
55 const auto& xs = mat_id_2_xs_map_.at(cell.material_id_);
56
57 //=========================================== Mark dirichlet nodes
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)
61 {
62 const auto& face = cell.faces_[f];
63 if (not face.has_neighbor_)
64 {
65 auto bc = DefaultBCDirichlet;
66 if (bcs_.count(face.neighbor_id_) > 0) bc = bcs_.at(face.neighbor_id_);
67
68 if (bc.type != BCType::DIRICHLET) continue;
69
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,
73 bc.values[0]};
74 }
75 }
76
77 for (size_t g = 0; g < num_groups; ++g)
78 {
79 //==================================== Get coefficient and nodal src
80 const double Dg = xs.Dg[g];
81 const double sigr_g = xs.sigR[g];
82
83 std::vector<double> qg(num_nodes, 0.0);
84 for (size_t j = 0; j < num_nodes; j++)
85 qg[j] = q_vector[sdm_.MapDOFLocal(cell, j, uk_man_, 0, g)];
86
87 //==================================== Assemble continuous terms
88 for (size_t i = 0; i < num_nodes; i++)
89 {
90 if (node_is_dirichlet[i].first) continue;
91 const int64_t imap = sdm_.MapDOF(cell, i, uk_man_, 0, g);
92 double entry_rhs_i = 0.0;
93 for (size_t j = 0; j < num_nodes; j++)
94 {
95 const int64_t jmap = sdm_.MapDOF(cell, j, uk_man_, 0, g);
96
97 const double entry_aij =
98 Dg * cell_K_matrix[i][j] + sigr_g * cell_M_matrix[i][j];
99
100 if (not node_is_dirichlet[j].first)
101 MatSetValue(A_, imap, jmap, entry_aij, ADD_VALUES);
102 else
103 {
104 const double bcvalue = node_is_dirichlet[j].second;
105 VecSetValue(rhs_, imap, -entry_aij * bcvalue, ADD_VALUES);
106 }
107
108 entry_rhs_i += qg[j] * cell_M_matrix[i][j];
109 } // for j
110
111 VecSetValue(rhs_, imap, entry_rhs_i, ADD_VALUES);
112 } // for i
113
114 //==================================== Assemble face terms
115 for (size_t f = 0; f < num_faces; ++f)
116 {
117 const auto& face = cell.faces_[f];
118 const size_t num_face_nodes = cell_mapping.NumFaceNodes(f);
119
120 const auto& face_M = unit_cell_matrices.face_M_matrices[f];
121 const auto& face_Si = unit_cell_matrices.face_Si_vectors[f];
122
123 if (not face.has_neighbor_)
124 {
125 auto bc = DefaultBCDirichlet;
126 if (bcs_.count(face.neighbor_id_) > 0)
127 bc = bcs_.at(face.neighbor_id_);
128
129 if (bc.type == BCType::DIRICHLET)
130 {
131 const double bc_value = bc.values[0];
132
133 for (size_t fi = 0; fi < num_face_nodes; ++fi)
134 {
135 const int i = cell_mapping.MapFaceNode(f, fi);
136 const int64_t imap = sdm_.MapDOF(cell, i, uk_man_, 0, g);
137
138 // MatSetValue(A_, imap, imap, cell_Vi[i], ADD_VALUES);
139 // VecSetValue(rhs_, imap, bc_value * cell_Vi[i], ADD_VALUES);
140 MatSetValue(A_, imap, imap, 1.0, ADD_VALUES);
141 VecSetValue(rhs_, imap, bc_value, ADD_VALUES);
142 } // for fi
143
144 } // Dirichlet BC
145 else if (bc.type == BCType::ROBIN)
146 {
147 const double aval = bc.values[0];
148 const double bval = bc.values[1];
149 const double fval = bc.values[2];
150
151 if (std::fabs(bval) < 1.0e-12) continue; // a and f assumed zero
152
153 for (size_t fi = 0; fi < num_face_nodes; fi++)
154 {
155 const int i = cell_mapping.MapFaceNode(f, fi);
156 const int64_t ir = sdm_.MapDOF(cell, i, uk_man_, 0, g);
157
158 if (std::fabs(aval) >= 1.0e-12)
159 {
160 for (size_t fj = 0; fj < num_face_nodes; fj++)
161 {
162 const int j = cell_mapping.MapFaceNode(f, fj);
163 const int64_t jr = sdm_.MapDOF(cell, j, uk_man_, 0, g);
164
165 const double aij = (aval / bval) * face_M[i][j];
166
167 MatSetValue(A_, ir, jr, aij, ADD_VALUES);
168 } // for fj
169 } // if a nonzero
170
171 if (std::fabs(fval) >= 1.0e-12)
172 {
173 const double rhs_val = (fval / bval) * face_Si[i];
174
175 VecSetValue(rhs_, ir, rhs_val, ADD_VALUES);
176 } // if f nonzero
177 } // for fi
178 } // Robin BC
179 } // boundary face
180 } // for face
181 } // for g
182 } // for cell
183
184 MatAssemblyBegin(A_, MAT_FINAL_ASSEMBLY);
185 MatAssemblyEnd(A_, MAT_FINAL_ASSEMBLY);
186 VecAssemblyBegin(rhs_);
187 VecAssemblyEnd(rhs_);
188
189 if (options.verbose)
190 {
191 MatInfo info;
192 MatGetInfo(A_, MAT_GLOBAL_SUM, &info);
193
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;
198 }
199
201 {
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");
206 }
207
208 KSPSetOperators(ksp_, A_, A_);
209
210 if (options.verbose)
212 << " Assembly completed";
213
214 PC pc;
215 KSPGetPC(ksp_, &pc);
216 PCSetUp(pc);
217
218 KSPSetUp(ksp_);
219}
220
221} // namespace lbs::acceleration
#define ChiInvalidArgumentIf(condition, message)
static chi::Timer program_timer
Definition: chi_runtime.h:79
static chi::ChiLog & log
Definition: chi_runtime.h:81
LogStream Log(LOG_LVL level=LOG_0)
Definition: chi_log.cc:35
std::string GetTimeString() const
Definition: chi_timer.cc:38
size_t NumNodes() const
Definition: CellMapping.cc:34
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_
Definition: diffusion.h:39
const chi_math::UnknownManager uk_man_
Definition: diffusion.h:40
const chi_mesh::MeshContinuum & grid_
Definition: diffusion.h:38
const std::map< uint64_t, BoundaryCondition > bcs_
Definition: diffusion.h:42
const MatID2XSMap mat_id_2_xs_map_
Definition: diffusion.h:44
const std::vector< UnitCellMatrices > & unit_cell_matrices_
Definition: diffusion.h:46
#define DefaultBCDirichlet
bool perform_symmetry_check
For debugging only (very expensive)
Definition: diffusion.h:63