SERiF 0.0.1a
3+1D Stellar Structure and Evolution
Loading...
Searching...
No Matches
utilities.cpp
Go to the documentation of this file.
1#include "utilities.h"
2#include "mfem.hpp"
3#include <memory>
4
5namespace serif::utilities {
6 mfem::SparseMatrix build_reduced_matrix(
7 const mfem::SparseMatrix& matrix,
8 const mfem::Array<int>& trialEssentialDofs,
9 const mfem::Array<int>& testEssentialDofs
10 ) {
11 int M_orig = matrix.Height();
12 int N_orig = matrix.Width();
13
14 // 1. Create boolean lookup tables for eliminated rows/columns
15 // These tables help quickly check if an original row/column index is eliminated.
16 mfem::Array<bool> row_is_eliminated(M_orig);
17 row_is_eliminated = false; // Initialize all to false (no rows eliminated yet)
18 for (int i = 0; i < testEssentialDofs.Size(); ++i) {
19 int r_idx = testEssentialDofs[i];
20 if (r_idx >= 0 && r_idx < M_orig) { // Check for valid index
21 row_is_eliminated[r_idx] = true;
22 }
23 }
24
25 mfem::Array<bool> col_is_eliminated(N_orig);
26 col_is_eliminated = false; // Initialize all to false (no columns eliminated yet)
27 for (int i = 0; i < trialEssentialDofs.Size(); ++i) {
28 int c_idx = trialEssentialDofs[i];
29 if (c_idx >= 0 && c_idx < N_orig) { // Check for valid index
30 col_is_eliminated[c_idx] = true;
31 }
32 }
33
34 // 2. Create mappings from old (original) indices to new indices
35 // Also, count the number of rows and columns in the new, reduced matrix.
36 mfem::Array<int> old_row_to_new_row(M_orig);
37 int num_new_rows = 0;
38 for (int i = 0; i < M_orig; ++i) {
39 if (row_is_eliminated[i]) {
40 old_row_to_new_row[i] = -1; // Mark as eliminated (no corresponding new row)
41 } else {
42 old_row_to_new_row[i] = num_new_rows++; // Assign the next available new row index
43 }
44 }
45
46 mfem::Array<int> old_col_to_new_col(N_orig);
47 int num_new_cols = 0;
48 for (int i = 0; i < N_orig; ++i) {
49 if (col_is_eliminated[i]) {
50 old_col_to_new_col[i] = -1; // Mark as eliminated (no corresponding new column)
51 } else {
52 old_col_to_new_col[i] = num_new_cols++; // Assign the next available new column index
53 }
54 }
55
56 // 3. Create the new sparse matrix with the calculated reduced dimensions.
57 // It's initially empty (no non-zero entries).
58 mfem::SparseMatrix A_new(num_new_rows, num_new_cols);
59
60 // 4. Iterate through the non-zero entries of the original matrix (A_orig).
61 // A_orig is typically stored in Compressed Sparse Row (CSR) format.
62 // GetI() returns row pointers, GetJ() returns column indices, GetData() returns values.
63 const int* I_orig = matrix.GetI();
64 const int* J_orig = matrix.GetJ();
65 const double* V_orig = matrix.GetData(); // Assuming double, can be templated if needed
66
67 for (int r_orig = 0; r_orig < M_orig; ++r_orig) {
68 // If the original row is marked for elimination, skip all its entries.
69 if (row_is_eliminated[r_orig]) {
70 continue;
71 }
72
73 // Get the new row index for the current original row.
74 int r_new = old_row_to_new_row[r_orig];
75
76 // Iterate through non-zero entries in the current original row r_orig.
77 // I_orig[r_orig] is the start index in J_orig and V_orig for row r_orig.
78 // I_orig[r_orig+1]-1 is the end index.
79 for (int k = I_orig[r_orig]; k < I_orig[r_orig + 1]; ++k) {
80 int c_orig = J_orig[k]; // Original column index of the non-zero entry
81 double val = V_orig[k]; // Value of the non-zero entry
82
83 if (col_is_eliminated[c_orig]) {
84 continue;
85 }
86
87 int c_new = old_col_to_new_col[c_orig];
88
89 A_new.Add(r_new, c_new, val);
90 }
91 }
92
93 A_new.Finalize();
94
95 return A_new;
96 }
97
98 mfem::Vector build_dof_identification_vector(const mfem::Array<int>& allDofs, const::mfem::Array<int>& highlightDofs) {
99 mfem::Vector v(allDofs.Size());
100 v = 0.0; // Initialize the vector to zero
101 v.SetSubVector(highlightDofs, 1.0); // Set the highlighted dofs to 1.0
102 return v;
103 }
104
105 mfem::GridFunction compute_curl(mfem::GridFunction& phi_gf) {
106 mfem::Mesh* mesh = phi_gf.FESpace()->GetMesh();
107 const int dim = mesh->Dimension();
108 // Match the polynomial order of the original RT space for consistency.
109 const int order = phi_gf.FESpace()->GetOrder(0);
110 mfem::Vector curl_mag_vec;
111
112 for (int ne = 0; ne < mesh->GetNE(); ++ne) {
113 if (mesh->GetElementType(ne) != mfem::Element::TRIANGLE &&
114 mesh->GetElementType(ne) != mfem::Element::QUADRILATERAL &&
115 mesh->GetElementType(ne) != mfem::Element::TETRAHEDRON &&
116 mesh->GetElementType(ne) != mfem::Element::HEXAHEDRON) {
117 throw std::invalid_argument("Mesh element type not supported for curl computation.");
118 }
119 mfem::IsoparametricTransformation T;
120 mesh->GetElementTransformation(ne, &T);
121 phi_gf.GetCurl(T, curl_mag_vec);
122 }
123 mfem::L2_FECollection fac(order, dim);
124 mfem::FiniteElementSpace fs(mesh, &fac);
125 mfem::GridFunction curl_gf(&fs);
126 curl_gf = 0.0;
127 return curl_gf;
128
129 }
130
131}
mfem::Vector build_dof_identification_vector(const mfem::Array< int > &allDofs, const::mfem::Array< int > &highlightDofs)
Generate a vector of 1s and 0s where 1 elemetns cooresponds to queried dofs. Useful for degugging.
Definition utilities.cpp:98
mfem::SparseMatrix build_reduced_matrix(const mfem::SparseMatrix &matrix, const mfem::Array< int > &trialEssentialDofs, const mfem::Array< int > &testEssentialDofs)
Definition utilities.cpp:6
mfem::GridFunction compute_curl(mfem::GridFunction &phi_gf)
Computes the curl of a given H(div) grid function (e.g., from an RT space).