libMesh
Classes | Functions
libMesh::RBDataDeserialization Namespace Reference

Classes

class  RBEIMEvaluationDeserialization
 This class de-serializes a RBEIMEvaluation object using the Cap'n Proto library. More...
 
class  RBEvaluationDeserialization
 This class de-serializes an RBEvaluation object using the Cap'n Proto library. More...
 
class  RBSCMEvaluationDeserialization
 This class de-serializes a RBSCMEvaluation object using the Cap'n Proto library. More...
 
class  TransientRBEvaluationDeserialization
 This class de-serializes a TransientRBEvaluation object using the Cap'n Proto library. More...
 

Functions

void load_parameter_ranges (RBParametrized &rb_evaluation, RBData::ParameterRanges::Reader &parameter_ranges, RBData::DiscreteParameterList::Reader &discrete_parameters_list)
 Load parameter ranges and discrete parameter values into an RBEvaluation from the corresponding structure in the buffer. More...
 
template<typename RBEvaluationReaderNumber >
void load_rb_evaluation_data (RBEvaluation &rb_evaluation, RBEvaluationReaderNumber &rb_evaluation_reader, bool read_error_bound_data)
 Load an RB evaluation from a corresponding reader structure in the buffer. More...
 
template<typename RBEvaluationReaderNumber , typename TransRBEvaluationReaderNumber >
void load_transient_rb_evaluation_data (TransientRBEvaluation &trans_rb_eval, RBEvaluationReaderNumber &rb_evaluation_reader, TransRBEvaluationReaderNumber &trans_rb_eval_reader, bool read_error_bound_data)
 Load an RB evaluation from a corresponding reader structure in the buffer. More...
 
template<typename RBEvaluationReaderNumber , typename RBEIMEvaluationReaderNumber >
void load_rb_eim_evaluation_data (RBEIMEvaluation &rb_eim_eval, RBEvaluationReaderNumber &rb_evaluation_reader, RBEIMEvaluationReaderNumber &rb_eim_eval_reader)
 Load an EIM RB evaluation from a corresponding reader structure in the buffer. More...
 
void load_rb_scm_evaluation_data (RBSCMEvaluation &rb_scm_eval, RBData::RBSCMEvaluation::Reader &rb_scm_eval_reader)
 Load an SCM RB evaluation from a corresponding reader structure in the buffer. More...
 
void load_point (RBData::Point3D::Reader point_reader, Point &point)
 Helper function that loads point data. More...
 
void load_elem_into_mesh (RBData::MeshElem::Reader mesh_elem_reader, libMesh::Elem *elem, libMesh::ReplicatedMesh &mesh)
 Helper function that loads element data. More...
 

Function Documentation

◆ load_elem_into_mesh()

void libMesh::RBDataDeserialization::load_elem_into_mesh ( RBData::MeshElem::Reader  mesh_elem_reader,
libMesh::Elem elem,
libMesh::ReplicatedMesh mesh 
)

Helper function that loads element data.

Definition at line 871 of file rb_data_deserialization.C.

874 {
875  auto mesh_elem_point_list = mesh_elem_reader.getPoints();
876  unsigned int n_points = mesh_elem_point_list.size();
877 
878  if (n_points != elem->n_nodes())
879  libmesh_error_msg("Wrong number of nodes for element type");
880 
881  for (unsigned int i=0; i < n_points; ++i)
882  {
883  libMesh::Node * node = new libMesh::Node(mesh_elem_point_list[i].getX(),
884  mesh_elem_point_list[i].getY(),
885  mesh_elem_point_list[i].getZ());
886 
887  mesh.add_node(node);
888 
889  elem->set_node(i) = node;
890  }
891 
892  elem->subdomain_id() = mesh_elem_reader.getSubdomainId();
893 
894  mesh.add_elem(elem);
895 }

References mesh, libMesh::Elem::n_nodes(), libMesh::Elem::set_node(), and libMesh::Elem::subdomain_id().

Referenced by load_rb_eim_evaluation_data().

◆ load_parameter_ranges()

void libMesh::RBDataDeserialization::load_parameter_ranges ( RBParametrized rb_evaluation,
RBData::ParameterRanges::Reader &  parameter_ranges,
RBData::DiscreteParameterList::Reader &  discrete_parameters_list 
)

Load parameter ranges and discrete parameter values into an RBEvaluation from the corresponding structure in the buffer.

Definition at line 287 of file rb_data_deserialization.C.

