Projects from \(H(\text{curl})\) to \(H(\text{div})\). More...
#include <project_Hcurl_to_Hdiv.hpp>
Public Member Functions | |
| ProjectHcurlToHdiv (unsigned int p, unsigned int mapping_degree, const Triangulation< 3 > &triangulation_Hcurl, const DoFHandler< 3 > &dof_handler_Hcurl, const Vector< double > &solution_Hcurl, std::string fname="Hdiv", const Function< 3 > *exact_solution=nullptr, bool print_time_tables=false, bool project_exact_solution=false, bool log_cg_convergence=false, bool write_higher_order_cells=false) | |
| The only constructor. More... | |
| double | get_L2_norm () |
| Returns \(L^2\) error norm. | |
| double | get_Linfty_norm () |
| Returns \(L^{\infty}\) error norm. | |
| unsigned int | get_n_cells () const |
| Returns the number of active cells in the mesh. | |
| unsigned int | get_n_dofs () const |
| Returns the total amount of the degrees of freedom. | |
| void | clear () |
| Releases computer memory associated with system matrix and right-hand side. | |
| const Triangulation< 3 > & | get_tria () const |
| Returns a reference to triangulation. | |
| const DoFHandler< 3 > & | get_dof_handler () const |
| Returns a reference to dof handler associated with the Raviart-Thomas finite elements. | |
| const Vector< double > & | get_solution () const |
| Returns a reference to solution, i.e., the result of the projection. | |
| void | save_matrix_and_rhs_to_csv (std::string fname) const |
| Saves the system matrix and the right-hand side into a csv file. More... | |
Projects from \(H(\text{curl})\) to \(H(\text{div})\).
This class template is not supposed to be used directly. Instead, one of the two wrap-around class templates, StaticVectorSolver::ProjectAtoB or StaticVectorSolver::ProjectTtoJ, must be used. The names of the wrap-around class templates are assumed to be more familiar to readers in electromagnetics.
Implements the following recipes:
This class template is supposed to be used in pair with StaticVectorSolver::Solver1 or StaticVectorSolver::Solver2. In all problems formulated in terms of the magnetic vector potential, \(\vec{A}\), the numerically calculated potential needs to be converted into magnetic field, \(\vec{B}\), as
\[ \vec{B} = \vec{\nabla} \times \vec{A}. \]
The magnetic vector potential belongs to the \(H(\text{curl})\) function space. The magnetic field belongs to the \(H(\text{div})\) function space. Therefore, one needs to compute the equation above such that the input, \(\vec{A}\), is in \(H(\text{curl})\) and the output, \(\vec{B}\), is in \(H(\text{div})\). Such computation can be envisioned as a some kind of projection from one function space, \(H(\text{curl})\), into another, \(H(\text{div})\).
An identical projection arises in problems in which a numerically computed current vector potential, \(\vec{T}\), needs to be converted back into the free-current density, \(\vec{J}_f \),
\[ \vec{J}_f = \vec{\nabla} \times \vec{T}, \]
for a comparison with a closed-form analytical expression of \(\vec{J}_f\) given by the formulation of the problem.
This class template does these two projections. The Bossavit's diagrams below illustrate them.
The user is supposed to derive an object from this class template. All usual computations, i.e., setup, assembling the linear system, etc., happen automatically.
Definition at line 106 of file project_Hcurl_to_Hdiv.hpp.
| StaticVectorSolver::ProjectHcurlToHdiv< stage >::ProjectHcurlToHdiv | ( | unsigned int | p, |
| unsigned int | mapping_degree, | ||
| const Triangulation< 3 > & | triangulation_Hcurl, | ||
| const DoFHandler< 3 > & | dof_handler_Hcurl, | ||
| const Vector< double > & | solution_Hcurl, | ||
| std::string | fname = "Hdiv", |
||
| const Function< 3 > * | exact_solution = nullptr, |
||
| bool | print_time_tables = false, |
||
| bool | project_exact_solution = false, |
||
| bool | log_cg_convergence = false, |
||
| bool | write_higher_order_cells = false |
||
| ) |
The only constructor.
| [in] | p | - The degree of the FE_RaviartThomas finite elements. |
| [in] | mapping_degree | - The degree of the interpolating polynomials used for mapping. Setting it to 1 will do in the most of the cases. Note, that it makes sense to attach a meaningful manifold to the triangulation if this parameter is greater than 1. |
| [in] | triangulation_Hcurl | - Reference to the triangulation inside the object of the class StaticVectorSolver::Solver1 or StaticVectorSolver::Solver2 that has yielded the vector potential in the \(H(\text{curl})\) function space. |
| [in] | dof_handler_Hcurl | - Reference to the dof handler inside the class StaticVectorSolver::Solver1 or StaticVectorSolver::Solver2 that has yielded the vector potential in the \(H(\text{curl})\) function space. |
| [in] | solution_Hcurl | - Vector filled with the degrees of freedom that together with the shape functions of the Nedelec finite elements model the vector potential in the \(H(\text{curl})\) function space. |
| [in] | fname | - The name of the output files without extension. |
| [in] | exact_solution | - Points to an object that describes the exact solution to the problem. It is needed for calculating error norms. |
| [in] | print_time_tables | - If true, prints time tables on the screen. |
| [in] | project_exact_solution | - If true, projects the exact solution onto the space spanned by the Raviart-Thomas finite elements and saves the result into the output file next to the solution. This may be useful for debugging purposes. |
| [in] | log_cg_convergence | - If true, logs convergence of the conjugate gradient solver into a file. The name of the file is generated by appending "_cg_convergence.csv" to fname. |
| [in] | write_higher_order_cells | - If false, the data is saved in the file fname.vtk. Higher-order cells are not saved. If true, the data is saved into fname.vtu file preserving the higher-order cells. In this case the file can be viewed with a help of Paraview version 5.5.0 or higher. |
Definition at line 327 of file project_Hcurl_to_Hdiv.hpp.
| void StaticVectorSolver::ProjectHcurlToHdiv< stage >::save_matrix_and_rhs_to_csv | ( | std::string | fname | ) | const |
Saves the system matrix and the right-hand side into a csv file.
All the zeros included into the csv files. This is a very dumb and inefficient way of saving sparse matrices. On the positive side - it is very easy and straightforward to read the csv files. This function may be useful for debugging. One can assemble the system on a coarse mesh (so there are a few mesh cells and the system matrix is small) and export the system matrix together with the right-hand side into another program such as Matlab or GNU Octave for an analysis.
| [in] | fname | - A stem of the names of the output files. The matrix will be saved into fname_matrix.csv file. The right-hand side will be save into fname_rhs.csv file. |
Definition at line 397 of file project_Hcurl_to_Hdiv.hpp.