Skip to content
Snippets Groups Projects
Commit 1ab0d79f authored by Adrien Crovato's avatar Adrien Crovato
Browse files

Improve MshDeform init. Private members

parent 6200ea36
No related branches found
No related tags found
No related merge requests found
......@@ -112,9 +112,11 @@ namespace std {
%shared_ptr(tbox::MshImport)
%nodefault tbox::MshImport;
%immutable tbox::MshImport::msh; // read-only variable (no setter)
%include "wMshImport.h"
%shared_ptr(tbox::MshExport)
%nodefault tbox::MshExport;
%immutable tbox::MshExport::msh; // read-only variable (no setter)
%include "wMshExport.h"
%shared_ptr(tbox::GmshImport)
%include "wGmshImport.h"
......@@ -140,6 +142,8 @@ namespace std {
%include "wGroups.h"
%warnfilter(+509);
%shared_ptr(tbox::MshDeform);
%immutable tbox::MshDeform::msh; // read-only variable (no setter)
%immutable tbox::MshDeform::nthreads;
%include "wMshDeform.h" // included after std templates because it needs them
%shared_ptr(tbox::MshCrack);
%include "wMshCrack.h"
......
......@@ -37,6 +37,8 @@ namespace tbox
class TBOX_API MshCrack : public fwk::wSharedObject
{
private:
std::shared_ptr<MshData> msh; ///< mesh
int nDim; ///< number of dimensions
bool hasCrack; ///< flag telling if mesh already has crack
std::map<Node *, Node *> nodMap; ///< map between nodes on both side of crack
Group *crckGrp; ///< physical group of crack
......@@ -50,9 +52,6 @@ private:
void swapNodes(Element *e);
public:
std::shared_ptr<MshData> msh;
int nDim;
MshCrack(std::shared_ptr<MshData> _msh, int _nDim);
virtual ~MshCrack();
......
......@@ -32,7 +32,10 @@
using namespace tbox;
MshDeform::MshDeform(std::shared_ptr<MshData> _msh, int _nDim) : wSharedObject(), nDim(_nDim), field(false), fixed(false), moving(false), msh(_msh)
MshDeform::MshDeform(std::shared_ptr<MshData> _msh,
int _nDim, int nthrds) : wSharedObject(), nDim(_nDim),
field(false), fixed(false), moving(false),
msh(_msh), nthreads(nthrds)
{
// Check problem dimension
if (nDim != 2 && nDim != 3)
......@@ -41,10 +44,13 @@ MshDeform::MshDeform(std::shared_ptr<MshData> _msh, int _nDim) : wSharedObject()
fldGrp = nullptr;
symBndGrp = nullptr;
// Initialize multithreading
nthreads = 1;
// Get mesh size
mshSize = msh->nodes.size();
// Display configuration
std::cout << "--- Mesh morpher (linear elasticity) ---\n"
<< "Number of threads: " << nthreads << "\n"
<< std::endl;
}
/**
......@@ -449,14 +455,13 @@ void MshDeform::deform()
Eigen::VectorXd f = Eigen::VectorXd::Zero(nDim * mshSize), q = Eigen::VectorXd::Zero(nDim * mshSize);
// Build stiffness matrix and RHS
std::cout << "--- Deforming mesh using linear elasticiy equations ---\n"
<< "Number of threads: " << nthreads << "\n"
<< "Assembling matrices... " << std::flush;
std::cout << "- Deforming mesh:\n"
<< "assembling matrices... " << std::flush;
this->build(K);
this->build(f);
std::cout << "done" << std::endl;
// Solve the SoE
std ::cout << "Solving equations... " << std::flush;
std ::cout << "solving equations... " << std::flush;
Gmres gmres;
// set the preconditioner
gmres.setFillFactor(2);
......@@ -471,7 +476,7 @@ void MshDeform::deform()
std::cout << "done (#it: " << gmres.getIterations() << ", error: " << std::scientific << gmres.getError() << std::fixed << ")" << std::endl;
// Update position
std ::cout << "Updating mesh... " << std::flush;
std ::cout << "updating mesh... " << std::flush;
for (auto n : msh->nodes)
{
for (int m = 0; m < nDim; m++)
......@@ -492,7 +497,7 @@ void MshDeform::deform()
e->update();
for (auto e : intBndElems)
e->update();
std::cout << "done" << std::endl;
std::cout << "done\n" << std::endl;
}
/**
......
......@@ -63,7 +63,7 @@ public:
std::shared_ptr<MshData> msh;
int nthreads;
MshDeform(std::shared_ptr<MshData> _msh, int _nDim);
MshDeform(std::shared_ptr<MshData> _msh, int _nDim, int nthrds = 1);
virtual ~MshDeform() { std::cout << "~MshDeform()\n"; }
void setField(std::string const &fldName);
......
......@@ -39,13 +39,11 @@ def main():
msh = gmsh.MeshLoader("models/beam.geo",__file__).execute(**pars)
gmshWriter = tbox.GmshExport(msh)
# create mesh deformation handler
mshDef = tbox.MshDeform(msh, 2)
mshDef = tbox.MshDeform(msh, 2, nthrds=parseargs().k)
mshDef.setField("internalField")
mshDef.addFixed(["clamp"])
mshDef.addMoving(["tip"])
mshDef.initialize()
args = parseargs()
mshDef.nthreads = args.k
# deform the mesh (dummy translation of "tip" nodes)
mshDef.savePos()
......
......@@ -38,13 +38,11 @@ def main():
msh = gmsh.MeshLoader("models/beam3.geo",__file__).execute(**pars)
gmshWriter = tbox.GmshExport(msh)
# create mesh deformation handler
mshDef = tbox.MshDeform(msh, 3)
mshDef = tbox.MshDeform(msh, 3, nthrds=parseargs().k)
mshDef.setField("field")
mshDef.addFixed(["clamp"])
mshDef.addMoving(["tip"])
mshDef.initialize()
args = parseargs()
mshDef.nthreads = args.k
# deform the mesh (dummy translation of "tip" nodes)
mshDef.savePos()
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment