libMesh
Classes | Functions
libMesh::RBDataSerialization Namespace Reference

Classes

class  RBEIMEvaluationSerialization
 This class serializes an RBEIMEvaluation object using the Cap'n Proto library. More...
 
class  RBEvaluationSerialization
 This class serializes an RBEvaluation object using the Cap'n Proto library. More...
 
class  RBSCMEvaluationSerialization
 This class serializes an RBSCMEvaluation object using the Cap'n Proto library. More...
 
class  TransientRBEvaluationSerialization
 This class serializes a TransientRBEvaluation object using the Cap'n Proto library. More...
 

Functions

void add_parameter_ranges_to_builder (const RBParametrized &rb_evaluation, RBData::ParameterRanges::Builder &parameter_ranges, RBData::DiscreteParameterList::Builder &discrete_parameters_list)
 Add parameter ranges for continuous and discrete parameters. More...
 
template<typename RBEvaluationBuilderNumber >
void add_rb_evaluation_data_to_builder (RBEvaluation &rb_eval, RBEvaluationBuilderNumber &rb_eval_builder)
 Add data for an RBEvaluation to the builder. More...
 
template<typename RBEvaluationBuilderNumber , typename TransRBEvaluationBuilderNumber >
void add_transient_rb_evaluation_data_to_builder (TransientRBEvaluation &trans_rb_eval, RBEvaluationBuilderNumber &rb_eval_builder, TransRBEvaluationBuilderNumber &trans_rb_eval_builder)
 Add data for a TransientRBEvaluation to the builder. More...
 
template<typename RBEvaluationBuilderNumber , typename RBEIMEvaluationBuilderNumber >
void add_rb_eim_evaluation_data_to_builder (RBEIMEvaluation &rb_eim_eval, RBEvaluationBuilderNumber &rb_eval_builder, RBEIMEvaluationBuilderNumber &rb_eim_eval_builder)
 Add data for an RBEIMEvaluation to the builder. More...
 
void add_rb_scm_evaluation_data_to_builder (RBSCMEvaluation &rb_scm_eval, RBData::RBSCMEvaluation::Builder &rb_scm_eval_builder)
 Add data for an RBSCMEvaluation to the builder. More...
 
void add_point_to_builder (const Point &point, RBData::Point3D::Builder point_builder)
 Helper function that adds point data. More...
 
void add_elem_to_builder (const Elem &elem, RBData::MeshElem::Builder mesh_elem_builder)
 Helper function that adds element data. More...
 

Function Documentation

◆ add_elem_to_builder()

void libMesh::RBDataSerialization::add_elem_to_builder ( const Elem elem,
RBData::MeshElem::Builder  mesh_elem_builder 
)

Helper function that adds element data.

Definition at line 749 of file rb_data_serialization.C.

750 {
751  std::string elem_type_string = libMesh::Utility::enum_to_string(elem.type());
752 
753  mesh_elem_builder.setType(elem_type_string.c_str());
754  mesh_elem_builder.setSubdomainId(elem.subdomain_id());
755 
756  unsigned int n_points = elem.n_nodes();
757  auto mesh_elem_point_list = mesh_elem_builder.initPoints(n_points);
758 
759  for (unsigned int j=0; j < n_points; ++j)
760  {
761  add_point_to_builder(elem.node_ref(j), mesh_elem_point_list[j]);
762  }
763 }

References add_point_to_builder(), libMesh::Utility::enum_to_string(), libMesh::Elem::n_nodes(), libMesh::Elem::node_ref(), libMesh::Elem::subdomain_id(), and libMesh::Elem::type().

Referenced by add_rb_eim_evaluation_data_to_builder().

◆ add_parameter_ranges_to_builder()

void libMesh::RBDataSerialization::add_parameter_ranges_to_builder ( const RBParametrized rb_evaluation,
RBData::ParameterRanges::Builder &  parameter_ranges,
RBData::DiscreteParameterList::Builder &  discrete_parameters_list 
)

Add parameter ranges for continuous and discrete parameters.

Definition at line 263 of file rb_data_serialization.C.

266 {
267  // Continuous parameters
268  {
269  unsigned int n_continuous_parameters = rb_evaluation.get_n_continuous_params();
270  auto names = parameter_ranges_list.initNames(n_continuous_parameters);
271  auto mins = parameter_ranges_list.initMinValues(n_continuous_parameters);
272  auto maxs = parameter_ranges_list.initMaxValues(n_continuous_parameters);
273 
274  std::set<std::string> parameter_names = rb_evaluation.get_parameter_names();
275  const RBParameters & parameters_min = rb_evaluation.get_parameters_min();
276  const RBParameters & parameters_max = rb_evaluation.get_parameters_max();
277 
278  unsigned int count = 0;
279  for (const auto & parameter_name : parameter_names)
280  {
281  if (!rb_evaluation.is_discrete_parameter(parameter_name))
282  {
283  names.set(count, parameter_name);
284  mins.set(count, parameters_min.get_value(parameter_name));
285  maxs.set(count, parameters_max.get_value(parameter_name));
286 
287  ++count;
288  }
289  }
290 
291  if (count != n_continuous_parameters)
292  libmesh_error_msg("Mismatch in number of continuous parameters");
293  }
294 
295  // Discrete parameters
296  {
297  unsigned int n_discrete_parameters = rb_evaluation.get_n_discrete_params();
298  auto names = discrete_parameters_list.initNames(n_discrete_parameters);
299  auto values_outer = discrete_parameters_list.initValues(n_discrete_parameters);
300 
301  const std::map<std::string, std::vector<Real>> & discrete_parameters =
302  rb_evaluation.get_discrete_parameter_values();
303 
304  unsigned int count = 0;
305  for (const auto & discrete_parameter : discrete_parameters)
306  {
307  names.set(count, discrete_parameter.first);
308 
309  const std::vector<Real> & values = discrete_parameter.second;
310  unsigned int n_values = values.size();
311 
312  values_outer.init(count, n_values);
313  auto values_inner = values_outer[count];
314  for (unsigned int i=0; i<n_values; ++i)
315  {
316  values_inner.set(i, values[i]);
317  }
318 
319  ++count;
320  }
321 
322  if (count != n_discrete_parameters)
323  libmesh_error_msg("Mismatch in number of discrete parameters");
324  }
325 }

References libMesh::RBParametrized::get_discrete_parameter_values(), libMesh::RBParametrized::get_n_continuous_params(), libMesh::RBParametrized::get_n_discrete_params(), libMesh::RBParametrized::get_parameter_names(), libMesh::RBParametrized::get_parameters_max(), libMesh::RBParametrized::get_parameters_min(), libMesh::RBParameters::get_value(), and libMesh::RBParametrized::is_discrete_parameter().

Referenced by add_rb_evaluation_data_to_builder(), and add_rb_scm_evaluation_data_to_builder().

◆ add_point_to_builder()

void libMesh::RBDataSerialization::add_point_to_builder ( const Point point,
RBData::Point3D::Builder  point_builder 
)

Helper function that adds point data.

Definition at line 738 of file rb_data_serialization.C.

739 {
740  point_builder.setX(point(0));
741 
742  if (LIBMESH_DIM >= 2)
743  point_builder.setY(point(1));
744 
745  if (LIBMESH_DIM >= 3)
746  point_builder.setZ(point(2));
747 }

Referenced by add_elem_to_builder(), and add_rb_eim_evaluation_data_to_builder().

◆ add_rb_eim_evaluation_data_to_builder()

template<typename RBEvaluationBuilderNumber , typename RBEIMEvaluationBuilderNumber >
void libMesh::RBDataSerialization::add_rb_eim_evaluation_data_to_builder ( RBEIMEvaluation rb_eim_eval,
RBEvaluationBuilderNumber &  rb_eval_builder,
RBEIMEvaluationBuilderNumber &  rb_eim_eval_builder 
)

Add data for an RBEIMEvaluation to the builder.

Templated to deal with both Real and Complex numbers.

