New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
modlinktomodel.F90 in branches/UKMO/r8727_WAVE-2_Clementi_add_coupling/NEMOGCM/EXTERNAL/AGRIF/AGRIF_FILES – NEMO

source: branches/UKMO/r8727_WAVE-2_Clementi_add_coupling/NEMOGCM/EXTERNAL/AGRIF/AGRIF_FILES/modlinktomodel.F90 @ 8749

Last change on this file since 8749 was 8749, checked in by jcastill, 6 years ago

Remove svn keywords

File size: 9.3 KB
Line 
1!     AGRIF (Adaptive Grid Refinement In Fortran)
2!
3!     Copyright (C) 2003 Laurent Debreu (Laurent.Debreu@imag.fr)
4!                        Christophe Vouland (Christophe.Vouland@imag.fr)
5!
6!     This program is free software; you can redistribute it and/or modify
7!     it under the terms of the GNU General Public License as published by
8!     the Free Software Foundation; either version 2 of the License, or
9!     (at your option) any later version.
10!
11!     This program is distributed in the hope that it will be useful,
12!     but WITHOUT ANY WARRANTY; without even the implied warranty of
13!     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
14!     GNU General Public License for more details.
15!
16!     You should have received a copy of the GNU General Public License
17!     along with this program; if not, write to the Free Software
18!     Foundation, Inc., 59 Temple Place -  Suite 330, Boston, MA 02111-1307, USA.
19!
20!
21!> Module Agrif_Link
22!>
23!> This module is used to link AGRIF files to the model.
24!
25module Agrif_Link
26!
27    interface
28!
29        subroutine Agrif_clustering_def ( )
30!
31        end subroutine Agrif_clustering_def
32!
33        subroutine Agrif_Set_numberofcells ( Agrif_Gr )
34            use Agrif_Grids, only : Agrif_Grid
35            type(Agrif_Grid), pointer :: Agrif_Gr   !< Pointer on the current grid
36        end subroutine Agrif_Set_numberofcells
37!
38        subroutine Agrif_Get_numberofcells ( Agrif_Gr )
39            use Agrif_Grids, only : Agrif_Grid
40            type(Agrif_Grid), pointer :: Agrif_Gr   !< Pointer on the current grid
41        end subroutine Agrif_Get_numberofcells
42!
43    end interface
44!
45    abstract interface
46!
47        subroutine alloc_proc ( Agrif_Gr )
48            use Agrif_Grids, only : Agrif_Grid
49            type(Agrif_Grid), pointer :: Agrif_Gr   !< Pointer on the current grid
50        end subroutine alloc_proc
51!
52        subroutine typedef_proc ( )
53            implicit none
54        end subroutine typedef_proc
55!
56    end interface
57   
58    procedure(alloc_proc)   :: Agrif_Allocationcalls
59    procedure(typedef_proc) :: Agrif_probdim_modtype_def
60!
61end module Agrif_Link
62!
63!===================================================================================================
64!  function Agrif_parent
65!        modify by conv. To use : un_parent = Agrif_Parent(un)
66!===================================================================================================
67!  function Agrif_Get_Coarse_Grid
68!        modify by conv. To use : un_Mygrid = Agrif_Get_Coarse_grid(un)
69!===================================================================================================
70!  function Agrif_Rhox
71!        modify by conv. To use : var = Agrif_Rhox()
72!                    REAL(Agrif_Curgrid % spaceref(1))
73!===================================================================================================
74!  function Agrif_Parent_Rhox
75!        modify by conv. To use : var = Agrif_Parent_Rhox()
76!                    REAL(Agrif_Curgrid % parent % spaceref(1))
77!===================================================================================================
78!  function Agrif_Irhox
79!        modify by conv. To use : var = Agrif_Parent_IRhox()
80!                    Agrif_Curgrid % spaceref(1)
81!===================================================================================================
82!  function Agrif_Parent_Irhox
83!        modify by conv. To use : var = Agrif_Parent_IRhox()
84!                    Agrif_Curgrid % parent % spaceref(1)
85!===================================================================================================
86!  function Agrif_Rhoy
87!        modify by conv. To use : var = Agrif_Rhoy()
88!                    REAL(Agrif_Curgrid % spaceref(2))
89!===================================================================================================
90!  function Agrif_Parent_Rhoy
91!        modify by conv. To use : var = Agrif_Parent_Rhoy()
92!                    REAL(Agrif_Curgrid % parent % spaceref(2))
93!===================================================================================================
94!  function Agrif_Irhoy
95!        modify by conv. To use : var = Agrif_Parent_IRhoy()
96!                    Agrif_Curgrid % spaceref(2)
97!===================================================================================================
98!  function Agrif_Parent_Irhoy
99!        modify by conv. To use : var = Agrif_Parent_IRhoy()
100!                    Agrif_Curgrid % parent % spaceref(2)
101!===================================================================================================
102!  function Agrif_Rhoz
103!        modify by conv. To use : var = Agrif_Rhoz()
104!                    REAL(Agrif_Curgrid % spaceref(3))
105!===================================================================================================
106!  function Agrif_Parent_Rhoz
107!        modify by conv. To use : var = Agrif_Parent_Rhoz()
108!                    REAL(Agrif_Curgrid % parent % spaceref(3))
109!===================================================================================================
110!  function Agrif_Irhoz
111!        modify by conv. To use : var = Agrif_Parent_IRhoz()
112!                    Agrif_Curgrid % spaceref(3)
113!===================================================================================================
114!  function Agrif_Parent_Irhoz
115!        modify by conv. To use : var = Agrif_Parent_IRhoz()
116!                    Agrif_Curgrid % parent % spaceref(3)
117!===================================================================================================
118!  function Agrif_NearCommonBorderX
119!        modify by conv. To use : var = Agrif_NearCommonBorderX()
120!                       Agrif_Curgrid % NearRootBorder(1)
121!===================================================================================================
122!  function Agrif_NearCommonBorderY
123!        modify by conv. To use : var = Agrif_NearCommonBorderY()
124!                       Agrif_Curgrid % NearRootBorder(2)
125!===================================================================================================
126!  function Agrif_NearCommonBorderZ
127!        modify by conv. To use : var = Agrif_NearCommonBorderZ()
128!                       Agrif_Curgrid % NearRootBorder(3)
129!===================================================================================================
130!  function Agrif_DistantCommonBorderX
131!        modify by conv. To use : var = Agrif_DistantCommonBorderX()
132!                       Agrif_Curgrid % DistantRootBorder(1)
133!===================================================================================================
134!  function Agrif_DistantCommonBorderY
135!        modify by conv. To use : var = Agrif_DistantCommonBorderY()
136!                       Agrif_Curgrid % DistantRootBorder(2)
137!===================================================================================================
138!  function Agrif_DistantCommonBorderZ
139!        modify by conv. To use : var = Agrif_DistantCommonBorderZ()
140!                       Agrif_Curgrid % DistantRootBorder(3)
141!===================================================================================================
142!  function Agrif_Nb_Step
143!        modify by conv. To use : var = Agrif_Nb_Step()
144!                          Agrif_Curgrid % ngridstep
145!===================================================================================================
146!  function Agrif_Nb_Fine_Grids
147!        modify by conv. To use : var = Agrif_Nb_Fine_Grids()
148!                         Agrif_nbfixedgrids
149!===================================================================================================
150!  function Agrif_Ix
151!        modify by conv. To use : var = Agrif_Ix()
152!                         Agrif_Curgrid % ix(1)
153!===================================================================================================
154!  function Agrif_Parent_Ix
155!        modify by conv. To use : var = Agrif_Parent_Ix()
156!                        Agrif_Curgrid % parent % ix(1)
157!===================================================================================================
158!  function Agrif_Iy
159!        modify by conv. To use : var = Agrif_Iy()
160!                        Agrif_Curgrid % ix(2)
161!===================================================================================================
162!  function Agrif_Parent_Iy
163!        modify by conv. To use : var = Agrif_Parent_Iy()
164!                       Agrif_Curgrid % parent % ix(2)
165!===================================================================================================
166!  function Agrif_Iz
167!        modify by conv. To use : var = Agrif_Iz()
168!                      Agrif_Curgrid % ix(3)
169!===================================================================================================
170!  function Agrif_Parent_Iz
171!        modify by conv. To use : var = Agrif_Parent_Iz()
172!                     Agrif_Curgrid % parent % ix(3)
173!===================================================================================================
174!  function Agrif_Get_grid_id
175!        modify by conv. To use : var = Agrif_Get_grid_id()
176!                    Agrif_Curgrid % grid_id
177!===================================================================================================
178!  function Agrif_Get_parent_id
179!        modify by conv. To use : var = Agrif_Get_parent_id()
180!                    Agrif_Curgrid % parent % grid_id
181!===================================================================================================
Note: See TracBrowser for help on using the repository browser.