HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
Loading...
Searching...
No Matches
Classes | Typedefs | Functions | Variables
DDR_magnetostatic

Implementation of the DDR scheme for the magnetostatic problem. More...

Classes

struct  HArDCore3D::Magnetostatics
 Assemble a magnetostatic problem. More...
 

Typedefs

typedef Eigen::SparseMatrix< doubleHArDCore3D::Magnetostatics::SystemMatrixType
 
typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Magnetostatics::ForcingTermType
 
typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Magnetostatics::SolutionPotentialType
 
typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Magnetostatics::SolutionCurlType
 
typedef IntegralWeight HArDCore3D::Magnetostatics::PermeabilityType
 
typedef Eigen::SparseMatrix< doubleHArDCore3D::Magnetostatics::SystemMatrixType
 
typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Magnetostatics::ForcingTermType
 
typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Magnetostatics::SolutionPotentialType
 
typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Magnetostatics::SolutionCurlType
 
typedef IntegralWeight HArDCore3D::Magnetostatics::PermeabilityType
 

Functions

 HArDCore3D::Magnetostatics::Magnetostatics (const DDRCore &ddrcore, bool use_threads, std::ostream &output=std::cout)
 Constructor.
 
void HArDCore3D::Magnetostatics::assembleLinearSystem (const ForcingTermType &f, const PermeabilityType &mu, const SolutionPotentialType &u)
 Assemble the global system

 
size_t HArDCore3D::Magnetostatics::dimensionSpace () const
 Returns the dimension of the magnetic field + potential space.
 
size_t HArDCore3D::Magnetostatics::nbSCDOFs () const
 Returns the number of statically condensed DOFs (here, the cell magnetic field DOFs)
 
size_t HArDCore3D::Magnetostatics::sizeSystem () const
 Returns the size of the statically condensed system.
 
const XCurlHArDCore3D::Magnetostatics::xCurl () const
 Returns the space XCurl.
 
const XDivHArDCore3D::Magnetostatics::xDiv () const
 Returns the space XDiv.
 
const SystemMatrixTypeHArDCore3D::Magnetostatics::systemMatrix () const
 Returns the linear system matrix.
 
SystemMatrixTypeHArDCore3D::Magnetostatics::systemMatrix ()
 Returns the linear system matrix.
 
const Eigen::VectorXd & HArDCore3D::Magnetostatics::systemVector () const
 Returns the linear system right-hand side vector.
 
Eigen::VectorXd & HArDCore3D::Magnetostatics::systemVector ()
 Returns the linear system right-hand side vector.
 
const doubleHArDCore3D::Magnetostatics::stabilizationParameter () const
 Returns the stabilization parameter.
 
doubleHArDCore3D::Magnetostatics::stabilizationParameter ()
 Returns the stabilization parameter.
 
const SystemMatrixTypeHArDCore3D::Magnetostatics::scMatrix () const
 Returns the static condensation recovery operator.
 
Eigen::VectorXd & HArDCore3D::Magnetostatics::scVector ()
 Returns the static condensation rhs.
 
std::vector< doubleHArDCore3D::Magnetostatics::computeNorms (const std::vector< Eigen::VectorXd > &list_dofs) const
 Compute the discrete Hcurl \times Hdiv norm of a family of vectors representing the dofs.
 
const SXCurlHArDCore3D::Magnetostatics::sxCurl () const
 Returns the space SXCurl.
 
const SXDivHArDCore3D::Magnetostatics::sxDiv () const
 Returns the space SXDiv.
 

Variables

static const double HArDCore3D::PI = boost::math::constants::pi<double>()
 
static Magnetostatics::SolutionPotentialType HArDCore3D::constant_u
 
static Magnetostatics::SolutionCurlType HArDCore3D::constant_sigma
 
static Magnetostatics::ForcingTermType HArDCore3D::constant_f
 
static Magnetostatics::PermeabilityType HArDCore3D::constant_mu = Magnetostatics::PermeabilityType(1.)
 
