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 : #ifdef MOOSE_MFEM_ENABLED 11 : 12 : #include "MultiApplibMeshToMFEMShapeEvaluationTransfer.h" 13 : #include "MFEMVectorUtils.h" 14 : #include "libmesh/mesh_function.h" 15 : #include "libmesh/parallel_algebra.h" // for communicator send and receive stuff 16 : 17 : registerMooseObject("MooseApp", MultiApplibMeshToMFEMShapeEvaluationTransfer); 18 : 19 : InputParameters 20 2632 : MultiApplibMeshToMFEMShapeEvaluationTransfer::validParams() 21 : { 22 2632 : InputParameters params = MFEMMultiAppTransfer::validParams(); 23 2632 : params.addClassDescription("Transfers variable values from a libMesh based application to an " 24 : "MFEM application, using shape function evaluations."); 25 2632 : return params; 26 0 : } 27 : 28 267 : MultiApplibMeshToMFEMShapeEvaluationTransfer::MultiApplibMeshToMFEMShapeEvaluationTransfer( 29 267 : InputParameters const & params) 30 267 : : MFEMMultiAppTransfer(params) 31 : { 32 267 : checkValidTransferProblemTypes<MFEMProblem, FEProblemBase>(); 33 267 : } 34 : 35 : MFEMProblem & 36 598 : MultiApplibMeshToMFEMShapeEvaluationTransfer::getActiveToProblem() 37 : { 38 598 : return static_cast<MFEMProblem &>(MFEMMultiAppTransfer::getActiveToProblem()); 39 : } 40 : 41 : void 42 303 : MultiApplibMeshToMFEMShapeEvaluationTransfer::transferVariables(bool is_target_local) 43 : { 44 606 : for (const auto var_index : make_range(numToVar())) 45 : { 46 : // Declare map of processor ID to corresponding vector of libMesh points 47 : // on that processor to interpolate source libMesh variable at 48 303 : std::map<processor_id_type, std::vector<Point>> outgoing_points; 49 303 : mfem::Vector interp_vals; 50 303 : if (is_target_local) 51 : { 52 : // Generate list of points where the grid function will be evaluated 53 : mfem::ParGridFunction & to_gf = 54 299 : *getActiveToProblem().getGridFunction(getToVarName(var_index)); 55 299 : mfem::ParFiniteElementSpace & to_pfespace = *to_gf.ParFESpace(); 56 299 : if (to_gf.VectorDim() > 1) 57 0 : mooseError("MultiApplibMeshToMFEMShapeEvaluationTransfer does not support transfers of " 58 : "vector variables from libMesh to MFEM-based subapps"); 59 299 : mfem::Vector vxyz; 60 : mfem::Ordering::Type point_ordering; 61 299 : _mfem_projector.extractNodePositions(to_pfespace, vxyz, point_ordering); 62 : 63 : // Populate outgoing point locations map between processor and points vector for libMesh to 64 : // use in interpolation 65 299 : const int dim = to_pfespace.GetParMesh()->Dimension(); 66 299 : const int nnodes = vxyz.Size() / dim; 67 154999 : for (const auto i : make_range(nnodes)) 68 : { 69 154700 : libMesh::Point point_in_target_frame; 70 544452 : for (const auto d : make_range(dim)) 71 389752 : point_in_target_frame(d) = vxyz[i + d * nnodes]; 72 : 73 154700 : const auto point_in_source_frame = mapPointToActiveSourceFrame(point_in_target_frame); 74 369916 : for (const auto i_proc : make_range(n_processors())) 75 215216 : outgoing_points[i_proc].push_back(point_in_source_frame); 76 : } 77 299 : interp_vals.SetSize(nnodes); 78 299 : } 79 : 80 : // Perform interpolation of libMesh variable at specified points 81 303 : interpolatelibMeshVariable(outgoing_points, var_index, interp_vals); 82 : 83 : // Project DoFs to MFEM GridFunction 84 303 : if (is_target_local) 85 : { 86 : mfem::ParGridFunction & to_gf = 87 299 : *getActiveToProblem().getGridFunction(getToVarName(var_index)); 88 299 : mfem::Ordering::Type libmesh_interp_ordering(mfem::Ordering::Type::byNODES); 89 299 : _mfem_projector.projectNodalValues(interp_vals, libmesh_interp_ordering, to_gf); 90 : } 91 303 : } 92 303 : } 93 : 94 : void 95 303 : MultiApplibMeshToMFEMShapeEvaluationTransfer::interpolatelibMeshVariable( 96 : std::map<processor_id_type, std::vector<Point>> & outgoing_points, 97 : const unsigned int var_index, 98 : mfem::Vector & interp_vals) 99 : { 100 303 : FEProblemBase & from_problem = getActiveFromProblem(); 101 909 : MooseVariableFieldBase & from_var = from_problem.getVariable(0, 102 303 : getFromVarName(var_index), 103 : Moose::VarKindType::VAR_ANY, 104 : Moose::VarFieldType::VAR_FIELD_ANY); 105 303 : System & from_sys = from_var.sys().system(); 106 303 : unsigned int from_var_num = from_sys.variable_number(getFromVarName(var_index)); 107 : // Construct a local mesh function for each origin problem 108 : libMesh::MeshFunction local_meshfuns( 109 303 : getlibMeshEquationSystem(from_problem, _displaced_source_mesh), 110 303 : *from_sys.current_local_solution, 111 303 : from_sys.get_dof_map(), 112 303 : from_var_num); 113 303 : local_meshfuns.init(); 114 303 : local_meshfuns.enable_out_of_mesh_mode(getMFEMOutOfMeshValue()); 115 : 116 : // Evaluate interpolated values at incoming points. 117 : auto gather_functor = 118 463 : [this, &local_meshfuns](processor_id_type /*pid*/, 119 : const std::vector<Point> & incoming_points, 120 : std::vector<mfem::real_t> & vals_for_incoming_points) 121 : { 122 463 : vals_for_incoming_points.resize(incoming_points.size(), getMFEMOutOfMeshValue()); 123 : // Compute interpolation values of the libMesh variable at all requested points 124 215679 : for (const auto i_pt : index_range(incoming_points)) 125 215216 : vals_for_incoming_points[i_pt] = local_meshfuns(incoming_points[i_pt]); 126 463 : }; 127 : // Copy data out to interp_vals 128 : auto action_functor = 129 463 : [&interp_vals, this](processor_id_type /*pid*/, 130 : const std::vector<Point> & /*my_outgoing_points*/, 131 : const std::vector<mfem::real_t> & vals_for_outgoing_points) 132 : { 133 215679 : for (const auto i : make_range(interp_vals.Size())) 134 : { 135 215216 : const auto val = vals_for_outgoing_points[i]; 136 215216 : if (val == getMFEMOutOfMeshValue()) 137 71431 : continue; 138 143785 : interp_vals(i) = val; 139 : } 140 463 : }; 141 : 142 : // Set interpolated field values at points on local processor 143 303 : interp_vals = getMFEMOutOfMeshValue(); // default to the out-of-mesh value 144 : // We assume incoming_vals is ordered in the same way as outgoing_points 145 303 : libMesh::Parallel::pull_parallel_vector_data( 146 303 : comm(), outgoing_points, gather_functor, action_functor, (mfem::real_t *)(nullptr)); 147 303 : } 148 : 149 : #endif