Definition at line 602 of file rb_data_serialization.C.

605 {
606  add_rb_evaluation_data_to_builder(rb_eim_evaluation, rb_evaluation_builder);
607 
608  unsigned int n_bfs = rb_eim_evaluation.get_n_basis_functions();
609 
610  // EIM interpolation matrix
611  {
612  // We store the lower triangular part of an NxN matrix, the size of which is given by
613  // (N(N + 1))/2
614  unsigned int half_matrix_size = n_bfs*(n_bfs+1)/2;
615 
616  auto interpolation_matrix_list =
617  rb_eim_evaluation_builder.initInterpolationMatrix(half_matrix_size);
618  for (unsigned int i=0; i < n_bfs; ++i)
619  for (unsigned int j=0; j <= i; ++j)
620  {
621  unsigned int offset = i*(i+1)/2 + j;
622  set_scalar_in_list(interpolation_matrix_list,
623  offset,
624  rb_eim_evaluation.interpolation_matrix(i,j));
625  }
626  }
627 
628  // Interpolation points
629  {
630  auto interpolation_points_list =
631  rb_eim_evaluation_builder.initInterpolationPoints(n_bfs);
632  for (unsigned int i=0; i < n_bfs; ++i)
633  add_point_to_builder(rb_eim_evaluation.interpolation_points[i],
634  interpolation_points_list[i]);
635  }
636 
637  // Interpolation points variables
638  {
639  auto interpolation_points_var_list =
640  rb_eim_evaluation_builder.initInterpolationPointsVar(n_bfs);
641  for (unsigned int i=0; i<n_bfs; ++i)
642  interpolation_points_var_list.set(i,
643  rb_eim_evaluation.interpolation_points_var[i]);
644  }
645 
646  // Interpolation elements
647  {
648  unsigned int n_interpolation_elems =
649  rb_eim_evaluation.interpolation_points_elem.size();
650  auto interpolation_points_elem_list =
651  rb_eim_evaluation_builder.initInterpolationPointsElems(n_interpolation_elems);
652 
653  if (n_interpolation_elems != n_bfs)
654  libmesh_error_msg("The number of elements should match the number of basis functions");
655 
656  for (unsigned int i=0; i<n_interpolation_elems; ++i)
657  {
658  const libMesh::Elem & elem = *rb_eim_evaluation.interpolation_points_elem[i];
659  auto mesh_elem_builder = interpolation_points_elem_list[i];
660  add_elem_to_builder(elem, mesh_elem_builder);
661  }
662  }
663 }

References add_elem_to_builder(), add_point_to_builder(), add_rb_evaluation_data_to_builder(), libMesh::RBEvaluation::get_n_basis_functions(), libMesh::RBEIMEvaluation::interpolation_matrix, libMesh::RBEIMEvaluation::interpolation_points, libMesh::RBEIMEvaluation::interpolation_points_elem, and libMesh::RBEIMEvaluation::interpolation_points_var.

Referenced by libMesh::RBDataSerialization::RBEIMEvaluationSerialization::write_to_file().

◆ add_rb_evaluation_data_to_builder()

template<typename RBEvaluationBuilderNumber >
void libMesh::RBDataSerialization::add_rb_evaluation_data_to_builder ( RBEvaluation rb_eval,
RBEvaluationBuilderNumber &  rb_eval_builder 
)

Add data for an RBEvaluation to the builder.

Definition at line 328 of file rb_data_serialization.C.

