source: codes/icosagcm/devel/Python/src/unstructured.pyx @ 631

Last change on this file since 631 was 631, checked in by dubos, 6 years ago

devel/Python : extract pure Python stuff from cython module unstructured.pyx

File size: 8.4 KB
Line 
1import time
2import math
3import numpy as np
4import dynamico.wrap as wrap
5from ctypes import c_void_p, c_int, c_double, c_bool
6
7#------------- direct Cython interface to DYNAMICO routines -------------#
8
9cdef extern from "functions.h":
10     cdef void dynamico_init_params()
11     cpdef void dynamico_setup_xios()
12     cpdef void dynamico_xios_set_timestep(double)
13     cpdef void dynamico_xios_update_calendar(int)
14
15#------------- import and wrap DYNAMICO routines -------------#
16
17ker=wrap.Struct() # store imported fun X as funs.X
18
19check_args = False # use True instead of False for debugging, probably with some overhead
20
21try:   
22    kernels = wrap.SharedLib(vars(ker), 'libicosa.so', check_args=check_args) 
23    setvar, setvars, getvar, getvars = kernels.setvar, kernels.setvars, kernels.getvar, kernels.getvars
24except OSError:
25    print """
26Unable to load shared library 'libkernels.so' !
27    """
28    raise
29
30# providing a full prototype enables type-checking when calling
31# if a number n is present in the prototype, the previous type is repeated n times
32kernels.import_funs([
33                     ['dynamico_setup_xios',None],
34                     ['dynamico_xios_set_timestep',c_double],
35                     ['dynamico_xios_update_calendar',c_int],
36                     ['dynamico_init_mesh',c_void_p,13],
37                     ['dynamico_init_metric', c_void_p,6],
38                     ['dynamico_init_hybrid', c_void_p,3],
39                     ['dynamico_caldyn_unstructured', c_double, c_void_p,20],
40                     ['dynamico_partition_graph', c_int,2, c_void_p,3, c_int, c_void_p],
41                     ])
42
43# set/get global variables
44eta_mass,eta_lag=(1,2)
45thermo_theta,thermo_entropy,thermo_moist,thermo_boussinesq=(1,2,3,4)
46
47kernels.addvars(
48    c_bool,'hydrostatic','debug_hevi_solver','rigid',
49    c_int,'llm','nqdyn','primal_num','max_primal_deg',
50    'dual_num','max_dual_deg','edge_num','max_trisk_deg',
51    'caldyn_thermo','caldyn_eta','nb_threads',
52    c_double,'elapsed','g', 'ptop', 'cpp', 'cppv',
53    'Rd', 'Rv', 'preff', 'Treff', 'pbot', 'rho_bot', 'Phi_bot')
54
55elapsed=0.
56
57#----------------------------- Base class for dynamics ------------------------
58
59class Caldyn:
60    def __init__(self,mesh):
61        self.mesh=mesh
62        fps, ftheta, fmass = mesh.field_ps, mesh.field_theta, mesh.field_mass
63        fw, fu, fz = mesh.field_w, mesh.field_u, mesh.field_z
64        self.ps, self.ms, self.dms       = fps(), fps(), fps()
65        self.s, self.hs, self.dhs        = ftheta(), ftheta(), ftheta()
66        self.pk, self.berni, self.geopot, self.hflux = fmass(),fmass(),fw(),fu()
67        self.qu, self.qv                 = fu(),fz()
68        self.fmass, self.ftheta, self.fu, self.fw = fmass, ftheta, fu, fw
69    def bwd_fast_slow(self, flow, tau):
70        global elapsed
71        time1=time.time()
72        flow,fast,slow = self._bwd_fast_slow_(flow,tau)
73        time2=time.time()
74        elapsed=elapsed+time2-time1
75        return flow,fast,slow
76
77# when calling caldyn_unstructured, arrays for tendencies must be re-created each time
78# to avoid overwriting in the same memory space when time scheme is multi-stage
79
80#-------------------------- Shallow-water dynamics ---------------------
81
82class Caldyn_RSW(Caldyn):
83    def __init__(self,mesh):
84        Caldyn.__init__(self,mesh)
85        setvars(('hydrostatic','caldyn_thermo','caldyn_eta'),
86                (True,thermo_boussinesq,eta_lag))
87        self.dhs = self.fmass()
88        dynamico_init_params()
89    def _bwd_fast_slow_(self, flow, tau):
90        h,u = flow
91        # h*s = h => uniform buoyancy s=1 => shallow-water
92        dh, du_slow, du_fast, hs, buf = self.fmass(), self.fu(), self.fu(), h.copy(), self.geopot
93        ker.dynamico_caldyn_unstructured(tau, self.ms, h, hs, u, self.geopot, buf,
94                  self.s, self.ps, self.pk, self.hflux, self.qv,
95                  self.dms, dh, self.dhs, du_fast, du_slow,
96                  buf, buf, buf, buf)
97        return (h,u), (0.,du_fast), (dh,du_slow)
98
99#----------------------------------- HPE ------------------------------------
100
101class Caldyn_HPE(Caldyn):
102    def __init__(self,caldyn_thermo,caldyn_eta, mesh,metric,thermo,BC,g):
103        Caldyn.__init__(self,mesh)
104        setvars(('hydrostatic','caldyn_thermo','caldyn_eta',
105                 'g','ptop','Rd','cpp','preff','Treff'),
106                (True,caldyn_thermo,caldyn_eta,
107                 g,BC.ptop,thermo.Rd,thermo.Cpd,thermo.p0,thermo.T0))
108        dynamico_init_params()
109    def _bwd_fast_slow_(self, flow, tau):
110        dm, dS, du_slow, du_fast, buf = self.fmass(), self.ftheta(), self.fu(), self.fu(), self.geopot
111        m,S,u = flow
112        ker.dynamico_caldyn_unstructured(tau, self.ms, m, S, u, self.geopot, buf,
113                  self.s, self.ps, self.pk, self.hflux, self.qv,
114                  self.dms, dm, dS, du_fast, du_slow,
115                  buf, buf, buf, buf)
116        return (m,S,u), (0.,0.,du_fast), (dm,dS,du_slow)
117
118#----------------------------------- NH ------------------------------------
119
120class Caldyn_NH(Caldyn):
121    def __init__(self,caldyn_thermo,caldyn_eta, mesh,metric,thermo,BC,g):
122        Caldyn.__init__(self,mesh)
123        setvars(('hydrostatic','caldyn_thermo','caldyn_eta',
124                 'g','ptop','Rd','cpp','preff','Treff',
125                 'pbot','rho_bot'),
126                (False,caldyn_thermo,caldyn_eta,
127                 g,BC.ptop,thermo.Rd,thermo.Cpd,thermo.p0,thermo.T0,
128                 BC.pbot.max(), BC.rho_bot.max()))
129        dynamico_init_params()
130    def bwd_fast_slow(self, flow, tau):
131        ftheta, fmass, fu, fw = self.ftheta, self.fmass, self.fu, self.fw
132        dm, dS, du_slow, du_fast = fmass(), ftheta(), fu(), fu()
133        dPhi_slow, dPhi_fast, dW_slow, dW_fast = fw(), fw(), fw(), fw()
134        m,S,u,Phi,W = flow
135        ker.dynamico_caldyn_unstructured(tau, self.ms, m, S, u, Phi, W,
136                  self.s, self.ps, self.pk, self.hflux, self.qv,
137                  self.dms, dm, dS, du_fast, du_slow,
138                  dPhi_fast, dPhi_slow, dW_fast, dW_slow)
139        return ((m,S,u,Phi,W), (0.,0.,du_fast,dPhi_fast,dW_fast), 
140                (dm,dS,du_slow,dPhi_slow,dW_slow))
141
142#------------------------ Copy mesh info to Fortran side -------------------
143
144def init_mesh(llm, nqdyn, edge_num, primal_num, dual_num,
145              max_trisk_deg, max_primal_deg, max_dual_deg,
146              primal_nb, primal_edge, primal_ne,
147              dual_nb,dual_edge,dual_ne,dual_vertex, 
148              left,right,down,up,trisk_deg,trisk,
149              Ai, Av, fv, le_de, Riv2, wee):
150    setvars( ('llm','nqdyn','edge_num','primal_num','dual_num',
151              'max_trisk_deg','max_primal_deg','max_dual_deg'),
152             (llm, nqdyn, edge_num, primal_num, dual_num, 
153              max_trisk_deg, max_primal_deg, max_dual_deg) )       
154    print('init_mesh ...')
155    ker.dynamico_init_mesh(primal_nb,primal_edge,primal_ne,
156              dual_nb,dual_edge,dual_ne,dual_vertex,
157              left,right,down,up,trisk_deg,trisk)
158    print ('...done')
159    print('init_metric ...')
160    ker.dynamico_init_metric(Ai,Av,fv,le_de,Riv2,wee)
161    print ('...done')
162
163#------------------------ Mesh partitioning ------------------------
164
165# Helper functions and interface to ParMETIS
166# list_stencil converts an adjacency graph from array format index[num_cells, MAX_EDGES] to compressed format
167# loc_stencil returns the start/end indices (vtxdist) expected by ParMETIS
168# i.e. index[start:end] with start=vtxdist[cell], end=vtxdist[cell+1] lists the edges of cell 'cell'
169
170def list_stencil(degree, stencil, cond=lambda x:True):
171    for i in range(degree.size):
172        for j in range(degree[i]):
173            s=stencil[i,j]
174            if cond(s): yield stencil[i,j]
175               
176def loc_stencil(degree, stencil):
177    loc=0
178    for i in range(degree.size):
179        yield loc
180        loc=loc+degree[i]
181    yield loc
182
183def partition_mesh(degree, stencil, nparts): 
184    # arguments : PArray1D and PArray2D describing mesh, number of desired partitions
185    dim_cell, degree, stencil = degree.dim, degree.data, stencil.data
186    comm, vtxdist, idx_start, idx_end = dim_cell.comm, dim_cell.vtxdist, dim_cell.start, dim_cell.end
187    mpi_rank, mpi_size = comm.Get_rank(), comm.Get_size()
188    adjncy_loc, xadj_loc = list_stencil(degree, stencil), loc_stencil(degree, stencil)
189    adjncy_loc, xadj_loc = [np.asarray(list(x), dtype=np.int32) for x in (adjncy_loc, xadj_loc)]
190    owner = np.zeros(idx_end-idx_start, dtype=np.int32);
191    ker.dynamico_partition_graph(mpi_rank, mpi_size, vtxdist, xadj_loc, adjncy_loc, nparts, owner)
192    return owner
Note: See TracBrowser for help on using the repository browser.