darcy_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-2024 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//====================================================================
26 #include "darcy_elements.h"
27 
28 namespace oomph
29 {
30  //=====================================================================
31  /// Performs a div-conserving transformation of the vector basis
32  /// functions from the reference element to the actual element
33  //=====================================================================
34  template<unsigned DIM>
36  const Shape& q_basis_local,
37  Shape& psi,
38  Shape& q_basis) const
39  {
40  // Get the number of nodes in the element
41  const unsigned n_node = this->nnode();
42 
43  // Storage for derivatives of the (geometric) shape functions, and call
44  // the shape functions
45  DShape dpsi(n_node, DIM);
46  this->dshape_local(s, psi, dpsi);
47 
48  // Storage for the (geometric) jacobian and its inverse
49  DenseMatrix<double> jacobian(DIM), inverse_jacobian(DIM);
50 
51  // Get the jacobian of the geometric mapping and its determinant
52  double det = local_to_eulerian_mapping(dpsi, jacobian, inverse_jacobian);
53 
54  // Get the number of computational basis vectors
55  const unsigned n_q_basis = this->nq_basis();
56 
57  // Loop over the basis vectors
58  for (unsigned l = 0; l < n_q_basis; l++)
59  {
60  // Loop over the spatial components
61  for (unsigned i = 0; i < DIM; i++)
62  {
63  // Zero the basis
64  q_basis(l, i) = 0.0;
65  }
66  }
67 
68  // Loop over the spatial components
69  for (unsigned i = 0; i < DIM; i++)
70  {
71  // And again
72  for (unsigned j = 0; j < DIM; j++)
73  {
74  // Get the element of the jacobian (must transpose it due to different
75  // conventions) and divide by the determinant
76  double jac_trans = jacobian(j, i) / det;
77 
78  // Loop over the computational basis vectors
79  for (unsigned l = 0; l < n_q_basis; l++)
80  {
81  // Apply the appropriate div-conserving mapping
82  q_basis(l, i) += jac_trans * q_basis_local(l, j);
83  }
84  }
85  }
86 
87  // Scale the basis by the ratio of the length of the edge to the length of
88  // the corresponding edge on the reference element
89  scale_basis(q_basis);
90 
91  return det;
92  }
93 
94  //=====================================================================
95  /// Output FE representation of soln: x,y,q1,q2,div_q,p,q \cdot n
96  //=====================================================================
97  template<unsigned DIM>
99  std::ostream& outfile,
100  const unsigned& nplot,
101  const Vector<double> unit_normal)
102  {
103  // Vector of local coordinates
104  Vector<double> s(DIM);
105 
106  // Tecplot header info
107  outfile << this->tecplot_zone_string(nplot);
108 
109  // Loop over plot points
110  unsigned num_plot_points = this->nplot_points(nplot);
111 
112  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
113  {
114  // Get local coordinates of plot point
115  this->get_s_plot(iplot, nplot, s);
116 
117  // Output the components of the position
118  for (unsigned i = 0; i < DIM; i++)
119  {
120  outfile << this->interpolated_x(s, i) << " ";
121  }
122 
123  // Output the components of the FE representation of q at s
124  for (unsigned i = 0; i < DIM; i++)
125  {
126  outfile << this->interpolated_q(s, i) << " ";
127  }
128 
129  // Output FE representation of div q at s
130  outfile << this->interpolated_div_q(s) << " ";
131 
132  // Output FE representation of p at s
133  outfile << this->interpolated_p(s) << " ";
134 
135  // Fluxes projected into the direction of the face normal
136  double flux = 0.0;
137  for (unsigned i = 0; i < 2; i++)
138  {
139  flux += this->interpolated_q(s, i) * unit_normal[i];
140  }
141  outfile << flux << " ";
142 
143  outfile << std::endl;
144  }
145 
146  // Write tecplot footer (e.g. FE connectivity lists)
147  this->write_tecplot_zone_footer(outfile, nplot);
148  }
149 
150 
151  //=====================================================================
152  /// Output FE representation of soln: x,y,q1,q2,div_q,p at
153  /// Nplot^DIM plot points
154  //=====================================================================
155  template<unsigned DIM>
156  void DarcyEquations<DIM>::output(std::ostream& outfile, const unsigned& nplot)
157  {
158  // Vector of local coordinates
159  Vector<double> s(DIM);
160 
161  // Tecplot header info
162  outfile << tecplot_zone_string(nplot);
163 
164  // Loop over plot points
165  unsigned num_plot_points = nplot_points(nplot);
166 
167  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
168  {
169  // Get local coordinates of plot point
170  get_s_plot(iplot, nplot, s);
171 
172  // Output the components of the position
173  for (unsigned i = 0; i < DIM; i++)
174  {
175  outfile << interpolated_x(s, i) << " ";
176  }
177 
178  // Output the components of the FE representation of q at s
179  for (unsigned i = 0; i < DIM; i++)
180  {
181  outfile << interpolated_q(s, i) << " ";
182  }
183 
184  // Output FE representation of div q at s
185  outfile << interpolated_div_q(s) << " ";
186 
187  // Output FE representation of p at s
188  outfile << interpolated_p(s) << " ";
189 
190  outfile << std::endl;
191  }
192 
193  // Write tecplot footer (e.g. FE connectivity lists)
194  this->write_tecplot_zone_footer(outfile, nplot);
195  }
196 
197  //=====================================================================
198  /// Output FE representation of exact soln: x,y,q1,q2,div_q,p at
199  /// Nplot^DIM plot points
200  //=====================================================================
201  template<unsigned DIM>
203  std::ostream& outfile,
204  const unsigned& nplot,
206  {
207  // Vector of local coordinates
208  Vector<double> s(DIM);
209 
210  // Vector for coordintes
211  Vector<double> x(DIM);
212 
213  // Tecplot header info
214  outfile << this->tecplot_zone_string(nplot);
215 
216  // Exact solution Vector
217  Vector<double> exact_soln(DIM + 2);
218 
219  // Loop over plot points
220  unsigned num_plot_points = this->nplot_points(nplot);
221 
222  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
223  {
224  // Get local coordinates of plot point
225  this->get_s_plot(iplot, nplot, s);
226 
227  // Get x position as Vector
228  this->interpolated_x(s, x);
229 
230  // Get exact solution at this point
231  (*exact_soln_pt)(x, exact_soln);
232 
233  // Output x,y,q_exact,p_exact,div_q_exact
234  for (unsigned i = 0; i < DIM; i++)
235  {
236  outfile << x[i] << " ";
237  }
238  for (unsigned i = 0; i < DIM + 2; i++)
239  {
240  outfile << exact_soln[i] << " ";
241  }
242  outfile << std::endl;
243  }
244 
245  // Write tecplot footer (e.g. FE connectivity lists)
246  this->write_tecplot_zone_footer(outfile, nplot);
247  }
248 
249  //=====================================================================
250  /// Compute the error between the FE solution and the exact solution
251  /// using the H(div) norm for q and L^2 norm for p
252  //=====================================================================
253  template<unsigned DIM>
255  std::ostream& outfile,
257  Vector<double>& error,
258  Vector<double>& norm)
259  {
260  for (unsigned i = 0; i < 2; i++)
261  {
262  error[i] = 0.0;
263  norm[i] = 0.0;
264  }
265 
266  // Vector of local coordinates
267  Vector<double> s(DIM);
268 
269  // Vector for coordinates
270  Vector<double> x(DIM);
271 
272  // Set the value of n_intpt
273  unsigned n_intpt = this->integral_pt()->nweight();
274 
275  outfile << "ZONE" << std::endl;
276 
277  // Exact solution Vector (u,v,[w])
278  Vector<double> exact_soln(DIM + 2);
279 
280  // Loop over the integration points
281  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
282  {
283  // Assign values of s
284  for (unsigned i = 0; i < DIM; i++)
285  {
286  s[i] = this->integral_pt()->knot(ipt, i);
287  }
288 
289  // Get the integral weight
290  double w = this->integral_pt()->weight(ipt);
291 
292  // Get jacobian of mapping
293  double J = this->J_eulerian(s);
294 
295  // Premultiply the weights and the Jacobian
296  double W = w * J;
297 
298  // Get x position as Vector
299  this->interpolated_x(s, x);
300 
301  // Get exact solution at this point
302  (*exact_soln_pt)(x, exact_soln);
303 
304  // Flux error
305  for (unsigned i = 0; i < DIM; i++)
306  {
307  norm[0] += exact_soln[i] * exact_soln[i] * W;
308  // Error due to q_i
309  error[0] += (exact_soln[i] - this->interpolated_q(s, i)) *
310  (exact_soln[i] - this->interpolated_q(s, i)) * W;
311  }
312 
313  // // Flux divergence error
314  // norm[0]+=exact_soln[DIM]*exact_soln[DIM]*W;
315  // error[0]+=(exact_soln[DIM]-interpolated_div_q(s))*
316  // (exact_soln[DIM]-interpolated_div_q(s))*W;
317 
318  // Pressure error
319  norm[1] += exact_soln[DIM + 1] * exact_soln[DIM + 1] * W;
320  error[1] += (exact_soln[DIM + 1] - this->interpolated_p(s)) *
321  (exact_soln[DIM + 1] - this->interpolated_p(s)) * W;
322 
323  // Output x,y,[z]
324  for (unsigned i = 0; i < DIM; i++)
325  {
326  outfile << x[i] << " ";
327  }
328 
329  // Output q_1_error,q_2_error,...
330  for (unsigned i = 0; i < DIM; i++)
331  {
332  outfile << exact_soln[i] - this->interpolated_q(s, i) << " ";
333  }
334 
335  // Output p_error
336  outfile << exact_soln[DIM + 1] - this->interpolated_p(s) << " ";
337 
338  outfile << std::endl;
339  }
340  }
341 
342  //=====================================================================
343  /// Fill in residuals and, if flag==true, jacobian
344  //=====================================================================
345  template<unsigned DIM>
347  Vector<double>& residuals, DenseMatrix<double>& jacobian, bool flag)
348  {
349  // Get the number of geometric nodes, total number of basis functions,
350  // and number of edges basis functions
351  const unsigned n_node = nnode();
352  const unsigned n_q_basis = nq_basis();
353  const unsigned n_q_basis_edge = nq_basis_edge();
354  const unsigned n_p_basis = np_basis();
355 
356  // Storage for the geometric and computational bases and the test functions
357  Shape psi(n_node), q_basis(n_q_basis, DIM), q_test(n_q_basis, DIM),
358  p_basis(n_p_basis), p_test(n_p_basis), div_q_basis_ds(n_q_basis),
359  div_q_test_ds(n_q_basis);
360 
361  // Get the number of integration points
362  unsigned n_intpt = integral_pt()->nweight();
363 
364  // Storage for the local coordinates
365  Vector<double> s(DIM);
366 
367  // Storage for the source function
368  Vector<double> f(DIM);
369 
370  // Storage for the mass source function
371  double mass_source_local;
372 
373  // Local equation and unknown numbers
374  int local_eqn = 0; //, local_unknown = 0;
375 
376  // Loop over the integration points
377  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
378  {
379  // Find the local coordinates at the integration point
380  for (unsigned i = 0; i < DIM; i++)
381  {
382  s[i] = integral_pt()->knot(ipt, i);
383  }
384 
385  // Get the weight of the intetgration point
386  double w = integral_pt()->weight(ipt);
387 
388  // Call the basis functions and test functions and get the
389  // (geometric) jacobian of the current element
390  double J = shape_basis_test_local_at_knot(ipt,
391  psi,
392  q_basis,
393  q_test,
394  p_basis,
395  p_test,
396  div_q_basis_ds,
397  div_q_test_ds);
398 
399  // Storage for interpolated values
400  Vector<double> interpolated_x(DIM, 0.0);
401  Vector<double> interpolated_q(DIM, 0.0);
402  double interpolated_div_q_ds = 0.0;
403  double interpolated_p = 0.0;
404 
405  // loop over the geometric basis functions to find interpolated x
406  for (unsigned l = 0; l < n_node; l++)
407  {
408  for (unsigned i = 0; i < DIM; i++)
409  {
410  interpolated_x[i] += nodal_position(l, i) * psi[l];
411  }
412  }
413 
414  // loop over the nodes and use the vector basis functions to find the
415  // interpolated flux
416  for (unsigned i = 0; i < DIM; i++)
417  {
418  // Loop over the edge basis vectors
419  for (unsigned l = 0; l < n_q_basis_edge; l++)
420  {
421  interpolated_q[i] += q_edge(l) * q_basis(l, i);
422  }
423  // Loop over the internal basis vectors
424  for (unsigned l = n_q_basis_edge; l < n_q_basis; l++)
425  {
426  interpolated_q[i] += q_internal(l - n_q_basis_edge) * q_basis(l, i);
427  }
428  }
429 
430  // loop over the pressure basis and find the interpolated pressure
431  for (unsigned l = 0; l < n_p_basis; l++)
432  {
433  interpolated_p += p_value(l) * p_basis(l);
434  }
435 
436  // loop over the q edge divergence basis and the q internal divergence
437  // basis to find interpolated div q
438  for (unsigned l = 0; l < n_q_basis_edge; l++)
439  {
440  interpolated_div_q_ds += q_edge(l) * div_q_basis_ds(l);
441  }
442  for (unsigned l = n_q_basis_edge; l < n_q_basis; l++)
443  {
444  interpolated_div_q_ds +=
445  q_internal(l - n_q_basis_edge) * div_q_basis_ds(l);
446  }
447 
448  // Get the source function
449  this->source(interpolated_x, f);
450 
451  // Get the mass source function
452  this->mass_source(interpolated_x, mass_source_local);
453 
454  // Loop over the test functions
455  for (unsigned l = 0; l < n_q_basis; l++)
456  {
457  if (l < n_q_basis_edge)
458  {
459  local_eqn = q_edge_local_eqn(l);
460  }
461  else // n_q_basis_edge <= l < n_basis
462  {
463  local_eqn = q_internal_local_eqn(l - n_q_basis_edge);
464  }
465 
466  // If it's not a boundary condition
467  if (local_eqn >= 0)
468  {
469  for (unsigned i = 0; i < DIM; i++)
470  {
471  residuals[local_eqn] +=
472  (interpolated_q[i] - f[i]) * q_test(l, i) * w * J;
473  }
474 
475  // deliberately no jacobian factor in this integral
476  residuals[local_eqn] -= (interpolated_p * div_q_test_ds(l)) * w;
477  }
478  } // End of loop over test functions
479 
480  // loop over pressure test functions
481  for (unsigned l = 0; l < n_p_basis; l++)
482  {
483  // get the local equation number
484  local_eqn = p_local_eqn(l);
485 
486  // If it's not a boundary condition
487  if (local_eqn >= 0)
488  {
489  // deliberately no jacobian factor in this integral
490  residuals[local_eqn] += interpolated_div_q_ds * p_test(l) * w;
491  residuals[local_eqn] -= mass_source_local * p_test(l) * w * J;
492  }
493  }
494  } // End of loop over integration points
495  }
496 
497 
498  //=====================================================================
499  // Force building of templates
500  //=====================================================================
501  template class DarcyEquations<2>;
502 
503 } // 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
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output FE representation of exact soln: x,y,q1,q2,div_q,p at Nplot^DIM plot points.
void output_with_projected_flux(std::ostream &outfile, const unsigned &nplot, const Vector< double > unit_normal)
Output incl. projection of fluxes into direction of the specified unit vector.
virtual void fill_in_generic_residual_contribution(Vector< double > &residuals, DenseMatrix< double > &jacobian, bool flag)
Fill in residuals and, if flag==true, jacobian.
double transform_basis(const Vector< double > &s, const Shape &q_basis_local, Shape &psi, Shape &q_basis) const
Performs a div-conserving transformation of the vector basis functions from the reference element to ...
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm)
Compute the error between the FE solution and the exact solution using the H(div) norm for q and L^2 ...
void output(std::ostream &outfile)
Output with default number of plot points.
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
Definition: elements.h:1763
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
//////////////////////////////////////////////////////////////////// ////////////////////////////////...