290 {
291  // Continuous parameters
292  RBParameters parameters_min;
293  RBParameters parameters_max;
294  {
295  unsigned int n_parameter_ranges = parameter_ranges.getNames().size();
296 
297  for (unsigned int i=0; i<n_parameter_ranges; ++i)
298  {
299  std::string parameter_name = parameter_ranges.getNames()[i];
300  Real min_value = parameter_ranges.getMinValues()[i];
301  Real max_value = parameter_ranges.getMaxValues()[i];
302 
303  parameters_min.set_value(parameter_name, min_value);
304  parameters_max.set_value(parameter_name, max_value);
305  }
306  }
307 
308  // Discrete parameters
309  std::map<std::string, std::vector<Real>> discrete_parameter_values;
310  {
311  unsigned int n_discrete_parameters = discrete_parameters_list.getNames().size();
312 
313  for (unsigned int i=0; i<n_discrete_parameters; ++i)
314  {
315  std::string parameter_name = discrete_parameters_list.getNames()[i];
316 
317  auto value_list = discrete_parameters_list.getValues()[i];
318  unsigned int n_values = value_list.size();
319  std::vector<Real> values(n_values);
320  for (unsigned int j=0; j<n_values; ++j)
321  values[j] = value_list[j];
322 
323  discrete_parameter_values[parameter_name] = values;
324  }
325  }
326 
327  rb_evaluation.initialize_parameters(parameters_min,
328  parameters_max,
329  discrete_parameter_values);
330 }

References libMesh::RBParametrized::initialize_parameters(), libMesh::Real, and libMesh::RBParameters::set_value().

Referenced by load_rb_evaluation_data(), and load_rb_scm_evaluation_data().

◆ load_point()

void libMesh::RBDataDeserialization::load_point ( RBData::Point3D::Reader  point_reader,
Point point 
)

Helper function that loads point data.

Definition at line 859 of file rb_data_deserialization.C.

860 {
861  point(0) = point_reader.getX();
862 
863  if (LIBMESH_DIM >= 2)
864  point(1) = point_reader.getY();
865 
866  if (LIBMESH_DIM >= 3)
867  point(2) = point_reader.getZ();
868 }

Referenced by load_rb_eim_evaluation_data().

◆ load_rb_eim_evaluation_data()

template<typename RBEvaluationReaderNumber , typename RBEIMEvaluationReaderNumber >
void libMesh::RBDataDeserialization::load_rb_eim_evaluation_data ( RBEIMEvaluation rb_eim_eval,
RBEvaluationReaderNumber &  rb_evaluation_reader,
RBEIMEvaluationReaderNumber &  rb_eim_eval_reader 
)

Load an EIM RB evaluation from a corresponding reader structure in the buffer.

Templated to deal with both Real and Complex numbers.

Definition at line 675 of file rb_data_deserialization.C.

678 {
679  // We use read_error_bound_data=false here, since the RBEvaluation error bound data
680  // is not relevant to EIM.
681  load_rb_evaluation_data(rb_eim_evaluation,
682  rb_evaluation_reader,
683  /*read_error_bound_data*/ false);
684 
685  unsigned int n_bfs = rb_eim_evaluation.get_n_basis_functions();
686 
687  // EIM interpolation matrix
688  {
689  auto interpolation_matrix_list = rb_eim_evaluation_reader.getInterpolationMatrix();
690 
691  if (interpolation_matrix_list.size() != n_bfs*(n_bfs+1)/2)
692  libmesh_error_msg("Size error while reading the eim inner product matrix.");
693 
694  for (unsigned int i=0; i<n_bfs; ++i)
695  for (unsigned int j=0; j<=i; ++j)
696  {
697  unsigned int offset = i*(i+1)/2 + j;
698  rb_eim_evaluation.interpolation_matrix(i,j) =
699  load_scalar_value(interpolation_matrix_list[offset]);
700  }
701  }
702 
703  // Interpolation points
704  {
705  auto interpolation_points_list =
706  rb_eim_evaluation_reader.getInterpolationPoints();
707 
708  if (interpolation_points_list.size() != n_bfs)
709  libmesh_error_msg("Size error while reading the eim interpolation points.");
710 
711  rb_eim_evaluation.interpolation_points.resize(n_bfs);
712  for (unsigned int i=0; i<n_bfs; ++i)
713  {
714  load_point(interpolation_points_list[i],
715  rb_eim_evaluation.interpolation_points[i]);
716  }
717  }
718 
719  // Interpolation points variables
720  {
721  auto interpolation_points_var_list =
722  rb_eim_evaluation_reader.getInterpolationPointsVar();
723  rb_eim_evaluation.interpolation_points_var.resize(n_bfs);
724 
725  if (interpolation_points_var_list.size() != n_bfs)
726  libmesh_error_msg("Size error while reading the eim interpolation variables.");
727 
728  for (unsigned int i=0; i<n_bfs; ++i)
729  {
730  rb_eim_evaluation.interpolation_points_var[i] =
731  interpolation_points_var_list[i];
732  }
733  }
734 
735  // Interpolation elements
736  {
737  libMesh::dof_id_type elem_id = 0;
738  libMesh::ReplicatedMesh & interpolation_points_mesh =
739  rb_eim_evaluation.get_interpolation_points_mesh();
740  interpolation_points_mesh.clear();
741 
742  auto interpolation_points_elem_list =
743  rb_eim_evaluation_reader.getInterpolationPointsElems();
744  unsigned int n_interpolation_elems = interpolation_points_elem_list.size();
745 
746  if (n_interpolation_elems != n_bfs)
747  libmesh_error_msg("The number of elements should match the number of basis functions");
748 
749  rb_eim_evaluation.interpolation_points_elem.resize(n_interpolation_elems);
750 
751  for (unsigned int i=0; i<n_interpolation_elems; ++i)
752  {
753  auto mesh_elem_reader = interpolation_points_elem_list[i];
754  std::string elem_type_name = mesh_elem_reader.getType().cStr();
755  libMesh::ElemType elem_type =
756  libMesh::Utility::string_to_enum<libMesh::ElemType>(elem_type_name);
757 
758  libMesh::Elem * elem = libMesh::Elem::build(elem_type).release();
759  elem->set_id(elem_id++);
760  load_elem_into_mesh(mesh_elem_reader, elem, interpolation_points_mesh);
761 
762  rb_eim_evaluation.interpolation_points_elem[i] = elem;
763  }
764  }
765 }

