libMesh
Public Member Functions | Static Public Member Functions | Private Attributes | List of all members
libMesh::InterMeshProjection Class Reference

This class implements inter mesh projection, i.e. More...

#include <inter_mesh_projection.h>

Public Member Functions

 InterMeshProjection (System &_from_system, System &_to_mesh)
 
void project_system_vectors ()
 

Static Public Member Functions

static Number fptr (const Point &p, const Parameters &, const std::string &libmesh_dbg_var(sys_name), const std::string &unknown_name)
 
static Gradient gptr (const Point &p, const Parameters &, const std::string &libmesh_dbg_var(sys_name), const std::string &unknown_name)
 

Private Attributes

Systemfrom_system
 
Systemto_system
 

Detailed Description

This class implements inter mesh projection, i.e.

projection of vectors defined on a given mesh (from_mesh associated with from_system) to another mesh (to_mesh of to_system).

Definition at line 48 of file inter_mesh_projection.h.

Constructor & Destructor Documentation

◆ InterMeshProjection()

libMesh::InterMeshProjection::InterMeshProjection ( System _from_system,
System _to_mesh 
)

Definition at line 36 of file inter_mesh_projection.C.

36  :
37  from_system(_from_system),
38  to_system(_to_system)
39 {
40  libmesh_experimental();
41 }

Member Function Documentation

◆ fptr()

static Number libMesh::InterMeshProjection::fptr ( const Point p,
const Parameters ,
const std::string &  libmesh_dbg_varsys_name,
const std::string &  unknown_name 
)
static

◆ gptr()

static Gradient libMesh::InterMeshProjection::gptr ( const Point p,
const Parameters ,
const std::string &  libmesh_dbg_varsys_name,
const std::string &  unknown_name 
)
static

Referenced by project_system_vectors().

◆ project_system_vectors()

void libMesh::InterMeshProjection::project_system_vectors ( )

Definition at line 43 of file inter_mesh_projection.C.

References libMesh::NumericVector< T >::build(), libMesh::ParallelObject::comm(), from_system, libMesh::System::get_dof_map(), libMesh::System::get_equation_systems(), libMesh::System::get_vector(), gptr(), libMesh::GradientMeshFunction::init(), libMesh::MeshFunction::init(), libMesh::NumericVector< T >::localize(), n_vars, libMesh::System::n_vars(), libMesh::System::n_vectors(), libMesh::System::project_vector(), libMesh::SERIAL, libMesh::NumericVector< T >::size(), libMesh::System::solution, to_system, libMesh::System::update_global_solution(), libMesh::System::variable_name(), libMesh::System::vector_is_adjoint(), libMesh::System::vectors_begin(), and libMesh::System::vectors_end().

