26 int main(
int argc,
char *argv[])
31 const char *source_mesh_file =
"../../data/inline-tri.mesh";
32 const char *destination_mesh_file =
"../../data/inline-quad.mesh";
34 int src_n_refinements = 0;
35 int dest_n_refinements = 0;
36 int source_fe_order = 1;
37 int dest_fe_order = 1;
38 bool visualization =
true;
39 bool use_vector_fe =
false;
43 args.
AddOption(&source_mesh_file,
"-s",
"--source_mesh",
44 "Mesh file to use for src.");
45 args.
AddOption(&destination_mesh_file,
"-d",
"--destination_mesh",
46 "Mesh file to use for dest.");
47 args.
AddOption(&src_n_refinements,
"-sr",
"--source_refinements",
48 "Number of src refinements");
49 args.
AddOption(&dest_n_refinements,
"-dr",
"--dest_refinements",
50 "Number of dest refinements");
51 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
53 "Enable or disable GLVis visualization.");
54 args.
AddOption(&source_fe_order,
"-so",
"--source_fe_order",
55 "Order of the src finite elements");
56 args.
AddOption(&dest_fe_order,
"-do",
"--dest_fe_order",
57 "Order of the dest finite elements");
58 args.
AddOption(&verbose,
"-verb",
"--verbose",
"--no-verb",
"--no-verbose",
59 "Enable/Disable verbose output");
60 args.
AddOption(&use_vector_fe,
"-vfe",
"--use_vector_fe",
"-no-vfe",
61 "--no-vector_fe",
"Use vector finite elements");
65 shared_ptr<Mesh> src_mesh, dest_mesh;
69 imesh.open(destination_mesh_file);
72 dest_mesh = make_shared<Mesh>(imesh, 1, 1);
77 mfem::err <<
"WARNING: Destination mesh file not found: " 78 << destination_mesh_file <<
"\n" 79 <<
"Using default 2D quad mesh.";
91 imesh.open(source_mesh_file);
95 src_mesh = make_shared<Mesh>(imesh, 1, 1);
100 mfem::err <<
"WARNING: Source mesh file not found: " << source_mesh_file
102 <<
"Using default box mesh.\n";
115 for (
int i = 0; i < src_mesh->
GetNV(); ++i)
119 for (
int d = 0; d <
dim; ++d)
126 for (
int i = 0; i < src_n_refinements; ++i)
131 for (
int i = 0; i < dest_n_refinements; ++i)
136 shared_ptr<FiniteElementCollection> src_fe_coll, dest_fe_coll;
141 make_shared<RT_FECollection>(source_fe_order, src_mesh->
Dimension());
143 make_shared<RT_FECollection>(dest_fe_order, dest_mesh->
Dimension());
148 make_shared<L2_FECollection>(source_fe_order, src_mesh->
Dimension());
150 make_shared<L2_FECollection>(dest_fe_order, dest_mesh->
Dimension());
154 make_shared<FiniteElementSpace>(src_mesh.get(), src_fe_coll.get());
157 make_shared<FiniteElementSpace>(dest_mesh.get(), dest_fe_coll.get());
171 src_fun.ProjectCoefficient(vector_coeff);
176 src_fun.ProjectCoefficient(coeff);
195 if (assembler.
Transfer(src_fun, dest_fun))
206 src_err = src_fun.ComputeL2Error(vector_coeff);
207 dest_err = dest_fun.ComputeL2Error(vector_coeff);
211 src_err = src_fun.ComputeL2Error(coeff);
212 dest_err = dest_fun.ComputeL2Error(coeff);
215 mfem::out <<
"l2 error: src: " << src_err <<
", dest: " << dest_err
218 plot(*src_mesh, src_fun,
"source");
219 plot(*dest_mesh, dest_fun,
"destination");
224 mfem::out <<
"No intersection -> no transfer!" << std::endl;
Class for grid function - Vector with associated FE space.
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
void SetVerbose(const bool verbose)
Expose process details with verbose output.
int Dimension() const
Dimension of the reference space used within the elements.
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Returns the minimum and maximum corners of the mesh bounding box.
This class implements the serial variational transfer between finite element spaces. Variational transfer has been shown to have better approximation properties than standard interpolation. This facilities can be used for supporting applications which require the handling of non matching meshes. For instance: General multi-physics problems, fluid structure interaction, or even visualization of average quantities within subvolumes This algorithm allows to perform quadrature in the intersection of elements of two separate and unrelated meshes. It generates quadrature rules in the intersection which allows us to integrate-with to machine precision using the mfem::MortarIntegrator interface. See https://doi.org/10.1137/15M1008361 for and in-depth explanation. At this time curved elements are not supported.
int main(int argc, char *argv[])
void check_options(mfem::OptionsParser &args)
void InitTransfer(int argc, char *argv[])
Initializes the par_moonolith library. It also calls MPI_Init.
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
bool Transfer(const GridFunction &src_fun, GridFunction &dest_fun)
transfer a function from source to destination. if the transfer is to be performed multiple times use...
const double * GetVertex(int i) const
Return pointer to vertex i's coordinates.
A general vector function coefficient.
void plot(mfem::Mesh &mesh, mfem::GridFunction &x, std::string title)
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Add a boolean option and set 'var' to receive the value. Enable/disable tags are used to set the bool...
double example_fun(const mfem::Vector &x)
void vector_fun(const mfem::Vector &x, mfem::Vector &f)
void AddMortarIntegrator(const std::shared_ptr< MortarIntegrator > &integrator)
This method must be called before Assemble or Transfer. It will assemble the operator in all intersec...
A general function coefficient.
int FinalizeTransfer()
Finalize the par_moonolith library.