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 RBEIMEvaluationReaderNumber >
void load_rb_eim_evaluation_data (RBEIMEvaluation &rb_eim_eval, 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...
 

Function Documentation

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

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

Referenced by load_rb_eim_evaluation_data(), load_rb_evaluation_data(), and load_rb_scm_evaluation_data().

277 {
278  // Continuous parameters
279  RBParameters parameters_min;
280  RBParameters parameters_max;
281 
282  // Discrete parameters
283  std::map<std::string, std::vector<Real>> discrete_parameter_values;
284 
285  // The RBData::ParameterRanges::Reader (CapnProto class) will throw
286  // an exception if there is a problem with the data file (e.g. if it
287  // doesn't exist). We don't use the libmesh_try/catch() macros here
288  // since they don't really allow you to do any processing of the
289  // exception object.
290 #ifdef LIBMESH_ENABLE_EXCEPTIONS
291  try
292 #endif
293  {
294  const auto & parameter_names = parameter_ranges.getNames();
295  for (auto i : make_range(parameter_names.size()))
296  {
297  parameters_min.set_value(parameter_names[i], parameter_ranges.getMinValues()[i]);
298  parameters_max.set_value(parameter_names[i], parameter_ranges.getMaxValues()[i]);
299  }
300 
301 
302  const auto & discrete_names = discrete_parameters_list.getNames();
303  for (auto i : make_range(discrete_names.size()))
304  {
305  const auto & value_list = discrete_parameters_list.getValues()[i];
306  unsigned int n_values = value_list.size();
307  std::vector<Real> values(n_values);
308  for (unsigned int j=0; j<n_values; ++j)
309  values[j] = value_list[j];
310 
311  discrete_parameter_values[discrete_names[i]] = values;
312  }
313  }
314 #ifdef LIBMESH_ENABLE_EXCEPTIONS
315  catch (std::exception & e)
316  {
317  libmesh_error_msg("Error loading parameter ranges from capnp reader.\n"
318  "This usually means that the training data either doesn't exist or is out of date.\n"
319  "Detailed information about the error is below:\n\n" << e.what() << "\n");
320  }
321 #endif
322 
323  rb_evaluation.initialize_parameters(parameters_min,
324  parameters_max,
325  discrete_parameter_values);
326 }
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:134

◆ load_point()

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

Helper function that loads point data.

Definition at line 1102 of file rb_data_deserialization.C.

Referenced by load_rb_eim_evaluation_data().

1103 {
1104  point(0) = point_reader.getX();
1105 
1106  if (LIBMESH_DIM >= 2)
1107  point(1) = point_reader.getY();
1108 
1109  if (LIBMESH_DIM >= 3)
1110  point(2) = point_reader.getZ();
1111 }

◆ load_rb_eim_evaluation_data()

template<typename RBEIMEvaluationReaderNumber >
void libMesh::RBDataDeserialization::load_rb_eim_evaluation_data ( RBEIMEvaluation rb_eim_eval,
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 671 of file rb_data_deserialization.C.

References libMesh::RBEIMEvaluation::add_interpolation_points_boundary_id(), libMesh::RBEIMEvaluation::add_interpolation_points_comp(), libMesh::RBEIMEvaluation::add_interpolation_points_elem_id(), libMesh::RBEIMEvaluation::add_interpolation_points_elem_type(), libMesh::RBEIMEvaluation::add_interpolation_points_JxW_all_qp(), libMesh::RBEIMEvaluation::add_interpolation_points_node_id(), libMesh::RBEIMEvaluation::add_interpolation_points_phi_i_all_qp(), libMesh::RBEIMEvaluation::add_interpolation_points_phi_i_qp(), libMesh::RBEIMEvaluation::add_interpolation_points_qp(), libMesh::RBEIMEvaluation::add_interpolation_points_side_index(), libMesh::RBEIMEvaluation::add_interpolation_points_spatial_indices(), libMesh::RBEIMEvaluation::add_interpolation_points_subdomain_id(), libMesh::RBEIMEvaluation::add_interpolation_points_xyz(), libMesh::RBEIMEvaluation::add_interpolation_points_xyz_perturbations(), libMesh::RBEIMEvaluation::get_eim_solutions_for_training_set(), libMesh::RBEIMEvaluation::get_parametrized_function(), libMesh::index_range(), libMesh::RBEIMEvaluation::initialize_param_fn_spatial_indices(), libMesh::RBParametrizedFunction::is_lookup_table, load_parameter_ranges(), load_point(), libMesh::make_range(), libMesh::RBParametrizedFunction::on_mesh_nodes(), libMesh::RBParametrizedFunction::on_mesh_sides(), libMesh::RBEIMEvaluation::resize_data_structures(), libMesh::RBEIMEvaluation::set_error_indicator_interpolation_row(), libMesh::RBEIMEvaluation::set_interpolation_matrix_entry(), and libMesh::RBEIMEvaluation::set_n_basis_functions().

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

673 {
674  // Set number of basis functions
675  unsigned int n_bfs = rb_eim_evaluation_reader.getNBfs();
676  rb_eim_evaluation.set_n_basis_functions(n_bfs);
677 
678  rb_eim_evaluation.resize_data_structures(n_bfs);
679 
680  // EIM interpolation matrix
681  {
682  auto interpolation_matrix_list = rb_eim_evaluation_reader.getInterpolationMatrix();
683 
684  libmesh_error_msg_if(interpolation_matrix_list.size() != n_bfs*(n_bfs+1)/2,
685  "Size error while reading the eim inner product matrix.");
686 
687  for (unsigned int i=0; i<n_bfs; ++i)
688  for (unsigned int j=0; j<=i; ++j)
689  {
690  unsigned int offset = i*(i+1)/2 + j;
691  rb_eim_evaluation.set_interpolation_matrix_entry(
692  i,j, load_scalar_value(interpolation_matrix_list[offset]));
693  }
694  }
695 
696  // Note that we do not check rb_eim_evaluation.use_eim_error_indicator()
697  // below because it can be false in some cases when we're reading in
698  // data (e.g. if we're using a base class RBEIMEvaluation to do the
699  // reading). The check below will still give the right behavior
700  // regardless of the value of use_eim_error_indicator().
701  if (rb_eim_evaluation_reader.getInterpolationXyz().size() > n_bfs)
702  {
703  auto error_indicator_data_list = rb_eim_evaluation_reader.getEimErrorIndicatorInterpData();
704 
705  libmesh_error_msg_if(error_indicator_data_list.size() != n_bfs,
706  "Size error while reading the eim error indicator data.");
707 
708  DenseVector<Number> eim_error_indicator_data(n_bfs);
709  for (unsigned int i=0; i<n_bfs; ++i)
710  eim_error_indicator_data(i) = load_scalar_value(error_indicator_data_list[i]);
711 
712  rb_eim_evaluation.set_error_indicator_interpolation_row(eim_error_indicator_data);
713  }
714 
715  // If we're using the EIM error indicator then we store one extra
716  // interpolation point and associated data, hence we increment n_bfs
717  // here so that we write out the extra data below.
718  //
719  // However the interpolation matrix and _error_indicator_interpolation_row
720  // does not include this extra point, which is why we read it in above
721  // before n_bfs is incremented.
722  if (rb_eim_evaluation_reader.getInterpolationXyz().size() > n_bfs)
723  n_bfs++;
724 
725  auto parameter_ranges =
726  rb_eim_evaluation_reader.getParameterRanges();
727  auto discrete_parameters_list =
728  rb_eim_evaluation_reader.getDiscreteParameters();
729 
730  load_parameter_ranges(rb_eim_evaluation,
731  parameter_ranges,
732  discrete_parameters_list);
733 
734  // Interpolation points
735  {
736  auto interpolation_points_list =
737  rb_eim_evaluation_reader.getInterpolationXyz();
738 
739  libmesh_error_msg_if(interpolation_points_list.size() != n_bfs,
740  "Size error while reading the eim interpolation points.");
741 
742  for (unsigned int i=0; i<n_bfs; ++i)
743  {
744  Point p;
745  load_point(interpolation_points_list[i], p);
746  rb_eim_evaluation.add_interpolation_points_xyz(p);
747  }
748  }
749 
750  // Interpolation points components
751  {
752  auto interpolation_points_comp_list =
753  rb_eim_evaluation_reader.getInterpolationComp();
754 
755  libmesh_error_msg_if(interpolation_points_comp_list.size() != n_bfs,
756  "Size error while reading the eim interpolation components.");
757 
758  for (unsigned int i=0; i<n_bfs; ++i)
759  {
760  rb_eim_evaluation.add_interpolation_points_comp(interpolation_points_comp_list[i]);
761  }
762  }
763 
764  // Interpolation points subdomain IDs
765  {
766  auto interpolation_points_subdomain_id_list =
767  rb_eim_evaluation_reader.getInterpolationSubdomainId();
768 
769  libmesh_error_msg_if(interpolation_points_subdomain_id_list.size() != n_bfs,
770  "Size error while reading the eim interpolation subdomain IDs.");
771 
772  for (unsigned int i=0; i<n_bfs; ++i)
773  {
774  rb_eim_evaluation.add_interpolation_points_subdomain_id(interpolation_points_subdomain_id_list[i]);
775  }
776  }
777 
778  // Interpolation points side indices, relevant if the parametrized function is defined
779  // on mesh sides or nodesets
780  if (rb_eim_evaluation.get_parametrized_function().on_mesh_sides() ||
781  rb_eim_evaluation.get_parametrized_function().on_mesh_nodes())
782  {
783  auto interpolation_points_boundary_id_list =
784  rb_eim_evaluation_reader.getInterpolationBoundaryId();
785 
786  libmesh_error_msg_if(interpolation_points_boundary_id_list.size() != n_bfs,
787  "Size error while reading the eim interpolation boundary IDs.");
788 
789  for (unsigned int i=0; i<n_bfs; ++i)
790  {
791  rb_eim_evaluation.add_interpolation_points_boundary_id(interpolation_points_boundary_id_list[i]);
792  }
793  }
794 
795  // Interpolation points element IDs
796  {
797  auto interpolation_points_elem_id_list =
798  rb_eim_evaluation_reader.getInterpolationElemId();
799 
800  libmesh_error_msg_if(interpolation_points_elem_id_list.size() != n_bfs,
801  "Size error while reading the eim interpolation element IDs.");
802 
803  for (unsigned int i=0; i<n_bfs; ++i)
804  {
805  rb_eim_evaluation.add_interpolation_points_elem_id(interpolation_points_elem_id_list[i]);
806  }
807  }
808 
809  // Interpolation points node IDs, relevant if the parametrized function is defined on mesh sides
810  if (rb_eim_evaluation.get_parametrized_function().on_mesh_nodes())
811  {
812  auto interpolation_points_node_id_list =
813  rb_eim_evaluation_reader.getInterpolationNodeId();
814 
815  libmesh_error_msg_if(interpolation_points_node_id_list.size() != n_bfs,
816  "Size error while reading the eim interpolation node IDs.");
817 
818  for (unsigned int i=0; i<n_bfs; ++i)
819  {
820  rb_eim_evaluation.add_interpolation_points_node_id(interpolation_points_node_id_list[i]);
821  }
822  }
823 
824  // Interpolation points side indices, relevant if the parametrized function is defined on mesh sides
825  if (rb_eim_evaluation.get_parametrized_function().on_mesh_sides())
826  {
827  auto interpolation_points_side_index_list =
828  rb_eim_evaluation_reader.getInterpolationSideIndex();
829 
830  libmesh_error_msg_if(interpolation_points_side_index_list.size() != n_bfs,
831  "Size error while reading the eim interpolation side indices.");
832 
833  for (unsigned int i=0; i<n_bfs; ++i)
834  {
835  rb_eim_evaluation.add_interpolation_points_side_index(interpolation_points_side_index_list[i]);
836  }
837  }
838 
839  // Interpolation points quad point indices
840  {
841  auto interpolation_points_qp_list =
842  rb_eim_evaluation_reader.getInterpolationQp();
843 
844  libmesh_error_msg_if(interpolation_points_qp_list.size() != n_bfs,
845  "Size error while reading the eim interpolation qps.");
846 
847  for (unsigned int i=0; i<n_bfs; ++i)
848  {
849  rb_eim_evaluation.add_interpolation_points_qp(interpolation_points_qp_list[i]);
850  }
851  }
852 
853  // Interpolation points JxW all qp values
854  {
855  auto interpolation_points_list_outer =
856  rb_eim_evaluation_reader.getInterpolationJxWAllQp();
857 
858  if (interpolation_points_list_outer.size() > 0)
859  {
860  libmesh_error_msg_if(interpolation_points_list_outer.size() != n_bfs,
861  "Size error while reading the eim JxW values.");
862 
863  for (unsigned int i=0; i<n_bfs; ++i)
864  {
865  auto interpolation_points_list_inner = interpolation_points_list_outer[i];
866 
867  std::vector<Real> JxW(interpolation_points_list_inner.size());
868  for (unsigned int j=0; j<JxW.size(); j++)
869  {
870  JxW[j] = interpolation_points_list_inner[j];
871  }
872  rb_eim_evaluation.add_interpolation_points_JxW_all_qp(JxW);
873  }
874  }
875  }
876 
877  // Interpolation points phi_i all qp values
878  {
879  auto interpolation_points_list_outer =
880  rb_eim_evaluation_reader.getInterpolationPhiValuesAllQp();
881 
882  if (interpolation_points_list_outer.size() > 0)
883  {
884  libmesh_error_msg_if(interpolation_points_list_outer.size() != n_bfs,
885  "Size error while reading the eim phi_i all qp values.");
886 
887  for (unsigned int i=0; i<n_bfs; ++i)
888  {
889  auto interpolation_points_list_middle = interpolation_points_list_outer[i];
890 
891  std::vector<std::vector<Real>> phi_i_all_qp(interpolation_points_list_middle.size());
892  for (auto j : index_range(phi_i_all_qp))
893  {
894  auto interpolation_points_list_inner = interpolation_points_list_middle[j];
895 
896  phi_i_all_qp[j].resize(interpolation_points_list_inner.size());
897  for (auto k : index_range(phi_i_all_qp[j]))
898  phi_i_all_qp[j][k] = interpolation_points_list_inner[k];
899  }
900  rb_eim_evaluation.add_interpolation_points_phi_i_all_qp(phi_i_all_qp);
901  }
902  }
903  }
904 
905  // Interpolation points element types
906  {
907  auto interpolation_points_elem_type_list =
908  rb_eim_evaluation_reader.getInterpolationElemType();
909 
910  if (interpolation_points_elem_type_list.size() > 0)
911  {
912  libmesh_error_msg_if(interpolation_points_elem_type_list.size() != n_bfs,
913  "Size error while reading the eim interpolation element types.");
914 
915  for (unsigned int i=0; i<n_bfs; ++i)
916  {
917  rb_eim_evaluation.add_interpolation_points_elem_type(static_cast<ElemType>(interpolation_points_elem_type_list[i]));
918  }
919  }
920  }
921 
922  // Interpolation points perturbations
923  {
924  auto interpolation_points_list_outer =
925  rb_eim_evaluation_reader.getInterpolationXyzPerturb();
926 
927  libmesh_error_msg_if(interpolation_points_list_outer.size() != n_bfs,
928  "Size error while reading the eim interpolation points.");
929 
930  for (unsigned int i=0; i<n_bfs; ++i)
931  {
932  auto interpolation_points_list_inner = interpolation_points_list_outer[i];
933 
934  std::vector<Point> perturbs(interpolation_points_list_inner.size());
935  for (unsigned int j=0; j<perturbs.size(); j++)
936  {
937  load_point(interpolation_points_list_inner[j], perturbs[j]);
938  }
939  rb_eim_evaluation.add_interpolation_points_xyz_perturbations(perturbs);
940  }
941  }
942 
943  // Optionally load EIM solutions for the training set
944  if (rb_eim_evaluation.get_parametrized_function().is_lookup_table)
945  {
946  auto eim_rhs_list_outer =
947  rb_eim_evaluation_reader.getEimSolutionsForTrainingSet();
948 
949  std::vector<DenseVector<Number>> & eim_solutions = rb_eim_evaluation.get_eim_solutions_for_training_set();
950  eim_solutions.clear();
951  eim_solutions.resize(eim_rhs_list_outer.size());
952 
953  for (auto i : make_range(eim_rhs_list_outer.size()))
954  {
955  auto eim_rhs_list_inner = eim_rhs_list_outer[i];
956 
957  DenseVector<Number> values(eim_rhs_list_inner.size());
958  for (auto j : index_range(values))
959  {
960  values(j) = load_scalar_value(eim_rhs_list_inner[j]);
961  }
962  eim_solutions[i] = values;
963  }
964  }
965 
966  // The shape function values at the interpolation points. This can be used to evaluate nodal data
967  // at EIM interpolation points, which are at quadrature points.
968  {
969  auto interpolation_points_list_outer =
970  rb_eim_evaluation_reader.getInterpolationPhiValues();
971 
972  libmesh_error_msg_if(interpolation_points_list_outer.size() != n_bfs,
973  "Size error while reading the eim interpolation points.");
974 
975  for (unsigned int i=0; i<n_bfs; ++i)
976  {
977  auto interpolation_points_list_inner = interpolation_points_list_outer[i];
978 
979  std::vector<Real> phi_i_qp(interpolation_points_list_inner.size());
980  for (unsigned int j=0; j<phi_i_qp.size(); j++)
981  {
982  phi_i_qp[j] = interpolation_points_list_inner[j];
983  }
984  rb_eim_evaluation.add_interpolation_points_phi_i_qp(phi_i_qp);
985  }
986  }
987 
988  // Read in the spatial indices at interpolation points. This data is optional
989  // so this may be empty.
990  {
991  auto interpolation_points_list_outer =
992  rb_eim_evaluation_reader.getInterpolationSpatialIndices();
993 
994  for (unsigned int i=0; i<interpolation_points_list_outer.size(); ++i)
995  {
996  auto interpolation_points_list_inner = interpolation_points_list_outer[i];
997 
998  std::vector<unsigned int> spatial_indices(interpolation_points_list_inner.size());
999  for (unsigned int j=0; j<spatial_indices.size(); j++)
1000  {
1001  spatial_indices[j] = interpolation_points_list_inner[j];
1002  }
1003  rb_eim_evaluation.add_interpolation_points_spatial_indices(spatial_indices);
1004  }
1005 
1006  rb_eim_evaluation.initialize_param_fn_spatial_indices();
1007  }
1008 }
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...
void load_point(RBData::Point3D::Reader point_reader, Point &point)
Helper function that loads point data.
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:134
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:111

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

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_transient_rb_evaluation_data(), and libMesh::RBDataDeserialization::RBEvaluationDeserialization::read_from_file().

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

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

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(), load_parameter_ranges(), libMesh::make_range(), libMesh::Real, and libMesh::RBSCMEvaluation::SCM_UB_vectors.

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

1013 {
1014  auto parameter_ranges =
1015  rb_scm_evaluation_reader.getParameterRanges();
1016  auto discrete_parameters_list =
1017  rb_scm_evaluation_reader.getDiscreteParameters();
1018  load_parameter_ranges(rb_scm_eval,
1019  parameter_ranges,
1020  discrete_parameters_list);
1021 
1022  unsigned int n_A_terms = rb_scm_eval.get_rb_theta_expansion().get_n_A_terms();
1023 
1024  {
1025  auto b_min_list = rb_scm_evaluation_reader.getBMin();
1026 
1027  libmesh_error_msg_if(b_min_list.size() != n_A_terms,
1028  "Size error while reading B_min");
1029 
1030  rb_scm_eval.B_min.clear();
1031  for (unsigned int i=0; i<n_A_terms; i++)
1032  rb_scm_eval.B_min.push_back(b_min_list[i]);
1033  }
1034 
1035  {
1036  auto b_max_list = rb_scm_evaluation_reader.getBMax();
1037 
1038  libmesh_error_msg_if(b_max_list.size() != n_A_terms,
1039  "Size error while reading B_max");
1040 
1041  rb_scm_eval.B_max.clear();
1042  for (unsigned int i=0; i<n_A_terms; i++)
1043  rb_scm_eval.B_max.push_back(b_max_list[i]);
1044  }
1045 
1046  {
1047  auto cJ_stability_vector =
1048  rb_scm_evaluation_reader.getCJStabilityVector();
1049 
1050  rb_scm_eval.C_J_stability_vector.clear();
1051  for (const auto & val : cJ_stability_vector)
1052  rb_scm_eval.C_J_stability_vector.push_back(val);
1053  }
1054 
1055  {
1056  auto cJ_parameters_outer =
1057  rb_scm_evaluation_reader.getCJ();
1058 
1059  rb_scm_eval.C_J.resize(cJ_parameters_outer.size());
1060  for (auto i : make_range(cJ_parameters_outer.size()))
1061  {
1062  auto cJ_parameters_inner =
1063  cJ_parameters_outer[i];
1064 
1065  for (auto j : make_range(cJ_parameters_inner.size()))
1066  {
1067  std::string param_name = cJ_parameters_inner[j].getName();
1068  Real param_value = cJ_parameters_inner[j].getValue();
1069  rb_scm_eval.C_J[i].set_value(param_name, param_value);
1070  }
1071  }
1072  }
1073 
1074  {
1075  auto scm_ub_vectors =
1076  rb_scm_evaluation_reader.getScmUbVectors();
1077 
1078  // The number of UB vectors is the same as the number of C_J values
1079  unsigned int n_C_J_values = rb_scm_evaluation_reader.getCJ().size();
1080 
1081  libmesh_error_msg_if(scm_ub_vectors.size() != n_C_J_values*n_A_terms,
1082  "Size mismatch in SCB UB vectors");
1083 
1084  rb_scm_eval.SCM_UB_vectors.resize( n_C_J_values );
1085  for (unsigned int i=0; i<n_C_J_values; i++)
1086  {
1087  rb_scm_eval.SCM_UB_vectors[i].resize(n_A_terms);
1088  for (unsigned int j=0; j<n_A_terms; j++)
1089  {
1090  unsigned int offset = i*n_A_terms + j;
1091  rb_scm_eval.SCM_UB_vectors[i][j] = scm_ub_vectors[offset];
1092 
1093  }
1094  }
1095  }
1096 
1097 }
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...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:134

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

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().

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