1937 lines
77 KiB
Python
1937 lines
77 KiB
Python
################################################################
|
|
#
|
|
# Workspace for MRI analysis of the 4Dflow data
|
|
#
|
|
# written by Jeremias Garay L: j.e.garay.labra@rug.nl
|
|
# Fernanda te amo
|
|
# for autoreload in ipython3
|
|
# %load_ext autoreload
|
|
# %autoreload 2
|
|
################################################################
|
|
import h5py
|
|
from dPdirectestim.dPdirectestim import *
|
|
from dolfin import *
|
|
import dolfin
|
|
import numpy as np
|
|
import sys, os
|
|
from mpi4py import MPI
|
|
from common import inout
|
|
import time
|
|
|
|
if '2017' in dolfin.__version__:
|
|
#class MPI(MPI):
|
|
# comm_world = mpi_comm_world()
|
|
set_log_active(False)
|
|
else:
|
|
PRINT_BARS = False
|
|
parameters["std_out_all_processes"] = False
|
|
|
|
def Etcetera(mode):
|
|
if mode=='happy':
|
|
print(r"""
|
|
/\ /\
|
|
//\\_//\\ ____
|
|
\_ _/ / /
|
|
/ ^ ^ \ /^^^]
|
|
\_\O/_/ [ ]
|
|
/ \_ [ /
|
|
\ \_ / /
|
|
[ [ / \/ _/
|
|
_[ [ \ /_/
|
|
""")
|
|
|
|
if mode=='serious':
|
|
print(r"""
|
|
|\__/|
|
|
/ \
|
|
/_.- -,_\
|
|
\@/
|
|
""")
|
|
|
|
def dealiased(VENC,vx,vy,vz):
|
|
|
|
vx2 = np.zeros(vx.shape)
|
|
vy2 = np.zeros(vy.shape)
|
|
vz2 = np.zeros(vz.shape)
|
|
|
|
badx_pos = np.where(vx>VENC)
|
|
badx_neg = np.where(vx<-VENC)
|
|
bady_pos = np.where(vy>VENC)
|
|
bady_neg = np.where(vy<-VENC)
|
|
badz_pos = np.where(vz>VENC)
|
|
badz_neg = np.where(vz<-VENC)
|
|
# More than VENC
|
|
vx2[badx_pos] = -(2*VENC - vx[badx_pos])
|
|
vy2[bady_pos] = -(2*VENC - vy[bady_pos])
|
|
vz2[badz_pos] = -(2*VENC - vz[badz_pos])
|
|
# Less than VENC
|
|
vx2[badx_neg] = 2*VENC + vx[badx_neg]
|
|
vy2[bady_neg] = 2*VENC + vy[bady_neg]
|
|
vz2[badz_neg] = 2*VENC + vz[badz_neg]
|
|
# The rest
|
|
otherx = np.where(vx2==0)
|
|
othery = np.where(vy2==0)
|
|
otherz = np.where(vz2==0)
|
|
vx2[otherx] = vx[otherx]
|
|
vy2[othery] = vy[othery]
|
|
vz2[otherz] = vz[otherz]
|
|
|
|
return [vx2,vy2,vz2]
|
|
|
|
def checkpoint_norm(MESH,fieldpath,mode):
|
|
V = MESH['FEM']
|
|
v = Function(V)
|
|
path = fieldpath + 'checkpoint/'
|
|
# Checkpoints folders
|
|
unsort_indexes = os.listdir(path)
|
|
indexes = [int(x) for x in unsort_indexes]
|
|
indexes.sort()
|
|
norm = np.zeros([len(indexes)])
|
|
j = 0
|
|
comm = MPI.COMM_WORLD
|
|
|
|
if mode=='ct':
|
|
# CT DATA
|
|
for l in indexes:
|
|
path1 = path + str(l) + '/u.h5'
|
|
hdf = HDF5File(MESH['mesh'].mpi_comm(),path1,'r')
|
|
hdf.read(v, 'u/vector_0')
|
|
hdf.close()
|
|
norm[j] = np.linalg.norm(v.vector().get_local())
|
|
j = j+1
|
|
|
|
if mode=='measurement':
|
|
# MEASUREMENT
|
|
path2 = fieldpath + 'measurements/'
|
|
j = 0
|
|
for l in indexes:
|
|
meapath = path2 + 'u' + str(l) + '.h5'
|
|
comm = MPI.COMM_WORLD
|
|
hdf = HDF5File(MESH['mesh'].mpi_comm(),meapath,'r')
|
|
hdf.read(v, 'u/vector_0')
|
|
hdf.close()
|
|
norm[j] = np.linalg.norm(v.vector().get_local())
|
|
j = j+1
|
|
|
|
|
|
if mode=='roukf':
|
|
indexes2 = indexes[1:-1]
|
|
j = 1
|
|
for l in indexes2:
|
|
path3 = fieldpath + 'roukf/checkpoint/' + str(l) + '/X0.h5'
|
|
hdf = HDF5File(MESH['mesh'].mpi_comm(),path3,'r')
|
|
hdf.read(v, 'X')
|
|
hdf.close()
|
|
norm[j] = np.linalg.norm(v.vector().get_local())
|
|
j = j+1
|
|
|
|
|
|
return norm
|
|
|
|
def LOADmesh(meshtype,resol=None,boxsize=None,meshpath=None, ranges=None):
|
|
|
|
if meshtype=='box_zeros':
|
|
# The size of the mesh is given
|
|
# The spacements are given
|
|
# The x range goes from 0
|
|
[Lx, Ly, Lz] = boxsize
|
|
[ddx, ddy, ddz] = resol
|
|
x0 = 0
|
|
y0 = 0
|
|
z0 = 0
|
|
x1 = (Lx-1)*ddx
|
|
y1 = (Ly-1)*ddy
|
|
z1 = (Lz-1)*ddz
|
|
mesh0 = BoxMesh(Point(x0, y0, z0), Point(x1, y1, z1), Lx-1, Ly-1, Lz-1)
|
|
P1 = FunctionSpace(mesh0, "Lagrange", 1)
|
|
# This need to be saved from Leo's construction
|
|
Xt = P1.tabulate_dof_coordinates()
|
|
#boundaries = MeshFunction('size_t', mesh0 , mesh0.topology().dim() - 1)
|
|
#XDMFFile(path2).write(boundaries)
|
|
if '2017' in dolfin.__version__:
|
|
Xt = np.reshape(Xt, (round(Xt.size/3), 3))
|
|
SqXt = np.zeros(Xt.shape)
|
|
for k in range(Xt.shape[0]):
|
|
xk = round((Xt[k][0]-x0)/ddx)
|
|
yk = round((Xt[k][1]-y0)/ddy)
|
|
zk = round((Xt[k][2]-z0)/ddz)
|
|
SqXt[k, 0] = int(xk)
|
|
SqXt[k, 1] = int(yk)
|
|
SqXt[k, 2] = int(zk)
|
|
if rank == 0:
|
|
print('BOX real created: H = ' + str(mesh0.hmin()) +
|
|
' H_max = ' + str(mesh0.hmax()))
|
|
|
|
BOX = {}
|
|
BOX['mesh'] = mesh0
|
|
BOX['FEM'] = P1
|
|
BOX['nodes'] = Xt
|
|
BOX['seq'] = SqXt.astype(int)
|
|
BOX['dimension'] = [Lx,Ly,Lz]
|
|
return BOX
|
|
|
|
if meshtype=='box_ranges':
|
|
# The size of the mesh is given
|
|
# The X ranges are given
|
|
# The spacements are calculated
|
|
[Lx,Ly,Lz] = boxsize
|
|
x0 = ranges['x'][0]
|
|
x1 = ranges['x'][1]
|
|
y0 = ranges['y'][0]
|
|
y1 = ranges['y'][1]
|
|
z0 = ranges['z'][0]
|
|
z1 = ranges['z'][1]
|
|
ddx = (x1-x0)/(Lx-1)
|
|
ddy = (y1-y0)/(Ly-1)
|
|
ddz = (z1-z0)/(Lz-1)
|
|
mesh0 = BoxMesh(Point(x0, y0, z0), Point(x1, y1, z1), Lx-1, Ly-1, Lz-1)
|
|
P1 = FunctionSpace(mesh0, "Lagrange", 1)
|
|
Xt = P1.tabulate_dof_coordinates()
|
|
if '2017' in dolfin.__version__:
|
|
Xt = np.reshape(Xt, (round(Xt.size/3), 3))
|
|
SqXt = np.zeros(Xt.shape)
|
|
for k in range(Xt.shape[0]):
|
|
xk = round((Xt[k][0]-x0)/ddx)
|
|
yk = round((Xt[k][1]-y0)/ddy)
|
|
zk = round((Xt[k][2]-z0)/ddz)
|
|
SqXt[k, 0] = int(xk)
|
|
SqXt[k, 1] = int(yk)
|
|
SqXt[k, 2] = int(zk)
|
|
if rank == 0:
|
|
print('BOX real created: H = ' + str(mesh0.hmin()) +
|
|
' H_max = ' + str(mesh0.hmax()))
|
|
|
|
BOX = {}
|
|
BOX['mesh'] = mesh0
|
|
BOX['FEM'] = P1
|
|
BOX['nodes'] = Xt
|
|
BOX['seq'] = SqXt.astype(int)
|
|
BOX['dimension'] = [Lx, Ly, Lz]
|
|
return BOX
|
|
|
|
if meshtype=='fix_resolution':
|
|
# The size of the mesh is given
|
|
# The spacements are given
|
|
# The X ranges are given
|
|
# The X ranges are redefined from the oldest
|
|
[Lx,Ly,Lz] = boxsize
|
|
[ddx,ddy,ddz] = resol
|
|
x0 = ranges['x'][0]
|
|
x1 = ranges['x'][1]
|
|
y0 = ranges['y'][0]
|
|
y1 = ranges['y'][1]
|
|
z0 = ranges['z'][0]
|
|
z1 = ranges['z'][1]
|
|
deltax = (ddx*(Lx-1) - (x1 - x0))/2
|
|
deltay = (ddy*(Ly-1) - (y1 - y0))/2
|
|
deltaz = (ddz*(Lz-1) - (z1 - z0))/2
|
|
x0_new = x0 - deltax
|
|
x1_new = x1 + deltax
|
|
y0_new = y0 - deltay
|
|
y1_new = y1 + deltay
|
|
z0_new = z0 - deltaz
|
|
z1_new = z1 + deltaz
|
|
|
|
mesh0 = BoxMesh(Point(x0_new,y0_new,z0_new),Point(x1_new,y1_new,z1_new),Lx-1,Ly-1,Lz-1)
|
|
P1 = FunctionSpace(mesh0, "Lagrange", 1)
|
|
Xt = P1.tabulate_dof_coordinates() # This need to be saved from Leo's construction
|
|
#boundaries = MeshFunction('size_t', mesh0 , mesh0.topology().dim() - 1)
|
|
#XDMFFile('/home/yeye/Desktop/bound.xdmf').write(boundaries)
|
|
if '2017' in dolfin.__version__:
|
|
Xt = np.reshape(Xt,( round(Xt.size/3) ,3))
|
|
SqXt = np.zeros(Xt.shape)
|
|
for k in range(Xt.shape[0]):
|
|
xk = round((Xt[k][0]-x0_new)/ddx)
|
|
yk = round((Xt[k][1]-y0_new)/ddy)
|
|
zk = round((Xt[k][2]-z0_new)/ddz)
|
|
SqXt[k, 0] = int(xk)
|
|
SqXt[k, 1] = int(yk)
|
|
SqXt[k, 2] = int(zk)
|
|
|
|
if rank==0:
|
|
print('BOX real created: H = ' + str(mesh0.hmin()) + ' H_max = ' + str(mesh0.hmax()) )
|
|
|
|
BOX = {}
|
|
BOX['mesh'] = mesh0
|
|
BOX['FEM'] = P1
|
|
BOX['nodes'] = Xt
|
|
BOX['seq'] = SqXt.astype(int)
|
|
BOX['dimension'] = [Lx, Ly, Lz]
|
|
|
|
return BOX
|
|
|
|
if meshtype=='cib':
|
|
ddx = 2.00032/1000; ddy = 2.00032/1000; ddz = 2.00526316/1000
|
|
x0 = ddx*17 ; y0 = ddy*17; z0 = ddz
|
|
Lx = 75; Ly = 85; Lz = 57
|
|
x1 = x0+ddx*(Lx-1); y1 = y0+ddy*(Ly-1); z1 = z0+ddz*(Lz-1)
|
|
|
|
|
|
mesh0 = BoxMesh(Point(x0,y0,z0),Point(x1,y1,z1),Lx-1,Ly-1,Lz-1)
|
|
P1 = FunctionSpace(mesh0, "Lagrange", 1)
|
|
Xt = P1.tabulate_dof_coordinates()
|
|
#XDMFFile('mesh_hex.xdmf').write(mesh0)
|
|
# MESH HEXA
|
|
#mesh0 = BoxMesh.create([Point(x0,y0,z0), Point(x1,y1,z1)], [Lx, Ly, Lz], CellType.Type.hexahedron)
|
|
#P1 = FunctionSpace(mesh0, "DG", 0)
|
|
#t = P1.tabulate_dof_coordinates()
|
|
|
|
if '2017' in dolfin.__version__:
|
|
Xt = np.reshape(Xt,( round(Xt.size/3) ,3))
|
|
SqXt = np.zeros(Xt.shape)
|
|
|
|
for k in range(Xt.shape[0]):
|
|
xk = round((Xt[k][0]-x0)/ddx)
|
|
yk = round((Xt[k][1]-y0)/ddy)
|
|
zk = round((Xt[k][2]-z0)/ddz)
|
|
SqXt[k,0] = int(xk)
|
|
SqXt[k,1] = int(yk)
|
|
SqXt[k,2] = int(zk)
|
|
if rank==0:
|
|
print('BOX cib created: H = ' + str(mesh0.hmin()) + ' H_max = ' + str(mesh0.hmax()) )
|
|
|
|
BOX = {}
|
|
BOX['mesh'] = mesh0
|
|
BOX['FEM'] = P1
|
|
BOX['nodes'] = Xt
|
|
BOX['seq'] = SqXt.astype(int)
|
|
BOX['name'] = 'cib'
|
|
BOX['dimension'] = [Lx,Ly,Lz]
|
|
return BOX
|
|
|
|
if meshtype=='cib_RT':
|
|
ddx = 1.8/1000; ddy = 1.8/1000; ddz = 1.8/1000
|
|
x0 = ddx*33 ; y0 = ddy*43; z0 = ddz*13
|
|
Lx = 103; Ly = 63; Lz = 47
|
|
x1 = x0+ddx*(Lx-1); y1 = y0+ddy*(Ly-1); z1 = z0+ddz*(Lz-1)
|
|
|
|
mesh0 = BoxMesh(Point(x0,y0,z0),Point(x1,y1,z1),Lx-1,Ly-1,Lz-1)
|
|
P1 = FunctionSpace(mesh0, "Lagrange", 1)
|
|
Xt = P1.tabulate_dof_coordinates()
|
|
#XDMFFile('mesh_hex.xdmf').write(mesh0)
|
|
# MESH HEXA
|
|
#mesh0 = BoxMesh.create([Point(x0,y0,z0), Point(x1,y1,z1)], [Lx, Ly, Lz], CellType.Type.hexahedron)
|
|
#P1 = FunctionSpace(mesh0, "DG", 0)
|
|
#t = P1.tabulate_dof_coordinates()
|
|
|
|
if '2017' in dolfin.__version__:
|
|
Xt = np.reshape(Xt,( round(Xt.size/3) ,3))
|
|
SqXt = np.zeros(Xt.shape)
|
|
|
|
for k in range(Xt.shape[0]):
|
|
xk = round((Xt[k][0]-x0)/ddx)
|
|
yk = round((Xt[k][1]-y0)/ddy)
|
|
zk = round((Xt[k][2]-z0)/ddz)
|
|
SqXt[k,0] = int(xk)
|
|
SqXt[k,1] = int(yk)
|
|
SqXt[k,2] = int(zk)
|
|
if rank==0:
|
|
print('BOX cib created: H = ' + str(mesh0.hmin()) + ' H_max = ' + str(mesh0.hmax()) )
|
|
|
|
BOX = {}
|
|
BOX['mesh'] = mesh0
|
|
BOX['FEM'] = P1
|
|
BOX['nodes'] = Xt
|
|
BOX['seq'] = SqXt.astype(int)
|
|
BOX['name'] = 'cib'
|
|
BOX['dimension'] = [Lx,Ly,Lz]
|
|
return BOX
|
|
|
|
if meshtype=='cib_GM':
|
|
ddx = 1.8/1000; ddy = 1.8/1000; ddz = 1.8/1000
|
|
x0 = ddx*20 ; y0 = ddy*43; z0 = -ddz*5
|
|
Lx = 110; Ly = 70; Lz = 43
|
|
x1 = x0+ddx*(Lx-1); y1 = y0+ddy*(Ly-1); z1 = z0+ddz*(Lz-1)
|
|
|
|
mesh0 = BoxMesh(Point(x0,y0,z0),Point(x1,y1,z1),Lx-1,Ly-1,Lz-1)
|
|
P1 = FunctionSpace(mesh0, "Lagrange", 1)
|
|
Xt = P1.tabulate_dof_coordinates()
|
|
#XDMFFile('mesh_hex.xdmf').write(mesh0)
|
|
# MESH HEXA
|
|
#mesh0 = BoxMesh.create([Point(x0,y0,z0), Point(x1,y1,z1)], [Lx, Ly, Lz], CellType.Type.hexahedron)
|
|
#P1 = FunctionSpace(mesh0, "DG", 0)
|
|
#t = P1.tabulate_dof_coordinates()
|
|
|
|
if '2017' in dolfin.__version__:
|
|
Xt = np.reshape(Xt,( round(Xt.size/3) ,3))
|
|
SqXt = np.zeros(Xt.shape)
|
|
|
|
for k in range(Xt.shape[0]):
|
|
xk = round((Xt[k][0]-x0)/ddx)
|
|
yk = round((Xt[k][1]-y0)/ddy)
|
|
zk = round((Xt[k][2]-z0)/ddz)
|
|
SqXt[k,0] = int(xk)
|
|
SqXt[k,1] = int(yk)
|
|
SqXt[k,2] = int(zk)
|
|
if rank==0:
|
|
print('BOX cib created: H = ' + str(mesh0.hmin()) + ' H_max = ' + str(mesh0.hmax()) )
|
|
|
|
BOX = {}
|
|
BOX['mesh'] = mesh0
|
|
BOX['FEM'] = P1
|
|
BOX['nodes'] = Xt
|
|
BOX['seq'] = SqXt.astype(int)
|
|
BOX['name'] = 'cib'
|
|
BOX['dimension'] = [Lx,Ly,Lz]
|
|
return BOX
|
|
|
|
if meshtype=='aorta':
|
|
|
|
if '.xml' in meshpath:
|
|
mesh1 = Mesh(meshpath)
|
|
Evel1 = VectorElement('Lagrange', mesh1.ufl_cell(), 1)
|
|
V1 = FunctionSpace(mesh1, Evel1)
|
|
boundaries1 = MeshFunction(
|
|
"size_t", mesh1, mesh1.topology().dim()-1)
|
|
if '.h5' in meshpath:
|
|
mesh1 = Mesh()
|
|
hdf1 = HDF5File(mesh1.mpi_comm(), meshpath, 'r')
|
|
hdf1.read(mesh1, '/mesh', False)
|
|
boundaries1 = MeshFunction(
|
|
'size_t', mesh1, mesh1.topology().dim() - 1)
|
|
hdf1.read(boundaries1, '/boundaries')
|
|
Evel1 = VectorElement('Lagrange', mesh1.ufl_cell(), 1)
|
|
V1 = FunctionSpace(mesh1, Evel1)
|
|
|
|
if rank == 0:
|
|
print('AORTA created: h_min = ' + str(mesh1.hmin()) +
|
|
' h_max = ' + str(mesh1.hmax()))
|
|
|
|
AORTA = {}
|
|
AORTA['mesh'] = mesh1
|
|
AORTA['boundaries'] = boundaries1
|
|
AORTA['FEM'] = V1
|
|
AORTA['name'] = 'aorta'
|
|
|
|
return AORTA
|
|
|
|
if meshtype == 'leo':
|
|
|
|
# LEO
|
|
mesh = Mesh(meshpath)
|
|
Evel = VectorElement('Lagrange', mesh.ufl_cell(), 1)
|
|
#Evel = VectorElement('DG',mesh2.ufl_cell(),0)
|
|
V = FunctionSpace(mesh, Evel)
|
|
boundaries = MeshFunction("size_t", mesh, mesh.topology().dim()-1)
|
|
marked = True
|
|
|
|
if marked:
|
|
class Outlet(SubDomain):
|
|
def inside(self, x, on_boundary):
|
|
return on_boundary and between(x[2], (11.18, 11.4))
|
|
|
|
class Inlet(SubDomain):
|
|
def inside(self, x, on_boundary):
|
|
return on_boundary and between(x[2], (7.56, 7.8)) and between(x[0], (10, 16))
|
|
|
|
outlet = Outlet()
|
|
inlet = Inlet()
|
|
boundaries.set_all(0)
|
|
outlet.mark(boundaries, 3)
|
|
inlet.mark(boundaries, 2)
|
|
|
|
if rank == 0:
|
|
print('LEO created: H = ' + str(mesh.hmin()) +
|
|
' H_max = ' + str(mesh.hmax()))
|
|
|
|
LEO = {}
|
|
LEO['mesh'] = mesh
|
|
LEO['FEM'] = V
|
|
LEO['boundaries'] = boundaries
|
|
|
|
return LEO
|
|
|
|
if meshtype=='other':
|
|
|
|
if '.xml' in meshpath:
|
|
mesh1 = Mesh(meshpath)
|
|
Evel1 = VectorElement('Lagrange', mesh1.ufl_cell(), 1)
|
|
V1 = FunctionSpace(mesh1, Evel1)
|
|
boundaries1 = MeshFunction(
|
|
"size_t", mesh1, mesh1.topology().dim()-1)
|
|
if '.h5' in meshpath:
|
|
mesh1 = Mesh()
|
|
hdf1 = HDF5File(mesh1.mpi_comm(), meshpath, 'r')
|
|
hdf1.read(mesh1, '/mesh', False)
|
|
boundaries1 = MeshFunction(
|
|
'size_t', mesh1, mesh1.topology().dim() - 1)
|
|
hdf1.read(boundaries1, '/boundaries')
|
|
Evel1 = VectorElement('Lagrange', mesh1.ufl_cell(), 1)
|
|
V1 = FunctionSpace(mesh1, Evel1)
|
|
|
|
if rank == 0:
|
|
print('MESH created: h_min = ' + str(mesh1.hmin()) +
|
|
' h_max = ' + str(mesh1.hmax()))
|
|
|
|
MESH = {}
|
|
MESH['mesh'] = mesh1
|
|
MESH['boundaries'] = boundaries1
|
|
MESH['FEM'] = V1
|
|
|
|
return MESH
|
|
|
|
def CREATEsequences(BOX,AORTA,options):
|
|
# Checkpoints folders
|
|
checkpoint_path = options['sequence']['checkpoint_path']+'checkpoint/'
|
|
unsort_indexes = os.listdir(checkpoint_path)
|
|
indexes = [int(x) for x in unsort_indexes]
|
|
indexes.sort()
|
|
[Lx, Ly, Lz] = BOX['dimension']
|
|
|
|
Lt = len(indexes)
|
|
Sqx = np.zeros([Lx, Ly, Lz, Lt])
|
|
Sqy = np.zeros([Lx, Ly, Lz, Lt])
|
|
Sqz = np.zeros([Lx, Ly, Lz, Lt])
|
|
|
|
v1x = Function(BOX['FEM'])
|
|
v1y = Function(BOX['FEM'])
|
|
v1z = Function(BOX['FEM'])
|
|
v2 = Function(AORTA['FEM'])
|
|
N_sp = int(options['sequence']['sampling_rate'])
|
|
|
|
if rank==0:
|
|
print('Total number of checkpoints folders: ',len(indexes))
|
|
|
|
if options['sequence']['sampling_rate']>1:
|
|
if rank==0:
|
|
print('Applying an undersampling of ' + str(N_sp))
|
|
print('list of indexes',indexes)
|
|
|
|
for k in range(0, len(indexes), N_sp):
|
|
if rank == 0:
|
|
print('READING CHECKPOINT FOLDER: ' + str(indexes[k]))
|
|
|
|
hdf = HDF5File(AORTA['mesh'].mpi_comm(
|
|
), checkpoint_path + str(indexes[k])+'/u.h5', 'r')
|
|
hdf.read(v2, 'u/vector_0')
|
|
hdf.close()
|
|
|
|
LagrangeInterpolator.interpolate(v1x, v2.sub(0))
|
|
LagrangeInterpolator.interpolate(v1y, v2.sub(1))
|
|
LagrangeInterpolator.interpolate(v1z, v2.sub(2))
|
|
|
|
v1_gather = v1x.vector().gather(
|
|
np.array(range(BOX['FEM'].dim()), 'intc'))
|
|
v2_gather = v1y.vector().gather(
|
|
np.array(range(BOX['FEM'].dim()), 'intc'))
|
|
v3_gather = v1z.vector().gather(
|
|
np.array(range(BOX['FEM'].dim()), 'intc'))
|
|
|
|
SqXt = BOX['seq']
|
|
#SqXt = SqXt.astype(int)
|
|
# To Gather the Sorting-Sequence in parallel with Gatherv
|
|
loc_X = SqXt[:, 0]
|
|
loc_Y = SqXt[:, 1]
|
|
loc_Z = SqXt[:, 2]
|
|
sendbuf_X = np.array(loc_X)
|
|
sendbuf_Y = np.array(loc_Y)
|
|
sendbuf_Z = np.array(loc_Z)
|
|
sendcounts_X = np.array(comm.gather(len(sendbuf_X), 0))
|
|
sendcounts_Y = np.array(comm.gather(len(sendbuf_Y), 0))
|
|
sendcounts_Z = np.array(comm.gather(len(sendbuf_Z), 0))
|
|
|
|
if rank == 0:
|
|
SStot_X = np.zeros([Lx*Ly*Lz], dtype=int)
|
|
SStot_Y = np.zeros([Lx*Ly*Lz], dtype=int)
|
|
SStot_Z = np.zeros([Lx*Ly*Lz], dtype=int)
|
|
else:
|
|
SStot_X = None
|
|
SStot_Y = None
|
|
SStot_Z = None
|
|
|
|
comm.Gatherv(sendbuf=sendbuf_X, recvbuf=(
|
|
SStot_X, sendcounts_X), root=0)
|
|
comm.Gatherv(sendbuf=sendbuf_Y, recvbuf=(
|
|
SStot_Y, sendcounts_Y), root=0)
|
|
comm.Gatherv(sendbuf=sendbuf_Z, recvbuf=(
|
|
SStot_Z, sendcounts_Z), root=0)
|
|
|
|
if rank == 0:
|
|
for l in range(Lx*Ly*Lz):
|
|
Sqx[SStot_X[l], SStot_Y[l], SStot_Z[l], k] = v1_gather[l]
|
|
Sqy[SStot_X[l], SStot_Y[l], SStot_Z[l], k] = v2_gather[l]
|
|
Sqz[SStot_X[l], SStot_Y[l], SStot_Z[l], k] = v3_gather[l]
|
|
|
|
if options['sequence']['sampling_rate']>1:
|
|
Sqx = Sqx[:, :, :, 0::N_sp]
|
|
Sqy = Sqy[:, :, :, 0::N_sp]
|
|
Sqz = Sqz[:, :, :, 0::N_sp]
|
|
|
|
if rank == 0:
|
|
print('saving the sequences' + options['sequence']['checkpoint_path'])
|
|
savepath = options['sequence']['checkpoint_path'] + options['sequence']['name'] + '.npz'
|
|
np.savez_compressed(savepath, x=Sqx, y=Sqy, z=Sqz)
|
|
|
|
|
|
return [Sqx,Sqy,Sqz]
|
|
|
|
def LOADsequences(loadpath):
|
|
# for loading existing sequences alocated in loadpath
|
|
if rank==0:
|
|
print('{reading} ' + loadpath)
|
|
|
|
if '.mat' in loadpath:
|
|
import scipy.io as sio
|
|
S = sio.loadmat(loadpath)
|
|
S = S['u_R1']
|
|
Sqx = S['x'][0][0]
|
|
Sqy = S['y'][0][0]
|
|
Sqz = S['z'][0][0]
|
|
else:
|
|
S = np.load(loadpath)
|
|
Sqx = S['x']
|
|
Sqy = S['y']
|
|
Sqz = S['z']
|
|
|
|
return [Sqx,Sqy,Sqz]
|
|
|
|
def SqtoH5(BOX,MESH,Sqx,Sqy,Sqz):
|
|
|
|
Xt = BOX['nodes']
|
|
SqXt = BOX['seq']
|
|
P1 = BOX['FEM']
|
|
V = MESH['FEM']
|
|
velocity = {}
|
|
Lt = Sqx.shape[3]
|
|
VX = V.sub(0).collapse()
|
|
VY = V.sub(1).collapse()
|
|
VZ = V.sub(2).collapse()
|
|
vv_x = Function(VX)
|
|
vv_y = Function(VY)
|
|
vv_z = Function(VZ)
|
|
|
|
if rank == 0:
|
|
print('{SQtoH5} total number of timesteps: ' + str(Lt))
|
|
|
|
nit = 0
|
|
nit_tot = (Lt-1)/20
|
|
|
|
for t in range(Lt):
|
|
if rank == 0:
|
|
if PRINT_BARS:
|
|
if np.mod(t, nit_tot) < 1:
|
|
sys.stdout.write('\r')
|
|
# Progress bar
|
|
sys.stdout.write(
|
|
"{SQtoH5} [%-40s] %d%%" % ('=='*nit, 100*t/Lt))
|
|
sys.stdout.flush()
|
|
nit = nit + 1
|
|
else:
|
|
print('iteration number:', t)
|
|
|
|
v1_x = Function(P1)
|
|
v1_y = Function(P1)
|
|
v1_z = Function(P1)
|
|
|
|
values_x = np.zeros(v1_x.vector().get_local().shape)
|
|
values_y = np.zeros(v1_y.vector().get_local().shape)
|
|
values_z = np.zeros(v1_z.vector().get_local().shape)
|
|
S0x = Sqx[:, :, :, t]
|
|
S0y = Sqy[:, :, :, t]
|
|
S0z = Sqz[:, :, :, t]
|
|
|
|
for k in range(Xt.shape[0]):
|
|
values_x[k] = S0x[SqXt[k, 0], SqXt[k, 1], SqXt[k, 2]]
|
|
values_y[k] = S0y[SqXt[k, 0], SqXt[k, 1], SqXt[k, 2]]
|
|
values_z[k] = S0z[SqXt[k, 0], SqXt[k, 1], SqXt[k, 2]]
|
|
|
|
v1_x.vector()[:] = values_x
|
|
v1_y.vector()[:] = values_y
|
|
v1_z.vector()[:] = values_z
|
|
|
|
LagrangeInterpolator.interpolate(vv_x, v1_x)
|
|
LagrangeInterpolator.interpolate(vv_y, v1_y)
|
|
LagrangeInterpolator.interpolate(vv_z, v1_z)
|
|
|
|
v = Function(V)
|
|
assign(v.sub(0), vv_x)
|
|
assign(v.sub(1), vv_y)
|
|
assign(v.sub(2), vv_z)
|
|
|
|
velocity[t] = v
|
|
|
|
|
|
return velocity
|
|
|
|
def MASKED(S_REF,S_RAW):
|
|
|
|
[Sq0x,Sq0y,Sq0z] = S_REF
|
|
[Sq1x,Sq1y,Sq1z] = S_RAW
|
|
|
|
if Sq0z.shape != Sq1z.shape:
|
|
print('Original',Sq0z.shape)
|
|
print('To mask',Sq1z.shape)
|
|
|
|
raise Exception('Dimension of masking sequences not match!')
|
|
|
|
# masked Sq1 according to Sq0. Both sequences must have the same sizes.
|
|
Sq2x = np.zeros(Sq0x.shape)
|
|
Sq2y = np.zeros(Sq0y.shape)
|
|
Sq2z = np.zeros(Sq0z.shape)
|
|
|
|
nonzerox = np.where(Sq0x!=0)
|
|
nonzeroy = np.where(Sq0y!=0)
|
|
nonzeroz = np.where(Sq0z!=0)
|
|
|
|
Sq2x[nonzerox] = Sq1x[nonzerox]
|
|
Sq2y[nonzeroy] = Sq1y[nonzeroy]
|
|
Sq2z[nonzeroz] = Sq1z[nonzeroz]
|
|
|
|
return [Sq2x,Sq2y,Sq2z]
|
|
|
|
def READcheckpoint(MESH,mode,ckpath,under_rate,name,zeros=True):
|
|
|
|
V = MESH['FEM']
|
|
W = MESH['FEM'].sub(0).collapse()
|
|
# Checkpoints folders
|
|
unsort_indexes = os.listdir(ckpath)
|
|
indexes = [int(x) for x in unsort_indexes]
|
|
indexes.sort()
|
|
|
|
numt = 0
|
|
nit = 0
|
|
nit_tot = (len(indexes)-1)/20
|
|
|
|
from dolfin import HDF5File
|
|
if mode == 'u':
|
|
velocity = {}
|
|
for k in range(0, len(indexes), under_rate):
|
|
path = ckpath + str(indexes[k]) + '/'+name+'.h5'
|
|
|
|
if zeros:
|
|
if indexes[k] < 10 and indexes[k] > 0:
|
|
path = ckpath + '0' + str(indexes[k]) + '/'+name+'.h5'
|
|
|
|
if rank == 0:
|
|
if PRINT_BARS:
|
|
if np.mod(numt, nit_tot) < 1:
|
|
sys.stdout.write('\r')
|
|
sys.stdout.write(
|
|
"{vel checkpoint} [%-40s] %d%%" % ('=='*nit, 100*numt/(len(indexes)-1)))
|
|
sys.stdout.flush()
|
|
nit = nit + 1
|
|
|
|
#if options['P_measurement']=='p2':
|
|
#v = Function(V)
|
|
#vmix = Function(V)
|
|
#v,p = vmix.split()
|
|
#from dolfin import HDF5File
|
|
#hdf = HDF5File(AORTA['mesh'].mpi_comm(), path, 'r')
|
|
#hdf.read(v, 'u')
|
|
#hdf.close()
|
|
#udat = h5py.File(path, 'r')
|
|
#v.vector()[:] = udat['u']['vector_0']
|
|
#else:
|
|
|
|
v = Function(V)
|
|
comm = MPI.COMM_WORLD
|
|
hdf = HDF5File(MESH['mesh'].mpi_comm(), path, 'r')
|
|
hdf.read(v, 'u/vector_0')
|
|
hdf.close()
|
|
velocity[numt] = v
|
|
numt = numt + 1
|
|
|
|
return velocity
|
|
|
|
if mode == 'p':
|
|
pressure = {}
|
|
for k in range(0, len(indexes), under_rate):
|
|
path = ckpath + str(indexes[k]) + '/'+name+'.h5'
|
|
if zeros:
|
|
if indexes[k] < 10 and indexes[k] > 0:
|
|
path = ckpath + '0' + str(indexes[k]) + '/'+name+'.h5'
|
|
|
|
if rank == 0:
|
|
if PRINT_BARS:
|
|
if np.mod(numt, nit_tot) < 1:
|
|
sys.stdout.write('\r')
|
|
sys.stdout.write(
|
|
"{press checkpoint} [%-40s] %d%%" % ('=='*nit, 100*numt/(len(indexes)-1)))
|
|
sys.stdout.flush()
|
|
nit = nit + 1
|
|
p = Function(W)
|
|
hdf = HDF5File(MESH['mesh'].mpi_comm(), path, 'r')
|
|
hdf.read(p, 'p/vector_0')
|
|
hdf.close()
|
|
pressure[numt] = p
|
|
numt = numt + 1
|
|
|
|
return pressure
|
|
|
|
def CHANGEmesh(MESH,field,mode):
|
|
# To change between meshes the same result
|
|
field2 = {}
|
|
if mode=='velocity':
|
|
V = MESH['FEM']
|
|
elif mode=='pressure':
|
|
V = MESH['FEM'].sub(0).collapse()
|
|
else:
|
|
raise Exception('You need to assign a mode for the changing!')
|
|
for k in range(len(list(field))):
|
|
if rank==0:
|
|
print('CHANGING MESH: index',k)
|
|
v2 = Function(V)
|
|
LagrangeInterpolator.interpolate(v2,field[k])
|
|
field2[k] = v2
|
|
|
|
|
|
return field2
|
|
|
|
def GradP(MESH,pressure,options):
|
|
|
|
W = MESH['FEM'].sub(0).collapse()
|
|
boundaries = MESH['boundaries']
|
|
barye2mmHg = 1/1333.22387415
|
|
|
|
Tf = Constant(options['Tf'])
|
|
dt = Constant(options['Tf']/len(list(pressure)))
|
|
NT = int(float(Tf)/float(dt))
|
|
|
|
ind = 0
|
|
t = float(dt)
|
|
ds = Measure("ds", subdomain_data=boundaries)
|
|
ones = interpolate(Constant(1),W)
|
|
DELTA_p = np.zeros([NT])
|
|
|
|
p0 = Function(W)
|
|
nit = 0
|
|
nit_tot = (NT-1)/20
|
|
|
|
zero_point = None
|
|
if options['reference']['fix_const']:
|
|
zero_point = options['reference']['zero_point']
|
|
zero_point = np.array(zero_point)
|
|
ndim = W.mesh().topology().dim()
|
|
zero_point = W.tabulate_dof_coordinates().reshape((-1, ndim))[0]
|
|
|
|
#################### THE TIME LOOP #########################
|
|
while ind < NT:
|
|
p0.assign(pressure[ind])
|
|
|
|
if rank == 0:
|
|
if PRINT_BARS:
|
|
if np.mod(ind, nit_tot) < 1:
|
|
sys.stdout.write('\r')
|
|
# Progress bar
|
|
sys.stdout.write("{ref} [%-40s] %f seg" %
|
|
('=='*nit, round(t, 2)))
|
|
sys.stdout.flush()
|
|
nit = nit + 1
|
|
|
|
if options['reference']['fix_const']:
|
|
p0.vector()[:] -= p0(zero_point)
|
|
|
|
i_pin = assemble(p0*ds(2))
|
|
i_pout = assemble(p0*ds(3))
|
|
i_s2 = assemble(ones*ds(2))
|
|
i_s3 = assemble(ones*ds(3))
|
|
|
|
DELTA_p[ind] = barye2mmHg*(i_pin/i_s2-i_pout/i_s3)
|
|
ind = ind + 1
|
|
t = t + float(dt)
|
|
|
|
if options['reference']['save']:
|
|
|
|
if rank==0:
|
|
np.savetxt(options['reference']['savepath'],DELTA_p)
|
|
print('\n Pressure drop save in ' + options['reference']['savepath'])
|
|
|
|
def Fluxes(MESH,velocity,options,bnd):
|
|
|
|
mesh = MESH['mesh']
|
|
W = MESH['FEM']
|
|
boundaries = MESH['boundaries']
|
|
n = FacetNormal(mesh)
|
|
Tf = Constant(options['Tf'])
|
|
rho = Constant(options['density'])
|
|
dt = Constant(options['Tf']/len(list(velocity)))
|
|
NT = int(float(Tf)/float(dt))
|
|
ind = 0
|
|
t = float(dt)
|
|
ds = Measure("ds", subdomain_data=boundaries)
|
|
QQ = {}
|
|
for k in bnd:
|
|
QQ[k] = np.zeros([NT])
|
|
|
|
u = Function(W)
|
|
nit = 0
|
|
nit_tot = (NT-1)/20
|
|
#################### THE TIME LOOP #########################
|
|
while ind<NT:
|
|
u.assign(velocity[ind])
|
|
if rank==0:
|
|
if PRINT_BARS:
|
|
if np.mod(ind,nit_tot)<1:
|
|
sys.stdout.write('\r')
|
|
# Progress bar
|
|
sys.stdout.write("{ref} [%-40s] %f seg" % ('=='*nit, round(t,2) ))
|
|
sys.stdout.flush()
|
|
nit = nit +1
|
|
|
|
for k in bnd:
|
|
QQ[k][ind] = assemble(dot(u,n)*ds(k))
|
|
|
|
ind = ind + 1
|
|
t = t + float(dt)
|
|
|
|
|
|
return QQ
|
|
|
|
def READ_CIB(user,path_to_cib):
|
|
|
|
ST0 = sio.loadmat(path_to_cib+'Segmented_Aorta.mat')
|
|
ST1 = sio.loadmat(path_to_cib+'VEL_AP.mat')
|
|
ST2 = sio.loadmat(path_to_cib+'VEL_FH.mat')
|
|
ST3 = sio.loadmat(path_to_cib+'VEL_RL.mat')
|
|
ST4 = sio.loadmat(path_to_cib+'voxel_size_n.mat')
|
|
|
|
mask = ST0['Segmented_Aorta']
|
|
v1 = ST1['VEL_AP']
|
|
v2 = ST2['VEL_FH']
|
|
v3 = ST3['VEL_RL']
|
|
sizes = ST4['voxel_size']
|
|
[dx,dy,dz] = sizes[0]
|
|
|
|
del ST0,ST1,ST2,ST3
|
|
|
|
v = np.zeros(v1.shape)
|
|
for k in range(v1.shape[3]):
|
|
v[:,:,:,k] = np.sqrt(v1[:,:,:,k]**2 + v2[:,:,:,k]**2 + v3[:,:,:,k]**2)*mask[:,:,:]
|
|
v1[:,:,:,k] = v1[:,:,:,k]*mask[:,:,:]
|
|
v2[:,:,:,k] = v2[:,:,:,k]*mask[:,:,:]
|
|
v3[:,:,:,k] = v3[:,:,:,k]*mask[:,:,:]
|
|
|
|
print('The maximum velocity is ',np.max(v))
|
|
|
|
|
|
# To order and clean up
|
|
[row,col,dep,numt] = v1.shape
|
|
|
|
#np.savez_compressed(path_to_cib + 'seq.npz', x=v1, y=v2,z=v3 , size = [dx,dy,dz])
|
|
np.savez_compressed(path_to_cib + 'seq0.npz', x=v1[:,20:80,:,:], y=v2[:,20:80,:,:],z=v3[:,20:80,:,:] , size = [dx,dy,dz])
|
|
|
|
def CIBtoH5(pathtocib,times,dt,outpath,interpolate=False,mesh_path=None,flip=False):
|
|
|
|
import csv
|
|
|
|
mesh = Mesh(pathtocib+'aorta.xml')
|
|
comm = mesh.mpi_comm()
|
|
Evel = VectorElement('Lagrange',mesh.ufl_cell(),1)
|
|
V = FunctionSpace(mesh,Evel)
|
|
|
|
VX = V.sub(0).collapse()
|
|
VY = V.sub(1).collapse()
|
|
VZ = V.sub(2).collapse()
|
|
vv_x = Function(VX)
|
|
vv_y = Function(VY)
|
|
vv_z = Function(VZ)
|
|
v = Function(V)
|
|
v.rename('velocity','meas')
|
|
|
|
mx = dof_to_vertex_map(VX)
|
|
my = dof_to_vertex_map(VY)
|
|
mz = dof_to_vertex_map(VZ)
|
|
|
|
if interpolate:
|
|
mesh2 = Mesh(mesh_path+'aorta.xml')
|
|
Evel2 = VectorElement('Lagrange',mesh2.ufl_cell(),1)
|
|
V2 = FunctionSpace(mesh2,Evel2)
|
|
v2 = Function(V2)
|
|
v2.rename('velocity','meas')
|
|
|
|
|
|
xdmf_v = XDMFFile(outpath+'meas.xdmf')
|
|
|
|
for ti in times:
|
|
print('reading aorta '+str(ti) + '.csv')
|
|
with open(pathtocib + 'aorta'+str(ti)+'.csv', 'r') as csvfile:
|
|
mylist = [row[0] for row in csv.reader(csvfile, delimiter=';')]
|
|
|
|
Values = np.array(mylist)
|
|
file_x = np.zeros([len(Values)])
|
|
file_y = np.zeros([len(Values)])
|
|
file_z = np.zeros([len(Values)])
|
|
|
|
for l in range(len(Values)):
|
|
row = Values[l].split(',')
|
|
file_x[l] = float(row[0])
|
|
file_y[l] = float(row[1])
|
|
file_z[l] = float(row[2])
|
|
|
|
vv_x.vector()[:] = file_x[mx]
|
|
vv_y.vector()[:] = file_y[my]
|
|
vv_z.vector()[:] = file_z[mz]
|
|
|
|
if flip:
|
|
# 102 for 13mm, Normal
|
|
assign(v.sub(1),vv_x)
|
|
assign(v.sub(0),vv_y)
|
|
assign(v.sub(2),vv_z)
|
|
else:
|
|
assign(v.sub(0),vv_x)
|
|
assign(v.sub(1),vv_y)
|
|
assign(v.sub(2),vv_z)
|
|
|
|
ck_path = outpath+'checkpoint/{i}/'.format(i=ti)
|
|
time_real = (ti-1)*dt
|
|
|
|
if interpolate:
|
|
LagrangeInterpolator.interpolate(v2,v)
|
|
inout.write_HDF5_data(comm, ck_path + '/u.h5', v2, '/u', t=time_real)
|
|
xdmf_v.write(v2, time_real)
|
|
else:
|
|
inout.write_HDF5_data(comm, ck_path + '/u.h5', v, '/u', t=time_real)
|
|
xdmf_v.write(v, time_real )
|
|
|
|
|
|
xdmf_v.close()
|
|
|
|
def CLOCK(rank,t1,t2):
|
|
if rank==0:
|
|
tot_time = np.round(t2 - t1,2)
|
|
if tot_time<60:
|
|
print('Total time: ' + str(tot_time) + ' s')
|
|
elif tot_time<3600:
|
|
time_min = tot_time//60
|
|
time_sec = tot_time - time_min*60
|
|
print('Total time: ' + str(time_min) + ' min, ' + str(time_sec) + ' s')
|
|
else:
|
|
time_hour = tot_time//3600
|
|
time_tot2 = tot_time - time_hour*3600
|
|
time_min = time_tot2//60
|
|
time_sec = time_tot2 - time_min*60
|
|
print('Total time: ' + str(time_hour) + ' hrs, ' + str(time_min) + ' min, ' + str(time_sec) + ' s')
|
|
|
|
|
|
|
|
def SCANNER(options):
|
|
|
|
########################################
|
|
#
|
|
# Basic Tools
|
|
#
|
|
########################################
|
|
if 'kspace_cib' in options:
|
|
if options['kspace_cib']['apply']:
|
|
print('--- kspace from CIB data ---')
|
|
import Graphics
|
|
import scipy.io as sio
|
|
S = sio.loadmat(options['kspace_cib']['kspace'])
|
|
VENC = options['kspace_cib']['VENC']
|
|
KS0 = S['k_space']['Flow0'][0][0]
|
|
KS1 = S['k_space']['Flow1'][0][0]
|
|
KS2 = S['k_space']['Flow2'][0][0]
|
|
KS3 = S['k_space']['Flow3'][0][0]
|
|
|
|
[Nx, Ny, Nz, Nt] = KS0.shape
|
|
|
|
M0 = np.zeros([Nx, Ny, Nz, Nt], dtype=complex)
|
|
M1 = np.zeros([Nx, Ny, Nz, Nt], dtype=complex)
|
|
M2 = np.zeros([Nx, Ny, Nz, Nt], dtype=complex)
|
|
M3 = np.zeros([Nx, Ny, Nz, Nt], dtype=complex)
|
|
|
|
for t in range(Nt):
|
|
KS0[:, :, :, t] = np.fft.fftshift(KS0[:, :, :, t])
|
|
KS1[:, :, :, t] = np.fft.fftshift(KS1[:, :, :, t])
|
|
KS2[:, :, :, t] = np.fft.fftshift(KS2[:, :, :, t])
|
|
KS3[:, :, :, t] = np.fft.fftshift(KS3[:, :, :, t])
|
|
M0[:, :, :, t] = np.fft.ifftn(KS0[:, :, :, t])
|
|
M1[:, :, :, t] = np.fft.ifftn(KS1[:, :, :, t])
|
|
M2[:, :, :, t] = np.fft.ifftn(KS2[:, :, :, t])
|
|
M3[:, :, :, t] = np.fft.ifftn(KS3[:, :, :, t])
|
|
|
|
M0 = np.fft.fftshift(M0, axes=0)
|
|
M1 = np.fft.fftshift(M1, axes=0)
|
|
M2 = np.fft.fftshift(M2, axes=0)
|
|
M3 = np.fft.fftshift(M3, axes=0)
|
|
|
|
FFE_PCA = {}
|
|
FFE_PCA['MR_FFE_FH'] = np.abs(M1)
|
|
FFE_PCA['MR_FFE_AP'] = np.abs(M2)
|
|
FFE_PCA['MR_FFE_RL'] = np.abs(M3)
|
|
|
|
phi0 = np.angle(M0)
|
|
vel_FH = VENC*(np.angle(M1)-np.angle(M0))/np.pi
|
|
vel_AP = VENC*(np.angle(M2)-np.angle(M0))/np.pi
|
|
vel_RL = VENC*(np.angle(M3)-np.angle(M0))/np.pi
|
|
|
|
[vel_FH, vel_AP, vel_RL] = dealiased(VENC, vel_FH, vel_AP, vel_RL)
|
|
|
|
FFE_PCA['MR_PCA_FH'] = vel_FH
|
|
FFE_PCA['MR_PCA_AP'] = vel_AP
|
|
FFE_PCA['MR_PCA_RL'] = vel_RL
|
|
|
|
#sio.savemat(options['kspace_cib']['output'],FFE_PCA)
|
|
Graphics.TakePhoto(FFE_PCA['MR_PCA_FH'][:, :, 22, 8])
|
|
Graphics.TakePhoto(FFE_PCA['MR_PCA_AP'][:, :, 22, 8])
|
|
Graphics.TakePhoto(FFE_PCA['MR_PCA_RL'][:, :, 22, 8])
|
|
|
|
print(FFE_PCA['MR_PCA_AP'][78, 20, 22, 8])
|
|
|
|
Graphics.PrintSequence(FFE_PCA['MR_PCA_FH'][:,:,:,8])
|
|
#Graphics.TakePhoto(FFE_PCA['MR_FFE_FH'][:,:,22,8])
|
|
|
|
if 'create_leo' in options:
|
|
if options['create_leo']['apply']:
|
|
print('Preparing files for create LEO mesh...')
|
|
# Creating the BOX
|
|
[vx, vy, vz] = LOADsequences(options['create_leo']['velseq'])
|
|
[Lx, Ly, Lz, Lt] = vx.shape
|
|
|
|
if options['create_leo']['masked']:
|
|
import Graphics
|
|
mask = Graphics.MATmedical(
|
|
options['create_leo']['segmentation'])
|
|
mask = mask.transpose((2, 0, 1))
|
|
raise Exception('a_crop and b_crop need to be checked!')
|
|
a_crop = 48
|
|
b_crop = 178
|
|
mask = mask[:, :, a_crop:b_crop]
|
|
for t in range(Lt):
|
|
vx[:, :, :, t] = vx[:, :, :, t]*mask
|
|
vy[:, :, :, t] = vy[:, :, :, t]*mask
|
|
vz[:, :, :, t] = vz[:, :, :, t]*mask
|
|
|
|
if '.mat' in options['create_leo']['velseq']:
|
|
import scipy.io as sio
|
|
S = sio.loadmat(options['create_leo']['velseq'])
|
|
[ddx,ddy,ddz] = S['u_R1']['VoxelSize'][0][0][0]
|
|
[ddx,ddy,ddz] = [ddx/10,ddy/10,ddz/10] #mm to cm
|
|
BOX = LOADmesh('box_zeros', resol=[ddx,ddy,ddz], boxsize=[Lx,Ly,Lz])
|
|
else:
|
|
[ddx,ddy,ddz] = options['create_leo']['resol']
|
|
ranges = options['create_leo']['ranges']
|
|
BOX = LOADmesh('fix_resolution', resol=[ddx,ddy,ddz], boxsize=[Lx,Ly,Lz], ranges=ranges)
|
|
|
|
Xt = BOX['nodes']
|
|
SqXt = BOX['seq']
|
|
|
|
# Saving the nodes...
|
|
print('Saving the nodes')
|
|
np.savetxt(options['create_leo']['outpath']+'nodes.txt', Xt)
|
|
|
|
# Interpolating the velocity
|
|
v1_x = Function(BOX['FEM'])
|
|
v1_y = Function(BOX['FEM'])
|
|
v1_z = Function(BOX['FEM'])
|
|
values_x = np.zeros(v1_x.vector().get_local().shape)
|
|
values_y = np.zeros(v1_y.vector().get_local().shape)
|
|
values_z = np.zeros(v1_z.vector().get_local().shape)
|
|
v = np.sqrt(vx**2 + vy**2 + vz**2)
|
|
tmax = np.where(v == np.max(v))[-1]
|
|
S0x = vx[:, :, :, tmax]
|
|
S0y = vy[:, :, :, tmax]
|
|
S0z = vz[:, :, :, tmax]
|
|
for k in range(Xt.shape[0]):
|
|
values_x[k] = S0x[SqXt[k, 0], SqXt[k, 1], SqXt[k, 2]]
|
|
values_y[k] = S0y[SqXt[k, 0], SqXt[k, 1], SqXt[k, 2]]
|
|
values_z[k] = S0z[SqXt[k, 0], SqXt[k, 1], SqXt[k, 2]]
|
|
|
|
print('Saving the interpolated velocities values...')
|
|
# Saving velocity values
|
|
np.savetxt(options['create_leo']['outpath'] +
|
|
'resol.txt', [ddx, ddy, ddz])
|
|
np.savetxt(options['create_leo']['outpath']+'ux.txt', values_x)
|
|
np.savetxt(options['create_leo']['outpath']+'uy.txt', values_y)
|
|
np.savetxt(options['create_leo']['outpath']+'uz.txt', values_z)
|
|
print('In order to continue with the procedure, MATLAB should be opened')
|
|
print('Then run the script CREATE_MESH.m in the LEO_matlab directory.')
|
|
|
|
if 'phase_contrast' in options:
|
|
if options['phase_contrast']['apply']:
|
|
if rank==0:
|
|
print('Computing Phase Contrast')
|
|
|
|
VENC = options['phase_contrast']['VENC']
|
|
M0S = np.load(options['phase_contrast']['meas0'])
|
|
MGS = np.load(options['phase_contrast']['measG'])
|
|
M0x = M0S['x']
|
|
MGx = MGS['x']
|
|
vx = VENC*(np.angle(MGx) - np.angle(M0x))/np.pi
|
|
del M0x,MGx
|
|
M0y = M0S['y']
|
|
MGy = MGS['y']
|
|
vy = VENC*(np.angle(MGy) - np.angle(M0y))/np.pi
|
|
del M0y,MGy
|
|
M0z = M0S['z']
|
|
MGz = MGS['z']
|
|
vz = VENC*(np.angle(MGz) - np.angle(M0z))/np.pi
|
|
del M0z,MGz
|
|
|
|
INDx = np.any(np.abs(vx)>VENC)
|
|
INDy = np.any(np.abs(vy)>VENC)
|
|
INDz = np.any(np.abs(vz)>VENC)
|
|
|
|
if INDx or INDy or INDz:
|
|
if not options['phase_contrast']['dealiased']:
|
|
if rank==0:
|
|
raise Exception('There is some velocity greater or lower than VENC. You should turn dealised true!')
|
|
|
|
|
|
if options['phase_contrast']['dealiased']:
|
|
[vx_d,vy_d,vz_d] = dealiased(VENC,vx,vy,vz)
|
|
else:
|
|
vx_d = vx; vy_d = vy; vz_d = vz
|
|
|
|
if rank==0:
|
|
print('saving velocity as ' + options['phase_contrast']['output'])
|
|
np.savez_compressed(options['phase_contrast']['output'], x=vx_d, y=vy_d, z=vz_d)
|
|
|
|
if 'resize' in options:
|
|
if options['resize']['apply']:
|
|
if rank == 0:
|
|
print('Changing the size of the sequence')
|
|
|
|
[Nx, Ny, Nz] = options['resize']['dim']
|
|
[Sqx, Sqy, Sqz] = LOADsequences(options['resize']['loadseq'])
|
|
[Nx0, Ny0, Nz0, Nt] = Sqx.shape
|
|
|
|
Sqx_new = np.zeros([Nx, Ny, Nz, Nt])
|
|
Sqy_new = np.zeros([Nx, Ny, Nz, Nt])
|
|
Sqz_new = np.zeros([Nx, Ny, Nz, Nt])
|
|
|
|
ax = int((Nx - Nx0)/2)
|
|
bx = ax + Nx0
|
|
ay = int((Ny - Ny0)/2)
|
|
by = ay + Ny0
|
|
az = int((Nz - Nz0)/2)
|
|
bz = az + Nz0
|
|
|
|
for t in range(Nt):
|
|
Sqx_new[ax:bx, ay:by, az:bz, t] = Sqx[:, :, :, t]
|
|
Sqy_new[ax:bx, ay:by, az:bz, t] = Sqy[:, :, :, t]
|
|
Sqz_new[ax:bx, ay:by, az:bz, t] = Sqz[:, :, :, t]
|
|
|
|
if options['resize']['save']:
|
|
if rank == 0:
|
|
print('saving the sequences' +
|
|
options['resize']['savepath'])
|
|
np.savez_compressed(
|
|
options['resize']['savepath'], x=Sqx_new, y=Sqy_new, z=Sqz_new)
|
|
|
|
if 'magnetization_CIB' in options:
|
|
if options['magnetization_CIB']['apply']:
|
|
if rank == 0:
|
|
print('-- Generating k-space for CIB phantom data-- ')
|
|
print('Using 0-G gradients')
|
|
|
|
import Graphics
|
|
masterpath = options['magnetization_CIB']['loadpath']
|
|
mask = Graphics.MATmedical(masterpath+'Segmented_Aorta.mat')
|
|
mask = mask.transpose((2, 0, 1))
|
|
VENC = Graphics.MATmedical(masterpath+'VENC.mat')
|
|
VENC = VENC[0][0]
|
|
if rank == 0:
|
|
print('Recognized VENC = ' + str(VENC) + ' cm/s')
|
|
|
|
# Magnitudes
|
|
magnitude_AP = Graphics.MATmedical(masterpath+'MR_FFE_AP.mat')
|
|
magnitude_FH = Graphics.MATmedical(masterpath+'MR_FFE_FH.mat')
|
|
magnitude_RL = Graphics.MATmedical(masterpath+'MR_FFE_RL.mat')
|
|
# Velocities
|
|
vel_AP = Graphics.MATmedical(masterpath+'MR_PCA_AP.mat')
|
|
vel_FH = Graphics.MATmedical(masterpath+'MR_PCA_FH.mat')
|
|
vel_RL = Graphics.MATmedical(masterpath+'MR_PCA_RL.mat')
|
|
|
|
# Cropping images
|
|
a_crop = 48
|
|
b_crop = 178
|
|
mask = mask[:, :, a_crop:b_crop]
|
|
magnitude_AP = magnitude_AP[:, :, a_crop:b_crop, :]
|
|
vel_AP = vel_AP[:, :, a_crop:b_crop, :]
|
|
magnitude_RL = magnitude_RL[:, :, a_crop:b_crop, :]
|
|
vel_RL = vel_RL[:, :, a_crop:b_crop, :]
|
|
magnitude_FH = magnitude_FH[:, :, a_crop:b_crop, :]
|
|
vel_FH = vel_FH[:, :, a_crop:b_crop, :]
|
|
|
|
# Reference phase
|
|
[Nx, Ny, Nz, Nt] = vel_AP.shape
|
|
gamma = 267.513e6 # rad/Tesla/sec Gyromagnetic ratio for H nuclei
|
|
B0 = 1.5 # Tesla Magnetic Field Strenght
|
|
TE = 5e-3 # Echo-time
|
|
[X, Y, Z] = np.meshgrid(np.linspace(
|
|
0, Ny, Ny), np.linspace(0, Nx, Nx), np.linspace(0, Nz, Nz))
|
|
frecX = 0.21
|
|
p_ref = np.zeros([Nx, Ny, Nz])
|
|
p_flow_AP = np.zeros([Nx, Ny, Nz])
|
|
p_flow_RL = np.zeros([Nx, Ny, Nz])
|
|
p_flow_FH = np.zeros([Nx, Ny, Nz])
|
|
#Sx2G = np.zeros([Nx,Ny,Nz,Nt],dtype=complex)
|
|
#Sy2G = np.zeros([Nx,Ny,Nz,Nt],dtype=complex)
|
|
#Sz2G = np.zeros([Nx,Ny,Nz,Nt],dtype=complex)
|
|
Sx20 = np.zeros([Nx, Ny, Nz, Nt], dtype=complex)
|
|
Sy20 = np.zeros([Nx, Ny, Nz, Nt], dtype=complex)
|
|
Sz20 = np.zeros([Nx, Ny, Nz, Nt], dtype=complex)
|
|
|
|
treshold = np.mean(magnitude_AP)
|
|
#p_ref = (gamma*B0*TE+frecX*X)*(np.abs(magnitude_AP[:,:,:,7])>treshold)
|
|
p_ref = (gamma*B0*TE+frecX*X)*(mask > 0.5)
|
|
|
|
for t in range(Nt):
|
|
p_flow_AP = p_ref + np.pi*vel_AP[:, :, :, t]/VENC
|
|
p_flow_RL = p_ref + np.pi*vel_RL[:, :, :, t]/VENC
|
|
p_flow_FH = p_ref + np.pi*vel_FH[:, :, :, t]/VENC
|
|
|
|
Sx20[:, :, :, t] = magnitude_AP[:, :, :, t] * \
|
|
np.cos(p_ref) + 1j*magnitude_AP[:, :, :, t]*np.sin(p_ref)
|
|
Sy20[:, :, :, t] = magnitude_RL[:, :, :, t] * \
|
|
np.cos(p_ref) + 1j*magnitude_RL[:, :, :, t]*np.sin(p_ref)
|
|
Sz20[:, :, :, t] = magnitude_FH[:, :, :, t] * \
|
|
np.cos(p_ref) + 1j*magnitude_FH[:, :, :, t]*np.sin(p_ref)
|
|
|
|
#Sx2G[:,:,:,t] = magnitude_AP[:,:,:,t]*np.cos(p_flow_AP) + 1j*magnitude_AP[:,:,:,t]*np.sin(p_flow_AP)
|
|
#Sy2G[:,:,:,t] = magnitude_RL[:,:,:,t]*np.cos(p_flow_RL) + 1j*magnitude_RL[:,:,:,t]*np.sin(p_flow_RL)
|
|
#Sz2G[:,:,:,t] = magnitude_FH[:,:,:,t]*np.cos(p_flow_FH) + 1j*magnitude_FH[:,:,:,t]*np.sin(p_flow_FH)
|
|
|
|
if rank == 0:
|
|
print('saving the sequences' +
|
|
options['magnetization_CIB']['savepath'])
|
|
savepath = options['magnetization_CIB']['savepath']
|
|
basepath = savepath.replace('.npz', '')
|
|
savepath0 = basepath + '0.npz'
|
|
savepathG = basepath + 'G.npz'
|
|
np.savez_compressed(savepath0, x=Sx20, y=Sy20, z=Sz20)
|
|
#np.savez_compressed(savepathG, x=Sx2G,y=Sy2G,z=Sz2G)
|
|
|
|
if 'magnetization' in options:
|
|
if options['magnetization']['apply']:
|
|
if rank==0:
|
|
print('-- Generating k-space -- ')
|
|
print('Using 0-G gradients')
|
|
|
|
[Sqx,Sqy,Sqz] = LOADsequences(options['magnetization']['loadseq'])
|
|
[Nx,Ny,Nz,Nt] = Sqx.shape
|
|
VENC = options['magnetization']['VENC']
|
|
gamma = 267.513e6 # rad/Tesla/sec Gyromagnetic ratio for H nuclei
|
|
B0 = 1.5 # Tesla Magnetic Field Strenght
|
|
TE = 5e-3 # Echo-time
|
|
[X,Y,Z] = np.meshgrid(np.linspace(0,Ny,Ny),np.linspace(0,Nx,Nx),np.linspace(0,Nz,Nz))
|
|
frecX = 0.21 # 0.21 without leaps: 0.01
|
|
|
|
Mx = np.zeros(Sqx.shape)
|
|
My = np.zeros(Sqx.shape)
|
|
Mz = np.zeros(Sqx.shape)
|
|
Sx2G = np.zeros(Sqx.shape,dtype=complex)
|
|
Sy2G = np.zeros(Sqx.shape,dtype=complex)
|
|
Sz2G = np.zeros(Sqx.shape,dtype=complex)
|
|
Sx20 = np.zeros(Sqx.shape,dtype=complex)
|
|
Sy20 = np.zeros(Sqx.shape,dtype=complex)
|
|
Sz20 = np.zeros(Sqx.shape,dtype=complex)
|
|
Fx0 = np.zeros(Sqx.shape)
|
|
Fy0 = np.zeros(Sqx.shape)
|
|
Fz0 = np.zeros(Sqx.shape)
|
|
FxG = np.zeros(Sqx.shape)
|
|
FyG = np.zeros(Sqx.shape)
|
|
FzG = np.zeros(Sqx.shape)
|
|
|
|
# Creating M
|
|
for t in range(Nt):
|
|
Mx[:,:,:,t] = 0.05 + 0.9*(np.abs(Sqx[:,:,:,t])>0.001) + 0.05*np.sqrt(np.abs(Sqx[:,:,:,t]))
|
|
My[:,:,:,t] = 0.05 + 0.9*(np.abs(Sqy[:,:,:,t])>0.001) + 0.05*np.sqrt(np.abs(Sqy[:,:,:,t]))
|
|
Mz[:,:,:,t] = 0.05 + 0.9*(np.abs(Sqz[:,:,:,t])>0.001) + 0.05*np.sqrt(np.abs(Sqz[:,:,:,t]))
|
|
# Organs
|
|
Mx[:,:,:,t] = Mx[:,:,:,t] + 1*( (X-Ny/2)**2 + 0.75*(Y-Nx/2)**2 < 35**2 ) * ( (X-Ny/2)**2 + 0.75*(Y-Nx/2)**2 > 33**2 )
|
|
Mx[:,:,:,t] = Mx[:,:,:,t] + 0.7*( 0.3*(X-Ny/2+9)**2 + 0.7*(Y-Nx/2-4)**2 < 3**2 )* ( np.abs(Z-Nz/2)<25 )
|
|
My[:,:,:,t] = My[:,:,:,t] + 1*( (X-Ny/2)**2 + 0.75*(Y-Nx/2)**2 < 35**2 ) * ( (X-Ny/2)**2 + 0.75*(Y-Nx/2)**2 > 33**2 )
|
|
My[:,:,:,t] = My[:,:,:,t] + 0.7*( 0.3*(X-Ny/2+9)**2 + 0.7*(Y-Nx/2-4)**2 < 3**2 )* ( np.abs(Z-Nz/2)<25 )
|
|
Mz[:,:,:,t] = Mz[:,:,:,t] + 1*( (X-Ny/2)**2 + 0.75*(Y-Nx/2)**2 < 35**2 ) * ( (X-Ny/2)**2 + 0.75*(Y-Nx/2)**2 > 33**2 )
|
|
Mz[:,:,:,t] = Mz[:,:,:,t] + 0.7*( 0.3*(X-Ny/2+9)**2 + 0.7*(Y-Nx/2-4)**2 < 3**2 )* ( np.abs(Z-Nz/2)<25 )
|
|
FxG[:,:,:,t] = (gamma*B0*TE+frecX*X)*(np.abs(Mz[:,:,:,t])>0.4) + np.pi*Sqx[:,:,:,t]/VENC
|
|
FyG[:,:,:,t] = (gamma*B0*TE+frecX*X)*(np.abs(Mz[:,:,:,t])>0.4) + np.pi*Sqy[:,:,:,t]/VENC
|
|
FzG[:,:,:,t] = (gamma*B0*TE+frecX*X)*(np.abs(Mz[:,:,:,t])>0.4) + np.pi*Sqz[:,:,:,t]/VENC
|
|
Fx0[:,:,:,t] = (gamma*B0*TE+frecX*X)*(np.abs(Mz[:,:,:,t])>0.4)
|
|
Fy0[:,:,:,t] = (gamma*B0*TE+frecX*X)*(np.abs(Mz[:,:,:,t])>0.4)
|
|
Fz0[:,:,:,t] = (gamma*B0*TE+frecX*X)*(np.abs(Mz[:,:,:,t])>0.4)
|
|
|
|
Sx20[:,:,:,t] = Mx[:,:,:,t]*np.cos(Fx0[:,:,:,t]) + 1j*Mx[:,:,:,t]*np.sin(Fx0[:,:,:,t])
|
|
Sy20[:,:,:,t] = My[:,:,:,t]*np.cos(Fy0[:,:,:,t]) + 1j*My[:,:,:,t]*np.sin(Fy0[:,:,:,t])
|
|
Sz20[:,:,:,t] = Mz[:,:,:,t]*np.cos(Fz0[:,:,:,t]) + 1j*Mz[:,:,:,t]*np.sin(Fz0[:,:,:,t])
|
|
Sx2G[:,:,:,t] = Mx[:,:,:,t]*np.cos(FxG[:,:,:,t]) + 1j*Mx[:,:,:,t]*np.sin(FxG[:,:,:,t])
|
|
Sy2G[:,:,:,t] = My[:,:,:,t]*np.cos(FyG[:,:,:,t]) + 1j*My[:,:,:,t]*np.sin(FyG[:,:,:,t])
|
|
Sz2G[:,:,:,t] = Mz[:,:,:,t]*np.cos(FzG[:,:,:,t]) + 1j*Mz[:,:,:,t]*np.sin(FzG[:,:,:,t])
|
|
|
|
|
|
# Adding noise by cartesian components
|
|
for t in range(Nt):
|
|
noise1 = np.random.normal(-0.0007,0.049,[Nx,Ny,Nz])
|
|
noise2 = np.random.normal(-0.0007,0.049,[Nx,Ny,Nz])
|
|
Sx20[:,:,:,t] += noise1 + 1j*noise2
|
|
Sy20[:,:,:,t] += noise1 + 1j*noise2
|
|
Sz20[:,:,:,t] += noise1 + 1j*noise2
|
|
|
|
noise1g = np.random.normal(-0.0007,0.049,[Nx,Ny,Nz])
|
|
noise2g = np.random.normal(-0.0007,0.049,[Nx,Ny,Nz])
|
|
Sx2G[:,:,:,t] += noise1g + 1j*noise2g
|
|
Sy2G[:,:,:,t] += noise1g + 1j*noise2g
|
|
Sz2G[:,:,:,t] += noise1g + 1j*noise2g
|
|
|
|
|
|
if rank==0:
|
|
print('saving the sequences' + options['magnetization']['savepath'])
|
|
savepath = options['magnetization']['savepath']
|
|
basepath = savepath.replace('.npz','')
|
|
savepath0 = basepath + '0.npz'
|
|
savepathG = basepath + 'G.npz'
|
|
np.savez_compressed(savepath0, x=Sx20,y=Sy20,z=Sz20)
|
|
np.savez_compressed(savepathG, x=Sx2G,y=Sy2G,z=Sz2G)
|
|
|
|
if 'sequence' in options:
|
|
if options['sequence']['apply']:
|
|
if rank == 0:
|
|
print('--- Creating Sequences ---')
|
|
|
|
resol = options['sequence']['resol']
|
|
boxsize = options['sequence']['boxsize']
|
|
ranges = options['sequence']['ranges']
|
|
BOX = LOADmesh(options['sequence']['boxtype'], resol=resol, boxsize=boxsize, ranges=ranges)
|
|
meshpath = options['sequence']['meshpath']
|
|
AORTA = LOADmesh('aorta', meshpath=meshpath)
|
|
[Sqx, Sqy, Sqz] = CREATEsequences(BOX, AORTA, options)
|
|
|
|
if 'norms' in options:
|
|
if options['norms']['apply']:
|
|
[BOX,AORTA,LEO] = CREATEmeshes(options)
|
|
if rank==0:
|
|
print('Computing Norms between fields')
|
|
ct_norm = checkpoint_norm(AORTA,options['norms']['field_path'],'ct')
|
|
kal_norm = checkpoint_norm(AORTA,options['norms']['field_path'],'roukf')
|
|
|
|
if options['norms']['plot']:
|
|
import matplotlib.pyplot as plt
|
|
plt.plot(ct_norm, color='blue',marker='.',label='ct')
|
|
#plt.plot(meas_norm,label='meas')
|
|
plt.plot(kal_norm,color='green',label='kalman')
|
|
#plt.ylim([0,8500])
|
|
plt.legend(fontsize=14)
|
|
plt.show()
|
|
|
|
if 'CIBtoH5' in options:
|
|
print('---Converting CIB files into H5---')
|
|
if options['CIBtoH5']['apply']:
|
|
path_to_cib = options['CIBtoH5']['data_path']
|
|
outpath = options['CIBtoH5']['outpath']
|
|
times = options['CIBtoH5']['times']
|
|
dt = options['CIBtoH5']['dt']
|
|
flip = options['CIBtoH5']['flip']
|
|
|
|
if options['CIBtoH5']['interpolate']:
|
|
mesh_path = options['CIBtoH5']['mesh_path']
|
|
CIBtoH5(path_to_cib,times,dt,outpath,interpolate=True,mesh_path=mesh_path,flip=flip)
|
|
else:
|
|
CIBtoH5(path_to_cib,times,dt,outpath,flip=flip)
|
|
|
|
########################################
|
|
#
|
|
# Undersampling
|
|
#
|
|
########################################
|
|
if 'cs' in options:
|
|
if options['cs']['apply']:
|
|
if rank==0:
|
|
print('Applying Compressed Sensing')
|
|
|
|
|
|
[Sqx, Sqy, Sqz] = LOADsequences(options['cs']['seqpath'])
|
|
|
|
|
|
import CS
|
|
if options['cs']['short']:
|
|
[Mx,My,Mz] = LOADsequences(options['cs']['Mpath'])
|
|
CS.undersampling_short(Mx,My,Mz,options)
|
|
else:
|
|
CS.undersampling(Sqx,Sqy,Sqz,options,options['cs']['savepath'])
|
|
|
|
if 'kt-BLAST' in options:
|
|
if options['kt-BLAST']['apply']:
|
|
if rank==0:
|
|
print('Applying kt-BLAST')
|
|
import ktBLAST
|
|
ktBLAST.undersampling(Sqx,Sqy,Sqz,options,options['kt-BLAST']['savepath'])
|
|
|
|
if 'SENSE' in options:
|
|
if options['SENSE']['apply']:
|
|
if rank==0:
|
|
print('-- SENSE reconstruction --')
|
|
import SENSE
|
|
|
|
[Mx,My,Mz] = LOADsequences(options['SENSE']['Mpath'])
|
|
[MxS,MyS,MzS] = SENSE.undersampling(Mx,My,Mz,options)
|
|
|
|
if rank==0:
|
|
print('saving the sequences' + options['SENSE']['savepath'])
|
|
np.savez_compressed(options['SENSE']['savepath'], x=MxS,y=MyS,z=MzS)
|
|
|
|
########################################
|
|
#
|
|
# Writing Checkpoint from Sequence
|
|
#
|
|
########################################
|
|
|
|
if 'create_checkpoint' in options:
|
|
if options['create_checkpoint']['apply']:
|
|
print('--- Create Checkpoint ---')
|
|
|
|
|
|
mesh_into = options['create_checkpoint']['mesh_into']
|
|
|
|
|
|
if options['create_checkpoint']['extension'] == '.mat':
|
|
import scipy.io as sio
|
|
seqpath = options['create_checkpoint']['loadseq']
|
|
seqpath = seqpath + options['create_checkpoint']['seqname']
|
|
seqpath = seqpath + '_R1.mat'
|
|
S = sio.loadmat(seqpath)
|
|
[ddx,ddy,ddz] = S['u_R1']['VoxelSize'][0][0][0]
|
|
[ddx,ddy,ddz] = [ddx/10,ddy/10,ddz/10] #mm to cm
|
|
S = S['u_R1']
|
|
Sqx = S['x'][0][0]
|
|
[Lx,Ly,Lz,Lt] = Sqx.shape
|
|
BOX = LOADmesh('box_zeros', resol=[ddx,ddy,ddz], boxsize=[Lx,Ly,Lz])
|
|
del S, Sqx, seqpath
|
|
else:
|
|
# Box mesh definition
|
|
boxsize = options['create_checkpoint']['boxsize']
|
|
resol = options['create_checkpoint']['resol']
|
|
ranges = options['create_checkpoint']['ranges']
|
|
boxtype = options['create_checkpoint']['boxtype']
|
|
if boxtype == 'fix_resolution':
|
|
BOX = LOADmesh(boxtype, resol=resol, boxsize=boxsize, ranges=ranges)
|
|
|
|
|
|
MESH = LOADmesh('other', meshpath=mesh_into)
|
|
|
|
|
|
for r in options['create_checkpoint']['Rseq']:
|
|
if rank == 0:
|
|
print('Writing checkpoint from sequence ' +
|
|
options['create_checkpoint']['seqname'] + '_R' + str(r) + '.npz')
|
|
|
|
seqpath = options['create_checkpoint']['loadseq'] + options['create_checkpoint']['seqname'] + \
|
|
'_R' + str(r) + options['create_checkpoint']['extension']
|
|
|
|
[Sx, Sy, Sz] = LOADsequences(seqpath)
|
|
|
|
if options['create_checkpoint']['under_rate']>1:
|
|
print('Undersampling with factor: ',
|
|
options['create_checkpoint']['under_rate'])
|
|
Sx = Sx[:, :, :, 0::options['create_checkpoint']['under_rate']]
|
|
Sy = Sy[:, :, :, 0::options['create_checkpoint']['under_rate']]
|
|
Sz = Sz[:, :, :, 0::options['create_checkpoint']['under_rate']]
|
|
|
|
if options['create_checkpoint']['masked']:
|
|
[S0x, S0y, S0z] = LOADsequences(
|
|
options['create_checkpoint']['masked_seq'])
|
|
[Sx, Sy, Sz] = MASKED([S0x, S0y, S0z], [Sx, Sy, Sz])
|
|
del S0x, S0y, S0z
|
|
|
|
vel_seq = SqtoH5(BOX, MESH, -Sx, Sy, -Sz)
|
|
#vel_seq = SqtoH5(BOX, MESH, Sx, Sy, Sz)
|
|
|
|
if rank == 0:
|
|
print(' ')
|
|
|
|
comm = MESH['mesh'].mpi_comm()
|
|
dt = options['create_checkpoint']['dt']
|
|
if options['create_checkpoint']['xdmf']:
|
|
xdmf_u = XDMFFile(options['create_checkpoint']['savepath']+'u.xdmf')
|
|
|
|
for l in range(len(vel_seq)):
|
|
if rank == 0:
|
|
print('saving checkpoint', l)
|
|
|
|
path = options['create_checkpoint']['savepath'] + \
|
|
'R' + str(r) + '/checkpoint/{i}/'.format(i=l)
|
|
|
|
if options['create_checkpoint']['xdmf']:
|
|
vel_seq[l].rename('velocity', 'u')
|
|
xdmf_u.write(vel_seq[l], l*dt)
|
|
|
|
inout.write_HDF5_data(
|
|
comm, path + '/u.h5', vel_seq[l], '/u', t=l*dt)
|
|
|
|
########################################
|
|
#
|
|
# Relative Pressure Estimators
|
|
#
|
|
########################################
|
|
if 'reference' in options:
|
|
if options['reference']['apply']:
|
|
if rank == 0:
|
|
print('---Applying Reference Pressure Gradient---')
|
|
|
|
ckpath = options['reference']['checkpoint_path']
|
|
under_rate = options['reference']['under_rate']
|
|
MESH = LOADmesh('other', meshpath=options['reference']['meshpath'])
|
|
zeros = not 'ct' in ckpath
|
|
pressure = READcheckpoint(
|
|
MESH, 'p', ckpath, under_rate, 'p', zeros=zeros)
|
|
GradP(MESH, pressure, options)
|
|
|
|
if 'ppe' in options:
|
|
|
|
if options['ppe']['apply']:
|
|
[BOX,AORTA,LEO] = CREATEmeshes(options)
|
|
if rank==0:
|
|
print('Applying PPE')
|
|
if options['ppe']['mesh_type']=='aorta':
|
|
velocity_model = READcheckpoint(AORTA,'v',options,options['checkpoint_path'],'u')
|
|
if rank==0:
|
|
print(' ')
|
|
PPE(AORTA,velocity_model,options,options['ppe']['name']+ '.txt')
|
|
if rank==0:
|
|
print(' ')
|
|
|
|
if options['ppe']['mesh_type']=='leo':
|
|
if not options['cs']['apply'] and not options['ppe']['read_check']:
|
|
raise Exception('You need first apply a CS sequence')
|
|
if options['ppe']['read_check']:
|
|
if options['meshname']=='leo':
|
|
velocity = READcheckpoint(LEO,'v',options,options['checkpoint_path'],'u')
|
|
MESH = LEO
|
|
if rank==0:
|
|
print(' ')
|
|
PPE(MESH,velocity,options,options['ppe']['name']+ '.txt')
|
|
if rank==0:
|
|
print(' ')
|
|
|
|
|
|
if not options['ppe']['read_check']:
|
|
for r in options['cs']['R']:
|
|
velocity[r] = SqtoH5(BOX,LEO,Sqx_cs[r],Sqy_cs[r],Sqz_cs[r])
|
|
if rank==0:
|
|
print(' ')
|
|
PPE(LEO,velocity[r],options,options['ppe']['name']+options['ppe']['under_type']+'_R'+str(r)+'.txt')
|
|
if rank==0:
|
|
print(' ')
|
|
|
|
if options['ppe']['mesh_type']=='leoct':
|
|
velocity_model = READcheckpoint(AORTA,'v',options,options['checkpoint_path'],'u')
|
|
if rank==0:
|
|
print(' ')
|
|
velocity_model_leo = CHANGEvel(LEO,velocity_model)
|
|
if rank==0:
|
|
print(' ')
|
|
PPE(LEO,velocity_model_leo,options,options['ppe']['name']+ '.txt')
|
|
if rank==0:
|
|
print(' ')
|
|
|
|
if 'ste' in options:
|
|
|
|
if options['ste']['apply']:
|
|
[BOX,AORTA,LEO] = CREATEmeshes(options)
|
|
if rank==0:
|
|
print('Applying STE')
|
|
|
|
|
|
if options['ste']['mesh_type']=='aorta':
|
|
if velocity_model:
|
|
STE(AORTA,velocity_model,options,options['ste']['name']+ '.txt')
|
|
else:
|
|
velocity_model = READcheckpoint(AORTA,'v',options,options['checkpoint_path'],'u')
|
|
|
|
if rank==0:
|
|
print('')
|
|
|
|
|
|
if options['ste']['int']:
|
|
STEint(AORTA,velocity_model,options,options['ste']['name']+ '.txt')
|
|
else:
|
|
STE(AORTA,velocity_model,options,options['ste']['name']+ '.txt')
|
|
|
|
if options['ste']['mesh_type']=='leo':
|
|
|
|
if not options['cs']['apply'] and not options['ste']['read_check']:
|
|
raise Exception('You need first apply a CS sequence')
|
|
if options['ste']['read_check']:
|
|
if options['meshname']=='leo':
|
|
velocity = READcheckpoint(LEO,'v',options,options['checkpoint_path'],'u')
|
|
MESH = LEO
|
|
if rank==0:
|
|
print(' ')
|
|
|
|
if options['ste']['int']:
|
|
STEint(MESH,velocity,options,options['ste']['name']+ '.txt')
|
|
else:
|
|
STE(MESH,velocity,options,options['ste']['name']+ '.txt')
|
|
if rank==0:
|
|
print(' ')
|
|
|
|
if not options['ste']['read_check']:
|
|
for r in options['cs']['R']:
|
|
if velocity[r]:
|
|
STE(LEO,velocity[r],options,options['ste']['name']+options['ste']['under_type']+'_R'+str(r)+'.txt')
|
|
if rank==0:
|
|
print(' ')
|
|
else:
|
|
velocity[r] = SqtoH5(BOX,LEO,Sqx_cs[r],Sqy_cs[r],Sqz_cs[r])
|
|
if rank==0:
|
|
print(' ')
|
|
if options['ste']['int']:
|
|
STEint(LEO,velocity[r],options,options['ste']['name']+options['ste']['under_type']+'_R'+str(r)+'.txt')
|
|
else:
|
|
STE(LEO,velocity[r],options,options['ste']['name']+options['ste']['under_type']+'_R'+str(r)+'.txt')
|
|
if rank==0:
|
|
print(' ')
|
|
|
|
if options['ste']['mesh_type']=='leoct':
|
|
if velocity_model_leo:
|
|
STE(LEO,velocity_model_leo,options,options['ste']['name']+ '.txt')
|
|
if rank==0:
|
|
print(' ')
|
|
else:
|
|
velocity_model = READcheckpoint(AORTA,'v',options,options['checkpoint_path'],'u')
|
|
if rank==0:
|
|
print(' ')
|
|
velocity_model_leo = CHANGEvel(LEO,velocity_model)
|
|
if rank==0:
|
|
print(' ')
|
|
if options['ste']['int']:
|
|
STEint(LEO,velocity_model_leo,options,options['ste']['name']+ '.txt')
|
|
else:
|
|
STE(LEO,velocity_model_leo,options,options['ste']['name']+ '.txt')
|
|
if rank==0:
|
|
print(' ')
|
|
|
|
if 'peak_pv' in options:
|
|
|
|
if options['peak_pv']['apply']:
|
|
import CS
|
|
import pickle
|
|
import sys
|
|
import logging
|
|
# DPestim
|
|
logging.getLogger().setLevel(logging.INFO)
|
|
parameters['form_compiler']['optimize'] = True
|
|
parameters['form_compiler']['cpp_optimize'] = True
|
|
parameters['form_compiler']['cpp_optimize_flags'] = '-O3 -ffast-math -march=native'
|
|
infile_dp = options['peak_pv']['infile_dp']
|
|
estimator = DPDirectEstim(infile_dp)
|
|
|
|
barye2mmHg = 1/1333.22387415
|
|
t_star = 0.185
|
|
if rank==0:
|
|
print('Computing Velocity and Pressure at Peak Systole')
|
|
print('The inlet max occurs at ' + str(t_star) + ' sec')
|
|
|
|
|
|
[BOX,AORTA,LEO] = CREATEmeshes(options)
|
|
|
|
if options['peak_pv']['mesh_into']=='leo':
|
|
MESH = LEO
|
|
del AORTA
|
|
if options['peak_pv']['mesh_into']=='aorta':
|
|
MESH = AORTA
|
|
del LEO
|
|
|
|
if options['peak_pv']['p_comp']=='error':
|
|
if rank==0:
|
|
print('Reading Pressure Reference')
|
|
P0_PPE = READcheckpoint(MESH,'p',options,options['checkpoint_path'],'p_PPE_impl_stan')
|
|
P0_STE = READcheckpoint(MESH,'p',options,options['checkpoint_path'],'p_STE_impl_stan')
|
|
|
|
[Sqx,Sqy,Sqz] = LOADsequences(options['peak_pv']['orig_seq'])
|
|
Nite = options['peak_pv']['N_CS']
|
|
[row,col,dep,numt] = Sqz.shape
|
|
frms_num = 3
|
|
peakslice = options['peak_pv']['peak_slice']
|
|
R = options['peak_pv']['R']
|
|
|
|
v0 = np.sqrt(Sqx[:,:,:,peakslice]**2 + Sqy[:,:,:,peakslice]**2 + Sqz[:,:,:,peakslice]**2)
|
|
max0 = np.where(v0==np.max(v0))
|
|
|
|
qqvec = {}
|
|
for ss in options['peak_pv']['flux_bnd']:
|
|
qqvec[ss] = np.zeros([Nite])
|
|
|
|
ppemax = np.zeros([Nite])
|
|
stemax = np.zeros([Nite])
|
|
vmax = np.zeros([Nite])
|
|
slrdmax = peakslice
|
|
# Selecting around the max only
|
|
slrdmax = 1
|
|
Sqx = Sqx[:,:,:,peakslice-1:peakslice+2]
|
|
Sqy = Sqy[:,:,:,peakslice-1:peakslice+2]
|
|
Sqz = Sqz[:,:,:,peakslice-1:peakslice+2]
|
|
|
|
for l in range(len(R)):
|
|
if rank==0:
|
|
print('Peak Systole velocity and pressure at R = ' + str(R[l]))
|
|
|
|
sx_cs = np.zeros([row,col,dep,frms_num,Nite])
|
|
sy_cs = np.zeros([row,col,dep,frms_num,Nite])
|
|
sz_cs = np.zeros([row,col,dep,frms_num,Nite])
|
|
|
|
for k in range(Nite):
|
|
if rank==0:
|
|
print('CS iteration number ' + str(k+1))
|
|
[xcs,ycs,zcs] = CS.undersampling_peakpv(Sqx,Sqy,Sqz,options,R[l])
|
|
sx_cs[:,:,:,:,k] = xcs
|
|
sy_cs[:,:,:,:,k] = ycs
|
|
sz_cs[:,:,:,:,k] = zcs
|
|
vk = np.sqrt(sx_cs[:,:,:,1,k]**2 + sy_cs[:,:,:,1,k]**2 + sz_cs[:,:,:,1,k]**2)
|
|
vmax[k] = vk[max0]
|
|
|
|
if rank==0:
|
|
print('\n CS done')
|
|
|
|
# To write the checkpoints
|
|
vel_seq = SqtoH5(BOX,MESH,sx_cs[:,:,:,:,k],sy_cs[:,:,:,:,k],sz_cs[:,:,:,:,k])
|
|
comm = MESH['mesh'].mpi_comm()
|
|
|
|
# Computing the Fluxes
|
|
if rank==0:
|
|
print('\n Computing the Flux')
|
|
QQ = Fluxes(MESH,vel_seq,options,options['peak_pv']['flux_bnd'])
|
|
|
|
for ss in options['peak_pv']['flux_bnd']:
|
|
qqvec[ss][k] = QQ[ss][slrdmax]
|
|
|
|
if rank==0:
|
|
print('\n Writing checkpoints')
|
|
|
|
for ns in range(len(vel_seq)):
|
|
pathss = options['peak_pv']['savepath'] + 'H5/checkpoint/{i}/'.format(i=ns)
|
|
if l<10 and l>0:
|
|
pathss = options['peak_pv']['savepath'] + 'H5/checkpoint/0{i}/'.format(i=ns)
|
|
inout.write_HDF5_data(comm, pathss + '/u.h5', vel_seq[ns], '/u', t=0)
|
|
if rank==0:
|
|
print('\n The checkpoints were wrote')
|
|
|
|
# Computing the Pressure Drop
|
|
estimator.estimate()
|
|
# Reading the results
|
|
if options['peak_pv']['p_comp']=='peak':
|
|
ppe_raw = open(options['peak_pv']['savepath'] + '/H5/pdrop_PPE_impl_stan.dat','rb')
|
|
ste_raw = open(options['peak_pv']['savepath'] + '/H5/pdrop_STE_impl_stan.dat','rb')
|
|
ppe = pickle.load(ppe_raw)['pdrop']*(-barye2mmHg)
|
|
ste = pickle.load(ste_raw)['pdrop']*(-barye2mmHg)
|
|
p1max[k] = ppe[slrdmax]
|
|
p2max[k] = ste[slrdmax]
|
|
elif options['peak_pv']['p_comp']=='error':
|
|
PPE = READcheckpoint(MESH,'p',options,options['peak_pv']['savepath']+'H5/checkpoint/','p_PPE_impl_stan')
|
|
STE = READcheckpoint(MESH,'p',options,options['peak_pv']['savepath']+'H5/checkpoint/','p_STE_impl_stan')
|
|
ppe_vec_0 = P0_PPE[peakslice].vector().get_local() - P0_PPE[peakslice].vector().get_local()[0]
|
|
ppe_vec = PPE[slrdmax].vector().get_local() - PPE[slrdmax].vector().get_local()[0]
|
|
ste_vec_0 = P0_STE[peakslice].vector().get_local() - P0_STE[peakslice].vector().get_local()[0]
|
|
ste_vec = STE[slrdmax].vector().get_local() - STE[slrdmax].vector().get_local()[0]
|
|
ppemax[k] = np.linalg.norm(ppe_vec_0 - ppe_vec)/np.linalg.norm(ppe_vec_0)
|
|
stemax[k] = np.linalg.norm(ste_vec_0 - ste_vec)/np.linalg.norm(ste_vec_0)
|
|
else:
|
|
raise Exception('Pressure computation not recognize!')
|
|
|
|
|
|
# VELOCITIES
|
|
vmean = np.mean(vmax)
|
|
vstd = np.std(vmax)
|
|
# PRESSURES
|
|
ppemean = np.mean(ppemax)
|
|
stemean = np.mean(stemax)
|
|
ppestd = np.std(ppemax)
|
|
stestd = np.std(stemax)
|
|
|
|
|
|
if options['peak_pv']['save']:
|
|
if rank==0:
|
|
print('\n saving the files in ' + options['peak_pv']['savepath'])
|
|
for ss in options['peak_pv']['flux_bnd']:
|
|
np.savetxt( options['peak_pv']['savepath'] + 'Q'+str(ss) +'_R'+str(R[l])+'.txt', [np.mean(qqvec[ss])])
|
|
np.savetxt( options['peak_pv']['savepath'] + 'Qstd'+str(ss) +'_R'+str(R[l])+'.txt', [np.std(qqvec[ss])])
|
|
np.savetxt( options['peak_pv']['savepath'] + 'ppemean_R'+str(R[l])+'.txt', [ppemean] )
|
|
np.savetxt( options['peak_pv']['savepath'] + 'stemean_R'+str(R[l])+'.txt', [stemean] )
|
|
np.savetxt( options['peak_pv']['savepath'] + 'vmean_R'+str(R[l])+'.txt', [vmean] )
|
|
np.savetxt( options['peak_pv']['savepath'] + 'ppestd_R'+str(R[l])+'.txt', [ppestd] )
|
|
np.savetxt( options['peak_pv']['savepath'] + 'stestd_R'+str(R[l])+'.txt', [stestd] )
|
|
np.savetxt( options['peak_pv']['savepath'] + 'vstd_R'+str(R[l])+'.txt', [vstd] )
|
|
|
|
if 'change_mesh' in options:
|
|
if options['change_mesh']['apply']:
|
|
if rank == 0:
|
|
print('Changing '+options['change_mesh']['mode'] + ' from: ')
|
|
print('Mesh input: ' + options['change_mesh']['mesh_in'])
|
|
print('Mesh output: ' + options['change_mesh']['mesh_out'])
|
|
|
|
MESH_in = LOADmesh(
|
|
'other', meshpath=options['change_mesh']['mesh_in'])
|
|
MESH_out = LOADmesh(
|
|
'other', meshpath=options['change_mesh']['mesh_out'])
|
|
|
|
ckpath = options['change_mesh']['checkpoint_path']
|
|
under_rate = options['change_mesh']['under_rate']
|
|
mode = options['change_mesh']['mode']
|
|
|
|
origin = READcheckpoint(MESH_in, mode, ckpath, under_rate, mode, zeros=False)
|
|
|
|
# To change between meshes the same result
|
|
changed = {}
|
|
#W = LEO['FEM'].sub(0).collapse()
|
|
comm = MESH_out['mesh'].mpi_comm()
|
|
v2 = Function(MESH_out['FEM'])
|
|
|
|
for k in range(len(list(origin))):
|
|
if rank == 0:
|
|
print('CHANGING: index', k)
|
|
LagrangeInterpolator.interpolate(v2, origin[k])
|
|
changed[k] = v2
|
|
print(max(origin[k].vector().get_local()),max(v2.vector().get_local()))
|
|
|
|
|
|
dt = options['change_mesh']['dt']
|
|
if options['create_checkpoint']['xdmf']:
|
|
xdmf_u = XDMFFile(options['change_mesh']['savepath']+'u.xdmf')
|
|
for l in range(len(changed)):
|
|
if rank == 0:
|
|
print('saving checkpoint', l)
|
|
path = options['change_mesh']['savepath'] + \
|
|
'R1/checkpoint/{i}/'.format(i=l)
|
|
writepath = path + '/'+options['change_mesh']['mode']+'.h5'
|
|
inout.write_HDF5_data(
|
|
comm, writepath, changed[l] ,'/'+options['change_mesh']['mode'], t=l*dt)
|
|
|
|
if options['create_checkpoint']['xdmf']:
|
|
changed[l].rename('velocity', 'u')
|
|
xdmf_u.write(changed[l], l*dt)
|
|
|
|
|
|
|
|
if __name__ == '__main__':
|
|
|
|
comm = MPI.COMM_WORLD
|
|
size = comm.Get_size()
|
|
rank = comm.Get_rank()
|
|
|
|
if len(sys.argv) > 1:
|
|
if os.path.exists(sys.argv[1]):
|
|
inputfile = sys.argv[1]
|
|
if rank==0:
|
|
print('Found input file ' + inputfile)
|
|
else:
|
|
Etcetera('serious')
|
|
raise Exception('Command line arg given but input file does not exist:'
|
|
' {}'.format(sys.argv[1]))
|
|
else:
|
|
Etcetera('serious')
|
|
raise Exception('An input file is required as argument!')
|
|
|
|
|
|
user = 'p283370' # Default's user
|
|
if 'Zion' in os.popen('hostname').read():
|
|
user = 'yeye'
|
|
np.set_printoptions(threshold=5)
|
|
|
|
if 'fwn-bborg-5-166' in os.popen('hostname').read():
|
|
user = 'p283370'
|
|
|
|
if rank==0:
|
|
print('Welcome user {uss}'.format(uss=user))
|
|
Etcetera('happy')
|
|
|
|
|
|
start_time = time.time()
|
|
|
|
options = inout.read_parameters(inputfile)
|
|
SCANNER(options)
|
|
|
|
end_time = time.time()
|
|
CLOCK(rank,start_time,end_time)
|
|
|