Graphite  0.5.0
GPU-accelerated graph optimization framework
Loading...
Searching...
No Matches
hessian.hpp
Go to the documentation of this file.
1
2#pragma once
4
5namespace graphite {
6
7namespace ops {
8
9template <typename T, typename S, size_t N, size_t E>
10__global__ void compute_hessian_block_kernel(
11 const size_t vi, const size_t vj, size_t dim_i, size_t dim_j,
12 const size_t *active_factors, const size_t num_active_factors,
13 const size_t *ids, const size_t *block_offsets, const uint8_t *vi_active,
14 const uint8_t *vj_active, const size_t *hessian_offset_i,
15 const size_t *hessian_offset_j, const S *jacobian_i, const S *jacobian_j,
16 const S *precision, const S *chi2_derivative, S *hessian) {
17 // TODO: simpify and optimize this kernel
18 const auto idx = get_thread_id();
19
20 const auto block_id = idx / (dim_i * dim_j);
21
22 if (block_id >= num_active_factors) {
23 return;
24 }
25
26 const auto factor_idx = active_factors[block_id];
27
28 const size_t vi_id = ids[factor_idx * N + vi];
29 const size_t vj_id = ids[factor_idx * N + vj];
30
31 if (is_vertex_active(vi_active, vi_id) &&
32 is_vertex_active(vj_active, vj_id)) {
33
34 const size_t block_size = dim_i * dim_j;
35 const size_t offset = idx % block_size;
36 // Hessian block may be rectangular
37 // output blocks are all column major
38
39 const bool transposed = hessian_offset_i[vi_id] > hessian_offset_j[vj_id];
40
41 auto ji = factor_idx * E * dim_i + jacobian_i;
42 auto jj = factor_idx * E * dim_j + jacobian_j;
43 const auto precision_offset = factor_idx * E * E;
44 const auto p = precision + precision_offset;
45
46 if (transposed) {
47 cuda::std::swap(ji, jj);
48 cuda::std::swap(dim_i, dim_j);
49 }
50
51 const size_t row = offset % dim_i;
52 const size_t col = offset / dim_i;
53
54 // Each thread computes one element of the Hessian block
55 using highp = T;
56 highp value = 0;
57
58 // computes J_i^T * P * J_j
59
60 const auto J = jj + col * E;
61 const auto Jt = ji + row * E;
62#pragma unroll
63 for (int i = 0; i < E; i++) { // p row
64 highp pj = 0;
65#pragma unroll
66 for (int j = 0; j < E; j++) { // p col
67 pj += (highp)p[i * E + j] * (highp)J[j];
68 }
69 value += (highp)Jt[i] * pj;
70 }
71
72 value *= (highp)chi2_derivative[factor_idx];
73
74 auto block = hessian + (block_offsets[block_id] + (row + col * dim_i));
75 S lp_value = static_cast<S>(value);
76 atomicAdd(block, lp_value);
77 }
78}
79
80template <typename S, int D>
81__global__ void augment_hessian_diagonal_kernel(S *diagonal_blocks,
82 S *scalar_diagonal, const S mu,
83 const bool use_identity,
84 const uint8_t *active_state,
85 const size_t num_threads) {
86 const size_t idx = get_thread_id();
87
88 if (idx >= num_threads) {
89 return;
90 }
91
92 constexpr auto block_size = D * D;
93
94 const auto vertex_id = idx;
95 if (!is_vertex_active(active_state, vertex_id)) {
96 return;
97 }
98
99 S *block = diagonal_blocks + vertex_id * block_size;
100 for (size_t i = 0; i < D; i++) {
101 const double diag = static_cast<double>(scalar_diagonal[vertex_id * D + i]);
102 double new_diag = diag;
103 if (use_identity) {
104 new_diag += static_cast<double>(mu);
105 } else {
106 new_diag += static_cast<double>(mu) * std::clamp(diag, 1.0e-6, 1.0e32);
107 }
108 block[i * D + i] = static_cast<S>(new_diag);
109 }
110}
111
112template <typename T, typename S, typename V>
113void augment_block_diagonal(V *v, InvP<T, S> *block_diagonal,
114 InvP<T, S> *scalar_diagonal, const T mu,
115 const bool use_identity, cudaStream_t stream) {
116 const size_t num_threads = v->count();
117 const auto threads_per_block = 256;
118 const auto num_blocks =
119 (num_threads + threads_per_block - 1) / threads_per_block;
120
121 augment_hessian_diagonal_kernel<InvP<T, S>, V::dim>
122 <<<num_blocks, threads_per_block, 0, stream>>>(
123 block_diagonal, scalar_diagonal, (InvP<T, S>)mu, use_identity,
124 v->get_active_state(), num_threads);
125}
126
127template <typename T, typename S, int D>
128__global__ void apply_block_jacobi_kernel(T *z, const T *r, S *block_diagonal,
129 const size_t *hessian_ids,
130 const uint8_t *active_state,
131 const size_t num_threads) {
132 const size_t idx = get_thread_id();
133 const auto local_vertex_id = idx / D;
134
135 if (idx >= num_threads || !is_vertex_active(active_state, local_vertex_id)) {
136 return;
137 }
138
139 constexpr auto block_size = D * D;
140
141 S *block = block_diagonal + local_vertex_id * block_size;
142 const auto hessian_offset = hessian_ids[local_vertex_id];
143 const auto offset = idx % D;
144 const auto row = offset;
145
146 T value = 0;
147#pragma unroll
148 for (size_t i = 0; i < D; i++) {
149 value += (T)block[row + i * D] * r[hessian_offset + i];
150 }
151 z[hessian_offset + row] = value;
152}
153
154template <typename T, typename S, typename V>
155void apply_block_jacobi(V *v, T *z, const T *r, InvP<T, S> *block_diagonal,
156 cudaStream_t stream) {
157 const size_t num_parameters = v->count() * v->dimension();
158 const size_t num_threads = num_parameters;
159 const auto threads_per_block = 256;
160 const auto num_blocks =
161 (num_threads + threads_per_block - 1) / threads_per_block;
162
163 apply_block_jacobi_kernel<T, InvP<T, S>, V::dim>
164 <<<num_blocks, threads_per_block, 0, stream>>>(
165 z, r, block_diagonal, v->get_hessian_ids(), v->get_active_state(),
166 num_threads);
167}
168
169template <typename highp, typename InvP, typename T, size_t I, size_t N,
170 size_t E, size_t D>
171__global__ void compute_hessian_diagonal_kernel(
172 InvP *diagonal_blocks, const T *jacs, const size_t *active_ids,
173 const size_t *ids, const uint8_t *active_state, const T *pmat,
174 const T *chi2_derivative, const size_t num_threads) {
175 const size_t idx = get_thread_id();
176
177 if (idx >= num_threads) {
178 return;
179 }
180
181 constexpr size_t jacobian_size = D * E;
182 constexpr size_t block_size = D * D;
183
184 // Stored as E x d col major, but we need to transpose it to d x E, where d is
185 // the vertex size
186 const size_t factor_id = active_ids[idx / block_size];
187 const size_t local_id =
188 ids[factor_id * N +
189 I]; // N is the number of vertices involved in the factor
190 if (!is_vertex_active(active_state, local_id)) {
191 return;
192 }
193 const auto jacobian_offset = factor_id * jacobian_size;
194
195 constexpr auto precision_matrix_size = E * E;
196 const auto precision_offset = factor_id * precision_matrix_size;
197
198 // Identify H block row and column (column major)
199 // const size_t row = idx % D;
200 // const size_t col = idx / D;
201
202 const size_t offset = idx % block_size;
203 const size_t row = offset % D;
204 const size_t col = offset / D;
205
206 // left[i]*pmat[i*E+j]*right[i] = h value
207 // where i goes from 0 to E
208 const T *Jt = jacs + jacobian_offset + row * E;
209 const T *J = jacs + jacobian_offset + col * E;
210
211 const T *precision_matrix = pmat + precision_offset;
212
213 highp value = 0;
214#pragma unroll
215 for (int i = 0; i < E; i++) { // pmat row
216 highp pj = 0;
217#pragma unroll
218 for (int j = 0; j < E; j++) { // pmat col
219 pj += (highp)precision_matrix[i * E + j] * (highp)J[j];
220 }
221 value += (highp)Jt[i] * pj;
222 // value += Jt[i]*J[i];
223 }
224
225 // if (row == col) {
226 // value = 1;
227 // // printf("Thread %d, row: %d, col: %d, value: %f\n", idx, row, col,
228 // value);
229 // // printf("D=%d\n", D);
230 // }
231 // else {
232 // value = 0;
233 // }
234 value *= (highp)chi2_derivative[factor_id];
235
236 // T* block = diagonal_blocks + local_id*block_size + (idx % block_size);
237 // printf("Thread %d, Hessian offset: %u\n", idx, local_id);
238 // if (value != 0) {
239 // printf("Thread %d, row: %d, col: %d, value: %f\n", idx, row, col, value);
240 // }
241 // if (D == 9 && local_id*block_size + row + col*D > 1701) {
242 // printf("Thread %d, vertex id: %u, row: %d, col: %d, value: %f, offset:
243 // %u \n", idx, local_id, row, col, value, local_id*block_size + row +
244 // col*D);
245
246 // }
247 // value = 5.0;
248 // T* block = diagonal_blocks + local_id*block_size + row + col*D;
249 InvP *block = diagonal_blocks + (local_id * block_size + row + col * D);
250 // if (row == col) {
251 // // local_id = 5;
252 // // row = 6; col = 7;
253 // // value = 8.5;
254 // printf("Thread %d, vertex id: %llu, row: %llu, col: %llu, value: %f,
255 // offset: %u \n", idx, local_id, row, col, value, local_id*block_size +
256 // row + col*D);
257 // }
258 // if (row == 8 && col == 8) {
259 // printf("Number of threads: %d\n", num_threads);
260 // printf("Thread %d, vertex id: %u, row: %zu, col: %zu, value: %f\n",
261 // idx, local_id, row, col, value);
262 // }
263 // T* block = diagonal_blocks + row*D + col;
264 InvP lp_value = static_cast<InvP>(value);
265 atomicAdd(block, lp_value);
266 // block[0] = value;
267 // *block = value;
268}
269
270template <typename highp, typename InvP, typename T, size_t I, size_t N,
271 size_t E, size_t D, typename VT, typename F>
272__global__ void compute_hessian_diagonal_dynamic_kernel(
273 InvP *diagonal_blocks, const size_t *active_ids, const size_t *ids,
274 const size_t *hessian_ids, const VT args,
275 const typename F::ObservationType *obs, const highp *jacobian_scales,
276 const typename F::ConstraintDataType *constraint_data,
277 const uint8_t *active_state, const T *pmat, const T *chi2_derivative,
278 const size_t num_threads) {
279 const size_t idx = get_thread_id();
280
281 if (idx >= num_threads) {
282 return;
283 }
284
285 constexpr size_t jacobian_size = D * E;
286 constexpr size_t block_size = D * D;
287
288 // Stored as E x d col major, but we need to transpose it to d x E, where d is
289 // the vertex size
290 const size_t factor_id = active_ids[idx / block_size];
291 const size_t local_id =
292 ids[factor_id * N +
293 I]; // N is the number of vertices involved in the factor
294 if (!is_vertex_active(active_state, local_id)) {
295 return;
296 }
297
298 constexpr auto precision_matrix_size = E * E;
299 const auto precision_offset = factor_id * precision_matrix_size;
300
301 // Identify H block row and column (column major)
302
303 const size_t offset = idx % block_size;
304 const size_t row = offset % D;
305 const size_t col = offset / D;
306
307 highp jacobian[jacobian_size];
308
309 compute_Jblock<highp, I, N, typename F::ObservationType, E, F, VT>(
310 jacobian, factor_id, local_id, obs, constraint_data, ids, hessian_ids,
311 args, std::make_index_sequence<N>{});
312
313 const auto hessian_offset = hessian_ids[local_id];
314
315 const highp *Jt = jacobian + row * E;
316 const highp *J = jacobian + col * E;
317
318 const T *precision_matrix = pmat + precision_offset;
319
320 const highp *jscale = jacobian_scales + hessian_offset;
321
322 highp value = 0;
323#pragma unroll
324 for (int i = 0; i < E; i++) { // pmat row
325 highp pj = 0;
326#pragma unroll
327 for (int j = 0; j < E; j++) { // pmat col
328 pj += (highp)precision_matrix[i * E + j] * (highp)J[j] * jscale[col];
329 }
330 value += (highp)Jt[i] * jscale[row] * pj;
331 }
332
333 value *= (highp)chi2_derivative[factor_id];
334
335 InvP *block = diagonal_blocks + (local_id * block_size + row + col * D);
336
337 InvP lp_value = static_cast<InvP>(value);
338 atomicAdd(block, lp_value);
339}
340
341template <typename T, typename S, typename F, std::size_t... Is>
342void launch_kernel_block_diagonal(
343 F *f, std::array<InvP<T, S> *, F::get_num_vertices()> &diagonal_blocks,
344 std::array<const size_t *, F::get_num_vertices()> &hessian_ids,
345 std::array<S *, F::get_num_vertices()> &jacs, const T *jacobian_scales,
346 const size_t num_factors, cudaStream_t stream, std::index_sequence<Is...>) {
347 (([&] {
348 constexpr size_t num_vertices = F::get_num_vertices();
349 constexpr size_t dimension = F::get_vertex_sizes()[Is];
350 const size_t num_threads = num_factors * dimension * dimension;
351 // std::cout << "Launching block diagonal kernel" << std::endl;
352 // std::cout << "Num threads: " << num_threads << std::endl;
353 // std::cout << "dimension: " << dimension << std::endl;
354 // std::cout << "num_factors: " << num_factors << std::endl;
355 size_t threads_per_block = 256;
356 size_t num_blocks =
357 (num_threads + threads_per_block - 1) / threads_per_block;
358
359 // std::cout << "Checking obs ptr: " << f->device_obs.data().get() <<
360 // std::endl; std::cout << "Checking residual ptr: " <<
361 // f->residuals.data().get() << std::endl; std::cout << "Checking ids
362 // ptr: " << f->device_ids.data().get() << std::endl;
363
364 if (f->store_jacobians() || !is_analytical<F>()) {
365 compute_hessian_diagonal_kernel<T, InvP<T, S>, S, Is, num_vertices,
366 F::error_dim, dimension>
367 <<<num_blocks, threads_per_block, 0, stream>>>(
368 diagonal_blocks[Is], jacs[Is], f->active_indices.data().get(),
369 f->device_ids.data().get(),
370 f->vertex_descriptors[Is]->get_active_state(),
371 f->precision_matrices.data().get(),
372 f->chi2_derivative.data().get(), num_threads);
373 } else {
374 if constexpr (is_analytical<F>()) {
375 compute_hessian_diagonal_dynamic_kernel<
376 T, InvP<T, S>, S, Is, num_vertices, F::error_dim, dimension,
377 typename F::VertexPointerPointerTuple, F>
378 <<<num_blocks, threads_per_block, 0, stream>>>(
379 diagonal_blocks[Is], f->active_indices.data().get(),
380 f->device_ids.data().get(), hessian_ids[Is], f->get_vertices(),
381 f->device_obs.data().get(), jacobian_scales,
382 f->data.data().get(),
383 f->vertex_descriptors[Is]->get_active_state(),
384 f->precision_matrices.data().get(),
385 f->chi2_derivative.data().get(), num_threads);
386 }
387 }
388 }()),
389 ...);
390}
391
392template <typename T, typename S, typename F>
393void compute_block_diagonal(
394 F *f, std::array<InvP<T, S> *, F::get_num_vertices()> &diagonal_blocks,
395 const T *jacobian_scales, cudaStream_t stream) {
396
397 // Then for each vertex, we need to compute the error
398 constexpr auto num_vertices = F::get_num_vertices();
399 constexpr auto vertex_sizes = F::get_vertex_sizes();
400
401 // At this point all necessary data should be on the GPU
402 // std::array<T*, num_vertices> verts;
403 auto verts = f->get_vertices();
404 std::array<S *, num_vertices> jacs;
405 std::array<const size_t *, num_vertices> hessian_ids;
406 for (int i = 0; i < num_vertices; i++) {
407 jacs[i] = f->jacobians[i].data.data().get();
408 hessian_ids[i] = f->vertex_descriptors[i]->get_hessian_ids();
409 }
410
411 const auto num_factors = f->active_count();
412
413 launch_kernel_block_diagonal<T, S>(f, diagonal_blocks, hessian_ids, jacs,
414 jacobian_scales, num_factors, stream,
415 std::make_index_sequence<num_vertices>{});
416}
417
418template <typename highp, typename T, size_t I, size_t N, size_t E, size_t D>
419__global__ void compute_hessian_scalar_diagonal_kernel(
420 highp *diagonal, const T *jacs, const size_t *active_ids, const size_t *ids,
421 const size_t *hessian_ids, const uint8_t *active_state, const T *pmat,
422 const T *chi2_derivative, const size_t num_threads) {
423 const size_t idx = get_thread_id();
424
425 if (idx >= num_threads) {
426 return;
427 }
428
429 constexpr size_t jacobian_size = D * E;
430
431 // Stored as E x d col major, but we need to transpose it to d x E, where d is
432 // the vertex size
433 const size_t factor_id = active_ids[idx / D];
434 const size_t local_id =
435 ids[factor_id * N +
436 I]; // N is the number of vertices involved in the factor
437 if (!is_vertex_active(active_state, local_id)) {
438 return;
439 }
440 const auto jacobian_offset = factor_id * jacobian_size;
441
442 constexpr auto precision_matrix_size = E * E;
443 const auto precision_offset = factor_id * precision_matrix_size;
444
445 // Identify H block row and column (column major)
446 const size_t row = idx % D;
447 const size_t col = row;
448
449 // left[i]*pmat[i*E+j]*right[i] = h value
450 // where i goes from 0 to E
451 const T *Jt = jacs + jacobian_offset + row * E;
452 const T *J = jacs + jacobian_offset + col * E;
453
454 const T *precision_matrix = pmat + precision_offset;
455
456 highp value = 0;
457#pragma unroll
458 for (int i = 0; i < E; i++) { // pmat row
459 highp pj = 0;
460#pragma unroll
461 for (int j = 0; j < E; j++) { // pmat col
462 pj += (highp)precision_matrix[i * E + j] * (highp)J[j];
463 }
464 value += (highp)Jt[i] * pj;
465 }
466
467 value *= (highp)chi2_derivative[factor_id];
468
469 // T* block = diagonal_blocks+(local_id*block_size + row + col*D);
470 const size_t hessian_offset = hessian_ids[local_id];
471 highp *location = diagonal + hessian_offset + row;
472 // T lp_value = static_cast<T>(value);
473 atomicAdd(location, value);
474}
475
476template <typename highp, typename T, size_t I, size_t N, size_t E, size_t D,
477 typename VT, typename F, bool use_scales>
478__global__ void compute_hessian_scalar_diagonal_dynamic_kernel(
479 highp *diagonal, const T *jacs, const size_t *active_ids, const size_t *ids,
480 const size_t *hessian_ids, const VT args,
481 const typename F::ObservationType *obs, const highp *jacobian_scales,
482 const typename F::ConstraintDataType *constraint_data,
483 const uint8_t *active_state, const T *pmat, const T *chi2_derivative,
484 const size_t num_threads) {
485 const size_t idx = get_thread_id();
486
487 if (idx >= num_threads) {
488 return;
489 }
490
491 constexpr size_t jacobian_size = D * E;
492
493 // Stored as E x d col major, but we need to transpose it to d x E, where d is
494 // the vertex size
495 const size_t factor_id = active_ids[idx / D];
496 const size_t local_id =
497 ids[factor_id * N +
498 I]; // N is the number of vertices involved in the factor
499 if (!is_vertex_active(active_state, local_id)) {
500 return;
501 }
502
503 constexpr auto precision_matrix_size = E * E;
504 const auto precision_offset = factor_id * precision_matrix_size;
505
506 // Identify H block row and column (column major)
507 const size_t row = idx % D;
508 const size_t col = row;
509
510 // using G = std::conditional_t<is_low_precision<T>::value, highp, T>;
511 // G jacobian[jacobian_size];
512 highp jacobian[jacobian_size];
513
514 compute_Jblock<highp, I, N, typename F::ObservationType, E, F, VT>(
515 jacobian, factor_id, local_id, obs, constraint_data, ids, hessian_ids,
516 args, std::make_index_sequence<N>{});
517 const auto hessian_offset = hessian_ids[local_id];
518 const highp *Jt = jacobian + row * E;
519 const highp *J = jacobian + col * E;
520
521 const T *precision_matrix = pmat + precision_offset;
522 highp value = 0;
523 if constexpr (use_scales) {
524 const highp *jscale = jacobian_scales + hessian_offset;
525
526#pragma unroll
527 for (int i = 0; i < E; i++) { // pmat row
528 highp pj = 0;
529#pragma unroll
530 for (int j = 0; j < E; j++) { // pmat col
531 pj += (highp)precision_matrix[i * E + j] * (highp)J[j] * jscale[col];
532 }
533 value += (highp)Jt[i] * jscale[row] * pj;
534 }
535 } else {
536
537#pragma unroll
538 for (int i = 0; i < E; i++) { // pmat row
539 highp pj = 0;
540#pragma unroll
541 for (int j = 0; j < E; j++) { // pmat col
542 pj += (highp)precision_matrix[i * E + j] * (highp)J[j];
543 }
544 value += (highp)Jt[i] * pj;
545 }
546 }
547
548 value *= (highp)chi2_derivative[factor_id];
549
550 // T* block = diagonal_blocks+(local_id*block_size + row + col*D);
551 highp *location = diagonal + hessian_offset + row;
552 // T lp_value = static_cast<T>(value);
553 atomicAdd(location, value);
554}
555
556template <typename T, typename S, typename F, std::size_t... Is>
557void launch_kernel_hessian_scalar_diagonal(
558 F *f, T *diagonal,
559 std::array<const size_t *, F::get_num_vertices()> &hessian_ids,
560 std::array<S *, F::get_num_vertices()> &jacs, const T *jacobian_scales,
561 const size_t num_factors, std::index_sequence<Is...>) {
562 (([&] {
563 constexpr size_t num_vertices = F::get_num_vertices();
564 constexpr size_t dimension = F::get_vertex_sizes()[Is];
565 const size_t num_threads = num_factors * dimension;
566
567 size_t threads_per_block = 256;
568 size_t num_blocks =
569 (num_threads + threads_per_block - 1) / threads_per_block;
570
571 if (f->store_jacobians() || !is_analytical<F>()) {
572 compute_hessian_scalar_diagonal_kernel<T, S, Is, num_vertices,
573 F::error_dim, dimension>
574 <<<num_blocks, threads_per_block>>>(
575 diagonal, jacs[Is], f->active_indices.data().get(),
576 f->device_ids.data().get(), hessian_ids[Is],
577 f->vertex_descriptors[Is]->get_active_state(),
578 f->precision_matrices.data().get(),
579 f->chi2_derivative.data().get(), num_threads);
580 } else {
581 if constexpr (is_analytical<F>()) {
582 if (jacobian_scales == nullptr) {
583 compute_hessian_scalar_diagonal_dynamic_kernel<
584 T, S, Is, num_vertices, F::error_dim, dimension,
585 typename F::VertexPointerPointerTuple, F, false>
586 <<<num_blocks, threads_per_block>>>(
587 diagonal, jacs[Is], f->active_indices.data().get(),
588 f->device_ids.data().get(), hessian_ids[Is],
589 f->get_vertices(), f->device_obs.data().get(), nullptr,
590 f->data.data().get(),
591 f->vertex_descriptors[Is]->get_active_state(),
592 f->precision_matrices.data().get(),
593 f->chi2_derivative.data().get(), num_threads);
594 } else {
595 compute_hessian_scalar_diagonal_dynamic_kernel<
596 T, S, Is, num_vertices, F::error_dim, dimension,
597 typename F::VertexPointerPointerTuple, F, true>
598 <<<num_blocks, threads_per_block>>>(
599 diagonal, jacs[Is], f->active_indices.data().get(),
600 f->device_ids.data().get(), hessian_ids[Is],
601 f->get_vertices(), f->device_obs.data().get(),
602 jacobian_scales, f->data.data().get(),
603 f->vertex_descriptors[Is]->get_active_state(),
604 f->precision_matrices.data().get(),
605 f->chi2_derivative.data().get(), num_threads);
606 }
607 }
608 }
609 }()),
610 ...);
611}
612
613template <typename T, typename S, typename F>
614void compute_hessian_scalar_diagonal(F *f, T *diagonal,
615 const T *jacobian_scales) {
616
617 // Then for each vertex, we need to compute the error
618 constexpr auto num_vertices = F::get_num_vertices();
619 constexpr auto vertex_sizes = F::get_vertex_sizes();
620
621 // At this point all necessary data should be on the GPU
622 // std::array<T*, num_vertices> verts;
623 auto verts = f->get_vertices();
624 std::array<S *, num_vertices> jacs;
625 std::array<const size_t *, num_vertices> hessian_ids;
626 for (int i = 0; i < num_vertices; i++) {
627 jacs[i] = f->jacobians[i].data.data().get();
628 hessian_ids[i] = f->vertex_descriptors[i]->get_hessian_ids();
629 }
630
631 const auto num_factors = f->active_count();
632
633 launch_kernel_hessian_scalar_diagonal<T, S>(
634 f, diagonal, hessian_ids, jacs, jacobian_scales, num_factors,
635 std::make_index_sequence<num_vertices>{});
636}
637
638} // namespace ops
639
640} // namespace graphite
The top-level namespace for Graphite.
Definition eigen_solver.cpp:4