libMesh
explicit_system.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 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 // C++ includes
21 
22 // Local includes
23 #include "libmesh/explicit_system.h"
24 #include "libmesh/numeric_vector.h"
25 
26 namespace libMesh
27 {
28 
29 
30 // ------------------------------------------------------------
31 // ExplicitSystem implementation
33  const std::string & name_in,
34  const unsigned int number_in) :
35  Parent (es, name_in, number_in),
36  rhs(nullptr)
37 
38 {
39  // Add the system RHS.
40  this->add_system_rhs ();
41 }
42 
43 
44 
46 {
47  // Clear the parent data
48  Parent::clear();
49 
50  // Restore us to a "basic" state
51  this->add_system_rhs ();
52 }
53 
54 
55 
56 void ExplicitSystem::assemble_qoi (const QoISet & qoi_indices)
57 {
58  // The user quantity of interest assembly gets to expect to
59  // accumulate on initially zero values
60  for (auto i : IntRange<unsigned int>(0, this->n_qois()))
61  if (qoi_indices.has_index(i))
62  qoi[i] = 0;
63 
64  Parent::assemble_qoi (qoi_indices);
65 }
66 
67 
68 
70  bool include_liftfunc,
71  bool apply_constraints)
72 {
73  // The user quantity of interest derivative assembly gets to expect
74  // to accumulate on initially zero vectors
75  for (auto i : IntRange<unsigned int>(0, this->n_qois()))
76  if (qoi_indices.has_index(i))
77  this->add_adjoint_rhs(i).zero();
78 
79  Parent::assemble_qoi_derivative (qoi_indices, include_liftfunc,
80  apply_constraints);
81 }
82 
83 
84 
86 {
87  // Assemble the linear system
88  this->assemble ();
89 
90  // Update the system after the solve
91  this->update();
92 }
93 
94 
95 
97 {
98  // Possible that we cleared the _vectors but
99  // forgot to update the rhs pointer?
100  if (this->n_vectors() == 0)
101  rhs = nullptr;
102 
103 
104  // Only need to add the rhs if it isn't there
105  // already!
106  if (rhs == nullptr)
107  rhs = &(this->add_vector ("RHS Vector", false));
108 
110 }
111 
112 } // namespace libMesh
libMesh::System
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:100
libMesh::NumericVector::zero
virtual void zero()=0
Set all entries to zero.
libMesh::ExplicitSystem::solve
virtual void solve() override
For explicit systems, just assemble the system which should directly compute A*x.
Definition: explicit_system.C:85
libMesh::ExplicitSystem::rhs
NumericVector< Number > * rhs
The system matrix.
Definition: explicit_system.h:114
libMesh::QoISet::has_index
bool has_index(std::size_t) const
Return whether or not this index is in the set to be calculated.
Definition: qoi_set.h:221
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::System::n_vectors
unsigned int n_vectors() const
Definition: system.h:2283
libMesh::ExplicitSystem::add_system_rhs
void add_system_rhs()
Add the system right-hand-side vector to the _vectors data structure.
Definition: explicit_system.C:96
libMesh::System::assemble
virtual void assemble()
Prepares matrix and _dof_map for matrix assembly.
Definition: system.C:462
libMesh::ExplicitSystem::ExplicitSystem
ExplicitSystem(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
Definition: explicit_system.C:32
libMesh::System::n_qois
unsigned int n_qois() const
Number of currently active quantities of interest.
Definition: system.h:2328
libMesh::System::add_adjoint_rhs
NumericVector< Number > & add_adjoint_rhs(unsigned int i=0)
Definition: system.C:1009
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::IntRange
The IntRange templated class is intended to make it easy to loop over integers which are indices of a...
Definition: int_range.h:53
libMesh::System::qoi
std::vector< Number > qoi
Values of the quantities of interest.
Definition: system.h:1574
libMesh::QoISet
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
Definition: qoi_set.h:45
libMesh::ExplicitSystem::assemble_qoi
virtual void assemble_qoi(const QoISet &qoi_indices=QoISet()) override
Prepares qoi for quantity of interest assembly, then calls user qoi function.
Definition: explicit_system.C:56
libMesh::System::assemble_qoi
virtual void assemble_qoi(const QoISet &qoi_indices=QoISet())
Calls user qoi function.
Definition: system.C:473
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::System::assemble_qoi_derivative
virtual void assemble_qoi_derivative(const QoISet &qoi_indices=QoISet(), bool include_liftfunc=true, bool apply_constraints=true)
Calls user qoi derivative function.
Definition: system.C:484
libMesh::ExplicitSystem::clear
virtual void clear() override
Clear all the data structures associated with the system.
Definition: explicit_system.C:45
libMesh::System::clear
virtual void clear()
Clear all the data structures associated with the system.
Definition: system.C:205
libMesh::System::add_vector
NumericVector< Number > & add_vector(const std::string &vec_name, const bool projections=true, const ParallelType type=PARALLEL)
Adds the additional vector vec_name to this system.
Definition: system.C:661
libMesh::System::update
virtual void update()
Update the local values to reflect the solution on neighboring processors.
Definition: system.C:408
libMesh::ExplicitSystem::assemble_qoi_derivative
virtual void assemble_qoi_derivative(const QoISet &qoi_indices=QoISet(), bool include_liftfunc=true, bool apply_constraints=true) override
Prepares adjoint_rhs for quantity of interest derivative assembly, then calls user qoi derivative fun...
Definition: explicit_system.C:69