19 #include "libmesh/dof_map.h"
20 #include "libmesh/linear_implicit_system.h"
21 #include "libmesh/mesh_function.h"
22 #include "libmesh/mesh_tools.h"
23 #include "libmesh/numeric_vector.h"
24 #include "libmesh/parallel_algebra.h"
25 #include "libmesh/quadrature_gauss.h"
26 #include "libmesh/sparse_matrix.h"
27 #include "libmesh/string_to_enum.h"
28 #include "libmesh/default_coupling.h"
31 assemble_l2(EquationSystems & es,
const std::string & system_name)
47 "Perform a projection between a master and sub-application mesh of a field variable.");
50 params.
addParam<
MooseEnum>(
"proj_type", proj_type,
"The type of the projection.");
52 params.
addParam<
bool>(
"fixed_meshes",
54 "Set to true when the meshes are not changing (ie, "
55 "no movement or adaptivity). This will cache some "
56 "information to speed up the transfer.");
67 _proj_type(getParam<
MooseEnum>(
"proj_type")),
68 _compute_matrix(true),
69 _fixed_meshes(getParam<bool>(
"fixed_meshes")),
73 paramError(
"variable",
" Support single to-variable only ");
76 paramError(
"source_variable",
" Support single from-variable only ");
88 for (
unsigned int i_to = 0; i_to <
_to_problems.size(); i_to++)
91 EquationSystems & to_es = to_problem.
es();
94 FEType fe_type = to_problem
101 LinearImplicitSystem & proj_sys = to_es.add_system<LinearImplicitSystem>(
"proj-sys-" +
name());
103 proj_sys.get_dof_map().add_coupling_functor(
104 proj_sys.get_dof_map().default_coupling(),
119 proj_sys.hide_output() =
true;
136 LinearImplicitSystem &
system = es.get_system<LinearImplicitSystem>(system_name);
137 MeshBase & to_mesh = es.get_mesh();
140 std::vector<Real> & final_evals = *es.parameters.get<std::vector<Real> *>(
"final_evals");
141 std::map<dof_id_type, unsigned int> & element_map =
142 *es.parameters.get<std::map<dof_id_type, unsigned int> *>(
"element_map");
145 FEType fe_type =
system.variable_type(0);
146 std::unique_ptr<FEBase> fe(FEBase::build(to_mesh.mesh_dimension(), fe_type));
147 QGauss qrule(to_mesh.mesh_dimension(), fe_type.default_quadrature_order());
148 fe->attach_quadrature_rule(&qrule);
150 DenseMatrix<Number> Ke;
151 DenseVector<Number> Fe;
152 std::vector<dof_id_type> dof_indices;
153 const std::vector<Real> & JxW = fe->get_JxW();
154 const std::vector<std::vector<Real>> & phi = fe->get_phi();
156 for (
const auto & elem : to_mesh.active_local_element_ptr_range())
160 dof_map.dof_indices(elem, dof_indices);
161 Ke.resize(dof_indices.size(), dof_indices.size());
162 Fe.resize(dof_indices.size());
164 for (
unsigned int qp = 0; qp < qrule.n_points(); qp++)
166 Real meshfun_eval = 0.;
167 if (element_map.find(elem->id()) != element_map.end())
170 meshfun_eval = final_evals[element_map[elem->id()] + qp];
174 for (
unsigned int i = 0; i < phi.size(); i++)
177 Fe(i) += JxW[qp] * (meshfun_eval * phi[i][qp]);
180 for (
unsigned int j = 0; j < phi.size(); j++)
183 Ke(i, j) += JxW[qp] * (phi[i][qp] * phi[j][qp]);
186 dof_map.constrain_element_matrix_and_vector(Ke, Fe, dof_indices);
189 system.matrix->add_matrix(Ke, dof_indices);
190 system.rhs->add_vector(Fe, dof_indices);
198 _console <<
"Beginning projection transfer " <<
name() << std::endl;
241 std::vector<std::vector<Point>> outgoing_qps(n_processors());
242 std::vector<std::map<std::pair<unsigned int, unsigned int>,
unsigned int>> element_index_map(
249 for (
unsigned int i_to = 0; i_to <
_to_problems.size(); i_to++)
251 MeshBase & to_mesh =
_to_meshes[i_to]->getMesh();
255 FEType fe_type =
system.variable_type(0);
256 std::unique_ptr<FEBase> fe(FEBase::build(to_mesh.mesh_dimension(), fe_type));
257 QGauss qrule(to_mesh.mesh_dimension(), fe_type.default_quadrature_order());
258 fe->attach_quadrature_rule(&qrule);
259 const std::vector<Point> & xyz = fe->get_xyz();
261 unsigned int from0 = 0;
262 for (processor_id_type i_proc = 0; i_proc < n_processors();
263 from0 += froms_per_proc[i_proc], i_proc++)
265 for (
const auto & elem :
266 as_range(to_mesh.local_elements_begin(), to_mesh.local_elements_end()))
271 for (
unsigned int i_from = 0; i_from < froms_per_proc[i_proc] && !qp_hit; i_from++)
273 for (
unsigned int qp = 0; qp < qrule.n_points() && !qp_hit; qp++)
276 if (bboxes[from0 + i_from].contains_point(qpt +
_to_positions[i_to]))
286 std::pair<unsigned int, unsigned int> key(i_to, elem->id());
287 element_index_map[i_proc][key] = outgoing_qps[i_proc].size();
288 for (
unsigned int qp = 0; qp < qrule.n_points(); qp++)
312 std::vector<Parallel::Request> send_qps(n_processors());
314 for (processor_id_type i_proc = 0; i_proc < n_processors(); i_proc++)
315 if (i_proc != processor_id())
316 _communicator.send(i_proc, outgoing_qps[i_proc], send_qps[i_proc]);
319 std::vector<BoundingBox> local_bboxes(froms_per_proc[processor_id()]);
322 unsigned int local_start = 0;
323 for (processor_id_type i_proc = 0; i_proc < n_processors() && i_proc != processor_id();
325 local_start += froms_per_proc[i_proc];
328 for (
unsigned int i_from = 0; i_from < froms_per_proc[processor_id()]; i_from++)
329 local_bboxes[i_from] = bboxes[local_start + i_from];
333 std::vector<MeshFunction *> local_meshfuns(froms_per_proc[processor_id()], NULL);
334 for (
unsigned int i_from = 0; i_from <
_from_problems.size(); i_from++)
339 System & from_sys = from_var.
sys().
system();
340 unsigned int from_var_num = from_sys.variable_number(from_var.
name());
342 MeshFunction * from_func =
new MeshFunction(
343 from_problem.
es(), *from_sys.current_local_solution, from_sys.get_dof_map(), from_var_num);
344 from_func->init(Trees::ELEMENTS);
346 local_meshfuns[i_from] = from_func;
351 std::vector<Parallel::Request> send_evals(n_processors());
352 std::vector<Parallel::Request> send_ids(n_processors());
353 std::vector<std::vector<Real>> outgoing_evals(n_processors());
354 std::vector<std::vector<unsigned int>> outgoing_ids(n_processors());
355 std::vector<std::vector<Real>> incoming_evals(n_processors());
356 std::vector<std::vector<unsigned int>> incoming_app_ids(n_processors());
357 for (processor_id_type i_proc = 0; i_proc < n_processors(); i_proc++)
360 std::vector<Point> incoming_qps;
363 if (i_proc == processor_id())
364 incoming_qps = outgoing_qps[i_proc];
366 _communicator.receive(i_proc, incoming_qps);
376 outgoing_evals[i_proc].resize(incoming_qps.size(),
OutOfMeshValue);
378 outgoing_ids[i_proc].resize(incoming_qps.size(), libMesh::invalid_uint);
379 for (
unsigned int qp = 0; qp < incoming_qps.size(); qp++)
381 Point qpt = incoming_qps[qp];
385 for (
unsigned int i_from = 0; i_from <
_from_problems.size(); i_from++)
387 if (local_bboxes[i_from].contains_point(qpt))
389 outgoing_evals[i_proc][qp] = (*local_meshfuns[i_from])(qpt -
_from_positions[i_from]);
396 if (i_proc == processor_id())
398 incoming_evals[i_proc] = outgoing_evals[i_proc];
400 incoming_app_ids[i_proc] = outgoing_ids[i_proc];
404 _communicator.send(i_proc, outgoing_evals[i_proc], send_evals[i_proc]);
406 _communicator.send(i_proc, outgoing_ids[i_proc], send_ids[i_proc]);
413 for (processor_id_type i_proc = 0; i_proc < n_processors(); i_proc++)
415 if (i_proc == processor_id())
417 _communicator.receive(i_proc, incoming_evals[i_proc]);
419 _communicator.receive(i_proc, incoming_app_ids[i_proc]);
422 std::vector<std::vector<Real>> final_evals(
_to_problems.size());
423 std::vector<std::map<dof_id_type, unsigned int>> trimmed_element_maps(
_to_problems.size());
425 for (
unsigned int i_to = 0; i_to <
_to_problems.size(); i_to++)
427 MeshBase & to_mesh =
_to_meshes[i_to]->getMesh();
430 FEType fe_type =
system.variable_type(0);
431 std::unique_ptr<FEBase> fe(FEBase::build(to_mesh.mesh_dimension(), fe_type));
432 QGauss qrule(to_mesh.mesh_dimension(), fe_type.default_quadrature_order());
433 fe->attach_quadrature_rule(&qrule);
434 const std::vector<Point> & xyz = fe->get_xyz();
436 for (
const auto & elem : to_mesh.active_local_element_ptr_range())
440 bool element_is_evaled =
false;
441 std::vector<Real> evals(qrule.n_points(), 0.);
443 for (
unsigned int qp = 0; qp < qrule.n_points(); qp++)
447 unsigned int lowest_app_rank = libMesh::invalid_uint;
448 for (
unsigned int i_proc = 0; i_proc < n_processors(); i_proc++)
452 std::map<std::pair<unsigned int, unsigned int>,
unsigned int> & map =
453 element_index_map[i_proc];
454 std::pair<unsigned int, unsigned int> key(i_to, elem->id());
455 if (map.find(key) == map.end())
457 unsigned int qp0 = map[key];
462 if (incoming_app_ids[i_proc][qp0 + qp] >= lowest_app_rank)
471 element_is_evaled =
true;
472 evals[qp] = incoming_evals[i_proc][qp0 + qp];
478 if (element_is_evaled)
480 trimmed_element_maps[i_to][elem->id()] = final_evals[i_to].size();
481 for (
unsigned int qp = 0; qp < qrule.n_points(); qp++)
482 final_evals[i_to].push_back(evals[qp]);
493 for (
unsigned int i_to = 0; i_to <
_to_problems.size(); i_to++)
495 _to_es[i_to]->parameters.set<std::vector<Real> *>(
"final_evals") = &final_evals[i_to];
496 _to_es[i_to]->parameters.set<std::map<dof_id_type, unsigned int> *>(
"element_map") =
497 &trimmed_element_maps[i_to];
499 _to_es[i_to]->parameters.set<std::vector<Real> *>(
"final_evals") = NULL;
500 _to_es[i_to]->parameters.set<std::map<dof_id_type, unsigned int> *>(
"element_map") = NULL;
504 delete local_meshfuns[i];
507 for (processor_id_type i_proc = 0; i_proc < n_processors(); i_proc++)
509 if (i_proc == processor_id())
512 send_qps[i_proc].wait();
513 send_evals[i_proc].wait();
515 send_ids[i_proc].wait();
521 _console <<
"Finished projection transfer " <<
name() << std::endl;
530 EquationSystems & proj_es = to_problem.
es();
531 LinearImplicitSystem & ls = *
_proj_sys[i_to];
537 Real tol = proj_es.
parameters.get<Real>(
"linear solver tolerance");
538 proj_es.parameters.
set<Real>(
"linear solver tolerance") = 1e-10;
541 proj_es.parameters.set<Real>(
"linear solver tolerance") = tol;
544 MeshBase & to_mesh = proj_es.get_mesh();
549 NumericVector<Number> * to_solution = to_sys.solution.get();
551 for (
const auto & node : to_mesh.local_node_ptr_range())
553 for (
unsigned int comp = 0; comp < node->n_comp(to_sys.number(), to_var.
number()); comp++)
555 const dof_id_type proj_index = node->dof_number(ls.number(),
_proj_var_num, comp);
556 const dof_id_type to_index = node->dof_number(to_sys.number(), to_var.
number(), comp);
557 to_solution->set(to_index, (*ls.solution)(proj_index));
560 for (
const auto & elem : to_mesh.active_local_element_ptr_range())
561 for (
unsigned int comp = 0; comp < elem->n_comp(to_sys.number(), to_var.
number()); comp++)
563 const dof_id_type proj_index = elem->dof_number(ls.number(),
_proj_var_num, comp);
564 const dof_id_type to_index = elem->dof_number(to_sys.number(), to_var.
number(), comp);
565 to_solution->set(to_index, (*ls.solution)(proj_index));
568 to_solution->close();