refineable_advection_diffusion_elements.cc
Go to the documentation of this file.
1 // LIC// ====================================================================
2 // LIC// This file forms part of oomph-lib, the object-oriented,
3 // LIC// multi-physics finite-element library, available
4 // LIC// at http://www.oomph-lib.org.
5 // LIC//
6 // LIC// Copyright (C) 2006-2023 Matthias Heil and Andrew Hazel
7 // LIC//
8 // LIC// This library is free software; you can redistribute it and/or
9 // LIC// modify it under the terms of the GNU Lesser General Public
10 // LIC// License as published by the Free Software Foundation; either
11 // LIC// version 2.1 of the License, or (at your option) any later version.
12 // LIC//
13 // LIC// This library is distributed in the hope that it will be useful,
14 // LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // LIC// Lesser General Public License for more details.
17 // LIC//
18 // LIC// You should have received a copy of the GNU Lesser General Public
19 // LIC// License along with this library; if not, write to the Free Software
20 // LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21 // LIC// 02110-1301 USA.
22 // LIC//
23 // LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24 // LIC//
25 // LIC//====================================================================
27 
28 namespace oomph
29 {
30  //==========================================================================
31  /// Add the element's contribution to the elemental residual vector
32  /// and/or elemental jacobian matrix.
33  /// This function overloads the standard version so that the possible
34  /// presence of hanging nodes is taken into account.
35  //=========================================================================
36  template<unsigned DIM>
39  Vector<double>& residuals,
40  DenseMatrix<double>& jacobian,
41  DenseMatrix<double>& mass_matrix,
42  unsigned flag)
43  {
44  // Find out how many nodes there are in the element
45  unsigned n_node = nnode();
46 
47  // Get the nodal index at which the unknown is stored
48  unsigned u_nodal_index = this->u_index_adv_diff();
49 
50  // Set up memory for the shape and test functions
51  Shape psi(n_node), test(n_node);
52  DShape dpsidx(n_node, DIM), dtestdx(n_node, DIM);
53 
54  // Set the value of n_intpt
55  unsigned n_intpt = integral_pt()->nweight();
56 
57  // Set the Vector to hold local coordinates
58  Vector<double> s(DIM);
59 
60  // Get Peclet number
61  double peclet = this->pe();
62 
63  // Get the Peclet multiplied by the Strouhal number
64  double peclet_st = this->pe_st();
65 
66  // Integers used to store the local equation number and local unknown
67  // indices for the residuals and jacobians
68  int local_eqn = 0, local_unknown = 0;
69 
70  // Local storage for pointers to hang_info objects
71  HangInfo *hang_info_pt = 0, *hang_info2_pt = 0;
72 
73  // Local variable to determine the ALE stuff
74  bool ALE_is_disabled_flag = this->ALE_is_disabled;
75 
76  // Loop over the integration points
77  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
78  {
79  // Assign values of s
80  for (unsigned i = 0; i < DIM; i++) s[i] = integral_pt()->knot(ipt, i);
81 
82  // Get the integral weight
83  double w = integral_pt()->weight(ipt);
84 
85  // Call the derivatives of the shape and test functions
86  double J = this->dshape_and_dtest_eulerian_at_knot_adv_diff(
87  ipt, psi, dpsidx, test, dtestdx);
88 
89  // Premultiply the weights and the Jacobian
90  double W = w * J;
91 
92  // Calculate local values of the function, initialise to zero
93  double dudt = 0.0;
94  double interpolated_u = 0.0;
95 
96  // These need to be a Vector to be ANSI C++, initialise to zero
97  Vector<double> interpolated_x(DIM, 0.0);
98  Vector<double> interpolated_dudx(DIM, 0.0);
99  Vector<double> mesh_velocity(DIM, 0.0);
100 
101  // Calculate function value and derivatives:
102  //-----------------------------------------
103 
104  // Loop over nodes
105  for (unsigned l = 0; l < n_node; l++)
106  {
107  // Get the value at the node
108  double u_value = this->nodal_value(l, u_nodal_index);
109  interpolated_u += u_value * psi(l);
110  dudt += this->du_dt_adv_diff(l) * psi(l);
111  // Loop over directions
112  for (unsigned j = 0; j < DIM; j++)
113  {
114  interpolated_x[j] += nodal_position(l, j) * psi(l);
115  interpolated_dudx[j] += u_value * dpsidx(l, j);
116  }
117  }
118 
119  // Get the mesh velocity, if required
120  if (!ALE_is_disabled_flag)
121  {
122  for (unsigned l = 0; l < n_node; l++)
123  {
124  // Loop over directions
125  for (unsigned j = 0; j < DIM; j++)
126  {
127  mesh_velocity[j] += dnodal_position_dt(l, j) * psi(l);
128  }
129  }
130  }
131 
132 
133  // Get body force
134  double source;
135  this->get_source_adv_diff(ipt, interpolated_x, source);
136 
137 
138  // Get wind
139  //--------
140  Vector<double> wind(DIM);
141  this->get_wind_adv_diff(ipt, s, interpolated_x, wind);
142 
143  // Assemble residuals and Jacobian
144  //================================
145 
146  // Loop over the nodes for the test functions
147  for (unsigned l = 0; l < n_node; l++)
148  {
149  // Local variables to store the number of master nodes and
150  // the weight associated with the shape function if the node is hanging
151  unsigned n_master = 1;
152  double hang_weight = 1.0;
153  // Local bool (is the node hanging)
154  bool is_node_hanging = this->node_pt(l)->is_hanging();
155 
156  // If the node is hanging, get the number of master nodes
157  if (is_node_hanging)
158  {
159  hang_info_pt = this->node_pt(l)->hanging_pt();
160  n_master = hang_info_pt->nmaster();
161  }
162  // Otherwise there is just one master node, the node itself
163  else
164  {
165  n_master = 1;
166  }
167 
168  // Loop over the number of master nodes
169  for (unsigned m = 0; m < n_master; m++)
170  {
171  // Get the local equation number and hang_weight
172  // If the node is hanging
173  if (is_node_hanging)
174  {
175  // Read out the local equation from the master node
176  local_eqn = this->local_hang_eqn(hang_info_pt->master_node_pt(m),
177  u_nodal_index);
178  // Read out the weight from the master node
179  hang_weight = hang_info_pt->master_weight(m);
180  }
181  // If the node is not hanging
182  else
183  {
184  // The local equation number comes from the node itself
185  local_eqn = this->nodal_local_eqn(l, u_nodal_index);
186  // The hang weight is one
187  hang_weight = 1.0;
188  }
189 
190  // If the nodal equation is not a boundary conditino
191  if (local_eqn >= 0)
192  {
193  // Add du/dt and body force/source term here
194  residuals[local_eqn] -=
195  (peclet_st * dudt + source) * test(l) * W * hang_weight;
196 
197  // The Mesh velocity and Advection--Diffusion bit
198  for (unsigned k = 0; k < DIM; k++)
199  {
200  // Terms that multiply the test function
201  double tmp = peclet * wind[k];
202  // If the mesh is moving need to subtract the mesh velocity
203  if (!ALE_is_disabled_flag) tmp -= peclet_st * mesh_velocity[k];
204  // Now combine the terms into the residual
205  residuals[local_eqn] -= interpolated_dudx[k] *
206  (tmp * test(l) + dtestdx(l, k)) * W *
207  hang_weight;
208  }
209 
210  // Calculate the Jacobian
211  if (flag)
212  {
213  // Local variables to store the number of master nodes
214  // and the weights associated with each hanging node
215  unsigned n_master2 = 1;
216  double hang_weight2 = 1.0;
217  // Loop over the nodes for the variables
218  for (unsigned l2 = 0; l2 < n_node; l2++)
219  {
220  // Local bool (is the node hanging)
221  bool is_node2_hanging = this->node_pt(l2)->is_hanging();
222  // If the node is hanging, get the number of master nodes
223  if (is_node2_hanging)
224  {
225  hang_info2_pt = this->node_pt(l2)->hanging_pt();
226  n_master2 = hang_info2_pt->nmaster();
227  }
228  // Otherwise there is one master node, the node itself
229  else
230  {
231  n_master2 = 1;
232  }
233 
234  // Loop over the master nodes
235  for (unsigned m2 = 0; m2 < n_master2; m2++)
236  {
237  // Get the local unknown and weight
238  // If the node is hanging
239  if (is_node2_hanging)
240  {
241  // Read out the local unknown from the master node
242  local_unknown = this->local_hang_eqn(
243  hang_info2_pt->master_node_pt(m2), u_nodal_index);
244  // Read out the hanging weight from the master node
245  hang_weight2 = hang_info2_pt->master_weight(m2);
246  }
247  // If the node is not hanging
248  else
249  {
250  // The local unknown number comes from the node
251  local_unknown = this->nodal_local_eqn(l2, u_nodal_index);
252  // The hang weight is one
253  hang_weight2 = 1.0;
254  }
255 
256  // If the unknown is not pinned
257  if (local_unknown >= 0)
258  {
259  // Add contribution to Elemental Matrix
260  // Mass matrix du/dt term
261  jacobian(local_eqn, local_unknown) -=
262  peclet_st * test(l) * psi(l2) *
263  this->node_pt(l2)->time_stepper_pt()->weight(1, 0) * W *
264  hang_weight * hang_weight2;
265 
266  // Add contribution to mass matrix
267  if (flag == 2)
268  {
269  mass_matrix(local_eqn, local_unknown) +=
270  peclet_st * test(l) * psi(l2) * W * hang_weight *
271  hang_weight2;
272  }
273 
274  // Add contribution to Elemental Matrix
275  for (unsigned k = 0; k < DIM; k++)
276  {
277  // Terms that multiply the test function
278  double tmp = peclet * wind[k];
279  // If the mesh is moving, subtract the mesh velocity
280  if (!ALE_is_disabled_flag)
281  {
282  tmp -= peclet_st * mesh_velocity[k];
283  }
284  // Construct the jacobian term
285  jacobian(local_eqn, local_unknown) -=
286  dpsidx(l2, k) * (tmp * test(l) + dtestdx(l, k)) * W *
287  hang_weight * hang_weight2;
288  }
289  }
290  } // End of loop over master nodes
291  } // End of loop over nodes
292  } // End of Jacobian calculation
293 
294  } // End of non-zero equation
295 
296  } // End of loop over the master nodes for residual
297  } // End of loop over nodes
298 
299  } // End of loop over integration points
300  }
301 
302 
303  //====================================================================
304  // Force build of templates
305  //====================================================================
309 
313 
314 } // namespace oomph
static char t char * s
Definition: cfortran.h:568
cstr elem_len * i
Definition: cfortran.h:603
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
Definition: shape.h:278
Class that contains data for hanging nodes.
Definition: nodes.h:742
double const & master_weight(const unsigned &i) const
Return weight for dofs on i-th master node.
Definition: nodes.h:808
Node *const & master_node_pt(const unsigned &i) const
Return a pointer to the i-th master node.
Definition: nodes.h:791
unsigned nmaster() const
Return the number of master nodes.
Definition: nodes.h:785
void fill_in_generic_residual_contribution_adv_diff(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Add the element's contribution to the elemental residual vector and/or Jacobian matrix flag=1: comput...
Refineable version of QAdvectionDiffusionElement. Inherit from the standard QAdvectionDiffusionElemen...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
//////////////////////////////////////////////////////////////////// ////////////////////////////////...