36 #ifndef GETFEM_ASSEMBLING_TENSORS_H__
37 #define GETFEM_ASSEMBLING_TENSORS_H__
47 #define ASM_THROW_PARSE_ERROR(x) \
48 GMM_ASSERT1(false, "parse error: " << x << endl << "found here:\n " \
49 << syntax_err_print());
50 #define ASM_THROW_TENSOR_ERROR(x) \
51 GMM_ASSERT1(false, "tensor error: " << x);
52 #define ASM_THROW_ERROR(x) GMM_ASSERT1(false, "error: " << x);
55 using bgeot::stride_type;
56 using bgeot::index_type;
57 using bgeot::index_set;
58 using bgeot::tensor_ranges;
59 using bgeot::tensor_strides;
60 using bgeot::tensor_mask;
61 using bgeot::tensor_shape;
62 using bgeot::tensor_ref;
63 using bgeot::multi_tensor_iterator;
73 std::deque< ATN_tensor* > childs_;
78 dim_type current_face;
80 ATN(
const std::string& n=std::string(
"unnamed")) :
81 name_(n), number_(unsigned(-1)), current_cv(
size_type(-1)),
82 current_face(dim_type(-1)) {}
85 void add_child(ATN_tensor& a) { childs_.push_back(&a); }
86 ATN_tensor& child(
size_type n) {
return *childs_[n]; }
87 size_type nchilds() {
return childs_.size(); }
90 void reinit() {
if (!is_zero_size()) reinit_(); }
93 if (cv != current_cv || face != current_face) {
100 const std::string& name() {
return name_; }
101 void set_name(
const std::string& n) { name_ = n; }
105 virtual void update_childs_required_shape();
107 virtual bool is_zero_size();
111 void set_number(
unsigned &gcnt);
112 unsigned number()
const {
return number_; }
114 virtual void reinit_() = 0;
115 virtual void exec_(
size_type , dim_type ) {}
118 class ATN_tensors_sum_scaled;
121 class ATN_tensor :
public ATN {
126 tensor_shape req_shape;
132 ATN_tensor() { shape_updated_ =
false; frozen_ =
false; }
133 bool is_shape_updated()
const {
return shape_updated_; }
134 void freeze() { frozen_ =
true; }
135 bool is_frozen()
const {
return frozen_; }
136 const tensor_ranges& ranges()
const {
return r_; }
137 const tensor_shape& required_shape()
const {
return req_shape; }
143 virtual void check_shape_update(
size_type , dim_type) {}
145 virtual void init_required_shape() { req_shape.set_empty(r_); }
148 virtual void update_childs_required_shape() {
149 for (dim_type i=0; i < nchilds(); ++i) {
150 child(i).merge_required_shape(req_shape);
156 tensor_ref& tensor() {
160 bool is_zero_size() {
return r_.is_zero_size(); }
162 void merge_required_shape(
const tensor_shape& shape_from_parent) {
163 req_shape.merge(shape_from_parent,
false);
168 virtual ATN_tensors_sum_scaled* is_tensors_sum_scaled() {
return 0; }
177 bool is_mf_ref()
const {
return (pmf != 0); }
178 vdim_specif() { dim =
size_type(-1); pmf = 0; }
179 vdim_specif(
size_type i) { dim = i; pmf = 0; }
180 vdim_specif(
const mesh_fem *pmf_) { dim = pmf_->nb_dof(); pmf = pmf_; }
182 class vdim_specif_list :
public std::vector< vdim_specif > {
184 vdim_specif_list() { reserve(8); }
187 void build_strides_for_cv(
size_type cv, tensor_ranges& r,
188 std::vector<tensor_strides >& str)
const;
194 template<
typename VEC >
class ATN_array_output :
public ATN {
196 vdim_specif_list vdim;
197 multi_tensor_iterator mti;
198 tensor_strides strides;
201 ATN_array_output(ATN_tensor& a, VEC& v_, vdim_specif_list &d)
204 strides.resize(vdim.size()+1);
208 for (
size_type i=0; i < vdim.size(); ++i) {
209 if (vdim[i].pmf) pmf = vdim[i].pmf;
210 strides[i+1] = strides[i]*int(vdim[i].dim);
212 if (gmm::vect_size(v) !=
size_type(strides[vdim.size()]))
213 ASM_THROW_TENSOR_ERROR(
"wrong size for output vector: supplied "
214 "vector size is " << gmm::vect_size(v)
215 <<
" while it should be "
216 << strides[vdim.size()]);
220 mti = multi_tensor_iterator(child(0).tensor(),
true);
225 std::vector< tensor_strides > str;
226 vdim.build_strides_for_cv(cv, r, str);
227 if (child(0).ranges() != r) {
228 ASM_THROW_TENSOR_ERROR(
"can't output a tensor of dimensions "
229 << child(0).ranges() <<
230 " into an output array of size " << r);
233 if (pmf && pmf->is_reduced()) {
234 if ( pmf->nb_dof() != 0)
238 dim_type qqdim = dim_type(gmm::vect_size(v) / nb_dof);
242 for (dim_type j=0; j < mti.ndim(); ++j) i += str[j][mti.index(j)];
243 gmm::add(gmm::scaled(gmm::mat_row(pmf->extension_matrix(), i),
247 GMM_ASSERT1(
false,
"To be verified ... ");
249 for (dim_type j=0; j < mti.ndim(); ++j) i += str[j][mti.index(j)];
250 gmm::add(gmm::scaled(gmm::mat_row(pmf->extension_matrix(),i/qqdim),
252 gmm::sub_vector(v, gmm::sub_slice(i%qqdim,nb_dof,qqdim)));
254 }
while (mti.qnext1());
259 typename gmm::linalg_traits<VEC>::iterator it = gmm::vect_begin(v);
260 for (dim_type j = 0; j < mti.ndim(); ++j) it += str[j][mti.index(j)];
262 }
while (mti.qnext1());
267 template <
typename MAT,
typename ROW,
typename COL>
268 void asmrankoneupdate(
const MAT &m_,
const ROW &row,
const COL &col,
270 MAT &m =
const_cast<MAT &
>(m_);
271 typename gmm::linalg_traits<ROW>::const_iterator itr = row.begin();
272 for (; itr != row.end(); ++itr) {
273 typename gmm::linalg_traits<COL>::const_iterator itc = col.begin();
274 for (; itc != col.end(); ++itc)
275 m(itr.index(), itc.index()) += (*itr) * (*itc) * r;
279 template <
typename MAT,
typename ROW>
280 void asmrankoneupdate(
const MAT &m_,
const ROW &row,
size_type j, scalar_type r) {
281 MAT &m =
const_cast<MAT &
>(m_);
282 typename gmm::linalg_traits<ROW>::const_iterator itr = row.begin();
283 for (; itr != row.end(); ++itr) m(itr.index(), j) += (*itr) * r;
286 template <
typename MAT,
typename COL>
287 void asmrankoneupdate(
const MAT &m_,
size_type j,
const COL &col, scalar_type r) {
288 MAT &m =
const_cast<MAT &
>(m_);
289 typename gmm::linalg_traits<COL>::const_iterator itc = col.begin();
290 for (; itc != col.end(); ++itc) m(j, itc.index()) += (*itc) * r;
294 template<
typename MAT >
class ATN_smatrix_output :
public ATN {
295 const mesh_fem &mf_r, &mf_c;
297 multi_tensor_iterator mti;
305 ATN_smatrix_output(ATN_tensor& a,
const mesh_fem& mf_r_,
306 const mesh_fem& mf_c_, MAT& m_)
307 : mf_r(mf_r_), mf_c(mf_c_), m(m_) {
313 mti = multi_tensor_iterator(child(0).tensor(),
true);
317 size_type nb_r = mf_r.nb_basic_dof_of_element(cv);
318 size_type nb_c = mf_c.nb_basic_dof_of_element(cv);
319 if (child(0).tensor().ndim() != 2)
320 ASM_THROW_TENSOR_ERROR(
"cannot write a " <<
321 int(child(0).tensor().ndim()) <<
322 "D-tensor into a matrix!");
323 if (child(0).tensor().dim(0) != nb_r ||
324 child(0).tensor().dim(1) != nb_c) {
325 ASM_THROW_TENSOR_ERROR(
"size mismatch for sparse matrix output:"
326 " tensor dimension is " << child(0).ranges()
327 <<
", while the elementary matrix for convex "
328 << cv <<
" should have " << nb_r <<
"x"
329 << nb_c <<
" elements");
331 std::vector<size_type> cvdof_r(mf_r.ind_basic_dof_of_element(cv).begin(),
332 mf_r.ind_basic_dof_of_element(cv).end());
333 std::vector<size_type> cvdof_c(mf_c.ind_basic_dof_of_element(cv).begin(),
334 mf_c.ind_basic_dof_of_element(cv).end());
336 if (it.size() == 0) {
344 }
while (mti.qnext1());
347 bool valid_mf_r = mf_r.nb_dof() > 0;
348 bool valid_mf_c = mf_c.nb_dof() > 0;
350 if (mf_r.is_reduced()) {
351 if (mf_c.is_reduced() && valid_mf_r && valid_mf_c) {
352 for (
unsigned i = 0; i < it.size(); ++i)
354 asmrankoneupdate(m, gmm::mat_row(mf_r.extension_matrix(),
356 gmm::mat_row(mf_c.extension_matrix(),
360 else if (valid_mf_r) {
361 for (
unsigned i = 0; i < it.size(); ++i)
363 asmrankoneupdate(m, gmm::mat_row(mf_r.extension_matrix(),
365 cvdof_c[it[i].j], *it[i].p);
369 if (mf_c.is_reduced() && valid_mf_c) {
370 for (
unsigned i = 0; i < it.size(); ++i)
372 asmrankoneupdate(m, cvdof_r[it[i].i],
373 gmm::mat_row(mf_c.extension_matrix(),
378 for (
unsigned i = 0; i < it.size(); ++i)
380 m(cvdof_r[it[i].i], cvdof_c[it[i].j]) += *it[i].p;
392 class base_asm_data {
395 virtual void copy_with_mti(
const std::vector<tensor_strides> &,
396 multi_tensor_iterator &,
397 const mesh_fem *)
const = 0;
398 virtual ~base_asm_data() {}
401 template<
typename VEC >
class asm_data :
public base_asm_data {
404 asm_data(
const VEC *v_) : v(*v_) {}
406 return gmm::vect_size(v);
410 void copy_with_mti(
const std::vector<tensor_strides> &str,
411 multi_tensor_iterator &mti,
const mesh_fem *pmf)
const {
413 if (pmf && pmf->is_reduced()) {
416 for (dim_type i = 0; i < mti.ndim(); ++i) ppos+=str[i][mti.index(i)];
418 =
gmm::vect_sp(gmm::mat_row(pmf->extension_matrix(), ppos), v);
419 }
while (mti.qnext1());
425 for (dim_type i = 0; i < mti.ndim(); ++i) ppos+=str[i][mti.index(i)];
427 }
while (mti.qnext1());
434 virtual std::unique_ptr<ATN> build_output_tensor(ATN_tensor &a,
435 vdim_specif_list& vdim)=0;
436 virtual ~base_asm_vec() {}
439 template<
typename VEC >
class asm_vec :
public base_asm_vec {
440 std::shared_ptr<VEC> v;
442 asm_vec(
const std::shared_ptr<VEC> &v_) : v(v_) {}
443 asm_vec(VEC *v_) : v(std::shared_ptr<VEC>(), v_) {}
444 virtual std::unique_ptr<ATN> build_output_tensor(ATN_tensor &a,
445 vdim_specif_list& vdim) {
446 return std::make_unique<ATN_array_output<VEC>>(a, *v, vdim);
448 VEC *vec() {
return v.get(); }
454 class base_vec_factory {
456 virtual base_asm_vec* create_vec(
const tensor_ranges& r) = 0;
457 virtual ~base_vec_factory() {}
460 template<
typename VEC >
class vec_factory
461 :
public base_vec_factory,
private std::deque<asm_vec<VEC> > {
463 base_asm_vec* create_vec(
const tensor_ranges& r) {
466 ASM_THROW_TENSOR_ERROR(
"can't create a vector of size " << r);
467 this->push_back(asm_vec<VEC>(std::make_shared<VEC>(sz)));
468 return &this->back();
476 virtual std::unique_ptr<ATN>
477 build_output_tensor(ATN_tensor& a,
const mesh_fem& mf1,
478 const mesh_fem& mf2) = 0;
479 virtual ~base_asm_mat() {}
482 template<
typename MAT >
class asm_mat :
public base_asm_mat {
483 std::shared_ptr<MAT> m;
485 asm_mat(
const std::shared_ptr<MAT> &m_) : m(m_) {}
486 asm_mat(MAT *m_) : m(std::shared_ptr<MAT>(), m_) {}
488 build_output_tensor(ATN_tensor& a,
const mesh_fem& mf1,
489 const mesh_fem& mf2) {
490 return std::make_unique<ATN_smatrix_output<MAT>>(a, mf1, mf2, *m);
492 MAT *mat() {
return m.get(); }
496 class base_mat_factory {
499 virtual ~base_mat_factory() {};
502 template<
typename MAT >
class mat_factory
503 :
public base_mat_factory,
private std::deque<asm_mat<MAT> > {
506 this->push_back(asm_mat<MAT>(std::make_shared<MAT>(m, n)));
507 return &this->back();
514 typedef enum { TNCONST, TNTENSOR, TNNONE } node_type;
520 tnode() : type_(TNNONE), x(1e300), t(NULL) {}
521 tnode(scalar_type x_) { assign(x_); }
522 tnode(ATN_tensor *t_) { assign(t_); }
523 void assign(scalar_type x_) { type_ = TNCONST; t = NULL; x = x_; }
524 void assign(ATN_tensor *t_) { type_ = TNTENSOR; t = t_; x = 1e300; }
525 ATN_tensor* tensor() { assert(type_ == TNTENSOR);
return t; }
526 scalar_type xval() { assert(type_ == TNCONST);
return x; }
527 node_type type() {
return type_; }
528 void check0() {
if (xval() == 0) ASM_THROW_ERROR(
"division by zero"); }
531 class asm_tokenizer {
533 typedef enum { OPEN_PAR=
'(', CLOSE_PAR=
')', COMMA=
',',
534 SEMICOLON=
';', COLON=
':', EQUAL=
'=', MFREF=
'#', IMREF=
'%',
535 PLUS=
'+',MINUS=
'-', PRODUCT=
'.',MULTIPLY=
'*',
536 DIVIDE=
'/', ARGNUM_SELECTOR=
'$',
537 OPEN_BRACE=
'{', CLOSE_BRACE=
'}',
538 END=0, IDENT=1, NUMBER=2 } tok_type_enum;
542 tok_type_enum curr_tok_type;
543 std::string curr_tok;
545 double curr_tok_dval;
547 std::deque<size_type> marks;
550 void set_str(
const std::string& s_) {
551 str = s_; tok_pos = 0; tok_len =
size_type(-1); curr_tok_type = END;
552 err_msg_mark = 0; get_tok();
554 std::string tok()
const {
return curr_tok; }
555 tok_type_enum tok_type()
const {
return curr_tok_type; }
558 {
return str.substr(i1, i2-i1); }
559 void err_set_mark() {
560 err_msg_mark = tok_pos;
562 void push_mark() { marks.push_back(tok_pos); }
563 void pop_mark() { assert(marks.size()); marks.pop_back(); }
564 std::string mark_txt() {
565 assert(marks.size());
566 return tok_substr(marks.back(),tok_pos);
570 std::string syntax_err_print();
571 void accept(tok_type_enum t,
const char *msg_=
"syntax error") {
572 if (tok_type() != t) ASM_THROW_PARSE_ERROR(msg_); advance();
574 void accept(tok_type_enum t, tok_type_enum t2,
575 const char *msg_=
"syntax error") {
576 if (tok_type() != t && tok_type() != t2)
577 ASM_THROW_PARSE_ERROR(msg_);
580 bool advance_if(tok_type_enum t) {
581 if (tok_type() == t) { advance();
return true; }
else return false;
583 void advance() { tok_pos += tok_len; get_tok(); }
585 double tok_number_dval()
586 { assert(tok_type()==NUMBER);
return curr_tok_dval; }
587 int tok_number_ival(
int maxval=10000000) {
588 int n=int(tok_number_dval());
589 if (n != curr_tok_dval) ASM_THROW_PARSE_ERROR(
"not an integer");
590 if (n > maxval) ASM_THROW_PARSE_ERROR(
"out of bound integer");
594 { assert(tok_type()==MFREF);
return curr_tok_ival; }
596 { assert(tok_type()==IMREF);
return curr_tok_ival; }
598 { assert(tok_type()==ARGNUM_SELECTOR);
return curr_tok_ival; }
607 std::vector<const mesh_fem *> mftab;
608 std::vector<const mesh_im *> imtab;
609 std::vector<pnonlinear_elem_term> innonlin;
611 std::vector<std::unique_ptr<base_asm_data>> indata;
612 std::vector<std::shared_ptr<base_asm_vec>> outvec;
614 std::vector<std::shared_ptr<base_asm_mat>> outmat;
617 base_vec_factory *vec_fact;
619 base_mat_factory *mat_fact;
622 std::vector<std::unique_ptr<ATN>> outvars;
625 std::map<std::string, ATN_tensor *> vars;
626 std::vector<std::unique_ptr<ATN_tensor>> atn_tensors;
637 vec_fact(0), mat_fact(0), parse_done(
false)
641 void set(
const std::string& s_) { set_str(s_); }
642 const std::vector<const mesh_fem*>& mf()
const {
return mftab; }
643 const std::vector<const mesh_im*>& im()
const {
return imtab; }
644 const std::vector<pnonlinear_elem_term> nonlin()
const {
return innonlin; }
645 const std::vector<std::unique_ptr<base_asm_data>>& data()
const {
return indata; }
646 const std::vector<std::shared_ptr<base_asm_vec>>& vec()
const {
return outvec; }
647 const std::vector<std::shared_ptr<base_asm_mat>>& mat()
const {
return outmat; }
652 void push_im(
const mesh_im& im_) { imtab.push_back(&im_); }
655 innonlin.push_back(net);
659 indata.push_back(std::make_unique<asm_data<VEC>>(&d));
663 outvec.push_back(std::make_shared<asm_vec<VEC>>(&(gmm::linalg_cast(v))));
666 template<
typename VEC >
void push_vec(
const VEC& v) {
667 outvec.push_back(std::make_shared<asm_vec<VEC>>(&(gmm::linalg_cast(v))));
670 template<
typename MAT >
void push_mat(
const MAT& m) {
671 outmat.push_back(std::make_shared<asm_mat<MAT>>(&(gmm::linalg_cast(m))));
675 outmat.push_back(std::make_shared<asm_mat<MAT>>(&(gmm::linalg_cast(m))));
678 template <
typename T>
void push_mat_or_vec(T &v) {
679 push_mat_or_vec(v,
typename gmm::linalg_traits<T>::linalg_type());
684 void set_mat_factory(base_mat_factory *fact) { mat_fact = fact; }
687 ATN_tensor* record(std::unique_ptr<ATN_tensor> &&t) {
688 t->set_name(mark_txt());
689 atn_tensors.push_back(std::move(t));
return atn_tensors.back().get();
691 ATN* record_out(std::unique_ptr<ATN> t) {
692 t->set_name(mark_txt());
693 outvars.push_back(std::move(t));
return outvars.back().get();
695 const mesh_fem& do_mf_arg_basic();
696 const mesh_fem& do_mf_arg(std::vector<const mesh_fem*> *multimf = 0);
697 void do_dim_spec(vdim_specif_list& lst);
698 std::string do_comp_red_ops();
699 ATN_tensor* do_comp();
700 ATN_tensor* do_data();
701 std::pair<ATN_tensor*, std::string> do_red_ops(ATN_tensor* t);
704 tnode do_prod_trans();
709 void consistency_check();
710 template <
typename T>
void push_mat_or_vec(T &v, gmm::abstract_vector) {
713 template <
typename T>
void push_mat_or_vec(T &v, gmm::abstract_matrix) {
721 void assembly(
const mesh_region ®ion =
Sparse tensors, used during the assembly.
Generic assembly of vectors, matrices.
void push_vec(const VEC &v)
Add a new output vector (fake const version..)
void push_mat(MAT &m)
Add a new output matrix.
void push_vec(VEC &v)
Add a new output vector.
void push_data(const VEC &d)
Add a new data (dense array)
void set_vec_factory(base_vec_factory *fact)
used by the getfem_interface..
void push_mat(const MAT &m)
Add a new output matrix (fake const version..)
void assembly(const mesh_region ®ion=mesh_region::all_convexes())
do the assembly on the specified region (boundary or set of convexes)
void push_nonlinear_term(pnonlinear_elem_term net)
Add a new non-linear term.
void push_mi(const mesh_im &im_)
Add a new mesh_im.
void push_mf(const mesh_fem &mf_)
Add a new mesh_fem.
Describe a finite element method linked to a mesh.
Describe an integration method linked to a mesh.
static mesh_region all_convexes()
provide a default value for the mesh_region parameters of assembly procedures etc.
abstract class for integration of non-linear terms into the mat_elem computations the nonlinear term ...
elementary computations (used by the generic assembly).
Build elementary tensors descriptors, used by generic assembly.
Define the getfem::mesh_fem class.
Define the getfem::mesh_im class (integration of getfem::mesh_fem).
strongest_value_type< V1, V2 >::value_type vect_sp(const V1 &v1, const V2 &v2)
*/
void add(const L1 &l1, L2 &l2)
*/
Include the base gmm files.
size_t size_type
used as the common size type in the library
GEneric Tool for Finite Element Methods.