Dhairya Malhotra 7 년 전
부모
커밋
165efcc066
1개의 변경된 파일84개의 추가작업 그리고 9개의 파일을 삭제
  1. 84 9
      include/sctl/sph_harm.hpp

+ 84 - 9
include/sctl/sph_harm.hpp

@@ -32,27 +32,99 @@ template <class Real> class SphericalHarmonics{
   public:
 
     // Scalar Spherical Harmonics
-    static void Grid2SHC(const Vector<Real>& X_in, Long Nt_in, Long Np_in, Long p_out, Vector<Real>& S_out, SHCArrange arrange_out);
 
-    static void SHC2Grid(const Vector<Real>& S_in, SHCArrange arrange_in, Long p_in, Long Nt_out, Long Np_out, Vector<Real>* X_out, Vector<Real>* X_theta_out=nullptr, Vector<Real>* X_phi_out=nullptr);
+    /**
+     * \brief Compute spherical harmonic coefficients from grid values.
+     * \param[in] X Grid values {X(t0,p0), X(t0,p1), ... , X(t1,p0), X(t1,p1), ... }, where, {cos(t0), cos(t1), ... } are the Gauss-Legendre nodes of order (Nt-1) in the interval [-1,1] and {p0, p1, ... } are equispaced in [0, 2*pi].
+     * \param[in] Nt Number of grid points \theta \in (1,pi).
+     * \param[in] Np Number of grid points \phi \in (1,2*pi).
+     * \param[in] p Order of spherical harmonic expansion.
+     * \param[in] arrange Arrangement of the coefficients.
+     * \param[out] S Spherical harmonic coefficients.
+     */
+    static void Grid2SHC(const Vector<Real>& X, Long Nt, Long Np, Long p, Vector<Real>& S, SHCArrange arrange);
+
+    /**
+     * \brief Evaluate grid values from spherical harmonic coefficients.
+     * \param[in] S Spherical harmonic coefficients.
+     * \param[in] arrange Arrangement of the coefficients.
+     * \param[in] p Order of spherical harmonic expansion.
+     * \param[in] Nt Number of grid points \theta \in (1,pi).
+     * \param[in] Np Number of grid points \phi \in (1,2*pi).
+     * \param[out] X Grid values {X(t0,p0), X(t0,p1), ... , X(t1,p0), X(t1,p1), ... }, where, {cos(t0), cos(t1), ... } are the Gauss-Legendre nodes of order (Nt-1) in the interval [-1,1] and {p0, p1, ... } are equispaced in [0, 2*pi].
+     * \param[out] X_theta \theta derivative of X evaluated at grid points.
+     * \param[out] X_phi \phi derivative of X evaluated at grid points.
+     */
+    static void SHC2Grid(const Vector<Real>& S, SHCArrange arrange, Long p, Long Nt, Long Np, Vector<Real>* X, Vector<Real>* X_theta=nullptr, Vector<Real>* X_phi=nullptr);
 
-    static void SHCEval(const Vector<Real>& S_in, SHCArrange arrange_in, Long p_in, const Vector<Real>& cos_theta_phi_in, Vector<Real>& X_out);
+    /**
+     * \brief Evaluate point values from spherical harmonic coefficients.
+     * \param[in] S Spherical harmonic coefficients.
+     * \param[in] arrange Arrangement of the coefficients.
+     * \param[in] p Order of spherical harmonic expansion.
+     * \param[in] cos_theta_phi Evaluation coordinates given as {cos(t0),p0, cos(t1),p1, ... }.
+     * \param[out] X Evaluated values {X0, X1, ... }.
+     */
+    static void SHCEval(const Vector<Real>& S, SHCArrange arrange, Long p, const Vector<Real>& cos_theta_phi, Vector<Real>& X);
 
-    static void SHC2Pole(const Vector<Real>& S_in, SHCArrange arrange_in, Long p_in, Vector<Real>& P_out);
+    static void SHC2Pole(const Vector<Real>& S, SHCArrange arrange, Long p, Vector<Real>& P);
 
     static void WriteVTK(const char* fname, const Vector<Real>* S, const Vector<Real>* f_val, SHCArrange arrange, Long p_in, Long p_out, Real period=0, const Comm& comm = Comm::World());
 
 
     // Vector Spherical Harmonics
-    static void Grid2VecSHC(const Vector<Real>& X_in, Long Nt_in, Long Np_in, Long p_out, Vector<Real>& S_out, SHCArrange arrange_out);
 
-    static void VecSHC2Grid(const Vector<Real>& S_in, SHCArrange arrange_in, Long p_in, Long Nt_out, Long Np_out, Vector<Real>& X_out);
+    /**
+     * \brief Compute vector spherical harmonic coefficients from grid values.
+     * \param[in] X Grid values {X(t0,p0), X(t0,p1), ... , X(t1,p0), ... , Y(t0,p0), ... , Z(t0,p0), ... }, where, {cos(t0), cos(t1), ... } are the Gauss-Legendre nodes of order (Nt-1) in the interval [-1,1] and {p0, p1, ... } are equispaced in [0, 2*pi].
+     * \param[in] Nt Number of grid points \theta \in (1,pi).
+     * \param[in] Np Number of grid points \phi \in (1,2*pi).
+     * \param[in] p Order of spherical harmonic expansion.
+     * \param[in] arrange Arrangement of the coefficients.
+     * \param[out] S Vector spherical harmonic coefficients.
+     */
+    static void Grid2VecSHC(const Vector<Real>& X, Long Nt, Long Np, Long p, Vector<Real>& S, SHCArrange arrange);
+
+    /**
+     * \brief Evaluate grid values from vector spherical harmonic coefficients.
+     * \param[in] S Vector spherical harmonic coefficients.
+     * \param[in] arrange Arrangement of the coefficients.
+     * \param[in] p Order of spherical harmonic expansion.
+     * \param[in] Nt Number of grid points \theta \in (1,pi).
+     * \param[in] Np Number of grid points \phi \in (1,2*pi).
+     * \param[out] X Grid values {X(t0,p0), X(t0,p1), ... , X(t1,p0), X(t1,p1), ... , Y(t0,p0), ... , Z(t0,p0), ... }, where, {cos(t0), cos(t1), ... } are the Gauss-Legendre nodes of order (Nt-1) in the interval [-1,1] and {p0, p1, ... } are equispaced in [0, 2*pi].
+     */
+    static void VecSHC2Grid(const Vector<Real>& S, SHCArrange arrange, Long p, Long Nt, Long Np, Vector<Real>& X);
 
-    static void VecSHCEval(const Vector<Real>& S_in, SHCArrange arrange_in, Long p_in, const Vector<Real>& cos_theta_phi_in, Vector<Real>& X_out);
+    /**
+     * \brief Evaluate point values from vector spherical harmonic coefficients.
+     * \param[in] S Vector spherical harmonic coefficients.
+     * \param[in] arrange Arrangement of the coefficients.
+     * \param[in] p Order of spherical harmonic expansion.
+     * \param[in] cos_theta_phi Evaluation coordinates given as {cos(t0),p0, cos(t1),p1, ... }.
+     * \param[out] X Evaluated values {X0,Y0,Z0, X1,Y1,Z1, ... }.
+     */
+    static void VecSHCEval(const Vector<Real>& S, SHCArrange arrange, Long p, const Vector<Real>& cos_theta_phi, Vector<Real>& X);
 
-    static void StokesEvalSL(const Vector<Real>& S_in, SHCArrange arrange_in, Long p_in, const Vector<Real>& coord_in, Vector<Real>& X_out);
+    /**
+     * \brief Evaluate Stokes single-layer operator at point values from the vector spherical harmonic coefficients for the density.
+     * \param[in] S Vector spherical harmonic coefficients.
+     * \param[in] arrange Arrangement of the coefficients.
+     * \param[in] p Order of spherical harmonic expansion.
+     * \param[in] Evaluation coordinates given as {x0,y0,z0, x1,y1,z1, ... }.
+     * \param[out] U Evaluated values {Ux0,Uy0,Uz0, Ux1,Uy1,Uz1, ... }.
+     */
+    static void StokesEvalSL(const Vector<Real>& S, SHCArrange arrange, Long p, const Vector<Real>& coord, Vector<Real>& U);
 
-    static void StokesEvalDL(const Vector<Real>& S_in, SHCArrange arrange_in, Long p_in, const Vector<Real>& coord_in, Vector<Real>& X_out);
+    /**
+     * \brief Evaluate Stokes double-layer operator at point values from the vector spherical harmonic coefficients for the density.
+     * \param[in] S Vector spherical harmonic coefficients.
+     * \param[in] arrange Arrangement of the coefficients.
+     * \param[in] p Order of spherical harmonic expansion.
+     * \param[in] Evaluation coordinates given as {x0,y0,z0, x1,y1,z1, ... }.
+     * \param[out] U Evaluated values {Ux0,Uy0,Uz0, Ux1,Uy1,Uz1, ... }.
+     */
+    static void StokesEvalDL(const Vector<Real>& S, SHCArrange arrange, Long p, const Vector<Real>& coord, Vector<Real>& U);
 
 
     static void test_stokes() {
@@ -154,6 +226,9 @@ template <class Real> class SphericalHarmonics{
       Clear();
     }
 
+    /**
+     * \brief Clear all precomputed data. This must be done before the program exits to avoid memory leaks.
+     */
     static void Clear() { MatrixStore().Resize(0); }
 
   private: