Chi-Tech
diffusion_mip_02c_assembleAand_b.cc
Go to the documentation of this file.
1#include "diffusion_mip.h"
2#include "acceleration.h"
3
5
8
10
12
13#include "chi_runtime.h"
14#include "chi_log.h"
15#include "utils/chi_timer.h"
16#include "console/chi_console.h"
17
18#define DefaultBCDirichlet BoundaryCondition{BCType::DIRICHLET,{0,0,0}}
19
20//###################################################################
21/**Assembles both the matrix and the RHS using unit cell-matrices. These are
22 * the routines used in the production versions.*/
24 AssembleAand_b(const std::vector<double>& q_vector)
25{
26 const std::string fname = "lbs::acceleration::DiffusionMIPSolver::"
27 "AssembleAand_b";
28 if (A_ == nullptr or rhs_ == nullptr or ksp_ == nullptr)
29 throw std::logic_error(fname + ": Some or all PETSc elements are null. "
30 "Check that Initialize has been called.");
31 if (options.verbose)
32 Chi::log.Log() << Chi::program_timer.GetTimeString() << " Starting assembly";
33
34 const size_t num_groups = uk_man_.unknowns_.front().num_components_;
35
36 VecSet(rhs_, 0.0);
37 for (const auto& cell : grid_.local_cells)
38 {
39 const size_t num_faces = cell.faces_.size();
40 const auto& cell_mapping = sdm_.GetCellMapping(cell);
41 const size_t num_nodes = cell_mapping.NumNodes();
42 const auto cc_nodes = cell_mapping.GetNodeLocations();
43 const auto& unit_cell_matrices = unit_cell_matrices_[cell.local_id_];
44
45 const auto& cell_K_matrix = unit_cell_matrices.K_matrix;
46 const auto& cell_M_matrix = unit_cell_matrices.M_matrix;
47
48 const auto& xs = mat_id_2_xs_map_.at(cell.material_id_);
49
50 for (size_t g=0; g<num_groups; ++g)
51 {
52 //==================================== Get coefficient and nodal src
53 const double Dg = xs.Dg[g];
54 const double sigr_g = xs.sigR[g];
55
56 std::vector<double> qg(num_nodes, 0.0);
57 for (size_t j=0; j<num_nodes; j++)
58 qg[j] = q_vector[sdm_.MapDOFLocal(cell, j, uk_man_, 0, g)];
59
60 //==================================== Assemble continuous terms
61 for (size_t i=0; i<num_nodes; i++)
62 {
63 const int64_t imap = sdm_.MapDOF(cell, i, uk_man_, 0, g);
64 double entry_rhs_i = 0.0;
65 for (size_t j=0; j<num_nodes; j++)
66 {
67 const int64_t jmap = sdm_.MapDOF(cell, j, uk_man_, 0, g);
68
69 const double entry_aij =
70 Dg * cell_K_matrix[i][j] +
71 sigr_g * cell_M_matrix[i][j];
72
73 entry_rhs_i += qg[j]* cell_M_matrix[i][j];
74
75 MatSetValue(A_, imap, jmap, entry_aij, ADD_VALUES);
76 }//for j
77
78 VecSetValue(rhs_, imap, entry_rhs_i, ADD_VALUES);
79 }//for i
80
81 //==================================== Assemble face terms
82 for (size_t f=0; f<num_faces; ++f)
83 {
84 const auto& face = cell.faces_[f];
85 const auto& n_f = face.normal_;
86 const size_t num_face_nodes = cell_mapping.NumFaceNodes(f);
87
88 const auto& face_M = unit_cell_matrices.face_M_matrices[f];
89 const auto& face_G = unit_cell_matrices.face_G_matrices[f];
90 const auto& face_Si = unit_cell_matrices.face_Si_vectors[f];
91
92 const double hm = HPerpendicular(cell, f);
93
94 typedef chi_mesh::MeshContinuum Grid;
95
96 if (face.has_neighbor_)
97 {
98 const auto& adj_cell = grid_.cells[face.neighbor_id_];
99 const auto& adj_cell_mapping = sdm_.GetCellMapping(adj_cell);
100 const auto ac_nodes = adj_cell_mapping.GetNodeLocations();
101 const size_t acf = Grid::MapCellFace(cell, adj_cell, f);
102 const double hp = HPerpendicular(adj_cell, acf);
103
104 const auto& adj_xs = mat_id_2_xs_map_.at(adj_cell.material_id_);
105 const double adj_Dg = adj_xs.Dg[g];
106
107 //========================= Compute kappa
108 double kappa = 1.0;
109 if (cell.Type() == chi_mesh::CellType::SLAB)
110 kappa = fmax(options.penalty_factor*(adj_Dg/hp + Dg/hm)*0.5,0.25);
111 if (cell.Type() == chi_mesh::CellType::POLYGON)
112 kappa = fmax(options.penalty_factor*(adj_Dg/hp + Dg/hm)*0.5,0.25);
113 if (cell.Type() == chi_mesh::CellType::POLYHEDRON)
114 kappa = fmax(options.penalty_factor*2.0*(adj_Dg/hp + Dg/hm)*0.5,0.25);
115
116 //========================= Assembly penalty terms
117 for (size_t fi=0; fi<num_face_nodes; ++fi)
118 {
119 const int i = cell_mapping.MapFaceNode(f,fi);
120 const int64_t imap = sdm_.MapDOF(cell, i, uk_man_, 0, g);
121
122 for (size_t fj=0; fj<num_face_nodes; ++fj)
123 {
124 const int jm = cell_mapping.MapFaceNode(f,fj); //j-minus
125 const int jp = MapFaceNodeDisc(cell, adj_cell, cc_nodes, ac_nodes,
126 f, acf, fj); //j-plus
127 const int64_t jmmap = sdm_.MapDOF(cell, jm, uk_man_, 0, g);
128 const int64_t jpmap = sdm_.MapDOF(adj_cell, jp, uk_man_, 0, g);
129
130 const double aij = kappa * face_M[i][jm];
131
132 MatSetValue(A_, imap, jmmap, aij, ADD_VALUES);
133 MatSetValue(A_, imap, jpmap, -aij, ADD_VALUES);
134 }//for fj
135 }//for fi
136
137 //========================= Assemble gradient terms
138 // For the following comments we use the notation:
139 // Dk = 0.5* n dot nabla bk
140
141 // 0.5*D* n dot (b_j^+ - b_j^-)*nabla b_i^-
142 for (int i=0; i<num_nodes; i++)
143 {
144 const int64_t imap = sdm_.MapDOF(cell, i, uk_man_, 0, g);
145
146 for (int fj=0; fj<num_face_nodes; fj++)
147 {
148 const int jm = cell_mapping.MapFaceNode(f,fj); //j-minus
149 const int jp = MapFaceNodeDisc(cell, adj_cell, cc_nodes, ac_nodes,
150 f, acf, fj); //j-plus
151 const int64_t jmmap = sdm_.MapDOF(cell, jm, uk_man_, 0, g);
152 const int64_t jpmap = sdm_.MapDOF(adj_cell, jp, uk_man_, 0, g);
153
154 const double aij = -0.5*Dg*n_f.Dot(face_G[jm][i]);
155
156 MatSetValue(A_, imap, jmmap, aij, ADD_VALUES);
157 MatSetValue(A_, imap, jpmap, -aij, ADD_VALUES);
158 }//for fj
159 }//for i
160
161 // 0.5*D* n dot (b_i^+ - b_i^-)*nabla b_j^-
162 for (int fi=0; fi<num_face_nodes; fi++)
163 {
164 const int im = cell_mapping.MapFaceNode(f,fi); //i-minus
165 const int ip = MapFaceNodeDisc(cell,adj_cell,cc_nodes,ac_nodes,
166 f,acf,fi); //i-plus
167 const int64_t immap = sdm_.MapDOF(cell, im, uk_man_, 0, g);
168 const int64_t ipmap = sdm_.MapDOF(adj_cell, ip, uk_man_, 0, g);
169
170 for (int j=0; j<num_nodes; j++)
171 {
172 const int64_t jmap = sdm_.MapDOF(cell, j, uk_man_, 0, g);
173
174 const double aij = -0.5*Dg*n_f.Dot(face_G[im][j]);
175
176 MatSetValue(A_, immap, jmap, aij, ADD_VALUES);
177 MatSetValue(A_, ipmap, jmap, -aij, ADD_VALUES);
178 }//for j
179 }//for fi
180
181 }//internal face
182 else
183 {
184 auto bc = DefaultBCDirichlet;
185 if (bcs_.count(face.neighbor_id_) > 0)
186 bc = bcs_.at(face.neighbor_id_);
187
188 if (bc.type == BCType::DIRICHLET)
189 {
190 const double bc_value = bc.values[0];
191
192 //========================= Compute kappa
193 double kappa = 1.0;
194 if (cell.Type() == chi_mesh::CellType::SLAB)
195 kappa = fmax(options.penalty_factor*Dg/hm,0.25);
196 if (cell.Type() == chi_mesh::CellType::POLYGON)
197 kappa = fmax(options.penalty_factor*Dg/hm,0.25);
198 if (cell.Type() == chi_mesh::CellType::POLYHEDRON)
199 kappa = fmax(options.penalty_factor*2.0*Dg/hm,0.25);
200
201 //========================= Assembly penalty terms
202 for (size_t fi=0; fi<num_face_nodes; ++fi)
203 {
204 const int i = cell_mapping.MapFaceNode(f,fi);
205 const int64_t imap = sdm_.MapDOF(cell, i, uk_man_, 0, g);
206
207 for (size_t fj=0; fj<num_face_nodes; ++fj)
208 {
209 const int jm = cell_mapping.MapFaceNode(f,fj);
210 const int64_t jmmap = sdm_.MapDOF(cell, jm, uk_man_, 0, g);
211
212 const double aij = kappa*face_M[i][jm];
213 const double aij_bc_value = aij*bc_value;
214
215 MatSetValue(A_ , imap, jmmap, aij, ADD_VALUES);
216 VecSetValue(rhs_, imap, aij_bc_value, ADD_VALUES);
217 }//for fj
218 }//for fi
219
220 //========================= Assemble gradient terms
221 // For the following comments we use the notation:
222 // Dk = n dot nabla bk
223
224 // D* n dot (b_j^+ - b_j^-)*nabla b_i^-
225 for (size_t i=0; i<num_nodes; i++)
226 {
227 const int64_t imap = sdm_.MapDOF(cell, i, uk_man_, 0, g);
228
229 for (size_t j=0; j<num_nodes; j++)
230 {
231 const int64_t jmap = sdm_.MapDOF(cell, j, uk_man_, 0, g);
232
233 const double aij = -Dg*n_f.Dot(face_G[j][i] + face_G[i][j]);
234 const double aij_bc_value = aij*bc_value;
235
236
237 MatSetValue(A_, imap, jmap, aij, ADD_VALUES);
238 VecSetValue(rhs_, imap, aij_bc_value, ADD_VALUES);
239 }//for fj
240 }//for i
241 }//Dirichlet BC
242 else if (bc.type == BCType::ROBIN)
243 {
244 const double aval = bc.values[0];
245 const double bval = bc.values[1];
246 const double fval = bc.values[2];
247
248 if (std::fabs(bval) < 1.0e-12) continue; //a and f assumed zero
249
250 for (size_t fi=0; fi<num_face_nodes; fi++)
251 {
252 const int i = cell_mapping.MapFaceNode(f,fi);
253 const int64_t ir = sdm_.MapDOF(cell, i, uk_man_, 0, g);
254
255 if (std::fabs(aval) >= 1.0e-12)
256 {
257 for (size_t fj=0; fj<num_face_nodes; fj++)
258 {
259 const int j = cell_mapping.MapFaceNode(f,fj);
260 const int64_t jr = sdm_.MapDOF(cell, j, uk_man_, 0, g);
261
262 const double aij = (aval/bval) * face_M[i][j];
263
264 MatSetValue(A_, ir , jr, aij, ADD_VALUES);
265 }//for fj
266 }//if a nonzero
267
268 if (std::fabs(fval) >= 1.0e-12)
269 {
270 const double rhs_val = (fval/bval) * face_Si[i];
271
272 VecSetValue(rhs_, ir, rhs_val, ADD_VALUES);
273 }//if f nonzero
274 }//for fi
275 }//Robin BC
276 }//boundary face
277 }//for face
278 }//for g
279 }//for cell
280
281 MatAssemblyBegin(A_, MAT_FINAL_ASSEMBLY);
282 MatAssemblyEnd(A_, MAT_FINAL_ASSEMBLY);
283 VecAssemblyBegin(rhs_);
284 VecAssemblyEnd(rhs_);
285
286 if (options.verbose)
287 {
288 MatInfo info;
289 MatGetInfo(A_, MAT_GLOBAL_SUM, &info);
290
291 Chi::log.Log() << "Number of mallocs used = " << info.mallocs
292 << "\nNumber of non-zeros allocated = "
293 << info.nz_allocated
294 << "\nNumber of non-zeros used = "
295 << info.nz_used
296 << "\nNumber of unneeded non-zeros = "
297 << info.nz_unneeded;
298 }
299
301 {
302 PetscBool symmetry = PETSC_FALSE;
303 MatIsSymmetric(A_, 1.0e-6, &symmetry);
304 if (symmetry == PETSC_FALSE)
305 throw std::logic_error(fname + ":Symmetry check failed");
306 }
307
308 KSPSetOperators(ksp_, A_, A_);
309
310 if (options.verbose)
311 Chi::log.Log() << Chi::program_timer.GetTimeString() << " Assembly completed";
312
313 PC pc;
314 KSPGetPC(ksp_, &pc);
315 PCSetUp(pc);
316
317 KSPSetUp(ksp_);
318}
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
const std::vector< chi_mesh::Vector3 > & GetNodeLocations() const
Definition: CellMapping.cc:156
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
GlobalCellHandler cells
void AssembleAand_b(const std::vector< double > &q_vector) override
int MapFaceNodeDisc(const chi_mesh::Cell &cur_cell, const chi_mesh::Cell &adj_cell, const std::vector< chi_mesh::Vector3 > &cc_node_locs, const std::vector< chi_mesh::Vector3 > &ac_node_locs, size_t ccf, size_t acf, size_t ccfi, double epsilon=1.0e-12)
double HPerpendicular(const chi_mesh::Cell &cell, unsigned int f)
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