330 {
331  const RBThetaExpansion & rb_theta_expansion = rb_eval.get_rb_theta_expansion();
332 
333  unsigned int n_F_terms = rb_theta_expansion.get_n_F_terms();
334  unsigned int n_A_terms = rb_theta_expansion.get_n_A_terms();
335 
336  // Number of basis functions
337  unsigned int n_bfs = rb_eval.get_n_basis_functions();
338  rb_evaluation_builder.setNBfs(n_bfs);
339 
340  // Fq representor inner-product data
341  {
342  unsigned int Q_f_hat = n_F_terms*(n_F_terms+1)/2;
343 
344  auto fq_innerprods_list = rb_evaluation_builder.initFqInnerprods(Q_f_hat);
345 
346  for (unsigned int i=0; i<Q_f_hat; i++)
347  set_scalar_in_list(fq_innerprods_list,
348  i,
349  rb_eval.Fq_representor_innerprods[i]);
350  }
351 
352  // FqAq representor inner-product data
353  {
354  auto fq_aq_innerprods_list =
355  rb_evaluation_builder.initFqAqInnerprods(n_F_terms*n_A_terms*n_bfs);
356 
357  for (unsigned int q_f=0; q_f < n_F_terms; ++q_f)
358  for (unsigned int q_a=0; q_a < n_A_terms; ++q_a)
359  for (unsigned int i=0; i < n_bfs; ++i)
360  {
361  unsigned int offset = q_f*n_A_terms*n_bfs + q_a*n_bfs + i;
362  set_scalar_in_list(
363  fq_aq_innerprods_list, offset,
364  rb_eval.Fq_Aq_representor_innerprods[q_f][q_a][i]);
365  }
366  }
367 
368  // AqAq representor inner-product data
369  {
370  unsigned int Q_a_hat = n_A_terms*(n_A_terms+1)/2;
371  auto aq_aq_innerprods_list =
372  rb_evaluation_builder.initAqAqInnerprods(Q_a_hat*n_bfs*n_bfs);
373 
374  for (unsigned int i=0; i < Q_a_hat; ++i)
375  for (unsigned int j=0; j < n_bfs; ++j)
376  for (unsigned int l=0; l < n_bfs; ++l)
377  {
378  unsigned int offset = i*n_bfs*n_bfs + j*n_bfs + l;
379  set_scalar_in_list(
380  aq_aq_innerprods_list,
381  offset,
382  rb_eval.Aq_Aq_representor_innerprods[i][j][l]);
383  }
384  }
385 
386  // Output dual inner-product data, and output vectors
387  {
388  unsigned int n_outputs = rb_theta_expansion.get_n_outputs();
389  auto output_innerprod_outer = rb_evaluation_builder.initOutputDualInnerprods(n_outputs);
390  auto output_vector_outer = rb_evaluation_builder.initOutputVectors(n_outputs);
391 
392  for (unsigned int output_id=0; output_id < n_outputs; ++output_id)
393  {
394  unsigned int n_output_terms = rb_theta_expansion.get_n_output_terms(output_id);
395 
396  {
397  unsigned int Q_l_hat = n_output_terms*(n_output_terms+1)/2;
398  auto output_innerprod_inner = output_innerprod_outer.init(output_id, Q_l_hat);
399  for (unsigned int q=0; q < Q_l_hat; ++q)
400  {
401  set_scalar_in_list(
402  output_innerprod_inner, q, rb_eval.output_dual_innerprods[output_id][q]);
403  }
404  }
405 
406  {
407  auto output_vector_middle = output_vector_outer.init(output_id, n_output_terms);
408  for (unsigned int q_l=0; q_l<n_output_terms; ++q_l)
409  {
410  auto output_vector_inner = output_vector_middle.init(q_l, n_bfs);
411  for (unsigned int j=0; j<n_bfs; ++j)
412  {
413  set_scalar_in_list(
414  output_vector_inner, j, rb_eval.RB_output_vectors[output_id][q_l](j));
415  }
416  }
417  }
418  }
419  }
420 
421  // Fq vectors and Aq matrices
422  {
423  unsigned int n_F_terms = rb_theta_expansion.get_n_F_terms();
424  unsigned int n_A_terms = rb_theta_expansion.get_n_A_terms();
425 
426  auto rb_fq_vectors_outer_list = rb_evaluation_builder.initRbFqVectors(n_F_terms);
427  for (unsigned int q_f=0; q_f < n_F_terms; ++q_f)
428  {
429  auto rb_fq_vectors_inner_list = rb_fq_vectors_outer_list.init(q_f, n_bfs);
430  for (unsigned int i=0; i<n_bfs; i++)
431  set_scalar_in_list(rb_fq_vectors_inner_list, i, rb_eval.RB_Fq_vector[q_f](i));
432  }
433 
434  auto rb_Aq_matrices_outer_list = rb_evaluation_builder.initRbAqMatrices(n_A_terms);
435  for (unsigned int q_a=0; q_a < n_A_terms; ++q_a)
436  {
437  auto rb_Aq_matrices_inner_list = rb_Aq_matrices_outer_list.init(q_a, n_bfs*n_bfs);
438  for (unsigned int i=0; i < n_bfs; ++i)
439  for (unsigned int j=0; j < n_bfs; ++j)
440  {
441  unsigned int offset = i*n_bfs+j;
442  set_scalar_in_list(rb_Aq_matrices_inner_list, offset, rb_eval.RB_Aq_vector[q_a](i,j));
443  }
444  }
445  }
446 
447  // Inner-product matrix
448  if (rb_eval.compute_RB_inner_product)
449  {
450  auto rb_inner_product_matrix_list =
451  rb_evaluation_builder.initRbInnerProductMatrix(n_bfs*n_bfs);
452 
453  for (unsigned int i=0; i < n_bfs; ++i)
454  for (unsigned int j=0; j < n_bfs; ++j)
455  {
456  unsigned int offset = i*n_bfs + j;
457  set_scalar_in_list(
458  rb_inner_product_matrix_list,
459  offset,
460  rb_eval.RB_inner_product_matrix(i,j) );
461  }
462  }
463 
464  auto parameter_ranges_list =
465  rb_evaluation_builder.initParameterRanges();
466  auto discrete_parameters_list =
467  rb_evaluation_builder.initDiscreteParameters();
469  parameter_ranges_list,
470  discrete_parameters_list);
471 }

