GetFEM  5.4.1
getfem_models.h
Go to the documentation of this file.
1 /* -*- c++ -*- (enables emacs c++ mode) */
2 /*===========================================================================
3 
4  Copyright (C) 2009-2020 Yves Renard
5 
6  This file is a part of GetFEM
7 
8  GetFEM is free software; you can redistribute it and/or modify it
9  under the terms of the GNU Lesser General Public License as published
10  by the Free Software Foundation; either version 3 of the License, or
11  (at your option) any later version along with the GCC Runtime Library
12  Exception either version 3.1 or (at your option) any later version.
13  This program is distributed in the hope that it will be useful, but
14  WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
15  or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
16  License and GCC Runtime Library Exception for more details.
17  You should have received a copy of the GNU Lesser General Public License
18  along with this program; if not, write to the Free Software Foundation,
19  Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301, USA.
20 
21  As a special exception, you may use this file as it is a part of a free
22  software library without restriction. Specifically, if other files
23  instantiate templates or use macros or inline functions from this file,
24  or you compile this file and link it with other files to produce an
25  executable, this file does not by itself cause the resulting executable
26  to be covered by the GNU Lesser General Public License. This exception
27  does not however invalidate any other reasons why the executable file
28  might be covered by the GNU Lesser General Public License.
29 
30 ===========================================================================*/
31 
32 /**
33  @file getfem_models.h
34  @author Yves Renard <Yves.Renard@insa-lyon.fr>
35  @date March 21, 2009.
36  @brief Model representation in Getfem.
37 */
38 
39 #ifndef GETFEM_MODELS_H__
40 #define GETFEM_MODELS_H__
41 
43 #include "getfem_assembling.h"
45 #include "getfem_im_data.h"
46 
47 namespace getfem {
48 
50  /** type of pointer on a brick */
51  typedef std::shared_ptr<const virtual_brick> pbrick;
52 
53  class virtual_dispatcher;
54  typedef std::shared_ptr<const virtual_dispatcher> pdispatcher;
55 
56  class virtual_time_scheme;
57  typedef std::shared_ptr<const virtual_time_scheme> ptime_scheme;
58 
59  // Event management : The model has to react when something has changed in
60  // the context and ask for corresponding (linear) bricks to recompute
61  // some terms.
62  // For the moment two events are taken into account
63  // - Change in a mesh_fem
64  // - Change in the data of a variable
65  // For this, a brick has to declare on which variable it depends and
66  // on which data. When a linear brick depend on a variable, the
67  // recomputation is done when the eventual corresponding mesh_fem
68  // is changed (or the size of the variable for a fixed size variable).
69  // When a linear brick depend on a data, the recomputation is done
70  // when the corresponding vector value is changed. If a variable is used
71  // as a data, it has to be declared as a data by the brick.
72  // A nonlinear brick is recomputed at each assembly of the tangent system.
73  // Remember this behavior when some changed are done on the variable
74  // and/or data.
75  // The change on a mesh_im is not taken into account for the moment.
76  // The different versions of the variables is not taken into account
77  // separately.
78 
79 
80 
81 
82  //=========================================================================
83  //
84  // Model object.
85  //
86  //=========================================================================
87 
88  // For backward compatibility with version 3.0
89  typedef model_real_plain_vector modeling_standard_plain_vector;
91  typedef model_real_sparse_matrix modeling_standard_sparse_matrix;
92  typedef model_complex_plain_vector modeling_standard_complex_plain_vector;
94  typedef model_complex_sparse_matrix modeling_standard_complex_sparse_matrix;
95 
96 
97  /** A prefix to refer to the previous version of a variable*/
98  const auto PREFIX_OLD = std::string{"Old_"};
99  const auto PREFIX_OLD_LENGTH = 4;
100 
101  /** Does the variable have Old_ prefix*/
102  bool is_old(const std::string &name);
103 
104  /** Strip the variable name from prefix Old_ if it has one*/
105  std::string no_old_prefix_name(const std::string &name);
106 
107  std::string sup_previous_and_dot_to_varname(std::string v);
108 
109  /** ``Model'' variables store the variables, the data and the
110  description of a model. This includes the global tangent matrix, the
111  right hand side and the constraints. There are two kinds of models, the
112  ``real'' and the ``complex'' models.
113  */
114  class APIDECL model : public context_dependencies,
115  virtual public dal::static_stored_object {
116 
117  protected:
118 
119  // State variables of the model
120  bool complex_version;
121  bool is_linear_;
122  bool is_symmetric_;
123  bool is_coercive_;
124  mutable model_real_sparse_matrix
125  rTM, // tangent matrix (only primary variables), real version
126  internal_rTM; // coupling matrix between internal and primary vars (no empty rows)
127  mutable model_complex_sparse_matrix cTM; // tangent matrix, complex version
128  mutable model_real_plain_vector
129  rrhs, // residual vector of primary variables (after condensation)
130  full_rrhs, // residual vector of primary and internal variables (pre-condensation)
131  internal_sol; // partial solution for internal variables (after condensation)
132  mutable model_complex_plain_vector crhs;
133  mutable bool act_size_to_be_done;
134  dim_type leading_dim;
135  getfem::lock_factory locks_;
136 
137  // Variables and parameters of the model
138 
139  enum var_description_filter {
140  VDESCRFILTER_NO = 0, // Variable being directly the dofs of a given fem
141  VDESCRFILTER_REGION = 1, /* Variable being the dofs of a fem on a mesh
142  * region (uses mf.dof_on_region). */
143  VDESCRFILTER_INFSUP = 2, /* Variable being the dofs of a fem on a mesh
144  * region with an additional filter on a mass
145  * matrix with respect to another fem. */
146  VDESCRFILTER_CTERM = 4, /* Variable being the dofs of a fem on a mesh
147  * region with an additional filter with the
148  * coupling term with respect to another variable.*/
149  VDESCRFILTER_REGION_CTERM = 5, /* both region and cterm. */
150  }; // INFSUP and CTERM are incompatible
151 
152  struct var_description {
153 
154  bool is_variable; // This is a variable or a parameter.
155  bool is_disabled; // For a variable, to be solved or not
156  bool is_complex; // The variable is complex numbers
157  bool is_affine_dependent; // The variable depends in an affine way
158  // to another variable.
159  bool is_internal; // An internal variable defined on integration
160  // points, condensed out of the global system.
161  bool is_fem_dofs; // The variable is the dofs of a fem
162  size_type n_iter; // Number of versions of the variable stored.
163  size_type n_temp_iter; // Number of additional temporary versions
164  size_type default_iter; // default iteration number.
165 
166  ptime_scheme ptsc; // For optional time integration scheme
167 
168  var_description_filter filter; // Version of an optional filter
169  // on the dofs
170  size_type filter_region; // Optional region for the filter
171  std::string filter_var; // Optional variable name for the filter
172  mesh_im const *filter_mim; // Optional integration method for the filter
173 
174  // fem or im_data description of the variable
175  mesh_fem const *mf; // Main fem of the variable.
176  ppartial_mesh_fem partial_mf; // Filter with respect to mf.
177  im_data const *imd; // im data description
178 
179  bgeot::multi_index qdims; // For data having a qdim != of the fem
180  // (dim per dof for dof data)
181  // and for constant variables.
182  gmm::uint64_type v_num;
183  std::vector<gmm::uint64_type> v_num_data;
184 
185  gmm::sub_interval I; // For a variable : indices on the whole system.
186  // For an affine dependent variable, should be the same as the
187  // original variable
188  std::vector<model_real_plain_vector> real_value;
189  std::vector<model_complex_plain_vector> complex_value;
190  std::vector<gmm::uint64_type> v_num_var_iter;
191  std::vector<gmm::uint64_type> v_num_iter;
192 
193  // For affine dependent variables
194  model_real_plain_vector affine_real_value;
195  model_complex_plain_vector affine_complex_value;
196  scalar_type alpha; // Factor for the affine dependent variables
197  std::string org_name; // Name of the original variable for affine
198  // dependent variables
199 
200  var_description(bool is_var = false, bool is_compl = false,
201  mesh_fem const *mf_ = 0, im_data const *imd_ = 0,
202  size_type n_it = 1,
203  var_description_filter filter_ = VDESCRFILTER_NO,
204  size_type filter_reg = size_type(-1),
205  const std::string &filter_var_ = std::string(""),
206  mesh_im const *filter_mim_ = 0)
207  : is_variable(is_var), is_disabled(false), is_complex(is_compl),
208  is_affine_dependent(false), is_internal(false),
209  is_fem_dofs(mf_ != 0),
210  n_iter(std::max(size_type(1), n_it)), n_temp_iter(0),
211  default_iter(0), ptsc(0),
212  filter(filter_), filter_region(filter_reg), filter_var(filter_var_),
213  filter_mim(filter_mim_), mf(mf_), imd(imd_), qdims(),
214  v_num(0), v_num_data(n_iter, act_counter()), I(0,0), alpha(1)
215  {
216  if (filter != VDESCRFILTER_NO && mf != 0)
217  partial_mf = std::make_shared<partial_mesh_fem>(*mf);
218  // v_num_data = v_num;
219  if (qdims.size() == 0) qdims.push_back(1);
220  GMM_ASSERT1(qdim(), "Attempt to create a null size variable");
221  }
222 
223  size_type qdim() const { return qdims.total_size(); }
224 
225  // add a temporary version for time integration schemes. Automatically
226  // set the default iter to it. id_num is an identifier. Do not add
227  // the version if a temporary already exist with this identifier.
228  size_type add_temporary(gmm::uint64_type id_num);
229 
230  void clear_temporaries();
231 
232  const mesh_fem &associated_mf() const {
233  GMM_ASSERT1(is_fem_dofs, "This variable is not linked to a fem");
234  return (filter == VDESCRFILTER_NO) ? *mf : *partial_mf;
235  }
236 
237  const mesh_fem *passociated_mf() const {
238  if (is_fem_dofs)
239  return (filter == VDESCRFILTER_NO || partial_mf.get() == 0)
240  ? mf : partial_mf.get();
241  return 0;
242  }
243 
244  size_type size() const // Should control that the variable is
245  // indeed initialized by actualize_sizes() ...
246  { return is_complex ? complex_value[0].size() : real_value[0].size(); }
247 
248  void set_size();
249  }; // struct var_description
250 
251  public:
252 
253  typedef std::vector<std::string> varnamelist;
254  typedef std::vector<const mesh_im *> mimlist;
255  typedef std::vector<model_real_sparse_matrix> real_matlist;
256  typedef std::vector<model_complex_sparse_matrix> complex_matlist;
257  typedef std::vector<model_real_plain_vector> real_veclist;
258  typedef std::vector<model_complex_plain_vector> complex_veclist;
259 
260  struct term_description {
261  bool is_matrix_term; // tangent matrix term or rhs term.
262  bool is_symmetric; // Term have to be symmetrized.
263  bool is_global; // Specific global term for highly coupling bricks
264  std::string var1, var2;
265 
266  term_description(const std::string &v)
267  : is_matrix_term(false), is_symmetric(false),
268  is_global(false), var1(sup_previous_and_dot_to_varname(v)) {}
269  term_description(const std::string &v1, const std::string &v2,
270  bool issym)
271  : is_matrix_term(true), is_symmetric(issym), is_global(false),
272  var1(sup_previous_and_dot_to_varname(v1)), var2(v2) {}
273  term_description(bool ism, bool issym)
274  : is_matrix_term(ism), is_symmetric(issym), is_global(true) {}
275  };
276 
277  typedef std::vector<term_description> termlist;
278 
279  enum build_version {
280  BUILD_RHS = 1,
281  BUILD_MATRIX = 2,
282  BUILD_ALL = 3,
283  BUILD_ON_DATA_CHANGE = 4,
284  BUILD_WITH_LIN = 8, // forced calculation of linear terms
285  BUILD_RHS_WITH_LIN = 9, // = BUILD_RHS | BUILD_WITH_LIN_RHS
286  BUILD_WITH_INTERNAL = 16,
287  BUILD_RHS_WITH_INTERNAL = 17, // = BUILD_RHS | BUILD_WITH_INTERNAL
288  BUILD_MATRIX_CONDENSED = 18, // = BUILD_MATRIX | BUILD_WITH_INTERNAL
289  BUILD_ALL_CONDENSED = 19, // = BUILD_ALL | BUILD_WITH_INTERNAL
290  };
291 
292  protected:
293 
294  // rmatlist and cmatlist could had been csc_matrix vectors to reduce the
295  // amount of memory (but this would require an additional copy).
296  struct brick_description {
297  mutable bool terms_to_be_computed;
298  mutable gmm::uint64_type v_num;
299  pbrick pbr; // brick pointer
300  pdispatcher pdispatch; // Optional dispatcher
301  size_type nbrhs; // Additional rhs for dispatcher.
302  varnamelist vlist; // List of variables used by the brick.
303  varnamelist dlist; // List of data used by the brick.
304  termlist tlist; // List of terms build by the brick
305  mimlist mims; // List of integration methods.
306  size_type region; // Optional region size_type(-1) for all.
307  bool is_update_brick; // Flag for declaring special type of brick
308  // with no contributions.
309  mutable scalar_type external_load; // External load computed in assembly
310 
311  mutable model_real_plain_vector coeffs;
312  mutable scalar_type matrix_coeff = scalar_type(0);
313  // Matrices the brick has to fill in (real version)
314  mutable real_matlist rmatlist;
315  // Rhs the brick has to fill in (real version)
316  mutable std::vector<real_veclist> rveclist;
317  // additional rhs for symmetric terms (real version)
318  mutable std::vector<real_veclist> rveclist_sym;
319  // Matrices the brick has to fill in (complex version)
320  mutable complex_matlist cmatlist;
321  // Rhs the brick has to fill in (complex version)
322  mutable std::vector<complex_veclist> cveclist;
323  // additional rhs for symmetric terms (complex version)
324  mutable std::vector<complex_veclist> cveclist_sym;
325 
326  brick_description() : v_num(0) {}
327 
328  brick_description(pbrick p, const varnamelist &vl,
329  const varnamelist &dl, const termlist &tl,
330  const mimlist &mms, size_type reg)
331  : terms_to_be_computed(true), v_num(0), pbr(p), pdispatch(0), nbrhs(1),
332  vlist(vl), dlist(dl), tlist(tl), mims(mms), region(reg),
333  is_update_brick(false), external_load(0),
334  rveclist(1), rveclist_sym(1), cveclist(1),
335  cveclist_sym(1) { }
336  };
337 
338  typedef std::map<std::string, var_description> VAR_SET;
339  mutable VAR_SET variables; // Variables list of the model
340  std::vector<brick_description> bricks; // Bricks list of the model
341  dal::bit_vector valid_bricks, active_bricks;
342  std::map<std::string, pinterpolate_transformation> transformations;
343  std::map<std::string, pelementary_transformation> elem_transformations;
344  std::map<std::string, psecondary_domain> secondary_domains;
345 
346  // Structure dealing with time integration scheme
347  int time_integration; // 0 : no, 1 : time step, 2 : init
348  bool init_step;
349  scalar_type time_step; // Time step (dt) for time integration schemes
350  scalar_type init_time_step; // Time step for initialization of derivatives
351 
352  // Structure dealing with simple dof constraints
353  typedef std::map<size_type, scalar_type> real_dof_constraints_var;
354  typedef std::map<size_type, complex_type> complex_dof_constraints_var;
355  mutable std::map<std::string, real_dof_constraints_var>
356  real_dof_constraints;
357  mutable std::map<std::string, complex_dof_constraints_var>
358  complex_dof_constraints;
359  void clear_dof_constraints()
360  { real_dof_constraints.clear(); complex_dof_constraints.clear(); }
361 
362  // Structure dealing with nonlinear expressions
363  struct gen_expr {
364  std::string expr;
365  const mesh_im &mim;
366  size_type region;
367  std::string secondary_domain;
368  gen_expr(const std::string &expr_, const mesh_im &mim_,
369  size_type region_, const std::string &secdom)
370  : expr(expr_), mim(mim_), region(region_), secondary_domain(secdom) {}
371  };
372 
373  // Structure for assignment in assembly
374  struct assignement_desc {
375  std::string varname;
376  std::string expr;
377  size_type region;
378  bool before;
379  size_type order;
380  };
381 
382  std::list<assignement_desc> assignments;
383 
384  mutable std::list<gen_expr> generic_expressions;
385 
386  // Groups of variables for interpolation on different meshes
387  // generic assembly
388  std::map<std::string, std::vector<std::string> > variable_groups;
389 
390  ga_macro_dictionary macro_dict;
391 
392 
393  virtual void actualize_sizes() const;
394  bool check_name_validity(const std::string &name, bool assert=true) const;
395  void brick_init(size_type ib, build_version version,
396  size_type rhs_ind = 0) const;
397 
398  void init() { complex_version = false; act_size_to_be_done = false; }
399 
400  void resize_global_system() const;
401 
402  //to be performed after to_variables is called.
403  virtual void post_to_variables_step();
404 
405  scalar_type approx_external_load_; // Computed by assembly procedure
406  // with BUILD_RHS option.
407 
408  VAR_SET::const_iterator find_variable(const std::string &name) const;
409  const var_description &variable_description(const std::string &name) const;
410 
411  public:
412 
413  void add_generic_expression(const std::string &expr, const mesh_im &mim,
414  size_type region,
415  const std::string &secondary_domain = "") const {
416  generic_expressions.push_back(gen_expr(expr, mim, region,
417  secondary_domain));
418  }
419  void add_external_load(size_type ib, scalar_type e) const
420  { bricks[ib].external_load = e; }
421  scalar_type approx_external_load() { return approx_external_load_; }
422  // call the brick if necessary
423  void update_brick(size_type ib, build_version version) const;
424  void linear_brick_add_to_rhs(size_type ib, size_type ind_data,
425  size_type n_iter) const;
426  void update_affine_dependent_variables();
427  void brick_call(size_type ib, build_version version,
428  size_type rhs_ind = 0) const;
429  model_real_plain_vector &rhs_coeffs_of_brick(size_type ib) const
430  { return bricks[ib].coeffs; }
431  scalar_type &matrix_coeff_of_brick(size_type ib) const
432  { return bricks[ib].matrix_coeff; }
433  bool is_var_newer_than_brick(const std::string &varname,
434  size_type ib, size_type niter = size_type(-1)) const;
435  bool is_var_mf_newer_than_brick(const std::string &varname,
436  size_type ib) const;
437  bool is_mim_newer_than_brick(const mesh_im &mim,
438  size_type ib) const;
439 
440  pbrick brick_pointer(size_type ib) const {
441  GMM_ASSERT1(valid_bricks[ib], "Inexistent brick");
442  return bricks[ib].pbr;
443  }
444 
445  void variable_list(varnamelist &vl) const
446  { for (const auto &v : variables) vl.push_back(v.first); }
447 
448  void define_variable_group(const std::string &group_name,
449  const std::vector<std::string> &nl);
450  bool variable_group_exists(const std::string &group_name) const
451  { return variable_groups.count(group_name) > 0; }
452 
453  const std::vector<std::string> &
454  variable_group(const std::string &group_name) const {
455  GMM_ASSERT1(variable_group_exists(group_name),
456  "Undefined variable group " << group_name);
457  return (variable_groups.find(group_name))->second;
458  }
459 
460  void clear_assembly_assignments(void) { assignments.clear(); }
461  void add_assembly_assignments(const std::string &dataname,
462  const std::string &expr,
463  size_type rg = size_type(-1),
464  size_type order = 1,
465  bool before = false);
466 
467  /* function to be called by Dirichlet bricks */
468  void add_real_dof_constraint(const std::string &varname, size_type dof,
469  scalar_type val) const
470  { real_dof_constraints[varname][dof] = val; }
471  /* function to be called by Dirichlet bricks */
472  void add_complex_dof_constraint(const std::string &varname, size_type dof,
473  complex_type val) const
474  { complex_dof_constraints[varname][dof] = val; }
475 
476 
477  void add_temporaries(const varnamelist &vl, gmm::uint64_type id_num) const;
478 
479  const mimlist &mimlist_of_brick(size_type ib) const
480  { return bricks[ib].mims; }
481 
482  const varnamelist &varnamelist_of_brick(size_type ib) const
483  { return bricks[ib].vlist; }
484 
485  const varnamelist &datanamelist_of_brick(size_type ib) const
486  { return bricks[ib].dlist; }
487 
488  size_type region_of_brick(size_type ib) const
489  { return bricks[ib].region; }
490 
491  bool temporary_uptodate(const std::string &varname,
492  gmm::uint64_type id_num, size_type &ind) const;
493 
494  size_type n_iter_of_variable(const std::string &name) const {
495  return variables.count(name) == 0 ? size_type(0)
496  : variables[name].n_iter;
497  }
498 
499  void set_default_iter_of_variable(const std::string &varname,
500  size_type ind) const;
501  void reset_default_iter_of_variables(const varnamelist &vl) const;
502 
503  void update_from_context() const { act_size_to_be_done = true; }
504 
505  const model_real_sparse_matrix &linear_real_matrix_term
506  (size_type ib, size_type iterm);
507 
508  const model_complex_sparse_matrix &linear_complex_matrix_term
509  (size_type ib, size_type iterm);
510 
511  /** Disable a brick. */
513  GMM_ASSERT1(valid_bricks[ib], "Inexistent brick");
514  active_bricks.del(ib);
515  }
516 
517  /** Enable a brick. */
519  GMM_ASSERT1(valid_bricks[ib], "Inexistent brick");
520  active_bricks.add(ib);
521  }
522 
523  /** Disable a variable (and its attached mutlipliers). */
524  void disable_variable(const std::string &name);
525 
526  /** Enable a variable (and its attached mutlipliers). */
527  void enable_variable(const std::string &name, bool enabled=true);
528 
529  /** States if a name corresponds to a declared variable. */
530  bool variable_exists(const std::string &name) const;
531 
532  /** States if a variable is disabled (treated as data). */
533  bool is_disabled_variable(const std::string &name) const;
534 
535  /** States if a name corresponds to a declared data or disabled variable. */
536  bool is_data(const std::string &name) const;
537 
538  /** States if a name corresponds to a declared data. */
539  bool is_true_data(const std::string &name) const;
540 
541  /** States if a variable is condensed out of the global system. */
542  bool is_internal_variable(const std::string &name) const;
543 
544  bool is_affine_dependent_variable(const std::string &name) const;
545 
546  const std::string &org_variable(const std::string &name) const;
547 
548  const scalar_type &factor_of_variable(const std::string &name) const;
549 
550  void set_factor_of_variable(const std::string &name, scalar_type a);
551 
552  bool is_im_data(const std::string &name) const;
553 
554  const im_data *pim_data_of_variable(const std::string &name) const;
555 
556  const gmm::uint64_type &
557  version_number_of_data_variable(const std::string &varname,
558  size_type niter = size_type(-1)) const;
559 
560  /** Boolean which says if the model deals with real or complex unknowns
561  and data. */
562  bool is_complex() const { return complex_version; }
563 
564  /** Return true if all the model terms do not affect the coercivity of
565  the whole tangent system. */
566  bool is_coercive() const { return is_coercive_; }
567 
568  /** Return true if the model has at least one internal variable. */
569  bool has_internal_variables() const {
570  for (const auto &v : variables)
571  if (v.second.is_internal && !v.second.is_disabled) return true;
572  return false;
573  }
574 
575  /** Return true if all the model terms do not affect the coercivity of
576  the whole tangent system. */
577  bool is_symmetric() const { return is_symmetric_; }
578 
579  /** Return true if all the model terms are linear. */
580  bool is_linear() const { return is_linear_; }
581 
582  /** Total number of degrees of freedom in the model. */
583  size_type nb_dof(bool with_internal=false) const;
584 
585  /** Number of internal degrees of freedom in the model. */
586  size_type nb_internal_dof() const { return nb_dof(true)-nb_dof(); }
587 
588  /** Number of primary degrees of freedom in the model. */
589  size_type nb_primary_dof() const { return nb_dof(); }
590 
591  /** Leading dimension of the meshes used in the model. */
592  dim_type leading_dimension() const { return leading_dim; }
593 
594  /** Gives a non already existing variable name begining by `name`. */
595  std::string new_name(const std::string &name);
596 
597  const gmm::sub_interval &
598  interval_of_variable(const std::string &name) const;
599 
600  /** Gives the access to the vector value of a variable. For the real
601  version. */
602  const model_real_plain_vector &
603  real_variable(const std::string &name, size_type niter) const;
604 
605  /**The same as above, but either accessing the latest variable version,
606  or the previous, if using "Old_" prefix*/
607  const model_real_plain_vector &
608  real_variable(const std::string &name) const;
609 
610  /** Gives the access to the vector value of a variable. For the complex
611  version. */
612  const model_complex_plain_vector &
613  complex_variable(const std::string &name, size_type niter) const;
614 
615  /**The same as above, but either accessing the latest variable version,
616  or the previous, if using "Old_" prefix*/
617  const model_complex_plain_vector &
618  complex_variable(const std::string &name) const;
619 
620  /** Gives the write access to the vector value of a variable. Make a
621  change flag of the variable set. For the real version. */
622  model_real_plain_vector &
623  set_real_variable(const std::string &name, size_type niter) const;
624 
625  /**The same as above, but for either latest variable, or
626  for the previous, if prefixed with "Old_".*/
627  model_real_plain_vector &
628  set_real_variable(const std::string &name) const;
629 
630  /** Gives the write access to the vector value of a variable. Make a
631  change flag of the variable set. For the complex version. */
632  model_complex_plain_vector &
633  set_complex_variable(const std::string &name, size_type niter) const;
634 
635  /**The same as above, but either accessing the latest variable version,
636  or the previous, if using "Old_" prefix*/
637  model_complex_plain_vector &
638  set_complex_variable(const std::string &name) const;
639 
640  model_real_plain_vector &
641  set_real_constant_part(const std::string &name) const;
642 
643  model_complex_plain_vector &
644  set_complex_constant_part(const std::string &name) const;
645 
646  private:
647  template<typename VECTOR, typename T>
648  void from_variables(VECTOR &V, bool with_internal, T) const {
649  for (const auto &v : variables)
650  if (v.second.is_variable && !v.second.is_affine_dependent
651  && !v.second.is_disabled
652  && (with_internal || !v.second.is_internal))
653  gmm::copy(v.second.real_value[0], gmm::sub_vector(V, v.second.I));
654  }
655 
656  template<typename VECTOR, typename T>
657  void from_variables(VECTOR &V, bool with_internal, std::complex<T>) const {
658  for (const auto &v : variables)
659  if (v.second.is_variable && !v.second.is_affine_dependent
660  && !v.second.is_disabled
661  && (with_internal || !v.second.is_internal))
662  gmm::copy(v.second.complex_value[0], gmm::sub_vector(V, v.second.I));
663  }
664 
665  public:
666  template<typename VECTOR>
667  void from_variables(VECTOR &V, bool with_internal=false) const {
668  typedef typename gmm::linalg_traits<VECTOR>::value_type T;
669  context_check(); if (act_size_to_be_done) actualize_sizes();
670  from_variables(V, with_internal, T());
671  }
672 
673  private:
674  template<typename VECTOR, typename T>
675  void to_variables(const VECTOR &V, bool with_internal, T) {
676  for (auto &&v : variables)
677  if (v.second.is_variable && !v.second.is_affine_dependent
678  && !v.second.is_disabled
679  && (with_internal || !v.second.is_internal)) {
680  gmm::copy(gmm::sub_vector(V, v.second.I), v.second.real_value[0]);
681  v.second.v_num_data[0] = act_counter();
682  }
683  update_affine_dependent_variables();
684  this->post_to_variables_step();
685  }
686 
687  template<typename VECTOR, typename T>
688  void to_variables(const VECTOR &V, bool with_internal, std::complex<T>) {
689  for (auto &&v : variables)
690  if (v.second.is_variable && !v.second.is_affine_dependent
691  && !v.second.is_disabled
692  && (with_internal || !v.second.is_internal)) {
693  gmm::copy(gmm::sub_vector(V, v.second.I), v.second.complex_value[0]);
694  v.second.v_num_data[0] = act_counter();
695  }
696  update_affine_dependent_variables();
697  this->post_to_variables_step();
698  }
699 
700  public:
701  template<typename VECTOR>
702  void to_variables(const VECTOR &V, bool with_internal=false) {
703  typedef typename gmm::linalg_traits<VECTOR>::value_type T;
704  context_check(); if (act_size_to_be_done) actualize_sizes();
705  to_variables(V, with_internal, T());
706  }
707 
708  /** Add a fixed size variable to the model assumed to be a vector.
709  niter is the number of version of the variable stored. */
710  void add_fixed_size_variable(const std::string &name, size_type size,
711  size_type niter = 1);
712 
713  /** Add a fixed size variable to the model whith given tensor dimensions.
714  niter is the number of version of the variable stored. */
715  void add_fixed_size_variable(const std::string &name,
716  const bgeot::multi_index &sizes,
717  size_type niter = 1);
718 
719  /** Add a fixed size data to the model. niter is the number of version
720  of the data stored, for time integration schemes. */
721  void add_fixed_size_data(const std::string &name, size_type size,
722  size_type niter = 1);
723 
724  /** Add a fixed size data to the model. niter is the number of version
725  of the data stored, for time integration schemes. */
726  void add_fixed_size_data(const std::string &name,
727  const bgeot::multi_index &sizes,
728  size_type niter = 1);
729 
730  /** Resize a fixed size variable (or data) of the model. */
731  void resize_fixed_size_variable(const std::string &name, size_type size);
732 
733  /** Resize a fixed size variable (or data) of the model. */
734  void resize_fixed_size_variable(const std::string &name,
735  const bgeot::multi_index &sizes);
736 
737  /** Add a fixed size data (assumed to be a vector) to the model and
738  initialized with v. */
739  template <typename VECT>
740  void add_initialized_fixed_size_data(const std::string &name,
741  const VECT &v) {
742  this->add_fixed_size_data(name, gmm::vect_size(v));
743  if (this->is_complex())
744  gmm::copy(v, this->set_complex_variable(name));
745  else
746  gmm::copy(gmm::real_part(v), this->set_real_variable(name));
747  }
748 
749  /** Add a fixed size data (assumed to be a vector) to the model and
750  initialized with v. */
751  template <typename VECT>
752  void add_initialized_fixed_size_data(const std::string &name,
753  const VECT &v,
754  const bgeot::multi_index &sizes) {
755  this->add_fixed_size_data(name, sizes);
756  if (this->is_complex())
757  gmm::copy(v, this->set_complex_variable(name));
758  else
759  gmm::copy(gmm::real_part(v), this->set_real_variable(name));
760  }
761 
762  /** Add a fixed size data (assumed to be a matrix) to the model and
763  initialized with M. */
764  void add_initialized_matrix_data(const std::string &name,
765  const base_matrix &M);
766  void add_initialized_matrix_data(const std::string &name,
767  const base_complex_matrix &M);
768 
769  /** Add a fixed size data (assumed to be a tensor) to the model and
770  initialized with t. */
771  void add_initialized_tensor_data(const std::string &name,
772  const base_tensor &t);
773  void add_initialized_tensor_data(const std::string &name,
774  const base_complex_tensor &t);
775 
776 
777  /** Add a scalar data (i.e. of size 1) to the model initialized with e. */
778  template <typename T>
779  void add_initialized_scalar_data(const std::string &name, T e) {
780  this->add_fixed_size_data(name, 1, 1);
781  if (this->is_complex())
782  this->set_complex_variable(name)[0] = e;
783  else
784  this->set_real_variable(name)[0] = gmm::real(e);
785  }
786 
787 
788  /** Add variable defined at integration points. */
789  void add_im_variable(const std::string &name, const im_data &imd,
790  size_type niter = 1);
791  /** Add internal variable, defined at integration points and condensed. */
792  void add_internal_im_variable(const std::string &name,
793  const im_data &imd);
794  /** Add data defined at integration points. */
795  void add_im_data(const std::string &name, const im_data &imd,
796  size_type niter = 1);
797 
798  /** Add a variable being the dofs of a finite element method to the model.
799  niter is the number of version of the variable stored, for time
800  integration schemes. */
801  void add_fem_variable(const std::string &name, const mesh_fem &mf,
802  size_type niter = 1);
803 
804  /** Add a variable linked to a fem with the dof filtered with respect
805  to a mesh region. Only the dof returned by the dof_on_region
806  method of `mf` will be kept. niter is the number of version
807  of the data stored, for time integration schemes. */
808  void add_filtered_fem_variable(const std::string &name, const mesh_fem &mf,
809  size_type region, size_type niter = 1);
810 
811 
812  /** Add a "virtual" variable be an affine depedent variable with respect
813  to another variable. Mainly used for time integration scheme for
814  instance to represent time derivative of variables.
815  `alpha` is the multiplicative scalar of the dependency. */
816  void add_affine_dependent_variable(const std::string &name,
817  const std::string &org_name,
818  scalar_type alpha = scalar_type(1));
819 
820  /** Add a data being the dofs of a finite element method to the model.*/
821  void add_fem_data(const std::string &name, const mesh_fem &mf,
822  dim_type qdim = 1, size_type niter = 1);
823 
824  /** Add a data being the dofs of a finite element method to the model.*/
825  void add_fem_data(const std::string &name, const mesh_fem &mf,
826  const bgeot::multi_index &sizes, size_type niter = 1);
827 
828  /** Add an initialized fixed size data to the model, assumed to be a
829  vector field if the size of the vector is a multiple of the dof
830  number. */
831  template <typename VECT>
832  void add_initialized_fem_data(const std::string &name, const mesh_fem &mf,
833  const VECT &v) {
834  this->add_fem_data(name, mf,
835  dim_type(gmm::vect_size(v) / mf.nb_dof()), 1);
836  if (this->is_complex())
837  gmm::copy(v, this->set_complex_variable(name));
838  else
839  gmm::copy(gmm::real_part(v), this->set_real_variable(name));
840  }
841 
842  /** Add a fixed size data to the model. The data is a tensor of given
843  sizes on each dof of the finite element method. */
844  template <typename VECT>
845  void add_initialized_fem_data(const std::string &name, const mesh_fem &mf,
846  const VECT &v,
847  const bgeot::multi_index &sizes) {
848  this->add_fem_data(name, mf, sizes, 1);
849  if (this->is_complex())
850  gmm::copy(v, this->set_complex_variable(name));
851  else
852  gmm::copy(gmm::real_part(v), this->set_real_variable(name));
853  }
854 
855  /** Add a particular variable linked to a fem being a multiplier with
856  respect to a primal variable. The dof will be filtered with the
857  gmm::range_basis function applied on the terms of the model which
858  link the multiplier and the primal variable. Optimized for boundary
859  multipliers. niter is the number of version of the data stored,
860  for time integration schemes. */
861  void add_multiplier(const std::string &name, const mesh_fem &mf,
862  const std::string &primal_name,
863  size_type niter = 1);
864 
865  /** Add a particular variable linked to a fem being a multiplier with
866  respect to a primal variable and a region. The dof will be filtered
867  both with the gmm::range_basis function applied on the terms of
868  the model which link the multiplier and the primal variable and on
869  the dof on the given region. Optimized for boundary
870  multipliers. niter is the number of version of the data stored,
871  for time integration schemes. */
872  void add_multiplier(const std::string &name, const mesh_fem &mf,
873  size_type region, const std::string &primal_name,
874  size_type niter = 1);
875 
876  /** Add a particular variable linked to a fem being a multiplier with
877  respect to a primal variable. The dof will be filtered with the
878  gmm::range_basis function applied on the mass matrix between the fem
879  of the multiplier and the one of the primal variable.
880  Optimized for boundary multipliers. niter is the number of version
881  of the data stored, for time integration schemes. */
882  void add_multiplier(const std::string &name, const mesh_fem &mf,
883  const std::string &primal_name, const mesh_im &mim,
884  size_type region, size_type niter = 1);
885 
886  /** Dictonnary of user defined macros. */
887  const ga_macro_dictionary &macro_dictionary() const { return macro_dict; }
888 
889  /** Add a macro definition for the high generic assembly language.
890  This macro can be used for the definition of generic assembly bricks.
891  The name of a macro cannot coincide with a variable name. */
892  void add_macro(const std::string &name, const std::string &expr);
893 
894  /** Delete a previously defined macro definition. */
895  void del_macro(const std::string &name);
896 
897  /** Says if a macro of that name has been defined. */
898  bool macro_exists(const std::string &name) const
899  { return macro_dict.macro_exists(name); }
900 
901  /** Delete a variable or data of the model. */
902  void delete_variable(const std::string &varname);
903 
904  /** Gives the access to the mesh_fem of a variable if any. Throw an
905  exception otherwise. */
906  const mesh_fem &mesh_fem_of_variable(const std::string &name) const;
907 
908  /** Gives a pointer to the mesh_fem of a variable if any. 0 otherwise.*/
909  const mesh_fem *pmesh_fem_of_variable(const std::string &name) const;
910 
911 
912  bgeot::multi_index qdims_of_variable(const std::string &name) const;
913  size_type qdim_of_variable(const std::string &name) const;
914 
915  /** Gives the access to the tangent matrix. For the real version. */
916  const model_real_sparse_matrix &
917  real_tangent_matrix(bool internal=false) const {
918  GMM_ASSERT1(!complex_version, "This model is a complex one");
919  context_check(); if (act_size_to_be_done) actualize_sizes();
920  return internal ? internal_rTM : rTM;
921  }
922 
923  /** Gives the access to the tangent matrix. For the complex version. */
924  const model_complex_sparse_matrix &complex_tangent_matrix() const {
925  GMM_ASSERT1(complex_version, "This model is a real one");
926  context_check(); if (act_size_to_be_done) actualize_sizes();
927  return cTM;
928  }
929 
930  /** Gives access to the right hand side of the tangent linear system.
931  For the real version. An assembly of the rhs has to be done first. */
932  const model_real_plain_vector &real_rhs(bool with_internal=false) const {
933  GMM_ASSERT1(!complex_version, "This model is a complex one");
934  context_check(); if (act_size_to_be_done) actualize_sizes();
935  return (with_internal && gmm::vect_size(full_rrhs)) ? full_rrhs : rrhs;
936  }
937 
938  /** Gives write access to the right hand side of the tangent linear system.
939  Some solvers need to manipulate the model rhs directly so that for
940  example internal condensed variables can be treated properly. */
941  model_real_plain_vector &set_real_rhs(bool with_internal=false) const {
942  GMM_ASSERT1(!complex_version, "This model is a complex one");
943  context_check(); if (act_size_to_be_done) actualize_sizes();
944  return (with_internal && gmm::vect_size(full_rrhs)) ? full_rrhs : rrhs;
945  }
946 
947  /** Gives access to the partial solution for condensed internal variables.
948  A matrix assembly with condensation has to be done first. */
949  const model_real_plain_vector &internal_solution() const {
950  GMM_ASSERT1(!complex_version, "This model is a complex one");
951  context_check(); if (act_size_to_be_done) actualize_sizes();
952  return internal_sol;
953  }
954 
955  /** Gives access to the part of the right hand side of a term of
956  a particular nonlinear brick. Does not account of the eventual time
957  dispatcher. An assembly of the rhs has to be done first.
958  For the real version. */
959  const model_real_plain_vector &real_brick_term_rhs
960  (size_type ib, size_type ind_term = 0, bool sym = false,
961  size_type ind_iter = 0) const
962  {
963  GMM_ASSERT1(!complex_version, "This model is a complex one");
964  context_check(); if (act_size_to_be_done) actualize_sizes();
965  GMM_ASSERT1(valid_bricks[ib], "Inexistent brick");
966  GMM_ASSERT1(ind_term < bricks[ib].tlist.size(), "Inexistent term");
967  GMM_ASSERT1(ind_iter < bricks[ib].nbrhs, "Inexistent iter");
968  GMM_ASSERT1(!sym || bricks[ib].tlist[ind_term].is_symmetric,
969  "Term is not symmetric");
970  if (sym)
971  return bricks[ib].rveclist_sym[ind_iter][ind_term];
972  else
973  return bricks[ib].rveclist[ind_iter][ind_term];
974  }
975 
976  /** Gives access to the right hand side of the tangent linear system.
977  For the complex version. */
978  const model_complex_plain_vector &complex_rhs() const {
979  GMM_ASSERT1(complex_version, "This model is a real one");
980  context_check(); if (act_size_to_be_done) actualize_sizes();
981  return crhs;
982  }
983 
984  /** Gives write access to the right hand side of the tangent linear system.
985  Some solvers need to manipulate the model rhs directly so that for
986  example internal condensed variables can be treated properly. */
987  model_complex_plain_vector &set_complex_rhs() const {
988  GMM_ASSERT1(complex_version, "This model is a real one");
989  context_check(); if (act_size_to_be_done) actualize_sizes();
990  return crhs;
991  }
992 
993  /** Gives access to the part of the right hand side of a term of a
994  particular nonlinear brick. Does not account of the eventual time
995  dispatcher. An assembly of the rhs has to be done first.
996  For the complex version. */
997  const model_complex_plain_vector &complex_brick_term_rhs
998  (size_type ib, size_type ind_term = 0, bool sym = false,
999  size_type ind_iter = 0) const
1000  {
1001  GMM_ASSERT1(!complex_version, "This model is a complex one");
1002  context_check(); if (act_size_to_be_done) actualize_sizes();
1003  GMM_ASSERT1(valid_bricks[ib], "Inexistent brick");
1004  GMM_ASSERT1(ind_term < bricks[ib].tlist.size(), "Inexistent term");
1005  GMM_ASSERT1(ind_iter < bricks[ib].nbrhs, "Inexistent iter");
1006  GMM_ASSERT1(!sym || bricks[ib].tlist[ind_term].is_symmetric,
1007  "Term is not symmetric");
1008  if (sym)
1009  return bricks[ib].cveclist_sym[ind_iter][ind_term];
1010  else
1011  return bricks[ib].cveclist[ind_iter][ind_term];
1012  }
1013 
1014  /** List the model variables and constant. */
1015  void listvar(std::ostream &ost) const;
1016 
1017  void listresiduals(std::ostream &ost) const;
1018 
1019  /** List the model bricks. */
1020  void listbricks(std::ostream &ost, size_type base_id = 0) const;
1021 
1022  /** Return the model brick ids. */
1023  const dal::bit_vector& get_active_bricks() const {
1024  return active_bricks;
1025  }
1026 
1027  /** Force the re-computation of a brick for the next assembly. */
1029  GMM_ASSERT1(valid_bricks[ib], "Inexistent brick");
1030  bricks[ib].terms_to_be_computed = true;
1031  }
1032 
1033  /** Add a brick to the model. varname is the list of variable used
1034  and datanames the data used. If a variable is used as a data, it
1035  should be declared in the datanames (it will depend on the value of
1036  the variable not only on the fem). Returns the brick index. */
1037  size_type add_brick(pbrick pbr, const varnamelist &varnames,
1038  const varnamelist &datanames,
1039  const termlist &terms, const mimlist &mims,
1040  size_type region);
1041 
1042  /** Delete the brick of index ib from the model. */
1043  void delete_brick(size_type ib);
1044 
1045  /** Add an integration method to a brick. */
1046  void add_mim_to_brick(size_type ib, const mesh_im &mim);
1047 
1048  /** Change the term list of a brick. Used for very special bricks only. */
1049  void change_terms_of_brick(size_type ib, const termlist &terms);
1050 
1051  /** Change the variable list of a brick. Used for very special bricks only.
1052  */
1053  void change_variables_of_brick(size_type ib, const varnamelist &vl);
1054 
1055  /** Change the data list of a brick. Used for very special bricks only.
1056  */
1057  void change_data_of_brick(size_type ib, const varnamelist &vl);
1058 
1059  /** Change the mim list of a brick. Used for very special bricks only.
1060  */
1061  void change_mims_of_brick(size_type ib, const mimlist &ml);
1062 
1063  /** Change the update flag of a brick. Used for very special bricks only.
1064  */
1065  void change_update_flag_of_brick(size_type ib, bool flag);
1066 
1067  void set_time(scalar_type t = scalar_type(0), bool to_init = true);
1068 
1069  scalar_type get_time();
1070 
1071  void set_time_step(scalar_type dt) { time_step = dt; }
1072  scalar_type get_time_step() const { return time_step; }
1073  scalar_type get_init_time_step() const { return init_time_step; }
1074  int is_time_integration() const { return time_integration; }
1075  void set_time_integration(int ti) { time_integration = ti; }
1076  bool is_init_step() const { return init_step; }
1077  void cancel_init_step() { init_step = false; }
1078  void call_init_affine_dependent_variables(int version);
1079  void shift_variables_for_time_integration();
1080  void copy_init_time_derivative();
1081  void add_time_integration_scheme(const std::string &varname,
1082  ptime_scheme ptsc);
1083  void perform_init_time_derivative(scalar_type ddt)
1084  { init_step = true; init_time_step = ddt; }
1085 
1086 
1087  /** Add a time dispacther to a brick. */
1088  void add_time_dispatcher(size_type ibrick, pdispatcher pdispatch);
1089 
1090  void set_dispatch_coeff();
1091 
1092  /** For transient problems. Initialisation of iterations. */
1093  virtual void first_iter();
1094 
1095  /** For transient problems. Prepare the next iterations. In particular
1096  shift the version of the variables.
1097  */
1098  virtual void next_iter();
1099 
1100  /** Add an interpolate transformation to the model to be used with the
1101  generic assembly.
1102  */
1103  void add_interpolate_transformation(const std::string &name,
1104  pinterpolate_transformation ptrans) {
1105  if (secondary_domain_exists(name))
1106  GMM_ASSERT1(false, "An secondary domain with the same "
1107  "name already exists");
1108  if (transformations.count(name) > 0)
1109  GMM_ASSERT1(name.compare("neighbor_element"), "neighbor_element is a "
1110  "reserved interpolate transformation name");
1111  transformations[name] = ptrans;
1112  }
1113 
1114  /** Get a pointer to the interpolate transformation `name`.
1115  */
1116  pinterpolate_transformation
1117  interpolate_transformation(const std::string &name) const {
1118  std::map<std::string, pinterpolate_transformation>::const_iterator
1119  it = transformations.find(name);
1120  GMM_ASSERT1(it != transformations.end(), "Inexistent transformation " << name);
1121  return it->second;
1122  }
1123 
1124  /** Tests if `name` corresponds to an interpolate transformation.
1125  */
1126  bool interpolate_transformation_exists(const std::string &name) const
1127  { return transformations.count(name) > 0; }
1128 
1129  /** Add an elementary transformation to the model to be used with the
1130  generic assembly.
1131  */
1132  void add_elementary_transformation(const std::string &name,
1133  pelementary_transformation ptrans) {
1134  elem_transformations[name] = ptrans;
1135  }
1136 
1137  /** Get a pointer to the elementary transformation `name`.
1138  */
1139  pelementary_transformation
1140  elementary_transformation(const std::string &name) const {
1141  std::map<std::string, pelementary_transformation>::const_iterator
1142  it = elem_transformations.find(name);
1143  GMM_ASSERT1(it != elem_transformations.end(),
1144  "Inexistent elementary transformation " << name);
1145  return it->second;
1146  }
1147 
1148  /** Tests if `name` corresponds to an elementary transformation.
1149  */
1150  bool elementary_transformation_exists(const std::string &name) const
1151  { return elem_transformations.count(name) > 0; }
1152 
1153 
1154  /** Add a secondary domain to the model to be used with the
1155  generic assembly.
1156  */
1157  void add_secondary_domain(const std::string &name,
1158  psecondary_domain ptrans) {
1159  if (interpolate_transformation_exists(name))
1160  GMM_ASSERT1(false, "An interpolate transformation with the same "
1161  "name already exists");secondary_domains[name] = ptrans;
1162  }
1163 
1164  /** Get a pointer to the interpolate transformation `name`.
1165  */
1166  psecondary_domain
1167  secondary_domain(const std::string &name) const {
1168  auto it = secondary_domains.find(name);
1169  GMM_ASSERT1(it != secondary_domains.end(),
1170  "Inexistent transformation " << name);
1171  return it->second;
1172  }
1173 
1174  /** Tests if `name` corresponds to an interpolate transformation.
1175  */
1176  bool secondary_domain_exists(const std::string &name) const
1177  { return secondary_domains.count(name) > 0; }
1178 
1179  /** Gives the name of the variable of index `ind_var` of the brick
1180  of index `ind_brick`. */
1181  const std::string &varname_of_brick(size_type ind_brick,
1182  size_type ind_var);
1183 
1184  /** Gives the name of the data of index `ind_data` of the brick
1185  of index `ind_brick`. */
1186  const std::string &dataname_of_brick(size_type ind_brick,
1187  size_type ind_data);
1188 
1189  /** Assembly of the tangent system taking into account the terms
1190  from all bricks. */
1191  virtual void assembly(build_version version);
1192 
1193  /** Gives the assembly string corresponding to the Neumann term of
1194  the fem variable `varname` on `region`. It is deduced from the
1195  assembly string declared by the model bricks.
1196  `region` should be the index of a boundary region
1197  on the mesh where `varname` is defined. Care to call this function
1198  only after all the volumic bricks have been declared.
1199  Complains, if a brick omit to declare an assembly string.
1200  */
1201  std::string Neumann_term(const std::string &varname, size_type region);
1202 
1203  virtual void clear();
1204 
1205  explicit model(bool comp_version = false);
1206 
1207  /** check consistency of RHS and Stiffness matrix for brick with
1208  @param ind_brick - index of the brick
1209  */
1210  void check_brick_stiffness_rhs(size_type ind_brick) const;
1211 
1212 
1213  };
1214 
1215  //=========================================================================
1216  //
1217  // Time integration scheme object.
1218  //
1219  //=========================================================================
1220 
1221  /** The time integration scheme object provides the necessary methods
1222  for the model object to apply a time integration scheme to an
1223  evolutionnary problem.
1224  **/
1225  class APIDECL virtual_time_scheme {
1226 
1227  public:
1228 
1229  virtual void init_affine_dependent_variables(model &md) const = 0;
1230  virtual void init_affine_dependent_variables_precomputation(model &md)
1231  const = 0;
1232  virtual void time_derivative_to_be_initialized
1233  (std::string &name_v, std::string &name_previous_v) const = 0;
1234  virtual void shift_variables(model &md) const = 0;
1235  virtual ~virtual_time_scheme() {}
1236  };
1237 
1238  void add_theta_method_for_first_order(model &md, const std::string &varname,
1239  scalar_type theta);
1240 
1241  void add_theta_method_for_second_order(model &md, const std::string &varname,
1242  scalar_type theta);
1243 
1244  void add_Newmark_scheme(model &md, const std::string &varname,
1245  scalar_type beta, scalar_type gamma);
1246 
1247  void add_Houbolt_scheme(model &md, const std::string &varname);
1248 
1249 
1250 
1251  //=========================================================================
1252  //
1253  // Time dispatcher object.
1254  //
1255  //=========================================================================
1256 
1257  /** The time dispatcher object modify the result of a brick in order to
1258  apply a time integration scheme.
1259  **/
1260  class APIDECL virtual_dispatcher {
1261 
1262  protected:
1263 
1264  size_type nbrhs_;
1265  std::vector<std::string> param_names;
1266 
1267  public:
1268 
1269  size_type nbrhs() const { return nbrhs_; }
1270 
1271  typedef model::build_version build_version;
1272 
1273  virtual void set_dispatch_coeff(const model &, size_type) const
1274  { GMM_ASSERT1(false, "Time dispatcher with not set_dispatch_coeff !"); }
1275 
1276  virtual void next_real_iter
1277  (const model &, size_type, const model::varnamelist &,
1278  const model::varnamelist &,
1279  model::real_matlist &,
1280  std::vector<model::real_veclist> &,
1281  std::vector<model::real_veclist> &,
1282  bool) const {
1283  GMM_ASSERT1(false, "Time dispatcher with not defined first real iter !");
1284  }
1285 
1286  virtual void next_complex_iter
1287  (const model &, size_type, const model::varnamelist &,
1288  const model::varnamelist &,
1289  model::complex_matlist &,
1290  std::vector<model::complex_veclist> &,
1291  std::vector<model::complex_veclist> &,
1292  bool) const{
1293  GMM_ASSERT1(false,"Time dispatcher with not defined first comples iter");
1294  }
1295 
1296  virtual void asm_real_tangent_terms
1297  (const model &, size_type,
1298  model::real_matlist &, std::vector<model::real_veclist> &,
1299  std::vector<model::real_veclist> &,
1300  build_version) const {
1301  GMM_ASSERT1(false, "Time dispatcher with not defined real tangent "
1302  "terms !");
1303  }
1304 
1305  virtual void asm_complex_tangent_terms
1306  (const model &, size_type,
1307  model::complex_matlist &, std::vector<model::complex_veclist> &,
1308  std::vector<model::complex_veclist> &,
1309  build_version) const {
1310  GMM_ASSERT1(false, "Time dispatcher with not defined complex tangent "
1311  "terms !");
1312  }
1313 
1314  virtual_dispatcher(size_type _nbrhs) : nbrhs_(_nbrhs) {
1315  GMM_ASSERT1(_nbrhs > 0, "Time dispatcher with no rhs");
1316  }
1317  virtual ~virtual_dispatcher() {}
1318 
1319  };
1320 
1321  // ----------------------------------------------------------------------
1322  //
1323  // theta-method dispatcher
1324  //
1325  // ----------------------------------------------------------------------
1326 
1327  class APIDECL theta_method_dispatcher : public virtual_dispatcher {
1328 
1329  public:
1330 
1331  typedef model::build_version build_version;
1332 
1333  void set_dispatch_coeff(const model &md, size_type ib) const;
1334 
1335  template <typename MATLIST, typename VECTLIST>
1336  inline void next_iter(const model &md, size_type ib,
1337  const model::varnamelist &/* vl */,
1338  const model::varnamelist &/* dl */,
1339  MATLIST &/* matl */,
1340  VECTLIST &vectl, VECTLIST &vectl_sym,
1341  bool first_iter) const {
1342  if (first_iter) md.update_brick(ib, model::BUILD_RHS);
1343 
1344  // shift the rhs
1345  for (size_type i = 0; i < vectl[0].size(); ++i)
1346  gmm::copy(vectl[0][i], vectl[1][i]);
1347  for (size_type i = 0; i < vectl_sym[0].size(); ++i)
1348  gmm::copy(vectl_sym[0][i], vectl_sym[1][i]);
1349 
1350  // add the component represented by the linear matrix terms to the
1351  // supplementary rhs
1352  md.linear_brick_add_to_rhs(ib, 1, 0);
1353  }
1354 
1355  void next_real_iter
1356  (const model &md, size_type ib, const model::varnamelist &vl,
1357  const model::varnamelist &dl, model::real_matlist &matl,
1358  std::vector<model::real_veclist> &vectl,
1359  std::vector<model::real_veclist> &vectl_sym, bool first_iter) const;
1360 
1361  void next_complex_iter
1362  (const model &md, size_type ib, const model::varnamelist &vl,
1363  const model::varnamelist &dl,
1364  model::complex_matlist &matl,
1365  std::vector<model::complex_veclist> &vectl,
1366  std::vector<model::complex_veclist> &vectl_sym,
1367  bool first_iter) const;
1368 
1369  void asm_real_tangent_terms
1370  (const model &md, size_type ib, model::real_matlist &/* matl */,
1371  std::vector<model::real_veclist> &/* vectl */,
1372  std::vector<model::real_veclist> &/* vectl_sym */,
1373  build_version version) const;
1374 
1375  virtual void asm_complex_tangent_terms
1376  (const model &md, size_type ib, model::complex_matlist &/* matl */,
1377  std::vector<model::complex_veclist> &/* vectl */,
1378  std::vector<model::complex_veclist> &/* vectl_sym */,
1379  build_version version) const;
1380 
1381  theta_method_dispatcher(const std::string &THETA);
1382  };
1383 
1384  //=========================================================================
1385  //
1386  // Functions adding standard time dispatchers.
1387  //
1388  //=========================================================================
1389 
1390  /** Add a theta-method time dispatcher to a list of bricks. For instance,
1391  a matrix term $K$ will be replaced by
1392  $\theta K U^{n+1} + (1-\theta) K U^{n}$.
1393  */
1394  void APIDECL add_theta_method_dispatcher(model &md, dal::bit_vector ibricks,
1395  const std::string &THETA);
1396 
1397  /** Function which udpate the velocity $v^{n+1}$ after the computation
1398  of the displacement $u^{n+1}$ and before the next iteration. Specific
1399  for theta-method and when the velocity is included in the data
1400  of the model.
1401  */
1403  (model &md, const std::string &U, const std::string &V,
1404  const std::string &pdt, const std::string &ptheta);
1405 
1406 
1407  /** Add a midpoint time dispatcher to a list of bricks. For instance,
1408  a nonlinear term $K(U)$ will be replaced by
1409  $K((U^{n+1} + U^{n})/2)$.
1410  */
1411  void APIDECL add_midpoint_dispatcher(model &md, dal::bit_vector ibricks);
1412 
1413  /** Function which udpate the velocity $v^{n+1}$ after the computation
1414  of the displacement $u^{n+1}$ and before the next iteration. Specific
1415  for Newmark scheme and when the velocity is included in the data
1416  of the model. This version inverts the mass matrix by a conjugate
1417  gradient.
1418  */
1420  (model &md, size_type id2dt2b, const std::string &U, const std::string &V,
1421  const std::string &pdt, const std::string &ptwobeta,
1422  const std::string &pgamma);
1423 
1424  //=========================================================================
1425  //
1426  // Brick object.
1427  //
1428  //=========================================================================
1429 
1430  /** The virtual brick has to be derived to describe real model bricks.
1431  The set_flags method has to be called by the derived class.
1432  The virtual methods asm_real_tangent_terms and/or
1433  asm_complex_tangent_terms have to be defined.
1434  The brick should not store data. The data have to be stored in the
1435  model object.
1436  **/
1437  class APIDECL virtual_brick {
1438  protected:
1439  bool islinear; // The brick add a linear term or not.
1440  bool issymmetric; // The brick add a symmetric term or not.
1441  bool iscoercive; // The brick add a potentially coercive terms or not.
1442  // (in particular, not a term involving a multiplier)
1443  bool isreal; // The brick admits a real version or not.
1444  bool iscomplex; // The brick admits a complex version or not.
1445  bool isinit; // internal flag.
1446  bool compute_each_time; // The brick is linear but needs to be computed
1447  // each time it is evaluated.
1448  bool isUpdateBrick; // The brick does not contribute any terms to the
1449  // system matrix or right-hand side, but only updates state variables.
1450  std::string name; // Name of the brick.
1451 
1452  public:
1453 
1454  typedef model::build_version build_version;
1455 
1456  virtual_brick() { isinit = false; }
1457  virtual ~virtual_brick() { }
1458  void set_flags(const std::string &bname, bool islin, bool issym,
1459  bool iscoer, bool ire, bool isco, bool each_time = false) {
1460  name = bname;
1461  islinear = islin; issymmetric = issym; iscoercive = iscoer;
1462  isreal = ire; iscomplex = isco; isinit = true;
1463  compute_each_time = each_time;
1464  }
1465 
1466 # define BRICK_NOT_INIT GMM_ASSERT1(isinit, "Set brick flags !")
1467  bool is_linear() const { BRICK_NOT_INIT; return islinear; }
1468  bool is_symmetric() const { BRICK_NOT_INIT; return issymmetric; }
1469  bool is_coercive() const { BRICK_NOT_INIT; return iscoercive; }
1470  bool is_real() const { BRICK_NOT_INIT; return isreal; }
1471  bool is_complex() const { BRICK_NOT_INIT; return iscomplex; }
1472  bool is_to_be_computed_each_time() const
1473  { BRICK_NOT_INIT; return compute_each_time; }
1474  const std::string &brick_name() const { BRICK_NOT_INIT; return name; }
1475 
1476 
1477  /** Assembly of bricks real tangent terms.
1478  In case of Getfem's compilation with OpenMP option,
1479  this method is executed on multiple threads.
1480  The parallelism is provided by distributing all tangent matrices and
1481  vectors and accumulating them later into the original. Additionally,
1482  by default, all mesh_region objects, participating in the assembly,
1483  are also partitioned. In order to avoid data race conditions, this
1484  method should not modify any data simultaneously accessible from
1485  multiple threads. In case this is unavoidable, the race can be
1486  prevented by distributing this data (of type T) between the threads
1487  via getfem::omp_distribute<T> (prefered method) or
1488  protected from concurrent access with mutexes (e.g. getfem::omp_lock)
1489  or OpenMP critical section. */
1490  virtual void asm_real_tangent_terms(const model &, size_type,
1491  const model::varnamelist &,
1492  const model::varnamelist &,
1493  const model::mimlist &,
1494  model::real_matlist &,
1495  model::real_veclist &,
1496  model::real_veclist &,
1497  size_type, build_version) const
1498  { /** doesn't have to be overriden if serial pre- post- assemblies are
1499  defined */
1500  }
1501 
1502 
1503  /** Assembly of bricks complex tangent terms.
1504  In case of Getfem's compilation with OpenMP option,
1505  this method is executed on multiple threads.
1506  The parallelism is provided by distributing all tangent matrices and
1507  vectors and accumulating them later into the original. Additionally,
1508  by default, all mesh_region objects, participating in the assembly,
1509  are also partitioned. In order to avoid data race conditions, this
1510  method should not modify any data simultaneously accessible from
1511  multiple threads. In case this is unavoidable, the race can be
1512  prevented by distributing this data (of type T) between the threads
1513  via getfem::omp_distribute<T> (prefered method) or
1514  protected from concurrent access with mutexes (e.g. getfem::omp_lock)
1515  or OpenMP critical section. */
1517  const model::varnamelist &,
1518  const model::varnamelist &,
1519  const model::mimlist &,
1520  model::complex_matlist &,
1521  model::complex_veclist &,
1522  model::complex_veclist &,
1523  size_type, build_version) const
1524  { /** doesn't have to be overriden if serial pre- post- assemblies are
1525  defined*/
1526  }
1527 
1528 
1529  /** Peform any pre assembly action for real term assembly. The purpose of
1530  this method is to do any action that cannot be peformed in the main
1531  assembly routines in parallel.
1532  Possible action can be modification of the model object, cashing
1533  some data that cannot be distributed etc. */
1535  const model::varnamelist &,
1536  const model::varnamelist &,
1537  const model::mimlist &,
1538  model::real_matlist &,
1539  model::real_veclist &,
1540  model::real_veclist &,
1541  size_type, build_version) const { };
1542 
1543  /** Peform any pre assembly action for complex term assembly. The purpose
1544  of this method is to do any action that cannot be peformed in the
1545  main assembly routines in parallel.
1546  Possible action can be modification of the model object, cashing
1547  some data that cannot be distributed etc. */
1549  const model::varnamelist &,
1550  const model::varnamelist &,
1551  const model::mimlist &,
1552  model::complex_matlist &,
1553  model::complex_veclist &,
1554  model::complex_veclist &,
1555  size_type, build_version) const { };
1556 
1557  /** Peform any post assembly action for real terms. The purpose of this
1558  method is to do any action that cannot be peformed in the main
1559  assembly routines in parallel.
1560  Possible action can be modification of the model object, cashing
1561  some data that cannot be distributed etc. */
1563  const model::varnamelist &,
1564  const model::varnamelist &,
1565  const model::mimlist &,
1566  model::real_matlist &,
1567  model::real_veclist &,
1568  model::real_veclist &,
1569  size_type, build_version) const { };
1570 
1571  /** Peform any post assembly action for complex terms. The purpose of this
1572  method is to do any action that cannot be peformed in the main
1573  assembly routines in parallel.
1574  Possible action can be modification of the model object, cashing
1575  some data that cannot be distributed etc. */
1577  const model::varnamelist &,
1578  const model::varnamelist &,
1579  const model::mimlist &,
1580  model::complex_matlist &,
1581  model::complex_veclist &,
1582  model::complex_veclist &,
1583  size_type, build_version) const { };
1584 
1585 
1586  /** check consistency of stiffness matrix and rhs */
1587  void check_stiffness_matrix_and_rhs(const model &, size_type,
1588  const model::termlist& tlist,
1589  const model::varnamelist &,
1590  const model::varnamelist &,
1591  const model::mimlist &,
1592  model::real_matlist &,
1593  model::real_veclist &,
1594  model::real_veclist &, size_type rg,
1595  const scalar_type delta = 1e-8) const;
1596  /** The brick may declare an assembly string for the computation of the
1597  Neumann terms (in order to prescribe boundary conditions with
1598  Nitche's method). */
1599  virtual std::string declare_volume_assembly_string
1600  (const model &, size_type, const model::varnamelist &,
1601  const model::varnamelist &) const {
1602  GMM_ASSERT1(false, "No assemby string declared, computation of Neumann "
1603  "term impossible for brick " << name);
1604  }
1605 
1606  private:
1607  /** simultaneous call to real_pre_assembly, real_assembly
1608  and real_post_assembly */
1609  void full_asm_real_tangent_terms_(const model &, size_type,
1610  const model::varnamelist &,
1611  const model::varnamelist &,
1612  const model::mimlist &,
1613  model::real_matlist &,
1614  model::real_veclist &,
1615  model::real_veclist &,
1616  size_type, build_version) const;
1617  };
1618 
1619  //=========================================================================
1620  //
1621  // Functions adding standard bricks to the model.
1622  //
1623  //=========================================================================
1624 
1625  /** Add a term given by the weak form language expression `expr` which will
1626  be assembled in region `region` and with the integration method `mim`.
1627  Only the matrix term will be taken into account, assuming that it is
1628  linear.
1629  The advantage of declaring a term linear instead of nonlinear is that
1630  it will be assembled only once and no assembly is necessary for the
1631  residual.
1632  Take care that if the expression contains some variables and if the
1633  expression is a potential or of first order (i.e. describe the weak
1634  form, not the derivative of the weak form), the expression will be
1635  derivated with respect to all variables.
1636  You can specify if the term is symmetric, coercive or not.
1637  If you are not sure, the better is to declare the term not symmetric
1638  and not coercive. But some solvers (conjugate gradient for instance)
1639  are not allowed for non-coercive problems.
1640  `brickname` is an optional name for the brick.
1641  */
1642  size_type APIDECL add_linear_term
1643  (model &md, const mesh_im &mim, const std::string &expr,
1644  size_type region = size_type(-1), bool is_sym = false,
1645  bool is_coercive = false, const std::string &brickname = "",
1646  bool return_if_nonlin = false);
1647 
1648  inline size_type APIDECL add_linear_generic_assembly_brick
1649  (model &md, const mesh_im &mim, const std::string &expr,
1650  size_type region = size_type(-1), bool is_sym = false,
1651  bool is_coercive = false, const std::string &brickname = "",
1652  bool return_if_nonlin = false) {
1653  return add_linear_term(md, mim, expr, region, is_sym,
1654  is_coercive, brickname, return_if_nonlin);
1655  }
1656 
1657  /** Add a nonlinear term given by the weak form language expression `expr`
1658  which will be assembled in region `region` and with the integration
1659  method `mim`.
1660  The expression can describe a potential or a weak form. Second order
1661  terms (i.e. containing second order test functions, Test2) are not
1662  allowed.
1663  You can specify if the term is symmetric, coercive or not.
1664  If you are not sure, the better is to declare the term not symmetric
1665  and not coercive. But some solvers (conjugate gradient for instance)
1666  are not allowed for non-coercive problems.
1667  `brickname` is an optional name for the brick.
1668  */
1670  (model &md, const mesh_im &mim, const std::string &expr,
1671  size_type region = size_type(-1), bool is_sym = false,
1672  bool is_coercive = false, const std::string &brickname = "");
1673 
1674  inline size_type APIDECL add_nonlinear_generic_assembly_brick
1675  (model &md, const mesh_im &mim, const std::string &expr,
1676  size_type region = size_type(-1), bool is_sym = false,
1677  bool is_coercive = false, const std::string &brickname = "") {
1678  return add_nonlinear_term(md, mim, expr, region,
1679  is_sym, is_coercive, brickname);
1680  }
1681 
1682 
1683  /** Add a source term given by the assembly string `expr` which will
1684  be assembled in region `region` and with the integration method `mim`.
1685  Only the residual term will be taken into account.
1686  Take care that if the expression contains some variables and if the
1687  expression is a potential, the expression will be
1688  derivated with respect to all variables.
1689  `brickname` is an optional name for the brick.
1690  */
1691  size_type APIDECL add_source_term
1692  (model &md, const mesh_im &mim, const std::string &expr,
1693  size_type region = size_type(-1),
1694  const std::string &brickname = std::string(),
1695  const std::string &directvarname = std::string(),
1696  const std::string &directdataname = std::string(),
1697  bool return_if_nonlin = false);
1698  inline size_type APIDECL add_source_term_generic_assembly_brick
1699  (model &md, const mesh_im &mim, const std::string &expr,
1700  size_type region = size_type(-1),
1701  const std::string &brickname = std::string(),
1702  const std::string &directvarname = std::string(),
1703  const std::string &directdataname = std::string(),
1704  bool return_if_nonlin = false) {
1705  return add_source_term(md, mim, expr, region, brickname,
1706  directvarname, directdataname, return_if_nonlin);
1707  }
1708 
1709  /** Adds a linear term given by a weak form language expression like
1710  ``add_linear_term`` function but for an integration on a direct
1711  product of two domains, a first specfied by ``mim`` and ``region``
1712  and a second one by ``secondary_domain`` which has to be declared
1713  first into the model.
1714  */
1716  (model &md, const mesh_im &mim, const std::string &expr,
1717  size_type region, const std::string &secondary_domain,
1718  bool is_sym = false, bool is_coercive = false,
1719  const std::string &brickname = "", bool return_if_nonlin = false);
1720 
1721  /** Adds a nonlinear term given by a weak form language expression like
1722  ``add_nonlinear_term`` function but for an integration on a direct
1723  product of two domains, a first specfied by ``mim`` and ``region``
1724  and a second one by ``secondary_domain`` which has to be declared
1725  first into the model.
1726  */
1728  (model &md, const mesh_im &mim, const std::string &expr,
1729  size_type region, const std::string &secondary_domain,
1730  bool is_sym = false, bool is_coercive = false,
1731  const std::string &brickname = "");
1732 
1733  /** Adds a source term given by a weak form language expression like
1734  ``add_source_term`` function but for an integration on a direct
1735  product of two domains, a first specfied by ``mim`` and ``region``
1736  and a second one by ``secondary_domain`` which has to be declared
1737  first into the model.
1738  */
1740  (model &md, const mesh_im &mim, const std::string &expr,
1741  size_type region, const std::string &secondary_domain,
1742  const std::string &brickname = std::string(),
1743  const std::string &directvarname = std::string(),
1744  const std::string &directdataname = std::string(),
1745  bool return_if_nonlin = false);
1746 
1747 
1748  /** Add a Laplacian term on the variable `varname` (in fact with a minus :
1749  :math:`-\text{div}(\nabla u)`). If it is a vector
1750  valued variable, the Laplacian term is componentwise. `region` is an
1751  optional mesh region on which the term is added. Return the brick index
1752  in the model.
1753  */
1755  (model &md, const mesh_im &mim, const std::string &varname,
1756  size_type region = size_type(-1));
1757 
1758 
1759  /** Add an elliptic term on the variable `varname`.
1760  The shape of the elliptic
1761  term depends both on the variable and the data. This corresponds to a
1762  term $-\text{div}(a\nabla u)$ where $a$ is the data and $u$ the variable.
1763  The data can be a scalar, a matrix or an order four tensor. The variable
1764  can be vector valued or not. If the data is a scalar or a matrix and
1765  the variable is vector valued then the term is added componentwise.
1766  An order four tensor data is allowed for vector valued variable only.
1767  The data can be constant or describbed on a fem. Of course, when
1768  the data is a tensor describe on a finite element method (a tensor
1769  field) the data can be a huge vector. The components of the
1770  matrix/tensor have to be stored with the fortran order (columnwise) in
1771  the data vector (compatibility with blas). The symmetry and coercivity
1772  of the given matrix/tensor is not verified (but assumed). `region` is an
1773  optional mesh region on which the term is added. Note that for the real
1774  version which uses the high-level generic assembly language, `dataexpr`
1775  can be any regular expression of the high-level generic assembly
1776  language (like "1", "sin(X[0])" or "Norm(u)" for instance) even
1777  depending on model variables.
1778  Return the brick index in the model.
1779  */
1781  (model &md, const mesh_im &mim, const std::string &varname,
1782  const std::string &dataexpr, size_type region = size_type(-1));
1783 
1784 
1785  /** Add a source term on the variable `varname`. The source term is
1786  represented by `dataexpr` which could be any regular expression of the
1787  high-level generic assembly language (except for the complex version
1788  where it has to be a declared data of the model). `region` is an
1789  optional mesh region on which the term is
1790  added. An additional optional data `directdataname` can be provided. The
1791  corresponding data vector will be directly added to the right hand
1792  side without assembly. Return the brick index in the model.
1793  */
1795  (model &md, const mesh_im &mim, const std::string &varname,
1796  const std::string &dataexpr, size_type region = size_type(-1),
1797  const std::string &directdataname = std::string());
1798 
1799  /** Add a source term on the variable `varname` on a boundary `region`.
1800  The source term is
1801  represented by the data `dataepxpr` which could be any regular
1802  expression of the high-level generic assembly language (except
1803  for the complex version where it has to be a declared data of
1804  the model). A scalar product with the outward normal unit vector to
1805  the boundary is performed. The main aim of this brick is to represent
1806  a Neumann condition with a vector data without performing the
1807  scalar product with the normal as a pre-processing.
1808  */
1810  (model &md, const mesh_im &mim, const std::string &varname,
1811  const std::string &dataexpr, size_type region);
1812 
1813 
1814  /** Add a (simple) Dirichlet condition on the variable `varname` and
1815  the mesh region `region`. The Dirichlet condition is prescribed by
1816  a simple post-treatment of the final linear system (tangent system
1817  for nonlinear problems) consisting of modifying the lines corresponding
1818  to the degree of freedom of the variable on `region` (0 outside the
1819  diagonal, 1 on the diagonal of the matrix and the expected value on
1820  the right hand side).
1821  The symmetry of the linear system is kept if all other bricks are
1822  symmetric.
1823  This brick is to be reserved for simple Dirichlet conditions (only dof
1824  declared on the correspodning boundary are prescribed). The application
1825  of this brick on reduced f.e.m. may be problematic. Intrinsic vectorial
1826  finite element method are not supported.
1827  `dataname` is the optional right hand side of the Dirichlet condition.
1828  It could be constant or (important) described on the same finite
1829  element method as `varname`.
1830  Returns the brick index in the model.
1831  */
1833  (model &md, const std::string &varname, size_type region,
1834  const std::string &dataname = std::string());
1835 
1836 
1837  /** Add a Dirichlet condition on the variable `varname` and the mesh
1838  region `region`. This region should be a boundary. The Dirichlet
1839  condition is prescribed with a multiplier variable `multname` which
1840  should be first declared as a multiplier
1841  variable on the mesh region in the model. `dataname` is the optional
1842  right hand side of the Dirichlet condition. It could be constant or
1843  described on a fem; scalar or vector valued, depending on the variable
1844  on which the Dirichlet condition is prescribed. Return the brick index
1845  in the model.
1846  */
1848  (model &md, const mesh_im &mim, const std::string &varname,
1849  const std::string &multname, size_type region,
1850  const std::string &dataname = std::string());
1851 
1852  /** Same function as the previous one but the multipliers variable will
1853  be declared to the brick by the function. `mf_mult` is the finite element
1854  method on which the multiplier will be build (it will be restricted to
1855  the mesh region `region` and eventually some conflicting dofs with some
1856  other multiplier variables will be suppressed).
1857  */
1859  (model &md, const mesh_im &mim, const std::string &varname,
1860  const mesh_fem &mf_mult, size_type region,
1861  const std::string &dataname = std::string());
1862 
1863  /** Same function as the previous one but the `mf_mult` parameter is
1864  replaced by `degree`. The multiplier will be described on a standard
1865  finite element method of the corresponding degree.
1866  */
1868  (model &md, const mesh_im &mim, const std::string &varname,
1869  dim_type degree, size_type region,
1870  const std::string &dataname = std::string());
1871 
1872 
1873  /** When `ind_brick` is the index of a Dirichlet brick with multiplier on
1874  the model `md`, the function return the name of the multiplier variable.
1875  Otherwise, it has an undefined behavior.
1876  */
1877  const APIDECL std::string &mult_varname_Dirichlet(model &md, size_type ind_brick);
1878 
1879  /** Add a Dirichlet condition on the variable `varname` and the mesh
1880  region `region`. This region should be a boundary. The Dirichlet
1881  condition is prescribed with penalization. The penalization coefficient
1882  is intially `penalization_coeff` and will be added to the data of
1883  the model. `dataname` is the optional
1884  right hand side of the Dirichlet condition. It could be constant or
1885  described on a fem; scalar or vector valued, depending on the variable
1886  on which the Dirichlet condition is prescribed.
1887  `mf_mult` is an optional parameter which allows to weaken the
1888  Dirichlet condition specifying a multiplier space.
1889  Returns the brick index in the model.
1890  */
1892  (model &md, const mesh_im &mim, const std::string &varname,
1893  scalar_type penalization_coeff, size_type region,
1894  const std::string &dataname = std::string(),
1895  const mesh_fem *mf_mult = 0);
1896 
1897  /** Add a Dirichlet condition on the variable `varname` and the mesh
1898  region `region`. This region should be a boundary. `Neumannterm`
1899  is the expression of the Neumann term (obtained by the Green formula)
1900  described as an expression of the high-level
1901  generic assembly language. This term can be obtained with
1902  md. Neumann_term(varname, region) once all volumic bricks have
1903  been added to the model. The Dirichlet
1904  condition is prescribed with Nitsche's method. `datag` is the optional
1905  right hand side of the Dirichlet condition. `datagamma0` is the
1906  Nitsche's method parameter. `theta` is a scalar value which can be
1907  positive or negative. `theta = 1` corresponds to the standard symmetric
1908  method which is conditionnaly coercive for `gamma0` small.
1909  `theta = -1` corresponds to the skew-symmetric method which is
1910  inconditionnaly coercive. `theta = 0` is the simplest method
1911  for which the second derivative of the Neumann term is not necessary
1912  even for nonlinear problems. Return the brick index in the model.
1913  */
1915  (model &md, const mesh_im &mim, const std::string &varname,
1916  const std::string &Neumannterm,
1917  const std::string &datagamma0, size_type region,
1918  scalar_type theta = scalar_type(0),
1919  const std::string &datag = std::string());
1920 
1921 
1922  /** Add a Dirichlet condition to the normal component of the vector
1923  (or tensor) valued variable `varname` and the mesh
1924  region `region`. This region should be a boundary. The Dirichlet
1925  condition is prescribed with a multiplier variable `multname` which
1926  should be first declared as a multiplier
1927  variable on the mesh region in the model. `dataname` is the optional
1928  right hand side of the normal Dirichlet condition.
1929  It could be constant or
1930  described on a fem; scalar or vector valued, depending on the variable
1931  on which the Dirichlet condition is prescribed (scalar if the variable
1932  is vector valued, vector if the variable is tensor valued).
1933  Returns the brick index in the model.
1934  */
1936  (model &md, const mesh_im &mim, const std::string &varname,
1937  const std::string &multname, size_type region,
1938  const std::string &dataname = std::string());
1939 
1940  /** Same function as the previous one but the multipliers variable will
1941  be declared to the brick by the function. `mf_mult` is the finite element
1942  method on which the multiplier will be build (it will be restricted to
1943  the mesh region `region` and eventually some conflicting dofs with some
1944  other multiplier variables will be suppressed).
1945  */
1947  (model &md, const mesh_im &mim, const std::string &varname,
1948  const mesh_fem &mf_mult, size_type region,
1949  const std::string &dataname = std::string());
1950 
1951  /** Same function as the previous one but the `mf_mult` parameter is
1952  replaced by `degree`. The multiplier will be described on a standard
1953  finite element method of the corresponding degree.
1954  */
1956  (model &md, const mesh_im &mim, const std::string &varname,
1957  dim_type degree, size_type region,
1958  const std::string &dataname = std::string());
1959 
1960  /** Add a Dirichlet condition to the normal component of the vector
1961  (or tensor) valued variable `varname` and the mesh
1962  region `region`. This region should be a boundary. The Dirichlet
1963  condition is prescribed with penalization. The penalization coefficient
1964  is intially `penalization_coeff` and will be added to the data of
1965  the model. `dataname` is the optional
1966  right hand side of the Dirichlet condition. It could be constant or
1967  described on a fem; scalar or vector valued, depending on the variable
1968  on which the Dirichlet condition is prescribed (scalar if the variable
1969  is vector valued, vector if the variable is tensor valued).
1970  `mf_mult` is an optional parameter which allows to weaken the
1971  Dirichlet condition specifying a multiplier space.
1972  Return the brick index in the model.
1973  */
1975  (model &md, const mesh_im &mim, const std::string &varname,
1976  scalar_type penalization_coeff, size_type region,
1977  const std::string &dataname = std::string(),
1978  const mesh_fem *mf_mult = 0);
1979 
1980 
1981 
1982  /** Add a Dirichlet condition on the normal component of the variable
1983  `varname` and the mesh
1984  region `region`. This region should be a boundary. `Neumannterm`
1985  is the expression of the Neumann term (obtained by the Green formula)
1986  described as an expression of the high-level
1987  generic assembly language. This term can be obtained with
1988  md.Neumann_term(varname, region) once all volumic bricks have
1989  been added to the model.The Dirichlet
1990  condition is prescribed with Nitsche's method. `datag` is the optional
1991  scalar right hand side of the Dirichlet condition. `datagamma0` is the
1992  Nitsche's method parameter. `theta` is a scalar value which can be
1993  positive or negative. `theta = 1` corresponds to the standard symmetric
1994  method which is conditionnaly coercive for `gamma0` small.
1995  `theta = -1` corresponds to the skew-symmetric method which is
1996  inconditionnaly coercive. `theta = 0` is the simplest method
1997  for which the second derivative of the Neumann term is not necessary
1998  even for nonlinear problems. Return the brick index in the model.
1999  */
2001  (model &md, const mesh_im &mim, const std::string &varname,
2002  const std::string &Neumannterm, const std::string &datagamma0,
2003  size_type region, scalar_type theta = scalar_type(0),
2004  const std::string &datag = std::string());
2005 
2006 
2007  /** Add some pointwise constraints on the variable `varname` thanks to
2008  a penalization. The penalization coefficient is initially
2009  `penalization_coeff` and will be added to the data of the model.
2010  The conditions are prescribed on a set of points given in the data
2011  `dataname_pt` whose dimension is the number of points times the dimension
2012  of the mesh. If the variable represents a vector field, the data
2013  `dataname_unitv` represents a vector of dimension the number of points
2014  times the dimension of the vector field which should store some
2015  unit vectors. In that case the prescribed constraint is the scalar
2016  product of the variable at the corresponding point with the corresponding
2017  unit vector.
2018  The optional data `dataname_val` is the vector of values to be prescribed
2019  at the different points.
2020  This brick is specifically designed to kill rigid displacement
2021  in a Neumann problem.
2022  */
2024  (model &md, const std::string &varname,
2025  scalar_type penalisation_coeff, const std::string &dataname_pt,
2026  const std::string &dataname_unitv = std::string(),
2027  const std::string &dataname_val = std::string());
2028 
2029 
2030  /** Add some pointwise constraints on the variable `varname` using a given
2031  multiplier `multname`.
2032  The conditions are prescribed on a set of points given in the data
2033  `dataname_pt` whose dimension is the number of points times the dimension
2034  of the mesh.
2035  The multiplier variable should be a fixed size variable of size the
2036  number of points.
2037  If the variable represents a vector field, the data
2038  `dataname_unitv` represents a vector of dimension the number of points
2039  times the dimension of the vector field which should store some
2040  unit vectors. In that case the prescribed constraint is the scalar
2041  product of the variable at the corresponding point with the corresponding
2042  unit vector.
2043  The optional data `dataname_val` is the vector of values to be prescribed
2044  at the different points.
2045  This brick is specifically designed to kill rigid displacement
2046  in a Neumann problem.
2047  */
2049  (model &md, const std::string &varname,
2050  const std::string &multname, const std::string &dataname_pt,
2051  const std::string &dataname_unitv = std::string(),
2052  const std::string &dataname_val = std::string());
2053 
2054  /** Add some pointwise constraints on the variable `varname` using
2055  multiplier. The multiplier variable is automatically added to the model.
2056  The conditions are prescribed on a set of points given in the data
2057  `dataname_pt` whose dimension is the number of points times the dimension
2058  of the mesh.
2059  If the variable represents a vector field, the data
2060  `dataname_unitv` represents a vector of dimension the number of points
2061  times the dimension of the vector field which should store some
2062  unit vectors. In that case the prescribed constraint is the scalar
2063  product of the variable at the corresponding point with the corresponding
2064  unit vector.
2065  The optional data `dataname_val` is the vector of values to be prescribed
2066  at the different points.
2067  This brick is specifically designed to kill rigid displacement
2068  in a Neumann problem.
2069  */
2071  (model &md, const std::string &varname, const std::string &dataname_pt,
2072  const std::string &dataname_unitv = std::string(),
2073  const std::string &dataname_val = std::string());
2074 
2075 
2076  /** Change the penalization coefficient of a Dirichlet condition with
2077  penalization brick. If the brick is not of this kind,
2078  this function has an undefined behavior.
2079  */
2080  void APIDECL change_penalization_coeff(model &md, size_type ind_brick,
2081  scalar_type penalisation_coeff);
2082 
2083  /** Add a generalized Dirichlet condition on the variable `varname` and
2084  the mesh region `region`. This version is for vector field.
2085  It prescribes a condition @f$ Hu = r @f$ where `H` is a matrix field.
2086  This region should be a boundary. The Dirichlet
2087  condition is prescribed with a multiplier variable `multname` which
2088  should be first declared as a multiplier
2089  variable on the mesh region in the model. `dataname` is the
2090  right hand side of the Dirichlet condition. It could be constant or
2091  described on a fem; scalar or vector valued, depending on the variable
2092  on which the Dirichlet condition is prescribed. `Hname' is the data
2093  corresponding to the matrix field `H`. It has to be a constant matrix
2094  or described on a scalar fem. Return the brick index in the model.
2095  */
2097  (model &md, const mesh_im &mim, const std::string &varname,
2098  const std::string &multname, size_type region,
2099  const std::string &dataname, const std::string &Hname);
2100 
2101  /** Same function as the preceeding one but the multipliers variable will
2102  be declared to the brick by the function. `mf_mult` is the finite element
2103  method on which the multiplier will be build (it will be restricted to
2104  the mesh region `region` and eventually some conflicting dofs with some
2105  other multiplier variables will be suppressed).
2106  */
2108  (model &md, const mesh_im &mim, const std::string &varname,
2109  const mesh_fem &mf_mult, size_type region,
2110  const std::string &dataname, const std::string &Hname);
2111 
2112  /** Same function as the preceeding one but the `mf_mult` parameter is
2113  replaced by `degree`. The multiplier will be described on a standard
2114  finite element method of the corresponding degree.
2115  */
2117  (model &md, const mesh_im &mim, const std::string &varname,
2118  dim_type degree, size_type region,
2119  const std::string &dataname, const std::string &Hname);
2120 
2121  /** Add a Dirichlet condition on the variable `varname` and the mesh
2122  region `region`. This version is for vector field.
2123  It prescribes a condition @f$ Hu = r @f$ where `H` is a matrix field.
2124  This region should be a boundary. This region should be a boundary.
2125  The Dirichlet
2126  condition is prescribed with penalization. The penalization coefficient
2127  is intially `penalization_coeff` and will be added to the data of
2128  the model. `dataname` is the
2129  right hand side of the Dirichlet condition. It could be constant or
2130  described on a fem; scalar or vector valued, depending on the variable
2131  on which the Dirichlet condition is prescribed. `Hname' is the data
2132  corresponding to the matrix field `H`. It has to be a constant matrix
2133  or described on a scalar fem. `mf_mult` is an optional parameter
2134  which allows to weaken the Dirichlet condition specifying a
2135  multiplier space. Return the brick index in the model.
2136  */
2138  (model &md, const mesh_im &mim, const std::string &varname,
2139  scalar_type penalization_coeff, size_type region,
2140  const std::string &dataname, const std::string &Hname,
2141  const mesh_fem *mf_mult = 0);
2142 
2143  /** Add a Dirichlet condition on the variable `varname` and the mesh
2144  region `region`. This region should be a boundary. This version
2145  is for vector field. It prescribes a condition
2146  @f$ Hu = r @f$ where `H` is a matrix field. `Neumannterm`
2147  is the expression of the Neumann term (obtained by the Green formula)
2148  described as an expression of the high-level
2149  generic assembly language. of the high-level
2150  generic assembly language. This term can be obtained with
2151  md.Neumann_term(varname, region) once all volumic bricks have
2152  been added to the model. The Dirichlet
2153  condition is prescribed with Nitsche's method. `datag` is the optional
2154  right hand side of the Dirichlet condition. `datagamma0` is the
2155  Nitsche's method parameter. `theta` is a scalar value which can be
2156  positive or negative. `theta = 1` corresponds to the standard symmetric
2157  method which is conditionnaly coercive for `gamma0` small.
2158  `theta = -1` corresponds to the skew-symmetric method which is
2159  inconditionnaly coercive. `theta = 0` is the simplest method
2160  for which the second derivative of the Neumann term is not necessary
2161  even for nonlinear problems. Return the brick index in the model.
2162  */
2164  (model &md, const mesh_im &mim, const std::string &varname,
2165  const std::string &Neumannterm, const std::string &datagamma0,
2166  size_type region, scalar_type theta,
2167  const std::string &datag, const std::string &dataH);
2168 
2169 
2170  /** Add a Helmoltz brick to the model. This corresponds to the scalar
2171  equation (@f$\Delta u + k^2u = 0@f$, with @f$K=k^2@f$).
2172  The weak formulation is (@f$\int k^2 u.v - \nabla u.\nabla v@f$)
2173 
2174  `dataexpr` should contain the wave number $k$. It can be real or
2175  complex.
2176  */
2177  size_type APIDECL add_Helmholtz_brick(model &md, const mesh_im &mim,
2178  const std::string &varname,
2179  const std::string &dataexpr,
2180  size_type region = size_type(-1));
2181 
2182 
2183  /** Add a Fourier-Robin brick to the model. This correspond to the weak term
2184  (@f$\int (qu).v @f$) on a boundary. It is used to represent a
2185  Fourier-Robin boundary condition.
2186 
2187  `dataexpr` is the parameter $q$ which should be a
2188  (@f$N\times N@f$) matrix term, where $N$ is the dimension of the
2189  variable `varname`. It can be an arbitrary valid expression of the
2190  high-level generic assembly language (except for the complex version
2191  for which it should be a data of the model). Note that an additional
2192  right hand side can be added with a source term brick.
2193  */
2194  size_type APIDECL add_Fourier_Robin_brick(model &md, const mesh_im &mim,
2195  const std::string &varname,
2196  const std::string &dataexpr,
2197  size_type region);
2198 
2199 
2200  // Constraint brick.
2201  model_real_sparse_matrix APIDECL &set_private_data_brick_real_matrix
2202  (model &md, size_type indbrick);
2203  model_real_plain_vector APIDECL &set_private_data_brick_real_rhs
2204  (model &md, size_type indbrick);
2205  model_complex_sparse_matrix APIDECL &set_private_data_brick_complex_matrix
2206  (model &md, size_type indbrick);
2207  model_complex_plain_vector APIDECL &set_private_data_brick_complex_rhs
2208  (model &md, size_type indbrick);
2209  size_type APIDECL add_constraint_with_penalization
2210  (model &md, const std::string &varname, scalar_type penalisation_coeff);
2211  size_type APIDECL add_constraint_with_multipliers
2212  (model &md, const std::string &varname, const std::string &multname);
2213 
2214  void set_private_data_rhs
2215  (model &md, size_type indbrick, const std::string &varname);
2216 
2217  template <typename VECT, typename T>
2218  void set_private_data_rhs(model &md, size_type ind,
2219  const VECT &L, T) {
2220  model_real_plain_vector &LL = set_private_data_brick_real_rhs(md, ind);
2221  gmm::resize(LL, gmm::vect_size(L));
2222  gmm::copy(L, LL);
2223  }
2224 
2225  template <typename VECT, typename T>
2226  void set_private_data_rhs(model &md, size_type ind, const VECT &L,
2227  std::complex<T>) {
2228  model_complex_plain_vector &LL = set_private_data_brick_complex_rhs(md, ind);
2229  gmm::resize(LL, gmm::vect_size(L));
2230  gmm::copy(L, LL);
2231  }
2232 
2233  /** For some specific bricks having an internal right hand side vector
2234  (explicit bricks: 'constraint brick' and 'explicit rhs brick'),
2235  set this rhs.
2236  */
2237  template <typename VECT>
2238  void set_private_data_rhs(model &md, size_type indbrick, const VECT &L) {
2239  typedef typename gmm::linalg_traits<VECT>::value_type T;
2240  set_private_data_rhs(md, indbrick, L, T());
2241  }
2242 
2243  template <typename MAT, typename T>
2244  void set_private_data_matrix(model &md, size_type ind,
2245  const MAT &B, T) {
2246  model_real_sparse_matrix &BB = set_private_data_brick_real_matrix(md, ind);
2247  gmm::resize(BB, gmm::mat_nrows(B), gmm::mat_ncols(B));
2248  gmm::copy(B, BB);
2249  }
2250 
2251  template <typename MAT, typename T>
2252  void set_private_data_matrix(model &md, size_type ind, const MAT &B,
2253  std::complex<T>) {
2254  model_complex_sparse_matrix &BB
2255  = set_private_data_brick_complex_matrix(md, ind);
2256  gmm::resize(BB, gmm::mat_nrows(B), gmm::mat_ncols(B));
2257  gmm::copy(B, BB);
2258  }
2259 
2260  /** For some specific bricks having an internal sparse matrix
2261  (explicit bricks: 'constraint brick' and 'explicit matrix brick'),
2262  set this matrix. @*/
2263  template <typename MAT>
2264  void set_private_data_matrix(model &md, size_type indbrick,
2265  const MAT &B) {
2266  typedef typename gmm::linalg_traits<MAT>::value_type T;
2267  set_private_data_matrix(md, indbrick, B, T());
2268  }
2269 
2270  /** Add an additional explicit penalized constraint on the variable
2271  `varname`. The constraint is $BU=L$ with `B` being a rectangular
2272  sparse matrix.
2273  Be aware that `B` should not contain a plain row, otherwise the whole
2274  tangent matrix will be plain. It is possible to change the constraint
2275  at any time with the methods set_private_matrix and set_private_rhs.
2276  The method change_penalization_coeff can also be used.
2277  */
2278  template <typename MAT, typename VECT>
2279  size_type add_constraint_with_penalization
2280  (model &md, const std::string &varname, scalar_type penalisation_coeff,
2281  const MAT &B, const VECT &L) {
2282  size_type ind
2283  = add_constraint_with_penalization(md, varname, penalisation_coeff);
2284  size_type n = gmm::mat_nrows(B), m = gmm::mat_ncols(B);
2285  set_private_data_rhs(md, ind, L);
2286  set_private_data_matrix(md, ind, B);
2287  return ind;
2288  }
2289 
2290  /** Add an additional explicit constraint on the variable `varname` thank to
2291  a multiplier `multname` peviously added to the model (should be a fixed
2292  size variable).
2293  The constraint is $BU=L$ with `B` being a rectangular sparse matrix.
2294  It is possible to change the constraint
2295  at any time with the methods set_private_matrix
2296  and set_private_rhs.
2297  */
2298  template <typename MAT, typename VECT>
2299  size_type add_constraint_with_multipliers
2300  (model &md, const std::string &varname, const std::string &multname,
2301  const MAT &B, const VECT &L) {
2302  size_type ind = add_constraint_with_multipliers(md, varname, multname);
2303  set_private_data_rhs(md, ind, L);
2304  set_private_data_matrix(md, ind, B);
2305  return ind;
2306  }
2307 
2308  template <typename MAT>
2309  size_type add_constraint_with_multipliers
2310  (model &md, const std::string &varname, const std::string &multname,
2311  const MAT &B, const std::string &Lname) {
2312  size_type ind = add_constraint_with_multipliers(md, varname, multname);
2313  set_private_data_rhs(md, ind, Lname);
2314  set_private_data_matrix(md, ind, B);
2315  return ind;
2316  }
2317 
2318  size_type APIDECL add_explicit_matrix(model &md, const std::string &varname1,
2319  const std::string &varname2,
2320  bool issymmetric, bool iscoercive);
2321  size_type APIDECL add_explicit_rhs(model &md, const std::string &varname);
2322 
2323  /** Add a brick reprenting an explicit matrix to be added to the tangent
2324  linear system relatively to the variables 'varname1' and 'varname2'.
2325  The given matrix should have as many rows as the dimension of
2326  'varname1' and as many columns as the dimension of 'varname2'.
2327  If the two variables are different and if `issymmetric' is set to true
2328  then the transpose of the matrix is also added to the tangent system
2329  (default is false). set `iscoercive` to true if the term does not
2330  affect the coercivity of the tangent system (default is false).
2331  The matrix can be changed by the command set_private_matrix.
2332  */
2333  template <typename MAT>
2334  size_type add_explicit_matrix(model &md, const std::string &varname1,
2335  const std::string &varname2, const MAT &B,
2336  bool issymmetric = false,
2337  bool iscoercive = false) {
2338  size_type ind = add_explicit_matrix(md, varname1, varname2,
2339  issymmetric, iscoercive);
2340  set_private_data_matrix(md, ind, B);
2341  return ind;
2342  }
2343 
2344  /** Add a brick representing an explicit right hand side to be added to
2345  the right hand side of the tangent
2346  linear system relatively to the variable 'varname'.
2347  The given rhs should have the same size than the dimension of
2348  'varname'. The rhs can be changed by the command set_private_rhs.
2349  */
2350  template <typename VECT>
2351  size_type add_explicit_rhs(model &md, const std::string &varname,
2352  const VECT &L) {
2353  size_type ind = add_explicit_rhs(md, varname);
2354  set_private_data_rhs(md, ind, L);
2355  return ind;
2356  }
2357 
2358 
2359  /** Linear elasticity brick ( @f$ \int \sigma(u):\varepsilon(v) @f$ ).
2360  for isotropic material. Parametrized by the Lamé coefficients
2361  lambda and mu.
2362  */
2364  (model &md, const mesh_im &mim, const std::string &varname,
2365  const std::string &dataname_lambda, const std::string &dataname_mu,
2366  size_type region = size_type(-1),
2367  const std::string &dataname_preconstraint = std::string());
2368 
2369  /** Linear elasticity brick ( @f$ \int \sigma(u):\varepsilon(v) @f$ ).
2370  for isotropic material. Parametrized by Young modulus and Poisson ratio
2371  For two-dimensional problems, corresponds to the plane strain
2372  approximation
2373  ( @f$ \lambda = E\nu/((1+\nu)(1-2\nu)), \mu = E/(2(1+\nu)) @f$ ).
2374  Corresponds to the standard model for three-dimensional problems.
2375  */
2377  (model &md, const mesh_im &mim, const std::string &varname,
2378  const std::string &data_E, const std::string &data_nu,
2379  size_type region);
2380 
2381  /**
2382  Linear elasticity brick ( @f$ \int \sigma(u):\varepsilon(v) @f$ ).
2383  for isotropic material. Parametrized by Young modulus and Poisson ratio.
2384  For two-dimensional problems, corresponds to the plane stress
2385  approximation
2386  ( @f$ \lambda^* = E\nu/(1-\nu^2), \mu = E/(2(1+\nu)) @f$ ).
2387  Corresponds to the standard model for three-dimensional problems.
2388  */
2390  (model &md, const mesh_im &mim, const std::string &varname,
2391  const std::string &data_E, const std::string &data_nu,
2392  size_type region);
2393 
2394  void APIDECL compute_isotropic_linearized_Von_Mises_or_Tresca
2395  (model &md, const std::string &varname, const std::string &dataname_lambda,
2396  const std::string &dataname_mu, const mesh_fem &mf_vm,
2397  model_real_plain_vector &VM, bool tresca);
2398 
2399  /**
2400  Compute the Von-Mises stress or the Tresca stress of a field
2401  (only valid for isotropic linearized elasticity in 3D)
2402  Parametrized by Lame coefficients.
2403  */
2404  template <class VECTVM>
2405  void compute_isotropic_linearized_Von_Mises_or_Tresca
2406  (model &md, const std::string &varname, const std::string &dataname_lambda,
2407  const std::string &dataname_mu, const mesh_fem &mf_vm,
2408  VECTVM &VM, bool tresca) {
2409  model_real_plain_vector VMM(mf_vm.nb_dof());
2410  compute_isotropic_linearized_Von_Mises_or_Tresca
2411  (md, varname, dataname_lambda, dataname_mu, mf_vm, VMM, tresca);
2412  gmm::copy(VMM, VM);
2413  }
2414 
2415  /**
2416  Compute the Von-Mises stress of a displacement field for isotropic
2417  linearized elasticity in 3D or in 2D with plane strain assumption.
2418  Parametrized by Young modulus and Poisson ratio.
2419  */
2421  (model &md, const std::string &varname, const std::string &data_E,
2422  const std::string &data_nu, const mesh_fem &mf_vm,
2423  model_real_plain_vector &VM);
2424 
2425  /**
2426  Compute the Von-Mises stress of a displacement field for isotropic
2427  linearized elasticity in 3D or in 2D with plane stress assumption.
2428  Parametrized by Young modulus and Poisson ratio.
2429  */
2431  (model &md, const std::string &varname, const std::string &data_E,
2432  const std::string &data_nu, const mesh_fem &mf_vm,
2433  model_real_plain_vector &VM);
2434 
2435 
2436  /**
2437  Mixed linear incompressibility condition brick.
2438 
2439  Update the tangent matrix with a pressure term:
2440  @f[
2441  T \longrightarrow
2442  \begin{array}{ll} T & B \\ B^t & M \end{array}
2443  @f]
2444  with @f$ B = - \int p.div u @f$ and
2445  @f$ M = \int \epsilon p.q @f$ ( @f$ \epsilon @f$ is an optional
2446  penalization coefficient).
2447 
2448  Be aware that an inf-sup condition between the finite element method
2449  describing the rpressure and the primal variable has to be satisfied.
2450 
2451  For nearly incompressible elasticity,
2452  @f[ p = -\lambda \textrm{div}~u @f]
2453  @f[ \sigma = 2 \mu \varepsilon(u) -p I @f]
2454  */
2456  (model &md, const mesh_im &mim, const std::string &varname,
2457  const std::string &multname_pressure, size_type region = size_type(-1),
2458  const std::string &dataexpr_penal_coeff = std::string());
2459 
2460  /** Mass brick ( @f$ \int \rho u.v @f$ ).
2461  Add a mass matix on a variable (eventually with a specified region).
2462  If the parameter $\rho$ is omitted it is assumed to be equal to 1.
2463  */
2464  size_type APIDECL add_mass_brick
2465  (model &md, const mesh_im &mim, const std::string &varname,
2466  const std::string &dataexpr_rho = std::string(),
2467  size_type region = size_type(-1));
2468 
2469  /** Lumped mass brick for first order.
2470  Add a lumped mass matix for first order on a variable (eventually with a specified region).
2471  If the parameter $\rho$ is omitted it is assumed to be equal to 1.
2472  */
2474  (model &md, const mesh_im &mim, const std::string &varname,
2475  const std::string &dataexpr_rho = std::string(),
2476  size_type region = size_type(-1));
2477 
2478  /** Basic d/dt brick ( @f$ \int \rho ((u^{n+1}-u^n)/dt).v @f$ ).
2479  Add the standard discretization of a first order time derivative. The
2480  parameter $rho$ is the density which could be omitted (the defaul value
2481  is 1). This brick should be used in addition to a time dispatcher for the
2482  other terms.
2483  */
2485  (model &md, const mesh_im &mim, const std::string &varname,
2486  const std::string &dataname_dt,
2487  const std::string &dataname_rho = std::string(),
2488  size_type region = size_type(-1));
2489 
2490  /** Basic d2/dt2 brick ( @f$ \int \rho ((u^{n+1}-u^n)/(\alpha dt^2) - v^n/(\alpha dt) ).w @f$ ).
2491  Add the standard discretization of a second order time derivative. The
2492  parameter $rho$ is the density which could be omitted (the defaul value
2493  is 1). This brick should be used in addition to a time dispatcher for the
2494  other terms. The time derivative $v$ of the variable $u$ is preferably
2495  computed as a post-traitement which depends on each scheme.
2496  */
2498  (model &md, const mesh_im &mim, const std::string &varnameU,
2499  const std::string &datanameV,
2500  const std::string &dataname_dt,
2501  const std::string &dataname_alpha,
2502  const std::string &dataname_rho = std::string(),
2503  size_type region = size_type(-1));
2504 
2505 
2506 } /* end of namespace getfem. */
2507 
2508 
2509 #endif /* GETFEM_MODELS_H_*/
getfem::model::is_complex
bool is_complex() const
Boolean which says if the model deals with real or complex unknowns and data.
Definition: getfem_models.h:562
getfem::change_penalization_coeff
void APIDECL change_penalization_coeff(model &md, size_type ind_brick, scalar_type penalisation_coeff)
Change the penalization coefficient of a Dirichlet condition with penalization brick.
Definition: getfem_models.cc:4827
getfem::model::add_elementary_transformation
void add_elementary_transformation(const std::string &name, pelementary_transformation ptrans)
Add an elementary transformation to the model to be used with the generic assembly.
Definition: getfem_models.h:1132
getfem::add_nonlinear_twodomain_term
size_type APIDECL add_nonlinear_twodomain_term(model &md, const mesh_im &mim, const std::string &expr, size_type region, const std::string &secondary_domain, bool is_sym=false, bool is_coercive=false, const std::string &brickname="")
Adds a nonlinear term given by a weak form language expression like add_nonlinear_term function but f...
Definition: getfem_models.cc:3670
getfem::mult_varname_Dirichlet
const APIDECL std::string & mult_varname_Dirichlet(model &md, size_type ind_brick)
When ind_brick is the index of a Dirichlet brick with multiplier on the model md, the function return...
Definition: getfem_models.cc:4690
getfem::virtual_dispatcher
The time dispatcher object modify the result of a brick in order to apply a time integration scheme.
Definition: getfem_models.h:1260
getfem::add_normal_Dirichlet_condition_with_multipliers
size_type APIDECL add_normal_Dirichlet_condition_with_multipliers(model &md, const mesh_im &mim, const std::string &varname, const std::string &multname, size_type region, const std::string &dataname=std::string())
Add a Dirichlet condition to the normal component of the vector (or tensor) valued variable varname a...
Definition: getfem_models.cc:4715
getfem::compute_isotropic_linearized_Von_Mises_pstrain
void APIDECL compute_isotropic_linearized_Von_Mises_pstrain(model &md, const std::string &varname, const std::string &data_E, const std::string &data_nu, const mesh_fem &mf_vm, model_real_plain_vector &VM)
Compute the Von-Mises stress of a displacement field for isotropic linearized elasticity in 3D or in ...
Definition: getfem_models.cc:6241
getfem::add_generalized_Dirichlet_condition_with_penalization
size_type APIDECL add_generalized_Dirichlet_condition_with_penalization(model &md, const mesh_im &mim, const std::string &varname, scalar_type penalization_coeff, size_type region, const std::string &dataname, const std::string &Hname, const mesh_fem *mf_mult=0)
Add a Dirichlet condition on the variable varname and the mesh region region.
Definition: getfem_models.cc:4806
getfem::add_linear_twodomain_term
size_type APIDECL add_linear_twodomain_term(model &md, const mesh_im &mim, const std::string &expr, size_type region, const std::string &secondary_domain, bool is_sym=false, bool is_coercive=false, const std::string &brickname="", bool return_if_nonlin=false)
Adds a linear term given by a weak form language expression like add_linear_term function but for an ...
Definition: getfem_models.cc:3579
getfem::model::macro_exists
bool macro_exists(const std::string &name) const
Says if a macro of that name has been defined.
Definition: getfem_models.h:898
getfem::model::elementary_transformation
pelementary_transformation elementary_transformation(const std::string &name) const
Get a pointer to the elementary transformation name.
Definition: getfem_models.h:1140
getfem::add_midpoint_dispatcher
void APIDECL add_midpoint_dispatcher(model &md, dal::bit_vector ibricks)
Add a midpoint time dispatcher to a list of bricks.
Definition: getfem_models.cc:7335
gmm::resize
void resize(M &v, size_type m, size_type n)
*‍/
Definition: gmm_blas.h:231
bgeot::size_type
size_t size_type
used as the common size type in the library
Definition: bgeot_poly.h:49
getfem::model::is_linear
bool is_linear() const
Return true if all the model terms are linear.
Definition: getfem_models.h:580
getfem::model::disable_brick
void disable_brick(size_type ib)
Disable a brick.
Definition: getfem_models.h:512
getfem::mesh_im
Describe an integration method linked to a mesh.
Definition: getfem_mesh_im.h:47
getfem::add_pointwise_constraints_with_penalization
size_type APIDECL add_pointwise_constraints_with_penalization(model &md, const std::string &varname, scalar_type penalisation_coeff, const std::string &dataname_pt, const std::string &dataname_unitv=std::string(), const std::string &dataname_val=std::string())
Add some pointwise constraints on the variable varname thanks to a penalization.
Definition: getfem_models.cc:5347
getfem::add_Dirichlet_condition_with_multipliers
size_type APIDECL add_Dirichlet_condition_with_multipliers(model &md, const mesh_im &mim, const std::string &varname, const std::string &multname, size_type region, const std::string &dataname=std::string())
Add a Dirichlet condition on the variable varname and the mesh region region.
Definition: getfem_models.cc:4656
getfem::add_linear_term
size_type APIDECL add_linear_term(model &md, const mesh_im &mim, const std::string &expr, size_type region=size_type(-1), bool is_sym=false, bool is_coercive=false, const std::string &brickname="", bool return_if_nonlin=false)
Add a term given by the weak form language expression expr which will be assembled in region region a...
Definition: getfem_models.cc:3571
getfem::model::add_initialized_fem_data
void add_initialized_fem_data(const std::string &name, const mesh_fem &mf, const VECT &v)
Add an initialized fixed size data to the model, assumed to be a vector field if the size of the vect...
Definition: getfem_models.h:832
getfem::model::nb_primary_dof
size_type nb_primary_dof() const
Number of primary degrees of freedom in the model.
Definition: getfem_models.h:589
getfem::add_nonlinear_term
size_type APIDECL add_nonlinear_term(model &md, const mesh_im &mim, const std::string &expr, size_type region=size_type(-1), bool is_sym=false, bool is_coercive=false, const std::string &brickname="")
Add a nonlinear term given by the weak form language expression expr which will be assembled in regio...
Definition: getfem_models.cc:3663
getfem::model::nb_internal_dof
size_type nb_internal_dof() const
Number of internal degrees of freedom in the model.
Definition: getfem_models.h:586
getfem::add_Dirichlet_condition_with_Nitsche_method
size_type APIDECL add_Dirichlet_condition_with_Nitsche_method(model &md, const mesh_im &mim, const std::string &varname, const std::string &Neumannterm, const std::string &datagamma0, size_type region, scalar_type theta=scalar_type(0), const std::string &datag=std::string())
Add a Dirichlet condition on the variable varname and the mesh region region.
Definition: getfem_models.cc:5025
getfem::model::enable_brick
void enable_brick(size_type ib)
Enable a brick.
Definition: getfem_models.h:518
getfem::add_generic_elliptic_brick
size_type APIDECL add_generic_elliptic_brick(model &md, const mesh_im &mim, const std::string &varname, const std::string &dataexpr, size_type region=size_type(-1))
Add an elliptic term on the variable varname.
Definition: getfem_models.cc:3931
getfem::add_theta_method_dispatcher
void APIDECL add_theta_method_dispatcher(model &md, dal::bit_vector ibricks, const std::string &THETA)
Add a theta-method time dispatcher to a list of bricks.
Definition: getfem_models.cc:7029
getfem::add_generalized_Dirichlet_condition_with_multipliers
size_type APIDECL add_generalized_Dirichlet_condition_with_multipliers(model &md, const mesh_im &mim, const std::string &varname, const std::string &multname, size_type region, const std::string &dataname, const std::string &Hname)
Add a generalized Dirichlet condition on the variable varname and the mesh region region.
Definition: getfem_models.cc:4770
getfem::add_mass_brick
size_type APIDECL add_mass_brick(model &md, const mesh_im &mim, const std::string &varname, const std::string &dataexpr_rho=std::string(), size_type region=size_type(-1))
Mass brick ( ).
Definition: getfem_models.cc:6511
getfem::add_linear_incompressibility
size_type APIDECL add_linear_incompressibility(model &md, const mesh_im &mim, const std::string &varname, const std::string &multname_pressure, size_type region=size_type(-1), const std::string &dataexpr_penal_coeff=std::string())
Mixed linear incompressibility condition brick.
Definition: getfem_models.cc:6362
getfem::velocity_update_for_order_two_theta_method
void APIDECL velocity_update_for_order_two_theta_method(model &md, const std::string &U, const std::string &V, const std::string &pdt, const std::string &ptheta)
Function which udpate the velocity $v^{n+1}$ after the computation of the displacement $u^{n+1}$ and ...
Definition: getfem_models.cc:7036
getfem::model::add_initialized_fixed_size_data
void add_initialized_fixed_size_data(const std::string &name, const VECT &v, const bgeot::multi_index &sizes)
Add a fixed size data (assumed to be a vector) to the model and initialized with v.
Definition: getfem_models.h:752
getfem::mesh_fem
Describe a finite element method linked to a mesh.
Definition: getfem_mesh_fem.h:148
getfem_partial_mesh_fem.h
a subclass of getfem::mesh_fem which allows to eliminate a number of dof of the original mesh_fem.
getfem::PREFIX_OLD
const auto PREFIX_OLD
A prefix to refer to the previous version of a variable.
Definition: getfem_models.h:98
getfem::model::add_initialized_fem_data
void add_initialized_fem_data(const std::string &name, const mesh_fem &mf, const VECT &v, const bgeot::multi_index &sizes)
Add a fixed size data to the model.
Definition: getfem_models.h:845
getfem::add_Helmholtz_brick
size_type APIDECL add_Helmholtz_brick(model &md, const mesh_im &mim, const std::string &varname, const std::string &dataexpr, size_type region=size_type(-1))
Add a Helmoltz brick to the model.
Definition: getfem_models.cc:5501
getfem::model
`‘Model’' variables store the variables, the data and the description of a model.
Definition: getfem_models.h:114
getfem
GEneric Tool for Finite Element Methods.
Definition: getfem_accumulated_distro.h:46
getfem::model::complex_rhs
const model_complex_plain_vector & complex_rhs() const
Gives access to the right hand side of the tangent linear system.
Definition: getfem_models.h:978
getfem::virtual_brick::complex_post_assembly_in_serial
virtual void complex_post_assembly_in_serial(const model &, size_type, const model::varnamelist &, const model::varnamelist &, const model::mimlist &, model::complex_matlist &, model::complex_veclist &, model::complex_veclist &, size_type, build_version) const
Peform any post assembly action for complex terms.
Definition: getfem_models.h:1576
getfem::virtual_brick::complex_pre_assembly_in_serial
virtual void complex_pre_assembly_in_serial(const model &, size_type, const model::varnamelist &, const model::varnamelist &, const model::mimlist &, model::complex_matlist &, model::complex_veclist &, model::complex_veclist &, size_type, build_version) const
Peform any pre assembly action for complex term assembly.
Definition: getfem_models.h:1548
getfem::model::interpolate_transformation_exists
bool interpolate_transformation_exists(const std::string &name) const
Tests if name corresponds to an interpolate transformation.
Definition: getfem_models.h:1126
getfem::add_normal_source_term_brick
size_type APIDECL add_normal_source_term_brick(model &md, const mesh_im &mim, const std::string &varname, const std::string &dataexpr, size_type region)
Add a source term on the variable varname on a boundary region.
Definition: getfem_models.cc:4264
getfem::add_basic_d_on_dt_brick
size_type APIDECL add_basic_d_on_dt_brick(model &md, const mesh_im &mim, const std::string &varname, const std::string &dataname_dt, const std::string &dataname_rho=std::string(), size_type region=size_type(-1))
Basic d/dt brick ( ).
Definition: getfem_models.cc:6766
getfem::model::add_initialized_fixed_size_data
void add_initialized_fixed_size_data(const std::string &name, const VECT &v)
Add a fixed size data (assumed to be a vector) to the model and initialized with v.
Definition: getfem_models.h:740
getfem_generic_assembly.h
A language for generic assembly of pde boundary value problems.
getfem::compute_isotropic_linearized_Von_Mises_pstress
void APIDECL compute_isotropic_linearized_Von_Mises_pstress(model &md, const std::string &varname, const std::string &data_E, const std::string &data_nu, const mesh_fem &mf_vm, model_real_plain_vector &VM)
Compute the Von-Mises stress of a displacement field for isotropic linearized elasticity in 3D or in ...
Definition: getfem_models.cc:6256
getfem::add_pointwise_constraints_with_given_multipliers
size_type APIDECL add_pointwise_constraints_with_given_multipliers(model &md, const std::string &varname, const std::string &multname, const std::string &dataname_pt, const std::string &dataname_unitv=std::string(), const std::string &dataname_val=std::string())
Add some pointwise constraints on the variable varname using a given multiplier multname.
Definition: getfem_models.cc:5369
getfem::model::touch_brick
void touch_brick(size_type ib)
Force the re-computation of a brick for the next assembly.
Definition: getfem_models.h:1028
getfem::velocity_update_for_Newmark_scheme
void APIDECL velocity_update_for_Newmark_scheme(model &md, size_type id2dt2b, const std::string &U, const std::string &V, const std::string &pdt, const std::string &ptwobeta, const std::string &pgamma)
Function which udpate the velocity $v^{n+1}$ after the computation of the displacement $u^{n+1}$ and ...
Definition: getfem_models.cc:7081
getfem::model::add_interpolate_transformation
void add_interpolate_transformation(const std::string &name, pinterpolate_transformation ptrans)
Add an interpolate transformation to the model to be used with the generic assembly.
Definition: getfem_models.h:1103
getfem::context_dependencies
Deal with interdependencies of objects.
Definition: getfem_context.h:81
getfem_im_data.h
Provides indexing of integration points for mesh_im.
getfem::model::has_internal_variables
bool has_internal_variables() const
Return true if the model has at least one internal variable.
Definition: getfem_models.h:569
getfem::model::complex_tangent_matrix
const model_complex_sparse_matrix & complex_tangent_matrix() const
Gives the access to the tangent matrix.
Definition: getfem_models.h:924
getfem::model::get_active_bricks
const dal::bit_vector & get_active_bricks() const
Return the model brick ids.
Definition: getfem_models.h:1023
getfem::model::update_from_context
void update_from_context() const
this function has to be defined and should update the object when the context is modified.
Definition: getfem_models.h:503
gmm::rsvector
sparse vector built upon std::vector.
Definition: gmm_def.h:488
getfem::model::is_coercive
bool is_coercive() const
Return true if all the model terms do not affect the coercivity of the whole tangent system.
Definition: getfem_models.h:566
getfem::virtual_time_scheme
The time integration scheme object provides the necessary methods for the model object to apply a tim...
Definition: getfem_models.h:1225
getfem::model::add_secondary_domain
void add_secondary_domain(const std::string &name, psecondary_domain ptrans)
Add a secondary domain to the model to be used with the generic assembly.
Definition: getfem_models.h:1157
getfem::model::macro_dictionary
const ga_macro_dictionary & macro_dictionary() const
Dictonnary of user defined macros.
Definition: getfem_models.h:887
getfem::model::set_real_rhs
model_real_plain_vector & set_real_rhs(bool with_internal=false) const
Gives write access to the right hand side of the tangent linear system.
Definition: getfem_models.h:941
getfem::model::real_rhs
const model_real_plain_vector & real_rhs(bool with_internal=false) const
Gives access to the right hand side of the tangent linear system.
Definition: getfem_models.h:932
getfem::model::add_initialized_scalar_data
void add_initialized_scalar_data(const std::string &name, T e)
Add a scalar data (i.e.
Definition: getfem_models.h:779
getfem::model::elementary_transformation_exists
bool elementary_transformation_exists(const std::string &name) const
Tests if name corresponds to an elementary transformation.
Definition: getfem_models.h:1150
getfem::virtual_brick::asm_real_tangent_terms
virtual void asm_real_tangent_terms(const model &, size_type, const model::varnamelist &, const model::varnamelist &, const model::mimlist &, model::real_matlist &, model::real_veclist &, model::real_veclist &, size_type, build_version) const
Assembly of bricks real tangent terms.
Definition: getfem_models.h:1490
getfem::add_Dirichlet_condition_with_penalization
size_type APIDECL add_Dirichlet_condition_with_penalization(model &md, const mesh_im &mim, const std::string &varname, scalar_type penalization_coeff, size_type region, const std::string &dataname=std::string(), const mesh_fem *mf_mult=0)
Add a Dirichlet condition on the variable varname and the mesh region region.
Definition: getfem_models.cc:4695
getfem::model::set_complex_rhs
model_complex_plain_vector & set_complex_rhs() const
Gives write access to the right hand side of the tangent linear system.
Definition: getfem_models.h:987
getfem::add_normal_Dirichlet_condition_with_Nitsche_method
size_type APIDECL add_normal_Dirichlet_condition_with_Nitsche_method(model &md, const mesh_im &mim, const std::string &varname, const std::string &Neumannterm, const std::string &datagamma0, size_type region, scalar_type theta=scalar_type(0), const std::string &datag=std::string())
Add a Dirichlet condition on the normal component of the variable varname and the mesh region region.
Definition: getfem_models.cc:5059
dal::static_stored_object
base class for static stored objects
Definition: dal_static_stored_objects.h:206
getfem::im_data
im_data provides indexing to the integration points of a mesh im object.
Definition: getfem_im_data.h:69
getfem::model::secondary_domain_exists
bool secondary_domain_exists(const std::string &name) const
Tests if name corresponds to an interpolate transformation.
Definition: getfem_models.h:1176
getfem::add_generalized_Dirichlet_condition_with_Nitsche_method
size_type APIDECL add_generalized_Dirichlet_condition_with_Nitsche_method(model &md, const mesh_im &mim, const std::string &varname, const std::string &Neumannterm, const std::string &datagamma0, size_type region, scalar_type theta, const std::string &datag, const std::string &dataH)
Add a Dirichlet condition on the variable varname and the mesh region region.
Definition: getfem_models.cc:5092
getfem::add_twodomain_source_term
size_type APIDECL add_twodomain_source_term(model &md, const mesh_im &mim, const std::string &expr, size_type region, const std::string &secondary_domain, const std::string &brickname=std::string(), const std::string &directvarname=std::string(), const std::string &directdataname=std::string(), bool return_if_nonlin=false)
Adds a source term given by a weak form language expression like add_source_term function but for an ...
Definition: getfem_models.cc:3419
getfem::virtual_brick::real_pre_assembly_in_serial
virtual void real_pre_assembly_in_serial(const model &, size_type, const model::varnamelist &, const model::varnamelist &, const model::mimlist &, model::real_matlist &, model::real_veclist &, model::real_veclist &, size_type, build_version) const
Peform any pre assembly action for real term assembly.
Definition: getfem_models.h:1534
getfem::virtual_brick
The virtual brick has to be derived to describe real model bricks.
Definition: getfem_models.h:1437
getfem::add_normal_Dirichlet_condition_with_penalization
size_type APIDECL add_normal_Dirichlet_condition_with_penalization(model &md, const mesh_im &mim, const std::string &varname, scalar_type penalization_coeff, size_type region, const std::string &dataname=std::string(), const mesh_fem *mf_mult=0)
Add a Dirichlet condition to the normal component of the vector (or tensor) valued variable varname a...
Definition: getfem_models.cc:4749
getfem::model::interpolate_transformation
pinterpolate_transformation interpolate_transformation(const std::string &name) const
Get a pointer to the interpolate transformation name.
Definition: getfem_models.h:1117
getfem::model::internal_solution
const model_real_plain_vector & internal_solution() const
Gives access to the partial solution for condensed internal variables.
Definition: getfem_models.h:949
getfem::add_source_term_brick
size_type APIDECL add_source_term_brick(model &md, const mesh_im &mim, const std::string &varname, const std::string &dataexpr, size_type region=size_type(-1), const std::string &directdataname=std::string())
Add a source term on the variable varname.
Definition: getfem_models.cc:4112
getfem::add_isotropic_linearized_elasticity_brick_pstrain
size_type APIDECL add_isotropic_linearized_elasticity_brick_pstrain(model &md, const mesh_im &mim, const std::string &varname, const std::string &data_E, const std::string &data_nu, size_type region)
Linear elasticity brick ( ).
Definition: getfem_models.cc:6127
getfem::virtual_brick::asm_complex_tangent_terms
virtual void asm_complex_tangent_terms(const model &, size_type, const model::varnamelist &, const model::varnamelist &, const model::mimlist &, model::complex_matlist &, model::complex_veclist &, model::complex_veclist &, size_type, build_version) const
Assembly of bricks complex tangent terms.
Definition: getfem_models.h:1516
getfem::add_isotropic_linearized_elasticity_brick
size_type APIDECL add_isotropic_linearized_elasticity_brick(model &md, const mesh_im &mim, const std::string &varname, const std::string &dataname_lambda, const std::string &dataname_mu, size_type region=size_type(-1), const std::string &dataname_preconstraint=std::string())
Linear elasticity brick ( ).
Definition: getfem_models.cc:6091
getfem::add_pointwise_constraints_with_multipliers
size_type APIDECL add_pointwise_constraints_with_multipliers(model &md, const std::string &varname, const std::string &dataname_pt, const std::string &dataname_unitv=std::string(), const std::string &dataname_val=std::string())
Add some pointwise constraints on the variable varname using multiplier.
Definition: getfem_models.cc:5385
getfem::model::is_symmetric
bool is_symmetric() const
Return true if all the model terms do not affect the coercivity of the whole tangent system.
Definition: getfem_models.h:577
getfem::model::real_tangent_matrix
const model_real_sparse_matrix & real_tangent_matrix(bool internal=false) const
Gives the access to the tangent matrix.
Definition: getfem_models.h:917
getfem::model::secondary_domain
psecondary_domain secondary_domain(const std::string &name) const
Get a pointer to the interpolate transformation name.
Definition: getfem_models.h:1167
getfem::virtual_brick::real_post_assembly_in_serial
virtual void real_post_assembly_in_serial(const model &, size_type, const model::varnamelist &, const model::varnamelist &, const model::mimlist &, model::real_matlist &, model::real_veclist &, model::real_veclist &, size_type, build_version) const
Peform any post assembly action for real terms.
Definition: getfem_models.h:1562
getfem::add_basic_d2_on_dt2_brick
size_type APIDECL add_basic_d2_on_dt2_brick(model &md, const mesh_im &mim, const std::string &varnameU, const std::string &datanameV, const std::string &dataname_dt, const std::string &dataname_alpha, const std::string &dataname_rho=std::string(), size_type region=size_type(-1))
Basic d2/dt2 brick ( ).
Definition: getfem_models.cc:6941
getfem::add_isotropic_linearized_elasticity_brick_pstress
size_type APIDECL add_isotropic_linearized_elasticity_brick_pstress(model &md, const mesh_im &mim, const std::string &varname, const std::string &data_E, const std::string &data_nu, size_type region)
Linear elasticity brick ( ).
Definition: getfem_models.cc:6156
getfem::add_Laplacian_brick
size_type APIDECL add_Laplacian_brick(model &md, const mesh_im &mim, const std::string &varname, size_type region=size_type(-1))
Add a Laplacian term on the variable varname (in fact with a minus : :math:-\text{div}(\nabla u)).
Definition: getfem_models.cc:3906
getfem::pbrick
std::shared_ptr< const virtual_brick > pbrick
type of pointer on a brick
Definition: getfem_models.h:49
getfem::add_Dirichlet_condition_with_simplification
size_type APIDECL add_Dirichlet_condition_with_simplification(model &md, const std::string &varname, size_type region, const std::string &dataname=std::string())
Add a (simple) Dirichlet condition on the variable varname and the mesh region region.
Definition: getfem_models.cc:5008
getfem::add_Fourier_Robin_brick
size_type APIDECL add_Fourier_Robin_brick(model &md, const mesh_im &mim, const std::string &varname, const std::string &dataexpr, size_type region)
Add a Fourier-Robin brick to the model.
Definition: getfem_models.cc:5622
getfem::add_lumped_mass_brick_for_first_order
size_type APIDECL add_lumped_mass_brick_for_first_order(model &md, const mesh_im &mim, const std::string &varname, const std::string &dataexpr_rho=std::string(), size_type region=size_type(-1))
Lumped mass brick for first order.
Definition: getfem_models.cc:6601
getfem::is_old
bool is_old(const std::string &name)
Does the variable have Old_ prefix.
Definition: getfem_models.cc:218
getfem::mesh_fem::nb_dof
virtual size_type nb_dof() const
Return the total number of degrees of freedom.
Definition: getfem_mesh_fem.h:562
getfem_assembling.h
Miscelleanous assembly routines for common terms. Use the low-level generic assembly....
getfem::add_source_term
size_type APIDECL add_source_term(model &md, const mesh_im &mim, const std::string &expr, size_type region=size_type(-1), const std::string &brickname=std::string(), const std::string &directvarname=std::string(), const std::string &directdataname=std::string(), bool return_if_nonlin=false)
Add a source term given by the assembly string expr which will be assembled in region region and with...
Definition: getfem_models.cc:3411
getfem::no_old_prefix_name
std::string no_old_prefix_name(const std::string &name)
Strip the variable name from prefix Old_ if it has one.
Definition: getfem_models.cc:222