Dendro-5.01

nlsm

For questions: Gitter

What is Dendro ?

“Dendro” in Greek language means tree. The Dendro library is a large scale (262K cores on ORNL’s Titan) distributed memory adaptive octree framework. The main goal of Dendro is to perform large scale multiphysics simulations efficeiently in mordern supercomputers. Dendro consists of efficient parallel data structures and algorithms to perform variational ( finite element) methods and finite difference mthods on 2:1 balanced arbitary adaptive octrees which enables the users to perform simulations raning from black holes (binary black hole mergers) to blood flow in human body, where applications ranging from relativity, astrophysics to biomedical engineering.


Get Dendro (git repo )

You can clone the repository using , git clone https://github.com/paralab/Dendro-5.01.git

Dendro-5.01 Documentation

Doxygen documentation can be found here

How to build Dendro-5.0 ?

You need CMake to build dendro. Create a build directory using ‘mkdir build’. Then go into the build directory by ‘cd build’ then execute ‘ccmake ..’ to generate the make files. You can build Dendro-5.0 with several options.

Simple simulation: Nonlinear Sigma Model (NLSigma)

NlSigma folder consists of simple, non lineat wave equation with adaptive mesh refinement (AMR). You can copy the parameter file from NLSigma/par folder and simply run mpirun -np 8 ./NLSigma/nlsmSolver nlsm.par.json, on your lattop to large supercomputer with higher resolution.

nlsm nlsm nlsm nlsm

You can write the equations in symbolic python which generate the C compute kernel. Look at nlsm.py

import dendro
from sympy import *
###############################################################
#  initialize
###############################################################
r = symbols('r')
# declare functions
chi = dendro.scalar("chi","[pp]")
phi = dendro.scalar("phi","[pp]")
d = dendro.set_first_derivative('grad')    # first argument is direction
d2s = dendro.set_second_derivative('grad2')  # first 2 arguments are directions
d2 = dendro.d2

###############################################################
#  evolution equations
###############################################################

phi_rhs = sum( d2(i,i,chi)  for i in dendro.e_i ) - sin(2*chi)/r**2
chi_rhs = phi

###############################################################
#  evolution equations
###############################################################
outs = [phi_rhs, chi_rhs]
vnames = ['phi_rhs', 'chi_rhs']
dendro.generate(outs, vnames, '[pp]')

Which generate the code,

// Dendro: original ops:  10
// Dendro: printing temp variables
// Dendro: printing variables
//--
phi_rhs[pp] = grad2_0_0_chi[pp] + grad2_1_1_chi[pp] + grad2_2_2_chi[pp] - sin(2*chi[pp])/pow(r, 2);
//--
chi_rhs[pp] = phi[pp];
// Dendro: reduced ops:  10

Parameters for NLSigma


Scalability on octree generation and partitioning.

We have performed octree generation and partitioning up to 262144 cores in ORNL’s titan super computer. We have managed to partition 1.3x10^12 octants among 262144 processors with in 4 seconds.

weak scaling on octrees strong scaling on octrees

Publications