References add_parameter_ranges_to_builder(), libMesh::RBEvaluation::Aq_Aq_representor_innerprods, libMesh::RBEvaluation::compute_RB_inner_product, libMesh::RBEvaluation::Fq_Aq_representor_innerprods, libMesh::RBEvaluation::Fq_representor_innerprods, libMesh::RBThetaExpansion::get_n_A_terms(), libMesh::RBEvaluation::get_n_basis_functions(), libMesh::RBThetaExpansion::get_n_F_terms(), libMesh::RBThetaExpansion::get_n_output_terms(), libMesh::RBThetaExpansion::get_n_outputs(), libMesh::RBEvaluation::get_rb_theta_expansion(), libMesh::RBEvaluation::output_dual_innerprods, libMesh::RBEvaluation::RB_Aq_vector, libMesh::RBEvaluation::RB_Fq_vector, libMesh::RBEvaluation::RB_inner_product_matrix, and libMesh::RBEvaluation::RB_output_vectors.

Referenced by add_rb_eim_evaluation_data_to_builder(), add_transient_rb_evaluation_data_to_builder(), and libMesh::RBDataSerialization::RBEvaluationSerialization::write_to_file().

◆ add_rb_scm_evaluation_data_to_builder()

void libMesh::RBDataSerialization::add_rb_scm_evaluation_data_to_builder ( RBSCMEvaluation rb_scm_eval,
RBData::RBSCMEvaluation::Builder &  rb_scm_eval_builder 
)

Add data for an RBSCMEvaluation to the builder.

Unlike the other functions above, this does not need to be templated because an RBSCMEvaluation only stores Real values, and hence doesn't depend on whether we're using complex numbers or not.

Definition at line 666 of file rb_data_serialization.C.

