MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
sbm_solver.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2021, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #include "marking.hpp"
13 #include "sbm_solver.hpp"
14 #include "mfem.hpp"
15 
16 namespace mfem
17 {
18 
20  const IntegrationPoint & ip,
21  const Vector &D)
22 {
23  Vector transip;
24  T.Transform(ip, transip);
25  for (int i = 0; i < D.Size(); i++)
26  {
27  transip(i) += D(i);
28  }
29 
30  return Function(transip);
31 }
32 
34  const FiniteElement &el1, const FiniteElement &el2,
36 {
37  int dim, ndof1, ndof2, ndof, ndoftotal;
38  double w;
39  DenseMatrix temp_elmat;
40 
41  dim = el1.GetDim();
42  ndof1 = el1.GetDof();
43  ndof2 = el2.GetDof();
44  ndoftotal = Trans.ElementType == ElementTransformation::BDR_FACE ?
45  ndof1 : ndof1 + ndof2;
46 
47  elmat.SetSize(ndoftotal);
48  elmat = 0.0;
49 
50  bool elem1f = true; // flag indicating whether Trans.Elem1No is part of the
51  // surrogate domain or not.
52  int elem1 = Trans.Elem1No,
53  elem2 = Trans.Elem2No,
54  marker1 = (*elem_marker)[elem1];
55 
56  int marker2;
57  if (Trans.Elem2No >= NEproc)
58  {
59  marker2 = (*elem_marker)[NEproc+par_shared_face_count];
61  }
63  {
64  marker2 = marker1;
65  }
66  else
67  {
68  marker2 = (*elem_marker)[elem2];
69  }
70 
71  if (!include_cut_cell)
72  {
73  // 1 is inside and 2 is cut or 1 is a boundary element.
74  if (marker1 == ShiftedFaceMarker::SBElementType::INSIDE &&
75  (marker2 == ShiftedFaceMarker::SBElementType::CUT ||
77  {
78  elem1f = true;
79  }
80  // 1 is cut, 2 is inside
81  else if (marker1 == ShiftedFaceMarker::SBElementType::CUT &&
82  marker2 == ShiftedFaceMarker::SBElementType::INSIDE)
83  {
84  if (Trans.Elem2No >= NEproc) { return; }
85  elem1f = false;
86  }
87  else
88  {
89  return;
90  }
91  }
92  else
93  {
94  // 1 is cut and 2 is outside or 1 is a boundary element.
95  if (marker1 == ShiftedFaceMarker::SBElementType::CUT &&
96  (marker2 == ShiftedFaceMarker::SBElementType::OUTSIDE ||
98  {
99  elem1f = true;
100  }
101  // 1 is outside, 2 is cut
102  else if (marker1 == ShiftedFaceMarker::SBElementType::OUTSIDE &&
103  marker2 == ShiftedFaceMarker::SBElementType::CUT)
104  {
105  if (Trans.Elem2No >= NEproc) { return; }
106  elem1f = false;
107  }
108  else
109  {
110  return;
111  }
112  }
113 
114  ndof = elem1f ? ndof1 : ndof2;
115 
116  temp_elmat.SetSize(ndof);
117  temp_elmat = 0.;
118 
119  nor.SetSize(dim);
120  nh.SetSize(dim);
121  ni.SetSize(dim);
122  adjJ.SetSize(dim);
123 
124  shape.SetSize(ndof);
125  dshape.SetSize(ndof, dim);
126  dshapephys.SetSize(ndof, dim);
127  Vector dshapephysdd(ndof);
128  dshapedn.SetSize(ndof);
129  Vector wrk = shape;
130 
131  const IntegrationRule *ir = IntRule;
132  if (ir == NULL)
133  {
134  int order = elem1f ? 4*el1.GetOrder() : 4*el2.GetOrder();
135  ir = &IntRules.Get(Trans.GetGeometryType(), order);
136  }
137 
138  Array<DenseMatrix *> dkphi_dxk;
139  DenseMatrix grad_phys;
140  Vector Factorial;
141  Array<DenseMatrix *> grad_phys_dir;
142 
143  if (nterms > 0)
144  {
145  if (elem1f)
146  {
147  el1.ProjectGrad(el1, *Trans.Elem1, grad_phys);
148  }
149  else
150  {
151  el2.ProjectGrad(el2, *Trans.Elem2, grad_phys);
152  }
153 
154  DenseMatrix grad_work;
155  grad_phys_dir.SetSize(dim); // NxN matrices for derivative in each direction
156  for (int i = 0; i < dim; i++)
157  {
158  grad_phys_dir[i] = new DenseMatrix(ndof, ndof);
159  grad_phys_dir[i]->CopyRows(grad_phys, i*ndof, (i+1)*ndof-1);
160  }
161 
162  DenseMatrix grad_phys_work = grad_phys;
163  grad_phys_work.SetSize(ndof, ndof*dim);
164 
165  dkphi_dxk.SetSize(nterms);
166 
167  for (int i = 0; i < nterms; i++)
168  {
169  int sz1 = pow(dim, i+1);
170  dkphi_dxk[i] = new DenseMatrix(ndof, ndof*sz1*dim);
171  int loc_col_per_dof = sz1;
172  int tot_col_per_dof = loc_col_per_dof*dim;
173  for (int k = 0; k < dim; k++)
174  {
175  grad_work.SetSize(ndof, ndof*sz1);
176  // grad_work[k] has derivative in kth direction for each DOF.
177  // grad_work[0] has d^2phi/dx^2 and d^2phi/dxdy terms and
178  // grad_work[1] has d^2phi/dydx and d^2phi/dy2 terms for each dof
179  if (i == 0)
180  {
181  Mult(*grad_phys_dir[k], grad_phys_work, grad_work);
182  }
183  else
184  {
185  Mult(*grad_phys_dir[k], *dkphi_dxk[i-1], grad_work);
186  }
187  // Now we must place columns for each dof together so that they are
188  // in order: d^2phi/dx^2, d^2phi/dxdy, d^2phi/dydx, d^2phi/dy2.
189  for (int j = 0; j < ndof; j++)
190  {
191  for (int d = 0; d < loc_col_per_dof; d++)
192  {
193  Vector col;
194  grad_work.GetColumn(j*loc_col_per_dof+d, col);
195  dkphi_dxk[i]->SetCol(j*tot_col_per_dof+k*loc_col_per_dof+d, col);
196  }
197  }
198  }
199  }
200 
201  for (int i = 0; i < grad_phys_dir.Size(); i++)
202  {
203  delete grad_phys_dir[i];
204  }
205 
206  Factorial.SetSize(nterms);
207  Factorial(0) = 2;
208  for (int i = 1; i < nterms; i++)
209  {
210  Factorial(i) = Factorial(i-1)*(i+2);
211  }
212  }
213 
214  DenseMatrix q_hess_dn(dim, ndof);
215  Vector q_hess_dn_work(q_hess_dn.GetData(), ndof*dim);
216  Vector q_hess_dot_d(ndof);
217 
218  Vector D(vD->GetVDim());
219  // Assemble: -< \nabla u.n, w >
220  // -< u + \nabla u.d + h.o.t, \nabla w.n>
221  // -<alpha h^{-1} (u + \nabla u.d + h.o.t), w + \nabla w.d + h.o.t>
222  for (int p = 0; p < ir->GetNPoints(); p++)
223  {
224  const IntegrationPoint &ip = ir->IntPoint(p);
225 
226  // Set the integration point in the face and the neighboring elements
227  Trans.SetAllIntPoints(&ip);
228 
229  // Access the neighboring elements' integration points
230  // Note: eip2 will only contain valid data if Elem2 exists
231  const IntegrationPoint &eip1 = Trans.GetElement1IntPoint();
232  const IntegrationPoint &eip2 = Trans.GetElement2IntPoint();
233 
234  if (dim == 1)
235  {
236  nor(0) = 2*eip1.x - 1.0;
237  }
238  else
239  {
240  // Note: this normal accounts for the weight of the surface transformation
241  // Jacobian i.e. nor = nhat*det(J)
242  CalcOrtho(Trans.Jacobian(), nor);
243  }
244  vD->Eval(D, Trans, ip);
245 
246  double nor_dot_d = nor*D;
247  // If we are clipping inside the domain, ntilde and d vector should be
248  // aligned.
249  if (!include_cut_cell && nor_dot_d < 0) { nor *= -1; }
250  if (include_cut_cell && nor_dot_d > 0) { nor *= -1; }
251 
252  if (elem1f)
253  {
254  el1.CalcShape(eip1, shape);
255  el1.CalcDShape(eip1, dshape);
256  w = ip.weight/Trans.Elem1->Weight();
257  CalcAdjugate(Trans.Elem1->Jacobian(), adjJ);
258  }
259  else
260  {
261  el1.CalcShape(eip2, shape);
262  el1.CalcDShape(eip2, dshape);
263  w = ip.weight/Trans.Elem2->Weight();
264  CalcAdjugate(Trans.Elem2->Jacobian(), adjJ);
265  }
266 
267  ni.Set(w, nor); // alpha_k*nor/det(J)
268  adjJ.Mult(ni, nh);
269  dshape.Mult(nh, dshapedn); // dphi/dn * Jinv * alpha_k * nor
270 
271  // <grad u.n, w> - Term 2
272  AddMult_a_VWt(-1., shape, dshapedn, temp_elmat);
273 
274  if (elem1f) { el1.CalcPhysDShape(*(Trans.Elem1), dshapephys); }
275  else { el1.CalcPhysDShape(*(Trans.Elem2), dshapephys); } // dphi/dx
276  dshapephys.Mult(D, dshapephysdd); // dphi/dx.D);
277 
278  q_hess_dot_d = 0.;
279  for (int i = 0; i < nterms; i++)
280  {
281  int sz1 = pow(dim, i+1);
282  DenseMatrix T1(dim, ndof*sz1);
283  Vector T1_wrk(T1.GetData(), dim*ndof*sz1);
284  dkphi_dxk[i]->MultTranspose(shape, T1_wrk);
285 
286  DenseMatrix T2;
287  Vector T2_wrk;
288  for (int j = 0; j < i+1; j++)
289  {
290  int sz2 = pow(dim, i-j);
291  T2.SetSize(dim, ndof*sz2);
292  T2_wrk.SetDataAndSize(T2.GetData(), dim*ndof*sz2);
293  T1.MultTranspose(D, T2_wrk);
294  T1 = T2;
295  }
296  Vector q_hess_dot_d_work(ndof);
297  T1.MultTranspose(D, q_hess_dot_d_work);
298  q_hess_dot_d_work *= 1./Factorial(i);
299  q_hess_dot_d += q_hess_dot_d_work;
300  }
301 
302  wrk = shape;
303  wrk += dshapephysdd;
304  wrk += q_hess_dot_d;
305  // <u + grad u.d + h.o.t, grad w.n> - Term 3
306  AddMult_a_VWt(-1., dshapedn, wrk, temp_elmat);
307 
308  double hinvdx;
309  if (elem1f) { hinvdx = nor*nor/Trans.Elem1->Weight(); }
310  else { hinvdx = nor*nor/Trans.Elem2->Weight(); }
311 
312  w = ip.weight*alpha*hinvdx;
313  // + <alpha * hinv * u + grad u.d + h.o.t, w + grad w.d + h.o.t> - Term 4
314  AddMult_a_VVt(w, wrk, temp_elmat);
315 
316  int offset = elem1f ? 0 : ndof1;
317  elmat.CopyMN(temp_elmat, offset, offset);
318  } // p < ir->GetNPoints()
319 
320  for (int i = 0; i < dkphi_dxk.Size(); i++)
321  {
322  delete dkphi_dxk[i];
323  }
324 }
325 
326 
328  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
329 {
330  mfem_error("SBM2DirichletLFIntegrator::AssembleRHSElementVect");
331 }
332 
334  const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
335 {
336  AssembleRHSElementVect(el, el, Tr, elvect);
337 }
338 
340  const FiniteElement &el1, const FiniteElement &el2,
341  FaceElementTransformations &Tr, Vector &elvect)
342 {
343  int dim, ndof1, ndof2, ndof, ndoftotal;
344  double w;
345  Vector temp_elvect;
346 
347  dim = el1.GetDim();
348  ndof1 = el1.GetDof();
349  ndof2 = el2.GetDof();
350  ndoftotal = ndof1 + ndof2;
351  if (Tr.Elem2No >= NEproc ||
353  {
354  ndoftotal = ndof1;
355  }
356 
357  elvect.SetSize(ndoftotal);
358  elvect = 0.0;
359 
360  bool elem1f = true;
361  int elem1 = Tr.Elem1No,
362  elem2 = Tr.Elem2No,
363  marker1 = (*elem_marker)[elem1];
364 
365  int marker2;
366  if (Tr.Elem2No >= NEproc)
367  {
368  marker2 = (*elem_marker)[NEproc+par_shared_face_count];
370  }
372  {
373  marker2 = marker1;
374  }
375  else
376  {
377  marker2 = (*elem_marker)[elem2];
378  }
379 
380  if (!include_cut_cell)
381  {
382  // 1 is inside and 2 is cut or 1 is a boundary element.
383  if ( marker1 == ShiftedFaceMarker::SBElementType::INSIDE &&
384  (marker2 == ShiftedFaceMarker::SBElementType::CUT ||
386  {
387  elem1f = true;
388  ndof = ndof1;
389  }
390  // 1 is cut, 2 is inside
391  else if (marker1 == ShiftedFaceMarker::SBElementType::CUT &&
392  marker2 == ShiftedFaceMarker::SBElementType::INSIDE)
393  {
394  if (Tr.Elem2No >= NEproc) { return; }
395  elem1f = false;
396  ndof = ndof2;
397  }
398  else
399  {
400  return;
401  }
402  }
403  else
404  {
405  // 1 is cut and 2 is outside or 1 is a boundary element.
406  if (marker1 == ShiftedFaceMarker::SBElementType::CUT &&
407  (marker2 == ShiftedFaceMarker::SBElementType::OUTSIDE ||
409  {
410  elem1f = true;
411  ndof = ndof1;
412  }
413  // 1 is outside, 2 is cut
414  else if (marker1 == ShiftedFaceMarker::SBElementType::OUTSIDE &&
415  marker2 == ShiftedFaceMarker::SBElementType::CUT)
416  {
417  if (Tr.Elem2No >= NEproc) { return; }
418  elem1f = false;
419  ndof = ndof2;
420  }
421  else
422  {
423  return;
424  }
425  }
426 
427  temp_elvect.SetSize(ndof);
428  temp_elvect = 0.0;
429 
430  nor.SetSize(dim);
431  nh.SetSize(dim);
432  ni.SetSize(dim);
433  adjJ.SetSize(dim);
434 
435  shape.SetSize(ndof);
436  dshape.SetSize(ndof, dim);
437  dshape_dd.SetSize(ndof);
438  dshape_dn.SetSize(ndof);
439 
440  const IntegrationRule *ir = IntRule;
441  if (ir == NULL)
442  {
443  // a simple choice for the integration order; is this OK?
444  int order = elem1f ? 4*el1.GetOrder() : 4*el2.GetOrder();
445  ir = &IntRules.Get(Tr.GetGeometryType(), order);
446  }
447 
448  Array<DenseMatrix *> dkphi_dxk;
449  DenseMatrix grad_phys;
450  Vector Factorial;
451  Array<DenseMatrix *> grad_phys_dir;
452 
453  if (nterms > 0)
454  {
455  if (elem1f)
456  {
457  el1.ProjectGrad(el1, *Tr.Elem1, grad_phys);
458  }
459  else
460  {
461  el2.ProjectGrad(el2, *Tr.Elem2, grad_phys);
462  }
463 
464  DenseMatrix grad_work;
465  grad_phys_dir.SetSize(dim); // NxN matrix for derivative in each direction
466  for (int i = 0; i < dim; i++)
467  {
468  grad_phys_dir[i] = new DenseMatrix(ndof, ndof);
469  grad_phys_dir[i]->CopyRows(grad_phys, i*ndof, (i+1)*ndof-1);
470  }
471 
472 
473  DenseMatrix grad_phys_work = grad_phys;
474  grad_phys_work.SetSize(ndof, ndof*dim);
475 
476  dkphi_dxk.SetSize(nterms);
477 
478  for (int i = 0; i < nterms; i++)
479  {
480  int sz1 = pow(dim, i+1);
481  dkphi_dxk[i] = new DenseMatrix(ndof, ndof*sz1*dim);
482  int loc_col_per_dof = sz1;
483  for (int k = 0; k < dim; k++)
484  {
485  grad_work.SetSize(ndof, ndof*sz1);
486  // grad_work[k] has derivative in kth direction for each DOF.
487  // grad_work[0] has d^2phi/dx^2 and d^2phi/dxdy terms and
488  // grad_work[1] has d^2phi/dydx and d^2phi/dy2 terms for each dof
489  if (i == 0)
490  {
491  Mult(*grad_phys_dir[k], grad_phys_work, grad_work);
492  }
493  else
494  {
495  Mult(*grad_phys_dir[k], *dkphi_dxk[i-1], grad_work);
496  }
497  // Now we must place columns for each dof together so that they are
498  // in order: d^2phi/dx^2, d^2phi/dxdy, d^2phi/dydx, d^2phi/dy2.
499  for (int j = 0; j < ndof; j++)
500  {
501  for (int d = 0; d < loc_col_per_dof; d++)
502  {
503  Vector col;
504  int tot_col_per_dof = loc_col_per_dof*dim;
505  grad_work.GetColumn(j*loc_col_per_dof+d, col);
506  dkphi_dxk[i]->SetCol(j*tot_col_per_dof+k*loc_col_per_dof+d, col);
507  }
508  }
509  }
510  }
511 
512  for (int i = 0; i < grad_phys_dir.Size(); i++)
513  {
514  delete grad_phys_dir[i];
515  }
516 
517  Factorial.SetSize(nterms);
518  Factorial(0) = 2;
519  for (int i = 1; i < nterms; i++)
520  {
521  Factorial(i) = Factorial(i-1)*(i+2);
522  }
523  }
524 
525  DenseMatrix q_hess_dn(dim, ndof);
526  Vector q_hess_dn_work(q_hess_dn.GetData(), ndof*dim);
527  Vector q_hess_dot_d(ndof);
528 
529  Vector D(vD->GetVDim());
530  Vector wrk = shape;
531  // Assemble: -< u_D, \nabla w.n >
532  // -<alpha h^{-1} u_D, w + \nabla w.d + h.o.t>
533  for (int p = 0; p < ir->GetNPoints(); p++)
534  {
535 
536  const IntegrationPoint &ip = ir->IntPoint(p);
537 
538  // Set the integration point in the face and the neighboring element
539  Tr.SetAllIntPoints(&ip);
540 
541  // Access the neighboring element's integration point
542  const IntegrationPoint &eip = Tr.GetElement1IntPoint();
543  const IntegrationPoint &eip1 = Tr.GetElement1IntPoint();
544  const IntegrationPoint &eip2 = Tr.GetElement2IntPoint();
545 
546  if (dim == 1)
547  {
548  nor(0) = 2*eip.x - 1.0;
549  }
550  else
551  {
552  CalcOrtho(Tr.Jacobian(), nor);
553  }
554  vD->Eval(D, Tr, ip);
555 
556  double nor_dot_d = nor*D;
557  if (!include_cut_cell && nor_dot_d < 0) { nor *= -1; }
558  if (include_cut_cell && nor_dot_d > 0) { nor *= -1; }
559  // note here that if we are clipping outside the domain, we will have to
560  // flip the sign if nor_dot_d is positive.
561 
562  double hinvdx;
563 
564  if (elem1f)
565  {
566  el1.CalcShape(eip1, shape);
567  el1.CalcDShape(eip1, dshape);
568  hinvdx =nor*nor/Tr.Elem1->Weight();
569  w = ip.weight * uD->Eval(Tr, ip, D) / Tr.Elem1->Weight();
571  }
572  else
573  {
574  el2.CalcShape(eip2, shape);
575  el2.CalcDShape(eip2, dshape);
576  hinvdx = nor*nor/Tr.Elem2->Weight();
577  w = ip.weight * uD->Eval(Tr, ip, D) / Tr.Elem2->Weight();
579  }
580 
581  ni.Set(w, nor);
582  adjJ.Mult(ni, nh);
583 
585  temp_elvect.Add(-1., dshape_dn); // T2
586 
587  double jinv;
588  if (elem1f)
589  {
590  w = ip.weight * uD->Eval(Tr, ip, D) * alpha * hinvdx;
591  jinv = 1./Tr.Elem1->Weight();
592  }
593  else
594  {
595  w = ip.weight * uD->Eval(Tr, ip, D) * alpha * hinvdx;
596  jinv = 1./Tr.Elem2->Weight();
597  }
598  adjJ.Mult(D, nh);
599  nh *= jinv;
601 
602  q_hess_dot_d = 0.;
603  for (int i = 0; i < nterms; i++)
604  {
605  int sz1 = pow(dim, i+1);
606  DenseMatrix T1(dim, ndof*sz1);
607  Vector T1_wrk(T1.GetData(), dim*ndof*sz1);
608  dkphi_dxk[i]->MultTranspose(shape, T1_wrk);
609 
610  DenseMatrix T2;
611  Vector T2_wrk;
612  for (int j = 0; j < i+1; j++)
613  {
614  int sz2 = pow(dim, i-j);
615  T2.SetSize(dim, ndof*sz2);
616  T2_wrk.SetDataAndSize(T2.GetData(), dim*ndof*sz2);
617  T1.MultTranspose(D, T2_wrk);
618  T1 = T2;
619  }
620  Vector q_hess_dot_d_work(ndof);
621  T1.MultTranspose(D, q_hess_dot_d_work);
622  q_hess_dot_d_work *= 1./Factorial(i);
623  q_hess_dot_d += q_hess_dot_d_work;
624  }
625 
626  wrk = shape;
627  wrk += dshape_dd; // \grad w .d
628  wrk += q_hess_dot_d;
629  temp_elvect.Add(w, wrk); // <u, gradw.d>
630 
631  int offset = elem1f ? 0 : ndof1;
632  for (int i = 0; i < temp_elvect.Size(); i++)
633  {
634  elvect(i+offset) = temp_elvect(i);
635  }
636  }
637 
638  for (int i = 0; i < dkphi_dxk.Size(); i++)
639  {
640  delete dkphi_dxk[i];
641  }
642 }
643 
644 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:247
Abstract class for all finite elements.
Definition: fe.hpp:243
int Size() const
Return the logical size of the array.
Definition: array.hpp:134
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe.hpp:317
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:920
Geometry::Type GetGeometryType() const
Return the Geometry::Type of the reference element.
Definition: eltrans.hpp:149
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the vector coefficient in the element described by T at the point ip, storing the result in ...
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:513
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:476
VectorCoefficient * vD
Definition: sbm_solver.hpp:60
A specialized ElementTransformation class representing a face and its two neighboring elements...
Definition: eltrans.hpp:467
void CalcAdjugate(const DenseMatrix &a, DenseMatrix &adja)
Definition: densemat.cpp:2089
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe.hpp:327
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:190
const IntegrationPoint & GetElement2IntPoint()
Get a const reference to the integration point in neighboring element 2 corresponding to the currentl...
Definition: eltrans.hpp:574
const IntegrationRule * IntRule
Definition: lininteg.hpp:25
double * GetData() const
Returns the matrix data array.
Definition: densemat.hpp:115
ElementTransformation * Elem2
Definition: eltrans.hpp:509
virtual void ProjectGrad(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Compute the discrete gradient matrix from the given FiniteElement onto &#39;this&#39; FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
Definition: fe.cpp:161
void CalcOrtho(const DenseMatrix &J, Vector &n)
Definition: densemat.cpp:2303
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:250
void MultTranspose(const double *x, double *y) const
Multiply a vector with the transpose matrix.
Definition: densemat.cpp:203
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:153
void AddMult_a_VWt(const double a, const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += a * v w^t.
Definition: densemat.cpp:2821
void SetAllIntPoints(const IntegrationPoint *face_ip)
Set the integration point in the Face and the two neighboring elements, if present.
Definition: eltrans.hpp:553
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().
Definition: eltrans.hpp:111
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: sbm_solver.cpp:327
double Weight()
Return the weight of the Jacobian matrix of the transformation at the currently set IntegrationPoint...
Definition: eltrans.hpp:123
void AddMult_a_VVt(const double a, const Vector &v, DenseMatrix &VVt)
VVt += a * v v^t.
Definition: densemat.cpp:2843
void CalcPhysDShape(ElementTransformation &Trans, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in physical space at the poi...
Definition: fe.cpp:192
int GetVDim()
Returns dimension of the vector.
void GetColumn(int c, Vector &col) const
Definition: densemat.cpp:1285
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe.hpp:323
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:674
const IntegrationPoint & GetElement1IntPoint()
Get a const reference to the integration point in neighboring element 1 corresponding to the currentl...
Definition: eltrans.hpp:564
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient in the element described by T at the point ip.
Definition: sbm_solver.hpp:31
void SetDataAndSize(double *d, int s)
Set the Vector data and size.
Definition: vector.hpp:147
const IntegrationRule * IntRule
Definition: nonlininteg.hpp:30
Vector & Set(const double a, const Vector &x)
(*this) = a * x
Definition: vector.cpp:258
std::function< double(const Vector &)> Function
Definition: sbm_solver.hpp:25
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Definition: vector.cpp:243
Class for integration point with weight.
Definition: intrules.hpp:25
int dim
Definition: ex24.cpp:53
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
Definition: sbm_solver.cpp:33
ElementTransformation * Elem1
Definition: eltrans.hpp:509
Vector data type.
Definition: vector.hpp:60
void Mult(const double *x, double *y) const
Matrix vector multiplication.
Definition: densemat.cpp:175
virtual void Transform(const IntegrationPoint &, Vector &)=0
Transform integration point from reference coordinates to physical coordinates and store them in the ...
void CopyMN(const DenseMatrix &A, int m, int n, int Aro, int Aco)
Copy the m x n submatrix of A at row/col offsets Aro/Aco to *this.
Definition: densemat.cpp:1532
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const =0
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:377
ShiftedFunctionCoefficient * uD
Definition: sbm_solver.hpp:115