Skip to content

Solenoid Stress & Strain

This package includes three complementary layers:

  • a 1D finite-difference radial stress solver for winding-pack models with zero rz shear,
  • 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
class SolenoidStress1D(NumpyModel):
    model_config = ConfigDict(validate_assignment=True, frozen=True, extra="forbid")

    rgrid: NpNDArray
    """[m] 1D grid of r-coordinates"""
    elasticity_modulus: float
    """[Pa] diagonal terms in material property matrix"""
    poisson_ratio: float
    """[dimensionless] factor determining off-diagonal terms in material property matrix"""
    order: Literal[2, 4] = 4
    """Finite-difference stencil polynomial order.
       Higher order operators produce excessive numerical error under typical use."""
    direct_inverse: bool = False
    """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."""

    @cached_property
    def operators(self) -> SolenoidStress1DOperators:
        """
        Linear operators for solving stress and strain in a pancake coil
        following Iwasa 2e section 3.6.
        """
        return solenoid_1d_structural_operators(
            np.array(self.rgrid), self.elasticity_modulus, self.poisson_ratio, self.order, self.direct_inverse
        )

    @cached_property
    def displacement_solver(self) -> Callable[[NDArray], NDArray]:
        """LU solver for load-displacement relation (A_ub)
        as an alternative to taking a direct inverse of A_bu"""
        return factorized(self.operators.a_bu)

direct_inverse class-attribute instance-attribute

direct_inverse: bool = False

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

displacement_solver: Callable[[NDArray], NDArray]

LU solver for load-displacement relation (A_ub) as an alternative to taking a direct inverse of A_bu

elasticity_modulus instance-attribute

elasticity_modulus: float

[Pa] diagonal terms in material property matrix

operators cached property

operators: SolenoidStress1DOperators

Linear operators for solving stress and strain in a pancake coil following Iwasa 2e section 3.6.

order class-attribute instance-attribute

order: Literal[2, 4] = 4

Finite-difference stencil polynomial order. Higher order operators produce excessive numerical error under typical use.

poisson_ratio instance-attribute

poisson_ratio: float

[dimensionless] factor determining off-diagonal terms in material property matrix

rgrid instance-attribute

rgrid: NpNDArray

[m] 1D grid of r-coordinates

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
@dataclass(frozen=True)
class SolenoidStress1DOperators:
    """
    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
    """

    a_bu: CSC
    """(n x n) sparse operator mapping displacement to the RHS like A @ u_r = -c * j * bz"""
    a_ub: NDArray | None
    """(n x n) fully-dense direct inverse of A_bu mapping RHS to displacement.
        Only generated if `direct_inverse` flag is set."""
    a_eu: CSR
    """(2n x n), sparse operator mapping displacement to strain; contains both radial and hoop components"""
    a_eu_radial: CSR
    """(n x n), sparse operators mapping displacement to strain; radial component only"""
    a_eu_hoop: CSR
    """(n x n), sparse operators mapping displacement to strain; hoop component only"""
    a_se: CSR
    """(2n x 2n), sparse operator mapping strain to stress"""

    def write_mat(self, dst: str | Path) -> str:
        """Write the collection of operators in .mat format.

        Args:
            dst: Target directory to place the file named "stress_operators.mat"

        Raises:
            IOError: If the directory does not exist
        """
        # Check directory
        dst = Path(dst).absolute()
        fpath = dst / "stress_operators.mat"
        getLogger("cfsem").info(f"Saving stress operator data to {fpath}")

        to_save = {
            "A_bu": self.a_bu,
            "A_ub": self.a_ub,
            "A_eu": self.a_eu,
            "A_eu_radial": self.a_eu_radial,
            "A_eu_hoop": self.a_eu_hoop,
            "A_se": self.a_se,
        }

        if self.a_ub is None:  # savemat fails on None value
            to_save.pop("A_ub")

        #    Note this will implicitly convert all CSR matrices to CSC, which is .mat's preferred I/O
        io.savemat(fpath, to_save)

        return f"{fpath}"

a_bu instance-attribute

a_bu: csc_matrix

(n x n) sparse operator mapping displacement to the RHS like A @ u_r = -c * j * bz

a_eu instance-attribute

a_eu: csr_matrix

(2n x n), sparse operator mapping displacement to strain; contains both radial and hoop components

a_eu_hoop instance-attribute

a_eu_hoop: csr_matrix

(n x n), sparse operators mapping displacement to strain; hoop component only

a_eu_radial instance-attribute

a_eu_radial: csr_matrix

(n x n), sparse operators mapping displacement to strain; radial component only

a_se instance-attribute

a_se: csr_matrix

(2n x 2n), sparse operator mapping strain to stress

a_ub instance-attribute

a_ub: NDArray | None