References libMesh::Elem::build(), libMesh::ReplicatedMesh::clear(), libMesh::RBEIMEvaluation::get_interpolation_points_mesh(), libMesh::RBEvaluation::get_n_basis_functions(), libMesh::RBEIMEvaluation::interpolation_matrix, libMesh::RBEIMEvaluation::interpolation_points, libMesh::RBEIMEvaluation::interpolation_points_elem, libMesh::RBEIMEvaluation::interpolation_points_var, load_elem_into_mesh(), load_point(), load_rb_evaluation_data(), and libMesh::DofObject::set_id().

Referenced by libMesh::RBDataDeserialization::RBEIMEvaluationDeserialization::read_from_file().

◆ load_rb_evaluation_data()

template<typename RBEvaluationReaderNumber >
void libMesh::RBDataDeserialization::load_rb_evaluation_data ( RBEvaluation rb_evaluation,
RBEvaluationReaderNumber &  rb_evaluation_reader,
bool  read_error_bound_data 
)

Load an RB evaluation from a corresponding reader structure in the buffer.

Definition at line 333 of file rb_data_deserialization.C.

336 {
337  // Set number of basis functions
338  unsigned int n_bfs = rb_evaluation_reader.getNBfs();
339  rb_evaluation.set_n_basis_functions(n_bfs);
340 
341  rb_evaluation.resize_data_structures(n_bfs, read_error_bound_data);
342 
343  auto parameter_ranges =
344  rb_evaluation_reader.getParameterRanges();
345  auto discrete_parameters_list =
346  rb_evaluation_reader.getDiscreteParameters();
347 
348  load_parameter_ranges(rb_evaluation,
349  parameter_ranges,
350  discrete_parameters_list);
351 
352  const RBThetaExpansion & rb_theta_expansion = rb_evaluation.get_rb_theta_expansion();
353 
354  unsigned int n_F_terms = rb_theta_expansion.get_n_F_terms();
355  unsigned int n_A_terms = rb_theta_expansion.get_n_A_terms();
356 
357  if (read_error_bound_data)
358  {
359 
360  // Fq representor inner-product data
361  {
362  unsigned int Q_f_hat = n_F_terms*(n_F_terms+1)/2;
363 
364  auto fq_innerprods_list = rb_evaluation_reader.getFqInnerprods();
365  if (fq_innerprods_list.size() != Q_f_hat)
366  libmesh_error_msg("Size error while reading Fq representor norm data from buffer.");
367 
368  for (unsigned int i=0; i < Q_f_hat; ++i)
369  rb_evaluation.Fq_representor_innerprods[i] = load_scalar_value(fq_innerprods_list[i]);
370  }
371 
372  // Fq_Aq representor inner-product data
373  {
374  auto fq_aq_innerprods_list = rb_evaluation_reader.getFqAqInnerprods();
375  if (fq_aq_innerprods_list.size() != n_F_terms*n_A_terms*n_bfs)
376  libmesh_error_msg("Size error while reading Fq-Aq representor norm data from buffer.");
377 
378  for (unsigned int q_f=0; q_f<n_F_terms; ++q_f)
379  for (unsigned int q_a=0; q_a<n_A_terms; ++q_a)
380  for (unsigned int i=0; i<n_bfs; ++i)
381  {
382  unsigned int offset = q_f*n_A_terms*n_bfs + q_a*n_bfs + i;
383  rb_evaluation.Fq_Aq_representor_innerprods[q_f][q_a][i] =
384  load_scalar_value(fq_aq_innerprods_list[offset]);
385  }
386  }
387 
388  // Aq_Aq representor inner-product data
389  {
390  unsigned int Q_a_hat = n_A_terms*(n_A_terms+1)/2;
391  auto aq_aq_innerprods_list = rb_evaluation_reader.getAqAqInnerprods();
392  if (aq_aq_innerprods_list.size() != Q_a_hat*n_bfs*n_bfs)
393  libmesh_error_msg("Size error while reading Aq-Aq representor norm data from buffer.");
394 
395  for (unsigned int i=0; i<Q_a_hat; ++i)
396  for (unsigned int j=0; j<n_bfs; ++j)
397  for (unsigned int l=0; l<n_bfs; ++l)
398  {
399  unsigned int offset = i*n_bfs*n_bfs + j*n_bfs + l;
400  rb_evaluation.Aq_Aq_representor_innerprods[i][j][l] =
401  load_scalar_value(aq_aq_innerprods_list[offset]);
402  }
403  }
404 
405  // Output dual inner-product data
406  {
407  unsigned int n_outputs = rb_theta_expansion.get_n_outputs();
408  auto output_innerprod_outer = rb_evaluation_reader.getOutputDualInnerprods();
409 
410  if (output_innerprod_outer.size() != n_outputs)
411  libmesh_error_msg("Incorrect number of outputs detected in the buffer");
412 
413  for (unsigned int output_id=0; output_id<n_outputs; ++output_id)
414  {
415  unsigned int n_output_terms = rb_theta_expansion.get_n_output_terms(output_id);
416 
417  unsigned int Q_l_hat = n_output_terms*(n_output_terms+1)/2;
418  auto output_innerprod_inner = output_innerprod_outer[output_id];
419 
420  if (output_innerprod_inner.size() != Q_l_hat)
421  libmesh_error_msg("Incorrect number of output terms detected in the buffer");
422 
423  for (unsigned int q=0; q<Q_l_hat; ++q)
424  {
425  rb_evaluation.output_dual_innerprods[output_id][q] =
426  load_scalar_value(output_innerprod_inner[q]);
427  }
428  }
429  }
430  }
431 
432  // Output vectors
433  {
434  unsigned int n_outputs = rb_theta_expansion.get_n_outputs();
435  auto output_vector_outer = rb_evaluation_reader.getOutputVectors();
436 
437  if (output_vector_outer.size() != n_outputs)
438  libmesh_error_msg("Incorrect number of outputs detected in the buffer");
439 
440  for (unsigned int output_id=0; output_id<n_outputs; ++output_id)
441  {
442  unsigned int n_output_terms = rb_theta_expansion.get_n_output_terms(output_id);
443 
444  auto output_vector_middle = output_vector_outer[output_id];
445  if (output_vector_middle.size() != n_output_terms)
446  libmesh_error_msg("Incorrect number of output terms detected in the buffer");
447 
448  for (unsigned int q_l=0; q_l<n_output_terms; ++q_l)
449  {
450  auto output_vectors_inner_list = output_vector_middle[q_l];
451 
452  if (output_vectors_inner_list.size() != n_bfs)
453  libmesh_error_msg("Incorrect number of output terms detected in the buffer");
454 
455  for (unsigned int j=0; j<n_bfs; ++j)
456  {
457  rb_evaluation.RB_output_vectors[output_id][q_l](j) =
458  load_scalar_value(output_vectors_inner_list[j]);
459  }
460  }
461  }
462  }
463 
464  // Fq vectors and Aq matrices
465  {
466  auto rb_fq_vectors_outer_list = rb_evaluation_reader.getRbFqVectors();
467  if (rb_fq_vectors_outer_list.size() != n_F_terms)
468  libmesh_error_msg("Incorrect number of Fq vectors detected in the buffer");
469 
470  for (unsigned int q_f=0; q_f<n_F_terms; ++q_f)
471  {
472  auto rb_fq_vectors_inner_list = rb_fq_vectors_outer_list[q_f];
473  if (rb_fq_vectors_inner_list.size() != n_bfs)
474  libmesh_error_msg("Incorrect Fq vector size detected in the buffer");
475 
476  for (unsigned int i=0; i < n_bfs; ++i)
477  {
478  rb_evaluation.RB_Fq_vector[q_f](i) =
479  load_scalar_value(rb_fq_vectors_inner_list[i]);
480  }
481  }
482 
483  auto rb_Aq_matrices_outer_list = rb_evaluation_reader.getRbAqMatrices();
484  if (rb_Aq_matrices_outer_list.size() != n_A_terms)
485  libmesh_error_msg("Incorrect number of Aq matrices detected in the buffer");
486 
487  for (unsigned int q_a=0; q_a<n_A_terms; ++q_a)
488  {
489  auto rb_Aq_matrices_inner_list = rb_Aq_matrices_outer_list[q_a];
490  if (rb_Aq_matrices_inner_list.size() != n_bfs*n_bfs)
491  libmesh_error_msg("Incorrect Aq matrix size detected in the buffer");
492 
493  for (unsigned int i=0; i<n_bfs; ++i)
494  for (unsigned int j=0; j<n_bfs; ++j)
495  {
496  unsigned int offset = i*n_bfs+j;
497  rb_evaluation.RB_Aq_vector[q_a](i,j) =
498  load_scalar_value(rb_Aq_matrices_inner_list[offset]);
499  }
500  }
501  }
502 
503  // Inner-product matrix
504  if (rb_evaluation.compute_RB_inner_product)
505  {
506  auto rb_inner_product_matrix_list =
507  rb_evaluation_reader.getRbInnerProductMatrix();
508 
509  if (rb_inner_product_matrix_list.size() != n_bfs*n_bfs)
510  libmesh_error_msg("Size error while reading the inner product matrix.");
511 
512  for (unsigned int i=0; i<n_bfs; ++i)
513  for (unsigned int j=0; j<n_bfs; ++j)
514  {
515  unsigned int offset = i*n_bfs + j;
516  rb_evaluation.RB_inner_product_matrix(i,j) =
517  load_scalar_value(rb_inner_product_matrix_list[offset]);
518  }
519  }
520 }

