Skip to content
Snippets Groups Projects
Forked from agnumpde / dune-tectonic
159 commits ahead of the upstream repository.
Code owners
Assign users and groups as approvers for specific file changes. Learn more.
levelcontactnetwork.cc 8.71 KiB
#ifdef HAVE_CONFIG_H
#include "config.h"
#endif

#include <dune/contact/projections/normalprojection.hh>

#include <dune/fufem/assemblers/localassemblers/neumannboundaryassembler.hh>
#include <dune/fufem/functions/constantfunction.hh>
#include <dune/fufem/facehierarchy.hh>

#include "../../assemblers.hh"
#include "../../utils/tobool.hh"
#include "../../utils/debugutils.hh"

#include "../friction/globalratestatefriction.hh"
#include "../friction/frictionpotential.hh"

#include "levelcontactnetwork.hh"

template <class GridType, class FrictionCouplingPair, class field_type>
LevelContactNetwork<GridType, FrictionCouplingPair, field_type>::LevelContactNetwork(
        int nBodies,
        int nCouplings,
        int level) :
    level_(level),
    bodies_(nBodies),
    couplings_(nCouplings),
    nonmortarPatches_(nCouplings),
    mortarPatches_(nCouplings),
    glues_(nCouplings)
{}


template <class GridType, class FrictionCouplingPair, class field_type>
void LevelContactNetwork<GridType, FrictionCouplingPair, field_type>::setBodies(const std::vector<std::shared_ptr<Body>> bodies) {
    assert(nBodies() == bodies.size());
    bodies_ = bodies;
}

template <class GridType, class FrictionCouplingPair, class field_type>
void LevelContactNetwork<GridType, FrictionCouplingPair, field_type>::setCouplings(const std::vector<std::shared_ptr<FrictionCouplingPair>> couplings) {
    assert(nCouplings() == couplings.size());
    couplings_ = couplings;
}

template <class GridType, class FrictionCouplingPair, class field_type>
void LevelContactNetwork<GridType, FrictionCouplingPair, field_type>::constructBody(
        const std::shared_ptr<BodyData<dimworld>>& bodyData,
        const std::shared_ptr<GridType>& grid,
        const size_t level,
        std::shared_ptr<Body>& body) const {

    body = std::make_shared<Body>(bodyData, grid, level);
}

// collects all leaf boundary nodes from the different bodies identified by "tag" in a single BitSetVector
template <class GridType, class FrictionCouplingPair, class field_type>
void LevelContactNetwork<GridType, FrictionCouplingPair, field_type>::totalNodes(
        const std::string& tag,
        Dune::BitSetVector<dimworld>& totalNodes) const {

    int totalSize = 0;
    for (size_t i=0; i<nBodies(); i++) {
        totalSize += this->body(i)->nVertices();
    }

    totalNodes.resize(totalSize);

    int idx=0;
    for (size_t i=0; i<nBodies(); i++) {
        const auto& body = this->body(i);
        std::vector<std::shared_ptr<typename Body::BoundaryCondition>> boundaryConditions;
        body->boundaryConditions(tag, boundaryConditions);

        const int idxBackup = idx;
        for (size_t bc = 0; bc<boundaryConditions.size(); bc++) {
            const auto& nodes = boundaryConditions[bc]->boundaryNodes();
            for (size_t j=0; j<nodes->size(); j++, idx++)
                for (int k=0; k<dimworld; k++)
                    totalNodes[idx][k] = (*nodes)[j][k];
            idx = (bc==boundaryConditions.size()-1 ? idx : idxBackup);
        }
    }
}

// collects all leaf boundary nodes from the different bodies identified by "tag" in a vector of BitSetVector
template <class GridType, class FrictionCouplingPair, class field_type>
void LevelContactNetwork<GridType, FrictionCouplingPair, field_type>::boundaryPatchNodes(
        const std::string& tag,
        BoundaryPatchNodes& nodes) const {

    nodes.resize(nBodies());

    for (size_t i=0; i<nBodies(); i++) {
        this->body(i)->boundaryPatchNodes(tag, nodes[i]);
    }
}

// collects all leaf boundary nodes from the different bodies identified by "tag" in a vector of BitSetVector
template <class GridType, class FrictionCouplingPair, class field_type>
void LevelContactNetwork<GridType, FrictionCouplingPair, field_type>::boundaryNodes(
        const std::string& tag,
        BoundaryNodes& nodes) const {

    nodes.resize(nBodies());

    for (size_t i=0; i<nBodies(); i++) {
        this->body(i)->boundaryNodes(tag, nodes[i]);
    }
}

// collects all leaf boundary functions from the different bodies identified by "tag"
template <class GridType, class FrictionCouplingPair, class field_type>
void LevelContactNetwork<GridType, FrictionCouplingPair, field_type>::boundaryFunctions(
        const std::string& tag,
        BoundaryFunctions& functions) const {

    functions.resize(nBodies());

    for (size_t i=0; i<nBodies(); i++) {
        this->body(i)->boundaryFunctions(tag, functions[i]);
    }
}

