ug4
dom_disc_embb.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2020: G-CSC, Goethe University Frankfurt
3  * Author: Dmitry Logashenko
4  *
5  * This file is part of UG4.
6  *
7  * UG4 is free software: you can redistribute it and/or modify it under the
8  * terms of the GNU Lesser General Public License version 3 (as published by the
9  * Free Software Foundation) with the following additional attribution
10  * requirements (according to LGPL/GPL v3 §7):
11  *
12  * (1) The following notice must be displayed in the Appropriate Legal Notices
13  * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14  *
15  * (2) The following notice must be displayed at a prominent place in the
16  * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17  *
18  * (3) The following bibliography is recommended for citation and must be
19  * preserved in all covered files:
20  * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21  * parallel geometric multigrid solver on hierarchically distributed grids.
22  * Computing and visualization in science 16, 4 (2013), 151-164"
23  * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24  * flexible software system for simulating pde based models on high performance
25  * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26  *
27  * This program is distributed in the hope that it will be useful,
28  * but WITHOUT ANY WARRANTY; without even the implied warranty of
29  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30  * GNU Lesser General Public License for more details.
31  */
32 
33 /*
34  * Global assembling of the problems with the embedded boundary.
35  */
36 #ifndef __H__UG__PLUGINS__D3F__EMBASS__
37 #define __H__UG__PLUGINS__D3F__EMBASS__
38 
39 #include <vector>
40 
41 // ug4 headers
42 #include "common/common.h"
44 #include "lib_disc/domain_traits.h"
48 #ifdef UG_FOR_LUA
50 #endif
51 
52 namespace ug {
53 
55 
62 template <typename TDomain, typename TAlgebra>
64 {
65 public:
66 
68  typedef TDomain domain_type;
69 
71  typedef TAlgebra algebra_type;
72 
74  typedef typename algebra_type::vector_type vector_type;
75 
77  typedef typename algebra_type::matrix_type matrix_type;
78 
80  static const int dim = domain_type::dim;
81 
84 
87 
89  virtual int check_elem_lsf
90  (
91  size_t n_co,
92  GridObject * pElem,
93  int si,
94  int g_level,
95  bool use_hanging,
96  const MathVector<dim> vCornerCoords [],
97  number time
98  )
99  {UG_THROW ("IInterfaceExtrapolation: Virtual functions are not implemented in the base class.");}
100 
102  virtual int check_elem_lsf
103  (
104  size_t n_co,
105  GridObject * pElem,
106  int si,
107  bool use_hanging,
108  const MathVector<dim> vCornerCoords [],
109  number time
110  )
111  {UG_THROW ("IInterfaceExtrapolation: Virtual functions are not implemented in the base class.");}
112 
114  virtual void extrapolate_by_lsf
115  (
116  size_t num_co,
117  size_t base_co,
118  number * u,
119  size_t fct
120  ) const
121  {UG_THROW ("IInterfaceExtrapolation: Virtual functions are not implemented in the base class.");}
122 
124  virtual void extrapolate_by_lsf
125  (
126  size_t num_co,
127  number * u,
128  size_t fct
129  ) const
130  {UG_THROW ("IInterfaceExtrapolation: Virtual functions are not implemented in the base class.");}
131 
133  virtual bool corner_inside
134  (
135  size_t co
136  ) const
137  {UG_THROW ("IInterfaceExtrapolation: Virtual functions are not implemented in the base class.");}
138 
140  virtual number lsf_at
141  (
142  size_t co
143  ) const
144  {UG_THROW ("IInterfaceExtrapolation: Virtual functions are not implemented in the base class.");}
145 
146 }; // class IInterfaceExtrapolation
147 
149 
157 template <typename TDomain, typename TAlgebra, typename TExtrapolation>
159 {
160 public:
161 
163  typedef TDomain domain_type;
164 
166  typedef TAlgebra algebra_type;
167 
170 
172  typedef typename algebra_type::vector_type vector_type;
173 
175  typedef typename algebra_type::matrix_type matrix_type;
176 
178  typedef TExtrapolation extrapolation_type;
179 
181  typedef typename extrapolation_type::ls_grid_func_type ls_grid_func_type;
182 
184  static const int dim = TDomain::dim;
185 
187 // Constructor/destructor
189 
190 public:
191 
194 
196  virtual ~LSGFGlobAssembler () {};
197 
199 // Assembling tools
201 
202 public:
203 
204  template <typename TElem, typename TIterator>
205  void
206  AssembleStiffnessMatrix( const std::vector<IElemDisc<domain_type>*>& vElemDisc,
209  TIterator iterBegin,
210  TIterator iterEnd,
211  int si, bool bNonRegularGrid,
212  matrix_type& A,
213  const vector_type& u,
215 
216  template <typename TElem, typename TIterator>
217  void
218  AssembleMassMatrix( const std::vector<IElemDisc<domain_type>*>& vElemDisc,
221  TIterator iterBegin,
222  TIterator iterEnd,
223  int si, bool bNonRegularGrid,
224  matrix_type& M,
225  const vector_type& u,
227 
228  template <typename TElem, typename TIterator>
229  void
230  AssembleJacobian( const std::vector<IElemDisc<domain_type>*>& vElemDisc,
233  TIterator iterBegin,
234  TIterator iterEnd,
235  int si, bool bNonRegularGrid,
236  matrix_type& J,
237  const vector_type& u,
239 
240  template <typename TElem, typename TIterator>
241  void
242  AssembleJacobian( const std::vector<IElemDisc<domain_type>*>& vElemDisc,
245  TIterator iterBegin,
246  TIterator iterEnd,
247  int si, bool bNonRegularGrid,
248  matrix_type& J,
250  number s_a0,
252 
253  template <typename TElem, typename TIterator>
254  void
255  AssembleDefect( const std::vector<IElemDisc<domain_type>*>& vElemDisc,
258  TIterator iterBegin,
259  TIterator iterEnd,
260  int si, bool bNonRegularGrid,
261  vector_type& d,
262  const vector_type& u,
264 
265  template <typename TElem, typename TIterator>
266  void
267  AssembleDefect( const std::vector<IElemDisc<domain_type>*>& vElemDisc,
270  TIterator iterBegin,
271  TIterator iterEnd,
272  int si, bool bNonRegularGrid,
273  vector_type& d,
275  const std::vector<number>& vScaleMass,
276  const std::vector<number>& vScaleStiff,
278 
279  template <typename TElem, typename TIterator>
280  void
281  AssembleLinear( const std::vector<IElemDisc<domain_type>*>& vElemDisc,
284  TIterator iterBegin,
285  TIterator iterEnd,
286  int si, bool bNonRegularGrid,
287  matrix_type& A,
288  vector_type& rhs,
290 
291  template <typename TElem, typename TIterator>
292  void
293  AssembleLinear( const std::vector<IElemDisc<domain_type>*>& vElemDisc,
296  TIterator iterBegin,
297  TIterator iterEnd,
298  int si, bool bNonRegularGrid,
299  matrix_type& A,
300  vector_type& rhs,
302  const std::vector<number>& vScaleMass,
303  const std::vector<number>& vScaleStiff,
305 
307 // Assemble Rhs: it cannot be done for the ghost-fluid method independently of the matrix
309 
310 public:
311 
312  template <typename TElem, typename TIterator>
313  static void
314  AssembleRhs( const std::vector<IElemDisc<domain_type>*>& vElemDisc,
317  TIterator iterBegin,
318  TIterator iterEnd,
319  int si, bool bNonRegularGrid,
320  vector_type& rhs,
321  const vector_type& u,
323  {
324  UG_THROW ("LSGFGlobAssembler::AssembleRhs: Cannot assemble the RHS in GF independently of the matrix");
325  }
326 
327  template <typename TElem, typename TIterator>
328  static void
329  AssembleRhs( const std::vector<IElemDisc<domain_type>*>& vElemDisc,
332  TIterator iterBegin,
333  TIterator iterEnd,
334  int si, bool bNonRegularGrid,
335  vector_type& rhs,
337  const std::vector<number>& vScaleMass,
338  const std::vector<number>& vScaleStiff,
340  {
341  UG_THROW ("LSGFGlobAssembler::AssembleRhs: Cannot assemble the RHS in GF independently of the matrix");
342  }
343 
345 // Prepare and Finish Timestep: these version merely skip the outer elements
347 
348 public:
360  void PrepareTimestep
361  (
362  const std::vector<IElemDisc<domain_type>*>& vElemDisc,
364  bool bNonRegularGrid,
366  number future_time,
368  );
369 
370  template <typename TElem, typename TIterator>
371  void
372  PrepareTimestepElem(const std::vector<IElemDisc<domain_type>*>& vElemDisc,
375  TIterator iterBegin,
376  TIterator iterEnd,
377  int si, bool bNonRegularGrid,
391  void FinishTimestep
392  (
393  const std::vector<IElemDisc<domain_type>*>& vElemDisc,
395  bool bNonRegularGrid,
398  );
399 
400  template <typename TElem, typename TIterator>
401  void
402  FinishTimestepElem(const std::vector<IElemDisc<domain_type>*>& vElemDisc,
405  TIterator iterBegin,
406  TIterator iterEnd,
407  int si, bool bNonRegularGrid,
410 
411  template <typename TElem, typename TIterator>
412  void
413  InitAllExports(const std::vector<IElemDisc<domain_type>*>& vElemDisc,
415  TIterator iterBegin,
416  TIterator iterEnd,
417  int si, bool bNonRegularGrid, bool bAsTimeDependent);
418 
420 // Error estimators: Not implemented for the ghost-fluid method
422 
423 public:
424 
425  template <typename TElem, typename TIterator>
426  static void
428  (
429  const std::vector<IElemError<domain_type>*>& vElemDisc,
432  TIterator iterBegin,
433  TIterator iterEnd,
434  int si,
435  bool bNonRegularGrid,
436  const vector_type& u
437  )
438  {
439  UG_THROW ("AssembleErrorEstimator: No error estimator implemented for the Ghost-Fluid method.");
440  }
441 
442  template <typename TElem, typename TIterator>
443  static void
445  (
446  const std::vector<IElemError<domain_type>*>& vElemDisc,
449  TIterator iterBegin,
450  TIterator iterEnd,
451  int si,
452  bool bNonRegularGrid,
453  std::vector<number> vScaleMass,
454  std::vector<number> vScaleStiff,
456  )
457  {
458  UG_THROW ("AssembleErrorEstimator: No error estimator implemented for the Ghost-Fluid method.");
459  }
460 
462 // Data
464 
465 private:
466 
468 
476 
477 public:
478 
480  void set_LSF
481  (
483  )
484  {
485  m_extrapol.set_LSF (spLSF);
486  }
487 
490  (
491  SmartPtr<approx_space_type> spApproxSpace
492  )
493  {
494  m_extrapol.prepare_interface_bc (spApproxSpace);
495  }
496 
499  (
500  SmartPtr<approx_space_type> spApproxSpace,
501  const char* fct_name,
502  number value
503  )
504  {
505  m_extrapol.set_Dirichlet_for (spApproxSpace, fct_name, value);
506  }
507 
510  (
511  SmartPtr<approx_space_type> spApproxSpace,
512  const char* fct_name,
514  )
515  {
516  m_extrapol.set_Dirichlet_for (spApproxSpace, fct_name, func);
517  }
518 
521  (
522  SmartPtr<approx_space_type> spApproxSpace,
523  const char* fct_name,
525  )
526  {
527  m_extrapol.set_plain_Dirichlet_for (spApproxSpace, fct_name, func);
528  }
529 
532  (
533  SmartPtr<approx_space_type> spApproxSpace,
534  const char* fct_name,
535  number value
536  )
537  {
538  m_extrapol.set_plain_Dirichlet_for (spApproxSpace, fct_name, value);
539  }
540 
543  (
544  SmartPtr<approx_space_type> spApproxSpace,
545  const char* fct_name
546  )
547  {
548  m_extrapol.set_Neumann0_for (spApproxSpace, fct_name);
549  }
550 
553  (
554  SmartPtr<approx_space_type> spApproxSpace,
555  const char* subset_names
556  )
557  {
558  m_extrapol.exclude_subsets (spApproxSpace, subset_names);
559  }
560 
562  void set_assemble_only_cut (bool b)
563  {
564  m_bAssembleOnlyCut = b;
565  }
566 
568  void project_LSF ()
569  {
570  m_extrapol.project_LSF ();
571  }
572 
575 
578  (
579  size_t n_co,
580  GridObject * pElem,
581  int si,
582  int g_level,
583  bool use_hanging,
584  const MathVector<dim> vCornerCoords [],
585  number time
586  )
587  {
588  return m_extrapol.check_elem_lsf (n_co, pElem, si, g_level, use_hanging, vCornerCoords, time);
589  }
590 
593  (
594  size_t num_co,
595  size_t base_co,
596  number * u,
597  size_t fct
598  ) const
599  {
600  m_extrapol.extrapolate_by_lsf (num_co, base_co, u, fct);
601  }
602 
605  (
606  size_t num_co,
607  number * u,
608  size_t fct
609  ) const
610  {
611  m_extrapol.extrapolate_by_lsf (num_co, u, fct);
612  }
613 
616  (
617  size_t co
618  ) const
619  {return m_extrapol.corner_inside (co);}
620 
623  (
624  size_t co
625  ) const
626  {return m_extrapol.lsf_at (co);}
627 
629  virtual void clear_outer_values
630  (
631  vector_type & d,
632  const DoFDistribution * dd
633  ) const
634  {m_extrapol.clear_outer_values (d, dd);}
635 
637  virtual void set_outer_values
638  (
639  vector_type & u,
640  const DoFDistribution * dd,
641  number time
642  )
643  {m_extrapol.set_outer_values (u, dd, time);}
644 
646  virtual void set_outer_matrices
647  (
648  matrix_type & A,
649  const DoFDistribution * dd
650  ) const
651  {m_extrapol.set_outer_matrices (A, dd);}
652 
653 }; // class LSGFGlobAssembler
654 
656 
664 template <typename TDomain, typename TAlgebra, typename TExtrapolation>
666 : public IDomainConstraint<TDomain, TAlgebra>
667 {
668 public:
669 
671  typedef TDomain domain_type;
672 
674  typedef TAlgebra algebra_type;
675 
677  typedef typename algebra_type::vector_type vector_type;
678 
680  typedef typename algebra_type::matrix_type matrix_type;
681 
683  typedef TExtrapolation extrapolation_type;
684 
687  (
688  extrapolation_type & rExtrapolation
689  )
690  : m_rExtrapolation (rExtrapolation)
691  {}
692 
694  virtual ~LSGFConstraint () {};
695 
698  (
699  matrix_type & J,
700  const vector_type & u,
702  int type,
703  number time = 0.0,
705  const number s_a0 = 1.0
706  )
707  {
708  m_rExtrapolation.set_outer_matrices (J, dd.get ());
709  }
710 
713  (
714  vector_type & d,
715  const vector_type & u,
717  int type,
718  number time = 0.0,
720  const std::vector<number> * vScaleMass = NULL,
721  const std::vector<number> * vScaleStiff = NULL
722  )
723  {
724  m_rExtrapolation.clear_outer_values (d, dd.get ());
725  }
726 
729  (
730  vector_type & u,
732  int type,
733  number time = 0.0
734  )
735  {
736  m_rExtrapolation.set_outer_values (u, dd.get (), time);
737  }
738 
741  (
742  matrix_type & A,
743  vector_type & b,
745  int type,
746  number time = 0.0
747  )
748  { // Note that this function is not really used, so it needs not to be optimal.
749  m_rExtrapolation.set_outer_matrices (A, dd.get ());
750  m_rExtrapolation.clear_outer_values (b, dd.get ());
751  }
752 
755  (
756  vector_type & b,
757  const vector_type & u,
759  int type,
760  number time = 0.0
761  )
762  {
763  m_rExtrapolation.clear_outer_values (b, dd.get ());
764  }
765 
767  int type () const {return CT_DIRICHLET;}
768 
769 private:
770 
773 };
774 
776 
783 template <typename TDomain, typename TAlgebra, typename TExtrapolation>
785 : public DomainDiscretizationBase<TDomain, TAlgebra, LSGFGlobAssembler<TDomain, TAlgebra, TExtrapolation> >,
786  public IInterfaceExtrapolation<TDomain, TAlgebra>
787 {
790 
793 
796 
797 public:
799  typedef TDomain domain_type;
800 
802  typedef TAlgebra algebra_type;
803 
805  typedef typename TDomain::grid_type grid_type;
806 
808  typedef typename algebra_type::matrix_type matrix_type;
809 
811  typedef typename algebra_type::vector_type vector_type;
812 
815 
818 
820  static const int dim = TDomain::dim;
821 
822 public:
827  {
828  // register the constraint
830  // set the default the boundary condtions at the interface
831  gass_type::prepare_interface_bc (pApproxSpace);
832  }
833 
836 
838  void set_LSF
839  (
841  )
842  {
843  gass_type::set_LSF (spLSF);
844  }
845 
848  (
849  const char* fct_name,
850  number value
851  )
852  {
854  }
855 
858  (
859  const char* fct_name,
861  )
862  {
864  }
865 
868  (
869  const char* fct_name,
870  number value
871  )
872  {
874  }
875 
878  (
879  const char* fct_name,
881  )
882  {
884  }
885 
886 #ifdef UG_FOR_LUA
887 
890  (
891  const char* fct_name,
892  const char* func_name
893  )
894  {
896  }
897 
900  (
901  const char* fct_name,
902  LuaFunctionHandle func
903  )
904  {
905  set_Dirichlet_on_if_for (fct_name, make_sp(new LuaUserData<number,dim>(func)));
906  }
907 
910  (
911  const char* fct_name,
912  const char* func_name
913  )
914  {
916  }
917 
920  (
921  const char* fct_name,
922  LuaFunctionHandle func
923  )
924  {
925  set_plain_Dirichlet_on_if_for (fct_name, make_sp(new LuaUserData<number,dim>(func)));
926  }
927 
928 #endif
929 
932  (
933  const char* fct_name
934  )
935  {
937  }
938 
941  (
942  const char* subset_names
943  )
944  {
946  }
947 
949  void set_assemble_only_cut (bool b)
950  {
952  }
953 
955  void project_LSF ()
956  {
958  }
959 
961  virtual int check_elem_lsf
962  (
963  size_t n_co,
964  GridObject * pElem,
965  int si,
966  int g_level,
967  bool use_hanging,
968  const MathVector<dim> vCornerCoords [],
969  number time
970  )
971  {
972  return gass_type::check_elem_lsf (n_co, pElem, si, g_level, use_hanging, vCornerCoords, time);
973  }
974 
976  virtual int check_elem_lsf
977  (
978  size_t n_co,
979  GridObject * pElem,
980  int si,
981  bool use_hanging,
982  const MathVector<dim> vCornerCoords [],
983  number time
984  )
985  {
987  int g_level = mg->get_level (pElem);
988  UG_ASSERT (g_level >= 0, "LSGFDomainDiscretization: Grid element without grid level.");
989  return gass_type::check_elem_lsf (n_co, pElem, si, g_level, use_hanging, vCornerCoords, time);
990  }
991 
993  virtual void extrapolate_by_lsf
994  (
995  size_t num_co,
996  size_t base_co,
997  number * u,
998  size_t fct
999  ) const
1000  {
1001  gass_type::extrapolate_by_lsf (num_co, base_co, u, fct);
1002  }
1003 
1005  virtual void extrapolate_by_lsf
1006  (
1007  size_t num_co,
1008  number * u,
1009  size_t fct
1010  ) const
1011  {
1012  gass_type::extrapolate_by_lsf (num_co, u, fct);
1013  }
1014 
1016  virtual bool corner_inside
1017  (
1018  size_t co
1019  ) const
1020  {return gass_type::corner_inside (co);}
1021 
1023  virtual number lsf_at
1024  (
1025  size_t co
1026  ) const
1027  {return gass_type::lsf_at (co);}
1028 
1029 protected:
1032 };
1033 
1034 } // end namespace ug
1035 
1036 #include "dom_disc_embb_impl.h"
1037 
1038 #endif // __H__UG__PLUGINS__D3F__EMBASS__
1039 
1040 /* End of File */
Definition: smart_pointer.h:108
base class for approximation spaces without type of algebra or dof distribution
Definition: approximation_space.h:279
The AssemblingTuner class combines tools to adapt the assembling routine.
Definition: ass_tuner.h:90
Type based UserData.
Definition: user_data.h:501
Definition: dof_distribution.h:51
generic domain discretization implementing the interface
Definition: domain_disc.h:81
void add(SmartPtr< IElemDisc< TDomain > > elem)
adds an element discretization to the assembling process
Definition: domain_disc.h:342
SmartPtr< approx_space_type > m_spApproxSpace
current approximation space
Definition: domain_disc.h:529
The base class for all geometric objects, such as vertices, edges, faces, volumes,...
Definition: grid_base_objects.h:157
TAlgebra algebra_type
Algebra type.
Definition: assemble_interface.h:113
algebra_type::vector_type vector_type
Type of algebra vector.
Definition: constraint_interface.h:72
algebra_type::matrix_type matrix_type
Type of algebra matrix.
Definition: constraint_interface.h:69
Definition: constraint_interface.h:162
ConstSmartPtr< DoFDistribution > dd(const GridLevel &gl) const
returns the level dof distribution
Definition: constraint_interface.h:242
Definition: elem_disc_interface.h:800
Definition: elem_disc_interface.h:760
Base class for the extrapolation over an embedded boundary.
Definition: dom_disc_embb.h:64
virtual int check_elem_lsf(size_t n_co, GridObject *pElem, int si, int g_level, bool use_hanging, const MathVector< dim > vCornerCoords[], number time)
checks whether the element is intersected by the interface, or what, and prepares the data
Definition: dom_disc_embb.h:90
IInterfaceExtrapolation()
Constructor.
Definition: dom_disc_embb.h:83
static const int dim
dimensionality (the World dimension)
Definition: dom_disc_embb.h:80
virtual ~IInterfaceExtrapolation()
Destructor.
Definition: dom_disc_embb.h:86
virtual bool corner_inside(size_t co) const
returns true if the corner is "inside" (use after check_elem_lsf)
Definition: dom_disc_embb.h:134
virtual int check_elem_lsf(size_t n_co, GridObject *pElem, int si, bool use_hanging, const MathVector< dim > vCornerCoords[], number time)
(slower version) checks whether the element is intersected by the interface, or what,...
Definition: dom_disc_embb.h:103
algebra_type::matrix_type matrix_type
matrix type
Definition: dom_disc_embb.h:77
TAlgebra algebra_type
algebra type for the functions to extrapolate
Definition: dom_disc_embb.h:71
virtual void extrapolate_by_lsf(size_t num_co, size_t base_co, number *u, size_t fct) const
extrapolates a component of the solution to the vertices behind the interface (w.r....
Definition: dom_disc_embb.h:115
virtual number lsf_at(size_t co) const
returns the effective value of the LSF at a corner (use after check_elem_lsf)
Definition: dom_disc_embb.h:141
TDomain domain_type
domain type
Definition: dom_disc_embb.h:68
virtual void extrapolate_by_lsf(size_t num_co, number *u, size_t fct) const
extrapolates a component of the solution to the vertices behind the interface (by averaging)
Definition: dom_disc_embb.h:125
algebra_type::vector_type vector_type
vector type (for the functions to extrapolate)
Definition: dom_disc_embb.h:74
a special constraint that sets functions and matrices in the outer subdomain to given values
Definition: dom_disc_embb.h:667
TExtrapolation extrapolation_type
Extrapolation type.
Definition: dom_disc_embb.h:683
algebra_type::matrix_type matrix_type
Matrix type in the algebra.
Definition: dom_disc_embb.h:680
TAlgebra algebra_type
Algebra type.
Definition: dom_disc_embb.h:674
TDomain domain_type
Domain type.
Definition: dom_disc_embb.h:671
int type() const
returns the type of the constraints
Definition: dom_disc_embb.h:767
algebra_type::vector_type vector_type
Vector type in the algebra.
Definition: dom_disc_embb.h:677
extrapolation_type & m_rExtrapolation
Extrapolation in the GF method.
Definition: dom_disc_embb.h:772
void adjust_jacobian(matrix_type &J, const vector_type &u, ConstSmartPtr< DoFDistribution > dd, int type, number time=0.0, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol=SPNULL, const number s_a0=1.0)
sets a unity row for all conductor indices
Definition: dom_disc_embb.h:698
void adjust_rhs(vector_type &b, const vector_type &u, ConstSmartPtr< DoFDistribution > dd, int type, number time=0.0)
sets the dirichlet value in the right-hand side
Definition: dom_disc_embb.h:755
LSGFConstraint(extrapolation_type &rExtrapolation)
class constructor
Definition: dom_disc_embb.h:687
void adjust_solution(vector_type &u, ConstSmartPtr< DoFDistribution > dd, int type, number time=0.0)
sets the value in the solution for all conductor indices
Definition: dom_disc_embb.h:729
void adjust_linear(matrix_type &A, vector_type &b, ConstSmartPtr< DoFDistribution > dd, int type, number time=0.0)
sets unity rows in A and dirichlet values in right-hand side b
Definition: dom_disc_embb.h:741
void adjust_defect(vector_type &d, const vector_type &u, ConstSmartPtr< DoFDistribution > dd, int type, number time=0.0, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol=SPNULL, const std::vector< number > *vScaleMass=NULL, const std::vector< number > *vScaleStiff=NULL)
sets a zero value in the defect for all conductor indices
Definition: dom_disc_embb.h:713
virtual ~LSGFConstraint()
virtual destructor
Definition: dom_disc_embb.h:694
domain discretization for the Level-Set Ghost-Fluid method
Definition: dom_disc_embb.h:787
virtual int check_elem_lsf(size_t n_co, GridObject *pElem, int si, bool use_hanging, const MathVector< dim > vCornerCoords[], number time)
(slow version) checks whether the element is intersected by the interface, or what,...
Definition: dom_disc_embb.h:977
TDomain::grid_type grid_type
Type of the grid.
Definition: dom_disc_embb.h:805
static const int dim
world dimension
Definition: dom_disc_embb.h:820
void set_plain_Dirichlet_on_if_for(const char *fct_name, number value)
sets the 'plain' Dirichlet boundary condition at the interface for a component of the solution
Definition: dom_disc_embb.h:868
void set_Neumann0_on_if_for(const char *fct_name)
sets the Neumann-0 boundary condition at the interface for a component of the solution
Definition: dom_disc_embb.h:932
void project_LSF()
project the level-set function to the coarse levels
Definition: dom_disc_embb.h:955
DomainDiscretizationBase< TDomain, TAlgebra, gass_type > base_type
Type of the base class.
Definition: dom_disc_embb.h:792
void set_Dirichlet_on_if_for(const char *fct_name, SmartPtr< CplUserData< number, dim > > func)
sets the Dirichlet boundary condition at the interface for a component of the solution
Definition: dom_disc_embb.h:858
virtual int check_elem_lsf(size_t n_co, GridObject *pElem, int si, int g_level, bool use_hanging, const MathVector< dim > vCornerCoords[], number time)
checks whether the element is intersected by the interface, or what, and prepares the data
Definition: dom_disc_embb.h:962
gass_type::ls_grid_func_type ls_grid_func_type
Type of the LSF grid functions.
Definition: dom_disc_embb.h:817
virtual void extrapolate_by_lsf(size_t num_co, number *u, size_t fct) const
extrapolates a component of the solution to the vertices behind the interface (by averaging)
Definition: dom_disc_embb.h:1006
TAlgebra algebra_type
Type of algebra.
Definition: dom_disc_embb.h:802
virtual number lsf_at(size_t co) const
returns the effective value of the LSF at a corner (use after check_elem_lsf)
Definition: dom_disc_embb.h:1024
LSGFDomainDiscretization(SmartPtr< approx_space_type > pApproxSpace)
default Constructor
Definition: dom_disc_embb.h:824
virtual ~LSGFDomainDiscretization()
virtual destructor
Definition: dom_disc_embb.h:835
void set_assemble_only_cut(bool b)
set the "assemble only in cut elements" flag
Definition: dom_disc_embb.h:949
virtual bool corner_inside(size_t co) const
returns true if the corner is "inside" (use after check_elem_lsf)
Definition: dom_disc_embb.h:1017
algebra_type::vector_type vector_type
Type of algebra vector.
Definition: dom_disc_embb.h:811
void set_plain_Dirichlet_on_if_for(const char *fct_name, SmartPtr< CplUserData< number, dim > > func)
sets the 'plain' Dirichlet boundary condition at the interface for a component of the solution
Definition: dom_disc_embb.h:878
LSGFConstraint< TDomain, TAlgebra, TExtrapolation > ls_constraint_type
Type of the constraint.
Definition: dom_disc_embb.h:795
virtual void extrapolate_by_lsf(size_t num_co, size_t base_co, number *u, size_t fct) const
extrapolates a component of the solution to the vertices behind the interface (w.r....
Definition: dom_disc_embb.h:994
TDomain domain_type
Type of Domain.
Definition: dom_disc_embb.h:799
void exclude_subsets(const char *subset_names)
excludes a (boundary) subsets from the extrapolation
Definition: dom_disc_embb.h:941
SmartPtr< ls_constraint_type > m_spLSFGFConstraint
the Level-Set Function constraint
Definition: dom_disc_embb.h:1031
ApproximationSpace< TDomain > approx_space_type
Type of approximation space.
Definition: dom_disc_embb.h:814
LSGFGlobAssembler< TDomain, TAlgebra, TExtrapolation > gass_type
Type of the global assembler.
Definition: dom_disc_embb.h:789
void set_LSF(SmartPtr< ls_grid_func_type > spLSF)
set the level-set function and check it
Definition: dom_disc_embb.h:839
void set_Dirichlet_on_if_for(const char *fct_name, number value)
sets the Dirichlet boundary condition at the interface for a component of the solution
Definition: dom_disc_embb.h:848
algebra_type::matrix_type matrix_type
Type of algebra matrix.
Definition: dom_disc_embb.h:808
Global assembler based on the ghost-fluid method with a level-set function.
Definition: dom_disc_embb.h:159
extrapolation_type::ls_grid_func_type ls_grid_func_type
Grid function type for the LSF.
Definition: dom_disc_embb.h:181
void FinishTimestep(const std::vector< IElemDisc< domain_type > * > &vElemDisc, ConstSmartPtr< DoFDistribution > dd, bool bNonRegularGrid, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, ConstSmartPtr< AssemblingTuner< TAlgebra > > spAssTuner)
Definition: dom_disc_embb_impl.h:1870
void set_Dirichlet_on_if_for(SmartPtr< approx_space_type > spApproxSpace, const char *fct_name, SmartPtr< CplUserData< number, dim > > func)
adds a Dirichlet BC with a given value on the interface
Definition: dom_disc_embb.h:510
void PrepareTimestep(const std::vector< IElemDisc< domain_type > * > &vElemDisc, ConstSmartPtr< DoFDistribution > dd, bool bNonRegularGrid, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, number future_time, ConstSmartPtr< AssemblingTuner< TAlgebra > > spAssTuner)
Definition: dom_disc_embb_impl.h:1720
static const int dim
world dimension
Definition: dom_disc_embb.h:184
TDomain domain_type
Domain type.
Definition: dom_disc_embb.h:163
virtual void clear_outer_values(vector_type &d, const DoFDistribution *dd) const
sets the values at the outer vertices to 0
Definition: dom_disc_embb.h:630
TAlgebra algebra_type
Algebra type.
Definition: dom_disc_embb.h:166
void set_Dirichlet_on_if_for(SmartPtr< approx_space_type > spApproxSpace, const char *fct_name, number value)
adds a Dirichlet BC with a given value on the interface
Definition: dom_disc_embb.h:499
void AssembleMassMatrix(const std::vector< IElemDisc< domain_type > * > &vElemDisc, ConstSmartPtr< domain_type > spDomain, ConstSmartPtr< DoFDistribution > dd, TIterator iterBegin, TIterator iterEnd, int si, bool bNonRegularGrid, matrix_type &M, const vector_type &u, ConstSmartPtr< AssemblingTuner< TAlgebra > > spAssTuner)
Definition: dom_disc_embb_impl.h:218
void extrapolate_by_lsf(size_t num_co, size_t base_co, number *u, size_t fct) const
extrapolates a component the solution to the vertices behind the interface (w.r.t....
Definition: dom_disc_embb.h:593
virtual ~LSGFGlobAssembler()
virtual destructor
Definition: dom_disc_embb.h:196
void set_plain_Dirichlet_on_if_for(SmartPtr< approx_space_type > spApproxSpace, const char *fct_name, number value)
adds a "plain" Dirichlet BC with a given value on the interface
Definition: dom_disc_embb.h:532
void project_LSF()
project the level-set function to the coarse levels
Definition: dom_disc_embb.h:568
number lsf_at(size_t co) const
returns the effective value of the LSF at a corner (use after check_elem_lsf)
Definition: dom_disc_embb.h:623
void AssembleDefect(const std::vector< IElemDisc< domain_type > * > &vElemDisc, ConstSmartPtr< domain_type > spDomain, ConstSmartPtr< DoFDistribution > dd, TIterator iterBegin, TIterator iterEnd, int si, bool bNonRegularGrid, vector_type &d, const vector_type &u, ConstSmartPtr< AssemblingTuner< TAlgebra > > spAssTuner)
Definition: dom_disc_embb_impl.h:747
void set_assemble_only_cut(bool b)
set the "assemble only in cut elements" flag
Definition: dom_disc_embb.h:562
static void AssembleRhs(const std::vector< IElemDisc< domain_type > * > &vElemDisc, ConstSmartPtr< domain_type > spDomain, ConstSmartPtr< DoFDistribution > dd, TIterator iterBegin, TIterator iterEnd, int si, bool bNonRegularGrid, vector_type &rhs, const vector_type &u, ConstSmartPtr< AssemblingTuner< TAlgebra > > spAssTuner)
Definition: dom_disc_embb.h:314
void set_Neumann0_on_if_for(SmartPtr< approx_space_type > spApproxSpace, const char *fct_name)
adds a Neumann-0 with on the interface
Definition: dom_disc_embb.h:543
static void AssembleErrorEstimator(const std::vector< IElemError< domain_type > * > &vElemDisc, ConstSmartPtr< domain_type > spDomain, ConstSmartPtr< DoFDistribution > dd, TIterator iterBegin, TIterator iterEnd, int si, bool bNonRegularGrid, const vector_type &u)
Definition: dom_disc_embb.h:428
LSGFGlobAssembler()
class constructor (may not have any arguments!)
Definition: dom_disc_embb.h:193
void prepare_interface_bc(SmartPtr< approx_space_type > spApproxSpace)
prepares the boundary conditions at the interface: sets all them to Dirichlet-0
Definition: dom_disc_embb.h:490
virtual void set_outer_values(vector_type &u, const DoFDistribution *dd, number time)
sets the values at the outer vertices to given values
Definition: dom_disc_embb.h:638
extrapolation_type m_extrapol
the extrapolation at the interface
Definition: dom_disc_embb.h:467
void exclude_subsets(SmartPtr< approx_space_type > spApproxSpace, const char *subset_names)
excludes a (boundary) subsets from the extrapolation
Definition: dom_disc_embb.h:553
void AssembleLinear(const std::vector< IElemDisc< domain_type > * > &vElemDisc, ConstSmartPtr< domain_type > spDomain, ConstSmartPtr< DoFDistribution > dd, TIterator iterBegin, TIterator iterEnd, int si, bool bNonRegularGrid, matrix_type &A, vector_type &rhs, ConstSmartPtr< AssemblingTuner< TAlgebra > > spAssTuner)
Definition: dom_disc_embb_impl.h:1264
virtual void set_outer_matrices(matrix_type &A, const DoFDistribution *dd) const
sets the matrices at outer vertices to identity
Definition: dom_disc_embb.h:647
void InitAllExports(const std::vector< IElemDisc< domain_type > * > &vElemDisc, ConstSmartPtr< DoFDistribution > dd, TIterator iterBegin, TIterator iterEnd, int si, bool bNonRegularGrid, bool bAsTimeDependent)
Definition: dom_disc_embb_impl.h:2036
void extrapolate_by_lsf(size_t num_co, number *u, size_t fct) const
extrapolates a component the solution to the vertices behind the interface (by averaging)
Definition: dom_disc_embb.h:605
void set_plain_Dirichlet_on_if_for(SmartPtr< approx_space_type > spApproxSpace, const char *fct_name, SmartPtr< CplUserData< number, dim > > func)
adds a "plain" Dirichlet BC with a given value on the interface
Definition: dom_disc_embb.h:521
static void AssembleRhs(const std::vector< IElemDisc< domain_type > * > &vElemDisc, ConstSmartPtr< domain_type > spDomain, ConstSmartPtr< DoFDistribution > dd, TIterator iterBegin, TIterator iterEnd, int si, bool bNonRegularGrid, vector_type &rhs, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, const std::vector< number > &vScaleMass, const std::vector< number > &vScaleStiff, ConstSmartPtr< AssemblingTuner< TAlgebra > > spAssTuner)
Definition: dom_disc_embb.h:329
int check_elem_lsf(size_t n_co, GridObject *pElem, int si, int g_level, bool use_hanging, const MathVector< dim > vCornerCoords[], number time)
checks whether the element is intersected by the interface, or what, and prepares the data
Definition: dom_disc_embb.h:578
static void AssembleErrorEstimator(const std::vector< IElemError< domain_type > * > &vElemDisc, ConstSmartPtr< domain_type > spDomain, ConstSmartPtr< DoFDistribution > dd, TIterator iterBegin, TIterator iterEnd, int si, bool bNonRegularGrid, std::vector< number > vScaleMass, std::vector< number > vScaleStiff, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol)
Definition: dom_disc_embb.h:445
extrapolation_type & extrapolation()
returns the extrapolation
Definition: dom_disc_embb.h:574
TExtrapolation extrapolation_type
Extrapolation type.
Definition: dom_disc_embb.h:178
bool m_bAssembleOnlyCut
Definition: dom_disc_embb.h:475
ApproximationSpace< domain_type > approx_space_type
type of approximation space
Definition: dom_disc_embb.h:169
void AssembleStiffnessMatrix(const std::vector< IElemDisc< domain_type > * > &vElemDisc, ConstSmartPtr< domain_type > spDomain, ConstSmartPtr< DoFDistribution > dd, TIterator iterBegin, TIterator iterEnd, int si, bool bNonRegularGrid, matrix_type &A, const vector_type &u, ConstSmartPtr< AssemblingTuner< TAlgebra > > spAssTuner)
Definition: dom_disc_embb_impl.h:60
void FinishTimestepElem(const std::vector< IElemDisc< domain_type > * > &vElemDisc, ConstSmartPtr< domain_type > spDomain, ConstSmartPtr< DoFDistribution > dd, TIterator iterBegin, TIterator iterEnd, int si, bool bNonRegularGrid, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, ConstSmartPtr< AssemblingTuner< TAlgebra > > spAssTuner)
Definition: dom_disc_embb_impl.h:1925
void PrepareTimestepElem(const std::vector< IElemDisc< domain_type > * > &vElemDisc, ConstSmartPtr< domain_type > spDomain, ConstSmartPtr< DoFDistribution > dd, TIterator iterBegin, TIterator iterEnd, int si, bool bNonRegularGrid, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, ConstSmartPtr< AssemblingTuner< TAlgebra > > spAssTuner)
Definition: dom_disc_embb_impl.h:1778
void set_LSF(SmartPtr< ls_grid_func_type > spLSF)
set the level-set function and check it
Definition: dom_disc_embb.h:481
bool corner_inside(size_t co) const
returns true if the corner is "inside" (use after check_elem_lsf)
Definition: dom_disc_embb.h:616
algebra_type::vector_type vector_type
Vector type in the algebra.
Definition: dom_disc_embb.h:172
algebra_type::matrix_type matrix_type
Matrix type in the algebra.
Definition: dom_disc_embb.h:175
void AssembleJacobian(const std::vector< IElemDisc< domain_type > * > &vElemDisc, ConstSmartPtr< domain_type > spDomain, ConstSmartPtr< DoFDistribution > dd, TIterator iterBegin, TIterator iterEnd, int si, bool bNonRegularGrid, matrix_type &J, const vector_type &u, ConstSmartPtr< AssemblingTuner< TAlgebra > > spAssTuner)
Definition: dom_disc_embb_impl.h:376
Factory providing LuaUserData.
Definition: lua_user_data.h:180
static SmartPtr< LuaUserData< TData, dim, TRet > > create(const std::string &name)
Definition: lua_user_data.h:223
Wrapper for sequential matrices to handle them in parallel.
Definition: parallel_matrix.h:65
time series of solutions and corresponding time point
Definition: solution_time_series.h:59
static const int dim
const NullSmartPtr SPNULL
The equivalent to NULL for smart pointers.
Definition: smart_pointer.h:90
#define UG_ASSERT(expr, msg)
Definition: assert.h:70
#define UG_THROW(msg)
Definition: error.h:57
double number
Definition: types.h:124
CPUAlgebra::vector_type vector_type
the ug namespace
@ CT_DIRICHLET
Definition: ass_tuner.h:59
function func(x, y, z, t, si)
SmartPtr< T, FreePolicy > make_sp(T *inst)
returns a SmartPtr for the passed raw pointer
Definition: smart_pointer.h:836
function ProblemDisc new(problemDesc, dom)