-
Notifications
You must be signed in to change notification settings - Fork 2
/
Copy pathsparse_triplet.c
170 lines (160 loc) · 8.33 KB
/
sparse_triplet.c
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
/** @file */
/**
* @brief Build a triplet-format sparse finite element matrix.
*
* @param[in] cf_local_matrix_f A function which computes local matrices.
* @param[in] mesh A struct describing the 2D finite element mesh.
* @param[in] ref_data A struct containing information about the finite element
* space on the reference element.
* @param[in] convection_function Function pointer for calculating the
* convection field.
* @param[out] matrix The triplet-format matrix.
*
* @retval status Integer of the status of the computation: zero for success,
* 1 for illegal array access, and 2 for incorrect parameter choices.
*/
int cf_build_triplet_matrix(cf_local_matrix_f local_matrix_function,
const cf_mesh_s* mesh,
const cf_ref_arrays_s* ref_data,
cf_convection_f convection_function,
cf_triplet_matrix_s* matrix)
{
cf_local_element_s element_data;
int element_num;
int local_index, global_index, node_index;
int left_basis_num, right_basis_num;
int i, j, status;
double xs[3], ys[3];
double local_matrix[ipow(mesh->num_basis_functions, 2)];
status = 0;
for (element_num = 0; element_num < mesh->num_elements; element_num++) {
cf_check(cf_get_corners(mesh, element_num, xs, ys));
element_data = cf_init_element_data(xs, ys,
ref_data->global_supg_constant);
local_matrix_function(ref_data, &element_data, convection_function,
local_matrix);
for (i = 0; i < mesh->num_basis_functions; i++) {
node_index = CF_INDEX(element_num, i,
mesh->num_elements,
mesh->num_basis_functions);
left_basis_num = mesh->elements[node_index];
for (j = 0; j < mesh->num_basis_functions; j++) {
node_index = CF_INDEX(element_num, j,
mesh->num_elements,
mesh->num_basis_functions);
global_index = element_num*
ipow(mesh->num_basis_functions, 2)
+ mesh->num_basis_functions*i + j;
if (global_index >= matrix->length) {
fprintf(stderr, CF_INDEX_TOO_LARGE,
global_index, matrix->length);
global_index = 0;
status = 1;
printf("index too large!\n");
cf_check(1);
}
local_index = CF_INDEX(i, j,
mesh->num_basis_functions,
mesh->num_basis_functions);
right_basis_num = mesh->elements[node_index];
matrix->values[global_index] =
local_matrix[local_index];
matrix->rows[global_index] = left_basis_num;
matrix->columns[global_index] = right_basis_num;
}
}
}
return status;
fail:
return status;
}
/**
* @brief Build a triplet-format sparse finite element matrix that depends on a
* previous solution.
*
* @param[in] cf_local_matrix_f A function which computes local matrices.
* @param[in] mesh A struct describing the 2D finite element mesh.
* @param[in] ref_data A struct containing information about the finite element
* space on the reference element.
* @param[in] convection_function Function pointer for calculating the
* convection field.
* @param[out] matrix The triplet-format matrix.
*
* @retval status Integer of the status of the computation: zero for success,
* 1 for illegal array access, and 2 for incorrect parameter choices.
*/
int cf_build_triplet_matrix_linearized(cf_local_matrix_f local_matrix_function,
const cf_mesh_s* mesh,
const cf_ref_arrays_s* ref_data,
cf_convection_f convection_function,
const cf_vector_s* solution,
cf_triplet_matrix_s* matrix)
{
cf_local_element_s element_data;
int element_num;
int local_index, global_index, node_index;
int left_basis_num, right_basis_num;
int i, j, status;
double xs[3], ys[3];
double local_matrix[ipow(mesh->num_basis_functions, 2)];
double local_values[ref_data->num_basis_functions];
status = 0;
for (element_num = 0; element_num < mesh->num_elements; element_num++) {
cf_check(cf_get_corners(mesh, element_num, xs, ys));
element_data = cf_init_element_data(xs, ys,
ref_data->global_supg_constant);
/*
* This is the only difference from the linear case: fill in the
* values of the element.
*/
for (i = 0; i < mesh->num_basis_functions; i++) {
node_index = CF_INDEX(element_num, i,
mesh->num_elements,
mesh->num_basis_functions);
left_basis_num = mesh->elements[node_index];
if (left_basis_num >= solution->length) {
fprintf(stderr, CF_INDEX_TOO_LARGE,
global_index, solution->length);
global_index = 0;
status = 1;
cf_check(1);
}
local_values[i] = solution->values[left_basis_num];
}
element_data.values = local_values;
local_matrix_function(ref_data, &element_data, convection_function,
local_matrix);
for (i = 0; i < mesh->num_basis_functions; i++) {
node_index = CF_INDEX(element_num, i,
mesh->num_elements,
mesh->num_basis_functions);
left_basis_num = mesh->elements[node_index];
for (j = 0; j < mesh->num_basis_functions; j++) {
node_index = CF_INDEX(element_num, j,
mesh->num_elements,
mesh->num_basis_functions);
global_index = element_num*
ipow(mesh->num_basis_functions, 2)
+ mesh->num_basis_functions*i + j;
if (global_index >= matrix->length) {
fprintf(stderr, CF_INDEX_TOO_LARGE,
global_index, matrix->length);
global_index = 0;
status = 1;
cf_check(1);
}
local_index = CF_INDEX(i, j,
mesh->num_basis_functions,
mesh->num_basis_functions);
right_basis_num = mesh->elements[node_index];
matrix->values[global_index] =
local_matrix[local_index];
matrix->rows[global_index] = left_basis_num;
matrix->columns[global_index] = right_basis_num;
}
}
}
return status;
fail:
return status;
}