Chi-Tech
diffusion_mip_02d_assemble_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 the RHS using unit cell-matrices. These are
22 * the routines used in the production versions.*/
24 Assemble_b(const std::vector<double>& q_vector)
25{
26 const std::string fname = "lbs::acceleration::DiffusionMIPSolver::"
27 "Assemble_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_M_matrix = unit_cell_matrices.M_matrix;
46
47 const auto& xs = mat_id_2_xs_map_.at(cell.material_id_);
48
49 for (size_t g=0; g<num_groups; ++g)
50 {
51 //==================================== Get coefficient and nodal src
52 const double Dg = xs.Dg[g];
53
54 std::vector<double> qg(num_nodes, 0.0);
55 for (size_t j=0; j<num_nodes; j++)
56 qg[j] = q_vector[sdm_.MapDOFLocal(cell, j, uk_man_, 0, g)];
57
58 //==================================== Assemble continuous terms
59 for (size_t i=0; i<num_nodes; i++)
60 {
61 const int64_t imap = sdm_.MapDOF(cell, i, uk_man_, 0, g);
62 double entry_rhs_i = 0.0;
63 for (size_t j=0; j<num_nodes; j++)
64 entry_rhs_i += qg[j]* cell_M_matrix[i][j];
65
66 VecSetValue(rhs_, imap, entry_rhs_i, ADD_VALUES);
67 }//for i
68
69 //==================================== Assemble face terms
70 for (size_t f=0; f<num_faces; ++f)
71 {
72 const auto& face = cell.faces_[f];
73 const auto& n_f = face.normal_;
74 const size_t num_face_nodes = cell_mapping.NumFaceNodes(f);
75
76 const auto& face_M = unit_cell_matrices.face_M_matrices[f];
77 const auto& face_G = unit_cell_matrices.face_G_matrices[f];
78 const auto& face_Si = unit_cell_matrices.face_Si_vectors[f];
79
80 const double hm = HPerpendicular(cell, f);
81
82 if (not face.has_neighbor_)
83 {
84 auto bc = DefaultBCDirichlet;
85 if (bcs_.count(face.neighbor_id_) > 0)
86 bc = bcs_.at(face.neighbor_id_);
87
88 if (bc.type == BCType::DIRICHLET)
89 {
90 const double bc_value = bc.values[0];
91
92 //========================= Compute kappa
93 double kappa = 1.0;
94 if (cell.Type() == chi_mesh::CellType::SLAB)
95 kappa = fmax(options.penalty_factor*Dg/hm,0.25);
96 if (cell.Type() == chi_mesh::CellType::POLYGON)
97 kappa = fmax(options.penalty_factor*Dg/hm,0.25);
98 if (cell.Type() == chi_mesh::CellType::POLYHEDRON)
99 kappa = fmax(options.penalty_factor*2.0*Dg/hm,0.25);
100
101 //========================= Assembly penalty terms
102 for (size_t fi=0; fi<num_face_nodes; ++fi)
103 {
104 const int i = cell_mapping.MapFaceNode(f,fi);
105 const int64_t imap = sdm_.MapDOF(cell, i, uk_man_, 0, g);
106
107 for (size_t fj=0; fj<num_face_nodes; ++fj)
108 {
109 const int jm = cell_mapping.MapFaceNode(f,fj);
110
111 const double aij = kappa*face_M[i][jm];
112 const double aij_bc_value = aij*bc_value;
113
114 VecSetValue(rhs_, imap, aij_bc_value, ADD_VALUES);
115 }//for fj
116 }//for fi
117
118 //========================= Assemble gradient terms
119 // For the following comments we use the notation:
120 // Dk = n dot nabla bk
121
122 // D* n dot (b_j^+ - b_j^-)*nabla b_i^-
123 for (size_t i=0; i<num_nodes; i++)
124 {
125 const int64_t imap = sdm_.MapDOF(cell, i, uk_man_, 0, g);
126
127 for (size_t j=0; j<num_nodes; j++)
128 {
129 const double aij = -Dg*n_f.Dot(face_G[j][i] + face_G[i][j]);
130 const double aij_bc_value = aij*bc_value;
131
132 VecSetValue(rhs_, imap, aij_bc_value, ADD_VALUES);
133 }//for fj
134 }//for i
135 }//Dirichlet BC
136 else if (bc.type == BCType::ROBIN)
137 {
138 const double bval = bc.values[1];
139 const double fval = bc.values[2];
140
141 if (std::fabs(bval) < 1.0e-12) continue; //a and f assumed zero
142
143 for (size_t fi=0; fi<num_face_nodes; fi++)
144 {
145 const int i = cell_mapping.MapFaceNode(f,fi);
146 const int64_t ir = sdm_.MapDOF(cell, i, uk_man_, 0, g);
147
148 if (std::fabs(fval) >= 1.0e-12)
149 {
150 const double rhs_val = (fval/bval) * face_Si[i];
151
152 VecSetValue(rhs_, ir, rhs_val, ADD_VALUES);
153 }//if f nonzero
154 }//for fi
155 }//Robin BC
156 }//boundary face
157 }//for face
158 }//for g
159 }//for cell
160
161 VecAssemblyBegin(rhs_);
162 VecAssemblyEnd(rhs_);
163
164 if (options.verbose)
165 Chi::log.Log() << Chi::program_timer.GetTimeString() << " Assembly completed";
166}
167
168//###################################################################
169/**Assembles the RHS using unit cell-matrices. These are
170 * the routines used in the production versions.*/
172Assemble_b(Vec petsc_q_vector)
173{
174 const std::string fname = "lbs::acceleration::DiffusionMIPSolver::"
175 "Assemble_b";
176 if (A_ == nullptr or rhs_ == nullptr or ksp_ == nullptr)
177 throw std::logic_error(fname + ": Some or all PETSc elements are null. "
178 "Check that Initialize has been called.");
179 if (options.verbose)
180 Chi::log.Log() << Chi::program_timer.GetTimeString() << " Starting assembly";
181
182 const size_t num_groups = uk_man_.unknowns_.front().num_components_;
183
184 const double* q_vector;
185 VecGetArrayRead(petsc_q_vector, &q_vector);
186
187 VecSet(rhs_, 0.0);
188 for (const auto& cell : grid_.local_cells)
189 {
190 const size_t num_faces = cell.faces_.size();
191 const auto& cell_mapping = sdm_.GetCellMapping(cell);
192 const size_t num_nodes = cell_mapping.NumNodes();
193 const auto cc_nodes = cell_mapping.GetNodeLocations();
194 const auto& unit_cell_matrices = unit_cell_matrices_[cell.local_id_];
195
196 const auto& cell_M_matrix = unit_cell_matrices.M_matrix;
197
198 const auto& xs = mat_id_2_xs_map_.at(cell.material_id_);
199
200 for (size_t g=0; g<num_groups; ++g)
201 {
202 //==================================== Get coefficient and nodal src
203 const double Dg = xs.Dg[g];
204
205 std::vector<double> qg(num_nodes, 0.0);
206 for (size_t j=0; j<num_nodes; j++)
207 qg[j] = q_vector[sdm_.MapDOFLocal(cell, j, uk_man_, 0, g)];
208
209 //==================================== Assemble continuous terms
210 for (size_t i=0; i<num_nodes; i++)
211 {
212 const int64_t imap = sdm_.MapDOF(cell, i, uk_man_, 0, g);
213 double entry_rhs_i = 0.0;
214 for (size_t j=0; j<num_nodes; j++)
215 entry_rhs_i += qg[j]* cell_M_matrix[i][j];
216
217 VecSetValue(rhs_, imap, entry_rhs_i, ADD_VALUES);
218 }//for i
219
220 //==================================== Assemble face terms
221 for (size_t f=0; f<num_faces; ++f)
222 {
223 const auto& face = cell.faces_[f];
224 const auto& n_f = face.normal_;
225 const size_t num_face_nodes = cell_mapping.NumFaceNodes(f);
226
227 const auto& face_M = unit_cell_matrices.face_M_matrices[f];
228 const auto& face_G = unit_cell_matrices.face_G_matrices[f];
229 const auto& face_Si = unit_cell_matrices.face_Si_vectors[f];
230
231 const double hm = HPerpendicular(cell, f);
232
233 if (not face.has_neighbor_)
234 {
235 auto bc = DefaultBCDirichlet;
236 if (bcs_.count(face.neighbor_id_) > 0)
237 bc = bcs_.at(face.neighbor_id_);
238
239 if (bc.type == BCType::DIRICHLET)
240 {
241 const double bc_value = bc.values[0];
242
243 //========================= Compute kappa
244 double kappa = 1.0;
245 if (cell.Type() == chi_mesh::CellType::SLAB)
246 kappa = fmax(options.penalty_factor*Dg/hm,0.25);
247 if (cell.Type() == chi_mesh::CellType::POLYGON)
248 kappa = fmax(options.penalty_factor*Dg/hm,0.25);
249 if (cell.Type() == chi_mesh::CellType::POLYHEDRON)
250 kappa = fmax(options.penalty_factor*2.0*Dg/hm,0.25);
251
252 //========================= Assembly penalty terms
253 for (size_t fi=0; fi<num_face_nodes; ++fi)
254 {
255 const int i = cell_mapping.MapFaceNode(f,fi);
256 const int64_t imap = sdm_.MapDOF(cell, i, uk_man_, 0, g);
257
258 for (size_t fj=0; fj<num_face_nodes; ++fj)
259 {
260 const int jm = cell_mapping.MapFaceNode(f,fj);
261
262 const double aij = kappa*face_M[i][jm];
263 const double aij_bc_value = aij*bc_value;
264
265 VecSetValue(rhs_, imap, aij_bc_value, ADD_VALUES);
266 }//for fj
267 }//for fi
268
269 //========================= Assemble gradient terms
270 // For the following comments we use the notation:
271 // Dk = n dot nabla bk
272
273 // D* n dot (b_j^+ - b_j^-)*nabla b_i^-
274 for (size_t i=0; i<num_nodes; i++)
275 {
276 const int64_t imap = sdm_.MapDOF(cell, i, uk_man_, 0, g);
277
278 for (size_t j=0; j<num_nodes; j++)
279 {
280 const double aij = -Dg*n_f.Dot(face_G[j][i] + face_G[i][j]);
281 const double aij_bc_value = aij*bc_value;
282
283 VecSetValue(rhs_, imap, aij_bc_value, ADD_VALUES);
284 }//for fj
285 }//for i
286 }//Dirichlet BC
287 else if (bc.type == BCType::ROBIN)
288 {
289 const double bval = bc.values[1];
290 const double fval = bc.values[2];
291
292 if (std::fabs(bval) < 1.0e-12) continue; //a and f assumed zero
293
294 for (size_t fi=0; fi<num_face_nodes; fi++)
295 {
296 const int i = cell_mapping.MapFaceNode(f,fi);
297 const int64_t ir = sdm_.MapDOF(cell, i, uk_man_, 0, g);
298
299 if (std::fabs(fval) >= 1.0e-12)
300 {
301 const double rhs_val = (fval/bval) * face_Si[i];
302
303 VecSetValue(rhs_, ir, rhs_val, ADD_VALUES);
304 }//if f nonzero
305 }//for fi
306 }//Robin BC
307 }//boundary face
308 }//for face
309 }//for g
310 }//for cell
311
312 VecRestoreArrayRead(petsc_q_vector, &q_vector);
313
314 VecAssemblyBegin(rhs_);
315 VecAssemblyEnd(rhs_);
316
317 if (options.verbose)
318 Chi::log.Log() << Chi::program_timer.GetTimeString() << " Assembly completed";
319}
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 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 Assemble_b(const std::vector< double > &q_vector) override
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
struct _p_Vec * Vec