All Classes Files Functions Variables Pages
partCubitMesh.c File Reference

This file contains C subroutines to partition a mesh with (METIS-5.1.0 library) and to initialize GLL nodes global numbering. Mesh file formats surported: ABAQUS *.inp and EXODUS_II *.ex2. More...

#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <math.h>
#include <assert.h>
#include <ctype.h>
#include "partCubitMesh.h"
#include "mpi.h"
#include "exodusII.h"
Include dependency graph for partCubitMesh.c:

Go to the source code of this file.

Functions

int count_elt (FILE *unit)
 
void setHexaWeight (mesh_t *mesh)
 
void readCubitMeshAbaqus (char *fileinp, mesh_t *mesh, int *number_of_elemnt_block)
 
void readCubitMeshExodus (char *fileinp, mesh_t *mesh, int *number_of_elemnt_block)
 
void read_part_cubit_mesh_ (char *fileinp, int *ig_ncpu, int *ig_myrank, int *size_real_t, int *ig_mesh_nnode, int *ig_ncpu_neighbor, int *ig_nhexa, int *ig_nhexa_outer, int *ig_nquad_parax, int *ig_nquad_fsurf, int *ig_hexa_nnode, int *ig_quad_nnode, float *rg_mesh_xmax, float *rg_mesh_ymax, float *rg_mesh_zmax, float *rg_mesh_xmin, float *rg_mesh_ymin, float *rg_mesh_zmin)
 
int fill_mesh_arrays_ (int *ig_ncpu, int *ig_myrank, int *ig_ngll_total, int *ig_nneighbor_all_kind, int *cpu_neighbor, int *ig_cpu_neighbor_info, int *ig_hexa_gnode_glonum, int *ig_quadp_gnode_glonum, int *ig_quadf_gnode_glonum, int *ig_hexa_gll_glonum, int *ig_quadp_gll_glonum, int *ig_quadf_gll_glonum, int *ig_quadp_neighbor_hexa, int *ig_quadp_neighbor_hexaface, int *ig_hexa_material_number)
 
void fill_efispec_arrays (info_t *info, mesh_t *mesh, int rank, int *ig_ngll_total, int *ig_nneighbor_all_kind, int *cpu_neighbor, int *ig_cpu_neighbor_info, int *ig_hexa_gnode_glonum, int *ig_quadp_gnode_glonum, int *ig_quadf_gnode_glonum, int *ig_hexa_gll_glonum, int *ig_quadp_gll_glonum, int *ig_quadf_gll_glonum, int *ig_quadp_neighbor_hexa, int *ig_quadp_neighbor_hexaface, int *ig_hexa_material_number)
 
info_t * getConnInfo (mesh_t *mesh, int rank, int *ig_mesh_nnode)
 
int getProcConn (int ielt, int iproc, proc_info_t *proc_tab, mesh_t *mesh)
 
void setEltConn (int ielt, int iproc, elt_info_t *elt_tab, mesh_t *mesh, proc_info_t *proc_tab)
 
void setHexaConnDetail (int elt_num, int neigh_num, int nb_conn, int *connex_nodes, elt_info_t *elt_tab, mesh_t *mesh)
 
void setQuadConnDetail (int elt_num, int neigh_num, int nb_conn, int *connex_nodes, elt_info_t *elt_tab, mesh_t *mesh)
 
int findFaceOrientation (int *connex_nodes, int num_conn_source, int num_conn_target)
 
int findEdgeOrientation (int *connex_nodes, int num_conn_source, int num_conn_target)
 
char * ltrim (char *s)
 
char * rtrim (char *s)
 
char * trim (char *s)
 
void printMemUsage (mesh_t *mesh)
 
char * rotateVect (char *ortn, char *ret_or, int cw_quarter_rot)
 
void add_conn (topo_t typecon, elt_info_t *elt_tab, int hexa_src, int numcon_src, idx_t num_neigh, int numcon_neigh, int orientation, elt_t type)
 
int getnbneigh (topo_t typecon, elt_info_t *elt_tab, int numcon_src)
 
void free_all_1 ()
 
void free_all_2 ()
 
void get_domain_min_max (mesh_t *mesh, float *rg_mesh_xmax, float *rg_mesh_ymax, float *rg_mesh_zmax, float *rg_mesh_xmin, float *rg_mesh_ymin, float *rg_mesh_zmin)
 
void mod_init_mesh_mp_init_gll_number_ (int *, int *)
 
void mod_init_mesh_mp_propagate_gll_nodes_face_ (int *, int *, int *, int *, int *)
 
