root/src/secOrd_op_initData_u.cpp

/* [<][>][^][v][top][bottom][index][help] */

DEFINITIONS

This source file includes following definitions.
  1. interpolate
  2. evaluate

   1 /*! \file secOrd_op_initData_u.cpp
   2 
   3     \brief <Program Name>
   4 
   5      Revision history:
   6 
   7           Revised by Christian Power dd.mm.yyyy
   8           Originally written by Christian Power
   9                (power22c@gmail.com) February 2016
  10 
  11      Implementation details for secOrd_op_initData_u.h
  12      Created by Christian Power on 04.02.2016
  13      Copyright (c) 2016 Christian Power. All rights reserved.
  14  */
  15 
  16 #include <stdexcept>
  17 #include <cmath>
  18 #include <config.h>
  19 #include <dune/fem/operator/lagrangeinterpolation.hh>
  20 #include "secOrd_op_initData_u.h"
  21 #include "grid.h"
  22 
  23 #ifdef DEBUG
  24 #include <iostream>
  25 #endif
  26 
  27 class Initial_data_u
  28   : public Dune::Fem::Function<Esfem::Grid::Grid_and_time::Function_space,
  29                                Initial_data_u>
  30 {
  31 public:
  32   using Base = Esfem::Grid::Grid_and_time::Function_space;
  33   using Domain = Base::DomainType;
  34   using Range = Base::RangeType;
  35   
  36   explicit Initial_data_u(const Dune::Fem::TimeProviderBase&);
  37   Initial_data_u(const Initial_data_u&) = delete;
  38   Initial_data_u& operator=(const Initial_data_u&) = delete;
  39 
  40   void evaluate(const Domain&, Range&) const;
  41 private:
  42   const Dune::Fem::TimeProviderBase& tp;
  43 };
  44 
  45 // ----------------------------------------------------------------------
  46 // Implementation esfem.h
  47 
  48 struct Esfem::SecOrd_op::Init_data_u::Data{
  49   Initial_data_u u0;
  50   Data(const Grid::Grid_and_time& gt)
  51     : u0 {gt.time_provider()}
  52   {}
  53 };
  54 
  55 Esfem::SecOrd_op::Init_data_u::Init_data_u(const Grid::Grid_and_time& gt){
  56   try{
  57     d_ptr = new Data {gt};
  58   }
  59   catch(const std::exception&){
  60     std::throw_with_nested(std::logic_error{"Error in constructor of Init_data_u."});
  61   }
  62   catch(...){
  63     throw std::logic_error{"Unknown error in constructor of Init_data_u."};
  64   }
  65 }
  66 Esfem::SecOrd_op::Init_data_u::~Init_data_u(){
  67   delete d_ptr;
  68   d_ptr = nullptr;
  69 #ifdef DEBUG
  70   std::cerr << "~Init_data_u(): delete d_ptr.\n";
  71 #endif
  72 }
  73 void Esfem::SecOrd_op::Init_data_u::interpolate(Grid::Scal_FEfun& fef) const{
  74   using FE_function = Esfem::Grid::Scal_FEfun::Dune_FEfun;
  75   Dune::LagrangeInterpolation<FE_function>::interpolateFunction(d_ptr -> u0, fef);  
  76 }
  77 
  78 // ----------------------------------------------------------------------
  79 // Implementation Initial_data_u
  80 
  81 Initial_data_u::Initial_data_u(const Dune::Fem::TimeProviderBase& tpb)
  82   : tp {tpb}
  83 {}
  84 void Initial_data_u::evaluate(const Domain& d, Range& r) const{
  85   static_assert(Domain::dimension == 3, "Bad domain dimension.");
  86   static_assert(Range::dimension == 1, "Bad range dimension.");
  87   const double x = d[0];
  88   const double y = d[1];
  89   const double t = tp.time();
  90   r = std::exp(-6.*t)*x*y;
  91 }
  92 
  93 /*! Log:
  94  */

/* [<][>][^][v][top][bottom][index][help] */