static Magnetostatics::SolutionPotentialType HArDCore3D::linear_u
 
static Magnetostatics::SolutionCurlType HArDCore3D::linear_sigma
 
static Magnetostatics::ForcingTermType HArDCore3D::linear_f
 
static Magnetostatics::PermeabilityType HArDCore3D::linear_mu = Magnetostatics::PermeabilityType(1.)
 
static Magnetostatics::SolutionPotentialType HArDCore3D::trigonometric_u
 
static Magnetostatics::SolutionCurlType HArDCore3D::trigonometric_sigma
 
static Magnetostatics::ForcingTermType HArDCore3D::trigonometric_f
 
static Magnetostatics::PermeabilityType HArDCore3D::trigonometric_mu = Magnetostatics::PermeabilityType(1.)
 
static Magnetostatics::SolutionPotentialType HArDCore3D::variable_permeability_u
 
static Magnetostatics::SolutionCurlType HArDCore3D::variable_permeability_sigma
 
static Magnetostatics::ForcingTermType HArDCore3D::variable_permeability_f
 
static Magnetostatics::PermeabilityType HArDCore3D::variable_permeability_mu
 
static const double HArDCore3D::PI = boost::math::constants::pi<double>()
 
static Magnetostatics::SolutionPotentialType HArDCore3D::constant_u
 
static Magnetostatics::SolutionCurlType HArDCore3D::constant_sigma
 
static Magnetostatics::ForcingTermType HArDCore3D::constant_f
 
static Magnetostatics::PermeabilityType HArDCore3D::constant_mu = Magnetostatics::PermeabilityType(1.)
 
static Magnetostatics::SolutionPotentialType HArDCore3D::linear_u
 
static Magnetostatics::SolutionCurlType HArDCore3D::linear_sigma
 
static Magnetostatics::ForcingTermType HArDCore3D::linear_f
 
static Magnetostatics::PermeabilityType HArDCore3D::linear_mu = Magnetostatics::PermeabilityType(1.)
 
static Magnetostatics::SolutionPotentialType HArDCore3D::trigonometric_u
 
static Magnetostatics::SolutionCurlType HArDCore3D::trigonometric_sigma
 
static Magnetostatics::ForcingTermType HArDCore3D::trigonometric_f
 
static Magnetostatics::PermeabilityType HArDCore3D::trigonometric_mu = Magnetostatics::PermeabilityType(1.)
 
static Magnetostatics::SolutionPotentialType HArDCore3D::variable_permeability_u
 
static Magnetostatics::SolutionCurlType HArDCore3D::variable_permeability_sigma
 
static Magnetostatics::ForcingTermType HArDCore3D::variable_permeability_f
 
static Magnetostatics::PermeabilityType HArDCore3D::variable_permeability_mu
 

Detailed Description

Implementation of the DDR scheme for the magnetostatic problem.

Implementation of the serendipity DDR scheme for the magnetostatic problem.

Typedef Documentation

◆ ForcingTermType [1/2]

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Magnetostatics::ForcingTermType

◆ ForcingTermType [2/2]

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Magnetostatics::ForcingTermType

◆ PermeabilityType [1/2]

◆ PermeabilityType [2/2]

◆ SolutionCurlType [1/2]

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Magnetostatics::SolutionCurlType

◆ SolutionCurlType [2/2]

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Magnetostatics::SolutionCurlType

◆ SolutionPotentialType [1/2]

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Magnetostatics::SolutionPotentialType

◆ SolutionPotentialType [2/2]

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Magnetostatics::SolutionPotentialType

◆ SystemMatrixType [1/2]

◆ SystemMatrixType [2/2]

Function Documentation

◆ assembleLinearSystem()

void Magnetostatics::assembleLinearSystem ( const ForcingTermType f,
const PermeabilityType mu,
const SolutionPotentialType u 
)

Assemble the global system

Parameters
fForcing term
muPermeability
uBoundary condition

◆ computeNorms()

std::vector< double > Magnetostatics::computeNorms ( const std::vector< Eigen::VectorXd > &  list_dofs) const

