yade-dev team mailing list archive
-
yade-dev team
-
Mailing list archive
-
Message #07842
[Branch ~yade-dev/yade/trunk] Rev 2908: 1. Changes in ConcretePM.cpp
------------------------------------------------------------
revno: 2908
committer: Jan Stransky <honzik.stransky@xxxxxxxxx>
branch nick: yade
timestamp: Sat 2011-08-27 19:40:30 +0200
message:
1. Changes in ConcretePM.cpp
2. added yade.ymport.unv function
3. added yade.export.VTKExporter class
4. added method outer to Vector3 in Python
added:
scripts/test/shell.unv
scripts/test/unvRead.py
scripts/test/vtkExporter.py
modified:
lib/base/Math.hpp
pkg/dem/ConcretePM.cpp
pkg/dem/ConcretePM.hpp
pkg/dem/VTKRecorder.cpp
pkg/dem/VTKRecorder.hpp
py/export.py
py/mathWrap/miniEigen.cpp
py/ymport.py
--
lp:yade
https://code.launchpad.net/~yade-dev/yade/trunk
Your team Yade developers is subscribed to branch lp:yade.
To unsubscribe from this branch go to https://code.launchpad.net/~yade-dev/yade/trunk/+edit-subscription
=== modified file 'lib/base/Math.hpp'
--- lib/base/Math.hpp 2011-03-24 21:38:33 +0000
+++ lib/base/Math.hpp 2011-08-27 17:40:30 +0000
@@ -177,6 +177,10 @@
VECTOR6_TEMPLATE(Scalar) ret; ret<<m(0,0),m(1,1),m(2,2),k*.5*(m(1,2)+m(2,1)),k*.5*(m(2,0)+m(0,2)),k*.5*(m(0,1)+m(1,0)); return ret;
}
+/* outer product m[i,j] = v1[i]*v2[j] */
+//Matrix3r outer(const Vector3r& v1, const Vector3r& v2) { return v1*v2.transpose(); }
+
+
__attribute__((unused))
const Real NaN(std::numeric_limits<Real>::signaling_NaN());
=== modified file 'pkg/dem/ConcretePM.cpp'
--- pkg/dem/ConcretePM.cpp 2011-04-08 06:31:21 +0000
+++ pkg/dem/ConcretePM.cpp 2011-08-27 17:40:30 +0000
@@ -25,23 +25,31 @@
// check unassigned values
assert(!isnan(mat1->G_over_E));
+ assert(!isnan(mat2->G_over_E));
if(!mat1->neverDamage) {
assert(!isnan(mat1->sigmaT));
assert(!isnan(mat1->epsCrackOnset));
- assert(!isnan(mat1->relDuctility));
+ assert(!isnan(mat1->crackOpening));
assert(!isnan(mat1->G_over_E));
}
+ if(!mat2->neverDamage) {
+ assert(!isnan(mat2->sigmaT));
+ assert(!isnan(mat2->epsCrackOnset));
+ assert(!isnan(mat2->crackOpening));
+ assert(!isnan(mat2->G_over_E));
+ }
+ cpmPhys->damLaw = mat1->damLaw;
// bodies sharing the same material; no averages necessary
if(mat1->id>=0 && mat1->id==mat2->id) {
cpmPhys->E=mat1->young;
cpmPhys->G=mat1->young*mat1->G_over_E;
cpmPhys->tanFrictionAngle=tan(mat1->frictionAngle);
cpmPhys->undamagedCohesion=mat1->sigmaT;
- cpmPhys->epsFracture=mat1->relDuctility*mat1->epsCrackOnset;
cpmPhys->isCohesive=(cohesiveThresholdIter<0 || scene->iter<cohesiveThresholdIter);
#define _CPATTR(a) cpmPhys->a=mat1->a
_CPATTR(epsCrackOnset);
+ _CPATTR(crackOpening);
_CPATTR(neverDamage);
_CPATTR(dmgTau);
_CPATTR(dmgRateExp);
@@ -56,7 +64,7 @@
cpmPhys->G=.5*(mat1->G_over_E+mat2->G_over_E)*.5*(mat1->young+mat2->young);
cpmPhys->tanFrictionAngle=tan(.5*(mat1->frictionAngle+mat2->frictionAngle));
cpmPhys->undamagedCohesion=.5*(mat1->sigmaT+mat2->sigmaT);
- cpmPhys->epsFracture=.5*(mat1->relDuctility+mat2->relDuctility)*.5*(mat1->epsCrackOnset+mat2->epsCrackOnset);
+ _AVGATTR(crackOpening);
cpmPhys->isCohesive=(cohesiveThresholdIter<0 || scene->iter<cohesiveThresholdIter);
_AVGATTR(epsCrackOnset);
cpmPhys->neverDamage=(mat1->neverDamage || mat2->neverDamage);
@@ -67,6 +75,7 @@
_AVGATTR(isoPrestress);
#undef _AVGATTR
}
+
// NOTE: some params are not assigned until in Law2_Dem3DofGeom_CpmPhys_Cpm, since they need geometry as well; those are:
// crossSection, kn, ks
}
@@ -152,11 +161,37 @@
// shorthands
Real& epsN(BC->epsN);
- Vector3r& epsT(BC->epsT); Real& kappaD(BC->kappaD); Real& epsPlSum(BC->epsPlSum); const Real& E(BC->E); const Real& undamagedCohesion(BC->undamagedCohesion); const Real& tanFrictionAngle(BC->tanFrictionAngle); const Real& G(BC->G); const Real& crossSection(BC->crossSection); const Real& omegaThreshold(Law2_Dem3DofGeom_CpmPhys_Cpm::omegaThreshold); const Real& epsCrackOnset(BC->epsCrackOnset); Real& relResidualStrength(BC->relResidualStrength); const Real& epsFracture(BC->epsFracture); const bool& neverDamage(BC->neverDamage); Real& omega(BC->omega); Real& sigmaN(BC->sigmaN); Vector3r& sigmaT(BC->sigmaT); Real& Fn(BC->Fn); Vector3r& Fs(BC->Fs); // for python access
- const bool& isCohesive(BC->isCohesive);
+ Vector3r& epsT(BC->epsT);
+ Real& kappaD(BC->kappaD);
+ Real& epsPlSum(BC->epsPlSum);
+ const Real& E(BC->E);
+ const Real& undamagedCohesion(BC->undamagedCohesion);
+ const Real& tanFrictionAngle(BC->tanFrictionAngle);
+ const Real& G(BC->G);
+ const Real& crossSection(BC->crossSection);
+ const Real& omegaThreshold(Law2_Dem3DofGeom_CpmPhys_Cpm::omegaThreshold);
+ const Real& epsCrackOnset(BC->epsCrackOnset);
+ Real& relResidualStrength(BC->relResidualStrength);
+ const Real& crackOpening(BC->crackOpening);
+ const int& damLaw(BC->damLaw);
+ const bool& neverDamage(BC->neverDamage);
+ Real& omega(BC->omega);
+ Real& sigmaN(BC->sigmaN);
+ Vector3r& sigmaT(BC->sigmaT);
+ Real& Fn(BC->Fn);
+ Vector3r& Fs(BC->Fs); // for python access
+ const bool& isCohesive(BC->isCohesive);
#ifdef CPM_MATERIAL_MODEL
- Real& epsNPl(BC->epsNPl); const Real& dt=scene->dt; const Real& dmgTau(BC->dmgTau); const Real& plTau(BC->plTau);const Real& yieldLogSpeed(Law2_Dem3DofGeom_CpmPhys_Cpm::yieldLogSpeed); const int& yieldSurfType(Law2_Dem3DofGeom_CpmPhys_Cpm::yieldSurfType); const Real& yieldEllipseShift(Law2_Dem3DofGeom_CpmPhys_Cpm::yieldEllipseShift); const Real& epsSoft(Law2_Dem3DofGeom_CpmPhys_Cpm::epsSoft); const Real& relKnSoft(Law2_Dem3DofGeom_CpmPhys_Cpm::relKnSoft);
+ Real& epsNPl(BC->epsNPl);
+ const Real& dt=scene->dt;
+ const Real& dmgTau(BC->dmgTau);
+ const Real& plTau(BC->plTau);
+ const Real& yieldLogSpeed(Law2_Dem3DofGeom_CpmPhys_Cpm::yieldLogSpeed);
+ const int& yieldSurfType(Law2_Dem3DofGeom_CpmPhys_Cpm::yieldSurfType);
+ const Real& yieldEllipseShift(Law2_Dem3DofGeom_CpmPhys_Cpm::yieldEllipseShift);
+ const Real& epsSoft(Law2_Dem3DofGeom_CpmPhys_Cpm::epsSoft);
+ const Real& relKnSoft(Law2_Dem3DofGeom_CpmPhys_Cpm::relKnSoft);
#endif
epsN=contGeom->strainN(); epsT=contGeom->strainT();
@@ -189,7 +224,8 @@
epsN+=BC->isoPrestress/E;
// very simplified version of the constitutive law
kappaD=max(max(0.,epsN),kappaD); // internal variable, max positive strain (non-decreasing)
- omega=isCohesive?funcG(kappaD,epsCrackOnset,epsFracture,neverDamage):1.; // damage variable (non-decreasing, as funcG is also non-decreasing)
+ Real epsFracture = crackOpening/contGeom->refLength;
+ omega=isCohesive?funcG(kappaD,epsCrackOnset,epsFracture,neverDamage,damLaw):1.; // damage variable (non-decreasing, as funcG is also non-decreasing)
sigmaN=(1-(epsN>0?omega:0))*E*epsN; // damage taken in account in tension only
sigmaT=G*epsT; // trial stress
Real yieldSigmaT=max((Real)0.,undamagedCohesion*(1-omega)-sigmaN*tanFrictionAngle); // Mohr-Coulomb law with damage
@@ -203,7 +239,7 @@
#endif
sigmaN-=BC->isoPrestress;
- NNAN(kappaD); NNAN(epsCrackOnset); NNAN(epsFracture); NNAN(omega);
+ NNAN(kappaD); NNAN(epsCrackOnset); NNAN(crackOpening); NNAN(omega);
NNAN(sigmaN); NNANV(sigmaT); NNAN(crossSection);
// handle broken contacts
@@ -399,18 +435,18 @@
shared_ptr<CpmPhys> phys=dynamic_pointer_cast<CpmPhys>(I->phys);
if(!phys) continue;
const Body::id_t id1=I->getId1(), id2=I->getId2();
- GenericSpheresContact* geom=YADE_CAST<GenericSpheresContact*>(I->geom.get());
+ Dem3DofGeom* geom=YADE_CAST<Dem3DofGeom*>(I->geom.get());
+
+ Matrix3r stress = Matrix3r::Zero();
+ const Vector3r& n= geom->normal;
+ const Real& Fn = phys->Fn;
+ const Vector3r& Fs = phys->Fs;
+ //stress[i,j] += geom->refLength*(Fn*n[i]*n[j]+0.5*(Fs[i]*n[j]+Fs[j]*n[i]));
+ //stress += geom->refLength*(Fn*outer(n,n)+.5*(outer(Fs,n)+outer(n,Fs)));
+ stress += geom->refLength*(Fn*n*n.transpose()+.5*(Fs*n.transpose()+n*Fs.transpose()));
- Vector3r normalStress=((1./phys->crossSection)*geom->normal.dot(phys->normalForce))*geom->normal;
- bodyStats[id1].sigma+=normalStress; bodyStats[id2].sigma+=normalStress;
- Vector3r shearStress;
- for(int i=0; i<3; i++){
- int ix1=(i+1)%3,ix2=(i+2)%3;
- shearStress[i]=geom->normal[ix1]*phys->shearForce[ix1]+geom->normal[ix2]*phys->shearForce[ix2];
- shearStress[i]/=phys->crossSection;
- }
- bodyStats[id1].tau+=shearStress;
- bodyStats[id2].tau+=shearStress;
+ bodyStats[id1].stress += stress;
+ bodyStats[id2].stress += stress;
bodyStats[id1].nLinks++; bodyStats[id2].nLinks++;
if(!phys->isCohesive) continue;
@@ -422,13 +458,13 @@
nAvgRelResidual+=1;
}
FOREACH(shared_ptr<Body> B, *scene->bodies){
- if (!B) continue;
+ //if (!B) continue;
const Body::id_t& id=B->getId();
// add damaged contacts that have already been deleted
CpmState* state=dynamic_cast<CpmState*>(B->state.get());
if(!state) continue;
- state->sigma=bodyStats[id].sigma;
- state->tau=bodyStats[id].tau;
+ //state->sigma=bodyStats[id].sigma;
+ //state->tau=bodyStats[id].tau;
int cohLinksWhenever=bodyStats[id].nCohLinks+state->numBrokenCohesive;
if(cohLinksWhenever>0){
state->normDmg=(bodyStats[id].dmgSum+state->numBrokenCohesive)/cohLinksWhenever;
=== modified file 'pkg/dem/ConcretePM.hpp'
--- pkg/dem/ConcretePM.hpp 2011-01-09 16:34:50 +0000
+++ pkg/dem/ConcretePM.hpp 2011-08-27 17:40:30 +0000
@@ -51,6 +51,7 @@
#include<yade/pkg/common/PeriodicEngines.hpp>
#include<yade/pkg/common/NormShearPhys.hpp>
#include<yade/pkg/dem/DemXDofGeom.hpp>
+#include<yade/pkg/dem/DemXDofGeom.hpp>
namespace py=boost::python;
@@ -65,8 +66,9 @@
((Real,normDmg,0,,"Average damage including already deleted contacts (it is really not damage, but 1-relResidualStrength now)"))
((Real,epsPlBroken,0,,"Plastic strain on contacts already deleted (bogus values)"))
((Real,normEpsPl,0,,"Sum of plastic strains normalized by number of contacts (bogus values)"))
- ((Vector3r,sigma,Vector3r::Zero(),,"Normal stresses on the particle"))
- ((Vector3r,tau,Vector3r::Zero(),,"Shear stresses on the particle.")),
+ //((Vector3r,sigma,Vector3r::Zero(),,"Normal stresses on the particle"))
+ //((Vector3r,tau,Vector3r::Zero(),,"Shear stresses on the particle."))
+ ((Matrix3r,stress,Matrix3r::Zero(),,"Stress tensor on the particle.")),
/*ctor*/ createIndex();
);
REGISTER_CLASS_INDEX(CpmState,State);
@@ -85,6 +87,8 @@
((bool,neverDamage,false,,"If true, no damage will occur (for testing only)."))
((Real,epsCrackOnset,NaN,,"Limit elastic strain [-]"))
((Real,relDuctility,NaN,,"Relative ductility, for damage evolution law peak right-tangent. [-]"))
+ ((Real,crackOpening,NaN,,"Crack opening when the crack is fully broken in tension. [m]"))
+ ((int,damLaw,1,,"Law for gamage evolution in uniaxial tension. 0 for linear stress-strain softening branch, 1 for exponential damage evolution law"))
((Real,dmgTau,((void)"deactivated if negative",-1),,"Characteristic time for normal viscosity. [s]"))
((Real,dmgRateExp,0,,"Exponent for normal viscosity function. [-]"))
((Real,plTau,((void)"deactivated if negative",-1),,"Characteristic time for visco-plasticity. [s]"))
@@ -122,7 +126,7 @@
((Real,undamagedCohesion,NaN,,"virgin material cohesion [Pa]"))
((Real,crossSection,NaN,,"equivalent cross-section associated with this contact [m²]"))
((Real,epsCrackOnset,NaN,,"strain at which the material starts to behave non-linearly"))
- ((Real,epsFracture,NaN,,"strain where the damage-evolution law tangent from the top (epsCrackOnset) touches the axis; since the softening law is exponential, this doesn't mean that the contact is fully damaged at this point, that happens only asymptotically"))
+ ((Real,crackOpening,NaN,,"Crack opening when the crack is fully broken in tension. [m]"))
((Real,dmgTau,-1,,"characteristic time for damage (if non-positive, the law without rate-dependence is used)"))
((Real,dmgRateExp,0,,"exponent in the rate-dependent damage evolution"))
((Real,dmgStrain,0,,"damage strain (at previous or current step)"))
@@ -133,6 +137,7 @@
((Real,kappaD,0,,"Up to now maximum normal strain (semi-norm), non-decreasing in time."))
((Real,epsNPl,0,,"normal plastic strain (initially zero)"))
((bool,neverDamage,false,,"the damage evolution function will always return virgin state"))
+ ((int,damLaw,-1,,"Law for softening part of uniaxial tension. 0 for linear, 1 for exponential"))
((Real,epsTrans,0,,"Transversal strain (perpendicular to the contact axis)"))
((Real,epsPlSum,0,,"cummulative shear plastic strain measure (scalar) on this contact"))
((bool,isCohesive,false,,"if not cohesive, interaction is deleted when distance is greater than zero."))
@@ -175,9 +180,14 @@
class Law2_Dem3DofGeom_CpmPhys_Cpm: public LawFunctor{
public:
/*! Damage evolution law */
- static Real funcG(const Real& kappaD, const Real& epsCrackOnset, const Real& epsFracture, const bool& neverDamage) {
+ static Real funcG(const Real& kappaD, const Real& epsCrackOnset, const Real& epsFracture, const bool& neverDamage, const int& damLaw) {
if(kappaD<epsCrackOnset || neverDamage) return 0;
- return 1.-(epsCrackOnset/kappaD)*exp(-(kappaD-epsCrackOnset)/epsFracture);
+ switch (damLaw) {
+ case 0: // linear
+ return (1.-epsCrackOnset/kappaD)/(1.-epsCrackOnset/epsFracture);
+ case 1: // exponential
+ return 1.-(epsCrackOnset/kappaD)*exp(-(kappaD-epsCrackOnset)/epsFracture);
+ }
}
//! return |sigmaT| at plastic surface for given sigmaN etc; not used by the law itself
Real yieldSigmaTMagnitude(Real sigmaN, Real omega, Real undamagedCohesion, Real tanFrictionAngle);
@@ -193,7 +203,7 @@
((Real,epsSoft,((void)"approximates confinement -20MPa precisely, -100MPa a little over, -200 and -400 are OK (secant)",-3e-3),,"Strain at which softening in compression starts (non-negative to deactivate)"))
((Real,relKnSoft,.3,,"Relative rigidity of the softening branch in compression (0=perfect elastic-plastic, <0 softening, >0 hardening)")),
/*ctor*/,
- .def("funcG",&Law2_Dem3DofGeom_CpmPhys_Cpm::funcG,(py::arg("kappaD"),py::arg("epsCrackOnset"),py::arg("epsFracture"),py::arg("neverDamage")=false),"Damage evolution law, evaluating the $\\omega$ parameter. $\\kappa_D$ is historically maximum strain, *epsCrackOnset* ($\\varepsilon_0$) = :yref:`CpmPhys.epsCrackOnset`, *epsFracture* = :yref:`CpmPhys.epsFracture`; if *neverDamage* is ``True``, the value returned will always be 0 (no damage).")
+ .def("funcG",&Law2_Dem3DofGeom_CpmPhys_Cpm::funcG,(py::arg("kappaD"),py::arg("epsCrackOnset"),py::arg("epsFracture"),py::arg("neverDamage")=false,py::arg("damLaw")=1),"Damage evolution law, evaluating the $\\omega$ parameter. $\\kappa_D$ is historically maximum strain, *epsCrackOnset* ($\\varepsilon_0$) = :yref:`CpmPhys.epsCrackOnset`, *epsFracture* = :yref:`CpmPhys.epsFracture`; if *neverDamage* is ``True``, the value returned will always be 0 (no damage). TODO")
.def("yieldSigmaTMagnitude",&Law2_Dem3DofGeom_CpmPhys_Cpm::yieldSigmaTMagnitude,(py::arg("sigmaN"),py::arg("omega"),py::arg("undamagedCohesion"),py::arg("tanFrictionAngle")),"Return radius of yield surface for given material and state parameters; uses attributes of the current instance (*yieldSurfType* etc), change them before calling if you need that.")
);
DECLARE_LOGGER;
@@ -237,7 +247,7 @@
#endif
class CpmStateUpdater: public PeriodicEngine {
- struct BodyStats{ int nCohLinks; int nLinks; Real dmgSum, epsPlSum; Vector3r sigma, tau; BodyStats(): nCohLinks(0), nLinks(0), dmgSum(0.), epsPlSum(0.), sigma(Vector3r::Zero()), tau(Vector3r::Zero()) {} };
+ struct BodyStats{ int nCohLinks; int nLinks; Real dmgSum, epsPlSum; Matrix3r stress; BodyStats(): nCohLinks(0), nLinks(0), dmgSum(0.), epsPlSum(0.), stress(Matrix3r::Zero()) {} };
public:
virtual void action(){ update(scene); }
void update(Scene* rb=NULL);
=== modified file 'pkg/dem/VTKRecorder.cpp'
--- pkg/dem/VTKRecorder.cpp 2011-04-18 08:00:23 +0000
+++ pkg/dem/VTKRecorder.cpp 2011-08-27 17:40:30 +0000
@@ -166,15 +166,9 @@
vtkSmartPointer<vtkFloatArray> cpmDamage = vtkSmartPointer<vtkFloatArray>::New();
cpmDamage->SetNumberOfComponents(1);
cpmDamage->SetName("cpmDamage");
- vtkSmartPointer<vtkFloatArray> cpmSigma = vtkSmartPointer<vtkFloatArray>::New();
- cpmSigma->SetNumberOfComponents(3);
- cpmSigma->SetName("cpmSigma");
- vtkSmartPointer<vtkFloatArray> cpmSigmaM = vtkSmartPointer<vtkFloatArray>::New();
- cpmSigmaM->SetNumberOfComponents(1);
- cpmSigmaM->SetName("cpmSigmaM");
- vtkSmartPointer<vtkFloatArray> cpmTau = vtkSmartPointer<vtkFloatArray>::New();
- cpmTau->SetNumberOfComponents(3);
- cpmTau->SetName("cpmTau");
+ vtkSmartPointer<vtkFloatArray> cpmStress = vtkSmartPointer<vtkFloatArray>::New();
+ cpmStress->SetNumberOfComponents(9);
+ cpmStress->SetName("cpmStress");
// extras for RPM
vtkSmartPointer<vtkFloatArray> rpmSpecNum = vtkSmartPointer<vtkFloatArray>::New();
@@ -312,13 +306,10 @@
if (recActive[REC_CPM]){
cpmDamage->InsertNextValue(YADE_PTR_CAST<CpmState>(b->state)->normDmg);
- const Vector3r& ss=YADE_PTR_CAST<CpmState>(b->state)->sigma;
- const Vector3r& tt=YADE_PTR_CAST<CpmState>(b->state)->tau;
- float s[3]={ss[0],ss[1],ss[2]};
- float t[3]={tt[0],tt[1],tt[2]};
- cpmSigma->InsertNextTupleValue(s);
- cpmSigmaM->InsertNextValue((ss[0]+ss[1]+ss[2])/3.);
- cpmTau->InsertNextTupleValue(t);
+ const Matrix3r& ss=YADE_PTR_CAST<CpmState>(b->state)->stress;
+ //float s[3]={ss[0],ss[1],ss[2]};
+ float s[9]={ss[0,0],ss[0,1],ss[0,2],ss[1,0],ss[1,1],ss[1,2],ss[2,0],ss[2,1],ss[2,2]};
+ cpmStress->InsertNextTupleValue(s);
}
if (recActive[REC_RPM]){
rpmSpecNum->InsertNextValue(YADE_PTR_CAST<RpmState>(b->state)->specimenNumber);
@@ -388,9 +379,7 @@
}
if (recActive[REC_CPM]){
spheresUg->GetPointData()->AddArray(cpmDamage);
- spheresUg->GetPointData()->AddArray(cpmSigma);
- spheresUg->GetPointData()->AddArray(cpmSigmaM);
- spheresUg->GetPointData()->AddArray(cpmTau);
+ spheresUg->GetPointData()->AddArray(cpmStress);
}
if (recActive[REC_RPM]){
spheresUg->GetPointData()->AddArray(rpmSpecNum);
=== modified file 'pkg/dem/VTKRecorder.hpp'
--- pkg/dem/VTKRecorder.hpp 2011-04-18 03:39:47 +0000
+++ pkg/dem/VTKRecorder.hpp 2011-08-27 17:40:30 +0000
@@ -20,7 +20,7 @@
((bool,multiblock,false,,"Use multi-block (``.vtm``) files to store data, rather than separate ``.vtu`` files."))
#endif
((string,fileName,"",,"Base file name; it will be appended with {spheres,intrs,facets}-243100.vtu (unless *multiblock* is ``True``) depending on active recorders and step number (243100 in this case). It can contain slashes, but the directory must exist already."))
- ((vector<string>,recorders,vector<string>(1,string("all")),,"List of active recorders (as strings). ``all`` (the default value) enables all base and generic recorders.\n\n.. admonition:: Base recorders\n\n\tBase recorders save the geometry (unstructured grids) on which other data is defined. They are implicitly activated by many of the other recorders. Each of them creates a new file (or a block, if :yref:`multiblock <VTKRecorder.multiblock>` is set).\n\n\t``spheres``\n\t\tSaves positions and radii (``radii``) of :yref:`spherical<Sphere>` particles.\n\t``facets``\n\t\tSave :yref:`facets<Facet>` positions (vertices).\n\t``intr``\n\t\tStore interactions as lines between nodes at respective particles positions. Additionally stores magnitude of normal (``forceN``) and shear (``absForceT``) forces on interactions (the :yref:`geom<Interaction.geom> must be of type :yref:`NormShearPhys`). \n\n.. admonition:: Generic recorders\n\n\tGeneric recorders do not depend on specific model being used and save commonly useful data.\n\n\t``id``\n\t\tSaves id's (field ``id``) of spheres; active only if ``spheres`` is active.\n\t``clumpId``\n\t\tSaves id's of clumps to which each sphere belongs (field ``clumpId``); active only if ``spheres`` is active.\n\t``colors``\n\t\tSaves colors of :yref:`spheres<Sphere>` and of :yref:`facets<Facet>` (field ``color``); only active if ``spheres`` or ``facets`` are activated.\n\t``mask``\n\t\tSaves groupMasks of :yref:`spheres<Sphere>` and of :yref:`facets<Facet>` (field ``mask``); only active if ``spheres`` or ``facets`` are activated.\n\t``materialId``\n\t\tSaves materialID of :yref:`spheres<Sphere>` and of :yref:`facets<Facet>`; only active if ``spheres`` or ``facets`` are activated.\n\t``velocity``\n\t\tSaves linear and angular velocities of spherical particles as Vector3 and length(fields ``linVelVec``, ``linVelLen`` and ``angVelVec``, ``angVelLen`` respectively``); only effective with ``spheres``.\n\t``stress``\n\t\tSaves stresses of :yref:`spheres<Sphere>` and of :yref:`facets<Facet>` as Vector3 and length; only active if ``spheres`` or ``facets`` are activated.\n\n.. admonition:: Specific recorders\n\n\tThe following should only be activated in appropriate cases, otherwise crashes can occur due to violation of type presuppositions.\n\n\t``cpm``\n\t\tSaves data pertaining to the :yref:`concrete model<Law2_Dem3DofGeom_CpmPhys_Cpm>`: ``cpmDamage`` (normalized residual strength averaged on particle), ``cpmSigma`` (stress on particle, normal components), ``cpmTau`` (shear components of stress on particle), ``cpmSigmaM`` (mean stress around particle); ``intr`` is activated automatically by ``cpm``\n\t``rpm``\n\t\tSaves data pertaining to the :yref:`rock particle model<Law2_Dem3DofGeom_RockPMPhys_Rpm>`: ``rpmSpecNum`` shows different pieces of separated stones, only ids. ``rpmSpecMass`` shows masses of separated stones.\n\t``wpm``\n\t\tSaves data pertaining to the :yref:`wire particle model<Law2_ScGeom_WirePhys_WirePM>`: ``wpmForceNFactor`` shows the loading factor for the wire, e.g. normal force devided by threshold normal force.\n\n"))
+ ((vector<string>,recorders,vector<string>(1,string("all")),,"List of active recorders (as strings). ``all`` (the default value) enables all base and generic recorders.\n\n.. admonition:: Base recorders\n\n\tBase recorders save the geometry (unstructured grids) on which other data is defined. They are implicitly activated by many of the other recorders. Each of them creates a new file (or a block, if :yref:`multiblock <VTKRecorder.multiblock>` is set).\n\n\t``spheres``\n\t\tSaves positions and radii (``radii``) of :yref:`spherical<Sphere>` particles.\n\t``facets``\n\t\tSave :yref:`facets<Facet>` positions (vertices).\n\t``intr``\n\t\tStore interactions as lines between nodes at respective particles positions. Additionally stores magnitude of normal (``forceN``) and shear (``absForceT``) forces on interactions (the :yref:`geom<Interaction.geom> must be of type :yref:`NormShearPhys`). \n\n.. admonition:: Generic recorders\n\n\tGeneric recorders do not depend on specific model being used and save commonly useful data.\n\n\t``id``\n\t\tSaves id's (field ``id``) of spheres; active only if ``spheres`` is active.\n\t``clumpId``\n\t\tSaves id's of clumps to which each sphere belongs (field ``clumpId``); active only if ``spheres`` is active.\n\t``colors``\n\t\tSaves colors of :yref:`spheres<Sphere>` and of :yref:`facets<Facet>` (field ``color``); only active if ``spheres`` or ``facets`` are activated.\n\t``mask``\n\t\tSaves groupMasks of :yref:`spheres<Sphere>` and of :yref:`facets<Facet>` (field ``mask``); only active if ``spheres`` or ``facets`` are activated.\n\t``materialId``\n\t\tSaves materialID of :yref:`spheres<Sphere>` and of :yref:`facets<Facet>`; only active if ``spheres`` or ``facets`` are activated.\n\t``velocity``\n\t\tSaves linear and angular velocities of spherical particles as Vector3 and length(fields ``linVelVec``, ``linVelLen`` and ``angVelVec``, ``angVelLen`` respectively``); only effective with ``spheres``.\n\t``stress``\n\t\tSaves stresses of :yref:`spheres<Sphere>` and of :yref:`facets<Facet>` as Vector3 and length; only active if ``spheres`` or ``facets`` are activated.\n\n.. admonition:: Specific recorders\n\n\tThe following should only be activated in appropriate cases, otherwise crashes can occur due to violation of type presuppositions.\n\n\t``cpm``\n\t\tSaves data pertaining to the :yref:`concrete model<Law2_Dem3DofGeom_CpmPhys_Cpm>`: ``cpmDamage`` (normalized residual strength averaged on particle), ``cpmStress`` (stress on particle); ``intr`` is activated automatically by ``cpm``\n\t``rpm``\n\t\tSaves data pertaining to the :yref:`rock particle model<Law2_Dem3DofGeom_RockPMPhys_Rpm>`: ``rpmSpecNum`` shows different pieces of separated stones, only ids. ``rpmSpecMass`` shows masses of separated stones.\n\t``wpm``\n\t\tSaves data pertaining to the :yref:`wire particle model<Law2_ScGeom_WirePhys_WirePM>`: ``wpmForceNFactor`` shows the loading factor for the wire, e.g. normal force devided by threshold normal force.\n\n"))
((int,mask,0,,"If mask defined, only bodies with corresponding groupMask will be exported. If 0, all bodies will be exported.")),
/*ctor*/
initRun=true;
=== modified file 'py/export.py'
--- py/export.py 2010-09-02 09:11:04 +0000
+++ py/export.py 2011-08-27 17:40:30 +0000
@@ -1,8 +1,9 @@
"""
-Export geometry to various formats.
+Export (not only) geometry to various formats.
"""
# encoding: utf-8
from yade.wrapper import *
+from yade import *
class VTKWriter:
"""
@@ -179,3 +180,129 @@
"""
return (textExt(filename=filename, format='x_y_z_r',consider=consider))
+
+
+
+
+class VTKExporter:
+ """Class for exporting data to VTK Simple Legacy File (for example if, for some reasin, you are not able to use VTKRecorder). Export of spheres and facets is supported.
+
+ USAGE:
+ create object vtkExporter = VTKExporter('baseFileName'),
+ add to engines PyRunner with command='vtkWriter.exportSomething(params)'
+ """
+ def __init__(self,baseName,startSnap=0):
+ self.spheresSnapCount = startSnap
+ self.facetsSnapCount = startSnap
+ self.baseName = baseName
+
+ def exportSpheres(self,ids='all',what=[],comment="comment",numLabel=None):
+ """
+ exports spheres (positions and radius) and defined properties.
+:parameters:
+ `ids`: list | "all"
+ if "all", then export all spheres, otherwise only spheres from integer list
+ `what`: [tuple(2)]
+ what other than then position and radius export. parameter is list of couple (name,command). Name is string under which it is save to vtk, command is string to evaluate. Node that the bodies are labeled as b in this function. Scalar, vector and tensor variables are supported. For example, to export velocity (with name particleVelocity) and the distance form point (0,0,0) (named as dist) you should write: ... what=[('particleVelocity','b.state.vel'),('dist','b.state.pos.norm()', ...
+ `comment`: string
+ comment to add to vtk file
+ `numLabel`: int
+ number of file (e.g. time step), if unspecified, the last used value + 1 will be used
+ """
+ allIds = False
+ if ids=='all':
+ ids=xrange(len(O.bodies))
+ allIds = True
+ bodies = []
+ for i in ids:
+ b = O.bodies[i]
+ if b.shape.__class__.__name__!="Sphere":
+ if not allIds: print "Warning: body %d is not Sphere"%(i)
+ continue
+ bodies.append(b)
+ n = len(bodies)
+ outFile = open(self.baseName+'-spheres-%04d'%self.spheresSnapCount+'.vtk', 'w')
+ outFile.write("# vtk DataFile Version 3.0.\n%s\nASCII\n\nDATASET POLYDATA\nPOINTS %d double\n"%(comment,n))
+ for b in bodies:
+ pos = b.state.pos
+ outFile.write("%g %g %g\n"%(pos[0],pos[1],pos[2]))
+ outFile.write("\nPOINT_DATA %d\nSCALARS radius double 1\nLOOKUP_TABLE default\n"%(n))
+ for b in bodies:
+ outFile.write("%g\n"%(b.shape.radius))
+ for name,command in what:
+ test = eval(command)
+ if isinstance(test,Matrix3):
+ outFile.write("\nTENSORS %s double\n"%(name))
+ for b in bodies:
+ t = eval(command)
+ outFile.write("%g %g %g\n%g %g %g\n%g %g %g\n\n"%(t[0,0],t[0,1],t[0,2],t[1,0],t[1,1],t[1,2],t[2,0],t[2,1],t[2,2]))
+ if isinstance(test,Vector3):
+ outFile.write("\nVECTORS %s double\n"%(name))
+ for b in bodies:
+ v = eval(command)
+ outFile.write("%g %g %g\n"%(v[0],v[1],v[2]))
+ else:
+ outFile.write("\nSCALARS %s double 1\nLOOKUP_TABLE default\n"%(name))
+ for b in bodies:
+ outFile.write("%g\n"%(eval(command)))
+ outFile.close()
+ self.spheresSnapCount += 1
+
+ def exportFacets(self,ids='all',what=[],comment="comment",numLabel=None):
+ """
+ exports facets (positions) and defined properties.
+:parameters:
+ `ids`: list | "all"
+ if "all", then export all spheres, otherwise only spheres from integer list
+ `what`: [tuple(2)]
+ see exportSpheres
+ `comment`: string
+ comment to add to vtk file
+ `numLabel`: int
+ number of file (e.g. time step), if unspecified, the last used value + 1 will be used
+ """
+ allIds = False
+ if ids=='all':
+ ids=xrange(len(O.bodies))
+ allIds = True
+ bodies = []
+ for i in ids:
+ b = O.bodies[i]
+ if b.shape.__class__.__name__!="Facet":
+ if not allIds: print "Warning: body %d is not Facet"%(i)
+ continue
+ bodies.append(b)
+ n = len(bodies)
+ outFile = open(self.baseName+'-facets-%04d'%self.facetsSnapCount+'.vtk', 'w')
+ outFile.write("# vtk DataFile Version 3.0.\n%s\nASCII\n\nDATASET POLYDATA\nPOINTS %d double\n"%(comment,3*n))
+ for b in bodies:
+ pt1 = b.state.pos + b.shape.vertices[0]
+ pt2 = b.state.pos + b.shape.vertices[1]
+ pt3 = b.state.pos + b.shape.vertices[2]
+ outFile.write("%g %g %g\n"%(pt1[0],pt1[1],pt1[2]))
+ outFile.write("%g %g %g\n"%(pt2[0],pt2[1],pt2[2]))
+ outFile.write("%g %g %g\n"%(pt3[0],pt3[1],pt3[2]))
+ outFile.write("\nPOLYGONS %d %d\n"%(n,4*n))
+ i = 0
+ for b in bodies:
+ outFile.write("3 %d %d %d\n"%(i,i+1,i+2))
+ i += 3
+ if what:
+ outFile.write("\nCELL_DATA %d"%(n))
+ for name,command in what:
+ test = eval(command)
+ if isinstance(test,Matrix3):
+ outFile.write("\nTENSORS %s double\n"%(name))
+ for b in bodies:
+ t = eval(command)
+ outFile.write("%g %g %g\n%g %g %g\n%g %g %g\n\n"%(t[0,0],t[0,1],t[0,2],t[1,0],t[1,1],t[1,2],t[2,0],t[2,1],t[2,2]))
+ if isinstance(test,Vector3):
+ outFile.write("\nVECTORS %s double\n"%(name))
+ for b in bodies:
+ v = eval(command)
+ outFile.write("%g %g %g\n"%(v[0],v[1],v[2]))
+ else:
+ outFile.write("\nSCALARS %s double 1\nLOOKUP_TABLE default\n"(name))
+ for b in bodies:
+ outFile.write("%g\n"%(eval(command)))
+ outFile.close()
=== modified file 'py/mathWrap/miniEigen.cpp'
--- py/mathWrap/miniEigen.cpp 2011-04-28 11:13:29 +0000
+++ py/mathWrap/miniEigen.cpp 2011-08-27 17:40:30 +0000
@@ -100,6 +100,7 @@
static bp::tuple Quaternionr_toAngleAxis(const Quaternionr& self){ AngleAxisr aa(self); return bp::make_tuple(aa.angle(),aa.axis());}
static Real Vector3r_dot(const Vector3r& self, const Vector3r& v){ return self.dot(v); }
+static Matrix3r Vector3r_outer(const Vector3r& self, const Vector3r& v){ return self*v.transpose(); }
static Real Vector3i_dot(const Vector3i& self, const Vector3i& v){ return self.dot(v); }
static Real Vector2r_dot(const Vector2r& self, const Vector2r& v){ return self.dot(v); }
static Real Vector2i_dot(const Vector2i& self, const Vector2i& v){ return self.dot(v); }
@@ -326,6 +327,7 @@
.add_static_property("UnitX",&Vector3r_UnitX).add_static_property("UnitY",&Vector3r_UnitY).add_static_property("UnitZ",&Vector3r_UnitZ)
// methods
.def("dot",&Vector3r_dot).def("cross",&Vector3r_cross)
+ .def("outer",&Vector3r_outer)
.def("norm",&Vector3r::norm).def("squaredNorm",&Vector3r::squaredNorm).def("normalize",&Vector3r::normalize).def("normalized",&Vector3r::normalized)
// operators
.def("__neg__",&Vector3r__neg__) // -v
=== modified file 'py/ymport.py'
--- py/ymport.py 2011-04-12 22:06:51 +0000
+++ py/ymport.py 2011-08-27 17:40:30 +0000
@@ -242,3 +242,70 @@
c=sphereList[i].Centre()
ret.append(utils.sphere([shift[0]+scale*float(c.X()),shift[1]+scale*float(c.Y()),shift[2]+scale*float(c.Z())],scale*float(r),**kw))
return ret
+
+
+class UNVReader:
+ # class used in ymport.unv function
+ # reads and evaluate given unv file and extracts all triangles
+ # can be extended to read tetrahedrons as well
+ def __init__(self,fileName,shift,scale):
+ self.unvFile = open(fileName,'r')
+ self.flag = 0
+ self.line = self.unvFile.readline()
+ self.lineSplit = self.line.split()
+ self.vertices = []
+ self.verticesTris = []
+ self.read()
+ def readLine(self):
+ self.line = self.unvFile.readline()
+ self.lineSplit = self.line.split()
+ def read(self):
+ while self.line:
+ self.evalLine()
+ self.line = self.unvFile.readline()
+ self.unvFile.close()
+ def evalLine(self):
+ self.lineSplit = self.line.split()
+ if len(self.lineSplit) <= 1: # eval special unv format
+ if self.lineSplit[0] == '-1': pass
+ elif self.lineSplit[0] == '2411': self.flag = 1; # nodes
+ elif self.lineSplit[0] == '2412': self.flag = 2; # edges (lines)
+ else: self.flag = 4; # volume elements or other, not interesting for us
+ elif self.flag == 1: self.evalVertices()
+ elif self.flag == 2: self.evalEdge()
+ elif self.flag == 3: self.evalFacet()
+ #elif self.flag == 4: self.evalGroup()
+ def evalVertices(self):
+ self.readLine()
+ self.vertices.append(Vector3(float(self.lineSplit[0]),float(self.lineSplit[1]),float(self.lineSplit[2])))
+ def evalEdge(self):
+ if self.lineSplit[1]=='41': self.flag = 3; self.evalFacet()
+ else: self.readLine(); self.readLine()
+ def evalFacet(self):
+ if self.lineSplit[1]=='41': # triangle
+ self.readLine()
+ self.verticesTris.append([self.vertices[int(self.lineSplit[0])-1],self.vertices[int(self.lineSplit[1])-1],self.vertices[int(self.lineSplit[2])-1]])
+ else: # is not triangle
+ self.readLine()
+ self.flag = 4
+ # can be added function to handle tetrahedrons
+
+def unv(fileName,shift=(0,0,0),scale=1.0,**kw):
+ """ Import geometry from unv file, return list of created facets.
+ :Parameters:
+ `fileName`: string
+ name of unv file
+ `shift`: [float,float,float] | Vector3
+ [X,Y,Z] parameter moves the specimen.
+ `scale`: float
+ factor scales the given data.
+ `**kw`: (unused keyword arguments)
+ is passed to :yref:`yade.utils.facet`
+
+ unv files are mainly used for FEM analyses (are used by OOFEM and Abakus TODO), but triangular elements can be imported as facets.
+ These files cen be created e.g. with open-source free software Salome (salome-platform.org TODO)
+
+ Example: :ysrc:`scripts/test/unvRead.py`."""
+ unvReader = UNVReader(fileName,shift,scale,**kw)
+ return [utils.facet(tri,**kw) for tri in unvReader.verticesTris]
+
=== added file 'scripts/test/shell.unv'
--- scripts/test/shell.unv 1970-01-01 00:00:00 +0000
+++ scripts/test/shell.unv 2011-08-27 17:40:30 +0000
@@ -0,0 +1,388 @@
+ -1
+ 2411
+ 1 0 0 11
+ 2.0000000000000000E+00 2.0000000000000000E+00 0.0000000000000000E+00
+ 2 0 0 11
+ -2.0000000000000000E+00 2.0000000000000000E+00 0.0000000000000000E+00
+ 3 0 0 11
+ -2.0000000000000000E+00 -2.0000000000000000E+00 0.0000000000000000E+00
+ 4 0 0 11
+ 2.0000000000000000E+00 -2.0000000000000000E+00 0.0000000000000000E+00
+ 5 0 0 11
+ -2.0000000000000000E+00 0.0000000000000000E+00 5.0000000000000000E+00
+ 6 0 0 11
+ 0.0000000000000000E+00 -2.0000000000000000E+00 5.0000000000000000E+00
+ 7 0 0 11
+ 2.0000000000000000E+00 0.0000000000000000E+00 5.0000000000000000E+00
+ 8 0 0 11
+ 0.0000000000000000E+00 2.0000000000000000E+00 4.0000000000000000E+00
+ 9 0 0 11
+ 0.0000000000000000E+00 2.0000000000000000E+00 0.0000000000000000E+00
+ 10 0 0 11
+ -2.0000000000000000E+00 0.0000000000000000E+00 0.0000000000000000E+00
+ 11 0 0 11
+ 0.0000000000000000E+00 -2.0000000000000000E+00 0.0000000000000000E+00
+ 12 0 0 11
+ 2.0000000000000000E+00 0.0000000000000000E+00 0.0000000000000000E+00
+ 13 0 0 11
+ -2.0000000000000000E+00 -1.3333333333333333E+00 1.6666666666666665E+00
+ 14 0 0 11
+ -2.0000000000000000E+00 -6.6666666666666674E-01 3.3333333333333330E+00
+ 15 0 0 11
+ -9.9999999999999989E-01 -1.0000000000000002E+00 5.0000000000000000E+00
+ 16 0 0 11
+ 6.6666666666666663E-01 -2.0000000000000000E+00 3.3333333333333335E+00
+ 17 0 0 11
+ 1.3333333333333333E+00 -2.0000000000000000E+00 1.6666666666666667E+00
+ 18 0 0 11
+ 1.0000000000000002E+00 -9.9999999999999989E-01 5.0000000000000000E+00
+ 19 0 0 11
+ 2.0000000000000000E+00 6.6666666666666663E-01 3.3333333333333335E+00
+ 20 0 0 11
+ 2.0000000000000000E+00 1.3333333333333333E+00 1.6666666666666667E+00
+ 21 0 0 11
+ -6.6666666666666663E-01 2.0000000000000000E+00 2.6666666666666665E+00
+ 22 0 0 11
+ -1.3333333333333333E+00 2.0000000000000000E+00 1.3333333333333335E+00
+ 23 0 0 11
+ 1.0000000000000000E+00 1.0000000000000000E+00 4.5000000000000000E+00
+ 24 0 0 11
+ -1.0000000000000000E+00 1.0000000000000000E+00 4.5000000000000000E+00
+ 25 0 0 11
+ 1.0204493184023047E+00 -7.7864862082500297E-01 0.0000000000000000E+00
+ 26 0 0 11
+ -1.0201187217404240E+00 -7.7812918103283835E-01 0.0000000000000000E+00
+ 27 0 0 11
+ 3.0453560848275885E-04 -4.8202294743143659E-01 0.0000000000000000E+00
+ 28 0 0 11
+ 7.3726756851835962E-01 5.7599033144961864E-01 0.0000000000000000E+00
+ 29 0 0 11
+ -7.4105532977775734E-01 5.7416265686003176E-01 0.0000000000000000E+00
+ 30 0 0 11
+ -4.5335881725944732E-01 -1.7106968539974061E+00 1.2783890157485245E+00
+ 31 0 0 11
+ -1.2134510285411346E+00 -9.2246231036720572E-01 3.9994005896714038E+00
+ 32 0 0 11
+ -6.8830573726249644E-01 -1.3864928868902922E+00 2.7549859841550672E+00
+ 33 0 0 11
+ 2.9354685262518415E-01 -1.7467305604786523E+00 2.1909364193444514E+00
+ 34 0 0 11
+ -3.9030144358386054E-01 -1.4566594712643419E+00 3.9933632500875595E+00
+ 35 0 0 11
+ 1.7106968523456150E+00 -4.5335882551308387E-01 1.2783890184401399E+00
+ 36 0 0 11
+ 9.2246208165153232E-01 -1.2134513531061526E+00 3.9994006917316649E+00
+ 37 0 0 11
+ 1.3864928554125591E+00 -6.8830581491128173E-01 2.7549860091085718E+00
+ 38 0 0 11
+ 1.7467305560005697E+00 2.9354683729707765E-01 2.1909364243104328E+00
+ 39 0 0 11
+ 1.4566592554516338E+00 -3.9030178030592577E-01 3.9933633564664008E+00
+ 40 0 0 11
+ -4.2824811618915049E-01 1.8389978287971704E+00 7.2413926942816254E-01
+ 41 0 0 11
+ 1.1661801462264525E+00 8.7072239860184431E-01 3.8450614971766779E+00
+ 42 0 0 11
+ 7.2901780443584108E-01 1.2662209225398113E+00 2.5182690674205359E+00
+ 43 0 0 11
+ -1.4336639743190543E-01 1.6442723473109442E+00 1.6085811466545532E+00
+ 44 0 0 11
+ 7.8730263538445089E-01 1.6659819261076851E+00 1.0233616742521123E+00
+ 45 0 0 11
+ -1.8760034743520608E+00 -3.4804595344083816E-01 8.4764344257365620E-01
+ 46 0 0 11
+ -1.2112412597444577E+00 8.1931532181887712E-01 3.5964020714959761E+00
+ 47 0 0 11
+ -1.5610645352608976E+00 5.0178591478731538E-01 2.0222341159297650E+00
+ 48 0 0 11
+ -1.7512252703454867E+00 8.8962091633202278E-01 8.1096264863768097E-01
+ 49 0 0 11
+ -6.1203710584758919E-01 6.2139416891930609E-03 4.7722723518183638E+00
+ 50 0 0 11
+ 3.5060158157754001E-01 -5.7446690939163791E-01 4.8816744009823605E+00
+ 51 0 0 11
+ 3.4509485944386387E-01 5.8653845581367858E-01 4.5890839104346739E+00
+ -1
+ -1
+ 2412
+ 1 11 2 1 7 2
+ 0 0 0
+ 1 9
+ 2 11 2 1 7 2
+ 0 0 0
+ 9 2
+ 3 11 2 1 7 2
+ 0 0 0
+ 2 10
+ 4 11 2 1 7 2
+ 0 0 0
+ 10 3
+ 5 11 2 1 7 2
+ 0 0 0
+ 3 11
+ 6 11 2 1 7 2
+ 0 0 0
+ 11 4
+ 7 11 2 1 7 2
+ 0 0 0
+ 4 12
+ 8 11 2 1 7 2
+ 0 0 0
+ 12 1
+ 9 11 2 1 7 2
+ 0 0 0
+ 3 13
+ 10 11 2 1 7 2
+ 0 0 0
+ 13 14
+ 11 11 2 1 7 2
+ 0 0 0
+ 14 5
+ 12 11 2 1 7 2
+ 0 0 0
+ 5 15
+ 13 11 2 1 7 2
+ 0 0 0
+ 15 6
+ 14 11 2 1 7 2
+ 0 0 0
+ 6 16
+ 15 11 2 1 7 2
+ 0 0 0
+ 16 17
+ 16 11 2 1 7 2
+ 0 0 0
+ 17 4
+ 17 11 2 1 7 2
+ 0 0 0
+ 6 18
+ 18 11 2 1 7 2
+ 0 0 0
+ 18 7
+ 19 11 2 1 7 2
+ 0 0 0
+ 7 19
+ 20 11 2 1 7 2
+ 0 0 0
+ 19 20
+ 21 11 2 1 7 2
+ 0 0 0
+ 20 1
+ 22 11 2 1 7 2
+ 0 0 0
+ 8 21
+ 23 11 2 1 7 2
+ 0 0 0
+ 21 22
+ 24 11 2 1 7 2
+ 0 0 0
+ 22 2
+ 25 11 2 1 7 2
+ 0 0 0
+ 7 23
+ 26 11 2 1 7 2
+ 0 0 0
+ 23 8
+ 27 11 2 1 7 2
+ 0 0 0
+ 8 24
+ 28 11 2 1 7 2
+ 0 0 0
+ 24 5
+ 29 41 2 1 7 3
+ 10 26 3
+ 30 41 2 1 7 3
+ 4 25 12
+ 31 41 2 1 7 3
+ 25 4 11
+ 32 41 2 1 7 3
+ 11 27 25
+ 33 41 2 1 7 3
+ 26 27 11
+ 34 41 2 1 7 3
+ 11 3 26
+ 35 41 2 1 7 3
+ 25 28 12
+ 36 41 2 1 7 3
+ 12 28 1
+ 37 41 2 1 7 3
+ 2 29 10
+ 38 41 2 1 7 3
+ 29 26 10
+ 39 41 2 1 7 3
+ 27 26 29
+ 40 41 2 1 7 3
+ 28 29 9
+ 41 41 2 1 7 3
+ 29 2 9
+ 42 41 2 1 7 3
+ 1 28 9
+ 43 41 2 1 7 3
+ 25 27 28
+ 44 41 2 1 7 3
+ 29 28 27
+ 45 41 2 1 7 3
+ 31 15 5
+ 46 41 2 1 7 3
+ 31 5 14
+ 47 41 2 1 7 3
+ 30 11 17
+ 48 41 2 1 7 3
+ 13 3 30
+ 49 41 2 1 7 3
+ 30 3 11
+ 50 41 2 1 7 3
+ 33 30 17
+ 51 41 2 1 7 3
+ 31 14 32
+ 52 41 2 1 7 3
+ 32 14 13
+ 53 41 2 1 7 3
+ 32 13 30
+ 54 41 2 1 7 3
+ 31 34 15
+ 55 41 2 1 7 3
+ 17 11 4
+ 56 41 2 1 7 3
+ 32 34 31
+ 57 41 2 1 7 3
+ 34 6 15
+ 58 41 2 1 7 3
+ 16 6 34
+ 59 41 2 1 7 3
+ 33 16 32
+ 60 41 2 1 7 3
+ 17 16 33
+ 61 41 2 1 7 3
+ 30 33 32
+ 62 41 2 1 7 3
+ 34 32 16
+ 63 41 2 1 7 3
+ 36 18 6
+ 64 41 2 1 7 3
+ 36 6 16
+ 65 41 2 1 7 3
+ 35 12 20
+ 66 41 2 1 7 3
+ 35 17 4
+ 67 41 2 1 7 3
+ 35 4 12
+ 68 41 2 1 7 3
+ 17 35 37
+ 69 41 2 1 7 3
+ 38 35 20
+ 70 41 2 1 7 3
+ 36 16 37
+ 71 41 2 1 7 3
+ 37 16 17
+ 72 41 2 1 7 3
+ 36 39 18
+ 73 41 2 1 7 3
+ 20 12 1
+ 74 41 2 1 7 3
+ 37 39 36
+ 75 41 2 1 7 3
+ 39 7 18
+ 76 41 2 1 7 3
+ 19 7 39
+ 77 41 2 1 7 3
+ 38 19 37
+ 78 41 2 1 7 3
+ 20 19 38
+ 79 41 2 1 7 3
+ 35 38 37
+ 80 41 2 1 7 3
+ 39 37 19
+ 81 41 2 1 7 3
+ 19 41 7
+ 82 41 2 1 7 3
+ 21 42 43
+ 83 41 2 1 7 3
+ 42 44 43
+ 84 41 2 1 7 3
+ 42 21 8
+ 85 41 2 1 7 3
+ 41 23 7
+ 86 41 2 1 7 3
+ 21 43 22
+ 87 41 2 1 7 3
+ 2 22 40
+ 88 41 2 1 7 3
+ 40 9 2
+ 89 41 2 1 7 3
+ 23 41 8
+ 90 41 2 1 7 3
+ 20 1 44
+ 91 41 2 1 7 3
+ 42 8 41
+ 92 41 2 1 7 3
+ 40 44 9
+ 93 41 2 1 7 3
+ 20 44 42
+ 94 41 2 1 7 3
+ 42 41 19
+ 95 41 2 1 7 3
+ 40 43 44
+ 96 41 2 1 7 3
+ 1 9 44
+ 97 41 2 1 7 3
+ 40 22 43
+ 98 41 2 1 7 3
+ 42 19 20
+ 99 41 2 1 7 3
+ 46 24 8
+ 100 41 2 1 7 3
+ 46 5 24
+ 101 41 2 1 7 3
+ 45 48 10
+ 102 41 2 1 7 3
+ 45 10 3
+ 103 41 2 1 7 3
+ 47 14 46
+ 104 41 2 1 7 3
+ 47 45 13
+ 105 41 2 1 7 3
+ 45 3 13
+ 106 41 2 1 7 3
+ 48 2 10
+ 107 41 2 1 7 3
+ 5 46 14
+ 108 41 2 1 7 3
+ 47 13 14
+ 109 41 2 1 7 3
+ 48 22 2
+ 110 41 2 1 7 3
+ 22 47 21
+ 111 41 2 1 7 3
+ 46 8 21
+ 112 41 2 1 7 3
+ 46 21 47
+ 113 41 2 1 7 3
+ 48 47 22
+ 114 41 2 1 7 3
+ 45 47 48
+ 115 41 2 1 7 3
+ 15 6 50
+ 116 41 2 1 7 3
+ 50 6 18
+ 117 41 2 1 7 3
+ 49 15 50
+ 118 41 2 1 7 3
+ 8 24 51
+ 119 41 2 1 7 3
+ 49 24 5
+ 120 41 2 1 7 3
+ 49 5 15
+ 121 41 2 1 7 3
+ 50 18 7
+ 122 41 2 1 7 3
+ 49 51 24
+ 123 41 2 1 7 3
+ 51 23 8
+ 124 41 2 1 7 3
+ 50 7 51
+ 125 41 2 1 7 3
+ 51 7 23
+ 126 41 2 1 7 3
+ 51 49 50
+ -1
=== added file 'scripts/test/unvRead.py'
--- scripts/test/unvRead.py 1970-01-01 00:00:00 +0000
+++ scripts/test/unvRead.py 2011-08-27 17:40:30 +0000
@@ -0,0 +1,11 @@
+from yade import *
+from yade import ymport
+
+facets = ymport.unv('shell.unv')
+O.bodies.append([f for f in facets])
+
+try:
+ import qt
+ qt.View()
+except:
+ pass
=== added file 'scripts/test/vtkExporter.py'
--- scripts/test/vtkExporter.py 1970-01-01 00:00:00 +0000
+++ scripts/test/vtkExporter.py 2011-08-27 17:40:30 +0000
@@ -0,0 +1,14 @@
+from yade import *
+from yade import utils,export
+
+O.bodies.append([
+ utils.sphere((0,0,0),1),
+ utils.sphere((0,2,0),1),
+ utils.sphere((0,2,3),2),
+ utils.facet([Vector3(0,-3,-1),Vector3(0,-2,5),Vector3(5,4,0)]),
+ utils.facet([Vector3(0,-3,-1),Vector3(0,-2,5),Vector3(-5,4,0)])
+])
+
+vtkExporter = export.VTKExporter('vtkExporterTesting')
+vtkExporter.exportSpheres(what=[('dist','b.state.pos.norm()')])
+vtkExporter.exportFacets(what=[('pos','b.state.pos')])