References 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::RBThetaExpansion::get_n_F_terms(), libMesh::RBThetaExpansion::get_n_output_terms(), libMesh::RBThetaExpansion::get_n_outputs(), libMesh::RBEvaluation::get_rb_theta_expansion(), load_parameter_ranges(), libMesh::RBEvaluation::output_dual_innerprods, libMesh::RBEvaluation::RB_Aq_vector, libMesh::RBEvaluation::RB_Fq_vector, libMesh::RBEvaluation::RB_inner_product_matrix, libMesh::RBEvaluation::RB_output_vectors, libMesh::RBEvaluation::resize_data_structures(), and libMesh::RBEvaluation::set_n_basis_functions().

Referenced by load_rb_eim_evaluation_data(), load_transient_rb_evaluation_data(), and libMesh::RBDataDeserialization::RBEvaluationDeserialization::read_from_file().

◆ load_rb_scm_evaluation_data()

void libMesh::RBDataDeserialization::load_rb_scm_evaluation_data ( RBSCMEvaluation rb_scm_eval,
RBData::RBSCMEvaluation::Reader &  rb_scm_eval_reader 
)

Load an SCM RB evaluation from a corresponding reader structure in the buffer.

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 768 of file rb_data_deserialization.C.

770 {
771  auto parameter_ranges =
772  rb_scm_evaluation_reader.getParameterRanges();
773  auto discrete_parameters_list =
774  rb_scm_evaluation_reader.getDiscreteParameters();
775  load_parameter_ranges(rb_scm_eval,
776  parameter_ranges,
777  discrete_parameters_list);
778 
779  unsigned int n_A_terms = rb_scm_eval.get_rb_theta_expansion().get_n_A_terms();
780 
781  {
782  auto b_min_list = rb_scm_evaluation_reader.getBMin();
783 
784  if (b_min_list.size() != n_A_terms)
785  libmesh_error_msg("Size error while reading B_min");
786 
787  rb_scm_eval.B_min.clear();
788  for (unsigned int i=0; i<n_A_terms; i++)
789  rb_scm_eval.B_min.push_back(b_min_list[i]);
790  }
791 
792  {
793  auto b_max_list = rb_scm_evaluation_reader.getBMax();
794 
795  if (b_max_list.size() != n_A_terms)
796  libmesh_error_msg("Size error while reading B_max");
797 
798  rb_scm_eval.B_max.clear();
799  for (unsigned int i=0; i<n_A_terms; i++)
800  rb_scm_eval.B_max.push_back(b_max_list[i]);
801  }
802 
803  {
804  auto cJ_stability_vector =
805  rb_scm_evaluation_reader.getCJStabilityVector();
806 
807  rb_scm_eval.C_J_stability_vector.clear();
808  for (const auto & val : cJ_stability_vector)
809  rb_scm_eval.C_J_stability_vector.push_back(val);
810  }
811 
812  {
813  auto cJ_parameters_outer =
814  rb_scm_evaluation_reader.getCJ();
815 
816  rb_scm_eval.C_J.resize(cJ_parameters_outer.size());
817  for (auto i : index_range(cJ_parameters_outer))
818  {
819  auto cJ_parameters_inner =
820  cJ_parameters_outer[i];
821 
822  for (auto j : index_range(cJ_parameters_inner))
823  {
824  std::string param_name = cJ_parameters_inner[j].getName();
825  Real param_value = cJ_parameters_inner[j].getValue();
826  rb_scm_eval.C_J[i].set_value(param_name, param_value);
827  }
828  }
829  }
830 
831  {
832  auto scm_ub_vectors =
833  rb_scm_evaluation_reader.getScmUbVectors();
834 
835  // The number of UB vectors is the same as the number of C_J values
836  unsigned int n_C_J_values = rb_scm_evaluation_reader.getCJ().size();
837 
838  if (scm_ub_vectors.size() != n_C_J_values*n_A_terms)
839  libmesh_error_msg("Size mismatch in SCB UB vectors");
840 
841  rb_scm_eval.SCM_UB_vectors.resize( n_C_J_values );
842  for (unsigned int i=0; i<n_C_J_values; i++)
843  {
844  rb_scm_eval.SCM_UB_vectors[i].resize(n_A_terms);
845  for (unsigned int j=0; j<n_A_terms; j++)
846  {
847  unsigned int offset = i*n_A_terms + j;
848  rb_scm_eval.SCM_UB_vectors[i][j] = scm_ub_vectors[offset];
849 
850  }
851  }
852  }
853 
854 }

