Back to home page

MITgcm

 
 

    


File indexing completed on 2023-09-21 05:10:51 UTC

view on githubraw file Latest commit 96b00645 on 2023-09-20 15:15:14 UTC
5fc0cd3689 Dani*0001 #include "STREAMICE_OPTIONS.h"
                0002 
                0003 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0004 
                0005 CBOP
6bd8ccb6c4 Dani*0006       SUBROUTINE STREAMICE_PETSC_NUMERATE (myThid)
5fc0cd3689 Dani*0007 
8a34959769 dngo*0008 C     *============================================================*
5fc0cd3689 Dani*0009 C     | SUBROUTINE                                                 |
                0010 C     | o                                                          |
8a34959769 dngo*0011 C     *============================================================*
5fc0cd3689 Dani*0012       IMPLICIT NONE
                0013 
                0014 C     === Global variables ===
                0015 #include "SIZE.h"
                0016 #include "GRID.h"
                0017 #include "EEPARAMS.h"
                0018 #include "PARAMS.h"
                0019 #include "STREAMICE.h"
                0020 #ifdef ALLOW_PETSC
acb2b21670 Dani*0021 #ifdef ALLOW_USE_MPI
                0022 #include "EESUPPORT.h"
                0023 #endif
5fc0cd3689 Dani*0024 #endif
                0025 
6bd8ccb6c4 Dani*0026       INTEGER myThid
5fc0cd3689 Dani*0027 
                0028 #ifdef ALLOW_STREAMICE
                0029 
96b006450c dngo*0030 #ifdef ALLOW_PETSC
                0031       INTEGER i, j, bi, bj
                0032 c     CHARACTER*(MAX_LEN_MBUF) msgBuf
5fc0cd3689 Dani*0033 #ifdef ALLOW_USE_MPI
8a34959769 dngo*0034       integer mpiRC
5fc0cd3689 Dani*0035 #endif
8a34959769 dngo*0036       integer mpimywid
5fc0cd3689 Dani*0037       _RS DoFCount
                0038       integer n_dofs_proc_loc (0:nPx*nPy-1)
                0039       integer n_dofs_cum_sum (0:nPx*nPy-1)
                0040 
                0041       DO bj = myByLo(myThid), myByHi(myThid)
                0042        DO bi = myBxLo(myThid), myBxHi(myThid)
                0043         DO j=1,sNy
                0044          DO i=1,sNx
                0045            streamice_petsc_dofs_u (i,j,bi,bj) = -2.0
                0046            streamice_petsc_dofs_v (i,j,bi,bj) = -2.0
                0047          ENDDO
                0048         ENDDO
                0049        ENDDO
                0050       ENDDO
                0051 
                0052       DoFCount = -1.0
                0053       DO bj = myByLo(myThid), myByHi(myThid)
                0054        DO bi = myBxLo(myThid), myBxHi(myThid)
                0055         DO j=1,sNy
                0056          DO i=1,sNx
                0057 
                0058 C   DOFS ARE NUMBERED AS FOLLOWS ON PROCESSOR DOMAIN:
                0059 C    grid is stepped through in order bj, bi, j, i
                0060 C    1) if umask(i,j,bi,bj)==1, the counter is updated by 1;
                0061 C        streamice_petsc_dofs_u is assigned the counter;
                0062 C        o/w  streamice_petsc_dofs_u is assigned -1
                0063 C    2) if vmask(i,j,bi,bj)==1, the counter is updated by 1;
                0064 C        streamice_petsc_dofs_v is assigned the counter;
                0065 C        o/w  streamice_petsc_dofs_v is assigned -1
                0066 C    NOTE THESE NUMBERING ARRAYS ARE USED TO CONSTRUCT PETSC VECTORS AND MATRIX
                0067 
                0068           if (STREAMICE_umask (i,j,bi,bj).eq.1.0) THEN
                0069            DoFCount = DoFCount + 1.0
                0070            streamice_petsc_dofs_u (i,j,bi,bj) = DoFCount
                0071           else
                0072            streamice_petsc_dofs_u (i,j,bi,bj) = -1.0
                0073           endif
                0074 
                0075           if (STREAMICE_vmask (i,j,bi,bj).eq.1.0) THEN
                0076            DoFCount = DoFCount + 1.0
                0077            streamice_petsc_dofs_v (i,j,bi,bj) = DoFCount
                0078           else
                0079            streamice_petsc_dofs_v (i,j,bi,bj) = -1.0
                0080           endif
                0081 
                0082          ENDDO
                0083         ENDDO
                0084        ENDDO
                0085       ENDDO
                0086 
                0087       print *, "DOF_COUNT", dofcount
                0088 
                0089 #ifdef ALLOW_USE_MPI
                0090 
                0091       DO i=0,nPx*nPy-1
                0092        n_dofs_proc_loc (i) = 0
                0093       ENDDO
                0094 
                0095       CALL MPI_COMM_RANK( MPI_COMM_WORLD, mpiMyWId, mpiRC )
                0096 
                0097       n_dofs_proc_loc (mpiMyWId) = INT(DoFCount)+1
                0098 
                0099       CALL MPI_Allreduce(n_dofs_proc_loc,n_dofs_process,nPx*nPy,
                0100      &       MPI_INTEGER, MPI_SUM,MPI_COMM_MODEL,mpiRC)
                0101 
                0102       n_dofs_cum_sum(0) = 0
                0103 
                0104       DO i=1,nPx*nPy-1
                0105        n_dofs_cum_sum(i) = n_dofs_cum_sum(i-1)+
                0106      &                     n_dofs_process(i-1)
                0107       ENDDO
                0108 
                0109 #else /* ALLOW_USE_MPI */
                0110 
                0111       n_dofs_process (0) = INT(DoFCount)+1
8a34959769 dngo*0112       n_dofs_cum_sum (0) = 0
                0113       mpimywid = 0
5fc0cd3689 Dani*0114 
                0115 #endif /* ALLOW_USE_MPI */
                0116 
                0117       DO bj = myByLo(myThid), myByHi(myThid)
                0118        DO bi = myBxLo(myThid), myBxHi(myThid)
                0119         DO j=1,sNy
                0120          DO i=1,sNx
                0121           IF (streamice_petsc_dofs_u(i,j,bi,bj).ge.0 ) THEN
                0122            streamice_petsc_dofs_u(i,j,bi,bj) =
                0123      &      streamice_petsc_dofs_u(i,j,bi,bj) +
                0124      &      n_dofs_cum_sum(mpimywid)
                0125           ENDIF
                0126           IF (streamice_petsc_dofs_v(i,j,bi,bj).ge.0 ) THEN
                0127            streamice_petsc_dofs_v(i,j,bi,bj) =
                0128      &      streamice_petsc_dofs_v(i,j,bi,bj) +
                0129      &      n_dofs_cum_sum(mpimywid)
                0130           ENDIF
                0131          ENDDO
                0132         ENDDO
                0133        ENDDO
                0134       ENDDO
                0135 
                0136       _EXCH_XY_RS(streamice_petsc_dofs_u,myThid)
                0137       _EXCH_XY_RS(streamice_petsc_dofs_v,myThid)
                0138 
                0139 #endif /* ALLOW_PETSC */
                0140 
                0141 #endif
                0142       RETURN
                0143       END