1 | MODULE lib_feti |
---|
2 | !!============================================================================== |
---|
3 | !! *** MODULE lib_feti *** |
---|
4 | !! Ocean solver : FETI library |
---|
5 | !!============================================================================== |
---|
6 | |
---|
7 | #if defined key_feti |
---|
8 | !! feti_creadr |
---|
9 | !! feti_prext |
---|
10 | !! feti_init |
---|
11 | !! feti_iclr |
---|
12 | !! feti_vclr |
---|
13 | !! feti_vmov |
---|
14 | !! feti_vneg (used in lib_feti) |
---|
15 | !! feti_vsub |
---|
16 | !! feti_vadd (used in lib_feti) |
---|
17 | !! feti_ldlt (used in lib_feti) |
---|
18 | !! feti_desremo (used in lib_feti) |
---|
19 | !! feti_mxv (used in lib_feti) |
---|
20 | !! feti_mxvadd (used in lib_feti) |
---|
21 | !! feti_saut (used in lib_feti) |
---|
22 | !! feti_moyen (used in lib_feti) |
---|
23 | !! feti_assemb (used in lib_feti) |
---|
24 | !! feti_poids |
---|
25 | !! feti_probit (used in lib_feti) |
---|
26 | !! feti_inisub |
---|
27 | !! feti_subound |
---|
28 | !! feti_subdir |
---|
29 | !! feti_listdir |
---|
30 | !! feti_blodir |
---|
31 | !! feti_numblo |
---|
32 | !! feti_blomat |
---|
33 | !! feti_blomat1 |
---|
34 | !! feti_nullsp |
---|
35 | !! feti_project (used in lib_feti) |
---|
36 | !! feti_etesolv (used in lib_feti) |
---|
37 | !! feti_osaxpy (used in lib_feti) |
---|
38 | !! feti_proe (used in lib_feti) |
---|
39 | !! feti_proet (used in lib_feti) |
---|
40 | !! feti_consrhs (used in lib_feti) |
---|
41 | !! feti_dualschur |
---|
42 | !! feti_extend (used in lib_feti) |
---|
43 | !! feti_restri (used in lib_feti) |
---|
44 | !! feti_halfint |
---|
45 | !! feti_front |
---|
46 | !! feti_calinv (used in lib_feti) |
---|
47 | !! feti_nbdia (used in lib_feti) |
---|
48 | !! feti_liblod (used in lib_feti) |
---|
49 | !! feti_liblos (used in lib_feti) |
---|
50 | !! feti_resloc |
---|
51 | !! feti_pblos (used in lib_feti) |
---|
52 | !! feti_pbloi (used in lib_feti) |
---|
53 | !! feti_proax |
---|
54 | !! |
---|
55 | !! EXTERNAL used (BLAS): sgemv, strsv, saxpy, irecv |
---|
56 | !! |
---|
57 | !! References : |
---|
58 | !! Guyon, M, Roux, F-X, Chartier, M and Fraunie, P, 1994 : |
---|
59 | !! A domain decomposition solver to compute the barotropic |
---|
60 | !! component of an OGCM in the parallel processing field. |
---|
61 | !! Ocean Modelling, issue 105, december 94. |
---|
62 | !! |
---|
63 | !! Modifications : |
---|
64 | !! original : 98-02 (M. Guyon) |
---|
65 | !!------------------------------------------------------------------------ |
---|
66 | !! * Modules used |
---|
67 | USE lib_mpp ! distribued memory computing |
---|
68 | |
---|
69 | CONTAINS |
---|
70 | |
---|
71 | subroutine feti_creadr(ial,ialmax,lmax,ltab,iatab,chtab) |
---|
72 | ! ial = derniere adresse libre dans le super-tableau de longueur lmax . |
---|
73 | ! ltab = taille de la structure de donnee tab , d'adresse iatab . |
---|
74 | implicit none |
---|
75 | external mynode |
---|
76 | |
---|
77 | integer ial,ialmax,lmax,ltab,iatab |
---|
78 | integer iialmax,iibidon,mynode |
---|
79 | character(len=12) chtab |
---|
80 | |
---|
81 | iatab=ial |
---|
82 | ial=iatab+ltab |
---|
83 | ialmax=max0(ial,ialmax) |
---|
84 | if(ial > lmax) then |
---|
85 | write(0,*)' creation de la structure de donnee ',chtab |
---|
86 | write(0,*)' depassement taille super-tableau :' |
---|
87 | write(0,*)' taille maximale , taille demandee :',lmax,ial |
---|
88 | stop |
---|
89 | endif |
---|
90 | |
---|
91 | !... test de la memoire |
---|
92 | |
---|
93 | iialmax=-ialmax |
---|
94 | call mpp_min(iialmax,1,iibidon) |
---|
95 | iialmax=-iialmax |
---|
96 | ! if(mynode() == 0)then |
---|
97 | ! write(0,*)' creation de la structure de donnee ',chtab |
---|
98 | ! write(0,*)'taille max vaut : ',iialmax |
---|
99 | ! endif |
---|
100 | |
---|
101 | end subroutine feti_creadr |
---|
102 | !*********************************************************************** |
---|
103 | subroutine feti_prext(n,x) |
---|
104 | implicit none |
---|
105 | integer n |
---|
106 | REAL(kind=8) x(n) |
---|
107 | REAL(kind=8) xmin,xmax |
---|
108 | integer i |
---|
109 | |
---|
110 | xmin=x(1) |
---|
111 | xmax=x(1) |
---|
112 | do i=2,n |
---|
113 | xmin=min(xmin,x(i)) |
---|
114 | xmax=max(xmax,x(i)) |
---|
115 | end do |
---|
116 | ! write(0,100)xmin,xmax |
---|
117 | 100 format(8d10.3) |
---|
118 | |
---|
119 | end subroutine feti_prext |
---|
120 | !********************************************************************* |
---|
121 | subroutine feti_init(noeuds,x) |
---|
122 | implicit none |
---|
123 | integer noeuds |
---|
124 | REAL(kind=8) x(noeuds) |
---|
125 | integer ji |
---|
126 | |
---|
127 | do ji=1,noeuds |
---|
128 | x(ji)=1. |
---|
129 | end do |
---|
130 | |
---|
131 | end subroutine feti_init |
---|
132 | !*********************************************************************** |
---|
133 | subroutine feti_iclr(n,x) |
---|
134 | implicit none |
---|
135 | integer n |
---|
136 | integer x(n) |
---|
137 | integer i |
---|
138 | |
---|
139 | do i=1,n |
---|
140 | x(i)=0 |
---|
141 | end do |
---|
142 | |
---|
143 | end subroutine feti_iclr |
---|
144 | !*********************************************************************** |
---|
145 | subroutine feti_vclr(n,x) |
---|
146 | implicit none |
---|
147 | integer n |
---|
148 | REAL(kind=8) x(n) |
---|
149 | integer i |
---|
150 | |
---|
151 | do i=1,n |
---|
152 | x(i)=0.e0 |
---|
153 | end do |
---|
154 | |
---|
155 | return |
---|
156 | end |
---|
157 | !*********************************************************************** |
---|
158 | subroutine feti_vmov(n,x,y) |
---|
159 | implicit none |
---|
160 | integer n |
---|
161 | REAL(kind=8) x(n),y(n) |
---|
162 | integer i |
---|
163 | do i=1,n |
---|
164 | y(i)=x(i) |
---|
165 | end do |
---|
166 | return |
---|
167 | end |
---|
168 | !*********************************************************************** |
---|
169 | subroutine feti_vneg(n,x,y) |
---|
170 | implicit none |
---|
171 | integer n |
---|
172 | REAL(kind=8) x(n),y(n) |
---|
173 | integer i |
---|
174 | |
---|
175 | do i=1,n |
---|
176 | y(i)=-x(i) |
---|
177 | end do |
---|
178 | |
---|
179 | return |
---|
180 | end |
---|
181 | !*********************************************************************** |
---|
182 | subroutine feti_vsub(n,x,y,z) |
---|
183 | implicit none |
---|
184 | integer n |
---|
185 | REAL(kind=8) x(n),y(n),z(n) |
---|
186 | integer i |
---|
187 | do i=1,n |
---|
188 | z(i)=x(i)-y(i) |
---|
189 | end do |
---|
190 | return |
---|
191 | end |
---|
192 | !*********************************************************************** |
---|
193 | subroutine feti_vadd(n,x,y,z) |
---|
194 | implicit none |
---|
195 | integer n |
---|
196 | REAL(kind=8) x(n),y(n),z(n) |
---|
197 | integer i |
---|
198 | |
---|
199 | do i=1,n |
---|
200 | z(i)=y(i)+x(i) |
---|
201 | end do |
---|
202 | |
---|
203 | return |
---|
204 | end |
---|
205 | !*********************************************************************** |
---|
206 | subroutine feti_ldlt(n,a,d,v,ndlblo,lisblo,nmdlblo,j0) |
---|
207 | |
---|
208 | !...ndlblo nombre de liberte bloques = nombre de pivots nuls = dim(Ker(Ep)) |
---|
209 | |
---|
210 | !...nmdlblo = Max (ndlblo(Ep)) <= superstructure management in f77 context |
---|
211 | ! p=1,Np |
---|
212 | |
---|
213 | implicit none |
---|
214 | integer n,ndlblo,nmdlblo |
---|
215 | integer j0 |
---|
216 | REAL(kind=8) a(n,n),v(n),d(n) |
---|
217 | integer lisblo(nmdlblo) |
---|
218 | integer i,j,k |
---|
219 | REAL(kind=8) piv,epsilo |
---|
220 | |
---|
221 | ! pivot de reference . |
---|
222 | piv=a(1,1) |
---|
223 | ! epsilon . |
---|
224 | epsilo=1.e-10 |
---|
225 | ! factorisation de Crout d'une matrice pleine . |
---|
226 | do j = 1 , n |
---|
227 | ! calcul du produit de la ligne numero j par la diagonale . |
---|
228 | do k = 1 , j - 1 |
---|
229 | v(k) = a(j,k) * a(k,k) |
---|
230 | end do |
---|
231 | ! calcul des termes de la partie inferieure de la colonne . |
---|
232 | call sgemv('N',n-j+1,j-1,-1.e0,a(j,1),n,v,1,1.e0,a(j,j),1) |
---|
233 | ! detection des pivots nuls . |
---|
234 | if(a(j,j) < (piv*epsilo)) then |
---|
235 | ndlblo=ndlblo+1 |
---|
236 | lisblo(ndlblo)=j0+j |
---|
237 | do i=1,j-1 |
---|
238 | a(i,j)=0.e0 |
---|
239 | end do |
---|
240 | do i=j+1,n |
---|
241 | a(j,i)=0.e0 |
---|
242 | end do |
---|
243 | a(j,j)=0.e0 |
---|
244 | d(j)=0.e0 |
---|
245 | else |
---|
246 | piv=a(j,j) |
---|
247 | d(j)=1.e0/a(j,j) |
---|
248 | do i = j+1 , n |
---|
249 | a(i,j) = a(i,j) * d(j) |
---|
250 | end do |
---|
251 | endif |
---|
252 | end do |
---|
253 | ! recopie de la partie triangulaire superieure . |
---|
254 | do i = 1, n |
---|
255 | do j = i+1, n |
---|
256 | a(i,j)=a(j,i) |
---|
257 | end do |
---|
258 | end do |
---|
259 | do i=1,n |
---|
260 | a(i,i)=d(i) |
---|
261 | end do |
---|
262 | |
---|
263 | return |
---|
264 | end |
---|
265 | !**********************************************************************c |
---|
266 | subroutine feti_desremo(n,j0,a,y,x) |
---|
267 | implicit none |
---|
268 | integer n,i,j,j0 |
---|
269 | REAL(kind=8) a(n,n),x(n),y(n) |
---|
270 | |
---|
271 | ! descente remontee du systeme . |
---|
272 | ! initialisation du vecteur resultat . |
---|
273 | do i=1,n |
---|
274 | x(i)=y(i) |
---|
275 | end do |
---|
276 | ! descente du systeme L x = x , L est a diagonale unitaire . |
---|
277 | call strsv('L','N','U',n-j0,a(j0+1,j0+1),n,x(j0+1),1) |
---|
278 | ! division par la diagonale . |
---|
279 | do j=j0+1,n |
---|
280 | x(j)=x(j)*a(j,j) |
---|
281 | end do |
---|
282 | ! remontee du systeme U x = x , U est a diagonale unitaire . |
---|
283 | call strsv('U','N','U',n-j0,a(j0+1,j0+1),n,x(j0+1),1) |
---|
284 | |
---|
285 | return |
---|
286 | end |
---|
287 | !*********************************************************************** |
---|
288 | subroutine feti_mxv(n,m,a,x,y) |
---|
289 | ! calcul du produit matrice-vecteur : y = A x |
---|
290 | implicit none |
---|
291 | integer n,m |
---|
292 | REAL(kind=8) a(n,m),x(m),y(n) |
---|
293 | REAL(kind=8) alpha,beta |
---|
294 | ! initialisation du vecteur y . |
---|
295 | alpha=1.e0 |
---|
296 | beta=0.e0 |
---|
297 | call sgemv('N',n,m,alpha,a,n,x,1,beta,y,1) |
---|
298 | return |
---|
299 | end |
---|
300 | !**********************************************************************c |
---|
301 | subroutine feti_mxvadd(n,m,a,x,y) |
---|
302 | ! calcul du produit matrice-vecteur : y = y + A x |
---|
303 | implicit none |
---|
304 | integer n,m |
---|
305 | REAL(kind=8) a(n,m),x(m),y(n) |
---|
306 | REAL(kind=8) alpha,beta |
---|
307 | ! initialisation du vecteur y . |
---|
308 | |
---|
309 | alpha=1.e0 |
---|
310 | beta=1.e0 |
---|
311 | if((n > 0).and.(m > 0))then |
---|
312 | call sgemv('N',n,m,alpha,a,n,x,1,beta,y,1) |
---|
313 | endif |
---|
314 | |
---|
315 | return |
---|
316 | end |
---|
317 | !**********************************************************************c |
---|
318 | subroutine feti_saut(numnes,listnes, & |
---|
319 | plistin,numins,listin, & |
---|
320 | numnos,v,w,bufin,bufout) |
---|
321 | !! narea,numnos,v,w,bufin,bufout) |
---|
322 | |
---|
323 | ! computation of the jump of a field on the interfaces, in the case |
---|
324 | ! where subdomain number i is allocated to processor number i-1 . |
---|
325 | |
---|
326 | ! INPUTS : numnes = number of neighbouring subdomains |
---|
327 | ! listnes = list of neighbouring subdomains |
---|
328 | ! lint = list of inner interfaces |
---|
329 | ! lext = list of extern interfaces |
---|
330 | ! plistin = pointer of sublists of interface nodes |
---|
331 | ! numins = number of interface nodes |
---|
332 | ! listin = list of interface nodes |
---|
333 | ! narea = subdomain number |
---|
334 | ! numnos = number of nodes in the subdomain |
---|
335 | ! v = local field |
---|
336 | |
---|
337 | ! OUTPUTS : w = jump of v |
---|
338 | |
---|
339 | ! WORKSPACE : bufin = contributions of local field on interface |
---|
340 | ! bufout = contributions of outer fields on interface |
---|
341 | ! work |
---|
342 | implicit none |
---|
343 | ! external irecv |
---|
344 | !! integer irecv, work(100) |
---|
345 | !! integer numnes,numins,numnos,narea |
---|
346 | integer numnes,numins,numnos |
---|
347 | integer listnes(numnes,3), plistin(numnes+1), listin(numins) |
---|
348 | REAL(kind=8) v(numnos),w(numins) |
---|
349 | REAL(kind=8) bufin(numins),bufout(numins) |
---|
350 | integer i,j,l |
---|
351 | integer mesg |
---|
352 | |
---|
353 | mesg=128 |
---|
354 | ! receiving the values of the outer fields on the interface |
---|
355 | ! do i=1,numnes |
---|
356 | ! l=(plistin(i+1)-plistin(i))*8 |
---|
357 | ! work(i)=irecv(mesg+listnes(i)-1,bufout(plistin(i)+1),l) |
---|
358 | ! end do |
---|
359 | ! gathering the values of the local field on the interface |
---|
360 | do i=1,numnes |
---|
361 | do j=plistin(i)+1,plistin(i+1) |
---|
362 | bufin(j)=v(listin(j)) |
---|
363 | end do |
---|
364 | end do |
---|
365 | ! sending the values of the local field on the interface |
---|
366 | do i=1,numnes |
---|
367 | l=(plistin(i+1)-plistin(i))*8 |
---|
368 | ! call mppsend(mesg+narea-1,bufin(plistin(i)+1),l,listnes(i)-1,0) |
---|
369 | call mppsend(mesg+listnes(i,3)-1,bufin(plistin(i)+1),l, & |
---|
370 | listnes(i,1)-1,0) |
---|
371 | end do |
---|
372 | ! computing the jump on each interface |
---|
373 | do i=1,numnes |
---|
374 | ! call msgwait(work(i)) |
---|
375 | l=(plistin(i+1)-plistin(i))*8 |
---|
376 | ! call mpprecv(mesg+listnes(i)-1,bufout(plistin(i)+1),l) |
---|
377 | call mpprecv(mesg+listnes(i,2)-1,bufout(plistin(i)+1),l) |
---|
378 | end do |
---|
379 | do j=1,numins |
---|
380 | w(j)=bufin(j)-bufout(j) |
---|
381 | end do |
---|
382 | |
---|
383 | return |
---|
384 | end |
---|
385 | !**********************************************************************c |
---|
386 | subroutine feti_moyen(numnes,listnes, & |
---|
387 | plistin,numins,listin, & |
---|
388 | narea,numnos,weight,v,w,bufin,bufout) |
---|
389 | |
---|
390 | ! averaging a field on the interfaces, in the case where subdomain |
---|
391 | ! number i is allocated to processor number i-1 . |
---|
392 | |
---|
393 | ! INPUTS : numnes = number of neighbouring subdomains |
---|
394 | ! listnes = list of neighbouring subdomains |
---|
395 | ! lint = list of inner interfaces |
---|
396 | ! lext = list of extern interfaces |
---|
397 | ! plistin = pointer of sublists of interface nodes |
---|
398 | ! numins = number of interface nodes |
---|
399 | ! listin = list of interface nodes |
---|
400 | ! narea = subdomain number |
---|
401 | ! numnos = number of nodes in the subdomain |
---|
402 | ! weight = weighting vector |
---|
403 | ! v = local field |
---|
404 | |
---|
405 | ! OUTPUTS : w = averaged field |
---|
406 | |
---|
407 | ! WORKSPACE : bufin = contributions of local field on interface |
---|
408 | ! bufout = contributions of outer fields on interface |
---|
409 | ! work |
---|
410 | implicit none |
---|
411 | ! external irecv |
---|
412 | !! integer irecv, work(100) |
---|
413 | integer numnes,numins,numnos,narea |
---|
414 | integer listnes(numnes,3), plistin(numnes+1), listin(numins) |
---|
415 | REAL(kind=8) weight(numnos),v(numnos),w(numnos) |
---|
416 | REAL(kind=8) bufin(numins),bufout(numins) |
---|
417 | integer i,j,l,n |
---|
418 | integer mesg |
---|
419 | |
---|
420 | mesg=128*2 |
---|
421 | ! initialization |
---|
422 | call feti_vmov(numnos,v,w) |
---|
423 | ! receiving the values of the outer fields on the interface |
---|
424 | ! do i=1,numnes |
---|
425 | ! l=(plistin(i+1)-plistin(i))*8 |
---|
426 | ! work(i)=irecv(mesg+listnes(i)-1,bufout(plistin(i)+1),l) |
---|
427 | ! end do |
---|
428 | ! gathering the values of the local field on the interface |
---|
429 | do i=1,numnes |
---|
430 | do j=plistin(i)+1,plistin(i+1) |
---|
431 | bufin(j)=v(listin(j)) |
---|
432 | end do |
---|
433 | end do |
---|
434 | ! sending the values of the local field on the interface |
---|
435 | do i=1,numnes |
---|
436 | l=(plistin(i+1)-plistin(i))*8 |
---|
437 | ! call mppsend(mesg+narea-1,bufin(plistin(i)+1),l,listnes(i)-1,0) |
---|
438 | call mppsend(mesg+listnes(i,3)-1,bufin(plistin(i)+1),l, & |
---|
439 | listnes(i,1)-1,0) |
---|
440 | end do |
---|
441 | ! computing the jump on each interface |
---|
442 | do i=1,numnes |
---|
443 | ! call msgwait(work(i)) |
---|
444 | l=(plistin(i+1)-plistin(i))*8 |
---|
445 | ! call mpprecv(mesg+listnes(i)-1,bufout(plistin(i)+1),l) |
---|
446 | call mpprecv(mesg+listnes(i,2)-1,bufout(plistin(i)+1),l) |
---|
447 | end do |
---|
448 | do j=1,numins |
---|
449 | n=listin(j) |
---|
450 | w(n)=w(n)+bufout(j) |
---|
451 | end do |
---|
452 | do n=1,numnos |
---|
453 | w(n)=w(n)*weight(n) |
---|
454 | end do |
---|
455 | |
---|
456 | return |
---|
457 | end |
---|
458 | !**********************************************************************c |
---|
459 | subroutine feti_assemb(numnes,listnes, & |
---|
460 | plistin,numins,listin, & |
---|
461 | narea,numnos,v,w,bufin,bufout) |
---|
462 | |
---|
463 | ! assembling a field on the interfaces, in the case where subdomain |
---|
464 | ! number i is allocated to processor number i-1 . |
---|
465 | |
---|
466 | ! INPUTS : numnes = number of neighbouring subdomains |
---|
467 | ! listnes = list of neighbouring subdomains |
---|
468 | ! lint = list of inner interfaces |
---|
469 | ! lext = list of extern interfaces |
---|
470 | ! plistin = pointer of sublists of interface nodes |
---|
471 | ! numins = number of interface nodes |
---|
472 | ! listin = list of interface nodes |
---|
473 | ! narea = subdomain number |
---|
474 | ! numnos = number of nodes in the subdomain |
---|
475 | ! v = unassembled field |
---|
476 | |
---|
477 | ! OUTPUTS : w = assembled field |
---|
478 | |
---|
479 | ! WORKSPACE : bufin = contributions of local field on interface |
---|
480 | ! bufout = contributions of outer fields on interface |
---|
481 | ! work |
---|
482 | implicit none |
---|
483 | ! external irecv |
---|
484 | !! integer irecv, work(100) |
---|
485 | integer numnes,numins,numnos,narea |
---|
486 | integer listnes(numnes,3), plistin(numnes+1), listin(numins) |
---|
487 | REAL(kind=8) v(numnos),w(numnos) |
---|
488 | REAL(kind=8) bufin(numins),bufout(numins) |
---|
489 | integer i,j,l,n |
---|
490 | integer mesg |
---|
491 | |
---|
492 | mesg=128*3 |
---|
493 | ! initialization |
---|
494 | call feti_vmov(numnos,v,w) |
---|
495 | ! receiving the values of the outer fields on the interface |
---|
496 | ! do i=1,numnes |
---|
497 | ! l=(plistin(i+1)-plistin(i))*8 |
---|
498 | ! work(i)=irecv(mesg+listnes(i)-1,bufout(plistin(i)+1),l) |
---|
499 | ! end do |
---|
500 | ! gathering the values of the local field on the interface |
---|
501 | do i=1,numnes |
---|
502 | do j=plistin(i)+1,plistin(i+1) |
---|
503 | bufin(j)=v(listin(j)) |
---|
504 | end do |
---|
505 | end do |
---|
506 | ! sending the values of the local field on the interface |
---|
507 | do i=1,numnes |
---|
508 | l=(plistin(i+1)-plistin(i))*8 |
---|
509 | ! write(90,*)'envoie a ',listnes(i),'de ',l,'octets' |
---|
510 | ! write(90,*)'tag est : ',listnes(i,3) |
---|
511 | ! call mppsend(mesg+narea-1,bufin(plistin(i)+1),l,listnes(i)-1,0) |
---|
512 | call mppsend(mesg+listnes(i,3)-1,bufin(plistin(i)+1),l, & |
---|
513 | listnes(i,1)-1,0) |
---|
514 | end do |
---|
515 | ! computing the jump on each interface |
---|
516 | do i=1,numnes |
---|
517 | ! call msgwait(work(i)) |
---|
518 | l=(plistin(i+1)-plistin(i))*8 |
---|
519 | ! write(90,*)'recoie de ',listnes(i),'de ',l,'octets' |
---|
520 | ! write(90,*)'recois avec tag :',listnes(i,2) |
---|
521 | ! call mpprecv(mesg+listnes(i)-1,bufout(plistin(i)+1),l) |
---|
522 | call mpprecv(mesg+listnes(i,2)-1,bufout(plistin(i)+1),l) |
---|
523 | end do |
---|
524 | do j=1,numins |
---|
525 | n=listin(j) |
---|
526 | w(n)=w(n)+bufout(j) |
---|
527 | end do |
---|
528 | |
---|
529 | return |
---|
530 | end |
---|
531 | !**********************************************************************c |
---|
532 | !! subroutine feti_poids(numnes,listnes,plistin,numins,listin, & |
---|
533 | !! narea,numnos,w,bufin,bufout) |
---|
534 | subroutine feti_poids(numnes, numins,listin, & |
---|
535 | numnos,w ) |
---|
536 | |
---|
537 | ! computation of the subdomain weighting vector, in the case where |
---|
538 | ! subdomain number i is allocated to processor number i-1 . |
---|
539 | |
---|
540 | ! INPUTS : numnes = number of neighbouring subdomains |
---|
541 | ! listnes = list of neighbouring subdomains |
---|
542 | ! plistin = pointer of sublists of interface nodes |
---|
543 | ! numins = number of interface nodes |
---|
544 | ! listin = list of interface nodes |
---|
545 | ! narea = subdomain number |
---|
546 | ! numnos = number of nodes in the subdomain |
---|
547 | |
---|
548 | ! OUTPUTS : w = weighting vector |
---|
549 | |
---|
550 | ! WORKSPACE : bufin = contributions of local field on interface |
---|
551 | ! bufout = contributions of outer fields on interface |
---|
552 | ! work |
---|
553 | implicit none |
---|
554 | !! integer numnes,numins,numnos,narea |
---|
555 | integer numnes,numins,numnos |
---|
556 | !! integer listnes(numnes), plistin(numnes+1), listin(numins) |
---|
557 | integer listin(numins) |
---|
558 | REAL(kind=8) w(numnos) |
---|
559 | !! REAL(kind=8) bufin(numins),bufout(numins) |
---|
560 | integer j,n |
---|
561 | |
---|
562 | ! initialization |
---|
563 | do n=1,numnos |
---|
564 | w(n)=1.e0 |
---|
565 | end do |
---|
566 | ! computing the weight on each interface |
---|
567 | do j=1,numins |
---|
568 | n=listin(j) |
---|
569 | w(n)=w(n)+1.e0 |
---|
570 | end do |
---|
571 | do n=1,numnos |
---|
572 | w(n)=1.e0/w(n) |
---|
573 | end do |
---|
574 | |
---|
575 | return |
---|
576 | end |
---|
577 | !**********************************************************************c |
---|
578 | subroutine feti_probit(numnes,listnes,plistin,numins,listin, & |
---|
579 | numnos,w,bitw) |
---|
580 | |
---|
581 | ! computing the local right-hand-side associated with a field |
---|
582 | ! on the interface . |
---|
583 | |
---|
584 | ! INPUTS : numnes = number of neighbouring subdomains |
---|
585 | ! listnes = list of neighbouring subdomains |
---|
586 | ! plistin = pointer of sublists of interface nodes |
---|
587 | ! numins = number of interface nodes |
---|
588 | ! listin = list of interface nodes |
---|
589 | ! numnos = number of nodes in the subdomain |
---|
590 | ! w = interface field |
---|
591 | |
---|
592 | ! OUTPUTS : bitw = local right-hand side |
---|
593 | |
---|
594 | implicit none |
---|
595 | integer numnes,numins,numnos |
---|
596 | integer listnes(numnes,3), plistin(numnes+1), listin(numins) |
---|
597 | REAL(kind=8) w(numins) |
---|
598 | REAL(kind=8) bitw(numnos) |
---|
599 | integer i,j,n |
---|
600 | |
---|
601 | ! initialization |
---|
602 | do n=1,numnos |
---|
603 | bitw(n)=0.e0 |
---|
604 | end do |
---|
605 | ! assembling the values on the interfaces |
---|
606 | do i=1,numnes |
---|
607 | do j=plistin(i)+1,plistin(i+1) |
---|
608 | n=listin(j) |
---|
609 | bitw(n)=bitw(n)+w(j) |
---|
610 | end do |
---|
611 | end do |
---|
612 | |
---|
613 | return |
---|
614 | end |
---|
615 | !**********************************************************************c |
---|
616 | subroutine feti_inisub(ni,nj,ibondi,ibondj,iperio, & |
---|
617 | ibsw,ibnw,ibse,ibne, & |
---|
618 | ninterf,ninterfc,nni,nnic) |
---|
619 | implicit none |
---|
620 | integer ni,nj,ibondi,ibondj,iperio |
---|
621 | integer ibsw,ibnw,ibse,ibne |
---|
622 | integer ibsw2,ibnw2,ibse2,ibne2 |
---|
623 | integer ninterf,ninterfc,nni,nnic |
---|
624 | ! integer i,j |
---|
625 | |
---|
626 | ! determination de la position du sous-domaine dans la grille |
---|
627 | ! cela depnd de ibondi et ibondj |
---|
628 | |
---|
629 | ! nombre de sous-domaines voisins et de la longueur du tableau des |
---|
630 | ! noeuds voisins |
---|
631 | |
---|
632 | ! initialisation |
---|
633 | |
---|
634 | ninterf=0 |
---|
635 | nni=0 |
---|
636 | |
---|
637 | !...first, the "segment-points" |
---|
638 | |
---|
639 | if((ibondi == 2).AND.(iperio /= 1))then |
---|
640 | |
---|
641 | ! la periodicite : pas d'interface east / west si nbondi = 2 |
---|
642 | ! ET nperio != 1 |
---|
643 | |
---|
644 | ! elseif(iperio == 1) then |
---|
645 | |
---|
646 | !...ca doit marcher!? |
---|
647 | |
---|
648 | else |
---|
649 | |
---|
650 | ! west |
---|
651 | |
---|
652 | if(ibondi /= -1)then |
---|
653 | ninterf=ninterf+1 |
---|
654 | nni=nni+nj+1 |
---|
655 | endif |
---|
656 | |
---|
657 | ! east |
---|
658 | |
---|
659 | if(ibondi /= 1)then |
---|
660 | ninterf=ninterf+1 |
---|
661 | nni=nni+nj+1 |
---|
662 | endif |
---|
663 | |
---|
664 | endif |
---|
665 | |
---|
666 | ! south |
---|
667 | |
---|
668 | if(ibondj /= -1.AND.ibondj /= 2)then |
---|
669 | ninterf=ninterf+1 |
---|
670 | nni=nni+ni+1 |
---|
671 | endif |
---|
672 | |
---|
673 | ! north |
---|
674 | |
---|
675 | if(ibondj /= 1.AND.ibondj /= 2)then |
---|
676 | ninterf=ninterf+1 |
---|
677 | nni=nni+ni+1 |
---|
678 | endif |
---|
679 | |
---|
680 | nnic=nni |
---|
681 | ninterfc=ninterf |
---|
682 | |
---|
683 | |
---|
684 | !...second, the "corner-points" |
---|
685 | |
---|
686 | !...determination of the neighboorings, they depends on iperio |
---|
687 | ! boundary condition => local description : ipXX becomes ipXX2... |
---|
688 | |
---|
689 | ibnw2 = ibnw |
---|
690 | ibne2 = ibne |
---|
691 | ibsw2 = ibsw |
---|
692 | ibse2 = ibse |
---|
693 | |
---|
694 | !...iperio boundary condition effect |
---|
695 | |
---|
696 | if(iperio == 1) then |
---|
697 | |
---|
698 | !...general case : Earth == infinite tube |
---|
699 | |
---|
700 | ibnw2 = 1 |
---|
701 | ibne2 = 1 |
---|
702 | ibsw2 = 1 |
---|
703 | ibse2 = 1 |
---|
704 | |
---|
705 | !...real boundary condition |
---|
706 | |
---|
707 | if(ibondj == -1.OR.ibondj == 2) then |
---|
708 | ibsw2 = 0 |
---|
709 | ibse2 = 0 |
---|
710 | endif |
---|
711 | |
---|
712 | if(ibondj == 1.OR.ibondj == 2) then |
---|
713 | ibnw2 = 0 |
---|
714 | ibne2 = 0 |
---|
715 | endif |
---|
716 | |
---|
717 | endif |
---|
718 | |
---|
719 | ! la periodicite : pas de coin si nbondi = 2 |
---|
720 | ! ET nperio != 1 |
---|
721 | |
---|
722 | ! north-west |
---|
723 | |
---|
724 | if(ibnw == 1)then |
---|
725 | ninterfc=ninterfc+1 |
---|
726 | nnic=nnic+1 |
---|
727 | endif |
---|
728 | |
---|
729 | ! north-east |
---|
730 | |
---|
731 | if(ibne == 1)then |
---|
732 | ninterfc=ninterfc+1 |
---|
733 | nnic=nnic+1 |
---|
734 | endif |
---|
735 | |
---|
736 | ! south-west |
---|
737 | |
---|
738 | if(ibsw == 1)then |
---|
739 | ninterfc=ninterfc+1 |
---|
740 | nnic=nnic+1 |
---|
741 | endif |
---|
742 | |
---|
743 | ! south-east |
---|
744 | |
---|
745 | if(ibse == 1)then |
---|
746 | ninterfc=ninterfc+1 |
---|
747 | nnic=nnic+1 |
---|
748 | endif |
---|
749 | |
---|
750 | ! write(0,*)' nbre de voisins =',ninterf |
---|
751 | ! write(0,*)' longueur des interfaces =',nni |
---|
752 | ! write(0,*)' longueur des interfaces =',nnic |
---|
753 | |
---|
754 | return |
---|
755 | end |
---|
756 | !**********************************************************************c |
---|
757 | subroutine feti_subound(ni,nj,ildi,ilei,ildj,ilej, & |
---|
758 | imoi,ibondi,ibondj,iperio, & |
---|
759 | ninterf,ninterfc, & |
---|
760 | iowe,ioea,ioso,iono, & |
---|
761 | ibsw,ibnw,ibse,ibne, & |
---|
762 | ipsw,ipnw,ipse,ipne, & |
---|
763 | nsdvois,nsdvoisc, & |
---|
764 | plistin,nni,listin) |
---|
765 | |
---|
766 | implicit none |
---|
767 | integer ni,nj,ildi,ilei,ildj,ilej |
---|
768 | integer imoi,ibondi,ibondj,iperio |
---|
769 | integer iowe,ioea,ioso,iono |
---|
770 | integer ibsw,ibnw,ibse,ibne |
---|
771 | integer ipsw,ipnw,ipse,ipne |
---|
772 | integer ninterf,ninterfc,nni |
---|
773 | integer nsdvois(ninterf,3),nsdvoisc(ninterfc,3) |
---|
774 | integer plistin(ninterfc+1),listin(nni) |
---|
775 | integer nint,nni0,ii,jj |
---|
776 | ! integer i,j,nint,nni0,ii,jj |
---|
777 | |
---|
778 | external mynode |
---|
779 | integer mynode |
---|
780 | |
---|
781 | ! determination de la position du sous-domaine dans la grille |
---|
782 | ! cela depnd de ibondi et ibondj |
---|
783 | |
---|
784 | ! liste des sous-domaines voisins et des noeuds interface |
---|
785 | |
---|
786 | ! be carefull!!! in the FETI algorithm, area number is considered |
---|
787 | ! instead of process number |
---|
788 | |
---|
789 | ! initialisation |
---|
790 | |
---|
791 | nint=0 |
---|
792 | nni0=0 |
---|
793 | plistin(1)=0 |
---|
794 | |
---|
795 | !...first, the "segment-points" |
---|
796 | |
---|
797 | if((ibondi == 2).AND.(iperio /= 1))then |
---|
798 | |
---|
799 | ! la periodicite : pas d'interface east / west si nbondi = 2 |
---|
800 | ! ET nperio != 1 |
---|
801 | |
---|
802 | elseif(iperio == 1) then |
---|
803 | |
---|
804 | ! west |
---|
805 | |
---|
806 | nint=nint+1 |
---|
807 | nsdvoisc(nint,1)=imoi |
---|
808 | do jj=1,nj+1 |
---|
809 | listin(nni0+jj)=(jj-1)*(ni+1)+1 |
---|
810 | end do |
---|
811 | nni0=nni0+nj+1 |
---|
812 | plistin(nint+1)=nni0 |
---|
813 | nsdvoisc(nint,2)=1 |
---|
814 | nsdvoisc(nint,3)=2 |
---|
815 | |
---|
816 | ! east |
---|
817 | |
---|
818 | nint=nint+1 |
---|
819 | nsdvoisc(nint,1)=imoi |
---|
820 | do jj=1,nj+1 |
---|
821 | listin(nni0+jj)=(jj-1)*(ni+1)+ilei |
---|
822 | end do |
---|
823 | nni0=nni0+nj+1 |
---|
824 | plistin(nint+1)=nni0 |
---|
825 | nsdvoisc(nint,2)=2 |
---|
826 | nsdvoisc(nint,3)=1 |
---|
827 | |
---|
828 | else |
---|
829 | |
---|
830 | ! west |
---|
831 | |
---|
832 | if(ibondi /= -1)then |
---|
833 | nint=nint+1 |
---|
834 | nsdvoisc(nint,1)=(iowe+1) |
---|
835 | do jj=1,nj+1 |
---|
836 | listin(nni0+jj)=(jj-1)*(ni+1)+1 |
---|
837 | end do |
---|
838 | nni0=nni0+nj+1 |
---|
839 | plistin(nint+1)=nni0 |
---|
840 | nsdvoisc(nint,2)=1 |
---|
841 | nsdvoisc(nint,3)=2 |
---|
842 | endif |
---|
843 | |
---|
844 | ! east |
---|
845 | |
---|
846 | if(ibondi /= 1)then |
---|
847 | nint=nint+1 |
---|
848 | nsdvoisc(nint,1)=(ioea+1) |
---|
849 | do jj=1,nj+1 |
---|
850 | listin(nni0+jj)=(jj-1)*(ni+1)+ilei |
---|
851 | end do |
---|
852 | nni0=nni0+nj+1 |
---|
853 | plistin(nint+1)=nni0 |
---|
854 | nsdvoisc(nint,2)=2 |
---|
855 | nsdvoisc(nint,3)=1 |
---|
856 | endif |
---|
857 | |
---|
858 | endif |
---|
859 | |
---|
860 | ! south |
---|
861 | |
---|
862 | if(ibondj /= -1.AND.ibondj /= 2)then |
---|
863 | nint=nint+1 |
---|
864 | nsdvoisc(nint,1)=(ioso+1) |
---|
865 | do ii=1,ni+1 |
---|
866 | listin(nni0+ii)=ii |
---|
867 | end do |
---|
868 | nni0=nni0+ni+1 |
---|
869 | plistin(nint+1)=nni0 |
---|
870 | nsdvoisc(nint,2)=3 |
---|
871 | nsdvoisc(nint,3)=4 |
---|
872 | endif |
---|
873 | |
---|
874 | ! north |
---|
875 | |
---|
876 | if(ibondj /= 1.AND.ibondj /= 2)then |
---|
877 | nint=nint+1 |
---|
878 | nsdvoisc(nint,1)=(iono+1) |
---|
879 | do ii=1,ni+1 |
---|
880 | listin(nni0+ii)=(ilej-1)*(ni+1)+ii |
---|
881 | end do |
---|
882 | nni0=nni0+ni+1 |
---|
883 | plistin(nint+1)=nni0 |
---|
884 | nsdvoisc(nint,2)=4 |
---|
885 | nsdvoisc(nint,3)=3 |
---|
886 | endif |
---|
887 | |
---|
888 | !...second, the "corner-points" |
---|
889 | |
---|
890 | !...determination of the neighboorings, they depends on iperio |
---|
891 | ! boundary condition => local description : inimpp or inimpp2... |
---|
892 | |
---|
893 | ! la periodicite : pas de coin si nbondi = 2 |
---|
894 | ! ET nperio != 1 |
---|
895 | |
---|
896 | !!!!!!!!!verifier la definition des voisins-coin |
---|
897 | |
---|
898 | ! north-west |
---|
899 | |
---|
900 | if(ibnw == 1)then |
---|
901 | nint=nint+1 |
---|
902 | nsdvoisc(nint,1)=(ipnw+1) |
---|
903 | listin(nni0+1)=(ilej-1)*(ni+1)+1 |
---|
904 | nni0=nni0+1 |
---|
905 | plistin(nint+1)=nni0 |
---|
906 | nsdvoisc(nint,2)=5 |
---|
907 | nsdvoisc(nint,3)=8 |
---|
908 | endif |
---|
909 | |
---|
910 | ! north-east |
---|
911 | |
---|
912 | if(ibne == 1)then |
---|
913 | nint=nint+1 |
---|
914 | nsdvoisc(nint,1)=(ipne+1) |
---|
915 | listin(nni0+1)=(ilej-1)*(ni+1)+ilei |
---|
916 | nni0=nni0+1 |
---|
917 | plistin(nint+1)=nni0 |
---|
918 | nsdvoisc(nint,2)=6 |
---|
919 | nsdvoisc(nint,3)=7 |
---|
920 | endif |
---|
921 | |
---|
922 | ! south-west |
---|
923 | |
---|
924 | if(ibsw == 1)then |
---|
925 | nint=nint+1 |
---|
926 | nsdvoisc(nint,1)=(ipsw+1) |
---|
927 | listin(nni0+1)=1 |
---|
928 | nni0=nni0+1 |
---|
929 | plistin(nint+1)=nni0 |
---|
930 | nsdvoisc(nint,2)=7 |
---|
931 | nsdvoisc(nint,3)=6 |
---|
932 | endif |
---|
933 | |
---|
934 | ! south-east |
---|
935 | |
---|
936 | if(ibse == 1)then |
---|
937 | nint=nint+1 |
---|
938 | nsdvoisc(nint,1)=(ipse+1) |
---|
939 | listin(nni0+1)=ilei |
---|
940 | nni0=nni0+1 |
---|
941 | plistin(nint+1)=nni0 |
---|
942 | nsdvoisc(nint,2)=8 |
---|
943 | nsdvoisc(nint,3)=5 |
---|
944 | endif |
---|
945 | |
---|
946 | !...je compte les coins |
---|
947 | |
---|
948 | do ii=1,ninterf |
---|
949 | nsdvois(ii,1) = nsdvoisc(ii,1) |
---|
950 | nsdvois(ii,2) = nsdvoisc(ii,2) |
---|
951 | nsdvois(ii,3) = nsdvoisc(ii,3) |
---|
952 | enddo |
---|
953 | |
---|
954 | |
---|
955 | ! impressions |
---|
956 | |
---|
957 | ! write(*,*)' sous-domaine :',imoi |
---|
958 | ! write(*,*)' nombre de sous-domaines voisins :',ninterfc |
---|
959 | ! do ii=1,ninterfc |
---|
960 | ! write(*,*)' sous-domaine voisin :',nsdvoisc(ii,1) |
---|
961 | ! write(*,*)' nombre de noeuds interface :', & |
---|
962 | ! plistin(ii+1)-plistin(ii) |
---|
963 | ! write(*,*)' noeuds interface :',(listin(jj),jj=plistin(ii)+1, & |
---|
964 | ! plistin(ii+1)) |
---|
965 | ! end do |
---|
966 | |
---|
967 | return |
---|
968 | end |
---|
969 | !**********************************************************************c |
---|
970 | subroutine feti_subdir(ni,nj,noeuds,ndir,mgcnum) |
---|
971 | implicit none |
---|
972 | integer ni,nj,noeuds,ndir |
---|
973 | integer mgcnum(ni+1,nj+1) |
---|
974 | integer i,j |
---|
975 | ! integer i,j,k |
---|
976 | |
---|
977 | ! determination de la position du sous-domaine dans la grille |
---|
978 | ! narea = (j-1)*iglo + i , avec 1<=i<=iglo et 1<=j<=jglo |
---|
979 | |
---|
980 | |
---|
981 | |
---|
982 | ! 1. NUMBER THE GRID-POINTS USING BMASK |
---|
983 | ! ------------------------------------- |
---|
984 | |
---|
985 | |
---|
986 | ndir=0 |
---|
987 | do j=1,nj+1 |
---|
988 | do i=1,ni+1 |
---|
989 | ndir=ndir+1-mgcnum(i,j) |
---|
990 | ! do k=1,noeuds |
---|
991 | ! ndir=ndir+1-mgcnum(k,1) |
---|
992 | enddo |
---|
993 | enddo |
---|
994 | |
---|
995 | return |
---|
996 | end |
---|
997 | !**********************************************************************c |
---|
998 | subroutine feti_listdir(ni,nj,logdir,ndir,lisdir) |
---|
999 | implicit none |
---|
1000 | integer ni,nj,ndir |
---|
1001 | integer logdir(ni,nj),lisdir(ndir) |
---|
1002 | integer i0,ji,jj |
---|
1003 | |
---|
1004 | ! creation de la liste des degres de liberte bloques |
---|
1005 | i0=0 |
---|
1006 | do ji=1,ni |
---|
1007 | do jj=1,nj |
---|
1008 | if(logdir(ji,jj) == 0) then |
---|
1009 | i0=i0+1 |
---|
1010 | lisdir(i0)=ji+(jj-1)*ni |
---|
1011 | endif |
---|
1012 | end do |
---|
1013 | end do |
---|
1014 | if(i0 /= ndir) then |
---|
1015 | ! write(*,*)' nombre de ddl bloques/prevus :',i0,ndir |
---|
1016 | ! write(*,*)'liste des ddl Dirichlet :',(lisdir(ji),ji=1,ndir) |
---|
1017 | stop |
---|
1018 | endif |
---|
1019 | |
---|
1020 | return |
---|
1021 | end |
---|
1022 | !*********************************************************************** |
---|
1023 | subroutine feti_blodir(n,x,ndlblo,list) |
---|
1024 | implicit none |
---|
1025 | integer n,ndlblo |
---|
1026 | integer list(ndlblo) |
---|
1027 | REAL(kind=8) x(n) |
---|
1028 | integer i |
---|
1029 | |
---|
1030 | ! remise a zero des ddl bloques |
---|
1031 | do i=1,ndlblo |
---|
1032 | x(list(i))=0.e0 |
---|
1033 | end do |
---|
1034 | |
---|
1035 | return |
---|
1036 | end |
---|
1037 | !*********************************************************************** |
---|
1038 | subroutine feti_numblo(ndlblo,lisblo) |
---|
1039 | implicit none |
---|
1040 | external mynode |
---|
1041 | !! integer mynode,indlblog,iibidon |
---|
1042 | integer mynode |
---|
1043 | integer ndlblo |
---|
1044 | integer lisblo(ndlblo) |
---|
1045 | !! integer n |
---|
1046 | |
---|
1047 | ! indlblog=-ndlblo |
---|
1048 | ! call mpp_min(indlblog,1,iibidon) |
---|
1049 | ! indlblog=-indlblog |
---|
1050 | ! if(mynode() == 0)THEN |
---|
1051 | ! write(0,*)' nombre de degres de liberte flottants :',ndlblo |
---|
1052 | ! write(0,*)' liste :',(lisblo(n),n=1,ndlblo) |
---|
1053 | ! endif |
---|
1054 | |
---|
1055 | return |
---|
1056 | end |
---|
1057 | !*********************************************************************** |
---|
1058 | subroutine feti_blomat(ni,nj,a,ndlblo,lisblo) |
---|
1059 | implicit none |
---|
1060 | integer ni,nj |
---|
1061 | REAL(kind=8) a(ni,nj,5) |
---|
1062 | integer ndlblo |
---|
1063 | integer lisblo(ndlblo) |
---|
1064 | !! integer i,j,k,l,n |
---|
1065 | integer i,j,k,n |
---|
1066 | |
---|
1067 | do n=1,ndlblo |
---|
1068 | ! degre de liberte bloque : lisblo(n) = (j-1)*ni + i |
---|
1069 | i=mod(lisblo(n)-1,ni)+1 |
---|
1070 | j=((lisblo(n)-1)/ni)+1 |
---|
1071 | ! annulation des coefficients sur la ligne |
---|
1072 | do k=1,5 |
---|
1073 | a(i,j,k)=0.e0 |
---|
1074 | end do |
---|
1075 | a(i,j,3)=1.e0 |
---|
1076 | ! annulation des coefficients dans les colonnes correspondantes |
---|
1077 | if(i > 1) then |
---|
1078 | a(i-1,j,4)=0.e0 |
---|
1079 | endif |
---|
1080 | if(i < ni) then |
---|
1081 | a(i+1,j,2)=0.e0 |
---|
1082 | endif |
---|
1083 | if(j > 1) then |
---|
1084 | a(i,j-1,5)=0.e0 |
---|
1085 | endif |
---|
1086 | if(j < nj) then |
---|
1087 | a(i,j+1,1)=0.e0 |
---|
1088 | endif |
---|
1089 | end do |
---|
1090 | |
---|
1091 | return |
---|
1092 | end |
---|
1093 | !*********************************************************************** |
---|
1094 | subroutine feti_blomat1(ni,nj,a,ndlblo,lisblo,nsp) |
---|
1095 | implicit none |
---|
1096 | integer ni,nj,ndlblo |
---|
1097 | REAL(kind=8) a(ni,nj,5),nsp(ni,nj,ndlblo) |
---|
1098 | integer lisblo(ndlblo) |
---|
1099 | !! integer i,j,k,l,n |
---|
1100 | integer i,j,n |
---|
1101 | |
---|
1102 | do n=1,ndlblo |
---|
1103 | ! degre de liberte bloque : lisblo(n) = (j-1)*ni + i |
---|
1104 | i=mod(lisblo(n)-1,ni)+1 |
---|
1105 | j=((lisblo(n)-1)/ni)+1 |
---|
1106 | ! annulation des coefficients dans les colonnes correspondantes |
---|
1107 | if(i > 1) then |
---|
1108 | nsp(i-1,j,n)=-a(i-1,j,4) |
---|
1109 | endif |
---|
1110 | if(i < ni) then |
---|
1111 | nsp(i+1,j,n)=-a(i+1,j,2) |
---|
1112 | endif |
---|
1113 | if(j > 1) then |
---|
1114 | nsp(i,j-1,n)=-a(i,j-1,5) |
---|
1115 | endif |
---|
1116 | if(j < nj) then |
---|
1117 | nsp(i,j+1,n)=-a(i,j+1,1) |
---|
1118 | endif |
---|
1119 | end do |
---|
1120 | |
---|
1121 | return |
---|
1122 | end |
---|
1123 | !*********************************************************************** |
---|
1124 | subroutine feti_nullsp(noeuds,ni,nj,lpblo,blo,a,ndlblo,lisblo,nsp,z) |
---|
1125 | ! calcul du noyau de la matrice |
---|
1126 | implicit none |
---|
1127 | integer noeuds,ni,nj,lpblo,ndlblo |
---|
1128 | REAL(kind=8) a(ni,nj,5) |
---|
1129 | integer lisblo(ndlblo) |
---|
1130 | REAL(kind=8) blo(lpblo),nsp(noeuds,ndlblo),z(noeuds) |
---|
1131 | external sdot |
---|
1132 | REAL(kind=8) sdot,res |
---|
1133 | integer j,i |
---|
1134 | |
---|
1135 | ! calcul de [Aii]-1 * Aib |
---|
1136 | do j=1,ndlblo |
---|
1137 | call feti_resloc(noeuds,ni,nj,a,lpblo,blo,nsp(1,j),nsp(1,j),z) |
---|
1138 | end do |
---|
1139 | ! remise a identite du bloc diagonal correspondant aux modes bloques |
---|
1140 | do i=1,ndlblo |
---|
1141 | do j=1,ndlblo |
---|
1142 | nsp(lisblo(i),j)=0.e0 |
---|
1143 | end do |
---|
1144 | nsp(lisblo(i),i)=1.e0 |
---|
1145 | end do |
---|
1146 | ! verifications |
---|
1147 | do j=1,ndlblo |
---|
1148 | call feti_proax(noeuds,ni,nj,a,nsp(1,j),z) |
---|
1149 | res=sqrt(sdot(noeuds,z,1,z,1)/ & |
---|
1150 | sdot(noeuds,nsp(1,j),1,nsp(1,j),1)) |
---|
1151 | ! write(0,*)' residu de la colonne ',j,' du noyau :',res |
---|
1152 | ! write(0,*)' vecteur colonne :' |
---|
1153 | ! call prvec(noeuds,nsp(1,j)) |
---|
1154 | end do |
---|
1155 | |
---|
1156 | return |
---|
1157 | end |
---|
1158 | !**********************************************************************c |
---|
1159 | subroutine feti_project(gint,pgint,nspdim,x,b,nmspdim,numit0 & |
---|
1160 | ,nitmax,gvec,agvec,d,ad,add,gamm,numnes,listnes,plistin,numins & |
---|
1161 | ,listin,narea,numnos,nsp,v,w,bufin,bufout,work) |
---|
1162 | |
---|
1163 | implicit none |
---|
1164 | |
---|
1165 | integer nitmax, nspdim |
---|
1166 | REAL(kind=8) x(nspdim), b(nspdim) |
---|
1167 | |
---|
1168 | !...the Krylov optimisation used in PCG associated to the Coarse Solver |
---|
1169 | ! involve to keep di & adi vectors in etesolve procedure |
---|
1170 | ! the nmspdim parameter is introduced to avoid message "out of bound" |
---|
1171 | ! arrising qhen the memory is checked! |
---|
1172 | |
---|
1173 | ! nmspdim = Max (nspdim(Ep)) <= memory management (f77 context) |
---|
1174 | ! p=1,Np |
---|
1175 | |
---|
1176 | ! REAL(kind=8) gvec(nspdim), agvec(nspdim), |
---|
1177 | ! & d(nspdim,0:nitmax-1), ad(nspdim,0:nitmax-1), |
---|
1178 | ! & add(0:nitmax-1), gamm(0:nitmax-1) |
---|
1179 | |
---|
1180 | integer nmspdim |
---|
1181 | REAL(kind=8) gvec(nspdim), agvec(nspdim), & |
---|
1182 | d(nmspdim,0:nitmax-1), ad(nmspdim,0:nitmax-1), & |
---|
1183 | add(0:nitmax-1), gamm(0:nitmax-1) |
---|
1184 | |
---|
1185 | integer numnes,numins,narea, numnos, numit0 |
---|
1186 | integer listnes(numnes,3), plistin(numnes+1), listin(numins) |
---|
1187 | REAL(kind=8) nsp(numnos,nspdim), v(numnos), & |
---|
1188 | w(numins), gint(numins), pgint(numins) |
---|
1189 | REAL(kind=8) bufin(numins), bufout(numins) |
---|
1190 | REAL(kind=8) work(0:nitmax-1) |
---|
1191 | |
---|
1192 | |
---|
1193 | call feti_proet(numnes,listnes,plistin,numins,listin,numnos, & |
---|
1194 | gint,nspdim,nsp,b,v) |
---|
1195 | call feti_etesolv(nspdim,x,b,nmspdim,numit0,nitmax,gvec,agvec,d,ad & |
---|
1196 | ,add,gamm,numnes,listnes,plistin,numins,listin,narea,numnos & |
---|
1197 | ,nsp,v,w,bufin,bufout,work) |
---|
1198 | call feti_proe(numnes,listnes,plistin,numins,listin,narea,numnos, & |
---|
1199 | nspdim,nsp,x,v,w,bufin,bufout) |
---|
1200 | call feti_vsub(numins,gint,w,pgint) |
---|
1201 | |
---|
1202 | return |
---|
1203 | end |
---|
1204 | !**********************************************************************c |
---|
1205 | subroutine feti_etesolv(nspdim,x,b,nmspdim,numit0,nitmax,gvec & |
---|
1206 | ,agvec,d,ad,add,gamm,numnes,listnes,plistin,numins,listin & |
---|
1207 | ,narea,numnos,nsp,v,w,bufin,bufout,work) |
---|
1208 | |
---|
1209 | ! INPUTS : nspdim = dimension of null space : LOCAL |
---|
1210 | ! x = initial solution |
---|
1211 | ! b = right hand side |
---|
1212 | ! nmspdim = Max (nspdim(Ep)) <= memory management |
---|
1213 | ! p=1,Np |
---|
1214 | ! numit0 = number of previous direction vectors, at most |
---|
1215 | ! nitmax |
---|
1216 | ! nitmax = number of maximum direction vectors for |
---|
1217 | ! reconjugation |
---|
1218 | ! gvec = gradient vector |
---|
1219 | ! agvec = A * g |
---|
1220 | ! d = direction vector |
---|
1221 | ! ad = A * d |
---|
1222 | ! add = dot products (Ad,d) |
---|
1223 | ! gamm = GAMMA coefficients |
---|
1224 | ! numnes = number of neighbouring subdomains |
---|
1225 | ! listnes = list of neighbouring subdomains |
---|
1226 | ! plistin = pointer of sublists of interface nodes |
---|
1227 | ! numins = number of interface nodes |
---|
1228 | ! listin = list of interface nodes |
---|
1229 | ! narea = subdomain number |
---|
1230 | ! numnos = number of nodes in the subdomain |
---|
1231 | |
---|
1232 | ! WORKSPACE : bufin = contributions of local fields on interface |
---|
1233 | ! bufout = contributions of outer fields on interface |
---|
1234 | ! work = workspace for mpp_sum calls |
---|
1235 | |
---|
1236 | implicit none |
---|
1237 | |
---|
1238 | integer nitmax, nspdim |
---|
1239 | REAL(kind=8) x(nspdim), b(nspdim) |
---|
1240 | |
---|
1241 | !...the Krylov optimisation involve to keep di & adi vectors |
---|
1242 | ! the nmspdim parameter is introduce to avoid message "out of bound" |
---|
1243 | ! arrising qhen the memory is checked! |
---|
1244 | |
---|
1245 | ! REAL(kind=8) gvec(nspdim), agvec(nspdim), |
---|
1246 | ! & d(nspdim,0:nitmax-1), ad(nspdim,0:nitmax-1), |
---|
1247 | ! & rho, temp1(2), temp2(2), add(0:nitmax-1), |
---|
1248 | ! & gamm(0:nitmax-1), facg0, facgn, facst, eps |
---|
1249 | |
---|
1250 | integer nmspdim |
---|
1251 | REAL(kind=8) gvec(nmspdim), agvec(nmspdim), & |
---|
1252 | d(nspdim,0:nitmax-1), ad(nspdim,0:nitmax-1), & |
---|
1253 | rho, temp1(2), temp2(2), add(0:nitmax-1), & |
---|
1254 | gamm(0:nitmax-1), facg0, facgn, facst, eps |
---|
1255 | |
---|
1256 | integer numnes,numins,narea, numnos |
---|
1257 | integer listnes(numnes,3), plistin(numnes+1), listin(numins) |
---|
1258 | REAL(kind=8) nsp(numnos,nspdim), v(numnos), w(numins) |
---|
1259 | REAL(kind=8) bufin(numins), bufout(numins) |
---|
1260 | REAL(kind=8) work(0:nitmax-1) |
---|
1261 | |
---|
1262 | integer i, j, n, numit0, nn, n1, numit1 |
---|
1263 | |
---|
1264 | external sdot |
---|
1265 | REAL(kind=8) sdot |
---|
1266 | |
---|
1267 | ! initialisation |
---|
1268 | |
---|
1269 | do 1 i=1,nspdim |
---|
1270 | x(i) = 0.e0 |
---|
1271 | gvec(i) = -b(i) |
---|
1272 | 1 continue |
---|
1273 | facg0 = sdot( nspdim, gvec, 1, gvec , 1) |
---|
1274 | call mpp_sum(facg0,1,work) |
---|
1275 | |
---|
1276 | if( facg0 == 0.e0 ) return |
---|
1277 | |
---|
1278 | eps=1.e-24 |
---|
1279 | facst=facg0*eps |
---|
1280 | |
---|
1281 | ! initial projection |
---|
1282 | |
---|
1283 | if(numit0 > 0) then |
---|
1284 | do 91 j = 0 , numit0 |
---|
1285 | gamm(j) = - sdot(nspdim,gvec,1,d(1,j),1) / add(j) |
---|
1286 | 91 continue |
---|
1287 | call mpp_sum(gamm,numit0+1,work) |
---|
1288 | call feti_mxvadd( nspdim, numit0+1, d, gamm, x ) |
---|
1289 | call feti_mxvadd( nspdim, numit0+1, ad, gamm, gvec ) |
---|
1290 | |
---|
1291 | ! test residual after initial projection |
---|
1292 | |
---|
1293 | facgn = sdot( nspdim, gvec, 1, gvec, 1 ) |
---|
1294 | call mpp_sum(facgn,1,work) |
---|
1295 | if( facgn <= facst ) then |
---|
1296 | ! write(0,*) 'residual after initial projection for feti_etesolv=' & |
---|
1297 | ! ,facgn/facg0 |
---|
1298 | return |
---|
1299 | endif |
---|
1300 | |
---|
1301 | endif |
---|
1302 | |
---|
1303 | ! compute (Et E) * gvec |
---|
1304 | |
---|
1305 | call feti_proe(numnes,listnes,plistin,numins,listin,narea,numnos, & |
---|
1306 | nspdim,nsp,gvec,v,w,bufin,bufout) |
---|
1307 | call feti_proet(numnes,listnes,plistin,numins,listin,numnos,w, & |
---|
1308 | nspdim,nsp,agvec,v) |
---|
1309 | |
---|
1310 | ! conjugate d . |
---|
1311 | |
---|
1312 | if(numit0 > 0) then |
---|
1313 | numit1=min0(numit0+1,nitmax-1) |
---|
1314 | do 94 i=0,numit0 |
---|
1315 | gamm(i) = -sdot(nspdim,gvec,1,ad(1,i),1)/ add(i) |
---|
1316 | 94 continue |
---|
1317 | call mpp_sum(gamm,numit0+1,work) |
---|
1318 | |
---|
1319 | call feti_osaxpy( nspdim,gamm(numit0),d(1,numit0), & |
---|
1320 | gvec,d(1,numit1)) |
---|
1321 | call feti_osaxpy( nspdim,gamm(numit0),ad(1,numit0), & |
---|
1322 | agvec,ad(1,numit1)) |
---|
1323 | |
---|
1324 | call feti_mxvadd( nspdim, numit1-1, d, gamm, d(1,numit1) ) |
---|
1325 | call feti_mxvadd( nspdim, numit1-1, ad, gamm, ad(1,numit1) ) |
---|
1326 | else |
---|
1327 | numit1=0 |
---|
1328 | |
---|
1329 | do 20 i = 1, nspdim |
---|
1330 | d(i,numit1) = gvec(i) |
---|
1331 | 20 continue |
---|
1332 | |
---|
1333 | do 30 i = 1, nspdim |
---|
1334 | ad(i,numit1) = agvec(i) |
---|
1335 | 30 continue |
---|
1336 | endif |
---|
1337 | |
---|
1338 | ! computing rho |
---|
1339 | |
---|
1340 | do 50 nn = numit1, nitmax |
---|
1341 | n=min0(nn,nitmax-1) |
---|
1342 | temp1(1) = sdot( nspdim, gvec, 1, d(1,n), 1 ) |
---|
1343 | temp1(2) = sdot( nspdim, ad(1,n), 1, d(1,n), 1 ) |
---|
1344 | call mpp_sum(temp1,2,temp2) |
---|
1345 | add(n) = temp1(2) |
---|
1346 | rho = - temp1(1) / add(n) |
---|
1347 | |
---|
1348 | call feti_osaxpy( nspdim, rho, d(1,n), x, x ) |
---|
1349 | call feti_osaxpy( nspdim, rho, ad(1,n), gvec, gvec ) |
---|
1350 | |
---|
1351 | ! test residual . |
---|
1352 | |
---|
1353 | facgn = sdot( nspdim, gvec, 1, gvec, 1 ) |
---|
1354 | call mpp_sum(facgn,1,work) |
---|
1355 | ! write(0,*) 'iteration =', nn+1, 'residual = ', facgn/facg0 |
---|
1356 | if( facgn <= facst ) then |
---|
1357 | goto 999 |
---|
1358 | endif |
---|
1359 | |
---|
1360 | ! reconjugation of gvec |
---|
1361 | |
---|
1362 | do 60 i=0,n |
---|
1363 | gamm(i) = -sdot(nspdim,gvec,1,ad(1,i),1)/ add(i) |
---|
1364 | 60 continue |
---|
1365 | call mpp_sum(gamm,n+1,work) |
---|
1366 | n1=min0(nn+1,nitmax-1) |
---|
1367 | |
---|
1368 | call feti_osaxpy( nspdim,gamm(n),d(1,n),gvec,d(1,n1)) |
---|
1369 | |
---|
1370 | call feti_mxvadd( nspdim, n1-1, d, gamm, d(1,n1) ) |
---|
1371 | |
---|
1372 | ! compute (Et E) * d |
---|
1373 | |
---|
1374 | call feti_proe(numnes,listnes,plistin,numins,listin,narea, & |
---|
1375 | numnos,nspdim,nsp,d(1,n1),v,w,bufin,bufout) |
---|
1376 | call feti_proet(numnes,listnes,plistin,numins,listin,numnos,w, & |
---|
1377 | nspdim,nsp,ad(1,n1),v) |
---|
1378 | |
---|
1379 | 50 continue |
---|
1380 | |
---|
1381 | if(narea == 1) then |
---|
1382 | write(0,*) 'No convergence for ete in ',nitmax,' iterations' |
---|
1383 | stop |
---|
1384 | endif |
---|
1385 | |
---|
1386 | |
---|
1387 | 999 continue |
---|
1388 | numit0=min0(n,nitmax-1) |
---|
1389 | ! write(0,*) 'number of iterations for ete :', nn+1, ' , residual =', facgn |
---|
1390 | |
---|
1391 | return |
---|
1392 | end |
---|
1393 | !**********************************************************************c |
---|
1394 | subroutine feti_osaxpy( n, a, x, y, z ) |
---|
1395 | |
---|
1396 | implicit none |
---|
1397 | |
---|
1398 | integer n |
---|
1399 | REAL(kind=8) a, x(n), y(n), z(n) |
---|
1400 | |
---|
1401 | integer i |
---|
1402 | |
---|
1403 | do 10 i = 1, n |
---|
1404 | z(i) = y(i) + a * x(i) |
---|
1405 | 10 continue |
---|
1406 | |
---|
1407 | return |
---|
1408 | end |
---|
1409 | !**********************************************************************c |
---|
1410 | subroutine feti_proe(numnes,listnes,plistin,numins,listin,narea, & |
---|
1411 | numnos,nspdim,nsp,alpha,v,w,bufin, & |
---|
1412 | bufout) |
---|
1413 | |
---|
1414 | ! jump of a velocity field on the interfaces, in the case where |
---|
1415 | ! subdomain number i is allocated to processor number i-1 . |
---|
1416 | ! the velocity field is Nsp * alpha, where Nsp is the nullspace . |
---|
1417 | |
---|
1418 | ! INPUTS : numnes = number of neighbouring subdomains |
---|
1419 | ! listnes = list of neighbouring subdomains |
---|
1420 | ! plistin = pointer of sublists of interface nodes |
---|
1421 | ! numins = number of interface nodes |
---|
1422 | ! listin = list of interface nodes |
---|
1423 | ! narea = subdomain number |
---|
1424 | ! numnos = number of nodes in the subdomain |
---|
1425 | ! nspdim = dimension of the null space |
---|
1426 | ! nsp = null space |
---|
1427 | ! alpha = components of the velocity field in the null |
---|
1428 | ! space |
---|
1429 | ! v = corresponding velocity field |
---|
1430 | |
---|
1431 | ! OUTPUTS : w = jump on the interface |
---|
1432 | |
---|
1433 | ! WORKSPACE : bufin = contributions of local field on interface |
---|
1434 | ! bufout = contributions of outer fields on interface |
---|
1435 | ! work |
---|
1436 | implicit none |
---|
1437 | integer numnes,numins,numnos,narea,nspdim |
---|
1438 | integer listnes(numnes,3), plistin(numnes+1), listin(numins) |
---|
1439 | REAL(kind=8) v(numnos),w(numins) |
---|
1440 | REAL(kind=8) bufin(numins),bufout(numins) |
---|
1441 | REAL(kind=8) alpha(nspdim), nsp(numnos,nspdim) |
---|
1442 | integer i |
---|
1443 | |
---|
1444 | ! computing the values of the local velocity field |
---|
1445 | do i=1,numnos |
---|
1446 | v(i)=0.e0 |
---|
1447 | end do |
---|
1448 | call feti_mxvadd( numnos, nspdim, nsp, alpha, v ) |
---|
1449 | ! computing the jump on the interface |
---|
1450 | call feti_saut(numnes,listnes,plistin,numins,listin, & |
---|
1451 | numnos,v,w,bufin,bufout) |
---|
1452 | !! narea,numnos,v,w,bufin,bufout) |
---|
1453 | |
---|
1454 | return |
---|
1455 | end |
---|
1456 | !**********************************************************************c |
---|
1457 | subroutine feti_proet(numnes,listnes,plistin,numins,listin, & |
---|
1458 | numnos,w,nspdim,nsp,alpha,v) |
---|
1459 | |
---|
1460 | ! projection in the null space of the right-hand-side associated with |
---|
1461 | ! a Lagrange multiplier w . |
---|
1462 | |
---|
1463 | ! INPUTS : numins = number of interface nodes |
---|
1464 | ! listin = list of interface nodes |
---|
1465 | ! numnos = number of nodes in the subdomain |
---|
1466 | ! w = Lagrange multiplier |
---|
1467 | ! nspdim = dimension of the null space |
---|
1468 | ! nsp = null space |
---|
1469 | |
---|
1470 | ! OUTPUTS : v = local right-hand-side |
---|
1471 | ! alpha = projection of v in the null space |
---|
1472 | |
---|
1473 | implicit none |
---|
1474 | external sdot |
---|
1475 | REAL(kind=8) sdot |
---|
1476 | integer numins,numnos,nspdim,numnes |
---|
1477 | integer listnes(numnes,3),plistin(numnes+1),listin(numins) |
---|
1478 | REAL(kind=8) v(numnos),w(numins) |
---|
1479 | REAL(kind=8) alpha(nspdim), nsp(numnos,nspdim) |
---|
1480 | integer i |
---|
1481 | |
---|
1482 | ! assembling the values of the Lagrange multiplier on the interface |
---|
1483 | call feti_probit(numnes,listnes,plistin,numins,listin,numnos,w,v) |
---|
1484 | ! computing the projection of v in the null space . |
---|
1485 | do 1 i=1,nspdim |
---|
1486 | alpha(i)=sdot(numnos,v,1,nsp(1,i),1) |
---|
1487 | 1 continue |
---|
1488 | |
---|
1489 | return |
---|
1490 | end |
---|
1491 | !**********************************************************************c |
---|
1492 | subroutine feti_consrhs( vdim,nspdim,rhs,nsp,cb ) |
---|
1493 | |
---|
1494 | |
---|
1495 | implicit none |
---|
1496 | |
---|
1497 | external sdot |
---|
1498 | REAL(kind=8) sdot |
---|
1499 | integer vdim, nspdim |
---|
1500 | |
---|
1501 | REAL(kind=8) rhs(vdim), nsp(vdim,nspdim), cb(nspdim) |
---|
1502 | |
---|
1503 | |
---|
1504 | integer i |
---|
1505 | |
---|
1506 | ! compute the right hand side for constraint |
---|
1507 | |
---|
1508 | do i=1,nspdim |
---|
1509 | cb(i)=-sdot(vdim,rhs,1,nsp(1,i),1) |
---|
1510 | end do |
---|
1511 | |
---|
1512 | return |
---|
1513 | end |
---|
1514 | !**********************************************************************c |
---|
1515 | subroutine feti_dualschur(noeuds,ni,nj,a,lpblo,blo, & |
---|
1516 | ninterf,ninterfc,nni,nnic, & |
---|
1517 | ndvois,ndvoisc,plistin,listin, & |
---|
1518 | poids,u,v,f,bitw,utilu, & |
---|
1519 | lambda,g,pg,mg,nitmax,nmaxd,j0,wj,dwj,dwwj, & |
---|
1520 | gamm,work,bufin,bufout,narea,epsilo,ndlblo, & |
---|
1521 | lisblo,ndkerep, & |
---|
1522 | xnul,ynul,numit0ete,nitmaxete,eteg, & |
---|
1523 | eteag,eted,etead,eteadd,etegamm,nsp,etev, & |
---|
1524 | etew,nnih,plistih,gh,w,dw, & |
---|
1525 | residu,indic,jn) |
---|
1526 | implicit none |
---|
1527 | ! nombre de ddl par noeuds, nombre de noeuds interne et interface |
---|
1528 | integer noeuds,ni,nj,nni,nnic |
---|
1529 | ! integer jpj,jpi |
---|
1530 | ! numero du sous-domaine = numero de processeur + 1 |
---|
1531 | integer narea |
---|
1532 | ! tableaux descripteurs de la matrice locale |
---|
1533 | REAL(kind=8) a(ni,nj,5) |
---|
1534 | ! tableaux descripteurs de l'inverse de la matrice du probleme local |
---|
1535 | integer lpblo |
---|
1536 | REAL(kind=8) blo(lpblo) |
---|
1537 | ! tableaux descripteurs de l'interface |
---|
1538 | integer ninterf,ninterfc |
---|
1539 | integer ndvois(ninterf,3),ndvoisc(ninterfc,3) |
---|
1540 | integer plistin(ninterfc+1),listin(nnic) |
---|
1541 | ! tableaux descripteurs de l'interface ajoutes pour nperio == 1 |
---|
1542 | REAL(kind=8) poids(noeuds) |
---|
1543 | ! utilitaires pour les assemblages aux interfaces |
---|
1544 | REAL(kind=8) bufin(nnic),bufout(nnic) |
---|
1545 | ! vecteurs locaux |
---|
1546 | REAL(kind=8) u(noeuds),v(noeuds),f(noeuds) |
---|
1547 | REAL(kind=8) bitw(noeuds),utilu(noeuds) |
---|
1548 | ! vecteurs inconnues aux interfaces |
---|
1549 | ! mg est le gradient preconditionne : M = sum(Bi Ai Bit) |
---|
1550 | REAL(kind=8) lambda(nni),g(nni),pg(nni),mg(nni) |
---|
1551 | ! tableaux servant a gerer le stockage de la moitie des coefficients |
---|
1552 | ! des tableaux des directions de descente sur chaque interface |
---|
1553 | integer nnih |
---|
1554 | integer plistih(ninterf+1) |
---|
1555 | REAL(kind=8) gh(nnih) |
---|
1556 | ! tableaux des directions de descente |
---|
1557 | ! on appelle D la matrice condensee a l'interface : sum(Bi [Ai]-1 Bit) |
---|
1558 | ! j0 est le nombre de vecteurs D orthogonaux deja stockes |
---|
1559 | ! j0 est egal au plus a nmaxd-1 |
---|
1560 | ! lorsque le nomtre total d'iterations atteint nmaxd, seule la derniere |
---|
1561 | ! direction calculee est stockee dans la colonne numero nmaxd |
---|
1562 | integer nitmax,nmaxd,j0 |
---|
1563 | REAL(kind=8) wj(nnih,nmaxd),dwj(nnih,nmaxd),dwwj(nmaxd) |
---|
1564 | REAL(kind=8) w(nni),dw(nni) |
---|
1565 | ! utilitaires pour les calculs de produits scalaires globaux |
---|
1566 | REAL(kind=8) work(nmaxd),gamm(nmaxd) |
---|
1567 | ! tableaux descripteurs du noyau de la matrice locale et tableaux |
---|
1568 | ! servant a la projection sur le noyau |
---|
1569 | |
---|
1570 | !...ndlblo is the dimension of the local nullspace .=<. the size of the |
---|
1571 | ! memory of the superstructure associated to the nullspace : ndkerep |
---|
1572 | ! indeed ndkerep = Max ndlblo = Max dim(Ker(Ep)) |
---|
1573 | ! p=1,Np |
---|
1574 | ! ndkerep is introduced to avoid messages "out of bounds" when memory |
---|
1575 | ! is checked |
---|
1576 | |
---|
1577 | integer ndkerep |
---|
1578 | integer ndlblo,numit0ete,nitmaxete |
---|
1579 | integer lisblo(ndlblo) |
---|
1580 | REAL(kind=8) xnul(ndlblo),ynul(ndlblo),eteg(ndlblo),eteag(ndlblo) |
---|
1581 | |
---|
1582 | !...the Krylov optimisation used in PCG associated to the Coarse Solver |
---|
1583 | ! involve to keep di & adi vectors in etesolve procedure |
---|
1584 | ! the ndkerep parameter is introduced to avoid message "out of bound" |
---|
1585 | ! arrising when the memory is checked! |
---|
1586 | |
---|
1587 | ! ndkerep = Max (ndlblo(Ep)) <= memory management (f77 context) |
---|
1588 | ! p=1,Np |
---|
1589 | |
---|
1590 | ! REAL(kind=8) eted(ndlblo,nitmaxete),etead(ndlblo,nitmaxete) |
---|
1591 | |
---|
1592 | REAL(kind=8) eted(ndkerep,nitmaxete),etead(ndkerep,nitmaxete) |
---|
1593 | REAL(kind=8) eteadd(nitmaxete),etegamm(nitmaxete) |
---|
1594 | REAL(kind=8) etev(noeuds),etew(nni),nsp(noeuds,ndlblo) |
---|
1595 | ! indicateur de convergence |
---|
1596 | |
---|
1597 | integer indic |
---|
1598 | REAL(kind=8) residu |
---|
1599 | integer jn |
---|
1600 | ! utilitaires |
---|
1601 | !! integer j,k,j1,ji,jj,jk |
---|
1602 | integer j,k,j1 |
---|
1603 | REAL(kind=8) scal(2),ff,roj,epsilo,epsilo2 |
---|
1604 | ! fonctions externes |
---|
1605 | external sdot |
---|
1606 | REAL(kind=8) sdot |
---|
1607 | ! booleen pour caracterisation |
---|
1608 | ! numerique si llconv=1 |
---|
1609 | integer llconv |
---|
1610 | |
---|
1611 | llconv=0 |
---|
1612 | |
---|
1613 | !...critere d'arret epsilo2 = epsilo **2 => evite sqrt |
---|
1614 | |
---|
1615 | epsilo2 = epsilo * epsilo |
---|
1616 | |
---|
1617 | ! calcul du carre scalaire du second membre global |
---|
1618 | |
---|
1619 | call feti_assemb(ninterfc,ndvoisc,plistin,nnic,listin,narea,noeuds & |
---|
1620 | ,f,utilu,bufin,bufout) |
---|
1621 | ff=sdot(noeuds,utilu,1,utilu,1) |
---|
1622 | call mpp_sum(ff,1,work) |
---|
1623 | ! write(90,*)'valeur du carre du second membre ds le Schur Dual :',ff |
---|
1624 | |
---|
1625 | ! calcul du second membre contraint |
---|
1626 | |
---|
1627 | call feti_consrhs(noeuds,ndlblo,f,nsp,ynul) |
---|
1628 | |
---|
1629 | ! calcul de lambda0 satisfaisant la contrainte |
---|
1630 | |
---|
1631 | ! call feti_vclr(nni,lambda) |
---|
1632 | call feti_etesolv(ndlblo,xnul,ynul,ndkerep,numit0ete,nitmaxete & |
---|
1633 | ,eteg,eteag,eted,etead,eteadd,etegamm,ninterf,ndvois,plistin & |
---|
1634 | ,nni,listin,narea,noeuds,nsp,etev,etew,bufin,bufout,work) |
---|
1635 | |
---|
1636 | call feti_proe(ninterf,ndvois,plistin,nni,listin,narea,noeuds, & |
---|
1637 | ndlblo,nsp,xnul,etev,lambda,bufin,bufout) |
---|
1638 | |
---|
1639 | ! calcul du nouveau second membre |
---|
1640 | |
---|
1641 | call feti_probit(ninterf,ndvois,plistin,nni,listin,noeuds,lambda,bitw) |
---|
1642 | |
---|
1643 | call feti_vadd(noeuds,bitw,f,bitw) |
---|
1644 | |
---|
1645 | |
---|
1646 | ! calcul du champ local initial u0 |
---|
1647 | ! call feti_blodir(noeuds,bitw,ndlblo,lisblo) |
---|
1648 | call feti_resloc(noeuds,ni,nj,a,lpblo,blo,bitw,u,utilu) |
---|
1649 | |
---|
1650 | ! calcul du gradient initial g0 = feti_saut(u0) |
---|
1651 | !! call feti_saut(ninterf,ndvois,plistin,nni,listin,narea,noeuds,u,g & |
---|
1652 | call feti_saut(ninterf,ndvois,plistin,nni,listin, noeuds,u,g & |
---|
1653 | ,bufin,bufout) |
---|
1654 | |
---|
1655 | ! reinitialisation de lambdaj0 par reconjugaison |
---|
1656 | |
---|
1657 | call feti_restri(ninterf,ndvois,plistin,nni,plistih,nnih,narea,g, gh) |
---|
1658 | do k=1,j0 |
---|
1659 | gamm(k)=sdot(nnih,gh,1,wj(1,k),1) |
---|
1660 | end do |
---|
1661 | if(j0 >= 1) call mpp_sum(gamm,j0,work) |
---|
1662 | do k=1,j0 |
---|
1663 | gamm(k)=-gamm(k)/dwwj(k) |
---|
1664 | end do |
---|
1665 | call feti_restri(ninterf,ndvois,plistin,nni,plistih,nnih,narea, lambda,gh) |
---|
1666 | call feti_mxvadd(nnih,j0,wj,gamm,gh) |
---|
1667 | call feti_extend(ninterf,ndvois,plistin,nni,plistih,nnih,narea,gh, lambda) |
---|
1668 | ! calcul du second membre associe a lambdaj0 |
---|
1669 | if(j0 > 0) then |
---|
1670 | call feti_probit(ninterf,ndvois,plistin,nni,listin,noeuds, lambda,bitw) |
---|
1671 | call feti_vadd(noeuds,bitw,f,bitw) |
---|
1672 | ! calcul du champ local associe uj0 |
---|
1673 | ! call feti_blodir(noeuds,bitw,ndlblo,lisblo) |
---|
1674 | call feti_resloc(noeuds,ni,nj,a,lpblo,blo,bitw,u,utilu) |
---|
1675 | ! calcul du gradient initial gj0= feti_saut(uj0) |
---|
1676 | !! call feti_saut(ninterf,ndvois,plistin,nni,listin,narea,noeuds & |
---|
1677 | call feti_saut(ninterf,ndvois,plistin,nni,listin, noeuds & |
---|
1678 | ,u,g,bufin,bufout) |
---|
1679 | endif |
---|
1680 | ! calcul du gradient projete |
---|
1681 | ! call feti_vmov(nni,g,pg) |
---|
1682 | call feti_project(g,pg,ndlblo,xnul,ynul,ndkerep,numit0ete & |
---|
1683 | ,nitmaxete,eteg,eteag,eted,etead,eteadd,etegamm,ninterf,ndvois & |
---|
1684 | ,plistin,nni,listin,narea,noeuds,nsp,etev,etew,bufin,bufout,work) |
---|
1685 | ! calcul du gradient preconditionne |
---|
1686 | ! call feti_vmov(nni,pg,mg) |
---|
1687 | |
---|
1688 | call feti_probit(ninterf,ndvois,plistin,nni,listin,noeuds, pg,bitw) |
---|
1689 | |
---|
1690 | call feti_proax(noeuds,ni,nj,a,bitw,v) |
---|
1691 | !! call feti_saut(ninterf,ndvois,plistin,nni,listin,narea,noeuds,v, mg,bufin,bufout) |
---|
1692 | call feti_saut(ninterf,ndvois,plistin,nni,listin, noeuds,v, mg,bufin,bufout) |
---|
1693 | ! verification du residu initial |
---|
1694 | residu=sdot(nni,mg,1,mg,1) |
---|
1695 | call mpp_sum(residu,1,work) |
---|
1696 | residu = (4. * residu) / ff |
---|
1697 | if(residu <= epsilo2)then |
---|
1698 | ! residu=10.*sqrt(residu/ff) |
---|
1699 | ! if(residu < epsilo)then |
---|
1700 | if(llconv == 1)then |
---|
1701 | if(narea == 1)then |
---|
1702 | write(0,*)' residu carre approche initial :',residu |
---|
1703 | endif |
---|
1704 | endif |
---|
1705 | ! calcul du champ global moyenne aux interfaces et du residu global |
---|
1706 | call feti_vsub(noeuds,u,etev,utilu) |
---|
1707 | call feti_moyen(ninterfc,ndvoisc,plistin,nnic,listin,narea, & |
---|
1708 | noeuds,poids,utilu,utilu,bufin,bufout) |
---|
1709 | call feti_vmov(noeuds,utilu,u) |
---|
1710 | call feti_proax(noeuds,ni,nj,a,utilu,v) |
---|
1711 | call feti_vsub(noeuds,v,f,v) |
---|
1712 | call feti_assemb(ninterfc,ndvoisc,plistin,nnic,listin,narea, & |
---|
1713 | noeuds,v,utilu,bufin,bufout) |
---|
1714 | residu=sdot(noeuds,utilu,1,utilu,1) |
---|
1715 | call mpp_sum(residu,1,work) |
---|
1716 | residu=sqrt(residu/ff) |
---|
1717 | ! residu=residu/ff |
---|
1718 | if(llconv == 1)then |
---|
1719 | if(narea == 1)then |
---|
1720 | write(0,*)' residu exact initial :',residu |
---|
1721 | endif |
---|
1722 | endif |
---|
1723 | |
---|
1724 | 101 format(10d9.2) |
---|
1725 | |
---|
1726 | return |
---|
1727 | |
---|
1728 | endif |
---|
1729 | |
---|
1730 | |
---|
1731 | ! calcul de la premiere direction de descente |
---|
1732 | ! calcul du gradient preconditionne projete |
---|
1733 | call feti_project(mg,mg,ndlblo,xnul,ynul,ndkerep,numit0ete & |
---|
1734 | ,nitmaxete,eteg,eteag,eted,etead,eteadd,etegamm,ninterf,ndvois & |
---|
1735 | ,plistin,nni,listin,narea,noeuds,nsp,etev,etew,bufin,bufout ,work) |
---|
1736 | call feti_restri(ninterf,ndvois,plistin,nni,plistih,nnih,narea, & |
---|
1737 | mg,wj(1,j0+1)) |
---|
1738 | do k=1,j0 |
---|
1739 | gamm(k)=sdot(nnih,wj(1,j0+1),1,dwj(1,k),1) |
---|
1740 | end do |
---|
1741 | if(j0 >= 1) call mpp_sum(gamm,j0,work) |
---|
1742 | do k=1,j0 |
---|
1743 | gamm(k)=-gamm(k)/dwwj(k) |
---|
1744 | end do |
---|
1745 | call feti_mxvadd(nnih,j0,wj,gamm,wj(1,j0+1)) |
---|
1746 | |
---|
1747 | ! iterations de gradient conjugue a l'interface |
---|
1748 | |
---|
1749 | |
---|
1750 | !debug |
---|
1751 | do j=j0+1,j0+nitmax |
---|
1752 | |
---|
1753 | if(llconv == 1.and.narea == 1) write(0,*)'etape numero',j-j0 |
---|
1754 | ! j1 est le numero de colonne de la derniere direction de descente |
---|
1755 | ! j1 = j, numero d'iteration total, tant que j < nmaxd |
---|
1756 | ! au-dela, j1 = nmaxd |
---|
1757 | j1=min0(j,nmaxd) |
---|
1758 | ! calcul de Dwi , D = sum ( Bi [Ai]-1 tBi ) |
---|
1759 | call feti_extend(ninterf,ndvois,plistin,nni,plistih,nnih,narea, wj(1,j1),w) |
---|
1760 | call feti_probit(ninterf,ndvois,plistin,nni,listin,noeuds, w,bitw) |
---|
1761 | ! call feti_blodir(noeuds,bitw,ndlblo,lisblo) |
---|
1762 | call feti_resloc(noeuds,ni,nj,a,lpblo,blo,bitw,v,utilu) |
---|
1763 | !! call feti_saut(ninterf,ndvois,plistin,nni,listin,narea,noeuds,v & |
---|
1764 | call feti_saut(ninterf,ndvois,plistin,nni,listin, noeuds,v & |
---|
1765 | ,dw,bufin,bufout) |
---|
1766 | call feti_restri(ninterf,ndvois,plistin,nni,plistih,nnih,narea, dw,dwj(1,j1)) |
---|
1767 | ! calcul du coefficient de descente |
---|
1768 | call feti_restri(ninterf,ndvois,plistin,nni,plistih,nnih,narea, g,gh) |
---|
1769 | scal(1)=sdot(nnih,gh,1,wj(1,j1),1) |
---|
1770 | scal(2)=sdot(nnih,dwj(1,j1),1,wj(1,j1),1) |
---|
1771 | call mpp_sum(scal,2,work) |
---|
1772 | roj=-scal(1)/scal(2) |
---|
1773 | dwwj(j1)=scal(2) |
---|
1774 | ! remise a jour de lambda, g et u |
---|
1775 | call saxpy(nni,roj,w,1,lambda,1) |
---|
1776 | call saxpy(nni,roj,dw,1,g,1) |
---|
1777 | call saxpy(noeuds,roj,v,1,u,1) |
---|
1778 | ! calcul du gradient projete |
---|
1779 | ! call feti_vmov(nni,g,pg) |
---|
1780 | call feti_project(g,pg,ndlblo,xnul,ynul,ndkerep,numit0ete & |
---|
1781 | ,nitmaxete,eteg,eteag,eted,etead,eteadd,etegamm,ninterf & |
---|
1782 | ,ndvois,plistin,nni,listin,narea,noeuds,nsp,etev,etew,bufin & |
---|
1783 | ,bufout,work) |
---|
1784 | ! calcul du gradient preconditionne Mg , M = sum ( Bi Ai tBi ) |
---|
1785 | ! call feti_vmov(nni,pg,mg) |
---|
1786 | call feti_probit(ninterf,ndvois,plistin,nni,listin,noeuds, pg,bitw) |
---|
1787 | call feti_proax(noeuds,ni,nj,a,bitw,v) |
---|
1788 | !! call feti_saut(ninterf,ndvois,plistin,nni,listin,narea,noeuds,v & |
---|
1789 | call feti_saut(ninterf,ndvois,plistin,nni,listin, noeuds,v & |
---|
1790 | ,mg,bufin,bufout) |
---|
1791 | ! verification du residu global approche |
---|
1792 | residu=sdot(nni,mg,1,mg,1) |
---|
1793 | call mpp_sum(residu,1,work) |
---|
1794 | ! residu=10.*sqrt(residu/ff) |
---|
1795 | residu = (4. * residu) / ff |
---|
1796 | if(llconv == 1.and.narea == 1)then |
---|
1797 | write(0,*)' residu carre global approche apres ', & |
---|
1798 | j-j0,' iterations :',residu |
---|
1799 | endif |
---|
1800 | |
---|
1801 | ! if(residu < epsilo) then |
---|
1802 | if(residu <= epsilo2)then |
---|
1803 | if(llconv == 1.and.narea == 1)then |
---|
1804 | write(0,*)'residu carre global approche apres ', & |
---|
1805 | j-j0,' iterations :',residu |
---|
1806 | endif |
---|
1807 | ! calcul du champ global moyenne aux interfaces et du residu global |
---|
1808 | call feti_vsub(noeuds,u,etev,utilu) |
---|
1809 | call feti_moyen(ninterfc,ndvoisc,plistin,nnic,listin,narea, & |
---|
1810 | noeuds,poids,utilu,utilu,bufin,bufout) |
---|
1811 | call feti_vmov(noeuds,utilu,u) |
---|
1812 | call feti_proax(noeuds,ni,nj,a,utilu,v) |
---|
1813 | call feti_vsub(noeuds,v,f,v) |
---|
1814 | call feti_assemb(ninterfc,ndvoisc,plistin,nnic,listin,narea, & |
---|
1815 | noeuds,v,utilu,bufin,bufout) |
---|
1816 | residu=sdot(noeuds,utilu,1,utilu,1) |
---|
1817 | call mpp_sum(residu,1,work) |
---|
1818 | ! residu=residu/ff |
---|
1819 | residu=sqrt(residu/ff) |
---|
1820 | if(llconv == 1.and.narea == 1)then |
---|
1821 | write(0,*)' residu global exact apres ', & |
---|
1822 | j-j0,' iterations :',residu |
---|
1823 | endif |
---|
1824 | |
---|
1825 | jn = j - j0 |
---|
1826 | j0=min0(j,nmaxd-1) |
---|
1827 | if(llconv == 1.and.narea == 1)then |
---|
1828 | write(0,*)' nombre de directions orthogonales conservees :',j0 |
---|
1829 | endif |
---|
1830 | |
---|
1831 | |
---|
1832 | return |
---|
1833 | endif |
---|
1834 | !################################################################## |
---|
1835 | ! calcul du champ global moyenne aux interfaces et du residu global |
---|
1836 | call feti_vsub(noeuds,u,etev,utilu) |
---|
1837 | call feti_moyen(ninterfc,ndvoisc,plistin,nnic,listin,narea, & |
---|
1838 | noeuds,poids,utilu,utilu,bufin,bufout) |
---|
1839 | call feti_proax(noeuds,ni,nj,a,utilu,v) |
---|
1840 | call feti_vsub(noeuds,v,f,v) |
---|
1841 | call feti_assemb(ninterfc,ndvoisc,plistin,nnic,listin,narea, & |
---|
1842 | noeuds,v,utilu,bufin,bufout) |
---|
1843 | residu=sdot(noeuds,utilu,1,utilu,1) |
---|
1844 | call mpp_sum(residu,1,work) |
---|
1845 | ! residu=sqrt(residu/ff) |
---|
1846 | residu=residu/ff |
---|
1847 | if(llconv == 1.and.narea == 1)then |
---|
1848 | write(0,*)' residu carre global exact apres ', & |
---|
1849 | j-j0,' iterations :',residu |
---|
1850 | endif |
---|
1851 | |
---|
1852 | !##################################################################### |
---|
1853 | ! calcul de la nouvelle direction de descente par reconjugaison |
---|
1854 | ! calcul du gradient preconditionne projete |
---|
1855 | call feti_project(mg,mg,ndlblo,xnul,ynul,ndkerep,numit0ete & |
---|
1856 | ,nitmaxete,eteg,eteag,eted,etead,eteadd,etegamm,ninterf & |
---|
1857 | ,ndvois,plistin,nni,listin,narea,noeuds,nsp,etev,etew,bufin & |
---|
1858 | ,bufout,work) |
---|
1859 | call feti_restri(ninterf,ndvois,plistin,nni,plistih,nnih,narea, mg,gh) |
---|
1860 | do k=1,j1 |
---|
1861 | gamm(k)=sdot(nnih,gh,1,dwj(1,k),1) |
---|
1862 | end do |
---|
1863 | call mpp_sum(gamm,j1,work) |
---|
1864 | do k=1,j1 |
---|
1865 | gamm(k)=-gamm(k)/dwwj(k) |
---|
1866 | end do |
---|
1867 | call feti_mxvadd(nnih,j1,wj,gamm,gh) |
---|
1868 | j1=min0(j1+1,nmaxd) |
---|
1869 | |
---|
1870 | !%%%%%%%%%%stabilisation numerique %%%%%%%%%%%%%%%%%%%%%%%%% |
---|
1871 | call feti_extend(ninterf,ndvois,plistin,nni,plistih,nnih,narea, gh,w) |
---|
1872 | call feti_project(w,w,ndlblo,xnul,ynul,ndkerep,numit0ete & |
---|
1873 | ,nitmaxete,eteg,eteag,eted,etead,eteadd,etegamm,ninterf & |
---|
1874 | ,ndvois,plistin,nni,listin,narea,noeuds,nsp,etev,etew,bufin & |
---|
1875 | ,bufout,work) |
---|
1876 | call feti_restri(ninterf,ndvois,plistin,nni,plistih,nnih,narea, w,gh) |
---|
1877 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
1878 | |
---|
1879 | call feti_vmov(nnih,gh,wj(1,j1)) |
---|
1880 | |
---|
1881 | end do |
---|
1882 | |
---|
1883 | residu = sqrt(residu) |
---|
1884 | if(narea == 1) then |
---|
1885 | write(0,*)' residu global approche apres ',nitmax, & |
---|
1886 | ' iterations :',residu |
---|
1887 | endif |
---|
1888 | indic=-2 |
---|
1889 | |
---|
1890 | return |
---|
1891 | end |
---|
1892 | !**********************************************************************c |
---|
1893 | subroutine feti_extend(numnes,listnes, & |
---|
1894 | plistin,numins,plistih, & |
---|
1895 | numinh,narea,wh,w) |
---|
1896 | |
---|
1897 | ! reassembling an half interface field |
---|
1898 | |
---|
1899 | ! INPUTS : numnes = number of neighbouring subdomains |
---|
1900 | ! listnes = list of neighbouring subdomains |
---|
1901 | ! lint = list of inner interfaces |
---|
1902 | ! lext = list of extern ibnterfaces |
---|
1903 | ! plistin = pointer of sublists of interface nodes |
---|
1904 | ! plistih = pointer of half-sublists of interface nodes |
---|
1905 | ! numins = number of interface nodes |
---|
1906 | ! numinh = half number of interface nodes |
---|
1907 | ! narea = subdomain number |
---|
1908 | ! wh = half interface vector |
---|
1909 | |
---|
1910 | ! OUTPUTS : w = complete interface vector |
---|
1911 | |
---|
1912 | ! WORKSPACE : work |
---|
1913 | implicit none |
---|
1914 | ! external irecv |
---|
1915 | !! integer irecv, work(100) |
---|
1916 | integer numnes,numins,numinh,narea |
---|
1917 | integer listnes(numnes,3), plistin(numnes+1), plistih(numnes+1) |
---|
1918 | REAL(kind=8) wh(numinh),w(numins) |
---|
1919 | integer i,l,lin,lout |
---|
1920 | integer mesg |
---|
1921 | |
---|
1922 | mesg=128*10 |
---|
1923 | ! if subdomains i < j are neighbours, subdomain i manage the first |
---|
1924 | ! half of the interface, and subdomain j the remaining part |
---|
1925 | |
---|
1926 | ! receiving the values of the outer fields on the interface |
---|
1927 | ! sending the values of the inner field on the interface |
---|
1928 | |
---|
1929 | do i=1,numnes |
---|
1930 | l=(plistin(i+1)-plistin(i)) |
---|
1931 | ! if(narea < listnes(i)) then |
---|
1932 | if(listnes(i,3) < listnes(i,2))then |
---|
1933 | lout=(l-(l/2)) |
---|
1934 | lin=(l/2) |
---|
1935 | ! work(i)=irecv(mesg+listnes(i)-1,w(plistin(i)+(l/2)+1), lout*8) |
---|
1936 | ! call mppsend(mesg+narea-1,wh(plistih(i)+1),lin*8,listnes(i) -1,0) |
---|
1937 | call mppsend(mesg+listnes(i,3)-1,wh(plistih(i)+1),lin*8, listnes(i,1)-1,0) |
---|
1938 | call feti_vmov(lin,wh(plistih(i)+1),w(plistin(i)+1)) |
---|
1939 | else |
---|
1940 | lout=(l/2) |
---|
1941 | lin=(l-(l/2)) |
---|
1942 | ! work(i)=irecv(mesg+listnes(i)-1,w(plistin(i)+1),lout*8) |
---|
1943 | ! call mppsend(mesg+narea-1,wh(plistih(i)+1),lin*8,listnes(i)-1,0) |
---|
1944 | call mppsend(mesg+listnes(i,3)-1,wh(plistih(i)+1),lin*8, listnes(i,1)-1,0) |
---|
1945 | call feti_vmov(lin,wh(plistih(i)+1),w(plistin(i)+(l/2)+1)) |
---|
1946 | endif |
---|
1947 | enddo |
---|
1948 | ! waiting for the completion on each interface |
---|
1949 | do i=1,numnes |
---|
1950 | l=(plistin(i+1)-plistin(i)) |
---|
1951 | ! if(narea < listnes(i)) then |
---|
1952 | if(listnes(i,3) < listnes(i,2))then |
---|
1953 | lout=(l-(l/2)) |
---|
1954 | ! call mpprecv(mesg+listnes(i)-1,w(plistin(i)+(l/2)+1), lout*8) |
---|
1955 | call mpprecv(mesg+listnes(i,2)-1,w(plistin(i)+(l/2)+1), lout*8) |
---|
1956 | call feti_vneg(lout,w(plistin(i)+(l/2)+1), & |
---|
1957 | w(plistin(i)+(l/2)+1)) |
---|
1958 | else |
---|
1959 | lout=(l/2) |
---|
1960 | ! call mpprecv(mesg+listnes(i)-1,w(plistin(i)+1),lout*8) |
---|
1961 | call mpprecv(mesg+listnes(i,2)-1,w(plistin(i)+1),lout*8) |
---|
1962 | call feti_vneg(lout,w(plistin(i)+1),w(plistin(i)+1)) |
---|
1963 | endif |
---|
1964 | end do |
---|
1965 | |
---|
1966 | return |
---|
1967 | end |
---|
1968 | !**********************************************************************c |
---|
1969 | subroutine feti_restri(numnes,listnes,plistin,numins,plistih, & |
---|
1970 | numinh,narea,w,wh) |
---|
1971 | |
---|
1972 | ! restriction of an interface field to one half on each interface |
---|
1973 | |
---|
1974 | ! INPUTS : numnes = number of neighbouring subdomains |
---|
1975 | ! listnes = list of neighbouring subdomains |
---|
1976 | ! plistin = pointer of sublists of interface nodes |
---|
1977 | ! plistih = pointer of half-sublists of interface nodes |
---|
1978 | ! numins = number of interface nodes |
---|
1979 | ! numinh = half number of interface nodes |
---|
1980 | ! narea = subdomain number |
---|
1981 | ! w = complete interface vector |
---|
1982 | |
---|
1983 | ! OUTPUTS : wh = half interface vector |
---|
1984 | |
---|
1985 | ! WORKSPACE : work |
---|
1986 | implicit none |
---|
1987 | integer numnes,numins,numinh,narea |
---|
1988 | integer listnes(numnes,3), plistin(numnes+1), plistih(numnes+1) |
---|
1989 | REAL(kind=8) wh(numinh),w(numins) |
---|
1990 | integer i,l,lin |
---|
1991 | |
---|
1992 | ! if subdomains i < j are neighbours, subdomain i manage the first |
---|
1993 | ! half of the interface, and subdomain j the remaining part |
---|
1994 | |
---|
1995 | do i=1,numnes |
---|
1996 | l=(plistin(i+1)-plistin(i)) |
---|
1997 | ! if(narea < listnes(i)) then |
---|
1998 | if(listnes(i,3) < listnes(i,2))then |
---|
1999 | lin=(l/2) |
---|
2000 | call feti_vmov(lin,w(plistin(i)+1),wh(plistih(i)+1)) |
---|
2001 | else |
---|
2002 | lin=(l-(l/2)) |
---|
2003 | call feti_vmov(lin,w(plistin(i)+(l/2)+1),wh(plistih(i)+1)) |
---|
2004 | endif |
---|
2005 | end do |
---|
2006 | |
---|
2007 | return |
---|
2008 | end |
---|
2009 | !**********************************************************************c |
---|
2010 | !! subroutine feti_halfint(numnes,listnes,plistin,numins,plistih & |
---|
2011 | !! ,numinh,narea) |
---|
2012 | subroutine feti_halfint(numnes,listnes,plistin, plistih & |
---|
2013 | ,numinh ) |
---|
2014 | |
---|
2015 | ! construction of the pointer of the restriction of an interface |
---|
2016 | ! field to one half on each interface |
---|
2017 | |
---|
2018 | ! INPUTS : numnes = number of neighbouring subdomains |
---|
2019 | ! listnes = list of neighbouring subdomains |
---|
2020 | ! plistin = pointer of sublists of interface nodes |
---|
2021 | ! numins = number of interface nodes |
---|
2022 | ! narea = subdomain number |
---|
2023 | |
---|
2024 | ! OUTPUTS : numinh = half number of interface nodes |
---|
2025 | ! plistih = pointer of half-sublists of interface nodes |
---|
2026 | |
---|
2027 | ! WORKSPACE : work |
---|
2028 | implicit none |
---|
2029 | !! integer numnes,numins,numinh,narea |
---|
2030 | integer numnes, numinh |
---|
2031 | integer listnes(numnes,3), plistin(numnes+1), plistih(numnes+1) |
---|
2032 | integer i,l |
---|
2033 | |
---|
2034 | ! if subdomains i < j are neighbours, subdomain i manage the first |
---|
2035 | ! half of the interface, and subdomain j the remaining part |
---|
2036 | |
---|
2037 | plistih(1)=0 |
---|
2038 | numinh=0 |
---|
2039 | do i=1,numnes |
---|
2040 | l=(plistin(i+1)-plistin(i)) |
---|
2041 | ! if(narea < listnes(i,1)) then |
---|
2042 | if(listnes(i,3) < listnes(i,2)) then |
---|
2043 | ! write(*,*)'je gere interface avec : ',listnes(i,1) |
---|
2044 | plistih(i+1)=plistih(i)+(l/2) |
---|
2045 | else |
---|
2046 | ! write(*,*)'il gere interface avec : ',listnes(i,1) |
---|
2047 | plistih(i+1)=plistih(i)+(l-(l/2)) |
---|
2048 | endif |
---|
2049 | end do |
---|
2050 | numinh=plistih(numnes+1) |
---|
2051 | |
---|
2052 | return |
---|
2053 | end |
---|
2054 | !*********************************************************************** |
---|
2055 | subroutine feti_front(n,ni,nj,a,lcmat,cmat,d,v,w,ndlblo,lisblo & |
---|
2056 | ,nmdlblo) |
---|
2057 | implicit REAL(kind=8) (a-h,o-z) |
---|
2058 | dimension a(n,5) |
---|
2059 | ! routine de calcul de la decomposition frontale par blocs de |
---|
2060 | ! la matrice pentadiagonale a . |
---|
2061 | ! cmat contient les termes des blocs diagonaux factorises . |
---|
2062 | dimension cmat(lcmat),d(n) |
---|
2063 | |
---|
2064 | !...ndlblo nombre de liberte bloques = nombre de pivots nuls = dim(Ker(Ep)) |
---|
2065 | |
---|
2066 | !...nmdlblo = Max (ndlblo(Ep)) <= superstructure management in f77 context |
---|
2067 | ! p=1,Np |
---|
2068 | |
---|
2069 | integer ndlblo,nmdlblo |
---|
2070 | integer lisblo(nmdlblo) |
---|
2071 | |
---|
2072 | ! d ne contient pas, en sortie, l'inverse de la diagonale de Crout . |
---|
2073 | ! v est un tableau servant a stocker la partie triangulaire inferieure |
---|
2074 | ! stricte d'un bloc diagonal , et w un tableau servant a stocker un |
---|
2075 | ! bloc sur-diagonal plein . |
---|
2076 | ! de dimension egale a la largeur maximale d'un front . |
---|
2077 | dimension v(ni,ni),w(ni) |
---|
2078 | !!! |
---|
2079 | ndimd=ni |
---|
2080 | lbd=ndimd*ndimd |
---|
2081 | |
---|
2082 | ! lecture du tout premier bloc diagonal |
---|
2083 | |
---|
2084 | nb=1 |
---|
2085 | ! adresse du bloc diagonal courant |
---|
2086 | iac=1 |
---|
2087 | call feti_liblod(n,ni,a,ndimd,nb,cmat(iac)) |
---|
2088 | ! inversion du premier bloc diagonal . |
---|
2089 | ! write(*,*)' inversion du premier bloc diagonal ' |
---|
2090 | j0=0 |
---|
2091 | call feti_ldlt(ndimd,cmat(iac),d(1),w,ndlblo,lisblo,nmdlblo,j0) |
---|
2092 | call feti_calinv(ndimd,cmat(iac),v) |
---|
2093 | |
---|
2094 | ! contraction puis inversion des blocs diagonaux de la la matrice |
---|
2095 | ! front par front . |
---|
2096 | |
---|
2097 | do 2 nb=2,nj |
---|
2098 | ! write(*,*)' front numero ',nb |
---|
2099 | ! adresse du bloc diagonal precedent |
---|
2100 | iap=iac |
---|
2101 | iac=(nb-1)*lbd+1 |
---|
2102 | |
---|
2103 | ! lecture du bloc sur-diagonal precedent: passage bande-vecteur |
---|
2104 | |
---|
2105 | call feti_liblos(n,a,ndimd,nb-1,w) |
---|
2106 | |
---|
2107 | ! elimination du bloc extra-diagonal . |
---|
2108 | |
---|
2109 | ! lecture du bloc diagonal courant: passage bande-plein |
---|
2110 | |
---|
2111 | call feti_liblod(n,ni,a,ndimd,nb,cmat(iac)) |
---|
2112 | |
---|
2113 | ! calcul du nouveau bloc diagonal courant . |
---|
2114 | |
---|
2115 | ! write(*,*)' calcul du nouveau bloc diagonal courant ' |
---|
2116 | call feti_nbdia(ndimd,w,cmat(iap),cmat(iac)) |
---|
2117 | |
---|
2118 | ! inversion du nouveau bloc diagonal courant . |
---|
2119 | |
---|
2120 | ! write(*,*)' factorisation du nouveau bloc diagonal courant ' |
---|
2121 | j0=(nb-1)*ndimd |
---|
2122 | call feti_ldlt(ndimd,cmat(iac),d((nb-1)*ndimd+1),w,ndlblo, & |
---|
2123 | lisblo,nmdlblo,j0) |
---|
2124 | call feti_calinv(ndimd,cmat(iac),v) |
---|
2125 | |
---|
2126 | 2 continue |
---|
2127 | |
---|
2128 | return |
---|
2129 | end |
---|
2130 | !**********************************************************************c |
---|
2131 | subroutine feti_calinv(dim,a,b) |
---|
2132 | ! calcul de la partie triangulaire inferieure du bloc : |
---|
2133 | ! b = [a]-1 puis recopie dans a complet . |
---|
2134 | implicit none |
---|
2135 | integer dim |
---|
2136 | REAL(kind=8) a(dim,dim),b(dim,dim) |
---|
2137 | integer i,j |
---|
2138 | |
---|
2139 | call feti_vclr(dim*dim,b) |
---|
2140 | do j=1,dim |
---|
2141 | b(j,j)=1.e0 |
---|
2142 | call feti_desremo(dim,j-1,a,b(1,j),b(1,j)) |
---|
2143 | end do |
---|
2144 | do j=1,dim |
---|
2145 | do i=j,dim |
---|
2146 | a(i,j)=b(i,j) |
---|
2147 | a(j,i)=b(i,j) |
---|
2148 | end do |
---|
2149 | end do |
---|
2150 | |
---|
2151 | return |
---|
2152 | end |
---|
2153 | !*********************************************************************** |
---|
2154 | subroutine feti_nbdia(n,d,w,v) |
---|
2155 | implicit REAL(kind=8) (a-h,o-z) |
---|
2156 | dimension d(n),v(n,n),w(n,n) |
---|
2157 | |
---|
2158 | do 1 j=1,n |
---|
2159 | do 2 i=1,n |
---|
2160 | v(i,j)=v(i,j)-d(i)*w(i,j)*d(j) |
---|
2161 | 2 continue |
---|
2162 | 1 continue |
---|
2163 | |
---|
2164 | return |
---|
2165 | end |
---|
2166 | !*********************************************************************** |
---|
2167 | subroutine feti_liblod(n,ni,a,ndimd,nb,db) |
---|
2168 | |
---|
2169 | implicit REAL(kind=8) (a-h,o-z) |
---|
2170 | dimension db(ndimd,ndimd) |
---|
2171 | ! tableaux descripteurs de la structure morse de la matrice . |
---|
2172 | dimension a(n,5) |
---|
2173 | |
---|
2174 | call feti_vclr(ndimd*ndimd,db) |
---|
2175 | ! do i=1,ndimd*ndimd |
---|
2176 | ! db(i,1)=0.e0 |
---|
2177 | !8 continue |
---|
2178 | n0=(nb-1)*ndimd |
---|
2179 | ! lecture de la diagonale db(i,i) . |
---|
2180 | do 3 i=1,ndimd |
---|
2181 | db(i,i)=a(n0+i,3) |
---|
2182 | 3 continue |
---|
2183 | ! lecture de la premiere sur-diagonale db(i,i+1) . |
---|
2184 | do 4 i=1,ndimd-1 |
---|
2185 | db(i,i+1)=a(n0+i,4) |
---|
2186 | 4 continue |
---|
2187 | ! lecture de la premiere sous-diagonale db(i,i-1) . |
---|
2188 | do 2 i=2,ndimd |
---|
2189 | db(i,i-1)=a(n0+i,2) |
---|
2190 | 2 continue |
---|
2191 | |
---|
2192 | return |
---|
2193 | end |
---|
2194 | !*********************************************************************** |
---|
2195 | subroutine feti_liblos(n,a,ndimd,nb,v) |
---|
2196 | |
---|
2197 | implicit REAL(kind=8) (a-h,o-z) |
---|
2198 | dimension v(ndimd) |
---|
2199 | ! tableaux descripteurs de la structure morse de la matrice . |
---|
2200 | dimension a(n,5) |
---|
2201 | |
---|
2202 | n0=(nb-1)*ndimd |
---|
2203 | ! lecture de la deuxieme sur-diagonale v(i) . |
---|
2204 | do 5 i=1,ndimd |
---|
2205 | v(i)=a(n0+i,5) |
---|
2206 | 5 continue |
---|
2207 | |
---|
2208 | return |
---|
2209 | end |
---|
2210 | !*********************************************************************** |
---|
2211 | subroutine feti_resloc(n,ni,nj,a,lcmat,cmat,y,x,z) |
---|
2212 | implicit REAL(kind=8) (a-h,o-z) |
---|
2213 | dimension a(n,5) |
---|
2214 | ! routine de resolution utilisant la decomposition frontale par blocs |
---|
2215 | ! de la matrice heptadiagonale a . |
---|
2216 | ! cmat contient les termes des blocs diagonaux factorises . |
---|
2217 | dimension cmat(lcmat) |
---|
2218 | ! d contient l'inverse de la diagonale de Crout . |
---|
2219 | dimension y(n),x(n),z(n) |
---|
2220 | ! |
---|
2221 | call feti_vmov(n,y,z) |
---|
2222 | ndimd=ni |
---|
2223 | lbd=ndimd*ndimd |
---|
2224 | ! descente du systeme . |
---|
2225 | do nb=1,nj-1 |
---|
2226 | ! calcul du second membre condense : z(i+1)=y(i+1)- FD-1 z(i) . |
---|
2227 | iacmat=(nb-1)*lbd+1 |
---|
2228 | n0=(nb-1)*ndimd |
---|
2229 | call feti_mxv(ndimd,ndimd,cmat(iacmat),z(n0+1),x(n0+1)) |
---|
2230 | call feti_pbloi(n,ni,nj,a,x,z,nb) |
---|
2231 | end do |
---|
2232 | ! inversion du dernier bloc diagonal . |
---|
2233 | iacmat=(nj-1)*lbd+1 |
---|
2234 | n0=(nj-1)*ndimd |
---|
2235 | call feti_mxv(ndimd,ndimd,cmat(iacmat),z(n0+1),x(n0+1)) |
---|
2236 | ! remontee du systeme par bloc . |
---|
2237 | do nb=nj-1,1,-1 |
---|
2238 | ! calcul du second membre condense : z(i)=z(i)- E x(i+1) . |
---|
2239 | call feti_pblos(n,ni,nj,a,x,z,nb) |
---|
2240 | ! calcul de la solution pour le bloc diagonal courant . |
---|
2241 | iacmat=(nb-1)*lbd+1 |
---|
2242 | n0=(nb-1)*ndimd |
---|
2243 | call feti_mxv(ndimd,ndimd,cmat(iacmat),z(n0+1),x(n0+1)) |
---|
2244 | end do |
---|
2245 | |
---|
2246 | end |
---|
2247 | !*********************************************************************** |
---|
2248 | subroutine feti_pblos(n,ni,nj,a,x,y,nb) |
---|
2249 | ! calcul du produit matrice-vecteur morse : y = y - A x . |
---|
2250 | ! produit par le bloc superieur numero nb . |
---|
2251 | implicit REAL(kind=8) (a-h,o-z) |
---|
2252 | dimension a(n,5) |
---|
2253 | dimension x(n),y(n) |
---|
2254 | |
---|
2255 | n0=(nb-1)*ni |
---|
2256 | ! produit par la deuxieme sur-diagonale a(i,i+ni) . |
---|
2257 | do i=n0+1,n0+ni |
---|
2258 | y(i)=y(i)-x(i+ni)*a(i,5) |
---|
2259 | end do |
---|
2260 | |
---|
2261 | end subroutine feti_pblos |
---|
2262 | !*********************************************************************** |
---|
2263 | subroutine feti_pbloi(n,ni,nj,a,x,y,nb) |
---|
2264 | ! calcul du produit matrice-vecteur morse : y = y - A x . |
---|
2265 | ! produit par le bloc inferieur numero nb . |
---|
2266 | implicit REAL(kind=8) (a-h,o-z) |
---|
2267 | dimension a(n,5) |
---|
2268 | dimension x(n),y(n) |
---|
2269 | |
---|
2270 | n0=nb*ni |
---|
2271 | ! produit par la deuxieme sous-diagonale a(i,i-ni) . |
---|
2272 | do i=n0+1,n0+ni |
---|
2273 | y(i)=y(i)-x(i-ni)*a(i,1) |
---|
2274 | end do |
---|
2275 | |
---|
2276 | end subroutine feti_pbloi |
---|
2277 | !*********************************************************************** |
---|
2278 | subroutine feti_proax(n,ni,nj,a,x,y) |
---|
2279 | ! calcul du produit matrice-vecteur morse par la matrice heptadiagonale |
---|
2280 | ! a stockee par bandes . |
---|
2281 | implicit REAL(kind=8) (a-h,o-z) |
---|
2282 | dimension a(n,5) |
---|
2283 | dimension x(n),y(n) |
---|
2284 | |
---|
2285 | ! produit par la diagonale . |
---|
2286 | do i=1,n |
---|
2287 | y(i)=x(i)*a(i,3) |
---|
2288 | end do |
---|
2289 | ! produit par la premiere sur-diagonale a(i,i+1) . |
---|
2290 | do i=1,n-1 |
---|
2291 | y(i)=y(i)+x(i+1)*a(i,4) |
---|
2292 | end do |
---|
2293 | ! produit par la deuxieme sur-diagonale a(i,i+ni) . |
---|
2294 | do i=1,n-ni |
---|
2295 | y(i)=y(i)+x(i+ni)*a(i,5) |
---|
2296 | end do |
---|
2297 | ! produit par la premiere sous-diagonale a(i,i-1) |
---|
2298 | do i=2,n |
---|
2299 | y(i)=y(i)+x(i-1)*a(i,2) |
---|
2300 | end do |
---|
2301 | ! produit par la deuxieme sous-diagonale a(i,i-ni) |
---|
2302 | do i=1+ni,n |
---|
2303 | y(i)=y(i)+x(i-ni)*a(i,1) |
---|
2304 | end do |
---|
2305 | |
---|
2306 | end subroutine feti_proax |
---|
2307 | |
---|
2308 | #else |
---|
2309 | !! no use of FETI librairy |
---|
2310 | #endif |
---|
2311 | |
---|
2312 | !!====================================================================== |
---|
2313 | END MODULE lib_feti |
---|