From f27fb702a4d976ee6b422a4b931d1bb5d74a3b96 Mon Sep 17 00:00:00 2001 From: jfa Date: Fri, 1 Feb 2008 08:07:48 +0000 Subject: [PATCH] NPAL17873: SMESH add a triangle instead of make quadrangles only. --- .../gui/SMESH/images/a-averagelength.png | Bin 10905 -> 13679 bytes .../gui/SMESH/input/1d_meshing_hypo.doc | 17 +++-- idl/SMESH_BasicHypothesis.idl | 24 ++++++- .../StdMeshers_CompositeSegment_1D.cxx | 7 +- src/StdMeshers/StdMeshers_LocalLength.cxx | 59 +++++++++++++--- src/StdMeshers/StdMeshers_LocalLength.hxx | 5 +- src/StdMeshers/StdMeshers_RadialPrism_3D.cxx | 4 +- src/StdMeshers/StdMeshers_Regular_1D.cxx | 64 ++++++++++++++---- src/StdMeshers/StdMeshers_Regular_1D.hxx | 18 ++--- .../StdMeshersGUI_StdHypothesisCreator.cxx | 11 ++- src/StdMeshersGUI/StdMeshers_msg_en.po | 3 + src/StdMeshers_I/StdMeshers_LocalLength_i.cxx | 52 +++++++++++--- src/StdMeshers_I/StdMeshers_LocalLength_i.hxx | 6 ++ 13 files changed, 217 insertions(+), 53 deletions(-) diff --git a/doc/salome/gui/SMESH/images/a-averagelength.png b/doc/salome/gui/SMESH/images/a-averagelength.png index dc007eb575faa7564fa27da30eba60d5de4c1242..70e2afd26dd6ed299365193b2073c28946ed0695 100755 GIT binary patch literal 13679 zcma*ObyQqU_cch6;1V=Ia3?sy-Q8U~xO?!Xkp!2<-JKx8HMj>04vo9J`}Fg?GvAuu z`rdEW{L|}pT{%^Cs?OPG?|Z{lm1WRTh)`f)V9@1cCDmbI-jP9%$B^Dbe*@>ft-!#L z!^laBY3c@_WO)K{`=18ae3_H5$dRe9sYUn6NrD3Ck>m+9Rhh2y#e?!Sixk3xDf7=I z5a8jN00^kLVJJvb0zWYF5 zi->}J7Cs$nqkn4vLvO})JG4a;0da>%uh}OSo}W8=&-CF91-J_)J}Qx8;W}x?zl#=+ zl;_AzjH0JT#n%hh(>6<7`I{kBU1+@C{I^j{$Nung_^DTohL$cN)?9X&Tu-BK6fCNG ztjb6;YIy|KpwK&?+IN|^R#TwWl#i_-;H0m~t>0JEy_SuSfp13~mPQ}PU_|L}G}5Vf zrK6=QP%qb{rY1Z&5C_5sZAMEoY@HSK4;q1`ZAjSm(S4%cQ(4v#28#zu7j%z}p-0Pu z^-cOrJ+?sc#P|xglxo(E8sGY9nQfC$lB6qXDHa2M{D@+zmm#Mgu+~%^o(Ss%6BdaF zPg-wwzdWTfRgv~|trIfm75*b#$3^>s(K7PK_;E0@=cruD z#r3rqYf`UL2iTbL?Au`9rQ!;^!tGp%BQs=6ZYu+ZPGR{)>4dM`x$s{pq zmFhPTot~c3P`i3X$K@p>G8 zJS9T6@cW}m)isp$^d38&c@d&4t2jE>D1Ne4@9^2?x}QHVky7}UkBm?+KGu8?mGyfl zVQd^T_M^ZWwmMB<@$AnHgEHoy5#an;8kfBR{8OWav8*dkHG16T4v*AYi}tWM+k`}W zauIW_)pBh!OFpULcSeNv>*g~#C`{FM-WlM2S4!o|N%!rOnCZ3x|6+$<#WKO^TB%g%d@131I zKORj0J716n^SAOKA939+-Q>6ui3eyiTXjBv7qa3Z)+6bBUQ;W@8VyD_P@fvTcG~Hm z8lwzc{*VZyLqd)4~cgg3Hh{XoZnCtwkN*8#uWt9iK1Dmox!t4HlE|lvanR=v7y1YtI#QAq}2PY zS9f7Al8`=oWcoKlZYn)b5w53|WxOW3u~mSCWH5lPZIZ(%qxO=fyu(g!!SsQ% zr)PA|FP%CnF9?)^WxVmU75=F*NiC*u=c5t3;@eOrU*Qk=>p4+*i; zoz4^;EH%Ps;;CGrq_P~a($)H0mKIDhree)b*_mcHWD=F)jJ5$S0s|C*{1ZE$l^{zP_ z#jLexqv_59rKONtAO7~{fmgr#inA6C-%sElk6VQ2@Fk5mg;Q=+-Zf@=v?7z)B2_$j zdO2~+yr48&iKh#zZG!rfzi3%#M7>=4H?J^}5X41Xm9G_&`P>ao*bUMpmy4wp+UL{3 z0zVp9a%}k5HHkwPA-KB>v{pd5B_67^=Xi~z(7(qACiKi@1+&*3G`7)bjs_uoD%WOD zMXg`X+F&tx!d)%i{+BP1FYkKBR?riMquSaKQG#Ki^V%joJhzioqg&A{A{WnQFlwF9 zaP4D;iAgHhWm@x1DM+T48-0as$?Wb7bnhQ_llyiK9#31~=mpvA-o#7D#DzQ& z)I?s!LTy!^x~B5u{~8;AS6&MXdwB9Ldt0ioTZ5fGX|spivmRm6@h!|9L@V#qx%7Xf z=u-;+V=tD{@9^FVlB$)Fd9yh13J(`e)X;DK#Y8;LhI9$U?1=z46;YA>HS zkzr}G)R^;6oMKG8k}R0ZowTkH@%sL*ZZT^7D=j-aVUm20X-)1t@n749<&GSC)^hc7 zwerfQp_e}f9M=x|OO|zuj!Oezh28^UA08OMn0;4!+g&U*a)W+WfkRi-=efY>ACeLG z``5)$X<4}nm0MpbGrQ*6->B4P=~cKh^IRtI(n#kIN3AEG2eonnuUmxzc!FhMhSajlm-p zvs2`wJ^Fh%nmMEA#WIvrC*2^SEytZVo7Mu_7Cbw-PbjwvRV&eAeKov8xyrtei%5F@ z1ldv)CN3X8NWHU1a$4qGFH^UK$nZ6=w{ctSDOmbZ!7l>+k=nL zM_)MO?nJWWAMUq+uABULL2xGz!uP;WCLDL52*XZ`xACJ*^PU{WzmZb z3<3bxLt^Od6E2p!B3hB=--GAc_4&8;&J~zO=9GKj;P?xI<-5ghT^AR}Rea?iJzg3~ zJu!K3Q4rpLfA_sRz}P8E@ouIt>IW0PCL@h)0}3nVQ30lLHp(bI>glPYtUTVlBQ%Lt z&f1>bqQC!cCu)r61vic5m!o$$KevRM>cvVpKC#g*-89U`xnJ%Bu4nGMUuxS zzN2Z)FsghD97JQ#{l3@fn_pcZ!ZjY&xyG@fUI3TJt=2;{QCzZ8=|j z40I`W!XFdv`#JHpjXmdfm|H}E^u;0mV$0RlNDv^4@XC`emIv;dd6WN6Gz8#b3ze_Yq@phzM@ zgI~YRX^oWzi3A*OUT4}ACj>7y_OApp6 za>Y{LVT9AK>4mcyZ^;%^Tx*jLXUh%d*UW>c8Dgn-@U81osD=IbgcclBNmg@ zmg2O8ADPB!HLlt{0p@mqnA-?&G9GCY_+Pk$$YlGa@0?)ZsCLqckq{K9>FiJnghk1D_1O4;Bu{Kv7wiNQQt9P?t#{B*OJ+83 zz3h{#GHeL~gE3cESD_xBw5DcZiSfNV-PqcasIwe9zh&%ACsV#0?9}0GaNy8$;3+*( zX_KjF1d|92JY4Nv+}so?rg6x|n>YFGjg$OIW>$#AXLha6!YeY}c6$KV0nqv_JcGxHa)cwEdeH$99$vLvn_1B^&8aD8uYeVE z*IIzNXt!2@+-m`;NX^@iS0~yBP+K3wBWeySkmRV>W&i@p3)&l%Vo^;;pNUlpEb>5$G z@HlPnX;GM)H{E`|7H8yMI?2!~RwMBj|7qbQ(1djH0STTZuhx-^Fh?8*2|kEYGz=_C z?u<<6*Z^1-p1}h_+tigG+r}ATp8ljJN_E;4q@uvXih+<((Ej9dYG0_yZ}wk6DGF7& zU!06EN$QF(&^vl)+Q$pSROG3}*sq}mf+c;jY9hemaoBYVKjiSC4HD3X@rP3CBRh;M zxtV-=7)na2O9urg6T&ie6a`x`Z^~HoLq5&(U4%Yv|2})%ZwTuYi+I`x+>eQS{8015 z{}?dzkYo5xu}AAyOSIty!(&!Ycv*t?$%Pm|?!_xt&WW!sEi50b6`dH7dFCU<9ASE$ zD@!?87~>49(Ij@?wpV{;}Ru-SdP+g0WS_WaRZEP_Z!OeLOVfi9g<9y$6;wAH~Q}6zVlIi9?|FiIe%+m#D&Rwn? z7#+j=PP(4%`KGf{pf=9EPI|Ic>|=x7Y|-7h-ZDpVM*ssoTUh%DrOj+LDLTXtzRQ5> z{fAdsINLe2BIOksFVL|)%fU2!FK2(PwY~pg$Pd+D!v&)d=P8%%ScxppOscFZxqce0r;kJ@RI8pdsA*)(Bh4JTbAr*BNoZVE!B98rHj-(B2Sxo6wd?K1zQc`;14EXL7gb2=h@*qH$TS``|qj> z-)vfTZdb12(*ohzTdP;_S5LfJzU-f~87?Ofz7n7!z{;Z$J}8J5x)w@|9|o#8FRVmr zS5IomsSIktnyp_q#FgE_+m;o{h&*{u`(Kq$YQ`=KwVkPULtbhXJ>?^P<j58R~@V!v7V#UPlOI;!%p#q!W8R-;Vjef(?l4%2q7 zv@bU-7`mIgSp+dnyN#@ph)$FiJWSG7v^@$)BN>wA1RSx-Q`!Qm&yLDKT zbRg6n(5YRTKJqPFNcg%wnMl?|0RsQq_mJm>;o-qdmNJB(9@_SCdKl^75a>pxy;M9K zbp28do#@p}ws)kuckcSuOFXFw0%jkDsN*NiD7~FE6pp2sTN6ztl(O6@`6)Rl2iuIB z?vyr*HG{MKPtBiVPVCoevPe4Ldhb;M!hUZ^FpPJqwDT96Hwe8g8Xho0BW+6~6Q5Rp z=+IWIt5@78CsVxQ`@h^(+#gkB3-(d?spq6CW>&|yz2~%iA<^#O zTA>hPbw9lE2_2)}8;hwE)}Ez*V13ua_TwVymSJk;Doar3l-kkHWTomi!W#*)tH3L13h zpTA{HW>4q5NUiXM2qZdTw}gkSNB(}+)ex5ux0&hS8EMKMZX9n=7`WKkcbq{hld2p4 z-mGe5A~{E^L84Hp%kMG_Xnu{{9F;?7vLe&dJ1 ztgq1*1XjYEa0H z;}%Pe*C=%mzf|i&K~y3<83HWvLpnu+k@BKsk3S4{&Oh>K*utpK!7QPYlF6#t42s89Jc{a*OZw#;jI zVP52I$Ckjgkouk;Yhftq-n@Z#JVM~KgH(`w>j@}dyZ;X?pwQV&SdqN$+H1pU81g|> z42u>e9}0i2j8F*0A_){-XQKR?@(=oestT}y3Rs&|SV)FJSY*&M*#AHt*c9?U8FH*< zcnV?Iz%Clxmixp;XWj#eDd?#^di4tz?oY~cC^l5$&<5`U*?zimXWm$E*^Er_xk5}} zitx!0L4u7|Cx$zb$iN=rsq1}tS-pc{iheZz`9rsF8N$`aTH3nUFDamie1(yaA3tXD z9SgL7H!ioN1R}#@5=+qw>xK5##0eKFDO1?FnnL=Ej|E-Rb9TR&%yp%B-PtODM7~8h zxMf_t#?Cb1?Q)qbSyI@%@4Wp$rR_k+|M+R}*3m*XBVgG*LG(N}_Lo)LKn!7q!()7r z!tokT6_bTkeu$p~tn>7@4L%VX4IZP~(B5}m{`&qvP3+Q&cRrR*5tgrkCh*M+fiNM= zKHxy0$zFPzT+msKlZN=?dDH&~#&MvK&!(-LL6&;CjpS^I(;T1d8`&tYF3DOs;mv z`j(@H`ONUQpSq{UGI*@glILnNdp1M|V!6GM5fhSQ6X~%~wD`!mG^50|#wj!vLD0#z zI)OYWKCPRzst`sP+*qPzJ+o7h{W4XVsOp3diWlz$y<9pJh-TYP$h~ADhr$s`K(=u) z1{BK3yO}hk8BFBEW8(-eacXJ(WYYt*Cu)H&bzdjl?wjzgiG;)A2s{?2udY(W$?27F zhcm&P8_pb5CC+R@pPM;-2WeuX?ceUl1fNoJV!5Cz1NhntSXST;?1_%RgA>ggF0L0n zKUb&T+_LA*yoMUxJ(-}brmICYxiPqioDlH7Wfh^Yn+h%Dt9B&KN)P#3?%FCf~5+AIc>gxUo z4c%*7O?mR9tB`M@ywZz}Q0TLo&3Cl-_hs92FIJphBVrFbmmyNEv!r)F6z;h)vvQMj zYAhiz*ZriE2VNOIMzw799XN-G@Sql9L>dF=Jt9&4-(EaB9A<3iQ;IS?2m*M5mN&64 z%2s_D;NZ^ArwlqR(og4qW9Gb^(-E*{@AY-(92k4>yuA0{iY+1U-E0n;R#@w4+7Ptt z(96jIs5mTzgO{)uJQ^fcMM8h5?G>#-fo9+Yov@0tC?62`XByNcX^K{6$|#5J&+<`z zECTQR7bkJRhY61YWi5X}!FC3NJQO1UFvHsDXX;4t)9YLlH%O`G zHe2`T4Hvp=v^08C0|0VEv5#}rR;30-f-sC=nD@$n{*(4Ifo&JGNVg$&rtvJLTOg+n zso4P{Rq1i1w7So<)n$p{5_gkxMhX>2LV@z-{*zFtkUgU(Sb6(vK9|09Wy-GWnJB!D zJ(bX+>^2Vpp>Q1^%u|V{`FQZa^9R3-X2cXi{Q4{miAFwMGJ=(KiTN|vqrN2=_C_j+pcvkNzX6bDxE_$AN?meIA zz8YN6Pv@i^CIZ64!*5lVJOhyD$+oV?5DD#(Bgmv)iDq-V$lG2je>G5Tew9&v70O_= ziCJ5G!i#g$cDq03y2eMA3%bAe^7YXt<~z~(Wk`s5#(^Qce(gtGF1Gnz{Nu&0oTsQ=^)Z3GmG%mR+cY@Rk){VrTk+3*gOdGiTYucn>x4Pmw?-^ zf$8F?oBnIs6>q_pM^jlr6kQqYr_Gysri7cwqp`>G<|Kz`IV(Q}^Uq{Fmc}t93Vb{; z*O;L0_UDA5T6Y%468>o_x%g=d8o5oy+BgU)qq*BB#wNy8C1mA?q{i$_`yELi15D2G z=WNtXuAeL(dg@6qFnEmrybEBov24*za6f-EntB!@@QLhs6F1xV7H)+9OP@hOmwmI? z?BUpzrVOZO6BoaYuA=MLq~OmMJ+sv44Aq>uv$0^cxiG?CqYhFdJUpUzJP&)7?lXin z!=&2#Iq+~RKl!D_w2C9y=P0wzmmjY*LVnjpy_B%zaT@e5SAgqE$mG9FP0o*Xob4`6 z=Kbb$w+!aAu=8Co%}?q#RVL}*-qx`+w@^@hc0S-IOY)@$jGm0n4e5}?(EU%iMLprZ zO~<5c@|AtXYILrAS>%Nq1OnYZJWQJ0s>M6{^+$ZPW9zXM5|;FLl+w=@KnyRhnYM+ajS^D@Sm*d4L_tJ_VxF>p3Sj`Uem7Hdh{15D<7G8 z@)ZhbSZQyU*eip=1vF`hQ&>6BNY5GkJ{!96HhK*6c6?fih5`-iId3+e4A#;7Z{z=V z4Xe(VX};+2lYn<(Zz!hQ{kciZC!xWa>Ptx_KA3<%R8e=tUYuSK7{|wsQ5X?S-v{tO zQ|{e4eg2t&f8t{^@$o{v;1d6q5=Y0*3;kAc7nWym(XTsq*&H6W#JFApHv^}EinjjE zxGj>YD|IWe#IfODM2_;eeY`k|*L#+~PF@HJgR6i=x`l0np z9;)jfFYPByIGgn9#%q|S-n_26#Q4{Lmrfw@cV;er8?V>jJ6LSyfU;JaD0UP1%yhMt zZc92EdcxHRq^lOCkY+?42yi!c)n%CF{ZCM#vNsv5t+4zkPC3bD!A{Y55Pe1h;PdLf zGWh4a`^QIq3sz1+LBaF$^PzJi)%(QKW+R%UYn<_Zl9d-SdIPifkN0p@9{quvhdMR| zrL>twphX&bQlRfX3QP#%FBxFvr-0CRNHoLZ zR6dZ?ZX*Omg~k0#9{r!;+o=cMJ4(~C+(L57|@+H|=UYswQ81IAwmX?%j;&-?grywzhK)C??T3bfwi$Xb!ESn? zjh%=WZrBq3f%&;oN~uONns{Ju34=1JY_!WWxz2uhYdNe%%vzQ{;N=SHewdhw?7nj; z121wu2qTJZ!_i1&+?DrlrC!Q%jq`I3VBdKrXG?%7=cNm6?SAIM_%R&5V}9~<|DEZ@GcEj{9! zn?Pw2bljWttlz5keS4z$abjbiSS4@|f<_S9mSB4xzOt7LfcUe`^D9+sa|NSSmOdSg z{7Dv7B^nyZp9}SfGXYUX&fTfk&~4!(0U+XQklk?(C+`IaAU3AFBIlV z<29l2M7JD$6sQ1>u=j{*MaREwq9d$g1m;p>L8-c(U3o)FVJ^FwZwr=y@ChCB$ws%v zzcLM3=hcue_KTrQV9-96p<{RT;KyfKiJrqBc-jfB!C~*0x$Z z=bSKv>+m7qCQn;wA&Qpws(xcylltvVO=s=l%FO-l3Ei*VKJz_$JOrBpQ%J{1tzjXQ zP-B6WD?!esY8f0Y`EK+)D zobiX95#$LD;YWP^f_Uky$2!o3S(S$zIlTui!A77uwJ_lAIi=hP3kkj@cjh0aR91?K z0~cnKIV32*g!(xLXHx?n_paX>>vTsUV zR!-NG!mHWa8bC{ng@ls4+0mgm`hWzF`JoCFsH^`^{ey|jy}`f&yB!}^lt+;_$6QaU zA}o+BQbjP1y6jNXio%#9s;Vx>%G)m%V90S<)Bl&>*On|U3mxr@Pe_6^lH%>mP?@-- zoct0G7lE`!le-3JNCrFNv!%ryE4(G#pSb&C4h>Tn86_Fm4ZgEMWd^Fp!}8XK%%Uvq z$Lggc=l8L})GN#Ix&mfrwX*tL6*AvF$n-derRTkR3%v5QmtoB;Dj9Aw!L3xko zwfc?v3IgfMzx>}0A8Lh|9rlbL86sGiiLtTOBU_{OFCdqkX}bwZEQ}7)$p*Tvoa&!E)_$8lHB1?aa0A?$Ku!;j+de29HwbAi2?YxSvg=Cd zXleQFH6%vV+V89ov^#UX+^BY?lQTZMl-=t_UD&?~V&7e=K)_0E&`^hwwZ6}C1WC=L zOe0NQ4DGl2zSz9?ac?~9>rVY1N2T4l9T~gntVe#qgp2)|c1S()(3iEf1mC({M+J}& zFkV9jR5aV(aGibKB#~V~S;DW-_(8z}E=L8Jy&m9tf1gCOBjdcyz%x48;`D0H z92r>>JmXR1tXqURzu5`l?#mM#*Poy7kQ`sB)1Z}}h78!1i!Jjf9^>Cam)(xBRo6go zL}~R&2?R~^ZJDbuTXy1Ir@FmiMps+0YJVS@UBmrhvPppYtj3^D9LWLOw%R!8Ml;!> zqa$HdL~CT>uDu*X%(uj-7(|lxgrA)x!$CRJ)&<1J73mx4t6jfNo7OAly3iIl-oqWf zw$6`>u~>1faQm;-UOdIxFg4DJpZZD5!w)lgpIHhO6KMgztF7-9?T?7Fe{}mj)__WbsKqNSioTI2=zEHJar#W(#s zSj^1{-Pi0Co0!;ctsXKnsRao)Bc6Ehz4V_~zr(4j>gg5PTApu>TT_2Iz|)+=eVS&W z7lj4{!FUbo>D}oLGxITxYa6v3NCPq3w_Tnz$p2oNa*WM;ta_rIFDmVYxFg&BOe ze}3}25dl(sOqmS&(t*$(D)Ye)fEgfs`d)yWd17w$oDf`l{Xl4_#BCzqI@&e?5)V!5 zCAs{#>A`-nZ$v44teWfa0oNbe4_s*YEw+hgT#_x0?@7q3)8_c5eS4ey>us{k2k)n= zB&a!i2s7)7ot2KQEv)EY#FJL=58_E7+gD8~lg&_GIe1o_yS43nCOIPkh2SRv<+?KH zjJ?55fS)_Dd-wSl&00h~sf?!i4TpLF!?;=JYUq251OU zDsuNCbbTdA^rw;S@US}MyZ^uFrT=JSV*d%>560%Y*vT%a(0X5{qpyUTX#*-_^(O$s z!Cv;F}5&h}#zq@+0 zGvkcj`=rYdo0`G?#SRknhVilotD`{L*g+&wV0c4$iK%0~1-2rBhK$ax4KImF{z$t=#*IyuRi7eErJh49 zea;{F;M?_;qBXo0vI7@ER#;-a8wJP5^1@q(aKZ;`r#HpD5?XAIb)(yy$lc|Z@1tH- z+-KKbQ4US@gs4Nb9bw1g*I6Lo+u++|pNz(~Fo8SZ;PfQdZXW69H865zP}${mldloT z@!_EhQR*;MbRCVxa=Xrbn>uK1Nm8eT@AD+*nJ@w zRS^&_s%nDJ{PR!;;s;tTE5F_&KZvSyB@C`+CE2*|+g+&8#D7z^#%`~nWitj7I9Kn5 z&6ek|tK{iyVcJ}Q%rKr2spD$5#MHQnDXP2FsS@dsEiNCG-PO1vj6q?d z$Bzo*T)XoC0-VjbGZ+%@;?gJEi*8>tDg zbF{=MUrpAc+w+%BChM^tldF?K%+cfC&AM*;Ww$7|J>_awh}&&_FY4pv#m=Mss&5o$ z$sK3ul;dN}n%=qt@Z+$>GBQ0+JZNmKT)D_w9;W(Y2+^BZ>K^7)#+U;;dOjJXzOXOQs7mh_)D%L@=6!r&r zMle>WStp1ma=hJ_WATc>^3GE*PD)PgPAe=lLStB($ec>l^A;qvojLnS*cZKBN<6;H zKFo6n=h*SZiQ-??=y(uJOLU++$E?xze~VBXsXg|ib{eFN$~k;l8@dp=?d98(0-N!; z3>DsiLLJnuk)V>kDBu7$%fzu0s|7FUpp7rlNpN68H@NuwgTwTf4)F0n)L8*%A_HSY z>n~$F{lFFx<1y}ub%K(=hYj?vd_EDJRw{0rsaGk|(Z6SpBc`pgcZ+rhsRsG{cI*Ax zyoTNV0O>;xYp-}YJ~7wG{souzs<@83>DZim(=Q2pnK`d$WR8mpwNlcJ)K_PoWYhi- z?_J=OkOb%rUt;0nCIk4Que190P>U+a^4M|M0Td=Exwaa<^OfTm4j7G|emh_NOmgCl zt+>8rK&?`=}-Ue`FzHvNsHh2}Kr09XBax2=LN89vFWc+c|!;5yU2vx8A2x z)bov46gkc{kOC`?Xxj6nlK?sfer~9oZ=B6XE^OtWefI1eeAr2mGJN%5aX1S7K=eue z*9W&aeW``33L}mfRvFI%OIhZ0chb?ZUDT|`H6X+@s=_-$nE6nn%>yoETrsQE@xLiV z41aHhDxhmBMaTmzTuG$f=I9H$+n-2IbBOqBu=j7^J>-%Qea~498Y5T~)@2+L+U%YR zk`dHR=MxJQOM3-HMjVqe`JD!*s1`4KwFB89{2R}jXF6i>lZ?0(fGerS0iz4ewz|)} zMvDh&YKeY3iez|AuigvZ3*&kzuM*BGdT#YhwzXP7^Jozc4wSo0f~^+Q@sh&qfc)d~ zc!&%pC*~5Lh7TeoUCC;Z28>~QT317=t_i2bO(empaOz?$?+rR}UGbL;*|r2o>Md@7kp=@f!?rY}LPg#u9zPs8 z*fmiy?e=LoMP33n)1H>~4^YY>;g~ivwBrr&?pObJ_$JKprtyd6j$cyo?|xP60UvLU zbgVH%D#C=Bv%{{`c?HAf-JU=k9CIb;S}>TV?n40vrn^TL0pd@euh{;dyoQxyy$-6t z2)HZ^18Y0Mc76R~rKd`cg>?D7+q%&XJuIvUs{2U&n@Xj=AK$BU{znKZ?c%>;GSSL! zu$#Yp8ecOyOcW)jH$qE)MK75+{3I>F2;AlGjEOOC1{r!7BDS}iDS4m zi8~jEpdS$Zo%v{|lu+4X-quw}uMXJ&b5W_MquiKoQ(*~yX|_vbNH-(i4A`nhe0nFB ztUkzvz0b;Qvc*HJ_D($lX4pa;&_tZ#@TM!HF;~v4-aoQYoQ7Pr(yUD%jSIyvm<^2X zhrI2y3pYK*WMu$3-ZyrrS?eq=SEwB`p$@Bj%KGz52t;9;=f!3nN|yL)hl0Kwg1aCZp=cXvr}cO7JJ-g|#+ z)!Ti$RTMQf)35a^B67cnIe2u>FGJcEn~{J&)+>;-{7 zfTYAkR5XH)SKYm2Z=d_m=ToRs#IReXkz-%y5X}Y^YW~z+E&$J{@y{!L%&AjVVJY~i z`NsOK&Yoddc?B~xEgP+Od}=~->BoYF)PQWRf#dGf&NkrA{=li z|L7@-l$Fl&J9~Tk@OMDx(J3t>e=4QjVX%*ZP4Rm_owPY6oTT|QsSu(Jf8|;Xs$3fv zSLB0Q$?Ya{Kb`|7G_(`=W zXTwmeGI`QmWRP!S)=eBazS-}FaHN|;36uioC*>gwv+PUd@D?={$Un!S2( zBMQ#+Tdhq@O>&qpuK7s#@Y^$g&N$)VsHhCekc)*mOp2^(zMr`n5i6aHjWyUA+8#F6 zWK=0pdR?v-v1qk(u{{i$fqp~d;X!b&H#bgJeOGBVScjwQ!i3488I@R9=TxE78ZNFH zq0_*GX1%zVJD;bF6LeVeeo=E>}}3y<8wxftMce{0Al!%@H|J8T|se;j&z;X2L(L(a6>_uyC4@u>#3Ab*y7h0!R*|5 zx|Qw1oNy7DRr;j!Qg;OJ1SGTA>a2*}@#@{C=sYRaL0OHU@t%xN6&SzQtr?J%VYs)B z_+XXKvx1zQ7Y@Cn*%aFns3j55qKQ9SI9Q z9?NTgX7u!$@~+$7*#=c_pB2x--$O2i=ur|mx6@JPi)9$^A!sQn^rhTn+hX0bhStlm zn?4C*@2Q@eBZdl6D*2tEt+D28j>iXFUS3{i6WMS~nctN6JWj{0z}+dYEM|&D+Ls8& zv{@Zb;z}lJyuH28&(B3vRZsRCQuXG8wlq41i_0I|8hChYRyIxpJ4Cz7ztrlTzm56|7PG9cw_Xy}|XzN>k;?OeF`C394$eVVydZN3}layjN5n8Y5W zebbjyqLekw_vQgzUG`;$QkjeYa&hO}-dNoudRcn|_LpGxR%lI&R5YNFL^}UyLQ3n~ zE#wwFyf;htgHZL&<6(*-CGomw?at`kU4SbA+c6Dv4y9?^nt!NWua~JO<>_7ShieNW zRs!@#$NP7P0e&schjE#iu}1GR~~HpZa$`b=DT+1=qztI3Y5F9M!qddCkD zqWvtB*I^Pf#TFaQq@xQ}PYUovYCU=3Xo_SR@z&d(@;cnu&w;%Skh*wu?E; zqvdmIgr2qCdpM?x8C-o&S!6OBKIt^PwR9oRC?3O`9WBgd64ZvLB~ong(Y$XqbON^j znGYkTK55?GjnVH-DTl!XLVkg%A!kZT=U?>}p+8Brk$YHuYOrvW08L&U9&)QTOZ8u~ z8{v!TVZSzlfh=amX=pcDbG(Y7!8n(#f)gYwL*1uBAxM4lXDER|gAs3vW9sb6u1J?e zCX;cX@@Qs8<@0u_3CSL(Lx>DM1{BI4`B&yDPf`Lq|QIOB|q=luZ zv|~CnJd>%CmPdbfeZ2-{Td%(|e(`A(vxc4HU~;#NMU#ZML9ekReaYi?Jw`S&c>3(F zA~QFTh%SaAun7p>_Ip;!>+4?adf{EtA7@$(r>hsS(*C6pTKV}abnO<4lvKXGQ0F_A zp`oG1!=(K4LZ=@b0^+-NBW#w@^?}AT{KtK^Fu$-j3`ImimP)YaA!O^66CZv__JiSGz(~Ye zYEgyjY#ri)|9(y}DHCtaE(IcnNA!5o$$)Z9jIi+U&^5B%g-WVR2}gm!sI|_Q7+bX) zYPYT;kXLODkBovI;F+ zZ3{oM2fntpb{Tm>raS~fU3ev1p{F2Tq+F%Jm|sv}#7Wp}$2;3mAzvw3+311#Y1P<321-dxdonm9xox!Io6aJEt_9LAP&Tzo&MrX&ZsIavXXY=)5eblpxOzhgh8 zC^pBAeyM}kJKZn)c8;n$3Ud_k#1xSubhH}xfJ)e1Vb8H zFCsi#&k9q&G+lIM4juMA1@C$=D!dmF5h?pNxw_wSUq|3vt5@G7YGP~wWT>XSYf3$Z zUuNtEgQ)drFZmnadY6v-J(4RNl3E!+CDgd$s*@J=O-(Ivndvq!8H1Hamedpy@nN!# zE{WHfZTL~}c({n|4#}8$)`;KLN||bo9-Crulg;<%qKd!K9mv8iq&nTG1#Vtj1cxm# zt2f)b6Y-pQEQN8|%3q$4?YL#7LLl>gz8V6y?-+QVUxkx?Z?ImDGd3IXmOHRJIO38W zH05>QA;esp54}8r<3uYS410Cxr0j7MV&qv4r4P@k@O}sjv~)VWqV10%%kz19aN3)I z+pK?;=~gUM$@~~AQ>6BHQzt`^U&!RJ;TUvvFdH#1e8({^u*^N@`-Ei6s4;Nm8XX|R zLe4j7WCw%dG=|WXxXXIGq|=}VTtJVU#n#Uot&J^v@Z8_i~u=rW%+m7FzK_C z_f$;XRj}35vvDT--_1>1s7-WMW%0_G?w34+r&Kde!l~)$^75{-ua2(OKQc4JQ1_fP zab90}SW=XPw+3J0p&qg>=;|RHWuX^h8)@ZjAh_gt;dl`lu`oB?g~zjWx>> zMh57o&!2YJGqT6070#cRSI^7lNx%CcAWc2`(G+Wz6TYog%c)dK!q5Ag-1z~u_HURc zE++M;cF?enuCcMP!-Ez3a1x8Dh?PqQ2YZczfPz+2Jf@I#(18+?tE($s@8iJg?cOKP zXFkhyAB?3-XC0^HMH=)3ziqz~EILe7^cP0THRQihM~xBAFVeg|prodXg|PqG*+KmL zl7xxdlsRY&&YS$Yd(x6$p!Lb)LNha~^l^xr0KFI*F5#!3YilrG!*!{}YW4t!UV&Ia zrN>ZXHINe2P#Yi(u)oI~YyJY+$g!P8KP{oNM&f4si}LXPd_9UBVH5kOI;J;E%v&?j zY_#FsI%_co^6lSyt@70-K4_Od#3y=%`~DRcnQY0{WG`OK7p}zi-KIatzv`7ZLPgym`~+p#&@b zu{pPn+e0Qk_uKzcJphJ&?T*b+$Mm<=cSdG)sK-l!Q=Yv$ad7%KKu)qbwhn0=Q6a(kIuHv(MV0E@kRAYo z$=|0o0BZi)_GNawA2NXju(7kVvl2$73;SjU2EGyy6f|FL;mRac{q{0zFSP8!*xxHn z6d_KvuS5^`20HnBUGt_F4Zqp$4?Lvyy|^N4+Gu`Zp~Lf|yL^^VR7eq$^->)jkpA3K zF;&{SUa_iN7g_WaO&*-qD8yIcm!5{76ra95ntYC`b}N z7^quWq?jjBXSMjQd@Lwx{c&))^wx1-dac#YA<=tO*99m}05GSdp#hYSEkzc%b?&Ui z`Q>GvYN?u;GxKksC!|R-GTOSM6@gJC`Ze-qtX#v3%M$7JF*)SHmd(+*npcP9wXLoE zgi*oDrWNOl4Mt@z5119&j?5D_riI8{Ft?iLq3ncJXYrOxu4&lvcMox`l#fA@dqhF zpOxnlu)v|=V4#R3PycCpR|iBqZG;Chxma_2t8F93XG_l{XQHo$N``c5Qpn^(WVS4) z0*hvEWXmrxSAi>sKDu%rfSfdW6X*i5nzzUL%DU==x&{V!5z)xS3V!etmHU-U8b?W; z{B~Ro#LA5u$cUPAOy=u-_~b$!h~?5;%LXIfVOpW548(G+Sy@bW{!(Gg#Ye>%gw@^0 zM#Qv$F~}u72f}mEU=GCMj$TYAN!CP>CJ!7KWh3ISvoZi4V}ShM2sj&ZCT)~$B$y)M zKL~t$ys`R>Ol}jfvRDi^VTOaVeRI_vsY%egVE+~!EvyR3)Tgbz&9-MQ2N+X*w@ zcj}d&fiM642BhzEZrVVsM#U2=75{($)hJ667`*Vwg57XPzam;50lxrz{2+3Z0Z8`q zS2Pxo&${EgL%5b$Ncq#5#Y}K%fo36Bp!?FT2ss#n1hm-Xa>JlO0G^8 zJ2Z%pV~|bFN5EwU=5-<4r?7Ys>!093vp2)5%kt71sN_TSZ!8n54My4fqCN0v|FC)I z;z_(^WB$m%xQcA?r?-(o;R5Db1TvdkioiR5x;tKuC$P)QSYLo2H(?Ib^;hsFW7cic zh*M=$Y=jLxr*c_JTWf(sa}&Pk1Mgwztru_i#R2IJrOkk=plJQ$yt-;aR2nzCeix|0 zYE9{#YS|C!AUXO&dSr599DY{oCq{~zXPxF$XMzrzuu}DLdp>*LwH}-JYk|TvI1er& znMwc#8YTq1?pTF{gd#}=Qx_IAU>4&-Fh>`HpE)%0w-c)T zG0(*^zJ#ZS|Fu7+M!-ipx>&k+8pjirchog7V?wzfPK0stX-Cs@m0AifK2PRb^T^+^ zY2-mz@7{IqFr~qTMn#!)eTNH_LpM4(bo|!nLvp#Etg!IAiun3)p@hPyFPhv#it}); zaa;lhS$}!v z$&-wszSm?2D3+MlA+RB)v7V)L*1~?K#o2OdYRc_&_4e%S?6Z9n52uh&rX;ERtM336 zF9@|4v!f}`v_ui{NGOI^i268GG6d4<^kT@Sy^f8QO%X^jPzUi?LK^nB`zLMVS>50e zn;galW+@LXW$ng{+X=SZ2A(%#gdJXeTsLr4LLR?ppei;&sFLvn5}<>!LlJDZfwx>dxf-~W3-X0tdzrR4nGv9RO^5DP4 z!wW7gWst~=rBg24Dw+R5E=s`W^-6b5cSUE{hVY*2BSaisF(Nma+nEp$k^&%b0ctH8 zuG^ohPkVCya{41p6Q8yU1LWi51H89ye0OMy^0dEx{W?fyGeSm2-q_XiG3mbCjF`3e z7BQwoJzH)#l0V?4#{N#!S1Pdr^*g8cvUf61x$f0I>x0*?wdpMPjilmyRwpt zS$D59(!@6)tW!VB0)qU2m|v{lg1AP7*vD9o`H-5JZ!qC5N$XeS>=8n?WFQE$B( zsOfpzuV38f`VHT$o1*)EYns4+&zGP3v<^v2aA54y?zqBQi8B`XQ<2!|RdsbeUJpfg z$>DKCk}Vf{Vfy96Em5Cqn9O2)lafsg4p;@5YO4F{=7M>vtq4w(i>Knpt$XsjpMcus zi^11cMTkAYoXIWzfIq`4AuCp|CORfJOMOKpKY7)-4X7-5C5)04P$8*z$zbF^lqj;~ zfp&4PfbTy^{fAak`T>F`AJ>HgMyB~sIsMP(1%%1C24(jB)m^~FiY2ch;_ofB`{41x zIQa{h>jc6!UOAU}Nn!}R{e)yZHpD%W=#4C!$J?XDig}+Lwvgv5_4<2y0wANQ)mDq_ zu=pRuWaLS(>r3p+UD#%;pruJ-G(kMQ>PwvlDNT3OQNbZs;;g##G_12SKs!>H9TQgVz*$Qc@j zP4jMptc#`d8GL;(i-L~4H7@)@3FIg^I5=W{m&o(YzS^g>qA^N^gP7lk~i$K7R%Qp^qta*H!X_3;rogPEJmxBp)WYR$jHmY*l8T2IMk+ zE{R4c)V*-FQcrk&o`+we2|z0KYNvfAr8wn@gSxBZ<$0mXCf|ZUEL`Oh0nFN}>ON~f zU~z{N=PjFYUfizx^0%+1kx;Bm0;OxRu&5#*55H@dcSWo z8k^ifQIo5eDVbMClkpfNe4qZ~LgJc!cwkbBJR48k)ULPsu8GLr+tMYcL-ay8o^z-!pe-~i2ij%n{2Nc5NG$|8ro$Kkb|q#O`g za!;zjVk>DX)(k2mCp7Y*S_j-@c8Nx;T!DBkivMTu(W|e952XqOeF+UOpv?U2#6ijG z7)$dOg&g%`$kzW5wDi|qz_Nf$_sQ|D0hA{Gm6~ZpK^TTpyQqzq)v{K<>FLP#9Klis@q2eFl6^5m3SM3UfXeOz;w$KRwKWh+L8T^oJUiIjY-M5xO};NaUf$j)n3$$O zsNU>n@Hzk8_!R&u(`ija!KSvmJ2T*O+WQQaj$AuwF=5;J)F+B5C>DZ20yJrCJV0gOXv}SRS+cA=-w42DcRSkvH}d*E zdx5%pdX5i6LPDl%Ok^@WPnr$AGr=n;LtFA>`1tr9S9>By4@a#QOLgWwo~Io^JEU;x z)0~ZRi~`S6wQ|kgvp+p3=;*ID7OIUn&8GqZQ=kSiN_f)w;-glnJDSexNWgB4T&iAK zXg*!|PjLS=0G`zI+iI(;Waa!cuvAFhHsC6>o54|pT*(}!7!(qbAK2OBZjP5PS6$bu z%_cv7|2XWVGE>v;acwx2F9k@{*`E+Hyr$6-P&fsG-9zNk3I$udi6sP;qQtv8U^*;k-5Qy6Y@vQB_qHRTBb` zFl#)nEmTyuD_T$f0uypHmVy&C5&S zbB22B##76MM$|^DG<)9Kz^-zK{>^7ILS?q~j%2q|s8zvQw-fPH!)dMeKFmrFRM^Qb z*zZs|=EJq_Yy+WP)vC9cr4L4bZ|Lj0?tQCRZ=rf~v?!X!ZB4*whSlzQ+tB!V-S@>i z7=t7nK!lQ#l95b7euWgEo~NRu3@d=|*|P{O}8XE6&s??S*_(8?P7 zY3|4|4FB2&JKKJmRFo<>L%~6qAZ`(pFeS6Fd*(0 zwXB*;pr~C$wr#OSEyrztfIk2EZ3@CqoWLsmH9lTBz<$ZlS$xdju?C=B$MUPg=B4gA z++6nfE=HpaepkreM2^$x8i}v3FEC4|1J<19Kp4Byk)+?P#?KkiTQ<{zu(aA0&w+xLY&HDi=m+vrmUOXouua_HC-s@Zh5Xn%;yw}U)c5?SmQw7 zMvFBjfB{zXUuA{iTxso5rR(XM9MJmP?yl542kOp z%b`>JcZIzN_}gHnSn(C&tBbk6jiJaf`h6qsUm+pYTFx^A6o=T=biN5J|7bx*xQFZ0 z{uy)S^V>^TK33S(zUb*%`{dy@RZA8)^l^3Q4Mjr3<^xbOnZew&V4me_X+}PqeNm2) zF(hPU;}kIh@LfHa!rEqw&GXI@C&y^yHl5)U_U+oQAV3JRv9WPpSc)krz2kR18DSH4 z{g{?K1p=X}{5uPRV|lTH~`N7zCUXSz^64+w2e+cn}SYS>ajRu)-2jRG}0`$fa7 z9oCuY<5~HQiBVnh^La(8zxuE0JpfF~%R^q+R)a%lr`hGb;R|mBw#B`=e7jKs6-) zV)6gPi~kqARQX6U&oE>y2S!#T?!_1hcaVw40^_JxBPiR@0BDQnfbq|;`M*JEZ-=yy z&&pj~aR|V%YQQ4tJ(*PodMsc8>>~k&k~B+ct5DSa7`niw=|SHAU$p&?q2|9kY&t;) zT%f6y6%qzEx@uq{iWR6p{ZO)@7Us4W(Qko@Huwz@s!m^Ng*x5r&KUN)ci^kIJ*Cw` zFTKVJ#*RuQOxYjHw6;9=UJQ8m9rPIg-&0efgjfiL{w@p6*;FRkUn*cA$vPN-F|0Sb zLX-I8(i?#VC3NlGUs;d70<2xL}LXA5bKxcFe6=!-Uq_3z=Vxc2uIcG`bv^F{ zTn4o*4V#|eH<83_aaUJYHN}u-SAaZ#y>Tdw7(Yp?2OTRVQ8SEdAkeuYj2pKNTm;+}HK&hPO-LTWnm|nwb7DqrvTbz$2 zPJgU99WF%yq6S}Oc#rz;E3>SxElwxhQw1_%F)=iC3bfuHw>hRUnE*HB%ifgMsSb4s zc-}Aq1tf3b>sPOSoVofxJpqQ!0#>VE!6{>JZ*M_esRntOUSsC-)5G*YrR zIg|R+b~hDK7Nav>XFd&J6q5d0n|s^+Tw0kBa%$29*iGZqD z9mr8o7N?)@HZe(k?hHVr?pJSNx3w16^oh(mLj!LqQaH`u{exftW9z)znAUXH>O5;m zfCCW)jT|P$N{BlEyrr-iA=1&&f%5b7S&jP8gE7a2knrni!vRiX_kDh>fkK*Z{>!(B zPyYGEQ0XJT%h3mJ?m3@7N4r^NA#<%R$KC7e>m9lx({P-GZvcP7sQVIo@7AhQX(@gv z=jHn|aO9xh$7BCCud~)>l@A!br5^Q%58w5CpK!lt`97y~b$tc_^-ka8{pGg~Ltd7S zyN?u-k(8{gKXh7MVv~Ok08>!S3Ux>qD*1`;ZG}&h>Av$``-&OZyM4RF!`Uy6Cay@p zXGJUMd1Ji3<}I*r@T&%D(_9C*NMIO*j2(#!4S?JLZ_n~wDkCGKu*gVCDypD?fiEzJ z*K&>fi?uYRLYm^c%o(;X@+uzzziaEdxpF$k8%aPwU}E374Yt;LJ_35(VPG#k(oLu9#4lmx*j!WihB11GC`IWZr4wid&SKR`1*?Ty;4y}81*V9Fv zy~oKsNe%$_waVjY10P3{BBk9YB>`Cv`rq6Hm4jVYQ7~BAGd>~Qs{=Q}OAlA31PP+p z^E>=6-LZkB`1!!CYV}VG{hQsT?Z!!06awvn ziF@8UA+$VhmZ2V9CRds?m;+{$P&<<|k@5gMn z?58Andb84&L`BsgtlXn544Gx#5t*bYV{NXbq-nTxURM$vyYaf8BLDd_eqzHdJ_}K^ z&pBf8T#o}u*CN`yB2IB(g02LxFfjMb%8WH1m;CXs4o?^u4K(rj=kFwRj|LWa+qIOHZgT%) zO8&LcU=xwsI0m3PvB)(*M*S+BjhIK8#IOxUu0C!;JMSOWZg)cfu}mY)MF=6bBYyp& zFuHLmYCL=4O3UHnoUG9@r#RIg;G_ETfQ3o1F`? zs1R`ck>DMk8VxOOI;h3?E%!}*a8l8{waU8tf!Z=YT1%QM z#b$jCm0XOixp^rQQUWn4fDacGQmpY;QmrKnz_>5muq%rBrJ7cs?`^}a`nuP^-)Imc zl4t9yF0!J}H^;HQ94fP{pL%oKjG{obDrv{sbQj$?rd%#la`SaWTCpEe=q!R0@}^J< z4LE)n{a(B78?U*gG>INkoYTQVM`-PumbjjIr5>8rMoA_6@$1B+?{EDS#BkOw2 znBRjKE5!KAqjuEiBi4&jAverage Length hypothesis can be applied for meshing of edges composing your geometrical object. Definition of this hypothesis consists of setting the \b length of segments, which will split these -edges. The points on the edges generated by these segments will -represent nodes of your mesh. Later these nodes will be used for -meshing of the faces abutting to these edges. +edges, and the \b precision of rounding. The points on the edges +generated by these segments will represent nodes of your mesh. +Later these nodes will be used for meshing of the faces abutting to +these edges. + +The \b precision parameter is used to allow rounding a number of +segments, calculated from the edge length and average length of +segment, to the lower integer, if this value outstands from it in +bounds of the precision. Otherwise, the number of segments is rounded +to the higher integer. Use value 0.5 to provide rounding to the +nearest integer, 1.0 for the lower integer, 0.0 for the higher +integer. Default value is 1e-07. \image html image41.gif @@ -155,4 +164,4 @@ minimum and maximum value of this parameter. \image html image148.gif -*/ \ No newline at end of file +*/ diff --git a/idl/SMESH_BasicHypothesis.idl b/idl/SMESH_BasicHypothesis.idl index 9fe9370ee..5978aa052 100644 --- a/idl/SMESH_BasicHypothesis.idl +++ b/idl/SMESH_BasicHypothesis.idl @@ -43,13 +43,33 @@ module StdMeshers /*! * Sets parameter value */ - void SetLength(in double length) + void SetLength(in double length) + raises (SALOME::SALOME_Exception); + + /*! + * Sets parameter value + * + * Precision parameter is used to allow rounding a number of segments, + * calculated from the edge length and average length of segment, + * to the lower integer, if this value outstands from it in bounds of the precision. + * Otherwise, the number of segments is rounded to the higher integer. + * Use value 0.5 to provide rounding to the nearest integer, + * 1.0 for the lower integer, 0.0 for the higher integer. + * Default value is 1e-07. In old studies, restored from file, + * this value will be set to zero, what corresponds to the old behaviour. + */ + void SetPrecision(in double precision) raises (SALOME::SALOME_Exception); /*! * Returns parameter value */ double GetLength(); + + /*! + * Returns parameter value + */ + double GetPrecision(); }; /*! @@ -85,7 +105,7 @@ module StdMeshers /*! * Sets parameter value */ - void SetNumberOfSegments(in long segmentsNumber) + void SetNumberOfSegments(in long segmentsNumber) raises (SALOME::SALOME_Exception); /*! diff --git a/src/StdMeshers/StdMeshers_CompositeSegment_1D.cxx b/src/StdMeshers/StdMeshers_CompositeSegment_1D.cxx index c60cbee61..381caade1 100644 --- a/src/StdMeshers/StdMeshers_CompositeSegment_1D.cxx +++ b/src/StdMeshers/StdMeshers_CompositeSegment_1D.cxx @@ -119,7 +119,7 @@ namespace { return; for ( int iE = 0; iE < side.NbEdges(); ++iE ) { - // set listener and its data + // set listener and its data EventListenerData * listenerData = new EventListenerData(true); const TopoDS_Edge& edge = side.Edge( iE ); SMESH_subMesh * sm = side.GetMesh()->GetSubMesh( edge ); @@ -333,7 +333,7 @@ bool StdMeshers_CompositeSegment_1D::Compute(SMESH_Mesh & aMesh, auto_ptr< BRepAdaptor_CompCurve > C3d ( side->GetCurve3d() ); double f = C3d->FirstParameter(), l = C3d->LastParameter(); list< double > params; - if ( !computeInternalParameters ( *C3d, side->Length(), f, l, params, false )) + if ( !computeInternalParameters ( aMesh, *C3d, side->Length(), f, l, params, false )) return false; // Redistribute parameters near ends @@ -349,7 +349,7 @@ bool StdMeshers_CompositeSegment_1D::Compute(SMESH_Mesh & aMesh, const SMDS_MeshNode * nFirst = SMESH_Algo::VertexNode( VFirst, meshDS ); const SMDS_MeshNode * nLast = SMESH_Algo::VertexNode( VLast, meshDS ); - if (!nFirst) + if (!nFirst) return error(COMPERR_BAD_INPUT_MESH, TComm("No node on vertex ") <ShapeToIndex(VFirst)); if (!nLast) @@ -414,4 +414,3 @@ bool StdMeshers_CompositeSegment_1D::Compute(SMESH_Mesh & aMesh, return true; } - diff --git a/src/StdMeshers/StdMeshers_LocalLength.cxx b/src/StdMeshers/StdMeshers_LocalLength.cxx index 6fd6bbcf7..e3fe6e898 100644 --- a/src/StdMeshers/StdMeshers_LocalLength.cxx +++ b/src/StdMeshers/StdMeshers_LocalLength.cxx @@ -43,6 +43,7 @@ #include #include #include +#include using namespace std; @@ -56,6 +57,7 @@ StdMeshers_LocalLength::StdMeshers_LocalLength(int hypId, int studyId, SMESH_Gen :SMESH_Hypothesis(hypId, studyId, gen) { _length = 1.; + _precision = Precision::Confusion(); _name = "LocalLength"; _param_algo_dim = 1; // is used by SMESH_Regular_1D } @@ -78,12 +80,13 @@ StdMeshers_LocalLength::~StdMeshers_LocalLength() void StdMeshers_LocalLength::SetLength(double length) throw(SALOME_Exception) { - double oldLength = _length; - if (length <= 0) - throw SALOME_Exception(LOCALIZED("length must be positive")); - _length = length; - if (oldLength != _length) - NotifySubMeshesHypothesisModification(); + double oldLength = _length; + if (length <= 0) + throw SALOME_Exception(LOCALIZED("length must be positive")); + _length = length; + const double precision = 1e-7; + if (fabs(oldLength - _length) > precision) + NotifySubMeshesHypothesisModification(); } //============================================================================= @@ -94,7 +97,33 @@ void StdMeshers_LocalLength::SetLength(double length) throw(SALOME_Exception) double StdMeshers_LocalLength::GetLength() const { - return _length; + return _length; +} + +//============================================================================= +/*! + * + */ +//============================================================================= +void StdMeshers_LocalLength::SetPrecision (double thePrecision) throw(SALOME_Exception) +{ + double oldPrecision = _precision; + if (_precision < 0) + throw SALOME_Exception(LOCALIZED("precision cannot be negative")); + _precision = thePrecision; + const double precision = 1e-8; + if (fabs(oldPrecision - _precision) > precision) + NotifySubMeshesHypothesisModification(); +} + +//============================================================================= +/*! + * + */ +//============================================================================= +double StdMeshers_LocalLength::GetPrecision() const +{ + return _precision; } //============================================================================= @@ -105,7 +134,7 @@ double StdMeshers_LocalLength::GetLength() const ostream & StdMeshers_LocalLength::SaveTo(ostream & save) { - save << this->_length; + save << this->_length << " " << this->_precision; return save; } @@ -119,11 +148,23 @@ istream & StdMeshers_LocalLength::LoadFrom(istream & load) { bool isOK = true; double a; + isOK = (load >> a); if (isOK) this->_length = a; else load.clear(ios::badbit | load.rdstate()); + + isOK = (load >> a); + if (isOK) + this->_precision = a; + else + { + load.clear(ios::badbit | load.rdstate()); + // old format, without precision + _precision = 0.; + } + return load; } @@ -190,5 +231,7 @@ bool StdMeshers_LocalLength::SetParametersByMesh(const SMESH_Mesh* theMesh, if ( nbEdges ) _length /= nbEdges; + _precision = Precision::Confusion(); + return nbEdges; } diff --git a/src/StdMeshers/StdMeshers_LocalLength.hxx b/src/StdMeshers/StdMeshers_LocalLength.hxx index 79a85df7f..93fc49a6f 100644 --- a/src/StdMeshers/StdMeshers_LocalLength.hxx +++ b/src/StdMeshers/StdMeshers_LocalLength.hxx @@ -35,15 +35,17 @@ #include "SMESH_Hypothesis.hxx" #include "Utils_SALOME_Exception.hxx" -class STDMESHERS_EXPORT StdMeshers_LocalLength:public SMESH_Hypothesis +class STDMESHERS_EXPORT StdMeshers_LocalLength: public SMESH_Hypothesis { public: StdMeshers_LocalLength(int hypId, int studyId, SMESH_Gen * gen); virtual ~ StdMeshers_LocalLength(); void SetLength(double length) throw(SALOME_Exception); + void SetPrecision(double precision) throw(SALOME_Exception); double GetLength() const; + double GetPrecision() const; virtual std::ostream & SaveTo(std::ostream & save); virtual std::istream & LoadFrom(std::istream & load); @@ -60,6 +62,7 @@ class STDMESHERS_EXPORT StdMeshers_LocalLength:public SMESH_Hypothesis protected: double _length; + double _precision; }; #endif diff --git a/src/StdMeshers/StdMeshers_RadialPrism_3D.cxx b/src/StdMeshers/StdMeshers_RadialPrism_3D.cxx index cc1b4f89e..4bcd5519a 100644 --- a/src/StdMeshers/StdMeshers_RadialPrism_3D.cxx +++ b/src/StdMeshers/StdMeshers_RadialPrism_3D.cxx @@ -160,7 +160,7 @@ bool StdMeshers_RadialPrism_3D::Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& a // to delete helper at exit from Compute() std::auto_ptr helperDeleter( myHelper ); - // get 2 shells + // get 2 shells TopoDS_Solid solid = TopoDS::Solid( aShape ); TopoDS_Shell outerShell = BRepTools::OuterShell( solid ); TopoDS_Shape innerShell; @@ -336,7 +336,7 @@ public: BRepAdaptor_Curve C3D(edge); double f = C3D.FirstParameter(), l = C3D.LastParameter(); list< double > params; - if ( !StdMeshers_Regular_1D::computeInternalParameters( C3D, len, f, l, params, false )) + if ( !StdMeshers_Regular_1D::computeInternalParameters( aMesh, C3D, len, f, l, params, false )) return error("StdMeshers_Regular_1D failed to compute layers distribution"); positions.clear(); diff --git a/src/StdMeshers/StdMeshers_Regular_1D.cxx b/src/StdMeshers/StdMeshers_Regular_1D.cxx index 9d051ae1a..b2a507433 100644 --- a/src/StdMeshers/StdMeshers_Regular_1D.cxx +++ b/src/StdMeshers/StdMeshers_Regular_1D.cxx @@ -146,7 +146,9 @@ bool StdMeshers_Regular_1D::CheckHypothesis const StdMeshers_LocalLength * hyp = dynamic_cast (theHyp); ASSERT(hyp); - _value[ BEG_LENGTH_IND ] = _value[ END_LENGTH_IND ] = hyp->GetLength(); + //_value[ BEG_LENGTH_IND ] = _value[ END_LENGTH_IND ] = hyp->GetLength(); + _value[ BEG_LENGTH_IND ] = hyp->GetLength(); + _value[ END_LENGTH_IND ] = hyp->GetPrecision(); ASSERT( _value[ BEG_LENGTH_IND ] > 0 ); _hypType = LOCAL_LENGTH; aStatus = SMESH_Hypothesis::HYP_OK; @@ -224,7 +226,9 @@ bool StdMeshers_Regular_1D::CheckHypothesis StdMeshers_AutomaticLength * hyp = const_cast (dynamic_cast (theHyp)); ASSERT(hyp); - _value[ BEG_LENGTH_IND ] = _value[ END_LENGTH_IND ] = hyp->GetLength( &aMesh, aShape ); + //_value[ BEG_LENGTH_IND ] = _value[ END_LENGTH_IND ] = hyp->GetLength( &aMesh, aShape ); + _value[ BEG_LENGTH_IND ] = hyp->GetLength( &aMesh, aShape ); + _value[ END_LENGTH_IND ] = Precision::Confusion(); // ?? or set to zero? ASSERT( _value[ BEG_LENGTH_IND ] > 0 ); _hypType = LOCAL_LENGTH; aStatus = SMESH_Hypothesis::HYP_OK; @@ -243,7 +247,7 @@ static bool computeParamByFunc(Adaptor3d_Curve& C3d, double first, double last, // never do this way //OSD::SetSignal( true ); - if( nbSeg<=0 ) + if (nbSeg <= 0) return false; MESSAGE( "computeParamByFunc" ); @@ -251,12 +255,12 @@ static bool computeParamByFunc(Adaptor3d_Curve& C3d, double first, double last, int nbPnt = 1 + nbSeg; vector x(nbPnt, 0.); - if( !buildDistribution( func, 0.0, 1.0, nbSeg, x, 1E-4 ) ) + if (!buildDistribution(func, 0.0, 1.0, nbSeg, x, 1E-4)) return false; MESSAGE( "Points:\n" ); char buf[1024]; - for( int i=0; i<=nbSeg; i++ ) + for ( int i=0; i<=nbSeg; i++ ) { sprintf( buf, "%f\n", float(x[i] ) ); MESSAGE( buf ); @@ -524,7 +528,7 @@ void StdMeshers_Regular_1D::redistributeNearVertices (SMESH_Mesh & theM std::swap( algo._value[ BEG_LENGTH_IND ], algo._value[ END_LENGTH_IND ]); } list params; - if ( algo.computeInternalParameters( theC3d, L, from, to, params, false )) + if ( algo.computeInternalParameters( theMesh, theC3d, L, from, to, params, false )) { if ( isEnd1 ) params.reverse(); while ( 1 + nHalf-- ) @@ -547,12 +551,14 @@ void StdMeshers_Regular_1D::redistributeNearVertices (SMESH_Mesh & theM * */ //============================================================================= -bool StdMeshers_Regular_1D::computeInternalParameters(Adaptor3d_Curve& theC3d, +bool StdMeshers_Regular_1D::computeInternalParameters(SMESH_Mesh & theMesh, + Adaptor3d_Curve& theC3d, double theLength, double theFirstU, double theLastU, list & theParams, - const bool theReverse) + const bool theReverse, + bool theConsiderPropagation) { theParams.clear(); @@ -568,6 +574,38 @@ bool StdMeshers_Regular_1D::computeInternalParameters(Adaptor3d_Curve& theC3d, { // Local Length hypothesis double nbseg = ceil(theLength / _value[ BEG_LENGTH_IND ]); // integer sup + + // NPAL17873: + bool isFound = false; + if (theConsiderPropagation && !_mainEdge.IsNull()) // propagated from some other edge + { + // Advanced processing to assure equal number of segments in case of Propagation + SMESH_subMesh* sm = theMesh.GetSubMeshContaining(_mainEdge); + if (sm) { + bool computed = sm->IsMeshComputed(); + if (!computed) { + if (sm->GetComputeState() == SMESH_subMesh::READY_TO_COMPUTE) { + sm->ComputeStateEngine(SMESH_subMesh::COMPUTE); + computed = sm->IsMeshComputed(); + } + } + if (computed) { + SMESHDS_SubMesh* smds = sm->GetSubMeshDS(); + int nb_segments = smds->NbElements(); + if (nbseg - 1 <= nb_segments && nb_segments <= nbseg + 1) { + isFound = true; + nbseg = nb_segments; + } + } + } + } + if (!isFound) // not found by meshed edge in the propagation chain, use precision + { + double aPrecision = _value[ END_LENGTH_IND ]; + double nbseg_prec = ceil((theLength / _value[ BEG_LENGTH_IND ]) - aPrecision); + if (nbseg_prec == (nbseg - 1)) nbseg--; + } + if (nbseg <= 0) nbseg = 1; // degenerated edge eltSize = theLength / nbseg; @@ -736,14 +774,14 @@ bool StdMeshers_Regular_1D::computeInternalParameters(Adaptor3d_Curve& theC3d, */ //============================================================================= -bool StdMeshers_Regular_1D::Compute(SMESH_Mesh & aMesh, const TopoDS_Shape & aShape) +bool StdMeshers_Regular_1D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape & theShape) { if ( _hypType == NONE ) return false; - SMESHDS_Mesh * meshDS = aMesh.GetMeshDS(); + SMESHDS_Mesh * meshDS = theMesh.GetMeshDS(); - const TopoDS_Edge & EE = TopoDS::Edge(aShape); + const TopoDS_Edge & EE = TopoDS::Edge(theShape); TopoDS_Edge E = TopoDS::Edge(EE.Oriented(TopAbs_FORWARD)); int shapeID = meshDS->ShapeToIndex( E ); @@ -769,10 +807,10 @@ bool StdMeshers_Regular_1D::Compute(SMESH_Mesh & aMesh, const TopoDS_Shape & aSh BRepAdaptor_Curve C3d( E ); double length = EdgeLength( E ); - if ( ! computeInternalParameters( C3d, length, f, l, params, reversed )) { + if ( ! computeInternalParameters( theMesh, C3d, length, f, l, params, reversed, true )) { return false; } - redistributeNearVertices( aMesh, C3d, length, params, VFirst, VLast ); + redistributeNearVertices( theMesh, C3d, length, params, VFirst, VLast ); // edge extrema (indexes : 1 & NbPoints) already in SMDS (TopoDS_Vertex) // only internal nodes receive an edge position with param on curve diff --git a/src/StdMeshers/StdMeshers_Regular_1D.hxx b/src/StdMeshers/StdMeshers_Regular_1D.hxx index 269a034bc..4a22e7253 100644 --- a/src/StdMeshers/StdMeshers_Regular_1D.hxx +++ b/src/StdMeshers/StdMeshers_Regular_1D.hxx @@ -73,12 +73,14 @@ public: protected: - virtual bool computeInternalParameters (Adaptor3d_Curve & theC3d, - double theLength, - double theFirstU, - double theLastU, - std::list< double > & theParameters, - const bool theReverse); + virtual bool computeInternalParameters (SMESH_Mesh & theMesh, + Adaptor3d_Curve & theC3d, + double theLength, + double theFirstU, + double theLastU, + std::list & theParameters, + const bool theReverse, + bool theConsiderPropagation = false); virtual void redistributeNearVertices (SMESH_Mesh & theMesh, Adaptor3d_Curve & theC3d, @@ -101,12 +103,12 @@ protected: BEG_LENGTH_IND = 0, END_LENGTH_IND = 1, DEFLECTION_IND = 0 - }; + }; enum IValueIndex { NB_SEGMENTS_IND = 0, DISTR_TYPE_IND = 1, - CONV_MODE_IND = 2 + CONV_MODE_IND = 2 }; enum VValueIndex { diff --git a/src/StdMeshersGUI/StdMeshersGUI_StdHypothesisCreator.cxx b/src/StdMeshersGUI/StdMeshersGUI_StdHypothesisCreator.cxx index 13f87e560..bfe80bdf1 100644 --- a/src/StdMeshersGUI/StdMeshersGUI_StdHypothesisCreator.cxx +++ b/src/StdMeshersGUI/StdMeshersGUI_StdHypothesisCreator.cxx @@ -388,6 +388,7 @@ QString StdMeshersGUI_StdHypothesisCreator::storeParams() const StdMeshers::StdMeshers_LocalLength::_narrow( hypothesis() ); h->SetLength( params[0].myValue.toDouble() ); + h->SetPrecision( params[1].myValue.toDouble() ); } else if( hypType()=="SegmentLengthAroundVertex" ) { @@ -535,6 +536,9 @@ bool StdMeshersGUI_StdHypothesisCreator::stdParams( ListOfStdParams& p ) const item.myName = tr("SMESH_LOCAL_LENGTH_PARAM"); item.myValue = h->GetLength(); p.append( item ); + item.myName = tr("SMESH_LOCAL_LENGTH_PRECISION"); + item.myValue = h->GetPrecision(); + p.append( item ); } else if( hypType()=="SegmentLengthAroundVertex" ) { @@ -704,12 +708,15 @@ bool StdMeshersGUI_StdHypothesisCreator::stdParams( ListOfStdParams& p ) const */ //================================================================================ -void StdMeshersGUI_StdHypothesisCreator::attuneStdWidget( QWidget* w, const int ) const +void StdMeshersGUI_StdHypothesisCreator::attuneStdWidget (QWidget* w, const int param) const { SMESHGUI_SpinBox* sb = w->inherits( "SMESHGUI_SpinBox" ) ? ( SMESHGUI_SpinBox* )w : 0; if( hypType()=="LocalLength" && sb ) { - sb->RangeStepAndValidator( VALUE_SMALL, VALUE_MAX, 1.0, 6 ); + if (param == 0) // Length + sb->RangeStepAndValidator( VALUE_SMALL, VALUE_MAX, 1.0, 6 ); + else // Precision + sb->RangeStepAndValidator( 0.0, 1.0, 0.05, 6 ); } else if( hypType()=="Arithmetic1D" && sb ) { diff --git a/src/StdMeshersGUI/StdMeshers_msg_en.po b/src/StdMeshersGUI/StdMeshers_msg_en.po index 922119e14..d8cfd1435 100644 --- a/src/StdMeshersGUI/StdMeshers_msg_en.po +++ b/src/StdMeshersGUI/StdMeshers_msg_en.po @@ -37,6 +37,9 @@ msgstr "Average Length" msgid "SMESH_LOCAL_LENGTH_PARAM" msgstr "Length" +msgid "SMESH_LOCAL_LENGTH_PRECISION" +msgstr "Precision" + msgid "SMESH_LOCAL_LENGTH_TITLE" msgstr "Hypothesis Construction" diff --git a/src/StdMeshers_I/StdMeshers_LocalLength_i.cxx b/src/StdMeshers_I/StdMeshers_LocalLength_i.cxx index 4190c1b8c..f360400ca 100644 --- a/src/StdMeshers_I/StdMeshers_LocalLength_i.cxx +++ b/src/StdMeshers_I/StdMeshers_LocalLength_i.cxx @@ -48,15 +48,15 @@ using namespace std; //============================================================================= StdMeshers_LocalLength_i::StdMeshers_LocalLength_i( PortableServer::POA_ptr thePOA, - int theStudyId, - ::SMESH_Gen* theGenImpl ) - : SALOME::GenericObj_i( thePOA ), + int theStudyId, + ::SMESH_Gen* theGenImpl ) + : SALOME::GenericObj_i( thePOA ), SMESH_Hypothesis_i( thePOA ) { MESSAGE( "StdMeshers_LocalLength_i::StdMeshers_LocalLength_i" ); myBaseImpl = new ::StdMeshers_LocalLength( theGenImpl->GetANewId(), - theStudyId, - theGenImpl ); + theStudyId, + theGenImpl ); } //============================================================================= @@ -79,7 +79,6 @@ StdMeshers_LocalLength_i::~StdMeshers_LocalLength_i() * Set length */ //============================================================================= - void StdMeshers_LocalLength_i::SetLength( CORBA::Double theLength ) throw ( SALOME::SALOME_Exception ) { @@ -97,6 +96,30 @@ void StdMeshers_LocalLength_i::SetLength( CORBA::Double theLength ) SMESH::TPythonDump() << _this() << ".SetLength( " << theLength << " )"; } +//============================================================================= +/*! + * StdMeshers_LocalLength_i::SetPrecision + * + * Set length + */ +//============================================================================= +void StdMeshers_LocalLength_i::SetPrecision( CORBA::Double thePrecision ) + throw ( SALOME::SALOME_Exception ) +{ + MESSAGE( "StdMeshers_LocalLength_i::SetPrecision" ); + ASSERT( myBaseImpl ); + try { + this->GetImpl()->SetPrecision( thePrecision ); + } + catch ( SALOME_Exception& S_ex ) { + THROW_SALOME_CORBA_EXCEPTION( S_ex.what(), + SALOME::BAD_PARAM ); + } + + // Update Python script + SMESH::TPythonDump() << _this() << ".SetPrecision( " << thePrecision << " )"; +} + //============================================================================= /*! * StdMeshers_LocalLength_i::GetLength @@ -104,7 +127,6 @@ void StdMeshers_LocalLength_i::SetLength( CORBA::Double theLength ) * Get length */ //============================================================================= - CORBA::Double StdMeshers_LocalLength_i::GetLength() { MESSAGE( "StdMeshers_LocalLength_i::GetLength" ); @@ -112,6 +134,20 @@ CORBA::Double StdMeshers_LocalLength_i::GetLength() return this->GetImpl()->GetLength(); } +//============================================================================= +/*! + * StdMeshers_LocalLength_i::GetPrecision + * + * Get precision + */ +//============================================================================= +CORBA::Double StdMeshers_LocalLength_i::GetPrecision() +{ + MESSAGE( "StdMeshers_LocalLength_i::GetPrecision" ); + ASSERT( myBaseImpl ); + return this->GetImpl()->GetPrecision(); +} + //============================================================================= /*! * StdMeshers_LocalLength_i::GetImpl @@ -119,7 +155,6 @@ CORBA::Double StdMeshers_LocalLength_i::GetLength() * Get implementation */ //============================================================================= - ::StdMeshers_LocalLength* StdMeshers_LocalLength_i::GetImpl() { MESSAGE( "StdMeshers_LocalLength_i::GetImpl" ); @@ -139,4 +174,3 @@ CORBA::Boolean StdMeshers_LocalLength_i::IsDimSupported( SMESH::Dimension type ) { return type == SMESH::DIM_1D; } - diff --git a/src/StdMeshers_I/StdMeshers_LocalLength_i.hxx b/src/StdMeshers_I/StdMeshers_LocalLength_i.hxx index 62c73cbc8..a51997803 100644 --- a/src/StdMeshers_I/StdMeshers_LocalLength_i.hxx +++ b/src/StdMeshers_I/StdMeshers_LocalLength_i.hxx @@ -58,8 +58,14 @@ public: // Set length void SetLength( CORBA::Double theLength ) throw ( SALOME::SALOME_Exception ); + // Set precision + void SetPrecision( CORBA::Double thePrecision ) + throw ( SALOME::SALOME_Exception ); + // Get length CORBA::Double GetLength(); + // Get precision + CORBA::Double GetPrecision(); // Get implementation ::StdMeshers_LocalLength* GetImpl(); -- 2.39.2