From: David Schneider Date: Wed, 2 Dec 2020 20:52:59 +0000 (+0100) Subject: Add a deal.II-preCICE example X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=81dd74d20978498d5b4c471f08be0f0b1d634199;p=code-gallery.git Add a deal.II-preCICE example --- diff --git a/.github/workflows/precice-ci.yml b/.github/workflows/precice-ci.yml new file mode 100644 index 0000000..6b2e20a --- /dev/null +++ b/.github/workflows/precice-ci.yml @@ -0,0 +1,35 @@ +name: PreCICE CI + +on: + push: + branches: [ master ] + pull_request: + branches: [ master ] + +jobs: + build: + + runs-on: ubuntu-latest + + steps: + - uses: actions/checkout@v2 + - name: setup + run: | + command="sudo apt-get -y update && \ + wget https://github.com/precice/precice/releases/download/v2.1.1/libprecice2_2.1.1_focal.deb && \ + sudo apt-get -y install ./libprecice2_2.1.1_focal.deb && \ + git clone https://github.com/${{ github.repository }} && \ + cd code-gallery && \ + git fetch origin ${{ github.ref }} && \ + git checkout FETCH_HEAD && \ + cd coupled_laplace_problem && \ + cmake . && \ + make && \ + (./coupled_laplace_problem 2>&1 & ./fancy_boundary_condition >fbc.log) && \ + sed -i '2d' solution-10.vtk && \ + numdiff solution-10.vtk test_data/reference-10.vtk"; + + echo $command + + docker pull dealii/dealii:v9.2.0-focal + docker run -t dealii/dealii:v9.2.0-focal /bin/sh -c "$command"; diff --git a/coupled_laplace_problem/Allclean b/coupled_laplace_problem/Allclean new file mode 100755 index 0000000..ad56970 --- /dev/null +++ b/coupled_laplace_problem/Allclean @@ -0,0 +1,25 @@ +#!/bin/sh +cd ${0%/*} || exit 1 # Run from this directory + +echo "Cleaning..." + +# Participant 1: coupled_laplace_problem +Participant1="coupled_laplace_problem" + +# Participant 2: fancy_boundary_condition +Participant2="fancy_boundary_condition" + +# Remove vtk result files +rm -fv solution-*.vtk + +# Remove the preCICE-related log files +echo "Deleting the preCICE log files..." +rm -fv \ + precice-*.log \ + precice-*-events.json + +rm -rfv precice-run +rm -fv .${Participant1}-${Participant2}.address + +echo "Cleaning complete!" +#------------------------------------------------------------------------------ diff --git a/coupled_laplace_problem/CMakeLists.txt b/coupled_laplace_problem/CMakeLists.txt new file mode 100644 index 0000000..ec2dc24 --- /dev/null +++ b/coupled_laplace_problem/CMakeLists.txt @@ -0,0 +1,51 @@ +## +# CMake script for the coupled laplace problem tutorial program: +## + +# Set the name of the project and target: +SET(TARGET "coupled_laplace_problem") + +# Declare all source files the target consists of. Here, this is only +# the one step-X.cc file, but as you expand your project you may wish +# to add other source files as well. If your project becomes much larger, +# you may want to either replace the following statement by something like +# FILE(GLOB_RECURSE TARGET_SRC "source/*.cc") +# FILE(GLOB_RECURSE TARGET_INC "include/*.h") +# SET(TARGET_SRC ${TARGET_SRC} ${TARGET_INC}) +# or switch altogether to the large project CMakeLists.txt file discussed +# in the "CMake in user projects" page accessible from the "User info" +# page of the documentation. +SET(TARGET_SRC + ${TARGET}.cc + ) + +# Usually, you will not need to modify anything beyond this point... + +CMAKE_MINIMUM_REQUIRED(VERSION 2.8.12) + +FIND_PACKAGE(deal.II 9.2.0 + HINTS ${deal.II_DIR} ${DEAL_II_DIR} ../ ../../ $ENV{DEAL_II_DIR} + ) +IF(NOT ${deal.II_FOUND}) + MESSAGE(FATAL_ERROR "\n" + "*** Could not locate a (sufficiently recent) version of deal.II. ***\n\n" + "You may want to either pass a flag -DDEAL_II_DIR=/path/to/deal.II to cmake\n" + "or set an environment variable \"DEAL_II_DIR\" that contains this path." + ) +ENDIF() + +DEAL_II_INITIALIZE_CACHED_VARIABLES() +PROJECT(${TARGET}) +DEAL_II_INVOKE_AUTOPILOT() + +FIND_PACKAGE(precice REQUIRED + HINTS ${PRECICE_DIR} + ) +MESSAGE(STATUS "Using the preCICE version found at ${precice_CONFIG}") + +TARGET_LINK_LIBRARIES(${TARGET} precice::precice) + +# In order to build the fancy_boundary_condition, the DEAL_II macros are not +# required since the code doesn't use deal.II at all. +ADD_EXECUTABLE(fancy_boundary_condition fancy_boundary_condition.cc) +TARGET_LINK_LIBRARIES(fancy_boundary_condition precice::precice) diff --git a/coupled_laplace_problem/README.md b/coupled_laplace_problem/README.md new file mode 100644 index 0000000..0f60b76 --- /dev/null +++ b/coupled_laplace_problem/README.md @@ -0,0 +1,69 @@ +Laplace equation coupled to an external simulation program +------------------------------------------ +## Overview + +preCICE allows to couple deal.II to external simulation software, such as OpenFOAM, SU2, or CalculiX. To keep dependencies of this example minimal, we couple deal.II to an external C++ program, which provides a time varying boundary condition. The deal.II code consists mainly of the [`step-4` tutorial program](https://www.dealii.org/developer/doxygen/deal.II/step_4.html), where a simple Laplace problem is solved. + +Coupling with preCICE is usually carried out along surfaces in order to apply a Dirichlet-Neumann coupling between two domains (volume coupling is also possible). For the sake of simplicity, we couple here an external C++ program in a unidirectional fashion to one side of our quadrilateral domain. The external C++ program generates a parabolic boundary profile with time varying amplitude. The boundary values are then used in the Laplace solver as a Dirichlet boundary condition. + +## Time discretization +Coupled simulations deal mostly with time-dependent problems. Hence, we make the stationary Laplace problem from step-4 time dependent. +@f{align*} + \\frac{\partial u}{\partial t}-\Delta u &= f \qquad\qquad & \text{in}\ \Omega, + \\ + u &= x^2+y^2 \qquad\qquad & \text{on}\ \partial\Omega_s, + \\ + u &= g(t) \qquad\qquad & \text{on}\ \partial\Omega_c. +@f} +with the fixed Dirichlet boundary \Omega_s, the coupling boundary \Omega_c and the time-dependent coupling data g(t). + +The system is consequently discretized by a first-order backward Euler method, resulting in +@f{align*} + \\frac{u^{n+1}-u^n}{\Delta t}-\Delta u^{n+1} &= f \qquad\qquad & \text{in}\ \Omega, +@f} +at the next time level n+1. + +## Requirements + +* `deal.II`, version `9.2` or greater. Older versions might work as well, but have not been tested. + +* [preCICE](https://github.com/precice/precice/wiki#1-get-precice), version `2.0` or greater. Have a look at the provided link for an installation guide. + +## Compiling and running + +Similar to the example programs, run +``` +cmake -DDEAL_II_DIR=/path/to/deal.II -Dprecice_DIR=/path/to/precice . +``` +in this directory to configure the problem. The explicit specification of the library locations can be omitted if they are installed globally or via environment variables. +You can then switch between debug and release mode by calling either +``` +make debug +``` +or +``` +make release +``` +This command will generate two executables: one for the `coupled_laplace_problem` and one for the `fancy_boundary_condition` participant. +``` +./coupled_laplace_problem +``` +executes the `coupled_laplace_problem`. In order to start the coupled simulation, execute +``` +./fancy_boundary_condition +``` +in the same directory from another terminal window. + +preCICE and the deal.II solver create several log files during the simulation. In order to remove all result-associated files and clean-up the simulation directory the Allclean script can be executed. + +## Results + +By default, the results are saved in each time step using `vtk` files, which can, for example, be visualized with ParaView. You should get something like: +![result](./doc/result.gif) + +## Further reading + +* A complete overview of the preCICE project can be found on the [preCICE webpage](https://www.precice.org/). +* The [source code](https://github.com/precice/precice/) of preCICE is hosted on Github. +* For more coupled deal.II codes, have a look at the [deal.II adapter repository](https://github.com/precice/dealii-adapter). In the [preCICE tutorials](https://github.com/precice/tutorials/tree/master/FSI/flap_perp_2D/OpenFOAM-deal.II), a fluid-structure interaction example is given coupling OpenFOAM (fluid dynamics) to deal.II (solid mechanics). +* [The preCICE reference paper](https://www.sciencedirect.com/science/article/abs/pii/S0045793016300974) diff --git a/coupled_laplace_problem/coupled_laplace_problem.cc b/coupled_laplace_problem/coupled_laplace_problem.cc new file mode 100644 index 0000000..072e67a --- /dev/null +++ b/coupled_laplace_problem/coupled_laplace_problem.cc @@ -0,0 +1,646 @@ +// The included deal.II header files are the same as in the other example +// programs: +#include +#include +#include + +#include +#include +#include + +#include +#include + +#include +#include +#include +#include + +#include +#include +#include +#include +#include +#include + +#include +#include +#include +// In addition to the deal.II header files, we include the preCICE API in order +// to obtain access to preCICE specific functionality +#include + +#include +#include + +using namespace dealii; + +// Configuration parameters +// +// We set up a simple hard-coded struct containing all names we need for +// external coupling. The struct includes the name of the preCICE +// configuration file as well as the name of the simulation participant, the +// name of the coupling mesh and the name of the exchanged data. The last three +// names you also find in the preCICE configuration file. For real application +// cases, these names are better handled by a parameter file. +struct CouplingParamters +{ + const std::string config_file = "precice-config.xml"; + const std::string participant_name = "laplace-solver"; + const std::string mesh_name = "dealii-mesh"; + const std::string read_data_name = "boundary-data"; +}; + + +// The Adapter class +// +// The Adapter class handles all functionalities to couple the deal.II solver +// code to other solvers with preCICE, i.e., data structures are set up and all +// relevant information is passed to preCICE. + +template +class Adapter +{ +public: + Adapter(const ParameterClass & parameters, + const types::boundary_id dealii_boundary_interface_id); + + double + initialize(const DoFHandler & dof_handler, + std::map &boundary_data, + const MappingQGeneric & mapping); + + double + advance(std::map &boundary_data, + const double computed_timestep_length); + + // public precCICE solver interface + precice::SolverInterface precice; + + // Boundary ID of the deal.II triangulation, associated with the coupling + // interface. The variable is defined in the constructor of this class and + // intentionally public so that it can be used during the grid generation and + // system assembly. The only thing, one needs to make sure is that this ID is + // unique for a particular triangulation. + const unsigned int dealii_boundary_interface_id; + +private: + // preCICE related initializations + // These variables are specified in and read from a parameter file, which is + // in this simple tutorial program the CouplingParameter struct already + // introduced in the beginning. + const std::string mesh_name; + const std::string read_data_name; + + // These IDs are filled by preCICE during the initialization. We set a default + // value of -1 in order to detect potential errors more easily. + int mesh_id; + int read_data_id; + int n_interface_nodes; + + // DoF IndexSet, containing relevant coupling DoF indices at the coupling + // boundary + IndexSet coupling_dofs; + + // Data containers which are passed to preCICE in an appropriate preCICE + // specific format + std::vector interface_nodes_ids; + std::vector read_data; + + // The MPI rank and total number of MPI ranks is required by preCICE when the + // SolverInterface is created. Since this tutorial runs only in serial mode we + // define the variables manually in this class instead of using the regular + // MPI interface. + static constexpr int this_mpi_process = 0; + static constexpr int n_mpi_processes = 1; + + // Function to transform the obtained data from preCICE into an appropriate + // map for Dirichlet boundary conditions + void + format_precice_to_dealii( + std::map &boundary_data) const; +}; + + + +// In the constructor of the Adapter class, we set up the preCICE +// SolverInterface. We need to tell preCICE our name as participant of the +// simulation and the name of the preCICE configuration file. Both have already +// been specified in the CouplingParameter class above. Thus, we pass the class +// directly to the constructor and read out all relevant information. As a +// second parameter, we need to specify the boundary ID of our triangulation, +// which is associated with the coupling interface. +template +Adapter::Adapter( + const ParameterClass & parameters, + const types::boundary_id deal_boundary_interface_id) + : precice(parameters.participant_name, + parameters.config_file, + this_mpi_process, + n_mpi_processes) + , dealii_boundary_interface_id(deal_boundary_interface_id) + , mesh_name(parameters.mesh_name) + , read_data_name(parameters.read_data_name) +{} + + + +// This function initializes preCICE (e.g. establishes communication channels +// and allocates memory) and passes all relevant data to preCICE. For surface +// coupling, relevant data is in particular the location of the data points at +// the associated interface(s). The `boundary_data` is an empty map, which is +// filled by preCICE, i.e., information of the other participant. Throughout +// the system assembly, the map can directly be used in order to apply the +// Dirichlet boundary conditions in the linear system. preCICE returns the +// maximum admissible time-step size during the initialization. +template +double +Adapter::initialize( + const DoFHandler & dof_handler, + std::map &boundary_data, + const MappingQGeneric & mapping) +{ + Assert(dim > 1, ExcNotImplemented()); + AssertDimension(dim, precice.getDimensions()); + + // In a first step, we get preCICE specific IDs from preCICE and store them in + // the respective variables. Later, they are used for data transfer. + mesh_id = precice.getMeshID(mesh_name); + read_data_id = precice.getDataID(read_data_name, mesh_id); + + + // Afterwards, we extract the number of interface nodes and the coupling DoFs + // at the coupling interface from our deal.II solver via + // `extract_boundary_dofs()` + std::set couplingBoundary; + couplingBoundary.insert(dealii_boundary_interface_id); + + // The `ComponentMask()` might be important in case we deal with vector valued + // problems, because vector valued problems have a DoF for each component. + DoFTools::extract_boundary_dofs(dof_handler, + ComponentMask(), + coupling_dofs, + couplingBoundary); + + // The coupling DoFs are used to set up the `boundary_data` map. At the end, + // we associate here each DoF with a respective boundary value. + for (const auto i : coupling_dofs) + boundary_data[i] = 0.0; + + // Since we deal with a scalar problem, the number of DoFs at the particular + // interface corresponds to the number of interface nodes. + n_interface_nodes = coupling_dofs.n_elements(); + + std::cout << "\t Number of coupling nodes: " << n_interface_nodes + << std::endl; + + // Now, we need to tell preCICE the coordinates of the interface nodes. Hence, + // we set up a std::vector to pass the node positions to preCICE. Each node is + // specified only once. + std::vector interface_nodes_positions; + interface_nodes_positions.reserve(dim * n_interface_nodes); + + // Set up the appropriate size of the data container needed for data + // exchange. Here, we deal with a scalar problem, so that only a scalar value + // is read/written per interface node. + read_data.resize(n_interface_nodes); + // The IDs are again filled by preCICE during the initializations. + interface_nodes_ids.resize(n_interface_nodes); + + // The node location is obtained using `map_dofs_to_support_points()`. + std::map> support_points; + DoFTools::map_dofs_to_support_points(mapping, dof_handler, support_points); + + // `support_points` contains now the coordinates of all DoFs. In the next + // step, the relevant coordinates are extracted using the IndexSet with the + // extracted coupling_dofs. + for (const auto element : coupling_dofs) + for (int i = 0; i < dim; ++i) + interface_nodes_positions.push_back(support_points[element][i]); + + // Now we have all information to define the coupling mesh and pass the + // information to preCICE. + precice.setMeshVertices(mesh_id, + n_interface_nodes, + interface_nodes_positions.data(), + interface_nodes_ids.data()); + + // Then, we initialize preCICE internally calling the API function + // `initialize()` + const double max_delta_t = precice.initialize(); + + + // read first coupling data from preCICE if available (i.e. deal.II is + // the second participant in a serial coupling scheme) + if (precice.isReadDataAvailable()) + { + precice.readBlockScalarData(read_data_id, + n_interface_nodes, + interface_nodes_ids.data(), + read_data.data()); + + // After receiving the coupling data in `read_data`, we convert it to + // the std::map `boundary_data` which is later needed in order to apply + // Dirichlet boundary conditions + format_precice_to_dealii(boundary_data); + } + + return max_delta_t; +} + + +// The function `advance()` is called in the main time loop after the +// computation in each time step. Here, +// coupling data is passed to and obtained from preCICE. +template +double +Adapter::advance( + std::map &boundary_data, + const double computed_timestep_length) +{ + // We specify the computed time-step length and pass it to preCICE. In + // return, preCICE tells us the maximum admissible time-step size our + // participant is allowed to compute in order to not exceed the next coupling + // time step. + const double max_delta_t = precice.advance(computed_timestep_length); + + // As a next step, we obtain data, i.e. the boundary condition, from another + // participant. We have already all IDs and just need to convert our obtained + // data to the deal.II compatible 'boundary map' , which is done in the + // format_deal_to_precice function. + precice.readBlockScalarData(read_data_id, + n_interface_nodes, + interface_nodes_ids.data(), + read_data.data()); + + format_precice_to_dealii(boundary_data); + + return max_delta_t; +} + + + +// This function takes the std::vector obtained by preCICE in `read_data` and +// inserts the values to the right position in the boundary map used throughout +// our deal.II solver for Dirichlet boundary conditions. The function is only +// used internally in the Adapter class and not called in the solver itself. The +// order, in which preCICE sorts the data in the `read_data` vector is exactly +// the same as the order of the initially passed vertices coordinates. +template +void +Adapter::format_precice_to_dealii( + std::map &boundary_data) const +{ + // We already stored the coupling DoF indices in the `boundary_data` map, so + // that we can simply iterate over all keys in the map. + auto dof_component = boundary_data.begin(); + for (int i = 0; i < n_interface_nodes; ++i) + { + AssertIndexRange(i, read_data.size()); + boundary_data[dof_component->first] = read_data[i]; + ++dof_component; + } +} + + +// The solver class is essentially the same as in step-4. We only extend the +// stationary problem to a time-dependent problem and introduced the coupling. +// Comments are added at any point, where the workflow differs from step-4. +template +class CoupledLaplaceProblem +{ +public: + CoupledLaplaceProblem(); + + void + run(); + +private: + void + make_grid(); + void + setup_system(); + void + assemble_system(); + void + solve(); + void + output_results() const; + + Triangulation triangulation; + FE_Q fe; + DoFHandler dof_handler; + MappingQ1 mapping; + + SparsityPattern sparsity_pattern; + SparseMatrix system_matrix; + + Vector solution; + Vector old_solution; + Vector system_rhs; + + // We allocate all structures required for the preCICE coupling: The map + // is used to apply Dirichlet boundary conditions and filled in the Adapter + // class with data from the other participant. The CouplingParameters hold the + // preCICE configuration as described above. The interface boundary ID is the + // ID associated to our coupling interface and needs to be specified, when we + // set up the Adapter class object, because we pass it directly to the + // Constructor of this class. + std::map boundary_data; + CouplingParamters parameters; + const types::boundary_id interface_boundary_id; + Adapter adapter; + + // The time-step size delta_t is the acutual time-step size used for all + // computations. The preCICE time-step size is obtained by preCICE in order to + // ensure a synchronization at all coupling time steps. The solver time + // step-size is the desired time-step size of our individual solver. In more + // sophisticated computations, it might be determined adaptively. The + // `time_step` counter is just used for the time-step number. + double delta_t; + double precice_delta_t; + const double solver_delta_t = 0.1; + unsigned int time_step = 0; +}; + + + +template +class RightHandSide : public Function +{ +public: + virtual double + value(const Point &p, const unsigned int component = 0) const override; +}; + + + +template +class BoundaryValues : public Function +{ +public: + virtual double + value(const Point &p, const unsigned int component = 0) const override; +}; + +template +double +RightHandSide::value(const Point &p, + const unsigned int /*component*/) const +{ + double return_value = 0.0; + for (unsigned int i = 0; i < dim; ++i) + return_value += 4.0 * std::pow(p(i), 4.0); + + return return_value; +} + + +template +double +BoundaryValues::value(const Point &p, + const unsigned int /*component*/) const +{ + return p.square(); +} + + + +template +CoupledLaplaceProblem::CoupledLaplaceProblem() + : fe(1) + , dof_handler(triangulation) + , interface_boundary_id(1) + , adapter(parameters, interface_boundary_id) +{} + + +template +void +CoupledLaplaceProblem::make_grid() +{ + GridGenerator::hyper_cube(triangulation, -1, 1); + triangulation.refine_global(4); + + for (const auto &cell : triangulation.active_cell_iterators()) + for (const auto face : cell->face_iterators()) + { + // We choose the boundary in positive x direction for the + // interface coupling. + if (face->at_boundary() && (face->center()[0] == 1)) + face->set_boundary_id(interface_boundary_id); + } + + std::cout << " Number of active cells: " << triangulation.n_active_cells() + << std::endl + << " Total number of cells: " << triangulation.n_cells() + << std::endl; +} + + +template +void +CoupledLaplaceProblem::setup_system() +{ + dof_handler.distribute_dofs(fe); + + std::cout << " Number of degrees of freedom: " << dof_handler.n_dofs() + << std::endl; + + DynamicSparsityPattern dsp(dof_handler.n_dofs()); + DoFTools::make_sparsity_pattern(dof_handler, dsp); + sparsity_pattern.copy_from(dsp); + + system_matrix.reinit(sparsity_pattern); + + solution.reinit(dof_handler.n_dofs()); + old_solution.reinit(dof_handler.n_dofs()); + system_rhs.reinit(dof_handler.n_dofs()); +} + + + +template +void +CoupledLaplaceProblem::assemble_system() +{ + // Reset global structures + system_rhs = 0; + system_matrix = 0; + // Update old solution values + old_solution = solution; + + QGauss quadrature_formula(fe.degree + 1); + + RightHandSide right_hand_side; + + FEValues fe_values(fe, + quadrature_formula, + update_values | update_gradients | + update_quadrature_points | update_JxW_values); + + const unsigned int dofs_per_cell = fe.n_dofs_per_cell(); + + FullMatrix cell_matrix(dofs_per_cell, dofs_per_cell); + Vector cell_rhs(dofs_per_cell); + std::vector local_dof_indices(dofs_per_cell); + // The solution values from previous time steps are stored for each quadrature + // point + std::vector local_values_old_solution(fe_values.n_quadrature_points); + + for (const auto &cell : dof_handler.active_cell_iterators()) + { + fe_values.reinit(cell); + cell_matrix = 0; + cell_rhs = 0; + // Get the local values from the `fe_values' object + fe_values.get_function_values(old_solution, local_values_old_solution); + + // The system matrix contains additionally a mass matrix due to the time + // discretization. The RHS has contributions from the old solution values. + for (const unsigned int q_index : fe_values.quadrature_point_indices()) + for (const unsigned int i : fe_values.dof_indices()) + { + for (const unsigned int j : fe_values.dof_indices()) + cell_matrix(i, j) += + ((fe_values.shape_value(i, q_index) * // phi_i(x_q) + fe_values.shape_value(j, q_index)) + // phi_j(x_q) + (delta_t * // delta t + fe_values.shape_grad(i, q_index) * // grad phi_i(x_q) + fe_values.shape_grad(j, q_index))) * // grad phi_j(x_q) + fe_values.JxW(q_index); // dx + + const auto x_q = fe_values.quadrature_point(q_index); + const auto &local_value = local_values_old_solution[q_index]; + cell_rhs(i) += ((delta_t * // delta t + fe_values.shape_value(i, q_index) * // phi_i(x_q) + right_hand_side.value(x_q)) + // f(x_q) + fe_values.shape_value(i, q_index) * + local_value) * // phi_i(x_q)*val + fe_values.JxW(q_index); // dx + } + + // Copy local to global + cell->get_dof_indices(local_dof_indices); + for (const unsigned int i : fe_values.dof_indices()) + { + for (const unsigned int j : fe_values.dof_indices()) + system_matrix.add(local_dof_indices[i], + local_dof_indices[j], + cell_matrix(i, j)); + + system_rhs(local_dof_indices[i]) += cell_rhs(i); + } + } + { + // At first, we apply the Dirichlet boundary condition from step-4, as + // usual. + std::map boundary_values; + VectorTools::interpolate_boundary_values(dof_handler, + 0, + BoundaryValues(), + boundary_values); + MatrixTools::apply_boundary_values(boundary_values, + system_matrix, + solution, + system_rhs); + } + { + // Afterwards, we apply the coupling boundary condition. The `boundary_data` + // has already been filled by preCICE. + MatrixTools::apply_boundary_values(boundary_data, + system_matrix, + solution, + system_rhs); + } +} + + + +template +void +CoupledLaplaceProblem::solve() +{ + SolverControl solver_control(1000, 1e-12); + SolverCG> solver(solver_control); + solver.solve(system_matrix, solution, system_rhs, PreconditionIdentity()); + + std::cout << " " << solver_control.last_step() + << " CG iterations needed to obtain convergence." << std::endl; +} + + + +template +void +CoupledLaplaceProblem::output_results() const +{ + DataOut data_out; + + data_out.attach_dof_handler(dof_handler); + data_out.add_data_vector(solution, "solution"); + + data_out.build_patches(mapping); + + std::ofstream output("solution-" + std::to_string(time_step) + ".vtk"); + data_out.write_vtk(output); +} + + + +template +void +CoupledLaplaceProblem::run() +{ + std::cout << "Solving problem in " << dim << " space dimensions." + << std::endl; + + make_grid(); + setup_system(); + + // After we set up the system, we initialize preCICE using the functionalities + // of the Adapter. preCICE returns the maximum admissible time-step size, + // which needs to be compared to our desired solver time-step size. + precice_delta_t = adapter.initialize(dof_handler, boundary_data, mapping); + delta_t = std::min(precice_delta_t, solver_delta_t); + + // preCICE steers the coupled simulation: `isCouplingOngoing` is + // used to synchronize the end of the simulation with the coupling partner + while (adapter.precice.isCouplingOngoing()) + { + // The time step number is solely used to generate unique output files + ++time_step; + // In the time loop, we assemble the coupled system and solve it as + // usual. + assemble_system(); + solve(); + + // After we solved the system, we advance the coupling to the next time + // level. In a bi-directional coupled simulation, we would pass our + // calculated data to and obtain new data from preCICE. Here, we simply + // obtain new data from preCICE, so from the other participant. As before, + // we obtain a maximum time-step size and compare it against the desired + // solver time-step size. + precice_delta_t = adapter.advance(boundary_data, delta_t); + delta_t = std::min(precice_delta_t, solver_delta_t); + + // Write an output file if the time step is completed. In case of an + // implicit coupling, where individual time steps are computed more than + // once, the function `isTimeWindowCompleted` prevents unnecessary result + // writing. For this simple tutorial configuration (explicit coupling), + // the function returns always `true`. + if (adapter.precice.isTimeWindowComplete()) + output_results(); + } +} + + + +int +main() +{ + CoupledLaplaceProblem<2> laplace_problem; + laplace_problem.run(); + + return 0; +} diff --git a/coupled_laplace_problem/doc/author b/coupled_laplace_problem/doc/author new file mode 100644 index 0000000..f236ba2 --- /dev/null +++ b/coupled_laplace_problem/doc/author @@ -0,0 +1,2 @@ +David Schneider +Benjamin Uekermann diff --git a/coupled_laplace_problem/doc/builds-on b/coupled_laplace_problem/doc/builds-on new file mode 100644 index 0000000..48a0f73 --- /dev/null +++ b/coupled_laplace_problem/doc/builds-on @@ -0,0 +1 @@ +step-4 diff --git a/coupled_laplace_problem/doc/dependencies b/coupled_laplace_problem/doc/dependencies new file mode 100644 index 0000000..8b501d8 --- /dev/null +++ b/coupled_laplace_problem/doc/dependencies @@ -0,0 +1 @@ +preCICE diff --git a/coupled_laplace_problem/doc/entry-name b/coupled_laplace_problem/doc/entry-name new file mode 100644 index 0000000..e28d1e6 --- /dev/null +++ b/coupled_laplace_problem/doc/entry-name @@ -0,0 +1 @@ +Laplace equation coupled to an external simulation program diff --git a/coupled_laplace_problem/doc/result.gif b/coupled_laplace_problem/doc/result.gif new file mode 100644 index 0000000..4f79cfd Binary files /dev/null and b/coupled_laplace_problem/doc/result.gif differ diff --git a/coupled_laplace_problem/doc/tooltip b/coupled_laplace_problem/doc/tooltip new file mode 100644 index 0000000..12df277 --- /dev/null +++ b/coupled_laplace_problem/doc/tooltip @@ -0,0 +1 @@ +Laplace equation surface coupled to an external simulation program (here simply a fancy boundary condition) using the coupling library preCICE. diff --git a/coupled_laplace_problem/fancy_boundary_condition.cc b/coupled_laplace_problem/fancy_boundary_condition.cc new file mode 100644 index 0000000..46fc24a --- /dev/null +++ b/coupled_laplace_problem/fancy_boundary_condition.cc @@ -0,0 +1,117 @@ +// This program does not use any deal.II functionality and depends only on +// preCICE and the standard libraries. +#include + +#include +#include + +// The program computes a time-varying parabolic boundary condition, which is +// passed to preCICE and serves as Dirichlet boundary condition for the other +// coupling participant. + +// Function to generate boundary values in each time step +void +define_boundary_values(std::vector &boundary_data, + const double time, + const double end_time) +{ + // Scale the current time value + const double relative_time = time / end_time; + // Define the amplitude. Values run from -0.5 to 0.5 + const double amplitude = (relative_time - 0.5); + // Specify the actual data we want to pass to the other participant. Here, we + // choose a parabola with boundary values 2 in order to enforce continuity + // to adjacent boundaries. + const double n_elements = boundary_data.size(); + const double right_zero = boundary_data.size() - 1; + const double left_zero = 0; + const double offset = 2; + for (uint i = 0; i < n_elements; ++i) + boundary_data[i] = + -amplitude * ((i - left_zero) * (i - right_zero)) + offset; +} + + +int +main() +{ + std::cout << "Boundary participant: starting... \n"; + + // Configuration + const std::string configFileName("precice-config.xml"); + const std::string solverName("boundary-participant"); + const std::string meshName("boundary-mesh"); + const std::string dataWriteName("boundary-data"); + + // Adjust to MPI rank and size for parallel computation + const int commRank = 0; + const int commSize = 1; + + precice::SolverInterface precice(solverName, + configFileName, + commRank, + commSize); + + const int meshID = precice.getMeshID(meshName); + const int dimensions = precice.getDimensions(); + const int numberOfVertices = 6; + + const int dataID = precice.getDataID(dataWriteName, meshID); + + // Set up data structures + std::vector writeData(numberOfVertices); + std::vector vertexIDs(numberOfVertices); + std::vector vertices(numberOfVertices * dimensions); + + // Define a boundary mesh + std::cout << "Boundary participant: defining boundary mesh \n"; + const double length = 2; + const double xCoord = 1; + const double deltaY = length / (numberOfVertices - 1); + for (int i = 0; i < numberOfVertices; ++i) + for (int j = 0; j < dimensions; ++j) + { + const unsigned int index = dimensions * i + j; + // The x-coordinate is always 1, i.e., the boundary is parallel to the + // y-axis. The y-coordinate is descending from 1 to -1. + if (j == 0) + vertices[index] = xCoord; + else + vertices[index] = 1 - deltaY * i; + } + + // Pass the vertices to preCICE + precice.setMeshVertices(meshID, + numberOfVertices, + vertices.data(), + vertexIDs.data()); + + // initialize the Solverinterface + double dt = precice.initialize(); + + // Start time loop + const double end_time = 1; + double time = 0; + while (precice.isCouplingOngoing()) + { + // Generate new boundary data + define_boundary_values(writeData, time, end_time); + + { + std::cout << "Boundary participant: writing coupling data \n"; + precice.writeBlockScalarData(dataID, + numberOfVertices, + vertexIDs.data(), + writeData.data()); + } + + dt = precice.advance(dt); + std::cout << "Boundary participant: advancing in time\n"; + + time += dt; + } + + std::cout << "Boundary participant: closing...\n"; + + return 0; +} diff --git a/coupled_laplace_problem/precice-config.xml b/coupled_laplace_problem/precice-config.xml new file mode 100644 index 0000000..d1e6722 --- /dev/null +++ b/coupled_laplace_problem/precice-config.xml @@ -0,0 +1,52 @@ + + + + + + + + + + + + + + + + + + + + + + + + + + + + + + + + + + + + + + + diff --git a/coupled_laplace_problem/test_data/reference-10.vtk b/coupled_laplace_problem/test_data/reference-10.vtk new file mode 100644 index 0000000..d8bc5e4 --- /dev/null +++ b/coupled_laplace_problem/test_data/reference-10.vtk @@ -0,0 +1,1294 @@ +# vtk DataFile Version 3.0 +ASCII +DATASET UNSTRUCTURED_GRID + +POINTS 1024 double +-1 -1 0 +-0.875 -1 0 +-1 -0.875 0 +-0.875 -0.875 0 +-0.875 -1 0 +-0.75 -1 0 +-0.875 -0.875 0 +-0.75 -0.875 0 +-1 -0.875 0 +-0.875 -0.875 0 +-1 -0.75 0 +-0.875 -0.75 0 +-0.875 -0.875 0 +-0.75 -0.875 0 +-0.875 -0.75 0 +-0.75 -0.75 0 +-0.75 -1 0 +-0.625 -1 0 +-0.75 -0.875 0 +-0.625 -0.875 0 +-0.625 -1 0 +-0.5 -1 0 +-0.625 -0.875 0 +-0.5 -0.875 0 +-0.75 -0.875 0 +-0.625 -0.875 0 +-0.75 -0.75 0 +-0.625 -0.75 0 +-0.625 -0.875 0 +-0.5 -0.875 0 +-0.625 -0.75 0 +-0.5 -0.75 0 +-1 -0.75 0 +-0.875 -0.75 0 +-1 -0.625 0 +-0.875 -0.625 0 +-0.875 -0.75 0 +-0.75 -0.75 0 +-0.875 -0.625 0 +-0.75 -0.625 0 +-1 -0.625 0 +-0.875 -0.625 0 +-1 -0.5 0 +-0.875 -0.5 0 +-0.875 -0.625 0 +-0.75 -0.625 0 +-0.875 -0.5 0 +-0.75 -0.5 0 +-0.75 -0.75 0 +-0.625 -0.75 0 +-0.75 -0.625 0 +-0.625 -0.625 0 +-0.625 -0.75 0 +-0.5 -0.75 0 +-0.625 -0.625 0 +-0.5 -0.625 0 +-0.75 -0.625 0 +-0.625 -0.625 0 +-0.75 -0.5 0 +-0.625 -0.5 0 +-0.625 -0.625 0 +-0.5 -0.625 0 +-0.625 -0.5 0 +-0.5 -0.5 0 +-0.5 -1 0 +-0.375 -1 0 +-0.5 -0.875 0 +-0.375 -0.875 0 +-0.375 -1 0 +-0.25 -1 0 +-0.375 -0.875 0 +-0.25 -0.875 0 +-0.5 -0.875 0 +-0.375 -0.875 0 +-0.5 -0.75 0 +-0.375 -0.75 0 +-0.375 -0.875 0 +-0.25 -0.875 0 +-0.375 -0.75 0 +-0.25 -0.75 0 +-0.25 -1 0 +-0.125 -1 0 +-0.25 -0.875 0 +-0.125 -0.875 0 +-0.125 -1 0 +0 -1 0 +-0.125 -0.875 0 +0 -0.875 0 +-0.25 -0.875 0 +-0.125 -0.875 0 +-0.25 -0.75 0 +-0.125 -0.75 0 +-0.125 -0.875 0 +0 -0.875 0 +-0.125 -0.75 0 +0 -0.75 0 +-0.5 -0.75 0 +-0.375 -0.75 0 +-0.5 -0.625 0 +-0.375 -0.625 0 +-0.375 -0.75 0 +-0.25 -0.75 0 +-0.375 -0.625 0 +-0.25 -0.625 0 +-0.5 -0.625 0 +-0.375 -0.625 0 +-0.5 -0.5 0 +-0.375 -0.5 0 +-0.375 -0.625 0 +-0.25 -0.625 0 +-0.375 -0.5 0 +-0.25 -0.5 0 +-0.25 -0.75 0 +-0.125 -0.75 0 +-0.25 -0.625 0 +-0.125 -0.625 0 +-0.125 -0.75 0 +0 -0.75 0 +-0.125 -0.625 0 +0 -0.625 0 +-0.25 -0.625 0 +-0.125 -0.625 0 +-0.25 -0.5 0 +-0.125 -0.5 0 +-0.125 -0.625 0 +0 -0.625 0 +-0.125 -0.5 0 +0 -0.5 0 +-1 -0.5 0 +-0.875 -0.5 0 +-1 -0.375 0 +-0.875 -0.375 0 +-0.875 -0.5 0 +-0.75 -0.5 0 +-0.875 -0.375 0 +-0.75 -0.375 0 +-1 -0.375 0 +-0.875 -0.375 0 +-1 -0.25 0 +-0.875 -0.25 0 +-0.875 -0.375 0 +-0.75 -0.375 0 +-0.875 -0.25 0 +-0.75 -0.25 0 +-0.75 -0.5 0 +-0.625 -0.5 0 +-0.75 -0.375 0 +-0.625 -0.375 0 +-0.625 -0.5 0 +-0.5 -0.5 0 +-0.625 -0.375 0 +-0.5 -0.375 0 +-0.75 -0.375 0 +-0.625 -0.375 0 +-0.75 -0.25 0 +-0.625 -0.25 0 +-0.625 -0.375 0 +-0.5 -0.375 0 +-0.625 -0.25 0 +-0.5 -0.25 0 +-1 -0.25 0 +-0.875 -0.25 0 +-1 -0.125 0 +-0.875 -0.125 0 +-0.875 -0.25 0 +-0.75 -0.25 0 +-0.875 -0.125 0 +-0.75 -0.125 0 +-1 -0.125 0 +-0.875 -0.125 0 +-1 0 0 +-0.875 0 0 +-0.875 -0.125 0 +-0.75 -0.125 0 +-0.875 0 0 +-0.75 0 0 +-0.75 -0.25 0 +-0.625 -0.25 0 +-0.75 -0.125 0 +-0.625 -0.125 0 +-0.625 -0.25 0 +-0.5 -0.25 0 +-0.625 -0.125 0 +-0.5 -0.125 0 +-0.75 -0.125 0 +-0.625 -0.125 0 +-0.75 0 0 +-0.625 0 0 +-0.625 -0.125 0 +-0.5 -0.125 0 +-0.625 0 0 +-0.5 0 0 +-0.5 -0.5 0 +-0.375 -0.5 0 +-0.5 -0.375 0 +-0.375 -0.375 0 +-0.375 -0.5 0 +-0.25 -0.5 0 +-0.375 -0.375 0 +-0.25 -0.375 0 +-0.5 -0.375 0 +-0.375 -0.375 0 +-0.5 -0.25 0 +-0.375 -0.25 0 +-0.375 -0.375 0 +-0.25 -0.375 0 +-0.375 -0.25 0 +-0.25 -0.25 0 +-0.25 -0.5 0 +-0.125 -0.5 0 +-0.25 -0.375 0 +-0.125 -0.375 0 +-0.125 -0.5 0 +0 -0.5 0 +-0.125 -0.375 0 +0 -0.375 0 +-0.25 -0.375 0 +-0.125 -0.375 0 +-0.25 -0.25 0 +-0.125 -0.25 0 +-0.125 -0.375 0 +0 -0.375 0 +-0.125 -0.25 0 +0 -0.25 0 +-0.5 -0.25 0 +-0.375 -0.25 0 +-0.5 -0.125 0 +-0.375 -0.125 0 +-0.375 -0.25 0 +-0.25 -0.25 0 +-0.375 -0.125 0 +-0.25 -0.125 0 +-0.5 -0.125 0 +-0.375 -0.125 0 +-0.5 0 0 +-0.375 0 0 +-0.375 -0.125 0 +-0.25 -0.125 0 +-0.375 0 0 +-0.25 0 0 +-0.25 -0.25 0 +-0.125 -0.25 0 +-0.25 -0.125 0 +-0.125 -0.125 0 +-0.125 -0.25 0 +0 -0.25 0 +-0.125 -0.125 0 +0 -0.125 0 +-0.25 -0.125 0 +-0.125 -0.125 0 +-0.25 0 0 +-0.125 0 0 +-0.125 -0.125 0 +0 -0.125 0 +-0.125 0 0 +0 0 0 +0 -1 0 +0.125 -1 0 +0 -0.875 0 +0.125 -0.875 0 +0.125 -1 0 +0.25 -1 0 +0.125 -0.875 0 +0.25 -0.875 0 +0 -0.875 0 +0.125 -0.875 0 +0 -0.75 0 +0.125 -0.75 0 +0.125 -0.875 0 +0.25 -0.875 0 +0.125 -0.75 0 +0.25 -0.75 0 +0.25 -1 0 +0.375 -1 0 +0.25 -0.875 0 +0.375 -0.875 0 +0.375 -1 0 +0.5 -1 0 +0.375 -0.875 0 +0.5 -0.875 0 +0.25 -0.875 0 +0.375 -0.875 0 +0.25 -0.75 0 +0.375 -0.75 0 +0.375 -0.875 0 +0.5 -0.875 0 +0.375 -0.75 0 +0.5 -0.75 0 +0 -0.75 0 +0.125 -0.75 0 +0 -0.625 0 +0.125 -0.625 0 +0.125 -0.75 0 +0.25 -0.75 0 +0.125 -0.625 0 +0.25 -0.625 0 +0 -0.625 0 +0.125 -0.625 0 +0 -0.5 0 +0.125 -0.5 0 +0.125 -0.625 0 +0.25 -0.625 0 +0.125 -0.5 0 +0.25 -0.5 0 +0.25 -0.75 0 +0.375 -0.75 0 +0.25 -0.625 0 +0.375 -0.625 0 +0.375 -0.75 0 +0.5 -0.75 0 +0.375 -0.625 0 +0.5 -0.625 0 +0.25 -0.625 0 +0.375 -0.625 0 +0.25 -0.5 0 +0.375 -0.5 0 +0.375 -0.625 0 +0.5 -0.625 0 +0.375 -0.5 0 +0.5 -0.5 0 +0.5 -1 0 +0.625 -1 0 +0.5 -0.875 0 +0.625 -0.875 0 +0.625 -1 0 +0.75 -1 0 +0.625 -0.875 0 +0.75 -0.875 0 +0.5 -0.875 0 +0.625 -0.875 0 +0.5 -0.75 0 +0.625 -0.75 0 +0.625 -0.875 0 +0.75 -0.875 0 +0.625 -0.75 0 +0.75 -0.75 0 +0.75 -1 0 +0.875 -1 0 +0.75 -0.875 0 +0.875 -0.875 0 +0.875 -1 0 +1 -1 0 +0.875 -0.875 0 +1 -0.875 0 +0.75 -0.875 0 +0.875 -0.875 0 +0.75 -0.75 0 +0.875 -0.75 0 +0.875 -0.875 0 +1 -0.875 0 +0.875 -0.75 0 +1 -0.75 0 +0.5 -0.75 0 +0.625 -0.75 0 +0.5 -0.625 0 +0.625 -0.625 0 +0.625 -0.75 0 +0.75 -0.75 0 +0.625 -0.625 0 +0.75 -0.625 0 +0.5 -0.625 0 +0.625 -0.625 0 +0.5 -0.5 0 +0.625 -0.5 0 +0.625 -0.625 0 +0.75 -0.625 0 +0.625 -0.5 0 +0.75 -0.5 0 +0.75 -0.75 0 +0.875 -0.75 0 +0.75 -0.625 0 +0.875 -0.625 0 +0.875 -0.75 0 +1 -0.75 0 +0.875 -0.625 0 +1 -0.625 0 +0.75 -0.625 0 +0.875 -0.625 0 +0.75 -0.5 0 +0.875 -0.5 0 +0.875 -0.625 0 +1 -0.625 0 +0.875 -0.5 0 +1 -0.5 0 +0 -0.5 0 +0.125 -0.5 0 +0 -0.375 0 +0.125 -0.375 0 +0.125 -0.5 0 +0.25 -0.5 0 +0.125 -0.375 0 +0.25 -0.375 0 +0 -0.375 0 +0.125 -0.375 0 +0 -0.25 0 +0.125 -0.25 0 +0.125 -0.375 0 +0.25 -0.375 0 +0.125 -0.25 0 +0.25 -0.25 0 +0.25 -0.5 0 +0.375 -0.5 0 +0.25 -0.375 0 +0.375 -0.375 0 +0.375 -0.5 0 +0.5 -0.5 0 +0.375 -0.375 0 +0.5 -0.375 0 +0.25 -0.375 0 +0.375 -0.375 0 +0.25 -0.25 0 +0.375 -0.25 0 +0.375 -0.375 0 +0.5 -0.375 0 +0.375 -0.25 0 +0.5 -0.25 0 +0 -0.25 0 +0.125 -0.25 0 +0 -0.125 0 +0.125 -0.125 0 +0.125 -0.25 0 +0.25 -0.25 0 +0.125 -0.125 0 +0.25 -0.125 0 +0 -0.125 0 +0.125 -0.125 0 +0 0 0 +0.125 0 0 +0.125 -0.125 0 +0.25 -0.125 0 +0.125 0 0 +0.25 0 0 +0.25 -0.25 0 +0.375 -0.25 0 +0.25 -0.125 0 +0.375 -0.125 0 +0.375 -0.25 0 +0.5 -0.25 0 +0.375 -0.125 0 +0.5 -0.125 0 +0.25 -0.125 0 +0.375 -0.125 0 +0.25 0 0 +0.375 0 0 +0.375 -0.125 0 +0.5 -0.125 0 +0.375 0 0 +0.5 0 0 +0.5 -0.5 0 +0.625 -0.5 0 +0.5 -0.375 0 +0.625 -0.375 0 +0.625 -0.5 0 +0.75 -0.5 0 +0.625 -0.375 0 +0.75 -0.375 0 +0.5 -0.375 0 +0.625 -0.375 0 +0.5 -0.25 0 +0.625 -0.25 0 +0.625 -0.375 0 +0.75 -0.375 0 +0.625 -0.25 0 +0.75 -0.25 0 +0.75 -0.5 0 +0.875 -0.5 0 +0.75 -0.375 0 +0.875 -0.375 0 +0.875 -0.5 0 +1 -0.5 0 +0.875 -0.375 0 +1 -0.375 0 +0.75 -0.375 0 +0.875 -0.375 0 +0.75 -0.25 0 +0.875 -0.25 0 +0.875 -0.375 0 +1 -0.375 0 +0.875 -0.25 0 +1 -0.25 0 +0.5 -0.25 0 +0.625 -0.25 0 +0.5 -0.125 0 +0.625 -0.125 0 +0.625 -0.25 0 +0.75 -0.25 0 +0.625 -0.125 0 +0.75 -0.125 0 +0.5 -0.125 0 +0.625 -0.125 0 +0.5 0 0 +0.625 0 0 +0.625 -0.125 0 +0.75 -0.125 0 +0.625 0 0 +0.75 0 0 +0.75 -0.25 0 +0.875 -0.25 0 +0.75 -0.125 0 +0.875 -0.125 0 +0.875 -0.25 0 +1 -0.25 0 +0.875 -0.125 0 +1 -0.125 0 +0.75 -0.125 0 +0.875 -0.125 0 +0.75 0 0 +0.875 0 0 +0.875 -0.125 0 +1 -0.125 0 +0.875 0 0 +1 0 0 +-1 0 0 +-0.875 0 0 +-1 0.125 0 +-0.875 0.125 0 +-0.875 0 0 +-0.75 0 0 +-0.875 0.125 0 +-0.75 0.125 0 +-1 0.125 0 +-0.875 0.125 0 +-1 0.25 0 +-0.875 0.25 0 +-0.875 0.125 0 +-0.75 0.125 0 +-0.875 0.25 0 +-0.75 0.25 0 +-0.75 0 0 +-0.625 0 0 +-0.75 0.125 0 +-0.625 0.125 0 +-0.625 0 0 +-0.5 0 0 +-0.625 0.125 0 +-0.5 0.125 0 +-0.75 0.125 0 +-0.625 0.125 0 +-0.75 0.25 0 +-0.625 0.25 0 +-0.625 0.125 0 +-0.5 0.125 0 +-0.625 0.25 0 +-0.5 0.25 0 +-1 0.25 0 +-0.875 0.25 0 +-1 0.375 0 +-0.875 0.375 0 +-0.875 0.25 0 +-0.75 0.25 0 +-0.875 0.375 0 +-0.75 0.375 0 +-1 0.375 0 +-0.875 0.375 0 +-1 0.5 0 +-0.875 0.5 0 +-0.875 0.375 0 +-0.75 0.375 0 +-0.875 0.5 0 +-0.75 0.5 0 +-0.75 0.25 0 +-0.625 0.25 0 +-0.75 0.375 0 +-0.625 0.375 0 +-0.625 0.25 0 +-0.5 0.25 0 +-0.625 0.375 0 +-0.5 0.375 0 +-0.75 0.375 0 +-0.625 0.375 0 +-0.75 0.5 0 +-0.625 0.5 0 +-0.625 0.375 0 +-0.5 0.375 0 +-0.625 0.5 0 +-0.5 0.5 0 +-0.5 0 0 +-0.375 0 0 +-0.5 0.125 0 +-0.375 0.125 0 +-0.375 0 0 +-0.25 0 0 +-0.375 0.125 0 +-0.25 0.125 0 +-0.5 0.125 0 +-0.375 0.125 0 +-0.5 0.25 0 +-0.375 0.25 0 +-0.375 0.125 0 +-0.25 0.125 0 +-0.375 0.25 0 +-0.25 0.25 0 +-0.25 0 0 +-0.125 0 0 +-0.25 0.125 0 +-0.125 0.125 0 +-0.125 0 0 +0 0 0 +-0.125 0.125 0 +0 0.125 0 +-0.25 0.125 0 +-0.125 0.125 0 +-0.25 0.25 0 +-0.125 0.25 0 +-0.125 0.125 0 +0 0.125 0 +-0.125 0.25 0 +0 0.25 0 +-0.5 0.25 0 +-0.375 0.25 0 +-0.5 0.375 0 +-0.375 0.375 0 +-0.375 0.25 0 +-0.25 0.25 0 +-0.375 0.375 0 +-0.25 0.375 0 +-0.5 0.375 0 +-0.375 0.375 0 +-0.5 0.5 0 +-0.375 0.5 0 +-0.375 0.375 0 +-0.25 0.375 0 +-0.375 0.5 0 +-0.25 0.5 0 +-0.25 0.25 0 +-0.125 0.25 0 +-0.25 0.375 0 +-0.125 0.375 0 +-0.125 0.25 0 +0 0.25 0 +-0.125 0.375 0 +0 0.375 0 +-0.25 0.375 0 +-0.125 0.375 0 +-0.25 0.5 0 +-0.125 0.5 0 +-0.125 0.375 0 +0 0.375 0 +-0.125 0.5 0 +0 0.5 0 +-1 0.5 0 +-0.875 0.5 0 +-1 0.625 0 +-0.875 0.625 0 +-0.875 0.5 0 +-0.75 0.5 0 +-0.875 0.625 0 +-0.75 0.625 0 +-1 0.625 0 +-0.875 0.625 0 +-1 0.75 0 +-0.875 0.75 0 +-0.875 0.625 0 +-0.75 0.625 0 +-0.875 0.75 0 +-0.75 0.75 0 +-0.75 0.5 0 +-0.625 0.5 0 +-0.75 0.625 0 +-0.625 0.625 0 +-0.625 0.5 0 +-0.5 0.5 0 +-0.625 0.625 0 +-0.5 0.625 0 +-0.75 0.625 0 +-0.625 0.625 0 +-0.75 0.75 0 +-0.625 0.75 0 +-0.625 0.625 0 +-0.5 0.625 0 +-0.625 0.75 0 +-0.5 0.75 0 +-1 0.75 0 +-0.875 0.75 0 +-1 0.875 0 +-0.875 0.875 0 +-0.875 0.75 0 +-0.75 0.75 0 +-0.875 0.875 0 +-0.75 0.875 0 +-1 0.875 0 +-0.875 0.875 0 +-1 1 0 +-0.875 1 0 +-0.875 0.875 0 +-0.75 0.875 0 +-0.875 1 0 +-0.75 1 0 +-0.75 0.75 0 +-0.625 0.75 0 +-0.75 0.875 0 +-0.625 0.875 0 +-0.625 0.75 0 +-0.5 0.75 0 +-0.625 0.875 0 +-0.5 0.875 0 +-0.75 0.875 0 +-0.625 0.875 0 +-0.75 1 0 +-0.625 1 0 +-0.625 0.875 0 +-0.5 0.875 0 +-0.625 1 0 +-0.5 1 0 +-0.5 0.5 0 +-0.375 0.5 0 +-0.5 0.625 0 +-0.375 0.625 0 +-0.375 0.5 0 +-0.25 0.5 0 +-0.375 0.625 0 +-0.25 0.625 0 +-0.5 0.625 0 +-0.375 0.625 0 +-0.5 0.75 0 +-0.375 0.75 0 +-0.375 0.625 0 +-0.25 0.625 0 +-0.375 0.75 0 +-0.25 0.75 0 +-0.25 0.5 0 +-0.125 0.5 0 +-0.25 0.625 0 +-0.125 0.625 0 +-0.125 0.5 0 +0 0.5 0 +-0.125 0.625 0 +0 0.625 0 +-0.25 0.625 0 +-0.125 0.625 0 +-0.25 0.75 0 +-0.125 0.75 0 +-0.125 0.625 0 +0 0.625 0 +-0.125 0.75 0 +0 0.75 0 +-0.5 0.75 0 +-0.375 0.75 0 +-0.5 0.875 0 +-0.375 0.875 0 +-0.375 0.75 0 +-0.25 0.75 0 +-0.375 0.875 0 +-0.25 0.875 0 +-0.5 0.875 0 +-0.375 0.875 0 +-0.5 1 0 +-0.375 1 0 +-0.375 0.875 0 +-0.25 0.875 0 +-0.375 1 0 +-0.25 1 0 +-0.25 0.75 0 +-0.125 0.75 0 +-0.25 0.875 0 +-0.125 0.875 0 +-0.125 0.75 0 +0 0.75 0 +-0.125 0.875 0 +0 0.875 0 +-0.25 0.875 0 +-0.125 0.875 0 +-0.25 1 0 +-0.125 1 0 +-0.125 0.875 0 +0 0.875 0 +-0.125 1 0 +0 1 0 +0 0 0 +0.125 0 0 +0 0.125 0 +0.125 0.125 0 +0.125 0 0 +0.25 0 0 +0.125 0.125 0 +0.25 0.125 0 +0 0.125 0 +0.125 0.125 0 +0 0.25 0 +0.125 0.25 0 +0.125 0.125 0 +0.25 0.125 0 +0.125 0.25 0 +0.25 0.25 0 +0.25 0 0 +0.375 0 0 +0.25 0.125 0 +0.375 0.125 0 +0.375 0 0 +0.5 0 0 +0.375 0.125 0 +0.5 0.125 0 +0.25 0.125 0 +0.375 0.125 0 +0.25 0.25 0 +0.375 0.25 0 +0.375 0.125 0 +0.5 0.125 0 +0.375 0.25 0 +0.5 0.25 0 +0 0.25 0 +0.125 0.25 0 +0 0.375 0 +0.125 0.375 0 +0.125 0.25 0 +0.25 0.25 0 +0.125 0.375 0 +0.25 0.375 0 +0 0.375 0 +0.125 0.375 0 +0 0.5 0 +0.125 0.5 0 +0.125 0.375 0 +0.25 0.375 0 +0.125 0.5 0 +0.25 0.5 0 +0.25 0.25 0 +0.375 0.25 0 +0.25 0.375 0 +0.375 0.375 0 +0.375 0.25 0 +0.5 0.25 0 +0.375 0.375 0 +0.5 0.375 0 +0.25 0.375 0 +0.375 0.375 0 +0.25 0.5 0 +0.375 0.5 0 +0.375 0.375 0 +0.5 0.375 0 +0.375 0.5 0 +0.5 0.5 0 +0.5 0 0 +0.625 0 0 +0.5 0.125 0 +0.625 0.125 0 +0.625 0 0 +0.75 0 0 +0.625 0.125 0 +0.75 0.125 0 +0.5 0.125 0 +0.625 0.125 0 +0.5 0.25 0 +0.625 0.25 0 +0.625 0.125 0 +0.75 0.125 0 +0.625 0.25 0 +0.75 0.25 0 +0.75 0 0 +0.875 0 0 +0.75 0.125 0 +0.875 0.125 0 +0.875 0 0 +1 0 0 +0.875 0.125 0 +1 0.125 0 +0.75 0.125 0 +0.875 0.125 0 +0.75 0.25 0 +0.875 0.25 0 +0.875 0.125 0 +1 0.125 0 +0.875 0.25 0 +1 0.25 0 +0.5 0.25 0 +0.625 0.25 0 +0.5 0.375 0 +0.625 0.375 0 +0.625 0.25 0 +0.75 0.25 0 +0.625 0.375 0 +0.75 0.375 0 +0.5 0.375 0 +0.625 0.375 0 +0.5 0.5 0 +0.625 0.5 0 +0.625 0.375 0 +0.75 0.375 0 +0.625 0.5 0 +0.75 0.5 0 +0.75 0.25 0 +0.875 0.25 0 +0.75 0.375 0 +0.875 0.375 0 +0.875 0.25 0 +1 0.25 0 +0.875 0.375 0 +1 0.375 0 +0.75 0.375 0 +0.875 0.375 0 +0.75 0.5 0 +0.875 0.5 0 +0.875 0.375 0 +1 0.375 0 +0.875 0.5 0 +1 0.5 0 +0 0.5 0 +0.125 0.5 0 +0 0.625 0 +0.125 0.625 0 +0.125 0.5 0 +0.25 0.5 0 +0.125 0.625 0 +0.25 0.625 0 +0 0.625 0 +0.125 0.625 0 +0 0.75 0 +0.125 0.75 0 +0.125 0.625 0 +0.25 0.625 0 +0.125 0.75 0 +0.25 0.75 0 +0.25 0.5 0 +0.375 0.5 0 +0.25 0.625 0 +0.375 0.625 0 +0.375 0.5 0 +0.5 0.5 0 +0.375 0.625 0 +0.5 0.625 0 +0.25 0.625 0 +0.375 0.625 0 +0.25 0.75 0 +0.375 0.75 0 +0.375 0.625 0 +0.5 0.625 0 +0.375 0.75 0 +0.5 0.75 0 +0 0.75 0 +0.125 0.75 0 +0 0.875 0 +0.125 0.875 0 +0.125 0.75 0 +0.25 0.75 0 +0.125 0.875 0 +0.25 0.875 0 +0 0.875 0 +0.125 0.875 0 +0 1 0 +0.125 1 0 +0.125 0.875 0 +0.25 0.875 0 +0.125 1 0 +0.25 1 0 +0.25 0.75 0 +0.375 0.75 0 +0.25 0.875 0 +0.375 0.875 0 +0.375 0.75 0 +0.5 0.75 0 +0.375 0.875 0 +0.5 0.875 0 +0.25 0.875 0 +0.375 0.875 0 +0.25 1 0 +0.375 1 0 +0.375 0.875 0 +0.5 0.875 0 +0.375 1 0 +0.5 1 0 +0.5 0.5 0 +0.625 0.5 0 +0.5 0.625 0 +0.625 0.625 0 +0.625 0.5 0 +0.75 0.5 0 +0.625 0.625 0 +0.75 0.625 0 +0.5 0.625 0 +0.625 0.625 0 +0.5 0.75 0 +0.625 0.75 0 +0.625 0.625 0 +0.75 0.625 0 +0.625 0.75 0 +0.75 0.75 0 +0.75 0.5 0 +0.875 0.5 0 +0.75 0.625 0 +0.875 0.625 0 +0.875 0.5 0 +1 0.5 0 +0.875 0.625 0 +1 0.625 0 +0.75 0.625 0 +0.875 0.625 0 +0.75 0.75 0 +0.875 0.75 0 +0.875 0.625 0 +1 0.625 0 +0.875 0.75 0 +1 0.75 0 +0.5 0.75 0 +0.625 0.75 0 +0.5 0.875 0 +0.625 0.875 0 +0.625 0.75 0 +0.75 0.75 0 +0.625 0.875 0 +0.75 0.875 0 +0.5 0.875 0 +0.625 0.875 0 +0.5 1 0 +0.625 1 0 +0.625 0.875 0 +0.75 0.875 0 +0.625 1 0 +0.75 1 0 +0.75 0.75 0 +0.875 0.75 0 +0.75 0.875 0 +0.875 0.875 0 +0.875 0.75 0 +1 0.75 0 +0.875 0.875 0 +1 0.875 0 +0.75 0.875 0 +0.875 0.875 0 +0.75 1 0 +0.875 1 0 +0.875 0.875 0 +1 0.875 0 +0.875 1 0 +1 1 0 + +CELLS 256 1280 +4 0 1 3 2 +4 4 5 7 6 +4 8 9 11 10 +4 12 13 15 14 +4 16 17 19 18 +4 20 21 23 22 +4 24 25 27 26 +4 28 29 31 30 +4 32 33 35 34 +4 36 37 39 38 +4 40 41 43 42 +4 44 45 47 46 +4 48 49 51 50 +4 52 53 55 54 +4 56 57 59 58 +4 60 61 63 62 +4 64 65 67 66 +4 68 69 71 70 +4 72 73 75 74 +4 76 77 79 78 +4 80 81 83 82 +4 84 85 87 86 +4 88 89 91 90 +4 92 93 95 94 +4 96 97 99 98 +4 100 101 103 102 +4 104 105 107 106 +4 108 109 111 110 +4 112 113 115 114 +4 116 117 119 118 +4 120 121 123 122 +4 124 125 127 126 +4 128 129 131 130 +4 132 133 135 134 +4 136 137 139 138 +4 140 141 143 142 +4 144 145 147 146 +4 148 149 151 150 +4 152 153 155 154 +4 156 157 159 158 +4 160 161 163 162 +4 164 165 167 166 +4 168 169 171 170 +4 172 173 175 174 +4 176 177 179 178 +4 180 181 183 182 +4 184 185 187 186 +4 188 189 191 190 +4 192 193 195 194 +4 196 197 199 198 +4 200 201 203 202 +4 204 205 207 206 +4 208 209 211 210 +4 212 213 215 214 +4 216 217 219 218 +4 220 221 223 222 +4 224 225 227 226 +4 228 229 231 230 +4 232 233 235 234 +4 236 237 239 238 +4 240 241 243 242 +4 244 245 247 246 +4 248 249 251 250 +4 252 253 255 254 +4 256 257 259 258 +4 260 261 263 262 +4 264 265 267 266 +4 268 269 271 270 +4 272 273 275 274 +4 276 277 279 278 +4 280 281 283 282 +4 284 285 287 286 +4 288 289 291 290 +4 292 293 295 294 +4 296 297 299 298 +4 300 301 303 302 +4 304 305 307 306 +4 308 309 311 310 +4 312 313 315 314 +4 316 317 319 318 +4 320 321 323 322 +4 324 325 327 326 +4 328 329 331 330 +4 332 333 335 334 +4 336 337 339 338 +4 340 341 343 342 +4 344 345 347 346 +4 348 349 351 350 +4 352 353 355 354 +4 356 357 359 358 +4 360 361 363 362 +4 364 365 367 366 +4 368 369 371 370 +4 372 373 375 374 +4 376 377 379 378 +4 380 381 383 382 +4 384 385 387 386 +4 388 389 391 390 +4 392 393 395 394 +4 396 397 399 398 +4 400 401 403 402 +4 404 405 407 406 +4 408 409 411 410 +4 412 413 415 414 +4 416 417 419 418 +4 420 421 423 422 +4 424 425 427 426 +4 428 429 431 430 +4 432 433 435 434 +4 436 437 439 438 +4 440 441 443 442 +4 444 445 447 446 +4 448 449 451 450 +4 452 453 455 454 +4 456 457 459 458 +4 460 461 463 462 +4 464 465 467 466 +4 468 469 471 470 +4 472 473 475 474 +4 476 477 479 478 +4 480 481 483 482 +4 484 485 487 486 +4 488 489 491 490 +4 492 493 495 494 +4 496 497 499 498 +4 500 501 503 502 +4 504 505 507 506 +4 508 509 511 510 +4 512 513 515 514 +4 516 517 519 518 +4 520 521 523 522 +4 524 525 527 526 +4 528 529 531 530 +4 532 533 535 534 +4 536 537 539 538 +4 540 541 543 542 +4 544 545 547 546 +4 548 549 551 550 +4 552 553 555 554 +4 556 557 559 558 +4 560 561 563 562 +4 564 565 567 566 +4 568 569 571 570 +4 572 573 575 574 +4 576 577 579 578 +4 580 581 583 582 +4 584 585 587 586 +4 588 589 591 590 +4 592 593 595 594 +4 596 597 599 598 +4 600 601 603 602 +4 604 605 607 606 +4 608 609 611 610 +4 612 613 615 614 +4 616 617 619 618 +4 620 621 623 622 +4 624 625 627 626 +4 628 629 631 630 +4 632 633 635 634 +4 636 637 639 638 +4 640 641 643 642 +4 644 645 647 646 +4 648 649 651 650 +4 652 653 655 654 +4 656 657 659 658 +4 660 661 663 662 +4 664 665 667 666 +4 668 669 671 670 +4 672 673 675 674 +4 676 677 679 678 +4 680 681 683 682 +4 684 685 687 686 +4 688 689 691 690 +4 692 693 695 694 +4 696 697 699 698 +4 700 701 703 702 +4 704 705 707 706 +4 708 709 711 710 +4 712 713 715 714 +4 716 717 719 718 +4 720 721 723 722 +4 724 725 727 726 +4 728 729 731 730 +4 732 733 735 734 +4 736 737 739 738 +4 740 741 743 742 +4 744 745 747 746 +4 748 749 751 750 +4 752 753 755 754 +4 756 757 759 758 +4 760 761 763 762 +4 764 765 767 766 +4 768 769 771 770 +4 772 773 775 774 +4 776 777 779 778 +4 780 781 783 782 +4 784 785 787 786 +4 788 789 791 790 +4 792 793 795 794 +4 796 797 799 798 +4 800 801 803 802 +4 804 805 807 806 +4 808 809 811 810 +4 812 813 815 814 +4 816 817 819 818 +4 820 821 823 822 +4 824 825 827 826 +4 828 829 831 830 +4 832 833 835 834 +4 836 837 839 838 +4 840 841 843 842 +4 844 845 847 846 +4 848 849 851 850 +4 852 853 855 854 +4 856 857 859 858 +4 860 861 863 862 +4 864 865 867 866 +4 868 869 871 870 +4 872 873 875 874 +4 876 877 879 878 +4 880 881 883 882 +4 884 885 887 886 +4 888 889 891 890 +4 892 893 895 894 +4 896 897 899 898 +4 900 901 903 902 +4 904 905 907 906 +4 908 909 911 910 +4 912 913 915 914 +4 916 917 919 918 +4 920 921 923 922 +4 924 925 927 926 +4 928 929 931 930 +4 932 933 935 934 +4 936 937 939 938 +4 940 941 943 942 +4 944 945 947 946 +4 948 949 951 950 +4 952 953 955 954 +4 956 957 959 958 +4 960 961 963 962 +4 964 965 967 966 +4 968 969 971 970 +4 972 973 975 974 +4 976 977 979 978 +4 980 981 983 982 +4 984 985 987 986 +4 988 989 991 990 +4 992 993 995 994 +4 996 997 999 998 +4 1000 1001 1003 1002 +4 1004 1005 1007 1006 +4 1008 1009 1011 1010 +4 1012 1013 1015 1014 +4 1016 1017 1019 1018 +4 1020 1021 1023 1022 + +CELL_TYPES 256 + 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 +POINT_DATA 1024 +SCALARS solution double 1 +LOOKUP_TABLE default +2 1.76562 1.76562 1.6963 1.76562 1.5625 1.6963 1.57867 1.76562 1.6963 1.5625 1.57812 1.6963 1.57867 1.57812 1.53421 1.5625 1.39062 1.57867 1.46147 1.39062 1.25 1.46147 1.36018 1.57867 1.46147 1.53421 1.46965 1.46147 1.36018 1.46965 1.40758 1.5625 1.57812 1.39062 1.45926 1.57812 1.53421 1.45926 1.46687 1.39062 1.45926 1.25 1.35455 1.45926 1.46687 1.35455 1.39852 1.53421 1.46965 1.46687 1.44589 1.46965 1.40758 1.44589 1.41861 1.46687 1.44589 1.39852 1.41061 1.44589 1.41861 1.41061 1.41165 1.25 1.14062 1.36018 1.28178 1.14062 1.0625 1.28178 1.22982 1.36018 1.28178 1.40758 1.35999 1.28178 1.22982 1.35999 1.33367 1.0625 1.01562 1.22982 1.2064 1.01562 1 1.2064 1.21316 1.22982 1.2064 1.33367 1.33295 1.2064 1.21316 1.33295 1.36123 1.40758 1.35999 1.41861 1.39919 1.35999 1.33367 1.39919 1.39652 1.41861 1.39919 1.41165 1.41606 1.39919 1.39652 1.41606 1.43377 1.33367 1.33295 1.39652 1.4167 1.33295 1.36123 1.4167 1.46473 1.39652 1.4167 1.43377 1.47212 1.4167 1.46473 1.47212 1.53733 1.25 1.35455 1.14062 1.27027 1.35455 1.39852 1.27027 1.33975 1.14062 1.27027 1.0625 1.20912 1.27027 1.33975 1.20912 1.29569 1.39852 1.41061 1.33975 1.37588 1.41061 1.41165 1.37588 1.39837 1.33975 1.37588 1.29569 1.34813 1.37588 1.39837 1.34813 1.38554 1.0625 1.20912 1.01562 1.17219 1.20912 1.29569 1.17219 1.26864 1.01562 1.17219 1 1.15985 1.17219 1.26864 1.15985 1.25955 1.29569 1.34813 1.26864 1.33057 1.34813 1.38554 1.33057 1.37675 1.26864 1.33057 1.25955 1.32459 1.33057 1.37675 1.32459 1.37366 1.41165 1.41606 1.39837 1.42145 1.41606 1.43377 1.42145 1.45548 1.39837 1.42145 1.38554 1.42191 1.42145 1.45548 1.42191 1.46776 1.43377 1.47212 1.45548 1.50854 1.47212 1.53733 1.50854 1.58772 1.45548 1.50854 1.46776 1.53159 1.50854 1.58772 1.53159 1.62106 1.38554 1.42191 1.37675 1.42104 1.42191 1.46776 1.42104 1.47403 1.37675 1.42104 1.37366 1.42059 1.42104 1.47403 1.42059 1.47596 1.46776 1.53159 1.47403 1.54442 1.53159 1.62106 1.54442 1.64017 1.47403 1.54442 1.47596 1.54855 1.54442 1.64017 1.54855 1.64641 1 1.01562 1.21316 1.2518 1.01562 1.0625 1.2518 1.32426 1.21316 1.2518 1.36123 1.42186 1.2518 1.32426 1.42186 1.51859 1.0625 1.14062 1.32426 1.43283 1.14062 1.25 1.43283 1.57989 1.32426 1.43283 1.51859 1.6556 1.43283 1.57989 1.6556 1.8373 1.36123 1.42186 1.46473 1.54546 1.42186 1.51859 1.54546 1.6642 1.46473 1.54546 1.53733 1.63552 1.54546 1.6642 1.63552 1.7733 1.51859 1.6556 1.6642 1.82684 1.6556 1.8373 1.82684 2.03968 1.6642 1.82684 1.7733 1.95803 1.82684 2.03968 1.95803 2.19774 1.25 1.39062 1.57989 1.76732 1.39062 1.5625 1.76732 1.99471 1.57989 1.76732 1.8373 2.0675 1.76732 1.99471 2.0675 2.34742 1.5625 1.76562 1.99471 2.25323 1.76562 2 2.25323 2.49069 1.99471 2.25323 2.34742 2.67101 2.25323 2.49069 2.67101 3.01884 1.8373 2.0675 2.03968 2.30883 2.0675 2.34742 2.30883 2.63931 2.03968 2.30883 2.19774 2.50081 2.30883 2.63931 2.50081 2.8758 2.34742 2.67101 2.63931 3.03526 2.67101 3.01884 3.03526 3.5131 2.63931 3.03526 2.8758 3.33261 3.03526 3.5131 3.33261 3.88415 1.53733 1.63552 1.58772 1.7001 1.63552 1.7733 1.7001 1.85333 1.58772 1.7001 1.62106 1.74389 1.7001 1.85333 1.74389 1.9085 1.7733 1.95803 1.85333 2.05595 1.95803 2.19774 2.05595 2.31738 1.85333 2.05595 1.9085 2.12427 2.05595 2.31738 2.12427 2.40156 1.62106 1.74389 1.64017 1.7694 1.74389 1.9085 1.7694 1.94097 1.64017 1.7694 1.64641 1.77779 1.7694 1.94097 1.77779 1.95171 1.9085 2.12427 1.94097 2.16476 2.12427 2.40156 2.16476 2.45166 1.94097 2.16476 1.95171 2.17819 2.16476 2.45166 2.17819 2.46831 2.19774 2.50081 2.31738 2.64769 2.50081 2.8758 2.64769 3.0574 2.31738 2.64769 2.40156 2.75147 2.64769 3.0574 2.75147 3.1853 2.8758 3.33261 3.0574 3.557 3.33261 3.88415 3.557 4.15312 3.0574 3.557 3.1853 3.71346 3.557 4.15312 3.71346 4.34499 2.40156 2.75147 2.45166 2.8133 2.75147 3.1853 2.8133 3.26129 2.45166 2.8133 2.46831 2.83384 2.8133 3.26129 2.83384 3.28651 3.1853 3.71346 3.26129 3.80619 3.71346 4.34499 3.80619 4.45638 3.26129 3.80619 3.28651 3.83687 3.80619 4.45638 3.83687 4.49129 1 1.15985 1.01562 1.17219 1.15985 1.25955 1.17219 1.26864 1.01562 1.17219 1.0625 1.20912 1.17219 1.26864 1.20912 1.29569 1.25955 1.32459 1.26864 1.33057 1.32459 1.37366 1.33057 1.37675 1.26864 1.33057 1.29569 1.34813 1.33057 1.37675 1.34813 1.38554 1.0625 1.20912 1.14062 1.27027 1.20912 1.29569 1.27027 1.33975 1.14062 1.27027 1.25 1.35455 1.27027 1.33975 1.35455 1.39852 1.29569 1.34813 1.33975 1.37588 1.34813 1.38554 1.37588 1.39837 1.33975 1.37588 1.39852 1.41061 1.37588 1.39837 1.41061 1.41165 1.37366 1.42059 1.37675 1.42104 1.42059 1.47596 1.42104 1.47403 1.37675 1.42104 1.38554 1.42191 1.42104 1.47403 1.42191 1.46776 1.47596 1.54855 1.47403 1.54442 1.54855 1.64641 1.54442 1.64017 1.47403 1.54442 1.46776 1.53159 1.54442 1.64017 1.53159 1.62106 1.38554 1.42191 1.39837 1.42145 1.42191 1.46776 1.42145 1.45548 1.39837 1.42145 1.41165 1.41606 1.42145 1.45548 1.41606 1.43377 1.46776 1.53159 1.45548 1.50854 1.53159 1.62106 1.50854 1.58772 1.45548 1.50854 1.43377 1.47212 1.50854 1.58772 1.47212 1.53733 1.25 1.35455 1.39062 1.45926 1.35455 1.39852 1.45926 1.46687 1.39062 1.45926 1.5625 1.57812 1.45926 1.46687 1.57812 1.53421 1.39852 1.41061 1.46687 1.44589 1.41061 1.41165 1.44589 1.41861 1.46687 1.44589 1.53421 1.46965 1.44589 1.41861 1.46965 1.40758 1.5625 1.57812 1.76562 1.6963 1.57812 1.53421 1.6963 1.57867 1.76562 1.6963 2 1.76562 1.6963 1.57867 1.76562 1.5625 1.53421 1.46965 1.57867 1.46147 1.46965 1.40758 1.46147 1.36018 1.57867 1.46147 1.5625 1.39062 1.46147 1.36018 1.39062 1.25 1.41165 1.41606 1.41861 1.39919 1.41606 1.43377 1.39919 1.39652 1.41861 1.39919 1.40758 1.35999 1.39919 1.39652 1.35999 1.33367 1.43377 1.47212 1.39652 1.4167 1.47212 1.53733 1.4167 1.46473 1.39652 1.4167 1.33367 1.33295 1.4167 1.46473 1.33295 1.36123 1.40758 1.35999 1.36018 1.28178 1.35999 1.33367 1.28178 1.22982 1.36018 1.28178 1.25 1.14062 1.28178 1.22982 1.14062 1.0625 1.33367 1.33295 1.22982 1.2064 1.33295 1.36123 1.2064 1.21316 1.22982 1.2064 1.0625 1.01562 1.2064 1.21316 1.01562 1 1.64641 1.77779 1.64017 1.7694 1.77779 1.95171 1.7694 1.94097 1.64017 1.7694 1.62106 1.74389 1.7694 1.94097 1.74389 1.9085 1.95171 2.17819 1.94097 2.16476 2.17819 2.46831 2.16476 2.45166 1.94097 2.16476 1.9085 2.12427 2.16476 2.45166 2.12427 2.40156 1.62106 1.74389 1.58772 1.7001 1.74389 1.9085 1.7001 1.85333 1.58772 1.7001 1.53733 1.63552 1.7001 1.85333 1.63552 1.7733 1.9085 2.12427 1.85333 2.05595 2.12427 2.40156 2.05595 2.31738 1.85333 2.05595 1.7733 1.95803 2.05595 2.31738 1.95803 2.19774 2.46831 2.83384 2.45166 2.8133 2.83384 3.28651 2.8133 3.26129 2.45166 2.8133 2.40156 2.75147 2.8133 3.26129 2.75147 3.1853 3.28651 3.83687 3.26129 3.80619 3.83687 4.49129 3.80619 4.45638 3.26129 3.80619 3.1853 3.71346 3.80619 4.45638 3.71346 4.34499 2.40156 2.75147 2.31738 2.64769 2.75147 3.1853 2.64769 3.0574 2.31738 2.64769 2.19774 2.50081 2.64769 3.0574 2.50081 2.8758 3.1853 3.71346 3.0574 3.557 3.71346 4.34499 3.557 4.15312 3.0574 3.557 2.8758 3.33261 3.557 4.15312 3.33261 3.88415 1.53733 1.63552 1.46473 1.54546 1.63552 1.7733 1.54546 1.6642 1.46473 1.54546 1.36123 1.42186 1.54546 1.6642 1.42186 1.51859 1.7733 1.95803 1.6642 1.82684 1.95803 2.19774 1.82684 2.03968 1.6642 1.82684 1.51859 1.6556 1.82684 2.03968 1.6556 1.8373 1.36123 1.42186 1.21316 1.2518 1.42186 1.51859 1.2518 1.32426 1.21316 1.2518 1 1.01562 1.2518 1.32426 1.01562 1.0625 1.51859 1.6556 1.32426 1.43283 1.6556 1.8373 1.43283 1.57989 1.32426 1.43283 1.0625 1.14062 1.43283 1.57989 1.14062 1.25 2.19774 2.50081 2.03968 2.30883 2.50081 2.8758 2.30883 2.63931 2.03968 2.30883 1.8373 2.0675 2.30883 2.63931 2.0675 2.34742 2.8758 3.33261 2.63931 3.03526 3.33261 3.88415 3.03526 3.5131 2.63931 3.03526 2.34742 2.67101 3.03526 3.5131 2.67101 3.01884 1.8373 2.0675 1.57989 1.76732 2.0675 2.34742 1.76732 1.99471 1.57989 1.76732 1.25 1.39062 1.76732 1.99471 1.39062 1.5625 2.34742 2.67101 1.99471 2.25323 2.67101 3.01884 2.25323 2.49069 1.99471 2.25323 1.5625 1.76562 2.25323 2.49069 1.76562 2