]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
fix Kovasnzay flow
authorkanschat <kanschat@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 26 Sep 2007 04:03:31 +0000 (04:03 +0000)
committerkanschat <kanschat@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 26 Sep 2007 04:03:31 +0000 (04:03 +0000)
git-svn-id: https://svn.dealii.org/trunk@15244 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/base/include/base/flow_function.h
deal.II/base/source/flow_function.cc

index ca49d3595599edd5a40d97a80fd78885e7e33da0..794c1f9dbd752e89e091a2f21b16f98c38b97662 100644 (file)
@@ -226,7 +226,6 @@ namespace Functions
                                     std::vector<std::vector<Tensor<1,2> > >& gradients) const;
       virtual void vector_laplacians (const std::vector<Point<2> > &points,
                                      std::vector<std::vector<double> >   &values) const;
-
     private:
                                       /// The auxiliary function Psi.
       double Psi(double phi) const;
@@ -286,9 +285,11 @@ namespace Functions
       virtual void vector_laplacians (const std::vector<Point<2> > &points,
                                      std::vector<std::vector<double> >   &values) const;
 
+                                      /// The value of lambda.
+      double lambda () const;
     private:
       const double Reynolds;
-      double lambda;
+      double lbda;
       double p_average;
       const bool stokes;
   };
index 8c219bb703a77d98b9c4dab765ee9e387677cf27..74412a12b0e6fa41f4d57b7e118cdeca953528f3 100644 (file)
@@ -613,7 +613,7 @@ namespace Functions
     long double r2 = Reynolds/2.;
     long double b = 4*M_PI*M_PI;
     long double l = -b/(r2+std::sqrt(r2*r2+b));
-    lambda = l;
+    lbda = l;
                                     // mean pressure for a domain
                                     // spreading from -.5 to 1.5 in
                                     // x-direction
@@ -640,11 +640,11 @@ namespace Functions
        const Point<2>& p = points[k];
        const double x = p(0);
        const double y = 2. * deal_II_numbers::PI * p(1);
-       const double elx = std::exp(lambda*x);
+       const double elx = std::exp(lbda*x);
        
        values[0][k] = 1. - elx * cos(y);
-       values[1][k] = .5 / deal_II_numbers::PI * lambda * elx * sin(y);
-       values[2][k] = .5 * elx * elx - p_average - this->mean_pressure;
+       values[1][k] = .5 / deal_II_numbers::PI * lbda * elx * sin(y);
+       values[2][k] = -.5 * elx * elx + p_average - this->mean_pressure;
       }
   }
   
@@ -653,29 +653,28 @@ namespace Functions
     const std::vector<Point<2> >& points,
     std::vector<std::vector<Tensor<1,2> > >& gradients) const
   {
-    Assert(false, ExcNotImplemented());
     unsigned int n = points.size();
     
-    Assert (gradients.size() == n, ExcDimensionMismatch(gradients.size(), n));
-    Assert (gradients[0].size() >= this->n_components,
-           ExcDimensionMismatch(gradients[0].size(), this->n_components));
+    Assert (gradients.size() == 3, ExcDimensionMismatch(gradients.size(), 3));
+    Assert (gradients[0].size() == n,
+           ExcDimensionMismatch(gradients[0].size(), n));
     
     for (unsigned int i=0;i<n;++i)
       {
        const double x = points[i](0);
        const double y = points[i](1);
        
-       const double elx = std::exp(lambda*x);
+       const double elx = std::exp(lbda*x);
        const double cy = cos(2*M_PI*y);
        const double sy = sin(2*M_PI*y);
        
                                         // u
-       gradients[0][i][0] = -lambda*elx*cy;
-       gradients[0][i][1] = 2*M_PI*elx*sy;
-       gradients[1][i][0] = lambda*lambda/(2*M_PI)*elx*sy;
-       gradients[1][i][1] =lambda*elx*cy;
+       gradients[0][i][0] = -lbda*elx*cy;
+       gradients[0][i][1] = 2. * deal_II_numbers::PI*elx*sy;
+       gradients[1][i][0] = lbda*lbda/(2*M_PI)*elx*sy;
+       gradients[1][i][1] =lbda*elx*cy;
                                         // p
-       gradients[2][i][0] = -lambda*elx*elx;
+       gradients[2][i][0] = -lbda*elx*elx;
        gradients[2][i][1] = 0.;
     }
   }
@@ -699,13 +698,13 @@ namespace Functions
            const Point<2>& p = points[k];
            const double x = p(0);
            const double y = zp * p(1);
-           const double elx = std::exp(lambda*x);
+           const double elx = std::exp(lbda*x);
            const double u  = 1. - elx * cos(y);
-           const double ux = -lambda * elx * cos(y);
+           const double ux = -lbda * elx * cos(y);
            const double uy = elx * zp * sin(y);
-           const double v  = lambda/zp * elx * sin(y);
-           const double vx = lambda*lambda/zp * elx * sin(y);
-           const double vy = zp*lambda/zp * elx * cos(y);
+           const double v  = lbda/zp * elx * sin(y);
+           const double vx = lbda*lbda/zp * elx * sin(y);
+           const double vy = zp*lbda/zp * elx * cos(y);
            
            values[0][k] = u*ux+v*uy;
            values[1][k] = u*vx+v*vy;
@@ -720,6 +719,11 @@ namespace Functions
       }
   }
   
+  double
+  Kovasznay::lambda () const
+  {
+    return lbda;
+  }
   
   
   

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.