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