1 /* Copyright 2002-2025 CS GROUP
2 * Licensed to CS GROUP (CS) under one or more
3 * contributor license agreements. See the NOTICE file distributed with
4 * this work for additional information regarding copyright ownership.
5 * CS licenses this file to You under the Apache License, Version 2.0
6 * (the "License"); you may not use this file except in compliance with
7 * the License. You may obtain a copy of the License at
8 *
9 * http://www.apache.org/licenses/LICENSE-2.0
10 *
11 * Unless required by applicable law or agreed to in writing, software
12 * distributed under the License is distributed on an "AS IS" BASIS,
13 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14 * See the License for the specific language governing permissions and
15 * limitations under the License.
16 */
17 package org.orekit.ssa.collision.shorttermencounter.probability.twod;
18
19 import org.hipparchus.CalculusFieldElement;
20 import org.hipparchus.linear.BlockFieldMatrix;
21 import org.hipparchus.linear.BlockRealMatrix;
22 import org.hipparchus.linear.FieldLUDecomposition;
23 import org.hipparchus.linear.FieldMatrix;
24 import org.hipparchus.linear.LUDecomposition;
25 import org.hipparchus.linear.RealMatrix;
26 import org.hipparchus.util.MathArrays;
27 import org.orekit.ssa.metrics.FieldProbabilityOfCollision;
28 import org.orekit.ssa.metrics.ProbabilityOfCollision;
29
30 /**
31 * Abstract class for Alfriend1999 normal and maximised methods as they share lots of similarities.
32 *
33 * @author Vincent Cucchietti
34 * @since 12.0
35 */
36 public abstract class AbstractAlfriend1999 extends AbstractShortTermEncounter2DPOCMethod {
37
38 /**
39 * Default constructor.
40 *
41 * @param name name of the method
42 */
43 protected AbstractAlfriend1999(final String name) {
44 super(name);
45 }
46
47 /** {@inheritDoc} */
48 @Override
49 public ProbabilityOfCollision compute(final double xm, final double ym, final double sigmaX, final double sigmaY,
50 final double radius) {
51 // Reconstruct necessary values from inputs
52 final double squaredMahalanobisDistance =
53 ShortTermEncounter2DDefinition.computeSquaredMahalanobisDistance(xm, ym, sigmaX, sigmaY);
54
55 // Compute covariance matrix determinant
56 final double covarianceMatrixDeterminant = computeCovarianceMatrixDeterminant(sigmaX, sigmaY);
57
58 // Compute probability of collision
59 final double value = computeValue(radius, squaredMahalanobisDistance, covarianceMatrixDeterminant);
60
61 return new ProbabilityOfCollision(value, getName(), this.isAMaximumProbabilityOfCollisionMethod());
62 }
63
64 /** {@inheritDoc} */
65 @Override
66 public <T extends CalculusFieldElement<T>> FieldProbabilityOfCollision<T> compute(final T xm, final T ym, final T sigmaX,
67 final T sigmaY,
68 final T radius) {
69 // Reconstruct necessary values from inputs
70 final T squaredMahalanobisDistance =
71 FieldShortTermEncounter2DDefinition.computeSquaredMahalanobisDistance(xm, ym, sigmaX, sigmaY);
72
73 // Compute covariance matrix determinant
74 final T covarianceMatrixDeterminant = computeCovarianceMatrixDeterminant(sigmaX, sigmaY);
75
76 // Compute probability of collision
77 final T value = computeValue(radius, squaredMahalanobisDistance, covarianceMatrixDeterminant);
78
79 return new FieldProbabilityOfCollision<>(value, getName(), this.isAMaximumProbabilityOfCollisionMethod());
80 }
81
82 /**
83 * Compute the covariance matrix determinant.
84 *
85 * @param sigmaX square root of the x-axis eigen value of the diagonalized combined covariance matrix projected onto the
86 * collision plane (m)
87 * @param sigmaY square root of the y-axis eigen value of the diagonalized combined covariance matrix projected onto the
88 * collision plane (m)
89 *
90 * @return covariance matrix determinant
91 */
92 private double computeCovarianceMatrixDeterminant(final double sigmaX, final double sigmaY) {
93 // Rebuild covariance matrix
94 final RealMatrix covarianceMatrix = new BlockRealMatrix(new double[][] {
95 { sigmaX * sigmaX, 0 },
96 { 0, sigmaY * sigmaY }
97 });
98
99 // Compute determinant
100 return new LUDecomposition(covarianceMatrix).getDeterminant();
101 }
102
103 /**
104 * Compute the covariance matrix determinant.
105 *
106 * @param sigmaX square root of the x-axis eigen value of the diagonalized combined covariance matrix projected onto the
107 * collision plane (m)
108 * @param sigmaY square root of the y-axis eigen value of the diagonalized combined covariance matrix projected onto the
109 * collision plane (m)
110 * @param <T> type of the field elements
111 *
112 * @return covariance matrix determinant
113 */
114 private <T extends CalculusFieldElement<T>> T computeCovarianceMatrixDeterminant(final T sigmaX, final T sigmaY) {
115 // Rebuild covariance matrix
116 final T[][] covarianceMatrixData = MathArrays.buildArray(sigmaX.getField(), 2, 2);
117 covarianceMatrixData[0][0] = sigmaX.multiply(sigmaX);
118 covarianceMatrixData[1][1] = sigmaY.multiply(sigmaY);
119
120 final FieldMatrix<T> covarianceMatrix = new BlockFieldMatrix<>(covarianceMatrixData);
121
122 // Compute determinant
123 return new FieldLUDecomposition<>(covarianceMatrix).getDeterminant();
124 }
125
126 /**
127 * Compute the value of the probability of collision.
128 *
129 * @param radius sum of primary and secondary collision object equivalent sphere radii (m)
130 * @param squaredMahalanobisDistance squared Mahalanobis distance
131 * @param covarianceMatrixDeterminant covariance matrix determinant
132 *
133 * @return value of the probability of collision
134 */
135 abstract double computeValue(double radius, double squaredMahalanobisDistance, double covarianceMatrixDeterminant);
136
137 /**
138 * Compute the value of the probability of collision.
139 *
140 * @param radius sum of primary and secondary collision object equivalent sphere radii (m)
141 * @param squaredMahalanobisDistance squared Mahalanobis distance
142 * @param covarianceMatrixDeterminant covariance matrix determinant
143 * @param <T> type of the field elements
144 *
145 * @return value of the probability of collision
146 */
147 abstract <T extends CalculusFieldElement<T>> T computeValue(T radius, T squaredMahalanobisDistance,
148 T covarianceMatrixDeterminant);
149 }