Skip to content
Snippets Groups Projects
Select Git revision
  • 9e5dad4ee6c3d99963a563a165a0598c5897bcdc
  • master default protected
  • porenetwork-bachelor
  • preconditioner
  • release_cleanup
  • dd-2f1s-comparision
  • 2p_cahnhilliard_equation_a
  • example_2p_cahnhilliard_A
  • porenetwork
  • release/1.0 protected
10 results

fff_chns.hh

Blame
  • fff_chns.hh 7.09 KiB
    #pragma once
    
    #include<dune/pdelab/constraints/common/constraints.hh>
    #include<dune/pdelab/constraints/conforming.hh>
    
    #include<dune/pdelab/gridfunctionspace/gridfunctionspace.hh>
    #include<dune/pdelab/gridfunctionspace/subspace.hh>
    #include<dune/pdelab/gridfunctionspace/vectorgridfunctionspace.hh>
    #include<dune/phasefield/utility/pf_3p_color_dgf.hh>
    
    
    //_________________________________________________________________________________________________________________________
    
    
    template<typename RF, typename ES, typename V_FEM, typename P_FEM, typename PHI_FEM>
    class GFS_fff_chns
    {
    public:
      typedef Dune::PDELab::ConformingDirichletConstraints ConstraintsAssembler;
      typedef Dune::PDELab::ISTL::VectorBackend<> VectorBackend;
      typedef Dune::PDELab::ISTL::VectorBackend<Dune::PDELab::ISTL::Blocking::fixed> VB_Block;
      typedef Dune::PDELab::LexicographicOrderingTag OrderingTagL;
      typedef Dune::PDELab::EntityBlockedOrderingTag OrderingTagB;
    
      ConstraintsAssembler conp, conphi, conmu;
    
      typedef Dune::PDELab::VectorGridFunctionSpace<ES,V_FEM,ES::GridView::dimension,
                                                    VectorBackend,VectorBackend,ConstraintsAssembler> V_GFS;
      V_GFS vGfs;
    
      typedef Dune::PDELab::GridFunctionSpace<ES, P_FEM, ConstraintsAssembler, VectorBackend> P_GFS;
      P_GFS pGfs;
    
      typedef Dune::PDELab::GridFunctionSpace<ES, PHI_FEM, ConstraintsAssembler, VectorBackend> PHI_GFS;
      PHI_GFS phiGfs;
    
      typedef Dune::PDELab::GridFunctionSpace<ES, PHI_FEM, ConstraintsAssembler, VectorBackend> MU_GFS;
      MU_GFS muGfs;
    
      typedef Dune::PDELab::CompositeGridFunctionSpace<VectorBackend, OrderingTagL, V_GFS, P_GFS> NS;
      NS ns;
    
      typedef Dune::PDELab::CompositeGridFunctionSpace<VB_Block, OrderingTagB, PHI_GFS, MU_GFS, PHI_GFS, MU_GFS> CH;
      CH ch;
    
      typedef typename NS::template ConstraintsContainer<RF>::Type NS_CC;
      NS_CC nsCC;
    
      typedef typename CH::template ConstraintsContainer<RF>::Type CH_CC;
      CH_CC chCC;
    
      GFS_fff_chns(const ES& es, const V_FEM& vFem, const P_FEM& pFem, const PHI_FEM& phiFem) :
          conp(),conphi(),conmu(),
          vGfs(es,vFem),
          pGfs(es,pFem, conp),
          phiGfs(es,phiFem, conphi),
          muGfs(es,phiFem, conmu),
          ns(vGfs, pGfs),
          ch(phiGfs, muGfs, phiGfs, muGfs),
          nsCC(),
          chCC()
      {
        vGfs.name("velocity");
        pGfs.name("pressure");
        phiGfs.name("phasefield");
        muGfs.name("potential");
      }
    
      template<typename Bdry>
      void updateConstraintsContainer(const Bdry& bdry)
      {
        Dune::PDELab::constraints(bdry.bdryNS,ns,nsCC, 0);
        Dune::PDELab::constraints(bdry.bdryCH,ch,chCC, 0);
      }
    };
    
    //_________________________________________________________________________________________________________________________
    
    
    template<typename GFS, typename NS_V, typename CH_V>
    class DGF_fff_chns
    {
    public:
    
      typedef Dune::PDELab::GridFunctionSubSpace<typename GFS::NS,Dune::TypeTree::TreePath<0> > VSubGFS;
      VSubGFS vSubGfs;
      typedef Dune::PDELab::VectorDiscreteGridFunction<VSubGFS,NS_V> VDGF;
      VDGF vDgf;
    
      typedef Dune::PDELab::GridFunctionSubSpace<typename GFS::NS,Dune::TypeTree::TreePath<1> > PSubGFS;
      PSubGFS pSubGfs;
      typedef Dune::PDELab::DiscreteGridFunction<PSubGFS,NS_V> PDGF;
      PDGF pDgf;
    
      typedef Dune::PDELab::GridFunctionSubSpace<typename GFS::CH,Dune::TypeTree::TreePath<0> > PhiSubGFS;
      PhiSubGFS phiSubGfs;
      typedef Dune::PDELab::DiscreteGridFunction<PhiSubGFS,CH_V> PhiDGF;
      PhiDGF phiDgf;
    
      typedef Dune::PDELab::GridFunctionSubSpace<typename GFS::CH,Dune::TypeTree::TreePath<1> > MuSubGFS;
      MuSubGFS muSubGfs;
      typedef Dune::PDELab::DiscreteGridFunction<MuSubGFS,CH_V> MuDGF;
      MuDGF muDgf;
    
      typedef Dune::PDELab::GridFunctionSubSpace<typename GFS::CH,Dune::TypeTree::TreePath<2> > Phi2SubGFS;
      Phi2SubGFS phi2SubGfs;
      typedef Dune::PDELab::DiscreteGridFunction<Phi2SubGFS,CH_V> Phi2DGF;
      Phi2DGF phi2Dgf;
    
      typedef Dune::PDELab::GridFunctionSubSpace<typename GFS::CH,Dune::TypeTree::TreePath<3> > Mu2SubGFS;
      Mu2SubGFS mu2SubGfs;
      typedef Dune::PDELab::DiscreteGridFunction<Mu2SubGFS,CH_V> Mu2DGF;
      Mu2DGF mu2Dgf;
    
      typedef PF_3P_Color_DGF<typename GFS::CH,CH_V,3> ColorDGF;
      ColorDGF colorDgf;
    
    
      DGF_fff_chns(const GFS& gfs, const NS_V& nsV, const CH_V& chV) :
          vSubGfs(gfs.ns),
          vDgf(vSubGfs,nsV),
          pSubGfs(gfs.ns),
          pDgf(pSubGfs,nsV),
    
          phiSubGfs(gfs.ch),
          phiDgf(phiSubGfs,chV),
          muSubGfs(gfs.ch),
          muDgf(muSubGfs,chV),
          phi2SubGfs(gfs.ch),
          phi2Dgf(phi2SubGfs,chV),
          mu2SubGfs(gfs.ch),
          mu2Dgf(mu2SubGfs,chV),
    
          colorDgf(gfs.ch,chV)
      {
      }
    
      template<typename VTKWRITER>
      void addToVTKWriter(VTKWRITER& vtkwriter)
      {
        vtkwriter.addVertexData(std::make_shared<Dune::PDELab::VTKGridFunctionAdapter<VDGF> >(vDgf,"v"));
        vtkwriter.addVertexData(std::make_shared<Dune::PDELab::VTKGridFunctionAdapter<PDGF> >(pDgf,"p"));
        vtkwriter.addVertexData(std::make_shared<Dune::PDELab::VTKGridFunctionAdapter<PhiDGF> >(phiDgf,"phi"));
        vtkwriter.addVertexData(std::make_shared<Dune::PDELab::VTKGridFunctionAdapter<MuDGF> >(muDgf,"mu"));
        vtkwriter.addVertexData(std::make_shared<Dune::PDELab::VTKGridFunctionAdapter<Phi2DGF> >(phi2Dgf,"phi2"));
        vtkwriter.addVertexData(std::make_shared<Dune::PDELab::VTKGridFunctionAdapter<Mu2DGF> >(mu2Dgf,"mu2"));
        vtkwriter.addVertexData(std::make_shared<Dune::PDELab::VTKGridFunctionAdapter<ColorDGF> >(colorDgf,"_coloring_"));
      }
    };
    
    //_________________________________________________________________________________________________________________________
    
    template<typename IniV, typename IniP, typename IniPhi, typename IniMu, typename IniPhi2, typename IniMu2>
    class Ini_fff_chns
    {
    public:
      IniV iniV;
      IniP iniP;
      IniPhi iniPhi;
      IniMu iniMu;
      IniPhi2 iniPhi2;
      IniMu2 iniMu2;
      typedef Dune::PDELab::CompositeGridFunction<IniV,IniP> IniNS;
      typedef Dune::PDELab::CompositeGridFunction<IniPhi,IniMu,IniPhi2,IniMu2> IniCH;
      IniNS iniNS;
      IniCH iniCH;
    
      template<typename GV, typename Params>
      Ini_fff_chns(const GV& gv, const Params& params) :
          iniV(gv, params),
          iniP(gv, params),
          iniPhi(gv, params),
          iniMu(gv, params),
          iniPhi2(gv, params),
          iniMu2(gv, params),
          iniNS(iniV,iniP),
          iniCH(iniPhi,iniMu,iniPhi2,iniMu2)
      {
      }
    };
    
    //_________________________________________________________________________________________________________________________
    
    template<typename BdryV1, typename BdryV2, typename BdryP,
              typename BdryPhi, typename BdryMu, typename BdryPhi2, typename BdryMu2>
    class Bdry_fff_chns_compositeV_2d
    {
    public:
      BdryV1 bdryV1;
      BdryV2 bdryV2;
      BdryP bdryP;
      BdryPhi bdryPhi;
      BdryMu bdryMu;
      BdryPhi2 bdryPhi2;
      BdryMu2 bdryMu2;
      typedef Dune::PDELab::CompositeConstraintsParameters<BdryV1,BdryV2>    BdryV;
      typedef Dune::PDELab::CompositeConstraintsParameters<BdryV, BdryP> BdryNS;
      typedef Dune::PDELab::CompositeConstraintsParameters<BdryPhi, BdryMu, BdryPhi2, BdryMu2> BdryCH;
      BdryV bdryV;
      BdryNS bdryNS;
      BdryCH bdryCH;
    
      Bdry_fff_chns_compositeV_2d() :
          bdryV1(),
          bdryV2(),
          bdryP(),
          bdryPhi(),
          bdryMu(),
          bdryPhi2(),
          bdryMu2(),
          bdryV(bdryV1,bdryV2),
          bdryNS(bdryV,bdryP),
          bdryCH(bdryPhi,bdryMu,bdryPhi2,bdryMu2)
      {
      }
    };