Compute the discrete Hcurl \times Hdiv norm of a family of vectors representing the dofs.

Parameters
list_dofsThe list of vectors representing the dofs

◆ dimensionSpace()

size_t HArDCore3D::Magnetostatics::dimensionSpace ( ) const
inline

Returns the dimension of the magnetic field + potential space.

◆ Magnetostatics()

Magnetostatics::Magnetostatics ( const DDRCore ddrcore,
bool  use_threads,
std::ostream &  output = std::cout 
)

Constructor.

Parameters
ddrcoreCore for the DDR space sequence
use_threadsTrue for parallel execution, false for sequential execution
outputOutput stream to print status messages

◆ nbSCDOFs()

size_t HArDCore3D::Magnetostatics::nbSCDOFs ( ) const
inline

Returns the number of statically condensed DOFs (here, the cell magnetic field DOFs)

◆ scMatrix()

const SystemMatrixType & HArDCore3D::Magnetostatics::scMatrix ( ) const
inline

Returns the static condensation recovery operator.

◆ scVector()

Eigen::VectorXd & HArDCore3D::Magnetostatics::scVector ( )
inline

Returns the static condensation rhs.

◆ sizeSystem()

size_t HArDCore3D::Magnetostatics::sizeSystem ( ) const
inline

Returns the size of the statically condensed system.

◆ stabilizationParameter() [1/2]

double & HArDCore3D::Magnetostatics::stabilizationParameter ( )
inline

Returns the stabilization parameter.

◆ stabilizationParameter() [2/2]

const double & HArDCore3D::Magnetostatics::stabilizationParameter ( ) const
inline

Returns the stabilization parameter.

◆ sxCurl()

const SXCurl & HArDCore3D::Magnetostatics::sxCurl ( ) const
inline

Returns the space SXCurl.

◆ sxDiv()

const SXDiv & HArDCore3D::Magnetostatics::sxDiv ( ) const
inline

Returns the space SXDiv.

◆ systemMatrix() [1/2]

SystemMatrixType & HArDCore3D::Magnetostatics::systemMatrix ( )
inline

Returns the linear system matrix.

◆ systemMatrix() [2/2]

const SystemMatrixType & HArDCore3D::Magnetostatics::systemMatrix ( ) const
inline

Returns the linear system matrix.

◆ systemVector() [1/2]

Eigen::VectorXd & HArDCore3D::Magnetostatics::systemVector ( )
inline

Returns the linear system right-hand side vector.

◆ systemVector() [2/2]

const Eigen::VectorXd & HArDCore3D::Magnetostatics::systemVector ( ) const
inline

Returns the linear system right-hand side vector.

◆ xCurl()

const XCurl & HArDCore3D::Magnetostatics::xCurl ( ) const
inline

Returns the space XCurl.

◆ xDiv()

const XDiv & HArDCore3D::Magnetostatics::xDiv ( ) const
inline

Returns the space XDiv.

Variable Documentation

◆ constant_f [1/2]

Magnetostatics::ForcingTermType HArDCore3D::constant_f
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d::Zero();
}

◆ constant_f [2/2]

Magnetostatics::ForcingTermType HArDCore3D::constant_f
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d::Zero();
}

◆ constant_mu [1/2]

◆ constant_mu [2/2]

◆ constant_sigma [1/2]

Magnetostatics::SolutionCurlType HArDCore3D::constant_sigma
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d::Zero();
}

◆ constant_sigma [2/2]

Magnetostatics::SolutionCurlType HArDCore3D::constant_sigma
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d::Zero();
}

◆ constant_u [1/2]

Magnetostatics::SolutionPotentialType HArDCore3D::constant_u
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(1., 2., 3.);
}

◆ constant_u [2/2]

Magnetostatics::SolutionPotentialType HArDCore3D::constant_u
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(1., 2., 3.);
}

◆ linear_f [1/2]

Magnetostatics::ForcingTermType HArDCore3D::linear_f
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d::Zero();
}

