LCOV - code coverage report
Current view: top level - src/transfers - MultiAppGeneralFieldTransfer.C (source / functions) Hit Total Coverage
Test: idaholab/moose framework: 2bf808 Lines: 801 912 87.8 %
Date: 2025-07-17 01:28:37 Functions: 40 41 97.6 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : //* This file is part of the MOOSE framework
       2             : //* https://mooseframework.inl.gov
       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             : 
      10             : #include "MultiAppGeneralFieldTransfer.h"
      11             : 
      12             : // MOOSE includes
      13             : #include "DisplacedProblem.h"
      14             : #include "FEProblem.h"
      15             : #include "MooseMesh.h"
      16             : #include "MooseTypes.h"
      17             : #include "MooseVariableFE.h"
      18             : #include "MeshDivision.h"
      19             : #include "Positions.h"
      20             : #include "Factory.h"
      21             : #include "MooseAppCoordTransform.h"
      22             : 
      23             : // libmesh includes
      24             : #include "libmesh/point_locator_base.h"
      25             : #include "libmesh/enum_point_locator_type.h"
      26             : 
      27             : // TIMPI includes
      28             : #include "timpi/communicator.h"
      29             : #include "timpi/parallel_sync.h"
      30             : 
      31             : using namespace libMesh;
      32             : 
      33             : namespace GeneralFieldTransfer
      34             : {
      35             : Number BetterOutOfMeshValue = std::numeric_limits<Real>::infinity();
      36             : }
      37             : 
      38             : InputParameters
      39       70392 : MultiAppGeneralFieldTransfer::validParams()
      40             : {
      41       70392 :   InputParameters params = MultiAppConservativeTransfer::validParams();
      42             :   // Expansion default to make a node in the target mesh overlapping with a node in the origin
      43             :   // mesh always register as being inside the origin application bounding box. The contains_point
      44             :   // bounding box checks uses exact comparisons
      45      211176 :   params.addRangeCheckedParam<Real>("bbox_factor",
      46      140784 :                                     1.00000001,
      47             :                                     "bbox_factor>0",
      48             :                                     "Factor to inflate or deflate the source app bounding boxes");
      49       70392 :   params.addRangeCheckedParam<std::vector<Real>>(
      50             :       "fixed_bounding_box_size",
      51             :       "fixed_bounding_box_size >= 0",
      52             :       "Override source app bounding box size(s) for searches. App bounding boxes will still be  "
      53             :       "centered on the same coordinates. Only non-zero components passed will override.");
      54      211176 :   params.addParam<Real>(
      55             :       "extrapolation_constant",
      56      140784 :       0,
      57             :       "Constant to use when no source app can provide a valid value for a target location.");
      58             : 
      59             :   // Block restrictions
      60       70392 :   params.addParam<std::vector<SubdomainName>>(
      61             :       "from_blocks",
      62             :       "Subdomain restriction to transfer from (defaults to all the origin app domain)");
      63       70392 :   params.addParam<std::vector<SubdomainName>>(
      64             :       "to_blocks", "Subdomain restriction to transfer to, (defaults to all the target app domain)");
      65             : 
      66             :   // Boundary restrictions
      67       70392 :   params.addParam<std::vector<BoundaryName>>(
      68             :       "from_boundaries",
      69             :       "The boundary we are transferring from (if not specified, whole domain is used).");
      70       70392 :   params.addParam<std::vector<BoundaryName>>(
      71             :       "to_boundaries",
      72             :       "The boundary we are transferring to (if not specified, whole domain is used).");
      73       70392 :   MooseEnum nodes_or_sides("nodes sides", "sides");
      74       70392 :   params.addParam<MooseEnum>("elemental_boundary_restriction",
      75             :                              nodes_or_sides,
      76             :                              "Whether elemental variable boundary restriction is considered by "
      77             :                              "element side or element nodes");
      78             : 
      79             :   // Mesh division restriction
      80       70392 :   params.addParam<MeshDivisionName>("from_mesh_division",
      81             :                                     "Mesh division object on the origin application");
      82       70392 :   params.addParam<MeshDivisionName>("to_mesh_division",
      83             :                                     "Mesh division object on the target application");
      84             :   MooseEnum mesh_division_uses("spatial_restriction matching_division matching_subapp_index none",
      85       70392 :                                "none");
      86       70392 :   params.addParam<MooseEnum>("from_mesh_division_usage",
      87             :                              mesh_division_uses,
      88             :                              "How to use the source mesh division in the transfer. See object "
      89             :                              "documentation for description of each option");
      90       70392 :   params.addParam<MooseEnum>("to_mesh_division_usage",
      91             :                              mesh_division_uses,
      92             :                              "How to use the target mesh division in the transfer. See object "
      93             :                              "documentation for description of each option");
      94             : 
      95             :   // Array and vector variables
      96      211176 :   params.addParam<std::vector<unsigned int>>("source_variable_components",
      97      140784 :                                              std::vector<unsigned int>(),
      98             :                                              "The source array or vector variable component(s).");
      99      211176 :   params.addParam<std::vector<unsigned int>>("target_variable_components",
     100      140784 :                                              std::vector<unsigned int>(),
     101             :                                              "The target array or vector variable component(s).");
     102             : 
     103             :   // Search options
     104      211176 :   params.addParam<bool>(
     105             :       "greedy_search",
     106      140784 :       false,
     107             :       "Whether or not to send a point to all the domains. If true, all the processors will be "
     108             :       "checked for a given point."
     109             :       "The code will be slow if this flag is on but it will give a better solution.");
     110      211176 :   params.addParam<bool>(
     111             :       "error_on_miss",
     112      140784 :       true,
     113             :       "Whether or not to error in the case that a target point is not found in the source domain.");
     114      211176 :   params.addParam<bool>("use_bounding_boxes",
     115      140784 :                         true,
     116             :                         "When set to false, bounding boxes will not be used to restrict the source "
     117             :                         "of the transfer. Either source applications must be set using the "
     118             :                         "from_mesh_division parameter, or a greedy search must be used.");
     119      211176 :   params.addParam<bool>(
     120             :       "use_nearest_app",
     121      140784 :       false,
     122             :       "When True, transfers from a child application will work by finding the nearest (using "
     123             :       "the `position` + mesh centroid) sub-app and query that app for the value to transfer.");
     124       70392 :   params.addParam<PositionsName>(
     125             :       "use_nearest_position",
     126             :       "Name of the the Positions object (in main app) such that transfers to/from a child "
     127             :       "application will work by finding the nearest position to a target and query only the "
     128             :       "app / points closer to this position than to any other position for the value to transfer.");
     129      211176 :   params.addParam<bool>(
     130             :       "from_app_must_contain_point",
     131      140784 :       false,
     132             :       "Wether on not the origin mesh must contain the point to evaluate data at. If false, this "
     133             :       "allows for interpolation between origin app meshes. Origin app bounding boxes are still "
     134             :       "considered so you may want to increase them with 'fixed_bounding_box_size'");
     135      211176 :   params.addParam<bool>("search_value_conflicts",
     136      140784 :                         false,
     137             :                         "Whether to look for potential conflicts between two valid and different "
     138             :                         "source values for any target point");
     139      211176 :   params.addParam<unsigned int>(
     140             :       "value_conflicts_output",
     141      140784 :       10,
     142             :       "Maximum number of conflicts to output if value-conflicts, from equidistant sources to a "
     143             :       "given transfer target location, search is turned on");
     144             : 
     145       70392 :   params.addParamNamesToGroup(
     146             :       "to_blocks from_blocks to_boundaries from_boundaries elemental_boundary_restriction "
     147             :       "from_mesh_division from_mesh_division_usage to_mesh_division to_mesh_division_usage",
     148             :       "Transfer spatial restriction");
     149       70392 :   params.addParamNamesToGroup(
     150             :       "greedy_search use_bounding_boxes use_nearest_app use_nearest_position "
     151             :       "search_value_conflicts",
     152             :       "Search algorithm");
     153       70392 :   params.addParamNamesToGroup("error_on_miss from_app_must_contain_point extrapolation_constant",
     154             :                               "Extrapolation behavior");
     155       70392 :   params.addParamNamesToGroup("bbox_factor fixed_bounding_box_size", "Source app bounding box");
     156      140784 :   return params;
     157       70392 : }
     158             : 
     159        6672 : MultiAppGeneralFieldTransfer::MultiAppGeneralFieldTransfer(const InputParameters & parameters)
     160             :   : MultiAppConservativeTransfer(parameters),
     161        6672 :     _from_var_components(getParam<std::vector<unsigned int>>("source_variable_components")),
     162        6672 :     _to_var_components(getParam<std::vector<unsigned int>>("target_variable_components")),
     163        6672 :     _use_bounding_boxes(getParam<bool>("use_bounding_boxes")),
     164        6672 :     _use_nearest_app(getParam<bool>("use_nearest_app")),
     165        6672 :     _nearest_positions_obj(
     166       13344 :         isParamValid("use_nearest_position")
     167        6672 :             ? &_fe_problem.getPositionsObject(getParam<PositionsName>("use_nearest_position"))
     168             :             : nullptr),
     169        6672 :     _source_app_must_contain_point(getParam<bool>("from_app_must_contain_point")),
     170        6672 :     _from_mesh_division_behavior(getParam<MooseEnum>("from_mesh_division_usage")),
     171        6672 :     _to_mesh_division_behavior(getParam<MooseEnum>("to_mesh_division_usage")),
     172        6672 :     _elemental_boundary_restriction_on_sides(
     173        6672 :         getParam<MooseEnum>("elemental_boundary_restriction") == "sides"),
     174        6672 :     _greedy_search(getParam<bool>("greedy_search")),
     175        6672 :     _search_value_conflicts(getParam<bool>("search_value_conflicts")),
     176        6672 :     _already_output_search_value_conflicts(false),
     177        6672 :     _search_value_conflicts_max_log(getParam<unsigned int>("value_conflicts_output")),
     178        6672 :     _error_on_miss(getParam<bool>("error_on_miss")),
     179        6672 :     _default_extrapolation_value(getParam<Real>("extrapolation_constant")),
     180        6672 :     _bbox_factor(getParam<Real>("bbox_factor")),
     181       13344 :     _fixed_bbox_size(isParamValid("fixed_bounding_box_size")
     182        6672 :                          ? getParam<std::vector<Real>>("fixed_bounding_box_size")
     183       33360 :                          : std::vector<Real>(3, 0))
     184             : {
     185        6672 :   _var_size = _to_var_names.size();
     186        6672 :   if (_to_var_names.size() != _from_var_names.size() && !parameters.isPrivate("source_variable"))
     187           4 :     paramError("variable", "The number of variables to transfer to and from should be equal");
     188             : 
     189             :   // Check the parameters of the components of the array / vector variable
     190        6668 :   if (_from_var_names.size() != _from_var_components.size() && _from_var_components.size() > 0)
     191           0 :     paramError("source_variable_components",
     192             :                "This parameter must be equal to the number of source variables");
     193        6668 :   if (_to_var_names.size() != _to_var_components.size() && _to_var_components.size() > 0)
     194           0 :     paramError("target_variable_components",
     195             :                "This parameter must be equal to the number of target variables");
     196             : 
     197             :   // Make simple 'use_nearest_app' parameter rely on Positions
     198        6668 :   if (_use_nearest_app)
     199             :   {
     200         148 :     if (_nearest_positions_obj)
     201           4 :       paramError("use_nearest_app", "Cannot use nearest-app and nearest-position together");
     202         144 :     if (!hasFromMultiApp())
     203           4 :       paramError("use_nearest_app",
     204             :                  "Should have a 'from_multiapp' when using the nearest-app informed search");
     205         140 :     auto pos_params = _app.getFactory().getValidParams("MultiAppPositions");
     206         420 :     pos_params.set<std::vector<MultiAppName>>("multiapps") = {getFromMultiApp()->name()};
     207         140 :     _fe_problem.addReporter("MultiAppPositions", "_created_for_" + name(), pos_params);
     208         140 :     _nearest_positions_obj = &_fe_problem.getPositionsObject("_created_for_" + name());
     209         140 :   }
     210             : 
     211             :   // Dont let users get wrecked by bounding boxes if it looks like they are trying to extrapolate
     212       13236 :   if (!_source_app_must_contain_point && _use_bounding_boxes &&
     213       13236 :       (_nearest_positions_obj || isParamSetByUser("from_app_must_contain_point")))
     214         518 :     if (!isParamSetByUser("bbox_factor") && !isParamSetByUser("fixed_bounding_box_size"))
     215           4 :       mooseWarning(
     216             :           "Extrapolation (nearest-source options, outside-app source) parameters have "
     217             :           "been passed, but no subapp bounding box expansion parameters have been passed.");
     218             : 
     219        6732 :   if (!_use_bounding_boxes &&
     220        6732 :       (isParamValid("fixed_bounding_box_size") || isParamSetByUser("bbox_factor")))
     221           0 :     paramError("use_bounding_boxes",
     222             :                "Cannot pass additional bounding box parameters (sizes, expansion, etc) if we are "
     223             :                "not using bounding boxes");
     224        6936 : }
     225             : 
     226             : void
     227        6484 : MultiAppGeneralFieldTransfer::initialSetup()
     228             : {
     229        6484 :   MultiAppConservativeTransfer::initialSetup();
     230             : 
     231             :   // Use IDs for block and boundary restriction
     232             :   // Loop over all source problems
     233       15211 :   for (const auto i_from : index_range(_from_problems))
     234             :   {
     235        8755 :     const auto & from_moose_mesh = _from_problems[i_from]->mesh(_displaced_source_mesh);
     236        8755 :     if (isParamValid("from_blocks"))
     237             :     {
     238        1326 :       auto & blocks = getParam<std::vector<SubdomainName>>("from_blocks");
     239        1326 :       std::vector<SubdomainID> ids = from_moose_mesh.getSubdomainIDs(blocks);
     240        1326 :       _from_blocks.insert(ids.begin(), ids.end());
     241        1326 :       if (_from_blocks.size() != blocks.size())
     242           0 :         paramError("from_blocks", "Some blocks were not found in the mesh");
     243        1326 :     }
     244             : 
     245        8755 :     if (isParamValid("from_boundaries"))
     246             :     {
     247         174 :       auto & boundary_names = getParam<std::vector<BoundaryName>>("from_boundaries");
     248         174 :       std::vector<BoundaryID> boundary_ids = from_moose_mesh.getBoundaryIDs(boundary_names);
     249         174 :       _from_boundaries.insert(boundary_ids.begin(), boundary_ids.end());
     250         174 :       if (_from_boundaries.size() != boundary_names.size())
     251           0 :         paramError("from_boundaries", "Some boundaries were not found in the mesh");
     252         174 :     }
     253             : 
     254        8755 :     if (isParamValid("from_mesh_division"))
     255             :     {
     256         328 :       const auto & mesh_div_name = getParam<MeshDivisionName>("from_mesh_division");
     257         328 :       _from_mesh_divisions.push_back(&_from_problems[i_from]->getMeshDivision(mesh_div_name));
     258             :       // Check that the behavior set makes sense
     259         328 :       if (_from_mesh_division_behavior == MeshDivisionTransferUse::RESTRICTION)
     260             :       {
     261         144 :         if (_from_mesh_divisions[i_from]->coversEntireMesh())
     262           0 :           mooseInfo("'from_mesh_division_usage' is set to use a spatial restriction but the "
     263           0 :                     "'from_mesh_division' for source app of global index " +
     264           0 :                     std::to_string(getGlobalSourceAppIndex(i_from)) +
     265             :                     " covers the entire mesh. Do not expect any restriction from a mesh "
     266             :                     "division that covers the entire mesh");
     267             :       }
     268         284 :       else if (_from_mesh_division_behavior == MeshDivisionTransferUse::MATCH_DIVISION_INDEX &&
     269         284 :                !isParamValid("to_mesh_division"))
     270           0 :         paramError("to_mesh_division_usage",
     271             :                    "Source mesh division cannot match target mesh division if no target mesh "
     272             :                    "division is specified");
     273         184 :       else if (_from_mesh_division_behavior == MeshDivisionTransferUse::MATCH_SUBAPP_INDEX)
     274             :       {
     275          84 :         if (!hasToMultiApp())
     276           4 :           paramError("from_mesh_division_usage",
     277             :                      "Cannot match source mesh division index to target subapp index if there is "
     278             :                      "only one target: the parent app (not a subapp)");
     279         160 :         else if (getToMultiApp()->numGlobalApps() !=
     280          80 :                  _from_mesh_divisions[i_from]->getNumDivisions())
     281           4 :           mooseWarning("Attempting to match target subapp index with the number of source mesh "
     282           4 :                        "divisions, which is " +
     283           8 :                        std::to_string(_from_mesh_divisions[i_from]->getNumDivisions()) +
     284           8 :                        " while there are " + std::to_string(getToMultiApp()->numGlobalApps()) +
     285             :                        " target subapps");
     286          76 :         if (_to_mesh_division_behavior == MeshDivisionTransferUse::MATCH_DIVISION_INDEX)
     287             :           // We do not support it because it would require sending the point + target app index +
     288             :           // target app division index, and we only send the Point + one number
     289           4 :           paramError("from_mesh_division_usage",
     290             :                      "We do not support using target subapp index for source division behavior and "
     291             :                      "matching the division index for the target mesh division behavior.");
     292             :       }
     293         100 :       else if (_from_mesh_division_behavior == "none")
     294           0 :         paramError("from_mesh_division_usage", "User must specify a 'from_mesh_division_usage'");
     295             :     }
     296        8427 :     else if (_from_mesh_division_behavior != "none")
     297           4 :       paramError("from_mesh_division",
     298             :                  "'from_mesh_division' must be specified if the usage method is specified");
     299             :   }
     300             : 
     301             :   // Loop over all target problems
     302       15174 :   for (const auto i_to : index_range(_to_problems))
     303             :   {
     304        8734 :     const auto & to_moose_mesh = _to_problems[i_to]->mesh(_displaced_target_mesh);
     305        8734 :     if (isParamValid("to_blocks"))
     306             :     {
     307        1140 :       auto & blocks = getParam<std::vector<SubdomainName>>("to_blocks");
     308        1140 :       std::vector<SubdomainID> ids = to_moose_mesh.getSubdomainIDs(blocks);
     309        1140 :       _to_blocks.insert(ids.begin(), ids.end());
     310        1140 :       if (_to_blocks.size() != blocks.size())
     311           0 :         paramError("to_blocks", "Some blocks were not found in the mesh");
     312        1140 :     }
     313             : 
     314        8734 :     if (isParamValid("to_boundaries"))
     315             :     {
     316        1288 :       auto & boundary_names = getParam<std::vector<BoundaryName>>("to_boundaries");
     317        1288 :       std::vector<BoundaryID> boundary_ids = to_moose_mesh.getBoundaryIDs(boundary_names);
     318        1288 :       _to_boundaries.insert(boundary_ids.begin(), boundary_ids.end());
     319        1288 :       if (_to_boundaries.size() != boundary_names.size())
     320           0 :         paramError("to_boundaries", "Some boundaries were not found in the mesh");
     321        1288 :     }
     322             : 
     323        8734 :     if (isParamValid("to_mesh_division"))
     324             :     {
     325         324 :       const auto & mesh_div_name = getParam<MeshDivisionName>("to_mesh_division");
     326         324 :       _to_mesh_divisions.push_back(&_to_problems[i_to]->getMeshDivision(mesh_div_name));
     327             :       // Check that the behavior set makes sense
     328         324 :       if (_to_mesh_division_behavior == MeshDivisionTransferUse::RESTRICTION)
     329             :       {
     330         144 :         if (_to_mesh_divisions[i_to]->coversEntireMesh())
     331           0 :           mooseInfo("'to_mesh_division_usage' is set to use a spatial restriction but the "
     332           0 :                     "'to_mesh_division' for target application of global index " +
     333           0 :                     std::to_string(getGlobalSourceAppIndex(i_to)) +
     334             :                     " covers the entire mesh. Do not expect any restriction from a mesh "
     335             :                     "division that covers the entire mesh");
     336             :       }
     337         180 :       else if (_to_mesh_division_behavior == MeshDivisionTransferUse::MATCH_DIVISION_INDEX)
     338             :       {
     339         100 :         if (!isParamValid("from_mesh_division"))
     340           0 :           paramError("to_mesh_division_usage",
     341             :                      "Target mesh division cannot match source mesh division if no source mesh "
     342             :                      "division is specified");
     343         200 :         else if ((*_from_mesh_divisions.begin())->getNumDivisions() !=
     344         100 :                  _to_mesh_divisions[i_to]->getNumDivisions())
     345           4 :           mooseWarning("Source and target mesh divisions do not have the same number of bins. If "
     346             :                        "this is what you expect, please reach out to a MOOSE or app developer to "
     347             :                        "ensure appropriate use");
     348             :       }
     349          80 :       else if (_to_mesh_division_behavior == MeshDivisionTransferUse::MATCH_SUBAPP_INDEX)
     350             :       {
     351          80 :         if (!hasFromMultiApp())
     352           4 :           paramError(
     353             :               "to_mesh_division_usage",
     354             :               "Cannot match target mesh division index to source subapp index if there is only one "
     355             :               "source: the parent app (not a subapp)");
     356          76 :         else if (getFromMultiApp()->numGlobalApps() != _to_mesh_divisions[i_to]->getNumDivisions())
     357           4 :           mooseWarning("Attempting to match source subapp index with the number of target mesh "
     358           4 :                        "divisions, which is " +
     359           8 :                        std::to_string(_to_mesh_divisions[i_to]->getNumDivisions()) +
     360           8 :                        " while there are " + std::to_string(getFromMultiApp()->numGlobalApps()) +
     361             :                        " source subapps");
     362          72 :         if (_from_mesh_division_behavior == MeshDivisionTransferUse::MATCH_DIVISION_INDEX)
     363           0 :           paramError(
     364             :               "from_mesh_division_usage",
     365             :               "We do not support using source subapp index for the target division behavior and "
     366             :               "matching the division index for the source mesh division behavior.");
     367             :       }
     368           0 :       else if (_to_mesh_division_behavior == "none")
     369           0 :         paramError("to_mesh_division_usage", "User must specify a 'to_mesh_division_usage'");
     370             :     }
     371        8410 :     else if (_to_mesh_division_behavior != "none")
     372           4 :       paramError("to_mesh_division",
     373           4 :                  "'to_mesh_division' must be specified if usage method '" +
     374           8 :                      Moose::stringify(_to_mesh_division_behavior) + "' is specified");
     375             :   }
     376             : 
     377             :   // Check if components are set correctly if using an array variable
     378       15139 :   for (const auto i_from : index_range(_from_problems))
     379             :   {
     380       16093 :     for (const auto var_index : make_range(_from_var_names.size()))
     381             :     {
     382             :       MooseVariableFieldBase & from_var =
     383       14788 :           _from_problems[i_from]->getVariable(0,
     384        7394 :                                               _from_var_names[var_index],
     385             :                                               Moose::VarKindType::VAR_ANY,
     386             :                                               Moose::VarFieldType::VAR_FIELD_ANY);
     387        7394 :       if (from_var.count() > 1 && _from_var_components.empty())
     388           4 :         paramError("source_variable_components", "Component must be passed for an array variable");
     389        7390 :       if (_from_var_components.size() && from_var.count() < _from_var_components[var_index])
     390           4 :         paramError("source_variable_components",
     391             :                    "Component passed is larger than size of variable");
     392             :     }
     393             :   }
     394       15134 :   for (const auto i_to : index_range(_to_problems))
     395             :   {
     396       17727 :     for (const auto var_index : make_range(_to_var_names.size()))
     397             :     {
     398             :       MooseVariableFieldBase & to_var =
     399       18050 :           _to_problems[i_to]->getVariable(0,
     400        9025 :                                           _to_var_names[var_index],
     401             :                                           Moose::VarKindType::VAR_ANY,
     402             :                                           Moose::VarFieldType::VAR_FIELD_ANY);
     403        9025 :       if (to_var.count() > 1 && _to_var_components.empty())
     404           4 :         paramError("target_variable_components", "Component must be passed for an array variable");
     405        9021 :       if (_to_var_components.size() && to_var.count() < _to_var_components[var_index])
     406           4 :         paramError("target_variable_components",
     407             :                    "Component passed is larger than size of variable");
     408             :     }
     409             :   }
     410             : 
     411             :   // Cache some quantities to avoid having to get them on every transferred point
     412        6424 :   if (_to_problems.size())
     413             :   {
     414        6409 :     _to_variables.resize(_to_var_names.size());
     415       12994 :     for (const auto i_var : index_range(_to_var_names))
     416       13170 :       _to_variables[i_var] = &_to_problems[0]->getVariable(
     417        6585 :           0, _to_var_names[i_var], Moose::VarKindType::VAR_ANY, Moose::VarFieldType::VAR_FIELD_ANY);
     418             :   }
     419        6424 : }
     420             : 
     421             : void
     422       64694 : MultiAppGeneralFieldTransfer::getAppInfo()
     423             : {
     424       64694 :   MultiAppFieldTransfer::getAppInfo();
     425             : 
     426             :   // Create the point locators to locate evaluation points in the origin mesh(es)
     427       64694 :   _from_point_locators.resize(_from_problems.size());
     428      134336 :   for (const auto i_from : index_range(_from_problems))
     429             :   {
     430       69642 :     const auto & from_moose_mesh = _from_problems[i_from]->mesh(_displaced_source_mesh);
     431       69642 :     _from_point_locators[i_from] =
     432      139284 :         PointLocatorBase::build(TREE_LOCAL_ELEMENTS, from_moose_mesh.getMesh());
     433       69642 :     _from_point_locators[i_from]->enable_out_of_mesh_mode();
     434             :   }
     435       64694 : }
     436             : 
     437             : void
     438       58210 : MultiAppGeneralFieldTransfer::execute()
     439             : {
     440       58210 :   TIME_SECTION(
     441             :       "MultiAppGeneralFieldTransfer::execute()_" + name(), 5, "Transfer execution " + name());
     442       58210 :   getAppInfo();
     443             : 
     444             :   // Set up bounding boxes, etc
     445       58210 :   prepareToTransfer();
     446             : 
     447             :   // loop over the vector of variables and make the transfer one by one
     448      116505 :   for (const auto i : make_range(_var_size))
     449       58342 :     transferVariable(i);
     450             : 
     451       58163 :   postExecute();
     452       58163 : }
     453             : 
     454             : void
     455       58210 : MultiAppGeneralFieldTransfer::prepareToTransfer()
     456             : {
     457             :   // Get the bounding boxes for the "from" domains.
     458             :   // Clean up _from_bboxes from the previous transfer execution
     459       58210 :   _from_bboxes.clear();
     460             : 
     461             :   // NOTE: This ignores the app's bounding box inflation and padding
     462       58210 :   _from_bboxes = getRestrictedFromBoundingBoxes();
     463             : 
     464             :   // Expand bounding boxes. Some desired points might be excluded
     465             :   // without an expansion
     466       58210 :   extendBoundingBoxes(_bbox_factor, _from_bboxes);
     467             : 
     468             :   // Figure out how many "from" domains each processor owns.
     469       58210 :   _froms_per_proc.clear();
     470       58210 :   _froms_per_proc = getFromsPerProc();
     471             : 
     472             :   // Get the index for the first source app every processor owns
     473       58210 :   _global_app_start_per_proc = getGlobalStartAppPerProc();
     474             : 
     475             :   // No need to keep searching for conflicts if the mesh has not changed
     476       58210 :   if (_already_output_search_value_conflicts && !_displaced_source_mesh && !_displaced_target_mesh)
     477         648 :     _search_value_conflicts = false;
     478       58210 : }
     479             : 
     480             : void
     481       58163 : MultiAppGeneralFieldTransfer::postExecute()
     482             : {
     483       58163 :   MultiAppConservativeTransfer::postExecute();
     484       58163 :   if (_search_value_conflicts)
     485         177 :     _already_output_search_value_conflicts = true;
     486       58163 : }
     487             : 
     488             : void
     489       58342 : MultiAppGeneralFieldTransfer::transferVariable(unsigned int i)
     490             : {
     491             :   mooseAssert(i < _var_size, "The variable of index " << i << " does not exist");
     492             : 
     493             :   // Find outgoing target points
     494             :   // We need to know what points we need to send which processors
     495             :   // One processor will receive many points from many processors
     496             :   // One point may go to different processors
     497       58342 :   ProcessorToPointVec outgoing_points;
     498       58342 :   extractOutgoingPoints(i, outgoing_points);
     499             : 
     500       58334 :   if (_from_var_names.size())
     501       56040 :     prepareEvaluationOfInterpValues(i);
     502             :   else
     503        2294 :     prepareEvaluationOfInterpValues(-1);
     504             : 
     505             :   // Fill values and app ids for incoming points
     506             :   // We are responsible to compute values for these incoming points
     507             :   auto gather_functor =
     508       89176 :       [this](processor_id_type /*pid*/,
     509             :              const std::vector<std::pair<Point, unsigned int>> & incoming_locations,
     510       89176 :              std::vector<std::pair<Real, Real>> & outgoing_vals)
     511             :   {
     512       89176 :     outgoing_vals.resize(
     513             :         incoming_locations.size(),
     514             :         {GeneralFieldTransfer::BetterOutOfMeshValue, GeneralFieldTransfer::BetterOutOfMeshValue});
     515             :     // Evaluate interpolation values for these incoming points
     516       89176 :     evaluateInterpValues(incoming_locations, outgoing_vals);
     517       89176 :   };
     518             : 
     519       58334 :   DofobjectToInterpValVec dofobject_to_valsvec(_to_problems.size());
     520       58334 :   InterpCaches interp_caches(_to_problems.size(), getMaxToProblemsBBoxDimensions());
     521       58334 :   InterpCaches distance_caches(_to_problems.size(), getMaxToProblemsBBoxDimensions());
     522             : 
     523             :   // Copy data out to incoming_vals_ids
     524       89176 :   auto action_functor = [this, &i, &dofobject_to_valsvec, &interp_caches, &distance_caches](
     525             :                             processor_id_type pid,
     526             :                             const std::vector<std::pair<Point, unsigned int>> & my_outgoing_points,
     527      178352 :                             const std::vector<std::pair<Real, Real>> & incoming_vals)
     528             :   {
     529       89176 :     auto & pointInfoVec = _processor_to_pointInfoVec[pid];
     530             : 
     531             :     // Cache interpolation values for each dof object / points
     532       89176 :     cacheIncomingInterpVals(pid,
     533             :                             i,
     534             :                             pointInfoVec,
     535             :                             my_outgoing_points,
     536             :                             incoming_vals,
     537             :                             dofobject_to_valsvec,
     538             :                             interp_caches,
     539             :                             distance_caches);
     540      147510 :   };
     541             : 
     542             :   // We assume incoming_vals_ids is ordered in the same way as outgoing_points
     543             :   // Hopefully, pull_parallel_vector_data will not mess up this
     544       58334 :   const std::pair<Real, Real> * ex = nullptr;
     545       58334 :   libMesh::Parallel::pull_parallel_vector_data(
     546       58334 :       comm(), outgoing_points, gather_functor, action_functor, ex);
     547             : 
     548             :   // Check for conflicts and overlaps from the maps that were built during the transfer
     549       58334 :   if (_search_value_conflicts)
     550         208 :     outputValueConflicts(i, dofobject_to_valsvec, distance_caches);
     551             : 
     552             :   // Set cached values into solution vector
     553       58303 :   setSolutionVectorValues(i, dofobject_to_valsvec, interp_caches);
     554       58295 : }
     555             : 
     556             : void
     557     5479250 : MultiAppGeneralFieldTransfer::locatePointReceivers(const Point point,
     558             :                                                    std::set<processor_id_type> & processors)
     559             : {
     560             :   // Check which processors have apps that may include or be near this point
     561             :   // A point may be close enough to several problems, hosted on several processes
     562     5479250 :   bool found = false;
     563             : 
     564             :   // Additional process-restriction techniques we could use (TODOs):
     565             :   // - create a heuristic for using nearest-positions
     566             :   // - from_mesh_divisions could be polled for which divisions they possess on each
     567             :   //   process, depending on the behavior chosen. This could limit potential senders.
     568             :   //   This should be done ahead of this function call, for all points at once
     569             : 
     570             :   // Determine the apps which will be receiving points (then sending values) using various
     571             :   // heuristics
     572     5479250 :   if (_use_nearest_app)
     573             :   {
     574             :     // Find the nearest position for the point
     575       17900 :     const bool initial = _fe_problem.getCurrentExecuteOnFlag() == EXEC_INITIAL;
     576             :     // The apps form the nearest positions here, this is the index of the nearest app
     577       17900 :     const auto nearest_index = _nearest_positions_obj->getNearestPositionIndex(point, initial);
     578             : 
     579             :     // Find the apps that are nearest to the same position
     580             :     // Global search over all applications
     581       42181 :     for (processor_id_type i_proc = 0; i_proc < n_processors(); ++i_proc)
     582             :     {
     583             :       // We need i_from to correspond to the global app index
     584       24281 :       unsigned int from0 = _global_app_start_per_proc[i_proc];
     585       63617 :       for (unsigned int i_from = from0; i_from < from0 + _froms_per_proc[i_proc]; ++i_from)
     586             :       {
     587       39336 :         if (_greedy_search || _search_value_conflicts || i_from == nearest_index)
     588             :         {
     589       20980 :           processors.insert(i_proc);
     590       20980 :           found = true;
     591             :         }
     592             :         mooseAssert(i_from < getFromMultiApp()->numGlobalApps(), "We should not reach this");
     593             :       }
     594             :     }
     595             :     mooseAssert((getFromMultiApp()->numGlobalApps() < n_processors() || processors.size() == 1) ||
     596             :                     _greedy_search || _search_value_conflicts,
     597             :                 "Should only be one source processor when using more processors than source apps");
     598             :   }
     599     5461350 :   else if (_use_bounding_boxes)
     600             :   {
     601             :     // We examine all (global) bounding boxes and find the minimum of the maximum distances within a
     602             :     // bounding box from the point. This creates a sphere around the point of interest. Any app
     603             :     // with a bounding box that intersects this sphere (with a bboxMinDistance <
     604             :     // nearest_max_distance) will be considered a potential source
     605             :     // NOTE: This is a heuristic. We could try others
     606             :     // NOTE: from_bboxes are in the reference space, as is the point.
     607     5459066 :     Real nearest_max_distance = std::numeric_limits<Real>::max();
     608    13535999 :     for (const auto & bbox : _from_bboxes)
     609             :     {
     610     8076933 :       Real distance = bboxMaxDistance(point, bbox);
     611     8076933 :       if (distance < nearest_max_distance)
     612     6689219 :         nearest_max_distance = distance;
     613             :     }
     614             : 
     615     5459066 :     unsigned int from0 = 0;
     616    13042372 :     for (processor_id_type i_proc = 0; i_proc < n_processors();
     617     7583306 :          from0 += _froms_per_proc[i_proc], ++i_proc)
     618             :       // i_from here is a hybrid index based on the cumulative sum of the apps per processor
     619    15660239 :       for (unsigned int i_from = from0; i_from < from0 + _froms_per_proc[i_proc]; ++i_from)
     620             :       {
     621     8076933 :         Real distance = bboxMinDistance(point, _from_bboxes[i_from]);
     622             :         // We will not break here because we want to send a point to all possible source domains
     623     8236372 :         if (_greedy_search || distance <= nearest_max_distance ||
     624      159439 :             _from_bboxes[i_from].contains_point(point))
     625             :         {
     626     7918212 :           processors.insert(i_proc);
     627     7918212 :           found = true;
     628             :         }
     629             :       }
     630             :   }
     631             :   // Greedy search will contact every single processor. It's not scalable, but if there's valid data
     632             :   // on any subapp on any process, it will find it
     633        2284 :   else if (_greedy_search)
     634             :   {
     635        2080 :     found = true;
     636        4940 :     for (const auto i_proc : make_range(n_processors()))
     637        2860 :       processors.insert(i_proc);
     638             :   }
     639             :   // Since we indicated that we only wanted values from a subapp with the same global index as the
     640             :   // target mesh division, we might as well only communicate with the process that owns this app
     641         404 :   else if (!_to_mesh_divisions.empty() &&
     642         200 :            _to_mesh_division_behavior == MeshDivisionTransferUse::MATCH_SUBAPP_INDEX)
     643             :   {
     644             :     // The target point could have a different index in each target mesh division. So on paper, we
     645             :     // would need to check all of them.
     646         200 :     auto saved_target_div = MooseMeshDivision::INVALID_DIVISION_INDEX;
     647         400 :     for (const auto i_to : index_range(_to_meshes))
     648             :     {
     649         400 :       const auto target_div = _to_mesh_divisions[i_to]->divisionIndex(
     650         200 :           _to_transforms[getGlobalTargetAppIndex(i_to)]->mapBack(point));
     651             :       // If it's the same division index, do not redo the search
     652         200 :       if (target_div == saved_target_div)
     653           0 :         continue;
     654             :       else
     655         200 :         saved_target_div = target_div;
     656             : 
     657             :       // Look for the processors owning a source-app with an index equal to the target mesh division
     658         475 :       for (const auto i_proc : make_range(n_processors()))
     659        1075 :         for (const auto i_from : make_range(_froms_per_proc[i_proc]))
     660         800 :           if (target_div == _global_app_start_per_proc[i_proc] + i_from)
     661             :           {
     662         200 :             processors.insert(i_proc);
     663         200 :             found = true;
     664             :           }
     665             :     }
     666             :   }
     667             :   else
     668           4 :     mooseError("No algorithm were selected to find which processes may send value data "
     669             :                "for a each target point. Please either specify using bounding boxes, "
     670             :                "greedy search, or to_mesh_division-based parameters");
     671             : 
     672             :   // Error out if we could not find this point when ask us to do so
     673     5479246 :   if (!found && _error_on_miss)
     674           0 :     mooseError(
     675             :         "Cannot find a source application to provide a value at point: ",
     676             :         point,
     677             :         " \n ",
     678             :         "It must be that mismatched meshes, between the source and target application, are being "
     679             :         "used.\nIf you are using the bounding boxes or nearest-app heuristics, or mesh-divisions, "
     680             :         "please consider using the greedy_search to confirm. Then consider choosing a different "
     681             :         "transfer type.\nThis check can be turned off by setting 'error_on_miss' to false. The "
     682             :         "'extrapolation_constant' parameter will be used to set the local value at missed points.");
     683     5479246 : }
     684             : 
     685             : void
     686     5479250 : MultiAppGeneralFieldTransfer::cacheOutgoingPointInfo(const Point point,
     687             :                                                      const dof_id_type dof_object_id,
     688             :                                                      const unsigned int problem_id,
     689             :                                                      ProcessorToPointVec & outgoing_points)
     690             : {
     691     5479250 :   std::set<processor_id_type> processors;
     692             :   // Find which processors will receive point data so they can send back value data
     693             :   // The list can be larger than needed, depending on the heuristic / algorithm used to make
     694             :   // the call on whether a processor (and the apps it runs) should be involved
     695     5479250 :   processors.clear();
     696     5479250 :   locatePointReceivers(point, processors);
     697             : 
     698             :   // We need to send this location data to these processors so they can send back values
     699    13043735 :   for (const auto pid : processors)
     700             :   {
     701             :     // Select which from_mesh_division the source data must come from for this point
     702     7564489 :     unsigned int required_source_division = 0;
     703     7564489 :     if (_from_mesh_division_behavior == MeshDivisionTransferUse::MATCH_SUBAPP_INDEX)
     704        8052 :       required_source_division = getGlobalTargetAppIndex(problem_id);
     705     7556437 :     else if (_from_mesh_division_behavior == MeshDivisionTransferUse::MATCH_DIVISION_INDEX ||
     706    15110014 :              _to_mesh_division_behavior == MeshDivisionTransferUse::MATCH_DIVISION_INDEX ||
     707     7553577 :              _to_mesh_division_behavior == MeshDivisionTransferUse::MATCH_SUBAPP_INDEX)
     708        7220 :       required_source_division = _to_mesh_divisions[problem_id]->divisionIndex(
     709        7220 :           _to_transforms[getGlobalTargetAppIndex(problem_id)]->mapBack(point));
     710             : 
     711             :     // Skip if we already know we don't want the point
     712     7564489 :     if (required_source_division == MooseMeshDivision::INVALID_DIVISION_INDEX)
     713           0 :       continue;
     714             : 
     715             :     // Store outgoing information for every source process
     716     7564489 :     outgoing_points[pid].push_back(std::pair<Point, unsigned int>(point, required_source_division));
     717             : 
     718             :     // Store point information locally for processing received data
     719             :     // We can use these information when inserting values into the solution vector
     720             :     PointInfo pointinfo;
     721     7564489 :     pointinfo.problem_id = problem_id;
     722     7564489 :     pointinfo.dof_object_id = dof_object_id;
     723     7564489 :     _processor_to_pointInfoVec[pid].push_back(pointinfo);
     724             :   }
     725     5479246 : }
     726             : 
     727             : void
     728       58342 : MultiAppGeneralFieldTransfer::extractOutgoingPoints(const unsigned int var_index,
     729             :                                                     ProcessorToPointVec & outgoing_points)
     730             : {
     731             :   // Get the variable name, with the accommodation for array/vector names
     732       58342 :   const auto & var_name = getToVarName(var_index);
     733             : 
     734             :   // Clean up the map from processor to pointInfo vector
     735             :   // This map should be consistent with outgoing_points
     736       58342 :   _processor_to_pointInfoVec.clear();
     737             : 
     738             :   // Loop over all problems
     739      119554 :   for (const auto i_to : index_range(_to_problems))
     740             :   {
     741       61220 :     const auto global_i_to = getGlobalTargetAppIndex(i_to);
     742             : 
     743             :     // libMesh EquationSystems
     744       61220 :     auto & es = getEquationSystem(*_to_problems[i_to], _displaced_target_mesh);
     745             :     // libMesh system that has this variable
     746       61220 :     System * to_sys = find_sys(es, var_name);
     747       61220 :     auto sys_num = to_sys->number();
     748       61220 :     auto var_num = _to_variables[var_index]->number();
     749       61220 :     auto & fe_type = _to_variables[var_index]->feType();
     750       61220 :     bool is_nodal = _to_variables[var_index]->isNodal();
     751             : 
     752             :     // Moose mesh
     753       61220 :     const auto & to_moose_mesh = _to_problems[i_to]->mesh(_displaced_target_mesh);
     754       61220 :     const auto & to_mesh = to_moose_mesh.getMesh();
     755             : 
     756             :     // We support more general variables via libMesh GenericProjector
     757       61220 :     if (fe_type.order > CONSTANT && !is_nodal)
     758             :     {
     759         335 :       GeneralFieldTransfer::RecordRequests<Number> f;
     760         335 :       GeneralFieldTransfer::RecordRequests<Gradient> g;
     761         335 :       GeneralFieldTransfer::NullAction<Number> nullsetter;
     762         335 :       const std::vector<unsigned int> varvec(1, var_num);
     763             : 
     764             :       libMesh::GenericProjector<GeneralFieldTransfer::RecordRequests<Number>,
     765             :                                 GeneralFieldTransfer::RecordRequests<Gradient>,
     766             :                                 Number,
     767             :                                 GeneralFieldTransfer::NullAction<Number>>
     768         335 :           request_gather(*to_sys, f, &g, nullsetter, varvec);
     769             : 
     770             :       // Defining only boundary values will not be enough to describe the variable, disallow it
     771         335 :       if (_to_boundaries.size() && (_to_variables[var_index]->getContinuity() == DISCONTINUOUS))
     772           4 :         mooseError("Higher order discontinuous elemental variables are not supported for "
     773             :                    "target-boundary "
     774             :                    "restricted transfers");
     775             : 
     776             :       // Not implemented as the target mesh division could similarly be cutting elements in an
     777             :       // arbitrary way with not enough requested points to describe the target variable
     778         331 :       if (!_to_mesh_divisions.empty() && !_to_mesh_divisions[i_to]->coversEntireMesh())
     779           0 :         mooseError("Higher order variable support not implemented for target mesh division "
     780             :                    "unless the mesh is fully covered / indexed in the mesh division. This must be "
     781             :                    "set programmatically in the MeshDivision object used.");
     782             : 
     783             :       // We dont look at boundary restriction, not supported for higher order target variables
     784             :       // Same for mesh divisions
     785         331 :       const auto & to_begin = _to_blocks.empty()
     786         557 :                                   ? to_mesh.active_local_elements_begin()
     787         331 :                                   : to_mesh.active_local_subdomain_set_elements_begin(_to_blocks);
     788             : 
     789         331 :       const auto & to_end = _to_blocks.empty()
     790         557 :                                 ? to_mesh.active_local_elements_end()
     791         331 :                                 : to_mesh.active_local_subdomain_set_elements_end(_to_blocks);
     792             : 
     793         331 :       ConstElemRange to_elem_range(to_begin, to_end);
     794             : 
     795         331 :       request_gather.project(to_elem_range);
     796             : 
     797         331 :       dof_id_type point_id = 0;
     798      265531 :       for (const Point & p : f.points_requested())
     799             :         // using the point number as a "dof_object_id" will serve to identify the point if we ever
     800             :         // rework interp/distance_cache into the dof_id_to_value maps
     801      265200 :         this->cacheOutgoingPointInfo(
     802      265200 :             (*_to_transforms[global_i_to])(p), point_id++, i_to, outgoing_points);
     803             : 
     804             :       // This is going to require more complicated transfer work
     805         331 :       if (!g.points_requested().empty())
     806           0 :         mooseError("We don't currently support variables with gradient degrees of freedom");
     807         331 :     }
     808       60885 :     else if (is_nodal)
     809             :     {
     810    10618662 :       for (const auto & node : to_mesh.local_node_ptr_range())
     811             :       {
     812             :         // Skip this node if the variable has no dofs at it.
     813     5281800 :         if (node->n_dofs(sys_num, var_num) < 1)
     814      267624 :           continue;
     815             : 
     816             :         // Skip if it is a block restricted transfer and current node does not have
     817             :         // specified blocks
     818     5014176 :         if (!_to_blocks.empty() && !inBlocks(_to_blocks, to_moose_mesh, node))
     819       13168 :           continue;
     820             : 
     821     5001008 :         if (!_to_boundaries.empty() && !onBoundaries(_to_boundaries, to_moose_mesh, node))
     822      169600 :           continue;
     823             : 
     824             :         // Skip if the node does not meet the target mesh division behavior
     825             :         // We cannot know from which app the data will come from so we cannot know
     826             :         // the source mesh division index and the source app global index
     827     4831408 :         if (!_to_mesh_divisions.empty() && _to_mesh_divisions[i_to]->divisionIndex(*node) ==
     828             :                                                MooseMeshDivision::INVALID_DIVISION_INDEX)
     829        4400 :           continue;
     830             : 
     831             :         // Cache point information
     832             :         // We will use this information later for setting values back to solution vectors
     833     4827008 :         cacheOutgoingPointInfo(
     834     4827008 :             (*_to_transforms[global_i_to])(*node), node->id(), i_to, outgoing_points);
     835       55066 :       }
     836             :     }
     837             :     else // Elemental, constant monomial
     838             :     {
     839        5815 :       for (const auto & elem :
     840      889922 :            as_range(to_mesh.local_elements_begin(), to_mesh.local_elements_end()))
     841             :       {
     842             :         // Skip this element if the variable has no dofs at it.
     843      439146 :         if (elem->n_dofs(sys_num, var_num) < 1)
     844        6976 :           continue;
     845             : 
     846             :         // Skip if the element is not inside the block restriction
     847      432170 :         if (!_to_blocks.empty() && !inBlocks(_to_blocks, elem))
     848        4544 :           continue;
     849             : 
     850             :         // Skip if the element does not have a side on the boundary
     851      427626 :         if (!_to_boundaries.empty() && !onBoundaries(_to_boundaries, to_moose_mesh, elem))
     852       37592 :           continue;
     853             : 
     854             :         // Skip if the element is not indexed within the mesh division
     855      390034 :         if (!_to_mesh_divisions.empty() && _to_mesh_divisions[i_to]->divisionIndex(*elem) ==
     856             :                                                MooseMeshDivision::INVALID_DIVISION_INDEX)
     857        2992 :           continue;
     858             : 
     859             :         // Cache point information
     860             :         // We will use this information later for setting values back to solution vectors
     861      387042 :         cacheOutgoingPointInfo((*_to_transforms[global_i_to])(elem->vertex_average()),
     862      387042 :                                elem->id(),
     863             :                                i_to,
     864             :                                outgoing_points);
     865        5815 :       } // for
     866             :     } // else
     867             :   } // for
     868       58334 : }
     869             : 
     870             : void
     871        6914 : MultiAppGeneralFieldTransfer::extractLocalFromBoundingBoxes(std::vector<BoundingBox> & local_bboxes)
     872             : {
     873        6914 :   local_bboxes.resize(_froms_per_proc[processor_id()]);
     874             :   // Find the index to the first of this processor's local bounding boxes.
     875        6914 :   unsigned int local_start = 0;
     876        8744 :   for (processor_id_type i_proc = 0; i_proc < n_processors() && i_proc != processor_id(); ++i_proc)
     877        1830 :     local_start += _froms_per_proc[i_proc];
     878             : 
     879             :   // Extract the local bounding boxes.
     880       15440 :   for (const auto i_from : make_range(_froms_per_proc[processor_id()]))
     881        8526 :     local_bboxes[i_from] = _from_bboxes[local_start + i_from];
     882        6914 : }
     883             : 
     884             : void
     885       89176 : MultiAppGeneralFieldTransfer::cacheIncomingInterpVals(
     886             :     processor_id_type pid,
     887             :     const unsigned int var_index,
     888             :     std::vector<PointInfo> & pointInfoVec,
     889             :     const std::vector<std::pair<Point, unsigned int>> & point_requests,
     890             :     const std::vector<std::pair<Real, Real>> & incoming_vals,
     891             :     DofobjectToInterpValVec & dofobject_to_valsvec,
     892             :     InterpCaches & interp_caches,
     893             :     InterpCaches & distance_caches)
     894             : {
     895             :   mooseAssert(pointInfoVec.size() == incoming_vals.size(),
     896             :               "Number of dof objects does not equal to the number of incoming values");
     897             : 
     898       89176 :   dof_id_type val_offset = 0;
     899     7653665 :   for (const auto & pointinfo : pointInfoVec)
     900             :   {
     901             :     // Retrieve target information from cached point infos
     902     7564489 :     const auto problem_id = pointinfo.problem_id;
     903     7564489 :     const auto dof_object_id = pointinfo.dof_object_id;
     904             : 
     905     7564489 :     auto & fe_type = _to_variables[var_index]->feType();
     906     7564489 :     bool is_nodal = _to_variables[var_index]->isNodal();
     907             : 
     908             :     // In the higher order elemental variable case, we receive point values, not nodal or
     909             :     // elemental. We use an InterpCache to store the values. The distance_cache is necessary to
     910             :     // choose between multiple origin problems sending values. This code could be unified with the
     911             :     // lower order order case by using the dofobject_to_valsvec
     912     7564489 :     if (fe_type.order > CONSTANT && !is_nodal)
     913             :     {
     914             :       // Cache solution on target mesh in its local frame of reference
     915      360675 :       InterpCache & value_cache = interp_caches[problem_id];
     916      360675 :       InterpCache & distance_cache = distance_caches[problem_id];
     917      360675 :       Point p = _to_transforms[getGlobalTargetAppIndex(problem_id)]->mapBack(
     918      360675 :           point_requests[val_offset].first);
     919      360675 :       const Number val = incoming_vals[val_offset].first;
     920             : 
     921             :       // Initialize distance to be able to compare
     922      360675 :       if (!distance_cache.hasKey(p))
     923      265200 :         distance_cache[p] = std::numeric_limits<Real>::max();
     924             : 
     925             :       // We should only have one closest value for each variable at any given point.
     926             :       // While there are shared Qps, on vertices for higher order variables usually,
     927             :       // the generic projector only queries each point once
     928      360675 :       if (_search_value_conflicts && !GeneralFieldTransfer::isBetterOutOfMeshValue(val) &&
     929      360675 :           value_cache.hasKey(p) != 0 && !MooseUtils::absoluteFuzzyEqual(value_cache[p], val) &&
     930           0 :           MooseUtils::absoluteFuzzyEqual(distance_cache[p], incoming_vals[val_offset].second))
     931           0 :         registerConflict(problem_id, dof_object_id, p, incoming_vals[val_offset].second, false);
     932             : 
     933             :       // if we use the nearest app, even if the value is bad we want to save the distance because
     934             :       // it's the distance to the app, if that's the closest app then so be it with the bad value
     935      428976 :       if ((!GeneralFieldTransfer::isBetterOutOfMeshValue(val) || _use_nearest_app) &&
     936       68301 :           MooseUtils::absoluteFuzzyGreaterThan(distance_cache[p], incoming_vals[val_offset].second))
     937             :       {
     938             :         // NOTE: We store the distance as well as the value. We really only need the
     939             :         // value to construct the variable, but the distance is used to make decisions in nearest
     940             :         // node schemes on which value to use
     941       64724 :         value_cache[p] = val;
     942       64724 :         distance_cache[p] = incoming_vals[val_offset].second;
     943             :       }
     944             :     }
     945             :     else
     946             :     {
     947             :       // Using the dof object pointer, so we can handle
     948             :       // both element and node using the same code
     949             : #ifndef NDEBUG
     950             :       auto var_num = _to_variables[var_index]->number();
     951             :       auto & to_sys = _to_variables[var_index]->sys();
     952             : 
     953             :       const MeshBase & to_mesh = _to_problems[problem_id]->mesh(_displaced_target_mesh).getMesh();
     954             :       const DofObject * dof_object_ptr = nullptr;
     955             :       const auto sys_num = to_sys.number();
     956             :       // It is a node
     957             :       if (is_nodal)
     958             :         dof_object_ptr = to_mesh.node_ptr(dof_object_id);
     959             :       // It is an element
     960             :       else
     961             :         dof_object_ptr = to_mesh.elem_ptr(dof_object_id);
     962             : 
     963             :       // We should only be supporting nodal and constant elemental
     964             :       // variables in this code path; if we see multiple DoFs on one
     965             :       // object we should have been using GenericProjector
     966             :       mooseAssert(dof_object_ptr->n_dofs(sys_num, var_num) == 1,
     967             :                   "Unexpectedly found " << dof_object_ptr->n_dofs(sys_num, var_num)
     968             :                                         << "dofs instead of 1");
     969             : #endif
     970             : 
     971     7203814 :       auto & dofobject_to_val = dofobject_to_valsvec[problem_id];
     972             : 
     973             :       // Check if we visited this dof object earlier
     974     7203814 :       auto values_ptr = dofobject_to_val.find(dof_object_id);
     975             :       // We did not visit this
     976     7203814 :       if (values_ptr == dofobject_to_val.end())
     977             :       {
     978             :         // Values for this dof object
     979     5214046 :         auto & val = dofobject_to_val[dof_object_id];
     980             :         // Interpolation value
     981     5214046 :         val.interp = incoming_vals[val_offset].first;
     982             :         // Where this value came from
     983     5214046 :         val.pid = pid;
     984             :         // Distance
     985     5214046 :         val.distance = incoming_vals[val_offset].second;
     986             :       }
     987             :       else
     988             :       {
     989     1989768 :         auto & val = values_ptr->second;
     990             : 
     991             :         // Look for value conflicts
     992     1989768 :         if (detectConflict(val.interp,
     993     1989768 :                            incoming_vals[val_offset].first,
     994             :                            val.distance,
     995     1989768 :                            incoming_vals[val_offset].second))
     996             :         {
     997             :           // Keep track of distance and value
     998             :           const Point p =
     999           0 :               getPointInTargetAppFrame(point_requests[val_offset].first,
    1000             :                                        problem_id,
    1001             :                                        "Registration of received equi-distant value conflict");
    1002           0 :           registerConflict(problem_id, dof_object_id, p, incoming_vals[val_offset].second, false);
    1003             :         }
    1004             : 
    1005             :         // We adopt values that are, in order of priority
    1006             :         // - valid (or from nearest app)
    1007             :         // - closest distance
    1008             :         // - the smallest rank with the same distance
    1009             :         // It is debatable whether we want invalid values from the nearest app. It could just be
    1010             :         // that the app position was closer but the extent of another child app was large enough
    1011     3979536 :         if ((!GeneralFieldTransfer::isBetterOutOfMeshValue(incoming_vals[val_offset].first) ||
    1012     3744829 :              _use_nearest_app) &&
    1013     1755061 :             (MooseUtils::absoluteFuzzyGreaterThan(val.distance, incoming_vals[val_offset].second) ||
    1014     1419540 :              ((val.pid > pid) &&
    1015      614163 :               MooseUtils::absoluteFuzzyEqual(val.distance, incoming_vals[val_offset].second))))
    1016             :         {
    1017      336049 :           val.interp = incoming_vals[val_offset].first;
    1018      336049 :           val.pid = pid;
    1019      336049 :           val.distance = incoming_vals[val_offset].second;
    1020             :         }
    1021             :       }
    1022             :     }
    1023             : 
    1024             :     // Move it to next position
    1025     7564489 :     val_offset++;
    1026             :   }
    1027       89176 : }
    1028             : 
    1029             : void
    1030         438 : MultiAppGeneralFieldTransfer::registerConflict(
    1031             :     unsigned int problem, dof_id_type dof_id, Point p, Real dist, bool local)
    1032             : {
    1033             :   // NOTE We could be registering the same conflict several times, we could count them instead
    1034         438 :   if (local)
    1035         438 :     _local_conflicts.push_back(std::make_tuple(problem, dof_id, p, dist));
    1036             :   else
    1037           0 :     _received_conflicts.push_back(std::make_tuple(problem, dof_id, p, dist));
    1038         438 : }
    1039             : 
    1040             : void
    1041         208 : MultiAppGeneralFieldTransfer::examineReceivedValueConflicts(
    1042             :     const unsigned int var_index,
    1043             :     const DofobjectToInterpValVec & dofobject_to_valsvec,
    1044             :     const InterpCaches & distance_caches)
    1045             : {
    1046         208 :   const auto var_name = getToVarName(var_index);
    1047             :   // We must check a posteriori because we could have two
    1048             :   // equidistant points with different values from two different problems, but a third point from
    1049             :   // another problem is actually closer, so there is no conflict because only that last one
    1050             :   // matters We check here whether the potential conflicts actually were the nearest points Loop
    1051             :   // over potential conflicts
    1052         208 :   for (auto conflict_it = _received_conflicts.begin(); conflict_it != _received_conflicts.end();)
    1053             :   {
    1054           0 :     const auto potential_conflict = *conflict_it;
    1055           0 :     bool overlap_found = false;
    1056             : 
    1057             :     // Extract info for the potential conflict
    1058           0 :     const unsigned int problem_id = std::get<0>(potential_conflict);
    1059           0 :     const dof_id_type dof_object_id = std::get<1>(potential_conflict);
    1060           0 :     const Point p = std::get<2>(potential_conflict);
    1061           0 :     const Real distance = std::get<3>(potential_conflict);
    1062             : 
    1063             :     // Extract target variable info
    1064           0 :     auto & es = getEquationSystem(*_to_problems[problem_id], _displaced_target_mesh);
    1065           0 :     System * to_sys = find_sys(es, var_name);
    1066           0 :     auto var_num = to_sys->variable_number(var_name);
    1067           0 :     auto & fe_type = to_sys->variable_type(var_num);
    1068           0 :     bool is_nodal = _to_variables[var_index]->isNodal();
    1069             : 
    1070             :     // Higher order elemental
    1071           0 :     if (fe_type.order > CONSTANT && !is_nodal)
    1072             :     {
    1073           0 :       auto cached_distance = distance_caches[problem_id].find(p);
    1074           0 :       if (cached_distance == distance_caches[problem_id].end())
    1075           0 :         mooseError("Conflict point was not found in the map of all origin-target distances");
    1076             :       // Distance is still the distance when we detected a potential overlap
    1077           0 :       if (MooseUtils::absoluteFuzzyEqual(cached_distance->second, distance))
    1078           0 :         overlap_found = true;
    1079             :     }
    1080             :     // Nodal and const monomial variable
    1081           0 :     else if (MooseUtils::absoluteFuzzyEqual(
    1082           0 :                  dofobject_to_valsvec[problem_id].find(dof_object_id)->second.distance, distance))
    1083           0 :       overlap_found = true;
    1084             : 
    1085             :     // Map will only keep the actual overlaps
    1086           0 :     if (!overlap_found)
    1087           0 :       _received_conflicts.erase(conflict_it);
    1088             :     else
    1089           0 :       ++conflict_it;
    1090             :   }
    1091         208 : }
    1092             : 
    1093             : void
    1094         208 : MultiAppGeneralFieldTransfer::examineLocalValueConflicts(
    1095             :     const unsigned int var_index,
    1096             :     const DofobjectToInterpValVec & dofobject_to_valsvec,
    1097             :     const InterpCaches & distance_caches)
    1098             : {
    1099         208 :   const auto var_name = getToVarName(var_index);
    1100             :   // We must check a posteriori because we could have:
    1101             :   // - two equidistant points with different values from two different problems
    1102             :   // - two (or more) equidistant points with different values from the same problem
    1103             :   // but a third point/value couple from another problem is actually closer, so there is no
    1104             :   // conflict because only that last one matters. We check here whether the potential conflicts
    1105             :   // actually were the nearest points. We use several global reductions. If there are not too many
    1106             :   // potential conflicts (and there should not be in a well-posed problem) it should be manageably
    1107             :   // expensive
    1108             : 
    1109             :   // Move relevant conflict info (location, distance) to a smaller data structure
    1110         208 :   std::vector<std::tuple<Point, Real>> potential_conflicts;
    1111         208 :   potential_conflicts.reserve(_local_conflicts.size());
    1112             : 
    1113             :   // Loop over potential conflicts to broadcast all the conflicts
    1114         646 :   for (auto conflict_it = _local_conflicts.begin(); conflict_it != _local_conflicts.end();
    1115         438 :        ++conflict_it)
    1116             :   {
    1117             :     // Extract info for the potential conflict
    1118         438 :     const auto potential_conflict = *conflict_it;
    1119         438 :     const unsigned int i_from = std::get<0>(potential_conflict);
    1120         438 :     Point p = std::get<2>(potential_conflict);
    1121         438 :     const Real distance = std::get<3>(potential_conflict);
    1122             :     // If not using nearest-positions: potential conflict was saved in the source frame
    1123             :     // If using nearest-positions: potential conflict was saved in the reference frame
    1124         438 :     if (!_nearest_positions_obj)
    1125             :     {
    1126         390 :       const auto from_global_num = getGlobalSourceAppIndex(i_from);
    1127         390 :       p = (*_from_transforms[from_global_num])(p);
    1128             :     }
    1129             : 
    1130             :     // Send data in the global frame of reference
    1131         438 :     potential_conflicts.push_back(std::make_tuple(p, distance));
    1132             :   }
    1133         208 :   _communicator.allgather(potential_conflicts, false);
    1134             :   // conflicts could have been reported multiple times within a tolerance
    1135         208 :   std::sort(potential_conflicts.begin(), potential_conflicts.end());
    1136         208 :   potential_conflicts.erase(unique(potential_conflicts.begin(),
    1137             :                                    potential_conflicts.end(),
    1138         407 :                                    [](auto l, auto r)
    1139             :                                    {
    1140         451 :                                      return std::get<0>(l).absolute_fuzzy_equals(std::get<0>(r)) &&
    1141         451 :                                             std::abs(std::get<1>(l) - std::get<1>(r)) < TOLERANCE;
    1142             :                                    }),
    1143         208 :                             potential_conflicts.end());
    1144             : 
    1145         208 :   std::vector<std::tuple<Point, Real>> real_conflicts;
    1146         208 :   real_conflicts.reserve(potential_conflicts.size());
    1147             : 
    1148             :   // For each potential conflict, we need to identify what problem asked for that value
    1149         602 :   for (auto conflict_it = potential_conflicts.begin(); conflict_it != potential_conflicts.end();
    1150         394 :        ++conflict_it)
    1151             :   {
    1152             :     // Extract info for the potential conflict
    1153         394 :     auto potential_conflict = *conflict_it;
    1154         394 :     const Point p = std::get<0>(potential_conflict);
    1155         394 :     const Real distance = std::get<1>(potential_conflict);
    1156             : 
    1157             :     // Check all the problems to try to find this requested point in the data structures filled
    1158             :     // with the received information
    1159         394 :     bool target_found = false;
    1160         394 :     bool conflict_real = false;
    1161         788 :     for (const auto i_to : index_range(_to_problems))
    1162             :     {
    1163             :       // Extract variable info
    1164         394 :       auto & es = getEquationSystem(*_to_problems[i_to], _displaced_target_mesh);
    1165         394 :       System * to_sys = find_sys(es, var_name);
    1166         394 :       auto var_num = to_sys->variable_number(var_name);
    1167         394 :       auto & fe_type = to_sys->variable_type(var_num);
    1168         394 :       bool is_nodal = _to_variables[var_index]->isNodal();
    1169             : 
    1170             :       // Move to the local frame of reference for the target problem
    1171             :       Point local_p =
    1172         394 :           getPointInTargetAppFrame(p, i_to, "Resolution of local value conflicts detected");
    1173             : 
    1174             :       // Higher order elemental
    1175         394 :       if (fe_type.order > CONSTANT && !is_nodal)
    1176             :       {
    1177             :         // distance_caches finds use a binned floating point search
    1178           0 :         auto cached_distance = distance_caches[i_to].find(local_p);
    1179           0 :         if (cached_distance != distance_caches[i_to].end())
    1180             :         {
    1181           0 :           target_found = true;
    1182             :           // Distance between source & target is still the distance we found in the sending
    1183             :           // process when we detected a potential overlap while gathering values to send
    1184           0 :           if (MooseUtils::absoluteFuzzyEqual(cached_distance->second, distance))
    1185           0 :             conflict_real = true;
    1186             :         }
    1187             :       }
    1188             :       // Nodal-value-dof-only and const monomial variable
    1189             :       else
    1190             :       {
    1191             :         // Find the dof id for the variable to be set
    1192         394 :         dof_id_type dof_object_id = std::numeric_limits<dof_id_type>::max();
    1193         394 :         auto pl = _to_problems[i_to]->mesh().getPointLocator();
    1194         394 :         pl->enable_out_of_mesh_mode();
    1195         394 :         if (is_nodal)
    1196             :         {
    1197         366 :           auto node = pl->locate_node(local_p);
    1198         366 :           if (node)
    1199             :             // this is not the dof_id for the variable, but the dof_object_id
    1200         366 :             dof_object_id = node->id();
    1201             :         }
    1202             :         else
    1203             :         {
    1204          28 :           auto elem = (*pl)(local_p);
    1205          28 :           if (elem)
    1206          28 :             dof_object_id = elem->id();
    1207             :         }
    1208         394 :         pl->disable_out_of_mesh_mode();
    1209             : 
    1210             :         // point isn't even in mesh
    1211         394 :         if (dof_object_id == std::numeric_limits<dof_id_type>::max())
    1212           0 :           continue;
    1213             : 
    1214             :         // this dof was not requested by this problem on this process
    1215         394 :         if (dofobject_to_valsvec[i_to].find(dof_object_id) == dofobject_to_valsvec[i_to].end())
    1216           0 :           continue;
    1217             : 
    1218         394 :         target_found = true;
    1219             :         // Check the saved distance in the vector of saved results. If the same, then the local
    1220             :         // conflict we detected with that distance is still an issue after receiving all values
    1221         788 :         if (MooseUtils::absoluteFuzzyEqual(
    1222         394 :                 dofobject_to_valsvec[i_to].find(dof_object_id)->second.distance, distance))
    1223         394 :           conflict_real = true;
    1224         394 :       }
    1225             :     }
    1226             :     // Only keep the actual conflicts / overlaps
    1227         394 :     if (target_found && conflict_real)
    1228         394 :       real_conflicts.push_back(potential_conflict);
    1229             :   }
    1230             : 
    1231             :   // Communicate real conflicts to all so they can be checked by every process
    1232         208 :   _communicator.allgather(real_conflicts, false);
    1233             : 
    1234             :   // Delete potential conflicts that were resolved
    1235             :   // Each local list of conflicts will now be updated. It's important to keep conflict lists local
    1236             :   // so we can give more context like the sending processor id (the domain of which can be
    1237             :   // inspected by the user)
    1238         646 :   for (auto conflict_it = _local_conflicts.begin(); conflict_it != _local_conflicts.end();)
    1239             :   {
    1240             :     // Extract info for the potential conflict
    1241         438 :     const auto potential_conflict = *conflict_it;
    1242         438 :     const unsigned int i_from = std::get<0>(potential_conflict);
    1243         438 :     Point p = std::get<2>(potential_conflict);
    1244         438 :     const Real distance = std::get<3>(potential_conflict);
    1245         438 :     if (!_nearest_positions_obj)
    1246             :     {
    1247         390 :       const auto from_global_num = getGlobalSourceAppIndex(i_from);
    1248         390 :       p = (*_from_transforms[from_global_num])(p);
    1249             :     }
    1250             : 
    1251             :     // If not in the vector of real conflicts, was not real so delete it
    1252         438 :     if (std::find_if(real_conflicts.begin(),
    1253             :                      real_conflicts.end(),
    1254        7037 :                      [p, distance](const auto & item)
    1255             :                      {
    1256        7037 :                        return std::get<0>(item).absolute_fuzzy_equals(p) &&
    1257        7037 :                               std::abs(std::get<1>(item) - distance) < TOLERANCE;
    1258         876 :                      }) == real_conflicts.end())
    1259           0 :       _local_conflicts.erase(conflict_it);
    1260             :     else
    1261         438 :       ++conflict_it;
    1262             :   }
    1263         208 : }
    1264             : 
    1265             : void
    1266         208 : MultiAppGeneralFieldTransfer::outputValueConflicts(
    1267             :     const unsigned int var_index,
    1268             :     const DofobjectToInterpValVec & dofobject_to_valsvec,
    1269             :     const InterpCaches & distance_caches)
    1270             : {
    1271             :   // Remove potential conflicts that did not materialize, the value did not end up being used
    1272         208 :   examineReceivedValueConflicts(var_index, dofobject_to_valsvec, distance_caches);
    1273         208 :   examineLocalValueConflicts(var_index, dofobject_to_valsvec, distance_caches);
    1274             : 
    1275             :   // Output the conflicts from the selection of local values (evaluateInterpValues-type routines)
    1276             :   // to send in response to value requests at target points
    1277         208 :   const std::string rank_str = std::to_string(_communicator.rank());
    1278         208 :   if (_local_conflicts.size())
    1279             :   {
    1280          31 :     unsigned int num_outputs = 0;
    1281          31 :     std::string local_conflicts_string = "";
    1282             :     std::string potential_reasons =
    1283             :         "Are some points in target mesh equidistant from the sources "
    1284          31 :         "(nodes/centroids/apps/positions, depending on transfer) in origin mesh(es)?\n";
    1285          31 :     if (hasFromMultiApp() && _from_problems.size() > 1)
    1286          23 :       potential_reasons += "Are multiple subapps overlapping?\n";
    1287         469 :     for (const auto & conflict : _local_conflicts)
    1288             :     {
    1289         438 :       const unsigned int problem_id = std::get<0>(conflict);
    1290         438 :       Point p = std::get<2>(conflict);
    1291         438 :       num_outputs++;
    1292             : 
    1293         438 :       std::string origin_domain_message;
    1294         438 :       if (hasFromMultiApp() && !_nearest_positions_obj)
    1295             :       {
    1296             :         // NOTES:
    1297             :         // - The origin app for a conflict may not be unique.
    1298             :         // - The conflicts vectors only store the conflictual points, not the original one
    1299             :         //   The original value found with a given distance could be retrieved from the main
    1300             :         //   caches
    1301         390 :         const auto app_id = _from_local2global_map[problem_id];
    1302         390 :         origin_domain_message = "In source child app " + std::to_string(app_id) + " mesh,";
    1303             :       }
    1304             :       // We can't locate the source app when considering nearest positions, so we saved the data
    1305             :       // in the reference space. So we return the conflict location in the target app (parent or
    1306             :       // sibling) instead
    1307          48 :       else if (hasFromMultiApp() && _nearest_positions_obj)
    1308             :       {
    1309          48 :         if (_to_problems.size() == 1 || _skip_coordinate_collapsing)
    1310             :         {
    1311          48 :           p = (*_to_transforms[0])(p);
    1312          48 :           origin_domain_message = "In target app mesh,";
    1313             :         }
    1314             :         else
    1315           0 :           origin_domain_message = "In reference (post-coordinate collapse) mesh,";
    1316             :       }
    1317             :       else
    1318           0 :         origin_domain_message = "In source parent app mesh,";
    1319             : 
    1320         438 :       if (num_outputs < _search_value_conflicts_max_log)
    1321         348 :         local_conflicts_string += origin_domain_message + " point: (" + std::to_string(p(0)) +
    1322         696 :                                   ", " + std::to_string(p(1)) + ", " + std::to_string(p(2)) +
    1323         696 :                                   "), equi-distance: " + std::to_string(std::get<3>(conflict)) +
    1324         174 :                                   "\n";
    1325         264 :       else if (num_outputs == _search_value_conflicts_max_log)
    1326             :         local_conflicts_string +=
    1327             :             "Maximum output of the search for value conflicts has been reached. Further conflicts "
    1328           8 :             "will not be output.\nIncrease 'search_value_conflicts_max_log' to output more.";
    1329         438 :     }
    1330             :     // Explicitly name source to give more context
    1331          31 :     std::string source_str = "unknown";
    1332          31 :     if (_from_var_names.size())
    1333          27 :       source_str = "variable '" + getFromVarName(var_index) + "'";
    1334           4 :     else if (isParamValid("source_user_object"))
    1335           4 :       source_str = "user object '" + getParam<UserObjectName>("source_user_object") + "'";
    1336             : 
    1337          62 :     mooseWarning("On rank " + rank_str +
    1338             :                  ", multiple valid values from equidistant points were "
    1339          31 :                  "found in the origin mesh for source " +
    1340          62 :                  source_str + " for " + std::to_string(_local_conflicts.size()) +
    1341          31 :                  " target points.\n" + potential_reasons + "Conflicts detected at :\n" +
    1342             :                  local_conflicts_string);
    1343           0 :   }
    1344             : 
    1345             :   // Output the conflicts discovered when receiving values from multiple origin problems
    1346         177 :   if (_received_conflicts.size())
    1347             :   {
    1348           0 :     unsigned int num_outputs = 0;
    1349           0 :     std::string received_conflicts_string = "";
    1350             :     std::string potential_reasons =
    1351             :         "Are some points in target mesh equidistant from the sources "
    1352           0 :         "(nodes/centroids/apps/positions, depending on transfer) in origin mesh(es)?\n";
    1353           0 :     if (hasToMultiApp() && _to_problems.size() > 1)
    1354           0 :       potential_reasons += "Are multiple subapps overlapping?\n";
    1355           0 :     for (const auto & conflict : _received_conflicts)
    1356             :     {
    1357             :       // Extract info for the potential overlap
    1358           0 :       const unsigned int problem_id = std::get<0>(conflict);
    1359           0 :       const Point p = std::get<2>(conflict);
    1360           0 :       num_outputs++;
    1361             : 
    1362           0 :       std::string target_domain_message;
    1363           0 :       if (hasToMultiApp())
    1364             :       {
    1365           0 :         const auto app_id = _to_local2global_map[problem_id];
    1366           0 :         target_domain_message = "In target child app " + std::to_string(app_id) + " mesh,";
    1367             :       }
    1368             :       else
    1369           0 :         target_domain_message = "In target parent app mesh,";
    1370             : 
    1371           0 :       if (num_outputs < _search_value_conflicts_max_log)
    1372           0 :         received_conflicts_string += target_domain_message + " point: (" + std::to_string(p(0)) +
    1373           0 :                                      ", " + std::to_string(p(1)) + ", " + std::to_string(p(2)) +
    1374           0 :                                      "), equi-distance: " + std::to_string(std::get<3>(conflict)) +
    1375           0 :                                      "\n";
    1376           0 :       else if (num_outputs == _search_value_conflicts_max_log)
    1377             :         received_conflicts_string +=
    1378             :             "Maximum output of the search for value conflicts has been reached. Further conflicts "
    1379           0 :             "will not be output.\nIncrease 'search_value_conflicts_max_log' to output more.";
    1380           0 :     }
    1381           0 :     mooseWarning("On rank " + rank_str +
    1382             :                  ", multiple valid values from equidistant points were "
    1383           0 :                  "received for target variable '" +
    1384           0 :                  getToVarName(var_index) + "' for " + std::to_string(_received_conflicts.size()) +
    1385           0 :                  " target points.\n" + potential_reasons + "Conflicts detected at :\n" +
    1386             :                  received_conflicts_string);
    1387           0 :   }
    1388             : 
    1389         177 :   if (_local_conflicts.empty() && _received_conflicts.empty())
    1390             :   {
    1391         177 :     if (isParamSetByUser("search_value_conflict"))
    1392           0 :       mooseInfo("Automated diagnosis did not detect floating point indetermination in transfer");
    1393         177 :     else if (_to_problems.size() > 10 || _from_problems.size() > 10 || _communicator.size() > 10)
    1394           0 :       mooseInfo(
    1395             :           "Automated diagnosis did not detect any floating point indetermination in "
    1396             :           "the transfer. You may consider turning it off using `search_value_conflicts=false` "
    1397             :           "to improve performance/scalability.");
    1398             :   }
    1399             : 
    1400             :   // Reset the conflicts vectors, to be used for checking conflicts when transferring the next
    1401             :   // variable
    1402         177 :   _local_conflicts.clear();
    1403         177 :   _received_conflicts.clear();
    1404         177 : }
    1405             : 
    1406             : void
    1407       58303 : MultiAppGeneralFieldTransfer::setSolutionVectorValues(
    1408             :     const unsigned int var_index,
    1409             :     const DofobjectToInterpValVec & dofobject_to_valsvec,
    1410             :     const InterpCaches & interp_caches)
    1411             : {
    1412             :   // Get the variable name, with the accommodation for array/vector names
    1413       58303 :   const auto & var_name = getToVarName(var_index);
    1414             : 
    1415      119476 :   for (const auto problem_id : index_range(_to_problems))
    1416             :   {
    1417       61181 :     auto & dofobject_to_val = dofobject_to_valsvec[problem_id];
    1418             : 
    1419             :     // libMesh EquationSystems
    1420             :     // NOTE: we would expect to set variables from the displaced equation system here
    1421       61181 :     auto & es = getEquationSystem(*_to_problems[problem_id], false);
    1422             : 
    1423             :     // libMesh system
    1424       61181 :     System * to_sys = find_sys(es, var_name);
    1425             : 
    1426             :     // libMesh mesh
    1427       61181 :     const MeshBase & to_mesh = _to_problems[problem_id]->mesh(_displaced_target_mesh).getMesh();
    1428       61181 :     auto var_num = to_sys->variable_number(var_name);
    1429       61181 :     auto sys_num = to_sys->number();
    1430             : 
    1431       61181 :     auto & fe_type = _to_variables[var_index]->feType();
    1432       61181 :     bool is_nodal = _to_variables[var_index]->isNodal();
    1433             : 
    1434       61181 :     if (fe_type.order > CONSTANT && !is_nodal)
    1435             :     {
    1436             :       // We may need to use existing data values in places where the
    1437             :       // from app domain doesn't overlap
    1438         331 :       MeshFunction to_func(es, *to_sys->current_local_solution, to_sys->get_dof_map(), var_num);
    1439         331 :       to_func.init();
    1440             : 
    1441             :       GeneralFieldTransfer::CachedData<Number> f(
    1442         331 :           interp_caches[problem_id], to_func, _default_extrapolation_value);
    1443         331 :       libMesh::VectorSetAction<Number> setter(*to_sys->solution);
    1444         331 :       const std::vector<unsigned int> varvec(1, var_num);
    1445             : 
    1446             :       libMesh::GenericProjector<GeneralFieldTransfer::CachedData<Number>,
    1447             :                                 GeneralFieldTransfer::CachedData<Gradient>,
    1448             :                                 Number,
    1449             :                                 libMesh::VectorSetAction<Number>>
    1450         331 :           set_solution(*to_sys, f, nullptr, setter, varvec);
    1451             : 
    1452             :       // We dont look at boundary restriction, not supported for higher order target variables
    1453         331 :       const auto & to_begin = _to_blocks.empty()
    1454         557 :                                   ? to_mesh.active_local_elements_begin()
    1455         331 :                                   : to_mesh.active_local_subdomain_set_elements_begin(_to_blocks);
    1456             : 
    1457         331 :       const auto & to_end = _to_blocks.empty()
    1458         557 :                                 ? to_mesh.active_local_elements_end()
    1459         331 :                                 : to_mesh.active_local_subdomain_set_elements_end(_to_blocks);
    1460             : 
    1461         331 :       ConstElemRange active_local_elem_range(to_begin, to_end);
    1462             : 
    1463         331 :       set_solution.project(active_local_elem_range);
    1464         331 :     }
    1465             :     else
    1466             :     {
    1467     5270105 :       for (const auto & val_pair : dofobject_to_val)
    1468             :       {
    1469     5209263 :         const auto dof_object_id = val_pair.first;
    1470             : 
    1471     5209263 :         const DofObject * dof_object = nullptr;
    1472     5209263 :         if (is_nodal)
    1473     4823257 :           dof_object = to_mesh.node_ptr(dof_object_id);
    1474             :         else
    1475      386006 :           dof_object = to_mesh.elem_ptr(dof_object_id);
    1476             : 
    1477     5209263 :         const auto dof = dof_object->dof_number(sys_num, var_num, 0);
    1478     5209263 :         const auto val = val_pair.second.interp;
    1479             : 
    1480             :         // This will happen if meshes are mismatched
    1481     5209263 :         if (_error_on_miss && GeneralFieldTransfer::isBetterOutOfMeshValue(val))
    1482             :         {
    1483             :           const auto target_location =
    1484           8 :               hasToMultiApp()
    1485           8 :                   ? " on target app " + std::to_string(getGlobalTargetAppIndex(problem_id))
    1486          12 :                   : " on parent app";
    1487           8 :           const auto info_msg = "\nThis check can be turned off by setting 'error_on_miss' to "
    1488             :                                 "false. The 'extrapolation_constant' parameter will be used to set "
    1489             :                                 "the local value at missed points.";
    1490           8 :           if (is_nodal)
    1491           4 :             mooseError("No source value for node ",
    1492             :                        dof_object_id,
    1493             :                        target_location,
    1494             :                        " could be located. Node details:\n",
    1495           4 :                        _to_meshes[problem_id]->nodePtr(dof_object_id)->get_info(),
    1496             :                        "\n",
    1497             :                        info_msg);
    1498             :           else
    1499           4 :             mooseError("No source value for element ",
    1500             :                        dof_object_id,
    1501             :                        target_location,
    1502             :                        " could be located. Element details:\n",
    1503           4 :                        _to_meshes[problem_id]->elemPtr(dof_object_id)->get_info(),
    1504             :                        "\n",
    1505             :                        info_msg);
    1506           0 :         }
    1507             : 
    1508             :         // We should not put garbage into our solution vector
    1509             :         // but it can be that we want to set it to a different value than what was already there
    1510             :         // for example: the source app has been displaced and was sending an indicator of its
    1511             :         // position
    1512     5209255 :         if (GeneralFieldTransfer::isBetterOutOfMeshValue(val))
    1513             :         {
    1514      523368 :           if (!GeneralFieldTransfer::isBetterOutOfMeshValue(_default_extrapolation_value))
    1515      522920 :             to_sys->solution->set(dof, _default_extrapolation_value);
    1516      523368 :           continue;
    1517             :         }
    1518             : 
    1519     4685887 :         to_sys->solution->set(dof, val);
    1520             :       }
    1521             :     }
    1522             : 
    1523       61173 :     to_sys->solution->close();
    1524             :     // Sync local solutions
    1525       61173 :     to_sys->update();
    1526             :   }
    1527       58295 : }
    1528             : 
    1529             : bool
    1530     1836484 : MultiAppGeneralFieldTransfer::acceptPointInOriginMesh(unsigned int i_from,
    1531             :                                                       const std::vector<BoundingBox> & local_bboxes,
    1532             :                                                       const Point & pt,
    1533             :                                                       const unsigned int only_from_mesh_div,
    1534             :                                                       Real & distance) const
    1535             : {
    1536     1836484 :   if (_use_bounding_boxes && !local_bboxes[i_from].contains_point(pt))
    1537     1448998 :     return false;
    1538             :   else
    1539             :   {
    1540      387486 :     auto * pl = _from_point_locators[i_from].get();
    1541      387486 :     const auto from_global_num = getGlobalSourceAppIndex(i_from);
    1542      387486 :     const auto transformed_pt = _from_transforms[from_global_num]->mapBack(pt);
    1543             : 
    1544             :     // Check point against source block restriction
    1545      387486 :     if (!_from_blocks.empty() && !inBlocks(_from_blocks, pl, transformed_pt))
    1546       16987 :       return false;
    1547             : 
    1548             :     // Check point against source boundary restriction. Block restriction will speed up the search
    1549      387290 :     if (!_from_boundaries.empty() &&
    1550           0 :         !onBoundaries(_from_boundaries, _from_blocks, *_from_meshes[i_from], pl, transformed_pt))
    1551           0 :       return false;
    1552             : 
    1553             :     // Check point against the source mesh division
    1554      394054 :     if ((!_from_mesh_divisions.empty() || !_to_mesh_divisions.empty()) &&
    1555        6764 :         !acceptPointMeshDivision(transformed_pt, i_from, only_from_mesh_div))
    1556        3796 :       return false;
    1557             : 
    1558             :     // Get nearest position (often a subapp position) for the target point
    1559             :     // We want values from the child app that is closest to the same position as the target
    1560      383494 :     Point nearest_position_source;
    1561      383494 :     if (_nearest_positions_obj)
    1562             :     {
    1563       26195 :       const bool initial = _fe_problem.getCurrentExecuteOnFlag() == EXEC_INITIAL;
    1564             :       // The search for the nearest position is done in the reference frame
    1565       26195 :       const Point nearest_position = _nearest_positions_obj->getNearestPosition(pt, initial);
    1566       26195 :       nearest_position_source = _nearest_positions_obj->getNearestPosition(
    1567       26195 :           (*_from_transforms[from_global_num])(Point(0, 0, 0)), initial);
    1568             : 
    1569       26195 :       if (!_skip_coordinate_collapsing &&
    1570           0 :           _from_transforms[from_global_num]->hasNonTranslationTransformation())
    1571           0 :         mooseError("Rotation and scaling currently unsupported with nearest positions transfer.");
    1572             : 
    1573             :       // Compute distance to nearest position and nearest position source
    1574       26195 :       const Real distance_to_position_nearest_source = (pt - nearest_position_source).norm();
    1575       26195 :       const Real distance_to_nearest_position = (pt - nearest_position).norm();
    1576             : 
    1577             :       // Source (usually app position) is not closest to the same positions as the target, dont
    1578             :       // send values. We check the distance instead of the positions because if they are the same
    1579             :       // that means there's two equidistant positions and we would want to capture that as a "value
    1580             :       // conflict"
    1581       26195 :       if (!MooseUtils::absoluteFuzzyEqual(distance_to_position_nearest_source,
    1582             :                                           distance_to_nearest_position))
    1583       12835 :         return false;
    1584             : 
    1585             :       // Set the distance as the distance from the nearest position to the target point
    1586       13360 :       distance = distance_to_position_nearest_source;
    1587             :     }
    1588             : 
    1589             :     // Check that the app actually contains the origin point
    1590             :     // We dont need to check if we already found it in a block or a boundary
    1591      370819 :     if (_from_blocks.empty() && _from_boundaries.empty() && _source_app_must_contain_point &&
    1592         160 :         !inMesh(pl, transformed_pt))
    1593         160 :       return false;
    1594             :   }
    1595      370499 :   return true;
    1596             : }
    1597             : 
    1598             : bool
    1599         160 : MultiAppGeneralFieldTransfer::inMesh(const PointLocatorBase * const pl, const Point & point) const
    1600             : {
    1601             :   // Note: we do not take advantage of a potential block restriction of the mesh here. This is
    1602             :   // because we can avoid this routine by calling inBlocks() instead
    1603         160 :   const Elem * elem = (*pl)(point);
    1604         160 :   return (elem != nullptr);
    1605             : }
    1606             : 
    1607             : bool
    1608       47328 : MultiAppGeneralFieldTransfer::inBlocks(const std::set<SubdomainID> & blocks,
    1609             :                                        const Elem * elem) const
    1610             : {
    1611       47328 :   return blocks.find(elem->subdomain_id()) != blocks.end();
    1612             : }
    1613             : 
    1614             : bool
    1615       30400 : MultiAppGeneralFieldTransfer::inBlocks(const std::set<SubdomainID> & blocks,
    1616             :                                        const MooseMesh & /* mesh */,
    1617             :                                        const Elem * elem) const
    1618             : {
    1619       30400 :   return inBlocks(blocks, elem);
    1620             : }
    1621             : 
    1622             : bool
    1623       37832 : MultiAppGeneralFieldTransfer::inBlocks(const std::set<SubdomainID> & blocks,
    1624             :                                        const MooseMesh & mesh,
    1625             :                                        const Node * node) const
    1626             : {
    1627       37832 :   const auto & node_blocks = mesh.getNodeBlockIds(*node);
    1628       37832 :   std::set<SubdomainID> u;
    1629       37832 :   std::set_intersection(blocks.begin(),
    1630             :                         blocks.end(),
    1631             :                         node_blocks.begin(),
    1632             :                         node_blocks.end(),
    1633             :                         std::inserter(u, u.begin()));
    1634       75664 :   return !u.empty();
    1635       37832 : }
    1636             : 
    1637             : bool
    1638        9324 : MultiAppGeneralFieldTransfer::inBlocks(const std::set<SubdomainID> & blocks,
    1639             :                                        const PointLocatorBase * const pl,
    1640             :                                        const Point & point) const
    1641             : {
    1642        9324 :   const Elem * elem = (*pl)(point, &blocks);
    1643        9324 :   return (elem != nullptr);
    1644             : }
    1645             : 
    1646             : bool
    1647      448792 : MultiAppGeneralFieldTransfer::onBoundaries(const std::set<BoundaryID> & boundaries,
    1648             :                                            const MooseMesh & mesh,
    1649             :                                            const Node * node) const
    1650             : {
    1651      448792 :   const BoundaryInfo & bnd_info = mesh.getMesh().get_boundary_info();
    1652      448792 :   std::vector<BoundaryID> vec_to_fill;
    1653      448792 :   bnd_info.boundary_ids(node, vec_to_fill);
    1654      448792 :   std::set<BoundaryID> vec_to_fill_set(vec_to_fill.begin(), vec_to_fill.end());
    1655      448792 :   std::set<BoundaryID> u;
    1656      448792 :   std::set_intersection(boundaries.begin(),
    1657             :                         boundaries.end(),
    1658             :                         vec_to_fill_set.begin(),
    1659             :                         vec_to_fill_set.end(),
    1660             :                         std::inserter(u, u.begin()));
    1661      897584 :   return !u.empty();
    1662      448792 : }
    1663             : 
    1664             : bool
    1665       68576 : MultiAppGeneralFieldTransfer::onBoundaries(const std::set<BoundaryID> & boundaries,
    1666             :                                            const MooseMesh & mesh,
    1667             :                                            const Elem * elem) const
    1668             : {
    1669             :   // Get all boundaries each side of the element is part of
    1670       68576 :   const BoundaryInfo & bnd_info = mesh.getMesh().get_boundary_info();
    1671       68576 :   std::vector<BoundaryID> vec_to_fill;
    1672       68576 :   std::vector<BoundaryID> vec_to_fill_temp;
    1673       68576 :   if (_elemental_boundary_restriction_on_sides)
    1674      424032 :     for (const auto side : make_range(elem->n_sides()))
    1675             :     {
    1676      363456 :       bnd_info.boundary_ids(elem, side, vec_to_fill_temp);
    1677      363456 :       vec_to_fill.insert(vec_to_fill.end(), vec_to_fill_temp.begin(), vec_to_fill_temp.end());
    1678             :     }
    1679             :   else
    1680       72000 :     for (const auto node_index : make_range(elem->n_nodes()))
    1681             :     {
    1682       64000 :       bnd_info.boundary_ids(elem->node_ptr(node_index), vec_to_fill_temp);
    1683       64000 :       vec_to_fill.insert(vec_to_fill.end(), vec_to_fill_temp.begin(), vec_to_fill_temp.end());
    1684             :     }
    1685       68576 :   std::set<BoundaryID> vec_to_fill_set(vec_to_fill.begin(), vec_to_fill.end());
    1686             : 
    1687             :   // Look for a match between the boundaries from the restriction and those near the element
    1688       68576 :   std::set<BoundaryID> u;
    1689       68576 :   std::set_intersection(boundaries.begin(),
    1690             :                         boundaries.end(),
    1691             :                         vec_to_fill_set.begin(),
    1692             :                         vec_to_fill_set.end(),
    1693             :                         std::inserter(u, u.begin()));
    1694      137152 :   return !u.empty();
    1695       68576 : }
    1696             : 
    1697             : bool
    1698           0 : MultiAppGeneralFieldTransfer::onBoundaries(const std::set<BoundaryID> & boundaries,
    1699             :                                            const std::set<SubdomainID> & block_restriction,
    1700             :                                            const MooseMesh & mesh,
    1701             :                                            const PointLocatorBase * const pl,
    1702             :                                            const Point & point) const
    1703             : {
    1704             :   // Find the element containing the point and use the block restriction if known for speed
    1705             :   const Elem * elem;
    1706           0 :   if (block_restriction.empty())
    1707           0 :     elem = (*pl)(point);
    1708             :   else
    1709           0 :     elem = (*pl)(point, &block_restriction);
    1710             : 
    1711           0 :   if (!elem)
    1712           0 :     return false;
    1713           0 :   return onBoundaries(boundaries, mesh, elem);
    1714             : }
    1715             : 
    1716             : bool
    1717        6764 : MultiAppGeneralFieldTransfer::acceptPointMeshDivision(
    1718             :     const Point & pt, const unsigned int i_local, const unsigned int only_from_this_mesh_div) const
    1719             : {
    1720             :   // This routine can also be called to examine if the to_mesh_division index matches the current
    1721             :   // source subapp index
    1722        6764 :   unsigned int source_mesh_div = MooseMeshDivision::INVALID_DIVISION_INDEX - 1;
    1723        6764 :   if (!_from_mesh_divisions.empty())
    1724        4824 :     source_mesh_div = _from_mesh_divisions[i_local]->divisionIndex(pt);
    1725             : 
    1726             :   // If the point is not indexed in the source division
    1727        6764 :   if (!_from_mesh_divisions.empty() && source_mesh_div == MooseMeshDivision::INVALID_DIVISION_INDEX)
    1728        2896 :     return false;
    1729             :   // If the point is not the at the same index in the target and the origin meshes, reject
    1730        3868 :   else if ((_from_mesh_division_behavior == MeshDivisionTransferUse::MATCH_DIVISION_INDEX ||
    1731        3868 :             _to_mesh_division_behavior == MeshDivisionTransferUse::MATCH_DIVISION_INDEX) &&
    1732             :            source_mesh_div != only_from_this_mesh_div)
    1733           0 :     return false;
    1734             :   // If the point is at a certain division index that is not the same as the index of the subapp
    1735             :   // we wanted the information to be from for that point, reject
    1736        3868 :   else if (_from_mesh_division_behavior == MeshDivisionTransferUse::MATCH_SUBAPP_INDEX &&
    1737             :            source_mesh_div != only_from_this_mesh_div)
    1738         462 :     return false;
    1739        4044 :   else if (_to_mesh_division_behavior == MeshDivisionTransferUse::MATCH_SUBAPP_INDEX &&
    1740         638 :            only_from_this_mesh_div != getGlobalSourceAppIndex(i_local))
    1741         438 :     return false;
    1742             :   else
    1743        2968 :     return true;
    1744             : }
    1745             : 
    1746             : bool
    1747      253973 : MultiAppGeneralFieldTransfer::closestToPosition(unsigned int pos_index, const Point & pt) const
    1748             : {
    1749             :   mooseAssert(_nearest_positions_obj, "Should not be here without a positions object");
    1750      253973 :   if (!_skip_coordinate_collapsing)
    1751           0 :     paramError("skip_coordinate_collapsing", "Coordinate collapsing not implemented");
    1752      253973 :   bool initial = _fe_problem.getCurrentExecuteOnFlag() == EXEC_INITIAL;
    1753      253973 :   if (!_search_value_conflicts)
    1754             :     // Faster to just compare the index
    1755      127549 :     return pos_index == _nearest_positions_obj->getNearestPositionIndex(pt, initial);
    1756             :   else
    1757             :   {
    1758             :     // Get the distance to the position and see if we are missing a value just because the position
    1759             :     // is not officially the closest, but it is actually at the same distance
    1760      126424 :     const auto nearest_position = _nearest_positions_obj->getNearestPosition(pt, initial);
    1761      126424 :     const auto nearest_position_at_index = _nearest_positions_obj->getPosition(pos_index, initial);
    1762      126424 :     Real distance_to_position_at_index = (pt - nearest_position_at_index).norm();
    1763      126424 :     const Real distance_to_nearest_position = (pt - nearest_position).norm();
    1764             : 
    1765      126424 :     if (!MooseUtils::absoluteFuzzyEqual(distance_to_position_at_index,
    1766             :                                         distance_to_nearest_position))
    1767       83685 :       return false;
    1768             :     // Actually the same position (point)
    1769       42739 :     else if (nearest_position == nearest_position_at_index)
    1770       42739 :       return true;
    1771             :     else
    1772             :     {
    1773           0 :       mooseWarning("Two equidistant positions ",
    1774             :                    nearest_position,
    1775             :                    " and ",
    1776             :                    nearest_position_at_index,
    1777             :                    " detected near point ",
    1778             :                    pt);
    1779           0 :       return true;
    1780             :     }
    1781             :   }
    1782             : }
    1783             : 
    1784             : Real
    1785     8076933 : MultiAppGeneralFieldTransfer::bboxMaxDistance(const Point & p, const BoundingBox & bbox) const
    1786             : {
    1787     8076933 :   std::array<Point, 2> source_points = {{bbox.first, bbox.second}};
    1788             : 
    1789     8076933 :   std::array<Point, 8> all_points;
    1790    24230799 :   for (unsigned int x = 0; x < 2; x++)
    1791    48461598 :     for (unsigned int y = 0; y < 2; y++)
    1792    96923196 :       for (unsigned int z = 0; z < 2; z++)
    1793    64615464 :         all_points[x + 2 * y + 4 * z] =
    1794   129230928 :             Point(source_points[x](0), source_points[y](1), source_points[z](2));
    1795             : 
    1796     8076933 :   Real max_distance = 0.;
    1797             : 
    1798    72692397 :   for (unsigned int i = 0; i < 8; i++)
    1799             :   {
    1800    64615464 :     Real distance = (p - all_points[i]).norm();
    1801    64615464 :     if (distance > max_distance)
    1802    17179135 :       max_distance = distance;
    1803             :   }
    1804             : 
    1805     8076933 :   return max_distance;
    1806             : }
    1807             : 
    1808             : Real
    1809     8076933 : MultiAppGeneralFieldTransfer::bboxMinDistance(const Point & p, const BoundingBox & bbox) const
    1810             : {
    1811     8076933 :   std::array<Point, 2> source_points = {{bbox.first, bbox.second}};
    1812             : 
    1813     8076933 :   std::array<Point, 8> all_points;
    1814    24230799 :   for (unsigned int x = 0; x < 2; x++)
    1815    48461598 :     for (unsigned int y = 0; y < 2; y++)
    1816    96923196 :       for (unsigned int z = 0; z < 2; z++)
    1817    64615464 :         all_points[x + 2 * y + 4 * z] =
    1818   129230928 :             Point(source_points[x](0), source_points[y](1), source_points[z](2));
    1819             : 
    1820     8076933 :   Real min_distance = std::numeric_limits<Real>::max();
    1821             : 
    1822    72692397 :   for (unsigned int i = 0; i < 8; i++)
    1823             :   {
    1824    64615464 :     Real distance = (p - all_points[i]).norm();
    1825    64615464 :     if (distance < min_distance)
    1826    17665903 :       min_distance = distance;
    1827             :   }
    1828             : 
    1829     8076933 :   return min_distance;
    1830             : }
    1831             : 
    1832             : std::vector<BoundingBox>
    1833       58210 : MultiAppGeneralFieldTransfer::getRestrictedFromBoundingBoxes() const
    1834             : {
    1835       58210 :   std::vector<std::pair<Point, Point>> bb_points(_from_meshes.size());
    1836       58210 :   const Real min_r = std::numeric_limits<Real>::lowest();
    1837       58210 :   const Real max_r = std::numeric_limits<Real>::max();
    1838             : 
    1839      119065 :   for (const auto j : make_range(_from_meshes.size()))
    1840             :   {
    1841       60855 :     Point min(max_r, max_r, max_r);
    1842       60855 :     Point max(min_r, min_r, min_r);
    1843       60855 :     bool at_least_one = false;
    1844       60855 :     const auto & from_mesh = _from_problems[j]->mesh(_displaced_source_mesh);
    1845             : 
    1846      121710 :     for (const auto & elem : as_range(from_mesh.getMesh().local_elements_begin(),
    1847     9673277 :                                       from_mesh.getMesh().local_elements_end()))
    1848             :     {
    1849     4745356 :       if (!_from_blocks.empty() && !inBlocks(_from_blocks, from_mesh, elem))
    1850       17152 :         continue;
    1851             : 
    1852    26622140 :       for (const auto & node : elem->node_ref_range())
    1853             :       {
    1854    21893936 :         if (!_from_boundaries.empty() && !onBoundaries(_from_boundaries, from_mesh, &node))
    1855      161408 :           continue;
    1856             : 
    1857    21732528 :         at_least_one = true;
    1858    86930112 :         for (const auto i : make_range(LIBMESH_DIM))
    1859             :         {
    1860    65197584 :           min(i) = std::min(min(i), node(i));
    1861    65197584 :           max(i) = std::max(max(i), node(i));
    1862             :         }
    1863             :       }
    1864       60855 :     }
    1865             : 
    1866             :     // For 2D RZ problems, we need to amend the bounding box to cover the whole XYZ projection
    1867             :     // - The XYZ-Y axis is assumed aligned with the RZ-Z axis
    1868             :     // - RZ systems also cover negative coordinates hence the use of the maximum R
    1869             :     // NOTE: We will only support the case where there is only one coordinate system
    1870       60855 :     if ((from_mesh.getUniqueCoordSystem() == Moose::COORD_RZ) && (LIBMESH_DIM == 3))
    1871             :     {
    1872         107 :       min(0) = -max(0);
    1873         107 :       min(2) = -max(0);
    1874         107 :       max(2) = max(0);
    1875             :     }
    1876             : 
    1877       60855 :     BoundingBox bbox(min, max);
    1878       60855 :     if (!at_least_one)
    1879          31 :       bbox.min() = max; // If we didn't hit any nodes, this will be _the_ minimum bbox
    1880             :     else
    1881             :     {
    1882             :       // Translate the bounding box to the from domain's position. We may have rotations so we
    1883             :       // must be careful in constructing the new min and max (first and second)
    1884       60824 :       const auto from_global_num = getGlobalSourceAppIndex(j);
    1885       60824 :       transformBoundingBox(bbox, *_from_transforms[from_global_num]);
    1886             :     }
    1887             : 
    1888             :     // Cast the bounding box into a pair of points (so it can be put through
    1889             :     // MPI communication).
    1890       60855 :     bb_points[j] = static_cast<std::pair<Point, Point>>(bbox);
    1891             :   }
    1892             : 
    1893             :   // Serialize the bounding box points.
    1894       58210 :   _communicator.allgather(bb_points);
    1895             : 
    1896             :   // Recast the points back into bounding boxes and return.
    1897       58210 :   std::vector<BoundingBox> bboxes(bb_points.size());
    1898      151168 :   for (const auto i : make_range(bb_points.size()))
    1899       92958 :     bboxes[i] = static_cast<BoundingBox>(bb_points[i]);
    1900             : 
    1901             :   // TODO move up
    1902             :   // Check for a user-set fixed bounding box size and modify the sizes as appropriate
    1903       58210 :   if (_fixed_bbox_size != std::vector<Real>(3, 0))
    1904        1216 :     for (const auto i : make_range(LIBMESH_DIM))
    1905         912 :       if (!MooseUtils::absoluteFuzzyEqual(_fixed_bbox_size[i], 0))
    1906        1932 :         for (const auto j : make_range(bboxes.size()))
    1907             :         {
    1908        1216 :           const auto current_width = (bboxes[j].second - bboxes[j].first)(i);
    1909        1216 :           bboxes[j].first(i) -= (_fixed_bbox_size[i] - current_width) / 2;
    1910        1216 :           bboxes[j].second(i) += (_fixed_bbox_size[i] - current_width) / 2;
    1911             :         }
    1912             : 
    1913      116420 :   return bboxes;
    1914       58210 : }
    1915             : 
    1916             : std::vector<unsigned int>
    1917       58210 : MultiAppGeneralFieldTransfer::getGlobalStartAppPerProc() const
    1918             : {
    1919       58210 :   std::vector<unsigned int> global_app_start_per_proc(1, -1);
    1920       58210 :   if (_from_local2global_map.size())
    1921       32573 :     global_app_start_per_proc[0] = _from_local2global_map[0];
    1922       58210 :   _communicator.allgather(global_app_start_per_proc, true);
    1923       58210 :   return global_app_start_per_proc;
    1924           0 : }
    1925             : 
    1926             : VariableName
    1927       61671 : MultiAppGeneralFieldTransfer::getFromVarName(unsigned int var_index)
    1928             : {
    1929             :   mooseAssert(var_index < _from_var_names.size(), "No source variable at this index");
    1930       61671 :   VariableName var_name = _from_var_names[var_index];
    1931       61671 :   if (_from_var_components.size())
    1932         280 :     var_name += "_" + std::to_string(_from_var_components[var_index]);
    1933       61671 :   return var_name;
    1934           0 : }
    1935             : 
    1936             : VariableName
    1937      117061 : MultiAppGeneralFieldTransfer::getToVarName(unsigned int var_index)
    1938             : {
    1939             :   mooseAssert(var_index < _to_var_names.size(), "No target variable at this index");
    1940      117061 :   VariableName var_name = _to_var_names[var_index];
    1941      117061 :   if (_to_var_components.size())
    1942         352 :     var_name += "_" + std::to_string(_to_var_components[var_index]);
    1943      117061 :   return var_name;
    1944           0 : }
    1945             : 
    1946             : Point
    1947      116668 : MultiAppGeneralFieldTransfer::getMaxToProblemsBBoxDimensions() const
    1948             : {
    1949             :   Point max_dimension = {std::numeric_limits<Real>::min(),
    1950             :                          std::numeric_limits<Real>::min(),
    1951      116668 :                          std::numeric_limits<Real>::min()};
    1952             : 
    1953      239092 :   for (const auto & to_mesh : _to_meshes)
    1954             :   {
    1955      122424 :     const auto bbox = to_mesh->getInflatedProcessorBoundingBox();
    1956      489696 :     for (const auto dim : make_range(LIBMESH_DIM))
    1957      367272 :       max_dimension(dim) = std::max(
    1958      734544 :           max_dimension(dim), std::max(std::abs(bbox.first(dim)), std::abs(bbox.second(dim))));
    1959             :   }
    1960             : 
    1961      116668 :   return max_dimension;
    1962             : }
    1963             : 
    1964             : bool
    1965     2360267 : MultiAppGeneralFieldTransfer::detectConflict(Real current_value,
    1966             :                                              Real new_value,
    1967             :                                              Real current_distance,
    1968             :                                              Real new_distance) const
    1969             : {
    1970             :   // No conflict if we're not looking for them
    1971     2360267 :   if (_search_value_conflicts)
    1972             :     // Only consider conflicts if the values are valid and different
    1973        3680 :     if (current_value != GeneralFieldTransfer::BetterOutOfMeshValue &&
    1974        5654 :         new_value != GeneralFieldTransfer::BetterOutOfMeshValue &&
    1975        1974 :         !MooseUtils::absoluteFuzzyEqual(current_value, new_value))
    1976             :       // Conflict only occurs if the origin points are equidistant
    1977        1767 :       if (MooseUtils::absoluteFuzzyEqual(current_distance, new_distance))
    1978          86 :         return true;
    1979     2360181 :   return false;
    1980             : }

Generated by: LCOV version 1.14