References libMesh::RBSCMEvaluation::B_max, libMesh::RBSCMEvaluation::B_min, libMesh::RBSCMEvaluation::C_J, libMesh::RBSCMEvaluation::C_J_stability_vector, libMesh::RBThetaExpansion::get_n_A_terms(), libMesh::RBSCMEvaluation::get_rb_theta_expansion(), libMesh::index_range(), load_parameter_ranges(), libMesh::Real, and libMesh::RBSCMEvaluation::SCM_UB_vectors.

Referenced by libMesh::RBDataDeserialization::RBSCMEvaluationDeserialization::read_from_file().

◆ load_transient_rb_evaluation_data()

template<typename RBEvaluationReaderNumber , typename TransRBEvaluationReaderNumber >
void libMesh::RBDataDeserialization::load_transient_rb_evaluation_data ( TransientRBEvaluation trans_rb_eval,
RBEvaluationReaderNumber &  rb_evaluation_reader,
TransRBEvaluationReaderNumber &  trans_rb_eval_reader,
bool  read_error_bound_data 
)

Load an RB evaluation from a corresponding reader structure in the buffer.

Templated to deal with both Real and Complex numbers.

Definition at line 523 of file rb_data_deserialization.C.

527 {
528  load_rb_evaluation_data(trans_rb_eval,
529  rb_eval_reader,
530  read_error_bound_data);
531 
532  trans_rb_eval.set_delta_t( trans_rb_eval_reader.getDeltaT() );
533  trans_rb_eval.set_euler_theta( trans_rb_eval_reader.getEulerTheta() );
534  trans_rb_eval.set_n_time_steps( trans_rb_eval_reader.getNTimeSteps() );
535  trans_rb_eval.set_time_step( trans_rb_eval_reader.getTimeStep() );
536 
537  unsigned int n_bfs = rb_eval_reader.getNBfs();
538  unsigned int n_F_terms = trans_rb_eval.get_rb_theta_expansion().get_n_F_terms();
539  unsigned int n_A_terms = trans_rb_eval.get_rb_theta_expansion().get_n_A_terms();
540 
541  TransientRBThetaExpansion & trans_theta_expansion =
542  cast_ref<TransientRBThetaExpansion &>(trans_rb_eval.get_rb_theta_expansion());
543  unsigned int n_M_terms = trans_theta_expansion.get_n_M_terms();
544 
545  // L2 matrix
546  {
547  auto rb_L2_matrix_list =
548  trans_rb_eval_reader.getRbL2Matrix();
549 
550  if (rb_L2_matrix_list.size() != n_bfs*n_bfs)
551  libmesh_error_msg("Size error while reading the L2 matrix.");
552 
553  for (unsigned int i=0; i<n_bfs; ++i)
554  for (unsigned int j=0; j<n_bfs; ++j)
555  {
556  unsigned int offset = i*n_bfs + j;
557  trans_rb_eval.RB_L2_matrix(i,j) =
558  load_scalar_value(rb_L2_matrix_list[offset]);
559  }
560  }
561 
562  // Mq matrices
563  {
564  auto rb_Mq_matrices_outer_list = trans_rb_eval_reader.getRbMqMatrices();
565 
566  if (rb_Mq_matrices_outer_list.size() != n_M_terms)
567  libmesh_error_msg("Incorrect number of Mq matrices detected in the buffer");
568 
569  for (unsigned int q_m=0; q_m < n_M_terms; ++q_m)
570  {
571  auto rb_Mq_matrices_inner_list = rb_Mq_matrices_outer_list[q_m];
572  if (rb_Mq_matrices_inner_list.size() != n_bfs*n_bfs)
573  libmesh_error_msg("Incorrect Mq matrix size detected in the buffer");
574 
575  for (unsigned int i=0; i<n_bfs; ++i)
576  for (unsigned int j=0; j<n_bfs; ++j)
577  {
578  unsigned int offset = i*n_bfs+j;
579  trans_rb_eval.RB_M_q_vector[q_m](i,j) =
580  load_scalar_value(rb_Mq_matrices_inner_list[offset]);
581  }
582  }
583  }
584 
585  // The initial condition and L2 error at t=0.
586  {
587  auto initial_l2_errors_reader =
588  trans_rb_eval_reader.getInitialL2Errors();
589  if (initial_l2_errors_reader.size() != n_bfs)
590  libmesh_error_msg("Incorrect number of initial L2 error terms detected in the buffer");
591 
592  auto initial_conditions_outer_list =
593  trans_rb_eval_reader.getInitialConditions();
594  if (initial_conditions_outer_list.size() != n_bfs)
595  libmesh_error_msg("Incorrect number of outer initial conditions detected in the buffer");
596 
597  for (unsigned int i=0; i<n_bfs; i++)
598  {
599  trans_rb_eval.initial_L2_error_all_N[i] =
600  initial_l2_errors_reader[i];
601 
602  auto initial_conditions_inner_list = initial_conditions_outer_list[i];
603  if (initial_conditions_inner_list.size() != (i+1))
604  libmesh_error_msg("Incorrect number of inner initial conditions detected in the buffer");
605 
606  for (unsigned int j=0; j<=i; j++)
607  {
608  trans_rb_eval.RB_initial_condition_all_N[i](j) =
609  load_scalar_value(initial_conditions_inner_list[j]);
610  }
611  }
612  }
613 
614 
615  if (read_error_bound_data)
616  {
617  // Fq_Mq data
618  {
619  auto fq_mq_innerprods_list = trans_rb_eval_reader.getFqMqInnerprods();
620  if (fq_mq_innerprods_list.size() != n_F_terms*n_M_terms*n_bfs)
621  libmesh_error_msg("Size error while reading Fq-Mq representor data from buffer.");
622 
623  for (unsigned int q_f=0; q_f<n_F_terms; ++q_f)
624  for (unsigned int q_m=0; q_m<n_M_terms; ++q_m)
625  for (unsigned int i=0; i<n_bfs; ++i)
626  {
627  unsigned int offset = q_f*n_M_terms*n_bfs + q_m*n_bfs + i;
628  trans_rb_eval.Fq_Mq_representor_innerprods[q_f][q_m][i] =
629  load_scalar_value(fq_mq_innerprods_list[offset]);
630  }
631  }
632 
633 
634  // Mq_Mq representor inner-product data
635  {
636  unsigned int Q_m_hat = n_M_terms*(n_M_terms+1)/2;
637  auto mq_mq_innerprods_list = trans_rb_eval_reader.getMqMqInnerprods();
638  if (mq_mq_innerprods_list.size() != Q_m_hat*n_bfs*n_bfs)
639  libmesh_error_msg("Size error while reading Mq-Mq representor data from buffer.");
640 
641  for (unsigned int i=0; i<Q_m_hat; ++i)
642  for (unsigned int j=0; j<n_bfs; ++j)
643  for (unsigned int l=0; l<n_bfs; ++l)
644  {
645  unsigned int offset = i*n_bfs*n_bfs + j*n_bfs + l;
646  trans_rb_eval.Mq_Mq_representor_innerprods[i][j][l] =
647  load_scalar_value(mq_mq_innerprods_list[offset]);
648  }
649  }
650 
651  // Aq_Mq representor inner-product data
652  {
653  auto aq_mq_innerprods_list =
654  trans_rb_eval_reader.getAqMqInnerprods();
655  if (aq_mq_innerprods_list.size() != n_A_terms*n_M_terms*n_bfs*n_bfs)
656  libmesh_error_msg("Size error while reading Aq-Mq representor data from buffer.");
657 
658  for (unsigned int q_a=0; q_a<n_A_terms; q_a++)
659  for (unsigned int q_m=0; q_m<n_M_terms; q_m++)
660  for (unsigned int i=0; i<n_bfs; i++)
661  for (unsigned int j=0; j<n_bfs; j++)
662  {
663  unsigned int offset =
664  q_a*(n_M_terms*n_bfs*n_bfs) + q_m*(n_bfs*n_bfs) + i*n_bfs + j;
665 
666  trans_rb_eval.Aq_Mq_representor_innerprods[q_a][q_m][i][j] =
667  load_scalar_value(aq_mq_innerprods_list[offset]);
668  }
669  }
670 
671  }
672 }

