Skip to content

cvxpy

CVXPy backend for lowering symbolic expressions to CVXPy format.

This module implements the CVXPy lowering backend that converts symbolic expression AST nodes into CVXPy expressions for convex optimization. The lowering uses a visitor pattern where each expression type has a corresponding visitor method.

Architecture

The CVXPy lowerer follows a visitor pattern with centralized registration:

  1. Visitor Registration: The @visitor decorator registers handler functions for each expression type in the _CVXPY_VISITORS dictionary
  2. Dispatch: The dispatch() function looks up and calls the appropriate visitor based on the expression's type
  3. Recursive Lowering: Each visitor recursively lowers child expressions and composes CVXPy operations
  4. Translation Only: This module only translates expressions; CVXPy itself validates DCP (Disciplined Convex Programming) rules when the problem is constructed/solved
Key Features
  • Expression Translation: Converts symbolic AST to CVXPy expression format
  • Variable Management: Maps symbolic States/Controls to CVXPy variables through a variable_map dictionary
  • Parameter Support: Handles both constant parameters and CVXPy Parameters for efficient parameter sweeps
  • Constraint Generation: Produces CVXPy constraint objects from symbolic equality and inequality expressions
Backend Usage

CVXPy lowering is used for convex constraints in the SCP subproblem. Unlike JAX lowering (which happens early during problem construction), CVXPy lowering occurs later during Problem.initialize() when CVXPy variables are available. See lower_symbolic_expressions() in symbolic/lower.py for details.

CVXPy Variable Mapping

The lowerer requires a variable_map dictionary that maps symbolic variable names to CVXPy expressions. For trajectory optimization::

variable_map = {
    "x": cvxpy.Variable((n_x,)),  # State vector
    "u": cvxpy.Variable((n_u,)),  # Control vector
    "param_name": cvxpy.Parameter((3,)),  # Runtime parameters
}

States and Controls use their slices (assigned during unification) to extract the correct portion of the unified x and u vectors.

Example

Basic usage::

import cvxpy as cp
from openscvx.symbolic.lowerers.cvxpy import CvxpyLowerer
import openscvx as ox

# Create symbolic expression
x = ox.State("x", shape=(3,))
u = ox.Control("u", shape=(2,))
expr = ox.Norm(x)**2 + 0.1 * ox.Norm(u)**2

# Create CVXPy variables
cvx_x = cp.Variable(3)
cvx_u = cp.Variable(2)

# Lower to CVXPy
lowerer = CvxpyLowerer(variable_map={"x": cvx_x, "u": cvx_u})
cvx_expr = lowerer.lower(expr)

# Use in optimization problem
prob = cp.Problem(cp.Minimize(cvx_expr), constraints=[...])
prob.solve()

Constraint lowering::

# Symbolic constraint
constraint = ox.Norm(x) <= 1.0

# Lower to CVXPy constraint
cvx_constraint = lowerer.lower(constraint)

# Add to problem
prob = cp.Problem(cp.Minimize(cost), constraints=[cvx_constraint])
For Contributors

Adding Support for New Expression Types

To add support for a new symbolic expression type to CVXPy lowering:

  1. Define the visitor method in CvxpyLowerer with the @visitor decorator::

    @visitor(MyNewExpr) def _visit_my_new_expr(self, node: MyNewExpr) -> cp.Expression: # Lower child expressions recursively operand = self.lower(node.operand)

    # Return CVXPy expression
    return cp.my_operation(operand)
    
  2. Key requirements:

    • Use the @visitor(ExprType) decorator to register the handler
    • Method name should be visit (private, lowercase, snake_case)
    • Recursively lower all child expressions using self.lower()
    • Return a cp.Expression or cp.Constraint object
    • Use cp.* operations for CVXPy atoms
  3. DCP considerations:

    • This module only translates; CVXPy validates DCP rules
    • Document the mathematical properties in the docstring (convex, concave, affine)
    • For non-DCP operations, raise NotImplementedError with helpful message
    • See _visit_sin, _visit_cos, _visit_ctcs for examples
  4. Example patterns:

    • Unary operation: return cp.my_func(self.lower(node.operand))
    • Binary operation: return self.lower(node.left) + self.lower(node.right)
    • Constraints: return self.lower(node.lhs) <= self.lower(node.rhs)
    • Not supported: Raise NotImplementedError with guidance
  5. Testing: Ensure your visitor works with:

    • Simple expressions: Direct lowering to cp.Expression
    • Constraint validation: CVXPy accepts the result
    • DCP checking: CVXPy's problem.solve() validates correctly
See Also
  • lower_to_cvxpy(): Convenience wrapper for single expression lowering
  • JaxLowerer: Alternative backend for non-convex constraints and dynamics
  • lower_symbolic_expressions(): Main orchestrator in symbolic/lower.py
  • CVXPy documentation: https://www.cvxpy.org/

_CVXPY_VISITORS: Dict[Type[Expr], Callable] = {} module-attribute

Registry mapping expression types to their visitor functions.

CvxpyLowerer

CVXPy backend for lowering symbolic expressions to disciplined convex programs.

This class implements the visitor pattern for converting symbolic expression AST nodes to CVXPy expressions and constraints. Each expression type has a corresponding visitor method decorated with @visitor that handles the lowering logic.

The lowering process is recursive: each visitor lowers its child expressions first, then composes them into a CVXPy operation. CVXPy will validate DCP (Disciplined Convex Programming) compliance when the problem is constructed.

Attributes:

Name Type Description
variable_map dict

Dictionary mapping variable names to CVXPy expressions. Must include "x" for states and "u" for controls. May include parameter names mapped to CVXPy Parameter objects or constants.

Example

Lower an expression to CVXPy:

import cvxpy as cp
lowerer = CvxpyLowerer(variable_map={
    "x": cp.Variable(3),
    "u": cp.Variable(2),
})
expr = ox.Norm(x)**2 + 0.1 * ox.Norm(u)**2
cvx_expr = lowerer.lower(expr)
Note

The lowerer is stateful (stores variable_map) unlike JaxLowerer which is stateless. Variables must be registered before lowering expressions that reference them.

Source code in openscvx/symbolic/lowerers/cvxpy.py
 245
 246
 247
 248
 249
 250
 251
 252
 253
 254
 255
 256
 257
 258
 259
 260
 261
 262
 263
 264
 265
 266
 267
 268
 269
 270
 271
 272
 273
 274
 275
 276
 277
 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
 730
 731
 732
 733
 734
 735
 736
 737
 738
 739
 740
 741
 742
 743
 744
 745
 746
 747
 748
 749
 750
 751
 752
 753
 754
 755
 756
 757
 758
 759
 760
 761
 762
 763
 764
 765
 766
 767
 768
 769
 770
 771
 772
 773
 774
 775
 776
 777
 778
 779
 780
 781
 782
 783
 784
 785
 786
 787
 788
 789
 790
 791
 792
 793
 794
 795
 796
 797
 798
 799
 800
 801
 802
 803
 804
 805
 806
 807
 808
 809
 810
 811
 812
 813
 814
 815
 816
 817
 818
 819
 820
 821
 822
 823
 824
 825
 826
 827
 828
 829
 830
 831
 832
 833
 834
 835
 836
 837
 838
 839
 840
 841
 842
 843
 844
 845
 846
 847
 848
 849
 850
 851
 852
 853
 854
 855
 856
 857
 858
 859
 860
 861
 862
 863
 864
 865
 866
 867
 868
 869
 870
 871
 872
 873
 874
 875
 876
 877
 878
 879
 880
 881
 882
 883
 884
 885
 886
 887
 888
 889
 890
 891
 892
 893
 894
 895
 896
 897
 898
 899
 900
 901
 902
 903
 904
 905
 906
 907
 908
 909
 910
 911
 912
 913
 914
 915
 916
 917
 918
 919
 920
 921
 922
 923
 924
 925
 926
 927
 928
 929
 930
 931
 932
 933
 934
 935
 936
 937
 938
 939
 940
 941
 942
 943
 944
 945
 946
 947
 948
 949
 950
 951
 952
 953
 954
 955
 956
 957
 958
 959
 960
 961
 962
 963
 964
 965
 966
 967
 968
 969
 970
 971
 972
 973
 974
 975
 976
 977
 978
 979
 980
 981
 982
 983
 984
 985
 986
 987
 988
 989
 990
 991
 992
 993
 994
 995
 996
 997
 998
 999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