void mod_init_mesh_mp_propagate_gll_nodes_edge_ (int *, int *, int *, int *, int *)
 
void mod_init_mesh_mp_propagate_gll_nodes_corner_ (int *, int *, int *, int *)
 
void mod_init_mesh_mp_propagate_gll_nodes_quad_ (int *, int *, int *, int *, int *)
 
void pass_address_back_ (float *, float *, float *)
 
void allocate_node_coord_arrays_ ()
 
void flush ()
 

Variables

float * xcoord
 
float * ycoord
 
float * zcoord
 

Function Documentation

void fill_efispec_arrays ( info_t *  info,
mesh_t *  mesh,
int  rank,
int *  ig_ngll_total,
int *  ig_nneighbor_all_kind,
int *  cpu_neighbor,
int *  ig_cpu_neighbor_info,
int *  ig_hexa_gnode_glonum,
int *  ig_quadp_gnode_glonum,
int *  ig_quadf_gnode_glonum,
int *  ig_hexa_gll_glonum,
int *  ig_quadp_gll_glonum,
int *  ig_quadf_gll_glonum,
int *  ig_quadp_neighbor_hexa,
int *  ig_quadp_neighbor_hexaface,
int *  ig_hexa_material_number 
)

!! attention iparax vaut +1 ŕ partir d'ici

!! attention ifsurf vaut +1 ŕ partir d'ici

Definition at line 606 of file partCubitMesh.c.

