Skip to content
Snippets Groups Projects
Forked from agnumpde / dune-tectonic
92 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++) {