class CvxpyLowerer:
    """CVXPy backend for lowering symbolic expressions to disciplined convex programs.

    This class implements the visitor pattern for converting symbolic expression
    AST nodes to CVXPy expressions and constraints. Each expression type has a
    corresponding visitor method decorated with @visitor that handles the lowering
    logic.

    The lowering process is recursive: each visitor lowers its child expressions
    first, then composes them into a CVXPy operation. CVXPy will validate DCP
    (Disciplined Convex Programming) compliance when the problem is constructed.

    Attributes:
        variable_map (dict): Dictionary mapping variable names to CVXPy expressions.
            Must include "x" for states and "u" for controls. May include parameter
            names mapped to CVXPy Parameter objects or constants.

    Example:
        Lower an expression to CVXPy:

            import cvxpy as cp
            lowerer = CvxpyLowerer(variable_map={
                "x": cp.Variable(3),
                "u": cp.Variable(2),
            })
            expr = ox.Norm(x)**2 + 0.1 * ox.Norm(u)**2
            cvx_expr = lowerer.lower(expr)

    Note:
        The lowerer is stateful (stores variable_map) unlike JaxLowerer which
        is stateless. Variables must be registered before lowering expressions
        that reference them.
    """

    def __init__(self, variable_map: Dict[str, cp.Expression] = None):
        """Initialize the CVXPy lowerer.

        Args:
            variable_map: Dictionary mapping variable names to CVXPy expressions.
                For State/Control objects, keys should be "x" and "u" respectively.
                For Parameter objects, keys should match their names. If None, an
                empty dictionary is created.

        Example:
            Initialize the CVXPy lowerer with the variable map:

                cvx_x = cp.Variable(3, name="x")
                cvx_u = cp.Variable(2, name="u")
                lowerer = CvxpyLowerer({"x": cvx_x, "u": cvx_u})
        """
        self.variable_map = variable_map or {}

    def lower(self, expr: Expr) -> cp.Expression:
        """Lower a symbolic expression to a CVXPy expression.

        Main entry point for lowering. Delegates to dispatch() which looks up
        the appropriate visitor method based on the expression type.

        Args:
            expr: Symbolic expression to lower (any Expr subclass)

        Returns:
            CVXPy expression or constraint object. For arithmetic expressions,
            returns cp.Expression. For Equality/Inequality, returns cp.Constraint.

        Raises:
            NotImplementedError: If no visitor exists for the expression type
            ValueError: If required variables are not in variable_map

        Example:
            Lower an expression to a CVXPy expression:

                lowerer = CvxpyLowerer(variable_map={"x": cvx_x, "u": cvx_u})
                x = ox.State("x", shape=(3,))
                expr = ox.Norm(x)
                cvx_expr = lowerer.lower(expr)
        """
        return dispatch(self, expr)

    def register_variable(self, name: str, cvx_expr: cp.Expression):
        """Register a CVXPy variable/expression for use in lowering.

        Adds or updates a variable in the variable_map. Useful for dynamically
        adding variables after the lowerer has been created.

        Args:
            name: Variable name (e.g., "x", "u", or parameter name)
            cvx_expr: CVXPy expression to associate with the name

        Example:
            Register a variable:

                lowerer = CvxpyLowerer()
                lowerer.register_variable("x", cp.Variable(3))
                lowerer.register_variable("obs_center", cp.Parameter(3))
        """
        self.variable_map[name] = cvx_expr

    @visitor(Constant)
    def _visit_constant(self, node: Constant) -> cp.Expression:
        """Lower a constant value to a CVXPy constant.

        Wraps the constant's numpy array value in a CVXPy Constant expression.

        Args:
            node: Constant expression node

        Returns:
            CVXPy constant expression wrapping the value
        """
        return cp.Constant(node.value)

    @visitor(State)
    def _visit_state(self, node: State) -> cp.Expression:
        """Lower a state variable to a CVXPy expression.

        Extracts the appropriate slice from the unified state vector "x" using
        the slice assigned during unification. The "x" variable must exist in
        the variable_map.

        Args:
            node: State expression node

        Returns:
            CVXPy expression representing the state slice: x[slice]

        Raises:
            ValueError: If "x" is not found in variable_map
        """
        if "x" not in self.variable_map:
            raise ValueError("State vector 'x' not found in variable_map.")

        cvx_var = self.variable_map["x"]

        # If the state has a slice assigned, apply it
        if node._slice is not None:
            return cvx_var[node._slice]
        return cvx_var

    @visitor(Control)
    def _visit_control(self, node: Control) -> cp.Expression:
        """Lower a control variable to a CVXPy expression.

        Extracts the appropriate slice from the unified control vector "u" using
        the slice assigned during unification. The "u" variable must exist in
        the variable_map.

        Args:
            node: Control expression node

        Returns:
            CVXPy expression representing the control slice: u[slice]

        Raises:
            ValueError: If "u" is not found in variable_map
        """
        if "u" not in self.variable_map:
            raise ValueError("Control vector 'u' not found in variable_map.")

        cvx_var = self.variable_map["u"]

        # If the control has a slice assigned, apply it
        if node._slice is not None:
            return cvx_var[node._slice]
        return cvx_var

    @visitor(NodeReference)
    def _visit_node_reference(self, node: "NodeReference") -> cp.Expression:
        """Lower NodeReference - extract value at a specific trajectory node.

        NodeReference enables cross-node constraints by referencing state/control
        values at specific discrete time points. This requires the variable_map to
        contain full trajectory arrays (N, n_x) or (N, n_u) rather than single-node
        vectors.

        Args:
            node: NodeReference expression with base and node_idx

        Returns:
            CVXPy expression representing the variable at the specified node:
            x[node_idx, slice] or u[node_idx, slice]

        Raises:
            ValueError: If the required trajectory variable is not in variable_map
            ValueError: If the base variable has no slice assigned
            NotImplementedError: If the base is a compound expression

        Example:
            For cross-node constraint: position.at(5) - position.at(4) <= 0.1

            variable_map = {
                "x": cp.vstack([x_nonscaled[k] for k in range(N)]),  # (N, n_x)
            }
            # position.at(5) lowers to x[5, position._slice]

        Note:
            The node_idx is already resolved to an absolute integer index during
            expression construction, so negative indices are already handled.
        """
        from openscvx.symbolic.expr.control import Control
        from openscvx.symbolic.expr.state import State

        idx = node.node_idx

        if isinstance(node.base, State):
            if "x" not in self.variable_map:
                raise ValueError(
                    "State vector 'x' not found in variable_map. "
                    "For cross-node constraints, 'x' must be the full trajectory (N, n_x)."
                )

            cvx_var = self.variable_map["x"]  # Should be (N, n_x) for cross-node constraints

            # Apply slice if state has one assigned
            if node.base._slice is not None:
                return cvx_var[idx, node.base._slice]
            else:
                # No slice means this is the entire unified state vector
                return cvx_var[idx, :]

        elif isinstance(node.base, Control):
            if "u" not in self.variable_map:
                raise ValueError(
                    "Control vector 'u' not found in variable_map. "
                    "For cross-node constraints, 'u' must be the full trajectory (N, n_u)."
                )

            cvx_var = self.variable_map["u"]  # Should be (N, n_u) for cross-node constraints

            # Apply slice if control has one assigned
            if node.base._slice is not None:
                return cvx_var[idx, node.base._slice]
            else:
                # No slice means this is the entire unified control vector
                return cvx_var[idx, :]

        else:
            # Compound expression (e.g., position[0].at(5))
            # This is more complex - would need to lower base in single-node context
            raise NotImplementedError(
                "Compound expressions in NodeReference are not yet supported for CVXPy lowering. "
                f"Base expression type: {type(node.base).__name__}. "
                "Only State and Control NodeReferences are currently supported."
            )

    @visitor(CrossNodeConstraint)
    def _visit_cross_node_constraint(self, node: CrossNodeConstraint) -> cp.Constraint:
        """Lower CrossNodeConstraint to CVXPy constraint.

        CrossNodeConstraint wraps constraints that reference multiple trajectory
        nodes via NodeReference (e.g., rate limits like x.at(k) - x.at(k-1) <= r).

        For CVXPy lowering, this simply lowers the inner constraint. The NodeReference
        nodes within the constraint will handle extracting values from the full
        trajectory arrays (which must be provided in variable_map as "x" and "u").

        Args:
            node: CrossNodeConstraint expression wrapping the inner constraint

        Returns:
            CVXPy constraint object

        Note:
            The variable_map must contain full trajectory arrays:
                - "x": (N, n_x) CVXPy expression (e.g., cp.vstack(x_nonscaled))
                - "u": (N, n_u) CVXPy expression (e.g., cp.vstack(u_nonscaled))

            NodeReference visitors will index into these arrays using the fixed
            node indices baked into the expression.

        Example:
            For constraint: position.at(5) - position.at(4) <= max_step

            With variable_map = {"x": cp.vstack([x[k] for k in range(N)])}

            The lowered constraint evaluates:
                x[5, pos_slice] - x[4, pos_slice] <= max_step
        """
        # Simply lower the inner constraint - NodeReference handles indexing
        return self.lower(node.constraint)

    @visitor(Parameter)
    def _visit_parameter(self, node: Parameter) -> cp.Expression:
        """Lower a parameter to a CVXPy expression.

        Parameters are looked up by name in the variable_map. They can be mapped
        to CVXPy Parameter objects (for efficient parameter sweeps) or constants.

        Args:
            node: Parameter expression node

        Returns:
            CVXPy expression from variable_map (Parameter or constant)

        Raises:
            ValueError: If parameter name is not found in variable_map

        Note:
            For parameter sweeps without recompilation, map to cp.Parameter.
            For fixed values, map to cp.Constant or numpy arrays.
        """
        param_name = node.name
        if param_name in self.variable_map:
            return self.variable_map[param_name]
        else:
            raise ValueError(
                f"Parameter '{param_name}' not found in variable_map. "
                f"Add it during CVXPy lowering or use cp.Parameter for parameter sweeps."
            )

    @visitor(Add)
    def _visit_add(self, node: Add) -> cp.Expression:
        """Lower addition to CVXPy expression.

        Recursively lowers all terms and composes them with element-wise addition.
        Addition is affine and always DCP-compliant.

        Args:
            node: Add expression node with multiple terms

        Returns:
            CVXPy expression representing the sum of all terms
        """
        terms = [self.lower(term) for term in node.terms]
        result = terms[0]
        for term in terms[1:]:
            result = result + term
        return result

    @visitor(Sub)
    def _visit_sub(self, node: Sub) -> cp.Expression:
        """Lower subtraction to CVXPy expression (element-wise left - right).

        Subtraction is affine and always DCP-compliant.

        Args:
            node: Sub expression node

        Returns:
            CVXPy expression representing left - right
        """
        left = self.lower(node.left)
        right = self.lower(node.right)
        return left - right

    @visitor(Mul)
    def _visit_mul(self, node: Mul) -> cp.Expression:
        """Lower element-wise multiplication to CVXPy expression.

        Element-wise multiplication is DCP-compliant when at least one operand
        is constant. For quadratic forms, use MatMul instead.

        Args:
            node: Mul expression node with multiple factors

        Returns:
            CVXPy expression representing element-wise product

        Note:
            For convex optimization, typically one factor should be constant.
            CVXPy will raise a DCP error if the composition violates DCP rules.
        """
        factors = [self.lower(factor) for factor in node.factors]
        result = factors[0]
        for factor in factors[1:]:
            result = result * factor
        return result

    @visitor(Div)
    def _visit_div(self, node: Div) -> cp.Expression:
        """Lower element-wise division to CVXPy expression.

        Division is DCP-compliant when the denominator is constant or when
        the numerator is constant and the denominator is concave.

        Args:
            node: Div expression node

        Returns:
            CVXPy expression representing left / right

        Note:
            CVXPy will raise a DCP error if the division violates DCP rules.
        """
        left = self.lower(node.left)
        right = self.lower(node.right)
        return left / right

    @visitor(MatMul)
    def _visit_matmul(self, node: MatMul) -> cp.Expression:
        """Lower matrix multiplication to CVXPy expression using @ operator.

        Matrix multiplication is DCP-compliant when at least one operand is
        constant. Used for quadratic forms like x.T @ Q @ x.

        Args:
            node: MatMul expression node

        Returns:
            CVXPy expression representing left @ right
        """
        left = self.lower(node.left)
        right = self.lower(node.right)
        return left @ right

    @visitor(Neg)
    def _visit_neg(self, node: Neg) -> cp.Expression:
        """Lower negation (unary minus) to CVXPy expression.

        Negation preserves DCP properties (negating convex gives concave).

        Args:
            node: Neg expression node

        Returns:
            CVXPy expression representing -operand
        """
        operand = self.lower(node.operand)
        return -operand

    @visitor(Sum)
    def _visit_sum(self, node: Sum) -> cp.Expression:
        """Lower sum reduction to CVXPy expression (sums all elements).

        Sum preserves DCP properties (sum of convex is convex).

        Args:
            node: Sum expression node

        Returns:
            CVXPy scalar expression representing the sum of all elements
        """
        operand = self.lower(node.operand)
        return cp.sum(operand)

    @visitor(Norm)
    def _visit_norm(self, node: Norm) -> cp.Expression:
        """Lower norm operation to CVXPy expression.

        Norms are convex functions and commonly used in convex optimization.
        Supports all CVXPy norm types (1, 2, inf, "fro", etc.).

        Args:
            node: Norm expression node with ord attribute

        Returns:
            CVXPy expression representing the norm of the operand

        Note:
            Common norms: ord=2 (Euclidean), ord=1 (Manhattan), ord="inf"
        """
        operand = self.lower(node.operand)
        return cp.norm(operand, node.ord)

    @visitor(Index)
    def _visit_index(self, node: Index) -> cp.Expression:
        """Lower indexing/slicing operation to CVXPy expression.

        Indexing preserves DCP properties (indexing into convex is convex).

        Args:
            node: Index expression node

        Returns:
            CVXPy expression representing base[index]
        """
        base = self.lower(node.base)
        return base[node.index]

    @visitor(Concat)
    def _visit_concat(self, node: Concat) -> cp.Expression:
        """Lower concatenation to CVXPy expression.

        Concatenates expressions horizontally along axis 0. Scalars are
        promoted to 1D arrays before concatenation. Preserves DCP properties.

        Args:
            node: Concat expression node

        Returns:
            CVXPy expression representing horizontal concatenation

        Note:
            Uses cp.hstack for concatenation. Scalars are reshaped to (1,).
        """
        exprs = [self.lower(child) for child in node.exprs]
        # Ensure all expressions are at least 1D for concatenation
        exprs_1d = []
        for expr in exprs:
            if expr.ndim == 0:  # scalar
                exprs_1d.append(cp.reshape(expr, (1,), order="C"))
            else:
                exprs_1d.append(expr)
        return cp.hstack(exprs_1d)

    @visitor(Sin)
    def _visit_sin(self, node: Sin) -> cp.Expression:
        """Raise NotImplementedError for sine function.

        Sine is not DCP-compliant in CVXPy as it is neither convex nor concave.

        Args:
            node: Sin expression node

        Raises:
            NotImplementedError: Always raised since sine is not DCP-compliant

        Note:
            For constraints involving trigonometric functions:
            - Use piecewise-linear approximations, or
            - Handle in the JAX dynamics/constraint layer instead of CVXPy
        """
        raise NotImplementedError(
            "Trigonometric functions like Sin are not DCP-compliant in CVXPy. "
            "Consider using piecewise-linear approximations or handle these constraints "
            "in the dynamics (JAX) layer instead."
        )

    @visitor(Cos)
    def _visit_cos(self, node: Cos) -> cp.Expression:
        """Raise NotImplementedError for cosine function.

        Cosine is not DCP-compliant in CVXPy as it is neither convex nor concave.

        Args:
            node: Cos expression node

        Raises:
            NotImplementedError: Always raised since cosine is not DCP-compliant

        Note:
            For constraints involving trigonometric functions:
            - Use piecewise-linear approximations, or
            - Handle in the JAX dynamics/constraint layer instead of CVXPy
        """
        raise NotImplementedError(
            "Trigonometric functions like Cos are not DCP-compliant in CVXPy. "
            "Consider using piecewise-linear approximations or handle these constraints "
            "in the dynamics (JAX) layer instead."
        )

    @visitor(Tan)
    def _visit_tan(self, node: Tan) -> cp.Expression:
        """Raise NotImplementedError for tangent function.

        Tangent is not DCP-compliant in CVXPy as it is neither convex nor concave.

        Args:
            node: Tan expression node

        Raises:
            NotImplementedError: Always raised since tangent is not DCP-compliant

        Note:
            For constraints involving trigonometric functions:
            - Use piecewise-linear approximations, or
            - Handle in the JAX dynamics/constraint layer instead of CVXPy
        """
        raise NotImplementedError(
            "Trigonometric functions like Tan are not DCP-compliant in CVXPy. "
            "Consider using piecewise-linear approximations or handle these constraints "
            "in the dynamics (JAX) layer instead."
        )

    @visitor(Exp)
    def _visit_exp(self, node: Exp) -> cp.Expression:
        """Lower exponential function to CVXPy expression.

        Exponential is a convex function and DCP-compliant when used in
        appropriate contexts (e.g., minimizing exp(x) or constraints like
        exp(x) <= c).

        Args:
            node: Exp expression node

        Returns:
            CVXPy expression representing exp(operand)

        Note:
            Exponential is convex increasing, so it's valid in:
            - Objective: minimize exp(x)
            - Constraints: exp(x) <= c (convex constraint)
        """
        operand = self.lower(node.operand)
        return cp.exp(operand)

    @visitor(Log)
    def _visit_log(self, node: Log) -> cp.Expression:
        """Lower natural logarithm to CVXPy expression.

        Logarithm is a concave function and DCP-compliant when used in
        appropriate contexts (e.g., maximizing log(x) or constraints like
        log(x) >= c).

        Args:
            node: Log expression node

        Returns:
            CVXPy expression representing log(operand)

        Note:
            Logarithm is concave increasing, so it's valid in:
            - Objective: maximize log(x)
            - Constraints: log(x) >= c (concave constraint, or equivalently c <= log(x))
        """
        operand = self.lower(node.operand)
        return cp.log(operand)

    @visitor(Abs)
    def _visit_abs(self, node: Abs) -> cp.Expression:
        """Lower absolute value to CVXPy expression.

        Absolute value is a convex function and DCP-compliant when used in
        appropriate contexts (e.g., minimizing |x| or constraints like |x| <= c).

        Args:
            node: Abs expression node

        Returns:
            CVXPy expression representing |operand|

        Note:
            Absolute value is convex, so it's valid in:
            - Objective: minimize abs(x)
            - Constraints: abs(x) <= c (convex constraint)
        """
        operand = self.lower(node.operand)
        return cp.abs(operand)

    @visitor(Equality)
    def _visit_equality(self, node: Equality) -> cp.Constraint:
        """Lower equality constraint to CVXPy constraint (lhs == rhs).

        Equality constraints require affine expressions on both sides for
        DCP compliance.

        Args:
            node: Equality constraint node

        Returns:
            CVXPy equality constraint object

        Note:
            For DCP compliance, both lhs and rhs must be affine. CVXPy will
            raise a DCP error if either side is non-affine.
        """
        left = self.lower(node.lhs)
        right = self.lower(node.rhs)
        return left == right

    @visitor(Inequality)
    def _visit_inequality(self, node: Inequality) -> cp.Constraint:
        """Lower inequality constraint to CVXPy constraint (lhs <= rhs).

        Inequality constraints must satisfy DCP rules: convex <= concave.

        Args:
            node: Inequality constraint node

        Returns:
            CVXPy inequality constraint object

        Note:
            For DCP compliance: lhs must be convex and rhs must be concave.
            Common form: convex_expr(x) <= constant
        """
        left = self.lower(node.lhs)
        right = self.lower(node.rhs)
        return left <= right

    @visitor(CTCS)
    def _visit_ctcs(self, node: CTCS) -> cp.Expression:
        """Raise NotImplementedError for CTCS constraints.

        CTCS (Continuous-Time Constraint Satisfaction) constraints are handled
        through dynamics augmentation using JAX, not CVXPy. They represent
        non-convex continuous-time constraints.

        Args:
            node: CTCS constraint node

        Raises:
            NotImplementedError: Always raised since CTCS uses JAX, not CVXPy

        Note:
            CTCS constraints are lowered to JAX during dynamics augmentation.
            They add virtual states and controls to enforce constraints over
            continuous time intervals. See JaxLowerer.visit_ctcs() instead.
        """
        raise NotImplementedError(
            "CTCS constraints are for continuous-time constraint satisfaction and "
            "should be handled through dynamics augmentation with JAX lowering, "
            "not CVXPy lowering. CTCS constraints represent non-convex dynamics "
            "augmentation."
        )

    @visitor(PositivePart)
    def _visit_pos(self, node: PositivePart) -> cp.Expression:
        """Lower positive part function to CVXPy.

        Computes max(x, 0), which is convex. Used in penalty methods for
        inequality constraints.

        Args:
            node: PositivePart expression node

        Returns:
            CVXPy expression representing max(operand, 0)

        Note:
            Positive part is convex and commonly used in hinge loss and
            penalty methods for inequality constraints.
        """
        operand = self.lower(node.x)
        return cp.maximum(operand, 0.0)

    @visitor(Square)
    def _visit_square(self, node: Square) -> cp.Expression:
        """Lower square function to CVXPy.

        Computes x^2, which is convex. Used in quadratic penalty methods
        and least-squares objectives.

        Args:
            node: Square expression node

        Returns:
            CVXPy expression representing operand^2

        Note:
            Square is convex increasing for x >= 0 and convex decreasing for
            x <= 0. It's always convex overall.
        """
        operand = self.lower(node.x)
        return cp.square(operand)

    @visitor(Huber)
    def _visit_huber(self, node: Huber) -> cp.Expression:
        """Lower Huber penalty function to CVXPy.

        Huber penalty is quadratic for small values and linear for large values,
        providing robustness to outliers. It is convex and DCP-compliant.

        The Huber function is defined as:
        - |x| <= delta: 0.5 * x^2
        - |x| > delta: delta * (|x| - 0.5 * delta)

        Args:
            node: Huber expression node with delta parameter

        Returns:
            CVXPy expression representing Huber penalty

        Note:
            Huber loss is convex and combines the benefits of squared error
            (smooth, differentiable) and absolute error (robust to outliers).
        """
        operand = self.lower(node.x)
        return cp.huber(operand, M=node.delta)

    @visitor(SmoothReLU)
    def _visit_srelu(self, node: SmoothReLU) -> cp.Expression:
        """Lower smooth ReLU penalty function to CVXPy.

        Smooth approximation to ReLU: sqrt(max(x, 0)^2 + c^2) - c
        Differentiable everywhere, approaches ReLU as c -> 0. Convex.

        Args:
            node: SmoothReLU expression node with smoothing parameter c

        Returns:
            CVXPy expression representing smooth ReLU penalty

        Note:
            This provides a smooth, convex approximation to the ReLU function
            max(x, 0). The parameter c controls the smoothness: smaller c gives
            a better approximation but less smoothness.
        """
        operand = self.lower(node.x)
        c = node.c
        # smooth_relu(x) = sqrt(max(x, 0)^2 + c^2) - c
        pos_part = cp.maximum(operand, 0.0)
        # For SmoothReLU, we use the 2-norm formulation
        return cp.sqrt(cp.sum_squares(pos_part) + c**2) - c

    @visitor(Sqrt)
    def _visit_sqrt(self, node: Sqrt) -> cp.Expression:
        """Lower square root to CVXPy expression.

        Square root is concave and DCP-compliant when used appropriately
        (e.g., maximizing sqrt(x) or constraints like sqrt(x) >= c).

        Args:
            node: Sqrt expression node

        Returns:
            CVXPy expression representing sqrt(operand)

        Note:
            Square root is concave increasing for x > 0. Valid in:
            - Objective: maximize sqrt(x)
            - Constraints: sqrt(x) >= c (concave constraint)
        """
        operand = self.lower(node.operand)
        return cp.sqrt(operand)

    @visitor(Max)
    def _visit_max(self, node: Max) -> cp.Expression:
        """Lower element-wise maximum to CVXPy expression.

        Maximum is convex (pointwise max of convex functions is convex).

        Args:
            node: Max expression node with multiple operands

        Returns:
            CVXPy expression representing element-wise maximum

        Note:
            For multiple operands, chains binary maximum operations.
            Maximum preserves convexity.
        """
        operands = [self.lower(op) for op in node.operands]
        # CVXPy's maximum can take multiple arguments
        if len(operands) == 2:
            return cp.maximum(operands[0], operands[1])
        else:
            # For more than 2 operands, chain maximum calls
            result = cp.maximum(operands[0], operands[1])
            for op in operands[2:]:
                result = cp.maximum(result, op)
            return result

    @visitor(LogSumExp)
    def _visit_logsumexp(self, node: LogSumExp) -> cp.Expression:
        """Lower log-sum-exp to CVXPy expression.

        Log-sum-exp is convex and is a smooth approximation to the maximum function.
        CVXPy's log_sum_exp atom computes log(sum(exp(x_i))) for stacked operands.

        Args:
            node: LogSumExp expression node with multiple operands

        Returns:
            CVXPy expression representing log-sum-exp

        Note:
            Log-sum-exp is convex and DCP-compliant. It satisfies:
            max(x₁, ..., xₙ) ≤ logsumexp(x₁, ..., xₙ) ≤ max(x₁, ..., xₙ) + log(n)
        """
        operands = [self.lower(op) for op in node.operands]

        # CVXPy's log_sum_exp expects a stacked expression with an axis parameter
        # For element-wise log-sum-exp, we stack along a new axis and reduce along it
        if len(operands) == 1:
            return operands[0]

        # Stack operands along a new axis (axis 0) and compute log_sum_exp along that axis
        stacked = cp.vstack(operands)
        return cp.log_sum_exp(stacked, axis=0)

    @visitor(Transpose)
    def _visit_transpose(self, node: Transpose) -> cp.Expression:
        """Lower matrix transpose to CVXPy expression.

        Transpose preserves DCP properties (transpose of convex is convex).

        Args:
            node: Transpose expression node

        Returns:
            CVXPy expression representing operand.T
        """
        operand = self.lower(node.operand)
        return operand.T

    @visitor(Power)
    def _visit_power(self, node: Power) -> cp.Expression:
        """Lower element-wise power (base**exponent) to CVXPy expression.

        Power is DCP-compliant for specific exponent values:
        - exponent >= 1: convex (when base >= 0)
        - 0 <= exponent <= 1: concave (when base >= 0)

        Args:
            node: Power expression node

        Returns:
            CVXPy expression representing base**exponent

        Note:
            CVXPy will verify DCP compliance at problem construction time.
            Common convex cases: x^2, x^3, x^4 (even powers)
        """
        base = self.lower(node.base)
        exponent = self.lower(node.exponent)
        return cp.power(base, exponent)

    @visitor(Stack)
    def _visit_stack(self, node: Stack) -> cp.Expression:
        """Lower vertical stacking to CVXPy expression.

        Stacks expressions vertically using cp.vstack. Preserves DCP properties.

        Args:
            node: Stack expression node with multiple rows

        Returns:
            CVXPy expression representing vertical stack of rows

        Note:
            Each row is stacked along axis 0 to create a 2D array.
        """
        rows = [self.lower(row) for row in node.rows]
        # Stack rows vertically
        return cp.vstack(rows)
