Dear all

I have a derived class NACA0012 which is shown below.

In the newer version of deal.II, I had to add a clone function. Could you tell 
me if this is correct and what is the purpose of this ?

I also have a project_to_manifold function which gives warnings while compiling


[ 10%] Building CXX object CMakeFiles/hyflo.dir/assemble.cc.o
In file included from 
/home/praveen/Applications/hyflo/examples/2d/naca0012/src/equation.h:22,
                 from 
/home/praveen/Applications/hyflo/examples/2d/naca0012/src/parameters.h:42,
                 from 
/home/praveen/Applications/hyflo/examples/2d/naca0012/src/claw.h:38,
                 from 
/home/praveen/Applications/hyflo/examples/2d/naca0012/src/assemble.cc:27:
/home/praveen/Applications/hyflo/examples/2d/naca0012/src/user.h:76:14: 
warning:   by ‘dealii::Point<2, double> NACA0012::project_to_manifold(const 
std::vector<dealii::Point<2, double>, std::allocator<dealii::Point<2, double> > 
>&, const dealii::Point<2, double>&) const’ [-Woverloaded-virtual]
     Point<2> project_to_manifold (const std::vector<Point<2>> 
&surrounding_points,

Can this be ignored or how to fix this ? Or am I doing something wrong in my 
code ?

Thanks
praveen

class NACA0012 : public Manifold<2,2>
{
private:
   // used in newton method for projection to manifold 
   void func(const double &m, const Point<2> &P, const double &x, double &f, 
double &d) const
   {  
      // set and mark independent variable
      Sacado::Fad::DFad<double> xd = x;
      xd.diff(0,1);
      
      Sacado::Fad::DFad<double> y = naca0012 (xd);
      if(P[1] < 0) y = -y;
      Sacado::Fad::DFad<double> F = m * (y - P[1]) + (xd - P[0]);
      f = F.val();
      d = F.fastAccessDx(0);
   }
   
public:
   virtual std::unique_ptr<Manifold<2,2>> clone() const override
   {  
      return std_cxx14::make_unique<NACA0012>();
   }
   
   virtual Point<2> project_to_manifold (const std::vector<Point<2>> 
&surrounding_points,
                                         const Point<2>              
&candidate) const
   {  
      const double GEOMTOL = 1.0e-13;
      Assert((surrounding_points[0][1] > -GEOMTOL && surrounding_points[1][1] > 
-GEOMTOL) ||
             (surrounding_points[0][1] <  GEOMTOL && surrounding_points[1][1] < 
 GEOMTOL),
             ExcMessage("End points not on same side of airfoil"));
      
      const double dx = surrounding_points[1][0] - surrounding_points[0][0];
      const double dy = surrounding_points[1][1] - surrounding_points[0][1];
      Assert (std::fabs(dx) > GEOMTOL, ExcMessage("dx is too small"));
      const double m  = dy/dx;
      
      // initial guess for newton
      double x = candidate[0];
      x = std::min(x, 1.0);
      x = std::max(x, 0.0);
      double f, d;
      func(m, candidate, x, f, d);
      unsigned int i = 0, ITMAX = 10;
      while (std::fabs(f) > GEOMTOL)
      {  
         double s = -f/d; 
         while(x+s < 0 || x+s > 1) s *= 0.5;
         x = x + s;
         func(m, candidate, x, f, d);
         ++i;
         AssertThrow(i < ITMAX, ExcMessage("Newton did not converge"));
      }
      
      double y = naca0012(x);
      if(candidate[1] < 0) y = -y;
      Point<2> p (x, y);
      return p;
   }
};

-- 
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].
For more options, visit https://groups.google.com/d/optout.

Reply via email to