]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
old exception replaced
authorguido <guido@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 30 May 2002 15:07:38 +0000 (15:07 +0000)
committerguido <guido@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 30 May 2002 15:07:38 +0000 (15:07 +0000)
git-svn-id: https://svn.dealii.org/trunk@5944 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/include/numerics/matrices.h
deal.II/deal.II/source/numerics/matrices.all_dimensions.cc
deal.II/deal.II/source/numerics/matrices.cc

index fffc309c4bf1c173c1c6b3a81367f80afd5a25a5..2c3a158c0db1741683715c18061037d935652870 100644 (file)
@@ -388,14 +388,6 @@ class MatrixCreator
                                      * Exception
                                      */
     DeclException0 (ExcComponentMismatch);    
-                                    /**
-                                     * Exception
-                                     */
-    DeclException2 (ExcDimensionsDontMatch,
-                   int, int,
-                   << "The dimensions " << arg1 << " and " << arg2
-                   << " don't match.");
-
   private:
                                     /**
                                      * Convenience abbreviation for
index df2171766f7cca84864be5fc5f6834b14a612b07..d9dd753f2a880b5467234c117791bb8070101c27 100644 (file)
@@ -31,11 +31,11 @@ MatrixTools::apply_boundary_values (const std::map<unsigned int,double> &boundar
                                    const bool        preserve_symmetry)
 {
   Assert (matrix.n() == matrix.m(),
-         ExcDimensionsDontMatch(matrix.n(), matrix.m()));
+         ExcDimensionMismatch(matrix.n(), matrix.m()));
   Assert (matrix.n() == right_hand_side.size(),
-         ExcDimensionsDontMatch(matrix.n(), right_hand_side.size()));
+         ExcDimensionMismatch(matrix.n(), right_hand_side.size()));
   Assert (matrix.n() == solution.size(),
-         ExcDimensionsDontMatch(matrix.n(), solution.size()));
+         ExcDimensionMismatch(matrix.n(), solution.size()));
                                   // if no boundary values are to be applied
                                   // simply return
   if (boundary_values.size() == 0)
@@ -199,11 +199,11 @@ MatrixTools::apply_boundary_values (const std::map<unsigned int,double> &boundar
   const unsigned int blocks = matrix.n_block_rows();
   
   Assert (matrix.n() == matrix.m(),
-         ExcDimensionsDontMatch(matrix.n(), matrix.m()));
+         ExcDimensionMismatch(matrix.n(), matrix.m()));
   Assert (matrix.n() == right_hand_side.size(),
-         ExcDimensionsDontMatch(matrix.n(), right_hand_side.size()));
+         ExcDimensionMismatch(matrix.n(), right_hand_side.size()));
   Assert (matrix.n() == solution.size(),
-         ExcDimensionsDontMatch(matrix.n(), solution.size()));
+         ExcDimensionMismatch(matrix.n(), solution.size()));
   Assert (matrix.n_block_rows() == matrix.n_block_cols(),
          ExcMatrixNotBlockSquare());
   Assert (matrix.get_sparsity_pattern().get_row_indices() == 
index 0d56e2c90096e53b1e0b8a0799690c5e768c988c..aef9a7cd0cff9e5c36065ad493a63a22c2b69641 100644 (file)
@@ -149,9 +149,6 @@ void MatrixCreator::create_mass_matrix_1 (const Mapping<dim>       &mapping,
       cell_matrix.clear ();
       cell->get_dof_indices (dof_indices);
       
-      const FullMatrix<double>  &values    = fe_values.get_shape_values ();
-      const std::vector<double> &weights   = fe_values.get_JxW_values ();
-      
       if (coefficient != 0)
        {
          if (coefficient->n_components==1)
@@ -159,46 +156,66 @@ void MatrixCreator::create_mass_matrix_1 (const Mapping<dim>       &mapping,
              coefficient->value_list (fe_values.get_quadrature_points(),
                                       coefficient_values);
              for (unsigned int point=0; point<n_q_points; ++point)
-               for (unsigned int i=0; i<dofs_per_cell; ++i) 
-                 for (unsigned int j=0; j<dofs_per_cell; ++j)
-                   if ((n_components==1) ||
-                       (fe.system_to_component_index(i).first ==
-                        fe.system_to_component_index(j).first))
-                     cell_matrix(i,j) += (values(i,point) *
-                                          values(j,point) *
-                                          weights[point] *
-                                          coefficient_values[point]);
+               {
+                 const double weight = fe_values.JxW(point);
+                 for (unsigned int i=0; i<dofs_per_cell; ++i)
+                   {
+                     const double v = fe_values.shape_value(i,point);
+                     for (unsigned int j=0; j<dofs_per_cell; ++j)
+                       {
+                         const double u = fe_values.shape_value(j,point);
+                         
+                         if ((n_components==1) ||
+                             (fe.system_to_component_index(i).first ==
+                              fe.system_to_component_index(j).first))
+                           cell_matrix(i,j) +=
+                             (u * v * weight * coefficient_values[point]);
+                       }
+                   }
+               }
            }
          else
            {
              coefficient->vector_value_list (fe_values.get_quadrature_points(),
                                              coefficient_vector_values);
              for (unsigned int point=0; point<n_q_points; ++point)
-               for (unsigned int i=0; i<dofs_per_cell; ++i)
-                 {
-                   const unsigned int component_i=
-                     fe.system_to_component_index(i).first;
-                   for (unsigned int j=0; j<dofs_per_cell; ++j)
-                     if ((n_components==1) ||
-                         (fe.system_to_component_index(j).first == component_i))
-                       cell_matrix(i,j) += (values(i,point) *
-                                            values(j,point) *
-                                            weights[point] *
-                                            coefficient_vector_values[point](component_i));    
-                 }
+               {
+                 const double weight = fe_values.JxW(point);
+                 for (unsigned int i=0; i<dofs_per_cell; ++i)
+                   {
+                     const double v = fe_values.shape_value(i,point);
+                     const unsigned int component_i=
+                       fe.system_to_component_index(i).first;
+                     for (unsigned int j=0; j<dofs_per_cell; ++j)
+                       {
+                         const double u = fe_values.shape_value(j,point);
+                         if ((n_components==1) ||
+                             (fe.system_to_component_index(j).first == component_i))
+                           cell_matrix(i,j) +=
+                             (u * v * weight *
+                              coefficient_vector_values[point](component_i));
+                       }
+                   }
+               }
            }
        }
       else
        for (unsigned int point=0; point<n_q_points; ++point)
-         for (unsigned int i=0; i<dofs_per_cell; ++i) 
-           for (unsigned int j=0; j<dofs_per_cell; ++j)
-             if ((n_components==1) ||
-                 (fe.system_to_component_index(i).first ==
-                  fe.system_to_component_index(j).first))
-               cell_matrix(i,j) += (values(i,point) *
-                                    values(j,point) *
-                                    weights[point]);
-
+         {
+           const double weight = fe_values.JxW(point);
+           for (unsigned int i=0; i<dofs_per_cell; ++i)
+             {
+               const double v = fe_values.shape_value(i,point); 
+               for (unsigned int j=0; j<dofs_per_cell; ++j)
+                 {
+                   const double u = fe_values.shape_value(j,point);
+                   if ((n_components==1) ||
+                       (fe.system_to_component_index(i).first ==
+                        fe.system_to_component_index(j).first))
+                     cell_matrix(i,j) += (u * v * weight);
+                 }
+             }
+         }
                                       // transfer everything into the
                                       // global object
       mutex.acquire ();

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.