21 #include "libmesh/transient_rb_evaluation.h"
22 #include "libmesh/transient_rb_theta_expansion.h"
25 #include "libmesh/getpot.h"
26 #include "libmesh/int_range.h"
27 #include "libmesh/libmesh_logging.h"
28 #include "libmesh/numeric_vector.h"
29 #include "libmesh/xdr_cxx.h"
41 _rb_solve_data_cached(false)
69 bool resize_error_bound_data)
71 LOG_SCOPE(
"resize_data_structures()",
"TransientRBEvaluation");
82 const unsigned int Q_m = trans_theta_expansion.
get_n_M_terms();
83 const unsigned int Q_a = trans_theta_expansion.
get_n_A_terms();
84 const unsigned int Q_f = trans_theta_expansion.
get_n_F_terms();
88 for (
unsigned int q=0; q<Q_m; q++)
105 if (resize_error_bound_data)
109 for (
unsigned int i=0; i<Q_f; i++)
112 for (
unsigned int j=0; j<Q_m; j++)
118 unsigned int Q_m_hat = Q_m*(Q_m+1)/2;
120 for (
unsigned int i=0; i<Q_m_hat; i++)
123 for (
unsigned int j=0; j<Nmax; j++)
130 for (
unsigned int i=0; i<Q_a; i++)
133 for (
unsigned int j=0; j<Q_m; j++)
136 for (
unsigned int k=0; k<Nmax; k++)
147 for (
unsigned int q_m=0; q_m<Q_m; q_m++)
156 LOG_SCOPE(
"rb_solve()",
"TransientRBEvaluation");
159 libmesh_error_msg(
"ERROR: N cannot be larger than the number of basis functions in rb_solve");
165 const unsigned int Q_m = trans_theta_expansion.
get_n_M_terms();
166 const unsigned int Q_a = trans_theta_expansion.
get_n_A_terms();
167 const unsigned int Q_f = trans_theta_expansion.
get_n_F_terms();
177 for (
unsigned int n=0; n<trans_theta_expansion.
get_n_outputs(); n++)
185 RB_mass_matrix_N.
zero();
187 for (
unsigned int q_m=0; q_m<Q_m; q_m++)
190 RB_mass_matrix_N.
add(trans_theta_expansion.
eval_M_theta(q_m, mu), RB_M_q_m);
203 for (
unsigned int q_a=0; q_a<Q_a; q_a++)
215 for (
unsigned int q_f=0; q_f<Q_f; q_f++)
242 for (
unsigned int time_level=0; time_level<=n_time_steps; time_level++)
252 for (
unsigned int n=0; n<trans_theta_expansion.
get_n_outputs(); n++)
264 Real error_bound_sum = 0.;
278 for (
unsigned int n=0; n<trans_theta_expansion.
get_n_outputs(); n++)
296 for (
unsigned int time_level=1; time_level<=n_time_steps; time_level++)
317 for (
unsigned int n=0; n<trans_theta_expansion.
get_n_outputs(); n++)
341 for (
unsigned int n=0; n<trans_theta_expansion.
get_n_outputs(); n++)
385 for (
unsigned int time_level=1; time_level<=n_time_steps; time_level++)
422 LOG_SCOPE(
"cache_online_residual_terms()",
"TransientRBEvaluation");
428 const unsigned int Q_m = trans_theta_expansion.
get_n_M_terms();
429 const unsigned int Q_a = trans_theta_expansion.
get_n_A_terms();
430 const unsigned int Q_f = trans_theta_expansion.
get_n_F_terms();
434 for (
unsigned int q_f1=0; q_f1<Q_f; q_f1++)
437 for (
unsigned int q_f2=q_f1; q_f2<Q_f; q_f2++)
439 Real delta = (q_f1==q_f2) ? 1. : 2.;
448 for (
unsigned int q_f=0; q_f<Q_f; q_f++)
451 for (
unsigned int q_a=0; q_a<Q_a; q_a++)
454 for (
unsigned int i=0; i<N; i++)
464 for (
unsigned int q_a1=0; q_a1<Q_a; q_a1++)
467 for (
unsigned int q_a2=q_a1; q_a2<Q_a; q_a2++)
470 Real delta = (q_a1==q_a2) ? 1. : 2.;
472 for (
unsigned int i=0; i<N; i++)
474 for (
unsigned int j=0; j<N; j++)
477 cached_theta_q_a1*cached_theta_q_a2*
486 for (
unsigned int q_f=0; q_f<Q_f; q_f++)
489 for (
unsigned int q_m=0; q_m<Q_m; q_m++)
492 for (
unsigned int i=0; i<N; i++)
500 for (
unsigned int q_a=0; q_a<Q_a; q_a++)
504 for (
unsigned int q_m=0; q_m<Q_m; q_m++)
508 for (
unsigned int i=0; i<N; i++)
510 for (
unsigned int j=0; j<N; j++)
520 for (
unsigned int q_m1=0; q_m1<Q_m; q_m1++)
523 for (
unsigned int q_m2=q_m1; q_m2<Q_m; q_m2++)
526 Real delta = (q_m1==q_m2) ? 1. : 2.;
528 for (
unsigned int i=0; i<N; i++)
530 for (
unsigned int j=0; j<N; j++)
533 cached_theta_q_m1*cached_theta_q_m2*
544 LOG_SCOPE(
"compute_residual_dual_norm()",
"TransientRBEvaluation");
555 for (
unsigned int i=0; i<N; i++)
567 for (
unsigned int i=0; i<N; i++)
568 for (
unsigned int j=0; j<N; j++)
578 libMesh::out <<
"Warning: Square of residual norm is negative "
579 <<
"in TransientRBEvaluation::compute_residual_dual_norm()" << std::endl;
584 residual_norm_sq =
std::abs(residual_norm_sq);
592 LOG_SCOPE(
"uncached_compute_residual_dual_norm()",
"TransientRBEvaluation");
601 const unsigned int Q_m = trans_theta_expansion.
get_n_M_terms();
602 const unsigned int Q_a = trans_theta_expansion.
get_n_A_terms();
603 const unsigned int Q_f = trans_theta_expansion.
get_n_F_terms();
608 std::vector<Number> RB_u_euler_theta(N);
609 std::vector<Number> mass_coeffs(N);
610 for (
unsigned int i=0; i<N; i++)
617 Number residual_norm_sq = 0.;
620 for (
unsigned int q_f1=0; q_f1<Q_f; q_f1++)
623 for (
unsigned int q_f2=q_f1; q_f2<Q_f; q_f2++)
625 Real delta = (q_f1==q_f2) ? 1. : 2.;
632 for (
unsigned int q_f=0; q_f<Q_f; q_f++)
635 for (
unsigned int q_a=0; q_a<Q_a; q_a++)
638 for (
unsigned int i=0; i<N; i++)
640 residual_norm_sq += 2.*RB_u_euler_theta[i]*cached_theta_q_f*cached_theta_q_a*
647 for (
unsigned int q_a1=0; q_a1<Q_a; q_a1++)
650 for (
unsigned int q_a2=q_a1; q_a2<Q_a; q_a2++)
653 Real delta = (q_a1==q_a2) ? 1. : 2.;
655 for (
unsigned int i=0; i<N; i++)
657 for (
unsigned int j=0; j<N; j++)
659 residual_norm_sq += delta*RB_u_euler_theta[i]*RB_u_euler_theta[j]*
660 cached_theta_q_a1*cached_theta_q_a2*
670 for (
unsigned int q_m1=0; q_m1<Q_m; q_m1++)
673 for (
unsigned int q_m2=q_m1; q_m2<Q_m; q_m2++)
676 Real delta = (q_m1==q_m2) ? 1. : 2.;
678 for (
unsigned int i=0; i<N; i++)
680 for (
unsigned int j=0; j<N; j++)
682 residual_norm_sq += delta*mass_coeffs[i]*mass_coeffs[j]*
683 cached_theta_q_m1*cached_theta_q_m2*
691 for (
unsigned int q_f=0; q_f<Q_f; q_f++)
694 for (
unsigned int q_m=0; q_m<Q_m; q_m++)
697 for (
unsigned int i=0; i<N; i++)
704 for (
unsigned int q_a=0; q_a<Q_a; q_a++)
708 for (
unsigned int q_m=0; q_m<Q_m; q_m++)
712 for (
unsigned int i=0; i<N; i++)
714 for (
unsigned int j=0; j<N; j++)
716 residual_norm_sq += 2.*RB_u_euler_theta[i]*mass_coeffs[j]*
717 cached_theta_q_a*cached_theta_q_m*
726 libMesh::out <<
"Warning: Square of residual norm is negative "
727 <<
"in TransientRBEvaluation::compute_residual_dual_norm()" << std::endl;
732 residual_norm_sq =
std::abs(residual_norm_sq);
742 const bool write_binary_data)
744 LOG_SCOPE(
"legacy_write_offline_data_to_files()",
"TransientRBEvaluation");
750 const unsigned int Q_m = trans_theta_expansion.
get_n_M_terms();
751 const unsigned int Q_a = trans_theta_expansion.
get_n_A_terms();
752 const unsigned int Q_f = trans_theta_expansion.
get_n_F_terms();
760 const std::string suffix = write_binary_data ?
".xdr" :
".dat";
764 std::ostringstream file_name;
768 file_name << directory_name <<
"/temporal_discretization_data" << suffix;
769 Xdr temporal_discretization_data_out(file_name.str(), mode);
771 Real real_value;
unsigned int int_value;
772 real_value =
get_delta_t(); temporal_discretization_data_out << real_value;
773 real_value =
get_euler_theta(); temporal_discretization_data_out << real_value;
774 int_value =
get_n_time_steps(); temporal_discretization_data_out << int_value;
775 int_value =
get_time_step(); temporal_discretization_data_out << int_value;
776 temporal_discretization_data_out.close();
781 file_name << directory_name <<
"/RB_L2_matrix" << suffix;
782 Xdr RB_L2_matrix_out(file_name.str(), mode);
784 for (
unsigned int i=0; i<n_bfs; i++)
786 for (
unsigned int j=0; j<n_bfs; j++)
791 RB_L2_matrix_out.close();
794 for (
unsigned int q_m=0; q_m<Q_m; q_m++)
797 file_name << directory_name <<
"/RB_M_";
798 file_name << std::setw(3)
799 << std::setprecision(0)
804 Xdr RB_M_q_m_out(file_name.str(), mode);
806 for (
unsigned int i=0; i<n_bfs; i++)
808 for (
unsigned int j=0; j<n_bfs; j++)
813 RB_M_q_m_out.close();
819 file_name << directory_name <<
"/initial_conditions" << suffix;
820 Xdr initial_conditions_out(file_name.str(), mode);
822 file_name << directory_name <<
"/initial_L2_error" << suffix;
823 Xdr initial_L2_error_out(file_name.str(), mode);
825 for (
unsigned int i=0; i<n_bfs; i++)
828 for (
unsigned int j=0; j<=i; j++)
833 initial_conditions_out.close();
834 initial_L2_error_out.close();
838 file_name << directory_name <<
"/Fq_Mq_terms" << suffix;
839 Xdr RB_Fq_Mq_terms_out(file_name.str(), mode);
841 for (
unsigned int q_f=0; q_f<Q_f; q_f++)
843 for (
unsigned int q_m=0; q_m<Q_m; q_m++)
845 for (
unsigned int i=0; i<n_bfs; i++)
851 RB_Fq_Mq_terms_out.close();
855 file_name << directory_name <<
"/Mq_Mq_terms" << suffix;
856 Xdr RB_Mq_Mq_terms_out(file_name.str(), mode);
858 unsigned int Q_m_hat = Q_m*(Q_m+1)/2;
859 for (
unsigned int q=0; q<Q_m_hat; q++)
861 for (
unsigned int i=0; i<n_bfs; i++)
863 for (
unsigned int j=0; j<n_bfs; j++)
869 RB_Mq_Mq_terms_out.close();
873 file_name << directory_name <<
"/Aq_Mq_terms" << suffix;
874 Xdr RB_Aq_Mq_terms_out(file_name.str(), mode);
876 for (
unsigned int q_a=0; q_a<Q_a; q_a++)
878 for (
unsigned int q_m=0; q_m<Q_m; q_m++)
880 for (
unsigned int i=0; i<n_bfs; i++)
882 for (
unsigned int j=0; j<n_bfs; j++)
889 RB_Aq_Mq_terms_out.close();
894 bool read_error_bound_data,
895 const bool read_binary_data)
897 LOG_SCOPE(
"legacy_read_offline_data_from_files()",
"TransientRBEvaluation");
903 const unsigned int Q_m = trans_theta_expansion.
get_n_M_terms();
904 const unsigned int Q_a = trans_theta_expansion.
get_n_A_terms();
905 const unsigned int Q_f = trans_theta_expansion.
get_n_F_terms();
915 const std::string suffix = read_binary_data ?
".xdr" :
".dat";
918 std::ostringstream file_name;
922 file_name << directory_name <<
"/temporal_discretization_data" << suffix;
925 Xdr temporal_discretization_data_in(file_name.str(), mode);
927 Real real_value;
unsigned int int_value;
928 temporal_discretization_data_in >> real_value;
set_delta_t(real_value);
929 temporal_discretization_data_in >> real_value;
set_euler_theta(real_value);
931 temporal_discretization_data_in >> int_value;
set_time_step(int_value);
932 temporal_discretization_data_in.close();
935 file_name << directory_name <<
"/RB_L2_matrix" << suffix;
938 Xdr RB_L2_matrix_in(file_name.str(), mode);
940 for (
unsigned int i=0; i<n_bfs; i++)
942 for (
unsigned int j=0; j<n_bfs; j++)
945 RB_L2_matrix_in >>
value;
949 RB_L2_matrix_in.close();
952 for (
unsigned int q_m=0; q_m<Q_m; q_m++)
955 file_name << directory_name <<
"/RB_M_";
956 file_name << std::setw(3)
957 << std::setprecision(0)
965 Xdr RB_M_q_m_in(file_name.str(), mode);
967 for (
unsigned int i=0; i<n_bfs; i++)
969 for (
unsigned int j=0; j<n_bfs; j++)
972 RB_M_q_m_in >>
value;
983 file_name << directory_name <<
"/initial_conditions" << suffix;
986 Xdr initial_conditions_in(file_name.str(), mode);
989 file_name << directory_name <<
"/initial_L2_error" << suffix;
992 Xdr initial_L2_error_in(file_name.str(), mode);
994 for (
unsigned int i=0; i<n_bfs; i++)
997 for (
unsigned int j=0; j<=i; j++)
1002 initial_conditions_in.close();
1003 initial_L2_error_in.close();
1006 if (read_error_bound_data)
1010 file_name << directory_name <<
"/Fq_Mq_terms" << suffix;
1013 Xdr RB_Fq_Mq_terms_in(file_name.str(), mode);
1015 for (
unsigned int q_f=0; q_f<Q_f; q_f++)
1017 for (
unsigned int q_m=0; q_m<Q_m; q_m++)
1019 for (
unsigned int i=0; i<n_bfs; i++)
1025 RB_Fq_Mq_terms_in.close();
1029 file_name << directory_name <<
"/Mq_Mq_terms" << suffix;
1032 Xdr RB_Mq_Mq_terms_in(file_name.str(), mode);
1034 unsigned int Q_m_hat = Q_m*(Q_m+1)/2;
1035 for (
unsigned int q=0; q<Q_m_hat; q++)
1037 for (
unsigned int i=0; i<n_bfs; i++)
1039 for (
unsigned int j=0; j<n_bfs; j++)
1045 RB_Mq_Mq_terms_in.close();
1049 file_name << directory_name <<
"/Aq_Mq_terms" << suffix;
1052 Xdr RB_Aq_Mq_terms_in(file_name.str(), mode);
1054 for (
unsigned int q_a=0; q_a<Q_a; q_a++)
1056 for (
unsigned int q_m=0; q_m<Q_m; q_m++)
1058 for (
unsigned int i=0; i<n_bfs; i++)
1060 for (
unsigned int j=0; j<n_bfs; j++)
1067 RB_Aq_Mq_terms_in.close();