References libMesh::TransientRBEvaluation::Aq_Mq_representor_innerprods, libMesh::TransientRBEvaluation::Fq_Mq_representor_innerprods, libMesh::RBThetaExpansion::get_n_A_terms(), libMesh::RBThetaExpansion::get_n_F_terms(), libMesh::TransientRBThetaExpansion::get_n_M_terms(), libMesh::RBEvaluation::get_rb_theta_expansion(), libMesh::TransientRBEvaluation::initial_L2_error_all_N, load_rb_evaluation_data(), libMesh::TransientRBEvaluation::Mq_Mq_representor_innerprods, libMesh::TransientRBEvaluation::RB_initial_condition_all_N, libMesh::TransientRBEvaluation::RB_L2_matrix, libMesh::TransientRBEvaluation::RB_M_q_vector, libMesh::RBTemporalDiscretization::set_delta_t(), libMesh::RBTemporalDiscretization::set_euler_theta(), libMesh::RBTemporalDiscretization::set_n_time_steps(), and libMesh::RBTemporalDiscretization::set_time_step().

Referenced by libMesh::RBDataDeserialization::TransientRBEvaluationDeserialization::read_from_file().

libMesh::dof_id_type
uint8_t dof_id_type
Definition: id_types.h:67
libMesh::Elem::n_nodes
virtual unsigned int n_nodes() const =0
libMesh::DofObject::set_id
dof_id_type & set_id()
Definition: dof_object.h:776
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
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::ReplicatedMesh::clear
virtual void clear() override
Clear all internal data.
Definition: replicated_mesh.C:593
libMesh::ReplicatedMesh
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
Definition: replicated_mesh.h:47
libMesh::Node
A Node is like a Point, but with more information.
Definition: node.h:52
libMesh::Elem::set_node
virtual Node *& set_node(const unsigned int i)
Definition: elem.h:2059
libMesh::RBDataDeserialization::load_point
void load_point(RBData::Point3D::Reader point_reader, Point &point)
Helper function that loads point data.
Definition: rb_data_deserialization.C:859
libMesh::RBDataDeserialization::load_elem_into_mesh
void load_elem_into_mesh(RBData::MeshElem::Reader mesh_elem_reader, libMesh::Elem *elem, libMesh::ReplicatedMesh &mesh)
Helper function that loads element data.
Definition: rb_data_deserialization.C:871
libMesh::Elem::subdomain_id
subdomain_id_type subdomain_id() const
Definition: elem.h:2069
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::Elem::build
static std::unique_ptr< Elem > build(const ElemType type, Elem *p=nullptr)
Definition: elem.C:246
libMesh::RBDataDeserialization::load_parameter_ranges
void load_parameter_ranges(RBParametrized &rb_evaluation, RBData::ParameterRanges::Reader &parameter_ranges, RBData::DiscreteParameterList::Reader &discrete_parameters_list)
Load parameter ranges and discrete parameter values into an RBEvaluation from the corresponding struc...
Definition: rb_data_deserialization.C:287
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33
libMesh::RBDataDeserialization::load_rb_evaluation_data
void load_rb_evaluation_data(RBEvaluation &rb_evaluation, RBEvaluationReaderNumber &rb_evaluation_reader, bool read_error_bound_data)
Load an RB evaluation from a corresponding reader structure in the buffer.
Definition: rb_data_deserialization.C:333