int main(int argc, char *argv[])
{
real_t refParameter;
bool plot = false;
bool manual = false;
bool HB = false;
bool dump;
RefineLoopMax = 2;
initUnifRef = 2;
degree = 2;
refCriterion = PUCA;
refParameter = 0.85;
dump = false;
gsCmdLine cmd(
"Solving a PDE with adaptive refinement using THB-splines.");
cmd.addSwitch("plot", "Plot resulting mesh in ParaView", plot);
cmd.addInt("r", "refine", "Maximum number of adaptive refinement steps to perform",
RefineLoopMax);
cmd.addInt("i", "initial-ref", "Initial number of uniform refinement steps to perform",
initUnifRef);
cmd.addInt("", "degree", "Spline degree of the THB basis", degree);
cmd.addInt("c", "criterion", "Criterion to be used for adaptive refinement (1-3, see documentation)",
refCriterion);
cmd.addInt("E", "extension", "Criterion to be used for adaptive refinement (1-3, see documentation)",
refExt);
cmd.addReal("p", "parameter", "Parameter for adaptive refinement", refParameter);
cmd.addSwitch("dump", "Write geometry and sequence of bases into XML files",
dump);
cmd.addSwitch("manual", "Uses the 'addLevel' feature of the THB basis, where levels can be specified manually", manual);
cmd.addSwitch("HB", "Uses a Hierarchical basis instead of a Truncated Hierarchical basis", HB);
try { cmd.getValues(argc,argv); } catch (int rv) { return rv; }
gsInfo<<
"Source function "<< f <<
"\n";
gsInfo<<
"Exact solution "<< g <<
".\n" <<
"\n";
gsInfo <<
" --- Geometry:\n" << *geo <<
"\n";
gsInfo <<
"Number of patches: " << patches.nPatches() <<
"\n";
if (dump)
gsInfo <<
"\nCoarse discretization basis:\n" << tbb <<
"\n";
for (int i = 0; i < initUnifRef; ++i)
{
if (manual)
tbb.uniformRefine(1,1);
else
tbb.uniformRefine();
}
if (manual)
{
if (HB)
else
for (
index_t k=0; k!=degree-2; k++)
{
tbb.uniformRefine(1,k+1);
tbb.reduceContinuity(1);
tbb.removeKnot(0.5,0,1);
}
{
tbb.uniformRefine(1,degree-1);
}
}
else
{
if (HB)
else
}
boxes.row(0)<<0.25,0.75;
boxes.row(1)<<0.00,1.00;
pa.options().setInt("DirichletValues", dirichlet::l2Projection);
if (dump)
gsWrite(bases[0],
"adapt_basis_0.xml");
for( int RefineLoop = 1; RefineLoop <= RefineLoopMax ; RefineLoop++ )
{
gsInfo <<
"\n============================== Loop " << RefineLoop <<
" of " << RefineLoopMax <<
" ==============================" <<
"\n" <<
"\n";
std::string unity_string = unity ? "has " : "does not have ";
gsInfo<<
" * Number of functions: "<<HTB_tmp->
size()<<
"\n";
gsInfo<<
" * Size per level: ";
for(
unsigned i = 0; i<= HTB_tmp->
maxLevel(); i++)
gsInfo << HTB_tmp->getXmatrix()[i].size()<<
" ";
gsInfo<<
" * Partition of unity: The basis "<<unity_string<<
"the partition of unity property\n";
gsInfo<<
"=========================================================================\n";
gsInfo <<
"Assembling... " << std::flush;
pa.assemble();
gsInfo <<
"Solving... " << std::flush;
gsMatrix<> solVector = gsSparseSolver<>::CGDiagonal(pa.matrix() ).solve( pa.rhs() );
pa.constructSolution(solVector, mpsol);
gsExprEvaluator<>::geometryMap Gm = ev.
getMap(patches);
const std::vector<real_t> & elErrEst = ev.
elementwise();
std::vector<bool> elMarked( elErrEst.size() );
gsInfo <<
"Marked "<< std::count(elMarked.begin(), elMarked.end(),
true)<<
"\n";
pa.refresh();
if (dump)
{
std::stringstream ss;
ss << "adapt_basis_" << RefineLoop << ".xml";
gsWrite(pa.multiBasis()[0], ss.str());
}
if ( (RefineLoop == RefineLoopMax) && plot)
{
gsInfo<<
"Plotting in Paraview...\n";
gsWriteParaview<>(sol, "p2d_adaRef_sol", 5001, true);
gsWriteParaview<>(pa.multiBasis()[0], "basis", 500, true);
}
}
gsInfo <<
"\nFinal basis: " << pa.multiBasis()[0] <<
"\n";
delete HTB;
return EXIT_SUCCESS;
}
Definition gsExpressions.h:928
void setDegree(short_t const &i)
Set the degree of the basis (either elevate or reduce) in order to have degree equal to i wrt to each...
Definition gsBasis.hpp:645
Class containing a set of boundary conditions.
Definition gsBoundaryConditions.h:342
void addCondition(int p, boxSide s, condition_type::type t, gsFunctionSet< T > *f, short_t unknown=0, bool parametric=false, int comp=-1)
Adds another boundary condition.
Definition gsBoundaryConditions.h:650
Class for command-line argument parsing.
Definition gsCmdLine.h:57
Generic evaluator of isogeometric expressions.
Definition gsExprEvaluator.h:39
void setIntegrationElements(const gsMultiBasis< T > &mesh)
Sets the domain of integration.
Definition gsExprEvaluator.h:110
const std::vector< T > & elementwise() const
Returns an std::vector containing the last computed values per element.
Definition gsExprEvaluator.h:99
T integralElWise(const expr::_expr< E > &expr)
Calculates the integral of the expression expr on each element.
Definition gsExprEvaluator.h:159
variable getVariable(const gsFunctionSet< T > &func, index_t dim=1)
Registers func as a variable and returns a handle to it.
Definition gsExprEvaluator.h:124
geometryMap getMap(const gsMultiPatch< T > &mp)
Registers mp as an isogeometric geometry map and return a handle to it.
Definition gsExprEvaluator.h:116
A scalar of vector field defined on a m_parametric geometry.
Definition gsField.h:55
static void open(const std::string &fn)
Opens the file fn using the preferred application of the OS.
Definition gsFileManager.cpp:688
Class defining a multivariate (real or vector) function given by a string mathematical expression.
Definition gsFunctionExpr.h:52
virtual const gsBasis< T > & basis() const =0
Returns a const reference to the basis of the geometry.
A hierarchical B-spline basis of parametric dimension d.
Definition gsHBSplineBasis.h:36
Class representing a (scalar) hierarchical tensor basis of functions .
Definition gsHTensorBasis.h:75
int treeSize() const
The number of nodes in the tree representation.
Definition gsHTensorBasis.h:647
virtual void refine(gsMatrix< T > const &boxes, int refExt)
Refine the basis to levels and in the areas defined by boxes with an extension.
Definition gsHTensorBasis.hpp:445
void printBases(std::ostream &os=gsInfo) const
Prints the spline-space hierarchy.
Definition gsHTensorBasis.h:542
void addLevel(const gsTensorBSplineBasis< d, T > &next_basis)
Adds a level, only if manual levels are activated.
Definition gsHTensorBasis.hpp:35
index_t size() const
The number of basis functions in this basis.
Definition gsHTensorBasis.hpp:298
bool testPartitionOfUnity(const index_t npts=100, const T tol=1e-12) const
Test the partition of unity.
Definition gsHTensorBasis.hpp:1972
unsigned maxLevel() const
Returns the level in which the indices are stored internally.
Definition gsHTensorBasis.h:811
size_t numElements(boxSide const &s=0) const
The number of elements on side s.
Definition gsHTensorBasis.h:972
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
Holds a set of patch-wise bases and their topology information.
Definition gsMultiBasis.h:37
Container class for a set of geometry patches and their topology, that is, the interface connections ...
Definition gsMultiPatch.h:100
Implementation of an (multiple right-hand side) Poisson assembler.
Definition gsPoissonAssembler.h:37
Truncated hierarchical B-spline basis.
Definition gsTHBSplineBasis.h:36
A tensor product B-spline basis.
Definition gsTensorBSplineBasis.h:37
A tensor product of d B-spline functions, with arbitrary target dimension.
Definition gsTensorBSpline.h:45
Main header to be included by clients using the G+Smo library.
void gsMarkElementsForRef(const std::vector< T > &elError, int refCriterion, T refParameter, std::vector< bool > &elMarked)
Marks elements/cells for refinement.
Definition gsAdaptiveRefUtils.h:196
void gsRefineMarkedElements(gsMultiBasis< T > &basis, const std::vector< bool > &elMarked, index_t refExtension=0)
Refine a gsMultiBasis, based on a vector of element-markings.
Definition gsAdaptiveRefUtils.h:294
Provides class for adaptive refinement.
#define index_t
Definition gsConfig.h:32
#define gsInfo
Definition gsDebug.h:43
EIGEN_STRONG_INLINE meas_expr< T > meas(const gsGeometryMap< T > &G)
The measure of a geometry map.
Definition gsExpressions.h:4555
The G+Smo namespace, containing all definitions for the library.
void gsWrite(const Object &obj, const std::string &fname)
Write an arbitrary Gismo object to an XML file with the given filename.
Definition gsReadFile.h:288
@ dirichlet
Dirichlet type.
Definition gsBoundaryConditions.h:31
Class gsNurbsCreator provides some simple examples of Nurbs Geometries.
Definition gsNurbsCreator.h:37