668 {
669  auto parameter_ranges_list =
670  rb_scm_eval_builder.initParameterRanges();
671  auto discrete_parameters_list =
672  rb_scm_eval_builder.initDiscreteParameters();
674  parameter_ranges_list,
675  discrete_parameters_list);
676 
677  {
678  if (rb_scm_eval.B_min.size() != rb_scm_eval.get_rb_theta_expansion().get_n_A_terms())
679  libmesh_error_msg("Size error while writing B_min");
680  auto b_min_list = rb_scm_eval_builder.initBMin( rb_scm_eval.B_min.size() );
681  for (auto i : index_range(rb_scm_eval.B_min))
682  b_min_list.set(i, rb_scm_eval.get_B_min(i));
683  }
684 
685  {
686  if (rb_scm_eval.B_max.size() != rb_scm_eval.get_rb_theta_expansion().get_n_A_terms())
687  libmesh_error_msg("Size error while writing B_max");
688 
689  auto b_max_list = rb_scm_eval_builder.initBMax( rb_scm_eval.B_max.size() );
690  for (auto i : index_range(rb_scm_eval.B_max))
691  b_max_list.set(i, rb_scm_eval.get_B_max(i));
692  }
693 
694  {
695  auto cj_stability_vector =
696  rb_scm_eval_builder.initCJStabilityVector( rb_scm_eval.C_J_stability_vector.size() );
697  for (auto i : index_range(rb_scm_eval.C_J_stability_vector))
698  cj_stability_vector.set(i, rb_scm_eval.get_C_J_stability_constraint(i));
699  }
700 
701  {
702  auto cj_parameters_outer =
703  rb_scm_eval_builder.initCJ( rb_scm_eval.C_J.size() );
704 
705  for (auto i : index_range(rb_scm_eval.C_J))
706  {
707  auto cj_parameters_inner =
708  cj_parameters_outer.init(i, rb_scm_eval.C_J[i].n_parameters());
709 
710  unsigned int count = 0;
711  for (const auto & pr : rb_scm_eval.C_J[i])
712  {
713  cj_parameters_inner[count].setName( pr.first );
714  cj_parameters_inner[count].setValue( pr.second );
715  count++;
716  }
717 
718  }
719  }
720 
721  {
722  unsigned int n_C_J_values = rb_scm_eval.C_J.size();
723  unsigned int n_A_terms = rb_scm_eval.get_rb_theta_expansion().get_n_A_terms();
724  unsigned int n_values = n_C_J_values*n_A_terms;
725  auto scm_ub_vectors =
726  rb_scm_eval_builder.initScmUbVectors( n_values );
727 
728  for (unsigned int i=0; i<n_C_J_values; i++)
729  for (unsigned int j=0; j<n_A_terms; j++)
730  {
731  unsigned int offset = i*n_A_terms + j;
732  scm_ub_vectors.set(offset, rb_scm_eval.get_SCM_UB_vector(i,j));
733  }
734  }
735 }

References add_parameter_ranges_to_builder(), libMesh::RBSCMEvaluation::B_max, libMesh::RBSCMEvaluation::B_min, libMesh::RBSCMEvaluation::C_J, libMesh::RBSCMEvaluation::C_J_stability_vector, libMesh::RBSCMEvaluation::get_B_max(), libMesh::RBSCMEvaluation::get_B_min(), libMesh::RBSCMEvaluation::get_C_J_stability_constraint(), libMesh::RBThetaExpansion::get_n_A_terms(), libMesh::RBSCMEvaluation::get_rb_theta_expansion(), libMesh::RBSCMEvaluation::get_SCM_UB_vector(), and libMesh::index_range().

Referenced by libMesh::RBDataSerialization::RBSCMEvaluationSerialization::write_to_file().

◆ add_transient_rb_evaluation_data_to_builder()

template<typename RBEvaluationBuilderNumber , typename TransRBEvaluationBuilderNumber >
void libMesh::RBDataSerialization::add_transient_rb_evaluation_data_to_builder ( TransientRBEvaluation trans_rb_eval,
RBEvaluationBuilderNumber &  rb_eval_builder,
TransRBEvaluationBuilderNumber &  trans_rb_eval_builder 
)

Add data for a TransientRBEvaluation to the builder.

Templated to deal with both Real and Complex numbers.

Definition at line 474 of file rb_data_serialization.C.