_visit_abs(node: Abs) -> cp.Expression

Lower absolute value to CVXPy expression.

Absolute value is a convex function and DCP-compliant when used in appropriate contexts (e.g., minimizing |x| or constraints like |x| <= c).

Parameters:

Name Type Description Default
node Abs

Abs expression node

required

Returns:

Type Description
Expression

CVXPy expression representing |operand|

Note

Absolute value is convex, so it's valid in: - Objective: minimize abs(x) - Constraints: abs(x) <= c (convex constraint)

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Abs)
def _visit_abs(self, node: Abs) -> cp.Expression:
    """Lower absolute value to CVXPy expression.

    Absolute value is a convex function and DCP-compliant when used in
    appropriate contexts (e.g., minimizing |x| or constraints like |x| <= c).

    Args:
        node: Abs expression node

    Returns:
        CVXPy expression representing |operand|

    Note:
        Absolute value is convex, so it's valid in:
        - Objective: minimize abs(x)
        - Constraints: abs(x) <= c (convex constraint)
    """
    operand = self.lower(node.operand)
    return cp.abs(operand)
_visit_add(node: Add) -> cp.Expression

Lower addition to CVXPy expression.

Recursively lowers all terms and composes them with element-wise addition. Addition is affine and always DCP-compliant.

Parameters:

