Jack:
You are creating a Function<dim*dim+dim+1>=Function<7>, but that is not right. The template argument denotes the number of components of the *input* vector 'x' for which you want to evaluate 'f(x)'. It does not denote the number of *output* vector components. It should only be ComponentSelectFunction<dim>.

Best
 W.

On 3/7/23 20:38, jack urombo wrote:
*** Caution: EXTERNAL Sender ***

Took a brief hiatus. and managed to pass a huddle following Marc's suggestion.

Now I have set to construct convergence tables for the polymer flow problem I am solving. It involves solving a Stokes Problem for velocity and pressure (u,p), (size dim+1)  and the Constitutive Equations for the stress tensor tau, (size dim*dim). The solution has 7 components for dim=2. To calculate the errors, I need the exact solution [(u,p), tau], ( size dim +1 + dim*dim). To pick the correct component I use the ComponentSelectFunction:

constComponentSelectFunction<dim*dim+dim+1>

velocity_mask(std::make_pair(0,dim),dim*dim+dim+1);

constComponentSelectFunction<dim*dim+dim+1>

pressure_mask(dim,dim+3);

constComponentSelectFunction<dim*dim+dim+1>

stress_mask(std::pair(dim+3,dim*dim+dim+1),dim*dim+dim+1);

Now the errors that I am getting are:

/home/jurombo/binaires/dealii/nnf/nnewton/oldroyd_b_stab/oldroyd_bn.cc:951:43:error:
 no matching fu
nction for call to ‘integrate_difference(dealii::MappingQ<2, 2>&, dealii::DoFHandler<2, 2>&, dealii: :BlockVector<double>&, Viscoelastic::ExactSolution<7>&, dealii::Vector<double>&, dealii::QGauss<7>&, dealii::VectorTools::NormType, const dealii::ComponentSelectFunction<7, double>*)’
  951 | VectorTools::integrate_difference (mapping, dof_handler,
      | ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~^~~~~~~~~~~~~~~~~~~~~~
  952 |                                            solution,
      | ~~~~~~~~~
  953 |                                            exact_solution,
      | ~~~~~~~~~~~~~~~
  954 |                                            cellwise_errors,
      | ~~~~~~~~~~~~~~~~
  955 |                                            quadrature,
      | ~~~~~~~~~~~
  956 |                                            VectorTools::H1_seminorm,
      | ~~~~~~~~~~~~~~~~~~~~~~~~~
  957 |                                            &velocity_mask);
      | ~~~~~~~~~~~~~~~
/usr/local/share/dealii/include/deal.II/numerics/vector_tools_integrate_difference.h:143:3:note:
 ca
