MFEM  v4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
linearform.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 // Implementation of class LinearForm
13 
14 #include "fem.hpp"
15 
16 namespace mfem
17 {
18 
20  : Vector(f->GetVSize())
21 {
22  // Linear forms are stored on the device
23  UseDevice(true);
24 
25  fes = f;
26  extern_lfs = 1;
27 
28  // Copy the pointers to the integrators
29  dlfi = lf->dlfi;
30 
31  dlfi_delta = lf->dlfi_delta;
32 
33  blfi = lf->blfi;
34 
35  flfi = lf->flfi;
37 }
38 
40 {
41  DeltaLFIntegrator *maybe_delta =
42  dynamic_cast<DeltaLFIntegrator *>(lfi);
43  if (!maybe_delta || !maybe_delta->IsDelta())
44  {
45  dlfi.Append(lfi);
46  }
47  else
48  {
49  dlfi_delta.Append(maybe_delta);
50  }
51 }
52 
54 {
55  blfi.Append (lfi);
56  blfi_marker.Append(NULL); // NULL -> all attributes are active
57 }
58 
60  Array<int> &bdr_attr_marker)
61 {
62  blfi.Append (lfi);
63  blfi_marker.Append(&bdr_attr_marker);
64 }
65 
67 {
68  flfi.Append(lfi);
69  flfi_marker.Append(NULL); // NULL -> all attributes are active
70 }
71 
73  Array<int> &bdr_attr_marker)
74 {
75  flfi.Append(lfi);
76  flfi_marker.Append(&bdr_attr_marker);
77 }
78 
80 {
81  Array<int> vdofs;
82  ElementTransformation *eltrans;
83  Vector elemvect;
84 
85  int i;
86 
87  Vector::operator=(0.0);
88 
89  // The above operation is executed on device because of UseDevice().
90  // The first use of AddElementVector() below will move it back to host
91  // because both 'vdofs' and 'elemvect' are on host.
92 
93  if (dlfi.Size())
94  {
95  for (i = 0; i < fes -> GetNE(); i++)
96  {
97  fes -> GetElementVDofs (i, vdofs);
98  eltrans = fes -> GetElementTransformation (i);
99  for (int k=0; k < dlfi.Size(); k++)
100  {
101  dlfi[k]->AssembleRHSElementVect(*fes->GetFE(i), *eltrans, elemvect);
102  AddElementVector (vdofs, elemvect);
103  }
104  }
105  }
106  AssembleDelta();
107 
108  if (blfi.Size())
109  {
110  Mesh *mesh = fes->GetMesh();
111 
112  // Which boundary attributes need to be processed?
113  Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
114  mesh->bdr_attributes.Max() : 0);
115  bdr_attr_marker = 0;
116  for (int k = 0; k < blfi.Size(); k++)
117  {
118  if (blfi_marker[k] == NULL)
119  {
120  bdr_attr_marker = 1;
121  break;
122  }
123  Array<int> &bdr_marker = *blfi_marker[k];
124  MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
125  "invalid boundary marker for boundary integrator #"
126  << k << ", counting from zero");
127  for (int i = 0; i < bdr_attr_marker.Size(); i++)
128  {
129  bdr_attr_marker[i] |= bdr_marker[i];
130  }
131  }
132 
133  for (i = 0; i < fes -> GetNBE(); i++)
134  {
135  const int bdr_attr = mesh->GetBdrAttribute(i);
136  if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
137  fes -> GetBdrElementVDofs (i, vdofs);
138  eltrans = fes -> GetBdrElementTransformation (i);
139  for (int k=0; k < blfi.Size(); k++)
140  {
141  blfi[k]->AssembleRHSElementVect(*fes->GetBE(i), *eltrans, elemvect);
142  AddElementVector (vdofs, elemvect);
143  }
144  }
145  }
146  if (flfi.Size())
147  {
149  Mesh *mesh = fes->GetMesh();
150 
151  // Which boundary attributes need to be processed?
152  Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
153  mesh->bdr_attributes.Max() : 0);
154  bdr_attr_marker = 0;
155  for (int k = 0; k < flfi.Size(); k++)
156  {
157  if (flfi_marker[k] == NULL)
158  {
159  bdr_attr_marker = 1;
160  break;
161  }
162  Array<int> &bdr_marker = *flfi_marker[k];
163  MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
164  "invalid boundary marker for boundary face integrator #"
165  << k << ", counting from zero");
166  for (int i = 0; i < bdr_attr_marker.Size(); i++)
167  {
168  bdr_attr_marker[i] |= bdr_marker[i];
169  }
170  }
171 
172  for (i = 0; i < mesh->GetNBE(); i++)
173  {
174  const int bdr_attr = mesh->GetBdrAttribute(i);
175  if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
176 
177  tr = mesh->GetBdrFaceTransformations(i);
178  if (tr != NULL)
179  {
180  fes -> GetElementVDofs (tr -> Elem1No, vdofs);
181  for (int k = 0; k < flfi.Size(); k++)
182  {
183  if (flfi_marker[k] &&
184  (*flfi_marker[k])[bdr_attr-1] == 0) { continue; }
185 
186  flfi[k] -> AssembleRHSElementVect (*fes->GetFE(tr -> Elem1No),
187  *tr, elemvect);
188  AddElementVector (vdofs, elemvect);
189  }
190  }
191  }
192  }
193 }
194 
195 void LinearForm::Update(FiniteElementSpace *f, Vector &v, int v_offset)
196 {
197  fes = f;
198  NewDataAndSize((double *)v + v_offset, fes->GetVSize());
200 }
201 
203 {
204  if (dlfi_delta.Size() == 0) { return; }
205 
206  if (!HaveDeltaLocations())
207  {
208  int sdim = fes->GetMesh()->SpaceDimension();
209  Vector center;
210  DenseMatrix centers(sdim, dlfi_delta.Size());
211  for (int i = 0; i < centers.Width(); i++)
212  {
213  centers.GetColumnReference(i, center);
214  dlfi_delta[i]->GetDeltaCenter(center);
215  MFEM_VERIFY(center.Size() == sdim,
216  "Point dim " << center.Size() <<
217  " does not match space dim " << sdim);
218  }
220  }
221 
222  Array<int> vdofs;
223  Vector elemvect;
224  for (int i = 0; i < dlfi_delta.Size(); i++)
225  {
226  int elem_id = dlfi_delta_elem_id[i];
227  // The delta center may be outside of this sub-domain, or
228  // (Par)Mesh::FindPoints() failed to find this point:
229  if (elem_id < 0) { continue; }
230 
231  const IntegrationPoint &ip = dlfi_delta_ip[i];
233  Trans.SetIntPoint(&ip);
234 
235  fes->GetElementVDofs(elem_id, vdofs);
236  dlfi_delta[i]->AssembleDeltaElementVect(*fes->GetFE(elem_id), Trans,
237  elemvect);
238  AddElementVector(vdofs, elemvect);
239  }
240 }
241 
243 {
244  Vector::operator=(value);
245  return *this;
246 }
247 
249 {
250  MFEM_ASSERT(fes && v.Size() == fes->GetVSize(), "");
252  return *this;
253 }
254 
256 {
257  if (!extern_lfs)
258  {
259  int k;
260  for (k=0; k < dlfi_delta.Size(); k++) { delete dlfi_delta[k]; }
261  for (k=0; k < dlfi.Size(); k++) { delete dlfi[k]; }
262  for (k=0; k < blfi.Size(); k++) { delete blfi[k]; }
263  for (k=0; k < flfi.Size(); k++) { delete flfi[k]; }
264  }
265 }
266 
267 
268 }
Array< LinearFormIntegrator * > blfi
Set of Boundary Integrators to be applied.
Definition: linearform.hpp:40
int Size() const
Logical size of the array.
Definition: array.hpp:118
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:347
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
Definition: mesh.hpp:976
void NewDataAndSize(double *d, int s)
Set the Vector data and size, deleting the old data, if owned.
Definition: vector.hpp:131
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:679
void GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:172
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:79
int extern_lfs
Indicates the LinerFormIntegrators stored in dlfi, dlfi_delta, blfi, and flfi are owned by another Li...
Definition: linearform.hpp:31
LinearForm & operator=(const LinearForm &rhs)
Copy assignment. Only the data of the base class Vector is copied.
Definition: linearform.hpp:90
void SetIntPoint(const IntegrationPoint *ip)
Definition: eltrans.hpp:53
Array< Array< int > * > blfi_marker
Entries are not owned.
Definition: linearform.hpp:41
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:150
void ResetDeltaLocations()
Force (re)computation of delta locations.
Definition: linearform.hpp:57
Array< DeltaLFIntegrator * > dlfi_delta
Separate array for integrators with delta function coefficients.
Definition: linearform.hpp:37
Array< IntegrationPoint > dlfi_delta_ip
The reference coordinates where the centers of the delta functions lie.
Definition: linearform.hpp:51
Abstract base class LinearFormIntegrator.
Definition: lininteg.hpp:22
Vector & operator=(const double *v)
Copy Size() entries from v.
Definition: vector.cpp:106
bool UseDevice() const
Return the device flag of the Memory object used by the Vector.
Definition: vector.hpp:89
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:272
Array< LinearFormIntegrator * > dlfi
Set of Domain Integrators to be applied.
Definition: linearform.hpp:34
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
FiniteElementSpace * fes
FE space on which the LinearForm lives. Not owned.
Definition: linearform.hpp:27
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Face Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:66
FaceElementTransformations * GetBdrFaceTransformations(int BdrElemNo)
Definition: mesh.cpp:915
void AddBoundaryIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:53
bool IsDelta() const
Returns true if the derived class instance uses a delta coefficient.
Definition: lininteg.hpp:70
void AddElementVector(const Array< int > &dofs, const Vector &elemvect)
Add (element) subvector to the vector.
Definition: vector.cpp:564
int SpaceDimension() const
Definition: mesh.hpp:714
LinearForm()
Create an empty LinearForm without an associated FiniteElementSpace.
Definition: linearform.hpp:82
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition: fespace.hpp:398
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:179
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:85
void GetColumnReference(int c, Vector &col)
Definition: densemat.hpp:262
Abstract class for integrators that support delta coefficients.
Definition: lininteg.hpp:47
void AssembleDelta()
Assembles delta functions of the linear form.
Definition: linearform.cpp:202
void Update()
Update the object according to the associated FE space fes.
Definition: linearform.hpp:158
Class for integration point with weight.
Definition: intrules.hpp:25
~LinearForm()
Destroys linear form.
Definition: linearform.cpp:255
Array< int > dlfi_delta_elem_id
The element ids where the centers of the delta functions lie.
Definition: linearform.hpp:48
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement associated with i&#39;th element.
Definition: fespace.cpp:1541
bool HaveDeltaLocations()
If true, the delta locations are not (re)computed during assembly.
Definition: linearform.hpp:54
virtual int FindPoints(DenseMatrix &point_mat, Array< int > &elem_ids, Array< IntegrationPoint > &ips, bool warn=true, InverseElementTransformation *inv_trans=NULL)
Find the ids of the elements that contain the given points, and their corresponding reference coordin...
Definition: mesh.cpp:9449
Vector data type.
Definition: vector.hpp:48
Class for linear form - Vector with associated FE space and LFIntegrators.
Definition: linearform.hpp:23
Array< LinearFormIntegrator * > flfi
Set of Boundary Face Integrators to be applied.
Definition: linearform.hpp:44
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement for the i&#39;th boundary element.
Definition: fespace.cpp:1781
Array< Array< int > * > flfi_marker
Entries are not owned.
Definition: linearform.hpp:45