Name Type Description Default
node Add

Add expression node with multiple terms

required

Returns:

Type Description
Expression

CVXPy expression representing the sum of all terms

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Add)
def _visit_add(self, node: Add) -> cp.Expression:
    """Lower addition to CVXPy expression.

    Recursively lowers all terms and composes them with element-wise addition.
    Addition is affine and always DCP-compliant.

    Args:
        node: Add expression node with multiple terms

    Returns:
        CVXPy expression representing the sum of all terms
    """
    terms = [self.lower(term) for term in node.terms]
    result = terms[0]
    for term in terms[1:]:
        result = result + term
    return result
_visit_concat(node: Concat) -> cp.Expression

Lower concatenation to CVXPy expression.

Concatenates expressions horizontally along axis 0. Scalars are promoted to 1D arrays before concatenation. Preserves DCP properties.

Parameters:

Name Type Description Default
node Concat

Concat expression node

required

Returns:

Type Description
Expression

CVXPy expression representing horizontal concatenation

Note

Uses cp.hstack for concatenation. Scalars are reshaped to (1,).

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Concat)
def _visit_concat(self, node: Concat) -> cp.Expression:
    """Lower concatenation to CVXPy expression.

    Concatenates expressions horizontally along axis 0. Scalars are
    promoted to 1D arrays before concatenation. Preserves DCP properties.

    Args:
        node: Concat expression node

    Returns:
        CVXPy expression representing horizontal concatenation

    Note:
        Uses cp.hstack for concatenation. Scalars are reshaped to (1,).
    """
    exprs = [self.lower(child) for child in node.exprs]
    # Ensure all expressions are at least 1D for concatenation
    exprs_1d = []
    for expr in exprs:
        if expr.ndim == 0:  # scalar
            exprs_1d.append(cp.reshape(expr, (1,), order="C"))
        else:
            exprs_1d.append(expr)
    return cp.hstack(exprs_1d)
_visit_constant(node: Constant) -> cp.Expression

Lower a constant value to a CVXPy constant.

Wraps the constant's numpy array value in a CVXPy Constant expression.

Parameters:

Name Type Description Default
node Constant

Constant expression node

required

Returns:

Type Description
Expression

CVXPy constant expression wrapping the value

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Constant)
def _visit_constant(self, node: Constant) -> cp.Expression:
    """Lower a constant value to a CVXPy constant.

    Wraps the constant's numpy array value in a CVXPy Constant expression.

    Args:
        node: Constant expression node

    Returns:
        CVXPy constant expression wrapping the value
    """
    return cp.Constant(node.value)
_visit_control(node: Control) -> cp.Expression

Lower a control variable to a CVXPy expression.

Extracts the appropriate slice from the unified control vector "u" using the slice assigned during unification. The "u" variable must exist in the variable_map.

Parameters:

Name Type Description Default
node Control

Control expression node

required

Returns:

Type Description
Expression

CVXPy expression representing the control slice: u[slice]

Raises:

Type Description
ValueError

If "u" is not found in variable_map

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Control)
def _visit_control(self, node: Control) -> cp.Expression:
    """Lower a control variable to a CVXPy expression.

    Extracts the appropriate slice from the unified control vector "u" using
    the slice assigned during unification. The "u" variable must exist in
    the variable_map.

    Args:
        node: Control expression node

    Returns:
        CVXPy expression representing the control slice: u[slice]

    Raises:
        ValueError: If "u" is not found in variable_map
    """
    if "u" not in self.variable_map:
        raise ValueError("Control vector 'u' not found in variable_map.")

    cvx_var = self.variable_map["u"]

    # If the control has a slice assigned, apply it
    if node._slice is not None:
        return cvx_var[node._slice]
    return cvx_var
_visit_cos(node: Cos) -> cp.Expression

Raise NotImplementedError for cosine function.

Cosine is not DCP-compliant in CVXPy as it is neither convex nor concave.

Parameters:

Name Type Description Default
node Cos

Cos expression node

required

Raises:

Type Description
NotImplementedError

Always raised since cosine is not DCP-compliant

Note

For constraints involving trigonometric functions: - Use piecewise-linear approximations, or - Handle in the JAX dynamics/constraint layer instead of CVXPy

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Cos)
def _visit_cos(self, node: Cos) -> cp.Expression:
    """Raise NotImplementedError for cosine function.

    Cosine is not DCP-compliant in CVXPy as it is neither convex nor concave.

    Args:
        node: Cos expression node

    Raises:
        NotImplementedError: Always raised since cosine is not DCP-compliant

    Note:
        For constraints involving trigonometric functions:
        - Use piecewise-linear approximations, or
        - Handle in the JAX dynamics/constraint layer instead of CVXPy
    """
    raise NotImplementedError(
        "Trigonometric functions like Cos are not DCP-compliant in CVXPy. "
        "Consider using piecewise-linear approximations or handle these constraints "
        "in the dynamics (JAX) layer instead."
    )
_visit_cross_node_constraint(node: CrossNodeConstraint) -> cp.Constraint

Lower CrossNodeConstraint to CVXPy constraint.

CrossNodeConstraint wraps constraints that reference multiple trajectory nodes via NodeReference (e.g., rate limits like x.at(k) - x.at(k-1) <= r).

For CVXPy lowering, this simply lowers the inner constraint. The NodeReference nodes within the constraint will handle extracting values from the full trajectory arrays (which must be provided in variable_map as "x" and "u").

Parameters:

Name Type Description Default
node CrossNodeConstraint

CrossNodeConstraint expression wrapping the inner constraint

required

Returns:

Type Description
Constraint

CVXPy constraint object

Note

The variable_map must contain full trajectory arrays: - "x": (N, n_x) CVXPy expression (e.g., cp.vstack(x_nonscaled)) - "u": (N, n_u) CVXPy expression (e.g., cp.vstack(u_nonscaled))

NodeReference visitors will index into these arrays using the fixed node indices baked into the expression.

Example

For constraint: position.at(5) - position.at(4) <= max_step

With variable_map = {"x": cp.vstack([x[k] for k in range(N)])}

The lowered constraint evaluates: x[5, pos_slice] - x[4, pos_slice] <= max_step

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(CrossNodeConstraint)
def _visit_cross_node_constraint(self, node: CrossNodeConstraint) -> cp.Constraint:
    """Lower CrossNodeConstraint to CVXPy constraint.

    CrossNodeConstraint wraps constraints that reference multiple trajectory
    nodes via NodeReference (e.g., rate limits like x.at(k) - x.at(k-1) <= r).

    For CVXPy lowering, this simply lowers the inner constraint. The NodeReference
    nodes within the constraint will handle extracting values from the full
    trajectory arrays (which must be provided in variable_map as "x" and "u").

    Args:
        node: CrossNodeConstraint expression wrapping the inner constraint

    Returns:
        CVXPy constraint object

    Note:
        The variable_map must contain full trajectory arrays:
            - "x": (N, n_x) CVXPy expression (e.g., cp.vstack(x_nonscaled))
            - "u": (N, n_u) CVXPy expression (e.g., cp.vstack(u_nonscaled))

        NodeReference visitors will index into these arrays using the fixed
        node indices baked into the expression.

    Example:
        For constraint: position.at(5) - position.at(4) <= max_step

        With variable_map = {"x": cp.vstack([x[k] for k in range(N)])}

        The lowered constraint evaluates:
            x[5, pos_slice] - x[4, pos_slice] <= max_step
    """
    # Simply lower the inner constraint - NodeReference handles indexing
    return self.lower(node.constraint)
