233 using Iterator =
typename DoFHandler<dim>::active_cell_iterator;
237 const auto cell_worker = [&](
const Iterator &t_cell, Scratch &scratch_data, CopyData ©_data) {
238 scratch_data.fe_values.reinit(t_cell);
239 const auto &fe_v = scratch_data.fe_values;
240 copy_data.cell_index = t_cell->active_cell_index();
243 const auto &JxW = fe_v.get_JxW_values();
244 const auto &q_points = fe_v.get_quadrature_points();
245 const auto &q_indices = fe_v.quadrature_point_indices();
247 auto &solution = scratch_data.solution;
248 std::vector<std::vector<Tensor<1, dim, NumberType>>> solution_grad;
249 std::vector<std::vector<Tensor<2, dim, NumberType>>> solution_hess;
250 solution_grad.resize(q_points.size(), std::vector<Tensor<1, dim, NumberType>>(
n_components));
251 solution_hess.resize(q_points.size(), std::vector<Tensor<2, dim, NumberType>>(
n_components));
252 fe_v.get_function_values(solution_global, solution);
253 fe_v.get_function_gradients(solution_global, solution_grad);
254 fe_v.get_function_hessians(solution_global, solution_hess);
256 double local_indicator = 0.;
257 for (
const auto &q_index : q_indices) {
258 const auto &x_q = q_points[q_index];
259 const std::vector<double>
nothing = {};
260 model.cell_indicator(local_indicator, x_q,
261 i_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index]));
262 copy_data.value += JxW[q_index] * local_indicator;
265 const auto face_worker = [&](
const Iterator &t_cell,
const uint &f,
const uint &sf,
const Iterator &t_ncell,
266 const uint &nf,
const unsigned int &nsf, Scratch &scratch_data,
267 CopyData ©_data) {
268 scratch_data.fe_interface_values.reinit(t_cell, f, sf, t_ncell, nf, nsf);
269 const auto &fe_iv = scratch_data.fe_interface_values;
270 const auto &fe_iv_s = scratch_data.fe_interface_values.get_fe_face_values(0);
271 const auto &fe_iv_n = scratch_data.fe_interface_values.get_fe_face_values(1);
273 copy_data.face_data.emplace_back();
274 auto ©_data_face = copy_data.face_data.emplace_back();
275 copy_data_face.cell_indices[0] = t_cell->active_cell_index();
276 copy_data_face.cell_indices[1] = t_ncell->active_cell_index();
277 copy_data_face.values[0] = 0;
278 copy_data_face.values[1] = 0;
280 const auto &JxW = fe_iv.get_JxW_values();
281 const auto &q_points = fe_iv.get_quadrature_points();
282 const auto &q_indices = fe_iv.quadrature_point_indices();
283 const std::vector<Tensor<1, dim>> &normals = fe_iv.get_normal_vectors();
285 auto &solution_s = scratch_data.solution_interface[0];
286 auto &solution_n = scratch_data.solution_interface[1];
287 std::vector<std::vector<Tensor<1, dim, NumberType>>> solution_grad_s;
288 std::vector<std::vector<Tensor<2, dim, NumberType>>> solution_hess_s;
289 solution_grad_s.resize(q_points.size(), std::vector<Tensor<1, dim, NumberType>>(
n_components));
290 solution_hess_s.resize(q_points.size(), std::vector<Tensor<2, dim, NumberType>>(
n_components));
291 std::vector<std::vector<Tensor<1, dim, NumberType>>> solution_grad_n;
292 std::vector<std::vector<Tensor<2, dim, NumberType>>> solution_hess_n;
293 solution_grad_n.resize(q_points.size(), std::vector<Tensor<1, dim, NumberType>>(
n_components));
294 solution_hess_n.resize(q_points.size(), std::vector<Tensor<2, dim, NumberType>>(
n_components));
295 fe_iv_s.get_function_values(solution_global, solution_s);
296 fe_iv_n.get_function_values(solution_global, solution_n);
297 fe_iv_s.get_function_gradients(solution_global, solution_grad_s);
298 fe_iv_n.get_function_gradients(solution_global, solution_grad_n);
299 fe_iv_s.get_function_hessians(solution_global, solution_hess_s);
300 fe_iv_n.get_function_hessians(solution_global, solution_hess_n);
302 array<double, 2> local_indicator{};
303 for (
const auto &q_index : q_indices) {
304 const auto &x_q = q_points[q_index];
305 model.face_indicator(local_indicator, normals[q_index], x_q,
306 i_tie(solution_s[q_index], solution_grad_s[q_index], solution_hess_s[q_index]),
307 i_tie(solution_n[q_index], solution_grad_n[q_index], solution_hess_n[q_index]));
309 copy_data_face.values[0] += JxW[q_index] * local_indicator[0] * (1. + t_cell->at_boundary());
310 copy_data_face.values[1] += JxW[q_index] * local_indicator[1] * (1. + t_ncell->at_boundary());
313 const auto copier = [&](
const CopyData &c) {
314 for (
auto &cdf : c.face_data)
315 for (
uint j = 0; j < 2; ++j)
316 indicator[cdf.cell_indices[j]] += cdf.values[j];
317 indicator[c.cell_index] += c.value;
319 const UpdateFlags update_flags =
320 update_values | update_gradients | update_quadrature_points | update_JxW_values | update_hessians;
321 const UpdateFlags interface_update_flags = update_values | update_gradients | update_quadrature_points |
322 update_JxW_values | update_normal_vectors | update_hessians;
326 MeshWorker::AssembleFlags assemble_flags = MeshWorker::assemble_own_cells;
385 using Iterator =
typename DoFHandler<dim>::active_cell_iterator;
391 std::array<
NumberType, Components::count_extractors()> __extracted_data{{}};
392 if constexpr (Components::count_extractors() > 0)
393 this->
extract(__extracted_data, solution_global, variables,
true,
false,
true);
394 const auto &extracted_data = __extracted_data;
396 const auto cell_worker = [&](
const Iterator &cell, Scratch &scratch_data, CopyData ©_data) {
397 scratch_data.fe_values.reinit(cell);
398 const auto &fe_v = scratch_data.fe_values;
399 const uint n_dofs = fe_v.get_fe().n_dofs_per_cell();
401 copy_data.reinit(cell, n_dofs);
402 const auto &JxW = fe_v.get_JxW_values();
403 const auto &q_points = fe_v.get_quadrature_points();
404 const auto &q_indices = fe_v.quadrature_point_indices();
406 auto &solution = scratch_data.solution;
407 auto &solution_dot = scratch_data.solution_dot;
408 fe_v.get_function_values(solution_global, solution);
409 fe_v.get_function_values(solution_global_dot, solution_dot);
411 array<NumberType, n_components>
mass{};
413 array<NumberType, n_components> source{};
414 for (
const auto &q_index : q_indices) {
415 const auto &x_q = q_points[q_index];
416 model.mass(
mass, x_q, solution[q_index], solution_dot[q_index]);
417 model.flux(flux, x_q,
fe_tie(solution[q_index], extracted_data, variables));
418 model.source(source, x_q,
fe_tie(solution[q_index], extracted_data, variables));
420 for (
uint i = 0; i < n_dofs; ++i) {
421 const auto component_i =
fe.system_to_component_index(i).first;
422 copy_data.cell_mass(i) += weight_mass * JxW[q_index] *
423 fe_v.shape_value_component(i, q_index, component_i) *
425 copy_data.cell_residual(i) += JxW[q_index] * weight *
426 (-scalar_product(fe_v.shape_grad_component(i, q_index, component_i),
428 + fe_v.shape_value_component(i, q_index, component_i) *
429 source[component_i]);
433 const auto boundary_worker = [&](
const Iterator &cell,
const uint &face_no, Scratch &scratch_data,
434 CopyData ©_data) {
435 scratch_data.fe_interface_values.reinit(cell, face_no);
436 const auto &fe_fv = scratch_data.fe_interface_values.get_fe_face_values(0);
437 const uint n_dofs = fe_fv.get_fe().n_dofs_per_cell();
439 const auto &JxW = fe_fv.get_JxW_values();
440 const auto &q_points = fe_fv.get_quadrature_points();
441 const auto &q_indices = fe_fv.quadrature_point_indices();
442 auto &solution = scratch_data.solution_interface[0];
444 array<Tensor<1, dim, NumberType>,
n_components> numflux{};
446 fe_fv.get_function_values(solution_global, solution);
447 const std::vector<Tensor<1, dim>> &normals = fe_fv.get_normal_vectors();
449 for (
const auto &q_index : q_indices) {
450 const auto &x_q = q_points[q_index];
451 model.boundary_numflux(numflux, normals[q_index], x_q,
452 fe_tie(solution[q_index], extracted_data, variables));
454 for (
uint i = 0; i < n_dofs; ++i) {
455 const auto component_i =
fe.system_to_component_index(i).first;
456 copy_data.cell_residual(i) +=
457 weight * JxW[q_index] *
458 (fe_fv.shape_value_component(i, q_index, component_i) *
459 scalar_product(numflux[component_i], normals[q_index]));
463 const auto face_worker = [&](
const Iterator &cell,
const uint &f,
const uint &sf,
const Iterator &ncell,
464 const uint &nf,
const uint &nsf, Scratch &scratch_data, CopyData ©_data) {
465 scratch_data.fe_interface_values.reinit(cell, f, sf, ncell, nf, nsf);
466 const auto &fe_iv = scratch_data.fe_interface_values;
467 const auto &fe_iv_s = scratch_data.fe_interface_values.get_fe_face_values(0);
468 const auto &fe_iv_n = scratch_data.fe_interface_values.get_fe_face_values(1);
469 const uint n_dofs = fe_iv.n_current_interface_dofs();
471 copy_data.face_data.emplace_back();
472 auto ©_data_face = copy_data.face_data.back();
473 copy_data_face.reinit(fe_iv);
475 const auto &JxW = fe_iv.get_JxW_values();
476 const auto &q_points = fe_iv.get_quadrature_points();
477 const auto &q_indices = fe_iv.quadrature_point_indices();
478 auto &solution_s = scratch_data.solution_interface[0];
479 auto &solution_n = scratch_data.solution_interface[1];
481 array<Tensor<1, dim, NumberType>,
n_components> numflux{};
483 fe_iv_s.get_function_values(solution_global, solution_s);
484 fe_iv_n.get_function_values(solution_global, solution_n);
485 const std::vector<Tensor<1, dim>> &normals = fe_iv.get_normal_vectors();
487 for (
const auto &q_index : q_indices) {
488 const auto &x_q = q_points[q_index];
489 model.numflux(numflux, normals[q_index], x_q,
fe_tie(solution_s[q_index], extracted_data, variables),
490 fe_tie(solution_n[q_index], extracted_data, variables));
492 for (
uint i = 0; i < n_dofs; ++i) {
493 const auto &cd_i = fe_iv.interface_dof_to_dof_indices(i);
494 const auto component_i = cd_i[0] == numbers::invalid_unsigned_int
495 ?
fe.system_to_component_index(cd_i[1]).first
496 :
fe.system_to_component_index(cd_i[0]).first;
497 copy_data_face.cell_residual(i) +=
498 weight * JxW[q_index] *
499 (fe_iv.jump_in_shape_values(i, q_index, component_i) *
500 scalar_product(numflux[component_i],
505 const auto copier = [&](
const CopyData &c) {
506 constraints.distribute_local_to_global(c.cell_residual, c.local_dof_indices,
residual);
507 constraints.distribute_local_to_global(c.cell_mass, c.local_dof_indices,
residual);
508 for (
auto &cdf : c.face_data)
509 constraints.distribute_local_to_global(cdf.cell_residual, cdf.joint_dof_indices,
residual);
514 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells | MeshWorker::assemble_boundary_faces |
515 MeshWorker::assemble_own_interior_faces_once;
584 using Iterator =
typename DoFHandler<dim>::active_cell_iterator;
590 std::array<
NumberType, Components::count_extractors()> extracted_data{{}};
591 if constexpr (Components::count_extractors() > 0) {
592 this->
extract(extracted_data, solution_global, variables,
true,
true,
true);
597 const auto cell_worker = [&](
const Iterator &cell, Scratch &scratch_data, CopyData ©_data) {
598 scratch_data.fe_values.reinit(cell);
599 const auto &fe_v = scratch_data.fe_values;
600 const uint n_dofs = fe_v.get_fe().n_dofs_per_cell();
602 copy_data.reinit(cell, n_dofs, Components::count_extractors());
603 const auto &JxW = fe_v.get_JxW_values();
604 const auto &q_points = fe_v.get_quadrature_points();
605 const auto &q_indices = fe_v.quadrature_point_indices();
607 auto &solution = scratch_data.solution;
608 auto &solution_dot = scratch_data.solution_dot;
609 fe_v.get_function_values(solution_global, solution);
610 fe_v.get_function_values(solution_global_dot, solution_dot);
618 for (
const auto &q_index : q_indices) {
619 const auto &x_q = q_points[q_index];
622 model.template jacobian_flux<0, 0>(j_flux, x_q,
fe_tie(solution[q_index], extracted_data, variables));
623 model.template jacobian_source<0, 0>(j_source, x_q,
fe_tie(solution[q_index], extracted_data, variables));
624 if constexpr (Components::count_extractors() > 0) {
625 model.template jacobian_flux_extr<1>(j_extr_flux, x_q,
626 fe_tie(solution[q_index], extracted_data, variables));
627 model.template jacobian_source_extr<1>(j_extr_source, x_q,
628 fe_tie(solution[q_index], extracted_data, variables));
631 for (
uint i = 0; i < n_dofs; ++i) {
632 const auto component_i =
fe.system_to_component_index(i).first;
633 for (
uint j = 0; j < n_dofs; ++j) {
634 const auto component_j =
fe.system_to_component_index(j).first;
635 copy_data.cell_jacobian(i, j) += weight * JxW[q_index] *
636 fe_v.shape_value_component(j, q_index, component_j) *
637 (-scalar_product(fe_v.shape_grad_component(i, q_index, component_i),
638 j_flux(component_i, component_j))
639 + fe_v.shape_value_component(i, q_index, component_i) *
640 j_source(component_i, component_j));
641 copy_data.cell_mass_jacobian(i, j) +=
642 JxW[q_index] * fe_v.shape_value_component(j, q_index, component_j) *
643 fe_v.shape_value_component(i, q_index, component_i) *
644 (alpha * j_mass_dot(component_i, component_j) + beta * j_mass(component_i, component_j));
647 if constexpr (Components::count_extractors() > 0)
648 for (
uint e = 0; e < Components::count_extractors(); ++e)
649 copy_data.extractor_cell_jacobian(i, e) +=
650 weight * JxW[q_index] *
651 (-scalar_product(fe_v.shape_grad_component(i, q_index, component_i),
652 j_extr_flux(component_i, e))
653 + fe_v.shape_value_component(i, q_index, component_i) *
654 j_extr_source(component_i, e));
658 const auto boundary_worker = [&](
const Iterator &cell,
const uint &face_no, Scratch &scratch_data,
659 CopyData ©_data) {
660 scratch_data.fe_interface_values.reinit(cell, face_no);
661 const auto &fe_fv = scratch_data.fe_interface_values.get_fe_face_values(0);
662 const uint n_dofs = fe_fv.get_fe().n_dofs_per_cell();
664 const auto &JxW = fe_fv.get_JxW_values();
665 const auto &q_points = fe_fv.get_quadrature_points();
666 const auto &q_indices = fe_fv.quadrature_point_indices();
667 auto &solution = scratch_data.solution_interface[0];
672 fe_fv.get_function_values(solution_global, solution);
673 const std::vector<Tensor<1, dim>> &normals = fe_fv.get_normal_vectors();
675 for (
const auto &q_index : q_indices) {
676 const auto &x_q = q_points[q_index];
677 model.template jacobian_boundary_numflux<0, 0>(j_boundary_numflux, normals[q_index], x_q,
678 fe_tie(solution[q_index], extracted_data, variables));
679 if constexpr (Components::count_extractors() > 0)
680 model.template jacobian_boundary_numflux_extr<1>(j_extr_boundary_numflux, normals[q_index], x_q,
681 fe_tie(solution[q_index], extracted_data, variables));
683 for (
uint i = 0; i < n_dofs; ++i) {
684 const auto component_i =
fe.system_to_component_index(i).first;
685 for (
uint j = 0; j < n_dofs; ++j) {
686 const auto component_j =
fe.system_to_component_index(j).first;
687 copy_data.cell_jacobian(i, j) +=
688 weight * JxW[q_index] * fe_fv.shape_value_component(j, q_index, component_j) *
689 (fe_fv.shape_value_component(i, q_index, component_i) *
690 scalar_product(j_boundary_numflux(component_i, component_j),
694 if constexpr (Components::count_extractors() > 0)
695 for (
uint e = 0; e < Components::count_extractors(); ++e)
696 copy_data.extractor_cell_jacobian(i, e) +=
697 weight * JxW[q_index] *
698 (fe_fv.shape_value_component(i, q_index, component_i) *
699 scalar_product(j_extr_boundary_numflux(component_i, e),
704 const auto face_worker = [&](
const Iterator &cell,
const uint &f,
const uint &sf,
const Iterator &ncell,
705 const uint &nf,
const uint &nsf, Scratch &scratch_data, CopyData ©_data) {
706 scratch_data.fe_interface_values.reinit(cell, f, sf, ncell, nf, nsf);
707 const auto &fe_iv = scratch_data.fe_interface_values;
708 const auto &fe_iv_s = scratch_data.fe_interface_values.get_fe_face_values(0);
709 const auto &fe_iv_n = scratch_data.fe_interface_values.get_fe_face_values(1);
710 const uint n_dofs = fe_iv.n_current_interface_dofs();
712 copy_data.face_data.emplace_back();
713 auto ©_data_face = copy_data.face_data.back();
714 copy_data_face.reinit(fe_iv, Components::count_extractors());
716 const auto &JxW = fe_iv.get_JxW_values();
717 const auto &q_points = fe_iv.get_quadrature_points();
718 const auto &q_indices = fe_iv.quadrature_point_indices();
719 auto &solution_s = scratch_data.solution_interface[0];
720 auto &solution_n = scratch_data.solution_interface[1];
722 fe_iv_s.get_function_values(solution_global, solution_s);
723 fe_iv_n.get_function_values(solution_global, solution_n);
724 const std::vector<Tensor<1, dim>> &normals = fe_iv.get_normal_vectors();
726 array<SimpleMatrix<Tensor<1, dim>,
n_components>, 2> j_numflux;
727 array<SimpleMatrix<Tensor<1, dim>,
n_components, Components::count_extractors()>, 2> j_extr_numflux;
728 for (
const auto &q_index : q_indices) {
729 const auto &x_q = q_points[q_index];
730 model.template jacobian_numflux<0, 0>(j_numflux, normals[q_index], x_q,
731 fe_tie(solution_s[q_index], extracted_data, variables),
732 fe_tie(solution_n[q_index], extracted_data, variables));
733 if constexpr (Components::count_extractors() > 0)
734 model.template jacobian_numflux_extr<1>(j_extr_numflux, normals[q_index], x_q,
735 fe_tie(solution_s[q_index], extracted_data, variables),
736 fe_tie(solution_n[q_index], extracted_data, variables));
738 for (
uint i = 0; i < n_dofs; ++i) {
739 const auto &cd_i = fe_iv.interface_dof_to_dof_indices(i);
740 const uint face_no_i = cd_i[0] == numbers::invalid_unsigned_int ? 1 : 0;
741 const auto &component_i =
fe.system_to_component_index(cd_i[face_no_i]).first;
742 for (
uint j = 0; j < n_dofs; ++j) {
743 const auto &cd_j = fe_iv.interface_dof_to_dof_indices(j);
744 const uint face_no_j = cd_j[0] == numbers::invalid_unsigned_int ? 1 : 0;
745 const auto &component_j =
fe.system_to_component_index(cd_j[face_no_j]).first;
747 copy_data_face.cell_jacobian(i, j) +=
748 weight * JxW[q_index] *
749 fe_iv.get_fe_face_values(face_no_j).shape_value_component(cd_j[face_no_j], q_index,
751 (fe_iv.jump_in_shape_values(i, q_index, component_i) *
752 scalar_product(j_numflux[face_no_j](component_i, component_j),
756 if constexpr (Components::count_extractors() > 0)
757 for (
uint e = 0; e < Components::count_extractors(); ++e)
758 copy_data_face.extractor_cell_jacobian(i, e) +=
759 weight * JxW[q_index] *
760 (fe_iv.jump_in_shape_values(i, q_index, component_i) *
761 scalar_product(j_extr_numflux[face_no_i](component_i, e),
766 const auto copier = [&](
const CopyData &c) {
767 constraints.distribute_local_to_global(c.cell_jacobian, c.local_dof_indices,
jacobian);
768 constraints.distribute_local_to_global(c.cell_mass_jacobian, c.local_dof_indices,
jacobian);
769 for (
auto &cdf : c.face_data) {
770 constraints.distribute_local_to_global(cdf.cell_jacobian, cdf.joint_dof_indices,
jacobian);
771 if constexpr (Components::count_extractors() > 0) {
772 FullMatrix<NumberType> extractor_dependence(cdf.joint_dof_indices.size(),
extractor_dof_indices.size());
774 constraints.distribute_local_to_global(extractor_dependence, cdf.joint_dof_indices,
extractor_dof_indices,
778 if constexpr (Components::count_extractors() > 0) {
779 FullMatrix<NumberType> extractor_dependence(c.local_dof_indices.size(),
extractor_dof_indices.size());
781 constraints.distribute_local_to_global(extractor_dependence, c.local_dof_indices,
extractor_dof_indices,
788 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells | MeshWorker::assemble_boundary_faces |
789 MeshWorker::assemble_own_interior_faces_once;