◆ linear_f [2/2]

Magnetostatics::ForcingTermType HArDCore3D::linear_f
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d::Zero();
}

◆ linear_mu [1/2]

◆ linear_mu [2/2]

◆ linear_sigma [1/2]

Magnetostatics::SolutionCurlType HArDCore3D::linear_sigma
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(-3., 2., 1.);
}

◆ linear_sigma [2/2]

Magnetostatics::SolutionCurlType HArDCore3D::linear_sigma
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(-3., 2., 1.);
}

◆ linear_u [1/2]

Magnetostatics::SolutionPotentialType HArDCore3D::linear_u
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
x(0) + x(1) + x(2),
2.*(x(0)-x(1)+x(2)),
-x(0) - x(1) + x(2)
);
}

◆ linear_u [2/2]

Magnetostatics::SolutionPotentialType HArDCore3D::linear_u
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
x(0) + x(1) + x(2),
2.*(x(0)-x(1)+x(2)),
-x(0) - x(1) + x(2)
);
}

◆ PI [1/2]

const double HArDCore3D::PI = boost::math::constants::pi<double>()
static

◆ PI [2/2]

const double HArDCore3D::PI = boost::math::constants::pi<double>()
static

◆ trigonometric_f [1/2]

Magnetostatics::ForcingTermType HArDCore3D::trigonometric_f
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return 3. * std::pow(PI, 2) * Eigen::Vector3d(
cos(PI*x(0)) * sin(PI*x(1)) * sin(PI*x(2)),
-2. * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2)),
sin(PI*x(0)) * sin(PI*x(1)) * cos(PI*x(2))
);
}

◆ trigonometric_f [2/2]

Magnetostatics::ForcingTermType HArDCore3D::trigonometric_f
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return 3. * std::pow(PI, 2) * Eigen::Vector3d(
cos(PI*x(0)) * sin(PI*x(1)) * sin(PI*x(2)),
-2. * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2)),
sin(PI*x(0)) * sin(PI*x(1)) * cos(PI*x(2))
);
}

◆ trigonometric_mu [1/2]

Magnetostatics::PermeabilityType HArDCore3D::trigonometric_mu = Magnetostatics::PermeabilityType(1.)
static

◆ trigonometric_mu [2/2]

Magnetostatics::PermeabilityType HArDCore3D::trigonometric_mu = Magnetostatics::PermeabilityType(1.)
static

◆ trigonometric_sigma [1/2]

Magnetostatics::SolutionCurlType HArDCore3D::trigonometric_sigma
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return 3. * PI * Eigen::Vector3d(
sin(PI*x(0)) * cos(PI*x(1)) * cos(PI*x(2)),
0.,
-cos(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2))
);
}

◆ trigonometric_sigma [2/2]

Magnetostatics::SolutionCurlType HArDCore3D::trigonometric_sigma
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return 3. * PI * Eigen::Vector3d(
sin(PI*x(0)) * cos(PI*x(1)) * cos(PI*x(2)),
0.,
-cos(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2))
);
}

◆ trigonometric_u [1/2]

Magnetostatics::SolutionPotentialType HArDCore3D::trigonometric_u
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
cos(PI*x(0)) * sin(PI*x(1)) * sin(PI*x(2)),
-2. * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2)),
sin(PI*x(0)) * sin(PI*x(1)) * cos(PI*x(2))
);
}

◆ trigonometric_u [2/2]

Magnetostatics::SolutionPotentialType HArDCore3D::trigonometric_u
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
cos(PI*x(0)) * sin(PI*x(1)) * sin(PI*x(2)),
-2. * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2)),
sin(PI*x(0)) * sin(PI*x(1)) * cos(PI*x(2))
);
}

◆ variable_permeability_f [1/2]

