DOpE
statpdeproblem.h
Go to the documentation of this file.
1 
24 #ifndef _STAT_PDE_PROBLEM_H_
25 #define _STAT_PDE_PROBLEM_H_
26 
27 #include "pdeprobleminterface.h"
28 #include "integrator.h"
29 #include "parameterreader.h"
30 #include "statevector.h"
31 #include "stateproblem.h"
32 
33 #include <lac/vector.h>
34 
35 #include <lac/block_sparsity_pattern.h>
36 #include <lac/block_sparse_matrix.h>
37 
38 #include "pdeproblemcontainer.h"
39 #include "pdeinterface.h"
40 #include "functionalinterface.h"
41 #include "dirichletdatainterface.h"
42 #include "dopeexception.h"
43 #include "newtonsolver.h"
44 #include "cglinearsolver.h"
45 #include "gmreslinearsolver.h"
46 #include "directlinearsolver.h"
47 #include "solutionextractor.h"
48 
49 #include <base/data_out_base.h>
50 #include <base/utilities.h>
51 #include <numerics/data_out.h>
52 #include <numerics/matrix_tools.h>
53 #include <numerics/vector_tools.h>
54 #include <base/function.h>
55 #include <lac/sparse_matrix.h>
56 #include <lac/compressed_simple_sparsity_pattern.h>
57 #include <lac/block_sparsity_pattern.h>
58 #include <lac/sparse_direct.h>
59 
60 #include <fstream>
61 namespace DOpE
62 {
77  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
78  typename VECTOR, int dealdim>
79  class StatPDEProblem : public PDEProblemInterface<PROBLEM, VECTOR, dealdim>
80  {
81  public:
94  template<typename INTEGRATORDATACONT>
95  StatPDEProblem(PROBLEM *OP, std::string state_behavior,
96  ParameterReader &param_reader, INTEGRATORDATACONT& idc,
97  int base_priority = 0);
98 
111  template<typename INTEGRATORDATACONT>
112  StatPDEProblem(PROBLEM *OP, std::string state_behavior,
113  ParameterReader &param_reader, INTEGRATORDATACONT& idc1,
114  INTEGRATORDATACONT& idc2, int base_priority = 0);
115 
116  virtual
117  ~StatPDEProblem();
118 
119  /******************************************************/
120 
126  static void
127  declare_params(ParameterReader &param_reader);
128 
129  /******************************************************/
130 
136  void
137  ReInit();
138 
139  /******************************************************/
140 
146  void
148 
149  /******************************************************/
150 
156  void
158 
159  /******************************************************/
160 
178  template<class DWRC, class PDE>
179  void
180  ComputeRefinementIndicators(DWRC& dwrc, PDE& pde);
181 
182  /******************************************************/
183 
189  void
190  StateSizeInfo(std::stringstream& out)
191  {
192  GetU().PrintInfos(out);
193  }
194 
195  /******************************************************/
196 
207  void
208  WriteToFile(const VECTOR &v, std::string name, std::string outfile,
209  std::string dof_type, std::string filetype);
210 
211  /******************************************************/
212 
226  void
227  WriteToFileElementwise(const Vector<double> &v, std::string name,
228  std::string outfile, std::string dof_type, std::string filetype);
229 
230  /******************************************************/
231 
242  void
243  WriteToFile(const ControlVector<VECTOR> &v, std::string name,
244  std::string outfile, std::string dof_type, std::string filetype);
245 
253  void
254  WriteToFile(const std::vector<double> &/*v*/,
255  std::string /*outfile*/)
256  {
257  abort();
258  }
259 
260  protected:
264  void
266 
267  /******************************************************/
268 
278  void
280 
281  /******************************************************/
285  const StateVector<VECTOR> &
286  GetU() const
287  {
288  return _u;
289  }
292  {
293  return _u;
294  }
298  const StateVector<VECTOR> &
299  GetZForEE() const
300  {
301  return _z_for_ee;
302  }
305  {
306  return _z_for_ee;
307  }
308 
309  NONLINEARSOLVER&
310  GetNonlinearSolver(std::string type);
311  INTEGRATOR&
313  {
314  return _integrator;
315  }
316 
317  const bool&
319  {
320  return _build_state_matrix;
321  }
322 
323  const bool&
325  {
326  return _build_adjoint_matrix;
327  }
328 
329  const bool&
331  {
332  return _state_reinit;
333  }
334 
335  const bool&
337  {
338  return _adjoint_reinit;
339  }
340 
341 
342  bool&
344  {
345  return _build_state_matrix;
346  }
347 
348  bool&
350  {
351  return _build_adjoint_matrix;
352  }
353 
354  bool&
356  {
357  return _state_reinit;
358  }
359 
360  bool&
362  {
363  return _adjoint_reinit;
364  }
365 
366 
367  private:
372  void
373  AddUDD()
374  {
375  for (auto it = this->GetUserDomainData().begin();
376  it != this->GetUserDomainData().end(); it++)
377  {
378  this->GetIntegrator().AddDomainData(it->first, it->second);
379  }
380  }
381 
386  void
387  DeleteUDD()
388  {
389  for (auto it = this->GetUserDomainData().begin();
390  it != this->GetUserDomainData().end(); it++)
391  {
392  this->GetIntegrator().DeleteDomainData(it->first);
393  }
394  }
395 
396  StateVector<VECTOR> _u;
397  StateVector<VECTOR> _z_for_ee;
398 
399  INTEGRATOR _integrator;
400  NONLINEARSOLVER _nonlinear_state_solver;
401  NONLINEARSOLVER _nonlinear_adjoint_solver;
402 
403  bool _build_state_matrix;
404  bool _build_adjoint_matrix;
405  bool _state_reinit;
406  bool _adjoint_reinit;
407 
408  int _n_patches;
409 
410  friend class SolutionExtractor<
411  StatPDEProblem<NONLINEARSOLVER, INTEGRATOR, PROBLEM, VECTOR, dealdim>,
412  VECTOR> ;
413  };
414 
415  /*************************************************************************/
416  /*****************************IMPLEMENTATION******************************/
417  /*************************************************************************/
418  using namespace dealii;
419 
420  /******************************************************/
421  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
422  typename VECTOR, int dealdim>
423  void
425  ParameterReader &param_reader)
426  {
427  NONLINEARSOLVER::declare_params(param_reader);
428  param_reader.SetSubsection("output parameters");
429  param_reader.declare_entry("number of patches", "0",
430  Patterns::Integer(0));
431 
432  }
433  /******************************************************/
434 
435  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
436  typename VECTOR, int dealdim>
437  template<typename INTEGRATORDATACONT>
439  PROBLEM *OP, std::string state_behavior,
440  ParameterReader &param_reader, INTEGRATORDATACONT& idc,
441  int base_priority)
442  : PDEProblemInterface<PROBLEM, VECTOR, dealdim>(OP, base_priority), _u(
443  OP->GetSpaceTimeHandler(), state_behavior, param_reader), _z_for_ee(
444  OP->GetSpaceTimeHandler(), state_behavior, param_reader), _integrator(
445  idc), _nonlinear_state_solver(_integrator, param_reader), _nonlinear_adjoint_solver(
446  _integrator, param_reader)
447  {
448  param_reader.SetSubsection("output parameters");
449  _n_patches = param_reader.get_integer("number of patches");
450  //PDEProblems should be ReInited
451  {
452  _state_reinit = true;
453  _adjoint_reinit = true;
454  }
455  }
456 
457  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
458  typename VECTOR, int dealdim>
459  template<typename INTEGRATORDATACONT>
461  PROBLEM *OP, std::string state_behavior,
462  ParameterReader &param_reader, INTEGRATORDATACONT& idc,
463  INTEGRATORDATACONT& idc2, int base_priority)
464  : PDEProblemInterface<PROBLEM, VECTOR, dealdim>(OP, base_priority), _u(
465  OP->GetSpaceTimeHandler(), state_behavior, param_reader), _z_for_ee(
466  OP->GetSpaceTimeHandler(), state_behavior, param_reader), _integrator(
467  idc, idc2), _nonlinear_state_solver(_integrator, param_reader), _nonlinear_adjoint_solver(
468  _integrator, param_reader)
469  {
470  //PDEProblems should be ReInited
471  {
472  _state_reinit = true;
473  _adjoint_reinit = true;
474  }
475  }
476 
477  /******************************************************/
478 
479  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
480  typename VECTOR, int dealdim>
482  {
483  }
484 
485  /******************************************************/
486 
487  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
488  typename VECTOR, int dealdim>
489  NONLINEARSOLVER&
491  std::string type)
492  {
493  if (type == "state")
494  {
495  return _nonlinear_state_solver;
496  }
497  else if (type == "adjoint_for_ee")
498  {
499  return _nonlinear_adjoint_solver;
500  }
501  else
502  {
503  throw DOpEException("No Solver for Problem type:`" + type + "' found",
504  "StatReducedProblem::GetNonlinearSolver");
505 
506  }
507  }
508 
509  /******************************************************/
510 
511  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
512  typename VECTOR, int dealdim>
513  void
515  {
517 
518  //Some Solvers must be reinited when called
519  // Better have subproblems, so that solver can be reinited here
520  {
521  _state_reinit = true;
522  _adjoint_reinit = true;
523  }
524 
525  _build_state_matrix = true;
526  _build_adjoint_matrix = true;
527 
528  GetU().ReInit();
529  GetZForEE().ReInit();
530  }
531 
532  /******************************************************/
533 
534  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
535  typename VECTOR, int dealdim>
536  void
538  {
539  this->InitializeFunctionalValues(this->GetProblem()->GetNFunctionals());
540  this->GetOutputHandler()->Write("Computing State Solution:",
541  4 + this->GetBasePriority());
542 
543  this->SetProblemType("state");
544  auto& problem = this->GetProblem()->GetStateProblem();
545  if (_state_reinit == true)
546  {
547  GetNonlinearSolver("state").ReInit(problem);
548  _state_reinit = false;
549  }
550 
551  this->GetProblem()->AddAuxiliaryToIntegrator(this->GetIntegrator());
552  AddUDD();
553  _build_state_matrix = this->GetNonlinearSolver("state").NonlinearSolve(
554  problem, (GetU().GetSpacialVector()), true, _build_state_matrix);
555  DeleteUDD();
556 
557  this->GetProblem()->DeleteAuxiliaryFromIntegrator(this->GetIntegrator());
558 
559  this->GetOutputHandler()->Write((GetU().GetSpacialVector()),
560  "State" + this->GetPostIndex(), problem.GetDoFType());
561 
562  }
563 
564  /******************************************************/
565 
566  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
567  typename VECTOR, int dealdim>
568  void
570  DOpEtypes::WeightComputation weight_comp)
571  {
572  this->GetOutputHandler()->Write("Computing Dual for Error Estimation:",
573  4 + this->GetBasePriority());
574  if (weight_comp == DOpEtypes::higher_order_interpolation)
575  {
576  this->SetProblemType("adjoint_for_ee");
577  }
578  else
579  {
580  throw DOpEException("Unknown WeightComputation",
581  "StatPDEProblem::ComputeDualForErrorEstimation");
582  }
583 
584  // auto& problem = this->GetProblem()->GetStateProblem();//Hier ist adjoint problem einzufuegen
585  auto& problem = *(this->GetProblem());
586  if (_adjoint_reinit == true)
587  {
588  GetNonlinearSolver("adjoint_for_ee").ReInit(problem);
589  _adjoint_reinit = false;
590  }
591 
592  this->GetProblem()->AddAuxiliaryToIntegrator(this->GetIntegrator());
593 
594  this->GetIntegrator().AddDomainData("state",
595  &(GetU().GetSpacialVector()));
596  AddUDD();
597 
598  _build_adjoint_matrix =
599  this->GetNonlinearSolver("adjoint_for_ee").NonlinearSolve(problem,
600  (GetZForEE().GetSpacialVector()), true, _build_adjoint_matrix);
601 
602  this->GetIntegrator().DeleteDomainData("state");
603  DeleteUDD();
604 
605  this->GetProblem()->DeleteAuxiliaryFromIntegrator(this->GetIntegrator());
606 
607  this->GetOutputHandler()->Write((GetZForEE().GetSpacialVector()),
608  "Adjoint_for_ee" + this->GetPostIndex(), problem.GetDoFType());
609 
610  }
611 
612  /******************************************************/
613 
614  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
615  typename VECTOR, int dealdim>
616  void
618  {
619  this->ComputeReducedState();
620 
621  this->GetOutputHandler()->Write("Computing Functionals:",
622  4 + this->GetBasePriority());
623 
624  this->GetProblem()->AddAuxiliaryToIntegrator(this->GetIntegrator());
625 
626  this->GetIntegrator().AddDomainData("state",
627  &(GetU().GetSpacialVector()));
628  AddUDD();
629 
630  for (unsigned int i = 0; i < this->GetProblem()->GetNFunctionals(); i++)
631  {
632  double ret = 0;
633  bool found = false;
634 
635  this->SetProblemType("aux_functional", i);
636  if (this->GetProblem()->GetFunctionalType().find("domain")
637  != std::string::npos)
638  {
639  found = true;
640  ret += this->GetIntegrator().ComputeDomainScalar(
641  *(this->GetProblem()));
642  }
643  if (this->GetProblem()->GetFunctionalType().find("point")
644  != std::string::npos)
645  {
646  found = true;
647  ret += this->GetIntegrator().ComputePointScalar(
648  *(this->GetProblem()));
649  }
650  if (this->GetProblem()->GetFunctionalType().find("boundary")
651  != std::string::npos)
652  {
653  found = true;
654  ret += this->GetIntegrator().ComputeBoundaryScalar(
655  *(this->GetProblem()));
656  }
657  if (this->GetProblem()->GetFunctionalType().find("face")
658  != std::string::npos)
659  {
660  found = true;
661  ret += this->GetIntegrator().ComputeFaceScalar(*(this->GetProblem()));
662  }
663 
664  if (!found)
665  {
666  throw DOpEException(
667  "Unknown Functional Type: "
668  + this->GetProblem()->GetFunctionalType(),
669  "StatPDEProblem::ComputeReducedFunctionals");
670  }
671  this->GetFunctionalValues()[i].push_back(ret);
672  std::stringstream out;
673  this->GetOutputHandler()->InitOut(out);
674  out << this->GetProblem()->GetFunctionalName() << ": " << ret;
675  this->GetOutputHandler()->Write(out, 2 + this->GetBasePriority());
676  }
677 
678  this->GetIntegrator().DeleteDomainData("state");
679  DeleteUDD();
680 
681  this->GetProblem()->DeleteAuxiliaryFromIntegrator(this->GetIntegrator());
682  }
683 
684  /******************************************************/
685 
686  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
687  typename VECTOR, int dealdim>
688  void
690  const StateVector<VECTOR>& st_u)
691  {
692  this->InitializeFunctionalValues(this->GetProblem()->GetNFunctionals());
693  GetU() = st_u;
694 
695  this->GetOutputHandler()->Write("Computing Functionals:",
696  4 + this->GetBasePriority());
697 
698  this->GetProblem()->AddAuxiliaryToIntegrator(this->GetIntegrator());
699 
700  this->GetIntegrator().AddDomainData("state",
701  &(GetU().GetSpacialVector()));
702  AddUDD();
703 
704  for (unsigned int i = 0; i < this->GetProblem()->GetNFunctionals(); i++)
705  {
706  double ret = 0;
707  bool found = false;
708 
709  this->SetProblemType("aux_functional", i);
710  if (this->GetProblem()->GetFunctionalType().find("domain")
711  != std::string::npos)
712  {
713  found = true;
714  ret += this->GetIntegrator().ComputeDomainScalar(
715  *(this->GetProblem()));
716  }
717  if (this->GetProblem()->GetFunctionalType().find("point")
718  != std::string::npos)
719  {
720  found = true;
721  ret += this->GetIntegrator().ComputePointScalar(
722  *(this->GetProblem()));
723  }
724  if (this->GetProblem()->GetFunctionalType().find("boundary")
725  != std::string::npos)
726  {
727  found = true;
728  ret += this->GetIntegrator().ComputeBoundaryScalar(
729  *(this->GetProblem()));
730  }
731  if (this->GetProblem()->GetFunctionalType().find("face")
732  != std::string::npos)
733  {
734  found = true;
735  ret += this->GetIntegrator().ComputeFaceScalar(*(this->GetProblem()));
736  }
737 
738  if (!found)
739  {
740  throw DOpEException(
741  "Unknown Functional Type: "
742  + this->GetProblem()->GetFunctionalType(),
743  "StatPDEProblem::ComputeReducedFunctionals");
744  }
745  this->GetFunctionalValues()[i].push_back(ret);
746  std::stringstream out;
747  this->GetOutputHandler()->InitOut(out);
748  out << this->GetProblem()->GetFunctionalName() << ": " << ret;
749  this->GetOutputHandler()->Write(out, 2 + this->GetBasePriority());
750  }
751 
752  this->GetIntegrator().DeleteDomainData("state");
753  DeleteUDD();
754 
755  this->GetProblem()->DeleteAuxiliaryFromIntegrator(this->GetIntegrator());
756  }
757 
758  /******************************************************/
759 
760  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
761  typename VECTOR, int dealdim>
762  template<class DWRC, class PDE>
763  void
765  DWRC& dwrc, PDE& pde)
766  {
767  //Attach the ResidualModifier to the PDE.
768  pde.ResidualModifier = boost::bind<void>(
769  boost::mem_fn(&DWRC::ResidualModifier), boost::ref(dwrc), _1);
770  pde.VectorResidualModifier = boost::bind<void>(
771  boost::mem_fn(&DWRC::VectorResidualModifier), boost::ref(dwrc), _1);
772  //first we reinit the dwrdatacontainer (this
773  //sets the weight-vectors to their correct length)
774  const unsigned int n_elements =
775  this->GetProblem()->GetSpaceTimeHandler()->GetStateDoFHandler().get_tria().n_active_cells();
776  dwrc.ReInit(n_elements);
777  //If we need the dual solution, compute it
778  if (dwrc.NeedDual())
779  this->ComputeDualForErrorEstimation(dwrc.GetWeightComputation());
780 
781  //some output
782  this->GetOutputHandler()->Write("Computing Error Indicators:",
783  4 + this->GetBasePriority());
784 
785  this->GetProblem()->AddAuxiliaryToIntegrator(this->GetIntegrator());
786 
787  //add the primal and (if needed) dual solution to the integrator
788  this->GetIntegrator().AddDomainData("state",
789  &(GetU().GetSpacialVector()));
790  if (dwrc.NeedDual())
791  this->GetIntegrator().AddDomainData("adjoint_for_ee",
792  &(GetZForEE().GetSpacialVector()));
793  AddUDD();
794 
795  this->SetProblemType("error_evaluation");
796 
797  //prepare the weights...
798  dwrc.PrepareWeights(GetU(), GetZForEE());
799 
800  //now we finally compute the refinement indicators
801  this->GetIntegrator().ComputeRefinementIndicators(*this->GetProblem(),
802  dwrc);
803 
804  // release the lock on the refinement indicators (see dwrcontainer.h)
805  dwrc.ReleaseLock();
806 
807  const float error = dwrc.GetError();
808 
809  // clear the data
810  dwrc.ClearWeightData();
811  this->GetIntegrator().DeleteDomainData("state");
812  if (dwrc.NeedDual())
813  this->GetIntegrator().DeleteDomainData("adjoint_for_ee");
814  DeleteUDD();
815  this->GetProblem()->DeleteAuxiliaryFromIntegrator(
816  this->GetIntegrator());
817 
818  std::stringstream out;
819  this->GetOutputHandler()->InitOut(out);
820  out << "Error estimate using " << dwrc.GetName();
821  if (dwrc.NeedDual())
822  out << " for the " << this->GetProblem()->GetFunctionalName();
823  out << ": " << error;
824  this->GetOutputHandler()->Write(out, 2 + this->GetBasePriority());
825  this->GetOutputHandler()->WriteElementwise(dwrc.GetErrorIndicators(),
826  "Error_Indicators" + this->GetPostIndex(),
827  this->GetProblem()->GetDoFType());
828 
829  }
830 
831  /******************************************************/
832 
833  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
834  typename VECTOR, int dealdim>
835  void
837  const Vector<double> &v, std::string name, std::string outfile,
838  std::string dof_type, std::string filetype)
839  {
840  if (dof_type == "state")
841  {
842  auto& data_out =
843  this->GetProblem()->GetSpaceTimeHandler()->GetDataOut();
844  data_out.attach_dof_handler(
845  this->GetProblem()->GetSpaceTimeHandler()->GetStateDoFHandler());
846 
847  data_out.add_data_vector(v, name);
848  data_out.build_patches(_n_patches);
849 
850  std::ofstream output(outfile.c_str());
851 
852  if (filetype == ".vtk")
853  {
854  data_out.write_vtk(output);
855  }
856  else
857  {
858  throw DOpEException(
859  "Don't know how to write filetype `" + filetype + "'!",
860  "StatPDEProblem::WriteToFileElementwise");
861  }
862  data_out.clear();
863  }
864  else
865  {
866  throw DOpEException("No such DoFHandler `" + dof_type + "'!",
867  "StatPDEProblem::WriteToFileElementwise");
868  }
869  }
870 
871  /******************************************************/
872 
873  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
874  typename VECTOR, int dealdim>
875  void
877  const VECTOR &v, std::string name, std::string outfile,
878  std::string dof_type, std::string filetype)
879  {
880  if (dof_type == "state")
881  {
882  auto& data_out =
883  this->GetProblem()->GetSpaceTimeHandler()->GetDataOut();
884  data_out.attach_dof_handler(
885  this->GetProblem()->GetSpaceTimeHandler()->GetStateDoFHandler());
886 
887  data_out.add_data_vector(v, name);
888  //TODO: mapping[0] is a workaround, as deal does not support interpolate
889  // boundary_values with a mapping collection at this point.
890  data_out.build_patches(
891  this->GetProblem()->GetSpaceTimeHandler()->GetMapping()[0]);
892 
893  std::ofstream output(outfile.c_str());
894 
895  if (filetype == ".vtk")
896  {
897  data_out.write_vtk(output);
898  }
899  else if (filetype == ".gpl")
900  {
901  data_out.write_gnuplot(output);
902  }
903  else
904  {
905  throw DOpEException(
906  "Don't know how to write filetype `" + filetype + "'!",
907  "StatPDEProblem::WriteToFile");
908  }
909  data_out.clear();
910  }
911  else
912  {
913  throw DOpEException("No such DoFHandler `" + dof_type + "'!",
914  "StatPDEProblem::WriteToFile");
915  }
916  }
917 
918  /******************************************************/
919 
920  template<typename NONLINEARSOLVER, typename INTEGRATOR, typename PROBLEM,
921  typename VECTOR, int dealdim>
922  void
924  const ControlVector<VECTOR> &v, std::string name, std::string outfile,
925  std::string dof_type, std::string filetype)
926  {
927  WriteToFile(v.GetSpacialVector(), name, outfile, dof_type, filetype);
928  }
929 }
930 #endif
StateVector< VECTOR > & GetU()
Definition: statpdeproblem.h:291
void StateSizeInfo(std::stringstream &out)
Definition: statpdeproblem.h:190
StatPDEProblem(PROBLEM *OP, std::string state_behavior, ParameterReader &param_reader, INTEGRATORDATACONT &idc, int base_priority=0)
Definition: statpdeproblem.h:438
void WriteToFileElementwise(const Vector< double > &v, std::string name, std::string outfile, std::string dof_type, std::string filetype)
Definition: statpdeproblem.h:836
Definition: statpdeproblem.h:79
void ComputeReducedFunctionals()
Definition: statpdeproblem.h:617
void WriteToFile(const std::vector< double > &, std::string)
Definition: statpdeproblem.h:254
Definition: parameterreader.h:36
const bool & GetBuildAdjointMatrix() const
Definition: statpdeproblem.h:324
INTEGRATOR & GetIntegrator()
Definition: statpdeproblem.h:312
const StateVector< VECTOR > & GetZForEE() const
Definition: statpdeproblem.h:299
const bool & GetStateReinit() const
Definition: statpdeproblem.h:330
void declare_entry(const std::string &entry, const std::string &default_value, const Patterns::PatternBase &pattern=Patterns::Anything(), const std::string &documentation=std::string())
Definition: parameterreader.h:98
const bool & GetBuildStateMatrix() const
Definition: statpdeproblem.h:318
NONLINEARSOLVER & GetNonlinearSolver(std::string type)
Definition: statpdeproblem.h:490
Definition: controlvector.h:48
Definition: pdeprobleminterface.h:53
StateVector< VECTOR > & GetZForEE()
Definition: statpdeproblem.h:304
bool & GetBuildAdjointMatrix()
Definition: statpdeproblem.h:349
static void declare_params(ParameterReader &param_reader)
Definition: statpdeproblem.h:424
Definition: statevector.h:49
VECTOR & GetSpacialVector()
Definition: controlvector.cc:118
bool & GetStateReinit()
Definition: statpdeproblem.h:355
const StateVector< VECTOR > & GetU() const
Definition: statpdeproblem.h:286
bool & GetAdjointReinit()
Definition: statpdeproblem.h:361
void WriteToFile(const VECTOR &v, std::string name, std::string outfile, std::string dof_type, std::string filetype)
Definition: statpdeproblem.h:876
virtual ~StatPDEProblem()
Definition: statpdeproblem.h:481
void ComputeRefinementIndicators(DWRC &dwrc, PDE &pde)
Definition: statpdeproblem.h:764
int get_integer(const std::string &entry_name)
Definition: parameterreader.h:126
void ComputeDualForErrorEstimation(DOpEtypes::WeightComputation)
Definition: statpdeproblem.h:569
void SetSubsection(const std::string subsection)
Definition: parameterreader.h:93
void ReInit()
Definition: statpdeproblem.h:514
WeightComputation
Definition: dopetypes.h:78
const std::map< std::string, const VECTOR * > & GetUserDomainData() const
Definition: reducedprobleminterface.h:303
bool & GetBuildStateMatrix()
Definition: statpdeproblem.h:343
const bool & GetAdjointReinit() const
Definition: statpdeproblem.h:336
void ComputeReducedState()
Definition: statpdeproblem.h:537
Definition: dopeexception.h:35