_visit_ctcs(node: CTCS) -> cp.Expression

Raise NotImplementedError for CTCS constraints.

CTCS (Continuous-Time Constraint Satisfaction) constraints are handled through dynamics augmentation using JAX, not CVXPy. They represent non-convex continuous-time constraints.

Parameters:

Name Type Description Default
node CTCS

CTCS constraint node

required

Raises:

Type Description
NotImplementedError

Always raised since CTCS uses JAX, not CVXPy

Note

CTCS constraints are lowered to JAX during dynamics augmentation. They add virtual states and controls to enforce constraints over continuous time intervals. See JaxLowerer.visit_ctcs() instead.

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(CTCS)
def _visit_ctcs(self, node: CTCS) -> cp.Expression:
    """Raise NotImplementedError for CTCS constraints.

    CTCS (Continuous-Time Constraint Satisfaction) constraints are handled
    through dynamics augmentation using JAX, not CVXPy. They represent
    non-convex continuous-time constraints.

    Args:
        node: CTCS constraint node

    Raises:
        NotImplementedError: Always raised since CTCS uses JAX, not CVXPy

    Note:
        CTCS constraints are lowered to JAX during dynamics augmentation.
        They add virtual states and controls to enforce constraints over
        continuous time intervals. See JaxLowerer.visit_ctcs() instead.
    """
    raise NotImplementedError(
        "CTCS constraints are for continuous-time constraint satisfaction and "
        "should be handled through dynamics augmentation with JAX lowering, "
        "not CVXPy lowering. CTCS constraints represent non-convex dynamics "
        "augmentation."
    )
_visit_div(node: Div) -> cp.Expression

Lower element-wise division to CVXPy expression.

Division is DCP-compliant when the denominator is constant or when the numerator is constant and the denominator is concave.

Parameters:

Name Type Description Default
node Div

Div expression node

required

Returns:

Type Description
Expression

CVXPy expression representing left / right

Note

CVXPy will raise a DCP error if the division violates DCP rules.

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Div)
def _visit_div(self, node: Div) -> cp.Expression:
    """Lower element-wise division to CVXPy expression.

    Division is DCP-compliant when the denominator is constant or when
    the numerator is constant and the denominator is concave.

    Args:
        node: Div expression node

    Returns:
        CVXPy expression representing left / right

    Note:
        CVXPy will raise a DCP error if the division violates DCP rules.
    """
    left = self.lower(node.left)
    right = self.lower(node.right)
    return left / right
_visit_equality(node: Equality) -> cp.Constraint

Lower equality constraint to CVXPy constraint (lhs == rhs).

Equality constraints require affine expressions on both sides for DCP compliance.

Parameters:

Name Type Description Default
node Equality

Equality constraint node

required

Returns:

Type Description
Constraint

CVXPy equality constraint object

Note

For DCP compliance, both lhs and rhs must be affine. CVXPy will raise a DCP error if either side is non-affine.

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Equality)
def _visit_equality(self, node: Equality) -> cp.Constraint:
    """Lower equality constraint to CVXPy constraint (lhs == rhs).

    Equality constraints require affine expressions on both sides for
    DCP compliance.

    Args:
        node: Equality constraint node

    Returns:
        CVXPy equality constraint object

    Note:
        For DCP compliance, both lhs and rhs must be affine. CVXPy will
        raise a DCP error if either side is non-affine.
    """
    left = self.lower(node.lhs)
    right = self.lower(node.rhs)
    return left == right
_visit_exp(node: Exp) -> cp.Expression

Lower exponential function to CVXPy expression.

Exponential is a convex function and DCP-compliant when used in appropriate contexts (e.g., minimizing exp(x) or constraints like exp(x) <= c).

Parameters:

Name Type Description Default
node Exp

Exp expression node

required

Returns:

Type Description
Expression

CVXPy expression representing exp(operand)

Note

Exponential is convex increasing, so it's valid in: - Objective: minimize exp(x) - Constraints: exp(x) <= c (convex constraint)

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Exp)
def _visit_exp(self, node: Exp) -> cp.Expression:
    """Lower exponential function to CVXPy expression.

    Exponential is a convex function and DCP-compliant when used in
    appropriate contexts (e.g., minimizing exp(x) or constraints like
    exp(x) <= c).

    Args:
        node: Exp expression node

    Returns:
        CVXPy expression representing exp(operand)

    Note:
        Exponential is convex increasing, so it's valid in:
        - Objective: minimize exp(x)
        - Constraints: exp(x) <= c (convex constraint)
    """
    operand = self.lower(node.operand)
    return cp.exp(operand)
_visit_huber(node: Huber) -> cp.Expression

Lower Huber penalty function to CVXPy.

Huber penalty is quadratic for small values and linear for large values, providing robustness to outliers. It is convex and DCP-compliant.

The Huber function is defined as: - |x| <= delta: 0.5 * x^2 - |x| > delta: delta * (|x| - 0.5 * delta)

Parameters:

Name Type Description Default
node Huber

Huber expression node with delta parameter

required

Returns:

Type Description
Expression

CVXPy expression representing Huber penalty

Note

Huber loss is convex and combines the benefits of squared error (smooth, differentiable) and absolute error (robust to outliers).

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Huber)
def _visit_huber(self, node: Huber) -> cp.Expression:
    """Lower Huber penalty function to CVXPy.

    Huber penalty is quadratic for small values and linear for large values,
    providing robustness to outliers. It is convex and DCP-compliant.

    The Huber function is defined as:
    - |x| <= delta: 0.5 * x^2
    - |x| > delta: delta * (|x| - 0.5 * delta)

    Args:
        node: Huber expression node with delta parameter

    Returns:
        CVXPy expression representing Huber penalty

    Note:
        Huber loss is convex and combines the benefits of squared error
        (smooth, differentiable) and absolute error (robust to outliers).
    """
    operand = self.lower(node.x)
    return cp.huber(operand, M=node.delta)
_visit_index(node: Index) -> cp.Expression

Lower indexing/slicing operation to CVXPy expression.

Indexing preserves DCP properties (indexing into convex is convex).

Parameters:

Name Type Description Default
node Index

Index expression node

required

Returns:

Type Description
Expression

CVXPy expression representing base[index]

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Index)
def _visit_index(self, node: Index) -> cp.Expression:
    """Lower indexing/slicing operation to CVXPy expression.

    Indexing preserves DCP properties (indexing into convex is convex).

    Args:
        node: Index expression node

    Returns:
        CVXPy expression representing base[index]
    """
    base = self.lower(node.base)
    return base[node.index]
_visit_inequality(node: Inequality) -> cp.Constraint

Lower inequality constraint to CVXPy constraint (lhs <= rhs).

Inequality constraints must satisfy DCP rules: convex <= concave.

Parameters:

Name Type Description Default
node Inequality

Inequality constraint node

required

Returns:

Type Description
Constraint

CVXPy inequality constraint object

Note

For DCP compliance: lhs must be convex and rhs must be concave. Common form: convex_expr(x) <= constant

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Inequality)
def _visit_inequality(self, node: Inequality) -> cp.Constraint:
    """Lower inequality constraint to CVXPy constraint (lhs <= rhs).

    Inequality constraints must satisfy DCP rules: convex <= concave.

    Args:
        node: Inequality constraint node

    Returns:
        CVXPy inequality constraint object

    Note:
        For DCP compliance: lhs must be convex and rhs must be concave.
        Common form: convex_expr(x) <= constant
    """
    left = self.lower(node.lhs)
    right = self.lower(node.rhs)
    return left <= right
_visit_log(node: Log) -> cp.Expression

Lower natural logarithm to CVXPy expression.

Logarithm is a concave function and DCP-compliant when used in appropriate contexts (e.g., maximizing log(x) or constraints like log(x) >= c).

Parameters:

Name Type Description Default
node Log

Log expression node

required

Returns:

Type Description
Expression

CVXPy expression representing log(operand)

Note

Logarithm is concave increasing, so it's valid in: - Objective: maximize log(x) - Constraints: log(x) >= c (concave constraint, or equivalently c <= log(x))

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Log)
def _visit_log(self, node: Log) -> cp.Expression:
    """Lower natural logarithm to CVXPy expression.

    Logarithm is a concave function and DCP-compliant when used in
    appropriate contexts (e.g., maximizing log(x) or constraints like
    log(x) >= c).

    Args:
        node: Log expression node

    Returns:
        CVXPy expression representing log(operand)

    Note:
        Logarithm is concave increasing, so it's valid in:
        - Objective: maximize log(x)
        - Constraints: log(x) >= c (concave constraint, or equivalently c <= log(x))
    """
    operand = self.lower(node.operand)
    return cp.log(operand)
_visit_logsumexp(node: LogSumExp) -> cp.Expression

Lower log-sum-exp to CVXPy expression.

Log-sum-exp is convex and is a smooth approximation to the maximum function. CVXPy's log_sum_exp atom computes log(sum(exp(x_i))) for stacked operands.

Parameters:

Name Type Description Default
node LogSumExp

LogSumExp expression node with multiple operands

required

Returns:

Type Description
Expression

CVXPy expression representing log-sum-exp

Note

Log-sum-exp is convex and DCP-compliant. It satisfies: max(x₁, ..., xₙ) ≤ logsumexp(x₁, ..., xₙ) ≤ max(x₁, ..., xₙ) + log(n)

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(LogSumExp)
def _visit_logsumexp(self, node: LogSumExp) -> cp.Expression:
    """Lower log-sum-exp to CVXPy expression.

    Log-sum-exp is convex and is a smooth approximation to the maximum function.
    CVXPy's log_sum_exp atom computes log(sum(exp(x_i))) for stacked operands.

    Args:
        node: LogSumExp expression node with multiple operands

    Returns:
        CVXPy expression representing log-sum-exp

    Note:
        Log-sum-exp is convex and DCP-compliant. It satisfies:
        max(x₁, ..., xₙ) ≤ logsumexp(x₁, ..., xₙ) ≤ max(x₁, ..., xₙ) + log(n)
    """
    operands = [self.lower(op) for op in node.operands]

    # CVXPy's log_sum_exp expects a stacked expression with an axis parameter
    # For element-wise log-sum-exp, we stack along a new axis and reduce along it
    if len(operands) == 1:
        return operands[0]

    # Stack operands along a new axis (axis 0) and compute log_sum_exp along that axis
    stacked = cp.vstack(operands)
    return cp.log_sum_exp(stacked, axis=0)
_visit_matmul(node: MatMul) -> cp.Expression

