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 "MultiAppMFEMTolibMeshShapeEvaluationTransfer.h" 13 : #include "MFEMVectorUtils.h" 14 : #include "libmesh/mesh_function.h" 15 : 16 : registerMooseObject("MooseApp", MultiAppMFEMTolibMeshShapeEvaluationTransfer); 17 : 18 : InputParameters 19 2626 : MultiAppMFEMTolibMeshShapeEvaluationTransfer::validParams() 20 : { 21 2626 : InputParameters params = MFEMMultiAppTransfer::validParams(); 22 2626 : params.addClassDescription("Transfers variable values from an MFEM based application to a " 23 : "libMesh application, using shape function evaluations."); 24 2626 : return params; 25 0 : } 26 : 27 264 : MultiAppMFEMTolibMeshShapeEvaluationTransfer::MultiAppMFEMTolibMeshShapeEvaluationTransfer( 28 264 : InputParameters const & params) 29 264 : : MFEMMultiAppTransfer(params), _mfem_interpolator(this->comm().get()) 30 : { 31 264 : checkValidTransferProblemTypes<FEProblemBase, MFEMProblem>(); 32 264 : } 33 : 34 : MFEMProblem & 35 336 : MultiAppMFEMTolibMeshShapeEvaluationTransfer::getActiveFromProblem() 36 : { 37 336 : return static_cast<MFEMProblem &>(MFEMMultiAppTransfer::getActiveFromProblem()); 38 : } 39 : 40 : /// Extract locations from libMesh-based MooseVariable at which projection will take place 41 : void 42 330 : MultiAppMFEMTolibMeshShapeEvaluationTransfer::extractlibMeshNodePositions( 43 : libMesh::System & to_sys, 44 : const MooseVariableFieldBase & to_var, 45 : std::vector<Point> & outgoing_libmesh_points) 46 : { 47 : // libMesh mesh 48 330 : const MeshBase & to_mesh = getActiveToProblem().mesh(_displaced_target_mesh).getMesh(); 49 330 : auto var_num = to_var.number(); 50 330 : auto sys_num = to_sys.number(); 51 330 : auto & fe_type = to_var.feType(); 52 330 : bool is_nodal = to_var.isNodal(); 53 : 54 : // Populate set of points 55 330 : if (fe_type.order > CONSTANT && !is_nodal) 56 0 : mooseError("Transfers of non-nodal FEs from MFEM to libMesh with order higher than " 57 : "CONSTANT are not supported."); 58 330 : else if (is_nodal) 59 : { 60 : // Consider nodes the variable has dofs at. 61 144048 : for (const auto & node : to_mesh.local_node_ptr_range()) 62 71904 : if (node->n_dofs(sys_num, var_num)) 63 43165 : outgoing_libmesh_points.push_back(*node); 64 : } 65 : else // Elemental, constant monomial 66 : { 67 : // Consider elements the variable has dofs at. 68 16630 : for (const auto & elem : as_range(to_mesh.local_elements_begin(), to_mesh.local_elements_end())) 69 8270 : if (elem->n_dofs(sys_num, var_num)) 70 8360 : outgoing_libmesh_points.push_back(elem->vertex_average()); 71 : } 72 330 : } 73 : 74 : /// Project interpolated values to set DoFs of libMesh-based MooseVariable 75 : void 76 330 : MultiAppMFEMTolibMeshShapeEvaluationTransfer::projectlibMeshNodalValues( 77 : libMesh::System & to_sys, const MooseVariableFieldBase & to_var, mfem::Vector & interp_vals) 78 : { 79 : // libMesh mesh 80 330 : const MeshBase & to_mesh = getActiveToProblem().mesh(_displaced_target_mesh).getMesh(); 81 330 : auto var_num = to_var.number(); 82 330 : auto sys_num = to_sys.number(); 83 330 : auto & fe_type = to_var.feType(); 84 330 : bool is_nodal = to_var.isNodal(); 85 : 86 330 : unsigned int mfem_point_index = 0; 87 330 : if (fe_type.order > CONSTANT && !is_nodal) 88 0 : mooseError("Transfers of non-nodal FEs from MFEM to libMesh with order higher than " 89 : "CONSTANT are not supported."); 90 330 : else if (is_nodal) 91 : { 92 : // Consider nodes the variable has dofs at. 93 144048 : for (const auto & node : to_mesh.local_node_ptr_range()) 94 71904 : if (node->n_dofs(sys_num, var_num)) 95 : { 96 42925 : const auto dof_object_id = node->id(); 97 42925 : const DofObject * dof_object = to_mesh.node_ptr(dof_object_id); 98 42925 : const auto dof = dof_object->dof_number(sys_num, var_num, 0); 99 42925 : const auto val = interp_vals[mfem_point_index]; 100 42925 : to_sys.solution->set(dof, val); 101 42925 : mfem_point_index++; 102 240 : } 103 : } 104 : else // Elemental, constant monomial 105 : { 106 : // Consider elements the variable has dofs at. 107 16630 : for (const auto & elem : as_range(to_mesh.local_elements_begin(), to_mesh.local_elements_end())) 108 8270 : if (elem->n_dofs(sys_num, var_num)) 109 : { 110 8270 : const auto dof_object_id = elem->id(); 111 8270 : const DofObject * dof_object = to_mesh.elem_ptr(dof_object_id); 112 8270 : const auto dof = dof_object->dof_number(sys_num, var_num, 0); 113 8270 : const auto val = interp_vals[mfem_point_index]; 114 8270 : to_sys.solution->set(dof, val); 115 8270 : mfem_point_index++; 116 90 : } 117 : } 118 330 : to_sys.solution->close(); 119 : // Sync local solutions 120 330 : to_sys.update(); 121 330 : } 122 : 123 : void 124 336 : MultiAppMFEMTolibMeshShapeEvaluationTransfer::transferVariables(bool is_target_local) 125 : { 126 : // Send from MFEM problem to libMesh problem 127 672 : for (const auto var_index : make_range(numToVar())) 128 : { 129 : // Declare source variable alias and construct vectors to store interpolation points and vals 130 336 : auto & from_var = *getActiveFromProblem().getGridFunction(getFromVarName(var_index)); 131 336 : if (from_var.VectorDim() > 1) 132 0 : mooseError("MultiAppMFEMTolibMeshShapeEvaluationTransfer does not support transfers of " 133 : "vector variables from MFEM to libMesh-based subapps"); 134 336 : from_var.ParFESpace()->GetParMesh()->EnsureNodes(); 135 : // Ordering of interpolated MFEM output, chosen to match libMesh ordering 136 336 : const mfem::Ordering::Type ordering = mfem::Ordering::byVDIM; 137 336 : mfem::Vector interp_vals; 138 336 : mfem::Vector outgoing_mfem_points; 139 : 140 : // If target libMesh app is on the current rank, update the set of interpolation 141 : // points with the libMesh variable's set 142 336 : if (is_target_local) 143 : { 144 : // Declare aliases for convenience 145 330 : FEProblemBase & to_problem = getActiveToProblem(); 146 : const MooseVariableFieldBase & to_var( 147 990 : to_problem.getVariable(0, 148 330 : getToVarName(var_index), 149 : Moose::VarKindType::VAR_ANY, 150 : Moose::VarFieldType::VAR_FIELD_ANY)); 151 330 : const auto & to_var_name = getToVarName(var_index); 152 330 : auto & es = getlibMeshEquationSystem(to_problem, _displaced_target_mesh); 153 330 : System & to_sys = *find_sys(es, to_var_name); 154 : // Extract set of target points in libMesh mesh to perform interpolation of MFEM variable at 155 330 : std::vector<Point> outgoing_libmesh_points; 156 330 : extractlibMeshNodePositions(to_sys, to_var, outgoing_libmesh_points); 157 330 : const MeshBase & to_mesh = to_problem.mesh(_displaced_target_mesh).getMesh(); 158 : // Perform interpolation of MFEM variable 159 51525 : for (auto & point : outgoing_libmesh_points) 160 51195 : point = mapPointToActiveSourceFrame(point); 161 660 : outgoing_mfem_points = Moose::MFEM::libMeshPointsToMFEMVector( 162 330 : outgoing_libmesh_points, to_mesh.mesh_dimension(), ordering); 163 330 : } 164 : 165 : // MFEM interpolation using FindPointsGSLib::Interpolate must be over all ranks 166 : // of the MFEM variable's MPI communicator 167 336 : _mfem_interpolator.SetDefaultInterpolationValue(getMFEMOutOfMeshValue()); 168 336 : _mfem_interpolator.Interpolate(*from_var.ParFESpace()->GetParMesh(), 169 : outgoing_mfem_points, 170 : from_var, 171 : interp_vals, 172 : ordering); 173 : 174 : // If target libMesh variable exists on the current rank, update the set of local 175 : // DoFs of the target libMesh variable using the interpolated MFEM values 176 336 : if (is_target_local) 177 : { 178 : // Declare aliases for convenience 179 330 : FEProblemBase & to_problem = getActiveToProblem(); 180 : const MooseVariableFieldBase & to_var( 181 990 : to_problem.getVariable(0, 182 330 : getToVarName(var_index), 183 : Moose::VarKindType::VAR_ANY, 184 : Moose::VarFieldType::VAR_FIELD_ANY)); 185 330 : const auto & to_var_name = getToVarName(var_index); 186 : // Get undisplaced system for writing transferred variable to 187 330 : auto & out_es = getlibMeshEquationSystem(to_problem, false); 188 330 : auto & out_sys = *find_sys(out_es, to_var_name); 189 : 190 : // Project interpolated values at destination nodes onto destination variables to set DoFs 191 330 : projectlibMeshNodalValues(out_sys, to_var, interp_vals); 192 : } 193 336 : } 194 336 : } 195 : 196 : #endif