From: Martin Kronbichler Date: Mon, 21 Apr 2008 15:42:19 +0000 (+0000) Subject: Implemented a faster matrix assembly routine, following the one in step-22. X-Git-Tag: v8.0.0~9207 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=f8cc05a619ee86ff86fda0e53d61f374359f38a6;p=dealii.git Implemented a faster matrix assembly routine, following the one in step-22. git-svn-id: https://svn.dealii.org/trunk@16005 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/deal.II/examples/step-31/step-31.cc b/deal.II/examples/step-31/step-31.cc index e4f3b1368c..c7b2bbb6cf 100644 --- a/deal.II/examples/step-31/step-31.cc +++ b/deal.II/examples/step-31/step-31.cc @@ -501,6 +501,13 @@ void BoussinesqFlowProblem::assemble_system () const double Raleigh_number = 10; + std::vector > phi_u (dofs_per_cell); + std::vector > phi_grads_u (dofs_per_cell); + std::vector div_phi_u (dofs_per_cell); + std::vector phi_p (dofs_per_cell); + std::vector phi_T (dofs_per_cell); + std::vector > grad_phi_T (dofs_per_cell); + const FEValuesExtractors::Vector velocities (0); const FEValuesExtractors::Scalar pressure (dim); const FEValuesExtractors::Scalar temperature (dim+1); @@ -520,47 +527,54 @@ void BoussinesqFlowProblem::assemble_system () { const double old_temperature = old_solution_values[q](dim+1); + // Extract the basis relevant + // terms in the inner products + // once in advance as shown + // in step-22. This accelerates + // the assembly process, + for (unsigned int k=0; k phi_i_u = fe_values[velocities].value (i, q); if (rebuild_matrices) - { - const SymmetricTensor<2,dim> - phi_i_grads_u = fe_values[velocities].symmetric_gradient (i, q); - const Tensor<1,dim> phi_i_u = fe_values[velocities].value (i, q); - const double div_phi_i_u = fe_values[velocities].divergence (i, q); - const double phi_i_p = fe_values[pressure].value (i, q); - const double phi_i_T = fe_values[temperature].value (i, q); - const Tensor<1,dim> grad_phi_i_T = fe_values[temperature].gradient(i, q); - - for (unsigned int j=0; j - phi_j_grads_u = fe_values[velocities].symmetric_gradient (j, q); - const double div_phi_j_u = fe_values[velocities].divergence (j, q); - const double phi_j_p = fe_values[pressure].value (j, q); - const double phi_j_T = fe_values[temperature].value (j, q); - - local_matrix(i,j) += (phi_i_grads_u * phi_j_grads_u - - div_phi_i_u * phi_j_p - - phi_i_p * div_phi_j_u - + phi_i_p * phi_j_p - + phi_i_T * phi_j_T) - * fe_values.JxW(q); - } - } + for (unsigned int j=0; j gravity (0,1); local_rhs(i) += (Raleigh_number * - gravity * phi_i_u * old_temperature)* + gravity * phi_u[i] * old_temperature)* fe_values.JxW(q); } } + // The assembly of the face + // cells which enters the + // right hand sides cannot + // be accelerated with the + // above technique, since + // all the basis functions are + // only evaluated once. for (unsigned int face_no=0; face_no::faces_per_cell; ++face_no)