1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2020 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17
18
19
20 // <h1>Adaptivity Example 2 - Solving a Transient System with Adaptive Mesh Refinement</h1>
21 // \author Benjamin S. Kirk
22 // \date 2003
23 //
24 // This example shows how a simple, linear transient
25 // system can be solved in parallel. The system is simple
26 // scalar convection-diffusion with a specified external
27 // velocity. The initial condition is given, and the
28 // solution is advanced in time with a standard Crank-Nicolson
29 // time-stepping strategy.
30 //
31 // Also, we use this example to demonstrate writing out and reading
32 // in of solutions. We do 25 time steps, then save the solution
33 // and do another 25 time steps starting from the saved solution.
34
35 // C++ include files that we need
36 #include <iostream>
37 #include <algorithm>
38 #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC
39 #include <cmath>
40
41 // Basic include file needed for the mesh functionality.
42 #include "libmesh/libmesh.h"
43 #include "libmesh/replicated_mesh.h"
44 #include "libmesh/mesh_refinement.h"
45 #include "libmesh/gmv_io.h"
46 #include "libmesh/equation_systems.h"
47 #include "libmesh/fe.h"
48 #include "libmesh/quadrature_gauss.h"
49 #include "libmesh/dof_map.h"
50 #include "libmesh/sparse_matrix.h"
51 #include "libmesh/numeric_vector.h"
52 #include "libmesh/dense_matrix.h"
53 #include "libmesh/dense_vector.h"
54 #include "libmesh/getpot.h"
55 #include "libmesh/enum_solver_package.h"
56 #include "libmesh/enum_xdr_mode.h"
57 #include "libmesh/enum_norm_type.h"
58
59 // This example will solve a linear transient system,
60 // so we need to include the TransientLinearImplicitSystem definition.
61 #include "libmesh/transient_system.h"
62 #include "libmesh/linear_implicit_system.h"
63 #include "libmesh/vector_value.h"
64
65 // To refine the mesh we need an ErrorEstimator
66 // object to figure out which elements to refine.
67 #include "libmesh/error_vector.h"
68 #include "libmesh/kelly_error_estimator.h"
69
70 // The definition of a geometric element
71 #include "libmesh/elem.h"
72
73 // Bring in everything from the libMesh namespace
74 using namespace libMesh;
75
76 // Function prototype. This function will assemble the system
77 // matrix and right-hand-side at each time step. Note that
78 // since the system is linear we technically do not need to
79 // assemble the matrix at each time step, but we will anyway.
80 // In subsequent examples we will employ adaptive mesh refinement,
81 // and with a changing mesh it will be necessary to rebuild the
82 // system matrix.
83 #ifdef LIBMESH_ENABLE_AMR
84 void assemble_cd (EquationSystems & es,
85 const std::string & system_name);
86 #endif
87
88 // Function prototype. This function will initialize the system.
89 // Initialization functions are optional for systems. They allow
90 // you to specify the initial values of the solution. If an
91 // initialization function is not provided then the default (0)
92 // solution is provided.
93 void init_cd (EquationSystems & es,
94 const std::string & system_name);
95
96 // Exact solution function prototype. This gives the exact
97 // solution as a function of space and time. In this case the
98 // initial condition will be taken as the exact solution at time 0,
99 // as will the Dirichlet boundary conditions at time t.
100 Real exact_solution (const Real x,
101 const Real y,
102 const Real t);
103
exact_value(const Point & p,const Parameters & parameters,const std::string &,const std::string &)104 Number exact_value (const Point & p,
105 const Parameters & parameters,
106 const std::string &,
107 const std::string &)
108 {
109 return exact_solution(p(0), p(1), parameters.get<Real> ("time"));
110 }
111
112
113
114 // Begin the main program. Note that the first
115 // statement in the program throws an error if
116 // you are in complex number mode, since this
117 // example is only intended to work with real
118 // numbers.
main(int argc,char ** argv)119 int main (int argc, char ** argv)
120 {
121 // Initialize libMesh.
122 LibMeshInit init (argc, argv);
123
124 // This example requires a linear solver package.
125 libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
126 "--enable-petsc, --enable-trilinos, or --enable-eigen");
127
128 #ifndef LIBMESH_ENABLE_AMR
129 libmesh_example_requires(false, "--enable-amr");
130 #else
131
132 // Our Trilinos interface does not yet support adaptive transient
133 // problems
134 libmesh_example_requires(libMesh::default_solver_package() != TRILINOS_SOLVERS, "--enable-petsc");
135
136 // Brief message to the user regarding the program name
137 // and command line arguments.
138
139 // Use commandline parameter to specify if we are to
140 // read in an initial solution or generate it ourself
141 libMesh::out << "Usage:\n"
142 <<"\t " << argv[0] << " -init_timestep 0 -n_timesteps 25 [-n_refinements 5]\n"
143 << "OR\n"
144 <<"\t " << argv[0] << " -read_solution -init_timestep 26 -n_timesteps 25\n"
145 << std::endl;
146
147 libMesh::out << "Running: " << argv[0];
148
149 for (int i=1; i<argc; i++)
150 libMesh::out << " " << argv[i];
151
152 libMesh::out << std::endl << std::endl;
153
154 // Create a GetPot object to parse the command line
155 GetPot command_line (argc, argv);
156
157
158 // This boolean value is obtained from the command line, it is true
159 // if the flag "-read_solution" is present, false otherwise.
160 // It indicates whether we are going to read in
161 // the mesh and solution files "saved_mesh.xda" and "saved_solution.xda"
162 // or whether we are going to start from scratch by just reading
163 // "mesh.xda"
164 const bool read_solution = command_line.search("-read_solution");
165
166 // This value is also obtained from the commandline and it specifies the
167 // initial value for the t_step looping variable. We must
168 // distinguish between the two cases here, whether we read in the
169 // solution or we started from scratch, so that we do not overwrite the
170 // gmv output files.
171 unsigned int init_timestep = 0;
172
173 // Search the command line for the "init_timestep" flag and if it is
174 // present, set init_timestep accordingly.
175 if (command_line.search("-init_timestep"))
176 init_timestep = command_line.next(0);
177 else
178 {
179 // This handy function will print the file name, line number,
180 // specified message, and then throw an exception.
181 libmesh_error_msg("ERROR: Initial timestep not specified!");
182 }
183
184 // This value is also obtained from the command line, and specifies
185 // the number of time steps to take.
186 unsigned int n_timesteps = 0;
187
188 // Again do a search on the command line for the argument
189 if (command_line.search("-n_timesteps"))
190 n_timesteps = command_line.next(0);
191 else
192 libmesh_error_msg("ERROR: Number of timesteps not specified");
193
194
195 // Skip this 2D example if libMesh was compiled as 1D-only.
196 libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
197
198 // Create a new mesh on the default MPI communicator.
199 // We still need some work on automatic parallel restarts with
200 // DistributedMesh
201 ReplicatedMesh mesh(init.comm());
202
203 // Create an equation systems object.
204 EquationSystems equation_systems (mesh);
205 MeshRefinement mesh_refinement (mesh);
206
207 // First we process the case where we do not read in the solution
208 if (!read_solution)
209 {
210 // Read the mesh from file.
211 mesh.read ("mesh.xda");
212
213 // Again do a search on the command line for an argument
214 unsigned int n_refinements = 5;
215 if (command_line.search("-n_refinements"))
216 n_refinements = command_line.next(0);
217
218 // Uniformly refine the mesh 5 times
219 if (!read_solution)
220 mesh_refinement.uniformly_refine (n_refinements);
221
222 // Print information about the mesh to the screen.
223 mesh.print_info();
224
225 // Declare the system and its variables.
226 // Begin by creating a transient system
227 // named "Convection-Diffusion".
228 TransientLinearImplicitSystem & system =
229 equation_systems.add_system<TransientLinearImplicitSystem>("Convection-Diffusion");
230
231 // Adds the variable "u" to "Convection-Diffusion". "u"
232 // will be approximated using first-order approximation.
233 system.add_variable ("u", FIRST);
234
235 // Give the system a pointer to the matrix assembly
236 // and initialization functions.
237 system.attach_assemble_function (assemble_cd);
238 system.attach_init_function (init_cd);
239
240 // Initialize the data structures for the equation system.
241 equation_systems.init ();
242 }
243 // Otherwise we read in the solution and mesh
244 else
245 {
246 // Read in the mesh stored in "saved_mesh.xda"
247 mesh.read("saved_mesh.xda");
248
249 // Print information about the mesh to the screen.
250 mesh.print_info();
251
252 // Read in the solution stored in "saved_solution.xda"
253 equation_systems.read("saved_solution.xda", READ);
254
255 // Get a reference to the system so that we can call update() on it
256 TransientLinearImplicitSystem & system =
257 equation_systems.get_system<TransientLinearImplicitSystem>("Convection-Diffusion");
258
259 // We need to call update to put system in a consistent state
260 // with the solution that was read in
261 system.update();
262
263 // Attach the same matrix assembly function as above. Note, we do not
264 // have to attach an init() function since we are initializing the
265 // system by reading in "saved_solution.xda"
266 system.attach_assemble_function (assemble_cd);
267
268 // Print out the H1 norm of the saved solution, for verification purposes:
269 Real H1norm = system.calculate_norm(*system.solution, SystemNorm(H1));
270
271 libMesh::out << "Initial H1 norm = " << H1norm << std::endl << std::endl;
272 }
273
274 // Prints information about the system to the screen.
275 equation_systems.print_info();
276
277 equation_systems.parameters.set<unsigned int>("linear solver maximum iterations") = 250;
278 equation_systems.parameters.set<Real>("linear solver tolerance") = TOLERANCE;
279
280 if (!read_solution)
281 // Write out the initial condition
282 GMVIO(mesh).write_equation_systems ("out.gmv.000", equation_systems);
283 else
284 // Write out the solution that was read in
285 GMVIO(mesh).write_equation_systems ("solution_read_in.gmv", equation_systems);
286
287 // The Convection-Diffusion system requires that we specify
288 // the flow velocity. We will specify it as a RealVectorValue
289 // data type and then use the Parameters object to pass it to
290 // the assemble function.
291 equation_systems.parameters.set<RealVectorValue>("velocity") =
292 RealVectorValue (0.8, 0.8);
293
294 // The Convection-Diffusion system also requires a specified
295 // diffusivity. We use an isotropic (hence Real) value.
296 equation_systems.parameters.set<Real>("diffusivity") = 0.01;
297
298 // Solve the system "Convection-Diffusion". This will be done by
299 // looping over the specified time interval and calling the
300 // solve() member at each time step. This will assemble the
301 // system and call the linear solver.
302
303 // Since only TransientLinearImplicitSystems (and systems
304 // derived from them) contain old solutions, to use the
305 // old_local_solution later we now need to specify the system
306 // type when we ask for it.
307 TransientLinearImplicitSystem & system =
308 equation_systems.get_system<TransientLinearImplicitSystem>("Convection-Diffusion");
309
310 const Real dt = 0.025;
311 system.time = init_timestep*dt;
312
313 // We do 25 timesteps both before and after writing out the
314 // intermediate solution
315 for (unsigned int t_step=init_timestep; t_step<(init_timestep+n_timesteps); t_step++)
316 {
317 // Increment the time counter, set the time and the
318 // time step size as parameters in the EquationSystem.
319 system.time += dt;
320
321 equation_systems.parameters.set<Real> ("time") = system.time;
322 equation_systems.parameters.set<Real> ("dt") = dt;
323
324 // A pretty update message
325 libMesh::out << " Solving time step ";
326
327 // Add a set of scope braces to enforce data locality.
328 {
329 std::ostringstream out;
330
331 out << std::setw(2)
332 << std::right
333 << t_step
334 << ", time="
335 << std::fixed
336 << std::setw(6)
337 << std::setprecision(3)
338 << std::setfill('0')
339 << std::left
340 << system.time
341 << "...";
342
343 libMesh::out << out.str() << std::endl;
344 }
345
346 // At this point we need to update the old
347 // solution vector. The old solution vector
348 // will be the current solution vector from the
349 // previous time step.
350
351 *system.old_local_solution = *system.current_local_solution;
352
353 // The number of refinement steps per time step.
354 const unsigned int max_r_steps = 2;
355
356 // A refinement loop.
357 for (unsigned int r_step=0; r_step<max_r_steps; r_step++)
358 {
359 // Assemble & solve the linear system
360 system.solve();
361
362 // Print out the H1 norm, for verification purposes:
363 Real H1norm = system.calculate_norm(*system.solution, SystemNorm(H1));
364
365 libMesh::out << "H1 norm = " << H1norm << std::endl;
366
367 // Possibly refine the mesh
368 if (r_step+1 != max_r_steps)
369 {
370 libMesh::out << " Refining the mesh..." << std::endl;
371
372 // The ErrorVector is a particular StatisticsVector
373 // for computing error information on a finite element mesh.
374 ErrorVector error;
375
376 // The ErrorEstimator class interrogates a finite element
377 // solution and assigns to each element a positive error value.
378 // This value is used for deciding which elements to refine
379 // and which to coarsen.
380 KellyErrorEstimator error_estimator;
381
382 // This is a subclass of JumpErrorEstimator, based on
383 // measuring discontinuities across sides between
384 // elements, and we can tell it to use a cheaper
385 // "unweighted" quadrature rule when numerically
386 // integrating those discontinuities.
387 error_estimator.use_unweighted_quadrature_rules = true;
388
389 // Compute the error for each active element using the provided
390 // flux_jump indicator. Note in general you will need to
391 // provide an error estimator specifically designed for your
392 // application.
393 error_estimator.estimate_error (system,
394 error);
395
396 // This takes the error in error and decides which elements
397 // will be coarsened or refined. Any element within 20% of the
398 // maximum error on any element will be refined, and any
399 // element within 7% of the minimum error on any element might
400 // be coarsened. Note that the elements flagged for refinement
401 // will be refined, but those flagged for coarsening _might_ be
402 // coarsened.
403 mesh_refinement.refine_fraction() = 0.80;
404 mesh_refinement.coarsen_fraction() = 0.07;
405 mesh_refinement.max_h_level() = 5;
406 mesh_refinement.flag_elements_by_error_fraction (error);
407
408 // This call actually refines and coarsens the flagged
409 // elements.
410 mesh_refinement.refine_and_coarsen_elements();
411
412 // This call reinitializes the EquationSystems object for
413 // the newly refined mesh. One of the steps in the
414 // reinitialization is projecting the solution,
415 // old_solution, etc... vectors from the old mesh to
416 // the current one.
417 equation_systems.reinit ();
418 }
419 }
420
421 // Again do a search on the command line for an argument
422 unsigned int output_freq = 10;
423 if (command_line.search("-output_freq"))
424 output_freq = command_line.next(0);
425
426 // Output every 10 timesteps to file.
427 if ((t_step+1)%output_freq == 0)
428 {
429 std::ostringstream file_name;
430
431 file_name << "out.gmv."
432 << std::setw(3)
433 << std::setfill('0')
434 << std::right
435 << t_step+1;
436
437 GMVIO(mesh).write_equation_systems (file_name.str(),
438 equation_systems);
439 }
440 }
441
442 if (!read_solution)
443 {
444 // Print out the H1 norm of the saved solution, for verification purposes:
445 Real H1norm = system.calculate_norm(*system.solution, SystemNorm(H1));
446
447 libMesh::out << "Final H1 norm = " << H1norm << std::endl << std::endl;
448
449 mesh.write("saved_mesh.xda");
450 equation_systems.write("saved_solution.xda", WRITE);
451 GMVIO(mesh).write_equation_systems ("saved_solution.gmv",
452 equation_systems);
453 }
454 #endif // #ifndef LIBMESH_ENABLE_AMR
455
456 return 0;
457 }
458
459 // Here we define the initialization routine for the
460 // Convection-Diffusion system. This routine is
461 // responsible for applying the initial conditions to
462 // the system.
init_cd(EquationSystems & es,const std::string & libmesh_dbg_var (system_name))463 void init_cd (EquationSystems & es,
464 const std::string & libmesh_dbg_var(system_name))
465 {
466 // It is a good idea to make sure we are initializing
467 // the proper system.
468 libmesh_assert_equal_to (system_name, "Convection-Diffusion");
469
470 // Get a reference to the Convection-Diffusion system object.
471 TransientLinearImplicitSystem & system =
472 es.get_system<TransientLinearImplicitSystem>("Convection-Diffusion");
473
474 // Project initial conditions at time 0
475 es.parameters.set<Real> ("time") = system.time = 0;
476
477 system.project_solution(exact_value, nullptr, es.parameters);
478 }
479
480
481
482 // This function defines the assembly routine which
483 // will be called at each time step. It is responsible
484 // for computing the proper matrix entries for the
485 // element stiffness matrices and right-hand sides.
486 #ifdef LIBMESH_ENABLE_AMR
assemble_cd(EquationSystems & es,const std::string & libmesh_dbg_var (system_name))487 void assemble_cd (EquationSystems & es,
488 const std::string & libmesh_dbg_var(system_name))
489 {
490 // It is a good idea to make sure we are assembling
491 // the proper system.
492 libmesh_assert_equal_to (system_name, "Convection-Diffusion");
493
494 // Get a constant reference to the mesh object.
495 const MeshBase & mesh = es.get_mesh();
496
497 // The dimension that we are running
498 const unsigned int dim = mesh.mesh_dimension();
499
500 // Get a reference to the Convection-Diffusion system object.
501 TransientLinearImplicitSystem & system =
502 es.get_system<TransientLinearImplicitSystem> ("Convection-Diffusion");
503
504 // Get the Finite Element type for the first (and only)
505 // variable in the system.
506 FEType fe_type = system.variable_type(0);
507
508 // Build a Finite Element object of the specified type. Since the
509 // FEBase::build() member dynamically creates memory we will
510 // store the object as a std::unique_ptr<FEBase>. This can be thought
511 // of as a pointer that will clean up after itself.
512 std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
513 std::unique_ptr<FEBase> fe_face (FEBase::build(dim, fe_type));
514
515 // A Gauss quadrature rule for numerical integration.
516 // Let the FEType object decide what order rule is appropriate.
517 QGauss qrule (dim, fe_type.default_quadrature_order());
518 QGauss qface (dim-1, fe_type.default_quadrature_order());
519
520 // Tell the finite element object to use our quadrature rule.
521 fe->attach_quadrature_rule (&qrule);
522 fe_face->attach_quadrature_rule (&qface);
523
524 // Here we define some references to cell-specific data that
525 // will be used to assemble the linear system. We will start
526 // with the element Jacobian * quadrature weight at each integration point.
527 const std::vector<Real> & JxW = fe->get_JxW();
528 const std::vector<Real> & JxW_face = fe_face->get_JxW();
529
530 // The element shape functions evaluated at the quadrature points.
531 const std::vector<std::vector<Real>> & phi = fe->get_phi();
532 const std::vector<std::vector<Real>> & psi = fe_face->get_phi();
533
534 // The element shape function gradients evaluated at the quadrature
535 // points.
536 const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
537
538 // The XY locations of the quadrature points used for face integration
539 const std::vector<Point> & qface_points = fe_face->get_xyz();
540
541 // A reference to the DofMap object for this system. The DofMap
542 // object handles the index translation from node and element numbers
543 // to degree of freedom numbers. We will talk more about the DofMap
544 // in future examples.
545 const DofMap & dof_map = system.get_dof_map();
546
547 // Define data structures to contain the element matrix
548 // and right-hand-side vector contribution. Following
549 // basic finite element terminology we will denote these
550 // "Ke" and "Fe".
551 DenseMatrix<Number> Ke;
552 DenseVector<Number> Fe;
553
554 // This vector will hold the degree of freedom indices for
555 // the element. These define where in the global system
556 // the element degrees of freedom get mapped.
557 std::vector<dof_id_type> dof_indices;
558
559 // Here we extract the velocity & parameters that we put in the
560 // EquationSystems object.
561 const RealVectorValue velocity =
562 es.parameters.get<RealVectorValue> ("velocity");
563
564 const Real diffusivity =
565 es.parameters.get<Real> ("diffusivity");
566
567 const Real dt = es.parameters.get<Real> ("dt");
568
569 // Now we will loop over all the elements in the mesh that
570 // live on the local processor. We will compute the element
571 // matrix and right-hand-side contribution. Since the mesh
572 // will be refined we want to only consider the ACTIVE elements,
573 // hence we use a variant of the active_elem_iterator.
574 for (const auto & elem : mesh.active_local_element_ptr_range())
575 {
576 // Get the degree of freedom indices for the
577 // current element. These define where in the global
578 // matrix and right-hand-side this element will
579 // contribute to.
580 dof_map.dof_indices (elem, dof_indices);
581
582 // Compute the element-specific data for the current
583 // element. This involves computing the location of the
584 // quadrature points (q_point) and the shape functions
585 // (phi, dphi) for the current element.
586 fe->reinit (elem);
587
588 const unsigned int n_dofs =
589 cast_int<unsigned int>(dof_indices.size());
590 libmesh_assert_equal_to (n_dofs, phi.size());
591
592 // Zero the element matrix and right-hand side before
593 // summing them. We use the resize member here because
594 // the number of degrees of freedom might have changed from
595 // the last element. Note that this will be the case if the
596 // element type is different (i.e. the last element was a
597 // triangle, now we are on a quadrilateral).
598 Ke.resize (n_dofs, n_dofs);
599
600 Fe.resize (n_dofs);
601
602 // Now we will build the element matrix and right-hand-side.
603 // Constructing the RHS requires the solution and its
604 // gradient from the previous timestep. This myst be
605 // calculated at each quadrature point by summing the
606 // solution degree-of-freedom values by the appropriate
607 // weight functions.
608 for (unsigned int qp=0; qp<qrule.n_points(); qp++)
609 {
610 // Values to hold the old solution & its gradient.
611 Number u_old = 0.;
612 Gradient grad_u_old;
613
614 // Compute the old solution & its gradient.
615 for (unsigned int l=0; l != n_dofs; l++)
616 {
617 u_old += phi[l][qp]*system.old_solution (dof_indices[l]);
618
619 // This will work,
620 // grad_u_old += dphi[l][qp]*system.old_solution (dof_indices[l]);
621 // but we can do it without creating a temporary like this:
622 grad_u_old.add_scaled (dphi[l][qp], system.old_solution (dof_indices[l]));
623 }
624
625 // Now compute the element matrix and RHS contributions.
626 for (unsigned int i=0; i != n_dofs; i++)
627 {
628 // The RHS contribution
629 Fe(i) += JxW[qp]*(
630 // Mass matrix term
631 u_old*phi[i][qp] +
632 -.5*dt*(
633 // Convection term
634 // (grad_u_old may be complex, so the
635 // order here is important!)
636 (grad_u_old*velocity)*phi[i][qp] +
637
638 // Diffusion term
639 diffusivity*(grad_u_old*dphi[i][qp]))
640 );
641
642 for (unsigned int j=0; j != n_dofs; j++)
643 {
644 // The matrix contribution
645 Ke(i,j) += JxW[qp]*(
646 // Mass-matrix
647 phi[i][qp]*phi[j][qp] +
648 .5*dt*(
649 // Convection term
650 (velocity*dphi[j][qp])*phi[i][qp] +
651 // Diffusion term
652 diffusivity*(dphi[i][qp]*dphi[j][qp]))
653 );
654 }
655 }
656 }
657
658 // At this point the interior element integration has
659 // been completed. However, we have not yet addressed
660 // boundary conditions. For this example we will only
661 // consider simple Dirichlet boundary conditions imposed
662 // via the penalty method.
663 //
664 // The following loops over the sides of the element.
665 // If the element has no neighbor on a side then that
666 // side MUST live on a boundary of the domain.
667 {
668 // The penalty value.
669 const Real penalty = 1.e10;
670
671 // The following loops over the sides of the element.
672 // If the element has no neighbor on a side then that
673 // side MUST live on a boundary of the domain.
674 for (auto s : elem->side_index_range())
675 if (elem->neighbor_ptr(s) == nullptr)
676 {
677 fe_face->reinit(elem, s);
678
679 libmesh_assert_equal_to (n_dofs, psi.size());
680
681 for (unsigned int qp=0; qp<qface.n_points(); qp++)
682 {
683 const Number value = exact_solution (qface_points[qp](0),
684 qface_points[qp](1),
685 system.time);
686
687 // RHS contribution
688 for (unsigned int i=0; i != n_dofs; i++)
689 Fe(i) += penalty*JxW_face[qp]*value*psi[i][qp];
690
691 // Matrix contribution
692 for (unsigned int i=0; i != n_dofs; i++)
693 for (unsigned int j=0; j != n_dofs; j++)
694 Ke(i,j) += penalty*JxW_face[qp]*psi[i][qp]*psi[j][qp];
695 }
696 }
697 }
698
699
700 // We have now built the element matrix and RHS vector in terms
701 // of the element degrees of freedom. However, it is possible
702 // that some of the element DOFs are constrained to enforce
703 // solution continuity, i.e. they are not really "free". We need
704 // to constrain those DOFs in terms of non-constrained DOFs to
705 // ensure a continuous solution. The
706 // DofMap::constrain_element_matrix_and_vector() method does
707 // just that.
708 dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
709
710 // The element matrix and right-hand-side are now built
711 // for this element. Add them to the global matrix and
712 // right-hand-side vector. The SparseMatrix::add_matrix()
713 // and NumericVector::add_vector() members do this for us.
714 system.matrix->add_matrix (Ke, dof_indices);
715 system.rhs->add_vector (Fe, dof_indices);
716
717 }
718 // Finished computing the system matrix and right-hand side.
719 }
720 #endif // #ifdef LIBMESH_ENABLE_AMR
721