v0.15.0
Loading...
Searching...
No Matches
LinearFormsIntegratorsImpl.hpp
Go to the documentation of this file.
1/** \file LinearFormsIntegratorsImpl.hpp
2 * \brief Linear forms integrators (implementation)
3 * \ingroup mofem_form
4
5*/
6
7#ifndef __LINEAR_FORMS_INTEGRATORS_IMPL_HPP__
8#define __LINEAR_FORMS_INTEGRATORS_IMPL_HPP__
9
10namespace MoFEM {
11
13 template <typename OpBase> struct S { S() = delete; };
15};
16
18 template <typename OpBase> struct S { S() = delete; };
20};
21
22template <int BASE_DIM, int FIELD_DIM, IntegrationType I, typename OpBase>
24
25/**
26 * @brief Integrate source
27 *
28 * @tparam OpBase
29 */
30template <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
60protected:
63};
64
65template <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 sourceFun(source_fun) {}
82
83 /**
84 * @brief Construct a new Op Source Impl object
85 *
86 * @param field_name
87 * @param source_fun
88 * @param ents_ptr
89 */
90 OpSourceImpl(const std::string field_name, VectorFun<FIELD_DIM> source_fun,
91 boost::shared_ptr<Range> ents_ptr = nullptr)
92 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
93 sourceFun(source_fun) {}
94
95protected:
98};
99
100template <int FIELD_DIM, typename OpBase>
102 SourceFunctionSpecialization::S<OpBase>> : public OpBase {
103
104 OpSourceImpl(const std::string field_name, TimeFun time_fun,
105 VectorFun<FIELD_DIM> source_fun,
106 boost::shared_ptr<Range> ents_ptr = nullptr)
107 : OpBase(field_name, field_name, OpBase::OPROW, time_fun, ents_ptr),
108 sourceFun(source_fun) {}
109
110 OpSourceImpl(const std::string field_name, VectorFun<FIELD_DIM> source_fun,
111 boost::shared_ptr<Range> ents_ptr = nullptr)
112 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
113 sourceFun(source_fun) {}
114
115protected:
118};
119
120template <int BASE_DIM, int S, IntegrationType I, typename OpBase>
122
123template <int S, typename OpBase>
124struct OpBaseTimesScalarImpl<1, S, GAUSS, OpBase> : public OpBase {
125
127 const std::string field_name, boost::shared_ptr<VectorDouble> vec,
128 ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; },
129 boost::shared_ptr<Range> ents_ptr = nullptr)
130 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr), sourceVec(vec),
131 betaCoeff(beta_coeff) {}
132
133protected:
135 boost::shared_ptr<VectorDouble> sourceVec;
137};
138
139template <int BASE_DIM, int FIELD_DIM, int S, IntegrationType I,
140 typename OpBase>
142
143template <int FIELD_DIM, int S, typename OpBase>
145
147 const std::string field_name, boost::shared_ptr<MatrixDouble> vec,
148 ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; },
149 boost::shared_ptr<Range> ents_ptr = nullptr)
150 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr), sourceVec(vec),
151 betaCoeff(beta_coeff) {}
152
153protected:
154 using OpBase::OpBase;
155
157 boost::shared_ptr<MatrixDouble> sourceVec;
164};
165
166template <int FIELD_DIM, int S, typename OpBase>
168
170 const std::string field_name, boost::shared_ptr<MatrixDouble> vec,
171 ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; },
172 boost::shared_ptr<Range> ents_ptr = nullptr)
173 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr), sourceVec(vec),
174 betaCoeff(beta_coeff) {}
175
176protected:
178 boost::shared_ptr<MatrixDouble> sourceVec;
181};
182
183template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S, IntegrationType I,
184 typename OpBase>
186
187template <int SPACE_DIM, int S, typename OpBase>
189 : public OpBase {
190
191 FTensor::Index<'i', SPACE_DIM> i; ///< summit Index
192
194 const std::string field_name, boost::shared_ptr<MatrixDouble> mat_vals,
195 ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; },
196 boost::shared_ptr<Range> ents_ptr = nullptr)
197 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
198 matVals(mat_vals), betaCoeff(beta_coeff) {}
199
200protected:
201 boost::shared_ptr<MatrixDouble> matVals;
204};
205
206template <int SPACE_DIM, int S, typename OpBase>
208 : public OpBase {
209
210 FTensor::Index<'i', SPACE_DIM> i; ///< summit Index
211 FTensor::Index<'j', SPACE_DIM> j; ///< summit Index
212
214 const std::string field_name, boost::shared_ptr<MatrixDouble> mat_vals,
215 ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; },
216 boost::shared_ptr<Range> ents_ptr = nullptr)
217 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
218 matVals(mat_vals), betaCoeff(beta_coeff) {}
219
220protected:
221 boost::shared_ptr<MatrixDouble> matVals;
224};
225
226template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S, IntegrationType I,
227 typename OpBase>
229
230template <int SPACE_DIM, int S, typename OpBase>
232 : public OpBase {
233
235 const std::string field_name, boost::shared_ptr<MatrixDouble> mat_vals,
236 ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; })
237 : OpBase(field_name, field_name, OpBase::OPROW), matVals(mat_vals),
238 betaCoeff(beta_coeff) {}
239
240protected:
241 boost::shared_ptr<MatrixDouble> matVals;
246};
247
248template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM, IntegrationType I,
249 typename OpBase, CoordinateTypes CoordSys>
251
252template <int FIELD_DIM, int SPACE_DIM, typename OpBase,
253 CoordinateTypes CoordSys>
255 : public OpBase {
257 const std::string field_name, boost::shared_ptr<MatrixDouble> mat_vals,
258 ScalarFun beta = [](double, double, double) { return 1; },
259 boost::shared_ptr<Range> ents_ptr = nullptr)
260 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
261 matVals(mat_vals), betaConst(beta) {}
262
263protected:
265 boost::shared_ptr<MatrixDouble> matVals;
269};
270
271template <int SPACE_DIM, typename OpBase, CoordinateTypes CoordSys>
272struct OpMixDivTimesUImpl<3, 1, SPACE_DIM, GAUSS, OpBase, CoordSys>
273 : public OpBase {
275 const std::string field_name, boost::shared_ptr<VectorDouble> vec_vals,
276 ScalarFun beta = [](double, double, double) constexpr { return 1; },
277 boost::shared_ptr<Range> ents_ptr = nullptr)
278 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
279 vecVals(vec_vals), betaConst(beta) {}
280
281protected:
283 boost::shared_ptr<VectorDouble> vecVals;
286};
287
288template <int FIELD_DIM, typename OpBase, CoordinateTypes CoordSys>
290 : public OpBase {
291
293 const std::string field_name, boost::shared_ptr<VectorDouble> vec,
294 ScalarFun beta = [](double, double, double) constexpr { return 1; },
295 boost::shared_ptr<Range> ents_ptr = nullptr)
296 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr), sourceVec(vec),
297 betaCoeff(beta) {}
298
299protected:
301 boost::shared_ptr<VectorDouble> sourceVec;
303};
304
305/**
306 * @brief Vector field time divergence of tensor
307 *
308 * \f[
309 * (v_i,\lambda_{ij,j})_\Omega
310 * \f]
311 *
312 * @tparam SPACE_DIM
313 * @tparam I
314 * @tparam OpBase
315 */
316template <int SPACE_DIM, IntegrationType I, typename OpBase>
318
319template <int SPACE_DIM, typename OpBase>
322 boost::shared_ptr<MatrixDouble> mat_vals)
323 : OpBase(field_name, field_name, OpBase::OPROW), matVals(mat_vals) {}
324
326 boost::shared_ptr<MatrixDouble> mat_vals,
327 ScalarFun beta_fun)
328 : OpBase(field_name, field_name, OpBase::OPROW), matVals(mat_vals),
329 betaCoeff(beta_fun) {}
330
331protected:
332 boost::shared_ptr<MatrixDouble> matVals;
333 ScalarFun betaCoeff = [](double, double, double) constexpr { return 1; };
334
337};
338
339/**
340 * @brief Tensor field time gradient of vector field
341 *
342 * \f[
343 * (u_i,\lambda_{ij,j})_\Omega
344 * \f]
345 *
346 * @tparam SPACE_DIM
347 * @tparam I
348 * @tparam OpBase
349 */
350template <int SPACE_DIM, IntegrationType I, typename OpBase>
352
353template <int SPACE_DIM, typename OpBase>
356 boost::shared_ptr<MatrixDouble> mat_vals)
357 : OpBase(field_name, field_name, OpBase::OPROW), matVals(mat_vals) {}
358
360 boost::shared_ptr<MatrixDouble> mat_vals,
361 ScalarFun beta_fun)
362 : OpBase(field_name, field_name, OpBase::OPROW), matVals(mat_vals),
363 betaCoeff(beta_fun) {}
364
365protected:
366 boost::shared_ptr<MatrixDouble> matVals;
367 ScalarFun betaCoeff = [](double, double, double) constexpr { return 1; };
368
372};
373
374/**
375 * @brief Multiply vector times normal on the face times scalar function
376 *
377 * This operator typically will be used to evaluate natural boundary conditions
378 * for mixed formulation.
379 *
380 * @tparam BASE_DIM
381 * @tparam SPACE_DIM
382 * @tparam OpBase
383 */
384template <int SPACE_DIM, IntegrationType I, typename OpBase>
386
387/**
388 * @brief This is specialisation for sources on boundary which depends on normal
389 *
390 * @tparam BASE_DIM
391 * @tparam OpBase
392 */
393template <int FIELD_DIM, IntegrationType I, typename OpBase>
396 : public OpNormalMixVecTimesScalarImpl<FIELD_DIM, I, OpBase> {
398 OpBase>::OpNormalMixVecTimesScalarImpl;
399};
400
401template <typename OpBase>
404 const std::string field_name,
405 ScalarFun source_fun = [](double, double, double) constexpr { return 1; },
406 boost::shared_ptr<Range> ents_ptr = nullptr)
407 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
408 sourceFun(source_fun) {}
409
410protected:
414};
415
416template <typename OpBase>
419 const std::string field_name,
420 ScalarFun source_fun = [](double, double, double) constexpr { return 1; },
421 boost::shared_ptr<Range> ents_ptr = nullptr)
422 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
423 sourceFun(source_fun) {}
424
425protected:
429};
430
431/**
432 * @brief Multiply vector times normal on the face times vector field
433 *
434 * This operator typically will be used to evaluate natural boundary conditions
435 * for mixed formulation.
436 *
437 * @tparam BASE_DIM
438 * @tparam SPACE_DIM
439 * @tparam OpBase
440 */
441template <int SPACE_DIM, IntegrationType I, typename OpBase>
443
444template <int SPACE_DIM, typename OpBase>
446 : public OpBase {
448 const std::string field_name, boost::shared_ptr<MatrixDouble> u_ptr,
449 ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; },
450 boost::shared_ptr<Range> ents_ptr = nullptr)
451 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr), uPtr(u_ptr),
452 betaCoeff(beta_coeff) {}
453
454protected:
455 boost::shared_ptr<MatrixDouble> uPtr;
458};
459
460template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM, IntegrationType I,
461 typename OpBase>
463
464template <int SPACE_DIM, typename OpBase>
467 const std::string field_name, boost::shared_ptr<MatrixDouble> u_ptr,
468 boost::shared_ptr<MatrixDouble> y_grad_ptr,
469 ConstantFun source_fun = []() constexpr { return 1; })
470 : OpBase(field_name, field_name, OpBase::OPROW), uPtr(u_ptr),
471 yGradPtr(y_grad_ptr), alphaConstant(source_fun) {}
472
473protected:
474 boost::shared_ptr<MatrixDouble> uPtr;
475 boost::shared_ptr<MatrixDouble> yGradPtr;
478};
479
480template <int FIELD_DIM, int SPACE_DIM, typename OpBase>
482 : public OpBase {
484 const std::string field_name, boost::shared_ptr<MatrixDouble> u_ptr,
485 boost::shared_ptr<MatrixDouble> y_grad_ptr,
486 ConstantFun source_fun = []() constexpr { return 1; })
487 : OpBase(field_name, field_name, OpBase::OPROW), uPtr(u_ptr),
488 yGradPtr(y_grad_ptr), alphaConstant(source_fun) {}
489
490protected:
491 boost::shared_ptr<MatrixDouble> uPtr;
492 boost::shared_ptr<MatrixDouble> yGradPtr;
495};
496
497template <typename OpBase>
500 EntitiesFieldData::EntData &row_data) {
502
503 // get element volume
504 const double vol = OpBase::getMeasure();
505 // get integration weights
506 auto t_w = OpBase::getFTensor0IntegrationWeight();
507 // get base function gradient on rows
508 auto t_row_base = row_data.getFTensor0N();
509 // get coordinate at integration points
510 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
511 // loop over integration points
512 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
513 // take into account Jacobian
514 const double alpha =
515 t_w * vol * sourceFun(t_coords(0), t_coords(1), t_coords(2));
516 // loop over rows base functions
517 int rr = 0;
518 for (; rr != OpBase::nbRows; ++rr) {
519 OpBase::locF[rr] += alpha * t_row_base;
520 ++t_row_base;
521 }
522 for (; rr < OpBase::nbRowBaseFunctions; ++rr)
523 ++t_row_base;
524 ++t_coords;
525 ++t_w; // move to another integration weight
526 }
528}
529
530template <int FIELD_DIM, typename OpBase>
533 iNtegrate(EntitiesFieldData::EntData &row_data) {
536
537 // get element volume
538 const double vol = OpBase::getMeasure();
539 // get integration weights
540 auto t_w = OpBase::getFTensor0IntegrationWeight();
541 // get base function gradient on rows
542 auto t_row_base = row_data.getFTensor0N();
543 // get coordinate at integration points
544 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
545 // loop over integration points
546 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
547 // source file
548 auto t_source = sourceFun(t_coords(0), t_coords(1), t_coords(2));
549 // take into account Jacobian
550 const double alpha = t_w * vol;
551 // loop over rows base functions
552 auto t_nf = getFTensor1FromArray<FIELD_DIM, FIELD_DIM>(OpBase::locF);
553 int rr = 0;
554 for (; rr != OpBase::nbRows / FIELD_DIM; ++rr) {
555 t_nf(i) += alpha * t_row_base * t_source(i);
556 ++t_row_base;
557 ++t_nf;
558 }
559 for (; rr < OpBase::nbRowBaseFunctions; ++rr)
560 ++t_row_base;
561 ++t_coords;
562 ++t_w; // move to another integration weight
563 }
565}
566
567template <int FIELD_DIM, typename OpBase>
570 iNtegrate(EntitiesFieldData::EntData &row_data) {
573
574 const size_t nb_base_functions = row_data.getN().size2() / 3;
575 // get element volume
576 const double vol = OpBase::getMeasure();
577 // get integration weights
578 auto t_w = OpBase::getFTensor0IntegrationWeight();
579 // get base function gradient on rows
580 auto t_row_base = row_data.getFTensor1N<3>();
581 // get coordinate at integration points
582 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
583 // loop over integration points
584 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
585 // source file
586 auto t_source = sourceFun(t_coords(0), t_coords(1), t_coords(2));
587 // take into account Jacobian
588 const double alpha = t_w * vol;
589 // loop over rows base functions
590 int rr = 0;
591 for (; rr != OpBase::nbRows; ++rr) {
592 OpBase::locF[rr] += alpha * t_row_base(i) * t_source(i);
593 ++t_row_base;
594 }
595 for (; rr < nb_base_functions; ++rr)
596 ++t_row_base;
597 ++t_coords;
598 ++t_w; // move to another integration weight
599 }
601}
602
603template <int S, typename OpBase>
605 EntitiesFieldData::EntData &row_data) {
607
608 // get element volume
609 const double vol = OpBase::getMeasure();
610 // get integration weights
611 auto t_w = OpBase::getFTensor0IntegrationWeight();
612 // get base function gradient on rows
613 auto t_row_base = row_data.getFTensor0N();
614 // get vector values
615 auto t_vec = getFTensor0FromVec<S>(*sourceVec);
616
617#ifndef NDEBUG
618 if (sourceVec->size() != OpBase::nbIntegrationPts) {
619 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
620 "Wrong number of integration points %d != %ld",
621 OpBase::nbIntegrationPts, sourceVec->size());
622 }
623#endif
624
625 // get coords
626 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
627 // loop over integration points
628 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
629 // take into account Jacobian
630 const double alpha =
631 t_w * vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
632 // loop over rows base functions
633 int rr = 0;
634 for (; rr != OpBase::nbRows; ++rr) {
635 OpBase::locF[rr] += alpha * t_row_base * t_vec;
636 ++t_row_base;
637 }
638 for (; rr < OpBase::nbRowBaseFunctions; ++rr)
639 ++t_row_base;
640 ++t_w; // move to another integration weight
641 ++t_vec;
642 ++t_coords;
643 }
645}
646
647template <int FIELD_DIM, int S, typename OpBase>
649 EntitiesFieldData::EntData &row_data) {
651
652 // get element volume
653 const double vol = OpBase::getMeasure();
654 // get integration weights
655 auto t_w = OpBase::getFTensor0IntegrationWeight();
656 // get base function gradient on rows
657 auto t_row_base = row_data.getFTensor0N();
658 // get coords
659 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
660 // get vector values
661 auto t_vec = getFTensor1FromMat<FIELD_DIM, S>(*sourceVec);
662 // loop over integration points
663 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
664 // take into account Jacobian
665 const double beta = betaCoeff(t_coords(0), t_coords(1), t_coords(2));
666 const double alpha = t_w * vol * beta;
667 // get loc vector tensor
668 auto t_nf = OpBase::template getNf<FIELD_DIM>();
669 // loop over rows base functions
670 int rr = 0;
671 for (; rr != OpBase::nbRows / FIELD_DIM; ++rr) {
672 t_nf(i) += alpha * t_row_base * t_vec(i);
673 ++t_row_base;
674 ++t_nf;
675 }
676 for (; rr < OpBase::nbRowBaseFunctions; ++rr)
677 ++t_row_base;
678 ++t_w; // move to another integration weight
679 ++t_vec;
680 ++t_coords;
681 }
683}
684
685template <int FIELD_DIM, int S, typename OpBase>
687 EntitiesFieldData::EntData &row_data) {
689
690 const size_t nb_base_functions = row_data.getN().size2() / 3;
691 // get element volume
692 const double vol = OpBase::getMeasure();
693 // get integration weights
694 auto t_w = OpBase::getFTensor0IntegrationWeight();
695 // get base function gradient on rows
696 auto t_row_base = row_data.getFTensor1N<3>();
697 // get coords
698 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
699 // get vector values
700 auto t_vec = getFTensor1FromMat<FIELD_DIM, S>(*sourceVec);
701 // loop over integration points
702 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
703 // take into account Jacobian
704 const double alpha =
705 t_w * vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
706 // loop over rows base functions
707 int rr = 0;
708 for (; rr != OpBase::nbRows; ++rr) {
709 OpBase::locF[rr] += alpha * t_row_base(i) * t_vec(i);
710 ++t_row_base;
711 }
712 for (; rr < nb_base_functions; ++rr)
713 ++t_row_base;
714 ++t_w; // move to another integration weight
715 ++t_vec;
716 ++t_coords;
717 }
719}
720
721template <int SPACE_DIM, int S, typename OpBase>
724 EntitiesFieldData::EntData &row_data) {
726
727 // get element volume
728 const double vol = OpBase::getMeasure();
729 // get integration weights
730 auto t_w = OpBase::getFTensor0IntegrationWeight();
731 // get base function gradient on rows
732 auto t_row_grad = row_data.getFTensor1DiffN<SPACE_DIM>();
733 // get field gradient values
734 auto t_val_grad = getFTensor1FromMat<SPACE_DIM, S>(*(matVals));
735 // get coordinate at integration points
736 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
737 // loop over integration points
738 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
739 const double beta = vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
740 // take into account Jacobian
741 const double alpha = t_w * beta;
742 // loop over rows base functions
743 int rr = 0;
744 for (; rr != OpBase::nbRows; rr++) {
745 // calculate element of local matrix
746 OpBase::locF[rr] += alpha * (t_row_grad(i) * t_val_grad(i));
747 ++t_row_grad; // move to another element of gradient of base
748 // function on row
749 }
750 for (; rr < OpBase::nbRowBaseFunctions; ++rr)
751 ++t_row_grad;
752
753 ++t_coords;
754 ++t_val_grad;
755 ++t_w; // move to another integration weight
756 }
758}
759
760template <int SPACE_DIM, int S, typename OpBase>
763 EntitiesFieldData::EntData &row_data) {
765
766 // get element volume
767 const double vol = OpBase::getMeasure();
768 // get integration weights
769 auto t_w = OpBase::getFTensor0IntegrationWeight();
770 // get base function gradient on rows
771 auto t_row_grad = row_data.getFTensor1DiffN<SPACE_DIM>();
772 // get field gradient values
773 auto t_val_grad = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(matVals));
774 // get coordinate at integration points
775 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
776 // loop over integration points
777 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
778 // take into account Jacobian
779 const double alpha =
780 t_w * vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
781 // get rhs vector
782 auto t_nf = OpBase::template getNf<SPACE_DIM>();
783 // loop over rows base functions
784 int rr = 0;
785 for (; rr != OpBase::nbRows / SPACE_DIM; rr++) {
786 // calculate element of local matrix
787 t_nf(i) += alpha * (t_row_grad(j) * t_val_grad(i, j));
788 ++t_row_grad; // move to another element of gradient of base
789 // function on row
790 ++t_nf;
791 }
792 for (; rr < OpBase::nbRowBaseFunctions; ++rr)
793 ++t_row_grad;
794
795 ++t_coords;
796 ++t_val_grad;
797 ++t_w; // move to another integration weight
798 }
800}
801
802template <int SPACE_DIM, int S, typename OpBase>
805 EntitiesFieldData::EntData &row_data) {
807 // get element volume
808 const double vol = OpBase::getMeasure();
809 // get coords
810 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
811 // get integration weights
812 auto t_w = OpBase::getFTensor0IntegrationWeight();
813 // get base function gradient on rows
814 auto t_row_grad = row_data.getFTensor1DiffN<SPACE_DIM>();
815 // get field gradient values
816 auto t_val_mat = getFTensor2SymmetricFromMat<SPACE_DIM, S>(*(matVals));
817 // loop over integration points
818 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
819 // take into account Jacobian
820 const double alpha =
821 t_w * vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
822 // get rhs vector
823 auto t_nf = OpBase::template getNf<SPACE_DIM>();
824 // loop over rows base functions
825 int rr = 0;
826 for (; rr != OpBase::nbRows / SPACE_DIM; rr++) {
827 // calculate element of local matrix
828 t_nf(j) += alpha * (t_row_grad(i) * t_val_mat(i, j));
829 ++t_row_grad; // move to another element of gradient of base
830 // function on row
831 ++t_nf;
832 }
833 for (; rr < OpBase::nbRowBaseFunctions; ++rr)
834 ++t_row_grad;
835 ++t_val_mat;
836 ++t_coords;
837 ++t_w; // move to another integration weight
838 }
840}
841
842template <int FIELD_DIM, int SPACE_DIM, typename OpBase,
843 CoordinateTypes CoordSys>
846 EntitiesFieldData::EntData &row_data) {
848
849 const size_t nb_base_functions = row_data.getN().size2() / 3;
850 auto t_w = this->getFTensor0IntegrationWeight();
851 // get coordinate at integration points
852 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
853 auto t_diff_base = row_data.getFTensor2DiffN<3, SPACE_DIM>();
854 auto t_base = row_data.getFTensor1N<3>();
855 auto t_u = getFTensor1FromMat<FIELD_DIM>(*(matVals));
856
857 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
858
859 const double alpha = this->getMeasure() * t_w *
860 betaConst(t_coords(0), t_coords(1), t_coords(2));
861 auto t_nf = OpBase::template getNf<FIELD_DIM>();
862
863 size_t bb = 0;
864 for (; bb != this->nbRows / FIELD_DIM; ++bb) {
865 const double t_div_base = t_diff_base(j, j);
866 t_nf(i) += alpha * t_div_base * t_u(i);
867 if constexpr (CoordSys == CYLINDRICAL) {
868 t_nf(i) += t_base(0) * (alpha / t_coords(0)) * t_u(i);
869 }
870 ++t_nf;
871 ++t_diff_base;
872 ++t_base;
873 }
874 for (; bb < nb_base_functions; ++bb) {
875 ++t_diff_base;
876 ++t_base;
877 }
878
879 ++t_u;
880 ++t_w;
881 ++t_coords;
882 }
883
885}
886
887template <int SPACE_DIM, typename OpBase, CoordinateTypes CoordSys>
890 EntitiesFieldData::EntData &row_data) {
892
893 const size_t nb_base_functions = row_data.getN().size2() / 3;
894 auto t_w = this->getFTensor0IntegrationWeight();
895 // get coordinate at integration points
896 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
897 auto t_diff_base = row_data.getFTensor2DiffN<3, SPACE_DIM>();
898 auto t_u = getFTensor0FromVec(*(vecVals));
899
900 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
901
902 const double alpha = this->getMeasure() * t_w *
903 betaConst(t_coords(0), t_coords(1), t_coords(2));
904 ;
905
906 size_t bb = 0;
907 for (; bb != this->nbRows; ++bb) {
908 const double t_div_base = t_diff_base(j, j);
909 OpBase::locF[bb] += alpha * t_div_base * t_u;
910 ++t_diff_base;
911 }
912 for (; bb < nb_base_functions; ++bb)
913 ++t_diff_base;
914
915 ++t_u;
916 ++t_w;
917 ++t_coords;
918 }
919
921}
922
923/**
924 * @brief div U times vector
925 *
926 * \f[
927 * \delta u_j = \phi^m\delta\overline{u}^m_j\\
928 * \delta u_{j,i} = \phi^m_{,i}\delta\overline{u}^m_j\\
929 * \textrm{tr}[\delta u_{j,i}] = \delta u_{j,i}\delta_{ji}\\
930 * (\textrm{tr}[\delta u_{j,i}], v) =\\
931 * (\delta u_{j,i} \delta_{ij}, v) =\\
932 * (\delta u_{j,i}, \delta_{ij} v) =\\
933 * (\phi^m_{,i}\delta\overline{u}^m_j, \delta_{ij} v) \\
934 * f_i^m=(\phi^m_{,i}, v)
935 * \f]
936 *
937 * @tparam FIELD_DIM
938 * @tparam SPACE_DIM
939 * @tparam OpBase
940 * @param row_data
941 * @return MoFEMErrorCode
942 */
943template <int FIELD_DIM, typename OpBase, CoordinateTypes CoordSys>
946 EntitiesFieldData::EntData &row_data) {
949
950 // get element volume
951 const double vol = OpBase::getMeasure();
952 // get integration weights
953 auto t_w = OpBase::getFTensor0IntegrationWeight();
954 // get coordinate at integration points
955 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
956 // get base function gradient on rows
957 auto t_row_grad = row_data.getFTensor1DiffN<FIELD_DIM>();
958 // get vector values
959 auto t_vec = getFTensor0FromVec(*sourceVec);
960 // loop over integration points
961 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
962 // take into account Jacobian
963 const double alpha =
964 t_w * vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
965 auto t_nf = getFTensor1FromArray<FIELD_DIM, FIELD_DIM>(OpBase::locF);
966 // loop over rows base functions
967 int rr = 0;
968 for (; rr != OpBase::nbRows / FIELD_DIM; ++rr) {
969 t_nf(i) += alpha * t_row_grad(i) * t_vec;
970 ++t_row_grad;
971 ++t_nf;
972 }
973 for (; rr < OpBase::nbRowBaseFunctions; ++rr)
974 ++t_row_grad;
975 ++t_w; // move to another integration weight
976 ++t_vec;
977 ++t_coords;
978 }
980}
981
982template <int SPACE_DIM, typename OpBase>
984 EntitiesFieldData::EntData &row_data) {
986
987 const size_t nb_base_functions = row_data.getN().size2() / 3;
988 auto t_w = this->getFTensor0IntegrationWeight();
989 auto t_coords = this->getFTensor1CoordsAtGaussPts();
990 auto t_base = row_data.getFTensor1N<3>();
991 auto t_grad = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(matVals));
992
993 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
994
995 const double alpha = this->getMeasure() * t_w;
996 auto t_nf = OpBase::template getNf<SPACE_DIM>();
997
998 size_t bb = 0;
999 for (; bb != this->nbRows / SPACE_DIM; ++bb) {
1000 t_nf(i) += alpha * betaCoeff(t_coords(0), t_coords(1), t_coords(2)) *
1001 t_base(j) * t_grad(i, j);
1002 ++t_nf;
1003 ++t_base;
1004 }
1005 for (; bb < nb_base_functions; ++bb)
1006 ++t_base;
1007
1008 ++t_grad;
1009 ++t_coords;
1010 ++t_w;
1011 }
1012
1014}
1015
1016template <int SPACE_DIM, typename OpBase>
1018 EntitiesFieldData::EntData &row_data) {
1020
1021 const size_t nb_base_functions = row_data.getN().size2();
1022 auto t_w = this->getFTensor0IntegrationWeight();
1023 auto t_coords = this->getFTensor1CoordsAtGaussPts();
1024 auto t_base = row_data.getFTensor0N();
1025 auto t_div = getFTensor1FromMat<SPACE_DIM>(*(matVals));
1026
1027 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
1028 const double alpha = this->getMeasure() * t_w;
1029 auto t_nf = OpBase::template getNf<SPACE_DIM>();
1030
1031 size_t bb = 0;
1032 for (; bb != this->nbRows / SPACE_DIM; ++bb) {
1033 t_nf(i) += alpha * t_base *
1034 betaCoeff(t_coords(0), t_coords(1), t_coords(2)) * t_div(i);
1035 ++t_nf;
1036 ++t_base;
1037 }
1038 for (; bb < nb_base_functions; ++bb)
1039 ++t_base;
1040
1041 ++t_coords;
1042 ++t_div;
1043 ++t_w;
1044 }
1045
1047}
1048
1049template <typename OpBase>
1051 EntitiesFieldData::EntData &row_data) {
1053
1054 const size_t nb_base_functions = row_data.getN().size2() / 3;
1055 // get element volume
1056 // get integration weights
1057 auto t_w = OpBase::getFTensor0IntegrationWeight();
1058 // get base function gradient on rows
1059 auto t_row_base = row_data.getFTensor1N<3>();
1060 // get coordinate at integration points
1061 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
1062 // get normal
1063 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
1064 double a = 1;
1065 if (this->getNumeredEntFiniteElementPtr()->getEntType() == MBTRI)
1066 a *= 2;
1067 // loop over integration points
1068 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
1069 // take into account Jacobian
1070 const double alpha =
1071 t_w * sourceFun(t_coords(0), t_coords(1), t_coords(2)) / a;
1072 // loop over rows base functions
1073 int rr = 0;
1074 for (; rr != OpBase::nbRows; ++rr) {
1075 OpBase::locF[rr] += alpha * t_row_base(i) * t_normal(i);
1076 ++t_row_base;
1077 }
1078 for (; rr < nb_base_functions; ++rr)
1079 ++t_row_base;
1080 ++t_coords;
1081 ++t_w; // move to another integration weight
1082 ++t_normal;
1083 }
1085}
1086
1087template <typename OpBase>
1089 EntitiesFieldData::EntData &row_data) {
1091
1092 const size_t nb_base_functions = row_data.getN().size2() / 3;
1093 FTensor::Tensor1<double, 3> t_z{0., 0., 1.};
1094 // get element volume
1095 // get integration weights
1096 auto t_w = OpBase::getFTensor0IntegrationWeight();
1097 // get base function gradient on rows
1098 auto t_row_base = row_data.getFTensor1N<3>();
1099 // get coordinate at integration points
1100 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
1101 // get normal
1102 auto t_tangent = OpBase::getFTensor1TangentAtGaussPts();
1103 // loop over integration points
1104 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
1105 // take into account Jacobian
1106 const double alpha = t_w * sourceFun(t_coords(0), t_coords(1), t_coords(2));
1108 FTensor::Index<'j', 3> j;
1109 FTensor::Index<'k', 3> k;
1110 t_normal(i) = FTensor::levi_civita(i, j, k) * t_tangent(j) * t_z(k);
1111 int rr = 0;
1112 for (; rr != OpBase::nbRows; ++rr) {
1113 OpBase::locF[rr] += alpha * t_row_base(i) * t_normal(i);
1114 ++t_row_base;
1115 }
1116 for (; rr < nb_base_functions; ++rr)
1117 ++t_row_base;
1118 ++t_coords;
1119 ++t_tangent;
1120 ++t_w; // move to another integration weight
1121 }
1123}
1124
1125template <int SPACE_DIM, typename OpBase>
1128 EntitiesFieldData::EntData &row_data) {
1130 FTensor::Index<'i', 3> i;
1132
1133 const size_t nb_base_functions = row_data.getN().size2() / 3;
1134 // get element volume
1135 // get integration weights
1136 auto t_w = OpBase::getFTensor0IntegrationWeight();
1137 // get base function gradient on rows
1138 auto t_row_base = row_data.getFTensor1N<3>();
1139 // get coordinate at integration points
1140 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
1141 // get normal
1142 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
1143 // get field
1144 auto t_u = getFTensor1FromMat<SPACE_DIM>(*uPtr);
1145 // loop over integration points
1146 auto a = OpBase::getMeasure();
1147 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
1148 // take into account Jacobian
1149 auto l2 = std::sqrt(t_normal(i) * t_normal(i));
1150 const double alpha =
1151 t_w * betaCoeff(t_coords(0), t_coords(1), t_coords(2)) * (a / l2);
1152 // get rhs vector
1153 auto t_nf = OpBase::template getNf<SPACE_DIM>();
1154 // loop over rows base functions
1155 int rr = 0;
1156 for (; rr != OpBase::nbRows / SPACE_DIM; ++rr) {
1157 t_nf(J) += alpha * (t_row_base(i) * t_normal(i)) * t_u(J);
1158 ++t_row_base;
1159 ++t_nf;
1160 }
1161 for (; rr < nb_base_functions; ++rr)
1162 ++t_row_base;
1163 ++t_coords;
1164 ++t_w; // move to another integration weight
1165 ++t_normal;
1166 ++t_u;
1167 }
1169}
1170
1171template <int SPACE_DIM, typename OpBase>
1174 EntitiesFieldData::EntData &row_data) {
1176
1177 auto t_w = this->getFTensor0IntegrationWeight();
1178 auto t_base = row_data.getFTensor0N();
1179
1180 auto t_u = getFTensor1FromMat<SPACE_DIM>(*uPtr);
1181 auto t_grad_y = getFTensor1FromMat<SPACE_DIM>(*yGradPtr);
1182
1184 const double alpha_constant = alphaConstant();
1185 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
1186
1187 // get element volume
1188 const double vol = OpBase::getMeasure();
1189 const double c = (t_grad_y(i) * t_u(i)) * (t_w * vol * alpha_constant);
1190
1191 // get element volume
1192 int rr = 0;
1193 for (; rr != OpBase::nbRows; ++rr) {
1194 OpBase::locF[rr] += c * t_base;
1195 ++t_base;
1196 }
1197 for (; rr < OpBase::nbRowBaseFunctions; ++rr)
1198 ++t_base;
1199
1200 ++t_w; // move to another integration weight
1201 ++t_u;
1202 ++t_grad_y;
1203 }
1204
1206}
1207
1208template <int FIELD_DIM, int SPACE_DIM, typename OpBase>
1211 EntitiesFieldData::EntData &row_data) {
1213
1214 auto t_w = this->getFTensor0IntegrationWeight();
1215 auto t_base = row_data.getFTensor0N();
1216
1217 auto t_u = getFTensor1FromMat<SPACE_DIM>(*uPtr);
1218 auto t_grad_y = getFTensor2FromMat<FIELD_DIM, SPACE_DIM>(*yGradPtr);
1219
1222 const double alpha_constant = alphaConstant();
1223 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
1224
1225 // get element volume
1226 const double vol = OpBase::getMeasure();
1227
1229 t_c(J) = (t_grad_y(J, i) * t_u(i)) * (t_w * vol * alpha_constant);
1230
1231 auto t_nf = getFTensor1FromArray<FIELD_DIM, FIELD_DIM>(OpBase::locF);
1232 int rr = 0;
1233 for (; rr != OpBase::nbRows / FIELD_DIM; ++rr) {
1234 t_nf(J) += t_base * t_c(J);
1235 ++t_base;
1236 ++t_nf;
1237 }
1238 for (; rr < OpBase::nbRowBaseFunctions; ++rr)
1239 ++t_base;
1240 ++t_w; // move to another integration weight
1241
1242 ++t_u;
1243 ++t_grad_y;
1244 }
1245
1247}
1248
1249} // namespace MoFEM
1250
1251#endif // __LINEAR_FORMS_INTEGRATORS_HPP__
constexpr double a
constexpr int SPACE_DIM
[Define dimension]
constexpr int FIELD_DIM
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
CoordinateTypes
Coodinate system.
@ CYLINDRICAL
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
constexpr int BASE_DIM
IntegrationType
Form integrator integration types.
boost::function< double(const double, const double, const double)> ScalarFun
Scalar function type.
boost::function< FTensor::Tensor1< double, DIM >(const double, const double, const double)> VectorFun
Vector function type.
@ GAUSS
Gaussian quadrature integration.
FTensor::Index< 'i', SPACE_DIM > i
const double c
speed of light (cm/ns)
FTensor::Index< 'J', DIM1 > J
Definition level_set.cpp:30
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
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
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
boost::function< double(double)> TimeFun
Lambda function used to scale with time.
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
boost::function< double()> ConstantFun
Constant function type.
constexpr IntegrationType I
constexpr auto field_name
OpBaseImpl< PETSC, EdgeEleOp > OpBase
Definition radiation.cpp:29
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2DiffN(FieldApproximationBase base)
Get derivatives of base functions for Hdiv space.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
VectorDouble locF
local entity vector
int nbRows
number of dofs on rows
int nbIntegrationPts
number of integration points
int nbRowBaseFunctions
number or row base functions
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)
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)
MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
Integrate grad-grad operator.
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)
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;})
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;})
OpGradTimesSymTensorImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals, ScalarFun beta_coeff=[](double, double, double) constexpr { return 1;})
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)
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)
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)
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)
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)
OpMixTensorTimesGradUImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals)
OpMixTensorTimesGradUImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals, ScalarFun beta_fun)
Tensor field time gradient of vector field.
OpMixVecTimesDivLambdaImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals, ScalarFun beta_fun)
OpMixVecTimesDivLambdaImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals)
Vector field time divergence of tensor.
OpNormalMixVecTimesScalarImpl(const std::string field_name, ScalarFun source_fun=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
OpNormalMixVecTimesScalarImpl(const std::string field_name, ScalarFun source_fun=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
Multiply vector times normal on the face times scalar function.
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)
Multiply vector times normal on the face times vector field.
OpSourceImpl(const std::string field_name, ScalarFun source_fun, boost::shared_ptr< Range > ents_ptr=nullptr)
Construct a new Op Source Impl object.
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.
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.
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.
OpSourceImpl(const std::string field_name, TimeFun time_fun, VectorFun< FIELD_DIM > source_fun, boost::shared_ptr< Range > ents_ptr=nullptr)
OpSourceImpl(const std::string field_name, VectorFun< FIELD_DIM > source_fun, boost::shared_ptr< Range > ents_ptr=nullptr)