You can not select more than 25 topics
Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.
347 lines
12 KiB
347 lines
12 KiB
#include <array>
|
|
#include <bitset>
|
|
#include <cassert>
|
|
#include <iostream>
|
|
#include <booluarray.hpp>
|
|
|
|
|
|
#include <cstddef>
|
|
#include <iostream>
|
|
#include <iomanip>
|
|
#include <fstream>
|
|
#include <vector>
|
|
#include "bernstein.hpp"
|
|
#include "multiloop.hpp"
|
|
#include "quadrature_multipoly.hpp"
|
|
#include "binomial.hpp"
|
|
|
|
#include "real.hpp"
|
|
#include "sparkstack.hpp"
|
|
#include "uvector.hpp"
|
|
#include "vector"
|
|
#include "xarray.hpp"
|
|
#include <chrono>
|
|
#include <cmath>
|
|
#include <memory>
|
|
|
|
#include "organizer/primitive.hpp"
|
|
|
|
namespace algoim::Organizer
|
|
{
|
|
|
|
namespace detail
|
|
{
|
|
void bernstein2PowerTensor(const tensor3& phiBernstein, tensor3& phiPower) {}
|
|
} // namespace detail
|
|
|
|
bool keepQuadraturePoint(std::vector<tensor3>& originTensors, const uvector3& originPt)
|
|
{
|
|
// TODO: using blobtree to filter tensors
|
|
for (auto& t : originTensors) {
|
|
if (evalPower(t, originPt) >= 0) { return false; }
|
|
}
|
|
return true;
|
|
}
|
|
|
|
// std::vector<std::shared_ptr<PrimitiveDesc>> primitives;
|
|
|
|
// BasicTask(std::vector<std::shared_ptr<PrimitiveDesc>> ps) {};
|
|
void restrictToInnerFace(const tensor3& phi, int k, real facePos, tensor2& res)
|
|
{
|
|
assert(0 <= k && k < 3);
|
|
assert(all(res.ext() == remove_component(phi.ext(), k)));
|
|
assert(0 < facePos && facePos < 1);
|
|
real* basisAlongK;
|
|
int P = phi.ext(k);
|
|
algoim_spark_alloc(real, &basisAlongK, P);
|
|
bernstein::evalBernsteinBasis(facePos, P, basisAlongK);
|
|
|
|
for (auto i = res.loop(); ~i; ++i) {
|
|
real evalAlongKAtFacePos = 0;
|
|
for (int j = 0; j < P; ++j) { evalAlongKAtFacePos += basisAlongK[j] * phi.m(add_component(i(), k, j)); }
|
|
res.l(i) = evalAlongKAtFacePos;
|
|
}
|
|
}
|
|
|
|
void basicTask(const std::shared_ptr<PrimitiveDesc>& p, int q = 20, real xmin = -1, real xmax = 1)
|
|
{
|
|
real volume = 0;
|
|
auto integrand = [](const uvector<real, 3>& x) { return 1.0; };
|
|
uvector3 range = xmax - xmin;
|
|
|
|
|
|
if (auto pt = std::dynamic_pointer_cast<SphereDesc>(p)) {
|
|
tensor3 tensor(nullptr, 3);
|
|
algoim_spark_alloc(real, tensor);
|
|
makeSphere(*pt, tensor);
|
|
detail::powerTransformation(range, xmin, tensor);
|
|
|
|
tensor3 phi(nullptr, 3);
|
|
algoim_spark_alloc(real, phi);
|
|
detail::power2BernsteinTensor(tensor, phi);
|
|
|
|
uvector<real, 3> testX(0., 0., 0.25);
|
|
real testEvalBernstein = bernstein::evalBernsteinPoly(phi, testX);
|
|
// auto vec1 = xarray2StdVector(phi);
|
|
std::cout << "eval bernstein without interpolation:" << testEvalBernstein << std::endl;
|
|
|
|
ImplicitPolyQuadrature<3> ipquad(phi);
|
|
ipquad.integrate(AutoMixed, q, [&](const uvector<real, 3>& x, real w) {
|
|
if (isInsideBernstein(phi, x)) volume += w * integrand(xmin + x * (xmax - xmin));
|
|
});
|
|
} else if (auto pt = std::dynamic_pointer_cast<MeshDesc>(p)) {
|
|
const int faceCount = pt->indexInclusiveScan.size();
|
|
assert(faceCount > 1);
|
|
std::vector<tensor3> planeTensors(faceCount, tensor3(nullptr, 2));
|
|
algoim_spark_alloc(real, planeTensors);
|
|
|
|
tensor3 compositeTensor(nullptr, 1 + faceCount);
|
|
algoim_spark_alloc(real, compositeTensor);
|
|
|
|
makeMesh(*pt, compositeTensor, planeTensors);
|
|
detail::powerTransformation(range, xmin, compositeTensor);
|
|
|
|
auto planeStdVector1 = xarray2StdVector(planeTensors[0]);
|
|
auto planeStdVector2 = xarray2StdVector(planeTensors[1]);
|
|
auto compositeTensorStdVector = xarray2StdVector(compositeTensor);
|
|
|
|
uvector<real, 3> testX(0., 0.75, 0.2);
|
|
real textEvalPower = evalPower(compositeTensor, testX);
|
|
|
|
tensor3 phi(nullptr, 1 + faceCount);
|
|
algoim_spark_alloc(real, phi);
|
|
detail::power2BernsteinTensor(compositeTensor, phi);
|
|
|
|
int quadraturePointCount = 0;
|
|
ImplicitPolyQuadrature<3> ipquad(phi);
|
|
|
|
real testEvalBernstein = bernstein::evalBernsteinPoly(phi, testX);
|
|
ipquad.integrate(AutoMixed, q, [&](const uvector<real, 3>& x, real w) {
|
|
quadraturePointCount++;
|
|
auto realX = x * range + xmin;
|
|
if (isInsidePowers(planeTensors, realX)) volume += w * integrand(realX);
|
|
});
|
|
std::cout << "textEvalPower: " << textEvalPower << std::endl;
|
|
std::cout << "quadraturePointCount: " << quadraturePointCount << std::endl;
|
|
}
|
|
volume *= pow(xmax - xmin, 3);
|
|
std::cout << "Volume xxx: " << volume << std::endl;
|
|
};
|
|
|
|
void basicTask(const std::vector<std::shared_ptr<PrimitiveDesc>>& primitives, int q = 20, real xmin = -1, real xmax = 1)
|
|
{
|
|
std::vector<SparkStack<real>*> phiStacks;
|
|
std::vector<tensor3> phis;
|
|
|
|
std::vector<SparkStack<real>*> originTensorStacks;
|
|
std::vector<tensor3> originTensors;
|
|
real volume;
|
|
auto integrand = [](const uvector<real, 3>& x) { return 1.0; };
|
|
uvector3 range = xmax - xmin;
|
|
uvector<real, 3> testX(0., 0.75, 0.2);
|
|
for (int i = 0; i < primitives.size(); i++) {
|
|
if (auto pt = std::dynamic_pointer_cast<SphereDesc>(primitives[i])) {
|
|
tensor3 originTensor(nullptr, 3), transformedTensor(nullptr, 3);
|
|
originTensorStacks.emplace_back(algoim_spark_alloc_heap(real, originTensor)); // 记录,用以最后bool
|
|
|
|
tensor3 phi(nullptr, 3);
|
|
phiStacks.emplace_back(
|
|
algoim_spark_alloc_heap(real, phi)); // 必须先于algoim_spark_alloc,使transformedTensor的内存以栈形式释放
|
|
|
|
algoim_spark_alloc(real, transformedTensor);
|
|
makeSphere(*pt, originTensor);
|
|
originTensors.emplace_back(originTensor);
|
|
detail::powerTransformation(range, xmin, originTensor, transformedTensor);
|
|
|
|
|
|
detail::power2BernsteinTensor(transformedTensor, phi);
|
|
phis.emplace_back(phi);
|
|
} else if (auto pt = std::dynamic_pointer_cast<MeshDesc>(primitives[i])) {
|
|
const int faceCount = pt->indexInclusiveScan.size();
|
|
assert(faceCount > 1);
|
|
std::vector<tensor3> planeTensors(faceCount, tensor3(nullptr, 2));
|
|
algoimSparkAllocHeapVector(originTensorStacks, planeTensors);
|
|
|
|
tensor3 phi(nullptr, 1 + faceCount);
|
|
phiStacks.emplace_back(
|
|
algoim_spark_alloc_heap(real, phi)); // 必须先于algoim_spark_alloc,使compositeTensor的内存以栈形式释放
|
|
|
|
tensor3 compositeTensor(nullptr, 1 + faceCount);
|
|
algoim_spark_alloc(real, compositeTensor);
|
|
|
|
makeMesh(*pt, compositeTensor, planeTensors);
|
|
detail::powerTransformation(range, xmin, compositeTensor);
|
|
|
|
real testEvalPower = evalPower(compositeTensor, testX);
|
|
|
|
originTensors.insert(originTensors.end(), planeTensors.begin(), planeTensors.end());
|
|
|
|
detail::power2BernsteinTensor(compositeTensor, phi);
|
|
real testEvalBernstein = bernstein::evalBernsteinPoly(phi, testX);
|
|
phis.emplace_back(phi);
|
|
}
|
|
}
|
|
real testEvalBernstein = bernstein::evalBernsteinPoly(phis[0], testX);
|
|
ImplicitPolyQuadrature<3> ipquad(phis);
|
|
|
|
ipquad.integrate(AutoMixed, q, [&](const uvector<real, 3>& x, real w) {
|
|
auto realX = x * range + xmin;
|
|
if (keepQuadraturePoint(originTensors, realX)) volume += w * integrand(realX);
|
|
});
|
|
volume *= pow(xmax - xmin, 3);
|
|
std::cout << "Volume xxx: " << volume << std::endl;
|
|
|
|
// free memory, thus deallocate memory of xarray
|
|
for (auto& p : phiStacks) delete p;
|
|
for (auto& p : originTensorStacks) delete p;
|
|
};
|
|
|
|
void quadratureTask(const Scene& scene) {}
|
|
|
|
void buildOctree(const Scene& scene, std::vector<Node>& nodes, const uvector3& min, const uvector3& max) {}
|
|
|
|
template <int N>
|
|
int numOfZero(const uvector<int, N>& v)
|
|
{
|
|
int res = 0;
|
|
for (int i = 0; i < N; ++i) {
|
|
if (v(i) == 0) { res++; }
|
|
}
|
|
return res;
|
|
}
|
|
|
|
template <int N>
|
|
int binaryToDecimal(const uvector<int, N>& v, int zeroRep = 0)
|
|
{
|
|
int res = 0;
|
|
int base = 1;
|
|
for (int i = N - 1; i >= 0; --i) {
|
|
if (v(i) != zeroRep) { res += base; }
|
|
base *= 2;
|
|
}
|
|
return res;
|
|
}
|
|
|
|
template <int N>
|
|
int replaceFirst(uvector<int, N>& v, int oldVal, int newVal)
|
|
{
|
|
for (int i = 0; i < N; ++i) {
|
|
if (v(i) == oldVal) {
|
|
v(i) = newVal;
|
|
return i;
|
|
}
|
|
}
|
|
return -1;
|
|
}
|
|
|
|
void build(const Scene& scene, std::vector<Node>& nodes, int nowNodeIdx, const uvector3& nowNodeMin, const uvector3& nowNodeMax)
|
|
{
|
|
const std::vector<int>& polyIntersectIndices = nodes[nowNodeIdx].polyIntersectIndices;
|
|
if (polyIntersectIndices.size() == 1) {
|
|
// TODO 相交很少时
|
|
}
|
|
int lastIdx = nodes.size();
|
|
nodes.resize(lastIdx + 8);
|
|
uvector3 nodeMid = (nowNodeMin + nowNodeMax) / 2;
|
|
for (int i = 0; i < polyIntersectIndices.size(); ++i) {
|
|
auto poly = scene.polys[polyIntersectIndices[i]];
|
|
uvector<int, 3> mark(0, 0, 0);
|
|
for (int faceAxis = 0; faceAxis < 3; ++faceAxis) {
|
|
real centerPlane = nodeMid(faceAxis);
|
|
tensor2 restrictToCenterPlane(nullptr, remove_component(poly.ext(), faceAxis));
|
|
algoim_spark_alloc(real, restrictToCenterPlane);
|
|
restrictToInnerFace(poly, faceAxis, centerPlane, restrictToCenterPlane);
|
|
if (bernstein::uniformSign(restrictToCenterPlane) == -1) {
|
|
// primitive contain the conterface
|
|
mark(faceAxis) = 2;
|
|
} else {
|
|
// primitive is on one side of the centerface
|
|
// deCasteljau to transformation
|
|
uvector3 halfCellMin = 0, halfCellMax = 1;
|
|
halfCellMax(faceAxis) = 0.5;
|
|
tensor3 negativeHalfCell(nullptr, poly.ext());
|
|
algoim_spark_alloc(real, negativeHalfCell);
|
|
bernstein::deCasteljau(poly, halfCellMin, halfCellMax, negativeHalfCell);
|
|
if (bernstein::uniformSign(negativeHalfCell) == -1) {
|
|
// 不在负半空间,此时一定在正半空间,可以assert一下正半空间的sign,不能为-1
|
|
mark(faceAxis) = 1;
|
|
} else {
|
|
mark(faceAxis) = -1;
|
|
}
|
|
}
|
|
}
|
|
if (any(mark == uvector3(2, 2, 2))) {
|
|
int subIdx = 0;
|
|
for (MultiLoop<3> j(0, 2); ~j; ++j, ++subIdx) {
|
|
nodes[lastIdx + subIdx].polyIntersectIndices.emplace_back(poly);
|
|
nodes[lastIdx + subIdx].min = nowNodeMin;
|
|
nodes[lastIdx + subIdx].max = nowNodeMax;
|
|
for (int k = 0; k < 3; ++k) {
|
|
if (j(k) == 0)
|
|
nodes[lastIdx + subIdx].max(k) = nodeMid(k);
|
|
else
|
|
nodes[lastIdx + subIdx].min(k) = nodeMid(k);
|
|
}
|
|
}
|
|
continue;
|
|
}
|
|
int zeroCount = numOfZero<3>(mark);
|
|
if (zeroCount == 0) {
|
|
// mark has -1 or 1 only
|
|
nodes[lastIdx + binaryToDecimal(mark, -1)].polyIntersectIndices.emplace_back(polyIntersectIndices[i]);
|
|
} else if (zeroCount == 1) {
|
|
// poly related to 2 subcells
|
|
uvector<int, 3> sideMark = mark;
|
|
int zeroIdx = replaceFirst(sideMark, 0, 1);
|
|
nodes[lastIdx + binaryToDecimal(sideMark, -1)].polyIntersectIndices.emplace_back(polyIntersectIndices[i]);
|
|
sideMark(zeroIdx) = -1;
|
|
nodes[lastIdx + binaryToDecimal(sideMark, -1)].polyIntersectIndices.emplace_back(polyIntersectIndices[i]);
|
|
} else {
|
|
// zeroCount == 2 or 3
|
|
}
|
|
}
|
|
// TODO: 考虑fully contain 问题
|
|
}
|
|
}; // namespace algoim::Organizer
|
|
|
|
// clang-format off
|
|
/**
|
|
最一开始的cell应当包含所有的primitive
|
|
func(given cell) {
|
|
得到8个subcell
|
|
for each primitive {
|
|
mark = [0,0,0]
|
|
for k in 0,1,2 {
|
|
测试primitive和k轴向的centerface有无交
|
|
if 无交 {
|
|
if primitive包裹centerface {
|
|
mark[k] = 2
|
|
} else {
|
|
判断primitive位于哪一侧
|
|
if k正方向 mark[k] = 1
|
|
else mark[k] = -1
|
|
}
|
|
}
|
|
}
|
|
if mark 含2 {
|
|
// TODO: 可能需要优化:在极端情况下,
|
|
// 例如primitive的边刚好与centerface重合,
|
|
// 此时不是所有subcell都关联primitive
|
|
将primitive放入所有subcell
|
|
} elif mark 没有0 {
|
|
将primitive放入对应的subcell
|
|
} elif mark 有1个0 {
|
|
将primitive放入两个对应的subcell
|
|
} else { // >= 2个0
|
|
测试每个可能有交subcell是否与primitive相交
|
|
相交的放入对于subcell
|
|
}
|
|
}
|
|
for each subcell {
|
|
func(subcell)
|
|
}
|
|
}
|
|
|
|
*/
|
|
|
|
// clang-format on
|