Lower matrix multiplication to CVXPy expression using @ operator.

Matrix multiplication is DCP-compliant when at least one operand is constant. Used for quadratic forms like x.T @ Q @ x.

Parameters:

Name Type Description Default
node MatMul

MatMul expression node

required

Returns:

Type Description
Expression

CVXPy expression representing left @ right

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(MatMul)
def _visit_matmul(self, node: MatMul) -> cp.Expression:
    """Lower matrix multiplication to CVXPy expression using @ operator.

    Matrix multiplication is DCP-compliant when at least one operand is
    constant. Used for quadratic forms like x.T @ Q @ x.

    Args:
        node: MatMul expression node

    Returns:
        CVXPy expression representing left @ right
    """
    left = self.lower(node.left)
    right = self.lower(node.right)
    return left @ right
_visit_max(node: Max) -> cp.Expression

Lower element-wise maximum to CVXPy expression.

Maximum is convex (pointwise max of convex functions is convex).

Parameters:

Name Type Description Default
node Max

Max expression node with multiple operands

required

Returns:

Type Description
Expression

CVXPy expression representing element-wise maximum

Note

For multiple operands, chains binary maximum operations. Maximum preserves convexity.

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Max)
def _visit_max(self, node: Max) -> cp.Expression:
    """Lower element-wise maximum to CVXPy expression.

    Maximum is convex (pointwise max of convex functions is convex).

    Args:
        node: Max expression node with multiple operands

    Returns:
        CVXPy expression representing element-wise maximum

    Note:
        For multiple operands, chains binary maximum operations.
        Maximum preserves convexity.
    """
    operands = [self.lower(op) for op in node.operands]
    # CVXPy's maximum can take multiple arguments
    if len(operands) == 2:
        return cp.maximum(operands[0], operands[1])
    else:
        # For more than 2 operands, chain maximum calls
        result = cp.maximum(operands[0], operands[1])
        for op in operands[2:]:
            result = cp.maximum(result, op)
        return result
_visit_mul(node: Mul) -> cp.Expression

Lower element-wise multiplication to CVXPy expression.

Element-wise multiplication is DCP-compliant when at least one operand is constant. For quadratic forms, use MatMul instead.

Parameters:

Name Type Description Default
node Mul

Mul expression node with multiple factors

required

Returns:

Type Description
Expression

CVXPy expression representing element-wise product

Note

For convex optimization, typically one factor should be constant. CVXPy will raise a DCP error if the composition violates DCP rules.

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Mul)
def _visit_mul(self, node: Mul) -> cp.Expression:
    """Lower element-wise multiplication to CVXPy expression.

    Element-wise multiplication is DCP-compliant when at least one operand
    is constant. For quadratic forms, use MatMul instead.

    Args:
        node: Mul expression node with multiple factors

    Returns:
        CVXPy expression representing element-wise product

    Note:
        For convex optimization, typically one factor should be constant.
        CVXPy will raise a DCP error if the composition violates DCP rules.
    """
    factors = [self.lower(factor) for factor in node.factors]
    result = factors[0]
    for factor in factors[1:]:
        result = result * factor
    return result
_visit_neg(node: Neg) -> cp.Expression

Lower negation (unary minus) to CVXPy expression.

Negation preserves DCP properties (negating convex gives concave).

Parameters:

Name Type Description Default
node Neg

Neg expression node

required

Returns:

Type Description
Expression

CVXPy expression representing -operand

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Neg)
def _visit_neg(self, node: Neg) -> cp.Expression:
    """Lower negation (unary minus) to CVXPy expression.

    Negation preserves DCP properties (negating convex gives concave).

    Args:
        node: Neg expression node

    Returns:
        CVXPy expression representing -operand
    """
    operand = self.lower(node.operand)
    return -operand
_visit_node_reference(node: NodeReference) -> cp.Expression

Lower NodeReference - extract value at a specific trajectory node.

NodeReference enables cross-node constraints by referencing state/control values at specific discrete time points. This requires the variable_map to contain full trajectory arrays (N, n_x) or (N, n_u) rather than single-node vectors.

Parameters:

Name Type Description Default
node NodeReference

NodeReference expression with base and node_idx

required

Returns:

Type Description
Expression

CVXPy expression representing the variable at the specified node:

Expression

x[node_idx, slice] or u[node_idx, slice]

Raises:

Type Description
ValueError

If the required trajectory variable is not in variable_map

ValueError

If the base variable has no slice assigned

NotImplementedError

If the base is a compound expression

Example

For cross-node constraint: position.at(5) - position.at(4) <= 0.1

variable_map = { "x": cp.vstack([x_nonscaled[k] for k in range(N)]), # (N, n_x) }

position.at(5) lowers to x[5, position._slice]
Note

The node_idx is already resolved to an absolute integer index during expression construction, so negative indices are already handled.

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(NodeReference)
def _visit_node_reference(self, node: "NodeReference") -> cp.Expression:
    """Lower NodeReference - extract value at a specific trajectory node.

    NodeReference enables cross-node constraints by referencing state/control
    values at specific discrete time points. This requires the variable_map to
    contain full trajectory arrays (N, n_x) or (N, n_u) rather than single-node
    vectors.

    Args:
        node: NodeReference expression with base and node_idx

    Returns:
        CVXPy expression representing the variable at the specified node:
        x[node_idx, slice] or u[node_idx, slice]

    Raises:
        ValueError: If the required trajectory variable is not in variable_map
        ValueError: If the base variable has no slice assigned
        NotImplementedError: If the base is a compound expression

    Example:
        For cross-node constraint: position.at(5) - position.at(4) <= 0.1

        variable_map = {
            "x": cp.vstack([x_nonscaled[k] for k in range(N)]),  # (N, n_x)
        }
        # position.at(5) lowers to x[5, position._slice]

    Note:
        The node_idx is already resolved to an absolute integer index during
        expression construction, so negative indices are already handled.
    """
    from openscvx.symbolic.expr.control import Control
    from openscvx.symbolic.expr.state import State

    idx = node.node_idx

    if isinstance(node.base, State):
        if "x" not in self.variable_map:
            raise ValueError(
                "State vector 'x' not found in variable_map. "
                "For cross-node constraints, 'x' must be the full trajectory (N, n_x)."
            )

        cvx_var = self.variable_map["x"]  # Should be (N, n_x) for cross-node constraints

        # Apply slice if state has one assigned
        if node.base._slice is not None:
            return cvx_var[idx, node.base._slice]
        else:
            # No slice means this is the entire unified state vector
            return cvx_var[idx, :]

    elif isinstance(node.base, Control):
        if "u" not in self.variable_map:
            raise ValueError(
                "Control vector 'u' not found in variable_map. "
                "For cross-node constraints, 'u' must be the full trajectory (N, n_u)."
            )

        cvx_var = self.variable_map["u"]  # Should be (N, n_u) for cross-node constraints

        # Apply slice if control has one assigned
        if node.base._slice is not None:
            return cvx_var[idx, node.base._slice]
        else:
            # No slice means this is the entire unified control vector
            return cvx_var[idx, :]

    else:
        # Compound expression (e.g., position[0].at(5))
        # This is more complex - would need to lower base in single-node context
        raise NotImplementedError(
            "Compound expressions in NodeReference are not yet supported for CVXPy lowering. "
            f"Base expression type: {type(node.base).__name__}. "
            "Only State and Control NodeReferences are currently supported."
        )
_visit_norm(node: Norm) -> cp.Expression

Lower norm operation to CVXPy expression.

Norms are convex functions and commonly used in convex optimization. Supports all CVXPy norm types (1, 2, inf, "fro", etc.).

Parameters:

Name Type Description Default
node Norm

Norm expression node with ord attribute

required

Returns:

Type Description
Expression

CVXPy expression representing the norm of the operand

Note

Common norms: ord=2 (Euclidean), ord=1 (Manhattan), ord="inf"

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Norm)
def _visit_norm(self, node: Norm) -> cp.Expression:
    """Lower norm operation to CVXPy expression.

    Norms are convex functions and commonly used in convex optimization.
    Supports all CVXPy norm types (1, 2, inf, "fro", etc.).

    Args:
        node: Norm expression node with ord attribute

    Returns:
        CVXPy expression representing the norm of the operand

    Note:
        Common norms: ord=2 (Euclidean), ord=1 (Manhattan), ord="inf"
    """
    operand = self.lower(node.operand)
    return cp.norm(operand, node.ord)
_visit_parameter(node: Parameter) -> cp.Expression

Lower a parameter to a CVXPy expression.

Parameters are looked up by name in the variable_map. They can be mapped to CVXPy Parameter objects (for efficient parameter sweeps) or constants.

Parameters:

Name Type Description Default
node Parameter

Parameter expression node

required

Returns:

Type Description
Expression

CVXPy expression from variable_map (Parameter or constant)

Raises:

Type Description
ValueError

If parameter name is not found in variable_map

Note

For parameter sweeps without recompilation, map to cp.Parameter. For fixed values, map to cp.Constant or numpy arrays.

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Parameter)
def _visit_parameter(self, node: Parameter) -> cp.Expression:
    """Lower a parameter to a CVXPy expression.

    Parameters are looked up by name in the variable_map. They can be mapped
    to CVXPy Parameter objects (for efficient parameter sweeps) or constants.

    Args:
        node: Parameter expression node

    Returns:
        CVXPy expression from variable_map (Parameter or constant)

    Raises:
        ValueError: If parameter name is not found in variable_map

    Note:
        For parameter sweeps without recompilation, map to cp.Parameter.
        For fixed values, map to cp.Constant or numpy arrays.
    """
    param_name = node.name
    if param_name in self.variable_map:
        return self.variable_map[param_name]
    else:
        raise ValueError(
            f"Parameter '{param_name}' not found in variable_map. "
            f"Add it during CVXPy lowering or use cp.Parameter for parameter sweeps."
        )
_visit_pos(node: PositivePart) -> cp.Expression

Lower positive part function to CVXPy.

Computes max(x, 0), which is convex. Used in penalty methods for inequality constraints.

Parameters:

Name Type Description Default
node PositivePart

PositivePart expression node

required

Returns:

Type Description
Expression

CVXPy expression representing max(operand, 0)

Note

Positive part is convex and commonly used in hinge loss and penalty methods for inequality constraints.

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(PositivePart)
def _visit_pos(self, node: PositivePart) -> cp.Expression:
    """Lower positive part function to CVXPy.

    Computes max(x, 0), which is convex. Used in penalty methods for
    inequality constraints.

    Args:
        node: PositivePart expression node

    Returns:
        CVXPy expression representing max(operand, 0)

    Note:
        Positive part is convex and commonly used in hinge loss and
        penalty methods for inequality constraints.
    """
    operand = self.lower(node.x)
    return cp.maximum(operand, 0.0)
_visit_power(node: Power) -> cp.Expression

Lower element-wise power (base**exponent) to CVXPy expression.

Power is DCP-compliant for specific exponent values: - exponent >= 1: convex (when base >= 0) - 0 <= exponent <= 1: concave (when base >= 0)

Parameters:

Name Type Description Default
node Power

