Chi-Tech
diffusion_mip_02b_assemble_b.cc
Go to the documentation of this file.
1#include "diffusion_mip.h"
2#include "acceleration.h"
3
5
8
10
11#include "chi_runtime.h" //TODO:Remove
12#include "chi_log.h" //TODO:Remove
13#include "utils/chi_timer.h"
14#include "console/chi_console.h"
15
16#define DefaultBCDirichlet BoundaryCondition{BCType::DIRICHLET,{0,0,0}}
17
18//###################################################################
19/**Assembles just the RHS using quadrature points. These
20 * routines exist for implementing MMS.*/
22 Assemble_b_wQpoints(const std::vector<double>& q_vector)
23{
24 const std::string fname = "lbs::acceleration::DiffusionMIPSolver::"
25 "AssembleAand_b_wQpoints";
26 if (A_ == nullptr or rhs_ == nullptr or ksp_ == nullptr)
27 throw std::logic_error(fname + ": Some or all PETSc elements are null. "
28 "Check that Initialize has been called.");
29 if (options.verbose)
30 Chi::log.Log() << Chi::program_timer.GetTimeString() << " Starting assembly";
31
32 lua_State* L = Chi::console.GetConsoleState();
33 const auto& source_function = options.source_lua_function;
34 const auto& solution_function = options.ref_solution_lua_function;
35
36 VecSet(rhs_, 0.0);
37
38 for (const auto& cell : grid_.local_cells)
39 {
40 const size_t num_faces = cell.faces_.size();
41 const auto& cell_mapping = sdm_.GetCellMapping(cell);
42 const size_t num_nodes = cell_mapping.NumNodes();
43 const auto qp_data = cell_mapping.MakeVolumetricQuadraturePointData();
44 const size_t num_groups = uk_man_.unknowns_.front().num_components_;
45
46 const auto& xs = mat_id_2_xs_map_.at(cell.material_id_);
47
48 //=========================================== For component/group
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; //entry may accumulate over j
63 if (source_function.empty())
64 for (size_t j=0; j<num_nodes; j++)
65 {
66 for (size_t qp : qp_data.QuadraturePointIndices())
67 {
68 entry_rhs_i +=
69 qg[j] *
70 qp_data.ShapeValue(i, qp) * qp_data.ShapeValue(j, qp) *
71 qp_data.JxW(qp);
72 }//for qp
73 }//for j
74 else
75 {
76 for (size_t qp : qp_data.QuadraturePointIndices())
77 entry_rhs_i +=
78 CallLuaXYZFunction(L, source_function, qp_data.QPointXYZ(qp)) *
79 qp_data.ShapeValue(i, qp) *
80 qp_data.JxW(qp);
81 }
82
83 VecSetValue(rhs_, imap, entry_rhs_i, ADD_VALUES);
84 }//for i
85
86 //==================================== Assemble face terms
87 for (size_t f=0; f<num_faces; ++f)
88 {
89 const auto& face = cell.faces_[f];
90 const auto& n_f = face.normal_;
91 const size_t num_face_nodes = cell_mapping.NumFaceNodes(f);
92 const auto fqp_data = cell_mapping.MakeSurfaceQuadraturePointData(f);
93
94 const double hm = HPerpendicular(cell, f);
95
96 if (not face.has_neighbor_)
97 {
98 auto bc = DefaultBCDirichlet;
99 if (bcs_.count(face.neighbor_id_) > 0)
100 bc = bcs_.at(face.neighbor_id_);
101
102 if (bc.type == BCType::DIRICHLET)
103 {
104 const double bc_value = bc.values[0];
105
106 //========================= Compute kappa
107 double kappa = 1.0;
108 if (cell.Type() == chi_mesh::CellType::SLAB)
109 kappa = fmax(options.penalty_factor*Dg/hm,0.25);
110 if (cell.Type() == chi_mesh::CellType::POLYGON)
111 kappa = fmax(options.penalty_factor*Dg/hm,0.25);
112 if (cell.Type() == chi_mesh::CellType::POLYHEDRON)
113 kappa = fmax(options.penalty_factor*2.0*Dg/hm,0.25);
114
115 //========================= Assembly penalty terms
116 for (size_t fi=0; fi<num_face_nodes; ++fi)
117 {
118 const int i = cell_mapping.MapFaceNode(f,fi);
119 const int64_t imap = sdm_.MapDOF(cell, i, uk_man_, 0, g);
120
121 for (size_t fj=0; fj<num_face_nodes; ++fj)
122 {
123 const int jm = cell_mapping.MapFaceNode(f,fj);
124
125 double aij = 0.0;
126 for (size_t qp : fqp_data.QuadraturePointIndices())
127 aij += kappa *
128 fqp_data.ShapeValue(i, qp) * fqp_data.ShapeValue(jm, qp) *
129 fqp_data.JxW(qp);
130 double aij_bc_value = aij*bc_value;
131
132 if (not solution_function.empty())
133 {
134 aij_bc_value = 0.0;
135 for (size_t qp : fqp_data.QuadraturePointIndices())
136 aij_bc_value +=
137 kappa *CallLuaXYZFunction(L, solution_function,
138 fqp_data.QPointXYZ(qp)) *
139 fqp_data.ShapeValue(i, qp) * fqp_data.ShapeValue(jm, qp) *
140 fqp_data.JxW(qp);
141 }
142
143 VecSetValue(rhs_, imap, aij_bc_value, ADD_VALUES);
144 }//for fj
145 }//for fi
146
147 //========================= Assemble gradient terms
148 // For the following comments we use the notation:
149 // Dk = 0.5* n dot nabla bk
150
151 // 0.5*D* n dot (b_j^+ - b_j^-)*nabla b_i^-
152 for (size_t i=0; i<num_nodes; i++)
153 {
154 const int64_t imap = sdm_.MapDOF(cell, i, uk_man_, 0, g);
155
156 for (size_t j=0; j<num_nodes; j++)
157 {
158 chi_mesh::Vector3 vec_aij;
159 for (size_t qp : fqp_data.QuadraturePointIndices())
160 vec_aij +=
161 fqp_data.ShapeValue(j, qp) * fqp_data.ShapeGrad(i, qp) *
162 fqp_data.JxW(qp) +
163 fqp_data.ShapeValue(i, qp) * fqp_data.ShapeGrad(j, qp) *
164 fqp_data.JxW(qp);
165 const double aij = -Dg*n_f.Dot(vec_aij);
166
167 double aij_bc_value = aij*bc_value;
168
169 if (not solution_function.empty())
170 {
171 chi_mesh::Vector3 vec_aij_mms;
172 for (size_t qp : fqp_data.QuadraturePointIndices())
173 vec_aij_mms +=
174 CallLuaXYZFunction(L, solution_function,
175 fqp_data.QPointXYZ(qp)) *
176 (fqp_data.ShapeValue(j, qp) * fqp_data.ShapeGrad(i, qp) *
177 fqp_data.JxW(qp) +
178 fqp_data.ShapeValue(i, qp) * fqp_data.ShapeGrad(j, qp) *
179 fqp_data.JxW(qp));
180 aij_bc_value = -Dg*n_f.Dot(vec_aij_mms);
181 }
182
183 VecSetValue(rhs_, imap, aij_bc_value, ADD_VALUES);
184 }//for fj
185 }//for i
186 }//Dirichlet BC
187 else if (bc.type == BCType::ROBIN)
188 {
189 const double bval = bc.values[1];
190 const double fval = bc.values[2];
191
192 if (std::fabs(bval) < 1.0e-12) continue; //a and f assumed zero
193
194 for (size_t fi=0; fi<num_face_nodes; fi++)
195 {
196 const int i = cell_mapping.MapFaceNode(f,fi);
197 const int64_t ir = sdm_.MapDOF(cell, i, uk_man_, 0, g);
198
199 if (std::fabs(fval) >= 1.0e-12)
200 {
201 double rhs_val = 0.0;
202 for (size_t qp : fqp_data.QuadraturePointIndices())
203 rhs_val += fqp_data.ShapeValue(i,qp) * fqp_data.JxW(qp);
204 rhs_val *= (fval/bval);
205
206 VecSetValue(rhs_, ir, rhs_val, ADD_VALUES);
207 }//if f nonzero
208 }//for fi
209 }//Robin BC
210 }//boundary face
211 }//for face
212 }//for g
213 }//for cell
214
215 VecAssemblyBegin(rhs_);
216 VecAssemblyEnd(rhs_);
217
218 KSPSetOperators(ksp_, A_, A_);
219
220 if (options.verbose)
221 Chi::log.Log() << Chi::program_timer.GetTimeString() << " Assembly completed";
222
223 PC pc;
224 KSPGetPC(ksp_, &pc);
225 PCSetUp(pc);
226
227 KSPSetUp(ksp_);
228}
static chi::Timer program_timer
Definition: chi_runtime.h:79
static chi::ChiLog & log
Definition: chi_runtime.h:81
static chi::Console & console
Definition: chi_runtime.h:80
LogStream Log(LOG_LVL level=LOG_0)
Definition: chi_log.cc:35
lua_State *& GetConsoleState()
Definition: chi_console.h:130
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
double HPerpendicular(const chi_mesh::Cell &cell, unsigned int f)
static double CallLuaXYZFunction(lua_State *L, const std::string &lua_func_name, const chi_mesh::Vector3 &xyz)
void Assemble_b_wQpoints(const std::vector< double > &q_vector)
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
#define DefaultBCDirichlet
Vector3 Dot(const chi_mesh::TensorRank2Dim3 &that) const
std::string ref_solution_lua_function
for mms
Definition: diffusion.h:66