Loading [MathJax]/extensions/tex2jax.js
www.mooseframework.org
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends
MultiAppProjectionTransfer.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
11 
12 // MOOSE includes
13 #include "AddVariableAction.h"
14 #include "FEProblem.h"
15 #include "MooseMesh.h"
16 #include "MooseVariableFE.h"
17 #include "SystemBase.h"
18 
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"
29 
30 void
31 assemble_l2(EquationSystems & es, const std::string & system_name)
32 {
33  MultiAppProjectionTransfer * transfer =
34  es.parameters.get<MultiAppProjectionTransfer *>("transfer");
35  transfer->assembleL2(es, system_name);
36 }
37 
39 
41 
44 {
46  params.addClassDescription(
47  "Perform a projection between a master and sub-application mesh of a field variable.");
48 
49  MooseEnum proj_type("l2", "l2");
50  params.addParam<MooseEnum>("proj_type", proj_type, "The type of the projection.");
51 
52  params.addParam<bool>("fixed_meshes",
53  false,
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.");
57 
58  // Need one layer of ghosting
59  params.addRelationshipManager("ElementSideNeighborLayers",
62  return params;
63 }
64 
66  : MultiAppConservativeTransfer(parameters),
67  _proj_type(getParam<MooseEnum>("proj_type")),
68  _compute_matrix(true),
69  _fixed_meshes(getParam<bool>("fixed_meshes")),
70  _qps_cached(false)
71 {
72  if (_to_var_names.size() != 1)
73  paramError("variable", " Support single to-variable only ");
74 
75  if (_from_var_names.size() != 1)
76  paramError("source_variable", " Support single from-variable only ");
77 }
78 
79 void
81 {
83 
84  getAppInfo();
85 
86  _proj_sys.resize(_to_problems.size(), NULL);
87 
88  for (unsigned int i_to = 0; i_to < _to_problems.size(); i_to++)
89  {
90  FEProblemBase & to_problem = *_to_problems[i_to];
91  EquationSystems & to_es = to_problem.es();
92 
93  // Add the projection system.
94  FEType fe_type = to_problem
95  .getVariable(0,
99  .feType();
100 
101  LinearImplicitSystem & proj_sys = to_es.add_system<LinearImplicitSystem>("proj-sys-" + name());
102 
103  proj_sys.get_dof_map().add_coupling_functor(
104  proj_sys.get_dof_map().default_coupling(),
105  false); // The false keeps it from getting added to the mesh
106 
107  _proj_var_num = proj_sys.add_variable("var", fe_type);
108  proj_sys.attach_assemble_function(assemble_l2);
109  _proj_sys[i_to] = &proj_sys;
110 
111  // Prevent the projection system from being written to checkpoint
112  // files. In the event of a recover or restart, we'll read the checkpoint
113  // before this initialSetup method is called. As a result, we'll find
114  // systems in the checkpoint (the projection systems) that we don't know
115  // what to do with, and there will be a crash. We could fix this by making
116  // the systems in the constructor, except we don't know how many sub apps
117  // there are at the time of construction. So instead, we'll just nuke the
118  // projection system and rebuild it from scratch every recover/restart.
119  proj_sys.hide_output() = true;
120 
121  // Reinitialize EquationSystems since we added a system.
122  to_es.reinit();
123  }
124 
125  if (_fixed_meshes)
126  {
127  _cached_qps.resize(n_processors());
128  _cached_index_map.resize(n_processors());
129  }
130 }
131 
132 void
133 MultiAppProjectionTransfer::assembleL2(EquationSystems & es, const std::string & system_name)
134 {
135  // Get the system and mesh from the input arguments.
136  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem>(system_name);
137  MeshBase & to_mesh = es.get_mesh();
138 
139  // Get the meshfunction evaluations and the map that was stashed in the es.
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");
143 
144  // Setup system vectors and matrices.
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);
149  const DofMap & dof_map = system.get_dof_map();
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();
155 
156  for (const auto & elem : to_mesh.active_local_element_ptr_range())
157  {
158  fe->reinit(elem);
159 
160  dof_map.dof_indices(elem, dof_indices);
161  Ke.resize(dof_indices.size(), dof_indices.size());
162  Fe.resize(dof_indices.size());
163 
164  for (unsigned int qp = 0; qp < qrule.n_points(); qp++)
165  {
166  Real meshfun_eval = 0.;
167  if (element_map.find(elem->id()) != element_map.end())
168  {
169  // We have evaluations for this element.
170  meshfun_eval = final_evals[element_map[elem->id()] + qp];
171  }
172 
173  // Now compute the element matrix and RHS contributions.
174  for (unsigned int i = 0; i < phi.size(); i++)
175  {
176  // RHS
177  Fe(i) += JxW[qp] * (meshfun_eval * phi[i][qp]);
178 
179  if (_compute_matrix)
180  for (unsigned int j = 0; j < phi.size(); j++)
181  {
182  // The matrix contribution
183  Ke(i, j) += JxW[qp] * (phi[i][qp] * phi[j][qp]);
184  }
185  }
186  dof_map.constrain_element_matrix_and_vector(Ke, Fe, dof_indices);
187 
188  if (_compute_matrix)
189  system.matrix->add_matrix(Ke, dof_indices);
190  system.rhs->add_vector(Fe, dof_indices);
191  }
192  }
193 }
194 
195 void
197 {
198  _console << "Beginning projection transfer " << name() << std::endl;
199 
200  getAppInfo();
201 
203  // We are going to project the solutions by solving some linear systems. In
204  // order to assemble the systems, we need to evaluate the "from" domain
205  // solutions at quadrature points in the "to" domain. Some parallel
206  // communication is necessary because each processor doesn't necessarily have
207  // all the "from" information it needs to set its "to" values. We don't want
208  // to use a bunch of big all-to-all broadcasts, so we'll use bounding boxes to
209  // figure out which processors have the information we need and only
210  // communicate with those processors.
211  //
212  // Each processor will
213  // 1. Check its local quadrature points in the "to" domains to see which
214  // "from" domains they might be in.
215  // 2. Send quadrature points to the processors with "from" domains that might
216  // contain those points.
217  // 3. Recieve quadrature points from other processors, evaluate its mesh
218  // functions at those points, and send the values back to the proper
219  // processor
220  // 4. Recieve mesh function evaluations from all relevant processors and
221  // decide which one to use at every quadrature point (the lowest global app
222  // index always wins)
223  // 5. And use the mesh function evaluations to assemble and solve an L2
224  // projection system on its local elements.
226 
228  // For every combination of global "from" problem and local "to" problem, find
229  // which "from" bounding boxes overlap with which "to" elements. Keep track
230  // of which processors own bounding boxes that overlap with which elements.
231  // Build vectors of quadrature points to send to other processors for mesh
232  // function evaluations.
234 
235  // Get the bounding boxes for the "from" domains.
236  std::vector<BoundingBox> bboxes = getFromBoundingBoxes();
237 
238  // Figure out how many "from" domains each processor owns.
239  std::vector<unsigned int> froms_per_proc = getFromsPerProc();
240 
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(
243  n_processors());
244  // element_index_map[i_to, element_id] = index
245  // outgoing_qps[index] is the first quadrature point in element
246 
247  if (!_qps_cached)
248  {
249  for (unsigned int i_to = 0; i_to < _to_problems.size(); i_to++)
250  {
251  MeshBase & to_mesh = _to_meshes[i_to]->getMesh();
252 
253  LinearImplicitSystem & system = *_proj_sys[i_to];
254 
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();
260 
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++)
264  {
265  for (const auto & elem :
266  as_range(to_mesh.local_elements_begin(), to_mesh.local_elements_end()))
267  {
268  fe->reinit(elem);
269 
270  bool qp_hit = false;
271  for (unsigned int i_from = 0; i_from < froms_per_proc[i_proc] && !qp_hit; i_from++)
272  {
273  for (unsigned int qp = 0; qp < qrule.n_points() && !qp_hit; qp++)
274  {
275  Point qpt = xyz[qp];
276  if (bboxes[from0 + i_from].contains_point(qpt + _to_positions[i_to]))
277  qp_hit = true;
278  }
279  }
280 
281  if (qp_hit)
282  {
283  // The selected processor's bounding box contains at least one
284  // quadrature point from this element. Send all qps from this element
285  // and remember where they are in the array using the map.
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++)
289  {
290  Point qpt = xyz[qp];
291  outgoing_qps[i_proc].push_back(qpt + _to_positions[i_to]);
292  }
293  }
294  }
295  }
296  }
297 
298  if (_fixed_meshes)
299  _cached_index_map = element_index_map;
300  }
301  else
302  {
303  element_index_map = _cached_index_map;
304  }
305 
307  // Request quadrature point evaluations from other processors and handle
308  // requests sent to this processor.
310 
311  // Non-blocking send quadrature points to other processors.
312  std::vector<Parallel::Request> send_qps(n_processors());
313  if (!_qps_cached)
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]);
317 
318  // Get the local bounding boxes.
319  std::vector<BoundingBox> local_bboxes(froms_per_proc[processor_id()]);
320  {
321  // Find the index to the first of this processor's local bounding boxes.
322  unsigned int local_start = 0;
323  for (processor_id_type i_proc = 0; i_proc < n_processors() && i_proc != processor_id();
324  i_proc++)
325  local_start += froms_per_proc[i_proc];
326 
327  // Extract the local bounding boxes.
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];
330  }
331 
332  // Setup the local mesh functions.
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++)
335  {
336  FEProblemBase & from_problem = *_from_problems[i_from];
337  MooseVariableFEBase & from_var = from_problem.getVariable(
339  System & from_sys = from_var.sys().system();
340  unsigned int from_var_num = from_sys.variable_number(from_var.name());
341 
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);
345  from_func->enable_out_of_mesh_mode(OutOfMeshValue);
346  local_meshfuns[i_from] = from_func;
347  }
348 
349  // Recieve quadrature points from other processors, evaluate mesh frunctions
350  // at those points, and send the values back.
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++)
358  {
359  // Use the cached qps if they're available.
360  std::vector<Point> incoming_qps;
361  if (!_qps_cached)
362  {
363  if (i_proc == processor_id())
364  incoming_qps = outgoing_qps[i_proc];
365  else
366  _communicator.receive(i_proc, incoming_qps);
367  // Cache these qps for later if _fixed_meshes
368  if (_fixed_meshes)
369  _cached_qps[i_proc] = incoming_qps;
370  }
371  else
372  {
373  incoming_qps = _cached_qps[i_proc];
374  }
375 
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++)
380  {
381  Point qpt = incoming_qps[qp];
382 
383  // Loop until we've found the lowest-ranked app that actually contains
384  // the quadrature point.
385  for (unsigned int i_from = 0; i_from < _from_problems.size(); i_from++)
386  {
387  if (local_bboxes[i_from].contains_point(qpt))
388  {
389  outgoing_evals[i_proc][qp] = (*local_meshfuns[i_from])(qpt - _from_positions[i_from]);
391  outgoing_ids[i_proc][qp] = _local2global_map[i_from];
392  }
393  }
394  }
395 
396  if (i_proc == processor_id())
397  {
398  incoming_evals[i_proc] = outgoing_evals[i_proc];
400  incoming_app_ids[i_proc] = outgoing_ids[i_proc];
401  }
402  else
403  {
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]);
407  }
408  }
409 
411  // Gather all of the qp evaluations and pick out the best ones for each qp.
413  for (processor_id_type i_proc = 0; i_proc < n_processors(); i_proc++)
414  {
415  if (i_proc == processor_id())
416  continue;
417  _communicator.receive(i_proc, incoming_evals[i_proc]);
419  _communicator.receive(i_proc, incoming_app_ids[i_proc]);
420  }
421 
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());
424 
425  for (unsigned int i_to = 0; i_to < _to_problems.size(); i_to++)
426  {
427  MeshBase & to_mesh = _to_meshes[i_to]->getMesh();
428  LinearImplicitSystem & system = *_proj_sys[i_to];
429 
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();
435 
436  for (const auto & elem : to_mesh.active_local_element_ptr_range())
437  {
438  fe->reinit(elem);
439 
440  bool element_is_evaled = false;
441  std::vector<Real> evals(qrule.n_points(), 0.);
442 
443  for (unsigned int qp = 0; qp < qrule.n_points(); qp++)
444  {
445  Point qpt = xyz[qp];
446 
447  unsigned int lowest_app_rank = libMesh::invalid_uint;
448  for (unsigned int i_proc = 0; i_proc < n_processors(); i_proc++)
449  {
450  // Ignore the selected processor if the element wasn't found in it's
451  // bounding box.
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())
456  continue;
457  unsigned int qp0 = map[key];
458 
459  // Ignore the selected processor if it's app has a higher rank than the
460  // previously found lowest app rank.
462  if (incoming_app_ids[i_proc][qp0 + qp] >= lowest_app_rank)
463  continue;
464 
465  // Ignore the selected processor if the qp was actually outside the
466  // processor's subapp's mesh.
467  if (incoming_evals[i_proc][qp0 + qp] == OutOfMeshValue)
468  continue;
469 
470  // This is the best meshfunction evaluation so far, save it.
471  element_is_evaled = true;
472  evals[qp] = incoming_evals[i_proc][qp0 + qp];
473  }
474  }
475 
476  // If we found good evaluations for any of the qps in this element, save
477  // those evaluations for later.
478  if (element_is_evaled)
479  {
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]);
483  }
484  }
485  }
486 
488  // We now have just one or zero mesh function values at all of our local
489  // quadrature points. Stash those values (and a map linking them to element
490  // ids) in the equation systems parameters and project the solution.
492 
493  for (unsigned int i_to = 0; i_to < _to_problems.size(); i_to++)
494  {
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];
498  projectSolution(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;
501  }
502 
503  for (unsigned int i = 0; i < _from_problems.size(); i++)
504  delete local_meshfuns[i];
505 
506  // Make sure all our sends succeeded.
507  for (processor_id_type i_proc = 0; i_proc < n_processors(); i_proc++)
508  {
509  if (i_proc == processor_id())
510  continue;
511  if (!_qps_cached)
512  send_qps[i_proc].wait();
513  send_evals[i_proc].wait();
515  send_ids[i_proc].wait();
516  }
517 
518  if (_fixed_meshes)
519  _qps_cached = true;
520 
521  _console << "Finished projection transfer " << name() << std::endl;
522 
523  postExecute();
524 }
525 
526 void
528 {
529  FEProblemBase & to_problem = *_to_problems[i_to];
530  EquationSystems & proj_es = to_problem.es();
531  LinearImplicitSystem & ls = *_proj_sys[i_to];
532  // activate the current transfer
533  proj_es.parameters.set<MultiAppProjectionTransfer *>("transfer") = this;
534 
535  // TODO: specify solver params in an input file
536  // solver tolerance
537  Real tol = proj_es.parameters.get<Real>("linear solver tolerance");
538  proj_es.parameters.set<Real>("linear solver tolerance") = 1e-10; // set our tolerance
539  // solve it
540  ls.solve();
541  proj_es.parameters.set<Real>("linear solver tolerance") = tol; // restore the original tolerance
542 
543  // copy projected solution into target es
544  MeshBase & to_mesh = proj_es.get_mesh();
545 
546  MooseVariableFEBase & to_var = to_problem.getVariable(
548  System & to_sys = to_var.sys().system();
549  NumericVector<Number> * to_solution = to_sys.solution.get();
550 
551  for (const auto & node : to_mesh.local_node_ptr_range())
552  {
553  for (unsigned int comp = 0; comp < node->n_comp(to_sys.number(), to_var.number()); comp++)
554  {
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));
558  }
559  }
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++)
562  {
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));
566  }
567 
568  to_solution->close();
569  to_sys.update();
570 }
MultiAppConservativeTransfer::_from_var_names
const std::vector< VariableName > _from_var_names
Name of variables transfering from.
Definition: MultiAppConservativeTransfer.h:44
MooseVariableFEBase
Definition: MooseVariableFEBase.h:27
MultiAppProjectionTransfer::initialSetup
virtual void initialSetup() override
Method called at the beginning of the simulation for checking integrity or doing one-time setup.
Definition: MultiAppProjectionTransfer.C:80
MultiAppTransfer::_from_positions
std::vector< Point > _from_positions
Definition: MultiAppTransfer.h:70
Transfer::FROM_MULTIAPP
Definition: Transfer.h:72
MultiAppProjectionTransfer::_proj_var_num
unsigned int _proj_var_num
Having one projection variable number seems weird, but there is always one variable in every system b...
Definition: MultiAppProjectionTransfer.h:54
MultiAppConservativeTransfer
Transfers variables on possibly differne meshes while conserving a user defined property of each vari...
Definition: MultiAppConservativeTransfer.h:24
MooseMesh.h
FEProblem.h
SystemBase.h
MooseVariableBase::feType
const FEType & feType() const
Get the type of finite element object.
Definition: MooseVariableBase.h:53
Moose::RelationshipManagerType::GEOMETRIC
MultiAppProjectionTransfer::projectSolution
void projectSolution(unsigned int to_problem)
Definition: MultiAppProjectionTransfer.C:527
InputParameters::addRelationshipManager
void addRelationshipManager(const std::string &name, Moose::RelationshipManagerType rm_type, Moose::RelationshipManagerInputParameterCallback input_parameter_callback=nullptr)
Tells MOOSE about a RelationshipManager that this object needs.
Definition: InputParameters.C:344
MultiAppTransfer::_to_meshes
std::vector< MooseMesh * > _to_meshes
Definition: MultiAppTransfer.h:67
MultiAppTransfer::_to_problems
std::vector< FEProblemBase * > _to_problems
Definition: MultiAppTransfer.h:63
MooseEnum
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:31
MultiAppProjectionTransfer::_cached_qps
std::vector< std::vector< Point > > _cached_qps
Definition: MultiAppProjectionTransfer.h:61
SystemBase::system
virtual System & system()=0
Get the reference to the libMesh system.
InputParameters::addParam
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an option parameter and a documentation string to the InputParameters object.
Definition: InputParameters.h:1198
MooseObject::parameters
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseObject.h:76
MultiAppProjectionTransfer.h
MultiAppConservativeTransfer::_from_var_name
VariableName _from_var_name
This values are used if a derived class only supports one variable.
Definition: MultiAppConservativeTransfer.h:49
assemble_l2
void assemble_l2(EquationSystems &es, const std::string &system_name)
Definition: MultiAppProjectionTransfer.C:31
MultiAppProjectionTransfer::MultiAppProjectionTransfer
MultiAppProjectionTransfer(const InputParameters &parameters)
Definition: MultiAppProjectionTransfer.C:65
MultiAppProjectionTransfer::_compute_matrix
bool _compute_matrix
True, if we need to recompute the projection matrix.
Definition: MultiAppProjectionTransfer.h:50
MultiAppTransfer::getFromBoundingBoxes
std::vector< BoundingBox > getFromBoundingBoxes()
Return the bounding boxes of all the "from" domains, including all the domains not local to this proc...
Definition: MultiAppTransfer.C:166
MultiAppProjectionTransfer::_cached_index_map
std::vector< std::map< std::pair< unsigned int, unsigned int >, unsigned int > > _cached_index_map
Definition: MultiAppProjectionTransfer.h:62
Moose::RelationshipManagerType::ALGEBRAIC
MultiAppConservativeTransfer::initialSetup
virtual void initialSetup() override
Method called at the beginning of the simulation for checking integrity or doing one-time setup.
Definition: MultiAppConservativeTransfer.C:104
Transfer::_current_direction
MooseEnum _current_direction
Definition: Transfer.h:106
MooseObject::paramError
void paramError(const std::string &param, Args... args) const
Emits an error prefixed with the file and line number of the given param (from the input file) along ...
Definition: MooseObject.h:215
ConsoleStreamInterface::_console
const ConsoleStream _console
An instance of helper class to write streams to the Console objects.
Definition: ConsoleStreamInterface.h:31
MultiAppTransfer::getFromsPerProc
std::vector< unsigned int > getFromsPerProc()
Return the number of "from" domains that each processor owns.
Definition: MultiAppTransfer.C:252
InputParameters
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system.
Definition: InputParameters.h:53
MultiAppConservativeTransfer::_to_var_names
const std::vector< AuxVariableName > _to_var_names
Name of variables transfering to.
Definition: MultiAppConservativeTransfer.h:46
InputParameters::set
T & set(const std::string &name, bool quiet_mode=false)
Returns a writable reference to the named parameters.
Definition: InputParameters.h:987
MultiAppProjectionTransfer::assemble_l2
friend void assemble_l2(EquationSystems &es, const std::string &system_name)
Definition: MultiAppProjectionTransfer.C:31
Moose::VAR_ANY
Definition: MooseTypes.h:610
defineLegacyParams
defineLegacyParams(MultiAppProjectionTransfer)
Transfer::OutOfMeshValue
const static Number OutOfMeshValue
Definition: Transfer.h:113
MultiAppProjectionTransfer::_proj_sys
std::vector< LinearImplicitSystem * > _proj_sys
Definition: MultiAppProjectionTransfer.h:51
MultiAppTransfer::_to_positions
std::vector< Point > _to_positions
Definition: MultiAppTransfer.h:69
MultiAppConservativeTransfer::_to_var_name
AuxVariableName _to_var_name
Definition: MultiAppConservativeTransfer.h:50
MultiAppConservativeTransfer::validParams
static InputParameters validParams()
Definition: MultiAppConservativeTransfer.C:23
InputParameters::addClassDescription
void addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump.
Definition: InputParameters.C:70
MooseVariableFE.h
MultiAppTransfer::_local2global_map
std::vector< unsigned int > _local2global_map
Definition: MultiAppTransfer.h:100
MultiAppTransfer::getAppInfo
void getAppInfo()
This method will fill information into the convenience member variables (_to_problems,...
Definition: MultiAppTransfer.C:99
FEProblemBase::getVariable
virtual MooseVariableFEBase & getVariable(THREAD_ID tid, const std::string &var_name, Moose::VarKindType expected_var_type=Moose::VarKindType::VAR_ANY, Moose::VarFieldType expected_var_field_type=Moose::VarFieldType::VAR_FIELD_ANY) override
Returns the variable reference for requested variable which must be of the expected_var_type (Nonline...
Definition: FEProblemBase.C:4207
MooseVariableBase::sys
SystemBase & sys()
Get the system this variable is part of.
Definition: MooseVariableBase.h:58
MultiAppTransfer::_to_es
std::vector< EquationSystems * > _to_es
Definition: MultiAppTransfer.h:65
MultiAppTransfer::_from_problems
std::vector< FEProblemBase * > _from_problems
Definition: MultiAppTransfer.h:64
system
nl system()
Definition: PetscDMMoose.C:1287
FEProblemBase::es
virtual EquationSystems & es() override
Definition: FEProblemBase.h:147
MultiAppProjectionTransfer::validParams
static InputParameters validParams()
Definition: MultiAppProjectionTransfer.C:43
MultiAppConservativeTransfer::postExecute
virtual void postExecute()
Add some extra work if necessary after execute().
Definition: MultiAppConservativeTransfer.C:148
registerMooseObject
registerMooseObject("MooseApp", MultiAppProjectionTransfer)
MultiAppProjectionTransfer::assembleL2
void assembleL2(EquationSystems &es, const std::string &system_name)
Definition: MultiAppProjectionTransfer.C:133
MultiAppProjectionTransfer::execute
virtual void execute() override
Execute the transfer.
Definition: MultiAppProjectionTransfer.C:196
dof_map
DofMap & dof_map
Definition: PetscDMMoose.C:1503
AddVariableAction.h
FEProblemBase
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
Definition: FEProblemBase.h:139
MultiAppProjectionTransfer
Project values from one domain to another.
Definition: MultiAppProjectionTransfer.h:28
Moose::VAR_FIELD_STANDARD
Definition: MooseTypes.h:615
MooseObject::name
virtual const std::string & name() const
Get the name of the object.
Definition: MooseObject.h:70
MultiAppProjectionTransfer::_qps_cached
bool _qps_cached
Definition: MultiAppProjectionTransfer.h:60
MooseVariableBase::name
const std::string & name() const override
Get the variable name.
Definition: MooseVariableBase.h:63
MultiAppProjectionTransfer::_fixed_meshes
bool _fixed_meshes
Definition: MultiAppProjectionTransfer.h:59
MooseVariableBase::number
unsigned int number() const
Get variable number coming from libMesh.
Definition: MooseVariableBase.h:48