Power expression node

required

Returns:

Type Description
Expression

CVXPy expression representing base**exponent

Note

CVXPy will verify DCP compliance at problem construction time. Common convex cases: x^2, x^3, x^4 (even powers)

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Power)
def _visit_power(self, node: Power) -> cp.Expression:
    """Lower element-wise power (base**exponent) to CVXPy expression.

    Power is DCP-compliant for specific exponent values:
    - exponent >= 1: convex (when base >= 0)
    - 0 <= exponent <= 1: concave (when base >= 0)

    Args:
        node: Power expression node

    Returns:
        CVXPy expression representing base**exponent

    Note:
        CVXPy will verify DCP compliance at problem construction time.
        Common convex cases: x^2, x^3, x^4 (even powers)
    """
    base = self.lower(node.base)
    exponent = self.lower(node.exponent)
    return cp.power(base, exponent)
_visit_sin(node: Sin) -> cp.Expression

Raise NotImplementedError for sine function.

Sine is not DCP-compliant in CVXPy as it is neither convex nor concave.

Parameters:

Name Type Description Default
node Sin

Sin expression node

required

Raises:

Type Description
NotImplementedError

Always raised since sine is not DCP-compliant

Note

For constraints involving trigonometric functions: - Use piecewise-linear approximations, or - Handle in the JAX dynamics/constraint layer instead of CVXPy

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Sin)
def _visit_sin(self, node: Sin) -> cp.Expression:
    """Raise NotImplementedError for sine function.

    Sine is not DCP-compliant in CVXPy as it is neither convex nor concave.

    Args:
        node: Sin expression node

    Raises:
        NotImplementedError: Always raised since sine is not DCP-compliant

    Note:
        For constraints involving trigonometric functions:
        - Use piecewise-linear approximations, or
        - Handle in the JAX dynamics/constraint layer instead of CVXPy
    """
    raise NotImplementedError(
        "Trigonometric functions like Sin are not DCP-compliant in CVXPy. "
        "Consider using piecewise-linear approximations or handle these constraints "
        "in the dynamics (JAX) layer instead."
    )
_visit_sqrt(node: Sqrt) -> cp.Expression

Lower square root to CVXPy expression.

Square root is concave and DCP-compliant when used appropriately (e.g., maximizing sqrt(x) or constraints like sqrt(x) >= c).

Parameters:

Name Type Description Default
node Sqrt

Sqrt expression node

required

Returns:

Type Description
Expression

CVXPy expression representing sqrt(operand)

Note

Square root is concave increasing for x > 0. Valid in: - Objective: maximize sqrt(x) - Constraints: sqrt(x) >= c (concave constraint)

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Sqrt)
def _visit_sqrt(self, node: Sqrt) -> cp.Expression:
    """Lower square root to CVXPy expression.

    Square root is concave and DCP-compliant when used appropriately
    (e.g., maximizing sqrt(x) or constraints like sqrt(x) >= c).

    Args:
        node: Sqrt expression node

    Returns:
        CVXPy expression representing sqrt(operand)

    Note:
        Square root is concave increasing for x > 0. Valid in:
        - Objective: maximize sqrt(x)
        - Constraints: sqrt(x) >= c (concave constraint)
    """
    operand = self.lower(node.operand)
    return cp.sqrt(operand)
_visit_square(node: Square) -> cp.Expression

Lower square function to CVXPy.

Computes x^2, which is convex. Used in quadratic penalty methods and least-squares objectives.

Parameters:

Name Type Description Default
node Square

Square expression node

required

Returns:

Type Description
Expression

CVXPy expression representing operand^2

Note

Square is convex increasing for x >= 0 and convex decreasing for x <= 0. It's always convex overall.

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Square)
def _visit_square(self, node: Square) -> cp.Expression:
    """Lower square function to CVXPy.

    Computes x^2, which is convex. Used in quadratic penalty methods
    and least-squares objectives.

    Args:
        node: Square expression node

    Returns:
        CVXPy expression representing operand^2

    Note:
        Square is convex increasing for x >= 0 and convex decreasing for
        x <= 0. It's always convex overall.
    """
    operand = self.lower(node.x)
    return cp.square(operand)
_visit_srelu(node: SmoothReLU) -> cp.Expression

Lower smooth ReLU penalty function to CVXPy.

Smooth approximation to ReLU: sqrt(max(x, 0)^2 + c^2) - c Differentiable everywhere, approaches ReLU as c -> 0. Convex.

Parameters:

Name Type Description Default
node SmoothReLU

SmoothReLU expression node with smoothing parameter c

required

Returns:

Type Description
Expression

CVXPy expression representing smooth ReLU penalty

Note

This provides a smooth, convex approximation to the ReLU function max(x, 0). The parameter c controls the smoothness: smaller c gives a better approximation but less smoothness.

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(SmoothReLU)
def _visit_srelu(self, node: SmoothReLU) -> cp.Expression:
    """Lower smooth ReLU penalty function to CVXPy.

    Smooth approximation to ReLU: sqrt(max(x, 0)^2 + c^2) - c
    Differentiable everywhere, approaches ReLU as c -> 0. Convex.

    Args:
        node: SmoothReLU expression node with smoothing parameter c

    Returns:
        CVXPy expression representing smooth ReLU penalty

    Note:
        This provides a smooth, convex approximation to the ReLU function
        max(x, 0). The parameter c controls the smoothness: smaller c gives
        a better approximation but less smoothness.
    """
    operand = self.lower(node.x)
    c = node.c
    # smooth_relu(x) = sqrt(max(x, 0)^2 + c^2) - c
    pos_part = cp.maximum(operand, 0.0)
    # For SmoothReLU, we use the 2-norm formulation
    return cp.sqrt(cp.sum_squares(pos_part) + c**2) - c
_visit_stack(node: Stack) -> cp.Expression

Lower vertical stacking to CVXPy expression.

Stacks expressions vertically using cp.vstack. Preserves DCP properties.

Parameters:

Name Type Description Default
node Stack

Stack expression node with multiple rows

required

Returns:

Type Description
Expression

CVXPy expression representing vertical stack of rows

Note

Each row is stacked along axis 0 to create a 2D array.

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Stack)
def _visit_stack(self, node: Stack) -> cp.Expression:
    """Lower vertical stacking to CVXPy expression.

    Stacks expressions vertically using cp.vstack. Preserves DCP properties.

    Args:
        node: Stack expression node with multiple rows

    Returns:
        CVXPy expression representing vertical stack of rows

    Note:
        Each row is stacked along axis 0 to create a 2D array.
    """
    rows = [self.lower(row) for row in node.rows]
    # Stack rows vertically
    return cp.vstack(rows)
_visit_state(node: State) -> cp.Expression

Lower a state variable to a CVXPy expression.

Extracts the appropriate slice from the unified state vector "x" using the slice assigned during unification. The "x" variable must exist in the variable_map.

Parameters:

Name Type Description Default
node State

State expression node

required

Returns:

Type Description
Expression

CVXPy expression representing the state slice: x[slice]

Raises:

Type Description
ValueError

If "x" is not found in variable_map

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(State)
def _visit_state(self, node: State) -> cp.Expression:
    """Lower a state variable to a CVXPy expression.

    Extracts the appropriate slice from the unified state vector "x" using
    the slice assigned during unification. The "x" variable must exist in
    the variable_map.

    Args:
        node: State expression node

    Returns:
        CVXPy expression representing the state slice: x[slice]

    Raises:
        ValueError: If "x" is not found in variable_map
    """
    if "x" not in self.variable_map:
        raise ValueError("State vector 'x' not found in variable_map.")

    cvx_var = self.variable_map["x"]

    # If the state has a slice assigned, apply it
    if node._slice is not None:
        return cvx_var[node._slice]
    return cvx_var
_visit_sub(node: Sub) -> cp.Expression

Lower subtraction to CVXPy expression (element-wise left - right).

Subtraction is affine and always DCP-compliant.

Parameters:

Name Type Description Default
node Sub

Sub expression node

required

Returns:

Type Description
Expression

CVXPy expression representing left - right

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Sub)
def _visit_sub(self, node: Sub) -> cp.Expression:
    """Lower subtraction to CVXPy expression (element-wise left - right).

    Subtraction is affine and always DCP-compliant.

    Args:
        node: Sub expression node

    Returns:
        CVXPy expression representing left - right
    """
    left = self.lower(node.left)
    right = self.lower(node.right)
    return left - right
_visit_sum(node: Sum) -> cp.Expression

Lower sum reduction to CVXPy expression (sums all elements).

Sum preserves DCP properties (sum of convex is convex).

Parameters:

Name Type Description Default
node Sum

Sum expression node

required

Returns:

Type Description
Expression

CVXPy scalar expression representing the sum of all elements

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Sum)
def _visit_sum(self, node: Sum) -> cp.Expression:
    """Lower sum reduction to CVXPy expression (sums all elements).

    Sum preserves DCP properties (sum of convex is convex).

    Args:
        node: Sum expression node

    Returns:
        CVXPy scalar expression representing the sum of all elements
    """
    operand = self.lower(node.operand)
    return cp.sum(operand)
_visit_tan(node: Tan) -> cp.Expression

Raise NotImplementedError for tangent function.

Tangent is not DCP-compliant in CVXPy as it is neither convex nor concave.

Parameters:

Name Type Description Default
node Tan

Tan expression node

required

Raises:

Type Description
NotImplementedError

Always raised since tangent is not DCP-compliant

Note

For constraints involving trigonometric functions: - Use piecewise-linear approximations, or - Handle in the JAX dynamics/constraint layer instead of CVXPy

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Tan)
def _visit_tan(self, node: Tan) -> cp.Expression:
    """Raise NotImplementedError for tangent function.

    Tangent is not DCP-compliant in CVXPy as it is neither convex nor concave.

    Args:
        node: Tan expression node

    Raises:
        NotImplementedError: Always raised since tangent is not DCP-compliant

    Note:
        For constraints involving trigonometric functions:
        - Use piecewise-linear approximations, or
        - Handle in the JAX dynamics/constraint layer instead of CVXPy
    """
    raise NotImplementedError(
        "Trigonometric functions like Tan are not DCP-compliant in CVXPy. "
        "Consider using piecewise-linear approximations or handle these constraints "
        "in the dynamics (JAX) layer instead."
    )
_visit_transpose(node: Transpose) -> cp.Expression

Lower matrix transpose to CVXPy expression.

Transpose preserves DCP properties (transpose of convex is convex).

Parameters:

Name Type Description Default
node Transpose

Transpose expression node

required

Returns:

Type Description
Expression

CVXPy expression representing operand.T

Source code in openscvx/symbolic/lowerers/cvxpy.py
@visitor(Transpose)
def _visit_transpose(self, node: Transpose) -> cp.Expression:
    """Lower matrix transpose to CVXPy expression.

    Transpose preserves DCP properties (transpose of convex is convex).

    Args:
        node: Transpose expression node

    Returns:
        CVXPy expression representing operand.T
    """
    operand = self.lower(node.operand)
    return operand.T
lower(expr: Expr) -> cp.Expression