(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_mat(dst: str | Path) -> str

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
def write_mat(self, dst: str | Path) -> str:
    """Write the collection of operators in .mat format.

    Args:
        dst: Target directory to place the file named "stress_operators.mat"

    Raises:
        IOError: If the directory does not exist
    """
    # Check directory
    dst = Path(dst).absolute()
    fpath = dst / "stress_operators.mat"
    getLogger("cfsem").info(f"Saving stress operator data to {fpath}")

    to_save = {
        "A_bu": self.a_bu,
        "A_ub": self.a_ub,
        "A_eu": self.a_eu,
        "A_eu_radial": self.a_eu_radial,
        "A_eu_hoop": self.a_eu_hoop,
        "A_se": self.a_se,
    }

    if self.a_ub is None:  # savemat fails on None value
        to_save.pop("A_ub")

    #    Note this will implicitly convert all CSR matrices to CSC, which is .mat's preferred I/O
    io.savemat(fpath, to_save)

    return f"{fpath}"

cfsem.solenoid_stress.solenoid_1d_structural_factor

solenoid_1d_structural_factor(
    elasticity_modulus: float, poisson_ratio: float
) -> float

Structural factor applied to RHS of solenoid stress solve

Source code in cfsem/solenoid_stress/solenoid_1d.py
def solenoid_1d_structural_factor(elasticity_modulus: float, poisson_ratio: float) -> float:
    """Structural factor applied to RHS of solenoid stress solve"""
    c = (1.0 - poisson_ratio**2) / elasticity_modulus  # [m/N]
    return c

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 solenoid_1d_structural_factor()

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
def 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.

    Args:
        c: [m/N] scalar structural factor; see `solenoid_1d_structural_factor()`
        j: [A/m^2] with shape (n x 1), current density at each point in the r-grid
        bz: [T] with shape (n x 1), Z-axis B-field at each point in the r-grid
        pi: [Pa] scalar pressure on inner wall, defined in +r direction
        po: [Pa] scalar pressure on outer wall, defined in -r direction

    Returns:
        -c * j * bz, [1/m^2] with shape (n x 1), the right-hand side of the solenoid stress PDE
    """
    # Guarantee arrays
    j = np.array(j)
    bz = np.array(bz)

    # RHS without BCs
    rhs = -c * j * bz

    # BC for r-stress at inner and outer radius
    # is the fluid or mechanical pressure, which is usually going to be set to zero
    # to represent an unsupported system, but could be set to a nonzero value
    # to represent a surface load.
    rhs[0] = -pi  # Sign convention: compression is negative stress
    rhs[-1] = -po

    return rhs

2D FEM

The FEM path supports:

  • axisymmetric and plane-strain structural formulations,
  • quad4, inferred quad9, and explicit quad9 elements,
  • gl3 and gl4 quadrature,
  • 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:

  1. call assemble_structural_2d(...) once with mesh, materials, load topology, and prescribed Dirichlet values,
  2. build each reduced load vector with model.build_rhs(...) or the exposed sparse operators,
  3. solve with model.solve(rhs), using the default cached LU factorization or model.solve(rhs, method="bicgstab", ...) for an iterative solve,
  4. 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
@dataclass(frozen=True, slots=True)
class ElementMeasures:
    """Per-element cross-section area and represented volume.

    `areas` and `volumes` both have shape `(nelem,)`.
    `areas` has units `[area]` and `volumes` has units `[volume]`.
    """

    areas: npt.NDArray[np.floating[Any]]
    volumes: npt.NDArray[np.floating[Any]]

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
@dataclass(frozen=True, slots=True)
class ElementQuadrature:
    """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]`.
    """

    points: npt.NDArray[np.floating[Any]]
    weights_area: npt.NDArray[np.floating[Any]]
    weights_volume: npt.NDArray[np.floating[Any]]
    nq_per_element: int

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
@dataclass(frozen=True, slots=True)
class ElevatedQuad9Mesh:
    """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]`.
    """

    input_nodes: npt.NDArray[np.floating[Any]]
    input_elements: npt.NDArray[np.uint64]
    analysis_nodes: npt.NDArray[np.floating[Any]]
    analysis_elements: npt.NDArray[np.uint64]
    corner_node_indices: npt.NDArray[np.int64]
    midside_node_indices: npt.NDArray[np.int64]
    center_node_indices: npt.NDArray[np.int64]

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
@dataclass(frozen=True, slots=True)
class QuadMeshInterpolation:
    """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.
    """

    values: npt.NDArray[np.floating[Any]]
    element_indices: npt.NDArray[np.int64]
    reference_points: npt.NDArray[np.floating[Any]]
    inside: npt.NDArray[np.bool_]

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
@dataclass(frozen=True, slots=True)
class QuadMeshQuery:
    """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.
    """

    nodes: npt.NDArray[np.floating[Any]]
    elements: npt.NDArray[np.uint64]
    points: npt.NDArray[np.floating[Any]]
    element_type: str
    nearest_node_indices: npt.NDArray[np.int64]
    nearest_node_points: npt.NDArray[np.floating[Any]]
    nearest_node_distances: npt.NDArray[np.floating[Any]]
    nearest_element_indices: npt.NDArray[np.int64]
    nearest_element_reference_points: npt.NDArray[np.floating[Any]]
    nearest_element_points: npt.NDArray[np.floating[Any]]
    nearest_element_distances: npt.NDArray[np.floating[Any]]
    nearest_face_element_indices: npt.NDArray[np.int64]
    nearest_face_local_faces: npt.NDArray[np.int64]
    nearest_face_reference_coordinates: npt.NDArray[np.floating[Any]]
    nearest_face_points: npt.NDArray[np.floating[Any]]
    nearest_face_distances: npt.NDArray[np.floating[Any]]

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
@dataclass(frozen=True, slots=True)
class QuadratureFieldSamples:
    """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]`.
    """

    points: npt.NDArray[np.floating[Any]]
    strain: npt.NDArray[np.floating[Any]]
    thermal_strain: npt.NDArray[np.floating[Any]]
    elastic_strain: npt.NDArray[np.floating[Any]]
    stress: npt.NDArray[np.floating[Any]]

    @property
    def total_strain(self) -> npt.NDArray[np.floating[Any]]:
        """Alias for `strain`, with shape `(nelem, nq_per_element, 4)` and units `[strain]`."""

        return self.strain
total_strain property
total_strain: NDArray[floating[Any]]

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
class 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.
    """

    def __init__(
        self,
        *,
        backend: Any,
        dtype: np.dtype[Any],
        input_nodes: npt.NDArray[np.floating[Any]],
        input_elements: npt.NDArray[np.uint64],
        analysis_nodes: npt.NDArray[np.floating[Any]],
        analysis_elements: npt.NDArray[np.uint64],
        elevated: ElevatedQuad9Mesh | None,
        pressure_faces: npt.NDArray[np.uint64],
        traction_faces: npt.NDArray[np.uint64],
        formulation: str,
        thickness: float,
        element_type: str,
        stiffness: sp.csc_matrix,
        body_force_to_rhs: sp.csr_matrix,
        pressure_to_rhs: sp.csr_matrix,
        traction_to_rhs: sp.csr_matrix,
        temperature_to_rhs: sp.csr_matrix,
        constant_rhs: npt.NDArray[np.floating[Any]],
        quadrature_points: npt.NDArray[np.floating[Any]],
        strain_operator: sp.csr_matrix,
        stress_operator: sp.csr_matrix,
        thermal_strain_operator: sp.csr_matrix,
        thermal_stress_operator: sp.csr_matrix,
        strain_constant: npt.NDArray[np.floating[Any]],
        stress_constant: npt.NDArray[np.floating[Any]],
        thermal_strain_constant: npt.NDArray[np.floating[Any]],
        thermal_stress_constant: npt.NDArray[np.floating[Any]],
        free_dofs: npt.NDArray[np.int64],
        fixed_dofs: npt.NDArray[np.int64],
        fixed_values: npt.NDArray[np.floating[Any]],
        ndof_full: int,
        ndof_reduced: int,
        nelem: int,
        nq_per_element: int,
        n_temperature_nodes: int,
    ) -> None:
        self._backend = backend
        self._dtype = dtype
        self._input_nodes = input_nodes
        self._input_elements = input_elements
        self._elevated = elevated
        self.stiffness = stiffness
        self.body_force_to_rhs = body_force_to_rhs
        self.pressure_to_rhs = pressure_to_rhs
        self.traction_to_rhs = traction_to_rhs
        self.temperature_to_rhs = temperature_to_rhs
        self.constant_rhs = constant_rhs
        self.pressure_faces = pressure_faces
        self.traction_faces = traction_faces
        self.analysis_nodes = analysis_nodes
        self.analysis_elements = analysis_elements
        self.quadrature_points = quadrature_points
        self.strain_operator = strain_operator
        self.stress_operator = stress_operator
        self.thermal_strain_operator = thermal_strain_operator
        self.thermal_stress_operator = thermal_stress_operator
        self.strain_constant = strain_constant
        self.stress_constant = stress_constant
        self.thermal_strain_constant = thermal_strain_constant
        self.thermal_stress_constant = thermal_stress_constant
        self.free_dofs = free_dofs
        self.fixed_dofs = fixed_dofs
        self.fixed_values = fixed_values
        self.formulation = formulation
        self.thickness = thickness
        self.element_type = element_type
        if formulation == "axisymmetric":
            self.coordinate_labels = ("r", "z")
            self.displacement_labels = ("u_r", "u_z")
            self.tensor_labels = ("rr", "zz", "tt", "rz")
            self.measure_label = "swept_volume"
        else:
            self.coordinate_labels = ("x", "y")
            self.displacement_labels = ("u_x", "u_y")
            self.tensor_labels = ("xx", "yy", "zz", "xy")
            self.measure_label = "volume"
        self.ndof_full = int(ndof_full)
        self.ndof_reduced = int(ndof_reduced)
        self.nelem = int(nelem)
        self.nq_per_element = int(nq_per_element)
        self.n_temperature_nodes = int(n_temperature_nodes)
        self._element_quadrature_cache: ElementQuadrature | None = None
        self._element_measures_cache: ElementMeasures | None = None

    @property
    def dtype(self) -> np.dtype[Any]:
        """Floating dtype used by all stored operators, arrays, and convenience-method outputs."""

        return self._dtype

    @property
    def ndof(self) -> int:
        """Compatibility alias for `ndof_full`, the full displacement-vector length `(ndof_full,)`."""

        return self.ndof_full

    @property
    def input_nodes(self) -> npt.NDArray[np.floating[Any]]:
        """Corner-node input mesh coordinates with shape `(nnode, 2)` and units `[length]`."""

        return self._input_nodes

    @property
    def input_elements(self) -> npt.NDArray[np.uint64]:
        """Input mesh connectivity with shape `(nelem, 4)`."""

        return self._input_elements

    def element_quadrature(self) -> ElementQuadrature:
        """Return physical quadrature points and mapped weights for each element.

        Returns:
            ElementQuadrature: Quadrature data with:
                `points` of shape `(nelem, nq_per_element, 2)` and units `[length]`,
                `weights_area` of shape `(nelem, nq_per_element)` and units `[area]`,
                `weights_volume` of shape `(nelem, nq_per_element)` and units `[volume]`.
        """

        cache = self._element_quadrature_cache
        if cache is not None:
            return cache
        points_flat, weights_area_flat, weights_volume_flat, nq = self._backend.element_quadrature()
        nelem = self.analysis_elements.shape[0]
        points = np.asarray(points_flat, dtype=self.dtype).reshape(nelem, int(nq), 2)
        weights_area = np.asarray(weights_area_flat, dtype=self.dtype).reshape(nelem, int(nq))
        weights_volume = np.asarray(weights_volume_flat, dtype=self.dtype).reshape(nelem, int(nq))
        cache = ElementQuadrature(
            points=points,
            weights_area=weights_area,
            weights_volume=weights_volume,
            nq_per_element=int(nq),
        )
        self._element_quadrature_cache = cache
        return cache

    def element_measures(self) -> ElementMeasures:
        """Return cross-section area and represented volume for each element.

        Returns:
            ElementMeasures: Per-element measures with:
                `areas` of shape `(nelem,)` and units `[area]`,
                `volumes` of shape `(nelem,)` and units `[volume]`.
        """

        cache = self._element_measures_cache
        if cache is not None:
            return cache
        areas, volumes = self._backend.element_measures()
        cache = ElementMeasures(
            areas=np.asarray(areas, dtype=self.dtype),
            volumes=np.asarray(volumes, dtype=self.dtype),
        )
        self._element_measures_cache = cache
        return cache

    def _normalize_temperature_for_backend(
        self,
        nodal_temperature: ArrayLike | None,
    ) -> npt.NDArray[np.floating[Any]] | None:
        if self.n_temperature_nodes == 0:
            values = (
                np.zeros((0,), dtype=self.dtype)
                if nodal_temperature is None
                else np.asarray(nodal_temperature, dtype=self.dtype).reshape(-1)
            )
            assert values.size == 0, "nodal_temperature was provided, but this model has no thermal operator"
            return None
        if nodal_temperature is None:
            raise ValueError("nodal_temperature is required because this model includes thermal materials")
        return _analysis_temperature_for_element_type(
            nodal_temperature,
            self._input_nodes.shape[0],
            self._elevated,
            self.dtype,
        )

    def build_rhs(
        self,
        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.

        Args:
            body_force: Elementwise body-force amplitudes with shape `(2,)` or `(nelem, 2)`.
                Components are `[b_r, b_z]` with units `[force / volume]`.
            pressure_values: Pressure amplitudes with shape `(n_pressure_faces,)` and units
                `[force / area]`. Positive values act in the inward normal direction.
            traction_values: Surface traction amplitudes with shape `(2,)` or
                `(n_traction_faces, 2)`. Components are `[t_r, t_z]` with units
                `[force / area]`.
            nodal_temperature: Input-node temperatures with shape `(n_input_nodes,)` and units
                `[temperature]`. Required only when the model includes thermal materials.

        Returns:
            NDArray: Reduced right-hand side with shape `(ndof_reduced,)` and units
            `[generalized force] = [energy / distance]`.

        Raises:
            ValueError: If thermal materials are present but `nodal_temperature` is omitted.
        """

        body_force_arr = (
            np.zeros((self.nelem, 2), dtype=self.dtype)
            if body_force is None
            else _normalize_body_force(body_force, self.nelem, self.dtype)
        )
        _, nload = _sparse_shape(self.pressure_to_rhs)
        pressure_arr = _normalize_pressure_values(pressure_values, nload, self.dtype)
        _, ntraction_cols = _sparse_shape(self.traction_to_rhs)
        traction_arr = _normalize_traction_values(
            traction_values,
            ntraction_cols // 2,
            self.dtype,
        )
        temperature_arr = self._normalize_temperature_for_backend(nodal_temperature)
        rhs = self._backend.build_rhs(
            body_force_arr.reshape(-1),
            pressure_arr if pressure_arr.size else None,
            traction_arr.reshape(-1) if traction_arr.size else None,
            temperature_arr,
        )
        return np.asarray(rhs, dtype=self.dtype)

    @overload
    def solve(
        self,
        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]]: ...

    @overload
    def solve(
        self,
        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: ...

    def solve(
        self,
        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.

        Args:
            rhs: Reduced right-hand side with shape `(ndof_reduced,)` and units
                `[generalized force] = [energy / distance]`.
            method: Linear solve method, either `"direct"` for cached sparse LU or `"bicgstab"`
                for an iterative solve.
            tolerance: 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.
            max_iterations: Optional maximum number of BiCGSTAB iterations.
            preconditioner: BiCGSTAB preconditioner, either `"diagonal"` or `"none"`.
            equilibration: BiCGSTAB equilibration mode, either `"ruiz"` or `"none"`.
            equilibration_max_iterations: Optional maximum number of equilibration iterations.
            equilibration_tolerance: Optional dimensionless equilibration tolerance.
            equilibration_norm_floor: Optional lower clamp for measured equilibration norms.
            equilibration_norm_ceil: Optional upper clamp for measured equilibration norms.
            initial_guess: BiCGSTAB initial guess, either `"zero"` or `"previous"`.
            return_diagnostics: If true, return `Structural2DSolveResult` with displacement and
                diagnostics. If false, return the displacement array directly.

        Returns:
            NDArray | Structural2DSolveResult: Full displacement vector with shape
            `(ndof_full,)` and component ordering `[u_r0, u_z0, u_r1, u_z1, ...]`, or the same
            displacement packaged with solver diagnostics. Units are `[length]`.
        """

        rhs_arr = np.ascontiguousarray(np.asarray(rhs, dtype=self.dtype).reshape(-1))
        assert (
            rhs_arr.shape[0] == self.ndof_reduced
        ), f"rhs must have length {self.ndof_reduced}; got {rhs_arr.shape}"
        norm_floor = _positive_float_or_none("equilibration_norm_floor", equilibration_norm_floor)
        norm_ceil = _positive_float_or_none("equilibration_norm_ceil", equilibration_norm_ceil)
        assert (
            norm_floor is None or norm_ceil is None or norm_ceil > norm_floor
        ), "equilibration_norm_ceil must be greater than equilibration_norm_floor"
        (
            displacement_raw,
            method_code,
            converged,
            iterations,
            residual_norm,
            equilibrated,
        ) = self._backend.solve(
            rhs_arr,
            _solve_method_code(method),
            _positive_float_or_none("tolerance", tolerance),
            _positive_int_or_none("max_iterations", max_iterations),
            _bicgstab_preconditioner_code(preconditioner),
            _bicgstab_equilibration_code(equilibration),
            _positive_int_or_none("equilibration_max_iterations", equilibration_max_iterations),
            _nonnegative_float_or_none("equilibration_tolerance", equilibration_tolerance),
            norm_floor,
            norm_ceil,
            _bicgstab_initial_guess_code(initial_guess),
        )
        displacement = np.asarray(displacement_raw, dtype=self.dtype)
        if not return_diagnostics:
            return displacement
        return Structural2DSolveResult(
            displacement=displacement,
            diagnostics=Structural2DSolveDiagnostics(
                method=_solve_method_name(int(method_code)),
                converged=bool(converged),
                iterations=int(iterations),
                residual_norm=float(residual_norm),
                equilibrated=bool(equilibrated),
            ),
        )

    def recover_full(self, reduced_solution: ArrayLike) -> npt.NDArray[np.floating[Any]]:
        """Reinsert prescribed Dirichlet values into a reduced displacement vector.

        Args:
            reduced_solution: Reduced displacement vector with shape `(ndof_reduced,)` and units
                `[length]`.

        Returns:
            NDArray: Full displacement vector with shape `(ndof_full,)` and component ordering
            `[u_r0, u_z0, u_r1, u_z1, ...]`. Units are `[length]`.
        """

        reduced_arr = np.asarray(reduced_solution, dtype=self.dtype).reshape(-1)
        assert (
            reduced_arr.shape[0] == self.ndof_reduced
        ), f"reduced_solution must have length {self.ndof_reduced}; got {reduced_arr.shape}"
        full = np.zeros((self.ndof_full,), dtype=self.dtype)
        full[self.fixed_dofs] = self.fixed_values
        full[self.free_dofs] = reduced_arr
        return full

    def evaluate_quadrature(
        self,
        displacements: ArrayLike,
        nodal_temperature: ArrayLike | None = None,
    ) -> QuadratureFieldSamples:
        """Evaluate quadrature-point strain and stress fields.

        Args:
            displacements: Either the reduced displacement solution with shape
                `(ndof_reduced,)`, or the full analysis displacement field with shape
                `(2 * n_analysis_nodes,)` or `(n_analysis_nodes, 2)`. Displacement units are
                `[length]`.
            nodal_temperature: Input-node temperatures with shape `(n_input_nodes,)` and units
                `[temperature]`. Required only when the model includes thermal materials.

        Returns:
            QuadratureFieldSamples: Recovered quadrature fields where:
                `points` has shape `(nelem, nq_per_element, 2)` and units `[length]`,
                `strain`, `thermal_strain`, and `elastic_strain` have shape
                `(nelem, nq_per_element, 4)` and units `[strain]`,
                `stress` has shape `(nelem, nq_per_element, 4)` and units `[stress]`.

        Raises:
            ValueError: If thermal materials are present but `nodal_temperature` is omitted.
        """

        arr = np.asarray(displacements, dtype=self.dtype)
        if arr.ndim == 1 and arr.shape == (self.ndof_reduced,):
            displacements_full = self.recover_full(arr)
        else:
            displacements_full = _normalize_displacements(
                displacements, self.analysis_nodes.shape[0], self.dtype
            ).reshape(-1)
        temperature_arr = self._normalize_temperature_for_backend(nodal_temperature)
        (
            points_flat,
            strain_flat,
            thermal_strain_flat,
            elastic_strain_flat,
            stress_flat,
            nq,
        ) = self._backend.evaluate_quadrature(displacements_full, temperature_arr)
        nelem = self.analysis_elements.shape[0]
        nq = int(nq)
        return QuadratureFieldSamples(
            points=np.asarray(points_flat, dtype=self.dtype).reshape(nelem, nq, 2),
            strain=np.asarray(strain_flat, dtype=self.dtype).reshape(nelem, nq, 4),
            thermal_strain=np.asarray(thermal_strain_flat, dtype=self.dtype).reshape(nelem, nq, 4),
            elastic_strain=np.asarray(elastic_strain_flat, dtype=self.dtype).reshape(nelem, nq, 4),
            stress=np.asarray(stress_flat, dtype=self.dtype).reshape(nelem, nq, 4),
        )
dtype property
dtype: dtype[Any]

Floating dtype used by all stored operators, arrays, and convenience-method outputs.

input_elements property
input_elements: NDArray[uint64]

Input mesh connectivity with shape (nelem, 4).

input_nodes property
input_nodes: NDArray[floating[Any]]

Corner-node input mesh coordinates with shape (nnode, 2) and units [length].

ndof property
ndof: int

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 (2,) or (nelem, 2). Components are [b_r, b_z] with units [force / volume].

None
pressure_values ArrayLike | None

Pressure amplitudes with shape (n_pressure_faces,) and units [force / area]. Positive values act in the inward normal direction.

None
traction_values ArrayLike | None

Surface traction amplitudes with shape (2,) or (n_traction_faces, 2). Components are [t_r, t_z] with units [force / area].

None
nodal_temperature ArrayLike | None

Input-node temperatures with shape (n_input_nodes,) and units [temperature]. Required only when the model includes thermal materials.

None

Returns:

Name Type Description
NDArray NDArray[floating[Any]]

Reduced right-hand side with shape (ndof_reduced,) and units

NDArray[floating[Any]]

[generalized force] = [energy / distance].

Raises:

Type Description
ValueError

If thermal materials are present but nodal_temperature is omitted.

Source code in cfsem/solenoid_stress/fem2d.py
def build_rhs(
    self,
    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.

    Args:
        body_force: Elementwise body-force amplitudes with shape `(2,)` or `(nelem, 2)`.
            Components are `[b_r, b_z]` with units `[force / volume]`.
        pressure_values: Pressure amplitudes with shape `(n_pressure_faces,)` and units
            `[force / area]`. Positive values act in the inward normal direction.
        traction_values: Surface traction amplitudes with shape `(2,)` or
            `(n_traction_faces, 2)`. Components are `[t_r, t_z]` with units
            `[force / area]`.
        nodal_temperature: Input-node temperatures with shape `(n_input_nodes,)` and units
            `[temperature]`. Required only when the model includes thermal materials.

    Returns:
        NDArray: Reduced right-hand side with shape `(ndof_reduced,)` and units
        `[generalized force] = [energy / distance]`.

    Raises:
        ValueError: If thermal materials are present but `nodal_temperature` is omitted.
    """

    body_force_arr = (
        np.zeros((self.nelem, 2), dtype=self.dtype)
        if body_force is None
        else _normalize_body_force(body_force, self.nelem, self.dtype)
    )
    _, nload = _sparse_shape(self.pressure_to_rhs)
    pressure_arr = _normalize_pressure_values(pressure_values, nload, self.dtype)
    _, ntraction_cols = _sparse_shape(self.traction_to_rhs)
    traction_arr = _normalize_traction_values(
        traction_values,
        ntraction_cols // 2,
        self.dtype,
    )
    temperature_arr = self._normalize_temperature_for_backend(nodal_temperature)
    rhs = self._backend.build_rhs(
        body_force_arr.reshape(-1),
        pressure_arr if pressure_arr.size else None,
        traction_arr.reshape(-1) if traction_arr.size else None,
        temperature_arr,
    )
    return np.asarray(rhs, dtype=self.dtype)
element_measures
element_measures() -> ElementMeasures

Return cross-section area and represented volume for each element.

Returns:

Name Type Description
ElementMeasures ElementMeasures

Per-element measures with: areas of shape (nelem,) and units [area], volumes of shape (nelem,) and units [volume].

Source code in cfsem/solenoid_stress/fem2d.py
def element_measures(self) -> ElementMeasures:
    """Return cross-section area and represented volume for each element.

    Returns:
        ElementMeasures: Per-element measures with:
            `areas` of shape `(nelem,)` and units `[area]`,
            `volumes` of shape `(nelem,)` and units `[volume]`.
    """

    cache = self._element_measures_cache
    if cache is not None:
        return cache
    areas, volumes = self._backend.element_measures()
    cache = ElementMeasures(
        areas=np.asarray(areas, dtype=self.dtype),
        volumes=np.asarray(volumes, dtype=self.dtype),
    )
    self._element_measures_cache = cache
    return cache
element_quadrature
element_quadrature() -> ElementQuadrature

Return physical quadrature points and mapped weights for each element.

Returns:

Name Type Description
ElementQuadrature ElementQuadrature

Quadrature data with: points of shape (nelem, nq_per_element, 2) and units [length], weights_area of shape (nelem, nq_per_element) and units [area], weights_volume of shape (nelem, nq_per_element) and units [volume].

Source code in cfsem/solenoid_stress/fem2d.py
def element_quadrature(self) -> ElementQuadrature:
    """Return physical quadrature points and mapped weights for each element.

    Returns:
        ElementQuadrature: Quadrature data with:
            `points` of shape `(nelem, nq_per_element, 2)` and units `[length]`,
            `weights_area` of shape `(nelem, nq_per_element)` and units `[area]`,
            `weights_volume` of shape `(nelem, nq_per_element)` and units `[volume]`.
    """

    cache = self._element_quadrature_cache
    if cache is not None:
        return cache
    points_flat, weights_area_flat, weights_volume_flat, nq = self._backend.element_quadrature()
    nelem = self.analysis_elements.shape[0]
    points = np.asarray(points_flat, dtype=self.dtype).reshape(nelem, int(nq), 2)
    weights_area = np.asarray(weights_area_flat, dtype=self.dtype).reshape(nelem, int(nq))
    weights_volume = np.asarray(weights_volume_flat, dtype=self.dtype).reshape(nelem, int(nq))
    cache = ElementQuadrature(
        points=points,
        weights_area=weights_area,
        weights_volume=weights_volume,
        nq_per_element=int(nq),
    )
    self._element_quadrature_cache = cache
    return cache
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 (ndof_reduced,), or the full analysis displacement field with shape (2 * n_analysis_nodes,) or (n_analysis_nodes, 2). Displacement units are [length].

required
nodal_temperature ArrayLike | None

Input-node temperatures with shape (n_input_nodes,) and units [temperature]. Required only when the model includes thermal materials.

None

Returns:

Name Type Description
QuadratureFieldSamples QuadratureFieldSamples

Recovered quadrature fields where: points has shape (nelem, nq_per_element, 2) and units [length], strain, thermal_strain, and elastic_strain have shape (nelem, nq_per_element, 4) and units [strain], stress has shape (nelem, nq_per_element, 4) and units [stress].

Raises:

Type Description
ValueError

If thermal materials are present but nodal_temperature is omitted.

Source code in cfsem/solenoid_stress/fem2d.py
def evaluate_quadrature(
    self,
    displacements: ArrayLike,
    nodal_temperature: ArrayLike | None = None,
) -> QuadratureFieldSamples:
    """Evaluate quadrature-point strain and stress fields.

    Args:
        displacements: Either the reduced displacement solution with shape
            `(ndof_reduced,)`, or the full analysis displacement field with shape
            `(2 * n_analysis_nodes,)` or `(n_analysis_nodes, 2)`. Displacement units are
            `[length]`.
        nodal_temperature: Input-node temperatures with shape `(n_input_nodes,)` and units
            `[temperature]`. Required only when the model includes thermal materials.

    Returns:
        QuadratureFieldSamples: Recovered quadrature fields where:
            `points` has shape `(nelem, nq_per_element, 2)` and units `[length]`,
            `strain`, `thermal_strain`, and `elastic_strain` have shape
            `(nelem, nq_per_element, 4)` and units `[strain]`,
            `stress` has shape `(nelem, nq_per_element, 4)` and units `[stress]`.

    Raises:
        ValueError: If thermal materials are present but `nodal_temperature` is omitted.
    """

    arr = np.asarray(displacements, dtype=self.dtype)
    if arr.ndim == 1 and arr.shape == (self.ndof_reduced,):
        displacements_full = self.recover_full(arr)
    else:
        displacements_full = _normalize_displacements(
            displacements, self.analysis_nodes.shape[0], self.dtype
        ).reshape(-1)
    temperature_arr = self._normalize_temperature_for_backend(nodal_temperature)
    (
        points_flat,
        strain_flat,
        thermal_strain_flat,
        elastic_strain_flat,
        stress_flat,
        nq,
    ) = self._backend.evaluate_quadrature(displacements_full, temperature_arr)
    nelem = self.analysis_elements.shape[0]
    nq = int(nq)
    return QuadratureFieldSamples(
        points=np.asarray(points_flat, dtype=self.dtype).reshape(nelem, nq, 2),
        strain=np.asarray(strain_flat, dtype=self.dtype).reshape(nelem, nq, 4),
        thermal_strain=np.asarray(thermal_strain_flat, dtype=self.dtype).reshape(nelem, nq, 4),
        elastic_strain=np.asarray(elastic_strain_flat, dtype=self.dtype).reshape(nelem, nq, 4),
        stress=np.asarray(stress_flat, dtype=self.dtype).reshape(nelem, nq, 4),
    )
recover_full
recover_full(
    reduced_solution: ArrayLike,
) -> npt.NDArray[np.floating[Any]]

Reinsert prescribed Dirichlet values into a reduced displacement vector.

Parameters:

Name Type Description Default
reduced_solution ArrayLike

Reduced displacement vector with shape (ndof_reduced,) and units [length].

required

Returns:

Name Type Description
NDArray NDArray[floating[Any]]

Full displacement vector with shape (ndof_full,) and component ordering

NDArray[floating[Any]]

[u_r0, u_z0, u_r1, u_z1, ...]. Units are [length].

Source code in cfsem/solenoid_stress/fem2d.py
def recover_full(self, reduced_solution: ArrayLike) -> npt.NDArray[np.floating[Any]]:
    """Reinsert prescribed Dirichlet values into a reduced displacement vector.

    Args:
        reduced_solution: Reduced displacement vector with shape `(ndof_reduced,)` and units
            `[length]`.

    Returns:
        NDArray: Full displacement vector with shape `(ndof_full,)` and component ordering
        `[u_r0, u_z0, u_r1, u_z1, ...]`. Units are `[length]`.
    """

    reduced_arr = np.asarray(reduced_solution, dtype=self.dtype).reshape(-1)
    assert (
        reduced_arr.shape[0] == self.ndof_reduced
    ), f"reduced_solution must have length {self.ndof_reduced}; got {reduced_arr.shape}"
    full = np.zeros((self.ndof_full,), dtype=self.dtype)
    full[self.fixed_dofs] = self.fixed_values
    full[self.free_dofs] = reduced_arr
    return full
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 (ndof_reduced,) and units [generalized force] = [energy / distance].

required
method Literal['direct', 'bicgstab']

Linear solve method, either "direct" for cached sparse LU or "bicgstab" for an iterative solve.

'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" or "none".

'diagonal'
equilibration Literal['none', 'ruiz']

BiCGSTAB equilibration mode, either "ruiz" or "none".

'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" or "previous".

'zero'
return_diagnostics bool

If true, return Structural2DSolveResult with displacement and diagnostics. If false, return the displacement array directly.

False

Returns:

Type Description
NDArray[floating[Any]] | Structural2DSolveResult

NDArray | Structural2DSolveResult: Full displacement vector with shape

NDArray[floating[Any]] | Structural2DSolveResult

(ndof_full,) and component ordering [u_r0, u_z0, u_r1, u_z1, ...], or the same

NDArray[floating[Any]] | Structural2DSolveResult

displacement packaged with solver diagnostics. Units are [length].

Source code in cfsem/solenoid_stress/fem2d.py
def solve(
    self,
    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.

    Args:
        rhs: Reduced right-hand side with shape `(ndof_reduced,)` and units
            `[generalized force] = [energy / distance]`.
        method: Linear solve method, either `"direct"` for cached sparse LU or `"bicgstab"`
            for an iterative solve.
        tolerance: 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.
        max_iterations: Optional maximum number of BiCGSTAB iterations.
        preconditioner: BiCGSTAB preconditioner, either `"diagonal"` or `"none"`.
        equilibration: BiCGSTAB equilibration mode, either `"ruiz"` or `"none"`.
        equilibration_max_iterations: Optional maximum number of equilibration iterations.
        equilibration_tolerance: Optional dimensionless equilibration tolerance.
        equilibration_norm_floor: Optional lower clamp for measured equilibration norms.
        equilibration_norm_ceil: Optional upper clamp for measured equilibration norms.
        initial_guess: BiCGSTAB initial guess, either `"zero"` or `"previous"`.
        return_diagnostics: If true, return `Structural2DSolveResult` with displacement and
            diagnostics. If false, return the displacement array directly.

    Returns:
        NDArray | Structural2DSolveResult: Full displacement vector with shape
        `(ndof_full,)` and component ordering `[u_r0, u_z0, u_r1, u_z1, ...]`, or the same
        displacement packaged with solver diagnostics. Units are `[length]`.
    """

    rhs_arr = np.ascontiguousarray(np.asarray(rhs, dtype=self.dtype).reshape(-1))
    assert (
        rhs_arr.shape[0] == self.ndof_reduced
    ), f"rhs must have length {self.ndof_reduced}; got {rhs_arr.shape}"
    norm_floor = _positive_float_or_none("equilibration_norm_floor", equilibration_norm_floor)
    norm_ceil = _positive_float_or_none("equilibration_norm_ceil", equilibration_norm_ceil)
    assert (
        norm_floor is None or norm_ceil is None or norm_ceil > norm_floor
    ), "equilibration_norm_ceil must be greater than equilibration_norm_floor"
    (
        displacement_raw,
        method_code,
        converged,
        iterations,
        residual_norm,
        equilibrated,
    ) = self._backend.solve(
        rhs_arr,
        _solve_method_code(method),
        _positive_float_or_none("tolerance", tolerance),
        _positive_int_or_none("max_iterations", max_iterations),
        _bicgstab_preconditioner_code(preconditioner),
        _bicgstab_equilibration_code(equilibration),
        _positive_int_or_none("equilibration_max_iterations", equilibration_max_iterations),
        _nonnegative_float_or_none("equilibration_tolerance", equilibration_tolerance),
        norm_floor,
        norm_ceil,
        _bicgstab_initial_guess_code(initial_guess),
    )
    displacement = np.asarray(displacement_raw, dtype=self.dtype)
    if not return_diagnostics:
        return displacement
    return Structural2DSolveResult(
        displacement=displacement,
        diagnostics=Structural2DSolveDiagnostics(
            method=_solve_method_name(int(method_code)),
            converged=bool(converged),
            iterations=int(iterations),
            residual_norm=float(residual_norm),
            equilibrated=bool(equilibrated),
        ),
    )

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
@dataclass(frozen=True, slots=True)
class Structural2DSolveDiagnostics:
    """Diagnostics from one 2D FEM linear solve.

    Args:
        method: Solver method used for the reduced structural system.
        converged: Whether the selected solver converged.
        iterations: Number of BiCGSTAB iterations, or zero for direct solves.
        residual_norm: Final residual norm in reduced-RHS units, or zero for direct solves.
        equilibrated: Whether row/column equilibration was applied before the solve.
    """

    method: Literal["direct", "bicgstab"]
    converged: bool
    iterations: int
    residual_norm: float
    equilibrated: bool

Structural2DSolveResult dataclass

Full displacement solution plus linear-solver diagnostics.

Parameters:

Name Type Description Default
displacement NDArray[floating[Any]]

Full displacement vector with shape (ndof_full,).

required
diagnostics Structural2DSolveDiagnostics

Solver diagnostics for the reduced structural solve.

required
Source code in cfsem/solenoid_stress/fem2d.py
@dataclass(frozen=True, slots=True)
class Structural2DSolveResult:
    """Full displacement solution plus linear-solver diagnostics.

    Args:
        displacement: Full displacement vector with shape `(ndof_full,)`.
        diagnostics: Solver diagnostics for the reduced structural solve.
    """

    displacement: npt.NDArray[np.floating[Any]]
    diagnostics: Structural2DSolveDiagnostics

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 (nnode, 2). Coordinates are (r, z) for formulation="axisymmetric" and (x, y) for formulation="plane_strain". Units are [length].

required
elements ArrayLike

Connectivity with shape (nelem, 4) for element_type="quad4". For element_type="quad9", pass either corner-only (nelem, 4) connectivity to infer a straight-sided quad9 mesh, or explicit (nelem, 9) connectivity in local order [corner0, corner1, corner2, corner3, face0_mid, face1_mid, face2_mid, face3_mid, center]. Corner nodes must be ordered counter-clockwise in the 2D analysis plane.

required
material_ids ArrayLike

Dense material row indices with shape (nelem,).

required
material_table ArrayLike

Elastic stress-strain matrices with shape (nmat, 4, 4). Matrix units are [stress / strain] = [pressure].

required
formulation str

Symmetry reduction, either "axisymmetric" or "plane_strain".

'axisymmetric'
thickness float | None

Plane-strain out-of-plane thickness. Required only for formulation="plane_strain".

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 (n_pressure_faces, 2). Each row is [element_index, local_face].

None
traction_faces ArrayLike | None

Optional traction-load topology with shape (n_traction_faces, 2). Each row is [element_index, local_face].

None
thermal_material_table ArrayLike | None

Optional thermal material rows with shape (nmat, 5) storing [alpha_r, alpha_z, alpha_t, alpha_rz, T_ref]. Thermal expansion coefficients have units [strain / temperature] and T_ref has units [temperature].

None
prescribed Mapping[int, float] | None

Optional mapping from full displacement DOF index to prescribed displacement value. Displacement units are [length].

None
quadrature str | int

Quadrature rule selector, either gl3, gl4, 3, or 4.

'gl3'
element_type str

Analysis element family, either quad4 or quad9.

'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 quadrature is unsupported.

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
def 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.

    Args:
        nodes: Corner-node coordinates with shape `(nnode, 2)`. Coordinates are `(r, z)` for
            `formulation="axisymmetric"` and `(x, y)` for `formulation="plane_strain"`.
            Units are `[length]`.
        elements: Connectivity with shape `(nelem, 4)` for `element_type="quad4"`. For
            `element_type="quad9"`, pass either corner-only `(nelem, 4)` connectivity to infer a
            straight-sided quad9 mesh, or explicit `(nelem, 9)` connectivity in local order
            `[corner0, corner1, corner2, corner3, face0_mid, face1_mid, face2_mid, face3_mid,
            center]`. Corner nodes must be ordered counter-clockwise in the 2D analysis plane.
        material_ids: Dense material row indices with shape `(nelem,)`.
        material_table: Elastic stress-strain matrices with shape `(nmat, 4, 4)`. Matrix units
            are `[stress / strain] = [pressure]`.
        formulation: Symmetry reduction, either `"axisymmetric"` or `"plane_strain"`.
        thickness: Plane-strain out-of-plane thickness. Required only for
            `formulation="plane_strain"`.
        material_orientation_angles: Optional scalar or per-element angles, in radians, rotating
            local material axes into the global 2D frame before assembly.
        pressure_faces: Optional pressure-load topology with shape `(n_pressure_faces, 2)`. Each
            row is `[element_index, local_face]`.
        traction_faces: Optional traction-load topology with shape `(n_traction_faces, 2)`. Each
            row is `[element_index, local_face]`.
        thermal_material_table: Optional thermal material rows with shape `(nmat, 5)` storing
            `[alpha_r, alpha_z, alpha_t, alpha_rz, T_ref]`. Thermal expansion coefficients have
            units `[strain / temperature]` and `T_ref` has units `[temperature]`.
        prescribed: Optional mapping from full displacement DOF index to prescribed displacement
            value. Displacement units are `[length]`.
        quadrature: Quadrature rule selector, either `gl3`, `gl4`, `3`, or `4`.
        element_type: Analysis element family, either `quad4` or `quad9`.
        par: Whether to assemble the stiffness matrix using threaded element batches.

    Returns:
        Structural2DFEMModel: Reusable model storing the reduced stiffness matrix, sparse load
        operators, sparse quadrature-recovery operators, and cached solve state.

    Raises:
        ValueError: If `quadrature` is unsupported.
        AssertionError: If array shapes are invalid or if mapping-style material inputs are passed
            instead of dense arrays.
    """

    dtype = _resolve_float_dtype(nodes, material_table, thermal_material_table)
    nodes_arr = _normalize_nodes(nodes, dtype)
    material_ids_arr, material_table_arr = _normalize_materials(material_ids, material_table, dtype)
    thermal_material_table_arr = _normalize_thermal_material_table(
        thermal_material_table,
        dtype,
    )
    pressure_faces_arr = _normalize_face_pairs("pressure_faces", pressure_faces)
    traction_faces_arr = _normalize_face_pairs("traction_faces", traction_faces)
    prescribed_dofs, prescribed_values = _normalize_prescribed_dirichlet(prescribed, dtype)
    quadrature_code = _quadrature_code(quadrature)
    normalized_formulation = _normalize_formulation(formulation)
    thickness_value = _normalize_thickness(normalized_formulation, thickness, dtype)
    normalized_element_type = _normalize_element_type(element_type)
    if normalized_element_type == "quad4":
        elements_arr = _normalize_elements(elements, 4)
        analysis_nodes, analysis_elements, elevated = nodes_arr, elements_arr, None
    else:
        elements_arr = _normalize_elements(elements, (4, 9))
        if elements_arr.shape[1] == 4:
            elevated = infer_quad9_mesh(nodes_arr, elements_arr)
            analysis_nodes, analysis_elements = elevated.analysis_nodes, elevated.analysis_elements
        else:
            analysis_nodes, analysis_elements, elevated = nodes_arr, elements_arr, None
    material_orientation_angles_arr = _normalize_material_orientation_angles(
        material_orientation_angles,
        int(analysis_elements.shape[0]),
        dtype,
    )
    low_level = _dispatch_pair(
        dtype,
        _assemble_model_2d_f32,
        _assemble_model_2d_f64,
    )
    backend = low_level(
        analysis_nodes,
        analysis_elements,
        material_ids_arr,
        material_table_arr,
        pressure_faces_arr,
        traction_faces_arr,
        (
            thermal_material_table_arr
            if thermal_material_table_arr is not None
            else np.zeros((0, 5), dtype=dtype)
        ),
        material_orientation_angles_arr,
        prescribed_dofs,
        prescribed_values,
        4 if normalized_element_type == "quad4" else 9,
        _formulation_code(normalized_formulation),
        thickness_value,
        quadrature_code,
        par,
    )

    stiffness = _csc_matrix_from_binding(backend.stiffness_csc(), dtype)
    body_force_to_rhs = _csr_matrix_from_binding(backend.body_force_to_rhs_csr(), dtype)
    pressure_to_rhs = _csr_matrix_from_binding(backend.pressure_to_rhs_csr(), dtype)
    traction_to_rhs = _csr_matrix_from_binding(backend.traction_to_rhs_csr(), dtype)
    analysis_temperature_to_rhs = _csr_matrix_from_binding(backend.temperature_to_rhs_csr(), dtype)
    strain_operator = _csr_matrix_from_binding(backend.strain_operator_csr(), dtype)
    stress_operator = _csr_matrix_from_binding(backend.stress_operator_csr(), dtype)
    analysis_thermal_strain_operator = _csr_matrix_from_binding(
        backend.thermal_strain_operator_csr(),
        dtype,
    )
    analysis_thermal_stress_operator = _csr_matrix_from_binding(
        backend.thermal_stress_operator_csr(),
        dtype,
    )
    if thermal_material_table_arr is None:
        temperature_to_rhs = sp.csr_matrix((int(backend.ndof_reduced), 0), dtype=dtype)
        thermal_strain_operator = sp.csr_matrix((_sparse_shape(strain_operator)[0], 0), dtype=dtype)
        thermal_stress_operator = sp.csr_matrix((_sparse_shape(stress_operator)[0], 0), dtype=dtype)
        n_temperature_nodes = 0
    else:
        if elevated is None:
            temperature_to_rhs = analysis_temperature_to_rhs
            thermal_strain_operator = analysis_thermal_strain_operator
            thermal_stress_operator = analysis_thermal_stress_operator
        else:
            temperature_elevation = _temperature_elevation_operator(elevated, dtype)
            temperature_to_rhs = _to_csr_matrix(analysis_temperature_to_rhs @ temperature_elevation)
            thermal_strain_operator = _to_csr_matrix(analysis_thermal_strain_operator @ temperature_elevation)
            thermal_stress_operator = _to_csr_matrix(analysis_thermal_stress_operator @ temperature_elevation)
        n_temperature_nodes = nodes_arr.shape[0]

    return Structural2DFEMModel(
        backend=backend,
        dtype=dtype,
        input_nodes=nodes_arr,
        input_elements=elements_arr,
        analysis_nodes=analysis_nodes,
        analysis_elements=analysis_elements,
        elevated=elevated,
        pressure_faces=pressure_faces_arr,
        traction_faces=traction_faces_arr,
        formulation=normalized_formulation,
        thickness=thickness_value,
        element_type=normalized_element_type,
        stiffness=stiffness,
        body_force_to_rhs=body_force_to_rhs,
        pressure_to_rhs=pressure_to_rhs,
        traction_to_rhs=traction_to_rhs,
        temperature_to_rhs=_to_csr_matrix(temperature_to_rhs),
        constant_rhs=np.asarray(backend.constant_rhs(), dtype=dtype),
        quadrature_points=np.asarray(
            backend.quadrature_points_flat(),
            dtype=dtype,
        ).reshape(analysis_elements.shape[0], int(backend.nq_per_element), 2),
        strain_operator=strain_operator,
        stress_operator=stress_operator,
        thermal_strain_operator=_to_csr_matrix(thermal_strain_operator),
        thermal_stress_operator=_to_csr_matrix(thermal_stress_operator),
        strain_constant=np.asarray(backend.strain_constant(), dtype=dtype),
        stress_constant=np.asarray(backend.stress_constant(), dtype=dtype),
        thermal_strain_constant=np.asarray(backend.thermal_strain_constant(), dtype=dtype),
        thermal_stress_constant=np.asarray(backend.thermal_stress_constant(), dtype=dtype),
        free_dofs=np.asarray(backend.free_dofs(), dtype=np.int64),
        fixed_dofs=np.asarray(backend.fixed_dofs(), dtype=np.int64),
        fixed_values=np.asarray(backend.fixed_values(), dtype=dtype),
        ndof_full=int(backend.ndof_full),
        ndof_reduced=int(backend.ndof_reduced),
        nelem=elements_arr.shape[0],
        nq_per_element=int(backend.nq_per_element),
        n_temperature_nodes=n_temperature_nodes,
    )

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 [pressure].

required
poisson_ratio float

Poisson ratio with units [dimensionless].

required
dtype DTypeLike

Output floating dtype.

float64

Returns:

Name Type Description
NDArray NDArray[floating[Any]]

Elastic stress-strain matrix with shape (4, 4) in component order

NDArray[floating[Any]]

[rr, zz, tt, rz]. Units are [stress / strain] = [pressure].

Source code in cfsem/solenoid_stress/fem2d.py
def cfsem_radial_material(
    youngs_modulus: float,
    poisson_ratio: float,
    dtype: npt.DTypeLike = np.float64,
) -> npt.NDArray[np.floating[Any]]:
    """Construct the reduced elastic matrix used by the 1D radial solver.

    Args:
        youngs_modulus: Young's modulus with units `[pressure]`.
        poisson_ratio: Poisson ratio with units `[dimensionless]`.
        dtype: Output floating dtype.

    Returns:
        NDArray: Elastic stress-strain matrix with shape `(4, 4)` in component order
        `[rr, zz, tt, rz]`. Units are `[stress / strain] = [pressure]`.
    """

    resolved_dtype = np.dtype(dtype)
    binding = _dispatch_pair(
        resolved_dtype,
        _cfsem_radial_material_f32,
        _cfsem_radial_material_f64,
    )
    return _as_float_array(binding(youngs_modulus, poisson_ratio), resolved_dtype).reshape(4, 4)

infer_quad9_mesh

infer_quad9_mesh(
    nodes: ArrayLike, elements: ArrayLike
) -> ElevatedQuad9Mesh

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 (nnode, 2) in (r, z) order. Units are [length].

required
elements ArrayLike

Quad4 connectivity with shape (nelem, 4). Corner nodes must be ordered counter-clockwise in the (r, z) plane.

required

Returns:

Name Type Description
ElevatedQuad9Mesh ElevatedQuad9Mesh

Elevated analysis mesh with: analysis_nodes of shape (n_analysis_nodes, 2) and units [length], analysis_elements of shape (nelem, 9), corner_node_indices, midside_node_indices, and center_node_indices as one-dimensional index arrays.

Source code in cfsem/solenoid_stress/fem2d.py
def infer_quad9_mesh(nodes: ArrayLike, elements: ArrayLike) -> ElevatedQuad9Mesh:
    """Elevate a corner-only quad mesh to an explicit 9-node Lagrange mesh.

    Args:
        nodes: Corner-node coordinates with shape `(nnode, 2)` in `(r, z)` order. Units are
            `[length]`.
        elements: Quad4 connectivity with shape `(nelem, 4)`. Corner nodes must be ordered
            counter-clockwise in the `(r, z)` plane.

    Returns:
        ElevatedQuad9Mesh: Elevated analysis mesh with:
            `analysis_nodes` of shape `(n_analysis_nodes, 2)` and units `[length]`,
            `analysis_elements` of shape `(nelem, 9)`,
            `corner_node_indices`, `midside_node_indices`, and `center_node_indices` as
            one-dimensional index arrays.
    """

    dtype = _resolve_float_dtype(nodes)
    nodes_arr = _normalize_nodes(nodes, dtype)
    elements_arr = _normalize_elements(elements)
    binding = _dispatch_pair(dtype, _infer_quad9_mesh_f32, _infer_quad9_mesh_f64)
    (
        analysis_nodes_flat,
        analysis_elements_flat,
        corner_node_indices,
        midside_node_indices,
        center_node_indices,
    ) = binding(nodes_arr, elements_arr)

    return ElevatedQuad9Mesh(
        input_nodes=nodes_arr,
        input_elements=elements_arr,
        analysis_nodes=np.asarray(analysis_nodes_flat, dtype=dtype).reshape(-1, 2),
        analysis_elements=np.asarray(analysis_elements_flat, dtype=np.uint64).reshape(-1, 9),
        corner_node_indices=np.asarray(corner_node_indices, dtype=np.int64),
        midside_node_indices=np.asarray(midside_node_indices, dtype=np.int64),
        center_node_indices=np.asarray(center_node_indices, dtype=np.int64),
    )

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 (nnode, 2).

required
elements ArrayLike

Quad connectivity with shape (nelem, 4) or (nelem, 9).

required
nodal_values ArrayLike

Values at mesh nodes with shape (nnode,) or (nnode, ...).

required
points ArrayLike

Query point coordinates with shape (npoint, 2).

required
element_type str

Element family, either "quad4" or "quad9".

'quad4'
outside str

Outside-mesh policy: "raise"/"error", "nan", or "nearest".

'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. values has shape (npoint,) or (npoint, ...) and the same units as

QuadMeshInterpolation

nodal_values; element_indices is unitless with shape (npoint,); reference_points

QuadMeshInterpolation

is unitless with shape (npoint, 2); inside has shape (npoint,).

Source code in cfsem/solenoid_stress/fem2d.py
def 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.

    Args:
        nodes: Mesh node coordinates with shape `(nnode, 2)`.
        elements: Quad connectivity with shape `(nelem, 4)` or `(nelem, 9)`.
        nodal_values: Values at mesh nodes with shape `(nnode,)` or `(nnode, ...)`.
        points: Query point coordinates with shape `(npoint, 2)`.
        element_type: Element family, either `"quad4"` or `"quad9"`.
        outside: Outside-mesh policy: `"raise"`/`"error"`, `"nan"`, or `"nearest"`.
        tolerance: Physical and reference-space tolerance for point containment.
        max_iterations: Maximum Newton/projection iterations per element.

    Returns:
        Interpolated values plus element indices, reference coordinates, and inside flags for the
        query points. `values` has shape `(npoint,)` or `(npoint, ...)` and the same units as
        `nodal_values`; `element_indices` is unitless with shape `(npoint,)`; `reference_points`
        is unitless with shape `(npoint, 2)`; `inside` has shape `(npoint,)`.
    """

    dtype = _resolve_float_dtype(nodes, nodal_values, points)
    query = query_quad_mesh(
        nodes,
        elements,
        points,
        element_type=element_type,
        max_iterations=max_iterations,
    )
    values_arr = np.asarray(nodal_values, dtype=dtype)
    assert (
        values_arr.ndim >= 1 and values_arr.shape[0] == query.nodes.shape[0]
    ), f"nodal_values must have shape (nnode,) or (nnode, ...); got {values_arr.shape}"
    values_shape = values_arr.shape[1:]
    values_2d = np.ascontiguousarray(values_arr.reshape(query.nodes.shape[0], -1), dtype=dtype)
    outside_policy = str(outside).strip().lower()
    tol = _normalize_query_tolerance(tolerance, dtype)
    inside = query.nearest_element_distances <= tol
    if outside_policy in {"raise", "error"} and not np.all(inside):
        first = int(np.flatnonzero(~inside)[0])
        raise ValueError(f"query point {first} is outside the quad mesh")
    if outside_policy not in {"nearest", "nan", "raise", "error"}:
        raise ValueError(f"unsupported outside policy {outside!r}; use 'raise', 'nan', or 'nearest'")
    operator = quad_mesh_interpolation_operator(query)
    values = np.asarray(operator @ values_2d, dtype=dtype).reshape((query.points.shape[0], *values_shape))
    if outside_policy == "nan" and np.any(~inside):
        values[~inside] = np.nan
    return QuadMeshInterpolation(
        values=values,
        element_indices=query.nearest_element_indices,
        reference_points=query.nearest_element_reference_points,
        inside=inside,
    )

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 [pressure].

required
poisson_ratio float

Poisson ratio with units [dimensionless].

required
dtype DTypeLike

Output floating dtype.

float64

Returns:

Name Type Description
NDArray NDArray[floating[Any]]

Elastic stress-strain matrix with shape (4, 4) in component order

NDArray[floating[Any]]

[rr, zz, tt, rz]. Units are [stress / strain] = [pressure].

Source code in cfsem/solenoid_stress/fem2d.py
def isotropic_axisymmetric_material(
    youngs_modulus: float,
    poisson_ratio: float,
    dtype: npt.DTypeLike = np.float64,
) -> npt.NDArray[np.floating[Any]]:
    """Construct the isotropic axisymmetric elastic stress-strain matrix.

    Args:
        youngs_modulus: Young's modulus with units `[pressure]`.
        poisson_ratio: Poisson ratio with units `[dimensionless]`.
        dtype: Output floating dtype.

    Returns:
        NDArray: Elastic stress-strain matrix with shape `(4, 4)` in component order
        `[rr, zz, tt, rz]`. Units are `[stress / strain] = [pressure]`.
    """

    resolved_dtype = np.dtype(dtype)
    binding = _dispatch_pair(
        resolved_dtype,
        _isotropic_axisymmetric_material_f32,
        _isotropic_axisymmetric_material_f64,
    )
    return _as_float_array(binding(youngs_modulus, poisson_ratio), resolved_dtype).reshape(4, 4)

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 [strain / temperature].

required
reference_temperature float

Stress-free reference temperature with units [temperature].

0.0
dtype DTypeLike

Output floating dtype.

float64

Returns:

Name Type Description
NDArray NDArray[floating[Any]]

Thermal material row with shape (5,) storing

NDArray[floating[Any]]

[alpha_r, alpha_z, alpha_t, alpha_rz, T_ref]. The first four entries have units

NDArray[floating[Any]]

[strain / temperature]; T_ref has units [temperature].

Source code in cfsem/solenoid_stress/fem2d.py
def isotropic_axisymmetric_thermal_material(
    alpha: float,
    reference_temperature: float = 0.0,
    dtype: npt.DTypeLike = np.float64,
) -> npt.NDArray[np.floating[Any]]:
    """Construct isotropic thermal-expansion data.

    Args:
        alpha: Isotropic thermal expansion coefficient with units `[strain / temperature]`.
        reference_temperature: Stress-free reference temperature with units `[temperature]`.
        dtype: Output floating dtype.

    Returns:
        NDArray: Thermal material row with shape `(5,)` storing
        `[alpha_r, alpha_z, alpha_t, alpha_rz, T_ref]`. The first four entries have units
        `[strain / temperature]`; `T_ref` has units `[temperature]`.
    """

    resolved_dtype = np.dtype(dtype)
    binding = _dispatch_pair(
        resolved_dtype,
        _isotropic_axisymmetric_thermal_material_f32,
        _isotropic_axisymmetric_thermal_material_f64,
    )
    return _as_float_array(binding(alpha, reference_temperature), resolved_dtype)

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
def isotropic_plane_strain_material(
    youngs_modulus: float,
    poisson_ratio: float,
    dtype: npt.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.
    """

    resolved_dtype = np.dtype(dtype)
    binding = _dispatch_pair(
        resolved_dtype,
        _isotropic_plane_strain_material_f32,
        _isotropic_plane_strain_material_f64,
    )
    return _as_float_array(binding(youngs_modulus, poisson_ratio), resolved_dtype).reshape(4, 4)

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
def isotropic_plane_strain_thermal_material(
    alpha: float,
    reference_temperature: float = 0.0,
    dtype: npt.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.
    """

    resolved_dtype = np.dtype(dtype)
    binding = _dispatch_pair(
        resolved_dtype,
        _isotropic_plane_strain_thermal_material_f32,
        _isotropic_plane_strain_thermal_material_f64,
    )
    return _as_float_array(binding(alpha, reference_temperature), resolved_dtype)

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 [strain / temperature].

required
alpha_z float

Axial thermal expansion coefficient with units [strain / temperature].

required
alpha_t float

Hoop thermal expansion coefficient with units [strain / temperature].

required
reference_temperature float

Stress-free reference temperature with units [temperature].

0.0
dtype DTypeLike

Output floating dtype.

float64

Returns:

Name Type Description
NDArray NDArray[floating[Any]]

Thermal material row with shape (5,) storing

NDArray[floating[Any]]

[alpha_r, alpha_z, alpha_t, alpha_rz, T_ref]. The first four entries have units

NDArray[floating[Any]]

[strain / temperature]; T_ref has units [temperature].

Source code in cfsem/solenoid_stress/fem2d.py
def orthotropic_axisymmetric_thermal_material(
    alpha_r: float,
    alpha_z: float,
    alpha_t: float,
    reference_temperature: float = 0.0,
    dtype: npt.DTypeLike = np.float64,
) -> npt.NDArray[np.floating[Any]]:
    """Construct orthotropic thermal-expansion data.

    Args:
        alpha_r: Radial thermal expansion coefficient with units `[strain / temperature]`.
        alpha_z: Axial thermal expansion coefficient with units `[strain / temperature]`.
        alpha_t: Hoop thermal expansion coefficient with units `[strain / temperature]`.
        reference_temperature: Stress-free reference temperature with units `[temperature]`.
        dtype: Output floating dtype.

    Returns:
        NDArray: Thermal material row with shape `(5,)` storing
        `[alpha_r, alpha_z, alpha_t, alpha_rz, T_ref]`. The first four entries have units
        `[strain / temperature]`; `T_ref` has units `[temperature]`.
    """

    resolved_dtype = np.dtype(dtype)
    binding = _dispatch_pair(
        resolved_dtype,
        _orthotropic_axisymmetric_thermal_material_f32,
        _orthotropic_axisymmetric_thermal_material_f64,
    )
    return _as_float_array(binding(alpha_r, alpha_z, alpha_t, reference_temperature), resolved_dtype)

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
def orthotropic_plane_strain_thermal_material(
    alpha_x: float,
    alpha_y: float,
    alpha_z: float,
    reference_temperature: float = 0.0,
    dtype: npt.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.
    """

    return orthotropic_axisymmetric_thermal_material(
        alpha_x,
        alpha_y,
        alpha_z,
        reference_temperature,
        dtype,
    )

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 (nelem,).

required
material_table_by_tag Mapping[int, ArrayLike]

Mapping from external material tag to elastic stress-strain matrix with shape (4, 4). Matrix units are [stress / strain] = [pressure].

required
thermal_material_table_by_tag Mapping[int, ArrayLike] | None

Optional mapping from external material tag to thermal row with shape (5,) storing [alpha_r, alpha_z, alpha_t, alpha_rz, T_ref]. Thermal expansion coefficients have units [strain / temperature] and T_ref has units [temperature].

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]

(packed_material_ids, packed_material_table, packed_thermal_material_table) where: packed_material_ids has shape (nelem,), packed_material_table has shape (nmat, 4, 4), packed_thermal_material_table has shape (nmat, 5) when provided, otherwise None.

Raises:

Type Description
ValueError

If an element tag is missing from material_table_by_tag, or if thermal tags do not match the elastic tags exactly.

Source code in cfsem/solenoid_stress/fem2d.py
def 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: npt.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.

    Args:
        material_ids: Element material tags with shape `(nelem,)`.
        material_table_by_tag: Mapping from external material tag to elastic stress-strain matrix
            with shape `(4, 4)`. Matrix units are `[stress / strain] = [pressure]`.
        thermal_material_table_by_tag: Optional mapping from external material tag to thermal row
            with shape `(5,)` storing `[alpha_r, alpha_z, alpha_t, alpha_rz, T_ref]`. Thermal
            expansion coefficients have units `[strain / temperature]` and `T_ref` has units
            `[temperature]`.
        dtype: Optional output floating dtype. Defaults to the resolved dtype of the provided
            material rows.

    Returns:
        tuple: `(packed_material_ids, packed_material_table, packed_thermal_material_table)` where:
            `packed_material_ids` has shape `(nelem,)`,
            `packed_material_table` has shape `(nmat, 4, 4)`,
            `packed_thermal_material_table` has shape `(nmat, 5)` when provided, otherwise `None`.

    Raises:
        ValueError: If an element tag is missing from `material_table_by_tag`, or if thermal tags
            do not match the elastic tags exactly.
    """

    assert material_table_by_tag, "material_table_by_tag cannot be empty"
    resolved_dtype = (
        np.dtype(dtype)
        if dtype is not None
        else _resolve_float_dtype(
            material_table_by_tag,
            thermal_material_table_by_tag,
        )
    )
    ids = np.asarray(material_ids, dtype=np.uint64)
    assert ids.ndim == 1, f"material_ids must have shape (nelem,); got {ids.shape}"

    material_tags = sorted(int(tag) for tag in material_table_by_tag)
    tag_to_index = {tag: index for index, tag in enumerate(material_tags)}
    try:
        packed_ids = np.asarray([tag_to_index[int(tag)] for tag in ids], dtype=np.uint64)
    except KeyError as exc:
        raise ValueError(
            f"material_ids contains tag {exc.args[0]} that is missing from material_table_by_tag"
        ) from exc

    material_rows = []
    for tag in material_tags:
        matrix = cast(
            npt.NDArray[np.floating[Any]],
            np.asarray(material_table_by_tag[tag], dtype=resolved_dtype),
        )
        assert matrix.shape == (
            4,
            4,
        ), f"material_table_by_tag[{tag}] must have shape (4, 4); got {matrix.shape}"
        material_rows.append(matrix)
    packed_material_table = cast(
        npt.NDArray[np.floating[Any]],
        np.ascontiguousarray(np.stack(material_rows, axis=0), dtype=resolved_dtype),
    )

    packed_thermal_table: npt.NDArray[np.floating[Any]] | None
    if thermal_material_table_by_tag is None:
        packed_thermal_table = None
    else:
        thermal_tags = {int(tag) for tag in thermal_material_table_by_tag}
        if thermal_tags != set(material_tags):
            raise ValueError(
                "thermal_material_table_by_tag must have exactly the same keys as material_table_by_tag"
            )
        thermal_rows = []
        for tag in material_tags:
            row = cast(
                npt.NDArray[np.floating[Any]],
                np.asarray(thermal_material_table_by_tag[tag], dtype=resolved_dtype),
            )
            assert row.shape == (
                5,
            ), f"thermal_material_table_by_tag[{tag}] must have shape (5,); got {row.shape}"
            thermal_rows.append(row)
        packed_thermal_table = cast(
            npt.NDArray[np.floating[Any]],
            np.ascontiguousarray(np.stack(thermal_rows, axis=0), dtype=resolved_dtype),
        )
        assert not np.any(
            packed_thermal_table[:, 3] != 0.0
        ), "shear thermal expansion (alpha_rz) is not yet supported"

    return packed_ids, packed_material_table, packed_thermal_table

quad_mesh_interpolation_operator

quad_mesh_interpolation_operator(
    query: QuadMeshQuery,
) -> sp.csr_matrix

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 query_quad_mesh(...).

required

Returns:

Type Description
csr_matrix

Sparse interpolation operator with shape (npoint, nnode). Entries are unitless 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
def quad_mesh_interpolation_operator(
    query: QuadMeshQuery,
) -> sp.csr_matrix:
    """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.

    Args:
        query: Mesh query data from `query_quad_mesh(...)`.

    Returns:
        Sparse interpolation operator with shape `(npoint, nnode)`. Entries are unitless shape
        function values, so output values have the same units as the nodal values supplied during
        matrix multiplication.
    """

    dtype = query.nodes.dtype
    binding = _dispatch_pair(
        dtype,
        _quad_mesh_interpolation_operator_f32,
        _quad_mesh_interpolation_operator_f64,
    )
    return _coo_operator_from_binding(
        binding(
            query.nodes,
            query.elements,
            np.asarray(query.nearest_element_indices, dtype=np.uint64),
            query.nearest_element_reference_points,
            query.element_type,
        ),
        dtype,
    )

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 query_quad_mesh(...).

required
formulation str

Symmetry reduction, either "axisymmetric" or "plane_strain".

required
thickness float | None

Plane-strain out-of-plane thickness with units [length]. Required only for formulation="plane_strain"; must be omitted for formulation="axisymmetric".

None

Returns:

Type Description
csr_matrix

Sparse strain-recovery operator with shape (4 * npoint, 2 * nnode). Entries have units

csr_matrix

[1 / length], so multiplying by full nodal displacements with shape (2 * nnode,) and

csr_matrix

units [length] returns unitless strain samples with shape (4 * npoint,).

Source code in cfsem/solenoid_stress/fem2d.py
def 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.

    Args:
        query: Mesh query data from `query_quad_mesh(...)`.
        formulation: Symmetry reduction, either `"axisymmetric"` or `"plane_strain"`.
        thickness: Plane-strain out-of-plane thickness with units `[length]`. Required only for
            `formulation="plane_strain"`; must be omitted for `formulation="axisymmetric"`.

    Returns:
        Sparse strain-recovery operator with shape `(4 * npoint, 2 * nnode)`. Entries have units
        `[1 / length]`, so multiplying by full nodal displacements with shape `(2 * nnode,)` and
        units `[length]` returns unitless strain samples with shape `(4 * npoint,)`.
    """

    normalized_formulation = _normalize_formulation(formulation)
    thickness_value = _normalize_thickness(normalized_formulation, thickness, query.nodes.dtype)
    dtype = query.nodes.dtype
    binding = _dispatch_pair(dtype, _quad_mesh_strain_operator_f32, _quad_mesh_strain_operator_f64)
    return _coo_operator_from_binding(
        binding(
            query.nodes,
            query.elements,
            np.asarray(query.nearest_element_indices, dtype=np.uint64),
            query.nearest_element_reference_points,
            query.element_type,
            _formulation_code(normalized_formulation),
            thickness_value,
        ),
        dtype,
    )

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 query_quad_mesh(...).

required
material_ids ArrayLike

Per-element material row indices with shape (nelem,). Entries are unitless indices into material_table.

required
material_table ArrayLike

Elastic stress-strain matrices with shape (nmat, 4, 4). Entries have units [stress / strain] = [pressure].

required
formulation str

Symmetry reduction, either "axisymmetric" or "plane_strain".

required
thickness float | None

Plane-strain out-of-plane thickness with units [length]. Required only for formulation="plane_strain"; must be omitted for formulation="axisymmetric".

None
material_orientation_angles ArrayLike | None

Optional scalar or per-element angles with shape (nelem,), in radians. Angles are unitless and rotate local material axes into the global 2D frame.

None

Returns:

Type Description
csr_matrix

Sparse stress-recovery operator with shape (4 * npoint, 2 * nnode). Entries have units

csr_matrix

[pressure / length], so multiplying by full nodal displacements with shape

csr_matrix

(2 * nnode,) and units [length] returns stress samples with shape (4 * npoint,) and

csr_matrix

units [pressure].

Source code in cfsem/solenoid_stress/fem2d.py
def 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.

    Args:
        query: Mesh query data from `query_quad_mesh(...)`.
        material_ids: Per-element material row indices with shape `(nelem,)`. Entries are unitless
            indices into `material_table`.
        material_table: Elastic stress-strain matrices with shape `(nmat, 4, 4)`. Entries have
            units `[stress / strain] = [pressure]`.
        formulation: Symmetry reduction, either `"axisymmetric"` or `"plane_strain"`.
        thickness: Plane-strain out-of-plane thickness with units `[length]`. Required only for
            `formulation="plane_strain"`; must be omitted for `formulation="axisymmetric"`.
        material_orientation_angles: Optional scalar or per-element angles with shape `(nelem,)`,
            in radians. Angles are unitless and rotate local material axes into the global 2D frame.

    Returns:
        Sparse stress-recovery operator with shape `(4 * npoint, 2 * nnode)`. Entries have units
        `[pressure / length]`, so multiplying by full nodal displacements with shape
        `(2 * nnode,)` and units `[length]` returns stress samples with shape `(4 * npoint,)` and
        units `[pressure]`.
    """

    dtype = _resolve_float_dtype(query.nodes, material_table, material_orientation_angles)
    material_ids_arr, material_table_arr = _normalize_materials(material_ids, material_table, dtype)
    assert material_ids_arr.shape == (
        query.elements.shape[0],
    ), f"material_ids must have shape ({query.elements.shape[0]},); got {material_ids_arr.shape}"
    material_orientation_angles_arr = _normalize_material_orientation_angles(
        material_orientation_angles,
        query.elements.shape[0],
        dtype,
    )
    normalized_formulation = _normalize_formulation(formulation)
    thickness_value = _normalize_thickness(normalized_formulation, thickness, dtype)
    binding = _dispatch_pair(dtype, _quad_mesh_stress_operator_f32, _quad_mesh_stress_operator_f64)
    return _coo_operator_from_binding(
        binding(
            np.asarray(query.nodes, dtype=dtype),
            query.elements,
            np.asarray(query.nearest_element_indices, dtype=np.uint64),
            np.asarray(query.nearest_element_reference_points, dtype=dtype),
            material_ids_arr,
            material_table_arr,
            material_orientation_angles_arr,
            query.element_type,
            _formulation_code(normalized_formulation),
            thickness_value,
        ),
        dtype,
    )

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 (nnode, 2) and units [length].

required
elements ArrayLike

Quad connectivity with shape (nelem, 4) for quad4 or (nelem, 9) for quad9. Entries are unitless node indices.

required
points ArrayLike

Query point coordinates with shape (npoint, 2) and units [length].

required
element_type str

Element family, either "quad4" or "quad9".

'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 [length], distances have units [length], reference coordinates are unitless,

QuadMeshQuery

and index arrays are unitless.

Source code in cfsem/solenoid_stress/fem2d.py
def 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.

    Args:
        nodes: Mesh node coordinates with shape `(nnode, 2)` and units `[length]`.
        elements: Quad connectivity with shape `(nelem, 4)` for `quad4` or `(nelem, 9)` for
            `quad9`. Entries are unitless node indices.
        points: Query point coordinates with shape `(npoint, 2)` and units `[length]`.
        element_type: Element family, either `"quad4"` or `"quad9"`.
        max_iterations: Maximum Newton/projection iterations per element. Unitless.

    Returns:
        Query data with nearest-node, nearest-element, and nearest-face arrays. Coordinate arrays
        have units `[length]`, distances have units `[length]`, reference coordinates are unitless,
        and index arrays are unitless.
    """

    dtype = _resolve_float_dtype(nodes, points)
    normalized_element_type = _normalize_element_type(element_type)
    nodes_arr = _normalize_nodes(nodes, dtype)
    elements_arr = _normalize_elements(
        elements,
        4 if normalized_element_type == "quad4" else 9,
    )
    points_arr = _normalize_query_points(points, dtype)
    binding = _dispatch_pair(dtype, _quad_mesh_query_f32, _quad_mesh_query_f64)
    data = binding(
        nodes_arr,
        elements_arr,
        points_arr,
        normalized_element_type,
        int(max_iterations),
    )

    return QuadMeshQuery(
        nodes=nodes_arr,
        elements=elements_arr,
        points=points_arr,
        element_type=normalized_element_type,
        nearest_node_indices=np.asarray(data["nearest_node_indices"], dtype=np.int64),
        nearest_node_points=np.asarray(data["nearest_node_points"], dtype=dtype).reshape(-1, 2),
        nearest_node_distances=np.asarray(data["nearest_node_distances"], dtype=dtype),
        nearest_element_indices=np.asarray(data["nearest_element_indices"], dtype=np.int64),
        nearest_element_reference_points=np.asarray(
            data["nearest_element_reference_points"],
            dtype=dtype,
        ).reshape(-1, 2),
        nearest_element_points=np.asarray(data["nearest_element_points"], dtype=dtype).reshape(-1, 2),
        nearest_element_distances=np.asarray(data["nearest_element_distances"], dtype=dtype),
        nearest_face_element_indices=np.asarray(data["nearest_face_element_indices"], dtype=np.int64),
        nearest_face_local_faces=np.asarray(data["nearest_face_local_faces"], dtype=np.int64),
        nearest_face_reference_coordinates=np.asarray(
            data["nearest_face_reference_coordinates"],
            dtype=dtype,
        ),
        nearest_face_points=np.asarray(data["nearest_face_points"], dtype=dtype).reshape(-1, 2),
        nearest_face_distances=np.asarray(data["nearest_face_distances"], dtype=dtype),
    )

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
def 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.

    Args:
        r: [m] (n x 1) array of radius points at which to evaluate the stress
        ri: [m] inner radius
        ro: [m] outer radius
        j: [A/m^2] current density
        bzi: [T] axial B-field at inner radius
        bzo: [T] axial B-field at outer radius
        poisson_ratio: [dimensionless] Material property; off-axis stress coupling term

    Returns:
        s_radial, s_hoop - each (n x 1) with units of [Pa]
    """

    # Terms shared between s_radial and s_hoop
    nu = poisson_ratio
    rho = r / ri
    alpha = ro / ri
    kappa = bzo / bzi

    jbr = j * bzi * ri  # [Pa]
    term1 = jbr / (alpha - 1.0)  # [Pa]
    term2 = (2.0 + nu) / 3.0
    term3 = (3.0 + nu) / 8.0
    term4 = alpha - kappa
    term5 = 1.0 - kappa

    # s_radial
    term6 = ((alpha**2 + alpha + 1.0 - alpha**2 / rho**2) / (alpha + 1.0)) - rho
    term7 = term3 * term5 * (alpha**2 + 1.0 - alpha**2 / rho**2 - rho**2)
    s_radial = term1 * (term2 * term4 * term6 - term7)  # [Pa]

    # s_hoop
    term8 = term2 * (alpha**2 + alpha + 1.0 + alpha**2 / rho**2) / (alpha + 1.0)
    term9 = rho * (1.0 + 2.0 * nu) / 3.0
    term10 = term4 * (term8 - term9)

    term11 = term3 * (alpha**2 + 1.0 + alpha**2 / rho**2)
    term12 = rho**2 * (1.0 + 3.0 * nu) / 8.0
    term13 = term5 * (term11 - term12)

    s_hoop = term1 * (term10 - term13)  # [Pa]

    return s_radial, s_hoop  # [Pa]

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
def 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

    Args:
        r: [m] radius at which to evaluate
        ri: [m] inner radius
        ro: [m] outer radius
        pin: [Pa] inside pressure
        pout: [Pa] outside pressure

    Returns:
        [Pa] radial stress
    """
    # Factors of pi cancel out
    # fmt: off
    s_radial = (
        ((pin * ri**2 - pout * ro**2) / (ro**2 - ri**2)) + \
        (ri**2 * ro**2 * (pout - pin) / (r**2 * (ro**2 - ri**2)))
    )
    # fmt: on

    return s_radial  # [Pa] radial stress

cfsem.solenoid_stress.s_hoop_thick_wall_cylinder

s_hoop_thick_wall_cylinder(
    r: NDArray,
    ri: float,
    ro: float,
    pin: float,
    pout: float,
) -> NDArray

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

Source code in cfsem/solenoid_stress/thick_wall_cylinder_handcalc.py
def s_hoop_thick_wall_cylinder(r: NDArray, ri: float, ro: float, pin: float, pout: float) -> NDArray:
    """
    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

    Args:
        r: [m] radius at which to evaluate
        ri: [m] inner radius
        ro: [m] outer radius
        pin: [Pa] inside pressure
        pout: [Pa] outside pressure

    Returns:
        [Pa] hoop stress
    """
    # Factors of pi cancel out
    # fmt: off
    s_hoop = (
        (pin * ri ** 2 - pout * ro ** 2) / (ro ** 2 - ri ** 2) -
        (ri ** 2 * ro ** 2 * (pout - pin) / (r ** 2 * (ro ** 2 - ri ** 2)))
    )
    # fmt: on

    return s_hoop  # [Pa] hoop stress