Dear Bangerth, Thank you for the reply. I modified the output results of the cook membrane problem to get tau(Kirchoff stress) in the output, initially, it was giving displacement as output. modified piece of code is able to run in version 9.4.0 but not able to run in version 9.4.2 which i am using. please see the modified output_results() as below and the error too. Thank you.
template <int dim,typename NumberType> void Solid<dim,NumberType>::output_results() const { DataOut<dim> data_out; data_out.attach_dof_handler(dof_handler_ref); std::vector<std::string> solution_names; solution_names.emplace_back("u_x"); solution_names.emplace_back("u_y"); solution_names.emplace_back("u_z"); data_out.add_data_vector(solution_n, solution_names); //// next part calculates Kirchhoff stress using scratch and pertaskdata FE_DGQ<dim> history_fe (1); DoFHandler<dim> history_dof_handler (triangulation); history_dof_handler.distribute_dofs (history_fe); std::vector< std::vector< Vector<double> > > history_field (dim, std::vector< Vector<double> >(dim)), local_history_values_at_qpoints (dim, std::vector< Vector<double> >(dim)), local_history_fe_values (dim, std::vector< Vector<double> >(dim)); for (unsigned int i=0; i<dim; i++) for (unsigned int j=0; j<dim; j++) { history_field[i][j].reinit(history_dof_handler.n_dofs()); local_history_values_at_qpoints[i][j].reinit(qf_cell.size()); local_history_fe_values[i][j].reinit(history_fe.dofs_per_cell); } FullMatrix<double> qpoint_to_dof_matrix (history_fe.dofs_per_cell, qf_cell.size()); FETools::compute_projection_from_quadrature_points_matrix (history_fe, qf_cell, qf_cell, qpoint_to_dof_matrix); typename DoFHandler<dim>::active_cell_iterator cell = dof_handler_ref. begin_active(), endc = dof_handler_ref.end(), dg_cell = history_dof_handler.begin_active(); for (; cell!=endc; ++cell, ++dg_cell) { /////////////////// added to define scratch and pertaskdata ////////////////////////////////////// const UpdateFlags uf_cell(update_gradients | update_JxW_values); const UpdateFlags uf_face(update_values | update_JxW_values); typename Assembler_Base<dim,NumberType>::PerTaskData_ASM data(this); typename Assembler_Base<dim,NumberType>::ScratchData_ASM scratch(fe, qf_cell, uf_cell, qf_face, uf_face, solution_n); const FEValuesExtractors::Vector &u_fe = data.solid->u_fe; data.reset(); scratch.reset(); scratch.fe_values_ref.reinit(cell); cell->get_dof_indices(data.local_dof_indices); const std::vector<std::shared_ptr<const PointHistory<dim,NumberType> > > lqph = const_cast<const Solid<dim,NumberType> *>(data.solid)-> quadrature_point_history.get_data(cell); Assert(lqph.size() == n_q_points, ExcInternalError()); scratch.fe_values_ref[u_fe].get_function_gradients(scratch.solution_total, scratch.solution_grads_u_total); /////////////////////////////////////////////////////////////////// for (unsigned int i=0; i<dim; i++) for (unsigned int j=0; j<dim; j++) { for (unsigned int q=0; q<qf_cell.size(); ++q) { //////////////// added to compute Kirchhoff stress tau ////////////////////////////////////// const Tensor<2,dim,NumberType> &grad_u = scratch.solution_grads_u_total[q]; const Tensor<2,dim,NumberType> F = Physics::Elasticity::Kinematics::F (grad_u); const NumberType det_F = determinant(F); const Tensor<2,dim,NumberType> F_bar = Physics::Elasticity::Kinematics:: F_iso(F); const SymmetricTensor<2,dim,NumberType> b_bar = Physics::Elasticity::Kinematics::b(F_bar); const SymmetricTensor<2,dim,NumberType> tau = lqph[q]->get_tau(det_F,b_bar); ////////////////////////////////////////////////////////////// local_history_values_at_qpoints[i][j](q) = tau[i][j]; // tau_locally stored } qpoint_to_dof_matrix.vmult (local_history_fe_values[i][j], local_history_values_at_qpoints[i][j]); // q_point to node extrapolation dg_cell->set_dof_values (local_history_fe_values[i][j], history_field[i][j]); // local to global dof transfer } } data_out.add_data_vector(history_dof_handler,history_field[0][0], "Kirchhoff_st_xx"); data_out.add_data_vector(history_dof_handler,history_field[1][1], "Kirchhoff_st_yy"); data_out.add_data_vector(history_dof_handler,history_field[2][2], "Kirchhoff_st_zz"); data_out.add_data_vector(history_dof_handler,history_field[0][1], "Kirchhoff_st_xy"); data_out.add_data_vector(history_dof_handler,history_field[0][2], "Kirchhoff_st_xz"); data_out.add_data_vector(history_dof_handler,history_field[1][2], "Kirchhoff_st_yz"); // Since we are dealing with a large deformation problem, it would be nice // to display the result on a displaced grid! The MappingQEulerian class // linked with the DataOut class provides an interface through which this // can be achieved without physically moving the grid points in the // Triangulation object ourselves. We first need to copy the solution to // a temporary vector and then create the Eulerian mapping. We also // specify the polynomial degree to the DataOut object in order to produce // a more refined output data set when higher order polynomials are used. Vector<double> soln(solution_n.size()); for (unsigned int i = 0; i < soln.size(); ++i) soln(i) = solution_n(i); MappingQEulerian<dim> q_mapping(degree, dof_handler_ref, soln); data_out.build_patches(q_mapping, degree); std::ostringstream filename; filename << "solution-" << time.get_timestep() << ".vtk"; std::ofstream output(filename.str().c_str()); data_out.write_vtk(output); } .............................................................................................. ................................................................................................... .................................................................................................... *error i am getting* make all Consolidate compiler generated dependencies of target cook_membrane [ 50%] Building CXX object CMakeFiles/cook_membrane.dir/cook_membrane.cc.o /home/newton/dealii-candi/deal.II-v9.4.2/examples/cook_membrane/cook_membrane.cc: In instantiation of ‘void Cook_Membrane::Solid<dim, NumberType>::output_results() const [with int dim = 3; NumberType = Sacado::Fad::Exp::GeneralFad<Sacado::Fad::Exp::DynamicStorage<double, double> >]’: /home/newton/dealii-candi/deal.II-v9.4.2/examples/cook_membrane/cook_membrane.cc:1006:5: required from ‘void Cook_Membrane::Solid<dim, NumberType>::run() [with int dim = 3; NumberType = Sacado::Fad::Exp::GeneralFad<Sacado::Fad::Exp::DynamicStorage<double, double> >]’ /home/newton/dealii-candi/deal.II-v9.4.2/examples/cook_membrane/cook_membrane.cc:2421:23: required from here /home/newton/dealii-candi/deal.II-v9.4.2/examples/cook_membrane/cook_membrane.cc:2288:62: error: cannot convert ‘std::vector<dealii::Tensor<2, 3, Sacado::Fad::Exp::GeneralFad<Sacado::Fad::Exp::DynamicStorage<double, double> > >, std::allocator<dealii::Tensor<2, 3, Sacado::Fad::Exp::GeneralFad<Sacado::Fad::Exp::DynamicStorage<double, double> > > > >’ to ‘std::vector<dealii::Tensor<2, 3, double>, std::allocator<dealii::Tensor<2, 3, double> > >&’ 2288 | scratch.solution_grads_u_total); | ~~~~~~~~^~~~~~~~~~~~~~~~~~~~~~ | | | std::vector<dealii::Tensor<2, 3, Sacado::Fad::Exp::GeneralFad<Sacado::Fad::Exp::DynamicStorage<double, double> > >, std::allocator<dealii::Tensor<2, 3, Sacado::Fad::Exp::GeneralFad<Sacado::Fad::Exp::DynamicStorage<double, double> > > > > In file included from /home/newton/dealii-candi/deal.II-v9.4.2/include/deal. II/grid/grid_tools.h:33, from /home/newton/dealii-candi/deal.II -v9.4.2/examples/cook_membrane/cook_membrane.cc:40: /home/newton/dealii-candi/deal.II-v9.4.2/include/deal.II/fe/fe_values.h:1180:10: note: initializing argument 2 of ‘void dealii::FEValuesViews::Vector<dim, spacedim>::get_function_gradients(const InputVector&, std::vector<typename dealii::ProductType<typename InputVector::value_type, dealii::Tensor<2, spacedim> >::type>&) const [with InputVector = dealii::BlockVector<double>; int dim = 3; int spacedim = 3; typename dealii::ProductType<typename InputVector::value_type, dealii::Tensor<2, spacedim> >::type = dealii::Tensor<2, 3, double>; typename InputVector::value_type = double]’ 1179 | std::vector<solution_gradient_type<typename InputVector::value_type>> | ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ 1180 | &gradients) const; | ~^~~~~~~~~ /home/newton/dealii-candi/deal.II-v9.4.2/examples/cook_membrane/cook_membrane.cc:2307:44: error: cannot convert ‘dealii::internal::SymmetricTensorAccessors::Accessor<2, 3, true, 1, Sacado::Fad::Exp::GeneralFad<Sacado::Fad::Exp::DynamicStorage<double, double> > >::reference’ {aka ‘Sacado::Fad::Exp::GeneralFad<Sacado::Fad::Exp::DynamicStorage<double, double> >’} to ‘double’ in assignment 2307 | local_history_values_at_qpoints[i][j](q) = tau[i][j]; // tau_locally stored | ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~^~~~~~~~~ /home/newton/dealii-candi/deal.II-v9.4.2/examples/cook_membrane/cook_membrane.cc: In instantiation of ‘void Cook_Membrane::Solid<dim, NumberType>::output_results() const [with int dim = 3; NumberType = Sacado::Rad::ADvar<Sacado::Fad::Exp::GeneralFad<Sacado::Fad::Exp::DynamicStorage<double, double> > >]’: /home/newton/dealii-candi/deal.II-v9.4.2/examples/cook_membrane/cook_membrane.cc:1006:5: required from ‘void Cook_Membrane::Solid<dim, NumberType>::run() [with int dim = 3; NumberType = Sacado::Rad::ADvar<Sacado::Fad::Exp::GeneralFad<Sacado::Fad::Exp::DynamicStorage<double, double> > >]’ /home/newton/dealii-candi/deal.II-v9.4.2/examples/cook_membrane/cook_membrane.cc:2434:23: required from here /home/newton/dealii-candi/deal.II-v9.4.2/examples/cook_membrane/cook_membrane.cc:2288:62: error: cannot convert ‘std::vector<dealii::Tensor<2, 3, Sacado::Rad::ADvar<Sacado::Fad::Exp::GeneralFad<Sacado::Fad::Exp::DynamicStorage<double, double> > > >, std::allocator<dealii::Tensor<2, 3, Sacado::Rad::ADvar<Sacado::Fad::Exp::GeneralFad<Sacado::Fad::Exp::DynamicStorage<double, double> > > > > >’ to ‘std::vector<dealii::Tensor<2, 3, double>, std::allocator<dealii::Tensor<2, 3, double> > >&’ 2288 | scratch.solution_grads_u_total); | ~~~~~~~~^~~~~~~~~~~~~~~~~~~~~~ | | | std::vector<dealii::Tensor<2, 3, Sacado::Rad::ADvar<Sacado::Fad::Exp::GeneralFad<Sacado::Fad::Exp::DynamicStorage<double, double> > > >, std::allocator<dealii::Tensor<2, 3, Sacado::Rad::ADvar<Sacado::Fad::Exp::GeneralFad<Sacado::Fad::Exp::DynamicStorage<double, double> > > > > > In file included from /home/newton/dealii-candi/deal.II-v9.4.2/include/deal. II/grid/grid_tools.h:33, from /home/newton/dealii-candi/deal.II -v9.4.2/examples/cook_membrane/cook_membrane.cc:40: /home/newton/dealii-candi/deal.II-v9.4.2/include/deal.II/fe/fe_values.h:1180:10: note: initializing argument 2 of ‘void dealii::FEValuesViews::Vector<dim, spacedim>::get_function_gradients(const InputVector&, std::vector<typename dealii::ProductType<typename InputVector::value_type, dealii::Tensor<2, spacedim> >::type>&) const [with InputVector = dealii::BlockVector<double>; int dim = 3; int spacedim = 3; typename dealii::ProductType<typename InputVector::value_type, dealii::Tensor<2, spacedim> >::type = dealii::Tensor<2, 3, double>; typename InputVector::value_type = double]’ 1179 | std::vector<solution_gradient_type<typename InputVector::value_type>> | ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ 1180 | &gradients) const; | ~^~~~~~~~~ /home/newton/dealii-candi/deal.II-v9.4.2/examples/cook_membrane/cook_membrane.cc:2307:44: error: cannot convert ‘dealii::internal::SymmetricTensorAccessors::Accessor<2, 3, true, 1, Sacado::Rad::ADvar<Sacado::Fad::Exp::GeneralFad<Sacado::Fad::Exp::DynamicStorage<double, double> > > >::reference’ {aka ‘Sacado::Rad::ADvar<Sacado::Fad::Exp::GeneralFad<Sacado::Fad::Exp::DynamicStorage<double, double> > >’} to ‘double’ in assignment 2307 | local_history_values_at_qpoints[i][j](q) = tau[i][j]; // tau_locally stored | ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~^~~~~~~~~ make[2]: *** [CMakeFiles/cook_membrane.dir/build.make:76: CMakeFiles/cook_membrane.dir/cook_membrane.cc.o] Error 1 make[1]: *** [CMakeFiles/Makefile2:90: CMakeFiles/cook_membrane.dir/all] Error 2 make: *** [Makefile:91: all] Error 2 "make all" terminated with exit code 2. Build might be incomplete. -- On Tue, Jun 20, 2023 at 9:43 AM Wolfgang Bangerth <bange...@colostate.edu> wrote: > On 6/18/23 11:43, ME20D503 NEWTON wrote: > > > /home/newton/dealii-candi/deal.II-v9.4.2/examples/cook_membrane/cook_membrane.cc:2288:62: > error: cannot convert ‘std::vector<dealii::Tensor<2, 3, > Sacado::Fad::Exp::GeneralFad<Sacado::Fad::Exp::DynamicStorage<double, > double> > >, std::allocator<dealii::Tensor<2, 3, > Sacado::Fad::Exp::GeneralFad<Sacado::Fad::Exp::DynamicStorage<double, > double> > > > >’ to ‘std::vector<dealii::Tensor<2, 3, double>, > std::allocator<dealii::Tensor<2, 3, double> > >&’ > > > > 2288 | scratch.solution_grads_u_total); > > > > | ~~~~~~~~^~~~~~~~~~~~~~~~~~~~~~ > > > I don't know what you have changed about this program, and since you don't > attach it, it's hard to tell *why* this error is happening. > > But the compiler is quite clear about *what* is happening: You are trying > to > call a function that expects an argument of type > std::vector<Tensor<2, 3, double>> > with an argument of type > std::vector<Tensor<2, 3, > Sacado::Fad::Exp::GeneralFad<Sacado::Fad::Exp::DynamicStorage<double, > double>>> > > The compiler is simply telling you that that won't work. Since I don't > know > the code, I can't tell you how to fix this beyond the description above. > > Best > W. > > -- > ------------------------------------------------------------------------ > Wolfgang Bangerth email: bange...@colostate.edu > www: http://www.math.colostate.edu/~bangerth/ > > > -- > 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 a topic in the > Google Groups "deal.II User Group" group. > To unsubscribe from this topic, visit > https://groups.google.com/d/topic/dealii/0wxyQezh-J4/unsubscribe. > To unsubscribe from this group and all its topics, send an email to > dealii+unsubscr...@googlegroups.com. > To view this discussion on the web visit > https://groups.google.com/d/msgid/dealii/c9275fa6-d6d4-d885-4887-17fa10e9177a%40colostate.edu > . > -- 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/CAM2td%3DWqAxm%3Dv1T_NxM51ir_rhqQx9Ng_zR1-wMAUQE8wEGENw%40mail.gmail.com.