Solenoid Stress & Strain¶
This package includes three complementary layers:
- a 1D finite-difference radial stress solver for winding-pack models with zero
rzshear, - a 2D quadrilateral FEM solver with axisymmetric and plane-strain formulations, reusable sparse load operators, cached Rust-side LU solves, and optional BiCGSTAB iterative solves,
- analytic reference formulas used for validation and convergence studies.
1D Finite-Difference Solver¶
cfsem.solenoid_stress.SolenoidStress1D ¶
Bases: NumpyModel
Source code in cfsem/solenoid_stress/solenoid_1d.py
direct_inverse
class-attribute
instance-attribute
¶
Whether to generate fully-dense direct inverse of the system, which can be useful as a linear operator. Alternatively, the system can be solved using an LU solver with reduced memory usage and better numerical conditioning.
displacement_solver
cached
property
¶
LU solver for load-displacement relation (A_ub) as an alternative to taking a direct inverse of A_bu
elasticity_modulus
instance-attribute
¶
[Pa] diagonal terms in material property matrix
operators
cached
property
¶
Linear operators for solving stress and strain in a pancake coil following Iwasa 2e section 3.6.
order
class-attribute
instance-attribute
¶
Finite-difference stencil polynomial order. Higher order operators produce excessive numerical error under typical use.
poisson_ratio
instance-attribute
¶
[dimensionless] factor determining off-diagonal terms in material property matrix
cfsem.solenoid_stress.SolenoidStress1DOperators
dataclass
¶
Linear operators for solving stress and strain in a pancake coil following Iwasa 2e section 3.6.
A_bu, (n x n) sparse operator mapping displacement to the RHS like A @ u_r = -c * j * bz A_ub, (n x n) fully-dense direct inverse of A_bu mapping RHS to displacement A_eu (2n x n), A_eu_radial (n x n), A_eu_hoop (n x n), sparse operators mapping displacement to strain * First entry is combined operator producing both strain components * Second and third entries are split operators, which are equivalent because they are fully decoupled A_se (2n x 2n), sparse operator mapping strain to stress
Source code in cfsem/solenoid_stress/solenoid_1d.py
a_bu
instance-attribute
¶
(n x n) sparse operator mapping displacement to the RHS like A @ u_r = -c * j * bz
a_eu
instance-attribute
¶
(2n x n), sparse operator mapping displacement to strain; contains both radial and hoop components
a_eu_hoop
instance-attribute
¶
(n x n), sparse operators mapping displacement to strain; hoop component only
a_eu_radial
instance-attribute
¶
(n x n), sparse operators mapping displacement to strain; radial component only
a_ub
instance-attribute
¶
(n x n) fully-dense direct inverse of A_bu mapping RHS to displacement.
Only generated if direct_inverse flag is set.
write_mat ¶
Write the collection of operators in .mat format.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
dst
|
str | Path
|
Target directory to place the file named "stress_operators.mat" |
required |
Raises:
| Type | Description |
|---|---|
IOError
|
If the directory does not exist |
Source code in cfsem/solenoid_stress/solenoid_1d.py
cfsem.solenoid_stress.solenoid_1d_structural_factor ¶
Structural factor applied to RHS of solenoid stress solve
cfsem.solenoid_stress.solenoid_1d_structural_rhs ¶
solenoid_1d_structural_rhs(
c: float,
j: NDArray | list[float],
bz: NDArray | list[float],
pi: float = 0.0,
po: float = 0.0,
) -> NDArray
Right-hand-side for solenoid stress solve, including zero values at the BCs.
From Iwasa 2e eqn. 3.64a
Recommend padding the grid with a dummy value at either end to make room for the BCs without losing accounting of nonzero current density at the inner/outer radius.
Padding for BCs can be done like:
rgrid = np.array([r0 - 1e-6] + rgrid.tolist() + [r1 + 1e-6])
Padding region is ultimately treated as structural material, so the padded region should be small to avoid introducing error, but not so small that it causes numerical error in the finite difference scheme.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
c
|
float
|
[m/N] scalar structural factor; see |
required |
j
|
NDArray | list[float]
|
[A/m^2] with shape (n x 1), current density at each point in the r-grid |
required |
bz
|
NDArray | list[float]
|
[T] with shape (n x 1), Z-axis B-field at each point in the r-grid |
required |
pi
|
float
|
[Pa] scalar pressure on inner wall, defined in +r direction |
0.0
|
po
|
float
|
[Pa] scalar pressure on outer wall, defined in -r direction |
0.0
|
Returns:
| Type | Description |
|---|---|
NDArray
|
-c * j * bz, [1/m^2] with shape (n x 1), the right-hand side of the solenoid stress PDE |
Source code in cfsem/solenoid_stress/solenoid_1d.py
2D FEM¶
The FEM path supports:
- axisymmetric and plane-strain structural formulations,
quad4, inferredquad9, and explicitquad9elements,gl3andgl4quadrature,- optional per-element in-plane material orientation angles,
- optional threaded stiffness assembly with
par=True, - reusable reduced-space operators for body force, pressure, traction, and nodal-temperature thermal strain,
- reduced quadrature-point recovery operators for strain and stress,
- direct sparse-LU or BiCGSTAB reduced-system solves,
- model-owned Dirichlet constraints applied during assembly.
The intended workflow is:
- call
assemble_structural_2d(...)once with mesh, materials, load topology, and prescribed Dirichlet values, - build each reduced load vector with
model.build_rhs(...)or the exposed sparse operators, - solve with
model.solve(rhs), using the default cached LU factorization ormodel.solve(rhs, method="bicgstab", ...)for an iterative solve, - recover quadrature strain and stress with
model.evaluate_quadrature(...).
By default, model.solve(rhs) uses the direct sparse-LU path and returns the full displacement
array. Passing method="bicgstab" selects the iterative solver; this path supports diagonal
preconditioning, Ruiz-style row/column equilibration, optional reuse of the previous converged
iterative solution as the initial guess, and optional diagnostics via return_diagnostics=True.
Equilibration scales and the scaled stiffness matrix are cached on the model for repeated
right-hand sides. The BiCGSTAB tolerance is passed unchanged to the system being solved. With
equilibration enabled, this means the tolerance applies to the equilibrated residual rather than
being rescaled to original reduced-RHS units.
Formulation Notes¶
The axisymmetric and plane-strain solvers share the same 2D quadrilateral mesh, two displacement unknowns per node, and four-component strain/stress storage. The difference is how that 2D mesh represents a 3D body.
For formulation="axisymmetric", the coordinates are interpreted as (r, z). Each quadrature
area sample represents a full ring, so stiffness and load integrals use the volume scale
2*pi*r*dA. The strain vector is [rr, zz, tt, rz]; the out-of-plane hoop strain is not an
independent displacement derivative, but is recovered from the radial displacement as
epsilon_tt = u_r / r.
For formulation="plane_strain", the coordinates are interpreted as (x, y). Each area sample
represents a prismatic slice with user-supplied thickness, so integrals use the volume scale
thickness*dA. The strain vector is [xx, yy, zz, xy]; the out-of-plane strain is constrained to
epsilon_zz = 0, while sigma_zz can still be nonzero through the constitutive matrix.
Plane strain and plane stress are different 2D reductions. Plane strain models a body that is long, periodic, or otherwise constrained in the out-of-plane direction, with zero out-of-plane strain and generally nonzero out-of-plane stress. Plane stress models a thin sheet or plate with traction-free faces through the thickness, with zero out-of-plane stress and generally nonzero out-of-plane strain. This FEM path currently implements axisymmetric and plane-strain reductions; it does not implement a plane-stress constitutive reduction.
cfsem.solenoid_stress.fem2d ¶
2D structural elasticity finite-element assembly.
This module provides a small displacement-based quadrilateral FEM solver for axisymmetric and plane-strain structural reductions. The backend stores sparse load operators, sparse quadrature-recovery operators, and a reduced stiffness matrix for repeated load solves.
The element formulation follows the standard small-strain Galerkin construction
K_e = integral(B^T D B c dA) where c is 2*pi*r for axisymmetric and the thickness of the
planar domain for plane strain.
with consistent body-force, surface-pressure, and surface-traction load vectors. The axisymmetric
engineering-strain vector is ordered as [e_rr, e_zz, e_tt, g_rz]. In Bower's terminology, the
underlying equations are the strain-displacement equation, the elastic stress-strain law, the
equation of static equilibrium for stresses, and the boundary conditions on displacement and
stress.
References
[1] Allan F. Bower, Applied Mechanics of Solids, CRC Press, 2009. See especially Section 8.1 and Table 8.3 for the general displacement-based finite-element construction and 2D interpolation functions.
[2] E. L. Wilson, "Structural Analysis of Axisymmetric Solids," AIAA Journal, 3(12), pp. 2269-2274, 1965.
[3] R. A. Mitchell, R. M. Woolley, and C. R. Fisher, "Formulation and experimental verification of an axisymmetric finite-element structural analysis," Journal of Research of the National Bureau of Standards Section C, 75C, 1971.
[4] I. Fried, "Notes on the finite element analysis of the axisymmetric elastic solid," International Journal of Solids and Structures, 10(3), 1974.
ElementMeasures
dataclass
¶
Per-element cross-section area and represented volume.
areas and volumes both have shape (nelem,).
areas has units [area] and volumes has units [volume].
Source code in cfsem/solenoid_stress/fem2d.py
ElementQuadrature
dataclass
¶
Per-element physical quadrature points and mapped weights.
points has shape (nelem, nq_per_element, 2).
weights_area and weights_volume have shape (nelem, nq_per_element).
points has units [length], weights_area has units [area], and
weights_volume has units [volume].
Source code in cfsem/solenoid_stress/fem2d.py
ElevatedQuad9Mesh
dataclass
¶
Explicit 9-node analysis mesh inferred from a corner-only quad4 mesh.
analysis_elements use the local quad9 ordering:
- corners 0..3 in counter-clockwise order [bottom-left, bottom-right, top-right, top-left]
- midsides 4..7 on faces [bottom, right, top, left]
- center node 8
input_nodes and analysis_nodes have units [length].
Source code in cfsem/solenoid_stress/fem2d.py
QuadMeshInterpolation
dataclass
¶
Interpolated nodal values and element-location metadata for query points.
values has shape (npoint, ...), where ... is the trailing shape of the nodal values.
element_indices stores the nearest element used for interpolation. inside reports whether
the nearest-element distance was within the containment tolerance.
Source code in cfsem/solenoid_stress/fem2d.py
QuadMeshQuery
dataclass
¶
One-pass geometric query results for points in a 2D quadrilateral mesh.
The query stores nearest-node, nearest-element, and nearest-face data for each query point.
Interpolation and recovery operators can reuse this object without repeating the mesh search.
For contained points, the nearest element is the containing element and
nearest_element_distances is zero to numerical tolerance.
Source code in cfsem/solenoid_stress/fem2d.py
QuadratureFieldSamples
dataclass
¶
Recovered strain and stress at element quadrature points.
Each field has shape (nelem, nq_per_element, 4) with component ordering
[rr, zz, tt, rz]. points has shape (nelem, nq_per_element, 2).
points has units [length], strain, thermal_strain, and elastic_strain
have units [strain], and stress has units [stress].
Source code in cfsem/solenoid_stress/fem2d.py
total_strain
property
¶
Alias for strain, with shape (nelem, nq_per_element, 4) and units [strain].
Structural2DFEMModel ¶
Reusable 2D structural FEM model with sparse operators and reduced solve state.
The stored sparse operators are the primary reusable objects:
- body_force_to_rhs, pressure_to_rhs, traction_to_rhs, and temperature_to_rhs
map load amplitudes to the reduced structural right-hand side,
- strain_operator, stress_operator, thermal_strain_operator, and
thermal_stress_operator map solved displacements or nodal temperatures to quadrature-point
fields.
build_rhs(...), solve(...), element_quadrature(), element_measures(), and
evaluate_quadrature(...) are convenience methods layered on top of those stored operators
and the reduced stiffness matrix.
Key public array shapes and units:
- stiffness has shape (ndof_reduced, ndof_reduced) with entry units
[generalized force / displacement] = [energy / distance^2],
- body_force_to_rhs has shape (ndof_reduced, 2 * nelem) with entry units [volume],
- pressure_to_rhs has shape (ndof_reduced, n_pressure_faces) with entry units [area],
- traction_to_rhs has shape (ndof_reduced, 2 * n_traction_faces) with entry units
[area],
- temperature_to_rhs has shape (ndof_reduced, n_temperature_nodes) with entry units
[generalized force / temperature] = [energy / (distance * temperature)].
input_nodes and analysis_nodes have shape (nnode, 2) and units [length].
input_elements and analysis_elements expose the original and analysis connectivity.
Source code in cfsem/solenoid_stress/fem2d.py
278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 | |
dtype
property
¶
Floating dtype used by all stored operators, arrays, and convenience-method outputs.
input_elements
property
¶
Input mesh connectivity with shape (nelem, 4).
input_nodes
property
¶
Corner-node input mesh coordinates with shape (nnode, 2) and units [length].
ndof
property
¶
Compatibility alias for ndof_full, the full displacement-vector length (ndof_full,).
build_rhs ¶
build_rhs(
body_force: ArrayLike | None = None,
pressure_values: ArrayLike | None = None,
traction_values: ArrayLike | None = None,
nodal_temperature: ArrayLike | None = None,
) -> npt.NDArray[np.floating[Any]]
Build one reduced structural right-hand side.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
body_force
|
ArrayLike | None
|
Elementwise body-force amplitudes with shape |
None
|
pressure_values
|
ArrayLike | None
|
Pressure amplitudes with shape |
None
|
traction_values
|
ArrayLike | None
|
Surface traction amplitudes with shape |
None
|
nodal_temperature
|
ArrayLike | None
|
Input-node temperatures with shape |
None
|
Returns:
| Name | Type | Description |
|---|---|---|
NDArray |
NDArray[floating[Any]]
|
Reduced right-hand side with shape |
NDArray[floating[Any]]
|
|
Raises:
| Type | Description |
|---|---|
ValueError
|
If thermal materials are present but |
Source code in cfsem/solenoid_stress/fem2d.py
element_measures ¶
Return cross-section area and represented volume for each element.
Returns:
| Name | Type | Description |
|---|---|---|
ElementMeasures |
ElementMeasures
|
Per-element measures with:
|
Source code in cfsem/solenoid_stress/fem2d.py
element_quadrature ¶
Return physical quadrature points and mapped weights for each element.
Returns:
| Name | Type | Description |
|---|---|---|
ElementQuadrature |
ElementQuadrature
|
Quadrature data with:
|
Source code in cfsem/solenoid_stress/fem2d.py
evaluate_quadrature ¶
evaluate_quadrature(
displacements: ArrayLike,
nodal_temperature: ArrayLike | None = None,
) -> QuadratureFieldSamples
Evaluate quadrature-point strain and stress fields.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
displacements
|
ArrayLike
|
Either the reduced displacement solution with shape
|
required |
nodal_temperature
|
ArrayLike | None
|
Input-node temperatures with shape |
None
|
Returns:
| Name | Type | Description |
|---|---|---|
QuadratureFieldSamples |
QuadratureFieldSamples
|
Recovered quadrature fields where:
|
Raises:
| Type | Description |
|---|---|
ValueError
|
If thermal materials are present but |
Source code in cfsem/solenoid_stress/fem2d.py
recover_full ¶
Reinsert prescribed Dirichlet values into a reduced displacement vector.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
reduced_solution
|
ArrayLike
|
Reduced displacement vector with shape |
required |
Returns:
| Name | Type | Description |
|---|---|---|
NDArray |
NDArray[floating[Any]]
|
Full displacement vector with shape |
NDArray[floating[Any]]
|
|
Source code in cfsem/solenoid_stress/fem2d.py
solve ¶
solve(
rhs: ArrayLike,
*,
method: Literal["direct", "bicgstab"] = "direct",
tolerance: float | None = None,
max_iterations: int | None = None,
preconditioner: Literal[
"none", "diagonal"
] = "diagonal",
equilibration: Literal["none", "ruiz"] = "ruiz",
equilibration_max_iterations: int | None = None,
equilibration_tolerance: float | None = None,
equilibration_norm_floor: float | None = None,
equilibration_norm_ceil: float | None = None,
initial_guess: Literal["zero", "previous"] = "zero",
return_diagnostics: Literal[False] = False,
) -> npt.NDArray[np.floating[Any]]
solve(
rhs: ArrayLike,
*,
method: Literal["direct", "bicgstab"] = "direct",
tolerance: float | None = None,
max_iterations: int | None = None,
preconditioner: Literal[
"none", "diagonal"
] = "diagonal",
equilibration: Literal["none", "ruiz"] = "ruiz",
equilibration_max_iterations: int | None = None,
equilibration_tolerance: float | None = None,
equilibration_norm_floor: float | None = None,
equilibration_norm_ceil: float | None = None,
initial_guess: Literal["zero", "previous"] = "zero",
return_diagnostics: Literal[True],
) -> Structural2DSolveResult
solve(
rhs: ArrayLike,
*,
method: Literal["direct", "bicgstab"] = "direct",
tolerance: float | None = None,
max_iterations: int | None = None,
preconditioner: Literal[
"none", "diagonal"
] = "diagonal",
equilibration: Literal["none", "ruiz"] = "ruiz",
equilibration_max_iterations: int | None = None,
equilibration_tolerance: float | None = None,
equilibration_norm_floor: float | None = None,
equilibration_norm_ceil: float | None = None,
initial_guess: Literal["zero", "previous"] = "zero",
return_diagnostics: bool = False,
) -> (
npt.NDArray[np.floating[Any]] | Structural2DSolveResult
)
Solve the reduced system and recover the full displacement field.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
rhs
|
ArrayLike
|
Reduced right-hand side with shape |
required |
method
|
Literal['direct', 'bicgstab']
|
Linear solve method, either |
'direct'
|
tolerance
|
float | None
|
Optional absolute BiCGSTAB residual tolerance for the solved system. If omitted, Rust chooses an RHS-scaled default. When Ruiz equilibration is enabled, this same value is applied to the equilibrated system rather than rescaled to original reduced-RHS units. |
None
|
max_iterations
|
int | None
|
Optional maximum number of BiCGSTAB iterations. |
None
|
preconditioner
|
Literal['none', 'diagonal']
|
BiCGSTAB preconditioner, either |
'diagonal'
|
equilibration
|
Literal['none', 'ruiz']
|
BiCGSTAB equilibration mode, either |
'ruiz'
|
equilibration_max_iterations
|
int | None
|
Optional maximum number of equilibration iterations. |
None
|
equilibration_tolerance
|
float | None
|
Optional dimensionless equilibration tolerance. |
None
|
equilibration_norm_floor
|
float | None
|
Optional lower clamp for measured equilibration norms. |
None
|
equilibration_norm_ceil
|
float | None
|
Optional upper clamp for measured equilibration norms. |
None
|
initial_guess
|
Literal['zero', 'previous']
|
BiCGSTAB initial guess, either |
'zero'
|
return_diagnostics
|
bool
|
If true, return |
False
|
Returns:
| Type | Description |
|---|---|
NDArray[floating[Any]] | Structural2DSolveResult
|
NDArray | Structural2DSolveResult: Full displacement vector with shape |
NDArray[floating[Any]] | Structural2DSolveResult
|
|
NDArray[floating[Any]] | Structural2DSolveResult
|
displacement packaged with solver diagnostics. Units are |
Source code in cfsem/solenoid_stress/fem2d.py
571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 | |
Structural2DSolveDiagnostics
dataclass
¶
Diagnostics from one 2D FEM linear solve.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
method
|
Literal['direct', 'bicgstab']
|
Solver method used for the reduced structural system. |
required |
converged
|
bool
|
Whether the selected solver converged. |
required |
iterations
|
int
|
Number of BiCGSTAB iterations, or zero for direct solves. |
required |
residual_norm
|
float
|
Final residual norm in reduced-RHS units, or zero for direct solves. |
required |
equilibrated
|
bool
|
Whether row/column equilibration was applied before the solve. |
required |
Source code in cfsem/solenoid_stress/fem2d.py
Structural2DSolveResult
dataclass
¶
Full displacement solution plus linear-solver diagnostics.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
displacement
|
NDArray[floating[Any]]
|
Full displacement vector with shape |
required |
diagnostics
|
Structural2DSolveDiagnostics
|
Solver diagnostics for the reduced structural solve. |
required |
Source code in cfsem/solenoid_stress/fem2d.py
assemble_structural_2d ¶
assemble_structural_2d(
nodes: ArrayLike,
elements: ArrayLike,
material_ids: ArrayLike,
material_table: ArrayLike,
*,
formulation: str = "axisymmetric",
thickness: float | None = None,
material_orientation_angles: ArrayLike | None = None,
pressure_faces: ArrayLike | None = None,
traction_faces: ArrayLike | None = None,
thermal_material_table: ArrayLike | None = None,
prescribed: Mapping[int, float] | None = None,
quadrature: str | int = "gl3",
element_type: str = "quad4",
par: bool = True,
) -> Structural2DFEMModel
Assemble the reusable 2D structural FEM model.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
nodes
|
ArrayLike
|
Corner-node coordinates with shape |
required |
elements
|
ArrayLike
|
Connectivity with shape |
required |
material_ids
|
ArrayLike
|
Dense material row indices with shape |
required |
material_table
|
ArrayLike
|
Elastic stress-strain matrices with shape |
required |
formulation
|
str
|
Symmetry reduction, either |
'axisymmetric'
|
thickness
|
float | None
|
Plane-strain out-of-plane thickness. Required only for
|
None
|
material_orientation_angles
|
ArrayLike | None
|
Optional scalar or per-element angles, in radians, rotating local material axes into the global 2D frame before assembly. |
None
|
pressure_faces
|
ArrayLike | None
|
Optional pressure-load topology with shape |
None
|
traction_faces
|
ArrayLike | None
|
Optional traction-load topology with shape |
None
|
thermal_material_table
|
ArrayLike | None
|
Optional thermal material rows with shape |
None
|
prescribed
|
Mapping[int, float] | None
|
Optional mapping from full displacement DOF index to prescribed displacement
value. Displacement units are |
None
|
quadrature
|
str | int
|
Quadrature rule selector, either |
'gl3'
|
element_type
|
str
|
Analysis element family, either |
'quad4'
|
par
|
bool
|
Whether to assemble the stiffness matrix using threaded element batches. |
True
|
Returns:
| Name | Type | Description |
|---|---|---|
Structural2DFEMModel |
Structural2DFEMModel
|
Reusable model storing the reduced stiffness matrix, sparse load |
Structural2DFEMModel
|
operators, sparse quadrature-recovery operators, and cached solve state. |
Raises:
| Type | Description |
|---|---|
ValueError
|
If |
AssertionError
|
If array shapes are invalid or if mapping-style material inputs are passed instead of dense arrays. |
Source code in cfsem/solenoid_stress/fem2d.py
1549 1550 1551 1552 1553 1554 1555 1556 1557 1558 1559 1560 1561 1562 1563 1564 1565 1566 1567 1568 1569 1570 1571 1572 1573 1574 1575 1576 1577 1578 1579 1580 1581 1582 1583 1584 1585 1586 1587 1588 1589 1590 1591 1592 1593 1594 1595 1596 1597 1598 1599 1600 1601 1602 1603 1604 1605 1606 1607 1608 1609 1610 1611 1612 1613 1614 1615 1616 1617 1618 1619 1620 1621 1622 1623 1624 1625 1626 1627 1628 1629 1630 1631 1632 1633 1634 1635 1636 1637 1638 1639 1640 1641 1642 1643 1644 1645 1646 1647 1648 1649 1650 1651 1652 1653 1654 1655 1656 1657 1658 1659 1660 1661 1662 1663 1664 1665 1666 1667 1668 1669 1670 1671 1672 1673 1674 1675 1676 1677 1678 1679 1680 1681 1682 1683 1684 1685 1686 1687 1688 1689 1690 1691 1692 1693 1694 1695 1696 1697 1698 1699 1700 1701 1702 1703 1704 1705 1706 1707 1708 1709 1710 1711 1712 1713 1714 1715 1716 1717 1718 1719 1720 1721 1722 1723 1724 1725 1726 1727 1728 1729 1730 1731 1732 1733 1734 1735 | |
cfsem_radial_material ¶
cfsem_radial_material(
youngs_modulus: float,
poisson_ratio: float,
dtype: DTypeLike = np.float64,
) -> npt.NDArray[np.floating[Any]]
Construct the reduced elastic matrix used by the 1D radial solver.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
youngs_modulus
|
float
|
Young's modulus with units |
required |
poisson_ratio
|
float
|
Poisson ratio with units |
required |
dtype
|
DTypeLike
|
Output floating dtype. |
float64
|
Returns:
| Name | Type | Description |
|---|---|---|
NDArray |
NDArray[floating[Any]]
|
Elastic stress-strain matrix with shape |
NDArray[floating[Any]]
|
|
Source code in cfsem/solenoid_stress/fem2d.py
infer_quad9_mesh ¶
Elevate a corner-only quad mesh to an explicit 9-node Lagrange mesh.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
nodes
|
ArrayLike
|
Corner-node coordinates with shape |
required |
elements
|
ArrayLike
|
Quad4 connectivity with shape |
required |
Returns:
| Name | Type | Description |
|---|---|---|
ElevatedQuad9Mesh |
ElevatedQuad9Mesh
|
Elevated analysis mesh with:
|
Source code in cfsem/solenoid_stress/fem2d.py
interpolate_quad_mesh_values ¶
interpolate_quad_mesh_values(
nodes: ArrayLike,
elements: ArrayLike,
nodal_values: ArrayLike,
points: ArrayLike,
*,
element_type: str = "quad4",
outside: str = "raise",
tolerance: float | None = None,
max_iterations: int = 20,
) -> QuadMeshInterpolation
Interpolate nodal values at arbitrary physical points in a quadrilateral mesh.
The interpolation uses the element's actual shape functions. nodal_values may have shape
(nnode,) or (nnode, ...); the returned values have shape (npoint,) or (npoint, ...).
Point location is Rust-backed but brute-force and scans all elements once per query point.
Outside policies are applied from the nearest-element distance: "raise" errors,
"nan" masks outside values, and "nearest" returns the nearest-element interpolation.
Complexity is O(npoint * nelem * max_iterations) for point location plus
O(npoint * nodes_per_element * ncomponent) for interpolation.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
nodes
|
ArrayLike
|
Mesh node coordinates with shape |
required |
elements
|
ArrayLike
|
Quad connectivity with shape |
required |
nodal_values
|
ArrayLike
|
Values at mesh nodes with shape |
required |
points
|
ArrayLike
|
Query point coordinates with shape |
required |
element_type
|
str
|
Element family, either |
'quad4'
|
outside
|
str
|
Outside-mesh policy: |
'raise'
|
tolerance
|
float | None
|
Physical and reference-space tolerance for point containment. |
None
|
max_iterations
|
int
|
Maximum Newton/projection iterations per element. |
20
|
Returns:
| Type | Description |
|---|---|
QuadMeshInterpolation
|
Interpolated values plus element indices, reference coordinates, and inside flags for the |
QuadMeshInterpolation
|
query points. |
QuadMeshInterpolation
|
|
QuadMeshInterpolation
|
is unitless with shape |
Source code in cfsem/solenoid_stress/fem2d.py
1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208 1209 1210 1211 1212 1213 1214 1215 1216 1217 1218 | |
isotropic_axisymmetric_material ¶
isotropic_axisymmetric_material(
youngs_modulus: float,
poisson_ratio: float,
dtype: DTypeLike = np.float64,
) -> npt.NDArray[np.floating[Any]]
Construct the isotropic axisymmetric elastic stress-strain matrix.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
youngs_modulus
|
float
|
Young's modulus with units |
required |
poisson_ratio
|
float
|
Poisson ratio with units |
required |
dtype
|
DTypeLike
|
Output floating dtype. |
float64
|
Returns:
| Name | Type | Description |
|---|---|---|
NDArray |
NDArray[floating[Any]]
|
Elastic stress-strain matrix with shape |
NDArray[floating[Any]]
|
|
Source code in cfsem/solenoid_stress/fem2d.py
isotropic_axisymmetric_thermal_material ¶
isotropic_axisymmetric_thermal_material(
alpha: float,
reference_temperature: float = 0.0,
dtype: DTypeLike = np.float64,
) -> npt.NDArray[np.floating[Any]]
Construct isotropic thermal-expansion data.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
alpha
|
float
|
Isotropic thermal expansion coefficient with units |
required |
reference_temperature
|
float
|
Stress-free reference temperature with units |
0.0
|
dtype
|
DTypeLike
|
Output floating dtype. |
float64
|
Returns:
| Name | Type | Description |
|---|---|---|
NDArray |
NDArray[floating[Any]]
|
Thermal material row with shape |
NDArray[floating[Any]]
|
|
|
NDArray[floating[Any]]
|
|
Source code in cfsem/solenoid_stress/fem2d.py
isotropic_plane_strain_material ¶
isotropic_plane_strain_material(
youngs_modulus: float,
poisson_ratio: float,
dtype: DTypeLike = np.float64,
) -> npt.NDArray[np.floating[Any]]
Construct the isotropic plane-strain elastic stress-strain matrix.
Returns a dense (4, 4) constitutive matrix in [xx, yy, zz, xy] order. The plane-strain
solver sets epsilon_zz = 0, but this matrix still recovers the nonzero sigma_zz implied
by the in-plane strains.
Source code in cfsem/solenoid_stress/fem2d.py
isotropic_plane_strain_thermal_material ¶
isotropic_plane_strain_thermal_material(
alpha: float,
reference_temperature: float = 0.0,
dtype: DTypeLike = np.float64,
) -> npt.NDArray[np.floating[Any]]
Construct isotropic plane-strain thermal-expansion data.
Returns a row [alpha_x, alpha_y, alpha_z, alpha_xy, T_ref] with equal normal expansion
coefficients and zero engineering shear expansion.
Source code in cfsem/solenoid_stress/fem2d.py
orthotropic_axisymmetric_thermal_material ¶
orthotropic_axisymmetric_thermal_material(
alpha_r: float,
alpha_z: float,
alpha_t: float,
reference_temperature: float = 0.0,
dtype: DTypeLike = np.float64,
) -> npt.NDArray[np.floating[Any]]
Construct orthotropic thermal-expansion data.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
alpha_r
|
float
|
Radial thermal expansion coefficient with units |
required |
alpha_z
|
float
|
Axial thermal expansion coefficient with units |
required |
alpha_t
|
float
|
Hoop thermal expansion coefficient with units |
required |
reference_temperature
|
float
|
Stress-free reference temperature with units |
0.0
|
dtype
|
DTypeLike
|
Output floating dtype. |
float64
|
Returns:
| Name | Type | Description |
|---|---|---|
NDArray |
NDArray[floating[Any]]
|
Thermal material row with shape |
NDArray[floating[Any]]
|
|
|
NDArray[floating[Any]]
|
|
Source code in cfsem/solenoid_stress/fem2d.py
orthotropic_plane_strain_thermal_material ¶
orthotropic_plane_strain_thermal_material(
alpha_x: float,
alpha_y: float,
alpha_z: float,
reference_temperature: float = 0.0,
dtype: DTypeLike = np.float64,
) -> npt.NDArray[np.floating[Any]]
Construct orthotropic plane-strain thermal-expansion data.
Returns a row [alpha_x, alpha_y, alpha_z, alpha_xy, T_ref] with zero engineering shear
expansion. Use material_orientation_angles during assembly to rotate local orthotropic axes.
Source code in cfsem/solenoid_stress/fem2d.py
pack_material_tables_from_tags ¶
pack_material_tables_from_tags(
material_ids: ArrayLike,
material_table_by_tag: Mapping[int, ArrayLike],
thermal_material_table_by_tag: Mapping[int, ArrayLike]
| None = None,
dtype: DTypeLike | None = None,
) -> tuple[
npt.NDArray[np.uint64],
npt.NDArray[np.floating[Any]],
npt.NDArray[np.floating[Any]] | None,
]
Pack tagged material definitions into the dense FEM input format.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
material_ids
|
ArrayLike
|
Element material tags with shape |
required |
material_table_by_tag
|
Mapping[int, ArrayLike]
|
Mapping from external material tag to elastic stress-strain matrix
with shape |
required |
thermal_material_table_by_tag
|
Mapping[int, ArrayLike] | None
|
Optional mapping from external material tag to thermal row
with shape |
None
|
dtype
|
DTypeLike | None
|
Optional output floating dtype. Defaults to the resolved dtype of the provided material rows. |
None
|
Returns:
| Name | Type | Description |
|---|---|---|
tuple |
tuple[NDArray[uint64], NDArray[floating[Any]], NDArray[floating[Any]] | None]
|
|
Raises:
| Type | Description |
|---|---|
ValueError
|
If an element tag is missing from |
Source code in cfsem/solenoid_stress/fem2d.py
1319 1320 1321 1322 1323 1324 1325 1326 1327 1328 1329 1330 1331 1332 1333 1334 1335 1336 1337 1338 1339 1340 1341 1342 1343 1344 1345 1346 1347 1348 1349 1350 1351 1352 1353 1354 1355 1356 1357 1358 1359 1360 1361 1362 1363 1364 1365 1366 1367 1368 1369 1370 1371 1372 1373 1374 1375 1376 1377 1378 1379 1380 1381 1382 1383 1384 1385 1386 1387 1388 1389 1390 1391 1392 1393 1394 1395 1396 1397 1398 1399 1400 1401 1402 1403 1404 1405 1406 1407 1408 1409 1410 1411 1412 1413 1414 1415 1416 1417 | |
quad_mesh_interpolation_operator ¶
Return a reusable sparse operator mapping nodal scalar values to query-point values.
The returned matrix has shape (npoint, nnode). Applying it to a dense (nnode,) vector gives
scalar values at the query points; applying it to (nnode, ncomponent) interpolates multiple
nodal fields with the same operator. The operator always uses the query's nearest element.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
query
|
QuadMeshQuery
|
Mesh query data from |
required |
Returns:
| Type | Description |
|---|---|
csr_matrix
|
Sparse interpolation operator with shape |
csr_matrix
|
function values, so output values have the same units as the nodal values supplied during |
csr_matrix
|
matrix multiplication. |
Source code in cfsem/solenoid_stress/fem2d.py
quad_mesh_strain_operator ¶
quad_mesh_strain_operator(
query: QuadMeshQuery,
*,
formulation: str,
thickness: float | None = None,
) -> sp.csr_matrix
Return a sparse operator mapping full nodal displacements to query-point strain.
The returned matrix has shape (4 * npoint, 2 * nnode). Rows are grouped by query point and
use the same four-component strain ordering as the structural FEM formulation.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
query
|
QuadMeshQuery
|
Mesh query data from |
required |
formulation
|
str
|
Symmetry reduction, either |
required |
thickness
|
float | None
|
Plane-strain out-of-plane thickness with units |
None
|
Returns:
| Type | Description |
|---|---|
csr_matrix
|
Sparse strain-recovery operator with shape |
csr_matrix
|
|
csr_matrix
|
units |
Source code in cfsem/solenoid_stress/fem2d.py
quad_mesh_stress_operator ¶
quad_mesh_stress_operator(
query: QuadMeshQuery,
material_ids: ArrayLike,
material_table: ArrayLike,
*,
formulation: str,
thickness: float | None = None,
material_orientation_angles: ArrayLike | None = None,
) -> sp.csr_matrix
Return a sparse operator mapping full nodal displacements to query-point stress.
The returned matrix has shape (4 * npoint, 2 * nnode). Rows are grouped by query point and
use the same four-component stress ordering as the structural FEM formulation. Each query point
uses the material row assigned to its nearest element; for points contained by the mesh, that is
the containing element. Optional material_orientation_angles follow assembly semantics and
rotate local anisotropic material axes into the global 2D frame per element.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
query
|
QuadMeshQuery
|
Mesh query data from |
required |
material_ids
|
ArrayLike
|
Per-element material row indices with shape |
required |
material_table
|
ArrayLike
|
Elastic stress-strain matrices with shape |
required |
formulation
|
str
|
Symmetry reduction, either |
required |
thickness
|
float | None
|
Plane-strain out-of-plane thickness with units |
None
|
material_orientation_angles
|
ArrayLike | None
|
Optional scalar or per-element angles with shape |
None
|
Returns:
| Type | Description |
|---|---|
csr_matrix
|
Sparse stress-recovery operator with shape |
csr_matrix
|
|
csr_matrix
|
|
csr_matrix
|
units |
Source code in cfsem/solenoid_stress/fem2d.py
1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 | |
query_quad_mesh ¶
query_quad_mesh(
nodes: ArrayLike,
elements: ArrayLike,
points: ArrayLike,
*,
element_type: str = "quad4",
max_iterations: int = 20,
) -> QuadMeshQuery
Query nearest node, nearest element, and nearest face in one pass.
The Rust backend scans all nodes once and all elements once per query point. The element scan
computes nearest-element and nearest-face metadata together so downstream interpolation and
recovery operators do not repeat point location. Complexity is
O(npoint * (nnode + nelem * max_iterations)). A point is contained when its
nearest-element distance is zero to the caller's tolerance.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
nodes
|
ArrayLike
|
Mesh node coordinates with shape |
required |
elements
|
ArrayLike
|
Quad connectivity with shape |
required |
points
|
ArrayLike
|
Query point coordinates with shape |
required |
element_type
|
str
|
Element family, either |
'quad4'
|
max_iterations
|
int
|
Maximum Newton/projection iterations per element. Unitless. |
20
|
Returns:
| Type | Description |
|---|---|
QuadMeshQuery
|
Query data with nearest-node, nearest-element, and nearest-face arrays. Coordinate arrays |
QuadMeshQuery
|
have units |
QuadMeshQuery
|
and index arrays are unitless. |
Source code in cfsem/solenoid_stress/fem2d.py
Analytic Reference Formulas¶
cfsem.solenoid_stress.s_long_solenoid ¶
s_long_solenoid(
r: NDArray,
ri: float,
ro: float,
j: float,
bzi: float,
bzo: float,
poisson_ratio: float,
) -> tuple[NDArray, NDArray]
Radial and hoop stress in an infinitely long solenoid under linearly-varying self field. The "infinite length" assumption is equivalent to assuming zero R-Z shear ("deck of cards") and assuming no B-field in the R-direction (no Z-load or stress).
The linearly varying B-field from inside to outside allows slightly extending this to partially account for the finite length of a real solenoid, which produces a region of negative Bz near the outer radius (as opposed to the true infinite solenoid, for which Bz trends to exactly zero at the outer radius).
Iwasa 2e pg 101 eqns 3.77a,b .
Assumes * Infinitely long solenoid (no R-field or Z-stress). * Linear B-field fall-off between inner and outer radius * "Very long" solenoid - allows some negative field at the OD, but always linearly varying * Uniform current density; no bulk regions of non-conducting structure * Radial stress at inner and outer radius is zero (BC due to no support) * Isotropic material * No thermal stress
Can acommodate a uniform or linearly-varying background field, but not general fields.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
r
|
NDArray
|
[m] (n x 1) array of radius points at which to evaluate the stress |
required |
ri
|
float
|
[m] inner radius |
required |
ro
|
float
|
[m] outer radius |
required |
j
|
float
|
[A/m^2] current density |
required |
bzi
|
float
|
[T] axial B-field at inner radius |
required |
bzo
|
float
|
[T] axial B-field at outer radius |
required |
poisson_ratio
|
float
|
[dimensionless] Material property; off-axis stress coupling term |
required |
Returns:
| Type | Description |
|---|---|
tuple[NDArray, NDArray]
|
s_radial, s_hoop - each (n x 1) with units of [Pa] |
Source code in cfsem/solenoid_stress/solenoid_handcalc.py
cfsem.solenoid_stress.s_radial_thick_wall_cylinder ¶
s_radial_thick_wall_cylinder(
r: NDArray,
ri: float,
ro: float,
pin: float,
pout: float,
) -> NDArray
Radial stress at a location in a thick walled cylinder under pressure load with ends "capped", although the capped constraint does not affect the hoop or radial stress compared to an infinite-length constraint.
https://www.engineeringtoolbox.com/stress-thick-walled-tube-d_949.html https://www.suncam.com/miva/downloads/docs/303.pdf
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
r
|
NDArray
|
[m] radius at which to evaluate |
required |
ri
|
float
|
[m] inner radius |
required |
ro
|
float
|
[m] outer radius |
required |
pin
|
float
|
[Pa] inside pressure |
required |
pout
|
float
|
[Pa] outside pressure |
required |
Returns:
| Type | Description |
|---|---|
NDArray
|
[Pa] radial stress |
Source code in cfsem/solenoid_stress/thick_wall_cylinder_handcalc.py
cfsem.solenoid_stress.s_hoop_thick_wall_cylinder ¶
Hoop stress at a location in a thick walled cylinder under pressure load with ends "capped", although the capped constraint does not affect the hoop or radial stress compared to an infinite-length constraint.
https://www.engineeringtoolbox.com/stress-thick-walled-tube-d_949.html https://www.suncam.com/miva/downloads/docs/303.pdf
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
r
|
NDArray
|
[m] radius at which to evaluate |
required |
ri
|
float
|
[m] inner radius |
required |
ro
|
float
|
[m] outer radius |
required |
pin
|
float
|
[Pa] inside pressure |
required |
pout
|
float
|
[Pa] outside pressure |
required |
Returns:
| Type | Description |
|---|---|
NDArray
|
[Pa] hoop stress |