335 using Iterator =
typename DoFHandler<dim>::active_cell_iterator;
341 std::array<
NumberType, Components::count_extractors()> __extracted_data{{}};
342 if constexpr (Components::count_extractors() > 0)
343 this->
extract(__extracted_data, solution_global, variables,
true,
false,
true);
344 const auto &extracted_data = __extracted_data;
346 const auto cell_worker = [&](
const Iterator &cell, Scratch &scratch_data, CopyData ©_data) {
347 scratch_data.fe_values.reinit(cell);
348 const auto &fe_v = scratch_data.fe_values;
349 const uint n_dofs = fe_v.get_fe().n_dofs_per_cell();
351 copy_data.reinit(cell, n_dofs);
352 const auto &JxW = fe_v.get_JxW_values();
353 const auto &q_points = fe_v.get_quadrature_points();
354 const auto &q_indices = fe_v.quadrature_point_indices();
356 auto &solution = scratch_data.solution;
357 auto &solution_grad = scratch_data.solution_grad;
358 auto &solution_hess = scratch_data.solution_hess;
359 auto &solution_dot = scratch_data.solution_dot;
360 fe_v.get_function_values(solution_global, solution);
361 fe_v.get_function_gradients(solution_global, solution_grad);
362 fe_v.get_function_hessians(solution_global, solution_hess);
363 fe_v.get_function_values(solution_global_dot, solution_dot);
365 array<Tensor<1, dim, NumberType>, Components::count_fe_functions()> flux{};
366 array<
NumberType, Components::count_fe_functions()> source{};
368 for (
const auto &q_index : q_indices) {
369 const auto &x_q = q_points[q_index];
372 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
375 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
376 model.mass(
mass, x_q, solution[q_index], solution_dot[q_index]);
378 for (
uint i = 0; i < n_dofs; ++i) {
379 const auto component_i =
fe.system_to_component_index(i).first;
380 copy_data.cell_residual(i) += JxW[q_index] * weight *
381 (-scalar_product(fe_v.shape_grad_component(i, q_index, component_i),
383 + fe_v.shape_value_component(i, q_index, component_i) *
384 source[component_i]);
385 copy_data.cell_residual(i) += weight_mass * JxW[q_index] *
386 fe_v.shape_value_component(i, q_index, component_i) *
391 const auto boundary_worker = [&](
const Iterator &cell,
const uint &face_no, Scratch &scratch_data,
392 CopyData ©_data) {
393 scratch_data.fe_interface_values.reinit(cell, face_no);
394 const auto &fe_fv = scratch_data.fe_interface_values.get_fe_face_values(0);
395 const uint n_dofs = fe_fv.get_fe().n_dofs_per_cell();
397 const auto &JxW = fe_fv.get_JxW_values();
398 const auto &q_points = fe_fv.get_quadrature_points();
399 const auto &q_indices = fe_fv.quadrature_point_indices();
400 const std::vector<Tensor<1, dim>> &normals = fe_fv.get_normal_vectors();
402 auto &solution = scratch_data.solution_interface[0];
403 auto &solution_grad = scratch_data.solution_grad_interface[0];
404 auto &solution_hess = scratch_data.solution_hess_interface[0];
405 fe_fv.get_function_values(solution_global, solution);
406 fe_fv.get_function_gradients(solution_global, solution_grad);
407 fe_fv.get_function_hessians(solution_global, solution_hess);
409 array<Tensor<1, dim, NumberType>, Components::count_fe_functions()> numflux{};
410 for (
const auto &q_index : q_indices) {
411 const auto &x_q = q_points[q_index];
412 model.boundary_numflux(
413 numflux, normals[q_index], x_q,
414 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
416 for (
uint i = 0; i < n_dofs; ++i) {
417 const auto component_i =
fe.system_to_component_index(i).first;
418 copy_data.cell_residual(i) +=
419 weight * JxW[q_index] *
420 (fe_fv.shape_value_component(i, q_index, component_i) *
421 scalar_product(numflux[component_i], normals[q_index]));
425 const auto copier = [&](
const CopyData &c) {
426 constraints.distribute_local_to_global(c.cell_residual, c.local_dof_indices,
residual);
431 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells | MeshWorker::assemble_boundary_faces;
499 using Iterator =
typename DoFHandler<dim>::active_cell_iterator;
505 std::array<
NumberType, Components::count_extractors()> extracted_data{{}};
506 if constexpr (Components::count_extractors() > 0) {
507 this->
extract(extracted_data, solution_global, variables,
true,
true,
true);
512 const auto cell_worker = [&](
const Iterator &cell, Scratch &scratch_data, CopyData ©_data) {
513 scratch_data.fe_values.reinit(cell);
514 const auto &fe_v = scratch_data.fe_values;
515 const uint n_dofs = fe_v.get_fe().n_dofs_per_cell();
517 copy_data.reinit(cell, n_dofs, Components::count_extractors());
518 const auto &JxW = fe_v.get_JxW_values();
519 const auto &q_points = fe_v.get_quadrature_points();
520 const auto &q_indices = fe_v.quadrature_point_indices();
522 auto &solution = scratch_data.solution;
523 auto &solution_dot = scratch_data.solution_dot;
524 auto &solution_grad = scratch_data.solution_grad;
525 auto &solution_hess = scratch_data.solution_hess;
527 fe_v.get_function_values(solution_global, solution);
528 fe_v.get_function_gradients(solution_global, solution_grad);
529 fe_v.get_function_hessians(solution_global, solution_hess);
530 fe_v.get_function_values(solution_global_dot, solution_dot);
544 for (
const auto &q_index : q_indices) {
545 const auto &x_q = q_points[q_index];
546 model.template jacobian_flux<0, 0>(
548 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
549 model.template jacobian_source<0, 0>(
551 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
552 model.template jacobian_flux_grad<1>(
554 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
555 model.template jacobian_source_grad<1>(
557 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
558 model.template jacobian_flux_hess<2>(
560 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
561 model.template jacobian_source_hess<2>(
563 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
564 if constexpr (Components::count_extractors() > 0) {
565 model.template jacobian_flux_extr<3>(
567 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
568 model.template jacobian_source_extr<3>(
570 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
576 tbb::blocked_range2d<uint>(0, n_dofs, 0, n_dofs), [&](
const tbb::blocked_range2d<uint> &range) {
577 for (
uint i = range.rows().begin(); i < range.rows().end(); ++i) {
578 const auto component_i = fe_v.get_fe().system_to_component_index(i).first;
579 for (
uint j = range.cols().begin(); j < range.cols().end(); ++j) {
580 const auto component_j =
fe.system_to_component_index(j).first;
582 copy_data.cell_jacobian(i, j) +=
583 weight * JxW[q_index] *
584 fe_v.shape_value_component(j, q_index, component_j) *
585 (-scalar_product(fe_v.shape_grad_component(i, q_index, component_i),
586 j_flux(component_i, component_j))
587 + fe_v.shape_value_component(i, q_index, component_i) *
588 j_source(component_i, component_j));
590 copy_data.cell_jacobian(i, j) +=
591 weight * JxW[q_index] *
592 scalar_product(fe_v.shape_grad_component(j, q_index, component_j),
593 -scalar_product(fe_v.shape_grad_component(i, q_index, component_i),
594 j_grad_flux(component_i, component_j))
595 + fe_v.shape_value_component(i, q_index, component_i) *
596 j_grad_source(component_i, component_j));
598 copy_data.cell_jacobian(i, j) +=
599 weight * JxW[q_index] *
600 scalar_product(fe_v.shape_hessian_component(j, q_index, component_j),
601 -scalar_product(fe_v.shape_grad_component(i, q_index, component_i),
602 j_hess_flux(component_i, component_j)) +
603 fe_v.shape_value_component(i, q_index, component_i) *
604 j_hess_source(component_i, component_j));
606 copy_data.cell_jacobian(i, j) +=
607 JxW[q_index] * fe_v.shape_value_component(j, q_index, component_j) *
608 fe_v.shape_value_component(i, q_index, component_i) *
609 (alpha * j_mass_dot(component_i, component_j) + beta * j_mass(component_i, component_j));
615 if constexpr (Components::count_extractors() > 0) {
616 for (
uint i = 0; i < n_dofs; ++i) {
617 const auto component_i =
fe.system_to_component_index(i).first;
618 for (
uint e = 0; e < Components::count_extractors(); ++e)
619 copy_data.extractor_cell_jacobian(i, e) +=
620 weight * JxW[q_index] *
621 (-scalar_product(fe_v.shape_grad_component(i, q_index, component_i),
622 j_extr_flux(component_i, e))
623 + fe_v.shape_value_component(i, q_index, component_i) *
624 j_extr_source(component_i, e));
630 const auto boundary_worker = [&](
const Iterator &cell,
const uint &face_no, Scratch &scratch_data,
631 CopyData ©_data) {
632 scratch_data.fe_interface_values.reinit(cell, face_no);
633 const auto &fe_fv = scratch_data.fe_interface_values.get_fe_face_values(0);
634 const uint n_dofs = fe_fv.get_fe().n_dofs_per_cell();
636 const auto &JxW = fe_fv.get_JxW_values();
637 const auto &q_points = fe_fv.get_quadrature_points();
638 const auto &q_indices = fe_fv.quadrature_point_indices();
639 const std::vector<Tensor<1, dim>> &normals = fe_fv.get_normal_vectors();
641 auto &solution = scratch_data.solution_interface[0];
642 auto &solution_grad = scratch_data.solution_grad_interface[0];
643 auto &solution_hess = scratch_data.solution_hess;
644 fe_fv.get_function_values(solution_global, solution);
645 fe_fv.get_function_gradients(solution_global, solution_grad);
646 fe_fv.get_function_hessians(solution_global, solution_hess);
650 j_grad_boundary_numflux;
652 j_hess_boundary_numflux;
654 j_extr_boundary_numflux;
655 for (
const auto &q_index : q_indices) {
656 const auto &x_q = q_points[q_index];
657 model.template jacobian_boundary_numflux<0, 0>(
658 j_boundary_numflux, normals[q_index], x_q,
659 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
660 model.template jacobian_boundary_numflux_grad<1>(
661 j_grad_boundary_numflux, normals[q_index], x_q,
662 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
663 model.template jacobian_boundary_numflux_hess<2>(
664 j_hess_boundary_numflux, normals[q_index], x_q,
665 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
666 if constexpr (Components::count_extractors() > 0) {
667 model.template jacobian_boundary_numflux_extr<3>(
668 j_extr_boundary_numflux, normals[q_index], x_q,
669 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
673 tbb::blocked_range2d<uint>(0, n_dofs, 0, n_dofs), [&](
const tbb::blocked_range2d<uint> &range) {
674 for (
uint i = range.rows().begin(); i < range.rows().end(); ++i) {
675 const auto component_i =
fe.system_to_component_index(i).first;
676 for (
uint j = range.cols().begin(); j < range.cols().end(); ++j) {
677 const auto component_j =
fe.system_to_component_index(j).first;
680 copy_data.cell_jacobian(i, j) +=
681 weight * JxW[q_index] *
682 fe_fv.shape_value_component(j, q_index, component_j) *
683 (fe_fv.shape_value_component(i, q_index, component_i) *
684 scalar_product(j_boundary_numflux(component_i, component_j),
687 copy_data.cell_jacobian(i, j) +=
688 weight * JxW[q_index] *
690 fe_fv.shape_grad_component(j, q_index, component_j),
691 fe_fv.shape_value_component(i, q_index, component_i) *
692 scalar_product(j_grad_boundary_numflux(component_i, component_j),
695 copy_data.cell_jacobian(i, j) +=
696 weight * JxW[q_index] *
698 fe_fv.shape_hessian_component(j, q_index, component_j),
699 fe_fv.shape_value_component(i, q_index, component_i) *
700 scalar_product(j_hess_boundary_numflux(component_i, component_j),
707 if constexpr (Components::count_extractors() > 0) {
708 for (
uint i = 0; i < n_dofs; ++i) {
709 const auto component_i =
fe.system_to_component_index(i).first;
710 for (
uint e = 0; e < Components::count_extractors(); ++e)
711 copy_data.extractor_cell_jacobian(i, e) +=
712 weight * JxW[q_index] *
713 (fe_fv.shape_value_component(i, q_index, component_i) *
714 scalar_product(j_extr_boundary_numflux(component_i, e),
721 const auto copier = [&](
const CopyData &c) {
722 constraints.distribute_local_to_global(c.cell_jacobian, c.local_dof_indices,
jacobian);
723 if constexpr (Components::count_extractors() > 0) {
724 FullMatrix<NumberType> extractor_dependence(c.local_dof_indices.size(),
extractor_dof_indices.size());
726 constraints.distribute_local_to_global(extractor_dependence, c.local_dof_indices,
extractor_dof_indices,
733 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells | MeshWorker::assemble_boundary_faces;