ndidate: ‘template<int dim, class InVector, class OutVector, int spacedim> void dealii::VectorTools: :integrate_difference(const dealii::Mapping<dim, spacedim>&, const dealii::DoFHandler<dim, spacedim> &, const InVector&, const dealii::Function<spacedim, typename InVector::value_type>&, OutVector&, co nst dealii::Quadrature<dim>&, const NormType&, const dealii::Function<spacedim, double>*, double)’
  143 | integrate_difference(
      | ^~~~~~~~~~~~~~~~~~~~
/usr/local/share/dealii/include/deal.II/numerics/vector_tools_integrate_difference.h:143:3:note:
template argument deduction/substitution failed:
/home/jurombo/binaires/dealii/nnf/nnewton/oldroyd_b_stab/oldroyd_bn.cc:951:43:note:
   deduced confl
icting values for non-type parameter ‘spacedim’ (‘2’ and ‘7’)
  951 | VectorTools::integrate_difference (mapping, dof_handler,
      | ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~^~~~~~~~~~~~~~~~~~~~~~
  952 |                                            solution,
      | ~~~~~~~~~
  953 |                                            exact_solution,
      | ~~~~~~~~~~~~~~~
  954 |                                            cellwise_errors,
      | ~~~~~~~~~~~~~~~~
  955 |                                            quadrature,
      | ~~~~~~~~~~~
  956 |                                            VectorTools::H1_seminorm,
      | ~~~~~~~~~~~~~~~~~~~~~~~~~
  957 |                                            &velocity_mask);
      | ~~~~~~~~~~~~~~~
/usr/local/share/dealii/include/deal.II/numerics/vector_tools_integrate_difference.h:160:3:note:
 ca
ndidate: ‘template<int dim, class InVector, class OutVector, int spacedim> void dealii::VectorTools: :integrate_difference(const dealii::DoFHandler<dim, spacedim>&, const InVector&, const dealii::Funct ion<spacedim, typename InVector::value_type>&, OutVector&, const dealii::Quadrature<dim>&, const Nor
mType&, const dealii::Function<spacedim, double>*, double)’

The issues seem to be emanating from the form of the *exact_solution *and how the components are selected.
Ant suggestions on resolving this?

----------------------------------------------------------------------------
On Monday, 6 February 2023 at 03:11:19 UTC+2 wrote:

    Hello Jack,

    the error message is quite descriptive and tells you that something is off
    with the number of components of your Function object:
    //The violated condition was:
        exact_solution.n_components == n_components//

    Check line 200 of your code:
             ExactSolution () : Function<dim>(dim+1) {}

    You call the constructor for Function with 3 components (since dim=2 in
    your case). However, your implementation returns 7 components. Make sure
    that you call the constructor with the right number of components.

    Marc

    On Saturday, February 4, 2023 at 8:54:47 AM UTC-7 jack urombo wrote:

        i am using the implimenation by Jaekwang Kim to model polymer flows.
        The major change that I want to do is to output convergence tables
        like in Step-7 and Step-22, which means I have to have an ExactSolution.

        The code is encountering an error because the dimension of the
        solutions do not match at a some point when generating output:
        /
        /
        //
        /An error occurred in line <820> of file
        </home/jurombo/Downloads/dealii-9.3.2/include/deal.II/n
        umerics/vector_tools_integrate_difference.templates.h> in function
            void dealii::VectorTools::internal::do_integrate_difference(const
        dealii::hp::MappingCollec
        tion<dim, spacedim>&, const dealii::DoFHandler<dim, spacedim>&, const
        InVector&, const dealii::
        Function<spacedim, typename InVector::value_type>&, OutVector&, const
        dealii::hp::QCollection<d
        im>&, const dealii::VectorTools::NormType&, const
        dealii::Function<spacedim>*, double) [with in
        t dim = 2; int spacedim = 2; InVector = dealii::BlockVector<double>;
        OutVector = dealii::Vector
        <double>; typename InVector::value_type = double]
        The violated condition was:
            exact_solution.n_components == n_components
        Additional information:
            Dimension 3 not equal to 7.

        Stacktrace:
        -----------
        #0  /usr/local/share/dealii/lib/libdeal_II.g.so.9.3.2:
        #1  /usr/local/share/dealii/lib/libdeal_II.g.so.9.3.2: void
        dealii::VectorTools::integrate_diff
        erence<2, dealii::BlockVector<double>, dealii::Vector<double>,
        2>(dealii::Mapping<2, 2> const&,
        dealii::DoFHandler<2, 2> const&, dealii::BlockVector<double> const&,
        dealii::Function<2, deali
        i::BlockVector<double>::value_type> const&, dealii::Vector<double>&,
        dealii::Quadrature<2> cons
        t&, dealii::VectorTools::NormType const&, dealii::Function<2, double>
        const*, double)
        #2  ./oldroyd_bn:
        Viscoelastic::StokesProblem<2>::output_results(unsigned int)
        #3  ./oldroyd_bn: Viscoelastic::StokesProblem<2>::run()
        #4  ./oldroyd_bn: main
        --------------------------------------------------------/
        I have checked the ExactSolution in the code to try and figure wher
        the error is occurring without success.
        Can anyone point me the direction to look.


The information in this message is confidential and legally privileged. It is intended solely for the addressee(s). Access to this message by anyone else is unauthorized. If received in error, please accept our apologies and notify the sender immediately. You must also delete the original message from your machine. If you are not the intended recipient, any use, disclosure, copying, distribution or action taken in reliance of it, is prohibited and may be unlawful. The information, attachments, opinions or advice contained in this email are not the views or opinions of Harare Institute of Technology, its subsidiaries or affiliates. Although this email and any attachments are believed to be free of any virus or other defects which might affect any computer or IT system into which they are received, no responsibility is accepted by Harare Institute of Technology and/or its subsidiaries for any loss or damage arising in any way from the receipt or use thereof.

--
The deal.II project is located at http://www.dealii.org/ <https://nam10.safelinks.protection.outlook.com/?url=http%3A%2F%2Fwww.dealii.org%2F&data=05%7C01%7CWolfgang.Bangerth%40colostate.edu%7Cee52a6bc9db54c2be9ef08db1f86a028%7Cafb58802ff7a4bb1ab21367ff2ecfc8b%7C0%7C0%7C638138435321352889%7CUnknown%7CTWFpbGZsb3d8eyJWIjoiMC4wLjAwMDAiLCJQIjoiV2luMzIiLCJBTiI6Ik1haWwiLCJXVCI6Mn0%3D%7C3000%7C%7C%7C&sdata=n3Kq7VQdSickKVWZtBCtcnQxHLtGIsnGqRoYcLrMaGQ%3D&reserved=0> For mailing list/forum options, see https://groups.google.com/d/forum/dealii?hl=en <https://nam10.safelinks.protection.outlook.com/?url=https%3A%2F%2Fgroups.google.com%2Fd%2Fforum%2Fdealii%3Fhl%3Den&data=05%7C01%7CWolfgang.Bangerth%40colostate.edu%7Cee52a6bc9db54c2be9ef08db1f86a028%7Cafb58802ff7a4bb1ab21367ff2ecfc8b%7C0%7C0%7C638138435321352889%7CUnknown%7CTWFpbGZsb3d8eyJWIjoiMC4wLjAwMDAiLCJQIjoiV2luMzIiLCJBTiI6Ik1haWwiLCJXVCI6Mn0%3D%7C3000%7C%7C%7C&sdata=IbbdzWELcQPM7rM7Lf%2F29rwszV9035i9E8bvJjPcfWs%3D&reserved=0>
---
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 [email protected] <mailto:[email protected]>. To view this discussion on the web visit https://groups.google.com/d/msgid/dealii/c77b68f3-c243-49f5-977e-15d5cd183937n%40googlegroups.com <https://nam10.safelinks.protection.outlook.com/?url=https%3A%2F%2Fgroups.google.com%2Fd%2Fmsgid%2Fdealii%2Fc77b68f3-c243-49f5-977e-15d5cd183937n%2540googlegroups.com%3Futm_medium%3Demail%26utm_source%3Dfooter&data=05%7C01%7CWolfgang.Bangerth%40colostate.edu%7Cee52a6bc9db54c2be9ef08db1f86a028%7Cafb58802ff7a4bb1ab21367ff2ecfc8b%7C0%7C0%7C638138435321352889%7CUnknown%7CTWFpbGZsb3d8eyJWIjoiMC4wLjAwMDAiLCJQIjoiV2luMzIiLCJBTiI6Ik1haWwiLCJXVCI6Mn0%3D%7C3000%7C%7C%7C&sdata=QbyKrxd0xzAk8Bi6jRoZXBU6BP%2F9GEqsJ3wyGDjFUxs%3D&reserved=0>.

--
------------------------------------------------------------------------
Wolfgang Bangerth          email:                 [email protected]
                           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 the Google Groups "deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email 
to [email protected].
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/819fbd9c-55e4-d1d9-f74d-8ddc9a85d336%40colostate.edu.

Reply via email to