v0.9.0
Classes | Public Member Functions | Public Attributes | List of all members
MagneticElement Struct Reference

Implementation of magneto-static problem (basic Implementation)Look for theory and details here: More...

#include <users_modules/basic_finite_elements/magnetostatic/src/MagneticElement.hpp>

Collaboration diagram for MagneticElement:
[legend]

Classes

struct  BlockData
 data structure storing material constants, model parameters, matrices, etc. More...
 
struct  OpCurlCurl
 calculate and assemble CurlCurl matrix

\[ \mathbf{A} = \int_\Omega \mu^{-1} \left( \nabla \times \mathbf{u} \cdot \nabla \times \mathbf{v} \right) \textrm{d}\Omega \]

where

\[ \mathbf{u} = \nabla \times \mathbf{B} \]

where \(\mathbf{B}\) is magnetic flux and \(\mu\) is magnetic permeability. More...

 
struct  OpNaturalBC
 calculate essential boundary conditions

\[ \mathbf{F} = \int_{\partial\Omega} \ \mathbf{u} \cdot \mathbf{j}_i \textrm{d}{\partial\Omega} \]

where \(\mathbf{j}_i\) is current density function. More...

 
struct  OpPostProcessCurl
 calculate and assemble CurlCurl matrix More...
 
struct  OpStab
 calculate and assemble stabilization matrix

\[ \mathbf{A} = \int_\Omega \epsilon \mathbf{u} \cdot \mathbf{v} \textrm{d}\Omega \]

where \(\epsilon\) is regularization parameter. More...

 
struct  TriFE
 define surface element More...
 
struct  VolumeFE
 definition of volume element More...
 

Public Member Functions

 MagneticElement (MoFEM::Interface &m_field)
 
virtual ~MagneticElement ()
 
MoFEMErrorCode getNaturalBc ()
 get natural boundary conditions More...
 
MoFEMErrorCode getEssentialBc ()
 get essential boundary conditions (only homogenous case is considered) More...
 
MoFEMErrorCode createFields ()
 build problem data structures More...
 
MoFEMErrorCode createElements ()
 create finite elementsCreate volume and surface element. Surface element is used to integrate natural boundary conditions. More...
 
MoFEMErrorCode createProblem ()
 create problemProblem is collection of finite elements. With the information on which fields finite elements operates the matrix and left and right hand side vector could be created. More...
 
MoFEMErrorCode destroyProblem ()
 destroy Distributed mesh manager More...
 
MoFEMErrorCode solveProblem (const bool regression_test=false)
 solve problemCreate matrices; integrate over elements; solve linear system of equations More...
 
MoFEMErrorCode postProcessResults ()
 post-process results, i.e. save solution on the mesh More...
 

Public Attributes

MoFEM::InterfacemField
 
BlockData blockData
 

Detailed Description

Implementation of magneto-static problem (basic Implementation)

Look for theory and details here:

[28] <www.hpfem.jku.at/publications/szthesis.pdf> https://pdfs.semanticscholar.org/0574/e5763d9b64bff16f908e9621f23af3b3dc86.pdf

Election file and all other problem related file are here maxwell_element.

Todo:

Extension for mix formulation

Use appropriate pre-conditioner for large problems

Examples
MagneticElement.hpp, and magnetostatic.cpp.

Definition at line 43 of file MagneticElement.hpp.

Constructor & Destructor Documentation

◆ MagneticElement()

MagneticElement::MagneticElement ( MoFEM::Interface m_field)
Examples
MagneticElement.hpp.

Definition at line 71 of file MagneticElement.hpp.

71 : mField(m_field) {}
MoFEM::Interface & mField

◆ ~MagneticElement()

virtual MagneticElement::~MagneticElement ( )
virtual
Examples
MagneticElement.hpp.

Definition at line 72 of file MagneticElement.hpp.

72 {}

Member Function Documentation

◆ createElements()

MoFEMErrorCode MagneticElement::createElements ( )

create finite elementsCreate volume and surface element. Surface element is used to integrate natural boundary conditions.

