OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
resol_init.F
Go to the documentation of this file.
1Copyright> OpenRadioss
2Copyright> Copyright (C) 1986-2025 Altair Engineering Inc.
3Copyright>
4Copyright> This program is free software: you can redistribute it and/or modify
5Copyright> it under the terms of the GNU Affero General Public License as published by
6Copyright> the Free Software Foundation, either version 3 of the License, or
7Copyright> (at your option) any later version.
8Copyright>
9Copyright> This program is distributed in the hope that it will be useful,
10Copyright> but WITHOUT ANY WARRANTY; without even the implied warranty of
11Copyright> MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12Copyright> GNU Affero General Public License for more details.
13Copyright>
14Copyright> You should have received a copy of the GNU Affero General Public License
15Copyright> along with this program. If not, see <https://www.gnu.org/licenses/>.
16Copyright>
17Copyright>
18Copyright> Commercial Alternative: Altair Radioss Software
19Copyright>
20Copyright> As an alternative to this open-source version, Altair also offers Altair Radioss
21Copyright> software under a commercial license. Contact Altair to discuss further if the
22Copyright> commercial version may interest you: https://www.altair.com/radioss/.
23!||====================================================================
24!|| resol_init ../engine/source/engine/resol_init.F
25!||--- called by ------------------------------------------------------
26!|| resol ../engine/source/engine/resol.F
27!||--- calls -----------------------------------------------------
28!|| admgvid ../engine/source/model/remesh/admgvid.F
29!|| admini ../engine/source/model/remesh/admini.F
30!|| admordr ../engine/source/model/remesh/admordr.F
31!|| anim_xfe_init ../engine/source/output/anim/generate/anim_crk_init.F
32!|| assadd2 ../engine/source/assembly/assadd2.f
33!|| chkinit ../engine/source/interfaces/interf/chkstfn3.F
34!|| cndmasi2_dim ../engine/source/elements/solid/solide10/s10cndf.F
35!|| cndmasi2_ini ../engine/source/elements/solid/solide10/s10cndf.F
36!|| cndordr ../engine/source/model/remesh/cndordr.F
37!|| dim_tshedg ../engine/source/elements/thickshell/solidec/dim_tshedg.f
38!|| fillipartl ../engine/source/engine/resol_init.F
39!|| findgroupc ../engine/source/elements/findgroup.F
40!|| findgroups ../engine/source/elements/findgroup.F
41!|| grpsplit ../engine/source/engine/resol_init.F
42!|| imp_init ../engine/source/implicit/imp_init.F
43!|| ind_tshedg ../engine/source/elements/thickshell/solidec/ind_tshedg.F
44!|| ini_tmax ../engine/source/output/ini_outmax.F
45!|| init_kyne ../engine/source/engine/resol_init.f
46!|| init_reac_nod ../engine/source/output/th/init_reac_nod.F
47!|| init_th_group ../engine/source/output/th/init_th_group.F
48!|| initimeg ../engine/source/system/timer.f
49!|| inivel_init ../engine/source/loads/general/inivel/inivel_init.f90
50!|| int_flushtime ../common_source/modules/interfaces/metric_mod.F
51!|| inter_sh_offset_ini ../engine/source/interfaces/shell_offset/inter_offset_ini.F90
52!|| kinini ../engine/source/constraints/general/kinini.F
53!|| mpp_init ../engine/source/mpi/interfaces/spmd_i7tool.F
54!|| my_barrier ../engine/source/system/machine.F
55!|| nloc_count_solnod ../engine/source/elements/solid/solide/nloc_count_solnod.F90
56!|| r2r_init ../engine/source/coupling/rad2rad/r2r_init.F
57!|| rbe2_init ../engine/source/constraints/general/rbe2/rbe2f.F
58!|| s10cnd_ini ../engine/source/elements/solid/solide10/s10cndf.F
59!|| s10cndi2_ini ../engine/source/elements/solid/solide10/s10cndf.F
60!|| s10cnds_dim ../engine/source/elements/solid/solide10/s10cndf.F
61!|| s10cnds_ini ../engine/source/elements/solid/solide10/s10cndf.f
62!|| section_init ../engine/source/tools/sect/section_init.F
63!|| spmd_anim_ply_init ../engine/source/mpi/anim/spmd_anim_ply_init.F
64!|| spmd_failwave_boundaries ../engine/source/mpi/output/spmd_exch_failwave.F
65!|| spmd_max_i ../engine/source/mpi/implicit/imp_spmd.f
66!|| spmd_sub_boundaries ../engine/source/mpi/spmd_exch_sub.F
67!|| tmax_ipart ../engine/source/output/tmax_ipart.F
68!|| tshcdcom_dim ../engine/source/elements/thickshell/solidec/tshcdcom_dim.F
69!|| tshcdcom_ini ../engine/source/elements/thickshell/solidec/tshcdcom_ini.F
70!|| zero1 ../engine/source/system/zero.F
71!|| zeror ../engine/source/system/zero.F
72!||--- uses -----------------------------------------------------
73!|| alemuscl_mod ../common_source/modules/ale/alemuscl_mod.F
74!|| crackxfem_mod ../engine/share/modules/crackxfem_mod.F
75!|| dtdc_mod ../engine/share/modules/dtdc_mod.F
76!|| ecnd_mod ../engine/share/modules/ecdn_mod.F
77!|| elbufdef_mod ../common_source/modules/mat_elem/elbufdef_mod.F90
78!|| failwave_mod ../common_source/modules/failwave_mod.F
79!|| glob_therm_mod ../common_source/modules/mat_elem/glob_therm_mod.F90
80!|| groupdef_mod ../common_source/modules/groupdef_mod.F
81!|| h3d_inc_mod ../engine/share/modules/h3d_inc_mod.F
82!|| h3d_mod ../engine/share/modules/h3d_mod.F
83!|| inivel_init_mod ../engine/source/loads/general/inivel/inivel_init.F90
84!|| intbufdef_mod ../common_source/modules/interfaces/intbufdef_mod.F90
85!|| inter_sh_offset_ini_mod ../engine/source/interfaces/shell_offset/inter_offset_ini.F90
86!|| inter_sh_offset_mod ../engine/source/modules/interfaces/sh_offset_mod.F90
87!|| inter_sorting_mod ../engine/share/modules/inter_sorting_mod.F
88!|| loads_mod ../common_source/modules/loads/loads_mod.F90
89!|| nloc_count_solnod_mod ../engine/source/elements/solid/solide/nloc_count_solnod.F90
90!|| nlocal_reg_mod ../common_source/modules/nlocal_reg_mod.F
91!|| outmax_mod ../common_source/modules/outmax_mod.F
92!|| output_mod ../common_source/modules/output/output_mod.F90
93!|| outputs_mod ../common_source/modules/outputs_mod.F
94!|| parith_on_mod ../common_source/modules/parith_on_mod.F90
95!|| pblast_mod ../common_source/modules/loads/pblast_mod.F90
96!|| pinchtype_mod ../common_source/modules/pinchtype_mod.F
97!|| plyxfem_mod ../engine/share/modules/plyxfem_mod.F
98!|| rbe3_mod ../common_source/modules/constraints/rbe3_mod.F90
99!|| sensor_mod ../common_source/modules/sensor_mod.F90
100!|| spmd_xv_inter_type1_mod ../engine/source/mpi/nodes/spmd_sd_xv_inter1.F90
101!|| stack_mod ../engine/share/modules/stack_mod.F
102!||====================================================================
103 SUBROUTINE resol_init(
104 1 ITASK ,FR_NBCC ,
105 2 ISENDTO ,IRCVFROM ,IAD_ELEM,FR_ELEM ,ITABM1 ,
106 3 IPARI ,IPARG ,ITAB ,IXS10 ,IXS20 ,
107 4 I13A ,I13B ,I13C ,I13D ,I13E ,
108 5 I13F ,I13G ,I13H ,I13I ,I15A ,
109 6 I15B ,I15C ,I15D ,I15E ,I15F ,
110 7 I15G ,I15H ,I15I ,I87A ,I87B ,
111 8 I87C ,I87D ,I87E ,I87F ,I87G ,
112 9 NFIA ,NFEA ,NFOA ,NDMA ,NDMA2 ,
113 A NODFT ,NODLT ,NDTASK ,NUMNTHREAD ,IXS16 ,
114 B IXS ,IXQ ,IXC ,IXT ,IXP ,
115 C IXR ,IXTG ,PON, IKINE ,
116 D A ,AR ,V ,VR ,
117 E X ,D ,MS ,IN ,STIFN ,
118 F STIFR ,DMAS ,DINER ,FANI ,ANIN ,
119 G WA ,UWA ,PM ,GEO ,
120 H PARTSAV ,PARTS0 ,MONVOL ,
121 I I87H ,I87I ,I87J ,I87K ,
122 J I15J ,KXX ,
123 K SECBUF ,SECFCUM ,NSTRF ,IGRNOD ,IEXLNK ,
124 L XFRAME ,
125 M IXTG1 ,IB ,VISCN ,DD_R2R ,
126 O ELBUF ,IPART ,MADPRT ,MADSH4 ,
127 P MADSH3 ,MADSOL ,MADNOD ,MADFAIL ,IGEO ,
128 Q INTLIST ,NBINTC ,PROCNE ,NISKYFI ,WEIGHT ,
129 R ISIZXV ,ILENXV ,ADDCNI2 ,PROCNI2 ,IAD_I2M ,
130 S FR_I2M ,FR_NBCCI2,I2SIZE ,FR_MAD ,LWIBEM ,
131 T LWRBEM ,FXBFP ,FXBEFW ,FXBEDP ,FXBGRP ,
132 U FXBGRW ,NDIN ,
133 V ISLEN7 ,IRLEN7 ,ISLEN11 ,IRLEN11 ,
134 W LWIFLOW ,LWRFLOW ,IFLOW ,ADDCNEL ,CNEL ,
135 X ADDTMPL ,IPARTL ,NPARTL ,NFNCA ,NFTCA ,
136 Y I15ATH ,I35ATH ,IPM ,SH4TREE ,IPADMESH ,
137 Z MSC ,INC ,SH3TREE ,MSTG ,INTG ,
138 a PTG ,FTHE ,FTHESKY ,FTHESKYI ,NME17 ,
139 b ISLEN17 ,IRLEN17 ,IRLEN7T ,ISLEN7T ,LINDIDEL ,
140 c LBUFIDEL,SH4TRIM ,SH3TRIM ,MSCND ,INCND ,
141 d IRLEN20 ,ISLEN20 ,IRLEN20T,ISLEN20T ,NBINT20 ,
142 e IRLEN20E,ISLEN20E ,NISKYFIE,
143 f MCP ,MS0 ,INOD_PXFEM,IEL_PXFEM,IADC_PXFEM,
144 g ADSKY_PXFEM,ICODT,ICODR ,IBFV ,ADMSMS ,
145 h NODREAC ,IGROUC ,NGROUC ,IGROUNC ,NGROUNC ,
146 i FR_RBY ,FR_RBY6 ,NPBY ,
147 j NOM_SECT ,MCPC ,MCPTG ,GRTH ,IGRTH ,
148 k NELEM ,LAG_SEC ,NPRW ,DIAG_SMS ,DMELC ,
149 l DMELTG ,NGRTH ,NFT2 ,DMELS ,DMELTR ,
150 m DMELP ,DMELRT ,RES_SMS ,I87L ,IRBE2 ,
151 n LRBE2 ,NMRBE2 ,IAD_RBE2 ,FR_RBE2 ,FR_RBE2M ,
152 o R2SIZE ,LPBY ,PROCNE_PXFEM ,ISENDP_PXFEM,IRECVP_PXFEM,
153 p IADSDP_PXFEM,IADRCP_PXFEM,FR_NBCC1,RBY,INT18KINE ,
154 q XDP ,I87M,INOD_CRKXFEM,IEL_CRKXFEM ,IADC_CRKXFEM,
155 r ADSKY_CRKXFEM,PROCNE_CRKXFEM,ISENDP_CRKXFEM,IRECVP_CRKXFEM,
156 s IADSDP_CRKXFEM,IADRCP_CRKXFEM ,INT24USE,NDAMA2 ,
157 t IGROUPC ,IGROUPTG ,IGROUPS ,IGROUPFLG ,DMINT2 ,IRBKIN_L,
158 u NRBYKIN_L,KINDRBY ,ELBUF_TAB ,SENSORS ,DD_R2R_ELEM,
159 v SDD_R2R_ELEM ,KINET, WEIGHT_MD,DMSPH ,IOLDSECT,LBUFIDEL24,
160 w INTBUF_TAB ,NUMSPH_GLO_R2R, FLG_SPHINOUT_R2R,I15K,
161 y CONDN ,CONDNSKY ,KXFENOD2ELC ,ELCUTC ,NODEDGE,
162 z IAD_EDGE ,CRKNODIAD,FR_EDGE ,FR_NBEDGE ,NODLEVXF,
163 x CRKEDGE ,XFEM_TAB ,ISENSINT , NISUBMAX,
164 1 INTLIST25 ,INT24E2EUSE ,TABMP_L ,
165 2 I87N ,TAB_MAT,H3D_DATA,TAGTRIMC,TAGTRIMTG,
166 3 IGRBRIC ,IGRQUAD ,IGRSH4N ,IGRSH3N ,IGRTRUSS ,
167 4 IGRBEAM ,IGRSPRING,IGRPART ,FORNEQS ,INT7ITIED ,
168 5 FXVEL_FGEO,FAILWAVE,NLOC_DMG,PINCH_DATA,SLLOADP ,
169 6 TAGSLV_RBY,NFNCA2 ,NFTCA2 ,IN0 ,SORT_COMM,STACK,OUTPUT,
170 7 THKE ,SFR_ELEM ,SH_OFFSET_TAB,
171 8 NEED_COMM_INT25_SOLID_EROSION,COMM_INT25_SOLID_EROSION,
172 9 ISKWN ,IFRAME, LOADS ,GLOB_THERM,PBLAST,RBE3)
173C-----------------------------------------------
174C M o d u l e s
175C-----------------------------------------------
176 USE plyxfem_mod
177 USE elbufdef_mod
178 USE intbufdef_mod
179 USE crackxfem_mod
180 USE ecnd_mod
181 USE h3d_mod
182 USE groupdef_mod
183 USE failwave_mod
185 USE pinchtype_mod
186 USE pblast_mod
187 USE dtdc_mod
189 USE stack_mod
190 USE outmax_mod
191 USE sensor_mod
192 USE h3d_inc_mod
193 USE outputs_mod
195 USE output_mod
196 USE nloc_count_solnod_mod
197 USE inter_sh_offset_ini_mod , only : inter_sh_offset_ini
198 USE inter_sh_offset_mod , only:sh_offset_
199 USE loads_mod
200 USE inivel_init_mod , only: inivel_init
201 use glob_therm_mod
202 use spmd_xv_inter_type1_mod , only : is_present_inter1
203 USE parith_on_mod, only: element_pon_
204 use rbe3_mod
205C-----------------------------------------------
206C I m p l i c i t T y p e s
207C-----------------------------------------------
208#include "implicit_f.inc"
209C-----------------------------------------------
210C C o m m o n B l o c k s
211C-----------------------------------------------
212#include "com01_c.inc"
213#include "com04_c.inc"
214#include "com08_c.inc"
215#include "com10_c.inc"
216#include "com_xfem1.inc"
217#include "param_c.inc"
218#include "scr02_c.inc"
219#include "scr03_c.inc"
220#include "scr07_c.inc"
221#include "scr12_c.inc"
222#include "scr14_c.inc"
223#include "scr16_c.inc"
224#include "scr17_c.inc"
225#include "scr23_c.inc"
226#include "units_c.inc"
227#include "cong2_c.inc"
228#include "task_c.inc"
229#include "parit_c.inc"
230#include "timerc_c.inc"
231#include "rad2r_c.inc"
232#include "scr18_c.inc"
233#include "spmd_c.inc"
234#include "fxbcom.inc"
235#include "flowcom.inc"
236#include "remesh_c.inc"
237#include "sms_c.inc"
238#include "lagmult.inc"
239#include "sphcom.inc"
240#include "intstamp_c.inc"
241C-----------------------------------------------------------------
242C D u m m y A r g u m e n t s
243C-----------------------------------------------
244 TYPE(element_pon_) :: PON
245 INTEGER ITASK, NBINTC, NODFT, NODLT, LINDIDEL, LBUFIDEL,
246 . numnthread, ndtask, nfia, nfea, nfoa ,ndma, nfnca, nftca,
247 . ndma2,ndin,n1,n2,n3,igtyp,npartl,ngrouc,ngrounc,
248 . i13a,i13b,i13c,i13d,i13e,i13f,i13g,i13h,i13i,
249 . i15a,i15b,i15c,i15d,i15e,i15f,i15g,i15h,i15i,i15j,i15k,
250 . i87a,i87b,i87c,i87d,i87e,i87f,i87g,i87h,i87i,i87j,
251 . i87k,i87l,i87m,i87n,nfnca2,nftca2,
252 . isizxv , ilenxv, i2size, islen7,irlen7 ,islen11 ,irlen11,
253 . i15ath, i35ath, nme17,islen17,irlen17,irlen7t,islen7t,
254 . irlen20,islen20,irlen20t,islen20t,nbint20,irlen20e,
255 . islen20e,nelem,lag_sec, ngrth, nft2,nmrbe2,
256 . int18kine,int24use,ndama2, nrbykin_l,ioldsect,lbufidel24,
257 . tabmp_l,tagtrimc(*),tagtrimtg(*), slloadp,sfr_elem
258 INTEGER
259 . ixs(nixs,*),ixs10(6,*) ,ixs20(12,*),
260 . ixs16(6,*) , igeo(npropgi,*),
261 . ixq(nixq,*),ixc(nixc,*), ixt(nixt,*), ixp(nixp,*),
262 . ixr(nixr,*), ixtg(nixtg,*), ixtg1(4,*),
263 . itab(*), iparg(nparg,*), ipari(npari,*),
264 . iexlnk(nr2r,*),
265 . weight(*), nstrf(*), ib(nibcld,*), itabm1(*),
266 . monvol(*),kxx(nixx,*),isendto(ninter+1,nspmd+1),
267 . fr_nbcc(2,nspmd+1), iad_elem(2,nspmd+1) ,fr_elem(*),
268 . ircvfrom(ninter+1,nspmd+1), intlist(ninter), procne(*),
269 . niskyfi(*),addcni2(*),procni2(*),iad_i2m(*),fr_i2m(*),
270 . fr_nbcci2(*), ipart(*),
271 . dd_r2r(nspmd+1,*),ipartl(*),
272 . madprt(*), madsh4(*), madsh3(*), madsol(*), madnod(*),
273 . madfail(*), fr_mad(5,*), lwibem, lwrbem, lwiflow, lwrflow,
274 . iflow(*), addcnel(0:*), cnel(0:*), addtmpl(0:*),
275 . ipm(npropmi,*), sh4tree(*), ipadmesh(*), sh3tree(*),
276 . sh4trim(*), sh3trim(*), niskyfie(*),
277 . icodt(*), icodr(*),ibfv(nifv,*),
278 . inod_pxfem(*),iel_pxfem(*) ,iadc_pxfem(4,*),elcutc(2,*),
279 . adsky_pxfem(*), kxfenod2elc(*),nodlevxf(*),crknodiad(*),
280 . nodedge(*),iad_edge(*),fr_edge(*),fr_nbedge(*), nodreac(*),
281 . igrouc(*),igrounc(*),fr_rby(*),fr_rby6(*),npby(*),
282 . nom_sect(*), grth(*),igrth(*), nprw(*),iad_rbe2(*),
283 . fr_rbe2(*),fr_rbe2m(*),r2size, irbe2(nrbe2l,*),lrbe2(*),
284 . ikine(numnod),lpby(*), procne_pxfem(*),
285 . isendp_pxfem(*),irecvp_pxfem(*),iadsdp_pxfem(*),
286 . iadrcp_pxfem(*),fr_nbcc1(2,*),inod_crkxfem(*),
287 . iel_crkxfem(*),iadc_crkxfem(*),adsky_crkxfem(0:*),
288 . procne_crkxfem(*),isendp_crkxfem(*),irecvp_crkxfem(*),
289 . iadsdp_crkxfem(*),iadrcp_crkxfem(*),
290 . igroupc(*),igrouptg(*),igroups(*),igroupflg(2),
291 . irbkin_l(*), kindrby(*), dd_r2r_elem(*),sdd_r2r_elem,
292 . kinet(*),weight_md(*),numsph_glo_r2r,flg_sphinout_r2r,
293 . isensint(nisubmax+1,ninter),nisubmax,
294 . intlist25(ninter25) ,int24e2euse ,fxvel_fgeo,
295 . tagslv_rby(numnod)
296 INTEGER, INTENT(IN ),DIMENSION(LISKN,NUMFRAM+1) :: IFRAME
297 INTEGER, INTENT(IN ),DIMENSION(LISKN,NUMSKW+1) :: ISKWN
298! INT7ITIED : check if an interface type 7 with ITIED /= 0 is used
299! in order to force the communication of a list of candidate nodes
300! INT7ITIED = 0 type 7 + ITIED/=0 not used
301! INT7ITIED = 1 type 7 + ITIED/=0 used
302 INTEGER, INTENT(INOUT) :: INT7ITIED
303 my_real
304 . X(3,*), D(3,*), V(3,*), VR(3,*),
305 . MS(*), IN(*), WA(*), A(3,*), AR(3,*),
306 . FANI(3,*), UWA(*), STIFN(*), STIFR(*),
307 . ANIN(*), PARTSAV(NPSAV,*),PARTS0(*),
308 . DMAS, DINER ,
309 . PM(NPROPM,*) , GEO(NPROPG,*),
310 . VISCN(*),
311 . SECBUF(*),SECFCUM(7,NUMNOD,NSECT),XFRAME(NXFRAME,*),
312 . ELBUF(*), MSC(*), INC(*), MSTG(*), INTG(*), PTG(*),
313 . mscnd(*), incnd(*), fthe(*), fthesky(*),ftheskyi(*), mcp(*),
314 . ms0(*), admsms(*), mcpc(*), mcptg(*), diag_sms(*),
315 . dmelc(*), dmeltg(*), dmels(*), dmeltr(*), dmelp(*), dmelrt(*),
316 . res_sms(3,*),rby(nrby,*), dmint2(4,i2nsn25),
317 . dmsph(*),condn(*),condnsky ,tab_mat(ngroup),forneqs(3,*)
318 my_real
319 . fxbfp(*), fxbefw(*), fxbedp(*), fxbgrp(*), fxbgrw(*),in0(*)
320 my_real
321 . thke(numelc+numeltg)
322c INTEGER*8
323c . I8A(3,3,*),I8AR(3,3,*),I8STIFN(3,*),I8STIFR(3,*),
324c . I8VISCN(3,*)
325C
326
327 LOGICAL, INTENT(inout) :: NEED_COMM_INT25_SOLID_EROSION !< boolean, true if the proc needs to comm some values related to interface type 25 with solid erosion
328 INTEGER, INTENT(inout) :: COMM_INT25_SOLID_EROSION !< integer, sub-communicator related to interface type 25 with solid erosion
329C
330 DOUBLE PRECISION XDP(3,*)
331 TYPE(INTBUF_STRUCT_) INTBUF_TAB(*)
332 TYPE (ELBUF_STRUCT_), DIMENSION(NGROUP) :: ELBUF_TAB
333 TYPE (ELBUF_STRUCT_), DIMENSION(NGROUP,NXEL) :: XFEM_TAB
334 TYPE (XFEM_EDGE_) , DIMENSION(*) :: CRKEDGE
335 TYPE(H3D_DATABASE) :: H3D_DATA
336 TYPE (PINCH) :: PINCH_DATA
337 TYPE (SENSORS_) :: SENSORS
338C-----------------------------------------------
339 TYPE (GROUP_) , DIMENSION(NGRBRIC) :: IGRBRIC
340 TYPE (GROUP_) , DIMENSION(NGRQUAD) :: IGRQUAD
341 TYPE (GROUP_) , DIMENSION(NGRSHEL) :: IGRSH4N
342 TYPE (GROUP_) , DIMENSION(NGRSH3N) :: IGRSH3N
343 TYPE (GROUP_) , DIMENSION(NGRTRUS) :: IGRTRUSS
344 TYPE (GROUP_) , DIMENSION(NGRBEAM) :: IGRBEAM
345 TYPE (GROUP_) , DIMENSION(NGRSPRI) :: IGRSPRING
346 TYPE (GROUP_) , DIMENSION(NGRPART) :: IGRPART
347 TYPE (GROUP_) , DIMENSION(NGRNOD) :: IGRNOD
348C-----------------------------------------------
349 TYPE (FAILWAVE_STR_) ,TARGET :: FAILWAVE
350 TYPE (NLOCAL_STR_) ,TARGET :: NLOC_DMG
351 TYPE(sorting_comm_type), DIMENSION(NINTER), INTENT(inout) :: SORT_COMM ! structure for interface sorting comm
352 TYPE (STACK_PLY) :: STACK
353C-----------------------------------------------
354 TYPE(OUTPUT_),INTENT(INOUT) :: OUTPUT
355 TYPE(sh_offset_) :: SH_OFFSET_TAB
356 TYPE (LOADS_) ,INTENT(INOUT) :: LOADS
357 type (glob_therm_) ,intent(inout) :: glob_therm
358 type (pblast_) ,intent(inout) :: pblast
359 type (rbe3_) ,intent(inout) :: rbe3
360C-----------------------------------------------
361C L o c a l V a r i a b l e s
362C-----------------------------------------------
363 INTEGER IMUEL, I, J, K, NG, NINT7,NNOD,K2S,K0,IAD1,IDUM,LLL,
364 . LRBUF, LIBUF, ITY, IAD, NNBEM, ITYP,IROTG,NS,LF,LT,LL,L,
365 . l1,l2,isectr,nfr,ic,icr,nisub, ni25,nbr,nsensor,inloc
366 INTEGER JD(50),KD(50),JFI,KFI,NMN,II,NINOUT,NNO,NEL,IFLGADM,
367 . N,JJ,KK, NFT, ISOLNOD,NBS
368 INTEGER, DIMENSION(SENSORS%NSENSOR) :: INDEX_SENSOR
369 INTEGER, DIMENSION(:), ALLOCATABLE :: ISEND,IRECV
370 INTEGER :: ITIED,NINIVELTG
371 my_real :: rdum
372 CHARACTER ZONE*5
373 INTEGER VALUES(8)
374C=======================================================================
375 IDUM = 0
376 rdum = zero
377 isectr = 0
378 nsensor = sensors%NSENSOR
379C-----------------------------------------------
380C //
381C-----------------------------------------------
382C
383C Sequential part
384C
385 IF (itask == 0)THEN
386C zeroing ITYPTS for DTIX
387C
388 itypts=0
389C
390C kinematic conditions : arrays init. (RBY & INT20)
391C
392 CALL init_kyne(ikine,npby,lpby,tagslv_rby)
393C
394C reaction force (node array)
395C
396 cptreac = 0
397 IF (ireac == 1 ) CALL init_reac_nod(cptreac,nodreac,nthgrp,output%TH%ITHGRP,output%TH%ITHBUF)
398C
399C TH init for group of elems
400C
401 ngrth = 0
402 IF (igrelem == 1 ) THEN
403 CALL init_th_group(grth ,igrth ,nelem ,ngrth ,iparg ,
404 . ipart ,igrbric ,igrquad ,igrsh4n ,igrsh3n,
405 . igrtruss ,igrbeam ,igrspring)
406 ENDIF
407C----- reset initial mass
408 IF (imassi /= 0) THEN
409 ms(1:numnod)=ms0(1:numnod)
410 IF (iroddl /=0) in(1:numnod)=in0(1:numnod)
411 END IF
412C
413C parallel structures init.
414C
415 irotg=0
416 DO i=1,nrbe3
417 irotg=max(irotg,rbe3%IRBE3(6,i))
418 ENDDO
419 CALL spmd_max_i(irotg)
420 rbe3%irotg = irotg
421 IF(irotg==0) THEN
422 rbe3%irotg_sz = 5
423 ELSE
424 rbe3%irotg_sz = 10
425 ENDIF
426
427C---------RBE2----
428 irotg=0
429 DO i=1,nrbe2
430 ic = irbe2(4,i)
431 icr=(ic-512*(ic/512))/64
432 irotg=max(irotg,icr)
433 IF (irbe2(11,i)==0) irotg =1
434 ENDDO
435 CALL spmd_max_i(irotg)
436 IF(irotg==0) THEN
437 r2size = 4
438 ELSE
439 r2size = 8
440 ENDIF
441 ns = nrbe2
442 CALL spmd_max_i(ns)
443 IF (ns==0) r2size = 0
444 nfr = iad_rbe2(nspmd+1)-iad_rbe2(1)
445 IF (nspmd==1) THEN
446 rbe3%irotg_sz = 0
447 r2size = 0
448 ENDIF
449
450c
451C IRBE2 init.
452 CALL rbe2_init(irbe2 ,lrbe2 ,nmrbe2 ,fr_rbe2 ,fr_rbe2m,nfr)
453C
454 CALL mpp_init(
455 1 ipari ,isendto ,ircvfrom,intlist ,nbintc ,
456 2 isizxv ,ilenxv ,iad_elem,i2size ,itask ,
457 3 islen7 ,irlen7 ,islen11 ,irlen11 ,igrbric ,
458 4 nme17 ,islen17 ,irlen17 ,irlen7t ,islen7t ,
459 5 lindidel,lbufidel,irlen20 ,islen20 ,irlen20t,
460 6 islen20t,nbint20 ,irlen20e,islen20e,fr_rby ,
461 7 fr_rby6 ,npby ,irbkin_l,nrbykin_l,kindrby,
462 8 nsensor ,sensors%SENSOR_TAB,lbufidel24, intbuf_tab,
463 9 sort_comm,need_comm_int25_solid_erosion,comm_int25_solid_erosion )
464C
465 IF(idel7ng>0.OR.irad2r>0.OR.alemuscl_param%IALEMUSCL>0.OR.pdel>0) THEN
466 CALL chkinit(
467 2 ixs ,ixq ,ixc ,ixt ,ixp ,
468 3 ixr ,ixtg ,ixs10 ,ixs20 ,
469 4 ixs16 ,ixtg1 ,geo ,addcnel ,cnel ,
470 5 addtmpl ,iparg )
471 ENDIF
472
473C
474 IF (irad2r /= 0) THEN
475 CALL r2r_init(iexlnk ,itab,igrnod,x ,
476 2 ms ,in ,dd_r2r,weight ,iad_elem,
477 3 fr_elem,addcnel,cnel,ixc,iparg,icodt,icodr,
478 4 ibfv,d,rby,npby,xdp,stifn,stifr,dd_r2r_elem,
479 5 sdd_r2r_elem,weight_md,ilenxv,numsph_glo_r2r,
480 6 flg_sphinout_r2r,ipari,nloc_dmg)
481 END IF
482C
483 nfia = numnod*min(1,anim_v(4)+outp_v(4)+h3d_data%N_VECT_CONT)
484 nfea = nfia + numnod*min(1,anim_v(5)+outp_v(5)+h3d_data%N_VECT_FINT)
485 nfnca= nfea + numnod*min(1,anim_v(6)+outp_v(6)+h3d_data%N_VECT_FEXT)
486 nftca= nfnca+ numnod*min(1,anim_v(12)+outp_v(12)+h3d_data%N_VECT_PCONT)
487 nfoa = nftca+ numnod*min(1,anim_v(12)+outp_v(12)+h3d_data%N_VECT_PCONT)
488 nft2 = nfoa+ 2*(nsect+nrbody+nrwall)
489 nfnca2= nft2 + numnod*min(1,anim_v(13)+h3d_data%N_VECT_CONT2)
490 nftca2= nfnca2+ numnod*min(1,anim_v(27)+h3d_data%N_VECT_PCONT2)
491 ndma = numnod*min(1,anim_n(1)+outp_n(1)+h3d_data%N_SCAL_DT)
492 ndin = ndma +numnod*min(1,anim_n(2)+outp_n(2)+h3d_data%N_SCAL_DMAS)
493 ndma2 = ndin+numnod*min(1,anim_n(12)+outp_n(3)+h3d_data%N_SCAL_DINER)
494 ndama2 = ndma2+numelr*(anim_fe(11)+anim_fe(12)+anim_fe(13))
495 IF(iroddl/=0)THEN
496 DO ng=1,ninter
497 ity = ipari(7,ng)
498 IF(ity==2) THEN
499 nmn=ipari(6,ng)
500 DO ii = 1, nmn
501 i = intbuf_tab(ng)%MSR(ii)
502 intbuf_tab(ng)%NMAS(nmn+ii) = in(i)
503C For multidomains inertia of main nodes on multidomains interface msut be non zero
504 IF (irad2r==1) in(i)=max(em20,in(i))
505 END DO
506 END IF
507 END DO
508 END IF
509 dmas = zero
510 diner = zero
511C
512 IF(mcheck==0)ncycle=0
513 i7kglo = 0
514 nabfwr = 0
515C
516 i13a=1+2*nsnod
517 i13b=i13a+nsels
518 i13c=i13b+nselq
519 i13d=i13c+nselc
520 i13e=i13d+nselt
521 i13f=i13e+nselp
522 i13g=i13f+nselr
523 i13h=i13g+nselu
524 i13i=i13h+nseltg
525 i15ath=1+lipart1*(npart+nthpart)
526 i15a=i15ath+2*9*(npart+nthpart)
527 i15b=i15a+numels
528 i15c=i15b+numelq
529 i15d=i15c+numelc
530 i15e=i15d+numelt
531 i15f=i15e+numelp
532 i15g=i15f+numelr
533 i15h=i15g
534 i15i=i15h+numeltg
535 i15j=i15i+numelx
536 i15k=i15j+numsph
537 i35ath=1+lisub1*nsubs
538C
539 i87a = 1
540 i87b = i87a + 8 * numels + 6 * numels10 + 12 * numels20 + 8 * numels16
541 i87c = i87b + 4 * numelq
542 i87d = i87c + 4 * numelc
543 i87e = i87d + 2 * numelt
544 i87f = i87e + 2 * numelp
545 i87g = i87f + 3 * numelr
546 i87h = i87g + 3 * numeltg
547 i87h = i87h + 3 * numeltg6
548 i87i = i87h
549 i87j = i87i + 4 * nskymv0
550 i87k = i87j + 4 * nconld
551 i87l = i87k + 4 * glob_therm%NUMCONV
552 i87m = i87l + 4 * glob_therm%NUMRADIA
553 i87n = i87m + slloadp
554C I87O = I87N + 4 * GLOB_THERM%NFXFLUX
555C
556C----------------------------
557 maxnx=0
558 DO i=1,numelx
559 IF (kxx(3,i)>maxnx) maxnx=kxx(3,i)
560 ENDDO
561C----------------------------
562 DO i=1,npart
563 partsav(8,i)=parts0(i)
564 ENDDO
565C----------------------------
566 IF (ispmd==0)THEN
567 CALL date_and_time(startdate, starttime, zone, values)
568 WRITE(istdo,'(A,I2.2,A,I2.2,A,I4.4)') ' ',values(3),'/',values(2),'/',values(1)
569 WRITE(iout,'(A,I2.2,A,I2.2,A,I4.4)') ' ',values(3),'/',values(2),'/',values(1)
570 END IF
571C
572 manim = 0
573 mrest = 0
574 mstop = 0
575 ictlstop = 0
576 h3d_data%MH3D = 0
577 IF(dtin/=0. .AND. mcheck==0)THEN !go on with previous time step in case of checkpoint restart (/CHKPT)
578 IF(dt2old==zero)THEN
579 dt2old=dtin/onep1
580 ELSE
581 dt2old= min(dt2old,dtin/onep1)
582 ENDIF
583 ENDIF
584 IF(anim_v(26)+h3d_data%N_VECT_CONT_MAX >0) ifcontmax=1
585 IF(h3d_data%N_VECT_PCONT_MAX >0) ifcontpmax=1
586 IF(h3d_data%N_VECT_CONT2_MAX >0) ifcont2max=1
587 IF(h3d_data%N_VECT_PCONT2_MAX >0) ifcontp2max=1
588 IF(h3d_data%N_VECT_CONT2_MIN >0) ifcont2min=1
589 IF(h3d_data%N_VECT_PCONT2_MIN >0) ifcontp2min=1
590 IF(h3d_data%N_SCAL_CSE_FRIC >0) THEN
591 s_efric = numnod
592 IF(nintstamp/=0) s_efricg = numnodg
593 ENDIF
594 IF(ninefric >0) s_efricint = numnod
595 IF(ninefric_stamp >0) s_efricintg = numnodg
596C------------------------
597C PARAL. ARITH.
598C------------------------
599 IF(iparit==3) THEN
600 write(6,*) 'Non supported /PARITH option'
601 ELSEIF(iparit/=0) THEN
602C
603C parith/on
604C
605 IF(ivector==1)THEN
606 iad1 = numnod+2
607 ELSE
608 iad1 = 1
609 ENDIF
610 CALL assadd2(
611 1 pon%ADSKY ,pon%ADSKY(iad1),pon%FSKY ,pon%FSKYM ,iad_elem ,
612 2 fr_elem ,fr_nbcc ,procne,niskyfi ,addcni2 ,
613 3 procni2 ,iad_i2m ,fr_i2m,fr_nbcci2,addcni2(iad1),
614 4 pon%IADSDP ,pon%IADRCP ,pon%ISENDP,pon%IRECVP ,fthesky ,
615 5 niskyfie,inod_pxfem ,adsky_pxfem,procne_pxfem,
616 6 isendp_pxfem,irecvp_pxfem ,iadsdp_pxfem,iadrcp_pxfem,
617 7 fr_nbcc1,inod_crkxfem,adsky_crkxfem,procne_crkxfem,
618 8 isendp_crkxfem,irecvp_crkxfem,iadsdp_crkxfem,iadrcp_crkxfem,
619 9 condnsky,glob_therm)
620 ENDIF
621C
622 CALL fillipartl(
623 1 ipartl ,ipart(i15a),ipart(i15b),ipart(i15c),ipart(i15d),
624 2 ipart(i15e),ipart(i15f),ipart(i15g),ipart(i15h),ipart(i15i),
625 3 ipart(i15j),ipart(i15k),npartl )
626C------------------------
627C SPLIT GROUP FOR OPTIMIZATION
628C------------------------
629 CALL grpsplit(
630 1 iparg, igrouc, ngrouc, igrounc, ngrounc,
631 2 ixc,ixs,ixtg,ipm,igeo,pm,geo,tabmp_l,tab_mat)
632C--------------------------
633C FIND GROUP FOR SHELLS
634C--------------------------
635 IF(igroupflg(1) == 1 ) CALL findgroupc(iparg, igrouc, ngrouc, igroupc, igrouptg)
636C--------------------------
637C FIND GROUP FOR BRICKS
638C--------------------------
639 IF(igroupflg(2) == 1 ) CALL findgroups(iparg, igroups)
640C----------------------------------------------------------
641C TAG : NODES FROM ALL SECTIONS
642C----------------------------------------------------------
643 IF(isecut/=0)THEN
644 k0=nstrf(25)
645 DO i=1,nsect
646 nnod=nstrf(k0+6)
647 k2s=k0+30+nstrf(k0+14)
648 DO j=1,nnod
649 secfcum(4,nstrf(k2s),i)=1.
650 k2s=k2s+1
651 ENDDO
652 IF (nstrf(k0) >= 100 ) isectr = i
653 k0=nstrf(k0+24)
654 ENDDO
655 CALL section_init(nstrf,secbuf,nom_sect,isectr,nsect,ioldsect)
656 ENDIF
657C-----------------------------------------------------
658C SQRT H1, H2, H3 for shell elements
659C-----------------------------------------------------
660 DO i = 1, numgeo
661 igtyp = igeo(11,i)
662 IF(igtyp==1.OR.(igtyp>=9 .AND. igtyp<=11).OR.igtyp==16) THEN
663 geo(18,i) = sqrt(geo(13,i))
664 geo(19,i) = sqrt(geo(14,i))
665 geo(20,i) = sqrt(geo(15,i))
666 ENDIF
667 ENDDO
668C-----------------------------------------------------
669C optional SQRT(G), SQRT(A11) SQRT(A12), SQRT(NU), SQRT(SHF) for former restart file
670C-----------------------------------------------------
671 IF(pminver<6)THEN
672 DO i = 1, numgeo
673 geo(100,i) = sqrt(geo(38,i)) ! SHFSR
674 END DO
675 DO i = 1, nummat
676 IF(ipm(2,i)==999)cycle !possible negative square root otherwise PM(25)=CPE(gas)
677 pm(12,i) = sqrt(abs(pm(22,i))) ! GSR
678 pm(13,i) = sqrt(abs(pm(24,i))) ! A11SR
679 pm(14,i) = sqrt(abs(pm(25,i))) ! A12SR
680 pm(190,i)= sqrt(abs(pm(21,i))) ! NUSR
681 END DO
682 END IF
683C----------------------------------------------------------
684C INIT FLEX BODY
685C----------------------------------------------------------
686 IF (nfxbody>0) THEN
687 DO i=1,lenvar
688 fxbfp(i)=zero
689 fxbgrp(i)=zero
690 ENDDO
691 DO i=1,nfxbody
692 fxbefw(i)=zero
693 fxbgrw(i)=zero
694 fxbedp(i)=zero
695 ENDDO
696 ENDIF
697C----------------------------------------------------------
698C LWORKING ARRAY SIZES - AIRBAG BEM
699C----------------------------------------------------------
700 iad=0
701 lwibem=0
702 lwrbem=0
703 DO i=1,nvolu
704 ityp=monvol(iad+2)
705 IF (ityp==7) THEN
706 nnbem=monvol(iad+32)
707 lwibem=lwibem+1+nnbem
708 lwrbem=lwrbem+nnbem**2
709 ENDIF
710 iad=iad+nimv
711 ENDDO
712C----------------------------------------------------------
713C WORKING ARRAY SIZES - FLOW BEM
714C----------------------------------------------------------
715 iad=0
716 lwiflow=0
717 lwrflow=0
718 DO i=1,nflow
719 ityp=iflow(iad+2)
720 IF (ityp == 1 .OR.ityp == 3) THEN
721 lwiflow=lwiflow+iflow(iad+8)
722 lwrflow=lwrflow+iflow(iad+9)
723 ENDIF
724 iad=iad+liflow
725 ENDDO
726C----------------------------------------------------------
727C Domain Decomposition Weight computation
728C----------------------------------------------------------
729 IF(iddw>0) CALL initimeg(ngroup)
730C----------------------------------------------------------
731C INIT ADAPTIVE MESHING (sequentielle)
732C----------------------------------------------------------
733 IF(nadmesh/=0)THEN
734 CALL admini(ixc ,ipart(i15c),ixtg ,ipart(i15h),ipart,
735 . igeo,ipm ,iparg ,x ,ms ,
736 . in ,elbuf_tab ,sh4tree,ipadmesh,msc ,
737 . inc ,sh3tree ,mstg ,intg ,ptg ,
738 . sh4trim ,sh3trim,mscnd ,incnd ,pm ,
739 . mcp ,mcpc ,mcptg ,tagtrimc ,tagtrimtg,
740 . glob_therm%ITHERM_FE)
741!
742 CALL admordr(sh4tree,sh3tree,ixc,ixtg)
743 iadmesh=0
744 ngdone=1
745 END IF
746 IF(istatcnd/=0)THEN
747C ADAPTIVE MESHING + STATIC CONDENSATION
748 CALL cndordr(ipart,ipart(i15c),ipart(i15h),sh4tree,sh3tree)
749 END IF
750C----------------------------------------------------------
751C INIT MULTIPLICATEURS LAGRANGES (sequentielle)
752C----------------------------------------------------------
753 IF(lag_ncf+lag_ncl > 0)THEN
754 lag_sec=0
755C numbering incompatible options if NSPMD > 1
756 DO i = 1, ninter
757 IF(ipari(33,i)/=0)lag_sec=1
758 END DO
759 DO i = 1, nrwall
760 IF(nprw(i+5*nrwall)==1)lag_sec=1
761 END DO
762 IF(nbcslag+ngjoint+nrbylag > 0)lag_sec=1
763C NUMMPC + NFVLAG : ok (parallele SPMD)
764 END IF
765
766C-----------------------
767C INTERFACE TYPE 1
768C-----------------------
769 is_present_inter1 = -1
770C-----------------------
771C INTERFACE TYPE 18 KINE
772C-----------------------
773 int18kine=0
774 DO i=1, ninter
775 IF(ipari(7,i) == 7 .AND. ipari(34,i) == -2 .AND. ipari(22,i) == 7)THEN
776 int18kine=1
777 ENDIF
778 ENDDO
779C-----------------------
780C INTERFACE TYPE 7 FLAG + ITIED /= 0
781C-----------------------
782 int7itied = 0
783 DO i=1, ninter
784 ityp = ipari(7,i)
785 itied = ipari(85,i)
786 IF(ityp==7.AND.itied/=0)THEN
787 int7itied = 1
788 ENDIF
789 IF(ityp==10) int7itied = 1
790 ENDDO
791C-----------------------
792C INTERFACE TYPE 24 FLAG
793C-----------------------
794 int24use = 0
795 DO i=1, ninter
796 IF(ipari(7,i)==24)THEN
797 int24use = 1
798C Check if type 24 has E2E , set INT24E2EUSE
799 IF(ipari(59,i) >0) int24e2euse=1
800 ENDIF
801 ENDDO
802C-----------------------
803C INTERFACE TYPE 25 LIST
804C-----------------------
805 ni25 = 0
806 DO i=1, ninter
807 IF(ipari(7,i)==25)THEN
808 ni25 = ni25 + 1
809 intlist25(ni25)=i
810 ENDIF
811 ENDDO
812C-----------------------
813C SENSOR INTERFACE
814C-----------------------
815 IF (sensors%STABSEN > 0) THEN
816 DO n=1,ninter
817 nisub =ipari(36,n)
818 isensint(1,n) = sensors%TABSENSOR(n+1 + nsect) - sensors%TABSENSOR(n + nsect)
819C
820 IF (ipari(71,n)>0) THEN
821C-- sensor associated to all interfaces of type19
822 isensint(1,n) = isensint(1,ipari(71,n))
823 ENDIF
824C
825 DO i=1,nisub
826 isensint(i+1,n) = sensors%TABSENSOR(i +1 + nsect + ninter) -
827 . sensors%TABSENSOR(i + nsect + ninter)
828 ENDDO
829 ENDDO
830 ENDIF
831C-----------------------
832C INTERFACE TYPE 2 PENALITE
833C-----------------------
834 int2pen=0
835 DO i=1, ninter
836 IF (ipari(7,i) == 2 .AND. ipari(20,i) == 25) THEN
837 int2pen=1
838 EXIT
839 ENDIF
840 ENDDO
841
842C-----------------------
843C /IMPDISP/FGEO
844C-----------------------
845 fxvel_fgeo=0
846 DO n=1,nfxvel
847 IF (ibfv(13,n) > 0 ) THEN
848 fxvel_fgeo = 1
849 EXIT
850 ENDIF
851 ENDDO
852
853
854 ENDIF ! ITASK==0
855C-----------------------------------------------------
856C END OF SEQUENTIAL PART
857C-----------------------------------------------------
858C
859 CALL my_barrier()
860C--- // --------------------------------------
861C FORCE & MOMENTUM INIT
862C---------------------------------------------
863 IF(ninter/=0.AND.anim_v(4)+outp_v(4)+h3d_data%N_VECT_CONT >0) CALL zeror(fani(1,nodft),numnthread)
864 IF(anim_v(12)+outp_v(12)+h3d_data%N_VECT_PCONT>0) THEN
865 CALL zeror(fani(1,nfnca+nodft),numnthread)
866 CALL zeror(fani(1,nftca+nodft),numnthread)
867 END IF
868 IF(anim_n(2)+outp_n(2)+h3d_data%N_SCAL_DMAS >0)THEN
869#include "vectorize.inc"
870 DO i=nodft,nodlt
871 anin(i+ndma) = zero
872 ENDDO
873 ENDIF
874 IF(anim_n(12)+outp_n(3)+h3d_data%N_SCAL_DINER >0)THEN
875#include "vectorize.inc"
876 DO i=nodft,nodlt
877 anin(i+ndin) = zero
878 ENDDO
879 END IF
880 IF(anim_n(15) == 1 .OR. anim_n(16) == 1 .OR. h3d_data%N_SCAL_DAMA2 == 1)THEN
881#include "vectorize.inc"
882 DO i=nodft,nodlt
883 anin(ndama2+2*(i-1)+1) = zero
884 anin(ndama2+2*(i-1)+2) = zero
885 ENDDO
886 ENDIF
887C-----------------------------------------------
888C RESTARTING RADIOSS ENGINE.
889 IF (iparit==0) THEN
890 CALL zeror(a(1,ndtask),numnod)
891 IF(iroddl/=0)CALL zeror(ar(1,ndtask),numnod)
892 DO i=ndtask,ndtask+numnod-1
893 stifn(i)=em20
894 ENDDO
895 IF(iroddl/=0)THEN
896 DO i=ndtask,ndtask+numnod-1
897 stifr(i)=em20
898 ENDDO
899 ENDIF
900C
901 IF(kdtint/=0)THEN
902 CALL zero1(viscn(ndtask),numnod)
903 ENDIF
904C
905 IF (glob_therm%ITHERM_FE > 0) THEN
906 CALL zero1(fthe(ndtask),numnod)
907 ENDIF
908C
909 IF(sol2sph_flag/=0)THEN
910 CALL zero1(dmsph(ndtask),numnod)
911 ENDIF
912
913 IF (glob_therm%NODADT_THERM > 0) THEN
914 CALL zero1(condn(ndtask),numnod)
915 ENDIF
916C
917 IF(npinch > 0) THEN
918 CALL zeror(pinch_data%APINCH(1,ndtask),npinch)
919 DO i=ndtask,ndtask+numnod-1
920 pinch_data%STIFPINCH(i)=em20
921 ENDDO
922 ENDIF
923 ELSE ! IPARIT>0
924 CALL zeror(a(1,nodft),numnthread)
925 IF(iroddl/=0)CALL zeror(ar(1,nodft),numnthread)
926 DO i=nodft,nodlt
927 stifn(i)=em20
928 ENDDO
929 IF(iroddl/=0)THEN
930 DO i=nodft,nodlt
931 stifr(i)=em20
932 ENDDO
933 ENDIF
934 IF(kdtint/=0)THEN
935 CALL zero1(viscn(nodft),numnthread)
936 ENDIF
937C
938 IF (glob_therm%ITHERM_FE > 0 ) THEN
939 CALL zero1(fthe(nodft),numnthread)
940 ENDIF
941C
942 IF(nsphsol/=0)THEN
943 CALL zero1(dmsph(nodft),numnthread)
944 ENDIF
945
946 IF (glob_therm%NODADT_THERM > 0) THEN
947 CALL zero1(condn(nodft),numnthread)
948 ENDIF
949C
950 IF(npinch > 0) THEN
951 CALL zeror(pinch_data%APINCH(1,nodft),numnthread)
952 DO i=nodft,nodlt
953 pinch_data%STIFPINCH(i)=em20
954 ENDDO
955 ENDIF
956 ENDIF
957
958C
959 IF(iparit==0) THEN
960 IF(iroddl==0) THEN
961 DO i = nodft, nodlt
962 stifn(i) = stifn(i)*weight(i)
963 ENDDO
964 ELSE
965 DO i = nodft, nodlt
966 stifn(i) = stifn(i)*weight(i)
967 stifr(i) = stifr(i)*weight(i)
968 ENDDO
969 ENDIF
970 ENDIF
971C-----------------------------------------------------
972C INIT IMPLICIT
973C----------------------------------------------------------
974C --default values----
975 IF (itask==0) CALL imp_init(v,vr,iparg,ipm,igeo,elbuf_tab)
976C----------------------------------------------------------
977C INIT ADAPTIVE MESHING //
978C----------------------------------------------------------
979 IF(nadmesh/=0)THEN
980 iflgadm=0
981 CALL admgvid(
982 1 iparg ,elbuf_tab ,pon%FSKY ,pon%FSKY ,fthesky,
983 2 pon%IADC,pon%IAD_TG,iflgadm,igrouc,ngrouc,
984 3 condnsky ,glob_therm%NODADT_THERM)
985 END IF
986C
987C----------------------------------------------------------
988 IF( itask == 0) CALL kinini()
989C----------------------------------------------------------
990C INIT SELECTIVE MASS SCALING
991C----------------------------------------------------------
992 IF(idtmins == 1 .AND. idtmins_old == 1)THEN
993 IF(dtfacs /= dtfacs_old .OR. dtmins /= dtmins_old)THEN
994C Forget about previous mass scaling (reversibility)
995 admsms(nodft:nodlt)=zero
996 res_sms(1:3,nodft:nodlt)=zero
997 ELSEIF(idtgrs_old/=0)THEN
998 IF( idtgrs < 0.AND.
999 . -idtgrs /= igrpart(idtgrs_old)%ID) THEN
1000C
1001C Forget about previous mass scaling (reversibility)
1002 admsms(nodft:nodlt)=zero
1003 res_sms(1:3,nodft:nodlt)=zero
1004 ELSE
1005C ..as if single run
1006 END IF
1007 ELSEIF(idtgrs_old==0.AND.idtgrs/=0)THEN
1008C
1009C Forget about previous mass scaling (reversibility)
1010 admsms(nodft:nodlt)=zero
1011 res_sms(1:3,nodft:nodlt)=zero
1012 ELSE
1013C ..as if single run
1014 END IF
1015C
1016 ELSEIF(idtmins == 2 .AND. idtmins_old == 2)THEN
1017 IF(dtfacs /= dtfacs_old .OR. dtmins /= dtmins_old)THEN
1018C ..keep non diagonal mass from previous run
1019 ELSEIF(idtgrs_old/=0)THEN
1020 IF( idtgrs < 0.AND.
1021 . -idtgrs/= igrpart(idtgrs_old)%ID) THEN
1022C
1023C Forget about previous mass scaling (reversibility)
1024 IF(itask==0)THEN
1025 dmelc(1:numelc )=zero
1026 dmeltg(1:numeltg)=zero
1027 dmels(1:numels )=zero
1028 dmeltr(1:numelt )=zero
1029 dmelp(1:numelp )=zero
1030 dmelrt(1:numelr )=zero
1031 dmint2(1:4,1:i2nsn25)=zero
1032 END IF
1033 res_sms(1:3,nodft:nodlt)=zero
1034 ELSE
1035C ..as if single run
1036 END IF
1037 ELSEIF(idtgrs_old==0.AND.idtgrs/=0)THEN
1038C
1039C Forget about previous mass scaling (reversibility)
1040 IF(itask==0)THEN
1041 dmelc(1:numelc )=zero
1042 dmeltg(1:numeltg)=zero
1043 dmels(1:numels )=zero
1044 dmeltr(1:numelt )=zero
1045 dmelp(1:numelp )=zero
1046 dmelrt(1:numelr )=zero
1047 dmint2(1:4,1:i2nsn25)=zero
1048 END IF
1049 res_sms(1:3,nodft:nodlt)=zero
1050 ELSE
1051C ..as if single run
1052 END IF
1053C
1054 ELSEIF(idtmins == 1 .AND. idtmins_old /= idtmins)THEN
1055C
1056 admsms(nodft:nodlt)=zero
1057 res_sms(1:3,nodft:nodlt)=zero
1058C
1059 ELSEIF(idtmins == 2 .AND. idtmins_old /= idtmins)THEN
1060C
1061 IF(itask==0)THEN
1062 dmelc(1:numelc )=zero
1063 dmeltg(1:numeltg)=zero
1064 dmels(1:numels )=zero
1065 dmeltr(1:numelt )=zero
1066 dmelp(1:numelp )=zero
1067 dmelrt(1:numelr )=zero
1068 dmint2(1:4,1:i2nsn25)=zero
1069 END IF
1070 res_sms(1:3,nodft:nodlt)=zero
1071C
1072 ELSEIF(idtmins_int /= 0 .AND. idtmins_int_old /= idtmins_int)THEN
1073C
1074 res_sms(1:3,nodft:nodlt)=zero
1075C
1076 END IF
1077C
1078 IF(itask == 0) THEN
1079 nisky_sms=0
1080C enforce sorting contacts
1081 kforsms=0
1082 IF((idtmins==2.AND.idtmins_old/=idtmins).OR.
1083 . (idtmins_int/=0.AND.idtmins_int_old/=idtmins_int))THEN
1084 kforsms=1
1085 END IF
1086 ENDIF
1087C
1088 IF(anim_ply > 0.AND. itask == 0) THEN
1089 CALL spmd_anim_ply_init ()
1090 ENDIF
1091C
1092 IF (icrack3d > 0 .AND. itask == 0)THEN
1093 CALL anim_xfe_init(ixc,ixtg,inod_crkxfem,iel_crkxfem,
1094 . iadc_crkxfem,iadc_crkxfem(1+4*ecrkxfec))
1095 ENDIF
1096C-----------------------
1097C ITET=2 OF S10
1098C-----------------------
1099 IF(ns10e > 0) THEN
1100 IF (itask == 0) THEN
1101 IF(nspmd>1) THEN
1102 CALL s10cnds_dim(icnds10,itagnd,fr_elem,iad_elem,nbs )
1103 ALLOCATE (iad_cnds(nspmd+1),fr_cnds(nbs))
1104 CALL s10cnds_ini(icnds10,itagnd,fr_elem,iad_elem,iad_cnds,fr_cnds )
1105 ELSE
1106 ALLOCATE (iad_cnds(0),fr_cnds(0))
1107 END IF
1108
1109 CALL cndmasi2_dim(ipari,intbuf_tab,icnds10,itagnd,weight,nkend,
1110 1 iad_cnds,fr_cnds,nbs,nspmd)
1111 IF(nkend>0) THEN
1112 ALLOCATE (imap2nd(nkend),masi2nd0(nkend))
1113 CALL cndmasi2_ini(ipari,intbuf_tab,icnds10,itagnd,
1114 . nkend,imap2nd,masi2nd0,ms0,weight, itab )
1115 IF(mcheck>0) nkend = -nkend
1116 END IF
1117 CALL s10cndi2_ini(ipari,intbuf_tab,icnds10,itagnd,weight,
1118 . fr_cnds,iad_cnds,itab )
1120 1 addcncnd,procncnd,vnd ,v ,itab ,
1122 END IF
1123 CALL my_barrier()
1124 ENDIF
1125C-----------------------
1126C TMAX OF H3D
1127C-----------------------
1128 IF (itask == 0)
1129 . CALL tmax_ipart(iparg ,ipart ,ipart(i15a),ipart(i15c),
1130 . ipart(i15i),h3d_data)
1131 CALL ini_tmax(elbuf_tab ,iparg ,geo ,pm ,
1132 . ixs ,ixs10 ,ixs16 ,ixs20 ,ixq ,
1133 . ixc ,ixtg ,ixt ,ixp ,ixr ,
1134 . x ,d ,v ,iad_elem ,fr_elem ,
1135 . weight ,ipm ,igeo ,stack ,itask )
1136!$OMP SINGLE
1137 IF (failwave%WAVE_MOD > 0) THEN
1138 CALL spmd_failwave_boundaries(failwave,iad_elem,fr_elem)
1139 ENDIF
1140 ! Non-local regularization is activated
1141 IF (nloc_dmg%IMOD > 0) THEN
1142 CALL spmd_sub_boundaries(nloc_dmg,iad_elem,fr_elem)
1143 CALL nloc_count_solnod(ngroup, nparg, iparg, elbuf_tab, ixs, nixs, numels)
1144 ENDIF
1145C-----------------------
1146C DT_DC OF TSH
1147C-----------------------
1148 ntsheg =0
1149 ntshegg =0
1150 IF (idttsh>0) CALL dim_tshedg(elbuf_tab ,ntsheg, ixs ,iparg )
1151 IF(nspmd>1) THEN
1152 ntshegg = ntsheg
1153 CALL spmd_max_i(ntshegg)
1154 END IF
1155 IF (ntsheg > 0) THEN
1156 ALLOCATE (ienunl(2*ntsheg),alpha_dc(numnod))
1157 ienunl=0
1158 alpha_dc=one
1159 CALL ind_tshedg(elbuf_tab ,ienunl, ixs ,iparg )
1160 IF(nspmd>1) THEN
1161 nbs = iad_elem(1,nspmd+1)-iad_elem(1,1)
1162 ALLOCATE (isend(nbs),irecv(nbs))
1163 isend=0
1164 CALL tshcdcom_dim(ienunl,fr_elem,iad_elem,nbs,nbr ,
1165 . isend ,irecv )
1166 ALLOCATE (iad_stsh(nspmd+1),fr_stsh(nbs))
1167 CALL tshcdcom_ini(isend,iad_elem,fr_elem,iad_stsh,fr_stsh)
1168 ALLOCATE (iad_rtsh(nspmd+1),fr_rtsh(nbr))
1169 CALL tshcdcom_ini(irecv,iad_elem,fr_elem,iad_rtsh,fr_rtsh)
1170 DEALLOCATE(isend,irecv)
1171 END IF
1172 END IF
1173C-----------------------
1174C offset for contact
1175C-----------------------
1176 CALL inter_sh_offset_ini(
1177 . ngroup, nparg, iparg, npropg,
1178 . numgeo, geo, numelc, nixc,
1179 . ixc, numeltg, nixtg, ixtg,
1180 . numnod, nspmd, iad_elem, fr_elem,
1181 . sfr_elem, thke, elbuf_tab, sh_offset_tab,
1182 . iparit )
1183! inivel w/ Tstart
1184 niniveltg = loads%NINIVELT
1185 IF (nspmd>1) CALL spmd_max_i(niniveltg)
1186 loads%NINIVELT_G = niniveltg
1187 IF (tt == zero .AND. loads%NINIVELT > 0) THEN
1188 CALL inivel_init(
1189 . ngrnod, ngrbric, ngrquad, ngrsh3n,
1190 . igrnod, igrbric, igrquad, igrsh3n,
1191 . numskw, liskn, iskwn, numfram,
1192 . iframe, loads%NINIVELT,loads%INIVELT,sensors)
1193 END IF
1194
1195 DO n = 1, ninter
1196 CALL int_flushtime(intbuf_tab(n)%METRIC)
1197 ENDDO
1198!$OMP END SINGLE
1199C-------------------------------------------
1200 RETURN
1201 END
1202C
1203C-----------------------------------------------
1204!||====================================================================
1205!|| grpsplit ../engine/source/engine/resol_init.F
1206!||--- called by ------------------------------------------------------
1207!|| resol_init ../engine/source/engine/resol_init.F
1208!||--- calls -----------------------------------------------------
1209!|| sort_mid_pid ../engine/source/system/sort_mid_pid.F
1210!||====================================================================
1211 SUBROUTINE grpsplit(
1212 1 IPARG, IGROUC, NGROUC, IGROUNC, NGROUNC,
1213 2 IXC,IXS,IXTG,IPM,IGEO,PM,GEO,TABMP_L,TAB_MAT)
1214C----6------------------------------------------------------------------
1215C I m p l i c i t T y p e s
1216C-----------------------------------------------
1217#include "implicit_f.inc"
1218C-----------------------------------------------
1219C C o m m o n B l o c k s
1220C-----------------------------------------------
1221#include "com01_c.inc"
1222#include "param_c.inc"
1223#include "com04_c.inc"
1224C-----------------------------------------------------------------
1225C D u m m y A r g u m e n t s
1226C-----------------------------------------------
1227 INTEGER IPARG(NPARG,*),IGROUC(*),IGROUNC(*),
1228 . ngrouc, ngrounc,tabmp_l
1229
1230 INTEGER IXC(NIXC,*),IXS(NIXS,*),IXTG(NIXTG,*),
1231 . ipm(npropmi,*),igeo(npropgi,*)
1232 my_real pm(npropm,*),geo(npropg,*)
1233 my_real tab_mat(ngroup)
1234! tab_mat_prop
1235! 1 : shell
1236! 2 : tri
1237! 3 --> 9 : solid
1238! 3 : ISOL=8
1239! 4 : ISOL=10
1240! 5 : ISOL=16
1241! 6 : ISOL=20
1242! 7 : ISOL=6
1243! 8 : ISOL=4
1244! 9 : ISOL=others
1245C-----------------------------------------------
1246C L o c a l V a r i a b l e s
1247C-----------------------------------------------
1248 INTEGER NG, ITY, N_SHELL, N_SOL(7),N_TRI,MARQUEUR,MARQUEUR_2,MARQUEUR_3
1249 INTEGER I,J,II,JJ,K,INDI
1250 INTEGER COMPTEUR_MAT_PROP_SHELL,COMPTEUR_MAT_PROP_SOL,COMPTEUR_MAT_PROP_TRI,
1251 . MID,PID,MTN,NEL,NFT,FIRST,LAST,SHIFT,ISOL,GR_ID,GR_ID2
1252 INTEGER, DIMENSION(:,:), ALLOCATABLE :: TAB_SHELL_LOC,TAB_TRI_LOC
1253 INTEGER, DIMENSION(:,:,:), ALLOCATABLE :: TAB_SOL_LOC
1254 INTEGER, DIMENSION(:,:), ALLOCATABLE :: TAB_LOC
1255 INTEGER, DIMENSION(:), ALLOCATABLE :: IGROUC_SHELL,IGROUC_TRI,MID_SHELL,MID_TRI
1256 INTEGER, DIMENSION(:), ALLOCATABLE :: POIN_GROUP_MID_SHELL,POIN_GROUP_MID_TRI
1257 INTEGER, DIMENSION(:), ALLOCATABLE :: POIN_GROUP_PID_SHELL,POIN_GROUP_PID_TRI
1258 INTEGER, DIMENSION(:,:), ALLOCATABLE :: POIN_GROUP_MID_SOL
1259 INTEGER, DIMENSION(:,:), ALLOCATABLE :: POIN_GROUP_PID_SOL
1260 INTEGER, DIMENSION(:,:), ALLOCATABLE :: IGROUC_SOL,MID_SOL
1261 my_real POIDS_J,POIDS_J1
1262C-----------------------------------------------
1263C
1264 n_shell = 0
1265 n_sol(1:7) = 0
1266 n_tri = 0
1267 ngrouc = 0
1268 ngrounc = 0
1269
1270 ALLOCATE(tab_shell_loc(ngroup,5))
1271 ALLOCATE(tab_tri_loc(ngroup,5))
1272 ALLOCATE(tab_sol_loc(ngroup,5,7))
1273 ALLOCATE( igrouc_shell(ngroup),igrouc_tri(ngroup) )
1274 ALLOCATE( igrouc_sol(ngroup,7) )
1275
1276 ALLOCATE( poin_group_mid_shell(ngroup),poin_group_mid_tri(ngroup) )
1277 ALLOCATE( poin_group_pid_shell(ngroup),poin_group_pid_tri(ngroup) )
1278 ALLOCATE( poin_group_mid_sol(ngroup,7),poin_group_pid_sol(ngroup,7) )
1279
1280 ALLOCATE(mid_shell(nummat))
1281 ALLOCATE(mid_tri(nummat))
1282 ALLOCATE(mid_sol(nummat,7))
1283
1284 compteur_mat_prop_shell = 0
1285 mid_shell(1:nummat) = 0
1286 mid_tri(1:nummat) = 0
1287 mid_sol(1:nummat,1:7) = 0
1288
1289 DO ng = 1, ngroup
1290 ity =iparg(5,ng)
1291 IF(ity==3.OR.ity==7)THEN
1292 ngrouc = ngrouc + 1
1293 igrouc(ngrouc)=ng
1294 IF(ity==3) THEN
1295 n_shell = n_shell + 1
1296 nft = iparg(3,ng)+1
1297 mid = ixc(1,nft)
1298 pid = ixc(6,nft)
1299 mtn = iparg(1,ng)
1300 mid_shell(mid) = mid_shell(mid) + 1
1301 poin_group_mid_shell(n_shell) = mid
1302 poin_group_pid_shell(n_shell) = pid
1303 igrouc_shell(n_shell) = ng
1304
1305 tab_shell_loc(n_shell,1) = iparg(2,ng)
1306 tab_shell_loc(n_shell,2) = ng
1307 tab_shell_loc(n_shell,3) = mid
1308 tab_shell_loc(n_shell,4) = pid
1309 tab_shell_loc(n_shell,5) = ngrouc
1310
1311 ELSEIF(ity==7) THEN
1312 n_tri = n_tri + 1
1313 nft = iparg(3,ng)+1
1314 mid = ixtg(1,nft)
1315 pid = ixtg(5,nft)
1316 mtn = iparg(1,ng)
1317 mid_tri(mid) = mid_tri(mid) + 1
1318 poin_group_mid_tri(n_tri) = mid
1319 poin_group_pid_tri(n_tri) = pid
1320 igrouc_tri(n_tri) = ng
1321
1322 tab_tri_loc(n_tri,1) = iparg(2,ng)
1323 tab_tri_loc(n_tri,2) = ng
1324 tab_tri_loc(n_tri,3) = mid
1325 tab_tri_loc(n_tri,4) = pid
1326 tab_tri_loc(n_tri,5) = ngrouc
1327
1328 ENDIF
1329 ELSE
1330 ngrounc = ngrounc + 1
1331 igrounc(ngrounc)=ng
1332 IF(ity==1) THEN
1333 nft = iparg(3,ng)+1
1334 mid = ixs(1,nft)
1335 pid = ixs(10,nft)
1336 mtn = iparg(1,ng)
1337 isol = iparg(28,ng)
1338 IF(isol==4) THEN
1339 indi = 6
1340 ELSEIF(isol==6) THEN
1341 indi = 5
1342 ELSEIF(isol==8) THEN
1343 indi = 1
1344 ELSEIF(isol==10) THEN
1345 indi = 2
1346 ELSEIF(isol==16) THEN
1347 indi = 3
1348 ELSEIF(isol==20) THEN
1349 indi = 4
1350 ELSE
1351 indi = 7
1352 ENDIF
1353
1354 n_sol(indi) = n_sol(indi) + 1
1355 igrouc_sol(n_sol(indi),indi) = ng
1356
1357 tab_sol_loc(n_sol(indi),1,indi) = iparg(2,ng)
1358 tab_sol_loc(n_sol(indi),2,indi) = ng
1359 tab_sol_loc(n_sol(indi),3,indi) = mid
1360 tab_sol_loc(n_sol(indi),4,indi) = pid
1361 tab_sol_loc(n_sol(indi),5,indi) = ngrounc
1362
1363 mid_sol(mid,indi) = mid_sol(mid,indi) + 1
1364 poin_group_mid_sol(n_sol(indi),indi) = mid
1365 poin_group_pid_sol(n_sol(indi),indi) = pid
1366
1367 ENDIF
1368 END IF
1369 END DO
1370! -------------------------
1371 IF(n_shell>0) THEN
1372
1373 ALLOCATE( tab_loc(n_shell,3) )
1374 tab_loc(1:n_shell,1:3) = -1
1375
1376
1377 CALL sort_mid_pid(n_shell,igrouc_shell,
1378 1 poin_group_mid_shell,poin_group_pid_shell,
1379 2 mid_shell,tab_loc,tab_shell_loc,tab_mat)
1380
1381
1382 DO i = 1,n_shell
1383 j = tab_loc(i,1)
1384 ii = tab_shell_loc(i,5)
1385 jj = tab_shell_loc(j,2)
1386 igrouc(ii) = jj
1387 ENDDO
1388
1389 DEALLOCATE( tab_loc )
1390 ENDIF ! N_SHELL>0
1391! -------------------------
1392 IF(n_tri>0) THEN
1393
1394 ALLOCATE( tab_loc(n_tri,3) )
1395 tab_loc(1:n_tri,1:3) = -1
1396
1397
1398 CALL sort_mid_pid(n_tri,igrouc_tri,
1399 1 poin_group_mid_tri,poin_group_pid_tri,
1400 2 mid_tri,tab_loc,tab_tri_loc,tab_mat)
1401
1402
1403 DO i = 1,n_tri
1404 j = tab_loc(i,1)
1405 ii = tab_tri_loc(i,5)
1406 jj = tab_tri_loc(j,2)
1407 igrouc(ii) = jj
1408 ENDDO
1409
1410 DEALLOCATE( tab_loc )
1411 ENDIF ! N_TRI>0
1412! -------------------------
1413 DO indi=1,7
1414 IF(n_sol(indi)>0) THEN
1415
1416 ALLOCATE( tab_loc(n_sol(indi),3) )
1417 tab_loc(1:n_sol(indi),1:3) = -1
1418
1419
1420 CALL sort_mid_pid(n_sol(indi),igrouc_sol(1,indi),
1421 1 poin_group_mid_sol(1,indi),poin_group_pid_sol(1,indi),
1422 2 mid_sol(1,indi),tab_loc,tab_sol_loc(1,1,indi),tab_mat)
1423
1424
1425 DO i = 1,n_sol(indi)
1426 j = tab_loc(i,1)
1427 ii = tab_sol_loc(i,5,indi)
1428 jj = tab_sol_loc(j,2,indi)
1429 igrounc(ii) = jj
1430 ENDDO
1431
1432 DEALLOCATE( tab_loc )
1433 ENDIF ! N_SOL>0
1434 ENDDO
1435! -------------------------
1436
1437 DEALLOCATE(mid_shell)
1438 DEALLOCATE(mid_tri)
1439 DEALLOCATE(mid_sol)
1440
1441 DEALLOCATE( poin_group_mid_shell,poin_group_mid_tri )
1442 DEALLOCATE( poin_group_pid_shell,poin_group_pid_tri )
1443 DEALLOCATE( poin_group_mid_sol,poin_group_pid_sol )
1444
1445
1446 DEALLOCATE(tab_shell_loc)
1447 DEALLOCATE(tab_tri_loc)
1448 DEALLOCATE(tab_sol_loc)
1449 DEALLOCATE( igrouc_shell,igrouc_tri )
1450 DEALLOCATE( igrouc_sol )
1451
1452 RETURN
1453 END
1454C
1455C-----------------------------------------------
1456!||====================================================================
1457!|| fillipartl ../engine/source/engine/resol_init.F
1458!||--- called by ------------------------------------------------------
1459!|| resol_init ../engine/source/engine/resol_init.F
1460!||====================================================================
1461 SUBROUTINE fillipartl(
1462 1 IPARTL ,IPARTS,IPARTQ ,IPARTC ,IPARTT,
1463 2 IPARTP ,IPARTR,IPARTUR,IPARTTG,IPARTX,
1464 3 IPARTSP,IPARTIG3D,NPARTL)
1465C----6---------------------------------------------------------------7---------8
1466C I m p l i c i t T y p e s
1467C-----------------------------------------------
1468#include "implicit_f.inc"
1469C-----------------------------------------------
1470C C o m m o n B l o c k s
1471C-----------------------------------------------
1472#include "com04_c.inc"
1473#include "sphcom.inc"
1474C-----------------------------------------------------------------
1475C D u m m y A r g u m e n t s
1476C-----------------------------------------------!$OMP+PRIVATE(
1477 INTEGER IPARTS(*),IPARTQ(*),IPARTC(*),IPARTT(*),IPARTSP(*),
1478 . ipartp(*),ipartr(*),ipartur(*),iparttg(*),ipartx(*),
1479 . ipartl(*),ipartig3d(*),
1480 . npartl
1481C-----------------------------------------------
1482C L o c a l V a r i a b l e s
1483C-----------------------------------------------
1484 INTEGER I
1485C-----------------------------------------------
1486C //
1487C-----------------------------------------------
1488C
1489 DO i = 1, npart
1490 ipartl(i) = 0
1491 END DO
1492C
1493 DO i = 1, numels
1494 ipartl(iparts(i))=1
1495 END DO
1496C
1497 DO i = 1, numelq
1498 ipartl(ipartq(i))=1
1499 END DO
1500C
1501 DO i = 1, numelc
1502 ipartl(ipartc(i))=1
1503 END DO
1504C
1505 DO i = 1, numelt
1506 ipartl(ipartt(i))=1
1507 END DO
1508C
1509 DO i = 1, numelp
1510 ipartl(ipartp(i))=1
1511 END DO
1512C
1513 DO i = 1, numelr
1514 ipartl(ipartr(i))=1
1515 END DO
1516C
1517 DO i = 1, numeltg
1518 ipartl(iparttg(i))=1
1519 END DO
1520C
1521 DO i = 1, numelx
1522 ipartl(ipartx(i))=1
1523 END DO
1524C
1525 DO i = 1, numels
1526 ipartl(iparts(i))=1
1527 END DO
1528C
1529 DO i = 1, numsph
1530 ipartl(ipartsp(i))=1
1531 END DO
1532C
1533 DO i = 1, numelig3d
1534 ipartl(ipartig3d(i))=1
1535 END DO
1536C
1537 npartl = 0
1538 DO i = 1, npart
1539 IF(ipartl(i)>0)THEN
1540 npartl = npartl + 1
1541 ipartl(npartl) = i
1542 END IF
1543 END DO
1544C
1545 RETURN
1546 END
1547C
1548!||====================================================================
1549!|| smp_init ../engine/source/engine/resol_init.F
1550!||--- called by ------------------------------------------------------
1551!|| resol ../engine/source/engine/resol.f
1552!||--- calls -----------------------------------------------------
1553!|| omp_get_thread_num ../engine/source/engine/openmp_stub.F90
1554!||====================================================================
1555 SUBROUTINE smp_init(
1556 1 ITSK ,NODFTSK ,NODLTSK ,NUMNTSK,NDTSK,
1557 2 IPMTSK,PARTFTSK,PARTLTSK,NWAFTSK,IGMTSK,
1558 3 GREFTSK,GRELTSK)
1559C-----------------------------------------------
1560C I m p l i c i t T y p e s
1561C-----------------------------------------------
1562#include "implicit_f.inc"
1563C-----------------------------------------------
1564C C o m m o n B l o c k s
1565C-----------------------------------------------
1566#include "com01_c.inc"
1567#include "com04_c.inc"
1568#include "param_c.inc"
1569#include "task_c.inc"
1570C-----------------------------------------------
1571C D u m m y A r g u m e n t s
1572C-----------------------------------------------
1573 INTEGER ITSK, NODFTSK, NODLTSK, NUMNTSK, NDTSK,
1574 1 ipmtsk, partftsk, partltsk, nwaftsk, igmtsk,
1575 2 greftsk,greltsk
1576C-----------------------------------------------
1577C L o c a l V a r i a b l e s
1578C-----------------------------------------------
1579 INTEGER LENWA_T, OMP_GET_THREAD_NUM
1580 EXTERNAL omp_get_thread_num
1581C-----------------------------------------------
1582C S o u r c e L i n e s
1583C-----------------------------------------------
1584C
1585C Initialisation // SMP
1586C
1587 itsk = omp_get_thread_num()
1588 nodftsk = 1+itsk*numnod/ nthread
1589 nodltsk = (itsk+1)*numnod/nthread
1590 numntsk = nodltsk - nodftsk + 1
1591 ndtsk = 1 + itsk*numnod
1592 ipmtsk = 1 + itsk*npsav*npart
1593 partftsk = 1+itsk*npsav*npart/ nthread
1594 partltsk = (itsk+1)*npsav*npart/nthread
1595 lenwa_t = lenwa / nthread
1596 nwaftsk = 1+itsk*lenwa_t
1597 igmtsk = 1 + itsk*npsav*ngpe
1598 greftsk = 1+itsk*npsav*ngpe/ nthread
1599 greltsk = (itsk+1)*npsav*ngpe/nthread
1600c NWALTSK = (ITSK+1)*LENWA_T
1601c LOUT = ISPMD==0.AND.ITSK==0
1602C
1603 RETURN
1604 END
1605!||====================================================================
1606!|| init_kyne ../engine/source/engine/resol_init.F
1607!||--- called by ------------------------------------------------------
1608!|| resol_init ../engine/source/engine/resol_init.F
1609!||====================================================================
1610 SUBROUTINE init_kyne(IKINE,NPBY,LPBY,TAGSLV_RBY)
1611C----6---------------------------------------------------------------7---------8
1612C I m p l i c i t T y p e s
1613C-----------------------------------------------
1614#include "implicit_f.inc"
1615#include "comlock.inc"
1616C-----------------------------------------------
1617C C o m m o n B l o c k s
1618C-----------------------------------------------
1619#include "com04_c.inc"
1620#include "lagmult.inc"
1621#include "param_c.inc"
1622C-----------------------------------------------------------------
1623C D u m m y A r g u m e n t s
1624C-----------------------------------------------
1625 INTEGER IKINE(NUMNOD),NPBY(NNPBY,*),LPBY(*),TAGSLV_RBY(*)
1626C-----------------------------------------------
1627C L o c a l V a r i a b l e s
1628C-----------------------------------------------
1629 INTEGER N,I,J,K,NSN
1630
1631 DO J=1,numnod
1632 ikine(j)=0
1633 ENDDO
1634C-------------------------------------
1635C Traitement Rigid Bodies
1636C-------------------------------------
1637 k = 0
1638 DO n=1,nrbykin
1639 nsn = npby(2,n)
1640 DO i=1,nsn
1641 j=lpby(k+i)
1642 ikine(j) = (ikine(j)/2)*2 + 1
1643 ENDDO
1644 k = k + nsn
1645 ENDDO
1646C-------------------------------------------
1647 tagslv_rby(1:numnod)=0
1648 k=0
1649 DO n=1,nrbykin
1650 nsn=npby(2,n)
1651 IF(npby(7,n)>=1)THEN
1652 DO i=1,nsn
1653 tagslv_rby(lpby(i+k))=n
1654 ENDDO
1655 ENDIF
1656 k=k+nsn
1657 ENDDO
1658C-------------------------------------------
1659 DO n=1,nrbylag
1660 nsn = npby(2,n)
1661 DO i=1,nsn
1662 j=lpby(k+i)
1663 ikine(j) = (ikine(j)/2)*2 + 1
1664 ENDDO
1665 k = k + 3*nsn
1666 ENDDO
1667 RETURN
1668 END
subroutine admgvid(iparg, elbuf_tab, fskyv, fsky, fthesky, iadc, iadtg, iflg, igrouc, ngrouc, condnsky, nodadt_therm)
Definition admgvid.F:35
subroutine admini(ixc, ipartc, ixtg, iparttg, ipart, igeo, ipm, iparg, x, ms, in, elbuf_tab, sh4tree, ipadmesh, msc, inc, sh3tree, mstg, intg, ptg, sh4trim, sh3trim, mscnd, incnd, pm, mcp, mcpc, mcptg, tagtrimc, tagtrimtg, itherm_fe)
Definition admini.F:44
subroutine admordr(sh4tree, sh3tree, ixc, ixtg)
Definition admordr.F:35
subroutine anim_xfe_init(ixc, ixtg, inod_crk, iel_crk, iadc_crk, iadtg_crk)
subroutine assadd2(addcne, indsky, fsky, fskym, iad_elem, fr_elem, fr_nbcc, procne, niskyfi, addcni2, procni2, iad_i2m, fr_i2m, fr_nbcci2, indskyi2, iadsdp, iadrcp, isendp, irecvp, fthesky, niskyfie, inod_pxfem, addcne_pxfem, procne_pxfem, isendp_pxfem, irecvp_pxfem, iadsdp_pxfem, iadrcp_pxfem, fr_nbcc1, inod_crkxfem, addcne_crkxfem, procne_crkxfem, isendp_crkxfem, irecvp_crkxfem, iadsdp_crkxfem, iadrcp_crkxfem, condnsky, glob_therm)
Definition assadd2.F:43
subroutine chkinit(ixs, ixq, ixc, ixt, ixp, ixr, ixtg, ixs10, ixs20, ixs16, ixtg1, geo, addcnel, cnel, adsky, iparg)
Definition chkstfn3.F:263
subroutine cndordr(ipart, ipartc, iparttg, sh4tree, sh3tree)
Definition cndordr.F:32
#define my_real
Definition cppsort.cpp:32
subroutine dim_tshedg(elbuf_str, nedg, ixs, iparg)
Definition dim_tshedg.F:31
subroutine initimeg(ng)
Definition timer.F:1432
subroutine findgroups(iparg, igroups)
Definition findgroup.F:74
subroutine findgroupc(iparg, igrouc, ngrouc, igroupc, igrouptg)
Definition findgroup.F:30
subroutine zero1(r, n)
subroutine imp_init(v, vr, iparg, ipm, igeo, elbuf_tab)
Definition imp_init.F:38
subroutine spmd_max_i(n)
Definition imp_spmd.F:1362
subroutine ind_tshedg(elbuf_str, ienunl, ixs, iparg)
Definition ind_tshedg.F:31
subroutine ini_tmax(elbuf_tab, iparg, geo, pm, ixs, ixs10, ixs16, ixs20, ixq, ixc, ixtg, ixt, ixp, ixr, x, d, v, iad_elem, fr_elem, weight, ipm, igeo, stack, itask)
Definition ini_outmax.F:43
subroutine init_reac_nod(cptreac, nodreac, nthgrp, ithgrp, ithbuf)
subroutine init_th_group(gr, igr, nelem, ngrth, iparg, ipart, igrbric, igrquad, igrsh4n, igrsh3n, igrtruss, igrbeam, igrspring)
subroutine inivel(v, vr, svr, itabm1)
Definition inivel.F:35
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21
type(alemuscl_param_) alemuscl_param
integer, dimension(:), pointer fr_stsh
Definition dtdc_mod.F:42
integer, dimension(:), pointer iad_stsh
Definition dtdc_mod.F:42
integer, dimension(:), pointer iad_rtsh
Definition dtdc_mod.F:43
integer ntshegg
Definition dtdc_mod.F:39
integer, dimension(:), pointer ienunl
Definition dtdc_mod.F:40
integer, dimension(:), pointer fr_rtsh
Definition dtdc_mod.F:43
integer, dimension(:), pointer iad_cndm1
Definition ecdn_mod.F:48
integer, dimension(:), pointer fr_nbcccnd1
Definition ecdn_mod.F:57
integer, dimension(:), pointer iad_cnds
Definition ecdn_mod.F:50
integer, dimension(:), allocatable imap2nd
Definition ecdn_mod.F:64
integer, dimension(:), pointer fr_cndm
Definition ecdn_mod.F:47
integer, dimension(:), pointer fr_cndm1
Definition ecdn_mod.F:48
integer, dimension(:), pointer itagnd
Definition ecdn_mod.F:54
integer, dimension(:), pointer procncnd
Definition ecdn_mod.F:49
integer, dimension(:), pointer icnds10
Definition ecdn_mod.F:42
integer, dimension(:), pointer fr_cnds
Definition ecdn_mod.F:50
integer, dimension(:), pointer addcncnd
Definition ecdn_mod.F:49
integer, dimension(:), pointer iad_cndm
Definition ecdn_mod.F:47
integer nkend
Definition ecdn_mod.F:63
integer, dimension(:), pointer fr_nbcccnd
Definition ecdn_mod.F:57
integer ifcontp2max
Definition outmax_mod.F:69
integer ifcontmax
Definition outmax_mod.F:69
integer ifcont2max
Definition outmax_mod.F:69
integer ifcontp2min
Definition outmax_mod.F:69
integer ifcontpmax
Definition outmax_mod.F:69
integer ifcont2min
Definition outmax_mod.F:69
integer ninefric
Definition outputs_mod.F:65
integer s_efricg
Definition outputs_mod.F:64
integer s_efricintg
Definition outputs_mod.F:64
integer s_efric
Definition outputs_mod.F:64
integer s_efricint
Definition outputs_mod.F:64
integer ninefric_stamp
Definition outputs_mod.F:65
subroutine r2r_init(iexlnk, itab, igrnod, x, ms, in, dd_r2r, weight, iad_elem, fr_elem, addcnel, cnel, ixc, iparg, icodt, icodr, ibfv, dx, rby, npby, xdp, stifn, stifr, dd_r2r_elem, sdd_r2r_elem, weight_md, ilenxv, numsph_glo_r2r, flg_sphinout_r2r, ipari, nloc_dmg)
Definition r2r_init.F:70
subroutine rbe2_init(irbe2, lrbe2, nmrbe2, fr_rbe2, fr_rbe2m, nfr)
Definition rbe2f.F:623
subroutine resol(timers, element, nodes, coupling, af, iaf, iskwn, neth, ipart, nom_opt, kxx, ixx, ixtg, ixs, ixq, ixt, ixp, ixr, ifill, mat_elem, ims, npc, ibcl, ibfv, idum, las, laccelm, nnlink, lnlink, iparg, dd_iad, igrv, iexlnk, kinet, ipari, nprw, iconx, npby, lpby, lrivet, nstrf, ljoint, nodpor, monvol, ilink, llink, linale, neflsw, nnflsw, icut, cluster, itask, inoise, thke, damp, pm, skews, geo, eani, bufmat, bufgeo, bufsf, w, veul, fill, dfill, alph, wb, dsave, asave, msnf, tf, forc, vel, fsav, fzero, xlas, accelm, agrv, fr_wave, failwave, parts0, elbuf, rwbuf, sensors, rwsav, rby, rivet, secbuf, volmon, lambda, wa, fv, partsav, uwa, val2, phi, segvar, r, crflsw, flsw, fani, xcut, anin, tani, secfcum, bufnois, idata, rdata, iframe, kxsp, ixsp, nod2sp, ispsym, ispcond, xframe, spbuf, xspsym, vspsym, pv, fsavd, ibvel, lbvel, wasph, w16, isphio, lprtsph, lonfsph, vsphio, fbvel, lagbuf, ibcslag, iactiv, dampr, gjbufi, gjbufr, rbmpc, ibmpc, sphveln, nbrcvois, nbsdvois, lnrcvois, lnsdvois, nercvois, nesdvois, lercvois, lesdvois, npsegcom, lsegcom, nporgeo, ixtg1, npbyl, lpbyl, rbyl, igeo, ipm, madprt, madsh4, madsh3, madsol, madnod, madfail, iad_rby, fr_rby, fr_wall, iad_rby2, fr_rby2, iad_i2m, fr_i2m, addcni2, procni2, iadi2, fr_mv, iadmv2, fr_ll, fr_rl, iadcj, fr_cj, fr_sec, iad_sec, iad_cut, fr_cut, rg_cut, newfront, fr_mad, fxbipm, fxbrpm, fxbnod, fxbmod, fxbglm, fxbcpm, fxbcps, fxblm, fxbfls, fxbdls, fxbdep, fxbvit, fxbacc, fxbelm, fxbsig, fxbgrvi, fxbgrvr, eigipm, eigibuf, eigrpm, lnodpor, fr_i18, graphe, iflow, rflow, lgrav, dd_r2r, fasolfr, fr_lagf, llagf, lprw, icontact, rcontact, sh4tree, sh3tree, ipadmesh, padmesh, msc, mstg, inc, intg, ptg, iskwp, nskwp, isensp, nsensp, iaccp, naccp, ipart_state, acontact, pcontact, factiv, sh4trim, sh3trim, mscnd, incnd, ibfflux, fbfflux, rbym, irbym, lnrbym, icodrbym, ibcv, fconv, ibftemp, fbftemp, iad_rbym, fr_rbym, weight_rm, ms_ply, zi_ply, inod_pxfem, iel_pxfem, iadc_pxfem, adsky_pxfem, icode_ply, icodt_ply, iskew_ply, admsms, madclnod, nom_sect, mcpc, mcptg, dmelc, dmeltg, mssa, dmels, mstr, dmeltr, msp, dmelp, msrt, dmelrt, ibcr, fradia, res_sms, table, irbe2, lrbe2, iad_rbe2, fr_rbe2, phie, msf, procne_pxfem, iadsdp_pxfem, iadrcp_pxfem, icfield, lcfield, cfield, msz2, diag_sms, iloadp, lloadp, loadp, inod_crk, iel_crk, iadc_crk, adsky_crk, cne_crk, procne_crk, iadsdp_crk, iadrcp_crk, ibufssg_io, ibc_ply, dmint2, ibordnode, elbuf_tab, por, nodedge, iad_edge, fr_edge, fr_nbedge, crknodiad, lgauge, gauge, igaup, ngaup, nodlevxf, dd_r2r_elem, nodglobxfe, sph2sol, sol2sph, irst, dmsph, wagap, xfem_tab, elcutc, nodenr, kxfenod2elc, enrtag, rthbu f, kxig3d, ixig3d, knot, wige, wsmcomp, stack, cputime_mp_glob, cputime_mp, tab_ump, poin_ump, sol2sph_typ, irunn_bis, addcsrect, iad_frnor, fr_nor, procnor, iad_fredg, fr_edg, drape_sh4n, drape_sh3n, tab_mat, nativ0_sms, multi_fvm, segquadfr, ms_2d, h3d_data, subsets, igrnod, igrbric, igrquad, igrsh4n, igrsh3n, igrtruss, igrbeam, igrspring, igrpart, igrsurf, forneqs, nloc_dmg, iskwp_l, knotlocpc, knotlocel, pinch_data, tag_skins6, ale_connectivity, xcell, xface, ne_nercvois, ne_nesdvois, ne_lercvois, ne_lesdvois, ibcscyc, lbcscyc, t_monvol, id_global_vois, face_vois, dynain_data, fcont_max, ebcs_tab, diffusion, kloadpinter, loadpinter, dgaploadint, drapeg, user_windows, output, interfaces, dt, loads, python, dpl0cld, vel0cld, ndamp_vrel, id_damp_vrel, fr_damp_vrel, ndamp_vrel_rbyg, names_and_titles, unitab, liflow, lrflow, glob_therm, pblast, rbe3)
Definition resol.F:633
subroutine fillipartl(ipartl, iparts, ipartq, ipartc, ipartt, ipartp, ipartr, ipartur, iparttg, ipartx, ipartsp, ipartig3d, npartl)
subroutine grpsplit(iparg, igrouc, ngrouc, igrounc, ngrounc, ixc, ixs, ixtg, ipm, igeo, pm, geo, tabmp_l, tab_mat)
subroutine smp_init(itsk, nodftsk, nodltsk, numntsk, ndtsk, ipmtsk, partftsk, partltsk, nwaftsk, igmtsk, greftsk, greltsk)
subroutine init_kyne(ikine, npby, lpby, tagslv_rby)
subroutine resol_init(itask, fr_nbcc, isendto, ircvfrom, iad_elem, fr_elem, itabm1, ipari, iparg, itab, ixs10, ixs20, i13a, i13b, i13c, i13d, i13e, i13f, i13g, i13h, i13i, i15a, i15b, i15c, i15d, i15e, i15f, i15g, i15h, i15i, i87a, i87b, i87c, i87d, i87e, i87f, i87g, nfia, nfea, nfoa, ndma, ndma2, nodft, nodlt, ndtask, numnthread, ixs16, ixs, ixq, ixc, ixt, ixp, ixr, ixtg, pon, ikine, a, ar, v, vr, x, d, ms, in, stifn, stifr, dmas, diner, fani, anin, wa, uwa, pm, geo, partsav, parts0, monvol, i87h, i87i, i87j, i87k, i15j, kxx, secbuf, secfcum, nstrf, igrnod, iexlnk, xframe, ixtg1, ib, viscn, dd_r2r, elbuf, ipart, madprt, madsh4, madsh3, madsol, madnod, madfail, igeo, intlist, nbintc, procne, niskyfi, weight, isizxv, ilenxv, addcni2, procni2, iad_i2m, fr_i2m, fr_nbcci2, i2size, fr_mad, lwibem, lwrbem, fxbfp, fxbefw, fxbedp, fxbgrp, fxbgrw, ndin, islen7, irlen7, islen11, irlen11, lwiflow, lwrflow, iflow, addcnel, cnel, addtmpl, ipartl, npartl, nfnca, nftca, i15ath, i35ath, ipm, sh4tree, ipadmesh, msc, inc, sh3tree, mstg, intg, ptg, fthe, fthesky, ftheskyi, nme17, islen17, irlen17, irlen7t, islen7t, lindidel, lbufidel, sh4trim, sh3trim, mscnd, incnd, irlen20, islen20, irlen20t, islen20t, nbint20, irlen20e, islen20e, niskyfie, mcp, ms0, inod_pxfem, iel_pxfem, iadc_pxfem, adsky_pxfem, icodt, icodr, ibfv, admsms, nodreac, igrouc, ngrouc, igrounc, ngrounc, fr_rby, fr_rby6, npby, nom_sect, mcpc, mcptg, grth, igrth, nelem, lag_sec, nprw, diag_sms, dmelc, dmeltg, ngrth, nft2, dmels, dmeltr, dmelp, dmelrt, res_sms, i87l, irbe2, lrbe2, nmrbe2, iad_rbe2, fr_rbe2, fr_rbe2m, r2size, lpby, procne_pxfem, isendp_pxfem, irecvp_pxfem, iadsdp_pxfem, iadrcp_pxfem, fr_nbcc1, rby, int18kine, xdp, i87m, inod_crkxfem, iel_crkxfem, iadc_crkxfem, adsky_crkxfem, procne_crkxfem, isendp_crkxfem, irecvp_crkxfem, iadsdp_crkxfem, iadrcp_crkxfem, int24use, ndama2, igroupc, igrouptg, igroups, igroupflg, dmint2, irbkin_l, nrbykin_l, kindrby, elbuf_tab, sensors, dd_r2r_elem, sdd_r2r_elem, kinet, weight_md, dmsph, ioldsect, lbufidel24, intbuf_tab, numsph_glo_r2r, flg_sphinout_r2r, i15k, condn, condnsky, kxfenod2elc, elcutc, nodedge, iad_edge, crknodiad, fr_edge, fr_nbedge, nodlevxf, crkedge, xfem_tab, isensint, nisubmax, intlist25, int24e2euse, tabmp_l, i87n, tab_mat, h3d_data, tagtrimc, tagtrimtg, igrbric, igrquad, igrsh4n, igrsh3n, igrtruss, igrbeam, igrspring, igrpart, forneqs, int7itied, fxvel_fgeo, failwave, nloc_dmg, pinch_data, slloadp, tagslv_rby, nfnca2, nftca2, in0, sort_comm, stack, output, thke, sfr_elem, sh_offset_tab, need_comm_int25_solid_erosion, comm_int25_solid_erosion, iskwn, iframe, loads, glob_therm, pblast, rbe3)
Definition resol_init.F:173
subroutine s10cnds_ini(icnds10, itags, fr_elem, iad_elem, iad_cdns, fr_cdns)
Definition s10cndf.F:861
subroutine s10cnd_ini(icnds10, itagnd, iad_cndm, fr_cndm, fr_nbcccnd, addcncnd, procncnd, vnd, v, itab, iad_cndm1, fr_cndm1, fr_nbcccnd1)
Definition s10cndf.F:399
subroutine s10cndi2_ini(ipari, intbuf_tab, icnds10, itagnd, weight, fr_cnds, iad_cnds, itab)
Definition s10cndf.F:513
subroutine cndmasi2_dim(ipari, intbuf_tab, icnds10, itagnd, weight, nkend, iad_cnds, fr_cnds, s_fr, nspmd)
Definition s10cndf.F:950
subroutine s10cnds_dim(icnds10, itags, fr_elem, iad_elem, nbdds)
Definition s10cndf.F:819
subroutine cndmasi2_ini(ipari, intbuf_tab, icnds10, itagnd, nkend, imap2nd, masi2nd0, ms, weight, itab)
Definition s10cndf.F:1049
subroutine section_init(nstrf, secbuf, nom_sect, isectr, nsect, ioldsect)
subroutine sort_mid_pid(n_shell, igrouc_shell, poin_group_mid_shell, poin_group_pid_shell, mid_shell, tab_loc, tab_shell_loc, tab_mat)
subroutine spmd_failwave_boundaries(failwave, iad_elem, fr_elem)
subroutine spmd_sub_boundaries(nloc_dmg, iad_elem, fr_elem)
subroutine mpp_init(ipari, isendto, ircvfrom, intlist, nbintc, isizxv, ilenxv, iad_elem, i2size, itask, islen7, irlen7, islen11, irlen11, igrbric, nme17, islen17, irlen17, irlen7t, islen7t, lindidel, lbufidel, irlen20, islen20, irlen20t, islen20t, nbint20, irlen20e, islen20e, fr_rby, fr_rby6, npby, irbkin_l, nrbykin_l, kindrby, nsensor, sensor_tab, lbufidel24, intbuf_tab, sort_comm, need_comm_int25_solid_erosion, comm_int25_solid_erosion)
subroutine kinini(ikine)
Definition kinini.F:29
subroutine spmd_anim_ply_init(igeo, geo, iparg, ixc, ixtg, ipartc, ipartq, iparttg, stack)
subroutine my_barrier
Definition machine.F:31
subroutine tmax_ipart(iparg, ipart, iparts, ipartc, ipartg, h3d_data)
Definition tmax_ipart.F:34
subroutine tshcdcom_dim(ienunl, fr_elem, iad_elem, nbdds, nbddr, isend, irecv)
subroutine tshcdcom_ini(isend, iad_elem, fr_elem, iad_stsh, fr_stsh)
subroutine zeror(a, n)
Definition zero.F:39