23 #include "getfem/getfem_arch_config.h" 28 size_type vdim_specif_list::nb_mf()
const {
29 return std::count_if(begin(),end(),
30 std::mem_fun_ref(&vdim_specif::is_mf_ref));
32 size_type vdim_specif_list::nbelt()
const {
34 for (const_iterator it = begin(); it != end(); ++it) sz *= (*it).dim;
37 void vdim_specif_list::build_strides_for_cv
38 (
size_type cv, tensor_ranges& r, std::vector<tensor_strides >& str)
const {
39 stride_type s = 1, cnt = 0;
42 for (const_iterator it = begin(); it != end(); ++it, ++cnt) {
43 if ((*it).is_mf_ref()) {
44 r[cnt] = unsigned((*it).pmf->nb_basic_dof_of_element(cv));
47 str[cnt].resize(r[cnt]);
48 for (index_type j=0; j < r[cnt]; ++j) {
49 str[cnt][j] = int((*it).pmf->ind_basic_dof_of_element(cv)[j]*s);
52 r[cnt] = int((*it).dim);
53 str[cnt].resize(r[cnt]);
54 for (index_type j=0; j < (*it).dim; ++j) {
58 s *= stride_type((*it).dim);
62 void ATN::update_childs_required_shape() {
63 for (dim_type i=0; i < nchilds(); ++i) {
64 child(i).merge_required_shape(tensor_shape(child(i).ranges()));
67 void ATN::set_number(
unsigned &gcnt) {
68 if (number_ ==
unsigned(-1)) {
69 for (
unsigned i=0; i < nchilds(); ++i) child(i).set_number(gcnt);
74 bool ATN::is_zero_size() {
75 return child(0).is_zero_size();
81 class ATN_tensor_w_data :
public ATN_tensor {
84 std::vector<scalar_type> data;
87 { ATN_tensor_w_data::reinit_(); std::fill(data.begin(), data.end(),0); }
91 void ATN_tensor_w_data::reinit_() {
92 tr.assign_shape(req_shape);
94 if (tr.card() > 10000000) {
95 cerr <<
"warning, a tensor of size " << tr.card()
96 <<
" will be created, it needs " 97 << tr.card()*
sizeof(scalar_type) <<
" bytes of memory\n";
100 cerr <<
"WARNING: tensor " << name()
101 <<
" will be created with a size of " 102 << ranges() <<
" and 0 non-null elements!" << endl;
104 data.resize(tr.card());
105 data_base = &data[0];
106 tr.set_base(data_base);
116 typedef std::vector<std::pair<ATN_tensor*,std::string> >
117 reduced_tensor_arg_type;
119 class ATN_reduced_tensor :
public ATN_tensor_w_data {
121 reduced_tensor_arg_type red;
122 bgeot::tensor_reduction tred;
124 void check_shape_update(
size_type , dim_type) {
125 shape_updated_ =
false;
126 for (dim_type i=0; i < nchilds(); ++i) {
127 if (child(i).is_shape_updated()) {
128 shape_updated_ =
true;
131 if (shape_updated_) {
133 for (dim_type i=0; i < nchilds(); ++i) {
134 std::string s = red_n(i);
135 if (s.size() != child(i).ranges().size()) {
136 ASM_THROW_TENSOR_ERROR(
"wrong number of indexes for the " 138 <<
"th argument of the reduction " 140 <<
" (ranges=" << child(i).ranges() <<
")");
142 for (
size_type j=0; j < s.length(); ++j) {
143 if (s[j] ==
' ') r_.push_back(child(i).ranges()[j]);
148 void update_childs_required_shape() {
151 for (dim_type n=0; n < nchilds(); ++n) {
152 tensor_shape ts(child(n).ranges());
153 tensor_ranges rn(child(n).ranges());
154 const std::string &s = red[n].second;
155 GMM_ASSERT1(rn.size() == s.size(),
"Wrong size !");
156 for (
unsigned i=0; i < rn.size(); ++i)
159 if (p !=
size_type(-1) && p < i && rn[p] != rn[i])
160 ASM_THROW_TENSOR_ERROR(
"can't reduce the diagonal of a tensor " 161 "of size " << rn <<
" with '" 166 bgeot::tensor_reduction::diag_shape(ts, red[n].second);
171 child(n).merge_required_shape(ts);
184 ATN_reduced_tensor(reduced_tensor_arg_type& r) : red(r) {
185 for (
size_type i=0; i < r.size(); ++i) add_child(*red[i].first);
189 std::string s = red[n].second;
191 s.append(red[n].first->ranges().size(),
' ');
198 for (dim_type i=0; i < red.size(); ++i) {
209 tred.insert(red[i].first->tensor(), red_n(i));
215 ATN_tensor_w_data::reinit0();
217 tred.prepare(&tensor());
221 std::fill(data.begin(), data.end(), 0.);
233 class ATN_sliced_tensor :
public ATN_tensor {
237 ATN_sliced_tensor(ATN_tensor& a, dim_type slice_dim_,
239 slice_dim(slice_dim_), slice_idx(slice_idx_) { add_child(a); }
240 void check_shape_update(
size_type , dim_type) {
241 if ((shape_updated_ = child(0).is_shape_updated())) {
242 if (slice_dim >= child(0).ranges().size() ||
243 slice_idx >= child(0).ranges()[slice_dim]) {
244 ASM_THROW_TENSOR_ERROR(
"can't slice tensor " << child(0).ranges()
245 <<
" at index " <<
int(slice_idx)
246 <<
" of dimension " <<
int(slice_dim));
248 r_ = child(0).ranges(); r_.erase(r_.begin()+slice_dim);
251 void update_childs_required_shape() {
252 tensor_shape ts = req_shape;
253 ts.set_ndim_noclean(dim_type(ts.ndim()+1));
254 ts.shift_dim_num_ge(slice_dim,+1);
255 ts.push_mask(tensor_mask(child(0).ranges()[slice_dim],
256 tensor_mask::Slice(slice_dim, index_type(slice_idx))));
257 child(0).merge_required_shape(ts);
261 tensor() = tensor_ref(child(0).tensor(),
262 tensor_mask::Slice(slice_dim, index_type(slice_idx)));
271 class ATN_permuted_tensor :
public ATN_tensor {
272 std::vector<dim_type> reorder;
275 ATN_permuted_tensor(ATN_tensor& a,
const std::vector<dim_type>& reorder_) :
276 reorder(reorder_) { add_child(a); }
278 void check_shape_update(
size_type , dim_type) {
279 if ((shape_updated_ = child(0).is_shape_updated())) {
280 if (reorder.size() != child(0).ranges().size())
281 ASM_THROW_TENSOR_ERROR(
"can't reorder tensor '" << name()
282 <<
"' of dimensions " << child(0).ranges()
283 <<
" with this permutation: " << vref(reorder));
284 r_.resize(reorder.size());
285 std::fill(r_.begin(), r_.end(), dim_type(-1));
290 for (
size_type i=0; i < reorder.size(); ++i)
291 r_[i] = child(0).ranges()[reorder[i]];
294 void update_childs_required_shape() {
295 tensor_shape ts = req_shape;
296 ts.permute(reorder,
true);
297 child(0).merge_required_shape(ts);
300 tensor() = child(0).tensor();
301 tensor().permute(reorder);
311 class ATN_diagonal_tensor :
public ATN_tensor {
314 ATN_diagonal_tensor(ATN_tensor& a, dim_type i1_, dim_type i2_) :
315 i1(i1_), i2(i2_) { add_child(a); }
317 void check_shape_update(
size_type , dim_type) {
318 if ((shape_updated_ = child(0).is_shape_updated())) {
319 if (i1 >= child(0).ranges().size() || i2 >= child(0).ranges().size() ||
320 i1 == i2 || child(0).ranges()[i1] != child(0).ranges()[i2])
321 ASM_THROW_TENSOR_ERROR(
"can't take the diagonal of a tensor of " 322 "sizes " << child(0).ranges() <<
323 " at indexes " <<
int(i1) <<
" and " 325 r_ = child(0).ranges();
328 void update_childs_required_shape() {
329 tensor_shape ts = req_shape.diag_shape(tensor_mask::Diagonal(i1,i2));
330 child(0).merge_required_shape(ts);
333 tensor() = tensor_ref(child(0).tensor(), tensor_mask::Diagonal(i1,i2));
340 struct computed_tensor_integration_callback
341 :
public mat_elem_integration_callback {
342 bgeot::tensor_reduction red;
344 std::vector<TDIter> tensor_bases;
346 virtual void exec(bgeot::base_tensor &t,
bool first, scalar_type c) {
349 std::fill(t.begin(), t.end(), 0.);
353 for (
unsigned k=0; k!=eltm.size(); ++k) {
354 tensor_bases[k] =
const_cast<TDIter
>(&(*eltm[k]->begin()));
357 long one = 1, n = int(red.out_data.size()); assert(n);
358 gmm::daxpy_(&n, &c, const_cast<double*>(&(red.out_data[0])),
359 &one, (
double*)&(t[0]), &one);
361 void resize_t(bgeot::base_tensor &t) {
362 bgeot::multi_index r;
363 if (red.reduced_range.size())
364 r.assign(red.reduced_range.begin(), red.reduced_range.end());
365 else { r.resize(1); r[0]=1; }
384 pnonlinear_elem_term nlt;
389 std::vector<const mesh_fem*> auxmf;
390 typedef enum { BASE=1, GRAD=2, HESS=3, NORMAL=4, GRADGT=5, GRADGTINV=6,
391 NONLIN=7, DATA=8 } op_type;
392 typedef enum { PLAIN_SHAPE = 0, VECTORIZED_SHAPE = 1,
393 MATRIXIZED_SHAPE = 2 } field_shape_type;
396 field_shape_type vshape;
402 std::string reduction;
414 mf_comp(mf_comp_vect *ow,
const mesh_fem* pmf_, op_type op_,
415 field_shape_type fst) :
416 nlt(0), pmf(pmf_), owner(ow), data(0), op(op_), vshape(fst) { }
417 mf_comp(mf_comp_vect *ow,
const std::vector<const mesh_fem*> vmf,
418 pnonlinear_elem_term nlt_) :
419 nlt(nlt_), pmf(vmf[0]), owner(ow), data(0),
420 auxmf(vmf.begin()+1, vmf.end()), op(NONLIN),
421 vshape(PLAIN_SHAPE) { }
422 mf_comp(mf_comp_vect *ow, ATN_tensor *t) :
423 nlt(0), pmf(0), owner(ow), data(t), op(DATA), vshape(PLAIN_SHAPE) {}
424 void push_back_dimensions(
size_type cv, tensor_ranges &rng,
425 bool only_reduced=
false)
const;
426 bool reduced(
unsigned i)
const {
427 if (i >= reduction.size())
return false;
428 else return reduction[i] !=
' ';
432 struct mf_comp_vect :
public std::vector<mf_comp> {
433 const mesh_im *main_im;
435 mf_comp_vect() :
std::vector<mf_comp>(), main_im(0) { }
436 mf_comp_vect(
const mf_comp_vect &other)
437 :
std::vector<mf_comp>(other), main_im(other.main_im) {
438 for (
size_type i=0; i < size(); ++i) (*
this)[i].owner =
this;
440 void set_im(
const mesh_im &mim) { main_im = &mim; }
441 const mesh_im& get_im()
const {
return *main_im; }
443 mf_comp_vect& operator=(
const mf_comp_vect &other);
446 void mf_comp::push_back_dimensions(
size_type cv, tensor_ranges &rng,
447 bool only_reduced)
const {
451 const bgeot::multi_index &sizes = nlt->sizes(cv);
452 for (
unsigned j=0; j < sizes.size(); ++j)
453 if (!only_reduced || !reduced(j))
458 for (
unsigned i=0; i < data->ranges().size(); ++i)
459 if (!only_reduced || !reduced(i))
460 rng.push_back(data->ranges()[i]);
464 assert(&owner->get_im());
465 assert(owner->get_im().linked_mesh().dim() != dim_type(-1));
466 rng.push_back(owner->get_im().linked_mesh().dim());
471 assert(&owner->get_im());
472 rng.push_back(owner->get_im().linked_mesh().dim());
473 rng.push_back(owner->get_im().linked_mesh().structure_of_convex(cv)->dim());
477 if (!only_reduced || !reduced(d))
478 rng.push_back(
unsigned(pmf->nb_basic_dof_of_element(cv)));
480 if (vshape == mf_comp::VECTORIZED_SHAPE) {
481 if (!only_reduced || !reduced(d)) rng.push_back(pmf->get_qdim());
484 if (vshape == mf_comp::MATRIXIZED_SHAPE) {
485 if (!only_reduced || !reduced(d)) {
486 GMM_ASSERT1(pmf->get_qdims().size() == 2,
"Non matrix field");
487 rng.push_back(dim_type(pmf->get_qdims()[0]));
490 if (!only_reduced || !reduced(d)) rng.push_back(dim_type(pmf->get_qdims()[1]));
494 if (op == GRAD || op == HESS) {
495 if (!only_reduced || !reduced(d))
496 rng.push_back(pmf->linked_mesh().dim());
500 if (!only_reduced || !reduced(d))
501 rng.push_back(pmf->linked_mesh().dim());
509 class ATN_computed_tensor :
public ATN_tensor {
512 pmat_elem_computation pmec;
514 pintegration_method pim;
517 std::vector<scalar_type> data;
520 dal::bit_vector req_bv;
523 bool has_inline_reduction;
525 computed_tensor_integration_callback icb;
530 bgeot::tensor_reduction fallback_red;
531 bool fallback_red_uptodate;
532 TDIter fallback_base;
537 ATN_computed_tensor(
const mf_comp_vect &mfcomp_) :
538 mfcomp(mfcomp_), pmec(0), pme(0), pim(0), pgt(0), data_base(0) {
539 has_inline_reduction =
false;
540 bool in_data =
false;
541 for (
size_type i=0; i < mfcomp.size(); ++i) {
542 if (mfcomp[i].reduction.size() || mfcomp[i].op == mf_comp::DATA) {
543 has_inline_reduction =
true;
544 if (mfcomp[i].op == mf_comp::DATA) { add_child(*mfcomp[i].data); in_data =
true; }
546 if (mfcomp[i].op != mf_comp::DATA && in_data) {
548 ASM_THROW_ERROR(
"data tensors inside comp() cannot be intermixed with Grad() and Base() etc., they must appear LAST");
560 stride_type add_dim(
const tensor_ranges& rng, dim_type d, stride_type s, tensor_ref &tref) {
561 assert(d < rng.size());
563 index_type r = rng[d];
564 tensor_mask m; m.set_full(d, r);
566 for (index_type i=0; i < r; ++i) v[i] = s*i;
567 assert(tref.masks().size() == tref.strides().size());
568 tref.set_ndim_noclean(dim_type(tref.ndim()+1));
570 tref.strides().push_back(v);
577 stride_type add_vdim(
const tensor_ranges& rng, dim_type d,
578 index_type target_dim, stride_type s,
580 assert(d < rng.size()-1);
581 index_type r = rng[d], q=rng[d+1];
582 index_type qmult = q/target_dim;
583 assert(r%qmult == 0); assert(q%qmult==0);
586 tensor_ranges trng(2); trng[0] = q; trng[1] = r;
587 index_set ti(2); ti[0] = dim_type(d+1); ti[1] = d;
588 tensor_mask m(trng,ti);
589 v.resize(r*target_dim);
590 tensor_ranges cnt(2);
591 for (cnt[1]=0; cnt[1] < r; cnt[1]++) {
592 for (index_type k=0; k < target_dim; ++k) {
593 cnt[0] = k*qmult + (cnt[1]%qmult);
594 m.set_mask_val(m.lpos(cnt),
true);
595 v[cnt[1]*target_dim+k] = s*( k * r/qmult + cnt[1]/qmult);
598 assert(tref.masks().size() == tref.strides().size());
599 tref.set_ndim_noclean(dim_type(tref.ndim()+2));
601 tref.strides().push_back(v);
602 return s*(r/qmult)*target_dim;
628 stride_type add_mdim(
const tensor_ranges& rng, dim_type d,
629 index_type target_dim, stride_type s, tensor_ref &tref) {
630 assert(d < rng.size()-2);
633 index_type r = rng[d], q=rng[d+1], p=rng[d+2];
634 index_type qmult = (q*p)/target_dim;
637 assert(p % target_dim == 0);
638 assert(r % (p/target_dim) == 0);
641 tensor_ranges trng(3); trng[0] = q; trng[1] = p; trng[2] = r;
642 index_set ti(3); ti[0] = dim_type(d+1); ti[1] = dim_type(d+2); ti[2] = d;
643 tensor_mask m(trng,ti);
644 v.resize(r*target_dim);
645 tensor_ranges cnt(3);
646 for (cnt[2]=0; cnt[2] < r; cnt[2]++) {
647 for (index_type k=0; k < target_dim; ++k) {
648 unsigned pos = (cnt[2]*target_dim+k) % (q*p);
650 unsigned ii = (pos/p), jj = (pos%p);
651 cnt[0] = ii; cnt[1] = jj;
653 m.set_mask_val(m.lpos(cnt),
true);
654 v[cnt[2]*target_dim+k] = s*(k * r/qmult + cnt[2]/qmult);
657 assert(tref.masks().size() == tref.strides().size());
658 tref.set_ndim_noclean(dim_type(tref.ndim()+3));
663 tref.strides().push_back(v);
664 return s*(r/qmult)*target_dim;
671 for (
size_type i=0; i < mfcomp.size(); ++i) {
672 if (mfcomp[i].op == mf_comp::DATA)
continue;
673 pfem fem = (mfcomp[i].pmf ? mfcomp[i].pmf->fem_of_element(cv)
675 pmat_elem_type pme2 = 0;
676 switch (mfcomp[i].op) {
681 case mf_comp::GRADGT:
682 case mf_comp::GRADGTINV:
685 case mf_comp::NONLIN: {
686 std::vector<pfem> ftab(1+mfcomp[i].auxmf.size());
688 for (
unsigned k=0; k < mfcomp[i].auxmf.size(); ++k)
689 ftab[k+1] = mfcomp[i].auxmf[k]->fem_of_element(cv);
692 case mf_comp::DATA: ;
694 if (pme == 0) pme = pme2;
698 if (pme == 0) pme = mat_elem_empty();
705 push_back_mfcomp_dimensions(
size_type cv,
const mf_comp& mc,
706 unsigned &d,
const bgeot::tensor_ranges &rng,
707 bgeot::tensor_ref &tref,
size_type tsz=1) {
708 if (mc.op == mf_comp::NONLIN) {
709 for (
size_type j=0; j < mc.nlt->sizes(cv).size(); ++j)
710 tsz = add_dim(rng, dim_type(d++), stride_type(tsz), tref);
711 }
else if (mc.op == mf_comp::DATA) {
713 tref = mc.data->tensor();
716 }
else if (mc.op == mf_comp::NORMAL) {
717 tsz = add_dim(rng, dim_type(d++), stride_type(tsz), tref);
718 }
else if (mc.op == mf_comp::GRADGT ||
719 mc.op == mf_comp::GRADGTINV) {
720 tsz = add_dim(rng, dim_type(d++), stride_type(tsz), tref);
721 tsz = add_dim(rng, dim_type(d++), stride_type(tsz), tref);
723 size_type target_dim = mc.pmf->fem_of_element(cv)->target_dim();
727 if (mc.vshape == mf_comp::VECTORIZED_SHAPE) {
728 if (target_dim == qdim) {
729 tsz = add_dim(rng, dim_type(d++), stride_type(tsz), tref);
730 tsz = add_dim(rng, dim_type(d++), stride_type(tsz), tref);
732 tsz = add_vdim(rng,dim_type(d),index_type(target_dim),
733 stride_type(tsz), tref);
736 }
else if (mc.vshape == mf_comp::MATRIXIZED_SHAPE) {
737 if (target_dim == qdim) {
738 tsz = add_dim(rng, dim_type(d++), stride_type(tsz), tref);
739 tsz = add_dim(rng, dim_type(d++), stride_type(tsz), tref);
740 tsz = add_dim(rng, dim_type(d++), stride_type(tsz), tref);
742 tsz = add_mdim(rng, dim_type(d), index_type(target_dim),
743 stride_type(tsz), tref);
746 }
else tsz = add_dim(rng, dim_type(d++), stride_type(tsz), tref);
747 if (mc.op == mf_comp::GRAD || mc.op == mf_comp::HESS) {
748 tsz = add_dim(rng, dim_type(d++), stride_type(tsz), tref);
750 if (mc.op == mf_comp::HESS) {
751 tsz = add_dim(rng, dim_type(d++), stride_type(tsz), tref);
757 void update_shape_with_inline_reduction(
size_type cv) {
758 fallback_red_uptodate =
false;
759 icb.tensor_bases.resize(mfcomp.size());
761 for (
size_type i=0; i < mfcomp.size(); ++i) {
765 mfcomp[i].push_back_dimensions(cv,rng);
766 push_back_mfcomp_dimensions(cv,mfcomp[i], d, rng, tref);
767 assert(tref.ndim() == rng.size() && d == rng.size());
768 if (mfcomp[i].reduction.size() == 0)
769 mfcomp[i].reduction.insert(
size_type(0), tref.ndim(),
' ');
770 if (mfcomp[i].op != mf_comp::DATA)
771 tref.set_base(icb.tensor_bases[i]);
772 tref.update_idx2mask();
773 if (mfcomp[i].reduction.size() != tref.ndim()) {
774 ASM_THROW_TENSOR_ERROR(
"wrong number of indices for the "<<
int(i+1)
775 <<
"th argument of the reduction "<< name()
776 <<
" (expected " <<
int(tref.ndim())
778 << mfcomp[i].reduction.size());
780 icb.red.insert(tref, mfcomp[i].reduction);
783 icb.red.result(tensor());
784 r_.resize(tensor().ndim());
785 for (dim_type i=0; i < tensor().ndim(); ++i) r_[i] = tensor().dim(i);
786 tsize = tensor().card();
791 void update_shape_with_expanded_tensor(
size_type cv) {
794 for (
size_type i=0; i < mfcomp.size(); ++i) {
795 tsize = stride_type(push_back_mfcomp_dimensions(cv, mfcomp[i], d, r_, tensor(), tsize));
797 assert(d == r_.size());
798 tensor().update_idx2mask();
801 void check_shape_update(
size_type cv, dim_type) {
802 const mesh_im& mi = mfcomp.get_im();
803 pintegration_method pim2;
805 bool fem_changed =
false;
806 pgt2 = mi.linked_mesh().trans_of_convex(cv);
807 pim2 = mi.int_method_of_element(cv);
812 cv_shape_update = cv;
813 shape_updated_ = (pme == 0);
815 shape_updated_ = shape_updated_ || child(i).is_shape_updated();
816 for (
size_type i=0; shape_updated_ ==
false && i < mfcomp.size(); ++i) {
817 if (mfcomp[i].pmf == 0)
continue;
819 shape_updated_ =
true; fem_changed =
true;
821 fem_changed = fem_changed ||
822 (mfcomp[i].pmf->fem_of_element(current_cv) !=
823 mfcomp[i].pmf->fem_of_element(cv));
825 shape_updated_ = shape_updated_ ||
826 (mfcomp[i].pmf->nb_basic_dof_of_element(current_cv) !=
827 mfcomp[i].pmf->nb_basic_dof_of_element(cv));
830 if (shape_updated_) {
833 for (
unsigned i=0; i < mfcomp.size(); ++i)
834 mfcomp[i].push_back_dimensions(cv, r_,
true);
836 if (fem_changed || shape_updated_) {
838 update_pmat_elem(cv);
840 if (shape_updated_ || fem_changed || pgt != pgt2 || pim != pim2) {
841 pgt = pgt2; pim = pim2;
842 pmec = mep(pme, pim, pgt, has_inline_reduction);
847 if (!shape_updated_)
return;
850 if (has_inline_reduction) {
851 update_shape_with_inline_reduction(cv_shape_update);
853 update_shape_with_expanded_tensor(cv_shape_update);
856 tensor().set_base(data_base);
858 void update_childs_required_shape() {
864 if (!fallback_red_uptodate) {
865 fallback_red_uptodate =
true;
868 tensor_ref tref; tensor_ranges r;
871 for (cnt=0; cnt < mfcomp.size() && mfcomp[cnt].op != mf_comp::DATA; ++cnt) {
872 mfcomp[cnt].push_back_dimensions(cv,r);
873 sz = push_back_mfcomp_dimensions(cv, mfcomp[cnt], d, r, tref, sz);
874 s += mfcomp[cnt].reduction;
876 fallback_red.clear();
877 tref.set_base(fallback_base);
878 fallback_red.insert(tref, s);
880 for ( ; cnt < mfcomp.size(); ++cnt) {
881 assert(mfcomp[cnt].op == mf_comp::DATA);
882 fallback_red.insert(mfcomp[cnt].data->tensor(), mfcomp[cnt].reduction);
884 fallback_red.prepare();
885 fallback_red.result(tensor());
886 assert(icb.red.reduced_range == fallback_red.reduced_range);
889 fallback_base = &(*t.begin());
890 fallback_red.do_reduction();
893 void exec_(
size_type cv, dim_type face) {
894 const mesh_im& mim = mfcomp.get_im();
895 for (
unsigned i=0; i < mfcomp.size(); ++i) {
896 if (mfcomp[i].op == mf_comp::DATA) {
898 for (
unsigned j=0; j < mfcomp[i].data->ranges().size(); ++j)
899 fullsz *= mfcomp[i].data->ranges()[j];
900 if (fullsz !=
size_type(mfcomp[i].data->tensor().card()))
901 ASM_THROW_TENSOR_ERROR(
"aaarg inline reduction will explode with non-full tensors. " 902 "Complain to the author, I was too lazy to do that properly");
905 icb.was_called =
false;
906 if (face == dim_type(-1)) {
907 pmec->gen_compute(t, mim.linked_mesh().points_of_convex(cv), cv,
908 has_inline_reduction ? &icb : 0);
910 pmec->gen_compute_on_face(t, mim.linked_mesh().points_of_convex(cv), face, cv,
911 has_inline_reduction ? &icb : 0);
914 if (has_inline_reduction && icb.was_called ==
false) {
915 do_post_reduction(cv);
916 data_base = &fallback_red.out_data[0];
917 }
else data_base = &(*t.begin());
918 GMM_ASSERT1(t.size() ==
size_type(tsize),
919 "Internal error: bad size " << t.size() <<
" should be " << tsize);
925 class ATN_tensor_from_dofs_data :
public ATN_tensor_w_data {
926 const base_asm_data *basm;
927 vdim_specif_list vdim;
928 multi_tensor_iterator mti;
930 std::vector< tensor_strides > e_str;
932 ATN_tensor_from_dofs_data(
const base_asm_data *basm_,
933 const vdim_specif_list& d) :
934 basm(basm_), vdim(d) {
936 void check_shape_update(
size_type cv, dim_type) {
937 shape_updated_ =
false;
938 r_.resize(vdim.size());
939 for (dim_type i=0; i < vdim.size(); ++i) {
940 if (vdim[i].is_mf_ref()) {
941 size_type nbde = vdim[i].pmf->nb_basic_dof_of_element(cv);
942 if (nbde != ranges()[i])
943 { r_[i] = unsigned(nbde); shape_updated_ =
true; }
944 }
else if (vdim[i].dim != ranges()[i]) {
945 r_[i] = unsigned(vdim[i].dim);
946 shape_updated_ =
true;
950 virtual void init_required_shape() { req_shape = tensor_shape(ranges()); }
954 ATN_tensor_w_data::reinit_();
955 mti.assign(tensor(),
true);
958 vdim.build_strides_for_cv(cv, e_r, e_str);
959 assert(e_r == ranges());
961 basm->copy_with_mti(e_str, mti, (vdim.nb_mf() >= 1) ? vdim[0].pmf : 0);
968 class ATN_symmetrized_tensor :
public ATN_tensor_w_data {
969 multi_tensor_iterator mti;
971 ATN_symmetrized_tensor(ATN_tensor& a) { add_child(a); }
972 void check_shape_update(
size_type , dim_type) {
973 if ((shape_updated_ = child(0).is_shape_updated())) {
974 if (child(0).ranges().size() != 2 ||
975 child(0).ranges()[0] != child(0).ranges()[1])
976 ASM_THROW_TENSOR_ERROR(
"can't symmetrize a non-square tensor " 977 "of sizes " << child(0).ranges());
978 r_ = child(0).ranges();
981 void update_childs_required_shape() {
982 tensor_shape ts = req_shape;
983 tensor_shape ts2 = req_shape;
984 index_set perm(2); perm[0] = 1; perm[1] = 0; ts2.permute(perm);
985 ts.merge(ts2,
false);
986 tensor_mask dm; dm.set_triangular(ranges()[0],0,1);
987 tensor_shape tsdm(2); tsdm.push_mask(dm);
988 ts.merge(tsdm,
true);
989 child(0).merge_required_shape(ts);
994 req_shape.set_full(ranges());
995 ATN_tensor_w_data::reinit0();
996 mti.assign(child(0).tensor(),
true);
999 std::fill(data.begin(), data.end(), 0.);
1001 index_type n = ranges()[0];
1003 index_type i=mti.index(0), j=mti.index(1);
1004 data[i*n+j]=data[j*n+i]=mti.p(0);
1005 }
while (mti.qnext1());
1010 template<
class UnaryOp>
class ATN_unary_op_tensor
1011 :
public ATN_tensor_w_data {
1012 multi_tensor_iterator mti;
1014 ATN_unary_op_tensor(ATN_tensor& a) { add_child(a); }
1015 void check_shape_update(
size_type , dim_type) {
1016 if ((shape_updated_ = (ranges() != child(0).ranges())))
1017 r_ = child(0).ranges();
1021 ATN_tensor_w_data::reinit0();
1022 mti.assign(tensor(), child(0).tensor(),
false);
1027 mti.p(0) = UnaryOp()(mti.p(1));
1028 }
while (mti.qnext2());
1033 class ATN_tensors_sum_scaled :
public ATN_tensor_w_data {
1034 std::vector<multi_tensor_iterator> mti;
1035 std::vector<scalar_type> scales;
1037 ATN_tensors_sum_scaled(ATN_tensor& t1, scalar_type s1) {
1039 scales.resize(1); scales[0]=s1;
1041 void push_scaled_tensor(ATN_tensor& t, scalar_type s) {
1042 add_child(t); scales.push_back(s);
1044 void check_shape_update(
size_type , dim_type) {
1045 if ((shape_updated_ = child(0).is_shape_updated()))
1046 r_ = child(0).ranges();
1048 if (ranges() != child(i).ranges())
1049 ASM_THROW_TENSOR_ERROR(
"can't add two tensors of sizes " <<
1050 ranges() <<
" and " << child(i).ranges());
1052 void apply_scale(scalar_type s) {
1053 for (
size_type i=0; i < scales.size(); ++i) scales[i] *= s;
1055 ATN_tensors_sum_scaled* is_tensors_sum_scaled() {
return this; }
1058 ATN_tensor_w_data::reinit0();
1059 mti.resize(nchilds());
1061 mti[i].assign(tensor(), child(i).tensor(),
false);
1068 std::fill(data.begin(), data.end(), 0.);
1071 mti[0].p(0) = mti[0].p(1)*scales[0];
1072 }
while (mti[0].qnext2());
1073 for (
size_type i=1; i < nchilds(); ++i) {
1076 mti[i].p(0) = mti[i].p(0)+mti[i].p(1)*scales[i];
1077 }
while (mti[i].qnext2());
1082 class ATN_tensor_scalar_add :
public ATN_tensor_w_data {
1084 multi_tensor_iterator mti;
1087 ATN_tensor_scalar_add(ATN_tensor& a, scalar_type v_,
int sgn_) :
1088 v(v_), sgn(sgn_) { add_child(a); }
1089 void check_shape_update(
size_type , dim_type) {
1090 if ((shape_updated_ = child(0).is_shape_updated()))
1091 r_ = child(0).ranges();
1095 ATN_tensor_w_data::reinit_();
1096 mti.assign(tensor(), child(0).tensor(),
false);
1099 std::fill(data.begin(), data.end(), v);
1103 mti.p(0) += mti.p(1);
1104 else mti.p(0) -= mti.p(1);
1105 }
while (mti.qnext2());
1109 class ATN_print_tensor :
public ATN {
1112 ATN_print_tensor(ATN_tensor& a, std::string n_)
1113 : name(n_) { add_child(a); }
1116 void exec_(
size_type cv, dim_type face) {
1117 multi_tensor_iterator mti(child(0).tensor(),
true);
1118 cout <<
"------- > evaluation of " << name <<
", at" << endl;
1119 cout <<
"convex " << cv;
1120 if (face != dim_type(-1)) cout <<
", face " << int(face);
1122 cout <<
" size = " << child(0).ranges() << endl;
1126 for (
size_type i=0; i < child(0).ranges().size(); ++i)
1127 cout <<(i>0 ?
"," :
"") << mti.index(dim_type(i));
1128 cout <<
"] = " << mti.p(0) << endl;
1129 }
while (mti.qnext1());
1140 std::string asm_tokenizer::syntax_err_print() {
1142 if (tok_pos - err_msg_mark > 80) err_msg_mark = tok_pos - 40;
1143 if (str.length() - err_msg_mark < 80) s = tok_substr(err_msg_mark, str.length());
1144 else { s = tok_substr(err_msg_mark,err_msg_mark+70); s.append(
" ... (truncated)"); }
1145 s +=
"\n" + std::string(std::max(
int(tok_pos - err_msg_mark),0),
'-') +
"^^";
1149 void asm_tokenizer::get_tok() {
1152 while (tok_pos < str.length() && isspace(str[tok_pos])) ++tok_pos;
1153 if (tok_pos == str.length()) {
1154 curr_tok_type = END; tok_len = 0;
1155 }
else if (strchr(
"{}(),;:=-.*/+", str[tok_pos])) {
1156 curr_tok_type = tok_type_enum(str[tok_pos]); tok_len = 1;
1157 }
else if (str[tok_pos] ==
'$' || str[tok_pos] ==
'#' || str[tok_pos] ==
'%') {
1158 curr_tok_type = str[tok_pos] ==
'$' ? ARGNUM_SELECTOR :
1159 (str[tok_pos] ==
'#' ? MFREF : IMREF);
1162 while (isdigit(str[tok_pos+tok_len])) {
1164 curr_tok_ival += str[tok_pos+tok_len] -
'0';
1168 }
else if (isalpha(str[tok_pos])) {
1169 curr_tok_type = IDENT;
1171 while (isalnum(str[tok_pos+tok_len]) || str[tok_pos+tok_len] ==
'_') ++tok_len;
1172 }
else if (isdigit(str[tok_pos])) {
1173 curr_tok_type = NUMBER;
1175 curr_tok_dval = strtod(&str[0]+tok_pos, &p);
1176 tok_len = p - &str[0] - tok_pos;
1178 if (tok_pos < str.length())
1179 curr_tok = str.substr(tok_pos, tok_len);
1184 const mesh_fem& generic_assembly::do_mf_arg_basic() {
1185 if (tok_type() != MFREF) ASM_THROW_PARSE_ERROR(
"expecting mesh_fem reference");
1186 if (tok_mfref_num() >= mftab.size())
1187 ASM_THROW_PARSE_ERROR(
"reference to a non-existant mesh_fem #" << tok_mfref_num()+1);
1188 const mesh_fem& mf_ = *mftab[tok_mfref_num()]; advance();
1192 const mesh_fem& generic_assembly::do_mf_arg(std::vector<const mesh_fem*> * multimf) {
1193 if (!multimf) advance();
1194 accept(OPEN_PAR,
"expecting '('");
1195 const mesh_fem &mf_ = do_mf_arg_basic();
1197 multimf->resize(1); (*multimf)[0] = &mf_;
1198 while (advance_if(COMMA)) {
1199 if (tok_type() != MFREF) ASM_THROW_PARSE_ERROR(
"expecting mesh_fem reference");
1200 if (tok_mfref_num() >= mftab.size())
1201 ASM_THROW_PARSE_ERROR(
"reference to a non-existant mesh_fem #" << tok_mfref_num()+1);
1202 multimf->push_back(mftab[tok_mfref_num()]); advance();
1205 accept(CLOSE_PAR,
"expecting ')'");
1210 std::string generic_assembly::do_comp_red_ops() {
1212 if (advance_if(OPEN_PAR)) {
1215 if (tok_type() == COLON) {
1216 s.push_back(
' '); advance(); j++;
1217 }
else if (tok_type() == IDENT) {
1218 if ((tok().length()==1 && isalpha(tok()[0])) || islower(tok()[0])) {
1219 s.push_back(tok()[0]); advance(); j++;
1220 }
else ASM_THROW_PARSE_ERROR(
"invalid reduction index '" << tok() <<
1221 "', only lower case characters allowed");
1223 }
while (advance_if(COMMA));
1224 accept(CLOSE_PAR,
"expecting ')'");
1229 static mf_comp::field_shape_type get_shape(
const std::string &s) {
1230 if (s[0] ==
'v')
return mf_comp::VECTORIZED_SHAPE;
1231 else if (s[0] ==
'm')
return mf_comp::MATRIXIZED_SHAPE;
1232 else return mf_comp::PLAIN_SHAPE;
1235 ATN_tensor* generic_assembly::do_comp() {
1236 accept(OPEN_PAR,
"expecting '('");
1238 bool in_data =
false;
1246 if (tok_type() == IMREF) {
1247 if (tok_imref_num() >= imtab.size())
1248 ASM_THROW_PARSE_ERROR(
"reference to a non-existant mesh_im %" << tok_imref_num()+1);
1249 what.set_im(*imtab[tok_imref_num()]); advance();
1250 accept(COMMA,
"expecting ','");
1252 what.set_im(*imtab[0]);
1255 if (tok_type() == CLOSE_PAR)
break;
1256 if (tok_type() != IDENT) ASM_THROW_PARSE_ERROR(
"expecting Base or Grad or Hess, Normal, etc..");
1257 std::string f = tok();
1258 const mesh_fem *pmf = 0;
1259 if (f.compare(
"Base")==0 || f.compare(
"vBase")==0 || f.compare(
"mBase")==0) {
1261 what.push_back(mf_comp(&what, pmf, mf_comp::BASE, get_shape(f)));
1262 }
else if (f.compare(
"Grad")==0 || f.compare(
"vGrad")==0 || f.compare(
"mGrad")==0) {
1264 what.push_back(mf_comp(&what, pmf, mf_comp::GRAD, get_shape(f)));
1265 }
else if (f.compare(
"Hess")==0 || f.compare(
"vHess")==0 || f.compare(
"mHess")==0) {
1267 what.push_back(mf_comp(&what, pmf, mf_comp::HESS, get_shape(f)));
1268 }
else if (f.compare(
"NonLin")==0) {
1271 if (tok_type() == ARGNUM_SELECTOR) { num = tok_argnum(); advance(); }
1272 if (num >= innonlin.size()) ASM_THROW_PARSE_ERROR(
"NonLin$" << num <<
" does not exist");
1273 std::vector<const mesh_fem*> allmf;
1274 pmf = &do_mf_arg(&allmf); what.push_back(mf_comp(&what, allmf, innonlin[num]));
1275 }
else if (f.compare(
"Normal") == 0) {
1277 accept(OPEN_PAR,
"expecting '('"); accept(CLOSE_PAR,
"expecting ')'");
1278 pmf = 0; what.push_back(mf_comp(&what, pmf, mf_comp::NORMAL, mf_comp::PLAIN_SHAPE));
1279 }
else if (f.compare(
"GradGT") == 0 ||
1280 f.compare(
"GradGTInv") == 0) {
1282 accept(OPEN_PAR,
"expecting '('"); accept(CLOSE_PAR,
"expecting ')'");
1284 what.push_back(mf_comp(&what, pmf,
1285 f.compare(
"GradGT") == 0 ?
1287 mf_comp::GRADGTINV, mf_comp::PLAIN_SHAPE));
1289 if (vars.find(f) != vars.end()) {
1290 what.push_back(mf_comp(&what, vars[f]));
1294 ASM_THROW_PARSE_ERROR(
"expecting Base, Grad, vBase, NonLin ...");
1298 if (!in_data && f[0] !=
'v' && f[0] !=
'm' && pmf && pmf->get_qdim() != 1 && f.compare(
"NonLin")!=0) {
1299 ASM_THROW_PARSE_ERROR(
"Attempt to use a vector mesh_fem as a scalar mesh_fem");
1301 what.back().reduction = do_comp_red_ops();
1302 }
while (advance_if(PRODUCT));
1303 accept(CLOSE_PAR,
"expecting ')'");
1305 return record(std::make_unique<ATN_computed_tensor>(what));
1308 void generic_assembly::do_dim_spec(vdim_specif_list& lst) {
1310 accept(OPEN_PAR,
"expecting '('");
1312 if (tok_type() == IDENT) {
1313 if (tok().compare(
"mdim")==0) lst.push_back(vdim_specif(do_mf_arg().linked_mesh().dim()));
1314 else if (tok().compare(
"qdim")==0) lst.push_back(vdim_specif(do_mf_arg().get_qdim()));
1315 else ASM_THROW_PARSE_ERROR(
"expecting mdim(#mf) or qdim(#mf) or a number or a mesh_fem #id");
1316 }
else if (tok_type() == NUMBER) {
1317 lst.push_back(vdim_specif(tok_number_ival()+1));
1319 }
else if (tok_type() == MFREF) {
1320 lst.push_back(vdim_specif(&do_mf_arg_basic()));
1321 }
else if (tok_type() != CLOSE_PAR) ASM_THROW_PARSE_ERROR(
"expecting mdim(#mf) or qdim(#mf) or a number or a mesh_fem #id");
1324 if (advance_if(CLOSE_PAR))
break;
1325 accept(COMMA,
"expecting ',' or ')'");
1330 ATN_tensor* generic_assembly::do_data() {
1333 if (tok_type() != OPEN_PAR) {
1334 if (tok_type() != ARGNUM_SELECTOR)
1335 ASM_THROW_PARSE_ERROR(
"expecting dataset number");
1336 datanum = tok_argnum();
1339 if (datanum >= indata.size())
1340 ASM_THROW_PARSE_ERROR(
"wrong dataset number: " << datanum);
1342 vdim_specif_list sz;
1345 if (sz.nbelt() != indata[datanum]->vect_size())
1346 ASM_THROW_PARSE_ERROR(
"invalid size for data argument " << datanum+1 <<
1347 " real size is " << indata[datanum]->vect_size()
1348 <<
" expected size is " << sz.nbelt());
1349 return record(std::make_unique<ATN_tensor_from_dofs_data>(indata[datanum].
get(), sz));
1352 std::pair<ATN_tensor*, std::string>
1353 generic_assembly::do_red_ops(ATN_tensor* t) {
1356 if (advance_if(OPEN_PAR)) {
1359 if (tok_type() == COLON) {
1360 s.push_back(
' '); advance(); j++;
1361 }
else if (tok_type() == NUMBER) {
1362 t = record(std::make_unique<ATN_sliced_tensor>(*t, dim_type(j),
1363 tok_number_ival())); advance();
1364 }
else if (tok_type() == IDENT) {
1365 if ((tok().length()==1 && isalpha(tok()[0])) || islower(tok()[0])) {
1366 s.push_back(tok()[0]); advance(); j++;
1367 }
else ASM_THROW_PARSE_ERROR(
"invalid reduction index '" << tok() <<
1368 "', only lower case chars allowed");
1370 }
while (advance_if(COMMA));
1371 accept(CLOSE_PAR,
"expecting ')'");
1373 return std::pair<ATN_tensor*,std::string>(t,s);
1382 tnode generic_assembly::do_tens() {
1385 if (advance_if(OPEN_PAR)) {
1387 accept(CLOSE_PAR,
"expecting ')'");
1388 }
else if (tok_type() == NUMBER) {
1389 t.assign(tok_number_dval()); advance();
1390 }
else if (tok_type() == IDENT) {
1391 if (vars.find(tok()) != vars.end()) {
1392 t.assign(vars[tok()]); advance();
1393 }
else if (tok().compare(
"comp")==0) {
1394 advance(); t.assign(do_comp());
1395 }
else if (tok().compare(
"data")==0) {
1396 advance(); t.assign(do_data());
1397 }
else if (tok().compare(
"sym")==0) {
1399 tnode t2 = do_expr();
1400 if (t2.type() != tnode::TNTENSOR)
1401 ASM_THROW_PARSE_ERROR(
"can't symmetrise a scalar!");
1402 t.assign(record(std::make_unique<ATN_symmetrized_tensor>(*t2.tensor())));
1403 }
else ASM_THROW_PARSE_ERROR(
"unknown identifier: " << tok());
1404 }
else ASM_THROW_PARSE_ERROR(
"unexpected token: " << tok());
1414 tnode generic_assembly::do_prod() {
1415 reduced_tensor_arg_type ttab;
1418 tnode t = do_tens();
1419 if (t.type() == tnode::TNCONST) {
1420 if (ttab.size() == 0)
return t;
1421 else ASM_THROW_PARSE_ERROR(
"can't mix tensor and scalar into a " 1422 "reduction expression");
1424 ttab.push_back(do_red_ops(t.tensor()));
1425 }
while (advance_if(PRODUCT));
1426 if (ttab.size() == 1 && ttab[0].second.length() == 0) {
1427 return tnode(ttab[0].first);
1429 return tnode(record(std::make_unique<ATN_reduced_tensor>(ttab)));
1435 tnode generic_assembly::do_prod_trans() {
1436 tnode t = do_prod();
1437 while (advance_if(OPEN_BRACE)) {
1440 dal::bit_vector check_permut;
1441 if (t.type() != tnode::TNTENSOR)
1442 ASM_THROW_PARSE_ERROR(
"can't use reorder braces on a constant!");
1445 if (tok_type() == COLON) i = j;
1446 else if (tok_type() == NUMBER) i = tok_number_ival(1000);
1447 else ASM_THROW_PARSE_ERROR(
"only numbers or colons allowed here");
1448 if (check_permut.is_in(i)) {
1449 t = tnode(record(std::make_unique<ATN_diagonal_tensor>(*t.tensor(), dim_type(i),
1451 check_permut.add(j);
1452 reorder.push_back(dim_type(j));
1454 check_permut.add(i);
1455 reorder.push_back(dim_type(i));
1458 if (advance_if(CLOSE_BRACE))
break;
1459 accept(COMMA,
"expecting ','");
1461 if (check_permut.first_false() != reorder.size()) {
1462 cerr << check_permut << endl;
1463 cerr << vref(reorder) << endl;
1464 ASM_THROW_PARSE_ERROR(
"you did not give a real permutation:" 1467 t = tnode(record(std::make_unique<ATN_permuted_tensor>(*t.tensor(), reorder)));
1475 tnode generic_assembly::do_term() {
1478 tnode t = do_prod_trans();
1481 if (advance_if(MULTIPLY)) mult =
true;
1482 else if (advance_if(DIVIDE)) mult =
false;
1484 tnode t2 = do_prod();
1485 if (mult ==
false && t.type() == tnode::TNCONST
1486 && t2.type() == tnode::TNTENSOR)
1487 ASM_THROW_PARSE_ERROR(
"can't divide a constant by a tensor");
1488 if (t.type() == tnode::TNTENSOR && t2.type() == tnode::TNTENSOR) {
1489 ASM_THROW_PARSE_ERROR(
"tensor term-by-term productor division not " 1490 "implemented yet! are you sure you need it ?");
1491 }
else if (t.type() == tnode::TNCONST && t2.type() == tnode::TNCONST) {
1493 t.assign(t.xval()*t2.xval());
1496 t.assign(t.xval()/t2.xval());
1499 if (t.type() != tnode::TNTENSOR) std::swap(t,t2);
1500 scalar_type v = t2.xval();
1502 if (v == 0.) { ASM_THROW_PARSE_ERROR(
"can't divide by zero"); }
1505 if (t.tensor()->is_tensors_sum_scaled() && !t.tensor()->is_frozen()) {
1506 t.tensor()->is_tensors_sum_scaled()->apply_scale(v);
1508 t.assign(record(std::make_unique<ATN_tensors_sum_scaled>(*t.tensor(), v)));
1521 tnode generic_assembly::do_expr() {
1524 if (advance_if(MINUS)) negt =
true;
1525 tnode t = do_term();
1527 if (t.type() == tnode::TNCONST) t.assign(-t.xval());
1528 else t.assign(record(std::make_unique<ATN_tensor_scalar_add>(*t.tensor(), 0., -1)));
1532 if (advance_if(PLUS)) plus = +1;
1533 else if (advance_if(MINUS)) plus = -1;
1535 tnode t2 = do_term();
1536 if (t.type() == tnode::TNTENSOR && t2.type() == tnode::TNTENSOR) {
1537 if (!t.tensor()->is_tensors_sum_scaled() || t.tensor()->is_frozen()) {
1538 t.assign(record(std::make_unique<ATN_tensors_sum_scaled>(*t.tensor(), +1)));
1540 t.tensor()->is_tensors_sum_scaled()
1541 ->push_scaled_tensor(*t2.tensor(), scalar_type(plus));
1542 }
else if (t.type() == tnode::TNCONST && t2.type() == tnode::TNCONST) {
1543 t.assign(t.xval()+t2.xval()*plus);
1546 if (t.type() != tnode::TNTENSOR)
1547 { std::swap(t,t2);
if (plus<0) tsgn = -1; }
1548 else if (plus<0) t2.assign(-t2.xval());
1549 t.assign(record(std::make_unique<ATN_tensor_scalar_add>(*t.tensor(), t2.xval(),
1561 void generic_assembly::do_instr() {
1562 enum { wALIAS, wOUTPUT_ARRAY, wOUTPUT_MATRIX, wPRINT, wERROR }
1567 if (tok_type() != IDENT) ASM_THROW_PARSE_ERROR(
"expecting identifier");
1568 if (vars.find(tok()) != vars.end())
1569 ASM_THROW_PARSE_ERROR(
"redefinition of identifier " << tok());
1578 vdim_specif_list vds;
1580 if (ident.compare(
"print") == 0) {
1581 print_mark = tok_mark();
1583 }
else if (tok_type() == ARGNUM_SELECTOR ||
1584 tok_type() == OPEN_PAR) {
1585 if (tok_type() == ARGNUM_SELECTOR) {
1586 arg_num = tok_argnum();
1588 }
else { arg_num = 0; }
1593 if (ident.compare(
"V")==0) {
1594 what = wOUTPUT_ARRAY;
1595 if (arg_num >= outvec.size())
1596 { outvec.resize(arg_num+1); outvec[arg_num] = 0; }
1598 if (outvec[arg_num] == 0) {
1599 if (vec_fact != 0) {
1600 tensor_ranges r(vds.size());
1601 for (
size_type i=0; i < vds.size(); ++i)
1602 r[i] =
unsigned(vds[i].dim);
1603 outvec[arg_num] = std::shared_ptr<base_asm_vec>(std::shared_ptr<base_asm_vec>(), vec_fact->create_vec(r));
1605 else ASM_THROW_PARSE_ERROR(
"output vector $" << arg_num+1
1606 <<
" does not exist");
1608 }
else if (vds.nb_mf()==2 && vds.size() == 2 && ident.compare(
"M")==0) {
1609 what = wOUTPUT_MATRIX;
1610 if (arg_num >= outmat.size())
1611 { outmat.resize(arg_num+1); outmat[arg_num] = 0; }
1613 if (outmat[arg_num] == 0) {
1615 outmat[arg_num] = std::shared_ptr<base_asm_mat>
1616 (std::shared_ptr<base_asm_mat>(),
1617 mat_fact->create_mat(vds[0].pmf->nb_dof(),
1618 vds[1].pmf->nb_dof()));
1619 else ASM_THROW_PARSE_ERROR(
"output matrix $" << arg_num+1
1620 <<
" does not exist");
1622 }
else ASM_THROW_PARSE_ERROR(
"not a valid output statement");
1626 }
else if (advance_if(EQUAL)) {
1628 }
else ASM_THROW_PARSE_ERROR(
"missing '=' or ':='");
1630 tnode t = do_expr();
1631 if (t.type() != tnode::TNTENSOR)
1632 ASM_THROW_PARSE_ERROR(
"left hand side is a constant, not a tensor!");
1636 record_out(std::make_unique<ATN_print_tensor>(*t.tensor(), tok_substr(print_mark,
1639 case wOUTPUT_ARRAY: {
1640 record_out(outvec[arg_num]->build_output_tensor(*t.tensor(), vds));
1642 case wOUTPUT_MATRIX: {
1643 record_out(outmat[arg_num]->build_output_tensor(*t.tensor(),
1648 vars[ident] = t.tensor(); t.tensor()->freeze();
1650 default: GMM_ASSERT3(
false,
"");
break;
1655 struct atn_number_compare {
1656 bool operator()(
const std::unique_ptr<ATN_tensor> &a,
1657 const std::unique_ptr<ATN_tensor> &b) {
1658 assert(a.get() && b.get());
1659 return (a->number() < b->number());
1663 struct outvar_compare {
1664 bool operator()(
const std::unique_ptr<ATN> &a,
1665 const std::unique_ptr<ATN> &b) {
1666 assert(a.get() && b.get());
1667 return (a->number() < b->number());
1671 void generic_assembly::parse() {
1672 if (parse_done)
return;
1674 if (tok_type() == END)
break;
1676 }
while (advance_if(SEMICOLON));
1677 if (tok_type() != END) ASM_THROW_PARSE_ERROR(
"unexpected token: '" 1679 if (outvars.size() == 0) cerr <<
"warning: assembly without output\n";
1683 for (
size_type i=0; i < outvars.size(); ++i)
1684 outvars[i]->set_number(gcnt);
1686 std::sort(atn_tensors.begin(), atn_tensors.end(), atn_number_compare());
1687 std::sort(outvars.begin(), outvars.end(), outvar_compare());
1690 while (atn_tensors.size()
1691 && atn_tensors.back()->number() == unsigned(-1)) {
1692 cerr <<
"warning: the expression " << atn_tensors.back()->name()
1693 <<
" won't be evaluated since it is not used!\n";
1694 atn_tensors.pop_back();
1700 void generic_assembly::exec(
size_type cv, dim_type face) {
1701 bool update_shapes =
false;
1702 for (
size_type i=0; i < atn_tensors.size(); ++i) {
1703 atn_tensors[i]->check_shape_update(cv,face);
1704 update_shapes = (update_shapes || atn_tensors[i]->is_shape_updated());
1714 if (update_shapes) {
1721 for (
auto &&t : atn_tensors)
1722 t->init_required_shape();
1724 for (
auto &&v : outvars)
1725 v->update_childs_required_shape();
1728 atn_tensors[i]->update_childs_required_shape();
1730 for (
auto &&t : atn_tensors)
1733 for (
auto &&v : outvars)
1736 for (
auto &&t : atn_tensors)
1738 for (
auto &&v : outvars)
1742 struct cv_fem_compare {
1743 const std::vector<const mesh_fem *> &mf;
1744 cv_fem_compare(
const std::vector<const mesh_fem *>& mf_) : mf(mf_) {}
1746 for (
size_type i=0; i < mf.size(); ++i) {
1747 pfem pfa(mf[i]->fem_of_element(a));
1748 pfem pfb(mf[i]->fem_of_element(b));
1750 unsigned nba = unsigned(pfa->nb_dof(a));
1751 unsigned nbb = unsigned(pfb->nb_dof(b));
1754 }
else if (nba > nbb) {
1756 }
else if (pfa < pfb) {
1767 static void get_convex_order(
const dal::bit_vector& cvlst,
1768 const std::vector<const mesh_im *>& imtab,
1769 const std::vector<const mesh_fem *>& mftab,
1770 const dal::bit_vector& candidates,
1771 std::vector<size_type>& cvorder) {
1772 cvorder.reserve(candidates.card()); cvorder.resize(0);
1774 for (dal::bv_visitor cv(candidates); !cv.finished(); ++cv) {
1775 if (cvlst.is_in(cv) &&
1776 imtab[0]->int_method_of_element(cv) !=
im_none()) {
1778 for (
size_type i=0; i < mftab.size(); ++i) {
1779 if (!mftab[i]->convex_index().is_in(cv)) {
1786 cvorder.push_back(cv);
1788 }
else if (!imtab[0]->linked_mesh().convex_index().is_in(cv)) {
1789 ASM_THROW_ERROR(
"the convex " << cv <<
" is not part of the mesh");
1797 void generic_assembly::consistency_check() {
1799 if (imtab.size() == 0)
1800 ASM_THROW_ERROR(
"no mesh_im (integration methods) given for assembly!");
1801 const mesh& m = imtab[0]->linked_mesh();
1802 for (
unsigned i=0; i < mftab.size(); ++i) {
1803 if (&mftab[i]->linked_mesh() != &m)
1804 ASM_THROW_ERROR(
"the mesh_fem/mesh_im live on different meshes!");
1806 for (
unsigned i=0; i < imtab.size(); ++i) {
1807 if (&imtab[i]->linked_mesh() != &m)
1808 ASM_THROW_ERROR(
"the mesh_fem/mesh_im live on different meshes!");
1810 if (imtab.size() == 0)
1811 ASM_THROW_ERROR(
"no integration method !");
1815 std::vector<size_type> cv;
1816 r.
from_mesh(imtab.at(0)->linked_mesh());
1817 r.error_if_not_homogeneous();
1820 consistency_check();
1821 get_convex_order(imtab.at(0)->convex_index(), imtab, mftab, r.
index(), cv);
1824 for (
size_type i=0; i < cv.size(); ++i) {
1825 mesh_region::face_bitset nf = r[cv[i]];
1826 dim_type f = dim_type(-1);
1829 if (nf[0]) exec(cv[i],f);
structure used to hold a set of convexes and/or convex faces.
pmat_elem_type mat_elem_nonlinear(pnonlinear_elem_term, std::vector< pfem > pfi)
Give a pointer to the structure describing the elementary matrix which compute the integral of a nonl...
const dal::bit_vector & index() const
Index of the region convexes, or the convexes from the partition on the current thread.
pmat_elem_type mat_elem_unit_normal(void)
Give a pointer to the structure describing the elementary matrix which compute the unit normal on the...
pmat_elem_type mat_elem_grad_geotrans(bool inverted)
Return the gradient of the geometrical transformation ("K" in the getfem++ kernel doc...
void assembly(const mesh_region ®ion=mesh_region::all_convexes())
do the assembly on the specified region (boundary or set of convexes)
const mesh_region & from_mesh(const mesh &m) const
For regions which have been built with just a number 'id', from_mesh(m) sets the current region to 'm...
size_t size_type
used as the common size type in the library
pmat_elem_type mat_elem_hessian(pfem pfi)
Give a pointer to the structure describing the elementary matrix which compute the integral of the he...
pmat_elem_type mat_elem_base(pfem pfi)
Give a pointer to the structure describing the elementary matrix which compute the integral of the ba...
pintegration_method im_none(void)
return IM_NONE
this is the above solutions for linux, but it still needs to be tested.
pmat_elem_type mat_elem_grad(pfem pfi)
Give a pointer to the structure describing the elementary matrix which compute the integral of the gr...
elementary computations (used by the generic assembly).
GEneric Tool for Finite Element Methods.
std::shared_ptr< const getfem::virtual_fem > pfem
type of pointer on a fem description
gmm::uint16_type short_type
used as the common short type integer in the library
gmm interface for fortran BLAS.
pmat_elem_type mat_elem_product(pmat_elem_type a, pmat_elem_type b)
Give a pointer to the structure describing the elementary matrix which computes the integral of produ...
Generic assembly implementation.
std::shared_ptr< const bgeot::geometric_trans > pgeometric_trans
pointer type for a geometric transformation