Returns
error code
Examples
MagneticElement.hpp, and magnetostatic.cpp.

Definition at line 243 of file MagneticElement.hpp.

243  {
245  // //Elements
254  "MESH_NODE_POSITIONS");
263  blockData.feNaturalBCName, "MESH_NODE_POSITIONS");
265  blockData.feName);
268  // build finite elemnts
270  // build adjacencies
273  }
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string &name_row)=0
set field col which finite element use
virtual MoFEMErrorCode add_ents_to_finite_element_by_type(const EntityHandle entities, const EntityType type, const std::string &name, const bool recursive=true)=0
add entities to finite element
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
MoFEM::Interface & mField
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
#define CHKERR
Inline error check.
Definition: definitions.h:596
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string &name_filed)=0
set finite element field data
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:51
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string &name_row)=0
set field row which finite element use
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elementsBuild finite element data structures. Have to be run before problem and adjacenc...

◆ createFields()

MoFEMErrorCode MagneticElement::createFields ( )

build problem data structures

Returns
error code
Examples
MagneticElement.hpp, and magnetostatic.cpp.

Definition at line 167 of file MagneticElement.hpp.

167  {
169 
170  // Set entities bit level. each entity has bit level depending for example
171  // on refinement level. In this case we do not refine mesh or not do
172  // topological changes, simply set refinement level to zero on all entities.
173 
174  CHKERR mField.getInterface<BitRefManager>()->setBitRefLevelByDim(
175  0, 3,BitRefLevel().set(0));
176 
177  // add fields
179  1);
180  CHKERR mField.add_field("MESH_NODE_POSITIONS", H1, AINSWORTH_LEGENDRE_BASE,
181  3);
182  // meshset consisting all entities in mesh
183  EntityHandle root_set = mField.get_moab().get_root_set();
184  // add entities to field
185  CHKERR mField.add_ents_to_field_by_type(root_set, MBTET,
187 
188  // // The higher-order gradients can be gauged by locally skipping the
189  // // corresponding degrees of freedom and basis functions in the
190  // higher-order
191  // // edge-based, face-based and cell-based finite element subspaces.
192  //
193  // Range tris,edges;
194  // CHKERR mField.get_moab().get_entities_by_type(root_set,MBTRI,tris,true);
195  // mField.get_moab().get_entities_by_type(root_set,MBEDGE,edges,true);
196  //
197  // // Set order in volume
198  // Range bc_ents = unite(blockData.naturalBc,blockData.essentialBc);
199  // Range vol_ents = subtract(unite(tris,edges),bc_ents);
200  // CHKERR
201  // mField.set_field_order(vol_ents,blockData.fieldName,blockData.oRder); int
202  // gauged_order = 1; CHKERR
203  // mField.set_field_order(bc_ents,blockData.fieldName,gauged_order);
204 
205  // Set order on tets
207  blockData.oRder);
209  blockData.oRder);
210  CHKERR mField.set_field_order(root_set, MBEDGE, blockData.fieldName,
211  blockData.oRder);
212 
213  // Set geometry approximation ordered
214  CHKERR mField.add_ents_to_field_by_type(root_set, MBTET,
215  "MESH_NODE_POSITIONS");
216  CHKERR mField.set_field_order(root_set, MBTET, "MESH_NODE_POSITIONS", 2);
217  CHKERR mField.set_field_order(root_set, MBTRI, "MESH_NODE_POSITIONS", 2);
218  CHKERR mField.set_field_order(root_set, MBEDGE, "MESH_NODE_POSITIONS", 2);
219  CHKERR mField.set_field_order(root_set, MBVERTEX, "MESH_NODE_POSITIONS", 1);
220 
221  // build field
223 
224  // get HO geometry for 10 node tets
225  // This method takes coordinates form edges mid nodes in 10 node tet and
226  // project values on 2nd order hierarchical basis used to approx. geometry.
227  Projection10NodeCoordsOnField ent_method_material(mField,
228  "MESH_NODE_POSITIONS");
229  CHKERR mField.loop_dofs("MESH_NODE_POSITIONS", ent_method_material);
230 
232  }
virtual moab::Interface & get_moab()=0
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
virtual MoFEMErrorCode add_field(const std::string &name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add field.
virtual MoFEMErrorCode set_field_order(const EntityHandle meshset, const EntityType type, const std::string &name, const ApproximationOrder order, int verb=DEFAULT_VERBOSITY)=0
Set order approximation of the entities in the field.
MoFEMErrorCode getInterface(const MOFEMuuid &uuid, IFACE *&iface) const
Get interface by uuid and return reference to pointer of interface.
int oRder
approximation order
virtual MoFEMErrorCode add_ents_to_field_by_type(const Range &ents, const EntityType type, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:146
virtual MoFEMErrorCode loop_dofs(const Problem *problem_ptr, const std::string &field_name, RowColData rc, DofMethod &method, int lower_rank, int upper_rank, int verb=DEFAULT_VERBOSITY)=0
Make a loop over dofs.
field with continuous tangents
Definition: definitions.h:172
MoFEM::Interface & mField
#define CHKERR
Inline error check.
Definition: definitions.h:596
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:51
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
continuous field
Definition: definitions.h:171

◆ createProblem()

MoFEMErrorCode MagneticElement::createProblem ( )

create problemProblem is collection of finite elements. With the information on which fields finite elements operates the matrix and left and right hand side vector could be created.

Here we use Distributed mesh manager from PETSc as a inteface.

Returns
error code
Examples
MagneticElement.hpp, and magnetostatic.cpp.

Definition at line 287 of file MagneticElement.hpp.

287  {
289  // set up DM
290  CHKERR DMRegister_MoFEM("DMMOFEM");
291  CHKERR DMCreate(PETSC_COMM_WORLD, &blockData.dM);
292  CHKERR DMSetType(blockData.dM, "DMMOFEM");
293  CHKERR DMMoFEMCreateMoFEM(blockData.dM, &mField, "MAGNETIC_PROBLEM",
294  BitRefLevel().set(0));
295  CHKERR DMSetFromOptions(blockData.dM);
296  // add elements to blockData.dM
299  CHKERR DMSetUp(blockData.dM);
300  // create matrices and vectors
301  CHKERR DMCreateGlobalVector(blockData.dM, &blockData.D);
303  }
PetscErrorCode DMMoFEMAddElement(DM dm, const char fe_name[])
add element to dm
Definition: DMMMoFEM.cpp:414
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMMoFEM.cpp:53
MoFEM::Interface & mField
PetscErrorCode DMMoFEMCreateMoFEM(DM dm, MoFEM::Interface *m_field_ptr, const char problem_name[], const MoFEM::BitRefLevel bit_level, const MoFEM::BitRefLevel bit_mask=MoFEM::BitRefLevel().set())
Must be called by user to set MoFEM data structures.
Definition: DMMMoFEM.cpp:109
#define CHKERR
Inline error check.
Definition: definitions.h:596
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:51
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407

◆ destroyProblem()

MoFEMErrorCode MagneticElement::destroyProblem ( )

destroy Distributed mesh manager

Returns
[description]
Examples
MagneticElement.hpp, and magnetostatic.cpp.

Definition at line 310 of file MagneticElement.hpp.

310  {
311 
313  CHKERR DMDestroy(&blockData.dM);
315  }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407

◆ getEssentialBc()

MoFEMErrorCode MagneticElement::getEssentialBc ( )

get essential boundary conditions (only homogenous case is considered)

Returns
error code
Examples
MagneticElement.hpp, and magnetostatic.cpp.

Definition at line 136 of file MagneticElement.hpp.

136  {
139  if (bit->getName().compare(0, 10, "ESSENTIALBC") == 0) {
140  Range faces;
141  CHKERR mField.get_moab().get_entities_by_type(bit->meshset, MBTRI,
142  faces, true);
143  CHKERR mField.get_moab().get_adjacencies(
144  faces, 1, true, blockData.essentialBc, moab::Interface::UNION);
145  blockData.essentialBc.merge(faces);
146  }
147  }
148  if (blockData.essentialBc.empty()) {
149  Range tets;
150  CHKERR mField.get_moab().get_entities_by_type(0, MBTET, tets);
151  Skinner skin(&mField.get_moab());
152  Range skin_faces; // skin faces from 3d ents
153  CHKERR skin.find_skin(0, tets, false, skin_faces);
154  skin_faces = subtract(skin_faces, blockData.naturalBc);
155  CHKERR mField.get_moab().get_adjacencies(
156  skin_faces, 1, true, blockData.essentialBc, moab::Interface::UNION);
157  blockData.essentialBc.merge(skin_faces);
158  }
160  }
virtual moab::Interface & get_moab()=0
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
MoFEM::Interface & mField
#define _IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet having a particular BC meshset in a moFEM field.
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407

◆ getNaturalBc()

MoFEMErrorCode MagneticElement::getNaturalBc ( )

get natural boundary conditions

Returns
error code
Examples
MagneticElement.hpp, and magnetostatic.cpp.

Definition at line 116 of file MagneticElement.hpp.

116  {
119  if (bit->getName().compare(0, 9, "NATURALBC") == 0) {
120  Range faces;
121  CHKERR mField.get_moab().get_entities_by_type(bit->meshset, MBTRI,
122  faces, true);
123  CHKERR mField.get_moab().get_adjacencies(
124  faces, 1, true, blockData.naturalBc, moab::Interface::UNION);
125  blockData.naturalBc.merge(faces);
126  }
127  }
129  }
virtual moab::Interface & get_moab()=0
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
MoFEM::Interface & mField
#define _IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet having a particular BC meshset in a moFEM field.
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407