Magnetostatics::ForcingTermType HArDCore3D::variable_permeability_f
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
double a = 3 * PI / std::pow(1. + x(0) + x(1) + x(2), 2);
double b = 3 * std::pow(PI, 2) / (1. + x(0) + x(1) + x(2));
double dy_sigmax = -a * sin(PI*x(0)) * cos(PI*x(1)) * cos(PI*x(2))
- b * sin(PI*x(0)) * sin(PI*x(1)) * cos(PI*x(2));
double dz_sigmax = -a * sin(PI*x(0)) * cos(PI*x(1)) * cos(PI*x(2))
- b * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2));
;
double dx_sigmaz = a * cos(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2))
+ b * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2));
double dy_sigmaz = a * cos(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2))
+ b * cos(PI*x(0)) * sin(PI*x(1)) * sin(PI*x(2));
return Eigen::Vector3d(
dy_sigmaz,
dz_sigmax - dx_sigmaz,
-dy_sigmax
);
}

◆ variable_permeability_f [2/2]

Magnetostatics::ForcingTermType HArDCore3D::variable_permeability_f
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
double a = 3 * PI / std::pow(1. + x(0) + x(1) + x(2), 2);
double b = 3 * std::pow(PI, 2) / (1. + x(0) + x(1) + x(2));
double dy_sigmax = -a * sin(PI*x(0)) * cos(PI*x(1)) * cos(PI*x(2))
- b * sin(PI*x(0)) * sin(PI*x(1)) * cos(PI*x(2));
double dz_sigmax = -a * sin(PI*x(0)) * cos(PI*x(1)) * cos(PI*x(2))
- b * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2));
;
double dx_sigmaz = a * cos(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2))
+ b * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2));
double dy_sigmaz = a * cos(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2))
+ b * cos(PI*x(0)) * sin(PI*x(1)) * sin(PI*x(2));
return Eigen::Vector3d(
dy_sigmaz,
dz_sigmax - dx_sigmaz,
-dy_sigmax
);
}

◆ variable_permeability_mu [1/2]

Magnetostatics::PermeabilityType HArDCore3D::variable_permeability_mu
static
Initial value:
= Magnetostatics::PermeabilityType(
[](const Cell & T, const Eigen::Vector3d & x) -> double { return 1. + x(0) + x(1) + x(2); },
[](const Cell & T) -> size_t { return 1; }
)

◆ variable_permeability_mu [2/2]

Magnetostatics::PermeabilityType HArDCore3D::variable_permeability_mu
static
Initial value:
= Magnetostatics::PermeabilityType(
[](const Cell & T, const Eigen::Vector3d & x) -> double { return 1. + x(0) + x(1) + x(2); },
[](const Cell & T) -> size_t { return 1; }
)

◆ variable_permeability_sigma [1/2]

Magnetostatics::SolutionCurlType HArDCore3D::variable_permeability_sigma
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return 3 * PI / (1. + x(0) + x(1) + x(2))
* Eigen::Vector3d(
sin(PI*x(0)) * cos(PI*x(1)) * cos(PI*x(2)),
0.,
-cos(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2))
);
}

◆ variable_permeability_sigma [2/2]

Magnetostatics::SolutionCurlType HArDCore3D::variable_permeability_sigma
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return 3 * PI / (1. + x(0) + x(1) + x(2))
* Eigen::Vector3d(
sin(PI*x(0)) * cos(PI*x(1)) * cos(PI*x(2)),
0.,
-cos(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2))
);
}

◆ variable_permeability_u [1/2]

Magnetostatics::SolutionPotentialType HArDCore3D::variable_permeability_u
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
cos(PI*x(0)) * sin(PI*x(1)) * sin(PI*x(2)),
-2. * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2)),
sin(PI*x(0)) * sin(PI*x(1)) * cos(PI*x(2))
);
}

◆ variable_permeability_u [2/2]

Magnetostatics::SolutionPotentialType HArDCore3D::variable_permeability_u
static
Initial value:
= [](const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
cos(PI*x(0)) * sin(PI*x(1)) * sin(PI*x(2)),
-2. * sin(PI*x(0)) * cos(PI*x(1)) * sin(PI*x(2)),
sin(PI*x(0)) * sin(PI*x(1)) * cos(PI*x(2))
);
}