MFEM  v3.2
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
coefficient.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 Coefficient class
13 
14 #include "fem.hpp"
15 
16 #include <cmath>
17 #include <limits>
18 
19 namespace mfem
20 {
21 
22 using namespace std;
23 
25  const IntegrationPoint & ip)
26 {
27  int att = T.Attribute;
28  return (constants(att-1));
29 }
30 
32  const IntegrationPoint & ip)
33 {
34  double x[3];
35  Vector transip(x, 3);
36 
37  T.Transform(ip, transip);
38 
39  if (Function)
40  {
41  return ((*Function)(transip));
42  }
43  else
44  {
45  return (*TDFunction)(transip, GetTime());
46  }
47 }
48 
50  const IntegrationPoint &ip)
51 {
52  return GridF -> GetValue (T.ElementNo, ip, Component);
53 }
54 
56  const IntegrationPoint &ip)
57 {
58  if (Q2)
59  {
60  return (*Transform2)(Q1->Eval(T, ip, GetTime()),
61  Q2->Eval(T, ip, GetTime()));
62  }
63  else
64  {
65  return (*Transform1)(Q1->Eval(T, ip, GetTime()));
66  }
67 }
68 
70  const IntegrationRule &ir)
71 {
72  Vector Mi;
73  M.SetSize(vdim, ir.GetNPoints());
74  for (int i = 0; i < ir.GetNPoints(); i++)
75  {
76  M.GetColumnReference(i, Mi);
77  const IntegrationPoint &ip = ir.IntPoint(i);
78  T.SetIntPoint(&ip);
79  Eval(Mi, T, ip);
80  }
81 }
82 
84  const IntegrationPoint &ip)
85 {
86  double x[3];
87  Vector transip(x, 3);
88 
89  T.Transform(ip, transip);
90 
91  V.SetSize(vdim);
92  if (Function)
93  {
94  (*Function)(transip, V);
95  }
96  else
97  {
98  (*TDFunction)(transip, GetTime(), V);
99  }
100  if (Q)
101  {
102  V *= Q->Eval(T, ip, GetTime());
103  }
104 }
105 
107  : VectorCoefficient(dim), Coeff(dim)
108 {
109  for (int i = 0; i < dim; i++)
110  {
111  Coeff[i] = NULL;
112  }
113 }
114 
116 {
117  for (int i = 0; i < vdim; i++)
118  {
119  delete Coeff[i];
120  }
121 }
122 
124  const IntegrationPoint &ip)
125 {
126  V.SetSize(vdim);
127  for (int i = 0; i < vdim; i++)
128  {
129  V(i) = Coeff[i]->Eval(T, ip, GetTime());
130  }
131 }
132 
134  GridFunction *gf) : VectorCoefficient (gf -> VectorDim())
135 {
136  GridFunc = gf;
137 }
138 
140  const IntegrationPoint &ip)
141 {
142  GridFunc->GetVectorValue(T.ElementNo, ip, V);
143 }
144 
147 {
148  GridFunc->GetVectorValues(T, ir, M);
149 }
150 
152  const IntegrationPoint &ip)
153 {
154  V.SetSize(vdim);
155  if (active_attr[T.Attribute-1])
156  {
157  c->SetTime(GetTime());
158  c->Eval(V, T, ip);
159  }
160  else
161  {
162  V = 0.0;
163  }
164 }
165 
168 {
169  if (active_attr[T.Attribute-1])
170  {
171  c->SetTime(GetTime());
172  c->Eval(M, T, ir);
173  }
174  else
175  {
176  M.SetSize(vdim);
177  M = 0.0;
178  }
179 }
180 
182  const IntegrationPoint &ip)
183 {
184  double x[3];
185  Vector transip(x, 3);
186 
187  T.Transform(ip, transip);
188 
189  K.SetSize(vdim);
190 
191  if (Function)
192  {
193  (*Function)(transip, K);
194  }
195  else
196  {
197  (*TDFunction)(transip, GetTime(), K);
198  }
199 }
200 
202  : MatrixCoefficient (dim)
203 {
204  Coeff.SetSize (vdim*vdim);
205 }
206 
208 {
209  for (int i=0; i< vdim*vdim; i++)
210  {
211  delete Coeff[i];
212  }
213 }
214 
216  const IntegrationPoint &ip)
217 {
218  int i, j;
219 
220  for (i = 0; i < vdim; i++)
221  for (j = 0; j < vdim; j++)
222  {
223  K(i,j) = Coeff[i*vdim+j] -> Eval(T, ip, GetTime());
224  }
225 }
226 
228  const IntegrationPoint &ip)
229 {
230  if (active_attr[T.Attribute-1])
231  {
232  c->SetTime(GetTime());
233  c->Eval(K, T, ip);
234  }
235  else
236  {
237  K.SetSize(vdim);
238  K = 0.0;
239  }
240 }
241 
242 double LpNormLoop(double p, Coefficient &coeff, Mesh &mesh,
243  const IntegrationRule *irs[])
244 {
245  double norm = 0.0;
247 
248  for (int i = 0; i < mesh.GetNE(); i++)
249  {
250  tr = mesh.GetElementTransformation(i);
251  const IntegrationRule &ir = *irs[mesh.GetElementType(i)];
252  for (int j = 0; j < ir.GetNPoints(); j++)
253  {
254  const IntegrationPoint &ip = ir.IntPoint(j);
255  tr->SetIntPoint(&ip);
256  double val = fabs(coeff.Eval(*tr, ip));
257  if (p < numeric_limits<double>::infinity())
258  {
259  norm += ip.weight * tr->Weight() * pow(val, p);
260  }
261  else
262  {
263  if (norm < val)
264  {
265  norm = val;
266  }
267  }
268  }
269  }
270  return norm;
271 }
272 
273 double LpNormLoop(double p, VectorCoefficient &coeff, Mesh &mesh,
274  const IntegrationRule *irs[])
275 {
276  double norm = 0.0;
278  int vdim = coeff.GetVDim();
279  Vector vval(vdim);
280  double val;
281 
282  for (int i = 0; i < mesh.GetNE(); i++)
283  {
284  tr = mesh.GetElementTransformation(i);
285  const IntegrationRule &ir = *irs[mesh.GetElementType(i)];
286  for (int j = 0; j < ir.GetNPoints(); j++)
287  {
288  const IntegrationPoint &ip = ir.IntPoint(j);
289  tr->SetIntPoint(&ip);
290  coeff.Eval(vval, *tr, ip);
291  if (p < numeric_limits<double>::infinity())
292  {
293  for (int idim(0); idim < vdim; ++idim)
294  {
295  norm += ip.weight * tr->Weight() * pow(fabs( vval(idim) ), p);
296  }
297  }
298  else
299  {
300  for (int idim(0); idim < vdim; ++idim)
301  {
302  val = fabs(vval(idim));
303  if (norm < val)
304  {
305  norm = val;
306  }
307  }
308  }
309  }
310  }
311 
312  return norm;
313 }
314 
315 double ComputeLpNorm(double p, Coefficient &coeff, Mesh &mesh,
316  const IntegrationRule *irs[])
317 {
318  double norm = LpNormLoop(p, coeff, mesh, irs);
319 
320  if (p < numeric_limits<double>::infinity())
321  {
322  // negative quadrature weights may cause norm to be negative
323  if (norm < 0.0)
324  {
325  norm = -pow(-norm, 1.0/p);
326  }
327  else
328  {
329  norm = pow(norm, 1.0/p);
330  }
331  }
332 
333  return norm;
334 }
335 
336 double ComputeLpNorm(double p, VectorCoefficient &coeff, Mesh &mesh,
337  const IntegrationRule *irs[])
338 {
339  double norm = LpNormLoop(p, coeff, mesh, irs);
340 
341  if (p < numeric_limits<double>::infinity())
342  {
343  // negative quadrature weights may cause norm to be negative
344  if (norm < 0.0)
345  {
346  norm = -pow(-norm, 1.0/p);
347  }
348  else
349  {
350  norm = pow(norm, 1.0/p);
351  }
352  }
353 
354  return norm;
355 }
356 
357 #ifdef MFEM_USE_MPI
358 double ComputeGlobalLpNorm(double p, Coefficient &coeff, ParMesh &pmesh,
359  const IntegrationRule *irs[])
360 {
361  double loc_norm = LpNormLoop(p, coeff, pmesh, irs);
362  double glob_norm = 0;
363 
364  MPI_Comm comm = pmesh.GetComm();
365 
366  if (p < numeric_limits<double>::infinity())
367  {
368  MPI_Allreduce(&loc_norm, &glob_norm, 1, MPI_DOUBLE, MPI_SUM, comm);
369 
370  // negative quadrature weights may cause norm to be negative
371  if (glob_norm < 0.0)
372  {
373  glob_norm = -pow(-glob_norm, 1.0/p);
374  }
375  else
376  {
377  glob_norm = pow(glob_norm, 1.0/p);
378  }
379  }
380  else
381  {
382  MPI_Allreduce(&loc_norm, &glob_norm, 1, MPI_DOUBLE, MPI_MAX, comm);
383  }
384 
385  return glob_norm;
386 }
387 
388 double ComputeGlobalLpNorm(double p, VectorCoefficient &coeff, ParMesh &pmesh,
389  const IntegrationRule *irs[])
390 {
391  double loc_norm = LpNormLoop(p, coeff, pmesh, irs);
392  double glob_norm = 0;
393 
394  MPI_Comm comm = pmesh.GetComm();
395 
396  if (p < numeric_limits<double>::infinity())
397  {
398  MPI_Allreduce(&loc_norm, &glob_norm, 1, MPI_DOUBLE, MPI_SUM, comm);
399 
400  // negative quadrature weights may cause norm to be negative
401  if (glob_norm < 0.0)
402  {
403  glob_norm = -pow(-glob_norm, 1.0/p);
404  }
405  else
406  {
407  glob_norm = pow(glob_norm, 1.0/p);
408  }
409  }
410  else
411  {
412  MPI_Allreduce(&loc_norm, &glob_norm, 1, MPI_DOUBLE, MPI_MAX, comm);
413  }
414 
415  return glob_norm;
416 }
417 #endif
418 
419 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:228
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Class for integration rule.
Definition: intrules.hpp:83
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
VectorGridFunctionCoefficient(GridFunction *gf)
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)=0
void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
Definition: gridfunc.cpp:316
void SetSize(int s)
Resize the vector if the new size is different.
Definition: vector.hpp:263
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)=0
void SetIntPoint(const IntegrationPoint *ip)
Definition: eltrans.hpp:35
Data type dense matrix using column-major storage.
Definition: densemat.hpp:22
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Definition: coefficient.cpp:83
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:496
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Definition: coefficient.cpp:55
void SetTime(double t)
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate coefficient.
Definition: coefficient.cpp:31
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:231
int dim
Definition: ex3.cpp:47
double ComputeLpNorm(double p, Coefficient &coeff, Mesh &mesh, const IntegrationRule *irs[])
double LpNormLoop(double p, Coefficient &coeff, Mesh &mesh, const IntegrationRule *irs[])
int GetElementType(int i) const
Returns the type of element i.
Definition: mesh.cpp:3680
double Eval(int i, int j, ElementTransformation &T, IntegrationPoint &ip)
virtual ~VectorArrayCoefficient()
Destroys vector coefficient.
virtual double Weight()=0
int GetVDim()
Returns dimension of the vector.
double Eval(int i, ElementTransformation &T, IntegrationPoint &ip)
Evaluates i&#39;th component of the vector.
MPI_Comm GetComm() const
Definition: pmesh.hpp:96
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient function.
Definition: coefficient.cpp:24
void GetColumnReference(int c, Vector &col)
Definition: densemat.hpp:201
void SetTime(double t)
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
double ComputeGlobalLpNorm(double p, Coefficient &coeff, ParMesh &pmesh, const IntegrationRule *irs[])
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Definition: coefficient.cpp:49
Class for integration point with weight.
Definition: intrules.hpp:25
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:237
void GetVectorValues(ElementTransformation &T, const IntegrationRule &ir, DenseMatrix &vals) const
Definition: gridfunc.cpp:429
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
VectorArrayCoefficient(int dim)
Construct vector of dim coefficients.
Vector data type.
Definition: vector.hpp:33
virtual void Transform(const IntegrationPoint &, Vector &)=0
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:71
Class for parallel meshes.
Definition: pmesh.hpp:28