477 {
478  add_rb_evaluation_data_to_builder(trans_rb_eval, rb_eval_builder);
479 
480  trans_rb_eval_builder.setDeltaT(trans_rb_eval.get_delta_t());
481  trans_rb_eval_builder.setEulerTheta(trans_rb_eval.get_euler_theta());
482  trans_rb_eval_builder.setNTimeSteps(trans_rb_eval.get_n_time_steps());
483  trans_rb_eval_builder.setTimeStep(trans_rb_eval.get_time_step());
484 
485  unsigned int n_bfs = trans_rb_eval.get_n_basis_functions();
486 
487  // L2-inner-product matrix
488  {
489  auto rb_L2_matrix_list =
490  trans_rb_eval_builder.initRbL2Matrix(n_bfs*n_bfs);
491 
492  for (unsigned int i=0; i<n_bfs; ++i)
493  for (unsigned int j=0; j<n_bfs; ++j)
494  {
495  unsigned int offset = i*n_bfs + j;
496  set_scalar_in_list(rb_L2_matrix_list,
497  offset,
498  trans_rb_eval.RB_L2_matrix(i,j));
499  }
500  }
501 
502  TransientRBThetaExpansion & trans_theta_expansion =
503  cast_ref<TransientRBThetaExpansion &>(trans_rb_eval.get_rb_theta_expansion());
504  unsigned int n_M_terms = trans_theta_expansion.get_n_M_terms();
505  // Mq matrices
506  {
507  auto rb_Mq_matrices_outer_list = trans_rb_eval_builder.initRbMqMatrices(n_M_terms);
508  for (unsigned int q_m=0; q_m < n_M_terms; ++q_m)
509  {
510  auto rb_Mq_matrices_inner_list = rb_Mq_matrices_outer_list.init(q_m, n_bfs*n_bfs);
511  for (unsigned int i=0; i < n_bfs; ++i)
512  for (unsigned int j=0; j < n_bfs; ++j)
513  {
514  unsigned int offset = i*n_bfs+j;
515  set_scalar_in_list(rb_Mq_matrices_inner_list,
516  offset,
517  trans_rb_eval.RB_M_q_vector[q_m](i,j));
518  }
519  }
520  }
521 
522  // The initial condition and L2 error at t=0.
523  // We store the values for each RB space of dimension (0,...,n_basis_functions).
524  {
525  auto initial_l2_errors_builder =
526  trans_rb_eval_builder.initInitialL2Errors(n_bfs);
527  auto initial_conditions_outer_list =
528  trans_rb_eval_builder.initInitialConditions(n_bfs);
529 
530  for (unsigned int i=0; i<n_bfs; i++)
531  {
532  initial_l2_errors_builder.set(i, trans_rb_eval.initial_L2_error_all_N[i]);
533 
534  auto initial_conditions_inner_list =
535  initial_conditions_outer_list.init(i, i+1);
536  for (unsigned int j=0; j<=i; j++)
537  {
538  set_scalar_in_list(initial_conditions_inner_list,
539  j,
540  trans_rb_eval.RB_initial_condition_all_N[i](j));
541  }
542  }
543  }
544 
545  // FqMq representor inner-product data
546  {
547  unsigned int n_F_terms = trans_theta_expansion.get_n_F_terms();
548  auto fq_mq_innerprods_list =
549  trans_rb_eval_builder.initFqMqInnerprods(n_F_terms*n_M_terms*n_bfs);
550 
551  for (unsigned int q_f=0; q_f<n_F_terms; ++q_f)
552  for (unsigned int q_m=0; q_m<n_M_terms; ++q_m)
553  for (unsigned int i=0; i<n_bfs; ++i)
554  {
555  unsigned int offset = q_f*n_M_terms*n_bfs + q_m*n_bfs + i;
556  set_scalar_in_list(fq_mq_innerprods_list,
557  offset,
558  trans_rb_eval.Fq_Mq_representor_innerprods[q_f][q_m][i]);
559  }
560  }
561 
562  // MqMq representor inner-product data
563  {
564  unsigned int Q_m_hat = n_M_terms*(n_M_terms+1)/2;
565  auto mq_mq_innerprods_list =
566  trans_rb_eval_builder.initMqMqInnerprods(Q_m_hat*n_bfs*n_bfs);
567 
568  for (unsigned int i=0; i < Q_m_hat; ++i)
569  for (unsigned int j=0; j < n_bfs; ++j)
570  for (unsigned int l=0; l < n_bfs; ++l)
571  {
572  unsigned int offset = i*n_bfs*n_bfs + j*n_bfs + l;
573  set_scalar_in_list(mq_mq_innerprods_list,
574  offset,
575  trans_rb_eval.Mq_Mq_representor_innerprods[i][j][l]);
576  }
577  }
578 
579  // AqMq representor inner-product data
580  {
581  unsigned int n_A_terms = trans_theta_expansion.get_n_A_terms();
582 
583  auto aq_mq_innerprods_list =
584  trans_rb_eval_builder.initAqMqInnerprods(n_A_terms*n_M_terms*n_bfs*n_bfs);
585 
586  for (unsigned int q_a=0; q_a<n_A_terms; q_a++)
587  for (unsigned int q_m=0; q_m<n_M_terms; q_m++)
588  for (unsigned int i=0; i<n_bfs; i++)
589  for (unsigned int j=0; j<n_bfs; j++)
590  {
591  unsigned int offset =
592  q_a*(n_M_terms*n_bfs*n_bfs) + q_m*(n_bfs*n_bfs) + i*n_bfs + j;
593  set_scalar_in_list(aq_mq_innerprods_list,
594  offset,
595  trans_rb_eval.Aq_Mq_representor_innerprods[q_a][q_m][i][j]);
596  }
597  }
598 
599 }

