import time import math from dynamico.dev import wrap from dynamico.dev.libs import libicosa from dynamico.util import list_stencil from ctypes import c_void_p, c_int, c_double, c_float, c_bool from numpy cimport ndarray cimport numpy as np import numpy as np from dynamico import getargs log_master, log_world = getargs.getLogger(__name__) INFO, DEBUG, ERROR = log_master.info, log_master.debug, log_world.error INFO_ALL, DEBUG_ALL = log_world.info, log_world.debug #-------------- choose precision of kernel computations ------------# # the compile-time constant CPP_MIXED_PREC is set by setup.py # based on the value of mixed_prec defined in module data_unstructured IF CPP_MIXED_PREC: c_num=c_float ctypedef float num np_num=np.float32 ELSE: c_num=c_double ctypedef double num np_num=np.float64 ctypedef num *num_ptr #------------- direct Cython interface to DYNAMICO routines -------------# cdef enum: max_nb_stage=5 cdef extern : cdef num tauj[max_nb_stage] cdef num cslj[max_nb_stage][max_nb_stage] cdef num cflj[max_nb_stage][max_nb_stage] cdef int nb_stage[1] cdef extern from "functions.h": cdef void dynamico_ARK_step(int nstep, num *mass_col, num *rhodz, num *theta_rhodz, num *u, num *geopot, num *w, num *theta, num *ps, num *pk, num *hflux, num *qv, num *dmass_col, num *drhodz, num *dtheta_rhodz, num *du_fast, num *du_slow, num *dPhi_fast, num *dPhi_slow, num *dW_fast, num *dW_slow) cdef void dynamico_remap(num *rhodz, num *theta_rhodz, num *u) cdef void dynamico_init_params() cpdef void dynamico_setup_xios() cpdef void dynamico_xios_set_timestep(double) cpdef void dynamico_xios_update_calendar(int) #------------- import and wrap DYNAMICO routines -------------# ker=wrap.Struct() # store imported fun X as funs.X check_args = False # use True instead of False for debugging, probably with some overhead try: kernels = wrap.SharedLib(vars(ker), libicosa, check_args=check_args) setvar, setvars, getvar, getvars = kernels.setvar, kernels.setvars, kernels.getvar, kernels.getvars except OSError: ERROR("Unable to load shared library 'libicosa.so' !") raise # providing a full prototype enables type-checking when calling # if a number n is present in the prototype, the previous type is repeated n times kernels.import_funs([ ['dynamico_setup_xios',None], ['dynamico_print_trace',None], ['dynamico_xios_set_timestep',c_double], ['dynamico_xios_update_calendar',c_int], ['dynamico_init_mesh',c_void_p,13], ['dynamico_init_metric', c_void_p,6], ['dynamico_init_hybrid', c_void_p,3], ['dynamico_caldyn_unstructured', c_num, c_void_p,20], ['dynamico_scalar_laplacian', c_void_p,2], ['dynamico_curl_laplacian', c_void_p,2], ['dynamico_init_transfer', c_int, c_int,2,c_void_p,3, c_int,2,c_void_p,3], ['dynamico_update_halo', c_int,3,c_void_p] ]) # set/get global variables eta_mass,eta_lag=(1,2) thermo_theta,thermo_entropy,thermo_moist,thermo_boussinesq=(1,2,3,4) kernels.addvars( c_bool,'hydrostatic','debug_hevi_solver', 'is_mpi_master', c_int,'llm','nqdyn','primal_num','max_primal_deg', 'dual_num','max_dual_deg','edge_num','max_trisk_deg', 'caldyn_thermo','caldyn_eta','nb_threads', c_double,'elapsed', c_num, 'g', 'ptop', 'cpp', 'cppv', 'Rd', 'Rv', 'preff', 'Treff', 'pbot', 'rho_bot', 'Phi_bot') elapsed=0. #------------------------ Extension type performing a full ARK time step ---------------------- cdef num_ptr ptr1(num[:] data) except *: return &data[0] cdef num_ptr ptr2(num[:,:] data) except *: return &data[0,0] cdef num_ptr ptr3(num[:,:,:] data) except *: return &data[0,0,0] cdef num_ptr ptr4(num[:,:,:,:] data) except *: return &data[0,0,0,0] cdef num_ptr ptr(data) except * : n=data.ndim if n==1 : return ptr1(data) if n==2 : return ptr2(data) if n==3 : return ptr3(data) if n==4 : return ptr4(data) if n>4: raise IndexError cdef alloc(num_ptr *p, allocator, n=1): data=allocator(n) p[0]=ptr(data) return data cdef class Caldyn_step: # number of time steps to do at each invocation of advance() cdef int nstep # pointer to allocated arrays cdef num_ptr p_mass, p_theta_rhodz, p_u, p_geopot, p_W # prognostic cdef num_ptr p_mass_col, p_dmass_col, p_ps, p_theta, p_pk, p_hflux, p_qv # diagnostic cdef num_ptr p_drhodz, p_dtheta_rhodz, p_du_fast, p_du_slow # tendencies cdef num_ptr p_dPhi_fast, p_dPhi_slow, p_dW_fast, p_dW_slow # tendencies # allocated arrays, must remain referenced or segfault cdef readonly ndarray mass, theta_rhodz, u, geopot, W cdef readonly ndarray mass_col, dmass_col, ps, theta, pk, hflux, qv cdef readonly ndarray drhodz, dtheta_rhodz, du_fast, du_slow cdef readonly ndarray dPhi_fast, dPhi_slow, dW_fast, dW_slow def __init__(self,mesh,time_scheme, nstep): self.nstep=nstep # self.mesh=mesh fps, ftheta, fmass = mesh.field_ps, mesh.field_theta, mesh.field_mass fw, fu, fz = mesh.field_w, mesh.field_u, mesh.field_z # collect coefficients of time scheme cdef double[:] tauj_ = time_scheme.tauj cdef double[:,:] cslj_ = time_scheme.csjl cdef double[:,:] cflj_ = time_scheme.cfjl ns = time_scheme.nstage nb_stage[0]=ns cdef int i,j for i in range(ns): tauj[i]=tauj_[i] for j in range(ns): cslj[i][j]=cslj_[i,j] cflj[i][j]=cflj_[i,j] # allocate arrays, store pointers to avoid overhead when calling dynamico # prognostic/diagnostic self.ps = alloc(&self.p_ps, fps) self.mass_col, self.dmass_col = alloc(&self.p_mass_col, fps), alloc(&self.p_dmass_col, fps,ns), self.mass, self.theta_rhodz = alloc(&self.p_mass, fmass), alloc(&self.p_theta_rhodz, fmass), self.theta, self.pk = alloc(&self.p_theta, fmass), alloc(&self.p_pk, fmass), self.geopot, self.W = alloc(&self.p_geopot, fw), alloc(&self.p_W, fw), self.hflux, self.u = alloc(&self.p_hflux, fu), alloc(&self.p_u, fu) self.qv = alloc(&self.p_qv,fz) # tendencies self.drhodz, self.dtheta_rhodz = alloc(&self.p_drhodz,fmass,ns), alloc(&self.p_dtheta_rhodz,fmass,ns) self.du_fast, self.du_slow = alloc(&self.p_du_fast,fu,ns), alloc(&self.p_du_slow,fu,ns) self.dPhi_fast, self.dPhi_slow = alloc(&self.p_dPhi_fast,fw,ns), alloc(&self.p_dPhi_slow,fw,ns) self.dW_fast, self.dW_slow = alloc(&self.p_dW_fast,fw,ns), alloc(&self.p_dW_slow,fw,ns) def next(self): # global elapsed # time1=time.time() dynamico_ARK_step(self.nstep, self.p_mass_col, self.p_mass, self.p_theta_rhodz, self.p_u, self.p_geopot, self.p_W, self.p_theta, self.p_ps, self.p_pk, self.p_hflux, self.p_qv, self.p_dmass_col, self.p_drhodz, self.p_dtheta_rhodz, self.p_du_fast, self.p_du_slow, self.p_dPhi_fast, self.p_dPhi_slow, self.p_dW_fast, self.p_dW_slow) #time2=time.time() #if time2>time1: elapsed=elapsed+time2-time1 def remap(self): dynamico_remap(self.p_mass, self.p_theta_rhodz, self.p_u) def caldyn_step_TRSW(mesh,time_scheme,nstep): setvars(('hydrostatic','caldyn_thermo','caldyn_eta'), (True,thermo_boussinesq,eta_lag)) dynamico_init_params() return Caldyn_step(mesh,time_scheme, nstep) def caldyn_step_HPE(mesh,time_scheme,nstep, caldyn_thermo,caldyn_eta, thermo,BC,g): setvars(('hydrostatic','caldyn_thermo','caldyn_eta', 'g','ptop','Rd','cpp','preff','Treff'), (True,caldyn_thermo,caldyn_eta, g,BC.ptop,thermo.Rd,thermo.Cpd,thermo.p0,thermo.T0)) dynamico_init_params() return Caldyn_step(mesh,time_scheme, nstep) def caldyn_step_NH(mesh,time_scheme,nstep, caldyn_thermo, caldyn_eta, thermo,BC,g): setvars(('hydrostatic','caldyn_thermo','caldyn_eta', 'g','ptop','Rd','cpp','preff','Treff','pbot','rho_bot'), (False,caldyn_thermo,caldyn_eta, g,BC.ptop,thermo.Rd,thermo.Cpd,thermo.p0,thermo.T0, BC.pbot.max(), BC.rho_bot.max())) dynamico_init_params() return Caldyn_step(mesh,time_scheme, nstep) #----------------------------- Base class for dynamics ------------------------ class Caldyn: def __init__(self,mesh): self.mesh=mesh fps, ftheta, fmass = mesh.field_ps, mesh.field_theta, mesh.field_mass fw, fu, fz = mesh.field_w, mesh.field_u, mesh.field_z self.ps, self.ms, self.dms = fps(), fps(), fps() self.s, self.hs, self.dhs = ftheta(), ftheta(), ftheta() self.pk, self.berni, self.geopot, self.hflux = fmass(),fmass(),fw(),fu() self.qu, self.qv = fu(),fz() self.fmass, self.ftheta, self.fu, self.fw = fmass, ftheta, fu, fw def bwd_fast_slow(self, flow, tau): global elapsed time1=time.time() flow,fast,slow = self._bwd_fast_slow_(flow,tau) time2=time.time() elapsed=elapsed+time2-time1 return flow,fast,slow # when calling caldyn_unstructured, arrays for tendencies must be re-created each time # to avoid overwriting in the same memory space when time scheme is multi-stage #-------------------------- Shallow-water dynamics --------------------- class Caldyn_RSW(Caldyn): def __init__(self,mesh): Caldyn.__init__(self,mesh) setvars(('hydrostatic','caldyn_thermo','caldyn_eta'), (True,thermo_boussinesq,eta_lag)) self.dhs = self.fmass() dynamico_init_params() def _bwd_fast_slow_(self, flow, tau): h,u = flow # h*s = h => uniform buoyancy s=1 => shallow-water dh, du_slow, du_fast, hs, buf = self.fmass(), self.fu(), self.fu(), h.copy(), self.geopot assert type(tau) is np_num, 'tau must be of type unstructured.np_num.' ker.dynamico_caldyn_unstructured(tau, self.ms, h, hs, u, self.geopot, buf, self.s, self.ps, self.pk, self.hflux, self.qv, self.dms, dh, self.dhs, du_fast, du_slow, buf, buf, buf, buf) return (h,u), (0.,du_fast), (dh,du_slow) #----------------------------------- HPE ------------------------------------ class Caldyn_HPE(Caldyn): def __init__(self,caldyn_thermo,caldyn_eta, mesh,thermo,BC,g): Caldyn.__init__(self,mesh) setvars(('hydrostatic','caldyn_thermo','caldyn_eta', 'g','ptop','Rd','cpp','preff','Treff'), (True,caldyn_thermo,caldyn_eta, g,BC.ptop,thermo.Rd,thermo.Cpd,thermo.p0,thermo.T0)) dynamico_init_params() def _bwd_fast_slow_(self, flow, tau): dm, dS, du_slow, du_fast, buf = self.fmass(), self.ftheta(), self.fu(), self.fu(), self.geopot m,S,u = flow ker.dynamico_caldyn_unstructured(tau, self.ms, m, S, u, self.geopot, buf, self.s, self.ps, self.pk, self.hflux, self.qv, self.dms, dm, dS, du_fast, du_slow, buf, buf, buf, buf) return (m,S,u), (0.,0.,du_fast), (dm,dS,du_slow) #----------------------------------- NH ------------------------------------ class Caldyn_NH(Caldyn): def __init__(self,caldyn_thermo,caldyn_eta, mesh,thermo,BC,g): Caldyn.__init__(self,mesh) setvars(('hydrostatic','caldyn_thermo','caldyn_eta', 'g','ptop','Rd','cpp','preff','Treff', 'pbot','rho_bot'), (False,caldyn_thermo,caldyn_eta, g,BC.ptop,thermo.Rd,thermo.Cpd,thermo.p0,thermo.T0, BC.pbot.max(), BC.rho_bot.max())) dynamico_init_params() def bwd_fast_slow(self, flow, tau): ftheta, fmass, fu, fw = self.ftheta, self.fmass, self.fu, self.fw dm, dS, du_slow, du_fast = fmass(), ftheta(), fu(), fu() dPhi_slow, dPhi_fast, dW_slow, dW_fast = fw(), fw(), fw(), fw() m,S,u,Phi,W = flow ker.dynamico_caldyn_unstructured(tau, self.ms, m, S, u, Phi, W, self.s, self.ps, self.pk, self.hflux, self.qv, self.dms, dm, dS, du_fast, du_slow, dPhi_fast, dPhi_slow, dW_fast, dW_slow) return ((m,S,u,Phi,W), (0.,0.,du_fast,dPhi_fast,dW_fast), (dm,dS,du_slow,dPhi_slow,dW_slow)) #------------------------ Copy mesh info to Fortran side ------------------- def init_mesh(llm, nqdyn, edge_num, primal_num, dual_num, max_trisk_deg, max_primal_deg, max_dual_deg, primal_nb, primal_edge, primal_ne, dual_nb,dual_edge,dual_ne,dual_vertex, left,right,down,up,trisk_deg,trisk, Ai, Av, fv, le_de, Riv2, wee): DEBUG('Types of Ai, Av, fv, le_de, Riv2, wee : %s' % ((Ai.dtype,Av.dtype,fv.dtype,le_de.dtype,Riv2.dtype,wee.dtype),) ) for var,varname in zip((Ai,Av,fv,le_de,Riv2,wee), ('Ai','Av','fv','le_de','Riv2','wee')): assert var.dtype == np.float64, '%s must be double precision'%varname setvars( ('llm','nqdyn','edge_num','primal_num','dual_num', 'max_trisk_deg','max_primal_deg','max_dual_deg'), (llm, nqdyn, edge_num, primal_num, dual_num, max_trisk_deg, max_primal_deg, max_dual_deg) ) INFO('init_mesh ...') ker.dynamico_init_mesh(primal_nb,primal_edge,primal_ne, dual_nb,dual_edge,dual_ne,dual_vertex, left,right,down,up,trisk_deg,trisk) INFO('...done') INFO('init_metric ...') ker.dynamico_init_metric(Ai,Av,fv,le_de,Riv2,wee) INFO('...done') #------------------------ Mesh partitioning ------------------------ # Helper functions and interface to ParMETIS # loc_stencil returns the start/end indices (vtxdist) expected by ParMETIS # i.e. index[start:end] with start=vtxdist[cell], end=vtxdist[cell+1] lists the edges of cell 'cell' def loc_stencil(degree, stencil): loc=0 for i in range(degree.size): yield loc loc=loc+degree[i] yield loc def partition_mesh(degree, stencil, nparts): # arguments : PArray1D and PArray2D describing mesh, number of desired partitions dim_cell, degree, stencil = degree.dim, degree.data, stencil.data comm, vtxdist, idx_start, idx_end = dim_cell.comm, dim_cell.vtxdist, dim_cell.start, dim_cell.end mpi_rank, mpi_size = comm.Get_rank(), comm.Get_size() adjncy_loc, xadj_loc = list_stencil(degree, stencil), loc_stencil(degree, stencil) adjncy_loc, xadj_loc = [np.asarray(list(x), dtype=np.int32) for x in (adjncy_loc, xadj_loc)] owner = np.zeros(idx_end-idx_start, dtype=np.int32); ker.dynamico_partition_graph(mpi_rank, mpi_size, vtxdist, xadj_loc, adjncy_loc, nparts, owner) return owner