Thank you very much. It works now. 

On Wednesday, 29 July 2020 13:49:15 UTC-7, Jean-Paul Pelteret wrote:
>
> Hi,
>
> The problem is most likely that you’ve got on (static) variable 
> “times_and_names” that is referencing two completely different of output 
> files. This means that for every tilmestep, you’ve now got two entries in 
> this file and Paraview probably doesn’t know what to do with it. You can 
> easily check this by commenting out one call to “output_results()” and then 
> the other. You’ll see that solution field because now the formatting of the 
> PVD file is what Paraview expects it to be.
>
> There are at least two solutions to the problem. You either (1) have to 
> move the “ times_and_names” variable out of this function and have one copy 
> for the “P” output files and one for the “T” output files, or (2) simply 
> add both solutions to the same output files. You could achieve this in the 
> following way (assuming that the DoFHandler is common to both fields, which 
> it looks like it must be):
>
> template <int dim>
> void CoupledTH<dim>::output_results(const 
> std::vector<std::pair<const Vector<double>*, 
> std::string>>& solutions_and_names) const {
>   DataOut<dim> data_out;
>
>   data_out.attach_dof_handler(dof_handler);
>
>   for (auto solution_and_name : solutions_and_names)
>   {
>     const Vector<double>*solution = solution_and_name.first;
>     const std::string var_name = solution_and_name.second;
>     data_out.add_data_vector(*solution, var_name);
>   }
>  ...
> }
>
> template <int dim>
> void CoupledTH<dim>::run() {
>     ...
>     do {
>         …
>
>         std::vector<std::pair<const Vector<double>*,std::string>> 
> solutions_and_names;
>         solutions_and_names.push_back(std::make_pair(&T_solution, "T"));
>         solutions_and_names.push_back(std::make_pair(&P_solution, "P"));
>         output_results(solutions_and_names );
>     } while (time < period);
> }
>
> Notice that you probably want to avoid copying the vector, so you probably 
> want to consider using pointers (or a std::reference) if you go this 
> (preferred) route. Its easiest to collect them up as a std::pair, since 
> they are logically bound to one another.
>
> I hope that this helps, either as a real solution to the problem or by 
> giving you a hint as to how you might be able to navigate your way around 
> the problem that you’re seeing.
>
> Best,
> Jean-Paul
>
> On 29 Jul 2020, at 19:39, 孙翔 <shya...@gmail.com <javascript:>> wrote:
>
>
> Hi, all,
>
> I followed step-18 and wrote an HPC code in which I used a set of DOFs to 
> output two variables T and P. I debugged the coed and the code runs well. 
> However, it only outputs one variable, for the other variable, it only 
> outputs the number of processors not the value of the variable. Parts of 
> the code is attached as follows. 
>
> The code looped to get the variables T and P solutions at each time step 
> using the function output_results. The output_resutlts function was almost 
> copied from step-18. I hope the code can output both T and P values, but 
> when I use ParaView to read the files, ParaView only displays the value of 
> P and can not display the other. Could someone help me check the code? 
> Thank you very much.
>
> template <int dim>
> void CoupledTH<dim>::run() {
>     make_grid();
>     setup_system();
>     do {
>         assemble_P_system();
>         linear_solve_P();
>         assemble_T_system();
>         linear_solve_T();
>         time += time_step;
>         output_results(T_solution, "T");
>         output_results(P_solution, "P");
>     } while (time < period);
> }
>
> template <int dim>
> void CoupledTH<dim>::output_results(Vector<double>& solution,
>                                     std::string var_name) const {
>   DataOut<dim> data_out;
>
>   data_out.attach_dof_handler(dof_handler);
>
>   data_out.add_data_vector(solution, var_name);
>
>   std::vector<types::subdomain_id> partition_int(
>       triangulation.n_active_cells());
>
>   GridTools::get_subdomain_association(triangulation, partition_int);
>
>   const Vector<double> partitioning(partition_int.begin(), partition_int.
> end());
>
>   data_out.add_data_vector(partitioning, "partitioning");
>
>   data_out.build_patches();
>
>   const std::string filename =
>       var_name + "-solution-" + Utilities::int_to_string(timestep_number, 
> 3);
>
>   const std::string pvtu_master_filename = data_out.
> write_vtu_with_pvtu_record(
>       "outputfiles/", filename, timestep_number, mpi_communicator);
>
>   if (this_mpi_process == 0) {
>     static std::vector<std::pair<double, std::string>> times_and_names;
>     times_and_names.push_back(
>         std::pair<double, std::string>(time, pvtu_master_filename));
>     std::ofstream pvd_output("outputfiles/"+ var_name + "_solution.pvd");
>     DataOutBase::write_pvd_record(pvd_output, times_and_names);
>   }
> }
>
>
> -- 
> The deal.II project is located at http://www.dealii.org/
> For mailing list/forum options, see 
> https://groups.google.com/d/forum/dealii?hl=en
> --- 
> You received this message because you are subscribed to the Google Groups 
> "deal.II User Group" group.
> To unsubscribe from this group and stop receiving emails from it, send an 
> email to dea...@googlegroups.com <javascript:>.
> To view this discussion on the web visit 
> https://groups.google.com/d/msgid/dealii/b9800f64-5131-403f-9371-a92c1359e413o%40googlegroups.com
>  
> <https://groups.google.com/d/msgid/dealii/b9800f64-5131-403f-9371-a92c1359e413o%40googlegroups.com?utm_medium=email&utm_source=footer>
> .
>
>
>

-- 
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see 
https://groups.google.com/d/forum/dealii?hl=en
--- 
You received this message because you are subscribed to the Google Groups 
"deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email 
to dealii+unsubscr...@googlegroups.com.
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/7b792d23-637a-4dbc-8048-92ee40dd146co%40googlegroups.com.

Reply via email to