|
Plugins
|
Classes | |
| class | CalcVolPowerElemHelperClass |
| Helper class for the integration of the power. More... | |
| class | CalcVolPowerHelperClass |
| Helper class for the computation of the power of the electromagnetic field. More... | |
| class | ComPol_MemCopy |
| class | ComPol_MemOp |
| class | ComputeElemFluxHelper |
| Helper class for computation of the flux in elements. (Helper for ComputeFlux.) More... | |
| class | ComputeElemFluxHelper< TGridFunc, RegularEdge > |
| Helper class for computation of the flux in elements. (Helper for ComputeFlux.) More... | |
| class | ComputeFluxHelper |
| Helper class for the loop over all the elements in the computation of the flux. (Helper for ComputeFlux.) More... | |
| class | EddyCurrent_E_Nedelec |
| FE-discretization of the time-harmonic E-based formulation of the eddy current model. More... | |
| class | EddyCurrentCurlEDependentCmpUserData |
| class | EddyCurrentHeat |
| class | EddyCurrentImBofEUserData |
| class | EddyCurrentReBofEUserData |
| class | EddyCurrentTraits |
| Auxiliary class defining some important constants. More... | |
| class | EMaterial |
| Class for subdomain-dependent data for the E-based formulated problems. More... | |
| class | EMDirichlet |
| Common interface to get the Dirichlet boundary conditions. More... | |
| struct | Functionality |
| class | NedelecCurlData |
| class | NedelecDirichletBC |
| Dirichlet BC for a rot-rot operator. More... | |
| class | NedelecGridFunctionData |
| class | NedelecInterpolation |
| Interpolation of the Nedelec dofs and their curls. More... | |
| class | NedelecInterpolation< TDomain, 2, 2 > |
| A specialization of NedelecInterpolation for 2d,. More... | |
| class | NedelecInterpolation< TDomain, 3, 3 > |
| A specialization of NedelecInterpolation for 3d,. More... | |
| class | NedelecLoopCurrent |
| Class for computation of loop currents. More... | |
| class | NedelecProject |
| Projection procedure class. More... | |
| class | NedelecProlongationMatrixHelper |
| class | NedelecProlongationMatrixHelper< TDomain, TAlgebra, RegularEdge > |
| class | NedelecSigmaEData |
| class | NedelecT1_LDisc |
| Tool kit for the Whitney-1 (Nedelec) based FE discretization of the rot-rot operators. More... | |
| class | NedelecT1_LDisc< TDomain, Tetrahedron > |
| Specialization of NedelecT1_LDisc for tetrahedra. More... | |
| class | NedelecT1_LDisc< TDomain, Triangle > |
| Specialization of NedelecT1_LDisc for triangles. More... | |
| class | NedelecT1_LDisc_forSimplex |
| Helper class for the specialization of NedelecT1_LDisc for simplices (triangles and tetrahedrons) More... | |
| class | NedelecTransfer |
| Class of the prolongation and the restriction of the Nedelec DoFs. More... | |
| class | TimeHarmonicNedelecHybridSmoother |
| The hybrid smoother by R. Hiptmair. More... | |
Functions | |
| template<typename TGridFunc > | |
| number | calc_EMF (const TGridFunc *gfE, const size_t fct[2], const SubsetGroup &ss_grp, const MathVector< TGridFunc::dim > &Normal, const typename TGridFunc::domain_type::position_type &base_pnt, const size_t n_pnt, const typename TGridFunc::domain_type::position_type &d_pnt, number emf[2]) |
| Computation of the magnetic flux through windings of a coil. | |
| template<typename TGridFunc > | |
| number | calc_magnetic_flux (const TGridFunc *gfE, const size_t fct[2], const SubsetGroup &ss_grp, const MathVector< TGridFunc::dim > &Normal, const typename TGridFunc::domain_type::position_type &base_pnt, const size_t n_pnt, const typename TGridFunc::domain_type::position_type &d_pnt, number flux[2]) |
| Computation of the magnetic flux through windings of a coil. | |
| template<typename TGridFunc > | |
| void | calc_power (TGridFunc *pJGGF, size_t JG_fct[], SubsetGroup &JG_ssg, TGridFunc *pEGF, size_t E_fct[], number pow[]) |
| Computes the power of the electromagnetic field (up to the contribution of the boundary) | |
| template<typename TGridFunc > | |
| void | CalcEMF (SmartPtr< TGridFunc > spGF, const char *cmps, const char *subsets, const std::vector< number > &Normal, const std::vector< number > &base_pnt, const size_t n_pnt, const std::vector< number > &d_pnt) |
| Prints of the electromotive force due to the time variation of the magnetic flux through windings of a coil. | |
| template<typename TGridFunc > | |
| void | CalcMagneticFlux (SmartPtr< TGridFunc > spGF, const char *cmps, const char *subsets, const std::vector< number > &Normal, const std::vector< number > &base_pnt, const size_t n_pnt, const std::vector< number > &d_pnt) |
| Prints of the magnetic flux through windings of a coil. | |
| template<typename TGridFunc > | |
| void | CalcPower (SmartPtr< TGridFunc > spJGGF, const char *JG_cmps, const char *JG_ss, SmartPtr< TGridFunc > spEGF, const char *E_cmps) |
| Prints the (complex-valued) power of the electromagnetic field. | |
| template<typename TGridFunc > | |
| void | ComputeFlux (SmartPtr< TGridFunc > spGF, const char *fct_names, const char *vol_subsets, const char *face_subsets) |
| prints the flux through a given surface in the shell | |
| template<typename TGridFunc > | |
| number | ComputeFlux (TGridFunc *pGF, size_t fct, SubsetGroup &volSSG, SubsetGroup &faceSSG) |
| returns the flux through a given surface | |
| template<typename TDomain > | |
| number | ComputeNedelecDoF (UserData< MathVector< TDomain::dim >, TDomain::dim > *function, Edge *edge, const typename TDomain::position_accessor_type &aaPos, const int si, number time) |
| template<typename TGridFunction > | |
| void | ComputeNedelecDoFs (SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > > spFunction, SmartPtr< TGridFunction > spGridFct, size_t fct, const SubsetGroup &ssGrp, number time) |
| template<int dim> | |
| void | get_ave_vector (const size_t n_vec, const MathVector< dim > *vec, MathVector< dim > &ave) |
| template<typename TArray , typename TLayout , typename TOp > | |
| void | MemAllReduce (TArray *pMem, const TLayout &masterLayout, const TLayout &slaveLayout, pcl::InterfaceCommunicator< TLayout > *pCom=NULL) |
| template<typename TPotGridFunc , typename TGridFunc > | |
| void | NedelecGradPotential (SmartPtr< TPotGridFunc > spPotGF, SmartPtr< TGridFunc > spGF, const char *cmp) |
| template<typename TPotGridFunc , typename TGridFunc > | |
| void | NedelecGradPotential (SmartPtr< TPotGridFunc > spPotGF, SmartPtr< TGridFunc > spGF, size_t fct) |
| template<typename TPotGridFunc > | |
| void | SetSubsetVertVal (SmartPtr< TPotGridFunc > spPotGF, const char *subsets, number value) |
| template<typename TPotGridFunc > | |
| void | SetSubsetVertVal (SmartPtr< TPotGridFunc > spPotGF, const SubsetGroup &ssGrp, number value) |
| template<typename TGridFunction > | |
| void | ComputeNedelecDoFs (SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > > spFunction, SmartPtr< TGridFunction > spGridFct, const char *cmp, const char *subsets, number time) |
| template<typename TGridFunction > | |
| void | ComputeNedelecDoFs (SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > > spFunction, SmartPtr< TGridFunction > spGridFct, const char *cmp, number time) |
| template<typename TGridFunction > | |
| void | ComputeNedelecDoFs (SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > > spFunction, SmartPtr< TGridFunction > spGridFct, const char *cmp, const char *subsets) |
| template<typename TGridFunction > | |
| void | ComputeNedelecDoFs (SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > > spFunction, SmartPtr< TGridFunction > spGridFct, const char *cmp) |
| number ug::Electromagnetism::calc_EMF | ( | const TGridFunc * | gfE, |
| const size_t | fct[2], | ||
| const SubsetGroup & | ss_grp, | ||
| const MathVector< TGridFunc::dim > & | Normal, | ||
| const typename TGridFunc::domain_type::position_type & | base_pnt, | ||
| const size_t | n_pnt, | ||
| const typename TGridFunc::domain_type::position_type & | d_pnt, | ||
| number | emf[2] | ||
| ) |
Computation of the magnetic flux through windings of a coil.
This function computes the magnetic flux through windings of a cylindric coil and returns the electomotive force due to the varying magnetic flux. The flux is a complex number (as well as the electric field that is used as the base of the computation). The function returns the total area of all the windings.
| [in] | gfE | grid function with the electric field |
| [in] | fct | function components (Re and Im) |
| [in] | ss_grp | subsets of the kerner of the coil |
| [in] | Normal | normal to the planes of the windings |
| [in] | base_pnt | point on the plane of the 1st winding |
| [in] | n_pnt | number of the windings |
| [in] | d_pnt | thickness of the windings |
| [out] | emf | the computed electromotive force (Re and Im) |
References ug::Electromagnetism::NedelecInterpolation< TDomain, refDim, WDim >::curl_flux(), ug::DoFRef(), ug::LFEID::NEDELEC, ug::SubsetGroup::size(), UG_THROW, and ug::VecLength().
Referenced by CalcEMF().
| number ug::Electromagnetism::calc_magnetic_flux | ( | const TGridFunc * | gfE, |
| const size_t | fct[2], | ||
| const SubsetGroup & | ss_grp, | ||
| const MathVector< TGridFunc::dim > & | Normal, | ||
| const typename TGridFunc::domain_type::position_type & | base_pnt, | ||
| const size_t | n_pnt, | ||
| const typename TGridFunc::domain_type::position_type & | d_pnt, | ||
| number | flux[2] | ||
| ) |
Computation of the magnetic flux through windings of a coil.
| [in] | gfE | grid function with the electric field |
| [in] | fct | function components (Re and Im) |
| [in] | ss_grp | subsets of the kerner of the coil |
| [in] | Normal | normal to the planes of the windings |
| [in] | base_pnt | point on the plane of the 1st winding |
| [in] | n_pnt | number of the windings |
| [in] | d_pnt | thickness of the windings |
| [out] | flux | the computed flux (Re and Im) |
| void ug::Electromagnetism::calc_power | ( | TGridFunc * | pJGGF, |
| size_t | JG_fct[], | ||
| SubsetGroup & | JG_ssg, | ||
| TGridFunc * | pEGF, | ||
| size_t | E_fct[], | ||
| number | pow[] | ||
| ) |
Computes the power of the electromagnetic field (up to the contribution of the boundary)
| pJGGF | grid function of the Nedelec-DoFs of the generator current \(\mathbf{J}_G\) |
| JG_fct | indices of the Re and Im parts in the grid function for \(\mathbf{J}_G\) |
| JG_ssg | (full-dim.) subsets where \(\mathbf{J}_G\) is defined (non-zero and in the kernel) |
| pEGF | grid function of the Nedelec-DoFs of the electric field \(\mathbf{E}\) |
| E_fct | indices of the Re and Im parts in the grid function for \(\mathbf{E}\) |
| pow | to add the integral |
Referenced by CalcPower().
| void ug::Electromagnetism::CalcEMF | ( | SmartPtr< TGridFunc > | spGF, |
| const char * | cmps, | ||
| const char * | subsets, | ||
| const std::vector< number > & | Normal, | ||
| const std::vector< number > & | base_pnt, | ||
| const size_t | n_pnt, | ||
| const std::vector< number > & | d_pnt | ||
| ) |
Prints of the electromotive force due to the time variation of the magnetic flux through windings of a coil.
| [in] | spGF | grid function with the electric field |
| [in] | cmps | names of the components of the grid function (for Re and Im) |
| [in] | subsets | subsets where to compute the flux |
| [in] | Normal | normal to the planes of the windings |
| [in] | base_pnt | point on the plane of the 1st winding |
| [in] | n_pnt | number of the windings |
| [in] | d_pnt | thickness of the windings |
References ug::SubsetGroup::add(), calc_EMF(), SmartPtr< typename T, template< class TT > class FreePolicy >::get(), ug::RemoveWhitespaceFromString(), ug::TokenizeString(), UG_CATCH_THROW, UG_LOG, and UG_THROW.
| void ug::Electromagnetism::CalcMagneticFlux | ( | SmartPtr< TGridFunc > | spGF, |
| const char * | cmps, | ||
| const char * | subsets, | ||
| const std::vector< number > & | Normal, | ||
| const std::vector< number > & | base_pnt, | ||
| const size_t | n_pnt, | ||
| const std::vector< number > & | d_pnt | ||
| ) |
Prints of the magnetic flux through windings of a coil.
| [in] | spGF | grid function with the electric field |
| [in] | cmps | names of the components of the grid function (for Re and Im) |
| [in] | subsets | subsets where to compute the flux |
| [in] | Normal | normal to the planes of the windings |
| [in] | base_pnt | point on the plane of the 1st winding |
| [in] | n_pnt | number of the windings |
| [in] | d_pnt | thickness of the windings |
| void ug::Electromagnetism::CalcPower | ( | SmartPtr< TGridFunc > | spJGGF, |
| const char * | JG_cmps, | ||
| const char * | JG_ss, | ||
| SmartPtr< TGridFunc > | spEGF, | ||
| const char * | E_cmps | ||
| ) |
Prints the (complex-valued) power of the electromagnetic field.
| [in] | spJGGF | grid function with the generator current |
| [in] | JG_cmps | names of the components of the grid function (for Re and Im) |
| JG_ss | (full-dim.) subsets where \(\mathbf{J}_G\) is defined (non-zero and in the kernel) | |
| [in] | spEGF | grid function with the electric field |
| [in] | E_cmps | names of the components of the grid function (for Re and Im) |
References ug::SubsetGroup::add(), calc_power(), SmartPtr< typename T, template< class TT > class FreePolicy >::get(), ug::RemoveWhitespaceFromString(), ug::SubsetGroup::set_subset_handler(), ug::TokenizeString(), UG_CATCH_THROW, UG_LOG, and UG_THROW.
| void ug::Electromagnetism::ComputeFlux | ( | SmartPtr< TGridFunc > | spGF, |
| const char * | fct_name, | ||
| const char * | vol_subsets, | ||
| const char * | face_subsets | ||
| ) |
prints the flux through a given surface in the shell
Computation of the flux of a vector field (given by the Whitney-1 elements) through a low-dimensional subset (surface). The function prints the value of the flux in the shell.
| spGF | grid function of the Nedelec-DoFs of the vector field |
| fct_name | name of the function in the grid function |
| vol_subsets | full-dim. subsets (adjacent to the surface) to indicate the negative direction |
| face_subsets | the surface (the low-dim. subsets) |
References ug::SubsetGroup::add(), ComputeFlux(), SmartPtr< typename T, template< class TT > class FreePolicy >::get(), ug::RemoveWhitespaceFromString(), ug::SubsetGroup::set_subset_handler(), ug::FunctionGroup::size(), ug::TokenizeString(), UG_CATCH_THROW, UG_LOG, and UG_THROW.
| number ug::Electromagnetism::ComputeFlux | ( | TGridFunc * | pGF, |
| size_t | fct, | ||
| SubsetGroup & | volSSG, | ||
| SubsetGroup & | faceSSG | ||
| ) |
returns the flux through a given surface
Computation of the flux of a vector field (given by the Whitney-1 elements) through a low-dimensional subset (surface). The function returns the value of the flux.
| pGF | grid function of the Nedelec-DoFs of the vector field |
| fct | index of the function in the grid function |
| volSSG | full-dim. subsets (adjacent to the surface) to indicate the negative direction |
| faceSSG | the surface (the low-dim. subsets) |
References ug::LFEID::NEDELEC, and UG_THROW.
Referenced by ComputeFlux().
|
inline |
Computes the value for the DoF associated with a given edge. Uses the trapezoid quadrature rule.
| function | the function to encode |
| edge | the edge to compute the DoF for |
| aaPos | coordinates of the vertices |
| si | subset index |
| time | the time argument |
References ug::EdgeVertices::vertex().
| void ug::Electromagnetism::ComputeNedelecDoFs | ( | SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > > | spFunction, |
| SmartPtr< TGridFunction > | spGridFct, | ||
| const char * | cmp | ||
| ) |
Computes the values of the DoFs of the Whitney-1 (Nedelec-type-1) basis functions for a given vector function. This function parses the string specifications of the function components and the subsets.
| spFunction | the function to encode |
| spGridFct | the grid function to store the result |
| cmp | name of the component |
References ComputeNedelecDoFs().
| void ug::Electromagnetism::ComputeNedelecDoFs | ( | SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > > | spFunction, |
| SmartPtr< TGridFunction > | spGridFct, | ||
| const char * | cmp, | ||
| const char * | subsets | ||
| ) |
Computes the values of the DoFs of the Whitney-1 (Nedelec-type-1) basis functions for a given vector function. This function parses the string specifications of the function components and the subsets.
| spFunction | the function to encode |
| spGridFct | the grid function to store the result |
| cmp | name of the component |
| subsets | names of the subsets (may be NULL meaning "all") |
References ComputeNedelecDoFs().
| void ug::Electromagnetism::ComputeNedelecDoFs | ( | SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > > | spFunction, |
| SmartPtr< TGridFunction > | spGridFct, | ||
| const char * | cmp, | ||
| const char * | subsets, | ||
| number | time | ||
| ) |
Computes the values of the DoFs of the Whitney-1 (Nedelec-type-1) basis functions for a given vector function. This function parses the string specifications of the function components and the subsets.
| spFunction | the function to encode |
| spGridFct | the grid function to store the result |
| cmp | name of the component |
| subsets | names of the subsets (may be NULL meaning "all") |
| time | the time argument |
References ug::SubsetGroup::add(), ug::SubsetGroup::add_all(), ComputeNedelecDoFs(), ug::LFEID::NEDELEC, ug::TokenizeString(), and UG_THROW.
| void ug::Electromagnetism::ComputeNedelecDoFs | ( | SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > > | spFunction, |
| SmartPtr< TGridFunction > | spGridFct, | ||
| const char * | cmp, | ||
| number | time | ||
| ) |
Computes the values of the DoFs of the Whitney-1 (Nedelec-type-1) basis functions for a given vector function. This function parses the string specifications of the function components and the subsets.
| spFunction | the function to encode |
| spGridFct | the grid function to store the result |
| cmp | name of the component |
| time | the time argument |
References ComputeNedelecDoFs().
| void ug::Electromagnetism::ComputeNedelecDoFs | ( | SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > > | spFunction, |
| SmartPtr< TGridFunction > | spGridFct, | ||
| size_t | fct, | ||
| const SubsetGroup & | ssGrp, | ||
| number | time | ||
| ) |
Computes the values of the DoFs of the Whitney-1 (Nedelec-type-1) basis functions for a given vector function.
| spFunction | the function to encode |
| spGridFct | the grid function to store the result |
| fct | index of the function component in the grid function |
| ssGrp | subset group where to encode |
| time | the time argument (for the function to encode) |
References ug::DoFRef(), ug::SubsetGroup::size(), UG_ASSERT, UG_LOG, and UG_THROW.
Referenced by ComputeNedelecDoFs(), ComputeNedelecDoFs(), ComputeNedelecDoFs(), and ComputeNedelecDoFs().
|
inline |
Computes the arithm. average of given MathVectors.
| [in] | n_vec | number of the vectors in the array |
| [in] | vec | array of the vectors |
| [out] | ave | the average |
| void ug::Electromagnetism::MemAllReduce | ( | TArray * | pMem, |
| const TLayout & | masterLayout, | ||
| const TLayout & | slaveLayout, | ||
| pcl::InterfaceCommunicator< TLayout > * | pCom = NULL |
||
| ) |
Syncronizes a memory array using a given reduction operation.
| TArray | type of the array |
| TLayout | type of the parallel layout |
| TOp | class of the reduction operation |
| pMem | the array to syncronize |
| masterLayout | the master layout |
| slaveLayout | the slave layout |
| pCom | (optionally) the communicator |
References pcl::InterfaceCommunicator< class TLayout >::communicate(), pcl::InterfaceCommunicator< class TLayout >::receive_data(), and pcl::InterfaceCommunicator< class TLayout >::send_data().
| void ug::Electromagnetism::NedelecGradPotential | ( | SmartPtr< TPotGridFunc > | spPotGF, |
| SmartPtr< TGridFunc > | spGF, | ||
| const char * | cmp | ||
| ) |
Computes the Whitney representation of the gradient of a vertex-centered scalar grid function. (This is the version with the symbolic function name.)
| spPotGF | the potential to compute the gradient for |
| spGF | to save the Whitney form of the gradient |
| cmp | name of the component |
References NedelecGradPotential(), and UG_THROW.
| void ug::Electromagnetism::NedelecGradPotential | ( | SmartPtr< TPotGridFunc > | spPotGF, |
| SmartPtr< TGridFunc > | spGF, | ||
| size_t | fct | ||
| ) |
Computes the Whitney representation of the gradient of a vertex-centered scalar grid function.
| spPotGF | the potential to compute the gradient for |
| spGF | to save the Whitney form of the gradient |
| fct | index of the function in spGF |
References ug::DoFRef(), ug::LFEID::NEDELEC, UG_THROW, and ug::EdgeVertices::vertex().
Referenced by NedelecGradPotential().
| void ug::Electromagnetism::SetSubsetVertVal | ( | SmartPtr< TPotGridFunc > | spPotGF, |
| const char * | subsets, | ||
| number | value | ||
| ) |
Sets a scalar potential field to a given value at vertices of a given subset, and to 0 for all other subsets.
| spPotGF | the potential to compute the gradient for |
| subsets | names of the subsets (may be NULL meaning "all") |
| value | value to set in the subset group |
References ug::SubsetGroup::add(), ug::SubsetGroup::add_all(), SetSubsetVertVal(), and ug::TokenizeString().
| void ug::Electromagnetism::SetSubsetVertVal | ( | SmartPtr< TPotGridFunc > | spPotGF, |
| const SubsetGroup & | ssGrp, | ||
| number | value | ||
| ) |
Sets a scalar potential field to a given value at vertices of a given subset, and to 0 for all other subsets.
| spPotGF | the potential to compute the gradient for |
| ssGrp | subset group where to set |
| value | value to set in the subset group |
References ug::SubsetGroup::size(), and UG_THROW.
Referenced by SetSubsetVertVal().