// collects all leaf boundary patches from the different bodies identified by "tag"
template <class GridType, class FrictionCouplingPair, class field_type>
void LevelContactNetwork<GridType, FrictionCouplingPair, field_type>::boundaryPatches(
        const std::string& tag,
        BoundaryPatches& patches) const {

    patches.resize(nBodies());

    for (size_t i=0; i<nBodies(); i++) {
        this->body(i)->boundaryPatches(tag, patches[i]);
    }
}

template <class GridType, class FrictionCouplingPair, class field_type>
auto LevelContactNetwork<GridType, FrictionCouplingPair, field_type>::level() const
-> int {

    return level_;
}

template <class GridType, class FrictionCouplingPair, class field_type>
auto LevelContactNetwork<GridType, FrictionCouplingPair, field_type>::nBodies() const
-> size_t {

    return bodies_.size();
}

template <class GridType, class FrictionCouplingPair, class field_type>
auto LevelContactNetwork<GridType, FrictionCouplingPair, field_type>::nCouplings() const
-> size_t {

    return couplings_.size();
}

template <class GridType, class FrictionCouplingPair, class field_type>
auto LevelContactNetwork<GridType, FrictionCouplingPair, field_type>::body(int i) const
-> const std::shared_ptr<Body>& {

    return bodies_[i];
}

template <class GridType, class FrictionCouplingPair, class field_type>
auto LevelContactNetwork<GridType, FrictionCouplingPair, field_type>::coupling(int i) const
-> const std::shared_ptr<FrictionCouplingPair>& {

    return couplings_[i];
}


template <class GridType, class FrictionCouplingPair, class field_type>
auto LevelContactNetwork<GridType, FrictionCouplingPair, field_type>::couplings() const
-> const std::vector<std::shared_ptr<FrictionCouplingPair>>& {

    return couplings_;
}

template <class GridType, class FrictionCouplingPair, class field_type>
void LevelContactNetwork<GridType, FrictionCouplingPair, field_type>::prolong(const BoundaryPatch& coarsePatch, BoundaryPatch& finePatch, const size_t fineLevel) {

    if (not(coarsePatch.isInitialized()))
        DUNE_THROW(Dune::Exception, "Coarse boundary patch has not been set up correctly!");

    const GridType& grid = coarsePatch.gridView().grid();

    // Set array sizes correctly
    finePatch.setup(grid.levelGridView(fineLevel));

    for (const auto& pIt : coarsePatch) {
        const auto& inside = pIt.inside();

        if (inside.level() == fineLevel)
            finePatch.addFace(inside, pIt.indexInInside());
        else {
            Face<GridType> face(grid, inside, pIt.indexInInside());

            typename Face<GridType>::HierarchicIterator it = face.hbegin(fineLevel);
            typename Face<GridType>::HierarchicIterator end = face.hend(fineLevel);

            for(; it!=end; ++it) {
                if (it->element().level() == fineLevel)
                    finePatch.addFace(it->element(), it->index());
            }
        }
    }
}
template <class GridType, class FrictionCouplingPair, class field_type>
void LevelContactNetwork<GridType, FrictionCouplingPair, field_type>::build(field_type overlap) {

    for (size_t i=0; i<couplings_.size(); i++) {
        auto& coupling = *couplings_[i];

        const size_t nmBody = coupling.gridIdx_[0];
        const size_t mBody = coupling.gridIdx_[1];

        this->prolong(*coupling.patch0(), nonmortarPatches_[i], body(nmBody)->level());
        this->prolong(*coupling.patch1(), mortarPatches_[i], body(mBody)->level());

        using Element = typename GridView::template Codim<0>::Entity;
        auto desc0 = [&] (const Element& e, unsigned int face) {
            return nonmortarPatches_[i].contains(e,face);
        };
        auto desc1 = [&] (const Element& e, unsigned int face) {
            return mortarPatches_[i].contains(e,face);
        };

        auto extract0 = std::make_shared<Extractor>(body(nmBody)->gridView(), desc0);
        auto extract1 = std::make_shared<Extractor>(body(mBody)->gridView(), desc1);

        gridGlueBackend_ = std::make_shared< Dune::GridGlue::ContactMerge<dim, ctype> >(overlap);
        glues_[i] = std::make_shared<Glue>(extract0, extract1, gridGlueBackend_);
        glues_[i]->build();
    }
}

template <class GridType, class FrictionCouplingPair, class field_type>
auto LevelContactNetwork<GridType, FrictionCouplingPair, field_type>::glue(int i) const
-> const std::shared_ptr<Glue>& {

    return glues_[i];
}

#include "levelcontactnetwork_tmpl.cc"