LCOV - code coverage report
Current view: top level - src/mfem/transfers - MultiAppMFEMShapeEvaluationTransfer.C (source / functions) Hit Total Coverage
Test: idaholab/moose framework: #32971 (54bef8) with base c6cf66 Lines: 66 69 95.7 %
Date: 2026-05-29 20:35:17 Functions: 7 7 100.0 %
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             : #ifdef MOOSE_MFEM_ENABLED
      11             : 
      12             : #include "MultiAppMFEMShapeEvaluationTransfer.h"
      13             : 
      14             : registerMooseObject("MooseApp", MultiAppMFEMShapeEvaluationTransfer);
      15             : 
      16             : InputParameters
      17        2270 : MultiAppMFEMShapeEvaluationTransfer::validParams()
      18             : {
      19        2270 :   InputParameters params = MFEMMultiAppTransfer::validParams();
      20        2270 :   params.addClassDescription("Transfers variable values from one MFEM application to another using "
      21             :                              "shape function evaluations.");
      22        2270 :   return params;
      23           0 : }
      24             : 
      25          86 : MultiAppMFEMShapeEvaluationTransfer::MultiAppMFEMShapeEvaluationTransfer(
      26          86 :     InputParameters const & params)
      27          86 :   : MFEMMultiAppTransfer(params)
      28             : {
      29          86 :   checkValidTransferProblemTypes<MFEMProblem, MFEMProblem>();
      30          86 : }
      31             : 
      32             : MFEMProblem &
      33         262 : MultiAppMFEMShapeEvaluationTransfer::getActiveFromProblem()
      34             : {
      35         262 :   return static_cast<MFEMProblem &>(MFEMMultiAppTransfer::getActiveFromProblem());
      36             : }
      37             : 
      38             : MFEMProblem &
      39         493 : MultiAppMFEMShapeEvaluationTransfer::getActiveToProblem()
      40             : {
      41         493 :   return static_cast<MFEMProblem &>(MFEMMultiAppTransfer::getActiveToProblem());
      42             : }
      43             : 
      44             : void
      45         138 : MultiAppMFEMShapeEvaluationTransfer::transferVariables(bool is_target_local)
      46             : {
      47             :   auto transformTargetPointsToSourceFrame =
      48         243 :       [this](mfem::Vector & point_coordinates, const unsigned int dimension)
      49             :   {
      50             :     mooseAssert(dimension == 2 || dimension == 3, "Target finite element space must be 2D or 3D");
      51             : 
      52         243 :     const auto n_points = point_coordinates.Size() / dimension;
      53      123715 :     for (const auto i : make_range(n_points))
      54             :     {
      55      123472 :       Point point_in_target_frame;
      56      467312 :       for (const auto d : make_range(dimension))
      57      343840 :         point_in_target_frame(d) = point_coordinates[i + d * n_points];
      58             : 
      59      123472 :       const auto point_in_source_frame = mapPointToActiveSourceFrame(point_in_target_frame);
      60      467312 :       for (const auto d : make_range(dimension))
      61      343840 :         point_coordinates[i + d * n_points] = point_in_source_frame(d);
      62             :     }
      63         243 :   };
      64             : 
      65             :   // Get GridFunction from problem by name. For complex variables, return the real component.
      66         741 :   auto getGridFunction = [&](MFEMProblem & problem,
      67             :                              const std::string & name,
      68             :                              bool & is_complex) -> mfem::ParGridFunction &
      69             :   {
      70         741 :     if (problem.getProblemData().gridfunctions.Has(name))
      71             :     {
      72         720 :       is_complex = false;
      73         720 :       return *problem.getGridFunction(name);
      74             :     }
      75          21 :     if (problem.getProblemData().cmplx_gridfunctions.Has(name))
      76             :     {
      77          21 :       is_complex = true;
      78          21 :       return problem.getComplexGridFunction(name)->real();
      79             :     }
      80           0 :     mooseError("No real or complex variable named '", name, "' found.");
      81         138 :   };
      82             : 
      83         393 :   for (const auto v : make_range(numToVar()))
      84             :   {
      85         255 :     bool is_from_complex{false};
      86             :     mfem::ParGridFunction & from_gf =
      87         255 :         getGridFunction(getActiveFromProblem(), getFromVarName(v), is_from_complex);
      88         255 :     mfem::ParFiniteElementSpace & from_pfespace = *from_gf.ParFESpace();
      89         255 :     from_pfespace.GetParMesh()->EnsureNodes();
      90             : 
      91         255 :     bool is_to_complex{false};
      92             :     // Vector to store position coords of nodes/interpolation pts
      93         255 :     mfem::Vector vxyz;
      94             :     // Ordering of node coordinates x,y,z
      95         255 :     mfem::Ordering::Type point_ordering = mfem::Ordering::Type::byNODES;
      96             :     // Ordering of (vector) GridFunction data
      97         255 :     mfem::Ordering::Type to_gf_ordering = mfem::Ordering::Type::byVDIM;
      98             :     // Vector to store GridFunction values at interp. pts
      99         255 :     mfem::Vector interp_vals;
     100             : 
     101             :     // Get points from target GridFunction on local rank
     102         255 :     if (is_target_local)
     103             :     {
     104             :       mfem::ParGridFunction & to_gf =
     105         243 :           getGridFunction(getActiveToProblem(), getToVarName(v), is_to_complex);
     106         243 :       mfem::ParFiniteElementSpace & to_pfespace = *to_gf.ParFESpace();
     107         243 :       to_pfespace.GetParMesh()->EnsureNodes();
     108             :       // Generate list of points where the grid function will be evaluated
     109         243 :       _mfem_projector.extractNodePositions(to_pfespace, vxyz, point_ordering);
     110             :       // Evaluate source grid function at target points
     111         243 :       const int dim = to_pfespace.GetParMesh()->Dimension();
     112         243 :       transformTargetPointsToSourceFrame(vxyz, dim);
     113             :       // Update ordering for interpolation
     114         243 :       to_gf_ordering = to_pfespace.GetOrdering();
     115             :     }
     116             : 
     117             :     // Evaluate source grid function at target points.
     118             :     // Interpolation must be performed over all source ranks.
     119         255 :     _mfem_interpolator.Setup(*from_gf.ParFESpace()->GetParMesh());
     120         255 :     _mfem_interpolator.SetDefaultInterpolationValue(getMFEMOutOfMeshValue());
     121         255 :     _mfem_interpolator.Interpolate(vxyz, from_gf, interp_vals, point_ordering, to_gf_ordering);
     122             : 
     123             :     // Project interpolated field onto target GridFunction data on the local rank
     124         255 :     if (is_target_local)
     125             :     {
     126             :       mfem::ParGridFunction & to_gf =
     127         243 :           getGridFunction(getActiveToProblem(), getToVarName(v), is_to_complex);
     128         243 :       mfem::ParFiniteElementSpace & to_pfespace = *to_gf.ParFESpace();
     129         243 :       _mfem_projector.projectNodalValues(interp_vals, to_gf_ordering, to_gf);
     130         243 :       if (is_to_complex)
     131             :       {
     132             :         // Get remaining imaginary component of destination GridFunction
     133             :         mfem::ParGridFunction & to_gf_im =
     134           7 :             getActiveToProblem().getComplexGridFunction(getToVarName(v))->imag();
     135           7 :         if (is_from_complex)
     136             :         {
     137             :           mfem::ParGridFunction & from_gf_im =
     138           7 :               getActiveFromProblem().getComplexGridFunction(getFromVarName(v))->imag();
     139           7 :           _mfem_interpolator.Interpolate(
     140           7 :               vxyz, from_gf_im, interp_vals, point_ordering, to_pfespace.GetOrdering());
     141           7 :           _mfem_projector.projectNodalValues(interp_vals, to_pfespace.GetOrdering(), to_gf_im);
     142             :         }
     143             :         else // Transfer from real variable to complex variable, so imag component zero
     144           0 :           to_gf_im = 0.0;
     145             :       }
     146             :     }
     147         255 :   }
     148         138 : }
     149             : 
     150             : #endif

Generated by: LCOV version 1.14