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