Logbook  (07-04-2025)
Static problems
StaticScalarSolver::ProjectAphiToHrz< stage > Class Template Reference

Computes the two-dimensional scaled auxiliary field \(\vec{H}'\) from scaled magnetic vector potential, \(A'\), i.e., \(\vec{H}' = - \mu^{-1} \vec{\nabla} \overset{V}{\times} A'\). More...

#include <project_Aphi_to_Hrz.hpp>

Inheritance diagram for StaticScalarSolver::ProjectAphiToHrz< stage >:
Collaboration diagram for StaticScalarSolver::ProjectAphiToHrz< stage >:

Public Member Functions

 ProjectAphiToHrz (unsigned int p, unsigned int mapping_degree, const Triangulation< 2 > &triangulation_A, const DoFHandler< 2 > &dof_handler_A, const Vector< double > &solution_A, std::string fname="Hrz", const Function< 2 > *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...
 
- Public Member Functions inherited from StaticScalarSolver::ProjectHgradToHcurl< 2, 1 >
 ProjectHgradToHcurl (unsigned int p, unsigned int mapping_degree, const Triangulation< dim > &triangulation_Hgrad, const DoFHandler< dim > &dof_handler_Hgrad, const Vector< double > &solution_Hgrad, const std::string fname="Hcurl", const Function< dim > *exact_solution=nullptr, bool axisymmetric=false, bool vector_potential=false, 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< dim > & get_tria () const
 Returns a reference to triangulation.
 
const DoFHandler< dim > & get_dof_handler () const
 Returns a reference to dof handler associated with the Nedelec 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...
 

Detailed Description

template<int stage = 1>
class StaticScalarSolver::ProjectAphiToHrz< stage >

Computes the two-dimensional scaled auxiliary field \(\vec{H}'\) from scaled magnetic vector potential, \(A'\), i.e., \(\vec{H}' = - \mu^{-1} \vec{\nabla} \overset{V}{\times} A'\).

The problem is assumed to be axisymmetric two-dimensional. The scaled auxiliary field \(\vec{H}'\) is an in-plane vector. The scaled magnetic vector potential, \(A'\), is an out-of-plane vector.

The scaled magnetic vector potential, \(A'\), is assumed to be an output of an object of the class StaticScalarSolver::Solver. That is, \(A'\) is expressed as a linear combination of the shape functions of the FE_Q finite elements. In other words, \(A'\) is in the \(H(\text{grad})\) function space. This class template computes scaled H-field such that \(\vec{H}\) is expressed as a linear combination of the shape functions of the FE_Nedelec finite elements. That is, \(\vec{H}'\) is in the \(H(\text{curl})\) function space. In other words, this class template takes a scalar function \(A'\) (an out-of-plane vector is, essentially a scalar) from the \(H(\text{grad})\) function space, calculates its vector curl, multiplies by \(-\mu^{-1}\), and projects the result into \(H(\text{curl})\) function space, which is a proper space for the scaled auxiliary H-field. The operation of projection is illustrated by means of the Bossavit's diagram below.

It is assumed that the object of the class StaticScalarSolver::Solver that yielded \(A'\) is still in the computer memory such that the triangulation, the degrees of freedom, and the handler of the degrees of freedom are accessible while \(\vec{H}'\) is computed. An object of the StaticScalarSolver::ProjectAzToHxy class template reuses the triangulation by creating an additional dof handler associated with the FE_Nedelec finite elements. That is, \(\vec{H}'\) and \(A'\) share the same triangulation but are modeled by different finite elements: \(\vec{H}'\) - by FE_Nedelec and \(A'\) - by FE_Q. In this disposition there are two DoFHandler objects associated with one Triangulation object. The algorithm walks synchronously through both DoFHandler objects and assembles the system matrix and the right-hand side. The algorithm utilizes the WorkStream technology of deal.II. The amount of threads used can be limited as the following

#include <deal.II/base/multithread_info.h>
MultithreadInfo::set_thread_limit(nr_threads_max);

The output data are saved into a vtk file. The following data are saved:

  • The calculated H'-field, \(\vec{H}' = - \mu^{-1} \vec{\nabla} \overset{V}{\times} A'\), under the name "VectorField".
  • The \(L^2\) and \(L^{\infty}\) error norms associated with the calculated H'-field under the names "L2norm" and "LinftyNorm". One value per mesh cell is saved.
  • The exact solution projected onto \(H(\text{curl})\) function space is saved under the name "VectorFieldExact". The "VectorField" and "VectorFieldExact" are modeled by exactly the same finite elements, i.e., FE_Nedelec.

The "L2norm", "LinftyNorm", and "VectorFieldExact" are saved only if an exact solution is passed to the constructor. Moreover, "VectorFieldExact" is calculated and saved only if the input parameter project_exact_solution passed to the constructor equals true.

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.

Note
Application examples:

Definition at line 101 of file project_Aphi_to_Hrz.hpp.

Constructor & Destructor Documentation

◆ ProjectAphiToHrz()

template<int stage = 1>
StaticScalarSolver::ProjectAphiToHrz< stage >::ProjectAphiToHrz ( unsigned int  p,
unsigned int  mapping_degree,
const Triangulation< 2 > &  triangulation_A,
const DoFHandler< 2 > &  dof_handler_A,
const Vector< double > &  solution_A,
std::string  fname = "Hrz",
const Function< 2 > *  exact_solution = nullptr,
bool  print_time_tables = false,
bool  project_exact_solution = false,
bool  log_cg_convergence = false,
bool  write_higher_order_cells = false 
)
inline

The only constructor.

Constructs the object and runs the calculations. That is, there is no need to call other functions such as run().

Parameters
[in]p- The degree of the FE_Nedelec 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_A- Reference to the triangulation inside an object of the class StaticScalarSolver::Solver that has yielded the scaled magnetic vector potential, \(A'\).
[in]dof_handler_A- Reference to the dof handler inside the class StaticScalarSolver::Solver that has yielded the scaled magnetic vector potential, \(A'\).
[in]solution_A- Vector filled with the degrees of freedom that together with the shape functions of the FE_Q finite elements model the scaled magnetic vector potential, \(A'\).
[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. The exact solution object must exist in the computer memory at the time of calling this constructor.
[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 Nedelec finite elements and saves the result into the output file next to \(\vec{H}'\).
[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 145 of file project_Aphi_to_Hrz.hpp.


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