References add_rb_evaluation_data_to_builder(), libMesh::TransientRBEvaluation::Aq_Mq_representor_innerprods, libMesh::TransientRBEvaluation::Fq_Mq_representor_innerprods, libMesh::RBTemporalDiscretization::get_delta_t(), libMesh::RBTemporalDiscretization::get_euler_theta(), libMesh::RBThetaExpansion::get_n_A_terms(), libMesh::RBEvaluation::get_n_basis_functions(), libMesh::RBThetaExpansion::get_n_F_terms(), libMesh::TransientRBThetaExpansion::get_n_M_terms(), libMesh::RBTemporalDiscretization::get_n_time_steps(), libMesh::RBEvaluation::get_rb_theta_expansion(), libMesh::RBTemporalDiscretization::get_time_step(), libMesh::TransientRBEvaluation::initial_L2_error_all_N, libMesh::TransientRBEvaluation::Mq_Mq_representor_innerprods, libMesh::TransientRBEvaluation::RB_initial_condition_all_N, libMesh::TransientRBEvaluation::RB_L2_matrix, and libMesh::TransientRBEvaluation::RB_M_q_vector.

Referenced by libMesh::RBDataSerialization::TransientRBEvaluationSerialization::write_to_file().

libMesh::RBDataSerialization::add_point_to_builder
void add_point_to_builder(const Point &point, RBData::Point3D::Builder point_builder)
Helper function that adds point data.
Definition: rb_data_serialization.C:738
libMesh::RBDataSerialization::add_elem_to_builder
void add_elem_to_builder(const Elem &elem, RBData::MeshElem::Builder mesh_elem_builder)
Helper function that adds element data.
Definition: rb_data_serialization.C:749
libMesh::index_range
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:106
libMesh::Utility::enum_to_string
std::string enum_to_string(const T e)
libMesh::RBDataSerialization::add_rb_evaluation_data_to_builder
void add_rb_evaluation_data_to_builder(RBEvaluation &rb_eval, RBEvaluationBuilderNumber &rb_eval_builder)
Add data for an RBEvaluation to the builder.
Definition: rb_data_serialization.C:328
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::RBDataSerialization::add_parameter_ranges_to_builder
void add_parameter_ranges_to_builder(const RBParametrized &rb_evaluation, RBData::ParameterRanges::Builder &parameter_ranges, RBData::DiscreteParameterList::Builder &discrete_parameters_list)
Add parameter ranges for continuous and discrete parameters.
Definition: rb_data_serialization.C:263