v0.14.0
LinearFormsIntegrators.hpp
Go to the documentation of this file.
1 /** \file LinearFormsIntegrators.hpp
2  * \brief Linear forms integrators
3  * \ingroup mofem_form
4 
5 */
6 
7 #ifndef __LINEAR_FORMS_INTEGRATORS_HPP__
8 #define __LINEAR_FORMS_INTEGRATORS_HPP__
9 
10 namespace MoFEM {
11 
13  template <typename OpBase> struct S { S() = delete; };
15 };
16 
18  template <typename OpBase> struct S { S() = delete; };
20 };
21 
22 template <int BASE_DIM, int FIELD_DIM, IntegrationType I, typename OpBase>
23 struct OpSourceImpl;
24 
25 /**
26  * @brief Integrate source
27  *
28  * @tparam OpBase
29  */
30 template <typename OpBase>
32  : public OpBase {
33 
34  /**
35  * @brief Construct a new Op Source Impl object
36  *
37  * @param field_name
38  * @param source_fun
39  * @param ents_ptr
40  */
41  OpSourceImpl(const std::string field_name, ScalarFun source_fun,
42  boost::shared_ptr<Range> ents_ptr = nullptr)
43  : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
44  sourceFun(source_fun) {}
45 
46  /**
47  * @brief Construct a new Op Source Impl object
48  *
49  * @param field_name
50  * @param time_fun
51  * @param source_fun
52  * @param ents_ptr
53  */
54  OpSourceImpl(const std::string field_name, TimeFun time_fun,
55  ScalarFun source_fun,
56  boost::shared_ptr<Range> ents_ptr = nullptr)
57  : OpBase(field_name, field_name, OpBase::OPROW, time_fun, ents_ptr),
58  sourceFun(source_fun) {}
59 
60 protected:
63 };
64 
65 template <int FIELD_DIM, typename OpBase>
67  SourceFunctionSpecialization::S<OpBase>> : public OpBase {
68 
69  /**
70  * @brief Construct a new Op Source Impl object
71  *
72  * @param field_name
73  * @param time_fun
74  * @param source_fun
75  * @param ents_ptr
76  */
77  OpSourceImpl(const std::string field_name, TimeFun time_fun,
78  VectorFun<FIELD_DIM> source_fun,
79  boost::shared_ptr<Range> ents_ptr = nullptr)
80  : OpBase(field_name, field_name, OpBase::OPROW, time_fun, ents_ptr) {}
81 
82  /**
83  * @brief Construct a new Op Source Impl object
84  *
85  * @param field_name
86  * @param source_fun
87  * @param ents_ptr
88  */
89  OpSourceImpl(const std::string field_name, VectorFun<FIELD_DIM> source_fun,
90  boost::shared_ptr<Range> ents_ptr = nullptr)
91  : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
92  sourceFun(source_fun) {}
93 
94 protected:
97 };
98 
99 template <int FIELD_DIM, typename OpBase>
101  SourceFunctionSpecialization::S<OpBase>> : public OpBase {
102 
103  OpSourceImpl(const std::string field_name, TimeFun time_fun,
104  VectorFun<FIELD_DIM> source_fun,
105  boost::shared_ptr<Range> ents_ptr = nullptr)
106  : OpBase(field_name, field_name, OpBase::OPROW, time_fun, ents_ptr),
107  sourceFun(source_fun) {}
108 
109  OpSourceImpl(const std::string field_name, VectorFun<FIELD_DIM> source_fun,
110  boost::shared_ptr<Range> ents_ptr = nullptr)
111  : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
112  sourceFun(source_fun) {}
113 
114 protected:
117 };
118 
119 template <int BASE_DIM, int S, IntegrationType I, typename OpBase>
121 
122 template <int S, typename OpBase>
123 struct OpBaseTimesScalarImpl<1, S, GAUSS, OpBase> : public OpBase {
124 
126  const std::string field_name, boost::shared_ptr<VectorDouble> vec,
127  ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; },
128  boost::shared_ptr<Range> ents_ptr = nullptr)
129  : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr), sourceVec(vec),
130  betaCoeff(beta_coeff) {}
131 
132 protected:
134  boost::shared_ptr<VectorDouble> sourceVec;
136 };
137 
138 template <int BASE_DIM, int FIELD_DIM, int S, IntegrationType I,
139  typename OpBase>
141 
142 template <int FIELD_DIM, int S, typename OpBase>
144 
146  const std::string field_name, boost::shared_ptr<MatrixDouble> vec,
147  ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; },
148  boost::shared_ptr<Range> ents_ptr = nullptr)
149  : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr), sourceVec(vec),
150  betaCoeff(beta_coeff) {}
151 
152 protected:
154  boost::shared_ptr<MatrixDouble> sourceVec;
157 };
158 
159 template <int FIELD_DIM, int S, typename OpBase>
161 
163  const std::string field_name, boost::shared_ptr<MatrixDouble> vec,
164  ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; },
165  boost::shared_ptr<Range> ents_ptr = nullptr)
166  : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr), sourceVec(vec),
167  betaCoeff(beta_coeff) {}
168 
169 protected:
171  boost::shared_ptr<MatrixDouble> sourceVec;
174 };
175 
176 template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S, IntegrationType I,
177  typename OpBase>
179 
180 template <int SPACE_DIM, int S, typename OpBase>
182  : public OpBase {
183 
185 
187  const std::string field_name, boost::shared_ptr<MatrixDouble> mat_vals,
188  ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; },
189  boost::shared_ptr<Range> ents_ptr = nullptr)
190  : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
191  matVals(mat_vals), betaCoeff(beta_coeff) {}
192 
193 protected:
194  boost::shared_ptr<MatrixDouble> matVals;
197 };
198 
199 template <int SPACE_DIM, int S, typename OpBase>
201  : public OpBase {
202 
205 
207  const std::string field_name, boost::shared_ptr<MatrixDouble> mat_vals,
208  ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; },
209  boost::shared_ptr<Range> ents_ptr = nullptr)
210  : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
211  matVals(mat_vals), betaCoeff(beta_coeff) {}
212 
213 protected:
214  boost::shared_ptr<MatrixDouble> matVals;
217 };
218 
219 template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S, IntegrationType I,
220  typename OpBase>
222 
223 template <int SPACE_DIM, int S, typename OpBase>
225  : public OpBase {
226 
228  const std::string field_name, boost::shared_ptr<MatrixDouble> mat_vals,
229  ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; })
230  : OpBase(field_name, field_name, OpBase::OPROW), matVals(mat_vals),
231  betaCoeff(beta_coeff) {}
232 
233 protected:
234  boost::shared_ptr<MatrixDouble> matVals;
239 };
240 
241 template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM, IntegrationType I,
242  typename OpBase, CoordinateTypes CoordSys>
244 
245 template <int FIELD_DIM, int SPACE_DIM, typename OpBase,
246  CoordinateTypes CoordSys>
248  : public OpBase {
250  const std::string field_name, boost::shared_ptr<MatrixDouble> mat_vals,
251  ScalarFun beta = [](double, double, double) { return 1; },
252  boost::shared_ptr<Range> ents_ptr = nullptr)
253  : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
254  matVals(mat_vals), betaConst(beta) {}
255 
256 protected:
258  boost::shared_ptr<MatrixDouble> matVals;
262 };
263 
264 template <int SPACE_DIM, typename OpBase, CoordinateTypes CoordSys>
265 struct OpMixDivTimesUImpl<3, 1, SPACE_DIM, GAUSS, OpBase, CoordSys>
266  : public OpBase {
268  const std::string field_name, boost::shared_ptr<VectorDouble> vec_vals,
269  ScalarFun beta = [](double, double, double) constexpr { return 1; },
270  boost::shared_ptr<Range> ents_ptr = nullptr)
271  : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
272  vecVals(vec_vals), betaConst(beta) {}
273 
274 protected:
276  boost::shared_ptr<VectorDouble> vecVals;
279 };
280 
281 template <int FIELD_DIM, typename OpBase, CoordinateTypes CoordSys>
283  : public OpBase {
284 
286  const std::string field_name, boost::shared_ptr<VectorDouble> vec,
287  ScalarFun beta = [](double, double, double) constexpr { return 1; },
288  boost::shared_ptr<Range> ents_ptr = nullptr)
289  : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr), sourceVec(vec),
290  betaCoeff(beta) {}
291 
292 protected:
294  boost::shared_ptr<VectorDouble> sourceVec;
296 };
297 
298 /**
299  * @brief Vector filed time divergence of tensor
300  *
301  * \f[
302  * (v_i,\lambda_{ij,j})_\Omega
303  * \f]
304  *
305  * @tparam SPACE_DIM
306  * @tparam I
307  * @tparam OpBase
308  */
309 template <int SPACE_DIM, IntegrationType I, typename OpBase>
311 
312 template <int SPACE_DIM, typename OpBase>
315  boost::shared_ptr<MatrixDouble> mat_vals)
316  : OpBase(field_name, field_name, OpBase::OPROW), matVals(mat_vals) {}
317 
319  boost::shared_ptr<MatrixDouble> mat_vals,
320  ScalarFun beta_fun)
321  : OpBase(field_name, field_name, OpBase::OPROW), matVals(mat_vals),
322  betaCoeff(beta_fun) {}
323 
324 protected:
325  boost::shared_ptr<MatrixDouble> matVals;
326  ScalarFun betaCoeff = [](double, double, double) constexpr { return 1; };
327 
330 };
331 
332 /**
333  * @brief Tensor field time gradient of vector field
334  *
335  * \f[
336  * (u_i,\lambda_{ij,j})_\Omega
337  * \f]
338  *
339  * @tparam SPACE_DIM
340  * @tparam I
341  * @tparam OpBase
342  */
343 template <int SPACE_DIM, IntegrationType I, typename OpBase>
345 
346 template <int SPACE_DIM, typename OpBase>
349  boost::shared_ptr<MatrixDouble> mat_vals)
350  : OpBase(field_name, field_name, OpBase::OPROW), matVals(mat_vals) {}
351 
353  boost::shared_ptr<MatrixDouble> mat_vals,
354  ScalarFun beta_fun)
355  : OpBase(field_name, field_name, OpBase::OPROW), matVals(mat_vals),
356  betaCoeff(beta_fun) {}
357 
358 protected:
359  boost::shared_ptr<MatrixDouble> matVals;
360  ScalarFun betaCoeff = [](double, double, double) constexpr { return 1; };
361 
365 };
366 
367 /**
368  * @brief Multiply vector times normal on the face times scalar function
369  *
370  * This operator typically will be used to evaluate natural boundary conditions
371  * for mixed formulation.
372  *
373  * @tparam BASE_DIM
374  * @tparam SPACE_DIM
375  * @tparam OpBase
376  */
377 template <int SPACE_DIM, IntegrationType I, typename OpBase>
379 
380 /**
381  * @brief This is specialisation for sources on boundary which depends on normal
382  *
383  * @tparam BASE_DIM
384  * @tparam OpBase
385  */
386 template <int FIELD_DIM, IntegrationType I, typename OpBase>
389  : public OpNormalMixVecTimesScalarImpl<FIELD_DIM, I, OpBase> {
392 };
393 
394 template <typename OpBase>
397  const std::string field_name,
398  ScalarFun source_fun = [](double, double, double) constexpr { return 1; },
399  boost::shared_ptr<Range> ents_ptr = nullptr)
400  : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
401  sourceFun(source_fun) {}
402 
403 protected:
407 };
408 
409 template <typename OpBase>
412  const std::string field_name,
413  ScalarFun source_fun = [](double, double, double) constexpr { return 1; },
414  boost::shared_ptr<Range> ents_ptr = nullptr)
415  : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
416  sourceFun(source_fun) {}
417 
418 protected:
422 };
423 
424 /**
425  * @brief Multiply vector times normal on the face times vector field
426  *
427  * This operator typically will be used to evaluate natural boundary conditions
428  * for mixed formulation.
429  *
430  * @tparam BASE_DIM
431  * @tparam SPACE_DIM
432  * @tparam OpBase
433  */
434 template <int SPACE_DIM, IntegrationType I, typename OpBase>
436 
437 template <int SPACE_DIM, typename OpBase>
439  : public OpBase {
441  const std::string field_name, boost::shared_ptr<MatrixDouble> u_ptr,
442  ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; },
443  boost::shared_ptr<Range> ents_ptr = nullptr)
444  : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr), uPtr(u_ptr),
445  betaCoeff(beta_coeff) {}
446 
447 protected:
448  boost::shared_ptr<MatrixDouble> uPtr;
451 };
452 
453 template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM, IntegrationType I,
454  typename OpBase>
456 
457 template <int SPACE_DIM, typename OpBase>
460  const std::string field_name, boost::shared_ptr<MatrixDouble> u_ptr,
461  boost::shared_ptr<MatrixDouble> y_grad_ptr,
462  ConstantFun source_fun = []() constexpr { return 1; })
463  : OpBase(field_name, field_name, OpBase::OPROW), uPtr(u_ptr),
464  yGradPtr(y_grad_ptr), alphaConstant(source_fun) {}
465 
466 protected:
467  boost::shared_ptr<MatrixDouble> uPtr;
468  boost::shared_ptr<MatrixDouble> yGradPtr;
471 };
472 
473 template <int FIELD_DIM, int SPACE_DIM, typename OpBase>
475  : public OpBase {
477  const std::string field_name, boost::shared_ptr<MatrixDouble> u_ptr,
478  boost::shared_ptr<MatrixDouble> y_grad_ptr,
479  ConstantFun source_fun = []() constexpr { return 1; })
480  : OpBase(field_name, field_name, OpBase::OPROW), uPtr(u_ptr),
481  yGradPtr(y_grad_ptr), alphaConstant(source_fun) {}
482 
483 protected:
484  boost::shared_ptr<MatrixDouble> uPtr;
485  boost::shared_ptr<MatrixDouble> yGradPtr;
488 };
489 
490 /**
491  * @brief Linear integrator form
492  * @ingroup mofem_forms
493  *
494  * @tparam EleOp
495  * @tparam A
496  * @tparam I
497  */
498 template <typename EleOp>
499 template <AssemblyType A>
500 template <IntegrationType I>
501 struct FormsIntegrators<EleOp>::Assembly<A>::LinearForm {
502 
503  /**
504  * @brief Integrate \f$(v,f(\mathbf{x}))_\Omega\f$, f is a scalar
505  * @ingroup mofem_forms
506  *
507  * @note \f$f(\mathbf{x})\$ is scalar function of coordinates
508  *
509  * @tparam BASE_DIM
510  * @tparam FIELD_DIM
511  */
512  template <int BASE_DIM, int FIELD_DIM,
513  typename S = SourceFunctionSpecialization>
514  using OpSource =
516 
517  /**
518  * @brief Vector field integrator \f$(v_i,f)_\Omega\f$, f is a vector
519  * @ingroup mofem_forms
520  *
521  * @tparam BASE_DIM
522  */
523  template <int BASE_DIM, int S = 1>
524  using OpBaseTimesScalar = OpBaseTimesScalarImpl<BASE_DIM, S, I, OpBase>;
525 
526  /** @deprecated use instead OpBaseTimesScalar
527  */
528  template <int BASE_DIM, int S = 1>
529  using OpBaseTimesScalarField = OpBaseTimesScalar<BASE_DIM, S>;
530 
531  /**
532  * @brief Vector field integrator \f$(v,f_i)_\Omega\f$, f is a vector
533  * @ingroup mofem_forms
534  *
535  * @tparam BASE_DIM
536  * @tparam FIELD_DIM
537  * @tparam 0
538  */
539  template <int BASE_DIM, int FIELD_DIM, int S>
540  using OpBaseTimesVector =
542  //! [Source operator]
543 
544  //! [Grad times tensor]
545 
546  /**
547  * @brief Integrate \f$(v_{,i},f_{ij})_\Omega\f$, f is a vector
548  * @ingroup mofem_forms
549  *
550  * @note \f$f_{ij}\$ is tensor at integration points
551  *
552  * @tparam BASE_DIM
553  * @tparam FIELD_DIM
554  * @tparam SPACE_DIM
555  */
556  template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S = 1>
557  using OpGradTimesTensor =
559 
560  /**
561  * @brief Integrate \f$(v_{,i},f_{ij})_\Omega\f$, f is symmetric tensor
562  * @ingroup mofem_forms
563  *
564  * @note \f$f_{ij}\$ is tensor at integration points
565  *
566  * @tparam BASE_DIM
567  * @tparam FIELD_DIM
568  * @tparam SPACE_DIM
569  */
570  template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S = 1>
571  using OpGradTimesSymTensor =
573 
574  /**
575  * @brief Integrate \f$(\lambda_{ij,j},u_{i})_\Omega\f$
576  *
577  * @tparam BASE_DIM
578  * @tparam FIELD_DIM
579  * @tparam SPACE_DIM
580  */
581  template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM,
582  CoordinateTypes CoordSys = CARTESIAN>
583  using OpMixDivTimesU =
585 
586  /**
587  * @brief Integrate \f$(\lambda_{ij},u_{i,j})_\Omega\f$
588  *
589  * @tparam SPACE_DIM
590  */
591  template <int SPACE_DIM>
592  using OpMixTensorTimesGradU = OpMixTensorTimesGradUImpl<SPACE_DIM, I, OpBase>;
593 
594  /**
595  * @brief Integrate \f$(u_{i},\lambda_{ij,j})_\Omega\f$
596  *
597  * @tparam SPACE_DIM
598  */
599  template <int SPACE_DIM>
600  using OpMixVecTimesDivLambda =
602 
603  /**
604  * @brief Multiply vector times normal on the face times scalar function
605  *
606  * This operator typically will be used to evaluate natural boundary
607  * conditions for mixed formulation.
608  *
609  * @tparam BASE_DIM
610  * @tparam SPACE_DIM
611  * @tparam OpBase
612  */
613  template <int SPACE_DIM>
614  using OpNormalMixVecTimesScalar =
616 
617  /**
618  * @brief Multiply vector times normal on the face times vector field
619  *
620  * This operator typically will be used to evaluate natural boundary
621  * conditions for mixed formulation.
622  *
623  * @tparam BASE_DIM
624  * @tparam SPACE_DIM
625  * @tparam OpBase
626  */
627  template <int SPACE_DIM>
628  using OpNormalMixVecTimesVectorField =
630 
631  /**
632  * @brief Convective term
633  *
634  * \f[
635  * (v, u_i \mathbf{y}_{,i})
636  * \f]
637  * where \f$\mathbf{y}\$ can be scalar or vector.
638  *
639  * @tparam BASE_DIM
640  * @tparam FIELD_DIM
641  * @tparam SPACE_DIM
642  */
643  template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM>
644  using OpConvectiveTermRhs =
646 };
647 
648 template <typename OpBase>
651  EntitiesFieldData::EntData &row_data) {
653 
654  // get element volume
655  const double vol = OpBase::getMeasure();
656  // get integration weights
657  auto t_w = OpBase::getFTensor0IntegrationWeight();
658  // get base function gradient on rows
659  auto t_row_base = row_data.getFTensor0N();
660  // get coordinate at integration points
661  auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
662  // loop over integration points
663  for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
664  // take into account Jacobian
665  const double alpha =
666  t_w * vol * sourceFun(t_coords(0), t_coords(1), t_coords(2));
667  // loop over rows base functions
668  int rr = 0;
669  for (; rr != OpBase::nbRows; ++rr) {
670  OpBase::locF[rr] += alpha * t_row_base;
671  ++t_row_base;
672  }
673  for (; rr < OpBase::nbRowBaseFunctions; ++rr)
674  ++t_row_base;
675  ++t_coords;
676  ++t_w; // move to another integration weight
677  }
679 }
680 
681 template <int FIELD_DIM, typename OpBase>
684  iNtegrate(EntitiesFieldData::EntData &row_data) {
687 
688  // get element volume
689  const double vol = OpBase::getMeasure();
690  // get integration weights
691  auto t_w = OpBase::getFTensor0IntegrationWeight();
692  // get base function gradient on rows
693  auto t_row_base = row_data.getFTensor0N();
694  // get coordinate at integration points
695  auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
696  // loop over integration points
697  for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
698  // source file
699  auto t_source = sourceFun(t_coords(0), t_coords(1), t_coords(2));
700  // take into account Jacobian
701  const double alpha = t_w * vol;
702  // loop over rows base functions
703  auto t_nf = getFTensor1FromArray<FIELD_DIM, FIELD_DIM>(OpBase::locF);
704  int rr = 0;
705  for (; rr != OpBase::nbRows / FIELD_DIM; ++rr) {
706  t_nf(i) += alpha * t_row_base * t_source(i);
707  ++t_row_base;
708  ++t_nf;
709  }
710  for (; rr < OpBase::nbRowBaseFunctions; ++rr)
711  ++t_row_base;
712  ++t_coords;
713  ++t_w; // move to another integration weight
714  }
716 }
717 
718 template <int FIELD_DIM, typename OpBase>
721  iNtegrate(EntitiesFieldData::EntData &row_data) {
724 
725  const size_t nb_base_functions = row_data.getN().size2() / 3;
726  // get element volume
727  const double vol = OpBase::getMeasure();
728  // get integration weights
729  auto t_w = OpBase::getFTensor0IntegrationWeight();
730  // get base function gradient on rows
731  auto t_row_base = row_data.getFTensor1N<3>();
732  // get coordinate at integration points
733  auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
734  // loop over integration points
735  for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
736  // source file
737  auto t_source = sourceFun(t_coords(0), t_coords(1), t_coords(2));
738  // take into account Jacobian
739  const double alpha = t_w * vol;
740  // loop over rows base functions
741  int rr = 0;
742  for (; rr != OpBase::nbRows; ++rr) {
743  OpBase::locF[rr] += alpha * t_row_base(i) * t_source(i);
744  ++t_row_base;
745  }
746  for (; rr < nb_base_functions; ++rr)
747  ++t_row_base;
748  ++t_coords;
749  ++t_w; // move to another integration weight
750  }
752 }
753 
754 template <int S, typename OpBase>
756  EntitiesFieldData::EntData &row_data) {
758 
759  // get element volume
760  const double vol = OpBase::getMeasure();
761  // get integration weights
762  auto t_w = OpBase::getFTensor0IntegrationWeight();
763  // get base function gradient on rows
764  auto t_row_base = row_data.getFTensor0N();
765  // get vector values
766  auto t_vec = getFTensor0FromVec<S>(*sourceVec);
767 
768 #ifndef NDEBUG
769  if (sourceVec->size() != OpBase::nbIntegrationPts) {
770  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
771  "Wrong number of integration points %d != %d",
772  OpBase::nbIntegrationPts, sourceVec->size());
773  }
774 #endif
775 
776  // get coords
777  auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
778  // loop over integration points
779  for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
780  // take into account Jacobian
781  const double alpha =
782  t_w * vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
783  // loop over rows base functions
784  int rr = 0;
785  for (; rr != OpBase::nbRows; ++rr) {
786  OpBase::locF[rr] += alpha * t_row_base * t_vec;
787  ++t_row_base;
788  }
789  for (; rr < OpBase::nbRowBaseFunctions; ++rr)
790  ++t_row_base;
791  ++t_w; // move to another integration weight
792  ++t_vec;
793  ++t_coords;
794  }
796 }
797 
798 template <int FIELD_DIM, int S, typename OpBase>
800  EntitiesFieldData::EntData &row_data) {
802 
803  // get element volume
804  const double vol = OpBase::getMeasure();
805  // get integration weights
806  auto t_w = OpBase::getFTensor0IntegrationWeight();
807  // get base function gradient on rows
808  auto t_row_base = row_data.getFTensor0N();
809  // get coords
810  auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
811  // get vector values
812  auto t_vec = getFTensor1FromMat<FIELD_DIM, S>(*sourceVec);
813  // loop over integration points
814  for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
815  // take into account Jacobian
816  const double alpha =
817  t_w * vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
818  // get loc vector tensor
819  auto t_nf = OpBase::template getNf<FIELD_DIM>();
820  // loop over rows base functions
821  int rr = 0;
822  for (; rr != OpBase::nbRows / FIELD_DIM; ++rr) {
823  t_nf(i) += alpha * t_row_base * t_vec(i);
824  ++t_row_base;
825  ++t_nf;
826  }
827  for (; rr < OpBase::nbRowBaseFunctions; ++rr)
828  ++t_row_base;
829  ++t_w; // move to another integration weight
830  ++t_vec;
831  ++t_coords;
832  }
834 }
835 
836 template <int FIELD_DIM, int S, typename OpBase>
838  EntitiesFieldData::EntData &row_data) {
840 
841  const size_t nb_base_functions = row_data.getN().size2() / 3;
842  // get element volume
843  const double vol = OpBase::getMeasure();
844  // get integration weights
845  auto t_w = OpBase::getFTensor0IntegrationWeight();
846  // get base function gradient on rows
847  auto t_row_base = row_data.getFTensor1N<3>();
848  // get coords
849  auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
850  // get vector values
851  auto t_vec = getFTensor1FromMat<FIELD_DIM, S>(*sourceVec);
852  // loop over integration points
853  for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
854  // take into account Jacobian
855  const double alpha =
856  t_w * vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
857  // loop over rows base functions
858  int rr = 0;
859  for (; rr != OpBase::nbRows; ++rr) {
860  OpBase::locF[rr] += alpha * t_row_base(i) * t_vec(i);
861  ++t_row_base;
862  }
863  for (; rr < nb_base_functions; ++rr)
864  ++t_row_base;
865  ++t_w; // move to another integration weight
866  ++t_vec;
867  ++t_coords;
868  }
870 }
871 
872 template <int SPACE_DIM, int S, typename OpBase>
875  EntitiesFieldData::EntData &row_data) {
877 
878  // get element volume
879  const double vol = OpBase::getMeasure();
880  // get integration weights
881  auto t_w = OpBase::getFTensor0IntegrationWeight();
882  // get base function gradient on rows
883  auto t_row_grad = row_data.getFTensor1DiffN<SPACE_DIM>();
884  // get filed gradient values
885  auto t_val_grad = getFTensor1FromMat<SPACE_DIM, S>(*(matVals));
886  // get coordinate at integration points
887  auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
888  // loop over integration points
889  for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
890  const double beta = vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
891  // take into account Jacobian
892  const double alpha = t_w * beta;
893  // loop over rows base functions
894  int rr = 0;
895  for (; rr != OpBase::nbRows; rr++) {
896  // calculate element of local matrix
897  OpBase::locF[rr] += alpha * (t_row_grad(i) * t_val_grad(i));
898  ++t_row_grad; // move to another element of gradient of base
899  // function on row
900  }
901  for (; rr < OpBase::nbRowBaseFunctions; ++rr)
902  ++t_row_grad;
903 
904  ++t_coords;
905  ++t_val_grad;
906  ++t_w; // move to another integration weight
907  }
909 }
910 
911 template <int SPACE_DIM, int S, typename OpBase>
914  EntitiesFieldData::EntData &row_data) {
916 
917  // get element volume
918  const double vol = OpBase::getMeasure();
919  // get integration weights
920  auto t_w = OpBase::getFTensor0IntegrationWeight();
921  // get base function gradient on rows
922  auto t_row_grad = row_data.getFTensor1DiffN<SPACE_DIM>();
923  // get filed gradient values
924  auto t_val_grad = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(matVals));
925  // get coordinate at integration points
926  auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
927  // loop over integration points
928  for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
929  // take into account Jacobian
930  const double alpha =
931  t_w * vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
932  // get rhs vector
933  auto t_nf = OpBase::template getNf<SPACE_DIM>();
934  // loop over rows base functions
935  int rr = 0;
936  for (; rr != OpBase::nbRows / SPACE_DIM; rr++) {
937  // calculate element of local matrix
938  t_nf(i) += alpha * (t_row_grad(j) * t_val_grad(i, j));
939  ++t_row_grad; // move to another element of gradient of base
940  // function on row
941  ++t_nf;
942  }
943  for (; rr < OpBase::nbRowBaseFunctions; ++rr)
944  ++t_row_grad;
945 
946  ++t_coords;
947  ++t_val_grad;
948  ++t_w; // move to another integration weight
949  }
951 }
952 
953 template <int SPACE_DIM, int S, typename OpBase>
956  EntitiesFieldData::EntData &row_data) {
958  // get element volume
959  const double vol = OpBase::getMeasure();
960  // get coords
961  auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
962  // get integration weights
963  auto t_w = OpBase::getFTensor0IntegrationWeight();
964  // get base function gradient on rows
965  auto t_row_grad = row_data.getFTensor1DiffN<SPACE_DIM>();
966  // get filed gradient values
967  auto t_val_mat = getFTensor2SymmetricFromMat<SPACE_DIM, S>(*(matVals));
968  // loop over integration points
969  for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
970  // take into account Jacobian
971  const double alpha =
972  t_w * vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
973  // get rhs vector
974  auto t_nf = OpBase::template getNf<SPACE_DIM>();
975  // loop over rows base functions
976  int rr = 0;
977  for (; rr != OpBase::nbRows / SPACE_DIM; rr++) {
978  // calculate element of local matrix
979  t_nf(j) += alpha * (t_row_grad(i) * t_val_mat(i, j));
980  ++t_row_grad; // move to another element of gradient of base
981  // function on row
982  ++t_nf;
983  }
984  for (; rr < OpBase::nbRowBaseFunctions; ++rr)
985  ++t_row_grad;
986  ++t_val_mat;
987  ++t_coords;
988  ++t_w; // move to another integration weight
989  }
991 }
992 
993 template <int FIELD_DIM, int SPACE_DIM, typename OpBase,
994  CoordinateTypes CoordSys>
997  EntitiesFieldData::EntData &row_data) {
999 
1000  const size_t nb_base_functions = row_data.getN().size2() / 3;
1001  auto t_w = this->getFTensor0IntegrationWeight();
1002  // get coordinate at integration points
1003  auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
1004  auto t_diff_base = row_data.getFTensor2DiffN<3, SPACE_DIM>();
1005  auto t_base = row_data.getFTensor1N<3>();
1006  auto t_u = getFTensor1FromMat<FIELD_DIM>(*(matVals));
1007 
1008  for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
1009 
1010  const double alpha = this->getMeasure() * t_w *
1011  betaConst(t_coords(0), t_coords(1), t_coords(2));
1012  auto t_nf = OpBase::template getNf<FIELD_DIM>();
1013 
1014  size_t bb = 0;
1015  for (; bb != this->nbRows / FIELD_DIM; ++bb) {
1016  const double t_div_base = t_diff_base(j, j);
1017  t_nf(i) += alpha * t_div_base * t_u(i);
1018  if constexpr (CoordSys == CYLINDRICAL) {
1019  t_nf(i) += t_base(0) * (alpha / t_coords(0)) * t_u(i);
1020  }
1021  ++t_nf;
1022  ++t_diff_base;
1023  ++t_base;
1024  }
1025  for (; bb < nb_base_functions; ++bb) {
1026  ++t_diff_base;
1027  ++t_base;
1028  }
1029 
1030  ++t_u;
1031  ++t_w;
1032  ++t_coords;
1033  }
1034 
1036 }
1037 
1038 template <int SPACE_DIM, typename OpBase, CoordinateTypes CoordSys>
1041  EntitiesFieldData::EntData &row_data) {
1043 
1044  const size_t nb_base_functions = row_data.getN().size2() / 3;
1045  auto t_w = this->getFTensor0IntegrationWeight();
1046  // get coordinate at integration points
1047  auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
1048  auto t_diff_base = row_data.getFTensor2DiffN<3, SPACE_DIM>();
1049  auto t_u = getFTensor0FromVec(*(vecVals));
1050 
1051  for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
1052 
1053  const double alpha = this->getMeasure() * t_w *
1054  betaConst(t_coords(0), t_coords(1), t_coords(2));
1055  ;
1056 
1057  size_t bb = 0;
1058  for (; bb != this->nbRows; ++bb) {
1059  const double t_div_base = t_diff_base(j, j);
1060  OpBase::locF[bb] += alpha * t_div_base * t_u;
1061  ++t_diff_base;
1062  }
1063  for (; bb < nb_base_functions; ++bb)
1064  ++t_diff_base;
1065 
1066  ++t_u;
1067  ++t_w;
1068  ++t_coords;
1069  }
1070 
1072 }
1073 
1074 /**
1075  * @brief div U times vector
1076  *
1077  * \f[
1078  * \delta u_j = \phi^m\delta\overline{u}^m_j\\
1079  * \delta u_{j,i} = \phi^m_{,i}\delta\overline{u}^m_j\\
1080  * \textrm{tr}[\delta u_{j,i}] = \delta u_{j,i}\delta_{ji}\\
1081  * (\textrm{tr}[\delta u_{j,i}], v) =\\
1082  * (\delta u_{j,i} \delta_{ij}, v) =\\
1083  * (\delta u_{j,i}, \delta_{ij} v) =\\
1084  * (\phi^m_{,i}\delta\overline{u}^m_j, \delta_{ij} v) \\
1085  * f_i^m=(\phi^m_{,i}, v)
1086  * \f]
1087  *
1088  * @tparam FIELD_DIM
1089  * @tparam SPACE_DIM
1090  * @tparam OpBase
1091  * @param row_data
1092  * @return MoFEMErrorCode
1093  */
1094 template <int FIELD_DIM, typename OpBase, CoordinateTypes CoordSys>
1097  EntitiesFieldData::EntData &row_data) {
1100 
1101  // get element volume
1102  const double vol = OpBase::getMeasure();
1103  // get integration weights
1104  auto t_w = OpBase::getFTensor0IntegrationWeight();
1105  // get coordinate at integration points
1106  auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
1107  // get base function gradient on rows
1108  auto t_row_grad = row_data.getFTensor1DiffN<FIELD_DIM>();
1109  // get vector values
1110  auto t_vec = getFTensor0FromVec(*sourceVec);
1111  // loop over integration points
1112  for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
1113  // take into account Jacobian
1114  const double alpha =
1115  t_w * vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
1116  auto t_nf = getFTensor1FromArray<FIELD_DIM, FIELD_DIM>(OpBase::locF);
1117  // loop over rows base functions
1118  int rr = 0;
1119  for (; rr != OpBase::nbRows / FIELD_DIM; ++rr) {
1120  t_nf(i) += alpha * t_row_grad(i) * t_vec;
1121  ++t_row_grad;
1122  ++t_nf;
1123  }
1124  for (; rr < OpBase::nbRowBaseFunctions; ++rr)
1125  ++t_row_grad;
1126  ++t_w; // move to another integration weight
1127  ++t_vec;
1128  ++t_coords;
1129  }
1131 }
1132 
1133 template <int SPACE_DIM, typename OpBase>
1135  EntitiesFieldData::EntData &row_data) {
1137 
1138  const size_t nb_base_functions = row_data.getN().size2() / 3;
1139  auto t_w = this->getFTensor0IntegrationWeight();
1140  auto t_coords = this->getFTensor1CoordsAtGaussPts();
1141  auto t_base = row_data.getFTensor1N<3>();
1142  auto t_grad = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(matVals));
1143 
1144  for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
1145 
1146  const double alpha = this->getMeasure() * t_w;
1147  auto t_nf = OpBase::template getNf<SPACE_DIM>();
1148 
1149  size_t bb = 0;
1150  for (; bb != this->nbRows / SPACE_DIM; ++bb) {
1151  t_nf(i) += alpha * betaCoeff(t_coords(0), t_coords(1), t_coords(2)) *
1152  t_base(j) * t_grad(i, j);
1153  ++t_nf;
1154  ++t_base;
1155  }
1156  for (; bb < nb_base_functions; ++bb)
1157  ++t_base;
1158 
1159  ++t_grad;
1160  ++t_coords;
1161  ++t_w;
1162  }
1163 
1165 }
1166 
1167 template <int SPACE_DIM, typename OpBase>
1169  EntitiesFieldData::EntData &row_data) {
1171 
1172  const size_t nb_base_functions = row_data.getN().size2();
1173  auto t_w = this->getFTensor0IntegrationWeight();
1174  auto t_coords = this->getFTensor1CoordsAtGaussPts();
1175  auto t_base = row_data.getFTensor0N();
1176  auto t_div = getFTensor1FromMat<SPACE_DIM>(*(matVals));
1177 
1178  for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
1179  const double alpha = this->getMeasure() * t_w;
1180  auto t_nf = OpBase::template getNf<SPACE_DIM>();
1181 
1182  size_t bb = 0;
1183  for (; bb != this->nbRows / SPACE_DIM; ++bb) {
1184  t_nf(i) += alpha * t_base *
1185  betaCoeff(t_coords(0), t_coords(1), t_coords(2)) * t_div(i);
1186  ++t_nf;
1187  ++t_base;
1188  }
1189  for (; bb < nb_base_functions; ++bb)
1190  ++t_base;
1191 
1192  ++t_coords;
1193  ++t_div;
1194  ++t_w;
1195  }
1196 
1198 }
1199 
1200 template <typename OpBase>
1202  EntitiesFieldData::EntData &row_data) {
1204 
1205  const size_t nb_base_functions = row_data.getN().size2() / 3;
1206  // get element volume
1207  // get integration weights
1208  auto t_w = OpBase::getFTensor0IntegrationWeight();
1209  // get base function gradient on rows
1210  auto t_row_base = row_data.getFTensor1N<3>();
1211  // get coordinate at integration points
1212  auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
1213  // get normal
1214  auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
1215  double a = 1;
1216  if (this->getNumeredEntFiniteElementPtr()->getEntType() == MBTRI)
1217  a *= 2;
1218  // loop over integration points
1219  for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
1220  // take into account Jacobian
1221  const double alpha =
1222  t_w * sourceFun(t_coords(0), t_coords(1), t_coords(2)) / a;
1223  // loop over rows base functions
1224  int rr = 0;
1225  for (; rr != OpBase::nbRows; ++rr) {
1226  OpBase::locF[rr] += alpha * t_row_base(i) * t_normal(i);
1227  ++t_row_base;
1228  }
1229  for (; rr < nb_base_functions; ++rr)
1230  ++t_row_base;
1231  ++t_coords;
1232  ++t_w; // move to another integration weight
1233  ++t_normal;
1234  }
1236 }
1237 
1238 template <typename OpBase>
1240  EntitiesFieldData::EntData &row_data) {
1242 
1243  const size_t nb_base_functions = row_data.getN().size2() / 3;
1244  FTensor::Tensor1<double, 3> t_z{0., 0., 1.};
1245  // get element volume
1246  // get integration weights
1247  auto t_w = OpBase::getFTensor0IntegrationWeight();
1248  // get base function gradient on rows
1249  auto t_row_base = row_data.getFTensor1N<3>();
1250  // get coordinate at integration points
1251  auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
1252  // get normal
1253  auto t_tangent = OpBase::getFTensor1TangentAtGaussPts();
1254  // loop over integration points
1255  for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
1256  // take into account Jacobian
1257  const double alpha = t_w * sourceFun(t_coords(0), t_coords(1), t_coords(2));
1258  FTensor::Tensor1<double, 3> t_normal;
1261  t_normal(i) = FTensor::levi_civita(i, j, k) * t_tangent(j) * t_z(k);
1262  int rr = 0;
1263  for (; rr != OpBase::nbRows; ++rr) {
1264  OpBase::locF[rr] += alpha * t_row_base(i) * t_normal(i);
1265  ++t_row_base;
1266  }
1267  for (; rr < nb_base_functions; ++rr)
1268  ++t_row_base;
1269  ++t_coords;
1270  ++t_tangent;
1271  ++t_w; // move to another integration weight
1272  }
1274 }
1275 
1276 template <int SPACE_DIM, typename OpBase>
1279  EntitiesFieldData::EntData &row_data) {
1283 
1284  const size_t nb_base_functions = row_data.getN().size2() / 3;
1285  // get element volume
1286  // get integration weights
1287  auto t_w = OpBase::getFTensor0IntegrationWeight();
1288  // get base function gradient on rows
1289  auto t_row_base = row_data.getFTensor1N<3>();
1290  // get coordinate at integration points
1291  auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
1292  // get normal
1293  auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
1294  // get field
1295  auto t_u = getFTensor1FromMat<SPACE_DIM>(*uPtr);
1296  // loop over integration points
1297  auto a = OpBase::getMeasure();
1298  for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
1299  // take into account Jacobian
1300  auto l2 = std::sqrt(t_normal(i) * t_normal(i));
1301  const double alpha =
1302  t_w * betaCoeff(t_coords(0), t_coords(1), t_coords(2)) * (a / l2);
1303  // get rhs vector
1304  auto t_nf = OpBase::template getNf<SPACE_DIM>();
1305  // loop over rows base functions
1306  int rr = 0;
1307  for (; rr != OpBase::nbRows / SPACE_DIM; ++rr) {
1308  t_nf(J) += alpha * (t_row_base(i) * t_normal(i)) * t_u(J);
1309  ++t_row_base;
1310  ++t_nf;
1311  }
1312  for (; rr < nb_base_functions; ++rr)
1313  ++t_row_base;
1314  ++t_coords;
1315  ++t_w; // move to another integration weight
1316  ++t_normal;
1317  ++t_u;
1318  }
1320 }
1321 
1322 template <int SPACE_DIM, typename OpBase>
1325  EntitiesFieldData::EntData &row_data) {
1327 
1328  auto t_w = this->getFTensor0IntegrationWeight();
1329  auto t_base = row_data.getFTensor0N();
1330 
1331  auto t_u = getFTensor1FromMat<SPACE_DIM>(*uPtr);
1332  auto t_grad_y = getFTensor1FromMat<SPACE_DIM>(*yGradPtr);
1333 
1335  const double alpha_constant = alphaConstant();
1336  for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
1337 
1338  // get element volume
1339  const double vol = OpBase::getMeasure();
1340  const double c = (t_grad_y(i) * t_u(i)) * (t_w * vol * alpha_constant);
1341 
1342  // get element volume
1343  int rr = 0;
1344  for (; rr != OpBase::nbRows; ++rr) {
1345  OpBase::locF[rr] += c * t_base;
1346  ++t_base;
1347  }
1348  for (; rr < OpBase::nbRowBaseFunctions; ++rr)
1349  ++t_base;
1350 
1351  ++t_w; // move to another integration weight
1352  ++t_u;
1353  ++t_grad_y;
1354  }
1355 
1357 }
1358 
1359 template <int FIELD_DIM, int SPACE_DIM, typename OpBase>
1362  EntitiesFieldData::EntData &row_data) {
1364 
1365  auto t_w = this->getFTensor0IntegrationWeight();
1366  auto t_base = row_data.getFTensor0N();
1367 
1368  auto t_u = getFTensor1FromMat<SPACE_DIM>(*uPtr);
1369  auto t_grad_y = getFTensor2FromMat<FIELD_DIM, SPACE_DIM>(*yGradPtr);
1370 
1373  const double alpha_constant = alphaConstant();
1374  for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
1375 
1376  // get element volume
1377  const double vol = OpBase::getMeasure();
1378 
1380  t_c(J) = (t_grad_y(J, i) * t_u(i)) * (t_w * vol * alpha_constant);
1381 
1382  auto t_nf = getFTensor1FromArray<FIELD_DIM, FIELD_DIM>(OpBase::locF);
1383  int rr = 0;
1384  for (; rr != OpBase::nbRows / FIELD_DIM; ++rr) {
1385  t_nf(J) += t_base * t_c(J);
1386  ++t_base;
1387  ++t_nf;
1388  }
1389  for (; rr < OpBase::nbRowBaseFunctions; ++rr)
1390  ++t_base;
1391  ++t_w; // move to another integration weight
1392 
1393  ++t_u;
1394  ++t_grad_y;
1395  }
1396 
1398 }
1399 
1400 } // namespace MoFEM
1401 
1402 #endif // __LINEAR_FORMS_INTEGRATORS_HPP__
MoFEM::OpConvectiveTermRhsImpl< 1, 1, SPACE_DIM, GAUSS, OpBase >::alphaConstant
ConstantFun alphaConstant
Definition: LinearFormsIntegrators.hpp:469
MoFEM::OpMixDivTimesUImpl< 3, 1, SPACE_DIM, GAUSS, OpBase, CoordSys >::vecVals
boost::shared_ptr< VectorDouble > vecVals
Definition: LinearFormsIntegrators.hpp:276
MoFEM::OpBaseImpl::nbIntegrationPts
int nbIntegrationPts
number of integration points
Definition: FormsIntegrators.hpp:228
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:127
MoFEM::OpMixDivTimesUImpl< 3, FIELD_DIM, SPACE_DIM, GAUSS, OpBase, CoordSys >::j
FTensor::Index< 'j', SPACE_DIM > j
Definition: LinearFormsIntegrators.hpp:260
MoFEM::SourceFunctionSpecialization::S
Definition: LinearFormsIntegrators.hpp:13
MoFEM::OpNormalMixVecTimesScalarImpl< 3, GAUSS, OpBase >::i
FTensor::Index< 'i', 3 > i
Definition: LinearFormsIntegrators.hpp:405
MoFEM::OpNormalMixVecTimesScalarImpl
Multiply vector times normal on the face times scalar function.
Definition: LinearFormsIntegrators.hpp:378
FTensor::Tensor1< double, 3 >
MoFEM::OpMixDivTimesUImpl< 3, FIELD_DIM, SPACE_DIM, GAUSS, OpBase, CoordSys >::i
FTensor::Index< 'i', FIELD_DIM > i
Definition: LinearFormsIntegrators.hpp:259
CARTESIAN
@ CARTESIAN
Definition: definitions.h:115
MoFEM::OpGradTimesTensorImpl< 1, 1, SPACE_DIM, S, GAUSS, OpBase >::betaCoeff
ScalarFun betaCoeff
Definition: LinearFormsIntegrators.hpp:196
MoFEM::OpGradTimesTensorImpl< 1, 1, SPACE_DIM, S, GAUSS, OpBase >::matVals
boost::shared_ptr< MatrixDouble > matVals
Definition: LinearFormsIntegrators.hpp:194
MoFEM::OpMixDivTimesUImpl< 3, 1, SPACE_DIM, GAUSS, OpBase, CoordSys >::OpMixDivTimesUImpl
OpMixDivTimesUImpl(const std::string field_name, boost::shared_ptr< VectorDouble > vec_vals, ScalarFun beta=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
Definition: LinearFormsIntegrators.hpp:267
MoFEM::OpMixVecTimesDivLambdaImpl< SPACE_DIM, GAUSS, OpBase >::OpMixVecTimesDivLambdaImpl
OpMixVecTimesDivLambdaImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals)
Definition: LinearFormsIntegrators.hpp:314
MoFEM::OpSourceImpl< 1, 1, GAUSS, SourceFunctionSpecialization::S< OpBase > >::sourceFun
ScalarFun sourceFun
Definition: LinearFormsIntegrators.hpp:61
MoFEM::OpNormalMixVecTimesScalarImpl< 3, GAUSS, OpBase >::sourceFun
ScalarFun sourceFun
Definition: LinearFormsIntegrators.hpp:404
MoFEM::OpNormalMixVecTimesScalarImpl< 3, GAUSS, OpBase >::OpNormalMixVecTimesScalarImpl
OpNormalMixVecTimesScalarImpl(const std::string field_name, ScalarFun source_fun=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
Definition: LinearFormsIntegrators.hpp:396
MoFEM::OpNormalMixVecTimesScalarImpl< 2, GAUSS, OpBase >::i
FTensor::Index< 'i', 3 > i
Definition: LinearFormsIntegrators.hpp:420
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
OpConvectiveTermRhs
FormsIntegrators< DomainEleOp >::Assembly< A >::LinearForm< I >::OpConvectiveTermRhs< 1, FIELD_DIM, SPACE_DIM > OpConvectiveTermRhs
Definition: operators_tests.cpp:52
MoFEM::OpBaseTimesVectorImpl< 3, FIELD_DIM, S, GAUSS, OpBase >::sourceVec
boost::shared_ptr< MatrixDouble > sourceVec
Definition: LinearFormsIntegrators.hpp:171
MoFEM::SourceFunctionSpecialization
Definition: LinearFormsIntegrators.hpp:12
MoFEM::SourceBoundaryNormalSpecialization
Definition: LinearFormsIntegrators.hpp:17
MoFEM::ConstantFun
boost::function< double()> ConstantFun
Constant function type.
Definition: FormsIntegrators.hpp:158
MoFEM::SourceBoundaryNormalSpecialization::S::S
S()=delete
MoFEM::OpNormalMixVecTimesVectorFieldImpl< SPACE_DIM, GAUSS, OpBase >::betaCoeff
ScalarFun betaCoeff
Definition: LinearFormsIntegrators.hpp:449
MoFEM::OpMixDivTimesUImpl< 3, FIELD_DIM, SPACE_DIM, GAUSS, OpBase, CoordSys >::matVals
boost::shared_ptr< MatrixDouble > matVals
Definition: LinearFormsIntegrators.hpp:258
MoFEM::OpGradTimesTensorImpl< 1, SPACE_DIM, SPACE_DIM, S, GAUSS, OpBase >::OpGradTimesTensorImpl
OpGradTimesTensorImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals, ScalarFun beta_coeff=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
Definition: LinearFormsIntegrators.hpp:206
J
FTensor::Index< 'J', DIM1 > J
Definition: level_set.cpp:30
MoFEM::OpNormalMixVecTimesVectorFieldImpl< SPACE_DIM, GAUSS, OpBase >::OpNormalMixVecTimesVectorFieldImpl
OpNormalMixVecTimesVectorFieldImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > u_ptr, ScalarFun beta_coeff=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
Definition: LinearFormsIntegrators.hpp:440
MoFEM::OpMixTensorTimesGradUImpl< SPACE_DIM, GAUSS, OpBase >::matVals
boost::shared_ptr< MatrixDouble > matVals
Definition: LinearFormsIntegrators.hpp:359
MoFEM::OpGradTimesTensorImpl
Definition: LinearFormsIntegrators.hpp:178
OpBase
OpBaseImpl< PETSC, EdgeEleOp > OpBase
Definition: radiation.cpp:29
FTensor::levi_civita
constexpr std::enable_if<(Dim0<=2 &&Dim1<=2), Tensor2_Expr< Levi_Civita< T >, T, Dim0, Dim1, i, j > >::type levi_civita(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
levi_civita functions to make for easy adhoc use
Definition: Levi_Civita.hpp:617
BASE_DIM
constexpr int BASE_DIM
Definition: dg_projection.cpp:15
MoFEM::OpBaseTimesVectorImpl< 3, FIELD_DIM, S, GAUSS, OpBase >::i
FTensor::Index< 'i', FIELD_DIM > i
Definition: LinearFormsIntegrators.hpp:172
CoordinateTypes
CoordinateTypes
Coodinate system.
Definition: definitions.h:114
MoFEM::OpSourceImpl< 1, FIELD_DIM, GAUSS, SourceFunctionSpecialization::S< OpBase > >::OpSourceImpl
OpSourceImpl(const std::string field_name, TimeFun time_fun, VectorFun< FIELD_DIM > source_fun, boost::shared_ptr< Range > ents_ptr=nullptr)
Construct a new Op Source Impl object.
Definition: LinearFormsIntegrators.hpp:77
MoFEM::SourceFunctionSpecialization::SourceFunctionSpecialization
SourceFunctionSpecialization()=delete
MoFEM::OpBaseTimesVectorImpl
Definition: LinearFormsIntegrators.hpp:140
FIELD_DIM
constexpr int FIELD_DIM
Definition: child_and_parent.cpp:15
MoFEM::SourceBoundaryNormalSpecialization::SourceBoundaryNormalSpecialization
SourceBoundaryNormalSpecialization()=delete
I
constexpr IntegrationType I
Definition: operators_tests.cpp:31
MoFEM::OpMixDivTimesUImpl< 3, FIELD_DIM, SPACE_DIM, GAUSS, OpBase, CoordSys >::betaConst
ScalarFun betaConst
Definition: LinearFormsIntegrators.hpp:257
MoFEM::OpMixDivTimesUImpl< 1, FIELD_DIM, FIELD_DIM, GAUSS, OpBase, CoordSys >::betaCoeff
ScalarFun betaCoeff
Definition: LinearFormsIntegrators.hpp:293
MoFEM::EntitiesFieldData::EntData::getFTensor0N
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
Definition: EntitiesFieldData.hpp:1489
EleOp
MoFEM::OpGradTimesTensorImpl< 1, 1, SPACE_DIM, S, GAUSS, OpBase >::i
FTensor::Index< 'i', SPACE_DIM > i
summit Index
Definition: LinearFormsIntegrators.hpp:184
MoFEM::OpBaseImpl::locF
VectorDouble locF
local entity vector
Definition: FormsIntegrators.hpp:241
c
const double c
speed of light (cm/ns)
Definition: initial_diffusion.cpp:39
MoFEM::OpMixDivTimesUImpl< 3, 1, SPACE_DIM, GAUSS, OpBase, CoordSys >::betaConst
ScalarFun betaConst
Definition: LinearFormsIntegrators.hpp:275
MoFEM::EntitiesFieldData::EntData::getFTensor2DiffN
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2DiffN(FieldApproximationBase base)
Get derivatives of base functions for Hdiv space.
Definition: EntitiesFieldData.cpp:660
MoFEM::OpBaseImpl
Definition: FormsIntegrators.hpp:170
MoFEM::SourceFunctionSpecialization::S::S
S()=delete
MoFEM::TimeFun
boost::function< double(double)> TimeFun
Lambda function used to scale with time.
Definition: FormsIntegrators.hpp:146
MoFEM::OpSourceImpl< 3, FIELD_DIM, GAUSS, SourceFunctionSpecialization::S< OpBase > >::OpSourceImpl
OpSourceImpl(const std::string field_name, VectorFun< FIELD_DIM > source_fun, boost::shared_ptr< Range > ents_ptr=nullptr)
Definition: LinearFormsIntegrators.hpp:109
MoFEM::OpMixTensorTimesGradUImpl< SPACE_DIM, GAUSS, OpBase >::j
FTensor::Index< 'j', SPACE_DIM > j
Definition: LinearFormsIntegrators.hpp:363
MoFEM::OpMixDivTimesUImpl< 3, 1, SPACE_DIM, GAUSS, OpBase, CoordSys >::j
FTensor::Index< 'j', SPACE_DIM > j
Definition: LinearFormsIntegrators.hpp:277
MoFEM::OpGradTimesSymTensorImpl
Definition: LinearFormsIntegrators.hpp:221
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
SPACE_DIM
constexpr int SPACE_DIM
Definition: child_and_parent.cpp:16
MoFEM::VectorFun
boost::function< FTensor::Tensor1< double, DIM >(const double, const double, const double)> VectorFun
Vector function type.
Definition: FormsIntegrators.hpp:168
a
constexpr double a
Definition: approx_sphere.cpp:30
MoFEM::OpConvectiveTermRhsImpl< 1, 1, SPACE_DIM, GAUSS, OpBase >::uPtr
boost::shared_ptr< MatrixDouble > uPtr
Definition: LinearFormsIntegrators.hpp:467
MoFEM::OpGradTimesTensorImpl< 1, SPACE_DIM, SPACE_DIM, S, GAUSS, OpBase >::i
FTensor::Index< 'i', SPACE_DIM > i
summit Index
Definition: LinearFormsIntegrators.hpp:203
MoFEM::OpGradTimesTensorImpl< 1, SPACE_DIM, SPACE_DIM, S, GAUSS, OpBase >::betaCoeff
ScalarFun betaCoeff
Definition: LinearFormsIntegrators.hpp:216
MoFEM::OpBaseTimesVectorImpl< 1, FIELD_DIM, S, GAUSS, OpBase >::i
FTensor::Index< 'i', FIELD_DIM > i
Definition: LinearFormsIntegrators.hpp:155
MoFEM::OpGradTimesSymTensorImpl< 1, SPACE_DIM, SPACE_DIM, S, GAUSS, OpBase >::betaCoeff
ScalarFun betaCoeff
Definition: LinearFormsIntegrators.hpp:238
double
MoFEM::OpBaseTimesVectorImpl< 3, FIELD_DIM, S, GAUSS, OpBase >::OpBaseTimesVectorImpl
OpBaseTimesVectorImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > vec, ScalarFun beta_coeff=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
Definition: LinearFormsIntegrators.hpp:162
MoFEM::OpBaseTimesScalarImpl< 1, S, GAUSS, OpBase >::sourceVec
boost::shared_ptr< VectorDouble > sourceVec
Definition: LinearFormsIntegrators.hpp:134
MoFEM::getFTensor0FromVec
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
MoFEM::OpNormalMixVecTimesScalarImpl< 2, GAUSS, OpBase >::sourceFun
ScalarFun sourceFun
Definition: LinearFormsIntegrators.hpp:419
MoFEM::OpConvectiveTermRhsImpl
Definition: LinearFormsIntegrators.hpp:455
MoFEM::GAUSS
@ GAUSS
Definition: FormsIntegrators.hpp:128
MoFEM::OpMixDivTimesUImpl
Definition: LinearFormsIntegrators.hpp:243
MoFEM::OpSourceImpl< 1, FIELD_DIM, GAUSS, SourceFunctionSpecialization::S< OpBase > >::sourceFun
VectorFun< FIELD_DIM > sourceFun
Definition: LinearFormsIntegrators.hpp:95
MoFEM::OpNormalMixVecTimesVectorFieldImpl< SPACE_DIM, GAUSS, OpBase >::uPtr
boost::shared_ptr< MatrixDouble > uPtr
Definition: LinearFormsIntegrators.hpp:448
MoFEM::ScalarFun
boost::function< double(const double, const double, const double)> ScalarFun
Scalar function type.
Definition: FormsIntegrators.hpp:136
MoFEM::OpConvectiveTermRhsImpl< 1, FIELD_DIM, SPACE_DIM, GAUSS, OpBase >::alphaConstant
ConstantFun alphaConstant
Definition: LinearFormsIntegrators.hpp:486
MoFEM::OpSourceImpl
Definition: LinearFormsIntegrators.hpp:23
MoFEM::OpBaseTimesVectorImpl< 1, FIELD_DIM, S, GAUSS, OpBase >::betaCoeff
ScalarFun betaCoeff
Definition: LinearFormsIntegrators.hpp:153
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
MoFEM::OpBaseTimesVectorImpl< 1, FIELD_DIM, S, GAUSS, OpBase >::sourceVec
boost::shared_ptr< MatrixDouble > sourceVec
Definition: LinearFormsIntegrators.hpp:154
OpGradTimesTensor
FormsIntegrators< DomainEleOp >::Assembly< A >::LinearForm< I >::OpGradTimesTensor< 1, FIELD_DIM, SPACE_DIM > OpGradTimesTensor
Definition: operators_tests.cpp:48
MoFEM::OpMixVecTimesDivLambdaImpl< SPACE_DIM, GAUSS, OpBase >::OpMixVecTimesDivLambdaImpl
OpMixVecTimesDivLambdaImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals, ScalarFun beta_fun)
Definition: LinearFormsIntegrators.hpp:318
MoFEM::OpMixDivTimesUImpl< 1, FIELD_DIM, FIELD_DIM, GAUSS, OpBase, CoordSys >::sourceVec
boost::shared_ptr< VectorDouble > sourceVec
Definition: LinearFormsIntegrators.hpp:294
MoFEM::OpConvectiveTermRhsImpl< 1, FIELD_DIM, SPACE_DIM, GAUSS, OpBase >::yGradPtr
boost::shared_ptr< MatrixDouble > yGradPtr
Definition: LinearFormsIntegrators.hpp:485
MoFEM::OpBaseTimesScalarImpl< 1, S, GAUSS, OpBase >::betaCoeff
ScalarFun betaCoeff
Definition: LinearFormsIntegrators.hpp:133
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
FTensor::Index< 'i', FIELD_DIM >
MoFEM::IntegrationType
IntegrationType
Form integrator integration types.
Definition: FormsIntegrators.hpp:128
MoFEM::OpSourceImpl< 3, FIELD_DIM, GAUSS, SourceFunctionSpecialization::S< OpBase > >::sourceFun
VectorFun< FIELD_DIM > sourceFun
Definition: LinearFormsIntegrators.hpp:115
MoFEM::OpGradTimesSymTensorImpl< 1, SPACE_DIM, SPACE_DIM, S, GAUSS, OpBase >::OpGradTimesSymTensorImpl
OpGradTimesSymTensorImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals, ScalarFun beta_coeff=[](double, double, double) constexpr { return 1;})
Definition: LinearFormsIntegrators.hpp:227
MoFEM::SourceBoundaryNormalSpecialization::S
Definition: LinearFormsIntegrators.hpp:18
MoFEM::OpConvectiveTermRhsImpl< 1, FIELD_DIM, SPACE_DIM, GAUSS, OpBase >::OpConvectiveTermRhsImpl
OpConvectiveTermRhsImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > u_ptr, boost::shared_ptr< MatrixDouble > y_grad_ptr, ConstantFun source_fun=[]() constexpr { return 1;})
Definition: LinearFormsIntegrators.hpp:476
MoFEM::OpMixTensorTimesGradUImpl
Tensor field time gradient of vector field.
Definition: LinearFormsIntegrators.hpp:344
MoFEM::OpBaseImpl::nbRowBaseFunctions
int nbRowBaseFunctions
number or row base functions
Definition: FormsIntegrators.hpp:229
MoFEM::OpGradTimesSymTensorImpl< 1, SPACE_DIM, SPACE_DIM, S, GAUSS, OpBase >::j
FTensor::Index< 'j', SPACE_DIM > j
Definition: LinearFormsIntegrators.hpp:236
MoFEM::OpMixTensorTimesGradUImpl< SPACE_DIM, GAUSS, OpBase >::i
FTensor::Index< 'i', SPACE_DIM > i
Definition: LinearFormsIntegrators.hpp:362
MoFEM::OpConvectiveTermRhsImpl< 1, 1, SPACE_DIM, GAUSS, OpBase >::yGradPtr
boost::shared_ptr< MatrixDouble > yGradPtr
Definition: LinearFormsIntegrators.hpp:468
MoFEM::EntitiesFieldData::EntData::getFTensor1DiffN
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
Definition: EntitiesFieldData.cpp:526
MoFEM::OpBaseImpl::nbRows
int nbRows
number of dofs on rows
Definition: FormsIntegrators.hpp:226
CYLINDRICAL
@ CYLINDRICAL
Definition: definitions.h:117
MoFEM::OpGradTimesTensorImpl< 1, SPACE_DIM, SPACE_DIM, S, GAUSS, OpBase >::j
FTensor::Index< 'j', SPACE_DIM > j
summit Index
Definition: LinearFormsIntegrators.hpp:204
MoFEM::EntitiesFieldData::EntData::getN
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
Definition: EntitiesFieldData.hpp:1305
MoFEM::OpNormalMixVecTimesScalarImpl< 2, GAUSS, OpBase >::OpNormalMixVecTimesScalarImpl
OpNormalMixVecTimesScalarImpl(const std::string field_name, ScalarFun source_fun=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
Definition: LinearFormsIntegrators.hpp:411
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
MoFEM::EntitiesFieldData::EntData::getFTensor1N
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
Definition: EntitiesFieldData.cpp:640
MoFEM::OpGradTimesSymTensorImpl< 1, SPACE_DIM, SPACE_DIM, S, GAUSS, OpBase >::matVals
boost::shared_ptr< MatrixDouble > matVals
Definition: LinearFormsIntegrators.hpp:234
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::OpMixTensorTimesGradUImpl< SPACE_DIM, GAUSS, OpBase >::OpMixTensorTimesGradUImpl
OpMixTensorTimesGradUImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals, ScalarFun beta_fun)
Definition: LinearFormsIntegrators.hpp:352
MoFEM::OpMixDivTimesUImpl< 1, FIELD_DIM, FIELD_DIM, GAUSS, OpBase, CoordSys >::OpMixDivTimesUImpl
OpMixDivTimesUImpl(const std::string field_name, boost::shared_ptr< VectorDouble > vec, ScalarFun beta=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
Definition: LinearFormsIntegrators.hpp:285
MoFEM::OpMixDivTimesUImpl< 3, FIELD_DIM, SPACE_DIM, GAUSS, OpBase, CoordSys >::OpMixDivTimesUImpl
OpMixDivTimesUImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals, ScalarFun beta=[](double, double, double) { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
Definition: LinearFormsIntegrators.hpp:249
MoFEM::OpMixVecTimesDivLambdaImpl< SPACE_DIM, GAUSS, OpBase >::i
FTensor::Index< 'i', SPACE_DIM > i
Definition: LinearFormsIntegrators.hpp:328
MoFEM::FormsIntegrators
Integrator forms.
Definition: FormsIntegrators.hpp:291
MoFEM::OpMixVecTimesDivLambdaImpl
Vector filed time divergence of tensor.
Definition: LinearFormsIntegrators.hpp:310
MoFEM::OpSourceImpl< 3, FIELD_DIM, GAUSS, SourceFunctionSpecialization::S< OpBase > >::OpSourceImpl
OpSourceImpl(const std::string field_name, TimeFun time_fun, VectorFun< FIELD_DIM > source_fun, boost::shared_ptr< Range > ents_ptr=nullptr)
Definition: LinearFormsIntegrators.hpp:103
MoFEM::OpSourceImpl< 1, 1, GAUSS, SourceFunctionSpecialization::S< OpBase > >::OpSourceImpl
OpSourceImpl(const std::string field_name, TimeFun time_fun, ScalarFun source_fun, boost::shared_ptr< Range > ents_ptr=nullptr)
Construct a new Op Source Impl object.
Definition: LinearFormsIntegrators.hpp:54
MoFEM::OpBaseTimesVectorImpl< 1, FIELD_DIM, S, GAUSS, OpBase >::OpBaseTimesVectorImpl
OpBaseTimesVectorImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > vec, ScalarFun beta_coeff=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
Definition: LinearFormsIntegrators.hpp:145
MoFEM::OpMixTensorTimesGradUImpl< SPACE_DIM, GAUSS, OpBase >::OpMixTensorTimesGradUImpl
OpMixTensorTimesGradUImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals)
Definition: LinearFormsIntegrators.hpp:348
MoFEM::OpBaseTimesVectorImpl< 3, FIELD_DIM, S, GAUSS, OpBase >::betaCoeff
ScalarFun betaCoeff
Definition: LinearFormsIntegrators.hpp:170
MoFEM::OpBaseTimesScalarImpl
Definition: LinearFormsIntegrators.hpp:120
MoFEM::OpGradTimesSymTensorImpl< 1, SPACE_DIM, SPACE_DIM, S, GAUSS, OpBase >::i
FTensor::Index< 'i', SPACE_DIM > i
Definition: LinearFormsIntegrators.hpp:235
k
FTensor::Index< 'k', 3 > k
Definition: matrix_function.cpp:20
MoFEM::OpGradTimesTensorImpl< 1, SPACE_DIM, SPACE_DIM, S, GAUSS, OpBase >::matVals
boost::shared_ptr< MatrixDouble > matVals
Definition: LinearFormsIntegrators.hpp:214
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEM::OpConvectiveTermRhsImpl< 1, 1, SPACE_DIM, GAUSS, OpBase >::OpConvectiveTermRhsImpl
OpConvectiveTermRhsImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > u_ptr, boost::shared_ptr< MatrixDouble > y_grad_ptr, ConstantFun source_fun=[]() constexpr { return 1;})
Definition: LinearFormsIntegrators.hpp:459
MoFEM::OpConvectiveTermRhsImpl< 1, FIELD_DIM, SPACE_DIM, GAUSS, OpBase >::uPtr
boost::shared_ptr< MatrixDouble > uPtr
Definition: LinearFormsIntegrators.hpp:484
MoFEM::OpSourceImpl< 1, FIELD_DIM, GAUSS, SourceFunctionSpecialization::S< OpBase > >::OpSourceImpl
OpSourceImpl(const std::string field_name, VectorFun< FIELD_DIM > source_fun, boost::shared_ptr< Range > ents_ptr=nullptr)
Construct a new Op Source Impl object.
Definition: LinearFormsIntegrators.hpp:89
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
MoFEM::OpNormalMixVecTimesVectorFieldImpl
Multiply vector times normal on the face times vector field.
Definition: LinearFormsIntegrators.hpp:435
MoFEM::OpBaseTimesScalarImpl< 1, S, GAUSS, OpBase >::OpBaseTimesScalarImpl
OpBaseTimesScalarImpl(const std::string field_name, boost::shared_ptr< VectorDouble > vec, ScalarFun beta_coeff=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
Definition: LinearFormsIntegrators.hpp:125
MoFEM::OpSourceImpl< 1, 1, GAUSS, SourceFunctionSpecialization::S< OpBase > >::OpSourceImpl
OpSourceImpl(const std::string field_name, ScalarFun source_fun, boost::shared_ptr< Range > ents_ptr=nullptr)
Construct a new Op Source Impl object.
Definition: LinearFormsIntegrators.hpp:41
MoFEM::OpGradTimesTensorImpl< 1, 1, SPACE_DIM, S, GAUSS, OpBase >::OpGradTimesTensorImpl
OpGradTimesTensorImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals, ScalarFun beta_coeff=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
Definition: LinearFormsIntegrators.hpp:186
MoFEM::OpMixVecTimesDivLambdaImpl< SPACE_DIM, GAUSS, OpBase >::matVals
boost::shared_ptr< MatrixDouble > matVals
Definition: LinearFormsIntegrators.hpp:325