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 
412 // Error estimators: Not implemented for the ghost-fluid method
414 
415 public:
416 
417  template <typename TElem, typename TIterator>
418  static void
420  (
421  const std::vector<IElemError<domain_type>*>& vElemDisc,
424  TIterator iterBegin,
425  TIterator iterEnd,
426  int si,
427  bool bNonRegularGrid,
428  const vector_type& u
429  )
430  {
431  UG_THROW ("AssembleErrorEstimator: No error estimator implemented for the Ghost-Fluid method.");
432  }
433 
434  template <typename TElem, typename TIterator>
435  static void
437  (
438  const std::vector<IElemError<domain_type>*>& vElemDisc,
441  TIterator iterBegin,
442  TIterator iterEnd,
443  int si,
444  bool bNonRegularGrid,
445  std::vector<number> vScaleMass,
446  std::vector<number> vScaleStiff,
448  )
449  {
450  UG_THROW ("AssembleErrorEstimator: No error estimator implemented for the Ghost-Fluid method.");
451  }
452 
454 // Data
456 
457 private:
458 
460 
468 
469 public:
470 
472  void set_LSF
473  (
475  )
476  {
477  m_extrapol.set_LSF (spLSF);
478  }
479 
482  (
483  SmartPtr<approx_space_type> spApproxSpace
484  )
485  {
486  m_extrapol.prepare_interface_bc (spApproxSpace);
487  }
488 
491  (
492  SmartPtr<approx_space_type> spApproxSpace,
493  const char* fct_name,
494  number value
495  )
496  {
497  m_extrapol.set_Dirichlet_for (spApproxSpace, fct_name, value);
498  }
499 
502  (
503  SmartPtr<approx_space_type> spApproxSpace,
504  const char* fct_name,
506  )
507  {
508  m_extrapol.set_Dirichlet_for (spApproxSpace, fct_name, func);
509  }
510 
513  (
514  SmartPtr<approx_space_type> spApproxSpace,
515  const char* fct_name,
517  )
518  {
519  m_extrapol.set_plain_Dirichlet_for (spApproxSpace, fct_name, func);
520  }
521 
524  (
525  SmartPtr<approx_space_type> spApproxSpace,
526  const char* fct_name,
527  number value
528  )
529  {
530  m_extrapol.set_plain_Dirichlet_for (spApproxSpace, fct_name, value);
531  }
532 
535  (
536  SmartPtr<approx_space_type> spApproxSpace,
537  const char* fct_name
538  )
539  {
540  m_extrapol.set_Neumann0_for (spApproxSpace, fct_name);
541  }
542 
545  (
546  SmartPtr<approx_space_type> spApproxSpace,
547  const char* subset_names
548  )
549  {
550  m_extrapol.exclude_subsets (spApproxSpace, subset_names);
551  }
552 
554  void set_assemble_only_cut (bool b)
555  {
556  m_bAssembleOnlyCut = b;
557  }
558 
560  void project_LSF ()
561  {
562  m_extrapol.project_LSF ();
563  }
564 
567 
570  (
571  size_t n_co,
572  GridObject * pElem,
573  int si,
574  int g_level,
575  bool use_hanging,
576  const MathVector<dim> vCornerCoords [],
577  number time
578  )
579  {
580  return m_extrapol.check_elem_lsf (n_co, pElem, si, g_level, use_hanging, vCornerCoords, time);
581  }
582 
585  (
586  size_t num_co,
587  size_t base_co,
588  number * u,
589  size_t fct
590  ) const
591  {
592  m_extrapol.extrapolate_by_lsf (num_co, base_co, u, fct);
593  }
594 
597  (
598  size_t num_co,
599  number * u,
600  size_t fct
601  ) const
602  {
603  m_extrapol.extrapolate_by_lsf (num_co, u, fct);
604  }
605 
608  (
609  size_t co
610  ) const
611  {return m_extrapol.corner_inside (co);}
612 
615  (
616  size_t co
617  ) const
618  {return m_extrapol.lsf_at (co);}
619 
621  virtual void clear_outer_values
622  (
623  vector_type & d,
624  const DoFDistribution * dd
625  ) const
626  {m_extrapol.clear_outer_values (d, dd);}
627 
629  virtual void set_outer_values
630  (
631  vector_type & u,
632  const DoFDistribution * dd,
633  number time
634  )
635  {m_extrapol.set_outer_values (u, dd, time);}
636 
638  virtual void set_outer_matrices
639  (
640  matrix_type & A,
641  const DoFDistribution * dd
642  ) const
643  {m_extrapol.set_outer_matrices (A, dd);}
644 
645 }; // class LSGFGlobAssembler
646 
648 
656 template <typename TDomain, typename TAlgebra, typename TExtrapolation>
658 : public IDomainConstraint<TDomain, TAlgebra>
659 {
660 public:
661 
663  typedef TDomain domain_type;
664 
666  typedef TAlgebra algebra_type;
667 
669  typedef typename algebra_type::vector_type vector_type;
670 
672  typedef typename algebra_type::matrix_type matrix_type;
673 
675  typedef TExtrapolation extrapolation_type;
676 
679  (
680  extrapolation_type & rExtrapolation
681  )
682  : m_rExtrapolation (rExtrapolation)
683  {}
684 
686  virtual ~LSGFConstraint () {};
687 
690  (
691  matrix_type & J,
692  const vector_type & u,
694  int type,
695  number time = 0.0,
697  const number s_a0 = 1.0
698  )
699  {
700  m_rExtrapolation.set_outer_matrices (J, dd.get ());
701  }
702 
705  (
706  vector_type & d,
707  const vector_type & u,
709  int type,
710  number time = 0.0,
712  const std::vector<number> * vScaleMass = NULL,
713  const std::vector<number> * vScaleStiff = NULL
714  )
715  {
716  m_rExtrapolation.clear_outer_values (d, dd.get ());
717  }
718 
721  (
722  vector_type & u,
724  int type,
725  number time = 0.0
726  )
727  {
728  m_rExtrapolation.set_outer_values (u, dd.get (), time);
729  }
730 
733  (
734  matrix_type & A,
735  vector_type & b,
737  int type,
738  number time = 0.0
739  )
740  { // Note that this function is not really used, so it needs not to be optimal.
741  m_rExtrapolation.set_outer_matrices (A, dd.get ());
742  m_rExtrapolation.clear_outer_values (b, dd.get ());
743  }
744 
747  (
748  vector_type & b,
749  const vector_type & u,
751  int type,
752  number time = 0.0
753  )
754  {
755  m_rExtrapolation.clear_outer_values (b, dd.get ());
756  }
757 
759  int type () const {return CT_DIRICHLET;}
760 
761 private:
762 
765 };
766 
768 
775 template <typename TDomain, typename TAlgebra, typename TExtrapolation>
777 : public DomainDiscretizationBase<TDomain, TAlgebra, LSGFGlobAssembler<TDomain, TAlgebra, TExtrapolation> >,
778  public IInterfaceExtrapolation<TDomain, TAlgebra>
779 {
782 
785 
788 
789 public:
791  typedef TDomain domain_type;
792 
794  typedef TAlgebra algebra_type;
795 
797  typedef typename TDomain::grid_type grid_type;
798 
800  typedef typename algebra_type::matrix_type matrix_type;
801 
803  typedef typename algebra_type::vector_type vector_type;
804 
807 
810 
812  static const int dim = TDomain::dim;
813 
814 public:
819  {
820  // register the constraint
822  // set the default the boundary condtions at the interface
823  gass_type::prepare_interface_bc (pApproxSpace);
824  }
825 
828 
830  void set_LSF
831  (
833  )
834  {
835  gass_type::set_LSF (spLSF);
836  }
837 
840  (
841  const char* fct_name,
842  number value
843  )
844  {
846  }
847 
850  (
851  const char* fct_name,
853  )
854  {
856  }
857 
860  (
861  const char* fct_name,
862  number value
863  )
864  {
866  }
867 
870  (
871  const char* fct_name,
873  )
874  {
876  }
877 
878 #ifdef UG_FOR_LUA
879 
882  (
883  const char* fct_name,
884  const char* func_name
885  )
886  {
888  }
889 
892  (
893  const char* fct_name,
894  LuaFunctionHandle func
895  )
896  {
897  set_Dirichlet_on_if_for (fct_name, make_sp(new LuaUserData<number,dim>(func)));
898  }
899 
902  (
903  const char* fct_name,
904  const char* func_name
905  )
906  {
908  }
909 
912  (
913  const char* fct_name,
914  LuaFunctionHandle func
915  )
916  {
917  set_plain_Dirichlet_on_if_for (fct_name, make_sp(new LuaUserData<number,dim>(func)));
918  }
919 
920 #endif
921 
924  (
925  const char* fct_name
926  )
927  {
929  }
930 
933  (
934  const char* subset_names
935  )
936  {
938  }
939 
941  void set_assemble_only_cut (bool b)
942  {
944  }
945 
947  void project_LSF ()
948  {
950  }
951 
953  virtual int check_elem_lsf
954  (
955  size_t n_co,
956  GridObject * pElem,
957  int si,
958  int g_level,
959  bool use_hanging,
960  const MathVector<dim> vCornerCoords [],
961  number time
962  )
963  {
964  return gass_type::check_elem_lsf (n_co, pElem, si, g_level, use_hanging, vCornerCoords, time);
965  }
966 
968  virtual int check_elem_lsf
969  (
970  size_t n_co,
971  GridObject * pElem,
972  int si,
973  bool use_hanging,
974  const MathVector<dim> vCornerCoords [],
975  number time
976  )
977  {
979  int g_level = mg->get_level (pElem);
980  UG_ASSERT (g_level >= 0, "LSGFDomainDiscretization: Grid element without grid level.");
981  return gass_type::check_elem_lsf (n_co, pElem, si, g_level, use_hanging, vCornerCoords, time);
982  }
983 
985  virtual void extrapolate_by_lsf
986  (
987  size_t num_co,
988  size_t base_co,
989  number * u,
990  size_t fct
991  ) const
992  {
993  gass_type::extrapolate_by_lsf (num_co, base_co, u, fct);
994  }
995 
997  virtual void extrapolate_by_lsf
998  (
999  size_t num_co,
1000  number * u,
1001  size_t fct
1002  ) const
1003  {
1004  gass_type::extrapolate_by_lsf (num_co, u, fct);
1005  }
1006 
1008  virtual bool corner_inside
1009  (
1010  size_t co
1011  ) const
1012  {return gass_type::corner_inside (co);}
1013 
1015  virtual number lsf_at
1016  (
1017  size_t co
1018  ) const
1019  {return gass_type::lsf_at (co);}
1020 
1021 protected:
1024 };
1025 
1026 } // end namespace ug
1027 
1028 #include "dom_disc_embb_impl.h"
1029 
1030 #endif // __H__UG__PLUGINS__D3F__EMBASS__
1031 
1032 /* 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:333
SmartPtr< approx_space_type > m_spApproxSpace
current approximation space
Definition: domain_disc.h:520
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:659
TExtrapolation extrapolation_type
Extrapolation type.
Definition: dom_disc_embb.h:675
algebra_type::matrix_type matrix_type
Matrix type in the algebra.
Definition: dom_disc_embb.h:672
TAlgebra algebra_type
Algebra type.
Definition: dom_disc_embb.h:666
TDomain domain_type
Domain type.
Definition: dom_disc_embb.h:663
int type() const
returns the type of the constraints
Definition: dom_disc_embb.h:759
algebra_type::vector_type vector_type
Vector type in the algebra.
Definition: dom_disc_embb.h:669
extrapolation_type & m_rExtrapolation
Extrapolation in the GF method.
Definition: dom_disc_embb.h:764
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:690
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:747
LSGFConstraint(extrapolation_type &rExtrapolation)
class constructor
Definition: dom_disc_embb.h:679
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:721
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:733
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:705
virtual ~LSGFConstraint()
virtual destructor
Definition: dom_disc_embb.h:686
domain discretization for the Level-Set Ghost-Fluid method
Definition: dom_disc_embb.h:779
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:969
TDomain::grid_type grid_type
Type of the grid.
Definition: dom_disc_embb.h:797
static const int dim
world dimension
Definition: dom_disc_embb.h:812
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:860
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:924
void project_LSF()
project the level-set function to the coarse levels
Definition: dom_disc_embb.h:947
DomainDiscretizationBase< TDomain, TAlgebra, gass_type > base_type
Type of the base class.
Definition: dom_disc_embb.h:784
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:850
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:954
gass_type::ls_grid_func_type ls_grid_func_type
Type of the LSF grid functions.
Definition: dom_disc_embb.h:809
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:998
TAlgebra algebra_type
Type of algebra.
Definition: dom_disc_embb.h:794
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:1016
LSGFDomainDiscretization(SmartPtr< approx_space_type > pApproxSpace)
default Constructor
Definition: dom_disc_embb.h:816
virtual ~LSGFDomainDiscretization()
virtual destructor
Definition: dom_disc_embb.h:827
void set_assemble_only_cut(bool b)
set the "assemble only in cut elements" flag
Definition: dom_disc_embb.h:941
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:1009
algebra_type::vector_type vector_type
Type of algebra vector.
Definition: dom_disc_embb.h:803
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:870
LSGFConstraint< TDomain, TAlgebra, TExtrapolation > ls_constraint_type
Type of the constraint.
Definition: dom_disc_embb.h:787
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:986
TDomain domain_type
Type of Domain.
Definition: dom_disc_embb.h:791
void exclude_subsets(const char *subset_names)
excludes a (boundary) subsets from the extrapolation
Definition: dom_disc_embb.h:933
SmartPtr< ls_constraint_type > m_spLSFGFConstraint
the Level-Set Function constraint
Definition: dom_disc_embb.h:1023
ApproximationSpace< TDomain > approx_space_type
Type of approximation space.
Definition: dom_disc_embb.h:806
LSGFGlobAssembler< TDomain, TAlgebra, TExtrapolation > gass_type
Type of the global assembler.
Definition: dom_disc_embb.h:781
void set_LSF(SmartPtr< ls_grid_func_type > spLSF)
set the level-set function and check it
Definition: dom_disc_embb.h:831
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:840
algebra_type::matrix_type matrix_type
Type of algebra matrix.
Definition: dom_disc_embb.h:800
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:502
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:622
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:491
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:585
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:524
void project_LSF()
project the level-set function to the coarse levels
Definition: dom_disc_embb.h:560
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:615
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:554
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:535
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:420
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:482
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:630
extrapolation_type m_extrapol
the extrapolation at the interface
Definition: dom_disc_embb.h:459
void exclude_subsets(SmartPtr< approx_space_type > spApproxSpace, const char *subset_names)
excludes a (boundary) subsets from the extrapolation
Definition: dom_disc_embb.h:545
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:639
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:597
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:513
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:570
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:437
extrapolation_type & extrapolation()
returns the extrapolation
Definition: dom_disc_embb.h:566
TExtrapolation extrapolation_type
Extrapolation type.
Definition: dom_disc_embb.h:178
bool m_bAssembleOnlyCut
Definition: dom_disc_embb.h:467
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:473
bool corner_inside(size_t co) const
returns true if the corner is "inside" (use after check_elem_lsf)
Definition: dom_disc_embb.h:608
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)