44 {
45  // Number of vectors to be projected
46  libmesh_assert_equal_to (to_system.n_vectors(), from_system.n_vectors());
47 
48  // Number of variable components in each vector
49  unsigned int n_vars = from_system.n_vars();
50  libmesh_assert_equal_to (to_system.n_vars(), n_vars);
51 
52  // We are going to use the multi-variable MeshFunction, so we can pass
53  // a single vector of variables rather than have a MeshFunction for each variable
54  std::vector<unsigned int> variables_vector;
55 
56  for (unsigned int j = 0; j != n_vars; ++j)
57  {
58  libmesh_assert_equal_to (to_system.variable_name(j), from_system.variable_name(j));
59  variables_vector.push_back(j);
60  }
61 
62  // Any system holds the solution along with the other vectors system.vectors
63  // We will first project the solution and then move to the system.vectors
64 
65  // Construct local version of the current system vector
66  // This has to be a serial vector
67  // Roy's FIXME: Technically it just has to be a ghosted vector whose algebraically
68  // ghosted values cover a domain which is a superset of the to-system's domain ...
69  // that's hard to do and we can skip it until the poor scalability bites someone.
70  std::unique_ptr<NumericVector<Number>> solution_vector_serial = NumericVector<Number>::build(from_system.comm());
71  solution_vector_serial->init(from_system.solution->size(), true, SERIAL);
72 
73  std::vector<Number> solution_vector;
74  from_system.update_global_solution(solution_vector);
75  (*solution_vector_serial) = solution_vector;
76 
77  // Construct a MeshFunction for the solution
78  MeshFunction mesh_func_solution(from_system.get_equation_systems(), *solution_vector_serial, from_system.get_dof_map(), variables_vector);
79 
80  mesh_func_solution.init();
81 
82  // For some element types (say C1) we also need to pass a gradient evaluation MeshFunction
83  // To do this evaluate, a new shim class GradientMeshFunction has been added which redirects
84  // gptr::operator evaluations inside projection methods into MeshFunction::gradient calls.
85  GradientMeshFunction gptr_solution(mesh_func_solution);
86  gptr_solution.init();
87 
88  to_system.project_vector(*to_system.solution, &mesh_func_solution, &gptr_solution);
89 
90  // Now loop over the vectors in system.vectors (includes old_nonlin_sol, rhs, adjoints, adjoint_rhs, sensitivity_rhs)
91  for (System::vectors_iterator vec = from_system.vectors_begin(), vec_end = from_system.vectors_end(); vec != vec_end; ++vec)
92  {
93  // The name of this vector
94  const std::string & vec_name = vec->first;
95 
96  // Construct local version of the current system vector
97  // This has to be a serial vector
98  // Roy's FIXME: Technically it just has to be a ghosted vector whose algebraically
99  // ghosted values cover a domain which is a superset of the to-system's domain ...
100  // that's hard to do and we can skip it until the poor scalability bites someone.
101  std::unique_ptr<NumericVector<Number>> current_vector_proxy = NumericVector<Number>::build(from_system.comm());
102  current_vector_proxy->init(from_system.get_vector(vec_name).size(), true, SERIAL);
103 
104  from_system.get_vector(vec_name).localize(*current_vector_proxy);
105 
106  // Construct a MeshFunction for the current component
107  MeshFunction mesh_func(from_system.get_equation_systems(), *current_vector_proxy, from_system.get_dof_map(), variables_vector);
108  mesh_func.init();
109 
110  // Project the current system vector, you need to check if this vector is an adjoint to pass
111  // the right options to project_vector
112  GradientMeshFunction gptr(mesh_func);
113  gptr.init();
114 
115  // The fourth argument here is whether the vector is an adjoint solution or not, we will be getting that information
116  // via the from_system instead of the to_system in case the user has not set the System::_vector_is_adjoint map to true.
117  to_system.project_vector(to_system.get_vector(vec_name), &mesh_func, &gptr, from_system.vector_is_adjoint(vec_name));
118 
119  }
120  // End loop over the vectors in the system
121 
122 }
std::map< std::string, std::unique_ptr< NumericVector< Number > >, std::less<> >::iterator vectors_iterator
Vector iterator typedefs.
Definition: system.h:757
vectors_iterator vectors_end()
End of vectors container.
Definition: system.h:2576
void update_global_solution(std::vector< Number > &global_soln) const
Fill the input vector global_soln so that it contains the global solution on all processors.
Definition: system.C:745
virtual numeric_index_type size() const =0
int vector_is_adjoint(std::string_view vec_name) const
Definition: system.C:1173
const EquationSystems & get_equation_systems() const
Definition: system.h:721
const Parallel::Communicator & comm() const
vectors_iterator vectors_begin()
Beginning of vectors container.
Definition: system.h:2564
unsigned int n_vars
static Gradient gptr(const Point &p, const Parameters &, const std::string &libmesh_dbg_var(sys_name), const std::string &unknown_name)
unsigned int n_vectors() const
Definition: system.h:2558
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1593
const std::string & variable_name(const unsigned int i) const
Definition: system.h:2478
static std::unique_ptr< NumericVector< T > > build(const Parallel::Communicator &comm, SolverPackage solver_package=libMesh::default_solver_package(), ParallelType parallel_type=AUTOMATIC)
Builds a NumericVector on the processors in communicator comm using the linear solver package specifi...
void project_vector(NumericVector< Number > &new_vector, FunctionBase< Number > *f, FunctionBase< Gradient > *g=nullptr, int is_adjoint=-1) const
Projects arbitrary functions onto a vector of degree of freedom values for the current system...
unsigned int n_vars() const
Definition: system.h:2430
const DofMap & get_dof_map() const
Definition: system.h:2374
const NumericVector< Number > & get_vector(std::string_view vec_name) const
Definition: system.C:943
virtual void localize(std::vector< T > &v_local) const =0
Creates a copy of the global vector in the local vector v_local.

Member Data Documentation

◆ from_system

System& libMesh::InterMeshProjection::from_system
private

Definition at line 66 of file inter_mesh_projection.h.

Referenced by project_system_vectors().

◆ to_system

System& libMesh::InterMeshProjection::to_system
private

Definition at line 69 of file inter_mesh_projection.h.

Referenced by project_system_vectors().


The documentation for this class was generated from the following files: