Chaste  Release::2018.1
HeartConfig.cpp
1 /*
2 
3 Copyright (c) 2005-2018, University of Oxford.
4 All rights reserved.
5 
6 University of Oxford means the Chancellor, Masters and Scholars of the
7 University of Oxford, having an administrative office at Wellington
8 Square, Oxford OX1 2JD, UK.
9 
10 This file is part of Chaste.
11 
12 Redistribution and use in source and binary forms, with or without
13 modification, are permitted provided that the following conditions are met:
14  * Redistributions of source code must retain the above copyright notice,
15  this list of conditions and the following disclaimer.
16  * Redistributions in binary form must reproduce the above copyright notice,
17  this list of conditions and the following disclaimer in the documentation
18  and/or other materials provided with the distribution.
19  * Neither the name of the University of Oxford nor the names of its
20  contributors may be used to endorse or promote products derived from this
21  software without specific prior written permission.
22 
23 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33 
34 */
35 
37 
38 #include "UblasCustomFunctions.hpp"
39 
40 #include "AbstractChasteRegion.hpp"
41 #include "ArchiveLocationInfo.hpp"
42 #include "ChastePoint.hpp"
43 #include "Exception.hpp"
44 #include "HeartConfig.hpp"
45 #include "HeartFileFinder.hpp"
46 #include "OutputFileHandler.hpp"
47 #include "Version.hpp"
48 #include "Warnings.hpp"
49 
50 #include "HeartRegionCodes.hpp"
51 
52 #include "RegularStimulus.hpp"
53 #include "SimpleStimulus.hpp"
54 
55 #include <cassert>
56 #include <fstream>
57 #include <istream>
58 #include <map>
59 #include <string>
60 
61 #include <xsd/cxx/tree/exceptions.hxx>
62 #include "XmlTools.hpp"
63 using namespace xsd::cxx::tree;
64 
65 // Coping with changes to XSD interface
66 #if (XSD_INT_VERSION >= 3000000L)
67 #define XSD_SEQUENCE_TYPE(base) base##_sequence
68 #define XSD_ITERATOR_TYPE(base) base##_iterator
69 #define XSD_NESTED_TYPE(t) t##_type
70 #define XSD_ANON_TYPE(t1, t2) \
71  t1::t2##_type
72 #else
73 #define XSD_SEQUENCE_TYPE(base) base::container
74 #define XSD_ITERATOR_TYPE(base) base::iterator
75 #define XSD_NESTED_TYPE(t) t::type
76 #define XSD_ANON_TYPE(t1, t2) \
77  t1::t2::_xsd_##t2##_::t2
78 #endif
79 
80 // These are for convenience
81 #define XSD_ANON_SEQUENCE_TYPE(t1, t2, t3) \
82  XSD_SEQUENCE_TYPE(XSD_ANON_TYPE(t1, t2)::t3)
83 #define XSD_ANON_ITERATOR_TYPE(t1, t2, t3) \
84  XSD_ITERATOR_TYPE(XSD_ANON_TYPE(t1, t2)::t3)
85 
86 // Newer versions don't allow you to set fixed attributes
87 #if (XSD_INT_VERSION >= 3020000L)
88 #define XSD_CREATE_WITH_FIXED_ATTR(type, name, attr) \
89  type name
90 #define XSD_CREATE_WITH_FIXED_ATTR1(type, name, arg1, attr) \
91  type name(arg1)
92 #define XSD_CREATE_WITH_FIXED_ATTR2(type, name, arg1, arg2, attr) \
93  type name(arg1, arg2)
94 #define XSD_CREATE_WITH_FIXED_ATTR3(type, name, arg1, arg2, arg3, attr) \
95  type name(arg1, arg2, arg3)
96 #else
97 #define XSD_CREATE_WITH_FIXED_ATTR(type, name, attr) \
98  type name(attr)
99 #define XSD_CREATE_WITH_FIXED_ATTR1(type, name, arg1, attr) \
100  type name(arg1, attr)
101 #define XSD_CREATE_WITH_FIXED_ATTR2(type, name, arg1, arg2, attr) \
102  type name(arg1, arg2, attr)
103 #define XSD_CREATE_WITH_FIXED_ATTR3(type, name, arg1, arg2, arg3, attr) \
104  type name(arg1, arg2, arg3, attr)
105 #endif
106 
110 #define ENSURE_SECTION_PRESENT(location, type) \
111  if (!location.present()) \
112  { \
113  type empty_item; \
114  location.set(empty_item); \
115  }
116 
117 #include <boost/current_function.hpp>
125 #define CHECK_EXISTS(test, path) \
126  do \
127  { \
128  if (!test) \
129  { \
130  EXCEPTION("No XML element " << path << " found in parameters when calling '" \
131  << BOOST_CURRENT_FUNCTION << "'"); \
132  } \
133  } while (false)
134 
140 {
141 public:
149  static void TransformIonicModelDefinitions(xercesc::DOMDocument* pDocument,
150  xercesc::DOMElement* pRootElement);
151 
160  static void TransformArchiveDirectory(xercesc::DOMDocument* pDocument,
161  xercesc::DOMElement* pRootElement);
162 
170  static void CheckForIluPreconditioner(xercesc::DOMDocument* pDocument,
171  xercesc::DOMElement* pRootElement);
172 
180  static void MoveConductivityHeterogeneities(xercesc::DOMDocument* pDocument,
181  xercesc::DOMElement* pRootElement);
182 
190  static void SetDefaultVisualizer(xercesc::DOMDocument* pDocument,
191  xercesc::DOMElement* pRootElement);
192 };
193 
194 //
195 // Default settings
196 //
197 #include "HeartConfigDefaults.hpp"
198 
199 //
200 // Definition of static member variables
201 //
202 boost::shared_ptr<HeartConfig> HeartConfig::mpInstance;
203 
204 //
205 // Methods
206 //
207 
209 {
210  if (mpInstance.get() == NULL)
211  {
212  mpInstance.reset(new HeartConfig);
213  }
214  return mpInstance.get();
215 }
216 
218  : mUseMassLumping(false),
219  mUseMassLumpingForPrecond(false),
220  mUseFixedNumberIterations(false),
221  mEvaluateNumItsEveryNSolves(UINT_MAX)
222 {
223  assert(mpInstance.get() == NULL);
226 
228  //CheckTimeSteps(); // necessity of this line of code is not tested -- remove with caution!
229 
230  //initialise the member variable of the layers
231  mEpiFraction = -1.0;
232  mEndoFraction = -1.0;
233  mMidFraction = -1.0;
235  // initialise to senseless values (these should be only 0, 1 and 2)
236  // note: the 'minus 3' is for checking purposes as we need to add 0, 1 or 2 to this initial value
237  // and UINT_MAX+1 seems to be 0
238  mIndexMid = UINT_MAX - 3u;
239  mIndexEpi = UINT_MAX - 3u;
240  mIndexEndo = UINT_MAX - 3u;
241 
243 
245  mTissueIdentifiers.insert(0);
246  mBathIdentifiers.insert(1);
247 }
248 
250 {
251 }
252 
253 void HeartConfig::Write(bool useArchiveLocationInfo, std::string subfolderName)
254 {
255  //Output file
256  std::string output_dirname;
257  if (useArchiveLocationInfo)
258  {
259  output_dirname = ArchiveLocationInfo::GetArchiveDirectory();
260  }
261  else
262  {
263  OutputFileHandler handler(GetOutputDirectory() + "/" + subfolderName, false);
264  output_dirname = handler.GetOutputDirectoryFullPath();
265  }
266 
267  // Sometimes this method is called collectively, sometimes not,
268  // in any case the below file operations only want to be performed by
269  // the master - so exit. Caller takes responsibility for nice
270  // exception handling.
271  if (!PetscTools::AmMaster())
272  {
273  return;
274  }
275 
276  out_stream p_parameters_file(new std::ofstream((output_dirname + "ChasteParameters.xml").c_str()));
277 
278  if (!p_parameters_file->is_open())
279  {
280  EXCEPTION("Could not open XML file in HeartConfig");
281  }
282 
283  //Schema map
284  //Note - this location is relative to where we are storing the xml
285  ::xml_schema::namespace_infomap map;
286  // Release 1.1 (and earlier) didn't use a namespace
287  map[""].schema = "ChasteParameters_1_1.xsd";
288  // Later releases use namespaces of the form https://chaste.comlab.ox.ac.uk/nss/parameters/N_M
289  map["cp20"].name = "https://chaste.comlab.ox.ac.uk/nss/parameters/2_0";
290  map["cp20"].schema = "ChasteParameters_2_0.xsd";
291  map["cp21"].name = "https://chaste.comlab.ox.ac.uk/nss/parameters/2_1";
292  map["cp21"].schema = "ChasteParameters_2_1.xsd";
293  map["cp22"].name = "https://chaste.comlab.ox.ac.uk/nss/parameters/2_2";
294  map["cp22"].schema = "ChasteParameters_2_2.xsd";
295  map["cp23"].name = "https://chaste.comlab.ox.ac.uk/nss/parameters/2_3";
296  map["cp23"].schema = "ChasteParameters_2_3.xsd";
297  map["cp30"].name = "https://chaste.comlab.ox.ac.uk/nss/parameters/3_0";
298  map["cp30"].schema = "ChasteParameters_3_0.xsd";
299  map["cp31"].name = "https://chaste.comlab.ox.ac.uk/nss/parameters/3_1";
300  map["cp31"].schema = "ChasteParameters_3_1.xsd";
301  map["cp33"].name = "https://chaste.comlab.ox.ac.uk/nss/parameters/3_3";
302  map["cp33"].schema = "ChasteParameters_3_3.xsd";
303  map["cp34"].name = "https://chaste.comlab.ox.ac.uk/nss/parameters/3_4";
304  map["cp34"].schema = "ChasteParameters_3_4.xsd";
305  // We use 'cp' as prefix for the latest version to avoid having to change saved
306  // versions for comparison at every release.
307  map["cp"].name = "https://chaste.comlab.ox.ac.uk/nss/parameters/2017_1";
308  map["cp"].schema = "ChasteParameters_2017_1.xsd";
309 
310  cp::ChasteParameters(*p_parameters_file, *mpParameters, map);
311 
312  // If we're archiving, try to save a copy of the latest schema too
313  if (useArchiveLocationInfo)
314  {
315  CopySchema(output_dirname);
316  }
317 }
318 
320 {
321  /*
322  * This method implements the logic required by HeartConfig to be able to handle resuming a simulation via the executable.
323  *
324  * When the control reaches the method mpParameters points to the file specified as resuming parameters.
325  * However SetParametersFile() will set this variable to point to the archived parameters.
326  *
327  * We make a temporary copy of mpParameters so we don't lose its content.
328  * At the end of the method we update the new mpParameters with the resuming parameters.
329  */
330  assert(mpParameters.use_count() > 0);
331  boost::shared_ptr<cp::chaste_parameters_type> p_new_parameters = mpParameters;
332 
333  /*
334  * When we unarchive a simulation, we load the old parameters file in order to inherit things such
335  * as default cell model, stimuli, heterogeneities, ... This has the side effect of inheriting the
336  * <CheckpointSimulation> element (if defined).
337  *
338  * We disable checkpointing definition coming from the unarchived config file. We will enable it again
339  * if defined in the resume config file.
340  */
341  std::string parameters_filename_xml = ArchiveLocationInfo::GetArchiveDirectory() + "ChasteParameters.xml";
342  mpParameters = ReadFile(parameters_filename_xml);
343  mParametersFilePath.SetPath(parameters_filename_xml, RelativeTo::AbsoluteOrCwd);
344 
345  // Release 3.0 and earlier wrote a separate defaults file in the checkpoint
346  std::string defaults_filename_xml = ArchiveLocationInfo::GetArchiveDirectory() + "ChasteDefaults.xml";
347  if (FileFinder(defaults_filename_xml).Exists())
348  {
349  boost::shared_ptr<cp::chaste_parameters_type> p_defaults = ReadFile(defaults_filename_xml);
350  MergeDefaults(mpParameters, p_defaults);
351  }
352 
354 
355  // If we are resuming a simulation, some parameters can be altered at this point.
356  if (p_new_parameters->ResumeSimulation().present())
357  {
358  UpdateParametersFromResumeSimulation(p_new_parameters);
359  }
360 
361  CheckTimeSteps(); // For consistency with SetParametersFile
362 }
363 
364 void HeartConfig::CopySchema(const std::string& rToDirectory)
365 {
366  // N.B. This method should only be called by the master process,
367  // in a situation where it can handle EXCEPTION()s nicely, e.g.
368  // TRY_IF_MASTER(CopySchema(...));
369 
370  std::string schema_name("ChasteParameters_2017_1.xsd");
371  FileFinder schema_location("heart/src/io/" + schema_name, RelativeTo::ChasteSourceRoot);
372  if (!schema_location.Exists())
373  {
374  // Try a relative path instead
375  schema_location.SetPath(schema_name, RelativeTo::CWD);
376  if (!schema_location.Exists())
377  {
378  // Warn the user
379  std::string message("Unable to locate schema file " + schema_name + ". You will need to ensure it is available when resuming from the checkpoint.");
380  WARN_ONCE_ONLY(message);
381  }
382  }
383  if (schema_location.Exists())
384  {
385  FileFinder output_directory(rToDirectory, RelativeTo::Absolute);
386  schema_location.CopyTo(output_directory);
387  }
388 }
389 
391 {
392  mSchemaLocations.clear();
393  // Location of schemas in the source tree
394  std::string root_dir = std::string(ChasteBuildInfo::GetRootDir()) + "/heart/src/io/";
395  // Release 1.1 (and earlier) didn't use a namespace
396  mSchemaLocations[""] = root_dir + "ChasteParameters_1_1.xsd";
397  // Later releases use namespaces of the form https://chaste.comlab.ox.ac.uk/nss/parameters/N_M
398  mSchemaLocations["https://chaste.comlab.ox.ac.uk/nss/parameters/2_0"] = root_dir + "ChasteParameters_2_0.xsd";
399  mSchemaLocations["https://chaste.comlab.ox.ac.uk/nss/parameters/2_1"] = root_dir + "ChasteParameters_2_1.xsd";
400  mSchemaLocations["https://chaste.comlab.ox.ac.uk/nss/parameters/2_2"] = root_dir + "ChasteParameters_2_2.xsd";
401  mSchemaLocations["https://chaste.comlab.ox.ac.uk/nss/parameters/2_3"] = root_dir + "ChasteParameters_2_3.xsd";
402  mSchemaLocations["https://chaste.comlab.ox.ac.uk/nss/parameters/3_0"] = root_dir + "ChasteParameters_3_0.xsd";
403  mSchemaLocations["https://chaste.comlab.ox.ac.uk/nss/parameters/3_1"] = root_dir + "ChasteParameters_3_1.xsd";
404  mSchemaLocations["https://chaste.comlab.ox.ac.uk/nss/parameters/3_3"] = root_dir + "ChasteParameters_3_3.xsd";
405  mSchemaLocations["https://chaste.comlab.ox.ac.uk/nss/parameters/3_4"] = root_dir + "ChasteParameters_3_4.xsd";
406  mSchemaLocations["https://chaste.comlab.ox.ac.uk/nss/parameters/2017_1"] = root_dir + "ChasteParameters_2017_1.xsd";
407 }
408 
409 unsigned HeartConfig::GetVersionFromNamespace(const std::string& rNamespaceUri)
410 {
411  unsigned version_major = 0;
412  unsigned version_minor = 0;
413  if (rNamespaceUri == "")
414  {
415  version_major = 1;
416  version_minor = 1;
417  }
418  else
419  {
420  std::string uri_base("https://chaste.comlab.ox.ac.uk/nss/parameters/");
421  if (rNamespaceUri.substr(0, uri_base.length()) == uri_base)
422  {
423  std::istringstream version_string(rNamespaceUri.substr(uri_base.length()));
424  version_string >> version_major;
425  version_string.ignore(1);
426  version_string >> version_minor;
427  if (version_string.fail())
428  {
429  version_major = 0;
430  version_minor = 0;
431  }
432  }
433  }
434 
435  unsigned version = version_major * 1000 + version_minor;
436  if (version == 0)
437  {
438  EXCEPTION(rNamespaceUri + " is not a recognised Chaste parameters namespace.");
439  }
440  return version;
441 }
442 
444 {
445  mSchemaLocations = rSchemaLocations;
447 }
448 
449 void HeartConfig::SetUseFixedSchemaLocation(bool useFixedSchemaLocation)
450 {
451  mUseFixedSchemaLocation = useFixedSchemaLocation;
452 }
453 
454 boost::shared_ptr<cp::chaste_parameters_type> HeartConfig::ReadFile(const std::string& rFileName)
455 {
456  // Determine whether to use the schema path given in the input XML, or our own schema
457  ::xml_schema::properties props;
459  {
460  for (SchemaLocationsMap::iterator it = mSchemaLocations.begin();
461  it != mSchemaLocations.end();
462  ++it)
463  {
464  if (it->first == "")
465  {
466  props.no_namespace_schema_location(XmlTools::EscapeSpaces(it->second));
467  }
468  else
469  {
470  props.schema_location(it->first, XmlTools::EscapeSpaces(it->second));
471  }
472  }
473  }
474 
475  // Get the parameters using the method 'ChasteParameters(rFileName)',
476  // which returns a std::auto_ptr. We convert to a shared_ptr for easier semantics.
477  try
478  {
479  // Make sure Xerces finalization happens
480  XmlTools::Finalizer finalizer(false);
481  // Parse XML to DOM
482  xsd::cxx::xml::dom::auto_ptr<xercesc::DOMDocument> p_doc = XmlTools::ReadXmlFile(rFileName, props);
483  // Test the namespace on the root element
484  xercesc::DOMElement* p_root_elt = p_doc->getDocumentElement();
485  std::string namespace_uri(X2C(p_root_elt->getNamespaceURI()));
486  const unsigned version = GetVersionFromNamespace(namespace_uri);
487  if (version < 2000) // Changes made in release 2.0
488  {
489  XmlTransforms::TransformIonicModelDefinitions(p_doc.get(), p_root_elt);
490  }
491  if (version < 2001) // Changes made in release 2.1
492  {
493  XmlTransforms::TransformArchiveDirectory(p_doc.get(), p_root_elt);
494  XmlTransforms::CheckForIluPreconditioner(p_doc.get(), p_root_elt);
495  }
496  if (version < 3001) // Changes made in release 3.1
497  {
498  XmlTransforms::MoveConductivityHeterogeneities(p_doc.get(), p_root_elt);
499  }
500  if (version < 3003) // Changes made in release 3.3
501  {
502  XmlTransforms::SetDefaultVisualizer(p_doc.get(), p_root_elt);
503  }
504  if (version < 3004) // Not the latest in release 3.4
505  {
506  XmlTools::SetNamespace(p_doc.get(), p_root_elt, "https://chaste.comlab.ox.ac.uk/nss/parameters/3_4");
507  }
508  if (version < 2017001) // Not the latest release
509  {
510  XmlTools::SetNamespace(p_doc.get(), p_root_elt, "https://chaste.comlab.ox.ac.uk/nss/parameters/2017_1");
511  }
512  // Parse DOM to object model
513  boost::shared_ptr<cp::chaste_parameters_type> p_params(cp::ChasteParameters(*p_doc, ::xml_schema::flags::dont_initialize, props));
514  // Get rid of the DOM stuff
515  p_doc.reset();
516 
517  return boost::shared_ptr<cp::chaste_parameters_type>(p_params);
518  }
519  catch (const xml_schema::exception& e)
520  {
521  std::cerr << e << std::endl;
522  // Make sure we don't store invalid parameters
523  mpParameters.reset();
524  EXCEPTION("XML parsing error in configuration file: " + rFileName);
525  }
526  catch (...)
527  {
528  // Make sure we don't store invalid parameters
529  mpParameters.reset();
530  throw;
531  }
532 }
533 
534 void HeartConfig::SetParametersFile(const std::string& rFileName)
535 {
536  mpParameters = ReadFile(rFileName);
539 
540  if (IsSimulationDefined())
541  {
542  CheckTimeSteps(); // Resume files might not have time steps defined
543  }
544 }
545 
547 {
548  return mParametersFilePath;
549 }
550 
551 void HeartConfig::UpdateParametersFromResumeSimulation(boost::shared_ptr<cp::chaste_parameters_type> pResumeParameters)
552 {
553  // Check for user foolishness
554  if ((pResumeParameters->ResumeSimulation()->SpaceDimension() != HeartConfig::Instance()->GetSpaceDimension())
555  || (pResumeParameters->ResumeSimulation()->Domain() != HeartConfig::Instance()->GetDomain()))
556  {
557  EXCEPTION("Problem type and space dimension should match when restarting a simulation.");
558  }
559 
560  // New simulation duration
561  HeartConfig::Instance()->SetSimulationDuration(pResumeParameters->ResumeSimulation()->SimulationDuration());
562 
563  // Stimulus definition. For these we always replace any previous definitions (at least for now...)
564  if (pResumeParameters->ResumeSimulation()->Stimuli().present())
565  {
566  mpParameters->Simulation()->Stimuli().set(pResumeParameters->ResumeSimulation()->Stimuli().get());
567  }
568 
569  // Cell heterogeneities. Note that while we copy the elements here, other code in CardiacSimulation actually updates
570  // the loaded simulation to take account of the new settings.
571  if (pResumeParameters->ResumeSimulation()->CellHeterogeneities().present())
572  {
573  if (!mpParameters->Simulation()->CellHeterogeneities().present())
574  {
575  // Original parameters had no heterogeneities, so just copy the whole element
576  mpParameters->Simulation()->CellHeterogeneities().set(pResumeParameters->ResumeSimulation()->CellHeterogeneities().get());
577  }
578  else
579  {
580  // Need to append the new heterogeneity defitions to the original sequence
581  XSD_SEQUENCE_TYPE(cp::cell_heterogeneities_type::CellHeterogeneity)& new_seq = pResumeParameters->ResumeSimulation()->CellHeterogeneities()->CellHeterogeneity();
582  XSD_SEQUENCE_TYPE(cp::cell_heterogeneities_type::CellHeterogeneity)& orig_seq = mpParameters->Simulation()->CellHeterogeneities()->CellHeterogeneity();
583  for (XSD_ITERATOR_TYPE(cp::cell_heterogeneities_type::CellHeterogeneity) i = new_seq.begin();
584  i != new_seq.end();
585  ++i)
586  {
587  orig_seq.push_back(*i);
588  }
589  }
590  }
591 
592  // Whether to checkpoint the resumed simulation
593  if (pResumeParameters->ResumeSimulation()->CheckpointSimulation().present())
594  {
596  pResumeParameters->ResumeSimulation()->CheckpointSimulation()->timestep(),
597  pResumeParameters->ResumeSimulation()->CheckpointSimulation()->max_checkpoints_on_disk());
598  }
599 
600  //Visualization parameters are no longer compulsory
601  if (pResumeParameters->ResumeSimulation()->OutputVisualizer().present())
602  {
603  HeartConfig::Instance()->SetVisualizeWithParallelVtk(pResumeParameters->ResumeSimulation()->OutputVisualizer()->parallel_vtk() == cp::yesno_type::yes);
604  HeartConfig::Instance()->SetVisualizeWithVtk(pResumeParameters->ResumeSimulation()->OutputVisualizer()->vtk() == cp::yesno_type::yes);
605  HeartConfig::Instance()->SetVisualizeWithCmgui(pResumeParameters->ResumeSimulation()->OutputVisualizer()->cmgui() == cp::yesno_type::yes);
606  HeartConfig::Instance()->SetVisualizeWithMeshalyzer(pResumeParameters->ResumeSimulation()->OutputVisualizer()->meshalyzer() == cp::yesno_type::yes);
607  }
608 
609  // Numerical parameters may be overridden
610  {
611  cp::numerical_type& r_resume = pResumeParameters->Numerical();
612  cp::numerical_type& r_user = mpParameters->Numerical();
613  if (r_resume.TimeSteps().present())
614  {
615  r_user.TimeSteps().set(r_resume.TimeSteps().get());
616  }
617  if (r_resume.KSPTolerances().present())
618  {
619  r_user.KSPTolerances().set(r_resume.KSPTolerances().get());
620  }
621  if (r_resume.KSPSolver().present())
622  {
623  r_user.KSPSolver().set(r_resume.KSPSolver().get());
624  }
625  if (r_resume.KSPPreconditioner().present())
626  {
627  r_user.KSPPreconditioner().set(r_resume.KSPPreconditioner().get());
628  }
629  if (r_resume.AdaptivityParameters().present())
630  {
631  r_user.AdaptivityParameters().set(r_resume.AdaptivityParameters().get());
632  }
633  }
634 
635  // Post-processing parameters may be overridden
636  if (pResumeParameters->PostProcessing().present())
637  {
639  cp::postprocessing_type& r_resume = pResumeParameters->PostProcessing().get();
640  cp::postprocessing_type& r_user = mpParameters->PostProcessing().get();
641  if (!r_resume.ActionPotentialDurationMap().empty())
642  {
643  r_user.ActionPotentialDurationMap() = r_resume.ActionPotentialDurationMap();
644  }
645  if (!r_resume.UpstrokeTimeMap().empty())
646  {
647  r_user.UpstrokeTimeMap() = r_resume.UpstrokeTimeMap();
648  }
649  if (!r_resume.MaxUpstrokeVelocityMap().empty())
650  {
651  r_user.MaxUpstrokeVelocityMap() = r_resume.MaxUpstrokeVelocityMap();
652  }
653  if (!r_resume.ConductionVelocityMap().empty())
654  {
655  r_user.ConductionVelocityMap() = r_resume.ConductionVelocityMap();
656  }
657  if (!r_resume.PseudoEcgElectrodePosition().empty())
658  {
659  r_user.PseudoEcgElectrodePosition() = r_resume.PseudoEcgElectrodePosition();
660  }
661  }
662 }
663 
665 {
666  // Throw it away first, so that mpInstance is NULL when we...
667  mpInstance.reset();
668  // ...make a new one
669  mpInstance.reset(new HeartConfig);
670 }
671 
673 {
674  return mpParameters->Simulation().present();
675 }
676 
678 {
679  return mpParameters->ResumeSimulation().present();
680 }
681 
682 void HeartConfig::CheckSimulationIsDefined(std::string callingMethod) const
683 {
684  if (IsSimulationResumed())
685  {
686  EXCEPTION(callingMethod + " information is not available in a resumed simulation.");
687  }
688 }
689 
690 void HeartConfig::CheckResumeSimulationIsDefined(std::string callingMethod) const
691 {
692  if (IsSimulationDefined())
693  {
694  EXCEPTION(callingMethod + " information is not available in a standard (non-resumed) simulation.");
695  }
696 }
697 
699 {
700  if (IsSimulationDefined())
701  {
702  CHECK_EXISTS(mpParameters->Simulation()->SpaceDimension().present(), "Simulation/SpaceDimension");
703  return mpParameters->Simulation()->SpaceDimension().get();
704  }
705  else
706  {
707  return mpParameters->ResumeSimulation()->SpaceDimension();
708  }
709 }
710 
712 {
713  if (IsSimulationDefined())
714  {
715  CHECK_EXISTS(mpParameters->Simulation()->SimulationDuration().present(), "Simulation/SimulationDuration");
716  return mpParameters->Simulation()->SimulationDuration().get();
717  }
718  else // IsSimulationResumed
719  {
720  return mpParameters->ResumeSimulation()->SimulationDuration();
721  }
722 }
723 
724 cp::domain_type HeartConfig::GetDomain() const
725 {
726  if (IsSimulationDefined())
727  {
728  CHECK_EXISTS(mpParameters->Simulation()->Domain().present(), "Simulation/Domain");
729  return mpParameters->Simulation()->Domain().get();
730  }
731  else
732  {
733  return mpParameters->ResumeSimulation()->Domain();
734  }
735 }
736 
737 cp::ionic_model_selection_type HeartConfig::GetDefaultIonicModel() const
738 {
739  CheckSimulationIsDefined("DefaultIonicModel");
740 
741  return mpParameters->Simulation()->IonicModels()->Default();
742 }
743 
744 template <unsigned DIM>
745 void HeartConfig::GetIonicModelRegions(std::vector<boost::shared_ptr<AbstractChasteRegion<DIM> > >& definedRegions,
746  std::vector<cp::ionic_model_selection_type>& ionicModels) const
747 {
748  CheckSimulationIsDefined("IonicModelRegions");
749  definedRegions.clear();
750  ionicModels.clear();
751 
752  XSD_SEQUENCE_TYPE(cp::ionic_models_type::Region)& regions = mpParameters->Simulation()->IonicModels()->Region();
753 
754  for (XSD_ITERATOR_TYPE(cp::ionic_models_type::Region) i = regions.begin();
755  i != regions.end();
756  ++i)
757  {
758  cp::ionic_model_region_type ionic_model_region(*i);
759 
760  if (ionic_model_region.Location().Cuboid().present() || ionic_model_region.Location().Ellipsoid().present())
761  {
762  if (ionic_model_region.Location().Cuboid().present())
763  {
764  cp::point_type point_a = ionic_model_region.Location().Cuboid()->LowerCoordinates();
765  cp::point_type point_b = ionic_model_region.Location().Cuboid()->UpperCoordinates();
766 
767  switch (DIM)
768  {
769  case 1:
770  {
771  ChastePoint<DIM> chaste_point_a(point_a.x());
772  ChastePoint<DIM> chaste_point_b(point_b.x());
773  definedRegions.push_back(boost::shared_ptr<AbstractChasteRegion<DIM> >(new ChasteCuboid<DIM>(chaste_point_a, chaste_point_b)));
774  break;
775  }
776  case 2:
777  {
778  ChastePoint<DIM> chaste_point_a(point_a.x(), point_a.y());
779  ChastePoint<DIM> chaste_point_b(point_b.x(), point_b.y());
780  definedRegions.push_back(boost::shared_ptr<AbstractChasteRegion<DIM> >(new ChasteCuboid<DIM>(chaste_point_a, chaste_point_b)));
781  break;
782  }
783  case 3:
784  {
785  ChastePoint<DIM> chaste_point_a(point_a.x(), point_a.y(), point_a.z());
786  ChastePoint<DIM> chaste_point_b(point_b.x(), point_b.y(), point_b.z());
787  definedRegions.push_back(boost::shared_ptr<AbstractChasteRegion<DIM> >(new ChasteCuboid<DIM>(chaste_point_a, chaste_point_b)));
788  break;
789  }
790  default:
792  break;
793  }
794  }
795  else if (ionic_model_region.Location().Ellipsoid().present())
796  {
797  cp::point_type centre = ionic_model_region.Location().Ellipsoid()->Centre();
798  cp::point_type radii = ionic_model_region.Location().Ellipsoid()->Radii();
799  switch (DIM)
800  {
801  case 1:
802  {
803  ChastePoint<DIM> chaste_point_a(centre.x());
804  ChastePoint<DIM> chaste_point_b(radii.x());
805  definedRegions.push_back(boost::shared_ptr<AbstractChasteRegion<DIM> >(new ChasteEllipsoid<DIM>(chaste_point_a, chaste_point_b)));
806  break;
807  }
808  case 2:
809  {
810  ChastePoint<DIM> chaste_point_a(centre.x(), centre.y());
811  ChastePoint<DIM> chaste_point_b(radii.x(), radii.y());
812  definedRegions.push_back(boost::shared_ptr<AbstractChasteRegion<DIM> >(new ChasteEllipsoid<DIM>(chaste_point_a, chaste_point_b)));
813  break;
814  }
815  case 3:
816  {
817  ChastePoint<DIM> chaste_point_a(centre.x(), centre.y(), centre.z());
818  ChastePoint<DIM> chaste_point_b(radii.x(), radii.y(), radii.z());
819  definedRegions.push_back(boost::shared_ptr<AbstractChasteRegion<DIM> >(new ChasteEllipsoid<DIM>(chaste_point_a, chaste_point_b)));
820  break;
821  }
822  default:
823  {
825  break;
826  }
827  }
828  }
829  else
830  {
832  }
833 
834  ionicModels.push_back(ionic_model_region.IonicModel());
835  }
836  else if (ionic_model_region.Location().EpiLayer().present() || ionic_model_region.Location().MidLayer().present() || ionic_model_region.Location().EndoLayer().present())
837  {
839  EXCEPTION("Definition of transmural layers is not yet supported for defining different ionic models, please use cuboids instead");
840  }
841  else
842  {
843  EXCEPTION("Invalid region type for ionic model definition");
844  }
845  }
846 }
847 
849 {
850  CheckSimulationIsDefined("Mesh");
851  return mpParameters->Simulation()->Mesh().present();
852 }
853 
855 {
856  CheckSimulationIsDefined("Mesh");
857  CHECK_EXISTS(IsMeshProvided(), "Simulation/Mesh");
858  cp::mesh_type mesh = mpParameters->Simulation()->Mesh().get();
859  return (mesh.Slab().present() || mesh.Sheet().present() || mesh.Fibre().present());
860 }
861 
863 {
864  CheckSimulationIsDefined("Mesh");
865  CHECK_EXISTS(IsMeshProvided(), "Simulation/Mesh");
866  cp::mesh_type mesh = mpParameters->Simulation()->Mesh().get();
867  return (mesh.Slab().present());
868 }
869 
871 {
872  CheckSimulationIsDefined("Mesh");
873  CHECK_EXISTS(IsMeshProvided(), "Simulation/Mesh");
874  cp::mesh_type mesh = mpParameters->Simulation()->Mesh().get();
875  return (mesh.Sheet().present());
876 }
877 
879 {
880  CheckSimulationIsDefined("Mesh");
881  CHECK_EXISTS(IsMeshProvided(), "Simulation/Mesh");
882  cp::mesh_type mesh = mpParameters->Simulation()->Mesh().get();
883  return (mesh.Fibre().present());
884 }
885 
887 {
888  CheckSimulationIsDefined("Mesh");
889  CHECK_EXISTS(IsMeshProvided(), "Simulation/Mesh");
890  return (mpParameters->Simulation()->Mesh()->LoadMesh().present());
891 }
892 
893 void HeartConfig::GetSlabDimensions(c_vector<double, 3>& slabDimensions) const
894 {
895  CheckSimulationIsDefined("Slab");
896 
897  if (GetSpaceDimension() != 3 || !GetCreateSlab())
898  {
899  EXCEPTION("Tissue slabs can only be defined in 3D");
900  }
901 
902  optional<cp::slab_type, false> slab_dimensions = mpParameters->Simulation()->Mesh()->Slab();
903 
904  slabDimensions[0] = slab_dimensions->x();
905  slabDimensions[1] = slab_dimensions->y();
906  slabDimensions[2] = slab_dimensions->z();
907 }
908 
909 void HeartConfig::GetSheetDimensions(c_vector<double, 2>& sheetDimensions) const
910 {
911  CheckSimulationIsDefined("Sheet");
912 
913  if (GetSpaceDimension() != 2 || !GetCreateSheet())
914  {
915  EXCEPTION("Tissue sheets can only be defined in 2D");
916  }
917 
918  optional<cp::sheet_type, false> sheet_dimensions = mpParameters->Simulation()->Mesh()->Sheet();
919 
920  sheetDimensions[0] = sheet_dimensions->x();
921  sheetDimensions[1] = sheet_dimensions->y();
922 }
923 
924 void HeartConfig::GetFibreLength(c_vector<double, 1>& fibreLength) const
925 {
926  CheckSimulationIsDefined("Fibre");
927 
928  if (GetSpaceDimension() != 1 || !GetCreateFibre())
929  {
930  EXCEPTION("Tissue fibres can only be defined in 1D");
931  }
932 
933  optional<cp::fibre_type, false> fibre_length = mpParameters->Simulation()->Mesh()->Fibre();
934 
935  fibreLength[0] = fibre_length->x();
936 }
937 
939 {
940  CheckSimulationIsDefined("InterNodeSpace");
941 
942  switch (GetSpaceDimension())
943  {
944  case 3:
945  CHECK_EXISTS(GetCreateSlab(), "Simulation/Mesh/Slab");
946  return mpParameters->Simulation()->Mesh()->Slab()->inter_node_space();
947  break;
948  case 2:
949  CHECK_EXISTS(GetCreateSheet(), "Simulation/Mesh/Sheet");
950  return mpParameters->Simulation()->Mesh()->Sheet()->inter_node_space();
951  break;
952  case 1:
953  CHECK_EXISTS(GetCreateFibre(), "Simulation/Mesh/Fibre");
954  return mpParameters->Simulation()->Mesh()->Fibre()->inter_node_space();
955  break;
956  default:
958  // LCOV_EXCL_START
959  return 0.0; //To fool the compiler
960  // LCOV_EXCL_STOP
961  }
962 }
963 
964 std::string HeartConfig::GetMeshName() const
965 {
966  CheckSimulationIsDefined("LoadMesh");
967  CHECK_EXISTS(GetLoadMesh(), "Mesh/LoadMesh");
968 
969  return mpParameters->Simulation()->Mesh()->LoadMesh()->name();
970 }
971 
972 cp::media_type HeartConfig::GetConductivityMedia() const
973 {
974  CheckSimulationIsDefined("LoadMesh");
975  CHECK_EXISTS(GetLoadMesh(), "Mesh/LoadMesh");
976 
977  return mpParameters->Simulation()->Mesh()->LoadMesh()->conductivity_media();
978 }
979 
980 template <unsigned DIM>
981 void HeartConfig::GetStimuli(std::vector<boost::shared_ptr<AbstractStimulusFunction> >& rStimuliApplied,
982  std::vector<boost::shared_ptr<AbstractChasteRegion<DIM> > >& rStimulatedAreas) const
983 {
984  CheckSimulationIsDefined("Stimuli");
985 
986  if (!mpParameters->Simulation()->Stimuli().present())
987  {
988  // Finding no stimuli defined is allowed (although HeartConfigRelatedFactory does
989  // throw an exception is no stimuli and no electrodes)
990  return;
991  }
992 
993  XSD_SEQUENCE_TYPE(cp::stimuli_type::Stimulus)
994  stimuli = mpParameters->Simulation()->Stimuli()->Stimulus();
995 
996  for (XSD_ITERATOR_TYPE(cp::stimuli_type::Stimulus) i = stimuli.begin();
997  i != stimuli.end();
998  ++i)
999  {
1000  cp::stimulus_type stimulus(*i);
1001  if (stimulus.Location().Cuboid().present() || stimulus.Location().Ellipsoid().present())
1002  {
1003  boost::shared_ptr<AbstractChasteRegion<DIM> > area_ptr;
1004  if (stimulus.Location().Cuboid().present())
1005  {
1006  cp::point_type point_a = stimulus.Location().Cuboid()->LowerCoordinates();
1007  cp::point_type point_b = stimulus.Location().Cuboid()->UpperCoordinates();
1008  switch (DIM)
1009  {
1010  case 1:
1011  {
1012  ChastePoint<DIM> chaste_point_a(point_a.x());
1013  ChastePoint<DIM> chaste_point_b(point_b.x());
1014  area_ptr.reset(new ChasteCuboid<DIM>(chaste_point_a, chaste_point_b));
1015  break;
1016  }
1017  case 2:
1018  {
1019  ChastePoint<DIM> chaste_point_a(point_a.x(), point_a.y());
1020  ChastePoint<DIM> chaste_point_b(point_b.x(), point_b.y());
1021  area_ptr.reset(new ChasteCuboid<DIM>(chaste_point_a, chaste_point_b));
1022  break;
1023  }
1024  case 3:
1025  {
1026  ChastePoint<DIM> chaste_point_a(point_a.x(), point_a.y(), point_a.z());
1027  ChastePoint<DIM> chaste_point_b(point_b.x(), point_b.y(), point_b.z());
1028  area_ptr.reset(new ChasteCuboid<DIM>(chaste_point_a, chaste_point_b));
1029  break;
1030  }
1031  default:
1032  NEVER_REACHED;
1033  break;
1034  }
1035  }
1036  else if (stimulus.Location().Ellipsoid().present())
1037  {
1038  cp::point_type centre = stimulus.Location().Ellipsoid()->Centre();
1039  cp::point_type radii = stimulus.Location().Ellipsoid()->Radii();
1040  switch (DIM)
1041  {
1042  case 1:
1043  {
1044  ChastePoint<DIM> chaste_point_a(centre.x());
1045  ChastePoint<DIM> chaste_point_b(radii.x());
1046  area_ptr.reset(new ChasteEllipsoid<DIM>(chaste_point_a, chaste_point_b));
1047  break;
1048  }
1049  case 2:
1050  {
1051  ChastePoint<DIM> chaste_point_a(centre.x(), centre.y());
1052  ChastePoint<DIM> chaste_point_b(radii.x(), radii.y());
1053  area_ptr.reset(new ChasteEllipsoid<DIM>(chaste_point_a, chaste_point_b));
1054  break;
1055  }
1056  case 3:
1057  {
1058  ChastePoint<DIM> chaste_point_a(centre.x(), centre.y(), centre.z());
1059  ChastePoint<DIM> chaste_point_b(radii.x(), radii.y(), radii.z());
1060  area_ptr.reset(new ChasteEllipsoid<DIM>(chaste_point_a, chaste_point_b));
1061  break;
1062  }
1063  default:
1064  {
1065  NEVER_REACHED;
1066  break;
1067  }
1068  }
1069  }
1070  rStimulatedAreas.push_back(area_ptr);
1071 
1072  boost::shared_ptr<AbstractStimulusFunction> stim;
1073 
1074  if (stimulus.Period().present())
1075  {
1076  if (stimulus.StopTime().present())
1077  {
1078  stim.reset(new RegularStimulus(stimulus.Strength(),
1079  stimulus.Duration(),
1080  stimulus.Period().get(),
1081  stimulus.Delay(),
1082  stimulus.StopTime().get()));
1083  }
1084  else
1085  {
1086  stim.reset(new RegularStimulus(stimulus.Strength(),
1087  stimulus.Duration(),
1088  stimulus.Period().get(),
1089  stimulus.Delay()));
1090  }
1091  }
1092  else
1093  {
1094  if (stimulus.StopTime().present())
1095  {
1096  EXCEPTION("Stop time can not be defined for SimpleStimulus. Use Duration instead.");
1097  }
1098 
1099  stim.reset(new SimpleStimulus(stimulus.Strength(),
1100  stimulus.Duration(),
1101  stimulus.Delay()));
1102  }
1103  rStimuliApplied.push_back(stim);
1104  }
1105  else if (stimulus.Location().EpiLayer().present() || stimulus.Location().MidLayer().present() || stimulus.Location().EndoLayer().present())
1106  {
1107  EXCEPTION("Definition of transmural layers is not yet supported for specifying stimulated areas, please use cuboids instead");
1108  }
1109  else
1110  {
1111  EXCEPTION("Invalid region type for stimulus definition");
1112  }
1113  }
1114 }
1115 
1116 template <unsigned DIM>
1117 void HeartConfig::GetCellHeterogeneities(std::vector<boost::shared_ptr<AbstractChasteRegion<DIM> > >& rCellHeterogeneityRegions,
1118  std::vector<double>& rScaleFactorGks,
1119  std::vector<double>& rScaleFactorIto,
1120  std::vector<double>& rScaleFactorGkr,
1121  std::vector<std::map<std::string, double> >* pParameterSettings)
1122 {
1123  CheckSimulationIsDefined("CellHeterogeneities");
1124 
1125  if (!mpParameters->Simulation()->CellHeterogeneities().present())
1126  {
1127  // finding no heterogeneities defined is allowed
1128  return;
1129  }
1130  XSD_SEQUENCE_TYPE(cp::cell_heterogeneities_type::CellHeterogeneity)
1131  cell_heterogeneity
1132  = mpParameters->Simulation()->CellHeterogeneities()->CellHeterogeneity();
1133 
1134  bool user_supplied_negative_value = false;
1135  mUserAskedForCellularTransmuralHeterogeneities = false; // overwritten with true below if necessary
1136  bool user_asked_for_cuboids_or_ellipsoids = false;
1137  unsigned counter_of_heterogeneities = 0;
1138 
1139  for (XSD_ITERATOR_TYPE(cp::cell_heterogeneities_type::CellHeterogeneity) i = cell_heterogeneity.begin();
1140  i != cell_heterogeneity.end();
1141  ++i)
1142  {
1143  cp::cell_heterogeneity_type ht(*i);
1144 
1145  if (ht.Location().Cuboid().present())
1146  {
1147  user_asked_for_cuboids_or_ellipsoids = true;
1148  cp::point_type point_a = ht.Location().Cuboid()->LowerCoordinates();
1149  cp::point_type point_b = ht.Location().Cuboid()->UpperCoordinates();
1150 
1151  ChastePoint<DIM> chaste_point_a(point_a.x(), point_a.y(), point_a.z());
1152  ChastePoint<DIM> chaste_point_b(point_b.x(), point_b.y(), point_b.z());
1153 
1154  rCellHeterogeneityRegions.push_back(boost::shared_ptr<AbstractChasteRegion<DIM> >(new ChasteCuboid<DIM>(chaste_point_a, chaste_point_b)));
1155  }
1156  else if (ht.Location().Ellipsoid().present())
1157  {
1158  user_asked_for_cuboids_or_ellipsoids = true;
1159  cp::point_type centre = ht.Location().Ellipsoid()->Centre();
1160  cp::point_type radii = ht.Location().Ellipsoid()->Radii();
1161 
1162  ChastePoint<DIM> chaste_point_a(centre.x(), centre.y(), centre.z());
1163  ChastePoint<DIM> chaste_point_b(radii.x(), radii.y(), radii.z());
1164  rCellHeterogeneityRegions.push_back(boost::shared_ptr<AbstractChasteRegion<DIM> >(new ChasteEllipsoid<DIM>(chaste_point_a, chaste_point_b)));
1165  }
1166  else if (ht.Location().EpiLayer().present())
1167  {
1168  mEpiFraction = ht.Location().EpiLayer().get();
1169 
1171  if (mEpiFraction < 0)
1172  {
1173  user_supplied_negative_value = true;
1174  }
1175  mIndexEpi = counter_of_heterogeneities;
1176  }
1177  else if (ht.Location().EndoLayer().present())
1178  {
1179  mEndoFraction = ht.Location().EndoLayer().get();
1180 
1182  if (mEndoFraction < 0)
1183  {
1184  user_supplied_negative_value = true;
1185  }
1186  mIndexEndo = counter_of_heterogeneities;
1187  }
1188  else if (ht.Location().MidLayer().present())
1189  {
1190  mMidFraction = ht.Location().MidLayer().get();
1191 
1193  if (mMidFraction < 0)
1194  {
1195  user_supplied_negative_value = true;
1196  }
1197  mIndexMid = counter_of_heterogeneities;
1198  }
1199  else
1200  {
1201  EXCEPTION("Invalid region type for cell heterogeneity definition");
1202  }
1203 
1204  // Old scale factors
1205  rScaleFactorGks.push_back(ht.ScaleFactorGks().present() ? (double)ht.ScaleFactorGks().get() : 1.0);
1206  rScaleFactorIto.push_back(ht.ScaleFactorIto().present() ? (double)ht.ScaleFactorIto().get() : 1.0);
1207  rScaleFactorGkr.push_back(ht.ScaleFactorGkr().present() ? (double)ht.ScaleFactorGkr().get() : 1.0);
1208 
1209  // Named parameters
1210  if (pParameterSettings)
1211  {
1212  std::map<std::string, double> param_settings;
1213  XSD_SEQUENCE_TYPE(cp::cell_heterogeneity_type::SetParameter)& params = ht.SetParameter();
1214  for (XSD_ITERATOR_TYPE(cp::cell_heterogeneity_type::SetParameter) param_it = params.begin();
1215  param_it != params.end();
1216  ++param_it)
1217  {
1218  cp::set_parameter_type param(*param_it);
1219  param_settings[param.name()] = param.value();
1220  }
1221  pParameterSettings->push_back(param_settings);
1222  }
1223 
1224  counter_of_heterogeneities++;
1225  }
1226 
1228  {
1229  // cuboids/ellipsoids and layers at the same time are not yet supported
1230  if (user_asked_for_cuboids_or_ellipsoids)
1231  {
1232  EXCEPTION("Specification of cellular heterogeneities by cuboids/ellipsoids and layers at the same time is not yet supported");
1233  }
1234 
1235  //check that the user supplied all three layers, the indexes should be 0, 1 and 2.
1236  // As they are initialised to a higher value, if their summation is higher than 3,
1237  // one (or more) is missing
1238  if ((mIndexMid + mIndexEndo + mIndexEpi) > 3)
1239  {
1240  EXCEPTION("Three specifications of layers must be supplied");
1241  }
1242  if (fabs((mEndoFraction + mMidFraction + mEpiFraction) - 1) > 1e-2)
1243  {
1244  EXCEPTION("Summation of epicardial, midmyocardial and endocardial fractions should be 1");
1245  }
1246  if (user_supplied_negative_value)
1247  {
1248  EXCEPTION("Fractions must be positive");
1249  }
1250  }
1251 }
1252 
1254 {
1256 }
1257 
1259 {
1260  return mEpiFraction;
1261 }
1262 
1264 {
1265  return mEndoFraction;
1266 }
1267 
1269 {
1270  return mMidFraction;
1271 }
1272 
1274 {
1275  return mIndexEpi;
1276 }
1277 
1279 {
1280  return mIndexEndo;
1281 }
1282 
1284 {
1285  return mIndexMid;
1286 }
1287 
1289 {
1290  CheckSimulationIsDefined("ConductivityHeterogeneities");
1291  return mpParameters->Physiological().ConductivityHeterogeneities().present();
1292 }
1293 
1294 template <unsigned DIM>
1296  std::vector<boost::shared_ptr<AbstractChasteRegion<DIM> > >& rConductivitiesHeterogeneityAreas,
1297  std::vector<c_vector<double, 3> >& rIntraConductivities,
1298  std::vector<c_vector<double, 3> >& rExtraConductivities) const
1299 {
1300  CheckSimulationIsDefined("ConductivityHeterogeneities");
1301  CHECK_EXISTS(GetConductivityHeterogeneitiesProvided(), "Physiological/ConductivityHeterogeneities");
1302  XSD_ANON_SEQUENCE_TYPE(cp::physiological_type, ConductivityHeterogeneities, ConductivityHeterogeneity)& conductivity_heterogeneity = mpParameters->Physiological().ConductivityHeterogeneities()->ConductivityHeterogeneity();
1303 
1304  for (XSD_ANON_ITERATOR_TYPE(cp::physiological_type, ConductivityHeterogeneities, ConductivityHeterogeneity) i = conductivity_heterogeneity.begin();
1305  i != conductivity_heterogeneity.end();
1306  ++i)
1307  {
1308  cp::conductivity_heterogeneity_type ht(*i);
1309 
1310  if (ht.Location().Cuboid().present())
1311  {
1312  cp::point_type point_a = ht.Location().Cuboid()->LowerCoordinates();
1313  cp::point_type point_b = ht.Location().Cuboid()->UpperCoordinates();
1314  ChastePoint<DIM> chaste_point_a(point_a.x(), point_a.y(), point_a.z());
1315  ChastePoint<DIM> chaste_point_b(point_b.x(), point_b.y(), point_b.z());
1316  rConductivitiesHeterogeneityAreas.push_back(boost::shared_ptr<AbstractChasteRegion<DIM> >(new ChasteCuboid<DIM>(chaste_point_a, chaste_point_b)));
1317  }
1318  else if (ht.Location().Ellipsoid().present())
1319  {
1320  cp::point_type centre = ht.Location().Ellipsoid()->Centre();
1321  cp::point_type radii = ht.Location().Ellipsoid()->Radii();
1322  ChastePoint<DIM> chaste_point_a(centre.x(), centre.y(), centre.z());
1323  ChastePoint<DIM> chaste_point_b(radii.x(), radii.y(), radii.z());
1324  rConductivitiesHeterogeneityAreas.push_back(boost::shared_ptr<AbstractChasteRegion<DIM> >(new ChasteEllipsoid<DIM>(chaste_point_a, chaste_point_b)));
1325  }
1326  else if (ht.Location().EpiLayer().present() || ht.Location().MidLayer().present() || ht.Location().EndoLayer().present())
1327  {
1329  EXCEPTION("Definition of transmural layers is not allowed for conductivities heterogeneities, you may use fibre orientation support instead");
1330  }
1331  else
1332  {
1333  EXCEPTION("Invalid region type for conductivity definition");
1334  }
1335 
1336  if (ht.IntracellularConductivities().present())
1337  {
1338  double intra_x = ht.IntracellularConductivities()->longi();
1339  double intra_y = ht.IntracellularConductivities()->trans();
1340  double intra_z = ht.IntracellularConductivities()->normal();
1341 
1342  rIntraConductivities.push_back(Create_c_vector(intra_x, intra_y, intra_z));
1343  }
1344  else
1345  {
1346  c_vector<double, 3> intra_conductivities;
1347  GetIntracellularConductivities(intra_conductivities);
1348  rIntraConductivities.push_back(intra_conductivities);
1349  }
1350 
1351  if (ht.ExtracellularConductivities().present())
1352  {
1353  double extra_x = ht.ExtracellularConductivities()->longi();
1354  double extra_y = ht.ExtracellularConductivities()->trans();
1355  double extra_z = ht.ExtracellularConductivities()->normal();
1356 
1357  rExtraConductivities.push_back(Create_c_vector(extra_x, extra_y, extra_z));
1358  }
1359  else
1360  {
1361  c_vector<double, 3> extra_conductivities;
1362  GetExtracellularConductivities(extra_conductivities);
1363  rExtraConductivities.push_back(extra_conductivities);
1364  }
1365  }
1366 }
1367 
1369 {
1370  CheckSimulationIsDefined("Simulation/OutputDirectory");
1371  CHECK_EXISTS(mpParameters->Simulation()->OutputDirectory().present(), "Simulation/OutputDirectory");
1372  return mpParameters->Simulation()->OutputDirectory().get();
1373 }
1374 
1376 {
1377  CheckSimulationIsDefined("Simulation/OutputFilenamePrefix");
1378  CHECK_EXISTS(mpParameters->Simulation()->OutputFilenamePrefix().present(), "Simulation/OutputFilenamePrefix");
1379  return mpParameters->Simulation()->OutputFilenamePrefix().get();
1380 }
1381 
1383 {
1384  CheckSimulationIsDefined("OutputVariables");
1385  return mpParameters->Simulation()->OutputVariables().present();
1386 }
1387 
1388 void HeartConfig::GetOutputVariables(std::vector<std::string>& rOutputVariables) const
1389 {
1390  CHECK_EXISTS(GetOutputVariablesProvided(), "Simulation/OutputVariables");
1391  XSD_SEQUENCE_TYPE(cp::output_variables_type::Var)& output_variables = mpParameters->Simulation()->OutputVariables()->Var();
1392  rOutputVariables.clear();
1393 
1394  for (XSD_ITERATOR_TYPE(cp::output_variables_type::Var) i = output_variables.begin();
1395  i != output_variables.end();
1396  ++i)
1397  {
1398  cp::var_type& r_var(*i);
1399 
1400  // Add to outputVariables the string returned by var.name()
1401  rOutputVariables.push_back(r_var.name());
1402  }
1403 }
1404 
1406 {
1407  CheckSimulationIsDefined("OutputUsingOriginalNodeOrdering");
1408  bool result = false;
1409  if (mpParameters->Simulation()->OutputUsingOriginalNodeOrdering().present())
1410  {
1411  result = (mpParameters->Simulation()->OutputUsingOriginalNodeOrdering().get() == cp::yesno_type::yes);
1412  }
1413  return result;
1414 }
1415 
1417 {
1418  return IsSimulationDefined() && mpParameters->Simulation()->CheckpointSimulation().present();
1419 }
1420 
1422 {
1423  CHECK_EXISTS(GetCheckpointSimulation(), "Simulation/CheckpointSimulation");
1424  return mpParameters->Simulation()->CheckpointSimulation()->timestep();
1425 }
1426 
1428 {
1429  CHECK_EXISTS(GetCheckpointSimulation(), "Simulation/CheckpointSimulation");
1430  return mpParameters->Simulation()->CheckpointSimulation()->max_checkpoints_on_disk();
1431 }
1432 
1434 {
1435  CheckResumeSimulationIsDefined("GetArchivedSimulationDir");
1436 
1437  return HeartFileFinder(mpParameters->ResumeSimulation()->ArchiveDirectory());
1438 }
1439 
1440 void HeartConfig::GetIntracellularConductivities(c_vector<double, 3>& rIntraConductivities) const
1441 {
1442  CHECK_EXISTS(mpParameters->Physiological().IntracellularConductivities().present(), "Physiological/IntracellularConductivities");
1443  cp::conductivities_type intra_conductivities
1444  = mpParameters->Physiological().IntracellularConductivities().get();
1445  double intra_x_cond = intra_conductivities.longi();
1446  double intra_y_cond = intra_conductivities.trans();
1447  double intra_z_cond = intra_conductivities.normal();
1448  ;
1449 
1450  assert(intra_y_cond != DBL_MAX);
1451  assert(intra_z_cond != DBL_MAX);
1452 
1453  rIntraConductivities[0] = intra_x_cond;
1454  rIntraConductivities[1] = intra_y_cond;
1455  rIntraConductivities[2] = intra_z_cond;
1456 }
1457 
1458 void HeartConfig::GetIntracellularConductivities(c_vector<double, 2>& rIntraConductivities) const
1459 {
1460  CHECK_EXISTS(mpParameters->Physiological().IntracellularConductivities().present(), "Physiological/IntracellularConductivities");
1461  cp::conductivities_type intra_conductivities
1462  = mpParameters->Physiological().IntracellularConductivities().get();
1463  double intra_x_cond = intra_conductivities.longi();
1464  double intra_y_cond = intra_conductivities.trans();
1465 
1466  assert(intra_y_cond != DBL_MAX);
1467 
1468  rIntraConductivities[0] = intra_x_cond;
1469  rIntraConductivities[1] = intra_y_cond;
1470 }
1471 
1472 void HeartConfig::GetIntracellularConductivities(c_vector<double, 1>& rIntraConductivities) const
1473 {
1474  CHECK_EXISTS(mpParameters->Physiological().IntracellularConductivities().present(), "Physiological/IntracellularConductivities");
1475  cp::conductivities_type intra_conductivities
1476  = mpParameters->Physiological().IntracellularConductivities().get();
1477  double intra_x_cond = intra_conductivities.longi();
1478 
1479  rIntraConductivities[0] = intra_x_cond;
1480 }
1481 
1482 void HeartConfig::GetExtracellularConductivities(c_vector<double, 3>& rExtraConductivities) const
1483 {
1484  CHECK_EXISTS(mpParameters->Physiological().ExtracellularConductivities().present(), "Physiological/ExtracellularConductivities");
1485  cp::conductivities_type extra_conductivities
1486  = mpParameters->Physiological().ExtracellularConductivities().get();
1487  double extra_x_cond = extra_conductivities.longi();
1488  double extra_y_cond = extra_conductivities.trans();
1489  double extra_z_cond = extra_conductivities.normal();
1490  ;
1491 
1492  assert(extra_y_cond != DBL_MAX);
1493  assert(extra_z_cond != DBL_MAX);
1494 
1495  rExtraConductivities[0] = extra_x_cond;
1496  rExtraConductivities[1] = extra_y_cond;
1497  rExtraConductivities[2] = extra_z_cond;
1498 }
1499 
1500 void HeartConfig::GetExtracellularConductivities(c_vector<double, 2>& rExtraConductivities) const
1501 {
1502  CHECK_EXISTS(mpParameters->Physiological().ExtracellularConductivities().present(), "Physiological/ExtracellularConductivities");
1503  cp::conductivities_type extra_conductivities
1504  = mpParameters->Physiological().ExtracellularConductivities().get();
1505  double extra_x_cond = extra_conductivities.longi();
1506  double extra_y_cond = extra_conductivities.trans();
1507 
1508  assert(extra_y_cond != DBL_MAX);
1509 
1510  rExtraConductivities[0] = extra_x_cond;
1511  rExtraConductivities[1] = extra_y_cond;
1512 }
1513 
1514 void HeartConfig::GetExtracellularConductivities(c_vector<double, 1>& rExtraConductivities) const
1515 {
1516  CHECK_EXISTS(mpParameters->Physiological().ExtracellularConductivities().present(), "Physiological/ExtracellularConductivities");
1517  cp::conductivities_type extra_conductivities
1518  = mpParameters->Physiological().ExtracellularConductivities().get();
1519  double extra_x_cond = extra_conductivities.longi();
1520 
1521  rExtraConductivities[0] = extra_x_cond;
1522 }
1523 
1524 double HeartConfig::GetBathConductivity(unsigned bathRegion) const
1525 {
1526  /*
1527  * We have to consider three cases: The user asks for ...
1528  * a) ... the default conductivity (bathRegion=UINT_MAX)
1529  * b) ... the conductivity of region defined to be heterogeneous
1530  * c) ... the conductivity of region NOT defined to be heterogeneous
1531  *
1532  * a) and c) should return the same
1533  */
1534 
1535  if (bathRegion == UINT_MAX)
1536  {
1537  /*bath conductivity mS/cm*/
1538  CHECK_EXISTS(mpParameters->Physiological().BathConductivity().present(), "Physiological/BathConductivity");
1539  return mpParameters->Physiological().BathConductivity().get();
1540  }
1541  else
1542  {
1543  assert(HeartRegionCode::IsRegionBath(bathRegion));
1544 
1545  std::map<unsigned, double>::const_iterator map_entry = mBathConductivities.find(bathRegion);
1546 
1547  if (map_entry != mBathConductivities.end())
1548  {
1549  return map_entry->second;
1550  }
1551  else
1552  {
1553  /*bath conductivity mS/cm*/
1554  CHECK_EXISTS(mpParameters->Physiological().BathConductivity().present(), "Physiological/BathConductivity");
1555  return mpParameters->Physiological().BathConductivity().get();
1556  }
1557  }
1558 }
1559 const std::set<unsigned>& HeartConfig::rGetTissueIdentifiers()
1560 {
1561  return mTissueIdentifiers;
1562 }
1563 
1564 const std::set<unsigned>& HeartConfig::rGetBathIdentifiers()
1565 {
1566  return mBathIdentifiers;
1567 }
1568 
1570 {
1571  CHECK_EXISTS(mpParameters->Physiological().SurfaceAreaToVolumeRatio().present(), "Physiological/SurfaceAreaToVolumeRatio");
1572  return mpParameters->Physiological().SurfaceAreaToVolumeRatio().get();
1573 }
1574 
1576 {
1577  CHECK_EXISTS(mpParameters->Physiological().Capacitance().present(), "Physiological/Capacitance");
1578  return mpParameters->Physiological().Capacitance().get();
1579 }
1580 
1582 {
1583  CHECK_EXISTS(mpParameters->Numerical().TimeSteps().present(), "Numerical/TimeSteps");
1584  return mpParameters->Numerical().TimeSteps()->ode();
1585 }
1586 
1588 {
1589  CHECK_EXISTS(mpParameters->Numerical().TimeSteps().present(), "Numerical/TimeSteps");
1590  return mpParameters->Numerical().TimeSteps()->pde();
1591 }
1592 
1594 {
1595  CHECK_EXISTS(mpParameters->Numerical().TimeSteps().present(), "Numerical/TimeSteps");
1596  return mpParameters->Numerical().TimeSteps()->printing();
1597 }
1598 
1600 {
1601  CHECK_EXISTS(mpParameters->Numerical().KSPTolerances().present(), "Numerical/KSPTolerances");
1602  return mpParameters->Numerical().KSPTolerances()->KSPAbsolute().present();
1603 }
1604 
1606 {
1607  CHECK_EXISTS(mpParameters->Numerical().KSPTolerances().present(), "Numerical/KSPTolerances");
1608  if (!GetUseAbsoluteTolerance())
1609  {
1610  EXCEPTION("Absolute tolerance is not set in Chaste parameters");
1611  }
1612  return mpParameters->Numerical().KSPTolerances()->KSPAbsolute().get();
1613 }
1614 
1616 {
1617  CHECK_EXISTS(mpParameters->Numerical().KSPTolerances().present(), "Numerical/KSPTolerances");
1618  return mpParameters->Numerical().KSPTolerances()->KSPRelative().present();
1619 }
1620 
1622 {
1623  CHECK_EXISTS(mpParameters->Numerical().KSPTolerances().present(), "Numerical/KSPTolerances");
1624  if (!GetUseRelativeTolerance())
1625  {
1626  EXCEPTION("Relative tolerance is not set in Chaste parameters");
1627  }
1628  return mpParameters->Numerical().KSPTolerances()->KSPRelative().get();
1629 }
1630 
1631 const char* HeartConfig::GetKSPSolver() const
1632 {
1633  CHECK_EXISTS(mpParameters->Numerical().KSPSolver().present(), "Numerical/KSPSolver");
1634  switch (mpParameters->Numerical().KSPSolver().get())
1635  {
1636  case cp::ksp_solver_type::gmres:
1637  return "gmres";
1638  case cp::ksp_solver_type::cg:
1639  return "cg";
1640  case cp::ksp_solver_type::symmlq:
1641  return "symmlq";
1642  case cp::ksp_solver_type::chebychev:
1643  return "chebychev";
1644  }
1645  // LCOV_EXCL_START
1646  EXCEPTION("Unknown ksp solver");
1647  // LCOV_EXCL_STOP
1648 }
1649 
1651 {
1652  CHECK_EXISTS(mpParameters->Numerical().KSPPreconditioner().present(), "Numerical/KSPPreconditioner");
1653  switch (mpParameters->Numerical().KSPPreconditioner().get())
1654  {
1655  case cp::ksp_preconditioner_type::jacobi:
1656  return "jacobi";
1657  case cp::ksp_preconditioner_type::bjacobi:
1658  return "bjacobi";
1659  case cp::ksp_preconditioner_type::hypre:
1660  return "hypre";
1661  case cp::ksp_preconditioner_type::ml:
1662  return "ml";
1663  case cp::ksp_preconditioner_type::spai:
1664  return "spai";
1665  case cp::ksp_preconditioner_type::blockdiagonal:
1666  return "blockdiagonal";
1667  case cp::ksp_preconditioner_type::ldufactorisation:
1668  return "ldufactorisation";
1669  case cp::ksp_preconditioner_type::twolevelsblockdiagonal:
1670  return "twolevelsblockdiagonal";
1671  case cp::ksp_preconditioner_type::none:
1672  return "none";
1673  }
1674  // LCOV_EXCL_START
1675  EXCEPTION("Unknown ksp preconditioner");
1676  // LCOV_EXCL_STOP
1677 }
1678 
1680 {
1681  CHECK_EXISTS(mpParameters->Numerical().MeshPartitioning().present(), "Numerical/MeshPartitioning");
1682  switch (mpParameters->Numerical().MeshPartitioning().get())
1683  {
1684  case cp::mesh_partitioning_type::dumb:
1685  return DistributedTetrahedralMeshPartitionType::DUMB;
1686  case cp::mesh_partitioning_type::metis:
1687  return DistributedTetrahedralMeshPartitionType::METIS_LIBRARY;
1688  case cp::mesh_partitioning_type::parmetis:
1689  return DistributedTetrahedralMeshPartitionType::PARMETIS_LIBRARY;
1690  case cp::mesh_partitioning_type::petsc:
1691  return DistributedTetrahedralMeshPartitionType::PETSC_MAT_PARTITION;
1692  }
1693  // LCOV_EXCL_START
1694  EXCEPTION("Unknown mesh partitioning type");
1695  // LCOV_EXCL_STOP
1696 }
1697 
1699 {
1700  bool IsAdaptivityParametersPresent = mpParameters->Numerical().AdaptivityParameters().present();
1701  if (IsAdaptivityParametersPresent)
1702  {
1703  WARNING("Use of the Adaptivity library is deprecated");
1704  }
1706 }
1707 
1708 /*
1709  * PostProcessing
1710  */
1711 
1713 {
1714  return mpParameters->PostProcessing().present();
1715 }
1716 
1718 {
1719  ENSURE_SECTION_PRESENT(mpParameters->PostProcessing(), cp::postprocessing_type);
1720 }
1721 
1723 {
1724  if (IsPostProcessingSectionPresent() == false)
1725  {
1726  return false;
1727  }
1728  else
1729  {
1731  }
1732 }
1734 {
1735  bool result = false;
1737  {
1738  XSD_SEQUENCE_TYPE(cp::postprocessing_type::ActionPotentialDurationMap)& apd_maps = mpParameters->PostProcessing()->ActionPotentialDurationMap();
1739  result = (apd_maps.begin() != apd_maps.end());
1740  }
1741  return result;
1742 }
1743 
1744 void HeartConfig::GetApdMaps(std::vector<std::pair<double, double> >& apd_maps) const
1745 {
1746  CHECK_EXISTS(IsApdMapsRequested(), "PostProcessing/ActionPotentialDurationMap");
1747  apd_maps.clear();
1748 
1749  XSD_SEQUENCE_TYPE(cp::postprocessing_type::ActionPotentialDurationMap)& apd_maps_sequence = mpParameters->PostProcessing()->ActionPotentialDurationMap();
1750 
1751  for (XSD_ITERATOR_TYPE(cp::postprocessing_type::ActionPotentialDurationMap) i = apd_maps_sequence.begin();
1752  i != apd_maps_sequence.end();
1753  ++i)
1754  {
1755  std::pair<double, double> map(i->repolarisation_percentage(), i->threshold());
1756 
1757  apd_maps.push_back(map);
1758  }
1759 }
1760 
1762 {
1763  bool result = false;
1765  {
1766  XSD_SEQUENCE_TYPE(cp::postprocessing_type::UpstrokeTimeMap)& upstroke_map = mpParameters->PostProcessing()->UpstrokeTimeMap();
1767  result = (upstroke_map.begin() != upstroke_map.end());
1768  }
1769  return result;
1770 }
1771 void HeartConfig::GetUpstrokeTimeMaps(std::vector<double>& upstroke_time_maps) const
1772 {
1773  CHECK_EXISTS(IsUpstrokeTimeMapsRequested(), "PostProcessing/UpstrokeTimeMap");
1774  assert(upstroke_time_maps.size() == 0);
1775 
1776  XSD_SEQUENCE_TYPE(cp::postprocessing_type::UpstrokeTimeMap)& upstroke_maps_sequence = mpParameters->PostProcessing()->UpstrokeTimeMap();
1777 
1778  for (XSD_ITERATOR_TYPE(cp::postprocessing_type::UpstrokeTimeMap) i = upstroke_maps_sequence.begin();
1779  i != upstroke_maps_sequence.end();
1780  ++i)
1781  {
1782  upstroke_time_maps.push_back(i->threshold());
1783  }
1784 }
1785 
1787 {
1788  bool result = false;
1790  {
1791  XSD_SEQUENCE_TYPE(cp::postprocessing_type::MaxUpstrokeVelocityMap)& max_upstroke_velocity_map = mpParameters->PostProcessing()->MaxUpstrokeVelocityMap();
1792  result = (max_upstroke_velocity_map.begin() != max_upstroke_velocity_map.end());
1793  }
1794  return result;
1795 }
1796 
1797 void HeartConfig::GetMaxUpstrokeVelocityMaps(std::vector<double>& upstroke_velocity_maps) const
1798 {
1799  CHECK_EXISTS(IsMaxUpstrokeVelocityMapRequested(), "PostProcessing/MaxUpstrokeVelocityMap");
1800  assert(upstroke_velocity_maps.size() == 0);
1801 
1802  XSD_SEQUENCE_TYPE(cp::postprocessing_type::MaxUpstrokeVelocityMap)& max_upstroke_velocity_maps_sequence = mpParameters->PostProcessing()->MaxUpstrokeVelocityMap();
1803 
1804  for (XSD_ITERATOR_TYPE(cp::postprocessing_type::MaxUpstrokeVelocityMap) i = max_upstroke_velocity_maps_sequence.begin();
1805  i != max_upstroke_velocity_maps_sequence.end();
1806  ++i)
1807  {
1808  upstroke_velocity_maps.push_back(i->threshold());
1809  }
1810 }
1811 
1813 {
1814  bool result = false;
1816  {
1817  XSD_SEQUENCE_TYPE(cp::postprocessing_type::ConductionVelocityMap)& cond_vel_maps = mpParameters->PostProcessing()->ConductionVelocityMap();
1818  result = (cond_vel_maps.begin() != cond_vel_maps.end());
1819  }
1820  return result;
1821 }
1822 
1823 void HeartConfig::GetConductionVelocityMaps(std::vector<unsigned>& conduction_velocity_maps) const
1824 {
1825  CHECK_EXISTS(IsConductionVelocityMapsRequested(), "PostProcessing/ConductionVelocityMap");
1826  assert(conduction_velocity_maps.size() == 0);
1827 
1828  XSD_SEQUENCE_TYPE(cp::postprocessing_type::ConductionVelocityMap)& cond_vel_maps_sequence = mpParameters->PostProcessing()->ConductionVelocityMap();
1829 
1830  for (XSD_ITERATOR_TYPE(cp::postprocessing_type::ConductionVelocityMap) i = cond_vel_maps_sequence.begin();
1831  i != cond_vel_maps_sequence.end();
1832  ++i)
1833  {
1834  conduction_velocity_maps.push_back(i->origin_node());
1835  }
1836 }
1837 
1839 {
1840  bool result = false;
1842  {
1843  XSD_SEQUENCE_TYPE(cp::postprocessing_type::TimeTraceAtNode)& requested_nodes = mpParameters->PostProcessing()->TimeTraceAtNode();
1844  result = (requested_nodes.begin() != requested_nodes.end());
1845  }
1846  return result;
1847 }
1848 
1849 void HeartConfig::GetNodalTimeTraceRequested(std::vector<unsigned>& rRequestedNodes) const
1850 {
1851  CHECK_EXISTS(IsAnyNodalTimeTraceRequested(), "PostProcessing/TimeTraceAtNode");
1852  assert(rRequestedNodes.size() == 0);
1853 
1854  XSD_SEQUENCE_TYPE(cp::postprocessing_type::TimeTraceAtNode)& req_nodes = mpParameters->PostProcessing()->TimeTraceAtNode();
1855 
1856  for (XSD_ITERATOR_TYPE(cp::postprocessing_type::TimeTraceAtNode) i = req_nodes.begin();
1857  i != req_nodes.end();
1858  ++i)
1859  {
1860  rRequestedNodes.push_back(i->node_number());
1861  }
1862 }
1863 
1865 {
1866  bool result = false;
1868  {
1869  XSD_SEQUENCE_TYPE(cp::postprocessing_type::PseudoEcgElectrodePosition)& electrodes = mpParameters->PostProcessing()->PseudoEcgElectrodePosition();
1870  result = (electrodes.begin() != electrodes.end());
1871  }
1872  return result;
1873 }
1874 
1875 template <unsigned SPACE_DIM>
1876 void HeartConfig::GetPseudoEcgElectrodePositions(std::vector<ChastePoint<SPACE_DIM> >& rPseudoEcgElectrodePositions) const
1877 {
1878  rPseudoEcgElectrodePositions.clear();
1879  XSD_SEQUENCE_TYPE(cp::postprocessing_type::PseudoEcgElectrodePosition)& electrodes = mpParameters->PostProcessing()->PseudoEcgElectrodePosition();
1880  for (XSD_ITERATOR_TYPE(cp::postprocessing_type::PseudoEcgElectrodePosition) i = electrodes.begin();
1881  i != electrodes.end();
1882  ++i)
1883  {
1884  rPseudoEcgElectrodePositions.push_back(ChastePoint<SPACE_DIM>(i->x(), i->y(), i->z()));
1885  }
1886 }
1887 
1888 /*
1889  * Output visualization
1890  */
1891 
1893 {
1894  CheckSimulationIsDefined("OutputVisualizer");
1895 
1896  return mpParameters->Simulation()->OutputVisualizer().present();
1897 }
1898 
1900 {
1902  {
1903  return false;
1904  }
1905  else
1906  {
1907  return mpParameters->Simulation()->OutputVisualizer()->meshalyzer() == cp::yesno_type::yes;
1908  }
1909 }
1910 
1912 {
1914  {
1915  return false;
1916  }
1917  else
1918  {
1919  return mpParameters->Simulation()->OutputVisualizer()->cmgui() == cp::yesno_type::yes;
1920  }
1921 }
1922 
1924 {
1926  {
1927  return false;
1928  }
1929  else
1930  {
1931  return mpParameters->Simulation()->OutputVisualizer()->parallel_vtk() == cp::yesno_type::yes;
1932  }
1933 }
1934 
1936 {
1938  {
1939  return false;
1940  }
1941  else
1942  {
1943  return mpParameters->Simulation()->OutputVisualizer()->vtk() == cp::yesno_type::yes;
1944  }
1945 }
1946 
1948 {
1950  {
1951  return 0u;
1952  }
1953  else
1954  {
1955  return mpParameters->Simulation()->OutputVisualizer()->precision();
1956  }
1957 }
1958 
1960 {
1961  return mpParameters->Simulation()->Electrodes().present();
1962 }
1963 
1964 /*
1965  * Set methods
1966  */
1967 void HeartConfig::SetSpaceDimension(unsigned spaceDimension)
1968 {
1969  mpParameters->Simulation()->SpaceDimension().set(spaceDimension);
1970 }
1971 
1972 void HeartConfig::SetSimulationDuration(double simulationDuration)
1973 {
1974  XSD_CREATE_WITH_FIXED_ATTR1(cp::time_type, time, simulationDuration, "ms");
1975  mpParameters->Simulation()->SimulationDuration().set(time);
1976 }
1977 
1978 void HeartConfig::SetDomain(const cp::domain_type& rDomain)
1979 {
1980  mpParameters->Simulation()->Domain().set(rDomain);
1981 }
1982 
1983 void HeartConfig::SetDefaultIonicModel(const cp::ionic_models_available_type& rIonicModel)
1984 {
1985  cp::ionic_model_selection_type ionic_model;
1986  ionic_model.Hardcoded(rIonicModel);
1987  cp::ionic_models_type container(ionic_model);
1988  mpParameters->Simulation()->IonicModels().set(container);
1989 }
1990 
1991 void HeartConfig::SetSlabDimensions(double x, double y, double z, double inter_node_space)
1992 {
1993  if (!mpParameters->Simulation()->Mesh().present())
1994  {
1995  XSD_CREATE_WITH_FIXED_ATTR(cp::mesh_type, mesh_to_load, "cm");
1996  mpParameters->Simulation()->Mesh().set(mesh_to_load);
1997  }
1998 
1999  cp::slab_type slab_definition(x, y, z, inter_node_space);
2000  mpParameters->Simulation()->Mesh()->Slab().set(slab_definition);
2001 }
2002 
2003 void HeartConfig::SetSheetDimensions(double x, double y, double inter_node_space)
2004 {
2005  if (!mpParameters->Simulation()->Mesh().present())
2006  {
2007  XSD_CREATE_WITH_FIXED_ATTR(cp::mesh_type, mesh_to_load, "cm");
2008  mpParameters->Simulation()->Mesh().set(mesh_to_load);
2009  }
2010 
2011  cp::sheet_type sheet_definition(x, y, inter_node_space);
2012  mpParameters->Simulation()->Mesh()->Sheet().set(sheet_definition);
2013 }
2014 
2015 void HeartConfig::SetFibreLength(double x, double inter_node_space)
2016 {
2017  if (!mpParameters->Simulation()->Mesh().present())
2018  {
2019  XSD_CREATE_WITH_FIXED_ATTR(cp::mesh_type, mesh_to_load, "cm");
2020  mpParameters->Simulation()->Mesh().set(mesh_to_load);
2021  }
2022 
2023  cp::fibre_type fibre_definition(x, inter_node_space);
2024  mpParameters->Simulation()->Mesh()->Fibre().set(fibre_definition);
2025 }
2026 
2027 void HeartConfig::SetMeshFileName(std::string meshPrefix, cp::media_type fibreDefinition)
2028 {
2029  if (!mpParameters->Simulation()->Mesh().present())
2030  {
2031  XSD_CREATE_WITH_FIXED_ATTR(cp::mesh_type, mesh_to_load, "cm");
2032  mpParameters->Simulation()->Mesh().set(mesh_to_load);
2033  }
2034 
2035  XSD_NESTED_TYPE(cp::mesh_type::LoadMesh)
2036  mesh_prefix(meshPrefix, fibreDefinition);
2037  mpParameters->Simulation()->Mesh()->LoadMesh().set(mesh_prefix);
2038 }
2039 
2040 void HeartConfig::SetIonicModelRegions(std::vector<ChasteCuboid<3> >& rDefinedRegions,
2041  std::vector<cp::ionic_model_selection_type>& rIonicModels) const
2042 {
2043  assert(rDefinedRegions.size() == rIonicModels.size());
2044  // You need to have defined a default model first...
2045  assert(mpParameters->Simulation()->IonicModels().present());
2046  XSD_SEQUENCE_TYPE(cp::ionic_models_type::Region)& regions = mpParameters->Simulation()->IonicModels()->Region();
2047  regions.clear();
2048  for (unsigned region_index = 0; region_index < rDefinedRegions.size(); region_index++)
2049  {
2050  cp::point_type point_a(rDefinedRegions[region_index].rGetLowerCorner()[0],
2051  rDefinedRegions[region_index].rGetLowerCorner()[1],
2052  rDefinedRegions[region_index].rGetLowerCorner()[2]);
2053 
2054  cp::point_type point_b(rDefinedRegions[region_index].rGetUpperCorner()[0],
2055  rDefinedRegions[region_index].rGetUpperCorner()[1],
2056  rDefinedRegions[region_index].rGetUpperCorner()[2]);
2057 
2058  XSD_CREATE_WITH_FIXED_ATTR(cp::location_type, locn, "cm");
2059  locn.Cuboid().set(cp::box_type(point_a, point_b));
2060 
2061  cp::ionic_model_region_type region(rIonicModels[region_index], locn);
2062  regions.push_back(region);
2063  }
2064 }
2065 
2066 void HeartConfig::SetConductivityHeterogeneities(std::vector<ChasteCuboid<3> >& rConductivityAreas,
2067  std::vector<c_vector<double, 3> >& rIntraConductivities,
2068  std::vector<c_vector<double, 3> >& rExtraConductivities)
2069 {
2070  assert(rConductivityAreas.size() == rIntraConductivities.size());
2071  assert(rIntraConductivities.size() == rExtraConductivities.size());
2072 
2073  XSD_ANON_SEQUENCE_TYPE(cp::physiological_type, ConductivityHeterogeneities, ConductivityHeterogeneity)
2074  heterogeneities_container;
2075 
2076  for (unsigned region_index = 0; region_index < rConductivityAreas.size(); region_index++)
2077  {
2078  cp::point_type point_a(rConductivityAreas[region_index].rGetLowerCorner()[0],
2079  rConductivityAreas[region_index].rGetLowerCorner()[1],
2080  rConductivityAreas[region_index].rGetLowerCorner()[2]);
2081 
2082  cp::point_type point_b(rConductivityAreas[region_index].rGetUpperCorner()[0],
2083  rConductivityAreas[region_index].rGetUpperCorner()[1],
2084  rConductivityAreas[region_index].rGetUpperCorner()[2]);
2085 
2086  XSD_CREATE_WITH_FIXED_ATTR(cp::location_type, locn, "cm");
2087  locn.Cuboid().set(cp::box_type(point_a, point_b));
2088  cp::conductivity_heterogeneity_type ht(locn);
2089 
2090  XSD_CREATE_WITH_FIXED_ATTR3(cp::conductivities_type, intra,
2091  rIntraConductivities[region_index][0],
2092  rIntraConductivities[region_index][1],
2093  rIntraConductivities[region_index][2],
2094  "mS/cm");
2095 
2096  ht.IntracellularConductivities(intra);
2097 
2098  XSD_CREATE_WITH_FIXED_ATTR3(cp::conductivities_type, extra,
2099  rExtraConductivities[region_index][0],
2100  rExtraConductivities[region_index][1],
2101  rExtraConductivities[region_index][2],
2102  "mS/cm");
2103 
2104  ht.ExtracellularConductivities(extra);
2105 
2106  heterogeneities_container.push_back(ht);
2107  }
2108 
2109  XSD_ANON_TYPE(cp::physiological_type, ConductivityHeterogeneities)
2110  heterogeneities_object;
2111  heterogeneities_object.ConductivityHeterogeneity(heterogeneities_container);
2112 
2113  mpParameters->Physiological().ConductivityHeterogeneities().set(heterogeneities_object);
2114 }
2115 
2117  std::vector<c_vector<double, 3> >& rIntraConductivities,
2118  std::vector<c_vector<double, 3> >& rExtraConductivities)
2119 {
2120  assert(rConductivityAreas.size() == rIntraConductivities.size());
2121  assert(rIntraConductivities.size() == rExtraConductivities.size());
2122 
2123  XSD_ANON_SEQUENCE_TYPE(cp::physiological_type, ConductivityHeterogeneities, ConductivityHeterogeneity)
2124  heterogeneities_container;
2125 
2126  for (unsigned region_index = 0; region_index < rConductivityAreas.size(); region_index++)
2127  {
2128  cp::point_type centre(rConductivityAreas[region_index].rGetCentre()[0],
2129  rConductivityAreas[region_index].rGetCentre()[1],
2130  rConductivityAreas[region_index].rGetCentre()[2]);
2131 
2132  cp::point_type radii(rConductivityAreas[region_index].rGetRadii()[0],
2133  rConductivityAreas[region_index].rGetRadii()[1],
2134  rConductivityAreas[region_index].rGetRadii()[2]);
2135 
2136  XSD_CREATE_WITH_FIXED_ATTR(cp::location_type, locn, "cm");
2137  locn.Ellipsoid().set(cp::ellipsoid_type(centre, radii));
2138  cp::conductivity_heterogeneity_type ht(locn);
2139 
2140  XSD_CREATE_WITH_FIXED_ATTR3(cp::conductivities_type, intra,
2141  rIntraConductivities[region_index][0],
2142  rIntraConductivities[region_index][1],
2143  rIntraConductivities[region_index][2],
2144  "mS/cm");
2145 
2146  ht.IntracellularConductivities(intra);
2147 
2148  XSD_CREATE_WITH_FIXED_ATTR3(cp::conductivities_type, extra,
2149  rExtraConductivities[region_index][0],
2150  rExtraConductivities[region_index][1],
2151  rExtraConductivities[region_index][2],
2152  "mS/cm");
2153 
2154  ht.ExtracellularConductivities(extra);
2155 
2156  heterogeneities_container.push_back(ht);
2157  }
2158 
2159  XSD_ANON_TYPE(cp::physiological_type, ConductivityHeterogeneities)
2160  heterogeneities_object;
2161  heterogeneities_object.ConductivityHeterogeneity(heterogeneities_container);
2162 
2163  mpParameters->Physiological().ConductivityHeterogeneities().set(heterogeneities_object);
2164 }
2165 
2166 void HeartConfig::SetOutputDirectory(const std::string& rOutputDirectory)
2167 {
2168  mpParameters->Simulation()->OutputDirectory().set(rOutputDirectory);
2169 }
2170 
2171 void HeartConfig::SetOutputFilenamePrefix(const std::string& rOutputFilenamePrefix)
2172 {
2173  mpParameters->Simulation()->OutputFilenamePrefix().set(rOutputFilenamePrefix);
2174 }
2175 
2176 void HeartConfig::SetOutputVariables(const std::vector<std::string>& rOutputVariables)
2177 {
2178  if (!mpParameters->Simulation()->OutputVariables().present())
2179  {
2180  cp::output_variables_type variables_requested;
2181  mpParameters->Simulation()->OutputVariables().set(variables_requested);
2182  }
2183 
2184  XSD_SEQUENCE_TYPE(cp::output_variables_type::Var)& var_type_sequence = mpParameters->Simulation()->OutputVariables()->Var();
2185  // Erase or create a sequence
2186  var_type_sequence.clear();
2187 
2188  for (unsigned i = 0; i < rOutputVariables.size(); i++)
2189  {
2190  cp::var_type temp(rOutputVariables[i]);
2191  var_type_sequence.push_back(temp);
2192  }
2193 }
2194 
2196 {
2197  //What if it doesn't exist?
2198  mpParameters->Simulation()->OutputUsingOriginalNodeOrdering().set(useOriginal ? cp::yesno_type::yes : cp::yesno_type::no);
2199 }
2200 
2201 void HeartConfig::SetCheckpointSimulation(bool saveSimulation, double checkpointTimestep, unsigned maxCheckpointsOnDisk)
2202 {
2203  if (saveSimulation)
2204  {
2205  // Make sure values for the optional parameters have been provided
2206  assert(checkpointTimestep != -1.0 && maxCheckpointsOnDisk != UINT_MAX);
2207 
2208  XSD_CREATE_WITH_FIXED_ATTR2(cp::simulation_type::XSD_NESTED_TYPE(CheckpointSimulation),
2209  cs,
2210  checkpointTimestep,
2211  maxCheckpointsOnDisk,
2212  "ms");
2213  mpParameters->Simulation()->CheckpointSimulation().set(cs);
2214  }
2215  else
2216  {
2217  mpParameters->Simulation()->CheckpointSimulation().reset();
2218  }
2219 
2220  CheckTimeSteps();
2221 }
2222 
2223 // Physiological
2224 
2225 void HeartConfig::SetIntracellularConductivities(const c_vector<double, 3>& rIntraConductivities)
2226 {
2227  XSD_CREATE_WITH_FIXED_ATTR3(cp::conductivities_type, intra,
2228  rIntraConductivities[0],
2229  rIntraConductivities[1],
2230  rIntraConductivities[2],
2231  "mS/cm");
2232 
2233  mpParameters->Physiological().IntracellularConductivities().set(intra);
2234 }
2235 
2236 void HeartConfig::SetIntracellularConductivities(const c_vector<double, 2>& rIntraConductivities)
2237 {
2238  XSD_CREATE_WITH_FIXED_ATTR3(cp::conductivities_type, intra,
2239  rIntraConductivities[0],
2240  rIntraConductivities[1],
2241  0.0, "mS/cm");
2242 
2243  mpParameters->Physiological().IntracellularConductivities().set(intra);
2244 }
2245 
2246 void HeartConfig::SetIntracellularConductivities(const c_vector<double, 1>& rIntraConductivities)
2247 {
2248  XSD_CREATE_WITH_FIXED_ATTR3(cp::conductivities_type, intra,
2249  rIntraConductivities[0],
2250  0.0, 0.0, "mS/cm");
2251 
2252  mpParameters->Physiological().IntracellularConductivities().set(intra);
2253 }
2254 
2255 void HeartConfig::SetExtracellularConductivities(const c_vector<double, 3>& rExtraConductivities)
2256 {
2257  XSD_CREATE_WITH_FIXED_ATTR3(cp::conductivities_type, extra,
2258  rExtraConductivities[0],
2259  rExtraConductivities[1],
2260  rExtraConductivities[2],
2261  "mS/cm");
2262 
2263  mpParameters->Physiological().ExtracellularConductivities().set(extra);
2264 }
2265 
2266 void HeartConfig::SetExtracellularConductivities(const c_vector<double, 2>& rExtraConductivities)
2267 {
2268  XSD_CREATE_WITH_FIXED_ATTR3(cp::conductivities_type, extra,
2269  rExtraConductivities[0],
2270  rExtraConductivities[1],
2271  0.0, "mS/cm");
2272 
2273  mpParameters->Physiological().ExtracellularConductivities().set(extra);
2274 }
2275 
2276 void HeartConfig::SetExtracellularConductivities(const c_vector<double, 1>& rExtraConductivities)
2277 {
2278  XSD_CREATE_WITH_FIXED_ATTR3(cp::conductivities_type, extra,
2279  rExtraConductivities[0],
2280  0.0, 0.0, "mS/cm");
2281 
2282  mpParameters->Physiological().ExtracellularConductivities().set(extra);
2283 }
2284 
2285 void HeartConfig::SetBathConductivity(double bathConductivity)
2286 {
2287  XSD_CREATE_WITH_FIXED_ATTR1(cp::conductivity_type, cond, bathConductivity, "mS/cm");
2288  mpParameters->Physiological().BathConductivity().set(cond);
2289 }
2290 
2291 void HeartConfig::SetBathMultipleConductivities(std::map<unsigned, double> bathConductivities)
2292 {
2294  mBathConductivities = bathConductivities;
2295 }
2296 
2297 //void HeartConfig::SetTissueIdentifiers(const std::set<unsigned>& tissueIds)
2298 //{
2299 // std::set<unsigned> empty_bath_identifiers; //Too dangerous (see GetValidBathId)
2300 // SetTissueAndBathIdentifiers(tissueIds, mBathIdentifiers);
2301 //}
2302 
2303 void HeartConfig::SetTissueAndBathIdentifiers(const std::set<unsigned>& tissueIds, const std::set<unsigned>& bathIds)
2304 {
2305  if (tissueIds.empty() || bathIds.empty())
2306  {
2307  EXCEPTION("Identifying set must be non-empty");
2308  }
2309  std::set<unsigned> shared_identifiers;
2310  std::set_intersection(tissueIds.begin(),
2311  tissueIds.end(),
2312  bathIds.begin(),
2313  bathIds.end(),
2314  std::inserter(shared_identifiers, shared_identifiers.begin()));
2315 
2316  if (!shared_identifiers.empty())
2317  {
2318  EXCEPTION("Tissue identifiers and bath identifiers overlap");
2319  }
2320  mTissueIdentifiers = tissueIds;
2321  mBathIdentifiers = bathIds;
2322 }
2323 
2325 {
2326  XSD_CREATE_WITH_FIXED_ATTR1(cp::inverse_length_type, ratio_object, ratio, "1/cm");
2327  mpParameters->Physiological().SurfaceAreaToVolumeRatio().set(ratio_object);
2328 }
2329 
2330 void HeartConfig::SetCapacitance(double capacitance)
2331 {
2332  XSD_CREATE_WITH_FIXED_ATTR1(cp::capacitance_type, capacitance_object, capacitance, "uF/cm^2");
2333  mpParameters->Physiological().Capacitance().set(capacitance_object);
2334 }
2335 
2336 // Numerical
2337 void HeartConfig::SetOdePdeAndPrintingTimeSteps(double odeTimeStep, double pdeTimeStep, double printingTimeStep)
2338 {
2339  XSD_CREATE_WITH_FIXED_ATTR3(cp::time_steps_type, time_steps,
2340  odeTimeStep, pdeTimeStep, printingTimeStep, "ms");
2341  mpParameters->Numerical().TimeSteps().set(time_steps);
2342  CheckTimeSteps();
2343 }
2344 
2345 void HeartConfig::SetOdeTimeStep(double odeTimeStep)
2346 {
2348 }
2349 
2350 void HeartConfig::SetPdeTimeStep(double pdeTimeStep)
2351 {
2353 }
2354 
2355 void HeartConfig::SetPrintingTimeStep(double printingTimeStep)
2356 {
2358 }
2359 
2361 {
2362  if (GetOdeTimeStep() <= 0)
2363  {
2364  EXCEPTION("Ode time-step should be positive");
2365  }
2366  if (GetPdeTimeStep() <= 0)
2367  {
2368  EXCEPTION("Pde time-step should be positive");
2369  }
2370  if (GetPrintingTimeStep() <= 0.0)
2371  {
2372  EXCEPTION("Printing time-step should be positive");
2373  }
2374 
2376  {
2377  EXCEPTION("Printing time-step should not be smaller than PDE time-step");
2378  }
2379 
2381  {
2382  EXCEPTION("Printing time-step should be a multiple of PDE time step");
2383  }
2384 
2385  if (GetOdeTimeStep() > GetPdeTimeStep())
2386  {
2387  EXCEPTION("Ode time-step should not be greater than PDE time-step");
2388  }
2389 
2391  {
2392  if (GetCheckpointTimestep() <= 0.0)
2393  {
2394  EXCEPTION("Checkpoint time-step should be positive");
2395  }
2396 
2398  {
2399  EXCEPTION("Checkpoint time-step should be a multiple of printing time-step");
2400  }
2401  }
2402 }
2403 
2404 void HeartConfig::SetUseRelativeTolerance(double relativeTolerance)
2405 {
2406  ENSURE_SECTION_PRESENT(mpParameters->Numerical().KSPTolerances(), cp::ksp_tolerances_type);
2407  //Remove any reference to tolerances is user parameters
2408  mpParameters->Numerical().KSPTolerances()->KSPAbsolute().reset();
2409  mpParameters->Numerical().KSPTolerances()->KSPRelative().set(relativeTolerance);
2410 }
2411 
2412 void HeartConfig::SetUseAbsoluteTolerance(double absoluteTolerance)
2413 {
2414  ENSURE_SECTION_PRESENT(mpParameters->Numerical().KSPTolerances(), cp::ksp_tolerances_type);
2415  //Remove any reference to tolerances is user parameters
2416  mpParameters->Numerical().KSPTolerances()->KSPRelative().reset();
2417  mpParameters->Numerical().KSPTolerances()->KSPAbsolute().set(absoluteTolerance);
2418 }
2419 
2420 void HeartConfig::SetKSPSolver(const char* kspSolver, bool warnOfChange)
2421 {
2422  if (warnOfChange && strcmp(GetKSPSolver(), kspSolver) != 0)
2423  {
2424  //Warn
2425  WARNING("Code has changed the KSP solver type from " << GetKSPSolver() << " to " << kspSolver);
2426  }
2427 
2428  /* Note that changes in these conditions need to be reflected in the Doxygen*/
2429  if (strcmp(kspSolver, "gmres") == 0)
2430  {
2431  mpParameters->Numerical().KSPSolver().set(cp::ksp_solver_type::gmres);
2432  return;
2433  }
2434  if (strcmp(kspSolver, "cg") == 0)
2435  {
2436  mpParameters->Numerical().KSPSolver().set(cp::ksp_solver_type::cg);
2437  return;
2438  }
2439  if (strcmp(kspSolver, "symmlq") == 0)
2440  {
2441  mpParameters->Numerical().KSPSolver().set(cp::ksp_solver_type::symmlq);
2442  return;
2443  }
2444  if (strcmp(kspSolver, "chebychev") == 0)
2445  {
2446  mpParameters->Numerical().KSPSolver().set(cp::ksp_solver_type::chebychev);
2447  return;
2448  }
2449 
2450  EXCEPTION("Unknown solver type provided");
2451 }
2452 
2453 void HeartConfig::SetKSPPreconditioner(const char* kspPreconditioner)
2454 {
2455  /* Note that changes in these conditions need to be reflected in the Doxygen*/
2456  if (strcmp(kspPreconditioner, "jacobi") == 0)
2457  {
2458  mpParameters->Numerical().KSPPreconditioner().set(cp::ksp_preconditioner_type::jacobi);
2459  return;
2460  }
2461  if (strcmp(kspPreconditioner, "bjacobi") == 0)
2462  {
2463  mpParameters->Numerical().KSPPreconditioner().set(cp::ksp_preconditioner_type::bjacobi);
2464  return;
2465  }
2466  if (strcmp(kspPreconditioner, "hypre") == 0)
2467  {
2468  mpParameters->Numerical().KSPPreconditioner().set(cp::ksp_preconditioner_type::hypre);
2469  return;
2470  }
2471  if (strcmp(kspPreconditioner, "ml") == 0)
2472  {
2473  mpParameters->Numerical().KSPPreconditioner().set(cp::ksp_preconditioner_type::ml);
2474  return;
2475  }
2476  if (strcmp(kspPreconditioner, "spai") == 0)
2477  {
2478  mpParameters->Numerical().KSPPreconditioner().set(cp::ksp_preconditioner_type::spai);
2479  return;
2480  }
2481  if (strcmp(kspPreconditioner, "twolevelsblockdiagonal") == 0)
2482  {
2483  mpParameters->Numerical().KSPPreconditioner().set(cp::ksp_preconditioner_type::twolevelsblockdiagonal);
2484  return;
2485  }
2486  if (strcmp(kspPreconditioner, "blockdiagonal") == 0)
2487  {
2488  mpParameters->Numerical().KSPPreconditioner().set(cp::ksp_preconditioner_type::blockdiagonal);
2489  return;
2490  }
2491  if (strcmp(kspPreconditioner, "ldufactorisation") == 0)
2492  {
2493  mpParameters->Numerical().KSPPreconditioner().set(cp::ksp_preconditioner_type::ldufactorisation);
2494  return;
2495  }
2496  if (strcmp(kspPreconditioner, "none") == 0)
2497  {
2498  mpParameters->Numerical().KSPPreconditioner().set(cp::ksp_preconditioner_type::none);
2499  return;
2500  }
2501 
2502  EXCEPTION("Unknown preconditioner type provided");
2503 }
2504 
2505 void HeartConfig::SetMeshPartitioning(const char* meshPartioningMethod)
2506 {
2507  /* Note that changes in these conditions need to be reflected in the Doxygen*/
2508  if (strcmp(meshPartioningMethod, "dumb") == 0)
2509  {
2510  mpParameters->Numerical().MeshPartitioning().set(cp::mesh_partitioning_type::dumb);
2511  return;
2512  }
2513  if (strcmp(meshPartioningMethod, "metis") == 0)
2514  {
2515  mpParameters->Numerical().MeshPartitioning().set(cp::mesh_partitioning_type::metis);
2516  return;
2517  }
2518  if (strcmp(meshPartioningMethod, "parmetis") == 0)
2519  {
2520  mpParameters->Numerical().MeshPartitioning().set(cp::mesh_partitioning_type::parmetis);
2521  return;
2522  }
2523  if (strcmp(meshPartioningMethod, "petsc") == 0)
2524  {
2525  mpParameters->Numerical().MeshPartitioning().set(cp::mesh_partitioning_type::petsc);
2526  return;
2527  }
2528 
2529  EXCEPTION("Unknown mesh partitioning method provided");
2530 }
2531 
2532 void HeartConfig::SetApdMaps(const std::vector<std::pair<double, double> >& apdMaps)
2533 {
2535  XSD_SEQUENCE_TYPE(cp::postprocessing_type::ActionPotentialDurationMap)& apd_maps_sequence
2536  = mpParameters->PostProcessing()->ActionPotentialDurationMap();
2537  //Erase or create a sequence
2538  apd_maps_sequence.clear();
2539 
2540  for (unsigned i = 0; i < apdMaps.size(); i++)
2541  {
2542  XSD_CREATE_WITH_FIXED_ATTR2(cp::apd_map_type, temp,
2543  apdMaps[i].first, apdMaps[i].second,
2544  "mV");
2545  apd_maps_sequence.push_back(temp);
2546  }
2547 }
2548 
2549 void HeartConfig::SetUpstrokeTimeMaps(std::vector<double>& upstrokeTimeMaps)
2550 {
2552  XSD_SEQUENCE_TYPE(cp::postprocessing_type::UpstrokeTimeMap)& var_type_sequence
2553  = mpParameters->PostProcessing()->UpstrokeTimeMap();
2554 
2555  //Erase or create a sequence
2556  var_type_sequence.clear();
2557 
2558  for (unsigned i = 0; i < upstrokeTimeMaps.size(); i++)
2559  {
2560  XSD_CREATE_WITH_FIXED_ATTR1(cp::upstrokes_map_type, temp,
2561  upstrokeTimeMaps[i],
2562  "mV");
2563  var_type_sequence.push_back(temp);
2564  }
2565 }
2566 
2567 void HeartConfig::SetMaxUpstrokeVelocityMaps(std::vector<double>& maxUpstrokeVelocityMaps)
2568 {
2570  XSD_SEQUENCE_TYPE(cp::postprocessing_type::MaxUpstrokeVelocityMap)& max_upstroke_velocity_maps_sequence
2571  = mpParameters->PostProcessing()->MaxUpstrokeVelocityMap();
2572 
2573  //Erase or create a sequence
2574  max_upstroke_velocity_maps_sequence.clear();
2575 
2576  for (unsigned i = 0; i < maxUpstrokeVelocityMaps.size(); i++)
2577  {
2578  XSD_CREATE_WITH_FIXED_ATTR1(cp::max_upstrokes_velocity_map_type, temp,
2579  maxUpstrokeVelocityMaps[i],
2580  "mV");
2581 
2582  max_upstroke_velocity_maps_sequence.push_back(temp);
2583  }
2584 }
2585 
2586 void HeartConfig::SetConductionVelocityMaps(std::vector<unsigned>& conductionVelocityMaps)
2587 {
2589  XSD_SEQUENCE_TYPE(cp::postprocessing_type::ConductionVelocityMap)& conduction_velocity_maps_sequence
2590  = mpParameters->PostProcessing()->ConductionVelocityMap();
2591 
2592  //Erase or create a sequence
2593  conduction_velocity_maps_sequence.clear();
2594 
2595  for (unsigned i = 0; i < conductionVelocityMaps.size(); i++)
2596  {
2597  cp::conduction_velocity_map_type temp(conductionVelocityMaps[i]);
2598  conduction_velocity_maps_sequence.push_back(temp);
2599  }
2600 }
2601 
2602 void HeartConfig::SetRequestedNodalTimeTraces(std::vector<unsigned>& requestedNodes)
2603 {
2605  XSD_SEQUENCE_TYPE(cp::postprocessing_type::TimeTraceAtNode)& requested_nodes_sequence
2606  = mpParameters->PostProcessing()->TimeTraceAtNode();
2607 
2608  //Erase or create a sequence
2609  requested_nodes_sequence.clear();
2610 
2611  for (unsigned i = 0; i < requestedNodes.size(); i++)
2612  {
2613  cp::node_number_type temp(requestedNodes[i]);
2614  requested_nodes_sequence.push_back(temp);
2615  }
2616 }
2617 
2618 template <unsigned SPACE_DIM>
2619 void HeartConfig::SetPseudoEcgElectrodePositions(const std::vector<ChastePoint<SPACE_DIM> >& rPseudoEcgElectrodePositions)
2620 {
2622  XSD_SEQUENCE_TYPE(cp::postprocessing_type::PseudoEcgElectrodePosition)& electrodes_sequence
2623  = mpParameters->PostProcessing()->PseudoEcgElectrodePosition();
2624 
2625  //Erase or create a sequence
2626  electrodes_sequence.clear();
2627 
2628  for (unsigned i = 0; i < rPseudoEcgElectrodePositions.size(); i++)
2629  {
2630  cp::point_type temp(rPseudoEcgElectrodePositions[i].GetWithDefault(0),
2631  rPseudoEcgElectrodePositions[i].GetWithDefault(1),
2632  rPseudoEcgElectrodePositions[i].GetWithDefault(2));
2633  electrodes_sequence.push_back(temp);
2634  }
2635 }
2636 
2637 /*
2638  * Output visualizer
2639  */
2640 
2642 {
2643  ENSURE_SECTION_PRESENT(mpParameters->Simulation()->OutputVisualizer(), cp::output_visualizer_type);
2644 }
2645 
2647 {
2649 
2650  mpParameters->Simulation()->OutputVisualizer()->meshalyzer(
2651  useMeshalyzer ? cp::yesno_type::yes : cp::yesno_type::no);
2652 }
2653 
2655 {
2657 
2658  mpParameters->Simulation()->OutputVisualizer()->cmgui(
2659  useCmgui ? cp::yesno_type::yes : cp::yesno_type::no);
2660 }
2661 
2663 {
2665 
2666  mpParameters->Simulation()->OutputVisualizer()->vtk(
2667  useVtk ? cp::yesno_type::yes : cp::yesno_type::no);
2668 }
2669 
2671 {
2673 
2674  mpParameters->Simulation()->OutputVisualizer()->parallel_vtk(
2675  useParallelVtk ? cp::yesno_type::yes : cp::yesno_type::no);
2676 }
2677 
2678 void HeartConfig::SetVisualizerOutputPrecision(unsigned numberOfDigits)
2679 {
2681 
2682  mpParameters->Simulation()->OutputVisualizer()->precision(numberOfDigits);
2683 }
2684 
2685 void HeartConfig::SetElectrodeParameters(bool groundSecondElectrode,
2686  unsigned index, double magnitude,
2687  double startTime, double duration)
2688 {
2689  assert(index < 3);
2690 
2691  cp::axis_type axis = cp::axis_type::x;
2692  if (index == 1)
2693  {
2694  axis = cp::axis_type::y;
2695  }
2696  else if (index == 2)
2697  {
2698  axis = cp::axis_type::z;
2699  }
2700 
2701  XSD_CREATE_WITH_FIXED_ATTR1(cp::surface_stimulus_strength_type, strength, magnitude, "uA/cm^2");
2702  XSD_CREATE_WITH_FIXED_ATTR1(cp::time_type, start_time, startTime, "ms");
2703  XSD_CREATE_WITH_FIXED_ATTR1(cp::time_type, duration_time, duration, "ms");
2704 
2705  if (!IsElectrodesPresent())
2706  {
2707  cp::electrodes_type element(groundSecondElectrode ? cp::yesno_type::yes : cp::yesno_type::no,
2708  axis,
2709  strength,
2710  start_time,
2711  duration_time);
2712  mpParameters->Simulation()->Electrodes().set(element);
2713  }
2714  else
2715  {
2716  mpParameters->Simulation()->Electrodes()->GroundSecondElectrode(groundSecondElectrode ? cp::yesno_type::yes : cp::yesno_type::no);
2717  mpParameters->Simulation()->Electrodes()->PerpendicularToAxis(axis);
2718  mpParameters->Simulation()->Electrodes()->Strength(strength);
2719  mpParameters->Simulation()->Electrodes()->StartTime(start_time);
2720  mpParameters->Simulation()->Electrodes()->Duration(duration_time);
2721  }
2722 }
2723 
2724 void HeartConfig::GetElectrodeParameters(bool& rGroundSecondElectrode,
2725  unsigned& rIndex, double& rMagnitude,
2726  double& rStartTime, double& rDuration)
2727 {
2728  if (!IsElectrodesPresent())
2729  {
2730  EXCEPTION("Attempted to get electrodes that have not been defined.");
2731  }
2732  else
2733  {
2734  rGroundSecondElectrode = (mpParameters->Simulation()->Electrodes()->GroundSecondElectrode() == cp::yesno_type::yes);
2735 
2736  cp::axis_type axis = mpParameters->Simulation()->Electrodes()->PerpendicularToAxis();
2737  if (axis == cp::axis_type::x)
2738  {
2739  rIndex = 0;
2740  }
2741  else if (axis == cp::axis_type::y)
2742  {
2743  rIndex = 1;
2744  }
2745  else
2746  {
2747  rIndex = 2;
2748  }
2749 
2750  rMagnitude = mpParameters->Simulation()->Electrodes()->Strength();
2751  rStartTime = mpParameters->Simulation()->Electrodes()->StartTime();
2752  rDuration = mpParameters->Simulation()->Electrodes()->Duration();
2753  }
2754 }
2755 
2757 {
2758  // If it's an older version parameters & defaults (we're loading a checkpoint) say 'no'
2759  bool result = false;
2760  if (mpParameters->Numerical().UseStateVariableInterpolation().present())
2761  {
2762  result = mpParameters->Numerical().UseStateVariableInterpolation().get() == cp::yesno_type::yes;
2763  }
2764  return result;
2765 }
2766 
2767 void HeartConfig::SetUseStateVariableInterpolation(bool useStateVariableInterpolation)
2768 {
2769  if (useStateVariableInterpolation)
2770  {
2771  mpParameters->Numerical().UseStateVariableInterpolation().set(cp::yesno_type::yes);
2772  }
2773  else
2774  {
2775  mpParameters->Numerical().UseStateVariableInterpolation().set(cp::yesno_type::no);
2776  }
2777 }
2778 
2780 {
2781  return mpParameters->Physiological().ApplyDrug().present();
2782 }
2783 
2785 {
2786  CHECK_EXISTS(HasDrugDose(), "Physiological/ApplyDrug");
2787  return mpParameters->Physiological().ApplyDrug()->concentration();
2788 }
2789 
2790 void HeartConfig::SetDrugDose(double drugDose)
2791 {
2792  if (!mpParameters->Physiological().ApplyDrug().present())
2793  {
2794  cp::apply_drug_type drug(drugDose);
2795  mpParameters->Physiological().ApplyDrug().set(drug);
2796  }
2797  else
2798  {
2799  mpParameters->Physiological().ApplyDrug()->concentration(drugDose);
2800  }
2801 }
2802 
2803 std::map<std::string, std::pair<double, double> > HeartConfig::GetIc50Values()
2804 {
2805  CHECK_EXISTS(HasDrugDose(), "Physiological/ApplyDrug");
2806  std::map<std::string, std::pair<double, double> > ic50s;
2807 
2808  XSD_SEQUENCE_TYPE(cp::apply_drug_type::IC50)& ic50_seq = mpParameters->Physiological().ApplyDrug()->IC50();
2809 
2810  for (XSD_ITERATOR_TYPE(cp::apply_drug_type::IC50) i = ic50_seq.begin();
2811  i != ic50_seq.end();
2812  ++i)
2813  {
2814  std::pair<double, double> ic50_hill(*i, i->hill());
2815  std::string current = i->current();
2816  ic50s[current] = ic50_hill;
2817  }
2818 
2819  return ic50s;
2820 }
2821 
2822 void HeartConfig::SetIc50Value(const std::string& rCurrentName, double ic50, double hill)
2823 {
2824  if (!mpParameters->Physiological().ApplyDrug().present())
2825  {
2826  SetDrugDose(0.0);
2827  }
2828  XSD_SEQUENCE_TYPE(cp::apply_drug_type::IC50)& ic50_seq = mpParameters->Physiological().ApplyDrug()->IC50();
2829  if (ic50_seq.empty())
2830  {
2831  // Erase or create a sequence
2832  ic50_seq.clear();
2833  }
2834  bool entry_exists = false;
2835  cp::ic50_type ic50_elt(ic50, rCurrentName);
2836  ic50_elt.hill(hill);
2837  for (XSD_ITERATOR_TYPE(cp::apply_drug_type::IC50) i = ic50_seq.begin();
2838  i != ic50_seq.end();
2839  ++i)
2840  {
2841  if (i->current() == rCurrentName)
2842  {
2843  entry_exists = true;
2844  *i = ic50_elt;
2845  break;
2846  }
2847  }
2848  if (!entry_exists)
2849  {
2850  ic50_seq.push_back(ic50_elt);
2851  }
2852 }
2853 
2854 void HeartConfig::SetUseMassLumping(bool useMassLumping)
2855 {
2856  mUseMassLumping = useMassLumping;
2857 }
2858 
2860 {
2861  return mUseMassLumping;
2862 }
2863 
2865 {
2866  mUseMassLumpingForPrecond = useMassLumping;
2867 }
2868 
2870 {
2872 }
2873 
2875 {
2876  mUseReactionDiffusionOperatorSplitting = useOperatorSplitting;
2877 }
2878 
2880 {
2882 }
2883 
2884 void HeartConfig::SetUseFixedNumberIterationsLinearSolver(bool useFixedNumberIterations, unsigned evaluateNumItsEveryNSolves)
2885 {
2886  mUseFixedNumberIterations = useFixedNumberIterations;
2887  mEvaluateNumItsEveryNSolves = evaluateNumItsEveryNSolves;
2888 }
2889 
2891 {
2893 }
2894 
2896 {
2898 }
2899 
2900 //
2901 // Purkinje methods
2902 //
2903 
2905 {
2906  CheckSimulationIsDefined("Purkinje");
2907  return mpParameters->Simulation()->Purkinje().present();
2908 }
2909 
2911 {
2912  CHECK_EXISTS(mpParameters->Physiological().Purkinje().present(), "Physiological/Purkinje");
2913  CHECK_EXISTS(mpParameters->Physiological().Purkinje()->Capacitance().present(),
2914  "Physiological/Purkinje/Capacitance");
2915  return mpParameters->Physiological().Purkinje()->Capacitance().get();
2916 }
2917 
2918 void HeartConfig::SetPurkinjeCapacitance(double capacitance)
2919 {
2920  ENSURE_SECTION_PRESENT(mpParameters->Physiological().Purkinje(), cp::purkinje_physiological_type);
2921  XSD_CREATE_WITH_FIXED_ATTR1(cp::capacitance_type, purk_Cm, capacitance, "uF/cm^2");
2922  mpParameters->Physiological().Purkinje()->Capacitance().set(purk_Cm);
2923 }
2924 
2926 {
2927  CHECK_EXISTS(mpParameters->Physiological().Purkinje().present(), "Physiological/Purkinje");
2928  CHECK_EXISTS(mpParameters->Physiological().Purkinje()->SurfaceAreaToVolumeRatio().present(),
2929  "Physiological/Purkinje/SurfaceAreaToVolumeRatio");
2930  return mpParameters->Physiological().Purkinje()->SurfaceAreaToVolumeRatio().get();
2931 }
2932 
2934 {
2935  ENSURE_SECTION_PRESENT(mpParameters->Physiological().Purkinje(), cp::purkinje_physiological_type);
2936  XSD_CREATE_WITH_FIXED_ATTR1(cp::inverse_length_type, purk_Am, ratio, "1/cm");
2937  mpParameters->Physiological().Purkinje()->SurfaceAreaToVolumeRatio().set(purk_Am);
2938 }
2939 
2941 {
2942  CHECK_EXISTS(mpParameters->Physiological().Purkinje().present(), "Physiological/Purkinje");
2943  CHECK_EXISTS(mpParameters->Physiological().Purkinje()->Conductivity().present(),
2944  "Physiological/Purkinje/Conductivity");
2945  return mpParameters->Physiological().Purkinje()->Conductivity().get();
2946 }
2947 
2948 void HeartConfig::SetPurkinjeConductivity(double conductivity)
2949 {
2950  ENSURE_SECTION_PRESENT(mpParameters->Physiological().Purkinje(), cp::purkinje_physiological_type);
2951  XSD_CREATE_WITH_FIXED_ATTR1(cp::conductivity_type, purkinje_conductivity, conductivity, "mS/cm");
2952  mpParameters->Physiological().Purkinje()->Conductivity().set(purkinje_conductivity);
2953 }
2954 
2955 /**********************************************************************
2956  * *
2957  * *
2958  * Utility methods for reading/transforming XML *
2959  * *
2960  * *
2961  **********************************************************************/
2962 
2963 void XmlTransforms::TransformArchiveDirectory(xercesc::DOMDocument* pDocument,
2964  xercesc::DOMElement* pRootElement)
2965 {
2966  using namespace xercesc;
2967  std::vector<xercesc::DOMElement*> elts = XmlTools::FindElements(
2968  pRootElement,
2969  "ResumeSimulation/ArchiveDirectory");
2970  if (elts.size() > 0)
2971  {
2972  // We have an ArchiveDirectory element, so add the relative_to='chaste_test_output' attribute
2973  DOMElement* p_dir_elt = elts[0];
2974  p_dir_elt->setAttribute(X("relative_to"), X("chaste_test_output"));
2975  }
2976 }
2977 
2978 void XmlTransforms::TransformIonicModelDefinitions(xercesc::DOMDocument* pDocument,
2979  xercesc::DOMElement* pRootElement)
2980 {
2981  // Default ionic model
2982  std::vector<xercesc::DOMElement*> p_elt_list = XmlTools::FindElements(
2983  pRootElement,
2984  "Simulation/IonicModels/Default");
2985  if (p_elt_list.size() > 0)
2986  {
2987  assert(p_elt_list.size() == 1); // Asserted by schema
2988  XmlTools::WrapContentInElement(pDocument, p_elt_list[0], X("Hardcoded"));
2989  // Now do any region-specific definitions
2990  p_elt_list = XmlTools::FindElements(pRootElement, "Simulation/IonicModels/Region/IonicModel");
2991  for (unsigned i = 0; i < p_elt_list.size(); i++)
2992  {
2993  XmlTools::WrapContentInElement(pDocument, p_elt_list[i], X("Hardcoded"));
2994  }
2995  }
2996 }
2997 
2998 void XmlTransforms::CheckForIluPreconditioner(xercesc::DOMDocument* pDocument,
2999  xercesc::DOMElement* pRootElement)
3000 {
3001  std::vector<xercesc::DOMElement*> p_elt_list = XmlTools::FindElements(
3002  pRootElement,
3003  "Numerical/KSPPreconditioner");
3004  if (p_elt_list.size() > 0)
3005  {
3006  assert(p_elt_list.size() == 1); // Asserted by schema
3007  std::string text_value = X2C(p_elt_list[0]->getTextContent());
3008  if (text_value == "ilu")
3009  {
3010  EXCEPTION("PETSc does not have a parallel implementation of ilu, so we no longer allow it as an option. Use bjacobi instead.");
3011  }
3012  }
3013 }
3014 
3015 void XmlTransforms::MoveConductivityHeterogeneities(xercesc::DOMDocument* pDocument,
3016  xercesc::DOMElement* pRootElement)
3017 {
3018  std::vector<xercesc::DOMElement*> p_elt_list = XmlTools::FindElements(
3019  pRootElement,
3020  "Simulation/ConductivityHeterogeneities");
3021  if (p_elt_list.size() > 0)
3022  {
3023  assert(p_elt_list.size() == 1); // Asserted by schema
3024  xercesc::DOMNode* p_parent = p_elt_list[0]->getParentNode();
3025  xercesc::DOMNode* p_child = p_parent->removeChild(p_elt_list[0]);
3026  std::vector<xercesc::DOMElement*> p_phys_list = XmlTools::FindElements(pRootElement, "Physiological");
3027  assert(p_phys_list.size() == 1); // Asserted by schema
3028  p_phys_list[0]->appendChild(p_child);
3029  }
3030 }
3031 
3032 void XmlTransforms::SetDefaultVisualizer(xercesc::DOMDocument* pDocument,
3033  xercesc::DOMElement* pRootElement)
3034 {
3035  std::vector<xercesc::DOMElement*> p_sim_list = XmlTools::FindElements(pRootElement, "Simulation");
3036  if (p_sim_list.size() > 0)
3037  {
3038  std::vector<xercesc::DOMElement*> p_viz_list = XmlTools::FindElements(p_sim_list[0], "OutputVisualizer");
3039  if (p_viz_list.empty())
3040  {
3041  // Create the element and set meshalyzer (only) to on
3042  xercesc::DOMElement* p_viz_elt = pDocument->createElementNS(X("https://chaste.comlab.ox.ac.uk/nss/parameters/3_3"), X("OutputVisualizer"));
3043  p_sim_list[0]->appendChild(p_viz_elt);
3044  p_viz_elt->setAttribute(X("meshalyzer"), X("yes"));
3045  }
3046  }
3047 }
3048 
3050 // Explicit instantiation of the templated functions
3052 // LCOV_EXCL_START //These methods are covered above with DIM=1,2,3 but the instantiations may fail spuriously
3057 template void HeartConfig::GetIonicModelRegions<3u>(std::vector<boost::shared_ptr<AbstractChasteRegion<3u> > >&, std::vector<cp::ionic_model_selection_type>&) const;
3058 template void HeartConfig::GetStimuli<3u>(std::vector<boost::shared_ptr<AbstractStimulusFunction> >&, std::vector<boost::shared_ptr<AbstractChasteRegion<3u> > >&) const;
3059 template void HeartConfig::GetCellHeterogeneities<3u>(std::vector<boost::shared_ptr<AbstractChasteRegion<3u> > >&, std::vector<double>&, std::vector<double>&, std::vector<double>&, std::vector<std::map<std::string, double> >*);
3060 template void HeartConfig::GetConductivityHeterogeneities<3u>(std::vector<boost::shared_ptr<AbstractChasteRegion<3u> > >&, std::vector<c_vector<double, 3> >&, std::vector<c_vector<double, 3> >&) const;
3061 
3062 template void HeartConfig::GetIonicModelRegions<2u>(std::vector<boost::shared_ptr<AbstractChasteRegion<2u> > >&, std::vector<cp::ionic_model_selection_type>&) const;
3063 template void HeartConfig::GetStimuli<2u>(std::vector<boost::shared_ptr<AbstractStimulusFunction> >&, std::vector<boost::shared_ptr<AbstractChasteRegion<2u> > >&) const;
3064 template void HeartConfig::GetCellHeterogeneities<2u>(std::vector<boost::shared_ptr<AbstractChasteRegion<2u> > >&, std::vector<double>&, std::vector<double>&, std::vector<double>&, std::vector<std::map<std::string, double> >*);
3065 template void HeartConfig::GetConductivityHeterogeneities<2u>(std::vector<boost::shared_ptr<AbstractChasteRegion<2u> > >&, std::vector<c_vector<double, 3> >&, std::vector<c_vector<double, 3> >&) const;
3066 
3067 template void HeartConfig::GetIonicModelRegions<1u>(std::vector<boost::shared_ptr<AbstractChasteRegion<1u> > >&, std::vector<cp::ionic_model_selection_type>&) const;
3068 template void HeartConfig::GetStimuli<1u>(std::vector<boost::shared_ptr<AbstractStimulusFunction> >&, std::vector<boost::shared_ptr<AbstractChasteRegion<1u> > >&) const;
3069 template void HeartConfig::GetCellHeterogeneities<1u>(std::vector<boost::shared_ptr<AbstractChasteRegion<1u> > >&, std::vector<double>&, std::vector<double>&, std::vector<double>&, std::vector<std::map<std::string, double> >*);
3070 template void HeartConfig::GetConductivityHeterogeneities<1u>(std::vector<boost::shared_ptr<AbstractChasteRegion<1u> > >&, std::vector<c_vector<double, 3> >&, std::vector<c_vector<double, 3> >&) const;
3071 
3072 template void HeartConfig::GetPseudoEcgElectrodePositions(std::vector<ChastePoint<1u> >& rPseudoEcgElectrodePositions) const;
3073 template void HeartConfig::GetPseudoEcgElectrodePositions(std::vector<ChastePoint<2u> >& rPseudoEcgElectrodePositions) const;
3074 template void HeartConfig::GetPseudoEcgElectrodePositions(std::vector<ChastePoint<3u> >& rPseudoEcgElectrodePositions) const;
3075 
3076 template void HeartConfig::SetPseudoEcgElectrodePositions(const std::vector<ChastePoint<1u> >& rPseudoEcgElectrodePositions);
3077 template void HeartConfig::SetPseudoEcgElectrodePositions(const std::vector<ChastePoint<2u> >& rPseudoEcgElectrodePositions);
3078 template void HeartConfig::SetPseudoEcgElectrodePositions(const std::vector<ChastePoint<3u> >& rPseudoEcgElectrodePositions);
3083 // LCOV_EXCL_STOP //These methods are covered above with DIM=1,2,3 but the instantiations may fail spuriously
3084 
3085 // Serialization for Boost >= 1.36
static std::string EscapeSpaces(const std::string &rPath)
Definition: XmlTools.cpp:442
void UpdateParametersFromResumeSimulation(boost::shared_ptr< cp::chaste_parameters_type > pResumeParameters)
double GetBathConductivity(unsigned bathRegion=UINT_MAX) const
void EnsureOutputVisualizerExists(void)
void SetIc50Value(const std::string &rCurrentName, double ic50, double hill=1.0)
void SetMeshFileName(std::string meshPrefix, cp::media_type fibreDefinition=cp::media_type::NoFibreOrientation)
void SetApdMaps(const std::vector< std::pair< double, double > > &rApdMaps)
void SetUseRelativeTolerance(double relativeTolerance)
std::map< std::string, std::string > SchemaLocationsMap
std::string GetOutputFilenamePrefix() const
bool GetUseAbsoluteTolerance() const
SchemaLocationsMap mSchemaLocations
void SetDomain(const cp::domain_type &rDomain)
double GetPurkinjeSurfaceAreaToVolumeRatio()
cp::media_type GetConductivityMedia() const
unsigned mEvaluateNumItsEveryNSolves
void SetRequestedNodalTimeTraces(std::vector< unsigned > &requestedNodes)
double GetSimulationDuration() const
void SetSpaceDimension(unsigned spaceDimension)
static bool IsRegionBath(HeartRegionType regionId)
bool GetCreateSheet() const
void GetIonicModelRegions(std::vector< boost::shared_ptr< AbstractChasteRegion< DIM > > > &rDefinedRegions, std::vector< cp::ionic_model_selection_type > &rIonicModels) const
const std::set< unsigned > & rGetTissueIdentifiers()
bool GetVisualizeWithVtk() const
bool IsAnyNodalTimeTraceRequested() const
void SetVisualizeWithMeshalyzer(bool useMeshalyzer=true)
void SetOutputDirectory(const std::string &rOutputDirectory)
void SetUseFixedNumberIterationsLinearSolver(bool useFixedNumberIterations=true, unsigned evaluateNumItsEveryNSolves=UINT_MAX)
bool IsConductionVelocityMapsRequested() const
void SetPseudoEcgElectrodePositions(const std::vector< ChastePoint< SPACE_DIM > > &rPseudoEcgElectrodePositions)
void SetSlabDimensions(double x, double y, double z, double inter_node_space)
unsigned GetVersionFromNamespace(const std::string &rNamespaceUri)
void SetBathMultipleConductivities(std::map< unsigned, double > bathConductivities)
unsigned GetEvaluateNumItsEveryNSolves()
void SetSheetDimensions(double x, double y, double inter_node_space)
bool HasPurkinje()
void GetCellHeterogeneities(std::vector< boost::shared_ptr< AbstractChasteRegion< DIM > > > &rCellHeterogeneityRegions, std::vector< double > &rScaleFactorGks, std::vector< double > &rScaleFactorIto, std::vector< double > &rScaleFactorGkr, std::vector< std::map< std::string, double > > *pParameterSettings)
void SetDefaultIonicModel(const cp::ionic_models_available_type &rIonicModel)
void GetIntracellularConductivities(c_vector< double, 3 > &rIntraConductivities) const
void GetExtracellularConductivities(c_vector< double, 3 > &rExtraConductivities) const
void CheckResumeSimulationIsDefined(std::string callingMethod="") const
void CheckSimulationIsDefined(std::string callingMethod="") const
bool GetUseMassLumpingForPrecond()
static xercesc::DOMElement * SetNamespace(xercesc::DOMDocument *pDocument, xercesc::DOMElement *pElement, const std::string &rNamespace)
Definition: XmlTools.cpp:335
unsigned mIndexMid
void GetConductionVelocityMaps(std::vector< unsigned > &rConductionVelocityMaps) const
void Write(bool useArchiveLocationInfo=false, std::string subfolderName="output")
void GetSlabDimensions(c_vector< double, 3 > &slabDimensions) const
void SetPdeTimeStep(double pdeTimeStep)
double GetCheckpointTimestep() const
bool mUseMassLumpingForPrecond
double mMidFraction
#define EXCEPTION(message)
Definition: Exception.hpp:143
boost::shared_ptr< cp::chaste_parameters_type > mpParameters
std::map< unsigned, double > mBathConductivities
bool mUseFixedSchemaLocation
void SetTissueAndBathIdentifiers(const std::set< unsigned > &rTissueIds, const std::set< unsigned > &rBathIds)
void EnsurePostProcessingSectionPresent()
bool Divides(double smallerNumber, double largerNumber)
static std::vector< xercesc::DOMElement * > FindElements(const xercesc::DOMElement *pContextElement, const std::string &rPath)
Definition: XmlTools.cpp:384
bool IsSimulationDefined() const
void SetVisualizeWithCmgui(bool useCmgui=true)
static void TransformArchiveDirectory(xercesc::DOMDocument *pDocument, xercesc::DOMElement *pRootElement)
void SetUseMassLumpingForPrecond(bool useMassLumping=true)
void SetUseStateVariableInterpolation(bool useStateVariableInterpolation=true)
bool GetOutputUsingOriginalNodeOrdering()
void SetKSPPreconditioner(const char *kspPreconditioner)
static bool AmMaster()
Definition: PetscTools.cpp:120
FileFinder CopyTo(const FileFinder &rDest) const
Definition: FileFinder.cpp:308
bool IsPseudoEcgCalculationRequested() const
unsigned mIndexEpi
unsigned GetEndoLayerIndex()
void SetElectrodeParameters(bool groundSecondElectrode, unsigned index, double magnitude, double startTime, double duration)
bool IsMaxUpstrokeVelocityMapRequested() const
void SetOutputUsingOriginalNodeOrdering(bool useOriginal)
bool GetVisualizeWithCmgui() const
void GetOutputVariables(std::vector< std::string > &rOutputVariables) const
DistributedTetrahedralMeshPartitionType::type GetMeshPartitioning() const
bool GetCreateFibre() const
double GetPrintingTimeStep() const
void SetVisualizeWithParallelVtk(bool useParallelVtk=true)
std::set< unsigned > mTissueIdentifiers
void SetPurkinjeSurfaceAreaToVolumeRatio(double ratio)
bool GetOutputVariablesProvided() const
void SetUseReactionDiffusionOperatorSplitting(bool useOperatorSplitting=true)
double GetAbsoluteTolerance() const
void SetDrugDose(double drugDose)
void SetCheckpointSimulation(bool checkpointSimulation, double checkpointTimestep=-1.0, unsigned maxCheckpointsOnDisk=UINT_MAX)
void SetUseFixedSchemaLocation(bool useFixedSchemaLocation)
std::string GetOutputDirectoryFullPath() const
double GetPdeTimeStep() const
void GetUpstrokeTimeMaps(std::vector< double > &rUpstrokeTimeMaps) const
bool IsApdMapsRequested() const
bool IsPostProcessingSectionPresent() const
#define NEVER_REACHED
Definition: Exception.hpp:206
void SetSurfaceAreaToVolumeRatio(double ratio)
double GetRelativeTolerance() const
bool mUseReactionDiffusionOperatorSplitting
std::string GetMeshName() const
void GetConductivityHeterogeneities(std::vector< boost::shared_ptr< AbstractChasteRegion< DIM > > > &conductivitiesHeterogeneityAreas, std::vector< c_vector< double, 3 > > &intraConductivities, std::vector< c_vector< double, 3 > > &extraConductivities) const
void GetFibreLength(c_vector< double, 1 > &fibreLength) const
bool GetUseFixedNumberIterationsLinearSolver()
bool IsAdaptivityParametersPresent() const
void MergeDefaults(boost::shared_ptr< cp::chaste_parameters_type > pParams, boost::shared_ptr< cp::chaste_parameters_type > pDefaults)
boost::shared_ptr< cp::chaste_parameters_type > CreateDefaultParameters()
void SetSimulationDuration(double simulationDuration)
void SetOdeTimeStep(double odeTimeStep)
cp::domain_type GetDomain() const
double GetEpiLayerFraction()
void SetExtracellularConductivities(const c_vector< double, 3 > &rExtraConductivities)
bool GetCreateSlab() const
static void TransformIonicModelDefinitions(xercesc::DOMDocument *pDocument, xercesc::DOMElement *pRootElement)
bool mUserAskedForCellularTransmuralHeterogeneities
static void WrapContentInElement(xercesc::DOMDocument *pDocument, xercesc::DOMElement *pElement, const XMLCh *pNewElementLocalName)
Definition: XmlTools.cpp:409
double GetDrugDose() const
void GetNodalTimeTraceRequested(std::vector< unsigned > &rRequestedNodes) const
void GetSheetDimensions(c_vector< double, 2 > &sheetDimensions) const
void SetConductivityHeterogeneities(std::vector< ChasteCuboid< 3 > > &rConductivityAreas, std::vector< c_vector< double, 3 > > &rIntraConductivities, std::vector< c_vector< double, 3 > > &rExtraConductivities)
double GetMidLayerFraction()
void SetUseMassLumping(bool useMassLumping=true)
bool GetUseStateVariableInterpolation() const
static void MoveConductivityHeterogeneities(xercesc::DOMDocument *pDocument, xercesc::DOMElement *pRootElement)
std::set< unsigned > mBathIdentifiers
bool IsPostProcessingRequested() const
unsigned GetVisualizerOutputPrecision()
void LoadFromCheckpoint()
void SetConductionVelocityMaps(std::vector< unsigned > &rConductionVelocityMaps)
void SetFibreLength(double x, double inter_node_space)
void CheckTimeSteps() const
static void CheckForIluPreconditioner(xercesc::DOMDocument *pDocument, xercesc::DOMElement *pRootElement)
void SetCapacitance(double capacitance)
double GetOdeTimeStep() const
unsigned GetSpaceDimension() const
void GetElectrodeParameters(bool &rGroundSecondElectrode, unsigned &rIndex, double &rMagnitude, double &rStartTime, double &rDuration)
unsigned GetEpiLayerIndex()
void GetPseudoEcgElectrodePositions(std::vector< ChastePoint< SPACE_DIM > > &rPseudoEcgElectrodePositions) const
bool GetUseReactionDiffusionOperatorSplitting()
void SetDefaultSchemaLocations()
bool IsElectrodesPresent() const
FileFinder GetParametersFilePath()
void SetMaxUpstrokeVelocityMaps(std::vector< double > &rMaxUpstrokeVelocityMaps)
bool GetUseMassLumping()
double mEndoFraction
c_vector< double, 1 > Create_c_vector(double x)
cp::ionic_model_selection_type GetDefaultIonicModel() const
void SetUpstrokeTimeMaps(std::vector< double > &rUpstrokeTimeMaps)
void SetMeshPartitioning(const char *meshPartioningMethod)
bool IsMeshProvided() const
static xsd::cxx::xml::dom::auto_ptr< xercesc::DOMDocument > ReadXmlFile(const std::string &rFileName, const ::xsd::cxx::tree::properties< char > &rProps, bool validate=true)
Definition: XmlTools.cpp:53
unsigned GetMaxCheckpointsOnDisk() const
boost::shared_ptr< cp::chaste_parameters_type > ReadFile(const std::string &rFileName)
void SetFixedSchemaLocations(const SchemaLocationsMap &rSchemaLocations)
double GetInterNodeSpace() const
double GetPurkinjeConductivity()
bool HasDrugDose() const
bool Exists() const
Definition: FileFinder.cpp:183
double GetEndoLayerFraction()
unsigned GetMidLayerIndex()
void GetApdMaps(std::vector< std::pair< double, double > > &rApdMaps) const
void SetParametersFile(const std::string &rFileName)
static void SetDefaultVisualizer(xercesc::DOMDocument *pDocument, xercesc::DOMElement *pRootElement)
bool GetCreateMesh() const
void SetIonicModelRegions(std::vector< ChasteCuboid< 3 > > &rDefinedRegions, std::vector< cp::ionic_model_selection_type > &rIonicModels) const
bool GetCheckpointSimulation() const
bool GetVisualizeWithParallelVtk() const
double GetPurkinjeCapacitance()
void SetKSPSolver(const char *kspSolver, bool warnOfChange=false)
bool GetLoadMesh() const
bool mUseFixedNumberIterations
const std::set< unsigned > & rGetBathIdentifiers()
double GetCapacitance() const
bool IsOutputVisualizerPresent() const
#define CHASTE_CLASS_EXPORT(T)
void SetIntracellularConductivities(const c_vector< double, 3 > &rIntraConductivities)
bool GetVisualizeWithMeshalyzer() const
FileFinder mParametersFilePath
void SetPurkinjeCapacitance(double capacitance)
double GetSurfaceAreaToVolumeRatio() const
std::string GetOutputDirectory() const
void SetOutputVariables(const std::vector< std::string > &rOutputVariables)
const char * GetKSPPreconditioner() const
bool IsSimulationResumed() const
virtual void SetPath(const std::string &rPath, RelativeTo::Value relativeTo)
Definition: FileFinder.cpp:102
void GetMaxUpstrokeVelocityMaps(std::vector< double > &rUpstrokeVelocityMaps) const
void SetVisualizerOutputPrecision(unsigned numberOfDigits)
static std::string GetArchiveDirectory()
HeartFileFinder GetArchivedSimulationDir() const
static void Reset()
void SetPrintingTimeStep(double printingTimeStep)
unsigned mIndexEndo
void SetPurkinjeConductivity(double conductivity)
bool AreCellularTransmuralHeterogeneitiesRequested()
void GetStimuli(std::vector< boost::shared_ptr< AbstractStimulusFunction > > &rStimuliApplied, std::vector< boost::shared_ptr< AbstractChasteRegion< DIM > > > &rStimulatedAreas) const
void SetUseAbsoluteTolerance(double absoluteTolerance)
void SetOutputFilenamePrefix(const std::string &rOutputFilenamePrefix)
bool mUseMassLumping
static HeartConfig * Instance()
static boost::shared_ptr< HeartConfig > mpInstance
static const char * GetRootDir()
bool GetUseRelativeTolerance() const
void SetBathConductivity(double bathConductivity)
bool IsUpstrokeTimeMapsRequested() const
void SetVisualizeWithVtk(bool useVtk=true)
const char * GetKSPSolver() const
bool GetConductivityHeterogeneitiesProvided() const
double mEpiFraction
void CopySchema(const std::string &rToDirectory)
void SetOdePdeAndPrintingTimeSteps(double odeTimeStep, double pdeTimeStep, double printingTimeStep)
void SetConductivityHeterogeneitiesEllipsoid(std::vector< ChasteEllipsoid< 3 > > &rConductivityAreas, std::vector< c_vector< double, 3 > > &rIntraConductivities, std::vector< c_vector< double, 3 > > &rExtraConductivities)
std::map< std::string, std::pair< double, double > > GetIc50Values()