Lower a symbolic expression to a CVXPy expression.

Main entry point for lowering. Delegates to dispatch() which looks up the appropriate visitor method based on the expression type.

Parameters:

Name Type Description Default
expr Expr

Symbolic expression to lower (any Expr subclass)

required

Returns:

Type Description
Expression

CVXPy expression or constraint object. For arithmetic expressions,

Expression

returns cp.Expression. For Equality/Inequality, returns cp.Constraint.

Raises:

Type Description
NotImplementedError

If no visitor exists for the expression type

ValueError

If required variables are not in variable_map

Example

Lower an expression to a CVXPy expression:

lowerer = CvxpyLowerer(variable_map={"x": cvx_x, "u": cvx_u})
x = ox.State("x", shape=(3,))
expr = ox.Norm(x)
cvx_expr = lowerer.lower(expr)
Source code in openscvx/symbolic/lowerers/cvxpy.py
def lower(self, expr: Expr) -> cp.Expression:
    """Lower a symbolic expression to a CVXPy expression.

    Main entry point for lowering. Delegates to dispatch() which looks up
    the appropriate visitor method based on the expression type.

    Args:
        expr: Symbolic expression to lower (any Expr subclass)

    Returns:
        CVXPy expression or constraint object. For arithmetic expressions,
        returns cp.Expression. For Equality/Inequality, returns cp.Constraint.

    Raises:
        NotImplementedError: If no visitor exists for the expression type
        ValueError: If required variables are not in variable_map

    Example:
        Lower an expression to a CVXPy expression:

            lowerer = CvxpyLowerer(variable_map={"x": cvx_x, "u": cvx_u})
            x = ox.State("x", shape=(3,))
            expr = ox.Norm(x)
            cvx_expr = lowerer.lower(expr)
    """
    return dispatch(self, expr)
register_variable(name: str, cvx_expr: cp.Expression)

Register a CVXPy variable/expression for use in lowering.

Adds or updates a variable in the variable_map. Useful for dynamically adding variables after the lowerer has been created.

Parameters:

Name Type Description Default
name str

Variable name (e.g., "x", "u", or parameter name)

required
cvx_expr Expression

CVXPy expression to associate with the name

required
Example

Register a variable:

lowerer = CvxpyLowerer()
lowerer.register_variable("x", cp.Variable(3))
lowerer.register_variable("obs_center", cp.Parameter(3))
Source code in openscvx/symbolic/lowerers/cvxpy.py
def register_variable(self, name: str, cvx_expr: cp.Expression):
    """Register a CVXPy variable/expression for use in lowering.

    Adds or updates a variable in the variable_map. Useful for dynamically
    adding variables after the lowerer has been created.

    Args:
        name: Variable name (e.g., "x", "u", or parameter name)
        cvx_expr: CVXPy expression to associate with the name

    Example:
        Register a variable:

            lowerer = CvxpyLowerer()
            lowerer.register_variable("x", cp.Variable(3))
            lowerer.register_variable("obs_center", cp.Parameter(3))
    """
    self.variable_map[name] = cvx_expr

dispatch(lowerer: Any, expr: Expr)

Dispatch an expression to its registered visitor function.

Looks up the visitor function for the expression's type and calls it. This is the core of the visitor pattern implementation.

Parameters:

Name Type Description Default
lowerer Any

The CvxpyLowerer instance (provides context for visitor methods)

required
expr Expr

The expression node to lower

required

Returns:

Type Description

The result of calling the visitor function (CVXPy expression or constraint)

Raises:

Type Description
NotImplementedError

If no visitor is registered for the expression type

Example

Dispatch an expression to lower it:

lowerer = CvxpyLowerer(variable_map={...})
expr = Add(x, y)
cvx_expr = dispatch(lowerer, expr)  # Calls visit_add
Source code in openscvx/symbolic/lowerers/cvxpy.py
def dispatch(lowerer: Any, expr: Expr):
    """Dispatch an expression to its registered visitor function.

    Looks up the visitor function for the expression's type and calls it.
    This is the core of the visitor pattern implementation.

    Args:
        lowerer: The CvxpyLowerer instance (provides context for visitor methods)
        expr: The expression node to lower

    Returns:
        The result of calling the visitor function (CVXPy expression or constraint)

    Raises:
        NotImplementedError: If no visitor is registered for the expression type

    Example:
        Dispatch an expression to lower it:

            lowerer = CvxpyLowerer(variable_map={...})
            expr = Add(x, y)
            cvx_expr = dispatch(lowerer, expr)  # Calls visit_add
    """
    fn = _CVXPY_VISITORS.get(type(expr))
    if fn is None:
        raise NotImplementedError(
            f"{lowerer.__class__.__name__!r} has no visitor for {type(expr).__name__}"
        )
    return fn(lowerer, expr)

lower_to_cvxpy(expr: Expr, variable_map: Dict[str, cp.Expression] = None) -> cp.Expression

Lower symbolic expression to CVXPy expression or constraint.

Convenience wrapper that creates a CvxpyLowerer and lowers a single symbolic expression to a CVXPy expression. The result can be used in CVXPy optimization problems.

Parameters:

Name Type Description Default
expr Expr

Symbolic expression to lower (any Expr subclass)

required
variable_map Dict[str, Expression]

Dictionary mapping variable names to CVXPy expressions. Must include "x" for states and "u" for controls. May include parameter names mapped to CVXPy Parameters or constants.

None

Returns:

Type Description
Expression

CVXPy expression for arithmetic expressions (Add, Mul, Norm, etc.)

Expression

or CVXPy constraint for constraint expressions (Equality, Inequality)

Raises:

Type Description
NotImplementedError

If the expression type is not supported (e.g., Sin, Cos, CTCS)

ValueError

If required variables are missing from variable_map

Example

Basic expression lowering::

import cvxpy as cp
import openscvx as ox

# Create CVXPy variables
cvx_x = cp.Variable(3, name="x")
cvx_u = cp.Variable(2, name="u")

# Create symbolic expression
x = ox.State("x", shape=(3,))
u = ox.Control("u", shape=(2,))
expr = ox.Norm(x)**2 + 0.1 * ox.Norm(u)**2

# Lower to CVXPy
cvx_expr = lower_to_cvxpy(expr, {"x": cvx_x, "u": cvx_u})

# Use in optimization problem
prob = cp.Problem(cp.Minimize(cvx_expr))
prob.solve()

Constraint lowering::

# Symbolic constraint
constraint = ox.Norm(x) <= 1.0

# Lower to CVXPy constraint
cvx_constraint = lower_to_cvxpy(constraint, {"x": cvx_x, "u": cvx_u})

# Use in problem
prob = cp.Problem(cp.Minimize(cost), constraints=[cvx_constraint])
See Also
  • CvxpyLowerer: The underlying lowerer class
  • lower_to_jax(): Convenience wrapper for JAX lowering
  • lower_symbolic_expressions(): Main orchestrator in symbolic/lower.py
Source code in openscvx/symbolic/lowerers/cvxpy.py
def lower_to_cvxpy(expr: Expr, variable_map: Dict[str, cp.Expression] = None) -> cp.Expression:
    """Lower symbolic expression to CVXPy expression or constraint.

    Convenience wrapper that creates a CvxpyLowerer and lowers a single
    symbolic expression to a CVXPy expression. The result can be used in
    CVXPy optimization problems.

    Args:
        expr: Symbolic expression to lower (any Expr subclass)
        variable_map: Dictionary mapping variable names to CVXPy expressions.
            Must include "x" for states and "u" for controls. May include
            parameter names mapped to CVXPy Parameters or constants.

    Returns:
        CVXPy expression for arithmetic expressions (Add, Mul, Norm, etc.)
        or CVXPy constraint for constraint expressions (Equality, Inequality)

    Raises:
        NotImplementedError: If the expression type is not supported (e.g., Sin, Cos, CTCS)
        ValueError: If required variables are missing from variable_map

    Example:
        Basic expression lowering::

            import cvxpy as cp
            import openscvx as ox

            # Create CVXPy variables
            cvx_x = cp.Variable(3, name="x")
            cvx_u = cp.Variable(2, name="u")

            # Create symbolic expression
            x = ox.State("x", shape=(3,))
            u = ox.Control("u", shape=(2,))
            expr = ox.Norm(x)**2 + 0.1 * ox.Norm(u)**2

            # Lower to CVXPy
            cvx_expr = lower_to_cvxpy(expr, {"x": cvx_x, "u": cvx_u})

            # Use in optimization problem
            prob = cp.Problem(cp.Minimize(cvx_expr))
            prob.solve()

        Constraint lowering::

            # Symbolic constraint
            constraint = ox.Norm(x) <= 1.0

            # Lower to CVXPy constraint
            cvx_constraint = lower_to_cvxpy(constraint, {"x": cvx_x, "u": cvx_u})

            # Use in problem
            prob = cp.Problem(cp.Minimize(cost), constraints=[cvx_constraint])

    See Also:
        - CvxpyLowerer: The underlying lowerer class
        - lower_to_jax(): Convenience wrapper for JAX lowering
        - lower_symbolic_expressions(): Main orchestrator in symbolic/lower.py
    """
    lowerer = CvxpyLowerer(variable_map)
    return lowerer.lower(expr)

visitor(expr_cls: Type[Expr])

Decorator to register a visitor function for an expression type.

This decorator registers a visitor method to handle a specific expression type during CVXPy lowering. The decorated function is stored in _CVXPY_VISITORS and will be called by dispatch() when lowering that expression type.

Parameters:

Name Type Description Default
expr_cls Type[Expr]

The Expr subclass this visitor handles (e.g., Add, Mul, Norm)

required

Returns:

Type Description

Decorator function that registers the visitor and returns it unchanged

Example

Register a function as the visitor for the Add expression:

@visitor(Add)
def _visit_add(self, node: Add):
    # Lower addition to CVXPy
    ...
Note

Multiple expression types can share a visitor by stacking decorators::

@visitor(Equality)
@visitor(Inequality)
def _visit_constraint(self, node: Constraint):
    # Handle both equality and inequality
    ...
Source code in openscvx/symbolic/lowerers/cvxpy.py
def visitor(expr_cls: Type[Expr]):
    """Decorator to register a visitor function for an expression type.

    This decorator registers a visitor method to handle a specific expression
    type during CVXPy lowering. The decorated function is stored in _CVXPY_VISITORS
    and will be called by dispatch() when lowering that expression type.

    Args:
        expr_cls: The Expr subclass this visitor handles (e.g., Add, Mul, Norm)

    Returns:
        Decorator function that registers the visitor and returns it unchanged

    Example:
        Register a function as the visitor for the Add expression:

            @visitor(Add)
            def _visit_add(self, node: Add):
                # Lower addition to CVXPy
                ...

    Note:
        Multiple expression types can share a visitor by stacking decorators::

            @visitor(Equality)
            @visitor(Inequality)
            def _visit_constraint(self, node: Constraint):
                # Handle both equality and inequality
                ...
    """

    def register(fn: Callable[[Any, Expr], cp.Expression]):
        _CVXPY_VISITORS[expr_cls] = fn
        return fn

    return register