◆ postProcessResults()

MoFEMErrorCode MagneticElement::postProcessResults ( )

post-process results, i.e. save solution on the mesh

Returns
[description]
Examples
MagneticElement.hpp, and magnetostatic.cpp.

Definition at line 412 of file MagneticElement.hpp.

412  {
413 
416  CHKERR post_proc.generateReferenceElementMesh();
417  CHKERR post_proc.addFieldValuesPostProc("MESH_NODE_POSITIONS");
418  CHKERR post_proc.addFieldValuesPostProc(blockData.fieldName);
419  post_proc.getOpPtrVector().push_back(new OpPostProcessCurl(
420  blockData, post_proc.postProcMesh, post_proc.mapGaussPts));
422  &post_proc);
423  CHKERR post_proc.writeFile("out_values.h5m");
425  }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method)
Executes FEMethod for finite elements in DM.
Definition: DMMMoFEM.cpp:496
Post processing.
MoFEM::Interface & mField
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407

◆ solveProblem()

MoFEMErrorCode MagneticElement::solveProblem ( const bool  regression_test = false)

solve problemCreate matrices; integrate over elements; solve linear system of equations

Examples
MagneticElement.hpp, and magnetostatic.cpp.

Definition at line 323 of file MagneticElement.hpp.

323  {
325  CHKERR DMCreateMatrix(blockData.dM, &blockData.A);
326  CHKERR DMCreateGlobalVector(blockData.dM, &blockData.F);
327  CHKERR VecDuplicate(blockData.F, &blockData.D);
328 
329  VolumeFE vol_fe(mField);
330  vol_fe.getOpPtrVector().push_back(new OpCurlCurl(blockData));
331  vol_fe.getOpPtrVector().push_back(new OpStab(blockData));
332  TriFE tri_fe(mField);
333  tri_fe.getOpPtrVector().push_back(new OpNaturalBC(blockData));
334 
335  CHKERR MatZeroEntries(blockData.A);
336  CHKERR VecZeroEntries(blockData.F);
337  CHKERR VecGhostUpdateBegin(blockData.F, ADD_VALUES, SCATTER_FORWARD);
338  CHKERR VecGhostUpdateEnd(blockData.F, ADD_VALUES, SCATTER_FORWARD);
339 
341  &vol_fe);
343  blockData.feNaturalBCName.c_str(), &tri_fe);
344 
345  CHKERR MatAssemblyBegin(blockData.A, MAT_FINAL_ASSEMBLY);
346  CHKERR MatAssemblyEnd(blockData.A, MAT_FINAL_ASSEMBLY);
347  CHKERR VecGhostUpdateBegin(blockData.F, ADD_VALUES, SCATTER_REVERSE);
348  CHKERR VecGhostUpdateEnd(blockData.F, ADD_VALUES, SCATTER_REVERSE);
349  CHKERR VecAssemblyBegin(blockData.F);
350  CHKERR VecAssemblyEnd(blockData.F);
351 
352  // Boundary conditions
353  std::vector<int> dofs_bc_indices;
354  const MoFEM::Problem *problem_ptr;
355  ierr = DMMoFEMGetProblemPtr(blockData.dM, &problem_ptr);
356  for (Range::iterator eit = blockData.essentialBc.begin();
357  eit != blockData.essentialBc.end(); eit++) {
359  problem_ptr, blockData.fieldName, *eit, mField.get_comm_rank(),
360  dof_ptr)) {
361  std::bitset<8> pstatus(dof_ptr->get()->getPStatus());
362  if (pstatus.test(0))
363  continue; // only local
364  dofs_bc_indices.push_back(dof_ptr->get()->getPetscGlobalDofIdx());
365  }
366  }
367 
368  const double diag = 1;
369  CHKERR MatZeroRowsColumns(
370  blockData.A, dofs_bc_indices.size(),
371  dofs_bc_indices.empty() ? PETSC_NULL : &*dofs_bc_indices.begin(), diag,
372  PETSC_NULL, PETSC_NULL);
373 
374  KSP solver;
375  CHKERR KSPCreate(PETSC_COMM_WORLD, &solver);
376  CHKERR KSPSetOperators(solver, blockData.A, blockData.A);
377  CHKERR KSPSetFromOptions(solver);
378  CHKERR KSPSetUp(solver);
379  CHKERR KSPSolve(solver, blockData.F, blockData.D);
380  CHKERR VecGhostUpdateBegin(blockData.D, INSERT_VALUES, SCATTER_FORWARD);
381  CHKERR VecGhostUpdateEnd(blockData.D, INSERT_VALUES, SCATTER_FORWARD);
382  CHKERR KSPDestroy(&solver);
383  CHKERR VecDestroy(&blockData.F);
384  CHKERR MatDestroy(&blockData.A);
385 
386  CHKERR VecGhostUpdateBegin(blockData.D, INSERT_VALUES, SCATTER_FORWARD);
387  CHKERR VecGhostUpdateEnd(blockData.D, INSERT_VALUES, SCATTER_FORWARD);
389  SCATTER_REVERSE);
390 
391  if (regression_test) {
392  // This test is for order = 1 only
393  double nrm2;
394  CHKERR VecNorm(blockData.D, NORM_2, &nrm2);
395  const double expected_value = 4.6772e+01;
396  const double tol = 1e-4;
397  if ((nrm2 - expected_value) / expected_value > tol) {
398  SETERRQ2(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
399  "Regression test failed %6.4e != %6.4e", nrm2, expected_value);
400  }
401  }
402 
403  CHKERR VecDestroy(&blockData.D);
405  }
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition: DMMMoFEM.cpp:434
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
keeps basic data about problemThis is low level structure with information about problem,...
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method)
Executes FEMethod for finite elements in DM.
Definition: DMMMoFEM.cpp:496
virtual int get_comm_rank() const =0
double tol
#define _IT_NUMEREDDOF_ROW_BY_NAME_ENT_PART_FOR_LOOP_(PROBLEMPTR, NAME, ENT, PART, IT)
use with loops to iterate row DOFs
MoFEM::Interface & mField
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr)
Get pointer to problem data structure.
Definition: DMMMoFEM.cpp:337

Member Data Documentation

◆ blockData

BlockData MagneticElement::blockData
Examples
MagneticElement.hpp, and magnetostatic.cpp.

Definition at line 109 of file MagneticElement.hpp.

◆ mField

MoFEM::Interface& MagneticElement::mField
Examples
MagneticElement.hpp.

Definition at line 45 of file MagneticElement.hpp.


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