622 {
623  int iproc, ineigh, ielt, i, j, k;
624  int* i_buf, *ibuf_hexa_material_number, *ibuf_hexa_gnode_glonum;
625  proc_info_t* proc_tab;
626  elt_info_t* elt_tab;
627  int* nneighbor_all_kind;
628  int* cpu_neighbor_info;
629  int* ibuf_quadp_neighbor_hexa;
630  int* ibuf_quadp_neighbor_hexaface;
631  int ibuf2[2];
632  MPI_Status status;
633  int icpu, iparax, ifsurf;
634 
635  idx_t *p_eptr, *p_eind;
636  if (mesh->hex27) {
637  p_eptr = mesh->eptr_27;
638  p_eind = mesh->eind_27;
639  } else {
640  p_eptr = mesh->eptr;
641  p_eind = mesh->eind;
642  }
643 
644  if (rank == 0) {
645  proc_tab = info->proc;
646  elt_tab = info->elt;
647 
648  for(iproc = 0; iproc<mesh->npart; iproc++) {
649  if (iproc == 0) {
650  icpu=0;//nb connected proc
651  for(ineigh = 0; ineigh<mesh->npart; ineigh++) {
652  if (proc_tab[iproc].connex[ineigh]) cpu_neighbor[icpu++] = ineigh; // list of connected proc
653  }
654  } else {
655  // envoi
656  i_buf = MALLOC(int, proc_tab[iproc].nb_conn);
657  icpu=0;//nb connected proc
658  for(ineigh = 0; ineigh<mesh->npart; ineigh++) {
659  if (proc_tab[iproc].connex[ineigh]) i_buf[icpu++] = ineigh; // list of connected proc
660  }
661  MPI_Send(i_buf, proc_tab[iproc].nb_conn, MPI_INT, iproc, 10, MPI_COMM_WORLD);
662  free(i_buf);
663  }
664  }
665  } else {
666  MPI_Recv(cpu_neighbor, ncpu_neighbor, MPI_INT, 0, 10, MPI_COMM_WORLD, &status);
667  }
668 
669  for(iproc = 0; iproc<mesh->npart; iproc++) {
670  int nneighbor_all_kind = 0;
671  int nb_items = 0;
672  if (rank == 0) {
673  if (iproc != 0) {
674  ibuf_hexa_material_number = MALLOC(int, proc_tab[iproc].nb_elt);
675  ibuf_hexa_gnode_glonum = MALLOC(int, proc_tab[iproc].nb_elt*mesh->num_nodes_hexa);
676 
677  ibuf_quadp_neighbor_hexa = MALLOC(int, proc_tab[iproc].nb_quad_p);
678  ibuf_quadp_neighbor_hexaface = MALLOC(int, proc_tab[iproc].nb_quad_p);
679 
680  i_buf = MALLOC(int, 7*26*info->proc[iproc].nb_elt); // 1 type + 5 params
681  cpu_neighbor_info = MALLOC(int, 3*26*info->proc[iproc].nb_ext);
682  }
683 
684  iparax = 0;
685  ifsurf = 0;
686 
687  for (ielt=0; ielt<proc_tab[iproc].nb_elt; ielt++) {// for each hexa of proc iproc
688 
689  int ihexa = ielt+1;
690  int ielt_num, ielt_number, ielt_face, ielt_edge, ielt_corner, ielt_coty, ielt_cpu;
691 
692  elt_info_t* curhex = proc_tab[iproc].local_elts[ielt];
693  int globnum = curhex->globalnum;
694 
695  // materiau
696  if (iproc != 0) ibuf_hexa_material_number[ielt] = mesh->layer[globnum];
697  else ig_hexa_material_number[ielt] = mesh->layer[globnum];
698 
699  // geometry
700  for (i=0; i<mesh->num_nodes_hexa; i++) {// geom coords of hexa geom nodes
701  // use efispec node order
702  // ig_hexa_gnode_glonum(inode,ihexa)
703  if (p_eind[p_eptr[globnum]+corner2efispec[i]-1] == 0) {
704  STOP("num node should not be 0 (fortran arrays start at 1)");
705  }
706  if (iproc != 0) ibuf_hexa_gnode_glonum[ielt * mesh->num_nodes_hexa + i] = p_eind[p_eptr[globnum]+corner2efispec[i]-1]; // numerotation starts from 1
707  else ig_hexa_gnode_glonum[ielt * mesh->num_nodes_hexa + i] = p_eind[p_eptr[globnum]+corner2efispec[i]-1]; // numerotation starts from 1
708  }
709  if (iproc == 0) mod_init_mesh_mp_init_gll_number_(&ihexa,ig_ngll_total);
710 
711  // faces
712  for (i=0; i<NFACE; i++) {
713  int iface = i+1;
714  // 6 faces connectivity
715  int type = curhex->faces[i].type;
716  switch (type) {
717  case HEXA : // hexa
718  ielt_num = elt_tab[curhex->faces[i].num_neigh].localnum + 1; // num commence ŕ 1
719  ielt_face = face2efispec[curhex->faces[i].num_conn];
720  ielt_coty = curhex->faces[i].orientation;
721  ielt_cpu = mesh->part[curhex->faces[i].num_neigh];
722  if (ielt_cpu == iproc) {
723  if (iproc == 0) mod_init_mesh_mp_propagate_gll_nodes_face_(&ihexa, &iface, &ielt_num, &ielt_face, &ielt_coty);
724  else {
725  i_buf[nb_items*7] = FACE*10+HEXA;
726  i_buf[nb_items*7+1] = ihexa;
727  i_buf[nb_items*7+2] = iface;
728  i_buf[nb_items*7+3] = ielt_num;
729  i_buf[nb_items*7+4] = ielt_face;
730  i_buf[nb_items*7+5] = ielt_coty;
731  i_buf[nb_items*7+6] = FACE;
732  nb_items++;
733  }
734  } else {
735  if (iproc == 0) {
736  if (*ig_nneighbor_all_kind >= 26*info->proc[iproc].nb_ext) STOP("ig_nneighbor_all_kind too large\n");
737  ig_cpu_neighbor_info[*ig_nneighbor_all_kind*3 + 0] = ihexa;
738  ig_cpu_neighbor_info[*ig_nneighbor_all_kind*3 + 1] = iface;
739  ig_cpu_neighbor_info[*ig_nneighbor_all_kind*3 + 2] = ielt_cpu;
740  (*ig_nneighbor_all_kind)++;
741  } else {
742  if (nneighbor_all_kind >= 26*info->proc[iproc].nb_ext) STOP("nneighbor_all_kind too large\n");
743  cpu_neighbor_info[nneighbor_all_kind*3 + 0] = ihexa;
744  cpu_neighbor_info[nneighbor_all_kind*3 + 1] = iface;
745  cpu_neighbor_info[nneighbor_all_kind*3 + 2] = ielt_cpu;
746  nneighbor_all_kind++;
747  }
748  }
749  break;
750  case QUAD_P : // quad p
751  if (iproc == 0) {
752  ig_quadp_neighbor_hexa[iparax] = ihexa;
753  ig_quadp_neighbor_hexaface[iparax] = iface;
754  } else {
755  ibuf_quadp_neighbor_hexa[iparax] = ihexa;
756  ibuf_quadp_neighbor_hexaface[iparax] = iface;
757  }
758  if (iproc == 0) {
759  if (mesh->num_node_per_dim*mesh->num_node_per_dim == 4) {
760  switch(iface) {
761  case 1 : ig_quadp_gnode_glonum[iparax*4+0] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+0];
762  ig_quadp_gnode_glonum[iparax*4+1] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+1];
763  ig_quadp_gnode_glonum[iparax*4+2] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+2];
764  ig_quadp_gnode_glonum[iparax*4+3] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+3];
765  break;
766  case 2 : ig_quadp_gnode_glonum[iparax*4+0] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+0];
767  ig_quadp_gnode_glonum[iparax*4+1] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+4];
768  ig_quadp_gnode_glonum[iparax*4+2] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+5];
769  ig_quadp_gnode_glonum[iparax*4+3] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+1];
770  break;
771  case 3 : ig_quadp_gnode_glonum[iparax*4+0] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+1];
772  ig_quadp_gnode_glonum[iparax*4+1] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+5];
773  ig_quadp_gnode_glonum[iparax*4+2] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+6];
774  ig_quadp_gnode_glonum[iparax*4+3] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+2];
775  break;
776  case 4 : ig_quadp_gnode_glonum[iparax*4+0] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+3];
777  ig_quadp_gnode_glonum[iparax*4+1] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+2];
778  ig_quadp_gnode_glonum[iparax*4+2] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+6];
779  ig_quadp_gnode_glonum[iparax*4+3] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+7];
780  break;
781  case 5 : ig_quadp_gnode_glonum[iparax*4+0] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+0];
782  ig_quadp_gnode_glonum[iparax*4+1] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+3];
783  ig_quadp_gnode_glonum[iparax*4+2] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+7];
784  ig_quadp_gnode_glonum[iparax*4+3] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+4];
785  break;
786  case 6 : ig_quadp_gnode_glonum[iparax*4+0] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+4];
787  ig_quadp_gnode_glonum[iparax*4+1] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+7];
788  ig_quadp_gnode_glonum[iparax*4+2] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+6];
789  ig_quadp_gnode_glonum[iparax*4+3] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+5];
790  break;
791  }
792  } else {
793  STOP("not yet ready for 27 nodes elements\n");
794  }
795  }
796  iparax++;
798  if (iproc == 0) mod_init_mesh_mp_propagate_gll_nodes_quad_(&ihexa, &iface, &iparax, ig_quadp_gll_glonum, &info->proc[iproc].nb_quad_p);
799  else {
800  i_buf[nb_items*7] = FACE*10+QUAD_P;
801  i_buf[nb_items*7+1] = ihexa;
802  i_buf[nb_items*7+2] = iface;
803  i_buf[nb_items*7+3] = iparax;
804  i_buf[nb_items*7+4] = 0;
805  i_buf[nb_items*7+5] = info->proc[iproc].nb_quad_p;
806  i_buf[nb_items*7+6] = FACE;
807  nb_items++;
808  // ŕ la construction de ig_quadp_gnode_glonum pour les autres proc, il faudra prendre ihexa-1 et iparax-1 !!!
809  }
810  break;
811  case QUAD_F : // quad f
812  if (iproc == 0) {
813  if (mesh->num_node_per_dim*mesh->num_node_per_dim == 4) {
814  switch(iface) {
815  case 1 : ig_quadf_gnode_glonum[ifsurf*4+0] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+0];
816  ig_quadf_gnode_glonum[ifsurf*4+1] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+1];
817  ig_quadf_gnode_glonum[ifsurf*4+2] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+2];
818  ig_quadf_gnode_glonum[ifsurf*4+3] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+3];
819  break;
820  case 2 : ig_quadf_gnode_glonum[ifsurf*4+0] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+0];
821  ig_quadf_gnode_glonum[ifsurf*4+1] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+4];
822  ig_quadf_gnode_glonum[ifsurf*4+2] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+5];
823  ig_quadf_gnode_glonum[ifsurf*4+3] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+1];
824  break;
825  case 3 : ig_quadf_gnode_glonum[ifsurf*4+0] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+1];
826  ig_quadf_gnode_glonum[ifsurf*4+1] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+5];
827  ig_quadf_gnode_glonum[ifsurf*4+2] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+6];
828  ig_quadf_gnode_glonum[ifsurf*4+3] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+2];
829  break;
830  case 4 : ig_quadf_gnode_glonum[ifsurf*4+0] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+3];
831  ig_quadf_gnode_glonum[ifsurf*4+1] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+2];
832  ig_quadf_gnode_glonum[ifsurf*4+2] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+6];
833  ig_quadf_gnode_glonum[ifsurf*4+3] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+7];
834  break;
835  case 5 : ig_quadf_gnode_glonum[ifsurf*4+0] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+0];
836  ig_quadf_gnode_glonum[ifsurf*4+1] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+3];
837  ig_quadf_gnode_glonum[ifsurf*4+2] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+7];
838  ig_quadf_gnode_glonum[ifsurf*4+3] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+4];
839  break;
840  case 6 : ig_quadf_gnode_glonum[ifsurf*4+0] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+4];
841  ig_quadf_gnode_glonum[ifsurf*4+1] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+7];
842  ig_quadf_gnode_glonum[ifsurf*4+2] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+6];
843  ig_quadf_gnode_glonum[ifsurf*4+3] = ig_hexa_gnode_glonum[ielt*mesh->num_nodes_hexa+5];
844  break;
845  }
846  } else {
847  STOP("not yet ready for 27 nodes elements\n");
848  }
849  }
850  ifsurf++;
852  if (iproc == 0) mod_init_mesh_mp_propagate_gll_nodes_quad_(&ihexa, &iface, &ifsurf, ig_quadf_gll_glonum, &info->proc->nb_quad_f);
853  else {
854  i_buf[nb_items*7] = FACE*10+QUAD_F;
855  i_buf[nb_items*7+1] = ihexa;
856  i_buf[nb_items*7+2] = iface;
857  i_buf[nb_items*7+3] = ifsurf;
858  i_buf[nb_items*7+4] = 0;
859  i_buf[nb_items*7+5] = info->proc[iproc].nb_quad_f;
860  i_buf[nb_items*7+6] = FACE;
861  nb_items++;
862  // ŕ la construction de ig_quadf_gnode_glonum pour les autres proc, il faudra prendre ihexa-1 et ifsurf-1 !!!
863  }
864  break;
865  }
866  }
867  // edges
868  for (i=0; i<NEDGE; i++) { // 12 edges connectivity
869  int iedge = i+1;
870  conn_info* ptr = &curhex->edges[i];
871  int nb_neighbor = getnbneigh(EDGE, curhex, i);
872  for (int ineigh = 0; ineigh < nb_neighbor; ineigh++) {
873  switch(ptr->type) {
874  case HEXA : ielt_number = elt_tab[ptr->num_neigh].localnum + 1; // start from 1
875  ielt_edge = edge2efispec[ptr->num_conn]; // num edge of neighbor
876  ielt_coty = ptr->orientation; // orientation
877  ielt_cpu = mesh->part[ptr->num_neigh]; // proc of neighbour
878  if (ielt_cpu == iproc) {
879  if (ielt_number > ihexa) {
880  if (iproc == 0 ) {
881  mod_init_mesh_mp_propagate_gll_nodes_edge_(&ihexa, &iedge, &ielt_number, &ielt_edge, &ielt_coty);
882  } else {
883  i_buf[nb_items*7] = EDGE*10+HEXA;
884  i_buf[nb_items*7+1] = ihexa;
885  i_buf[nb_items*7+2] = iedge;
886  i_buf[nb_items*7+3] = ielt_number;
887  i_buf[nb_items*7+4] = ielt_edge;
888  i_buf[nb_items*7+5] = ielt_coty;
889  i_buf[nb_items*7+6] = nb_neighbor;
890  nb_items++;
891  }
892  }
893  } else {
894  if (iproc == 0) {
895  ig_cpu_neighbor_info[*ig_nneighbor_all_kind*3 + 0] = ihexa;
896  ig_cpu_neighbor_info[*ig_nneighbor_all_kind*3 + 1] = NFACE + iedge;
897  ig_cpu_neighbor_info[*ig_nneighbor_all_kind*3 + 2] = ielt_cpu;
898  (*ig_nneighbor_all_kind)++;
899  } else {
900  cpu_neighbor_info[nneighbor_all_kind*3 + 0] = ihexa;
901  cpu_neighbor_info[nneighbor_all_kind*3 + 1] = NFACE + iedge;
902  cpu_neighbor_info[nneighbor_all_kind*3 + 2] = ielt_cpu;
903  nneighbor_all_kind++;
904  }
905  }
906  case NONE : break;
907  default : STOP("edge neighbor should be an hexa");
908  }
909  ptr = ptr->next;
910  }
911  }
912  // corner
913  for (j=0; j<NCORNER; j++) { // 8 corners connectivity
914  int inode = j+1;
915  i = cornerefispec2cubit[j]; // written in efispec order
916  conn_info* ptr = &curhex->corners[i];
917  int nb_neighbor = getnbneigh(CORNER, curhex, i); // nb of neighbors connected to corner i
918  for (k=0; k<nb_neighbor; k++) {
919  switch(ptr->type) {
920  case HEXA : ielt_number = elt_tab[ptr->num_neigh].localnum + 1; // start from 1
921  ielt_corner = corner2efispec[ptr->num_conn]; // num edge of neighbor
922  ielt_cpu = mesh->part[ptr->num_neigh]; // proc of neighbour
923  if (ielt_cpu == iproc) {
924  if (ielt_number > ihexa) {
925  if (iproc == 0 ) {
926  mod_init_mesh_mp_propagate_gll_nodes_corner_(&ihexa, &inode, &ielt_number, &ielt_corner);
927  } else {
928  i_buf[nb_items*7] = CORNER*10+HEXA;
929  i_buf[nb_items*7+1] = ihexa;
930  i_buf[nb_items*7+2] = inode;
931  i_buf[nb_items*7+3] = ielt_number;
932  i_buf[nb_items*7+4] = ielt_corner;
933  i_buf[nb_items*7+5] = 0;
934  i_buf[nb_items*7+6] = nb_neighbor;
935  nb_items++;
936  }
937  }
938  } else {
939  if (iproc == 0) {
940  ig_cpu_neighbor_info[*ig_nneighbor_all_kind*3 + 0] = ihexa;
941  ig_cpu_neighbor_info[*ig_nneighbor_all_kind*3 + 1] = NFACE + NEDGE + inode;
942  ig_cpu_neighbor_info[*ig_nneighbor_all_kind*3 + 2] = ielt_cpu;
943  (*ig_nneighbor_all_kind)++;
944  } else {
945  cpu_neighbor_info[nneighbor_all_kind*3 + 0] = ihexa;
946  cpu_neighbor_info[nneighbor_all_kind*3 + 1] = NFACE + NEDGE + inode;
947  cpu_neighbor_info[nneighbor_all_kind*3 + 2] = ielt_cpu;
948  nneighbor_all_kind++;
949  }
950  }
951  case NONE : break;
952  default : STOP("corner neighbor should be an hexa");
953  }
954  ptr = ptr->next;
955  }
956  }
957 
958  }
959  if (iparax != info->proc[iproc].nb_quad_p || ifsurf != info->proc[iproc].nb_quad_f) STOP("problem with quad number");
960 
961  if (iproc != 0) {
962  MPI_Send(ibuf_hexa_material_number, proc_tab[iproc].nb_elt, MPI_INT, iproc, 21, MPI_COMM_WORLD);
963  MPI_Send(ibuf_hexa_gnode_glonum, proc_tab[iproc].nb_elt*mesh->num_nodes_hexa, MPI_INT, iproc, 22, MPI_COMM_WORLD);
964  MPI_Send(ibuf_quadp_neighbor_hexa, proc_tab[iproc].nb_quad_p, MPI_INT, iproc, 23, MPI_COMM_WORLD);
965  MPI_Send(ibuf_quadp_neighbor_hexaface, proc_tab[iproc].nb_quad_p, MPI_INT, iproc, 24, MPI_COMM_WORLD);
966 
967  ibuf2[0] = nb_items;
968  ibuf2[1] = nneighbor_all_kind;
969  MPI_Send(ibuf2, 2, MPI_INT, iproc, 25, MPI_COMM_WORLD);
970 
971  MPI_Send(i_buf, nb_items*7, MPI_INT, iproc, 26, MPI_COMM_WORLD);
972  MPI_Send(cpu_neighbor_info, 3*nneighbor_all_kind, MPI_INT, iproc, 27, MPI_COMM_WORLD);
973 
974  free(ibuf_hexa_material_number);
975  free(ibuf_hexa_gnode_glonum);
976  free(ibuf_quadp_neighbor_hexa);
977  free(ibuf_quadp_neighbor_hexaface);
978  free(i_buf);
979  free(cpu_neighbor_info);
980  }
981  } else {
982  if (rank == iproc) {
983  i_buf = MALLOC(int, 6*26*nhexa_local);
984  MPI_Recv(ig_hexa_material_number, nhexa_local, MPI_INT, 0, 21, MPI_COMM_WORLD, &status);
985  MPI_Recv(ig_hexa_gnode_glonum, nhexa_local*num_nodes_hexa, MPI_INT, 0, 22, MPI_COMM_WORLD, &status);
986  MPI_Recv(ig_quadp_neighbor_hexa, nquad_parax, MPI_INT, 0, 23, MPI_COMM_WORLD, &status);
987  MPI_Recv(ig_quadp_neighbor_hexaface, nquad_parax, MPI_INT, 0, 24, MPI_COMM_WORLD, &status);
988 
989  MPI_Recv(ibuf2, 2, MPI_INT, 0, 25, MPI_COMM_WORLD, &status);
990  nb_items = ibuf2[0];
991  *ig_nneighbor_all_kind = nneighbor_all_kind = ibuf2[1];
992 
993  MPI_Recv(i_buf, nb_items*7, MPI_INT, 0, 26, MPI_COMM_WORLD, &status);
994  MPI_Recv(ig_cpu_neighbor_info, 3*nneighbor_all_kind, MPI_INT, 0, 27, MPI_COMM_WORLD, &status);
995  int iibuf = 0;
996  for (ielt=0; ielt < nhexa_local; ielt++) {
997 
998  int ihexa = ielt+1;
999  mod_init_mesh_mp_init_gll_number_(&ihexa,ig_ngll_total);
1000  // besoin de savoir combien d'items -> tests avec ihexa & iface
1001  for (i=0; i<NFACE; i++) { // a revérifier ...
1002  int iface = i+1;
1003  if (i_buf[iibuf*7+6] == FACE && i_buf[iibuf*7+1] == ihexa && i_buf[iibuf*7+2] == iface) {
1004  int type = i_buf[iibuf*7];
1005  switch (type) {
1006  case FACE*10+HEXA : mod_init_mesh_mp_propagate_gll_nodes_face_(&i_buf[iibuf*7+1], &i_buf[iibuf*7+2], &i_buf[iibuf*7+3], &i_buf[iibuf*7+4], &i_buf[iibuf*7+5]);
1007  iibuf++;
1008  break;
1009  case FACE*10+QUAD_P : iparax = i_buf[iibuf*7+3]-1;
1010  assert(i_buf[iibuf*7+1]-1 == ielt);
1011  assert(i_buf[iibuf*7+2] == iface);
1012  switch(iface) {
1013  case 1 : ig_quadp_gnode_glonum[iparax*4+0] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+0];
1014  ig_quadp_gnode_glonum[iparax*4+1] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+1];
1015  ig_quadp_gnode_glonum[iparax*4+2] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+2];
1016  ig_quadp_gnode_glonum[iparax*4+3] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+3];
1017  break;
1018  case 2 : ig_quadp_gnode_glonum[iparax*4+0] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+0];
1019  ig_quadp_gnode_glonum[iparax*4+1] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+4];
1020  ig_quadp_gnode_glonum[iparax*4+2] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+5];
1021  ig_quadp_gnode_glonum[iparax*4+3] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+1];
1022  break;
1023  case 3 : ig_quadp_gnode_glonum[iparax*4+0] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+1];
1024  ig_quadp_gnode_glonum[iparax*4+1] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+5];
1025  ig_quadp_gnode_glonum[iparax*4+2] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+6];
1026  ig_quadp_gnode_glonum[iparax*4+3] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+2];
1027  break;
1028  case 4 : ig_quadp_gnode_glonum[iparax*4+0] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+3];
1029  ig_quadp_gnode_glonum[iparax*4+1] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+2];
1030  ig_quadp_gnode_glonum[iparax*4+2] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+6];
1031  ig_quadp_gnode_glonum[iparax*4+3] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+7];
1032  break;
1033  case 5 : ig_quadp_gnode_glonum[iparax*4+0] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+0];
1034  ig_quadp_gnode_glonum[iparax*4+1] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+3];
1035  ig_quadp_gnode_glonum[iparax*4+2] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+7];
1036  ig_quadp_gnode_glonum[iparax*4+3] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+4];
1037  break;
1038  case 6 : ig_quadp_gnode_glonum[iparax*4+0] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+4];
1039  ig_quadp_gnode_glonum[iparax*4+1] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+7];
1040  ig_quadp_gnode_glonum[iparax*4+2] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+6];
1041  ig_quadp_gnode_glonum[iparax*4+3] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+5];
1042  break;
1043  }
1044  mod_init_mesh_mp_propagate_gll_nodes_quad_(&i_buf[iibuf*7+1], &i_buf[iibuf*7+2], &i_buf[iibuf*7+3], ig_quadp_gll_glonum, &i_buf[iibuf*7+5]);
1045  iibuf++;
1046  break;
1047  case FACE*10+QUAD_F : ifsurf = i_buf[iibuf*7+3]-1;
1048  assert(i_buf[iibuf*7+1]-1 == ielt);
1049  assert(i_buf[iibuf*7+2] == iface);
1050  switch(iface) {
1051  case 1 : ig_quadf_gnode_glonum[ifsurf*4+0] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+0];
1052  ig_quadf_gnode_glonum[ifsurf*4+1] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+1];
1053  ig_quadf_gnode_glonum[ifsurf*4+2] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+2];
1054  ig_quadf_gnode_glonum[ifsurf*4+3] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+3];
1055  break;
1056  case 2 : ig_quadf_gnode_glonum[ifsurf*4+0] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+0];
1057  ig_quadf_gnode_glonum[ifsurf*4+1] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+4];
1058  ig_quadf_gnode_glonum[ifsurf*4+2] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+5];
1059  ig_quadf_gnode_glonum[ifsurf*4+3] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+1];
1060  break;
1061  case 3 : ig_quadf_gnode_glonum[ifsurf*4+0] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+1];
1062  ig_quadf_gnode_glonum[ifsurf*4+1] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+5];
1063  ig_quadf_gnode_glonum[ifsurf*4+2] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+6];
1064  ig_quadf_gnode_glonum[ifsurf*4+3] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+2];
1065  break;
1066  case 4 : ig_quadf_gnode_glonum[ifsurf*4+0] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+3];
1067  ig_quadf_gnode_glonum[ifsurf*4+1] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+2];
1068  ig_quadf_gnode_glonum[ifsurf*4+2] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+6];
1069  ig_quadf_gnode_glonum[ifsurf*4+3] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+7];
1070  break;
1071  case 5 : ig_quadf_gnode_glonum[ifsurf*4+0] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+0];
1072  ig_quadf_gnode_glonum[ifsurf*4+1] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+3];
1073  ig_quadf_gnode_glonum[ifsurf*4+2] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+7];
1074  ig_quadf_gnode_glonum[ifsurf*4+3] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+4];
1075  break;
1076  case 6 : ig_quadf_gnode_glonum[ifsurf*4+0] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+4];
1077  ig_quadf_gnode_glonum[ifsurf*4+1] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+7];
1078  ig_quadf_gnode_glonum[ifsurf*4+2] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+6];
1079  ig_quadf_gnode_glonum[ifsurf*4+3] = ig_hexa_gnode_glonum[ielt*num_nodes_hexa+5];
1080  break;
1081  }
1082  mod_init_mesh_mp_propagate_gll_nodes_quad_(&i_buf[iibuf*7+1], &i_buf[iibuf*7+2], &i_buf[iibuf*7+3], ig_quadf_gll_glonum, &i_buf[iibuf*7+5]);
1083  iibuf++;
1084  break;
1085  }
1086  }
1087  }
1088  for (i=0; i<NEDGE; i++) {
1089  int iedge = i+1;
1090  if (i_buf[iibuf*7] == EDGE*10+HEXA && i_buf[iibuf*7+1] == ihexa && i_buf[iibuf*7+2] == iedge) {
1091  int nb_neighbor = i_buf[iibuf*7+6];
1092  for (int ineigh = 0; ineigh < nb_neighbor; ineigh++) {
1093  if (i_buf[iibuf*7] == EDGE*10+HEXA && i_buf[iibuf*7+1] == ihexa && i_buf[iibuf*7+2] == iedge) {
1094  mod_init_mesh_mp_propagate_gll_nodes_edge_(&i_buf[iibuf*7+1], &i_buf[iibuf*7+2], &i_buf[iibuf*7+3], &i_buf[iibuf*7+4], &i_buf[iibuf*7+5]);
1095  iibuf++;
1096  }
1097  }
1098  }
1099  }
1100  for (j=0; j<NCORNER; j++) {
1101  int inode = j+1;
1102  if (i_buf[iibuf*7] == CORNER*10+HEXA && i_buf[iibuf*7+1] == ihexa && i_buf[iibuf*7+2] == inode) {
1103  int nb_neighbor = i_buf[iibuf*7+6];
1104  for (int ineigh = 0; ineigh < nb_neighbor; ineigh++) {
1105  if (i_buf[iibuf*7] == CORNER*10+HEXA && i_buf[iibuf*7+1] == ihexa && i_buf[iibuf*7+2] == inode) {
1106  mod_init_mesh_mp_propagate_gll_nodes_corner_(&i_buf[iibuf*7+1], &i_buf[iibuf*7+2], &i_buf[iibuf*7+3], &i_buf[iibuf*7+4]);
1107  iibuf++;
1108  }
1109  }
1110  }
1111  }
1112  }
1113  free(i_buf);
1114  assert(iibuf == nb_items);
1115  }
1116  }
1117  } // endfor iproc
1118 
1119  return;
1120 }