[3] | 1 | MODULE solmat |
---|
| 2 | !!====================================================================== |
---|
| 3 | !! *** MODULE solmat *** |
---|
| 4 | !! solver : construction of the matrix |
---|
| 5 | !!====================================================================== |
---|
| 6 | |
---|
| 7 | !!---------------------------------------------------------------------- |
---|
| 8 | !! sol_mat : Construction of the matrix of used by the elliptic solvers |
---|
| 9 | !! fetsch : |
---|
| 10 | !! fetmat : |
---|
| 11 | !! fetstr : |
---|
| 12 | !!---------------------------------------------------------------------- |
---|
| 13 | !! * Modules used |
---|
| 14 | USE oce ! ocean dynamics and active tracers |
---|
| 15 | USE dom_oce ! ocean space and time domain |
---|
| 16 | USE sol_oce ! ocean solver |
---|
| 17 | USE phycst ! physical constants |
---|
| 18 | USE obc_oce ! ocean open boundary conditions |
---|
| 19 | USE lib_mpp ! distributed memory computing |
---|
| 20 | |
---|
| 21 | IMPLICIT NONE |
---|
| 22 | PRIVATE |
---|
| 23 | |
---|
| 24 | !! * Routine accessibility |
---|
| 25 | PUBLIC sol_mat ! routine called by inisol.F90 |
---|
| 26 | !!---------------------------------------------------------------------- |
---|
| 27 | !! OPA 9.0 , LODYC-IPSL (2003) |
---|
| 28 | !!---------------------------------------------------------------------- |
---|
| 29 | |
---|
| 30 | CONTAINS |
---|
| 31 | |
---|
| 32 | SUBROUTINE sol_mat |
---|
| 33 | !!---------------------------------------------------------------------- |
---|
| 34 | !! *** ROUTINE sol_mat *** |
---|
| 35 | !! |
---|
| 36 | !! ** Purpose : Construction of the matrix of used by the elliptic |
---|
| 37 | !! solvers (either sor, pcg or feti methods). |
---|
| 38 | !! |
---|
| 39 | !! ** Method : The matrix depends on the type of free surface: |
---|
| 40 | !! * default option: rigid lid and bsf |
---|
| 41 | !! The matrix is built for the barotropic stream function system. |
---|
| 42 | !! a diagonal preconditioning matrix is also defined. |
---|
| 43 | !! * 'key_dynspg_fsc' defined: free surface |
---|
| 44 | !! The matrix is built for the divergence of the transport system |
---|
| 45 | !! a diagonal preconditioning matrix is also defined. |
---|
| 46 | !! Note that for feti solver (nsolv=3) a specific initialization |
---|
| 47 | !! is required (call to fetstr.F) for memory allocation and inter- |
---|
| 48 | !! face definition. |
---|
| 49 | !! |
---|
| 50 | !! ** Action : - gcp : extra-diagonal elements of the matrix |
---|
| 51 | !! - gcdmat : preconditioning matrix (diagonal elements) |
---|
| 52 | !! - gcdprc : inverse of the preconditioning matrix |
---|
| 53 | !! |
---|
| 54 | !! History : |
---|
| 55 | !! 1.0 ! 88-04 (G. Madec) Original code |
---|
| 56 | !! ! 91-11 (G. Madec) |
---|
| 57 | !! ! 93-03 (M. Guyon) symetrical conditions |
---|
| 58 | !! ! 93-06 (M. Guyon) suppress pointers |
---|
| 59 | !! ! 96-05 (G. Madec) merge sor and pcg formulations |
---|
| 60 | !! ! 96-11 (A. Weaver) correction to preconditioning |
---|
| 61 | !! ! 98-02 (M. Guyon) FETI method |
---|
| 62 | !! 8.5 ! 02-08 (G. Madec) F90: Free form |
---|
| 63 | !! ! 02-11 (C. Talandier, A-M. Treguier) Free surface & Open boundaries |
---|
| 64 | !!---------------------------------------------------------------------- |
---|
| 65 | !! * Local declarations |
---|
| 66 | INTEGER :: ji, jj ! dummy loop indices |
---|
| 67 | INTEGER :: ii, ij, iiend, ijend ! temporary integers |
---|
| 68 | REAL(wp) :: zcoefs, zcoefw, zcoefe, zcoefn ! temporary scalars |
---|
| 69 | REAL(wp) :: z2dt |
---|
| 70 | #if defined key_dynspg_fsc |
---|
| 71 | REAL(wp) :: zcoef |
---|
| 72 | #endif |
---|
| 73 | !!---------------------------------------------------------------------- |
---|
| 74 | |
---|
| 75 | ! FETI method ( nsolv = 3) |
---|
| 76 | ! memory allocation and interface definition for the solver |
---|
| 77 | |
---|
| 78 | IF( nsolv == 3 ) CALL fetstr |
---|
| 79 | |
---|
| 80 | |
---|
| 81 | ! 1. Construction of the matrix |
---|
| 82 | ! ----------------------------- |
---|
| 83 | |
---|
| 84 | ! initialize to zero |
---|
| 85 | gcp(:,:,1) = 0.e0 |
---|
| 86 | gcp(:,:,2) = 0.e0 |
---|
| 87 | gcp(:,:,3) = 0.e0 |
---|
| 88 | gcp(:,:,4) = 0.e0 |
---|
| 89 | |
---|
| 90 | gcdprc(:,:) = 0.e0 |
---|
| 91 | gcdmat(:,:) = 0.e0 |
---|
| 92 | |
---|
| 93 | z2dt = 2. * rdt |
---|
| 94 | |
---|
| 95 | #if defined key_dynspg_fsc && ! defined key_obc |
---|
| 96 | |
---|
| 97 | ! defined the coefficients for free surface elliptic system |
---|
| 98 | |
---|
| 99 | DO jj = 2, jpjm1 |
---|
| 100 | DO ji = 2, jpim1 |
---|
| 101 | zcoef = z2dt * z2dt * g * rnu * bmask(ji,jj) |
---|
| 102 | zcoefs = -zcoef * hv(ji ,jj-1) * e1v(ji ,jj-1) / e2v(ji ,jj-1) ! south coefficient |
---|
| 103 | zcoefw = -zcoef * hu(ji-1,jj ) * e2u(ji-1,jj ) / e1u(ji-1,jj ) ! west coefficient |
---|
| 104 | zcoefe = -zcoef * hu(ji ,jj ) * e2u(ji ,jj ) / e1u(ji ,jj ) ! east coefficient |
---|
| 105 | zcoefn = -zcoef * hv(ji ,jj ) * e1v(ji ,jj ) / e2v(ji ,jj ) ! north coefficient |
---|
| 106 | gcp(ji,jj,1) = zcoefs |
---|
| 107 | gcp(ji,jj,2) = zcoefw |
---|
| 108 | gcp(ji,jj,3) = zcoefe |
---|
| 109 | gcp(ji,jj,4) = zcoefn |
---|
| 110 | gcdmat(ji,jj) = e1t(ji,jj) * e2t(ji,jj) * bmask(ji,jj) & ! diagonal coefficient |
---|
| 111 | - zcoefs -zcoefw -zcoefe -zcoefn |
---|
| 112 | END DO |
---|
| 113 | END DO |
---|
| 114 | |
---|
| 115 | # elif defined key_dynspg_fsc && defined key_obc |
---|
| 116 | |
---|
| 117 | ! defined gcdmat in the case of open boundaries |
---|
| 118 | |
---|
| 119 | DO jj = 2, jpjm1 |
---|
| 120 | DO ji = 2, jpim1 |
---|
| 121 | zcoef = z2dt * z2dt * g * rnu * bmask(ji,jj) |
---|
| 122 | ! south coefficient |
---|
| 123 | IF( ( lpsouthobc ) .AND. ( jj == njs0p1 ) ) THEN |
---|
| 124 | zcoefs = -zcoef * hv(ji,jj-1) * e1v(ji,jj-1)/e2v(ji,jj-1)*(1.-vsmsk(ji,1)) |
---|
| 125 | ELSE |
---|
| 126 | zcoefs = -zcoef * hv(ji,jj-1) * e1v(ji,jj-1)/e2v(ji,jj-1) |
---|
| 127 | END IF |
---|
| 128 | gcp(ji,jj,1) = zcoefs |
---|
| 129 | |
---|
| 130 | ! west coefficient |
---|
| 131 | IF( ( lpwestobc ) .AND. ( ji == niw0p1 ) ) THEN |
---|
| 132 | zcoefw = -zcoef * hu(ji-1,jj) * e2u(ji-1,jj)/e1u(ji-1,jj)*(1.-uwmsk(jj,1)) |
---|
| 133 | ELSE |
---|
| 134 | zcoefw = -zcoef * hu(ji-1,jj) * e2u(ji-1,jj)/e1u(ji-1,jj) |
---|
| 135 | END IF |
---|
| 136 | gcp(ji,jj,2) = zcoefw |
---|
| 137 | |
---|
| 138 | ! east coefficient |
---|
| 139 | IF( ( lpeastobc ) .AND. ( ji == nie0 ) ) THEN |
---|
| 140 | zcoefe = -zcoef * hu(ji,jj) * e2u(ji,jj)/e1u(ji,jj)*(1.-uemsk(jj,1)) |
---|
| 141 | ELSE |
---|
| 142 | zcoefe = -zcoef * hu(ji,jj) * e2u(ji,jj)/e1u(ji,jj) |
---|
| 143 | END IF |
---|
| 144 | gcp(ji,jj,3) = zcoefe |
---|
| 145 | |
---|
| 146 | ! north coefficient |
---|
| 147 | IF( ( lpnorthobc ) .AND. ( jj == njn0 ) ) THEN |
---|
| 148 | zcoefn = -zcoef * hv(ji,jj) * e1v(ji,jj)/e2v(ji,jj)*(1.-vnmsk(ji,1)) |
---|
| 149 | ELSE |
---|
| 150 | zcoefn = -zcoef * hv(ji,jj) * e1v(ji,jj)/e2v(ji,jj) |
---|
| 151 | END IF |
---|
| 152 | gcp(ji,jj,4) = zcoefn |
---|
| 153 | |
---|
| 154 | ! diagonal coefficient |
---|
| 155 | gcdmat(ji,jj) = e1t(ji,jj)*e2t(ji,jj)*bmask(ji,jj) & |
---|
| 156 | - zcoefs -zcoefw -zcoefe -zcoefn |
---|
| 157 | END DO |
---|
| 158 | END DO |
---|
| 159 | |
---|
| 160 | # else |
---|
| 161 | |
---|
| 162 | ! defined the coefficients for bsf elliptic system |
---|
| 163 | |
---|
| 164 | DO jj = 2, jpjm1 |
---|
| 165 | DO ji = 2, jpim1 |
---|
| 166 | zcoefs = -hur(ji ,jj ) * e1u(ji ,jj ) / e2u(ji ,jj ) * bmask(ji,jj) ! south coefficient |
---|
| 167 | zcoefw = -hvr(ji ,jj ) * e2v(ji ,jj ) / e1v(ji ,jj ) * bmask(ji,jj) ! west coefficient |
---|
| 168 | zcoefe = -hvr(ji+1,jj ) * e2v(ji+1,jj ) / e1v(ji+1,jj ) * bmask(ji,jj) ! east coefficient |
---|
| 169 | zcoefn = -hur(ji ,jj+1) * e1u(ji ,jj+1) / e2u(ji ,jj+1) * bmask(ji,jj) ! north coefficient |
---|
| 170 | gcp(ji,jj,1) = zcoefs |
---|
| 171 | gcp(ji,jj,2) = zcoefw |
---|
| 172 | gcp(ji,jj,3) = zcoefe |
---|
| 173 | gcp(ji,jj,4) = zcoefn |
---|
| 174 | gcdmat(ji,jj) = -zcoefs -zcoefw -zcoefe -zcoefn ! diagonal coefficient |
---|
| 175 | |
---|
| 176 | END DO |
---|
| 177 | END DO |
---|
| 178 | |
---|
| 179 | #endif |
---|
| 180 | |
---|
| 181 | ! 2. Boundary conditions |
---|
| 182 | ! ---------------------- |
---|
| 183 | |
---|
| 184 | ! Cyclic east-west boundary conditions |
---|
| 185 | ! ji=2 is the column east of ji=jpim1 and reciprocally, |
---|
| 186 | ! ji=jpim1 is the column west of ji=2 |
---|
| 187 | ! all the coef are already set to zero as bmask is initialized to |
---|
| 188 | ! zero for ji=1 and ji=jpj in dommsk. |
---|
| 189 | |
---|
| 190 | ! Symetrical conditions |
---|
| 191 | ! free surface: no specific action |
---|
| 192 | ! bsf system: n-s gradient of bsf = 0 along j=2 (perhaps a bug !!!!!!) |
---|
| 193 | ! the diagonal coefficient of the southern grid points must be modify to |
---|
| 194 | ! account for the existence of the south symmetric bassin. |
---|
| 195 | |
---|
| 196 | #if ! defined key_dynspg_fsc |
---|
| 197 | IF( nperio == 2 ) THEN |
---|
| 198 | DO ji = 1, jpi |
---|
| 199 | IF( bmask(ji,2) /= 0.e0 ) THEN |
---|
| 200 | zcoefs = - hur(ji,2)*e1u(ji,2)/e2u(ji,2) |
---|
| 201 | gcdmat(ji,2) = gcdmat(ji,2) - zcoefs |
---|
| 202 | ENDIF |
---|
| 203 | END DO |
---|
| 204 | ENDIF |
---|
| 205 | #endif |
---|
| 206 | |
---|
| 207 | ! North fold boundary condition |
---|
| 208 | ! all the coef are already set to zero as bmask is initialized to |
---|
| 209 | ! zero on duplicated lignes and portion of lignes |
---|
| 210 | |
---|
| 211 | ! 3. Preconditioned matrix |
---|
| 212 | ! ------------------------ |
---|
| 213 | |
---|
| 214 | IF( nsolv /= 3 ) THEN |
---|
| 215 | |
---|
| 216 | ! SOR and PCG solvers |
---|
| 217 | DO jj = 1, jpj |
---|
| 218 | DO ji = 1, jpi |
---|
| 219 | IF( bmask(ji,jj) /= 0. ) gcdprc(ji,jj) = 1.e0 / gcdmat(ji,jj) |
---|
| 220 | END DO |
---|
| 221 | END DO |
---|
| 222 | |
---|
| 223 | gcp(:,:,1) = gcp(:,:,1) * gcdprc(:,:) |
---|
| 224 | gcp(:,:,2) = gcp(:,:,2) * gcdprc(:,:) |
---|
| 225 | gcp(:,:,3) = gcp(:,:,3) * gcdprc(:,:) |
---|
| 226 | gcp(:,:,4) = gcp(:,:,4) * gcdprc(:,:) |
---|
| 227 | |
---|
| 228 | ELSE |
---|
| 229 | |
---|
| 230 | ! FETI method |
---|
| 231 | ! if feti solver : gcdprc is a mask for the non-overlapping |
---|
| 232 | ! data structuring |
---|
| 233 | |
---|
| 234 | DO jj = 1, jpj |
---|
| 235 | DO ji = 1, jpi |
---|
| 236 | IF( bmask(ji,jj) /= 0. ) THEN |
---|
| 237 | gcdprc(ji,jj) = 1. |
---|
| 238 | ELSE |
---|
| 239 | gcdprc(ji,jj) = 0. |
---|
| 240 | ENDIF |
---|
| 241 | END DO |
---|
| 242 | END DO |
---|
| 243 | |
---|
| 244 | ! so "common" line & "common" column have to be !=0 except on global |
---|
| 245 | ! domain boundaries |
---|
| 246 | ! pbs with nbondi if nperio != 2 ? |
---|
| 247 | ! ii = nldi-1 |
---|
| 248 | ! pb with nldi value if jperio==1 : nbondi modifyed at the end |
---|
| 249 | ! of inimpp.F => pb |
---|
| 250 | ! pb with periodicity conditions : iiend, ijend |
---|
| 251 | |
---|
| 252 | ijend = nlej |
---|
| 253 | iiend = nlei |
---|
| 254 | IF( jperio == 1 .OR. jperio == 4 .OR. jperio == 6 ) iiend = nlci - jpreci |
---|
| 255 | ii = jpreci |
---|
| 256 | |
---|
| 257 | ! case number 1 |
---|
| 258 | |
---|
| 259 | IF( nbondi /= -1 .AND. nbondi /= 2 ) THEN |
---|
| 260 | DO jj = 1, ijend |
---|
| 261 | IF( fmask(ii,jj,1) == 1. ) gcdprc(ii,jj) = 1. |
---|
| 262 | END DO |
---|
| 263 | ENDIF |
---|
| 264 | |
---|
| 265 | ! case number 2 |
---|
| 266 | |
---|
| 267 | IF( nperio == 1 .OR. nperio == 4 .OR. nperio == 6 ) THEN |
---|
| 268 | DO jj = 1, ijend |
---|
| 269 | IF( fmask(ii,jj,1) == 1. ) gcdprc(ii,jj) = 1. |
---|
| 270 | END DO |
---|
| 271 | ENDIF |
---|
| 272 | |
---|
| 273 | ! ij = nldj-1 |
---|
| 274 | ! pb with nldi value if jperio==1 : nbondi modifyed at the end |
---|
| 275 | ! of inimpp.F => pb, here homogeneisation... |
---|
| 276 | |
---|
| 277 | ij = jprecj |
---|
| 278 | IF( nbondj /= -1 .AND. nbondj /= 2 ) THEN |
---|
| 279 | DO ji = 1, iiend |
---|
| 280 | IF( fmask(ji,ij,1) == 1. ) gcdprc(ji,ij) = 1. |
---|
| 281 | END DO |
---|
| 282 | ENDIF |
---|
| 283 | ENDIF |
---|
| 284 | |
---|
| 285 | |
---|
| 286 | ! 4. Initialization the arrays used in pcg |
---|
| 287 | ! ---------------------------------------- |
---|
| 288 | gcx (:,:) = 0.e0 |
---|
| 289 | gcxb (:,:) = 0.e0 |
---|
| 290 | gcb (:,:) = 0.e0 |
---|
| 291 | gcr (:,:) = 0.e0 |
---|
| 292 | gcdes(:,:) = 0.e0 |
---|
| 293 | gccd (:,:) = 0.e0 |
---|
| 294 | |
---|
| 295 | ! FETI method |
---|
| 296 | IF( nsolv == 3 ) THEN |
---|
| 297 | CALL fetmat ! Matrix treatment : Neumann condition, inverse computation |
---|
| 298 | CALL fetsch ! data framework for the Schur Dual solver |
---|
| 299 | ENDIF |
---|
| 300 | |
---|
| 301 | END SUBROUTINE sol_mat |
---|
| 302 | |
---|
| 303 | #if defined key_feti |
---|
| 304 | |
---|
| 305 | SUBROUTINE fetstr |
---|
| 306 | !!--------------------------------------------------------------------- |
---|
| 307 | !! *** ROUTINE fetstr *** |
---|
| 308 | !! |
---|
| 309 | !! ** Purpose : Construction of the matrix of the barotropic stream |
---|
| 310 | !! function system. |
---|
| 311 | !! Finite Elements Tearing & Interconnecting (FETI) approach |
---|
| 312 | !! Memory allocation and interface definition for the solver |
---|
| 313 | !! |
---|
| 314 | !! ** Method : |
---|
| 315 | !! |
---|
| 316 | !! References : |
---|
| 317 | !! Guyon, M, Roux, F-X, Chartier, M and Fraunie, P, 1994 : |
---|
| 318 | !! A domain decomposition solver to compute the barotropic |
---|
| 319 | !! component of an OGCM in the parallel processing field. |
---|
| 320 | !! Ocean Modelling, issue 105, december 94. |
---|
| 321 | !! |
---|
| 322 | !! History : |
---|
| 323 | !! ! 98-02 (M. Guyon) Original code |
---|
| 324 | !! 8.5 ! 02-09 (G. Madec) F90: Free form and module |
---|
| 325 | !!---------------------------------------------------------------------- |
---|
| 326 | !! * Modules used |
---|
| 327 | USE lib_feti ! feti librairy |
---|
| 328 | !! * Local declarations |
---|
| 329 | INTEGER :: iiend, ijend, iperio ! temporary integers |
---|
| 330 | !!--------------------------------------------------------------------- |
---|
| 331 | |
---|
| 332 | |
---|
| 333 | ! Preconditioning technics of the Dual Schur Operator |
---|
| 334 | ! <= definition of the Coarse Grid solver |
---|
| 335 | ! <= dimension of the nullspace of the local operators |
---|
| 336 | ! <= Neumann boundaries conditions |
---|
| 337 | |
---|
| 338 | ! 0. Initializations |
---|
| 339 | ! ------------------ |
---|
| 340 | |
---|
| 341 | ndkerep = 1 |
---|
| 342 | |
---|
| 343 | ! initialization of the superstructures management |
---|
| 344 | |
---|
| 345 | malxm = 1 |
---|
| 346 | malim = 1 |
---|
| 347 | |
---|
| 348 | ! memory space for the pcpg associated with the FETI dual formulation |
---|
| 349 | ! ndkerep is associated to the list of rigid modes, |
---|
| 350 | ! ndkerep == 1 because the Dual Operator |
---|
| 351 | ! is a first order operator due to SPG elliptic Operator is a |
---|
| 352 | ! second order operator |
---|
| 353 | |
---|
| 354 | nim = 50 |
---|
| 355 | nim = nim + ndkerep |
---|
| 356 | nim = nim + 2*jpi + 2*jpj |
---|
| 357 | nim = nim + jpi*jpj |
---|
| 358 | |
---|
| 359 | nxm = 33 |
---|
| 360 | nxm = nxm + 4*jpnij |
---|
| 361 | nxm = nxm + 19*(jpi+jpj) |
---|
| 362 | nxm = nxm + 13*jpi*jpj |
---|
| 363 | nxm = nxm + jpi*jpi*jpj |
---|
| 364 | |
---|
| 365 | ! krylov space memory |
---|
| 366 | |
---|
| 367 | iperio = 0 |
---|
| 368 | IF( jperio == 1 .OR. jperio == 4 .OR. jperio == 6) iperio = 1 |
---|
| 369 | nxm = nxm + 3*(jpnij-jpni)*jpi |
---|
| 370 | nxm = nxm + 3*(jpnij-jpnj+iperio)*jpj |
---|
| 371 | nxm = nxm + 2*(jpi+jpj)*(jpnij-jpni)*jpi |
---|
| 372 | nxm = nxm + 2*(jpi+jpj)*(jpnij-jpnj+iperio)*jpj |
---|
| 373 | |
---|
| 374 | ! Resolution with the Schur dual Method ( frontal and local solver by |
---|
| 375 | ! blocks |
---|
| 376 | ! Case with a local symetrical matrix |
---|
| 377 | ! The local matrix is stored in a multi-column form |
---|
| 378 | ! The total number of nodes for this subdomain is named "noeuds" |
---|
| 379 | |
---|
| 380 | noeuds = jpi*jpj |
---|
| 381 | nifmat = jpi-1 |
---|
| 382 | njfmat = jpj-1 |
---|
| 383 | nelem = nifmat*njfmat |
---|
| 384 | npe = 4 |
---|
| 385 | nmorse = 5*noeuds |
---|
| 386 | |
---|
| 387 | ! 1. mesh building |
---|
| 388 | ! ---------------- |
---|
| 389 | |
---|
| 390 | ! definition of specific information for a subdomain |
---|
| 391 | ! narea : subdomain number = processor number +1 |
---|
| 392 | ! ninterf : neighbour subdomain number |
---|
| 393 | ! nni : interface point number |
---|
| 394 | ! ndvois array : neighbour subdomain list |
---|
| 395 | ! maplistin array : node pointer at each interface |
---|
| 396 | ! maplistin array : concatened list of interface nodes |
---|
| 397 | |
---|
| 398 | ! messag coding is necessary by interface type for avoid collision |
---|
| 399 | ! if nperio == 1 |
---|
| 400 | |
---|
| 401 | ! lint array : indoor interface list / type |
---|
| 402 | ! lext array : outdoor interface list / type |
---|
| 403 | |
---|
| 404 | ! domain with jpniXjpnj subdomains |
---|
| 405 | |
---|
| 406 | CALL feti_inisub(nifmat,njfmat,nbondi,nbondj,nperio, & |
---|
| 407 | nbsw,nbnw,nbse,nbne,ninterf,ninterfc,nni,nnic) |
---|
| 408 | |
---|
| 409 | CALL feti_creadr(malim,malimax,nim,3*ninterf ,mandvois ,'ndvois' ) |
---|
| 410 | CALL feti_creadr(malim,malimax,nim,3*ninterfc,mandvoisc,'ndvoisc') |
---|
| 411 | CALL feti_creadr(malim,malimax,nim,ninterfc+1,maplistin,'plistin') |
---|
| 412 | CALL feti_creadr(malim,malimax,nim,nnic ,malistin ,'listin' ) |
---|
| 413 | |
---|
| 414 | ! pb with periodicity conditions : iiend, ijend |
---|
| 415 | |
---|
| 416 | ijend = nlej |
---|
| 417 | iiend = nlei |
---|
| 418 | IF (jperio == 1) iiend = nlci - jpreci |
---|
| 419 | |
---|
| 420 | CALL feti_subound(nifmat,njfmat,nldi,iiend,nldj,ijend, & |
---|
| 421 | narea,nbondi,nbondj,nperio, & |
---|
| 422 | ninterf,ninterfc, & |
---|
| 423 | nowe,noea,noso,nono, & |
---|
| 424 | nbsw,nbnw,nbse,nbne, & |
---|
| 425 | npsw,npnw,npse,npne, & |
---|
| 426 | mfet(mandvois),mfet(mandvoisc), & |
---|
| 427 | mfet(maplistin),nnic,mfet(malistin) ) |
---|
| 428 | |
---|
| 429 | END SUBROUTINE fetstr |
---|
| 430 | |
---|
| 431 | |
---|
| 432 | SUBROUTINE fetmat |
---|
| 433 | !!--------------------------------------------------------------------- |
---|
| 434 | !! *** ROUTINE fetmat *** |
---|
| 435 | !! |
---|
| 436 | !! ** Purpose : Construction of the matrix of the barotropic stream |
---|
| 437 | !! function system. |
---|
| 438 | !! Finite Elements Tearing & Interconnecting (FETI) approach |
---|
| 439 | !! Matrix treatment : Neumann condition, inverse computation |
---|
| 440 | !! |
---|
| 441 | !! ** Method : |
---|
| 442 | !! |
---|
| 443 | !! References : |
---|
| 444 | !! Guyon, M, Roux, F-X, Chartier, M and Fraunie, P, 1994 : |
---|
| 445 | !! A domain decomposition solver to compute the barotropic |
---|
| 446 | !! component of an OGCM in the parallel processing field. |
---|
| 447 | !! Ocean Modelling, issue 105, december 94. |
---|
| 448 | !! |
---|
| 449 | !! History : |
---|
| 450 | !! ! 98-02 (M. Guyon) Original code |
---|
| 451 | !! 8.5 ! 02-09 (G. Madec) F90: Free form and module |
---|
| 452 | !!---------------------------------------------------------------------- |
---|
| 453 | !! * Modules used |
---|
| 454 | USE lib_feti ! feti librairy |
---|
| 455 | !! * Local declarations |
---|
| 456 | INTEGER :: ji, jj, jk, jl |
---|
| 457 | INTEGER :: iimask(jpi,jpj) |
---|
| 458 | INTEGER :: iiend, ijend |
---|
| 459 | REAL(wp) :: zres, zres2, zdemi |
---|
| 460 | !!--------------------------------------------------------------------- |
---|
| 461 | |
---|
| 462 | ! Matrix computation |
---|
| 463 | ! ------------------ |
---|
| 464 | |
---|
| 465 | CALL feti_creadr(malxm,malxmax,nxm,nmorse,maan,'matrice a') |
---|
| 466 | |
---|
| 467 | nnitot = nni |
---|
| 468 | |
---|
| 469 | CALL mpp_sum(nnitot,1,numit0ete) |
---|
| 470 | CALL feti_creadr(malxm,malxmax,nxm,npe*npe,maae,'ae') |
---|
| 471 | |
---|
| 472 | ! initialisation of the local barotropic matrix |
---|
| 473 | ! local boundary conditions on the halo |
---|
| 474 | |
---|
| 475 | CALL lbc_lnk( gcp(:,:,1), 'F', 1) |
---|
| 476 | CALL lbc_lnk( gcp(:,:,2), 'F', 1) |
---|
| 477 | CALL lbc_lnk( gcp(:,:,3), 'F', 1) |
---|
| 478 | CALL lbc_lnk( gcp(:,:,4), 'F', 1) |
---|
| 479 | CALL lbc_lnk( gcdmat , 'T', 1) |
---|
| 480 | |
---|
| 481 | ! Neumann conditions |
---|
| 482 | ! initialisation of the integer Neumann Mask |
---|
| 483 | |
---|
| 484 | CALL feti_iclr(jpi*jpj,iimask) |
---|
| 485 | DO jj = 1, jpj |
---|
| 486 | DO ji = 1, jpi |
---|
| 487 | iimask(ji,jj) = INT( gcdprc(ji,jj) ) |
---|
| 488 | END DO |
---|
| 489 | END DO |
---|
| 490 | |
---|
| 491 | ! regularization of the local matrix |
---|
| 492 | |
---|
| 493 | DO jj = 1, jpj |
---|
| 494 | DO ji = 1, jpi |
---|
| 495 | gcdmat(ji,jj) = gcdmat(ji,jj) * gcdprc(ji,jj) + 1. - gcdprc(ji,jj) |
---|
| 496 | END DO |
---|
| 497 | END DO |
---|
| 498 | |
---|
| 499 | DO jk = 1, 4 |
---|
| 500 | DO jj = 1, jpj |
---|
| 501 | DO ji = 1, jpi |
---|
| 502 | gcp(ji,jj,jk) = gcp(ji,jj,jk) * gcdprc(ji,jj) |
---|
| 503 | END DO |
---|
| 504 | END DO |
---|
| 505 | END DO |
---|
| 506 | |
---|
| 507 | ! implementation of the west, east, north & south Neumann conditions |
---|
| 508 | |
---|
| 509 | zdemi = 0.5 |
---|
| 510 | |
---|
| 511 | ! pb with periodicity conditions : iiend, ijend |
---|
| 512 | |
---|
| 513 | ijend = nlej |
---|
| 514 | iiend = nlei |
---|
| 515 | IF( jperio == 1 .OR. jperio == 4 .OR. jperio == 6 ) iiend = nlci - jpreci |
---|
| 516 | |
---|
| 517 | IF( nbondi == 2 .AND. (nperio /= 1 .OR. nperio /= 4 .OR. nperio == 6) ) THEN |
---|
| 518 | |
---|
| 519 | ! with the periodicity : no east/west interface if nbondi = 2 |
---|
| 520 | ! and nperio != 1 |
---|
| 521 | |
---|
| 522 | ELSE |
---|
| 523 | ! west |
---|
| 524 | IF( nbondi /= -1 ) THEN |
---|
| 525 | DO jj = 1, jpj |
---|
| 526 | IF( iimask(1,jj) /= 0 ) THEN |
---|
| 527 | gcp(1,jj,2) = 0. |
---|
| 528 | gcp(1,jj,1) = zdemi * gcp(1,jj,1) |
---|
| 529 | gcp(1,jj,4) = zdemi * gcp(1,jj,4) |
---|
| 530 | ENDIF |
---|
| 531 | END DO |
---|
| 532 | DO jj = 1, jpj |
---|
| 533 | IF( iimask(1,jj) /= 0 ) THEN |
---|
| 534 | gcdmat(1,jj) = - ( gcp(1,jj,1) + gcp(1,jj,2) + gcp(1,jj,3) + gcp(1,jj,4) ) |
---|
| 535 | ENDIF |
---|
| 536 | END DO |
---|
| 537 | ENDIF |
---|
| 538 | ! east |
---|
| 539 | IF( nbondi /= 1 ) THEN |
---|
| 540 | DO jj = 1, jpj |
---|
| 541 | IF( iimask(iiend,jj) /= 0 ) THEN |
---|
| 542 | gcp(iiend,jj,3) = 0. |
---|
| 543 | gcp(iiend,jj,1) = zdemi * gcp(iiend,jj,1) |
---|
| 544 | gcp(iiend,jj,4) = zdemi * gcp(iiend,jj,4) |
---|
| 545 | ENDIF |
---|
| 546 | END DO |
---|
| 547 | DO jj = 1, jpj |
---|
| 548 | IF( iimask(iiend,jj) /= 0 ) THEN |
---|
| 549 | gcdmat(iiend,jj) = - ( gcp(iiend,jj,1) + gcp(iiend,jj,2) & |
---|
| 550 | + gcp(iiend,jj,3) + gcp(iiend,jj,4) ) |
---|
| 551 | ENDIF |
---|
| 552 | END DO |
---|
| 553 | ENDIF |
---|
| 554 | ENDIF |
---|
| 555 | |
---|
| 556 | ! south |
---|
| 557 | IF( nbondj /= -1 .AND. nbondj /= 2 ) THEN |
---|
| 558 | DO ji = 1, jpi |
---|
| 559 | IF( iimask(ji,1) /= 0 ) THEN |
---|
| 560 | gcp(ji,1,1) = 0. |
---|
| 561 | gcp(ji,1,2) = zdemi * gcp(ji,1,2) |
---|
| 562 | gcp(ji,1,3) = zdemi * gcp(ji,1,3) |
---|
| 563 | ENDIF |
---|
| 564 | END DO |
---|
| 565 | DO ji = 1, jpi |
---|
| 566 | IF( iimask(ji,1) /= 0 ) THEN |
---|
| 567 | gcdmat(ji,1) = - ( gcp(ji,1,1) + gcp(ji,1,2) + gcp(ji,1,3) + gcp(ji,1,4) ) |
---|
| 568 | ENDIF |
---|
| 569 | END DO |
---|
| 570 | ENDIF |
---|
| 571 | |
---|
| 572 | ! north |
---|
| 573 | IF( nbondj /= 1 .AND. nbondj /= 2 ) THEN |
---|
| 574 | DO ji = 1, jpi |
---|
| 575 | IF( iimask(ji,ijend) /= 0 ) THEN |
---|
| 576 | gcp(ji,ijend,4) = 0. |
---|
| 577 | gcp(ji,ijend,2) = zdemi * gcp(ji,ijend,2) |
---|
| 578 | gcp(ji,ijend,3) = zdemi * gcp(ji,ijend,3) |
---|
| 579 | ENDIF |
---|
| 580 | END DO |
---|
| 581 | DO ji = 1, jpi |
---|
| 582 | IF( iimask(ji,ijend) /= 0 ) THEN |
---|
| 583 | gcdmat(ji,ijend) = - ( gcp(ji,ijend,1) + gcp(ji,ijend,2) & |
---|
| 584 | + gcp(ji,ijend,3) + gcp(ji,ijend,4) ) |
---|
| 585 | ENDIF |
---|
| 586 | END DO |
---|
| 587 | ENDIF |
---|
| 588 | |
---|
| 589 | ! matrix terms are saved in FETI solver arrays |
---|
| 590 | CALL feti_vmov(noeuds,gcp(1,1,1),wfeti(maan)) |
---|
| 591 | CALL feti_vmov(noeuds,gcp(1,1,2),wfeti(maan+noeuds)) |
---|
| 592 | CALL feti_vmov(noeuds,gcdmat,wfeti(maan+2*noeuds)) |
---|
| 593 | CALL feti_vmov(noeuds,gcp(1,1,3),wfeti(maan+3*noeuds)) |
---|
| 594 | CALL feti_vmov(noeuds,gcp(1,1,4),wfeti(maan+4*noeuds)) |
---|
| 595 | |
---|
| 596 | ! construction of Dirichlet liberty degrees array |
---|
| 597 | CALL feti_subdir(nifmat,njfmat,noeuds,ndir,iimask) |
---|
| 598 | CALL feti_creadr(malim,malimax,nim,ndir,malisdir,'lisdir') |
---|
| 599 | CALL feti_listdir(jpi,jpj,iimask,ndir,mfet(malisdir)) |
---|
| 600 | |
---|
| 601 | ! stop onto matrix term for Dirichlet conditions |
---|
| 602 | CALL feti_blomat(nifmat+1,njfmat+1,wfeti(maan),ndir,mfet(malisdir)) |
---|
| 603 | |
---|
| 604 | ! reservation of factorized diagonal blocs and temporary array for |
---|
| 605 | ! factorization |
---|
| 606 | npblo = (njfmat+1) * (nifmat+1) * (nifmat+1) |
---|
| 607 | ndimax = nifmat+1 |
---|
| 608 | |
---|
| 609 | CALL feti_creadr(malxm,malxmax,nxm,npblo,mablo,'blo') |
---|
| 610 | CALL feti_creadr(malxm,malxmax,nxm,noeuds,madia,'dia') |
---|
| 611 | CALL feti_creadr(malxm,malxmax,nxm,noeuds,mav,'v') |
---|
| 612 | CALL feti_creadr(malxm,malxmax,nxm,ndimax*ndimax,mautil,'util') |
---|
| 613 | |
---|
| 614 | ! stoping the rigid modes |
---|
| 615 | |
---|
| 616 | ! the number of rigid modes =< Max [dim(Ker(Ep))] |
---|
| 617 | ! p=1,Np |
---|
| 618 | |
---|
| 619 | CALL feti_creadr(malim,malimax,nim,ndkerep,malisblo,'lisblo') |
---|
| 620 | |
---|
| 621 | ! Matrix factorization |
---|
| 622 | |
---|
| 623 | CALL feti_front(noeuds,nifmat+1,njfmat+1,wfeti(maan),npblo, & |
---|
| 624 | wfeti(mablo),wfeti(madia), & |
---|
| 625 | wfeti(mautil),wfeti(mav),ndlblo,mfet(malisblo),ndkerep) |
---|
| 626 | CALL feti_prext(noeuds,wfeti(madia)) |
---|
| 627 | |
---|
| 628 | ! virtual dealloc => we have to see for a light f90 version |
---|
| 629 | ! the super structure is removed to clean the coarse grid |
---|
| 630 | ! solver structure |
---|
| 631 | |
---|
| 632 | malxm = madia |
---|
| 633 | CALL feti_vclr(noeuds,wfeti(madia)) |
---|
| 634 | CALL feti_vclr(noeuds,wfeti(mav)) |
---|
| 635 | CALL feti_vclr(ndimax*ndimax,wfeti(mautil)) |
---|
| 636 | |
---|
| 637 | ! ndlblo is the dimension of the local nullspace .=<. the size of the |
---|
| 638 | ! memory of the superstructure associated to the nullspace : ndkerep |
---|
| 639 | ! ndkerep is introduced to avoid messages "out of bounds" when memory |
---|
| 640 | ! is checked |
---|
| 641 | |
---|
| 642 | ! copy matrix for Dirichlet condition |
---|
| 643 | |
---|
| 644 | CALL feti_creadr(malxm,malxmax,nxm,noeuds,miax,'x') |
---|
| 645 | CALL feti_creadr(malxm,malxmax,nxm,noeuds,may,'y') |
---|
| 646 | CALL feti_creadr(malxm,malxmax,nxm,noeuds,maz,'z') |
---|
| 647 | |
---|
| 648 | ! stoping the rigid modes |
---|
| 649 | |
---|
| 650 | ! ndlblo is the dimension of the local nullspace .=<. the size of the |
---|
| 651 | ! memory of the superstructure associated to the nullspace : ndkerep |
---|
| 652 | ! ndkerep is introduced to avoid messages "out of bounds" when memory |
---|
| 653 | ! is checked |
---|
| 654 | |
---|
| 655 | CALL feti_creadr(malxm,malxmax,nxm,ndkerep*noeuds,mansp,'nsp') |
---|
| 656 | CALL feti_blomat1(nifmat+1,njfmat+1,wfeti(maan),ndlblo, & |
---|
| 657 | mfet(malisblo),wfeti(mansp)) |
---|
| 658 | |
---|
| 659 | ! computation of operator kernel |
---|
| 660 | |
---|
| 661 | CALL feti_nullsp(noeuds,nifmat+1,njfmat+1,npblo,wfeti(mablo), & |
---|
| 662 | wfeti(maan),ndlblo,mfet(malisblo),wfeti(mansp), & |
---|
| 663 | wfeti(maz)) |
---|
| 664 | |
---|
| 665 | ! test of the factorisation onto each sub domain |
---|
| 666 | |
---|
| 667 | CALL feti_init(noeuds,wfeti(may)) |
---|
| 668 | CALL feti_blodir(noeuds,wfeti(may),ndir,mfet(malisdir)) |
---|
| 669 | CALL feti_blodir(noeuds,wfeti(may),ndlblo,mfet(malisblo)) |
---|
| 670 | CALL feti_vclr(noeuds,wfeti(miax)) |
---|
| 671 | CALL feti_resloc(noeuds,nifmat+1,njfmat+1,wfeti(maan),npblo, & |
---|
| 672 | wfeti(mablo),wfeti(may),wfeti(miax),wfeti(maz)) |
---|
| 673 | CALL feti_proax(noeuds,nifmat+1,njfmat+1,wfeti(maan),wfeti(miax), & |
---|
| 674 | wfeti(maz)) |
---|
| 675 | CALL feti_blodir(noeuds,wfeti(maz),ndlblo,mfet(malisblo)) |
---|
| 676 | CALL feti_vsub(noeuds,wfeti(may),wfeti(maz),wfeti(maz)) |
---|
| 677 | |
---|
| 678 | zres2 = 0. |
---|
| 679 | DO jl = 1, noeuds |
---|
| 680 | zres2 = zres2 + wfeti(may+jl-1) * wfeti(may+jl-1) |
---|
| 681 | END DO |
---|
| 682 | CALL mpp_sum(zres2,1,zres) |
---|
| 683 | |
---|
| 684 | res2 = 0. |
---|
| 685 | DO jl = 1, noeuds |
---|
| 686 | res2 = res2 + wfeti(maz+jl-1) * wfeti(maz+jl-1) |
---|
| 687 | END DO |
---|
| 688 | res2 = res2 / zres2 |
---|
| 689 | CALL mpp_sum(res2,1,zres) |
---|
| 690 | |
---|
| 691 | res2 = SQRT(res2) |
---|
| 692 | IF(lwp) WRITE(numout,*) 'global residu : sqrt((Ax-b,Ax-b)/(b.b)) =', res2 |
---|
| 693 | |
---|
| 694 | IF( res2 > (eps/100.) ) THEN |
---|
| 695 | IF(lwp) WRITE (numout,*) 'eps is :',eps |
---|
| 696 | IF(lwp) WRITE (numout,*) 'factorized matrix precision :',res2 |
---|
| 697 | STOP |
---|
| 698 | ENDIF |
---|
| 699 | |
---|
| 700 | END SUBROUTINE fetmat |
---|
| 701 | |
---|
| 702 | |
---|
| 703 | SUBROUTINE fetsch |
---|
| 704 | !!--------------------------------------------------------------------- |
---|
| 705 | !! *** ROUTINE fetsch *** |
---|
| 706 | !! |
---|
| 707 | !! ** Purpose : |
---|
| 708 | !! Construction of the matrix of the barotropic stream function |
---|
| 709 | !! system. |
---|
| 710 | !! Finite Elements Tearing & Interconnecting (FETI) approach |
---|
| 711 | !! Data framework for the Schur Dual solve |
---|
| 712 | !! |
---|
| 713 | !! ** Method : |
---|
| 714 | !! |
---|
| 715 | !! References : |
---|
| 716 | !! Guyon, M, Roux, F-X, Chartier, M and Fraunie, P, 1994 : |
---|
| 717 | !! A domain decomposition solver to compute the barotropic |
---|
| 718 | !! component of an OGCM in the parallel processing field. |
---|
| 719 | !! Ocean Modelling, issue 105, december 94. |
---|
| 720 | !! |
---|
| 721 | !! History : |
---|
| 722 | !! ! 98-02 (M. Guyon) Original code |
---|
| 723 | !! 8.5 ! 02-09 (G. Madec) F90: Free form and module |
---|
| 724 | !!---------------------------------------------------------------------- |
---|
| 725 | !! * Modules used |
---|
| 726 | USE lib_feti ! feti librairy |
---|
| 727 | !! * Local declarations |
---|
| 728 | !!--------------------------------------------------------------------- |
---|
| 729 | |
---|
| 730 | ! computing weights for the conform construction |
---|
| 731 | |
---|
| 732 | CALL feti_creadr(malxm,malxmax,nxm,noeuds,mapoids ,'poids' ) |
---|
| 733 | CALL feti_creadr(malxm,malxmax,nxm,nnic ,mabufin ,'bufin' ) |
---|
| 734 | CALL feti_creadr(malxm,malxmax,nxm,nnic ,mabufout,'bufout') |
---|
| 735 | |
---|
| 736 | !! CALL feti_poids(ninterfc,mfet(mandvoisc),mfet(maplistin),nnic, & |
---|
| 737 | !! mfet(malistin),narea,noeuds,wfeti(mapoids),wfeti(mabufin), & |
---|
| 738 | !! wfeti(mabufout) ) |
---|
| 739 | CALL feti_poids(ninterfc, nnic, & |
---|
| 740 | mfet(malistin), noeuds,wfeti(mapoids) ) |
---|
| 741 | |
---|
| 742 | |
---|
| 743 | ! Schur dual arrays |
---|
| 744 | |
---|
| 745 | CALL feti_creadr(malxm,malxmax,nxm,noeuds,mabitw,'bitw') |
---|
| 746 | CALL feti_creadr(malxm,malxmax,nxm,noeuds,mautilu,'utilu') |
---|
| 747 | CALL feti_creadr(malxm,malxmax,nxm,nni,malambda,'lambda') |
---|
| 748 | CALL feti_creadr(malxm,malxmax,nxm,nni,mag,'g') |
---|
| 749 | CALL feti_creadr(malxm,malxmax,nxm,nni,mapg,'pg') |
---|
| 750 | CALL feti_creadr(malxm,malxmax,nxm,nni,mamg,'mg') |
---|
| 751 | CALL feti_creadr(malxm,malxmax,nxm,nni,maw,'w') |
---|
| 752 | CALL feti_creadr(malxm,malxmax,nxm,nni,madw,'dw') |
---|
| 753 | |
---|
| 754 | ! coarse grid solver dimension and arrays |
---|
| 755 | |
---|
| 756 | nitmaxete = ndlblo |
---|
| 757 | CALL mpp_sum(nitmaxete,1,numit0ete) |
---|
| 758 | |
---|
| 759 | nitmaxete = nitmaxete + 1 |
---|
| 760 | CALL feti_creadr(malxm,malxmax,nxm,ndkerep,maxnul,'xnul') |
---|
| 761 | CALL feti_creadr(malxm,malxmax,nxm,ndkerep,maynul,'ynul') |
---|
| 762 | CALL feti_creadr(malxm,malxmax,nxm,ndkerep,maeteg,'eteg') |
---|
| 763 | CALL feti_creadr(malxm,malxmax,nxm,ndkerep,maeteag,'eteag') |
---|
| 764 | CALL feti_creadr(malxm,malxmax,nxm,ndkerep*nitmaxete,maeted,'eted') |
---|
| 765 | CALL feti_creadr(malxm,malxmax,nxm,ndkerep*nitmaxete,maetead,'etead') |
---|
| 766 | CALL feti_creadr(malxm,malxmax,nxm,nitmaxete,maeteadd,'eteadd') |
---|
| 767 | CALL feti_creadr(malxm,malxmax,nxm,nitmaxete,maetegamm,'etegamm') |
---|
| 768 | CALL feti_creadr(malxm,malxmax,nxm,nni,maetew,'etew') |
---|
| 769 | CALL feti_creadr(malxm,malxmax,nxm,noeuds,maetev,'etev') |
---|
| 770 | |
---|
| 771 | ! construction of semi interface arrays |
---|
| 772 | |
---|
| 773 | CALL feti_creadr(malim,malimax,nim,ninterf+1,maplistih,'plistih') |
---|
| 774 | !! CALL feti_halfint(ninterf,mfet(mandvois),mfet(maplistin),nni, & |
---|
| 775 | !! mfet(maplistih),nnih,narea) |
---|
| 776 | CALL feti_halfint(ninterf,mfet(mandvois),mfet(maplistin), & |
---|
| 777 | mfet(maplistih),nnih ) |
---|
| 778 | |
---|
| 779 | CALL feti_creadr(malxm,malxmax,nxm,nnih,magh,'gh') |
---|
| 780 | |
---|
| 781 | ! Schur Dual Method |
---|
| 782 | |
---|
| 783 | nmaxd = nnitot / 2 |
---|
| 784 | |
---|
| 785 | ! computation of the remain array for descent directions |
---|
| 786 | |
---|
| 787 | nmaxd = min0(nmaxd,(nxm-nitmaxete-malxm)/(2*nnih+3)) |
---|
| 788 | CALL mpp_min(nmaxd,1,numit0ete) |
---|
| 789 | |
---|
| 790 | nitmax = nnitot/2 |
---|
| 791 | epsilo = eps |
---|
| 792 | ntest = 0 |
---|
| 793 | |
---|
| 794 | ! Krylov space construction |
---|
| 795 | |
---|
| 796 | CALL feti_creadr(malxm,malxmax,nxm,nnih*nmaxd,mawj,'wj') |
---|
| 797 | CALL feti_creadr(malxm,malxmax,nxm,nnih*nmaxd,madwj,'dwj') |
---|
| 798 | CALL feti_creadr(malxm,malxmax,nxm,nmaxd,madwwj,'dwwj') |
---|
| 799 | CALL feti_creadr(malxm,malxmax,nxm,nmaxd,magamm,'gamm') |
---|
| 800 | CALL feti_creadr(malxm,malxmax,nxm,max0(nmaxd,nitmaxete),mawork,'work') |
---|
| 801 | mjj0 = 0 |
---|
| 802 | numit0ete = 0 |
---|
| 803 | |
---|
| 804 | END SUBROUTINE fetsch |
---|
| 805 | |
---|
| 806 | #else |
---|
| 807 | SUBROUTINE fetstr ! Empty routine |
---|
| 808 | END SUBROUTINE fetstr |
---|
| 809 | SUBROUTINE fetmat ! Empty routine |
---|
| 810 | END SUBROUTINE fetmat |
---|
| 811 | SUBROUTINE fetsch ! Empty routine |
---|
| 812 | END SUBROUTINE fetsch |
---|
| 813 | #endif |
---|
| 814 | |
---|
| 815 | !!====================================================================== |
---|
| 816 | END MODULE solmat |
---|