ug4
grid_function_util.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2010-2015: G-CSC, Goethe University Frankfurt
3  * Author: Andreas Vogel
4  *
5  * This file is part of UG4.
6  *
7  * UG4 is free software: you can redistribute it and/or modify it under the
8  * terms of the GNU Lesser General Public License version 3 (as published by the
9  * Free Software Foundation) with the following additional attribution
10  * requirements (according to LGPL/GPL v3 §7):
11  *
12  * (1) The following notice must be displayed in the Appropriate Legal Notices
13  * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14  *
15  * (2) The following notice must be displayed at a prominent place in the
16  * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17  *
18  * (3) The following bibliography is recommended for citation and must be
19  * preserved in all covered files:
20  * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21  * parallel geometric multigrid solver on hierarchically distributed grids.
22  * Computing and visualization in science 16, 4 (2013), 151-164"
23  * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24  * flexible software system for simulating pde based models on high performance
25  * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26  *
27  * This program is distributed in the hope that it will be useful,
28  * but WITHOUT ANY WARRANTY; without even the implied warranty of
29  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30  * GNU Lesser General Public License for more details.
31  */
32 
33 #ifndef __H__UG__LIB_DISC__FUNCTION_SPACE__GRID_FUNCTION_UTIL__
34 #define __H__UG__LIB_DISC__FUNCTION_SPACE__GRID_FUNCTION_UTIL__
35 
36 #include <vector>
37 #include <string>
38 #include <cmath> // for isinf, isnan
39 #include <boost/function.hpp>
40 
41 
42 #include "common/util/file_util.h"
43 
51 #include "lib_disc/io/vtkoutput.h"
58 #include "lib_grid/algorithms/debug_util.h" // for ElementDebugInfo
60 
61 #include "grid_function.h"
62 #include "dof_position_util.h"
63 
64 #ifdef UG_PARALLEL
65 #include "pcl/pcl.h"
66 #endif
67 
68 namespace ug {
69 
70 #ifndef isnan
71 using boost::math::isnan;
72 #endif
73 
74 #ifndef isinf
75 using boost::math::isinf;
76 #endif
77 
78 
79 template <typename TBaseElem, typename TGridFunction>
80 static void ScaleGFOnElems
81 (
85  const std::vector<number>& vScale
86 )
87 {
89  std::vector<DoFIndex> vInd;
90 
91  try
92  {
93  // iterate all elements (including SHADOW_RIM_COPY!)
94  iter = dd->template begin<TBaseElem>(SurfaceView::ALL);
95  iterEnd = dd->template end<TBaseElem>(SurfaceView::ALL);
96  for (; iter != iterEnd; ++iter)
97  {
98  for (size_t fi = 0; fi < dd->num_fct(); ++fi)
99  {
100  size_t nInd = dd->inner_dof_indices(*iter, fi, vInd);
101 
102  // remember multi indices
103  for (size_t dof = 0; dof < nInd; ++dof)
104  DoFRef(*vecOut, vInd[dof]) = vScale[fi] * DoFRef(*vecIn, vInd[dof]);
105  }
106  }
107  }
108  UG_CATCH_THROW("Error while scaling vector.")
109 }
110 
111 
121 template <typename TGridFunction>
123 (
124  SmartPtr<TGridFunction> scaledVecOut,
126  const std::vector<number>& scalingFactors
127 )
128 {
129  // check that the correct numbers of scaling factors are given
130  size_t n = scalingFactors.size();
131  UG_COND_THROW(n != vecIn->num_fct(), "Number of scaling factors (" << n << ") "
132  "does not match number of functions given in dimless vector (" << vecIn->num_fct() << ").");
133 
134  // check that input and output vectors have the same number of components and dofs
135  UG_COND_THROW(n != scaledVecOut->num_fct(), "Input and output vectors do not have "
136  "the same number of functions (" << n << " vs. " << scaledVecOut->num_fct() << ").");
137  for (size_t fct = 0; fct < n; ++fct)
138  {
139  UG_COND_THROW(vecIn->num_dofs(fct) != scaledVecOut->num_dofs(fct),
140  "Input and output vectors do not have the same number of DoFs for function " << fct
141  << " (" << vecIn->num_dofs(fct) << " vs. " << scaledVecOut->num_dofs(fct) << ").");
142  }
143 
144  ConstSmartPtr<DoFDistribution> dd = vecIn->dof_distribution();
145 
146  if (dd->max_dofs(VERTEX))
147  ScaleGFOnElems<Vertex, TGridFunction>(dd, scaledVecOut, vecIn, scalingFactors);
148  if (dd->max_dofs(EDGE))
149  ScaleGFOnElems<Edge, TGridFunction>(dd, scaledVecOut, vecIn, scalingFactors);
150  if (dd->max_dofs(FACE))
151  ScaleGFOnElems<Face, TGridFunction>(dd, scaledVecOut, vecIn, scalingFactors);
152  if (dd->max_dofs(VOLUME))
153  ScaleGFOnElems<Volume, TGridFunction>(dd, scaledVecOut, vecIn, scalingFactors);
154 }
155 
156 
157 
158 
160 // AverageComponent
162 
170 template<typename TGridFunction, typename TBaseElem>
172 {
173  typedef TGridFunction GF;
174  typedef typename GF::template traits<TBaseElem>::const_iterator iter_type;
175 
176  iter_type iter = spGF->template begin<TBaseElem>();
177  iter_type iterEnd = spGF->template end<TBaseElem>();
178 
179  PeriodicBoundaryManager* pbm = spGF->domain()->grid()->periodic_boundary_manager();
180 
181 // loop elems
182  std::vector<DoFIndex> vMultInd;
183  for(; iter != iterEnd; ++iter)
184  {
185  // get element
186  TBaseElem* elem = *iter;
187 
188  // skip periodic ghosts
189  if (pbm && pbm->is_slave(elem)) continue;
190 
191  // get global indices
192  spGF->inner_dof_indices(elem, fct, vMultInd);
193 
194  // sum up value
195  for(size_t i = 0; i < vMultInd.size(); ++i)
196  {
197  DoFRef(*spGF, vMultInd[i]) -= sub;
198  }
199  }
200 }
201 
212 template<typename TGridFunction>
214 (
216  const std::vector<std::string>& vCmp,
217  number mean
218 )
219 {
220  typedef TGridFunction GF;
221  PROFILE_FUNC_GROUP("gmg");
222 
223  if(vCmp.empty())
224  return;
225 
226  if(spGF.invalid())
227  UG_THROW("AverageComponent: expects a valid GridFunction.");
228 
230  spGF->approx_space()->dof_distribution_info();
231 
232 // compute integral of components
233  const number area = Integral(1.0, spGF);
234  std::vector<number> vIntegral(vCmp.size(), 0.0);
235  for(size_t f = 0; f < vCmp.size(); f++){
236  const size_t fct = spGF->fct_id_by_name(vCmp[f].c_str());
237  vIntegral[f] = Integral(spGF, vCmp[f].c_str(), NULL, ddinfo->lfeid(fct).order());
238  }
239 
240 // subtract value
241  for(size_t f = 0; f < vCmp.size(); f++)
242  {
243  const number sub = (vIntegral[f] - mean) / area;
244  const size_t fct = spGF->fct_id_by_name(vCmp[f].c_str());
245 
246  if(ddinfo->max_fct_dofs(fct, VERTEX)) SubtractValueFromComponent<GF, Vertex>(spGF, fct, sub);
247  if(ddinfo->max_fct_dofs(fct, EDGE)) SubtractValueFromComponent<GF, Edge>(spGF, fct, sub);
248  if(ddinfo->max_fct_dofs(fct, FACE)) SubtractValueFromComponent<GF, Face>(spGF, fct, sub);
249  if(ddinfo->max_fct_dofs(fct, VOLUME)) SubtractValueFromComponent<GF, Volume>(spGF, fct, sub);
250  }
251 }
252 
260 template<typename TGridFunction>
261 void AdjustMeanValue(SmartPtr<TGridFunction> spGF, const std::vector<std::string>& vCmp)
262 {
263  AdjustMeanValue(spGF, vCmp, 0.0);
264 }
265 
273 template<typename TGridFunction>
274 void AdjustMeanValue(SmartPtr<TGridFunction> spGF, const std::string& fcts)
275 {
276  AdjustMeanValue(spGF, TokenizeTrimString(fcts), 0.0);
277 }
278 
289 template<typename TGridFunction>
290 void AdjustMeanValue(SmartPtr<TGridFunction> spGF, const std::string& fcts, number mean)
291 {
292  AdjustMeanValue(spGF, TokenizeTrimString(fcts), mean);
293 }
294 
296 // Summing up components
298 
305 template <typename TGridFunction, typename TBaseElem>
307 (
308  TGridFunction * u,
309  size_t fct
310 )
311 {
312  typedef typename TGridFunction::template traits<TBaseElem>::const_iterator t_elem_iter;
313 
314  std::vector<DoFIndex> ind;
315  number sum (1);
316 
317  sum = 0;
318 // Loop the elements in the subset
319  for (t_elem_iter vi = u->template begin<TBaseElem> ();
320  vi != u->template end<TBaseElem> (); ++vi)
321  {
322  TBaseElem * vert = *vi;
323 
324  // indices at this element
325  u->inner_dof_indices (vert, fct, ind);
326 
327  if (ind.size () != 1)
328  UG_THROW ("SumGFValuesAt: The function must be scalar!");
329 
330  // add the contribution
331  sum += DoFRef (*u, ind [0]);
332  }
333 
334 #ifdef UG_PARALLEL
335  pcl::ProcessCommunicator procComm;
336  sum = procComm.allreduce (sum, PCL_RO_SUM);
337 #endif
338 
339  return sum;
340 }
341 
350 template <typename TGridFunction, typename TBaseElem>
352 (
353  TGridFunction * u,
354  const char * fct_names
355 )
356 {
357 // Get the function index
358  std::vector<std::string> vfctNames;
359  TokenizeString (fct_names, vfctNames);
360  if (vfctNames.size () != 1)
361  UG_THROW ("SumGFValuesAt: Exactly one function name must be specified.");
362  FunctionGroup fctGroup (u->function_pattern ());
363  fctGroup.add (vfctNames [0]);
364 
365 // Compute the sum
366  return SumGFValuesAt<TGridFunction, TBaseElem> (u, fctGroup[0]);
367 }
368 
375 template <typename TGridFunction, typename TBaseElem>
377 (
378  TGridFunction * u,
379  size_t fct,
380  SubsetGroup & ssGroup
381 )
382 {
383  typedef typename TGridFunction::template traits<TBaseElem>::const_iterator t_elem_iter;
384 
385  std::vector<DoFIndex> ind;
386  number sum (1);
387 
388 // Loop the selected subsets
389  sum = 0;
390  for (size_t i = 0; i < ssGroup.size (); i++)
391  {
392  int ssi = ssGroup [i];
393 
394  // Loop the elements in the subset
395  for (t_elem_iter vi = u->template begin<TBaseElem> (ssi);
396  vi != u->template end<TBaseElem> (ssi); ++vi)
397  {
398  TBaseElem * vert = *vi;
399 
400  // indices at this element
401  u->inner_dof_indices (vert, fct, ind);
402 
403  if (ind.size () != 1)
404  UG_THROW ("SumGFValuesAt: The function must be scalar!");
405 
406  // add the contribution
407  sum += DoFRef (*u, ind [0]);
408  }
409  }
410 
411 #ifdef UG_PARALLEL
412  pcl::ProcessCommunicator procComm;
413  sum = procComm.allreduce (sum, PCL_RO_SUM);
414 #endif
415 
416  return sum;
417 }
418 
427 template <typename TGridFunction, typename TBaseElem>
429 (
430  TGridFunction * u,
431  const char * fct_names,
432  const char * subset_names
433 )
434 {
435 // Get the function index
436  std::vector<std::string> vfctNames;
437  TokenizeString (fct_names, vfctNames);
438  if (vfctNames.size () != 1)
439  UG_THROW ("SumGFValuesAt: Exactly one function name must be specified.");
440  FunctionGroup fctGroup (u->function_pattern ());
441  fctGroup.add (vfctNames [0]);
442 
443 // Get the subset indices
444  std::vector<std::string> vssNames;
445  TokenizeString (subset_names, vssNames);
446  SubsetGroup ssGroup (u->domain()->subset_handler ());
447  ssGroup.add (vssNames);
448 
449 // Compute the sum
450  return SumGFValuesAt<TGridFunction, TBaseElem> (u, fctGroup[0], ssGroup);
451 }
452 
454 // Checking components for nan and inf
456 
463 template <typename TGridFunction, typename TBaseElem>
465 (
466  const TGridFunction * u,
467  size_t fct
468 )
469 {
470  typedef typename TGridFunction::template traits<TBaseElem>::const_iterator t_elem_iter;
471 
472  std::vector<DoFIndex> ind;
473 
474 // Loop the elements in the subset
475  for (t_elem_iter vi = u->template begin<TBaseElem> ();
476  vi != u->template end<TBaseElem> (); ++vi)
477  {
478  TBaseElem * vert = *vi;
479 
480  // indices at this element
481  u->inner_dof_indices (vert, fct, ind);
482 
483  if (ind.size () != 1)
484  UG_THROW ("CheckGFforNaN: The function must be scalar!");
485 
486  // check the value
487  number value = DoFRef (*u, ind [0]);
488  if (isnan (value))
489  {
490  int si = u->domain()->subset_handler()->get_subset_index (vert);
491  UG_LOG ("nan at index " << ind [0] << ", grid data idx " << vert->grid_data_index ());
492  if (si >= 0)
493  UG_LOG (", subset " << u->domain()->subset_handler()->get_subset_name (si))
494  UG_LOG ('\n');
495  return true;
496  }
497  if (isinf (value))
498  {
499  int si = u->domain()->subset_handler()->get_subset_index (vert);
500  UG_LOG ("inf at index " << ind [0] << ", grid data idx " << vert->grid_data_index ());
501  if (si >= 0)
502  UG_LOG (", subset " << u->domain()->subset_handler()->get_subset_name (si))
503  UG_LOG ('\n');
504  return true;
505  }
506  }
507 
508  return false;
509 }
510 
519 template <typename TGridFunction, typename TBaseElem>
521 (
522  const TGridFunction * u,
523  const char * fct_names
524 )
525 {
526 // Get the function index
527  std::vector<std::string> vfctNames;
528  TokenizeString (fct_names, vfctNames);
529  FunctionGroup fctGroup (u->function_pattern ());
530  for (size_t i = 0; i < vfctNames.size (); i++)
531  fctGroup.add (vfctNames [i]);
532 
533 // Check the functions
534  bool result = false;
535  for (size_t i = 0; i < vfctNames.size (); i++)
536  {
537  UG_LOG ("Checking " << vfctNames[i] << " ... ");
538  if (CheckGFforNaN<TGridFunction, TBaseElem> (u, fctGroup[i]))
539  result = true;
540  else
541  UG_LOG ("OK\n");
542  }
543  return result;
544 }
545 
554 template <typename TGridFunction, typename TBaseElem>
556 (
557  const TGridFunction * u,
558  const FunctionGroup & fctGroup
559 )
560 {
561 // Check the functions
562  bool result = false;
563  for (size_t i = 0; i < fctGroup.size (); i++)
564  {
565  UG_LOG ("Checking fct #" << i << " ... ");
566  if (CheckGFforNaN<TGridFunction, TBaseElem> (u, fctGroup[i]))
567  result = true;
568  else
569  UG_LOG ("OK\n");
570  }
571  return result;
572 }
573 
574 // //////////////////////////////////////////////////////////////////////////////
575 // Check values are within bounds
576 // //////////////////////////////////////////////////////////////////////////////
577 template <typename TGridFunction, typename TBaseElem>
579 (
581  size_t cmp,
582  number lowerBnd,
583  number upperBnd
584 )
585 {
586  typedef typename TGridFunction::template traits<TBaseElem>::const_iterator elem_it;
587 
588  bool ret = true;
589 
590  // loop the elements in the subset
591  std::vector<DoFIndex> vDI;
592  elem_it it = u->template begin<TBaseElem>();
593  elem_it itEnd = u->template end<TBaseElem>();
594  for (; it != itEnd; ++it)
595  {
596  TBaseElem* elem = *it;
597 
598  // loop indices at this element
599  const size_t nInd = u->inner_dof_indices(elem, cmp, vDI, true);
600  for (size_t i = 0; i < nInd; ++i)
601  {
602  const number& val = DoFRef(*u, vDI[i]);
603  if (val < lowerBnd || val > upperBnd)
604  {
605  UG_LOG_ALL_PROCS("Function value for component " << cmp << " (" << val << ") "
606  "is outside the specified range [" << lowerBnd << ", " << upperBnd << "] "
607  "at " << ElementDebugInfo(*u->domain()->grid(), elem) << std::endl);
608  ret = false;
609  }
610  }
611  }
612 
613 #ifdef UG_PARALLEL
614  if (pcl::NumProcs() > 1)
615  {
617  int retInt = ret;
618  ret = pc.allreduce(retInt, PCL_RO_BAND);
619  }
620 #endif
621 
622  return ret;
623 }
624 
625 template <typename TGridFunction>
627 (
629  size_t cmp,
630  number lowerBnd,
631  number upperBnd
632 )
633 {
634  bool ret = true;
635  if (u->max_fct_dofs(cmp, 0))
636  ret = ret && CheckGFValuesWithinBounds<TGridFunction, Vertex>(u, cmp, lowerBnd, upperBnd);
637  if (u->max_fct_dofs(cmp, 1))
638  ret = ret && CheckGFValuesWithinBounds<TGridFunction, Edge>(u, cmp, lowerBnd, upperBnd);
639  if (u->max_fct_dofs(cmp, 2))
640  ret = ret && CheckGFValuesWithinBounds<TGridFunction, Face>(u, cmp, lowerBnd, upperBnd);
641  if (u->max_fct_dofs(cmp, 3))
642  ret = ret && CheckGFValuesWithinBounds<TGridFunction, Volume>(u, cmp, lowerBnd, upperBnd);
643 
644  return ret;
645 }
646 
647 template <typename TGridFunction>
649 (
651  const char* fctNames,
652  number lowerBnd,
653  number upperBnd
654 )
655 {
656  std::vector<std::string> vFctNames;
657  TokenizeTrimString (fctNames, vFctNames);
658  FunctionGroup fctGroup (u->function_pattern());
659  const size_t nFct = vFctNames.size();
660  for (size_t f = 0; f < nFct; ++f)
661  {
662  try {fctGroup.add(vFctNames[f]);}
663  UG_CATCH_THROW("Could not add function " << vFctNames[f] << " to function group.");
664  }
665 
666  bool ret = true;
667  const size_t fctGrpSz = fctGroup.size();
668  for (size_t f = 0; f < fctGrpSz; ++f)
669  ret = ret && CheckGFValuesWithinBounds(u, fctGroup[f], lowerBnd, upperBnd);
670 
671  return ret;
672 }
673 
674 
675 
676 
678 // Writing algebra
680 
681 template<typename TDomain>
683 {
684  PROFILE_FUNC_GROUP("debug");
685 
686 #ifdef UG_PARALLEL
687  name = ConnectionViewer::GetParallelName(name, dd->layouts()->proc_comm());
688 #endif
689 
690  std::vector<size_t> fctIndex;
691  std::vector<std::string> fctNames;
692 
693  ExtractAlgebraIndices<TDomain>(domain, dd, fctIndex);
694 
695  size_t numFct = dd->num_fct();
696  if(numFct <= 1) return;
697 
698  fctNames.resize(numFct);
699  for(size_t i=0; i<numFct; i++)
700  fctNames[i] = dd->name(i);
701 
702  name.append(".indices");
703  std::fstream file(name.c_str(), std::ios::out);
704 
705  //std::cout << "name is " << name << "\n";
706  file << "NUMDOF " << fctNames.size() << "\n";
707  for(size_t i=0; i<numFct; i++)
708  file << fctNames[i] << "\n";
709 
710  for(size_t i=0; i<fctIndex.size(); i++)
711  file << fctIndex[i] << "\n";
712 }
713 
714 template<class TFunction>
715 void WriteMatrixToConnectionViewer(const char *filename,
716  const typename TFunction::algebra_type::matrix_type &A,
717  const TFunction &u) {
718 
719  PROFILE_FUNC_GROUP("debug");
720 // check name
721  if ( !FileTypeIs( filename, ".mat") ) {
722  UG_THROW( "Only '.mat' format supported for matrices, but"
723  " filename is '" << filename << "'." );
724  }
725 
726 // position array
727  const static int dim = TFunction::domain_type::dim;
728  std::vector<MathVector<dim> > vPos;
729  ExtractPositions(u, vPos);
730 
731 // write matrix
732  if(vPos.empty())
733  ConnectionViewer::WriteMatrixPar( filename, A, (MathVector<dim>*)NULL, dim );
734  else
735  ConnectionViewer::WriteMatrixPar( filename, A, &vPos[0], dim );
736 
737  WriteAlgebraIndices(filename, u.domain(),u.dof_distribution());
738 
739 }
740 
741 template<typename TGridFunction>
743  TGridFunction& u,
744  MatrixOperator<typename TGridFunction::algebra_type::matrix_type,
745  typename TGridFunction::vector_type>& A,
746  const char* filename) {
747  PROFILE_FUNC_GROUP("debug");
748 // forward
749  WriteMatrixToConnectionViewer(filename, A.get_matrix(), u);
750 }
751 
761 // TODO extend MatrixIO to support other than CPUAlgebra
762 inline void SaveMatrixToMTX( const char *filename,
764  std::string comment="%Generated with ug4." ) {
765  PROFILE_FUNC_GROUP("debug");
766 
767  // check name
768  if ( !FileTypeIs( filename, ".mtx" ) ) {
769  UG_THROW( "Please use '.mtx' as file extension for MatrixMarket exchange files."
770  " (Filename is '" << filename << "')" );
771  }
772 
773  // automatically detect the file mode to use
774  MatrixIO::OpenMode openMode = (FileExists( filename )) ? MatrixIO::EXISTING : MatrixIO::NEW;
775 
776  // create MatrixIO object for handling the export
777  MatrixIOMtx mtx( filename, openMode );
778  mtx.write_from( A.get_matrix(), comment );
779 }
780 
781 template<class TFunction>
782 void WriteVectorToConnectionViewer(const char *filename,
783  const typename TFunction::algebra_type::vector_type &b,
784  const TFunction &u,
785  const typename TFunction::algebra_type::vector_type *pCompareVec = NULL) {
786  PROFILE_FUNC_GROUP("debug");
787 // check name
788  if ( !FileTypeIs( filename, ".vec") ) {
789  UG_THROW( "Only '.vec' format supported for vectors, but"
790  " filename is '" << filename << "'." );
791  }
792 
793 // get positions of vertices
794  const static int dim = TFunction::domain_type::dim;
795  std::vector<MathVector<dim> > vPos;
796  ExtractPositions(u, vPos);
797 
798 // write vector
799  ConnectionViewer::WriteVectorPar( filename, b, &vPos[0], dim, pCompareVec);
800  WriteAlgebraIndices(filename, u.domain(),u.dof_distribution());
801 }
802 
803 template<class TFunction>
805  const char *filename,
806  const typename TFunction::algebra_type::matrix_type &A,
807  const typename TFunction::algebra_type::vector_type &b,
808  const TFunction &u,
809  const typename TFunction::algebra_type::vector_type *pCompareVec = NULL) {
810  PROFILE_FUNC_GROUP("debug");
811 // get dimension
812  const static int dim = TFunction::domain_type::dim;
813 
814 // check name
815  if ( !FileTypeIs( filename, ".vec") ) {
816  UG_THROW( "Only '.vec' format supported for vectors." );
817  }
818 
819 // get positions of vertices
820  std::vector<MathVector<dim> > positions;
821  ExtractPositions(u, positions);
822 
823 // write vector
824  ConnectionViewer::WriteVectorPar( filename, A, b, &positions[0], dim, pCompareVec );
825 }
826 
827 template<typename TGridFunction>
828 void SaveVectorForConnectionViewer(TGridFunction& b, const char* filename) {
829  WriteVectorToConnectionViewer(filename, b, b);
830 }
831 
832 template<typename TGridFunction>
833 void SaveVectorDiffForConnectionViewer(TGridFunction& b, TGridFunction& bCompare, const char* filename) {
834  WriteVectorToConnectionViewer(filename, b, b, &bCompare);
835 }
836 
837 template<typename TGridFunction>
839  TGridFunction& u,
840  MatrixOperator<typename TGridFunction::algebra_type::matrix_type,
841  typename TGridFunction::vector_type>& A,
842  const char* filename) {
843  WriteVectorToConnectionViewer(filename, A.get_matrix(), u, u);
844 }
845 
846 template<typename TGridFunction>
848  TGridFunction& u,
849  TGridFunction& compareVec,
850  MatrixOperator<typename TGridFunction::algebra_type::matrix_type,
851  typename TGridFunction::vector_type>& A,
852  const char* filename) {
853 // forward
854  WriteVectorToConnectionViewer(filename, A.get_matrix(), u, u, &compareVec);
855 }
856 
857 // from connection_viewer_input.h
858 // with additional checks
859 template<typename vector_type>
860 bool ReadVector(std::string filename, vector_type &vec,int dim)
861 {
862  Progress p;
863  std::cout << " Reading std::vector from " << filename << "... ";
864  std::fstream matfile(filename.c_str(), std::ios::in);
865  if(matfile.is_open() == false) { std::cout << "failed.\n"; return false; }
866 
867  int version=-1, dimension=-1, gridsize;
868 
869  matfile >> version;
870  matfile >> dimension;
871  matfile >> gridsize;
872 
873  assert(version == 1);
874  assert(dimension == dim);
875  // todo check positions and not just size
876  assert(gridsize == (int)vec.size());
877 
878 
879  PROGRESS_START(prog, gridsize*2, "ReadVector " << dimension << "d from " << filename << " , " << gridsize << " x " << gridsize);
880  for(int i=0; i<gridsize; i++)
881  {
882  if(i%100) { PROGRESS_UPDATE(prog, i); }
883  if(matfile.eof())
884  {
885  std::cout << " failed.\n";
886  assert(0);
887  return false;
888  }
889  double x, y, z;
890  matfile >> x >> y;
891  if(dimension==3) matfile >> z;
892  }
893 
894  int printStringsInWindow;
895  matfile >> printStringsInWindow;
896 
897  // vec.resize(gridsize);
898  bool bEOF = matfile.eof();
899  while(!bEOF)
900  {
901  int from, to; double value;
902  char c = matfile.peek();
903  if(c == -1 || c == 'c' || c == 'v' || matfile.eof())
904  break;
905 
906  matfile >> from >> to >> value;
907  assert(from == to);
908  vec[from] = value;
909  if(from%100) { PROGRESS_UPDATE(prog, from); }
910  bEOF = matfile.eof();
911  }
912  return true;
913 }
914 
915 // load vector that has been saved in connection viewer format and write it
916 // into grid function
917 template<typename TGridFunction>
918 void LoadVector(TGridFunction& u,const char* filename){
919  PROFILE_FUNC();
920  typename TGridFunction::algebra_type::vector_type b;
921  b.resize(u.num_indices());
922  ReadVector(filename,b,TGridFunction::dim);
923  u.assign(b);
924 }
925 
926 // Same as before, but for comma separated value (CSV)
927 template<class TFunction>
928 void WriteVectorCSV(const char *filename,
929  const typename TFunction::algebra_type::vector_type &b,
930  const TFunction &u) {
931  PROFILE_FUNC_GROUP("debug");
932 // get dimension
933  const static int dim = TFunction::domain_type::dim;
934 
935 // check name
936  if ( !FileTypeIs( filename, ".csv") ) {
937  UG_THROW( "Only '.csv' format supported for vectors, but"
938  " filename is '" << filename << "'." );
939  }
940 
941 // extended filename
942 // add p000X extension in parallel
943 #ifdef UG_PARALLEL
944  std::string name(filename);
945  size_t iExtPos = name.find_last_of(".");
946  name.resize(iExtPos);
947  int rank = pcl::ProcRank();
948  char ext[20];
949  sprintf(ext, "_p%05d.csv", rank);
950  name.append(ext);
951 #endif
952 
953 // get positions of vertices
954  std::vector<MathVector<dim> > positions;
955  ExtractPositions(u, positions);
956 
957 // write vector
958  WriteVectorCSV( filename, b, &positions[0], dim );
959 }
960 
961 template<typename TGridFunction>
962 void SaveVectorCSV(TGridFunction& b, const char* filename) {
963  PROFILE_FUNC();
964  WriteVectorCSV(filename, b, b);
965 }
966 
980 template<typename TDomain, typename TAlgebra>
983  std::string subset, std::string fct1, std::string fct2 )
984 {
985  PROFILE_FUNC();
986  // get subset index
987  size_t subSetID = spGridFct->subset_id_by_name( subset.c_str() );
988 
989  // get function indices
990  size_t fct1ID = spGridFct->fct_id_by_name( fct1.c_str() );
991  size_t fct2ID = spGridFct->fct_id_by_name( fct2.c_str() );
992 
993  // create space for sum of difference
994  number sum = 0.0;
995  size_t numElements = 0;
996 
997  // loop over all vertices in given subset and compare values of fct1 and fct2
998  typedef typename GridFunction<TDomain, TAlgebra>::template traits<Vertex>::const_iterator gridFctIterator;
999  for( gridFctIterator iter = spGridFct->template begin<Vertex>((int)subSetID);
1000  iter != spGridFct->template end<Vertex>((int)subSetID); ++iter ) {
1001  // get dof_indices for the two functions on given subset
1002  std::vector< DoFIndex > indFct1, indFct2;
1003  spGridFct->template dof_indices<Vertex>( *iter, fct1ID, indFct1 );
1004  spGridFct->template dof_indices<Vertex>( *iter, fct2ID, indFct2 );
1005 
1006  // calculate the difference between the two functions at this grid point
1007  sum += DoFRef( *spGridFct, indFct1[0] ) - DoFRef( *spGridFct, indFct2[0] );
1008  numElements++;
1009  }
1010 
1011  // return overal arithmetic average of differences
1012  return sum / numElements;
1013 }
1014 
1015 
1017 // Grid Function Debug Writer
1019 
1020 template<typename TDomain, typename TAlgebra>
1021 class GridFunctionDebugWriter: public IDebugWriter<TAlgebra>
1022 {
1024  static const int dim = TDomain::dim;
1025 
1026 public:
1028  typedef TAlgebra algebra_type;
1029 
1031  typedef typename algebra_type::vector_type vector_type;
1032 
1034  typedef typename algebra_type::matrix_type matrix_type;
1035 
1038 
1040  using super::get_base_dir;
1041 
1042 public:
1045  SmartPtr<ApproximationSpace<TDomain> > spApproxSpace) :
1046  m_spApproxSpace(spApproxSpace), bConnViewerOut(
1047  true), bConnViewerIndices(false), bVTKOut(true), m_printConsistent(true)
1048  {
1050  reset();
1051  }
1052 
1054 
1056  void set_grid_level(const GridLevel& gridLevel) {
1057  m_glFrom = gridLevel; m_glTo = gridLevel;
1058  }
1059 
1061  void set_grid_levels(const GridLevel& glFrom, const GridLevel& glTo) {
1062  m_glFrom = glFrom; m_glTo = glTo;
1063  }
1064 
1066  GridLevel grid_level() const {return m_glFrom;}
1067 
1069  void reset() {
1071  }
1072 
1073 
1075  void set_vtk_output(bool b) {bVTKOut = b;}
1076 
1079 
1082 
1085 
1087  virtual void write_vector(const vector_type& vec, const char* filename) {
1088  // write to conn viewer
1089  if (bConnViewerOut)
1090  write_vector_to_conn_viewer(vec, filename);
1091 
1092  // write to vtk
1093  if (bVTKOut)
1094  write_vector_to_vtk(vec, filename);
1095  }
1096 
1097  virtual void update_positions()
1098  {
1100  }
1101 
1103  virtual void write_matrix(const matrix_type& mat, const char* filename) {
1104  PROFILE_FUNC_GROUP("debug");
1105  // something to do ?
1106  if (!bConnViewerOut)
1107  return;
1108 
1109  update_positions();
1110 
1111  std::string name;
1112  this->compose_file_path (name);
1113  name += filename;
1114 
1115  if ( !FileTypeIs( filename, ".mat" ) ) {
1116  UG_THROW( "Only '.mat' format supported for matrices, but"
1117  " filename is '" << filename << "'." );
1118  }
1119 
1120  // write to file
1122  if(m_glFrom == m_glTo){
1123  if(mat.num_rows() != mat.num_cols())
1124  UG_THROW("DebugWriter: grid level the same, but non-square matrix.");
1125 
1126  const std::vector<MathVector<dim> >& vPos = this->template get_positions<dim>();
1127  if(vPos.empty())
1129  else
1130  ConnectionViewer::WriteMatrixPar(name.c_str(), mat, &vPos[0], dim);
1131  }
1132  else{
1133  const std::vector<MathVector<dim> >& vFromPos = this->template get_positions<dim>();
1134 
1135  std::vector<MathVector<dim> > vToPos;
1136  ExtractPositions<TDomain>(m_spApproxSpace->domain(),
1137  m_spApproxSpace->dof_distribution(m_glTo),
1138  vToPos);
1139 
1140  if(mat.num_cols() == vFromPos.size() && mat.num_rows() == vToPos.size())
1141  {
1142  ConnectionViewer::WriteMatrixPar(name, mat, vFromPos, vToPos, dim);
1143  }
1144  else{
1145  UG_THROW("GridFunctionDebugWriter: Wrong size of matrix for writing"
1146  "matrix ("<<m_glTo<<" x "<<m_glFrom<<"), that would be a ("
1147  <<vToPos.size()<<" x "<<vFromPos.size()<<") matrix. But "
1148  "passed matrix of size ("<<mat.num_rows()<<" x "
1149  <<mat.num_cols()<<").");
1150  }
1151  }
1152  }
1153 
1154 protected:
1155 
1156  void write_algebra_indices_CV(std::string name)
1157  {
1158  if(bConnViewerIndices)
1159  WriteAlgebraIndices<TDomain>(name, m_spApproxSpace->domain(), m_spApproxSpace->dof_distribution(m_glFrom));
1160  }
1161 
1163  void extract_positions(const GridLevel& gridLevel) {
1164  // extract positions for this grid function
1165  std::vector<MathVector<dim> >& vPos =
1166  this->template positions<dim>();
1167 
1168  vPos.clear();
1169  ExtractPositions<TDomain>(
1170  m_spApproxSpace->domain(),
1171  m_spApproxSpace->dof_distribution(gridLevel),
1172  vPos);
1173  }
1174 
1176  virtual void write_vector_to_conn_viewer(const vector_type& vec,
1177  const char* filename)
1178  {
1179  PROFILE_FUNC_GROUP("debug");
1180  update_positions();
1181 
1182  std::string name;
1183  this->compose_file_path (name);
1184  name += filename;
1185 
1186  if ( !FileTypeIs( filename, ".vec" ) ) {
1187  UG_THROW( "Only '.vec' format supported for vectors, but"
1188  " filename is '" << name << "'.");
1189  }
1190 
1191  // write
1192  write_algebra_indices_CV(filename);
1193  const std::vector<MathVector<dim> >& vPos =
1194  this->template get_positions<dim>();
1195  if(vPos.empty())
1197  else
1198  ConnectionViewer::WriteVectorPar(name.c_str(), vec, &vPos[0], dim);
1199  }
1200 
1201  void write_vector_to_vtk(const vector_type& vec, const char* filename) {
1202  PROFILE_FUNC_GROUP("debug");
1203  // check name
1204 
1205  std::string name;
1206  this->compose_file_path (name);
1207  name += filename;
1208 
1209  typedef GridFunction<TDomain, TAlgebra> TGridFunction;
1210  TGridFunction vtkFunc(
1212  m_spApproxSpace->dof_distribution(m_glFrom));
1213  vtkFunc.resize_values(vec.size());
1214  vtkFunc.assign(vec);
1215  VTKOutput<dim> out;
1216  out.print(name.c_str(), vtkFunc, m_printConsistent);
1217  }
1218 
1219 protected:
1220  // underlying approximation space
1222 
1223  // flag if write to conn viewer
1226 
1227  // flag if write to vtk
1228  bool bVTKOut;
1229 
1230  // flag if data shall be made consistent before printing
1232 
1233  // current grid level
1235 };
1236 
1238 // Grid Function Position Provider
1240 
1241 template<typename TGridFunction>
1244 public:
1247  m_pGridFunc(NULL) {
1248  }
1249  GridFunctionPositionProvider(const TGridFunction& u) :
1250  m_pGridFunc(&u) {
1251  }
1252  void set_reference_grid_function(const TGridFunction& u) {
1253  m_pGridFunc = &u;
1254  }
1255 
1256  virtual bool get_positions(
1257  std::vector<MathVector<TGridFunction::domain_type::dim> >&vec) {
1258  UG_ASSERT(m_pGridFunc != NULL,
1259  "provide a grid function with set_reference_grid_function");
1260  ExtractPositions(*m_pGridFunc, vec);
1261  return true;
1262  }
1263 
1264 protected:
1265  const TGridFunction *m_pGridFunc;
1266 };
1267 
1269 // Grid Function Vector Writer
1271 
1272 template<typename TGridFunction, typename TVector>
1273 class GridFunctionVectorWriter: public IVectorWriter<TVector> {
1274 public:
1276  typedef typename TGridFunction::domain_type domain_type;
1277  typedef TVector vector_type;
1278 
1279 public:
1282  m_pGridFunc(NULL) {
1283  }
1284 
1286  m_userData = userData;
1287  }
1288 
1289  void set_reference_grid_function(const TGridFunction& u) {
1290  m_pGridFunc = &u;
1291  }
1292 
1293  /*virtual double calculate(MathVector<3> &v, double time)
1294  {
1295  }*/
1296 
1297  virtual bool update(vector_type &vec) {
1298  PROFILE_FUNC_GROUP("debug");
1299  UG_ASSERT(m_pGridFunc != NULL,
1300  "provide a grid function with set_reference_grid_function");
1301  UG_ASSERT(m_userData.valid(), "provide user data with set_user_data");
1302 
1303  const TGridFunction &u = *m_pGridFunc;
1304 
1305  // get position accessor
1306 
1307  const typename domain_type::position_accessor_type& aaPos =
1308  u.domain()->position_accessor();
1309 
1310  // number of total dofs
1311  int nr = u.num_indices();
1312 
1313  // resize positions
1314  vec.resize(nr);
1315 
1316  typedef typename TGridFunction::template traits<Vertex>::const_iterator const_iterator;
1317 
1318  // loop all subsets
1319  for (int si = 0; si < u.num_subsets(); ++si) {
1320  // loop all vertices
1321  const_iterator iter = u.template begin<Vertex>(si);
1322  const_iterator iterEnd = u.template end<Vertex>(si);
1323 
1324  for (; iter != iterEnd; ++iter) {
1325  // get vertex
1326  Vertex* v = *iter;
1327 
1328  // algebra indices vector
1329  std::vector<size_t> ind;
1330 
1331  // load indices associated with vertex
1332  u.inner_algebra_indices(v, ind);
1333 
1334  number t = 0.0;
1335 
1336  number d;
1337  (*m_userData)(d, aaPos[v], t, si);
1338 
1339  // write
1340  for (size_t i = 0; i < ind.size(); ++i) {
1341  const size_t index = ind[i];
1342  for (size_t alpha = 0; alpha < GetSize(vec[index]); ++alpha)
1343  BlockRef(vec[index], alpha) = d;
1344  }
1345  }
1346  }
1347  return true;
1348  }
1349 
1350 protected:
1351  const TGridFunction *m_pGridFunc;
1353 };
1354 
1356 // Grid Function Vector Writer Dirichlet 0
1358 
1359 template<typename TGridFunction>
1361  typename TGridFunction::algebra_type::vector_type> {
1362 public:
1363  typedef typename TGridFunction::domain_type domain_type;
1364 
1365  typedef typename TGridFunction::approximation_space_type approximation_space_type;
1366 
1367  typedef typename TGridFunction::algebra_type algebra_type;
1368  typedef typename algebra_type::vector_type vector_type;
1370 
1371 public:
1374  m_pApproxSpace(NULL), m_spPostProcess(NULL), m_level(-1) {
1375  }
1376 
1377  void set_level(size_t level) {
1378  m_level = level;
1379  }
1380 
1382  approximation_space_type& approxSpace) {
1383  m_spPostProcess = pp;
1384  m_pApproxSpace = &approxSpace;
1385  }
1386 
1387  /*virtual double calculate(MathVector<3> &v, double time)
1388  {
1389  }*/
1390 
1391  virtual bool update(vector_type &vec) {
1392  PROFILE_FUNC_GROUP("debug");
1393  UG_ASSERT(m_spPostProcess.valid(), "provide a post process with init");
1394  UG_ASSERT(m_pApproxSpace != NULL, "provide approximation space init");
1395 
1396  size_t numDoFs = 0;
1397  GridLevel gl;
1398  if (m_level == (size_t) -1) {
1399  numDoFs = m_pApproxSpace->dof_distribution(GridLevel(GridLevel::TOP, GridLevel::SURFACE))->num_indices();
1400  gl = GridLevel();
1401  } else {
1402  numDoFs =
1403  m_pApproxSpace->dof_distribution(GridLevel(m_level, GridLevel::LEVEL, true))->num_indices();
1405  }
1406  vec.resize(numDoFs);
1407  vec.set(1.0);
1408 
1409  m_spPostProcess->adjust_defect(vec, vec, m_pApproxSpace->dof_distribution(gl), CT_DIRICHLET);
1410 
1411  return true;
1412  }
1413 
1414 protected:
1417  size_t m_level;
1418 };
1419 
1420 } // end namespace ug
1421 
1422 #endif /* __H__UG__LIB_DISC__FUNCTION_SPACE__GRID_FUNCTION_UTIL__ */
parameterString p
location name
Definition: checkpoint_util.lua:128
Definition: smart_pointer.h:296
bool invalid() const
returns true if the pointer is invalid, false if not.
Definition: smart_pointer.h:212
Definition: pcl_process_communicator.h:70
void allreduce(const void *sendBuf, void *recBuf, int count, DataType type, ReduceOperation op) const
performs MPI_Allreduce on the processes of the communicator.
Definition: pcl_process_communicator.cpp:318
base class for approximation spaces without type of algebra or dof distribution
Definition: approximation_space.h:279
Definition: function_group.h:52
void add(size_t fct)
adds a function by id to this group
Definition: function_group.cpp:79
size_t size() const
number of functions in this group
Definition: function_group.h:116
Definition: grid_function_util.h:1022
bool bConnViewerOut
Definition: grid_function_util.h:1224
ApproximationSpace< TDomain > approximation_space_type
type of approximation space
Definition: grid_function_util.h:1037
GridFunctionDebugWriter(SmartPtr< ApproximationSpace< TDomain > > spApproxSpace)
Constructor.
Definition: grid_function_util.h:1044
static const int dim
dimension
Definition: grid_function_util.h:1024
algebra_type::vector_type vector_type
type of vector
Definition: grid_function_util.h:1031
void set_print_consistent(bool b)
sets if data shall be made consistent before printing
Definition: grid_function_util.h:1084
virtual void write_matrix(const matrix_type &mat, const char *filename)
write matrix
Definition: grid_function_util.h:1103
GridLevel m_glTo
Definition: grid_function_util.h:1234
bool bConnViewerIndices
Definition: grid_function_util.h:1225
void set_grid_levels(const GridLevel &glFrom, const GridLevel &glTo)
sets the grid level
Definition: grid_function_util.h:1061
void set_grid_level(const GridLevel &gridLevel)
sets the grid level
Definition: grid_function_util.h:1056
virtual void write_vector(const vector_type &vec, const char *filename)
write vector
Definition: grid_function_util.h:1087
GridLevel grid_level() const
returns current grid level
Definition: grid_function_util.h:1066
void set_vtk_output(bool b)
sets if writing to vtk
Definition: grid_function_util.h:1075
SmartPtr< approximation_space_type > m_spApproxSpace
Definition: grid_function_util.h:1221
bool bVTKOut
Definition: grid_function_util.h:1228
IDebugWriter< TAlgebra > super
Definition: grid_function_util.h:1039
GridLevel m_glFrom
Definition: grid_function_util.h:1234
void set_conn_viewer_indices(bool b)
sets if .indices file is written or conn viewer
Definition: grid_function_util.h:1081
bool m_printConsistent
Definition: grid_function_util.h:1231
virtual ~GridFunctionDebugWriter()
Definition: grid_function_util.h:1053
void write_algebra_indices_CV(std::string name)
Definition: grid_function_util.h:1156
virtual void write_vector_to_conn_viewer(const vector_type &vec, const char *filename)
write vector
Definition: grid_function_util.h:1176
void reset()
sets to toplevel on surface
Definition: grid_function_util.h:1069
TAlgebra algebra_type
type of matrix
Definition: grid_function_util.h:1028
void set_conn_viewer_output(bool b)
sets if writing to conn viewer
Definition: grid_function_util.h:1078
algebra_type::matrix_type matrix_type
type of matrix
Definition: grid_function_util.h:1034
virtual void update_positions()
Definition: grid_function_util.h:1097
void extract_positions(const GridLevel &gridLevel)
reads the positions
Definition: grid_function_util.h:1163
void write_vector_to_vtk(const vector_type &vec, const char *filename)
Definition: grid_function_util.h:1201
represents numerical solutions on a grid using an algebraic vector
Definition: grid_function.h:121
Definition: grid_function_util.h:1243
GridFunctionPositionProvider()
Constructor.
Definition: grid_function_util.h:1246
GridFunctionPositionProvider(const TGridFunction &u)
Definition: grid_function_util.h:1249
void set_reference_grid_function(const TGridFunction &u)
Definition: grid_function_util.h:1252
virtual bool get_positions(std::vector< MathVector< TGridFunction::domain_type::dim > > &vec)
Definition: grid_function_util.h:1256
const TGridFunction * m_pGridFunc
Definition: grid_function_util.h:1265
Definition: grid_function_util.h:1361
TGridFunction::domain_type domain_type
Definition: grid_function_util.h:1363
algebra_type::vector_type vector_type
Definition: grid_function_util.h:1368
GridFunctionVectorWriterDirichlet0()
Constructor.
Definition: grid_function_util.h:1373
vector_type::value_type value_type
Definition: grid_function_util.h:1369
virtual bool update(vector_type &vec)
Definition: grid_function_util.h:1391
void set_level(size_t level)
Definition: grid_function_util.h:1377
TGridFunction::approximation_space_type approximation_space_type
Definition: grid_function_util.h:1365
size_t m_level
Definition: grid_function_util.h:1417
SmartPtr< IConstraint< algebra_type > > m_spPostProcess
Definition: grid_function_util.h:1416
TGridFunction::algebra_type algebra_type
Definition: grid_function_util.h:1367
approximation_space_type * m_pApproxSpace
Definition: grid_function_util.h:1415
void init(SmartPtr< IConstraint< algebra_type > > pp, approximation_space_type &approxSpace)
Definition: grid_function_util.h:1381
Definition: grid_function_util.h:1273
void set_user_data(SmartPtr< UserData< number, domain_type::dim > > userData)
Definition: grid_function_util.h:1285
virtual bool update(vector_type &vec)
write vector
Definition: grid_function_util.h:1297
GridFunctionVectorWriter()
Constructor.
Definition: grid_function_util.h:1281
void set_reference_grid_function(const TGridFunction &u)
Definition: grid_function_util.h:1289
TVector vector_type
Definition: grid_function_util.h:1277
const TGridFunction * m_pGridFunc
Definition: grid_function_util.h:1351
TGridFunction::domain_type domain_type
Definition: grid_function_util.h:1276
SmartPtr< UserData< number, domain_type::dim > > m_userData
Definition: grid_function_util.h:1352
TVector::value_type value_type
Definition: grid_function_util.h:1275
Definition: grid_level.h:42
@ TOP
Definition: grid_level.h:45
@ SURFACE
Definition: grid_level.h:48
@ LEVEL
Definition: grid_level.h:48
interface for adjustment of constraints
Definition: constraint_interface.h:63
base class for all debug writer
Definition: debug_writer.h:244
TAlgebra::vector_type vector_type
type of vector
Definition: debug_writer.h:250
TAlgebra::matrix_type matrix_type
type of matrix
Definition: debug_writer.h:253
Interface for providing vertex positions.
Definition: debug_writer.h:46
base class for all vector debug writer
Definition: debug_writer.h:118
std::string get_base_dir()
Definition: debug_writer.h:192
void compose_file_path(std::string &path)
composes the path for the files and creates the intermediate directories (up to the base one):
Definition: debug_writer.h:210
Interface for modifying a vector (e.g, setting Dirichlet values, ...)
Definition: vector_writer.h:42
OpenMode
Specifies how file names pointing to non-existing files should be handeld.
Definition: matrix_io.h:129
@ EXISTING
Only existing files are valid files (i.e. non-existing file result in error)
Definition: matrix_io.h:131
@ NEW
Non-existing files are created with the specified file name.
Definition: matrix_io.h:133
Provides I/O functionality for MatrixMarket exchange file format.
Definition: matrix_io_mtx.h:84
void write_from(matrix_type &matrix, std::string comment="%Generated with ug4.")
Writes a ug matrix into the associated MatrixMarket exchange file.
Definition: matrix_io_mtx.h:287
Definition: matrix_operator.h:49
virtual M & get_matrix()
Definition: matrix_operator.h:74
Definition: periodic_boundary_manager.h:99
bool is_slave(TElem *) const
Definition: periodic_boundary_manager_impl.hpp:498
Definition: progress.h:50
Group of subsets.
Definition: subset_group.h:51
size_t size() const
number of subsets in this group
Definition: subset_group.h:122
void add(int si)
adds a subset by number to this group
Definition: subset_group.cpp:64
Const iterator to traverse the surface of a multi-grid hierarchy.
Definition: surface_view.h:237
@ ALL
Definition: surface_view.h:102
Type based UserData.
Definition: user_data.h:117
output writer to the VTK file format
Definition: vtkoutput.h:150
void print(const char *filename, TFunction &u, int step, number time, bool makeConsistent)
Definition: vtkoutput_impl.h:184
Base-class for all vertex-types.
Definition: grid_base_objects.h:231
File utility functions.
std::string ElementDebugInfo(const Grid &grid, GridObject *e)
Returns a string containing information on the given element.
Definition: debug_util.cpp:991
#define PCL_RO_SUM
Definition: pcl_methods.h:63
int ProcRank()
returns the rank of the process
Definition: pcl_base.cpp:83
int NumProcs()
returns the number of processes
Definition: pcl_base.cpp:91
#define PCL_RO_BAND
Definition: pcl_methods.h:66
function common sprintf(s,...)
static const int dim
number alpha
base_type::TBaseElem TBaseElem
UG_API bool FileExists(const char *filename)
!!! Serial i/o version !!!
Definition: file_util.cpp:53
bool FileTypeIs(const char *filename, const char *extension)
!!! Serial i/o version !!!
Definition: file_util.cpp:86
#define UG_LOG_ALL_PROCS(msg)
Definition: log.h:371
#define UG_ASSERT(expr, msg)
Definition: assert.h:70
#define UG_CATCH_THROW(msg)
Definition: error.h:64
#define UG_THROW(msg)
Definition: error.h:57
#define UG_LOG(msg)
Definition: log.h:367
#define UG_COND_THROW(cond, msg)
UG_COND_THROW(cond, msg) : performs a UG_THROW(msg) if cond == true.
Definition: error.h:61
double number
Definition: types.h:124
void WriteVectorPar(std::string filename, const Vector_type &b, const postype *positions, int dimensions, const Vector_type *compareVec=NULL)
Definition: connection_viewer_output.h:421
void WriteMatrixPar(std::string name, const Matrix_type &A, const postype *positions, int dimensions)
Definition: connection_viewer_output.h:150
string GetParallelName(string name, const pcl::ProcessCommunicator &pc, bool bWriteHeader)
Definition: connection_viewer_output.cpp:51
CPUAlgebra::vector_type vector_type
the ug namespace
bool CheckGFforNaN(const TGridFunction *u, size_t fct)
Definition: grid_function_util.h:465
bool ReadVector(std::string filename, vector_type &vec, int dim)
Definition: grid_function_util.h:860
@ CT_DIRICHLET
Definition: ass_tuner.h:59
number SumGFValuesAt(TGridFunction *u, size_t fct)
Definition: grid_function_util.h:307
void SaveMatrixForConnectionViewer(TGridFunction &u, MatrixOperator< typename TGridFunction::algebra_type::matrix_type, typename TGridFunction::vector_type > &A, const char *filename)
Definition: grid_function_util.h:742
void ScaleGF(SmartPtr< TGridFunction > scaledVecOut, ConstSmartPtr< TGridFunction > vecIn, const std::vector< number > &scalingFactors)
Scales all functions contained in a grid function.
Definition: grid_function_util.h:123
number & DoFRef(TMatrix &mat, const DoFIndex &iInd, const DoFIndex &jInd)
Definition: multi_index.h:276
void SaveVectorForConnectionViewer(TGridFunction &b, const char *filename)
Definition: grid_function_util.h:828
void ExtractPositions(ConstSmartPtr< TDomain > domain, ConstSmartPtr< DoFDistribution > dd, std::vector< MathVector< TDomain::dim > > &vPos)
Definition: dof_position_util.cpp:424
void TokenizeTrimString(const string &str, vector< string > &vToken, const char delimiter)
Definition: string_util.cpp:83
void SaveMatrixToMTX(const char *filename, MatrixOperator< CPUAlgebra::matrix_type, CPUAlgebra::vector_type > &A, std::string comment="%Generated with ug4.")
Save the assembled matrix of a matrix operator to MatrixMarket format.
Definition: grid_function_util.h:762
void WriteVectorCSV(const char *filename, const Vector_type &b, postype *positions, int dimensions)
writes to a file in somewhat SparseMatrix-market format (for connection viewer)
Definition: csv_gnuplot_output.h:59
void TokenizeString(const string &str, vector< string > &vToken, const char delimiter)
Definition: string_util.cpp:56
void WriteAlgebraIndices(std::string name, ConstSmartPtr< TDomain > domain, ConstSmartPtr< DoFDistribution > dd)
Definition: grid_function_util.h:682
void SaveVectorDiffForConnectionViewer(TGridFunction &b, TGridFunction &bCompare, const char *filename)
Definition: grid_function_util.h:833
void LoadVector(TGridFunction &u, const char *filename)
Definition: grid_function_util.h:918
void SubtractValueFromComponent(SmartPtr< TGridFunction > spGF, size_t fct, number sub)
Definition: grid_function_util.h:171
void WriteVectorToConnectionViewer(const char *filename, const typename TFunction::algebra_type::vector_type &b, const TFunction &u, const typename TFunction::algebra_type::vector_type *pCompareVec=NULL)
Definition: grid_function_util.h:782
double & BlockRef(T &vec, size_t i)
Definition: blocks.h:66
number AverageFunctionDifference(SmartPtr< GridFunction< TDomain, TAlgebra > > spGridFct, std::string subset, std::string fct1, std::string fct2)
Calculates the average of the pointwise difference of two functions on given subset.
Definition: grid_function_util.h:981
void AdjustMeanValue(SmartPtr< TGridFunction > spGF, const std::vector< std::string > &vCmp, number mean)
Definition: grid_function_util.h:214
@ VOLUME
Definition: grid_base_objects.h:63
@ VERTEX
Definition: grid_base_objects.h:60
@ EDGE
Definition: grid_base_objects.h:61
@ FACE
Definition: grid_base_objects.h:62
void WriteMatrixToConnectionViewer(const char *filename, const typename TFunction::algebra_type::matrix_type &A, const TFunction &u)
Definition: grid_function_util.h:715
void SaveVectorCSV(TGridFunction &b, const char *filename)
Definition: grid_function_util.h:962
static void ScaleGFOnElems(ConstSmartPtr< DoFDistribution > dd, SmartPtr< TGridFunction > vecOut, ConstSmartPtr< TGridFunction > vecIn, const std::vector< number > &vScale)
Definition: grid_function_util.h:81
bool CheckGFValuesWithinBounds(ConstSmartPtr< TGridFunction > u, size_t cmp, number lowerBnd, number upperBnd)
Definition: grid_function_util.h:579
size_t GetSize(const T &t)
number Integral(SmartPtr< UserData< number, TGridFunction::dim > > spData, TGridFunction &spGridFct, const char *subsets, number time, int quadOrder, std::string quadType)
Definition: integrate.h:856
#define PROFILE_FUNC_GROUP(groups)
Definition: profiler.h:258
#define PROFILE_FUNC()
Definition: profiler.h:257
#define PROGRESS_START(progVarName, dSize, msg)
Definition: progress.h:111
#define PROGRESS_UPDATE(progVarName, d)
Definition: progress.h:117
T value_type
Definition: sparsematrix_interface.h:2